MAST
Multidisciplinary-design Adaptation and Sensitivity Toolkit (MAST)
structural_buckling_eigenproblem_elem_operations.cpp
Go to the documentation of this file.
1 /*
2  * MAST: Multidisciplinary-design Adaptation and Sensitivity Toolkit
3  * Copyright (C) 2013-2020 Manav Bhatia and MAST authors
4  *
5  * This library is free software; you can redistribute it and/or
6  * modify it under the terms of the GNU Lesser General Public
7  * License as published by the Free Software Foundation; either
8  * version 2.1 of the License, or (at your option) any later version.
9  *
10  * This library is distributed in the hope that it will be useful,
11  * but WITHOUT ANY WARRANTY; without even the implied warranty of
12  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
13  * Lesser General Public License for more details.
14  *
15  * You should have received a copy of the GNU Lesser General Public
16  * License along with this library; if not, write to the Free Software
17  * Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA
18  */
19 
20 
21 // MAST includes
27 #include "base/assembly_base.h"
28 #include "mesh/geom_elem.h"
29 
33 
34 }
35 
36 
37 
40 
41 }
42 
43 
44 void
46 set_elem_data(unsigned int dim,
47  const libMesh::Elem& ref_elem,
48  MAST::GeomElem& elem) const {
49 
50  libmesh_assert(!_physics_elem);
51 
52  if (dim == 1) {
53 
55  dynamic_cast<const MAST::ElementPropertyCard1D&>(_discipline->get_property_card(ref_elem));
56 
57  elem.set_local_y_vector(p.y_vector());
58  }
59 }
60 
61 
62 void
64 
65  libmesh_assert(!_physics_elem);
66  libmesh_assert(_system);
67  libmesh_assert(_assembly);
68 
70  dynamic_cast<const MAST::ElementPropertyCardBase&>
72 
74  MAST::build_structural_element(*_system, elem, p).release();
75 }
76 
77 
78 
79 void
82  RealMatrixX& mat_B) {
83 
85  dynamic_cast<MAST::StructuralElementBase&>(*_physics_elem);
86 
88  vec = RealVectorX::Zero(mat_A.rows()); // dummy vector
90  dummy = RealMatrixX::Zero(vec.size(), vec.size());
91 
92  mat_A.setZero();
93 
94  // calculate the Jacobian components
95  e.internal_residual(true, vec, mat_A);
96  e.prestress_residual(true, vec, mat_A);
98  vec,
99  dummy,
100  mat_A,
103  vec,
104  dummy,
105  mat_A,
107 }
108 
109 
110 
111 
112 void
115  bool base_sol,
116  RealMatrixX& mat_A,
117  RealMatrixX& mat_B) {
118 
120  dynamic_cast<MAST::StructuralElementBase&>(*_physics_elem);
121 
123  vec = RealVectorX::Zero(mat_A.rows()); // dummy vector
125  dummy = RealMatrixX::Zero(vec.size(), vec.size());
126 
127  mat_A.setZero();
128 
129  // calculate the Jacobian components
130  e.internal_residual_sensitivity(f, true, vec, mat_A);
131  e.prestress_residual_sensitivity(f, true, vec, mat_A);
133  vec,
134  dummy,
135  mat_A,
138  vec,
139  dummy,
140  mat_A,
142 }
143 
const MAST::ElementPropertyCardBase & get_property_card(const libMesh::Elem &elem) const
get property card for the specified element
RealVectorX & y_vector()
returns value of the property val.
virtual void init(const MAST::GeomElem &elem)
initializes the object for the geometric element elem.
void set_local_y_vector(const RealVectorX &y_vec)
for 1D elements the transformed coordinate system attached to the element defines the local x-axis al...
Definition: geom_elem.cpp:119
virtual bool internal_residual_sensitivity(const MAST::FunctionBase &p, bool request_jacobian, RealVectorX &f, RealMatrixX &jac)=0
sensitivity of the internal force contribution to system residual
MAST::PhysicsDisciplineBase * _discipline
virtual void set_elem_data(unsigned int dim, const libMesh::Elem &ref_elem, MAST::GeomElem &elem) const
sets the structural element y-vector if 1D element is used.
const MAST::SideBCMapType & side_loads() const
Matrix< Real, Dynamic, Dynamic > RealMatrixX
virtual bool prestress_residual_sensitivity(const MAST::FunctionBase &p, bool request_jacobian, RealVectorX &f, RealMatrixX &jac)=0
sensitivity of the prestress force contribution to system residual
bool side_external_residual(bool request_jacobian, RealVectorX &f, RealMatrixX &jac_xdot, RealMatrixX &jac, std::multimap< libMesh::boundary_id_type, MAST::BoundaryConditionBase *> &bc)
side external force contribution to system residual.
Matrix< Real, Dynamic, 1 > RealVectorX
bool volume_external_residual(bool request_jacobian, RealVectorX &f, RealMatrixX &jac_xdot, RealMatrixX &jac, std::multimap< libMesh::subdomain_id_type, MAST::BoundaryConditionBase *> &bc)
volume external force contribution to system residual.
std::unique_ptr< MAST::StructuralElementBase > build_structural_element(MAST::SystemInitialization &sys, const MAST::GeomElem &elem, const MAST::ElementPropertyCardBase &p)
builds the structural element for the specified element type
This class acts as a wrapper around libMesh::Elem for the purpose of providing a uniform interface fo...
Definition: geom_elem.h:59
const MAST::VolumeBCMapType & volume_loads() const
bool side_external_residual_sensitivity(const MAST::FunctionBase &p, bool request_jacobian, RealVectorX &f, RealMatrixX &jac_xdot, RealMatrixX &jac, std::multimap< libMesh::boundary_id_type, MAST::BoundaryConditionBase *> &bc)
sensitivity of the side external force contribution to system residual
virtual bool prestress_residual(bool request_jacobian, RealVectorX &f, RealMatrixX &jac)=0
prestress force contribution to system residual
bool volume_external_residual_sensitivity(const MAST::FunctionBase &p, bool request_jacobian, RealVectorX &f, RealMatrixX &jac_xdot, RealMatrixX &jac, std::multimap< libMesh::subdomain_id_type, MAST::BoundaryConditionBase *> &bc)
sensitivity of the volume external force contribution to system residual
virtual bool internal_residual(bool request_jacobian, RealVectorX &f, RealMatrixX &jac)=0
internal force contribution to system residual
virtual void elem_calculations(RealMatrixX &mat_A, RealMatrixX &mat_B)
performs the element calculations over elem, and returns the element matrices for the eigenproblem ...
virtual void elem_sensitivity_calculations(const MAST::FunctionBase &f, bool base_sol, RealMatrixX &mat_A, RealMatrixX &mat_B)
performs the element sensitivity calculations over elem, and returns the element matrices for the eig...
MAST::SystemInitialization * _system