MAST
stress_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
25 #include "base/nonlinear_system.h"
26 #include "mesh/geom_elem.h"
27 
28 // libMesh includes
29 #include "libmesh/system.h"
30 #include "libmesh/dof_map.h"
31 #include "libmesh/numeric_vector.h"
32 
34 MAST::AssemblyBase() {
35 
36 }
37 
38 
40 
41 }
42 
43 
44 
45 void
46 get_max_stress_strain_values(const std::vector<MAST::StressStrainOutputBase::Data*>& data,
47  RealVectorX& max_strain,
48  RealVectorX& max_stress,
49  Real& max_vm,
50  const MAST::FunctionBase* p) {
51 
52  max_strain = RealVectorX::Zero(6);
53  max_stress = RealVectorX::Zero(6);
54  max_vm = 0.;
55 
56  // if there is only one data point, then simply copy the value to the output
57  // routines
58  if (data.size() == 1) {
59  if (p == nullptr) {
60  max_strain = data[0]->strain();
61  max_stress = data[0]->stress();
62  max_vm = data[0]->von_Mises_stress();
63  }
64  else {
65  max_strain = data[0]->get_strain_sensitivity(*p);
66  max_stress = data[0]->get_stress_sensitivity(*p);
67  max_vm = data[0]->dvon_Mises_stress_dp (*p);
68  }
69 
70  return;
71  }
72 
73  // if multiple values are provided for an element, then we need to compare
74  std::vector<MAST::StressStrainOutputBase::Data*>::const_iterator
75  it = data.begin(),
76  end = data.end();
77 
78  Real
79  vm = 0.;
80 
81  for ( ; it != end; it++) {
82 
83  // get the strain value at this point
84  const RealVectorX& strain = (*it)->strain();
85  const RealVectorX& stress = (*it)->stress();
86  vm = (*it)->von_Mises_stress();
87 
88  // now compare
89  if (vm > max_vm) max_vm = vm;
90 
91  for ( unsigned int i=0; i<6; i++) {
92  if (fabs(strain(i)) > fabs(max_strain(i))) max_strain(i) = strain(i);
93  if (fabs(stress(i)) > fabs(max_stress(i))) max_stress(i) = stress(i);
94  }
95  }
96 }
97 
98 
99 
100 
103  const libMesh::NumericVector<Real>& X) {
104 
105  // make sure it has been initialized
106  libmesh_assert(_system);
107  libmesh_assert(_discipline);
108  libmesh_assert(!_elem_ops);
109 
110  this->set_elem_operation_object(ops);
111 
112  RealVectorX sol;
113 
114  const MAST::NonlinearSystem&
115  structural_sys = _system->system();
116 
117  libMesh::System&
118  stress_sys = dynamic_cast<MAST::StructuralSystemInitialization*>(_system)->get_stress_sys();
119  const std::vector<unsigned int>&
120  stress_vars = dynamic_cast<MAST::StructuralSystemInitialization*>(_system)->get_stress_var_ids();
121 
122  stress_sys.solution->zero();
123 
124  std::vector<libMesh::dof_id_type> dof_indices;
125  const libMesh::DofMap& dof_map = structural_sys.get_dof_map();
126 
127  std::unique_ptr<libMesh::NumericVector<Real> > localized_solution;
128  localized_solution.reset(this->build_localized_vector(structural_sys,
129  X).release());
130 
131 
132  libMesh::MeshBase::const_element_iterator el =
133  structural_sys.get_mesh().active_local_elements_begin();
134  const libMesh::MeshBase::const_element_iterator end_el =
135  structural_sys.get_mesh().active_local_elements_end();
136 
137  libMesh::dof_id_type
138  dof_id = 0;
139 
140  unsigned int
141  sys_num = stress_sys.number();
142 
144  max_strain_vals = RealVectorX::Zero(6),
145  max_stress_vals = RealVectorX::Zero(6);
146 
147  Real
148  max_vm_stress = 0.;
149 
150  for ( ; el != end_el; el++) {
151 
152  const libMesh::Elem* elem = *el;
153 
154  dof_map.dof_indices (elem, dof_indices);
155 
156  unsigned int ndofs = (unsigned int)dof_indices.size();
157  sol.setZero(ndofs);
158 
159  for (unsigned int i=0; i<dof_indices.size(); i++)
160  sol(i) = (*localized_solution)(dof_indices[i]);
161 
162  // clear before calculating the data
163  ops.clear();
164  MAST::GeomElem geom_elem;
165  ops.set_elem_data(elem->dim(), *elem, geom_elem);
166  geom_elem.init(*elem, *_system);
167 
168  ops.init(geom_elem);
169  ops.set_elem_solution(sol);
170  ops.evaluate();
171  ops.clear_elem();
172 
173  // get the stress-strain data map from the object
174  const std::map<const libMesh::dof_id_type,
175  std::vector<MAST::StressStrainOutputBase::Data*> >& output_map =
177 
178  // make sure that only one element has been added to this data,
179  // and that the element id is the same as the one being computed
180  libmesh_assert_equal_to(output_map.size(), 1);
181  libmesh_assert_equal_to(output_map.begin()->first, elem->id());
182 
183  // now iterate over all the elements and set the value in the
184  // new system used for output
185  std::map<const libMesh::dof_id_type,
186  std::vector<MAST::StressStrainOutputBase::Data*> >::const_iterator
187  e_it = output_map.begin(),
188  e_end = output_map.end();
189 
190  for ( ; e_it != e_end; e_it++) {
191 
192  get_max_stress_strain_values(e_it->second,
193  max_strain_vals,
194  max_stress_vals,
195  max_vm_stress,
196  nullptr);
197 
198  // set the values in the system
199  // stress value
200  dof_id = elem->dof_number(sys_num, stress_vars[12], 0);
201  stress_sys.solution->set(dof_id, max_vm_stress);
202 
203  for (unsigned int i=0; i<6; i++) {
204  // strain value
205  dof_id = elem->dof_number(sys_num, stress_vars[i], 0);
206  stress_sys.solution->set(dof_id, max_strain_vals(i));
207 
208  // stress value
209  dof_id = elem->dof_number(sys_num, stress_vars[i+6], 0);
210  stress_sys.solution->set(dof_id, max_stress_vals(i));
211  }
212  }
213  }
214 
215  stress_sys.solution->close();
217 }
218 
219 
220 
221 
222 
225  const libMesh::NumericVector<Real>& X,
226  const libMesh::NumericVector<Real>& dXdp,
227  const MAST::FunctionBase& p,
228  libMesh::NumericVector<Real>& dsigmadp) {
229 
230  // make sure it has been initialized
231  libmesh_assert(_system);
232  libmesh_assert(_discipline);
233  libmesh_assert(!_elem_ops);
234 
235  this->set_elem_operation_object(ops);
236 
238  sol,
239  dsol;
240 
241  const MAST::NonlinearSystem&
242  structural_sys = _system->system();
243 
244  libMesh::System&
245  stress_sys = dynamic_cast<MAST::StructuralSystemInitialization*>(_system)->get_stress_sys();
246  const std::vector<unsigned int>&
247  stress_vars = dynamic_cast<MAST::StructuralSystemInitialization*>(_system)->get_stress_var_ids();
248 
249  dsigmadp.zero();
250 
251  std::vector<libMesh::dof_id_type> dof_indices;
252  const libMesh::DofMap& dof_map = structural_sys.get_dof_map();
253 
254  std::unique_ptr<libMesh::NumericVector<Real> >
255  localized_solution,
256  localized_solution_sens;
257  localized_solution.reset(this->build_localized_vector(structural_sys,
258  X).release());
259  localized_solution_sens.reset(this->build_localized_vector(structural_sys,
260  dXdp).release());
261 
262 
263  libMesh::MeshBase::const_element_iterator el =
264  structural_sys.get_mesh().active_local_elements_begin();
265  const libMesh::MeshBase::const_element_iterator end_el =
266  structural_sys.get_mesh().active_local_elements_end();
267 
268  libMesh::dof_id_type
269  dof_id = 0;
270 
271  unsigned int
272  sys_num = stress_sys.number();
273 
275  max_strain_vals = RealVectorX::Zero(6),
276  max_stress_vals = RealVectorX::Zero(6);
277 
278  Real
279  max_vm_stress = 0.;
280 
281  for ( ; el != end_el; el++) {
282 
283  const libMesh::Elem* elem = *el;
284 
285  dof_map.dof_indices (elem, dof_indices);
286 
287  unsigned int ndofs = (unsigned int)dof_indices.size();
288  sol.setZero(ndofs);
289  dsol.setZero(ndofs);
290 
291  for (unsigned int i=0; i<dof_indices.size(); i++) {
292  sol (i) = (*localized_solution) (dof_indices[i]);
293  dsol(i) = (*localized_solution_sens)(dof_indices[i]);
294  }
295 
296  // clear before calculating the data. We have to calculate the stress
297  // before calculating sensitivity since the sensitivity assumes the
298  // presence of stress data, which is cleared after each element.
299  ops.clear();
300  MAST::GeomElem geom_elem;
301  ops.set_elem_data(elem->dim(), *elem, geom_elem);
302  geom_elem.init(*elem, *_system);
303 
304  ops.init(geom_elem);
305  ops.set_stress_plot_mode(true);
306  ops.set_elem_solution(sol);
308  ops.evaluate();
309  ops.evaluate_sensitivity(p);
310  ops.clear_elem();
311 
312  // get the stress-strain data map from the object
313  const std::map<const libMesh::dof_id_type,
314  std::vector<MAST::StressStrainOutputBase::Data*> >& output_map =
316 
317  // make sure that only one element has been added to this data,
318  // and that the element id is the same as the one being computed
319  libmesh_assert_equal_to(output_map.size(), 1);
320  libmesh_assert_equal_to(output_map.begin()->first, elem->id());
321 
322  // now iterate over all the elements and set the value in the
323  // new system used for output
324  std::map<const libMesh::dof_id_type,
325  std::vector<MAST::StressStrainOutputBase::Data*> >::const_iterator
326  e_it = output_map.begin(),
327  e_end = output_map.end();
328 
329  for ( ; e_it != e_end; e_it++) {
330 
331  get_max_stress_strain_values(e_it->second,
332  max_strain_vals,
333  max_stress_vals,
334  max_vm_stress,
335  &p);
336 
337  // set the values in the system
338  // stress value
339  dof_id = elem->dof_number(sys_num, stress_vars[12], 0);
340  dsigmadp.set(dof_id, max_vm_stress);
341 
342  for (unsigned int i=0; i<6; i++) {
343  // strain value
344  dof_id = elem->dof_number(sys_num, stress_vars[i], 0);
345  dsigmadp.set(dof_id, max_strain_vals(i));
346 
347  // stress value
348  dof_id = elem->dof_number(sys_num, stress_vars[i+6], 0);
349  dsigmadp.set(dof_id, max_stress_vals(i));
350  }
351  }
352  }
353 
354  dsigmadp.close();
355  ops.set_stress_plot_mode(false);
357 }
358 
359 
360 
MAST::AssemblyElemOperations * _elem_ops
provides assembly elem operations for use by this class
MAST::NonlinearSystem & system()
virtual void evaluate_sensitivity(const MAST::FunctionBase &f)
this evaluates all relevant stress sensitivity components on the element to evaluate the p-averaged q...
StressAssembly()
constructor associates this assembly object with the system
Data structure provides the mechanism to store stress and strain output from a structural analysis...
virtual ~StressAssembly()
destructor resets the association of this assembly object with the system
This class implements a system for solution of nonlinear systems.
virtual void clear_elem_operation_object()
clears the association of this object with the assembly element operation object. ...
void clear()
clears the data structure of any stored values so that it can be used for another element...
virtual void set_elem_solution_sensitivity(const RealVectorX &sol)
sets the element solution sensitivity
void set_stress_plot_mode(bool f)
tells the object that the calculation is for stress to be output for plotting.
virtual void set_elem_operation_object(MAST::AssemblyElemOperations &elem_ops)
attaches a element operation to this object, and associated this with the element operation object...
virtual void set_elem_solution(const RealVectorX &sol)
sets the element solution
void get_max_stress_strain_values(const std::vector< MAST::StressStrainOutputBase::Data * > &data, RealVectorX &max_strain, RealVectorX &max_stress, Real &max_vm, const MAST::FunctionBase *p)
libMesh::Real Real
MAST::SystemInitialization * _system
System for which this assembly is performed.
virtual void update_stress_strain_data(MAST::StressStrainOutputBase &ops, const libMesh::NumericVector< Real > &X)
updates the stresses and strains for the specified solution vector X.
MAST::PhysicsDisciplineBase * _discipline
PhysicsDisciplineBase object for which this class is assembling.
virtual const std::map< const libMesh::dof_id_type, std::vector< MAST::StressStrainOutputBase::Data * > > & get_stress_strain_data() const
virtual void init(const MAST::GeomElem &elem)
initialize for the element.
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.
Matrix< Real, Dynamic, 1 > RealVectorX
std::unique_ptr< libMesh::NumericVector< Real > > build_localized_vector(const libMesh::System &sys, const libMesh::NumericVector< Real > &global) const
localizes the parallel vector so that the local copy stores all values necessary for calculation of t...
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 evaluate()
this evaluates all relevant stress components on the element to evaluate the p-averaged quantity...
virtual void clear_elem()
clears the element initialization
virtual void update_stress_strain_sensitivity_data(MAST::StressStrainOutputBase &ops, const libMesh::NumericVector< Real > &X, const libMesh::NumericVector< Real > &dXdp, const MAST::FunctionBase &p, libMesh::NumericVector< Real > &dsigmadp)
updates the sensitivity of stresses and strains for the specified solution vector X and its sensitivi...
virtual void init(const libMesh::Elem &elem, const MAST::SystemInitialization &sys_init)
initialize the object for the specified reference elem.
Definition: geom_elem.cpp:128