9 void ConvergenceStudy::Reset()
22 DGFaceErrors.SetSize(0);
23 DGFaceRates.SetSize(0);
30 bool norm_set =
false;
33 int order_quad = std::max(2, 2*order+1);
35 for (
int i=0; i < Geometry::NumGeom; ++i)
40 ParGridFunction *pgf =
dynamic_cast<ParGridFunction *
>(gf);
43 ParMesh *pmesh = pgf->ParFESpace()->GetParMesh();
70 void ConvergenceStudy::AddL2Error(GridFunction *gf,
71 Coefficient *scalar_u, VectorCoefficient *vector_u)
74 int dim = gf->FESpace()->GetMesh()->Dimension();
77 ParGridFunction *pgf =
dynamic_cast<ParGridFunction *
>(gf);
80 MPI_Comm comm = pgf->ParFESpace()->GetComm();
82 MPI_Comm_rank(comm, &rank);
84 if (rank==0) { print_flag = 1; }
85 tdofs = pgf->ParFESpace()->GlobalTrueVSize();
88 if (!tdofs) { tdofs = gf->FESpace()->GetTrueVSize(); }
93 L2Err = gf->ComputeL2Error(*scalar_u);
94 CoeffNorm = GetNorm(gf,scalar_u,
nullptr);
98 L2Err = gf->ComputeL2Error(*vector_u);
99 CoeffNorm = GetNorm(gf,
nullptr,vector_u);
103 MFEM_ABORT(
"Exact Solution Coefficient pointer is NULL");
105 L2Errors.Append(L2Err);
111 double num =
log(L2Errors[counter-1]/L2Err);
112 double den =
log((
double)ndofs[counter]/ndofs[counter-1]);
119 void ConvergenceStudy::AddGf(GridFunction *gf, Coefficient *scalar_u,
120 VectorCoefficient *grad,
121 Coefficient *ell_coeff,
122 JumpScaling jump_scaling)
124 cont_type = gf->FESpace()->FEColl()->GetContType();
128 "This constructor is intended for H1 or L2 Elements")
130 AddL2Error(gf,scalar_u,
nullptr);
131 int dim = gf->FESpace()->
GetMesh()->Dimension();
135 double GradErr = gf->ComputeGradError(grad);
136 DErrors.Append(GradErr);
138 sqrt(L2Errors[counter-1]*L2Errors[counter-1]+GradErr*GradErr);
139 EnErrors.Append(error);
146 double num =
log(DErrors[dcounter-1]/GradErr);
147 double den =
log((
double)ndofs[dcounter]/ndofs[dcounter-1]);
149 num =
log(EnErrors[dcounter-1]/error);
150 eval = dim * num/den;
153 EnRates.Append(eval);
154 CoeffDNorm = GetNorm(gf,
nullptr,grad);
156 MFEM_VERIFY(counter == dcounter,
157 "Number of added solutions and derivatives do not match")
162 double DGErr = gf->ComputeDGFaceJumpError(scalar_u,ell_coeff,jump_scaling);
163 DGFaceErrors.Append(DGErr);
169 double num =
log(DGFaceErrors[fcounter-1]/DGErr);
170 double den =
log((
double)ndofs[fcounter]/ndofs[fcounter-1]);
173 DGFaceRates.Append(val);
175 MFEM_VERIFY(fcounter == counter,
"Number of added solutions mismatch");
179 void ConvergenceStudy::AddGf(GridFunction *gf, VectorCoefficient *vector_u,
180 VectorCoefficient *curl, Coefficient *div)
182 cont_type = gf->FESpace()->FEColl()->GetContType();
184 AddL2Error(gf,
nullptr,vector_u);
185 int dim = gf->FESpace()->GetMesh()->Dimension();
187 bool derivative =
false;
190 DErr = gf->ComputeCurlError(curl);
191 CoeffDNorm = GetNorm(gf,
nullptr,curl);
196 DErr = gf->ComputeDivError(div);
198 CoeffDNorm = GetNorm(gf,div,
nullptr);
203 double error =
sqrt(L2Errors[counter-1]*L2Errors[counter-1] + DErr*DErr);
204 DErrors.Append(DErr);
205 EnErrors.Append(error);
212 double num =
log(DErrors[dcounter-1]/DErr);
213 double den =
log((
double)ndofs[dcounter]/ndofs[dcounter-1]);
215 num =
log(EnErrors[dcounter-1]/error);
216 eval = dim * num/den;
219 EnRates.Append(eval);
221 MFEM_VERIFY(counter == dcounter,
222 "Number of added solutions and derivatives do not match")
226 void ConvergenceStudy::Print(
bool relative, std::ostream &os)
230 std::string title = (relative) ?
"Relative " :
"Absolute ";
232 os <<
" -------------------------------------------" <<
"\n";
233 os << std::setw(21) << title <<
"L2 Error " <<
"\n";
234 os <<
" -------------------------------------------"
236 os << std::right<< std::setw(11)<<
"DOFs "<< std::setw(13) <<
"Error ";
237 os << std::setw(15) <<
"Rate " <<
"\n";
238 os <<
" -------------------------------------------"
240 os << std::setprecision(4);
241 double d = (relative) ? CoeffNorm : 1.0;
242 for (
int i =0; i<counter; i++)
244 os << std::right << std::setw(10)<< ndofs[i] << std::setw(16)
245 << std::scientific << L2Errors[i]/d << std::setw(13)
246 << std::fixed << L2Rates[i] <<
"\n";
249 if (dcounter == counter)
254 case 0: dname =
"Grad";
break;
255 case 1: dname =
"Curl";
break;
256 case 2: dname =
"Div";
break;
257 case 3: dname =
"DG Grad";
break;
260 os <<
" -------------------------------------------" <<
"\n";
261 os << std::setw(21) << title << dname <<
" Error " <<
"\n";
262 os <<
" -------------------------------------------" <<
"\n";
263 os << std::right<<std::setw(11)<<
"DOFs "<< std::setw(13) <<
"Error";
264 os << std::setw(15) <<
"Rate " <<
"\n";
265 os <<
" -------------------------------------------"
267 os << std::setprecision(4);
268 d = (relative) ? CoeffDNorm : 1.0;
269 for (
int i =0; i<dcounter; i++)
271 os << std::right << std::setw(10)<< ndofs[i] << std::setw(16)
272 << std::scientific << DErrors[i]/d << std::setw(13)
273 << std::fixed << DRates[i] <<
"\n";
278 case 0: dname =
"H1";
break;
279 case 1: dname =
"H(Curl)";
break;
280 case 2: dname =
"H(Div)";
break;
281 case 3: dname =
"DG H1";
break;
288 sqrt(CoeffNorm*CoeffNorm + CoeffDNorm*CoeffDNorm):1.0;
290 os <<
" -------------------------------------------" <<
"\n";
291 os << std::setw(21) << title << dname <<
" Error " <<
"\n";
292 os <<
" -------------------------------------------" <<
"\n";
293 os << std::right<< std::setw(11)<<
"DOFs "<< std::setw(13);
295 os << std::setw(15) <<
"Rate " <<
"\n";
296 os <<
" -------------------------------------------"
298 os << std::setprecision(4);
299 for (
int i =0; i<dcounter; i++)
301 os << std::right << std::setw(10)<< ndofs[i] << std::setw(16)
302 << std::scientific << EnErrors[i]/d << std::setw(13)
303 << std::fixed << EnRates[i] <<
"\n";
308 if (cont_type == 3 && fcounter)
310 os <<
" -------------------------------------------" <<
"\n";
311 os <<
" DG Face Jump Error " <<
"\n";
312 os <<
" -------------------------------------------"
314 os << std::right<< std::setw(11)<<
"DOFs "<< std::setw(13);
316 os << std::setw(15) <<
"Rate " <<
"\n";
317 os <<
" -------------------------------------------"
319 os << std::setprecision(4);
320 for (
int i =0; i<fcounter; i++)
322 os << std::right << std::setw(10)<< ndofs[i] << std::setw(16)
323 << std::scientific << DGFaceErrors[i] << std::setw(13)
324 << std::fixed << DGFaceRates[i] <<
"\n";
Field is discontinuous across element interfaces.
Class for an integration rule - an Array of IntegrationPoint.
Class for grid function - Vector with associated FE space.
const IntegrationRule & Get(int GeomType, int Order)
Returns an integration rule for given GeomType and Order.
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.
int GetMaxElementOrder() const
Return the maximum polynomial order.
FiniteElementSpace * FESpace()
Base class Coefficients that optionally depend on space and time. These are used by the BilinearFormI...
Field is continuous across element interfaces.
double ComputeGlobalLpNorm(double p, Coefficient &coeff, ParMesh &pmesh, const IntegrationRule *irs[])
Compute the global Lp norm of a function f. .
FDualNumber< tbase > log(const FDualNumber< tbase > &f)
log([dual number])
FDualNumber< tbase > sqrt(const FDualNumber< tbase > &f)
sqrt([dual number])
IntegrationRules IntRules(0, Quadrature1D::GaussLegendre)
A global object with all integration rules (defined in intrules.cpp)