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

Solving the nonlinear advection diffusion equation. More...

#include <IBVP_bundle.h>

Go to the source code of this file.

Classes

class  CppNoddy::Example::Nlin_adv_equation
 
class  CppNoddy::Example::Nlin_adv_left_BC
 
class  CppNoddy::Example::Nlin_adv_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  { U , Ud }
 

Functions

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

Detailed Description

Solving the nonlinear advection diffusion equation.

\[ U_{yy} - Re U U_t = S( y, t ) \]

subject to $ U(y=0,t) = 1 $ and $ U(y=1,t) = 2 $ with initial condition $ U(y,t=0) = 1 + y + \epsilon y(1-y) $ and a source term

\[ S(y,t) = ( 1 + y + \epsilon y(1-y)e^{-t} ) \epsilon y(1-y ) e^{-t} - 2\epsilon e^{-t} / Re \]

for some constant parameters $\epsilon $ and $ Re $. This source term is chosen to allow an exact solution

\[ U(y,t) = 1 + y + \epsilon y(1-y) e^{-t} \]

against which the numerical computation is compared, by computing the maximum absolute deviation over all spatial and temporal nodes.

Definition in file IBVPNonlinearAdvDiffusion.cpp.

Enumeration Type Documentation

◆ anonymous enum

anonymous enum
Enumerator
Ud 

Definition at line 19 of file IBVPNonlinearAdvDiffusion.cpp.

Function Documentation

◆ main()

int main ( )

Definition at line 107 of file IBVPNonlinearAdvDiffusion.cpp.

108{
109 cout << "\n";
110 cout << "=== IBVP: Nonlinear advection diffusion =============\n";
111 cout << "\n";
112
113 // instantiate the problem
114 Example::Nlin_adv_equation problem;
115 Example::Nlin_adv_left_BC BC_left;
116 Example::Nlin_adv_right_BC BC_right;
117
118 // domain definition
119 double left = 0.0;
120 double right = 1.0;
121 // number of (spatial) nodal points
122 unsigned ny = 400;
123 // number of time steps
124 unsigned max_steps = 4000;
125 // time step
126 double dt = 1.0 / ( max_steps - 1 );
127
128 // initial rotation frequency
129 Example::Re = 1.0;
130
131 // construct our IBVP
132 PDE_IBVP<double> advection( &problem, Utility::uniform_node_vector( left, right, ny ), &BC_left, &BC_right );
133
134 // initial conditions
135 for ( unsigned i = 0; i < ny; ++i )
136 {
137 double y = advection.solution().coord( i );
138 advection.solution()( i, U ) = 1.0 + y + Example::eps * y * ( 1 - y );
139 advection.solution()( i, Ud ) = 1.0;
140 }
141
142 // set up output details
143 std::string dirname("./DATA");
144 mkdir( dirname.c_str(), S_IRWXU );
145 TrackerFile my_file( "./DATA/IBVP_NlinAdvection.dat" );
146 my_file.push_ptr( &advection.coord(), "time" );
147 my_file.push_ptr( &advection.solution(), "Solution profile" );
148 my_file.header();
149
150 double max_err( 0.0 );
151 // time step
152 for ( unsigned i = 1; i < max_steps; ++i )
153 {
154 // take a time step
155 try
156 {
157 advection.step2( dt );
158 }
159 catch ( const std::runtime_error &error )
160 {
161 cout << " \033[1;31;48m * FAILED THROUGH EXCEPTION BEING RAISED \033[0m\n";
162 return 1;
163 }
164 if ( i % 100 == 0 )
165 {
166 my_file.update();
167 my_file.newline();
168 }
169 // generate the exact solution to compare to
170 for ( unsigned i = 0; i < ny; ++i )
171 {
172 double y = advection.solution().coord( i );
173 const double exact_solution = 1.0 + y + Example::eps * y * ( 1 - y ) * std::exp( -advection.coord() );
174 max_err = std::max( max_err, std::abs( advection.solution()( i, 0 ) - exact_solution ) );
175 }
176 }
177
178 if ( max_err > 1.e-6 )
179 {
180 cout << "\033[1;31;48m * FAILED \033[0m\n";
181 cout.precision( 14 );
182 cout << 1. / ( ny - 1 ) << " " << dt << " " << max_err << "\n";
183 return 1;
184 }
185 else
186 {
187 cout << "\033[1;32;48m * PASSED \033[0m\n";
188 return 0;
189 }
190}
A templated object for real/complex vector system of unsteady equations.
Definition: PDE_IBVP.h:37
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_IBVP< _Type >::coord(), CppNoddy::TrackerFile::header(), CppNoddy::TrackerFile::newline(), CppNoddy::TrackerFile::push_ptr(), CppNoddy::PDE_IBVP< _Type >::solution(), CppNoddy::PDE_IBVP< _Type >::step2(), U, Ud, CppNoddy::Utility::uniform_node_vector(), and CppNoddy::TrackerFile::update().

© 2012

R.E. Hewitt