MAST
stress_output_base.h
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 #ifndef __mast__stress_output_base__
21 #define __mast__stress_output_base__
22 
23 // C++ includes
24 #include <map>
25 #include <vector>
26 
27 // MAST includes
28 #include "base/mast_data_types.h"
31 
32 
33 // libMesh includes
34 #include "libmesh/elem.h"
35 
36 namespace MAST {
37 
38 
39  // Forward declerations
40  class FunctionBase;
41 
42 
63 
64  public:
65 
66 
72  class Data {
73 
74  public:
75  Data(const RealVectorX& stress,
76  const RealVectorX& strain,
77  const libMesh::Point& qp,
78  const libMesh::Point& xyz,
79  Real JxW);
80 
81 
83 
84 
89  const libMesh::Point&
91 
95  const RealVectorX& stress() const;
96 
97 
101  const RealVectorX& strain() const;
102 
103 
107  Real quadrature_point_JxW() const;
108 
112  Real von_Mises_stress() const;
113 
114 
119 
120 
125 
126 
130  void set_derivatives(const RealMatrixX& dstress_dX,
131  const RealMatrixX& dstrain_dX);
132 
133 
137  const RealMatrixX& get_dstress_dX() const;
138 
139 
143  const RealMatrixX& get_dstrain_dX() const;
144 
145 
149  void set_sensitivity(const MAST::FunctionBase& f,
150  const RealVectorX& dstress_df,
151  const RealVectorX& dstrain_df);
152 
153 
158  bool
160 
165  const RealVectorX&
167 
168 
173  const RealVectorX&
175 
176 
177  protected:
178 
183 
184 
189 
190 
194  std::map<const MAST::FunctionBase*, RealVectorX> _stress_sensitivity;
195 
196 
200  std::map<const MAST::FunctionBase*, RealVectorX> _strain_sensitivity;
201 
202 
207 
212 
213 
217  libMesh::Point _qp;
218 
219 
223  libMesh::Point _xyz;
224 
225 
231  };
232 
233 
234 
239 
240  virtual ~StressStrainOutputBase();
241 
242 
246  void set_aggregation_coefficients(Real p1, Real p2, Real rho, Real sigma0) {
247 
248  _p_norm_stress = p1;
249  _p_norm_weight = p2;
250  _rho = rho;
251  _sigma0 = sigma0;
252  }
253 
258 
259  return _p_norm_stress;
260  }
261 
262 
267  void set_stress_plot_mode(bool f) {
268 
270  }
271 
272 
276  virtual void
277  set_elem_data(unsigned int dim,
278  const libMesh::Elem& ref_elem,
279  MAST::GeomElem& elem) const;
280 
284  virtual void init(const MAST::GeomElem& elem);
285 
294  virtual void zero_for_analysis();
295 
296 
302  virtual void zero_for_sensitivity();
303 
310  virtual void evaluate();
311 
318  virtual void evaluate_sensitivity(const MAST::FunctionBase& f);
319 
327 
328  libmesh_assert(false); // to be implemented
329  }
330 
342  virtual void
345 
349  virtual Real output_for_elem() {
350  //
351  libmesh_error();
352  }
353 
357  virtual Real output_total();
358 
365 
374 
375 
382  virtual void output_derivative_for_elem(RealVectorX& dq_dX);
383 
384 
385 
386  bool stress_plot_mode() const {
387  return _if_stress_plot_mode;
388  }
389 
390  bool primal_data_initialized() const {
392  }
393 
398  void clear();
399 
403  virtual void clear_sensitivity_data();
404 
405 
413 
414 
419  unsigned int
421 
422 
427  unsigned int
429 
430 
431 
438  const unsigned int qp,
439  const libMesh::Point& quadrature_pt,
440  const libMesh::Point& physical_pt,
441  const RealVectorX& stress,
442  const RealVectorX& strain,
443  Real JxW);
444 
445 
452  const unsigned int s,
453  const unsigned int qp,
454  const libMesh::Point& quadrature_pt,
455  const libMesh::Point& physical_pt,
456  const RealVectorX& stress,
457  const RealVectorX& strain,
458  Real JxW_Vn);
459 
460 
461 
465  virtual const std::map<const libMesh::dof_id_type,
466  std::vector<MAST::StressStrainOutputBase::Data*> >&
467  get_stress_strain_data() const;
468 
469 
474 
478  virtual const std::vector<MAST::StressStrainOutputBase::Data*>&
480 
481 
488  const unsigned int qp);
489 
490 
491 
499  virtual void functional_for_all_elems();
500 
501 
513  (const MAST::FunctionBase& f,
514  Real& dsigma_vm_val_df) const;
515 
516 
530  (const MAST::FunctionBase& f,
531  Real& dsigma_vm_val_df) const;
532 
533 
539  (const MAST::FunctionBase& f,
540  const libMesh::dof_id_type e_id,
541  Real& dsigma_vm_val_df) const;
542 
543 
549  (const MAST::FunctionBase& f,
550  const libMesh::dof_id_type e_id,
551  Real& dsigma_vm_val_df) const;
552 
553 
554 
567  (const libMesh::dof_id_type e_id,
568  RealVectorX& dq_dX) const;
569 
570 
571 
572  protected:
573 
579 
584 
589 
591 
597 
605 
609  std::map<const libMesh::dof_id_type, std::vector<MAST::StressStrainOutputBase::Data*>>
611 
612 
616  std::map<const libMesh::dof_id_type, std::vector<MAST::StressStrainOutputBase::Data*>>
618  };
619 }
620 
621 #endif // __mast__stress_output_base__
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)
virtual Real output_for_elem()
should not get called for this output.
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
libMesh::Point _qp
quadrature point location in element coordinates
Data structure provides the mechanism to store stress and strain output from a structural analysis...
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
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
const RealMatrixX & get_dstrain_dX() const
This provides the base class for definitin of element level contribution of output quantity in an ana...
Real _sigma0
reference stress value used in scaling volume.
void set_stress_plot_mode(bool f)
tells the object that the calculation is for stress to be output for plotting.
virtual void zero_for_sensitivity()
zeroes the output quantity values stored inside this object so that assembly process can begin...
libMesh::Real Real
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
libMesh::Point _xyz
quadrature point location in physical coordinates
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.
void set_aggregation_coefficients(Real p1, Real p2, Real rho, Real sigma0)
sets the norm for calculation of stress functional
const libMesh::Point & point_location_in_element_coordinate() const
virtual void evaluate_shape_sensitivity(const MAST::FunctionBase &f)
this evaluates all relevant shape sensitivity components on the element.
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 ...
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...
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
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 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
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.