MFEM  v4.3.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-2021, Lawrence Livermore National Security, LLC. Produced
2 // at the Lawrence Livermore National Laboratory. All Rights reserved. See files
3 // LICENSE and NOTICE for details. LLNL-CODE-806117.
4 //
5 // This file is part of the MFEM library. For more information and source code
6 // availability visit https://mfem.org.
7 //
8 // MFEM is free software; you can redistribute it and/or modify it under the
9 // terms of the BSD-3 license. We welcome feedback and contributions, see file
10 // CONTRIBUTING.md for details.
11 
12 #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  // Check if the starting mesh (given by x) is inverted. Note that x hasn't
402  // been modified by the Newton update yet.
403  const double min_detT_in = ComputeMinDet(x_out_loc, *fes);
404  const bool untangling = (min_detT_in <= 0.0) ? true : false;
405  const double untangle_factor = 1.5;
406  if (untangling)
407  {
408  // Needed for the line search below. The untangling metrics see this
409  // reference to detect deteriorations.
410  *min_det_ptr = untangle_factor * min_detT_in;
411  }
412 
413  const bool have_b = (b.Size() == Height());
414 
415  Vector x_out(x.Size());
416  bool x_out_ok = false;
417  double scale = 1.0, energy_out = 0.0, min_detT_out;
418  const double norm_in = Norm(r);
419 
420  const double detJ_factor = (solver_type == 1) ? 0.25 : 0.5;
421 
422  // Perform the line search.
423  for (int i = 0; i < 12; i++)
424  {
425  // Update the mesh and get the L-vector in x_out_loc.
426  add(x, -scale, c, x_out);
427  if (serial)
428  {
429  const SparseMatrix *cP = fes->GetConformingProlongation();
430  if (!cP) { x_out_loc = x_out; }
431  else { cP->Mult(x_out, x_out_loc); }
432  }
433 #ifdef MFEM_USE_MPI
434  else
435  {
436  fes->GetProlongationMatrix()->Mult(x_out, x_out_loc);
437  }
438 #endif
439 
440  // Check the changes in detJ.
441  min_detT_out = ComputeMinDet(x_out_loc, *fes);
442  if (untangling == false && min_detT_out < 0.0)
443  {
444  // No untangling, and detJ got negative -- no good.
445  if (print_level >= 0)
446  { mfem::out << "Scale = " << scale << " Neg det(J) found.\n"; }
447  scale *= detJ_factor; continue;
448  }
449  if (untangling == true && min_detT_out < *min_det_ptr)
450  {
451  // Untangling, and detJ got even more negative -- no good.
452  if (print_level >= 0)
453  { mfem::out << "Scale = " << scale << " Neg det(J) decreased.\n"; }
454  scale *= detJ_factor; continue;
455  }
456 
457  // Skip the energy and residual checks when we're untangling. The
458  // untangling metrics change their denominators, which can affect the
459  // energy and residual, so their increase/decrease is not relevant.
460  if (untangling) { x_out_ok = true; break; }
461 
462  // Check the changes in total energy.
463  ProcessNewState(x_out);
464 
465  if (serial)
466  {
467  energy_out = nlf->GetGridFunctionEnergy(x_out_loc);
468  }
469 #ifdef MFEM_USE_MPI
470  else
471  {
472  energy_out = p_nlf->GetParGridFunctionEnergy(x_out_loc);
473  }
474 #endif
475  if (energy_out > energy_in + 0.2*fabs(energy_in) ||
476  std::isnan(energy_out) != 0)
477  {
478  if (print_level >= 0)
479  {
480  mfem::out << "Scale = " << scale << " Increasing energy: "
481  << energy_in << " --> " << energy_out << '\n';
482  }
483  scale *= 0.5; continue;
484  }
485 
486  // Check the changes in the Newton residual.
487  oper->Mult(x_out, r);
488  if (have_b) { r -= b; }
489  double norm_out = Norm(r);
490 
491  if (norm_out > 1.2*norm_in)
492  {
493  if (print_level >= 0)
494  {
495  mfem::out << "Scale = " << scale << " Norm increased: "
496  << norm_in << " --> " << norm_out << '\n';
497  }
498  scale *= 0.5; continue;
499  }
500  else { x_out_ok = true; break; }
501  } // end line search
502 
503  if (untangling)
504  {
505  // Update the global min detJ. Untangling metrics see this min_det_ptr.
506  if (min_detT_out > 0.0)
507  {
508  *min_det_ptr = 0.0;
509  if (print_level >= 0)
510  { mfem::out << "The mesh has been untangled at the used points!\n"; }
511  }
512  else { *min_det_ptr = untangle_factor * min_detT_out; }
513  }
514 
515  if (print_level >= 0)
516  {
517  if (untangling)
518  {
519  mfem::out << "Min det(T) change: "
520  << min_detT_in << " -> " << min_detT_out
521  << " with " << scale << " scaling.\n";
522  }
523  else
524  {
525  mfem::out << "Energy decrease: "
526  << energy_in << " --> " << energy_out << " or "
527  << (energy_in - energy_out) / energy_in * 100.0
528  << "% with " << scale << " scaling.\n";
529  }
530  }
531 
532  if (x_out_ok == false) { scale = 0.0; }
533  return scale;
534 }
535 
537 {
538  const NonlinearForm *nlf = dynamic_cast<const NonlinearForm *>(oper);
539  const Array<NonlinearFormIntegrator*> &integs = *nlf->GetDNFI();
540 
541  // Reset the update flags of all TargetConstructors. This is done to avoid
542  // repeated updates of shared TargetConstructors.
543  TMOP_Integrator *ti = NULL;
544  TMOPComboIntegrator *co = NULL;
545  DiscreteAdaptTC *dtc = NULL;
546  for (int i = 0; i < integs.Size(); i++)
547  {
548  ti = dynamic_cast<TMOP_Integrator *>(integs[i]);
549  if (ti)
550  {
551  dtc = ti->GetDiscreteAdaptTC();
552  if (dtc) { dtc->ResetUpdateFlags(); }
553  }
554  co = dynamic_cast<TMOPComboIntegrator *>(integs[i]);
555  if (co)
556  {
558  for (int j = 0; j < ati.Size(); j++)
559  {
560  dtc = ati[j]->GetDiscreteAdaptTC();
561  if (dtc) { dtc->ResetUpdateFlags(); }
562  }
563  }
564  }
565 
566  if (parallel)
567  {
568 #ifdef MFEM_USE_MPI
569  const ParNonlinearForm *nlf =
570  dynamic_cast<const ParNonlinearForm *>(oper);
571  const ParFiniteElementSpace *pfesc = nlf->ParFESpace();
572  Vector x_loc(pfesc->GetVSize());
573  pfesc->GetProlongationMatrix()->Mult(x, x_loc);
574  for (int i = 0; i < integs.Size(); i++)
575  {
576  ti = dynamic_cast<TMOP_Integrator *>(integs[i]);
577  if (ti)
578  {
579  ti->UpdateAfterMeshChange(x_loc);
580  ti->ComputeFDh(x_loc, *pfesc);
581  UpdateDiscreteTC(*ti, x_loc);
582  }
583  co = dynamic_cast<TMOPComboIntegrator *>(integs[i]);
584  if (co)
585  {
587  for (int j = 0; j < ati.Size(); j++)
588  {
589  ati[j]->UpdateAfterMeshChange(x_loc);
590  ati[j]->ComputeFDh(x_loc, *pfesc);
591  UpdateDiscreteTC(*ati[j], x_loc);
592  }
593  }
594  }
595 #endif
596  }
597  else
598  {
599  const FiniteElementSpace *fesc = nlf->FESpace();
600  const Operator *P = nlf->GetProlongation();
601  Vector x_loc;
602  if (P)
603  {
604  x_loc.SetSize(P->Height());
605  P->Mult(x,x_loc);
606  }
607  else
608  {
609  x_loc = x;
610  }
611  for (int i = 0; i < integs.Size(); i++)
612  {
613  ti = dynamic_cast<TMOP_Integrator *>(integs[i]);
614  if (ti)
615  {
616  ti->UpdateAfterMeshChange(x_loc);
617  ti->ComputeFDh(x_loc, *fesc);
618  UpdateDiscreteTC(*ti, x_loc);
619  }
620  co = dynamic_cast<TMOPComboIntegrator *>(integs[i]);
621  if (co)
622  {
624  for (int j = 0; j < ati.Size(); j++)
625  {
626  ati[j]->UpdateAfterMeshChange(x_loc);
627  ati[j]->ComputeFDh(x_loc, *fesc);
628  UpdateDiscreteTC(*ati[j], x_loc);
629  }
630  }
631  }
632  }
633 }
634 
636  const Vector &x_new) const
637 {
638  const bool update_flag = true;
639  DiscreteAdaptTC *discrtc = ti.GetDiscreteAdaptTC();
640  if (discrtc)
641  {
642  discrtc->UpdateTargetSpecification(x_new, update_flag);
643  if (ti.GetFDFlag())
644  {
645  double dx = ti.GetFDh();
646  discrtc->UpdateGradientTargetSpecification(x_new, dx, update_flag);
647  discrtc->UpdateHessianTargetSpecification(x_new, dx, update_flag);
648  }
649  }
650 }
651 
653  const FiniteElementSpace &fes) const
654 {
655  double min_detJ = infinity();
656  const int NE = fes.GetNE(), dim = fes.GetMesh()->Dimension();
657  Array<int> xdofs;
658  DenseMatrix Jpr(dim);
659  const bool mixed_mesh = fes.GetMesh()->GetNumGeometries(dim) > 1;
660  if (dim == 1 || mixed_mesh || UsesTensorBasis(fes) == false)
661  {
662  for (int i = 0; i < NE; i++)
663  {
664  const int dof = fes.GetFE(i)->GetDof();
665  DenseMatrix dshape(dof, dim), pos(dof, dim);
666  Vector posV(pos.Data(), dof * dim);
667 
668  fes.GetElementVDofs(i, xdofs);
669  x_loc.GetSubVector(xdofs, posV);
670 
671  const IntegrationRule &irule = GetIntegrationRule(*fes.GetFE(i));
672  const int nsp = irule.GetNPoints();
673  for (int j = 0; j < nsp; j++)
674  {
675  fes.GetFE(i)->CalcDShape(irule.IntPoint(j), dshape);
676  MultAtB(pos, dshape, Jpr);
677  min_detJ = std::min(min_detJ, Jpr.Det());
678  }
679  }
680  }
681  else
682  {
683  min_detJ = dim == 2 ? MinDetJpr_2D(&fes, x_loc) :
684  dim == 3 ? MinDetJpr_3D(&fes, x_loc) : 0.0;
685  }
686  double min_detT_all = min_detJ;
687 #ifdef MFEM_USE_MPI
688  if (parallel)
689  {
690  auto p_nlf = dynamic_cast<const ParNonlinearForm *>(oper);
691  MPI_Allreduce(&min_detJ, &min_detT_all, 1, MPI_DOUBLE, MPI_MIN,
692  p_nlf->ParFESpace()->GetComm());
693  }
694 #endif
695  const DenseMatrix &Wideal =
697  min_detT_all /= Wideal.Det();
698 
699  return min_detT_all;
700 }
701 
702 #ifdef MFEM_USE_MPI
703 // Metric values are visualized by creating an L2 finite element functions and
704 // computing the metric values at the nodes.
706  const TargetConstructor &tc, ParMesh &pmesh,
707  char *title, int position)
708 {
710  ParFiniteElementSpace fes(&pmesh, &fec, 1);
711  ParGridFunction metric(&fes);
712  InterpolateTMOP_QualityMetric(qm, tc, pmesh, metric);
713  socketstream sock;
714  if (pmesh.GetMyRank() == 0)
715  {
716  sock.open("localhost", 19916);
717  sock << "solution\n";
718  }
719  pmesh.PrintAsOne(sock);
720  metric.SaveAsOne(sock);
721  if (pmesh.GetMyRank() == 0)
722  {
723  sock << "window_title '"<< title << "'\n"
724  << "window_geometry "
725  << position << " " << 0 << " " << 600 << " " << 600 << "\n"
726  << "keys jRmclA\n";
727  }
728 }
729 #endif
730 
731 // Metric values are visualized by creating an L2 finite element functions and
732 // computing the metric values at the nodes.
734  const TargetConstructor &tc, Mesh &mesh,
735  char *title, int position)
736 {
738  FiniteElementSpace fes(&mesh, &fec, 1);
739  GridFunction metric(&fes);
740  InterpolateTMOP_QualityMetric(qm, tc, mesh, metric);
741  osockstream sock(19916, "localhost");
742  sock << "solution\n";
743  mesh.Print(sock);
744  metric.Save(sock);
745  sock.send();
746  sock << "window_title '"<< title << "'\n"
747  << "window_geometry "
748  << position << " " << 0 << " " << 600 << " " << 600 << "\n"
749  << "keys jRmclA\n";
750 }
751 
752 }
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:705
virtual void Interpolate(const GridFunction &field_in, Vector &field_out)
Definition: gslib.cpp:582
virtual const Operator * GetProlongation() const
Get the finite element space prolongation matrix.
virtual double GetTime() const
Read the currently set time.
Definition: operator.hpp:323
int Size() const
Return the logical size of the array.
Definition: array.hpp:134
int GetVSize() const
Return the number of vector dofs, i.e. GetNDofs() x GetVDim().
Definition: fespace.hpp:540
virtual void Print(std::ostream &out=mfem::out) const
Definition: mesh.hpp:1387
Conjugate gradient method.
Definition: solvers.hpp:316
int GetDim() const
Returns the reference space dimension for the finite element.
Definition: fe.hpp:317
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
void UpdateAfterMeshChange(const Vector &new_x)
Definition: tmop.cpp:3042
const SparseMatrix * GetConformingProlongation() const
The returned SparseMatrix is owned by the FiniteElementSpace.
Definition: fespace.cpp:1168
const Operator * oper
Definition: solvers.hpp:75
Data type for scaled Jacobi-type smoother of sparse matrix.
ParMesh * GetParMesh() const
Definition: pfespace.hpp:267
void UpdateDiscreteTC(const TMOP_Integrator &ti, const Vector &x_new) const
Definition: tmop_tools.cpp:635
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:616
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:910
int GetNumGeometries(int dim) const
Return the number of geometries of the given dimension present in the mesh.
Definition: mesh.cpp:5428
void SetSize(int s)
Resize the vector to size s.
Definition: vector.hpp:513
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:3254
void GetElementVDofs(int i, Array< int > &vdofs) const
Returns indexes of degrees of freedom in array dofs for i&#39;th element.
Definition: fespace.cpp:262
double MinDetJpr_2D(const FiniteElementSpace *, const Vector &) const
double Det() const
Definition: densemat.cpp:451
Base abstract class for first order time dependent operators.
Definition: operator.hpp:282
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:525
ParFiniteElementSpace * pfes
Definition: tmop.hpp:823
Data type dense matrix using column-major storage.
Definition: densemat.hpp:23
int Size() const
Returns the size of the vector.
Definition: vector.hpp:190
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:142
int GetNE() const
Returns number of elements.
Definition: mesh.hpp:846
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:733
void UpdateGradientTargetSpecification(const Vector &x, double dx, bool use_flag=false)
Definition: tmop.cpp:2068
bool UsesTensorBasis(const FiniteElementSpace &fes)
Definition: fespace.hpp:957
double GetFDh() const
Definition: tmop.hpp:1520
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:291
const DenseMatrix & GetGeomToPerfGeomJac(int GeomType) const
Definition: geom.hpp:97
double Norm(const Vector &x) const
Definition: solvers.hpp:87
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:1486
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:3484
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.hpp:320
int GetNE() const
Returns number of elements in the mesh.
Definition: fespace.hpp:567
Abstract class for local mesh quality metrics in the target-matrix optimization paradigm (TMOP) by P...
Definition: tmop.hpp:23
Geometry Geometries
Definition: fe.cpp:13507
IntegrationPoint & IntPoint(int i)
Returns a reference to the i-th integration point.
Definition: intrules.hpp:250
void SetPrintLevel(int print_lvl)
Definition: solvers.cpp:71
MPI_Comm GetComm() const
Definition: pfespace.hpp:263
double f(const Vector &xvec)
Definition: lor_mms.hpp:32
double GetElementSize(ElementTransformation *T, int type=0)
Definition: mesh.cpp:76
Data type sparse matrix.
Definition: sparsemat.hpp:41
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
Definition: operator.hpp:66
virtual void FreeData()
Definition: gslib.cpp:212
Mesh * GetMesh() const
Returns the mesh.
Definition: fespace.hpp:410
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:130
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:100
virtual const Operator * GetProlongationMatrix() const
The returned Operator is owned by the FiniteElementSpace.
Definition: fespace.hpp:448
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
Parallel smoothers in hypre.
Definition: hypre.hpp:840
FiniteElementSpace * FESpace()
Definition: gridfunc.hpp:629
int Dimension() const
Definition: mesh.hpp:911
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:99
double Min() const
Returns the minimal element of the vector.
Definition: vector.cpp:1173
virtual void ProcessNewState(const Vector &x) const
This method can be overloaded in derived classes to perform computations that need knowledge of the n...
Definition: tmop_tools.cpp:536
double * Data() const
Returns the matrix data array.
Definition: densemat.hpp:111
void SetAbsTol(double atol)
Definition: solvers.hpp:99
int GetMyRank() const
Definition: pmesh.hpp:278
void SetRelTol(double rtol)
Definition: solvers.hpp:98
void UpdateHessianTargetSpecification(const Vector &x, double dx, bool use_flag=false)
Definition: tmop.cpp:2094
MPI_Comm GetComm() const
Definition: pmesh.hpp:276
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:323
virtual void Mult(const Vector &x, Vector &y) const
Matrix vector multiplication.
Definition: sparsemat.cpp:613
const AssemblyLevel al
Definition: tmop_tools.hpp:109
void subtract(const Vector &x, const Vector &y, Vector &z)
Definition: vector.cpp:438
MemoryType
Memory types supported by MFEM.
Definition: mem_manager.hpp:31
Performs a single remap advection step in parallel.
Definition: tmop_tools.hpp:101
DiscreteAdaptTC * GetDiscreteAdaptTC() const
Definition: tmop.hpp:1503
void ResetUpdateFlags()
Used in combination with the Update methods to avoid extra computations.
Definition: tmop.hpp:1128
void SaveAsOne(const char *fname, int precision=16) const
Definition: pgridfunc.cpp:867
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:3052
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:873
double MinDetJpr_3D(const FiniteElementSpace *, const Vector &) const
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:225
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:2388
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:330
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:2648
FindPointsGSLIB can robustly evaluate a GridFunction on an arbitrary collection of points...
Definition: gslib.hpp:47
bool GetFDFlag() const
Definition: tmop.hpp:1519
Vector data type.
Definition: vector.hpp:60
const FiniteElementCollection * FEColl() const
Definition: fespace.hpp:554
void MakeRef(Vector &base, int offset, int size)
Reset the Vector to be a reference to a sub-vector of base.
Definition: vector.hpp:577
void GetNodes(Vector &node_coord) const
Definition: mesh.cpp:7343
virtual double ComputeScalingFactor(const Vector &x, const Vector &b) const
This method can be overloaded in derived classes to implement line search algorithms.
Definition: tmop_tools.cpp:362
void PrintAsOne(std::ostream &out=mfem::out) const
Definition: pmesh.cpp:4527
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
Definition: solvers.cpp:92
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:825
RefCoord s[3]
const SparseMatrix & SpMat() const
Returns a const reference to the sparse matrix.
double u(const Vector &xvec)
Definition: lor_mms.hpp:24
void SetNodes(const Vector &node_coord)
Definition: mesh.cpp:7355
Base class for solvers.
Definition: operator.hpp:648
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:868
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:652
const Array< TMOP_Integrator * > & GetTMOPIntegrators() const
Definition: tmop.hpp:1543
virtual double * HostReadWrite()
Shortcut for mfem::ReadWrite(vec.GetMemory(), vec.Size(), false).
Definition: vector.hpp:446
void Reset(OpType *A, bool own_A=true)
Reset the OperatorHandle to the given OpType pointer, A.
Definition: handle.hpp:140
virtual void SetInitialField(const Vector &init_nodes, const Vector &init_field)
Definition: tmop_tools.cpp:315
FiniteElementSpace * fes
Definition: tmop.hpp:818
Arbitrary order &quot;L2-conforming&quot; discontinuous finite elements.
Definition: fe_coll.hpp:285
alpha (q . grad u, v)
A TMOP integrator class based on any given TMOP_QualityMetric and TargetConstructor.
Definition: tmop.hpp:1198