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