MFEM  v4.2.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-2020, 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 -ls 2 -li 100 -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 -ls 2 -li 100 -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 -ls 2 -li 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 7 -tid 5 -ni 200 -ls 2 -li 100 -bnd -qt 1 -qo 8
44 // Adapted discrete size; explicit combo of metrics; mixed tri/quad mesh:
45 // 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
46 // Adapted discrete size+aspect_ratio:
47 // mpirun -np 4 pmesh-optimizer -m square01.mesh -o 2 -rs 2 -mid 7 -tid 6 -ni 100
48 // mpirun -np 4 pmesh-optimizer -m square01.mesh -o 2 -rs 2 -mid 7 -tid 6 -ni 100 -qo 6 -ex -st 1 -nor
49 // Adapted discrete size+orientation (requires GSLIB):
50 // * mpirun -np 4 pmesh-optimizer -m square01.mesh -o 2 -rs 2 -mid 14 -tid 8 -ni 100 -ls 2 -li 100 -bnd -qt 1 -qo 8 -fd -ae 1
51 // Adapted discrete aspect_ratio+orientation (requires GSLIB):
52 // * mpirun -np 4 pmesh-optimizer -m square01.mesh -o 2 -rs 2 -mid 85 -tid 8 -ni 10 -ls 2 -li 100 -bnd -qt 1 -qo 8 -fd -ae 1
53 // Adapted discrete aspect ratio (3D):
54 // mpirun -np 4 pmesh-optimizer -m cube.mesh -o 2 -rs 2 -mid 302 -tid 7 -ni 20 -ls 2 -li 100 -bnd -qt 1 -qo 8
55 //
56 // Adaptive limiting:
57 // mpirun -np 4 pmesh-optimizer -m stretched2D.mesh -o 2 -mid 2 -tid 1 -ni 50 -qo 5 -nor -vl 1 -alc 0.5
58 // Adaptive limiting through the L-BFGS solver:
59 // 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
60 // Adaptive limiting through FD (requires GSLIB):
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 -fd -ae 1
62 //
63 // Blade shape:
64 // mpirun -np 4 pmesh-optimizer -m blade.mesh -o 4 -rs 0 -mid 2 -tid 1 -ni 200 -ls 2 -li 100 -bnd -qt 1 -qo 8
65 // Blade shape with FD-based solver:
66 // mpirun -np 4 pmesh-optimizer -m blade.mesh -o 4 -rs 0 -mid 2 -tid 1 -ni 200 -ls 2 -li 100 -bnd -qt 1 -qo 8 -fd
67 // Blade limited shape:
68 // mpirun -np 4 pmesh-optimizer -m blade.mesh -o 4 -rs 0 -mid 2 -tid 1 -ni 200 -ls 2 -li 100 -bnd -qt 1 -qo 8 -lc 5000
69 // ICF shape and equal size:
70 // mpirun -np 4 pmesh-optimizer -o 3 -rs 0 -mid 9 -tid 2 -ni 200 -ls 2 -li 100 -bnd -qt 1 -qo 8
71 // ICF shape and initial size:
72 // mpirun -np 4 pmesh-optimizer -o 3 -rs 0 -mid 9 -tid 3 -ni 100 -ls 2 -li 100 -bnd -qt 1 -qo 8
73 // ICF shape:
74 // mpirun -np 4 pmesh-optimizer -o 3 -rs 0 -mid 1 -tid 1 -ni 100 -ls 2 -li 100 -bnd -qt 1 -qo 8
75 // ICF limited shape:
76 // mpirun -np 4 pmesh-optimizer -o 3 -rs 0 -mid 1 -tid 1 -ni 100 -ls 2 -li 100 -bnd -qt 1 -qo 8 -lc 10
77 // ICF combo shape + size (rings, slow convergence):
78 // mpirun -np 4 pmesh-optimizer -o 3 -rs 0 -mid 1 -tid 1 -ni 1000 -ls 2 -li 100 -bnd -qt 1 -qo 8 -cmb 1
79 // Mixed tet / cube / hex mesh with limiting:
80 // 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
81 // 3D pinched sphere shape (the mesh is in the mfem/data GitHub repository):
82 // * mpirun -np 4 pmesh-optimizer -m ../../../mfem_data/ball-pert.mesh -o 4 -rs 0 -mid 303 -tid 1 -ni 20 -ls 2 -li 500 -fix-bnd
83 // 2D non-conforming shape and equal size:
84 // mpirun -np 4 pmesh-optimizer -m ./amr-quad-q2.mesh -o 2 -rs 1 -mid 9 -tid 2 -ni 200 -ls 2 -li 100 -bnd -qt 1 -qo 8
85 
86 #include "mfem.hpp"
87 #include "../common/mfem-common.hpp"
88 #include <iostream>
89 #include <fstream>
90 #include "mesh-optimizer.hpp"
91 
92 using namespace mfem;
93 using namespace std;
94 
95 int main (int argc, char *argv[])
96 {
97  // 0. Initialize MPI.
98  int num_procs, myid;
99  MPI_Init(&argc, &argv);
100  MPI_Comm_size(MPI_COMM_WORLD, &num_procs);
101  MPI_Comm_rank(MPI_COMM_WORLD, &myid);
102 
103  // 1. Set the method's default parameters.
104  const char *mesh_file = "icf.mesh";
105  int mesh_poly_deg = 1;
106  int rs_levels = 0;
107  int rp_levels = 0;
108  double jitter = 0.0;
109  int metric_id = 1;
110  int target_id = 1;
111  double lim_const = 0.0;
112  double adapt_lim_const = 0.0;
113  int quad_type = 1;
114  int quad_order = 8;
115  int solver_type = 0;
116  int solver_iter = 10;
117  double solver_rtol = 1e-10;
118  int lin_solver = 2;
119  int max_lin_iter = 100;
120  bool move_bnd = true;
121  int combomet = 0;
122  bool normalization = false;
123  bool visualization = true;
124  int verbosity_level = 0;
125  bool fdscheme = false;
126  int adapt_eval = 0;
127  bool exactaction = false;
128 
129  // 2. Parse command-line options.
130  OptionsParser args(argc, argv);
131  args.AddOption(&mesh_file, "-m", "--mesh",
132  "Mesh file to use.");
133  args.AddOption(&mesh_poly_deg, "-o", "--order",
134  "Polynomial degree of mesh finite element space.");
135  args.AddOption(&rs_levels, "-rs", "--refine-serial",
136  "Number of times to refine the mesh uniformly in serial.");
137  args.AddOption(&rp_levels, "-rp", "--refine-parallel",
138  "Number of times to refine the mesh uniformly in parallel.");
139  args.AddOption(&jitter, "-ji", "--jitter",
140  "Random perturbation scaling factor.");
141  args.AddOption(&metric_id, "-mid", "--metric-id",
142  "Mesh optimization metric:\n\t"
143  "1 : |T|^2 -- 2D shape\n\t"
144  "2 : 0.5|T|^2/tau-1 -- 2D shape (condition number)\n\t"
145  "7 : |T-T^-t|^2 -- 2D shape+size\n\t"
146  "9 : tau*|T-T^-t|^2 -- 2D shape+size\n\t"
147  "22 : 0.5(|T|^2-2*tau)/(tau-tau_0) -- 2D untangling\n\t"
148  "50 : 0.5|T^tT|^2/tau^2-1 -- 2D shape\n\t"
149  "55 : (tau-1)^2 -- 2D size\n\t"
150  "56 : 0.5(sqrt(tau)-1/sqrt(tau))^2 -- 2D size\n\t"
151  "58 : |T^tT|^2/(tau^2)-2*|T|^2/tau+2 -- 2D shape\n\t"
152  "77 : 0.5(tau-1/tau)^2 -- 2D size\n\t"
153  "211: (tau-1)^2-tau+sqrt(tau^2) -- 2D untangling\n\t"
154  "252: 0.5(tau-1)^2/(tau-tau_0) -- 2D untangling\n\t"
155  "301: (|T||T^-1|)/3-1 -- 3D shape\n\t"
156  "302: (|T|^2|T^-1|^2)/9-1 -- 3D shape\n\t"
157  "303: (|T|^2)/3*tau^(2/3)-1 -- 3D shape\n\t"
158  "315: (tau-1)^2 -- 3D size\n\t"
159  "316: 0.5(sqrt(tau)-1/sqrt(tau))^2 -- 3D size\n\t"
160  "321: |T-T^-t|^2 -- 3D shape+size\n\t"
161  "352: 0.5(tau-1)^2/(tau-tau_0) -- 3D untangling");
162  args.AddOption(&target_id, "-tid", "--target-id",
163  "Target (ideal element) type:\n\t"
164  "1: Ideal shape, unit size\n\t"
165  "2: Ideal shape, equal size\n\t"
166  "3: Ideal shape, initial size\n\t"
167  "4: Given full analytic Jacobian (in physical space)\n\t"
168  "5: Ideal shape, given size (in physical space)");
169  args.AddOption(&lim_const, "-lc", "--limit-const", "Limiting constant.");
170  args.AddOption(&adapt_lim_const, "-alc", "--adapt-limit-const",
171  "Adaptive limiting coefficient constant.");
172  args.AddOption(&quad_type, "-qt", "--quad-type",
173  "Quadrature rule type:\n\t"
174  "1: Gauss-Lobatto\n\t"
175  "2: Gauss-Legendre\n\t"
176  "3: Closed uniform points");
177  args.AddOption(&quad_order, "-qo", "--quad_order",
178  "Order of the quadrature rule.");
179  args.AddOption(&solver_type, "-st", "--solver-type",
180  " Type of solver: (default) 0: Newton, 1: LBFGS");
181  args.AddOption(&solver_iter, "-ni", "--newton-iters",
182  "Maximum number of Newton iterations.");
183  args.AddOption(&solver_rtol, "-rtol", "--newton-rel-tolerance",
184  "Relative tolerance for the Newton solver.");
185  args.AddOption(&lin_solver, "-ls", "--lin-solver",
186  "Linear solver: 0 - l1-Jacobi, 1 - CG, 2 - MINRES.");
187  args.AddOption(&max_lin_iter, "-li", "--lin-iter",
188  "Maximum number of iterations in the linear solve.");
189  args.AddOption(&move_bnd, "-bnd", "--move-boundary", "-fix-bnd",
190  "--fix-boundary",
191  "Enable motion along horizontal and vertical boundaries.");
192  args.AddOption(&combomet, "-cmb", "--combo-type",
193  "Combination of metrics options:"
194  "0: Use single metric\n\t"
195  "1: Shape + space-dependent size given analytically\n\t"
196  "2: Shape + adapted size given discretely; shared target");
197  args.AddOption(&normalization, "-nor", "--normalization", "-no-nor",
198  "--no-normalization",
199  "Make all terms in the optimization functional unitless.");
200  args.AddOption(&fdscheme, "-fd", "--fd_approximation",
201  "-no-fd", "--no-fd-approx",
202  "Enable finite difference based derivative computations.");
203  args.AddOption(&exactaction, "-ex", "--exact_action",
204  "-no-ex", "--no-exact-action",
205  "Enable exact action of TMOP_Integrator.");
206  args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
207  "--no-visualization",
208  "Enable or disable GLVis visualization.");
209  args.AddOption(&verbosity_level, "-vl", "--verbosity-level",
210  "Set the verbosity level - 0, 1, or 2.");
211  args.AddOption(&adapt_eval, "-ae", "--adaptivity-evaluator",
212  "0 - Advection based (DEFAULT), 1 - GSLIB.");
213  args.Parse();
214  if (!args.Good())
215  {
216  if (myid == 0) { args.PrintUsage(cout); }
217  return 1;
218  }
219  if (myid == 0) { args.PrintOptions(cout); }
220 
221  // 3. Initialize and refine the starting mesh.
222  Mesh *mesh = new Mesh(mesh_file, 1, 1, false);
223  for (int lev = 0; lev < rs_levels; lev++)
224  {
225  mesh->UniformRefinement();
226  }
227  const int dim = mesh->Dimension();
228  if (myid == 0)
229  {
230  cout << "Mesh curvature: ";
231  if (mesh->GetNodes()) { cout << mesh->GetNodes()->OwnFEC()->Name(); }
232  else { cout << "(NONE)"; }
233  cout << endl;
234  }
235 
236  ParMesh *pmesh = new ParMesh(MPI_COMM_WORLD, *mesh);
237 
238  delete mesh;
239  for (int lev = 0; lev < rp_levels; lev++)
240  {
241  pmesh->UniformRefinement();
242  }
243 
244  // 4. Define a finite element space on the mesh. Here we use vector finite
245  // elements which are tensor products of quadratic finite elements. The
246  // number of components in the vector finite element space is specified by
247  // the last parameter of the FiniteElementSpace constructor.
249  if (mesh_poly_deg <= 0)
250  {
251  fec = new QuadraticPosFECollection;
252  mesh_poly_deg = 2;
253  }
254  else { fec = new H1_FECollection(mesh_poly_deg, dim); }
255  ParFiniteElementSpace *pfespace = new ParFiniteElementSpace(pmesh, fec, dim);
256 
257  // 5. Make the mesh curved based on the above finite element space. This
258  // means that we define the mesh elements through a fespace-based
259  // transformation of the reference element.
260  pmesh->SetNodalFESpace(pfespace);
261 
262  // 6. Set up an empty right-hand side vector b, which is equivalent to b=0.
263  Vector b(0);
264 
265  // 7. Get the mesh nodes (vertices and other degrees of freedom in the finite
266  // element space) as a finite element grid function in fespace. Note that
267  // changing x automatically changes the shapes of the mesh elements.
268  ParGridFunction x(pfespace);
269  pmesh->SetNodalGridFunction(&x);
270 
271  // 8. Define a vector representing the minimal local mesh size in the mesh
272  // nodes. We index the nodes using the scalar version of the degrees of
273  // freedom in pfespace. Note: this is partition-dependent.
274  //
275  // In addition, compute average mesh size and total volume.
276  Vector h0(pfespace->GetNDofs());
277  h0 = infinity();
278  double vol_loc = 0.0;
279  Array<int> dofs;
280  for (int i = 0; i < pmesh->GetNE(); i++)
281  {
282  // Get the local scalar element degrees of freedom in dofs.
283  pfespace->GetElementDofs(i, dofs);
284  // Adjust the value of h0 in dofs based on the local mesh size.
285  const double hi = pmesh->GetElementSize(i);
286  for (int j = 0; j < dofs.Size(); j++)
287  {
288  h0(dofs[j]) = min(h0(dofs[j]), hi);
289  }
290  vol_loc += pmesh->GetElementVolume(i);
291  }
292  double volume;
293  MPI_Allreduce(&vol_loc, &volume, 1, MPI_DOUBLE, MPI_SUM, MPI_COMM_WORLD);
294  const double small_phys_size = pow(volume, 1.0 / dim) / 100.0;
295 
296  // 9. Add a random perturbation to the nodes in the interior of the domain.
297  // We define a random grid function of fespace and make sure that it is
298  // zero on the boundary and its values are locally of the order of h0.
299  // The latter is based on the DofToVDof() method which maps the scalar to
300  // the vector degrees of freedom in pfespace.
301  ParGridFunction rdm(pfespace);
302  rdm.Randomize();
303  rdm -= 0.25; // Shift to random values in [-0.5,0.5].
304  rdm *= jitter;
305  // Scale the random values to be of order of the local mesh size.
306  for (int i = 0; i < pfespace->GetNDofs(); i++)
307  {
308  for (int d = 0; d < dim; d++)
309  {
310  rdm(pfespace->DofToVDof(i,d)) *= h0(i);
311  }
312  }
313  Array<int> vdofs;
314  for (int i = 0; i < pfespace->GetNBE(); i++)
315  {
316  // Get the vector degrees of freedom in the boundary element.
317  pfespace->GetBdrElementVDofs(i, vdofs);
318  // Set the boundary values to zero.
319  for (int j = 0; j < vdofs.Size(); j++) { rdm(vdofs[j]) = 0.0; }
320  }
321  x -= rdm;
322  // Set the perturbation of all nodes from the true nodes.
323  x.SetTrueVector();
324  x.SetFromTrueVector();
325 
326  // 10. Save the starting (prior to the optimization) mesh to a file. This
327  // output can be viewed later using GLVis: "glvis -m perturbed -np
328  // num_mpi_tasks".
329  {
330  ostringstream mesh_name;
331  mesh_name << "perturbed.mesh";
332  ofstream mesh_ofs(mesh_name.str().c_str());
333  mesh_ofs.precision(8);
334  pmesh->PrintAsOne(mesh_ofs);
335  }
336 
337  // 11. Store the starting (prior to the optimization) positions.
338  ParGridFunction x0(pfespace);
339  x0 = x;
340 
341  // 12. Form the integrator that uses the chosen metric and target.
342  double tauval = -0.1;
343  TMOP_QualityMetric *metric = NULL;
344  switch (metric_id)
345  {
346  case 1: metric = new TMOP_Metric_001; break;
347  case 2: metric = new TMOP_Metric_002; break;
348  case 7: metric = new TMOP_Metric_007; break;
349  case 9: metric = new TMOP_Metric_009; break;
350  case 14: metric = new TMOP_Metric_SSA2D; break;
351  case 22: metric = new TMOP_Metric_022(tauval); break;
352  case 50: metric = new TMOP_Metric_050; break;
353  case 55: metric = new TMOP_Metric_055; break;
354  case 56: metric = new TMOP_Metric_056; break;
355  case 58: metric = new TMOP_Metric_058; break;
356  case 77: metric = new TMOP_Metric_077; break;
357  case 85: metric = new TMOP_Metric_085; break;
358  case 211: metric = new TMOP_Metric_211; break;
359  case 252: metric = new TMOP_Metric_252(tauval); break;
360  case 301: metric = new TMOP_Metric_301; break;
361  case 302: metric = new TMOP_Metric_302; break;
362  case 303: metric = new TMOP_Metric_303; break;
363  case 315: metric = new TMOP_Metric_315; break;
364  case 316: metric = new TMOP_Metric_316; break;
365  case 321: metric = new TMOP_Metric_321; break;
366  case 352: metric = new TMOP_Metric_352(tauval); break;
367  default:
368  if (myid == 0) { cout << "Unknown metric_id: " << metric_id << endl; }
369  return 3;
370  }
372  TargetConstructor *target_c = NULL;
373  HessianCoefficient *adapt_coeff = NULL;
374  H1_FECollection ind_fec(mesh_poly_deg, dim);
375  ParFiniteElementSpace ind_fes(pmesh, &ind_fec);
376  ParFiniteElementSpace ind_fesv(pmesh, &ind_fec, dim);
377  ParGridFunction size(&ind_fes), aspr(&ind_fes), disc(&ind_fes), ori(&ind_fes);
378  ParGridFunction aspr3d(&ind_fesv), size3d(&ind_fesv);
379 
380  switch (target_id)
381  {
382  case 1: target_t = TargetConstructor::IDEAL_SHAPE_UNIT_SIZE; break;
383  case 2: target_t = TargetConstructor::IDEAL_SHAPE_EQUAL_SIZE; break;
384  case 3: target_t = TargetConstructor::IDEAL_SHAPE_GIVEN_SIZE; break;
385  case 4:
386  {
388  AnalyticAdaptTC *tc = new AnalyticAdaptTC(target_t);
389  adapt_coeff = new HessianCoefficient(dim, metric_id);
390  tc->SetAnalyticTargetSpec(NULL, NULL, adapt_coeff);
391  target_c = tc;
392  break;
393  }
394  case 5:
395  {
397  DiscreteAdaptTC *tc = new DiscreteAdaptTC(target_t);
398  if (adapt_eval == 0)
399  {
401  }
402  else
403  {
404 #ifdef MFEM_USE_GSLIB
406 #else
407  MFEM_ABORT("MFEM is not built with GSLIB.");
408 #endif
409  }
411  size.ProjectCoefficient(ind_coeff);
412  tc->SetParDiscreteTargetSize(size);
413  target_c = tc;
414  break;
415  }
416  case 6: //material indicator 2D
417  {
418  ParGridFunction d_x(&ind_fes), d_y(&ind_fes);
419 
421  DiscreteAdaptTC *tc = new DiscreteAdaptTC(target_t);
423  disc.ProjectCoefficient(ind_coeff);
424  if (adapt_eval == 0)
425  {
427  }
428  else
429  {
430 #ifdef MFEM_USE_GSLIB
432 #else
433  MFEM_ABORT("MFEM is not built with GSLIB.");
434 #endif
435  }
436  //Diffuse the interface
437  DiffuseField(disc,2);
438 
439  //Get partials with respect to x and y of the grid function
440  disc.GetDerivative(1,0,d_x);
441  disc.GetDerivative(1,1,d_y);
442 
443  //Compute the squared magnitude of the gradient
444  for (int i = 0; i < size.Size(); i++)
445  {
446  size(i) = std::pow(d_x(i),2)+std::pow(d_y(i),2);
447  }
448  const double max = size.Max();
449  double max_all;
450  MPI_Allreduce(&max, &max_all, 1, MPI_DOUBLE, MPI_MAX, MPI_COMM_WORLD);
451 
452  for (int i = 0; i < d_x.Size(); i++)
453  {
454  d_x(i) = std::abs(d_x(i));
455  d_y(i) = std::abs(d_y(i));
456  }
457  const double eps = 0.01;
458  const double aspr_ratio = 20.0;
459  const double size_ratio = 40.0;
460 
461  for (int i = 0; i < size.Size(); i++)
462  {
463  size(i) = (size(i)/max_all);
464  aspr(i) = (d_x(i)+eps)/(d_y(i)+eps);
465  aspr(i) = 0.1 + 0.9*(1-size(i))*(1-size(i));
466  if (aspr(i) > aspr_ratio) {aspr(i) = aspr_ratio;}
467  if (aspr(i) < 1.0/aspr_ratio) {aspr(i) = 1.0/aspr_ratio;}
468  }
469  Vector vals;
470  const int NE = pmesh->GetNE();
471  double volume = 0.0, volume_ind = 0.0;
472 
473  for (int i = 0; i < NE; i++)
474  {
476  const IntegrationRule &ir =
477  IntRules.Get(pmesh->GetElementBaseGeometry(i), Tr->OrderJ());
478  size.GetValues(i, ir, vals);
479  for (int j = 0; j < ir.GetNPoints(); j++)
480  {
481  const IntegrationPoint &ip = ir.IntPoint(j);
482  Tr->SetIntPoint(&ip);
483  volume += ip.weight * Tr->Weight();
484  volume_ind += vals(j) * ip.weight * Tr->Weight();
485  }
486  }
487  double volume_all, volume_ind_all;
488  MPI_Allreduce(&volume, &volume_all, 1, MPI_DOUBLE, MPI_SUM, MPI_COMM_WORLD);
489  MPI_Allreduce(&volume_ind, &volume_ind_all, 1, MPI_DOUBLE, MPI_SUM,
490  MPI_COMM_WORLD);
491  const int NE_ALL = pmesh->GetGlobalNE();
492 
493  const double avg_zone_size = volume_all / NE_ALL;
494 
495  const double small_avg_ratio =
496  (volume_ind_all + (volume_all - volume_ind_all) / size_ratio)
497  / volume_all;
498 
499  const double small_zone_size = small_avg_ratio * avg_zone_size;
500  const double big_zone_size = size_ratio * small_zone_size;
501 
502  for (int i = 0; i < size.Size(); i++)
503  {
504  const double val = size(i);
505  const double a = (big_zone_size - small_zone_size) / small_zone_size;
506  size(i) = big_zone_size / (1.0+a*val);
507  }
508 
509  DiffuseField(size, 2);
510  DiffuseField(aspr, 2);
511 
512  tc->SetParDiscreteTargetSize(size);
514  target_c = tc;
515  break;
516  }
517  case 7: // aspect-ratio 3D
518  {
520  DiscreteAdaptTC *tc = new DiscreteAdaptTC(target_t);
521  if (adapt_eval == 0)
522  {
524  }
525  else
526  {
527 #ifdef MFEM_USE_GSLIB
529 #else
530  MFEM_ABORT("MFEM is not built with GSLIB.");
531 #endif
532  }
534  aspr3d.ProjectCoefficient(fd_aspr3d);
536  target_c = tc;
537  break;
538  }
539  case 8: // shape/size + orientation 2D
540  {
542  DiscreteAdaptTC *tc = new DiscreteAdaptTC(target_t);
543  if (adapt_eval == 0)
544  {
546  }
547  else
548  {
549 #ifdef MFEM_USE_GSLIB
551 #else
552  MFEM_ABORT("MFEM is not built with GSLIB.");
553 #endif
554  }
555 
556  if (metric_id == 14)
557  {
558  ConstantCoefficient ind_coeff(0.1*0.1);
559  size.ProjectCoefficient(ind_coeff);
560  tc->SetParDiscreteTargetSize(size);
561  }
562 
563  if (metric_id == 85)
564  {
566  aspr.ProjectCoefficient(aspr_coeff);
567  DiffuseField(aspr,2);
569  }
570 
572  ori.ProjectCoefficient(ori_coeff);
574  target_c = tc;
575  break;
576  }
577  default:
578  if (myid == 0) { cout << "Unknown target_id: " << target_id << endl; }
579  return 3;
580  }
581 
582  if (target_c == NULL)
583  {
584  target_c = new TargetConstructor(target_t, MPI_COMM_WORLD);
585  }
586  target_c->SetNodes(x0);
587  TMOP_Integrator *he_nlf_integ= new TMOP_Integrator(metric, target_c);
588  if (fdscheme) { he_nlf_integ->EnableFiniteDifferences(x); }
589  he_nlf_integ->SetExactActionFlag(exactaction);
590 
591  // Setup the quadrature rules for the TMOP integrator.
592  IntegrationRules *irules = NULL;
593  switch (quad_type)
594  {
595  case 1: irules = &IntRulesLo; break;
596  case 2: irules = &IntRules; break;
597  case 3: irules = &IntRulesCU; break;
598  default:
599  if (myid == 0) { cout << "Unknown quad_type: " << quad_type << endl; }
600  return 3;
601  }
602  he_nlf_integ->SetIntegrationRules(*irules, quad_order);
603  if (myid == 0 && dim == 2)
604  {
605  cout << "Triangle quadrature points: "
606  << irules->Get(Geometry::TRIANGLE, quad_order).GetNPoints()
607  << "\nQuadrilateral quadrature points: "
608  << irules->Get(Geometry::SQUARE, quad_order).GetNPoints() << endl;
609  }
610  if (myid == 0 && dim == 3)
611  {
612  cout << "Tetrahedron quadrature points: "
613  << irules->Get(Geometry::TETRAHEDRON, quad_order).GetNPoints()
614  << "\nHexahedron quadrature points: "
615  << irules->Get(Geometry::CUBE, quad_order).GetNPoints()
616  << "\nPrism quadrature points: "
617  << irules->Get(Geometry::PRISM, quad_order).GetNPoints() << endl;
618  }
619 
620  if (normalization) { he_nlf_integ->ParEnableNormalization(x0); }
621 
622  // Limit the node movement.
623  // The limiting distances can be given by a general function of space.
624  ParGridFunction dist(pfespace);
625  dist = 1.0;
626  // The small_phys_size is relevant only with proper normalization.
627  if (normalization) { dist = small_phys_size; }
628  ConstantCoefficient lim_coeff(lim_const);
629  if (lim_const != 0.0) { he_nlf_integ->EnableLimiting(x0, dist, lim_coeff); }
630 
631  // Adaptive limiting.
632  ParGridFunction zeta_0(&ind_fes);
633  ConstantCoefficient coef_zeta(adapt_lim_const);
634  AdaptivityEvaluator *adapt_evaluator = NULL;
635  if (adapt_lim_const > 0.0)
636  {
638  zeta_0.ProjectCoefficient(alim_coeff);
639 
640  if (adapt_eval == 0) { adapt_evaluator = new AdvectorCG; }
641  else if (adapt_eval == 1)
642  {
643 #ifdef MFEM_USE_GSLIB
644  adapt_evaluator = new InterpolatorFP;
645 #else
646  MFEM_ABORT("MFEM is not built with GSLIB support!");
647 #endif
648  }
649  else { MFEM_ABORT("Bad interpolation option."); }
650 
651  he_nlf_integ->EnableAdaptiveLimiting(zeta_0, coef_zeta, *adapt_evaluator);
652  if (visualization)
653  {
654  socketstream vis1;
655  common::VisualizeField(vis1, "localhost", 19916, zeta_0, "Zeta 0",
656  300, 600, 300, 300);
657  }
658  }
659 
660  // 13. Setup the final NonlinearForm (which defines the integral of interest,
661  // its first and second derivatives). Here we can use a combination of
662  // metrics, i.e., optimize the sum of two integrals, where both are
663  // scaled by used-defined space-dependent weights. Note that there are
664  // no command-line options for the weights and the type of the second
665  // metric; one should update those in the code.
666  ParNonlinearForm a(pfespace);
667  ConstantCoefficient *coeff1 = NULL;
668  TMOP_QualityMetric *metric2 = NULL;
669  TargetConstructor *target_c2 = NULL;
671 
672  // Explicit combination of metrics.
673  if (combomet > 0)
674  {
675  // First metric.
676  coeff1 = new ConstantCoefficient(1.0);
677  he_nlf_integ->SetCoefficient(*coeff1);
678 
679  // Second metric.
680  metric2 = new TMOP_Metric_077;
681  TMOP_Integrator *he_nlf_integ2 = NULL;
682  if (combomet == 1)
683  {
684  target_c2 = new TargetConstructor(
686  target_c2->SetVolumeScale(0.01);
687  target_c2->SetNodes(x0);
688  he_nlf_integ2 = new TMOP_Integrator(metric2, target_c2);
689  he_nlf_integ2->SetCoefficient(coeff2);
690  }
691  else { he_nlf_integ2 = new TMOP_Integrator(metric2, target_c); }
692  he_nlf_integ2->SetIntegrationRules(*irules, quad_order);
693  if (fdscheme) { he_nlf_integ2->EnableFiniteDifferences(x); }
694  he_nlf_integ2->SetExactActionFlag(exactaction);
695 
697  combo->AddTMOPIntegrator(he_nlf_integ);
698  combo->AddTMOPIntegrator(he_nlf_integ2);
699  if (normalization) { combo->ParEnableNormalization(x0); }
700  if (lim_const != 0.0) { combo->EnableLimiting(x0, dist, lim_coeff); }
701 
702  a.AddDomainIntegrator(combo);
703  }
704  else { a.AddDomainIntegrator(he_nlf_integ); }
705 
706  const double init_energy = a.GetParGridFunctionEnergy(x);
707 
708  // Visualize the starting mesh and metric values.
709  // Note that for combinations of metrics, this only shows the first metric.
710  if (visualization)
711  {
712  char title[] = "Initial metric values";
713  vis_tmop_metric_p(mesh_poly_deg, *metric, *target_c, *pmesh, title, 0);
714  }
715 
716  // 14. Fix all boundary nodes, or fix only a given component depending on the
717  // boundary attributes of the given mesh. Attributes 1/2/3 correspond to
718  // fixed x/y/z components of the node. Attribute 4 corresponds to an
719  // entirely fixed node. Other boundary attributes do not affect the node
720  // movement boundary conditions.
721  if (move_bnd == false)
722  {
723  Array<int> ess_bdr(pmesh->bdr_attributes.Max());
724  ess_bdr = 1;
725  a.SetEssentialBC(ess_bdr);
726  }
727  else
728  {
729  int n = 0;
730  for (int i = 0; i < pmesh->GetNBE(); i++)
731  {
732  const int nd = pfespace->GetBE(i)->GetDof();
733  const int attr = pmesh->GetBdrElement(i)->GetAttribute();
734  MFEM_VERIFY(!(dim == 2 && attr == 3),
735  "Boundary attribute 3 must be used only for 3D meshes. "
736  "Adjust the attributes (1/2/3/4 for fixed x/y/z/all "
737  "components, rest for free nodes), or use -fix-bnd.");
738  if (attr == 1 || attr == 2 || attr == 3) { n += nd; }
739  if (attr == 4) { n += nd * dim; }
740  }
741  Array<int> ess_vdofs(n), vdofs;
742  n = 0;
743  for (int i = 0; i < pmesh->GetNBE(); i++)
744  {
745  const int nd = pfespace->GetBE(i)->GetDof();
746  const int attr = pmesh->GetBdrElement(i)->GetAttribute();
747  pfespace->GetBdrElementVDofs(i, vdofs);
748  if (attr == 1) // Fix x components.
749  {
750  for (int j = 0; j < nd; j++)
751  { ess_vdofs[n++] = vdofs[j]; }
752  }
753  else if (attr == 2) // Fix y components.
754  {
755  for (int j = 0; j < nd; j++)
756  { ess_vdofs[n++] = vdofs[j+nd]; }
757  }
758  else if (attr == 3) // Fix z components.
759  {
760  for (int j = 0; j < nd; j++)
761  { ess_vdofs[n++] = vdofs[j+2*nd]; }
762  }
763  else if (attr == 4) // Fix all components.
764  {
765  for (int j = 0; j < vdofs.Size(); j++)
766  { ess_vdofs[n++] = vdofs[j]; }
767  }
768  }
769  a.SetEssentialVDofs(ess_vdofs);
770  }
771 
772  // 15. As we use the Newton method to solve the resulting nonlinear system,
773  // here we setup the linear solver for the system's Jacobian.
774  Solver *S = NULL;
775  const double linsol_rtol = 1e-12;
776  if (lin_solver == 0)
777  {
778  S = new DSmoother(1, 1.0, max_lin_iter);
779  }
780  else if (lin_solver == 1)
781  {
782  CGSolver *cg = new CGSolver(MPI_COMM_WORLD);
783  cg->SetMaxIter(max_lin_iter);
784  cg->SetRelTol(linsol_rtol);
785  cg->SetAbsTol(0.0);
786  cg->SetPrintLevel(verbosity_level >= 2 ? 3 : -1);
787  S = cg;
788  }
789  else
790  {
791  MINRESSolver *minres = new MINRESSolver(MPI_COMM_WORLD);
792  minres->SetMaxIter(max_lin_iter);
793  minres->SetRelTol(linsol_rtol);
794  minres->SetAbsTol(0.0);
795  minres->SetPrintLevel(verbosity_level >= 2 ? 3 : -1);
796  S = minres;
797  }
798 
799  // Compute the minimum det(J) of the starting mesh.
800  tauval = infinity();
801  const int NE = pmesh->GetNE();
802  for (int i = 0; i < NE; i++)
803  {
804  const IntegrationRule &ir =
805  irules->Get(pfespace->GetFE(i)->GetGeomType(), quad_order);
807  for (int j = 0; j < ir.GetNPoints(); j++)
808  {
809  transf->SetIntPoint(&ir.IntPoint(j));
810  tauval = min(tauval, transf->Jacobian().Det());
811  }
812  }
813  double minJ0;
814  MPI_Allreduce(&tauval, &minJ0, 1, MPI_DOUBLE, MPI_MIN, MPI_COMM_WORLD);
815  tauval = minJ0;
816  if (myid == 0)
817  { cout << "Minimum det(J) of the original mesh is " << tauval << endl; }
818  double h0min = h0.Min(), h0min_all;
819  MPI_Allreduce(&h0min, &h0min_all, 1, MPI_DOUBLE, MPI_MIN, MPI_COMM_WORLD);
820  tauval -= 0.01 * h0min_all; // Slightly below minJ0 to avoid div by 0.
821 
822  // Perform the nonlinear optimization.
823  const IntegrationRule &ir =
824  irules->Get(pfespace->GetFE(0)->GetGeomType(), quad_order);
825  TMOPNewtonSolver solver(pfespace->GetComm(), ir, solver_type);
826  // Provide all integration rules in case of a mixed mesh.
827  solver.SetIntegrationRules(*irules, quad_order);
828  if (solver_type == 0)
829  {
830  // Specify linear solver when we use a Newton-based solver.
831  solver.SetPreconditioner(*S);
832  }
833  solver.SetMaxIter(solver_iter);
834  solver.SetRelTol(solver_rtol);
835  solver.SetAbsTol(0.0);
836  solver.SetPrintLevel(verbosity_level >= 1 ? 1 : -1);
837  solver.SetOperator(a);
838  solver.Mult(b, x.GetTrueVector());
839  x.SetFromTrueVector();
840  if (myid == 0 && solver.GetConverged() == false)
841  {
842  cout << "Nonlinear solver: rtol = " << solver_rtol << " not achieved.\n";
843  }
844 
845  // 16. Save the optimized mesh to a file. This output can be viewed later
846  // using GLVis: "glvis -m optimized -np num_mpi_tasks".
847  {
848  ostringstream mesh_name;
849  mesh_name << "optimized.mesh";
850  ofstream mesh_ofs(mesh_name.str().c_str());
851  mesh_ofs.precision(8);
852  pmesh->PrintAsOne(mesh_ofs);
853  }
854 
855  // 17. Compute the amount of energy decrease.
856  const double fin_energy = a.GetParGridFunctionEnergy(x);
857  double metric_part = fin_energy;
858  if (lim_const > 0.0 || adapt_lim_const > 0.0)
859  {
860  lim_coeff.constant = 0.0;
861  coef_zeta.constant = 0.0;
862  metric_part = a.GetParGridFunctionEnergy(x);
863  lim_coeff.constant = lim_const;
864  coef_zeta.constant = adapt_lim_const;
865  }
866  if (myid == 0)
867  {
868  cout << "Initial strain energy: " << init_energy
869  << " = metrics: " << init_energy
870  << " + limiting term: " << 0.0 << endl;
871  cout << " Final strain energy: " << fin_energy
872  << " = metrics: " << metric_part
873  << " + limiting term: " << fin_energy - metric_part << endl;
874  cout << "The strain energy decreased by: " << setprecision(12)
875  << (init_energy - fin_energy) * 100.0 / init_energy << " %." << endl;
876  }
877 
878  // 18. Visualize the final mesh and metric values.
879  if (visualization)
880  {
881  char title[] = "Final metric values";
882  vis_tmop_metric_p(mesh_poly_deg, *metric, *target_c, *pmesh, title, 600);
883  }
884 
885  if (adapt_lim_const > 0.0 && visualization)
886  {
887  socketstream vis0;
888  common::VisualizeField(vis0, "localhost", 19916, zeta_0, "Xi 0",
889  600, 600, 300, 300);
890  }
891 
892  // 19. Visualize the mesh displacement.
893  if (visualization)
894  {
895  x0 -= x;
896  socketstream sock;
897  if (myid == 0)
898  {
899  sock.open("localhost", 19916);
900  sock << "solution\n";
901  }
902  pmesh->PrintAsOne(sock);
903  x0.SaveAsOne(sock);
904  if (myid == 0)
905  {
906  sock << "window_title 'Displacements'\n"
907  << "window_geometry "
908  << 1200 << " " << 0 << " " << 600 << " " << 600 << "\n"
909  << "keys jRmclA" << endl;
910  }
911  }
912 
913  // 20. Free the used memory.
914  delete S;
915  delete target_c2;
916  delete metric2;
917  delete coeff1;
918  delete adapt_evaluator;
919  delete target_c;
920  delete adapt_coeff;
921  delete metric;
922  delete pfespace;
923  delete fec;
924  delete pmesh;
925 
926  MPI_Finalize();
927  return 0;
928 }
int GetNPoints() const
Returns the number of the points in the integration rule.
Definition: intrules.hpp:245
void vis_tmop_metric_p(int order, TMOP_QualityMetric &qm, const TargetConstructor &tc, ParMesh &pmesh, char *title, int position)
Definition: tmop_tools.cpp:623
void discrete_aspr_3d(const Vector &x, Vector &v)
virtual void SetParDiscreteTargetSize(const ParGridFunction &tspec_)
Definition: tmop.cpp:1049
int Size() const
Return the logical size of the array.
Definition: array.hpp:124
Shifted barrier form of metric 56 (area, ideal barrier metric), 2D.
Definition: tmop.hpp:354
Shifted barrier form of 3D metric 16 (volume, ideal barrier metric), 3D.
Definition: tmop.hpp:472
Conjugate gradient method.
Definition: solvers.hpp:258
Definition: ex25.cpp:144
int GetNDofs() const
Returns number of degrees of freedom.
Definition: fespace.hpp:397
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:2747
int DofToVDof(int dof, int vd, int ndofs=-1) const
Definition: fespace.cpp:144
double discrete_size_2d(const Vector &x)
Shape &amp; volume, ideal barrier metric, 3D.
Definition: tmop.hpp:456
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:915
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:147
void SetFromTrueVector()
Shortcut for calling SetFromTrueDofs() with GetTrueVector() as argument.
Definition: gridfunc.hpp:142
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:763
int GetNBE() const
Returns number of boundary elements.
Definition: mesh.hpp:740
Shape, ideal barrier metric, 3D.
Definition: tmop.hpp:406
void SetAdaptivityEvaluator(AdaptivityEvaluator *ae)
Definition: tmop.hpp:849
double Det() const
Definition: densemat.cpp:451
void DiffuseField(GridFunction &field, int smooth_steps)
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
void SaveAsOne(std::ostream &out=mfem::out)
Merge the local grid functions.
Definition: pgridfunc.cpp:855
Container class for integration rules.
Definition: intrules.hpp:309
Shape, ideal barrier metric, 3D.
Definition: tmop.hpp:390
void SetVolumeScale(double vol_scale)
Used by target type IDEAL_SHAPE_EQUAL_SIZE. The default volume scale is 1.
Definition: tmop.hpp:668
Parallel non-linear operator on the true dofs.
Shape, ideal barrier metric, 2D.
Definition: tmop.hpp:166
void EnableAdaptiveLimiting(const GridFunction &z0, Coefficient &coeff, AdaptivityEvaluator &ae)
Restriction of the node positions to certain regions.
Definition: tmop.cpp:1921
Volume metric, 3D.
Definition: tmop.hpp:422
Area, ideal barrier metric, 2D.
Definition: tmop.hpp:303
void SetIntegrationRules(IntegrationRules &irules, int order)
Prescribe a set of integration rules; relevant for mixed meshes.
Definition: tmop.hpp:1032
int GetNE() const
Returns number of elements.
Definition: mesh.hpp:737
Abstract parallel finite element space.
Definition: pfespace.hpp:28
virtual void ProjectCoefficient(Coefficient &coeff)
Definition: pgridfunc.cpp:474
void SetExactActionFlag(bool flag_)
Flag to control if exact action of Integration is effected.
Definition: tmop.hpp:1123
MINRES method.
Definition: solvers.hpp:368
void Randomize(int seed=0)
Set random values in the vector.
Definition: vector.cpp:725
int main(int argc, char *argv[])
Definition: ex1.cpp:66
virtual int OrderJ() const =0
Return the order of the elements of the Jacobian of the transformation.
Shape &amp; area, ideal barrier metric, 2D.
Definition: tmop.hpp:182
Geometry::Type GetElementBaseGeometry(int i) const
Definition: mesh.hpp:834
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:2775
Shifted barrier form of metric 2 (shape, ideal barrier metric), 2D.
Definition: tmop.hpp:214
Shape, ideal barrier metric, 2D.
Definition: tmop.hpp:285
Geometry::Type GetGeomType() const
Returns the Geometry::Type of the reference element.
Definition: fe.hpp:312
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:248
void SetTrueVector()
Shortcut for calling GetTrueDofs() with GetTrueVector() as argument.
Definition: gridfunc.hpp:136
void SetPrintLevel(int print_lvl)
Definition: solvers.cpp:70
MPI_Comm GetComm() const
Definition: pfespace.hpp:239
void EnableLimiting(const GridFunction &n0, const GridFunction &dist, Coefficient &w0, TMOP_LimiterFunction *lfunc=NULL)
Limiting of the mesh displacements (general version).
Definition: tmop.cpp:1896
int GetNBE() const
Returns number of boundary elements in the mesh.
Definition: fespace.hpp:436
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:1044
Area metric, 2D.
Definition: tmop.hpp:249
Volume, ideal barrier metric, 3D.
Definition: tmop.hpp:438
double b
Definition: lissajous.cpp:42
void PrintAsOne(std::ostream &out=mfem::out)
Definition: pmesh.cpp:4305
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:8382
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:665
void SetMaxIter(int max_it)
Definition: solvers.hpp:98
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:483
virtual void SetParDiscreteTargetAspectRatio(const ParGridFunction &tspec_)
Definition: tmop.cpp:1065
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
void ParEnableNormalization(const ParGridFunction &x)
Definition: tmop.cpp:2864
void SetNodalFESpace(FiniteElementSpace *nfes)
Definition: mesh.cpp:4126
int Dimension() const
Definition: mesh.hpp:788
void PrintUsage(std::ostream &out) const
Print the usage message.
Definition: optparser.cpp:434
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:124
A general vector function coefficient.
double GetElementSize(int i, int type=0)
Get the size of the i-th element relative to the perfect reference element.
Definition: mesh.cpp:74
void SetAbsTol(double atol)
Definition: solvers.hpp:97
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
Definition: mesh.hpp:201
void SetRelTol(double rtol)
Definition: solvers.hpp:96
Shape, ideal barrier metric, 2D.
Definition: tmop.hpp:233
Untangling metric, 2D.
Definition: tmop.hpp:335
int GetDof() const
Returns the number of degrees of freedom in the finite element.
Definition: fe.hpp:315
double weight_fun(const Vector &x)
Area, ideal barrier metric, 2D.
Definition: tmop.hpp:266
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:1141
double material_indicator_2d(const Vector &x)
Shape &amp; area metric, 2D.
Definition: tmop.hpp:198
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:460
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:336
int dim
Definition: ex24.cpp:53
IntegrationRules IntRulesCU(0, Quadrature1D::ClosedUniform)
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:1798
void PrintOptions(std::ostream &out) const
Print the options.
Definition: optparser.cpp:304
double infinity()
Define a shortcut for std::numeric_limits&lt;double&gt;::infinity()
Definition: vector.hpp:45
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:51
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:6603
Shape+Size+Orientation metric, 2D.
Definition: tmop.hpp:151
Arbitrary order H1-conforming (continuous) finite elements.
Definition: fe_coll.hpp:159
TargetType
Target-matrix construction algorithms supported by this class.
Definition: tmop.hpp:607
void SetEssentialVDofs(const Array< int > &ess_vdofs_list)
Specify essential boundary conditions.
virtual void SetParDiscreteTargetOrientation(const ParGridFunction &tspec_)
Definition: tmop.cpp:1074
Base class for solvers.
Definition: operator.hpp:634
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:2047
void GetBdrElementVDofs(int i, Array< int > &vdofs) const
Returns indexes of degrees of freedom for i&#39;th boundary element.
Definition: fespace.cpp:179
void SetNodalGridFunction(GridFunction *nodes, bool make_owner=false)
Definition: mesh.cpp:4154
Base class representing target-matrix construction algorithms for mesh optimization via the target-ma...
Definition: tmop.hpp:603
virtual void SetAnalyticTargetSpec(Coefficient *sspec, VectorCoefficient *vspec, TMOPMatrixCoefficient *mspec)
Definition: tmop.cpp:941
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:378
Shape &amp; orientation metric, 2D.
Definition: tmop.hpp:320
void ParEnableNormalization(const ParGridFunction &x)
Definition: tmop.cpp:2627
double GetElementVolume(int i)
Definition: mesh.cpp:101
const Element * GetBdrElement(int i) const
Definition: mesh.hpp:823
Shape, ideal barrier metric, 3D.
Definition: tmop.hpp:374
Metric without a type, 2D.
Definition: tmop.hpp:75
bool Good() const
Return true if the command line options were parsed successfully.
Definition: optparser.hpp:145
A TMOP integrator class based on any given TMOP_QualityMetric and TargetConstructor.
Definition: tmop.hpp:884