MFEM  v3.2
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 = new Mesh(mesh_file, 1, 1);
80  int dim = mesh->Dimension();
81 
82  // 4. Refine the serial mesh on all processors to increase the resolution. In
83  // this example we do 'ref_levels' of uniform refinement. We choose
84  // 'ref_levels' to be the largest number that gives a final mesh with no
85  // more than 10,000 elements.
86  {
87  int ref_levels =
88  (int)floor(log(10000./mesh->GetNE())/log(2.)/dim);
89  for (int l = 0; l < ref_levels; l++)
90  {
91  mesh->UniformRefinement();
92  }
93  }
94 
95  // 5. Define a parallel mesh by a partitioning of the serial mesh. Refine
96  // this mesh further in parallel to increase the resolution. Once the
97  // parallel mesh is defined, the serial mesh can be deleted.
98  ParMesh *pmesh = new ParMesh(MPI_COMM_WORLD, *mesh);
99  delete mesh;
100  {
101  int par_ref_levels = 1;
102  for (int l = 0; l < par_ref_levels; l++)
103  {
104  pmesh->UniformRefinement();
105  }
106  }
107  pmesh->ReorientTetMesh();
108 
109  // 6. Define the trial, interfacial (trace) and test DPG spaces:
110  // - The trial space, x0_space, contains the non-interfacial unknowns and
111  // has the essential BC.
112  // - The interfacial space, xhat_space, contains the interfacial unknowns
113  // and does not have essential BC.
114  // - The test space, test_space, is an enriched space where the enrichment
115  // degree may depend on the spatial dimension of the domain, the type of
116  // the mesh and the trial space order.
117  int trial_order = order;
118  int trace_order = order - 1;
119  int test_order = order; // reduced order, full order is (order + dim - 1)
120  if (dim == 2 && (order%2 == 0 || (pmesh->MeshGenerator() & 2 && order > 1)))
121  {
122  test_order++;
123  }
124  if (test_order < trial_order)
125  if (myid == 0)
126  cerr << "Warning, test space not enriched enough to handle primal"
127  << " trial space\n";
128 
129  FiniteElementCollection *x0_fec, *xhat_fec, *test_fec;
130 
131  x0_fec = new H1_FECollection(trial_order, dim);
132  xhat_fec = new RT_Trace_FECollection(trace_order, dim);
133  test_fec = new L2_FECollection(test_order, dim);
134 
135  ParFiniteElementSpace *x0_space, *xhat_space, *test_space;
136 
137  x0_space = new ParFiniteElementSpace(pmesh, x0_fec);
138  xhat_space = new ParFiniteElementSpace(pmesh, xhat_fec);
139  test_space = new ParFiniteElementSpace(pmesh, test_fec);
140 
141  HYPRE_Int glob_true_s0 = x0_space->GlobalTrueVSize();
142  HYPRE_Int glob_true_s1 = xhat_space->GlobalTrueVSize();
143  HYPRE_Int glob_true_s_test = test_space->GlobalTrueVSize();
144  if (myid == 0)
145  {
146  cout << "\nNumber of Unknowns:\n"
147  << " Trial space, X0 : " << glob_true_s0
148  << " (order " << trial_order << ")\n"
149  << " Interface space, Xhat : " << glob_true_s1
150  << " (order " << trace_order << ")\n"
151  << " Test space, Y : " << glob_true_s_test
152  << " (order " << test_order << ")\n\n";
153  }
154 
155  // 7. Set up the linear form F(.) which corresponds to the right-hand side of
156  // the FEM linear system, which in this case is (f,phi_i) where f=1.0 and
157  // phi_i are the basis functions in the test finite element fespace.
158  ConstantCoefficient one(1.0);
159  ParLinearForm * F = new ParLinearForm(test_space);
161  F->Assemble();
162 
163  ParGridFunction * x0 = new ParGridFunction(x0_space);
164  *x0 = 0.;
165 
166  // 8. Set up the mixed bilinear form for the primal trial unknowns, B0,
167  // the mixed bilinear form for the interfacial unknowns, Bhat,
168  // the inverse stiffness matrix on the discontinuous test space, Sinv,
169  // and the stiffness matrix on the continuous trial space, S0.
170  Array<int> ess_bdr(pmesh->bdr_attributes.Max());
171  ess_bdr = 1;
172  Array<int> ess_dof;
173  x0_space->GetEssentialVDofs(ess_bdr, ess_dof);
174 
175  ParMixedBilinearForm *B0 = new ParMixedBilinearForm(x0_space,test_space);
177  B0->Assemble();
178  B0->EliminateEssentialBCFromTrialDofs(ess_dof, *x0, *F);
179  B0->Finalize();
180 
181  ParMixedBilinearForm *Bhat = new ParMixedBilinearForm(xhat_space,test_space);
183  Bhat->Assemble();
184  Bhat->Finalize();
185 
186  ParBilinearForm *Sinv = new ParBilinearForm(test_space);
187  SumIntegrator *Sum = new SumIntegrator;
188  Sum->AddIntegrator(new DiffusionIntegrator(one));
189  Sum->AddIntegrator(new MassIntegrator(one));
190  Sinv->AddDomainIntegrator(new InverseIntegrator(Sum));
191  Sinv->Assemble();
192  Sinv->Finalize();
193 
194  ParBilinearForm *S0 = new ParBilinearForm(x0_space);
196  S0->Assemble();
197  S0->EliminateEssentialBC(ess_bdr);
198  S0->Finalize();
199 
200  HypreParMatrix * matB0 = B0->ParallelAssemble(); delete B0;
201  HypreParMatrix * matBhat = Bhat->ParallelAssemble(); delete Bhat;
202  HypreParMatrix * matSinv = Sinv->ParallelAssemble(); delete Sinv;
203  HypreParMatrix * matS0 = S0->ParallelAssemble(); delete S0;
204 
205  // 9. Define the block structure of the problem, by creating the offset
206  // variables. Also allocate two BlockVector objects to store the solution
207  // and rhs.
208  enum {x0_var, xhat_var, NVAR};
209 
210  int true_s0 = x0_space->TrueVSize();
211  int true_s1 = xhat_space->TrueVSize();
212  int true_s_test = test_space->TrueVSize();
213 
214  Array<int> true_offsets(NVAR+1);
215  true_offsets[0] = 0;
216  true_offsets[1] = true_s0;
217  true_offsets[2] = true_s0+true_s1;
218 
219  Array<int> true_offsets_test(2);
220  true_offsets_test[0] = 0;
221  true_offsets_test[1] = true_s_test;
222 
223  BlockVector x(true_offsets), b(true_offsets);
224  x = 0.0;
225  b = 0.0;
226 
227  // 10. Set up the 1x2 block Least Squares DPG operator, B = [B0 Bhat],
228  // the normal equation operator, A = B^t Sinv B, and
229  // the normal equation right-hand-size, b = B^t Sinv F.
230  BlockOperator B(true_offsets_test, true_offsets);
231  B.SetBlock(0, 0, matB0);
232  B.SetBlock(0, 1, matBhat);
233 
234  RAPOperator A(B, *matSinv, B);
235 
236  HypreParVector *trueF = F->ParallelAssemble();
237  {
238  HypreParVector SinvF(test_space);
239  matSinv->Mult(*trueF, SinvF);
240  B.MultTranspose(SinvF, b);
241  }
242 
243  // 11. Set up a block-diagonal preconditioner for the 2x2 normal equation
244  //
245  // [ S0^{-1} 0 ]
246  // [ 0 Shat^{-1} ] Shat = (Bhat^T Sinv Bhat)
247  //
248  // corresponding to the primal (x0) and interfacial (xhat) unknowns.
249  // Since the Shat operator is equivalent to an H(div) matrix reduced to
250  // the interfacial skeleton, we approximate its inverse with one V-cycle
251  // of the ADS preconditioner from the hypre library (in 2D we use AMS for
252  // the rotated H(curl) problem).
253  HypreBoomerAMG *S0inv = new HypreBoomerAMG(*matS0);
254  S0inv->SetPrintLevel(0);
255 
256  HypreParMatrix *Shat = RAP(matSinv, matBhat);
257  HypreSolver *Shatinv;
258  if (dim == 2) { Shatinv = new HypreAMS(*Shat, xhat_space); }
259  else { Shatinv = new HypreADS(*Shat, xhat_space); }
260 
261  BlockDiagonalPreconditioner P(true_offsets);
262  P.SetDiagonalBlock(0, S0inv);
263  P.SetDiagonalBlock(1, Shatinv);
264 
265  // 12. Solve the normal equation system using the PCG iterative solver.
266  // Check the weighted norm of residual for the DPG least square problem.
267  // Wrap the primal variable in a GridFunction for visualization purposes.
268  CGSolver pcg(MPI_COMM_WORLD);
269  pcg.SetOperator(A);
270  pcg.SetPreconditioner(P);
271  pcg.SetRelTol(1e-6);
272  pcg.SetMaxIter(200);
273  pcg.SetPrintLevel(1);
274  pcg.Mult(b, x);
275 
276  {
277  HypreParVector LSres(test_space), tmp(test_space);
278  B.Mult(x, LSres);
279  LSres -= *trueF;
280  matSinv->Mult(LSres, tmp);
281  double res = sqrt(InnerProduct(LSres, tmp));
282  if (myid == 0)
283  {
284  cout << "\n|| B0*x0 + Bhat*xhat - F ||_{S^-1} = " << res << endl;
285  }
286  }
287 
288  x0->Distribute(x.GetBlock(x0_var));
289 
290  // 13. Save the refined mesh and the solution in parallel. This output can
291  // be viewed later using GLVis: "glvis -np <np> -m mesh -g sol".
292  {
293  ostringstream mesh_name, sol_name;
294  mesh_name << "mesh." << setfill('0') << setw(6) << myid;
295  sol_name << "sol." << setfill('0') << setw(6) << myid;
296 
297  ofstream mesh_ofs(mesh_name.str().c_str());
298  mesh_ofs.precision(8);
299  pmesh->Print(mesh_ofs);
300 
301  ofstream sol_ofs(sol_name.str().c_str());
302  sol_ofs.precision(8);
303  x0->Save(sol_ofs);
304  }
305 
306  // 14. Send the solution by socket to a GLVis server.
307  if (visualization)
308  {
309  char vishost[] = "localhost";
310  int visport = 19916;
311  socketstream sol_sock(vishost, visport);
312  sol_sock << "parallel " << num_procs << " " << myid << "\n";
313  sol_sock.precision(8);
314  sol_sock << "solution\n" << *pmesh << *x0 << flush;
315  }
316 
317  // 15. Free the used memory.
318  delete trueF;
319  delete Shatinv;
320  delete S0inv;
321  delete Shat;
322  delete matB0;
323  delete matBhat;
324  delete matSinv;
325  delete matS0;
326  delete x0;
327  delete F;
328  delete test_space;
329  delete xhat_space;
330  delete x0_space;
331  delete test_fec;
332  delete xhat_fec;
333  delete x0_fec;
334  delete pmesh;
335 
336  MPI_Finalize();
337 
338  return 0;
339 }
Class for domain integration L(v) := (f, v)
Definition: lininteg.hpp:47
Conjugate gradient method.
Definition: solvers.hpp:110
The Auxiliary-space Maxwell Solver in hypre.
Definition: hypre.hpp:768
Integrator defining a sum of multiple Integrators.
Definition: bilininteg.hpp:150
The Auxiliary-space Divergence Solver in hypre.
Definition: hypre.hpp:800
Subclass constant coefficient.
Definition: coefficient.hpp:57
virtual void ReorientTetMesh()
See the remarks for the serial version in mesh.hpp.
Definition: pmesh.cpp:1588
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:496
virtual void Save(std::ostream &out) const
Definition: pgridfunc.cpp:312
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.
int main(int argc, char *argv[])
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:940
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:705
Class for parallel linear form.
Definition: plinearform.hpp:26
int dim
Definition: ex3.cpp:47
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:6402
void SetPrintLevel(int print_level)
Definition: hypre.hpp:746
void SetMaxIter(int max_it)
Definition: solvers.hpp:62
HypreParMatrix * RAP(HypreParMatrix *A, HypreParMatrix *P)
Returns the matrix P^t * A * P.
Definition: hypre.cpp:1365
int MeshGenerator()
Definition: mesh.hpp:489
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:523
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:141
void SetRelTol(double rtol)
Definition: solvers.hpp:60
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:85
double InnerProduct(HypreParVector *x, HypreParVector *y)
Definition: hypre.cpp:196
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:529
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
Definition: solvers.hpp:124
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:77
virtual void GetEssentialVDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_dofs) const
Determine the boundary degrees of freedom.
Definition: pfespace.cpp:466
Class for parallel grid function.
Definition: pgridfunc.hpp:31
Wrapper for hypre&#39;s ParCSR matrix class.
Definition: hypre.hpp:150
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:130
virtual void Print(std::ostream &out=std::cout) const
Definition: pmesh.cpp:2927
bool Good() const
Definition: optparser.hpp:120