Cells

Cells represent the basic building blocks of the MPM grid. They are used to compute the material point kernels and to interpolate the field solutions back to the material points. Cells are used to compute the residual and stiffness matrices for the global system of equations.

Implementing your own cells can be done easily by subclassing from the abstract base class CellBase.

class edelweissmeshfree.cells.base.cell.CellBase(cellType, cellNumber, nodes)[source]

MPM cells in EdelweissMPM should be derived from this base class in order to follow the general interface.

EdelweissMPM expects the layout of the internal and external load vectors, P, PExt, (and the stiffness) to be of the form

[ node 1 - dofs field 1,
  node 1 - dofs field 2,
  node 1 - ... ,
  node 1 - dofs field n,
  node 2 - dofs field 1,
  ... ,
  node N - dofs field n].
Parameters:
  • cellType (str) – A string identifying the requested element formulation.

  • cellNumber (int) – A unique integer label used for all kinds of purposes.

  • nodes (list[Node]) – The list of Nodes assigned to this cell.

Attributes:
assignedMaterialPoints

The shape of the element in Ensight Gold notation.

cellNumber

The unique number of this cell

dofIndicesPermutation

The permutation pattern for the residual vector and the stiffness matrix to aggregate all entries in order to resemble the defined fields nodewise.

ensightType

The shape of the element in Ensight Gold notation.

fields

The list of fields per grid nodes.

nDof

The total number of degrees of freedom this cell has

nNodes

The list of nodes this cell holds

Methods

assignMaterialPoints(materialPoints)

Assign a list of material points which are currently residing within the cell.

computeBodyLoad(loadType, load, P, K, ...)

Compute bulk loads (body loads) for given time based for all assigned MaterialPoints.

computeDistributedLoad(loadType, surfaceID, ...)

Compute distributed (surface) a for given time based for a specific assigned material point.

computeMaterialPointKernels(dU, P, K, ...)

Evaluate the kernel residual and stiffness for given time based for all assigned MaterialPoints.

getBoundingBox()

Get the bounding box min.

getCoordinatesAtCenter()

Compute the underlying MarmotElement centroid coordinates.

getInterpolationVector(coordinate)

Get the interpolation vector for a given global coordinate.

getVIJContributionSize()

Return the number of entries this entity contributes to the VIJ (COO) system matrix.

initializeVIJContribution(idcs, I_, J_, offset)

Fill the global I and J index arrays for this constraint's VIJ contribution.

interpolateSolutionContributionToMaterialPoints(dU)

Interpolate field solutions to the assigned MaterialPoints.

isCoordinateInCell(coordinate)

Check if a given coordinate is located within this cell.

abstractmethod assignMaterialPoints(materialPoints)[source]

Assign a list of material points which are currently residing within the cell.

Parameters:

materialPoints (list[MaterialPointBase]) – The list of material points to be assigned.

abstractmethod computeBodyLoad(loadType, load, P, K, timeTotal, dTime)[source]

Compute bulk loads (body loads) for given time based for all assigned MaterialPoints.

Parameters:
  • loadType (str) – The type of load to be computed (e.g., ‘bodyforce’)

  • load (ndarray) – The float (vector) describing the load.

  • P (ndarray) – The external load vector to be defined.

  • K (ndarray) – The stiffness matrix to be defined.

  • timeTotal (float) – The current total time.

  • dTime (float) – The time increment.

abstractmethod computeDistributedLoad(loadType, surfaceID, materialPoint, load, P, K, timeTotal, dTime)[source]

Compute distributed (surface) a for given time based for a specific assigned material point.

Parameters:
  • loadType (str) – The type of load to be computed (e.g., ‘pressure’)

  • surfaceID (int) – The ID describing the surface of the material point.

  • materialPoint (MaterialPointBase) – The specific (already assigned) mateiral point for computing the surface load.

  • load (ndarray) – The float (vector) describing the load.

  • P (ndarray) – The external load vector to be defined.

  • K (ndarray) – The stiffness matrix to be defined.

  • timeTotal (float) – The current total time.

  • dTime (float) – The time increment.

abstractmethod computeMaterialPointKernels(dU, P, K, timeTotal, dTime)[source]

Evaluate the kernel residual and stiffness for given time based for all assigned MaterialPoints.

Parameters:
  • dU (ndarray) – The current solution increment.

  • P (ndarray) – The external load vector to be defined.

  • K (ndarray) – The stiffness matrix to be defined.

  • timeTotal (float) – The current total time.

  • dTime (float) – The time increment.

abstractmethod getBoundingBox()[source]

Get the bounding box min. and max. of a cell

Returns:

The tuple containing the min. and max. coordinates of the bounding box.

Return type:

tuple

abstractmethod getCoordinatesAtCenter()[source]

Compute the underlying MarmotElement centroid coordinates.

Returns:

The cells’s central coordinates.

Return type:

np.ndarray

abstractmethod getInterpolationVector(coordinate)[source]

Get the interpolation vector for a given global coordinate.

Returns:

The interpolation vector for all nodes.

Return type:

np.ndarray

Parameters:

coordinate (ndarray)

abstractmethod interpolateSolutionContributionToMaterialPoints(dU)[source]

Interpolate field solutions to the assigned MaterialPoints.

Parameters:

dU (ndarray) – The current solution vector contribution increment for all fields.

abstractmethod isCoordinateInCell(coordinate)[source]

Check if a given coordinate is located within this cell.

Returns:

The truth value if this coordinate is located in the cell.

Return type:

bool

Parameters:

coordinate (ndarray)

abstract property assignedMaterialPoints: list[MaterialPointBase]

The shape of the element in Ensight Gold notation.

abstract property cellNumber: int

The unique number of this cell

abstract property dofIndicesPermutation: ndarray

The permutation pattern for the residual vector and the stiffness matrix to aggregate all entries in order to resemble the defined fields nodewise.

abstract property ensightType: str

The shape of the element in Ensight Gold notation.

abstract property fields: list[list[str]]

The list of fields per grid nodes.

abstract property nDof: int

The total number of degrees of freedom this cell has

abstract property nNodes: int

The list of nodes this cell holds

MarmotCell class

class edelweissmeshfree.cells.marmotcell.marmotcell.MarmotCellWrapper

This cell as a wrapper for MarmotCells. It is an abstract class and cannot be used directly. Rather, it is used as a base class for the specific cell types, such as LagrangianMarmotCell and BSplineMarmotCell.

For the documentation of MarmotCells, please refer to Marmot.

Parameters:
  • cellType – The Marmot element which should be represented, e.g., CPE4.

  • cellNumber – The (unique) number of this cell.

  • nodes – The list of nodes for this Cell.

Attributes:
assignedMaterialPoints
cellNumber
cellType
dofIndicesPermutation
ensightType
fields
nDof
nNodes
nodes

Methods

computeMaterialPointKernels(dUc, Pc, Kc, ...)

Evaluate residual and stiffness for given time, field, and field increment.

getVIJContributionSize()

Return the number of entries this entity contributes to the VIJ (COO) system matrix.

initializeVIJContribution(idcs, I_, J_, offset)

Initialize the I and J arrays for the VIJ (COO) system matrix assembly.

assignMaterialPoints

computeBodyLoad

computeConsistentInertia

computeDistributedLoad

computeLumpedInertia

getBoundingBox

getInterpolationVector

interpolateFieldsToMaterialPoints

isCoordinateInCell

computeMaterialPointKernels(dUc, Pc, Kc, timeNew, dTime)

Evaluate residual and stiffness for given time, field, and field increment.

getVIJContributionSize()

Return the number of entries this entity contributes to the VIJ (COO) system matrix.

Return type:

int

initializeVIJContribution(idcs, I_, J_, offset)

Initialize the I and J arrays for the VIJ (COO) system matrix assembly.

Parameters:
  • idcs (ndarray)

  • I_ (ndarray)

  • J_ (ndarray)

  • offset (int)

Return type:

None

LagrangianMarmotCell class

class edelweissmeshfree.cells.marmotcell.lagrangianmarmotcell.LagrangianMarmotCellWrapper

This class is a wrapper for Lagrangian Marmot cells. It is responsible for creating the Marmot cell and for managing the material points associated with the cell.

Parameters:
  • cellType – The type of the cell, e.g., CPE4.

  • cellNumber – The number of the cell.

  • nodes – The nodes of the cell.

Attributes:
assignedMaterialPoints
cellNumber
cellType
dofIndicesPermutation
ensightType
fields
nDof
nNodes
nodes

Methods

computeMaterialPointKernels(dUc, Pc, Kc, ...)

Evaluate residual and stiffness for given time, field, and field increment.

getVIJContributionSize()

Return the number of entries this entity contributes to the VIJ (COO) system matrix.

initializeVIJContribution(idcs, I_, J_, offset)

Initialize the I and J arrays for the VIJ (COO) system matrix assembly.

assignMaterialPoints

computeBodyLoad

computeConsistentInertia

computeDistributedLoad

computeLumpedInertia

getBoundingBox

getInterpolationVector

interpolateFieldsToMaterialPoints

isCoordinateInCell

BSplineMarmotCell class

class edelweissmeshfree.cells.marmotcell.bsplinemarmotcell.BSplineMarmotCellWrapper

This class is a wrapper for the MarmotCell class. It is used to create a BSpline-MarmotCell object from a set of nodes and knot vectors. The MarmotCell object is then used to evaluate the shape functions and their derivatives at a given point.

Attributes:
assignedMaterialPoints
cellNumber
cellType
dofIndicesPermutation
ensightType
fields
nDof
nNodes
nodes

Methods

computeMaterialPointKernels(dUc, Pc, Kc, ...)

Evaluate residual and stiffness for given time, field, and field increment.

getVIJContributionSize()

Return the number of entries this entity contributes to the VIJ (COO) system matrix.

initializeVIJContribution(idcs, I_, J_, offset)

Initialize the I and J arrays for the VIJ (COO) system matrix assembly.

assignMaterialPoints

computeBodyLoad

computeConsistentInertia

computeDistributedLoad

computeLumpedInertia

getBoundingBox

getInterpolationVector

interpolateFieldsToMaterialPoints

isCoordinateInCell

PythonCell class

A pure Python implementation of a bilinear Quad4 cell for 2D displacement problems.

This cell uses standard bilinear shape functions on a quadrilateral to interpolate displacements from grid nodes to material points, and to assemble residual vectors and stiffness matrices from material point contributions.

class edelweissmeshfree.cells.pythoncell.cell.PythonCell(cellType, cellNumber, nodes)[source]

MPM cells in EdelweissMPM should be derived from this base class in order to follow the general interface.

EdelweissMPM expects the layout of the internal and external load vectors, P, PExt, (and the stiffness) to be of the form

[ node 1 - dofs field 1,
  node 1 - dofs field 2,
  node 1 - ... ,
  node 1 - dofs field n,
  node 2 - dofs field 1,
  ... ,
  node N - dofs field n].
Parameters:
  • cellType (str) – A string identifying the requested element formulation.

  • cellNumber (int) – A unique integer label used for all kinds of purposes.

  • nodes (list[Node]) – The list of Nodes assigned to this cell.

Attributes:
assignedMaterialPoints

The shape of the element in Ensight Gold notation.

cellNumber

The unique number of this cell

dofIndicesPermutation

The permutation pattern for the residual vector and the stiffness matrix to aggregate all entries in order to resemble the defined fields nodewise.

ensightType

The shape of the element in Ensight Gold notation.

fields

The list of fields per grid nodes.

nDof

The total number of degrees of freedom this cell has

nNodes

The list of nodes this cell holds

Methods

assignMaterialPoints(materialPoints)

Assign a list of material points which are currently residing within the cell.

computeBodyLoad(loadType, load, P, K, ...)

Compute body load contribution (e.g., gravity).

computeDistributedLoad(loadType, surfaceID, ...)

Compute distributed (surface) load for a specific material point.

computeMaterialPointKernels(dU, P, K, ...)

Compute residual and stiffness contributions from all assigned material points.

getBoundingBox()

Get the bounding box min.

getCoordinatesAtCenter()

Compute the underlying MarmotElement centroid coordinates.

getInterpolationVector(coordinate)

Get the interpolation vector for a given global coordinate.

getVIJContributionSize()

Return the number of entries this entity contributes to the VIJ (COO) system matrix.

initializeVIJContribution(idcs, I_, J_, offset)

Fill the global I and J index arrays for this constraint's VIJ contribution.

interpolateFieldsToMaterialPoints(dU)

Interpolate nodal displacement increments to material points.

interpolateSolutionContributionToMaterialPoints(dU)

Alias for interpolateFieldsToMaterialPoints (base class interface).

isCoordinateInCell(coordinate)

Check if a given coordinate is located within this cell.

_get_B_matrix(dNdx)[source]

Compute the strain-displacement matrix B for plane strain/stress.

Parameters:

dNdx – Shape function gradients (4, 2).

Returns:

(4, 8) strain-displacement matrix for Voigt notation [eps_xx, eps_yy, eps_zz, gamma_xy].

Return type:

B

assignMaterialPoints(materialPoints)[source]

Assign a list of material points which are currently residing within the cell.

Parameters:

materialPoints (list) – The list of material points to be assigned.

computeBodyLoad(loadType, load, P, K, timeTotal, dTime)[source]

Compute body load contribution (e.g., gravity).

Parameters:
  • loadType (str)

  • load (ndarray)

  • P (ndarray)

  • K (ndarray)

  • timeTotal (float)

  • dTime (float)

computeDistributedLoad(loadType, surfaceID, materialPoint, load, P, K, timeTotal, dTime)[source]

Compute distributed (surface) load for a specific material point.

Parameters:
  • loadType (str)

  • surfaceID (int)

  • materialPoint (MaterialPointBase)

  • load (ndarray)

  • P (ndarray)

  • K (ndarray)

  • timeTotal (float)

  • dTime (float)

computeMaterialPointKernels(dU, P, K, timeTotal, dTime)[source]

Compute residual and stiffness contributions from all assigned material points.

Parameters:
  • dU (ndarray) – The current solution increment (8,).

  • P (ndarray) – The internal load vector to be filled (8,).

  • K (ndarray) – The stiffness matrix to be filled (64,) in row-major order.

  • timeTotal (float) – The current total time.

  • dTime (float) – The time increment.

getBoundingBox()[source]

Get the bounding box min. and max. of a cell

Returns:

The tuple containing the min. and max. coordinates of the bounding box.

Return type:

tuple

getCoordinatesAtCenter()[source]

Compute the underlying MarmotElement centroid coordinates.

Returns:

The cells’s central coordinates.

Return type:

np.ndarray

getInterpolationVector(coordinate)[source]

Get the interpolation vector for a given global coordinate.

Returns:

The interpolation vector for all nodes.

Return type:

np.ndarray

Parameters:

coordinate (ndarray)

interpolateFieldsToMaterialPoints(dU)[source]

Interpolate nodal displacement increments to material points.

Parameters:

dU (ndarray) – The nodal displacement increment vector (8,) = [u1x, u1y, u2x, u2y, u3x, u3y, u4x, u4y].

interpolateSolutionContributionToMaterialPoints(dU)[source]

Alias for interpolateFieldsToMaterialPoints (base class interface).

Parameters:

dU (ndarray)

isCoordinateInCell(coordinate)[source]

Check if a given coordinate is located within this cell.

Returns:

The truth value if this coordinate is located in the cell.

Return type:

bool

Parameters:

coordinate (ndarray)

property assignedMaterialPoints: list

The shape of the element in Ensight Gold notation.

property cellNumber: int

The unique number of this cell

property dofIndicesPermutation: ndarray

The permutation pattern for the residual vector and the stiffness matrix to aggregate all entries in order to resemble the defined fields nodewise.

property ensightType: str

The shape of the element in Ensight Gold notation.

property fields: list[list[str]]

The list of fields per grid nodes.

property nDof: int

The total number of degrees of freedom this cell has

property nNodes: int

The list of nodes this cell holds

edelweissmeshfree.cells.pythoncell.cell._compute_shape_functions_and_gradients(xi, eta, node_coords)[source]

Compute shape functions and their spatial gradients at a given parametric point.

Parameters:
  • xi – Parametric coordinate xi in [-1, 1].

  • eta – Parametric coordinate eta in [-1, 1].

  • node_coords – (4, 2) array of node coordinates.

Returns:

  • N – Shape function values (4,).

  • dNdx – Shape function gradients w.r.t. physical coordinates (4, 2).

  • detJ – Determinant of the Jacobian.

edelweissmeshfree.cells.pythoncell.cell._global_to_parametric(coordinate, node_coords)[source]

Map a global coordinate to parametric coordinates using Newton iteration.

Parameters:
  • coordinate – Global (x, y) coordinate.

  • node_coords – (4, 2) array of node coordinates.

Returns:

Parametric coordinates.

Return type:

xi, eta

Theory (Quad4 interpolation and cell kernels)

The Python cell implementation uses a bilinear Quad4 interpolation in the parent space \((\xi,\eta)\in[-1,1]^2\):

\[N_1 = \frac{1}{4}(1-\xi)(1-\eta),\quad N_2 = \frac{1}{4}(1+\xi)(1-\eta),\quad N_3 = \frac{1}{4}(1+\xi)(1+\eta),\quad N_4 = \frac{1}{4}(1-\xi)(1+\eta).\]

With nodal displacement vector \(\mathbf{u}_n=[u_{1x},u_{1y},u_{2x},u_{2y},u_{3x},u_{3y},u_{4x},u_{4y}]^\mathsf{T}\), the material point displacement increment is interpolated as:

\[\Delta\mathbf{u}_{mp} = \sum_{a=1}^{4} N_a\,\Delta\mathbf{u}_a.\]

The strain increment follows from the standard strain-displacement matrix:

\[\Delta\boldsymbol{\varepsilon}_{mp} = \mathbf{B}\,\Delta\mathbf{u}_n,\]

where \(\mathbf{B}\) is assembled from the spatial gradients \(\nabla N_a\) (obtained using the Jacobian mapping from parent to physical space).

Cell residual and tangent contributions are assembled from material point stress \(\boldsymbol{\sigma}\) and consistent tangent \(\mathbf{C}\):

\[\mathbf{P}_{int} \mathrel{+}= \mathbf{B}^\mathsf{T}\boldsymbol{\sigma}\,V_{mp}, \qquad \mathbf{K} \mathrel{+}= \mathbf{B}^\mathsf{T}\mathbf{C}\mathbf{B}\,V_{mp}.\]