GemaCoreLib
The GeMA Core library
gmWedgeShape.h
Go to the documentation of this file.
1 /************************************************************************
2 **
3 ** Copyright (C) 2014 by Carlos Augusto Teixera Mendes
4 ** All rights reserved.
5 **
6 ** This file is part of the "GeMA" software. It's use should respect
7 ** the terms in the license agreement that can be found together
8 ** with this source code.
9 ** It is provided AS IS, with NO WARRANTY OF ANY KIND,
10 ** INCLUDING THE WARRANTY OF DESIGN, MERCHANTABILITY AND FITNESS FOR
11 ** A PARTICULAR PURPOSE.
12 **
13 ************************************************************************/
14 
25 #ifndef _GEMA_WEDGE_SHAPE_H_
26 #define _GEMA_WEDGE_SHAPE_H_
27 
28 #include "gmShape.h"
29 
30 
31 //--------------------------------------------------------------------------------------------------------------------------
32 // Wedge Elements
33 //--------------------------------------------------------------------------------------------------------------------------
34 
43 {
44 public:
45  // See comments on the base class
46  virtual int numNaturalCoord() const { return 3; }
47 
48  // See comments on the base class
49  virtual int numCartesianCoord() const { return 3; }
50 
51  // See comments on the base class
52  virtual void naturalCoordLimits(int coord, double* min, double* max) const
53  {
54  *min = 0.0; *max = 1.0;
55  if (coord == 2) { *min = -1.0; }
56  }
57 
58  // See comments on the base class
59  virtual void naturalCenter(GmVector& coord) const { coord.set_size(3); coord[0] = coord[1] = 1/3.0; coord[2] = 0.0; } // Center = (1/3, 1/3, 0) for all wedges
60 
61  virtual bool translateEdgePoint(int edge, const GmVector& srcEdgeCoord, GmVector& elementCoord) const;
62  virtual bool translateFacePoint(int face, const GmVector& srcFaceCoord, GmVector& elementCoord) const;
63 
64  virtual double shapeCartesianPartialsFromPJ(const GmMatrix& P, const GmMatrix& J, GmMatrix& dN, bool transposed = false) const;
65 
66  // See comments on the base class
67  virtual double scaledJacobianDet(const GmMatrix& J) const { return 0.5* arma::det(J); }
68 
69  // See comments on the base class
70  virtual void jacobianAndPartials(const GmVector& ncoord, const GmMatrix& X, GmMatrix& J, GmMatrix& P, bool transposed = false) const
71  {
72  jacIndependentNatCoord(ncoord, X, J, P, transposed);
73  }
74 
75  // See comments on the base class
76  virtual double borderScalingFactor(int border, const GmVector& borderCoord, const GmVector& elementCoord,
77  const GmMatrix& X, bool transposed = false) const
78  {
79  return faceScalingFactor(border, borderCoord, elementCoord, X, transposed);
80  }
81 
82  virtual double edgeScalingFactor(int border, const GmVector& borderCoord, const GmVector& elementCoord,
83  const GmMatrix& X, bool transposed = false) const;
84 
85  virtual double faceScalingFactor(int border, const GmVector& borderCoord, const GmVector& elementCoord,
86  const GmMatrix& X, bool transposed = false) const;
87 
88 protected:
89  const int* fixedEdgeCoord(int border) const;
90  const double* fixedEdgeValue(int border) const;
91  int fixedFaceCoord(int border) const;
92  double fixedFaceValue(int border) const;
93 
94  // See comments on the base class. Calls the default algorithm for Wedge elements
95  virtual bool gradientBasedCartesianToNatural(const GmVector& coord, const GmMatrix& X, double tol,
96  int maxIter, double natTol, GmVector& ncoord, bool* inside) const
97  {
98  return gradientBasedCartesianToNaturalWedge(coord, X, tol, maxIter, natTol, ncoord, inside);
99  }
100 };
101 
111 {
112 public:
113  // See comments on the base class
114  virtual GmCellType elemType() const { return GM_WEDGE6; }
115 
116  // See comments on the base class
117  virtual int numFunctions() const { return 6; }
118 
119  virtual void nodeNaturalCoord(int node, GmVector& coord) const;
120 
121  virtual void shapeValues(const GmVector& ncoord, GmVector& N) const;
122  virtual void shapePartials(const GmVector& ncoord, GmMatrix& dN, bool transposed = false) const;
123 };
124 
125 
135 {
136 public:
137  // See comments on the base class
138  virtual GmCellType elemType() const { return GM_WEDGE15; }
139 
140  // See comments on the base class
141  virtual int numFunctions() const { return 15; }
142 
143  virtual void nodeNaturalCoord(int node, GmVector& coord) const;
144 
145  virtual void shapeValues(const GmVector& ncoord, GmVector& N) const;
146  virtual void shapePartials(const GmVector& ncoord, GmMatrix& dN, bool transposed = false) const;
147 };
148 
149 #endif
150 
A wedge with 6 nodes - linear.
Definition: gmCellType.h:62
virtual bool translateFacePoint(int face, const GmVector &srcFaceCoord, GmVector &elementCoord) const =0
Given a "quad like" or "tri like" face coordinate (-1 to 1 pair for quad faces and 0 to 1 barycentric...
virtual double borderScalingFactor(int border, const GmVector &borderCoord, const GmVector &elementCoord, const GmMatrix &X, bool transposed=false) const
Returns the scaling factor needed when calculating integrals over borders (edges or faces)....
Definition: gmWedgeShape.h:76
virtual GmCellType elemType() const
Returns the type of this element.
Definition: gmWedgeShape.h:138
virtual GmCellType elemType() const
Returns the type of this element.
Definition: gmWedgeShape.h:114
GmShape specialization for a Wedge with 6 nodes object.
Definition: gmWedgeShape.h:134
virtual double scaledJacobianDet(const GmMatrix &J) const
Returns the Jacobian determinant multiplied by the scaling factor needed for transforming the differe...
Definition: gmWedgeShape.h:67
virtual void shapePartials(const GmVector &ncoord, GmMatrix &dN, bool transposed=false) const =0
Function used to calculate shape function partial derivatives with respect to its natural coordinates...
virtual int numFunctions() const
Returns the number of shape functions of this element type (equal to the number of nodes)
Definition: gmWedgeShape.h:141
virtual bool translateEdgePoint(int edge, const GmVector &srcEdgeCoord, GmVector &elementCoord) const =0
Given a "bar like" edge coordinate (from -1 to 1), at the given edge, fills elementCoord with the equ...
virtual int numCartesianCoord() const
Returns the number of cartesian coordinates expected by this element type.
Definition: gmWedgeShape.h:49
Shape function handling base classe.
Definition: gmShape.h:37
bool gradientBasedCartesianToNaturalWedge(const GmVector &coord, const GmMatrix &X, double tol, int maxIter, double natTol, GmVector &ncoord, bool *inside) const
An implementation of the gradient based cartesian to natural algorithm specific to Wedge elements....
Definition: gmShape.cpp:367
A wedge with 15 nodes - quadratic.
Definition: gmCellType.h:63
virtual void shapeValues(const GmVector &ncoord, GmVector &N) const =0
Function used to evaluate the set of shape functions over a point defined by its natural coordinates.
virtual void jacIndependentNatCoord(const GmVector &ncoord, const GmMatrix &X, GmMatrix &J, GmMatrix &P, bool transposed) const
A generic implementation of jacobianAndPartials() for the case where all natural coordinates are inde...
Definition: gmShape.cpp:800
virtual int numFunctions() const
Returns the number of shape functions of this element type (equal to the number of nodes)
Definition: gmWedgeShape.h:117
virtual double faceScalingFactor(int border, const GmVector &borderCoord, const GmVector &elementCoord, const GmMatrix &X, bool transposed=false) const =0
Returns the scaling factor needed when calculating integrals over element faces. This has the same ef...
virtual void naturalCoordLimits(int coord, double *min, double *max) const
Fills min and max with the domain limits for the given natural coordinate (between 0 and numNaturalCo...
Definition: gmWedgeShape.h:52
GmShape specialization for a Wedge with 6 nodes object.
Definition: gmWedgeShape.h:110
#define GMC_API_EXPORT
Macro for controling if the class is being exported (GEMA_CORE_LIB defined) or imported (GEMA_CORE_LI...
Definition: gmCoreConfig.h:35
virtual void jacobianAndPartials(const GmVector &ncoord, const GmMatrix &X, GmMatrix &J, GmMatrix &P, bool transposed=false) const
This function does the same calculations as the jacobian() call, but also filling the extra parameter...
Definition: gmWedgeShape.h:70
GmCellType
Mesh Cell types. Don't change type orders or add types without reading comments below.
Definition: gmCellType.h:30
virtual double shapeCartesianPartialsFromPJ(const GmMatrix &P, const GmMatrix &J, GmMatrix &dN, bool transposed=false) const =0
Alternative version of shapeCartesianPartialsXxxx() to calculate shape function partial derivatives w...
virtual int numNaturalCoord() const
Returns the number of natural coordinates used by this element type.
Definition: gmWedgeShape.h:46
GmShape specialization for Wedges, containing common functions used by Wedge elements,...
Definition: gmWedgeShape.h:42
virtual void nodeNaturalCoord(int node, GmVector &coord) const =0
Fills the coord vector with the set of natural coordinates for the reference shape function node....
arma::vec GmVector
The basic type for a GeMA vector object. Currently based on an Armadillo vector.
Definition: gmVector.h:34
arma::mat GmMatrix
The basic type for a GeMA matrix object. Currently based on an Armadillo matrix.
Definition: gmMatrix.h:38
Declaration of the GmShape base class.
virtual void naturalCenter(GmVector &coord) const
Fills the coord vector with the set of natural coordinates for the element center.
Definition: gmWedgeShape.h:59
virtual double edgeScalingFactor(int border, const GmVector &borderCoord, const GmVector &elementCoord, const GmMatrix &X, bool transposed=false) const =0
Returns the scaling factor needed when calculating integrals over element edges. This has the same ef...