MAST
MAST::OutputAssemblyElemOperations Class Referenceabstract

This provides the base class for definitin of element level contribution of output quantity in an analysis. More...

#include <output_assembly_elem_operations.h>

Inheritance diagram for MAST::OutputAssemblyElemOperations:
Collaboration diagram for MAST::OutputAssemblyElemOperations:

Public Member Functions

 OutputAssemblyElemOperations ()
 
virtual ~OutputAssemblyElemOperations ()
 virtual destructor More...
 
virtual void zero_for_analysis ()=0
 zeroes the output quantity values stored inside this object so that assembly process can begin. More...
 
virtual void zero_for_sensitivity ()=0
 zeroes the output quantity values stored inside this object so that assembly process can begin. More...
 
virtual Real output_for_elem ()=0
 
virtual Real output_total ()=0
 
virtual Real output_sensitivity_for_elem (const MAST::FunctionBase &p)=0
 
virtual Real output_sensitivity_total (const MAST::FunctionBase &p)=0
 
virtual void output_derivative_for_elem (RealVectorX &dq_dX)=0
 returns the output quantity derivative with respect to state vector in dq_dX. More...
 
virtual void evaluate ()=0
 this is the abstract interface to be implemented by derived classes. More...
 
virtual void evaluate_sensitivity (const MAST::FunctionBase &f)=0
 this evaluates all relevant sensitivity components on the element. More...
 
virtual void evaluate_shape_sensitivity (const MAST::FunctionBase &f)=0
 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)=0
 this evaluates all relevant topological sensitivity components on the element. 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 set_elem_data (unsigned int dim, const libMesh::Elem &ref_elem, MAST::GeomElem &elem) const =0
 some analyses may want to set additional element data before initialization of the GeomElem. More...
 
virtual void init (const MAST::GeomElem &elem)=0
 initializes the object for calculation of element quantities for the specified elem. 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

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

This provides the base class for definitin of element level contribution of output quantity in an analysis.

Apart from providing evaluation of the quantity $ q(X, p) $, the implementations are expected to evaluate the sensitivity $ \frac{\partial q}{\partial p} $ and the derivative $ \frac{\partial q}{\partial X} $. The latter quantity is used in an adjoint analysis.

Definition at line 54 of file output_assembly_elem_operations.h.

Constructor & Destructor Documentation

MAST::OutputAssemblyElemOperations::OutputAssemblyElemOperations ( )

Definition at line 31 of file output_assembly_elem_operations.cpp.

MAST::OutputAssemblyElemOperations::~OutputAssemblyElemOperations ( )
virtual

virtual destructor

Definition at line 38 of file output_assembly_elem_operations.cpp.

Here is the call graph for this function:

Member Function Documentation

virtual void MAST::OutputAssemblyElemOperations::evaluate ( )
pure virtual

this is the abstract interface to be implemented by derived classes.

This method calculates the quantity $ q(X, p) $. This is provided as an interface since some quantities requires evaluation of all inputs before the final quantity can be calculated. So, the user can call this routine on each element before requesting the output values using either output_for_elem() or output_total().

Implemented in MAST::StressStrainOutputBase, MAST::LevelSetPerimeter, MAST::LevelSetVolume, MAST::IntegratedForceOutput, and MAST::ComplianceOutput.

Here is the caller graph for this function:

virtual void MAST::OutputAssemblyElemOperations::evaluate_sensitivity ( const MAST::FunctionBase f)
pure virtual

this evaluates all relevant sensitivity components on the element.

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

Implemented in MAST::StressStrainOutputBase, MAST::LevelSetPerimeter, MAST::LevelSetVolume, MAST::IntegratedForceOutput, and MAST::ComplianceOutput.

Here is the caller graph for this function:

virtual void MAST::OutputAssemblyElemOperations::evaluate_shape_sensitivity ( const MAST::FunctionBase f)
pure virtual

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.

Implemented in MAST::StressStrainOutputBase, MAST::LevelSetPerimeter, MAST::LevelSetVolume, MAST::IntegratedForceOutput, and MAST::ComplianceOutput.

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

this evaluates all relevant topological sensitivity components on the element.

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

Implemented in MAST::StressStrainOutputBase, MAST::LevelSetPerimeter, MAST::LevelSetVolume, MAST::IntegratedForceOutput, and MAST::ComplianceOutput.

Here is the caller graph for this function:

const std::set< libMesh::boundary_id_type > & MAST::OutputAssemblyElemOperations::get_participating_boundaries ( )
Returns
the set of boundary ids for which data will be stored. This is set using the set_participating_boundaries() method.

Definition at line 113 of file output_assembly_elem_operations.cpp.

Here is the call graph for this function:

const std::set< const libMesh::Elem * > & MAST::OutputAssemblyElemOperations::get_participating_elements ( ) const
Returns
the set of elements for which data will be stored. This is set using the set_participating_elements() method.

Definition at line 97 of file output_assembly_elem_operations.cpp.

const std::set< libMesh::subdomain_id_type > & MAST::OutputAssemblyElemOperations::get_participating_subdomains ( )
Returns
the set of subdomain ids for which data will be stored. This is set using the set_participating_subdomains() method.

Definition at line 105 of file output_assembly_elem_operations.cpp.

bool MAST::OutputAssemblyElemOperations::if_evaluate_for_boundary ( const MAST::GeomElem elem,
const unsigned int  s 
) const
virtual

checks to see if the specified side of the element needs evaluation of the output contribution.

Definition at line 155 of file output_assembly_elem_operations.cpp.

Here is the call graph for this function:

Here is the caller graph for this function:

bool MAST::OutputAssemblyElemOperations::if_evaluate_for_element ( const MAST::GeomElem elem) const
virtual

checks to see if the object has been told about the subset of elements and if the specified element is in the subset.

Reimplemented in MAST::StressTemperatureAdjoint.

Definition at line 121 of file output_assembly_elem_operations.cpp.

Here is the call graph for this function:

Here is the caller graph for this function:

virtual void MAST::OutputAssemblyElemOperations::output_derivative_for_elem ( RealVectorX dq_dX)
pure virtual

returns the output quantity derivative with respect to state vector in dq_dX.

This method calculates the quantity

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

for this output function. This is returned for the element for which this

Implemented in MAST::StressStrainOutputBase, MAST::ComplianceOutput, MAST::LevelSetPerimeter, MAST::LevelSetVolume, MAST::IntegratedForceOutput, and MAST::StressTemperatureAdjoint.

Here is the caller graph for this function:

virtual Real MAST::OutputAssemblyElemOperations::output_for_elem ( )
pure virtual
Returns
the output quantity value contribution for the present element for which this object has been initialized. This does not make sense for some output quantities that are defined as averaged over a domain, example the p-norm of the stress (see MAST::StressStrainOutputBase). For such cases, the output_total must be called.

Implemented in MAST::StressStrainOutputBase, MAST::ComplianceOutput, MAST::LevelSetPerimeter, MAST::LevelSetVolume, and MAST::IntegratedForceOutput.

virtual Real MAST::OutputAssemblyElemOperations::output_sensitivity_for_elem ( const MAST::FunctionBase p)
pure 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 evaluated for on element for which this object is initialized.

Implemented in MAST::StressStrainOutputBase, MAST::ComplianceOutput, MAST::LevelSetPerimeter, MAST::LevelSetVolume, and MAST::IntegratedForceOutput.

virtual Real MAST::OutputAssemblyElemOperations::output_sensitivity_total ( const MAST::FunctionBase p)
pure 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.

Implemented in MAST::StressStrainOutputBase, MAST::ComplianceOutput, MAST::LevelSetPerimeter, MAST::LevelSetVolume, and MAST::IntegratedForceOutput.

Here is the caller graph for this function:

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

Implemented in MAST::StressStrainOutputBase, MAST::ComplianceOutput, MAST::LevelSetPerimeter, MAST::LevelSetVolume, and MAST::IntegratedForceOutput.

void MAST::OutputAssemblyElemOperations::set_participating_boundaries ( const std::set< libMesh::boundary_id_type > &  bids)

The assembly will integration over boudnaries with ids specified in bids.

Definition at line 88 of file output_assembly_elem_operations.cpp.

Here is the caller graph for this function:

void MAST::OutputAssemblyElemOperations::set_participating_elements ( const std::set< const libMesh::Elem * > &  elems)

sets the elements for which this object will evaluate and store the output data.

This allows the user to specify a smaller subset of elements that will be grouped together in the output for constraint evaluation. If this method is not called, then the object will store data for all elements in the subdomain.

Definition at line 73 of file output_assembly_elem_operations.cpp.

Here is the call graph for this function:

Here is the caller graph for this function:

void MAST::OutputAssemblyElemOperations::set_participating_elements_to_all ( )

This will allow volume contribution from all elements.

Definition at line 61 of file output_assembly_elem_operations.cpp.

Here is the call graph for this function:

void MAST::OutputAssemblyElemOperations::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.

The user must identify the subdomains on which the quantity will be defined using this method, or the set_participating_elements() method.

Definition at line 47 of file output_assembly_elem_operations.cpp.

Here is the caller graph for this function:

virtual void MAST::OutputAssemblyElemOperations::zero_for_analysis ( )
pure 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.

Implemented in MAST::StressStrainOutputBase, MAST::ComplianceOutput, MAST::LevelSetPerimeter, MAST::LevelSetVolume, and MAST::IntegratedForceOutput.

Here is the caller graph for this function:

virtual void MAST::OutputAssemblyElemOperations::zero_for_sensitivity ( )
pure 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.

Implemented in MAST::StressStrainOutputBase, MAST::ComplianceOutput, MAST::LevelSetPerimeter, MAST::LevelSetVolume, and MAST::IntegratedForceOutput.

Here is the caller graph for this function:

Member Data Documentation

std::set<libMesh::boundary_id_type> MAST::OutputAssemblyElemOperations::_bids
protected

set of bids for which data will be processed

Definition at line 259 of file output_assembly_elem_operations.h.

std::set<const libMesh::Elem*> MAST::OutputAssemblyElemOperations::_elem_subset
protected

set of elements for which the data will be stored.

If this is empty, then data for all elements will be stored.

Definition at line 247 of file output_assembly_elem_operations.h.

bool MAST::OutputAssemblyElemOperations::_if_evaluate_on_all_elems
protected

if true, evaluates on all elements.

Definition at line 241 of file output_assembly_elem_operations.h.

std::set<libMesh::subdomain_id_type> MAST::OutputAssemblyElemOperations::_sub_domain_ids
protected

set of subdomain ids for which data will be processed.

Definition at line 253 of file output_assembly_elem_operations.h.


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