MFEM  v4.6.0
Finite element discretization library
pmesh-fitting.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 // Boundary and Interface Fitting Miniapp
14 // --------------------------------------------------------------
15 //
16 // This miniapp performs mesh optimization for controlling mesh quality and
17 // aligning a selected set of nodes to boundary and/or interface of interest
18 // defined using a level-set function. The mesh quality aspect is based on a
19 // variational formulation of the Target-Matrix Optimization Paradigm (TMOP).
20 // Boundary/interface alignment is weakly enforced using a penalization term
21 // that moves a selected set of nodes towards the zero level set of a signed
22 // smooth discrete function.
23 //
24 // See the following papers for more details:
25 //
26 // [1] "Adaptive Surface Fitting and Tangential Relaxation for High-Order Mesh Optimization" by
27 // Knupp, Kolev, Mittal, Tomov.
28 // [2] "High-Order Mesh Morphing for Boundary and Interface Fitting to Implicit Geometries" by
29 // Barrera, Kolev, Mittal, Tomov.
30 // [3] "The target-matrix optimization paradigm for high-order meshes" by
31 // Dobrev, Knupp, Kolev, Mittal, Tomov.
32 //
33 // Compile with: make pmesh-fitting
34 //
35 // Sample runs:
36 // Interface fitting:
37 // mpirun -np 4 pmesh-fitting -o 3 -mid 58 -tid 1 -ni 200 -vl 1 -sfc 5e4 -rtol 1e-5
38 // mpirun -np 4 pmesh-fitting -m square01-tri.mesh -o 3 -rs 0 -mid 58 -tid 1 -ni 200 -vl 1 -sfc 1e4 -rtol 1e-5
39 // Surface fitting with weight adaptation and termination based on fitting error:
40 // mpirun -np 4 pmesh-fitting -o 2 -mid 2 -tid 1 -ni 100 -vl 2 -sfc 10 -rtol 1e-20 -st 0 -sfa 10.0 -sft 1e-5
41 // Fitting to Fischer-Tropsch reactor like domain (requires GSLIB):
42 // * mpirun -np 6 pmesh-fitting -m ../../data/inline-tri.mesh -o 2 -rs 4 -mid 2 -tid 1 -vl 2 -sfc 100 -rtol 1e-12 -ni 100 -li 40 -ae 1 -bnd -sbgmesh -slstype 2 -smtype 0 -sfa 10.0 -sft 1e-4 -amriter 5 -dist -mod-bndr-attr
43 
44 #include "mesh-fitting.hpp"
45 
46 using namespace mfem;
47 using namespace std;
48 
49 int main (int argc, char *argv[])
50 {
51  // 0. Initialize MPI and HYPRE.
52  Mpi::Init(argc, argv);
53  int myid = Mpi::WorldRank();
54  Hypre::Init();
55 
56  // 1. Set the method's default parameters.
57  const char *mesh_file = "square01.mesh";
58  int mesh_poly_deg = 1;
59  int rs_levels = 1;
60  int rp_levels = 0;
61  int metric_id = 2;
62  int target_id = 1;
63  double surface_fit_const = 100.0;
64  int quad_type = 1;
65  int quad_order = 8;
66  int solver_type = 0;
67  int solver_iter = 20;
68  double solver_rtol = 1e-10;
69  int solver_art_type = 0;
70  int lin_solver = 2;
71  int max_lin_iter = 100;
72  bool move_bnd = true;
73  bool visualization = true;
74  int verbosity_level = 0;
75  int adapt_eval = 0;
76  const char *devopt = "cpu";
77  double surface_fit_adapt = 0.0;
78  double surface_fit_threshold = -10;
79  bool adapt_marking = false;
80  bool surf_bg_mesh = false;
81  bool comp_dist = false;
82  int surf_ls_type = 1;
83  int marking_type = 0;
84  bool mod_bndr_attr = false;
85  bool material = false;
86  int mesh_node_ordering = 0;
87  int amr_iters = 0;
88 
89  // 2. Parse command-line options.
90  OptionsParser args(argc, argv);
91  args.AddOption(&mesh_file, "-m", "--mesh",
92  "Mesh file to use.");
93  args.AddOption(&mesh_poly_deg, "-o", "--order",
94  "Polynomial degree of mesh finite element space.");
95  args.AddOption(&rs_levels, "-rs", "--refine-serial",
96  "Number of times to refine the mesh uniformly in serial.");
97  args.AddOption(&rp_levels, "-rp", "--refine-parallel",
98  "Number of times to refine the mesh uniformly in parallel.");
99  args.AddOption(&metric_id, "-mid", "--metric-id",
100  "Mesh optimization metric. See list in mesh-optimizer.");
101  args.AddOption(&target_id, "-tid", "--target-id",
102  "Target (ideal element) type:\n\t"
103  "1: Ideal shape, unit size\n\t"
104  "2: Ideal shape, equal size\n\t"
105  "3: Ideal shape, initial size\n\t"
106  "4: Given full analytic Jacobian (in physical space)\n\t"
107  "5: Ideal shape, given size (in physical space)");
108  args.AddOption(&surface_fit_const, "-sfc", "--surface-fit-const",
109  "Surface preservation constant.");
110  args.AddOption(&quad_type, "-qt", "--quad-type",
111  "Quadrature rule type:\n\t"
112  "1: Gauss-Lobatto\n\t"
113  "2: Gauss-Legendre\n\t"
114  "3: Closed uniform points");
115  args.AddOption(&quad_order, "-qo", "--quad_order",
116  "Order of the quadrature rule.");
117  args.AddOption(&solver_type, "-st", "--solver-type",
118  " Type of solver: (default) 0: Newton, 1: LBFGS");
119  args.AddOption(&solver_iter, "-ni", "--newton-iters",
120  "Maximum number of Newton iterations.");
121  args.AddOption(&solver_rtol, "-rtol", "--newton-rel-tolerance",
122  "Relative tolerance for the Newton solver.");
123  args.AddOption(&solver_art_type, "-art", "--adaptive-rel-tol",
124  "Type of adaptive relative linear solver tolerance:\n\t"
125  "0: None (default)\n\t"
126  "1: Eisenstat-Walker type 1\n\t"
127  "2: Eisenstat-Walker type 2");
128  args.AddOption(&lin_solver, "-ls", "--lin-solver",
129  "Linear solver:\n\t"
130  "0: l1-Jacobi\n\t"
131  "1: CG\n\t"
132  "2: MINRES\n\t"
133  "3: MINRES + Jacobi preconditioner\n\t"
134  "4: MINRES + l1-Jacobi preconditioner");
135  args.AddOption(&max_lin_iter, "-li", "--lin-iter",
136  "Maximum number of iterations in the linear solve.");
137  args.AddOption(&move_bnd, "-bnd", "--move-boundary", "-fix-bnd",
138  "--fix-boundary",
139  "Enable motion along horizontal and vertical boundaries.");
140  args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
141  "--no-visualization",
142  "Enable or disable GLVis visualization.");
143  args.AddOption(&verbosity_level, "-vl", "--verbosity-level",
144  "Set the verbosity level - 0, 1, or 2.");
145  args.AddOption(&adapt_eval, "-ae", "--adaptivity-evaluator",
146  "0 - Advection based (DEFAULT), 1 - GSLIB.");
147  args.AddOption(&devopt, "-d", "--device",
148  "Device configuration string, see Device::Configure().");
149  args.AddOption(&surface_fit_adapt, "-sfa", "--adaptive-surface-fit",
150  "Enable or disable adaptive surface fitting.");
151  args.AddOption(&surface_fit_threshold, "-sft", "--surf-fit-threshold",
152  "Set threshold for surface fitting. TMOP solver will"
153  "terminate when max surface fitting error is below this limit");
154  args.AddOption(&adapt_marking, "-marking", "--adaptive-marking", "-no-amarking",
155  "--no-adaptive-marking",
156  "Enable or disable adaptive marking surface fitting.");
157  args.AddOption(&surf_bg_mesh, "-sbgmesh", "--surf-bg-mesh",
158  "-no-sbgmesh","--no-surf-bg-mesh",
159  "Use background mesh for surface fitting.");
160  args.AddOption(&comp_dist, "-dist", "--comp-dist",
161  "-no-dist","--no-comp-dist",
162  "Compute distance from 0 level set or not.");
163  args.AddOption(&surf_ls_type, "-slstype", "--surf-ls-type",
164  "1 - Circle (DEFAULT), 2 - Squircle, 3 - Butterfly.");
165  args.AddOption(&marking_type, "-smtype", "--surf-marking-type",
166  "1 - Interface (DEFAULT), 2 - Boundary attribute.");
167  args.AddOption(&mod_bndr_attr, "-mod-bndr-attr", "--modify-boundary-attribute",
168  "-fix-bndr-attr", "--fix-boundary-attribute",
169  "Change boundary attribute based on alignment with Cartesian axes.");
170  args.AddOption(&material, "-mat", "--mat",
171  "-no-mat","--no-mat", "Use default material attributes.");
172  args.AddOption(&mesh_node_ordering, "-mno", "--mesh_node_ordering",
173  "Ordering of mesh nodes."
174  "0 (default): byNodes, 1: byVDIM");
175  args.AddOption(&amr_iters, "-amriter", "--amr-iter",
176  "Number of amr iterations on background mesh");
177  args.Parse();
178  if (!args.Good())
179  {
180  if (myid == 0) { args.PrintUsage(cout); }
181  return 1;
182  }
183  if (myid == 0) { args.PrintOptions(cout); }
184 
185  Device device(devopt);
186  if (myid == 0) { device.Print();}
187 
188  // 3. Initialize and refine the starting mesh.
189  Mesh *mesh = new Mesh(mesh_file, 1, 1, false);
190  for (int lev = 0; lev < rs_levels; lev++)
191  {
192  mesh->UniformRefinement();
193  }
194  const int dim = mesh->Dimension();
195 
196  // Define level-set coefficient
197  FunctionCoefficient *ls_coeff = NULL;
198  if (surf_ls_type == 1) //Circle
199  {
200  ls_coeff = new FunctionCoefficient(circle_level_set);
201  }
202  else if (surf_ls_type == 2) // reactor
203  {
204  ls_coeff = new FunctionCoefficient(reactor);
205  }
206  else if (surf_ls_type == 6) // 3D shape
207  {
208  ls_coeff = new FunctionCoefficient(csg_cubecylsph);
209  }
210  else
211  {
212  MFEM_ABORT("Surface fitting level set type not implemented yet.")
213  }
214 
215  ParMesh *pmesh = new ParMesh(MPI_COMM_WORLD, *mesh);
216  delete mesh;
217  for (int lev = 0; lev < rp_levels; lev++) { pmesh->UniformRefinement(); }
218 
219  // 4. Setup background mesh for surface fitting
220  ParMesh *pmesh_surf_fit_bg = NULL;
221  if (surf_bg_mesh)
222  {
223  Mesh *mesh_surf_fit_bg = NULL;
224  if (dim == 2)
225  {
226  mesh_surf_fit_bg =
228  }
229  else if (dim == 3)
230  {
231  mesh_surf_fit_bg =
232  new Mesh(Mesh::MakeCartesian3D(4, 4, 4, Element::HEXAHEDRON, true));
233  }
234  mesh_surf_fit_bg->EnsureNCMesh();
235  pmesh_surf_fit_bg = new ParMesh(MPI_COMM_WORLD, *mesh_surf_fit_bg);
236  delete mesh_surf_fit_bg;
237  }
238 
239  // 5. Define a finite element space on the mesh. Here we use vector finite
240  // elements which are tensor products of quadratic finite elements. The
241  // number of components in the vector finite element space is specified by
242  // the last parameter of the FiniteElementSpace constructor.
244  if (mesh_poly_deg <= 0)
245  {
246  fec = new QuadraticPosFECollection;
247  mesh_poly_deg = 2;
248  }
249  else { fec = new H1_FECollection(mesh_poly_deg, dim); }
250  ParFiniteElementSpace *pfespace =
251  new ParFiniteElementSpace(pmesh, fec, dim, mesh_node_ordering);
252 
253  // 6. Make the mesh curved based on the above finite element space. This
254  // means that we define the mesh elements through a fespace-based
255  // transformation of the reference element.
256  pmesh->SetNodalFESpace(pfespace);
257 
258  // 7. Get the mesh nodes (vertices and other degrees of freedom in the finite
259  // element space) as a finite element grid function in fespace. Note that
260  // changing x automatically changes the shapes of the mesh elements.
261  ParGridFunction x(pfespace);
262  pmesh->SetNodalGridFunction(&x);
263  x.SetTrueVector();
264 
265  // 10. Save the starting (prior to the optimization) mesh to a file. This
266  // output can be viewed later using GLVis: "glvis -m perturbed -np
267  // num_mpi_tasks".
268  {
269  ostringstream mesh_name;
270  mesh_name << "perturbed.mesh";
271  ofstream mesh_ofs(mesh_name.str().c_str());
272  mesh_ofs.precision(8);
273  pmesh->PrintAsSerial(mesh_ofs);
274  }
275 
276  // 11. Store the starting (prior to the optimization) positions.
277  ParGridFunction x0(pfespace);
278  x0 = x;
279 
280  // 12. Form the integrator that uses the chosen metric and target.
281  TMOP_QualityMetric *metric = NULL;
282  switch (metric_id)
283  {
284  // T-metrics
285  case 2: metric = new TMOP_Metric_002; break;
286  case 58: metric = new TMOP_Metric_058; break;
287  case 80: metric = new TMOP_Metric_080(0.5); break;
288  case 303: metric = new TMOP_Metric_303; break;
289  case 328: metric = new TMOP_Metric_328; break;
290  default:
291  if (myid == 0) { cout << "Unknown metric_id: " << metric_id << endl; }
292  return 3;
293  }
294 
295  if (metric_id < 300)
296  {
297  MFEM_VERIFY(dim == 2, "Incompatible metric for 3D meshes");
298  }
299  if (metric_id >= 300)
300  {
301  MFEM_VERIFY(dim == 3, "Incompatible metric for 2D meshes");
302  }
303 
305  TargetConstructor *target_c = NULL;
306  switch (target_id)
307  {
308  case 1: target_t = TargetConstructor::IDEAL_SHAPE_UNIT_SIZE; break;
309  case 2: target_t = TargetConstructor::IDEAL_SHAPE_EQUAL_SIZE; break;
310  case 3: target_t = TargetConstructor::IDEAL_SHAPE_GIVEN_SIZE; break;
311  case 4: target_t = TargetConstructor::GIVEN_SHAPE_AND_SIZE; break;
312  default:
313  if (myid == 0) { cout << "Unknown target_id: " << target_id << endl; }
314  return 3;
315  }
316 
317  if (target_c == NULL)
318  {
319  target_c = new TargetConstructor(target_t, MPI_COMM_WORLD);
320  }
321  target_c->SetNodes(x0);
322  TMOP_Integrator *tmop_integ = new TMOP_Integrator(metric, target_c);
323 
324  // Setup the quadrature rules for the TMOP integrator.
325  IntegrationRules *irules = NULL;
326  switch (quad_type)
327  {
328  case 1: irules = &IntRulesLo; break;
329  case 2: irules = &IntRules; break;
330  case 3: irules = &IntRulesCU; break;
331  default:
332  if (myid == 0) { cout << "Unknown quad_type: " << quad_type << endl; }
333  return 3;
334  }
335  tmop_integ->SetIntegrationRules(*irules, quad_order);
336  if (myid == 0 && dim == 2)
337  {
338  cout << "Triangle quadrature points: "
339  << irules->Get(Geometry::TRIANGLE, quad_order).GetNPoints()
340  << "\nQuadrilateral quadrature points: "
341  << irules->Get(Geometry::SQUARE, quad_order).GetNPoints() << endl;
342  }
343  if (myid == 0 && dim == 3)
344  {
345  cout << "Tetrahedron quadrature points: "
346  << irules->Get(Geometry::TETRAHEDRON, quad_order).GetNPoints()
347  << "\nHexahedron quadrature points: "
348  << irules->Get(Geometry::CUBE, quad_order).GetNPoints()
349  << "\nPrism quadrature points: "
350  << irules->Get(Geometry::PRISM, quad_order).GetNPoints() << endl;
351  }
352 
353  // Modify boundary attribute for surface node movement
354  // Sets attributes of a boundary element to 1/2/3 if it is parallel to x/y/z.
355  if (mod_bndr_attr)
356  {
358  pmesh->SetAttributes();
359  }
360  pmesh->ExchangeFaceNbrData();
361 
362 
363  // Surface fitting.
364  L2_FECollection mat_coll(0, dim);
365  H1_FECollection surf_fit_fec(mesh_poly_deg, dim);
366  ParFiniteElementSpace surf_fit_fes(pmesh, &surf_fit_fec);
367  ParFiniteElementSpace mat_fes(pmesh, &mat_coll);
368  ParGridFunction mat(&mat_fes);
369  ParGridFunction surf_fit_mat_gf(&surf_fit_fes);
370  ParGridFunction surf_fit_gf0(&surf_fit_fes);
371  Array<bool> surf_fit_marker(surf_fit_gf0.Size());
372  ConstantCoefficient surf_fit_coeff(surface_fit_const);
373  AdaptivityEvaluator *adapt_surface = NULL;
374  AdaptivityEvaluator *adapt_grad_surface = NULL;
375  AdaptivityEvaluator *adapt_hess_surface = NULL;
376 
377  // Background mesh FECollection, FESpace, and GridFunction
378  H1_FECollection *surf_fit_bg_fec = NULL;
379  ParFiniteElementSpace *surf_fit_bg_fes = NULL;
380  ParGridFunction *surf_fit_bg_gf0 = NULL;
381  ParFiniteElementSpace *surf_fit_bg_grad_fes = NULL;
382  ParGridFunction *surf_fit_bg_grad = NULL;
383  ParFiniteElementSpace *surf_fit_bg_hess_fes = NULL;
384  ParGridFunction *surf_fit_bg_hess = NULL;
385 
386  // If a background mesh is used, we interpolate the Gradient and Hessian
387  // from that mesh to the current mesh being optimized.
388  ParFiniteElementSpace *surf_fit_grad_fes = NULL;
389  ParGridFunction *surf_fit_grad = NULL;
390  ParFiniteElementSpace *surf_fit_hess_fes = NULL;
391  ParGridFunction *surf_fit_hess = NULL;
392 
393  if (surf_bg_mesh)
394  {
395  pmesh_surf_fit_bg->SetCurvature(mesh_poly_deg);
396 
397  Vector p_min(dim), p_max(dim);
398  pmesh->GetBoundingBox(p_min, p_max);
399  GridFunction &x_bg = *pmesh_surf_fit_bg->GetNodes();
400  const int num_nodes = x_bg.Size() / dim;
401  for (int i = 0; i < num_nodes; i++)
402  {
403  for (int d = 0; d < dim; d++)
404  {
405  double length_d = p_max(d) - p_min(d),
406  extra_d = 0.2 * length_d;
407  x_bg(i + d*num_nodes) = p_min(d) - extra_d +
408  x_bg(i + d*num_nodes) * (length_d + 2*extra_d);
409  }
410  }
411  surf_fit_bg_fec = new H1_FECollection(mesh_poly_deg+1, dim);
412  surf_fit_bg_fes = new ParFiniteElementSpace(pmesh_surf_fit_bg, surf_fit_bg_fec);
413  surf_fit_bg_gf0 = new ParGridFunction(surf_fit_bg_fes);
414  }
415 
416  Array<int> vdofs;
417  if (surface_fit_const > 0.0)
418  {
419  surf_fit_gf0.ProjectCoefficient(*ls_coeff);
420  if (surf_bg_mesh)
421  {
422  OptimizeMeshWithAMRAroundZeroLevelSet(*pmesh_surf_fit_bg, *ls_coeff,
423  amr_iters, *surf_fit_bg_gf0);
424  pmesh_surf_fit_bg->Rebalance();
425  surf_fit_bg_fes->Update();
426  surf_fit_bg_gf0->Update();
427 
428  if (comp_dist)
429  {
430  ComputeScalarDistanceFromLevelSet(*pmesh_surf_fit_bg, *ls_coeff,
431  *surf_fit_bg_gf0);
432  }
433  else { surf_fit_bg_gf0->ProjectCoefficient(*ls_coeff); }
434 
435 
436  surf_fit_bg_grad_fes =
437  new ParFiniteElementSpace(pmesh_surf_fit_bg, surf_fit_bg_fec, dim);
438  surf_fit_bg_grad = new ParGridFunction(surf_fit_bg_grad_fes);
439 
440  surf_fit_grad_fes =
441  new ParFiniteElementSpace(pmesh, &surf_fit_fec, dim);
442  surf_fit_grad = new ParGridFunction(surf_fit_grad_fes);
443 
444  surf_fit_bg_hess_fes =
445  new ParFiniteElementSpace(pmesh_surf_fit_bg, surf_fit_bg_fec, dim * dim);
446  surf_fit_bg_hess = new ParGridFunction(surf_fit_bg_hess_fes);
447 
448  surf_fit_hess_fes =
449  new ParFiniteElementSpace(pmesh, &surf_fit_fec, dim * dim);
450  surf_fit_hess = new ParGridFunction(surf_fit_hess_fes);
451 
452  //Setup gradient of the background mesh
453  const int size_bg = surf_fit_bg_gf0->Size();
454  for (int d = 0; d < pmesh_surf_fit_bg->Dimension(); d++)
455  {
456  ParGridFunction surf_fit_bg_grad_comp(
457  surf_fit_bg_fes, surf_fit_bg_grad->GetData() + d * size_bg);
458  surf_fit_bg_gf0->GetDerivative(1, d, surf_fit_bg_grad_comp);
459  }
460 
461  //Setup Hessian on background mesh
462  int id = 0;
463  for (int d = 0; d < pmesh_surf_fit_bg->Dimension(); d++)
464  {
465  for (int idir = 0; idir < pmesh_surf_fit_bg->Dimension(); idir++)
466  {
467  ParGridFunction surf_fit_bg_grad_comp(
468  surf_fit_bg_fes, surf_fit_bg_grad->GetData() + d * size_bg);
469  ParGridFunction surf_fit_bg_hess_comp(
470  surf_fit_bg_fes, surf_fit_bg_hess->GetData()+ id * size_bg);
471  surf_fit_bg_grad_comp.GetDerivative(1, idir,
472  surf_fit_bg_hess_comp);
473  id++;
474  }
475  }
476  }
477  else // !surf_bg_mesh
478  {
479  if (comp_dist)
480  {
481  ComputeScalarDistanceFromLevelSet(*pmesh, *ls_coeff, surf_fit_gf0);
482  }
483  }
484 
485  // Set material gridfunction
486  for (int i = 0; i < pmesh->GetNE(); i++)
487  {
488  if (material)
489  {
490  mat(i) = pmesh->GetAttribute(i)-1;
491  }
492  else
493  {
494  mat(i) = material_id(i, surf_fit_gf0);
495  pmesh->SetAttribute(i, mat(i) + 1);
496  }
497  }
498 
499  // Adapt attributes for marking such that if all but 1 face of an element
500  // are marked, the element attribute is switched.
501  if (adapt_marking)
502  {
503  ModifyAttributeForMarkingDOFS(pmesh, mat, 0);
504  ModifyAttributeForMarkingDOFS(pmesh, mat, 1);
505  }
506 
507  GridFunctionCoefficient coeff_mat(&mat);
508  surf_fit_mat_gf.ProjectDiscCoefficient(coeff_mat,
510  surf_fit_mat_gf.SetTrueVector();
511  surf_fit_mat_gf.SetFromTrueVector();
512 
513  // Set DOFs for fitting
514  // Strategy 1: Choose face between elements of different attributes.
515  if (marking_type == 0)
516  {
517  mat.ExchangeFaceNbrData();
518  const Vector &FaceNbrData = mat.FaceNbrData();
519  for (int j = 0; j < surf_fit_marker.Size(); j++)
520  {
521  surf_fit_marker[j] = false;
522  }
523  surf_fit_mat_gf = 0.0;
524 
525  Array<int> dof_list;
526  Array<int> dofs;
527  for (int i = 0; i < pmesh->GetNumFaces(); i++)
528  {
529  auto tr = pmesh->GetInteriorFaceTransformations(i);
530  if (tr != NULL)
531  {
532  int mat1 = mat(tr->Elem1No);
533  int mat2 = mat(tr->Elem2No);
534  if (mat1 != mat2)
535  {
536  surf_fit_gf0.ParFESpace()->GetFaceDofs(i, dofs);
537  dof_list.Append(dofs);
538  }
539  }
540  }
541  for (int i = 0; i < pmesh->GetNSharedFaces(); i++)
542  {
543  auto tr = pmesh->GetSharedFaceTransformations(i);
544  if (tr != NULL)
545  {
546  int faceno = pmesh->GetSharedFace(i);
547  int mat1 = mat(tr->Elem1No);
548  int mat2 = FaceNbrData(tr->Elem2No-pmesh->GetNE());
549  if (mat1 != mat2)
550  {
551  surf_fit_gf0.ParFESpace()->GetFaceDofs(faceno, dofs);
552  dof_list.Append(dofs);
553  }
554  }
555  }
556  for (int i = 0; i < dof_list.Size(); i++)
557  {
558  surf_fit_marker[dof_list[i]] = true;
559  surf_fit_mat_gf(dof_list[i]) = 1.0;
560  }
561  }
562  // Strategy 2: Mark all boundaries with attribute marking_type
563  else if (marking_type > 0)
564  {
565  for (int i = 0; i < pmesh->GetNBE(); i++)
566  {
567  const int attr = pmesh->GetBdrElement(i)->GetAttribute();
568  if (attr == marking_type)
569  {
570  surf_fit_fes.GetBdrElementVDofs(i, vdofs);
571  for (int j = 0; j < vdofs.Size(); j++)
572  {
573  surf_fit_marker[vdofs[j]] = true;
574  surf_fit_mat_gf(vdofs[j]) = 1.0;
575  }
576  }
577  }
578  }
579 
580  // Set AdaptivityEvaluators for transferring information from initial
581  // mesh to current mesh as it moves during adaptivity.
582  if (adapt_eval == 0)
583  {
584  adapt_surface = new AdvectorCG;
585  MFEM_VERIFY(!surf_bg_mesh, "Background meshes require GSLIB.");
586  }
587  else if (adapt_eval == 1)
588  {
589 #ifdef MFEM_USE_GSLIB
590  adapt_surface = new InterpolatorFP;
591  if (surf_bg_mesh)
592  {
593  adapt_grad_surface = new InterpolatorFP;
594  adapt_hess_surface = new InterpolatorFP;
595  }
596 #else
597  MFEM_ABORT("MFEM is not built with GSLIB support!");
598 #endif
599  }
600  else { MFEM_ABORT("Bad interpolation option."); }
601 
602  if (!surf_bg_mesh)
603  {
604  tmop_integ->EnableSurfaceFitting(surf_fit_gf0, surf_fit_marker,
605  surf_fit_coeff, *adapt_surface);
606  }
607  else
608  {
609  tmop_integ->EnableSurfaceFittingFromSource(
610  *surf_fit_bg_gf0, surf_fit_gf0,
611  surf_fit_marker, surf_fit_coeff, *adapt_surface,
612  *surf_fit_bg_grad, *surf_fit_grad, *adapt_grad_surface,
613  *surf_fit_bg_hess, *surf_fit_hess, *adapt_hess_surface);
614  }
615 
616  if (visualization)
617  {
618  socketstream vis1, vis2, vis3, vis4, vis5;
619  common::VisualizeField(vis1, "localhost", 19916, surf_fit_gf0,
620  "Level Set", 0, 0, 300, 300);
621  common::VisualizeField(vis2, "localhost", 19916, mat,
622  "Materials", 300, 0, 300, 300);
623  common::VisualizeField(vis3, "localhost", 19916, surf_fit_mat_gf,
624  "Surface DOFs", 600, 0, 300, 300);
625  if (surf_bg_mesh)
626  {
627  common::VisualizeField(vis4, "localhost", 19916, *surf_fit_bg_gf0,
628  "Level Set - Background",
629  0, 400, 300, 300);
630  }
631  }
632  }
633 
634  // 13. Setup the final NonlinearForm (which defines the integral of interest,
635  // its first and second derivatives). Here we can use a combination of
636  // metrics, i.e., optimize the sum of two integrals, where both are
637  // scaled by used-defined space-dependent weights. Note that there are
638  // no command-line options for the weights and the type of the second
639  // metric; one should update those in the code.
640  ParNonlinearForm a(pfespace);
641  ConstantCoefficient *metric_coeff1 = NULL;
642  a.AddDomainIntegrator(tmop_integ);
643 
644  // Compute the minimum det(J) of the starting mesh.
645  double min_detJ = infinity();
646  const int NE = pmesh->GetNE();
647  for (int i = 0; i < NE; i++)
648  {
649  const IntegrationRule &ir =
650  irules->Get(pfespace->GetFE(i)->GetGeomType(), quad_order);
652  for (int j = 0; j < ir.GetNPoints(); j++)
653  {
654  transf->SetIntPoint(&ir.IntPoint(j));
655  min_detJ = min(min_detJ, transf->Jacobian().Det());
656  }
657  }
658  MPI_Allreduce(MPI_IN_PLACE, &min_detJ, 1,
659  MPI_DOUBLE, MPI_MIN, MPI_COMM_WORLD);
660  if (myid == 0)
661  { cout << "Minimum det(J) of the original mesh is " << min_detJ << endl; }
662 
663  MFEM_VERIFY(min_detJ > 0, "The input mesh is inverted, use mesh-optimizer.");
664 
665  const double init_energy = a.GetParGridFunctionEnergy(x);
666  double init_metric_energy = init_energy;
667  if (surface_fit_const > 0.0)
668  {
669  surf_fit_coeff.constant = 0.0;
670  init_metric_energy = a.GetParGridFunctionEnergy(x);
671  surf_fit_coeff.constant = surface_fit_const;
672  }
673 
674  // 14. Fix all boundary nodes, or fix only a given component depending on the
675  // boundary attributes of the given mesh. Attributes 1/2/3 correspond to
676  // fixed x/y/z components of the node. Attribute dim+1 corresponds to
677  // an entirely fixed node.
678  if (move_bnd == false)
679  {
680  Array<int> ess_bdr(pmesh->bdr_attributes.Max());
681  ess_bdr = 1;
682  if (marking_type > 0)
683  {
684  ess_bdr[marking_type-1] = 0;
685  }
686  a.SetEssentialBC(ess_bdr);
687  }
688  else
689  {
690  int n = 0;
691  for (int i = 0; i < pmesh->GetNBE(); i++)
692  {
693  const int nd = pfespace->GetBE(i)->GetDof();
694  const int attr = pmesh->GetBdrElement(i)->GetAttribute();
695  MFEM_VERIFY(!(dim == 2 && attr == 3),
696  "Boundary attribute 3 must be used only for 3D meshes. "
697  "Adjust the attributes (1/2/3/4 for fixed x/y/z/all "
698  "components, rest for free nodes), or use -fix-bnd.");
699  if (attr == 1 || attr == 2 || attr == 3) { n += nd; }
700  if (attr == 4) { n += nd * dim; }
701  }
702  Array<int> ess_vdofs(n);
703  n = 0;
704  for (int i = 0; i < pmesh->GetNBE(); i++)
705  {
706  const int nd = pfespace->GetBE(i)->GetDof();
707  const int attr = pmesh->GetBdrElement(i)->GetAttribute();
708  pfespace->GetBdrElementVDofs(i, vdofs);
709  if (attr == 1) // Fix x components.
710  {
711  for (int j = 0; j < nd; j++)
712  { ess_vdofs[n++] = vdofs[j]; }
713  }
714  else if (attr == 2) // Fix y components.
715  {
716  for (int j = 0; j < nd; j++)
717  { ess_vdofs[n++] = vdofs[j+nd]; }
718  }
719  else if (attr == 3) // Fix z components.
720  {
721  for (int j = 0; j < nd; j++)
722  { ess_vdofs[n++] = vdofs[j+2*nd]; }
723  }
724  else if (attr == 4) // Fix all components.
725  {
726  for (int j = 0; j < vdofs.Size(); j++)
727  { ess_vdofs[n++] = vdofs[j]; }
728  }
729  }
730  a.SetEssentialVDofs(ess_vdofs);
731  }
732 
733  // 15. As we use the Newton method to solve the resulting nonlinear system,
734  // here we setup the linear solver for the system's Jacobian.
735  Solver *S = NULL, *S_prec = NULL;
736  const double linsol_rtol = 1e-12;
737  if (lin_solver == 0)
738  {
739  S = new DSmoother(1, 1.0, max_lin_iter);
740  }
741  else if (lin_solver == 1)
742  {
743  CGSolver *cg = new CGSolver(MPI_COMM_WORLD);
744  cg->SetMaxIter(max_lin_iter);
745  cg->SetRelTol(linsol_rtol);
746  cg->SetAbsTol(0.0);
747  cg->SetPrintLevel(verbosity_level >= 2 ? 3 : -1);
748  S = cg;
749  }
750  else
751  {
752  MINRESSolver *minres = new MINRESSolver(MPI_COMM_WORLD);
753  minres->SetMaxIter(max_lin_iter);
754  minres->SetRelTol(linsol_rtol);
755  minres->SetAbsTol(0.0);
756  if (verbosity_level > 2) { minres->SetPrintLevel(1); }
757  else { minres->SetPrintLevel(verbosity_level == 2 ? 3 : -1); }
758  if (lin_solver == 3 || lin_solver == 4)
759  {
760  auto hs = new HypreSmoother;
761  hs->SetType((lin_solver == 3) ? HypreSmoother::Jacobi
762  /* */ : HypreSmoother::l1Jacobi, 1);
763  hs->SetPositiveDiagonal(true);
764  S_prec = hs;
765  minres->SetPreconditioner(*S_prec);
766  }
767  S = minres;
768  }
769 
770  // Perform the nonlinear optimization.
771  const IntegrationRule &ir =
772  irules->Get(pfespace->GetFE(0)->GetGeomType(), quad_order);
773  TMOPNewtonSolver solver(pfespace->GetComm(), ir, solver_type);
774  if (surface_fit_adapt > 0.0)
775  {
776  solver.SetAdaptiveSurfaceFittingScalingFactor(surface_fit_adapt);
777  }
778  if (surface_fit_threshold > 0)
779  {
780  solver.SetTerminationWithMaxSurfaceFittingError(surface_fit_threshold);
781  }
782  // Provide all integration rules in case of a mixed mesh.
783  solver.SetIntegrationRules(*irules, quad_order);
784  if (solver_type == 0)
785  {
786  // Specify linear solver when we use a Newton-based solver.
787  solver.SetPreconditioner(*S);
788  }
789  solver.SetMaxIter(solver_iter);
790  solver.SetRelTol(solver_rtol);
791  solver.SetAbsTol(0.0);
792  solver.SetMinimumDeterminantThreshold(0.001*min_detJ);
793  if (solver_art_type > 0)
794  {
795  solver.SetAdaptiveLinRtol(solver_art_type, 0.5, 0.9);
796  }
797  solver.SetPrintLevel(verbosity_level >= 1 ? 1 : -1);
798  solver.SetOperator(a);
799  Vector b(0);
800  solver.Mult(b, x.GetTrueVector());
801  x.SetFromTrueVector();
802 
803  // 16. Save the optimized mesh to a file. This output can be viewed later
804  // using GLVis: "glvis -m optimized -np num_mpi_tasks".
805  {
806  ostringstream mesh_name;
807  mesh_name << "optimized.mesh";
808  ofstream mesh_ofs(mesh_name.str().c_str());
809  mesh_ofs.precision(8);
810  pmesh->PrintAsSerial(mesh_ofs);
811  }
812 
813  // Compute the final energy of the functional.
814  const double fin_energy = a.GetParGridFunctionEnergy(x);
815  double fin_metric_energy = fin_energy;
816  if (surface_fit_const > 0.0)
817  {
818  surf_fit_coeff.constant = 0.0;
819  fin_metric_energy = a.GetParGridFunctionEnergy(x);
820  surf_fit_coeff.constant = surface_fit_const;
821  }
822 
823  if (myid == 0)
824  {
825  std::cout << std::scientific << std::setprecision(4);
826  cout << "Initial strain energy: " << init_energy
827  << " = metrics: " << init_metric_energy
828  << " + extra terms: " << init_energy - init_metric_energy << endl;
829  cout << " Final strain energy: " << fin_energy
830  << " = metrics: " << fin_metric_energy
831  << " + extra terms: " << fin_energy - fin_metric_energy << endl;
832  cout << "The strain energy decreased by: "
833  << (init_energy - fin_energy) * 100.0 / init_energy << " %." << endl;
834  }
835 
836  if (surface_fit_const > 0.0)
837  {
838  if (visualization)
839  {
840  socketstream vis2, vis3;
841  common::VisualizeField(vis2, "localhost", 19916, mat,
842  "Materials", 300, 400, 300, 300);
843  common::VisualizeField(vis3, "localhost", 19916, surf_fit_mat_gf,
844  "Surface DOFs", 600, 400, 300, 300);
845  }
846  double err_avg, err_max;
847  tmop_integ->GetSurfaceFittingErrors(x, err_avg, err_max);
848  if (myid == 0)
849  {
850  std::cout << "Avg fitting error: " << err_avg << std::endl
851  << "Max fitting error: " << err_max << std::endl;
852  }
853  }
854 
855  // 18. Visualize the mesh displacement.
856  if (visualization)
857  {
858  x0 -= x;
859  socketstream vis;
860  common::VisualizeField(vis, "localhost", 19916, x0,
861  "Displacements", 900, 400, 300, 300, "jRmclA");
862  }
863 
864  delete S;
865  delete S_prec;
866  delete metric_coeff1;
867  delete adapt_surface;
868  delete adapt_grad_surface;
869  delete adapt_hess_surface;
870  delete ls_coeff;
871  delete surf_fit_hess;
872  delete surf_fit_hess_fes;
873  delete surf_fit_bg_hess;
874  delete surf_fit_bg_hess_fes;
875  delete surf_fit_grad;
876  delete surf_fit_grad_fes;
877  delete surf_fit_bg_grad;
878  delete surf_fit_bg_grad_fes;
879  delete surf_fit_bg_gf0;
880  delete surf_fit_bg_fes;
881  delete surf_fit_bg_fec;
882  delete target_c;
883  delete metric;
884  delete pfespace;
885  delete fec;
886  delete pmesh_surf_fit_bg;
887  delete pmesh;
888 
889  return 0;
890 }
static void Init()
Initialize hypre by calling HYPRE_Init() and set default options. After calling Hypre::Init(), hypre will be finalized automatically at program exit.
Definition: hypre.hpp:80
Conjugate gradient method.
Definition: solvers.hpp:493
int GetNPoints() const
Returns the number of the points in the integration rule.
Definition: intrules.hpp:253
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
void GetDerivative(int comp, int der_comp, ParGridFunction &der)
Parallel version of GridFunction::GetDerivative(); see its documentation.
Definition: pgridfunc.cpp:521
Data type for scaled Jacobi-type smoother of sparse matrix.
const IntegrationRule & Get(int GeomType, int Order)
Returns an integration rule for given GeomType and Order.
Definition: intrules.cpp:980
IntegrationRules IntRulesLo(0, Quadrature1D::GaussLobatto)
void SetFromTrueVector()
Shortcut for calling SetFromTrueDofs() with GetTrueVector() as argument.
Definition: gridfunc.hpp:150
A coefficient that is constant across space and time.
Definition: coefficient.hpp:84
void PrintOptions(std::ostream &out) const
Print the options.
Definition: optparser.cpp:331
int Dimension() const
Dimension of the reference space used within the elements.
Definition: mesh.hpp:1020
int GetNumFaces() const
Return the number of faces (3D), edges (2D) or vertices (1D).
Definition: mesh.cpp:5668
3D barrier Shape (S) metric, well-posed (polyconvex & invex).
Definition: tmop.hpp:703
void PrintUsage(std::ostream &out) const
Print the usage message.
Definition: optparser.cpp:462
void ModifyAttributeForMarkingDOFS(ParMesh *pmesh, ParGridFunction &mat, int attr_to_switch)
static Mesh MakeCartesian3D(int nx, int ny, int nz, Element::Type type, double sx=1.0, double sy=1.0, double sz=1.0, bool sfc_ordering=true)
Definition: mesh.cpp:3783
virtual void Update(bool want_transform=true)
Definition: pfespace.cpp:3323
void SetIntPoint(const IntegrationPoint *ip)
Set the integration point ip that weights and Jacobians will be evaluated at.
Definition: eltrans.hpp:93
int GetAttribute() const
Return element&#39;s attribute.
Definition: element.hpp:55
void ExchangeFaceNbrData(Table *gr_sface, int *s2l_face)
Definition: pmesh.cpp:2207
double Det() const
Definition: densemat.cpp:487
int Size() const
Returns the size of the vector.
Definition: vector.hpp:197
Coefficient defined by a GridFunction. This coefficient is mesh dependent.
T Max() const
Find the maximal element in the array, using the comparison operator < for class T.
Definition: array.cpp:68
Container class for integration rules.
Definition: intrules.hpp:412
Parallel non-linear operator on the true dofs.
void SetIntegrationRules(IntegrationRules &irules, int order)
Prescribe a set of integration rules; relevant for mixed meshes.
Definition: tmop.hpp:2004
void Print(std::ostream &out=mfem::out)
Print the configuration of the MFEM virtual device object.
Definition: device.cpp:279
virtual void ProjectDiscCoefficient(VectorCoefficient &coeff)
Project a discontinuous vector coefficient as a grid function on a continuous finite element space...
Definition: pgridfunc.cpp:582
bool Good() const
Return true if the command line options were parsed successfully.
Definition: optparser.hpp:159
Abstract parallel finite element space.
Definition: pfespace.hpp:28
virtual void ProjectCoefficient(Coefficient &coeff)
Project coeff Coefficient to this GridFunction. The projection computation depends on the choice of t...
Definition: pgridfunc.cpp:561
STL namespace.
MINRES method.
Definition: solvers.hpp:603
FaceElementTransformations * GetInteriorFaceTransformations(int FaceNo)
Definition: mesh.hpp:1517
int GetNBE() const
Returns number of boundary elements.
Definition: mesh.hpp:1089
IntegrationRules IntRules(0, Quadrature1D::GaussLegendre)
A global object with all integration rules (defined in intrules.cpp)
Definition: intrules.hpp:480
2D barrier shape (S) metric (not polyconvex).
Definition: tmop.hpp:466
Geometry::Type GetGeomType() const
Returns the Geometry::Type of the reference element.
Definition: fe_base.hpp:320
int GetAttribute(int i) const
Return the attribute of element i.
Definition: mesh.hpp:1190
Abstract class for local mesh quality metrics in the target-matrix optimization paradigm (TMOP) by P...
Definition: tmop.hpp:23
IntegrationPoint & IntPoint(int i)
Returns a reference to the i-th integration point.
Definition: intrules.hpp:256
int material(Vector &x, Vector &xmin, Vector &xmax)
Definition: shaper.cpp:53
void Rebalance()
Definition: pmesh.cpp:4044
void SetTrueVector()
Shortcut for calling GetTrueDofs() with GetTrueVector() as argument.
Definition: gridfunc.hpp:144
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
Definition: solvers.hpp:616
virtual void SetPrintLevel(int print_lvl)
Legacy method to set the level of verbosity of the solver output.
Definition: solvers.cpp:71
void GetBoundingBox(Vector &p_min, Vector &p_max, int ref=2)
Definition: pmesh.cpp:6145
double reactor(const Vector &x)
void Parse()
Parse the command-line options. Note that this function expects all the options provided through the ...
Definition: optparser.cpp:151
int Append(const T &el)
Append element &#39;el&#39; to array, resize if necessary.
Definition: array.hpp:759
void EnsureNCMesh(bool simplices_nonconforming=false)
Definition: mesh.cpp:9888
virtual void Update()
Transform by the Space UpdateMatrix (e.g., on Mesh change).
Definition: pgridfunc.cpp:90
double b
Definition: lissajous.cpp:42
void OptimizeMeshWithAMRAroundZeroLevelSet(ParMesh &pmesh, FunctionCoefficient &ls_coeff, int amr_iter, ParGridFunction &distance_s, const int quad_order=5, Array< ParGridFunction *> *pgf_to_update=NULL)
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:10232
const DenseMatrix & Jacobian()
Return the Jacobian matrix of the transformation at the currently set IntegrationPoint, using the method SetIntPoint().
Definition: eltrans.hpp:119
void SetNodes(const GridFunction &n)
Set the nodes to be used in the target-matrix construction.
Definition: tmop.hpp:1409
FaceElementTransformations * GetSharedFaceTransformations(int sf, bool fill2=true)
Definition: pmesh.cpp:3080
void SetMaxIter(int max_it)
Definition: solvers.hpp:201
double csg_cubecylsph(const Vector &x)
Version of QuadraticFECollection with positive basis functions.
Definition: fe_coll.hpp:764
Parallel smoothers in hypre.
Definition: hypre.hpp:971
static void Init()
Singleton creation with Mpi::Init();.
int GetSharedFace(int sface) const
Return the local face index for the given shared face.
Definition: pmesh.cpp:3215
double * GetData() const
Return a pointer to the beginning of the Vector data.
Definition: vector.hpp:206
Vector & FaceNbrData()
Definition: pgridfunc.hpp:208
3D compound barrier Shape+Size (VS) metric (polyconvex, balanced).
Definition: tmop.hpp:906
virtual int GetFaceDofs(int i, Array< int > &dofs, int variant=0) const
Definition: pfespace.cpp:518
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
int GetNSharedFaces() const
Return the number of shared faces (3D), edges (2D), vertices (1D)
Definition: pmesh.cpp:3196
virtual const FiniteElement * GetFE(int i) const
Definition: pfespace.cpp:534
Collection of finite elements from the same family in multiple dimensions. This class is used to matc...
Definition: fe_coll.hpp:26
ParFiniteElementSpace * ParFESpace() const
Definition: pgridfunc.hpp:113
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
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
void SetAttributes() override
Determine the sets of unique attribute values in domain and boundary elements.
Definition: pmesh.cpp:1580
const Vector & GetTrueVector() const
Read only access to the (optional) internal true-dof Vector.
Definition: gridfunc.hpp:130
int GetNE() const
Returns number of elements.
Definition: mesh.hpp:1086
double a
Definition: lissajous.cpp:41
void ModifyBoundaryAttributesForNodeMovement(ParMesh *pmesh, ParGridFunction &x)
double circle_level_set(const Vector &x)
void GetElementTransformation(int i, IsoparametricTransformation *ElTr)
Definition: mesh.cpp:355
DofTransformation * GetBdrElementVDofs(int i, Array< int > &vdofs) const
Returns indices of degrees of freedom for i&#39;th boundary element. The returned indices are offsets int...
Definition: fespace.cpp:297
void EnableSurfaceFittingFromSource(const ParGridFunction &s_bg, ParGridFunction &s0, const Array< bool > &smarker, Coefficient &coeff, AdaptivityEvaluator &ae, const ParGridFunction &s_bg_grad, ParGridFunction &s0_grad, AdaptivityEvaluator &age, const ParGridFunction &s_bg_hess, ParGridFunction &s0_hess, AdaptivityEvaluator &ahe)
Fitting of certain DOFs in the current mesh to the zero level set of a function defined on another (f...
Definition: tmop.cpp:2978
int dim
Definition: ex24.cpp:53
void PrintAsSerial(std::ostream &out=mfem::out) const
Definition: pmesh.cpp:5278
void SetCurvature(int order, bool discont=false, int space_dim=-1, int ordering=1) override
Set the curvature of the mesh nodes using the given polynomial degree.
Definition: pmesh.cpp:2056
int main(int argc, char *argv[])
IntegrationRules IntRulesCU(0, Quadrature1D::ClosedUniform)
void GetSurfaceFittingErrors(const Vector &pos, double &err_avg, double &err_max)
Definition: tmop.cpp:3047
double infinity()
Define a shortcut for std::numeric_limits<double>::infinity()
Definition: vector.hpp:44
static int WorldRank()
Return the MPI rank in MPI_COMM_WORLD.
int Size() const
Return the logical size of the array.
Definition: array.hpp:141
MPI_Comm GetComm() const
Definition: pfespace.hpp:269
A general function coefficient.
void SetAttribute(int i, int attr)
Set the attribute of element i.
Definition: mesh.hpp:1193
Vector data type.
Definition: vector.hpp:58
void VisualizeField(socketstream &sock, const char *vishost, int visport, GridFunction &gf, const char *title, int x, int y, int w, int h, const char *keys, bool vec)
Definition: fem_extras.cpp:92
void SetNodalFESpace(FiniteElementSpace *nfes) override
Definition: pmesh.cpp:2076
Arbitrary order H1-conforming (continuous) finite elements.
Definition: fe_coll.hpp:259
TargetType
Target-matrix construction algorithms supported by this class.
Definition: tmop.hpp:1333
void GetNodes(Vector &node_coord) const
Definition: mesh.cpp:8302
void EnableSurfaceFitting(const GridFunction &s0, const Array< bool > &smarker, Coefficient &coeff, AdaptivityEvaluator &ae)
Fitting of certain DOFs to the zero level set of a function.
Definition: tmop.cpp:2911
void ComputeScalarDistanceFromLevelSet(ParMesh &pmesh, FunctionCoefficient &ls_coeff, ParGridFunction &distance_s, const int nDiffuse=2, const int pLapOrder=5, const int pLapNewton=50)
const Element * GetBdrElement(int i) const
Return pointer to the i&#39;th boundary element object.
Definition: mesh.hpp:1158
Base class for solvers.
Definition: operator.hpp:682
Class for parallel grid function.
Definition: pgridfunc.hpp:32
void SetNodalGridFunction(GridFunction *nodes, bool make_owner=false)
Definition: mesh.cpp:5624
Base class representing target-matrix construction algorithms for mesh optimization via the target-ma...
Definition: tmop.hpp:1329
The MFEM Device class abstracts hardware devices such as GPUs, as well as programming models such as ...
Definition: device.hpp:121
void SetAdaptiveSurfaceFittingScalingFactor(double factor)
Definition: tmop_tools.hpp:240
int material_id(int el_id, const GridFunction &g)
Class for parallel meshes.
Definition: pmesh.hpp:32
const FiniteElement * GetBE(int i) const
Returns pointer to the FiniteElement in the FiniteElementCollection associated with i&#39;th boundary fac...
Definition: fespace.cpp:3166
void SetType(HypreSmoother::Type type, int relax_times=1)
Set the relaxation type and number of sweeps.
Definition: hypre.cpp:3429
Arbitrary order "L2-conforming" discontinuous finite elements.
Definition: fe_coll.hpp:327
A TMOP integrator class based on any given TMOP_QualityMetric and TargetConstructor.
Definition: tmop.hpp:1733