MAST
stress_output_base.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 
36  const RealVectorX& strain,
37  const libMesh::Point& qp,
38  const libMesh::Point& xyz,
39  Real JxW):
40 _stress(stress),
41 _strain(strain),
42 _qp(qp),
43 _xyz(xyz),
44 _JxW(JxW) {
45 
46  // make sure that both the stress and strain are for a 3D configuration,
47  // which is the default for this data structure
48  libmesh_assert_equal_to(stress.size(), 6);
49  libmesh_assert_equal_to(strain.size(), 6);
50 
51 }
52 
53 
54 
55 void
57 
58  _stress_sensitivity.clear();
59  _strain_sensitivity.clear();
60 }
61 
62 
63 const libMesh::Point&
66 
67  return _qp;
68 }
69 
70 
71 const RealVectorX&
73 
74  return _stress;
75 }
76 
77 
78 
79 const RealVectorX&
81  return _strain;
82 }
83 
84 
85 
86 void
88  const RealMatrixX& dstrain_dX) {
89 
90  // make sure that the number of rows is 6.
91  libmesh_assert_equal_to(dstress_dX.rows(), 6);
92  libmesh_assert_equal_to(dstrain_dX.rows(), 6);
93 
94  _dstress_dX = dstress_dX;
95  _dstrain_dX = dstrain_dX;
96 }
97 
98 
99 
100 const RealMatrixX&
102 
103  return _dstress_dX;
104 }
105 
106 
107 const RealMatrixX&
109 
110  return _dstrain_dX;
111 }
112 
113 
114 Real
116 
117  return _JxW;
118 }
119 
120 
121 void
123  const RealVectorX& dstress_df,
124  const RealVectorX& dstrain_df) {
125 
126  // make sure that both the stress and strain are for a 3D configuration,
127  // which is the default for this data structure
128  libmesh_assert_equal_to(dstress_df.size(), 6);
129  libmesh_assert_equal_to(dstrain_df.size(), 6);
130 
131  _stress_sensitivity[&f] = dstress_df;
132  _strain_sensitivity[&f] = dstrain_df;
133 }
134 
135 
136 bool
139 
140  // make sure that the data exists
141  std::map<const MAST::FunctionBase*, RealVectorX>::const_iterator
142  it = _stress_sensitivity.find(&f);
143 
144  return it != _stress_sensitivity.end();
145 }
146 
147 
148 const RealVectorX&
151 
152  // make sure that the data exists
153  std::map<const MAST::FunctionBase*, RealVectorX>::const_iterator
154  it = _stress_sensitivity.find(&f);
155 
156  libmesh_assert(it != _stress_sensitivity.end());
157 
158 
159  return it->second;
160 }
161 
162 
163 
164 const RealVectorX&
167 
168  // make sure that the data exists
169  std::map<const MAST::FunctionBase*, RealVectorX>::const_iterator
170  it = _strain_sensitivity.find(&f);
171 
172  libmesh_assert(it != _strain_sensitivity.end());
173 
174 
175  return it->second;
176 }
177 
178 
179 
180 Real
182 
183  // make sure that the data is available
184  libmesh_assert_equal_to(_stress.size(), 6);
185 
186  return
187  pow(0.5 * (pow(_stress(0)-_stress(1),2) + //(((sigma_xx - sigma_yy)^2 +
188  pow(_stress(1)-_stress(2),2) + // (sigma_yy - sigma_zz)^2 +
189  pow(_stress(2)-_stress(0),2)) + // (sigma_zz - sigma_xx)^2)/2 +
190  3.0 * (pow(_stress(3), 2) + // 3* (tau_xx^2 +
191  pow(_stress(4), 2) + // tau_yy^2 +
192  pow(_stress(5), 2)), 0.5); // tau_zz^2))^.5
193 }
194 
195 
196 
199 
200  // make sure that the data is available
201  libmesh_assert_equal_to(_stress.size(), 6);
202  libmesh_assert_equal_to(_dstress_dX.rows(), 6);
203 
204  Real
205  p =
206  0.5 * (pow(_stress(0)-_stress(1),2) + //((sigma_xx - sigma_yy)^2 +
207  pow(_stress(1)-_stress(2),2) + // (sigma_yy - sigma_zz)^2 +
208  pow(_stress(2)-_stress(0),2)) + // (sigma_zz - sigma_xx)^2)/2 +
209  3.0 * (pow(_stress(3), 2) + // 3* (tau_xx^2 +
210  pow(_stress(4), 2) + // tau_yy^2 +
211  pow(_stress(5), 2)); // tau_zz^2)
212 
214  dp = RealVectorX::Zero(_dstress_dX.cols());
215 
216  // if p == 0, then the sensitivity returns nan
217  // Hence, we are avoiding this by setting it to zero whenever p = 0.
218  if (fabs(p) > 0.)
219  dp =
220  (((_dstress_dX.row(0) - _dstress_dX.row(1)) * (_stress(0) - _stress(1)) +
221  (_dstress_dX.row(1) - _dstress_dX.row(2)) * (_stress(1) - _stress(2)) +
222  (_dstress_dX.row(2) - _dstress_dX.row(0)) * (_stress(2) - _stress(0))) +
223  6.0 * (_dstress_dX.row(3) * _stress(3)+
224  _dstress_dX.row(4) * _stress(4)+
225  _dstress_dX.row(5) * _stress(5))) * 0.5 * pow(p, -0.5);
226 
227  return dp;
228 }
229 
230 
231 
232 
233 Real
236 
237  // make sure that the data is available
238  libmesh_assert_equal_to(_stress.size(), 6);
239 
240  if (!this->has_stress_sensitivity(f))
241  return 0.;
242 
243  // get the stress sensitivity data
244  const RealVectorX dstress_dp = this->get_stress_sensitivity(f);
245 
246 
247  Real
248  p =
249  0.5 * (pow(_stress(0)-_stress(1),2) + //((sigma_xx - sigma_yy)^2 +
250  pow(_stress(1)-_stress(2),2) + // (sigma_yy - sigma_zz)^2 +
251  pow(_stress(2)-_stress(0),2)) + // (sigma_zz - sigma_xx)^2)/2 +
252  3.0 * (pow(_stress(3), 2) + // 3* (tau_xx^2 +
253  pow(_stress(4), 2) + // tau_yy^2 +
254  pow(_stress(5), 2)), // tau_zz^2)
255  dp = 0.;
256 
257  // if p == 0, then the sensitivity returns nan
258  // Hence, we are avoiding this by setting it to zero whenever p = 0.
259  if (fabs(p) > 0.)
260  dp =
261  (((dstress_dp(0) - dstress_dp(1)) * (_stress(0) - _stress(1)) +
262  (dstress_dp(1) - dstress_dp(2)) * (_stress(1) - _stress(2)) +
263  (dstress_dp(2) - dstress_dp(0)) * (_stress(2) - _stress(0))) +
264  6.0 * (dstress_dp(3) * _stress(3)+
265  dstress_dp(4) * _stress(4)+
266  dstress_dp(5) * _stress(5))) * 0.5 * pow(p, -0.5);
267 
268  return dp;
269 }
270 
271 
272 
275 _p_norm_stress (2.),
276 _p_norm_weight (1),
277 _rho (1.),
278 _sigma0 (0.),
279 _exp_arg_lim (1.e2),
281 _JxW_val (0.),
282 _sigma_vm_int (0.),
283 _sigma_vm_p_norm (0.),
284 _if_stress_plot_mode (false) {
285 
286 }
287 
288 
289 
290 
292 
293  this->clear();
294 }
295 
296 
297 
298 
299 void
301 
302  _JxW_val = 0.;
303  _sigma_vm_int = 0.;
304  _sigma_vm_p_norm = 0.;
305 }
306 
307 
308 
309 void
311 
312  // nothign to be done here.
313 }
314 
315 
316 void
318 
319  // make sure that this has not been initialized ana calculated for all elems
320  libmesh_assert(_physics_elem);
321  libmesh_assert(!_primal_data_initialized);
322 
324  // ask for the values
325  dynamic_cast<MAST::StructuralElementBase*>
326  (_physics_elem)->calculate_stress(false,
327  nullptr,
328  *this);
329 }
330 
331 
332 
333 void
335 
336  // the primal data should have been calculated
337  libmesh_assert(_physics_elem);
339  libmesh_assert(_primal_data_initialized);
340 
342 
343  // ask for the values
344  dynamic_cast<MAST::StructuralElementBase*>
345  (_physics_elem)->calculate_stress(false,
346  &f,
347  *this);
348  }
349 }
350 
351 
352 
353 void
357 
358  // the primal data should have been calculated
359  libmesh_assert(_physics_elem);
360  libmesh_assert(f.is_topology_parameter());
362  libmesh_assert(_primal_data_initialized);
363 
365  &elem = dynamic_cast<const MAST::LevelSetIntersectedElem&>(_physics_elem->elem());
366 
367  // sensitivity only exists at the boundary. So, we proceed with calculation
368  // only if this element has an intersection in the interior, or with a side.
369 
370  if (this->if_evaluate_for_element(elem) &&
373 
374  // ask for the values
375  dynamic_cast<MAST::StructuralElementBase*>
376  (_physics_elem)->calculate_stress_boundary_velocity(f, *this,
378  vel);
379  }
380 }
381 
382 
383 
384 Real
386 
387  libmesh_assert(!_if_stress_plot_mode);
388 
389  // if this has not been initialized, then we should do so now
391  this->functional_for_all_elems();
392 
393  return _sigma_vm_p_norm;
394 }
395 
396 
397 Real
399 
400  libmesh_assert(!_if_stress_plot_mode);
401  libmesh_assert(_primal_data_initialized);
402 
403  Real val = 0.;
404 
406  (p, _physics_elem->elem().get_quadrature_elem().id(), val);
407 
408  return val;
409 }
410 
411 
412 
413 Real
415 
416  libmesh_assert(!_if_stress_plot_mode);
417  libmesh_assert(_primal_data_initialized);
418 
419  Real
420  val = 0.,
421  val_b = 0.;
422 
424  if (p.is_topology_parameter())
426  return val+val_b;
427 }
428 
429 
430 
431 void
433 
434  libmesh_assert(_physics_elem);
435  libmesh_assert(!_if_stress_plot_mode);
436  libmesh_assert(_primal_data_initialized);
437 
439 
440  dq_dX.setZero();
441 
442  dynamic_cast<MAST::StructuralElementBase*>
443  (_physics_elem)->calculate_stress(true,
444  nullptr,
445  *this);
446 
448  (_physics_elem->elem().get_quadrature_elem().id(), dq_dX);
449  }
450 }
451 
452 
453 
454 void
456 set_elem_data(unsigned int dim,
457  const libMesh::Elem& ref_elem,
458  MAST::GeomElem& elem) const {
459 
460  libmesh_assert(!_physics_elem);
461 
462  if (dim == 1) {
463 
464  const MAST::ElementPropertyCard1D& p =
465  dynamic_cast<const MAST::ElementPropertyCard1D&>(_discipline->get_property_card(ref_elem));
466 
467  elem.set_local_y_vector(p.y_vector());
468  }
469 }
470 
471 
472 void
474 
475  libmesh_assert(!_physics_elem);
476  libmesh_assert(_assembly);
477  libmesh_assert(_system);
478 
480  dynamic_cast<const MAST::ElementPropertyCardBase&>(_discipline->get_property_card(elem));
481 
482  _physics_elem =
483  MAST::build_structural_element(*_system, *_assembly, elem, p).release();
484 }
485 
486 
487 
488 void
490 
491  std::map<const libMesh::dof_id_type, std::vector<MAST::StressStrainOutputBase::Data*> >::iterator
492  map_it = _stress_data.begin(),
493  map_end = _stress_data.end();
494 
495  for ( ; map_it != map_end; map_it++) {
496 
497  // iterate over all the data and delete them
498  std::vector<MAST::StressStrainOutputBase::Data*>::iterator
499  it = map_it->second.begin(),
500  end = map_it->second.end();
501 
502  for ( ; it != end; it++)
503  delete *it;
504 
505  map_it->second.clear();
506  }
507 
508  _stress_data.clear();
509 
510 
511  map_it = _boundary_stress_data.begin();
512  map_end = _boundary_stress_data.end();
513 
514  for ( ; map_it != map_end; map_it++) {
515 
516  // iterate over all the data and delete them
517  std::vector<MAST::StressStrainOutputBase::Data*>::iterator
518  it = map_it->second.begin(),
519  end = map_it->second.end();
520 
521  for ( ; it != end; it++)
522  delete *it;
523 
524  map_it->second.clear();
525  }
526 
527  _boundary_stress_data.clear();
528 
529  this->clear_elem();
530 
531  _primal_data_initialized = false;
532  _JxW_val = 0.;
533  _sigma_vm_int = 0.;
534  _sigma_vm_p_norm = 0.;
535  _if_stress_plot_mode = false;
536 }
537 
538 
539 
540 
541 void
543 
544  std::map<const libMesh::dof_id_type, std::vector<MAST::StressStrainOutputBase::Data*> >::iterator
545  map_it = _stress_data.begin(),
546  map_end = _stress_data.end();
547 
548  for ( ; map_it != map_end; map_it++) {
549 
550  // iterate over all the data and delete them
551  std::vector<MAST::StressStrainOutputBase::Data*>::iterator
552  it = map_it->second.begin(),
553  end = map_it->second.end();
554 
555  for ( ; it != end; it++)
556  (*it)->clear_sensitivity_data();
557  }
558 
559 
560  map_it = _boundary_stress_data.begin();
561  map_end = _boundary_stress_data.end();
562 
563  for ( ; map_it != map_end; map_it++) {
564 
565  // iterate over all the data and delete them
566  std::vector<MAST::StressStrainOutputBase::Data*>::iterator
567  it = map_it->second.begin(),
568  end = map_it->second.end();
569 
570  for ( ; it != end; it++)
571  delete *it;
572  }
573  _boundary_stress_data.clear();
574 }
575 
576 
577 
578 
579 
580 
584  const unsigned int qp,
585  const libMesh::Point& quadrature_pt,
586  const libMesh::Point& physical_pt,
587  const RealVectorX& stress,
588  const RealVectorX& strain,
589  Real JxW) {
590 
591  // if the element subset has been provided, then make sure that this
592  // element exists in the subdomain.
593  if (_elem_subset.size())
594  libmesh_assert(_elem_subset.count(&e.get_reference_elem()));
595 
596 
599  strain,
600  quadrature_pt,
601  physical_pt,
602  JxW);
603 
604 
605  // check if the specified element exists in the map. If not, add it
606  std::map<const libMesh::dof_id_type, std::vector<MAST::StressStrainOutputBase::Data*> >::iterator
607  it = _stress_data.find(e.get_quadrature_elem().id());
608 
609  // if the element does not exist in the map, add it to the map.
610  if (it == _stress_data.end())
611  it =
612  _stress_data.insert(std::pair<const libMesh::dof_id_type, std::vector<MAST::StressStrainOutputBase::Data*> >
613  (e.get_quadrature_elem().id(),
614  std::vector<MAST::StressStrainOutputBase::Data*>())).first;
615  else
616  // this assumes that the previous qp data is provided and
617  // therefore, this qp number should be == size of the vector.
618  libmesh_assert_equal_to(qp, it->second.size());
619 
620  it->second.push_back(d);
621 
622  return *d;
623 }
624 
625 
626 
627 
631  const unsigned int s,
632  const unsigned int qp,
633  const libMesh::Point& quadrature_pt,
634  const libMesh::Point& physical_pt,
635  const RealVectorX& stress,
636  const RealVectorX& strain,
637  Real JxW_Vn) {
638 
639  // if the element subset has been provided, then make sure that this
640  // element exists in the subdomain.
641  if (_elem_subset.size())
642  libmesh_assert(_elem_subset.count(&e.get_reference_elem()));
643 
644 
647  strain,
648  quadrature_pt,
649  physical_pt,
650  JxW_Vn);
651 
652 
653  // check if the specified element exists in the map. If not, add it
654  std::map<const libMesh::dof_id_type, std::vector<MAST::StressStrainOutputBase::Data*> >::iterator
655  it = _boundary_stress_data.find(e.get_quadrature_elem().id());
656 
657  // if the element does not exist in the map, add it to the map.
658  if (it == _boundary_stress_data.end())
659  it =
660  _boundary_stress_data.insert
661  (std::pair<const libMesh::dof_id_type, std::vector<MAST::StressStrainOutputBase::Data*> >
662  (e.get_quadrature_elem().id(),
663  std::vector<MAST::StressStrainOutputBase::Data*>())).first;
664  else
665  // this assumes that the previous qp data is provided and
666  // therefore, this qp number should be == size of the vector.
667  libmesh_assert_equal_to(qp, it->second.size());
668 
669  it->second.push_back(d);
670 
671  return *d;
672 }
673 
674 
675 
676 
677 const std::map<const libMesh::dof_id_type, std::vector<MAST::StressStrainOutputBase::Data*> >&
679 
680  return _stress_data;
681 }
682 
683 
684 
685 Real
687 
688  std::map<const libMesh::dof_id_type, std::vector<MAST::StressStrainOutputBase::Data*>>::const_iterator
689  it = _stress_data.begin(),
690  end = _stress_data.end();
691 
692  Real
693  vm = 0.,
694  max_vm = 0.;
695 
696  for ( ; it != end; it++) {
697 
698  std::vector<MAST::StressStrainOutputBase::Data*>::const_iterator
699  s_it = it->second.begin(),
700  s_end = it->second.end();
701 
702  for ( ; s_it != s_end; s_it++) {
703  vm = (*s_it)->von_Mises_stress();
704  max_vm = vm>max_vm?vm:max_vm;
705  }
706  }
707 
708  // now, identify the max stress on all ranks.
709  _system->system().comm().max(max_vm);
710 
711  return max_vm;
712 }
713 
714 
715 
716 unsigned int
719 
720  unsigned int n = 0;
721 
722  std::map<const libMesh::dof_id_type, std::vector<MAST::StressStrainOutputBase::Data*> >::const_iterator
723  it = _stress_data.find(e.get_quadrature_elem().id());
724 
725  if ( it != _stress_data.end())
726  n = (unsigned int)it->second.size();
727 
728  return n;
729 }
730 
731 
732 
733 unsigned int
736 
737  unsigned int n = 0;
738 
739  std::map<const libMesh::dof_id_type, std::vector<MAST::StressStrainOutputBase::Data*> >::const_iterator
740  it = _boundary_stress_data.find(e.get_quadrature_elem().id());
741 
742  if ( it != _boundary_stress_data.end())
743  n = (unsigned int)it->second.size();
744 
745  return n;
746 }
747 
748 
749 
750 const std::vector<MAST::StressStrainOutputBase::Data*>&
753 
754  // check if the specified element exists in the map. If not, add it
755  std::map<const libMesh::dof_id_type, std::vector<MAST::StressStrainOutputBase::Data*> >::const_iterator
756  it = _stress_data.find(e.get_quadrature_elem().id());
757 
758  // make sure that the specified elem exists in the map
759  libmesh_assert(it != _stress_data.end());
760 
761  return it->second;
762 }
763 
764 
765 
769  const unsigned int qp) {
770 
771 
772  // check if the specified element exists in the map. If not, add it
773  std::map<const libMesh::dof_id_type, std::vector<MAST::StressStrainOutputBase::Data*> >::const_iterator
774  it = _stress_data.find(e.get_quadrature_elem().id());
775 
776  // make sure that the specified elem exists in the map
777  libmesh_assert(it != _stress_data.end());
778 
779  libmesh_assert_less(qp, it->second.size());
780 
781  return *(it->second[qp]);
782 }
783 
784 
785 
788 
789  MAST::BoundaryConditionBase *rval = nullptr;
790 
791  // this should only be called if the user has specified evaluation of
792  // stress for this
793  libmesh_assert(this->if_evaluate_for_element(elem));
794 
795 
797  vol_loads = _discipline->volume_loads();
798 
799  std::pair<std::multimap<libMesh::subdomain_id_type, MAST::BoundaryConditionBase*>::const_iterator,
800  std::multimap<libMesh::subdomain_id_type, MAST::BoundaryConditionBase*>::const_iterator> it;
801 
802  it = vol_loads.equal_range(elem.get_reference_elem().subdomain_id());
803 
804  // FIXME: that this assumes that only one temperature boundary condition
805  // is applied for an element.
806  for ( ; it.first != it.second; it.first++)
807  if (it.first->second->type() == MAST::TEMPERATURE) {
808 
809  // make sure that only one thermal load exists for an element
810  libmesh_assert(!rval);
811 
812  rval = it.first->second;
813  }
814 
815  return rval;
816 }
817 
818 
819 
820 
821 void
823 
824  libmesh_assert(!_if_stress_plot_mode);
825  libmesh_assert(!_primal_data_initialized);
826  libmesh_assert_greater(_sigma0, 0.);
827 
828  Real
829  sp = 0.,
830  exp_sp = 0.,
831  e_val = 0.,
832  JxW = 0.;
833 
834  _JxW_val = 0.;
835  _sigma_vm_int = 0.;
836  _sigma_vm_p_norm = 0.;
837 
838  // first find the data with the maximum value, to be used for scaling
839 
840  // iterate over all element data
841  std::map<const libMesh::dof_id_type, std::vector<MAST::StressStrainOutputBase::Data*> >::const_iterator
842  map_it = _stress_data.begin(),
843  map_end = _stress_data.end();
844 
845  for ( ; map_it != map_end; map_it++) {
846 
847  std::vector<MAST::StressStrainOutputBase::Data*>::const_iterator
848  vec_it = map_it->second.begin(),
849  vec_end = map_it->second.end();
850 
851  for ( ; vec_it != vec_end; vec_it++) {
852 
853  // ask this data point for the von Mises stress value
854  e_val = (*vec_it)->von_Mises_stress();
855  JxW = (*vec_it)->quadrature_point_JxW();
856 
857  // we do not use absolute value here, since von Mises stress
858  // is >= 0.
859  sp = pow((e_val-_sigma0)/_sigma0, _p_norm_weight);
860  if (_rho * sp > _exp_arg_lim)
861  exp_sp = exp(_exp_arg_lim);
862  else
863  exp_sp = exp(_rho * sp);
864  _sigma_vm_int += pow(e_val/_sigma0, _p_norm_stress) * exp_sp * JxW;
865  _JxW_val += exp_sp * JxW;
866  }
867  }
868 
869  // sum over all processors, since part of the mesh will exist on the
870  // other processors.
871  _system->system().comm().sum(_sigma_vm_int);
872  _system->system().comm().sum(_JxW_val);
873 
876 }
877 
878 
879 
880 void
883  Real& dsigma_vm_val_df) const {
884 
885  libmesh_assert(!_if_stress_plot_mode);
886  libmesh_assert(_primal_data_initialized);
887 
888  Real
889  val = 0.;
890 
891  dsigma_vm_val_df = 0.;
892 
893  // iterate over all element data
894  std::map<const libMesh::dof_id_type, std::vector<MAST::StressStrainOutputBase::Data*> >::const_iterator
895  map_it = _stress_data.begin(),
896  map_end = _stress_data.end();
897 
898  for ( ; map_it != map_end; map_it++) {
899 
900  this->functional_sensitivity_for_elem(f, map_it->first, val);
901  dsigma_vm_val_df += val;
902  }
903 
904  // sum over all processors, since part of the mesh will exist on the
905  // other processors
906  _system->system().comm().sum(dsigma_vm_val_df);
907 }
908 
909 
910 
911 
912 void
915  Real& dsigma_vm_val_df) const {
916 
917  libmesh_assert(!_if_stress_plot_mode);
918  libmesh_assert(_primal_data_initialized);
919 
920  Real
921  val = 0.;
922 
923  dsigma_vm_val_df = 0.;
924 
925  // iterate over all element data
926  std::map<const libMesh::dof_id_type, std::vector<MAST::StressStrainOutputBase::Data*> >::const_iterator
927  map_it = _boundary_stress_data.begin(),
928  map_end = _boundary_stress_data.end();
929 
930  for ( ; map_it != map_end; map_it++) {
931 
932  this->functional_boundary_sensitivity_for_elem(f, map_it->first, val);
933  dsigma_vm_val_df += val;
934  }
935 
936  // sum over all processors, since part of the mesh will exist on the
937  // other processors
938  _system->system().comm().sum(dsigma_vm_val_df);
939 }
940 
941 
942 
943 void
946  const libMesh::dof_id_type e_id,
947  Real& dsigma_vm_val_df) const {
948 
949  libmesh_assert(!_if_stress_plot_mode);
950  libmesh_assert(_primal_data_initialized);
951  libmesh_assert_greater(_sigma0, 0.);
952 
953  Real
954  sp_stress = 0.,
955  sp1_stress = 0.,
956  sp_weight = 0.,
957  sp1_weight = 0.,
958  exp_sp = 0.,
959  e_val = 0.,
960  de_val = 0.,
961  num_sens = 0.,
962  denom_sens = 0.,
963  JxW = 0.;
964 
965  dsigma_vm_val_df = 0.;
966 
967  // iterate over all element data
968  std::map<const libMesh::dof_id_type, std::vector<MAST::StressStrainOutputBase::Data*> >::const_iterator
969  map_it = _stress_data.find(e_id),
970  map_end = _stress_data.end();
971 
972  libmesh_assert(map_it != map_end);
973 
974  std::vector<MAST::StressStrainOutputBase::Data*>::const_iterator
975  vec_it = map_it->second.begin(),
976  vec_end = map_it->second.end();
977 
978  for ( ; vec_it != vec_end; vec_it++) {
979 
980  // ask this data point for the von Mises stress value
981  e_val = (*vec_it)->von_Mises_stress();
982  de_val = (*vec_it)->dvon_Mises_stress_dp(f);
983  JxW = (*vec_it)->quadrature_point_JxW();
984 
985  // we do not use absolute value here, since von Mises stress
986  // is >= 0.
987  sp_stress = pow(e_val/_sigma0, _p_norm_stress);
988  sp1_stress = pow(e_val/_sigma0, _p_norm_stress-1.);
989  sp_weight = pow((e_val-_sigma0)/_sigma0, _p_norm_weight);
990  sp1_weight = pow((e_val-_sigma0)/_sigma0, _p_norm_weight-1.);
991  if (_rho*sp_weight > _exp_arg_lim) {
992  exp_sp = exp(_exp_arg_lim);
993  denom_sens += 0.;
994  num_sens += _p_norm_stress * sp1_stress * exp_sp * de_val/_sigma0 * JxW;
995  }
996  else {
997  exp_sp = exp(_rho * sp_weight);
998  denom_sens += exp_sp * _rho * _p_norm_weight * sp1_weight * de_val/_sigma0 * JxW;
999  num_sens += (de_val/_sigma0 * JxW * exp_sp *
1000  (_p_norm_stress * sp1_stress +
1001  sp_stress * _rho * _p_norm_weight * sp1_weight));
1002  }
1003 
1004  }
1005 
1006  dsigma_vm_val_df = _sigma0/_p_norm_stress * pow(_sigma_vm_int/_JxW_val, 1./_p_norm_stress - 1.) *
1007  (num_sens / _JxW_val - _sigma_vm_int / pow(_JxW_val, 2.) * denom_sens);
1008 }
1009 
1010 
1011 
1012 
1013 
1014 void
1017  const libMesh::dof_id_type e_id,
1018  Real& dsigma_vm_val_df) const {
1019 
1020  libmesh_assert(!_if_stress_plot_mode);
1021  libmesh_assert(_primal_data_initialized);
1022  libmesh_assert_greater(_sigma0, 0.);
1023 
1024  Real
1025  sp = 0.,
1026  exp_sp = 0.,
1027  e_val = 0.,
1028  JxW_Vn = 0.,
1029  num_sens = 0.,
1030  denom_sens = 0.;
1031 
1032  dsigma_vm_val_df = 0.;
1033 
1034  // iterate over all element data
1035  std::map<const libMesh::dof_id_type, std::vector<MAST::StressStrainOutputBase::Data*> >::const_iterator
1036  map_it = _boundary_stress_data.find(e_id),
1037  map_end = _boundary_stress_data.end();
1038 
1039  libmesh_assert(map_it != map_end);
1040 
1041  std::vector<MAST::StressStrainOutputBase::Data*>::const_iterator
1042  vec_it = map_it->second.begin(),
1043  vec_end = map_it->second.end();
1044 
1045  for ( ; vec_it != vec_end; vec_it++) {
1046 
1047  // ask this data point for the von Mises stress value
1048  e_val = (*vec_it)->von_Mises_stress();
1049  JxW_Vn = (*vec_it)->quadrature_point_JxW();
1050 
1051  // we do not use absolute value here, since von Mises stress
1052  // is >= 0.
1053  sp = pow((e_val-_sigma0)/_sigma0, _p_norm_weight);
1054  if (_rho*sp > _exp_arg_lim)
1055  exp_sp = exp(_exp_arg_lim);
1056  else
1057  exp_sp = exp(_rho * sp);
1058  denom_sens += exp_sp * JxW_Vn;
1059  num_sens += pow(e_val/_sigma0, _p_norm_stress) * exp_sp * JxW_Vn;
1060  }
1061 
1062  dsigma_vm_val_df = _sigma0/_p_norm_stress * pow(_sigma_vm_int/_JxW_val, 1./_p_norm_stress - 1.) *
1063  (num_sens / _JxW_val - _sigma_vm_int / pow(_JxW_val, 2.) * denom_sens);
1064 }
1065 
1066 
1067 
1068 
1069 void
1071  RealVectorX& dq_dX) const {
1072  libmesh_assert(!_if_stress_plot_mode);
1073  libmesh_assert(_primal_data_initialized);
1074  libmesh_assert_greater(_sigma0, 0.);
1075 
1076  Real
1077  sp_stress = 0.,
1078  sp1_stress = 0.,
1079  sp_weight = 0.,
1080  sp1_weight = 0.,
1081  exp_sp = 0.,
1082  e_val = 0.,
1083  JxW = 0.;
1084 
1085  RealVectorX
1086  denom_sens = RealVectorX::Zero(dq_dX.size()),
1087  num_sens = RealVectorX::Zero(dq_dX.size()),
1088  de_val = RealVectorX::Zero(dq_dX.size());
1089  dq_dX.setZero();
1090 
1091  // first find the data with the maximum value, to be used for scaling
1092 
1093  // iterate over all element data
1094  std::map<const libMesh::dof_id_type, std::vector<MAST::StressStrainOutputBase::Data*> >::const_iterator
1095  map_it = _stress_data.find(e_id),
1096  map_end = _stress_data.end();
1097 
1098  // make sure that the data exists
1099  libmesh_assert(map_it != map_end);
1100 
1101  std::vector<MAST::StressStrainOutputBase::Data*>::const_iterator
1102  vec_it = map_it->second.begin(),
1103  vec_end = map_it->second.end();
1104 
1105 
1106  for ( ; vec_it != vec_end; vec_it++) {
1107 
1108  // ask this data point for the von Mises stress value
1109  e_val = (*vec_it)->von_Mises_stress();
1110  de_val = (*vec_it)->dvon_Mises_stress_dX();
1111  JxW = (*vec_it)->quadrature_point_JxW();
1112 
1113  // we do not use absolute value here, since von Mises stress
1114  // is >= 0.
1115  sp_stress = pow(e_val/_sigma0, _p_norm_stress);
1116  sp1_stress = pow(e_val/_sigma0, _p_norm_stress-1.);
1117  sp_weight = pow((e_val-_sigma0)/_sigma0, _p_norm_weight);
1118  sp1_weight = pow((e_val-_sigma0)/_sigma0, _p_norm_weight-1.);
1119  if (_rho*sp_weight > _exp_arg_lim) {
1120 
1121  exp_sp = exp(_exp_arg_lim);
1122  denom_sens += 0. * de_val;
1123  num_sens += 1. * _p_norm_stress * sp1_stress * exp_sp/_sigma0 * JxW * de_val;
1124  }
1125  else {
1126 
1127  exp_sp = exp(_rho * sp_weight);
1128  denom_sens += exp_sp * _rho * _p_norm_weight * sp1_weight/_sigma0 * JxW * de_val;
1129  num_sens += (de_val/_sigma0 * JxW * exp_sp *
1130  (_p_norm_stress * sp1_stress +
1131  sp_stress * _rho * _p_norm_weight * sp1_weight));
1132  }
1133  }
1134 
1135  dq_dX = _sigma0/_p_norm_stress * pow(_sigma_vm_int/_JxW_val, 1./_p_norm_stress - 1.) *
1136  (num_sens / _JxW_val - _sigma_vm_int / pow(_JxW_val, 2.) * denom_sens);
1137 }
1138 
1139 
1140 
1141 
virtual Real output_sensitivity_for_elem(const MAST::FunctionBase &p)
Data(const RealVectorX &stress, const RealVectorX &strain, const libMesh::Point &qp, const libMesh::Point &xyz, Real JxW)
const MAST::VolumeBCMapType & volume_loads() const
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...
std::map< const MAST::FunctionBase *, RealVectorX > _stress_sensitivity
map of sensitivity of the stress with respect to a parameter
std::multimap< libMesh::subdomain_id_type, MAST::BoundaryConditionBase * > VolumeBCMapType
libMesh::Point _qp
quadrature point location in element coordinates
RealMatrixX _dstress_dX
derivative of stress wrt state vector
virtual void functional_for_all_elems()
calculates and returns the von Mises p-norm functional for all the elements that this object currentl...
const RealVectorX & stress() const
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...
std::unique_ptr< MAST::StructuralElementBase > build_structural_element(MAST::SystemInitialization &sys, MAST::AssemblyBase &assembly, const MAST::GeomElem &elem, const MAST::ElementPropertyCardBase &p)
builds the structural element for the specified element type
unsigned int n_boundary_stress_strain_data_for_elem(const GeomElem &e) const
void clear()
clears the data structure of any stored values so that it can be used for another element...
virtual void functional_sensitivity_for_all_elems(const MAST::FunctionBase &f, Real &dsigma_vm_val_df) const
calculates and returns the sensitivity of von Mises p-norm functional for all the elements that this ...
std::map< const libMesh::dof_id_type, std::vector< MAST::StressStrainOutputBase::Data * > > _boundary_stress_data
vector of stress with the associated location details
RealVectorX & y_vector()
returns value of the property val.
const MAST::GeomElem & elem() const
Definition: elem_base.h:117
const RealMatrixX & get_dstrain_dX() const
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
Real _sigma0
reference stress value used in scaling volume.
virtual void zero_for_sensitivity()
zeroes the output quantity values stored inside this object so that assembly process can begin...
libMesh::Real Real
MAST::PhysicsDisciplineBase * _discipline
Real _rho
exponent used in scaling volume based on stress value.
const RealVectorX & get_strain_sensitivity(const MAST::FunctionBase &f) const
@ returns the sensitivity of the data with respect to a function
virtual const std::vector< MAST::StressStrainOutputBase::Data * > & get_stress_strain_data_for_elem(const MAST::GeomElem &e) const
bool _primal_data_initialized
primal data, needed for sensitivity and adjoints
virtual void evaluate_topology_sensitivity(const MAST::FunctionBase &f, const MAST::FieldFunction< RealVectorX > &vel)
This evaluates the contribution to the topology sensitivity on the boundary.
Real _p_norm_stress
norm to be used for calculation of output stress function.
virtual const libMesh::Elem & get_reference_elem() const
Definition: geom_elem.cpp:54
const libMesh::Point & point_location_in_element_coordinate() const
virtual Real output_sensitivity_total(const MAST::FunctionBase &p)
Real dvon_Mises_stress_dp(const MAST::FunctionBase &f) const
unsigned int n_stress_strain_data_for_elem(const MAST::GeomElem &e) const
std::map< const libMesh::dof_id_type, std::vector< MAST::StressStrainOutputBase::Data * > > _stress_data
vector of stress with the associated location details
RealMatrixX _dstrain_dX
derivative of strain data wrt state vector
virtual const std::map< const libMesh::dof_id_type, std::vector< MAST::StressStrainOutputBase::Data * > > & get_stress_strain_data() const
Matrix< Real, Dynamic, Dynamic > RealMatrixX
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()...
MAST::BoundaryConditionBase * get_thermal_load_for_elem(const MAST::GeomElem &elem)
virtual void functional_boundary_sensitivity_for_all_elems(const MAST::FunctionBase &f, Real &dsigma_vm_val_df) const
calculates and returns the sensitivity of von Mises p-norm functional for all the elements that this ...
This class inherits from MAST::GeomElem and provides an interface to initialize FE objects on sub-ele...
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
virtual void zero_for_analysis()
zeroes the output quantity values stored inside this object so that assembly process can begin...
const MAST::ElementPropertyCardBase & get_property_card(const libMesh::Elem &elem) const
get property card for the specified element
virtual void clear_sensitivity_data()
clears the data stored for sensitivity analysis.
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...
This class provides a mechanism to store stress/strain values, their derivatives and sensitivity valu...
const RealVectorX & strain() const
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 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...
StressStrainOutputBase()
default constructor
Real _JxW
quadrature point JxW (product of transformation Jacobian and quadrature weight) for use in definition...
void set_derivatives(const RealMatrixX &dstress_dX, const RealMatrixX &dstrain_dX)
adds the derivative data
virtual const libMesh::Elem & get_quadrature_elem() const
Definition: geom_elem.cpp:72
const RealVectorX & get_stress_sensitivity(const MAST::FunctionBase &f) const
@ returns the sensitivity of the data with respect to a function
virtual void evaluate()
this evaluates all relevant stress components on the element to evaluate the p-averaged quantity...
std::map< const MAST::FunctionBase *, RealVectorX > _strain_sensitivity
map of sensitivity of the strain with respect to a parameter
virtual MAST::StressStrainOutputBase::Data & get_stress_strain_data_for_elem_at_qp(const MAST::GeomElem &e, const unsigned int qp)
virtual void clear_elem()
clears the element initialization
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 has_stress_sensitivity(const MAST::FunctionBase &f) const
@ returns true if sensitivity data is available for function f .
bool _if_stress_plot_mode
identifies the mode in which evaluation is peformed.
virtual MAST::StressStrainOutputBase::Data & add_stress_strain_at_qp_location(const MAST::GeomElem &e, const unsigned int qp, const libMesh::Point &quadrature_pt, const libMesh::Point &physical_pt, const RealVectorX &stress, const RealVectorX &strain, Real JxW)
add the stress tensor associated with the qp.
const RealMatrixX & get_dstress_dX() const
std::set< const libMesh::Elem * > _elem_subset
set of elements for which the data will be stored.
virtual bool is_topology_parameter() const
Definition: function_base.h:97
void set_sensitivity(const MAST::FunctionBase &f, const RealVectorX &dstress_df, const RealVectorX &dstrain_df)
sets the sensitivity of the data with respect to a function
virtual MAST::StressStrainOutputBase::Data & add_stress_strain_at_boundary_qp_location(const MAST::GeomElem &e, const unsigned int s, const unsigned int qp, const libMesh::Point &quadrature_pt, const libMesh::Point &physical_pt, const RealVectorX &stress, const RealVectorX &strain, Real JxW_Vn)
add the stress tensor associated with the qp on side s of element e.
MAST::SystemInitialization * _system