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. More...
#include <assembly_handler.h>
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 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 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... | |
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) |
Calculate the elemental Jacobian matrix "d equation
/ d variable" for elem_pt. More... | |
virtual 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... | |
virtual void | get_dresiduals_dparameter (GeneralisedElement *const &elem_pt, double *const ¶meter_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 ¶meter_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... | |
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.
Definition at line 62 of file assembly_handler.h.
|
inline |
Empty constructor.
Definition at line 66 of file assembly_handler.h.
|
inlinevirtual |
Empty virtual destructor.
Definition at line 169 of file assembly_handler.h.
|
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.
Definition at line 168 of file assembly_handler.cc.
Referenced by oomph::Problem::bifurcation_parameter_pt().
|
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.
Definition at line 134 of file assembly_handler.h.
Referenced by oomph::Problem::adapt(), oomph::Problem::doc_errors(), and oomph::Problem::p_adapt().
|
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.
Definition at line 62 of file assembly_handler.cc.
References oomph::GeneralisedElement::dof_pt_vector().
|
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.
Definition at line 51 of file assembly_handler.cc.
References oomph::GeneralisedElement::dof_vector(), and t.
|
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 in oomph::PeriodicOrbitAssemblyHandler< NNODE_1D >, oomph::HopfHandler, oomph::PitchForkHandler, oomph::FoldHandler, oomph::ParameterDerivativeHandler, oomph::ParallelResidualsHandler, oomph::EigenProblemHandler, and oomph::ExplicitTimeStepHandler.
Definition at line 82 of file assembly_handler.cc.
References oomph::GeneralisedElement::eqn_number().
Referenced by oomph::ParallelResidualsHandler::eqn_number(), oomph::ParameterDerivativeHandler::eqn_number(), oomph::Problem::get_hessian_vector_products(), oomph::Problem::get_jacobian(), oomph::PitchForkHandler::get_jacobian(), oomph::Problem::get_my_eqns(), oomph::Problem::get_residuals(), oomph::Problem::parallel_sparse_assemble(), oomph::PitchForkHandler::PitchForkHandler(), oomph::HSL_MA42::reorder_elements(), oomph::Problem::setup_element_count_per_dof(), oomph::HSL_MA42::solve(), 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().
|
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.
Definition at line 113 of file assembly_handler.cc.
References get_jacobian().
Referenced by oomph::Problem::parallel_sparse_assemble(), 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().
|
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.
Definition at line 138 of file assembly_handler.cc.
References oomph::GeneralisedElement::get_djacobian_dparameter().
Referenced by oomph::ParameterDerivativeHandler::get_jacobian().
|
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.
Definition at line 125 of file assembly_handler.cc.
References oomph::GeneralisedElement::get_dresiduals_dparameter().
Referenced by oomph::ParameterDerivativeHandler::get_residuals().
|
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.
Definition at line 188 of file assembly_handler.cc.
Referenced by oomph::Problem::get_bifurcation_eigenfunction().
|
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.
Definition at line 153 of file assembly_handler.cc.
References oomph::GeneralisedElement::get_hessian_vector_products().
Referenced by oomph::Problem::get_hessian_vector_products().
|
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.
Definition at line 221 of file assembly_handler.cc.
References oomph::GeneralisedElement::get_inner_product_vectors().
|
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.
Definition at line 206 of file assembly_handler.cc.
References oomph::GeneralisedElement::get_inner_products().
|
virtual |
Calculate the elemental Jacobian matrix "d equation / d variable" for elem_pt.
Calculate the elemental Jacobian matrix "d equation / d variable" by calling the element's get_jacobian function.
Reimplemented in oomph::FpPreconditionerAssemblyHandler, oomph::PeriodicOrbitAssemblyHandler< NNODE_1D >, oomph::HopfHandler, oomph::PitchForkHandler, oomph::FoldHandler, oomph::ParameterDerivativeHandler, oomph::ParallelResidualsHandler, oomph::EigenProblemHandler, and oomph::ExplicitTimeStepHandler.
Definition at line 102 of file assembly_handler.cc.
References oomph::GeneralisedElement::get_jacobian().
Referenced by get_all_vectors_and_matrices(), oomph::Problem::get_hessian_vector_products(), oomph::Problem::get_jacobian(), oomph::ParallelResidualsHandler::get_jacobian(), oomph::HSL_MA42::solve(), and oomph::HSL_MA42::solve_for_one_dof().
|
virtual |
Return the contribution to the residuals of the element elem_pt.
Get the residuals by calling the underlying element's residuals directly.
Reimplemented in oomph::FpPreconditionerAssemblyHandler, oomph::PeriodicOrbitAssemblyHandler< NNODE_1D >, oomph::HopfHandler, oomph::PitchForkHandler, oomph::FoldHandler, oomph::ParameterDerivativeHandler, oomph::ParallelResidualsHandler, oomph::EigenProblemHandler, and oomph::ExplicitTimeStepHandler.
Definition at line 92 of file assembly_handler.cc.
References oomph::GeneralisedElement::get_residuals().
Referenced by oomph::ParallelResidualsHandler::get_all_vectors_and_matrices(), oomph::Problem::get_residuals(), and oomph::ParallelResidualsHandler::get_residuals().
|
virtual |
Return the t-th level of storage associated with the i-th (local) dof stored in the problem.
Definition at line 70 of file assembly_handler.cc.
References oomph::Problem::dof_pt(), i, and t.
|
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 in oomph::PeriodicOrbitAssemblyHandler< NNODE_1D >, oomph::HopfHandler, oomph::PitchForkHandler, oomph::FoldHandler, oomph::ParameterDerivativeHandler, oomph::ParallelResidualsHandler, oomph::EigenProblemHandler, and oomph::ExplicitTimeStepHandler.
Definition at line 42 of file assembly_handler.cc.
References oomph::GeneralisedElement::ndof().
Referenced by oomph::Problem::get_hessian_vector_products(), oomph::Problem::get_jacobian(), oomph::Problem::get_my_eqns(), oomph::Problem::get_residuals(), oomph::ParallelResidualsHandler::ndof(), oomph::ParameterDerivativeHandler::ndof(), oomph::Problem::parallel_sparse_assemble(), oomph::PitchForkHandler::PitchForkHandler(), oomph::HSL_MA42::reorder_elements(), oomph::Problem::setup_element_count_per_dof(), oomph::HSL_MA42::solve(), oomph::HSL_MA42::solve_for_one_dof(), 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().
|
inlinevirtual |
Function that is used to perform any synchronisation required during the solution.
Reimplemented in oomph::PitchForkHandler.
Definition at line 165 of file assembly_handler.h.
Referenced by oomph::Problem::synchronise_all_dofs().