MAST
|
#include <transient_solver_base.h>
Public Member Functions | |
TransientSolverBase (unsigned int o, unsigned int n) | |
constructor requires the number of iterations to store for the derived solver. More... | |
virtual | ~TransientSolverBase () |
virtual void | set_assembly (MAST::AssemblyBase &assembly) |
sets the assembly object More... | |
virtual void | clear_assembly () |
clears the assembly object More... | |
virtual void | set_elem_operation_object (MAST::TransientAssemblyElemOperations &elem_ops) |
Attaches the assembly elem operations object that provides the x_dot, M and J quantities for the element. More... | |
virtual MAST::TransientAssemblyElemOperations & | get_elem_operation_object () |
virtual void | clear_elem_operation_object () |
Clears the assembly elem operations object. More... | |
libMesh::NumericVector< Real > & | solution (unsigned int prev_iter=0) const |
libMesh::NumericVector< Real > & | solution_sensitivity (unsigned int prev_iter=0) const |
libMesh::NumericVector< Real > & | velocity (unsigned int prev_iter=0) const |
libMesh::NumericVector< Real > & | velocity_sensitivity (unsigned int prev_iter=0) const |
libMesh::NumericVector< Real > & | acceleration (unsigned int prev_iter=0) const |
libMesh::NumericVector< Real > & | acceleration_sensitivity (unsigned int prev_iter=0) const |
virtual void | update_velocity (libMesh::NumericVector< Real > &vel, const libMesh::NumericVector< Real > &sol)=0 |
update the transient velocity based on the current solution More... | |
virtual void | update_acceleration (libMesh::NumericVector< Real > &acc, const libMesh::NumericVector< Real > &sol)=0 |
update the transient acceleration based on the current solution More... | |
virtual void | update_sensitivity_velocity (libMesh::NumericVector< Real > &vel, const libMesh::NumericVector< Real > &sol)=0 |
update the transient sensitivity velocity based on the current sensitivity solution More... | |
virtual void | update_sensitivity_acceleration (libMesh::NumericVector< Real > &acc, const libMesh::NumericVector< Real > &sol)=0 |
update the transient sensitivity acceleration based on the current sensitivity solution More... | |
virtual void | update_delta_velocity (libMesh::NumericVector< Real > &vel, const libMesh::NumericVector< Real > &sol)=0 |
update the perturbation in transient velocity based on the current perturbed solution More... | |
virtual void | update_delta_acceleration (libMesh::NumericVector< Real > &acc, const libMesh::NumericVector< Real > &sol)=0 |
update the perturbation in transient acceleration based on the current perturbed solution More... | |
virtual void | solve (MAST::AssemblyBase &assembly) |
solves the current time step for solution and velocity More... | |
virtual void | sensitivity_solve (MAST::AssemblyBase &assembly, const MAST::FunctionBase &f) |
solvers the current time step for sensitivity wrt f More... | |
void | solve_highest_derivative_and_advance_time_step (MAST::AssemblyBase &assembly) |
To be used only for initial conditions. More... | |
void | solve_highest_derivative_and_advance_time_step_with_sensitivity (MAST::AssemblyBase &assembly, const MAST::FunctionBase &f) |
solves for the sensitivity of highest derivative and advances the time-step. More... | |
virtual void | advance_time_step () |
advances the time step and copies the current solution to old solution, and so on. More... | |
virtual void | advance_time_step_with_sensitivity () |
advances the time step and copies the current sensitivity solution to old sensitivity solution, and so on. More... | |
virtual void | build_local_quantities (const libMesh::NumericVector< Real > ¤t_sol, std::vector< libMesh::NumericVector< Real > * > &qtys) |
localizes the relevant solutions for system assembly. More... | |
virtual void | build_sensitivity_local_quantities (unsigned int prev_iter, std::vector< libMesh::NumericVector< Real > * > &qtys) |
localizes the relevant solutions for system assembly. More... | |
Public Member Functions inherited from MAST::NonlinearImplicitAssemblyElemOperations | |
NonlinearImplicitAssemblyElemOperations () | |
virtual | ~NonlinearImplicitAssemblyElemOperations () |
virtual void | elem_calculations (bool if_jac, RealVectorX &vec, RealMatrixX &mat)=0 |
performs the element calculations over elem , and returns the element vector and matrix quantities in mat and vec , respectively. More... | |
virtual void | elem_linearized_jacobian_solution_product (RealVectorX &vec)=0 |
performs the element calculations over elem , and returns the element vector quantity in vec . More... | |
virtual void | elem_sensitivity_calculations (const MAST::FunctionBase &f, RealVectorX &vec)=0 |
performs the element sensitivity calculations over elem , and returns the element residual sensitivity in vec . More... | |
virtual void | elem_shape_sensitivity_calculations (const MAST::FunctionBase &f, RealVectorX &vec)=0 |
performs the element shape sensitivity calculations over elem , and returns the element residual sensitivity in vec . More... | |
virtual void | elem_topology_sensitivity_calculations (const MAST::FunctionBase &f, const MAST::FieldFunction< RealVectorX > &vel, RealVectorX &vec)=0 |
performs the element topology sensitivity calculations over elem , and returns the element residual sensitivity in vec . More... | |
virtual void | elem_second_derivative_dot_solution_assembly (RealMatrixX &mat)=0 |
calculates over elem , and returns the matrix in vec . More... | |
void | check_element_numerical_jacobian (RealVectorX &sol) |
a helper function to evaluate the numerical Jacobian and compare it with the analytical Jacobian. More... | |
Public Member Functions inherited from MAST::AssemblyElemOperations | |
AssemblyElemOperations () | |
virtual | ~AssemblyElemOperations () |
MAST::SystemInitialization & | get_system_initialization () |
MAST::PhysicsDisciplineBase & | get_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 MAST::AssemblyBase & | get_assembly () |
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::ElementBase & | get_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... | |
Public Attributes | |
Real | dt |
time step More... | |
Additional Inherited Members | |
Protected Attributes inherited from MAST::AssemblyElemOperations | |
MAST::SystemInitialization * | _system |
MAST::PhysicsDisciplineBase * | _discipline |
MAST::AssemblyBase * | _assembly |
MAST::ElementBase * | _physics_elem |
Definition at line 40 of file transient_solver_base.h.
MAST::TransientSolverBase::TransientSolverBase | ( | unsigned int | o, |
unsigned int | n | ||
) |
constructor requires the number of iterations to store for the derived solver.
Definition at line 36 of file transient_solver_base.cpp.
|
virtual |
libMesh::NumericVector< Real > & MAST::TransientSolverBase::acceleration | ( | unsigned int | prev_iter = 0 | ) | const |
prev_iter
= 0 gives the current acceleration estimate. Note that prev_iter
cannot be greater than the total number of iterations that this solver stores solutions for. Definition at line 282 of file transient_solver_base.cpp.
libMesh::NumericVector< Real > & MAST::TransientSolverBase::acceleration_sensitivity | ( | unsigned int | prev_iter = 0 | ) | const |
prev_iter
parameter is similar to the acceleration() method. Definition at line 301 of file transient_solver_base.cpp.
|
virtual |
advances the time step and copies the current solution to old solution, and so on.
Definition at line 752 of file transient_solver_base.cpp.
|
virtual |
advances the time step and copies the current sensitivity solution to old sensitivity solution, and so on.
Does not influence the primary solution.
Definition at line 794 of file transient_solver_base.cpp.
|
virtual |
localizes the relevant solutions for system assembly.
The calling function has to delete the pointers to these vectors
Definition at line 534 of file transient_solver_base.cpp.
|
virtual |
localizes the relevant solutions for system assembly.
The calling function has to delete the pointers to these vectors. prev_iter
= 0 implies the current iterate, while increasing values will identify decreasing iterations for which data is stored.
Definition at line 612 of file transient_solver_base.cpp.
|
virtual |
clears the assembly object
Reimplemented from MAST::AssemblyElemOperations.
Definition at line 69 of file transient_solver_base.cpp.
|
virtual |
Clears the assembly elem operations object.
Definition at line 145 of file transient_solver_base.cpp.
|
virtual |
Definition at line 78 of file transient_solver_base.cpp.
|
virtual |
solvers the current time step for sensitivity wrt f
Reimplemented in MAST::StabilizedFirstOrderNewmarkTransientSensitivitySolver.
Definition at line 333 of file transient_solver_base.cpp.
|
virtual |
sets the assembly object
Reimplemented from MAST::AssemblyElemOperations.
Definition at line 60 of file transient_solver_base.cpp.
|
virtual |
Attaches the assembly elem operations object that provides the x_dot, M and J quantities for the element.
Definition at line 88 of file transient_solver_base.cpp.
libMesh::NumericVector< Real > & MAST::TransientSolverBase::solution | ( | unsigned int | prev_iter = 0 | ) | const |
prev_iter
= 0 gives the current solution estimate. Note that prev_iter
cannot be greater than the total number of iterations that this solver stores solutions for. Definition at line 197 of file transient_solver_base.cpp.
libMesh::NumericVector< Real > & MAST::TransientSolverBase::solution_sensitivity | ( | unsigned int | prev_iter = 0 | ) | const |
prev_iter
parameter is similar to the solution() method. Definition at line 220 of file transient_solver_base.cpp.
|
virtual |
solves the current time step for solution and velocity
Definition at line 320 of file transient_solver_base.cpp.
void MAST::TransientSolverBase::solve_highest_derivative_and_advance_time_step | ( | MAST::AssemblyBase & | assembly | ) |
To be used only for initial conditions.
Initializes the highest derivative solution using the solution and its low order time derivatives at the specified time step. Then advances the time step so that the solver is ready for time integration.
Definition at line 347 of file transient_solver_base.cpp.
void MAST::TransientSolverBase::solve_highest_derivative_and_advance_time_step_with_sensitivity | ( | MAST::AssemblyBase & | assembly, |
const MAST::FunctionBase & | f | ||
) |
solves for the sensitivity of highest derivative and advances the time-step.
Definition at line 440 of file transient_solver_base.cpp.
|
pure virtual |
update the transient acceleration based on the current solution
Implemented in MAST::StabilizedFirstOrderNewmarkTransientSensitivitySolver, and MAST::SecondOrderNewmarkTransientSolver.
|
pure virtual |
update the perturbation in transient acceleration based on the current perturbed solution
Implemented in MAST::StabilizedFirstOrderNewmarkTransientSensitivitySolver, and MAST::SecondOrderNewmarkTransientSolver.
|
pure virtual |
update the perturbation in transient velocity based on the current perturbed solution
Implemented in MAST::StabilizedFirstOrderNewmarkTransientSensitivitySolver, and MAST::SecondOrderNewmarkTransientSolver.
|
pure virtual |
update the transient sensitivity acceleration based on the current sensitivity solution
Implemented in MAST::StabilizedFirstOrderNewmarkTransientSensitivitySolver, and MAST::SecondOrderNewmarkTransientSolver.
|
pure virtual |
update the transient sensitivity velocity based on the current sensitivity solution
Implemented in MAST::StabilizedFirstOrderNewmarkTransientSensitivitySolver, and MAST::SecondOrderNewmarkTransientSolver.
|
pure virtual |
update the transient velocity based on the current solution
Implemented in MAST::StabilizedFirstOrderNewmarkTransientSensitivitySolver, and MAST::SecondOrderNewmarkTransientSolver.
libMesh::NumericVector< Real > & MAST::TransientSolverBase::velocity | ( | unsigned int | prev_iter = 0 | ) | const |
prev_iter
= 0 gives the current velocity estimate. Note that prev_iter
cannot be greater than the total number of iterations that this solver stores solutions for. Definition at line 243 of file transient_solver_base.cpp.
libMesh::NumericVector< Real > & MAST::TransientSolverBase::velocity_sensitivity | ( | unsigned int | prev_iter = 0 | ) | const |
prev_iter
parameter is similar to the velocity() method. Definition at line 262 of file transient_solver_base.cpp.
Real MAST::TransientSolverBase::dt |
time step
Definition at line 77 of file transient_solver_base.h.