MFEM  v4.2.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-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 #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 
37  for (int i = 0; i < ncomp; i++)
38  {
39  Vector new_field_temp(new_field.GetData()+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);
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);
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  for (int i = 0; i < s; i++)
98  {
99  double vel = 0.;
100  for (int j = 0; j < dim; j++)
101  {
102  vel += u(i+j*s)*u(i+j*s);
103  }
104  v_max = std::max(v_max, vel);
105  }
106 
107 #ifdef MFEM_USE_MPI
108  if (pfes)
109  {
110  double v_loc = v_max, h_loc = h_min;
111  MPI_Allreduce(&v_loc, &v_max, 1, MPI_DOUBLE, MPI_MAX, pfes->GetComm());
112  MPI_Allreduce(&h_loc, &h_min, 1, MPI_DOUBLE, MPI_MIN, pfes->GetComm());
113  }
114 #endif
115 
116  if (v_max == 0.0) // No need to change the field.
117  {
118  delete oper;
119  delete fess;
120 #ifdef MFEM_USE_MPI
121  delete pfess;
122 #endif
123  return;
124  }
125 
126  v_max = std::sqrt(v_max);
127  double dt = dt_scale * h_min / v_max;
128 
129  double t = 0.0;
130  bool last_step = false;
131  for (int ti = 1; !last_step; ti++)
132  {
133  if (t + dt >= 1.0)
134  {
135  dt = 1.0 - t;
136  last_step = true;
137  }
138  ode_solver.Step(new_field, t, dt);
139  }
140 
141  double glob_minv = minv,
142  glob_maxv = maxv;
143 #ifdef MFEM_USE_MPI
144  if (pfes)
145  {
146  MPI_Allreduce(&minv, &glob_minv, 1, MPI_DOUBLE, MPI_MIN, pfes->GetComm());
147  MPI_Allreduce(&maxv, &glob_maxv, 1, MPI_DOUBLE, MPI_MAX, pfes->GetComm());
148  }
149 #endif
150 
151  // Trim the overshoots and undershoots.
152  for (int i = 0; i < s; i++)
153  {
154  if (new_field(i) < glob_minv) { new_field(i) = glob_minv; }
155  if (new_field(i) > glob_maxv) { new_field(i) = glob_maxv; }
156  }
157 
158  delete oper;
159  delete fess;
160 #ifdef MFEM_USE_MPI
161  delete pfess;
162 #endif
163 }
164 
166  GridFunction &vel,
167  FiniteElementSpace &fes)
168  : TimeDependentOperator(fes.GetVSize()),
169  x0(x_start), x_now(*fes.GetMesh()->GetNodes()),
170  u(vel), u_coeff(&u), M(&fes), K(&fes)
171 {
173  K.AddDomainIntegrator(Kinteg);
174  K.Assemble(0);
175  K.Finalize(0);
176 
177  MassIntegrator *Minteg = new MassIntegrator;
178  M.AddDomainIntegrator(Minteg);
179  M.Assemble(0);
180  M.Finalize(0);
181 }
182 
183 void SerialAdvectorCGOper::Mult(const Vector &ind, Vector &di_dt) const
184 {
185  // Move the mesh.
186  const double t = GetTime();
187  add(x0, t, u, x_now);
188 
189  // Assemble on the new mesh.
190  K.BilinearForm::operator=(0.0);
191  K.Assemble();
192  Vector rhs(K.Size());
193  K.Mult(ind, rhs);
194  M.BilinearForm::operator=(0.0);
195  M.Assemble();
196 
197  di_dt = 0.0;
198  CGSolver lin_solver;
199  DSmoother prec;
200  lin_solver.SetPreconditioner(prec);
201  lin_solver.SetOperator(M.SpMat());
202  lin_solver.SetRelTol(1e-12); lin_solver.SetAbsTol(0.0);
203  lin_solver.SetMaxIter(100);
204  lin_solver.SetPrintLevel(0);
205  lin_solver.Mult(rhs, di_dt);
206 }
207 
208 #ifdef MFEM_USE_MPI
210  GridFunction &vel,
211  ParFiniteElementSpace &pfes)
212  : TimeDependentOperator(pfes.GetVSize()),
213  x0(x_start), x_now(*pfes.GetMesh()->GetNodes()),
214  u(vel), u_coeff(&u), M(&pfes), K(&pfes)
215 {
217  K.AddDomainIntegrator(Kinteg);
218  K.Assemble(0);
219  K.Finalize(0);
220 
221  MassIntegrator *Minteg = new MassIntegrator;
222  M.AddDomainIntegrator(Minteg);
223  M.Assemble(0);
224  M.Finalize(0);
225 }
226 
227 void ParAdvectorCGOper::Mult(const Vector &ind, Vector &di_dt) const
228 {
229  // Move the mesh.
230  const double t = GetTime();
231  add(x0, t, u, x_now);
232 
233  // Assemble on the new mesh.
234  K.BilinearForm::operator=(0.0);
235  K.Assemble();
237  K.Mult(ind, rhs);
238  M.BilinearForm::operator=(0.0);
239  M.Assemble();
240 
241  HypreParVector *RHS = rhs.ParallelAssemble();
243  X = 0.0;
245 
246  CGSolver lin_solver(M.ParFESpace()->GetParMesh()->GetComm());
247  HypreSmoother prec;
249  lin_solver.SetPreconditioner(prec);
250  lin_solver.SetOperator(*Mh);
251  lin_solver.SetRelTol(1e-8);
252  lin_solver.SetAbsTol(0.0);
253  lin_solver.SetMaxIter(100);
254  lin_solver.SetPrintLevel(0);
255  lin_solver.Mult(*RHS, X);
256  K.ParFESpace()->GetProlongationMatrix()->Mult(X, di_dt);
257 
258  delete Mh;
259  delete RHS;
260 }
261 #endif
262 
263 #ifdef MFEM_USE_GSLIB
265  const Vector &init_field)
266 {
267  nodes0 = init_nodes;
268  Mesh *m = mesh;
269 #ifdef MFEM_USE_MPI
270  if (pmesh) { m = pmesh; }
271 #endif
272  m->SetNodes(nodes0);
273 
274  const double rel_bbox_el = 0.1;
275  const double newton_tol = 1.0e-12;
276  const int npts_at_once = 256;
277 
278  if (finder)
279  {
280  finder->FreeData();
281  delete finder;
282  }
283 
285 #ifdef MFEM_USE_MPI
286  if (pfes)
287  {
288  f = pfes;
289  finder = new FindPointsGSLIB(pfes->GetComm());
290  }
291  else { finder = new FindPointsGSLIB(); }
292 #else
293  finder = new FindPointsGSLIB();
294 #endif
295  finder->Setup(*m, rel_bbox_el, newton_tol, npts_at_once);
296 
297  field0_gf.SetSpace(f);
298  field0_gf = init_field;
299 
300  dim = f->GetFE(0)->GetDim();
301 }
302 
304  Vector &new_field)
305 {
306  finder->Interpolate(new_nodes, field0_gf, new_field);
307 }
308 
309 #endif
310 
312  const Vector &b) const
313 {
314  const FiniteElementSpace *fes = NULL;
315  double energy_in = 0.0;
316 #ifdef MFEM_USE_MPI
317  const ParNonlinearForm *p_nlf = dynamic_cast<const ParNonlinearForm *>(oper);
318  MFEM_VERIFY(!(parallel && p_nlf == NULL), "Invalid Operator subclass.");
319  if (parallel)
320  {
321  fes = p_nlf->FESpace();
322  energy_in = p_nlf->GetEnergy(x);
323  }
324 #endif
325  const bool serial = !parallel;
326  const NonlinearForm *nlf = dynamic_cast<const NonlinearForm *>(oper);
327  MFEM_VERIFY(!(serial && nlf == NULL), "Invalid Operator subclass.");
328  if (serial)
329  {
330  fes = nlf->FESpace();
331  energy_in = nlf->GetEnergy(x);
332  }
333 
334  const int NE = fes->GetMesh()->GetNE(), dim = fes->GetMesh()->Dimension();
335  Array<int> xdofs;
336  DenseMatrix Jpr(dim);
337 
338  // Get the local prolongation of the solution vector.
339  Vector x_out_loc(fes->GetVSize());
340  if (serial)
341  {
342  const SparseMatrix *cP = fes->GetConformingProlongation();
343  if (!cP) { x_out_loc = x; }
344  else { cP->Mult(x, x_out_loc); }
345  }
346 #ifdef MFEM_USE_MPI
347  else
348  {
349  fes->GetProlongationMatrix()->Mult(x, x_out_loc);
350  }
351 #endif
352 
353  // Check if the starting mesh (given by x) is inverted.
354  // Note that x hasn't been modified by the Newton update yet.
355  double min_detJ = infinity();
356  for (int i = 0; i < NE; i++)
357  {
358  const int dof = fes->GetFE(i)->GetDof();
359  DenseMatrix dshape(dof, dim), pos(dof, dim);
360  Vector posV(pos.Data(), dof * dim);
361 
362  fes->GetElementVDofs(i, xdofs);
363  x_out_loc.GetSubVector(xdofs, posV);
364 
365  const IntegrationRule &irule = GetIntegrationRule(*fes->GetFE(i));
366  const int nsp = irule.GetNPoints();
367  for (int j = 0; j < nsp; j++)
368  {
369  fes->GetFE(i)->CalcDShape(irule.IntPoint(j), dshape);
370  MultAtB(pos, dshape, Jpr);
371  min_detJ = std::min(min_detJ, Jpr.Det());
372  }
373  }
374  double min_detJ_all = min_detJ;
375 #ifdef MFEM_USE_MPI
376  if (parallel)
377  {
378  MPI_Allreduce(&min_detJ, &min_detJ_all, 1, MPI_DOUBLE, MPI_MIN,
379  p_nlf->ParFESpace()->GetComm());
380  }
381 #endif
382  const bool untangling = (min_detJ_all <= 0) ? true : false;
383 
384  const bool have_b = (b.Size() == Height());
385 
386  Vector x_out(x.Size());
387  bool x_out_ok = false;
388  double scale = 1.0, energy_out = 0.0;
389  const double norm0 = Norm(r);
390 
391  const double detJ_factor = (solver_type == 1) ? 0.25 : 0.5;
392 
393  // Perform the line search.
394  for (int i = 0; i < 12; i++)
395  {
396  add(x, -scale, c, x_out);
397 
398  if (serial)
399  {
400  const SparseMatrix *cP = fes->GetConformingProlongation();
401  if (!cP) { x_out_loc = x_out; }
402  else { cP->Mult(x_out, x_out_loc); }
403  }
404 #ifdef MFEM_USE_MPI
405  else
406  {
407  fes->GetProlongationMatrix()->Mult(x_out, x_out_loc);
408  }
409 #endif
410 
411  // Check det(Jpr) > 0.
412  if (!untangling)
413  {
414  int jac_ok = 1;
415  for (int i = 0; i < NE; i++)
416  {
417  const int dof = fes->GetFE(i)->GetDof();
418  DenseMatrix dshape(dof, dim), pos(dof, dim);
419  Vector posV(pos.Data(), dof * dim);
420 
421  fes->GetElementVDofs(i, xdofs);
422  x_out_loc.GetSubVector(xdofs, posV);
423 
424  const IntegrationRule &irule = GetIntegrationRule(*fes->GetFE(i));
425  const int nsp = irule.GetNPoints();
426  for (int j = 0; j < nsp; j++)
427  {
428  fes->GetFE(i)->CalcDShape(irule.IntPoint(j), dshape);
429  MultAtB(pos, dshape, Jpr);
430  if (Jpr.Det() <= 0.0) { jac_ok = 0; goto break2; }
431  }
432  }
433 
434  break2:
435  int jac_ok_all = jac_ok;
436 #ifdef MFEM_USE_MPI
437  if (parallel)
438  {
439  MPI_Allreduce(&jac_ok, &jac_ok_all, 1, MPI_INT, MPI_LAND,
440  p_nlf->ParFESpace()->GetComm());
441  }
442 #endif
443 
444  if (jac_ok_all == 0)
445  {
446  if (print_level >= 0)
447  { mfem::out << "Scale = " << scale << " Neg det(J) found.\n"; }
448  scale *= detJ_factor; continue;
449  }
450  } // endif(!untangling)
451 
452  ProcessNewState(x_out);
453  if (serial)
454  {
455  energy_out = nlf->GetGridFunctionEnergy(x_out_loc);
456  }
457 #ifdef MFEM_USE_MPI
458  else
459  {
460  energy_out = p_nlf->GetParGridFunctionEnergy(x_out_loc);
461  }
462 #endif
463 
464  if (untangling)
465  {
466  if (energy_out > energy_in || std::isnan(energy_out) != 0)
467  {
468  scale *= 0.5;
469  }
470  else { x_out_ok = true; break; }
471  }
472  else
473  {
474  if (energy_out > 1.2*energy_in || std::isnan(energy_out) != 0)
475  {
476  if (print_level >= 0)
477  { mfem::out << "Scale = " << scale << " Increasing energy.\n"; }
478  scale *= 0.5; continue;
479  }
480 
481  oper->Mult(x_out, r);
482  if (have_b) { r -= b; }
483  double norm = Norm(r);
484 
485  if (norm > 1.2*norm0)
486  {
487  if (print_level >= 0)
488  { mfem::out << "Scale = " << scale << " Norm increased.\n"; }
489  scale *= 0.5; continue;
490  }
491  else { x_out_ok = true; break; }
492  } // endif (untangling)
493  } // enddo (i)
494 
495  if (print_level >= 0)
496  {
497  mfem::out << "Energy decrease: "
498  << (energy_in - energy_out) / energy_in * 100.0
499  << "% with " << scale << " scaling.\n";
500  }
501 
502  if (x_out_ok == false) { scale = 0.0; }
503  return scale;
504 }
505 
507 {
508  const NonlinearForm *nlf = dynamic_cast<const NonlinearForm *>(oper);
509  const Array<NonlinearFormIntegrator*> &integs = *nlf->GetDNFI();
510 
511  // Reset the update flags of all TargetConstructors.
512  // This is done to avoid repeated updates of shared TargetConstructors.
513  TMOP_Integrator *ti = NULL;
514  TMOPComboIntegrator *co = NULL;
515  DiscreteAdaptTC *dtc = NULL;
516  for (int i = 0; i < integs.Size(); i++)
517  {
518  ti = dynamic_cast<TMOP_Integrator *>(integs[i]);
519  if (ti)
520  {
521  dtc = ti->GetDiscreteAdaptTC();
522  if (dtc) { dtc->ResetUpdateFlags(); }
523  }
524  co = dynamic_cast<TMOPComboIntegrator *>(integs[i]);
525  if (co)
526  {
528  for (int j = 0; j < ati.Size(); j++)
529  {
530  dtc = ati[j]->GetDiscreteAdaptTC();
531  if (dtc) { dtc->ResetUpdateFlags(); }
532  }
533  }
534  }
535 
536  if (parallel)
537  {
538 #ifdef MFEM_USE_MPI
539  const ParNonlinearForm *nlf =
540  dynamic_cast<const ParNonlinearForm *>(oper);
541  const ParFiniteElementSpace *pfesc = nlf->ParFESpace();
542  Vector x_loc(pfesc->GetVSize());
543  pfesc->GetProlongationMatrix()->Mult(x, x_loc);
544  for (int i = 0; i < integs.Size(); i++)
545  {
546  ti = dynamic_cast<TMOP_Integrator *>(integs[i]);
547  if (ti)
548  {
549  ti->UpdateAfterMeshChange(x_loc);
550  ti->ComputeFDh(x_loc, *pfesc);
551  UpdateDiscreteTC(*ti, x_loc);
552  }
553  co = dynamic_cast<TMOPComboIntegrator *>(integs[i]);
554  if (co)
555  {
557  for (int j = 0; j < ati.Size(); j++)
558  {
559  ati[j]->ComputeFDh(x_loc, *pfesc);
560  UpdateDiscreteTC(*ati[j], x_loc);
561  }
562  }
563  }
564 #endif
565  }
566  else
567  {
568  const FiniteElementSpace *fesc = nlf->FESpace();
569  const Operator *P = nlf->GetProlongation();
570  Vector x_loc;
571  if (P)
572  {
573  x_loc.SetSize(P->Height());
574  P->Mult(x,x_loc);
575  }
576  else
577  {
578  x_loc = x;
579  }
580  for (int i = 0; i < integs.Size(); i++)
581  {
582  ti = dynamic_cast<TMOP_Integrator *>(integs[i]);
583  if (ti)
584  {
585  ti->UpdateAfterMeshChange(x_loc);
586  ti->ComputeFDh(x_loc, *fesc);
587  UpdateDiscreteTC(*ti, x_loc);
588  }
589  co = dynamic_cast<TMOPComboIntegrator *>(integs[i]);
590  if (co)
591  {
593  for (int j = 0; j < ati.Size(); j++)
594  {
595  ati[j]->ComputeFDh(x_loc, *fesc);
596  UpdateDiscreteTC(*ati[j], x_loc);
597  }
598  }
599  }
600  }
601 }
602 
604  const Vector &x_new) const
605 {
606  const bool update_flag = true;
607  DiscreteAdaptTC *discrtc = ti.GetDiscreteAdaptTC();
608  if (discrtc)
609  {
610  discrtc->UpdateTargetSpecification(x_new, update_flag);
611  if (ti.GetFDFlag())
612  {
613  double dx = ti.GetFDh();
614  discrtc->UpdateGradientTargetSpecification(x_new, dx, update_flag);
615  discrtc->UpdateHessianTargetSpecification(x_new, dx, update_flag);
616  }
617  }
618 }
619 
620 #ifdef MFEM_USE_MPI
621 // Metric values are visualized by creating an L2 finite element functions and
622 // computing the metric values at the nodes.
624  const TargetConstructor &tc, ParMesh &pmesh,
625  char *title, int position)
626 {
628  ParFiniteElementSpace fes(&pmesh, &fec, 1);
629  ParGridFunction metric(&fes);
630  InterpolateTMOP_QualityMetric(qm, tc, pmesh, metric);
631  socketstream sock;
632  if (pmesh.GetMyRank() == 0)
633  {
634  sock.open("localhost", 19916);
635  sock << "solution\n";
636  }
637  pmesh.PrintAsOne(sock);
638  metric.SaveAsOne(sock);
639  if (pmesh.GetMyRank() == 0)
640  {
641  sock << "window_title '"<< title << "'\n"
642  << "window_geometry "
643  << position << " " << 0 << " " << 600 << " " << 600 << "\n"
644  << "keys jRmclA\n";
645  }
646 }
647 #endif
648 
649 // Metric values are visualized by creating an L2 finite element functions and
650 // computing the metric values at the nodes.
652  const TargetConstructor &tc, Mesh &mesh,
653  char *title, int position)
654 {
656  FiniteElementSpace fes(&mesh, &fec, 1);
657  GridFunction metric(&fes);
658  InterpolateTMOP_QualityMetric(qm, tc, mesh, metric);
659  osockstream sock(19916, "localhost");
660  sock << "solution\n";
661  mesh.Print(sock);
662  metric.Save(sock);
663  sock.send();
664  sock << "window_title '"<< title << "'\n"
665  << "window_geometry "
666  << position << " " << 0 << " " << 600 << " " << 600 << "\n"
667  << "keys jRmclA\n";
668 }
669 
670 }
int GetNPoints() const
Returns the number of the points in the integration rule.
Definition: intrules.hpp:245
VectorGridFunctionCoefficient u_coeff
Definition: tmop_tools.hpp:77
void vis_tmop_metric_p(int order, TMOP_QualityMetric &qm, const TargetConstructor &tc, ParMesh &pmesh, char *title, int position)
Definition: tmop_tools.cpp:623
void Interpolate(const GridFunction &field_in, Vector &field_out)
Definition: gslib.cpp:575
virtual const Operator * GetProlongation() const
Get the finite element space prolongation matrix.
virtual double GetTime() const
Read the currently set time.
Definition: operator.hpp:308
int Size() const
Return the logical size of the array.
Definition: array.hpp:124
int GetVSize() const
Return the number of vector dofs, i.e. GetNDofs() x GetVDim().
Definition: fespace.hpp:400
virtual void Print(std::ostream &out=mfem::out) const
Definition: mesh.hpp:1234
Conjugate gradient method.
Definition: solvers.hpp:258
int GetDim() const
Returns the reference space dimension for the finite element.
Definition: fe.hpp:309
ParFiniteElementSpace * ParFESpace() const
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:2723
const SparseMatrix * GetConformingProlongation() const
The returned SparseMatrix is owned by the FiniteElementSpace.
Definition: fespace.cpp:874
const Operator * oper
Definition: solvers.hpp:73
Data type for scaled Jacobi-type smoother of sparse matrix.
ParMesh * GetParMesh() const
Definition: pfespace.hpp:243
void UpdateDiscreteTC(const TMOP_Integrator &ti, const Vector &x_new) const
Definition: tmop_tools.cpp:603
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:535
virtual const Operator * GetProlongationMatrix() const
The returned Operator is owned by the FiniteElementSpace.
Definition: pfespace.cpp:895
void SetSize(int s)
Resize the vector to size s.
Definition: vector.hpp:459
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:2881
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:173
double Det() const
Definition: densemat.cpp:451
Base abstract class for first order time dependent operators.
Definition: operator.hpp:267
void SaveAsOne(std::ostream &out=mfem::out)
Merge the local grid functions.
Definition: pgridfunc.cpp:855
ParFiniteElementSpace * pfes
Definition: tmop.hpp:560
Data type dense matrix using column-major storage.
Definition: densemat.hpp:23
int Size() const
Returns the size of the vector.
Definition: vector.hpp:160
Parallel non-linear operator on the true dofs.
const IntegrationRule & GetIntegrationRule(const FiniteElement &el) const
Definition: tmop_tools.hpp:123
int GetNE() const
Returns number of elements.
Definition: mesh.hpp:737
virtual void Mult(const Vector &x, Vector &y) const =0
Operator application: y=A(x).
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:651
void UpdateGradientTargetSpecification(const Vector &x, double dx, bool use_flag=false)
Definition: tmop.cpp:1769
double * GetData() const
Return a pointer to the beginning of the Vector data.
Definition: vector.hpp:169
ParBilinearForm M
Definition: tmop_tools.hpp:98
double GetFDh() const
Definition: tmop.hpp:1120
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:261
double Norm(const Vector &x) const
Definition: solvers.hpp:85
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:1188
virtual void ComputeAtNewPosition(const Vector &new_nodes, Vector &new_field)
Definition: tmop_tools.cpp:29
void Init(TimeDependentOperator &_f) override
Associate a TimeDependentOperator with the ODE solver.
Definition: ode.cpp:100
virtual void Save(std::ostream &out) const
Save the GridFunction to an output stream.
Definition: gridfunc.cpp:3417
double GetEnergy(const ParGridFunction &x) const
Compute the energy of a ParGridFunction.
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 SetPrintLevel(int print_lvl)
Definition: solvers.cpp:70
MPI_Comm GetComm() const
Definition: pfespace.hpp:239
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:65
ParBilinearForm K
Definition: tmop_tools.hpp:98
Mesh * GetMesh() const
Returns the mesh.
Definition: fespace.hpp:314
ParFiniteElementSpace * ParFESpace() const
Return the parallel FE space associated with the ParBilinearForm.
SerialAdvectorCGOper(const Vector &x_start, GridFunction &vel, FiniteElementSpace &fes)
Definition: tmop_tools.cpp:165
Array< NonlinearFormIntegrator * > * GetDNFI()
Access all integrators added with AddDomainIntegrator().
double b
Definition: lissajous.cpp:42
void PrintAsOne(std::ostream &out=mfem::out)
Definition: pmesh.cpp:4305
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:183
void SetMaxIter(int max_it)
Definition: solvers.hpp:98
virtual const Operator * GetProlongationMatrix() const
The returned Operator is owned by the FiniteElementSpace.
Definition: fespace.hpp:330
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:71
Performs a single remap advection step in serial.
Definition: tmop_tools.hpp:71
Parallel smoothers in hypre.
Definition: hypre.hpp:592
FiniteElementSpace * FESpace()
Definition: gridfunc.hpp:595
int Dimension() const
Definition: mesh.hpp:788
double GetParGridFunctionEnergy(const Vector &x) const
Compute the energy corresponding to the state x.
FiniteElementSpace * FESpace()
double GetElementSize(int i, int type=0)
Get the size of the i-th element relative to the perfect reference element.
Definition: mesh.cpp:74
Wrapper for hypre&#39;s parallel vector class.
Definition: hypre.hpp:70
double Min() const
Returns the minimal element of the vector.
Definition: vector.cpp:1105
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:506
double * Data() const
Returns the matrix data array.
Definition: densemat.hpp:96
void SetAbsTol(double atol)
Definition: solvers.hpp:97
int GetMyRank() const
Definition: pmesh.hpp:238
void SetRelTol(double rtol)
Definition: solvers.hpp:96
void UpdateHessianTargetSpecification(const Vector &x, double dx, bool use_flag=false)
Definition: tmop.cpp:1795
MPI_Comm GetComm() const
Definition: pmesh.hpp:236
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:315
virtual void Mult(const Vector &x, Vector &y) const
Matrix vector multiplication.
Definition: sparsemat.cpp:594
ParAdvectorCGOper(const Vector &x_start, GridFunction &vel, ParFiniteElementSpace &pfes)
Definition: tmop_tools.cpp:209
void subtract(const Vector &x, const Vector &y, Vector &z)
Definition: vector.cpp:408
Performs a single remap advection step in parallel.
Definition: tmop_tools.hpp:91
DiscreteAdaptTC * GetDiscreteAdaptTC() const
Definition: tmop.hpp:1103
void ResetUpdateFlags()
Used in combination with the Update methods to avoid extra computations.
Definition: tmop.hpp:819
virtual void ComputeAtNewPosition(const Vector &new_nodes, Vector &new_field)
Definition: tmop_tools.cpp:303
void ComputeFDh(const Vector &x, const FiniteElementSpace &fes)
Determines the perturbation, h, for FD-based approximation.
Definition: tmop.cpp:2729
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:832
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:97
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.
const FiniteElement * GetFE(int i) const
Returns pointer to the FiniteElement in the FiniteElementCollection associated with i&#39;th element in t...
Definition: fespace.cpp:1798
double infinity()
Define a shortcut for std::numeric_limits&lt;double&gt;::infinity()
Definition: vector.hpp:45
virtual void SetSpace(FiniteElementSpace *f)
Associate a new FiniteElementSpace with the GridFunction.
Definition: gridfunc.cpp:180
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:227
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
Definition: solvers.hpp:272
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:2650
FindPointsGSLIB can robustly evaluate a GridFunction on an arbitrary collection of points...
Definition: gslib.hpp:45
bool GetFDFlag() const
Definition: tmop.hpp:1119
Vector data type.
Definition: vector.hpp:51
const FiniteElementCollection * FEColl() const
Definition: fespace.hpp:414
void GetNodes(Vector &node_coord) const
Definition: mesh.cpp:6603
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:311
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
Definition: solvers.cpp:91
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
const SparseMatrix & SpMat() const
Returns a const reference to the sparse matrix.
void SetNodes(const Vector &node_coord)
Definition: mesh.cpp:6615
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:603
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
Wrapper for hypre&#39;s ParCSR matrix class.
Definition: hypre.hpp:181
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
const Array< TMOP_Integrator * > & GetTMOPIntegrators() const
Definition: tmop.hpp:1143
void SetType(HypreSmoother::Type type, int relax_times=1)
Set the relaxation type and number of sweeps.
Definition: hypre.cpp:2187
virtual void SetInitialField(const Vector &init_nodes, const Vector &init_field)
Definition: tmop_tools.cpp:264
FiniteElementSpace * fes
Definition: tmop.hpp:555
Arbitrary order &quot;L2-conforming&quot; discontinuous finite elements.
Definition: fe_coll.hpp:221
double f(const Vector &p)
alpha (q . grad u, v)
A TMOP integrator class based on any given TMOP_QualityMetric and TargetConstructor.
Definition: tmop.hpp:884