MAST
MAST::LevelSetNonlinearImplicitAssembly Class Reference

#include <level_set_nonlinear_implicit_assembly.h>

Inheritance diagram for MAST::LevelSetNonlinearImplicitAssembly:
Collaboration diagram for MAST::LevelSetNonlinearImplicitAssembly:

Public Member Functions

 LevelSetNonlinearImplicitAssembly (bool enable_dof_handler)
 constructor associates this assembly object with the system More...
 
virtual ~LevelSetNonlinearImplicitAssembly ()
 destructor resets the association of this assembly object with the system More...
 
void set_evaluate_output_on_negative_phi (bool f)
 sets the flag on whether or not to evaluate the output on negative level set function More...
 
bool if_use_dof_handler () const
 
virtual void set_level_set_function (MAST::FieldFunction< Real > &level_set, const MAST::FilterBase &filter)
 attaches level set function to this More...
 
virtual void set_indicator_function (MAST::FieldFunction< RealVectorX > &indicator)
 attaches indicator function to this. More...
 
virtual void clear_level_set_function ()
 clears association with level set function More...
 
virtual void set_level_set_velocity_function (MAST::FieldFunction< RealVectorX > &velocity)
 the velocity function used to calculate topology sensitivity More...
 
virtual void clear_level_set_velocity_function ()
 clears the velocity function More...
 
MAST::LevelSetIntersectionget_intersection ()
 
MAST::LevelSetInterfaceDofHandlerget_dof_handler ()
 
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_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 (const libMesh::NumericVector< Real > &X, MAST::OutputAssemblyElemOperations &output)
 calculates the value of quantity $ q(X,p) $. 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...
 
- Public Member Functions inherited from MAST::NonlinearImplicitAssembly
 NonlinearImplicitAssembly ()
 constructor associates this assembly object with the system More...
 
virtual ~NonlinearImplicitAssembly ()
 destructor resets the association of this assembly object with the system More...
 
Real res_l2_norm () const
 L2 norm of the last-assembled residual. More...
 
Real first_iter_res_l2_norm () const
 
void reset_residual_norm_history ()
 reset L2 norm of the last-assembled residual More...
 
void set_post_assembly_operation (MAST::NonlinearImplicitAssembly::PostAssemblyOperation &post)
 sets the PostAssemblyOperation object for use after assembly. More...
 
virtual void linearized_jacobian_solution_product (const libMesh::NumericVector< Real > &X, const libMesh::NumericVector< Real > &dX, libMesh::NumericVector< Real > &JdX, libMesh::NonlinearImplicitSystem &S)
 calculates the product of the Jacobian and a perturbation in solution vector $ [J] \{\Delta X\} $. More...
 
virtual void second_derivative_dot_solution_assembly (const libMesh::NumericVector< Real > &X, const libMesh::NumericVector< Real > &dX, libMesh::SparseMatrix< Real > &d_JdX_dX, libMesh::NonlinearImplicitSystem &S)
 calculates $ d ([J] \{\Delta X\})/ dX $. More...
 
- Public Member Functions inherited from MAST::AssemblyBase
 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 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...
 

Protected Attributes

bool _enable_dof_handler
 Evaluates the total sensitivity of output wrt p using the adjoint solution provided in dq_dX for a linearization about solution X. More...
 
bool _evaluate_output_on_negative_phi
 
MAST::FieldFunction< Real > * _level_set
 
MAST::FieldFunction< RealVectorX > * _indicator
 
MAST::LevelSetIntersection_intersection
 
MAST::LevelSetInterfaceDofHandler_dof_handler
 
MAST::LevelSetVoidSolution_void_solution_monitor
 
MAST::FieldFunction< RealVectorX > * _velocity
 
const MAST::FilterBase_filter
 
- Protected Attributes inherited from MAST::NonlinearImplicitAssembly
MAST::NonlinearImplicitAssembly::PostAssemblyOperation_post_assembly
 this object, if non-NULL is user-provided to perform actions after assembly and before returning to the solver More...
 
Real _res_l2_norm
 L2 norm of the last-assembled residual. More...
 
Real _first_iter_res_l2_norm
 
- Protected Attributes inherited from MAST::AssemblyBase
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...
 

Additional Inherited Members

- Public Attributes inherited from MAST::AssemblyBase
bool close_matrix
 flag to control the closing fo the Jacobian after assembly More...
 

Detailed Description

Definition at line 37 of file level_set_nonlinear_implicit_assembly.h.

Constructor & Destructor Documentation

MAST::LevelSetNonlinearImplicitAssembly::LevelSetNonlinearImplicitAssembly ( bool  enable_dof_handler)

constructor associates this assembly object with the system

Definition at line 47 of file level_set_nonlinear_implicit_assembly.cpp.

MAST::LevelSetNonlinearImplicitAssembly::~LevelSetNonlinearImplicitAssembly ( )
virtual

destructor resets the association of this assembly object with the system

Definition at line 63 of file level_set_nonlinear_implicit_assembly.cpp.

Member Function Documentation

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

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

Reimplemented from MAST::AssemblyBase.

Definition at line 761 of file level_set_nonlinear_implicit_assembly.cpp.

Here is the call graph for this function:

Here is the caller graph for this function:

void MAST::LevelSetNonlinearImplicitAssembly::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 from MAST::AssemblyBase.

Definition at line 916 of file level_set_nonlinear_implicit_assembly.cpp.

Here is the call graph for this function:

Here is the caller graph for this function:

void MAST::LevelSetNonlinearImplicitAssembly::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 nullptr, the calculated sensitivity will be the partial derivarive of output wrt p.

Reimplemented from MAST::AssemblyBase.

Definition at line 1139 of file level_set_nonlinear_implicit_assembly.cpp.

Here is the call graph for this function:

Here is the caller graph for this function:

void MAST::LevelSetNonlinearImplicitAssembly::clear_level_set_function ( )
virtual

clears association with level set function

Definition at line 137 of file level_set_nonlinear_implicit_assembly.cpp.

Here is the call graph for this function:

void MAST::LevelSetNonlinearImplicitAssembly::clear_level_set_velocity_function ( )
virtual

clears the velocity function

Definition at line 171 of file level_set_nonlinear_implicit_assembly.cpp.

Here is the call graph for this function:

MAST::LevelSetInterfaceDofHandler & MAST::LevelSetNonlinearImplicitAssembly::get_dof_handler ( )
Returns
a reference to the LevelSetInterfaceDofHandler object

Definition at line 94 of file level_set_nonlinear_implicit_assembly.cpp.

Here is the call graph for this function:

MAST::LevelSetIntersection & MAST::LevelSetNonlinearImplicitAssembly::get_intersection ( )
Returns
a reference to the level set function

Definition at line 72 of file level_set_nonlinear_implicit_assembly.cpp.

bool MAST::LevelSetNonlinearImplicitAssembly::if_use_dof_handler ( ) const
Returns
flag if using dof_handler or not

Definition at line 88 of file level_set_nonlinear_implicit_assembly.cpp.

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

function that assembles the matrices and vectors quantities for nonlinear solution

Reimplemented from MAST::NonlinearImplicitAssembly.

Definition at line 322 of file level_set_nonlinear_implicit_assembly.cpp.

Here is the call graph for this function:

Here is the caller graph for this function:

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

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 from MAST::NonlinearImplicitAssembly.

Definition at line 553 of file level_set_nonlinear_implicit_assembly.cpp.

Here is the call graph for this function:

Here is the caller graph for this function:

void MAST::LevelSetNonlinearImplicitAssembly::set_evaluate_output_on_negative_phi ( bool  f)

sets the flag on whether or not to evaluate the output on negative level set function

Definition at line 81 of file level_set_nonlinear_implicit_assembly.cpp.

void MAST::LevelSetNonlinearImplicitAssembly::set_indicator_function ( MAST::FieldFunction< RealVectorX > &  indicator)
virtual

attaches indicator function to this.

Definition at line 126 of file level_set_nonlinear_implicit_assembly.cpp.

Here is the caller graph for this function:

void MAST::LevelSetNonlinearImplicitAssembly::set_level_set_function ( MAST::FieldFunction< Real > &  level_set,
const MAST::FilterBase filter 
)
virtual

attaches level set function to this

Definition at line 104 of file level_set_nonlinear_implicit_assembly.cpp.

Here is the call graph for this function:

Here is the caller graph for this function:

void MAST::LevelSetNonlinearImplicitAssembly::set_level_set_velocity_function ( MAST::FieldFunction< RealVectorX > &  velocity)
virtual

the velocity function used to calculate topology sensitivity

Definition at line 158 of file level_set_nonlinear_implicit_assembly.cpp.

Here is the caller graph for this function:

Member Data Documentation

MAST::LevelSetInterfaceDofHandler* MAST::LevelSetNonlinearImplicitAssembly::_dof_handler
protected

Definition at line 200 of file level_set_nonlinear_implicit_assembly.h.

bool MAST::LevelSetNonlinearImplicitAssembly::_enable_dof_handler
protected

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 190 of file level_set_nonlinear_implicit_assembly.h.

bool MAST::LevelSetNonlinearImplicitAssembly::_evaluate_output_on_negative_phi
protected

Definition at line 192 of file level_set_nonlinear_implicit_assembly.h.

const MAST::FilterBase* MAST::LevelSetNonlinearImplicitAssembly::_filter
protected

Definition at line 206 of file level_set_nonlinear_implicit_assembly.h.

MAST::FieldFunction<RealVectorX>* MAST::LevelSetNonlinearImplicitAssembly::_indicator
protected

Definition at line 196 of file level_set_nonlinear_implicit_assembly.h.

MAST::LevelSetIntersection* MAST::LevelSetNonlinearImplicitAssembly::_intersection
protected

Definition at line 198 of file level_set_nonlinear_implicit_assembly.h.

MAST::FieldFunction<Real>* MAST::LevelSetNonlinearImplicitAssembly::_level_set
protected

Definition at line 194 of file level_set_nonlinear_implicit_assembly.h.

MAST::FieldFunction<RealVectorX>* MAST::LevelSetNonlinearImplicitAssembly::_velocity
protected

Definition at line 204 of file level_set_nonlinear_implicit_assembly.h.

MAST::LevelSetVoidSolution* MAST::LevelSetNonlinearImplicitAssembly::_void_solution_monitor
protected

Definition at line 202 of file level_set_nonlinear_implicit_assembly.h.


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