MFEM  v4.1.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 
19 namespace mfem
20 {
21 
22 // Performs the full remap advection loop.
24 {
25 private:
26  RK4Solver ode_solver;
27  Vector nodes0;
28  Vector field0;
29 
30 public:
31  AdvectorCG() : AdaptivityEvaluator(), ode_solver(), nodes0(), field0() { }
32 
33  virtual void SetInitialField(const Vector &init_nodes,
34  const Vector &init_field);
35 
36  virtual void ComputeAtNewPosition(const Vector &new_nodes,
37  Vector &new_field);
38 };
39 
40 /// Performs a single remap advection step in serial.
42 {
43 protected:
44  const Vector &x0;
48  mutable BilinearForm M, K;
49 
50 public:
51  /** Here @a fes is the FESpace of the function that will be moved. Note
52  that Mult() moves the nodes of the mesh corresponding to @a fes. */
53  SerialAdvectorCGOper(const Vector &x_start, GridFunction &vel,
54  FiniteElementSpace &fes);
55 
56  virtual void Mult(const Vector &ind, Vector &di_dt) const;
57 };
58 
59 #ifdef MFEM_USE_MPI
60 /// Performs a single remap advection step in parallel.
62 {
63 protected:
64  const Vector &x0;
68  mutable ParBilinearForm M, K;
69 
70 public:
71  /** Here @a pfes is the ParFESpace of the function that will be moved. Note
72  that Mult() moves the nodes of the mesh corresponding to @a pfes. */
73  ParAdvectorCGOper(const Vector &x_start, GridFunction &vel,
74  ParFiniteElementSpace &pfes);
75 
76  virtual void Mult(const Vector &ind, Vector &di_dt) const;
77 };
78 #endif
79 
81 {
82 private:
83  bool parallel;
84 
85  // Quadrature points that are checked for negative Jacobians etc.
86  const IntegrationRule &ir;
87 
88  mutable DiscreteAdaptTC *discr_tc;
89 
90 public:
91 #ifdef MFEM_USE_MPI
92  TMOPNewtonSolver(MPI_Comm comm, const IntegrationRule &irule)
93  : NewtonSolver(comm), parallel(true), ir(irule), discr_tc(NULL) { }
94 #endif
96  : NewtonSolver(), parallel(false), ir(irule), discr_tc(NULL) { }
97 
98  void SetDiscreteAdaptTC(DiscreteAdaptTC *tc) { discr_tc = tc; }
99 
100  virtual double ComputeScalingFactor(const Vector &x, const Vector &b) const;
101 
102  virtual void ProcessNewState(const Vector &x) const;
103 };
104 
105 /// Allows negative Jacobians. Used for untangling.
107 {
108 private:
109  bool parallel;
110 
111  // Quadrature points that are checked for negative Jacobians etc.
112  const IntegrationRule &ir;
113 
114  mutable DiscreteAdaptTC *discr_tc;
115 
116 public:
117 #ifdef MFEM_USE_MPI
118  TMOPDescentNewtonSolver(MPI_Comm comm, const IntegrationRule &irule)
119  : NewtonSolver(comm), parallel(true), ir(irule), discr_tc(NULL) { }
120 #endif
122  : NewtonSolver(), parallel(false), ir(irule), discr_tc(NULL) { }
123 
124  virtual double ComputeScalingFactor(const Vector &x, const Vector &b) const;
125 
126  virtual void ProcessNewState(const Vector &x) const;
127 };
128 
129 void vis_tmop_metric_s(int order, TMOP_QualityMetric &qm,
130  const TargetConstructor &tc, Mesh &pmesh,
131  char *title, int position);
132 #ifdef MFEM_USE_MPI
133 void vis_tmop_metric_p(int order, TMOP_QualityMetric &qm,
134  const TargetConstructor &tc, ParMesh &pmesh,
135  char *title, int position);
136 #endif
137 
138 }
139 
140 #endif
VectorGridFunctionCoefficient u_coeff
Definition: tmop_tools.hpp:47
void vis_tmop_metric_p(int order, TMOP_QualityMetric &qm, const TargetConstructor &tc, ParMesh &pmesh, char *title, int position)
Definition: tmop_tools.cpp:473
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:27
Base abstract class for first order time dependent operators.
Definition: operator.hpp:259
TMOPNewtonSolver(MPI_Comm comm, const IntegrationRule &irule)
Definition: tmop_tools.hpp:92
Allows negative Jacobians. Used for untangling.
Definition: tmop_tools.hpp:106
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:501
ParBilinearForm M
Definition: tmop_tools.hpp:68
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:452
virtual void ComputeAtNewPosition(const Vector &new_nodes, Vector &new_field)
Definition: tmop_tools.cpp:29
ParBilinearForm K
Definition: tmop_tools.hpp:68
SerialAdvectorCGOper(const Vector &x_start, GridFunction &vel, FiniteElementSpace &fes)
Definition: tmop_tools.cpp:123
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:141
Newton's method for solving F(x)=b for a given operator F.
Definition: solvers.hpp:300
Performs a single remap advection step in serial.
Definition: tmop_tools.hpp:41
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:337
TMOPDescentNewtonSolver(MPI_Comm comm, const IntegrationRule &irule)
Definition: tmop_tools.hpp:118
The classical explicit forth-order Runge-Kutta method, RK4.
Definition: ode.hpp:145
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition: fespace.hpp:87
ParAdvectorCGOper(const Vector &x_start, GridFunction &vel, ParFiniteElementSpace &pfes)
Definition: tmop_tools.cpp:167
Performs a single remap advection step in parallel.
Definition: tmop_tools.hpp:61
TMOPDescentNewtonSolver(const IntegrationRule &irule)
Definition: tmop_tools.hpp:121
TMOPNewtonSolver(const IntegrationRule &irule)
Definition: tmop_tools.hpp:95
virtual void SetInitialField(const Vector &init_nodes, const Vector &init_field)
Definition: tmop_tools.cpp:22
VectorGridFunctionCoefficient u_coeff
Definition: tmop_tools.hpp:67
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:185
Vector data type.
Definition: vector.hpp:48
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:221
Vector coefficient defined by a vector GridFunction.
void SetDiscreteAdaptTC(DiscreteAdaptTC *tc)
Definition: tmop_tools.hpp:98
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:355