MFEM  v3.4
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 = new Mesh(mesh_file, 1, 1);
66  int dim = mesh->Dimension();
67 
68  // 3. Refine the mesh to increase the resolution. In this example we do
69  // 'ref_levels' of uniform refinement. We choose 'ref_levels' to be the
70  // largest number that gives a final mesh with no more than 10,000
71  // elements.
72  {
73  int ref_levels =
74  (int)floor(log(10000./mesh->GetNE())/log(2.)/dim);
75  for (int l = 0; l < ref_levels; l++)
76  {
77  mesh->UniformRefinement();
78  }
79  }
80 
81  // 4. Define the trial, interfacial (trace) and test DPG spaces:
82  // - The trial space, x0_space, contains the non-interfacial unknowns and
83  // has the essential BC.
84  // - The interfacial space, xhat_space, contains the interfacial unknowns
85  // and does not have essential BC.
86  // - The test space, test_space, is an enriched space where the enrichment
87  // degree may depend on the spatial dimension of the domain, the type of
88  // the mesh and the trial space order.
89  unsigned int trial_order = order;
90  unsigned int trace_order = order - 1;
91  unsigned int test_order = order; /* reduced order, full order is
92  (order + dim - 1) */
93  if (dim == 2 && (order%2 == 0 || (mesh->MeshGenerator() & 2 && order > 1)))
94  {
95  test_order++;
96  }
97  if (test_order < trial_order)
98  cerr << "Warning, test space not enriched enough to handle primal"
99  << " trial space\n";
100 
101  FiniteElementCollection *x0_fec, *xhat_fec, *test_fec;
102 
103  x0_fec = new H1_FECollection(trial_order, dim);
104  xhat_fec = new RT_Trace_FECollection(trace_order, dim);
105  test_fec = new L2_FECollection(test_order, dim);
106 
107  FiniteElementSpace *x0_space = new FiniteElementSpace(mesh, x0_fec);
108  FiniteElementSpace *xhat_space = new FiniteElementSpace(mesh, xhat_fec);
109  FiniteElementSpace *test_space = new FiniteElementSpace(mesh, test_fec);
110 
111  // 5. Define the block structure of the problem, by creating the offset
112  // variables. Also allocate two BlockVector objects to store the solution
113  // and rhs.
114  enum {x0_var, xhat_var, NVAR};
115 
116  int s0 = x0_space->GetVSize();
117  int s1 = xhat_space->GetVSize();
118  int s_test = test_space->GetVSize();
119 
120  Array<int> offsets(NVAR+1);
121  offsets[0] = 0;
122  offsets[1] = s0;
123  offsets[2] = s0+s1;
124 
125  Array<int> offsets_test(2);
126  offsets_test[0] = 0;
127  offsets_test[1] = s_test;
128 
129  std::cout << "\nNumber of Unknowns:\n"
130  << " Trial space, X0 : " << s0
131  << " (order " << trial_order << ")\n"
132  << " Interface space, Xhat : " << s1
133  << " (order " << trace_order << ")\n"
134  << " Test space, Y : " << s_test
135  << " (order " << test_order << ")\n\n";
136 
137  BlockVector x(offsets), b(offsets);
138  x = 0.;
139 
140  // 6. Set up the linear form F(.) which corresponds to the right-hand side of
141  // the FEM linear system, which in this case is (f,phi_i) where f=1.0 and
142  // phi_i are the basis functions in the test finite element fespace.
143  ConstantCoefficient one(1.0);
144  LinearForm F(test_space);
146  F.Assemble();
147 
148  // 7. Set up the mixed bilinear form for the primal trial unknowns, B0,
149  // the mixed bilinear form for the interfacial unknowns, Bhat,
150  // the inverse stiffness matrix on the discontinuous test space, Sinv,
151  // and the stiffness matrix on the continuous trial space, S0.
152  Array<int> ess_bdr(mesh->bdr_attributes.Max());
153  ess_bdr = 1;
154 
155  MixedBilinearForm *B0 = new MixedBilinearForm(x0_space,test_space);
157  B0->Assemble();
158  B0->EliminateTrialDofs(ess_bdr, x.GetBlock(x0_var), F);
159  B0->Finalize();
160 
161  MixedBilinearForm *Bhat = new MixedBilinearForm(xhat_space,test_space);
163  Bhat->Assemble();
164  Bhat->Finalize();
165 
166  BilinearForm *Sinv = new BilinearForm(test_space);
167  SumIntegrator *Sum = new SumIntegrator;
168  Sum->AddIntegrator(new DiffusionIntegrator(one));
169  Sum->AddIntegrator(new MassIntegrator(one));
170  Sinv->AddDomainIntegrator(new InverseIntegrator(Sum));
171  Sinv->Assemble();
172  Sinv->Finalize();
173 
174  BilinearForm *S0 = new BilinearForm(x0_space);
176  S0->Assemble();
177  S0->EliminateEssentialBC(ess_bdr);
178  S0->Finalize();
179 
180  SparseMatrix &matB0 = B0->SpMat();
181  SparseMatrix &matBhat = Bhat->SpMat();
182  SparseMatrix &matSinv = Sinv->SpMat();
183  SparseMatrix &matS0 = S0->SpMat();
184 
185  // 8. Set up the 1x2 block Least Squares DPG operator, B = [B0 Bhat],
186  // the normal equation operator, A = B^t Sinv B, and
187  // the normal equation right-hand-size, b = B^t Sinv F.
188  BlockOperator B(offsets_test, offsets);
189  B.SetBlock(0,0,&matB0);
190  B.SetBlock(0,1,&matBhat);
191  RAPOperator A(B, matSinv, B);
192  {
193  Vector SinvF(s_test);
194  matSinv.Mult(F,SinvF);
195  B.MultTranspose(SinvF, b);
196  }
197 
198  // 9. Set up a block-diagonal preconditioner for the 2x2 normal equation
199  //
200  // [ S0^{-1} 0 ]
201  // [ 0 Shat^{-1} ] Shat = (Bhat^T Sinv Bhat)
202  //
203  // corresponding to the primal (x0) and interfacial (xhat) unknowns.
204  SparseMatrix * Shat = RAP(matBhat, matSinv, matBhat);
205 
206 #ifndef MFEM_USE_SUITESPARSE
207  const double prec_rtol = 1e-3;
208  const int prec_maxit = 200;
209  CGSolver *S0inv = new CGSolver;
210  S0inv->SetOperator(matS0);
211  S0inv->SetPrintLevel(-1);
212  S0inv->SetRelTol(prec_rtol);
213  S0inv->SetMaxIter(prec_maxit);
214  CGSolver *Shatinv = new CGSolver;
215  Shatinv->SetOperator(*Shat);
216  Shatinv->SetPrintLevel(-1);
217  Shatinv->SetRelTol(prec_rtol);
218  Shatinv->SetMaxIter(prec_maxit);
219  // Disable 'iterative_mode' when using CGSolver (or any IterativeSolver) as
220  // a preconditioner:
221  S0inv->iterative_mode = false;
222  Shatinv->iterative_mode = false;
223 #else
224  Operator *S0inv = new UMFPackSolver(matS0);
225  Operator *Shatinv = new UMFPackSolver(*Shat);
226 #endif
227 
228  BlockDiagonalPreconditioner P(offsets);
229  P.SetDiagonalBlock(0, S0inv);
230  P.SetDiagonalBlock(1, Shatinv);
231 
232  // 10. Solve the normal equation system using the PCG iterative solver.
233  // Check the weighted norm of residual for the DPG least square problem.
234  // Wrap the primal variable in a GridFunction for visualization purposes.
235  PCG(A, P, b, x, 1, 200, 1e-12, 0.0);
236 
237  {
238  Vector LSres(s_test);
239  B.Mult(x, LSres);
240  LSres -= F;
241  double res = sqrt(matSinv.InnerProduct(LSres, LSres));
242  cout << "\n|| B0*x0 + Bhat*xhat - F ||_{S^-1} = " << res << endl;
243  }
244 
245  GridFunction x0;
246  x0.MakeRef(x0_space, x.GetBlock(x0_var), 0);
247 
248  // 11. Save the refined mesh and the solution. This output can be viewed
249  // later using GLVis: "glvis -m refined.mesh -g sol.gf".
250  {
251  ofstream mesh_ofs("refined.mesh");
252  mesh_ofs.precision(8);
253  mesh->Print(mesh_ofs);
254  ofstream sol_ofs("sol.gf");
255  sol_ofs.precision(8);
256  x0.Save(sol_ofs);
257  }
258 
259  // 12. Send the solution by socket to a GLVis server.
260  if (visualization)
261  {
262  char vishost[] = "localhost";
263  int visport = 19916;
264  socketstream sol_sock(vishost, visport);
265  sol_sock.precision(8);
266  sol_sock << "solution\n" << *mesh << x0 << flush;
267  }
268 
269  // 13. Free the used memory.
270  delete S0inv;
271  delete Shatinv;
272  delete Shat;
273  delete Bhat;
274  delete B0;
275  delete S0;
276  delete Sinv;
277  delete test_space;
278  delete test_fec;
279  delete xhat_space;
280  delete xhat_fec;
281  delete x0_space;
282  delete x0_fec;
283  delete mesh;
284 
285  return 0;
286 }
Class for domain integration L(v) := (f, v)
Definition: lininteg.hpp:93
int GetVSize() const
Return the number of vector dofs, i.e. GetNDofs() x GetVDim().
Definition: fespace.hpp:253
virtual void Print(std::ostream &out=mfem::out) const
Definition: mesh.hpp:1034
Conjugate gradient method.
Definition: solvers.hpp:111
Class for grid function - Vector with associated FE space.
Definition: gridfunc.hpp:27
Integrator defining a sum of multiple Integrators.
Definition: bilininteg.hpp:152
void Assemble(int skip_zeros=1)
Assembles the form i.e. sums over all domain/bdr integrators.
Subclass constant coefficient.
Definition: coefficient.hpp:57
HypreParMatrix * RAP(const HypreParMatrix *A, const HypreParMatrix *P)
Returns the matrix P^t * A * P.
Definition: hypre.cpp:1558
void Assemble()
Assembles the linear form i.e. sums over all domain/bdr integrators.
Definition: linearform.cpp:51
void EliminateTrialDofs(Array< int > &bdr_attr_is_ess, const Vector &sol, Vector &rhs)
virtual void MultTranspose(const Vector &x, Vector &y) const
Action of the transpose operator.
int GetNE() const
Returns number of elements.
Definition: mesh.hpp:618
void AddIntegrator(BilinearFormIntegrator *integ)
Definition: bilininteg.hpp:162
bool iterative_mode
If true, use the second argument of Mult() as an initial guess.
Definition: operator.hpp:272
virtual void Finalize(int skip_zeros=1)
Finalizes the matrix initialization.
int main(int argc, char *argv[])
Definition: ex1.cpp:45
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
Direct sparse solver using UMFPACK.
Definition: solvers.hpp:344
Integrator that inverts the matrix assembled by another integrator.
Definition: bilininteg.hpp:134
virtual void Save(std::ostream &out) const
Save the GridFunction to an output stream.
Definition: gridfunc.cpp:2349
int dim
Definition: ex3.cpp:47
void SetPrintLevel(int print_lvl)
Definition: solvers.cpp:72
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:6741
void EliminateEssentialBC(const Array< int > &bdr_attr_is_ess, const Vector &sol, Vector &rhs, DiagonalPolicy dpolicy=DIAG_ONE)
Eliminate essential boundary DOFs from the system.
void SetMaxIter(int max_it)
Definition: solvers.hpp:63
int MeshGenerator()
Get the mesh generator/type.
Definition: mesh.hpp:611
T Max() const
Find the maximal element in the array, using the comparison operator &lt; for class T.
Definition: array.cpp:109
void Assemble(int skip_zeros=1)
The operator x -&gt; R*A*P*x constructed through the actions of R^T, A and P.
Definition: operator.hpp:349
virtual void MakeRef(FiniteElementSpace *f, double *v)
Make the GridFunction reference external data on a new FiniteElementSpace.
Definition: gridfunc.cpp:182
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:457
int Dimension() const
Definition: mesh.hpp:645
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.
void AddDomainIntegrator(LinearFormIntegrator *lfi)
Adds new Domain Integrator.
Definition: linearform.cpp:19
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
Definition: mesh.hpp:174
void SetRelTol(double rtol)
Definition: solvers.hpp:61
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition: fespace.hpp:66
virtual void Mult(const Vector &x, Vector &y) const
Matrix vector multiplication.
Definition: sparsemat.cpp:486
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:125
Vector data type.
Definition: vector.hpp:48
Arbitrary order H1-conforming (continuous) finite elements.
Definition: fe_coll.hpp:79
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:690
void SetBlock(int iRow, int iCol, Operator *op, double c=1.0)
Add a block op in the block-entry (iblock, jblock).
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:128
bool Good() const
Definition: optparser.hpp:120