MAST
Multidisciplinary-design Adaptation and Sensitivity Toolkit (MAST)
structural_modal_eigenproblem_assembly.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 // MAST includes
27 #include "base/function_base.h"
28 #include "base/assembly_base.h"
30 
31 
32 
36 _incompatible_sol_assembly(nullptr) { }
37 
38 
39 
40 
43 { }
44 
45 
46 
47 
48 void
51 
52  unsigned int
53  n = (unsigned int)sol.size();
54 
56  zero = RealVectorX::Zero(n);
57 
59  _physics_elem->set_velocity (zero); // set to zero vector for a quasi-steady analysis
60  _physics_elem->set_acceleration(zero); // set to zero vector for a quasi-steady analysis
61 
62 
64 
65  // set the incompatible mode solution if required by the
66  // element
67 
69  dynamic_cast<MAST::StructuralElementBase&>(*_physics_elem);
70 
71  if (s_elem.if_incompatible_modes())
73  }
74 }
75 
76 
77 
78 
79 void
82 
83  unsigned int
84  n = (unsigned int)sol.size();
85 
87  zero = RealVectorX::Zero(n);
88 
89  _physics_elem->set_solution (sol, true);
90 }
91 
92 
93 
94 void
97  RealMatrixX& mat_B) {
98 
99  libmesh_assert(_physics_elem);
100 
102  dynamic_cast<MAST::StructuralElementBase&>(*_physics_elem);
103 
104  RealVectorX vec = RealVectorX::Zero(mat_A.rows()); // dummy vector
105  RealMatrixX mat = RealMatrixX::Zero(mat_A.rows(), mat_A.cols()); // dummy matrix
106  mat_A.setZero();
107  mat_B.setZero();
108 
109  // calculate the Jacobian components
110  e.internal_residual(true, vec, mat_A);
111  e.side_external_residual(true, vec, mat, mat_A, _discipline->side_loads());
112  e.volume_external_residual(true, vec, mat, mat_A, _discipline->volume_loads());
113 
114  // calculate the mass matrix components
115  e.inertial_residual(true, vec, mat_B, mat, mat_A);
116 }
117 
118 
119 
120 
121 void
124  bool base_sol,
125  RealMatrixX& mat_A,
126  RealMatrixX& mat_B) {
127 
128  libmesh_assert(_physics_elem);
129 
131  dynamic_cast<MAST::StructuralElementBase&>(*_physics_elem);
132 
133  RealVectorX vec = RealVectorX::Zero(mat_A.rows()); // dummy vector
134  RealMatrixX mat = RealMatrixX::Zero(mat_A.rows(), mat_A.cols()); // dummy matrix
135  mat_A.setZero();
136  mat_B.setZero();
137 
138  // calculate the Jacobian components
139  e.internal_residual_sensitivity(f, true, vec, mat_A);
140  e.side_external_residual_sensitivity(f, true, vec, mat, mat_A, _discipline->side_loads());
141  e.volume_external_residual_sensitivity(f, true, vec, mat, mat_A, _discipline->volume_loads());
142 
143  // calculate the mass matrix components
144  e.inertial_residual_sensitivity(f, true, vec, mat_B, mat, mat_A);
145 
146  // if the linearization is about a base state, then the sensitivity of
147  // the base state will influence the sensitivity of the Jacobian
148  if (base_sol)
150 }
151 
152 
153 
154 
155 void
158  bool base_sol,
159  RealMatrixX& mat_A,
160  RealMatrixX& mat_B) {
161 
162  libmesh_assert(_physics_elem);
163  libmesh_assert(f.is_topology_parameter());
164 
165  RealVectorX vec = RealVectorX::Zero(mat_A.rows()); // dummy vector
166  RealMatrixX mat = RealMatrixX::Zero(mat_A.rows(), mat_A.cols()); // dummy matrix
167  mat_A.setZero();
168  mat_B.setZero();
169 
170  std::pair<const MAST::FieldFunction<RealVectorX>*, unsigned int>
171  val = this->get_elem_boundary_velocity_data();
172 
173  if (val.first) {
174 
176  dynamic_cast<MAST::StructuralElementBase&>(*_physics_elem);
177 
179  val.second,
180  *val.first,
181  true,
182  vec,
183  mat_A);
185  val.second,
186  *val.first,
188  true,
189  vec,
190  mat_A);
191 
193  val.second,
194  *val.first,
195  true,
196  vec,
197  mat_B,
198  mat,
199  mat_A);
200 
201  // if the linearization is about a base state, then the sensitivity of
202  // the base state will influence the sensitivity of the Jacobian
203  if (base_sol)
204  libmesh_assert(false); // to be implemented
205  //e.internal_residual_jac_dot_state_sensitivity(mat_A);
206  }
207 }
208 
209 
210 
211 void
214  bool base_sol,
216  RealMatrixX& mat_A,
217  RealMatrixX& mat_B) {
218 
219  libmesh_assert(_physics_elem);
220  libmesh_assert(f.is_topology_parameter());
221 
223  &elem = dynamic_cast<const MAST::LevelSetIntersectedElem&>(_physics_elem->elem());
224 
225  RealVectorX vec = RealVectorX::Zero(mat_A.rows()); // dummy vector
226  RealMatrixX mat = RealMatrixX::Zero(mat_A.rows(), mat_A.cols()); // dummy matrix
227  mat_A.setZero();
228  mat_B.setZero();
229 
230  // sensitivity only exists at the boundary. So, we proceed with calculation
231  // only if this element has an intersection in the interior, or with a side.
232  if (elem.if_elem_has_level_set_boundary() &&
233  elem.if_subelem_has_side_on_level_set_boundary()) {
234 
236  dynamic_cast<MAST::StructuralElementBase&>(*_physics_elem);
237 
239  elem.get_subelem_side_on_level_set_boundary(),
240  vel,
241  true,
242  vec,
243  mat_A);
245  elem.get_subelem_side_on_level_set_boundary(),
246  vel,
248  true,
249  vec,
250  mat_A);
251 
253  elem.get_subelem_side_on_level_set_boundary(),
254  vel,
255  true,
256  vec,
257  mat_B,
258  mat,
259  mat_A);
260 
261  // if the linearization is about a base state, then the sensitivity of
262  // the base state will influence the sensitivity of the Jacobian
263  if (base_sol)
264  libmesh_assert(false); // to be implemented
265  //e.internal_residual_jac_dot_state_sensitivity(mat_A);
266  }
267 }
268 
269 
270 
271 void
273 set_elem_data(unsigned int dim,
274  const libMesh::Elem& ref_elem,
275  MAST::GeomElem& elem) const {
276 
277  libmesh_assert(!_physics_elem);
278 
279  if (dim == 1) {
280 
281  const MAST::ElementPropertyCard1D& p =
282  dynamic_cast<const MAST::ElementPropertyCard1D&>(_discipline->get_property_card(ref_elem));
283 
284  elem.set_local_y_vector(p.y_vector());
285  }
286 }
287 
288 
289 void
291 init(const MAST::GeomElem& elem) {
292 
293  libmesh_assert(!_physics_elem);
294  libmesh_assert(_system);
295  libmesh_assert(_assembly);
296 
298  dynamic_cast<const MAST::ElementPropertyCardBase&>
300 
301  _physics_elem =
302  MAST::build_structural_element(*_system, elem, p).release();
303 }
304 
virtual void internal_residual_boundary_velocity(const MAST::FunctionBase &p, const unsigned int s, const MAST::FieldFunction< RealVectorX > &vel_f, bool request_jacobian, RealVectorX &f, RealMatrixX &jac)=0
calculates the term on side s: .
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...
const MAST::ElementPropertyCardBase & get_property_card(const libMesh::Elem &elem) const
get property card for the specified element
virtual bool inertial_residual_sensitivity(const MAST::FunctionBase &p, bool request_jacobian, RealVectorX &f, RealMatrixX &jac_xddot, RealMatrixX &jac_xdot, RealMatrixX &jac)
sensitivity of the inertial force contribution to system residual
virtual void set_velocity(const RealVectorX &vec, bool if_sens=false)
stores vec as velocity for element level calculations, or its sensitivity if if_sens is true...
Definition: elem_base.cpp:98
RealVectorX & y_vector()
returns value of the property val.
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 inertial_residual(bool request_jacobian, RealVectorX &f, RealMatrixX &jac_xddot, RealMatrixX &jac_xdot, RealMatrixX &jac)
inertial force contribution to system residual
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
void set_elem_incompatible_sol(MAST::StructuralElementBase &elem)
virtual std::pair< const MAST::FieldFunction< RealVectorX > *, unsigned int > get_elem_boundary_velocity_data()
searches through the side load data and populates the data with the boundary id and velocity function...
virtual bool if_incompatible_modes() const =0
virtual bool is_topology_parameter() const
Definition: function_base.h:97
const MAST::SideBCMapType & side_loads() const
virtual void elem_calculations(RealMatrixX &mat_A, RealMatrixX &mat_B)
performs the element calculations over elem, and returns the element matrices for the eigenproblem ...
Matrix< Real, Dynamic, Dynamic > RealMatrixX
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.
This class inherits from MAST::GeomElem and provides an interface to initialize FE objects on sub-ele...
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
virtual void elem_topology_sensitivity_calculations(const MAST::FunctionBase &f, bool base_sol, RealMatrixX &mat_A, RealMatrixX &mat_B)
performs the element topology sensitivity calculations over elem.
StructuralModalEigenproblemAssemblyElemOperations()
constructor associates the eigen system with this assembly object
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
virtual void set_elem_solution_sensitivity(const RealVectorX &sol)
sets the element solution sensitivity before calculations
const MAST::GeomElem & elem() const
Definition: elem_base.h:109
virtual void set_solution(const RealVectorX &vec, bool if_sens=false)
stores vec as solution for element level calculations, or its sensitivity if if_sens is true...
Definition: elem_base.cpp:60
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 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.
virtual void init(const MAST::GeomElem &elem)
initializes the object for the geometric element elem.
virtual bool internal_residual_jac_dot_state_sensitivity(RealMatrixX &jac)=0
calculates d[J]/d{x} .
virtual ~StructuralModalEigenproblemAssemblyElemOperations()
destructor resets the association with the eigen system from this assembly object ...
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 inertial_residual_boundary_velocity(const MAST::FunctionBase &p, const unsigned int s, const MAST::FieldFunction< RealVectorX > &vel_f, bool request_jacobian, RealVectorX &f, RealMatrixX &jac_xddot, RealMatrixX &jac_xdot, RealMatrixX &jac)
sensitivity of the inertial force contribution to system residual
void volume_external_residual_boundary_velocity(const MAST::FunctionBase &p, const unsigned int s, const MAST::FieldFunction< RealVectorX > &vel_f, std::multimap< libMesh::subdomain_id_type, MAST::BoundaryConditionBase *> &bc, bool request_jacobian, RealVectorX &f, RealMatrixX &jac)
boundary velocity contribution of volume external force.
virtual void set_elem_solution(const RealVectorX &sol)
sets the element solution(s) before calculations
virtual void set_acceleration(const RealVectorX &vec, bool if_sens=false)
stores vec as acceleration for element level calculations, or its sensitivity if if_sens is true...
Definition: elem_base.cpp:122
MAST::SystemInitialization * _system