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-2019 Manav Bhatia
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 
37 
38 
39 
42 { }
43 
44 
45 
46 
47 void
50 
51  unsigned int
52  n = (unsigned int)sol.size();
53 
55  zero = RealVectorX::Zero(n);
56 
58  _physics_elem->set_velocity (zero); // set to zero vector for a quasi-steady analysis
59  _physics_elem->set_acceleration(zero); // set to zero vector for a quasi-steady analysis
60 
61 
63 
64  // set the incompatible mode solution if required by the
65  // element
66 
68  dynamic_cast<MAST::StructuralElementBase&>(*_physics_elem);
69 
70  if (s_elem.if_incompatible_modes())
72  }
73 }
74 
75 
76 
77 
78 void
81 
82  unsigned int
83  n = (unsigned int)sol.size();
84 
86  zero = RealVectorX::Zero(n);
87 
88  _physics_elem->set_solution (sol, true);
89 }
90 
91 
92 
93 void
96  RealMatrixX& mat_B) {
97 
98  libmesh_assert(_physics_elem);
99 
101  dynamic_cast<MAST::StructuralElementBase&>(*_physics_elem);
102 
103  RealVectorX vec = RealVectorX::Zero(mat_A.rows()); // dummy vector
104  RealMatrixX mat = RealMatrixX::Zero(mat_A.rows(), mat_A.cols()); // dummy matrix
105  mat_A.setZero();
106  mat_B.setZero();
107 
108  // calculate the Jacobian components
109  e.internal_residual(true, vec, mat_A);
110  e.side_external_residual(true, vec, mat, mat_A, _discipline->side_loads());
111  e.volume_external_residual(true, vec, mat, mat_A, _discipline->volume_loads());
112 
113  // calculate the mass matrix components
114  e.inertial_residual(true, vec, mat_B, mat, mat_A);
115 }
116 
117 
118 
119 
120 void
123  bool base_sol,
124  RealMatrixX& mat_A,
125  RealMatrixX& mat_B) {
126 
127  libmesh_assert(_physics_elem);
128 
130  dynamic_cast<MAST::StructuralElementBase&>(*_physics_elem);
131 
132  RealVectorX vec = RealVectorX::Zero(mat_A.rows()); // dummy vector
133  RealMatrixX mat = RealMatrixX::Zero(mat_A.rows(), mat_A.cols()); // dummy matrix
134  mat_A.setZero();
135  mat_B.setZero();
136 
137  // calculate the Jacobian components
138  e.internal_residual_sensitivity(f, true, vec, mat_A);
139  e.side_external_residual_sensitivity(f, true, vec, mat, mat_A, _discipline->side_loads());
140  e.volume_external_residual_sensitivity(f, true, vec, mat, mat_A, _discipline->volume_loads());
141 
142  // calculate the mass matrix components
143  e.inertial_residual_sensitivity(f, true, vec, mat_B, mat, mat_A);
144 
145  // if the linearization is about a base state, then the sensitivity of
146  // the base state will influence the sensitivity of the Jacobian
147  if (base_sol)
149 }
150 
151 
152 
153 
154 void
157  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 
166  &elem = dynamic_cast<const MAST::LevelSetIntersectedElem&>(_physics_elem->elem());
167 
168  RealVectorX vec = RealVectorX::Zero(mat_A.rows()); // dummy vector
169  RealMatrixX mat = RealMatrixX::Zero(mat_A.rows(), mat_A.cols()); // dummy matrix
170  mat_A.setZero();
171  mat_B.setZero();
172 
173  // sensitivity only exists at the boundary. So, we proceed with calculation
174  // only if this element has an intersection in the interior, or with a side.
175  if (elem.if_elem_has_level_set_boundary() &&
176  elem.if_subelem_has_side_on_level_set_boundary()) {
177 
179  dynamic_cast<MAST::StructuralElementBase&>(*_physics_elem);
180 
182  elem.get_subelem_side_on_level_set_boundary(),
183  vel,
184  true,
185  vec,
186  mat_A);
188  elem.get_subelem_side_on_level_set_boundary(),
189  vel,
191  true,
192  vec,
193  mat_A);
194 
196  elem.get_subelem_side_on_level_set_boundary(),
197  vel,
198  true,
199  vec,
200  mat_B,
201  mat,
202  mat_A);
203 
204  // if the linearization is about a base state, then the sensitivity of
205  // the base state will influence the sensitivity of the Jacobian
206  if (base_sol)
207  libmesh_assert(false); // to be implemented
208  //e.internal_residual_jac_dot_state_sensitivity(mat_A);
209  }
210 }
211 
212 
213 
214 void
216 set_elem_data(unsigned int dim,
217  const libMesh::Elem& ref_elem,
218  MAST::GeomElem& elem) const {
219 
220  libmesh_assert(!_physics_elem);
221 
222  if (dim == 1) {
223 
224  const MAST::ElementPropertyCard1D& p =
225  dynamic_cast<const MAST::ElementPropertyCard1D&>(_discipline->get_property_card(ref_elem));
226 
227  elem.set_local_y_vector(p.y_vector());
228  }
229 }
230 
231 
232 void
234 init(const MAST::GeomElem& elem) {
235 
236  libmesh_assert(!_physics_elem);
237  libmesh_assert(_system);
238  libmesh_assert(_assembly);
239 
241  dynamic_cast<const MAST::ElementPropertyCardBase&>
243 
244  _physics_elem =
245  MAST::build_structural_element(*_system, *_assembly, elem, p).release();
246 }
247 
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::VolumeBCMapType & volume_loads() const
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.
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
std::unique_ptr< MAST::StructuralElementBase > build_structural_element(MAST::SystemInitialization &sys, MAST::AssemblyBase &assembly, const MAST::GeomElem &elem, const MAST::ElementPropertyCardBase &p)
builds the structural element for the specified element type
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:100
RealVectorX & y_vector()
returns value of the property val.
const MAST::GeomElem & elem() const
Definition: elem_base.h:117
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
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.
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 bool if_incompatible_modes() const =0
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 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
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.
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.
This class inherits from MAST::GeomElem and provides an interface to initialize FE objects on sub-ele...
Matrix< Real, Dynamic, 1 > RealVectorX
const MAST::ElementPropertyCardBase & get_property_card(const libMesh::Elem &elem) const
get property card for the specified element
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
virtual void set_elem_solution_sensitivity(const RealVectorX &sol)
sets the element solution sensitivity before calculations
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:62
virtual void elem_topology_sensitivity_calculations(const MAST::FunctionBase &f, bool base_sol, const MAST::FieldFunction< RealVectorX > &vel, RealMatrixX &mat_A, RealMatrixX &mat_B)
performs the element topology sensitivity calculations over elem.
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} .
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 ~StructuralModalEigenproblemAssemblyElemOperations()
destructor resets the association with the eigen system from this assembly object ...
const MAST::SideBCMapType & side_loads() const
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
virtual bool is_topology_parameter() const
Definition: function_base.h:97
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:124
MAST::SystemInitialization * _system