MFEM  v4.3.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
pmesh-optimizer.cpp
Go to the documentation of this file.
1 // Copyright (c) 2010-2021, 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 discrete size:
43 // mpirun -np 4 pmesh-optimizer -m square01.mesh -o 2 -rs 2 -mid 80 -tid 5 -ni 50 -qo 4 -nor
44 // Adapted discrete size 3D with PA:
45 // mpirun -np 4 pmesh-optimizer -m cube.mesh -o 2 -rs 2 -mid 321 -tid 5 -ls 3 -nor -pa
46 // Adapted discrete size 3D with PA on device (requires CUDA).
47 // * 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
48 // Adapted discrete size; explicit combo of metrics; mixed tri/quad mesh:
49 // 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
50 // Adapted discrete size+aspect_ratio:
51 // mpirun -np 4 pmesh-optimizer -m square01.mesh -o 2 -rs 2 -mid 7 -tid 6 -ni 100
52 // mpirun -np 4 pmesh-optimizer -m square01.mesh -o 2 -rs 2 -mid 7 -tid 6 -ni 100 -qo 6 -ex -st 1 -nor
53 // Adapted discrete size+orientation (requires GSLIB):
54 // * mpirun -np 4 pmesh-optimizer -m square01.mesh -o 2 -rs 2 -mid 36 -tid 8 -qo 4 -fd -ae 1 -nor
55 // Adapted discrete aspect_ratio+orientation (requires GSLIB):
56 // * 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
57 // Adapted discrete aspect ratio (3D):
58 // mpirun -np 4 pmesh-optimizer -m cube.mesh -o 2 -rs 2 -mid 302 -tid 7 -ni 20 -bnd -qt 1 -qo 8
59 //
60 // Adaptive limiting:
61 // mpirun -np 4 pmesh-optimizer -m stretched2D.mesh -o 2 -mid 2 -tid 1 -ni 50 -qo 5 -nor -vl 1 -alc 0.5
62 // Adaptive limiting through the L-BFGS solver:
63 // 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
64 // Adaptive limiting through FD (requires GSLIB):
65 // * 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
66 //
67 // Blade shape:
68 // 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
69 // Blade shape with FD-based solver:
70 // mpirun -np 4 pmesh-optimizer -m blade.mesh -o 4 -mid 2 -tid 1 -ni 30 -ls 4 -bnd -qt 1 -qo 8 -fd
71 // Blade limited shape:
72 // mpirun -np 4 pmesh-optimizer -m blade.mesh -o 4 -mid 2 -tid 1 -bnd -qt 1 -qo 8 -lc 5000
73 // ICF shape and equal size:
74 // mpirun -np 4 pmesh-optimizer -o 3 -mid 9 -tid 2 -ni 25 -ls 3 -art 2 -qo 5
75 // ICF shape and initial size:
76 // mpirun -np 4 pmesh-optimizer -o 3 -mid 9 -tid 3 -ni 30 -ls 3 -bnd -qt 1 -qo 8
77 // ICF shape:
78 // mpirun -np 4 pmesh-optimizer -o 3 -mid 1 -tid 1 -ni 100 -bnd -qt 1 -qo 8
79 // ICF limited shape:
80 // mpirun -np 4 pmesh-optimizer -o 3 -mid 1 -tid 1 -ni 100 -bnd -qt 1 -qo 8 -lc 10
81 // ICF combo shape + size (rings, slow convergence):
82 // mpirun -np 4 pmesh-optimizer -o 3 -mid 1 -tid 1 -ni 1000 -bnd -qt 1 -qo 8 -cmb 1
83 // Mixed tet / cube / hex mesh with limiting:
84 // 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
85 // 3D pinched sphere shape (the mesh is in the mfem/data GitHub repository):
86 // * mpirun -np 4 pmesh-optimizer -m ../../../mfem_data/ball-pert.mesh -o 4 -mid 303 -tid 1 -ni 20 -li 500 -fix-bnd
87 // 2D non-conforming shape and equal size:
88 // 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
89 //
90 // 2D untangling:
91 // mpirun -np 4 pmesh-optimizer -m jagged.mesh -o 2 -mid 22 -tid 1 -ni 50 -li 50 -qo 4 -fd -vl 1
92 // 3D untangling (the mesh is in the mfem/data GitHub repository):
93 // * 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
94 //
95 
96 #include "mfem.hpp"
97 #include "../common/mfem-common.hpp"
98 #include <iostream>
99 #include <fstream>
100 #include "mesh-optimizer.hpp"
101 
102 using namespace mfem;
103 using namespace std;
104 
105 int main (int argc, char *argv[])
106 {
107  // 0. Initialize MPI.
108  int num_procs, myid;
109  MPI_Init(&argc, &argv);
110  MPI_Comm_size(MPI_COMM_WORLD, &num_procs);
111  MPI_Comm_rank(MPI_COMM_WORLD, &myid);
112 
113  // 1. Set the method's default parameters.
114  const char *mesh_file = "icf.mesh";
115  int mesh_poly_deg = 1;
116  int rs_levels = 0;
117  int rp_levels = 0;
118  double jitter = 0.0;
119  int metric_id = 1;
120  int target_id = 1;
121  double lim_const = 0.0;
122  double adapt_lim_const = 0.0;
123  int quad_type = 1;
124  int quad_order = 8;
125  int solver_type = 0;
126  int solver_iter = 20;
127  double solver_rtol = 1e-10;
128  int solver_art_type = 0;
129  int lin_solver = 2;
130  int max_lin_iter = 100;
131  bool move_bnd = true;
132  int combomet = 0;
133  bool normalization = false;
134  bool visualization = true;
135  int verbosity_level = 0;
136  bool fdscheme = false;
137  int adapt_eval = 0;
138  bool exactaction = false;
139  const char *devopt = "cpu";
140  bool pa = false;
141 
142  // 2. Parse command-line options.
143  OptionsParser args(argc, argv);
144  args.AddOption(&mesh_file, "-m", "--mesh",
145  "Mesh file to use.");
146  args.AddOption(&mesh_poly_deg, "-o", "--order",
147  "Polynomial degree of mesh finite element space.");
148  args.AddOption(&rs_levels, "-rs", "--refine-serial",
149  "Number of times to refine the mesh uniformly in serial.");
150  args.AddOption(&rp_levels, "-rp", "--refine-parallel",
151  "Number of times to refine the mesh uniformly in parallel.");
152  args.AddOption(&jitter, "-ji", "--jitter",
153  "Random perturbation scaling factor.");
154  args.AddOption(&metric_id, "-mid", "--metric-id",
155  "Mesh optimization metric:\n\t"
156  "T-metrics\n\t"
157  "1 : |T|^2 -- 2D shape\n\t"
158  "2 : 0.5|T|^2/tau-1 -- 2D shape (condition number)\n\t"
159  "7 : |T-T^-t|^2 -- 2D shape+size\n\t"
160  "9 : tau*|T-T^-t|^2 -- 2D shape+size\n\t"
161  "14 : |T-I|^2 -- 2D shape+size+orientation\n\t"
162  "22 : 0.5(|T|^2-2*tau)/(tau-tau_0) -- 2D untangling\n\t"
163  "50 : 0.5|T^tT|^2/tau^2-1 -- 2D shape\n\t"
164  "55 : (tau-1)^2 -- 2D size\n\t"
165  "56 : 0.5(sqrt(tau)-1/sqrt(tau))^2 -- 2D size\n\t"
166  "58 : |T^tT|^2/(tau^2)-2*|T|^2/tau+2 -- 2D shape\n\t"
167  "77 : 0.5(tau-1/tau)^2 -- 2D size\n\t"
168  "80 : (1-gamma)mu_2 + gamma mu_77 -- 2D shape+size\n\t"
169  "85 : |T-|T|/sqrt(2)I|^2 -- 2D shape+orientation\n\t"
170  "98 : (1/tau)|T-I|^2 -- 2D shape+size+orientation\n\t"
171  // "211: (tau-1)^2-tau+sqrt(tau^2+eps) -- 2D untangling\n\t"
172  // "252: 0.5(tau-1)^2/(tau-tau_0) -- 2D untangling\n\t"
173  "301: (|T||T^-1|)/3-1 -- 3D shape\n\t"
174  "302: (|T|^2|T^-1|^2)/9-1 -- 3D shape\n\t"
175  "303: (|T|^2)/3*tau^(2/3)-1 -- 3D shape\n\t"
176  // "311: (tau-1)^2-tau+sqrt(tau^2+eps)-- 3D untangling\n\t"
177  "313: (|T|^2)(tau-tau0)^(-2/3)/3 -- 3D untangling\n\t"
178  "315: (tau-1)^2 -- 3D size\n\t"
179  "316: 0.5(sqrt(tau)-1/sqrt(tau))^2 -- 3D size\n\t"
180  "321: |T-T^-t|^2 -- 3D shape+size\n\t"
181  // "352: 0.5(tau-1)^2/(tau-tau_0) -- 3D untangling\n\t"
182  "A-metrics\n\t"
183  "11 : (1/4*alpha)|A-(adjA)^T(W^TW)/omega|^2 -- 2D shape\n\t"
184  "36 : (1/alpha)|A-W|^2 -- 2D shape+size+orientation\n\t"
185  "107: (1/2*alpha)|A-|A|/|W|W|^2 -- 2D shape+orientation\n\t"
186  "126: (1-gamma)nu_11 + gamma*nu_14a -- 2D shape+size\n\t"
187  );
188  args.AddOption(&target_id, "-tid", "--target-id",
189  "Target (ideal element) type:\n\t"
190  "1: Ideal shape, unit size\n\t"
191  "2: Ideal shape, equal size\n\t"
192  "3: Ideal shape, initial size\n\t"
193  "4: Given full analytic Jacobian (in physical space)\n\t"
194  "5: Ideal shape, given size (in physical space)");
195  args.AddOption(&lim_const, "-lc", "--limit-const", "Limiting constant.");
196  args.AddOption(&adapt_lim_const, "-alc", "--adapt-limit-const",
197  "Adaptive limiting coefficient constant.");
198  args.AddOption(&quad_type, "-qt", "--quad-type",
199  "Quadrature rule type:\n\t"
200  "1: Gauss-Lobatto\n\t"
201  "2: Gauss-Legendre\n\t"
202  "3: Closed uniform points");
203  args.AddOption(&quad_order, "-qo", "--quad_order",
204  "Order of the quadrature rule.");
205  args.AddOption(&solver_type, "-st", "--solver-type",
206  " Type of solver: (default) 0: Newton, 1: LBFGS");
207  args.AddOption(&solver_iter, "-ni", "--newton-iters",
208  "Maximum number of Newton iterations.");
209  args.AddOption(&solver_rtol, "-rtol", "--newton-rel-tolerance",
210  "Relative tolerance for the Newton solver.");
211  args.AddOption(&solver_art_type, "-art", "--adaptive-rel-tol",
212  "Type of adaptive relative linear solver tolerance:\n\t"
213  "0: None (default)\n\t"
214  "1: Eisenstat-Walker type 1\n\t"
215  "2: Eisenstat-Walker type 2");
216  args.AddOption(&lin_solver, "-ls", "--lin-solver",
217  "Linear solver:\n\t"
218  "0: l1-Jacobi\n\t"
219  "1: CG\n\t"
220  "2: MINRES\n\t"
221  "3: MINRES + Jacobi preconditioner\n\t"
222  "4: MINRES + l1-Jacobi preconditioner");
223  args.AddOption(&max_lin_iter, "-li", "--lin-iter",
224  "Maximum number of iterations in the linear solve.");
225  args.AddOption(&move_bnd, "-bnd", "--move-boundary", "-fix-bnd",
226  "--fix-boundary",
227  "Enable motion along horizontal and vertical boundaries.");
228  args.AddOption(&combomet, "-cmb", "--combo-type",
229  "Combination of metrics options:\n\t"
230  "0: Use single metric\n\t"
231  "1: Shape + space-dependent size given analytically\n\t"
232  "2: Shape + adapted size given discretely; shared target");
233  args.AddOption(&normalization, "-nor", "--normalization", "-no-nor",
234  "--no-normalization",
235  "Make all terms in the optimization functional unitless.");
236  args.AddOption(&fdscheme, "-fd", "--fd_approximation",
237  "-no-fd", "--no-fd-approx",
238  "Enable finite difference based derivative computations.");
239  args.AddOption(&exactaction, "-ex", "--exact_action",
240  "-no-ex", "--no-exact-action",
241  "Enable exact action of TMOP_Integrator.");
242  args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
243  "--no-visualization",
244  "Enable or disable GLVis visualization.");
245  args.AddOption(&verbosity_level, "-vl", "--verbosity-level",
246  "Set the verbosity level - 0, 1, or 2.");
247  args.AddOption(&adapt_eval, "-ae", "--adaptivity-evaluator",
248  "0 - Advection based (DEFAULT), 1 - GSLIB.");
249  args.AddOption(&devopt, "-d", "--device",
250  "Device configuration string, see Device::Configure().");
251  args.AddOption(&pa, "-pa", "--partial-assembly", "-no-pa",
252  "--no-partial-assembly", "Enable Partial Assembly.");
253  args.Parse();
254  if (!args.Good())
255  {
256  if (myid == 0) { args.PrintUsage(cout); }
257  return 1;
258  }
259  if (myid == 0) { args.PrintOptions(cout); }
260 
261  Device device(devopt);
262  if (myid == 0) { device.Print();}
263 
264  // 3. Initialize and refine the starting mesh.
265  Mesh *mesh = new Mesh(mesh_file, 1, 1, false);
266  for (int lev = 0; lev < rs_levels; lev++)
267  {
268  mesh->UniformRefinement();
269  }
270  const int dim = mesh->Dimension();
271  if (myid == 0)
272  {
273  cout << "Mesh curvature: ";
274  if (mesh->GetNodes()) { cout << mesh->GetNodes()->OwnFEC()->Name(); }
275  else { cout << "(NONE)"; }
276  cout << endl;
277  }
278 
279  ParMesh *pmesh = new ParMesh(MPI_COMM_WORLD, *mesh);
280 
281  delete mesh;
282  for (int lev = 0; lev < rp_levels; lev++)
283  {
284  pmesh->UniformRefinement();
285  }
286 
287  // 4. Define a finite element space on the mesh. Here we use vector finite
288  // elements which are tensor products of quadratic finite elements. The
289  // number of components in the vector finite element space is specified by
290  // the last parameter of the FiniteElementSpace constructor.
292  if (mesh_poly_deg <= 0)
293  {
294  fec = new QuadraticPosFECollection;
295  mesh_poly_deg = 2;
296  }
297  else { fec = new H1_FECollection(mesh_poly_deg, dim); }
298  ParFiniteElementSpace *pfespace = new ParFiniteElementSpace(pmesh, fec, dim);
299 
300  // 5. Make the mesh curved based on the above finite element space. This
301  // means that we define the mesh elements through a fespace-based
302  // transformation of the reference element.
303  pmesh->SetNodalFESpace(pfespace);
304 
305  // 6. Set up an empty right-hand side vector b, which is equivalent to b=0.
306  Vector b(0);
307 
308  // 7. Get the mesh nodes (vertices and other degrees of freedom in the finite
309  // element space) as a finite element grid function in fespace. Note that
310  // changing x automatically changes the shapes of the mesh elements.
311  ParGridFunction x(pfespace);
312  pmesh->SetNodalGridFunction(&x);
313 
314  // 8. Define a vector representing the minimal local mesh size in the mesh
315  // nodes. We index the nodes using the scalar version of the degrees of
316  // freedom in pfespace. Note: this is partition-dependent.
317  //
318  // In addition, compute average mesh size and total volume.
319  Vector h0(pfespace->GetNDofs());
320  h0 = infinity();
321  double vol_loc = 0.0;
322  Array<int> dofs;
323  for (int i = 0; i < pmesh->GetNE(); i++)
324  {
325  // Get the local scalar element degrees of freedom in dofs.
326  pfespace->GetElementDofs(i, dofs);
327  // Adjust the value of h0 in dofs based on the local mesh size.
328  const double hi = pmesh->GetElementSize(i);
329  for (int j = 0; j < dofs.Size(); j++)
330  {
331  h0(dofs[j]) = min(h0(dofs[j]), hi);
332  }
333  vol_loc += pmesh->GetElementVolume(i);
334  }
335  double volume;
336  MPI_Allreduce(&vol_loc, &volume, 1, MPI_DOUBLE, MPI_SUM, MPI_COMM_WORLD);
337  const double small_phys_size = pow(volume, 1.0 / dim) / 100.0;
338 
339  // 9. Add a random perturbation to the nodes in the interior of the domain.
340  // We define a random grid function of fespace and make sure that it is
341  // zero on the boundary and its values are locally of the order of h0.
342  // The latter is based on the DofToVDof() method which maps the scalar to
343  // the vector degrees of freedom in pfespace.
344  ParGridFunction rdm(pfespace);
345  rdm.Randomize();
346  rdm -= 0.25; // Shift to random values in [-0.5,0.5].
347  rdm *= jitter;
348  rdm.HostReadWrite();
349  // Scale the random values to be of order of the local mesh size.
350  for (int i = 0; i < pfespace->GetNDofs(); i++)
351  {
352  for (int d = 0; d < dim; d++)
353  {
354  rdm(pfespace->DofToVDof(i,d)) *= h0(i);
355  }
356  }
357  Array<int> vdofs;
358  for (int i = 0; i < pfespace->GetNBE(); i++)
359  {
360  // Get the vector degrees of freedom in the boundary element.
361  pfespace->GetBdrElementVDofs(i, vdofs);
362  // Set the boundary values to zero.
363  for (int j = 0; j < vdofs.Size(); j++) { rdm(vdofs[j]) = 0.0; }
364  }
365  x -= rdm;
366  // Set the perturbation of all nodes from the true nodes.
367  x.SetTrueVector();
368  x.SetFromTrueVector();
369 
370  // 10. Save the starting (prior to the optimization) mesh to a file. This
371  // output can be viewed later using GLVis: "glvis -m perturbed -np
372  // num_mpi_tasks".
373  {
374  ostringstream mesh_name;
375  mesh_name << "perturbed.mesh";
376  ofstream mesh_ofs(mesh_name.str().c_str());
377  mesh_ofs.precision(8);
378  pmesh->PrintAsOne(mesh_ofs);
379  }
380 
381  // 11. Store the starting (prior to the optimization) positions.
382  ParGridFunction x0(pfespace);
383  x0 = x;
384 
385  // 12. Form the integrator that uses the chosen metric and target.
386  double tauval = -0.1;
387  TMOP_QualityMetric *metric = NULL;
388  switch (metric_id)
389  {
390  // T-metrics
391  case 1: metric = new TMOP_Metric_001; break;
392  case 2: metric = new TMOP_Metric_002; break;
393  case 7: metric = new TMOP_Metric_007; break;
394  case 9: metric = new TMOP_Metric_009; break;
395  case 14: metric = new TMOP_Metric_014; break;
396  case 22: metric = new TMOP_Metric_022(tauval); break;
397  case 50: metric = new TMOP_Metric_050; break;
398  case 55: metric = new TMOP_Metric_055; break;
399  case 56: metric = new TMOP_Metric_056; break;
400  case 58: metric = new TMOP_Metric_058; break;
401  case 77: metric = new TMOP_Metric_077; break;
402  case 80: metric = new TMOP_Metric_080(0.5); break;
403  case 85: metric = new TMOP_Metric_085; break;
404  case 98: metric = new TMOP_Metric_098; break;
405  // case 211: metric = new TMOP_Metric_211; break;
406  // case 252: metric = new TMOP_Metric_252(tauval); break;
407  case 301: metric = new TMOP_Metric_301; break;
408  case 302: metric = new TMOP_Metric_302; break;
409  case 303: metric = new TMOP_Metric_303; break;
410  // case 311: metric = new TMOP_Metric_311; break;
411  case 313: metric = new TMOP_Metric_313(tauval); break;
412  case 315: metric = new TMOP_Metric_315; break;
413  case 316: metric = new TMOP_Metric_316; break;
414  case 321: metric = new TMOP_Metric_321; break;
415  // case 352: metric = new TMOP_Metric_352(tauval); break;
416  // A-metrics
417  case 11: metric = new TMOP_AMetric_011; break;
418  case 36: metric = new TMOP_AMetric_036; break;
419  case 107: metric = new TMOP_AMetric_107a; break;
420  case 126: metric = new TMOP_AMetric_126(0.9); break;
421  default:
422  if (myid == 0) { cout << "Unknown metric_id: " << metric_id << endl; }
423  return 3;
424  }
426  TargetConstructor *target_c = NULL;
427  HessianCoefficient *adapt_coeff = NULL;
428  H1_FECollection ind_fec(mesh_poly_deg, dim);
429  ParFiniteElementSpace ind_fes(pmesh, &ind_fec);
430  ParFiniteElementSpace ind_fesv(pmesh, &ind_fec, dim);
431  ParGridFunction size(&ind_fes), aspr(&ind_fes), disc(&ind_fes), ori(&ind_fes);
432  ParGridFunction aspr3d(&ind_fesv);
433 
434  const AssemblyLevel al =
436 
437  switch (target_id)
438  {
439  case 1: target_t = TargetConstructor::IDEAL_SHAPE_UNIT_SIZE; break;
440  case 2: target_t = TargetConstructor::IDEAL_SHAPE_EQUAL_SIZE; break;
441  case 3: target_t = TargetConstructor::IDEAL_SHAPE_GIVEN_SIZE; break;
442  case 4:
443  {
445  AnalyticAdaptTC *tc = new AnalyticAdaptTC(target_t);
446  adapt_coeff = new HessianCoefficient(dim, metric_id);
447  tc->SetAnalyticTargetSpec(NULL, NULL, adapt_coeff);
448  target_c = tc;
449  break;
450  }
451  case 5: // Discrete size 2D or 3D
452  {
454  DiscreteAdaptTC *tc = new DiscreteAdaptTC(target_t);
455  if (adapt_eval == 0)
456  {
457  tc->SetAdaptivityEvaluator(new AdvectorCG(al));
458  }
459  else
460  {
461 #ifdef MFEM_USE_GSLIB
463 #else
464  MFEM_ABORT("MFEM is not built with GSLIB.");
465 #endif
466  }
467  if (dim == 2)
468  {
470  size.ProjectCoefficient(ind_coeff);
471  }
472  else if (dim == 3)
473  {
475  size.ProjectCoefficient(ind_coeff);
476  }
477  tc->SetParDiscreteTargetSize(size);
478  target_c = tc;
479  break;
480  }
481  case 6: // material indicator 2D
482  {
483  ParGridFunction d_x(&ind_fes), d_y(&ind_fes);
484 
486  DiscreteAdaptTC *tc = new DiscreteAdaptTC(target_t);
488  disc.ProjectCoefficient(ind_coeff);
489  if (adapt_eval == 0)
490  {
491  tc->SetAdaptivityEvaluator(new AdvectorCG(al));
492  }
493  else
494  {
495 #ifdef MFEM_USE_GSLIB
497 #else
498  MFEM_ABORT("MFEM is not built with GSLIB.");
499 #endif
500  }
501  // Diffuse the interface
502  DiffuseField(disc,2);
503 
504  // Get partials with respect to x and y of the grid function
505  disc.GetDerivative(1,0,d_x);
506  disc.GetDerivative(1,1,d_y);
507 
508  // Compute the squared magnitude of the gradient
509  for (int i = 0; i < size.Size(); i++)
510  {
511  size(i) = std::pow(d_x(i),2)+std::pow(d_y(i),2);
512  }
513  const double max = size.Max();
514  double max_all;
515  MPI_Allreduce(&max, &max_all, 1, MPI_DOUBLE, MPI_MAX, MPI_COMM_WORLD);
516 
517  for (int i = 0; i < d_x.Size(); i++)
518  {
519  d_x(i) = std::abs(d_x(i));
520  d_y(i) = std::abs(d_y(i));
521  }
522  const double eps = 0.01;
523  const double aspr_ratio = 20.0;
524  const double size_ratio = 40.0;
525 
526  for (int i = 0; i < size.Size(); i++)
527  {
528  size(i) = (size(i)/max_all);
529  aspr(i) = (d_x(i)+eps)/(d_y(i)+eps);
530  aspr(i) = 0.1 + 0.9*(1-size(i))*(1-size(i));
531  if (aspr(i) > aspr_ratio) {aspr(i) = aspr_ratio;}
532  if (aspr(i) < 1.0/aspr_ratio) {aspr(i) = 1.0/aspr_ratio;}
533  }
534  Vector vals;
535  const int NE = pmesh->GetNE();
536  double volume = 0.0, volume_ind = 0.0;
537 
538  for (int i = 0; i < NE; i++)
539  {
541  const IntegrationRule &ir =
542  IntRules.Get(pmesh->GetElementBaseGeometry(i), Tr->OrderJ());
543  size.GetValues(i, ir, vals);
544  for (int j = 0; j < ir.GetNPoints(); j++)
545  {
546  const IntegrationPoint &ip = ir.IntPoint(j);
547  Tr->SetIntPoint(&ip);
548  volume += ip.weight * Tr->Weight();
549  volume_ind += vals(j) * ip.weight * Tr->Weight();
550  }
551  }
552  double volume_all, volume_ind_all;
553  MPI_Allreduce(&volume, &volume_all, 1, MPI_DOUBLE, MPI_SUM, MPI_COMM_WORLD);
554  MPI_Allreduce(&volume_ind, &volume_ind_all, 1, MPI_DOUBLE, MPI_SUM,
555  MPI_COMM_WORLD);
556  const int NE_ALL = pmesh->GetGlobalNE();
557 
558  const double avg_zone_size = volume_all / NE_ALL;
559 
560  const double small_avg_ratio =
561  (volume_ind_all + (volume_all - volume_ind_all) / size_ratio)
562  / volume_all;
563 
564  const double small_zone_size = small_avg_ratio * avg_zone_size;
565  const double big_zone_size = size_ratio * small_zone_size;
566 
567  for (int i = 0; i < size.Size(); i++)
568  {
569  const double val = size(i);
570  const double a = (big_zone_size - small_zone_size) / small_zone_size;
571  size(i) = big_zone_size / (1.0+a*val);
572  }
573 
574  DiffuseField(size, 2);
575  DiffuseField(aspr, 2);
576 
577  tc->SetParDiscreteTargetSize(size);
579  target_c = tc;
580  break;
581  }
582  case 7: // Discrete aspect ratio 3D
583  {
585  DiscreteAdaptTC *tc = new DiscreteAdaptTC(target_t);
586  if (adapt_eval == 0)
587  {
588  tc->SetAdaptivityEvaluator(new AdvectorCG(al));
589  }
590  else
591  {
592 #ifdef MFEM_USE_GSLIB
594 #else
595  MFEM_ABORT("MFEM is not built with GSLIB.");
596 #endif
597  }
599  aspr3d.ProjectCoefficient(fd_aspr3d);
601  target_c = tc;
602  break;
603  }
604  case 8: // shape/size + orientation 2D
605  {
607  DiscreteAdaptTC *tc = new DiscreteAdaptTC(target_t);
608  if (adapt_eval == 0)
609  {
610  tc->SetAdaptivityEvaluator(new AdvectorCG(al));
611  }
612  else
613  {
614 #ifdef MFEM_USE_GSLIB
616 #else
617  MFEM_ABORT("MFEM is not built with GSLIB.");
618 #endif
619  }
620 
621  if (metric_id == 14 || metric_id == 36)
622  {
623  ConstantCoefficient ind_coeff(0.1*0.1);
624  size.ProjectCoefficient(ind_coeff);
625  tc->SetParDiscreteTargetSize(size);
626  }
627 
628  if (metric_id == 85)
629  {
631  aspr.ProjectCoefficient(aspr_coeff);
632  DiffuseField(aspr,2);
634  }
635 
637  ori.ProjectCoefficient(ori_coeff);
639  target_c = tc;
640  break;
641  }
642  default:
643  if (myid == 0) { cout << "Unknown target_id: " << target_id << endl; }
644  return 3;
645  }
646 
647  if (target_c == NULL)
648  {
649  target_c = new TargetConstructor(target_t, MPI_COMM_WORLD);
650  }
651  target_c->SetNodes(x0);
652  TMOP_Integrator *he_nlf_integ= new TMOP_Integrator(metric, target_c);
653 
654  // Finite differences for computations of derivatives.
655  if (fdscheme)
656  {
657  MFEM_VERIFY(pa == false, "PA for finite differences is not implemented.");
658 
659  he_nlf_integ->EnableFiniteDifferences(x);
660  }
661  he_nlf_integ->SetExactActionFlag(exactaction);
662 
663  // Setup the quadrature rules for the TMOP integrator.
664  IntegrationRules *irules = NULL;
665  switch (quad_type)
666  {
667  case 1: irules = &IntRulesLo; break;
668  case 2: irules = &IntRules; break;
669  case 3: irules = &IntRulesCU; break;
670  default:
671  if (myid == 0) { cout << "Unknown quad_type: " << quad_type << endl; }
672  return 3;
673  }
674  he_nlf_integ->SetIntegrationRules(*irules, quad_order);
675  if (myid == 0 && dim == 2)
676  {
677  cout << "Triangle quadrature points: "
678  << irules->Get(Geometry::TRIANGLE, quad_order).GetNPoints()
679  << "\nQuadrilateral quadrature points: "
680  << irules->Get(Geometry::SQUARE, quad_order).GetNPoints() << endl;
681  }
682  if (myid == 0 && dim == 3)
683  {
684  cout << "Tetrahedron quadrature points: "
685  << irules->Get(Geometry::TETRAHEDRON, quad_order).GetNPoints()
686  << "\nHexahedron quadrature points: "
687  << irules->Get(Geometry::CUBE, quad_order).GetNPoints()
688  << "\nPrism quadrature points: "
689  << irules->Get(Geometry::PRISM, quad_order).GetNPoints() << endl;
690  }
691 
692  if (normalization) { he_nlf_integ->ParEnableNormalization(x0); }
693 
694  // Limit the node movement.
695  // The limiting distances can be given by a general function of space.
696  ParFiniteElementSpace dist_pfespace(pmesh, fec); // scalar space
697  ParGridFunction dist(&dist_pfespace);
698  dist = 1.0;
699  // The small_phys_size is relevant only with proper normalization.
700  if (normalization) { dist = small_phys_size; }
701  ConstantCoefficient lim_coeff(lim_const);
702  if (lim_const != 0.0) { he_nlf_integ->EnableLimiting(x0, dist, lim_coeff); }
703 
704  // Adaptive limiting.
705  ParGridFunction zeta_0(&ind_fes);
706  ConstantCoefficient coef_zeta(adapt_lim_const);
707  AdaptivityEvaluator *adapt_evaluator = NULL;
708  if (adapt_lim_const > 0.0)
709  {
710  MFEM_VERIFY(pa == false, "PA is not implemented for adaptive limiting");
711 
713  zeta_0.ProjectCoefficient(alim_coeff);
714 
715  if (adapt_eval == 0) { adapt_evaluator = new AdvectorCG(al); }
716  else if (adapt_eval == 1)
717  {
718 #ifdef MFEM_USE_GSLIB
719  adapt_evaluator = new InterpolatorFP;
720 #else
721  MFEM_ABORT("MFEM is not built with GSLIB support!");
722 #endif
723  }
724  else { MFEM_ABORT("Bad interpolation option."); }
725 
726  he_nlf_integ->EnableAdaptiveLimiting(zeta_0, coef_zeta, *adapt_evaluator);
727  if (visualization)
728  {
729  socketstream vis1;
730  common::VisualizeField(vis1, "localhost", 19916, zeta_0, "Zeta 0",
731  300, 600, 300, 300);
732  }
733  }
734 
735  // 13. Setup the final NonlinearForm (which defines the integral of interest,
736  // its first and second derivatives). Here we can use a combination of
737  // metrics, i.e., optimize the sum of two integrals, where both are
738  // scaled by used-defined space-dependent weights. Note that there are
739  // no command-line options for the weights and the type of the second
740  // metric; one should update those in the code.
741  ParNonlinearForm a(pfespace);
743  ConstantCoefficient *coeff1 = NULL;
744  TMOP_QualityMetric *metric2 = NULL;
745  TargetConstructor *target_c2 = NULL;
747 
748  // Explicit combination of metrics.
749  if (combomet > 0)
750  {
751  // First metric.
752  coeff1 = new ConstantCoefficient(1.0);
753  he_nlf_integ->SetCoefficient(*coeff1);
754 
755  // Second metric.
756  metric2 = new TMOP_Metric_077;
757  TMOP_Integrator *he_nlf_integ2 = NULL;
758  if (combomet == 1)
759  {
760  target_c2 = new TargetConstructor(
762  target_c2->SetVolumeScale(0.01);
763  target_c2->SetNodes(x0);
764  he_nlf_integ2 = new TMOP_Integrator(metric2, target_c2);
765  he_nlf_integ2->SetCoefficient(coeff2);
766  }
767  else { he_nlf_integ2 = new TMOP_Integrator(metric2, target_c); }
768  he_nlf_integ2->SetIntegrationRules(*irules, quad_order);
769  if (fdscheme) { he_nlf_integ2->EnableFiniteDifferences(x); }
770  he_nlf_integ2->SetExactActionFlag(exactaction);
771 
773  combo->AddTMOPIntegrator(he_nlf_integ);
774  combo->AddTMOPIntegrator(he_nlf_integ2);
775  if (normalization) { combo->ParEnableNormalization(x0); }
776  if (lim_const != 0.0) { combo->EnableLimiting(x0, dist, lim_coeff); }
777 
778  a.AddDomainIntegrator(combo);
779  }
780  else { a.AddDomainIntegrator(he_nlf_integ); }
781 
782  if (pa) { a.Setup(); }
783 
784  // Compute the minimum det(J) of the starting mesh.
785  tauval = infinity();
786  const int NE = pmesh->GetNE();
787  for (int i = 0; i < NE; i++)
788  {
789  const IntegrationRule &ir =
790  irules->Get(pfespace->GetFE(i)->GetGeomType(), quad_order);
792  for (int j = 0; j < ir.GetNPoints(); j++)
793  {
794  transf->SetIntPoint(&ir.IntPoint(j));
795  tauval = min(tauval, transf->Jacobian().Det());
796  }
797  }
798  double minJ0;
799  MPI_Allreduce(&tauval, &minJ0, 1, MPI_DOUBLE, MPI_MIN, MPI_COMM_WORLD);
800  tauval = minJ0;
801  if (myid == 0)
802  { cout << "Minimum det(J) of the original mesh is " << tauval << endl; }
803 
804  if (tauval < 0.0 && metric_id != 22 && metric_id != 211 && metric_id != 252
805  && metric_id != 311 && metric_id != 313 && metric_id != 352)
806  {
807  MFEM_ABORT("The input mesh is inverted! Try an untangling metric.");
808  }
809  if (tauval < 0.0)
810  {
811  MFEM_VERIFY(target_t == TargetConstructor::IDEAL_SHAPE_UNIT_SIZE,
812  "Untangling is supported only for ideal targets.");
813 
814  const DenseMatrix &Wideal =
816  tauval /= Wideal.Det();
817 
818  double h0min = h0.Min(), h0min_all;
819  MPI_Allreduce(&h0min, &h0min_all, 1, MPI_DOUBLE, MPI_MIN, MPI_COMM_WORLD);
820  // Slightly below minJ0 to avoid div by 0.
821  tauval -= 0.01 * h0min_all;
822  }
823 
824  const double init_energy = a.GetParGridFunctionEnergy(x);
825 
826  // Visualize the starting mesh and metric values.
827  // Note that for combinations of metrics, this only shows the first metric.
828  if (visualization)
829  {
830  char title[] = "Initial metric values";
831  vis_tmop_metric_p(mesh_poly_deg, *metric, *target_c, *pmesh, title, 0);
832  }
833 
834  // 14. Fix all boundary nodes, or fix only a given component depending on the
835  // boundary attributes of the given mesh. Attributes 1/2/3 correspond to
836  // fixed x/y/z components of the node. Attribute 4 corresponds to an
837  // entirely fixed node. Other boundary attributes do not affect the node
838  // movement boundary conditions.
839  if (move_bnd == false)
840  {
841  Array<int> ess_bdr(pmesh->bdr_attributes.Max());
842  ess_bdr = 1;
843  a.SetEssentialBC(ess_bdr);
844  }
845  else
846  {
847  int n = 0;
848  for (int i = 0; i < pmesh->GetNBE(); i++)
849  {
850  const int nd = pfespace->GetBE(i)->GetDof();
851  const int attr = pmesh->GetBdrElement(i)->GetAttribute();
852  MFEM_VERIFY(!(dim == 2 && attr == 3),
853  "Boundary attribute 3 must be used only for 3D meshes. "
854  "Adjust the attributes (1/2/3/4 for fixed x/y/z/all "
855  "components, rest for free nodes), or use -fix-bnd.");
856  if (attr == 1 || attr == 2 || attr == 3) { n += nd; }
857  if (attr == 4) { n += nd * dim; }
858  }
859  Array<int> ess_vdofs(n), vdofs;
860  n = 0;
861  for (int i = 0; i < pmesh->GetNBE(); i++)
862  {
863  const int nd = pfespace->GetBE(i)->GetDof();
864  const int attr = pmesh->GetBdrElement(i)->GetAttribute();
865  pfespace->GetBdrElementVDofs(i, vdofs);
866  if (attr == 1) // Fix x components.
867  {
868  for (int j = 0; j < nd; j++)
869  { ess_vdofs[n++] = vdofs[j]; }
870  }
871  else if (attr == 2) // Fix y components.
872  {
873  for (int j = 0; j < nd; j++)
874  { ess_vdofs[n++] = vdofs[j+nd]; }
875  }
876  else if (attr == 3) // Fix z components.
877  {
878  for (int j = 0; j < nd; j++)
879  { ess_vdofs[n++] = vdofs[j+2*nd]; }
880  }
881  else if (attr == 4) // Fix all components.
882  {
883  for (int j = 0; j < vdofs.Size(); j++)
884  { ess_vdofs[n++] = vdofs[j]; }
885  }
886  }
887  a.SetEssentialVDofs(ess_vdofs);
888  }
889 
890  // 15. As we use the Newton method to solve the resulting nonlinear system,
891  // here we setup the linear solver for the system's Jacobian.
892  Solver *S = NULL, *S_prec = NULL;
893  const double linsol_rtol = 1e-12;
894  if (lin_solver == 0)
895  {
896  S = new DSmoother(1, 1.0, max_lin_iter);
897  }
898  else if (lin_solver == 1)
899  {
900  CGSolver *cg = new CGSolver(MPI_COMM_WORLD);
901  cg->SetMaxIter(max_lin_iter);
902  cg->SetRelTol(linsol_rtol);
903  cg->SetAbsTol(0.0);
904  cg->SetPrintLevel(verbosity_level >= 2 ? 3 : -1);
905  S = cg;
906  }
907  else
908  {
909  MINRESSolver *minres = new MINRESSolver(MPI_COMM_WORLD);
910  minres->SetMaxIter(max_lin_iter);
911  minres->SetRelTol(linsol_rtol);
912  minres->SetAbsTol(0.0);
913  if (verbosity_level > 2) { minres->SetPrintLevel(1); }
914  else { minres->SetPrintLevel(verbosity_level == 2 ? 3 : -1); }
915  if (lin_solver == 3 || lin_solver == 4)
916  {
917  if (pa)
918  {
919  MFEM_VERIFY(lin_solver != 4, "PA l1-Jacobi is not implemented");
920  S_prec = new OperatorJacobiSmoother;
921  }
922  else
923  {
924  HypreSmoother *hs = new HypreSmoother;
925  hs->SetType((lin_solver == 3) ? HypreSmoother::Jacobi
927  S_prec = hs;
928  }
929  minres->SetPreconditioner(*S_prec);
930  }
931  S = minres;
932  }
933 
934  // Perform the nonlinear optimization.
935  const IntegrationRule &ir =
936  irules->Get(pfespace->GetFE(0)->GetGeomType(), quad_order);
937  TMOPNewtonSolver solver(pfespace->GetComm(), ir, solver_type);
938  // Provide all integration rules in case of a mixed mesh.
939  solver.SetIntegrationRules(*irules, quad_order);
940  if (solver_type == 0)
941  {
942  // Specify linear solver when we use a Newton-based solver.
943  solver.SetPreconditioner(*S);
944  }
945  // For untangling, the solver will update the min det(T) values.
946  if (tauval < 0.0) { solver.SetMinDetPtr(&tauval); }
947  solver.SetMaxIter(solver_iter);
948  solver.SetRelTol(solver_rtol);
949  solver.SetAbsTol(0.0);
950  if (solver_art_type > 0)
951  {
952  solver.SetAdaptiveLinRtol(solver_art_type, 0.5, 0.9);
953  }
954  solver.SetPrintLevel(verbosity_level >= 1 ? 1 : -1);
955  solver.SetOperator(a);
956  solver.Mult(b, x.GetTrueVector());
957  x.SetFromTrueVector();
958  if (myid == 0 && solver.GetConverged() == false)
959  {
960  cout << "Nonlinear solver: rtol = " << solver_rtol << " not achieved.\n";
961  }
962 
963  // 16. Save the optimized mesh to a file. This output can be viewed later
964  // using GLVis: "glvis -m optimized -np num_mpi_tasks".
965  {
966  ostringstream mesh_name;
967  mesh_name << "optimized.mesh";
968  ofstream mesh_ofs(mesh_name.str().c_str());
969  mesh_ofs.precision(8);
970  pmesh->PrintAsOne(mesh_ofs);
971  }
972 
973  // 17. Compute the amount of energy decrease.
974  const double fin_energy = a.GetParGridFunctionEnergy(x);
975  double metric_part = fin_energy;
976  if (lim_const > 0.0 || adapt_lim_const > 0.0)
977  {
978  lim_coeff.constant = 0.0;
979  coef_zeta.constant = 0.0;
980  metric_part = a.GetParGridFunctionEnergy(x);
981  lim_coeff.constant = lim_const;
982  coef_zeta.constant = adapt_lim_const;
983  }
984  if (myid == 0)
985  {
986  cout << "Initial strain energy: " << init_energy
987  << " = metrics: " << init_energy
988  << " + limiting term: " << 0.0 << endl;
989  cout << " Final strain energy: " << fin_energy
990  << " = metrics: " << metric_part
991  << " + limiting term: " << fin_energy - metric_part << endl;
992  cout << "The strain energy decreased by: " << setprecision(12)
993  << (init_energy - fin_energy) * 100.0 / init_energy << " %." << endl;
994  }
995 
996  // 18. Visualize the final mesh and metric values.
997  if (visualization)
998  {
999  char title[] = "Final metric values";
1000  vis_tmop_metric_p(mesh_poly_deg, *metric, *target_c, *pmesh, title, 600);
1001  }
1002 
1003  if (adapt_lim_const > 0.0 && visualization)
1004  {
1005  socketstream vis0;
1006  common::VisualizeField(vis0, "localhost", 19916, zeta_0, "Xi 0",
1007  600, 600, 300, 300);
1008  }
1009 
1010  // 19. Visualize the mesh displacement.
1011  if (visualization)
1012  {
1013  x0 -= x;
1014  socketstream sock;
1015  if (myid == 0)
1016  {
1017  sock.open("localhost", 19916);
1018  sock << "solution\n";
1019  }
1020  pmesh->PrintAsOne(sock);
1021  x0.SaveAsOne(sock);
1022  if (myid == 0)
1023  {
1024  sock << "window_title 'Displacements'\n"
1025  << "window_geometry "
1026  << 1200 << " " << 0 << " " << 600 << " " << 600 << "\n"
1027  << "keys jRmclA" << endl;
1028  }
1029  }
1030 
1031  // 20. Free the used memory.
1032  delete S;
1033  delete S_prec;
1034  delete target_c2;
1035  delete metric2;
1036  delete coeff1;
1037  delete adapt_evaluator;
1038  delete target_c;
1039  delete adapt_coeff;
1040  delete metric;
1041  delete pfespace;
1042  delete fec;
1043  delete pmesh;
1044 
1045  MPI_Finalize();
1046  return 0;
1047 }
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:705
void discrete_aspr_3d(const Vector &x, Vector &v)
virtual void SetParDiscreteTargetSize(const ParGridFunction &tspec_)
Definition: tmop.cpp:1347
int Size() const
Return the logical size of the array.
Definition: array.hpp:134
Conjugate gradient method.
Definition: solvers.hpp:316
Definition: ex25.cpp:148
int GetNDofs() const
Returns number of degrees of freedom.
Definition: fespace.hpp:537
AssemblyLevel
Enumeration defining the assembly level for bilinear and nonlinear form classes derived from Operator...
Class for an integration rule - an Array of IntegrationPoint.
Definition: intrules.hpp:90
void EnableFiniteDifferences(const GridFunction &x)
Enables FD-based approximation and computes dx.
Definition: tmop.cpp:3070
int DofToVDof(int dof, int vd, int ndofs=-1) const
Definition: fespace.cpp:233
double discrete_size_2d(const Vector &x)
3D barrier Shape+Size (VS) metric.
Definition: tmop.hpp:578
Data type for scaled Jacobi-type smoother of sparse matrix.
const IntegrationRule & Get(int GeomType, int Order)
Returns an integration rule for given GeomType and Order.
Definition: intrules.cpp:920
IntegrationRules IntRulesLo(0, Quadrature1D::GaussLobatto)
void SetIntegrationRules(IntegrationRules &irules, int order)
Prescribe a set of integration rules; relevant for mixed meshes.
Definition: tmop_tools.hpp:172
void SetFromTrueVector()
Shortcut for calling SetFromTrueDofs() with GetTrueVector() as argument.
Definition: gridfunc.hpp:143
A coefficient that is constant across space and time.
Definition: coefficient.hpp:78
long GetGlobalNE() const
Return the total (global) number of elements.
Definition: mesh.hpp:872
int GetNBE() const
Returns number of boundary elements.
Definition: mesh.hpp:849
3D barrier Shape (S) metric.
Definition: tmop.hpp:484
void SetAdaptivityEvaluator(AdaptivityEvaluator *ae)
Definition: tmop.hpp:1158
virtual void Setup()
Setup the NonlinearForm: based on the current AssemblyLevel and the current mesh, optionally...
double Det() const
Definition: densemat.cpp:451
2D barrier Shape+Size (VS) metric (polyconvex).
Definition: tmop.hpp:358
void AddDomainIntegrator(NonlinearFormIntegrator *nlfi)
Adds new Domain Integrator.
void SetIntPoint(const IntegrationPoint *ip)
Set the integration point ip that weights and Jacobians will be evaluated at.
Definition: eltrans.hpp:85
Container class for integration rules.
Definition: intrules.hpp:311
Data type dense matrix using column-major storage.
Definition: densemat.hpp:23
3D barrier Shape (S) metric.
Definition: tmop.hpp:466
void SetVolumeScale(double vol_scale)
Used by target type IDEAL_SHAPE_EQUAL_SIZE. The default volume scale is 1.
Definition: tmop.hpp:950
Parallel non-linear operator on the true dofs.
2D barrier shape (S) metric (polyconvex).
Definition: tmop.hpp:186
void EnableAdaptiveLimiting(const GridFunction &z0, Coefficient &coeff, AdaptivityEvaluator &ae)
Restriction of the node positions to certain regions.
Definition: tmop.cpp:2239
3D non-barrier Size (V) metric.
Definition: tmop.hpp:542
2D barrier size (V) metric (polyconvex).
Definition: tmop.hpp:340
void SetIntegrationRules(IntegrationRules &irules, int order)
Prescribe a set of integration rules; relevant for mixed meshes.
Definition: tmop.hpp:1418
int GetNE() const
Returns number of elements.
Definition: mesh.hpp:846
2D barrier Shape+Size (VS) metric (polyconvex).
Definition: tmop.hpp:734
void Print(std::ostream &out=mfem::out)
Print the configuration of the MFEM virtual device object.
Definition: device.cpp:279
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:493
void SetExactActionFlag(bool flag_)
Flag to control if exact action of Integration is effected.
Definition: tmop.hpp:1523
MINRES method.
Definition: solvers.hpp:426
void Randomize(int seed=0)
Set random values in the vector.
Definition: vector.cpp:763
virtual int OrderJ() const =0
Return the order of the elements of the Jacobian of the transformation.
2D barrier Shape+Size (VS) metric (not polyconvex).
Definition: tmop.hpp:204
Geometry::Type GetElementBaseGeometry(int i) const
Definition: mesh.hpp:971
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:3098
2D Shifted barrier form of shape metric (mu_2).
Definition: tmop.hpp:253
const DenseMatrix & GetGeomToPerfGeomJac(int GeomType) const
Definition: geom.hpp:97
2D barrier shape (S) metric (not polyconvex).
Definition: tmop.hpp:323
Geometry::Type GetGeomType() const
Returns the Geometry::Type of the reference element.
Definition: fe.hpp:320
Abstract class for local mesh quality metrics in the target-matrix optimization paradigm (TMOP) by P...
Definition: tmop.hpp:23
Geometry Geometries
Definition: fe.cpp:13507
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:137
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
Definition: solvers.hpp:439
void SetPrintLevel(int print_lvl)
Definition: solvers.cpp:71
MPI_Comm GetComm() const
Definition: pfespace.hpp:263
double GetElementSize(ElementTransformation *T, int type=0)
Definition: mesh.cpp:76
virtual const FiniteElement * GetFE(int i) const
Definition: pfespace.cpp:511
void EnableLimiting(const GridFunction &n0, const GridFunction &dist, Coefficient &w0, TMOP_LimiterFunction *lfunc=NULL)
Limiting of the mesh displacements (general version).
Definition: tmop.cpp:2214
int GetNBE() const
Returns number of boundary elements in the mesh.
Definition: fespace.hpp:576
void Parse()
Parse the command-line options. Note that this function expects all the options provided through the ...
Definition: optparser.cpp:150
void SetCoefficient(Coefficient &w1)
Sets a scaling Coefficient for the quality metric term of the integrator.
Definition: tmop.hpp:1430
2D non-barrier size (V) metric (not polyconvex).
Definition: tmop.hpp:288
3D barrier Size (V) metric.
Definition: tmop.hpp:560
Jacobi smoothing for a given bilinear form (no matrix necessary).
Definition: solvers.hpp:130
double b
Definition: lissajous.cpp:42
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:9143
const DenseMatrix & Jacobian()
Return the Jacobian matrix of the transformation at the currently set IntegrationPoint, using the method SetIntPoint().
Definition: eltrans.hpp:111
void SetNodes(const GridFunction &n)
Set the nodes to be used in the target-matrix construction.
Definition: tmop.hpp:944
void SetMaxIter(int max_it)
Definition: solvers.hpp:100
T Max() const
Find the maximal element in the array, using the comparison operator &lt; for class T.
Definition: array.cpp:68
Version of QuadraticFECollection with positive basis functions.
Definition: fe_coll.hpp:579
virtual void SetParDiscreteTargetAspectRatio(const ParGridFunction &tspec_)
Definition: tmop.cpp:1363
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:123
Parallel smoothers in hypre.
Definition: hypre.hpp:840
void ParEnableNormalization(const ParGridFunction &x)
Definition: tmop.cpp:3187
2D barrier Shape+Size+Orientation (VOS) metric (polyconvex).
Definition: tmop.hpp:396
void SetNodalFESpace(FiniteElementSpace *nfes)
Definition: mesh.cpp:4843
int Dimension() const
Definition: mesh.hpp:911
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:457
void DiffuseField(ParGridFunction &field, int smooth_steps)
Definition: dist_solver.cpp:17
2D barrier Shape+Orientation (OS) metric (polyconvex).
Definition: tmop.hpp:716
double GetParGridFunctionEnergy(const Vector &x) const
Compute the energy corresponding to the state x.
const Vector & GetTrueVector() const
Read only access to the (optional) internal true-dof Vector.
Definition: gridfunc.hpp:127
A general vector function coefficient.
void SetAbsTol(double atol)
Definition: solvers.hpp:99
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
Definition: mesh.hpp:204
void SetRelTol(double rtol)
Definition: solvers.hpp:98
2D barrier (not a shape) metric (polyconvex).
Definition: tmop.hpp:272
int GetDof() const
Returns the number of degrees of freedom in the finite element.
Definition: fe.hpp:323
double weight_fun(const Vector &x)
2D barrier size (V) metric (polyconvex).
Definition: tmop.hpp:305
Collection of finite elements from the same family in multiple dimensions. This class is used to matc...
Definition: fe_coll.hpp:26
void AddOption(bool *var, const char *enable_short_name, const char *enable_long_name, const char *disable_short_name, const char *disable_long_name, const char *description, bool required=false)
Add a boolean option and set &#39;var&#39; to receive the value. Enable/disable tags are used to set the bool...
Definition: optparser.hpp:82
double adapt_lim_fun(const Vector &x)
void AddTMOPIntegrator(TMOP_Integrator *ti)
Adds a new TMOP_Integrator to the combination.
Definition: tmop.hpp:1541
void SaveAsOne(const char *fname, int precision=16) const
Definition: pgridfunc.cpp:867
double material_indicator_2d(const Vector &x)
2D barrier Shape+Size (VS) metric (not polyconvex).
Definition: tmop.hpp:222
double a
Definition: lissajous.cpp:41
virtual void GetElementDofs(int i, Array< int > &dofs) const
Returns indexes of degrees of freedom in array dofs for i&#39;th element.
Definition: pfespace.cpp:467
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:347
2D barrier Shape+Size+Orientation (VOS) metric (polyconvex).
Definition: tmop.hpp:698
3D Shape (S) metric, untangling version of 303.
Definition: tmop.hpp:521
int dim
Definition: ex24.cpp:53
IntegrationRules IntRulesCU(0, Quadrature1D::ClosedUniform)
void PrintOptions(std::ostream &out) const
Print the options.
Definition: optparser.cpp:327
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;.
A general function coefficient.
Vector data type.
Definition: vector.hpp:60
void VisualizeField(socketstream &sock, const char *vishost, int visport, GridFunction &gf, const char *title, int x, int y, int w, int h, const char *keys, bool vec)
Definition: fem_extras.cpp:94
void GetNodes(Vector &node_coord) const
Definition: mesh.cpp:7343
void PrintAsOne(std::ostream &out=mfem::out) const
Definition: pmesh.cpp:4527
Arbitrary order H1-conforming (continuous) finite elements.
Definition: fe_coll.hpp:216
TargetType
Target-matrix construction algorithms supported by this class.
Definition: tmop.hpp:872
void SetEssentialVDofs(const Array< int > &ess_vdofs_list)
Specify essential boundary conditions.
virtual void SetParDiscreteTargetOrientation(const ParGridFunction &tspec_)
Definition: tmop.cpp:1372
Base class for solvers.
Definition: operator.hpp:648
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:2686
void GetBdrElementVDofs(int i, Array< int > &vdofs) const
Returns indexes of degrees of freedom for i&#39;th boundary element.
Definition: fespace.cpp:268
void SetNodalGridFunction(GridFunction *nodes, bool make_owner=false)
Definition: mesh.cpp:4871
Base class representing target-matrix construction algorithms for mesh optimization via the target-ma...
Definition: tmop.hpp:868
The MFEM Device class abstracts hardware devices such as GPUs, as well as programming models such as ...
Definition: device.hpp:121
virtual void SetAnalyticTargetSpec(Coefficient *sspec, VectorCoefficient *vspec, TMOPMatrixCoefficient *mspec)
Definition: tmop.cpp:1227
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:377
void SetType(HypreSmoother::Type type, int relax_times=1)
Set the relaxation type and number of sweeps.
Definition: hypre.cpp:3028
2D barrier Shape+Orientation (OS) metric (polyconvex).
Definition: tmop.hpp:381
virtual double * HostReadWrite()
Shortcut for mfem::ReadWrite(vec.GetMemory(), vec.Size(), false).
Definition: vector.hpp:446
void ParEnableNormalization(const ParGridFunction &x)
Definition: tmop.cpp:2946
int main()
double GetElementVolume(int i)
Definition: mesh.cpp:112
const Element * GetBdrElement(int i) const
Definition: mesh.hpp:946
3D barrier Shape (S) metric.
Definition: tmop.hpp:450
2D non-barrier metric without a type.
Definition: tmop.hpp:108
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:1198
2D non-barrier Shape+Size+Orientation (VOS) metric (polyconvex).
Definition: tmop.hpp:238