CppNoddy  0.92
Loading...
Searching...
No Matches
Classes | Namespaces | Enumerations | Functions | Variables
EVPOrrSommerfeldEasy_lapack.cpp File Reference

Solves the following linear eigenvalue problem for values $ c $ that satisfy : More...

#include <BVP_bundle.h>
#include <EVP_bundle.h>
#include "../Utils_Fill.h"

Go to the source code of this file.

Classes

class  CppNoddy::Example::OS_evp_equation
 Define the OS equation for the global QZ EVP. More...
 
class  CppNoddy::Example::OS_bvp_equation
 Define the OSE for the local refinement procedure. More...
 
class  CppNoddy::Example::OS_evp_both_BC
 
class  CppNoddy::Example::OS_bvp_left_BC
 
class  CppNoddy::Example::OS_bvp_right_BC
 

Namespaces

namespace  CppNoddy
 A collection of OO numerical routines aimed at simple (typical) applied problems in continuum mechanics.
 
namespace  CppNoddy::Example
 

Enumerations

enum  {
  phi , phid , psi , psid ,
  eval
}
 

Functions

double CppNoddy::Example::U (double y)
 Globally define the base flow. More...
 
double CppNoddy::Example::Udd (double y)
 Globally define the base flow curvature. More...
 
int main ()
 

Variables

double CppNoddy::Example::Re
 Globally define the Reynolds number and wavenumber. More...
 
double CppNoddy::Example::alpha
 

Detailed Description

Solves the following linear eigenvalue problem for values $ c $ that satisfy :

\[ \phi''(y) - \alpha^2 \phi(y) - \psi(y) = 0\,, \]

\[ \psi''(y) - \alpha^2 \psi(y) - i \alpha Re \left \{ ( U(y) - c ) \psi(y) - U''(y) \phi \right \} = 0\,, \]

subject to $ \phi(\pm 1) = \phi'(\pm 1) = 0 $ where $ \alpha = 1.02 $, $ Re = 5772.2 $ and $ U(y) = 1 - y^2 $. These values approximately correspond to the first neutral temporal mode in plane Poiseuille flow, therefore the test to be satisfied is that an eigenvalue exists with $ \Im ( c ) \approx 0 $. Here we use the ODE_EVP class with a very low number of mesh points to get an approximation to the temporal spectrum of modes. The 'most dangerous' mode is then extracted from this set, interpolated onto a MUCH finer mesh, and used as the initial guess in a nonlinear boundary-value problem using the ODE_BVP class. This leads to a local refinement of the selected mode that could equally be applied to other eigenvalues very simply.

Definition in file EVPOrrSommerfeldEasy_lapack.cpp.

Enumeration Type Documentation

◆ anonymous enum

anonymous enum
Enumerator
phi 
phid 
psi 
psid 
eval 

Definition at line 26 of file EVPOrrSommerfeldEasy_lapack.cpp.

Function Documentation

◆ main()

int main ( )

Definition at line 161 of file EVPOrrSommerfeldEasy_lapack.cpp.

162{
163 cout << "\n";
164 cout << "=== EVP: Temporal OSE via the equation interface ====\n";
165 cout << "=== followed by local refinement of the eigenvalue.\n";
166 cout << "\n";
167
168
169 // instantiate the eigenproblem
170 Example::OS_evp_equation evp;
171 // instantiate the boundary conditions -- same applies at both boundaries
172 Example::OS_evp_both_BC BC_both;
173 Example::Re = 5772.2;
174 Example::alpha = 1.02;
175 // domain is -1 to 1
176 double left = -1.0;
177 double right = 1.0;
178 // number of nodal points
179 int N = 32;
180 // mesh
181 DenseVector<double> nodes( Utility::uniform_node_vector( left, right, N ) );
182
183 // pass it to the ode
184 ODE_EVP<D_complex> ode_global( &evp, nodes, &BC_both, &BC_both );
185 try
186 {
187 // solve the global eigenvalue problem on ropey mesh
188 ode_global.eigensolve();
189 }
190 catch (const std::runtime_error &error )
191 {
192 cout << " \033[1;31;48m * FAILED THROUGH EXCEPTION BEING RAISED \033[0m\n";
193 return 1;
194 }
195 // get the eigenvalues near the origin
196 ode_global.p_eigensystem() -> tag_eigenvalues_disc( + 1, 0.6 );
197 // make a mesh of each tagged eigenmode in the ODE_EVP class
198 ode_global.add_tagged_to_mesh();
199
200 // instantiate the nonlinear BVP formulation
201 Example::OS_bvp_equation bvp;
202 // instantiate the BCs
203 Example::OS_bvp_left_BC BC_left;
204 Example::OS_bvp_right_BC BC_right;
205 // pass it to the ode BVP solver
206 ODE_BVP<D_complex> ode_local( &bvp, nodes, &BC_left, &BC_right );
207 // use our not-great initial guess from QZ code to start the solver off
208 // here we just use the first (0) tagged eigenvalue
209 ode_local.solution() = ode_global.get_mesh( 0 );
210 // interpolate onto a finer mesh
211 ode_local.solution().remesh1( Utility::uniform_node_vector( left, right, 1024 ) );
212 ode_local.set_monitor_det( false );
213 // solve the nonlinear BVP
214 ode_local.solve2();
215 // check the eigenvalue
216 const double tol = 1.e-4;
217 if ( abs( ode_local.solution()( 0, eval ).imag() ) < tol )
218 {
219 cout << "\033[1;32;48m * PASSED \033[0m\n";
220 return 0;
221 }
222 cout << "\033[1;31;48m * FAILED \033[0m\n";
223 return 1;
224}
An DenseVector class – a dense vector object.
Definition: DenseVector.h:34
A templated object for real/complex vector system of first-order ordinary differential equations.
Definition: ODE_BVP.h:37
A templated object for real/complex vector system of first-order ordinary differential equations.
Definition: ODE_EVP.h:35
DenseVector< double > uniform_node_vector(const double &lower, const double &upper, const std::size_t &N)
Return a DENSE vector with the nodal points of a uniform mesh distributed between the upper/lower bou...
Definition: Utility.cpp:113

References CppNoddy::ODE_EVP< _Type >::add_tagged_to_mesh(), CppNoddy::ODE_EVP< _Type >::eigensolve(), eval, CppNoddy::ODE_EVP< _Type >::get_mesh(), CppNoddy::ODE_EVP< _Type >::p_eigensystem(), CppNoddy::ODE_BVP< _Type, _Xtype >::set_monitor_det(), CppNoddy::ODE_BVP< _Type, _Xtype >::solution(), CppNoddy::ODE_BVP< _Type, _Xtype >::solve2(), and CppNoddy::Utility::uniform_node_vector().

© 2012

R.E. Hewitt