MFEM  v4.2.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
convergence.cpp
Go to the documentation of this file.
1 #include "convergence.hpp"
2 
3 using namespace std;
4 
5 
6 namespace mfem
7 {
8 
9 void ConvergenceStudy::Reset()
10 {
11  counter=0;
12  dcounter=0;
13  fcounter=0;
14  cont_type=-1;
15  print_flag=1;
16  L2Errors.SetSize(0);
17  L2Rates.SetSize(0);
18  DErrors.SetSize(0);
19  DRates.SetSize(0);
20  EnErrors.SetSize(0);
21  EnRates.SetSize(0);
22  DGFaceErrors.SetSize(0);
23  DGFaceRates.SetSize(0);
24  ndofs.SetSize(0);
25 }
26 
27 double ConvergenceStudy::GetNorm(GridFunction *gf, Coefficient *scalar_u,
28  VectorCoefficient *vector_u)
29 {
30  bool norm_set = false;
31  double norm=0.0;
32  int order = gf->FESpace()->GetOrder(0);
33  int order_quad = std::max(2, 2*order+1);
34  const IntegrationRule *irs[Geometry::NumGeom];
35  for (int i=0; i < Geometry::NumGeom; ++i)
36  {
37  irs[i] = &(IntRules.Get(i, order_quad));
38  }
39 #ifdef MFEM_USE_MPI
40  ParGridFunction *pgf = dynamic_cast<ParGridFunction *>(gf);
41  if (pgf)
42  {
43  ParMesh *pmesh = pgf->ParFESpace()->GetParMesh();
44  if (scalar_u)
45  {
46  norm = ComputeGlobalLpNorm(2.0,*scalar_u,*pmesh,irs);
47  }
48  else if (vector_u)
49  {
50  norm = ComputeGlobalLpNorm(2.0,*vector_u,*pmesh,irs);
51  }
52  norm_set = true;
53  }
54 #endif
55  if (!norm_set)
56  {
57  Mesh *mesh = gf->FESpace()->GetMesh();
58  if (scalar_u)
59  {
60  norm = ComputeLpNorm(2.0,*scalar_u,*mesh,irs);
61  }
62  else if (vector_u)
63  {
64  norm = ComputeLpNorm(2.0,*vector_u,*mesh,irs);
65  }
66  }
67  return norm;
68 }
69 
70 void ConvergenceStudy::AddL2Error(GridFunction *gf,
71  Coefficient *scalar_u, VectorCoefficient *vector_u)
72 {
73  int tdofs=0;
74 #ifdef MFEM_USE_MPI
75  ParGridFunction *pgf = dynamic_cast<ParGridFunction *>(gf);
76  if (pgf)
77  {
78  MPI_Comm comm = pgf->ParFESpace()->GetComm();
79  int rank;
80  MPI_Comm_rank(comm, &rank);
81  print_flag = 0;
82  if (rank==0) { print_flag = 1; }
83  tdofs = pgf->ParFESpace()->GlobalTrueVSize();
84  }
85 #endif
86  if (!tdofs) { tdofs = gf->FESpace()->GetTrueVSize(); }
87  ndofs.Append(tdofs);
88  double L2Err;
89  if (scalar_u)
90  {
91  L2Err = gf->ComputeL2Error(*scalar_u);
92  CoeffNorm = GetNorm(gf,scalar_u,nullptr);
93  }
94  else if (vector_u)
95  {
96  L2Err = gf->ComputeL2Error(*vector_u);
97  CoeffNorm = GetNorm(gf,nullptr,vector_u);
98  }
99  else
100  {
101  MFEM_ABORT("Exact Solution Coefficient pointer is NULL");
102  }
103  L2Errors.Append(L2Err);
104  // Compute the rate of convergence by:
105  // rate = log (||u - u_h|| / ||u - u_{h/2}||)/log(2)
106  double val = (counter) ? log(L2Errors[counter-1]/L2Err)/log(2.0) : 0.0;
107  L2Rates.Append(val);
108  counter++;
109 }
110 
111 void ConvergenceStudy::AddGf(GridFunction *gf, Coefficient *scalar_u,
112  VectorCoefficient *grad,
113  Coefficient *ell_coeff, double Nu)
114 {
115  cont_type = gf->FESpace()->FEColl()->GetContType();
116 
117  MFEM_VERIFY((cont_type == mfem::FiniteElementCollection::CONTINUOUS) ||
119  "This constructor is intended for H1 or L2 Elements")
120 
121  AddL2Error(gf,scalar_u, nullptr);
122 
123  if (grad)
124  {
125  double GradErr = gf->ComputeGradError(grad);
126  DErrors.Append(GradErr);
127  double err = sqrt(L2Errors[counter-1]*L2Errors[counter-1]+GradErr*GradErr);
128  EnErrors.Append(err);
129  // Compute the rate of convergence by:
130  // rate = log (||u - u_h|| / ||u - u_{h/2}||)/log(2)
131  double val = (dcounter) ? log(DErrors[dcounter-1]/GradErr)/log(2.0) : 0.0;
132  double eval = (dcounter) ? log(EnErrors[dcounter-1]/err)/log(2.0) : 0.0;
133  DRates.Append(val);
134  EnRates.Append(eval);
135  CoeffDNorm = GetNorm(gf,nullptr,grad);
136  dcounter++;
137  MFEM_VERIFY(counter == dcounter,
138  "Number of added solutions and derivatives do not match")
139  }
140 
141  if (cont_type == mfem::FiniteElementCollection::DISCONTINUOUS && ell_coeff)
142  {
143  double DGErr = gf->ComputeDGFaceJumpError(scalar_u,ell_coeff,Nu);
144  DGFaceErrors.Append(DGErr);
145  // Compute the rate of convergence by:
146  // rate = log (||u - u_h|| / ||u - u_{h/2}||)/log(2)
147  double val=(fcounter) ? log(DGFaceErrors[fcounter-1]/DGErr)/log(2.0):0.;
148  DGFaceRates.Append(val);
149  fcounter++;
150  MFEM_VERIFY(fcounter == counter, "Number of added solutions mismatch");
151  }
152 }
153 
154 void ConvergenceStudy::AddGf(GridFunction *gf, VectorCoefficient *vector_u,
155  VectorCoefficient *curl, Coefficient *div)
156 {
157  cont_type = gf->FESpace()->FEColl()->GetContType();
158 
159  AddL2Error(gf,nullptr,vector_u);
160  double DErr = 0.0;
161  bool derivative = false;
162  if (curl)
163  {
164  DErr = gf->ComputeCurlError(curl);
165  CoeffDNorm = GetNorm(gf,nullptr,curl);
166  derivative = true;
167  }
168  else if (div)
169  {
170  DErr = gf->ComputeDivError(div);
171  // update coefficient norm
172  CoeffDNorm = GetNorm(gf,div,nullptr);
173  derivative = true;
174  }
175  if (derivative)
176  {
177  double err = sqrt(L2Errors[counter-1]*L2Errors[counter-1] + DErr*DErr);
178  DErrors.Append(DErr);
179  EnErrors.Append(err);
180  // Compute the rate of convergence by:
181  // rate = log (||u - u_h|| / ||u - u_{h/2}||)/log(2)
182  double val = (dcounter) ? log(DErrors[dcounter-1]/DErr)/log(2.0) : 0.0;
183  double eval = (dcounter) ? log(EnErrors[dcounter-1]/err)/log(2.0) : 0.0;
184  DRates.Append(val);
185  EnRates.Append(eval);
186  dcounter++;
187  MFEM_VERIFY(counter == dcounter,
188  "Number of added solutions and derivatives do not match")
189  }
190 }
191 
192 void ConvergenceStudy::Print(bool relative, std::ostream &out)
193 {
194  if (print_flag)
195  {
196  std::string title = (relative) ? "Relative " : "Absolute ";
197  out << "\n";
198  out << " -------------------------------------------" << "\n";
199  out << std::setw(21) << title << "L2 Error " << "\n";
200  out << " -------------------------------------------"
201  << "\n";
202  out << std::right<< std::setw(11)<< "DOFs "<< std::setw(13) << "Error ";
203  out << std::setw(15) << "Rate " << "\n";
204  out << " -------------------------------------------"
205  << "\n";
206  out << std::setprecision(4);
207  double d = (relative) ? CoeffNorm : 1.0;
208  for (int i =0; i<counter; i++)
209  {
210  out << std::right << std::setw(10)<< ndofs[i] << std::setw(16)
211  << std::scientific << L2Errors[i]/d << std::setw(13)
212  << std::fixed << L2Rates[i] << "\n";
213  }
214  out << "\n";
215  if (dcounter == counter)
216  {
217  std::string dname;
218  switch (cont_type)
219  {
220  case 0: dname = "Grad"; break;
221  case 1: dname = "Curl"; break;
222  case 2: dname = "Div"; break;
223  case 3: dname = "DG Grad"; break;
224  default: break;
225  }
226  out << " -------------------------------------------" << "\n";
227  out << std::setw(21) << title << dname << " Error " << "\n";
228  out << " -------------------------------------------" << "\n";
229  out << std::right<<std::setw(11)<< "DOFs "<< std::setw(13) << "Error";
230  out << std::setw(15) << "Rate " << "\n";
231  out << " -------------------------------------------"
232  << "\n";
233  out << std::setprecision(4);
234  d = (relative) ? CoeffDNorm : 1.0;
235  for (int i =0; i<dcounter; i++)
236  {
237  out << std::right << std::setw(10)<< ndofs[i] << std::setw(16)
238  << std::scientific << DErrors[i]/d << std::setw(13)
239  << std::fixed << DRates[i] << "\n";
240  }
241  out << "\n";
242  switch (cont_type)
243  {
244  case 0: dname = "H1"; break;
245  case 1: dname = "H(Curl)"; break;
246  case 2: dname = "H(Div)"; break;
247  case 3: dname = "DG H1"; break;
248  default: break;
249  }
250 
251  if (dcounter)
252  {
253  d = (relative) ?
254  sqrt(CoeffNorm*CoeffNorm + CoeffDNorm*CoeffDNorm):1.0;
255 
256  out << " -------------------------------------------" << "\n";
257  out << std::setw(21) << title << dname << " Error " << "\n";
258  out << " -------------------------------------------" << "\n";
259  out << std::right<< std::setw(11)<< "DOFs "<< std::setw(13);
260  out << "Error ";
261  out << std::setw(15) << "Rate " << "\n";
262  out << " -------------------------------------------"
263  << "\n";
264  out << std::setprecision(4);
265  for (int i =0; i<dcounter; i++)
266  {
267  out << std::right << std::setw(10)<< ndofs[i] << std::setw(16)
268  << std::scientific << EnErrors[i]/d << std::setw(13)
269  << std::fixed << EnRates[i] << "\n";
270  }
271  out << "\n";
272  }
273  if (cont_type == 3 && fcounter)
274  {
275  out << " -------------------------------------------" << "\n";
276  out << " DG Face Jump Error " << "\n";
277  out << " -------------------------------------------"
278  << "\n";
279  out << std::right<< std::setw(11)<< "DOFs "<< std::setw(13);
280  out << "Error ";
281  out << std::setw(15) << "Rate " << "\n";
282  out << " -------------------------------------------"
283  << "\n";
284  out << std::setprecision(4);
285  for (int i =0; i<fcounter; i++)
286  {
287  out << std::right << std::setw(10)<< ndofs[i] << std::setw(16)
288  << std::scientific << DGFaceErrors[i] << std::setw(13)
289  << std::fixed << DGFaceRates[i] << "\n";
290  }
291  out << "\n";
292  }
293  }
294  }
295 }
296 
297 } // namespace mfem
Field is discontinuous across element interfaces.
Definition: fe_coll.hpp:48
Class for an integration rule - an Array of IntegrationPoint.
Definition: intrules.hpp:90
Class for grid function - Vector with associated FE space.
Definition: gridfunc.hpp:30
const IntegrationRule & Get(int GeomType, int Order)
Returns an integration rule for given GeomType and Order.
Definition: intrules.cpp:915
Base class for vector Coefficients that optionally depend on time and space.
double ComputeLpNorm(double p, Coefficient &coeff, Mesh &mesh, const IntegrationRule *irs[])
Compute the Lp norm of a function f. .
Mesh * GetMesh() const
Returns the mesh.
Definition: fespace.hpp:314
FiniteElementSpace * FESpace()
Definition: gridfunc.hpp:595
OutStream err(std::cerr)
Global stream used by the library for standard error output. Initially it uses the same std::streambu...
Definition: globals.hpp:71
Base class Coefficients that optionally depend on space and time. These are used by the BilinearFormI...
Definition: coefficient.hpp:39
int GetOrder(int i) const
Returns the order of the i&#39;th finite element.
Definition: fespace.cpp:96
Field is continuous across element interfaces.
Definition: fe_coll.hpp:45
double ComputeGlobalLpNorm(double p, Coefficient &coeff, ParMesh &pmesh, const IntegrationRule *irs[])
Compute the global Lp norm of a function f. .
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
IntegrationRules IntRules(0, Quadrature1D::GaussLegendre)
A global object with all integration rules (defined in intrules.cpp)
Definition: intrules.hpp:378