CppNoddy  0.92
All Classes Namespaces Files Functions Variables Typedefs Enumerator Friends Macros Modules Pages
Classes | Namespaces | Enumerations | Functions
IBVPLinear.cpp File Reference

Solving the linear equation. More...

#include <IBVP_double_bundle.h>
#include "../Utils_Fill.h"

Go to the source code of this file.

Classes

class  CppNoddy::Example::diffusion_double
 
class  CppNoddy::Example::dd_BC
 

Namespaces

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

Enumerations

enum  { U , Ud }
 

Functions

double CppNoddy::Example::source (const double &x, const double &y, const double &t)
 
int main ()
 

Detailed Description

Solving the linear equation.

\[ - u_x - u_t + u_{yy} = ( ( 1 - y^2 )( x + t ) - 2 ) e^{-xt}\,, \]

subject to $ u(x=0,y,t) = 1-y^2 $ $ u(x,t=0,y) = 1 - y^2 $ and $ f(x, y = \pm 1,t ) = 0 $. The solution is computed over a range of $ t $ for $ (x,y) \in [0,10]\times [-1,1] $ and the maximum deviation away from the exact solution $ u = (1-y^2) e^{-xt} $ is found. The test fails if this deviation is larger than a set tolerance $ \tau $. The example is solved using the PDE_double_IBVP for problems that are parabolic in 2 coordinates $ (x,t) $ with a BVP in $ y $.

Definition in file IBVPLinear.cpp.

Enumeration Type Documentation

◆ anonymous enum

anonymous enum
Enumerator
Ud 

Definition at line 18 of file IBVPLinear.cpp.

18{ U, Ud };
@ U
Definition: IBVPLinear.cpp:18
@ Ud
Definition: IBVPLinear.cpp:18

Function Documentation

◆ main()

int main ( )

Definition at line 104 of file IBVPLinear.cpp.

105{
106 cout << "\n";
107 cout << "=== double_IBVP: linear diffusive problem ===========\n";
108 cout << "\n";
109
110 // instantiate the problem
111 Example::diffusion_double problem;
112 Example::dd_BC BC_bottom;
113 Example::dd_BC BC_top;
114
115 // domain definition
116 double top = 1.0;
117 double bottom = -1.0;
118 double left = 0.0;
119 double right = 10.0;
120 // number of (spatial) nodal points
121 unsigned ny = 301;
122 unsigned nx = 301;
123 // time and time step
124 double t_end = 1.0;
125 double dt = 0.01;
126
127 // construct our IBVP
128 PDE_double_IBVP<double> diffusion( &problem,
129 Utility::uniform_node_vector( left, right, nx ),
130 Utility::uniform_node_vector( bottom, top, ny ),
131 &BC_bottom, &BC_top );
132
133 // initial conditions
134 for ( unsigned i = 0; i < nx; ++i )
135 {
136 for ( unsigned j = 0; j < ny; ++j )
137 {
138 double y = diffusion.solution().coord( i, j ).second;
139 diffusion.solution()( i, j, U ) = ( 1.0 - y * y );
140 diffusion.solution()( i, j, Ud ) = - 2 * y ;
141 }
142 }
143
144 double max_error( 0.0 );
145 do
146 {
147 diffusion.update_previous_solution();
148 diffusion.step2( dt );
149 for ( unsigned i = 0; i < nx; ++i )
150 {
151 for ( unsigned j = 0; j < ny; ++j )
152 {
153 double x = diffusion.solution().coord( i, j ).first;
154 double y = diffusion.solution().coord( i, j ).second;
155 double exact_u = ( 1 - y * y ) * exp( - x * diffusion.t() );
156 max_error = max( max_error, abs( exact_u - diffusion.solution()( i, j, U ) ) );
157 }
158 }
159 }
160 while ( diffusion.t() < t_end );
161
162 const double tol( 1.e-4 );
163 // check the BL transpiration vs the known solution
164 if ( max_error > tol )
165 {
166 cout << "\033[1;31;48m * FAILED \033[0m\n";
167 cout << " Difference = " << max_error << "\n";
168 return 1;
169 }
170 else
171 {
172 cout << "\033[1;32;48m * PASSED \033[0m\n";
173 return 0;
174 }
175
176}
A templated object for real/complex vector system of unsteady equations.
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::PDE_double_IBVP< _Type >::solution(), CppNoddy::PDE_double_IBVP< _Type >::step2(), CppNoddy::PDE_double_IBVP< _Type >::t(), U, Ud, CppNoddy::Utility::uniform_node_vector(), and CppNoddy::PDE_double_IBVP< _Type >::update_previous_solution().

© 2012

R.E. Hewitt