MAST
Multidisciplinary-design Adaptation and Sensitivity Toolkit (MAST)
compliance_output.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
24 #include "base/assembly_base.h"
26 #include "base/nonlinear_system.h"
32 #include "mesh/geom_elem.h"
33 
34 // libMesh includes
35 #include "libmesh/parallel.h"
36 
39 _compliance (0.),
40 _dcompliance_dp (0.) {
41 
42 }
43 
44 
45 
46 
48 
49 }
50 
51 
52 
53 
54 void
56 
57  _compliance = 0.;
58  _dcompliance_dp = 0.;
59 }
60 
61 
62 
63 void
65 
66  _compliance = 0.;
67  _dcompliance_dp = 0.;
68 }
69 
70 
71 void
73 
74  // make sure that this has not been initialized ana calculated for all elems
75  libmesh_assert(_physics_elem);
76 
78 
80  dynamic_cast<MAST::StructuralElementBase&>(*_physics_elem);
81 
83  vec = RealVectorX::Zero(e.sol().size());
84 
86  dummy = RealMatrixX::Zero(vec.size(), vec.size());
87 
88  e.side_external_residual(false,
89  vec,
90  dummy,
91  dummy,
94  vec,
95  dummy,
96  dummy,
98 
99  // compute the contribution of this element to compliance
100  _compliance -= vec.dot(e.sol());
101  }
102 }
103 
104 
105 
106 void
108 
109  // make sure that this has not been initialized ana calculated for all elems
110  libmesh_assert(_physics_elem);
111 
113 
115  dynamic_cast<MAST::StructuralElementBase&>(*_physics_elem);
116 
118  vec = RealVectorX::Zero(e.sol().size());
119 
121  dummy = RealMatrixX::Zero(vec.size(), vec.size());
122 
124  vec,
125  dummy,
126  dummy,
129  vec,
130  dummy,
131  dummy,
133 
134  // ask for the values
135  _dcompliance_dp -= vec.dot(e.sol());
136 
137  vec.setZero();
138  e.side_external_residual(false,
139  vec,
140  dummy,
141  dummy,
143  e.volume_external_residual(false,
144  vec,
145  dummy,
146  dummy,
148 
149  // ask for the values
150  _dcompliance_dp -= vec.dot(e.sol(true));
151  }
152 }
153 
154 
155 
156 void
159 
160  // the primal data should have been calculated
161  libmesh_assert(_physics_elem);
162  libmesh_assert(f.is_topology_parameter());
163 
165 
166  std::pair<const MAST::FieldFunction<RealVectorX>*, unsigned int>
167  val = this->get_elem_boundary_velocity_data();
168 
169  if (val.first) {
170 
172  dynamic_cast<MAST::StructuralElementBase&>(*_physics_elem);
173 
175  vec = RealVectorX::Zero(e.sol().size());
176 
178  dummy = RealMatrixX::Zero(vec.size(), vec.size());
179 
181  val.second,
182  *val.first,
184  false,
185  vec,
186  dummy);
187 
188  // compute the contribution of this element to compliance
189  _dcompliance_dp -= vec.dot(e.sol());
190  }
191  }
192 }
193 
194 
195 
196 void
200 
201  // the primal data should have been calculated
202  libmesh_assert(_physics_elem);
203  libmesh_assert(f.is_topology_parameter());
204 
206  &elem = dynamic_cast<const MAST::LevelSetIntersectedElem&>(_physics_elem->elem());
207 
208  // sensitivity only exists at the boundary. So, we proceed with calculation
209  // only if this element has an intersection in the interior, or with a side.
210 
211  if (this->if_evaluate_for_element(elem) &&
212  elem.if_elem_has_level_set_boundary() &&
213  elem.if_subelem_has_side_on_level_set_boundary()) {
214 
216  dynamic_cast<MAST::StructuralElementBase&>(*_physics_elem);
217 
219  vec = RealVectorX::Zero(e.sol().size());
220 
222  dummy = RealMatrixX::Zero(vec.size(), vec.size());
223 
225  elem.get_subelem_side_on_level_set_boundary(),
226  vel,
228  false,
229  vec,
230  dummy);
231 
232  // compute the contribution of this element to compliance
233  _dcompliance_dp -= vec.dot(e.sol());
234  }
235 }
236 
237 
238 
239 Real
241 
242  Real val = _compliance;
243 
244  if (!_skip_comm_sum)
245  _system->system().comm().sum(val);
246 
247  return val;
248 }
249 
250 
251 
252 Real
254 
255  Real val = _dcompliance_dp;
256 
257  if (!_skip_comm_sum)
258  _system->system().comm().sum(val);
259 
260  return val;
261 }
262 
263 
264 
265 void
267 
268  // make sure that this has not been initialized ana calculated for all elems
269  libmesh_assert(_physics_elem);
270 
271  // since compliance = -X^T F, derivative wrt X = -F.
272 
274 
276  dynamic_cast<MAST::StructuralElementBase&>(*_physics_elem);
277 
278  dq_dX.setZero();
279 
281  dummy = RealMatrixX::Zero(dq_dX.size(), dq_dX.size());
282 
283  e.side_external_residual(false,
284  dq_dX,
285  dummy,
286  dummy,
288  e.volume_external_residual(false,
289  dq_dX,
290  dummy,
291  dummy,
293 
294  dq_dX *= -1.;
295  }
296 }
297 
298 
299 
300 void
302 set_elem_data(unsigned int dim,
303  const libMesh::Elem& ref_elem,
304  MAST::GeomElem& elem) const {
305 
306  libmesh_assert(!_physics_elem);
307 
308  if (dim == 1) {
309 
310  const MAST::ElementPropertyCard1D& p =
311  dynamic_cast<const MAST::ElementPropertyCard1D&>(_discipline->get_property_card(ref_elem));
312 
313  elem.set_local_y_vector(p.y_vector());
314  }
315 }
316 
317 
318 void
320 
321  libmesh_assert(!_physics_elem);
322  libmesh_assert(_assembly);
323  libmesh_assert(_system);
324 
326  dynamic_cast<const MAST::ElementPropertyCardBase&>(_discipline->get_property_card(elem));
327 
328  _physics_elem =
329  MAST::build_structural_element(*_system, elem, p).release();
330 }
331 
332 
virtual void zero_for_sensitivity()
zeroes the output quantity values stored inside this object so that assembly process can begin...
MAST::NonlinearSystem & system()
const MAST::ElementPropertyCardBase & get_property_card(const libMesh::Elem &elem) const
get property card for the specified element
virtual void init(const MAST::GeomElem &elem)
initialize for the element.
virtual bool if_evaluate_for_element(const MAST::GeomElem &elem) const
checks to see if the object has been told about the subset of elements and if the specified element i...
virtual void zero_for_analysis()
zeroes the output quantity values stored inside this object so that assembly process can begin...
const RealVectorX & sol(bool if_sens=false) const
Definition: elem_base.cpp:51
RealVectorX & y_vector()
returns value of the property val.
This provides the base class for definitin of element level contribution of output quantity in an ana...
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
libMesh::Real Real
MAST::PhysicsDisciplineBase * _discipline
virtual Real output_sensitivity_total(const MAST::FunctionBase &p)
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 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 bool is_topology_parameter() const
Definition: function_base.h:97
const MAST::SideBCMapType & side_loads() const
bool _skip_comm_sum
If an output has contrinutions only from local processor then the user can request that the global co...
virtual void evaluate_sensitivity(const MAST::FunctionBase &f)
this evaluates all relevant stress sensitivity components on the element to evaluate the p-averaged q...
Matrix< Real, Dynamic, Dynamic > RealMatrixX
ComplianceOutput()
default constructor
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
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 evaluate()
this evaluates all relevant stress components on the element to evaluate the p-averaged quantity...
const MAST::GeomElem & elem() const
Definition: elem_base.h:109
virtual void output_derivative_for_elem(RealVectorX &dq_dX)
calculates the derivative of p-norm von Mises stress for the norm identified using set_p_val()...
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
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 void evaluate_topology_sensitivity(const MAST::FunctionBase &f)
this evaluates all relevant topological sensitivity components on the element.
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.
MAST::SystemInitialization * _system