#include <assembly_handler.h>
Public Member Functions | |
HopfHandler (Problem *const &problem_pt, double *const ¶meter_pt) | |
Constructor. More... | |
HopfHandler (Problem *const &problem_pt, double *const ¶mter_pt, const double &omega, const DoubleVector &phi, const DoubleVector &psi) | |
~HopfHandler () | |
Destructor return the problem to its original (non-augmented) state. 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 Hopf bifurcation by returning 3. 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... | |
const double & | omega () const |
Return the frequency of the bifurcation. More... | |
void | solve_standard_system () |
Set to solve the standard system. More... | |
void | solve_complex_system () |
Set to solve the complex system. More... | |
void | solve_full_system () |
Solve non-block system. More... | |
![]() | |
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 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), and complex system (2), where the matrix is a combination of the jacobian and mass matrices. More... | |
Problem * | Problem_pt |
Pointer to the problem. More... | |
double * | Parameter_pt |
Pointer to the parameter. More... | |
unsigned | Ndof |
Store the number of degrees of freedom in the non-augmented problem. More... | |
double | Omega |
The critical frequency of the bifurcation. More... | |
Vector< double > | Phi |
The real part of the null vector. More... | |
Vector< double > | Psi |
The imaginary part of the null vector. More... | |
Vector< double > | C |
A constant vector used to ensure that the null vector is not trivial. 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... | |
Friends | |
class | BlockHopfLinearSolver |
A class that is used to assemble the augmented system that defines a Hopf bifurcation. The "standard" problem must be a function of a global parameter and a solution is
, where
are the unknowns in the problem. A Hopf bifurcation may be specified by the augmented system of size
.
In the above is the usual Jacobian matrix,
and
is the mass matrix that multiplies the time derivative terms.
is the (complex) null vector of the complex matrix
, where
is the critical frequency.
is a constant vector that is used to ensure that the null vector is non-trivial.
Definition at line 1011 of file assembly_handler.h.
oomph::HopfHandler::HopfHandler | ( | Problem *const & | problem_pt, |
double *const & | parameter_pt | ||
) |
Constructor.
Constructor: Initialise the hopf handler, by setting initial guesses for Phi, Psi and calculating Count. If the system changes, a new handler must be constructed.
Definition at line 4227 of file assembly_handler.cc.
References oomph::LinearAlgebraDistribution::build(), C, 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(), oomph::GeneralisedElement::ndof(), Ndof, oomph::Problem::ndof(), oomph::Mesh::nelement(), Omega, Phi, Problem_pt, Psi, oomph::LinearSolver::resolve(), oomph::LinearSolver::solve(), and oomph::Problem::Sparse_assemble_with_arrays_previous_allocation.
oomph::HopfHandler::HopfHandler | ( | Problem *const & | problem_pt, |
double *const & | parameter_pt, | ||
const double & | omega, | ||
const DoubleVector & | phi, | ||
const DoubleVector & | psi | ||
) |
Constructor with initial guesses for the frequency and null vectors, such as might be provided by an eigensolver
Constructor: Initialise the hopf handler, by setting initial guesses for Phi, Psi, Omega and calculating Count. If the system changes, a new handler must be constructed.
Definition at line 4358 of file assembly_handler.cc.
References oomph::LinearAlgebraDistribution::build(), C, 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(), oomph::GeneralisedElement::ndof(), Ndof, oomph::Problem::ndof(), oomph::Mesh::nelement(), Omega, Phi, Problem_pt, Psi, and oomph::Problem::Sparse_assemble_with_arrays_previous_allocation.
oomph::HopfHandler::~HopfHandler | ( | ) |
Destructor return the problem to its original (non-augmented) state.
Destructor, return the problem to its original state, before the augmented system was added
Definition at line 4427 of file assembly_handler.cc.
References oomph::LinearAlgebraDistribution::build(), oomph::Problem::communicator_pt(), oomph::Problem::Dof_distribution_pt, oomph::Problem::Dof_pt, oomph::BlockHopfLinearSolver::linear_solver_pt(), oomph::Problem::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 1108 of file assembly_handler.h.
References oomph::AssemblyHandler::get_eigenfunction().
|
inlinevirtual |
Indicate that we are tracking a Hopf bifurcation by returning 3.
Reimplemented from oomph::AssemblyHandler.
Definition at line 1104 of file assembly_handler.h.
|
virtual |
Get the global equation number of the local unknown.
Reimplemented from oomph::AssemblyHandler.
Definition at line 4480 of file assembly_handler.cc.
References oomph::GeneralisedElement::eqn_number(), oomph::GeneralisedElement::ndof(), and Ndof.
Referenced by get_jacobian().
|
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 4780 of file assembly_handler.cc.
|
virtual |
Overload the derivatives of the residuals with respect to a parameter to apply to the augmented system.
Get the derivatives of the augmented residuals with respect to a parameter
Reimplemented from oomph::AssemblyHandler.
Definition at line 4724 of file assembly_handler.cc.
References oomph::GeneralisedElement::eqn_number(), oomph::GeneralisedElement::get_djacobian_and_dmass_matrix_dparameter(), i, oomph::GeneralisedElement::ndof(), Omega, Phi, Psi, and Solve_which_system.
|
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
Reimplemented from oomph::AssemblyHandler.
Definition at line 4826 of file assembly_handler.cc.
References oomph::Problem::communicator_pt(), Ndof, Phi, Problem_pt, and Psi.
|
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 4803 of file assembly_handler.cc.
|
virtual |
Calculate the elemental Jacobian matrix "d equation / d variable".
Reimplemented from oomph::AssemblyHandler.
Definition at line 4570 of file assembly_handler.cc.
References oomph::Problem::actions_after_change_in_bifurcation_parameter(), C, Count, oomph::Problem::Dof_pt, oomph::GeneralisedElement::eqn_number(), eqn_number(), oomph::BlackBoxFDNewtonSolver::FD_step, oomph::GeneralisedElement::get_jacobian(), oomph::GeneralisedElement::get_jacobian_and_mass_matrix(), get_residuals(), oomph::GeneralisedElement::ndof(), Ndof, ndof(), Omega, Phi, Problem_pt, Psi, and Solve_which_system.
|
virtual |
Get the residuals.
Reimplemented from oomph::AssemblyHandler.
Definition at line 4512 of file assembly_handler.cc.
References C, Count, oomph::GeneralisedElement::eqn_number(), oomph::GeneralisedElement::get_jacobian_and_mass_matrix(), i, oomph::Problem::mesh_pt(), oomph::GeneralisedElement::ndof(), oomph::Mesh::nelement(), Omega, Phi, Problem_pt, Psi, and Solve_which_system.
Referenced by get_jacobian().
|
virtual |
Get the number of elemental degrees of freedom.
Reimplemented from oomph::AssemblyHandler.
Definition at line 4452 of file assembly_handler.cc.
References oomph::GeneralisedElement::ndof(), and Solve_which_system.
Referenced by get_jacobian().
|
inline |
Return the frequency of the bifurcation.
Definition at line 1116 of file assembly_handler.h.
void oomph::HopfHandler::solve_complex_system | ( | ) |
Set to solve the complex system.
Set to solve the complex (jacobian and mass matrix) system.
Definition at line 4864 of file assembly_handler.cc.
References oomph::LinearAlgebraDistribution::build(), oomph::Problem::communicator_pt(), oomph::Problem::Dof_distribution_pt, oomph::Problem::Dof_pt, Ndof, Phi, Problem_pt, Solve_which_system, and oomph::Problem::Sparse_assemble_with_arrays_previous_allocation.
Referenced by oomph::BlockHopfLinearSolver::solve(), and oomph::BlockHopfLinearSolver::solve_for_two_rhs().
void oomph::HopfHandler::solve_full_system | ( | ) |
Solve non-block system.
Set to Solve full system system.
Definition at line 4890 of file assembly_handler.cc.
References oomph::LinearAlgebraDistribution::build(), oomph::Problem::communicator_pt(), oomph::Problem::Dof_distribution_pt, oomph::Problem::Dof_pt, Ndof, Omega, Parameter_pt, Phi, Problem_pt, Psi, Solve_which_system, and oomph::Problem::Sparse_assemble_with_arrays_previous_allocation.
Referenced by oomph::BlockHopfLinearSolver::solve(), and oomph::BlockHopfLinearSolver::solve_for_two_rhs().
void oomph::HopfHandler::solve_standard_system | ( | ) |
Set to solve the standard system.
Set to solve the standard (underlying jacobian) system.
Definition at line 4847 of file assembly_handler.cc.
References 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.
Referenced by oomph::BlockHopfLinearSolver::solve(), and oomph::BlockHopfLinearSolver::solve_for_two_rhs().
|
friend |
Definition at line 1013 of file assembly_handler.h.
|
private |
A constant vector used to ensure that the null vector is not trivial.
Definition at line 1043 of file assembly_handler.h.
Referenced by get_jacobian(), get_residuals(), HopfHandler(), oomph::BlockHopfLinearSolver::solve(), and oomph::BlockHopfLinearSolver::solve_for_two_rhs().
|
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 1048 of file assembly_handler.h.
Referenced by get_jacobian(), get_residuals(), and HopfHandler().
|
private |
Store the number of degrees of freedom in the non-augmented problem.
Definition at line 1030 of file assembly_handler.h.
Referenced by eqn_number(), get_eigenfunction(), get_jacobian(), HopfHandler(), solve_complex_system(), solve_full_system(), solve_standard_system(), and ~HopfHandler().
|
private |
The critical frequency of the bifurcation.
Definition at line 1033 of file assembly_handler.h.
Referenced by get_dresiduals_dparameter(), get_jacobian(), get_residuals(), HopfHandler(), oomph::BlockHopfLinearSolver::solve(), oomph::BlockHopfLinearSolver::solve_for_two_rhs(), and solve_full_system().
|
private |
Pointer to the parameter.
Definition at line 1026 of file assembly_handler.h.
Referenced by solve_full_system().
|
private |
The real part of the null vector.
Definition at line 1036 of file assembly_handler.h.
Referenced by get_dresiduals_dparameter(), get_eigenfunction(), get_jacobian(), get_residuals(), HopfHandler(), oomph::BlockHopfLinearSolver::solve(), solve_complex_system(), oomph::BlockHopfLinearSolver::solve_for_two_rhs(), and solve_full_system().
|
private |
Pointer to the problem.
Definition at line 1023 of file assembly_handler.h.
Referenced by get_eigenfunction(), get_jacobian(), get_residuals(), HopfHandler(), solve_complex_system(), solve_full_system(), solve_standard_system(), and ~HopfHandler().
|
private |
The imaginary part of the null vector.
Definition at line 1039 of file assembly_handler.h.
Referenced by get_dresiduals_dparameter(), get_eigenfunction(), get_jacobian(), get_residuals(), HopfHandler(), oomph::BlockHopfLinearSolver::solve(), oomph::BlockHopfLinearSolver::solve_for_two_rhs(), and solve_full_system().
|
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), and complex system (2), where the matrix is a combination of the jacobian and mass matrices.
Definition at line 1020 of file assembly_handler.h.
Referenced by get_dresiduals_dparameter(), get_jacobian(), get_residuals(), ndof(), solve_complex_system(), solve_full_system(), and solve_standard_system().