Public Member Functions | Private Attributes | Friends | List of all members
oomph::HopfHandler Class Reference

#include <assembly_handler.h>

+ Inheritance diagram for oomph::HopfHandler:

Public Member Functions

 HopfHandler (Problem *const &problem_pt, double *const &parameter_pt)
 Constructor. More...
 
 HopfHandler (Problem *const &problem_pt, double *const &paramter_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 &parameter_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 &parameter_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...
 
- 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 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...
 
ProblemProblem_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
 

Detailed Description

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 $ \lambda $ and a solution is $ R(u,\lambda) = 0 $, where $ u $ are the unknowns in the problem. A Hopf bifurcation may be specified by the augmented system of size $ 3N+2 $.

\[ R(u,\lambda) = 0, \]

\[ J\phi + \omega M \psi = 0, \]

\[ J\psi - \omega M \phi = 0, \]

\[ c \cdot \phi = 1. \]

\[ c \cdot \psi = 0. \]

In the above $ J $ is the usual Jacobian matrix, $ dR/du $ and $ M $ is the mass matrix that multiplies the time derivative terms. $ \phi + i\psi $ is the (complex) null vector of the complex matrix $ J - i\omega M $, where $ \omega $ is the critical frequency. $ c $ 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.

Constructor & Destructor Documentation

◆ HopfHandler() [1/2]

oomph::HopfHandler::HopfHandler ( Problem *const &  problem_pt,
double *const &  parameter_pt 
)

◆ HopfHandler() [2/2]

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.

◆ ~HopfHandler()

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.

Member Function Documentation

◆ bifurcation_parameter_pt()

double* oomph::HopfHandler::bifurcation_parameter_pt ( ) const
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().

◆ bifurcation_type()

int oomph::HopfHandler::bifurcation_type ( ) const
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.

◆ eqn_number()

unsigned long oomph::HopfHandler::eqn_number ( GeneralisedElement *const &  elem_pt,
const unsigned &  ieqn_local 
)
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().

◆ get_djacobian_dparameter()

void oomph::HopfHandler::get_djacobian_dparameter ( GeneralisedElement *const &  elem_pt,
double *const &  parameter_pt,
Vector< double > &  dres_dparam,
DenseMatrix< double > &  djac_dparam 
)
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.

◆ get_dresiduals_dparameter()

void oomph::HopfHandler::get_dresiduals_dparameter ( GeneralisedElement *const &  elem_pt,
double *const &  parameter_pt,
Vector< double > &  dres_dparam 
)
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.

◆ get_eigenfunction()

void oomph::HopfHandler::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

Reimplemented from oomph::AssemblyHandler.

Definition at line 4826 of file assembly_handler.cc.

References oomph::Problem::communicator_pt(), Ndof, Phi, Problem_pt, and Psi.

◆ get_hessian_vector_products()

void oomph::HopfHandler::get_hessian_vector_products ( GeneralisedElement *const &  elem_pt,
Vector< double > const &  Y,
DenseMatrix< double > const &  C,
DenseMatrix< double > &  product 
)
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.

◆ get_jacobian()

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

◆ get_residuals()

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

◆ ndof()

unsigned oomph::HopfHandler::ndof ( GeneralisedElement *const &  elem_pt)
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().

◆ omega()

const double& oomph::HopfHandler::omega ( ) const
inline

Return the frequency of the bifurcation.

Definition at line 1116 of file assembly_handler.h.

◆ solve_complex_system()

void oomph::HopfHandler::solve_complex_system ( )

◆ solve_full_system()

void oomph::HopfHandler::solve_full_system ( )

◆ solve_standard_system()

void oomph::HopfHandler::solve_standard_system ( )

Friends And Related Function Documentation

◆ BlockHopfLinearSolver

friend class BlockHopfLinearSolver
friend

Definition at line 1013 of file assembly_handler.h.

Member Data Documentation

◆ C

Vector<double> oomph::HopfHandler::C
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().

◆ Count

Vector<int> oomph::HopfHandler::Count
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().

◆ Ndof

unsigned oomph::HopfHandler::Ndof
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().

◆ Omega

double oomph::HopfHandler::Omega
private

◆ Parameter_pt

double* oomph::HopfHandler::Parameter_pt
private

Pointer to the parameter.

Definition at line 1026 of file assembly_handler.h.

Referenced by solve_full_system().

◆ Phi

Vector<double> oomph::HopfHandler::Phi
private

◆ Problem_pt

Problem* oomph::HopfHandler::Problem_pt
private

◆ Psi

Vector<double> oomph::HopfHandler::Psi
private

◆ Solve_which_system

unsigned oomph::HopfHandler::Solve_which_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().


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