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