MFEM  v4.6.0
Finite element discretization library
convection-diffusion.cpp
Go to the documentation of this file.
1 // Copyright (c) 2010-2023, 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 // MFEM Ultraweak DPG example for convection-diffusion
13 //
14 // Compile with: make convection-diffusion
15 //
16 // sample runs
17 // convection-diffusion -m ../../data/star.mesh -o 2 -ref 2 -theta 0.0 -eps 1e-1 -beta '2 3'
18 // convection-diffusion -m ../../data/beam-hex.mesh -o 2 -ref 2 -theta 0.0 -eps 1e0 -beta '1 0 2'
19 // convection-diffusion -m ../../data/inline-tri.mesh -o 3 -ref 2 -theta 0.0 -eps 1e-2 -beta '4 2' -sc
20 
21 // AMR runs
22 // convection-diffusion -o 3 -ref 5 -prob 1 -eps 1e-1 -theta 0.75
23 // convection-diffusion -o 2 -ref 9 -prob 1 -eps 1e-2 -theta 0.75
24 // convection-diffusion -o 3 -ref 9 -prob 1 -eps 1e-3 -theta 0.75 -sc
25 
26 // Description:
27 // This example code demonstrates the use of MFEM to define and solve
28 // the "ultraweak" (UW) DPG formulation for the convection-diffusion problem
29 
30 // - εΔu + ∇⋅(βu) = f, in Ω
31 // u = u_0, on ∂Ω
32 
33 // It solves the following kinds of problems
34 // (a) A manufactured solution where u_exact = sin(π * (x + y + z)).
35 // (b) The 2D Erickson-Johnson problem
36 
37 // The DPG UW deals with the First Order System
38 // - ∇⋅σ + ∇⋅(βu) = f, in Ω
39 // 1/ε σ - ∇u = 0, in Ω
40 // u = u_0, on ∂Ω
41 
42 // Ultraweak-DPG is obtained by integration by parts of both equations and the
43 // introduction of trace unknowns on the mesh skeleton
44 //
45 // u ∈ L²(Ω), σ ∈ (L²(Ω))ᵈⁱᵐ
46 // û ∈ H^1/2, σ̂ ∈ H^-1/2
47 // -(βu , ∇v) + (σ , ∇v) + < f̂ , v > = (f,v), ∀ v ∈ H¹(Ω)
48 // (u , ∇⋅τ) + 1/ε (σ , τ) + < û , τ⋅n > = 0, ∀ τ ∈ H(div,Ω)
49 // û = u_0 on ∂Ω
50 
51 // Note:
52 // f̂ := βu - σ, û := -u on the mesh skeleton
53 
54 // -------------------------------------------------------------
55 // | | u | σ | û | f̂ | RHS |
56 // -------------------------------------------------------------
57 // | v |-(βu , ∇v) | (σ , ∇v) | | < f̂ ,v > | (f,v) |
58 // | | | | | | |
59 // | τ | (u ,∇⋅τ) | 1/ε(σ , τ)| <û,τ⋅n> | | 0 |
60 
61 // where (v,τ) ∈ H¹(Ωₕ) × H(div,Ωₕ)
62 
63 // For more information see https://doi.org/10.1016/j.camwa.2013.06.010
64 
65 #include "mfem.hpp"
66 #include "util/weakform.hpp"
67 #include "../common/mfem-common.hpp"
68 #include <fstream>
69 #include <iostream>
70 
71 using namespace std;
72 using namespace mfem;
73 using namespace mfem::common;
74 
76 {
78  EJ // see https://doi.org/10.1016/j.camwa.2013.06.010
79 };
80 
83 double epsilon;
84 
85 double exact_u(const Vector & X);
86 void exact_gradu(const Vector & X, Vector & du);
87 double exact_laplacian_u(const Vector & X);
88 void exact_sigma(const Vector & X, Vector & sigma);
89 double exact_hatu(const Vector & X);
90 void exact_hatf(const Vector & X, Vector & hatf);
91 double f_exact(const Vector & X);
92 void setup_test_norm_coeffs(GridFunction & c1_gf, GridFunction & c2_gf);
93 
94 int main(int argc, char *argv[])
95 {
96  // 1. Parse command-line options.
97  const char *mesh_file = "../../data/inline-quad.mesh";
98  int order = 1;
99  int delta_order = 1;
100  int ref = 1;
101  bool visualization = true;
102  int iprob = 0;
103  double theta = 0.0;
104  bool static_cond = false;
105  epsilon = 1e0;
106 
107  OptionsParser args(argc, argv);
108  args.AddOption(&mesh_file, "-m", "--mesh",
109  "Mesh file to use.");
110  args.AddOption(&order, "-o", "--order",
111  "Finite element order (polynomial degree).");
112  args.AddOption(&delta_order, "-do", "--delta-order",
113  "Order enrichment for DPG test space.");
114  args.AddOption(&epsilon, "-eps", "--epsilon",
115  "Epsilon coefficient");
116  args.AddOption(&ref, "-ref", "--num-refinements",
117  "Number of uniform refinements");
118  args.AddOption(&theta, "-theta", "--theta",
119  "Theta parameter for AMR");
120  args.AddOption(&iprob, "-prob", "--problem", "Problem case"
121  " 0: manufactured, 1: Erickson-Johnson ");
122  args.AddOption(&beta, "-beta", "--beta",
123  "Vector Coefficient beta");
124  args.AddOption(&static_cond, "-sc", "--static-condensation", "-no-sc",
125  "--no-static-condensation", "Enable static condensation.");
126  args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
127  "--no-visualization",
128  "Enable or disable GLVis visualization.");
129  args.Parse();
130  if (!args.Good())
131  {
132  args.PrintUsage(cout);
133  return 1;
134  }
135 
136  if (iprob > 1) { iprob = 1; }
137  prob = (prob_type)iprob;
138 
139  if (prob == prob_type::EJ)
140  {
141  mesh_file = "../../data/inline-quad.mesh";
142  }
143 
144  Mesh mesh(mesh_file, 1, 1);
145  int dim = mesh.Dimension();
146  MFEM_VERIFY(dim > 1, "Dimension = 1 is not supported in this example");
147 
148  if (beta.Size() == 0)
149  {
150  beta.SetSize(dim);
151  beta = 0.0;
152  beta[0] = 1.;
153  }
154 
155  args.PrintOptions(cout);
156 
157  // Define spaces
158  enum TrialSpace
159  {
160  u_space = 0,
161  sigma_space = 1,
162  hatu_space = 2,
163  hatf_space = 3
164  };
165  enum TestSpace
166  {
167  v_space = 0,
168  tau_space = 1
169  };
170  // L2 space for u
171  FiniteElementCollection *u_fec = new L2_FECollection(order-1,dim);
172  FiniteElementSpace *u_fes = new FiniteElementSpace(&mesh,u_fec);
173 
174  // Vector L2 space for σ
175  FiniteElementCollection *sigma_fec = new L2_FECollection(order-1,dim);
176  FiniteElementSpace *sigma_fes = new FiniteElementSpace(&mesh,sigma_fec, dim);
177 
178  // H^1/2 space for û
179  FiniteElementCollection * hatu_fec = new H1_Trace_FECollection(order,dim);
180  FiniteElementSpace *hatu_fes = new FiniteElementSpace(&mesh,hatu_fec);
181 
182  // H^-1/2 space for σ̂
183  FiniteElementCollection * hatf_fec = new RT_Trace_FECollection(order-1,dim);
184  FiniteElementSpace *hatf_fes = new FiniteElementSpace(&mesh,hatf_fec);
185 
186  // testspace fe collections
187  int test_order = order+delta_order;
188  FiniteElementCollection * v_fec = new H1_FECollection(test_order, dim);
189  FiniteElementCollection * tau_fec = new RT_FECollection(test_order-1, dim);
190 
191  // Coefficients
192  ConstantCoefficient one(1.0);
193  ConstantCoefficient negone(-1.0);
195  ConstantCoefficient eps1(1./epsilon);
196  ConstantCoefficient negeps1(-1./epsilon);
198 
199  ConstantCoefficient negeps(-epsilon);
200  VectorConstantCoefficient betacoeff(beta);
201  Vector negbeta = beta; negbeta.Neg();
202  DenseMatrix bbt(beta.Size());
203  MultVVt(beta, bbt);
204  MatrixConstantCoefficient bbtcoeff(bbt);
205  VectorConstantCoefficient negbetacoeff(negbeta);
206 
209 
210  trial_fes.Append(u_fes);
211  trial_fes.Append(sigma_fes);
212  trial_fes.Append(hatu_fes);
213  trial_fes.Append(hatf_fes);
214  test_fec.Append(v_fec);
215  test_fec.Append(tau_fec);
216 
217  FiniteElementCollection *coeff_fec = new L2_FECollection(0,dim);
218  FiniteElementSpace *coeff_fes = new FiniteElementSpace(&mesh,coeff_fec);
219  GridFunction c1_gf, c2_gf;
220  GridFunctionCoefficient c1_coeff(&c1_gf);
221  GridFunctionCoefficient c2_coeff(&c2_gf);
222 
223  DPGWeakForm * a = new DPGWeakForm(trial_fes,test_fec);
224  a->StoreMatrices(true); // needed for residual calculation
225 
226  //-(βu , ∇v)
227  a->AddTrialIntegrator(new MixedScalarWeakDivergenceIntegrator(betacoeff),
228  TrialSpace::u_space, TestSpace::v_space);
229 
230  // (σ,∇ v)
231  a->AddTrialIntegrator(new TransposeIntegrator(new GradientIntegrator(one)),
232  TrialSpace::sigma_space, TestSpace::v_space);
233 
234  // (u ,∇⋅τ)
235  a->AddTrialIntegrator(new MixedScalarWeakGradientIntegrator(negone),
236  TrialSpace::u_space, TestSpace::tau_space);
237 
238  // 1/ε (σ,τ)
239  a->AddTrialIntegrator(new TransposeIntegrator(new VectorFEMassIntegrator(eps1)),
240  TrialSpace::sigma_space, TestSpace::tau_space);
241 
242  // <û,τ⋅n>
243  a->AddTrialIntegrator(new NormalTraceIntegrator,
244  TrialSpace::hatu_space, TestSpace::tau_space);
245 
246  // <f̂ ,v>
247  a->AddTrialIntegrator(new TraceIntegrator,
248  TrialSpace::hatf_space, TestSpace::v_space);
249 
250  // mesh dependent test norm
251  c1_gf.SetSpace(coeff_fes);
252  c2_gf.SetSpace(coeff_fes);
253  setup_test_norm_coeffs(c1_gf,c2_gf);
254 
255  // c1 (v,δv)
256  a->AddTestIntegrator(new MassIntegrator(c1_coeff),
257  TestSpace::v_space, TestSpace::v_space);
258  // ε (∇v,∇δv)
259  a->AddTestIntegrator(new DiffusionIntegrator(eps),
260  TestSpace::v_space, TestSpace::v_space);
261  // (β⋅∇v, β⋅∇δv)
262  a->AddTestIntegrator(new DiffusionIntegrator(bbtcoeff),
263  TestSpace::v_space, TestSpace::v_space);
264  // c2 (τ,δτ)
265  a->AddTestIntegrator(new VectorFEMassIntegrator(c2_coeff),
266  TestSpace::tau_space, TestSpace::tau_space);
267  // (∇⋅τ,∇⋅δτ)
268  a->AddTestIntegrator(new DivDivIntegrator(one),
269  TestSpace::tau_space, TestSpace::tau_space);
270 
272  a->AddDomainLFIntegrator(new DomainLFIntegrator(f),TestSpace::v_space);
273 
276  Array<int> elements_to_refine;
279 
280  GridFunction hatu_gf, hatf_gf;
281 
282  socketstream u_out;
283  socketstream sigma_out;
284 
285  double res0 = 0.;
286  double err0 = 0.;
287  int dof0 = 0; // init to suppress gcc warning
288  std::cout << "\n Ref |"
289  << " Dofs |"
290  << " L2 Error |"
291  << " Rate |"
292  << " Residual |"
293  << " Rate |" << endl;
294  std::cout << std::string(64,'-') << endl;
295 
296  if (static_cond) { a->EnableStaticCondensation(); }
297  for (int it = 0; it<=ref; it++)
298  {
299  a->Assemble();
300 
301  Array<int> ess_tdof_list_uhat;
302  Array<int> ess_tdof_list_fhat;
303  Array<int> ess_bdr_uhat;
304  Array<int> ess_bdr_fhat;
305  if (mesh.bdr_attributes.Size())
306  {
307  ess_bdr_uhat.SetSize(mesh.bdr_attributes.Max());
308  ess_bdr_fhat.SetSize(mesh.bdr_attributes.Max());
309  ess_bdr_uhat = 1; ess_bdr_fhat = 0;
310  if (prob == prob_type::EJ)
311  {
312  ess_bdr_uhat = 0;
313  ess_bdr_fhat = 1;
314  ess_bdr_uhat[1] = 1;
315  ess_bdr_fhat[1] = 0;
316  }
317  hatu_fes->GetEssentialTrueDofs(ess_bdr_uhat, ess_tdof_list_uhat);
318  hatf_fes->GetEssentialTrueDofs(ess_bdr_fhat, ess_tdof_list_fhat);
319  }
320 
321  // shift the ess_tdofs
322  int n = ess_tdof_list_uhat.Size();
323  int m = ess_tdof_list_fhat.Size();
324  Array<int> ess_tdof_list(n+m);
325  for (int j = 0; j < n; j++)
326  {
327  ess_tdof_list[j] = ess_tdof_list_uhat[j]
328  + u_fes->GetTrueVSize()
329  + sigma_fes->GetTrueVSize();
330  }
331  for (int j = 0; j < m; j++)
332  {
333  ess_tdof_list[j+n] = ess_tdof_list_fhat[j]
334  + u_fes->GetTrueVSize()
335  + sigma_fes->GetTrueVSize()
336  + hatu_fes->GetTrueVSize();
337  }
338 
339  Array<int> offsets(5);
340  offsets[0] = 0;
341  int dofs = 0;
342  for (int i = 0; i<trial_fes.Size(); i++)
343  {
344  offsets[i+1] = trial_fes[i]->GetVSize();
345  dofs += trial_fes[i]->GetTrueVSize();
346  }
347  offsets.PartialSum();
348 
349  BlockVector x(offsets); x = 0.0;
350  hatu_gf.MakeRef(hatu_fes,x.GetBlock(2),0);
351  hatf_gf.MakeRef(hatf_fes,x.GetBlock(3),0);
352  hatu_gf.ProjectBdrCoefficient(hatuex,ess_bdr_uhat);
353  hatf_gf.ProjectBdrCoefficientNormal(hatfex,ess_bdr_fhat);
354 
355  OperatorPtr Ah;
356  Vector X,B;
357  a->FormLinearSystem(ess_tdof_list,x,Ah,X,B);
358 
359  BlockMatrix * A = Ah.As<BlockMatrix>();
360 
362  M.owns_blocks = 1;
363  for (int i = 0 ; i < A->NumRowBlocks(); i++)
364  {
365  M.SetDiagonalBlock(i,new DSmoother(A->GetBlock(i,i)));
366  }
367 
368  CGSolver cg;
369  cg.SetRelTol(1e-8);
370  cg.SetMaxIter(20000);
371  cg.SetPrintLevel(0);
372  cg.SetPreconditioner(M);
373  cg.SetOperator(*A);
374  cg.Mult(B, X);
375 
376  a->RecoverFEMSolution(X,x);
377 
378  GridFunction u_gf, sigma_gf;
379  u_gf.MakeRef(u_fes,x.GetBlock(0),0);
380  sigma_gf.MakeRef(sigma_fes,x.GetBlock(1),0);
381 
382  double u_err = u_gf.ComputeL2Error(uex);
383  double sigma_err = sigma_gf.ComputeL2Error(sigmaex);
384  double L2Error = sqrt(u_err*u_err + sigma_err*sigma_err);
385 
386  Vector & residuals = a->ComputeResidual(x);
387  double residual = residuals.Norml2();
388 
389  double rate_err = (it) ? dim*log(err0/L2Error)/log((double)dof0/dofs) : 0.0;
390  double rate_res = (it) ? dim*log(res0/residual)/log((double)dof0/dofs) : 0.0;
391 
392  err0 = L2Error;
393  res0 = residual;
394  dof0 = dofs;
395 
396  std::ios oldState(nullptr);
397  oldState.copyfmt(std::cout);
398  std::cout << std::right << std::setw(5) << it << " | "
399  << std::setw(10) << dof0 << " | "
400  << std::setprecision(3)
401  << std::setw(10) << std::scientific << err0 << " | "
402  << std::setprecision(2)
403  << std::setw(6) << std::fixed << rate_err << " | "
404  << std::setprecision(3)
405  << std::setw(10) << std::scientific << res0 << " | "
406  << std::setprecision(2)
407  << std::setw(6) << std::fixed << rate_res << " | "
408  << std::endl;
409  std::cout.copyfmt(oldState);
410 
411  if (visualization)
412  {
413  const char * keys = (it == 0 && dim == 2) ? "jRcm\n" : nullptr;
414  char vishost[] = "localhost";
415  int visport = 19916;
416  VisualizeField(u_out,vishost, visport, u_gf,
417  "Numerical u", 0,0, 500, 500, keys);
418  VisualizeField(sigma_out,vishost, visport, sigma_gf,
419  "Numerical flux", 501,0,500, 500, keys);
420  }
421 
422  if (it == ref)
423  {
424  break;
425  }
426 
427  elements_to_refine.SetSize(0);
428  double max_resid = residuals.Max();
429  for (int iel = 0; iel<mesh.GetNE(); iel++)
430  {
431  if (residuals[iel] > theta * max_resid)
432  {
433  elements_to_refine.Append(iel);
434  }
435  }
436 
437  mesh.GeneralRefinement(elements_to_refine,1,1);
438  for (int i =0; i<trial_fes.Size(); i++)
439  {
440  trial_fes[i]->Update(false);
441  }
442  a->Update();
443 
444  coeff_fes->Update();
445  c1_gf.Update();
446  c2_gf.Update();
447  setup_test_norm_coeffs(c1_gf,c2_gf);
448  }
449 
450  delete coeff_fes;
451  delete coeff_fec;
452  delete a;
453  delete tau_fec;
454  delete v_fec;
455  delete hatf_fes;
456  delete hatf_fec;
457  delete hatu_fes;
458  delete hatu_fec;
459  delete sigma_fes;
460  delete sigma_fec;
461  delete u_fec;
462  delete u_fes;
463 
464  return 0;
465 }
466 
467 double exact_u(const Vector & X)
468 {
469  double x = X[0];
470  double y = X[1];
471  double z = 0.;
472  if (X.Size() == 3) { z = X[2]; }
473  switch (prob)
474  {
475  case EJ:
476  {
477  double alpha = sqrt(1. + 4. * epsilon * epsilon * M_PI * M_PI);
478  double r1 = (1. + alpha) / (2.*epsilon);
479  double r2 = (1. - alpha) / (2.*epsilon);
480  double denom = exp(-r2) - exp(-r1);
481 
482  double g1 = exp(r2*(x-1.));
483  double g2 = exp(r1*(x-1.));
484  double g = g1-g2;
485 
486  return g * cos(M_PI * y)/denom;
487  }
488  break;
489  default:
490  {
491  double alpha = M_PI * (x + y + z);
492  return sin(alpha);
493  }
494  break;
495  }
496 }
497 
498 void exact_gradu(const Vector & X, Vector & du)
499 {
500  double x = X[0];
501  double y = X[1];
502  double z = 0.;
503  if (X.Size() == 3) { z = X[2]; }
504  du.SetSize(X.Size());
505  du = 0.;
506  switch (prob)
507  {
508  case EJ:
509  {
510  double alpha = sqrt(1. + 4. * epsilon * epsilon * M_PI * M_PI);
511  double r1 = (1. + alpha) / (2.*epsilon);
512  double r2 = (1. - alpha) / (2.*epsilon);
513  double denom = exp(-r2) - exp(-r1);
514 
515  double g1 = exp(r2*(x-1.));
516  double g1_x = r2*g1;
517  double g2 = exp(r1*(x-1.));
518  double g2_x = r1*g2;
519  double g = g1-g2;
520  double g_x = g1_x - g2_x;
521 
522  du[0] = g_x * cos(M_PI * y)/denom;
523  du[1] = -M_PI * g * sin(M_PI*y)/denom;
524  }
525  break;
526  default:
527  {
528  double alpha = M_PI * (x + y + z);
529  du.SetSize(X.Size());
530  for (int i = 0; i<du.Size(); i++)
531  {
532  du[i] = M_PI * cos(alpha);
533  }
534  }
535  break;
536  }
537 }
538 
539 double exact_laplacian_u(const Vector & X)
540 {
541  double x = X[0];
542  double y = X[1];
543  double z = 0.;
544  if (X.Size() == 3) { z = X[2]; }
545 
546  switch (prob)
547  {
548  case EJ:
549  {
550  double alpha = sqrt(1. + 4. * epsilon * epsilon * M_PI * M_PI);
551  double r1 = (1. + alpha) / (2.*epsilon);
552  double r2 = (1. - alpha) / (2.*epsilon);
553  double denom = exp(-r2) - exp(-r1);
554 
555  double g1 = exp(r2*(x-1.));
556  double g1_x = r2*g1;
557  double g1_xx = r2*g1_x;
558  double g2 = exp(r1*(x-1.));
559  double g2_x = r1*g2;
560  double g2_xx = r1*g2_x;
561  double g = g1-g2;
562  double g_xx = g1_xx - g2_xx;
563 
564  double u = g * cos(M_PI * y)/denom;
565  double u_xx = g_xx * cos(M_PI * y)/denom;
566  double u_yy = -M_PI * M_PI * u;
567  return u_xx + u_yy;
568  }
569  break;
570  default:
571  {
572  double alpha = M_PI * (x + y + z);
573  double u = sin(alpha);
574  return -M_PI*M_PI * u * X.Size();
575  }
576  break;
577  }
578 }
579 
580 void exact_sigma(const Vector & X, Vector & sigma)
581 {
582  // σ = ε ∇ u
583  exact_gradu(X,sigma);
584  sigma *= epsilon;
585 }
586 
587 double exact_hatu(const Vector & X)
588 {
589  return -exact_u(X);
590 }
591 
592 void exact_hatf(const Vector & X, Vector & hatf)
593 {
594  Vector sigma;
595  exact_sigma(X,sigma);
596  double u = exact_u(X);
597  hatf.SetSize(X.Size());
598  for (int i = 0; i<hatf.Size(); i++)
599  {
600  hatf[i] = beta[i] * u - sigma[i];
601  }
602 }
603 
604 double f_exact(const Vector & X)
605 {
606  // f = - εΔu + ∇⋅(βu)
607  Vector du;
608  exact_gradu(X,du);
609  double d2u = exact_laplacian_u(X);
610 
611  double s = 0;
612  for (int i = 0; i<du.Size(); i++)
613  {
614  s += beta[i] * du[i];
615  }
616  return -epsilon * d2u + s;
617 }
618 
620 {
621  Array<int> vdofs;
622  FiniteElementSpace * fes = c1_gf.FESpace();
623  Mesh * mesh = fes->GetMesh();
624  for (int i = 0; i < mesh->GetNE(); i++)
625  {
626  double volume = mesh->GetElementVolume(i);
627  double c1 = min(epsilon/volume, 1.);
628  double c2 = min(1./epsilon, 1./volume);
629  fes->GetElementDofs(i,vdofs);
630  c1_gf.SetSubVector(vdofs,c1);
631  c2_gf.SetSubVector(vdofs,c2);
632  }
633 }
virtual double ComputeL2Error(Coefficient *exsol[], const IntegrationRule *irs[]=NULL, const Array< int > *elems=NULL) const
Definition: gridfunc.cpp:2786
A matrix coefficient that is constant in space and time.
void SetSubVector(const Array< int > &dofs, const double value)
Set the entries listed in dofs to the given value.
Definition: vector.cpp:605
Class for domain integration L(v) := (f, v)
Definition: lininteg.hpp:108
int visport
Conjugate gradient method.
Definition: solvers.hpp:493
Class for grid function - Vector with associated FE space.
Definition: gridfunc.hpp:30
virtual void Update(bool want_transform=true)
Reflect changes in the mesh: update number of DOFs, etc. Also, calculate GridFunction transformation ...
Definition: fespace.cpp:3402
Data type for scaled Jacobi-type smoother of sparse matrix.
A class to handle Vectors in a block fashion.
Definition: blockvector.hpp:30
A coefficient that is constant across space and time.
Definition: coefficient.hpp:84
void PrintOptions(std::ostream &out) const
Print the options.
Definition: optparser.cpp:331
int Dimension() const
Dimension of the reference space used within the elements.
Definition: mesh.hpp:1020
double exact_u(const Vector &X)
void SetSize(int s)
Resize the vector to size s.
Definition: vector.hpp:517
void PrintUsage(std::ostream &out) const
Print the usage message.
Definition: optparser.cpp:462
Vector coefficient that is constant in space and time.
Pointer to an Operator of a specified type.
Definition: handle.hpp:33
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
Definition: solvers.cpp:718
int Size() const
Returns the size of the vector.
Definition: vector.hpp:197
Coefficient defined by a GridFunction. This coefficient is mesh dependent.
T Max() const
Find the maximal element in the array, using the comparison operator < for class T.
Definition: array.cpp:68
Data type dense matrix using column-major storage.
Definition: densemat.hpp:23
double exact_laplacian_u(const Vector &X)
virtual void GetEssentialTrueDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list, int component=-1)
Get a list of essential true dofs, ess_tdof_list, corresponding to the boundary attributes marked in ...
Definition: fespace.cpp:587
int main(int argc, char *argv[])
bool Good() const
Return true if the command line options were parsed successfully.
Definition: optparser.hpp:159
Vector beta
STL namespace.
(Q div u, div v) for RT elements
prob_type prob
void exact_gradu(const Vector &X, Vector &du)
void exact_hatf(const Vector &X, Vector &hatf)
void ProjectBdrCoefficientNormal(VectorCoefficient &vcoeff, Array< int > &bdr_attr)
Definition: gridfunc.cpp:2694
int NumRowBlocks() const
Return the number of row blocks.
Definition: blockmatrix.hpp:46
virtual void SetPrintLevel(int print_lvl)
Legacy method to set the level of verbosity of the solver output.
Definition: solvers.cpp:71
void Parse()
Parse the command-line options. Note that this function expects all the options provided through the ...
Definition: optparser.cpp:151
A class to handle Block diagonal preconditioners in a matrix-free implementation. ...
char vishost[]
int Append(const T &el)
Append element &#39;el&#39; to array, resize if necessary.
Definition: array.hpp:759
Arbitrary order "H^{1/2}-conforming" trace finite elements defined on the interface between mesh elem...
Definition: fe_coll.hpp:319
void SetMaxIter(int max_it)
Definition: solvers.hpp:201
double f_exact(const Vector &X)
virtual void MakeRef(FiniteElementSpace *f, double *v)
Make the GridFunction reference external data on a new FiniteElementSpace.
Definition: gridfunc.cpp:203
Array< int > & RowOffsets()
Return the row offsets for block starts.
Definition: blockmatrix.hpp:57
FiniteElementSpace * FESpace()
Definition: gridfunc.hpp:691
void setup_test_norm_coeffs(GridFunction &c1_gf, GridFunction &c2_gf)
Arbitrary order H(div)-conforming Raviart-Thomas finite elements.
Definition: fe_coll.hpp:380
void exact_sigma(const Vector &X, Vector &sigma)
A general vector function coefficient.
virtual int GetTrueVSize() const
Return the number of vector true (conforming) dofs.
Definition: fespace.hpp:712
Mesh * GetMesh() const
Returns the mesh.
Definition: fespace.hpp:555
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
Definition: mesh.hpp:275
void SetRelTol(double rtol)
Definition: solvers.hpp:199
void MultVVt(const Vector &v, DenseMatrix &vvt)
Make a matrix from a vector V.Vt.
Definition: densemat.cpp:3116
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition: fespace.hpp:219
Collection of finite elements from the same family in multiple dimensions. This class is used to matc...
Definition: fe_coll.hpp:26
void AddOption(bool *var, const char *enable_short_name, const char *enable_long_name, const char *disable_short_name, const char *disable_long_name, const char *description, bool required=false)
Add a boolean option and set &#39;var&#39; to receive the value. Enable/disable tags are used to set the bool...
Definition: optparser.hpp:82
void SetSize(int nsize)
Change the logical size of the array, keep existing entries.
Definition: array.hpp:687
void PartialSum()
Fill the entries of the array with the cumulative sum of the entries.
Definition: array.cpp:103
virtual DofTransformation * GetElementDofs(int elem, Array< int > &dofs) const
Returns indices of degrees of freedom of element &#39;elem&#39;. The returned indices are offsets into an ldo...
Definition: fespace.cpp:2731
double Max() const
Returns the maximal element of the vector.
Definition: vector.cpp:925
virtual void Update()
Transform by the Space UpdateMatrix (e.g., on Mesh change).
Definition: gridfunc.cpp:164
int GetNE() const
Returns number of elements.
Definition: mesh.hpp:1086
double a
Definition: lissajous.cpp:41
OpType * As() const
Return the Operator pointer statically cast to a specified OpType. Similar to the method Get()...
Definition: handle.hpp:104
int dim
Definition: ex24.cpp:53
SparseMatrix & GetBlock(int i, int j)
Return a reference to block (i,j). Reference may be invalid if Aij(i,j) == NULL.
Definition: blockmatrix.cpp:82
double Norml2() const
Returns the l2 norm of the vector.
Definition: vector.cpp:835
virtual void SetSpace(FiniteElementSpace *f)
Associate a new FiniteElementSpace with the GridFunction.
Definition: gridfunc.cpp:195
int Size() const
Return the logical size of the array.
Definition: array.hpp:141
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
Definition: solvers.hpp:507
const double alpha
Definition: ex15.cpp:369
A general function coefficient.
Vector data type.
Definition: vector.hpp:58
Class representing the DPG weak formulation. Given the variational formulation a(u,v) = b(v), (or A u = b, where <Au,v> = a(u,v)) this class forms the DPG linear system A^T G^-1 A u = A^T G^-1 b This system results from the minimum residual formulation u = argmin_w ||G^-1(b - Aw)||. Here G is a symmetric positive definite matrix resulting from the discretization of the Riesz operator on the test space. Since the test space is broken (discontinuous), G is defined and inverted element-wise and the assembly of the global system is performed in the same manner as the standard FEM method. Note that DPGWeakForm can handle multiple Finite Element spaces.
Definition: weakform.hpp:33
double epsilon
Arbitrary order "H^{-1/2}-conforming" face finite elements defined on the interface between mesh elem...
Definition: fe_coll.hpp:434
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
Definition: solvers.cpp:173
Arbitrary order H1-conforming (continuous) finite elements.
Definition: fe_coll.hpp:259
double sigma(const Vector &x)
Definition: maxwell.cpp:102
RefCoord s[3]
double u(const Vector &xvec)
Definition: lor_mms.hpp:22
prob_type
Definition: ex25.cpp:146
void GeneralRefinement(const Array< Refinement > &refinements, int nonconforming=-1, int nc_limit=0)
Definition: mesh.cpp:9820
double GetElementVolume(int i)
Definition: mesh.cpp:119
double exact_hatu(const Vector &X)
void ProjectBdrCoefficient(Coefficient &coeff, Array< int > &attr)
Project a Coefficient on the GridFunction, modifying only DOFs on the boundary associated with the bo...
Definition: gridfunc.hpp:468
Vector & GetBlock(int i)
Get the i-th vector in the block.
Definition: blockvector.hpp:93
Arbitrary order "L2-conforming" discontinuous finite elements.
Definition: fe_coll.hpp:327
double f(const Vector &p)
void Neg()
(*this) = -(*this)
Definition: vector.cpp:301