oomph::ParallelResidualsHandler Class Reference

#include <assembly_handler.h>

+ Inheritance diagram for oomph::ParallelResidualsHandler:

Public Member Functions

 ParallelResidualsHandler (AssemblyHandler *const &assembly_handler_pt)
 Constructor, set the original assembly handler. More...
 
unsigned ndof (GeneralisedElement *const &elem_pt)
 
unsigned long eqn_number (GeneralisedElement *const &elem_pt, const unsigned &ieqn_local)
 
void get_residuals (GeneralisedElement *const &elem_pt, Vector< double > &residuals)
 
void get_jacobian (GeneralisedElement *const &elem_pt, Vector< double > &residuals, DenseMatrix< double > &jacobian)
 
void get_all_vectors_and_matrices (GeneralisedElement *const &elem_pt, Vector< Vector< double >> &vec, Vector< DenseMatrix< double >> &matrix)
 
 ~ParallelResidualsHandler ()
 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 doublelocal_problem_dof (Problem *const &problem_pt, const unsigned &t, const unsigned &i)
 
virtual void get_dresiduals_dparameter (GeneralisedElement *const &elem_pt, double *const &parameter_pt, Vector< double > &dres_dparam)
 
virtual void get_djacobian_dparameter (GeneralisedElement *const &elem_pt, double *const &parameter_pt, Vector< double > &dres_dparam, DenseMatrix< double > &djac_dparam)
 
virtual void get_hessian_vector_products (GeneralisedElement *const &elem_pt, Vector< double > const &Y, DenseMatrix< double > const &C, DenseMatrix< double > &product)
 
virtual int bifurcation_type () const
 
virtual doublebifurcation_parameter_pt () const
 
virtual void get_eigenfunction (Vector< DoubleVector > &eigenfunction)
 
virtual void get_inner_products (GeneralisedElement *const &elem_pt, Vector< std::pair< unsigned, unsigned >> const &history_index, Vector< double > &inner_product)
 
virtual void get_inner_product_vectors (GeneralisedElement *const &elem_pt, Vector< unsigned > const &history_index, Vector< Vector< double >> &inner_product_vector)
 
virtual ~AssemblyHandler ()
 Empty virtual destructor. More...
 

Private Attributes

AssemblyHandlerAssembly_handler_pt
 The original assembly handler. More...
 

Detailed Description

A class that is used to assemble the residuals in parallel by overloading the get_all_vectors_and_matrices, so that only the residuals are returned. This ensures that the (moderately complex) distributed parallel assembly loops are only in one place.

Constructor & Destructor Documentation

◆ ParallelResidualsHandler()

oomph::ParallelResidualsHandler::ParallelResidualsHandler ( AssemblyHandler *const &  assembly_handler_pt)
inline

Constructor, set the original assembly handler.

276  : Assembly_handler_pt(assembly_handler_pt)
277  {
278  }
AssemblyHandler * Assembly_handler_pt
The original assembly handler.
Definition: assembly_handler.h:271

◆ ~ParallelResidualsHandler()

oomph::ParallelResidualsHandler::~ParallelResidualsHandler ( )
inline

Empty virtual destructor.

327 {}

Member Function Documentation

◆ eqn_number()

unsigned long oomph::ParallelResidualsHandler::eqn_number ( GeneralisedElement *const &  elem_pt,
const unsigned ieqn_local 
)
inlinevirtual

Use underlying AssemblyHandler to return the global equation number of the local unknown ieqn_local in elem_pt.

Reimplemented from oomph::AssemblyHandler.

291  {
292  return Assembly_handler_pt->eqn_number(elem_pt, ieqn_local);
293  }
virtual unsigned long eqn_number(GeneralisedElement *const &elem_pt, const unsigned &ieqn_local)
Definition: assembly_handler.cc:82

References Assembly_handler_pt, and oomph::AssemblyHandler::eqn_number().

◆ get_all_vectors_and_matrices()

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

Calculate all desired vectors and matrices provided by the element elem_pt This function calls only the get_residuals function associated with the original assembly handler

Reimplemented from oomph::AssemblyHandler.

322  {
323  Assembly_handler_pt->get_residuals(elem_pt, vec[0]);
324  }
virtual void get_residuals(GeneralisedElement *const &elem_pt, Vector< double > &residuals)
Return the contribution to the residuals of the element elem_pt.
Definition: assembly_handler.cc:92

References Assembly_handler_pt, and oomph::AssemblyHandler::get_residuals().

◆ get_jacobian()

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

Use underlying AssemblyHandler to Calculate the elemental Jacobian matrix "d equation / d variable" for elem_pt.

Reimplemented from oomph::AssemblyHandler.

310  {
311  Assembly_handler_pt->get_jacobian(elem_pt, residuals, jacobian);
312  }
virtual void get_jacobian(GeneralisedElement *const &elem_pt, Vector< double > &residuals, DenseMatrix< double > &jacobian)
Definition: assembly_handler.cc:102

References Assembly_handler_pt, and oomph::AssemblyHandler::get_jacobian().

◆ get_residuals()

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

Use underlying AssemblyHandler to return the contribution to the residuals of the element elem_pt

Reimplemented from oomph::AssemblyHandler.

299  {
300  Assembly_handler_pt->get_residuals(elem_pt, residuals);
301  }

References Assembly_handler_pt, and oomph::AssemblyHandler::get_residuals().

◆ ndof()

unsigned oomph::ParallelResidualsHandler::ndof ( GeneralisedElement *const &  elem_pt)
inlinevirtual

Use underlying assembly handler to return the number of degrees of freedom in the element elem_pt

Reimplemented from oomph::AssemblyHandler.

283  {
284  return Assembly_handler_pt->ndof(elem_pt);
285  }
virtual unsigned ndof(GeneralisedElement *const &elem_pt)
Return the number of degrees of freedom in the element elem_pt.
Definition: assembly_handler.cc:42

References Assembly_handler_pt, and oomph::AssemblyHandler::ndof().

Member Data Documentation

◆ Assembly_handler_pt

AssemblyHandler* oomph::ParallelResidualsHandler::Assembly_handler_pt
private

The original assembly handler.

Referenced by eqn_number(), get_all_vectors_and_matrices(), get_jacobian(), get_residuals(), and ndof().


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