MFEM  v4.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
eltrans.hpp
Go to the documentation of this file.
1 // Copyright (c) 2010, Lawrence Livermore National Security, LLC. Produced at
2 // the Lawrence Livermore National Laboratory. LLNL-CODE-443211. All Rights
3 // reserved. See file COPYRIGHT for details.
4 //
5 // This file is part of the MFEM library. For more information and source code
6 // availability see http://mfem.org.
7 //
8 // MFEM is free software; you can redistribute it and/or modify it under the
9 // terms of the GNU Lesser General Public License (as published by the Free
10 // Software Foundation) version 2.1 dated February 1999.
11 
12 #ifndef MFEM_ELEMENTTRANSFORM
13 #define MFEM_ELEMENTTRANSFORM
14 
15 #include "../config/config.hpp"
16 #include "../linalg/linalg.hpp"
17 #include "intrules.hpp"
18 #include "fe.hpp"
19 
20 namespace mfem
21 {
22 
24 {
25 protected:
28  double Wght;
29  int EvalState;
31  {
36  };
38  int space_dim;
39 
40  // Evaluate the Jacobian of the transformation at the IntPoint and store it
41  // in dFdx.
42  virtual const DenseMatrix &EvalJacobian() = 0;
43 
44  double EvalWeight();
45  const DenseMatrix &EvalAdjugateJ();
46  const DenseMatrix &EvalInverseJ();
47 
48 public:
50 
52 
53  void SetIntPoint(const IntegrationPoint *ip)
54  { IntPoint = ip; EvalState = 0; }
55  const IntegrationPoint &GetIntPoint() { return *IntPoint; }
56 
57  virtual void Transform(const IntegrationPoint &, Vector &) = 0;
58  virtual void Transform(const IntegrationRule &, DenseMatrix &) = 0;
59 
60  /// Transform columns of 'matrix', store result in 'result'.
61  virtual void Transform(const DenseMatrix &matrix, DenseMatrix &result) = 0;
62 
63  /** @brief Return the Jacobian matrix of the transformation at the currently
64  set IntegrationPoint, using the method SetIntPoint(). */
65  /** The dimensions of the Jacobian matrix are physical-space-dim by
66  reference-space-dim. The first column contains the x derivatives of the
67  transformation, the second -- the y derivatives, etc. */
69  { return (EvalState & JACOBIAN_MASK) ? dFdx : EvalJacobian(); }
70 
71  double Weight() { return (EvalState & WEIGHT_MASK) ? Wght : EvalWeight(); }
72 
74  { return (EvalState & ADJUGATE_MASK) ? adjJ : EvalAdjugateJ(); }
75 
77  { return (EvalState & INVERSE_MASK) ? invJ : EvalInverseJ(); }
78 
79  virtual int Order() = 0;
80  virtual int OrderJ() = 0;
81  virtual int OrderW() = 0;
82  /// Order of adj(J)^t.grad(fi)
83  virtual int OrderGrad(const FiniteElement *fe) = 0;
84 
85  /// Return the Geometry::Type of the reference element.
86  Geometry::Type GetGeometryType() const { return geom; }
87 
88  /// Return the dimension of the reference element.
89  int GetDimension() const { return Geometry::Dimension[geom]; }
90 
91  /// Get the dimension of the target (physical) space.
92  /** We support 2D meshes embedded in 3D; in this case the function will
93  return "3". */
94  int GetSpaceDim() const { return space_dim; }
95 
96  /** @brief Transform a point @a pt from physical space to a point @a ip in
97  reference space. */
98  /** Attempt to find the IntegrationPoint that is transformed into the given
99  point in physical space. If the inversion fails a non-zero value is
100  returned. This method is not 100 percent reliable for non-linear
101  transformations. */
102  virtual int TransformBack(const Vector &pt, IntegrationPoint &ip) = 0;
103 
105 };
106 
107 
108 /// The inverse transformation of a given ElementTransformation.
110 {
111 public:
112  /// Algorithms for selecting an initial guess.
114  {
115  Center = 0, ///< Use the center of the reference element.
116  ClosestPhysNode = 1, /**<
117  Use the point returned by FindClosestPhysPoint() from a reference-space
118  grid of type and size controlled by SetInitGuessPointsType() and
119  SetInitGuessRelOrder(), respectively. */
120  ClosestRefNode = 2, /**<
121  Use the point returned by FindClosestRefPoint() from a reference-space
122  grid of type and size controlled by SetInitGuessPointsType() and
123  SetInitGuessRelOrder(), respectively. */
124  GivenPoint = 3 ///< Use a specific point, set with SetInitialGuess().
125  };
126 
127  /// Solution strategy.
129  {
130  Newton = 0, /**<
131  Use Newton's algorithm, without restricting the reference-space points
132  (iterates) to the reference element. */
134  Use Newton's algorithm, restricting the reference-space points to the
135  reference element by scaling back the Newton increments, i.e.
136  projecting new iterates, x_new, lying outside the element, to the
137  intersection of the line segment [x_old, x_new] with the boundary. */
139  Use Newton's algorithm, restricting the reference-space points to the
140  reference element by projecting new iterates, x_new, lying outside the
141  element, to the point on the boundary closest (in reference-space) to
142  x_new. */
143  };
144 
145  /// Values returned by Transform().
147  {
148  Inside = 0, ///< The point is inside the element
149  Outside = 1, ///< The point is _probably_ outside the element
150  Unknown = 2 ///< The algorithm failed to determine where the point is
151  };
152 
153 protected:
154  // Pointer to the forward transformation. Not owned.
156 
157  // Parameters of the inversion algorithms:
159  int init_guess_type; // algorithm to use
160  int qpts_type; // Quadrature1D type for the initial guess type
161  int rel_qpts_order; // num_1D_qpts = max(trans_order+rel_qpts_order,0)+1
162  int solver_type; // solution strategy to use
163  int max_iter; // max. number of Newton iterations
164  double ref_tol; // reference space tolerance
165  double phys_rtol; // physical space tolerance (relative)
166  double ip_tol; // tolerance for checking if a point is inside the ref. elem.
168 
169  void NewtonPrint(int mode, double val);
170  void NewtonPrintPoint(const char *prefix, const Vector &pt,
171  const char *suffix);
172  int NewtonSolve(const Vector &pt, IntegrationPoint &ip);
173 
174 public:
175  /// Construct the InverseElementTransformation with default parameters.
176  /** Some practical considerations regarding the choice of initial guess type
177  and solver type:
178  1. The combination of #Center and #NewtonSegmentProject provides the
179  fastest way to estimate if a point lies inside an element, assuming
180  that most queried elements are not very deformed, e.g. if most
181  elements are convex.
182  2. [Default] The combination of #Center and #NewtonElementProject
183  provides a somewhat slower alternative to 1 with the benefit of being
184  more reliable in the case when the query point is inside the element
185  but potentially slower in the case when the query point is outside the
186  element.
187  3. The combination of #ClosestPhysNode and #NewtonElementProject is
188  slower than 1 and 2 but much more reliable, especially in the case of
189  highly distorted elements which do not have very high aspect ratios.
190  4. The combination of #ClosestRefNode and #NewtonElementProject should
191  generally be the most reliable, coming at a bit higher computational
192  cost than 3 while performing better on distorted meshes with elements
193  having high aspect ratios.
194  @note None of these choices provide a guarantee that if a query point is
195  inside the element then it will be found. The only guarantee is that if
196  the Transform() method returns #Inside then the point lies inside the
197  element up to one of the specified physical- or reference-space
198  tolerances. */
200  : T(Trans),
201  ip0(NULL),
203  qpts_type(Quadrature1D::OpenHalfUniform),
204  rel_qpts_order(-1),
206  max_iter(16),
207  ref_tol(1e-15),
208  phys_rtol(1e-15),
209  ip_tol(1e-8),
210  print_level(-1)
211  { }
212 
214 
215  /// Set a new forward ElementTransformation, @a Trans.
217 
218  /** @brief Choose how the initial guesses for subsequent calls to Transform()
219  will be selected. */
221 
222  /** @brief Set the initial guess for subsequent calls to Transform(),
223  switching to the #GivenPoint #InitGuessType at the same time. */
224  void SetInitialGuess(const IntegrationPoint &init_ip)
225  { ip0 = &init_ip; SetInitialGuessType(GivenPoint); }
226 
227  /// Set the Quadrature1D type used for the `Closest*` initial guess types.
228  void SetInitGuessPointsType(int q_type) { qpts_type = q_type; }
229 
230  /// Set the relative order used for the `Closest*` initial guess types.
231  /** The number of points in each spatial direction is given by the formula
232  max(trans_order+order,0)+1, where trans_order is the order of the current
233  ElementTransformation. */
234  void SetInitGuessRelOrder(int order) { rel_qpts_order = order; }
235 
236  /** @brief Specify which algorithm to use for solving the transformation
237  equation, i.e. when calling the Transform() method. */
238  void SetSolverType(SolverType stype) { solver_type = stype; }
239 
240  /// Set the maximum number of iterations when solving for a reference point.
241  void SetMaxIter(int max_it) { max_iter = max_it; }
242 
243  /// Set the reference-space convergence tolerance.
244  void SetReferenceTol(double ref_sp_tol) { ref_tol = ref_sp_tol; }
245 
246  /// Set the relative physical-space convergence tolerance.
247  void SetPhysicalRelTol(double phys_rel_tol) { phys_rtol = phys_rel_tol; }
248 
249  /** @brief Set the tolerance used to determine if a point lies inside or
250  outside of the reference element. */
251  /** This tolerance is used only with the pure #Newton solver. */
252  void SetElementTol(double el_tol) { ip_tol = el_tol; }
253 
254  /// Set the desired print level, useful for debugging.
255  /** The valid options are: -1 - never print (default); 0 - print only errors;
256  1 - print the first and last last iterations; 2 - print every iteration;
257  and 3 - print every iteration including point coordinates. */
258  void SetPrintLevel(int pr_level) { print_level = pr_level; }
259 
260  /** @brief Find the IntegrationPoint mapped closest to @a pt. */
261  /** This function uses the given IntegrationRule, @a ir, maps its points to
262  physical space and finds the one that is closest to the point @a pt.
263 
264  @param pt The query point.
265  @param ir The IntegrationRule, i.e. the set of reference points to map
266  to physical space and check.
267  @return The index of the IntegrationPoint in @a ir whose mapped point is
268  closest to @a pt.
269  @see FindClosestRefPoint(). */
270  int FindClosestPhysPoint(const Vector& pt, const IntegrationRule &ir);
271 
272  /** @brief Find the IntegrationPoint mapped closest to @a pt, using a norm
273  that approximates the (unknown) distance in reference coordinates. */
274  /** @see FindClosestPhysPoint(). */
275  int FindClosestRefPoint(const Vector& pt, const IntegrationRule &ir);
276 
277  /** @brief Given a point, @a pt, in physical space, find its reference
278  coordinates, @a ip.
279 
280  @returns A value of type #TransformResult. */
281  virtual int Transform(const Vector &pt, IntegrationPoint &ip);
282 };
283 
284 
286 {
287 private:
288  DenseMatrix dshape;
289  Vector shape;
290 
291  const FiniteElement *FElem;
292  DenseMatrix PointMat; // dim x dof
293 
294  // Evaluate the Jacobian of the transformation at the IntPoint and store it
295  // in dFdx.
296  virtual const DenseMatrix &EvalJacobian();
297 
298 public:
299  void SetFE(const FiniteElement *FE) { FElem = FE; geom = FE->GetGeomType(); }
300  const FiniteElement* GetFE() const { return FElem; }
301 
302  /** @brief Read and write access to the underlying point matrix describing
303  the transformation. */
304  /** The dimensions of the matrix are space-dim x dof. The transformation is
305  defined as
306 
307  x=F(xh)=P.phi(xh),
308 
309  where xh (x hat) is the reference point, x is the corresponding physical
310  point, P is the point matrix, and phi(xh) is the column-vector of all
311  basis functions evaluated at xh. The columns of P represent the control
312  points in physical space defining the transformation. */
313  DenseMatrix &GetPointMat() { return PointMat; }
314  void FinalizeTransformation() { space_dim = PointMat.Height(); }
315 
317 
318  virtual void Transform(const IntegrationPoint &, Vector &);
319  virtual void Transform(const IntegrationRule &, DenseMatrix &);
320  virtual void Transform(const DenseMatrix &matrix, DenseMatrix &result);
321 
322  virtual int Order() { return FElem->GetOrder(); }
323  virtual int OrderJ();
324  virtual int OrderW();
325  virtual int OrderGrad(const FiniteElement *fe);
326 
327  virtual int TransformBack(const Vector & v, IntegrationPoint & ip)
328  {
329  InverseElementTransformation inv_tr(this);
330  return inv_tr.Transform(v, ip);
331  }
332 
334 };
335 
337 {
338 public:
340  void Transform (const IntegrationPoint &, IntegrationPoint &);
341  void Transform (const IntegrationRule &, IntegrationRule &);
342 };
343 
345 {
346 public:
350 };
351 
352 /* Elem1(Loc1(x)) = Face(x) = Elem2(Loc2(x))
353 
354 
355  Physical Space
356 
357  *--------* ^ *--------*
358  Elem1No / / \ / \ / \ \ Elem2No
359  / / \ / \ / \ \
360  / / n \ / \ / \ \
361  *--------* ==> * ( ) * *--------*
362  \ \ / \ / \ / /
363  \ \ / \ / \ / /
364  \ \ / \ / \ / /
365  *--------* v *--------*
366 
367  ^ ^
368  | ^ |
369  Elem1 | | | Elem2
370  | | Face |
371  |
372  *--------* *--------*
373  / /| / /|
374  1 *--------* | 1 *--------* 1 *--------* |
375  | | | Loc1 | | Loc2 | | |
376  | | * <----- | x | -----> | | *
377  | |/ | | | |/
378  *--------* *--------* *--------*
379  0 1 0 1 0 1
380 
381  Reference Space
382 */
383 
384 }
385 
386 #endif
Abstract class for Finite Elements.
Definition: fe.hpp:229
virtual ~ElementTransformation()
Definition: eltrans.hpp:104
void SetPrintLevel(int pr_level)
Set the desired print level, useful for debugging.
Definition: eltrans.hpp:258
int NewtonSolve(const Vector &pt, IntegrationPoint &ip)
Definition: eltrans.cpp:159
const DenseMatrix & AdjugateJacobian()
Definition: eltrans.hpp:73
ElementTransformation * Face
Definition: eltrans.hpp:348
Class for an integration rule - an Array of IntegrationPoint.
Definition: intrules.hpp:85
DenseMatrix & GetPointMat()
Read and write access to the underlying point matrix describing the transformation.
Definition: eltrans.hpp:313
Use a specific point, set with SetInitialGuess().
Definition: eltrans.hpp:124
void SetPhysicalRelTol(double phys_rel_tol)
Set the relative physical-space convergence tolerance.
Definition: eltrans.hpp:247
void SetInitialGuess(const IntegrationPoint &init_ip)
Set the initial guess for subsequent calls to Transform(), switching to the GivenPoint InitGuessType ...
Definition: eltrans.hpp:224
Geometry::Type GetGeometryType() const
Return the Geometry::Type of the reference element.
Definition: eltrans.hpp:86
const IntegrationPoint * IntPoint
Definition: eltrans.hpp:26
ElementTransformation * T
Definition: eltrans.hpp:155
Use the center of the reference element.
Definition: eltrans.hpp:115
virtual int Transform(const Vector &pt, IntegrationPoint &ip)
Given a point, pt, in physical space, find its reference coordinates, ip.
Definition: eltrans.cpp:323
virtual int TransformBack(const Vector &v, IntegrationPoint &ip)
Transform a point pt from physical space to a point ip in reference space.
Definition: eltrans.hpp:327
int GetDimension() const
Return the dimension of the reference element.
Definition: eltrans.hpp:89
void SetIntPoint(const IntegrationPoint *ip)
Definition: eltrans.hpp:53
void NewtonPrint(int mode, double val)
Definition: eltrans.cpp:116
int GetOrder() const
Returns the order of the finite element. In the case of anisotropic orders, returns the maximum order...
Definition: fe.hpp:315
Data type dense matrix using column-major storage.
Definition: densemat.hpp:23
void SetInitialGuessType(InitGuessType itype)
Choose how the initial guesses for subsequent calls to Transform() will be selected.
Definition: eltrans.hpp:220
IntegrationPointTransformation Loc2
Definition: eltrans.hpp:349
const DenseMatrix & EvalAdjugateJ()
Definition: eltrans.cpp:34
const DenseMatrix & InverseJacobian()
Definition: eltrans.hpp:76
const FiniteElement * GetFE() const
Definition: eltrans.hpp:300
A class container for 1D quadrature type constants.
Definition: intrules.hpp:278
virtual void Transform(const IntegrationPoint &, Vector &)
Definition: eltrans.cpp:463
const IntegrationPoint & GetIntPoint()
Definition: eltrans.hpp:55
The inverse transformation of a given ElementTransformation.
Definition: eltrans.hpp:109
ElementTransformation * Elem2
Definition: eltrans.hpp:348
void SetSolverType(SolverType stype)
Specify which algorithm to use for solving the transformation equation, i.e. when calling the Transfo...
Definition: eltrans.hpp:238
void NewtonPrintPoint(const char *prefix, const Vector &pt, const char *suffix)
Definition: eltrans.cpp:145
void SetTransformation(ElementTransformation &Trans)
Set a new forward ElementTransformation, Trans.
Definition: eltrans.hpp:216
Geometry::Type GetGeomType() const
Returns the Geometry::Type of the reference element.
Definition: fe.hpp:308
int FindClosestRefPoint(const Vector &pt, const IntegrationRule &ir)
Find the IntegrationPoint mapped closest to pt, using a norm that approximates the (unknown) distance...
Definition: eltrans.cpp:84
The point is inside the element.
Definition: eltrans.hpp:148
static const int Dimension[NumGeom]
Definition: geom.hpp:46
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
Definition: operator.hpp:36
IntegrationPointTransformation Loc1
Definition: eltrans.hpp:349
The point is probably outside the element.
Definition: eltrans.hpp:149
virtual int OrderGrad(const FiniteElement *fe)
Order of adj(J)^t.grad(fi)
Definition: eltrans.cpp:444
const DenseMatrix & Jacobian()
Return the Jacobian matrix of the transformation at the currently set IntegrationPoint, using the method SetIntPoint().
Definition: eltrans.hpp:68
int GetSpaceDim() const
Get the dimension of the target (physical) space.
Definition: eltrans.hpp:94
void SetReferenceTol(double ref_sp_tol)
Set the reference-space convergence tolerance.
Definition: eltrans.hpp:244
void SetInitGuessPointsType(int q_type)
Set the Quadrature1D type used for the Closest* initial guess types.
Definition: eltrans.hpp:228
InverseElementTransformation(ElementTransformation *Trans=NULL)
Construct the InverseElementTransformation with default parameters.
Definition: eltrans.hpp:199
void SetElementTol(double el_tol)
Set the tolerance used to determine if a point lies inside or outside of the reference element...
Definition: eltrans.hpp:252
const DenseMatrix & EvalInverseJ()
Definition: eltrans.cpp:44
void SetIdentityTransformation(Geometry::Type GeomType)
Definition: eltrans.cpp:370
void Transform(const IntegrationPoint &, IntegrationPoint &)
Definition: eltrans.cpp:517
Class for integration point with weight.
Definition: intrules.hpp:25
The algorithm failed to determine where the point is.
Definition: eltrans.hpp:150
SolverType
Solution strategy.
Definition: eltrans.hpp:128
virtual int OrderGrad(const FiniteElement *fe)=0
Order of adj(J)^t.grad(fi)
TransformResult
Values returned by Transform().
Definition: eltrans.hpp:146
int FindClosestPhysPoint(const Vector &pt, const IntegrationRule &ir)
Find the IntegrationPoint mapped closest to pt.
Definition: eltrans.cpp:57
void SetFE(const FiniteElement *FE)
Definition: eltrans.hpp:299
IsoparametricTransformation Transf
Definition: eltrans.hpp:339
ElementTransformation * Elem1
Definition: eltrans.hpp:348
virtual const DenseMatrix & EvalJacobian()=0
InitGuessType
Algorithms for selecting an initial guess.
Definition: eltrans.hpp:113
Vector data type.
Definition: vector.hpp:48
virtual void Transform(const IntegrationPoint &, Vector &)=0
const IntegrationPoint * ip0
Definition: eltrans.hpp:158
void SetMaxIter(int max_it)
Set the maximum number of iterations when solving for a reference point.
Definition: eltrans.hpp:241
void SetInitGuessRelOrder(int order)
Set the relative order used for the Closest* initial guess types.
Definition: eltrans.hpp:234
virtual int TransformBack(const Vector &pt, IntegrationPoint &ip)=0
Transform a point pt from physical space to a point ip in reference space.