Public Member Functions | List of all members
oomph::ExplicitTimeStepHandler Class Reference

A class that is used to define the functions used to assemble and invert the mass matrix when taking an explicit timestep. The idea is simply to replace the jacobian matrix with the mass matrix and then our standard linear solvers will solve the required system. More...

#include <assembly_handler.h>

+ Inheritance diagram for oomph::ExplicitTimeStepHandler:

Public Member Functions

 ExplicitTimeStepHandler ()
 Empty Constructor. More...
 
unsigned ndof (GeneralisedElement *const &elem_pt)
 Return the number of degrees of freedom in the element elem_pt. More...
 
unsigned long eqn_number (GeneralisedElement *const &elem_pt, const unsigned &ieqn_local)
 Return the global equation number of the local unknown ieqn_local in elem_pt. More...
 
void get_residuals (GeneralisedElement *const &elem_pt, Vector< double > &residuals)
 Return the contribution to the residuals of the element elem_pt This is deliberately broken in our eigenproblem. More...
 
void get_jacobian (GeneralisedElement *const &elem_pt, Vector< double > &residuals, DenseMatrix< double > &jacobian)
 Calculate the elemental Jacobian matrix "d equation / d variable" for elem_pt. Again deliberately broken in the eigenproblem. More...
 
void get_all_vectors_and_matrices (GeneralisedElement *const &elem_pt, Vector< Vector< double >> &vec, Vector< DenseMatrix< double >> &matrix)
 Calculate all desired vectors and matrices provided by the element elem_pt. More...
 
 ~ExplicitTimeStepHandler ()
 Empty virtual destructor. More...
 
- Public Member Functions inherited from oomph::AssemblyHandler
 AssemblyHandler ()
 Empty constructor. More...
 
virtual void dof_vector (GeneralisedElement *const &elem_pt, const unsigned &t, Vector< double > &dof)
 Return vector of dofs at time level t in the element elem_pt. More...
 
virtual void dof_pt_vector (GeneralisedElement *const &elem_pt, Vector< double * > &dof_pt)
 Return vector of pointers to dofs in the element elem_pt. More...
 
virtual double & local_problem_dof (Problem *const &problem_pt, const unsigned &t, const unsigned &i)
 Return the t-th level of storage associated with the i-th (local) dof stored in the problem. More...
 
virtual void get_dresiduals_dparameter (GeneralisedElement *const &elem_pt, double *const &parameter_pt, Vector< double > &dres_dparam)
 Calculate the derivative of the residuals with respect to a parameter. More...
 
virtual void get_djacobian_dparameter (GeneralisedElement *const &elem_pt, double *const &parameter_pt, Vector< double > &dres_dparam, DenseMatrix< double > &djac_dparam)
 Calculate the derivative of the residuals and jacobian with respect to a parameter. More...
 
virtual void get_hessian_vector_products (GeneralisedElement *const &elem_pt, Vector< double > const &Y, DenseMatrix< double > const &C, DenseMatrix< double > &product)
 Calculate the product of the Hessian (derivative of Jacobian with respect to all variables) an eigenvector, Y, and other specified vectors, C (d(J_{ij})/d u_{k}) Y_{j} C_{k}. More...
 
virtual int bifurcation_type () const
 Return an unsigned integer to indicate whether the handler is a bifurcation tracking handler. The default is zero (not) More...
 
virtual double * bifurcation_parameter_pt () const
 Return a pointer to the bifurcation parameter in bifurcation tracking problems. More...
 
virtual void get_eigenfunction (Vector< DoubleVector > &eigenfunction)
 Return the eigenfunction(s) associated with the bifurcation that has been detected in bifurcation tracking problems. More...
 
virtual void get_inner_products (GeneralisedElement *const &elem_pt, Vector< std::pair< unsigned, unsigned >> const &history_index, Vector< double > &inner_product)
 Compute the inner products of the given vector of pairs of history values over the element. More...
 
virtual void get_inner_product_vectors (GeneralisedElement *const &elem_pt, Vector< unsigned > const &history_index, Vector< Vector< double >> &inner_product_vector)
 Compute the vectors that when taken as a dot product with other history values give the inner product over the element. More...
 
virtual void synchronise ()
 Function that is used to perform any synchronisation required during the solution. More...
 
virtual ~AssemblyHandler ()
 Empty virtual destructor. More...
 

Detailed Description

A class that is used to define the functions used to assemble and invert the mass matrix when taking an explicit timestep. The idea is simply to replace the jacobian matrix with the mass matrix and then our standard linear solvers will solve the required system.

Definition at line 180 of file assembly_handler.h.

Constructor & Destructor Documentation

◆ ExplicitTimeStepHandler()

oomph::ExplicitTimeStepHandler::ExplicitTimeStepHandler ( )
inline

Empty Constructor.

Definition at line 184 of file assembly_handler.h.

◆ ~ExplicitTimeStepHandler()

oomph::ExplicitTimeStepHandler::~ExplicitTimeStepHandler ( )
inline

Empty virtual destructor.

Definition at line 212 of file assembly_handler.h.

Member Function Documentation

◆ eqn_number()

unsigned long oomph::ExplicitTimeStepHandler::eqn_number ( GeneralisedElement *const &  elem_pt,
const unsigned &  ieqn_local 
)
virtual

Return the global equation number of the local unknown ieqn_local in elem_pt.

Get the global equation number of the local unknown. Direct call to the function in the element.

Reimplemented from oomph::AssemblyHandler.

Definition at line 248 of file assembly_handler.cc.

References oomph::GeneralisedElement::eqn_number().

◆ get_all_vectors_and_matrices()

void oomph::ExplicitTimeStepHandler::get_all_vectors_and_matrices ( GeneralisedElement *const &  elem_pt,
Vector< Vector< double >> &  vec,
Vector< DenseMatrix< double >> &  matrix 
)
virtual

Calculate all desired vectors and matrices provided by the element elem_pt.

Calculate all desired vectors and matrices that are required by the problem by calling those of the underlying element.

Reimplemented from oomph::AssemblyHandler.

Definition at line 278 of file assembly_handler.cc.

References oomph::GeneralisedElement::get_mass_matrix().

◆ get_jacobian()

void oomph::ExplicitTimeStepHandler::get_jacobian ( GeneralisedElement *const &  elem_pt,
Vector< double > &  residuals,
DenseMatrix< double > &  jacobian 
)
virtual

Calculate the elemental Jacobian matrix "d equation / d variable" for elem_pt. Again deliberately broken in the eigenproblem.

Replace get jacobian with the call to get the mass matrix.

Reimplemented from oomph::AssemblyHandler.

Definition at line 266 of file assembly_handler.cc.

References oomph::GeneralisedElement::get_mass_matrix().

◆ get_residuals()

void oomph::ExplicitTimeStepHandler::get_residuals ( GeneralisedElement *const &  elem_pt,
Vector< double > &  residuals 
)
virtual

Return the contribution to the residuals of the element elem_pt This is deliberately broken in our eigenproblem.

Call the element's residuals.

Reimplemented from oomph::AssemblyHandler.

Definition at line 257 of file assembly_handler.cc.

References oomph::GeneralisedElement::get_residuals().

◆ ndof()

unsigned oomph::ExplicitTimeStepHandler::ndof ( GeneralisedElement *const &  elem_pt)
virtual

Return the number of degrees of freedom in the element elem_pt.

Get the number of elemental degrees of freedom. Direct call to the function in the element.

Reimplemented from oomph::AssemblyHandler.

Definition at line 239 of file assembly_handler.cc.

References oomph::GeneralisedElement::ndof().


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