MFEM  v4.3.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-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 #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  // Minimum determinant over the whole mesh. Used for mesh untangling.
132  double *min_det_ptr = nullptr;
133 
134  // Quadrature points that are checked for negative Jacobians etc.
136  // These fields are relevant for mixed meshes.
139 
141 
143  {
144  if (IntegRules)
145  {
146  return IntegRules->Get(el.GetGeomType(), integ_order);
147  }
148  return ir;
149  }
150 
151  void UpdateDiscreteTC(const TMOP_Integrator &ti, const Vector &x_new) const;
152 
153  double ComputeMinDet(const Vector &x_loc,
154  const FiniteElementSpace &fes) const;
155 
156  double MinDetJpr_2D(const FiniteElementSpace*, const Vector&) const;
157  double MinDetJpr_3D(const FiniteElementSpace*, const Vector&) const;
158 
159 public:
160 #ifdef MFEM_USE_MPI
161  TMOPNewtonSolver(MPI_Comm comm, const IntegrationRule &irule, int type = 0)
162  : LBFGSSolver(comm), solver_type(type), parallel(true),
163  ir(irule), IntegRules(NULL), integ_order(-1) { }
164 #endif
165  TMOPNewtonSolver(const IntegrationRule &irule, int type = 0)
166  : LBFGSSolver(), solver_type(type), parallel(false),
167  ir(irule), IntegRules(NULL), integ_order(-1) { }
168 
169  /// Prescribe a set of integration rules; relevant for mixed meshes.
170  /** If called, this function has priority over the IntegrationRule given to
171  the constructor of the class. */
172  void SetIntegrationRules(IntegrationRules &irules, int order)
173  {
174  IntegRules = &irules;
175  integ_order = order;
176  }
177 
178  void SetMinDetPtr(double *md_ptr) { min_det_ptr = md_ptr; }
179 
180  // Set the memory type for temporary memory allocations.
182 
183  virtual double ComputeScalingFactor(const Vector &x, const Vector &b) const;
184 
185  virtual void ProcessNewState(const Vector &x) const;
186 
187  virtual void Mult(const Vector &b, Vector &x) const
188  {
189  if (solver_type == 0)
190  {
191  NewtonSolver::Mult(b, x);
192  }
193  else if (solver_type == 1)
194  {
195  LBFGSSolver::Mult(b, x);
196  }
197  else { MFEM_ABORT("Invalid type"); }
198  }
199 
200  virtual void SetSolver(Solver &solver)
201  {
202  if (solver_type == 0)
203  {
204  NewtonSolver::SetSolver(solver);
205  }
206  else if (solver_type == 1)
207  {
208  LBFGSSolver::SetSolver(solver);
209  }
210  else { MFEM_ABORT("Invalid type"); }
211  }
212  virtual void SetPreconditioner(Solver &pr) { SetSolver(pr); }
213 };
214 
215 void vis_tmop_metric_s(int order, TMOP_QualityMetric &qm,
216  const TargetConstructor &tc, Mesh &pmesh,
217  char *title, int position);
218 #ifdef MFEM_USE_MPI
219 void vis_tmop_metric_p(int order, TMOP_QualityMetric &qm,
220  const TargetConstructor &tc, ParMesh &pmesh,
221  char *title, int position);
222 #endif
223 
224 }
225 
226 #endif
Abstract class for all finite elements.
Definition: fe.hpp:243
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
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:920
void UpdateDiscreteTC(const TMOP_Integrator &ti, const Vector &x_new) const
Definition: tmop_tools.cpp:635
void SetIntegrationRules(IntegrationRules &irules, int order)
Prescribe a set of integration rules; relevant for mixed meshes.
Definition: tmop_tools.hpp:172
const IntegrationRule & ir
Definition: tmop_tools.hpp:135
double MinDetJpr_2D(const FiniteElementSpace *, const Vector &) const
Base abstract class for first order time dependent operators.
Definition: operator.hpp:282
Container class for integration rules.
Definition: intrules.hpp:311
const IntegrationRule & GetIntegrationRule(const FiniteElement &el) const
Definition: tmop_tools.hpp:142
virtual void SetSolver(Solver &solver)
Set the linear solver for inverting the Jacobian.
Definition: solvers.hpp:512
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 SetMinDetPtr(double *md_ptr)
Definition: tmop_tools.hpp:178
virtual void ComputeAtNewPosition(const Vector &new_nodes, Vector &new_field)
Definition: tmop_tools.cpp:29
Geometry::Type GetGeomType() const
Returns the Geometry::Type of the reference element.
Definition: fe.hpp:320
virtual void FreeData()
Definition: gslib.cpp:212
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
TMOPNewtonSolver(const IntegrationRule &irule, int type=0)
Definition: tmop_tools.hpp:165
virtual void Mult(const Vector &b, Vector &x) const
Solve the nonlinear system with right-hand side b.
Definition: solvers.cpp:1653
IntegrationRules * IntegRules
Definition: tmop_tools.hpp:137
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
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:212
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition: fespace.hpp:87
virtual void Mult(const Vector &b, Vector &x) const
Solve the nonlinear system with right-hand side b.
Definition: solvers.cpp:1825
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:161
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:569
virtual void ComputeAtNewPosition(const Vector &new_nodes, Vector &new_field)
Definition: tmop_tools.cpp:354
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: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
SerialAdvectorCGOper(const Vector &x_start, GridFunction &vel, FiniteElementSpace &fes, AssemblyLevel al=AssemblyLevel::LEGACY)
Definition: tmop_tools.cpp:167
void SetTempMemoryType(MemoryType mt)
Definition: tmop_tools.hpp:181
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:187
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
This method can be overloaded in derived classes to implement line search algorithms.
Definition: tmop_tools.cpp:362
Vector coefficient defined by a vector GridFunction.
Base class for solvers.
Definition: operator.hpp:648
double ComputeMinDet(const Vector &x_loc, const FiniteElementSpace &fes) const
Definition: tmop_tools.cpp:652
virtual void SetSolver(Solver &solver)
Set the linear solver for inverting the Jacobian.
Definition: tmop_tools.hpp:200
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:1198