proteus.FemTools module

Class hierarchies for constructing and working with finite element spaces

Inheritance diagram of proteus.FemTools

class proteus.FemTools.ReferenceElement(dim=0, nNodes=0, nElementBoundaries=0)[source]

Bases: object

A base class for simple domains on which local function spaces will be built

onElement(xi)[source]
class proteus.FemTools.ReferenceSimplex(nd=3)[source]

Bases: proteus.FemTools.ReferenceElement

The unit simplex in \(R^n, n<=3\)

onElement(xi)[source]
class proteus.FemTools.ReferenceCube(nd=3)[source]

Bases: proteus.FemTools.ReferenceElement

The unit cube in \(R^n, n<=3\)

onElement(xi)[source]
class proteus.FemTools.LocalFunctionSpace(dim=0, referenceElement=None)[source]

Bases: object

Base class for low-dimensional spaces of functions on a reference element.

For example, linear functions on a reference triangle

defineTraceFunctions()[source]
class proteus.FemTools.LinearOnSimplexWithNodalBasis(nd=3)[source]

Bases: proteus.FemTools.LocalFunctionSpace

First order polynomials on the unit nd-simplex with the nodal basis.

Nodal basis functions on the reference nd-simplex (nd <= 3) with coordinates xi[0],xi[1],and xi[2]. The basis functions are numbered according to the nodes.

class proteus.FemTools.LinearOnCubeWithNodalBasis(nd=3)[source]

Bases: proteus.FemTools.LocalFunctionSpace

First order polynomials on the unit nd-cube with the nodal basis.

Nodal basis functions on the reference nd-cube (nd <=3) with coordinates xi[0],xi[1],and xi[2]. The basis functions are numbered according to the nodes.

class proteus.FemTools.LagrangeOnCubeWithNodalBasis(nd=3, order=2)[source]

Bases: proteus.FemTools.LocalFunctionSpace

Lagrange polynomials on the unit nd-cube with the nodal basis.

Nodal basis functions on the reference nd-cube (nd <=3) with coordinates xi[0],xi[1],and xi[2]. The basis functions are numbered according to the nodes.

class proteus.FemTools.BernsteinOnCube(nd=3, order=2)[source]

Bases: proteus.FemTools.LocalFunctionSpace

Bernstein polynomials on the unit nd-cube

Bernstein basis functions on the reference nd-cube (nd <=3) with coordinates xi[0],xi[1],and xi[2]. The basis functions are numbered as follows: ############## # * 2D * # ############## ########### # Order=1 # ########### 3–2 | | 0–1 ########### # Order=2 # ########### 3-6-2 | | | 7-8-5 | | | 0-4-1 ########### # Order=3 # ########### 3–9–8–2 | | | | 10-14-15-7 | | | | 11-12-13-6 | | | | 0–4–5–1 ########### # Order=4 # ########### 3–12-11-10-2 | | | | | 13-22-23-24-9 | | | | | 14-19-20-21-8 | | | | | 15-16-17-18-7 | | | | | 0–4–5–6–1

NOTE (mql): *** The polynomials are defined for arbitrary order in 1D.

This definition is carried to multi-D via tensor products.

* The numbering is defined for arbitrary order (via funMap) in 2D. * TODO: Generalize numbering for 1D and 3D. * TODO: The quad rule Lobatto edge alt must be generalized. * TODO: define hessians and higher order deriavtives

factorial()[source]

Find x!.

Raise a ValueError if x is negative or non-integral.

nChooseK(n, k)[source]
class proteus.FemTools.QuadraticOnSimplexWithNodalBasis(nd=3)[source]

Bases: proteus.FemTools.LocalFunctionSpace

Quadratic polynomials on the unit nd-simplex with the nodal basis.

Nodal basis functions on the reference nd-simplex (nd <=3) with coordinates xi[0],xi[1],and xi[2]. The basis functions are numbered according to

\[\]

psi &= lambda_i(2lambda_i-1) 0<= i<= d psi &= 4lambda_jlambda_k 0<= j < k <= d

where \(\lambda_i\) is the barycentric coordinate associated with node i (i.e., it’s 1 at node i and zero elsewhere)

Gradients of shape functions are

\[\]

abla psi_i &= (4lambda_i-1) ablalambda_i 0<= i <= d

abla psi_i &= 4lambda_k ablalambda_j + 4lambda_j ablalambda_k mbox{for} 0 <= j < k <= d

In 2d we have

\[\]

psi_i &= lambda_i(2lambda_i-1) 0<= i<= 2 psi_3 &= 4lambda_0lambda_1 psi_4 &= 4lambda_1lambda_2 psi_5 &= 4lambda_0lambda_2

2d numberings for \(\psi\)

2 | | | 5 4 | | 0—3–1

Note that mesh numbers edges according to the node they are across from, so that local dof 3 corresponds to edge 2, local dof 4 corresponds to edge 0, local dof 5 corresponds to edge 1,

3d should be

\[\]

psi_i &= lambda_i(2lambda_i-1) 0<= i<= 3 psi_4 &= 4lambda_0lambda_1 psi_5 &= 4lambda_1lambda_2 psi_6 &= 4lambda_2lambda_3 psi_7 &= 4lambda_0lambda_2 psi_8 &= 4lambda_1lambda_3 psi_9 &= 4lambda_0lambda_3

class proteus.FemTools.BernsteinOnSimplex(nd=3)[source]

Bases: proteus.FemTools.LocalFunctionSpace

Quadratic Bernstein polynomials on the unit nd-simplex.

The basis functions are numbered according to

\[\]

psi &= lambda_i^2 0<= i<= d psi &= 2lambda_jlambda_k 0<= j < k <= d

where \(\lambda_i\) is the barycentric coordinate associated with node i (i.e., it’s 1 at node i and zero elsewhere). These can be generalized to:

psi &= p!/(i!j!k!)lambda_0^ilambda_1^jlambda_2^k, 0<=i,j,k<=order, i+j+k=order

Gradients of shape functions are

\[\]

abla psi_i &= 2lambda_i ablalambda_i 0<= i <= d

abla psi_i &= 2lambda_k ablalambda_j + 2lambda_j ablalambda_k mbox{for} 0 <= j < k <= d

which again can be generalized via the general formula above.

In 2d we have

\[\]

psi_i &= lambda_i^2 0<= i<= 2 psi_3 &= 2lambda_0lambda_1 psi_4 &= 2lambda_1lambda_2 psi_5 &= 2lambda_0lambda_2

2d numberings for \(\psi\)

2 | | | 5 4 | | 0—3–1

Note that mesh numbers edges according to the node they are across from, so that local dof 3 corresponds to edge 2, local dof 4 corresponds to edge 0, local dof 5 corresponds to edge 1,

3d should be

\[\]

psi_i &= lambda_i^2 0<= i<= 3 psi_4 &= 2lambda_0lambda_1 psi_5 &= 2lambda_1lambda_2 psi_6 &= 2lambda_2lambda_3 psi_7 &= 2lambda_0lambda_2 psi_8 &= 2lambda_1lambda_3 psi_9 &= 2lambda_0lambda_3

class proteus.FemTools.CrouzeixRaviartWithNodalBasis(nd=3)[source]

Bases: proteus.FemTools.LocalFunctionSpace

Crouzeix-Raviart element implemented as a Lagrange type element. That is, the degrees of freedom are the values at barycenters of faces. Defined for the reference nd-simplex (nd <=3)

class proteus.FemTools.p0(nd=3)[source]

Bases: proteus.FemTools.LocalFunctionSpace

class proteus.FemTools.Monomials(nd=3, k=0)[source]

Bases: proteus.FemTools.LocalFunctionSpace

class proteus.FemTools.P1BubblesWithNodalBasis(nd=3)[source]

Bases: proteus.FemTools.LocalFunctionSpace

First order polynomials on the unit nd-simplex with the nodal basis + bubble

\[\]

b = (n_d+1)^{n_d+1} Pi_{i=0}^{n_d} lambda_i

First nd+1 basis functions are nodal ones on the reference nd-simplex (nd <=3) with coordinates xi[0],xi[1],and xi[2]. The basis functions are numbered according to the nodes. The last shape function is the bubble, b

class proteus.FemTools.P1P0BubblesWithNodalBasis(nd=3)[source]

Bases: proteus.FemTools.LocalFunctionSpace

First order polynomials on the unit nd-simplex with the nodal basis + dg bubble

\[\]

b_e = 1 for x in Omega_e

First nd+1 basis functions are nodal ones on the reference nd-simplex (nd <=3) with coordinates xi[0],xi[1],and xi[2]. The basis functions are numbered according to the nodes. The last shape function is the bubble, b

class proteus.FemTools.InterpolationConditions(dim=0, referenceElement=None)[source]

Bases: object

Base class for generalized interpolation conditions for function spaces.

For example, a function’s values at a set of points that is large enough to uniquely specify an element of the local function space that will be paired with the interpolation conditions.

quadrature2DOF_element(k)[source]
definedOnlocalElementBoundary(k, ebN_local)[source]
quadrature2Node_element(k)[source]
projectFiniteElementFunctionFromInterpolationConditions_opt = None[source]
class proteus.FemTools.NodalInterpolationConditions(referenceElement)[source]

Bases: proteus.FemTools.InterpolationConditions

Obtains the DOF from the function values at the nodes

definedOnLocalElementBoundary(k, ebN_local)[source]
quadrature2Node_element(k)[source]
projectFiniteElementFunctionFromInterpolationConditions_opt(finiteElementFunction, interpolationValues)[source]

Allow the interpolation conditions to control projection of a (global) finite element function from an array of interpolation values in order to take advantage of specific structure, otherwise can just use functionals interface

class proteus.FemTools.CubeNodalInterpolationConditions(referenceElement)[source]

Bases: proteus.FemTools.NodalInterpolationConditions

quadrilateralLocalBoundaryLookup = {0: [3, 0], 1: [0, 1], 2: [1, 2], 3: [2, 3]}[source]
hexahedronLocalBoundaryLookup = {0: [0, 1, 4], 1: [0, 1, 2], 2: [0, 2, 3], 3: [0, 3, 4], 4: [1, 4, 5], 5: [1, 2, 5], 6: [2, 3, 5], 7: [3, 4, 5]}[source]
definedOnLocalElementBoundary(k, ebN_local)[source]
class proteus.FemTools.QuadraticLagrangeNodalInterpolationConditions(referenceElement)[source]

Bases: proteus.FemTools.InterpolationConditions

Obtains the DOF from the function values at vertices and midpoints of edges (whole element is considered an edge in 1d)

p2tetrahedronLocalBoundaryLookup = {0: [1, 2, 3], 1: [0, 2, 3], 2: [0, 1, 3], 3: [0, 1, 2], 4: [2, 3], 5: [0, 3], 6: [0, 1], 7: [1, 3], 8: [0, 2], 9: [1, 2]}[source]
fmod(y, /)[source]

Return fmod(x, y), according to platform C.

x % y may differ.

definedOnLocalElementBoundary(k, ebN_local)[source]
quadrature2Node_element(k)[source]
projectFiniteElementFunctionFromInterpolationConditions_opt(finiteElementFunction, interpolationValues)[source]

Allow the interpolation conditions to control projection of a (global) finite element function from an array of interpolation values in order to take advantage of specific structure, otherwise can just use functionals interface

class proteus.FemTools.QuadraticLagrangeCubeNodalInterpolationConditions(referenceElement)[source]

Bases: proteus.FemTools.InterpolationConditions

Obtains the DOF from the function values at vertices and midpoints of edges (whole element is considered an edge in 1d)

q2quadrilateralLocalBoundaryLookup = {0: [3, 0], 1: [0, 1], 2: [1, 2], 3: [2, 3], 4: [0], 5: [1], 6: [2], 7: [3], 8: []}[source]
q2hexahedronLocalBoundaryLookup = {0: [0, 1, 4], 1: [0, 1, 2], 2: [0, 2, 3], 3: [0, 3, 4], 4: [1, 4, 5], 5: [1, 2, 5], 6: [2, 3, 5], 7: [3, 4, 5], 8: [3, 4, 5], 9: [3, 4, 5], 10: [3, 4, 5], 11: [3, 4, 5], 12: [3, 4, 5], 13: [3, 4, 5], 14: [3, 4, 5], 15: [3, 4, 5], 16: [3, 4, 5], 17: [3, 4, 5], 18: [3, 4, 5], 19: [3, 4, 5], 20: [3, 4, 5], 21: [3, 4, 5], 22: [3, 4, 5], 23: [3, 4, 5], 24: [3, 4, 5], 25: [3, 4, 5], 26: []}[source]
fmod(y, /)[source]

Return fmod(x, y), according to platform C.

x % y may differ.

definedOnLocalElementBoundary(k, ebN_local)[source]
quadrature2Node_element(k)[source]
projectFiniteElementFunctionFromInterpolationConditions_opt(finiteElementFunction, interpolationValues)[source]

Allow the interpolation conditions to control projection of a (global) finite element function from an array of interpolation values in order to take advantage of specific structure, otherwise can just use functionals interface

class proteus.FemTools.FaceBarycenterInterpolationConditions(referenceElement)[source]

Bases: proteus.FemTools.InterpolationConditions

Obtains the DOF from the function values at the barycenter of faces

definedOnLocalElementBoundary(k, ebN_local)[source]
projectFiniteElementFunctionFromInterpolationConditions_opt(finiteElementFunction, interpolationValues)[source]

Allow the interpolation conditions to control projection of a (global) finite element function from an array of interpolation values in order to take advantage of specific structure, otherwise can just use functionals interface

class proteus.FemTools.p0InterpolationConditions(referenceElement)[source]

Bases: proteus.FemTools.InterpolationConditions

Quadrature = <module 'proteus.Quadrature' from '/home/travis/build/erdc/proteus/proteus/Quadrature.py'>[source]
quadrature2DOF_element(k)[source]
definedOnLocalElementBoundary(k, ebN_local)[source]
class proteus.FemTools.MonomialInterpolationConditions(referenceElement, monomialSpace)[source]

Bases: proteus.FemTools.InterpolationConditions

Quadrature = <module 'proteus.Quadrature' from '/home/travis/build/erdc/proteus/proteus/Quadrature.py'>[source]
quadrature2DOF_element(k)[source]
definedOnLocalElementBoundary(k, ebN_local)[source]
class proteus.FemTools.P1BubbleInterpolationConditions(referenceElement)[source]

Bases: proteus.FemTools.InterpolationConditions

Interpolation conditions for space P1 enriched with bubbles

definedOnLocalElementBoundary(k, ebN_local)[source]
quadrature2Node_element(k)[source]
class proteus.FemTools.P1P0BubbleInterpolationConditions(referenceElement, monomialSpace)[source]

Bases: proteus.FemTools.InterpolationConditions

Quadrature = <module 'proteus.Quadrature' from '/home/travis/build/erdc/proteus/proteus/Quadrature.py'>[source]
quadrature2DOF_element(k)[source]
definedOnLocalElementBoundary(k, ebN_local)[source]
class proteus.FemTools.ReferenceFiniteElement(localFunctionSpace, interpolationConditions)[source]

Bases: object

The geometric element, local function space, and interpolation conditions.

I distinguish between the “elements” (patches) of the mesh and the “finite elements” which also include information about the function space. This is sort of consistent with the mathematical theory where a finite element is a triple: (T,PI,SIGMA)

class proteus.FemTools.DOFMap(nDOF=0)[source]

Bases: object

Base class for integer mappings between local degrees of freedom and global degrees of freedom.

updateAfterParallelPartitioning(mesh)[source]
class proteus.FemTools.NodalDOFMap(mesh)[source]

Bases: proteus.FemTools.DOFMap

The mapping that associates a local degree of freedom number with the global node numbers of the element.

updateAfterParallelPartitioning(mesh)[source]
class proteus.FemTools.DiscontinuousGalerkinDOFMap(mesh, localFunctionSpace)[source]

Bases: proteus.FemTools.DOFMap

A DOF map to use with discontinuous Galerkin spaces, which have unique degrees of freedom on each element.

updateAfterParallelPartitioning(globalMesh)[source]

Fix nDOF_all_processes, nDOF_subdomain, and max_dof_neighbors

class proteus.FemTools.ElementBoundaryDOFMap(mesh)[source]

Bases: proteus.FemTools.DOFMap

The mapping that associates a local degree of freedom number with the global edge numbers of the element.

updateAfterParallelPartitioning(mesh)[source]
class proteus.FemTools.QuadraticLagrangeCubeDOFMap(mesh, localFunctionSpace, nd)[source]

Bases: proteus.FemTools.DOFMap

DOF mapping for quadratic lagrange finite element functions on unit cubes

The mapping associates local degree of freedom with global vertex number for iloc 0<= iloc<= space dim global edge number for spacedim < iloc

total dimension is number of vertices + number of edges

updateAfterParallelPartitioning(globalMesh)[source]

Fix self.nDOF_all_processes,self.nDOF_subdomain, self.max_dof_neighbors

class proteus.FemTools.QuadraticLagrangeDOFMap(mesh, localFunctionSpace, nd)[source]

Bases: proteus.FemTools.DOFMap

DOF mapping for quadratic lagrange finite element functions on unit simplexes

The mapping associates local degree of freedom with global vertex number for iloc 0<= iloc<= space dim global edge number for spacedim < iloc

total dimension is number of vertices + number of edges

updateAfterParallelPartitioning(globalMesh)[source]

Fix self.nDOF_all_processes,self.nDOF_subdomain, self.max_dof_neighbors

class proteus.FemTools.p0DOFMap(mesh)[source]

Bases: proteus.FemTools.DOFMap

updateAfterParallelPartitioning(globalMesh)[source]

Fix self.nDOF_all_processes,self.nDOF_subdomain, self.max_dof_neighbors

class proteus.FemTools.P1BubbleDOFMap(mesh, localFunctionSpace, nd)[source]

Bases: proteus.FemTools.DOFMap

DOF mapping for Lagrange P1 + bubble finite element functions on unit simplexes

The mapping associates local degree of freedom with

global vertex number for iloc 0<= iloc<= space dim global element number for iloc = space dim + 1

total dimension is number of vertices + number of elements TODO: implement for parallel

class proteus.FemTools.ElementMaps(mesh)[source]

Bases: object

Base class for a set of real number vector fields that map a reference element into the physical domain.

getBasisValuesRef(xiArray)[source]

Evaluate the basis of the map at a set of points, xiArray, given on the reference element. Store in member array self.psi

getValues(xiArray, xArray)[source]

Evaluate the set of nElements maps at a set of points, xiArray, given on the reference element.

getBasisGradientValuesRef(xiArray)[source]

Evaluate the basis gradients of the map at a set of points, xiArray. Store in member self.grad_psi

getJacobianValues(xiArray, jacobianArray, jacobianDeterminantArray, inverseJacobianArray)[source]

Evaluate the jacobian, jacobian determinant, and inverse jacobian of the set of nElements maps at a set of points, xiArray.

getInverseValues(inverseJacobianArray, xArray, xiArray)[source]

Evaluate the set of nElements inverse maps at a set of points, xArray, given on the physical elements. Return the value of the inverse map.

getInverseValue(eN, x)[source]

Evaluate the element inverse map at a point, x, given in physical coorindates. Return the value of the inverse map, xi, in reference coordinates.

Note: The mapped point may not lie on the reference element

getBasisValuesTraceRef(xiArray)[source]

Evaluate the basis of the maps at a set of points, xiArray, on the boundaries of the reference element. Store in member self.psi_trace

getValuesTrace(xiArray, xArray)[source]

Evaluate the nElements x nElementBoundaries_element maps at a set of points, xiArray, on the element boundary reference element.

getJacobianValuesTrace(xiArray, jacobianInverseArray, metricTensorArray, metricTensorDeterminantSqrtArray, unitNormalArray)[source]

Evaluate the metric tensor, square root of the determinant of the metric tensor, and the unit normal at the nElements x nElementBoundaries_element element boundaries.

getBasisGradientValuesTraceRef(xiArray)[source]

Evaluate the basis gradients of the map at a set of points, xiArray, on the reference element boundaries.

getJacobianValuesTrace_movingDomain(xiArray, xt, jacobianInverseArray, metricTensorArray, metricTensorDeterminantSqrtArray, unitNormalArray)[source]

Evaluate the metric tensor, square root of the determinant of the metric tensor, and the unit normal at the nElements x nElementBoundaries_element element boundaries.

getValuesGlobalExteriorTrace(xiArray, xArray)[source]

Evaluate the nExteriorElementBoundaries_global maps at a set of points, xiArray, on the element boundary reference element.

getJacobianValuesGlobalExteriorTrace(xiArray, jacobianInverseArray, metricTensorArray, metricTensorDeterminantSqrtArray, unitNormalArray)[source]

Evaluate the metric tensor, square root of the determinant of the metric tensor, and the unit normal at the nExteriorElementBoundaries_global element boundaries.

getJacobianValuesGlobalExteriorTrace_movingDomain(xiArray, xt, jacobianInverseArray, metricTensorArray, metricTensorDeterminantSqrtArray, unitNormalArray)[source]

Evaluate the metric tensor, square root of the determinant of the metric tensor, and the unit normal at the nExteriorElementBoundaries_global element boundaries.

class proteus.FemTools.ParametricMaps(mesh, referenceElement, localFunctionSpace)[source]

Bases: proteus.FemTools.ElementMaps

A class that calculates the element maps using the nodes of a mesh as the degrees of freedom and a local function space for which the nodes are the local DOF.

getBasisValuesIP(interpolationPoints)[source]
getBasisValuesRef(xiArray)[source]

Evaluate the basis of the map at a set of points, xiArray, given on the reference element. Store in member array self.psi

getValues(xiArray, xArray)[source]

Evaluate the set of nElements maps at a set of points, xiArray, given on the reference element.

getBasisGradientValuesRef(xiArray)[source]

Evaluate the basis gradients of the map at a set of points, xiArray. Store in member self.grad_psi

getJacobianValues(xiArray, jacobianArray, jacobianInverseArray, jacobianDeterminantArray)[source]

Evaluate the jacobian, jacobian determinant, and inverse jacobian of the set of nElements maps at a set of points, xiArray.

getBasisValuesTraceRef(xiArray)[source]

Evaluate the basis of the maps at a set of points, xiArray, on the boundaries of the reference element. Store in member self.psi_trace

getValuesTrace(xiArray, xArray)[source]

Evaluate the nElements x nElementBoundaries_element maps at a set of points, xiArray, on the element boundary reference element.

getBasisGradientValuesTraceRef(xiArray)[source]

Evaluate the basis gradients of the map at a set of points, xiArray, on the reference element boundaries.

getJacobianValuesTrace(xiArray, jacobianInverseArray, metricTensorArray, metricTensorDeterminantSqrtArray, unitNormalArray)[source]

Evaluate the metric tensor, square root of the determinant of the metric tensor, and the unit normal at the nElements x nElementBoundaries_element element boundaries.

getJacobianValuesTrace_movingDomain(xiArray, xt, jacobianInverseArray, metricTensorArray, metricTensorDeterminantSqrtArray, unitNormalArray)[source]

Evaluate the metric tensor, square root of the determinant of the metric tensor, and the unit normal at the nElements x nElementBoundaries_element element boundaries.

getPermutations(xiArray)[source]
getValuesGlobalExteriorTrace(xiArray, xArray)[source]

treat this like element boundary version but only store values on exterior elements

getJacobianValuesGlobalExteriorTrace(xiArray, jacobianInverseArray, metricTensorArray, metricTensorDeterminantSqrtArray, unitNormalArray)[source]

treat this like element boundary version but only store values on exterior elements

getJacobianValuesGlobalExteriorTrace_movingDomain(xiArray, xt, jacobianInverseArray, metricTensorArray, metricTensorDeterminantSqrtArray, unitNormalArray)[source]

treat this like element boundary version but only store values on exterior elements

class proteus.FemTools.AffineMaps(mesh, referenceElement, localFunctionSpace)[source]

Bases: proteus.FemTools.ParametricMaps

getInverseValues(inverseJacobian, xArray, xiArray)[source]

Evaluate the set of nElements inverse maps at a set of points, xArray, given on the physical elements. Return the value of the inverse map.

getInverseValue(eN, x)[source]

Evaluate the element inverse map at a point, x, given in physical coorindates. Return the value of the inverse map, xi, in reference coordinates.

Note: The mapped point may not lie on the reference element

getInverseValuesTrace(inverseJacobian, xArray, xiArray)[source]
getInverseValuesGlobalExteriorTrace(inverseJacobian, xArray, xiArray)[source]
class proteus.FemTools.ParametricFiniteElementSpace(referenceFiniteElement, elementMaps, dofMap)[source]

Bases: object

Base class for spaces of functions defined by a set of finite elements that are each related to the same reference finite element.

dim – global dimension of the space mesh – a Mesh object, which is a partition of the domain finiteElements – a dictionary of the FiniteElement objects indexed by element number dofMap – a DOF

getBasisValuesRef(xiArray)[source]
getBasisValues(xiArray, vArray)[source]
getBasisValuesAtArray(xiArrayArray, vArray)[source]
getBasisGradientValuesRef(xiArray)[source]
getBasisGradientValues(xiArray, inverseJacobianArray, grad_vArray)[source]

This function calculates the BasisGradientValues for calculations on the reference element.

Parameters
  • xiArray (input, list of tuple) – A list of quadrature points (x,y,z) in 2D case, z = 0.

  • inverseJacobianArray (input, numpy.array) – Values of the inverseJacobian matrix used in the affine transformation from the physical domain to the reference element.

  • grad_vArray (output, numpy.array) – Gradient values of basis functions on reference triangle, adjusted for transformation from physical domain.

getBasisHessianValuesRef(xiArray)[source]
getBasisHessianValues(xiArray, inverseJacobianArray, Hessian_vArray)[source]
getBasisGradientValuesAtArray(xiArrayArray, inverseJacobianArray, grad_vArray)[source]
getBasisValuesTrace(permutations, xiArray, vArray)[source]

This function calculates the basis function values on the trace of the element boundaries. permutations (input) - xiArray (input) - a list of the element boundary quarature points mapped from the physical domain to the reference triangle. vArray (output) - the vector to store the basis function trace values on the reference triangle

getBasisValuesTraceAtArray(xiArrayArray, vArray)[source]
getBasisGradientValuesTrace(permutations, xiArray, inverseJacobianTraceArray, grad_vArray)[source]
getBasisGradientValuesTraceAtArray(xiArrayArray, inverseJacobianTraceArray, grad_vArray)[source]
getBasisValuesTraceRef(xiArray)[source]
getBasisValuesGlobalExteriorTrace(xiArray, vArray)[source]
getBasisValuesGlobalExteriorTraceAtArray(xiArrayArray, vArray)[source]
getBasisGradientValuesTraceRef(xiArray)[source]
getBasisGradientValuesGlobalExteriorTrace(xiArray, inverseJacobianTraceArray, grad_vArray)[source]
getBasisGradientValuesGlobalExteriorTraceAtArray(xiArrayArray, inverseJacobianTraceArray, grad_vArray)[source]
updateInterpolationPoints()[source]
endTimeSeriesEnsight(timeValues, filename, description, ts=1)[source]
getValuesAtMeshNodes(dof, nodalValues, isVector, dim_dof)[source]

we would like to be able to get references to existing values for values at nodes for some calculations like vtkVisualization

class proteus.FemTools.C0_AffineLinearOnSimplexWithNodalBasis(mesh, nd=3)[source]

Bases: proteus.FemTools.ParametricFiniteElementSpace

The standard linear CG space.

Globally C0 Each geometric element is the image of the reference simplex under a linear affine mapping. The nodal basis is used on the reference simplex.

writeMeshEnsight(filename, description=None)[source]
writeMeshXdmf(ar, name, t=0.0, init=False, meshChanged=False, arGrid=None, tCount=0)[source]
writeFunctionGnuplot(u, filename)[source]
writeFunctionXdmf(ar, u, tCount=0, init=True)[source]
readFunctionXdmf(ar, u, tCount=0)[source]
writeVectorFunctionXdmf(ar, uList, components, vectorName, tCount=0, init=True)[source]
writeFunctionEnsight(u, filename, append=False, firstVariable=True, case_filename=None)[source]
writeE2VectorFunctionEnsight(u, v, filename, nOutput, append=False, firstVariable=True, case_filename=None)[source]
writeE2VectorFunctionHeaderEnsight(u, v, filename, nOutput, append=False, firstVariable=True, case_filename=None)[source]
writeE3VectorFunctionEnsight(u, v, w, filename, nOutput, append=False, firstVariable=True, case_filename=None)[source]
writeE3VectorFunctionHeaderEnsight(u, v, w, filename, nOutput, append=False, firstVariable=True, case_filename=None)[source]
writeFunctionHeaderEnsight(u, filename, append=False, firstVariable=True, case_filename=None)[source]
writeFunctionMatlab(u, output, append=True, storeMeshData=True, figureOffset=1)[source]

save a scalar finite element function to matlab format for viewing returns number of function representations written

getValuesAtMeshNodes(dof, nodalValues, isVector, dim_dof)[source]
class proteus.FemTools.C0_LinearOnCubeWithNodalBasis(mesh, nd=3)[source]

Bases: proteus.FemTools.C0_AffineLinearOnSimplexWithNodalBasis

The standard linear CG space.

Globally C0 Each geometric element is the image of the reference cube under a n-linear(non-affine) mapping. The nodal basis is used on the reference cube.

class proteus.FemTools.C0_AffineLinearOnCubeWithNodalBasis(mesh, nd=3)[source]

Bases: proteus.FemTools.ParametricFiniteElementSpace

The standard linear CG space.

Globally C0 Each geometric element is the image of the reference simplex under a linear affine mapping. The nodal basis is used on the reference simplex.

writeMeshXdmf(ar, name, t=0.0, init=False, meshChanged=False, arGrid=None, tCount=0)[source]
writeFunctionXdmf(ar, u, tCount=0, init=True)[source]
writeVectorFunctionXdmf(ar, uList, components, vectorName, tCount=0, init=True)[source]
proteus.FemTools.P1[source]

alias of proteus.FemTools.C0_AffineLinearOnSimplexWithNodalBasis

class proteus.FemTools.C0_LagrangeOnCubeWithNodalBasis(mesh, nd=3, order=2)[source]

Bases: proteus.FemTools.C0_AffineLinearOnSimplexWithNodalBasis

The standard linear CG space.

Globally C0 Each geometric element is the image of the reference cube under a n-linear(non-affine) mapping. The nodal basis is used on the reference cube.

class proteus.FemTools.C0_AffineLagrangeOnCubeWithNodalBasis(mesh, nd=3, order=2)[source]

Bases: proteus.FemTools.ParametricFiniteElementSpace

The standard linear CG space.

Globally C0 Each geometric element is the image of the reference simplex under a linear affine mapping. The nodal basis is used on the reference simplex.

writeMeshXdmf(ar, name, t=0.0, init=False, meshChanged=False, arGrid=None, tCount=0)[source]
writeFunctionXdmf(ar, u, tCount=0, init=True)[source]
writeVectorFunctionXdmf(ar, uList, components, vectorName, tCount=0, init=True)[source]
proteus.FemTools.LagrangeCubeFactory(OrderIn)[source]
proteus.FemTools.Q1[source]

alias of proteus.FemTools.C0_AffineLinearOnCubeWithNodalBasis

proteus.FemTools.Q2[source]

alias of proteus.FemTools.LagrangeCubeFactory.<locals>.LagrangeCubeOrderN

proteus.FemTools.C0_AffineQuadraticOnCubeWithNodalBasis[source]

alias of proteus.FemTools.LagrangeCubeFactory.<locals>.LagrangeCubeOrderN

class proteus.FemTools.C0_BernsteinOnCube(mesh, nd=3, order=2)[source]

Bases: proteus.FemTools.C0_AffineLinearOnSimplexWithNodalBasis

The standard linear CG space. Globally C0 Each geometric element is the image of the reference cube under a n-linear(non-affine) mapping. The Bernstein basis is used on the reference cube.

class proteus.FemTools.C0_AffineBernsteinOnCube(mesh, nd=3, order=2)[source]

Bases: proteus.FemTools.ParametricFiniteElementSpace

Bernstein CG space. Globally C0 Each geometric element is the image of the reference simplex under a linear affine mapping. The Bernstein basis is used on the reference simplex.

writeMeshXdmf(ar, name, t=0.0, init=False, meshChanged=False, arGrid=None, tCount=0)[source]
writeFunctionXdmf(ar, u, tCount=0, init=True)[source]
writeVectorFunctionXdmf(ar, uList, components, vectorName, tCount=0, init=True)[source]
class proteus.FemTools.DG_AffinePolynomialsOnSimplexWithMonomialBasis(mesh, nd=3, k=0)[source]

Bases: proteus.FemTools.ParametricFiniteElementSpace

writeMeshEnsight(filename, description=None)[source]
writeFunctionGnuplot(u, filename)[source]
writeFunctionEnsight(u, filename, append=False, firstVariable=True, case_filename=None)[source]
writeFunctionHeaderEnsight(u, filename, append=False, firstVariable=True, case_filename=None)[source]
writeFunctionMatlab(u, output, append=True, storeMeshData=True, figureOffset=1)[source]

save a scalar finite element function to matlab format for viewing returns number of function representations written

Warning, currently has to compute values at interpolation points! tries to use basis values at interpolation points if in u

writeMeshXdmf(ar, name, t=0.0, init=False, meshChanged=False, arGrid=None, tCount=0)[source]
writeFunctionXdmf(ar, u, tCount=0, init=True)[source]

not much choice except to get u values at interpolation points?

writeVectorFunctionXdmf(ar, uList, components, vectorName, tCount=0, init=True)[source]
class proteus.FemTools.DG_AffineP0_OnSimplexWithMonomialBasis(mesh, nd=3)[source]

Bases: proteus.FemTools.DG_AffinePolynomialsOnSimplexWithMonomialBasis

writeMeshEnsight(filename, description=None)[source]
writeFunctionHeaderEnsight(u, filename, append=False, firstVariable=True, case_filename=None)[source]
writeFunctionEnsight(u, filename, append=False, firstVariable=True, case_filename=None)[source]
writeMeshXdmf(ar, name, t=0.0, init=False, meshChanged=False, arGrid=None, tCount=0)[source]
writeFunctionXdmf(ar, u, tCount=0, init=True)[source]

not much choice except to get u values at interpolation points?

writeVectorFunctionXdmf(ar, uList, components, vectorName, tCount=0, init=True)[source]
getValuesAtMeshNodes(dof, nodalValues, isVector, dim_dof)[source]
class proteus.FemTools.DG_AffineP1_OnSimplexWithMonomialBasis(mesh, nd=3)[source]

Bases: proteus.FemTools.DG_AffinePolynomialsOnSimplexWithMonomialBasis

class proteus.FemTools.DG_AffineP2_OnSimplexWithMonomialBasis(mesh, nd=3)[source]

Bases: proteus.FemTools.DG_AffinePolynomialsOnSimplexWithMonomialBasis

class proteus.FemTools.DG_AffineP3_OnSimplexWithMonomialBasis(mesh, nd=3)[source]

Bases: proteus.FemTools.DG_AffinePolynomialsOnSimplexWithMonomialBasis

class proteus.FemTools.DG_AffineP4_OnSimplexWithMonomialBasis(mesh, nd=3)[source]

Bases: proteus.FemTools.DG_AffinePolynomialsOnSimplexWithMonomialBasis

class proteus.FemTools.DG_AffineP5_OnSimplexWithMonomialBasis(mesh, nd=3)[source]

Bases: proteus.FemTools.DG_AffinePolynomialsOnSimplexWithMonomialBasis

class proteus.FemTools.DG_AffineP6_OnSimplexWithMonomialBasis(mesh, nd=3)[source]

Bases: proteus.FemTools.DG_AffinePolynomialsOnSimplexWithMonomialBasis

class proteus.FemTools.DG_AffineLinearOnSimplexWithNodalBasis(mesh, nd=3)[source]

Bases: proteus.FemTools.ParametricFiniteElementSpace

A linear DG space with the nodal basis.

Globally piecewise continuous. Each geometric element is the image of the reference simplex under a piecewise linear, continuous, affine mapping. The nodal basis is used on the reference simplex. mwf/cek

writeFunctionGnuplot(u, filename)[source]
writeMeshEnsight(filename, description=None)[source]
writeFunctionHeaderEnsight(u, filename, append=False, firstVariable=True, case_filename=None)[source]
writeFunctionEnsight(u, filename, append=False, firstVariable=True, case_filename=None)[source]
writeFunctionMatlab(u, output, append=True, storeMeshData=True, figureOffset=1)[source]

save a scalar finite element function to matlab format for viewing returns number of function representations written

writeMeshXdmf(ar, name, t=0.0, init=False, meshChanged=False, arGrid=None, tCount=0)[source]
writeFunctionXdmf(ar, u, tCount=0, init=True)[source]
writeVectorFunctionXdmf(ar, uList, components, vectorName, tCount=0, init=True)[source]
getValuesAtMeshNodes(dof, nodalValues, isVector, dim_dof)[source]

Calculate function at mesh nodes from degrees of freedom defined elsewhere

class proteus.FemTools.C0_AffineQuadraticOnSimplexWithNodalBasis(mesh, nd=3)[source]

Bases: proteus.FemTools.ParametricFiniteElementSpace

A quadratic C0 space with the nodal basis.

Globally piecewise continuous. Each geometric element is the image of the reference simplex under a piecewise linear, continuous, affine mapping. The nodal basis is used on the reference simplex.

writeFunctionGnuplot(u, filename)[source]

for now, just print out nodal dofs for vertices

writeMeshXdmf(ar, name, t=0.0, init=False, meshChanged=False, arGrid=None, tCount=0)[source]
writeMeshEnsight(filename, description=None)[source]
writeFunctionHeaderEnsight(u, filename, append=False, firstVariable=True, case_filename=None)[source]
writeFunctionXdmf(ar, u, tCount=0, init=True)[source]
readFunctionXdmf(ar, u, tCount=0)[source]
writeVectorFunctionXdmf(ar, uList, components, vectorName, tCount=0, init=True)[source]
writeFunctionEnsight(u, filename, append=False, firstVariable=True, case_filename=None)[source]
writeE2VectorFunctionEnsight(u, v, filename, nOutput, append=False, firstVariable=True, case_filename=None)[source]

old ensight printing

writeE2VectorFunctionHeaderEnsight(u, v, filename, nOutput, append=False, firstVariable=True, case_filename=None)[source]
writeE3VectorFunctionEnsight(u, v, w, filename, nOutput, append=False, firstVariable=True, case_filename=None)[source]
writeE3VectorFunctionHeaderEnsight(u, v, w, filename, nOutput, append=False, firstVariable=True, case_filename=None)[source]
writeFunctionMatlab(u, output, append=True, storeMeshData=True, figureOffset=1)[source]

save a scalar finite element function to matlab format for viewing returns number of function representations written

proteus.FemTools.P2[source]

alias of proteus.FemTools.C0_AffineQuadraticOnSimplexWithNodalBasis

class proteus.FemTools.C0_AffineBernsteinOnSimplex(mesh, nd=3)[source]

Bases: proteus.FemTools.C0_AffineQuadraticOnSimplexWithNodalBasis

A quadratic C0 space with Bernstein basis.

Globally piecewise continuous. Each geometric element is the image of the reference simplex under a piecewise linear, continuous, affine mapping. The nodal basis is used on the reference simplex.

class proteus.FemTools.DG_AffineQuadraticOnSimplexWithNodalBasis(mesh, nd=3)[source]

Bases: proteus.FemTools.ParametricFiniteElementSpace

A quadratic DG space with the nodal basis.

Globally piecewise continuous. Each geometric element is the image of the reference simplex under a piecewise linear, continuous, affine mapping. The nodal basis is used on the reference simplex.

writeMeshEnsight(filename, description=None)[source]
writeFunctionGnuplot(u, filename)[source]
writeFunctionHeaderEnsight(u, filename, append=False, firstVariable=True, case_filename=None)[source]
writeFunctionEnsight(u, filename, append=False, firstVariable=True, case_filename=None)[source]

For now only works for triangles

writeFunctionMatlab(u, output, append=True, storeMeshData=True, figureOffset=1)[source]

save a scalar finite element function to matlab format for viewing returns number of function representations written

writeMeshXdmf(ar, name, t=0.0, init=False, meshChanged=False, arGrid=None, tCount=0)[source]
writeFunctionXdmf(ar, u, tCount=0, init=True)[source]
writeVectorFunctionXdmf(ar, uList, components, vectorName, tCount=0, init=True)[source]
getValuesAtMeshNodes(dof, nodalValues, isVector, dim_dof)[source]

Calculate solution at mesh nodes from degrees of freedom defind elsewhere.

class proteus.FemTools.NC_AffineLinearOnSimplexWithNodalBasis(mesh, nd=3)[source]

Bases: proteus.FemTools.ParametricFiniteElementSpace

Space of functions that are P^1 on elements and continuous at element boundary barycenters.

Each geometric element is the image of the reference simplex under a linear affine mapping. The nodal basis is used on the reference simplex.

writeFunctionGnuplot(u, filename)[source]
writeMeshEnsight(filename, description=None)[source]
writeFunctionHeaderEnsight(u, filename, append=False, firstVariable=True, case_filename=None)[source]
writeFunctionEnsight(u, filename, append=False, firstVariable=True, case_filename=None)[source]
writeMeshXdmf(ar, name, t=0.0, init=False, meshChanged=False, arGrid=None, tCount=0)[source]
writeFunctionXdmf(ar, u, tCount=0, init=True)[source]
writeVectorFunctionXdmf(ar, uList, components, vectorName, tCount=0, init=True)[source]
writeFunctionMatlab(u, output, append=True, storeMeshData=True, figureOffset=1)[source]

save a scalar finite element function to matlab format for viewing returns number of function representations written

getValuesAtMeshNodes(dof, nodalValues, isVector, dim_dof)[source]

Calculate function at mesh nodes from degrees of freedom defined elsewhere

class proteus.FemTools.DG_Constants(mesh, nd=3)[source]

Bases: proteus.FemTools.ParametricFiniteElementSpace

The constant DG space

writeMeshEnsight(filename, description=None)[source]
writeFunctionGnuplot(u, filename)[source]
writeFunctionEnsight(u, filename, append=False, firstVariable=True, case_filename=None)[source]
writeFunctionHeaderEnsight(u, filename, append=False, firstVariable=True, case_filename=None)[source]
writeFunctionMatlab(u, output, append=True, storeMeshData=True, figureOffset=1)[source]

save a scalar finite element function to matlab format for viewing returns number of function representations written

writeMeshXdmf(ar, name, t=0.0, init=False, meshChanged=False, arGrid=None, tCount=0)[source]
writeFunctionXdmf(ar, u, tCount=0, init=True)[source]
writeVectorFunctionXdmf(ar, uList, components, vectorName, tCount=0, init=True)[source]
getValuesAtMeshNodes(dof, nodalValues, isVector, dim_dof)[source]

Compute function at mesh nodes from degrees of freedom defined elsewhere

class proteus.FemTools.C0_AffineP1BubbleOnSimplexWithNodalBasis(mesh, nd=3)[source]

Bases: proteus.FemTools.ParametricFiniteElementSpace

C0P1 enriched with bubble functions

writeMeshEnsight(filename, description=None)[source]
writeFunctionGnuplot(u, filename)[source]
writeFunctionEnsight(u, filename, append=False, firstVariable=True, case_filename=None)[source]
writeFunctionHeaderEnsight(u, filename, append=False, firstVariable=True, case_filename=None)[source]
writeFunctionMatlab(u, output, append=True, storeMeshData=True, figureOffset=1)[source]

save a scalar finite element function to matlab format for viewing returns number of function representations written

Warning, currently has to compute values at interpolation points! tries to use basis values at interpolation points if in u

writeMeshXdmf(ar, name, t=0.0, init=False, meshChanged=False, arGrid=None, tCount=0)[source]
writeFunctionXdmf(ar, u, tCount=0, init=True)[source]

not much choice except to get u values at interpolation points?

writeVectorFunctionXdmf(ar, uList, components, vectorName, tCount=0, init=True)[source]

Write function to XDFM rep

getValuesAtMeshNodes(dof, nodalValues, isVector, dim_dof)[source]

Calculate function values at mesh nodes from degrees of freedom defined elsewhere

class proteus.FemTools.C0_AffineP1P0BubbleOnSimplexWithNodalBasis(mesh, nd=3)[source]

Bases: proteus.FemTools.ParametricFiniteElementSpace

TODO set output to take advantage of nodal information

writeMeshEnsight(filename, description=None)[source]
writeFunctionGnuplot(u, filename)[source]
writeFunctionEnsight(u, filename, append=False, firstVariable=True, case_filename=None)[source]
writeFunctionHeaderEnsight(u, filename, append=False, firstVariable=True, case_filename=None)[source]
writeFunctionMatlab(u, output, append=True, storeMeshData=True, figureOffset=1)[source]

save a scalar finite element function to matlab format for viewing returns number of function representations written

Warning, currently has to compute values at interpolation points! tries to use basis values at interpolation points if in u

writeMeshXdmf(ar, name, t=0.0, init=False, meshChanged=False, arGrid=None, tCount=0)[source]
writeFunctionXdmf(ar, u, tCount=0, init=True)[source]

not much choice except to get u values at interpolation points?

writeVectorFunctionXdmf(ar, uList, components, vectorName, tCount=0, init=True)[source]

Write function to XDMF representation

class proteus.FemTools.FiniteElementFunction(finiteElementSpace, dof=None, dof_last=None, dof_last_last=None, dim_dof=1, name='no_name', isVector=False)[source]

Bases: object

A member of a finite element space of scalar functions.

This class manages the calculation of function, gradient, tensor, and trace values associated with a finite element function. Typically, this function will be a scalar function (For example, in a 2D stokes problem with a pressure, and two velocity components (say u or v), this class will manage the calculation of one of the scalar functions for p, u or v.

Parameters

finiteElementSpace (proteus.FemTools.ParametricFiniteElementSpace) –

Notes

femSpace and dof are the two main attributes of this class. The dof attributes provide the scalar multiples for the underlying basis functions that comes from femSpace’s ParametricFiniteElementSpace object

name[source]

Name assigned to the FiniteElementFunction instance

isVector[source]

boolean variable identifying whether the finiteElementFunction is vector quantitied.

femSpace[source]

proteus.FemTools.ParametricFiniteElementSpace provides information on the underlying finite element space.

nDOF_global[source]

integer variable with the global number of degrees of freedom

dim_dof[source]

dimension of the degree of freedom

dof[source]

numpy array of degree-of-freedom values used to calculate the finite element function.

useC[source]

boolean variable indentifying whether calculations should use C-routines.

par_dof[source]

for parallel FiniteElementFunciton capability

copy()[source]
projectFromInterpolationConditions(interpolationValues)[source]
getValue(eN, xi)[source]

Calculate the function value at a point on an element.

Parameters
  • eN (int) – Global element number.

  • xi (point) – Evaluation coordinate.

Returns

value – The finite element functions value at xi.

Return type

float

getValues(v, u)[source]
getGradientValues(grad_v, grad_u)[source]
getHessianValues(Hessian_v, Hessian_u)[source]
getGradientTensorValues(grad_v_x_grad_w, grad_u_x_grad_w)[source]
getValuesTrace(v, u)[source]
getGradientValuesTrace(grad_v, grad_u)[source]
getValuesGlobalExteriorTrace(v, u)[source]
getGradientValuesGlobalExteriorTrace(grad_v, grad_u)[source]
writeFunctionGnuplot(filename, append=False)[source]
writeFunctionEnsight(filename, append=False, case_filename=None)[source]
writeFunctionMatlab(output, append=True, storeMeshData=True, figureOffset=1)[source]
calculateValuesAtMeshNodes()[source]
setupParallelCommunication()[source]

build data structure for parallel communication

class proteus.FemTools.DOFBoundaryConditions(femSpace, getPointwiseBoundaryConditions=None, weakDirichletConditions=False, getPeriodicBoundaryConditions=None, allowNodalMaterialBoundaryTypes=True)[source]

Bases: object

A class for generating the set of DOF that are replaced by Dirichlet conditions and the values that are to be assigned to those DOF. For now I will ignore the ability to specify different boundary conditions for interpolatory DOF that correspond to the same physical point (i.e. discontinuous approximations where left and right limits can be specified)

DOFBoundaryConditionDict – a dictionary of boundary condition functions accessed by global DOF number DOFBoundryConditionPointDict – a dictionary of physical points associated with the boundary condition freeDOFSet – the DOF not specified by boundary conditions global2freeGlobal – an integer mapping from global DOF numbers to the free DOF numbers.

The constructor requires a function that takes a point as input and, if the point is on a part of the boundary where values are specified, then it returns the function of t,x, and the unknown that computes the boundary condition.

TODO

For now allow for flag to specify type of dirichlet condition to allow say nonlinear function of solution to be specified initially, only weak bc’s will allow this functionality though

class proteus.FemTools.DOFBoundaryConditions_alt(femSpace, getPointwiseBoundaryConditions=None, weakDirichletConditions=False, getPeriodicBoundaryConditions=None, allowNodalMaterialBoundaryTypes=True)[source]

Bases: object

A class for generating the set of DOF that are replaced by Dirichlet conditions and the values that are to be assigned to those DOF. For now I will ignore the ability to specify different boundary conditions for interpolatory DOF that correspond to the same physical point (i.e. discontinuous approximations where left and right limits can be specified)

DOFBoundaryConditionDict – a dictionary of boundary condition functions accessed by global DOF number DOFBoundryConditionPointDict – a dictionary of physical points associated with the boundary condition freeDOFSet – the DOF not specified by boundary conditions global2freeGlobal – an integer mapping from global DOF numbers to the free DOF numbers.

The constructor requires a function that takes a point as input and, if the point is on a part of the boundary where values are specified, then it returns the function of t,x, and the unknown that computes the boundary condition.

TODO

For now allow for flag to specify type of dirichlet condition to allow say nonlinear function of solution to be specified initially, only weak bc’s will allow this functionality though

class proteus.FemTools.FluxBoundaryConditions(mesh, nElementBoundaryQuadraturePoints_elementBoundary, x, getAdvectiveFluxBoundaryConditions=None, getDiffusiveFluxBoundaryConditions={}, getStressFluxBoundaryConditions=None)[source]

Bases: object

Generating class for the fluxBoundaryConditions dictionaries.

This class manages the generation of the advective, stress and diffusive fluxBoundaryCondition dictionaries. These dictionaries store lists of functions that describe the flux conditions along the finite element boundaries.

Variables
  • advectiveFluxBoundaryConditionsDict (dict) – Stores advective flux boundary conditions functions.

  • stressFluxBoundaryConditionsDict (dict) – Stores stress flux boundary conditions functions.

  • diffusiveFluxBoundaryConditionsDictDict (dict) – Stores diffusive flux boundary conditions functions.

Notes

The boundary condition dictionaries all use the key convention (ebNE,k) where ebNE is the global edge number and k is the local quadrature point number.

class proteus.FemTools.FluxBoundaryConditionsGlobalElementBoundaries(mesh, nElementBoundaryQuadraturePoints_elementBoundary, x, getAdvectiveFluxBoundaryConditions=None, getDiffusiveFluxBoundaryConditions={})[source]

Bases: object

mwf original version that sets indeces based on all element boundaries A class for generating the list of element boundaries where values are specified.

class proteus.FemTools.StressBoundaryConditions(mesh, nElementBoundaryQuadraturePoints_elementBoundary, x, getStressTraceBoundaryConditions=None)[source]

Bases: object

A class for generating the list of element boundaries where values are specified.

class proteus.FemTools.MultilevelProjectionOperators(multilevelMesh, femSpaceDictList, offsetListList, strideListList, dofBoundaryConditionsDictList)[source]

Bases: object

A class that takes a hierarchical (multiLevel) mesh and generates the interpolation and restriction operators.

restrictList/prolongList – projection and restriction at only the free nodes restrict_bcList/prolong_bcList – includes dirichlet boundaries as well

By default this is set up for conforming spaces. Since the spaces are conforming the coarse basis functions are in the fine space so we need only find the coefficients of the fine space basis functions that yield the coarse space basis functions. This is the matrix of the trivial injection from coarse to fine and it is used as the projection operator. Restriction is taken as the matrix of the adjoint of the injection, which is simply the transpose of the projection matrix. These operators fall out if you try to solve for the error on the coarse grid: Starting with u_f we have a(u_f+e_f,w_f) = <f,w_f>, and we want to solve instead a(e_c,w_c) = <f - a(u_f,w_f),w_c> on the coarse grid Using the injection i we can express this in the fine space as a(i e_c, i w_c) = <f - a(u_f,w_f),i w_c> writing this in matrix form yields p^t A_f p E = p^t R_f — P1 nonconforming space —- Try to set up now for nonconforming P1 approximation following Chen_96b. Setup prolongation by evaluating coarse grid basis functions at fine grid interpolation condition points (face barycenters).

Then, just need to know if fine grid interpolationCondition point falls on interface of coarse grid elements or not. If so, use average value of coarse grid quantity on fine grid. Otherwise just evaluate it

Use this simple interpolation from coarse to fine as the projection operator. Restriction is taken as the matrix of the adjoint of the injection, which is simply the transpose of the projection matrix.

I don’t think these fall out as nicely since they’re nonconforming.