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