21 #ifndef __mast__transient_solver_base__    22 #define __mast__transient_solver_base__    29 #include "libmesh/numeric_vector.h"    35     class TransientAssemblyElemOperations;
    37     class NonlinearSystem;
    87         libMesh::NumericVector<Real>&
    88         solution(
unsigned int prev_iter = 0) 
const;
    94         libMesh::NumericVector<Real>&
   104         libMesh::NumericVector<Real>&
   105         velocity(
unsigned int prev_iter = 0) 
const;
   112         libMesh::NumericVector<Real>&
   123         libMesh::NumericVector<Real>&
   130         libMesh::NumericVector<Real>&
   138                                      const libMesh::NumericVector<Real>& sol) = 0;
   144                                          const libMesh::NumericVector<Real>& sol) = 0;
   151                                                  const libMesh::NumericVector<Real>& sol) = 0;
   158                                                      const libMesh::NumericVector<Real>& sol) = 0;
   167                               const libMesh::NumericVector<Real>& sol) = 0;
   175                                   const libMesh::NumericVector<Real>& sol) = 0;
   200                                                        bool increment_time = 
true);
   234                                std::vector<libMesh::NumericVector<Real>*>& qtys);
   245                                            std::vector<libMesh::NumericVector<Real>*>& qtys);
   256         build_perturbed_local_quantities
   257         (
const libMesh::NumericVector<Real>& current_sol,
   258          std::vector<libMesh::NumericVector<Real>*>& qtys);
   265         set_element_data(
const std::vector<libMesh::dof_id_type>& dof_indices,
   266                          const std::vector<libMesh::NumericVector<Real>*>& sols) = 0;
   274         extract_element_sensitivity_data(
const std::vector<libMesh::dof_id_type>& dof_indices,
   275                                          const std::vector<libMesh::NumericVector<Real>*>& sols,
   276                                          std::vector<RealVectorX>& local_sols) = 0;
   290         set_element_perturbed_data
   291         (
const std::vector<libMesh::dof_id_type>& dof_indices,
   292          const std::vector<libMesh::NumericVector<Real>*>& sols) = 0;
   300                       const libMesh::Elem& ref_elem,
   319         bool  _first_step, _first_sensitivity_step;
   326         const unsigned int _ode_order;
   332         const unsigned int _n_iters_to_store;
   345         bool   _if_highest_derivative_solution;
   352 #endif // __mast__transient_solver_base__ 
virtual void build_local_quantities(const libMesh::NumericVector< Real > ¤t_sol, std::vector< libMesh::NumericVector< Real > *> &qtys)
localizes the relevant solutions for system assembly. 
virtual void advance_time_step_with_sensitivity()
advances the time step and copies the current sensitivity solution to old sensitivity solution...
virtual void advance_time_step(bool increment_time=true)
advances the time step and copies the current solution to old solution, and so on. 
virtual void elem_sensitivity_contribution_previous_timestep(const std::vector< RealVectorX > &prev_sols, RealVectorX &vec)
computes the contribution for this element from previous time step 
virtual void sensitivity_solve(MAST::AssemblyBase &assembly, const MAST::FunctionBase &f)
solvers the current time step for sensitivity wrt f 
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 ...
virtual void set_assembly(MAST::AssemblyBase &assembly)
sets the assembly object 
virtual MAST::TransientAssemblyElemOperations & get_elem_operation_object()
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 ...
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 ...
virtual ~TransientSolverBase()
libMesh::NumericVector< Real > & solution(unsigned int prev_iter=0) const
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 elem...
libMesh::NumericVector< Real > & solution_sensitivity(unsigned int prev_iter=0) const
libMesh::NumericVector< Real > & acceleration_sensitivity(unsigned int prev_iter=0) const
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...
virtual void update_velocity(libMesh::NumericVector< Real > &vel, const libMesh::NumericVector< Real > &sol)=0
update the transient velocity based on the current solution 
libMesh::NumericVector< Real > & velocity(unsigned int prev_iter=0) const
libMesh::NumericVector< Real > & acceleration(unsigned int prev_iter=0) const
Matrix< Real, Dynamic, 1 > RealVectorX
virtual void init(const MAST::GeomElem &elem)=0
initializes the object for calculation of element quantities for the specified elem. 
virtual void update_acceleration(libMesh::NumericVector< Real > &acc, const libMesh::NumericVector< Real > &sol)=0
update the transient acceleration based on the current solution 
void solve_highest_derivative_and_advance_time_step(MAST::AssemblyBase &assembly, bool increment_time=true)
To be used only for initial conditions. 
virtual void clear_assembly()
clears the assembly object 
This class acts as a wrapper around libMesh::Elem for the purpose of providing a uniform interface fo...
virtual void solve(MAST::AssemblyBase &assembly)
solves the current time step for solution and velocity 
virtual void clear_elem()
clears the element initialization 
virtual void build_sensitivity_local_quantities(unsigned int prev_iter, std::vector< libMesh::NumericVector< Real > *> &qtys)
localizes the relevant solutions for system assembly. 
libMesh::NumericVector< Real > & velocity_sensitivity(unsigned int prev_iter=0) const
TransientSolverBase(unsigned int o, unsigned int n)
constructor requires the number of iterations to store for the derived solver. 
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 ...
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. 
virtual void clear_elem_operation_object()
Clears the assembly elem operations object.