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-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 += pow(1. + pow(e_val/_sigma0, _p_norm_stress), 1./_p_norm_stress) * 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  JxW = 0.;
118 
119  dsigma_vm_val_df = 0.;
120 
121  // iterate over all element data
122  std::map<const libMesh::dof_id_type, std::vector<MAST::StressStrainOutputBase::Data*> >::const_iterator
123  map_it = _stress_data.find(e_id),
124  map_end = _stress_data.end();
125 
126  libmesh_assert(map_it != map_end);
127 
128  std::vector<MAST::StressStrainOutputBase::Data*>::const_iterator
129  vec_it = map_it->second.begin(),
130  vec_end = map_it->second.end();
131 
132  for ( ; vec_it != vec_end; vec_it++) {
133 
134  // ask this data point for the von Mises stress value
135  e_val = (*vec_it)->von_Mises_stress();
136  de_val = (*vec_it)->dvon_Mises_stress_dp(f);
137  JxW = (*vec_it)->quadrature_point_JxW();
138 
139  dsigma_vm_val_df +=
140  pow(1. + pow(e_val/_sigma0, _p_norm_stress), 1./_p_norm_stress-1.) *
141  pow(e_val/_sigma0, _p_norm_stress-1.) * de_val/_sigma0 * JxW;
142  }
143 }
144 
145 
146 
147 
148 
149 void
152  const libMesh::dof_id_type e_id,
153  Real& dsigma_vm_val_df) const {
154 
155  libmesh_assert(!_if_stress_plot_mode);
156  libmesh_assert(_primal_data_initialized);
157  libmesh_assert_greater(_sigma0, 0.);
158 
159  Real
160  e_val = 0.,
161  JxW_Vn = 0.;
162 
163  dsigma_vm_val_df = 0.;
164 
165  // iterate over all element data
166  std::map<const libMesh::dof_id_type, std::vector<MAST::StressStrainOutputBase::Data*> >::const_iterator
167  map_it = _boundary_stress_data.find(e_id),
168  map_end = _boundary_stress_data.end();
169 
170  libmesh_assert(map_it != map_end);
171 
172  std::vector<MAST::StressStrainOutputBase::Data*>::const_iterator
173  vec_it = map_it->second.begin(),
174  vec_end = map_it->second.end();
175 
176  for ( ; vec_it != vec_end; vec_it++) {
177 
178  // ask this data point for the von Mises stress value
179  e_val = (*vec_it)->von_Mises_stress();
180  JxW_Vn = (*vec_it)->quadrature_point_JxW();
181 
182  dsigma_vm_val_df +=
183  (pow(1. + pow(e_val/_sigma0, _p_norm_stress), 1./_p_norm_stress)-1.) * JxW_Vn;
184  }
185 }
186 
187 
188 
189 
190 void
192  RealVectorX& dq_dX) const {
193  libmesh_assert(!_if_stress_plot_mode);
194  libmesh_assert(_primal_data_initialized);
195  libmesh_assert_greater(_sigma0, 0.);
196 
197  Real
198  e_val = 0.,
199  JxW = 0.;
200 
202  de_val = RealVectorX::Zero(dq_dX.size());
203  dq_dX.setZero();
204 
205  // first find the data with the maximum value, to be used for scaling
206 
207  // iterate over all element data
208  std::map<const libMesh::dof_id_type, std::vector<MAST::StressStrainOutputBase::Data*> >::const_iterator
209  map_it = _stress_data.find(e_id),
210  map_end = _stress_data.end();
211 
212  // make sure that the data exists
213  libmesh_assert(map_it != map_end);
214 
215  std::vector<MAST::StressStrainOutputBase::Data*>::const_iterator
216  vec_it = map_it->second.begin(),
217  vec_end = map_it->second.end();
218 
219 
220  for ( ; vec_it != vec_end; vec_it++) {
221 
222  // ask this data point for the von Mises stress value
223  e_val = (*vec_it)->von_Mises_stress();
224  de_val = (*vec_it)->dvon_Mises_stress_dX();
225  JxW = (*vec_it)->quadrature_point_JxW();
226 
227  dq_dX +=
228  pow(1. + pow(e_val/_sigma0, _p_norm_stress), 1./_p_norm_stress-1.) *
229  pow(e_val/_sigma0, _p_norm_stress-1.) * de_val/_sigma0 * JxW;
230  }
231 }
232 
233 
234 
235 
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.
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_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...
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...
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...
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