CppNoddy  0.92
Loading...
Searching...
No Matches
Public Member Functions | List of all members
CppNoddy::ODE_EVP< _Type > Class Template Reference

A templated object for real/complex vector system of first-order ordinary differential equations. More...

#include <ODE_EVP.h>

Inheritance diagram for CppNoddy::ODE_EVP< _Type >:
CppNoddy::Uncopyable

Public Member Functions

 ODE_EVP (Equation_2matrix< _Type > *equation_ptr, const DenseVector< double > &nodes, Residual< _Type > *ptr_to_left_residual, Residual< _Type > *ptr_to_right_residual)
 The class is defined by a vector function for the system. More...
 
 ~ODE_EVP ()
 Destructor. More...
 
void eigensolve ()
 Formulate and solve the global eigenvalue problem for a linear system. More...
 
LinearEigenSystem_basep_eigensystem ()
 Allow access to the underlying dense linear eigensystem through a pointer to the private member data. More...
 
void add_tagged_to_mesh ()
 
OneD_Node_Mesh< D_complexget_mesh (const unsigned &i) const
 

Detailed Description

template<typename _Type>
class CppNoddy::ODE_EVP< _Type >

A templated object for real/complex vector system of first-order ordinary differential equations.

Definition at line 35 of file ODE_EVP.h.

Constructor & Destructor Documentation

◆ ODE_EVP()

template<typename _Type >
CppNoddy::ODE_EVP< _Type >::ODE_EVP ( Equation_2matrix< _Type > *  equation_ptr,
const DenseVector< double > &  nodes,
Residual< _Type > *  ptr_to_left_residual,
Residual< _Type > *  ptr_to_right_residual 
)

The class is defined by a vector function for the system.

Parameters
equation_ptrA pointer to an equation with 2 associated matrices; matrix1 will define the eigenvalue problem.
nodesA vector of nodal points.
ptr_to_left_residualA pointer to a residual object that defines the LHS boundary conditions.
ptr_to_right_residualA pointer to a residual object that defines the RHS boundary conditions.

Definition at line 26 of file ODE_EVP.cpp.

29 :
30 p_EQUATION(ptr_to_equation),
31 p_LEFT_RESIDUAL(ptr_to_left_residual),
32 p_RIGHT_RESIDUAL(ptr_to_right_residual),
33 NODES(nodes) {
34 unsigned n(nodes.size());
35 unsigned order(p_EQUATION -> get_order());
36 // first make the dense matrix equivalent of the problem
37 // using the Dense( banded ) constructor
38 p_A_DENSE = new DenseMatrix<_Type>(n * order, n * order, 0.0);
39 p_B_DENSE = new DenseMatrix<_Type>(n * order, n * order, 0.0);
40 // point the system pointer to the dense problem
41 p_SYSTEM = new DenseLinearEigenSystem<_Type>(p_A_DENSE, p_B_DENSE);
42 }
std::size_t size() const
A pass-thru definition to get the size of the vector.
Definition: DenseVector.h:330

References CppNoddy::DenseVector< _Type >::size().

◆ ~ODE_EVP()

template<typename _Type >
CppNoddy::ODE_EVP< _Type >::~ODE_EVP

Destructor.

Definition at line 45 of file ODE_EVP.cpp.

45 {
46 // clean up the matrices & evp system
47 delete p_SYSTEM;
48 delete p_A_DENSE;
49 delete p_B_DENSE;
50 }

Member Function Documentation

◆ add_tagged_to_mesh()

template<typename _Type >
void CppNoddy::ODE_EVP< _Type >::add_tagged_to_mesh ( )
inline

Definition at line 60 of file ODE_EVP.h.

60 {
61 // clear the existing data if any
62 MESHES.clear();
63 // order of the equation
64 unsigned order = p_EQUATION -> get_order();
65 // get the eigenvalues
66 DenseVector<D_complex> vals(p_SYSTEM -> get_tagged_eigenvalues());
67 // get the eigenvectors
68 DenseMatrix<D_complex> vecs(p_SYSTEM -> get_tagged_eigenvectors());
69 // loop through the eigenvectors
70 for(unsigned ivec = 0; ivec < vals.size(); ++ivec) {
71 // make a mesh with the right node distribution
72 // we'll increase the order by 1 to allow the eigenvalue to be
73 // stored at each nodal point -- this is wasteful but very useful
74 // for feeding this into a BVP for local refinement.
75 OneD_Node_Mesh<D_complex> eigfn(NODES, order + 1);
76 // loop through all nodes
77 for(unsigned node = 0; node < NODES.size(); ++node) {
78 // complex vector of the dof at this node ( + 1 for the eigenvalue)
79 DenseVector<D_complex> vars_at_node(order + 1, 0.0);
80 // get the dof from the eigenvector
81 for(unsigned var = 0; var < order; ++var) {
82 vars_at_node[ var ] = vecs[ ivec ][ node * order + var ];
83 }
84 // the last variable at each node is the corresponding eigenvalue
85 vars_at_node[ order ] = vals[ ivec ];
86 // set the first 'order' dof to be those from the eigenvector
87 eigfn.set_nodes_vars(node, vars_at_node);
88 }
89 //// store the eigenvalue in the mesh at each node ... wasteful, but useful
90 //// a complex vector filled with the same value N times
91 //DenseVector<D_complex> c( NODES.size(), vals[ ivec ] );
92 //// add it to the mesh -- for use in nonlinear local refinement via ODE_BVP
93 //eigfn.push_var( c );
94 // add the eigenfunction to the vector of meshes
95 MESHES.push_back(eigfn);
96 }
97 }

References CppNoddy::OneD_Node_Mesh< _Type, _Xtype >::set_nodes_vars(), and CppNoddy::DenseVector< _Type >::size().

Referenced by main(), and CppNoddy::Example::Neutral_residual::Neutral_residual().

◆ eigensolve()

template<typename _Type >
void CppNoddy::ODE_EVP< _Type >::eigensolve

Formulate and solve the global eigenvalue problem for a linear system.

Definition at line 58 of file ODE_EVP.cpp.

58 {
59 // NOTE: moving construct_banded_system to the constructor is a bad idea, because
60 // sometimes we will want to construct an EVP that depends on an
61 // underlying baseflow solution that is unknown at the time of construction.
62 //
63 // construct the banded matrix eigenvalue problem that corresponds
64 // to the equation & boundary conditions specified in the constructor.
65 assemble_dense_problem();
66 // solve the system
67 p_SYSTEM -> eigensolve();
68 }
void eigensolve()
Formulate and solve the global eigenvalue problem for a linear system.
Definition: ODE_EVP.cpp:58

Referenced by main(), and CppNoddy::Example::Neutral_residual::Neutral_residual().

◆ get_mesh()

template<typename _Type >
OneD_Node_Mesh< D_complex > CppNoddy::ODE_EVP< _Type >::get_mesh ( const unsigned &  i) const
inline

Definition at line 99 of file ODE_EVP.h.

99 {
100#ifdef PARANOID
101 if(i > MESHES.size()) {
102 std::string problem;
103 problem = "You have tried to extract an eigenfunction from the ODE_EVP class\n";
104 problem += "whose index is outside the range of stored meshes.\n";
105 throw ExceptionRange(problem, MESHES.size(), i);
106 }
107#endif
108 return MESHES[ i ];
109 }

Referenced by main(), and CppNoddy::Example::Neutral_residual::Neutral_residual().

◆ p_eigensystem()

template<typename _Type >
LinearEigenSystem_base * CppNoddy::ODE_EVP< _Type >::p_eigensystem

Allow access to the underlying dense linear eigensystem through a pointer to the private member data.

Definition at line 53 of file ODE_EVP.cpp.

53 {
54 return p_SYSTEM;
55 }

Referenced by main(), and CppNoddy::Example::Neutral_residual::Neutral_residual().


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

© 2012

R.E. Hewitt