MFEM  v3.1
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
ex8p.cpp
Go to the documentation of this file.
1 // MFEM Example 8 - Parallel Version
2 //
3 // Compile with: make ex8p
4 //
5 // Sample runs: mpirun -np 4 ex8p -m ../data/square-disc.mesh
6 // mpirun -np 4 ex8p -m ../data/star.mesh
7 // mpirun -np 4 ex8p -m ../data/escher.mesh
8 // mpirun -np 4 ex8p -m ../data/fichera.mesh
9 // mpirun -np 4 ex8p -m ../data/square-disc-p2.vtk
10 // mpirun -np 4 ex8p -m ../data/square-disc-p3.mesh
11 // mpirun -np 4 ex8p -m ../data/star-surf.mesh -o 2
12 //
13 // Description: This example code demonstrates the use of the Discontinuous
14 // Petrov-Galerkin (DPG) method in its primal 2x2 block form as a
15 // simple finite element discretization of the Laplace problem
16 // -Delta u = f with homogeneous Dirichlet boundary conditions. We
17 // use high-order continuous trial space, a high-order interfacial
18 // (trace) space, and a high-order discontinuous test space
19 // defining a local dual (H^{-1}) norm.
20 //
21 // We use the primal form of DPG, see "A primal DPG method without
22 // a first-order reformulation", Demkowicz and Gopalakrishnan, CAM
23 // 2013, DOI:10.1016/j.camwa.2013.06.029.
24 //
25 // The example highlights the use of interfacial (trace) finite
26 // elements and spaces, trace face integrators and the definition
27 // of block operators and preconditioners. The use of the ADS
28 // preconditioner from hypre for interfacially-reduced H(div)
29 // problems is also illustrated.
30 //
31 // We recommend viewing examples 1-5 before viewing this example.
32 
33 #include "mfem.hpp"
34 #include <fstream>
35 #include <iostream>
36 
37 using namespace std;
38 using namespace mfem;
39 
40 int main(int argc, char *argv[])
41 {
42  // 1. Initialize MPI.
43  int num_procs, myid;
44  MPI_Init(&argc, &argv);
45  MPI_Comm_size(MPI_COMM_WORLD, &num_procs);
46  MPI_Comm_rank(MPI_COMM_WORLD, &myid);
47 
48  // 2. Parse command-line options.
49  const char *mesh_file = "../data/star.mesh";
50  int order = 1;
51  bool visualization = 1;
52 
53  OptionsParser args(argc, argv);
54  args.AddOption(&mesh_file, "-m", "--mesh",
55  "Mesh file to use.");
56  args.AddOption(&order, "-o", "--order",
57  "Finite element order (polynomial degree).");
58  args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
59  "--no-visualization",
60  "Enable or disable GLVis visualization.");
61  args.Parse();
62  if (!args.Good())
63  {
64  if (myid == 0)
65  {
66  args.PrintUsage(cout);
67  }
68  MPI_Finalize();
69  return 1;
70  }
71  if (myid == 0)
72  {
73  args.PrintOptions(cout);
74  }
75 
76  // 3. Read the (serial) mesh from the given mesh file on all processors. We
77  // can handle triangular, quadrilateral, tetrahedral, hexahedral, surface
78  // and volume meshes with the same code.
79  Mesh *mesh;
80  ifstream imesh(mesh_file);
81  if (!imesh)
82  {
83  if (myid == 0)
84  {
85  cerr << "\nCan not open mesh file: " << mesh_file << '\n' << endl;
86  }
87  MPI_Finalize();
88  return 2;
89  }
90  mesh = new Mesh(imesh, 1, 1);
91  imesh.close();
92  int dim = mesh->Dimension();
93 
94  // 4. Refine the serial mesh on all processors to increase the resolution. In
95  // this example we do 'ref_levels' of uniform refinement. We choose
96  // 'ref_levels' to be the largest number that gives a final mesh with no
97  // more than 10,000 elements.
98  {
99  int ref_levels =
100  (int)floor(log(10000./mesh->GetNE())/log(2.)/dim);
101  for (int l = 0; l < ref_levels; l++)
102  {
103  mesh->UniformRefinement();
104  }
105  }
106 
107  // 5. Define a parallel mesh by a partitioning of the serial mesh. Refine
108  // this mesh further in parallel to increase the resolution. Once the
109  // parallel mesh is defined, the serial mesh can be deleted.
110  ParMesh *pmesh = new ParMesh(MPI_COMM_WORLD, *mesh);
111  delete mesh;
112  {
113  int par_ref_levels = 1;
114  for (int l = 0; l < par_ref_levels; l++)
115  {
116  pmesh->UniformRefinement();
117  }
118  }
119  pmesh->ReorientTetMesh();
120 
121  // 6. Define the trial, interfacial (trace) and test DPG spaces:
122  // - The trial space, x0_space, contains the non-interfacial unknowns and
123  // has the essential BC.
124  // - The interfacial space, xhat_space, contains the interfacial unknowns
125  // and does not have essential BC.
126  // - The test space, test_space, is an enriched space where the enrichment
127  // degree may depend on the spatial dimension of the domain, the type of
128  // the mesh and the trial space order.
129  int trial_order = order;
130  int trace_order = order - 1;
131  int test_order = order; // reduced order, full order is (order + dim - 1)
132  if (dim == 2 && (order%2 == 0 || (pmesh->MeshGenerator() & 2 && order > 1)))
133  {
134  test_order++;
135  }
136  if (test_order < trial_order)
137  if (myid == 0)
138  cerr << "Warning, test space not enriched enough to handle primal"
139  << " trial space\n";
140 
141  FiniteElementCollection *x0_fec, *xhat_fec, *test_fec;
142 
143  x0_fec = new H1_FECollection(trial_order, dim);
144  xhat_fec = new RT_Trace_FECollection(trace_order, dim);
145  test_fec = new L2_FECollection(test_order, dim);
146 
147  ParFiniteElementSpace *x0_space, *xhat_space, *test_space;
148 
149  x0_space = new ParFiniteElementSpace(pmesh, x0_fec);
150  xhat_space = new ParFiniteElementSpace(pmesh, xhat_fec);
151  test_space = new ParFiniteElementSpace(pmesh, test_fec);
152 
153  HYPRE_Int glob_true_s0 = x0_space->GlobalTrueVSize();
154  HYPRE_Int glob_true_s1 = xhat_space->GlobalTrueVSize();
155  HYPRE_Int glob_true_s_test = test_space->GlobalTrueVSize();
156  if (myid == 0)
157  {
158  cout << "\nNumber of Unknowns:\n"
159  << " Trial space, X0 : " << glob_true_s0
160  << " (order " << trial_order << ")\n"
161  << " Interface space, Xhat : " << glob_true_s1
162  << " (order " << trace_order << ")\n"
163  << " Test space, Y : " << glob_true_s_test
164  << " (order " << test_order << ")\n\n";
165  }
166 
167  // 7. Set up the linear form F(.) which corresponds to the right-hand side of
168  // the FEM linear system, which in this case is (f,phi_i) where f=1.0 and
169  // phi_i are the basis functions in the test finite element fespace.
170  ConstantCoefficient one(1.0);
171  ParLinearForm * F = new ParLinearForm(test_space);
173  F->Assemble();
174 
175  ParGridFunction * x0 = new ParGridFunction(x0_space);
176  *x0 = 0.;
177 
178  // 8. Set up the mixed bilinear form for the primal trial unknowns, B0,
179  // the mixed bilinear form for the interfacial unknowns, Bhat,
180  // the inverse stiffness matrix on the discontinuous test space, Sinv,
181  // and the stiffness matrix on the continuous trial space, S0.
182  Array<int> ess_bdr(pmesh->bdr_attributes.Max());
183  ess_bdr = 1;
184  Array<int> ess_dof;
185  x0_space->GetEssentialVDofs(ess_bdr, ess_dof);
186 
187  ParMixedBilinearForm *B0 = new ParMixedBilinearForm(x0_space,test_space);
189  B0->Assemble();
190  B0->EliminateEssentialBCFromTrialDofs(ess_dof, *x0, *F);
191  B0->Finalize();
192 
193  ParMixedBilinearForm *Bhat = new ParMixedBilinearForm(xhat_space,test_space);
195  Bhat->Assemble();
196  Bhat->Finalize();
197 
198  ParBilinearForm *Sinv = new ParBilinearForm(test_space);
199  SumIntegrator *Sum = new SumIntegrator;
200  Sum->AddIntegrator(new DiffusionIntegrator(one));
201  Sum->AddIntegrator(new MassIntegrator(one));
202  Sinv->AddDomainIntegrator(new InverseIntegrator(Sum));
203  Sinv->Assemble();
204  Sinv->Finalize();
205 
206  ParBilinearForm *S0 = new ParBilinearForm(x0_space);
208  S0->Assemble();
209  S0->EliminateEssentialBC(ess_bdr);
210  S0->Finalize();
211 
212  HypreParMatrix * matB0 = B0->ParallelAssemble(); delete B0;
213  HypreParMatrix * matBhat = Bhat->ParallelAssemble(); delete Bhat;
214  HypreParMatrix * matSinv = Sinv->ParallelAssemble(); delete Sinv;
215  HypreParMatrix * matS0 = S0->ParallelAssemble(); delete S0;
216 
217  // 9. Define the block structure of the problem, by creating the offset
218  // variables. Also allocate two BlockVector objects to store the solution
219  // and rhs.
220  enum {x0_var, xhat_var, NVAR};
221 
222  int true_s0 = x0_space->TrueVSize();
223  int true_s1 = xhat_space->TrueVSize();
224  int true_s_test = test_space->TrueVSize();
225 
226  Array<int> true_offsets(NVAR+1);
227  true_offsets[0] = 0;
228  true_offsets[1] = true_s0;
229  true_offsets[2] = true_s0+true_s1;
230 
231  Array<int> true_offsets_test(2);
232  true_offsets_test[0] = 0;
233  true_offsets_test[1] = true_s_test;
234 
235  BlockVector x(true_offsets), b(true_offsets);
236  x = 0.0;
237  b = 0.0;
238 
239  // 10. Set up the 1x2 block Least Squares DPG operator, B = [B0 Bhat],
240  // the normal equation operator, A = B^t Sinv B, and
241  // the normal equation right-hand-size, b = B^t Sinv F.
242  BlockOperator B(true_offsets_test, true_offsets);
243  B.SetBlock(0, 0, matB0);
244  B.SetBlock(0, 1, matBhat);
245 
246  RAPOperator A(B, *matSinv, B);
247 
248  HypreParVector *trueF = F->ParallelAssemble();
249  {
250  HypreParVector SinvF(test_space);
251  matSinv->Mult(*trueF, SinvF);
252  B.MultTranspose(SinvF, b);
253  }
254 
255  // 11. Set up a block-diagonal preconditioner for the 2x2 normal equation
256  //
257  // [ S0^{-1} 0 ]
258  // [ 0 Shat^{-1} ] Shat = (Bhat^T Sinv Bhat)
259  //
260  // corresponding to the primal (x0) and interfacial (xhat) unknowns.
261  // Since the Shat operator is equivalent to an H(div) matrix reduced to
262  // the interfacial skeleton, we approximate its inverse with one V-cycle
263  // of the ADS preconditioner from the hypre library (in 2D we use AMS for
264  // the rotated H(curl) problem).
265  HypreBoomerAMG *S0inv = new HypreBoomerAMG(*matS0);
266  S0inv->SetPrintLevel(0);
267 
268  HypreParMatrix *Shat = RAP(matSinv, matBhat);
269  HypreSolver *Shatinv;
270  if (dim == 2) { Shatinv = new HypreAMS(*Shat, xhat_space); }
271  else { Shatinv = new HypreADS(*Shat, xhat_space); }
272 
273  BlockDiagonalPreconditioner P(true_offsets);
274  P.SetDiagonalBlock(0, S0inv);
275  P.SetDiagonalBlock(1, Shatinv);
276 
277  // 12. Solve the normal equation system using the PCG iterative solver.
278  // Check the weighted norm of residual for the DPG least square problem.
279  // Wrap the primal variable in a GridFunction for visualization purposes.
280  CGSolver pcg(MPI_COMM_WORLD);
281  pcg.SetOperator(A);
282  pcg.SetPreconditioner(P);
283  pcg.SetRelTol(1e-6);
284  pcg.SetMaxIter(200);
285  pcg.SetPrintLevel(1);
286  pcg.Mult(b, x);
287 
288  {
289  HypreParVector LSres(test_space), tmp(test_space);
290  B.Mult(x, LSres);
291  LSres -= *trueF;
292  matSinv->Mult(LSres, tmp);
293  double res = sqrt(InnerProduct(LSres, tmp));
294  if (myid == 0)
295  {
296  cout << "\n|| B0*x0 + Bhat*xhat - F ||_{S^-1} = " << res << endl;
297  }
298  }
299 
300  x0->Distribute(x.GetBlock(x0_var));
301 
302  // 13. Save the refined mesh and the solution in parallel. This output can
303  // be viewed later using GLVis: "glvis -np <np> -m mesh -g sol".
304  {
305  ostringstream mesh_name, sol_name;
306  mesh_name << "mesh." << setfill('0') << setw(6) << myid;
307  sol_name << "sol." << setfill('0') << setw(6) << myid;
308 
309  ofstream mesh_ofs(mesh_name.str().c_str());
310  mesh_ofs.precision(8);
311  pmesh->Print(mesh_ofs);
312 
313  ofstream sol_ofs(sol_name.str().c_str());
314  sol_ofs.precision(8);
315  x0->Save(sol_ofs);
316  }
317 
318  // 14. Send the solution by socket to a GLVis server.
319  if (visualization)
320  {
321  char vishost[] = "localhost";
322  int visport = 19916;
323  socketstream sol_sock(vishost, visport);
324  sol_sock << "parallel " << num_procs << " " << myid << "\n";
325  sol_sock.precision(8);
326  sol_sock << "solution\n" << *pmesh << *x0 << flush;
327  }
328 
329  // 15. Free the used memory.
330  delete trueF;
331  delete Shatinv;
332  delete S0inv;
333  delete Shat;
334  delete matB0;
335  delete matBhat;
336  delete matSinv;
337  delete matS0;
338  delete x0;
339  delete F;
340  delete test_space;
341  delete xhat_space;
342  delete x0_space;
343  delete test_fec;
344  delete xhat_fec;
345  delete x0_fec;
346  delete pmesh;
347 
348  MPI_Finalize();
349 
350  return 0;
351 }
Class for domain integration L(v) := (f, v)
Definition: lininteg.hpp:47
Conjugate gradient method.
Definition: solvers.hpp:109
The Auxiliary-space Maxwell Solver in hypre.
Definition: hypre.hpp:761
Integrator defining a sum of multiple Integrators.
Definition: bilininteg.hpp:150
The Auxiliary-space Divergence Solver in hypre.
Definition: hypre.hpp:793
Subclass constant coefficient.
Definition: coefficient.hpp:57
virtual void ReorientTetMesh()
See the remarks for the serial version in mesh.hpp.
Definition: pmesh.cpp:1432
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
virtual void Save(std::ostream &out) const
Definition: pgridfunc.cpp:306
Abstract parallel finite element space.
Definition: pfespace.hpp:28
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).
virtual void Mult(const Vector &x, Vector &y) const
Operator application.
Definition: solvers.cpp:294
virtual void Finalize(int skip_zeros=1)
Finalizes the matrix initialization.
HYPRE_Int Mult(HypreParVector &x, HypreParVector &y, double alpha=1.0, double beta=0.0)
Computes y = alpha * A * x + beta * y.
Definition: hypre.cpp:903
HypreParMatrix * ParallelAssemble()
Returns the matrix assembled on the true dofs, i.e. P^t A P.
HypreParMatrix * ParallelAssemble()
Returns the matrix assembled on the true dofs, i.e. P^t A P.
void EliminateEssentialBCFromTrialDofs(Array< int > &marked_vdofs, Vector &sol, Vector &rhs)
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
Integrator that inverts the matrix assembled by another integrator.
Definition: bilininteg.hpp:132
The BoomerAMG solver in hypre.
Definition: hypre.hpp:698
Class for parallel linear form.
Definition: plinearform.hpp:26
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. ...
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:7316
void SetPrintLevel(int print_level)
Definition: hypre.hpp:739
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)
void Assemble(int skip_zeros=1)
Assemble the local matrix.
The operator x -&gt; R*A*P*x.
Definition: operator.hpp:164
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.
Wrapper for hypre&#39;s parallel vector class.
Definition: hypre.hpp:58
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
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
void Distribute(const Vector *tv)
Definition: pgridfunc.cpp:79
double InnerProduct(HypreParVector *x, HypreParVector *y)
Definition: hypre.cpp:192
Class for parallel bilinear form.
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
Class for parallel bilinear form.
Abstract class for hypre&#39;s solvers and preconditioners.
Definition: hypre.hpp:522
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
Definition: solvers.hpp:123
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
Definition: solvers.cpp:92
Arbitrary order H1-conforming (continuous) finite elements.
Definition: fe_coll.hpp:54
virtual void GetEssentialVDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_dofs) const
Determine the boundary degrees of freedom.
Definition: pfespace.cpp:532
Class for parallel grid function.
Definition: pgridfunc.hpp:31
Wrapper for hypre&#39;s ParCSR matrix class.
Definition: hypre.hpp:143
A class to handle Block systems in a matrix-free implementation.
Class for parallel meshes.
Definition: pmesh.hpp:28
void EliminateEssentialBC(Array< int > &bdr_attr_is_ess, Vector &sol, Vector &rhs, int d=0)
void ParallelAssemble(Vector &tv)
Assemble the vector on the true dofs, i.e. P^t v.
Definition: plinearform.cpp:34
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
virtual void Print(std::ostream &out=std::cout) const
Definition: pmesh.cpp:2850
bool Good() const
Definition: optparser.hpp:120