MFEM  v3.1
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
ex8.cpp
Go to the documentation of this file.
1 // MFEM Example 8
2 //
3 // Compile with: make ex8
4 //
5 // Sample runs: ex8 -m ../data/square-disc.mesh
6 // ex8 -m ../data/star.mesh
7 // ex8 -m ../data/escher.mesh
8 // ex8 -m ../data/fichera.mesh
9 // ex8 -m ../data/square-disc-p2.vtk
10 // ex8 -m ../data/square-disc-p3.mesh
11 // ex8 -m ../data/star-surf.mesh -o 2
12 // ex8 -m ../data/mobius-strip.mesh
13 //
14 // Description: This example code demonstrates the use of the Discontinuous
15 // Petrov-Galerkin (DPG) method in its primal 2x2 block form as a
16 // simple finite element discretization of the Laplace problem
17 // -Delta u = f with homogeneous Dirichlet boundary conditions. We
18 // use high-order continuous trial space, a high-order interfacial
19 // (trace) space, and a high-order discontinuous test space
20 // defining a local dual (H^{-1}) norm.
21 //
22 // We use the primal form of DPG, see "A primal DPG method without
23 // a first-order reformulation", Demkowicz and Gopalakrishnan, CAM
24 // 2013, DOI:10.1016/j.camwa.2013.06.029.
25 //
26 // The example highlights the use of interfacial (trace) finite
27 // elements and spaces, trace face integrators and the definition
28 // of block operators and preconditioners.
29 //
30 // We recommend viewing examples 1-5 before viewing this example.
31 
32 #include "mfem.hpp"
33 #include <fstream>
34 #include <iostream>
35 
36 using namespace std;
37 using namespace mfem;
38 
39 int main(int argc, char *argv[])
40 {
41  // 1. Parse command-line options.
42  const char *mesh_file = "../data/star.mesh";
43  int order = 1;
44  bool visualization = 1;
45 
46  OptionsParser args(argc, argv);
47  args.AddOption(&mesh_file, "-m", "--mesh",
48  "Mesh file to use.");
49  args.AddOption(&order, "-o", "--order",
50  "Finite element order (polynomial degree).");
51  args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
52  "--no-visualization",
53  "Enable or disable GLVis visualization.");
54  args.Parse();
55  if (!args.Good())
56  {
57  args.PrintUsage(cout);
58  return 1;
59  }
60  args.PrintOptions(cout);
61 
62  // 2. Read the mesh from the given mesh file. We can handle triangular,
63  // quadrilateral, tetrahedral, hexahedral, surface and volume meshes with
64  // the same code.
65  Mesh *mesh;
66  ifstream imesh(mesh_file);
67  if (!imesh)
68  {
69  cerr << "\nCan not open mesh file: " << mesh_file << '\n' << endl;
70  return 2;
71  }
72  mesh = new Mesh(imesh, 1, 1);
73  imesh.close();
74  int dim = mesh->Dimension();
75 
76  // 3. Refine the mesh to increase the resolution. In this example we do
77  // 'ref_levels' of uniform refinement. We choose 'ref_levels' to be the
78  // largest number that gives a final mesh with no more than 10,000
79  // elements.
80  {
81  int ref_levels =
82  (int)floor(log(10000./mesh->GetNE())/log(2.)/dim);
83  for (int l = 0; l < ref_levels; l++)
84  {
85  mesh->UniformRefinement();
86  }
87  }
88 
89  // 4. Define the trial, interfacial (trace) and test DPG spaces:
90  // - The trial space, x0_space, contains the non-interfacial unknowns and
91  // has the essential BC.
92  // - The interfacial space, xhat_space, contains the interfacial unknowns
93  // and does not have essential BC.
94  // - The test space, test_space, is an enriched space where the enrichment
95  // degree may depend on the spatial dimension of the domain, the type of
96  // the mesh and the trial space order.
97  int trial_order = order;
98  int trace_order = order - 1;
99  int test_order = order; // reduced order, full order is (order + dim - 1)
100  if (dim == 2 && (order%2 == 0 || (mesh->MeshGenerator() & 2 && order > 1)))
101  {
102  test_order++;
103  }
104  if (test_order < trial_order)
105  cerr << "Warning, test space not enriched enough to handle primal"
106  << " trial space\n";
107 
108  FiniteElementCollection *x0_fec, *xhat_fec, *test_fec;
109 
110  x0_fec = new H1_FECollection(trial_order, dim);
111  xhat_fec = new RT_Trace_FECollection(trace_order, dim);
112  test_fec = new L2_FECollection(test_order, dim);
113 
114  FiniteElementSpace *x0_space = new FiniteElementSpace(mesh, x0_fec);
115  FiniteElementSpace *xhat_space = new FiniteElementSpace(mesh, xhat_fec);
116  FiniteElementSpace *test_space = new FiniteElementSpace(mesh, test_fec);
117 
118  // 5. Define the block structure of the problem, by creating the offset
119  // variables. Also allocate two BlockVector objects to store the solution
120  // and rhs.
121  enum {x0_var, xhat_var, NVAR};
122 
123  int s0 = x0_space->GetVSize();
124  int s1 = xhat_space->GetVSize();
125  int s_test = test_space->GetVSize();
126 
127  Array<int> offsets(NVAR+1);
128  offsets[0] = 0;
129  offsets[1] = s0;
130  offsets[2] = s0+s1;
131 
132  Array<int> offsets_test(2);
133  offsets_test[0] = 0;
134  offsets_test[1] = s_test;
135 
136  std::cout << "\nNumber of Unknowns:\n"
137  << " Trial space, X0 : " << s0
138  << " (order " << trial_order << ")\n"
139  << " Interface space, Xhat : " << s1
140  << " (order " << trace_order << ")\n"
141  << " Test space, Y : " << s_test
142  << " (order " << test_order << ")\n\n";
143 
144  BlockVector x(offsets), b(offsets);
145  x = 0.;
146 
147  // 6. Set up the linear form F(.) which corresponds to the right-hand side of
148  // the FEM linear system, which in this case is (f,phi_i) where f=1.0 and
149  // phi_i are the basis functions in the test finite element fespace.
150  ConstantCoefficient one(1.0);
151  LinearForm F(test_space);
153  F.Assemble();
154 
155  // 7. Set up the mixed bilinear form for the primal trial unknowns, B0,
156  // the mixed bilinear form for the interfacial unknowns, Bhat,
157  // the inverse stiffness matrix on the discontinuous test space, Sinv,
158  // and the stiffness matrix on the continuous trial space, S0.
159  Array<int> ess_bdr(mesh->bdr_attributes.Max());
160  ess_bdr = 1;
161 
162  MixedBilinearForm *B0 = new MixedBilinearForm(x0_space,test_space);
164  B0->Assemble();
165  B0->EliminateTrialDofs(ess_bdr, x.GetBlock(x0_var), F);
166  B0->Finalize();
167 
168  MixedBilinearForm *Bhat = new MixedBilinearForm(xhat_space,test_space);
170  Bhat->Assemble();
171  Bhat->Finalize();
172 
173  BilinearForm *Sinv = new BilinearForm(test_space);
174  SumIntegrator *Sum = new SumIntegrator;
175  Sum->AddIntegrator(new DiffusionIntegrator(one));
176  Sum->AddIntegrator(new MassIntegrator(one));
177  Sinv->AddDomainIntegrator(new InverseIntegrator(Sum));
178  Sinv->Assemble();
179  Sinv->Finalize();
180 
181  BilinearForm *S0 = new BilinearForm(x0_space);
183  S0->Assemble();
184  S0->EliminateEssentialBC(ess_bdr);
185  S0->Finalize();
186 
187  SparseMatrix &matB0 = B0->SpMat();
188  SparseMatrix &matBhat = Bhat->SpMat();
189  SparseMatrix &matSinv = Sinv->SpMat();
190  SparseMatrix &matS0 = S0->SpMat();
191 
192  // 8. Set up the 1x2 block Least Squares DPG operator, B = [B0 Bhat],
193  // the normal equation operator, A = B^t Sinv B, and
194  // the normal equation right-hand-size, b = B^t Sinv F.
195  BlockOperator B(offsets_test, offsets);
196  B.SetBlock(0,0,&matB0);
197  B.SetBlock(0,1,&matBhat);
198  RAPOperator A(B, matSinv, B);
199  {
200  Vector SinvF(s_test);
201  matSinv.Mult(F,SinvF);
202  B.MultTranspose(SinvF, b);
203  }
204 
205  // 9. Set up a block-diagonal preconditioner for the 2x2 normal equation
206  //
207  // [ S0^{-1} 0 ]
208  // [ 0 Shat^{-1} ] Shat = (Bhat^T Sinv Bhat)
209  //
210  // corresponding to the primal (x0) and interfacial (xhat) unknowns.
211  SparseMatrix * Shat = RAP(matBhat, matSinv, matBhat);
212 
213 #ifndef MFEM_USE_SUITESPARSE
214  const double prec_rtol = 1e-3;
215  const int prec_maxit = 200;
216  CGSolver *S0inv = new CGSolver;
217  S0inv->SetOperator(matS0);
218  S0inv->SetPrintLevel(-1);
219  S0inv->SetRelTol(prec_rtol);
220  S0inv->SetMaxIter(prec_maxit);
221  CGSolver *Shatinv = new CGSolver;
222  Shatinv->SetOperator(*Shat);
223  Shatinv->SetPrintLevel(-1);
224  Shatinv->SetRelTol(prec_rtol);
225  Shatinv->SetMaxIter(prec_maxit);
226  // Disable 'iterative_mode' when using CGSolver (or any IterativeSolver) as
227  // a preconditioner:
228  S0inv->iterative_mode = false;
229  Shatinv->iterative_mode = false;
230 #else
231  Operator *S0inv = new UMFPackSolver(matS0);
232  Operator *Shatinv = new UMFPackSolver(*Shat);
233 #endif
234 
235  BlockDiagonalPreconditioner P(offsets);
236  P.SetDiagonalBlock(0, S0inv);
237  P.SetDiagonalBlock(1, Shatinv);
238 
239  // 10. Solve the normal equation system using the PCG iterative solver.
240  // Check the weighted norm of residual for the DPG least square problem.
241  // Wrap the primal variable in a GridFunction for visualization purposes.
242  PCG(A, P, b, x, 1, 200, 1e-12, 0.0);
243 
244  {
245  Vector LSres(s_test);
246  B.Mult(x, LSres);
247  LSres -= F;
248  double res = sqrt(matSinv.InnerProduct(LSres, LSres));
249  cout << "\n|| B0*x0 + Bhat*xhat - F ||_{S^-1} = " << res << endl;
250  }
251 
252  GridFunction x0;
253  x0.Update(x0_space, x.GetBlock(x0_var), 0);
254 
255  // 11. Save the refined mesh and the solution. This output can be viewed
256  // later using GLVis: "glvis -m refined.mesh -g sol.gf".
257  {
258  ofstream mesh_ofs("refined.mesh");
259  mesh_ofs.precision(8);
260  mesh->Print(mesh_ofs);
261  ofstream sol_ofs("sol.gf");
262  sol_ofs.precision(8);
263  x0.Save(sol_ofs);
264  }
265 
266  // 12. Send the solution by socket to a GLVis server.
267  if (visualization)
268  {
269  char vishost[] = "localhost";
270  int visport = 19916;
271  socketstream sol_sock(vishost, visport);
272  sol_sock.precision(8);
273  sol_sock << "solution\n" << *mesh << x0 << flush;
274  }
275 
276  // 13. Free the used memory.
277  delete S0inv;
278  delete Shatinv;
279  delete Shat;
280  delete Bhat;
281  delete B0;
282  delete S0;
283  delete Sinv;
284  delete test_space;
285  delete test_fec;
286  delete xhat_space;
287  delete xhat_fec;
288  delete x0_space;
289  delete x0_fec;
290  delete mesh;
291 
292  return 0;
293 }
Class for domain integration L(v) := (f, v)
Definition: lininteg.hpp:47
int GetVSize() const
Definition: fespace.hpp:164
Conjugate gradient method.
Definition: solvers.hpp:109
Class for grid function - Vector with associated FE space.
Definition: gridfunc.hpp:27
Integrator defining a sum of multiple Integrators.
Definition: bilininteg.hpp:150
void Assemble(int skip_zeros=1)
Assembles the form i.e. sums over all domain/bdr integrators.
Subclass constant coefficient.
Definition: coefficient.hpp:57
void Assemble()
Assembles the linear form i.e. sums over all domain/bdr integrators.
Definition: linearform.cpp:34
virtual void MultTranspose(const Vector &x, Vector &y) const
Action of the transpose operator.
int GetNE() const
Returns number of elements.
Definition: mesh.hpp:454
void AddIntegrator(BilinearFormIntegrator *integ)
Definition: bilininteg.hpp:160
void SetBlock(int iRow, int iCol, Operator *op)
Add a block op in the block-entry (iblock, jblock).
bool iterative_mode
If true, use the second argument of Mult as an initial guess.
Definition: operator.hpp:106
virtual void Finalize(int skip_zeros=1)
Finalizes the matrix initialization.
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
Direct sparse solver using UMFPACK.
Definition: solvers.hpp:332
Integrator that inverts the matrix assembled by another integrator.
Definition: bilininteg.hpp:132
virtual void Save(std::ostream &out) const
Save the GridFunction to an output stream.
Definition: gridfunc.cpp:2138
int dim
Definition: ex3.cpp:48
void SetPrintLevel(int print_lvl)
Definition: solvers.cpp:71
A class to handle Block diagonal preconditioners in a matrix-free implementation. ...
Data type sparse matrix.
Definition: sparsemat.hpp:38
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:7316
void SetMaxIter(int max_it)
Definition: solvers.hpp:61
HypreParMatrix * RAP(HypreParMatrix *A, HypreParMatrix *P)
Returns the matrix P^t * A * P.
Definition: hypre.cpp:1318
int MeshGenerator()
Truegrid or NetGen?
Definition: mesh.hpp:447
T Max() const
Definition: array.cpp:90
void Assemble(int skip_zeros=1)
The operator x -&gt; R*A*P*x.
Definition: operator.hpp:164
void EliminateTrialDofs(Array< int > &bdr_attr_is_ess, Vector &sol, Vector &rhs)
void PCG(const Operator &A, Solver &B, const Vector &b, Vector &x, int print_iter, int max_num_iter, double RTOLERANCE, double ATOLERANCE)
Preconditioned conjugate gradient method. (tolerances are squared)
Definition: solvers.cpp:440
int Dimension() const
Definition: mesh.hpp:475
void PrintUsage(std::ostream &out) const
Definition: optparser.cpp:434
void AddTraceFaceIntegrator(BilinearFormIntegrator *bfi)
virtual void Mult(const Vector &x, Vector &y) const
Operator application.
virtual void Print(std::ostream &out=std::cout) const
Print the mesh to the given stream using the default MFEM mesh format.
Definition: mesh.cpp:8332
void AddDomainIntegrator(LinearFormIntegrator *lfi)
Adds new Domain Integrator.
Definition: linearform.cpp:19
Array< int > bdr_attributes
Definition: mesh.hpp:140
int main(int argc, char *argv[])
Definition: ex1.cpp:45
void SetRelTol(double rtol)
Definition: solvers.hpp:59
Abstract finite element space.
Definition: fespace.hpp:62
virtual void Mult(const Vector &x, Vector &y) const
Matrix vector multiplication.
Definition: sparsemat.cpp:417
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)
Definition: optparser.hpp:74
const SparseMatrix & SpMat() const
virtual void Finalize(int skip_zeros=1)
Finalizes the matrix initialization.
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
Adds new Domain Integrator.
void PrintOptions(std::ostream &out) const
Definition: optparser.cpp:304
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
Definition: solvers.hpp:123
Vector data type.
Definition: vector.hpp:33
Arbitrary order H1-conforming (continuous) finite elements.
Definition: fe_coll.hpp:54
Class for linear form - Vector with associated FE space and LFIntegrators.
Definition: linearform.hpp:23
const SparseMatrix & SpMat() const
Returns a reference to the sparse matrix.
Abstract operator.
Definition: operator.hpp:21
A class to handle Block systems in a matrix-free implementation.
double InnerProduct(const Vector &x, const Vector &y) const
Compute y^t A x.
Definition: sparsemat.cpp:619
void EliminateEssentialBC(Array< int > &bdr_attr_is_ess, Vector &sol, Vector &rhs, int d=0)
void SetDiagonalBlock(int iblock, Operator *op)
Add a square block op in the block-entry (iblock, iblock).
Arbitrary order &quot;L2-conforming&quot; discontinuous finite elements.
Definition: fe_coll.hpp:95
bool Good() const
Definition: optparser.hpp:120