Rheolef  7.2
an efficient C++ finite element environment
 
Loading...
Searching...
No Matches
navier_stokes_taylor_newton_dg.cc
Go to the documentation of this file.
1
25#include "rheolef.h"
26using namespace rheolef;
27using namespace std;
28#include "taylor.h"
30#include "inertia.h"
31#include "navier_stokes_dg.h"
32int main(int argc, char**argv) {
33 environment rheolef (argc, argv);
34 Float eps = numeric_limits<Float>::epsilon();
35 geo omega (argv[1]);
36 string approx = (argc > 2) ? argv[2] : "P1d";
37 Float Re = (argc > 3) ? atof(argv[3]) : 100;
38 Float tol = (argc > 4) ? atof(argv[4]) : eps;
39 size_t max_iter = (argc > 5) ? atoi(argv[5]) : 100;
40 string restart = (argc > 6) ? argv[6] : "";
41 navier_stokes_dg F (Re, omega, approx);
43 int status = damped_newton (F, xh, tol, max_iter, &derr);
44 dout << catchmark("Re") << Re << endl
45 << catchmark("u") << xh[0]
46 << catchmark("p") << xh[1];
47 return status;
48}
see the Float page for the full documentation
see the geo page for the full documentation
see the catchmark page for the full documentation
Definition catchmark.h:67
see the environment page for the full documentation
int main()
Definition field2bb.cc:58
The inertia term of the Navier-Stokes equation with the discontinuous Galerkin method – di Pietro & E...
This file is part of Rheolef.
int damped_newton(const Problem &P, const Preconditioner &T, Field &u, Real &tol, Size &max_iter, odiststream *p_derr=0)
see the damped_newton page for the full documentation
STL namespace.
The Navier-Stokes equations with the discontinuous Galerkin method – class header.
rheolef - reference manual
The Stokes problem with Dirichlet boundary condition by the discontinuous Galerkin method – solver fu...
value_type initial(string restart) const
Eigen::Matrix< field, 2, 1 > value_type
The Taylor benchmark – right-hand-side and boundary condition.