FunctionSpace

A function space is a mesh combined with a basis function. This together defines the shape of elements and their nodes. A function space has two template parameters, the mesh and the basis function. The following is an example of a 3D quadratic function space.

FunctionSpace::FunctionSpace<Mesh::StructuredDeformableOfDimension<3>, BasisFunction::LagrangeOfOrder<2>>

The function space also knows the location of the nodes in the physical space. Depending on the Mesh class, it stores all node positions (for Mesh::StructuredDeformableOfDimension<D> and Mesh::UnstructuredDeformableOfDimension<D>) or only the mesh width (for Mesh::RegularFixedOfDimension<D>).

A function space is partioned over multiple processes, on the local process, only the portion of the local subdomain is stored. The case where there is only one process is just a special case of this general concept. Information about the local portion is given by the meshPartition object.

A function space object appears usually as a shared pointer, since multiple object have to know the function space and with a shared pointer it can be provided to every class without copy.

// define the function space class
typedef FunctionSpace::FunctionSpace<Mesh::StructuredDeformableOfDimension<2>, BasisFunction::LagrangeOfOrder<1>> FunctionSpaceType;  // this is an example

// this is how a function space is usually stored, as shared_ptr
std::shared_ptr<FunctionSpaceType> functionSpace;

// retrieve the mesh partition via
// functionSpace->meshPartition()

To create a new function space, use the meshManager object. This is a singleton class, the mesh manager is available from the DihuContext context object. It stores all existent function spaces and takes care that the same function space won’t be created twice.

Create a function space from settings (see Mesh for details). If this function space was already created, it returns a pointer to the already created function space. If it was not yet created, it constructs a new object. This works for all mesh types.

// assuming `context_` is the object of type DihuContext and `specificSettings_` is the python settings object of type `PythonConfig`
// (can be retrived by `specificSettings_ = context_.getPythonConfig();`)
std::shared_ptr<FunctionSpaceType> functionSpace
  = context_.meshManager()->functionSpace<FunctionSpaceType>(specificSettings_);

Create a function space with a Mesh::RegularFixedOfDimension<D> mesh directly, nElements is the local number of elements, in the subdomain, physicalExtent is the physical “size” of the subdomain. You can give the function space a name in the first argument, instead of “name”. Later the function space can be retrieved from the mesh manager with the functionSpace<> method, as above.

std::array<element_no_t, D> nElements;
std::array<double, D> physicalExtent;

std::shared_ptr<FunctionSpaceType> functionSpace
  = context_.meshManager()->createFunctionSpace<FunctionSpaceType>("name", nElements, physicalExtent);

Create a function space with a Mesh::StructuredDeformableOfDimension<D> mesh directly. nodePositions is a vector of all local node positions without ghosts, in the local ordering. nElementsPerCoordinateDirection is the local number of elements in the subdomain and nRanksPerCoordinateDirection is the global number of ranks in each coordinate direction.

const int D = 3;   // or 2 or 1
std::vector<Vec3> nodePositions;
const std::array<element_no_t,D> nElementsPerCoordinateDirection;
const std::array<int,D> nRanksPerCoordinateDirection;

std::shared_ptr<FunctionSpaceType> functionSpace
  = context_.meshManager()->createFunctionSpace<FunctionSpaceType>("name", nodePositions, nElementsPerCoordinateDirection, nRanksPerCoordinateDirection);

MeshPartition

A meshPartition object has the following methods:

int nRanks(int coordinateDirection) const

Number of ranks in a coordinate direction.

Parameters:

int coordinateDirection

element_no_t nElementsLocal() const

Number of elements in the current partition.

global_no_t nElementsGlobal() const

Number of elements in total.

dof_no_t nDofsLocalWithGhosts() const

Number of dofs in the local partition.

dof_no_t nDofsLocalWithoutGhosts() const

Number of dofs in the local partition, without ghosts.

global_no_t nDofsGlobal() const

Number of dofs in total.

node_no_t nNodesLocalWithGhosts() const

Number of nodes in the local partition.

node_no_t nNodesLocalWithoutGhosts() const

Number of nodes in the local partition.

node_no_t nNodesLocalWithGhosts(int coordinateDirection, int partitionIndex = -1) const

Number of nodes in the local partition specified by partitionIndex or the current partition if partitionIndex == -1.

Parameters:
  • int coordinateDirection

  • int partitionIndex = -1

node_no_t nNodesLocalWithoutGhosts(int coordinateDirection, int partitionIndex = -1) const

Number of nodes in the partition specified by partitionIndex or the current partition if partitionIndex == -1.

Parameters:
  • int coordinateDirection

  • int partitionIndex = -1

element_no_t nElementsLocal(int coordinateDirection) const

Number of elments in the local partition.

Parameters:

int coordinateDirection

element_no_t nElementsGlobal(int coordinateDirection) const

Number of elments in total.

Parameters:

int coordinateDirection

int beginElementGlobal(int coordinateDirection) const

Global no of first local element.

Parameters:

int coordinateDirection

global_no_t nNodesGlobal() const

Number of nodes in total.

global_no_t beginNodeGlobalNatural(int coordinateDirection, int partitionIndex = -1) const

Global no of first local node in the partition specified by partitionIndex or the current partition if partitionIndex == -1.

Parameters:
  • int coordinateDirection

  • int partitionIndex = -1

global_no_t nNodesGlobal(int coordinateDirection) const

Number of nodes in total.

Parameters:

int coordinateDirection

bool hasFullNumberOfNodes(int coordinateDirection, int partitionIndex = -1) const

Get if there are nodes on both boundaries in the given coordinate direction this is the case if the partition touches the right/top/back boundary. Consider the partition specified by partitionIndex or the current partition if partitionIndex == -1.

Parameters:
  • int coordinateDirection

  • int partitionIndex = -1

const std::vector<element_no_t> &localSizesOnRanks(int coordinateDirection) const

Get a vector with the local sizes on every rank.

Parameters:

int coordinateDirection

ISLocalToGlobalMapping localToGlobalMappingDofs()

Get the local to global mapping for the current partition, for the dof numbering.

global_no_t getElementNoGlobalNatural(element_no_t elementNoLocal) const

Get the global natural element no for a local element no.

Parameters:

element_no_t elementNoLocal

global_no_t getNodeNoGlobalNatural(std::array<global_no_t, MeshType::dim()> coordinatesGlobal) const

Get the global natural node no for the global coordinates of this node, this can be combined with getCoordinatesGlobal.

Parameters:

std::array<global_no_t,MeshType::dim()>

global_no_t getNodeNoGlobalPetsc(node_no_t nodeNoLocal) const

Get the node no in global petsc ordering from a local node no.

Parameters:

node_no_t nodeNoLocal

void getDofNoGlobalPetsc(const std::vector<dof_no_t> &dofNosLocal, std::vector<PetscInt> &dofNosGlobalPetsc) const

Transfer the local nos in global dof nos, using the PETSc localToGlobal mapping for the dofs.

Parameters:
  • const std::vector<dof_no_t> &dofNosLocal

  • std::vector<PetscInt> &dofNosGlobalPetsc

global_no_t getDofNoGlobalPetsc(dof_no_t dofNoLocal) const

Get the global petsc dof no for the local no, using the PETSc localToGlobal mapping for the dofs.

Parameters:

dof_no_t dofNoLocal

std::array<global_no_t, MeshType::dim()> getCoordinatesGlobal(node_no_t nodeNoLocal) const

Get the global node coordinates (x,y,z) of the node given by its local node no. This also works for ghost nodes.

std::array<int, MeshType::dim()> getCoordinatesLocal(node_no_t nodeNoLocal) const

Get the local coordinates for a local node no, also for ghost nodes. With this method and functionSpace->getNodeNo(coordinatesLocal) it is possible to implement a global-to-local mapping.

std::array<int, MeshType::dim()> getCoordinatesLocal(std::array<global_no_t, MeshType::dim()> coordinatesGlobal, bool &isOnLocalDomain) const

From global natural coordinates compute the local coordinates, set isOnLocalDomain to true if the node with global coordinates is in the local domain.

Parameters:
  • )> getCoordinatesLocal(std::array<global_no_t

  • MeshType::dim()> coordinatesGlobal

  • bool &isOnLocalDomain

std::array<int, MeshType::dim()> getElementCoordinatesLocal(element_no_t elementNoLocal) const

Get the local coordinates for a local element no.

element_no_t getElementNoLocal(std::array<int, MeshType::dim()> elementCoordinates) const

Get the local element no. from coordinates.

Parameters:

std::array<int,MeshType::dim()>

node_no_t getNodeNoLocal(global_no_t nodeNoGlobalPetsc) const

Get the local node no for a global petsc node no, does not work for ghost nodes.

Parameters:

global_no_t nodeNoGlobalPetsc

dof_no_t getDofNoLocal(global_no_t dofNoGlobalPetsc) const

Get the local dof no for a global petsc dof no, does not work for ghost nodes.

Parameters:

global_no_t dofNoGlobalPetsc

template<typename T>
void extractLocalNodesWithoutGhosts(std::vector<T> &vector, int nComponents = 1) const

From a vector of values of global/natural node numbers remove all that are non-local, nComponents consecutive values for each dof are assumed.

Parameters:
  • std::vector<T> &vector

  • int nComponents=1

template<typename T>
void extractLocalDofsWithoutGhosts(std::vector<T> &values) const

From a vector of values of global/natural dofs remove all that are non-local.

Parameters:

std::vector<T> &values

void extractLocalDofsWithoutGhosts(std::vector<double> &values) const

From a vector of values of global/natural dofs remove all that are non-local.

Parameters:

std::vector<double> &values

int convertRankNoToPartitionIndex(int coordinateDirection, int rankNo)

Get the partition index in a given coordinate direction from the rankNo.

Parameters:
  • int coordinateDirection

  • int rankNo

void output(std::ostream &stream)

Output to stream for debugging.

Parameters:

std::ostream &stream

const std::vector<PetscInt> &dofNosLocal(bool onlyNodalValues = false) const

Get a vector of local dof nos, range [0,nDofsLocalWithoutGhosts] are the dofs without ghost dofs, the whole vector are the dofs with ghost dofs @param onlyNodalValues: if for Hermite only get every second dof such that derivatives are not returned.

Parameters:

bool onlyNodalValues=false

void getDofNosGlobalNatural(std::vector<global_no_t> &dofNosGlobalNatural) const

Get a vector of global natural dof nos of the locally stored non-ghost dofs, needed for setParameters callback function in cellml adapter.

Parameters:

std::vector<global_no_t> &dofNosGlobalNatural

const std::vector<PetscInt> &ghostDofNosGlobalPetsc() const

Get the global dof nos of the ghost dofs in the local partition.

void initializeDofNosLocalNaturalOrdering(std::shared_ptr<FunctionSpace::FunctionSpace<MeshType, BasisFunctionType>> functionSpace)

Initialize the vector dofNosLocalNaturalOrdering_, this needs the functionSpace and has to be called before dofNosLocalNaturalOrdering() can be used. If the vector is already initialized by a previous call to this method, it has no effect.

Parameters:

std::shared_ptr<FunctionSpace::FunctionSpace<MeshType,BasisFunctionType>> functionSpace

const std::vector<dof_no_t> &dofNosLocalNaturalOrdering() const

Get a vector of local dof nos in local natural ordering, initializeDofNosLocalNaturalOrdering has to be called beforehand.

bool isNonGhost(node_no_t nodeNoLocal, int &neighbourRankNo) const

Check if the given dof is owned by the own rank, then return true, if not, neighbourRankNo is set to the rank by which the dof is owned.

Parameters:
  • node_no_t nodeNoLocal

  • int &neighbourRankNo

void getBoundaryElements(Mesh::face_t face, int &neighbourRankNo, std::array<element_no_t, MeshType::dim()> &nBoundaryElements, std::vector<dof_no_t> &dofNos)

Get information about neighbouring rank and boundary elements for specified face, @param neighbourRankNo: the rank of the neighbouring process that shares the face, @param nElements: Size of one-layer mesh that contains boundary elements that touch the neighbouring process.

Parameters:
  • Mesh::face_t face

  • int &neighbourRankNo

  • std::array<element_no_t,MeshType::dim()> &nBoundaryElements

  • std::vector<dof_no_t> &dofNos

int neighbourRank(Mesh::face_t face)

Get the rank no of the neighbour in direction face, -1 if there is no such neighbour.

Parameters:

Mesh::face_t face

int ownRankPartitioningIndex(int coordinateDirection)

Get the partitioning index in the coordinate direction, i.e. the no. of this rank in this direction, the total number of ranks in each direction can be retrieved by nRanks.

Parameters:

int coordinateDirection

FunctionSpace

The function space object has the following methods:

static constexpr int nDofsPerElement()

Number of degrees of freedom of this basis.

static constexpr int nNodesPerElement()

Number of nodes per element.

static constexpr int nDofsPerNode()

Number of dofs per node.

static constexpr int averageNDofsPerElement()

If one assigns every dof to an element it is contained in, the number of degrees of freedom per element (not considering boundary elements).

static constexpr int averageNNodesPerElement()

If one assigns every node to an element it is contained in, the number of nodes per element (not considering boundary elements).

static double phi(int dofIndex, std::array<double, MeshType::dim()> xi)

Evaluate the basis function corresponding to element-local dof dofIndex at xi, xi lives in [0,1]^D.

Parameters:
  • int dofIndex

  • std::array<double,MeshType::dim()>

static double dphi_dxi(int dofIndex, int derivativeIdx, std::array<double, MeshType::dim()> xi)

Evaluate the derivative of Phi(xi) w.r.t xi_i, where i is given by derivativeIdx, i.e. Phi_{dofIndex,derivativeIdx}(xi).

Parameters:
  • int dofIndex

  • int derivativeIdx

  • std::array<double,MeshType::dim()>

static std::array<double, MeshType::dim()> gradPhi(int dofIndex, std::array<double, MeshType::dim()> xi)

Evaluate the first derivative of the basis function corresponding to element-local dof dofIndex at xi, interval for xi is [0,1]^D.

Parameters:
  • int dofIndex

  • std::array<double,MeshType::dim()>

void setMeshPartition(std::shared_ptr<Partition::MeshPartition<FunctionSpace<MeshType, BasisFunctionType>, MeshType>> meshPartition)

Set the partition, call this prior to initialize to not initialize the partition from settings but use the given meshPartition.

Parameters:

std::shared_ptr<Partition::MeshPartition<FunctionSpace<MeshType,BasisFunctionType>,MeshType>> meshPartition

std::shared_ptr<Partition::MeshPartition<FunctionSpace<MeshType, BasisFunctionType>, MeshType>> meshPartition() const

Get the partition.

std::shared_ptr<Partition::MeshPartitionBase> meshPartitionBase()

Get the partition as pointer of type meshPartitionBase, this is in the itnerface in mesh.

dof_no_t getDofNo(element_no_t elementNoLocal, int dofIndex) const

Return the local dof number of element-local dof dofIndex of element elementNoLocal.

Parameters:
  • element_no_t elementNoLocal

  • int dofIndex

node_no_t getNodeNo(element_no_t elementNoLocal, int nodeIndex) const

Return the local node number of element-local node nodeIndex of element with local no elementNoLocal.

Parameters:
  • element_no_t elementNoLocal

  • int nodeIndex

global_no_t getNodeNoGlobalNatural(global_no_t elementNoLocalGlobalNatural, int nodeIndex) const

Return the global/natural node number of element-local node nodeIndex of element with global no elementNoLocalGlobal.

Parameters:
  • global_no_t elementNoLocalGlobalNatural

  • int nodeIndex

void getNodeDofs(node_no_t nodeGlobalNo, std::vector<dof_no_t> &dofGlobalNos) const

Get all dofs of a specific node, as vector.

Parameters:
  • node_no_t nodeGlobalNo

  • std::vector<dof_no_t> &dofGlobalNos

void getNodeDofs(node_no_t nodeGlobalNo, std::array<dof_no_t, FunctionSpaceBaseDim<1, BasisFunctionType>::nDofsPerNode()> &dofGlobalNos) const

Get all dofs of a specific node, as array.

Parameters:
  • node_no_t nodeGlobalNo

  • std::array<dof_no_t,FunctionSpaceBaseDim<1,BasisFunctionType>::nDofsPerNode(

dof_no_t getNodeDofNo(node_no_t nodeGlobalNo, int dofIndex) const

Get the dof no of the specified dof at the node.

Parameters:
  • node_no_t nodeGlobalNo

  • int dofIndex

node_no_t getNeighbourNodeNoLocal(node_no_t nodeNoLocal, Mesh::face_t direction) const

Get neighbouring node to nodeNoLocal or -1 if there is no such node, nodeNoLocal has to be a non-ghost local node.

Parameters:
  • node_no_t nodeNoLocal

  • Mesh::face_t direction

node_no_t getNodeNo(std::array<int, MeshType::dim()> coordinateLocal) const

Get node local no from the local coordinate in natural local numbering.

Parameters:

std::array<int,MeshType::dim()>

std::shared_ptr<FieldVariableBaseFunctionSpaceType> fieldVariable(std::string name)

Return a field variable with given name, this is not implemented for structured meshes since there are no extra stored field variables, only for unstructured meshes is it implemented and then stores field variables that were present in parsed exfiles.

Parameters:

std::string name

dof_no_t getDofNoLocal(std::array<global_no_t, MeshType::dim()> coordinatesGlobal, int nodalDofIndex, bool &isOnLocalDomain)

Get the local dof no. for the global coordinates.

Parameters:
  • std::array<global_no_t,MeshType::dim()> coordinatesGlobal

  • int nodalDofIndex

  • bool &isOnLocalDomain

double meshWidth() const

Get mesh width (=distance between nodes) of the given coordinate direction.

node_no_t nNodesLocalWithGhosts() const

Return number of nodes including ghost nodes, i.e. these nodes are known locally but some of them are owned by other ranks.

node_no_t nNodesLocalWithGhosts(int dimension) const

Return number of nodes in specified coordinate direction.

Parameters:

int dimension

node_no_t nNodesLocalWithoutGhosts() const

Return number of nodes that are owned by this partition.

node_no_t nNodesLocalWithoutGhosts(int dimension) const

Return number of nodes in specified coordinate direction that are owned by this partition.

Parameters:

int dimension

dof_no_t nDofsLocalWithGhosts() const

Return number of dofs.

dof_no_t nDofsLocalWithoutGhosts() const

Return number of dofs.

global_no_t nNodesGlobal(int dimension) const

Return number of nodes in specified coordinate direction for the whole global domain.

Parameters:

int dimension

global_no_t nNodesGlobal() const

Return global number of nodes.

global_no_t nDofsGlobal() const

Return global number of dofs.

void getNodePositions(std::vector<double> &nodes) const

Fill a vector with the node position entries, nodes will contain consecutively the (x,y,z) values of just all nodes, i.e. for Hermite not the derivatives.

Parameters:

std::vector<double> &nodes

static void getFaceDofs(Mesh::face_t face, std::array<dof_no_t, FunctionSpaceBaseDim<1, BasisFunctionType>::nDofsPerNode()> &dofIndices)

Get all dof indices of a face, note: dimension in FunctionSpaceBaseDim is current-1 (=0), in this case the dofIndices array has exactly so many entries as there are dofs for a node.

Parameters:
  • Mesh::face_t face

  • std::array<dof_no_t,FunctionSpaceBaseDim<1,BasisFunctionType>::nDofsPerNode(

static int getNeighbourNodeIndex(int nodeIndex, Mesh::face_t face)

Get the neighbouring elemental node index in given direction inside one element or -1 if there is no such node in the element in that direction.

Parameters:
  • int nodeIndex

  • Mesh::face_t face

std::array<dof_no_t, FunctionSpaceFunction<MeshType, BasisFunctionType>::nNodesPerElement()> getElementNodeNos(element_no_t elementNo) const

Return an array of all node nos. of the element.

bool findPosition(Vec3 point, element_no_t &elementNo, int &ghostMeshNo, std::array<double, D> &xi, bool startSearchInCurrentElement, double xiTolerance = 1e-4)

Get the element no and the xi value of the point, return true if the point is inside the mesh or false otherwise. Start search at given elementNo ghostMeshNo: -1 means main mesh, 0-5 means ghost Mesh with respecitve Mesh::face_t.

Parameters:
  • Vec3 point

  • element_no_t &elementNo

  • int &ghostMeshNo

  • std::array<double,D> &xi

  • bool startSearchInCurrentElement

  • double xiTolerance = 1e-4

bool pointIsInElement(Vec3 point, element_no_t elementNo, std::array<double, D> &xi, double xiTolerance)

Check if the point lies inside the element, if yes, return true and set xi to the value of the point, defined in 11_function_space_xi.h.

Parameters:
  • Vec3 point

  • element_no_t elementNo

  • std::array<double,D> &xi

  • double xiTolerance

void setGhostMesh(Mesh::face_t face, const std::shared_ptr<FunctionSpace<MeshType, BasisFunctionType>> ghostMesh)

Store a ghost mesh which is a neighouring mesh with only one layer of elements, this will be used by pointIsInElement and findPosition.

Parameters:
  • Mesh::face_t face

  • const std::shared_ptr<FunctionSpace<MeshType,BasisFunctionType>> ghostMesh

bool findPosition(Vec3 point, element_no_t &elementNo, int &ghostMeshNo, std::array<double, MeshType::dim()> &xi, bool startSearchInCurrentElement, double xiTolerance = 1e-4)

Get the element no and the xi value of the point, return true if the point is inside the mesh or false otherwise. Start search at given elementNo ghostMeshNo: -1 means main mesh, 0-5 means ghost Mesh with respecitve Mesh::face_t.

Parameters:
  • Vec3 point

  • element_no_t &elementNo

  • int &ghostMeshNo

  • std::array<double,MeshType::dim()> &xi

  • bool startSearchInCurrentElement

  • double xiTolerance = 1e-4

bool checkNeighbouringElements(const Vec3 &point, element_no_t &elementNo, int &ghostMeshNo, std::array<double, MeshType::dim()> &xi)

Check if the point is in a neighbouring element to elementNo on ghostMeshNo (-1=main mesh, 0-5=ghost mesh on respective face, 0=face0Minus, 1=face0Plus, etc.), return true if the element was found amoung the neighbours set elementNo, ghostMeshNo and xi appropriately.

Parameters:
  • const Vec3 &point

  • element_no_t &elementNo

  • int &ghostMeshNo

  • std::array<double,MeshType::dim()>

std::shared_ptr<FieldVariable::FieldVariableBaseFunctionSpace<FunctionSpace<MeshType, BasisFunctionType>>> createFieldVariable(std::string name, std::vector<std::string> componentNames)

Create a non-geometry field field variable with no values being set, with given component names.

Parameters:
  • std::string name

  • std::vector<std::string> componentNames

std::shared_ptr<FieldVariable::FieldVariableBaseFunctionSpace<FunctionSpace<MeshType, BasisFunctionType>>> createFieldVariable(std::string name, int nComponents = 1)

Create a non-geometry field field variable with no values being set, with given number of components, the component names will be the numbers.

Parameters:
  • std::string name

  • int nComponents=1

template<int nComponents>
std::shared_ptr<FieldVariable::FieldVariable<FunctionSpace<MeshType, BasisFunctionType>, nComponents>> createFieldVariable(std::string name)

Create a non-geometry field field variable with no values being set, with given number of components, the component names will be the numbers.

Parameters:

std::string name

template<int nComponents>
std::shared_ptr<FieldVariable::FieldVariable<FunctionSpace<MeshType, BasisFunctionType>, nComponents>> createFieldVariable(std::string name, std::vector<std::string> componentNames)

Create a non-geometry field field variable with no values being set, with given number of components and component names.

Parameters:
  • std::string name

  • std::vector<std::string> componentNames

int getNumberScaleFactors(element_no_t elementGlobalNo)

Get the number of scale factors that are stored for the global element no.

Parameters:

element_no_t elementGlobalNo

std::array<std::array<double, MeshType::dim()>, FunctionSpaceFunction<MeshType, BasisFunctionType>::nDofsPerElement()> getGradPhi(std::array<double, MeshType::dim()> xi) const

Return an array of the gradients of all nodal basis functions, evaluated at xi.

Parameters:
  • )>

  • FunctionSpaceFunction<MeshType,BasisFunctionType>::nDofsPerElement()>getGradPhi(std::array<double

  • MeshType::dim()>

template<int nComponents>
std::array<double, nComponents> interpolateValueInElement(std::array<std::array<double, nComponents>, FunctionSpaceFunction<MeshType, BasisFunctionType>::nDofsPerElement()> &elementalDofValues, std::array<double, MeshType::dim()> xi) const

Interpolate the nComponents values within an element at the given xi position using the basis functions.

Parameters:
  • std::array<std::array<double,nComponents>,FunctionSpaceFunction<MeshType,BasisFunctionType>::nDofsPerElement()> &elementalDofValues

  • std::array<double,MeshType::dim()>

double interpolateValueInElement(std::array<double, FunctionSpaceFunction<MeshType, BasisFunctionType>::nDofsPerElement()> &elementalDofValues, std::array<double, MeshType::dim()> xi) const

Interpolate the value within an element at the given xi position using the basis functions.

Parameters:
  • std::array<double,FunctionSpaceFunction<MeshType,BasisFunctionType>::nDofsPerElement()> &elementalDofValues

  • std::array<double,MeshType::dim()>

std::array<double, MeshType::dim()> interpolateGradientInElement(std::array<double, FunctionSpaceFunction<MeshType, BasisFunctionType>::nDofsPerElement()> &elementalDofValues, Tensor2<MeshType::dim()> inverseJacobianParameterSpace, std::array<double, MeshType::dim()> xi) const

Interpolate the gradient of a scalar field within an element at the given xi position using the basis functions the inverseJacobianParameterSpace can be computed by getInverseJacobian.

Parameters:
  • )> interpolateGradientInElement(std::array<double

  • FunctionSpaceFunction<MeshType,BasisFunctionType>::nDofsPerElement()> &elementalDofValues

  • Tensor2<MeshType::dim()> inverseJacobianParameterSpace

  • std::array<double,MeshType::dim()>

Vec3 getNormal(Mesh::face_t face, std::array<Vec3, FunctionSpaceFunction<MeshType, BasisFunctionType>::nDofsPerElement()> geometryValues, std::array<double, MeshType::dim()> xi)

Compute the normal in world space, normal to face at xi, use the given geometry values, that can by obtained by fieldVariable->getElementValues(elementNo, geometryValues) or mesh->getElementGeometry(elementNo, geometryValues).

Parameters:
  • Mesh::face_t face

  • std::array<Vec3,FunctionSpaceFunction<MeshType,BasisFunctionType>::nDofsPerElement()> geometryValues

  • std::array<double,MeshType::dim()>

Vec3 getNormal(Mesh::face_t face, element_no_t elementNoLocal, std::array<double, MeshType::dim()> xi)

Compute the normal in world space, normal to face at xi.

Parameters:
  • Mesh::face_t face

  • element_no_t elementNoLocal

  • std::array<double,MeshType::dim()>

Tensor2<MeshType::dim()> getInverseJacobian(std::array<Vec3, FunctionSpaceFunction<MeshType, BasisFunctionType>::nDofsPerElement()> &geometryValues, element_no_t elementNo, std::array<double, MeshType::dim()> xi)

Compute the inverseJacobian that is needed to transform a gradient vector from parameter space to world space, for an element at a xi position. This version of the method needs the values of the geometry field, if the jacobian is needed at multiple positions in the same element, these values can be retrieved once and used for all computations of the jacobians. There is also the convienience method which does not need the geometryValues but gets them itself. The following properties of the jacobian hold: jacobianParameterSpace[columnIdx][rowIdx] = dX_rowIdx/dxi_columnIdx inverseJacobianParameterSpace[columnIdx][rowIdx] = dxi_rowIdx/dX_columnIdx because of inverse function theorem.

Parameters:
  • )> getInverseJacobian(std::array<Vec3

  • FunctionSpaceFunction<MeshType,BasisFunctionType>::nDofsPerElement()> &geometryValues

  • element_no_t elementNo

  • std::array<double,MeshType::dim()>

Tensor2<MeshType::dim()> getInverseJacobian(element_no_t elementNo, std::array<double, MeshType::dim()> xi)

Compute the inverseJacobian that is needed to transform a gradient vector from parameter space to world space, for an element at a xi position. The following properties of the jacobian hold: jacobianParameterSpace[columnIdx][rowIdx] = dX_rowIdx/dxi_columnIdx inverseJacobianParameterSpace[columnIdx][rowIdx] = dxi_rowIdx/dX_columnIdx because of inverse function theorem.

Parameters:
  • )> getInverseJacobian(element_no_t elementNo

  • std::array<double,MeshType::dim()>

bool pointIsInElement(Vec3 point, element_no_t elementNo, std::array<double, 1> &xi, double xiTolerance = 1e-4)

Check if the point lies inside the element, if yes, return true and set xi to the value of the point.

Parameters:
  • Vec3 point

  • element_no_t elementNo

  • std::array<double,1> &xi

  • double xiTolerance = 1e-4

std::array<dof_no_t, FunctionSpaceFunction<MeshType, BasisFunctionType>::nDofsPerElement()> getElementDofNosLocal(element_no_t elementNo) const

Return an array of all dof nos. of the element, including ghost dofs (local dof nos).

void getElementDofNosLocalWithoutGhosts(element_no_t elementNo, std::vector<dof_no_t> &dofNosLocal) const

Fill a vector of all local dof nos. of the element, without ghost dofs.

Parameters:
  • element_no_t elementNo

  • std::vector<dof_no_t> &dofNosLocal

Vec3 getGeometry(node_no_t dofGlobalNo) const

Return the geometry field entry (node position for Lagrange elements) of a specific dof.

Parameters:

node_no_t dofGlobalNo

void getElementGeometry(element_no_t elementNoLocal, std::array<Vec3, FunctionSpaceBaseDim<MeshType::dim(), BasisFunctionType>::nDofsPerElement()> &values)

Get all geometry entries for an element.

Parameters:
  • element_no_t elementNoLocal

  • std::array<Vec3,FunctionSpaceBaseDim<MeshType::dim(),BasisFunctionType>::nDofsPerElement(

void extractSurfaceGeometry(const std::array<Vec3, FunctionSpaceBaseDim<MeshType::dim(), BasisFunctionType>::nDofsPerElement()> &geometryVolume, Mesh::face_t face, std::array<Vec3, FunctionSpaceBaseDim<MeshType::dim() - 1, BasisFunctionType>::nNodesPerElement()> &geometrySurface)

From the function space geometry, extract geometry data for a surface with has one lower dimensionality, only the nodal dofs are extracted, also for Hermite.

Parameters:
  • const std::array<Vec3,FunctionSpaceBaseDim<MeshType::dim(),BasisFunctionType>::nDofsPerElement()> &geometryVolume

  • Mesh::face_t face

  • std::array<Vec3,FunctionSpaceBaseDim<MeshType::dim()-1,BasisFunctionType>::nNodesPerElement(

GeometryFieldType &geometryField()

Return the internal geometry field variable.