MFEM  v4.4.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
pmesh-optimizer.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 // Mesh Optimizer Miniapp: Optimize high-order meshes - Parallel Version
14 // ---------------------------------------------------------------------
15 //
16 // This miniapp performs mesh optimization using the Target-Matrix Optimization
17 // Paradigm (TMOP) by P.Knupp et al., and a global variational minimization
18 // approach. It minimizes the quantity sum_T int_T mu(J(x)), where T are the
19 // target (ideal) elements, J is the Jacobian of the transformation from the
20 // target to the physical element, and mu is the mesh quality metric. This
21 // metric can measure shape, size or alignment of the region around each
22 // quadrature point. The combination of targets & quality metrics is used to
23 // optimize the physical node positions, i.e., they must be as close as possible
24 // to the shape / size / alignment of their targets. This code also demonstrates
25 // a possible use of nonlinear operators (the class TMOP_QualityMetric, defining
26 // mu(J), and the class TMOP_Integrator, defining int mu(J)), as well as their
27 // coupling to Newton methods for solving minimization problems. Note that the
28 // utilized Newton methods are oriented towards avoiding invalid meshes with
29 // negative Jacobian determinants. Each Newton step requires the inversion of a
30 // Jacobian matrix, which is done through an inner linear solver.
31 //
32 // Compile with: make pmesh-optimizer
33 //
34 // Sample runs:
35 // Adapted analytic shape:
36 // mpirun -np 4 pmesh-optimizer -m square01.mesh -o 2 -rs 2 -mid 2 -tid 4 -ni 200 -bnd -qt 1 -qo 8
37 // Adapted analytic size+orientation:
38 // mpirun -np 4 pmesh-optimizer -m square01.mesh -o 2 -rs 2 -mid 14 -tid 4 -ni 200 -bnd -qt 1 -qo 8 -fd
39 // Adapted analytic shape+orientation:
40 // mpirun -np 4 pmesh-optimizer -m square01.mesh -o 2 -rs 2 -mid 85 -tid 4 -ni 100 -bnd -qt 1 -qo 8 -fd
41 //
42 // Adapted analytic shape and/or size with hr-adaptivity:
43 // mpirun -np 4 pmesh-optimizer -m square01.mesh -o 2 -tid 9 -ni 50 -li 20 -hmid 55 -mid 7 -hr
44 // mpirun -np 4 pmesh-optimizer -m square01.mesh -o 2 -tid 10 -ni 50 -li 20 -hmid 55 -mid 7 -hr
45 // mpirun -np 4 pmesh-optimizer -m square01.mesh -o 2 -tid 11 -ni 50 -li 20 -hmid 58 -mid 7 -hr
46 //
47 // Adapted discrete size:
48 // mpirun -np 4 pmesh-optimizer -m square01.mesh -o 2 -rs 2 -mid 80 -tid 5 -ni 50 -qo 4 -nor
49 // Adapted discrete size 3D with PA:
50 // mpirun -np 4 pmesh-optimizer -m cube.mesh -o 2 -rs 2 -mid 321 -tid 5 -ls 3 -nor -pa
51 // Adapted discrete size 3D with PA on device (requires CUDA).
52 // * mpirun -n 4 pmesh-optimizer -m cube.mesh -o 3 -rs 3 -mid 321 -tid 5 -ls 3 -nor -lc 0.1 -pa -d cuda
53 // Adapted discrete size; explicit combo of metrics; mixed tri/quad mesh:
54 // mpirun -np 4 pmesh-optimizer -m ../../data/square-mixed.mesh -o 2 -rs 2 -mid 2 -tid 5 -ni 200 -bnd -qo 6 -cmb 2 -nor
55 // Adapted discrete size+aspect_ratio:
56 // mpirun -np 4 pmesh-optimizer -m square01.mesh -o 2 -rs 2 -mid 7 -tid 6 -ni 100
57 // mpirun -np 4 pmesh-optimizer -m square01.mesh -o 2 -rs 2 -mid 7 -tid 6 -ni 100 -qo 6 -ex -st 1 -nor
58 // Adapted discrete size+orientation (requires GSLIB):
59 // * mpirun -np 4 pmesh-optimizer -m square01.mesh -o 2 -rs 2 -mid 36 -tid 8 -qo 4 -fd -ae 1 -nor
60 // Adapted discrete aspect_ratio+orientation (requires GSLIB):
61 // * mpirun -np 4 pmesh-optimizer -m square01.mesh -o 2 -rs 2 -mid 85 -tid 8 -ni 10 -bnd -qt 1 -qo 8 -fd -ae 1
62 // Adapted discrete aspect ratio (3D):
63 // mpirun -np 4 pmesh-optimizer -m cube.mesh -o 2 -rs 2 -mid 302 -tid 7 -ni 20 -bnd -qt 1 -qo 8
64 //
65 // Adaptive limiting:
66 // mpirun -np 4 pmesh-optimizer -m stretched2D.mesh -o 2 -mid 2 -tid 1 -ni 50 -qo 5 -nor -vl 1 -alc 0.5
67 // Adaptive limiting through the L-BFGS solver:
68 // mpirun -np 4 pmesh-optimizer -m stretched2D.mesh -o 2 -mid 2 -tid 1 -ni 400 -qo 5 -nor -vl 1 -alc 0.5 -st 1
69 // Adaptive limiting through FD (requires GSLIB):
70 // * mpirun -np 4 pmesh-optimizer -m stretched2D.mesh -o 2 -mid 2 -tid 1 -ni 50 -qo 5 -nor -vl 1 -alc 0.5 -fd -ae 1
71 //
72 // Adaptive surface fitting:
73 // mpirun -np 4 pmesh-optimizer -m square01.mesh -o 3 -rs 1 -mid 58 -tid 1 -ni 200 -vl 1 -sfc 5e4 -rtol 1e-5 -nor
74 // mpirun -np 4 pmesh-optimizer -m square01-tri.mesh -o 3 -rs 0 -mid 58 -tid 1 -ni 200 -vl 1 -sfc 1e4 -rtol 1e-5 -nor
75 // Surface fitting with weight adaptation and termination based on fitting error
76 // mpirun -np 4 pmesh-optimizer -m square01.mesh -o 2 -rs 1 -mid 2 -tid 1 -ni 100 -vl 2 -sfc 10 -rtol 1e-20 -st 0 -sfa -sft 1e-5
77 //
78 // Blade shape:
79 // mpirun -np 4 pmesh-optimizer -m blade.mesh -o 4 -mid 2 -tid 1 -ni 30 -ls 3 -art 1 -bnd -qt 1 -qo 8
80 // * mpirun -np 4 pmesh-optimizer -m blade.mesh -o 4 -mid 2 -tid 1 -ni 30 -ls 3 -art 1 -bnd -qt 1 -qo 8 -d cuda
81 // Blade shape with FD-based solver:
82 // mpirun -np 4 pmesh-optimizer -m blade.mesh -o 4 -mid 2 -tid 1 -ni 30 -ls 4 -bnd -qt 1 -qo 8 -fd
83 // Blade limited shape:
84 // mpirun -np 4 pmesh-optimizer -m blade.mesh -o 4 -mid 2 -tid 1 -bnd -qt 1 -qo 8 -lc 5000
85 // ICF shape and equal size:
86 // mpirun -np 4 pmesh-optimizer -o 3 -mid 9 -tid 2 -ni 25 -ls 3 -art 2 -qo 5
87 // ICF shape and initial size:
88 // mpirun -np 4 pmesh-optimizer -o 3 -mid 9 -tid 3 -ni 30 -ls 3 -bnd -qt 1 -qo 8
89 // ICF shape:
90 // mpirun -np 4 pmesh-optimizer -o 3 -mid 1 -tid 1 -ni 100 -bnd -qt 1 -qo 8
91 // ICF limited shape:
92 // mpirun -np 4 pmesh-optimizer -o 3 -mid 1 -tid 1 -ni 100 -bnd -qt 1 -qo 8 -lc 10
93 // ICF combo shape + size (rings, slow convergence):
94 // mpirun -np 4 pmesh-optimizer -o 3 -mid 1 -tid 1 -ni 1000 -bnd -qt 1 -qo 8 -cmb 1
95 // Mixed tet / cube / hex mesh with limiting:
96 // mpirun -np 4 pmesh-optimizer -m ../../data/fichera-mixed-p2.mesh -o 4 -rs 1 -mid 301 -tid 1 -fix-bnd -qo 6 -nor -lc 0.25
97 // 3D pinched sphere shape (the mesh is in the mfem/data GitHub repository):
98 // * mpirun -np 4 pmesh-optimizer -m ../../../mfem_data/ball-pert.mesh -o 4 -mid 303 -tid 1 -ni 20 -li 500 -fix-bnd
99 // 2D non-conforming shape and equal size:
100 // mpirun -np 4 pmesh-optimizer -m ./amr-quad-q2.mesh -o 2 -rs 1 -mid 9 -tid 2 -ni 200 -bnd -qt 1 -qo 8
101 //
102 // 2D untangling:
103 // mpirun -np 4 pmesh-optimizer -m jagged.mesh -o 2 -mid 22 -tid 1 -ni 50 -li 50 -qo 4 -fd -vl 1
104 // 3D untangling (the mesh is in the mfem/data GitHub repository):
105 // * mpirun -np 4 pmesh-optimizer -m ../../../mfem_data/cube-holes-inv.mesh -o 3 -mid 313 -tid 1 -rtol 1e-5 -li 50 -qo 4 -fd -vl 1
106 
107 #include "mfem.hpp"
108 #include "../common/mfem-common.hpp"
109 #include <iostream>
110 #include <fstream>
111 #include "mesh-optimizer.hpp"
112 
113 using namespace mfem;
114 using namespace std;
115 
116 int main (int argc, char *argv[])
117 {
118  // 0. Initialize MPI and HYPRE.
119  Mpi::Init(argc, argv);
120  int myid = Mpi::WorldRank();
121  Hypre::Init();
122 
123  // 1. Set the method's default parameters.
124  const char *mesh_file = "icf.mesh";
125  int mesh_poly_deg = 1;
126  int rs_levels = 0;
127  int rp_levels = 0;
128  double jitter = 0.0;
129  int metric_id = 1;
130  int target_id = 1;
131  double lim_const = 0.0;
132  double adapt_lim_const = 0.0;
133  double surface_fit_const = 0.0;
134  int quad_type = 1;
135  int quad_order = 8;
136  int solver_type = 0;
137  int solver_iter = 20;
138  double solver_rtol = 1e-10;
139  int solver_art_type = 0;
140  int lin_solver = 2;
141  int max_lin_iter = 100;
142  bool move_bnd = true;
143  int combomet = 0;
144  bool hradaptivity = false;
145  int h_metric_id = -1;
146  bool normalization = false;
147  bool visualization = true;
148  int verbosity_level = 0;
149  bool fdscheme = false;
150  int adapt_eval = 0;
151  bool exactaction = false;
152  const char *devopt = "cpu";
153  bool pa = false;
154  int n_hr_iter = 5;
155  int n_h_iter = 1;
156  bool surface_fit_adapt = false;
157  double surface_fit_threshold = -10;
158 
159  // 2. Parse command-line options.
160  OptionsParser args(argc, argv);
161  args.AddOption(&mesh_file, "-m", "--mesh",
162  "Mesh file to use.");
163  args.AddOption(&mesh_poly_deg, "-o", "--order",
164  "Polynomial degree of mesh finite element space.");
165  args.AddOption(&rs_levels, "-rs", "--refine-serial",
166  "Number of times to refine the mesh uniformly in serial.");
167  args.AddOption(&rp_levels, "-rp", "--refine-parallel",
168  "Number of times to refine the mesh uniformly in parallel.");
169  args.AddOption(&jitter, "-ji", "--jitter",
170  "Random perturbation scaling factor.");
171  args.AddOption(&metric_id, "-mid", "--metric-id",
172  "Mesh optimization metric:\n\t"
173  "T-metrics\n\t"
174  "1 : |T|^2 -- 2D shape\n\t"
175  "2 : 0.5|T|^2/tau-1 -- 2D shape (condition number)\n\t"
176  "7 : |T-T^-t|^2 -- 2D shape+size\n\t"
177  "9 : tau*|T-T^-t|^2 -- 2D shape+size\n\t"
178  "14 : |T-I|^2 -- 2D shape+size+orientation\n\t"
179  "22 : 0.5(|T|^2-2*tau)/(tau-tau_0) -- 2D untangling\n\t"
180  "50 : 0.5|T^tT|^2/tau^2-1 -- 2D shape\n\t"
181  "55 : (tau-1)^2 -- 2D size\n\t"
182  "56 : 0.5(sqrt(tau)-1/sqrt(tau))^2 -- 2D size\n\t"
183  "58 : |T^tT|^2/(tau^2)-2*|T|^2/tau+2 -- 2D shape\n\t"
184  "77 : 0.5(tau-1/tau)^2 -- 2D size\n\t"
185  "80 : (1-gamma)mu_2 + gamma mu_77 -- 2D shape+size\n\t"
186  "85 : |T-|T|/sqrt(2)I|^2 -- 2D shape+orientation\n\t"
187  "98 : (1/tau)|T-I|^2 -- 2D shape+size+orientation\n\t"
188  // "211: (tau-1)^2-tau+sqrt(tau^2+eps) -- 2D untangling\n\t"
189  // "252: 0.5(tau-1)^2/(tau-tau_0) -- 2D untangling\n\t"
190  "301: (|T||T^-1|)/3-1 -- 3D shape\n\t"
191  "302: (|T|^2|T^-1|^2)/9-1 -- 3D shape\n\t"
192  "303: (|T|^2)/3*tau^(2/3)-1 -- 3D shape\n\t"
193  // "311: (tau-1)^2-tau+sqrt(tau^2+eps)-- 3D untangling\n\t"
194  "313: (|T|^2)(tau-tau0)^(-2/3)/3 -- 3D untangling\n\t"
195  "315: (tau-1)^2 -- 3D size\n\t"
196  "316: 0.5(sqrt(tau)-1/sqrt(tau))^2 -- 3D size\n\t"
197  "321: |T-T^-t|^2 -- 3D shape+size\n\t"
198  // "352: 0.5(tau-1)^2/(tau-tau_0) -- 3D untangling\n\t"
199  "A-metrics\n\t"
200  "11 : (1/4*alpha)|A-(adjA)^T(W^TW)/omega|^2 -- 2D shape\n\t"
201  "36 : (1/alpha)|A-W|^2 -- 2D shape+size+orientation\n\t"
202  "107: (1/2*alpha)|A-|A|/|W|W|^2 -- 2D shape+orientation\n\t"
203  "126: (1-gamma)nu_11 + gamma*nu_14a -- 2D shape+size\n\t"
204  );
205  args.AddOption(&target_id, "-tid", "--target-id",
206  "Target (ideal element) type:\n\t"
207  "1: Ideal shape, unit size\n\t"
208  "2: Ideal shape, equal size\n\t"
209  "3: Ideal shape, initial size\n\t"
210  "4: Given full analytic Jacobian (in physical space)\n\t"
211  "5: Ideal shape, given size (in physical space)");
212  args.AddOption(&lim_const, "-lc", "--limit-const", "Limiting constant.");
213  args.AddOption(&adapt_lim_const, "-alc", "--adapt-limit-const",
214  "Adaptive limiting coefficient constant.");
215  args.AddOption(&surface_fit_const, "-sfc", "--surface-fit-const",
216  "Surface preservation constant.");
217  args.AddOption(&quad_type, "-qt", "--quad-type",
218  "Quadrature rule type:\n\t"
219  "1: Gauss-Lobatto\n\t"
220  "2: Gauss-Legendre\n\t"
221  "3: Closed uniform points");
222  args.AddOption(&quad_order, "-qo", "--quad_order",
223  "Order of the quadrature rule.");
224  args.AddOption(&solver_type, "-st", "--solver-type",
225  " Type of solver: (default) 0: Newton, 1: LBFGS");
226  args.AddOption(&solver_iter, "-ni", "--newton-iters",
227  "Maximum number of Newton iterations.");
228  args.AddOption(&solver_rtol, "-rtol", "--newton-rel-tolerance",
229  "Relative tolerance for the Newton solver.");
230  args.AddOption(&solver_art_type, "-art", "--adaptive-rel-tol",
231  "Type of adaptive relative linear solver tolerance:\n\t"
232  "0: None (default)\n\t"
233  "1: Eisenstat-Walker type 1\n\t"
234  "2: Eisenstat-Walker type 2");
235  args.AddOption(&lin_solver, "-ls", "--lin-solver",
236  "Linear solver:\n\t"
237  "0: l1-Jacobi\n\t"
238  "1: CG\n\t"
239  "2: MINRES\n\t"
240  "3: MINRES + Jacobi preconditioner\n\t"
241  "4: MINRES + l1-Jacobi preconditioner");
242  args.AddOption(&max_lin_iter, "-li", "--lin-iter",
243  "Maximum number of iterations in the linear solve.");
244  args.AddOption(&move_bnd, "-bnd", "--move-boundary", "-fix-bnd",
245  "--fix-boundary",
246  "Enable motion along horizontal and vertical boundaries.");
247  args.AddOption(&combomet, "-cmb", "--combo-type",
248  "Combination of metrics options:\n\t"
249  "0: Use single metric\n\t"
250  "1: Shape + space-dependent size given analytically\n\t"
251  "2: Shape + adapted size given discretely; shared target");
252  args.AddOption(&hradaptivity, "-hr", "--hr-adaptivity", "-no-hr",
253  "--no-hr-adaptivity",
254  "Enable hr-adaptivity.");
255  args.AddOption(&h_metric_id, "-hmid", "--h-metric",
256  "Same options as metric_id. Used to determine refinement"
257  " type for each element if h-adaptivity is enabled.");
258  args.AddOption(&normalization, "-nor", "--normalization", "-no-nor",
259  "--no-normalization",
260  "Make all terms in the optimization functional unitless.");
261  args.AddOption(&fdscheme, "-fd", "--fd_approximation",
262  "-no-fd", "--no-fd-approx",
263  "Enable finite difference based derivative computations.");
264  args.AddOption(&exactaction, "-ex", "--exact_action",
265  "-no-ex", "--no-exact-action",
266  "Enable exact action of TMOP_Integrator.");
267  args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
268  "--no-visualization",
269  "Enable or disable GLVis visualization.");
270  args.AddOption(&verbosity_level, "-vl", "--verbosity-level",
271  "Set the verbosity level - 0, 1, or 2.");
272  args.AddOption(&adapt_eval, "-ae", "--adaptivity-evaluator",
273  "0 - Advection based (DEFAULT), 1 - GSLIB.");
274  args.AddOption(&devopt, "-d", "--device",
275  "Device configuration string, see Device::Configure().");
276  args.AddOption(&pa, "-pa", "--partial-assembly", "-no-pa",
277  "--no-partial-assembly", "Enable Partial Assembly.");
278  args.AddOption(&n_hr_iter, "-nhr", "--n_hr_iter",
279  "Number of hr-adaptivity iterations.");
280  args.AddOption(&n_h_iter, "-nh", "--n_h_iter",
281  "Number of h-adaptivity iterations per r-adaptivity"
282  "iteration.");
283  args.AddOption(&surface_fit_adapt, "-sfa", "--adaptive-surface-fit", "-no-sfa",
284  "--no-adaptive-surface-fit",
285  "Enable or disable adaptive surface fitting.");
286  args.AddOption(&surface_fit_threshold, "-sft", "--surf-fit-threshold",
287  "Set threshold for surface fitting. TMOP solver will"
288  "terminate when max surface fitting error is below this limit");
289  args.Parse();
290  if (!args.Good())
291  {
292  if (myid == 0) { args.PrintUsage(cout); }
293  return 1;
294  }
295  if (myid == 0) { args.PrintOptions(cout); }
296  if (h_metric_id < 0) { h_metric_id = metric_id; }
297 
298  if (hradaptivity)
299  {
300  MFEM_VERIFY(strcmp(devopt,"cpu")==0, "HR-adaptivity is currently only"
301  " supported on cpus.");
302  }
303  Device device(devopt);
304  if (myid == 0) { device.Print();}
305 
306  // 3. Initialize and refine the starting mesh.
307  Mesh *mesh = new Mesh(mesh_file, 1, 1, false);
308  for (int lev = 0; lev < rs_levels; lev++)
309  {
310  mesh->UniformRefinement();
311  }
312  const int dim = mesh->Dimension();
313 
314  if (hradaptivity) { mesh->EnsureNCMesh(); }
315  ParMesh *pmesh = new ParMesh(MPI_COMM_WORLD, *mesh);
316 
317  delete mesh;
318  for (int lev = 0; lev < rp_levels; lev++)
319  {
320  pmesh->UniformRefinement();
321  }
322 
323  // 4. Define a finite element space on the mesh. Here we use vector finite
324  // elements which are tensor products of quadratic finite elements. The
325  // number of components in the vector finite element space is specified by
326  // the last parameter of the FiniteElementSpace constructor.
328  if (mesh_poly_deg <= 0)
329  {
330  fec = new QuadraticPosFECollection;
331  mesh_poly_deg = 2;
332  }
333  else { fec = new H1_FECollection(mesh_poly_deg, dim); }
334  ParFiniteElementSpace *pfespace = new ParFiniteElementSpace(pmesh, fec, dim);
335 
336  // 5. Make the mesh curved based on the above finite element space. This
337  // means that we define the mesh elements through a fespace-based
338  // transformation of the reference element.
339  pmesh->SetNodalFESpace(pfespace);
340 
341  // 6. Set up an empty right-hand side vector b, which is equivalent to b=0.
342  Vector b(0);
343 
344  // 7. Get the mesh nodes (vertices and other degrees of freedom in the finite
345  // element space) as a finite element grid function in fespace. Note that
346  // changing x automatically changes the shapes of the mesh elements.
347  ParGridFunction x(pfespace);
348  pmesh->SetNodalGridFunction(&x);
349 
350  // 8. Define a vector representing the minimal local mesh size in the mesh
351  // nodes. We index the nodes using the scalar version of the degrees of
352  // freedom in pfespace. Note: this is partition-dependent.
353  //
354  // In addition, compute average mesh size and total volume.
355  Vector h0(pfespace->GetNDofs());
356  h0 = infinity();
357  double vol_loc = 0.0;
358  Array<int> dofs;
359  for (int i = 0; i < pmesh->GetNE(); i++)
360  {
361  // Get the local scalar element degrees of freedom in dofs.
362  pfespace->GetElementDofs(i, dofs);
363  // Adjust the value of h0 in dofs based on the local mesh size.
364  const double hi = pmesh->GetElementSize(i);
365  for (int j = 0; j < dofs.Size(); j++)
366  {
367  h0(dofs[j]) = min(h0(dofs[j]), hi);
368  }
369  vol_loc += pmesh->GetElementVolume(i);
370  }
371  double vol_glb;
372  MPI_Allreduce(&vol_loc, &vol_glb, 1, MPI_DOUBLE, MPI_SUM, MPI_COMM_WORLD);
373  const double small_phys_size = pow(vol_glb, 1.0 / dim) / 100.0;
374 
375  // 9. Add a random perturbation to the nodes in the interior of the domain.
376  // We define a random grid function of fespace and make sure that it is
377  // zero on the boundary and its values are locally of the order of h0.
378  // The latter is based on the DofToVDof() method which maps the scalar to
379  // the vector degrees of freedom in pfespace.
380  ParGridFunction rdm(pfespace);
381  rdm.Randomize();
382  rdm -= 0.25; // Shift to random values in [-0.5,0.5].
383  rdm *= jitter;
384  rdm.HostReadWrite();
385  // Scale the random values to be of order of the local mesh size.
386  for (int i = 0; i < pfespace->GetNDofs(); i++)
387  {
388  for (int d = 0; d < dim; d++)
389  {
390  rdm(pfespace->DofToVDof(i,d)) *= h0(i);
391  }
392  }
393  Array<int> vdofs;
394  for (int i = 0; i < pfespace->GetNBE(); i++)
395  {
396  // Get the vector degrees of freedom in the boundary element.
397  pfespace->GetBdrElementVDofs(i, vdofs);
398  // Set the boundary values to zero.
399  for (int j = 0; j < vdofs.Size(); j++) { rdm(vdofs[j]) = 0.0; }
400  }
401  x -= rdm;
402  // Set the perturbation of all nodes from the true nodes.
403  x.SetTrueVector();
404  x.SetFromTrueVector();
405 
406  // 10. Save the starting (prior to the optimization) mesh to a file. This
407  // output can be viewed later using GLVis: "glvis -m perturbed -np
408  // num_mpi_tasks".
409  {
410  ostringstream mesh_name;
411  mesh_name << "perturbed.mesh";
412  ofstream mesh_ofs(mesh_name.str().c_str());
413  mesh_ofs.precision(8);
414  pmesh->PrintAsOne(mesh_ofs);
415  }
416 
417  // 11. Store the starting (prior to the optimization) positions.
418  ParGridFunction x0(pfespace);
419  x0 = x;
420 
421  // 12. Form the integrator that uses the chosen metric and target.
422  double tauval = -0.1;
423  TMOP_QualityMetric *metric = NULL;
424  switch (metric_id)
425  {
426  // T-metrics
427  case 1: metric = new TMOP_Metric_001; break;
428  case 2: metric = new TMOP_Metric_002; break;
429  case 7: metric = new TMOP_Metric_007; break;
430  case 9: metric = new TMOP_Metric_009; break;
431  case 14: metric = new TMOP_Metric_014; break;
432  case 22: metric = new TMOP_Metric_022(tauval); break;
433  case 50: metric = new TMOP_Metric_050; break;
434  case 55: metric = new TMOP_Metric_055; break;
435  case 56: metric = new TMOP_Metric_056; break;
436  case 58: metric = new TMOP_Metric_058; break;
437  case 77: metric = new TMOP_Metric_077; break;
438  case 80: metric = new TMOP_Metric_080(0.5); break;
439  case 85: metric = new TMOP_Metric_085; break;
440  case 98: metric = new TMOP_Metric_098; break;
441  // case 211: metric = new TMOP_Metric_211; break;
442  // case 252: metric = new TMOP_Metric_252(tauval); break;
443  case 301: metric = new TMOP_Metric_301; break;
444  case 302: metric = new TMOP_Metric_302; break;
445  case 303: metric = new TMOP_Metric_303; break;
446  // case 311: metric = new TMOP_Metric_311; break;
447  case 313: metric = new TMOP_Metric_313(tauval); break;
448  case 315: metric = new TMOP_Metric_315; break;
449  case 316: metric = new TMOP_Metric_316; break;
450  case 321: metric = new TMOP_Metric_321; break;
451  case 328: metric = new TMOP_Metric_328(0.5); break;
452  case 332: metric = new TMOP_Metric_332(0.5); break;
453  case 333: metric = new TMOP_Metric_333(0.5); break;
454  case 334: metric = new TMOP_Metric_334(0.5); break;
455  // case 352: metric = new TMOP_Metric_352(tauval); break;
456  // A-metrics
457  case 11: metric = new TMOP_AMetric_011; break;
458  case 36: metric = new TMOP_AMetric_036; break;
459  case 107: metric = new TMOP_AMetric_107a; break;
460  case 126: metric = new TMOP_AMetric_126(0.9); break;
461  default:
462  if (myid == 0) { cout << "Unknown metric_id: " << metric_id << endl; }
463  return 3;
464  }
465  TMOP_QualityMetric *h_metric = NULL;
466  if (hradaptivity)
467  {
468  switch (h_metric_id)
469  {
470  case 1: h_metric = new TMOP_Metric_001; break;
471  case 2: h_metric = new TMOP_Metric_002; break;
472  case 7: h_metric = new TMOP_Metric_007; break;
473  case 9: h_metric = new TMOP_Metric_009; break;
474  case 55: h_metric = new TMOP_Metric_055; break;
475  case 56: h_metric = new TMOP_Metric_056; break;
476  case 58: h_metric = new TMOP_Metric_058; break;
477  case 77: h_metric = new TMOP_Metric_077; break;
478  case 315: h_metric = new TMOP_Metric_315; break;
479  case 316: h_metric = new TMOP_Metric_316; break;
480  case 321: h_metric = new TMOP_Metric_321; break;
481  default: cout << "Metric_id not supported for h-adaptivity: " << h_metric_id <<
482  endl;
483  return 3;
484  }
485  }
486 
487  if (metric_id < 300 || h_metric_id < 300)
488  {
489  MFEM_VERIFY(dim == 2, "Incompatible metric for 3D meshes");
490  }
491  if (metric_id >= 300 || h_metric_id >= 300)
492  {
493  MFEM_VERIFY(dim == 3, "Incompatible metric for 2D meshes");
494  }
495 
497  TargetConstructor *target_c = NULL;
498  HessianCoefficient *adapt_coeff = NULL;
499  HRHessianCoefficient *hr_adapt_coeff = NULL;
500  H1_FECollection ind_fec(mesh_poly_deg, dim);
501  ParFiniteElementSpace ind_fes(pmesh, &ind_fec);
502  ParFiniteElementSpace ind_fesv(pmesh, &ind_fec, dim);
503  ParGridFunction size(&ind_fes), aspr(&ind_fes), ori(&ind_fes);
504  ParGridFunction aspr3d(&ind_fesv);
505 
506  const AssemblyLevel al =
508 
509  switch (target_id)
510  {
511  case 1: target_t = TargetConstructor::IDEAL_SHAPE_UNIT_SIZE; break;
512  case 2: target_t = TargetConstructor::IDEAL_SHAPE_EQUAL_SIZE; break;
513  case 3: target_t = TargetConstructor::IDEAL_SHAPE_GIVEN_SIZE; break;
514  case 4:
515  {
517  AnalyticAdaptTC *tc = new AnalyticAdaptTC(target_t);
518  adapt_coeff = new HessianCoefficient(dim, metric_id);
519  tc->SetAnalyticTargetSpec(NULL, NULL, adapt_coeff);
520  target_c = tc;
521  break;
522  }
523  case 5: // Discrete size 2D or 3D
524  {
526  DiscreteAdaptTC *tc = new DiscreteAdaptTC(target_t);
527  if (adapt_eval == 0)
528  {
529  tc->SetAdaptivityEvaluator(new AdvectorCG(al));
530  }
531  else
532  {
533 #ifdef MFEM_USE_GSLIB
535 #else
536  MFEM_ABORT("MFEM is not built with GSLIB.");
537 #endif
538  }
539  if (dim == 2)
540  {
542  size.ProjectCoefficient(size_coeff);
543  }
544  else if (dim == 3)
545  {
547  size.ProjectCoefficient(size_coeff);
548  }
549  tc->SetParDiscreteTargetSize(size);
550  target_c = tc;
551  break;
552  }
553  case 6: // material indicator 2D
554  {
555  ParGridFunction d_x(&ind_fes), d_y(&ind_fes), disc(&ind_fes);
556 
558  DiscreteAdaptTC *tc = new DiscreteAdaptTC(target_t);
560  disc.ProjectCoefficient(mat_coeff);
561  if (adapt_eval == 0)
562  {
563  tc->SetAdaptivityEvaluator(new AdvectorCG(al));
564  }
565  else
566  {
567 #ifdef MFEM_USE_GSLIB
569 #else
570  MFEM_ABORT("MFEM is not built with GSLIB.");
571 #endif
572  }
573  // Diffuse the interface
574  DiffuseField(disc,2);
575 
576  // Get partials with respect to x and y of the grid function
577  disc.GetDerivative(1,0,d_x);
578  disc.GetDerivative(1,1,d_y);
579 
580  // Compute the squared magnitude of the gradient
581  for (int i = 0; i < size.Size(); i++)
582  {
583  size(i) = std::pow(d_x(i),2)+std::pow(d_y(i),2);
584  }
585  const double max = size.Max();
586  double max_all;
587  MPI_Allreduce(&max, &max_all, 1, MPI_DOUBLE, MPI_MAX, MPI_COMM_WORLD);
588 
589  for (int i = 0; i < d_x.Size(); i++)
590  {
591  d_x(i) = std::abs(d_x(i));
592  d_y(i) = std::abs(d_y(i));
593  }
594  const double eps = 0.01;
595  const double aspr_ratio = 20.0;
596  const double size_ratio = 40.0;
597 
598  for (int i = 0; i < size.Size(); i++)
599  {
600  size(i) = (size(i)/max_all);
601  aspr(i) = (d_x(i)+eps)/(d_y(i)+eps);
602  aspr(i) = 0.1 + 0.9*(1-size(i))*(1-size(i));
603  if (aspr(i) > aspr_ratio) {aspr(i) = aspr_ratio;}
604  if (aspr(i) < 1.0/aspr_ratio) {aspr(i) = 1.0/aspr_ratio;}
605  }
606  Vector vals;
607  const int NE = pmesh->GetNE();
608  double volume = 0.0, volume_ind = 0.0;
609 
610  for (int i = 0; i < NE; i++)
611  {
613  const IntegrationRule &ir =
614  IntRules.Get(pmesh->GetElementBaseGeometry(i), Tr->OrderJ());
615  size.GetValues(i, ir, vals);
616  for (int j = 0; j < ir.GetNPoints(); j++)
617  {
618  const IntegrationPoint &ip = ir.IntPoint(j);
619  Tr->SetIntPoint(&ip);
620  volume += ip.weight * Tr->Weight();
621  volume_ind += vals(j) * ip.weight * Tr->Weight();
622  }
623  }
624  double volume_all, volume_ind_all;
625  MPI_Allreduce(&volume, &volume_all, 1, MPI_DOUBLE, MPI_SUM, MPI_COMM_WORLD);
626  MPI_Allreduce(&volume_ind, &volume_ind_all, 1, MPI_DOUBLE, MPI_SUM,
627  MPI_COMM_WORLD);
628  const int NE_ALL = pmesh->GetGlobalNE();
629 
630  const double avg_zone_size = volume_all / NE_ALL;
631 
632  const double small_avg_ratio =
633  (volume_ind_all + (volume_all - volume_ind_all) / size_ratio)
634  / volume_all;
635 
636  const double small_zone_size = small_avg_ratio * avg_zone_size;
637  const double big_zone_size = size_ratio * small_zone_size;
638 
639  for (int i = 0; i < size.Size(); i++)
640  {
641  const double val = size(i);
642  const double a = (big_zone_size - small_zone_size) / small_zone_size;
643  size(i) = big_zone_size / (1.0+a*val);
644  }
645 
646  DiffuseField(size, 2);
647  DiffuseField(aspr, 2);
648 
649  tc->SetParDiscreteTargetSize(size);
651  target_c = tc;
652  break;
653  }
654  case 7: // Discrete aspect ratio 3D
655  {
657  DiscreteAdaptTC *tc = new DiscreteAdaptTC(target_t);
658  if (adapt_eval == 0)
659  {
660  tc->SetAdaptivityEvaluator(new AdvectorCG(al));
661  }
662  else
663  {
664 #ifdef MFEM_USE_GSLIB
666 #else
667  MFEM_ABORT("MFEM is not built with GSLIB.");
668 #endif
669  }
671  aspr3d.ProjectCoefficient(fd_aspr3d);
673  target_c = tc;
674  break;
675  }
676  case 8: // shape/size + orientation 2D
677  {
679  DiscreteAdaptTC *tc = new DiscreteAdaptTC(target_t);
680  if (adapt_eval == 0)
681  {
682  tc->SetAdaptivityEvaluator(new AdvectorCG(al));
683  }
684  else
685  {
686 #ifdef MFEM_USE_GSLIB
688 #else
689  MFEM_ABORT("MFEM is not built with GSLIB.");
690 #endif
691  }
692 
693  if (metric_id == 14 || metric_id == 36)
694  {
695  ConstantCoefficient size_coeff(0.1*0.1);
696  size.ProjectCoefficient(size_coeff);
697  tc->SetParDiscreteTargetSize(size);
698  }
699 
700  if (metric_id == 85)
701  {
703  aspr.ProjectCoefficient(aspr_coeff);
704  DiffuseField(aspr,2);
706  }
707 
709  ori.ProjectCoefficient(ori_coeff);
711  target_c = tc;
712  break;
713  }
714  // Targets used for hr-adaptivity tests.
715  case 9: // size target in an annular region.
716  case 10: // size+aspect-ratio in an annular region.
717  case 11: // size+aspect-ratio target for a rotate sine wave
718  {
720  AnalyticAdaptTC *tc = new AnalyticAdaptTC(target_t);
721  hr_adapt_coeff = new HRHessianCoefficient(dim, target_id - 9);
722  tc->SetAnalyticTargetSpec(NULL, NULL, hr_adapt_coeff);
723  target_c = tc;
724  break;
725  }
726  default:
727  if (myid == 0) { cout << "Unknown target_id: " << target_id << endl; }
728  return 3;
729  }
730 
731  if (target_c == NULL)
732  {
733  target_c = new TargetConstructor(target_t, MPI_COMM_WORLD);
734  }
735  target_c->SetNodes(x0);
736  TMOP_Integrator *tmop_integ = new TMOP_Integrator(metric, target_c,
737  h_metric);
738 
739  // Finite differences for computations of derivatives.
740  if (fdscheme)
741  {
742  MFEM_VERIFY(pa == false, "PA for finite differences is not implemented.");
743  tmop_integ->EnableFiniteDifferences(x);
744  }
745  tmop_integ->SetExactActionFlag(exactaction);
746 
747  // Setup the quadrature rules for the TMOP integrator.
748  IntegrationRules *irules = NULL;
749  switch (quad_type)
750  {
751  case 1: irules = &IntRulesLo; break;
752  case 2: irules = &IntRules; break;
753  case 3: irules = &IntRulesCU; break;
754  default:
755  if (myid == 0) { cout << "Unknown quad_type: " << quad_type << endl; }
756  return 3;
757  }
758  tmop_integ->SetIntegrationRules(*irules, quad_order);
759  if (myid == 0 && dim == 2)
760  {
761  cout << "Triangle quadrature points: "
762  << irules->Get(Geometry::TRIANGLE, quad_order).GetNPoints()
763  << "\nQuadrilateral quadrature points: "
764  << irules->Get(Geometry::SQUARE, quad_order).GetNPoints() << endl;
765  }
766  if (myid == 0 && dim == 3)
767  {
768  cout << "Tetrahedron quadrature points: "
769  << irules->Get(Geometry::TETRAHEDRON, quad_order).GetNPoints()
770  << "\nHexahedron quadrature points: "
771  << irules->Get(Geometry::CUBE, quad_order).GetNPoints()
772  << "\nPrism quadrature points: "
773  << irules->Get(Geometry::PRISM, quad_order).GetNPoints() << endl;
774  }
775 
776  // Limit the node movement.
777  // The limiting distances can be given by a general function of space.
778  ParFiniteElementSpace dist_pfespace(pmesh, fec); // scalar space
779  ParGridFunction dist(&dist_pfespace);
780  dist = 1.0;
781  // The small_phys_size is relevant only with proper normalization.
782  if (normalization) { dist = small_phys_size; }
783  ConstantCoefficient lim_coeff(lim_const);
784  if (lim_const != 0.0) { tmop_integ->EnableLimiting(x0, dist, lim_coeff); }
785 
786  // Adaptive limiting.
787  ParGridFunction adapt_lim_gf0(&ind_fes);
788  ConstantCoefficient adapt_lim_coeff(adapt_lim_const);
789  AdaptivityEvaluator *adapt_lim_eval = NULL;
790  if (adapt_lim_const > 0.0)
791  {
792  MFEM_VERIFY(pa == false, "PA is not implemented for adaptive limiting");
793 
794  FunctionCoefficient adapt_lim_gf0_coeff(adapt_lim_fun);
795  adapt_lim_gf0.ProjectCoefficient(adapt_lim_gf0_coeff);
796 
797  if (adapt_eval == 0) { adapt_lim_eval = new AdvectorCG(al); }
798  else if (adapt_eval == 1)
799  {
800 #ifdef MFEM_USE_GSLIB
801  adapt_lim_eval = new InterpolatorFP;
802 #else
803  MFEM_ABORT("MFEM is not built with GSLIB support!");
804 #endif
805  }
806  else { MFEM_ABORT("Bad interpolation option."); }
807 
808  tmop_integ->EnableAdaptiveLimiting(adapt_lim_gf0, adapt_lim_coeff,
809  *adapt_lim_eval);
810  if (visualization)
811  {
812  socketstream vis1;
813  common::VisualizeField(vis1, "localhost", 19916, adapt_lim_gf0, "Zeta 0",
814  300, 600, 300, 300);
815  }
816  }
817 
818  // Surface fitting.
819  L2_FECollection mat_coll(0, dim);
820  H1_FECollection surf_fit_fec(mesh_poly_deg, dim);
821  ParFiniteElementSpace surf_fit_fes(pmesh, &surf_fit_fec);
822  ParFiniteElementSpace mat_fes(pmesh, &mat_coll);
823  ParGridFunction mat(&mat_fes);
824  ParGridFunction surf_fit_mat_gf(&surf_fit_fes);
825  ParGridFunction surf_fit_gf0(&surf_fit_fes);
826  Array<bool> surf_fit_marker(surf_fit_gf0.Size());
827  ConstantCoefficient surf_fit_coeff(surface_fit_const);
828  AdaptivityEvaluator *adapt_surface = NULL;
829  if (surface_fit_const > 0.0)
830  {
831  MFEM_VERIFY(hradaptivity == false,
832  "Surface fitting with HR is not implemented yet.");
833  MFEM_VERIFY(pa == false,
834  "Surface fitting with PA is not implemented yet.");
835 
837  surf_fit_gf0.ProjectCoefficient(ls_coeff);
838 
839  for (int i = 0; i < pmesh->GetNE(); i++)
840  {
841  mat(i) = material_id(i, surf_fit_gf0);
842  pmesh->SetAttribute(i, mat(i) + 1);
843  }
844 
845  GridFunctionCoefficient coeff_mat(&mat);
846  surf_fit_mat_gf.ProjectDiscCoefficient(coeff_mat, GridFunction::ARITHMETIC);
847  for (int j = 0; j < surf_fit_marker.Size(); j++)
848  {
849  if (surf_fit_mat_gf(j) > 0.1 && surf_fit_mat_gf(j) < 0.9)
850  {
851  surf_fit_marker[j] = true;
852  surf_fit_mat_gf(j) = 1.0;
853  }
854  else
855  {
856  surf_fit_marker[j] = false;
857  surf_fit_mat_gf(j) = 0.0;
858  }
859  }
860 
861  if (adapt_eval == 0) { adapt_surface = new AdvectorCG; }
862  else if (adapt_eval == 1)
863  {
864 #ifdef MFEM_USE_GSLIB
865  adapt_surface = new InterpolatorFP;
866 #else
867  MFEM_ABORT("MFEM is not built with GSLIB support!");
868 #endif
869  }
870  else { MFEM_ABORT("Bad interpolation option."); }
871 
872  tmop_integ->EnableSurfaceFitting(surf_fit_gf0, surf_fit_marker, surf_fit_coeff,
873  *adapt_surface);
874  if (visualization)
875  {
876  socketstream vis1, vis2, vis3;
877  common::VisualizeField(vis1, "localhost", 19916, surf_fit_gf0, "Level Set 0",
878  300, 600, 300, 300);
879  common::VisualizeField(vis2, "localhost", 19916, mat, "Materials",
880  600, 600, 300, 300);
881  common::VisualizeField(vis3, "localhost", 19916, surf_fit_mat_gf,
882  "Dofs to Move",
883  900, 600, 300, 300);
884  }
885  }
886 
887  // Has to be after the enabling of the limiting / alignment, as it computes
888  // normalization factors for these terms as well.
889  if (normalization) { tmop_integ->ParEnableNormalization(x0); }
890 
891  // 13. Setup the final NonlinearForm (which defines the integral of interest,
892  // its first and second derivatives). Here we can use a combination of
893  // metrics, i.e., optimize the sum of two integrals, where both are
894  // scaled by used-defined space-dependent weights. Note that there are
895  // no command-line options for the weights and the type of the second
896  // metric; one should update those in the code.
897  ParNonlinearForm a(pfespace);
899  ConstantCoefficient *metric_coeff1 = NULL;
900  TMOP_QualityMetric *metric2 = NULL;
901  TargetConstructor *target_c2 = NULL;
902  FunctionCoefficient metric_coeff2(weight_fun);
903 
904  // Explicit combination of metrics.
905  if (combomet > 0)
906  {
907  // First metric.
908  metric_coeff1 = new ConstantCoefficient(1.0);
909  tmop_integ->SetCoefficient(*metric_coeff1);
910 
911  // Second metric.
912  if (dim == 2) { metric2 = new TMOP_Metric_077; }
913  else { metric2 = new TMOP_Metric_315; }
914  TMOP_Integrator *tmop_integ2 = NULL;
915  if (combomet == 1)
916  {
917  target_c2 = new TargetConstructor(
919  target_c2->SetVolumeScale(0.01);
920  target_c2->SetNodes(x0);
921  tmop_integ2 = new TMOP_Integrator(metric2, target_c2, h_metric);
922  tmop_integ2->SetCoefficient(metric_coeff2);
923  }
924  else { tmop_integ2 = new TMOP_Integrator(metric2, target_c, h_metric); }
925  tmop_integ2->SetIntegrationRules(*irules, quad_order);
926  if (fdscheme) { tmop_integ2->EnableFiniteDifferences(x); }
927  tmop_integ2->SetExactActionFlag(exactaction);
928 
930  combo->AddTMOPIntegrator(tmop_integ);
931  combo->AddTMOPIntegrator(tmop_integ2);
932  if (normalization) { combo->ParEnableNormalization(x0); }
933  if (lim_const != 0.0) { combo->EnableLimiting(x0, dist, lim_coeff); }
934 
935  a.AddDomainIntegrator(combo);
936  }
937  else
938  {
939  a.AddDomainIntegrator(tmop_integ);
940  }
941 
942  if (pa) { a.Setup(); }
943 
944  // Compute the minimum det(J) of the starting mesh.
945  tauval = infinity();
946  const int NE = pmesh->GetNE();
947  for (int i = 0; i < NE; i++)
948  {
949  const IntegrationRule &ir =
950  irules->Get(pfespace->GetFE(i)->GetGeomType(), quad_order);
952  for (int j = 0; j < ir.GetNPoints(); j++)
953  {
954  transf->SetIntPoint(&ir.IntPoint(j));
955  tauval = min(tauval, transf->Jacobian().Det());
956  }
957  }
958  double minJ0;
959  MPI_Allreduce(&tauval, &minJ0, 1, MPI_DOUBLE, MPI_MIN, MPI_COMM_WORLD);
960  tauval = minJ0;
961  if (myid == 0)
962  { cout << "Minimum det(J) of the original mesh is " << tauval << endl; }
963 
964  if (tauval < 0.0 && metric_id != 22 && metric_id != 211 && metric_id != 252
965  && metric_id != 311 && metric_id != 313 && metric_id != 352)
966  {
967  MFEM_ABORT("The input mesh is inverted! Try an untangling metric.");
968  }
969  if (tauval < 0.0)
970  {
971  MFEM_VERIFY(target_t == TargetConstructor::IDEAL_SHAPE_UNIT_SIZE,
972  "Untangling is supported only for ideal targets.");
973 
974  const DenseMatrix &Wideal =
976  tauval /= Wideal.Det();
977 
978  double h0min = h0.Min(), h0min_all;
979  MPI_Allreduce(&h0min, &h0min_all, 1, MPI_DOUBLE, MPI_MIN, MPI_COMM_WORLD);
980  // Slightly below minJ0 to avoid div by 0.
981  tauval -= 0.01 * h0min_all;
982  }
983 
984  // For HR tests, the energy is normalized by the number of elements.
985  const double init_energy = a.GetParGridFunctionEnergy(x) /
986  (hradaptivity ? pmesh->GetGlobalNE() : 1);
987  double init_metric_energy = init_energy;
988  if (lim_const > 0.0 || adapt_lim_const > 0.0 || surface_fit_const > 0.0)
989  {
990  lim_coeff.constant = 0.0;
991  adapt_lim_coeff.constant = 0.0;
992  surf_fit_coeff.constant = 0.0;
993  init_metric_energy = a.GetParGridFunctionEnergy(x) /
994  (hradaptivity ? pmesh->GetGlobalNE() : 1);
995  lim_coeff.constant = lim_const;
996  adapt_lim_coeff.constant = adapt_lim_const;
997  surf_fit_coeff.constant = surface_fit_const;
998  }
999 
1000  // Visualize the starting mesh and metric values.
1001  // Note that for combinations of metrics, this only shows the first metric.
1002  if (visualization)
1003  {
1004  char title[] = "Initial metric values";
1005  vis_tmop_metric_p(mesh_poly_deg, *metric, *target_c, *pmesh, title, 0);
1006  }
1007 
1008  // 14. Fix all boundary nodes, or fix only a given component depending on the
1009  // boundary attributes of the given mesh. Attributes 1/2/3 correspond to
1010  // fixed x/y/z components of the node. Attribute dim+1 corresponds to
1011  // an entirely fixed node.
1012  if (move_bnd == false)
1013  {
1014  Array<int> ess_bdr(pmesh->bdr_attributes.Max());
1015  ess_bdr = 1;
1016  a.SetEssentialBC(ess_bdr);
1017  }
1018  else
1019  {
1020  int n = 0;
1021  for (int i = 0; i < pmesh->GetNBE(); i++)
1022  {
1023  const int nd = pfespace->GetBE(i)->GetDof();
1024  const int attr = pmesh->GetBdrElement(i)->GetAttribute();
1025  MFEM_VERIFY(!(dim == 2 && attr == 3),
1026  "Boundary attribute 3 must be used only for 3D meshes. "
1027  "Adjust the attributes (1/2/3/4 for fixed x/y/z/all "
1028  "components, rest for free nodes), or use -fix-bnd.");
1029  if (attr == 1 || attr == 2 || attr == 3) { n += nd; }
1030  if (attr == 4) { n += nd * dim; }
1031  }
1032  Array<int> ess_vdofs(n);
1033  n = 0;
1034  for (int i = 0; i < pmesh->GetNBE(); i++)
1035  {
1036  const int nd = pfespace->GetBE(i)->GetDof();
1037  const int attr = pmesh->GetBdrElement(i)->GetAttribute();
1038  pfespace->GetBdrElementVDofs(i, vdofs);
1039  if (attr == 1) // Fix x components.
1040  {
1041  for (int j = 0; j < nd; j++)
1042  { ess_vdofs[n++] = vdofs[j]; }
1043  }
1044  else if (attr == 2) // Fix y components.
1045  {
1046  for (int j = 0; j < nd; j++)
1047  { ess_vdofs[n++] = vdofs[j+nd]; }
1048  }
1049  else if (attr == 3) // Fix z components.
1050  {
1051  for (int j = 0; j < nd; j++)
1052  { ess_vdofs[n++] = vdofs[j+2*nd]; }
1053  }
1054  else if (attr == 4) // Fix all components.
1055  {
1056  for (int j = 0; j < vdofs.Size(); j++)
1057  { ess_vdofs[n++] = vdofs[j]; }
1058  }
1059  }
1060  a.SetEssentialVDofs(ess_vdofs);
1061  }
1062 
1063  // 15. As we use the Newton method to solve the resulting nonlinear system,
1064  // here we setup the linear solver for the system's Jacobian.
1065  Solver *S = NULL, *S_prec = NULL;
1066  const double linsol_rtol = 1e-12;
1067  if (lin_solver == 0)
1068  {
1069  S = new DSmoother(1, 1.0, max_lin_iter);
1070  }
1071  else if (lin_solver == 1)
1072  {
1073  CGSolver *cg = new CGSolver(MPI_COMM_WORLD);
1074  cg->SetMaxIter(max_lin_iter);
1075  cg->SetRelTol(linsol_rtol);
1076  cg->SetAbsTol(0.0);
1077  cg->SetPrintLevel(verbosity_level >= 2 ? 3 : -1);
1078  S = cg;
1079  }
1080  else
1081  {
1082  MINRESSolver *minres = new MINRESSolver(MPI_COMM_WORLD);
1083  minres->SetMaxIter(max_lin_iter);
1084  minres->SetRelTol(linsol_rtol);
1085  minres->SetAbsTol(0.0);
1086  if (verbosity_level > 2) { minres->SetPrintLevel(1); }
1087  else { minres->SetPrintLevel(verbosity_level == 2 ? 3 : -1); }
1088  if (lin_solver == 3 || lin_solver == 4)
1089  {
1090  if (pa)
1091  {
1092  MFEM_VERIFY(lin_solver != 4, "PA l1-Jacobi is not implemented");
1093  auto js = new OperatorJacobiSmoother;
1094  js->SetPositiveDiagonal(true);
1095  S_prec = js;
1096  }
1097  else
1098  {
1099  auto hs = new HypreSmoother;
1100  hs->SetType((lin_solver == 3) ? HypreSmoother::Jacobi
1101  /* */ : HypreSmoother::l1Jacobi, 1);
1102  hs->SetPositiveDiagonal(true);
1103  S_prec = hs;
1104  }
1105  minres->SetPreconditioner(*S_prec);
1106  }
1107  S = minres;
1108  }
1109 
1110  // Perform the nonlinear optimization.
1111  const IntegrationRule &ir =
1112  irules->Get(pfespace->GetFE(0)->GetGeomType(), quad_order);
1113  TMOPNewtonSolver solver(pfespace->GetComm(), ir, solver_type);
1114  if (surface_fit_adapt) { solver.EnableAdaptiveSurfaceFitting(); }
1115  if (surface_fit_threshold > 0)
1116  {
1117  solver.SetTerminationWithMaxSurfaceFittingError(surface_fit_threshold);
1118  }
1119  // Provide all integration rules in case of a mixed mesh.
1120  solver.SetIntegrationRules(*irules, quad_order);
1121  if (solver_type == 0)
1122  {
1123  // Specify linear solver when we use a Newton-based solver.
1124  solver.SetPreconditioner(*S);
1125  }
1126  // For untangling, the solver will update the min det(T) values.
1127  if (tauval < 0.0) { solver.SetMinDetPtr(&tauval); }
1128  solver.SetMaxIter(solver_iter);
1129  solver.SetRelTol(solver_rtol);
1130  solver.SetAbsTol(0.0);
1131  if (solver_art_type > 0)
1132  {
1133  solver.SetAdaptiveLinRtol(solver_art_type, 0.5, 0.9);
1134  }
1135  solver.SetPrintLevel(verbosity_level >= 1 ? 1 : -1);
1136 
1137  // hr-adaptivity solver.
1138  // If hr-adaptivity is disabled, r-adaptivity is done once using the
1139  // TMOPNewtonSolver.
1140  // Otherwise, "hr_iter" iterations of r-adaptivity are done followed by
1141  // "h_per_r_iter" iterations of h-adaptivity after each r-adaptivity.
1142  // The solver terminates if an h-adaptivity iteration does not modify
1143  // any element in the mesh.
1144  TMOPHRSolver hr_solver(*pmesh, a, solver,
1145  x, move_bnd, hradaptivity,
1146  mesh_poly_deg, h_metric_id,
1147  n_hr_iter, n_h_iter);
1148  hr_solver.AddGridFunctionForUpdate(&x0);
1149  if (adapt_lim_const > 0.)
1150  {
1151  hr_solver.AddGridFunctionForUpdate(&adapt_lim_gf0);
1152  hr_solver.AddFESpaceForUpdate(&ind_fes);
1153  }
1154  hr_solver.Mult();
1155 
1156  // 16. Save the optimized mesh to a file. This output can be viewed later
1157  // using GLVis: "glvis -m optimized -np num_mpi_tasks".
1158  {
1159  ostringstream mesh_name;
1160  mesh_name << "optimized.mesh";
1161  ofstream mesh_ofs(mesh_name.str().c_str());
1162  mesh_ofs.precision(8);
1163  pmesh->PrintAsOne(mesh_ofs);
1164  }
1165 
1166  // Compute the final energy of the functional.
1167  const double fin_energy = a.GetParGridFunctionEnergy(x) /
1168  (hradaptivity ? pmesh->GetGlobalNE() : 1);
1169  double fin_metric_energy = fin_energy;
1170  if (lim_const > 0.0 || adapt_lim_const > 0.0 || surface_fit_const > 0.0)
1171  {
1172  lim_coeff.constant = 0.0;
1173  adapt_lim_coeff.constant = 0.0;
1174  surf_fit_coeff.constant = 0.0;
1175  fin_metric_energy = a.GetParGridFunctionEnergy(x) /
1176  (hradaptivity ? pmesh->GetGlobalNE() : 1);
1177  lim_coeff.constant = lim_const;
1178  adapt_lim_coeff.constant = adapt_lim_const;
1179  surf_fit_coeff.constant = surface_fit_const;
1180  }
1181  if (myid == 0)
1182  {
1183  std::cout << std::scientific << std::setprecision(4);
1184  cout << "Initial strain energy: " << init_energy
1185  << " = metrics: " << init_metric_energy
1186  << " + extra terms: " << init_energy - init_metric_energy << endl;
1187  cout << " Final strain energy: " << fin_energy
1188  << " = metrics: " << fin_metric_energy
1189  << " + extra terms: " << fin_energy - fin_metric_energy << endl;
1190  cout << "The strain energy decreased by: "
1191  << (init_energy - fin_energy) * 100.0 / init_energy << " %." << endl;
1192  }
1193 
1194  // 18. Visualize the final mesh and metric values.
1195  if (visualization)
1196  {
1197  char title[] = "Final metric values";
1198  vis_tmop_metric_p(mesh_poly_deg, *metric, *target_c, *pmesh, title, 600);
1199  }
1200 
1201  if (adapt_lim_const > 0.0 && visualization)
1202  {
1203  socketstream vis0;
1204  common::VisualizeField(vis0, "localhost", 19916, adapt_lim_gf0, "Xi 0",
1205  600, 600, 300, 300);
1206  }
1207 
1208  if (surface_fit_const > 0.0)
1209  {
1210  if (visualization)
1211  {
1212  socketstream vis2, vis3;
1213  common::VisualizeField(vis2, "localhost", 19916, mat,
1214  "Materials", 600, 900, 300, 300);
1215  common::VisualizeField(vis3, "localhost", 19916, surf_fit_mat_gf,
1216  "Surface dof", 900, 900, 300, 300);
1217  }
1218  double err_avg, err_max;
1219  tmop_integ->GetSurfaceFittingErrors(err_avg, err_max);
1220  if (myid == 0)
1221  {
1222  std::cout << "Avg fitting error: " << err_avg << std::endl
1223  << "Max fitting error: " << err_max << std::endl;
1224  }
1225  }
1226 
1227  // 19. Visualize the mesh displacement.
1228  if (visualization)
1229  {
1230  x0 -= x;
1231  socketstream sock;
1232  if (myid == 0)
1233  {
1234  sock.open("localhost", 19916);
1235  sock << "solution\n";
1236  }
1237  pmesh->PrintAsOne(sock);
1238  x0.SaveAsOne(sock);
1239  if (myid == 0)
1240  {
1241  sock << "window_title 'Displacements'\n"
1242  << "window_geometry "
1243  << 1200 << " " << 0 << " " << 600 << " " << 600 << "\n"
1244  << "keys jRmclA" << endl;
1245  }
1246  }
1247 
1248  // 20. Free the used memory.
1249  delete S;
1250  delete S_prec;
1251  delete target_c2;
1252  delete metric2;
1253  delete metric_coeff1;
1254  delete adapt_lim_eval;
1255  delete adapt_surface;
1256  delete target_c;
1257  delete hr_adapt_coeff;
1258  delete adapt_coeff;
1259  delete h_metric;
1260  delete metric;
1261  delete pfespace;
1262  delete fec;
1263  delete pmesh;
1264 
1265  return 0;
1266 }
int GetNPoints() const
Returns the number of the points in the integration rule.
Definition: intrules.hpp:247
void vis_tmop_metric_p(int order, TMOP_QualityMetric &qm, const TargetConstructor &tc, ParMesh &pmesh, char *title, int position)
Definition: tmop_tools.cpp:863
void discrete_aspr_3d(const Vector &x, Vector &v)
virtual void SetParDiscreteTargetSize(const ParGridFunction &tspec_)
Definition: tmop.cpp:1375
int Size() const
Return the logical size of the array.
Definition: array.hpp:138
static void Init()
Initialize hypre by calling HYPRE_Init() and set default options. After calling Hypre::Init(), hypre will be finalized automatically at program exit.
Definition: hypre.hpp:67
Conjugate gradient method.
Definition: solvers.hpp:465
Definition: ex25.cpp:148
int GetNDofs() const
Returns number of degrees of freedom.
Definition: fespace.hpp:560
AssemblyLevel
Enumeration defining the assembly level for bilinear and nonlinear form classes derived from Operator...
Class for an integration rule - an Array of IntegrationPoint.
Definition: intrules.hpp:90
void EnableFiniteDifferences(const GridFunction &x)
Enables FD-based approximation and computes dx.
Definition: tmop.cpp:3642
int DofToVDof(int dof, int vd, int ndofs=-1) const
Definition: fespace.cpp:251
void GetDerivative(int comp, int der_comp, ParGridFunction &der)
Parallel version of GridFunction::GetDerivative(); see its documentation.
Definition: pgridfunc.cpp:485
double discrete_size_2d(const Vector &x)
3D barrier Shape+Size (VS) metric.
Definition: tmop.hpp:578
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:923
IntegrationRules IntRulesLo(0, Quadrature1D::GaussLobatto)
void SetPositiveDiagonal(bool pos_diag=true)
Replace diagonal entries with their absolute values.
Definition: solvers.hpp:311
void SetFromTrueVector()
Shortcut for calling SetFromTrueDofs() with GetTrueVector() as argument.
Definition: gridfunc.hpp:150
3D barrier Shape+Size (VS) metric (polyconvex).
Definition: tmop.hpp:661
A coefficient that is constant across space and time.
Definition: coefficient.hpp:78
long GetGlobalNE() const
Return the total (global) number of elements.
Definition: mesh.hpp:960
int GetNBE() const
Returns number of boundary elements.
Definition: mesh.hpp:931
3D barrier Shape (S) metric.
Definition: tmop.hpp:484
void SetAdaptivityEvaluator(AdaptivityEvaluator *ae)
Definition: tmop.hpp:1233
virtual void Setup()
Setup the NonlinearForm: based on the current AssemblyLevel and the current mesh, optionally...
double Det() const
Definition: densemat.cpp:436
2D barrier Shape+Size (VS) metric (polyconvex).
Definition: tmop.hpp:358
void AddDomainIntegrator(NonlinearFormIntegrator *nlfi)
Adds new Domain Integrator.
void SetIntPoint(const IntegrationPoint *ip)
Set the integration point ip that weights and Jacobians will be evaluated at.
Definition: eltrans.hpp:93
double surface_level_set(const Vector &x)
Coefficient defined by a GridFunction. This coefficient is mesh dependent.
Container class for integration rules.
Definition: intrules.hpp:311
Data type dense matrix using column-major storage.
Definition: densemat.hpp:23
3D barrier Shape (S) metric.
Definition: tmop.hpp:466
void GetSurfaceFittingErrors(double &err_avg, double &err_max)
Definition: tmop.cpp:2436
int Size() const
Returns the size of the vector.
Definition: vector.hpp:199
void SetVolumeScale(double vol_scale)
Used by target type IDEAL_SHAPE_EQUAL_SIZE. The default volume scale is 1.
Definition: tmop.hpp:994
Parallel non-linear operator on the true dofs.
2D barrier shape (S) metric (polyconvex).
Definition: tmop.hpp:186
void EnableAdaptiveLimiting(const GridFunction &z0, Coefficient &coeff, AdaptivityEvaluator &ae)
Restriction of the node positions to certain regions.
Definition: tmop.cpp:2365
3D non-barrier Size (V) metric.
Definition: tmop.hpp:542
2D barrier size (V) metric (polyconvex).
Definition: tmop.hpp:340
void SetIntegrationRules(IntegrationRules &irules, int order)
Prescribe a set of integration rules; relevant for mixed meshes.
Definition: tmop.hpp:1556
int GetNE() const
Returns number of elements.
Definition: mesh.hpp:928
2D barrier Shape+Size (VS) metric (polyconvex).
Definition: tmop.hpp:774
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:546
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:525
void SetExactActionFlag(bool flag_)
Flag to control if exact action of Integration is effected.
Definition: tmop.hpp:1711
MINRES method.
Definition: solvers.hpp:575
void Randomize(int seed=0)
Set random values in the vector.
Definition: vector.cpp:772
virtual int OrderJ() const =0
Return the order of the elements of the Jacobian of the transformation.
2D barrier Shape+Size (VS) metric (not polyconvex).
Definition: tmop.hpp:204
Geometry::Type GetElementBaseGeometry(int i) const
Definition: mesh.hpp:1062
void EnableLimiting(const GridFunction &n0, const GridFunction &dist, Coefficient &w0, TMOP_LimiterFunction *lfunc=NULL)
Adds the limiting term to the first integrator. Disables it for the rest.
Definition: tmop.cpp:3670
2D Shifted barrier form of shape metric (mu_2).
Definition: tmop.hpp:253
const DenseMatrix & GetGeomToPerfGeomJac(int GeomType) const
Definition: geom.hpp:98
DofTransformation * GetBdrElementVDofs(int i, Array< int > &vdofs) const
Returns indexes of degrees of freedom for i&#39;th boundary element.
Definition: fespace.cpp:297
2D barrier shape (S) metric (not polyconvex).
Definition: tmop.hpp:323
Geometry::Type GetGeomType() const
Returns the Geometry::Type of the reference element.
Definition: fe_base.hpp:320
Abstract class for local mesh quality metrics in the target-matrix optimization paradigm (TMOP) by P...
Definition: tmop.hpp:23
Geometry Geometries
Definition: fe.cpp:49
IntegrationPoint & IntPoint(int i)
Returns a reference to the i-th integration point.
Definition: intrules.hpp:250
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:588
virtual void SetPrintLevel(int print_lvl)
Legacy method to set the level of verbosity of the solver output.
Definition: solvers.cpp:71
MPI_Comm GetComm() const
Definition: pfespace.hpp:273
double GetElementSize(ElementTransformation *T, int type=0)
Definition: mesh.cpp:76
virtual const FiniteElement * GetFE(int i) const
Definition: pfespace.cpp:535
void EnableLimiting(const GridFunction &n0, const GridFunction &dist, Coefficient &w0, TMOP_LimiterFunction *lfunc=NULL)
Limiting of the mesh displacements (general version).
Definition: tmop.cpp:2340
int GetNBE() const
Returns number of boundary elements in the mesh.
Definition: fespace.hpp:599
void Parse()
Parse the command-line options. Note that this function expects all the options provided through the ...
Definition: optparser.cpp:151
void SetCoefficient(Coefficient &w1)
Sets a scaling Coefficient for the quality metric term of the integrator.
Definition: tmop.hpp:1568
2D non-barrier size (V) metric (not polyconvex).
Definition: tmop.hpp:288
void EnsureNCMesh(bool simplices_nonconforming=false)
Definition: mesh.cpp:9394
3D barrier Size (V) metric.
Definition: tmop.hpp:560
Jacobi smoothing for a given bilinear form (no matrix necessary).
Definition: solvers.hpp:274
double b
Definition: lissajous.cpp:42
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:9737
const DenseMatrix & Jacobian()
Return the Jacobian matrix of the transformation at the currently set IntegrationPoint, using the method SetIntPoint().
Definition: eltrans.hpp:119
void SetNodes(const GridFunction &n)
Set the nodes to be used in the target-matrix construction.
Definition: tmop.hpp:988
void SetMaxIter(int max_it)
Definition: solvers.hpp:200
T Max() const
Find the maximal element in the array, using the comparison operator &lt; for class T.
Definition: array.cpp:68
Version of QuadraticFECollection with positive basis functions.
Definition: fe_coll.hpp:716
virtual void SetParDiscreteTargetAspectRatio(const ParGridFunction &tspec_)
Definition: tmop.cpp:1391
int GetAttribute() const
Return element&#39;s attribute.
Definition: element.hpp:55
double discrete_ori_2d(const Vector &x)
double Weight()
Return the weight of the Jacobian matrix of the transformation at the currently set IntegrationPoint...
Definition: eltrans.hpp:131
Parallel smoothers in hypre.
Definition: hypre.hpp:896
static void Init()
Singleton creation with Mpi::Init();.
void ParEnableNormalization(const ParGridFunction &x)
Definition: tmop.cpp:3785
2D barrier Shape+Size+Orientation (VOS) metric (polyconvex).
Definition: tmop.hpp:396
void SetNodalFESpace(FiniteElementSpace *nfes)
Definition: mesh.cpp:5273
void AddFESpaceForUpdate(FiniteElementSpace *fes)
Definition: tmop_amr.hpp:253
int Dimension() const
Definition: mesh.hpp:999
void SetAssemblyLevel(AssemblyLevel assembly_level)
Set the desired assembly level. The default is AssemblyLevel::LEGACY.
void PrintUsage(std::ostream &out) const
Print the usage message.
Definition: optparser.cpp:454
void DiffuseField(ParGridFunction &field, int smooth_steps)
Definition: dist_solver.cpp:17
void EnableAdaptiveSurfaceFitting()
Definition: tmop_tools.hpp:215
virtual DofTransformation * GetElementDofs(int i, Array< int > &dofs) const
Returns indexes of degrees of freedom in array dofs for i&#39;th element.
Definition: pfespace.cpp:471
2D barrier Shape+Orientation (OS) metric (polyconvex).
Definition: tmop.hpp:756
double GetParGridFunctionEnergy(const Vector &x) const
Compute the energy corresponding to the state x.
FDualNumber< tbase > pow(const FDualNumber< tbase > &a, const FDualNumber< tbase > &b)
pow([dual number],[dual number])
Definition: fdual.hpp:543
A general vector function coefficient.
3D barrier Shape+Size (VS) metric (polyconvex).
Definition: tmop.hpp:596
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 AddGridFunctionForUpdate(GridFunction *gf)
Definition: tmop_amr.hpp:252
void SetRelTol(double rtol)
Definition: solvers.hpp:198
2D barrier (not a shape) metric (polyconvex).
Definition: tmop.hpp:272
int GetDof() const
Returns the number of degrees of freedom in the finite element.
Definition: fe_base.hpp:323
double weight_fun(const Vector &x)
2D barrier size (V) metric (polyconvex).
Definition: tmop.hpp:305
Collection of finite elements from the same family in multiple dimensions. This class is used to matc...
Definition: fe_coll.hpp:26
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
double adapt_lim_fun(const Vector &x)
void AddTMOPIntegrator(TMOP_Integrator *ti)
Adds a new TMOP_Integrator to the combination.
Definition: tmop.hpp:1735
void SaveAsOne(const char *fname, int precision=16) const
Definition: pgridfunc.cpp:899
double material_indicator_2d(const Vector &x)
2D barrier Shape+Size (VS) metric (not polyconvex).
Definition: tmop.hpp:222
double a
Definition: lissajous.cpp:41
double discrete_size_3d(const Vector &x)
double discrete_aspr_2d(const Vector &x)
Class for integration point with weight.
Definition: intrules.hpp:25
void GetElementTransformation(int i, IsoparametricTransformation *ElTr)
Definition: mesh.cpp:348
2D barrier Shape+Size+Orientation (VOS) metric (polyconvex).
Definition: tmop.hpp:738
3D Shape (S) metric, untangling version of 303.
Definition: tmop.hpp:521
int dim
Definition: ex24.cpp:53
IntegrationRules IntRulesCU(0, Quadrature1D::ClosedUniform)
void PrintOptions(std::ostream &out) const
Print the options.
Definition: optparser.cpp:324
double infinity()
Define a shortcut for std::numeric_limits&lt;double&gt;::infinity()
Definition: vector.hpp:46
void SetEssentialBC(const Array< int > &bdr_attr_is_ess, Vector *rhs=NULL)
Specify essential boundary conditions.
int open(const char hostname[], int port)
Open the socket stream on &#39;port&#39; at &#39;hostname&#39;.
static int WorldRank()
Return the MPI rank in MPI_COMM_WORLD.
A general function coefficient.
void SetAttribute(int i, int attr)
Set the attribute of element i.
Definition: mesh.hpp:1458
Vector data type.
Definition: vector.hpp:60
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:94
3D barrier Shape+Size (VS) metric (polyconvex).
Definition: tmop.hpp:640
void PrintAsOne(std::ostream &out=mfem::out) const
Definition: pmesh.cpp:4829
Arbitrary order H1-conforming (continuous) finite elements.
Definition: fe_coll.hpp:216
TargetType
Target-matrix construction algorithms supported by this class.
Definition: tmop.hpp:912
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:2400
void SetEssentialVDofs(const Array< int > &ess_vdofs_list)
Specify essential boundary conditions.
virtual void SetParDiscreteTargetOrientation(const ParGridFunction &tspec_)
Definition: tmop.cpp:1399
Base class for solvers.
Definition: operator.hpp:651
Class for parallel grid function.
Definition: pgridfunc.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:3102
void SetNodalGridFunction(GridFunction *nodes, bool make_owner=false)
Definition: mesh.cpp:5301
Base class representing target-matrix construction algorithms for mesh optimization via the target-ma...
Definition: tmop.hpp:908
The MFEM Device class abstracts hardware devices such as GPUs, as well as programming models such as ...
Definition: device.hpp:121
virtual void SetAnalyticTargetSpec(Coefficient *sspec, VectorCoefficient *vspec, TMOPMatrixCoefficient *mspec)
Definition: tmop.cpp:1227
3D barrier Shape+Size (VS) metric (polyconvex).
Definition: tmop.hpp:617
int material_id(int el_id, const GridFunction &g)
Class for parallel meshes.
Definition: pmesh.hpp:32
IntegrationRules IntRules(0, Quadrature1D::GaussLegendre)
A global object with all integration rules (defined in intrules.cpp)
Definition: intrules.hpp:379
void SetType(HypreSmoother::Type type, int relax_times=1)
Set the relaxation type and number of sweeps.
Definition: hypre.cpp:3224
2D barrier Shape+Orientation (OS) metric (polyconvex).
Definition: tmop.hpp:381
virtual double * HostReadWrite()
Shortcut for mfem::ReadWrite(vec.GetMemory(), vec.Size(), false).
Definition: vector.hpp:458
void ParEnableNormalization(const ParGridFunction &x)
Definition: tmop.cpp:3490
int main()
double GetElementVolume(int i)
Definition: mesh.cpp:112
const Element * GetBdrElement(int i) const
Definition: mesh.hpp:1037
3D barrier Shape (S) metric.
Definition: tmop.hpp:450
2D non-barrier metric without a type.
Definition: tmop.hpp:108
Arbitrary order &quot;L2-conforming&quot; discontinuous finite elements.
Definition: fe_coll.hpp:284
bool Good() const
Return true if the command line options were parsed successfully.
Definition: optparser.hpp:150
A TMOP integrator class based on any given TMOP_QualityMetric and TargetConstructor.
Definition: tmop.hpp:1303
2D non-barrier Shape+Size+Orientation (VOS) metric (polyconvex).
Definition: tmop.hpp:238