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

A templated object for real/complex vector system of unsteady equations. More...

#include <PDE_double_IBVP.h>

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

Public Member Functions

 PDE_double_IBVP (Equation_3matrix< _Type > *equation_ptr, const DenseVector< double > &xnodes, const DenseVector< double > &ynodes, Residual_with_coords< _Type > *ptr_to_bottom_residual, Residual_with_coords< _Type > *ptr_to_top_residual)
 The class is defined by a vector function for the system. More...
 
 ~PDE_double_IBVP ()
 Destructor. More...
 
void update_previous_solution ()
 Copy the current solution to the previous solution. More...
 
void step2 (const double &dt)
 A Crank-Nicolson 'time' stepper. More...
 
double & t ()
 Return a reference to the current value of the 'timelike' coordinate. More...
 
double & tolerance ()
 Return a reference to the convergence tolerance. More...
 
TwoD_Node_Mesh< _Type > & solution ()
 

Detailed Description

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

A templated object for real/complex vector system of unsteady equations.

Definition at line 44 of file PDE_double_IBVP.h.

Constructor & Destructor Documentation

◆ PDE_double_IBVP()

template<typename _Type >
CppNoddy::PDE_double_IBVP< _Type >::PDE_double_IBVP ( Equation_3matrix< _Type > *  equation_ptr,
const DenseVector< double > &  xnodes,
const DenseVector< double > &  ynodes,
Residual_with_coords< _Type > *  ptr_to_bottom_residual,
Residual_with_coords< _Type > *  ptr_to_top_residual 
)

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

Parameters
equation_ptrA pointer to an inherited Equation object.
xnodesA vector that defines the nodal x-positions.
ynodesA vector that defines the nodal y-positions.
ptr_to_bottom_residualA pointer to a residual object that defines the y=y1 boundary conditions.
ptr_to_top_residualA pointer to a residual object that defines the y=y2 boundary conditions.

Definition at line 18 of file PDE_double_IBVP.cpp.

23 :
24 TOL(1.e-8),
25 T(0.0),
26 MAX_ITERATIONS(20),
27 p_EQUATION(ptr_to_equation),
28 p_BOTTOM_RESIDUAL(ptr_to_bottom_residual),
29 p_TOP_RESIDUAL(ptr_to_top_residual),
30 UPDATED(false) {
31 // create the 2D mesh for the current time level and previous time level
32 SOLN = TwoD_Node_Mesh<_Type>(xnodes, ynodes, p_EQUATION -> get_order());
33 // copy construct the previous solution's storage
34 PREV_SOLN = SOLN;
35 // initialise the eqn object
36 // "x"
37 p_EQUATION -> coord(2) = xnodes[0];
38#ifdef TIME
39 // timers
40 T_ASSEMBLE = Timer("Assembling of the matrix (incl. equation updates):");
41 T_SOLVE = Timer("Solving of the matrix:");
42#endif
43 }

◆ ~PDE_double_IBVP()

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

Destructor.

Definition at line 46 of file PDE_double_IBVP.cpp.

46 {
47#ifdef TIME
48 std::cout << "\n";
49 T_ASSEMBLE.stop();
50 T_ASSEMBLE.print();
51 T_SOLVE.stop();
52 T_SOLVE.print();
53#endif
54 }

Member Function Documentation

◆ solution()

template<typename _Type >
TwoD_Node_Mesh< _Type > & CppNoddy::PDE_double_IBVP< _Type >::solution
inline
Returns
A handle to the solution mesh

Definition at line 123 of file PDE_double_IBVP.h.

123 {
124 return SOLN;
125 }

Referenced by main().

◆ step2()

template<typename _Type >
void CppNoddy::PDE_double_IBVP< _Type >::step2 ( const double &  dt)

A Crank-Nicolson 'time' stepper.

Definition at line 57 of file PDE_double_IBVP.cpp.

57 {
58 DT = dt;
59 if(!UPDATED) {
60 std::string problem;
61 problem = " You have called the PDE_double_IBVP::step2 method without calling\n";
62 problem += " the PDE_double_IBVP::update_previous_solution method first. This\n";
63 problem += " method is required to copy/store the previous time level's solution\n";
64 problem += " and then allow you to specify/alter the x-initial condition by\n";
65 problem += " for the current time level.";
66 throw ExceptionRuntime(problem);
67 }
68 // the order of the problem
69 unsigned order(p_EQUATION -> get_order());
70 // get the number of nodes in the mesh
71 // -- this may have been refined by the user since the last call.
72 unsigned nx(SOLN.get_nnodes().first);
73 unsigned ny(SOLN.get_nnodes().second);
74 // measure of maximum residual
75 double max_residual(1.0);
76 // the current solution moves to the previous solution
77 // ANY LARGE STORAGE USED IN THE MAIN LOOP IS
78 // DEFINED HERE TO AVOID REPEATED CONSTRUCTION.
79 // Note we blank the A matrix after every iteration.
80 //
81 // Banded LHS matrix - max obove diagonal band wTidth is
82 // from first variable at node i to last variable at node i+1
83 BandedMatrix<_Type> a(ny * order, 2 * order - 1, 0.0);
84 // RHS
85 DenseVector<_Type> b(ny * order, 0.0);
86 // linear solver definition
87#ifdef LAPACK
88 BandedLinearSystem<_Type> system(&a, &b, "lapack");
89#else
90 BandedLinearSystem<_Type> system(&a, &b, "native");
91#endif
92 for(std::size_t j = 0; j < nx - 1; ++j) {
93 //// next x-soln guess is the previous x-soln
94 //for ( std::size_t var = 0; var < order; ++var )
95 //{
96 //for ( std::size_t i = 0; i < ny; ++i )
97 //{
98 //SOLN( j + 1, i, var ) = SOLN( j, i, var );
99 //}
100 //}
101 //
102 // iteration counter
103 int counter = 0;
104 // loop until converged or too many iterations
105 do {
106 // iteration counter
107 ++counter;
108 // time the assemble phase
109#ifdef TIME
110 T_ASSEMBLE.start();
111#endif
112 assemble_matrix_problem(a, b, j);
113 max_residual = b.inf_norm();
114#ifdef DEBUG
115 std::cout.precision(12);
116 std::cout << " PDE_double_IBVP.step2 : x_j = " << SOLN.coord(j,0).first << " Residual_max = " << max_residual << " tol = " << TOL << " counter = " << counter << "\n";
117#endif
118#ifdef TIME
119 T_ASSEMBLE.stop();
120 T_SOLVE.start();
121#endif
122 // linear solver
123 system.solve();
124 // keep the solution in a OneD_GenMesh object
125 for(std::size_t i = 0; i < ny; ++i) {
126 for(std::size_t var = 0; var < order; ++var) {
127 SOLN(j + 1, i, var) += b[ i * order + var ];
128 }
129 }
130#ifdef TIME
131 T_SOLVE.stop();
132#endif
133 } while((max_residual > TOL) && (counter < MAX_ITERATIONS));
134 if(counter >= MAX_ITERATIONS) {
135 std::string problem("\n The PDE_double_IBVP.step2 method took too many iterations. \n");
136 throw ExceptionItn(problem, counter, max_residual);
137 }
138 }
139 // set the time to the updated level
140 T += DT;
141#ifdef DEBUG
142 std::cout << "[DEBUG] solved at t = " << T << "\n";
143#endif
144 UPDATED = false;
145 }

References CppNoddy::DenseVector< _Type >::inf_norm(), and CppNoddy::BandedLinearSystem< _Type >::solve().

Referenced by main().

◆ t()

template<typename _Type >
double & CppNoddy::PDE_double_IBVP< _Type >::t
inline

Return a reference to the current value of the 'timelike' coordinate.

Returns
A handle to the current timelinke coordinate stored in the object

Definition at line 128 of file PDE_double_IBVP.h.

128 {
129 return T;
130 }

Referenced by main().

◆ tolerance()

template<typename _Type >
double & CppNoddy::PDE_double_IBVP< _Type >::tolerance ( )
inline

Return a reference to the convergence tolerance.

Definition at line 76 of file PDE_double_IBVP.h.

76 {
77 return TOL;
78 }

◆ update_previous_solution()

template<typename _Type >
void CppNoddy::PDE_double_IBVP< _Type >::update_previous_solution ( )
inline

Copy the current solution to the previous solution.

Definition at line 63 of file PDE_double_IBVP.h.

63 {
64 PREV_SOLN = SOLN;
65 UPDATED = true;
66 }

Referenced by main().


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

© 2012

R.E. Hewitt