MFEM  v4.4.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
tmop_tools.cpp
Go to the documentation of this file.
1 // Copyright (c) 2010-2022, Lawrence Livermore National Security, LLC. Produced
2 // at the Lawrence Livermore National Laboratory. All Rights reserved. See files
3 // LICENSE and NOTICE for details. LLNL-CODE-806117.
4 //
5 // This file is part of the MFEM library. For more information and source code
6 // availability visit https://mfem.org.
7 //
8 // MFEM is free software; you can redistribute it and/or modify it under the
9 // terms of the BSD-3 license. We welcome feedback and contributions, see file
10 // CONTRIBUTING.md for details.
11 
12 #include "tmop_tools.hpp"
13 #include "nonlinearform.hpp"
14 #include "pnonlinearform.hpp"
15 #include "../general/osockstream.hpp"
16 
17 namespace mfem
18 {
19 
20 using namespace mfem;
21 
22 void AdvectorCG::SetInitialField(const Vector &init_nodes,
23  const Vector &init_field)
24 {
25  nodes0 = init_nodes;
26  field0 = init_field;
27 }
28 
30  Vector &new_field)
31 {
32  // TODO: Implement for AMR meshes.
33  const int pnt_cnt = new_field.Size()/ncomp;
34 
35  new_field = field0;
36  Vector new_field_temp;
37  for (int i = 0; i < ncomp; i++)
38  {
39  new_field_temp.MakeRef(new_field, i*pnt_cnt, pnt_cnt);
40  ComputeAtNewPositionScalar(new_nodes, new_field_temp);
41  }
42 
43  field0 = new_field;
44  nodes0 = new_nodes;
45 }
46 
47 void AdvectorCG::ComputeAtNewPositionScalar(const Vector &new_nodes,
48  Vector &new_field)
49 {
50  Mesh *m = mesh;
51 #ifdef MFEM_USE_MPI
52  if (pmesh) { m = pmesh; }
53 #endif
54 
55  MFEM_VERIFY(m != NULL, "No mesh has been given to the AdaptivityEvaluator.");
56 
57  // This will be used to move the positions.
58  GridFunction *mesh_nodes = m->GetNodes();
59  *mesh_nodes = nodes0;
60  double minv = new_field.Min(), maxv = new_field.Max();
61 
62  // Velocity of the positions.
63  GridFunction u(mesh_nodes->FESpace());
64  subtract(new_nodes, nodes0, u);
65 
66  // Define a scalar FE space for the solution, and the advection operator.
67  TimeDependentOperator *oper = NULL;
68  FiniteElementSpace *fess = NULL;
69 #ifdef MFEM_USE_MPI
70  ParFiniteElementSpace *pfess = NULL;
71 #endif
72  if (fes)
73  {
74  fess = new FiniteElementSpace(fes->GetMesh(), fes->FEColl(), 1);
75  oper = new SerialAdvectorCGOper(nodes0, u, *fess, al);
76  }
77 #ifdef MFEM_USE_MPI
78  else if (pfes)
79  {
80  pfess = new ParFiniteElementSpace(pfes->GetParMesh(), pfes->FEColl(), 1);
81  oper = new ParAdvectorCGOper(nodes0, u, *pfess, al, opt_mt);
82  }
83 #endif
84  MFEM_VERIFY(oper != NULL,
85  "No FE space has been given to the AdaptivityEvaluator.");
86  ode_solver.Init(*oper);
87 
88  // Compute some time step [mesh_size / speed].
90  for (int i = 0; i < m->GetNE(); i++)
91  {
92  h_min = std::min(h_min, m->GetElementSize(i));
93  }
94  double v_max = 0.0;
95  const int s = new_field.Size();
96 
97  u.HostReadWrite();
98  for (int i = 0; i < s; i++)
99  {
100  double vel = 0.;
101  for (int j = 0; j < dim; j++)
102  {
103  vel += u(i+j*s)*u(i+j*s);
104  }
105  v_max = std::max(v_max, vel);
106  }
107 
108 #ifdef MFEM_USE_MPI
109  if (pfes)
110  {
111  double v_loc = v_max, h_loc = h_min;
112  MPI_Allreduce(&v_loc, &v_max, 1, MPI_DOUBLE, MPI_MAX, pfes->GetComm());
113  MPI_Allreduce(&h_loc, &h_min, 1, MPI_DOUBLE, MPI_MIN, pfes->GetComm());
114  }
115 #endif
116 
117  if (v_max == 0.0) // No need to change the field.
118  {
119  delete oper;
120  delete fess;
121 #ifdef MFEM_USE_MPI
122  delete pfess;
123 #endif
124  return;
125  }
126 
127  v_max = std::sqrt(v_max);
128  double dt = dt_scale * h_min / v_max;
129 
130  double t = 0.0;
131  bool last_step = false;
132  for (int ti = 1; !last_step; ti++)
133  {
134  if (t + dt >= 1.0)
135  {
136  dt = 1.0 - t;
137  last_step = true;
138  }
139  ode_solver.Step(new_field, t, dt);
140  }
141 
142  double glob_minv = minv,
143  glob_maxv = maxv;
144 #ifdef MFEM_USE_MPI
145  if (pfes)
146  {
147  MPI_Allreduce(&minv, &glob_minv, 1, MPI_DOUBLE, MPI_MIN, pfes->GetComm());
148  MPI_Allreduce(&maxv, &glob_maxv, 1, MPI_DOUBLE, MPI_MAX, pfes->GetComm());
149  }
150 #endif
151 
152  // Trim the overshoots and undershoots.
153  new_field.HostReadWrite();
154  for (int i = 0; i < s; i++)
155  {
156  if (new_field(i) < glob_minv) { new_field(i) = glob_minv; }
157  if (new_field(i) > glob_maxv) { new_field(i) = glob_maxv; }
158  }
159 
160  delete oper;
161  delete fess;
162 #ifdef MFEM_USE_MPI
163  delete pfess;
164 #endif
165 }
166 
168  GridFunction &vel,
169  FiniteElementSpace &fes,
170  AssemblyLevel al)
171  : TimeDependentOperator(fes.GetVSize()),
172  x0(x_start), x_now(*fes.GetMesh()->GetNodes()),
173  u(vel), u_coeff(&u), M(&fes), K(&fes), al(al)
174 {
176  K.AddDomainIntegrator(Kinteg);
177  K.SetAssemblyLevel(al);
178  K.Assemble(0);
179  K.Finalize(0);
180 
181  MassIntegrator *Minteg = new MassIntegrator;
182  M.AddDomainIntegrator(Minteg);
183  M.SetAssemblyLevel(al);
184  M.Assemble(0);
185  M.Finalize(0);
186 }
187 
188 void SerialAdvectorCGOper::Mult(const Vector &ind, Vector &di_dt) const
189 {
190  // Move the mesh.
191  const double t = GetTime();
192  add(x0, t, u, x_now);
193 
194  if (al == AssemblyLevel::PARTIAL)
195  {
197  }
198 
199  // Assemble on the new mesh.
200  K.BilinearForm::operator=(0.0);
201  K.Assemble();
202  Vector rhs(K.Size());
203  K.Mult(ind, rhs);
204  M.BilinearForm::operator=(0.0);
205  M.Assemble();
206 
207  di_dt = 0.0;
208  CGSolver lin_solver;
209  Solver *prec = nullptr;
210  Array<int> ess_tdof_list;
211  if (al == AssemblyLevel::PARTIAL)
212  {
213  prec = new OperatorJacobiSmoother(M, ess_tdof_list);
214  lin_solver.SetOperator(M);
215  }
216  else
217  {
218  prec = new DSmoother(M.SpMat());
219  lin_solver.SetOperator(M.SpMat());
220  }
221  lin_solver.SetPreconditioner(*prec);
222  lin_solver.SetRelTol(1e-12); lin_solver.SetAbsTol(0.0);
223  lin_solver.SetMaxIter(100);
224  lin_solver.SetPrintLevel(0);
225  lin_solver.Mult(rhs, di_dt);
226 
227  delete prec;
228 }
229 
230 #ifdef MFEM_USE_MPI
232  GridFunction &vel,
233  ParFiniteElementSpace &pfes,
234  AssemblyLevel al,
235  MemoryType mt)
236  : TimeDependentOperator(pfes.GetVSize()),
237  x0(x_start), x_now(*pfes.GetMesh()->GetNodes()),
238  u(vel), u_coeff(&u), M(&pfes), K(&pfes), al(al)
239 {
241  if (al == AssemblyLevel::PARTIAL)
242  {
243  Kinteg->SetPAMemoryType(mt);
244  }
245  K.AddDomainIntegrator(Kinteg);
246  K.SetAssemblyLevel(al);
247  K.Assemble(0);
248  K.Finalize(0);
249 
250  MassIntegrator *Minteg = new MassIntegrator;
251  if (al == AssemblyLevel::PARTIAL)
252  {
253  Minteg->SetPAMemoryType(mt);
254  }
255  M.AddDomainIntegrator(Minteg);
256  M.SetAssemblyLevel(al);
257  M.Assemble(0);
258  M.Finalize(0);
259 }
260 
261 void ParAdvectorCGOper::Mult(const Vector &ind, Vector &di_dt) const
262 {
263  // Move the mesh.
264  const double t = GetTime();
265  add(x0, t, u, x_now);
266 
267  if (al == AssemblyLevel::PARTIAL)
268  {
270  }
271 
272  // Assemble on the new mesh.
273  K.BilinearForm::operator=(0.0);
274  K.Assemble();
276  K.Mult(ind, rhs);
277  M.BilinearForm::operator=(0.0);
278  M.Assemble();
279 
280  HypreParVector *RHS = rhs.ParallelAssemble();
282  X = 0.0;
283 
284  OperatorHandle Mop;
285  Solver *prec = nullptr;
286  Array<int> ess_tdof_list;
287  if (al == AssemblyLevel::PARTIAL)
288  {
289  M.FormSystemMatrix(ess_tdof_list, Mop);
290  prec = new OperatorJacobiSmoother(M, ess_tdof_list);
291  }
292  else
293  {
294  Mop.Reset(M.ParallelAssemble());
295  prec = new HypreSmoother;
296  static_cast<HypreSmoother*>(prec)->SetType(HypreSmoother::Jacobi, 1);
297  }
298 
299  CGSolver lin_solver(M.ParFESpace()->GetParMesh()->GetComm());
300  lin_solver.SetPreconditioner(*prec);
301  lin_solver.SetOperator(*Mop);
302  lin_solver.SetRelTol(1e-8);
303  lin_solver.SetAbsTol(0.0);
304  lin_solver.SetMaxIter(100);
305  lin_solver.SetPrintLevel(0);
306  lin_solver.Mult(*RHS, X);
307  K.ParFESpace()->GetProlongationMatrix()->Mult(X, di_dt);
308 
309  delete RHS;
310  delete prec;
311 }
312 #endif
313 
314 #ifdef MFEM_USE_GSLIB
316  const Vector &init_field)
317 {
318  nodes0 = init_nodes;
319  Mesh *m = mesh;
320 #ifdef MFEM_USE_MPI
321  if (pmesh) { m = pmesh; }
322 #endif
323  m->SetNodes(nodes0);
324 
325  const double rel_bbox_el = 0.1;
326  const double newton_tol = 1.0e-12;
327  const int npts_at_once = 256;
328 
329  if (finder)
330  {
331  finder->FreeData();
332  delete finder;
333  }
334 
336 #ifdef MFEM_USE_MPI
337  if (pfes)
338  {
339  f = pfes;
340  finder = new FindPointsGSLIB(pfes->GetComm());
341  }
342  else { finder = new FindPointsGSLIB(); }
343 #else
344  finder = new FindPointsGSLIB();
345 #endif
346  finder->Setup(*m, rel_bbox_el, newton_tol, npts_at_once);
347 
348  field0_gf.SetSpace(f);
349  field0_gf = init_field;
350 
351  dim = f->GetFE(0)->GetDim();
352 }
353 
355  Vector &new_field)
356 {
357  finder->Interpolate(new_nodes, field0_gf, new_field);
358 }
359 
360 #endif
361 
363  const Vector &b) const
364 {
365  const FiniteElementSpace *fes = NULL;
366  double energy_in = 0.0;
367 #ifdef MFEM_USE_MPI
368  const ParNonlinearForm *p_nlf = dynamic_cast<const ParNonlinearForm *>(oper);
369  MFEM_VERIFY(!(parallel && p_nlf == NULL), "Invalid Operator subclass.");
370  if (parallel)
371  {
372  fes = p_nlf->FESpace();
373  energy_in = p_nlf->GetEnergy(x);
374  }
375 #endif
376  const bool serial = !parallel;
377  const NonlinearForm *nlf = dynamic_cast<const NonlinearForm *>(oper);
378  MFEM_VERIFY(!(serial && nlf == NULL), "Invalid Operator subclass.");
379  if (serial)
380  {
381  fes = nlf->FESpace();
382  energy_in = nlf->GetEnergy(x);
383  }
384 
385  // Get the local prolongation of the solution vector.
386  Vector x_out_loc(fes->GetVSize(),
388  if (serial)
389  {
390  const SparseMatrix *cP = fes->GetConformingProlongation();
391  if (!cP) { x_out_loc = x; }
392  else { cP->Mult(x, x_out_loc); }
393  }
394 #ifdef MFEM_USE_MPI
395  else
396  {
397  fes->GetProlongationMatrix()->Mult(x, x_out_loc);
398  }
399 #endif
400 
401  double scale = 1.0;
402  if (surf_fit_max_threshold > 0.0)
403  {
404  double avg_err, max_err;
405  GetSurfaceFittingError(avg_err, max_err);
406  if (max_err < surf_fit_max_threshold)
407  {
409  {
410  mfem::out << "TMOPNewtonSolver converged "
411  "based on the surface fitting error.\n";
412  }
413  scale = 0.0;
414  return scale;
415  }
416  }
417 
418  // Check if the starting mesh (given by x) is inverted. Note that x hasn't
419  // been modified by the Newton update yet.
420  const double min_detT_in = ComputeMinDet(x_out_loc, *fes);
421  const bool untangling = (min_detT_in <= 0.0) ? true : false;
422  const double untangle_factor = 1.5;
423  if (untangling)
424  {
425  // Needed for the line search below. The untangling metrics see this
426  // reference to detect deteriorations.
427  MFEM_VERIFY(min_det_ptr != NULL, " Initial mesh was valid, but"
428  " intermediate mesh is invalid. Contact TMOP Developers.");
429  *min_det_ptr = untangle_factor * min_detT_in;
430  }
431 
432  const bool have_b = (b.Size() == Height());
433 
434  Vector x_out(x.Size());
435  bool x_out_ok = false;
436  double energy_out = 0.0, min_detT_out;
437  const double norm_in = Norm(r);
438 
439  const double detJ_factor = (solver_type == 1) ? 0.25 : 0.5;
440 
441  // Perform the line search.
442  for (int i = 0; i < 12; i++)
443  {
444  // Update the mesh and get the L-vector in x_out_loc.
445  add(x, -scale, c, x_out);
446  if (serial)
447  {
448  const SparseMatrix *cP = fes->GetConformingProlongation();
449  if (!cP) { x_out_loc = x_out; }
450  else { cP->Mult(x_out, x_out_loc); }
451  }
452 #ifdef MFEM_USE_MPI
453  else
454  {
455  fes->GetProlongationMatrix()->Mult(x_out, x_out_loc);
456  }
457 #endif
458 
459  // Check the changes in detJ.
460  min_detT_out = ComputeMinDet(x_out_loc, *fes);
461  if (untangling == false && min_detT_out < 0.0)
462  {
463  // No untangling, and detJ got negative -- no good.
465  {
466  mfem::out << "Scale = " << scale << " Neg det(J) found.\n";
467  }
468  scale *= detJ_factor; continue;
469  }
470  if (untangling == true && min_detT_out < *min_det_ptr)
471  {
472  // Untangling, and detJ got even more negative -- no good.
474  {
475  mfem::out << "Scale = " << scale << " Neg det(J) decreased.\n";
476  }
477  scale *= detJ_factor; continue;
478  }
479 
480  // Skip the energy and residual checks when we're untangling. The
481  // untangling metrics change their denominators, which can affect the
482  // energy and residual, so their increase/decrease is not relevant.
483  if (untangling) { x_out_ok = true; break; }
484 
485  // Check the changes in total energy.
486  ProcessNewState(x_out);
487 
488  if (serial)
489  {
490  energy_out = nlf->GetGridFunctionEnergy(x_out_loc);
491  }
492 #ifdef MFEM_USE_MPI
493  else
494  {
495  energy_out = p_nlf->GetParGridFunctionEnergy(x_out_loc);
496  }
497 #endif
498  if (energy_out > energy_in + 0.2*fabs(energy_in) ||
499  std::isnan(energy_out) != 0)
500  {
502  {
503  mfem::out << "Scale = " << scale << " Increasing energy: "
504  << energy_in << " --> " << energy_out << '\n';
505  }
506  scale *= 0.5; continue;
507  }
508 
509  // Check the changes in the Newton residual.
510  oper->Mult(x_out, r);
511  if (have_b) { r -= b; }
512  double norm_out = Norm(r);
513 
514  if (norm_out > 1.2*norm_in)
515  {
517  {
518  mfem::out << "Scale = " << scale << " Norm increased: "
519  << norm_in << " --> " << norm_out << '\n';
520  }
521  scale *= 0.5; continue;
522  }
523  else { x_out_ok = true; break; }
524  } // end line search
525 
526  if (untangling)
527  {
528  // Update the global min detJ. Untangling metrics see this min_det_ptr.
529  if (min_detT_out > 0.0)
530  {
531  *min_det_ptr = 0.0;
534  { mfem::out << "The mesh has been untangled at the used points!\n"; }
535  }
536  else { *min_det_ptr = untangle_factor * min_detT_out; }
537  }
538 
541  {
542  if (untangling)
543  {
544  mfem::out << "Min det(T) change: "
545  << min_detT_in << " -> " << min_detT_out
546  << " with " << scale << " scaling.\n";
547  }
548  else
549  {
550  mfem::out << "Energy decrease: "
551  << energy_in << " --> " << energy_out << " or "
552  << (energy_in - energy_out) / energy_in * 100.0
553  << "% with " << scale << " scaling.\n";
554  }
555  }
556 
557  if (x_out_ok == false) { scale = 0.0; }
558 
560 
561  return scale;
562 }
563 
565 {
566  const NonlinearForm *nlf = dynamic_cast<const NonlinearForm *>(oper);
567  const Array<NonlinearFormIntegrator*> &integs = *nlf->GetDNFI();
568  TMOP_Integrator *ti = NULL;
569  TMOPComboIntegrator *co = NULL;
570  for (int i = 0; i < integs.Size(); i++)
571  {
572  ti = dynamic_cast<TMOP_Integrator *>(integs[i]);
573  if (ti)
574  {
575  ti->UpdateSurfaceFittingWeight(factor);
576  }
577  co = dynamic_cast<TMOPComboIntegrator *>(integs[i]);
578  if (co)
579  {
581  for (int j = 0; j < ati.Size(); j++)
582  {
583  ati[j]->UpdateSurfaceFittingWeight(factor);
584  }
585  }
586  }
587 }
588 
590 {
591  const NonlinearForm *nlf = dynamic_cast<const NonlinearForm *>(oper);
592  const Array<NonlinearFormIntegrator*> &integs = *nlf->GetDNFI();
593  TMOP_Integrator *ti = NULL;
594  TMOPComboIntegrator *co = NULL;
595  weights.SetSize(0);
596  double weight;
597 
598  for (int i = 0; i < integs.Size(); i++)
599  {
600  ti = dynamic_cast<TMOP_Integrator *>(integs[i]);
601  if (ti)
602  {
603  weight = ti->GetSurfaceFittingWeight();
604  weights.Append(weight);
605  }
606  co = dynamic_cast<TMOPComboIntegrator *>(integs[i]);
607  if (co)
608  {
610  for (int j = 0; j < ati.Size(); j++)
611  {
612  weight = ati[j]->GetSurfaceFittingWeight();
613  weights.Append(weight);
614  }
615  }
616  }
617 }
618 
620  double &err_max) const
621 {
622  const NonlinearForm *nlf = dynamic_cast<const NonlinearForm *>(oper);
623  const Array<NonlinearFormIntegrator*> &integs = *nlf->GetDNFI();
624  TMOP_Integrator *ti = NULL;
625  TMOPComboIntegrator *co = NULL;
626 
627  err_avg = 0.0;
628  err_max = 0.0;
629  double err_avg_loc, err_max_loc;
630  for (int i = 0; i < integs.Size(); i++)
631  {
632  ti = dynamic_cast<TMOP_Integrator *>(integs[i]);
633  if (ti)
634  {
635  if (ti->IsSurfaceFittingEnabled())
636  {
637  ti->GetSurfaceFittingErrors(err_avg_loc, err_max_loc);
638  err_avg = std::fmax(err_avg_loc, err_avg);
639  err_max = std::fmax(err_max_loc, err_max);
640  }
641  }
642  co = dynamic_cast<TMOPComboIntegrator *>(integs[i]);
643  if (co)
644  {
646  for (int j = 0; j < ati.Size(); j++)
647  {
648  if (ati[j]->IsSurfaceFittingEnabled())
649  {
650  ati[j]->GetSurfaceFittingErrors(err_avg_loc, err_max_loc);
651  err_avg = std::fmax(err_avg_loc, err_avg);
652  err_max = std::fmax(err_max_loc, err_max);
653  }
654  }
655  }
656  }
657 }
658 
660 {
661  const NonlinearForm *nlf = dynamic_cast<const NonlinearForm *>(oper);
662  const Array<NonlinearFormIntegrator*> &integs = *nlf->GetDNFI();
663 
664  // Reset the update flags of all TargetConstructors. This is done to avoid
665  // repeated updates of shared TargetConstructors.
666  TMOP_Integrator *ti = NULL;
667  TMOPComboIntegrator *co = NULL;
668  DiscreteAdaptTC *dtc = NULL;
669  for (int i = 0; i < integs.Size(); i++)
670  {
671  ti = dynamic_cast<TMOP_Integrator *>(integs[i]);
672  if (ti)
673  {
674  dtc = ti->GetDiscreteAdaptTC();
675  if (dtc) { dtc->ResetUpdateFlags(); }
676  }
677  co = dynamic_cast<TMOPComboIntegrator *>(integs[i]);
678  if (co)
679  {
681  for (int j = 0; j < ati.Size(); j++)
682  {
683  dtc = ati[j]->GetDiscreteAdaptTC();
684  if (dtc) { dtc->ResetUpdateFlags(); }
685  }
686  }
687  }
688 
689  if (parallel)
690  {
691 #ifdef MFEM_USE_MPI
692  const ParNonlinearForm *pnlf =
693  dynamic_cast<const ParNonlinearForm *>(oper);
694  const ParFiniteElementSpace *pfesc = pnlf->ParFESpace();
695  Vector x_loc(pfesc->GetVSize());
696  pfesc->GetProlongationMatrix()->Mult(x, x_loc);
697  for (int i = 0; i < integs.Size(); i++)
698  {
699  ti = dynamic_cast<TMOP_Integrator *>(integs[i]);
700  if (ti)
701  {
703  ti->ComputeFDh(x_loc, *pfesc);
704  UpdateDiscreteTC(*ti, x_loc);
705  }
706  co = dynamic_cast<TMOPComboIntegrator *>(integs[i]);
707  if (co)
708  {
710  for (int j = 0; j < ati.Size(); j++)
711  {
712  ati[j]->UpdateAfterMeshPositionChange(x_loc);
713  ati[j]->ComputeFDh(x_loc, *pfesc);
714  UpdateDiscreteTC(*ati[j], x_loc);
715  }
716  }
717  }
718 #endif
719  }
720  else
721  {
722  const FiniteElementSpace *fesc = nlf->FESpace();
723  const Operator *P = nlf->GetProlongation();
724  Vector x_loc;
725  if (P)
726  {
727  x_loc.SetSize(P->Height());
728  P->Mult(x,x_loc);
729  }
730  else
731  {
732  x_loc = x;
733  }
734  for (int i = 0; i < integs.Size(); i++)
735  {
736  ti = dynamic_cast<TMOP_Integrator *>(integs[i]);
737  if (ti)
738  {
740  ti->ComputeFDh(x_loc, *fesc);
741  UpdateDiscreteTC(*ti, x_loc);
742  }
743  co = dynamic_cast<TMOPComboIntegrator *>(integs[i]);
744  if (co)
745  {
747  for (int j = 0; j < ati.Size(); j++)
748  {
749  ati[j]->UpdateAfterMeshPositionChange(x_loc);
750  ati[j]->ComputeFDh(x_loc, *fesc);
751  UpdateDiscreteTC(*ati[j], x_loc);
752  }
753  }
754  }
755  }
756 
757  // Constant coefficient associated with the surface fitting terms if
758  // adaptive surface fitting is enabled. The idea is to increase the
759  // coefficient if the surface fitting error does not sufficiently
760  // decrease between subsequent TMOPNewtonSolver iterations.
762  {
763  double surf_fit_err_max = -10;
764  double surf_fit_err_avg = -10;
765  // Get surface fitting errors.
766  GetSurfaceFittingError(surf_fit_err_avg, surf_fit_err_max);
767  // Get array with surface fitting weights.
768  Array<double> weights;
769  GetSurfaceFittingWeight(weights);
770 
772  {
773  mfem::out << "Avg/Max surface fitting error: " <<
774  surf_fit_err_avg << " " <<
775  surf_fit_err_max << "\n";
776  mfem::out << "Min/Max surface fitting weight: " <<
777  weights.Min() << " " << weights.Max() << "\n";
778  }
779 
780  double change_surf_fit_err = surf_fit_err_avg_prvs-surf_fit_err_avg;
781  double rel_change_surf_fit_err = change_surf_fit_err/surf_fit_err_avg_prvs;
782  // Increase the surface fitting coefficient if the surface fitting error
783  // does not decrease sufficiently.
784  if (rel_change_surf_fit_err < 1.e-2)
785  {
787  }
788  surf_fit_err_avg_prvs = surf_fit_err_avg;
789  update_surf_fit_coeff = false;
790  }
791 }
792 
794  const Vector &x_new) const
795 {
796  const bool update_flag = true;
797  DiscreteAdaptTC *discrtc = ti.GetDiscreteAdaptTC();
798  if (discrtc)
799  {
800  discrtc->UpdateTargetSpecification(x_new, update_flag);
801  if (ti.GetFDFlag())
802  {
803  double dx = ti.GetFDh();
804  discrtc->UpdateGradientTargetSpecification(x_new, dx, update_flag);
805  discrtc->UpdateHessianTargetSpecification(x_new, dx, update_flag);
806  }
807  }
808 }
809 
811  const FiniteElementSpace &fes) const
812 {
813  double min_detJ = infinity();
814  const int NE = fes.GetNE(), dim = fes.GetMesh()->Dimension();
815  Array<int> xdofs;
816  DenseMatrix Jpr(dim);
817  const bool mixed_mesh = fes.GetMesh()->GetNumGeometries(dim) > 1;
818  if (dim == 1 || mixed_mesh || UsesTensorBasis(fes) == false)
819  {
820  for (int i = 0; i < NE; i++)
821  {
822  const int dof = fes.GetFE(i)->GetDof();
823  DenseMatrix dshape(dof, dim), pos(dof, dim);
824  Vector posV(pos.Data(), dof * dim);
825 
826  fes.GetElementVDofs(i, xdofs);
827  x_loc.GetSubVector(xdofs, posV);
828 
829  const IntegrationRule &irule = GetIntegrationRule(*fes.GetFE(i));
830  const int nsp = irule.GetNPoints();
831  for (int j = 0; j < nsp; j++)
832  {
833  fes.GetFE(i)->CalcDShape(irule.IntPoint(j), dshape);
834  MultAtB(pos, dshape, Jpr);
835  min_detJ = std::min(min_detJ, Jpr.Det());
836  }
837  }
838  }
839  else
840  {
841  min_detJ = dim == 2 ? MinDetJpr_2D(&fes, x_loc) :
842  dim == 3 ? MinDetJpr_3D(&fes, x_loc) : 0.0;
843  }
844  double min_detT_all = min_detJ;
845 #ifdef MFEM_USE_MPI
846  if (parallel)
847  {
848  auto p_nlf = dynamic_cast<const ParNonlinearForm *>(oper);
849  MPI_Allreduce(&min_detJ, &min_detT_all, 1, MPI_DOUBLE, MPI_MIN,
850  p_nlf->ParFESpace()->GetComm());
851  }
852 #endif
853  const DenseMatrix &Wideal =
855  min_detT_all /= Wideal.Det();
856 
857  return min_detT_all;
858 }
859 
860 #ifdef MFEM_USE_MPI
861 // Metric values are visualized by creating an L2 finite element functions and
862 // computing the metric values at the nodes.
864  const TargetConstructor &tc, ParMesh &pmesh,
865  char *title, int position)
866 {
868  ParFiniteElementSpace fes(&pmesh, &fec, 1);
869  ParGridFunction metric(&fes);
870  InterpolateTMOP_QualityMetric(qm, tc, pmesh, metric);
871  socketstream sock;
872  if (pmesh.GetMyRank() == 0)
873  {
874  sock.open("localhost", 19916);
875  sock << "solution\n";
876  }
877  pmesh.PrintAsOne(sock);
878  metric.SaveAsOne(sock);
879  if (pmesh.GetMyRank() == 0)
880  {
881  sock << "window_title '"<< title << "'\n"
882  << "window_geometry "
883  << position << " " << 0 << " " << 600 << " " << 600 << "\n"
884  << "keys jRmclA\n";
885  }
886 }
887 #endif
888 
889 // Metric values are visualized by creating an L2 finite element functions and
890 // computing the metric values at the nodes.
892  const TargetConstructor &tc, Mesh &mesh,
893  char *title, int position)
894 {
896  FiniteElementSpace fes(&mesh, &fec, 1);
897  GridFunction metric(&fes);
898  InterpolateTMOP_QualityMetric(qm, tc, mesh, metric);
899  osockstream sock(19916, "localhost");
900  sock << "solution\n";
901  mesh.Print(sock);
902  metric.Save(sock);
903  sock.send();
904  sock << "window_title '"<< title << "'\n"
905  << "window_geometry "
906  << position << " " << 0 << " " << 600 << " " << 600 << "\n"
907  << "keys jRmclA\n";
908 }
909 
910 }
int GetNPoints() const
Returns the number of the points in the integration rule.
Definition: intrules.hpp:247
VectorGridFunctionCoefficient u_coeff
Definition: tmop_tools.hpp:85
void vis_tmop_metric_p(int order, TMOP_QualityMetric &qm, const TargetConstructor &tc, ParMesh &pmesh, char *title, int position)
Definition: tmop_tools.cpp:863
virtual void Interpolate(const GridFunction &field_in, Vector &field_out)
Definition: gslib.cpp:584
virtual const Operator * GetProlongation() const
Get the finite element space prolongation matrix.
virtual double GetTime() const
Read the currently set time.
Definition: operator.hpp:326
int Size() const
Return the logical size of the array.
Definition: array.hpp:138
int GetVSize() const
Return the number of vector dofs, i.e. GetNDofs() x GetVDim().
Definition: fespace.hpp:563
Conjugate gradient method.
Definition: solvers.hpp:465
int GetDim() const
Returns the reference space dimension for the finite element.
Definition: fe_base.hpp:311
ParFiniteElementSpace * ParFESpace() const
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
virtual void Mult(const Vector &x, Vector &y) const
Matrix vector multiplication: .
Class for grid function - Vector with associated FE space.
Definition: gridfunc.hpp:30
const SparseMatrix * GetConformingProlongation() const
The returned SparseMatrix is owned by the FiniteElementSpace.
Definition: fespace.cpp:1234
const Operator * oper
Definition: solvers.hpp:120
Data type for scaled Jacobi-type smoother of sparse matrix.
double GetSurfaceFittingWeight()
Get the surface fitting weight.
Definition: tmop.cpp:3469
ParMesh * GetParMesh() const
Definition: pfespace.hpp:277
void UpdateDiscreteTC(const TMOP_Integrator &ti, const Vector &x_new) const
Definition: tmop_tools.cpp:793
void Assemble(int skip_zeros=1)
Assembles the form i.e. sums over all domain/bdr integrators.
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
Definition: solvers.cpp:711
FiniteElementSpace * FESpace()
Return the FE space associated with the BilinearForm.
virtual const Operator * GetProlongationMatrix() const
The returned Operator is owned by the FiniteElementSpace.
Definition: pfespace.cpp:1153
int GetNumGeometries(int dim) const
Return the number of geometries of the given dimension present in the mesh.
Definition: mesh.cpp:5877
void SetSize(int s)
Resize the vector to size s.
Definition: vector.hpp:521
void InterpolateTMOP_QualityMetric(TMOP_QualityMetric &metric, const TargetConstructor &tc, const Mesh &mesh, GridFunction &metric_gf)
Interpolates the metric&#39;s values at the nodes of metric_gf.
Definition: tmop.cpp:3852
double MinDetJpr_2D(const FiniteElementSpace *, const Vector &) const
double Det() const
Definition: densemat.cpp:436
Base abstract class for first order time dependent operators.
Definition: operator.hpp:285
bool summary
A summary of the solver process will be reported after the last iteration to mfem::out.
Definition: solvers.hpp:91
bool first_and_last
Information about the first and last iteration will be printed to mfem::out.
Definition: solvers.hpp:94
Pointer to an Operator of a specified type.
Definition: handle.hpp:33
void GetSubVector(const Array< int > &dofs, Vector &elemvect) const
Extract entries listed in dofs to the output Vector elemvect.
Definition: vector.cpp:534
ParFiniteElementSpace * pfes
Definition: tmop.hpp:863
Data type dense matrix using column-major storage.
Definition: densemat.hpp:23
void GetSurfaceFittingErrors(double &err_avg, double &err_max)
Definition: tmop.cpp:2436
int Size() const
Returns the size of the vector.
Definition: vector.hpp:199
Parallel non-linear operator on the true dofs.
void Init(TimeDependentOperator &f_) override
Associate a TimeDependentOperator with the ODE solver.
Definition: ode.cpp:100
const IntegrationRule & GetIntegrationRule(const FiniteElement &el) const
Definition: tmop_tools.hpp:148
int GetNE() const
Returns number of elements.
Definition: mesh.hpp:928
bool IsSurfaceFittingEnabled()
Definition: tmop.hpp:1633
virtual void Mult(const Vector &x, Vector &y) const =0
Operator application: y=A(x).
void SetAssemblyLevel(AssemblyLevel assembly_level)
Set the desired assembly level.
Abstract parallel finite element space.
Definition: pfespace.hpp:28
void vis_tmop_metric_s(int order, TMOP_QualityMetric &qm, const TargetConstructor &tc, Mesh &mesh, char *title, int position)
Definition: tmop_tools.cpp:891
void UpdateSurfaceFittingWeight(double factor)
Update the surface fitting weight as surf_fit_coeff *= factor;.
Definition: tmop.cpp:3457
void UpdateGradientTargetSpecification(const Vector &x, double dx, bool use_flag=false)
Definition: tmop.cpp:2179
bool UsesTensorBasis(const FiniteElementSpace &fes)
Definition: fespace.hpp:983
T Min() const
Find the minimal element in the array, using the comparison operator &lt; for class T.
Definition: array.cpp:85
double GetFDh() const
Definition: tmop.hpp:1708
virtual double GetEnergy(const Vector &x) const
Compute the enery corresponding to the state x.
void add(const Vector &v1, const Vector &v2, Vector &v)
Definition: vector.cpp:300
const DenseMatrix & GetGeomToPerfGeomJac(int GeomType) const
Definition: geom.hpp:98
double Norm(const Vector &x) const
Definition: solvers.hpp:171
HypreParMatrix * ParallelAssemble()
Returns the matrix assembled on the true dofs, i.e. P^t A P.
void UpdateTargetSpecification(const Vector &new_x, bool use_flag=false)
Definition: tmop.cpp:1534
virtual void ComputeAtNewPosition(const Vector &new_nodes, Vector &new_field)
Definition: tmop_tools.cpp:29
virtual void Save(std::ostream &out) const
Save the GridFunction to an output stream.
Definition: gridfunc.cpp:3619
double GetEnergy(const ParGridFunction &x) const
Compute the energy of a ParGridFunction.
Geometry::Type GetGeomType() const
Returns the Geometry::Type of the reference element.
Definition: fe_base.hpp:320
int GetNE() const
Returns number of elements in the mesh.
Definition: fespace.hpp:590
Abstract class for local mesh quality metrics in the target-matrix optimization paradigm (TMOP) by P...
Definition: tmop.hpp:23
Geometry Geometries
Definition: fe.cpp:49
IntegrationPoint & IntPoint(int i)
Returns a reference to the i-th integration point.
Definition: intrules.hpp:250
void UpdateAfterMeshPositionChange(const Vector &new_x)
Definition: tmop.cpp:3608
DofTransformation * GetElementVDofs(int i, Array< int > &vdofs) const
Returns indexes of degrees of freedom in array dofs for i&#39;th element.
Definition: fespace.cpp:281
virtual void SetPrintLevel(int print_lvl)
Legacy method to set the level of verbosity of the solver output.
Definition: solvers.cpp:71
MPI_Comm GetComm() const
Definition: pfespace.hpp:273
double f(const Vector &xvec)
Definition: lor_mms.hpp:32
double GetElementSize(ElementTransformation *T, int type=0)
Definition: mesh.cpp:76
void GetSurfaceFittingWeight(Array< double > &weights) const
Get the surface fitting weight for all the TMOP integrators.
Definition: tmop_tools.cpp:589
Data type sparse matrix.
Definition: sparsemat.hpp:46
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
Definition: operator.hpp:66
int Append(const T &el)
Append element &#39;el&#39; to array, resize if necessary.
Definition: array.hpp:751
virtual void FreeData()
Definition: gslib.cpp:214
Mesh * GetMesh() const
Returns the mesh.
Definition: fespace.hpp:433
ParFiniteElementSpace * ParFESpace() const
Return the parallel FE space associated with the ParBilinearForm.
Array< NonlinearFormIntegrator * > * GetDNFI()
Access all integrators added with AddDomainIntegrator().
Jacobi smoothing for a given bilinear form (no matrix necessary).
Definition: solvers.hpp:274
double b
Definition: lissajous.cpp:42
const AssemblyLevel al
Definition: tmop_tools.hpp:87
static MemoryType GetDeviceMemoryType()
Get the current Device MemoryType. This is the MemoryType used by most MFEM classes when allocating m...
Definition: device.hpp:273
int Size() const
Get the size of the BilinearForm as a square matrix.
virtual void Mult(const Vector &ind, Vector &di_dt) const
Perform the action of the operator: y = k = f(x, t), where k solves the algebraic equation F(x...
Definition: tmop_tools.cpp:188
void SetMaxIter(int max_it)
Definition: solvers.hpp:200
T Max() const
Find the maximal element in the array, using the comparison operator &lt; for class T.
Definition: array.cpp:68
virtual const Operator * GetProlongationMatrix() const
The returned Operator is owned by the FiniteElementSpace.
Definition: fespace.hpp:471
void Assemble(int skip_zeros=1)
Assemble the local matrix.
void Setup(Mesh &m, const double bb_t=0.1, const double newt_tol=1.0e-12, const int npt_max=256)
Definition: gslib.cpp:75
Performs a single remap advection step in serial.
Definition: tmop_tools.hpp:79
void UpdateSurfaceFittingWeight(double factor) const
Update surface fitting weight as surf_fit_weight *= factor.
Definition: tmop_tools.cpp:564
Parallel smoothers in hypre.
Definition: hypre.hpp:896
FiniteElementSpace * FESpace()
Definition: gridfunc.hpp:652
int Dimension() const
Definition: mesh.hpp:999
void SetPAMemoryType(MemoryType mt)
Definition: nonlininteg.hpp:50
double GetParGridFunctionEnergy(const Vector &x) const
Compute the energy corresponding to the state x.
FiniteElementSpace * FESpace()
Wrapper for hypre&#39;s parallel vector class.
Definition: hypre.hpp:148
double Min() const
Returns the minimal element of the vector.
Definition: vector.cpp:1182
virtual void ProcessNewState(const Vector &x) const
Definition: tmop_tools.cpp:659
double * Data() const
Returns the matrix data array.
Definition: densemat.hpp:111
void SetAbsTol(double atol)
Definition: solvers.hpp:199
PrintLevel print_options
Output behavior for the iterative solver.
Definition: solvers.hpp:137
int GetMyRank() const
Definition: pmesh.hpp:290
void SetRelTol(double rtol)
Definition: solvers.hpp:198
void UpdateHessianTargetSpecification(const Vector &x, double dx, bool use_flag=false)
Definition: tmop.cpp:2205
MPI_Comm GetComm() const
Definition: pmesh.hpp:288
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition: fespace.hpp:88
int GetDof() const
Returns the number of degrees of freedom in the finite element.
Definition: fe_base.hpp:323
virtual void Mult(const Vector &x, Vector &y) const
Matrix vector multiplication.
Definition: sparsemat.cpp:638
const AssemblyLevel al
Definition: tmop_tools.hpp:109
void subtract(const Vector &x, const Vector &y, Vector &z)
Definition: vector.cpp:447
MemoryType
Memory types supported by MFEM.
Definition: mem_manager.hpp:31
Performs a single remap advection step in parallel.
Definition: tmop_tools.hpp:101
void SetSize(int nsize)
Change the logical size of the array, keep existing entries.
Definition: array.hpp:679
DiscreteAdaptTC * GetDiscreteAdaptTC() const
Definition: tmop.hpp:1691
void ResetUpdateFlags()
Used in combination with the Update methods to avoid extra computations.
Definition: tmop.hpp:1189
void SaveAsOne(const char *fname, int precision=16) const
Definition: pgridfunc.cpp:899
virtual void ComputeAtNewPosition(const Vector &new_nodes, Vector &new_field)
Definition: tmop_tools.cpp:354
virtual void FormSystemMatrix(const Array< int > &ess_tdof_list, OperatorHandle &A)
Form the linear system matrix A, see FormLinearSystem() for details.
void ComputeFDh(const Vector &x, const FiniteElementSpace &fes)
Determines the perturbation, h, for FD-based approximation.
Definition: tmop.cpp:3624
double GetGridFunctionEnergy(const Vector &x) const
Compute the enery corresponding to the state x.
virtual void Finalize(int skip_zeros=1)
Finalizes the matrix initialization.
double Max() const
Returns the maximal element of the vector.
Definition: vector.cpp:882
double MinDetJpr_3D(const FiniteElementSpace *, const Vector &) const
virtual void Print(std::ostream &os=mfem::out) const
Definition: mesh.hpp:1646
FDualNumber< tbase > sqrt(const FDualNumber< tbase > &f)
sqrt([dual number])
Definition: fdual.hpp:600
Mesh * GetMesh(int type)
Definition: ex29.cpp:218
int dim
Definition: ex24.cpp:53
virtual void SetInitialField(const Vector &init_nodes, const Vector &init_field)
Definition: tmop_tools.cpp:22
VectorGridFunctionCoefficient u_coeff
Definition: tmop_tools.hpp:107
virtual void Mult(const Vector &x, Vector &y) const
Perform the action of the operator: y = k = f(x, t), where k solves the algebraic equation F(x...
Definition: operator.cpp:230
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
Adds new Domain Integrator. Assumes ownership of bfi.
virtual const FiniteElement * GetFE(int i) const
Returns pointer to the FiniteElement in the FiniteElementCollection associated with i&#39;th element in t...
Definition: fespace.cpp:2783
virtual void GetSurfaceFittingError(double &err_avg, double &err_max) const
Definition: tmop_tools.cpp:619
SerialAdvectorCGOper(const Vector &x_start, GridFunction &vel, FiniteElementSpace &fes, AssemblyLevel al=AssemblyLevel::LEGACY)
Definition: tmop_tools.cpp:167
double infinity()
Define a shortcut for std::numeric_limits&lt;double&gt;::infinity()
Definition: vector.hpp:46
virtual void SetSpace(FiniteElementSpace *f)
Associate a new FiniteElementSpace with the GridFunction.
Definition: gridfunc.cpp:193
int open(const char hostname[], int port)
Open the socket stream on &#39;port&#39; at &#39;hostname&#39;.
virtual void Mult(const Vector &ind, Vector &di_dt) const
Perform the action of the operator: y = k = f(x, t), where k solves the algebraic equation F(x...
Definition: tmop_tools.cpp:261
RefCoord t[3]
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
Definition: solvers.hpp:479
ParAdvectorCGOper(const Vector &x_start, GridFunction &vel, ParFiniteElementSpace &pfes, AssemblyLevel al=AssemblyLevel::LEGACY, MemoryType mt=MemoryType::DEFAULT)
Definition: tmop_tools.cpp:231
void MultAtB(const DenseMatrix &A, const DenseMatrix &B, DenseMatrix &AtB)
Multiply the transpose of a matrix A with a matrix B: At*B.
Definition: densemat.cpp:2633
FindPointsGSLIB can robustly evaluate a GridFunction on an arbitrary collection of points...
Definition: gslib.hpp:47
bool GetFDFlag() const
Definition: tmop.hpp:1707
Vector data type.
Definition: vector.hpp:60
const FiniteElementCollection * FEColl() const
Definition: fespace.hpp:577
void MakeRef(Vector &base, int offset, int size)
Reset the Vector to be a reference to a sub-vector of base.
Definition: vector.hpp:585
void GetNodes(Vector &node_coord) const
Definition: mesh.cpp:7841
virtual double ComputeScalingFactor(const Vector &x, const Vector &b) const
Definition: tmop_tools.cpp:362
void PrintAsOne(std::ostream &out=mfem::out) const
Definition: pmesh.cpp:4829
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
Definition: solvers.cpp:173
void Step(Vector &x, double &t, double &dt) override
Perform a time step from time t [in] to time t [out] based on the requested step size dt [in]...
Definition: ode.cpp:109
void DeleteGeometricFactors()
Destroy all GeometricFactors stored by the Mesh.
Definition: mesh.cpp:881
RefCoord s[3]
const SparseMatrix & SpMat() const
Returns a const reference to the sparse matrix.
bool iterations
Detailed information about each iteration will be reported to mfem::out.
Definition: solvers.hpp:88
double u(const Vector &xvec)
Definition: lor_mms.hpp:24
void SetNodes(const Vector &node_coord)
Definition: mesh.cpp:7853
Base class for solvers.
Definition: operator.hpp:651
Class for parallel grid function.
Definition: pgridfunc.hpp:32
Base class representing target-matrix construction algorithms for mesh optimization via the target-ma...
Definition: tmop.hpp:908
OutStream out(std::cout)
Global stream used by the library for standard output. Initially it uses the same std::streambuf as s...
Definition: globals.hpp:66
Abstract operator.
Definition: operator.hpp:24
virtual void CalcDShape(const IntegrationPoint &ip, DenseMatrix &dshape) const =0
Evaluate the gradients of all shape functions of a scalar finite element in reference space at the gi...
Class for parallel meshes.
Definition: pmesh.hpp:32
double ComputeMinDet(const Vector &x_loc, const FiniteElementSpace &fes) const
Definition: tmop_tools.cpp:810
const Array< TMOP_Integrator * > & GetTMOPIntegrators() const
Definition: tmop.hpp:1737
virtual double * HostReadWrite()
Shortcut for mfem::ReadWrite(vec.GetMemory(), vec.Size(), false).
Definition: vector.hpp:458
void Reset(OpType *A, bool own_A=true)
Reset the OperatorHandle to the given OpType pointer, A.
Definition: handle.hpp:145
virtual void SetInitialField(const Vector &init_nodes, const Vector &init_field)
Definition: tmop_tools.cpp:315
FiniteElementSpace * fes
Definition: tmop.hpp:858
Arbitrary order &quot;L2-conforming&quot; discontinuous finite elements.
Definition: fe_coll.hpp:284
alpha (q . grad u, v)
A TMOP integrator class based on any given TMOP_QualityMetric and TargetConstructor.
Definition: tmop.hpp:1303