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