Skip to content

FE_transformation_2d

Module: fe_transformation_2d.py

This module provides the abstract base class for all 2D finite element transformations. It defines the interface for mapping between reference and physical coordinates in two-dimensional finite element analysis.

Classes:

Name Description
FETransformation2D

Abstract base class for 2D finite element transformations

Dependencies
  • abc: For abstract base class functionality
  • quad_affine: For affine transformation implementations
  • quad_bilinear: For bilinear transformation implementations
Key Features
  • Abstract interface for coordinate transformations
  • Reference to physical domain mapping
  • Jacobian matrix computation
  • Support for different element geometries
  • Cell geometry specification interface
  • Systematic transformation validation
Authors

Thivin Anandh D (https://thivinanandh.github.io)

Version Info

27/Dec/2024: Initial version - Thivin Anandh D

References

None

FETransforamtion2D

A base class for 2D finite element transformations.

This abstract class defines the interface for mapping between reference and physical coordinates in 2D finite element analysis. Implementations must provide specific transformation rules for different element types.

Methods:

Name Description
set_cell

Sets the physical coordinates of the element vertices. Must be implemented by derived classes.

get_original_from_ref

Maps coordinates from reference to physical domain. Must be implemented by derived classes.

get_jacobian

Computes the Jacobian matrix of the transformation. Must be implemented by derived classes.

Example

class QuadTransform(FETransformation2D): ... def set_cell(self, vertices): ... self.vertices = vertices ... def get_original_from_ref(self, xi:np.ndarray, eta:np.ndarray) -> np.ndarray: ... # Implementation for quad element ... pass ... def get_jacobian(self, xi: np.ndarray, eta:np.ndarray) -> np.ndarray: ... # Implementation for quad element ... pass

Notes
  • Reference domain is typically [-1,1] × [-1,1]
  • Transformations must be invertible
  • Implementations should handle element distortion
  • Jacobian is used for both mapping and integration
Source code in scirex/core/sciml/fe/fe_transformation_2d.py
class FETransforamtion2D:
    """
    A base class for 2D finite element transformations.

    This abstract class defines the interface for mapping between reference and physical
    coordinates in 2D finite element analysis. Implementations must provide specific
    transformation rules for different element types.

    Attributes:
        None

    Methods:
        set_cell():
            Sets the physical coordinates of the element vertices.
            Must be implemented by derived classes.

        get_original_from_ref(xi, eta):
            Maps coordinates from reference to physical domain.
            Must be implemented by derived classes.

        get_jacobian(xi, eta):
            Computes the Jacobian matrix of the transformation.
            Must be implemented by derived classes.

    Example:
        >>> class QuadTransform(FETransformation2D):
        ...     def set_cell(self, vertices):
        ...         self.vertices = vertices
        ...     def get_original_from_ref(self, xi:np.ndarray, eta:np.ndarray) -> np.ndarray:
        ...         # Implementation for quad element
        ...         pass
        ...     def get_jacobian(self, xi: np.ndarray, eta:np.ndarray) -> np.ndarray:
        ...         # Implementation for quad element
        ...         pass

    Notes:
        - Reference domain is typically [-1,1] × [-1,1]
        - Transformations must be invertible
        - Implementations should handle element distortion
        - Jacobian is used for both mapping and integration
    """

    def __init__(self):
        """
        Constructor for the FETransforamtion2D class.
        """

    @abstractmethod
    def set_cell(self):
        """
        Set the cell coordinates, which will be used to calculate the Jacobian and actual values.

        :return: None
        """

    @abstractmethod
    def get_original_from_ref(self, xi: np.ndarray, eta: np.ndarray) -> np.ndarray:
        """
        This method returns the original coordinates from the reference coordinates.

        Args:
            xi (np.ndarray): The xi coordinate.
            eta (np.ndarray): The eta coordinate.

        Returns:
            np.ndarray: Returns the transformed original coordinates from the reference coordinates.
        """

    @abstractmethod
    def get_jacobian(self, xi: np.ndarray, eta: np.ndarray) -> np.ndarray:
        """
        This method returns the Jacobian of the transformation.

        Args:
            xi (np.ndarray): The xi coordinate.
            eta (np.ndarray): The eta coordinate.

        Returns:
            np.ndarray: Returns the Jacobian of the transformation.
        """

__init__()

Constructor for the FETransforamtion2D class.

Source code in scirex/core/sciml/fe/fe_transformation_2d.py
def __init__(self):
    """
    Constructor for the FETransforamtion2D class.
    """

get_jacobian(xi, eta) abstractmethod

This method returns the Jacobian of the transformation.

Parameters:

Name Type Description Default
xi ndarray

The xi coordinate.

required
eta ndarray

The eta coordinate.

required

Returns:

Type Description
ndarray

np.ndarray: Returns the Jacobian of the transformation.

Source code in scirex/core/sciml/fe/fe_transformation_2d.py
@abstractmethod
def get_jacobian(self, xi: np.ndarray, eta: np.ndarray) -> np.ndarray:
    """
    This method returns the Jacobian of the transformation.

    Args:
        xi (np.ndarray): The xi coordinate.
        eta (np.ndarray): The eta coordinate.

    Returns:
        np.ndarray: Returns the Jacobian of the transformation.
    """

get_original_from_ref(xi, eta) abstractmethod

This method returns the original coordinates from the reference coordinates.

Parameters:

Name Type Description Default
xi ndarray

The xi coordinate.

required
eta ndarray

The eta coordinate.

required

Returns:

Type Description
ndarray

np.ndarray: Returns the transformed original coordinates from the reference coordinates.

Source code in scirex/core/sciml/fe/fe_transformation_2d.py
@abstractmethod
def get_original_from_ref(self, xi: np.ndarray, eta: np.ndarray) -> np.ndarray:
    """
    This method returns the original coordinates from the reference coordinates.

    Args:
        xi (np.ndarray): The xi coordinate.
        eta (np.ndarray): The eta coordinate.

    Returns:
        np.ndarray: Returns the transformed original coordinates from the reference coordinates.
    """

set_cell() abstractmethod

Set the cell coordinates, which will be used to calculate the Jacobian and actual values.

:return: None

Source code in scirex/core/sciml/fe/fe_transformation_2d.py
@abstractmethod
def set_cell(self):
    """
    Set the cell coordinates, which will be used to calculate the Jacobian and actual values.

    :return: None
    """

QuadAffin

Bases: FETransforamtion2D

Implements affine transformation for quadrilateral elements.

This class provides methods to transform between reference and physical quadrilateral elements using affine mapping. It handles coordinate transformations, Jacobian computations, and derivative mappings.

Attributes:

Name Type Description
co_ordinates

Array of physical element vertex coordinates Shape: (4, 2) for 2D quadrilateral

x0, (x1, x2, x3)

x-coordinates of vertices

y0, (y1, y2, y3)

y-coordinates of vertices

xc0, (xc1, xc2)

x-coordinate transformation coefficients

yc0, (yc1, yc2)

y-coordinate transformation coefficients

detjk (yc1, yc2)

Determinant of the Jacobian

rec_detjk (yc1, yc2)

Reciprocal of Jacobian determinant

Example

coords = np.array([[0,0], [1,0], [1,1], [0,1]]) quad = QuadAffin(coords) ref_point = np.array([0.5, 0.5]) physical_point = quad.get_original_from_ref(*ref_point)

Note

The implementation assumes counterclockwise vertex ordering and non-degenerate quadrilateral elements.

References

[1] ParMooN Project: QuadAffine.C implementation

Source code in scirex/core/sciml/fe/quad_affine.py
class QuadAffin(FETransforamtion2D):
    """
    Implements affine transformation for quadrilateral elements.

    This class provides methods to transform between reference and physical
    quadrilateral elements using affine mapping. It handles coordinate
    transformations, Jacobian computations, and derivative mappings.

    Attributes:
        co_ordinates: Array of physical element vertex coordinates
            Shape: (4, 2) for 2D quadrilateral
        x0, x1, x2, x3: x-coordinates of vertices
        y0, y1, y2, y3: y-coordinates of vertices
        xc0, xc1, xc2: x-coordinate transformation coefficients
        yc0, yc1, yc2: y-coordinate transformation coefficients
        detjk: Determinant of the Jacobian
        rec_detjk: Reciprocal of Jacobian determinant

    Example:
        >>> coords = np.array([[0,0], [1,0], [1,1], [0,1]])
        >>> quad = QuadAffin(coords)
        >>> ref_point = np.array([0.5, 0.5])
        >>> physical_point = quad.get_original_from_ref(*ref_point)

    Note:
        The implementation assumes counterclockwise vertex ordering and
        non-degenerate quadrilateral elements.

    References:
        [1] ParMooN Project: QuadAffine.C implementation
    """

    def __init__(self, co_ordinates: np.ndarray) -> None:
        """
        Constructor for the QuadAffin class.

        Args:
            co_ordinates: Array of physical element vertex coordinates
                Shape: (4, 2) for 2D quadrilateral

        Returns:
            None
        """
        self.co_ordinates = co_ordinates
        self.set_cell()
        self.get_jacobian(
            0, 0
        )  # 0,0 is just a dummy value # this sets the jacobian and the inverse of the jacobian

    def set_cell(self):
        """
        Set the cell coordinates, which will be used to calculate the Jacobian and actual values.

        Args:
            None

        Returns:
            None
        """

        self.x0 = self.co_ordinates[0][0]
        self.x1 = self.co_ordinates[1][0]
        self.x2 = self.co_ordinates[2][0]
        self.x3 = self.co_ordinates[3][0]

        # get the y-coordinates of the cell
        self.y0 = self.co_ordinates[0][1]
        self.y1 = self.co_ordinates[1][1]
        self.y2 = self.co_ordinates[2][1]
        self.y3 = self.co_ordinates[3][1]

        self.xc0 = (self.x1 + self.x3) * 0.5
        self.xc1 = (self.x1 - self.x0) * 0.5
        self.xc2 = (self.x3 - self.x0) * 0.5

        self.yc0 = (self.y1 + self.y3) * 0.5
        self.yc1 = (self.y1 - self.y0) * 0.5
        self.yc2 = (self.y3 - self.y0) * 0.5

    def get_original_from_ref(self, xi: np.ndarray, eta: np.ndarray) -> np.ndarray:
        """
        Returns the original coordinates from the reference coordinates.

        Args:
            xi (np.ndarray): The xi coordinate.
            eta (np.ndarray): The eta coordinate.

        Returns:
            np.ndarray: The transformed original coordinates from the reference coordinates.
        """
        x = self.xc0 + self.xc1 * xi + self.xc2 * eta
        y = self.yc0 + self.yc1 * xi + self.yc2 * eta

        return np.array([x, y])

    def get_jacobian(self, xi: np.ndarray, eta: np.ndarray) -> np.ndarray:
        """
        Returns the Jacobian of the transformation.

        Args:
            xi (np.ndarray): The xi coordinate.
            eta (np.ndarray): The eta coordinate.

        Returns:
            np.ndarray: The Jacobian of the transformation.
        """
        self.detjk = self.xc1 * self.yc2 - self.xc2 * self.yc1
        self.rec_detjk = 1 / self.detjk

        return abs(self.detjk)

    def get_orig_from_ref_derivative(self, ref_gradx, ref_grady, xi, eta):
        """
        Returns the derivatives of the original coordinates with respect to the reference coordinates.

        Args:
            ref_gradx (np.ndarray): The reference gradient in the x-direction.
            ref_grady (np.ndarray): The reference gradient in the y-direction.
            xi (np.ndarray): The xi coordinate.
            eta (np.ndarray): The eta coordinate.

        Returns:
            tuple: The derivatives of the original coordinates with respect to the reference coordinates.
        """
        gradx_orig = np.zeros(ref_gradx.shape)
        grady_orig = np.zeros(ref_grady.shape)

        for i in range(ref_gradx.shape[0]):
            gradx_orig[i] = (
                self.yc2 * ref_gradx[i] - self.yc1 * ref_grady[i]
            ) * self.rec_detjk
            grady_orig[i] = (
                -self.xc2 * ref_gradx[i] + self.xc1 * ref_grady[i]
            ) * self.rec_detjk

        return gradx_orig, grady_orig

    def get_orig_from_ref_second_derivative(
        self, grad_xx_ref, grad_xy_ref, grad_yy_ref, xi, eta
    ):
        """
        Returns the second derivatives (xx, xy, yy) of the original coordinates with respect to the reference coordinates.

        Args:
            grad_xx_ref (np.ndarray): The reference second derivative in the x-direction.
            grad_xy_ref (np.ndarray): The reference second derivative in the xy-direction.
            grad_yy_ref (np.ndarray): The reference second derivative in the y-direction.
            xi (np.ndarray): The xi coordinate.
            eta (np.ndarray): The eta coordinate.

        Returns:
            tuple: The second derivatives (xx, xy, yy) of the original coordinates with respect to the reference coordinates.
        """
        GeoData = np.zeros((3, 3))
        Eye = np.identity(3)

        # Populate GeoData (assuming xc1, xc2, yc1, yc2 are defined)
        GeoData[0, 0] = self.xc1 * self.xc1
        GeoData[0, 1] = 2 * self.xc1 * self.yc1
        GeoData[0, 2] = self.yc1 * self.yc1
        GeoData[1, 0] = self.xc1 * self.xc2
        GeoData[1, 1] = self.yc1 * self.xc2 + self.xc1 * self.yc2
        GeoData[1, 2] = self.yc1 * self.yc2
        GeoData[2, 0] = self.xc2 * self.xc2
        GeoData[2, 1] = 2 * self.xc2 * self.yc2
        GeoData[2, 2] = self.yc2 * self.yc2

        # solve the linear system
        solution = np.linalg.solve(GeoData, Eye)

        # generate empty arrays for the original second derivatives
        grad_xx_orig = np.zeros(grad_xx_ref.shape)
        grad_xy_orig = np.zeros(grad_xy_ref.shape)
        grad_yy_orig = np.zeros(grad_yy_ref.shape)

        for j in range(grad_xx_ref.shape[0]):
            r20 = grad_xx_ref[j]
            r11 = grad_xy_ref[j]
            r02 = grad_yy_ref[j]

            grad_xx_orig[j] = (
                solution[0, 0] * r20 + solution[0, 1] * r11 + solution[0, 2] * r02
            )
            grad_xy_orig[j] = (
                solution[1, 0] * r20 + solution[1, 1] * r11 + solution[1, 2] * r02
            )
            grad_yy_orig[j] = (
                solution[2, 0] * r20 + solution[2, 1] * r11 + solution[2, 2] * r02
            )

        return grad_xx_orig, grad_xy_orig, grad_yy_orig

__init__(co_ordinates)

Constructor for the QuadAffin class.

Parameters:

Name Type Description Default
co_ordinates ndarray

Array of physical element vertex coordinates Shape: (4, 2) for 2D quadrilateral

required

Returns:

Type Description
None

None

Source code in scirex/core/sciml/fe/quad_affine.py
def __init__(self, co_ordinates: np.ndarray) -> None:
    """
    Constructor for the QuadAffin class.

    Args:
        co_ordinates: Array of physical element vertex coordinates
            Shape: (4, 2) for 2D quadrilateral

    Returns:
        None
    """
    self.co_ordinates = co_ordinates
    self.set_cell()
    self.get_jacobian(
        0, 0
    )  # 0,0 is just a dummy value # this sets the jacobian and the inverse of the jacobian

get_jacobian(xi, eta)

Returns the Jacobian of the transformation.

Parameters:

Name Type Description Default
xi ndarray

The xi coordinate.

required
eta ndarray

The eta coordinate.

required

Returns:

Type Description
ndarray

np.ndarray: The Jacobian of the transformation.

Source code in scirex/core/sciml/fe/quad_affine.py
def get_jacobian(self, xi: np.ndarray, eta: np.ndarray) -> np.ndarray:
    """
    Returns the Jacobian of the transformation.

    Args:
        xi (np.ndarray): The xi coordinate.
        eta (np.ndarray): The eta coordinate.

    Returns:
        np.ndarray: The Jacobian of the transformation.
    """
    self.detjk = self.xc1 * self.yc2 - self.xc2 * self.yc1
    self.rec_detjk = 1 / self.detjk

    return abs(self.detjk)

get_orig_from_ref_derivative(ref_gradx, ref_grady, xi, eta)

Returns the derivatives of the original coordinates with respect to the reference coordinates.

Parameters:

Name Type Description Default
ref_gradx ndarray

The reference gradient in the x-direction.

required
ref_grady ndarray

The reference gradient in the y-direction.

required
xi ndarray

The xi coordinate.

required
eta ndarray

The eta coordinate.

required

Returns:

Name Type Description
tuple

The derivatives of the original coordinates with respect to the reference coordinates.

Source code in scirex/core/sciml/fe/quad_affine.py
def get_orig_from_ref_derivative(self, ref_gradx, ref_grady, xi, eta):
    """
    Returns the derivatives of the original coordinates with respect to the reference coordinates.

    Args:
        ref_gradx (np.ndarray): The reference gradient in the x-direction.
        ref_grady (np.ndarray): The reference gradient in the y-direction.
        xi (np.ndarray): The xi coordinate.
        eta (np.ndarray): The eta coordinate.

    Returns:
        tuple: The derivatives of the original coordinates with respect to the reference coordinates.
    """
    gradx_orig = np.zeros(ref_gradx.shape)
    grady_orig = np.zeros(ref_grady.shape)

    for i in range(ref_gradx.shape[0]):
        gradx_orig[i] = (
            self.yc2 * ref_gradx[i] - self.yc1 * ref_grady[i]
        ) * self.rec_detjk
        grady_orig[i] = (
            -self.xc2 * ref_gradx[i] + self.xc1 * ref_grady[i]
        ) * self.rec_detjk

    return gradx_orig, grady_orig

get_orig_from_ref_second_derivative(grad_xx_ref, grad_xy_ref, grad_yy_ref, xi, eta)

Returns the second derivatives (xx, xy, yy) of the original coordinates with respect to the reference coordinates.

Parameters:

Name Type Description Default
grad_xx_ref ndarray

The reference second derivative in the x-direction.

required
grad_xy_ref ndarray

The reference second derivative in the xy-direction.

required
grad_yy_ref ndarray

The reference second derivative in the y-direction.

required
xi ndarray

The xi coordinate.

required
eta ndarray

The eta coordinate.

required

Returns:

Name Type Description
tuple

The second derivatives (xx, xy, yy) of the original coordinates with respect to the reference coordinates.

Source code in scirex/core/sciml/fe/quad_affine.py
def get_orig_from_ref_second_derivative(
    self, grad_xx_ref, grad_xy_ref, grad_yy_ref, xi, eta
):
    """
    Returns the second derivatives (xx, xy, yy) of the original coordinates with respect to the reference coordinates.

    Args:
        grad_xx_ref (np.ndarray): The reference second derivative in the x-direction.
        grad_xy_ref (np.ndarray): The reference second derivative in the xy-direction.
        grad_yy_ref (np.ndarray): The reference second derivative in the y-direction.
        xi (np.ndarray): The xi coordinate.
        eta (np.ndarray): The eta coordinate.

    Returns:
        tuple: The second derivatives (xx, xy, yy) of the original coordinates with respect to the reference coordinates.
    """
    GeoData = np.zeros((3, 3))
    Eye = np.identity(3)

    # Populate GeoData (assuming xc1, xc2, yc1, yc2 are defined)
    GeoData[0, 0] = self.xc1 * self.xc1
    GeoData[0, 1] = 2 * self.xc1 * self.yc1
    GeoData[0, 2] = self.yc1 * self.yc1
    GeoData[1, 0] = self.xc1 * self.xc2
    GeoData[1, 1] = self.yc1 * self.xc2 + self.xc1 * self.yc2
    GeoData[1, 2] = self.yc1 * self.yc2
    GeoData[2, 0] = self.xc2 * self.xc2
    GeoData[2, 1] = 2 * self.xc2 * self.yc2
    GeoData[2, 2] = self.yc2 * self.yc2

    # solve the linear system
    solution = np.linalg.solve(GeoData, Eye)

    # generate empty arrays for the original second derivatives
    grad_xx_orig = np.zeros(grad_xx_ref.shape)
    grad_xy_orig = np.zeros(grad_xy_ref.shape)
    grad_yy_orig = np.zeros(grad_yy_ref.shape)

    for j in range(grad_xx_ref.shape[0]):
        r20 = grad_xx_ref[j]
        r11 = grad_xy_ref[j]
        r02 = grad_yy_ref[j]

        grad_xx_orig[j] = (
            solution[0, 0] * r20 + solution[0, 1] * r11 + solution[0, 2] * r02
        )
        grad_xy_orig[j] = (
            solution[1, 0] * r20 + solution[1, 1] * r11 + solution[1, 2] * r02
        )
        grad_yy_orig[j] = (
            solution[2, 0] * r20 + solution[2, 1] * r11 + solution[2, 2] * r02
        )

    return grad_xx_orig, grad_xy_orig, grad_yy_orig

get_original_from_ref(xi, eta)

Returns the original coordinates from the reference coordinates.

Parameters:

Name Type Description Default
xi ndarray

The xi coordinate.

required
eta ndarray

The eta coordinate.

required

Returns:

Type Description
ndarray

np.ndarray: The transformed original coordinates from the reference coordinates.

Source code in scirex/core/sciml/fe/quad_affine.py
def get_original_from_ref(self, xi: np.ndarray, eta: np.ndarray) -> np.ndarray:
    """
    Returns the original coordinates from the reference coordinates.

    Args:
        xi (np.ndarray): The xi coordinate.
        eta (np.ndarray): The eta coordinate.

    Returns:
        np.ndarray: The transformed original coordinates from the reference coordinates.
    """
    x = self.xc0 + self.xc1 * xi + self.xc2 * eta
    y = self.yc0 + self.yc1 * xi + self.yc2 * eta

    return np.array([x, y])

set_cell()

Set the cell coordinates, which will be used to calculate the Jacobian and actual values.

Returns:

Type Description

None

Source code in scirex/core/sciml/fe/quad_affine.py
def set_cell(self):
    """
    Set the cell coordinates, which will be used to calculate the Jacobian and actual values.

    Args:
        None

    Returns:
        None
    """

    self.x0 = self.co_ordinates[0][0]
    self.x1 = self.co_ordinates[1][0]
    self.x2 = self.co_ordinates[2][0]
    self.x3 = self.co_ordinates[3][0]

    # get the y-coordinates of the cell
    self.y0 = self.co_ordinates[0][1]
    self.y1 = self.co_ordinates[1][1]
    self.y2 = self.co_ordinates[2][1]
    self.y3 = self.co_ordinates[3][1]

    self.xc0 = (self.x1 + self.x3) * 0.5
    self.xc1 = (self.x1 - self.x0) * 0.5
    self.xc2 = (self.x3 - self.x0) * 0.5

    self.yc0 = (self.y1 + self.y3) * 0.5
    self.yc1 = (self.y1 - self.y0) * 0.5
    self.yc2 = (self.y3 - self.y0) * 0.5

QuadBilinear

Bases: FETransforamtion2D

Implements bilinear transformation for quadrilateral elements.

This class provides methods to transform between reference and physical quadrilateral elements using bilinear mapping. It handles coordinate transformations, Jacobian computations, and derivative mappings for more general quadrilateral elements than affine transformations.

Attributes:

Name Type Description
co_ordinates

Array of physical element vertex coordinates Shape: (4, 2) for 2D quadrilateral

x0, (x1, x2, x3)

x-coordinates of vertices

y0, (y1, y2, y3)

y-coordinates of vertices

xc0, (xc1, xc2, xc3)

x-coordinate transformation coefficients

yc0, (yc1, yc2, yc3)

y-coordinate transformation coefficients

detjk

Determinant of the Jacobian matrix

Example

coords = np.array([[0,0], [1,0], [1.2,1], [0.2,1.1]]) quad = QuadBilinear(coords) ref_point = np.array([0.5, 0.5]) physical_point = quad.get_original_from_ref(*ref_point)

Note
  • Implementation assumes counterclockwise vertex ordering
  • Second derivatives computation is not fully implemented
  • Jacobian is computed point-wise due to non-constant nature of bilinear transformation
References

[1] ParMooN Project: QuadBilineare.C implementation

Source code in scirex/core/sciml/fe/quad_bilinear.py
class QuadBilinear(FETransforamtion2D):
    """
    Implements bilinear transformation for quadrilateral elements.

    This class provides methods to transform between reference and physical
    quadrilateral elements using bilinear mapping. It handles coordinate
    transformations, Jacobian computations, and derivative mappings for more
    general quadrilateral elements than affine transformations.

    Attributes:
        co_ordinates: Array of physical element vertex coordinates
            Shape: (4, 2) for 2D quadrilateral
        x0, x1, x2, x3: x-coordinates of vertices
        y0, y1, y2, y3: y-coordinates of vertices
        xc0, xc1, xc2, xc3: x-coordinate transformation coefficients
        yc0, yc1, yc2, yc3: y-coordinate transformation coefficients
        detjk: Determinant of the Jacobian matrix

    Example:
        >>> coords = np.array([[0,0], [1,0], [1.2,1], [0.2,1.1]])
        >>> quad = QuadBilinear(coords)
        >>> ref_point = np.array([0.5, 0.5])
        >>> physical_point = quad.get_original_from_ref(*ref_point)

    Note:
        - Implementation assumes counterclockwise vertex ordering
        - Second derivatives computation is not fully implemented
        - Jacobian is computed point-wise due to non-constant nature
        of bilinear transformation

    References:
        [1] ParMooN Project: QuadBilineare.C implementation
    """

    def __init__(self, co_ordinates: np.ndarray) -> None:
        """
        Constructor for the QuadBilinear class.

        Args:
            co_ordinates: Array of physical element vertex coordinates
                Shape: (4, 2) for 2D quadrilateral

        Returns:
            None
        """
        self.co_ordinates = co_ordinates
        self.set_cell()
        self.detjk = None  # Jacobian of the transformation

    def set_cell(self):
        """
        Set the cell coordinates, which will be used as intermediate values to calculate the Jacobian and actual values.

        Args:
            None

        Returns:
            None
        """
        self.x0 = self.co_ordinates[0][0]
        self.x1 = self.co_ordinates[1][0]
        self.x2 = self.co_ordinates[2][0]
        self.x3 = self.co_ordinates[3][0]

        # get the y-coordinates of the cell
        self.y0 = self.co_ordinates[0][1]
        self.y1 = self.co_ordinates[1][1]
        self.y2 = self.co_ordinates[2][1]
        self.y3 = self.co_ordinates[3][1]

        self.xc0 = (self.x0 + self.x1 + self.x2 + self.x3) * 0.25
        self.xc1 = (-self.x0 + self.x1 + self.x2 - self.x3) * 0.25
        self.xc2 = (-self.x0 - self.x1 + self.x2 + self.x3) * 0.25
        self.xc3 = (self.x0 - self.x1 + self.x2 - self.x3) * 0.25

        self.yc0 = (self.y0 + self.y1 + self.y2 + self.y3) * 0.25
        self.yc1 = (-self.y0 + self.y1 + self.y2 - self.y3) * 0.25
        self.yc2 = (-self.y0 - self.y1 + self.y2 + self.y3) * 0.25
        self.yc3 = (self.y0 - self.y1 + self.y2 - self.y3) * 0.25

    def get_original_from_ref(self, xi: np.ndarray, eta: np.ndarray) -> np.ndarray:
        """
        This method returns the original coordinates from the reference coordinates.

        Args:
            xi (np.ndarray): The xi coordinate.
            eta (np.ndarray): The eta coordinate.

        Returns:
            np.ndarray: Returns the transformed original coordinates from the reference coordinates.
        """
        x = self.xc0 + self.xc1 * xi + self.xc2 * eta + self.xc3 * xi * eta
        y = self.yc0 + self.yc1 * xi + self.yc2 * eta + self.yc3 * xi * eta

        return np.array([x, y], dtype=np.float64)

    def get_jacobian(self, xi: np.ndarray, eta: np.ndarray) -> np.ndarray:
        """
        This method returns the Jacobian of the transformation.

        Args:
            xi (np.ndarray): The xi coordinate.
            eta (np.ndarray): The eta coordinate.

        Returns:
            np.ndarray: Returns the Jacobian of the transformation.
        """
        self.detjk = abs(
            (self.xc1 + self.xc3 * eta) * (self.yc2 + self.yc3 * xi)
            - (self.xc2 + self.xc3 * xi) * (self.yc1 + self.yc3 * eta)
        )
        return self.detjk

    def get_orig_from_ref_derivative(
        self,
        ref_gradx: np.ndarray,
        ref_grady: np.ndarray,
        xi: np.ndarray,
        eta: np.ndarray,
    ) -> np.ndarray:
        """
        This method returns the derivatives of the original coordinates with respect to the reference coordinates.

        Args:
            ref_gradx (np.ndarray): The derivative of the xi coordinate in the reference element.
            ref_grady (np.ndarray): The derivative of the eta coordinate in the reference element.
            xi (np.ndarray): The xi coordinate.
            eta (np.ndarray): The eta coordinate.

        Returns:
            np.ndarray: The derivatives of the original coordinates [x, y] with respect to the reference coordinates.

        """
        n_test = ref_gradx.shape[0]
        gradx_orig = np.zeros(ref_gradx.shape, dtype=np.float64)
        grady_orig = np.zeros(ref_grady.shape, dtype=np.float64)

        for j in range(n_test):
            Xi = xi
            Eta = eta
            rec_detjk = 1 / (
                (self.xc1 + self.xc3 * Eta) * (self.yc2 + self.yc3 * Xi)
                - (self.xc2 + self.xc3 * Xi) * (self.yc1 + self.yc3 * Eta)
            )
            gradx_orig[j] = (
                (self.yc2 + self.yc3 * Xi) * ref_gradx[j]
                - (self.yc1 + self.yc3 * Eta) * ref_grady[j]
            ) * rec_detjk
            grady_orig[j] = (
                -(self.xc2 + self.xc3 * Xi) * ref_gradx[j]
                + (self.xc1 + self.xc3 * Eta) * ref_grady[j]
            ) * rec_detjk

        return gradx_orig, grady_orig

    def get_orig_from_ref_second_derivative(
        self,
        grad_xx_ref: np.ndarray,
        grad_xy_ref: np.ndarray,
        grad_yy_ref: np.ndarray,
        xi: np.ndarray,
        eta: np.ndarray,
    ):
        """
        This method returns the second derivatives of the original coordinates with respect to the reference coordinates.

        Args:
            grad_xx_ref (np.ndarray): The second derivative of the xi coordinate in the reference element.
            grad_xy_ref (np.ndarray): The second derivative of the xi and eta coordinates in the reference element.
            grad_yy_ref (np.ndarray): The second derivative of the eta coordinate in the reference element.
            xi (np.ndarray): The xi coordinate.
            eta (np.ndarray): The eta coordinate.

        Note:
            Second derivative calculations are not fully implemented in this method. Needs further development.
        """
        # print(" Error : Second Derivative not implemented -- Ignore this error, if second derivative is not required ")
        return grad_xx_ref, grad_xy_ref, grad_yy_ref

__init__(co_ordinates)

Constructor for the QuadBilinear class.

Parameters:

Name Type Description Default
co_ordinates ndarray

Array of physical element vertex coordinates Shape: (4, 2) for 2D quadrilateral

required

Returns:

Type Description
None

None

Source code in scirex/core/sciml/fe/quad_bilinear.py
def __init__(self, co_ordinates: np.ndarray) -> None:
    """
    Constructor for the QuadBilinear class.

    Args:
        co_ordinates: Array of physical element vertex coordinates
            Shape: (4, 2) for 2D quadrilateral

    Returns:
        None
    """
    self.co_ordinates = co_ordinates
    self.set_cell()
    self.detjk = None  # Jacobian of the transformation

get_jacobian(xi, eta)

This method returns the Jacobian of the transformation.

Parameters:

Name Type Description Default
xi ndarray

The xi coordinate.

required
eta ndarray

The eta coordinate.

required

Returns:

Type Description
ndarray

np.ndarray: Returns the Jacobian of the transformation.

Source code in scirex/core/sciml/fe/quad_bilinear.py
def get_jacobian(self, xi: np.ndarray, eta: np.ndarray) -> np.ndarray:
    """
    This method returns the Jacobian of the transformation.

    Args:
        xi (np.ndarray): The xi coordinate.
        eta (np.ndarray): The eta coordinate.

    Returns:
        np.ndarray: Returns the Jacobian of the transformation.
    """
    self.detjk = abs(
        (self.xc1 + self.xc3 * eta) * (self.yc2 + self.yc3 * xi)
        - (self.xc2 + self.xc3 * xi) * (self.yc1 + self.yc3 * eta)
    )
    return self.detjk

get_orig_from_ref_derivative(ref_gradx, ref_grady, xi, eta)

This method returns the derivatives of the original coordinates with respect to the reference coordinates.

Parameters:

Name Type Description Default
ref_gradx ndarray

The derivative of the xi coordinate in the reference element.

required
ref_grady ndarray

The derivative of the eta coordinate in the reference element.

required
xi ndarray

The xi coordinate.

required
eta ndarray

The eta coordinate.

required

Returns:

Type Description
ndarray

np.ndarray: The derivatives of the original coordinates [x, y] with respect to the reference coordinates.

Source code in scirex/core/sciml/fe/quad_bilinear.py
def get_orig_from_ref_derivative(
    self,
    ref_gradx: np.ndarray,
    ref_grady: np.ndarray,
    xi: np.ndarray,
    eta: np.ndarray,
) -> np.ndarray:
    """
    This method returns the derivatives of the original coordinates with respect to the reference coordinates.

    Args:
        ref_gradx (np.ndarray): The derivative of the xi coordinate in the reference element.
        ref_grady (np.ndarray): The derivative of the eta coordinate in the reference element.
        xi (np.ndarray): The xi coordinate.
        eta (np.ndarray): The eta coordinate.

    Returns:
        np.ndarray: The derivatives of the original coordinates [x, y] with respect to the reference coordinates.

    """
    n_test = ref_gradx.shape[0]
    gradx_orig = np.zeros(ref_gradx.shape, dtype=np.float64)
    grady_orig = np.zeros(ref_grady.shape, dtype=np.float64)

    for j in range(n_test):
        Xi = xi
        Eta = eta
        rec_detjk = 1 / (
            (self.xc1 + self.xc3 * Eta) * (self.yc2 + self.yc3 * Xi)
            - (self.xc2 + self.xc3 * Xi) * (self.yc1 + self.yc3 * Eta)
        )
        gradx_orig[j] = (
            (self.yc2 + self.yc3 * Xi) * ref_gradx[j]
            - (self.yc1 + self.yc3 * Eta) * ref_grady[j]
        ) * rec_detjk
        grady_orig[j] = (
            -(self.xc2 + self.xc3 * Xi) * ref_gradx[j]
            + (self.xc1 + self.xc3 * Eta) * ref_grady[j]
        ) * rec_detjk

    return gradx_orig, grady_orig

get_orig_from_ref_second_derivative(grad_xx_ref, grad_xy_ref, grad_yy_ref, xi, eta)

This method returns the second derivatives of the original coordinates with respect to the reference coordinates.

Parameters:

Name Type Description Default
grad_xx_ref ndarray

The second derivative of the xi coordinate in the reference element.

required
grad_xy_ref ndarray

The second derivative of the xi and eta coordinates in the reference element.

required
grad_yy_ref ndarray

The second derivative of the eta coordinate in the reference element.

required
xi ndarray

The xi coordinate.

required
eta ndarray

The eta coordinate.

required
Note

Second derivative calculations are not fully implemented in this method. Needs further development.

Source code in scirex/core/sciml/fe/quad_bilinear.py
def get_orig_from_ref_second_derivative(
    self,
    grad_xx_ref: np.ndarray,
    grad_xy_ref: np.ndarray,
    grad_yy_ref: np.ndarray,
    xi: np.ndarray,
    eta: np.ndarray,
):
    """
    This method returns the second derivatives of the original coordinates with respect to the reference coordinates.

    Args:
        grad_xx_ref (np.ndarray): The second derivative of the xi coordinate in the reference element.
        grad_xy_ref (np.ndarray): The second derivative of the xi and eta coordinates in the reference element.
        grad_yy_ref (np.ndarray): The second derivative of the eta coordinate in the reference element.
        xi (np.ndarray): The xi coordinate.
        eta (np.ndarray): The eta coordinate.

    Note:
        Second derivative calculations are not fully implemented in this method. Needs further development.
    """
    # print(" Error : Second Derivative not implemented -- Ignore this error, if second derivative is not required ")
    return grad_xx_ref, grad_xy_ref, grad_yy_ref

get_original_from_ref(xi, eta)

This method returns the original coordinates from the reference coordinates.

Parameters:

Name Type Description Default
xi ndarray

The xi coordinate.

required
eta ndarray

The eta coordinate.

required

Returns:

Type Description
ndarray

np.ndarray: Returns the transformed original coordinates from the reference coordinates.

Source code in scirex/core/sciml/fe/quad_bilinear.py
def get_original_from_ref(self, xi: np.ndarray, eta: np.ndarray) -> np.ndarray:
    """
    This method returns the original coordinates from the reference coordinates.

    Args:
        xi (np.ndarray): The xi coordinate.
        eta (np.ndarray): The eta coordinate.

    Returns:
        np.ndarray: Returns the transformed original coordinates from the reference coordinates.
    """
    x = self.xc0 + self.xc1 * xi + self.xc2 * eta + self.xc3 * xi * eta
    y = self.yc0 + self.yc1 * xi + self.yc2 * eta + self.yc3 * xi * eta

    return np.array([x, y], dtype=np.float64)

set_cell()

Set the cell coordinates, which will be used as intermediate values to calculate the Jacobian and actual values.

Returns:

Type Description

None

Source code in scirex/core/sciml/fe/quad_bilinear.py
def set_cell(self):
    """
    Set the cell coordinates, which will be used as intermediate values to calculate the Jacobian and actual values.

    Args:
        None

    Returns:
        None
    """
    self.x0 = self.co_ordinates[0][0]
    self.x1 = self.co_ordinates[1][0]
    self.x2 = self.co_ordinates[2][0]
    self.x3 = self.co_ordinates[3][0]

    # get the y-coordinates of the cell
    self.y0 = self.co_ordinates[0][1]
    self.y1 = self.co_ordinates[1][1]
    self.y2 = self.co_ordinates[2][1]
    self.y3 = self.co_ordinates[3][1]

    self.xc0 = (self.x0 + self.x1 + self.x2 + self.x3) * 0.25
    self.xc1 = (-self.x0 + self.x1 + self.x2 - self.x3) * 0.25
    self.xc2 = (-self.x0 - self.x1 + self.x2 + self.x3) * 0.25
    self.xc3 = (self.x0 - self.x1 + self.x2 - self.x3) * 0.25

    self.yc0 = (self.y0 + self.y1 + self.y2 + self.y3) * 0.25
    self.yc1 = (-self.y0 + self.y1 + self.y2 - self.y3) * 0.25
    self.yc2 = (-self.y0 - self.y1 + self.y2 + self.y3) * 0.25
    self.yc3 = (self.y0 - self.y1 + self.y2 - self.y3) * 0.25