MAST
ks_stress_output.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
24 #include "base/assembly_base.h"
26 #include "base/nonlinear_system.h"
32 #include "mesh/geom_elem.h"
33 
34 
35 
38 
39 }
40 
41 
42 
43 
45 
46 }
47 
48 
49 
50 
51 void
53 
54  libmesh_assert(!_if_stress_plot_mode);
55  libmesh_assert(!_primal_data_initialized);
56  libmesh_assert_greater(_sigma0, 0.);
57 
58  Real
59  e_val = 0.,
60  JxW = 0.;
61 
62  _JxW_val = 0.;
63  _sigma_vm_int = 0.;
64  _sigma_vm_p_norm = 0.;
65 
66  // first find the data with the maximum value, to be used for scaling
67 
68  // iterate over all element data
69  std::map<const libMesh::dof_id_type, std::vector<MAST::StressStrainOutputBase::Data*> >::const_iterator
70  map_it = _stress_data.begin(),
71  map_end = _stress_data.end();
72 
73  for ( ; map_it != map_end; map_it++) {
74 
75  std::vector<MAST::StressStrainOutputBase::Data*>::const_iterator
76  vec_it = map_it->second.begin(),
77  vec_end = map_it->second.end();
78 
79  for ( ; vec_it != vec_end; vec_it++) {
80 
81  // ask this data point for the von Mises stress value
82  e_val = (*vec_it)->von_Mises_stress();
83  JxW = (*vec_it)->quadrature_point_JxW();
84 
85  // we do not use absolute value here, since von Mises stress
86  // is >= 0.
87  _sigma_vm_int += exp(_p_norm_stress * (e_val-_sigma0)/_sigma0) * JxW;
88  _JxW_val += JxW;
89  }
90  }
91 
92  // sum over all processors, since part of the mesh will exist on the
93  // other processors.
94  _system->system().comm().sum(_sigma_vm_int);
95  _system->system().comm().sum(_JxW_val);
96 
99 }
100 
101 
102 
103 
104 void
107  const libMesh::dof_id_type e_id,
108  Real& dsigma_vm_val_df) const {
109 
110  libmesh_assert(!_if_stress_plot_mode);
111  libmesh_assert(_primal_data_initialized);
112  libmesh_assert_greater(_sigma0, 0.);
113 
114  Real
115  e_val = 0.,
116  de_val = 0.,
117  num_sens = 0.,
118  JxW = 0.;
119 
120  dsigma_vm_val_df = 0.;
121 
122  // iterate over all element data
123  std::map<const libMesh::dof_id_type, std::vector<MAST::StressStrainOutputBase::Data*> >::const_iterator
124  map_it = _stress_data.find(e_id),
125  map_end = _stress_data.end();
126 
127  libmesh_assert(map_it != map_end);
128 
129  std::vector<MAST::StressStrainOutputBase::Data*>::const_iterator
130  vec_it = map_it->second.begin(),
131  vec_end = map_it->second.end();
132 
133  for ( ; vec_it != vec_end; vec_it++) {
134 
135  // ask this data point for the von Mises stress value
136  e_val = (*vec_it)->von_Mises_stress();
137  de_val = (*vec_it)->dvon_Mises_stress_dp(f);
138  JxW = (*vec_it)->quadrature_point_JxW();
139 
140  num_sens += _p_norm_stress * de_val/_sigma0 * exp(_p_norm_stress * (e_val-_sigma0)/_sigma0) * JxW;
141  }
142 
143  dsigma_vm_val_df = 1./_p_norm_stress / (_sigma_vm_int/_JxW_val) * num_sens / _JxW_val;
144 }
145 
146 
147 
148 
149 
150 void
153  const libMesh::dof_id_type e_id,
154  Real& dsigma_vm_val_df) const {
155 
156  libmesh_assert(!_if_stress_plot_mode);
157  libmesh_assert(_primal_data_initialized);
158  libmesh_assert_greater(_sigma0, 0.);
159 
160  Real
161  e_val = 0.,
162  JxW_Vn = 0.,
163  num_sens = 0.,
164  denom_sens = 0.;
165 
166  dsigma_vm_val_df = 0.;
167 
168  // iterate over all element data
169  std::map<const libMesh::dof_id_type, std::vector<MAST::StressStrainOutputBase::Data*> >::const_iterator
170  map_it = _boundary_stress_data.find(e_id),
171  map_end = _boundary_stress_data.end();
172 
173  libmesh_assert(map_it != map_end);
174 
175  std::vector<MAST::StressStrainOutputBase::Data*>::const_iterator
176  vec_it = map_it->second.begin(),
177  vec_end = map_it->second.end();
178 
179  for ( ; vec_it != vec_end; vec_it++) {
180 
181  // ask this data point for the von Mises stress value
182  e_val = (*vec_it)->von_Mises_stress();
183  JxW_Vn = (*vec_it)->quadrature_point_JxW();
184 
185  denom_sens += JxW_Vn;
186  num_sens += exp(_p_norm_stress * (e_val-_sigma0)/_sigma0) * JxW_Vn;
187  }
188 
189  dsigma_vm_val_df = 1./_p_norm_stress / (_sigma_vm_int/_JxW_val) *
190  (num_sens / _JxW_val - _sigma_vm_int / pow(_JxW_val, 2.) * denom_sens);
191 }
192 
193 
194 
195 
196 void
198  RealVectorX& dq_dX) const {
199  libmesh_assert(!_if_stress_plot_mode);
200  libmesh_assert(_primal_data_initialized);
201  libmesh_assert_greater(_sigma0, 0.);
202 
203  Real
204  e_val = 0.,
205  JxW = 0.;
206 
208  num_sens = RealVectorX::Zero(dq_dX.size()),
209  de_val = RealVectorX::Zero(dq_dX.size());
210  dq_dX.setZero();
211 
212  // first find the data with the maximum value, to be used for scaling
213 
214  // iterate over all element data
215  std::map<const libMesh::dof_id_type, std::vector<MAST::StressStrainOutputBase::Data*> >::const_iterator
216  map_it = _stress_data.find(e_id),
217  map_end = _stress_data.end();
218 
219  // make sure that the data exists
220  libmesh_assert(map_it != map_end);
221 
222  std::vector<MAST::StressStrainOutputBase::Data*>::const_iterator
223  vec_it = map_it->second.begin(),
224  vec_end = map_it->second.end();
225 
226 
227  for ( ; vec_it != vec_end; vec_it++) {
228 
229  // ask this data point for the von Mises stress value
230  e_val = (*vec_it)->von_Mises_stress();
231  de_val = (*vec_it)->dvon_Mises_stress_dX();
232  JxW = (*vec_it)->quadrature_point_JxW();
233 
234  num_sens += _p_norm_stress * de_val/_sigma0 * exp(_p_norm_stress * (e_val-_sigma0)/_sigma0) * JxW;
235  }
236 
237  dq_dX = 1./_p_norm_stress / (_sigma_vm_int/_JxW_val) * num_sens / _JxW_val;
238 }
239 
240 
241 
242 
MAST::NonlinearSystem & system()
Data structure provides the mechanism to store stress and strain output from a structural analysis...
std::map< const libMesh::dof_id_type, std::vector< MAST::StressStrainOutputBase::Data * > > _boundary_stress_data
vector of stress with the associated location details
virtual void functional_boundary_sensitivity_for_elem(const MAST::FunctionBase &f, const libMesh::dof_id_type e_id, Real &dsigma_vm_val_df) const
calculates and returns the boundary sensitivity of von Mises p-norm functional for the element e...
Real _sigma0
reference stress value used in scaling volume.
libMesh::Real Real
bool _primal_data_initialized
primal data, needed for sensitivity and adjoints
Real _p_norm_stress
norm to be used for calculation of output stress function.
KSStressStrainOutput()
default constructor
std::map< const libMesh::dof_id_type, std::vector< MAST::StressStrainOutputBase::Data * > > _stress_data
vector of stress with the associated location details
virtual void functional_for_all_elems()
calculates and returns the von Mises p-norm functional for all the elements that this object currentl...
Matrix< Real, Dynamic, 1 > RealVectorX
virtual void functional_sensitivity_for_elem(const MAST::FunctionBase &f, const libMesh::dof_id_type e_id, Real &dsigma_vm_val_df) const
calculates and returns the sensitivity of von Mises p-norm functional for the element e...
bool _if_stress_plot_mode
identifies the mode in which evaluation is peformed.
virtual void functional_state_derivartive_for_elem(const libMesh::dof_id_type e_id, RealVectorX &dq_dX) const
calculates and returns the derivative of von Mises p-norm functional wrt state vector for the specifi...
MAST::SystemInitialization * _system