MAST
MAST::AssemblyBase Class Reference

#include <assembly_base.h>

Inheritance diagram for MAST::AssemblyBase:
Collaboration diagram for MAST::AssemblyBase:

Classes

class  ElemParameterDependence
 Inherited objects from this class can be provided by the user provide assessment of whether or not an element is influenced by a give parameter. More...
 
class  SolverMonitor
 

Public Member Functions

 AssemblyBase ()
 constructor takes a reference to the discipline that provides the boundary conditions, volume loads, properties, etc. More...
 
virtual ~AssemblyBase ()
 virtual destructor More...
 
const MAST::PhysicsDisciplineBasediscipline () const
 
MAST::PhysicsDisciplineBasediscipline ()
 
MAST::AssemblyElemOperationsget_elem_ops ()
 
virtual void set_discipline_and_system (MAST::PhysicsDisciplineBase &discipline, MAST::SystemInitialization &system)
 attaches a system to this discipline More...
 
void attach_elem_parameter_dependence_object (MAST::AssemblyBase::ElemParameterDependence &dep)
 This object, if provided by user, will be used to reduce unnecessary computations in sensitivity analysis assembly operations. More...
 
void clear_elem_parameter_dependence_object ()
 
virtual void clear_discipline_and_system ()
 clears association with a system to this discipline More...
 
virtual void set_elem_operation_object (MAST::AssemblyElemOperations &elem_ops)
 attaches a element operation to this object, and associated this with the element operation object. More...
 
virtual void clear_elem_operation_object ()
 clears the association of this object with the assembly element operation object. More...
 
const MAST::NonlinearSystemsystem () const
 
MAST::NonlinearSystemsystem ()
 
MAST::SystemInitializationsystem_init ()
 
void set_solver_monitor (MAST::AssemblyBase::SolverMonitor &monitor)
 attaches the solver monitor, which is a user provided routine that is called each time More...
 
MAST::AssemblyBase::SolverMonitorget_solver_monitor ()
 
void clear_solver_monitor ()
 clears the monitor object More...
 
void attach_solution_function (MAST::MeshFieldFunction &f)
 tells the assembly object that this function is will need to be initialized before each residual evaluation More...
 
void detach_solution_function ()
 removes the attachment of the solution function More...
 
virtual void residual_and_jacobian (const libMesh::NumericVector< Real > &X, libMesh::NumericVector< Real > *R, libMesh::SparseMatrix< Real > *J, libMesh::NonlinearImplicitSystem &S)
 function that assembles the matrices and vectors quantities for nonlinear solution More...
 
virtual bool sensitivity_assemble (const MAST::FunctionBase &f, libMesh::NumericVector< Real > &sensitivity_rhs)
 Assembly function. More...
 
virtual void calculate_output (const libMesh::NumericVector< Real > &X, MAST::OutputAssemblyElemOperations &output)
 calculates the value of quantity $ q(X,p) $. More...
 
virtual void calculate_output_derivative (const libMesh::NumericVector< Real > &X, MAST::OutputAssemblyElemOperations &output, libMesh::NumericVector< Real > &dq_dX)
 calculates $ \frac{\partial q(X, p)}{\partial X} $ More...
 
virtual void calculate_output_direct_sensitivity (const libMesh::NumericVector< Real > &X, const libMesh::NumericVector< Real > *dXdp, const MAST::FunctionBase &p, MAST::OutputAssemblyElemOperations &output)
 evaluates the sensitivity of the outputs in the attached discipline with respect to the parametrs in params. More...
 
virtual Real calculate_output_adjoint_sensitivity (const libMesh::NumericVector< Real > &X, const libMesh::NumericVector< Real > &dq_dX, const MAST::FunctionBase &p, MAST::AssemblyElemOperations &elem_ops, MAST::OutputAssemblyElemOperations &output, const bool include_partial_sens=true)
 Evaluates the total sensitivity of output wrt p using the adjoint solution provided in dq_dX for a linearization about solution X. More...
 
std::unique_ptr< libMesh::NumericVector< Real > > build_localized_vector (const libMesh::System &sys, const libMesh::NumericVector< Real > &global) const
 localizes the parallel vector so that the local copy stores all values necessary for calculation of the element quantities More...
 

Public Attributes

bool close_matrix
 flag to control the closing fo the Jacobian after assembly More...
 

Protected Attributes

MAST::AssemblyElemOperations_elem_ops
 provides assembly elem operations for use by this class More...
 
MAST::PhysicsDisciplineBase_discipline
 PhysicsDisciplineBase object for which this class is assembling. More...
 
MAST::SystemInitialization_system
 System for which this assembly is performed. More...
 
MAST::MeshFieldFunction_sol_function
 system solution that will be initialized before each solution More...
 
MAST::AssemblyBase::SolverMonitor_solver_monitor
 User provided solver monitor is attached to the linear nonlinear solvers, if provided. More...
 
MAST::AssemblyBase::ElemParameterDependence_param_dependence
 If provided by user, this object is used by sensitiivty analysis to check for whether or the current design parameter influences an element. More...
 

Detailed Description

Definition at line 51 of file assembly_base.h.

Constructor & Destructor Documentation

MAST::AssemblyBase::AssemblyBase ( )

constructor takes a reference to the discipline that provides the boundary conditions, volume loads, properties, etc.

Definition at line 39 of file assembly_base.cpp.

MAST::AssemblyBase::~AssemblyBase ( )
virtual

virtual destructor

Definition at line 53 of file assembly_base.cpp.

Here is the call graph for this function:

Member Function Documentation

void MAST::AssemblyBase::attach_elem_parameter_dependence_object ( MAST::AssemblyBase::ElemParameterDependence dep)

This object, if provided by user, will be used to reduce unnecessary computations in sensitivity analysis assembly operations.

This association is cleard when clear_discipline_and_system() is called.

Definition at line 129 of file assembly_base.cpp.

Here is the caller graph for this function:

void MAST::AssemblyBase::attach_solution_function ( MAST::MeshFieldFunction f)

tells the assembly object that this function is will need to be initialized before each residual evaluation

Definition at line 234 of file assembly_base.cpp.

std::unique_ptr< libMesh::NumericVector< Real > > MAST::AssemblyBase::build_localized_vector ( const libMesh::System &  sys,
const libMesh::NumericVector< Real > &  global 
) const

localizes the parallel vector so that the local copy stores all values necessary for calculation of the element quantities

Definition at line 211 of file assembly_base.cpp.

Here is the caller graph for this function:

void MAST::AssemblyBase::calculate_output ( const libMesh::NumericVector< Real > &  X,
MAST::OutputAssemblyElemOperations output 
)
virtual

calculates the value of quantity $ q(X,p) $.

Reimplemented in MAST::LevelSetNonlinearImplicitAssembly.

Definition at line 254 of file assembly_base.cpp.

Here is the call graph for this function:

Here is the caller graph for this function:

Real MAST::AssemblyBase::calculate_output_adjoint_sensitivity ( const libMesh::NumericVector< Real > &  X,
const libMesh::NumericVector< Real > &  dq_dX,
const MAST::FunctionBase p,
MAST::AssemblyElemOperations elem_ops,
MAST::OutputAssemblyElemOperations output,
const bool  include_partial_sens = true 
)
virtual

Evaluates the total sensitivity of output wrt p using the adjoint solution provided in dq_dX for a linearization about solution X.

Definition at line 511 of file assembly_base.cpp.

Here is the call graph for this function:

Here is the caller graph for this function:

void MAST::AssemblyBase::calculate_output_derivative ( const libMesh::NumericVector< Real > &  X,
MAST::OutputAssemblyElemOperations output,
libMesh::NumericVector< Real > &  dq_dX 
)
virtual

calculates $ \frac{\partial q(X, p)}{\partial X} $

Reimplemented in MAST::LevelSetNonlinearImplicitAssembly.

Definition at line 328 of file assembly_base.cpp.

Here is the call graph for this function:

Here is the caller graph for this function:

void MAST::AssemblyBase::calculate_output_direct_sensitivity ( const libMesh::NumericVector< Real > &  X,
const libMesh::NumericVector< Real > *  dXdp,
const MAST::FunctionBase p,
MAST::OutputAssemblyElemOperations output 
)
virtual

evaluates the sensitivity of the outputs in the attached discipline with respect to the parametrs in params.

The base solution should be provided in X. If total sensitivity is desired, then dXdp should contain the sensitivity of solution wrt the parameter p, otherwise it can be set to nullptr. If dXdp is zero, the calculated sensitivity will be the partial derivarive of output wrt p.

Reimplemented in MAST::LevelSetNonlinearImplicitAssembly.

Definition at line 412 of file assembly_base.cpp.

Here is the call graph for this function:

Here is the caller graph for this function:

void MAST::AssemblyBase::clear_discipline_and_system ( )
virtual

clears association with a system to this discipline

Reimplemented in MAST::StructuralBucklingEigenproblemAssembly, MAST::FSIGeneralizedAeroForceAssembly, and MAST::StructuralFluidInteractionAssembly.

Definition at line 176 of file assembly_base.cpp.

Here is the caller graph for this function:

void MAST::AssemblyBase::clear_elem_operation_object ( )
virtual

clears the association of this object with the assembly element operation object.

Definition at line 199 of file assembly_base.cpp.

Here is the call graph for this function:

Here is the caller graph for this function:

void MAST::AssemblyBase::clear_elem_parameter_dependence_object ( )

Definition at line 139 of file assembly_base.cpp.

void MAST::AssemblyBase::clear_solver_monitor ( )

clears the monitor object

Definition at line 154 of file assembly_base.cpp.

Here is the call graph for this function:

void MAST::AssemblyBase::detach_solution_function ( )

removes the attachment of the solution function

Definition at line 246 of file assembly_base.cpp.

const MAST::PhysicsDisciplineBase & MAST::AssemblyBase::discipline ( ) const
Returns
a const reference to the PhysicsDisciplineBase object associated with this object

Definition at line 62 of file assembly_base.cpp.

Here is the caller graph for this function:

MAST::PhysicsDisciplineBase & MAST::AssemblyBase::discipline ( )
Returns
a non-const reference to the PhysicsDisciplineBase object associated with this object

Definition at line 72 of file assembly_base.cpp.

MAST::AssemblyElemOperations & MAST::AssemblyBase::get_elem_ops ( )
Returns
a reference to the element operations object

Definition at line 101 of file assembly_base.cpp.

Here is the caller graph for this function:

MAST::AssemblyBase::SolverMonitor * MAST::AssemblyBase::get_solver_monitor ( )
Returns
the solver monitor, which is a user provided routine that is called each time

Definition at line 147 of file assembly_base.cpp.

Here is the caller graph for this function:

virtual void MAST::AssemblyBase::residual_and_jacobian ( const libMesh::NumericVector< Real > &  X,
libMesh::NumericVector< Real > *  R,
libMesh::SparseMatrix< Real > *  J,
libMesh::NonlinearImplicitSystem &  S 
)
inlinevirtual

function that assembles the matrices and vectors quantities for nonlinear solution

Reimplemented in MAST::LevelSetNonlinearImplicitAssembly, MAST::NonlinearImplicitAssembly, MAST::TransientAssembly, and MAST::LevelSetReinitializationTransientAssembly.

Definition at line 220 of file assembly_base.h.

Here is the caller graph for this function:

virtual bool MAST::AssemblyBase::sensitivity_assemble ( const MAST::FunctionBase f,
libMesh::NumericVector< Real > &  sensitivity_rhs 
)
inlinevirtual

Assembly function.

This function will be called to assemble the RHS of the sensitivity equations (which is -1 times sensitivity of system residual) prior to a solve and must be provided by the user in a derived class. The method provides dR/dp for f parameter.

If the routine is not able to provide sensitivity for this parameter, then it should return false, and the system will attempt to use finite differencing.

Reimplemented in MAST::NonlinearImplicitAssembly, MAST::ComplexAssemblyBase, MAST::LevelSetNonlinearImplicitAssembly, and MAST::TransientAssembly.

Definition at line 239 of file assembly_base.h.

Here is the call graph for this function:

Here is the caller graph for this function:

void MAST::AssemblyBase::set_discipline_and_system ( MAST::PhysicsDisciplineBase discipline,
MAST::SystemInitialization system 
)
virtual

attaches a system to this discipline

Definition at line 162 of file assembly_base.cpp.

Here is the call graph for this function:

Here is the caller graph for this function:

void MAST::AssemblyBase::set_elem_operation_object ( MAST::AssemblyElemOperations elem_ops)
virtual

attaches a element operation to this object, and associated this with the element operation object.

Definition at line 187 of file assembly_base.cpp.

Here is the call graph for this function:

Here is the caller graph for this function:

void MAST::AssemblyBase::set_solver_monitor ( MAST::AssemblyBase::SolverMonitor monitor)

attaches the solver monitor, which is a user provided routine that is called each time

Definition at line 119 of file assembly_base.cpp.

Here is the call graph for this function:

const MAST::NonlinearSystem & MAST::AssemblyBase::system ( ) const
Returns
a const reference to the MAST::NonlinearSystem object associated with this object

Definition at line 83 of file assembly_base.cpp.

Here is the call graph for this function:

Here is the caller graph for this function:

MAST::NonlinearSystem & MAST::AssemblyBase::system ( )
Returns
a non-const reference to the MAST::NonlinearSystem object associated with this object

Definition at line 92 of file assembly_base.cpp.

Here is the call graph for this function:

MAST::SystemInitialization & MAST::AssemblyBase::system_init ( )
Returns
a non-const reference to the MAST::SystemInitialization

Definition at line 110 of file assembly_base.cpp.

Here is the caller graph for this function:

Member Data Documentation

MAST::PhysicsDisciplineBase* MAST::AssemblyBase::_discipline
protected

PhysicsDisciplineBase object for which this class is assembling.

Definition at line 312 of file assembly_base.h.

MAST::AssemblyElemOperations* MAST::AssemblyBase::_elem_ops
protected

provides assembly elem operations for use by this class

Definition at line 306 of file assembly_base.h.

MAST::AssemblyBase::ElemParameterDependence* MAST::AssemblyBase::_param_dependence
protected

If provided by user, this object is used by sensitiivty analysis to check for whether or the current design parameter influences an element.

This can be used to enhance computational efficiency.

Definition at line 336 of file assembly_base.h.

MAST::MeshFieldFunction* MAST::AssemblyBase::_sol_function
protected

system solution that will be initialized before each solution

Definition at line 323 of file assembly_base.h.

MAST::AssemblyBase::SolverMonitor* MAST::AssemblyBase::_solver_monitor
protected

User provided solver monitor is attached to the linear nonlinear solvers, if provided.

Definition at line 329 of file assembly_base.h.

MAST::SystemInitialization* MAST::AssemblyBase::_system
protected

System for which this assembly is performed.

Definition at line 318 of file assembly_base.h.

bool MAST::AssemblyBase::close_matrix

flag to control the closing fo the Jacobian after assembly

Definition at line 103 of file assembly_base.h.


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