MAST
MAST::StressStrainOutputBase Class Reference

Data structure provides the mechanism to store stress and strain output from a structural analysis. More...

#include <stress_output_base.h>

Inheritance diagram for MAST::StressStrainOutputBase:
Collaboration diagram for MAST::StressStrainOutputBase:

Classes

class  Data
 This class provides a mechanism to store stress/strain values, their derivatives and sensitivity values corresponding to a specific quadrature point on the element. More...
 

Public Member Functions

 StressStrainOutputBase ()
 default constructor More...
 
virtual ~StressStrainOutputBase ()
 
void set_aggregation_coefficients (Real p1, Real p2, Real rho, Real sigma0)
 sets the $p-$norm for calculation of stress functional More...
 
Real get_p_stress_val ()
 
void set_stress_plot_mode (bool f)
 tells the object that the calculation is for stress to be output for plotting. More...
 
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. More...
 
virtual void init (const MAST::GeomElem &elem)
 initialize for the element. More...
 
virtual void zero_for_analysis ()
 zeroes the output quantity values stored inside this object so that assembly process can begin. More...
 
virtual void zero_for_sensitivity ()
 zeroes the output quantity values stored inside this object so that assembly process can begin. More...
 
virtual void evaluate ()
 this evaluates all relevant stress components on the element to evaluate the p-averaged quantity. More...
 
virtual void evaluate_sensitivity (const MAST::FunctionBase &f)
 this evaluates all relevant stress sensitivity components on the element to evaluate the p-averaged quantity sensitivity. More...
 
virtual void evaluate_shape_sensitivity (const MAST::FunctionBase &f)
 this evaluates all relevant shape sensitivity components on the element. More...
 
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. More...
 
virtual Real output_for_elem ()
 should not get called for this output. More...
 
virtual Real output_total ()
 
virtual Real output_sensitivity_for_elem (const MAST::FunctionBase &p)
 
virtual Real output_sensitivity_total (const MAST::FunctionBase &p)
 
virtual void output_derivative_for_elem (RealVectorX &dq_dX)
 calculates the derivative of p-norm von Mises stress for the $p-$norm identified using set_p_val(). More...
 
bool stress_plot_mode () const
 
bool primal_data_initialized () const
 
void clear ()
 clears the data structure of any stored values so that it can be used for another element. More...
 
virtual void clear_sensitivity_data ()
 clears the data stored for sensitivity analysis. More...
 
MAST::BoundaryConditionBaseget_thermal_load_for_elem (const MAST::GeomElem &elem)
 
unsigned int n_stress_strain_data_for_elem (const MAST::GeomElem &e) const
 
unsigned int n_boundary_stress_strain_data_for_elem (const GeomElem &e) const
 
virtual MAST::StressStrainOutputBase::Dataadd_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. More...
 
virtual MAST::StressStrainOutputBase::Dataadd_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. More...
 
virtual const std::map< const libMesh::dof_id_type, std::vector< MAST::StressStrainOutputBase::Data * > > & get_stress_strain_data () const
 
Real get_maximum_von_mises_stress () const
 
virtual const std::vector< MAST::StressStrainOutputBase::Data * > & get_stress_strain_data_for_elem (const MAST::GeomElem &e) const
 
virtual MAST::StressStrainOutputBase::Dataget_stress_strain_data_for_elem_at_qp (const MAST::GeomElem &e, const unsigned int qp)
 
virtual void functional_for_all_elems ()
 calculates and returns the von Mises p-norm functional for all the elements that this object currently stores data for. More...
 
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 object currently stores data for. More...
 
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 object currently stores data for. More...
 
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. More...
 
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. More...
 
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 specified element. More...
 
- Public Member Functions inherited from MAST::OutputAssemblyElemOperations
 OutputAssemblyElemOperations ()
 
virtual ~OutputAssemblyElemOperations ()
 virtual destructor More...
 
void set_participating_subdomains (const std::set< libMesh::subdomain_id_type > &sids)
 The output function can be a boundary integrated quantity, volume integrated quantity or a combination of these two. More...
 
void set_participating_elements (const std::set< const libMesh::Elem * > &elems)
 sets the elements for which this object will evaluate and store the output data. More...
 
void set_participating_elements_to_all ()
 This will allow volume contribution from all elements. More...
 
void set_participating_boundaries (const std::set< libMesh::boundary_id_type > &bids)
 The assembly will integration over boudnaries with ids specified in bids. More...
 
const std::set< const libMesh::Elem * > & get_participating_elements () const
 
const std::set< libMesh::subdomain_id_type > & get_participating_subdomains ()
 
const std::set< libMesh::boundary_id_type > & get_participating_boundaries ()
 
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 is in the subset. More...
 
virtual bool if_evaluate_for_boundary (const MAST::GeomElem &elem, const unsigned int s) const
 checks to see if the specified side of the element needs evaluation of the output contribution. More...
 
- Public Member Functions inherited from MAST::AssemblyElemOperations
 AssemblyElemOperations ()
 
virtual ~AssemblyElemOperations ()
 
MAST::SystemInitializationget_system_initialization ()
 
MAST::PhysicsDisciplineBaseget_discipline ()
 
virtual void set_discipline_and_system (MAST::PhysicsDisciplineBase &discipline, MAST::SystemInitialization &system)
 attaches a system to this discipline More...
 
virtual void clear_discipline_and_system ()
 clears association with a system to this discipline More...
 
virtual void set_assembly (MAST::AssemblyBase &assembly)
 sets the assembly object More...
 
virtual MAST::AssemblyBaseget_assembly ()
 
virtual void clear_assembly ()
 clears the assembly object More...
 
virtual void clear_elem ()
 clears the element initialization More...
 
MAST::ElementBaseget_physics_elem ()
 
virtual void set_elem_solution (const RealVectorX &sol)
 sets the element solution More...
 
virtual void set_elem_solution_sensitivity (const RealVectorX &sol)
 sets the element solution sensitivity More...
 
virtual void set_elem_perturbed_solution (const RealVectorX &sol)
 sets the element perturbed solution More...
 
virtual void set_elem_velocity (const RealVectorX &vel)
 sets the element velocity More...
 
virtual void set_elem_velocity_sensitivity (const RealVectorX &vel)
 sets the element velocity sensitivity More...
 
virtual void set_elem_perturbed_velocity (const RealVectorX &vel)
 sets the element perturbed velocity More...
 
virtual void set_elem_acceleration (const RealVectorX &accel)
 sets the element acceleration More...
 
virtual void set_elem_acceleration_sensitivity (const RealVectorX &accel)
 sets the element acceleration More...
 
virtual void set_elem_perturbed_acceleration (const RealVectorX &accel)
 sets the element perturbed acceleration More...
 

Protected Attributes

Real _p_norm_stress
 $ p-$norm to be used for calculation of output stress function. More...
 
Real _p_norm_weight
 
Real _rho
 exponent used in scaling volume based on stress value. More...
 
Real _sigma0
 reference stress value used in scaling volume. More...
 
Real _exp_arg_lim
 
bool _primal_data_initialized
 primal data, needed for sensitivity and adjoints More...
 
Real _JxW_val
 
Real _sigma_vm_int
 
Real _sigma_vm_p_norm
 
bool _if_stress_plot_mode
 identifies the mode in which evaluation is peformed. More...
 
std::map< const libMesh::dof_id_type, std::vector< MAST::StressStrainOutputBase::Data * > > _stress_data
 vector of stress with the associated location details More...
 
std::map< const libMesh::dof_id_type, std::vector< MAST::StressStrainOutputBase::Data * > > _boundary_stress_data
 vector of stress with the associated location details More...
 
- Protected Attributes inherited from MAST::OutputAssemblyElemOperations
bool _if_evaluate_on_all_elems
 if true, evaluates on all elements. More...
 
std::set< const libMesh::Elem * > _elem_subset
 set of elements for which the data will be stored. More...
 
std::set< libMesh::subdomain_id_type > _sub_domain_ids
 set of subdomain ids for which data will be processed. More...
 
std::set< libMesh::boundary_id_type > _bids
 set of bids for which data will be processed More...
 
- Protected Attributes inherited from MAST::AssemblyElemOperations
MAST::SystemInitialization_system
 
MAST::PhysicsDisciplineBase_discipline
 
MAST::AssemblyBase_assembly
 
MAST::ElementBase_physics_elem
 

Detailed Description

Data structure provides the mechanism to store stress and strain output from a structural analysis.

The user may specify one of three evaluation modes: centroid, element quadrature points, or user specified element local points.

This is, by default, a 3D tensor that should be initialized with 6 components. The first three components of stress are the direct stresses in the $ x, y, z$ directions, $ \sigma_{xx}, \sigma_{yy}, \sigma_{zz}$, and the next three components are the shear stresses $ \tau_{xy}, \tau_{yz}, \tau_{xz} $.

The first three components of strain are the direct strains in the $x, y, z$ directions, $ \epsilon_{xx}, \epsilon_{yy}, \epsilon_{zz} $, and the next three components are the engineering shear strains $ \gamma_{xy}, \gamma_{yz}, \gamma_{xz} $.

Definition at line 61 of file stress_output_base.h.

Constructor & Destructor Documentation

MAST::StressStrainOutputBase::StressStrainOutputBase ( )

default constructor

Definition at line 273 of file stress_output_base.cpp.

MAST::StressStrainOutputBase::~StressStrainOutputBase ( )
virtual

Definition at line 291 of file stress_output_base.cpp.

Here is the call graph for this function:

Member Function Documentation

MAST::StressStrainOutputBase::Data & MAST::StressStrainOutputBase::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 
)
virtual

add the stress tensor associated with the qp on side s of element e.

Returns
a reference to Data.

Reimplemented in MAST::StressTemperatureAdjoint.

Definition at line 630 of file stress_output_base.cpp.

Here is the call graph for this function:

Here is the caller graph for this function:

MAST::StressStrainOutputBase::Data & MAST::StressStrainOutputBase::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 
)
virtual

add the stress tensor associated with the qp.

Returns
a reference to Data.

Reimplemented in MAST::StressTemperatureAdjoint.

Definition at line 583 of file stress_output_base.cpp.

Here is the call graph for this function:

Here is the caller graph for this function:

void MAST::StressStrainOutputBase::clear ( )

clears the data structure of any stored values so that it can be used for another element.

Definition at line 489 of file stress_output_base.cpp.

Here is the call graph for this function:

Here is the caller graph for this function:

void MAST::StressStrainOutputBase::clear_sensitivity_data ( )
virtual

clears the data stored for sensitivity analysis.

Definition at line 542 of file stress_output_base.cpp.

Here is the call graph for this function:

void MAST::StressStrainOutputBase::evaluate ( )
virtual

this evaluates all relevant stress components on the element to evaluate the p-averaged quantity.

This is only done on the current element for which this object has been initialized.

Implements MAST::OutputAssemblyElemOperations.

Definition at line 317 of file stress_output_base.cpp.

Here is the call graph for this function:

Here is the caller graph for this function:

void MAST::StressStrainOutputBase::evaluate_sensitivity ( const MAST::FunctionBase f)
virtual

this evaluates all relevant stress sensitivity components on the element to evaluate the p-averaged quantity sensitivity.

This is only done on the current element for which this object has been initialized.

Implements MAST::OutputAssemblyElemOperations.

Definition at line 334 of file stress_output_base.cpp.

Here is the call graph for this function:

Here is the caller graph for this function:

virtual void MAST::StressStrainOutputBase::evaluate_shape_sensitivity ( const MAST::FunctionBase f)
inlinevirtual

this evaluates all relevant shape sensitivity components on the element.

This is only done on the current element for which this object has been initialized.

Implements MAST::OutputAssemblyElemOperations.

Definition at line 326 of file stress_output_base.h.

Here is the call graph for this function:

void MAST::StressStrainOutputBase::evaluate_topology_sensitivity ( const MAST::FunctionBase f,
const MAST::FieldFunction< RealVectorX > &  vel 
)
virtual

This evaluates the contribution to the topology sensitivity on the boundary.

Given that the integral is nonlinear due to the $p-$norm, the expression is quite involved:

\[ \frac{ \frac{1}{p} \left( \int_\Omega (\sigma_{VM}(\Omega))^p ~ d\Omega \right)^{\frac{1}{p}-1}}{\left( \int_\Omega ~ d\Omega \right)^{\frac{1}{p}}} \int_\Gamma V_n \sigma_{VM}^p ~d\Gamma + \frac{ \frac{-1}{p} \left( \int_\Omega (\sigma_{VM}(\Omega))^p ~ d\Omega \right)^{\frac{1}{p}}}{\left( \int_\Omega ~ d\Omega \right)^{\frac{1+p}{p}}} \int_\Gamma V_n ~d\Gamma \]

Implements MAST::OutputAssemblyElemOperations.

Definition at line 355 of file stress_output_base.cpp.

Here is the call graph for this function:

Here is the caller graph for this function:

void MAST::StressStrainOutputBase::functional_boundary_sensitivity_for_all_elems ( const MAST::FunctionBase f,
Real dsigma_vm_val_df 
) const
virtual

calculates and returns the sensitivity of von Mises p-norm functional for all the elements that this object currently stores data for.

This is defined as

\[ \frac{ \frac{1}{p} \left( \int_\Omega (\sigma_{VM}(\Omega))^p ~ d\Omega \right)^{\frac{1}{p}-1}}{\left( \int_\Omega ~ d\Omega \right)^{\frac{1}{p}}} \int_\Gamma V_n \sigma_{VM}^p ~d\Gamma + \frac{ \frac{-1}{p} \left( \int_\Omega (\sigma_{VM}(\Omega))^p ~ d\Omega \right)^{\frac{1}{p}}}{\left( \int_\Omega ~ d\Omega \right)^{\frac{1+p}{p}}} \int_\Gamma V_n ~d\Gamma \]

Definition at line 914 of file stress_output_base.cpp.

Here is the call graph for this function:

Here is the caller graph for this function:

void MAST::StressStrainOutputBase::functional_boundary_sensitivity_for_elem ( const MAST::FunctionBase f,
const libMesh::dof_id_type  e_id,
Real dsigma_vm_val_df 
) const
virtual

calculates and returns the boundary sensitivity of von Mises p-norm functional for the element e.

Reimplemented in MAST::KSStressStrainOutput, and MAST::SmoothRampStressStrainOutput.

Definition at line 1016 of file stress_output_base.cpp.

Here is the caller graph for this function:

void MAST::StressStrainOutputBase::functional_for_all_elems ( )
virtual

calculates and returns the von Mises p-norm functional for all the elements that this object currently stores data for.

This is defined as

\[ \left( \frac{\int_\Omega (\sigma_{VM}(\Omega))^p ~ d\Omega}{\int_\Omega ~ d\Omega} \right)^{\frac{1}{p}} \]

Reimplemented in MAST::KSStressStrainOutput, and MAST::SmoothRampStressStrainOutput.

Definition at line 822 of file stress_output_base.cpp.

Here is the call graph for this function:

Here is the caller graph for this function:

void MAST::StressStrainOutputBase::functional_sensitivity_for_all_elems ( const MAST::FunctionBase f,
Real dsigma_vm_val_df 
) const
virtual

calculates and returns the sensitivity of von Mises p-norm functional for all the elements that this object currently stores data for.

This is defined as

\[ \frac{ \frac{1}{p} \left( \int_\Omega (\sigma_{VM}(\Omega))^p ~ d\Omega \right)^{\frac{1}{p}-1}}{\left( \int_\Omega ~ d\Omega \right)^{\frac{1}{p}}} \int_\Omega p (\sigma_{VM}(\Omega))^{p-1} \frac{d \sigma_{VM}(\Omega)}{d\alpha} ~ d\Omega \]

Definition at line 882 of file stress_output_base.cpp.

Here is the call graph for this function:

Here is the caller graph for this function:

void MAST::StressStrainOutputBase::functional_sensitivity_for_elem ( const MAST::FunctionBase f,
const libMesh::dof_id_type  e_id,
Real dsigma_vm_val_df 
) const
virtual

calculates and returns the sensitivity of von Mises p-norm functional for the element e.

Reimplemented in MAST::KSStressStrainOutput, and MAST::SmoothRampStressStrainOutput.

Definition at line 945 of file stress_output_base.cpp.

Here is the call graph for this function:

Here is the caller graph for this function:

void MAST::StressStrainOutputBase::functional_state_derivartive_for_elem ( const libMesh::dof_id_type  e_id,
RealVectorX dq_dX 
) const
virtual

calculates and returns the derivative of von Mises p-norm functional wrt state vector for the specified element.

This assumes that the von_Mises_p_norm_functional_for_all_elems() has been called to calculate the primal data. This is defined as

\[ \frac{ \frac{1}{p} \left( \int_\Omega (\sigma_{VM}(\Omega))^p ~ d\Omega \right)^{\frac{1}{p}-1}}{\left( \int_\Omega ~ d\Omega \right)^{\frac{1}{p}}} \int_\Omega p (\sigma_{VM}(\Omega))^{p-1} \frac{d \sigma_{VM}(\Omega)}{dX} ~ d\Omega \]

Reimplemented in MAST::KSStressStrainOutput, and MAST::SmoothRampStressStrainOutput.

Definition at line 1070 of file stress_output_base.cpp.

Here is the caller graph for this function:

Real MAST::StressStrainOutputBase::get_maximum_von_mises_stress ( ) const
Returns
the maximum von Mises stress of all stored components

Definition at line 686 of file stress_output_base.cpp.

Here is the call graph for this function:

Here is the caller graph for this function:

Real MAST::StressStrainOutputBase::get_p_stress_val ( )
inline
Returns
the $p-$norm for calculation of stress functional

Definition at line 257 of file stress_output_base.h.

const std::map< const libMesh::dof_id_type, std::vector< MAST::StressStrainOutputBase::Data * > > & MAST::StressStrainOutputBase::get_stress_strain_data ( ) const
virtual
Returns
the map of stress/strain data for all elems

Reimplemented in MAST::StressTemperatureAdjoint.

Definition at line 678 of file stress_output_base.cpp.

Here is the caller graph for this function:

const std::vector< MAST::StressStrainOutputBase::Data * > & MAST::StressStrainOutputBase::get_stress_strain_data_for_elem ( const MAST::GeomElem e) const
virtual
Returns
the vector of stress/strain data for specified elem.

Reimplemented in MAST::StressTemperatureAdjoint.

Definition at line 752 of file stress_output_base.cpp.

Here is the call graph for this function:

Here is the caller graph for this function:

MAST::StressStrainOutputBase::Data & MAST::StressStrainOutputBase::get_stress_strain_data_for_elem_at_qp ( const MAST::GeomElem e,
const unsigned int  qp 
)
virtual
Returns
the vector of stress/strain data for specified elem at the specified quadrature point.

Reimplemented in MAST::StressTemperatureAdjoint.

Definition at line 768 of file stress_output_base.cpp.

Here is the call graph for this function:

Here is the caller graph for this function:

MAST::BoundaryConditionBase * MAST::StressStrainOutputBase::get_thermal_load_for_elem ( const MAST::GeomElem elem)
Returns
the thermal load for this element, if present in the the volume loads. A nullptr pointer is returned if no load is present.

Definition at line 787 of file stress_output_base.cpp.

Here is the call graph for this function:

Here is the caller graph for this function:

void MAST::StressStrainOutputBase::init ( const MAST::GeomElem elem)
virtual

initialize for the element.

Implements MAST::AssemblyElemOperations.

Definition at line 473 of file stress_output_base.cpp.

Here is the call graph for this function:

Here is the caller graph for this function:

unsigned int MAST::StressStrainOutputBase::n_boundary_stress_strain_data_for_elem ( const GeomElem e) const
Returns
the number of points for which stress-strain data is stored for the boundary of element.

Definition at line 735 of file stress_output_base.cpp.

Here is the call graph for this function:

Here is the caller graph for this function:

unsigned int MAST::StressStrainOutputBase::n_stress_strain_data_for_elem ( const MAST::GeomElem e) const
Returns
the number of points for which stress-strain data is stored for the given element.

Definition at line 718 of file stress_output_base.cpp.

Here is the call graph for this function:

Here is the caller graph for this function:

void MAST::StressStrainOutputBase::output_derivative_for_elem ( RealVectorX dq_dX)
virtual

calculates the derivative of p-norm von Mises stress for the $p-$norm identified using set_p_val().

The quantity is evaluated over the current element for which this object is initialized.

Implements MAST::OutputAssemblyElemOperations.

Reimplemented in MAST::StressTemperatureAdjoint.

Definition at line 432 of file stress_output_base.cpp.

Here is the call graph for this function:

Here is the caller graph for this function:

virtual Real MAST::StressStrainOutputBase::output_for_elem ( )
inlinevirtual

should not get called for this output.

Use output_total() instead.

Implements MAST::OutputAssemblyElemOperations.

Definition at line 349 of file stress_output_base.h.

Here is the call graph for this function:

Real MAST::StressStrainOutputBase::output_sensitivity_for_elem ( const MAST::FunctionBase p)
virtual
Returns
the sensitivity of p-norm von Mises stress for the $p-$norm identified using set_p_val(). The returned quantity is evaluated for the element for which this object is initialized.

Implements MAST::OutputAssemblyElemOperations.

Definition at line 398 of file stress_output_base.cpp.

Here is the call graph for this function:

Here is the caller graph for this function:

Real MAST::StressStrainOutputBase::output_sensitivity_total ( const MAST::FunctionBase p)
virtual
Returns
the output quantity sensitivity for parameter. This method calculates the partial derivative of quantity

\[ \frac{\partial q(X, p)}{\partial p} \]

with respect to parameter $ p $. This returns the quantity accumulated over all elements.

Implements MAST::OutputAssemblyElemOperations.

Definition at line 414 of file stress_output_base.cpp.

Here is the call graph for this function:

Here is the caller graph for this function:

Real MAST::StressStrainOutputBase::output_total ( )
virtual
Returns
the output quantity value accumulated over all elements

Implements MAST::OutputAssemblyElemOperations.

Definition at line 385 of file stress_output_base.cpp.

Here is the call graph for this function:

Here is the caller graph for this function:

bool MAST::StressStrainOutputBase::primal_data_initialized ( ) const
inline

Definition at line 390 of file stress_output_base.h.

Here is the call graph for this function:

Here is the caller graph for this function:

void MAST::StressStrainOutputBase::set_aggregation_coefficients ( Real  p1,
Real  p2,
Real  rho,
Real  sigma0 
)
inline

sets the $p-$norm for calculation of stress functional

Definition at line 246 of file stress_output_base.h.

void MAST::StressStrainOutputBase::set_elem_data ( unsigned int  dim,
const libMesh::Elem &  ref_elem,
MAST::GeomElem elem 
) const
virtual

sets the structural element y-vector if 1D element is used.

Implements MAST::AssemblyElemOperations.

Definition at line 456 of file stress_output_base.cpp.

Here is the call graph for this function:

Here is the caller graph for this function:

void MAST::StressStrainOutputBase::set_stress_plot_mode ( bool  f)
inline

tells the object that the calculation is for stress to be output for plotting.

Definition at line 267 of file stress_output_base.h.

Here is the call graph for this function:

Here is the caller graph for this function:

bool MAST::StressStrainOutputBase::stress_plot_mode ( ) const
inline

Definition at line 386 of file stress_output_base.h.

Here is the caller graph for this function:

void MAST::StressStrainOutputBase::zero_for_analysis ( )
virtual

zeroes the output quantity values stored inside this object so that assembly process can begin.

This will zero out data so that it is ready for a new evaluation. Before sensitivity analysis, call the other method, since some nonlinear functionals need the forward quantities for sensitivity analysis, eg., stress output.

Implements MAST::OutputAssemblyElemOperations.

Definition at line 300 of file stress_output_base.cpp.

Here is the caller graph for this function:

void MAST::StressStrainOutputBase::zero_for_sensitivity ( )
virtual

zeroes the output quantity values stored inside this object so that assembly process can begin.

This will only zero the data to compute new sensitivity analysis.

Implements MAST::OutputAssemblyElemOperations.

Definition at line 310 of file stress_output_base.cpp.

Here is the caller graph for this function:

Member Data Documentation

std::map<const libMesh::dof_id_type, std::vector<MAST::StressStrainOutputBase::Data*> > MAST::StressStrainOutputBase::_boundary_stress_data
protected

vector of stress with the associated location details

Definition at line 617 of file stress_output_base.h.

Real MAST::StressStrainOutputBase::_exp_arg_lim
protected

Definition at line 590 of file stress_output_base.h.

bool MAST::StressStrainOutputBase::_if_stress_plot_mode
protected

identifies the mode in which evaluation is peformed.

if p-norm functional is being evaluated then certain requirements are not enforced. This is to be used when stress is being calculated per element for plotting.

Definition at line 604 of file stress_output_base.h.

Real MAST::StressStrainOutputBase::_JxW_val
protected

Definition at line 596 of file stress_output_base.h.

Real MAST::StressStrainOutputBase::_p_norm_stress
protected

$ p-$norm to be used for calculation of output stress function.

Default value is 2.0.

Definition at line 578 of file stress_output_base.h.

Real MAST::StressStrainOutputBase::_p_norm_weight
protected

Definition at line 578 of file stress_output_base.h.

bool MAST::StressStrainOutputBase::_primal_data_initialized
protected

primal data, needed for sensitivity and adjoints

Definition at line 595 of file stress_output_base.h.

Real MAST::StressStrainOutputBase::_rho
protected

exponent used in scaling volume based on stress value.

Definition at line 583 of file stress_output_base.h.

Real MAST::StressStrainOutputBase::_sigma0
protected

reference stress value used in scaling volume.

Definition at line 588 of file stress_output_base.h.

Real MAST::StressStrainOutputBase::_sigma_vm_int
protected

Definition at line 596 of file stress_output_base.h.

Real MAST::StressStrainOutputBase::_sigma_vm_p_norm
protected

Definition at line 596 of file stress_output_base.h.

std::map<const libMesh::dof_id_type, std::vector<MAST::StressStrainOutputBase::Data*> > MAST::StressStrainOutputBase::_stress_data
protected

vector of stress with the associated location details

Definition at line 610 of file stress_output_base.h.


The documentation for this class was generated from the following files: