The Stokes problem for the Taylor benchmark by the discontinuous Galerkin method – error analysis.
The Stokes problem for the Taylor benchmark by the discontinuous Galerkin method – error analysis
int main(
int argc,
char**argv) {
Float err_linf_expected = (argc > 1) ? atof(argv[1]) : 1e+38;
din >> uh;
space Xh = uh.get_space();
geo omega = Xh.get_geo();
space Qh (omega, Xh.get_approx());
din >> ph;
size_t k = Xh.degree();
size_t d = omega.dimension();
Float meas_omega = integrate (omega);
Float moy_ph = integrate (omega, ph, iopt);
ph = ph - moy_ph/meas_omega;
string high_approx = "P"+to_string(k+1)+"d";
space Xh1 (omega, high_approx,
"vector"),
Qh1 (omega, high_approx);
Float err_u_l2 = sqrt(integrate (omega, norm2(uh-
u_exact()), iopt));
Float err_u_linf = euh.max_abs();
Float err_u_h1 = sqrt(integrate (omega, norm2(grad_h(euh)), iopt)
+ integrate (omega.sides(), (1/h_local())*norm2(jump(euh)), iopt));
Float err_p_l2 = sqrt(integrate (omega, sqr(ph-
p_exact()), iopt));
Float err_p_linf = eph.max_abs();
derr << "err_u_l2 = " << err_u_l2 << endl
<< "err_u_linf = " << err_u_linf << endl
<< "err_u_h1 = " << err_u_h1 << endl
<< "err_p_l2 = " << err_p_l2 << endl
<< "err_p_linf = " << err_p_linf << endl;
}
return (max(err_u_linf,err_p_linf) <= err_linf_expected) ? 0 : 1;
}
see the Float page for the full documentation
see the field page for the full documentation
see the geo page for the full documentation
see the catchmark page for the full documentation
see the environment page for the full documentation
see the integrate_option page for the full documentation
void set_family(family_type type)
see the space page for the full documentation
verbose clean transpose logscale grid shrink ball stereo iso volume skipvtk deformation fastfieldload lattice reader_on_stdin color format format format format format format format format format format format format format format format format format format dump
This file is part of Rheolef.
rheolef - reference manual
The Taylor benchmark – the exact solution of the Stokes problem.