MAST
structural_transient_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 
21 // MAST includes
26 #include "base/assembly_base.h"
27 #include "mesh/geom_elem.h"
28 
29 
33 
34 }
35 
36 
37 
38 
41 
42 }
43 
44 
45 
46 void
48 elem_calculations(bool if_jac,
49  RealVectorX& f_m,
50  RealVectorX& f_x,
51  RealMatrixX& f_m_jac_xdot,
52  RealMatrixX& f_m_jac,
53  RealMatrixX& f_x_jac) {
54 
55  libmesh_assert(_physics_elem);
56 
58  dynamic_cast<MAST::StructuralElementBase&>(*_physics_elem);
59 
60  f_m.setZero();
61  f_x.setZero();
62  f_m_jac_xdot.setZero();
63  f_m_jac.setZero();
64  f_x_jac.setZero();
65 
66  // assembly of the flux terms
67  e.internal_residual(if_jac, f_x, f_x_jac);
68  e.side_external_residual(if_jac,
69  f_x,
70  f_m_jac_xdot,
71  f_x_jac,
73  e.volume_external_residual(if_jac,
74  f_x,
75  f_m_jac_xdot,
76  f_x_jac,
78 
79  //assembly of the capacitance term
80  e.damping_residual(if_jac, f_m, f_m_jac_xdot, f_m_jac);
81 }
82 
83 
84 
85 
86 void
88 elem_calculations(bool if_jac,
89  RealVectorX& f_m,
90  RealVectorX& f_x,
91  RealMatrixX& f_m_jac_xddot,
92  RealMatrixX& f_m_jac_xdot,
93  RealMatrixX& f_m_jac,
94  RealMatrixX& f_x_jac_xdot,
95  RealMatrixX& f_x_jac) {
96 
97  libmesh_assert(_physics_elem);
98 
100  dynamic_cast<MAST::StructuralElementBase&>(*_physics_elem);
101 
102  f_m.setZero();
103  f_x.setZero();
104  f_m_jac_xddot.setZero();
105  f_m_jac_xdot.setZero();
106  f_m_jac.setZero();
107  f_x_jac_xdot.setZero();
108  f_x_jac.setZero();
109 
110  // assembly of the flux terms
111  e.internal_residual(if_jac, f_x, f_x_jac);
112  e.damping_residual(if_jac, f_x, f_x_jac_xdot, f_x_jac);
113  e.side_external_residual(if_jac,
114  f_x,
115  f_m_jac_xdot,
116  f_x_jac,
118  e.volume_external_residual(if_jac,
119  f_x,
120  f_m_jac_xdot,
121  f_x_jac,
123 
124  //assembly of the capacitance term
125  e.inertial_residual(if_jac, f_m, f_m_jac_xddot, f_m_jac_xdot, f_m_jac);
126 }
127 
128 
129 
130 void
133 
134  libmesh_assert(_physics_elem);
135 
137  dynamic_cast<MAST::StructuralElementBase&>(*_physics_elem);
138 
139  const unsigned int
140  n = (unsigned int) f.size();
141 
143  dummy = RealMatrixX::Zero(n, n);
144 
145  f.setZero();
146 
147  // assembly of the flux terms
148  e.linearized_internal_residual(false, f, dummy);
149  //e.linearized_damping_residual(false, f, dummy, dummy);
151  f,
152  dummy,
153  dummy,
156  f,
157  dummy,
158  dummy,
160 
161  //assembly of the capacitance term
162  e.linearized_inertial_residual(false, f, dummy, dummy, dummy);
163 }
164 
165 
166 
167 
168 void
171  RealVectorX& f_m,
172  RealVectorX& f_x) {
173 
174  libmesh_assert(_physics_elem);
175 
177  dynamic_cast<MAST::StructuralElementBase&>(*_physics_elem);
178 
179  unsigned int
180  n = (unsigned int)f_m.size();
181 
183  dummy = RealMatrixX::Zero(n, n);
184 
185  f_m.setZero();
186  f_x.setZero();
187 
188  // assembly of the flux terms
189  e.internal_residual_sensitivity(f, false, f_x, dummy);
190  e.damping_residual_sensitivity(f, false, f_x, dummy, dummy);
192  f_x,
193  dummy,
194  dummy,
197  f_x,
198  dummy,
199  dummy,
201 
202  //assembly of the capacitance term
203  e.inertial_residual_sensitivity(f, false, f_m, dummy, dummy, dummy);
204 }
205 
206 
207 
208 void
211 
212  libmesh_error(); // to be implemented
213 }
214 
215 
216 void
218 set_elem_data(unsigned int dim,
219  const libMesh::Elem& ref_elem,
220  MAST::GeomElem& elem) const {
221 
222  libmesh_assert(!_physics_elem);
223 
224  if (dim == 1) {
225 
226  const MAST::ElementPropertyCard1D& p =
227  dynamic_cast<const MAST::ElementPropertyCard1D&>(_discipline->get_property_card(ref_elem));
228 
229  elem.set_local_y_vector(p.y_vector());
230  }
231 }
232 
233 
234 void
236 
237  libmesh_assert(!_physics_elem);
238  libmesh_assert(_system);
239  libmesh_assert(_assembly);
240 
242  dynamic_cast<const MAST::ElementPropertyCardBase&>(_discipline->get_property_card(elem));
243 
244  _physics_elem =
245  MAST::build_structural_element(*_system, *_assembly, elem, p).release();
246 }
247 
const MAST::VolumeBCMapType & volume_loads() const
virtual void linearized_jacobian_solution_product(RealVectorX &f)
Calculates the product of Jacobian-solution, and Jacobian-velocity over the element for a system of t...
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
virtual bool linearized_inertial_residual(bool request_jacobian, RealVectorX &f, RealMatrixX &jac_xddot, RealMatrixX &jac_xdot, RealMatrixX &jac)
inertial force contribution to system residual of linerized problem
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
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 damping_residual(bool request_jacobian, RealVectorX &f, RealMatrixX &jac_xdot, RealMatrixX &jac)
damping 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 ~StructuralTransientAssemblyElemOperations()
destructor resets the association of this assembly object with the system
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 elem_calculations(bool if_jac, RealVectorX &f_m, RealVectorX &f_x, RealMatrixX &f_m_jac_x_dot, RealMatrixX &f_m_jac, RealMatrixX &f_x_jac)
This call for first order ode should not be used for this transient assembly.
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.
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_sensitivity_calculations(const MAST::FunctionBase &f, RealVectorX &f_m, RealVectorX &f_x)
performs the element sensitivity calculations over elem, and returns the component of element residua...
Matrix< Real, Dynamic, Dynamic > RealMatrixX
Matrix< Real, Dynamic, 1 > RealVectorX
const MAST::ElementPropertyCardBase & get_property_card(const libMesh::Elem &elem) const
get property card for the specified element
bool linearized_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.
This class acts as a wrapper around libMesh::Elem for the purpose of providing a uniform interface fo...
Definition: geom_elem.h:59
virtual bool damping_residual_sensitivity(const MAST::FunctionBase &p, bool request_jacobian, RealVectorX &f, RealMatrixX &jac_xdot, RealMatrixX &jac)
sensitivity of the damping force contribution to system residual
virtual void elem_second_derivative_dot_solution_assembly(RealMatrixX &mat)
calculates over elem, and returns the matrix in vec .
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 linearized_internal_residual(bool request_jacobian, RealVectorX &f, RealMatrixX &jac)
internal force contribution to system residual of the linearized problem
StructuralTransientAssemblyElemOperations()
constructor associates this assembly object with the system
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 init(const MAST::GeomElem &elem)
initializes the object for the geometric element elem.
bool linearized_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.
MAST::SystemInitialization * _system