A class that is used to assemble the augmented system that defines a fold (saddle-node) or limit point. The "standard" problem must be a function of a global paramter , and a solution is , where are the unknowns in the problem. A limit point is formally specified by the augmented system of size . More...
#include <assembly_handler.h>
Public Member Functions | |
FoldHandler (Problem *const &problem_pt, double *const ¶meter_pt) | |
Constructor: initialise the fold handler, by setting initial guesses for Y, Phi and calculating count. If the system changes, a new fold handler must be constructed. More... | |
FoldHandler (Problem *const &problem_pt, double *const ¶meter_pt, const DoubleVector &eigenvector) | |
Constructor in which initial eigenvector can be passed. More... | |
FoldHandler (Problem *const &problem_pt, double *const ¶meter_pt, const DoubleVector &eigenvector, const DoubleVector &normalisation) | |
Constructor in which initial eigenvector and normalisation can be passed. More... | |
~FoldHandler () | |
Destructor, return the problem to its original state before the augmented system was added. More... | |
unsigned | ndof (GeneralisedElement *const &elem_pt) |
Get the number of elemental degrees of freedom. More... | |
unsigned long | eqn_number (GeneralisedElement *const &elem_pt, const unsigned &ieqn_local) |
Get the global equation number of the local unknown. More... | |
void | get_residuals (GeneralisedElement *const &elem_pt, Vector< double > &residuals) |
Get the residuals. More... | |
void | get_jacobian (GeneralisedElement *const &elem_pt, Vector< double > &residuals, DenseMatrix< double > &jacobian) |
Calculate the elemental Jacobian matrix "d equation
/ d variable". More... | |
void | get_dresiduals_dparameter (GeneralisedElement *const &elem_pt, double *const ¶meter_pt, Vector< double > &dres_dparam) |
Overload the derivatives of the residuals with respect to a parameter to apply to the augmented system. More... | |
void | get_djacobian_dparameter (GeneralisedElement *const &elem_pt, double *const ¶meter_pt, Vector< double > &dres_dparam, DenseMatrix< double > &djac_dparam) |
Overload the derivative of the residuals and jacobian with respect to a parameter so that it breaks. More... | |
void | get_hessian_vector_products (GeneralisedElement *const &elem_pt, Vector< double > const &Y, DenseMatrix< double > const &C, DenseMatrix< double > &product) |
Overload the hessian vector product function so that it breaks. More... | |
int | bifurcation_type () const |
Indicate that we are tracking a fold bifurcation by returning 1. More... | |
double * | bifurcation_parameter_pt () const |
Return a pointer to the bifurcation parameter in bifurcation tracking problems. More... | |
void | get_eigenfunction (Vector< DoubleVector > &eigenfunction) |
Return the eigenfunction(s) associated with the bifurcation that has been detected in bifurcation tracking problems. More... | |
void | solve_augmented_block_system () |
Set to solve the augmented block system. More... | |
void | solve_block_system () |
Set to solve the block system. More... | |
void | solve_full_system () |
Solve non-block system. 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_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_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... | |
Private Types | |
enum | { Full_augmented , Block_J , Block_augmented_J } |
A little private enum to determine whether we are solving the block system or not. More... | |
Private Attributes | |
unsigned | Solve_which_system |
Integer flag to indicate which system should be assembled. There are three possibilities. The full augmented system (0), the non-augmented jacobian system (1), a system in which the jacobian is augmented by 1 row and column to ensure that it is non-singular (2). See the enum above. More... | |
Problem * | Problem_pt |
Pointer to the problem. More... | |
unsigned | Ndof |
Store the number of degrees of freedom in the non-augmented problem. More... | |
Vector< double > | Phi |
A constant vector used to ensure that the null vector is not trivial. More... | |
Vector< double > | Y |
Storage for the null vector. More... | |
Vector< int > | Count |
A vector that is used to determine how many elements contribute to a particular equation. It is used to ensure that the global system is correctly formulated. More... | |
double * | Parameter_pt |
Storage for the pointer to the parameter. More... | |
Friends | |
class | AugmentedBlockFoldLinearSolver |
A class that is used to assemble the augmented system that defines a fold (saddle-node) or limit point. The "standard" problem must be a function of a global paramter , and a solution is , where are the unknowns in the problem. A limit point is formally specified by the augmented system of size .
In the above is the usual Jacobian matrix, , and is a constant vector that is chosen to ensure that the null vector, , is not trivial.
Definition at line 472 of file assembly_handler.h.
|
private |
A little private enum to determine whether we are solving the block system or not.
Enumerator | |
---|---|
Full_augmented | |
Block_J | |
Block_augmented_J |
Definition at line 478 of file assembly_handler.h.
oomph::FoldHandler::FoldHandler | ( | Problem *const & | problem_pt, |
double *const & | parameter_pt | ||
) |
Constructor: initialise the fold handler, by setting initial guesses for Y, Phi and calculating count. If the system changes, a new fold handler must be constructed.
Constructor: Initialise the fold handler, by setting initial guesses for Y, Phi and calculating Count. If the system changes, a new handler must be constructed.
Definition at line 864 of file assembly_handler.cc.
References oomph::LinearAlgebraDistribution::build(), oomph::Problem::communicator_pt(), Count, oomph::LinearSolver::disable_resolve(), oomph::Problem::Dof_distribution_pt, oomph::Problem::Dof_pt, e, oomph::Mesh::element_pt(), oomph::LinearSolver::enable_resolve(), oomph::GeneralisedElement::eqn_number(), oomph::Problem::get_derivative_wrt_global_parameter(), oomph::LinearSolver::is_resolve_enabled(), oomph::Problem::linear_solver_pt(), oomph::Problem::mesh_pt(), Ndof, oomph::GeneralisedElement::ndof(), oomph::Problem::ndof(), oomph::Mesh::nelement(), Phi, Problem_pt, oomph::LinearSolver::resolve(), oomph::LinearSolver::solve(), oomph::Problem::Sparse_assemble_with_arrays_previous_allocation, and Y.
oomph::FoldHandler::FoldHandler | ( | Problem *const & | problem_pt, |
double *const & | parameter_pt, | ||
const DoubleVector & | eigenvector | ||
) |
Constructor in which initial eigenvector can be passed.
Constructor in which the eigenvector is passed as an initial guess.
Definition at line 970 of file assembly_handler.cc.
References oomph::LinearAlgebraDistribution::build(), oomph::Problem::communicator_pt(), Count, oomph::Problem::Dof_distribution_pt, oomph::Problem::Dof_pt, e, oomph::Mesh::element_pt(), oomph::GeneralisedElement::eqn_number(), oomph::Problem::mesh_pt(), Ndof, oomph::GeneralisedElement::ndof(), oomph::Problem::ndof(), oomph::Mesh::nelement(), Phi, Problem_pt, oomph::Problem::Sparse_assemble_with_arrays_previous_allocation, and Y.
oomph::FoldHandler::FoldHandler | ( | Problem *const & | problem_pt, |
double *const & | parameter_pt, | ||
const DoubleVector & | eigenvector, | ||
const DoubleVector & | normalisation | ||
) |
Constructor in which initial eigenvector and normalisation can be passed.
Constructor in which the eigenvector and normalisation vector are passed as an initial guess.
Definition at line 1036 of file assembly_handler.cc.
References oomph::LinearAlgebraDistribution::build(), oomph::Problem::communicator_pt(), Count, oomph::Problem::Dof_distribution_pt, oomph::Problem::Dof_pt, e, oomph::Mesh::element_pt(), oomph::GeneralisedElement::eqn_number(), oomph::Problem::mesh_pt(), Ndof, oomph::GeneralisedElement::ndof(), oomph::Problem::ndof(), oomph::Mesh::nelement(), Phi, Problem_pt, oomph::Problem::Sparse_assemble_with_arrays_previous_allocation, and Y.
oomph::FoldHandler::~FoldHandler | ( | ) |
Destructor, return the problem to its original state before the augmented system was added.
Destructor return the problem to its original (non-augmented) state.
Definition at line 1573 of file assembly_handler.cc.
References oomph::LinearAlgebraDistribution::build(), oomph::Problem::communicator_pt(), oomph::Problem::Dof_distribution_pt, oomph::Problem::Dof_pt, oomph::Problem::linear_solver_pt(), oomph::AugmentedBlockFoldLinearSolver::linear_solver_pt(), Ndof, Problem_pt, and oomph::Problem::Sparse_assemble_with_arrays_previous_allocation.
|
inlinevirtual |
Return a pointer to the bifurcation parameter in bifurcation tracking problems.
Reimplemented from oomph::AssemblyHandler.
Definition at line 584 of file assembly_handler.h.
References Parameter_pt.
|
inlinevirtual |
Indicate that we are tracking a fold bifurcation by returning 1.
Reimplemented from oomph::AssemblyHandler.
Definition at line 577 of file assembly_handler.h.
|
virtual |
Get the global equation number of the local unknown.
Return the global equation number associated with local equation number ieqn_local. We choose to number the unknowns according to the augmented system.
Reimplemented from oomph::AssemblyHandler.
Definition at line 1138 of file assembly_handler.cc.
References oomph::GeneralisedElement::eqn_number(), Ndof, and oomph::GeneralisedElement::ndof().
Referenced by get_jacobian(), oomph::AugmentedBlockFoldLinearSolver::resolve(), and oomph::AugmentedBlockFoldLinearSolver::solve().
|
virtual |
Overload the derivative of the residuals and jacobian with respect to a parameter so that it breaks.
Overload the derivative of the residuals and jacobian with respect to a parameter so that it breaks because it should not be required.
Reimplemented from oomph::AssemblyHandler.
Definition at line 1512 of file assembly_handler.cc.
|
virtual |
Overload the derivatives of the residuals with respect to a parameter to apply to the augmented system.
Formulate the derivatives of the augmented system with respect to a parameter.
Reimplemented from oomph::AssemblyHandler.
Definition at line 1442 of file assembly_handler.cc.
References Block_augmented_J, Block_J, oomph::GeneralisedElement::eqn_number(), Full_augmented, oomph::GeneralisedElement::get_djacobian_dparameter(), oomph::GeneralisedElement::get_dresiduals_dparameter(), i, oomph::GeneralisedElement::ndof(), Solve_which_system, and Y.
|
virtual |
Return the eigenfunction(s) associated with the bifurcation that has been detected in bifurcation tracking problems.
Reimplemented from oomph::AssemblyHandler.
Definition at line 1555 of file assembly_handler.cc.
References oomph::Problem::communicator_pt(), Ndof, Problem_pt, and Y.
|
virtual |
Overload the hessian vector product function so that it breaks.
Overload the hessian vector product function so that it breaks because it should not be required.
Reimplemented from oomph::AssemblyHandler.
Definition at line 1533 of file assembly_handler.cc.
|
virtual |
Calculate the elemental Jacobian matrix "d equation / d variable".
Calculate the elemental Jacobian matrix "d equation / d variable" for the augmented system.
ALICE
ALICE
Reimplemented from oomph::AssemblyHandler.
Definition at line 1244 of file assembly_handler.cc.
References oomph::Problem::actions_after_change_in_bifurcation_parameter(), oomph::Problem::actions_before_newton_convergence_check(), Block_augmented_J, Block_J, Count, oomph::Problem::Dof_pt, oomph::GeneralisedElement::eqn_number(), eqn_number(), oomph::BlackBoxFDNewtonSolver::FD_step, Full_augmented, oomph::GeneralisedElement::get_jacobian(), get_residuals(), Ndof, oomph::GeneralisedElement::ndof(), ndof(), Phi, Problem_pt, Solve_which_system, and Y.
Referenced by oomph::AugmentedBlockFoldLinearSolver::resolve(), and oomph::AugmentedBlockFoldLinearSolver::solve().
|
virtual |
Get the residuals.
Formulate the augmented system.
Reimplemented from oomph::AssemblyHandler.
Definition at line 1171 of file assembly_handler.cc.
References Block_augmented_J, Block_J, Count, oomph::GeneralisedElement::eqn_number(), Full_augmented, oomph::GeneralisedElement::get_jacobian(), oomph::GeneralisedElement::get_residuals(), i, oomph::Problem::mesh_pt(), oomph::GeneralisedElement::ndof(), oomph::Mesh::nelement(), Phi, Problem_pt, Solve_which_system, and Y.
Referenced by get_jacobian().
|
virtual |
Get the number of elemental degrees of freedom.
The number of unknowns is 2n+1.
Reimplemented from oomph::AssemblyHandler.
Definition at line 1105 of file assembly_handler.cc.
References Block_augmented_J, Block_J, Full_augmented, oomph::GeneralisedElement::ndof(), and Solve_which_system.
Referenced by get_jacobian(), oomph::AugmentedBlockFoldLinearSolver::resolve(), and oomph::AugmentedBlockFoldLinearSolver::solve().
void oomph::FoldHandler::solve_augmented_block_system | ( | ) |
Set to solve the augmented block system.
Set to solve the augmented-by-one block system.
Definition at line 1600 of file assembly_handler.cc.
References Block_augmented_J, Block_J, oomph::LinearAlgebraDistribution::build(), oomph::Problem::communicator_pt(), oomph::Problem::Dof_distribution_pt, oomph::Problem::Dof_pt, Ndof, Parameter_pt, Problem_pt, Solve_which_system, and oomph::Problem::Sparse_assemble_with_arrays_previous_allocation.
Referenced by oomph::AugmentedBlockFoldLinearSolver::resolve(), and oomph::AugmentedBlockFoldLinearSolver::solve().
void oomph::FoldHandler::solve_block_system | ( | ) |
Set to solve the block system.
Set to solve the block system. The boolean flag specifies whether the block decomposition should use exactly the same jacobian.
Definition at line 1629 of file assembly_handler.cc.
References Block_J, oomph::LinearAlgebraDistribution::build(), oomph::Problem::communicator_pt(), oomph::Problem::Dof_distribution_pt, oomph::Problem::Dof_pt, Ndof, Problem_pt, Solve_which_system, and oomph::Problem::Sparse_assemble_with_arrays_previous_allocation.
void oomph::FoldHandler::solve_full_system | ( | ) |
Solve non-block system.
Set to Solve non-block system.
Definition at line 1649 of file assembly_handler.cc.
References Block_J, oomph::LinearAlgebraDistribution::build(), oomph::Problem::communicator_pt(), oomph::Problem::Dof_distribution_pt, oomph::Problem::Dof_pt, Full_augmented, Ndof, Parameter_pt, Problem_pt, Solve_which_system, oomph::Problem::Sparse_assemble_with_arrays_previous_allocation, and Y.
Referenced by oomph::AugmentedBlockFoldLinearSolver::resolve(), and oomph::AugmentedBlockFoldLinearSolver::solve().
|
friend |
Definition at line 474 of file assembly_handler.h.
|
private |
A vector that is used to determine how many elements contribute to a particular equation. It is used to ensure that the global system is correctly formulated.
Definition at line 509 of file assembly_handler.h.
Referenced by FoldHandler(), get_jacobian(), and get_residuals().
|
private |
Store the number of degrees of freedom in the non-augmented problem.
Definition at line 497 of file assembly_handler.h.
Referenced by eqn_number(), FoldHandler(), get_eigenfunction(), get_jacobian(), solve_augmented_block_system(), solve_block_system(), solve_full_system(), and ~FoldHandler().
|
private |
Storage for the pointer to the parameter.
Definition at line 512 of file assembly_handler.h.
Referenced by bifurcation_parameter_pt(), solve_augmented_block_system(), and solve_full_system().
|
private |
A constant vector used to ensure that the null vector is not trivial.
Definition at line 501 of file assembly_handler.h.
Referenced by FoldHandler(), get_jacobian(), and get_residuals().
|
private |
Pointer to the problem.
Definition at line 493 of file assembly_handler.h.
Referenced by FoldHandler(), get_eigenfunction(), get_jacobian(), get_residuals(), solve_augmented_block_system(), solve_block_system(), solve_full_system(), and ~FoldHandler().
|
private |
Integer flag to indicate which system should be assembled. There are three possibilities. The full augmented system (0), the non-augmented jacobian system (1), a system in which the jacobian is augmented by 1 row and column to ensure that it is non-singular (2). See the enum above.
Definition at line 490 of file assembly_handler.h.
Referenced by get_dresiduals_dparameter(), get_jacobian(), get_residuals(), ndof(), solve_augmented_block_system(), solve_block_system(), and solve_full_system().
|
private |
Storage for the null vector.
Definition at line 504 of file assembly_handler.h.
Referenced by FoldHandler(), get_dresiduals_dparameter(), get_eigenfunction(), get_jacobian(), get_residuals(), oomph::AugmentedBlockFoldLinearSolver::resolve(), oomph::AugmentedBlockFoldLinearSolver::solve(), and solve_full_system().