CppNoddy  0.92
Loading...
Searching...
No Matches
IBVPNonlinearSlower.cpp
Go to the documentation of this file.
1/// \file IBVPNonlinearSlower.cpp
2/// \ingroup Tests
3/// \ingroup IBVP_double
4/// Solving the nonlinear problem
5/// \f[ -u u_x - u_t + u_{yy} = y^2 t^2 e^{-2x} - y e^{-x}\,, \f]
6/// subject to \f$ u(x=0,y,t) = ty \f$, \f$ u(x,t=0,y) = 0 \f$
7/// and \f$ f(x, y = 0,t ) = 0 \f$, \f$ f(x, y = 1,t ) = t \f$.
8/// The solution is computed over a range of \f$ t \f$ for \f$ (x,y) \in [0,10]\times [0,1] \f$
9/// and the maximum deviation away from the exact solution \f$ u = yte^{-x} \f$ is found.
10/// The test fails if this deviation is larger than a set tolerance \f$ \tau \f$.
11/// The example is solved using the PDE_double_IBVP for problems that are
12/// parabolic in 2 coordinates \f$ (x,t) \f$ with a BVP in \f$ y \f$.
13
14#include <IBVP_double_bundle.h>
15
16#include "../Utils_Fill.h"
17
18enum { U, Ud };
19
20namespace CppNoddy
21{
22 namespace Example
23 {
24
25 // A source term rigged up to give a nice exact solution
26 double source( const double& x, const double& y, const double& t )
27 {
28 return -y * std::exp( -x ) + std::pow( y * t * std::exp( -x ), 2 );
29 }
30
31 class nonlinear : public Equation_3matrix<double>
32 {
33 public:
34 /// The problem is 2nd order and real
35 nonlinear() : Equation_3matrix<double> ( 2 ) {}
36
37 /// Define a nonlinear advection diffusion problem
39 {
40 // The system
41 f[ U ] = z[ Ud ];
42 f[ Ud ] = source( coord(2), coord(0), coord(1) );
43 }
44
45 /// Define the unsteady terms by providing the mass matrix for t evolution
47 {
49 }
50
51 /// Providing the matrix for t evolution
53 {
54 // eqn 1 variable 0
55 m( 1, 0 ) = -1.0;
56 }
57
58 /// Providing the matrix for x evolution
60 {
61 // eqn 1 variable 0
62 m( 1, 0 ) = -z[ U ];
63 }
64
65 };
66
67 // BOUNDARY CONDITIONS
68 class BC_lower : public Residual_with_coords<double>
69 {
70 public:
71 // 1 constraint, 2nd order system, 2 coordinates (x & t)
72 BC_lower() : Residual_with_coords<double> ( 1, 2, 2 ) {}
73
75 {
76 b[ 0 ] = z[ U ];
77 }
78 };
79
80 class BC_upper : public Residual_with_coords<double>
81 {
82 public:
83 // 1 constraint, 2nd order system, 2 coordinates (x & t)
84 BC_upper() : Residual_with_coords<double> ( 1, 2, 2 ) {}
85
87 {
88 const double x = coord( 0 );
89 const double t = coord( 1 );
90 b[ 0 ] = z[ U ] - t * std::exp( - x );
91 }
92 };
93
94 } // end Test namespace
95} // end CppNoddy namespace
96
97using namespace CppNoddy;
98using namespace std;
99
100int main()
101{
102 cout << "\n";
103 cout << "=== double_IBVP: non-constant mass matrix example ===\n";
104 cout << "\n";
105
106 // instantiate the problem
107 Example::nonlinear problem;
108 Example::BC_lower BC_bottom;
109 Example::BC_upper BC_top;
110
111 // domain definition
112 double top = 1.0;
113 double bottom = 0.0;
114 double left = 0.0;
115 double right = 10.0;
116 // number of (spatial) nodal points
117 unsigned ny = 11;
118 unsigned nx = 1001;
119 // time and time step
120 double t_end = 1.;
121 double dt = 0.01;
122
123 // construct our IBVP
124 PDE_double_IBVP<double> nlin( &problem,
125 Utility::uniform_node_vector( left, right, nx ),
126 Utility::uniform_node_vector( bottom, top, ny ),
127 &BC_bottom, &BC_top );
128
129 // initial conditions
130 for ( unsigned i = 0; i < nx; ++i )
131 {
132 for ( unsigned j = 0; j < ny; ++j )
133 {
134 nlin.solution()( i, j, U ) = 0.0;
135 nlin.solution()( i, j, Ud ) = 0.0;
136 }
137 }
138
139 double max_error( 0.0 );
140 int counter( 0 );
141 do
142 {
143 // inlet profile is time dependent, so we set it here for the next time level
145 // now we can alter the solution stored in the double_IBVP class to define
146 // the desired x=0 'initial' conditions at the next time step
147 double t_next = nlin.t() + dt;
148 for ( unsigned j = 0; j < ny; ++j )
149 {
150 double y = nlin.solution().coord( 0, j ).second;
151 nlin.solution()( 0, j, U ) = t_next * y;
152 nlin.solution()( 0, j, Ud ) = t_next;
153 }
154 nlin.step2( dt );
155 for ( unsigned i = 0; i < nx; ++i )
156 {
157 for ( unsigned j = 0; j < ny; ++j )
158 {
159 double x = nlin.solution().coord( i, j ).first;
160 double y = nlin.solution().coord( i, j ).second;
161 double exact_u = y * nlin.t() * exp( - x );
162 max_error = max( max_error, abs( exact_u - nlin.solution()( i, j, U ) ) );
163 }
164 }
165 ++counter;
166 }
167 while ( nlin.t() < t_end );
168
169 const double tol( 1.e-4 );
170 // check the error from the exact solution
171 if ( max_error > tol )
172 {
173 cout << "\033[1;31;48m * FAILED \033[0m\n";
174 cout << " Difference = " << max_error << "\n";
175 return 1;
176 }
177 else
178 {
179 cout << "\033[1;32;48m * PASSED \033[0m\n";
180 return 0;
181 }
182
183}
@ f
Definition: BVPBerman.cpp:15
@ U
Definition: BVPKarman.cpp:20
@ Ud
Definition: BVPKarman.cpp:20
int main()
A shorter bundled include file for double initial boundary value problems.
A matrix class that constructs a DENSE matrix as a row major std::vector of DenseVectors.
Definition: DenseMatrix.h:25
An DenseVector class – a dense vector object.
Definition: DenseVector.h:34
An equation object base class used in the PDE_double_IBVP class.
void residual_fn(const DenseVector< double > &z, DenseVector< double > &b) const
void residual_fn(const DenseVector< double > &z, DenseVector< double > &b) const
void residual_fn(const DenseVector< double > &z, DenseVector< double > &f) const
Define a nonlinear advection diffusion problem.
void matrix2(const DenseVector< double > &z, DenseMatrix< double > &m) const
Providing the matrix for x evolution.
void matrix1(const DenseVector< double > &z, DenseMatrix< double > &m) const
Providing the matrix for t evolution.
nonlinear()
The problem is 2nd order and real.
void matrix0(const DenseVector< double > &z, DenseMatrix< double > &m) const
Define the unsteady terms by providing the mass matrix for t evolution.
A templated object for real/complex vector system of unsteady equations.
void update_previous_solution()
Copy the current solution to the previous solution.
double & t()
Return a reference to the current value of the 'timelike' coordinate.
void step2(const double &dt)
A Crank-Nicolson 'time' stepper.
TwoD_Node_Mesh< _Type > & solution()
A base class to be inherited by objects that define residuals.
_Xtype & coord(const unsigned &i)
General handle access to the coordinates.
double source(const double &x, const double &y, const double &t)
Definition: IBVPLinear.cpp:26
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
A collection of OO numerical routines aimed at simple (typical) applied problems in continuum mechani...
void fill_identity(CppNoddy::Sequential_Matrix_base< _Type > &A)
Fill diagonal with unit values.
Definition: Utils_Fill.h:22

© 2012

R.E. Hewitt