MFEM v4.7.0
Finite element discretization library
|
Class to compute error and convergence rates. It supports H1, H(curl) (ND elements), H(div) (RT elements) and L2 (DG). More...
#include <convergence.hpp>
Public Member Functions | |
void | Reset () |
Clear any internal data. | |
void | AddL2GridFunction (GridFunction *gf, Coefficient *scalar_u, VectorCoefficient *grad=nullptr, Coefficient *ell_coeff=nullptr, JumpScaling jump_scaling={1.0, JumpScaling::ONE_OVER_H}) |
void | AddH1GridFunction (GridFunction *gf, Coefficient *scalar_u, VectorCoefficient *grad=nullptr) |
Add H1 GridFunction, the exact solution and possibly its gradient. | |
void | AddHcurlGridFunction (GridFunction *gf, VectorCoefficient *vector_u, VectorCoefficient *curl=nullptr) |
Add H(curl) GridFunction, the exact solution and possibly its curl. | |
void | AddHdivGridFunction (GridFunction *gf, VectorCoefficient *vector_u, Coefficient *div=nullptr) |
Add H(div) GridFunction, the exact solution and possibly its div. | |
real_t | GetL2Error (int n) |
Get the L2 error at step n. | |
void | GetL2Errors (Array< real_t > &L2Errors_) |
Get all L2 errors. | |
real_t | GetDError (int n) |
Get the Grad/Curl/Div error at step n. | |
void | GetDErrors (Array< real_t > &DErrors_) |
Get all Grad/Curl/Div errors. | |
real_t | GetDGFaceJumpsError (int n) |
Get the DGFaceJumps error at step n. | |
void | GetDGFaceJumpsErrors (Array< real_t > &DGFaceErrors_) |
Get all DGFaceJumps errors. | |
void | Print (bool relative=false, std::ostream &out=mfem::out) |
Print rates and errors. | |
Class to compute error and convergence rates. It supports H1, H(curl) (ND elements), H(div) (RT elements) and L2 (DG).
For "smooth enough" solutions the Galerkin error measured in the appropriate norm satisfies || u - u_h || ~ h^k
Here, k is called the asymptotic rate of convergence
For successive h-refinements the rate can be estimated by k = log(||u - u_h|| / ||u - u_{h/2}||)/(1/dim * log(N_{h/2}/N_{h})
Definition at line 35 of file convergence.hpp.
|
inline |
Add H1 GridFunction, the exact solution and possibly its gradient.
Definition at line 86 of file convergence.hpp.
|
inline |
Add H(curl) GridFunction, the exact solution and possibly its curl.
Definition at line 93 of file convergence.hpp.
|
inline |
Add H(div) GridFunction, the exact solution and possibly its div.
Definition at line 100 of file convergence.hpp.
|
inline |
Add L2 GridFunction, the exact solution and possibly its gradient and/or DG face jumps parameters
Definition at line 77 of file convergence.hpp.
|
inline |
Get the Grad/Curl/Div error at step n.
Definition at line 120 of file convergence.hpp.
Get all Grad/Curl/Div errors.
Definition at line 127 of file convergence.hpp.
|
inline |
Get the DGFaceJumps error at step n.
Definition at line 133 of file convergence.hpp.
Get all DGFaceJumps errors.
Definition at line 140 of file convergence.hpp.
|
inline |
Get the L2 error at step n.
Definition at line 107 of file convergence.hpp.
Get all L2 errors.
Definition at line 114 of file convergence.hpp.
void mfem::ConvergenceStudy::Print | ( | bool | relative = false, |
std::ostream & | out = mfem::out ) |
Print rates and errors.
Definition at line 236 of file convergence.cpp.
void mfem::ConvergenceStudy::Reset | ( | ) |
Clear any internal data.
Definition at line 19 of file convergence.cpp.