FIAT package¶
Submodules¶
FIAT.P0 module¶
- class FIAT.P0.P0(ref_el)[source]¶
Bases:
CiarletElement
FIAT.Sminus module¶
- class FIAT.Sminus.TrimmedSerendipity(ref_el, degree, mapping)[source]¶
Bases:
FiniteElement
- entity_closure_dofs()[source]¶
Return the map of topological entities to degrees of freedom on the closure of those entities for the finite element.
- entity_dofs()[source]¶
Return the map of topological entities to degrees of freedom for the finite element.
- tabulate(order, points, entity=None)[source]¶
Return tabulated values of derivatives up to given order of basis functions at given points.
- Parameters:
order – The maximum order of derivative.
points – An iterable of points.
entity – Optional (dimension, entity number) pair indicating which topological entity of the reference element to tabulate on. If
None
, default cell-wise tabulation is performed.
- class FIAT.Sminus.TrimmedSerendipityEdge(ref_el, degree)[source]¶
Bases:
TrimmedSerendipity
- class FIAT.Sminus.TrimmedSerendipityFace(ref_el, degree)[source]¶
Bases:
TrimmedSerendipity
FIAT.SminusCurl module¶
- class FIAT.SminusCurl.TrimmedSerendipity(ref_el, degree, mapping)[source]¶
Bases:
FiniteElement
- entity_closure_dofs()[source]¶
Return the map of topological entities to degrees of freedom on the closure of those entities for the finite element.
- entity_dofs()[source]¶
Return the map of topological entities to degrees of freedom for the finite element.
- tabulate(order, points, entity=None)[source]¶
Return tabulated values of derivatives up to given order of basis functions at given points.
- Parameters:
order – The maximum order of derivative.
points – An iterable of points.
entity – Optional (dimension, entity number) pair indicating which topological entity of the reference element to tabulate on. If
None
, default cell-wise tabulation is performed.
- class FIAT.SminusCurl.TrimmedSerendipityCurl(ref_el, degree)[source]¶
Bases:
TrimmedSerendipity
FIAT.SminusDiv module¶
- class FIAT.SminusDiv.TrimmedSerendipity(ref_el, degree, mapping)[source]¶
Bases:
FiniteElement
- entity_closure_dofs()[source]¶
Return the map of topological entities to degrees of freedom on the closure of those entities for the finite element.
- entity_dofs()[source]¶
Return the map of topological entities to degrees of freedom for the finite element.
- tabulate(order, points, entity=None)[source]¶
Return tabulated values of derivatives up to given order of basis functions at given points.
- Parameters:
order – The maximum order of derivative.
points – An iterable of points.
entity – Optional (dimension, entity number) pair indicating which topological entity of the reference element to tabulate on. If
None
, default cell-wise tabulation is performed.
- class FIAT.SminusDiv.TrimmedSerendipityDiv(ref_el, degree)[source]¶
Bases:
TrimmedSerendipity
FIAT.alfeld_sorokina module¶
- class FIAT.alfeld_sorokina.AlfeldSorokina(ref_el, degree=2)[source]¶
Bases:
CiarletElement
The Alfeld-Sorokina C0 quadratic macroelement with C0 divergence.
This element belongs to a Stokes complex, and is paired with CG1(Alfeld).
FIAT.argyris module¶
- class FIAT.argyris.Argyris(ref_el, degree=5, variant=None)[source]¶
Bases:
CiarletElement
The Argyris finite element.
- Parameters:
ref_el – The reference element.
degree – The degree.
variant – optional variant specifying the types of nodes.
variant can be chosen from [“point”, “integral”, “integral(q)”] “point” -> dofs are evaluated by point evaluation. “integral” -> dofs are evaluated by quadrature rules with the minimum degree required for unisolvence. “integral(q)” -> dofs are evaluated by quadrature rules with the minimum degree required for unisolvence plus q.
FIAT.arnold_qin module¶
- class FIAT.arnold_qin.ArnoldQin(ref_el, degree=2, reduced=False)[source]¶
Bases:
CiarletElement
The Arnold-Qin C0(Alfeld) quadratic macroelement with divergence in P0. This element belongs to a Stokes complex, and is paired with unsplit DG0.
FIAT.arnold_winther module¶
Implementation of the Arnold-Winther finite elements.
- class FIAT.arnold_winther.ArnoldWinther(ref_el, degree=3)[source]¶
Bases:
CiarletElement
The definition of the conforming Arnold-Winther element.
- class FIAT.arnold_winther.ArnoldWintherNC(ref_el, degree=2)[source]¶
Bases:
CiarletElement
The definition of the nonconforming Arnold-Winther element.
FIAT.barycentric_interpolation module¶
- class FIAT.barycentric_interpolation.LagrangeLineExpansionSet(*args, **kwargs)[source]¶
Bases:
LineExpansionSet
Lagrange polynomial expansion set for given points the line.
- class FIAT.barycentric_interpolation.LagrangePolynomialSet(ref_el, pts, shape=())[source]¶
Bases:
PolynomialSet
- FIAT.barycentric_interpolation.barycentric_interpolation(nodes, wts, dmat, pts, order=0)[source]¶
Evaluates a Lagrange basis on a line reference element via the second barycentric interpolation formula. See Berrut and Trefethen (2004) https://doi.org/10.1137/S0036144502417715 Eq. (4.2) & (9.4)
FIAT.bell module¶
- class FIAT.bell.Bell(ref_el)[source]¶
Bases:
CiarletElement
The Bell finite element.
FIAT.bernardi_raugel module¶
- class FIAT.bernardi_raugel.BernardiRaugel(ref_el, order=1)[source]¶
Bases:
CiarletElement
The Bernardi-Raugel (extended) element.
This element does not belong to a Stokes complex, but can be paired with DG_{k-1}. This pair is inf-sup stable, but only weakly divergence-free.
FIAT.bernstein module¶
- class FIAT.bernstein.Bernstein(ref_el, degree)[source]¶
Bases:
FiniteElement
A finite element with Bernstein polynomials as basis functions.
- tabulate(order, points, entity=None)[source]¶
Return tabulated values of derivatives up to given order of basis functions at given points.
- Parameters:
order – The maximum order of derivative.
points – An iterable of points.
entity – Optional (dimension, entity number) pair indicating which topological entity of the reference element to tabulate on. If
None
, default cell-wise tabulation is performed.
- class FIAT.bernstein.BernsteinDualSet(ref_el, degree)[source]¶
Bases:
DualSet
The dual basis for Bernstein elements.
- FIAT.bernstein.bernstein_Dx(points, ks, order, R2B)[source]¶
Evaluates Bernstein polynomials or its derivatives according to reference coordinates.
- Parameters:
points – array of points in BARYCENTRIC COORDINATES
ks – exponents defining the Bernstein polynomial
alpha – derivative order (returns all derivatives of this specified order)
R2B – linear mapping from reference to barycentric coordinates
- Returns:
dictionary mapping from derivative tuples to arrays of Bernstein polynomial values at given points.
- FIAT.bernstein.bernstein_db(points, ks, alpha=None)[source]¶
Evaluates Bernstein polynomials or its derivative at barycentric points.
- Parameters:
points – array of points in barycentric coordinates
ks – exponents defining the Bernstein polynomial
alpha – derivative tuple
- Returns:
array of Bernstein polynomial values at given points.
FIAT.brezzi_douglas_fortin_marini module¶
- class FIAT.brezzi_douglas_fortin_marini.BrezziDouglasFortinMarini(ref_el, degree)[source]¶
Bases:
CiarletElement
The BDFM element
FIAT.brezzi_douglas_marini module¶
- class FIAT.brezzi_douglas_marini.BDMDualSet(ref_el, degree, variant, interpolant_deg)[source]¶
Bases:
DualSet
- class FIAT.brezzi_douglas_marini.BrezziDouglasMarini(ref_el, degree, variant=None)[source]¶
Bases:
CiarletElement
The BDM element
- Parameters:
ref_el – The reference element.
degree – The degree.
variant – optional variant specifying the types of nodes.
variant can be chosen from [“point”, “integral”, “integral(q)”] “point” -> dofs are evaluated by point evaluation. Note that this variant has suboptimal convergence order in the H(div)-norm “integral” -> dofs are evaluated by quadrature rules with the minimum degree required for unisolvence. “integral(q)” -> dofs are evaluated by quadrature rules with the minimum degree required for unisolvence plus q. You might want to choose a high quadrature degree to make sure that expressions will be interpolated exactly. This is important when you want to have (nearly) div-preserving interpolation.
FIAT.brezzi_douglas_marini_cube module¶
- class FIAT.brezzi_douglas_marini_cube.BrezziDouglasMariniCube(ref_el, degree, mapping)[source]¶
Bases:
FiniteElement
The Brezzi-Douglas-Marini element on quadrilateral cells.
- Parameters:
ref_el – The reference element.
k – The degree.
mapping – A string giving the Piola mapping. Either ‘contravariant Piola’ or ‘covariant Piola’.
- entity_closure_dofs()[source]¶
Return the map of topological entities to degrees of freedom on the closure of those entities for the finite element.
- entity_dofs()[source]¶
Return the map of topological entities to degrees of freedom for the finite element.
- tabulate(order, points, entity=None)[source]¶
Return tabulated values of derivatives up to a given order of basis functions at given points.
- Parameters:
order – The maximum order of derivative.
points – An iterable of points.
entity – Optional (dimension, entity number) pair indicating which topological entity of the reference element to tabulate on. If
None
, tabulated values are computed by geometrically approximating which facet the points are on.
- class FIAT.brezzi_douglas_marini_cube.BrezziDouglasMariniCubeEdge(ref_el, degree)[source]¶
Bases:
BrezziDouglasMariniCube
The Brezzi-Douglas-Marini HCurl element on quadrilateral cells.
- Parameters:
ref_el – The reference element.
k – The degree.
- class FIAT.brezzi_douglas_marini_cube.BrezziDouglasMariniCubeFace(ref_el, degree)[source]¶
Bases:
BrezziDouglasMariniCube
The Brezzi-Douglas-Marini HDiv element on quadrilateral cells.
- Parameters:
ref_el – The reference element.
k – The degree.
- FIAT.brezzi_douglas_marini_cube.bdmce_edge_basis(deg, dx, dy, x_mid, y_mid)[source]¶
Returns the basis functions associated with DoFs on the edges of elements for the HCurl Brezz-Douglas-Marini element on quadrilateral cells.
These were introduced by Brezzi, Douglas, Marini (1985) “Two families of mixed finite elements for Second Order Elliptic Problems”
Following, e.g. Brezzi, Douglas, Fortin, Marini (1987) “Efficient rectangular mixed finite elements in two and three space variables” For rectangle K and degree j: BDM_j(K) = [P_j(K)^2 + Span(curl(xy^{j+1}, x^{j+1}y))] x P_{j-1}(K)
The resulting basis functions all have a curl whose polynomials are of degree (j - 1).
- Parameters:
deg – The element degree.
dx – A tuple of sympy expressions, expanding the interval in the x direction. Probably (1-x, x).
dy – A tuple of sympy expressions, expanding the interval in the y direction. Probably (1-y, y).
x_mid – A sympy expression, probably 2*x-1.
y_mid – A sympy expression, probably 2*y-1.
- FIAT.brezzi_douglas_marini_cube.bdmce_face_basis(deg, dx, dy, x_mid, y_mid)[source]¶
Returns the basis functions associated with DoFs on the faces of elements for the HCurl Brezz-Douglas-Marini element on quadrilateral cells.
These were introduced by Brezzi, Douglas, Marini (1985) “Two families of mixed finite elements for Second Order Elliptic Problems”
Following, e.g. Brezzi, Douglas, Fortin, Marini (1987) “Efficient rectangular mixed finite elements in two and three space variables” For rectangle K and degree j: BDM_j(K) = [P_j(K)^2 + Span(curl(xy^{j+1}, x^{j+1}y))] x P_{j-1}(K)
The resulting basis functions all have a curl whose polynomials are of degree (j - 1).
- Parameters:
deg – The element degree.
dx – A tuple of sympy expressions, expanding the interval in the x direction. Probably (1-x, x).
dy – A tuple of sympy expressions, expanding the interval in the y direction. Probably (1-y, y).
x_mid – A sympy expression, probably 2*x-1.
y_mid – A sympy expression, probably 2*y-1.
- FIAT.brezzi_douglas_marini_cube.construct_bdmce_basis(ref_el, degree)[source]¶
Return the basis functions for a particular BDMCE space as a list.
- Parameters:
ref_el – The reference element.
k – The degree.
- FIAT.brezzi_douglas_marini_cube.numpy_lambdify(X, F, modules='numpy', dummify=False)[source]¶
Unfortunately, SymPy’s own lambdify() doesn’t work well with NumPy in that simple functions like
lambda x: 1.0,
when evaluated with NumPy arrays, return just “1.0” instead of an array of 1s with the same shape as x. This function does that.
FIAT.bubble module¶
- class FIAT.bubble.Bubble(ref_el, degree)[source]¶
Bases:
CodimBubble
The bubble finite element: the dofs of the Lagrange FE in the interior of the cell
- class FIAT.bubble.CodimBubble(ref_el, degree, codim)[source]¶
Bases:
RestrictedElement
Bubbles of a certain codimension.
- class FIAT.bubble.FacetBubble(ref_el, degree)[source]¶
Bases:
CodimBubble
The facet bubble finite element: the dofs of the Lagrange FE in the interior of the facets
FIAT.check_format_variant module¶
- FIAT.check_format_variant.parse_lagrange_variant(variant, discontinuous=False, integral=False)[source]¶
Parses variant options for Lagrange elements.
variant may be a single option or comma-separated pair indicating the dof type (integral, equispaced, spectral, etc) and the type of splitting to give a macro-element (Alfeld, Powell-Sabin, iso)
FIAT.christiansen_hu module¶
- class FIAT.christiansen_hu.ChristiansenHu(ref_el, degree=1)[source]¶
Bases:
CiarletElement
The Christiansen-Hu C^0(Worsey-Farin) linear macroelement with divergence in P0. This element belongs to a Stokes complex, and is paired with unsplit DG0.
FIAT.crouzeix_raviart module¶
- class FIAT.crouzeix_raviart.CrouzeixRaviart(ref_el, degree, variant=None)[source]¶
Bases:
CiarletElement
The Crouzeix-Raviart finite element:
K: Triangle/Tetrahedron Polynomial space: P_k Dual basis: Evaluation at points or integral moments
FIAT.discontinuous module¶
- class FIAT.discontinuous.DiscontinuousElement(element)[source]¶
Bases:
CiarletElement
A copy of an existing element where all dofs are associated with the cell
- get_nodal_basis()[source]¶
Return the nodal basis, encoded as a PolynomialSet object, for the finite element.
- mapping()[source]¶
Return a list of appropriate mappings from the reference element to a physical element for each basis function of the finite element.
FIAT.discontinuous_lagrange module¶
- class FIAT.discontinuous_lagrange.BrokenLagrangeDualSet(ref_el, degree, point_variant='equispaced')[source]¶
Bases:
DualSet
The dual basis for Lagrange elements. This class works for simplices of any dimension. Nodes are point evaluation at equispaced points. This is the broken version where all nodes are topologically associated with the cell itself
- class FIAT.discontinuous_lagrange.DiscontinuousLagrange(ref_el, degree, variant='equispaced')[source]¶
Bases:
CiarletElement
The discontinuous Lagrange finite element.
- Parameters:
ref_el – The reference element, which could be a standard FIAT simplex or a split complex
degree – The polynomial degree
variant –
A comma-separated string that may specify the type of point distribution and the splitting strategy if a macro element is desired. Either option may be omitted. The default point type is equispaced and the default splitting strategy is None. Example: variant=’gl’ gives a standard unsplit point distribution with
spectral points.
variant=’equispaced,Iso(2)’ with degree=1 gives the P2:P1 iso element. variant=’Alfeld’ can be used to obtain a barycentrically refined
macroelement for Scott-Vogelius.
FIAT.discontinuous_pc module¶
- class FIAT.discontinuous_pc.DPC0(ref_el)[source]¶
Bases:
CiarletElement
- class FIAT.discontinuous_pc.DPCDualSet(ref_el, flat_el, degree)[source]¶
Bases:
DualSet
The dual basis for DPC elements. This class works for hypercubes of any dimension. Nodes are point evaluation at equispaced points. This is the discontinuous version where all nodes are topologically associated with the cell itself
- class FIAT.discontinuous_pc.HigherOrderDPC(ref_el, degree)[source]¶
Bases:
CiarletElement
The DPC finite element. It is what it is.
FIAT.discontinuous_raviart_thomas module¶
- class FIAT.discontinuous_raviart_thomas.DRTDualSet(ref_el, degree)[source]¶
Bases:
DualSet
Dual basis for Raviart-Thomas elements consisting of point evaluation of normals on facets of codimension 1 and internal moments against polynomials. This is the discontinuous version where all nodes are topologically associated with the cell itself
- class FIAT.discontinuous_raviart_thomas.DiscontinuousRaviartThomas(ref_el, degree)[source]¶
Bases:
CiarletElement
The discontinuous Raviart-Thomas finite element
FIAT.discontinuous_taylor module¶
- class FIAT.discontinuous_taylor.DiscontinuousTaylorDualSet(ref_el, degree)[source]¶
Bases:
DualSet
The dual basis for Taylor elements. This class works for intervals. Nodes are function and derivative evaluation at the midpoint.
- class FIAT.discontinuous_taylor.HigherOrderDiscontinuousTaylor(ref_el, degree)[source]¶
Bases:
CiarletElement
The discontinuous Taylor finite element. Use a Taylor basis for DG.
FIAT.dual_set module¶
- class FIAT.dual_set.DualSet(nodes, ref_el, entity_ids, entity_permutations=None)[source]¶
Bases:
object
- get_entity_permutations()[source]¶
This method returns a nested dictionary that gives, for each dimension, for each entity, and for each possible entity orientation, the DoF permutation array that maps the entity local DoF ordering to the canonical global DoF ordering; see
Simplex
andUFCQuadrilateral
for how we define entity orientations for standard cells.The entity permutations dict for the degree 4 Lagrange finite element on the interval, for instance, is given by:
{0: {0: {0: [0]}, 1: {0: [0]}}, 1: {0: {0: [0, 1, 2], 1: [2, 1, 0]}}}
Note that there are two entities on dimension
0
(vertices), each of which has only one possible orientation, while there is a single entity on dimension1
(interval), which has two possible orientations representing non-reflected and reflected intervals.
- get_indices(restriction_domain, take_closure=True)[source]¶
Returns the list of dofs with support on a given restriction domain.
- Parameters:
restriction_domain – can be ‘interior’, ‘vertex’, ‘edge’, ‘face’ or ‘facet’
take_closure – Are we taking the closure of the restriction domain?
- to_riesz(poly_set)[source]¶
This method gives the action of the entire dual set on each member of the expansion set underlying poly_set. Then, applying the linear functionals of the dual set to an arbitrary polynomial in poly_set is accomplished by (generalized) matrix multiplication.
For scalar-valued spaces, this produces a matrix :math:R_{i, j} such that :math:ell_i(f) = sum_{j} a_j ell_i(phi_j) for :math:f=sum_{j} a_j phi_j.
More generally, it will have shape concatenating the number of functionals in the dual set, the value shape of functions it takes, and the number of members of the expansion set.
- FIAT.dual_set.lexsort_nodes(ref_el, nodes, entity=None, offset=0)[source]¶
Sort PointEvaluation nodes in lexicographical ordering.
FIAT.enriched module¶
- class FIAT.enriched.EnrichedElement(*elements)[source]¶
Bases:
FiniteElement
Class implementing a finite element that combined the degrees of freedom of two existing finite elements.
This is an implementation which does not care about orthogonality of primal and dual basis.
- get_nodal_basis()[source]¶
Return the nodal basis, encoded as a PolynomialSet object, for the finite element.
FIAT.expansions module¶
Principal orthogonal expansion functions as defined by Karniadakis and Sherwin. These are parametrized over a reference element so as to allow users to get coordinates that they want.
- FIAT.expansions.C0_basis(dim, n, tabulations)[source]¶
Modify a tabulation of a hierarchical basis to enforce C0-continuity.
- Parameters:
dim – The spatial dimension of the simplex.
n – The polynomial degree.
tabulations – An iterable tabulations of the hierarchical basis.
- Returns:
A tuple of tabulations of the C0 basis.
- class FIAT.expansions.ExpansionSet(*args, **kwargs)[source]¶
Bases:
object
- get_dmats(degree, cell=0)[source]¶
Returns a numpy array with the expansion coefficients dmat[k, j, i] of the gradient of each member of the expansion set:
d/dx_k phi_j = sum_i dmat[k, j, i] phi_i.
- tabulate_jumps(n, points, order=0)[source]¶
Tabulates derivative jumps on given points.
- Parameters:
n – the polynomial degree.
points – an iterable of points on the cell complex.
order – the order of differentiation.
- Returns:
a dictionary of tabulations of derivative jumps across interior facets.
- tabulate_normal_jumps(n, ref_pts, facet, order=0)[source]¶
Tabulates the normal derivative jumps on reference points on a facet.
- Parameters:
n – the polynomial degree.
ref_pts – an iterable of points on the reference facet.
facet – the facet id.
order – the order of differentiation.
- Returns:
a numpy array of tabulations of normal derivative jumps.
- class FIAT.expansions.LineExpansionSet(*args, **kwargs)[source]¶
Bases:
ExpansionSet
Evaluates the Legendre basis on a line reference element.
- class FIAT.expansions.PointExpansionSet(*args, **kwargs)[source]¶
Bases:
ExpansionSet
Evaluates the point basis on a point reference element.
- class FIAT.expansions.TetrahedronExpansionSet(*args, **kwargs)[source]¶
Bases:
ExpansionSet
Collapsed orthonormal polynomial expansion on a tetrahedron.
- class FIAT.expansions.TriangleExpansionSet(*args, **kwargs)[source]¶
Bases:
ExpansionSet
Evaluates the orthonormal Dubiner basis on a triangular reference element.
- FIAT.expansions.compute_cell_point_map(ref_el, pts, unique=True, tol=1e-12)[source]¶
Maps cells on a simplicial complex to points. Points outside the complex are binned to the nearest cell.
- Parameters:
ref_el – a SimplicialComplex.
pts – an iterable of physical points on the complex.
unique – Are we assigning a unique cell to points on facets?
tol – the absolute tolerance.
- Returns:
a dict mapping cell id to the point ids nearest to that cell.
- FIAT.expansions.compute_partition_of_unity(ref_el, pt, unique=True, tol=1e-12)[source]¶
Computes the partition of unity functions for each subcell.
- Parameters:
ref_el – a SimplicialComplex.
pt – a physical point on the complex.
unique – Are we assigning a unique cell to points on facets?
tol – the absolute tolerance.
- Returns:
a list of (weighted) characteristic functions for each subcell.
- FIAT.expansions.dubiner_recurrence(dim, n, order, ref_pts, Jinv, scale, variant=None)[source]¶
Tabulate a Dubiner expansion set using the recurrence from (Kirby 2010).
- Parameters:
dim – The spatial dimension of the simplex.
n – The polynomial degree.
order – The maximum order of differentiation.
ref_pts – An
ndarray
with the coordinates on the default (-1, 1)^d simplex.Jinv – The inverse of the Jacobian of the coordinate mapping from the default simplex.
scale – A scale factor that sets the first member of expansion set.
variant – Choose between the default (None) orthogonal basis, ‘bubble’ for integrated Jacobi polynomials, or ‘dual’ for the L2-duals of the integrated Jacobi polynomials.
- Returns:
A tuple with tabulations of the expansion set and its derivatives.
- FIAT.expansions.pad_coordinates(ref_pts, embedded_dim)[source]¶
Pad reference coordinates by appending -1.0.
- FIAT.expansions.pad_jacobian(A, embedded_dim)[source]¶
Pad coordinate mapping Jacobian by appending zero rows.
- FIAT.expansions.polynomial_cell_node_map(ref_el, n, continuity=None)[source]¶
Maps cells on a simplicial complex to members of a polynomial basis.
- Parameters:
ref_el – a SimplicialComplex.
n – the polynomial degree of the expansion set.
continuity – the continuity of the expansion set.
- Returns:
a numpy array mapping cell id to basis functions supported on that cell.
- FIAT.expansions.polynomial_dimension(ref_el, n, continuity=None)[source]¶
Returns the dimension of the space of polynomials of degree no greater than n on the reference complex.
- FIAT.expansions.polynomial_entity_ids(ref_el, n, continuity=None)[source]¶
Maps entites of a cell complex to members of a polynomial basis.
- Parameters:
ref_el – a SimplicialComplex.
n – the polynomial degree of the expansion set.
continuity – the continuity of the expansion set.
- Returns:
a dict of dicts mapping dimension and entity id to basis functions.
FIAT.fdm_element module¶
- class FIAT.fdm_element.FDMBrokenH1(ref_el, degree)[source]¶
Bases:
FDMFiniteElement
1D DG element with shape functions that diagonalize the Laplacian.
- class FIAT.fdm_element.FDMBrokenL2(ref_el, degree)[source]¶
Bases:
FDMFiniteElement
1D DG element with the derivates of DG FDM shape functions.
- class FIAT.fdm_element.FDMDiscontinuousLagrange(ref_el, degree)[source]¶
Bases:
FDMFiniteElement
1D DG element with derivatives of interior CG FDM shape functions.
- class FIAT.fdm_element.FDMDual(ref_el, degree, bc_order=1, formdegree=0, orthogonalize=False)[source]¶
Bases:
DualSet
The dual basis for 1D elements with FDM shape functions.
- class FIAT.fdm_element.FDMFiniteElement(ref_el, degree)[source]¶
Bases:
CiarletElement
1D element that diagonalizes bilinear forms with BCs.
- class FIAT.fdm_element.FDMHermite(ref_el, degree)[source]¶
Bases:
FDMFiniteElement
1D CG element with interior shape functions that diagonalize the biharmonic operator.
- class FIAT.fdm_element.FDMLagrange(ref_el, degree)[source]¶
Bases:
FDMFiniteElement
1D CG element with interior shape functions that diagonalize the Laplacian.
- class FIAT.fdm_element.FDMQuadrature(ref_el, degree)[source]¶
Bases:
FDMFiniteElement
1D DG element with interior CG FDM shape functions and orthogonalized vertex modes.
FIAT.finite_element module¶
- class FIAT.finite_element.CiarletElement(poly_set, dual, order, formdegree=None, mapping='affine', ref_complex=None)[source]¶
Bases:
FiniteElement
Class implementing Ciarlet’s abstraction of a finite element being a domain, function space, and set of nodes.
Elements derived from this class are nodal finite elements, with a nodal basis generated from polynomials encoded in a PolynomialSet.
- get_nodal_basis()[source]¶
Return the nodal basis, encoded as a PolynomialSet object, for the finite element.
- static is_nodal()[source]¶
True if primal and dual bases are orthogonal. If false, dual basis is not implemented or is undefined.
All implementations/subclasses are nodal including this one.
- tabulate(order, points, entity=None)[source]¶
Return tabulated values of derivatives up to given order of basis functions at given points.
- Parameters:
order – The maximum order of derivative.
points – An iterable of points.
entity – Optional (dimension, entity number) pair indicating which topological entity of the reference element to tabulate on. If
None
, default cell-wise tabulation is performed.
- class FIAT.finite_element.FiniteElement(ref_el, dual, order, formdegree=None, mapping='affine', ref_complex=None)[source]¶
Bases:
object
Class implementing a basic abstraction template for general finite element families. Finite elements which inherit from this class are non-nodal unless they are CiarletElement subclasses.
- entity_closure_dofs()[source]¶
Return the map of topological entities to degrees of freedom on the closure of those entities for the finite element.
- entity_dofs()[source]¶
Return the map of topological entities to degrees of freedom for the finite element.
- entity_permutations()[source]¶
This method returns a nested dictionary that gives, for each dimension, for each entity, and for each possible entity orientation, the DoF permutation array that maps the entity local DoF ordering to the canonical global DoF ordering; see
Simplex
andUFCQuadrilateral
for how we define entity orientations for standard cells.The entity permutations dict for the degree 4 Lagrange finite element on the interval, for instance, is given by:
{0: {0: {0: [0]}, 1: {0: [0]}}, 1: {0: {0: [0, 1, 2], 1: [2, 1, 0]}}}
Note that there are two entities on dimension
0
(vertices), each of which has only one possible orientation, while there is a single entity on dimension1
(interval), which has two possible orientations representing non-reflected and reflected intervals.
- get_reference_complex()[source]¶
Return the reference element complex, which is either the reference element or its subdivision in the case of a macro element.
- static is_nodal()[source]¶
True if primal and dual bases are orthogonal. If false, dual basis is not implemented or is undefined.
Subclasses may not necessarily be nodal, unless it is a CiarletElement.
- mapping()[source]¶
Return a list of appropriate mappings from the reference element to a physical element for each basis function of the finite element.
- tabulate(order, points, entity=None)[source]¶
Return tabulated values of derivatives up to given order of basis functions at given points.
- Parameters:
order – The maximum order of derivative.
points – An iterable of points.
entity – Optional (dimension, entity number) pair indicating which topological entity of the reference element to tabulate on. If
None
, default cell-wise tabulation is performed.
FIAT.functional module¶
- class FIAT.functional.ComponentPointEvaluation(ref_el, comp, shp, x)[source]¶
Bases:
Functional
Class representing point evaluation of a particular component of a vector/tensor function at a particular point x.
- class FIAT.functional.FrobeniusIntegralMoment(ref_el, Q, f_at_qpts, nm=None)[source]¶
Bases:
IntegralMoment
- class FIAT.functional.Functional(ref_el, target_shape, pt_dict, deriv_dict, functional_type)[source]¶
Bases:
object
Abstract class representing a linear functional. All FIAT functionals are discrete in the sense that they are written as a weighted sum of (derivatives of components of) their argument evaluated at particular points.
- Parameters:
ref_el – a
Cell
target_shape – a tuple indicating the value shape of functions on the functional operates (e.g. if the function eats 2-vectors then target_shape is (2,) and if it eats scalars then target_shape is ()
pt_dict – A dict mapping points to lists of information about how the functional is evaluated. Each entry in the list takes the form of a tuple (wt, comp) so that (at least if the deriv_dict argument is empty), the functional takes the form \(\ell(f) = \sum_{q=1}^{N_q} \sum_{k=1}^{K_q} w^q_k f_{c_k}(x_q)\) where \(f_{c_k}\) indicates a particular vector or tensor component
deriv_dict – A dict that is similar to pt_dict, although the entries of each list are tuples (wt, alpha, comp) with alpha a tuple of nonnegative integers corresponding to the order of partial differentiation in each spatial direction.
functional_type – a string labeling the kind of functional this is.
- evaluate(f)[source]¶
Obsolete and broken functional evaluation.
To evaluate the functional, call it on the target function:
functional(function)
- get_point_dict()[source]¶
Returns the functional information, which is a dictionary mapping each point in the support of the functional to a list of pairs containing the weight and component.
- get_type_tag()[source]¶
Returns the type of function (e.g. point evaluation or normal component, which is probably handy for clients of FIAT
- to_riesz(poly_set)[source]¶
Constructs an array representation of the functional so that the functional may be applied to a function expressed in in terms of the expansion set underlying poly_set by means of contracting coefficients.
That is, poly_set will have members all expressed in the form \(p = \sum_{i} \alpha^i \phi_i\) where \(\{\phi_i\}_{i}\) is some orthonormal expansion set and \(\alpha^i\) are coefficients. Note: the orthonormal expansion set is always scalar-valued but if the members of poly_set are vector or tensor valued the \(\alpha^i\) will be scalars or vectors.
This function constructs a tensor \(R\) such that the contraction of \(R\) with the array of coefficients \(\alpha\) produces the effect of \(\ell(f)\)
In the case of scalar-value functions, \(R\) is just a vector of the same length as the expansion set, and \(R_i = \ell(\phi_i)\). For vector-valued spaces, \(R_{ij}\) will be \(\ell(e^i \phi_j)\) where \(e^i\) is the canonical unit vector nonzero only in one entry \(i\).
- class FIAT.functional.IntegralLegendreBidirectionalMoment(cell, s1, s2, entity, mom_deg, comp_deg, nm='')[source]¶
Bases:
IntegralLegendreDirectionalMoment
Moment of dot(s1, dot(tau, s2)) against Legendre on entity, multiplied by the size of the reference facet
- class FIAT.functional.IntegralLegendreDirectionalMoment(cell, s, entity, mom_deg, quad_deg, nm='')[source]¶
Bases:
FrobeniusIntegralMoment
Moment of v.s against a Legendre polynomial over an edge
- class FIAT.functional.IntegralLegendreNormalMoment(cell, entity, mom_deg, comp_deg)[source]¶
Bases:
IntegralLegendreDirectionalMoment
Moment of v.n against a Legendre polynomial over an edge
- class FIAT.functional.IntegralLegendreNormalNormalMoment(cell, entity, mom_deg, comp_deg)[source]¶
Bases:
IntegralLegendreBidirectionalMoment
Moment of dot(n, dot(tau, n)) against Legendre on entity.
- class FIAT.functional.IntegralLegendreNormalTangentialMoment(cell, entity, mom_deg, comp_deg)[source]¶
Bases:
IntegralLegendreBidirectionalMoment
Moment of dot(n, dot(tau, t)) against Legendre on entity.
- class FIAT.functional.IntegralLegendreTangentialMoment(cell, entity, mom_deg, comp_deg)[source]¶
Bases:
IntegralLegendreDirectionalMoment
Moment of v.t against a Legendre polynomial over an edge
- class FIAT.functional.IntegralLegendreTangentialTangentialMoment(cell, entity, mom_deg, comp_deg)[source]¶
Bases:
IntegralLegendreBidirectionalMoment
Moment of dot(t, dot(tau, t)) against Legendre on entity.
- class FIAT.functional.IntegralMoment(ref_el, Q, f_at_qpts, comp=(), shp=())[source]¶
Bases:
Functional
Functional representing integral of the input against some tabulated function f.
- Parameters:
ref_el – a
Cell
.Q – a
QuadratureRule
.f_at_qpts – an array tabulating the function f at the quadrature points.
comp – Optional argument indicating that only a particular component of the input function should be integrated against f
shp – Optional argument giving the value shape of input functions.
- class FIAT.functional.IntegralMomentOfDerivative(ref_el, s, Q, f_at_qpts, comp=(), shp=())[source]¶
Bases:
Functional
Functional giving directional derivative integrated against some function on a facet.
- class FIAT.functional.IntegralMomentOfDivergence(ref_el, Q, f_at_qpts)[source]¶
Bases:
Functional
Functional representing integral of the divergence of the input against some tabulated function f.
- class FIAT.functional.IntegralMomentOfEdgeTangentEvaluation(ref_el, Q, P_at_qpts, edge)[source]¶
Bases:
Functional
int_e vcdot t p ds
p in Polynomials
- Parameters:
ref_el – reference element for which e is a dim-1 entity
Q – quadrature rule on the face
P_at_qpts – polynomials evaluated at quad points
edge – which edge.
- class FIAT.functional.IntegralMomentOfFaceTangentEvaluation(ref_el, Q, P_at_qpts, facet)[source]¶
Bases:
Functional
int_F v times n cdot p ds
p in Polynomials
- Parameters:
ref_el – reference element for which F is a codim-1 entity
Q – quadrature rule on the face
P_at_qpts – polynomials evaluated at quad points
facet – which facet.
- class FIAT.functional.IntegralMomentOfNormalDerivative(ref_el, facet_no, Q, f_at_qpts)[source]¶
Bases:
Functional
Functional giving normal derivative integrated against some function on a facet.
- class FIAT.functional.IntegralMomentOfNormalEvaluation(ref_el, Q, P_at_qpts, facet)[source]¶
Bases:
Functional
int_F vcdot n p ds p in Polynomials :arg ref_el: reference element for which F is a codim-1 entity :arg Q: quadrature rule on the face :arg P_at_qpts: polynomials evaluated at quad points :arg facet: which facet.
- class FIAT.functional.IntegralMomentOfScaledNormalEvaluation(ref_el, Q, P_at_qpts, facet)[source]¶
Bases:
Functional
int_F vcdot n p ds
p in Polynomials
- Parameters:
ref_el – reference element for which F is a codim-1 entity
Q – quadrature rule on the face
P_at_qpts – polynomials evaluated at quad points
facet – which facet.
- class FIAT.functional.IntegralMomentOfTangentialEvaluation(ref_el, Q, P_at_qpts, facet)[source]¶
Bases:
Functional
int_F vcdot n p ds p in Polynomials :arg ref_el: reference element for which F is a codim-1 entity :arg Q: quadrature rule on the face :arg P_at_qpts: polynomials evaluated at quad points :arg facet: which facet.
- class FIAT.functional.IntegralMomentOfTensorDivergence(ref_el, Q, f_at_qpts)[source]¶
Bases:
Functional
Like IntegralMomentOfDivergence, but on symmetric tensors.
- class FIAT.functional.PointDerivative(ref_el, x, alpha)[source]¶
Bases:
Functional
Class representing point partial differentiation of scalar functions at a particular point x.
- class FIAT.functional.PointDirectionalDerivative(ref_el, s, pt, comp=(), shp=(), nm=None)[source]¶
Bases:
Functional
Represents d/ds at a point.
- class FIAT.functional.PointDivergence(ref_el, x)[source]¶
Bases:
Functional
Class representing point divergence of vector functions at a particular point x.
- class FIAT.functional.PointEdgeTangentEvaluation(ref_el, edge_no, pt)[source]¶
Bases:
Functional
Implements the evaluation of the tangential component of a vector at a point on a facet of dimension 1.
- class FIAT.functional.PointEvaluation(ref_el, x)[source]¶
Bases:
Functional
Class representing point evaluation of scalar functions at a particular point x.
- class FIAT.functional.PointFaceTangentEvaluation(ref_el, face_no, tno, pt)[source]¶
Bases:
Functional
Implements the evaluation of a tangential component of a vector at a point on a facet of codimension 1.
- class FIAT.functional.PointNormalDerivative(ref_el, facet_no, pt, comp=(), shp=())[source]¶
Bases:
PointDirectionalDerivative
Represents d/dn at a point on a facet.
- class FIAT.functional.PointNormalEvaluation(ref_el, facet_no, pt)[source]¶
Bases:
Functional
Implements the evaluation of the normal component of a vector at a point on a facet of codimension 1.
- class FIAT.functional.PointNormalSecondDerivative(ref_el, facet_no, pt, comp=(), shp=())[source]¶
Bases:
PointSecondDerivative
Represents d^/dn^2 at a point on a facet.
- class FIAT.functional.PointScaledNormalEvaluation(ref_el, facet_no, pt)[source]¶
Bases:
Functional
Implements the evaluation of the normal component of a vector at a point on a facet of codimension 1, where the normal is scaled by the volume of that facet.
- class FIAT.functional.PointSecondDerivative(ref_el, s1, s2, pt, comp=(), shp=(), nm=None)[source]¶
Bases:
Functional
Represents d/ds1 d/ds2 at a point.
- class FIAT.functional.PointTangentialDerivative(ref_el, edge_no, pt, comp=(), shp=())[source]¶
Bases:
PointDirectionalDerivative
Represents d/dt at a point on an edge.
- class FIAT.functional.PointTangentialSecondDerivative(ref_el, edge_no, pt, comp=(), shp=())[source]¶
Bases:
PointSecondDerivative
Represents d^/dt^2 at a point on an edge.
- class FIAT.functional.PointwiseInnerProductEvaluation(ref_el, v, w, pt)[source]¶
Bases:
Functional
This is a functional on symmetric 2-tensor fields. Let u be such a field, p be a point, and v,w be vectors. This implements the evaluation v^T u(p) w.
Clearly v^iu_{ij}w^j = u_{ij}v^iw^j. Thus the value can be computed from the Frobenius inner product of u with wv^T. This gives the correct weights.
- class FIAT.functional.TensorBidirectionalIntegralMoment(ref_el, v, w, Q, f_at_qpts)[source]¶
Bases:
FrobeniusIntegralMoment
This is a functional on symmetric 2-tensor fields. Let u be such a field, f a function tabulated at points, and v,w be vectors. This implements the evaluation int v^T u(x) w f(x). Clearly v^iu_{ij}w^j = u_{ij}v^iw^j. Thus the value can be computed from the Frobenius inner product of u with vw^T. This gives the correct weights.
FIAT.gauss_legendre module¶
- class FIAT.gauss_legendre.GaussLegendre(ref_el, degree, variant='equispaced')[source]¶
Bases:
DiscontinuousLagrange
Simplicial discontinuous element with nodes at the (recursive) Gauss-Legendre points.
FIAT.gauss_lobatto_legendre module¶
FIAT.gauss_radau module¶
- class FIAT.gauss_radau.GaussRadau(ref_el, degree)[source]¶
Bases:
CiarletElement
1D discontinuous element with nodes at the Gauss-Radau points.
FIAT.gopalakrishnan_lederer_schoberl module¶
- FIAT.gopalakrishnan_lederer_schoberl.GopalakrishnanLedererSchoberlFirstKind(ref_el, degree)[source]¶
The GLS element used for the Mass-Conserving mixed Stress (MCS) formulation for Stokes flow.
GLS^1(k) is the space of trace-free polynomials of degree k with continuous normal-tangential components of degree k-1.
Reference: https://doi.org/10.1093/imanum/drz022
- class FIAT.gopalakrishnan_lederer_schoberl.GopalakrishnanLedererSchoberlSecondKind(ref_el, degree)[source]¶
Bases:
CiarletElement
The GLS element used for the Mass-Conserving mixed Stress (MCS) formulation for Stokes flow with weakly imposed stress symmetry.
GLS^2(k) is the space of trace-free polynomials of degree k with continuous normal-tangential components.
Reference: https://doi.org/10.1137/19M1248960
Notes¶
This element does not include the bubbles required for inf-sup stability of the weak symmetry constraint.
FIAT.guzman_neilan module¶
- class FIAT.guzman_neilan.GuzmanNeilanFirstKindH1(ref_el, order=1)[source]¶
Bases:
GuzmanNeilanH1
The Guzman-Neilan H1-conforming (extended) macroelement of the first kind.
Reference element: a simplex of any dimension. Function space: Pk^d + normal facet bubbles with div in P0, with 1 <= k < dim. Degrees of freedom: evaluation at Pk lattice points, and normal moments on faces.
This element belongs to a Stokes complex, and is paired with unsplit DG_{k-1}.
- class FIAT.guzman_neilan.GuzmanNeilanH1(ref_el, order=1, kind=1)[source]¶
Bases:
CiarletElement
The Guzman-Neilan H1-conforming (extended) macroelement.
- FIAT.guzman_neilan.GuzmanNeilanH1div(ref_el, degree=2, reduced=False)[source]¶
The Guzman-Neilan H1(div)-conforming (extended) macroelement.
Reference element: a simplex of any dimension. Function space: C0 P2^d(Alfeld) with C0 P1 divergence + normal facet bubbles with div in P0. Degrees of freedom: evaluation at P2(Alfeld) lattice points, divergence at P1 lattice points,
and normal moments on faces.
This element belongs to a Stokes complex, and is paired with CG1(Alfeld).
- class FIAT.guzman_neilan.GuzmanNeilanSecondKindH1(ref_el, order=1)[source]¶
Bases:
GuzmanNeilanH1
The Guzman-Neilan H1-conforming (extended) macroelement of the second kind.
Reference element: a simplex of any dimension. Function space: C0 Pk^d(Alfeld) + normal facet bubbles with div in P0, with 1 <= k < dim. Degrees of freedom: evaluation at Pk(Alfeld) lattice points, and normal moments on faces.
This element belongs to a Stokes complex, and is paired with DG_{k-1}(Alfeld).
- FIAT.guzman_neilan.GuzmanNeilanSpace(ref_el, order, kind=1, reduced=False)[source]¶
Return a basis for the (extended) Guzman-Neilan H1 space.
Project the extended Bernardi-Raugel space (Pk + FacetBubble)^d into C0 Pk(Alfeld)^d with P_{k-1} divergence, preserving its trace.
- Parameters:
ref_el – a simplex
order – the maximal polynomial degree
kind – kind = 1 gives Pk^d + GN bubbles, kind = 2 gives C0 Pk(Alfeld)^d + GN bubbles.
reduced – Include tangential bubbles if reduced = False.
- Returns:
a PolynomialSet basis for the Guzman-Neilan H1 space.
- FIAT.guzman_neilan.constant_div_projection(BR, C0, M, num_bubbles)[source]¶
Project the BR space into C0 Pk(Alfeld)^d with P_{k-1} divergence.
- FIAT.guzman_neilan.inner(v, u, qwts)[source]¶
Compute the L2 inner product from tabulation arrays and quadrature weights
FIAT.hct module¶
- class FIAT.hct.HsiehCloughTocher(ref_el, degree=3, reduced=False)[source]¶
Bases:
CiarletElement
The HCT macroelement. For degree higher than 3, we implement the super-smooth C^1 space from Groselj and Knez (2022) on a barycentric split, although there the basis functions are positive on an incenter split.
FIAT.hdiv_trace module¶
- class FIAT.hdiv_trace.HDivTrace(ref_el, degree, variant=None)[source]¶
Bases:
FiniteElement
Class implementing the trace of hdiv elements. This class is a stand-alone element family that produces a DG-facet field. This element is what’s produced after performing the trace operation on an existing H(Div) element.
This element is also known as the discontinuous trace field that arises in several DG formulations.
- get_nodal_basis()[source]¶
Return the nodal basis, encoded as a PolynomialSet object, for the finite element.
- static is_nodal()[source]¶
True if primal and dual bases are orthogonal. If false, dual basis is not implemented or is undefined.
Subclasses may not necessarily be nodal, unless it is a CiarletElement.
- tabulate(order, points, entity=None)[source]¶
Return tabulated values of derivatives up to a given order of basis functions at given points.
- Parameters:
order – The maximum order of derivative.
points – An iterable of points.
entity – Optional (dimension, entity number) pair indicating which topological entity of the reference element to tabulate on. If
None
, tabulated values are computed by geometrically approximating which facet the points are on.
Note
Performing illegal tabulations on this element will result in either a tabulation table of numpy.nan arrays (entity=None case), or insertions of the TraceError exception class. This is due to the fact that performing cell-wise tabulations, or asking for any order of derivative evaluations, are not mathematically well-defined.
- exception FIAT.hdiv_trace.TraceError(msg)[source]¶
Bases:
Exception
Exception caused by tabulating a trace element on the interior of a cell, or the gradient of a trace element.
- FIAT.hdiv_trace.barycentric_coordinates(points, vertices)[source]¶
Computes the barycentric coordinates for a set of points relative to a simplex defined by a set of vertices.
- Parameters:
points – A set of points.
vertices – A set of vertices that define the simplex.
- FIAT.hdiv_trace.construct_dg_element(ref_el, degree, variant)[source]¶
Constructs a discontinuous galerkin element of a given degree on a particular reference cell.
- FIAT.hdiv_trace.extract_unique_facet(coordinates, tolerance=1e-10)[source]¶
Determines whether a set of points (described in its barycentric coordinates) are all on one of the facet sub-entities, and return the particular facet and whether the search has been successful.
- Parameters:
coordinates – A set of points described in barycentric coordinates.
tolerance – A fixed tolerance for geometric identifications.
- FIAT.hdiv_trace.map_from_reference_facet(point, vertices)[source]¶
Evaluates the physical coordinate of a point using barycentric coordinates.
- Parameters:
point – The reference points to be mapped to the facet.
vertices – The vertices defining the physical element.
- FIAT.hdiv_trace.map_to_reference_facet(points, vertices, facet)[source]¶
Given a set of points and vertices describing a facet of a simplex in n-dimensional coordinates (where the points lie on the facet), map the points to the reference simplex of dimension (n-1).
- Parameters:
points – A set of points in n-D.
vertices – A set of vertices describing a facet of a simplex in n-D.
facet – Integer representing the facet number.
FIAT.hdivcurl module¶
FIAT.hellan_herrmann_johnson module¶
Implementation of the Hellan-Herrmann-Johnson finite elements.
- class FIAT.hellan_herrmann_johnson.HellanHerrmannJohnson(ref_el, degree=0, variant=None)[source]¶
Bases:
CiarletElement
The definition of Hellan-Herrmann-Johnson element. HHJ(k) is the space of symmetric-matrix-valued polynomials of degree k or less with normal-normal continuity.
FIAT.hermite module¶
- class FIAT.hermite.CubicHermite(ref_el, deg=3)[source]¶
Bases:
CiarletElement
The cubic Hermite finite element. It is what it is.
FIAT.hierarchical module¶
- class FIAT.hierarchical.IntegratedLegendre(ref_el, degree, variant=None)[source]¶
Bases:
CiarletElement
Simplicial continuous element with integrated Legendre polynomials.
- class FIAT.hierarchical.IntegratedLegendreDual(ref_el, degree)[source]¶
Bases:
DualSet
The dual basis for integrated Legendre elements.
- class FIAT.hierarchical.Legendre(ref_el, degree, variant=None)[source]¶
Bases:
CiarletElement
Simplicial discontinuous element with Legendre polynomials.
FIAT.hu_zhang module¶
Implementation of the Hu-Zhang finite elements.
- class FIAT.hu_zhang.HuZhang(ref_el, degree=3, variant=None)[source]¶
Bases:
CiarletElement
The definition of the Hu-Zhang element.
FIAT.jacobi module¶
Several functions related to the one-dimensional jacobi polynomials: Evaluation, evaluation of derivatives, plus computation of the roots via Newton’s method. These mainly are used in defining the expansion functions over the simplices and in defining quadrature rules over each domain.
- FIAT.jacobi.eval_jacobi(a, b, n, x)[source]¶
Evaluates the nth jacobi polynomial with weight parameters a,b at a point x. Recurrence relations implemented from the pseudocode given in Karniadakis and Sherwin, Appendix B
- FIAT.jacobi.eval_jacobi_batch(a, b, n, xs)[source]¶
Evaluates all jacobi polynomials with weights a,b up to degree n. xs is a numpy.array of points. Returns a two-dimensional array of tabulations, where the rows correspond to the Jacobi polynomials and the columns correspond to the points.
- FIAT.jacobi.eval_jacobi_deriv(a, b, n, x)[source]¶
Evaluates the first derivative of P_{n}^{a,b} at a point x.
- FIAT.jacobi.eval_jacobi_deriv_batch(a, b, n, xs)[source]¶
Evaluates the first derivatives of all jacobi polynomials with weights a,b up to degree n. xs is a numpy.array of points. Returns a two-dimensional array of points, where the rows correspond to the Jacobi polynomials and the columns correspond to the points.
FIAT.johnson_mercier module¶
- class FIAT.johnson_mercier.JohnsonMercier(ref_el, degree=1, variant=None)[source]¶
Bases:
CiarletElement
The Johnson-Mercier finite element.
FIAT.kong_mulder_veldhuizen module¶
- class FIAT.kong_mulder_veldhuizen.KongMulderVeldhuizen(ref_el, degree)[source]¶
Bases:
CiarletElement
The “lumped” simplical finite element (NB: requires custom quad. “KMV” points to achieve a diagonal mass matrix).
References¶
Higher-order triangular and tetrahedral finite elements with mass lumping for solving the wave equation M. J. S. CHIN-JOE-KONG, W. A. MULDER and M. VAN VELDHUIZEN
HIGHER-ORDER MASS-LUMPED FINITE ELEMENTS FOR THE WAVE EQUATION W.A. MULDER
NEW HIGHER-ORDER MASS-LUMPED TETRAHEDRAL ELEMENTS S. GEEVERS, W.A. MULDER, AND J.J.W. VAN DER VEGT
More Continuous Mass-Lumped Triangular Finite Elements W. A. MULDER
FIAT.lagrange module¶
- class FIAT.lagrange.Lagrange(ref_el, degree, variant='equispaced', sort_entities=False)[source]¶
Bases:
CiarletElement
The Lagrange finite element.
- Parameters:
ref_el – The reference element, which could be a standard FIAT simplex or a split complex
degree – The polynomial degree
variant –
A comma-separated string that may specify the type of point distribution and the splitting strategy if a macro element is desired. Either option may be omitted. The default point type is equispaced and the default splitting strategy is None. Example: variant=’gll’ gives a standard unsplit point distribution with
spectral points.
variant=’equispaced,Iso(2)’ with degree=1 gives the P2:P1 iso element. variant=’Alfeld’ can be used to obtain a barycentrically refined
macroelement for Scott-Vogelius.
sort_entities – A flag to sort entities by support vertex ids. If false then entities are sorted first by dimension and then by entity id. The DOFs are always sorted by the entity ordering and then lexicographically by lattice multiindex.
- class FIAT.lagrange.LagrangeDualSet(ref_el, degree, point_variant='equispaced', sort_entities=False)[source]¶
Bases:
DualSet
The dual basis for Lagrange elements. This class works for simplicial complexes of any dimension. Nodes are point evaluation at recursively-defined points.
- Parameters:
ref_el – The simplicial complex.
degree – The polynomial degree.
point_variant – The point distribution variant passed on to recursivenodes.
sort_entities – A flag to sort entities by support vertex ids. If false then entities are sorted first by dimension and then by entity id. The DOFs are always sorted by the entity ordering and then lexicographically by lattice multiindex.
FIAT.macro module¶
- class FIAT.macro.AlfeldSplit(ref_el)[source]¶
Bases:
PowellSabinSplit
Splits a simplicial complex by connecting cell vertices to their barycenter.
- class FIAT.macro.CkPolynomialSet(ref_el, degree, order=1, vorder=None, shape=(), **kwargs)[source]¶
Bases:
PolynomialSet
Constructs a C^k-continuous PolynomialSet on a simplicial complex.
- Parameters:
ref_el – The simplicial complex.
degree – The polynomial degree.
order – The order of continuity across facets or a dict of dicts mapping dimension and entity_id to the order of continuity at each entity.
vorder – The order of super-smoothness at interior vertices.
shape – The value shape.
variant – The variant for the underlying ExpansionSet.
scale – The scale for the underlying ExpansionSet.
- class FIAT.macro.HDivSymPolynomialSet(ref_el, degree, order=0, **kwargs)[source]¶
Bases:
PolynomialSet
- Constructs a symmetric tensor-valued PolynomialSet with continuous
normal components on a simplicial complex.
- Parameters:
ref_el – The simplicial complex.
degree – The polynomial degree.
order – The order of continuity across subcells.
variant – The variant for the underlying ExpansionSet.
scale – The scale for the underlying ExpansionSet.
- class FIAT.macro.IsoSplit(ref_el, degree=2, variant=None)[source]¶
Bases:
SplitSimplicialComplex
Splits simplex into the simplicial complex obtained by connecting points on a regular lattice.
- Parameters:
ref_el – The parent Simplex to split.
degree – The number of subdivisions along each edge of the simplex.
variant – The point distribution variant.
- class FIAT.macro.MacroQuadratureRule(ref_el, Q_ref, parent_facets=None)[source]¶
Bases:
QuadratureRule
Composite quadrature rule on parent facets that respects the splitting.
- Parameters:
ref_el – A simplicial complex.
Q_ref – A QuadratureRule on the reference simplex.
- Args parent_facets:
An iterable of facets of the same dimension as Q_ref, defaults to all facets.
- Returns:
A quadrature rule on the sub entities of the simplicial complex.
- class FIAT.macro.PowellSabin12Split(ref_el)[source]¶
Bases:
SplitSimplicialComplex
Splits a triangle (only!) by connecting each vertex to the opposite edge midpoint and edge midpoints to each other.
- class FIAT.macro.PowellSabinSplit(ref_el, dimension=1)[source]¶
Bases:
SplitSimplicialComplex
Splits a simplicial complex by connecting barycenters of subentities to the barycenter of the superentities of the given dimension or higher.
- class FIAT.macro.SplitSimplicialComplex(parent, vertices, topology)[source]¶
Bases:
SimplicialComplex
Abstract class to implement a split on a Simplex.
- Parameters:
parent – The parent Simplex to split.
vertices – The vertices of the simplicial complex.
topology – The topology of the simplicial complex.
- construct_subelement(dimension)[source]¶
Constructs the reference element of a cell subentity specified by subelement dimension.
- Parameters:
dimension – subentity dimension (integer)
- get_cell_connectivity()[source]¶
Connectitivity from cell in a complex to global facet ids and respects the entity numbering on the reference cell.
N.B. cell_connectivity[cell][dim] has the same contents as self.connectivity[(sd, dim)][cell], except those are sorted.
- get_interior_facets(dimension)[source]¶
Returns the list of entities of the given dimension that are supported on the parent’s interior.
- Parameters:
dimension – subentity dimension (integer)
- class FIAT.macro.WorseyFarinSplit(ref_el)[source]¶
Bases:
PowellSabinSplit
Splits a simplicial complex by connecting cell and facet vertices to their barycenter. This reduces to Powell-Sabin on the triangle, and Alfeld on the interval.
- FIAT.macro.bary_to_xy(verts, bary, result=None)[source]¶
Maps barycentric coordinates to physical points.
- Parameters:
verts – A tuple of points.
bary – A row-stacked numpy array of barycentric coordinates.
result – A row-stacked numpy array of physical points.
- Returns:
result
- FIAT.macro.facet_support(facet_coords, tol=1e-12)[source]¶
Returns the support of a facet.
- Parameters:
facet_coords – An iterable of tuples (barycentric coordinates) describing the facet.
- Returns:
A tuple of vertex ids where some coordinate is nonzero.
FIAT.mardal_tai_winther module¶
Implementation of the Mardal-Tai-Winther finite elements.
- class FIAT.mardal_tai_winther.MardalTaiWinther(ref_el, degree=3)[source]¶
Bases:
CiarletElement
The definition of the Mardal-Tai-Winther element.
FIAT.mixed module¶
- class FIAT.mixed.MixedElement(elements, ref_el=None)[source]¶
Bases:
FiniteElement
A FIAT-like representation of a mixed element.
- Parameters:
elements – An iterable of FIAT elements.
ref_el – The reference element (optional).
This object offers tabulation of the concatenated basis function tables along with an entity_dofs dict.
- mapping()[source]¶
Return a list of appropriate mappings from the reference element to a physical element for each basis function of the finite element.
FIAT.morley module¶
- class FIAT.morley.Morley(ref_el)[source]¶
Bases:
CiarletElement
The Morley finite element.
FIAT.nedelec module¶
- class FIAT.nedelec.Nedelec(ref_el, degree, variant=None)[source]¶
Bases:
CiarletElement
Nedelec finite element
- Parameters:
ref_el – The reference element.
degree – The degree.
variant – optional variant specifying the types of nodes.
variant can be chosen from [“point”, “integral”, “integral(q)”] “point” -> dofs are evaluated by point evaluation. Note that this variant has suboptimal convergence order in the H(curl)-norm “integral” -> dofs are evaluated by quadrature rules with the minimum degree required for unisolvence. “integral(q)” -> dofs are evaluated by quadrature rules with the minimum degree required for unisolvence plus q. You might want to choose a high quadrature degree to make sure that expressions will be interpolated exactly. This is important when you want to have (nearly) curl-preserving interpolation.
- class FIAT.nedelec.NedelecDual(ref_el, degree, variant, interpolant_deg)[source]¶
Bases:
DualSet
Dual basis for first-kind Nedelec.
FIAT.nedelec_second_kind module¶
- class FIAT.nedelec_second_kind.NedelecSecondKind(cell, degree, variant=None)[source]¶
Bases:
CiarletElement
The H(curl) Nedelec elements of the second kind on triangles and tetrahedra: the polynomial space described by the full polynomials of degree k, with a suitable set of degrees of freedom to ensure H(curl) conformity.
- Parameters:
ref_el – The reference element.
degree – The degree.
variant – optional variant specifying the types of nodes.
variant can be chosen from [“point”, “integral”, “integral(q)”] “point” -> dofs are evaluated by point evaluation. Note that this variant has suboptimal convergence order in the H(curl)-norm “integral” -> dofs are evaluated by quadrature rules with the minimum degree required for unisolvence. “integral(q)” -> dofs are evaluated by quadrature rules with the minimum degree required for unisolvence plus q. You might want to choose a high quadrature degree to make sure that expressions will be interpolated exactly. This is important when you want to have (nearly) curl-preserving interpolation.
- class FIAT.nedelec_second_kind.NedelecSecondKindDual(cell, degree, variant, interpolant_deg)[source]¶
Bases:
DualSet
This class represents the dual basis for the Nedelec H(curl) elements of the second kind. The degrees of freedom (L) for the elements of the k’th degree are
d = 2:
vertices: None
edges: L(f) = f (x_i) * t for (k+1) points x_i on each edge
cell: L(f) = int f * g * dx for g in RT_{k-1}
d = 3:
vertices: None
edges: L(f) = f(x_i) * t for (k+1) points x_i on each edge
faces: L(f) = int_F f * g * ds for g in RT_{k-1}(F) for each face F
cell: L(f) = int f * g * dx for g in RT_{k-2}
Higher spatial dimensions are not yet implemented. (For d = 1, these elements coincide with the CG_k elements.)
FIAT.nodal_enriched module¶
- class FIAT.nodal_enriched.NodalEnrichedElement(*elements)[source]¶
Bases:
CiarletElement
NodalEnriched element is a direct sum of a sequence of finite elements. Primal basis is reorthogonalized to the dual basis for nodality.
- The following is equivalent:
the constructor is well-defined,
the resulting element is unisolvent and its basis is nodal,
the supplied elements are unisolvent with nodal basis and their primal bases are mutually linearly independent,
the supplied elements are unisolvent with nodal basis and their dual bases are mutually linearly independent.
FIAT.orientation_utils module¶
- class FIAT.orientation_utils.Orientation[source]¶
Bases:
object
Base class representing unsigned integer orientations.
Orientations represented by this class are to be consistent with those used in make_entity_permutations_simplex and make_entity_permutations_tensorproduct.
- FIAT.orientation_utils.check_permutation_even_or_odd(_l)[source]¶
Check if the given permutation is even or odd relative to
[0, 1, ...)
.- Parameters:
_l – The permutation.
- Returns:
0
if even and1
if odd.
- FIAT.orientation_utils.make_entity_permutations_simplex(dim, npoints)[source]¶
Make orientation-permutation map for the given simplex dimension, dim, and the number of points along each axis
As an example, we first compute the orientation of a triangular cell:
|1 0 47 42 | | +–2—+ +–43–+
FIAT canonical Mapped example physical cell
Suppose that the facets of the physical cell are canonically ordered as:
C = [43, 42, 47]
FIAT facet to Physical facet map is given by:
M = [42, 47, 43]
Then the orientation of the cell is computed as:
C.index(M[0]) = 1; C.remove(M[0]) C.index(M[1]) = 1; C.remove(M[1]) C.index(M[2]) = 0; C.remove(M[2])
o = (1 * 2!) + (1 * 1!) + (0 * 0!) = 3
For npoints = 3, there are 6 DoFs:
5 0 3 4 1 3 0 1 2 2 4 5
FIAT canonical Physical cell canonical
The permutation associated with o = 3 then is:
[2, 4, 5, 1, 3, 0]
The output of this function contains one such permutation for each orientation for the given simplex dimension and the number of points along each axis.
- FIAT.orientation_utils.make_entity_permutations_tensorproduct(cells, dim, o_p_maps)[source]¶
Make orientation-permutation map for an entity of a tensor product cell.
- Parameters:
cells – List of cells composing the tensor product cell.
dim – List of (sub)dimensions of the component cells that makes the tensor product (sub)cell.
o_p_maps – List of orientation-permutation maps of the component (sub)cells.
- Returns:
The orientation-permutation map of the tensor product (sub)cell.
Example¶
from FIAT.reference_element import UFCInterval from FIAT.orientation_utils import make_entity_permutations_tensorproduct cells = [UFCInterval(), UFCInterval()] m = make_entity_permutations_tensorproduct(cells, [1, 0], [{0: [0, 1], 1: [1, 0]}, {0: [0]}]) print(m) # prints: # {(0, 0, 0): [0, 1], # (0, 1, 0): [1, 0]} m = make_entity_permutations_tensorproduct(cells, [1, 1], [{0: [0, 1], 1: [1, 0]}, {0: [0, 1], 1: [1, 0]}]) print(m) # prints: # {(0, 0, 0): [0, 1, 2, 3], # (0, 0, 1): [1, 0, 3, 2], # (0, 1, 0): [2, 3, 0, 1], # (0, 1, 1): [3, 2, 1, 0], # (1, 0, 0): [0, 2, 1, 3], # (1, 0, 1): [2, 0, 3, 1], # (1, 1, 0): [1, 3, 0, 2], # (1, 1, 1): [3, 1, 2, 0]}
FIAT.orthopoly module¶
orthopoly.py - A suite of functions for generating orthogonal polynomials and quadrature rules.
Copyright (c) 2014 Greg von Winckel All rights reserved.
Permission is hereby granted, free of charge, to any person obtaining a copy of this software and associated documentation files (the “Software”), to deal in the Software without restriction, including without limitation the rights to use, copy, modify, merge, publish, distribute, sublicense, and/or sell copies of the Software, and to permit persons to whom the Software is furnished to do so, subject to the following conditions:
The above copyright notice and this permission notice shall be included in all copies or substantial portions of the Software.
THE SOFTWARE IS PROVIDED “AS IS”, WITHOUT WARRANTY OF ANY KIND, EXPRESS OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY, FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER DEALINGS IN THE SOFTWARE.
Last updated on Wed Jan 1 14:29:25 MST 2014
Modified by David A. Ham (david.ham@imperial.ac.uk), 2016
- FIAT.orthopoly.gauss(alpha, beta)[source]¶
Compute the Gauss nodes and weights from the recursion coefficients associated with a set of orthogonal polynomials
Inputs: alpha - recursion coefficients beta - recursion coefficients
Outputs: x - quadrature nodes w - quadrature weights
Adapted from the MATLAB code by Walter Gautschi http://www.cs.purdue.edu/archives/2002/wxg/codes/gauss.m
- FIAT.orthopoly.jacobi(N, a, b, x, NOPT=1)[source]¶
JACOBI computes the Jacobi polynomials which are orthogonal on [-1,1] with respect to the weight w(x)=[(1-x)^a]*[(1+x)^b] and evaluate them on the given grid up to P_N(x). Setting NOPT=2 returns the L2-normalized polynomials
- FIAT.orthopoly.jacobiD(N, a, b, x, NOPT=1)[source]¶
JACOBID computes the first derivatives of the normalized Jacobi polynomials which are orthogonal on [-1,1] with respect to the weight w(x)=[(1-x)^a]*[(1+x)^b] and evaluate them on the given grid up to P_N(x). Setting NOPT=2 returns the derivatives of the L2-normalized polynomials
- FIAT.orthopoly.lobatto(alpha, beta, xl1, xl2)[source]¶
Compute the Lobatto nodes and weights with the preassigned nodea xl1,xl2
Inputs: alpha - recursion coefficients beta - recursion coefficients xl1 - assigned node location xl2 - assigned node location
Outputs: x - quadrature nodes w - quadrature weights
Based on the section 7 of the paper “Some modified matrix eigenvalue problems” by Gene Golub, SIAM Review Vol 15, No. 2, April 1973, pp.318–334
- FIAT.orthopoly.mm_log(N, a)[source]¶
MM_LOG Modified moments for a logarithmic weight function.
The call mm=MM_LOG(n,a) computes the first n modified moments of the logarithmic weight function w(t)=t^a log(1/t) on [0,1] relative to shifted Legendre polynomials.
- REFERENCE: Walter Gautschi,``On the preceding paper `A Legendre
polynomial integral’ by James L. Blue’’, Math. Comp. 33 (1979), 742-743.
Adapted from the MATLAB implementation: https://www.cs.purdue.edu/archives/2002/wxg/codes/mm_log.m
- FIAT.orthopoly.mod_chebyshev(N, mom, alpham, betam)[source]¶
Calcuate the recursion coefficients for the orthogonal polynomials which are are orthogonal with respect to a weight function which is represented in terms of its modifed moments which are obtained by integrating the monic polynomials against the weight function.
References¶
John C. Wheeler, “Modified moments and Gaussian quadratures” Rocky Mountain Journal of Mathematics, Vol. 4, Num. 2 (1974), 287–296
Walter Gautschi, “Orthogonal Polynomials (in Matlab) Journal of Computational and Applied Mathematics, Vol. 178 (2005) 215–234
Adapted from the MATLAB implementation: https://www.cs.purdue.edu/archives/2002/wxg/codes/chebyshev.m
- FIAT.orthopoly.polyval(alpha, beta, x)[source]¶
Evaluate polynomials on x given the recursion coefficients alpha and beta
- FIAT.orthopoly.rec_jaclog(N, a)[source]¶
Generate the recursion coefficients alpha_k, beta_k
P_{k+1}(x) = (x-alpha_k)*P_{k}(x) - beta_k P_{k-1}(x)
for the monic polynomials which are orthogonal on [0,1] with respect to the weight w(x)=x^a*log(1/x)
Inputs: N - polynomial order a - weight parameter
Outputs: alpha - recursion coefficients beta - recursion coefficients
Adated from the MATLAB code: https://www.cs.purdue.edu/archives/2002/wxg/codes/r_jaclog.m
- FIAT.orthopoly.rec_jacobi(N, a, b)[source]¶
Generate the recursion coefficients alpha_k, beta_k
P_{k+1}(x) = (x-alpha_k)*P_{k}(x) - beta_k P_{k-1}(x)
for the Jacobi polynomials which are orthogonal on [-1,1] with respect to the weight w(x)=[(1-x)^a]*[(1+x)^b]
Inputs: N - polynomial order a - weight parameter b - weight parameter
Outputs: alpha - recursion coefficients beta - recursion coefficients
Adapted from the MATLAB code by Dirk Laurie and Walter Gautschi http://www.cs.purdue.edu/archives/2002/wxg/codes/r_jacobi.m
- FIAT.orthopoly.rec_jacobi01(N, a, b)[source]¶
Generate the recursion coefficients alpha_k, beta_k for the Jacobi polynomials which are orthogonal on [0,1]
See rec_jacobi for the recursion coefficients on [-1,1]
Inputs: N - polynomial order a - weight parameter b - weight parameter
Outputs: alpha - recursion coefficients beta - recursion coefficients
Adapted from the MATLAB implementation: https://www.cs.purdue.edu/archives/2002/wxg/codes/r_jacobi01.m
FIAT.pointwise_dual module¶
- FIAT.pointwise_dual.compute_pointwise_dual(el, pts)[source]¶
Constructs a dual basis to the basis for el as a linear combination of a set of pointwise evaluations. This is useful when the prescribed finite element isn’t Ciarlet (e.g. the basis functions are provided explicitly as formulae). Alternately, the element’s given dual basis may involve differentiation, making run-time interpolation difficult in FIAT clients. The pointwise dual, consisting only of pointwise evaluations, will effectively replace these derivatives with (automatically determined) finite differences. This is exact on the polynomial space, but is an approximation if applied to functions outside the space.
- Parameters:
el – a
FiniteElement
.pts – an iterable of points with the same length as el’s dimension. These points must be unisolvent for the polynomial space
- Returns:
a :class DualSet
FIAT.polynomial_set module¶
- class FIAT.polynomial_set.ONPolynomialSet(ref_el, degree, shape=(), **kwargs)[source]¶
Bases:
PolynomialSet
Constructs an orthonormal basis out of expansion set by having an identity matrix of coefficients. Can be used to specify ON bases for vector- and tensor-valued sets as well.
- class FIAT.polynomial_set.ONSymTensorPolynomialSet(ref_el, degree, size=None, **kwargs)[source]¶
Bases:
PolynomialSet
Constructs an orthonormal basis for symmetric-tensor-valued polynomials on a reference element.
- class FIAT.polynomial_set.PolynomialSet(ref_el, degree, embedded_degree, expansion_set, coeffs)[source]¶
Bases:
object
Implements a set of polynomials as linear combinations of an expansion set over a reference element. ref_el: the reference element degree: an order labeling the space embedded degree: the degree of polynomial expansion basis that
must be used to evaluate this space
- coeffs: A numpy array containing the coefficients of the expansion
basis for each member of the set. Coeffs is ordered by coeffs[i,j,k] where i is the label of the member, k is the label of the expansion function, and j is a (possibly empty) tuple giving the index for a vector- or tensor-valued function.
- class FIAT.polynomial_set.TracelessTensorPolynomialSet(ref_el, degree, size=None, **kwargs)[source]¶
Bases:
PolynomialSet
Constructs an orthonormal basis for traceless-tensor-valued polynomials on a reference element.
- FIAT.polynomial_set.construct_new_coeffs(ref_el, A, B)[source]¶
Constructs new coefficients for the set union of A and B If A and B are discontinuous and do not have the same degree the smaller one is extended to match the larger.
This does not handle the case that A and B have continuity but not the same degree.
- FIAT.polynomial_set.make_bubbles(ref_el, degree, codim=0, shape=(), scale='L2 piola')[source]¶
Construct a polynomial set with codim bubbles up to the given degree.
- FIAT.polynomial_set.mis(m, n)[source]¶
Returns all m-tuples of nonnegative integers that sum up to n.
- FIAT.polynomial_set.polynomial_set_union_normalized(A, B)[source]¶
Given polynomial sets A and B, constructs a new polynomial set whose span is the same as that of span(A) union span(B). It may not contain any of the same members of the set, as we construct a span via SVD.
FIAT.powell_sabin module¶
- class FIAT.powell_sabin.QuadraticPowellSabin12(ref_el, degree=2)[source]¶
Bases:
CiarletElement
The PS12 macroelement is a C^1 quadratic macroelement defined on the 12-way Powell-Sabin split of a triangle.
- class FIAT.powell_sabin.QuadraticPowellSabin12DualSet(ref_complex, degree=2)[source]¶
Bases:
DualSet
- class FIAT.powell_sabin.QuadraticPowellSabin6(ref_el, degree=2)[source]¶
Bases:
CiarletElement
The PS6 macroelement is a C^1 quadratic macroelement defined on the 6-way Powell-Sabin split of a triangle.
FIAT.quadrature module¶
- class FIAT.quadrature.CollapsedQuadratureSimplexRule(ref_el, m)[source]¶
Bases:
QuadratureRule
Implements the collapsed quadrature rules defined in Karniadakis & Sherwin by mapping products of Gauss-Jacobi rules from the hypercube to the simplex.
- class FIAT.quadrature.CollapsedQuadratureTetrahedronRule(ref_el, m)[source]¶
Bases:
CollapsedQuadratureSimplexRule
Implements the collapsed quadrature rules defined in Karniadakis & Sherwin by mapping products of Gauss-Jacobi rules from the cube to the tetrahedron.
- class FIAT.quadrature.CollapsedQuadratureTriangleRule(ref_el, m)[source]¶
Bases:
CollapsedQuadratureSimplexRule
Implements the collapsed quadrature rules defined in Karniadakis & Sherwin by mapping products of Gauss-Jacobi rules from the square to the triangle.
- class FIAT.quadrature.FacetQuadratureRule(ref_el, entity_dim, entity_id, Q_ref)[source]¶
Bases:
QuadratureRule
A quadrature rule on a facet mapped from a reference quadrature rule.
- class FIAT.quadrature.GaussJacobiQuadratureLineRule(ref_el, m, a=0, b=0)[source]¶
Bases:
QuadratureRule
Gauss-Jacobi quadature rule determined by Jacobi weights a and b using m roots of m:th order Jacobi polynomial.
- class FIAT.quadrature.GaussLegendreQuadratureLineRule(ref_el, m)[source]¶
Bases:
GaussJacobiQuadratureLineRule
Gauss–Legendre quadrature rule on the interval.
The quadrature rule uses m points for a degree of precision of 2m-1.
- class FIAT.quadrature.GaussLobattoLegendreQuadratureLineRule(ref_el, m)[source]¶
Bases:
QuadratureRule
Gauss-Lobatto-Legendre quadrature rule on the interval.
The quadrature rule uses m points for a degree of precision of 2m-3.
- class FIAT.quadrature.QuadratureRule(ref_el, pts, wts)[source]¶
Bases:
object
General class that models integration over a reference element as the weighted sum of a function evaluated at a set of points.
- class FIAT.quadrature.RadauQuadratureLineRule(ref_el, m, right=True)[source]¶
Bases:
QuadratureRule
Gauss–Radau quadrature rule on the interval.
The quadrature rule uses m points for a degree of precision of 2m-1.
- FIAT.quadrature.make_quadrature(ref_el, m)[source]¶
Returns the collapsed quadrature rule using m points per direction on the given reference element. In the tensor product case, m is a tuple.
- FIAT.quadrature.make_tensor_product_quadrature(*quad_rules)[source]¶
Returns the quadrature rule for a TensorProduct cell, by combining the quadrature rules of the components.
FIAT.quadrature_element module¶
- class FIAT.quadrature_element.QuadratureElement(ref_el, points, weights=None)[source]¶
Bases:
FiniteElement
A set of quadrature points pretending to be a finite element.
- static is_nodal()[source]¶
True if primal and dual bases are orthogonal. If false, dual basis is not implemented or is undefined.
Subclasses may not necessarily be nodal, unless it is a CiarletElement.
FIAT.quadrature_schemes module¶
Quadrature schemes on cells
This module generates quadrature schemes on reference cells that integrate exactly a polynomial of a given degree using a specified scheme.
Scheme options are:
scheme=”default”
scheme=”canonical” (collapsed Gauss scheme)
Background on the schemes:
- Keast rules for tetrahedra:
Keast, P. Moderate-degree tetrahedral quadrature formulas, Computer Methods in Applied Mechanics and Engineering 55(3):339-348, 1986. http://dx.doi.org/10.1016/0045-7825(86)90059-9
- Xiao-Gimbutas rules for simplices:
Xiao, H., and Gimbutas, Z. A numerical algorithm for the construction of efficient quadrature rules in two and higher dimensions, Computers & mathematics with applications 59(2): 663-676, 2010.
- FIAT.quadrature_schemes.create_quadrature(ref_el, degree, scheme='default', entity=None)[source]¶
Generate quadrature rule for given reference element that will integrate a polynomial of order ‘degree’ exactly.
For low-degree polynomials on triangles (<=50) and tetrahedra (<=15), this uses hard-coded rules, otherwise it falls back to a collapsed Gauss scheme on simplices. On tensor-product cells, it is a tensor-product quadrature rule of the subcells.
- Parameters:
ref_el – The FIAT cell to create the quadrature for.
degree – The degree of polynomial that the rule should integrate exactly.
scheme – The quadrature scheme, can be choosen from [“default”, “canonical”, “KMV”] “default” -> hard-coded scheme for low degree and collapsed Gauss scheme for high degree, “canonical” -> collapsed Gauss scheme, “KMV” -> spectral lumped scheme for low degree (<=5 on triangles, <=3 on tetrahedra).
entity – A tuple of entity dimension and entity id specifying the integration domain. If not provided, the domain is the entire cell.
- FIAT.quadrature_schemes.xg_scheme(ref_el, degree)[source]¶
A (nearly) Gaussian simplicial quadrature with very few quadrature nodes, available for low-to-moderate orders.
Raises ValueError if no quadrature rule for the requested parameters is available.
See
H. Xiao and Z. Gimbutas, “A numerical algorithm for the construction of efficient quadrature rules in two and higher dimensions,” Computers & Mathematics with Applications, vol. 59, no. 2, pp. 663-676, 2010. http://dx.doi.org/10.1016/j.camwa.2009.10.027
FIAT.raviart_thomas module¶
- class FIAT.raviart_thomas.RTDualSet(ref_el, degree, variant, interpolant_deg)[source]¶
Bases:
DualSet
Dual basis for Raviart-Thomas elements consisting of point evaluation of normals on facets of codimension 1 and internal moments against polynomials
- FIAT.raviart_thomas.RTSpace(ref_el, degree)[source]¶
Constructs a basis for the Raviart-Thomas space (P_{degree-1})^d + P_{degree-1} x
- class FIAT.raviart_thomas.RaviartThomas(ref_el, degree, variant=None)[source]¶
Bases:
CiarletElement
The Raviart Thomas element
- Parameters:
ref_el – The reference element.
degree – The degree.
variant – optional variant specifying the types of nodes.
variant can be chosen from [“point”, “integral”, “integral(q)”] “point” -> dofs are evaluated by point evaluation. Note that this variant has suboptimal convergence order in the H(div)-norm “integral” -> dofs are evaluated by quadrature rules with the minimum degree required for unisolvence. “integral(q)” -> dofs are evaluated by quadrature rules with the minimum degree required for unisolvence plus q. You might want to choose a high quadrature degree to make sure that expressions will be interpolated exactly. This is important when you want to have (nearly) div-preserving interpolation.
FIAT.reference_element module¶
Abstract class and particular implementations of finite element reference simplex geometry/topology.
Provides an abstract base class and particular implementations for the reference simplex geometry and topology. The rest of FIAT is abstracted over this module so that different reference element geometry (e.g. a vertex at (0,0) versus at (-1,-1)) and orderings of entities have a single point of entry.
Currently implemented are UFC and Default Line, Triangle and Tetrahedron.
- class FIAT.reference_element.Cell(shape, vertices, topology)[source]¶
Bases:
object
Abstract class for a reference cell. Provides accessors for geometry (vertex coordinates) as well as topology (orderings of vertices that make up edges, faces, etc.
- cell_orientation_reflection_map()[source]¶
Return the map indicating whether each possible cell orientation causes reflection (
1
) or not (0
).
- construct_subcomplex(dimension)[source]¶
Constructs the reference subcomplex of the parent cell subentity specified by subcomplex dimension.
- Parameters:
dimension – tuple for tensor product cells, int otherwise
- construct_subelement(dimension)[source]¶
Constructs the reference element of a cell subentity specified by subelement dimension.
- Parameters:
dimension – tuple for tensor product cells, int otherwise
- extract_extrinsic_orientation(o)[source]¶
Extract extrinsic orientation.
Parameters¶
- oOrientation
Total orientation.
Returns¶
- Orientation
Extrinsic orientation.
- extract_intrinsic_orientation(o, axis)[source]¶
Extract intrinsic orientation.
Parameters¶
- oOrientation
Total orientation.
- axisint
Reference cell axis for which intrinsic orientation is computed.
Returns¶
- Orientation
Intrinsic orientation.
- property extrinsic_orientation_permutation_map¶
A map from extrinsic orientations to corresponding axis permutation matrices.
Notes¶
result[eo] gives the physical axis-reference axis permutation matrix corresponding to eo (extrinsic orientation).
- get_connectivity()[source]¶
Returns a dictionary encoding the connectivity of the element.
The dictionary’s keys are the spatial dimensions pairs ((1, 0), (2, 0), (2, 1), …) and each value is a list with entities of second dimension ordered by local dim0-dim1 numbering.
- get_dimension()[source]¶
Returns the subelement dimension of the cell. For tensor product cells, this a tuple of dimensions for each cell in the product. For all other cells, this is the same as the spatial dimension.
- get_entity_transform(dim, entity_i)[source]¶
Returns a mapping of point coordinates from the entity_i-th subentity of dimension dim to the cell.
- Parameters:
dim – tuple for tensor product cells, int otherwise
entity_i – entity number (integer)
- get_topology()[source]¶
Returns a dictionary encoding the topology of the element.
The dictionary’s keys are the spatial dimensions (0, 1, …) and each value is a dictionary mapping.
- class FIAT.reference_element.DefaultLine[source]¶
Bases:
DefaultSimplex
This is the reference line with vertices (-1.0,) and (1.0,).
- class FIAT.reference_element.DefaultTetrahedron[source]¶
Bases:
DefaultSimplex
This is the reference tetrahedron with vertices (-1,-1,-1), (1,-1,-1),(-1,1,-1), and (-1,-1,1).
- class FIAT.reference_element.DefaultTriangle[source]¶
Bases:
DefaultSimplex
This is the reference triangle with vertices (-1.0,-1.0), (1.0,-1.0), and (-1.0,1.0).
- class FIAT.reference_element.Hypercube(dimension, product)[source]¶
Bases:
Cell
Abstract class for a reference hypercube
- cell_orientation_reflection_map()[source]¶
Return the map indicating whether each possible cell orientation causes reflection (
1
) or not (0
).
- compute_reference_normal(facet_dim, facet_i)[source]¶
Returns the unit normal in infinity norm to facet_i.
- construct_subelement(dimension)[source]¶
Constructs the reference element of a cell subentity specified by subelement dimension.
- Parameters:
dimension – subentity dimension (integer)
- contains_point(point, epsilon=0)[source]¶
Checks if reference cell contains given point (with numerical tolerance as given by the L1 distance (aka ‘manhattan’, ‘taxicab’ or rectilinear distance) to the cell.
Parameters¶
- pointnumpy.ndarray, list or symbolic expression
The coordinates of the point.
- epsilonfloat
The tolerance for the check.
Returns¶
bool : True if the point is inside the cell, False otherwise.
- distance_to_point_l1(point, rescale=False)[source]¶
Get the L1 distance (aka ‘manhattan’, ‘taxicab’ or rectilinear distance) to a point with 0.0 if the point is inside the cell.
For more information see the docstring for the UFCSimplex method.
- get_dimension()[source]¶
Returns the subelement dimension of the cell. Same as the spatial dimension.
- class FIAT.reference_element.IntrepidTetrahedron[source]¶
Bases:
Simplex
This is the reference tetrahedron with vertices (0,0,0), (1,0,0),(0,1,0), and (0,0,1) used in the Intrepid project.
- class FIAT.reference_element.IntrepidTriangle[source]¶
Bases:
Simplex
This is the Intrepid triangle with vertices (0,0),(1,0),(0,1)
- class FIAT.reference_element.Simplex(shape, vertices, topology)[source]¶
Bases:
SimplicialComplex
Abstract class for a reference simplex.
Orientation of a physical cell is computed systematically by comparing the canonical orderings of its facets and the facets in the FIAT reference cell.
As an example, we compute the orientation of a triangular cell:
|1 0 47 42 | | +–2—+ +–43–+
FIAT canonical Mapped example physical cell
Suppose that the facets of the physical cell are canonically ordered as:
C = [43, 42, 47]
FIAT facet to Physical facet map is given by:
M = [42, 47, 43]
Then the orientation of the cell is computed as:
C.index(M[0]) = 1; C.remove(M[0]) C.index(M[1]) = 1; C.remove(M[1]) C.index(M[2]) = 0; C.remove(M[2])
o = (1 * 2!) + (1 * 1!) + (0 * 0!) = 3
- class FIAT.reference_element.SimplicialComplex(shape, vertices, topology)[source]¶
Bases:
Cell
Abstract class for a simplicial complex.
This consists of list of vertex locations and a topology map defining facets.
- compute_barycentric_coordinates(points, entity=None, rescale=False)[source]¶
Returns the barycentric coordinates of a list of points on the complex.
- compute_bubble(points, entity=None)[source]¶
Returns the lowest-order bubble on an entity evaluated at the given points on the cell.
- compute_edge_tangent(edge_i)[source]¶
Computes the nonnormalized tangent to a 1-dimensional facet. returns a single vector.
- compute_face_edge_tangents(dim, entity_id)[source]¶
Computes all the edge tangents of any k-face with k>=1. The result is a array of binom(dim+1,2) vectors. This agrees with compute_edge_tangent when dim=1.
- compute_face_tangents(face_i)[source]¶
Computes the two tangents to a face. Only implemented for a tetrahedron.
- compute_normal(facet_i, cell=None)[source]¶
Returns the unit normal vector to facet i of codimension 1.
- compute_normalized_edge_tangent(edge_i)[source]¶
Computes the unit tangent vector to a 1-dimensional facet
- compute_normalized_tangents(dim, i)[source]¶
Computes tangents in any dimension based on differences between vertices and the first vertex of the i:th facet of dimension dim. Returns a (possibly empty) list. These tangents are normalized to have unit length.
- compute_reference_normal(facet_dim, facet_i)[source]¶
Returns the unit normal in infinity norm to facet_i.
- compute_scaled_normal(facet_i)[source]¶
Returns the unit normal to facet_i of scaled by the volume of that facet.
- compute_tangents(dim, i)[source]¶
Computes tangents in any dimension based on differences between vertices and the first vertex of the i:th facet of dimension dim. Returns a (possibly empty) list. These tangents are NOT normalized to have unit length.
- contains_point(point, epsilon=0.0, entity=None)[source]¶
Checks if reference cell contains given point (with numerical tolerance as given by the L1 distance (aka ‘manhatten’, ‘taxicab’ or rectilinear distance) to the cell.
Parameters¶
- pointnumpy.ndarray, list or symbolic expression
The coordinates of the point.
- epsilonfloat
The tolerance for the check.
- entitytuple or None
A tuple of entity dimension and entity id.
Returns¶
bool : True if the point is inside the cell, False otherwise.
- distance_to_point_l1(points, entity=None, rescale=False)[source]¶
Get the L1 distance (aka ‘manhatten’, ‘taxicab’ or rectilinear distance) from an entity to a point with 0.0 if the point is inside the entity.
Parameters¶
- pointsnumpy.ndarray or list
The coordinates of the points.
- entitytuple or None
A tuple of entity dimension and entity id.
- rescalebool
If true, the L1 distance is measured with respect to rescaled barycentric coordinates, such that the L1 and L2 distances agree for points opposite to a single facet.
Returns¶
- numpy.float64 or numpy.ndarray
The L1 distance, also known as taxicab, manhatten or rectilinear distance, of the cell to the point. If 0.0 the point is inside the cell.
Notes¶
This is done with the help of barycentric coordinates where the general algorithm is to compute the most negative (i.e. minimum) barycentric coordinate then return its negative. For implementation reasons we return the sum of all the negative barycentric coordinates. In each of the below examples the point coordinate is X with appropriate dimensions.
Consider, for example, a UFCInterval. We have two vertices which make the interval,
P0 = [0] and P1 = [1].
- Our point is
X = [x].
- Barycentric coordinates are defined as
X = alpha * P0 + beta * P1 where alpha + beta = 1.0.
- The solution is
alpha = 1 - X[0] = 1 - x and beta = X[0] = x.
If both alpha and beta are positive, the point is inside the reference interval.
—regionA—P0=0——P1=1—regionB—
If we are in regionA, alpha is negative and -alpha = X[0] - 1.0 is the (positive) distance from P0. If we are in regionB, beta is negative and -beta = -X[0] is the exact (positive) distance from P1. Since we are in 1D the L1 distance is the same as the L2 distance. If we are in the interval we can just return 0.0.
Things get more complicated when we consider higher dimensions. Consider a UFCTriangle. We have three vertices which make the reference triangle,
P0 = (0, 0), P1 = (1, 0) and P2 = (0, 1).
- Our point is
X = [x, y].
Below is a diagram of the cell (which may not render correctly in sphinx):
- Barycentric coordinates are defined as
X = alpha * P0 + beta * P1 + gamma * P2 where alpha + beta + gamma = 1.0.
- The solution is
alpha = 1 - X[0] - X[1] = 1 - x - y, beta = X[0] = x and gamma = X[1] = y.
If all three are positive, the point is inside the reference cell. If any are negative, we are outside it. The absolute sum of any negative barycentric coordinates usefully gives the L1 distance from the cell to the point. For example the point (1,1) has L1 distance 1 from the cell: on this case alpha = -1, beta = 1 and gamma = 1. -alpha = 1 is the L1 distance. For comparison the L2 distance (the length of the vector from the nearest point on the cell to the point) is sqrt(0.5^2 + 0.5^2) = 0.707. Similarly the point (-1.0, -1.0) has alpha = 3, beta = -1 and gamma = -1. The absolute sum of beta and gamma 2 which is again the L1 distance. The L2 distance in this case is sqrt(1^2 + 1^2) = 1.414.
- For a UFCTetrahedron we have four vertices
P0 = (0,0,0), P1 = (1,0,0), P2 = (0,1,0) and P3 = (0,0,1).
- Our point is
X = [x, y, z].
- The barycentric coordinates are defined as
X = alpha * P0 + beta * P1 + gamma * P2 + delta * P3 where alpha + beta + gamma + delta = 1.0.
- The solution is
alpha = 1 - X[0] - X[1] - X[2] = 1 - x - y - z, beta = X[0] = x, gamma = X[1] = y and delta = X[2] = z.
The rules are the same as for the tetrahedron but with one extra barycentric coordinate. Our approximate distance, the absolute sum of the negative barycentric coordinates, is at worse around 4 times the actual distance to the tetrahedron.
- extract_extrinsic_orientation(o)[source]¶
Extract extrinsic orientation.
Parameters¶
- oOrientation
Total orientation.
Returns¶
- Orientation
Extrinsic orientation.
- extract_intrinsic_orientation(o, axis)[source]¶
Extract intrinsic orientation.
Parameters¶
- oOrientation
Total orientation.
- axisint
Reference cell axis for which intrinsic orientation is computed.
Returns¶
- Orientation
Intrinsic orientation.
- property extrinsic_orientation_permutation_map¶
A map from extrinsic orientations to corresponding axis permutation matrices.
Notes¶
result[eo] gives the physical axis-reference axis permutation matrix corresponding to eo (extrinsic orientation).
- get_dimension()[source]¶
Returns the subelement dimension of the cell. Same as the spatial dimension.
- get_entity_transform(dim, entity)[source]¶
Returns a mapping of point coordinates from the entity-th subentity of dimension dim to the cell.
- Parameters:
dim – subentity dimension (integer)
entity – entity number (integer)
- make_points(dim, entity_id, order, variant=None, interior=1)[source]¶
Constructs a lattice of points on the entity_id:th facet of dimension dim. Order indicates how many points to include in each direction.
- class FIAT.reference_element.TensorProductCell(*cells)[source]¶
Bases:
Cell
A cell that is the product of FIAT cells.
- cell_orientation_reflection_map()[source]¶
Return the map indicating whether each possible cell orientation causes reflection (
1
) or not (0
).
- compare(op, other)[source]¶
Parent-based comparison between simplicial complexes. This is done dimension by dimension.
- compute_reference_normal(facet_dim, facet_i)[source]¶
Returns the unit normal in infinity norm to facet_i of subelement dimension facet_dim.
- construct_subcomplex(dimension)[source]¶
Constructs the reference subcomplex of the parent cell subentity specified by subcomplex dimension.
- Parameters:
dimension – dimension in each “direction” (tuple)
- construct_subelement(dimension)[source]¶
Constructs the reference element of a cell subentity specified by subelement dimension.
- Parameters:
dimension – dimension in each “direction” (tuple)
- contains_point(point, epsilon=0.0)[source]¶
Checks if reference cell contains given point (with numerical tolerance as given by the L1 distance (aka ‘manhattan’, ‘taxicab’ or rectilinear distance) to the cell.
Parameters¶
- pointnumpy.ndarray, list or symbolic expression
The coordinates of the point.
- epsilonfloat
The tolerance for the check.
Returns¶
bool : True if the point is inside the cell, False otherwise.
- distance_to_point_l1(point, rescale=False)[source]¶
Get the L1 distance (aka ‘manhatten’, ‘taxicab’ or rectilinear distance) to a point with 0.0 if the point is inside the cell.
For more information see the docstring for the UFCSimplex method.
- extract_extrinsic_orientation(o)[source]¶
Extract extrinsic orientation.
Parameters¶
- oOrientation
Total orientation.
Returns¶
- Orientation
Extrinsic orientation.
Notes¶
The difinition of orientations used here must be consistent with that used in make_entity_permutations_tensorproduct.
- extract_intrinsic_orientation(o, axis)[source]¶
Extract intrinsic orientation.
Parameters¶
- oOrientation
Total orientation.
//
and%
must be overloaded in type(o).- axisint
Reference cell axis for which intrinsic orientation is computed.
Returns¶
- Orientation
Intrinsic orientation.
Notes¶
Must be consistent with make_entity_permutations_tensorproduct.
- property extrinsic_orientation_permutation_map¶
A map from extrinsic orientations to corresponding axis permutation matrices.
Notes¶
result[eo] gives the physical axis-reference axis permutation matrix corresponding to eo (extrinsic orientation).
- get_dimension()[source]¶
Returns the subelement dimension of the cell, a tuple of dimensions for each cell in the product.
- get_entity_transform(dim, entity_i)[source]¶
Returns a mapping of point coordinates from the entity_i-th subentity of dimension dim to the cell.
- Parameters:
dim – subelement dimension (tuple)
entity_i – entity number (integer)
- class FIAT.reference_element.UFCHexahedron[source]¶
Bases:
UFCHypercube
This is the reference hexahedron with vertices (0.0, 0.0, 0.0), (0.0, 0.0, 1.0), (0.0, 1.0, 0.0), (0.0, 1.0, 1.0), (1.0, 0.0, 0.0), (1.0, 0.0, 1.0), (1.0, 1.0, 0.0) and (1.0, 1.0, 1.0).
- class FIAT.reference_element.UFCHypercube(dim)[source]¶
Bases:
Hypercube
Reference UFC Hypercube
UFCHypercube: [0, 1]^d with vertices in lexicographical order.
- class FIAT.reference_element.UFCInterval[source]¶
Bases:
UFCSimplex
This is the reference interval with vertices (0.0,) and (1.0,).
- class FIAT.reference_element.UFCQuadrilateral[source]¶
Bases:
UFCHypercube
This is the reference quadrilateral with vertices (0.0, 0.0), (0.0, 1.0), (1.0, 0.0) and (1.0, 1.0).
Orientation of a physical cell is computed systematically by comparing the canonical orderings of its facets and the facets in the FIAT reference cell.
As an example, we compute the orientation of a quadrilateral cell:
+—3—+ +–57—+ | | | | 0 1 43 55 | | | | +—2—+ +–42—+
FIAT canonical Mapped example physical cell
Suppose that the facets of the physical cell are canonically ordered as:
C = [55, 42, 43, 57]
FIAT index to Physical index map must be such that C[0] = 55 is mapped to a vertical facet; in this example it is:
M = [43, 55, 42, 57]
C and M are decomposed into “vertical” and “horizontal” parts, keeping the relative orders of numbers:
C -> C0 = [55, 43], C1 = [42, 57] M -> M0 = [43, 55], M1 = [42, 57]
Then the orientation of the cell is computed as the following:
C0.index(M0[0]) = 1; C0.remove(M0[0]) C0.index(M0[1]) = 0; C0.remove(M0[1]) C1.index(M1[0]) = 0; C1.remove(M1[0]) C1.index(M1[1]) = 0; C1.remove(M1[1])
o = 2 * 1 + 0 = 2
- class FIAT.reference_element.UFCTetrahedron[source]¶
Bases:
UFCSimplex
This is the reference tetrahedron with vertices (0,0,0), (1,0,0),(0,1,0), and (0,0,1).
- class FIAT.reference_element.UFCTriangle[source]¶
Bases:
UFCSimplex
This is the reference triangle with vertices (0.0,0.0), (1.0,0.0), and (0.0,1.0).
- FIAT.reference_element.compute_unflattening_map(topology_dict)[source]¶
This function returns unflattening map for the given tensor product topology dict.
- FIAT.reference_element.default_simplex(spatial_dim)[source]¶
Factory function that maps spatial dimension to an instance of the default reference simplex of that dimension.
- FIAT.reference_element.flatten_entities(topology_dict)[source]¶
This function flattens topology dict of TensorProductCell and entity_dofs dict of TensorProductElement
- FIAT.reference_element.flatten_permutations(perm_dict)[source]¶
This function flattens permutation dict of TensorProductElement
- FIAT.reference_element.flatten_reference_cube(ref_el)[source]¶
This function flattens a Tensor Product hypercube to the corresponding UFC hypercube
- FIAT.reference_element.lattice_iter(start, finish, depth)[source]¶
Generator iterating over the depth-dimensional lattice of integers between start and (finish-1). This works on simplices in 0d, 1d, 2d, 3d, and beyond
- FIAT.reference_element.linalg_subspace_intersection(A, B)[source]¶
Computes the intersection of the subspaces spanned by the columns of 2-dimensional arrays A,B using the algorithm found in Golub and van Loan (3rd ed) p. 604. A should be in R^{m,p} and B should be in R^{m,q}. Returns an orthonormal basis for the intersection of the spaces, stored in the columns of the result.
- FIAT.reference_element.make_affine_mapping(xs, ys)[source]¶
Constructs (A,b) such that x –> A * x + b is the affine mapping from the simplex defined by xs to the simplex defined by ys.
- FIAT.reference_element.make_lattice(verts, n, interior=0, variant=None)[source]¶
Constructs a lattice of points on the simplex defined by verts. For example, the 1:st order lattice will be just the vertices. The optional argument interior specifies how many points from the boundary to omit. For example, on a line with n = 2, and interior = 0, this function will return the vertices and midpoint, but with interior = 1, it will only return the midpoint.
- FIAT.reference_element.multiindex_equal(d, isum, imin=0)[source]¶
A generator for d-tuple multi-indices whose sum is isum and minimum is imin.
- FIAT.reference_element.tuple_sum(tree)[source]¶
This function calculates the sum of elements in a tuple, it is needed to handle nested tuples in TensorProductCell. Example: tuple_sum(((1, 0), 1)) returns 2 If input argument is not the tuple, returns input.
- FIAT.reference_element.ufc_hypercube(spatial_dim)[source]¶
Factory function that maps spatial dimension to an instance of the UFC reference hypercube of that dimension.
FIAT.regge module¶
Implementation of the generalized Regge finite elements.
- class FIAT.regge.Regge(ref_el, degree=0, variant=None)[source]¶
Bases:
CiarletElement
The generalized Regge elements for symmetric-matrix-valued functions. REG(k) is the space of symmetric-matrix-valued polynomials of degree k or less with tangential-tangential continuity.
FIAT.restricted module¶
- class FIAT.restricted.RestrictedDualSet(dual, indices)[source]¶
Bases:
DualSet
Restrict the given DualSet to the specified list of dofs.
- class FIAT.restricted.RestrictedElement(element, indices=None, restriction_domain=None, take_closure=True)[source]¶
Bases:
CiarletElement
Restrict the given element to the specified list of dofs.
FIAT.serendipity module¶
- class FIAT.serendipity.Serendipity(ref_el, degree)[source]¶
Bases:
FiniteElement
- entity_closure_dofs()[source]¶
Return the map of topological entities to degrees of freedom on the closure of those entities for the finite element.
- entity_dofs()[source]¶
Return the map of topological entities to degrees of freedom for the finite element.
- tabulate(order, points, entity=None)[source]¶
Return tabulated values of derivatives up to given order of basis functions at given points.
- Parameters:
order – The maximum order of derivative.
points – An iterable of points.
entity – Optional (dimension, entity number) pair indicating which topological entity of the reference element to tabulate on. If
None
, default cell-wise tabulation is performed.
- FIAT.serendipity.unisolvent_pts_hex(K, deg)[source]¶
Gives a set of unisolvent points for the hex serendipity space of order deg. The S element is not dual to these nodes, but a dual basis can be constructed from them.
FIAT.tensor_product module¶
- class FIAT.tensor_product.FlattenedDimensions(element)[source]¶
Bases:
FiniteElement
A wrapper class that flattens entity dimensions of a FIAT element defined on a TensorProductCell to one with quadrilateral/hexahedron entities. TensorProductCell has dimension defined as a tuple of factor element dimensions (i, j) in 2D and (i, j, k) in 3D. Flattened dimension is a sum of the tuple elements.
- get_nodal_basis()[source]¶
Return the nodal basis, encoded as a PolynomialSet object, for the finite element.
- is_nodal()[source]¶
True if primal and dual bases are orthogonal. If false, dual basis is not implemented or is undefined.
Subclasses may not necessarily be nodal, unless it is a CiarletElement.
- class FIAT.tensor_product.TensorProductElement(A, B)[source]¶
Bases:
FiniteElement
Class implementing a finite element that is the tensor product of two existing finite elements.
- get_nodal_basis()[source]¶
Return the nodal basis, encoded as a PolynomialSet object, for the finite element.
- is_nodal()[source]¶
True if primal and dual bases are orthogonal. If false, dual basis is not implemented or is undefined.
Subclasses may not necessarily be nodal, unless it is a CiarletElement.
FIAT.xg_quad_data module¶
Module contents¶
FInite element Automatic Tabulator – supports constructing and evaluating arbitrary order Lagrange and many other elements. Simplices in one, two, and three dimensions are supported.