MFEM v4.7.0
Finite element discretization library
Loading...
Searching...
No Matches
eltrans.hpp
Go to the documentation of this file.
1// Copyright (c) 2010-2024, 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_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
20namespace mfem
21{
22
24{
25protected:
41
42 /** @brief Evaluate the Jacobian of the transformation at the IntPoint and
43 store it in dFdx. */
44 virtual const DenseMatrix &EvalJacobian() = 0;
45
46 /** @brief Evaluate the Hessian of the transformation at the IntPoint and
47 store it in d2Fdx2. */
48 virtual const DenseMatrix &EvalHessian() = 0;
49
54
55public:
56
57 /** This enumeration declares the values stored in
58 ElementTransformation::ElementType and indicates which group of objects
59 the index stored in ElementTransformation::ElementNo refers:
60
61 | ElementType | Range of ElementNo
62 +-------------+-------------------------
63 | ELEMENT | [0, Mesh::GetNE() )
64 | BDR_ELEMENT | [0, Mesh::GetNBE() )
65 | EDGE | [0, Mesh::GetNEdges() )
66 | FACE | [0, Mesh::GetNFaces() )
67 | BDR_FACE | [0, Mesh::GetNBE() )
68 */
69 enum
70 {
73 EDGE = 3,
74 FACE = 4,
75 BDR_FACE = 5
76 };
77
79
80 /// The Mesh object containing the element.
81 /** If the element transformation belongs to a mesh, this will point to the
82 containing Mesh object. ElementNo will be the number of the element in
83 this Mesh. This will be NULL if the element does not belong to a mesh. */
84 const Mesh *mesh;
85
87
88 /** @brief Force the reevaluation of the Jacobian in the next call. */
89 void Reset() { EvalState = 0; }
90
91 /** @brief Set the integration point @a ip that weights and Jacobians will
92 be evaluated at. */
94 { IntPoint = ip; EvalState = 0; }
95
96 /** @brief Get a const reference to the currently set integration point. This
97 will return NULL if no integration point is set. */
99
100 /** @brief Transform integration point from reference coordinates to
101 physical coordinates and store them in the vector. */
102 virtual void Transform(const IntegrationPoint &, Vector &) = 0;
103
104 /** @brief Transform all the integration points from the integration rule
105 from reference coordinates to physical
106 coordinates and store them as column vectors in the matrix. */
107 virtual void Transform(const IntegrationRule &, DenseMatrix &) = 0;
108
109 /** @brief Transform all the integration points from the column vectors
110 of @a matrix from reference coordinates to physical
111 coordinates and store them as column vectors in @a result. */
112 virtual void Transform(const DenseMatrix &matrix, DenseMatrix &result) = 0;
113
114 /** @brief Return the Jacobian matrix of the transformation at the currently
115 set IntegrationPoint, using the method SetIntPoint(). */
116 /** The dimensions of the Jacobian matrix are physical-space-dim by
117 reference-space-dim. The first column contains the x derivatives of the
118 transformation, the second -- the y derivatives, etc. */
120 { return (EvalState & JACOBIAN_MASK) ? dFdx : EvalJacobian(); }
121
122
123 /** @brief Return the Hessian matrix of the transformation at the currently
124 set IntegrationPoint, using the method SetIntPoint(). */
126 { return (EvalState & HESSIAN_MASK) ? d2Fdx2 : EvalHessian(); }
127
128 /** @brief Return the weight of the Jacobian matrix of the transformation
129 at the currently set IntegrationPoint.
130 The Weight evaluates to $ \sqrt{\lvert J^T J \rvert} $. */
132
133 /** @brief Return the adjugate of the Jacobian matrix of the transformation
134 at the currently set IntegrationPoint. */
137
138 /** @brief Return the transpose of the adjugate of the Jacobian matrix of
139 the transformation at the currently set IntegrationPoint. */
142
143 /** @brief Return the inverse of the Jacobian matrix of the transformation
144 at the currently set IntegrationPoint. */
147
148 /// Return the order of the current element we are using for the transformation.
149 virtual int Order() const = 0;
150
151 /// Return the order of the elements of the Jacobian of the transformation.
152 virtual int OrderJ() const = 0;
153
154 /** @brief Return the order of the determinant of the Jacobian (weight)
155 of the transformation. */
156 virtual int OrderW() const = 0;
157
158 /// Return the order of $ adj(J)^T \nabla fi $
159 virtual int OrderGrad(const FiniteElement *fe) const = 0;
160
161 /// Return the Geometry::Type of the reference element.
163
164 /// Return the topological dimension of the reference element.
165 int GetDimension() const { return Geometry::Dimension[geom]; }
166
167 /// Get the dimension of the target (physical) space.
168 /** We support 2D meshes embedded in 3D; in this case the function will
169 return "3". */
170 virtual int GetSpaceDim() const = 0;
171
172 /** @brief Transform a point @a pt from physical space to a point @a ip in
173 reference space and optionally can set a solver tolerance using @a phys_tol. */
174 /** Attempt to find the IntegrationPoint that is transformed into the given
175 point in physical space. If the inversion fails a non-zero value is
176 returned. This method is not 100 percent reliable for non-linear
177 transformations. */
178 virtual int TransformBack(const Vector &pt, IntegrationPoint &ip,
179 const real_t phys_tol = 1e-15) = 0;
180
182};
183
184
185/// The inverse transformation of a given ElementTransformation.
187{
188public:
189 /// Algorithms for selecting an initial guess.
191 {
192 Center = 0, ///< Use the center of the reference element.
193 ClosestPhysNode = 1, /**<
194 Use the point returned by FindClosestPhysPoint() from a reference-space
195 grid of type and size controlled by SetInitGuessPointsType() and
196 SetInitGuessRelOrder(), respectively. */
197 ClosestRefNode = 2, /**<
198 Use the point returned by FindClosestRefPoint() from a reference-space
199 grid of type and size controlled by SetInitGuessPointsType() and
200 SetInitGuessRelOrder(), respectively. */
201 GivenPoint = 3 ///< Use a specific point, set with SetInitialGuess().
202 };
203
204 /// Solution strategy.
206 {
207 Newton = 0, /**<
208 Use Newton's algorithm, without restricting the reference-space points
209 (iterates) to the reference element. */
210 NewtonSegmentProject = 1, /**<
211 Use Newton's algorithm, restricting the reference-space points to the
212 reference element by scaling back the Newton increments, i.e.
213 projecting new iterates, x_new, lying outside the element, to the
214 intersection of the line segment [x_old, x_new] with the boundary. */
215 NewtonElementProject = 2 /**<
216 Use Newton's algorithm, restricting the reference-space points to the
217 reference element by projecting new iterates, x_new, lying outside the
218 element, to the point on the boundary closest (in reference-space) to
219 x_new. */
220 };
221
222 /// Values returned by Transform().
224 {
225 Inside = 0, ///< The point is inside the element
226 Outside = 1, ///< The point is _probably_ outside the element
227 Unknown = 2 ///< The algorithm failed to determine where the point is
228 };
229
230protected:
231 // Pointer to the forward transformation. Not owned.
233
234 // Parameters of the inversion algorithms:
236 int init_guess_type; // algorithm to use
237 GeometryRefiner refiner; // geometry refiner for initial guess
238 int rel_qpts_order; // num_1D_qpts = max(trans_order+rel_qpts_order,0)+1
239 int solver_type; // solution strategy to use
240 int max_iter; // max. number of Newton iterations
241 real_t ref_tol; // reference space tolerance
242 real_t phys_rtol; // physical space tolerance (relative)
243 real_t ip_tol; // tolerance for checking if a point is inside the ref. elem.
245
246 void NewtonPrint(int mode, real_t val);
247 void NewtonPrintPoint(const char *prefix, const Vector &pt,
248 const char *suffix);
249 int NewtonSolve(const Vector &pt, IntegrationPoint &ip);
250
251public:
252 /// Construct the InverseElementTransformation with default parameters.
253 /** Some practical considerations regarding the choice of initial guess type
254 and solver type:
255 1. The combination of #Center and #NewtonSegmentProject provides the
256 fastest way to estimate if a point lies inside an element, assuming
257 that most queried elements are not very deformed, e.g. if most
258 elements are convex.
259 2. [Default] The combination of #Center and #NewtonElementProject
260 provides a somewhat slower alternative to 1 with the benefit of being
261 more reliable in the case when the query point is inside the element
262 but potentially slower in the case when the query point is outside the
263 element.
264 3. The combination of #ClosestPhysNode and #NewtonElementProject is
265 slower than 1 and 2 but much more reliable, especially in the case of
266 highly distorted elements which do not have very high aspect ratios.
267 4. The combination of #ClosestRefNode and #NewtonElementProject should
268 generally be the most reliable, coming at a bit higher computational
269 cost than 3 while performing better on distorted meshes with elements
270 having high aspect ratios.
271 @note None of these choices provide a guarantee that if a query point is
272 inside the element then it will be found. The only guarantee is that if
273 the Transform() method returns #Inside then the point lies inside the
274 element up to one of the specified physical- or reference-space
275 tolerances. */
277 : T(Trans),
278 ip0(NULL),
280 refiner(Quadrature1D::OpenHalfUniform),
281 rel_qpts_order(-1),
283 max_iter(16),
284 ref_tol(1e-15),
285 phys_rtol(1e-15),
286 ip_tol(1e-8),
287 print_level(-1)
288 { }
289
291
292 /// Set a new forward ElementTransformation, @a Trans.
293 void SetTransformation(ElementTransformation &Trans) { T = &Trans; }
294
295 /** @brief Choose how the initial guesses for subsequent calls to Transform()
296 will be selected. */
298
299 /** @brief Set the initial guess for subsequent calls to Transform(),
300 switching to the #GivenPoint #InitGuessType at the same time. */
302 { ip0 = &init_ip; SetInitialGuessType(GivenPoint); }
303
304 /// Set the Quadrature1D type used for the `Closest*` initial guess types.
305 void SetInitGuessPointsType(int q_type) { refiner.SetType(q_type); }
306
307 /// Set the relative order used for the `Closest*` initial guess types.
308 /** The number of points in each spatial direction is given by the formula
309 max(trans_order+order,0)+1, where trans_order is the order of the current
310 ElementTransformation. */
311 void SetInitGuessRelOrder(int order) { rel_qpts_order = order; }
312
313 /** @brief Specify which algorithm to use for solving the transformation
314 equation, i.e. when calling the Transform() method. */
315 void SetSolverType(SolverType stype) { solver_type = stype; }
316
317 /// Set the maximum number of iterations when solving for a reference point.
318 void SetMaxIter(int max_it) { max_iter = max_it; }
319
320 /// Set the reference-space convergence tolerance.
321 void SetReferenceTol(real_t ref_sp_tol) { ref_tol = ref_sp_tol; }
322
323 /// Set the relative physical-space convergence tolerance.
324 void SetPhysicalRelTol(real_t phys_rel_tol) { phys_rtol = phys_rel_tol; }
325
326 /** @brief Set the tolerance used to determine if a point lies inside or
327 outside of the reference element. */
328 /** This tolerance is used only with the pure #Newton solver. */
329 void SetElementTol(real_t el_tol) { ip_tol = el_tol; }
330
331 /// Set the desired print level, useful for debugging.
332 /** The valid options are: -1 - never print (default); 0 - print only errors;
333 1 - print the first and last iterations; 2 - print every iteration;
334 and 3 - print every iteration including point coordinates. */
335 void SetPrintLevel(int pr_level) { print_level = pr_level; }
336
337 /** @brief Find the IntegrationPoint mapped closest to @a pt. */
338 /** This function uses the given IntegrationRule, @a ir, maps its points to
339 physical space and finds the one that is closest to the point @a pt.
340
341 @param pt The query point.
342 @param ir The IntegrationRule, i.e. the set of reference points to map
343 to physical space and check.
344 @return The index of the IntegrationPoint in @a ir whose mapped point is
345 closest to @a pt.
346 @see FindClosestRefPoint(). */
347 int FindClosestPhysPoint(const Vector& pt, const IntegrationRule &ir);
348
349 /** @brief Find the IntegrationPoint mapped closest to @a pt, using a norm
350 that approximates the (unknown) distance in reference coordinates. */
351 /** @see FindClosestPhysPoint(). */
352 int FindClosestRefPoint(const Vector& pt, const IntegrationRule &ir);
353
354 /** @brief Given a point, @a pt, in physical space, find its reference
355 coordinates, @a ip.
356
357 @returns A value of type #TransformResult. */
358 virtual int Transform(const Vector &pt, IntegrationPoint &ip);
359};
360
361/// A standard isoparametric element transformation
363{
364private:
365 DenseMatrix dshape, d2shape;
366 Vector shape;
367
368 const FiniteElement *FElem;
369 DenseMatrix PointMat; // dim x dof
370
371 /** @brief Evaluate the Jacobian of the transformation at the IntPoint and
372 store it in dFdx. */
373 virtual const DenseMatrix &EvalJacobian();
374 // Evaluate the Hessian of the transformation at the IntPoint and store it
375 // in d2Fdx2.
376 virtual const DenseMatrix &EvalHessian();
377
378public:
380
381 /// Set the element that will be used to compute the transformations
382 void SetFE(const FiniteElement *FE)
383 {
384 MFEM_ASSERT(FE != NULL, "Must provide a valid FiniteElement object!");
385 EvalState = (FE != FElem) ? 0 : EvalState;
386 FElem = FE; geom = FE->GetGeomType();
387 }
388
389 /// Get the current element used to compute the transformations
390 const FiniteElement* GetFE() const { return FElem; }
391
392 /// @brief Set the underlying point matrix describing the transformation.
393 /** The dimensions of the matrix are space-dim x dof. The transformation is
394 defined as
395 $ x = F( \hat x ) = P \phi( \hat x ) $
396
397 where $ \hat x $ is the reference point, @a x is the corresponding
398 physical point, @a P is the point matrix, and $ \phi( \hat x ) $ is
399 the column-vector of all basis functions evaluated at $ \hat x $ .
400 The columns of @a P represent the control points in physical space
401 defining the transformation. */
402 void SetPointMat(const DenseMatrix &pm) { PointMat = pm; EvalState = 0; }
403
404 /// Return the stored point matrix.
405 const DenseMatrix &GetPointMat() const { return PointMat; }
406
407 /// @brief Write access to the stored point matrix. Use with caution.
408 /** If the point matrix is altered using this member function the Reset
409 function should also be called to force the reevaluation of the
410 Jacobian, etc.. */
411 DenseMatrix &GetPointMat() { return PointMat; }
412
413 /// Set the FiniteElement Geometry for the reference elements being used.
415
416 /** @brief Transform integration point from reference coordinates to
417 physical coordinates and store them in the vector. */
418 virtual void Transform(const IntegrationPoint &, Vector &);
419
420 /** @brief Transform all the integration points from the integration rule
421 from reference coordinates to physical
422 coordinates and store them as column vectors in the matrix. */
423 virtual void Transform(const IntegrationRule &, DenseMatrix &);
424
425 /** @brief Transform all the integration points from the column vectors
426 of @a matrix from reference coordinates to physical
427 coordinates and store them as column vectors in @a result. */
428 virtual void Transform(const DenseMatrix &matrix, DenseMatrix &result);
429
430 /// Return the order of the current element we are using for the transformation.
431 virtual int Order() const { return FElem->GetOrder(); }
432
433 /// Return the order of the elements of the Jacobian of the transformation.
434 virtual int OrderJ() const;
435
436 /** @brief Return the order of the determinant of the Jacobian (weight)
437 of the transformation. */
438 virtual int OrderW() const;
439
440 /// Return the order of $ adj(J)^T \nabla fi $
441 virtual int OrderGrad(const FiniteElement *fe) const;
442
443 virtual int GetSpaceDim() const { return PointMat.Height(); }
444
445 /** @brief Transform a point @a pt from physical space to a point @a ip in
446 reference space and optionally can set a solver tolerance using @a phys_tol. */
447 /** Attempt to find the IntegrationPoint that is transformed into the given
448 point in physical space. If the inversion fails a non-zero value is
449 returned. This method is not 100 percent reliable for non-linear
450 transformations. */
451 virtual int TransformBack(const Vector & v, IntegrationPoint & ip,
452 const real_t phys_rel_tol = 1e-15)
453 {
454 InverseElementTransformation inv_tr(this);
455 inv_tr.SetPhysicalRelTol(phys_rel_tol);
456 return inv_tr.Transform(v, ip);
457 }
458
460
461 MFEM_DEPRECATED void FinalizeTransformation() {}
462};
463
471
472/** @brief A specialized ElementTransformation class representing a face and
473 its two neighboring elements.
474
475 This class can be used as a container for the element transformation data
476 needed for integrating discontinuous fields on element interfaces in a
477 Discontinuous Galerkin (DG) context.
478
479 The secondary purpose of this class is to enable the
480 GridFunction::GetValue function, and various related functions, to properly
481 evaluate fields with limited continuity on boundary elements.
482*/
484{
485private:
486
487 // Bitwise OR of ConfigMasks
488 int mask;
489
490 IntegrationPoint eip1, eip2;
491
492protected: // interface for Mesh to be able to configure this object.
493
494 friend class Mesh;
495#ifdef MFEM_USE_MPI
496 friend class ParMesh;
497#endif
498
499 /// Set the mask indicating which portions of the object have been setup
500 /** The argument @a m is a bitmask used in
501 Mesh::GetFaceElementTransformations to indicate which portions of the
502 FaceElementTransformations object have been configured.
503
504 mask & 1: Elem1 is configured
505 mask & 2: Elem2 is configured
506 mask & 4: Loc1 is configured
507 mask & 8: Loc2 is configured
508 mask & 16: The Face transformation itself is configured
509 */
510 void SetConfigurationMask(int m) { mask = m; }
511
512public:
513
515 {
516 HAVE_ELEM1 = 1, ///< Element on side 1 is configured
517 HAVE_ELEM2 = 2, ///< Element on side 2 is configured
518 HAVE_LOC1 = 4, ///< Point transformation for side 1 is configured
519 HAVE_LOC2 = 8, ///< Point transformation for side 2 is configured
520 HAVE_FACE = 16 ///< Face transformation is configured
521 };
522
524 Geometry::Type &FaceGeom; ///< @deprecated Use GetGeometryType instead
526 ElementTransformation *Face; ///< @deprecated No longer necessary
528
530
531 /** @brief Method to set the geometry type of the face.
532
533 @note This method is designed to be used when
534 [Par]Mesh::GetFaceTransformation will not be called i.e. when the face
535 transformation will not be needed but the neighboring element
536 transformations will be. Using this method to override the GeometryType
537 should only be done with great care.
538 */
540
541 /** @brief Return the mask defining the configuration state.
542
543 The mask value indicates which portions of FaceElementTransformations
544 object have been configured.
545
546 mask & 1: Elem1 is configured
547 mask & 2: Elem2 is configured
548 mask & 4: Loc1 is configured
549 mask & 8: Loc2 is configured
550 mask & 16: The Face transformation itself is configured
551 */
552 int GetConfigurationMask() const { return mask; }
553
554 /** @brief Set the integration point in the Face and the two neighboring
555 elements, if present.
556
557 The point @a face_ip must be in the reference coordinate system of the
558 face.
559 */
560 void SetIntPoint(const IntegrationPoint *face_ip);
561
562 /** @brief Set the integration point in the Face and the two neighboring
563 elements, if present.
564
565 This is a more expressive member function name than SetIntPoint, which
566 in this special case, does the same thing. This function can be used for
567 greater code clarity.
568 */
569 inline void SetAllIntPoints(const IntegrationPoint *face_ip)
571
572 /** @brief Get a const reference to the integration point in neighboring
573 element 1 corresponding to the currently set integration point on the
574 face.
575
576 This IntegrationPoint object will only contain up-to-date data if
577 SetIntPoint or SetAllIntPoints has been called with the latest
578 integration point for the face and the appropriate point transformation
579 has been configured. */
580 const IntegrationPoint &GetElement1IntPoint() { return eip1; }
581
582 /** @brief Get a const reference to the integration point in neighboring
583 element 2 corresponding to the currently set integration point on the
584 face.
585
586 This IntegrationPoint object will only contain up-to-date data if
587 SetIntPoint or SetAllIntPoints has been called with the latest
588 integration point for the face and the appropriate point transformation
589 has been configured. */
590 const IntegrationPoint &GetElement2IntPoint() { return eip2; }
591
592 virtual void Transform(const IntegrationPoint &, Vector &);
593 virtual void Transform(const IntegrationRule &, DenseMatrix &);
594 virtual void Transform(const DenseMatrix &matrix, DenseMatrix &result);
595
600
601 /** @brief Check for self-consistency: compares the result of mapping the
602 reference face vertices to physical coordinates using the three
603 transformations: face, element 1, and element 2.
604
605 @param[in] print_level If set to a positive number, print the physical
606 coordinates of the face vertices computed through
607 all available transformations: face, element 1,
608 and/or element 2.
609 @param[in,out] out The output stream to use for printing.
610
611 @returns A maximal distance between physical coordinates of face vertices
612 that should coincide. A successful check should return a small
613 number relative to the mesh extents. If less than 2 of the three
614 transformations are set, returns 0.
615
616 @warning This check will generally fail on periodic boundary faces.
617 */
618 real_t CheckConsistency(int print_level = 0,
619 std::ostream &out = mfem::out);
620};
621
622/** Elem1(Loc1(x)) = Face(x) = Elem2(Loc2(x))
623
624
625 Physical Space
626
627 *--------* ^ *--------*
628 Elem1No / / \ / \ / \ \ Elem2No
629 / / \ / \ / \ \
630 / / n \ / \ / \ \
631 *--------* ==> * ( ) * *--------*
632 \ \ / \ / \ / /
633 \ \ / \ / \ / /
634 \ \ / \ / \ / /
635 *--------* v *--------*
636
637 ^ ^
638 | ^ |
639 Elem1 | | | Elem2
640 | | Face |
641 |
642 *--------* *--------*
643 / /| / /|
644 1 *--------* | 1 *--------* 1 *--------* |
645 | | | Loc1 | | Loc2 | | |
646 | | * <----- | x | -----> | | *
647 | |/ | | | |/
648 *--------* *--------* *--------*
649 0 1 0 1 0 1
650
651 Reference Space
652*/
653
654}
655
656#endif
Data type dense matrix using column-major storage.
Definition densemat.hpp:24
const Mesh * mesh
The Mesh object containing the element.
Definition eltrans.hpp:84
virtual const DenseMatrix & EvalJacobian()=0
Evaluate the Jacobian of the transformation at the IntPoint and store it in dFdx.
virtual int Order() const =0
Return the order of the current element we are using for the transformation.
virtual int TransformBack(const Vector &pt, IntegrationPoint &ip, const real_t phys_tol=1e-15)=0
Transform a point pt from physical space to a point ip in reference space and optionally can set a so...
const DenseMatrix & TransposeAdjugateJacobian()
Return the transpose of the adjugate of the Jacobian matrix of the transformation at the currently se...
Definition eltrans.hpp:140
Geometry::Type GetGeometryType() const
Return the Geometry::Type of the reference element.
Definition eltrans.hpp:162
virtual void Transform(const IntegrationRule &, DenseMatrix &)=0
Transform all the integration points from the integration rule from reference coordinates to physical...
virtual int OrderGrad(const FiniteElement *fe) const =0
Return the order of .
const DenseMatrix & Hessian()
Return the Hessian matrix of the transformation at the currently set IntegrationPoint,...
Definition eltrans.hpp:125
const DenseMatrix & InverseJacobian()
Return the inverse of the Jacobian matrix of the transformation at the currently set IntegrationPoint...
Definition eltrans.hpp:145
const DenseMatrix & AdjugateJacobian()
Return the adjugate of the Jacobian matrix of the transformation at the currently set IntegrationPoin...
Definition eltrans.hpp:135
const IntegrationPoint & GetIntPoint()
Get a const reference to the currently set integration point. This will return NULL if no integration...
Definition eltrans.hpp:98
const DenseMatrix & EvalAdjugateJ()
Definition eltrans.cpp:36
virtual int GetSpaceDim() const =0
Get the dimension of the target (physical) space.
real_t Weight()
Return the weight of the Jacobian matrix of the transformation at the currently set IntegrationPoint....
Definition eltrans.hpp:131
virtual int OrderJ() const =0
Return the order of the elements of the Jacobian of the transformation.
const IntegrationPoint * IntPoint
Definition eltrans.hpp:26
const DenseMatrix & Jacobian()
Return the Jacobian matrix of the transformation at the currently set IntegrationPoint,...
Definition eltrans.hpp:119
int GetDimension() const
Return the topological dimension of the reference element.
Definition eltrans.hpp:165
virtual const DenseMatrix & EvalHessian()=0
Evaluate the Hessian of the transformation at the IntPoint and store it in d2Fdx2.
const DenseMatrix & EvalInverseJ()
Definition eltrans.cpp:57
virtual int OrderW() const =0
Return the order of the determinant of the Jacobian (weight) of the transformation.
void SetIntPoint(const IntegrationPoint *ip)
Set the integration point ip that weights and Jacobians will be evaluated at.
Definition eltrans.hpp:93
const DenseMatrix & EvalTransAdjugateJ()
Definition eltrans.cpp:46
virtual void Transform(const IntegrationPoint &, Vector &)=0
Transform integration point from reference coordinates to physical coordinates and store them in the ...
virtual void Transform(const DenseMatrix &matrix, DenseMatrix &result)=0
Transform all the integration points from the column vectors of matrix from reference coordinates to ...
void Reset()
Force the reevaluation of the Jacobian in the next call.
Definition eltrans.hpp:89
A specialized ElementTransformation class representing a face and its two neighboring elements.
Definition eltrans.hpp:484
ElementTransformation * Elem2
Definition eltrans.hpp:525
void SetIntPoint(const IntegrationPoint *face_ip)
Set the integration point in the Face and the two neighboring elements, if present.
Definition eltrans.cpp:565
ElementTransformation * Elem1
Definition eltrans.hpp:525
const IntegrationPoint & GetElement1IntPoint()
Get a const reference to the integration point in neighboring element 1 corresponding to the currentl...
Definition eltrans.hpp:580
ElementTransformation * Face
Definition eltrans.hpp:526
@ HAVE_ELEM2
Element on side 2 is configured.
Definition eltrans.hpp:517
@ HAVE_LOC1
Point transformation for side 1 is configured.
Definition eltrans.hpp:518
@ HAVE_ELEM1
Element on side 1 is configured.
Definition eltrans.hpp:516
@ HAVE_FACE
Face transformation is configured.
Definition eltrans.hpp:520
@ HAVE_LOC2
Point transformation for side 2 is configured.
Definition eltrans.hpp:519
int GetConfigurationMask() const
Return the mask defining the configuration state.
Definition eltrans.hpp:552
const IntegrationPoint & GetElement2IntPoint()
Get a const reference to the integration point in neighboring element 2 corresponding to the currentl...
Definition eltrans.hpp:590
IntegrationPointTransformation Loc1
Definition eltrans.hpp:527
IntegrationPointTransformation & GetIntPoint1Transformation()
Definition eltrans.cpp:604
void SetGeometryType(Geometry::Type g)
Method to set the geometry type of the face.
Definition eltrans.hpp:539
void SetAllIntPoints(const IntegrationPoint *face_ip)
Set the integration point in the Face and the two neighboring elements, if present.
Definition eltrans.hpp:569
ElementTransformation & GetElement1Transformation()
Definition eltrans.cpp:588
virtual void Transform(const IntegrationPoint &, Vector &)
Transform integration point from reference coordinates to physical coordinates and store them in the ...
Definition eltrans.cpp:619
void SetConfigurationMask(int m)
Set the mask indicating which portions of the object have been setup.
Definition eltrans.hpp:510
IntegrationPointTransformation Loc2
Definition eltrans.hpp:527
IntegrationPointTransformation & GetIntPoint2Transformation()
Definition eltrans.cpp:612
real_t CheckConsistency(int print_level=0, std::ostream &out=mfem::out)
Check for self-consistency: compares the result of mapping the reference face vertices to physical co...
Definition eltrans.cpp:643
ElementTransformation & GetElement2Transformation()
Definition eltrans.cpp:596
Abstract class for all finite elements.
Definition fe_base.hpp:239
int GetOrder() const
Returns the order of the finite element. In the case of anisotropic orders, returns the maximum order...
Definition fe_base.hpp:333
Geometry::Type GetGeomType() const
Returns the Geometry::Type of the reference element.
Definition fe_base.hpp:326
void SetType(int t)
Set the Quadrature1D type of points to use for subdivision.
Definition geom.hpp:341
static const int Dimension[NumGeom]
Definition geom.hpp:47
IsoparametricTransformation Transf
Definition eltrans.hpp:467
void Transform(const IntegrationPoint &, IntegrationPoint &)
Definition eltrans.cpp:543
Class for integration point with weight.
Definition intrules.hpp:35
Class for an integration rule - an Array of IntegrationPoint.
Definition intrules.hpp:100
The inverse transformation of a given ElementTransformation.
Definition eltrans.hpp:187
virtual int Transform(const Vector &pt, IntegrationPoint &ip)
Given a point, pt, in physical space, find its reference coordinates, ip.
Definition eltrans.cpp:336
void SetPrintLevel(int pr_level)
Set the desired print level, useful for debugging.
Definition eltrans.hpp:335
InitGuessType
Algorithms for selecting an initial guess.
Definition eltrans.hpp:191
@ GivenPoint
Use a specific point, set with SetInitialGuess().
Definition eltrans.hpp:201
@ Center
Use the center of the reference element.
Definition eltrans.hpp:192
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:97
const IntegrationPoint * ip0
Definition eltrans.hpp:235
InverseElementTransformation(ElementTransformation *Trans=NULL)
Construct the InverseElementTransformation with default parameters.
Definition eltrans.hpp:276
void SetInitGuessPointsType(int q_type)
Set the Quadrature1D type used for the Closest* initial guess types.
Definition eltrans.hpp:305
void NewtonPrintPoint(const char *prefix, const Vector &pt, const char *suffix)
Definition eltrans.cpp:158
void SetElementTol(real_t el_tol)
Set the tolerance used to determine if a point lies inside or outside of the reference element.
Definition eltrans.hpp:329
void SetSolverType(SolverType stype)
Specify which algorithm to use for solving the transformation equation, i.e. when calling the Transfo...
Definition eltrans.hpp:315
void SetInitialGuess(const IntegrationPoint &init_ip)
Set the initial guess for subsequent calls to Transform(), switching to the GivenPoint InitGuessType ...
Definition eltrans.hpp:301
TransformResult
Values returned by Transform().
Definition eltrans.hpp:224
@ Outside
The point is probably outside the element.
Definition eltrans.hpp:226
@ Unknown
The algorithm failed to determine where the point is.
Definition eltrans.hpp:227
@ Inside
The point is inside the element.
Definition eltrans.hpp:225
void SetPhysicalRelTol(real_t phys_rel_tol)
Set the relative physical-space convergence tolerance.
Definition eltrans.hpp:324
int FindClosestPhysPoint(const Vector &pt, const IntegrationRule &ir)
Find the IntegrationPoint mapped closest to pt.
Definition eltrans.cpp:70
int NewtonSolve(const Vector &pt, IntegrationPoint &ip)
Definition eltrans.cpp:172
void SetMaxIter(int max_it)
Set the maximum number of iterations when solving for a reference point.
Definition eltrans.hpp:318
void SetReferenceTol(real_t ref_sp_tol)
Set the reference-space convergence tolerance.
Definition eltrans.hpp:321
void SetInitGuessRelOrder(int order)
Set the relative order used for the Closest* initial guess types.
Definition eltrans.hpp:311
void SetInitialGuessType(InitGuessType itype)
Choose how the initial guesses for subsequent calls to Transform() will be selected.
Definition eltrans.hpp:297
void NewtonPrint(int mode, real_t val)
Definition eltrans.cpp:129
void SetTransformation(ElementTransformation &Trans)
Set a new forward ElementTransformation, Trans.
Definition eltrans.hpp:293
ElementTransformation * T
Definition eltrans.hpp:232
A standard isoparametric element transformation.
Definition eltrans.hpp:363
void SetPointMat(const DenseMatrix &pm)
Set the underlying point matrix describing the transformation.
Definition eltrans.hpp:402
void SetFE(const FiniteElement *FE)
Set the element that will be used to compute the transformations.
Definition eltrans.hpp:382
virtual int TransformBack(const Vector &v, IntegrationPoint &ip, const real_t phys_rel_tol=1e-15)
Transform a point pt from physical space to a point ip in reference space and optionally can set a so...
Definition eltrans.hpp:451
virtual int OrderJ() const
Return the order of the elements of the Jacobian of the transformation.
Definition eltrans.cpp:439
const FiniteElement * GetFE() const
Get the current element used to compute the transformations.
Definition eltrans.hpp:390
virtual int OrderGrad(const FiniteElement *fe) const
Return the order of .
Definition eltrans.cpp:467
virtual int Order() const
Return the order of the current element we are using for the transformation.
Definition eltrans.hpp:431
DenseMatrix & GetPointMat()
Write access to the stored point matrix. Use with caution.
Definition eltrans.hpp:411
virtual int GetSpaceDim() const
Get the dimension of the target (physical) space.
Definition eltrans.hpp:443
void SetIdentityTransformation(Geometry::Type GeomType)
Set the FiniteElement Geometry for the reference elements being used.
Definition eltrans.cpp:379
const DenseMatrix & GetPointMat() const
Return the stored point matrix.
Definition eltrans.hpp:405
virtual int OrderW() const
Return the order of the determinant of the Jacobian (weight) of the transformation.
Definition eltrans.cpp:453
MFEM_DEPRECATED void FinalizeTransformation()
Definition eltrans.hpp:461
virtual void Transform(const IntegrationPoint &, Vector &)
Transform integration point from reference coordinates to physical coordinates and store them in the ...
Definition eltrans.cpp:488
Mesh data type.
Definition mesh.hpp:56
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
Definition operator.hpp:66
Class for parallel meshes.
Definition pmesh.hpp:34
A class container for 1D quadrature type constants.
Definition intrules.hpp:394
Vector data type.
Definition vector.hpp:80
IntegrationType itype
Definition ex38.cpp:46
OutStream out(std::cout)
Global stream used by the library for standard output. Initially it uses the same std::streambuf as s...
Definition globals.hpp:66
float real_t
Definition config.hpp:43