MAST
Multidisciplinary-design Adaptation and Sensitivity Toolkit (MAST)
smooth_ramp_stress_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 
37 
38 
41 
42 }
43 
44 
45 
46 
48 
49 }
50 
51 
52 
53 
54 void
56 
57  libmesh_assert(!_if_stress_plot_mode);
58  libmesh_assert(!_primal_data_initialized);
59  libmesh_assert_greater(_sigma0, 0.);
60 
61  Real
62  e_val = 0.,
63  JxW = 0.;
64 
65  _JxW_val = 0.;
66  _sigma_vm_int = 0.;
67  _sigma_vm_p_norm = 0.;
68 
69  // first find the data with the maximum value, to be used for scaling
70 
71  // iterate over all element data
72  std::map<const libMesh::dof_id_type, std::vector<MAST::StressStrainOutputBase::Data*> >::const_iterator
73  map_it = _stress_data.begin(),
74  map_end = _stress_data.end();
75 
76  for ( ; map_it != map_end; map_it++) {
77 
78  std::vector<MAST::StressStrainOutputBase::Data*>::const_iterator
79  vec_it = map_it->second.begin(),
80  vec_end = map_it->second.end();
81 
82  for ( ; vec_it != vec_end; vec_it++) {
83 
84  // ask this data point for the von Mises stress value
85  e_val = (*vec_it)->von_Mises_stress();
86  JxW = (*vec_it)->quadrature_point_JxW();
87 
88  // we do not use absolute value here, since von Mises stress
89  // is >= 0.
90  _sigma_vm_int += pow(1. + pow(e_val/_sigma0, _p_norm_stress), 1./_p_norm_stress) * JxW;
91  _JxW_val += JxW;
92  }
93  }
94 
95  // sum over all processors, since part of the mesh will exist on the
96  // other processors.
97  if (!_skip_comm_sum) {
98 
99  _system->system().comm().sum(_sigma_vm_int);
100  _system->system().comm().sum(_JxW_val);
101  }
102 
105 }
106 
107 
108 
109 
110 void
113  const libMesh::dof_id_type e_id,
114  Real& dsigma_vm_val_df) const {
115 
116  libmesh_assert(!_if_stress_plot_mode);
117  libmesh_assert(_primal_data_initialized);
118  libmesh_assert_greater(_sigma0, 0.);
119 
120  Real
121  e_val = 0.,
122  de_val = 0.,
123  JxW = 0.;
124 
125  dsigma_vm_val_df = 0.;
126 
127  // iterate over all element data
128  std::map<const libMesh::dof_id_type, std::vector<MAST::StressStrainOutputBase::Data*> >::const_iterator
129  map_it = _stress_data.find(e_id),
130  map_end = _stress_data.end();
131 
132  libmesh_assert(map_it != map_end);
133 
134  std::vector<MAST::StressStrainOutputBase::Data*>::const_iterator
135  vec_it = map_it->second.begin(),
136  vec_end = map_it->second.end();
137 
138  for ( ; vec_it != vec_end; vec_it++) {
139 
140  // ask this data point for the von Mises stress value
141  e_val = (*vec_it)->von_Mises_stress();
142  de_val = (*vec_it)->dvon_Mises_stress_dp(f);
143  JxW = (*vec_it)->quadrature_point_JxW();
144 
145  dsigma_vm_val_df +=
146  pow(1. + pow(e_val/_sigma0, _p_norm_stress), 1./_p_norm_stress-1.) *
147  pow(e_val/_sigma0, _p_norm_stress-1.) * de_val/_sigma0 * JxW;
148  }
149 }
150 
151 
152 
153 
154 
155 void
158  const libMesh::dof_id_type e_id,
159  Real& dsigma_vm_val_df) const {
160 
161  libmesh_assert(!_if_stress_plot_mode);
162  libmesh_assert(_primal_data_initialized);
163  libmesh_assert_greater(_sigma0, 0.);
164 
165  Real
166  e_val = 0.,
167  JxW_Vn = 0.;
168 
169  dsigma_vm_val_df = 0.;
170 
171  // iterate over all element data
172  std::map<const libMesh::dof_id_type, std::vector<MAST::StressStrainOutputBase::Data*> >::const_iterator
173  map_it = _boundary_stress_data.find(e_id),
174  map_end = _boundary_stress_data.end();
175 
176  libmesh_assert(map_it != map_end);
177 
178  std::vector<MAST::StressStrainOutputBase::Data*>::const_iterator
179  vec_it = map_it->second.begin(),
180  vec_end = map_it->second.end();
181 
182  for ( ; vec_it != vec_end; vec_it++) {
183 
184  // ask this data point for the von Mises stress value
185  e_val = (*vec_it)->von_Mises_stress();
186  JxW_Vn = (*vec_it)->quadrature_point_JxW();
187 
188  dsigma_vm_val_df +=
189  (pow(1. + pow(e_val/_sigma0, _p_norm_stress), 1./_p_norm_stress)-1.) * JxW_Vn;
190  }
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  de_val = RealVectorX::Zero(dq_dX.size());
209  dq_dX.setZero();
210 
211  // first find the data with the maximum value, to be used for scaling
212 
213  // iterate over all element data
214  std::map<const libMesh::dof_id_type, std::vector<MAST::StressStrainOutputBase::Data*> >::const_iterator
215  map_it = _stress_data.find(e_id),
216  map_end = _stress_data.end();
217 
218  // make sure that the data exists
219  libmesh_assert(map_it != map_end);
220 
221  std::vector<MAST::StressStrainOutputBase::Data*>::const_iterator
222  vec_it = map_it->second.begin(),
223  vec_end = map_it->second.end();
224 
225 
226  for ( ; vec_it != vec_end; vec_it++) {
227 
228  // ask this data point for the von Mises stress value
229  e_val = (*vec_it)->von_Mises_stress();
230  de_val = (*vec_it)->dvon_Mises_stress_dX();
231  JxW = (*vec_it)->quadrature_point_JxW();
232 
233  dq_dX +=
234  pow(1. + pow(e_val/_sigma0, _p_norm_stress), 1./_p_norm_stress-1.) *
235  pow(e_val/_sigma0, _p_norm_stress-1.) * de_val/_sigma0 * JxW;
236  }
237 }
238 
239 
240 
241 
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
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.
bool _skip_comm_sum
If an output has contrinutions only from local processor then the user can request that the global co...
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_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...
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...
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...
Matrix< Real, Dynamic, 1 > RealVectorX
bool _if_stress_plot_mode
identifies the mode in which evaluation is peformed.
virtual void functional_for_all_elems()
calculates and returns the von Mises p-norm functional for all the elements that this object currentl...
MAST::SystemInitialization * _system