MFEM v4.7.0
Finite element discretization library
Loading...
Searching...
No Matches
schwarz_ex1.cpp
Go to the documentation of this file.
1// Copyright (c) 2010-2024, 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_ex1
28//
29// Sample runs: schwarz_ex1
30// schwarz_ex1 -m1 ../../data/star.mesh -m2 ../../data/beam-quad.mesh
31
32#include "mfem.hpp"
33#include <fstream>
34#include <iostream>
35
36using namespace std;
37using namespace mfem;
38
39// Method to use FindPointsGSLIB to determine the boundary points of a mesh that
40// are interior to another mesh.
42 FindPointsGSLIB &finder2,
43 Vector &mesh_nodes_1, Vector &mesh_nodes_2,
44 Array<int> ess_tdof_list1,
45 Array<int> ess_tdof_list2,
46 Array<int> &ess_tdof_list1_int,
47 Array<int> &ess_tdof_list2_int,
48 const int dim);
49
50int main(int argc, char *argv[])
51{
52 // Parse command-line options.
53 const char *mesh_file_1 = "../../data/square-disc.mesh";
54 const char *mesh_file_2 = "../../data/inline-quad.mesh";
55 int order = 2;
56 bool visualization = true;
57 int r1_levels = 0;
58 int r2_levels = 0;
59 double rel_tol = 1.e-8;
60
61 OptionsParser args(argc, argv);
62 args.AddOption(&mesh_file_1, "-m1", "--mesh",
63 "Mesh file to use.");
64 args.AddOption(&mesh_file_2, "-m2", "--mesh",
65 "Mesh file to use.");
66 args.AddOption(&order, "-o", "--order",
67 "Finite element order (polynomial degree) or -1 for"
68 " isoparametric space.");
69 args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
70 "--no-visualization",
71 "Enable or disable GLVis visualization.");
72 args.AddOption(&r1_levels, "-r1", "--refine-serial",
73 "Number of times to refine the mesh uniformly in serial.");
74 args.AddOption(&r2_levels, "-r2", "--refine-serial",
75 "Number of times to refine the mesh uniformly in serial.");
76 args.AddOption(&rel_tol, "-rt", "--relative tolerance",
77 "Tolerance for Schwarz iteration convergence criterion.");
78
79 args.Parse();
80 if (!args.Good())
81 {
82 args.PrintUsage(cout);
83 return 1;
84 }
85 args.PrintOptions(cout);
86
87 const int nmeshes = 2;
88
89 // Read the mesh from the given mesh file. We can handle triangular,
90 // quadrilateral, tetrahedral, hexahedral, surface and volume meshes with
91 // the same code.
92 Array <Mesh*> mesharr(2);
93 mesharr[0] = new Mesh(mesh_file_1, 1, 1);
94 mesharr[1] = new Mesh(mesh_file_2, 1, 1);
95 int dim = mesharr[0]->Dimension();
96
97
98 // Refine the mesh to increase the resolution. In this example we do
99 // 'ref_levels' of uniform refinement. We choose 'ref_levels' to be the
100 // largest number that gives a final mesh with no more than 50,000 elements.
101 for (int lev = 0; lev < r1_levels; lev++) { mesharr[0]->UniformRefinement(); }
102 for (int lev = 0; lev < r2_levels; lev++) { mesharr[1]->UniformRefinement(); }
103
104 // Define a finite element space on the mesh. Here we use continuous
105 // Lagrange finite elements of the specified order. If order < 1, we
106 // instead use an isoparametric/isogeometric space.
107 Array <FiniteElementCollection*> fecarr(nmeshes);
108 Array <FiniteElementSpace*> fespacearr(nmeshes);
109 for (int i = 0; i< nmeshes; i ++)
110 {
111 if (order > 0)
112 {
113 fecarr[i] = new H1_FECollection(order, dim);
114 }
115 else if (mesharr[i]->GetNodes())
116 {
117 fecarr[i] = mesharr[i]->GetNodes()->OwnFEC();
118 cout << "Using isoparametric FEs: " << fecarr[0]->Name() << endl;
119 }
120 else
121 {
122 fecarr[i] = new H1_FECollection(order = 1, dim);
123 }
124 fespacearr[i] = new FiniteElementSpace(mesharr[i], fecarr[i]);
125 }
126
127 // Determine the list of true (i.e. conforming) essential boundary dofs.
128 // In this example, the boundary conditions are defined by marking all
129 // the boundary attributes from the mesh as essential (Dirichlet) and
130 // converting them to a list of true dofs.
131 Array<int> ess_tdof_list1, ess_tdof_list2;
132 if (mesharr[0]->bdr_attributes.Size())
133 {
134 Array<int> ess_bdr(mesharr[0]->bdr_attributes.Max());
135 ess_bdr = 1;
136 fespacearr[0]->GetEssentialTrueDofs(ess_bdr, ess_tdof_list1);
137 }
138
139 if (mesharr[1]->bdr_attributes.Size())
140 {
141 Array<int> ess_bdr(mesharr[1]->bdr_attributes.Max());
142 ess_bdr = 1;
143 fespacearr[1]->GetEssentialTrueDofs(ess_bdr, ess_tdof_list2);
144 }
145
146 // Set up the linear form b(.) which corresponds to the right-hand side of
147 // the FEM linear system, which in this case is (1,phi_i) where phi_i are
148 // the basis functions in the finite element fespace1.
149 ConstantCoefficient one(1.0);
150 Array<LinearForm*> b_ar(nmeshes);
151 for (int i = 0; i< nmeshes; i ++)
152 {
153 b_ar[i] = new LinearForm(fespacearr[i]);
154 b_ar[i]->AddDomainIntegrator(new DomainLFIntegrator(one));
155 b_ar[i]->Assemble();
156 }
157
158 // Define the solution vector x as a finite element grid function
159 // corresponding to fespace1. Initialize x with initial guess of zero,
160 // which satisfies the boundary conditions.
161 GridFunction x1(fespacearr[0]), x2(fespacearr[1]);
162 x1 = 0;
163 x2 = 0;
164
165 // Setup FindPointsGSLIB and determine points on each mesh's boundary that
166 // are interior to another mesh.
167 mesharr[0]->SetCurvature(order, false, dim, Ordering::byNODES);
168 Vector mesh_nodes_1 = *mesharr[0]->GetNodes();
169 mesharr[1]->SetCurvature(order, false, dim, Ordering::byNODES);
170 Vector mesh_nodes_2 = *mesharr[1]->GetNodes();
171
172 // For the default mesh inputs, we need to rescale inline-quad.mesh
173 // such that it does not cover the entire domain [0, 1]^2 and still has a
174 // non-trivial overlap with the other mesh.
175 if (strcmp(mesh_file_1, "../../data/square-disc.mesh") == 0 &&
176 strcmp(mesh_file_2, "../../data/inline-quad.mesh") == 0 )
177 {
178 for (int i = 0; i < mesh_nodes_2.Size(); i++)
179 {
180 mesh_nodes_2(i) = 0.5 + 0.5*(mesh_nodes_2(i)-0.5);
181 }
182 mesharr[1]->SetNodes(mesh_nodes_2);
183 }
184
185 FindPointsGSLIB finder1, finder2;
186 finder1.Setup(*mesharr[0]);
187 finder2.Setup(*mesharr[1]);
188
189 Array<int> ess_tdof_list1_int, ess_tdof_list2_int;
190 GetInterdomainBoundaryPoints(finder1, finder2, mesh_nodes_1, mesh_nodes_2,
191 ess_tdof_list1, ess_tdof_list2,
192 ess_tdof_list1_int, ess_tdof_list2_int, dim);
193
194 // Use FindPointsGSLIB to interpolate the solution at interdomain boundary
195 // points.
196 const int number_boundary_1 = ess_tdof_list1_int.Size(),
197 number_boundary_2 = ess_tdof_list2_int.Size(),
198 number_true_1 = mesh_nodes_1.Size()/dim,
199 number_true_2 = mesh_nodes_2.Size()/dim;
200
201 MFEM_VERIFY(number_boundary_1 != 0 ||
202 number_boundary_2 != 0, " Please use overlapping grids.");
203
204 Vector bnd1(number_boundary_1*dim);
205 for (int i = 0; i < number_boundary_1; i++)
206 {
207 int idx = ess_tdof_list1_int[i];
208 for (int d = 0; d < dim; d++)
209 {
210 bnd1(i+d*number_boundary_1) = mesh_nodes_1(idx + d*number_true_1);
211 }
212 }
213
214 Vector bnd2(number_boundary_2*dim);
215 for (int i = 0; i < number_boundary_2; i++)
216 {
217 int idx = ess_tdof_list2_int[i];
218 for (int d = 0; d < dim; d++)
219 {
220 bnd2(i+d*number_boundary_2) = mesh_nodes_2(idx + d*number_true_2);
221 }
222 }
223
224 Vector interp_vals1(number_boundary_1), interp_vals2(number_boundary_2);
225 finder1.Interpolate(bnd2, x1, interp_vals2);
226 finder2.Interpolate(bnd1, x2, interp_vals1);
227
228 // Set up the bilinear form a(.,.) on the finite element space corresponding
229 // to the Laplacian operator -Delta, by adding a Diffusion integrator.
230 Array <BilinearForm*> a_ar(2);
231 a_ar[0] = new BilinearForm(fespacearr[0]);
232 a_ar[1] = new BilinearForm(fespacearr[1]);
233 a_ar[0]->AddDomainIntegrator(new DiffusionIntegrator(one));
234 a_ar[1]->AddDomainIntegrator(new DiffusionIntegrator(one));
235
236 // Assemble the bilinear form and the corresponding linear system,
237 // applying any necessary transformations such as: eliminating boundary
238 // conditions, applying conforming constraints for non-conforming AMR,
239 // static condensation, etc.
240 a_ar[0]->Assemble();
241 a_ar[1]->Assemble();
242
243 delete b_ar[0];
244 delete b_ar[1];
245
246 // Use simultaneous Schwarz iterations to iteratively solve the PDE and
247 // interpolate interdomain boundary data to impose Dirichlet boundary
248 // conditions.
249
250 int NiterSchwarz = 100;
251 for (int schwarz = 0; schwarz < NiterSchwarz; schwarz++)
252 {
253 for (int i = 0; i < nmeshes; i ++)
254 {
255 b_ar[i] = new LinearForm(fespacearr[i]);
256 b_ar[i]->AddDomainIntegrator(new DomainLFIntegrator(one));
257 b_ar[i]->Assemble();
258 }
259
260 OperatorPtr A1, A2;
261 Vector B1, X1, B2, X2;
262
263 a_ar[0]->FormLinearSystem(ess_tdof_list1, x1, *b_ar[0], A1, X1, B1);
264 a_ar[1]->FormLinearSystem(ess_tdof_list2, x2, *b_ar[1], A2, X2, B2);
265
266 // Solve the linear system A X = B.
267 // Use a simple symmetric Gauss-Seidel preconditioner with PCG.
268 GSSmoother M1((SparseMatrix&)(*A1));
269 PCG(*A1, M1, B1, X1, 0, 200, 1e-12, 0.0);
270 GSSmoother M2((SparseMatrix&)(*A2));
271 PCG(*A2, M2, B2, X2, 0, 200, 1e-12, 0.0);
272
273 // Recover the solution as a finite element grid function.
274 a_ar[0]->RecoverFEMSolution(X1, *b_ar[0], x1);
275 a_ar[1]->RecoverFEMSolution(X2, *b_ar[1], x2);
276
277 // Interpolate boundary condition
278 finder1.Interpolate(x1, interp_vals2);
279 finder2.Interpolate(x2, interp_vals1);
280
281 double dxmax = std::numeric_limits<float>::min();
282 double x1inf = x1.Normlinf();
283 double x2inf = x2.Normlinf();
284 for (int i = 0; i < number_boundary_1; i++)
285 {
286 int idx = ess_tdof_list1_int[i];
287 double dx = std::abs(x1(idx)-interp_vals1(i))/x1inf;
288 if (dx > dxmax) { dxmax = dx; }
289 x1(idx) = interp_vals1(i);
290 }
291
292 for (int i = 0; i < number_boundary_2; i++)
293 {
294 int idx = ess_tdof_list2_int[i];
295 double dx = std::abs(x2(idx)-interp_vals2(i))/x2inf;
296 if (dx > dxmax) { dxmax = dx; }
297 x2(idx) = interp_vals2(i);
298 }
299
300 delete b_ar[0];
301 delete b_ar[1];
302
303 std::cout << std::setprecision(8) <<
304 "Iteration: " << schwarz <<
305 ", Relative residual: " << dxmax << endl;
306 if (dxmax < rel_tol) { break; }
307 }
308
309 // Send the solution by socket to a GLVis server.
310 if (visualization)
311 {
312 char vishost[] = "localhost";
313 int visport = 19916;
314 for (int ip = 0; ip<mesharr.Size(); ++ip)
315 {
316 socketstream sol_sock(vishost, visport);
317 sol_sock.precision(8);
318 sol_sock << "parallel " << mesharr.Size() << " " << ip << "\n";
319 if (ip==0) { sol_sock << "solution\n" << *mesharr[ip] << x1 << flush; }
320 if (ip==1) { sol_sock << "solution\n" << *mesharr[ip] << x2 << flush; }
321 sol_sock << "window_title 'Overlapping grid solution'\n"
322 << "window_geometry "
323 << 0 << " " << 0 << " " << 450 << " " << 350 << "\n"
324 << "keys jmcA]]]" << endl;
325 }
326 }
327
328 // Free the used memory.
329 finder1.FreeData();
330 finder2.FreeData();
331 for (int i = 0; i < nmeshes; i++)
332 {
333 delete a_ar[i];
334 delete fespacearr[i];
335 if (order > 0) { delete fecarr[i]; }
336 delete mesharr[i];
337 }
338
339 return 0;
340}
341
343 FindPointsGSLIB &finder2,
344 Vector &mesh_nodes_1, Vector &mesh_nodes_2,
345 Array<int> ess_tdof_list1,
346 Array<int> ess_tdof_list2,
347 Array<int> &ess_tdof_list1_int,
348 Array<int> &ess_tdof_list2_int,
349 const int dim)
350{
351 int number_boundary_1 = ess_tdof_list1.Size(),
352 number_boundary_2 = ess_tdof_list2.Size(),
353 number_true_1 = mesh_nodes_1.Size()/dim,
354 number_true_2 = mesh_nodes_2.Size()/dim;
355
356 Vector bnd1(number_boundary_1*dim);
357 for (int i = 0; i < number_boundary_1; i++)
358 {
359 int idx = ess_tdof_list1[i];
360 for (int d = 0; d < dim; d++)
361 {
362 bnd1(i+d*number_boundary_1) = mesh_nodes_1(idx + d*number_true_1);
363 }
364 }
365
366 Vector bnd2(number_boundary_2*dim);
367 for (int i = 0; i < number_boundary_2; i++)
368 {
369 int idx = ess_tdof_list2[i];
370 for (int d = 0; d < dim; d++)
371 {
372 bnd2(i+d*number_boundary_2) = mesh_nodes_2(idx + d*number_true_2);
373 }
374 }
375
376 finder1.FindPoints(bnd2);
377 finder2.FindPoints(bnd1);
378
379 const Array<unsigned int> &code_out1 = finder1.GetCode();
380 const Array<unsigned int> &code_out2 = finder2.GetCode();
381
382 // Setup ess_tdof_list_int
383 for (int i = 0; i < number_boundary_1; i++)
384 {
385 int idx = ess_tdof_list1[i];
386 if (code_out2[i] != 2) { ess_tdof_list1_int.Append(idx); }
387 }
388
389 for (int i = 0; i < number_boundary_2; i++)
390 {
391 int idx = ess_tdof_list2[i];
392 if (code_out1[i] != 2) { ess_tdof_list2_int.Append(idx); }
393 }
394}
T Max() const
Find the maximal element in the array, using the comparison operator < for class T.
Definition array.cpp:68
int Size() const
Return the logical size of the array.
Definition array.hpp:144
int Append(const T &el)
Append element 'el' to array, resize if necessary.
Definition array.hpp:769
A "square matrix" operator for the associated FE space and BLFIntegrators The sum of all the BLFInteg...
A coefficient that is constant across space and time.
Class for domain integration .
Definition lininteg.hpp:109
FindPointsGSLIB can robustly evaluate a GridFunction on an arbitrary collection of points....
Definition gslib.hpp:53
void FindPoints(const Vector &point_pos, int point_pos_ordering=Ordering::byNODES)
Definition gslib.cpp:186
void Setup(Mesh &m, const double bb_t=0.1, const double newt_tol=1.0e-12, const int npt_max=256)
Definition gslib.cpp:109
virtual void Interpolate(const GridFunction &field_in, Vector &field_out)
Definition gslib.cpp:856
virtual const Array< unsigned int > & GetCode() const
Definition gslib.hpp:208
virtual void FreeData()
Definition gslib.cpp:283
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition fespace.hpp:220
Data type for Gauss-Seidel smoother of sparse matrix.
Class for grid function - Vector with associated FE space.
Definition gridfunc.hpp:31
Arbitrary order H1-conforming (continuous) finite elements.
Definition fe_coll.hpp:260
Vector with associated FE space and LinearFormIntegrators.
Mesh data type.
Definition mesh.hpp:56
Pointer to an Operator of a specified type.
Definition handle.hpp:34
void Parse()
Parse the command-line options. Note that this function expects all the options provided through the ...
void PrintUsage(std::ostream &out) const
Print the usage message.
void PrintOptions(std::ostream &out) const
Print the options.
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 'var' to receive the value. Enable/disable tags are used to set the bool...
Definition optparser.hpp:82
bool Good() const
Return true if the command line options were parsed successfully.
Data type sparse matrix.
Definition sparsemat.hpp:51
Vector data type.
Definition vector.hpp:80
real_t Normlinf() const
Returns the l_infinity norm of the vector.
Definition vector.cpp:850
int Size() const
Returns the size of the vector.
Definition vector.hpp:218
int dim
Definition ex24.cpp:53
int main()
const int visport
void PCG(const Operator &A, Solver &B, const Vector &b, Vector &x, int print_iter, int max_num_iter, real_t RTOLERANCE, real_t ATOLERANCE)
Preconditioned conjugate gradient method. (tolerances are squared)
Definition solvers.cpp:913
const char vishost[]
struct schwarz_common schwarz
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)