MFEM  v4.6.0
Finite element discretization library
mesh-optimizer.cpp
Go to the documentation of this file.
1 // Copyright (c) 2010-2023, Lawrence Livermore National Security, LLC. Produced
2 // at the Lawrence Livermore National Laboratory. All Rights reserved. See files
3 // LICENSE and NOTICE for details. LLNL-CODE-806117.
4 //
5 // This file is part of the MFEM library. For more information and source code
6 // availability visit https://mfem.org.
7 //
8 // MFEM is free software; you can redistribute it and/or modify it under the
9 // terms of the BSD-3 license. We welcome feedback and contributions, see file
10 // CONTRIBUTING.md for details.
11 //
12 // --------------------------------------------------
13 // Mesh Optimizer Miniapp: Optimize high-order meshes
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 mesh-optimizer
33 //
34 // Sample runs:
35 // Adapted analytic shape:
36 // mesh-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 // mesh-optimizer -m square01.mesh -o 2 -rs 2 -mid 14 -tid 4 -ni 100 -bnd -qt 1 -qo 8 -fd
39 // Adapted analytic shape+orientation:
40 // mesh-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 // mesh-optimizer -m square01.mesh -o 2 -tid 9 -ni 50 -li 20 -hmid 55 -mid 7 -hr
44 // mesh-optimizer -m square01.mesh -o 2 -tid 10 -ni 50 -li 20 -hmid 55 -mid 7 -hr
45 // mesh-optimizer -m square01.mesh -o 2 -tid 11 -ni 50 -li 20 -hmid 58 -mid 7 -hr
46 //
47 // Adapted discrete size:
48 // mesh-optimizer -m square01.mesh -o 2 -rs 2 -mid 94 -tid 5 -ni 50 -qo 4 -nor
49 // (requires GSLIB):
50 // * mesh-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 NC mesh;
52 // mesh-optimizer -m amr-quad-q2.mesh -o 2 -rs 2 -mid 94 -tid 5 -ni 50 -qo 4 -nor
53 // Adapted discrete size 3D with PA:
54 // mesh-optimizer -m cube.mesh -o 2 -rs 2 -mid 321 -tid 5 -ls 3 -nor -pa
55 // Adapted discrete size 3D with PA on device (requires CUDA):
56 // * mesh-optimizer -m cube.mesh -o 3 -rs 3 -mid 321 -tid 5 -ls 3 -nor -lc 0.1 -pa -d cuda
57 // Adapted discrete size; explicit combo of metrics; mixed tri/quad mesh:
58 // mesh-optimizer -m ../../data/square-mixed.mesh -o 2 -rs 2 -mid 2 -tid 5 -ni 200 -bnd -qo 6 -cmb 2 -nor
59 // Adapted discrete size+aspect_ratio:
60 // mesh-optimizer -m square01.mesh -o 2 -rs 2 -mid 7 -tid 6 -ni 100
61 // mesh-optimizer -m square01.mesh -o 2 -rs 2 -mid 7 -tid 6 -ni 100 -qo 6 -ex -st 1 -nor
62 // Adapted discrete size+orientation:
63 // mesh-optimizer -m square01.mesh -o 2 -rs 2 -mid 36 -tid 8 -qo 4 -fd -nor
64 // Adapted discrete aspect ratio (3D):
65 // mesh-optimizer -m cube.mesh -o 2 -rs 2 -mid 302 -tid 7 -ni 20 -bnd -qt 1 -qo 8
66 //
67 // Adaptive limiting:
68 // mesh-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 // mesh-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 // * mesh-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 // Blade shape:
75 // mesh-optimizer -m blade.mesh -o 4 -mid 2 -tid 1 -ni 30 -ls 3 -art 1 -bnd -qt 1 -qo 8
76 // (requires CUDA):
77 // * mesh-optimizer -m blade.mesh -o 4 -mid 2 -tid 1 -ni 30 -ls 3 -art 1 -bnd -qt 1 -qo 8 -d cuda
78 // Blade shape with FD-based solver:
79 // mesh-optimizer -m blade.mesh -o 4 -mid 2 -tid 1 -ni 30 -ls 4 -bnd -qt 1 -qo 8 -fd
80 // Blade limited shape:
81 // mesh-optimizer -m blade.mesh -o 4 -mid 2 -tid 1 -bnd -qt 1 -qo 8 -lc 5000
82 // ICF shape and equal size:
83 // mesh-optimizer -o 3 -mid 80 -bec -tid 2 -ni 25 -ls 3 -art 2 -qo 5
84 // ICF shape and initial size:
85 // mesh-optimizer -o 3 -mid 9 -tid 3 -ni 30 -ls 3 -bnd -qt 1 -qo 8
86 // ICF shape:
87 // mesh-optimizer -o 3 -mid 1 -tid 1 -ni 100 -bnd -qt 1 -qo 8
88 // ICF limited shape:
89 // mesh-optimizer -o 3 -mid 1 -tid 1 -ni 100 -bnd -qt 1 -qo 8 -lc 10
90 // ICF combo shape + size (rings, slow convergence):
91 // mesh-optimizer -o 3 -mid 1 -tid 1 -ni 1000 -bnd -qt 1 -qo 8 -cmb 1
92 // Mixed tet / cube / hex mesh with limiting:
93 // mesh-optimizer -m ../../data/fichera-mixed-p2.mesh -o 4 -rs 1 -mid 301 -tid 1 -fix-bnd -qo 6 -nor -lc 0.25
94 // 3D pinched sphere shape (the mesh is in the mfem/data GitHub repository):
95 // * mesh-optimizer -m ../../../mfem_data/ball-pert.mesh -o 4 -mid 303 -tid 1 -ni 20 -li 500 -fix-bnd
96 // 2D non-conforming shape and equal size:
97 // mesh-optimizer -m ./amr-quad-q2.mesh -o 2 -rs 1 -mid 9 -tid 2 -ni 200 -bnd -qt 1 -qo 8
98 //
99 // 2D untangling:
100 // mesh-optimizer -m jagged.mesh -o 2 -mid 22 -tid 1 -ni 50 -li 50 -qo 4 -fd -vl 1
101 // 2D untangling with shifted barrier metric:
102 // mesh-optimizer -m jagged.mesh -o 2 -mid 4 -tid 1 -ni 50 -qo 4 -fd -vl 1 -btype 1
103 // 3D untangling (the mesh is in the mfem/data GitHub repository):
104 // * mesh-optimizer -m ../../../mfem_data/cube-holes-inv.mesh -o 3 -mid 313 -tid 1 -rtol 1e-5 -li 50 -qo 4 -fd -vl 1
105 
106 #include "mfem.hpp"
107 #include "../common/mfem-common.hpp"
108 #include <fstream>
109 #include <iostream>
110 #include "mesh-optimizer.hpp"
111 
112 using namespace mfem;
113 using namespace std;
114 
115 int main(int argc, char *argv[])
116 {
117  // 0. Set the method's default parameters.
118  const char *mesh_file = "icf.mesh";
119  int mesh_poly_deg = 1;
120  int rs_levels = 0;
121  double jitter = 0.0;
122  int metric_id = 1;
123  int target_id = 1;
124  double lim_const = 0.0;
125  double adapt_lim_const = 0.0;
126  int quad_type = 1;
127  int quad_order = 8;
128  int solver_type = 0;
129  int solver_iter = 20;
130  double solver_rtol = 1e-10;
131  int solver_art_type = 0;
132  int lin_solver = 2;
133  int max_lin_iter = 100;
134  bool move_bnd = true;
135  int combomet = 0;
136  bool bal_expl_combo = false;
137  bool hradaptivity = false;
138  int h_metric_id = -1;
139  bool normalization = false;
140  bool visualization = true;
141  int verbosity_level = 0;
142  bool fdscheme = false;
143  int adapt_eval = 0;
144  bool exactaction = false;
145  bool integ_over_targ = true;
146  const char *devopt = "cpu";
147  bool pa = false;
148  int n_hr_iter = 5;
149  int n_h_iter = 1;
150  int mesh_node_ordering = 0;
151  int barrier_type = 0;
152  int worst_case_type = 0;
153 
154  // 1. Parse command-line options.
155  OptionsParser args(argc, argv);
156  args.AddOption(&mesh_file, "-m", "--mesh",
157  "Mesh file to use.");
158  args.AddOption(&mesh_poly_deg, "-o", "--order",
159  "Polynomial degree of mesh finite element space.");
160  args.AddOption(&rs_levels, "-rs", "--refine-serial",
161  "Number of times to refine the mesh uniformly in serial.");
162  args.AddOption(&jitter, "-ji", "--jitter",
163  "Random perturbation scaling factor.");
164  args.AddOption(&metric_id, "-mid", "--metric-id",
165  "Mesh optimization metric:\n\t"
166  "T-metrics\n\t"
167  "1 : |T|^2 -- 2D no type\n\t"
168  "2 : 0.5|T|^2/tau-1 -- 2D shape (condition number)\n\t"
169  "7 : |T-T^-t|^2 -- 2D shape+size\n\t"
170  "9 : tau*|T-T^-t|^2 -- 2D shape+size\n\t"
171  "14 : |T-I|^2 -- 2D shape+size+orientation\n\t"
172  "22 : 0.5(|T|^2-2*tau)/(tau-tau_0) -- 2D untangling\n\t"
173  "50 : 0.5|T^tT|^2/tau^2-1 -- 2D shape\n\t"
174  "55 : (tau-1)^2 -- 2D size\n\t"
175  "56 : 0.5(sqrt(tau)-1/sqrt(tau))^2 -- 2D size\n\t"
176  "58 : |T^tT|^2/(tau^2)-2*|T|^2/tau+2 -- 2D shape\n\t"
177  "77 : 0.5(tau-1/tau)^2 -- 2D size\n\t"
178  "80 : (1-gamma)mu_2 + gamma mu_77 -- 2D shape+size\n\t"
179  "85 : |T-|T|/sqrt(2)I|^2 -- 2D shape+orientation\n\t"
180  "90 : balanced combo mu_50 & mu_77 -- 2D shape+size\n\t"
181  "94 : balanced combo mu_2 & mu_56 -- 2D shape+size\n\t"
182  "98 : (1/tau)|T-I|^2 -- 2D shape+size+orientation\n\t"
183  // "211: (tau-1)^2-tau+sqrt(tau^2+eps) -- 2D untangling\n\t"
184  // "211: (tau-1)^2-tau+sqrt(tau^2) -- 2D untangling\n\t"
185  // "252: 0.5(tau-1)^2/(tau-tau_0) -- 2D untangling\n\t"
186  "301: (|T||T^-1|)/3-1 -- 3D shape\n\t"
187  "302: (|T|^2|T^-1|^2)/9-1 -- 3D shape\n\t"
188  "303: (|T|^2)/3/tau^(2/3)-1 -- 3D shape\n\t"
189  "304: (|T|^3)/3^{3/2}/tau-1 -- 3D shape\n\t"
190  //"311: (tau-1)^2-tau+sqrt(tau^2+eps)-- 3D untangling\n\t"
191  "313: (|T|^2)(tau-tau0)^(-2/3)/3 -- 3D untangling\n\t"
192  "315: (tau-1)^2 -- 3D no type\n\t"
193  "316: 0.5(sqrt(tau)-1/sqrt(tau))^2 -- 3D no type\n\t"
194  "321: |T-T^-t|^2 -- 3D shape+size\n\t"
195  "322: |T-adjT^-t|^2 -- 3D shape+size\n\t"
196  "323: |J|^3-3sqrt(3)ln(det(J))-3sqrt(3) -- 3D shape+size\n\t"
197  "328: balanced combo mu_301 & mu_316 -- 3D shape+size\n\t"
198  "332: (1-gamma) mu_302 + gamma mu_315 -- 3D shape+size\n\t"
199  "333: (1-gamma) mu_302 + gamma mu_316 -- 3D shape+size\n\t"
200  "334: (1-gamma) mu_303 + gamma mu_316 -- 3D shape+size\n\t"
201  "328: balanced combo mu_302 & mu_318 -- 3D shape+size\n\t"
202  "347: (1-gamma) mu_304 + gamma mu_316 -- 3D shape+size\n\t"
203  // "352: 0.5(tau-1)^2/(tau-tau_0) -- 3D untangling\n\t"
204  "360: (|T|^3)/3^{3/2}-tau -- 3D shape\n\t"
205  "A-metrics\n\t"
206  "11 : (1/4*alpha)|A-(adjA)^T(W^TW)/omega|^2 -- 2D shape\n\t"
207  "36 : (1/alpha)|A-W|^2 -- 2D shape+size+orientation\n\t"
208  "107: (1/2*alpha)|A-|A|/|W|W|^2 -- 2D shape+orientation\n\t"
209  "126: (1-gamma)nu_11 + gamma*nu_14a -- 2D shape+size\n\t"
210  );
211  args.AddOption(&target_id, "-tid", "--target-id",
212  "Target (ideal element) type:\n\t"
213  "1: Ideal shape, unit size\n\t"
214  "2: Ideal shape, equal size\n\t"
215  "3: Ideal shape, initial size\n\t"
216  "4: Given full analytic Jacobian (in physical space)\n\t"
217  "5: Ideal shape, given size (in physical space)");
218  args.AddOption(&lim_const, "-lc", "--limit-const", "Limiting constant.");
219  args.AddOption(&adapt_lim_const, "-alc", "--adapt-limit-const",
220  "Adaptive limiting coefficient constant.");
221  args.AddOption(&quad_type, "-qt", "--quad-type",
222  "Quadrature rule type:\n\t"
223  "1: Gauss-Lobatto\n\t"
224  "2: Gauss-Legendre\n\t"
225  "3: Closed uniform points");
226  args.AddOption(&quad_order, "-qo", "--quad_order",
227  "Order of the quadrature rule.");
228  args.AddOption(&solver_type, "-st", "--solver-type",
229  " Type of solver: (default) 0: Newton, 1: LBFGS");
230  args.AddOption(&solver_iter, "-ni", "--newton-iters",
231  "Maximum number of Newton iterations.");
232  args.AddOption(&solver_rtol, "-rtol", "--newton-rel-tolerance",
233  "Relative tolerance for the Newton solver.");
234  args.AddOption(&solver_art_type, "-art", "--adaptive-rel-tol",
235  "Type of adaptive relative linear solver tolerance:\n\t"
236  "0: None (default)\n\t"
237  "1: Eisenstat-Walker type 1\n\t"
238  "2: Eisenstat-Walker type 2");
239  args.AddOption(&lin_solver, "-ls", "--lin-solver",
240  "Linear solver:\n\t"
241  "0: l1-Jacobi\n\t"
242  "1: CG\n\t"
243  "2: MINRES\n\t"
244  "3: MINRES + Jacobi preconditioner\n\t"
245  "4: MINRES + l1-Jacobi preconditioner");
246  args.AddOption(&max_lin_iter, "-li", "--lin-iter",
247  "Maximum number of iterations in the linear solve.");
248  args.AddOption(&move_bnd, "-bnd", "--move-boundary", "-fix-bnd",
249  "--fix-boundary",
250  "Enable motion along horizontal and vertical boundaries.");
251  args.AddOption(&combomet, "-cmb", "--combo-type",
252  "Combination of metrics options:\n\t"
253  "0: Use single metric\n\t"
254  "1: Shape + space-dependent size given analytically\n\t"
255  "2: Shape + adapted size given discretely; shared target");
256  args.AddOption(&bal_expl_combo, "-bec", "--balance-explicit-combo",
257  "-no-bec", "--balance-explicit-combo",
258  "Automatic balancing of explicit combo metrics.");
259  args.AddOption(&hradaptivity, "-hr", "--hr-adaptivity", "-no-hr",
260  "--no-hr-adaptivity",
261  "Enable hr-adaptivity.");
262  args.AddOption(&h_metric_id, "-hmid", "--h-metric",
263  "Same options as metric_id. Used to determine refinement"
264  " type for each element if h-adaptivity is enabled.");
265  args.AddOption(&normalization, "-nor", "--normalization", "-no-nor",
266  "--no-normalization",
267  "Make all terms in the optimization functional unitless.");
268  args.AddOption(&fdscheme, "-fd", "--fd_approximation",
269  "-no-fd", "--no-fd-approx",
270  "Enable finite difference based derivative computations.");
271  args.AddOption(&exactaction, "-ex", "--exact_action",
272  "-no-ex", "--no-exact-action",
273  "Enable exact action of TMOP_Integrator.");
274  args.AddOption(&integ_over_targ, "-it", "--integrate-target",
275  "-ir", "--integrate-reference",
276  "Integrate over target (-it) or reference (-ir) element.");
277  args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
278  "--no-visualization",
279  "Enable or disable GLVis visualization.");
280  args.AddOption(&verbosity_level, "-vl", "--verbosity-level",
281  "Verbosity level for the involved iterative solvers:\n\t"
282  "0: no output\n\t"
283  "1: Newton iterations\n\t"
284  "2: Newton iterations + linear solver summaries\n\t"
285  "3: newton iterations + linear solver iterations");
286  args.AddOption(&adapt_eval, "-ae", "--adaptivity-evaluator",
287  "0 - Advection based (DEFAULT), 1 - GSLIB.");
288  args.AddOption(&devopt, "-d", "--device",
289  "Device configuration string, see Device::Configure().");
290  args.AddOption(&pa, "-pa", "--partial-assembly", "-no-pa",
291  "--no-partial-assembly", "Enable Partial Assembly.");
292  args.AddOption(&n_hr_iter, "-nhr", "--n_hr_iter",
293  "Number of hr-adaptivity iterations.");
294  args.AddOption(&n_h_iter, "-nh", "--n_h_iter",
295  "Number of h-adaptivity iterations per r-adaptivity"
296  "iteration.");
297  args.AddOption(&mesh_node_ordering, "-mno", "--mesh_node_ordering",
298  "Ordering of mesh nodes."
299  "0 (default): byNodes, 1: byVDIM");
300  args.AddOption(&barrier_type, "-btype", "--barrier-type",
301  "0 - None,"
302  "1 - Shifted Barrier,"
303  "2 - Pseudo Barrier.");
304  args.AddOption(&worst_case_type, "-wctype", "--worst-case-type",
305  "0 - None,"
306  "1 - Beta,"
307  "2 - PMean.");
308  args.Parse();
309  if (!args.Good())
310  {
311  args.PrintUsage(cout);
312  return 1;
313  }
314  args.PrintOptions(cout);
315 
316  if (h_metric_id < 0) { h_metric_id = metric_id; }
317 
318  if (hradaptivity)
319  {
320  MFEM_VERIFY(strcmp(devopt,"cpu")==0, "HR-adaptivity is currently only"
321  " supported on cpus.");
322  }
323  Device device(devopt);
324  device.Print();
325 
326  // 2. Initialize and refine the starting mesh.
327  Mesh *mesh = new Mesh(mesh_file, 1, 1, false);
328  for (int lev = 0; lev < rs_levels; lev++) { mesh->UniformRefinement(); }
329  const int dim = mesh->Dimension();
330 
331  if (hradaptivity) { mesh->EnsureNCMesh(); }
332 
333  // 3. Define a finite element space on the mesh-> Here we use vector finite
334  // elements which are tensor products of quadratic finite elements. The
335  // number of components in the vector finite element space is specified by
336  // the last parameter of the FiniteElementSpace constructor.
338  if (mesh_poly_deg <= 0)
339  {
340  fec = new QuadraticPosFECollection;
341  mesh_poly_deg = 2;
342  }
343  else { fec = new H1_FECollection(mesh_poly_deg, dim); }
344  FiniteElementSpace *fespace = new FiniteElementSpace(mesh, fec, dim,
345  mesh_node_ordering);
346 
347  // 4. Make the mesh curved based on the above finite element space. This
348  // means that we define the mesh elements through a fespace-based
349  // transformation of the reference element.
350  mesh->SetNodalFESpace(fespace);
351 
352  // 5. Set up an empty right-hand side vector b, which is equivalent to b=0.
353  Vector b(0);
354 
355  // 6. Get the mesh nodes (vertices and other degrees of freedom in the finite
356  // element space) as a finite element grid function in fespace. Note that
357  // changing x automatically changes the shapes of the mesh elements.
358  GridFunction x(fespace);
359  mesh->SetNodalGridFunction(&x);
360 
361  // 7. Define a vector representing the minimal local mesh size in the mesh
362  // nodes. We index the nodes using the scalar version of the degrees of
363  // freedom in fespace. Note: this is partition-dependent.
364  //
365  // In addition, compute average mesh size and total volume.
366  Vector h0(fespace->GetNDofs());
367  h0 = infinity();
368  double mesh_volume = 0.0;
369  Array<int> dofs;
370  for (int i = 0; i < mesh->GetNE(); i++)
371  {
372  // Get the local scalar element degrees of freedom in dofs.
373  fespace->GetElementDofs(i, dofs);
374  // Adjust the value of h0 in dofs based on the local mesh size.
375  const double hi = mesh->GetElementSize(i);
376  for (int j = 0; j < dofs.Size(); j++)
377  {
378  h0(dofs[j]) = min(h0(dofs[j]), hi);
379  }
380  mesh_volume += mesh->GetElementVolume(i);
381  }
382  const double small_phys_size = pow(mesh_volume, 1.0 / dim) / 100.0;
383 
384  // 8. Add a random perturbation to the nodes in the interior of the domain.
385  // We define a random grid function of fespace and make sure that it is
386  // zero on the boundary and its values are locally of the order of h0.
387  // The latter is based on the DofToVDof() method which maps the scalar to
388  // the vector degrees of freedom in fespace.
389  GridFunction rdm(fespace);
390  rdm.Randomize();
391  rdm -= 0.25; // Shift to random values in [-0.5,0.5].
392  rdm *= jitter;
393  rdm.HostReadWrite();
394  // Scale the random values to be of order of the local mesh size.
395  for (int i = 0; i < fespace->GetNDofs(); i++)
396  {
397  for (int d = 0; d < dim; d++)
398  {
399  rdm(fespace->DofToVDof(i,d)) *= h0(i);
400  }
401  }
402  Array<int> vdofs;
403  for (int i = 0; i < fespace->GetNBE(); i++)
404  {
405  // Get the vector degrees of freedom in the boundary element.
406  fespace->GetBdrElementVDofs(i, vdofs);
407  // Set the boundary values to zero.
408  for (int j = 0; j < vdofs.Size(); j++) { rdm(vdofs[j]) = 0.0; }
409  }
410  x -= rdm;
411  x.SetTrueVector();
412  x.SetFromTrueVector();
413 
414  // 9. Save the starting (prior to the optimization) mesh to a file. This
415  // output can be viewed later using GLVis: "glvis -m perturbed.mesh".
416  {
417  ofstream mesh_ofs("perturbed.mesh");
418  mesh->Print(mesh_ofs);
419  }
420 
421  // 10. Store the starting (prior to the optimization) positions.
422  GridFunction x0(fespace);
423  x0 = x;
424 
425  // 11. Form the integrator that uses the chosen metric and target.
426  double min_detJ = -0.1;
427  TMOP_QualityMetric *metric = NULL;
428  switch (metric_id)
429  {
430  // T-metrics
431  case 1: metric = new TMOP_Metric_001; break;
432  case 2: metric = new TMOP_Metric_002; break;
433  case 4: metric = new TMOP_Metric_004; break;
434  case 7: metric = new TMOP_Metric_007; break;
435  case 9: metric = new TMOP_Metric_009; break;
436  case 14: metric = new TMOP_Metric_014; break;
437  case 22: metric = new TMOP_Metric_022(min_detJ); break;
438  case 50: metric = new TMOP_Metric_050; break;
439  case 55: metric = new TMOP_Metric_055; break;
440  case 56: metric = new TMOP_Metric_056; break;
441  case 58: metric = new TMOP_Metric_058; break;
442  case 66: metric = new TMOP_Metric_066(0.5); break;
443  case 77: metric = new TMOP_Metric_077; break;
444  case 80: metric = new TMOP_Metric_080(0.5); break;
445  case 85: metric = new TMOP_Metric_085; break;
446  case 90: metric = new TMOP_Metric_090; break;
447  case 94: metric = new TMOP_Metric_094; break;
448  case 98: metric = new TMOP_Metric_098; break;
449  // case 211: metric = new TMOP_Metric_211; break;
450  // case 252: metric = new TMOP_Metric_252(min_detJ); break;
451  case 301: metric = new TMOP_Metric_301; break;
452  case 302: metric = new TMOP_Metric_302; break;
453  case 303: metric = new TMOP_Metric_303; break;
454  case 304: metric = new TMOP_Metric_304; break;
455  // case 311: metric = new TMOP_Metric_311; break;
456  case 313: metric = new TMOP_Metric_313(min_detJ); break;
457  case 315: metric = new TMOP_Metric_315; break;
458  case 316: metric = new TMOP_Metric_316; break;
459  case 321: metric = new TMOP_Metric_321; break;
460  case 322: metric = new TMOP_Metric_322; break;
461  case 323: metric = new TMOP_Metric_323; break;
462  case 328: metric = new TMOP_Metric_328; break;
463  case 332: metric = new TMOP_Metric_332(0.5); break;
464  case 333: metric = new TMOP_Metric_333(0.5); break;
465  case 334: metric = new TMOP_Metric_334(0.5); break;
466  case 338: metric = new TMOP_Metric_338; break;
467  case 347: metric = new TMOP_Metric_347(0.5); break;
468  // case 352: metric = new TMOP_Metric_352(min_detJ); break;
469  case 360: metric = new TMOP_Metric_360; break;
470  // A-metrics
471  case 11: metric = new TMOP_AMetric_011; break;
472  case 36: metric = new TMOP_AMetric_036; break;
473  case 107: metric = new TMOP_AMetric_107a; break;
474  case 126: metric = new TMOP_AMetric_126(0.9); break;
475  default:
476  cout << "Unknown metric_id: " << metric_id << endl;
477  return 3;
478  }
479  TMOP_QualityMetric *h_metric = NULL;
480  if (hradaptivity)
481  {
482  switch (h_metric_id)
483  {
484  case 1: h_metric = new TMOP_Metric_001; break;
485  case 2: h_metric = new TMOP_Metric_002; break;
486  case 7: h_metric = new TMOP_Metric_007; break;
487  case 9: h_metric = new TMOP_Metric_009; break;
488  case 55: h_metric = new TMOP_Metric_055; break;
489  case 56: h_metric = new TMOP_Metric_056; break;
490  case 58: h_metric = new TMOP_Metric_058; break;
491  case 77: h_metric = new TMOP_Metric_077; break;
492  case 315: h_metric = new TMOP_Metric_315; break;
493  case 316: h_metric = new TMOP_Metric_316; break;
494  case 321: h_metric = new TMOP_Metric_321; break;
495  default: cout << "Metric_id not supported for h-adaptivity: " << h_metric_id <<
496  endl;
497  return 3;
498  }
499  }
500 
502  switch (barrier_type)
503  {
505  break;
507  break;
509  break;
510  default: cout << "barrier_type not supported: " << barrier_type << endl;
511  return 3;
512  }
513 
515  switch (worst_case_type)
516  {
518  break;
520  break;
522  break;
523  default: cout << "worst_case_type not supported: " << worst_case_type << endl;
524  return 3;
525  }
526 
527  TMOP_QualityMetric *untangler_metric = NULL;
528  if (barrier_type > 0 || worst_case_type > 0)
529  {
530  if (barrier_type > 0)
531  {
532  MFEM_VERIFY(metric_id == 4 || metric_id == 14 || metric_id == 66,
533  "Metric not supported for shifted/pseudo barriers.");
534  }
535  untangler_metric = new TMOP_WorstCaseUntangleOptimizer_Metric(*metric,
536  2,
537  1.5,
538  0.001, //0.01 for pseudo barrier
539  0.001,
540  btype,
541  wctype);
542  }
543 
544  if (metric_id < 300 || h_metric_id < 300)
545  {
546  MFEM_VERIFY(dim == 2, "Incompatible metric for 3D meshes");
547  }
548  if (metric_id >= 300 || h_metric_id >= 300)
549  {
550  MFEM_VERIFY(dim == 3, "Incompatible metric for 2D meshes");
551  }
552 
554  TargetConstructor *target_c = NULL;
555  HessianCoefficient *adapt_coeff = NULL;
556  HRHessianCoefficient *hr_adapt_coeff = NULL;
557  H1_FECollection ind_fec(mesh_poly_deg, dim);
558  FiniteElementSpace ind_fes(mesh, &ind_fec);
559  FiniteElementSpace ind_fesv(mesh, &ind_fec, dim);
560  GridFunction size(&ind_fes), aspr(&ind_fes), ori(&ind_fes);
561  GridFunction aspr3d(&ind_fesv);
562 
563  const AssemblyLevel al =
565 
566  switch (target_id)
567  {
568  case 1: target_t = TargetConstructor::IDEAL_SHAPE_UNIT_SIZE; break;
569  case 2: target_t = TargetConstructor::IDEAL_SHAPE_EQUAL_SIZE; break;
570  case 3: target_t = TargetConstructor::IDEAL_SHAPE_GIVEN_SIZE; break;
571  case 4: // Analytic
572  {
574  AnalyticAdaptTC *tc = new AnalyticAdaptTC(target_t);
575  adapt_coeff = new HessianCoefficient(dim, metric_id);
576  tc->SetAnalyticTargetSpec(NULL, NULL, adapt_coeff);
577  target_c = tc;
578  break;
579  }
580  case 5: // Discrete size 2D or 3D
581  {
583  DiscreteAdaptTC *tc = new DiscreteAdaptTC(target_t);
584  if (adapt_eval == 0)
585  {
586  tc->SetAdaptivityEvaluator(new AdvectorCG(al));
587  }
588  else
589  {
590 #ifdef MFEM_USE_GSLIB
592 #else
593  MFEM_ABORT("MFEM is not built with GSLIB.");
594 #endif
595  }
596  ConstructSizeGF(size);
597  tc->SetSerialDiscreteTargetSize(size);
598  target_c = tc;
599  break;
600  }
601  case 6: // Discrete size + aspect ratio - 2D
602  {
603  GridFunction d_x(&ind_fes), d_y(&ind_fes), disc(&ind_fes);
604 
606  DiscreteAdaptTC *tc = new DiscreteAdaptTC(target_t);
608  disc.ProjectCoefficient(mat_coeff);
609  if (adapt_eval == 0)
610  {
611  tc->SetAdaptivityEvaluator(new AdvectorCG(al));
612  }
613  else
614  {
615 #ifdef MFEM_USE_GSLIB
617 #else
618  MFEM_ABORT("MFEM is not built with GSLIB.");
619 #endif
620  }
621 
622  // Diffuse the interface
623  DiffuseField(disc,2);
624 
625  // Get partials with respect to x and y of the grid function
626  disc.GetDerivative(1,0,d_x);
627  disc.GetDerivative(1,1,d_y);
628 
629  // Compute the squared magnitude of the gradient
630  for (int i = 0; i < size.Size(); i++)
631  {
632  size(i) = std::pow(d_x(i),2)+std::pow(d_y(i),2);
633  }
634  const double max = size.Max();
635 
636  for (int i = 0; i < d_x.Size(); i++)
637  {
638  d_x(i) = std::abs(d_x(i));
639  d_y(i) = std::abs(d_y(i));
640  }
641  const double eps = 0.01;
642  const double aspr_ratio = 20.0;
643  const double size_ratio = 40.0;
644 
645  for (int i = 0; i < size.Size(); i++)
646  {
647  size(i) = (size(i)/max);
648  aspr(i) = (d_x(i)+eps)/(d_y(i)+eps);
649  aspr(i) = 0.1 + 0.9*(1-size(i))*(1-size(i));
650  if (aspr(i) > aspr_ratio) {aspr(i) = aspr_ratio;}
651  if (aspr(i) < 1.0/aspr_ratio) {aspr(i) = 1.0/aspr_ratio;}
652  }
653  Vector vals;
654  const int NE = mesh->GetNE();
655  double volume = 0.0, volume_ind = 0.0;
656 
657  for (int i = 0; i < NE; i++)
658  {
660  const IntegrationRule &ir =
661  IntRules.Get(mesh->GetElementBaseGeometry(i), Tr->OrderJ());
662  size.GetValues(i, ir, vals);
663  for (int j = 0; j < ir.GetNPoints(); j++)
664  {
665  const IntegrationPoint &ip = ir.IntPoint(j);
666  Tr->SetIntPoint(&ip);
667  volume += ip.weight * Tr->Weight();
668  volume_ind += vals(j) * ip.weight * Tr->Weight();
669  }
670  }
671 
672  const double avg_zone_size = volume / NE;
673 
674  const double small_avg_ratio = (volume_ind + (volume - volume_ind) /
675  size_ratio) /
676  volume;
677 
678  const double small_zone_size = small_avg_ratio * avg_zone_size;
679  const double big_zone_size = size_ratio * small_zone_size;
680 
681  for (int i = 0; i < size.Size(); i++)
682  {
683  const double val = size(i);
684  const double a = (big_zone_size - small_zone_size) / small_zone_size;
685  size(i) = big_zone_size / (1.0+a*val);
686  }
687 
688  DiffuseField(size, 2);
689  DiffuseField(aspr, 2);
690 
691  tc->SetSerialDiscreteTargetSize(size);
693  target_c = tc;
694  break;
695  }
696  case 7: // Discrete aspect ratio 3D
697  {
699  DiscreteAdaptTC *tc = new DiscreteAdaptTC(target_t);
700  if (adapt_eval == 0)
701  {
702  tc->SetAdaptivityEvaluator(new AdvectorCG(al));
703  }
704  else
705  {
706 #ifdef MFEM_USE_GSLIB
708 #else
709  MFEM_ABORT("MFEM is not built with GSLIB.");
710 #endif
711  }
713  aspr3d.ProjectCoefficient(fd_aspr3d);
714 
716  target_c = tc;
717  break;
718  }
719  case 8: // shape/size + orientation 2D
720  {
722  DiscreteAdaptTC *tc = new DiscreteAdaptTC(target_t);
723  if (adapt_eval == 0)
724  {
725  tc->SetAdaptivityEvaluator(new AdvectorCG(al));
726  }
727  else
728  {
729 #ifdef MFEM_USE_GSLIB
731 #else
732  MFEM_ABORT("MFEM is not built with GSLIB.");
733 #endif
734  }
735 
736  ConstantCoefficient size_coeff(0.1*0.1);
737  size.ProjectCoefficient(size_coeff);
738  tc->SetSerialDiscreteTargetSize(size);
739 
741  ori.ProjectCoefficient(ori_coeff);
743  target_c = tc;
744  break;
745  }
746  // Targets used for hr-adaptivity tests.
747  case 9: // size target in an annular region.
748  case 10: // size+aspect-ratio in an annular region.
749  case 11: // size+aspect-ratio target for a rotate sine wave
750  {
752  AnalyticAdaptTC *tc = new AnalyticAdaptTC(target_t);
753  hr_adapt_coeff = new HRHessianCoefficient(dim, target_id - 9);
754  tc->SetAnalyticTargetSpec(NULL, NULL, hr_adapt_coeff);
755  target_c = tc;
756  break;
757  }
758  default: cout << "Unknown target_id: " << target_id << endl; return 3;
759  }
760  if (target_c == NULL)
761  {
762  target_c = new TargetConstructor(target_t);
763  }
764  target_c->SetNodes(x0);
765 
766  // Automatically balanced gamma in composite metrics.
767  auto metric_combo = dynamic_cast<TMOP_Combo_QualityMetric *>(metric);
768  if (metric_combo && bal_expl_combo)
769  {
770  Vector bal_weights;
771  metric_combo->ComputeBalancedWeights(x, *target_c, bal_weights);
772  metric_combo->SetWeights(bal_weights);
773  }
774 
775  TMOP_QualityMetric *metric_to_use = barrier_type > 0 || worst_case_type > 0
776  ? untangler_metric
777  : metric;
778  auto tmop_integ = new TMOP_Integrator(metric_to_use, target_c, h_metric);
779  tmop_integ->IntegrateOverTarget(integ_over_targ);
780  if (barrier_type > 0 || worst_case_type > 0)
781  {
782  tmop_integ->ComputeUntangleMetricQuantiles(x, *fespace);
783  }
784 
785  // Finite differences for computations of derivatives.
786  if (fdscheme)
787  {
788  MFEM_VERIFY(pa == false, "PA for finite differences is not implemented.");
789  tmop_integ->EnableFiniteDifferences(x);
790  }
791  tmop_integ->SetExactActionFlag(exactaction);
792 
793  // Setup the quadrature rules for the TMOP integrator.
794  IntegrationRules *irules = NULL;
795  switch (quad_type)
796  {
797  case 1: irules = &IntRulesLo; break;
798  case 2: irules = &IntRules; break;
799  case 3: irules = &IntRulesCU; break;
800  default: cout << "Unknown quad_type: " << quad_type << endl; return 3;
801  }
802  tmop_integ->SetIntegrationRules(*irules, quad_order);
803  if (dim == 2)
804  {
805  cout << "Triangle quadrature points: "
806  << irules->Get(Geometry::TRIANGLE, quad_order).GetNPoints()
807  << "\nQuadrilateral quadrature points: "
808  << irules->Get(Geometry::SQUARE, quad_order).GetNPoints() << endl;
809  }
810  if (dim == 3)
811  {
812  cout << "Tetrahedron quadrature points: "
813  << irules->Get(Geometry::TETRAHEDRON, quad_order).GetNPoints()
814  << "\nHexahedron quadrature points: "
815  << irules->Get(Geometry::CUBE, quad_order).GetNPoints()
816  << "\nPrism quadrature points: "
817  << irules->Get(Geometry::PRISM, quad_order).GetNPoints() << endl;
818  }
819 
820  // Limit the node movement.
821  // The limiting distances can be given by a general function of space.
822  FiniteElementSpace dist_fespace(mesh, fec); // scalar space
823  GridFunction dist(&dist_fespace);
824  dist = 1.0;
825  // The small_phys_size is relevant only with proper normalization.
826  if (normalization) { dist = small_phys_size; }
827  ConstantCoefficient lim_coeff(lim_const);
828  if (lim_const != 0.0) { tmop_integ->EnableLimiting(x0, dist, lim_coeff); }
829 
830  // Adaptive limiting.
831  GridFunction adapt_lim_gf0(&ind_fes);
832  ConstantCoefficient adapt_lim_coeff(adapt_lim_const);
833  AdaptivityEvaluator *adapt_lim_eval = NULL;
834  if (adapt_lim_const > 0.0)
835  {
836  MFEM_VERIFY(pa == false, "PA is not implemented for adaptive limiting");
837 
838  FunctionCoefficient adapt_lim_gf0_coeff(adapt_lim_fun);
839  adapt_lim_gf0.ProjectCoefficient(adapt_lim_gf0_coeff);
840 
841  if (adapt_eval == 0) { adapt_lim_eval = new AdvectorCG(al); }
842  else if (adapt_eval == 1)
843  {
844 #ifdef MFEM_USE_GSLIB
845  adapt_lim_eval = new InterpolatorFP;
846 #else
847  MFEM_ABORT("MFEM is not built with GSLIB support!");
848 #endif
849  }
850  else { MFEM_ABORT("Bad interpolation option."); }
851 
852  tmop_integ->EnableAdaptiveLimiting(adapt_lim_gf0, adapt_lim_coeff,
853  *adapt_lim_eval);
854  if (visualization)
855  {
856  socketstream vis1;
857  common::VisualizeField(vis1, "localhost", 19916, adapt_lim_gf0, "Zeta 0",
858  300, 600, 300, 300);
859  }
860  }
861 
862  // Has to be after the enabling of the limiting / alignment, as it computes
863  // normalization factors for these terms as well.
864  if (normalization) { tmop_integ->EnableNormalization(x0); }
865 
866  // 12. Setup the final NonlinearForm (which defines the integral of interest,
867  // its first and second derivatives). Here we can use a combination of
868  // metrics, i.e., optimize the sum of two integrals, where both are
869  // scaled by used-defined space-dependent weights. Note that there are no
870  // command-line options for the weights and the type of the second
871  // metric; one should update those in the code.
872  NonlinearForm a(fespace);
873  if (pa) { a.SetAssemblyLevel(AssemblyLevel::PARTIAL); }
874  ConstantCoefficient *metric_coeff1 = NULL;
875  TMOP_QualityMetric *metric2 = NULL;
876  TargetConstructor *target_c2 = NULL;
877  FunctionCoefficient metric_coeff2(weight_fun);
878 
879  // Explicit combination of metrics.
880  if (combomet > 0)
881  {
882  // First metric.
883  metric_coeff1 = new ConstantCoefficient(1.0);
884  tmop_integ->SetCoefficient(*metric_coeff1);
885 
886  // Second metric.
887  if (dim == 2) { metric2 = new TMOP_Metric_077; }
888  else { metric2 = new TMOP_Metric_315; }
889  TMOP_Integrator *tmop_integ2 = NULL;
890  if (combomet == 1)
891  {
892  target_c2 = new TargetConstructor(
894  target_c2->SetVolumeScale(0.01);
895  target_c2->SetNodes(x0);
896  tmop_integ2 = new TMOP_Integrator(metric2, target_c2, h_metric);
897  tmop_integ2->SetCoefficient(metric_coeff2);
898  }
899  else { tmop_integ2 = new TMOP_Integrator(metric2, target_c, h_metric); }
900  tmop_integ2->IntegrateOverTarget(integ_over_targ);
901  tmop_integ2->SetIntegrationRules(*irules, quad_order);
902  if (fdscheme) { tmop_integ2->EnableFiniteDifferences(x); }
903  tmop_integ2->SetExactActionFlag(exactaction);
904 
906  combo->AddTMOPIntegrator(tmop_integ);
907  combo->AddTMOPIntegrator(tmop_integ2);
908  if (normalization) { combo->EnableNormalization(x0); }
909  if (lim_const != 0.0) { combo->EnableLimiting(x0, dist, lim_coeff); }
910 
911  a.AddDomainIntegrator(combo);
912  }
913  else
914  {
915  a.AddDomainIntegrator(tmop_integ);
916  }
917 
918  if (pa) { a.Setup(); }
919 
920  // Compute the minimum det(J) of the starting mesh.
921  min_detJ = infinity();
922  const int NE = mesh->GetNE();
923  for (int i = 0; i < NE; i++)
924  {
925  const IntegrationRule &ir =
926  irules->Get(fespace->GetFE(i)->GetGeomType(), quad_order);
928  for (int j = 0; j < ir.GetNPoints(); j++)
929  {
930  transf->SetIntPoint(&ir.IntPoint(j));
931  min_detJ = min(min_detJ, transf->Jacobian().Det());
932  }
933  }
934  cout << "Minimum det(J) of the original mesh is " << min_detJ << endl;
935 
936  if (min_detJ < 0.0 && barrier_type == 0
937  && metric_id != 22 && metric_id != 211 && metric_id != 252
938  && metric_id != 311 && metric_id != 313 && metric_id != 352)
939  {
940  MFEM_ABORT("The input mesh is inverted! Try an untangling metric.");
941  }
942  if (min_detJ < 0.0)
943  {
944  MFEM_VERIFY(target_t == TargetConstructor::IDEAL_SHAPE_UNIT_SIZE,
945  "Untangling is supported only for ideal targets.");
946 
947  const DenseMatrix &Wideal =
949  min_detJ /= Wideal.Det();
950 
951  // Slightly below minJ0 to avoid div by 0.
952  min_detJ -= 0.01 * h0.Min();
953  }
954 
955  // For HR tests, the energy is normalized by the number of elements.
956  const double init_energy = a.GetGridFunctionEnergy(x) /
957  (hradaptivity ? mesh->GetNE() : 1);
958  double init_metric_energy = init_energy;
959  if (lim_const > 0.0 || adapt_lim_const > 0.0)
960  {
961  lim_coeff.constant = 0.0;
962  adapt_lim_coeff.constant = 0.0;
963  init_metric_energy = a.GetGridFunctionEnergy(x) /
964  (hradaptivity ? mesh->GetNE() : 1);
965  lim_coeff.constant = lim_const;
966  adapt_lim_coeff.constant = adapt_lim_const;
967  }
968 
969  // Visualize the starting mesh and metric values.
970  // Note that for combinations of metrics, this only shows the first metric.
971  if (visualization)
972  {
973  char title[] = "Initial metric values";
974  vis_tmop_metric_s(mesh_poly_deg, *metric, *target_c, *mesh, title, 0);
975  }
976 
977  // 13. Fix all boundary nodes, or fix only a given component depending on the
978  // boundary attributes of the given mesh. Attributes 1/2/3 correspond to
979  // fixed x/y/z components of the node. Attribute 4 corresponds to an
980  // entirely fixed node. Other boundary attributes do not affect the node
981  // movement boundary conditions.
982  if (move_bnd == false)
983  {
984  Array<int> ess_bdr(mesh->bdr_attributes.Max());
985  ess_bdr = 1;
986  a.SetEssentialBC(ess_bdr);
987  }
988  else
989  {
990  int n = 0;
991  for (int i = 0; i < mesh->GetNBE(); i++)
992  {
993  const int nd = fespace->GetBE(i)->GetDof();
994  const int attr = mesh->GetBdrElement(i)->GetAttribute();
995  MFEM_VERIFY(!(dim == 2 && attr == 3),
996  "Boundary attribute 3 must be used only for 3D meshes. "
997  "Adjust the attributes (1/2/3/4 for fixed x/y/z/all "
998  "components, rest for free nodes), or use -fix-bnd.");
999  if (attr == 1 || attr == 2 || attr == 3) { n += nd; }
1000  if (attr == 4) { n += nd * dim; }
1001  }
1002  Array<int> ess_vdofs(n);
1003  n = 0;
1004  for (int i = 0; i < mesh->GetNBE(); i++)
1005  {
1006  const int nd = fespace->GetBE(i)->GetDof();
1007  const int attr = mesh->GetBdrElement(i)->GetAttribute();
1008  fespace->GetBdrElementVDofs(i, vdofs);
1009  if (attr == 1) // Fix x components.
1010  {
1011  for (int j = 0; j < nd; j++)
1012  { ess_vdofs[n++] = vdofs[j]; }
1013  }
1014  else if (attr == 2) // Fix y components.
1015  {
1016  for (int j = 0; j < nd; j++)
1017  { ess_vdofs[n++] = vdofs[j+nd]; }
1018  }
1019  else if (attr == 3) // Fix z components.
1020  {
1021  for (int j = 0; j < nd; j++)
1022  { ess_vdofs[n++] = vdofs[j+2*nd]; }
1023  }
1024  else if (attr == 4) // Fix all components.
1025  {
1026  for (int j = 0; j < vdofs.Size(); j++)
1027  { ess_vdofs[n++] = vdofs[j]; }
1028  }
1029  }
1030  a.SetEssentialVDofs(ess_vdofs);
1031  }
1032 
1033  // As we use the inexact Newton method to solve the resulting nonlinear
1034  // system, here we setup the linear solver for the system's Jacobian.
1035  Solver *S = NULL, *S_prec = NULL;
1036  const double linsol_rtol = 1e-12;
1037  // Level of output.
1038  IterativeSolver::PrintLevel linsolver_print;
1039  if (verbosity_level == 2)
1040  { linsolver_print.Errors().Warnings().FirstAndLast(); }
1041  if (verbosity_level > 2)
1042  { linsolver_print.Errors().Warnings().Iterations(); }
1043  if (lin_solver == 0)
1044  {
1045  S = new DSmoother(1, 1.0, max_lin_iter);
1046  }
1047  else if (lin_solver == 1)
1048  {
1049  CGSolver *cg = new CGSolver;
1050  cg->SetMaxIter(max_lin_iter);
1051  cg->SetRelTol(linsol_rtol);
1052  cg->SetAbsTol(0.0);
1053  cg->SetPrintLevel(linsolver_print);
1054  S = cg;
1055  }
1056  else
1057  {
1058  MINRESSolver *minres = new MINRESSolver;
1059  minres->SetMaxIter(max_lin_iter);
1060  minres->SetRelTol(linsol_rtol);
1061  minres->SetAbsTol(0.0);
1062  minres->SetPrintLevel(linsolver_print);
1063  if (lin_solver == 3 || lin_solver == 4)
1064  {
1065  if (pa)
1066  {
1067  MFEM_VERIFY(lin_solver != 4, "PA l1-Jacobi is not implemented");
1068  auto js = new OperatorJacobiSmoother;
1069  js->SetPositiveDiagonal(true);
1070  S_prec = js;
1071  }
1072  else
1073  {
1074  auto ds = new DSmoother((lin_solver == 3) ? 0 : 1, 1.0, 1);
1075  ds->SetPositiveDiagonal(true);
1076  S_prec = ds;
1077  }
1078  minres->SetPreconditioner(*S_prec);
1079  }
1080  S = minres;
1081  }
1082 
1083  //
1084  // Perform the nonlinear optimization.
1085  //
1086  const IntegrationRule &ir =
1087  irules->Get(fespace->GetFE(0)->GetGeomType(), quad_order);
1088  TMOPNewtonSolver solver(ir, solver_type);
1089  // Provide all integration rules in case of a mixed mesh.
1090  solver.SetIntegrationRules(*irules, quad_order);
1091  // Specify linear solver when we use a Newton-based solver.
1092  if (solver_type == 0) { solver.SetPreconditioner(*S); }
1093  // For untangling, the solver will update the min det(T) values.
1094  solver.SetMinDetPtr(&min_detJ);
1095  solver.SetMaxIter(solver_iter);
1096  solver.SetRelTol(solver_rtol);
1097  solver.SetAbsTol(0.0);
1098  if (solver_art_type > 0)
1099  {
1100  solver.SetAdaptiveLinRtol(solver_art_type, 0.5, 0.9);
1101  }
1102  // Level of output.
1103  IterativeSolver::PrintLevel newton_print;
1104  if (verbosity_level > 0)
1105  { newton_print.Errors().Warnings().Iterations(); }
1106  solver.SetPrintLevel(newton_print);
1107  // hr-adaptivity solver.
1108  // If hr-adaptivity is disabled, r-adaptivity is done once using the
1109  // TMOPNewtonSolver.
1110  // Otherwise, "hr_iter" iterations of r-adaptivity are done followed by
1111  // "h_per_r_iter" iterations of h-adaptivity after each r-adaptivity.
1112  // The solver terminates if an h-adaptivity iteration does not modify
1113  // any element in the mesh.
1114  TMOPHRSolver hr_solver(*mesh, a, solver,
1115  x, move_bnd, hradaptivity,
1116  mesh_poly_deg, h_metric_id,
1117  n_hr_iter, n_h_iter);
1118  hr_solver.AddGridFunctionForUpdate(&x0);
1119  if (adapt_lim_const > 0.)
1120  {
1121  hr_solver.AddGridFunctionForUpdate(&adapt_lim_gf0);
1122  hr_solver.AddFESpaceForUpdate(&ind_fes);
1123  }
1124  hr_solver.Mult();
1125 
1126  // 15. Save the optimized mesh to a file. This output can be viewed later
1127  // using GLVis: "glvis -m optimized.mesh".
1128  {
1129  ofstream mesh_ofs("optimized.mesh");
1130  mesh_ofs.precision(14);
1131  mesh->Print(mesh_ofs);
1132  }
1133 
1134  // Report the final energy of the functional.
1135  const double fin_energy = a.GetGridFunctionEnergy(x) /
1136  (hradaptivity ? mesh->GetNE() : 1);
1137  double fin_metric_energy = fin_energy;
1138  if (lim_const > 0.0 || adapt_lim_const > 0.0)
1139  {
1140  lim_coeff.constant = 0.0;
1141  adapt_lim_coeff.constant = 0.0;
1142  fin_metric_energy = a.GetGridFunctionEnergy(x) /
1143  (hradaptivity ? mesh->GetNE() : 1);
1144  lim_coeff.constant = lim_const;
1145  adapt_lim_coeff.constant = adapt_lim_const;
1146  }
1147  std::cout << std::scientific << std::setprecision(4);
1148  cout << "Initial strain energy: " << init_energy
1149  << " = metrics: " << init_metric_energy
1150  << " + extra terms: " << init_energy - init_metric_energy << endl;
1151  cout << " Final strain energy: " << fin_energy
1152  << " = metrics: " << fin_metric_energy
1153  << " + extra terms: " << fin_energy - fin_metric_energy << endl;
1154  cout << "The strain energy decreased by: "
1155  << (init_energy - fin_energy) * 100.0 / init_energy << " %." << endl;
1156 
1157  // Visualize the final mesh and metric values.
1158  if (visualization)
1159  {
1160  char title[] = "Final metric values";
1161  vis_tmop_metric_s(mesh_poly_deg, *metric, *target_c, *mesh, title, 600);
1162  }
1163 
1164  if (adapt_lim_const > 0.0 && visualization)
1165  {
1166  socketstream vis0;
1167  common::VisualizeField(vis0, "localhost", 19916, adapt_lim_gf0, "Xi 0",
1168  600, 600, 300, 300);
1169  }
1170 
1171  // Visualize the mesh displacement.
1172  if (visualization)
1173  {
1174  osockstream sock(19916, "localhost");
1175  sock << "solution\n";
1176  mesh->Print(sock);
1177  x0 -= x;
1178  x0.Save(sock);
1179  sock.send();
1180  sock << "window_title 'Displacements'\n"
1181  << "window_geometry "
1182  << 1200 << " " << 0 << " " << 600 << " " << 600 << "\n"
1183  << "keys jRmclA" << endl;
1184  }
1185 
1186  delete S;
1187  delete S_prec;
1188  delete target_c2;
1189  delete metric2;
1190  delete metric_coeff1;
1191  delete adapt_lim_eval;
1192  delete target_c;
1193  delete hr_adapt_coeff;
1194  delete adapt_coeff;
1195  delete h_metric;
1196  delete metric;
1197  delete untangler_metric;
1198  delete fespace;
1199  delete fec;
1200  delete mesh;
1201 
1202  return 0;
1203 }
void discrete_aspr_3d(const Vector &x, Vector &v)
Conjugate gradient method.
Definition: solvers.hpp:493
int GetNPoints() const
Returns the number of the points in the integration rule.
Definition: intrules.hpp:253
Definition: ex25.cpp:149
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:96
3D non-barrier Shape (S) metric.
Definition: tmop.hpp:1051
void EnableFiniteDifferences(const GridFunction &x)
Enables FD-based approximation and computes dx.
Definition: tmop.cpp:4497
Class for grid function - Vector with associated FE space.
Definition: gridfunc.hpp:30
3D barrier Shape+Size (VS) metric, well-posed (invex).
Definition: tmop.hpp:843
Data type for scaled Jacobi-type smoother of sparse matrix.
const IntegrationRule & Get(int GeomType, int Order)
Returns an integration rule for given GeomType and Order.
Definition: intrules.cpp:980
IntegrationRules IntRulesLo(0, Quadrature1D::GaussLobatto)
void SetPositiveDiagonal(bool pos_diag=true)
Replace diagonal entries with their absolute values.
Definition: solvers.hpp:339
void SetIntegrationRules(IntegrationRules &irules, int order)
Prescribe a set of integration rules; relevant for mixed meshes.
Definition: tmop_tools.hpp:207
void SetFromTrueVector()
Shortcut for calling SetFromTrueDofs() with GetTrueVector() as argument.
Definition: gridfunc.hpp:150
3D barrier Shape+Size (VS) metric, well-posed (polyconvex).
Definition: tmop.hpp:967
A coefficient that is constant across space and time.
Definition: coefficient.hpp:84
void PrintOptions(std::ostream &out) const
Print the options.
Definition: optparser.cpp:331
Geometry::Type GetElementBaseGeometry(int i) const
Definition: mesh.hpp:1242
virtual void SetSerialDiscreteTargetOrientation(const GridFunction &tspec_)
Definition: tmop.cpp:1960
int main(int argc, char *argv[])
int Dimension() const
Dimension of the reference space used within the elements.
Definition: mesh.hpp:1020
3D barrier Shape (S) metric, well-posed (polyconvex & invex).
Definition: tmop.hpp:703
void SetAdaptivityEvaluator(AdaptivityEvaluator *ae)
Definition: tmop.hpp:1658
void PrintUsage(std::ostream &out) const
Print the usage message.
Definition: optparser.cpp:462
void SetIntPoint(const IntegrationPoint *ip)
Set the integration point ip that weights and Jacobians will be evaluated at.
Definition: eltrans.hpp:93
int GetAttribute() const
Return element&#39;s attribute.
Definition: element.hpp:55
void ConstructSizeGF(GridFunction &size)
double Det() const
Definition: densemat.cpp:487
T Max() const
Find the maximal element in the array, using the comparison operator < for class T.
Definition: array.cpp:68
Container class for integration rules.
Definition: intrules.hpp:412
Data type dense matrix using column-major storage.
Definition: densemat.hpp:23
3D barrier Shape (S) metric, well-posed (polyconvex & invex).
Definition: tmop.hpp:682
void SetVolumeScale(double vol_scale)
Used by target type IDEAL_SHAPE_EQUAL_SIZE. The default volume scale is 1.
Definition: tmop.hpp:1415
int GetNDofs() const
Returns number of degrees of freedom. This is the number of Local Degrees of Freedom.
Definition: fespace.hpp:706
3D Size (V) metric.
Definition: tmop.hpp:785
void SetIntegrationRules(IntegrationRules &irules, int order)
Prescribe a set of integration rules; relevant for mixed meshes.
Definition: tmop.hpp:2004
void EnableNormalization(const GridFunction &x)
Normalization factor that considers all integrators in the combination.
Definition: tmop.cpp:4784
3D compound barrier Shape+Size (VS) metric (polyconvex, balanced).
Definition: tmop.hpp:989
2D barrier Shape+Size (VS) metric (polyconvex).
Definition: tmop.hpp:1145
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:885
bool Good() const
Return true if the command line options were parsed successfully.
Definition: optparser.hpp:159
void vis_tmop_metric_s(int order, TMOP_QualityMetric &qm, const TargetConstructor &tc, Mesh &mesh, char *title, int position)
Definition: tmop_tools.cpp:907
void SetExactActionFlag(bool flag_)
Flag to control if exact action of Integration is effected.
Definition: tmop.hpp:2224
STL namespace.
MINRES method.
Definition: solvers.hpp:603
void Randomize(int seed=0)
Set random values in the vector.
Definition: vector.cpp:817
void IntegrateOverTarget(bool integ_over_target_)
Definition: tmop.hpp:2013
virtual int OrderJ() const =0
Return the order of the elements of the Jacobian of the transformation.
void SetMinDetPtr(double *md_ptr)
Definition: tmop_tools.hpp:213
2D barrier Shape+Size (VS) metric (not polyconvex).
Definition: tmop.hpp:340
3D barrier Shape+Size (VS) metric, well-posed (invex).
Definition: tmop.hpp:864
int GetNBE() const
Returns number of boundary elements.
Definition: mesh.hpp:1089
Geometry Geometries
Definition: fe.cpp:49
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:4685
2D Shifted barrier form of shape metric (mu_2).
Definition: tmop.hpp:389
virtual const FiniteElement * GetFE(int i) const
Returns pointer to the FiniteElement in the FiniteElementCollection associated with i&#39;th element in t...
Definition: fespace.cpp:2841
virtual void SetSerialDiscreteTargetAspectRatio(const GridFunction &tspec_)
Definition: tmop.cpp:1950
IntegrationRules IntRules(0, Quadrature1D::GaussLegendre)
A global object with all integration rules (defined in intrules.cpp)
Definition: intrules.hpp:480
2D barrier shape (S) metric (not polyconvex).
Definition: tmop.hpp:466
3D barrier Shape (S) metric, well-posed (polyconvex & invex).
Definition: tmop.hpp:724
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
IntegrationPoint & IntPoint(int i)
Returns a reference to the i-th integration point.
Definition: intrules.hpp:256
void SetTrueVector()
Shortcut for calling GetTrueDofs() with GetTrueVector() as argument.
Definition: gridfunc.hpp:144
2D compound barrier Shape+Size (VS) metric (balanced).
Definition: tmop.hpp:588
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
Definition: solvers.hpp:616
virtual void SetPrintLevel(int print_lvl)
Legacy method to set the level of verbosity of the solver output.
Definition: solvers.cpp:71
double GetElementSize(ElementTransformation *T, int type=0)
Definition: mesh.cpp:83
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:2028
void EnsureNCMesh(bool simplices_nonconforming=false)
Definition: mesh.cpp:9888
3D Size (V) metric.
Definition: tmop.hpp:803
Jacobi smoothing for a given bilinear form (no matrix necessary).
Definition: solvers.hpp:302
double b
Definition: lissajous.cpp:42
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:10232
const DenseMatrix & Jacobian()
Return the Jacobian matrix of the transformation at the currently set IntegrationPoint, using the method SetIntPoint().
Definition: eltrans.hpp:119
void SetNodes(const GridFunction &n)
Set the nodes to be used in the target-matrix construction.
Definition: tmop.hpp:1409
void SetMaxIter(int max_it)
Definition: solvers.hpp:201
3D barrier Shape+Size (VS) metric, well-posed (polyconvex).
Definition: tmop.hpp:1010
Version of QuadraticFECollection with positive basis functions.
Definition: fe_coll.hpp:764
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
2D barrier Shape+Size+Orientation (VOS) metric (polyconvex).
Definition: tmop.hpp:609
virtual void SetNodalFESpace(FiniteElementSpace *nfes)
Definition: mesh.cpp:5577
void AddFESpaceForUpdate(FiniteElementSpace *fes)
Definition: tmop_amr.hpp:253
2D barrier Shape+Orientation (OS) metric (polyconvex).
Definition: tmop.hpp:1127
A general vector function coefficient.
3D compound barrier Shape+Size (VS) metric (polyconvex, balanced).
Definition: tmop.hpp:906
void SetAbsTol(double atol)
Definition: solvers.hpp:200
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
Definition: mesh.hpp:275
void AddGridFunctionForUpdate(GridFunction *gf)
Definition: tmop_amr.hpp:252
void SetRelTol(double rtol)
Definition: solvers.hpp:199
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
Definition: tmop_tools.hpp:286
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition: fespace.hpp:219
double weight_fun(const Vector &x)
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
void SetAdaptiveLinRtol(const int type=2, const double rtol0=0.5, const double rtol_max=0.9, const double alpha=0.5 *(1.0+sqrt(5.0)), const double gamma=1.0)
Enable adaptive linear solver relative tolerance algorithm.
Definition: solvers.cpp:1927
double adapt_lim_fun(const Vector &x)
int GetDof() const
Returns the number of degrees of freedom in the finite element.
Definition: fe_base.hpp:323
virtual DofTransformation * GetElementDofs(int elem, Array< int > &dofs) const
Returns indices of degrees of freedom of element &#39;elem&#39;. The returned indices are offsets into an ldo...
Definition: fespace.cpp:2731
void AddTMOPIntegrator(TMOP_Integrator *ti)
Adds a new TMOP_Integrator to the combination.
Definition: tmop.hpp:2254
double material_indicator_2d(const Vector &x)
2D barrier Shape+Size (VS) metric (not polyconvex).
Definition: tmop.hpp:358
int GetNE() const
Returns number of elements.
Definition: mesh.hpp:1086
double a
Definition: lissajous.cpp:41
Class for integration point with weight.
Definition: intrules.hpp:31
void GetElementTransformation(int i, IsoparametricTransformation *ElTr)
Definition: mesh.cpp:355
void DiffuseField(ParGridFunction &field, int smooth_steps)
Definition: dist_solver.cpp:22
DofTransformation * GetBdrElementVDofs(int i, Array< int > &vdofs) const
Returns indices of degrees of freedom for i&#39;th boundary element. The returned indices are offsets int...
Definition: fespace.cpp:297
2D barrier Shape+Size+Orientation (VOS) metric (polyconvex).
Definition: tmop.hpp:1109
3D Shape (S) metric, untangling version of 303.
Definition: tmop.hpp:764
int dim
Definition: ex24.cpp:53
IntegrationRules IntRulesCU(0, Quadrature1D::ClosedUniform)
2D compound barrier Shape+Size (VS) metric (balanced).
Definition: tmop.hpp:567
int GetNBE() const
Returns number of boundary elements in the mesh.
Definition: fespace.hpp:745
double infinity()
Define a shortcut for std::numeric_limits<double>::infinity()
Definition: vector.hpp:44
const DenseMatrix & GetGeomToPerfGeomJac(int GeomType) const
Definition: geom.hpp:98
virtual void ProjectCoefficient(Coefficient &coeff)
Project coeff Coefficient to this GridFunction. The projection computation depends on the choice of t...
Definition: gridfunc.cpp:2415
int Size() const
Return the logical size of the array.
Definition: array.hpp:141
A general function coefficient.
int DofToVDof(int dof, int vd, int ndofs=-1) const
Compute a single vdof corresponding to the index dof and the vector index vd.
Definition: fespace.cpp:251
Vector data type.
Definition: vector.hpp:58
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
virtual void Print(std::ostream &os=mfem::out) const
Definition: mesh.hpp:2011
3D barrier Shape+Size (VS) metric, well-posed (polyconvex).
Definition: tmop.hpp:948
Arbitrary order H1-conforming (continuous) finite elements.
Definition: fe_coll.hpp:259
TargetType
Target-matrix construction algorithms supported by this class.
Definition: tmop.hpp:1333
const Element * GetBdrElement(int i) const
Return pointer to the i&#39;th boundary element object.
Definition: mesh.hpp:1158
Base class for solvers.
Definition: operator.hpp:682
void SetNodalGridFunction(GridFunction *nodes, bool make_owner=false)
Definition: mesh.cpp:5624
Base class representing target-matrix construction algorithms for mesh optimization via the target-ma...
Definition: tmop.hpp:1329
virtual void Save(std::ostream &out) const
Save the GridFunction to an output stream.
Definition: gridfunc.cpp:3696
The MFEM Device class abstracts hardware devices such as GPUs, as well as programming models such as ...
Definition: device.hpp:121
virtual void SetAnalyticTargetSpec(Coefficient *sspec, VectorCoefficient *vspec, TMOPMatrixCoefficient *mspec)
Definition: tmop.cpp:1701
3D compound barrier Shape+Size (VS) metric (polyconvex).
Definition: tmop.hpp:927
virtual void SetSerialDiscreteTargetSize(const GridFunction &tspec_)
Definition: tmop.cpp:1930
const FiniteElement * GetBE(int i) const
Returns pointer to the FiniteElement in the FiniteElementCollection associated with i&#39;th boundary fac...
Definition: fespace.cpp:3166
2D barrier Shape+Orientation (OS) metric (polyconvex).
Definition: tmop.hpp:552
virtual double * HostReadWrite()
Shortcut for mfem::ReadWrite(vec.GetMemory(), vec.Size(), false).
Definition: vector.hpp:473
double GetElementVolume(int i)
Definition: mesh.cpp:119
Settings for the output behavior of the IterativeSolver.
Definition: solvers.hpp:78
3D barrier Shape (S) metric, well-posed (polyconvex & invex).
Definition: tmop.hpp:663
2D non-barrier metric without a type.
Definition: tmop.hpp:221
A TMOP integrator class based on any given TMOP_QualityMetric and TargetConstructor.
Definition: tmop.hpp:1733
2D non-barrier Shape+Size+Orientation (VOS) metric (polyconvex).
Definition: tmop.hpp:374