50int main(
int argc,
char *argv[])
53 const char *mesh_file_1 =
"../../data/square-disc.mesh";
54 const char *mesh_file_2 =
"../../data/inline-quad.mesh";
56 bool visualization =
true;
60 double rel_tol = 1.e-8;
63 args.
AddOption(&mesh_file_1,
"-m1",
"--mesh",
65 args.
AddOption(&mesh_file_2,
"-m2",
"--mesh",
68 "Finite element order (polynomial degree) or -1 for"
69 " isoparametric space.");
70 args.
AddOption(&visualization,
"-vis",
"--visualization",
"-no-vis",
72 "Enable or disable GLVis visualization.");
73 args.
AddOption(&r1_levels,
"-r1",
"--refine-serial",
74 "Number of times to refine the mesh uniformly in serial.");
75 args.
AddOption(&r2_levels,
"-r2",
"--refine-serial",
76 "Number of times to refine the mesh uniformly in serial.");
77 args.
AddOption(&rel_tol,
"-rt",
"--relative tolerance",
78 "Tolerance for Schwarz iteration convergence criterion.");
79 args.
AddOption(&visport,
"-p",
"--send-port",
"Socket for GLVis.");
88 const int nmeshes = 2;
93 Array <Mesh*> mesharr(2);
94 mesharr[0] =
new Mesh(mesh_file_1, 1, 1);
95 mesharr[1] =
new Mesh(mesh_file_2, 1, 1);
96 int dim = mesharr[0]->Dimension();
102 for (
int lev = 0; lev < r1_levels; lev++) { mesharr[0]->UniformRefinement(); }
103 for (
int lev = 0; lev < r2_levels; lev++) { mesharr[1]->UniformRefinement(); }
108 Array <FiniteElementCollection*> fecarr(nmeshes);
109 Array <FiniteElementSpace*> fespacearr(nmeshes);
110 for (
int i = 0; i< nmeshes; i ++)
116 else if (mesharr[i]->GetNodes())
118 fecarr[i] = mesharr[i]->GetNodes()->OwnFEC();
119 cout <<
"Using isoparametric FEs: " << fecarr[0]->Name() << endl;
133 if (mesharr[0]->bdr_attributes.
Size())
137 fespacearr[0]->GetEssentialTrueDofs(ess_bdr, ess_tdof_list1);
140 if (mesharr[1]->bdr_attributes.
Size())
144 fespacearr[1]->GetEssentialTrueDofs(ess_bdr, ess_tdof_list2);
152 for (
int i = 0; i< nmeshes; i ++)
169 Vector mesh_nodes_1 = *mesharr[0]->GetNodes();
171 Vector mesh_nodes_2 = *mesharr[1]->GetNodes();
176 if (strcmp(mesh_file_1,
"../../data/square-disc.mesh") == 0 &&
177 strcmp(mesh_file_2,
"../../data/inline-quad.mesh") == 0 )
179 for (
int i = 0; i < mesh_nodes_2.
Size(); i++)
181 mesh_nodes_2(i) = 0.5 + 0.5*(mesh_nodes_2(i)-0.5);
183 mesharr[1]->SetNodes(mesh_nodes_2);
187 finder1.
Setup(*mesharr[0]);
188 finder2.
Setup(*mesharr[1]);
190 Array<int> ess_tdof_list1_int, ess_tdof_list2_int;
192 ess_tdof_list1, ess_tdof_list2,
193 ess_tdof_list1_int, ess_tdof_list2_int,
dim);
197 const int number_boundary_1 = ess_tdof_list1_int.
Size(),
198 number_boundary_2 = ess_tdof_list2_int.
Size(),
199 number_true_1 = mesh_nodes_1.
Size()/
dim,
200 number_true_2 = mesh_nodes_2.
Size()/
dim;
202 MFEM_VERIFY(number_boundary_1 != 0 ||
203 number_boundary_2 != 0,
" Please use overlapping grids.");
206 for (
int i = 0; i < number_boundary_1; i++)
208 int idx = ess_tdof_list1_int[i];
209 for (
int d = 0; d <
dim; d++)
211 bnd1(i+d*number_boundary_1) = mesh_nodes_1(idx + d*number_true_1);
216 for (
int i = 0; i < number_boundary_2; i++)
218 int idx = ess_tdof_list2_int[i];
219 for (
int d = 0; d <
dim; d++)
221 bnd2(i+d*number_boundary_2) = mesh_nodes_2(idx + d*number_true_2);
225 Vector interp_vals1(number_boundary_1), interp_vals2(number_boundary_2);
231 Array <BilinearForm*> a_ar(2);
251 int NiterSchwarz = 100;
254 for (
int i = 0; i < nmeshes; i ++)
264 a_ar[0]->FormLinearSystem(ess_tdof_list1, x1, *b_ar[0], A1, X1, B1);
265 a_ar[1]->FormLinearSystem(ess_tdof_list2, x2, *b_ar[1], A2, X2, B2);
270 PCG(*A1, M1, B1, X1, 0, 200, 1e-12, 0.0);
272 PCG(*A2, M2, B2, X2, 0, 200, 1e-12, 0.0);
275 a_ar[0]->RecoverFEMSolution(X1, *b_ar[0], x1);
276 a_ar[1]->RecoverFEMSolution(X2, *b_ar[1], x2);
282 double dxmax = std::numeric_limits<float>::min();
283 double x1inf = x1.Normlinf();
285 for (
int i = 0; i < number_boundary_1; i++)
287 int idx = ess_tdof_list1_int[i];
288 double dx = std::abs(x1(idx)-interp_vals1(i))/x1inf;
289 if (dx > dxmax) { dxmax = dx; }
290 x1(idx) = interp_vals1(i);
293 for (
int i = 0; i < number_boundary_2; i++)
295 int idx = ess_tdof_list2_int[i];
296 double dx = std::abs(x2(idx)-interp_vals2(i))/x2inf;
297 if (dx > dxmax) { dxmax = dx; }
298 x2(idx) = interp_vals2(i);
304 std::cout << std::setprecision(8) <<
306 ", Relative residual: " << dxmax << endl;
307 if (dxmax < rel_tol) {
break; }
314 for (
int ip = 0; ip<mesharr.
Size(); ++ip)
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;
331 for (
int i = 0; i < nmeshes; i++)
334 delete fespacearr[i];
335 if (order > 0) {
delete fecarr[i]; }