MFEM  v4.5.2
Finite element discretization library
tmop_tools.hpp
Go to the documentation of this file.
1 // Copyright (c) 2010-2023, 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 #ifndef MFEM_TMOP_TOOLS_HPP
13 #define MFEM_TMOP_TOOLS_HPP
14 
15 #include "bilinearform.hpp"
16 #include "pbilinearform.hpp"
17 #include "tmop.hpp"
18 #include "gslib.hpp"
19 
20 namespace mfem
21 {
22 
23 // Performs the full remap advection loop.
25 {
26 private:
27  RK4Solver ode_solver;
28  Vector nodes0;
29  Vector field0;
30  const double dt_scale;
31  const AssemblyLevel al;
33 
34  void ComputeAtNewPositionScalar(const Vector &new_nodes, Vector &new_field);
35 public:
37  double timestep_scale = 0.5)
39  ode_solver(), nodes0(), field0(), dt_scale(timestep_scale), al(al) { }
40 
41  virtual void SetInitialField(const Vector &init_nodes,
42  const Vector &init_field);
43 
44  virtual void ComputeAtNewPosition(const Vector &new_nodes,
45  Vector &new_field,
46  int new_nodes_ordering = Ordering::byNODES);
47 
48  /// Set the memory type used for large memory allocations. This memory type
49  /// is used when constructing the AdvectorCGOper but currently only for the
50  /// parallel variant.
51  void SetMemoryType(MemoryType mt) { opt_mt = mt; }
52 };
53 
54 #ifdef MFEM_USE_GSLIB
56 {
57 private:
58  Vector nodes0;
59  GridFunction field0_gf;
60  FindPointsGSLIB *finder;
61 public:
62  InterpolatorFP() : finder(NULL) { }
63 
64  virtual void SetInitialField(const Vector &init_nodes,
65  const Vector &init_field);
66 
67  virtual void ComputeAtNewPosition(const Vector &new_nodes,
68  Vector &new_field,
69  int new_nodes_ordering = Ordering::byNODES);
70 
72  {
73  return finder;
74  }
75 
77  {
78  finder->FreeData();
79  delete finder;
80  }
81 };
82 #endif
83 
84 /// Performs a single remap advection step in serial.
86 {
87 protected:
88  const Vector &x0;
92  mutable BilinearForm M, K;
94 
95 public:
96  /** Here @a fes is the FESpace of the function that will be moved. Note
97  that Mult() moves the nodes of the mesh corresponding to @a fes. */
99  FiniteElementSpace &fes,
101 
102  virtual void Mult(const Vector &ind, Vector &di_dt) const;
103 };
104 
105 #ifdef MFEM_USE_MPI
106 /// Performs a single remap advection step in parallel.
108 {
109 protected:
110  const Vector &x0;
114  mutable ParBilinearForm M, K;
116 
117 public:
118  /** Here @a pfes is the ParFESpace of the function that will be moved. Note
119  that Mult() moves the nodes of the mesh corresponding to @a pfes.
120  @a mt is used to set the memory type of the integrators. */
121  ParAdvectorCGOper(const Vector &x_start, GridFunction &vel,
122  ParFiniteElementSpace &pfes,
125 
126  virtual void Mult(const Vector &ind, Vector &di_dt) const;
127 };
128 #endif
129 
131 {
132 protected:
133  // 0 - Newton, 1 - LBFGS.
135  bool parallel;
136 
137  // Surface fitting variables.
138  bool adaptive_surf_fit = false;
139  mutable double surf_fit_err_avg_prvs = 10000.0;
140  mutable bool update_surf_fit_coeff = false;
141  double surf_fit_max_threshold = -1.0;
142 
143  // Minimum determinant over the whole mesh. Used for mesh untangling.
144  double *min_det_ptr = nullptr;
145  // Flag to compute minimum determinant and maximum metric in ProcessNewState,
146  // which is required for TMOP_WorstCaseUntangleOptimizer_Metric.
147  mutable bool compute_metric_quantile_flag = true;
148 
149  // Quadrature points that are checked for negative Jacobians etc.
151  // These fields are relevant for mixed meshes.
154 
156 
158  {
159  if (IntegRules)
160  {
161  return IntegRules->Get(el.GetGeomType(), integ_order);
162  }
163  return ir;
164  }
165 
166  void UpdateDiscreteTC(const TMOP_Integrator &ti, const Vector &x_new,
167  int x_ordering = Ordering::byNODES) const;
168 
169  double ComputeMinDet(const Vector &x_loc,
170  const FiniteElementSpace &fes) const;
171 
172  double MinDetJpr_2D(const FiniteElementSpace*, const Vector&) const;
173  double MinDetJpr_3D(const FiniteElementSpace*, const Vector&) const;
174 
175  /** @name Methods for adaptive surface fitting weight. */
176  ///@{
177  /// Get the average and maximum surface fitting error at the marked nodes.
178  /// If there is more than 1 TMOP integrator, we get the maximum of the
179  /// average and maximum error over all integrators.
180  virtual void GetSurfaceFittingError(double &err_avg, double &err_max) const;
181 
182  /// Update surface fitting weight as surf_fit_weight *= factor.
183  void UpdateSurfaceFittingWeight(double factor) const;
184 
185  /// Get the surface fitting weight for all the TMOP integrators.
186  void GetSurfaceFittingWeight(Array<double> &weights) const;
187  ///@}
188 
189 public:
190 #ifdef MFEM_USE_MPI
191  TMOPNewtonSolver(MPI_Comm comm, const IntegrationRule &irule, int type = 0)
192  : LBFGSSolver(comm), solver_type(type), parallel(true),
193  ir(irule), IntegRules(NULL), integ_order(-1) { }
194 #endif
195  TMOPNewtonSolver(const IntegrationRule &irule, int type = 0)
196  : LBFGSSolver(), solver_type(type), parallel(false),
197  ir(irule), IntegRules(NULL), integ_order(-1) { }
198 
199  /// Prescribe a set of integration rules; relevant for mixed meshes.
200  /** If called, this function has priority over the IntegrationRule given to
201  the constructor of the class. */
202  void SetIntegrationRules(IntegrationRules &irules, int order)
203  {
204  IntegRules = &irules;
205  integ_order = order;
206  }
207 
208  void SetMinDetPtr(double *md_ptr) { min_det_ptr = md_ptr; }
209 
210  /// Set the memory type for temporary memory allocations.
212 
213  /// Compute scaling factor for the node movement direction using line-search.
214  /// We impose constraints on TMOP energy, gradient, minimum Jacobian of
215  /// the mesh, and (optionally) on the surface fitting error.
216  virtual double ComputeScalingFactor(const Vector &x, const Vector &b) const;
217 
218  /// Update (i) discrete functions at new nodal positions, and
219  /// (ii) surface fitting weight.
220  virtual void ProcessNewState(const Vector &x) const;
221 
222  /** @name Methods for adaptive surface fitting weight. (Experimental) */
223  /// Enable adaptive surface fitting weight.
224  /// The weight is modified after each TMOPNewtonSolver iteration.
226 
227  /// Set the termination criterion for mesh optimization based on
228  /// the maximum surface fitting error.
230  {
231  surf_fit_max_threshold = max_error;
232  }
233 
234  virtual void Mult(const Vector &b, Vector &x) const
235  {
236  if (solver_type == 0)
237  {
238  NewtonSolver::Mult(b, x);
239  }
240  else if (solver_type == 1)
241  {
242  LBFGSSolver::Mult(b, x);
243  }
244  else { MFEM_ABORT("Invalid type"); }
245  }
246 
247  virtual void SetSolver(Solver &solver)
248  {
249  if (solver_type == 0)
250  {
251  NewtonSolver::SetSolver(solver);
252  }
253  else if (solver_type == 1)
254  {
255  LBFGSSolver::SetSolver(solver);
256  }
257  else { MFEM_ABORT("Invalid type"); }
258  }
259  virtual void SetPreconditioner(Solver &pr) { SetSolver(pr); }
260 };
261 
262 void vis_tmop_metric_s(int order, TMOP_QualityMetric &qm,
263  const TargetConstructor &tc, Mesh &pmesh,
264  char *title, int position);
265 #ifdef MFEM_USE_MPI
266 void vis_tmop_metric_p(int order, TMOP_QualityMetric &qm,
267  const TargetConstructor &tc, ParMesh &pmesh,
268  char *title, int position);
269 #endif
270 
271 }
272 
273 #endif
Abstract class for all finite elements.
Definition: fe_base.hpp:232
VectorGridFunctionCoefficient u_coeff
Definition: tmop_tools.hpp:91
void vis_tmop_metric_p(int order, TMOP_QualityMetric &qm, const TargetConstructor &tc, ParMesh &pmesh, char *title, int position)
Definition: tmop_tools.cpp:899
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
Class for grid function - Vector with associated FE space.
Definition: gridfunc.hpp:30
virtual void ComputeAtNewPosition(const Vector &new_nodes, Vector &new_field, int new_nodes_ordering=Ordering::byNODES)
Definition: tmop_tools.cpp:370
const IntegrationRule & Get(int GeomType, int Order)
Returns an integration rule for given GeomType and Order.
Definition: intrules.cpp:923
void SetIntegrationRules(IntegrationRules &irules, int order)
Prescribe a set of integration rules; relevant for mixed meshes.
Definition: tmop_tools.hpp:202
const IntegrationRule & ir
Definition: tmop_tools.hpp:150
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:214
void GetSurfaceFittingWeight(Array< double > &weights) const
Get the surface fitting weight for all the TMOP integrators.
Definition: tmop_tools.cpp:608
double ComputeMinDet(const Vector &x_loc, const FiniteElementSpace &fes) const
Definition: tmop_tools.cpp:846
Base abstract class for first order time dependent operators.
Definition: operator.hpp:316
virtual void Mult(const Vector &b, Vector &x) const
Solve the nonlinear system with right-hand side b.
Definition: solvers.cpp:2005
Container class for integration rules.
Definition: intrules.hpp:311
double MinDetJpr_2D(const FiniteElementSpace *, const Vector &) const
Definition: tmop_pa_jp2.cpp:75
virtual void SetSolver(Solver &solver)
Set the linear solver for inverting the Jacobian.
Definition: solvers.hpp:689
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:927
void SetMinDetPtr(double *md_ptr)
Definition: tmop_tools.hpp:208
virtual void ComputeAtNewPosition(const Vector &new_nodes, Vector &new_field, int new_nodes_ordering=Ordering::byNODES)
Definition: tmop_tools.cpp:29
void SetTerminationWithMaxSurfaceFittingError(double max_error)
Definition: tmop_tools.hpp:229
Geometry::Type GetGeomType() const
Returns the Geometry::Type of the reference element.
Definition: fe_base.hpp:319
virtual void FreeData()
Definition: gslib.cpp:267
double b
Definition: lissajous.cpp:42
const AssemblyLevel al
Definition: tmop_tools.hpp:93
virtual double ComputeScalingFactor(const Vector &x, const Vector &b) const
Definition: tmop_tools.cpp:379
Performs a single remap advection step in serial.
Definition: tmop_tools.hpp:85
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:283
void EnableAdaptiveSurfaceFitting()
Definition: tmop_tools.hpp:225
TMOPNewtonSolver(const IntegrationRule &irule, int type=0)
Definition: tmop_tools.hpp:195
IntegrationRules * IntegRules
Definition: tmop_tools.hpp:152
const FindPointsGSLIB * GetFindPointsGSLIB() const
Definition: tmop_tools.hpp:71
virtual void Mult(const Vector &b, Vector &x) const
Solve the nonlinear system with right-hand side b.
Definition: solvers.cpp:1817
The classical explicit forth-order Runge-Kutta method, RK4.
Definition: ode.hpp:162
virtual void ProcessNewState(const Vector &x) const
Definition: tmop_tools.cpp:678
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
Definition: tmop_tools.hpp:259
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition: fespace.hpp:96
AdvectorCG(AssemblyLevel al=AssemblyLevel::LEGACY, double timestep_scale=0.5)
Definition: tmop_tools.hpp:36
const AssemblyLevel al
Definition: tmop_tools.hpp:115
MemoryType
Memory types supported by MFEM.
Definition: mem_manager.hpp:31
TMOPNewtonSolver(MPI_Comm comm, const IntegrationRule &irule, int type=0)
Definition: tmop_tools.hpp:191
Performs a single remap advection step in parallel.
Definition: tmop_tools.hpp:107
double MinDetJpr_3D(const FiniteElementSpace *, const Vector &) const
Definition: tmop_pa_jp3.cpp:77
virtual void SetSolver(Solver &solver)
Set the linear solver for inverting the Jacobian.
Definition: solvers.hpp:780
const IntegrationRule & GetIntegrationRule(const FiniteElement &el) const
Definition: tmop_tools.hpp:157
A "square matrix" operator for the associated FE space and BLFIntegrators The sum of all the BLFInteg...
void SetMemoryType(MemoryType mt)
Definition: tmop_tools.hpp:51
virtual void SetInitialField(const Vector &init_nodes, const Vector &init_field)
Definition: tmop_tools.cpp:22
VectorGridFunctionCoefficient u_coeff
Definition: tmop_tools.hpp:113
SerialAdvectorCGOper(const Vector &x_start, GridFunction &vel, FiniteElementSpace &fes, AssemblyLevel al=AssemblyLevel::LEGACY)
Definition: tmop_tools.cpp:193
void SetTempMemoryType(MemoryType mt)
Set the memory type for temporary memory allocations.
Definition: tmop_tools.hpp:211
Class for parallel bilinear form.
virtual void GetSurfaceFittingError(double &err_avg, double &err_max) const
Definition: tmop_tools.cpp:638
ParAdvectorCGOper(const Vector &x_start, GridFunction &vel, ParFiniteElementSpace &pfes, AssemblyLevel al=AssemblyLevel::LEGACY, MemoryType mt=MemoryType::DEFAULT)
Definition: tmop_tools.cpp:253
FindPointsGSLIB can robustly evaluate a GridFunction on an arbitrary collection of points...
Definition: gslib.hpp:52
Vector data type.
Definition: vector.hpp:60
Vector coefficient defined by a vector GridFunction.
Base class for solvers.
Definition: operator.hpp:682
virtual void Mult(const Vector &b, Vector &x) const
Solve the nonlinear system with right-hand side b.
Definition: tmop_tools.hpp:234
virtual void SetSolver(Solver &solver)
Set the linear solver for inverting the Jacobian.
Definition: tmop_tools.hpp:247
void UpdateSurfaceFittingWeight(double factor) const
Update surface fitting weight as surf_fit_weight *= factor.
Definition: tmop_tools.cpp:583
virtual void SetInitialField(const Vector &init_nodes, const Vector &init_field)
Definition: tmop_tools.cpp:333
A TMOP integrator class based on any given TMOP_QualityMetric and TargetConstructor.
Definition: tmop.hpp:1638
void UpdateDiscreteTC(const TMOP_Integrator &ti, const Vector &x_new, int x_ordering=Ordering::byNODES) const
Definition: tmop_tools.cpp:828