MFEM  v4.2.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-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 #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 
32  void ComputeAtNewPositionScalar(const Vector &new_nodes, Vector &new_field);
33 public:
34  AdvectorCG(double timestep_scale = 0.5)
36  ode_solver(), nodes0(), field0(), dt_scale(timestep_scale) { }
37 
38  virtual void SetInitialField(const Vector &init_nodes,
39  const Vector &init_field);
40 
41  virtual void ComputeAtNewPosition(const Vector &new_nodes,
42  Vector &new_field);
43 };
44 
45 #ifdef MFEM_USE_GSLIB
47 {
48 private:
49  Vector nodes0;
50  GridFunction field0_gf;
51  FindPointsGSLIB *finder;
52  int dim;
53 public:
54  InterpolatorFP() : finder(NULL) { }
55 
56  virtual void SetInitialField(const Vector &init_nodes,
57  const Vector &init_field);
58 
59  virtual void ComputeAtNewPosition(const Vector &new_nodes,
60  Vector &new_field);
61 
63  {
64  finder->FreeData();
65  delete finder;
66  }
67 };
68 #endif
69 
70 /// Performs a single remap advection step in serial.
72 {
73 protected:
74  const Vector &x0;
78  mutable BilinearForm M, K;
79 
80 public:
81  /** Here @a fes is the FESpace of the function that will be moved. Note
82  that Mult() moves the nodes of the mesh corresponding to @a fes. */
84  FiniteElementSpace &fes);
85 
86  virtual void Mult(const Vector &ind, Vector &di_dt) const;
87 };
88 
89 #ifdef MFEM_USE_MPI
90 /// Performs a single remap advection step in parallel.
92 {
93 protected:
94  const Vector &x0;
98  mutable ParBilinearForm M, K;
99 
100 public:
101  /** Here @a pfes is the ParFESpace of the function that will be moved. Note
102  that Mult() moves the nodes of the mesh corresponding to @a pfes. */
103  ParAdvectorCGOper(const Vector &x_start, GridFunction &vel,
104  ParFiniteElementSpace &pfes);
105 
106  virtual void Mult(const Vector &ind, Vector &di_dt) const;
107 };
108 #endif
109 
111 {
112 protected:
113  // 0 - Newton, 1 - LBFGS.
115  bool parallel;
116 
117  // Quadrature points that are checked for negative Jacobians etc.
119  // These fields are relevant for mixed meshes.
122 
124  {
125  if (IntegRules)
126  {
127  return IntegRules->Get(el.GetGeomType(), integ_order);
128  }
129  return ir;
130  }
131 
132  void UpdateDiscreteTC(const TMOP_Integrator &ti, const Vector &x_new) const;
133 
134 public:
135 #ifdef MFEM_USE_MPI
136  TMOPNewtonSolver(MPI_Comm comm, const IntegrationRule &irule, int type = 0)
137  : LBFGSSolver(comm), solver_type(type), parallel(true),
138  ir(irule), IntegRules(NULL), integ_order(-1) { }
139 #endif
140  TMOPNewtonSolver(const IntegrationRule &irule, int type = 0)
141  : LBFGSSolver(), solver_type(type), parallel(false),
142  ir(irule), IntegRules(NULL), integ_order(-1) { }
143 
144  /// Prescribe a set of integration rules; relevant for mixed meshes.
145  /** If called, this function has priority over the IntegrationRule given to
146  the constructor of the class. */
147  void SetIntegrationRules(IntegrationRules &irules, int order)
148  {
149  IntegRules = &irules;
150  integ_order = order;
151  }
152 
153  virtual double ComputeScalingFactor(const Vector &x, const Vector &b) const;
154 
155  virtual void ProcessNewState(const Vector &x) const;
156 
157  virtual void Mult(const Vector &b, Vector &x) const
158  {
159  if (solver_type == 0)
160  {
161  NewtonSolver::Mult(b, x);
162  }
163  else if (solver_type == 1)
164  {
165  LBFGSSolver::Mult(b, x);
166  }
167  else { MFEM_ABORT("Invalid type"); }
168  }
169 
170  virtual void SetSolver(Solver &solver)
171  {
172  if (solver_type == 0)
173  {
174  NewtonSolver::SetSolver(solver);
175  }
176  else if (solver_type == 1)
177  {
178  LBFGSSolver::SetSolver(solver);
179  }
180  else { MFEM_ABORT("Invalid type"); }
181  }
182  virtual void SetPreconditioner(Solver &pr) { SetSolver(pr); }
183 };
184 
185 void vis_tmop_metric_s(int order, TMOP_QualityMetric &qm,
186  const TargetConstructor &tc, Mesh &pmesh,
187  char *title, int position);
188 #ifdef MFEM_USE_MPI
189 void vis_tmop_metric_p(int order, TMOP_QualityMetric &qm,
190  const TargetConstructor &tc, ParMesh &pmesh,
191  char *title, int position);
192 #endif
193 
194 }
195 
196 #endif
Abstract class for all finite elements.
Definition: fe.hpp:235
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
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:915
void UpdateDiscreteTC(const TMOP_Integrator &ti, const Vector &x_new) const
Definition: tmop_tools.cpp:603
void SetIntegrationRules(IntegrationRules &irules, int order)
Prescribe a set of integration rules; relevant for mixed meshes.
Definition: tmop_tools.hpp:147
const IntegrationRule & ir
Definition: tmop_tools.hpp:118
Base abstract class for first order time dependent operators.
Definition: operator.hpp:267
Container class for integration rules.
Definition: intrules.hpp:309
const IntegrationRule & GetIntegrationRule(const FiniteElement &el) const
Definition: tmop_tools.hpp:123
virtual void SetSolver(Solver &solver)
Set the linear solver for inverting the Jacobian.
Definition: solvers.hpp:421
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
ParBilinearForm M
Definition: tmop_tools.hpp:98
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:312
ParBilinearForm K
Definition: tmop_tools.hpp:98
AdvectorCG(double timestep_scale=0.5)
Definition: tmop_tools.hpp:34
SerialAdvectorCGOper(const Vector &x_start, GridFunction &vel, FiniteElementSpace &fes)
Definition: tmop_tools.cpp:165
double b
Definition: lissajous.cpp:42
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
Performs a single remap advection step in serial.
Definition: tmop_tools.hpp:71
TMOPNewtonSolver(const IntegrationRule &irule, int type=0)
Definition: tmop_tools.hpp:140
virtual void Mult(const Vector &b, Vector &x) const
Solve the nonlinear system with right-hand side b.
Definition: solvers.cpp:1563
IntegrationRules * IntegRules
Definition: tmop_tools.hpp:120
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
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:182
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:1644
ParAdvectorCGOper(const Vector &x_start, GridFunction &vel, ParFiniteElementSpace &pfes)
Definition: tmop_tools.cpp:209
TMOPNewtonSolver(MPI_Comm comm, const IntegrationRule &irule, int type=0)
Definition: tmop_tools.hpp:136
Performs a single remap advection step in parallel.
Definition: tmop_tools.hpp:91
virtual void SetSolver(Solver &solver)
Set the linear solver for inverting the Jacobian.
Definition: solvers.hpp:461
virtual void ComputeAtNewPosition(const Vector &new_nodes, Vector &new_field)
Definition: tmop_tools.cpp:303
A "square matrix" operator for the associated FE space and BLFIntegrators The sum of all the BLFInteg...
virtual void SetInitialField(const Vector &init_nodes, const Vector &init_field)
Definition: tmop_tools.cpp:22
VectorGridFunctionCoefficient u_coeff
Definition: tmop_tools.hpp:97
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:227
virtual void Mult(const Vector &b, Vector &x) const
Solve the nonlinear system with right-hand side b.
Definition: tmop_tools.hpp:157
FindPointsGSLIB can robustly evaluate a GridFunction on an arbitrary collection of points...
Definition: gslib.hpp:45
Vector data type.
Definition: vector.hpp:51
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
Vector coefficient defined by a vector GridFunction.
Base class for solvers.
Definition: operator.hpp:634
virtual void SetSolver(Solver &solver)
Set the linear solver for inverting the Jacobian.
Definition: tmop_tools.hpp:170
virtual void SetInitialField(const Vector &init_nodes, const Vector &init_field)
Definition: tmop_tools.cpp:264
A TMOP integrator class based on any given TMOP_QualityMetric and TargetConstructor.
Definition: tmop.hpp:884