MFEM  v4.4.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
tmop_tools.hpp
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 #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 
47  /// Set the memory type used for large memory allocations. This memory type
48  /// is used when constructing the AdvectorCGOper but currently only for the
49  /// parallel variant.
50  void SetMemoryType(MemoryType mt) { opt_mt = mt; }
51 };
52 
53 #ifdef MFEM_USE_GSLIB
55 {
56 private:
57  Vector nodes0;
58  GridFunction field0_gf;
59  FindPointsGSLIB *finder;
60  int dim;
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 
71  {
72  finder->FreeData();
73  delete finder;
74  }
75 };
76 #endif
77 
78 /// Performs a single remap advection step in serial.
80 {
81 protected:
82  const Vector &x0;
86  mutable BilinearForm M, K;
88 
89 public:
90  /** Here @a fes is the FESpace of the function that will be moved. Note
91  that Mult() moves the nodes of the mesh corresponding to @a fes. */
93  FiniteElementSpace &fes,
95 
96  virtual void Mult(const Vector &ind, Vector &di_dt) const;
97 };
98 
99 #ifdef MFEM_USE_MPI
100 /// Performs a single remap advection step in parallel.
102 {
103 protected:
104  const Vector &x0;
108  mutable ParBilinearForm M, K;
110 
111 public:
112  /** Here @a pfes is the ParFESpace of the function that will be moved. Note
113  that Mult() moves the nodes of the mesh corresponding to @a pfes.
114  @a mt is used to set the memory type of the integrators. */
115  ParAdvectorCGOper(const Vector &x_start, GridFunction &vel,
116  ParFiniteElementSpace &pfes,
119 
120  virtual void Mult(const Vector &ind, Vector &di_dt) const;
121 };
122 #endif
123 
125 {
126 protected:
127  // 0 - Newton, 1 - LBFGS.
129  bool parallel;
130 
131  // Surface fitting variables.
132  bool adaptive_surf_fit = false;
133  mutable double surf_fit_err_avg_prvs = 10000.0;
134  mutable bool update_surf_fit_coeff = false;
135  double surf_fit_max_threshold = -1.0;
136 
137  // Minimum determinant over the whole mesh. Used for mesh untangling.
138  double *min_det_ptr = nullptr;
139 
140  // Quadrature points that are checked for negative Jacobians etc.
142  // These fields are relevant for mixed meshes.
145 
147 
149  {
150  if (IntegRules)
151  {
152  return IntegRules->Get(el.GetGeomType(), integ_order);
153  }
154  return ir;
155  }
156 
157  void UpdateDiscreteTC(const TMOP_Integrator &ti, const Vector &x_new) const;
158 
159  double ComputeMinDet(const Vector &x_loc,
160  const FiniteElementSpace &fes) const;
161 
162  double MinDetJpr_2D(const FiniteElementSpace*, const Vector&) const;
163  double MinDetJpr_3D(const FiniteElementSpace*, const Vector&) const;
164 
165  /** @name Methods for adaptive surface fitting weight. */
166  ///@{
167  /// Get the average and maximum surface fitting error at the marked nodes.
168  /// If there is more than 1 TMOP integrator, we get the maximum of the
169  /// average and maximum error over all integrators.
170  virtual void GetSurfaceFittingError(double &err_avg, double &err_max) const;
171 
172  /// Update surface fitting weight as surf_fit_weight *= factor.
173  void UpdateSurfaceFittingWeight(double factor) const;
174 
175  /// Get the surface fitting weight for all the TMOP integrators.
176  void GetSurfaceFittingWeight(Array<double> &weights) const;
177  ///@}
178 
179 public:
180 #ifdef MFEM_USE_MPI
181  TMOPNewtonSolver(MPI_Comm comm, const IntegrationRule &irule, int type = 0)
182  : LBFGSSolver(comm), solver_type(type), parallel(true),
183  ir(irule), IntegRules(NULL), integ_order(-1) { }
184 #endif
185  TMOPNewtonSolver(const IntegrationRule &irule, int type = 0)
186  : LBFGSSolver(), solver_type(type), parallel(false),
187  ir(irule), IntegRules(NULL), integ_order(-1) { }
188 
189  /// Prescribe a set of integration rules; relevant for mixed meshes.
190  /** If called, this function has priority over the IntegrationRule given to
191  the constructor of the class. */
192  void SetIntegrationRules(IntegrationRules &irules, int order)
193  {
194  IntegRules = &irules;
195  integ_order = order;
196  }
197 
198  void SetMinDetPtr(double *md_ptr) { min_det_ptr = md_ptr; }
199 
200  /// Set the memory type for temporary memory allocations.
202 
203  /// Compute scaling factor for the node movement direction using line-search.
204  /// We impose constraints on TMOP energy, gradient, minimum Jacobian of
205  /// the mesh, and (optionally) on the surface fitting error.
206  virtual double ComputeScalingFactor(const Vector &x, const Vector &b) const;
207 
208  /// Update (i) discrete functions at new nodal positions, and
209  /// (ii) surface fitting weight.
210  virtual void ProcessNewState(const Vector &x) const;
211 
212  /** @name Methods for adaptive surface fitting weight. (Experimental) */
213  /// Enable adaptive surface fitting weight.
214  /// The weight is modified after each TMOPNewtonSolver iteration.
216 
217  /// Set the termination criterion for mesh optimization based on
218  /// the maximum surface fitting error.
220  {
221  surf_fit_max_threshold = max_error;
222  }
223 
224  virtual void Mult(const Vector &b, Vector &x) const
225  {
226  if (solver_type == 0)
227  {
228  NewtonSolver::Mult(b, x);
229  }
230  else if (solver_type == 1)
231  {
232  LBFGSSolver::Mult(b, x);
233  }
234  else { MFEM_ABORT("Invalid type"); }
235  }
236 
237  virtual void SetSolver(Solver &solver)
238  {
239  if (solver_type == 0)
240  {
241  NewtonSolver::SetSolver(solver);
242  }
243  else if (solver_type == 1)
244  {
245  LBFGSSolver::SetSolver(solver);
246  }
247  else { MFEM_ABORT("Invalid type"); }
248  }
249  virtual void SetPreconditioner(Solver &pr) { SetSolver(pr); }
250 };
251 
252 void vis_tmop_metric_s(int order, TMOP_QualityMetric &qm,
253  const TargetConstructor &tc, Mesh &pmesh,
254  char *title, int position);
255 #ifdef MFEM_USE_MPI
256 void vis_tmop_metric_p(int order, TMOP_QualityMetric &qm,
257  const TargetConstructor &tc, ParMesh &pmesh,
258  char *title, int position);
259 #endif
260 
261 }
262 
263 #endif
Abstract class for all finite elements.
Definition: fe_base.hpp:235
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
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
const IntegrationRule & Get(int GeomType, int Order)
Returns an integration rule for given GeomType and Order.
Definition: intrules.cpp:923
void UpdateDiscreteTC(const TMOP_Integrator &ti, const Vector &x_new) const
Definition: tmop_tools.cpp:793
void SetIntegrationRules(IntegrationRules &irules, int order)
Prescribe a set of integration rules; relevant for mixed meshes.
Definition: tmop_tools.hpp:192
const IntegrationRule & ir
Definition: tmop_tools.hpp:141
double MinDetJpr_2D(const FiniteElementSpace *, const Vector &) const
Base abstract class for first order time dependent operators.
Definition: operator.hpp:285
Container class for integration rules.
Definition: intrules.hpp:311
const IntegrationRule & GetIntegrationRule(const FiniteElement &el) const
Definition: tmop_tools.hpp:148
virtual void SetSolver(Solver &solver)
Set the linear solver for inverting the Jacobian.
Definition: solvers.hpp:661
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 SetMinDetPtr(double *md_ptr)
Definition: tmop_tools.hpp:198
virtual void ComputeAtNewPosition(const Vector &new_nodes, Vector &new_field)
Definition: tmop_tools.cpp:29
void SetTerminationWithMaxSurfaceFittingError(double max_error)
Definition: tmop_tools.hpp:219
Geometry::Type GetGeomType() const
Returns the Geometry::Type of the reference element.
Definition: fe_base.hpp:320
void GetSurfaceFittingWeight(Array< double > &weights) const
Get the surface fitting weight for all the TMOP integrators.
Definition: tmop_tools.cpp:589
virtual void FreeData()
Definition: gslib.cpp:214
double b
Definition: lissajous.cpp:42
const AssemblyLevel al
Definition: tmop_tools.hpp:87
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
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
void EnableAdaptiveSurfaceFitting()
Definition: tmop_tools.hpp:215
TMOPNewtonSolver(const IntegrationRule &irule, int type=0)
Definition: tmop_tools.hpp:185
virtual void Mult(const Vector &b, Vector &x) const
Solve the nonlinear system with right-hand side b.
Definition: solvers.cpp:1804
IntegrationRules * IntegRules
Definition: tmop_tools.hpp:143
virtual void ProcessNewState(const Vector &x) const
Definition: tmop_tools.cpp:659
The classical explicit forth-order Runge-Kutta method, RK4.
Definition: ode.hpp:162
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
Definition: tmop_tools.hpp:249
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition: fespace.hpp:88
virtual void Mult(const Vector &b, Vector &x) const
Solve the nonlinear system with right-hand side b.
Definition: solvers.cpp:1992
AdvectorCG(AssemblyLevel al=AssemblyLevel::LEGACY, double timestep_scale=0.5)
Definition: tmop_tools.hpp:36
const AssemblyLevel al
Definition: tmop_tools.hpp:109
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:181
Performs a single remap advection step in parallel.
Definition: tmop_tools.hpp:101
virtual void SetSolver(Solver &solver)
Set the linear solver for inverting the Jacobian.
Definition: solvers.hpp:752
virtual void ComputeAtNewPosition(const Vector &new_nodes, Vector &new_field)
Definition: tmop_tools.cpp:354
A &quot;square matrix&quot; operator for the associated FE space and BLFIntegrators The sum of all the BLFInteg...
void SetMemoryType(MemoryType mt)
Definition: tmop_tools.hpp:50
double MinDetJpr_3D(const FiniteElementSpace *, const Vector &) const
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 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
void SetTempMemoryType(MemoryType mt)
Set the memory type for temporary memory allocations.
Definition: tmop_tools.hpp:201
Class for parallel bilinear form.
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
virtual void Mult(const Vector &b, Vector &x) const
Solve the nonlinear system with right-hand side b.
Definition: tmop_tools.hpp:224
ParAdvectorCGOper(const Vector &x_start, GridFunction &vel, ParFiniteElementSpace &pfes, AssemblyLevel al=AssemblyLevel::LEGACY, MemoryType mt=MemoryType::DEFAULT)
Definition: tmop_tools.cpp:231
FindPointsGSLIB can robustly evaluate a GridFunction on an arbitrary collection of points...
Definition: gslib.hpp:47
Vector data type.
Definition: vector.hpp:60
virtual double ComputeScalingFactor(const Vector &x, const Vector &b) const
Definition: tmop_tools.cpp:362
Vector coefficient defined by a vector GridFunction.
Base class for solvers.
Definition: operator.hpp:651
double ComputeMinDet(const Vector &x_loc, const FiniteElementSpace &fes) const
Definition: tmop_tools.cpp:810
virtual void SetSolver(Solver &solver)
Set the linear solver for inverting the Jacobian.
Definition: tmop_tools.hpp:237
virtual void SetInitialField(const Vector &init_nodes, const Vector &init_field)
Definition: tmop_tools.cpp:315
A TMOP integrator class based on any given TMOP_QualityMetric and TargetConstructor.
Definition: tmop.hpp:1303