hvfem.py¶
Define functions for high-order vector finite element method.
-
petgem.hvfem.
AffineTetrahedron
(X)[source]¶ Compute affine coordinates and their gradients.
- Parameters
X (ndarray) – point coordinates
- Returns
affine coordinates and gradients of affine coordinates
- Return type
ndarray
Note
References:n Fuentes, F., Keith, B., Demkowicz, L., & Nagaraj, S. (2015). Orientation embedded high order shape functions for the exact sequence elements of all shapes. Computers & Mathematics with applications, 70(4), 353-458.
-
petgem.hvfem.
AncEE
(S, DS, Nord, Idec, N)[source]¶ Compute edge Hcurl ancillary functions and their curls.
- Parameters
S (ndarray) – affine coordinates associated to edge
DS (ndarray) – derivatives of S in R^N
Nord (int) – polynomial order
Idec (bool) – Binary flag
N (int) – spatial dimension
- Returns
edge Hcurl ancillary functions, curls of edge Hcurl ancillary functions
- Return type
ndarray
Note
References:n Idec: = FALSE s0+s1 != 1
= TRUE s0+s1 = 1
-
petgem.hvfem.
AncETri
(S, DS, Nord, Idec, N)[source]¶ Compute triangle face Hcurl ancillary functions and their curls.
- Parameters
S (ndarray) – (s0,s1,s2) affine coordinates associated to triangle face
DS (ndarray) – derivatives of S0,S1,S2
Nord (int) – polynomial order
Idec (boll) – Binary flag:
N (int) – spatial dimension
- Returns
triangle Hcurl ancillary functions and curls of triangle Hcurl ancillary functions
- Return type
ndarray
-
petgem.hvfem.
HomIJacobi
(S, DS, Nord, Minalpha, Idec, N)[source]¶ Compute values of integrated homogenized Jacobi polynomials and their gradients. Result is half of a matrix with each row associated to a fixed alpha. Alpha grows by 2 in each row.
- Parameters
S (ndarray) – (s0,s1) affine(like) coordinates
DS (ndarray) – gradients of S in R(N)
Nord (int) – max polynomial order
Minalpha (int) – first row value of alpha (integer)
Idec (bool) – decision flag to compute
- Returns
polynomial values and derivatives in x (Jacobi polynomials)
- Return type
ndarray
-
petgem.hvfem.
HomLegendre
(S, Nord)[source]¶ Compute values of homogenized Legendre polynomials.
- Parameters
S (ndarray) – affine(like) coordinates
Nord (int) – polynomial order
- Returns
polynomial values
- Return type
ndarray
-
petgem.hvfem.
OrientE
(S, DS, Nori, N)[source]¶ Compute the local to global transformations of edges.
- Parameters
S (ndarray) – projection of affine coordinates on edges
DS (ndarray) – projection of gradients of affine coordinates on edges
Nori (ndarray) – edge orientation
N (int) – number of dimensions
- Returns
global transformation of edges and global transformation of gradients of edges
- Return type
ndarray
-
petgem.hvfem.
OrientTri
(S, DS, Nori, N)[source]¶ Compute the local to global transformations of edges.
- Parameters
S (ndarray) – projection of affine coordinates on faces
DS (ndarray) – projection of gradients of affine coordinates on faces
Nori (ndarray) – face orientation
N (int) – number of dimensions
- Returns
global transformation of faces and global transformation of gradients of faces
- Return type
ndarray
-
petgem.hvfem.
PolyIJacobi
(X, T, Nord, Minalpha, Idec)[source]¶ Compute values of integrated shifted scaled Jacobi polynomials and their derivatives starting with p=1.
Result is ‘half’ of a matrix with each row associated to a fixed alpha. Alpha grows by 2 in each row.
- Parameters
X (ndarray) – coordinate from [0,1]
T (ndarray) – scaling parameter
Nord (int) – max polynomial order
Minalpha (int) – first row value of alpha
Idec (bool) – decision flag to compute (= FALSE polynomials with x and t derivatives, = TRUE polynomials with x derivatives only)
- Returns
polynomial values, derivatives in x (Jacobi polynomials), derivatives in t
-
petgem.hvfem.
PolyJacobi
(X, T, Nord, Minalpha)[source]¶ Compute values of shifted scaled Jacobi polynomials P**alpha-i.
Result is a half of a matrix with each row associated to a fixed alpha. Alpha grows by 2 in each row.
- Parameters
X (ndarray) – coordinate from [0,1]
T (float) – scaling parameter
Nord (int) – max polynomial order
Minalpha (int) – first row value of alpha (integer)
- Returns
polynomial values
- Return type
ndarray
-
petgem.hvfem.
PolyLegendre
(X, T, Nord)[source]¶ Compute values of shifted scaled Legendre polynomials.
- Parameters
X (ndarray) – coordinate from [0,1]
T (float) – scaling parameter
Nord (int) – polynomial order
- Returns
polynomial values
- Return type
ndarray
-
petgem.hvfem.
ProjectTetE
(Lam, DLam)[source]¶ Projection of tetrahedral edges in concordance with numbering of topological entities (vertices, edges, faces).
- Parameters
Lam (ndarray) – affine coordinates
DLam (ndarray) – gradients of affine coordinates
- Returns
projection of affine coordinates on edges, projection of gradients of affine coordinates on edges
- Return type
ndarray
Note
References:n Fuentes, F., Keith, B., Demkowicz, L., & Nagaraj, S. (2015). Orientation embedded high order shape functions for the exact sequence elements of all shapes. Computers & Mathematics with applications, 70(4), 353-458.
-
petgem.hvfem.
ProjectTetF
(Lam, DLam)[source]¶ Projection of tetrahedral faces in concordance with numbering of topological entities (vertices, edges, faces).
- Parameters
Lam (ndarray) – affine coordinates
DLam (ndarray) – gradients of affine coordinates
- Returns
projection of affine coordinates on faces, projection of gradients of affine coordinates on faces
- Return type
ndarray
-
petgem.hvfem.
compute2DGaussPoints
(Nord)[source]¶ Compute 2D gauss points for high-order nédélec elements.
- Parameters
Nord (int) – polynomial order of nedelec basis functions
- Returns
coordinates of gauss points and its weights
- Return type
ndarray.
-
petgem.hvfem.
compute3DGaussPoints
(Nord)[source]¶ Compute 3D gauss points for high-order nédélec elements.
- Parameters
Nord (int) – polynomial order of nedelec basis functions
- Returns
coordinates of gauss points and its weights
- Return type
ndarray.
-
petgem.hvfem.
computeBasisFunctions
(edge_orientation, face_orientation, jacobian, invjacob, Nord, points)[source]¶ Compute the basis function for a given element.
- Parameters
edges_orientation (ndarray) – orientation for edges
faces_orientation (ndarray) – orientation for faces
jacobian (ndarray) – jacobian matrix
invjacob (ndarray) – inverse of jacobian matrix
Nord (int) – polynomial order of nedelec basis functions
points (ndarray) – spatial points at which basis functions will be computed
- Returns
basis functions and its curl for p-order=Nord
- Return type
ndarray
-
petgem.hvfem.
computeBasisFunctionsReferenceElement
(edge_orientation, face_orientation, Nord, points)[source]¶ Compute the basis function for the reference element.
- Parameters
edges_orientation (ndarray) – orientation for edges
faces_orientation (ndarray) – orientation for faces
Nord (int) – polynomial order of nedelec basis functions
points (ndarray) – spatial points at which basis functions will be computed
- Returns
basis functions on reference element
- Return type
ndarray
-
petgem.hvfem.
computeConnectivityDOFS
(elemsE, elemsF, Nord)[source]¶ Compute the degrees of freedom connectivity for a given list of edges, faces and elements.
- Parameters
elemsE (ndarray) – elements-edge connectivity with dimensions = (number_elements, 6)
elemsF (ndarray) – element/faces connectivity with dimensions = (number_elements, 4)
Nord (int) – polynomial order of nedelec basis functions
- Returns
local/global dofs list for elements, dofs index on edges, dofs index on faces, dofs index on volumes, total number of dofs
- Return type
ndarray and int
Note
References:n Amor-Martin, A., Garcia-Castillo, L. E., & Garcia-Doñoro, D. D. (2016). Second-order Nédélec curl-conforming prismatic element for computational electromagnetics. IEEE Transactions on Antennas and Propagation, 64(10), 4384-4395.
-
petgem.hvfem.
computeElementOrientation
(edgesEle, nodesEle, edgesNodesEle, globalEdgesInFace)[source]¶ Compute the orientation for the computation of hierarchical basis functions of high-order (High-order nédélec basis functions).
:param ndarray edgesEle:list of element’s edges :param ndarray nodesEle: list of element’s nodes :param ndarray edgesNodesEle: list of nodes for each edge in edgesEle :param ndarray globalEdgesInFace: list of edges for each face :return: orientation for edges and orientation for faces :rtype: ndarray.
Note
References:n Amor-Martin, A., Garcia-Castillo, L. E., & Garcia-Doñoro, D. D. (2016). Second-order Nédélec curl-conforming prismatic element for computational electromagnetics. IEEE Transactions on Antennas and Propagation, 64(10), 4384-4395.
-
petgem.hvfem.
computeElementalMatrices
(edge_orientation, face_orientation, jacobian, invjacob, Nord, sigmaEle)[source]¶ Compute the elemental mass matrix and stiffness matrix based ons high-order vector finite element.
- Parameters
edges_orientation (ndarray) – orientation for edges
faces_orientation (ndarray) – orientation for faces
jacobian (ndarray) – jacobian matrix
invjacob (ndarray) – inverse of jacobian matrix
Nord (int) – polynomial order of vector basis functions
sigmaEle (ndarray) – element conductivity with dimensions (1, 2), (horizontal and vertical)
- Returns
elemental mass matrix and elemental stiffness matrix
- Return type
ndarray
Note
References:n Fuentes, F., Keith, B., Demkowicz, L., & Nagaraj, S. (2015). Orientation embedded high order shape functions for the exact sequence elements of all shapes. Computers & Mathematics with applications, 70(4), 353-458.
-
petgem.hvfem.
computeJacobian
(eleNodes)[source]¶ Compute the jacobian and its inverse.
- Parameters
eleNodes (ndarray) – spatial coordinates of the nodes with dimensions = (4,3)
- Returns
jacobian matrix and its inverse.
- Return type
ndarray
-
petgem.hvfem.
computeSourceVectorRotation
(azimuth, dip)[source]¶ Compute the weigths vector for source rotation in the xyz plane.
- Parameters
azimuth (float) – degrees for x-y plane rotation
dip (float) – degrees for x-z plane rotation
- Returns
weigths for source rotation
- Return type
ndarray.
-
petgem.hvfem.
get2DJacobDet
(coordEle, faceNumber)[source]¶ Compute the determinant of the jacobian for 2D integrals (when 3D basis functions are used)
- Parameters
coordEle (ndarray) – coordinates of the tetrahedron
faceNumber (int) – local face number
- Returns
determinant of the 2D jacobian
- Return type
ndarray
-
petgem.hvfem.
getFaceByLocalNodes
(faceNumber)[source]¶ Get local nodes ordering for a given face
- Parameters
faceNumber (int) – local face number
- Returns
list of local nodes for faceNumber
- Return type
ndarray
-
petgem.hvfem.
getNeumannBCface
(face_flag, polarization, ud)[source]¶ Get Neumann boundary condition for boundary face.
- Parameters
face_flag (int) – face flag (side on which face belongs)
polarization (int) – polarization mode (x-mode or y-mode)
ud (ndarray) – magnetic field vector
- Returns
value of Neumann boundary condition for boundary face
- Return type
ndarray
-
petgem.hvfem.
getNormalVector
(faceNumber, invJacobMatrix)[source]¶ This function computes the normal vector for a given tetrahedral face.
- Parameters
faceNumber (int) – local face number
invJacobMatrix (int) – inverse of jacobian matrix
- Returns
face normal vector
- Return type
ndarray
-
petgem.hvfem.
getRealFromReference
(rRef, verticesReal)[source]¶ Translate a point defined in the reference element (rRef) to the real element (rReal) defined by verticesReal
- Parameters
rRef (ndarray) – reference coordinates
verticesReal (ndarray) – real coordinates of element (4x3) = (4 nodes x 3 coordinates)
- Returns
points in real element
- Return type
ndarray
-
petgem.hvfem.
shape3DETet
(X, Nord, NoriE, NoriF)[source]¶ Compute values of 3D tetrahedron element H(curl) shape functions and their derivatives.
- Parameters
X (ndarray) – master tetrahedron coordinates from (0,1)^3
Nord (int) – polynomial order
NoriE (ndarray) – edge orientation
NoriF (ndarray) – face orientation
- Returns
number of dof, values of the shape functions at the point, curl of the shape functions
- Return type
ndarray.
Note
References:n Fuentes, F., Keith, B., Demkowicz, L., & Nagaraj, S. (2015). Orientation embedded high order shape functions for the exact sequence elements of all shapes. Computers & Mathematics with applications, 70(4), 353-458.
-
petgem.hvfem.
tetrahedronXYZToXiEtaZeta
(eleNodes, points)[source]¶ Compute the reference tetrahedron coordinates from xyz global tetrahedron coordinates.
- Parameters
eleNodes (ndarray) – spatial coordinates of the nodes with dimensions = (4,3)
points (ndarray) – xyz points coordinates to be transformed
- Returns
xietazeta points coordinates
- Return type
ndarray