MFEM  v4.1.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
mesh-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
14 // --------------------------------------------------
15 //
16 // This miniapp performs mesh optimization using the Target-Matrix Optimization
17 // Paradigm (TMOP) by P.Knupp et al., and a global variational minimization
18 // approach. It minimizes the quantity sum_T int_T mu(J(x)), where T are the
19 // target (ideal) elements, J is the Jacobian of the transformation from the
20 // target to the physical element, and mu is the mesh quality metric. This
21 // metric can measure shape, size or alignment of the region around each
22 // quadrature point. The combination of targets & quality metrics is used to
23 // optimize the physical node positions, i.e., they must be as close as possible
24 // to the shape / size / alignment of their targets. This code also demonstrates
25 // a possible use of nonlinear operators (the class TMOP_QualityMetric, defining
26 // mu(J), and the class TMOP_Integrator, defining int mu(J)), as well as their
27 // coupling to Newton methods for solving minimization problems. Note that the
28 // utilized Newton methods are oriented towards avoiding invalid meshes with
29 // negative Jacobian determinants. Each Newton step requires the inversion of a
30 // Jacobian matrix, which is done through an inner linear solver.
31 //
32 // Compile with: make mesh-optimizer
33 //
34 // Sample runs:
35 // Adapted analytic Hessian:
36 // mesh-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 discrete size:
38 // mesh-optimizer -m square01.mesh -o 2 -rs 2 -mid 7 -tid 5 -ni 200 -ls 2 -li 100 -bnd -qt 1 -qo 8
39 //
40 // Blade shape:
41 // mesh-optimizer -m blade.mesh -o 4 -rs 0 -mid 2 -tid 1 -ni 200 -ls 2 -li 100 -bnd -qt 1 -qo 8
42 // Blade limited shape:
43 // mesh-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
44 // ICF shape and equal size:
45 // mesh-optimizer -o 3 -rs 0 -mid 9 -tid 2 -ni 200 -ls 2 -li 100 -bnd -qt 1 -qo 8
46 // ICF shape and initial size:
47 // mesh-optimizer -o 3 -rs 0 -mid 9 -tid 3 -ni 100 -ls 2 -li 100 -bnd -qt 1 -qo 8
48 // ICF shape:
49 // mesh-optimizer -o 3 -rs 0 -mid 1 -tid 1 -ni 100 -ls 2 -li 100 -bnd -qt 1 -qo 8
50 // ICF limited shape:
51 // mesh-optimizer -o 3 -rs 0 -mid 1 -tid 1 -ni 100 -ls 2 -li 100 -bnd -qt 1 -qo 8 -lc 10
52 // ICF combo shape + size (rings, slow convergence):
53 // mesh-optimizer -o 3 -rs 0 -mid 1 -tid 1 -ni 1000 -ls 2 -li 100 -bnd -qt 1 -qo 8 -cmb
54 // 3D pinched sphere shape (the mesh is in the mfem/data GitHub repository):
55 // * mesh-optimizer -m ../../../mfem_data/ball-pert.mesh -o 4 -rs 0 -mid 303 -tid 1 -ni 20 -ls 2 -li 500 -fix-bnd
56 // 2D non-conforming shape and equal size:
57 // mesh-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
58 
59 
60 #include "mfem.hpp"
61 #include <fstream>
62 #include <iostream>
63 
64 using namespace mfem;
65 using namespace std;
66 
67 double weight_fun(const Vector &x);
68 
69 double ind_values(const Vector &x)
70 {
71  const int opt = 6;
72  const double small = 0.001, big = 0.01;
73 
74  // Sine wave.
75  if (opt==1)
76  {
77  const double X = x(0), Y = x(1);
78  const double ind = std::tanh((10*(Y-0.5) + std::sin(4.0*M_PI*X)) + 1) -
79  std::tanh((10*(Y-0.5) + std::sin(4.0*M_PI*X)) - 1);
80 
81  return ind * small + (1.0 - ind) * big;
82  }
83 
84  if (opt==2)
85  {
86  // Circle in the middle.
87  double val = 0.;
88  const double xc = x(0) - 0.5, yc = x(1) - 0.5;
89  const double r = sqrt(xc*xc + yc*yc);
90  double r1 = 0.15; double r2 = 0.35; double sf=30.0;
91  val = 0.5*(std::tanh(sf*(r-r1)) - std::tanh(sf*(r-r2)));
92  if (val > 1.) {val = 1;}
93 
94  return val * small + (1.0 - val) * big;
95  }
96 
97  if (opt == 3)
98  {
99  // cross
100  const double X = x(0), Y = x(1);
101  const double r1 = 0.45, r2 = 0.55;
102  const double sf = 40.0;
103 
104  double val = 0.5 * ( std::tanh(sf*(X-r1)) - std::tanh(sf*(X-r2)) +
105  std::tanh(sf*(Y-r1)) - std::tanh(sf*(Y-r2)) );
106  if (val > 1.) { val = 1.0; }
107 
108  return val * small + (1.0 - val) * big;
109  }
110 
111  if (opt==4)
112  {
113  // Multiple circles
114  double r1,r2,val,rval;
115  double sf = 10;
116  val = 0.;
117  // circle 1
118  r1= 0.25; r2 = 0.25; rval = 0.1;
119  double xc = x(0) - r1, yc = x(1) - r2;
120  double r = sqrt(xc*xc+yc*yc);
121  val = 0.5*(1+std::tanh(sf*(r+rval))) - 0.5*(1+std::tanh(sf*
122  (r-rval)));// std::exp(val1);
123  // circle 2
124  r1= 0.75; r2 = 0.75;
125  xc = x(0) - r1, yc = x(1) - r2;
126  r = sqrt(xc*xc+yc*yc);
127  val += (0.5*(1+std::tanh(sf*(r+rval))) - 0.5*(1+std::tanh(sf*
128  (r-rval))));// std::exp(val1);
129  // circle 3
130  r1= 0.75; r2 = 0.25;
131  xc = x(0) - r1, yc = x(1) - r2;
132  r = sqrt(xc*xc+yc*yc);
133  val += 0.5*(1+std::tanh(sf*(r+rval))) - 0.5*(1+std::tanh(sf*
134  (r-rval)));// std::exp(val1);
135  // circle 4
136  r1= 0.25; r2 = 0.75;
137  xc = x(0) - r1, yc = x(1) - r2;
138  r = sqrt(xc*xc+yc*yc);
139  val += 0.5*(1+std::tanh(sf*(r+rval))) - 0.5*(1+std::tanh(sf*(r-rval)));
140  if (val > 1.0) {val = 1.;}
141  if (val < 0.0) {val = 0.;}
142 
143  return val * small + (1.0 - val) * big;
144  }
145 
146  if (opt==5)
147  {
148  // cross
149  double val = 0.;
150  double X = x(0)-0.5, Y = x(1)-0.5;
151  double rval = std::sqrt(X*X + Y*Y);
152  double thval = 60.*M_PI/180.;
153  double Xmod,Ymod;
154  Xmod = X*std::cos(thval) + Y*std::sin(thval);
155  Ymod= -X*std::sin(thval) + Y*std::cos(thval);
156  X = Xmod+0.5; Y = Ymod+0.5;
157  double r1 = 0.45; double r2 = 0.55; double sf=30.0;
158  val = ( 0.5*(1+std::tanh(sf*(X-r1))) - 0.5*(1+std::tanh(sf*(X-r2)))
159  + 0.5*(1+std::tanh(sf*(Y-r1))) - 0.5*(1+std::tanh(sf*(Y-r2))) );
160  if (rval > 0.4) {val = 0.;}
161  if (val > 1.0) {val = 1.;}
162  if (val < 0.0) {val = 0.;}
163 
164  return val * small + (1.0 - val) * big;
165  }
166 
167  if (opt==6)
168  {
169  double val = 0.;
170  const double xc = x(0) - 0.0, yc = x(1) - 0.5;
171  const double r = sqrt(xc*xc + yc*yc);
172  double r1 = 0.45; double r2 = 0.55; double sf=30.0;
173  val = 0.5*(1+std::tanh(sf*(r-r1))) - 0.5*(1+std::tanh(sf*(r-r2)));
174  if (val > 1.) {val = 1;}
175  if (val < 0.) {val = 0;}
176 
177  return val * small + (1.0 - val) * big;
178  }
179 
180  return 0.0;
181 }
182 
183 class HessianCoefficient : public MatrixCoefficient
184 {
185 private:
186  int type;
187 
188 public:
189  HessianCoefficient(int dim, int type_)
190  : MatrixCoefficient(dim), type(type_) { }
191 
192  virtual void Eval(DenseMatrix &K, ElementTransformation &T,
193  const IntegrationPoint &ip)
194  {
195  Vector pos(3);
196  T.Transform(ip, pos);
197 
198  if (type == 0)
199  {
200  K(0, 0) = 1.0 + 3.0 * std::sin(M_PI*pos(0));
201  K(0, 1) = 0.0;
202  K(1, 0) = 0.0;
203  K(1, 1) = 1.0;
204  }
205  else
206  {
207  const double xc = pos(0) - 0.5, yc = pos(1) - 0.5;
208  const double r = sqrt(xc*xc + yc*yc);
209  double r1 = 0.15; double r2 = 0.35; double sf=30.0;
210  const double eps = 0.5;
211 
212  const double tan1 = std::tanh(sf*(r-r1)),
213  tan2 = std::tanh(sf*(r-r2));
214 
215  K(0, 0) = eps + 1.0 * (tan1 - tan2);
216  K(0, 1) = 0.0;
217  K(1, 0) = 0.0;
218  K(1, 1) = 1.0;
219  }
220  }
221 };
222 
223 // Additional IntegrationRules that can be used with the --quad-type option.
226 
227 
228 int main (int argc, char *argv[])
229 {
230  // 0. Set the method's default parameters.
231  const char *mesh_file = "icf.mesh";
232  int mesh_poly_deg = 1;
233  int rs_levels = 0;
234  double jitter = 0.0;
235  int metric_id = 1;
236  int target_id = 1;
237  double lim_const = 0.0;
238  int quad_type = 1;
239  int quad_order = 8;
240  int newton_iter = 10;
241  double newton_rtol = 1e-12;
242  int lin_solver = 2;
243  int max_lin_iter = 100;
244  bool move_bnd = true;
245  bool combomet = 0;
246  bool normalization = false;
247  bool visualization = true;
248  int verbosity_level = 0;
249 
250  // 1. Parse command-line options.
251  OptionsParser args(argc, argv);
252  args.AddOption(&mesh_file, "-m", "--mesh",
253  "Mesh file to use.");
254  args.AddOption(&mesh_poly_deg, "-o", "--order",
255  "Polynomial degree of mesh finite element space.");
256  args.AddOption(&rs_levels, "-rs", "--refine-serial",
257  "Number of times to refine the mesh uniformly in serial.");
258  args.AddOption(&jitter, "-ji", "--jitter",
259  "Random perturbation scaling factor.");
260  args.AddOption(&metric_id, "-mid", "--metric-id",
261  "Mesh optimization metric:\n\t"
262  "1 : |T|^2 -- 2D shape\n\t"
263  "2 : 0.5|T|^2/tau-1 -- 2D shape (condition number)\n\t"
264  "7 : |T-T^-t|^2 -- 2D shape+size\n\t"
265  "9 : tau*|T-T^-t|^2 -- 2D shape+size\n\t"
266  "22 : 0.5(|T|^2-2*tau)/(tau-tau_0) -- 2D untangling\n\t"
267  "50 : 0.5|T^tT|^2/tau^2-1 -- 2D shape\n\t"
268  "55 : (tau-1)^2 -- 2D size\n\t"
269  "56 : 0.5(sqrt(tau)-1/sqrt(tau))^2 -- 2D size\n\t"
270  "58 : |T^tT|^2/(tau^2)-2*|T|^2/tau+2 -- 2D shape\n\t"
271  "77 : 0.5(tau-1/tau)^2 -- 2D size\n\t"
272  "211: (tau-1)^2-tau+sqrt(tau^2) -- 2D untangling\n\t"
273  "252: 0.5(tau-1)^2/(tau-tau_0) -- 2D untangling\n\t"
274  "301: (|T||T^-1|)/3-1 -- 3D shape\n\t"
275  "302: (|T|^2|T^-1|^2)/9-1 -- 3D shape\n\t"
276  "303: (|T|^2)/3*tau^(2/3)-1 -- 3D shape\n\t"
277  "315: (tau-1)^2 -- 3D size\n\t"
278  "316: 0.5(sqrt(tau)-1/sqrt(tau))^2 -- 3D size\n\t"
279  "321: |T-T^-t|^2 -- 3D shape+size\n\t"
280  "352: 0.5(tau-1)^2/(tau-tau_0) -- 3D untangling");
281  args.AddOption(&target_id, "-tid", "--target-id",
282  "Target (ideal element) type:\n\t"
283  "1: Ideal shape, unit size\n\t"
284  "2: Ideal shape, equal size\n\t"
285  "3: Ideal shape, initial size\n\t"
286  "4: Given full analytic Jacobian (in physical space)\n\t"
287  "5: Ideal shape, given size (in physical space)");
288  args.AddOption(&lim_const, "-lc", "--limit-const", "Limiting constant.");
289  args.AddOption(&quad_type, "-qt", "--quad-type",
290  "Quadrature rule type:\n\t"
291  "1: Gauss-Lobatto\n\t"
292  "2: Gauss-Legendre\n\t"
293  "3: Closed uniform points");
294  args.AddOption(&quad_order, "-qo", "--quad_order",
295  "Order of the quadrature rule.");
296  args.AddOption(&newton_iter, "-ni", "--newton-iters",
297  "Maximum number of Newton iterations.");
298  args.AddOption(&newton_rtol, "-rtol", "--newton-rel-tolerance",
299  "Relative tolerance for the Newton solver.");
300  args.AddOption(&lin_solver, "-ls", "--lin-solver",
301  "Linear solver: 0 - l1-Jacobi, 1 - CG, 2 - MINRES.");
302  args.AddOption(&max_lin_iter, "-li", "--lin-iter",
303  "Maximum number of iterations in the linear solve.");
304  args.AddOption(&move_bnd, "-bnd", "--move-boundary", "-fix-bnd",
305  "--fix-boundary",
306  "Enable motion along horizontal and vertical boundaries.");
307  args.AddOption(&combomet, "-cmb", "--combo-met", "-no-cmb", "--no-combo-met",
308  "Combination of metrics.");
309  args.AddOption(&normalization, "-nor", "--normalization", "-no-nor",
310  "--no-normalization",
311  "Make all terms in the optimization functional unitless.");
312  args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
313  "--no-visualization",
314  "Enable or disable GLVis visualization.");
315  args.AddOption(&verbosity_level, "-vl", "--verbosity-level",
316  "Set the verbosity level - 0, 1, or 2.");
317  args.Parse();
318  if (!args.Good())
319  {
320  args.PrintUsage(cout);
321  return 1;
322  }
323  args.PrintOptions(cout);
324 
325  // 2. Initialize and refine the starting mesh.
326  Mesh *mesh = new Mesh(mesh_file, 1, 1, false);
327  for (int lev = 0; lev < rs_levels; lev++) { mesh->UniformRefinement(); }
328  const int dim = mesh->Dimension();
329  cout << "Mesh curvature: ";
330  if (mesh->GetNodes()) { cout << mesh->GetNodes()->OwnFEC()->Name(); }
331  else { cout << "(NONE)"; }
332  cout << endl;
333 
334  // 3. Define a finite element space on the mesh. Here we use vector finite
335  // elements which are tensor products of quadratic finite elements. The
336  // number of components in the vector finite element space is specified by
337  // the last parameter of the FiniteElementSpace constructor.
339  if (mesh_poly_deg <= 0)
340  {
341  fec = new QuadraticPosFECollection;
342  mesh_poly_deg = 2;
343  }
344  else { fec = new H1_FECollection(mesh_poly_deg, dim); }
345  FiniteElementSpace *fespace = new FiniteElementSpace(mesh, fec, dim);
346 
347  // 4. Make the mesh curved based on the above finite element space. This
348  // means that we define the mesh elements through a fespace-based
349  // transformation of the reference element.
350  mesh->SetNodalFESpace(fespace);
351 
352  // 5. Set up an empty right-hand side vector b, which is equivalent to b=0.
353  Vector b(0);
354 
355  // 6. Get the mesh nodes (vertices and other degrees of freedom in the finite
356  // element space) as a finite element grid function in fespace. Note that
357  // changing x automatically changes the shapes of the mesh elements.
358  GridFunction x(fespace);
359  mesh->SetNodalGridFunction(&x);
360 
361  // 7. Define a vector representing the minimal local mesh size in the mesh
362  // nodes. We index the nodes using the scalar version of the degrees of
363  // freedom in pfespace. Note: this is partition-dependent.
364  //
365  // In addition, compute average mesh size and total volume.
366  Vector h0(fespace->GetNDofs());
367  h0 = infinity();
368  double volume = 0.0;
369  Array<int> dofs;
370  for (int i = 0; i < mesh->GetNE(); i++)
371  {
372  // Get the local scalar element degrees of freedom in dofs.
373  fespace->GetElementDofs(i, dofs);
374  // Adjust the value of h0 in dofs based on the local mesh size.
375  const double hi = mesh->GetElementSize(i);
376  for (int j = 0; j < dofs.Size(); j++)
377  {
378  h0(dofs[j]) = min(h0(dofs[j]), hi);
379  }
380  volume += mesh->GetElementVolume(i);
381  }
382  const double small_phys_size = pow(volume, 1.0 / dim) / 100.0;
383 
384  // 8. Add a random perturbation to the nodes in the interior of the domain.
385  // We define a random grid function of fespace and make sure that it is
386  // zero on the boundary and its values are locally of the order of h0.
387  // The latter is based on the DofToVDof() method which maps the scalar to
388  // the vector degrees of freedom in fespace.
389  GridFunction rdm(fespace);
390  rdm.Randomize();
391  rdm -= 0.25; // Shift to random values in [-0.5,0.5].
392  rdm *= jitter;
393  // Scale the random values to be of order of the local mesh size.
394  for (int i = 0; i < fespace->GetNDofs(); i++)
395  {
396  for (int d = 0; d < dim; d++)
397  {
398  rdm(fespace->DofToVDof(i,d)) *= h0(i);
399  }
400  }
401  Array<int> vdofs;
402  for (int i = 0; i < fespace->GetNBE(); i++)
403  {
404  // Get the vector degrees of freedom in the boundary element.
405  fespace->GetBdrElementVDofs(i, vdofs);
406  // Set the boundary values to zero.
407  for (int j = 0; j < vdofs.Size(); j++) { rdm(vdofs[j]) = 0.0; }
408  }
409  x -= rdm;
410  x.SetTrueVector();
411  x.SetFromTrueVector();
412 
413  // 9. Save the starting (prior to the optimization) mesh to a file. This
414  // output can be viewed later using GLVis: "glvis -m perturbed.mesh".
415  {
416  ofstream mesh_ofs("perturbed.mesh");
417  mesh->Print(mesh_ofs);
418  }
419 
420  // 10. Store the starting (prior to the optimization) positions.
421  GridFunction x0(fespace);
422  x0 = x;
423 
424  // 11. Form the integrator that uses the chosen metric and target.
425  double tauval = -0.1;
426  TMOP_QualityMetric *metric = NULL;
427  switch (metric_id)
428  {
429  case 1: metric = new TMOP_Metric_001; break;
430  case 2: metric = new TMOP_Metric_002; break;
431  case 7: metric = new TMOP_Metric_007; break;
432  case 9: metric = new TMOP_Metric_009; break;
433  case 22: metric = new TMOP_Metric_022(tauval); break;
434  case 50: metric = new TMOP_Metric_050; break;
435  case 55: metric = new TMOP_Metric_055; break;
436  case 56: metric = new TMOP_Metric_056; break;
437  case 58: metric = new TMOP_Metric_058; break;
438  case 77: metric = new TMOP_Metric_077; break;
439  case 211: metric = new TMOP_Metric_211; break;
440  case 252: metric = new TMOP_Metric_252(tauval); break;
441  case 301: metric = new TMOP_Metric_301; break;
442  case 302: metric = new TMOP_Metric_302; break;
443  case 303: metric = new TMOP_Metric_303; break;
444  case 315: metric = new TMOP_Metric_315; break;
445  case 316: metric = new TMOP_Metric_316; break;
446  case 321: metric = new TMOP_Metric_321; break;
447  case 352: metric = new TMOP_Metric_352(tauval); break;
448  default: cout << "Unknown metric_id: " << metric_id << endl; return 3;
449  }
451  TargetConstructor *target_c = NULL;
452  HessianCoefficient *adapt_coeff = NULL;
453  H1_FECollection ind_fec(3, dim);
454  FiniteElementSpace ind_fes(mesh, &ind_fec);
455  GridFunction size;
456  switch (target_id)
457  {
458  case 1: target_t = TargetConstructor::IDEAL_SHAPE_UNIT_SIZE; break;
459  case 2: target_t = TargetConstructor::IDEAL_SHAPE_EQUAL_SIZE; break;
460  case 3: target_t = TargetConstructor::IDEAL_SHAPE_GIVEN_SIZE; break;
461  case 4:
462  {
464  AnalyticAdaptTC *tc = new AnalyticAdaptTC(target_t);
465  adapt_coeff = new HessianCoefficient(dim, 1);
466  tc->SetAnalyticTargetSpec(NULL, NULL, adapt_coeff);
467  target_c = tc;
468  break;
469  }
470  case 5:
471  {
473  DiscreteAdaptTC *tc = new DiscreteAdaptTC(target_t);
474  size.SetSpace(&ind_fes);
475  FunctionCoefficient ind_coeff(ind_values);
476  size.ProjectCoefficient(ind_coeff);
477  tc->SetSerialDiscreteTargetSpec(size);
478  target_c = tc;
479  break;
480  }
481  default: cout << "Unknown target_id: " << target_id << endl; return 3;
482  }
483 
484  if (target_c == NULL)
485  {
486  target_c = new TargetConstructor(target_t);
487  }
488  target_c->SetNodes(x0);
489  TMOP_Integrator *he_nlf_integ = new TMOP_Integrator(metric, target_c);
490 
491  // 12. Setup the quadrature rule for the non-linear form integrator.
492  const IntegrationRule *ir = NULL;
493  const int geom_type = fespace->GetFE(0)->GetGeomType();
494  switch (quad_type)
495  {
496  case 1: ir = &IntRulesLo.Get(geom_type, quad_order); break;
497  case 2: ir = &IntRules.Get(geom_type, quad_order); break;
498  case 3: ir = &IntRulesCU.Get(geom_type, quad_order); break;
499  default: cout << "Unknown quad_type: " << quad_type << endl;
500  delete he_nlf_integ; return 3;
501  }
502  cout << "Quadrature points per cell: " << ir->GetNPoints() << endl;
503  he_nlf_integ->SetIntegrationRule(*ir);
504 
505  if (normalization) { he_nlf_integ->EnableNormalization(x0); }
506 
507  // 13. Limit the node movement.
508  // The limiting distances can be given by a general function of space.
509  GridFunction dist(fespace);
510  dist = 1.0;
511  // The small_phys_size is relevant only with proper normalization.
512  if (normalization) { dist = small_phys_size; }
513  ConstantCoefficient lim_coeff(lim_const);
514  if (lim_const != 0.0) { he_nlf_integ->EnableLimiting(x0, dist, lim_coeff); }
515 
516  // 14. Setup the final NonlinearForm (which defines the integral of interest,
517  // its first and second derivatives). Here we can use a combination of
518  // metrics, i.e., optimize the sum of two integrals, where both are
519  // scaled by used-defined space-dependent weights. Note that there are no
520  // command-line options for the weights and the type of the second
521  // metric; one should update those in the code.
522  NonlinearForm a(fespace);
523  ConstantCoefficient *coeff1 = NULL;
524  TMOP_QualityMetric *metric2 = NULL;
525  TargetConstructor *target_c2 = NULL;
527 
528  if (combomet == 1)
529  {
530  // TODO normalization of combinations.
531  // We will probably drop this example and replace it with adaptivity.
532  if (normalization) { MFEM_ABORT("Not implemented."); }
533 
534  // Weight of the original metric.
535  coeff1 = new ConstantCoefficient(1.0);
536  he_nlf_integ->SetCoefficient(*coeff1);
537  a.AddDomainIntegrator(he_nlf_integ);
538 
539  metric2 = new TMOP_Metric_077;
540  target_c2 = new TargetConstructor(
542  target_c2->SetVolumeScale(0.01);
543  target_c2->SetNodes(x0);
544  TMOP_Integrator *he_nlf_integ2 = new TMOP_Integrator(metric2, target_c2);
545  he_nlf_integ2->SetIntegrationRule(*ir);
546 
547  // Weight of metric2.
548  he_nlf_integ2->SetCoefficient(coeff2);
549  a.AddDomainIntegrator(he_nlf_integ2);
550  }
551  else { a.AddDomainIntegrator(he_nlf_integ); }
552 
553  const double init_energy = a.GetGridFunctionEnergy(x);
554 
555  // 15. Visualize the starting mesh and metric values.
556  if (visualization)
557  {
558  char title[] = "Initial metric values";
559  vis_tmop_metric_s(mesh_poly_deg, *metric, *target_c, *mesh, title, 0);
560  }
561 
562  // 16. Fix all boundary nodes, or fix only a given component depending on the
563  // boundary attributes of the given mesh. Attributes 1/2/3 correspond to
564  // fixed x/y/z components of the node. Attribute 4 corresponds to an
565  // entirely fixed node. Other boundary attributes do not affect the node
566  // movement boundary conditions.
567  if (move_bnd == false)
568  {
569  Array<int> ess_bdr(mesh->bdr_attributes.Max());
570  ess_bdr = 1;
571  a.SetEssentialBC(ess_bdr);
572  }
573  else
574  {
575  const int nd = fespace->GetBE(0)->GetDof();
576  int n = 0;
577  for (int i = 0; i < mesh->GetNBE(); i++)
578  {
579  const int attr = mesh->GetBdrElement(i)->GetAttribute();
580  MFEM_VERIFY(!(dim == 2 && attr == 3),
581  "Boundary attribute 3 must be used only for 3D meshes. "
582  "Adjust the attributes (1/2/3/4 for fixed x/y/z/all "
583  "components, rest for free nodes), or use -fix-bnd.");
584  if (attr == 1 || attr == 2 || attr == 3) { n += nd; }
585  if (attr == 4) { n += nd * dim; }
586  }
587  Array<int> ess_vdofs(n), vdofs;
588  n = 0;
589  for (int i = 0; i < mesh->GetNBE(); i++)
590  {
591  const int attr = mesh->GetBdrElement(i)->GetAttribute();
592  fespace->GetBdrElementVDofs(i, vdofs);
593  if (attr == 1) // Fix x components.
594  {
595  for (int j = 0; j < nd; j++)
596  { ess_vdofs[n++] = vdofs[j]; }
597  }
598  else if (attr == 2) // Fix y components.
599  {
600  for (int j = 0; j < nd; j++)
601  { ess_vdofs[n++] = vdofs[j+nd]; }
602  }
603  else if (attr == 3) // Fix z components.
604  {
605  for (int j = 0; j < nd; j++)
606  { ess_vdofs[n++] = vdofs[j+2*nd]; }
607  }
608  else if (attr == 4) // Fix all components.
609  {
610  for (int j = 0; j < vdofs.Size(); j++)
611  { ess_vdofs[n++] = vdofs[j]; }
612  }
613  }
614  a.SetEssentialVDofs(ess_vdofs);
615  }
616 
617  // 17. As we use the Newton method to solve the resulting nonlinear system,
618  // here we setup the linear solver for the system's Jacobian.
619  Solver *S = NULL;
620  const double linsol_rtol = 1e-12;
621  if (lin_solver == 0)
622  {
623  S = new DSmoother(1, 1.0, max_lin_iter);
624  }
625  else if (lin_solver == 1)
626  {
627  CGSolver *cg = new CGSolver;
628  cg->SetMaxIter(max_lin_iter);
629  cg->SetRelTol(linsol_rtol);
630  cg->SetAbsTol(0.0);
631  cg->SetPrintLevel(verbosity_level >= 2 ? 3 : -1);
632  S = cg;
633  }
634  else
635  {
636  MINRESSolver *minres = new MINRESSolver;
637  minres->SetMaxIter(max_lin_iter);
638  minres->SetRelTol(linsol_rtol);
639  minres->SetAbsTol(0.0);
640  minres->SetPrintLevel(verbosity_level >= 2 ? 3 : -1);
641  S = minres;
642  }
643 
644  // 18. Compute the minimum det(J) of the starting mesh.
645  tauval = infinity();
646  const int NE = mesh->GetNE();
647  for (int i = 0; i < NE; i++)
648  {
650  for (int j = 0; j < ir->GetNPoints(); j++)
651  {
652  transf->SetIntPoint(&ir->IntPoint(j));
653  tauval = min(tauval, transf->Jacobian().Det());
654  }
655  }
656  cout << "Minimum det(J) of the original mesh is " << tauval << endl;
657 
658  // 19. Finally, perform the nonlinear optimization.
659  NewtonSolver *newton = NULL;
660  if (tauval > 0.0)
661  {
662  tauval = 0.0;
663  TMOPNewtonSolver *tns = new TMOPNewtonSolver(*ir);
664  if (target_id == 5)
665  {
666  tns->SetDiscreteAdaptTC(dynamic_cast<DiscreteAdaptTC *>(target_c));
667  }
668  newton = tns;
669  cout << "TMOPNewtonSolver is used (as all det(J) > 0).\n";
670  }
671  else
672  {
673  if ( (dim == 2 && metric_id != 22 && metric_id != 252) ||
674  (dim == 3 && metric_id != 352) )
675  {
676  cout << "The mesh is inverted. Use an untangling metric." << endl;
677  return 3;
678  }
679  tauval -= 0.01 * h0.Min(); // Slightly below minJ0 to avoid div by 0.
680  newton = new TMOPDescentNewtonSolver(*ir);
681  cout << "The TMOPDescentNewtonSolver is used (as some det(J) < 0).\n";
682  }
683  newton->SetPreconditioner(*S);
684  newton->SetMaxIter(newton_iter);
685  newton->SetRelTol(newton_rtol);
686  newton->SetAbsTol(0.0);
687  newton->SetPrintLevel(verbosity_level >= 1 ? 1 : -1);
688  newton->SetOperator(a);
689  newton->Mult(b, x.GetTrueVector());
690  x.SetFromTrueVector();
691 
692  if (newton->GetConverged() == false)
693  {
694  cout << "NewtonIteration: rtol = " << newton_rtol << " not achieved."
695  << endl;
696  }
697  delete newton;
698 
699  // 20. Save the optimized mesh to a file. This output can be viewed later
700  // using GLVis: "glvis -m optimized.mesh".
701  {
702  ofstream mesh_ofs("optimized.mesh");
703  mesh_ofs.precision(14);
704  mesh->Print(mesh_ofs);
705  }
706 
707  // 21. Compute the amount of energy decrease.
708  const double fin_energy = a.GetGridFunctionEnergy(x);
709  double metric_part = fin_energy;
710  if (lim_const != 0.0)
711  {
712  lim_coeff.constant = 0.0;
713  metric_part = a.GetGridFunctionEnergy(x);
714  lim_coeff.constant = lim_const;
715  }
716  cout << "Initial strain energy: " << init_energy
717  << " = metrics: " << init_energy
718  << " + limiting term: " << 0.0 << endl;
719  cout << " Final strain energy: " << fin_energy
720  << " = metrics: " << metric_part
721  << " + limiting term: " << fin_energy - metric_part << endl;
722  cout << "The strain energy decreased by: " << setprecision(12)
723  << (init_energy - fin_energy) * 100.0 / init_energy << " %." << endl;
724 
725  // 22. Visualize the final mesh and metric values.
726  if (visualization)
727  {
728  char title[] = "Final metric values";
729  vis_tmop_metric_s(mesh_poly_deg, *metric, *target_c, *mesh, title, 600);
730  }
731 
732  // 23. Visualize the mesh displacement.
733  if (visualization)
734  {
735  x0 -= x;
736  osockstream sock(19916, "localhost");
737  sock << "solution\n";
738  mesh->Print(sock);
739  x0.Save(sock);
740  sock.send();
741  sock << "window_title 'Displacements'\n"
742  << "window_geometry "
743  << 1200 << " " << 0 << " " << 600 << " " << 600 << "\n"
744  << "keys jRmclA" << endl;
745  }
746 
747  // 24. Free the used memory.
748  delete S;
749  delete target_c2;
750  delete metric2;
751  delete coeff1;
752  delete target_c;
753  delete metric;
754  delete fespace;
755  delete fec;
756  delete mesh;
757 
758  return 0;
759 }
760 
761 // Defined with respect to the icf mesh.
762 double weight_fun(const Vector &x)
763 {
764  const double r = sqrt(x(0)*x(0) + x(1)*x(1) + 1e-12);
765  const double den = 0.002;
766  double l2 = 0.2 + 0.5*std::tanh((r-0.16)/den) - 0.5*std::tanh((r-0.17)/den)
767  + 0.5*std::tanh((r-0.23)/den) - 0.5*std::tanh((r-0.24)/den);
768  return l2;
769 }
int GetNPoints() const
Returns the number of the points in the integration rule.
Definition: intrules.hpp:245
int Size() const
Logical size of the array.
Definition: array.hpp:124
Shifted barrier form of metric 56 (area, ideal barrier metric), 2D.
Definition: tmop.hpp:324
Shifted barrier form of 3D metric 16 (volume, ideal barrier metric), 3D.
Definition: tmop.hpp:442
virtual void Print(std::ostream &out=mfem::out) const
Definition: mesh.hpp:1188
Conjugate gradient method.
Definition: solvers.hpp:152
int GetNDofs() const
Returns number of degrees of freedom.
Definition: fespace.hpp:378
Class for an integration rule - an Array of IntegrationPoint.
Definition: intrules.hpp:90
Class for grid function - Vector with associated FE space.
Definition: gridfunc.hpp:27
int DofToVDof(int dof, int vd, int ndofs=-1) const
Definition: fespace.cpp:143
Shape &amp; volume, ideal barrier metric, 3D.
Definition: tmop.hpp:426
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:890
void SetFromTrueVector()
Shortcut for calling SetFromTrueDofs() with GetTrueVector() as argument.
Definition: gridfunc.hpp:139
Subclass constant coefficient.
Definition: coefficient.hpp:67
void EnableNormalization(const GridFunction &x)
Computes the normalization factors of the metric and limiting integrals using the mesh position given...
Definition: tmop.cpp:1366
int GetNBE() const
Returns number of boundary elements.
Definition: mesh.hpp:696
Shape, ideal barrier metric, 3D.
Definition: tmop.hpp:376
double Det() const
Definition: densemat.cpp:449
void AddDomainIntegrator(NonlinearFormIntegrator *nlfi)
Adds new Domain Integrator.
Allows negative Jacobians. Used for untangling.
Definition: tmop_tools.hpp:106
void SetIntPoint(const IntegrationPoint *ip)
Definition: eltrans.hpp:56
Container class for integration rules.
Definition: intrules.hpp:307
Data type dense matrix using column-major storage.
Definition: densemat.hpp:23
Shape, ideal barrier metric, 3D.
Definition: tmop.hpp:360
void SetVolumeScale(double vol_scale)
Used by target type IDEAL_SHAPE_EQUAL_SIZE. The default volume scale is 1.
Definition: tmop.hpp:635
Shape, ideal barrier metric, 2D.
Definition: tmop.hpp:151
Volume metric, 3D.
Definition: tmop.hpp:392
Area, ideal barrier metric, 2D.
Definition: tmop.hpp:288
int GetNE() const
Returns number of elements.
Definition: mesh.hpp:693
virtual void SetSerialDiscreteTargetSpec(GridFunction &tspec)
Definition: tmop.cpp:943
void vis_tmop_metric_s(int order, TMOP_QualityMetric &qm, const TargetConstructor &tc, Mesh &mesh, char *title, int position)
Definition: tmop_tools.cpp:501
MINRES method.
Definition: solvers.hpp:262
void Randomize(int seed=0)
Set random values in the vector.
Definition: vector.cpp:692
int main(int argc, char *argv[])
Definition: ex1.cpp:62
double weight_fun(const Vector &x)
Shape &amp; area, ideal barrier metric, 2D.
Definition: tmop.hpp:167
Shifted barrier form of metric 2 (shape, ideal barrier metric), 2D.
Definition: tmop.hpp:199
Shape, ideal barrier metric, 2D.
Definition: tmop.hpp:270
virtual void Save(std::ostream &out) const
Save the GridFunction to an output stream.
Definition: gridfunc.cpp:2731
Geometry::Type GetGeomType() const
Returns the Geometry::Type of the reference element.
Definition: fe.hpp:311
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:133
void SetPrintLevel(int print_lvl)
Definition: solvers.cpp:70
void EnableLimiting(const GridFunction &n0, const GridFunction &dist, Coefficient &w0, TMOP_LimiterFunction *lfunc=NULL)
Adds a limiting term to the integrator (general version).
Definition: tmop.cpp:1038
int GetNBE() const
Returns number of boundary elements in the mesh.
Definition: fespace.hpp:417
void SetCoefficient(Coefficient &w1)
Sets a scaling Coefficient for the quality metric term of the integrator.
Definition: tmop.hpp:784
Area metric, 2D.
Definition: tmop.hpp:234
Volume, ideal barrier metric, 3D.
Definition: tmop.hpp:408
double b
Definition: lissajous.cpp:42
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:7982
IntegrationRules IntRulesLo(0, Quadrature1D::GaussLobatto)
const DenseMatrix & Jacobian()
Return the Jacobian matrix of the transformation at the currently set IntegrationPoint, using the method SetIntPoint().
Definition: eltrans.hpp:71
void SetNodes(const GridFunction &n)
Set the nodes to be used in the target-matrix construction.
Definition: tmop.hpp:632
void SetMaxIter(int max_it)
Definition: solvers.hpp:65
T Max() const
Find the maximal element in the array, using the comparison operator &lt; for class T.
Definition: array.cpp:68
virtual void GetElementDofs(int i, Array< int > &dofs) const
Returns indexes of degrees of freedom in array dofs for i&#39;th element.
Definition: fespace.cpp:1540
Newton&#39;s method for solving F(x)=b for a given operator F.
Definition: solvers.hpp:300
Version of QuadraticFECollection with positive basis functions.
Definition: fe_coll.hpp:394
int GetAttribute() const
Return element&#39;s attribute.
Definition: element.hpp:55
double ind_values(const Vector &x)
void SetNodalFESpace(FiniteElementSpace *nfes)
Definition: mesh.cpp:3885
int GetConverged() const
Definition: solvers.hpp:69
int Dimension() const
Definition: mesh.hpp:744
void PrintUsage(std::ostream &out) const
Definition: optparser.cpp:434
virtual void Mult(const Vector &b, Vector &x) const
Solve the nonlinear system with right-hand side b.
Definition: solvers.cpp:1339
const Vector & GetTrueVector() const
Read only access to the (optional) internal true-dof Vector.
Definition: gridfunc.hpp:121
double GetElementSize(int i, int type=0)
Get the size of the i-th element relative to the perfect reference element.
Definition: mesh.cpp:76
void SetAbsTol(double atol)
Definition: solvers.hpp:64
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
Definition: mesh.hpp:189
void SetRelTol(double rtol)
Definition: solvers.hpp:63
Shape, ideal barrier metric, 2D.
Definition: tmop.hpp:218
Untangling metric, 2D.
Definition: tmop.hpp:305
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition: fespace.hpp:87
int GetDof() const
Returns the number of degrees of freedom in the finite element.
Definition: fe.hpp:314
Area, ideal barrier metric, 2D.
Definition: tmop.hpp:251
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)
Definition: optparser.hpp:76
virtual void SetAnalyticTargetSpec(Coefficient *sspec, VectorCoefficient *vspec, MatrixCoefficient *mspec)
Definition: tmop.cpp:882
Shape &amp; area metric, 2D.
Definition: tmop.hpp:183
double a
Definition: lissajous.cpp:41
Class for integration point with weight.
Definition: intrules.hpp:25
double GetGridFunctionEnergy(const Vector &x) const
Compute the enery corresponding to the state x.
void GetElementTransformation(int i, IsoparametricTransformation *ElTr)
Definition: mesh.cpp:338
void SetIntegrationRule(const IntegrationRule &irule)
Prescribe a fixed IntegrationRule to use.
Definition: nonlininteg.hpp:41
int dim
Definition: ex24.cpp:43
const FiniteElement * GetFE(int i) const
Returns pointer to the FiniteElement associated with i&#39;th element.
Definition: fespace.cpp:1642
void PrintOptions(std::ostream &out) const
Definition: optparser.cpp:304
double infinity()
Define a shortcut for std::numeric_limits&lt;double&gt;::infinity()
Definition: vector.hpp:42
virtual void ProjectCoefficient(Coefficient &coeff)
Definition: gridfunc.cpp:1685
virtual void SetSpace(FiniteElementSpace *f)
Associate a new FiniteElementSpace with the GridFunction.
Definition: gridfunc.cpp:180
void SetEssentialBC(const Array< int > &bdr_attr_is_ess, Vector *rhs=NULL)
Specify essential boundary conditions.
class for C-function coefficient
Vector data type.
Definition: vector.hpp:48
virtual void Transform(const IntegrationPoint &, Vector &)=0
IntegrationRules IntRulesCU(0, Quadrature1D::ClosedUniform)
void GetNodes(Vector &node_coord) const
Definition: mesh.cpp:6203
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
Definition: solvers.cpp:91
Arbitrary order H1-conforming (continuous) finite elements.
Definition: fe_coll.hpp:83
TargetType
Target-matrix construction algorithms supported by this class.
Definition: tmop.hpp:574
void SetDiscreteAdaptTC(DiscreteAdaptTC *tc)
Definition: tmop_tools.hpp:98
void SetEssentialVDofs(const Array< int > &ess_vdofs_list)
(DEPRECATED) Specify essential boundary conditions.
Base class for solvers.
Definition: operator.hpp:500
const FiniteElement * GetBE(int i) const
Returns pointer to the FiniteElement for the i&#39;th boundary element.
Definition: fespace.cpp:1882
void GetBdrElementVDofs(int i, Array< int > &vdofs) const
Returns indexes of degrees of freedom for i&#39;th boundary element.
Definition: fespace.cpp:178
void SetNodalGridFunction(GridFunction *nodes, bool make_owner=false)
Definition: mesh.cpp:3913
Base class representing target-matrix construction algorithms for mesh optimization via the target-ma...
Definition: tmop.hpp:570
IntegrationRules IntRules(0, Quadrature1D::GaussLegendre)
A global object with all integration rules (defined in intrules.cpp)
Definition: intrules.hpp:376
double GetElementVolume(int i)
Definition: mesh.cpp:103
const Element * GetBdrElement(int i) const
Definition: mesh.hpp:779
aka closed Newton-Cotes
Definition: intrules.hpp:295
Shape, ideal barrier metric, 3D.
Definition: tmop.hpp:344
Metric without a type, 2D.
Definition: tmop.hpp:75
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
Definition: solvers.cpp:1328
bool Good() const
Definition: optparser.hpp:125
A TMOP integrator class based on any given TMOP_QualityMetric and TargetConstructor.
Definition: tmop.hpp:728