MFEM v4.8.0
Finite element discretization library
Loading...
Searching...
No Matches
pminimal-surface.cpp
Go to the documentation of this file.
1// Copyright (c) 2010-2025, 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// Minimal Surface Miniapp: Compute minimal surfaces - Parallel Version
14// --------------------------------------------------------------------
15//
16// This miniapp solves Plateau's problem: the Dirichlet problem for the minimal
17// surface equation.
18//
19// Two problems can be run:
20//
21// - Problem 0 solves the minimal surface equation of parametric surfaces.
22// The surface (-s) option allow the selection of different
23// parametrization:
24// s=0: Uses the given mesh from command line options
25// s=1: Catenoid
26// s=2: Helicoid
27// s=3: Enneper
28// s=4: Hold
29// s=5: Costa
30// s=6: Shell
31// s=7: Scherk
32// s=8: FullPeach
33// s=9: QuarterPeach
34// s=10: SlottedSphere
35//
36// - Problem 1 solves the minimal surface equation of the form z=f(x,y),
37// for the Dirichlet problem, using Picard iterations:
38// -div( q grad u) = 0, with q(u) = (1 + |∇u|²)^{-1/2}
39//
40// Compile with: make pminimal-surface
41//
42// Sample runs: mpirun -np 4 pminimal-surface
43// mpirun -np 4 pminimal-surface -a
44// mpirun -np 4 pminimal-surface -c
45// mpirun -np 4 pminimal-surface -c -a
46// mpirun -np 4 pminimal-surface -no-pa
47// mpirun -np 4 pminimal-surface -no-pa -a
48// mpirun -np 4 pminimal-surface -no-pa -a -c
49// mpirun -np 4 pminimal-surface -p 1
50//
51// Device sample runs:
52// mpirun -np 4 pminimal-surface -d debug
53// mpirun -np 4 pminimal-surface -d debug -a
54// mpirun -np 4 pminimal-surface -d debug -c
55// mpirun -np 4 pminimal-surface -d debug -c -a
56// mpirun -np 4 pminimal-surface -d cuda
57// mpirun -np 4 pminimal-surface -d cuda -a
58// mpirun -np 4 pminimal-surface -d cuda -c
59// mpirun -np 4 pminimal-surface -d cuda -c -a
60// mpirun -np 4 pminimal-surface -d cuda -no-pa
61// mpirun -np 4 pminimal-surface -d cuda -no-pa -a
62// mpirun -np 4 pminimal-surface -d cuda -no-pa -c
63// mpirun -np 4 pminimal-surface -d cuda -no-pa -c -a
64
65#include "mfem.hpp"
67
68using namespace mfem;
69using namespace std;
70
71// Constant variables
72constexpr int DIM = 2;
73constexpr int SDIM = 3;
74constexpr real_t PI = M_PI;
75constexpr real_t NRM = 1.e-4;
76constexpr real_t EPS = 1.e-14;
78constexpr real_t NL_DMAX = std::numeric_limits<real_t>::max();
79
80// Static variables for GLVis
81constexpr int GLVIZ_W = 1024;
82constexpr int GLVIZ_H = 1024;
83constexpr char vishost[] = "localhost";
84
85// Context/Options for the solver
86struct Opt
87{
88 int sz, id;
89 int pb = 0;
90 int nx = 6;
91 int ny = 6;
92 int order = 3;
93 int refine = 2;
94 int niters = 8;
95 int surface = 5;
96 // Socket to send visualization data
97 int visport = 19916;
98 bool pa = true;
99 bool vis = true;
100 bool amr = false;
101 bool wait = false;
102 bool print = false;
103 bool radial = false;
104 bool by_vdim = false;
105 bool snapshot = false;
106 // bool vis_mesh = false; // Not supported by GLVis
107 real_t tau = 1.0;
108 real_t lambda = 0.1;
109 real_t amr_threshold = 0.6;
110 const char *keys = "gAm";
111 const char *device_config = "cpu";
112 const char *mesh_file = "../../data/mobius-strip.mesh";
113 void (*Tptr)(const Vector&, Vector&) = nullptr;
114};
115
116class Surface: public Mesh
117{
118protected:
119 Opt &opt;
120 ParMesh *mesh;
121 Array<int> bc;
122 socketstream glvis;
123 H1_FECollection *fec;
125public:
126 // Reading from mesh file
127 Surface(Opt &opt, const char *file): Mesh(file, true), opt(opt) { }
128
129 // Generate 2D empty surface mesh
130 Surface(Opt &opt, bool): Mesh(), opt(opt) { }
131
132 // Generate 2D quad surface mesh
133 Surface(Opt &opt)
134 : Mesh(Mesh::MakeCartesian2D(opt.nx, opt.ny, QUAD, true)), opt(opt) { }
135
136 // Generate 2D generic surface mesh
137 Surface(Opt &opt, int nv, int ne, int nbe):
138 Mesh(DIM, nv, ne, nbe, SDIM), opt(opt) { }
139
140 void Create()
141 {
142 if (opt.surface > 0)
143 {
144 Prefix();
145 Transform();
146 }
147 Postfix();
148 Refine();
149 Snap();
150 fec = new H1_FECollection(opt.order, DIM);
151 if (opt.amr) { EnsureNCMesh(); }
152 mesh = new ParMesh(MPI_COMM_WORLD, *this);
153 fes = new ParFiniteElementSpace(mesh, fec, opt.by_vdim ? 1 : SDIM);
154 BoundaryConditions();
155 }
156
157 int Solve()
158 {
159 // Initialize GLVis server if 'visualization' is set
160 if (opt.vis) { opt.vis = glvis.open(vishost, opt.visport) == 0; }
161 // Send to GLVis the first mesh
162 if (opt.vis) { Visualize(glvis, opt, mesh, GLVIZ_W, GLVIZ_H); }
163 // Create and launch the surface solver
164 if (opt.by_vdim)
165 {
166 ByVDim(*this, opt).Solve();
167 }
168 else
169 {
170 ByNodes(*this, opt).Solve();
171 }
172 if (opt.vis && opt.snapshot)
173 {
174 opt.keys = "Sq";
175 Visualize(glvis, opt, mesh, mesh->GetNodes());
176 }
177 return 0;
178 }
179
180 ~Surface() override
181 {
182 if (opt.vis) { glvis.close(); }
183 delete mesh; delete fec; delete fes;
184 }
185
186 virtual void Prefix()
187 {
188 SetCurvature(opt.order, false, SDIM, Ordering::byNODES);
189 }
190
191 virtual void Transform() { if (opt.Tptr) { Mesh::Transform(opt.Tptr); } }
192
193 virtual void Postfix()
194 {
195 SetCurvature(opt.order, false, SDIM, Ordering::byNODES);
196 }
197
198 virtual void Refine()
199 {
200 for (int l = 0; l < opt.refine; l++)
201 {
203 }
204 }
205
206 virtual void Snap()
207 {
209 for (int i = 0; i < nodes.Size(); i++)
210 {
211 if (std::abs(nodes(i)) < EPS)
212 {
213 nodes(i) = 0.0;
214 }
215 }
216 }
217
218 void SnapNodesToUnitSphere()
219 {
220 Vector node(SDIM);
222 for (int i = 0; i < nodes.FESpace()->GetNDofs(); i++)
223 {
224 for (int d = 0; d < SDIM; d++)
225 {
226 node(d) = nodes(nodes.FESpace()->DofToVDof(i, d));
227 }
228 node /= node.Norml2();
229 for (int d = 0; d < SDIM; d++)
230 {
231 nodes(nodes.FESpace()->DofToVDof(i, d)) = node(d);
232 }
233 }
234 }
235
236 virtual void BoundaryConditions()
237 {
238 if (bdr_attributes.Size())
239 {
240 Array<int> ess_bdr(bdr_attributes.Max());
241 ess_bdr = 1;
242 bc.HostReadWrite();
243 fes->GetEssentialTrueDofs(ess_bdr, bc);
244 }
245 }
246
247 // Initialize visualization of some given mesh
248 static void Visualize(socketstream &glvis,
249 Opt &opt, const Mesh *mesh,
250 const int w, const int h,
251 const GridFunction *sol = nullptr)
252 {
253 const GridFunction &solution = sol ? *sol : *mesh->GetNodes();
254 glvis << "parallel " << opt.sz << " " << opt.id << "\n";
255 glvis << "solution\n" << *mesh << solution;
256 glvis.precision(8);
257 glvis << "window_size " << w << " " << h << "\n";
258 glvis << "keys " << opt.keys << "\n";
259 opt.keys = nullptr;
260 if (opt.wait) { glvis << "pause\n"; }
261 glvis << std::flush;
262 }
263
264 // Visualize some solution on the given mesh
265 static void Visualize(socketstream &glvis,
266 const Opt &opt, const Mesh *mesh,
267 const GridFunction *sol = nullptr)
268 {
269 glvis << "parallel " << opt.sz << " " << opt.id << "\n";
270 const GridFunction &solution = sol ? *sol : *mesh->GetNodes();
271 glvis << "solution\n" << *mesh << solution;
272 if (opt.wait) { glvis << "pause\n"; }
273 if (opt.snapshot && opt.keys) { glvis << "keys " << opt.keys << "\n"; }
274 glvis << std::flush;
275 }
276
277 using Mesh::Print;
278 static void Print(const Opt &opt, ParMesh *mesh, const GridFunction *sol)
279 {
280 const char *mesh_file = "surface.mesh";
281 const char *sol_file = "sol.gf";
282 if (!opt.id)
283 {
284 mfem::out << "Printing " << mesh_file << ", " << sol_file << std::endl;
285 }
286
287 std::ostringstream mesh_name;
288 mesh_name << mesh_file << "." << std::setfill('0') << std::setw(6) << opt.id;
289 std::ofstream mesh_ofs(mesh_name.str().c_str());
290 mesh_ofs.precision(8);
291 mesh->Print(mesh_ofs);
292 mesh_ofs.close();
293
294 std::ostringstream sol_name;
295 sol_name << sol_file << "." << std::setfill('0') << std::setw(6) << opt.id;
296 std::ofstream sol_ofs(sol_name.str().c_str());
297 sol_ofs.precision(8);
298 sol->Save(sol_ofs);
299 sol_ofs.close();
300 }
301
302 // Surface Solver class
303 class Solver
304 {
305 protected:
306 Opt &opt;
307 Surface &S;
308 CGSolver cg;
309 OperatorPtr A;
311 ParGridFunction x, x0, b;
313 mfem::Solver *M = nullptr;
314 const int print_iter = -1, max_num_iter = 2000;
315 const real_t RTOLERANCE = EPS, ATOLERANCE = EPS*EPS;
316 public:
317 Solver(Surface &S, Opt &opt): opt(opt), S(S), cg(MPI_COMM_WORLD),
318 a(S.fes), x(S.fes), x0(S.fes), b(S.fes), one(1.0)
319 {
320 cg.SetRelTol(RTOLERANCE);
321 cg.SetAbsTol(ATOLERANCE);
322 cg.SetMaxIter(max_num_iter);
323 cg.SetPrintLevel(print_iter);
324 }
325
326 ~Solver() { delete M; }
327
328 void Solve()
329 {
330 constexpr bool converged = true;
331 if (opt.pa) { a.SetAssemblyLevel(AssemblyLevel::PARTIAL);}
332 for (int i=0; i < opt.niters; ++i)
333 {
334 if (opt.amr) { Amr(); }
335 if (opt.vis) { Surface::Visualize(S.glvis, opt, S.mesh); }
336 if (!opt.id) { mfem::out << "Iteration " << i << ": "; }
337 S.mesh->NodesUpdated();
338 a.Update();
339 a.Assemble();
340 if (Step() == converged) { break; }
341 }
342 if (opt.print) { Surface::Print(opt, S.mesh, S.mesh->GetNodes()); }
343 }
344
345 virtual bool Step() = 0;
346
347 protected:
348 bool Converged(const real_t rnorm) { return rnorm < NRM; }
349
350 bool ParAXeqB()
351 {
352 b = 0.0;
353 Vector X, B;
354 a.FormLinearSystem(S.bc, x, b, A, X, B);
355 if (!opt.pa) { M = new HypreBoomerAMG; }
356 if (M) { cg.SetPreconditioner(*M); }
357 cg.SetOperator(*A);
358 cg.Mult(B, X);
359 a.RecoverFEMSolution(X, b, x);
360 const bool by_vdim = opt.by_vdim;
361 GridFunction *nodes = by_vdim ? &x0 : S.fes->GetMesh()->GetNodes();
362 x.HostReadWrite();
363 nodes->HostRead();
364 real_t rnorm = nodes->DistanceTo(x) / nodes->Norml2();
365 real_t glob_norm;
366 MPI_Allreduce(&rnorm, &glob_norm, 1, MPITypeMap<real_t>::mpi_type,
367 MPI_MAX, MPI_COMM_WORLD);
368 rnorm = glob_norm;
369 if (!opt.id) { mfem::out << "rnorm = " << rnorm << std::endl; }
370 const real_t lambda = opt.lambda;
371 if (by_vdim)
372 {
373 MFEM_VERIFY(!opt.radial,"'VDim solver can't use radial option!");
374 return Converged(rnorm);
375 }
376 if (opt.radial)
377 {
378 GridFunction delta(S.fes);
379 subtract(x, *nodes, delta); // delta = x - nodes
380 // position and Δ vectors at point i
381 Vector ni(SDIM), di(SDIM);
382 for (int i = 0; i < delta.Size()/SDIM; i++)
383 {
384 // extract local vectors
385 const int ndof = S.fes->GetNDofs();
386 for (int d = 0; d < SDIM; d++)
387 {
388 ni(d) = (*nodes)(d*ndof + i);
389 di(d) = delta(d*ndof + i);
390 }
391 // project the delta vector in radial direction
392 const real_t ndotd = (ni*di) / (ni*ni);
393 di.Set(ndotd,ni);
394 // set global vectors
395 for (int d = 0; d < SDIM; d++) { delta(d*ndof + i) = di(d); }
396 }
397 add(*nodes, delta, *nodes);
398 }
399 // x = lambda*nodes + (1-lambda)*x
400 add(lambda, *nodes, (1.0-lambda), x, x);
401 return Converged(rnorm);
402 }
403
404 void Amr()
405 {
406 MFEM_VERIFY(opt.amr_threshold >= 0.0 && opt.amr_threshold <= 1.0, "");
407 Mesh *smesh = S.mesh;
409 const int NE = smesh->GetNE();
410 DenseMatrix Jadjt, Jadj(DIM, SDIM);
411 for (int e = 0; e < NE; e++)
412 {
413 real_t minW = +NL_DMAX;
414 real_t maxW = -NL_DMAX;
416 const Geometry::Type &type =
417 smesh->GetElement(e)->GetGeometryType();
418
419 const IntegrationRule *ir = &IntRules.Get(type, opt.order);
420 const int NQ = ir->GetNPoints();
421 for (int q = 0; q < NQ; q++)
422 {
423 eTr->SetIntPoint(&ir->IntPoint(q));
424 const DenseMatrix &J = eTr->Jacobian();
425 CalcAdjugate(J, Jadj);
426 Jadjt = Jadj;
427 Jadjt.Transpose();
428 const real_t w = Jadjt.Weight();
429 minW = std::min(minW, w);
430 maxW = std::max(maxW, w);
431 }
432 if (std::abs(maxW) != 0.0)
433 {
434 const real_t rho = minW / maxW;
435 MFEM_VERIFY(rho <= 1.0, "");
436 if (rho < opt.amr_threshold) { amr.Append(Refinement(e)); }
437 }
438 }
439 if (amr.Size()>0)
440 {
441 smesh->GetNodes()->HostReadWrite();
442 smesh->GeneralRefinement(amr);
443 S.fes->Update();
444 x.HostReadWrite();
445 x.Update();
446 a.Update();
447 b.HostReadWrite();
448 b.Update();
449 S.BoundaryConditions();
450 }
451 }
452 };
453
454 // Surface solver 'by vector'
455 class ByNodes: public Solver
456 {
457 public:
458 ByNodes(Surface &S, Opt &opt): Solver(S, opt)
459 { a.AddDomainIntegrator(new VectorDiffusionIntegrator(one)); }
460
461 bool Step() override
462 {
463 x = *S.fes->GetMesh()->GetNodes();
464 bool converge = ParAXeqB();
465 S.mesh->SetNodes(x);
466 return converge ? true : false;
467 }
468 };
469
470 // Surface solver 'by ByVDim'
471 class ByVDim: public Solver
472 {
473 public:
474 void SetNodes(const GridFunction &Xi, const int c)
475 {
476 auto d_Xi = Xi.Read();
477 auto d_nodes = S.fes->GetMesh()->GetNodes()->Write();
478 const int ndof = S.fes->GetNDofs();
479 mfem::forall(ndof, [=] MFEM_HOST_DEVICE (int i)
480 {
481 d_nodes[c*ndof + i] = d_Xi[i];
482 });
483 }
484
485 void GetNodes(GridFunction &Xi, const int c)
486 {
487 auto d_Xi = Xi.Write();
488 const int ndof = S.fes->GetNDofs();
489 auto d_nodes = S.fes->GetMesh()->GetNodes()->Read();
490 mfem::forall(ndof, [=] MFEM_HOST_DEVICE (int i)
491 {
492 d_Xi[i] = d_nodes[c*ndof + i];
493 });
494 }
495
496 ByVDim(Surface &S, Opt &opt): Solver(S, opt)
497 { a.AddDomainIntegrator(new DiffusionIntegrator(one)); }
498
499 bool Step() override
500 {
501 bool cvg[SDIM] {false};
502 for (int c=0; c < SDIM; ++c)
503 {
504 GetNodes(x, c);
505 x0 = x;
506 cvg[c] = ParAXeqB();
507 SetNodes(x, c);
508 }
509 const bool converged = cvg[0] && cvg[1] && cvg[2];
510 return converged ? true : false;
511 }
512 };
513};
514
515// #0: Default surface mesh file
516struct MeshFromFile: public Surface
517{
518 MeshFromFile(Opt &opt): Surface(opt, opt.mesh_file) { }
519};
520
521// #1: Catenoid surface
522struct Catenoid: public Surface
523{
524 Catenoid(Opt &opt): Surface((opt.Tptr = Parametrization, opt)) { }
525
526 void Prefix() override
527 {
528 SetCurvature(opt.order, false, SDIM, Ordering::byNODES);
529 Array<int> v2v(GetNV());
530 for (int i = 0; i < v2v.Size(); i++) { v2v[i] = i; }
531 // identify vertices on vertical lines
532 for (int j = 0; j <= opt.ny; j++)
533 {
534 const int v_old = opt.nx + j * (opt.nx + 1);
535 const int v_new = j * (opt.nx + 1);
536 v2v[v_old] = v_new;
537 }
538 // renumber elements
539 for (int i = 0; i < GetNE(); i++)
540 {
541 Element *el = GetElement(i);
542 int *v = el->GetVertices();
543 const int nv = el->GetNVertices();
544 for (int j = 0; j < nv; j++)
545 {
546 v[j] = v2v[v[j]];
547 }
548 }
549 // renumber boundary elements
550 for (int i = 0; i < GetNBE(); i++)
551 {
552 Element *el = GetBdrElement(i);
553 int *v = el->GetVertices();
554 const int nv = el->GetNVertices();
555 for (int j = 0; j < nv; j++)
556 {
557 v[j] = v2v[v[j]];
558 }
559 }
562 }
563
564 static void Parametrization(const Vector &x, Vector &p)
565 {
566 p.SetSize(SDIM);
567 // u in [0,2π] and v in [-π/6,π/6]
568 const real_t u = 2.0*PI*x[0];
569 const real_t v = PI*(x[1]-0.5)/3.;
570 p[0] = cos(u);
571 p[1] = sin(u);
572 p[2] = v;
573 }
574};
575
576// #2: Helicoid surface
577struct Helicoid: public Surface
578{
579 Helicoid(Opt &opt): Surface((opt.Tptr = Parametrization, opt)) { }
580
581 static void Parametrization(const Vector &x, Vector &p)
582 {
583 p.SetSize(SDIM);
584 // u in [0,2π] and v in [-2π/3,2π/3]
585 const real_t u = 2.0*PI*x[0];
586 const real_t v = 2.0*PI*(2.0*x[1]-1.0)/3.0;
587 p(0) = sin(u)*v;
588 p(1) = cos(u)*v;
589 p(2) = u;
590 }
591};
592
593// #3: Enneper's surface
594struct Enneper: public Surface
595{
596 Enneper(Opt &opt): Surface((opt.Tptr = Parametrization, opt)) { }
597
598 static void Parametrization(const Vector &x, Vector &p)
599 {
600 p.SetSize(SDIM);
601 // (u,v) in [-2, +2]
602 const real_t u = 4.0*(x[0]-0.5);
603 const real_t v = 4.0*(x[1]-0.5);
604 p[0] = +u - u*u*u/3.0 + u*v*v;
605 p[1] = -v - u*u*v + v*v*v/3.0;
606 p[2] = u*u - v*v;
607 }
608};
609
610// #4: Hold surface
611struct Hold: public Surface
612{
613 Hold(Opt &opt): Surface((opt.Tptr = Parametrization, opt)) { }
614
615 void Prefix() override
616 {
617 SetCurvature(opt.order, false, SDIM, Ordering::byNODES);
618 Array<int> v2v(GetNV());
619 for (int i = 0; i < v2v.Size(); i++) { v2v[i] = i; }
620 // identify vertices on vertical lines
621 for (int j = 0; j <= opt.ny; j++)
622 {
623 const int v_old = opt.nx + j * (opt.nx + 1);
624 const int v_new = j * (opt.nx + 1);
625 v2v[v_old] = v_new;
626 }
627 // renumber elements
628 for (int i = 0; i < GetNE(); i++)
629 {
630 Element *el = GetElement(i);
631 int *v = el->GetVertices();
632 const int nv = el->GetNVertices();
633 for (int j = 0; j < nv; j++)
634 {
635 v[j] = v2v[v[j]];
636 }
637 }
638 // renumber boundary elements
639 for (int i = 0; i < GetNBE(); i++)
640 {
641 Element *el = GetBdrElement(i);
642 int *v = el->GetVertices();
643 const int nv = el->GetNVertices();
644 for (int j = 0; j < nv; j++)
645 {
646 v[j] = v2v[v[j]];
647 }
648 }
651 }
652
653 static void Parametrization(const Vector &x, Vector &p)
654 {
655 p.SetSize(SDIM);
656 // u in [0,2π] and v in [0,1]
657 const real_t u = 2.0*PI*x[0];
658 const real_t v = x[1];
659 p[0] = cos(u)*(1.0 + 0.3*sin(3.*u + PI*v));
660 p[1] = sin(u)*(1.0 + 0.3*sin(3.*u + PI*v));
661 p[2] = v;
662 }
663};
664
665// #5: Costa minimal surface
666#include <complex>
667using cdouble = std::complex<real_t>;
668#define I cdouble(0.0, 1.0)
669
670// https://dlmf.nist.gov/20.2
671cdouble EllipticTheta(const int a, const cdouble u, const cdouble q)
672{
673 cdouble J = 0.0;
674 real_t delta = std::numeric_limits<real_t>::max();
675 switch (a)
676 {
677 case 1:
678 for (int n=0; delta > EPS; n+=1)
679 {
680 const cdouble j(pow(-real_t(1),real_t(n))*pow(q,
681 real_t(n*(n+1)))*sin(real_t(2*n+1)*u));
682 delta = abs(j);
683 J += j;
684 }
685 return cdouble(real_t(2)*pow(q,real_t(0.25))*J);
686
687 case 2:
688 for (int n=0; delta > EPS; n+=1)
689 {
690 const cdouble j(pow(q,real_t(n*(n+1)))*cos(real_t(2*n+1)*u));
691 delta = abs(j);
692 J += j;
693 }
694 return cdouble(real_t(2)*pow(q,real_t(0.25))*J);
695 case 3:
696 for (int n=1; delta > EPS; n+=1)
697 {
698 const cdouble j = pow(q,real_t(n*n))*cos(real_t(2*n)*u);
699 delta = abs(j);
700 J += j;
701 }
702 return real_t(1) + real_t(2)*J;
703 case 4:
704 for (int n=1; delta > EPS; n+=1)
705 {
706 const cdouble j = pow(-real_t(1),real_t(n))*pow(q,
707 real_t(n*n))*cos(real_t(2*n)*u);
708 delta = abs(j);
709 J += j;
710 }
711 return real_t(1) + real_t(2)*J;
712 }
713 return J;
714}
715
716// https://dlmf.nist.gov/23.6#E5
718 const cdouble w1 = 0.5,
719 const cdouble w3 = real_t(0.5)*I)
720{
721 const cdouble tau = w3/w1;
722 const cdouble q = exp(I*PI*tau);
723 const cdouble e1 = PI*PI/(real_t(12)*w1*w1)*
724 (pow(EllipticTheta(2,0,q),real_t(4)) +
725 real_t(2)*pow(EllipticTheta(4,0,q),real_t(4)));
726 const cdouble u = PI*z / (real_t(2)*w1);
727 const cdouble P = PI * EllipticTheta(3,0,q)*EllipticTheta(4,0,q) *
728 EllipticTheta(2,u,q)/(real_t(2)*w1*EllipticTheta(1,u,q));
729 return P*P + e1;
730}
731
732cdouble EllipticTheta1Prime(const int k, const cdouble u, const cdouble q)
733{
734 cdouble J = 0.0;
735 real_t delta = std::numeric_limits<real_t>::max();
736 for (int n=0; delta > EPS; n+=1)
737 {
738 const real_t alpha = 2.0*n+1.0;
739 const cdouble Dcosine = pow(alpha,real_t(k))*sin(k*PI/real_t(2) + alpha*u);
740 const cdouble j(pow(-real_t(1),real_t(n))*pow(q,real_t(n*(n+1)))*Dcosine);
741 delta = abs(j);
742 J += j;
743 }
744 return cdouble(real_t(2)*pow(q,real_t(0.25))*J);
745}
746
747// Logarithmic Derivative of Theta Function 1
749{
750 cdouble J = 0.0;
751 real_t delta = std::numeric_limits<real_t>::max();
752 for (int n=1; delta > EPS; n+=1)
753 {
754 cdouble q2n = pow(q, real_t(2*n));
755 if (abs(q2n) < EPS) { q2n = 0.0; }
756 const cdouble j = q2n/(real_t(1)-q2n)*sin(real_t(2*n)*u);
757 delta = abs(j);
758 J += j;
759 }
760 return real_t(1)/tan(u) + real_t(4)*J;
761}
762
763// https://dlmf.nist.gov/23.6#E13
765 const cdouble w1 = 0.5,
766 const cdouble w3 = real_t(0.5)*I)
767{
768 const cdouble tau = w3/w1;
769 const cdouble q = exp(I*PI*tau);
770 const cdouble n1 = -PI*PI/(real_t(12)*w1) *
771 (EllipticTheta1Prime(3,0,q)/
772 EllipticTheta1Prime(1,0,q));
773 const cdouble u = PI*z / (real_t(2)*w1);
774 return z*n1/w1 + PI/(real_t(2)*w1)*LogEllipticTheta1Prime(u,q);
775}
776
777// https://www.mathcurve.com/surfaces.gb/costa/costa.shtml
778static real_t ALPHA[4] {0.0};
779struct Costa: public Surface
780{
781 Costa(Opt &opt): Surface((opt.Tptr = Parametrization, opt), false) { }
782
783 void Prefix() override
784 {
785 ALPHA[3] = opt.tau;
786 const int nx = opt.nx, ny = opt.ny;
787 MFEM_VERIFY(nx>2 && ny>2, "");
788 const int nXhalf = (nx%2)==0 ? 4 : 2;
789 const int nYhalf = (ny%2)==0 ? 4 : 2;
790 const int nxh = nXhalf + nYhalf;
791 const int NVert = (nx+1) * (ny+1);
792 const int NElem = nx*ny - 4 - nxh;
793 const int NBdrElem = 0;
794 InitMesh(DIM, SDIM, NVert, NElem, NBdrElem);
795 // Sets vertices and the corresponding coordinates
796 for (int j = 0; j <= ny; j++)
797 {
798 const real_t cy = ((real_t) j / ny) ;
799 for (int i = 0; i <= nx; i++)
800 {
801 const real_t cx = ((real_t) i / nx);
802 const real_t coords[SDIM] = {cx, cy, 0.0};
803 AddVertex(coords);
804 }
805 }
806 // Sets elements and the corresponding indices of vertices
807 for (int j = 0; j < ny; j++)
808 {
809 for (int i = 0; i < nx; i++)
810 {
811 if (i == 0 && j == 0) { continue; }
812 if (i+1 == nx && j == 0) { continue; }
813 if (i == 0 && j+1 == ny) { continue; }
814 if (i+1 == nx && j+1 == ny) { continue; }
815 if ((j == 0 || j+1 == ny) && (abs(nx-(i<<1)-1)<=1)) { continue; }
816 if ((i == 0 || i+1 == nx) && (abs(ny-(j<<1)-1)<=1)) { continue; }
817 const int i0 = i + j*(nx+1);
818 const int i1 = i+1 + j*(nx+1);
819 const int i2 = i+1 + (j+1)*(nx+1);
820 const int i3 = i + (j+1)*(nx+1);
821 const int ind[4] = {i0, i1, i2, i3};
822 AddQuad(ind);
823 }
824 }
826 FinalizeQuadMesh(false, 0, true);
828 SetCurvature(opt.order, false, SDIM, Ordering::byNODES);
829 }
830
831 static void Parametrization(const Vector &X, Vector &p)
832 {
833 const real_t tau = ALPHA[3];
834 Vector x = X;
835 x -= +0.5;
836 x *= tau;
837 x -= -0.5;
838
839 p.SetSize(3);
840 const bool y_top = x[1] > 0.5;
841 const bool x_top = x[0] > 0.5;
842 real_t u = x[0];
843 real_t v = x[1];
844 if (y_top) { v = 1.0 - x[1]; }
845 if (x_top) { u = 1.0 - x[0]; }
846 const cdouble w = u + I*v;
847 const cdouble w3 = I/real_t(2);
848 const cdouble w1 = 1./2.;
849 const cdouble pw = WeierstrassP(w);
850 const cdouble e1 = WeierstrassP(0.5);
851 const cdouble zw = WeierstrassZeta(w);
852 const cdouble dw = WeierstrassZeta(w-w1) - WeierstrassZeta(w-w3);
853 p[0] = real(PI*(u+PI/(real_t(4)*e1))- zw +PI/(real_t(2)*e1)*(dw));
854 p[1] = real(PI*(v+PI/(real_t(4)*e1))-I*zw-PI*I/(real_t(2)*e1)*(dw));
855 p[2] = sqrt(PI/2.)*log(abs((pw-e1)/(pw+e1)));
856 if (y_top) { p[1] *= -1.0; }
857 if (x_top) { p[0] *= -1.0; }
858 const bool nan = std::isnan(p[0]) || std::isnan(p[1]) || std::isnan(p[2]);
859 MFEM_VERIFY(!nan, "nan");
860 ALPHA[0] = std::fmax(p[0], ALPHA[0]);
861 ALPHA[1] = std::fmax(p[1], ALPHA[1]);
862 ALPHA[2] = std::fmax(p[2], ALPHA[2]);
863 }
864
865 void Snap() override
866 {
867 Vector node(SDIM);
868 MFEM_VERIFY(ALPHA[0] > 0.0,"");
869 MFEM_VERIFY(ALPHA[1] > 0.0,"");
870 MFEM_VERIFY(ALPHA[2] > 0.0,"");
872 const real_t phi = (1.0 + sqrt(5.0)) / 2.0;
873 for (int i = 0; i < nodes.FESpace()->GetNDofs(); i++)
874 {
875 for (int d = 0; d < SDIM; d++)
876 {
877 const real_t alpha = d==2 ? phi : 1.0;
878 const int vdof = nodes.FESpace()->DofToVDof(i, d);
879 nodes(vdof) /= alpha * ALPHA[d];
880 }
881 }
882 }
883};
884
885// #6: Shell surface model
886struct Shell: public Surface
887{
888 Shell(Opt &opt): Surface((opt.niters = 1, opt.Tptr = Parametrization, opt)) { }
889
890 static void Parametrization(const Vector &x, Vector &p)
891 {
892 p.SetSize(3);
893 // u in [0,2π] and v in [-15, 6]
894 const real_t u = 2.0*PI*x[0];
895 const real_t v = 21.0*x[1]-15.0;
896 p[0] = +1.0*pow(1.16,v)*cos(v)*(1.0+cos(u));
897 p[1] = -1.0*pow(1.16,v)*sin(v)*(1.0+cos(u));
898 p[2] = -2.0*pow(1.16,v)*(1.0+sin(u));
899 }
900};
901
902// #7: Scherk's doubly periodic surface
903struct Scherk: public Surface
904{
905 static void Parametrization(const Vector &x, Vector &p)
906 {
907 p.SetSize(SDIM);
908 const real_t alpha = 0.49;
909 // (u,v) in [-απ, +απ]
910 const real_t u = alpha*PI*(2.0*x[0]-1.0);
911 const real_t v = alpha*PI*(2.0*x[1]-1.0);
912 p[0] = u;
913 p[1] = v;
914 p[2] = log(cos(v)/cos(u));
915 }
916
917 Scherk(Opt &opt): Surface((opt.Tptr = Parametrization, opt)) { }
918};
919
920// #8: Full Peach street model
921struct FullPeach: public Surface
922{
923 static constexpr int NV = 8;
924 static constexpr int NE = 6;
925
926 FullPeach(Opt &opt):
927 Surface((opt.niters = std::min(4, opt.niters), opt), NV, NE, 0) { }
928
929 void Prefix() override
930 {
931 const real_t quad_v[NV][SDIM] =
932 {
933 {-1, -1, -1}, {+1, -1, -1}, {+1, +1, -1}, {-1, +1, -1},
934 {-1, -1, +1}, {+1, -1, +1}, {+1, +1, +1}, {-1, +1, +1}
935 };
936 const int quad_e[NE][4] =
937 {
938 {3, 2, 1, 0}, {0, 1, 5, 4}, {1, 2, 6, 5},
939 {2, 3, 7, 6}, {3, 0, 4, 7}, {4, 5, 6, 7}
940
941 };
942 for (int j = 0; j < NV; j++) { AddVertex(quad_v[j]); }
943 for (int j = 0; j < NE; j++) { AddQuad(quad_e[j], j+1); }
944
945 FinalizeQuadMesh(false, 0, true);
946 FinalizeTopology(false);
948 SetCurvature(opt.order, false, SDIM, Ordering::byNODES);
949 }
950
951 void Snap() override { SnapNodesToUnitSphere(); }
952
953 void BoundaryConditions() override
954 {
955 Vector X(SDIM);
956 Array<int> dofs;
957 Array<int> ess_vdofs, ess_tdofs;
958 ess_vdofs.SetSize(fes->GetVSize());
959 MFEM_VERIFY(fes->GetVSize() >= fes->GetTrueVSize(),"");
960 ess_vdofs = 0;
961 DenseMatrix PointMat;
962 mesh->GetNodes()->HostRead();
963 for (int e = 0; e < fes->GetNE(); e++)
964 {
965 fes->GetElementDofs(e, dofs);
966 const IntegrationRule &ir = fes->GetFE(e)->GetNodes();
968 eTr->Transform(ir, PointMat);
969 Vector one(dofs.Size());
970 for (int dof = 0; dof < dofs.Size(); dof++)
971 {
972 one = 0.0;
973 one[dof] = 1.0;
974 const int k = dofs[dof];
975 MFEM_ASSERT(k >= 0, "");
976 PointMat.Mult(one, X);
977 const bool halfX = std::abs(X[0]) < EPS && X[1] <= 0.0;
978 const bool halfY = std::abs(X[2]) < EPS && X[1] >= 0.0;
979 const bool is_on_bc = halfX || halfY;
980 for (int c = 0; c < SDIM; c++)
981 { ess_vdofs[fes->DofToVDof(k, c)] = is_on_bc; }
982 }
983 }
984 const SparseMatrix *R = fes->GetRestrictionMatrix();
985 if (!R)
986 {
987 ess_tdofs.MakeRef(ess_vdofs);
988 }
989 else
990 {
991 R->BooleanMult(ess_vdofs, ess_tdofs);
992 }
993 bc.HostReadWrite();
994 ParFiniteElementSpace::MarkerToList(ess_tdofs, bc);
995 }
996};
997
998// #9: 1/4th Peach street model
999struct QuarterPeach: public Surface
1000{
1001 QuarterPeach(Opt &opt): Surface((opt.Tptr = Parametrization, opt)) { }
1002
1003 static void Parametrization(const Vector &X, Vector &p)
1004 {
1005 p = X;
1006 const real_t x = 2.0*X[0]-1.0;
1007 const real_t y = X[1];
1008 const real_t r = sqrt(x*x + y*y);
1009 const real_t t = (x==0.0) ? PI/2.0 :
1010 (y==0.0 && x>0.0) ? 0. :
1011 (y==0.0 && x<0.0) ? PI : acos(x/r);
1012 const real_t sqrtx = sqrt(1.0 + x*x);
1013 const real_t sqrty = sqrt(1.0 + y*y);
1014 const bool yaxis = PI/4.0<t && t < 3.0*PI/4.0;
1015 const real_t R = yaxis?sqrtx:sqrty;
1016 const real_t gamma = r/R;
1017 p[0] = gamma * cos(t);
1018 p[1] = gamma * sin(t);
1019 p[2] = 1.0 - gamma;
1020 }
1021
1022 void Postfix() override
1023 {
1024 for (int i = 0; i < GetNBE(); i++)
1025 {
1026 Element *el = GetBdrElement(i);
1027 const int fn = GetBdrElementFaceIndex(i);
1028 MFEM_VERIFY(!FaceIsTrueInterior(fn),"");
1031 const GridFunction *nodes = GetNodes();
1032 Vector nval;
1033 real_t R[2], X[2][SDIM];
1034 for (int v = 0; v < 2; v++)
1035 {
1036 R[v] = 0.0;
1037 const int iv = vertices[v];
1038 for (int d = 0; d < SDIM; d++)
1039 {
1040 nodes->GetNodalValues(nval, d+1);
1041 const real_t x = X[v][d] = nval[iv];
1042 if (d < 2) { R[v] += x*x; }
1043 }
1044 }
1045 if (std::abs(X[0][1])<=EPS && std::abs(X[1][1])<=EPS &&
1046 (R[0]>0.1 || R[1]>0.1))
1047 { el->SetAttribute(1); }
1048 else { el->SetAttribute(2); }
1049 }
1050 }
1051};
1052
1053// #10: Slotted sphere mesh
1054struct SlottedSphere: public Surface
1055{
1056 SlottedSphere(Opt &opt): Surface((opt.niters = 4, opt), 64, 40, 0) { }
1057
1058 void Prefix() override
1059 {
1060 constexpr real_t delta = 0.15;
1061 constexpr int NV1D = 4;
1062 constexpr int NV = NV1D*NV1D*NV1D;
1063 constexpr int NEPF = (NV1D-1)*(NV1D-1);
1064 constexpr int NE = NEPF*6;
1065 const real_t V1D[NV1D] = {-1.0, -delta, delta, 1.0};
1066 real_t QV[NV][3];
1067 for (int iv=0; iv<NV; ++iv)
1068 {
1069 int ix = iv % NV1D;
1070 int iy = (iv / NV1D) % NV1D;
1071 int iz = (iv / NV1D) / NV1D;
1072
1073 QV[iv][0] = V1D[ix];
1074 QV[iv][1] = V1D[iy];
1075 QV[iv][2] = V1D[iz];
1076 }
1077 int QE[NE][4];
1078 for (int ix=0; ix<NV1D-1; ++ix)
1079 {
1080 for (int iy=0; iy<NV1D-1; ++iy)
1081 {
1082 int el_offset = ix + iy*(NV1D-1);
1083 // x = 0
1084 QE[0*NEPF + el_offset][0] = NV1D*ix + NV1D*NV1D*iy;
1085 QE[0*NEPF + el_offset][1] = NV1D*(ix+1) + NV1D*NV1D*iy;
1086 QE[0*NEPF + el_offset][2] = NV1D*(ix+1) + NV1D*NV1D*(iy+1);
1087 QE[0*NEPF + el_offset][3] = NV1D*ix + NV1D*NV1D*(iy+1);
1088 // x = 1
1089 int x_off = NV1D-1;
1090 QE[1*NEPF + el_offset][3] = x_off + NV1D*ix + NV1D*NV1D*iy;
1091 QE[1*NEPF + el_offset][2] = x_off + NV1D*(ix+1) + NV1D*NV1D*iy;
1092 QE[1*NEPF + el_offset][1] = x_off + NV1D*(ix+1) + NV1D*NV1D*(iy+1);
1093 QE[1*NEPF + el_offset][0] = x_off + NV1D*ix + NV1D*NV1D*(iy+1);
1094 // y = 0
1095 QE[2*NEPF + el_offset][0] = NV1D*NV1D*iy + ix;
1096 QE[2*NEPF + el_offset][1] = NV1D*NV1D*iy + ix + 1;
1097 QE[2*NEPF + el_offset][2] = NV1D*NV1D*(iy+1) + ix + 1;
1098 QE[2*NEPF + el_offset][3] = NV1D*NV1D*(iy+1) + ix;
1099 // y = 1
1100 int y_off = NV1D*(NV1D-1);
1101 QE[3*NEPF + el_offset][0] = y_off + NV1D*NV1D*iy + ix;
1102 QE[3*NEPF + el_offset][1] = y_off + NV1D*NV1D*iy + ix + 1;
1103 QE[3*NEPF + el_offset][2] = y_off + NV1D*NV1D*(iy+1) + ix + 1;
1104 QE[3*NEPF + el_offset][3] = y_off + NV1D*NV1D*(iy+1) + ix;
1105 // z = 0
1106 QE[4*NEPF + el_offset][0] = NV1D*iy + ix;
1107 QE[4*NEPF + el_offset][1] = NV1D*iy + ix + 1;
1108 QE[4*NEPF + el_offset][2] = NV1D*(iy+1) + ix + 1;
1109 QE[4*NEPF + el_offset][3] = NV1D*(iy+1) + ix;
1110 // z = 1
1111 int z_off = NV1D*NV1D*(NV1D-1);
1112 QE[5*NEPF + el_offset][0] = z_off + NV1D*iy + ix;
1113 QE[5*NEPF + el_offset][1] = z_off + NV1D*iy + ix + 1;
1114 QE[5*NEPF + el_offset][2] = z_off + NV1D*(iy+1) + ix + 1;
1115 QE[5*NEPF + el_offset][3] = z_off + NV1D*(iy+1) + ix;
1116 }
1117 }
1118 // Delete on x = 0 face
1119 QE[0*NEPF + 1 + 2*(NV1D-1)][0] = -1;
1120 QE[0*NEPF + 1 + 1*(NV1D-1)][0] = -1;
1121 // Delete on x = 1 face
1122 QE[1*NEPF + 1 + 2*(NV1D-1)][0] = -1;
1123 QE[1*NEPF + 1 + 1*(NV1D-1)][0] = -1;
1124 // Delete on y = 1 face
1125 QE[3*NEPF + 1 + 0*(NV1D-1)][0] = -1;
1126 QE[3*NEPF + 1 + 1*(NV1D-1)][0] = -1;
1127 // Delete on z = 1 face
1128 QE[5*NEPF + 0 + 1*(NV1D-1)][0] = -1;
1129 QE[5*NEPF + 1 + 1*(NV1D-1)][0] = -1;
1130 QE[5*NEPF + 2 + 1*(NV1D-1)][0] = -1;
1131 // Delete on z = 0 face
1132 QE[4*NEPF + 1 + 0*(NV1D-1)][0] = -1;
1133 QE[4*NEPF + 1 + 1*(NV1D-1)][0] = -1;
1134 QE[4*NEPF + 1 + 2*(NV1D-1)][0] = -1;
1135 // Delete on y = 0 face
1136 QE[2*NEPF + 1 + 0*(NV1D-1)][0] = -1;
1137 QE[2*NEPF + 1 + 1*(NV1D-1)][0] = -1;
1138 for (int j = 0; j < NV; j++) { AddVertex(QV[j]); }
1139 for (int j = 0; j < NE; j++)
1140 {
1141 if (QE[j][0] < 0) { continue; }
1142 AddQuad(QE[j], j+1);
1143 }
1145 FinalizeQuadMesh(false, 0, true);
1146 EnsureNodes();
1148 }
1149
1150 void Snap() override { SnapNodesToUnitSphere(); }
1151};
1152
1153static int Problem0(Opt &opt)
1154{
1155 // Create our surface mesh from command line options
1156 Surface *S = nullptr;
1157 switch (opt.surface)
1158 {
1159 case 0: S = new MeshFromFile(opt); break;
1160 case 1: S = new Catenoid(opt); break;
1161 case 2: S = new Helicoid(opt); break;
1162 case 3: S = new Enneper(opt); break;
1163 case 4: S = new Hold(opt); break;
1164 case 5: S = new Costa(opt); break;
1165 case 6: S = new Shell(opt); break;
1166 case 7: S = new Scherk(opt); break;
1167 case 8: S = new FullPeach(opt); break;
1168 case 9: S = new QuarterPeach(opt); break;
1169 case 10: S = new SlottedSphere(opt); break;
1170 default: MFEM_ABORT("Unknown surface (surface <= 10)!");
1171 }
1172 S->Create();
1173 S->Solve();
1174 delete S;
1175 return 0;
1176}
1177
1178// Problem 1: solve the Dirichlet problem for the minimal surface equation
1179// of the form z=f(x,y), using Picard iterations.
1180static real_t u0(const Vector &x) { return sin(3.0 * PI * (x[1] + x[0])); }
1181
1182enum {NORM, AREA};
1183
1184static real_t qf(const int order, const int ker, Mesh &m,
1186{
1188 const IntegrationRule &ir(IntRules.Get(type, order));
1190
1191 const int NE(m.GetNE());
1192 const int ND(fes.GetTypicalFE()->GetDof());
1193 const int NQ(ir.GetNPoints());
1195 const GeometricFactors *geom = m.GetGeometricFactors(ir, flags);
1196
1197 const int D1D = fes.GetTypicalFE()->GetOrder() + 1;
1198 const int Q1D = IntRules.Get(Geometry::SEGMENT, ir.GetOrder()).GetNPoints();
1199 MFEM_VERIFY(ND == D1D*D1D, "");
1200 MFEM_VERIFY(NQ == Q1D*Q1D, "");
1201
1202 Vector Eu(ND*NE), grad_u(DIM*NQ*NE), sum(NE*NQ), one(NE*NQ);
1203 qi->SetOutputLayout(QVectorLayout::byVDIM);
1204 const ElementDofOrdering e_ordering = ElementDofOrdering::LEXICOGRAPHIC;
1205 const Operator *G(fes.GetElementRestriction(e_ordering));
1206 G->Mult(u, Eu);
1207 qi->Derivatives(Eu, grad_u);
1208
1209 auto W = Reshape(ir.GetWeights().Read(), Q1D, Q1D);
1210 auto J = Reshape(geom->J.Read(), Q1D, Q1D, DIM, DIM, NE);
1211 auto detJ = Reshape(geom->detJ.Read(), Q1D, Q1D, NE);
1212 auto grdU = Reshape(grad_u.Read(), DIM, Q1D, Q1D, NE);
1213 auto S = Reshape(sum.Write(), Q1D, Q1D, NE);
1214
1215 mfem::forall_2D(NE, Q1D, Q1D, [=] MFEM_HOST_DEVICE (int e)
1216 {
1217 MFEM_FOREACH_THREAD(qy,y,Q1D)
1218 {
1219 MFEM_FOREACH_THREAD(qx,x,Q1D)
1220 {
1221 const real_t w = W(qx, qy);
1222 const real_t J11 = J(qx, qy, 0, 0, e);
1223 const real_t J12 = J(qx, qy, 1, 0, e);
1224 const real_t J21 = J(qx, qy, 0, 1, e);
1225 const real_t J22 = J(qx, qy, 1, 1, e);
1226 const real_t det = detJ(qx, qy, e);
1227 const real_t area = w * det;
1228 const real_t gu0 = grdU(0, qx, qy, e);
1229 const real_t gu1 = grdU(1, qx, qy, e);
1230 const real_t tgu0 = (J22*gu0 - J12*gu1)/det;
1231 const real_t tgu1 = (J11*gu1 - J21*gu0)/det;
1232 const real_t ngu = tgu0*tgu0 + tgu1*tgu1;
1233 const real_t s = (ker == AREA) ? sqrt(1.0 + ngu) :
1234 (ker == NORM) ? ngu : 0.0;
1235 S(qx, qy, e) = area * s;
1236 }
1237 }
1238 });
1239 one = 1.0;
1240 return sum * one;
1241}
1242
1243static int Problem1(Opt &opt)
1244{
1245 const int order = opt.order;
1246 Mesh smesh = Mesh::MakeCartesian2D(opt.nx, opt.ny, QUAD);
1247 smesh.SetCurvature(opt.order, false, DIM, Ordering::byNODES);
1248 for (int l = 0; l < opt.refine; l++) { smesh.UniformRefinement(); }
1249 ParMesh mesh(MPI_COMM_WORLD, smesh);
1250 const H1_FECollection fec(order, DIM);
1251 ParFiniteElementSpace fes(&mesh, &fec);
1252 Array<int> ess_tdof_list;
1253 if (mesh.bdr_attributes.Size())
1254 {
1255 Array<int> ess_bdr(mesh.bdr_attributes.Max());
1256 ess_bdr = 1;
1257 fes.GetEssentialTrueDofs(ess_bdr, ess_tdof_list);
1258 }
1259 ParGridFunction uold(&fes), u(&fes), b(&fes);
1260 FunctionCoefficient u0_fc(u0);
1261 u.ProjectCoefficient(u0_fc);
1262 socketstream glvis;
1263 if (opt.vis) { opt.vis = glvis.open(vishost, opt.visport) == 0; }
1264 if (opt.vis) { Surface::Visualize(glvis, opt, &mesh, GLVIZ_W, GLVIZ_H, &u); }
1265 Vector B, X;
1266 OperatorPtr A;
1267 CGSolver cg(MPI_COMM_WORLD);
1268 cg.SetRelTol(EPS);
1269 cg.SetAbsTol(EPS*EPS);
1270 cg.SetMaxIter(400);
1271 cg.SetPrintLevel(0);
1272 ParGridFunction eps(&fes);
1273 for (int i = 0; i < opt.niters; i++)
1274 {
1275 b = 0.0;
1276 uold = u;
1277 ParBilinearForm a(&fes);
1278 if (opt.pa) { a.SetAssemblyLevel(AssemblyLevel::PARTIAL); }
1279 const real_t q_uold = qf(order, AREA, mesh, fes, uold);
1280 MFEM_VERIFY(std::abs(q_uold) > EPS,"");
1281 ConstantCoefficient q_uold_cc(1.0/sqrt(q_uold));
1282 a.AddDomainIntegrator(new DiffusionIntegrator(q_uold_cc));
1283 a.Assemble();
1284 a.FormLinearSystem(ess_tdof_list, u, b, A, X, B);
1285 cg.SetOperator(*A);
1286 cg.Mult(B, X);
1287 a.RecoverFEMSolution(X, b, u);
1288 subtract(u, uold, eps);
1289 const real_t norm = sqrt(std::abs(qf(order, NORM, mesh, fes, eps)));
1290 const real_t area = qf(order, AREA, mesh, fes, u);
1291 if (!opt.id)
1292 {
1293 mfem::out << "Iteration " << i << ", norm: " << norm
1294 << ", area: " << area << std::endl;
1295 }
1296 if (opt.vis) { Surface::Visualize(glvis, opt, &mesh, &u); }
1297 if (opt.print) { Surface::Print(opt, &mesh, &u); }
1298 if (norm < NRM) { break; }
1299 }
1300 return 0;
1301}
1302
1303int main(int argc, char *argv[])
1304{
1305 Opt opt;
1306 Mpi::Init(argc, argv);
1307 opt.id = Mpi::WorldRank();
1308 opt.sz = Mpi::WorldSize();
1309 Hypre::Init();
1310
1311 // Parse command-line options.
1312 OptionsParser args(argc, argv);
1313 args.AddOption(&opt.visport, "-p", "--send-port", "Socket for GLVis.");
1314 args.AddOption(&opt.pb, "-p", "--problem", "Problem to solve.");
1315 args.AddOption(&opt.mesh_file, "-m", "--mesh", "Mesh file to use.");
1316 args.AddOption(&opt.wait, "-w", "--wait", "-no-w", "--no-wait",
1317 "Enable or disable a GLVis pause.");
1318 args.AddOption(&opt.radial, "-rad", "--radial", "-no-rad", "--no-radial",
1319 "Enable or disable radial constraints in solver.");
1320 args.AddOption(&opt.nx, "-x", "--num-elements-x",
1321 "Number of elements in x-direction.");
1322 args.AddOption(&opt.ny, "-y", "--num-elements-y",
1323 "Number of elements in y-direction.");
1324 args.AddOption(&opt.order, "-o", "--order", "Finite element order.");
1325 args.AddOption(&opt.refine, "-r", "--ref-levels", "Refinement");
1326 args.AddOption(&opt.niters, "-n", "--niter-max", "Max number of iterations");
1327 args.AddOption(&opt.surface, "-s", "--surface", "Choice of the surface.");
1328 args.AddOption(&opt.pa, "-pa", "--partial-assembly", "-no-pa",
1329 "--no-partial-assembly", "Enable Partial Assembly.");
1330 args.AddOption(&opt.tau, "-t", "--tau", "Costa scale factor.");
1331 args.AddOption(&opt.lambda, "-l", "--lambda", "Lambda step toward solution.");
1332 args.AddOption(&opt.amr, "-a", "--amr", "-no-a", "--no-amr", "Enable AMR.");
1333 args.AddOption(&opt.amr_threshold, "-at", "--amr-threshold", "AMR threshold.");
1334 args.AddOption(&opt.device_config, "-d", "--device",
1335 "Device configuration string, see Device::Configure().");
1336 args.AddOption(&opt.keys, "-k", "--keys", "GLVis configuration keys.");
1337 args.AddOption(&opt.vis, "-vis", "--visualization", "-no-vis",
1338 "--no-visualization", "Enable or disable visualization.");
1339 args.AddOption(&opt.by_vdim, "-c", "--solve-byvdim",
1340 "-no-c", "--solve-bynodes",
1341 "Enable or disable the 'ByVdim' solver");
1342 args.AddOption(&opt.print, "-print", "--print", "-no-print", "--no-print",
1343 "Enable or disable result output (files in mfem format).");
1344 args.AddOption(&opt.snapshot, "-ss", "--snapshot", "-no-ss", "--no-snapshot",
1345 "Enable or disable GLVis snapshot.");
1346 args.AddOption(&opt.visport, "-p", "--send-port", "Socket for GLVis.");
1347 args.Parse();
1348 if (!args.Good()) { args.PrintUsage(mfem::out); return 1; }
1349 MFEM_VERIFY(opt.lambda >= 0.0 && opt.lambda <= 1.0,"");
1350 if (!opt.id) { args.PrintOptions(mfem::out); }
1351
1352 // Initialize hardware devices
1353 Device device(opt.device_config);
1354 if (!opt.id) { device.Print(); }
1355
1356 if (opt.pb == 0) { Problem0(opt); }
1357
1358 if (opt.pb == 1) { Problem1(opt); }
1359
1360 return 0;
1361}
T Max() const
Find the maximal element in the array, using the comparison operator < for class T.
Definition array.cpp:68
void SetSize(int nsize)
Change the logical size of the array, keep existing entries.
Definition array.hpp:758
int Size() const
Return the logical size of the array.
Definition array.hpp:147
void MakeRef(T *data_, int size_, bool own_data=false)
Make this Array a reference to a pointer.
Definition array.hpp:943
int Append(const T &el)
Append element 'el' to array, resize if necessary.
Definition array.hpp:830
const T * Read(bool on_dev=true) const
Shortcut for mfem::Read(a.GetMemory(), a.Size(), on_dev).
Definition array.hpp:337
T * HostReadWrite()
Shortcut for mfem::ReadWrite(a.GetMemory(), a.Size(), false).
Definition array.hpp:357
Conjugate gradient method.
Definition solvers.hpp:538
void Mult(const Vector &b, Vector &x) const override
Iterative solution of the linear system using the Conjugate Gradient method.
Definition solvers.cpp:751
void SetOperator(const Operator &op) override
Set/update the solver for the given operator.
Definition solvers.hpp:551
A coefficient that is constant across space and time.
Data type dense matrix using column-major storage.
Definition densemat.hpp:24
void Mult(const real_t *x, real_t *y) const
Matrix vector multiplication.
Definition densemat.cpp:125
void Transpose()
(*this) = (*this)^t
real_t Weight() const
Definition densemat.cpp:573
The MFEM Device class abstracts hardware devices such as GPUs, as well as programming models such as ...
Definition device.hpp:123
void Print(std::ostream &out=mfem::out)
Print the configuration of the MFEM virtual device object.
Definition device.cpp:297
const DenseMatrix & Jacobian()
Return the Jacobian matrix of the transformation at the currently set IntegrationPoint,...
Definition eltrans.hpp:132
void SetIntPoint(const IntegrationPoint *ip)
Set the integration point ip that weights and Jacobians will be evaluated at.
Definition eltrans.hpp:106
virtual void Transform(const IntegrationPoint &, Vector &)=0
Transform integration point from reference coordinates to physical coordinates and store them in the ...
Abstract data type element.
Definition element.hpp:29
Geometry::Type GetGeometryType() const
Definition element.hpp:55
virtual void GetVertices(Array< int > &v) const =0
Get the indices defining the vertices.
void SetAttribute(const int attr)
Set element's attribute.
Definition element.hpp:61
Type
Constants for the classes derived from Element.
Definition element.hpp:41
virtual int GetNVertices() const =0
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition fespace.hpp:244
int GetNDofs() const
Returns number of degrees of freedom. This is the number of Local Degrees of Freedom.
Definition fespace.hpp:845
virtual void GetEssentialTrueDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list, int component=-1) const
Get a list of essential true dofs, ess_tdof_list, corresponding to the boundary attributes marked in ...
Definition fespace.cpp:658
const QuadratureInterpolator * GetQuadratureInterpolator(const IntegrationRule &ir) const
Return a QuadratureInterpolator that interpolates E-vectors to quadrature point values and/or derivat...
Definition fespace.cpp:1582
int GetNE() const
Returns number of elements in the mesh.
Definition fespace.hpp:891
const ElementRestrictionOperator * GetElementRestriction(ElementDofOrdering e_ordering) const
Return an Operator that converts L-vectors to E-vectors.
Definition fespace.cpp:1510
virtual void Update(bool want_transform=true)
Reflect changes in the mesh: update number of DOFs, etc. Also, calculate GridFunction transformation ...
Definition fespace.cpp:4157
Mesh * GetMesh() const
Returns the mesh.
Definition fespace.hpp:679
int GetVSize() const
Return the number of vector dofs, i.e. GetNDofs() x GetVDim().
Definition fespace.hpp:848
const FiniteElement * GetTypicalFE() const
Return GetFE(0) if the local mesh is not empty; otherwise return a typical FE based on the Geometry t...
Definition fespace.cpp:3871
int DofToVDof(int dof, int vd, int ndofs=-1) const
Compute a single vdof corresponding to the index dof and the vector index vd.
Definition fespace.cpp:294
int GetOrder() const
Returns the order of the finite element. In the case of anisotropic orders, returns the maximum order...
Definition fe_base.hpp:338
const IntegrationRule & GetNodes() const
Get a const reference to the nodes of the element.
Definition fe_base.hpp:400
int GetDof() const
Returns the number of degrees of freedom in the finite element.
Definition fe_base.hpp:334
A general function coefficient.
Structure for storing mesh geometric factors: coordinates, Jacobians, and determinants of the Jacobia...
Definition mesh.hpp:2937
Vector detJ
Determinants of the Jacobians at all quadrature points.
Definition mesh.hpp:2982
Vector J
Jacobians of the element transformations at all quadrature points.
Definition mesh.hpp:2976
Class for grid function - Vector with associated FE space.
Definition gridfunc.hpp:31
virtual void Save(std::ostream &out) const
Save the GridFunction to an output stream.
Arbitrary order H1-conforming (continuous) finite elements.
Definition fe_coll.hpp:275
The BoomerAMG solver in hypre.
Definition hypre.hpp:1717
static void Init()
Initialize hypre by calling HYPRE_Init() and set default options. After calling Hypre::Init(),...
Definition hypre.cpp:33
Class for an integration rule - an Array of IntegrationPoint.
Definition intrules.hpp:100
int GetOrder() const
Returns the order of the integration rule.
Definition intrules.hpp:249
int GetNPoints() const
Returns the number of the points in the integration rule.
Definition intrules.hpp:256
const Array< real_t > & GetWeights() const
Return the quadrature weights in a contiguous array.
Definition intrules.cpp:86
IntegrationPoint & IntPoint(int i)
Returns a reference to the i-th integration point.
Definition intrules.hpp:259
const IntegrationRule & Get(int GeomType, int Order)
Returns an integration rule for given GeomType and Order.
void SetRelTol(real_t rtol)
Definition solvers.hpp:229
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
Definition solvers.cpp:174
virtual void SetPrintLevel(int print_lvl)
Legacy method to set the level of verbosity of the solver output.
Definition solvers.cpp:72
void SetMaxIter(int max_it)
Definition solvers.hpp:231
void SetAbsTol(real_t atol)
Definition solvers.hpp:230
Mesh data type.
Definition mesh.hpp:64
Array< Vertex > vertices
Definition mesh.hpp:105
void InitMesh(int Dim_, int spaceDim_, int NVert, int NElem, int NBdrElem)
Begin construction of a mesh.
Definition mesh.cpp:1850
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
Definition mesh.hpp:290
void GeneralRefinement(const Array< Refinement > &refinements, int nonconforming=-1, int nc_limit=0)
Definition mesh.cpp:10883
void NodesUpdated()
This function should be called after the mesh node coordinates have been updated externally,...
Definition mesh.hpp:2211
void EnsureNodes()
Make sure that the mesh has valid nodes, i.e. its geometry is described by a vector finite element gr...
Definition mesh.cpp:6432
int AddQuad(int v1, int v2, int v3, int v4, int attr=1)
Adds a quadrilateral to the mesh given by 4 vertices v1 through v4.
Definition mesh.cpp:1958
Geometry::Type GetTypicalElementGeometry() const
If the local mesh is not empty, return GetElementGeometry(0); otherwise, return a typical Geometry pr...
Definition mesh.cpp:1528
int GetBdrElementFaceIndex(int be_idx) const
Return the local face (codimension-1) index for the given boundary element index.
Definition mesh.hpp:1588
const Element * GetElement(int i) const
Return pointer to the i'th element object.
Definition mesh.hpp:1339
bool FaceIsTrueInterior(int FaceNo) const
Definition mesh.hpp:571
void FinalizeTopology(bool generate_bdr=true)
Finalize the construction of the secondary topology (connectivity) data of a Mesh.
Definition mesh.cpp:3362
virtual void Print(std::ostream &os=mfem::out, const std::string &comments="") const
Definition mesh.hpp:2433
int AddVertex(real_t x, real_t y=0.0, real_t z=0.0)
Definition mesh.cpp:1873
int GetNE() const
Returns number of elements.
Definition mesh.hpp:1282
const Element * GetBdrElement(int i) const
Return pointer to the i'th boundary element object.
Definition mesh.hpp:1354
void GetElementTransformation(int i, IsoparametricTransformation *ElTr) const
Builds the transformation defining the i-th element in ElTr. ElTr must be allocated in advance and wi...
Definition mesh.cpp:357
void FinalizeQuadMesh(int generate_edges=0, int refine=0, bool fix_orientation=true)
Finalize the construction of a quadrilateral Mesh.
Definition mesh.cpp:2404
const GeometricFactors * GetGeometricFactors(const IntegrationRule &ir, const int flags, MemoryType d_mt=MemoryType::DEFAULT)
Return the mesh geometric factors corresponding to the given integration rule.
Definition mesh.cpp:880
void GetNodes(Vector &node_coord) const
Definition mesh.cpp:9294
int GetNV() const
Returns number of vertices. Vertices are only at the corners of elements, where you would expect them...
Definition mesh.hpp:1279
void GetFaceVertices(int i, Array< int > &vert) const
Returns the indices of the vertices of face i.
Definition mesh.hpp:1526
friend class ParMesh
Definition mesh.hpp:68
int GetNBE() const
Returns number of boundary elements.
Definition mesh.hpp:1285
GridFunction * GetNodes()
Return a pointer to the internal node GridFunction (may be NULL).
Definition mesh.hpp:2229
void EnsureNCMesh(bool simplices_nonconforming=false)
Definition mesh.cpp:10951
static Mesh MakeCartesian2D(int nx, int ny, Element::Type type, bool generate_edges=false, real_t sx=1.0, real_t sy=1.0, bool sfc_ordering=true)
Creates mesh for the rectangle [0,sx]x[0,sy], divided into nx*ny quadrilaterals if type = QUADRILATER...
Definition mesh.cpp:4471
void SetNodes(const Vector &node_coord)
Updates the vertex/node locations. Invokes NodesUpdated().
Definition mesh.cpp:9306
virtual void SetCurvature(int order, bool discont=false, int space_dim=-1, int ordering=1)
Set the curvature of the mesh nodes using the given polynomial degree.
Definition mesh.cpp:6484
void Transform(void(*f)(const Vector &, Vector &))
Definition mesh.cpp:13146
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition mesh.cpp:11295
void RemoveInternalBoundaries()
Definition mesh.cpp:13304
void RemoveUnusedVertices()
Remove unused vertices and rebuild mesh connectivity.
Definition mesh.cpp:13197
static int WorldRank()
Return the MPI rank in MPI_COMM_WORLD.
static int WorldSize()
Return the size of MPI_COMM_WORLD.
static void Init(int &argc, char **&argv, int required=default_thread_required, int *provided=nullptr)
Singleton creation with Mpi::Init(argc, argv).
Pointer to an Operator of a specified type.
Definition handle.hpp:34
Abstract operator.
Definition operator.hpp:25
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.
Class for parallel bilinear form.
Abstract parallel finite element space.
Definition pfespace.hpp:29
void GetEssentialTrueDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list, int component=-1) const override
int GetTrueVSize() const override
Return the number of local vector true dofs.
Definition pfespace.hpp:350
void GetElementDofs(int i, Array< int > &dofs, DofTransformation &doftrans) const override
The same as GetElementDofs(), but with a user-allocated DofTransformation object. doftrans must be al...
Definition pfespace.cpp:561
const SparseMatrix * GetRestrictionMatrix() const override
Get the R matrix which restricts a local dof vector to true dof vector.
Definition pfespace.hpp:467
const FiniteElement * GetFE(int i) const override
Definition pfespace.cpp:627
Class for parallel grid function.
Definition pgridfunc.hpp:50
void Update() override
Transform by the Space UpdateMatrix (e.g., on Mesh change).
Definition pgridfunc.cpp:91
Class for parallel meshes.
Definition pmesh.hpp:34
void Print(std::ostream &out=mfem::out, const std::string &comments="") const override
Definition pmesh.cpp:4800
A class that performs interpolation from an E-vector to quadrature point values and/or derivatives (Q...
Base class for solvers.
Definition operator.hpp:780
Data type sparse matrix.
Definition sparsemat.hpp:51
void BooleanMult(const Array< int > &x, Array< int > &y) const
y = A * x, treating all entries as booleans (zero=false, nonzero=true).
Vector data type.
Definition vector.hpp:82
virtual const real_t * Read(bool on_dev=true) const
Shortcut for mfem::Read(vec.GetMemory(), vec.Size(), on_dev).
Definition vector.hpp:494
void SetSize(int s)
Resize the vector to size s.
Definition vector.hpp:558
virtual real_t * HostReadWrite()
Shortcut for mfem::ReadWrite(vec.GetMemory(), vec.Size(), false).
Definition vector.hpp:514
virtual real_t * Write(bool on_dev=true)
Shortcut for mfem::Write(vec.GetMemory(), vec.Size(), on_dev).
Definition vector.hpp:502
int open(const char hostname[], int port)
Open the socket stream on 'port' at 'hostname'.
int close()
Close the socketstream.
const real_t alpha
Definition ex15.cpp:369
real_t Surface()
Analytic surface integral.
Definition ex38.cpp:89
int main()
real_t b
Definition lissajous.cpp:42
real_t delta
Definition lissajous.cpp:43
real_t a
Definition lissajous.cpp:41
std::complex< real_t > cdouble
real_t u(const Vector &xvec)
Definition lor_mms.hpp:22
OutStream out(std::cout)
Global stream used by the library for standard output. Initially it uses the same std::streambuf as s...
Definition globals.hpp:66
void add(const Vector &v1, const Vector &v2, Vector &v)
Definition vector.cpp:391
MFEM_HOST_DEVICE DeviceTensor< sizeof...(Dims), T > Reshape(T *ptr, Dims... dims)
Wrap a pointer as a DeviceTensor with automatically deduced template parameters.
Definition dtensor.hpp:131
void CalcAdjugate(const DenseMatrix &a, DenseMatrix &adja)
void forall_2D(int N, int X, int Y, lambda &&body)
Definition forall.hpp:762
void subtract(const Vector &x, const Vector &y, Vector &z)
Definition vector.cpp:547
float real_t
Definition config.hpp:43
ElementDofOrdering
Constants describing the possible orderings of the DOFs in one element.
Definition fespace.hpp:83
void forall(int N, lambda &&body)
Definition forall.hpp:753
IntegrationRules IntRules(0, Quadrature1D::GaussLegendre)
A global object with all integration rules (defined in intrules.cpp)
Definition intrules.hpp:492
real_t p(const Vector &x, real_t t)
constexpr int GLVIZ_W
constexpr real_t NL_DMAX
constexpr real_t NRM
constexpr int SDIM
constexpr char vishost[]
cdouble LogEllipticTheta1Prime(const cdouble u, const cdouble q)
constexpr int GLVIZ_H
constexpr int DIM
constexpr Element::Type QUAD
constexpr real_t EPS
std::complex< real_t > cdouble
cdouble WeierstrassZeta(const cdouble z, const cdouble w1=0.5, const cdouble w3=real_t(0.5) *I)
cdouble EllipticTheta1Prime(const int k, const cdouble u, const cdouble q)
cdouble WeierstrassP(const cdouble z, const cdouble w1=0.5, const cdouble w3=real_t(0.5) *I)
constexpr real_t PI
cdouble EllipticTheta(const int a, const cdouble u, const cdouble q)
Helper struct to convert a C++ type to an MPI type.
std::array< int, NCMesh::MaxFaceNodes > nodes