MFEM  v4.2.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
convergence.hpp
Go to the documentation of this file.
1 // Copyright (c) 2010-2020, Lawrence Livermore National Security, LLC. Produced
2 // at the Lawrence Livermore National Laboratory. All Rights reserved. See files
3 // LICENSE and NOTICE for details. LLNL-CODE-806117.
4 //
5 // This file is part of the MFEM library. For more information and source code
6 // availability visit https://mfem.org.
7 //
8 // MFEM is free software; you can redistribute it and/or modify it under the
9 // terms of the BSD-3 license. We welcome feedback and contributions, see file
10 // CONTRIBUTING.md for details.
11 
12 #ifndef MFEM_CONVERGENCE
13 #define MFEM_CONVERGENCE
14 
15 #include "../linalg/linalg.hpp"
16 #include "gridfunc.hpp"
17 #ifdef MFEM_USE_MPI
18 #include "pgridfunc.hpp"
19 #endif
20 
21 namespace mfem
22 {
23 
24 /** @brief Class to compute error and convergence rates.
25  It supports H1, H(curl) (ND elements), H(div) (RT elements) and L2 (DG).
26 
27  For "smooth enough" solutions the Galerkin error measured in the appropriate
28  norm satisfies || u - u_h || ~ h^k
29 
30  Here, k is called the asymptotic rate of convergence
31 
32  For successive uniform h-refinements the rate can be estimated by
33  k = log(||u - u_h|| / ||u - u_{h/2}||)/log(2)
34 */
36 {
37 private:
38  // counters for solutions/derivatives
39  int counter=0;
40  int dcounter=0;
41  int fcounter=0;
42 
43  // space continuity type
44  int cont_type=-1;
45 
46  // printing flag for helpful for MPI calls
47  int print_flag=1;
48 
49  // exact solution and derivatives
50  double CoeffNorm;
51  double CoeffDNorm;
52 
53  // Arrays to store error/rates
54  Array<double> L2Errors, DGFaceErrors, DErrors, EnErrors;
55  Array<double> L2Rates, DGFaceRates, DRates, EnRates;
56  Array<int> ndofs;
57 
58  void AddL2Error(GridFunction *gf, Coefficient *scalar_u,
59  VectorCoefficient *vector_u);
60  void AddGf(GridFunction *gf, Coefficient *scalar_u,
61  VectorCoefficient *grad=nullptr,
62  Coefficient *ell_coeff=nullptr, double Nu=1.0);
63  void AddGf(GridFunction *gf, VectorCoefficient *vector_u,
64  VectorCoefficient *curl, Coefficient *div);
65  // returns the L2-norm of scalar_u or vector_u
66  double GetNorm(GridFunction *gf, Coefficient *scalar_u,
67  VectorCoefficient *vector_u);
68 
69 public:
70 
71  /// Clear any internal data
72  void Reset();
73 
74  /// Add L2 GridFunction, the exact solution and possibly its gradient and/or
75  /// DG face jumps parameters
77  VectorCoefficient *grad=nullptr,
78  Coefficient *ell_coeff=nullptr, double Nu=1.0)
79  {
80  AddGf(gf, scalar_u, grad, ell_coeff, Nu);
81  }
82 
83  /// Add H1 GridFunction, the exact solution and possibly its gradient
85  VectorCoefficient *grad=nullptr)
86  {
87  AddGf(gf, scalar_u, grad);
88  }
89 
90  /// Add H(curl) GridFunction, the exact solution and possibly its curl
92  VectorCoefficient *curl=nullptr)
93  {
94  AddGf(gf, vector_u, curl, nullptr);
95  }
96 
97  /// Add H(div) GridFunction, the exact solution and possibly its div
99  Coefficient *div=nullptr)
100  {
101  AddGf(gf,vector_u, nullptr, div);
102  }
103 
104  /// Get the L2 error at step n
105  double GetL2Error(int n)
106  {
107  MFEM_VERIFY( n <= counter,"Step out of bounds")
108  return L2Errors[n];
109  }
110 
111  /// Get all L2 errors
112  void GetL2Errors(Array<double> & L2Errors_)
113  {
114  L2Errors_ = L2Errors;
115  }
116 
117  /// Get the Grad/Curl/Div error at step n
118  double GetDError(int n)
119  {
120  MFEM_VERIFY(n <= dcounter,"Step out of bounds")
121  return DErrors[n];
122  }
123 
124  /// Get all Grad/Curl/Div errors
125  void GetDErrors(Array<double> & DErrors_)
126  {
127  DErrors_ = DErrors;
128  }
129 
130  /// Get the DGFaceJumps error at step n
131  double GetDGFaceJumpsError(int n)
132  {
133  MFEM_VERIFY(n<= fcounter,"Step out of bounds")
134  return DGFaceErrors[n];
135  }
136 
137  /// Get all DGFaceJumps errors
138  void GetDGFaceJumpsErrors(Array<double> & DGFaceErrors_)
139  {
140  DGFaceErrors_ = DGFaceErrors;
141  }
142 
143  /// Print rates and errors
144  void Print(bool relative = false, std::ostream &out = mfem::out);
145 };
146 
147 } // namespace mfem
148 
149 #endif // MFEM_CONVERGENCE
Class for grid function - Vector with associated FE space.
Definition: gridfunc.hpp:30
void Reset()
Clear any internal data.
Definition: convergence.cpp:9
Base class for vector Coefficients that optionally depend on time and space.
void AddHcurlGridFunction(GridFunction *gf, VectorCoefficient *vector_u, VectorCoefficient *curl=nullptr)
Add H(curl) GridFunction, the exact solution and possibly its curl.
Definition: convergence.hpp:91
void AddHdivGridFunction(GridFunction *gf, VectorCoefficient *vector_u, Coefficient *div=nullptr)
Add H(div) GridFunction, the exact solution and possibly its div.
Definition: convergence.hpp:98
void GetDGFaceJumpsErrors(Array< double > &DGFaceErrors_)
Get all DGFaceJumps errors.
double GetL2Error(int n)
Get the L2 error at step n.
void GetDErrors(Array< double > &DErrors_)
Get all Grad/Curl/Div errors.
void AddL2GridFunction(GridFunction *gf, Coefficient *scalar_u, VectorCoefficient *grad=nullptr, Coefficient *ell_coeff=nullptr, double Nu=1.0)
Definition: convergence.hpp:76
double GetDGFaceJumpsError(int n)
Get the DGFaceJumps error at step n.
void GetL2Errors(Array< double > &L2Errors_)
Get all L2 errors.
Base class Coefficients that optionally depend on space and time. These are used by the BilinearFormI...
Definition: coefficient.hpp:39
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).
Definition: convergence.hpp:35
OutStream out(std::cout)
Global stream used by the library for standard output. Initially it uses the same std::streambuf as s...
Definition: globals.hpp:66
void AddH1GridFunction(GridFunction *gf, Coefficient *scalar_u, VectorCoefficient *grad=nullptr)
Add H1 GridFunction, the exact solution and possibly its gradient.
Definition: convergence.hpp:84
double GetDError(int n)
Get the Grad/Curl/Div error at step n.