oomph::AssemblyHandler Class Reference

#include <assembly_handler.h>

+ Inheritance diagram for oomph::AssemblyHandler:

Public Member Functions

 AssemblyHandler ()
 Empty constructor. More...
 
virtual unsigned ndof (GeneralisedElement *const &elem_pt)
 Return the number of degrees of freedom in the element elem_pt. 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 unsigned long eqn_number (GeneralisedElement *const &elem_pt, const unsigned &ieqn_local)
 
virtual void get_residuals (GeneralisedElement *const &elem_pt, Vector< double > &residuals)
 Return the contribution to the residuals of the element elem_pt. More...
 
virtual void get_jacobian (GeneralisedElement *const &elem_pt, Vector< double > &residuals, DenseMatrix< double > &jacobian)
 
virtual void get_all_vectors_and_matrices (GeneralisedElement *const &elem_pt, Vector< Vector< double >> &vec, Vector< DenseMatrix< double >> &matrix)
 
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...
 

Detailed Description

A class that is used to define the functions used to assemble the elemental contributions to the residuals vector and Jacobian matrix that define the problem being solved. The main use of this class is to assemble and solve the augmented systems used in bifurcation detection and tracking. The default implementation merely calls the underlying elemental functions with no augmentation.

Constructor & Destructor Documentation

◆ AssemblyHandler()

oomph::AssemblyHandler::AssemblyHandler ( )
inline

Empty constructor.

66 {}

◆ ~AssemblyHandler()

virtual oomph::AssemblyHandler::~AssemblyHandler ( )
inlinevirtual

Empty virtual destructor.

169 {}

Member Function Documentation

◆ bifurcation_parameter_pt()

double * oomph::AssemblyHandler::bifurcation_parameter_pt ( ) const
virtual

Return a pointer to the bifurcation parameter in bifurcation tracking problems

Return the eigenfunction(s) associated with the bifurcation that has been detected in bifurcation tracking problems. Default Broken implementation

Reimplemented in oomph::HopfHandler, oomph::PitchForkHandler, and oomph::FoldHandler.

169  {
170  std::ostringstream error_stream;
171  error_stream << "There is no bifurcation parameter associated with the "
172  "current assembly handler.\n"
173  << "Eigenfunction are only calculated by the Fold, PitchFork "
174  "and Hopf Handlers"
175  << "\n";
176 
177  throw OomphLibError(
178  error_stream.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
179  return 0;
180  }
#define OOMPH_EXCEPTION_LOCATION
Definition: oomph_definitions.h:61
#define OOMPH_CURRENT_FUNCTION
Definition: oomph_definitions.h:86

References OOMPH_CURRENT_FUNCTION, and OOMPH_EXCEPTION_LOCATION.

Referenced by oomph::Problem::bifurcation_parameter_pt().

◆ bifurcation_type()

virtual int oomph::AssemblyHandler::bifurcation_type ( ) const
inlinevirtual

Return an unsigned integer to indicate whether the handler is a bifurcation tracking handler. The default is zero (not)

Reimplemented in oomph::HopfHandler, oomph::PitchForkHandler, and oomph::FoldHandler.

135  {
136  return 0;
137  }

Referenced by oomph::Problem::adapt(), oomph::Problem::doc_errors(), and oomph::Problem::p_adapt().

◆ dof_pt_vector()

void oomph::AssemblyHandler::dof_pt_vector ( GeneralisedElement *const &  elem_pt,
Vector< double * > &  dof_pt 
)
virtual

Return vector of pointers to dofs in the element elem_pt.

Get the vector of pointers to dofs in the element elem_pt Direct call to the function in the element.

64  {
65  return elem_pt->dof_pt_vector(dof_pt);
66  }

References oomph::GeneralisedElement::dof_pt_vector().

◆ dof_vector()

void oomph::AssemblyHandler::dof_vector ( GeneralisedElement *const &  elem_pt,
const unsigned t,
Vector< double > &  dof 
)
virtual

Return vector of dofs at time level t in the element elem_pt.

Get the vector of dofs in the element elem_pt at time level t Direct call to the function in the element.

54  {
55  return elem_pt->dof_vector(t, dof);
56  }
t
Definition: plotPSD.py:36

References oomph::GeneralisedElement::dof_vector(), and plotPSD::t.

◆ eqn_number()

◆ get_all_vectors_and_matrices()

void oomph::AssemblyHandler::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 element by calling the get_jacobian function

Reimplemented in oomph::ParallelResidualsHandler, oomph::EigenProblemHandler, and oomph::ExplicitTimeStepHandler.

117  {
118  get_jacobian(elem_pt, vec[0], matrix[0]);
119  }
virtual void get_jacobian(GeneralisedElement *const &elem_pt, Vector< double > &residuals, DenseMatrix< double > &jacobian)
Definition: assembly_handler.cc:102
Eigen::Map< Eigen::Matrix< T, Eigen::Dynamic, Eigen::Dynamic, Eigen::ColMajor >, 0, Eigen::OuterStride<> > matrix(T *data, int rows, int cols, int stride)
Definition: common.h:85

References get_jacobian(), and matrix().

Referenced by oomph::Problem::sparse_assemble_row_or_column_compressed_with_lists(), oomph::Problem::sparse_assemble_row_or_column_compressed_with_maps(), oomph::Problem::sparse_assemble_row_or_column_compressed_with_two_arrays(), oomph::Problem::sparse_assemble_row_or_column_compressed_with_two_vectors(), and oomph::Problem::sparse_assemble_row_or_column_compressed_with_vectors_of_pairs().

◆ get_djacobian_dparameter()

void oomph::AssemblyHandler::get_djacobian_dparameter ( GeneralisedElement *const &  elem_pt,
double *const &  parameter_pt,
Vector< double > &  dres_dparam,
DenseMatrix< double > &  djac_dparam 
)
virtual

Calculate the derivative of the residuals and jacobian with respect to a parameter

Calculate the derivative of the residuals and jacobian with respect to a parameter by calling the elemental function

Reimplemented in oomph::HopfHandler, oomph::PitchForkHandler, and oomph::FoldHandler.

143  {
144  elem_pt->get_djacobian_dparameter(parameter_pt, dres_dparam, djac_dparam);
145  }

References oomph::GeneralisedElement::get_djacobian_dparameter().

Referenced by oomph::ParameterDerivativeHandler::get_jacobian().

◆ get_dresiduals_dparameter()

void oomph::AssemblyHandler::get_dresiduals_dparameter ( GeneralisedElement *const &  elem_pt,
double *const &  parameter_pt,
Vector< double > &  dres_dparam 
)
virtual

Calculate the derivative of the residuals with respect to a parameter

Calculate the derivative of the residuals with respect to a parameter, by calling the elemental function

Reimplemented in oomph::HopfHandler, oomph::PitchForkHandler, and oomph::FoldHandler.

129  {
130  elem_pt->get_dresiduals_dparameter(parameter_pt, dres_dparam);
131  }

References oomph::GeneralisedElement::get_dresiduals_dparameter().

Referenced by oomph::ParameterDerivativeHandler::get_residuals().

◆ get_eigenfunction()

void oomph::AssemblyHandler::get_eigenfunction ( Vector< DoubleVector > &  eigenfunction)
virtual

Return the eigenfunction(s) associated with the bifurcation that has been detected in bifurcation tracking problems

Return the eigenfunction(s) associated with the bifurcation that has been detected in bifurcation tracking problems. Default Broken implementation

Reimplemented in oomph::HopfHandler, oomph::PitchForkHandler, and oomph::FoldHandler.

189  {
190  std::ostringstream error_stream;
191  error_stream << "There is no eigenfunction associated with the current "
192  "assembly handler.\n"
193  << "Eigenfunction are only calculated by the Fold, PitchFork "
194  "and Hopf Handlers"
195  << "\n";
196 
197  throw OomphLibError(
198  error_stream.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
199  }

References OOMPH_CURRENT_FUNCTION, and OOMPH_EXCEPTION_LOCATION.

Referenced by oomph::Problem::get_bifurcation_eigenfunction().

◆ get_hessian_vector_products()

void oomph::AssemblyHandler::get_hessian_vector_products ( GeneralisedElement *const &  elem_pt,
Vector< double > const &  Y,
DenseMatrix< double > const &  C,
DenseMatrix< double > &  product 
)
virtual

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}

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} At the moment the dof pointer is passed in to enable easy calculation of finite difference default

Reimplemented in oomph::HopfHandler, oomph::PitchForkHandler, and oomph::FoldHandler.

158  {
159  elem_pt->get_hessian_vector_products(Y, C, product);
160  }
Definition: matrices.h:74
void product(const MatrixType &m)
Definition: product.h:42
const char Y
Definition: test/EulerAngles.cpp:32

References oomph::GeneralisedElement::get_hessian_vector_products(), product(), and Y.

Referenced by oomph::Problem::get_hessian_vector_products().

◆ get_inner_product_vectors()

void oomph::AssemblyHandler::get_inner_product_vectors ( GeneralisedElement *const &  elem_pt,
Vector< unsigned > const &  history_index,
Vector< Vector< double >> &  inner_product_vector 
)
virtual

Compute the vectors that when taken as a dot product with other history values give the inner product over the element

Compute the vectors that when taken as a dot product with other history values give the inner product over the element. In other words if we call get_inner_product_vectors({0,1},output) output[0] will be a vector such that dofs.output[0] gives the inner product of current dofs with themselves.

225  {
226  elem_pt->get_inner_product_vectors(history_index, inner_product_vector);
227  }

References oomph::GeneralisedElement::get_inner_product_vectors().

◆ get_inner_products()

void oomph::AssemblyHandler::get_inner_products ( GeneralisedElement *const &  elem_pt,
Vector< std::pair< unsigned, unsigned >> const &  history_index,
Vector< double > &  inner_product 
)
virtual

Compute the inner products of the given vector of pairs of history values over the element.

========================================================================== Compute the inner products of the given vector of pairs of history values over the element. The values of the index in the pair may be the same.

210  {
211  elem_pt->get_inner_products(history_index, inner_product);
212  }

References oomph::GeneralisedElement::get_inner_products().

◆ get_jacobian()

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

◆ get_residuals()

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

◆ local_problem_dof()

double & oomph::AssemblyHandler::local_problem_dof ( Problem *const &  problem_pt,
const unsigned t,
const unsigned i 
)
virtual

Return the t-th level of storage associated with the i-th (local) dof stored in the problem

73  {
74  return *(problem_pt->dof_pt(i) + t);
75  }
int i
Definition: BiCGSTAB_step_by_step.cpp:9

References oomph::Problem::dof_pt(), i, and plotPSD::t.

◆ ndof()


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