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

Solving the nonlinear problem. More...

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

Go to the source code of this file.

Classes

class  CppNoddy::Example::nonlinear
 
class  CppNoddy::Example::BC_lower
 
class  CppNoddy::Example::BC_upper
 

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 nonlinear problem.

\[ -u u_x - u_t + u_{yy} = y^2 t^2 e^{-2x} - y e^{-x}\,, \]

subject to $ u(x=0,y,t) = ty $, $ u(x,t=0,y) = 0 $ and $ f(x, y = 0,t ) = 0 $, $ f(x, y = 1,t ) = t $. The solution is computed over a range of $ t $ for $ (x,y) \in [0,10]\times [0,1] $ and the maximum deviation away from the exact solution $ u = yte^{-x} $ 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 $. This is the ‘fast’ version in which a number of default finite-differenced methods are overloaded with their analytic versions.

Definition in file IBVPNonlinearFast.cpp.

Enumeration Type Documentation

◆ anonymous enum

anonymous enum
Enumerator
Ud 

Definition at line 20 of file IBVPNonlinearFast.cpp.

20{ U, Ud };

Function Documentation

◆ main()

int main ( )

Definition at line 123 of file IBVPNonlinearFast.cpp.

124{
125 cout << "\n";
126 cout << "=== double_IBVP: non-constant mass matrix example ===\n";
127 cout << "\n";
128
129 // instantiate the problem
130 Example::nonlinear problem;
131 Example::BC_lower BC_bottom;
132 Example::BC_upper BC_top;
133
134 // domain definition
135 double top = 1.0;
136 double bottom = 0.0;
137 double left = 0.0;
138 double right = 10.0;
139 // number of (spatial) nodal points
140 unsigned ny = 11;
141 unsigned nx = 1001;
142 // time and time step
143 double t_end = 1.;
144 double dt = 0.01;
145
146 // construct our IBVP
147 PDE_double_IBVP<double> nlin( &problem,
148 Utility::uniform_node_vector( left, right, nx ),
149 Utility::uniform_node_vector( bottom, top, ny ),
150 &BC_bottom, &BC_top );
151
152 // initial conditions
153 for ( unsigned i = 0; i < nx; ++i )
154 {
155 for ( unsigned j = 0; j < ny; ++j )
156 {
157 nlin.solution()( i, j, U ) = 0.0;
158 nlin.solution()( i, j, Ud ) = 0.0;
159 }
160 }
161
162 double max_error( 0.0 );
163 int counter( 0 );
164 do
165 {
166 // inlet profile is time dependent, so we set it here for the next time level
167 nlin.update_previous_solution();
168 // now we can alter the solution stored in the double_IBVP class to define
169 // the desired x=0 'initial' conditions at the next time step
170 double t_next = nlin.t() + dt;
171 for ( unsigned j = 0; j < ny; ++j )
172 {
173 double y = nlin.solution().coord( 0, j ).second;
174 nlin.solution()( 0, j, U ) = t_next * y;
175 nlin.solution()( 0, j, Ud ) = t_next;
176 }
177 nlin.step2( dt );
178 for ( unsigned i = 0; i < nx; ++i )
179 {
180 for ( unsigned j = 0; j < ny; ++j )
181 {
182 double x = nlin.solution().coord( i, j ).first;
183 double y = nlin.solution().coord( i, j ).second;
184 double exact_u = y * nlin.t() * exp( - x );
185 max_error = max( max_error, abs( exact_u - nlin.solution()( i, j, U ) ) );
186 }
187 }
188 ++counter;
189 }
190 while ( nlin.t() < t_end );
191
192 const double tol( 1.e-4 );
193 // check the error from the exact solution
194 if ( max_error > tol )
195 {
196 cout << "\033[1;31;48m * FAILED \033[0m\n";
197 cout << " Difference = " << max_error << "\n";
198 return 1;
199 }
200 else
201 {
202 cout << "\033[1;32;48m * PASSED \033[0m\n";
203 return 0;
204 }
205
206}
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