MFEM  v4.3.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
schwarz_ex1p.cpp
Go to the documentation of this file.
1 // Copyright (c) 2010-2021, 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 // -------------------------------------------------
13 // Overlapping Grids Miniapp: Poisson problem (ex1p)
14 // -------------------------------------------------
15 //
16 // This example code demonstrates use of MFEM to solve the Poisson problem:
17 //
18 // -Delta u = 1 \in [0, 1]^2, u_b = 0 \in \dO
19 //
20 // on two overlapping grids. Using simultaneous Schwarz iterations, the Poisson
21 // equation is solved iteratively, with boundary data interpolated between the
22 // overlapping boundaries for each grid. The overlapping Schwarz method was
23 // introduced by H. A. Schwarz in 1870, see also Section 2.2 of "Stability
24 // analysis of a singlerate and multirate predictor-corrector scheme for
25 // overlapping grids" by Mittal, Dutta and Fischer, arXiv:2010.00118.
26 //
27 // Compile with: make schwarz_ex1p
28 //
29 // Sample runs: mpirun -np 4 schwarz_ex1p -nm 3 -np1 2 -np2 1 -np3 1
30 // mpirun -np 4 schwarz_ex1p -nm 2 -np1 2 -np2 2
31 // mpirun -np 4 schwarz_ex1p -np1 2 -np2 2 -m1 ../../data/star.mesh -m2 ../../data/beam-quad.mesh
32 
33 #include "mfem.hpp"
34 #include <fstream>
35 #include <iostream>
36 
37 using namespace std;
38 using namespace mfem;
39 
40 // Method to use FindPointsGSLIB to determine the boundary points of a mesh that
41 // are interior to another mesh.
43  Vector &vxyz, int color,
44  Array<int> ess_tdof_list,
45  Array<int> &ess_tdof_list_int, int dim)
46 {
47  int number_boundary = ess_tdof_list.Size(),
48  number_true = vxyz.Size()/dim;
49 
50  Vector bnd(number_boundary*dim);
51  Array<unsigned int> colorv(number_boundary);
52  for (int i = 0; i < number_boundary; i++)
53  {
54  int idx = ess_tdof_list[i];
55  for (int d = 0; d < dim; d++)
56  {
57  bnd(i+d*number_boundary) = vxyz(idx + d*number_true);
58  }
59  colorv[i] = (unsigned int)color;
60  }
61 
62  finder.FindPoints(bnd, colorv);
63 
64  const Array<unsigned int> &code_out = finder.GetCode();
65 
66  // Setup ess_tdof_list_int
67  for (int i = 0; i < number_boundary; i++)
68  {
69  int idx = ess_tdof_list[i];
70  if (code_out[i] != 2) { ess_tdof_list_int.Append(idx); }
71  }
72 }
73 
74 int main(int argc, char *argv[])
75 {
76  // Initialize MPI.
77  int num_procs, myid;
78  MPI_Init(&argc, &argv);
79  MPI_Comm_size(MPI_COMM_WORLD, &num_procs);
80  MPI_Comm_rank(MPI_COMM_WORLD, &myid);
81 
82  // Parse command-line options.
83  int lim_meshes = 3; // should be greater than nmeshes
84  Array <const char *> mesh_file_list(lim_meshes);
85  Array <int> np_list(lim_meshes), rs_levels(lim_meshes),
86  rp_levels(lim_meshes);
87  mesh_file_list[0] = "../../data/square-disc.mesh";
88  mesh_file_list[1] = "../../data/inline-quad.mesh";
89  mesh_file_list[2] = "../../data/inline-quad.mesh";
90  int order = 2;
91  bool visualization = true;
92  rs_levels = 0;
93  rp_levels = 0;
94  np_list = 0;
95  double rel_tol = 1.e-8;
96  int nmeshes = 2;
97 
98  OptionsParser args(argc, argv);
99  args.AddOption(&mesh_file_list[0], "-m1", "--mesh",
100  "Mesh file to use.");
101  args.AddOption(&mesh_file_list[1], "-m2", "--mesh",
102  "Mesh file to use.");
103  args.AddOption(&mesh_file_list[2], "-m3", "--mesh",
104  "Mesh file to use.");
105  args.AddOption(&order, "-o", "--order",
106  "Finite element order (polynomial degree) or -1 for"
107  " isoparametric space.");
108  args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
109  "--no-visualization",
110  "Enable or disable GLVis visualization.");
111  args.AddOption(&rs_levels[0], "-r1", "--refine-serial",
112  "Number of times to refine the mesh 1 uniformly in serial.");
113  args.AddOption(&rs_levels[1], "-r2", "--refine-serial",
114  "Number of times to refine the mesh 2 uniformly in serial.");
115  args.AddOption(&rs_levels[2], "-r3", "--refine-serial",
116  "Number of times to refine the mesh 3 uniformly in serial.");
117  args.AddOption(&rp_levels[0], "-rp1", "--refine-parallel",
118  "Number of times to refine the mesh 1 uniformly in parallel.");
119  args.AddOption(&rp_levels[1], "-rp2", "--refine-parallel",
120  "Number of times to refine the mesh 2 uniformly in parallel.");
121  args.AddOption(&rp_levels[2], "-rp3", "--refine-parallel",
122  "Number of times to refine the mesh 3 uniformly in parallel.");
123  args.AddOption(&np_list[0], "-np1", "--np1",
124  "number of MPI ranks for mesh 1");
125  args.AddOption(&np_list[1], "-np2", "--np2",
126  "number of MPI ranks for mesh 2");
127  args.AddOption(&np_list[2], "-np3", "--np3",
128  "number of MPI ranks for mesh 3");
129  args.AddOption(&nmeshes, "-nm", "--nm",
130  "number of meshes");
131  args.AddOption(&rel_tol, "-rt", "--relative tolerance",
132  "Tolerance for Schwarz iteration convergence criterion.");
133 
134  args.Parse();
135  if (!args.Good())
136  {
137  args.PrintUsage(cout);
138  return 1;
139  }
140  if (myid == 0)
141  {
142  args.PrintOptions(cout);
143  }
144 
145  // Setup MPI communicator for each mesh by splitting MPI_COMM_WORLD.
146  MPI_Comm *comml = new MPI_Comm;
147  int color = 0;
148  int npsum = 0;
149  for (int i = 0; i < nmeshes; i++)
150  {
151  npsum += np_list[i];
152  if (myid < npsum) { color = i; break; }
153  }
154 
155  MPI_Comm_split(MPI_COMM_WORLD, color, myid, comml);
156  int myidlocal, numproclocal;
157  MPI_Comm_rank(*comml, &myidlocal);
158  MPI_Comm_size(*comml, &numproclocal);
159 
160  // Check number of mpi ranks specified for each mesh.
161  MFEM_VERIFY(np_list.Sum() == num_procs, " The individual mpi ranks for each"
162  " of the meshes do not add up to"
163  " the total ranks specified.");
164 
165  // Read the mesh from the given mesh file. We can handle triangular,
166  // quadrilateral, tetrahedral, hexahedral, surface and volume meshes with
167  // the same code.
168  Mesh *mesh = new Mesh(mesh_file_list[color], 1, 1);
169  int dim = mesh->Dimension();
170 
171  // Refine the mesh to increase the resolution. In this example we do
172  // 'ref_levels' of uniform refinement. We choose 'ref_levels' to be the
173  // largest number that gives a final mesh with no more than 50,000 elements.
174  for (int lev = 0; lev < rs_levels[color]; lev++) { mesh->UniformRefinement(); }
175 
176  // Define a parallel mesh by a partitioning of the serial mesh. Refine this
177  // mesh further in parallel to increase the resolution. Once the parallel
178  // mesh is defined, the serial mesh can be deleted.
179  ParMesh *pmesh;
180  pmesh = new ParMesh(*comml, *mesh);
181  for (int l = 0; l < rp_levels[color]; l++)
182  {
183  pmesh->UniformRefinement();
184  }
185  delete mesh;
186 
187  // Define a finite element space on the mesh. Here we use continuous
188  // Lagrange finite elements of the specified order. If order < 1, we
189  // instead use an isoparametric/isogeometric space.
191  if (order > 0)
192  {
193  fec = new H1_FECollection(order, dim);
194  }
195  else if (pmesh->GetNodes())
196  {
197  fec = pmesh->GetNodes()->OwnFEC();
198  if (myid == 0)
199  {
200  cout << "Using isoparametric FEs: " << fec->Name() << endl;
201  }
202  }
203  else
204  {
205  fec = new H1_FECollection(order = 1, dim);
206  }
207  ParFiniteElementSpace *fespace = new ParFiniteElementSpace(pmesh, fec);
208  HYPRE_BigInt size = fespace->GlobalTrueVSize();
209  if (myid == 0)
210  {
211  cout << "Number of finite element unknowns: " << size << endl;
212  }
213 
214  // Determine the list of true (i.e. conforming) essential boundary dofs.
215  // In this example, the boundary conditions are defined by marking all
216  // the boundary attributes from the mesh as essential (Dirichlet) and
217  // converting them to a list of true dofs.
218  Array<int> ess_tdof_list;
219  if (pmesh->bdr_attributes.Size())
220  {
221  Array<int> ess_bdr(pmesh->bdr_attributes.Max());
222  ess_bdr = 1;
223  fespace->GetEssentialTrueDofs(ess_bdr, ess_tdof_list);
224  }
225 
226  // Set up the linear form b(.) which corresponds to the right-hand side of
227  // the FEM linear system, which in this case is (1,phi_i) where phi_i are
228  // the basis functions in the finite element fespace1.
229  ParLinearForm *b = new ParLinearForm(fespace);
230  ConstantCoefficient one(1.0);
232  b->Assemble();
233 
234  // Define the solution vector x as a finite element grid function
235  // corresponding to fespace1. Initialize x with initial guess of zero,
236  // which satisfies the boundary conditions.
237  ParGridFunction x(fespace);
238  x = 0.0;
239  x.SetTrueVector();
240 
241  // Setup FindPointsGSLIB and determine points on each mesh's boundary that
242  // are interior to another mesh.
243  pmesh->SetCurvature(order, false, dim, Ordering::byNODES);
244  {
245  Vector vxyz = *pmesh->GetNodes();
246 
247  // For the default mesh inputs, we need to rescale inline-quad.mesh such
248  // that it does not cover the entire domain [0, 1]^2 and still has a
249  // non-trivial overlap with the other mesh.
250  if (strcmp(mesh_file_list[0], "../../data/square-disc.mesh") == 0 &&
251  strcmp(mesh_file_list[1], "../../data/inline-quad.mesh") == 0 &&
252  strcmp(mesh_file_list[2], "../../data/inline-quad.mesh") == 0 )
253  {
254  if (nmeshes == 2)
255  {
256  if (color == 1) // rescale from [0, 1]^2 to [0.25, 0.75]^2
257  {
258  for (int i = 0; i < vxyz.Size(); i++)
259  {
260  vxyz(i) = 0.5 + 0.5*(vxyz(i)-0.5);
261  }
262  }
263  }
264  else if (nmeshes == 3)
265  {
266  if (color == 1)
267  {
268  // rescale from [0, 1]^2 to [0.21, 0.61] in x and [0.25, 0.75] in y
269  const int pts_cnt = vxyz.Size()/dim;
270  for (int i = 0; i < pts_cnt; i++)
271  {
272  vxyz(i) = 0.41 + 0.4*(vxyz(i)-0.5);
273  }
274  for (int i = 0; i < pts_cnt; i++)
275  {
276  vxyz(i+pts_cnt) = 0.5 + 0.5*(vxyz(i+pts_cnt)-0.5);
277  }
278  }
279  else if (color == 2)
280  {
281  // rescale from [0, 1]^2 to [0.4, 0.8] in x and [0.2, 0.8] in y
282  const int pts_cnt = vxyz.Size()/dim;
283  for (int i = 0; i < pts_cnt; i++)
284  {
285  vxyz(i) = 0.6 + 0.4*(vxyz(i)-0.5);
286  }
287  for (int i = 0; i < pts_cnt; i++)
288  {
289  vxyz(i+pts_cnt) = 0.5 + 0.6*(vxyz(i+pts_cnt)-0.5);
290  }
291  }
292  }
293  }
294  pmesh->SetNodes(vxyz);
295  }
296 
297  pmesh->GetNodes()->SetTrueVector();
298  Vector vxyz = pmesh->GetNodes()->GetTrueVector();
299 
300 
301  OversetFindPointsGSLIB finder(MPI_COMM_WORLD);
302  finder.Setup(*pmesh, color);
303 
304  Array<int> ess_tdof_list_int;
305  GetInterdomainBoundaryPoints(finder, vxyz, color,
306  ess_tdof_list, ess_tdof_list_int, dim);
307 
308  // Use FindPointsGSLIB to interpolate the solution at interdomain boundary
309  // points.
310  const int number_boundary = ess_tdof_list_int.Size(),
311  number_true = vxyz.Size()/dim;
312 
313  int number_boundary_g = number_boundary;
314  MPI_Allreduce(&number_boundary, &number_boundary_g, 1, MPI_INT, MPI_SUM,
315  *comml);
316  MFEM_VERIFY(number_boundary_g != 0, " Please use overlapping grids.");
317 
318  Array<unsigned int> colorv;
319  colorv.SetSize(number_boundary);
320 
321  MPI_Barrier(MPI_COMM_WORLD);
322  Vector bnd(number_boundary*dim);
323  for (int i = 0; i < number_boundary; i++)
324  {
325  int idx = ess_tdof_list_int[i];
326  for (int d = 0; d < dim; d++)
327  {
328  bnd(i+d*number_boundary) = vxyz(idx + d*number_true);
329  }
330  colorv[i] = (unsigned int)color;
331  }
332  Vector interp_vals1(number_boundary);
333  finder.Interpolate(bnd, colorv, x, interp_vals1);
334 
335  // Set up the bilinear form a(.,.) on the finite element space corresponding
336  // to the Laplacian operator -Delta, by adding a Diffusion integrator.
337  ParBilinearForm *a = new ParBilinearForm(fespace);
339 
340  // Assemble the bilinear form and the corresponding linear system,
341  // applying any necessary transformations such as: eliminating boundary
342  // conditions, applying conforming constraints for non-conforming AMR,
343  // static condensation, etc.
344  a->Assemble();
345 
346  delete b;
347 
348  // Use simultaneous Schwarz iterations to iteratively solve the PDE and
349  // interpolate interdomain boundary data to impose Dirichlet boundary
350  // conditions.
351 
352  int NiterSchwarz = 100;
353  for (int schwarz = 0; schwarz < NiterSchwarz; schwarz++)
354  {
355  ParLinearForm *b = new ParLinearForm(fespace);
357  b->Assemble();
358 
359  OperatorPtr A;
360  Vector B, X;
361  a->FormLinearSystem(ess_tdof_list, x, *b, A, X, B);
362 
363  // Solve the linear system A X = B.
364  // Use a simple symmetric Gauss-Seidel preconditioner with PCG.
365  Solver *prec = NULL;
366  prec = new HypreBoomerAMG;
367  dynamic_cast<HypreBoomerAMG *>(prec)->SetPrintLevel(-1);
368  CGSolver cg(*comml);
369  cg.SetRelTol(1e-12);
370  cg.SetMaxIter(2000);
371  cg.SetPrintLevel(0);
372  cg.SetPreconditioner(*prec);
373  cg.SetOperator(*A);
374  cg.Mult(B, X);
375  delete prec;
376 
377  // Recover the solution as a finite element grid function.
378  a->RecoverFEMSolution(X, *b, x);
379 
380  // Interpolate boundary condition
381  finder.Interpolate(x, interp_vals1);
382 
383  double dxmax = std::numeric_limits<float>::min();
384  double xinf = x.Normlinf();
385  double xinfg = xinf;
386  MPI_Allreduce(&xinf, &xinfg, 1, MPI_DOUBLE, MPI_MAX, MPI_COMM_WORLD);
387  x.SetTrueVector();
388  Vector xt = x.GetTrueVector();
389  for (int i = 0; i < number_boundary; i++)
390  {
391  int idx = ess_tdof_list_int[i];
392  double dx = std::abs(xt(idx)-interp_vals1(i))/xinfg;
393  if (dx > dxmax) { dxmax = dx; }
394  xt(idx) = interp_vals1(i);
395  }
396  x.SetFromTrueDofs(xt);
397  double dxmaxg = dxmax;
398  MPI_Allreduce(&dxmax, &dxmaxg, 1, MPI_DOUBLE, MPI_MAX, MPI_COMM_WORLD);
399 
400  delete b;
401 
402  if (myid == 0)
403  {
404  std::cout << std::setprecision(8) <<
405  "Iteration: " << schwarz <<
406  ", Relative residual: " << dxmaxg << endl;
407  }
408 
409  if (dxmaxg < rel_tol) { break; }
410  }
411 
412  // Send the solution by socket to a GLVis server.
413  if (visualization)
414  {
415  char vishost[] = "localhost";
416  int visport = 19916;
417  socketstream sol_sock(vishost, visport);
418  sol_sock << "parallel " << num_procs << " " << myid << "\n";
419  sol_sock.precision(8);
420  sol_sock << "solution\n" << *pmesh << x << flush;
421  }
422 
423  // 15. Free the used memory.
424  finder.FreeData();
425  delete a;
426  delete fespace;
427  if (order > 0) { delete fec; }
428  delete pmesh;
429  delete comml;
430 
431 
432  MPI_Finalize();
433 
434  return 0;
435 }
int Size() const
Return the logical size of the array.
Definition: array.hpp:134
Class for domain integration L(v) := (f, v)
Definition: lininteg.hpp:97
virtual void GetEssentialTrueDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list, int component=-1)
Definition: pfespace.cpp:790
Conjugate gradient method.
Definition: solvers.hpp:316
A coefficient that is constant across space and time.
Definition: coefficient.hpp:78
virtual void FormLinearSystem(const Array< int > &ess_tdof_list, Vector &x, Vector &b, OperatorHandle &A, Vector &X, Vector &B, int copy_interior=0)
Form the linear system A X = B, corresponding to this bilinear form and the linear form b(...
virtual void SetCurvature(int order, bool discont=false, int space_dim=-1, int ordering=1)
Definition: pmesh.cpp:1958
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
Definition: solvers.cpp:616
Pointer to an Operator of a specified type.
Definition: handle.hpp:33
HYPRE_BigInt GlobalTrueVSize() const
Definition: pfespace.hpp:275
int Size() const
Returns the size of the vector.
Definition: vector.hpp:190
Abstract parallel finite element space.
Definition: pfespace.hpp:28
double Normlinf() const
Returns the l_infinity norm of the vector.
Definition: vector.cpp:801
The BoomerAMG solver in hypre.
Definition: hypre.hpp:1387
Class for parallel linear form.
Definition: plinearform.hpp:26
void SetTrueVector()
Shortcut for calling GetTrueDofs() with GetTrueVector() as argument.
Definition: gridfunc.hpp:137
void SetPrintLevel(int print_lvl)
Definition: solvers.cpp:71
virtual void SetFromTrueDofs(const Vector &tv)
Set the GridFunction from the given true-dof vector.
Definition: pgridfunc.hpp:164
void Parse()
Parse the command-line options. Note that this function expects all the options provided through the ...
Definition: optparser.cpp:150
OversetFindPointsGSLIB enables use of findpts for arbitrary number of overlapping grids...
Definition: gslib.hpp:192
int Append(const T &el)
Append element &#39;el&#39; to array, resize if necessary.
Definition: array.hpp:746
constexpr char vishost[]
virtual void FreeData()
Definition: gslib.cpp:212
double b
Definition: lissajous.cpp:42
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:9143
constexpr int visport
void SetMaxIter(int max_it)
Definition: solvers.hpp:100
T Max() const
Find the maximal element in the array, using the comparison operator &lt; for class T.
Definition: array.cpp:68
void Assemble(int skip_zeros=1)
Assemble the local matrix.
void PrintUsage(std::ostream &out) const
Print the usage message.
Definition: optparser.cpp:457
const Vector & GetTrueVector() const
Read only access to the (optional) internal true-dof Vector.
Definition: gridfunc.hpp:127
void FindPoints(const Vector &point_pos, Array< unsigned int > &point_id)
Definition: gslib.cpp:919
void AddDomainIntegrator(LinearFormIntegrator *lfi)
Adds new Domain Integrator. Assumes ownership of lfi.
Definition: linearform.cpp:39
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
Definition: mesh.hpp:204
void SetRelTol(double rtol)
Definition: solvers.hpp:98
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:674
HYPRE_Int HYPRE_BigInt
virtual void RecoverFEMSolution(const Vector &X, const Vector &b, Vector &x)
double a
Definition: lissajous.cpp:41
virtual const char * Name() const
Definition: fe_coll.hpp:61
int dim
Definition: ex24.cpp:53
virtual const Array< unsigned int > & GetCode() const
Definition: gslib.hpp:169
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
Adds new Domain Integrator. Assumes ownership of bfi.
void PrintOptions(std::ostream &out) const
Print the options.
Definition: optparser.cpp:327
void Interpolate(const Vector &point_pos, Array< unsigned int > &point_id, const GridFunction &field_in, Vector &field_out)
Definition: gslib.cpp:986
Class for parallel bilinear form.
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
Definition: solvers.hpp:330
Vector data type.
Definition: vector.hpp:60
void GetNodes(Vector &node_coord) const
Definition: mesh.cpp:7343
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:216
void GetInterdomainBoundaryPoints(FindPointsGSLIB &finder1, FindPointsGSLIB &finder2, Vector &mesh_nodes_1, Vector &mesh_nodes_2, Array< int > ess_tdof_list1, Array< int > ess_tdof_list2, Array< int > &ess_tdof_list1_int, Array< int > &ess_tdof_list2_int, const int dim)
void SetNodes(const Vector &node_coord)
Definition: mesh.cpp:7355
Base class for solvers.
Definition: operator.hpp:648
Class for parallel grid function.
Definition: pgridfunc.hpp:32
Class for parallel meshes.
Definition: pmesh.hpp:32
void Setup(Mesh &m, const int meshid, GridFunction *gfmax=NULL, const double bb_t=0.1, const double newt_tol=1.0e-12, const int npt_max=256)
Definition: gslib.cpp:847
int main()
bool Good() const
Return true if the command line options were parsed successfully.
Definition: optparser.hpp:150