Skip to content

jacobi

Module: basis_2d_QN_Jacobi.py

This module implements a specialized basis function class for 2D QN finite elements using Jacobi polynomials. It provides functionality for computing basis functions and their derivatives in two dimensions, primarily used in variational physics-informed neural networks (VPINNs) with domain decomposition.

Classes: Basis2DQNJacobi: Main class implementing 2D basis functions using Jacobi polynomials

Dependencies: - numpy: For numerical computations - scipy.special: For Jacobi polynomial calculations - .basis_function_2d: For base class implementation

Key Features: - Implementation of 2D QN element basis functions using Jacobi polynomials - Computation of function values and derivatives up to second order - Tensor product construction of 2D basis functions from 1D components - Specialized handling of Jacobi polynomials for test functions - Support for variable number of shape functions

Authors

Thivin Anandh (http://thivinanandh.github.io/)

Version Info: 27/Dec/2024: Initial version: Thivin Anandh D

References: - hp-VPINNs: Variational Physics-Informed Neural Networks With Domain Decomposition: https://github.com/ehsankharazmi/hp-VPINNs/

Basis2DQNJacobi

Bases: BasisFunction2D

A specialized implementation of two-dimensional basis functions using Jacobi polynomials for QN elements.

This class provides a complete implementation for computing basis functions and their derivatives in two dimensions, specifically designed for use in variational physics-informed neural networks (VPINNs) with domain decomposition. The basis functions are constructed using Jacobi polynomials with efficient derivative computations.

The class inherits from BasisFunction2D and implements all required methods for computing function values and first/second order derivatives. The implementation follows the methodology described in hp-VPINNs research.

Attributes:

Name Type Description
num_shape_functions int

Total number of shape functions in the 2D element. Must be a perfect square as it represents tensor product of 1D functions.

Methods:

Name Description
jacobi_wrapper

Evaluates Jacobi polynomial at given points

djacobi

Computes kth derivative of Jacobi polynomial

test_fcnx

Computes x-component test functions

test_fcny

Computes y-component test functions

dtest_fcn

Computes first derivatives of test functions

ddtest_fcn

Computes second derivatives of test functions

value

Computes values of all basis functions

gradx

Computes x-derivatives of all basis functions

grady

Computes y-derivatives of all basis functions

gradxx

Computes second x-derivatives of all basis functions

gradyy

Computes second y-derivatives of all basis functions

gradxy

Computes mixed xy-derivatives of all basis functions

Implementation Details
  • Basis functions are constructed as tensor products of 1D test functions
  • Test functions are derived from Jacobi polynomials with parameters (0,0)
  • All computations maintain double precision (float64)
  • Efficient vectorized operations using numpy arrays
Example
basis = Basis2DQNJacobi(num_shape_functions=16)  # Creates 4x4 basis functions
xi = np.linspace(-1, 1, 100)
eta = np.linspace(-1, 1, 100)
values = basis.value(xi, eta)
x_derivatives = basis.gradx(xi, eta)
Source code in scirex/core/sciml/fe/basis_2d_qn_jacobi.py
class Basis2DQNJacobi(BasisFunction2D):
    """
    A specialized implementation of two-dimensional basis functions using Jacobi polynomials for QN elements.

    This class provides a complete implementation for computing basis functions and their derivatives
    in two dimensions, specifically designed for use in variational physics-informed neural networks
    (VPINNs) with domain decomposition. The basis functions are constructed using Jacobi polynomials
    with efficient derivative computations.

    The class inherits from BasisFunction2D and implements all required methods for computing
    function values and first/second order derivatives. The implementation follows the methodology
    described in hp-VPINNs research.

    Attributes:
        num_shape_functions (int): Total number of shape functions in the 2D element.
            Must be a perfect square as it represents tensor product of 1D functions.

    Methods:
        jacobi_wrapper(n, a, b, x): Evaluates Jacobi polynomial at given points
        djacobi(n, a, b, x, k): Computes kth derivative of Jacobi polynomial
        test_fcnx(n_test, x): Computes x-component test functions
        test_fcny(n_test, y): Computes y-component test functions
        dtest_fcn(n_test, x): Computes first derivatives of test functions
        ddtest_fcn(n_test, x): Computes second derivatives of test functions
        value(xi, eta): Computes values of all basis functions
        gradx(xi, eta): Computes x-derivatives of all basis functions
        grady(xi, eta): Computes y-derivatives of all basis functions
        gradxx(xi, eta): Computes second x-derivatives of all basis functions
        gradyy(xi, eta): Computes second y-derivatives of all basis functions
        gradxy(xi, eta): Computes mixed xy-derivatives of all basis functions

    Implementation Details:
        - Basis functions are constructed as tensor products of 1D test functions
        - Test functions are derived from Jacobi polynomials with parameters (0,0)
        - All computations maintain double precision (float64)
        - Efficient vectorized operations using numpy arrays

    Example:
        ```python
        basis = Basis2DQNJacobi(num_shape_functions=16)  # Creates 4x4 basis functions
        xi = np.linspace(-1, 1, 100)
        eta = np.linspace(-1, 1, 100)
        values = basis.value(xi, eta)
        x_derivatives = basis.gradx(xi, eta)
        ```
    """

    def __init__(self, num_shape_functions: int):
        super().__init__(num_shape_functions)

    def jacobi_wrapper(self, n: int, a: int, b: int, x: np.ndarray) -> np.ndarray:
        """
        Evaluate the Jacobi polynomial of degree `n` with parameters `a` and `b` at the given points `x`.

        Args:
            n (int): Degree of the Jacobi polynomial.
            a (float): First parameter of the Jacobi polynomial.
            b (float): Second parameter of the Jacobi polynomial.
            x (np.ndarray): Points at which to evaluate the Jacobi polynomial.

        Returns:
            np.ndarray: Values of the Jacobi polynomial at the given points `x`.
        """

        x = np.array(x, dtype=np.float64)
        return jacobi(n, a, b)(x)

    # Derivative of the Jacobi polynomials
    def djacobi(self, n: int, a: int, b: int, x: np.ndarray, k: int) -> np.ndarray:
        """
        Evaluate the k-th derivative of the Jacobi polynomial of degree n with parameters a and b at the given points x.

        Args:
            n (int): Degree of the Jacobi polynomial.
            a (float): First parameter of the Jacobi polynomial.
            b (float): Second parameter of the Jacobi polynomial.
            x (np.ndarray): Points at which to evaluate the Jacobi polynomial.
            k (int): Order of the derivative.

        Returns:
            np.ndarray: Values of the k-th derivative of the Jacobi polynomial at the given points x.

        Raises:
            ValueError: If the derivative order is not 1 or 2
        """
        x = np.array(x, dtype=np.float64)
        if k == 1:
            return jacobi(n, a, b).deriv()(x)
        if k == 2:
            return jacobi(n, a, b).deriv(2)(x)
        else:
            print(f"Invalid derivative order {k} in {__name__}.")
            raise ValueError("Derivative order should be 1 or 2.")

    ## Helper Function
    def test_fcnx(self, n_test: int, x: np.ndarray) -> np.ndarray:
        """
        Compute the x-component of the test functions for a given number of test functions and x-coordinates.

        Args:
            n_test (int): Number of test functions.
            x (np.ndarray): x-coordinates at which to evaluate the test functions.

        Returns:
            np.ndarray: Values of the x-component of the test functions.

        """
        test_total = []
        for n in range(1, n_test + 1):
            test = self.jacobi_wrapper(n - 1, 0, 0, x)
            test_total.append(test)
        return np.asarray(test_total, np.float64)

    def test_fcny(self, n_test: int, y: np.ndarray) -> np.ndarray:
        """
        Compute the y-component of the test functions for a given number of test functions and y-coordinates.

        Args:
            n_test (int): Number of test functions.
            y (np.ndarray): y-coordinates at which to evaluate the test functions.

        Returns:
            np.ndarray: Values of the y-component of the test functions.

        """
        test_total = []
        for n in range(1, n_test + 1):
            test = self.jacobi_wrapper(n - 1, 0, 0, y)
            test_total.append(test)
        return np.asarray(test_total, np.float64)

    def dtest_fcn(self, n_test: int, x: np.ndarray) -> np.ndarray:
        """
        Compute the x-derivatives of the test functions for a given number of test functions and x-coordinates.

        Args:
            n_test (int): Number of test functions.
            x (np.ndarray): x-coordinates at which to evaluate the test functions.

        Returns:
            np.ndarray: Values of the x-derivatives of the test functions.
        """
        d1test_total = []
        for n in range(1, n_test + 1):
            d1test = self.djacobi(n - 1, 0, 0, x, 1)
            d1test_total.append(d1test)
        return np.asarray(d1test_total)

    def ddtest_fcn(self, n_test: int, x: np.ndarray) -> np.ndarray:
        """
        Compute the x-derivatives of the test functions for a given number of test functions and x-coordinates.

        Args:
            n_test (int): Number of test functions.
            x (np.ndarray): x-coordinates at which to evaluate the test functions.

        Returns:
            np.ndarray: Values of the x-derivatives of the test functions.
        """
        d1test_total = []
        for n in range(1, n_test + 1):
            d1test = self.djacobi(n - 1, 0, 0, x, 2)
            d1test_total.append(d1test)
        return np.asarray(d1test_total)

    def value(self, xi: np.ndarray, eta: np.ndarray) -> np.ndarray:
        """
        This method returns the values of the basis functions at the given (xi, eta) coordinates.

        Args:
            xi (np.ndarray): x-coordinates at which to evaluate the basis functions.
            eta (np.ndarray): y-coordinates at which to evaluate the basis functions.

        Returns:
            np.ndarray: Values of the basis functions.
        """
        num_shape_func_in_1d = int(np.sqrt(self.num_shape_functions))
        test_x = self.test_fcnx(num_shape_func_in_1d, xi)
        test_y = self.test_fcny(num_shape_func_in_1d, eta)
        values = np.zeros((self.num_shape_functions, len(xi)), dtype=np.float64)

        for i in range(num_shape_func_in_1d):
            values[num_shape_func_in_1d * i : num_shape_func_in_1d * (i + 1), :] = (
                test_x[i, :] * test_y
            )

        return values

    def gradx(self, xi: np.ndarray, eta: np.ndarray) -> np.ndarray:
        """
        This method returns the x-derivatives of the basis functions at the given (xi, eta) coordinates.

        Args:
            xi (np.ndarray): x-coordinates at which to evaluate the basis functions.
            eta (np.ndarray): y-coordinates at which to evaluate the basis functions.

        Returns:
            np.ndarray: Values of the x-derivatives of the basis functions.
        """
        num_shape_func_in_1d = int(np.sqrt(self.num_shape_functions))
        grad_test_x = self.dtest_fcn(num_shape_func_in_1d, xi)
        test_y = self.test_fcny(num_shape_func_in_1d, eta)
        values = np.zeros((self.num_shape_functions, len(xi)), dtype=np.float64)

        for i in range(num_shape_func_in_1d):
            values[num_shape_func_in_1d * i : num_shape_func_in_1d * (i + 1), :] = (
                grad_test_x[i, :] * test_y
            )

        return values

    def grady(self, xi: np.ndarray, eta: np.ndarray) -> np.ndarray:
        """
        This method returns the y-derivatives of the basis functions at the given (xi, eta) coordinates.

        Args:
            xi (np.ndarray): x-coordinates at which to evaluate the basis functions.
            eta (np.ndarray): y-coordinates at which to evaluate the basis functions.

        Returns:
            np.ndarray: Values of the y-derivatives of the basis functions.
        """
        num_shape_func_in_1d = int(np.sqrt(self.num_shape_functions))
        test_x = self.test_fcnx(num_shape_func_in_1d, xi)
        grad_test_y = self.dtest_fcn(num_shape_func_in_1d, eta)
        values = np.zeros((self.num_shape_functions, len(xi)), dtype=np.float64)

        for i in range(num_shape_func_in_1d):
            values[num_shape_func_in_1d * i : num_shape_func_in_1d * (i + 1), :] = (
                test_x[i, :] * grad_test_y
            )

        return values

    def gradxx(self, xi: np.ndarray, eta: np.ndarray) -> np.ndarray:
        """
        This method returns the xx-derivatives of the basis functions at the given (xi, eta) coordinates.

        Args:
            xi (np.ndarray): x-coordinates at which to evaluate the basis functions.
            eta (np.ndarray): y-coordinates at which to evaluate the basis functions.

        Returns:
            np.ndarray: Values of the xx-derivatives of the basis functions.
        """
        num_shape_func_in_1d = int(np.sqrt(self.num_shape_functions))
        grad_grad_x = self.ddtest_fcn(num_shape_func_in_1d, xi)
        test_y = self.test_fcny(num_shape_func_in_1d, eta)
        values = np.zeros((self.num_shape_functions, len(xi)), dtype=np.float64)

        for i in range(num_shape_func_in_1d):
            values[num_shape_func_in_1d * i : num_shape_func_in_1d * (i + 1), :] = (
                grad_grad_x[i, :] * test_y
            )

        return values

    def gradxy(self, xi: np.ndarray, eta: np.ndarray) -> np.ndarray:
        """
        This method returns the xy-derivatives of the basis functions at the given (xi, eta) coordinates.

        Args:
            xi (np.ndarray): x-coordinates at which to evaluate the basis functions.
            eta (np.ndarray): y-coordinates at which to evaluate the basis functions.

        Returns:
            np.ndarray: Values of the xy-derivatives of the basis functions.
        """
        num_shape_func_in_1d = int(np.sqrt(self.num_shape_functions))
        grad_test_x = self.dtest_fcn(num_shape_func_in_1d, xi)
        grad_test_y = self.dtest_fcn(num_shape_func_in_1d, eta)
        values = np.zeros((self.num_shape_functions, len(xi)), dtype=np.float64)

        for i in range(num_shape_func_in_1d):
            values[num_shape_func_in_1d * i : num_shape_func_in_1d * (i + 1), :] = (
                grad_test_x[i, :] * grad_test_y
            )

        return values

    def gradyy(self, xi: np.ndarray, eta: np.ndarray) -> np.ndarray:
        """
        This method returns the yy-derivatives of the basis functions at the given (xi, eta) coordinates.

        Args:
            xi (np.ndarray): x-coordinates at which to evaluate the basis functions.
            eta (np.ndarray): y-coordinates at which to evaluate the basis functions.

        Returns:
            np.ndarray: Values of the yy-derivatives of the basis functions.
        """
        num_shape_func_in_1d = int(np.sqrt(self.num_shape_functions))
        test_x = self.test_fcnx(num_shape_func_in_1d, xi)
        grad_grad_y = self.ddtest_fcn(num_shape_func_in_1d, eta)
        values = np.zeros((self.num_shape_functions, len(xi)), dtype=np.float64)

        for i in range(num_shape_func_in_1d):
            values[num_shape_func_in_1d * i : num_shape_func_in_1d * (i + 1), :] = (
                test_x[i, :] * grad_grad_y
            )

        return values

ddtest_fcn(n_test, x)

Compute the x-derivatives of the test functions for a given number of test functions and x-coordinates.

Parameters:

Name Type Description Default
n_test int

Number of test functions.

required
x ndarray

x-coordinates at which to evaluate the test functions.

required

Returns:

Type Description
ndarray

np.ndarray: Values of the x-derivatives of the test functions.

Source code in scirex/core/sciml/fe/basis_2d_qn_jacobi.py
def ddtest_fcn(self, n_test: int, x: np.ndarray) -> np.ndarray:
    """
    Compute the x-derivatives of the test functions for a given number of test functions and x-coordinates.

    Args:
        n_test (int): Number of test functions.
        x (np.ndarray): x-coordinates at which to evaluate the test functions.

    Returns:
        np.ndarray: Values of the x-derivatives of the test functions.
    """
    d1test_total = []
    for n in range(1, n_test + 1):
        d1test = self.djacobi(n - 1, 0, 0, x, 2)
        d1test_total.append(d1test)
    return np.asarray(d1test_total)

djacobi(n, a, b, x, k)

Evaluate the k-th derivative of the Jacobi polynomial of degree n with parameters a and b at the given points x.

Parameters:

Name Type Description Default
n int

Degree of the Jacobi polynomial.

required
a float

First parameter of the Jacobi polynomial.

required
b float

Second parameter of the Jacobi polynomial.

required
x ndarray

Points at which to evaluate the Jacobi polynomial.

required
k int

Order of the derivative.

required

Returns:

Type Description
ndarray

np.ndarray: Values of the k-th derivative of the Jacobi polynomial at the given points x.

Raises:

Type Description
ValueError

If the derivative order is not 1 or 2

Source code in scirex/core/sciml/fe/basis_2d_qn_jacobi.py
def djacobi(self, n: int, a: int, b: int, x: np.ndarray, k: int) -> np.ndarray:
    """
    Evaluate the k-th derivative of the Jacobi polynomial of degree n with parameters a and b at the given points x.

    Args:
        n (int): Degree of the Jacobi polynomial.
        a (float): First parameter of the Jacobi polynomial.
        b (float): Second parameter of the Jacobi polynomial.
        x (np.ndarray): Points at which to evaluate the Jacobi polynomial.
        k (int): Order of the derivative.

    Returns:
        np.ndarray: Values of the k-th derivative of the Jacobi polynomial at the given points x.

    Raises:
        ValueError: If the derivative order is not 1 or 2
    """
    x = np.array(x, dtype=np.float64)
    if k == 1:
        return jacobi(n, a, b).deriv()(x)
    if k == 2:
        return jacobi(n, a, b).deriv(2)(x)
    else:
        print(f"Invalid derivative order {k} in {__name__}.")
        raise ValueError("Derivative order should be 1 or 2.")

dtest_fcn(n_test, x)

Compute the x-derivatives of the test functions for a given number of test functions and x-coordinates.

Parameters:

Name Type Description Default
n_test int

Number of test functions.

required
x ndarray

x-coordinates at which to evaluate the test functions.

required

Returns:

Type Description
ndarray

np.ndarray: Values of the x-derivatives of the test functions.

Source code in scirex/core/sciml/fe/basis_2d_qn_jacobi.py
def dtest_fcn(self, n_test: int, x: np.ndarray) -> np.ndarray:
    """
    Compute the x-derivatives of the test functions for a given number of test functions and x-coordinates.

    Args:
        n_test (int): Number of test functions.
        x (np.ndarray): x-coordinates at which to evaluate the test functions.

    Returns:
        np.ndarray: Values of the x-derivatives of the test functions.
    """
    d1test_total = []
    for n in range(1, n_test + 1):
        d1test = self.djacobi(n - 1, 0, 0, x, 1)
        d1test_total.append(d1test)
    return np.asarray(d1test_total)

gradx(xi, eta)

This method returns the x-derivatives of the basis functions at the given (xi, eta) coordinates.

Parameters:

Name Type Description Default
xi ndarray

x-coordinates at which to evaluate the basis functions.

required
eta ndarray

y-coordinates at which to evaluate the basis functions.

required

Returns:

Type Description
ndarray

np.ndarray: Values of the x-derivatives of the basis functions.

Source code in scirex/core/sciml/fe/basis_2d_qn_jacobi.py
def gradx(self, xi: np.ndarray, eta: np.ndarray) -> np.ndarray:
    """
    This method returns the x-derivatives of the basis functions at the given (xi, eta) coordinates.

    Args:
        xi (np.ndarray): x-coordinates at which to evaluate the basis functions.
        eta (np.ndarray): y-coordinates at which to evaluate the basis functions.

    Returns:
        np.ndarray: Values of the x-derivatives of the basis functions.
    """
    num_shape_func_in_1d = int(np.sqrt(self.num_shape_functions))
    grad_test_x = self.dtest_fcn(num_shape_func_in_1d, xi)
    test_y = self.test_fcny(num_shape_func_in_1d, eta)
    values = np.zeros((self.num_shape_functions, len(xi)), dtype=np.float64)

    for i in range(num_shape_func_in_1d):
        values[num_shape_func_in_1d * i : num_shape_func_in_1d * (i + 1), :] = (
            grad_test_x[i, :] * test_y
        )

    return values

gradxx(xi, eta)

This method returns the xx-derivatives of the basis functions at the given (xi, eta) coordinates.

Parameters:

Name Type Description Default
xi ndarray

x-coordinates at which to evaluate the basis functions.

required
eta ndarray

y-coordinates at which to evaluate the basis functions.

required

Returns:

Type Description
ndarray

np.ndarray: Values of the xx-derivatives of the basis functions.

Source code in scirex/core/sciml/fe/basis_2d_qn_jacobi.py
def gradxx(self, xi: np.ndarray, eta: np.ndarray) -> np.ndarray:
    """
    This method returns the xx-derivatives of the basis functions at the given (xi, eta) coordinates.

    Args:
        xi (np.ndarray): x-coordinates at which to evaluate the basis functions.
        eta (np.ndarray): y-coordinates at which to evaluate the basis functions.

    Returns:
        np.ndarray: Values of the xx-derivatives of the basis functions.
    """
    num_shape_func_in_1d = int(np.sqrt(self.num_shape_functions))
    grad_grad_x = self.ddtest_fcn(num_shape_func_in_1d, xi)
    test_y = self.test_fcny(num_shape_func_in_1d, eta)
    values = np.zeros((self.num_shape_functions, len(xi)), dtype=np.float64)

    for i in range(num_shape_func_in_1d):
        values[num_shape_func_in_1d * i : num_shape_func_in_1d * (i + 1), :] = (
            grad_grad_x[i, :] * test_y
        )

    return values

gradxy(xi, eta)

This method returns the xy-derivatives of the basis functions at the given (xi, eta) coordinates.

Parameters:

Name Type Description Default
xi ndarray

x-coordinates at which to evaluate the basis functions.

required
eta ndarray

y-coordinates at which to evaluate the basis functions.

required

Returns:

Type Description
ndarray

np.ndarray: Values of the xy-derivatives of the basis functions.

Source code in scirex/core/sciml/fe/basis_2d_qn_jacobi.py
def gradxy(self, xi: np.ndarray, eta: np.ndarray) -> np.ndarray:
    """
    This method returns the xy-derivatives of the basis functions at the given (xi, eta) coordinates.

    Args:
        xi (np.ndarray): x-coordinates at which to evaluate the basis functions.
        eta (np.ndarray): y-coordinates at which to evaluate the basis functions.

    Returns:
        np.ndarray: Values of the xy-derivatives of the basis functions.
    """
    num_shape_func_in_1d = int(np.sqrt(self.num_shape_functions))
    grad_test_x = self.dtest_fcn(num_shape_func_in_1d, xi)
    grad_test_y = self.dtest_fcn(num_shape_func_in_1d, eta)
    values = np.zeros((self.num_shape_functions, len(xi)), dtype=np.float64)

    for i in range(num_shape_func_in_1d):
        values[num_shape_func_in_1d * i : num_shape_func_in_1d * (i + 1), :] = (
            grad_test_x[i, :] * grad_test_y
        )

    return values

grady(xi, eta)

This method returns the y-derivatives of the basis functions at the given (xi, eta) coordinates.

Parameters:

Name Type Description Default
xi ndarray

x-coordinates at which to evaluate the basis functions.

required
eta ndarray

y-coordinates at which to evaluate the basis functions.

required

Returns:

Type Description
ndarray

np.ndarray: Values of the y-derivatives of the basis functions.

Source code in scirex/core/sciml/fe/basis_2d_qn_jacobi.py
def grady(self, xi: np.ndarray, eta: np.ndarray) -> np.ndarray:
    """
    This method returns the y-derivatives of the basis functions at the given (xi, eta) coordinates.

    Args:
        xi (np.ndarray): x-coordinates at which to evaluate the basis functions.
        eta (np.ndarray): y-coordinates at which to evaluate the basis functions.

    Returns:
        np.ndarray: Values of the y-derivatives of the basis functions.
    """
    num_shape_func_in_1d = int(np.sqrt(self.num_shape_functions))
    test_x = self.test_fcnx(num_shape_func_in_1d, xi)
    grad_test_y = self.dtest_fcn(num_shape_func_in_1d, eta)
    values = np.zeros((self.num_shape_functions, len(xi)), dtype=np.float64)

    for i in range(num_shape_func_in_1d):
        values[num_shape_func_in_1d * i : num_shape_func_in_1d * (i + 1), :] = (
            test_x[i, :] * grad_test_y
        )

    return values

gradyy(xi, eta)

This method returns the yy-derivatives of the basis functions at the given (xi, eta) coordinates.

Parameters:

Name Type Description Default
xi ndarray

x-coordinates at which to evaluate the basis functions.

required
eta ndarray

y-coordinates at which to evaluate the basis functions.

required

Returns:

Type Description
ndarray

np.ndarray: Values of the yy-derivatives of the basis functions.

Source code in scirex/core/sciml/fe/basis_2d_qn_jacobi.py
def gradyy(self, xi: np.ndarray, eta: np.ndarray) -> np.ndarray:
    """
    This method returns the yy-derivatives of the basis functions at the given (xi, eta) coordinates.

    Args:
        xi (np.ndarray): x-coordinates at which to evaluate the basis functions.
        eta (np.ndarray): y-coordinates at which to evaluate the basis functions.

    Returns:
        np.ndarray: Values of the yy-derivatives of the basis functions.
    """
    num_shape_func_in_1d = int(np.sqrt(self.num_shape_functions))
    test_x = self.test_fcnx(num_shape_func_in_1d, xi)
    grad_grad_y = self.ddtest_fcn(num_shape_func_in_1d, eta)
    values = np.zeros((self.num_shape_functions, len(xi)), dtype=np.float64)

    for i in range(num_shape_func_in_1d):
        values[num_shape_func_in_1d * i : num_shape_func_in_1d * (i + 1), :] = (
            test_x[i, :] * grad_grad_y
        )

    return values

jacobi_wrapper(n, a, b, x)

Evaluate the Jacobi polynomial of degree n with parameters a and b at the given points x.

Parameters:

Name Type Description Default
n int

Degree of the Jacobi polynomial.

required
a float

First parameter of the Jacobi polynomial.

required
b float

Second parameter of the Jacobi polynomial.

required
x ndarray

Points at which to evaluate the Jacobi polynomial.

required

Returns:

Type Description
ndarray

np.ndarray: Values of the Jacobi polynomial at the given points x.

Source code in scirex/core/sciml/fe/basis_2d_qn_jacobi.py
def jacobi_wrapper(self, n: int, a: int, b: int, x: np.ndarray) -> np.ndarray:
    """
    Evaluate the Jacobi polynomial of degree `n` with parameters `a` and `b` at the given points `x`.

    Args:
        n (int): Degree of the Jacobi polynomial.
        a (float): First parameter of the Jacobi polynomial.
        b (float): Second parameter of the Jacobi polynomial.
        x (np.ndarray): Points at which to evaluate the Jacobi polynomial.

    Returns:
        np.ndarray: Values of the Jacobi polynomial at the given points `x`.
    """

    x = np.array(x, dtype=np.float64)
    return jacobi(n, a, b)(x)

test_fcnx(n_test, x)

Compute the x-component of the test functions for a given number of test functions and x-coordinates.

Parameters:

Name Type Description Default
n_test int

Number of test functions.

required
x ndarray

x-coordinates at which to evaluate the test functions.

required

Returns:

Type Description
ndarray

np.ndarray: Values of the x-component of the test functions.

Source code in scirex/core/sciml/fe/basis_2d_qn_jacobi.py
def test_fcnx(self, n_test: int, x: np.ndarray) -> np.ndarray:
    """
    Compute the x-component of the test functions for a given number of test functions and x-coordinates.

    Args:
        n_test (int): Number of test functions.
        x (np.ndarray): x-coordinates at which to evaluate the test functions.

    Returns:
        np.ndarray: Values of the x-component of the test functions.

    """
    test_total = []
    for n in range(1, n_test + 1):
        test = self.jacobi_wrapper(n - 1, 0, 0, x)
        test_total.append(test)
    return np.asarray(test_total, np.float64)

test_fcny(n_test, y)

Compute the y-component of the test functions for a given number of test functions and y-coordinates.

Parameters:

Name Type Description Default
n_test int

Number of test functions.

required
y ndarray

y-coordinates at which to evaluate the test functions.

required

Returns:

Type Description
ndarray

np.ndarray: Values of the y-component of the test functions.

Source code in scirex/core/sciml/fe/basis_2d_qn_jacobi.py
def test_fcny(self, n_test: int, y: np.ndarray) -> np.ndarray:
    """
    Compute the y-component of the test functions for a given number of test functions and y-coordinates.

    Args:
        n_test (int): Number of test functions.
        y (np.ndarray): y-coordinates at which to evaluate the test functions.

    Returns:
        np.ndarray: Values of the y-component of the test functions.

    """
    test_total = []
    for n in range(1, n_test + 1):
        test = self.jacobi_wrapper(n - 1, 0, 0, y)
        test_total.append(test)
    return np.asarray(test_total, np.float64)

value(xi, eta)

This method returns the values of the basis functions at the given (xi, eta) coordinates.

Parameters:

Name Type Description Default
xi ndarray

x-coordinates at which to evaluate the basis functions.

required
eta ndarray

y-coordinates at which to evaluate the basis functions.

required

Returns:

Type Description
ndarray

np.ndarray: Values of the basis functions.

Source code in scirex/core/sciml/fe/basis_2d_qn_jacobi.py
def value(self, xi: np.ndarray, eta: np.ndarray) -> np.ndarray:
    """
    This method returns the values of the basis functions at the given (xi, eta) coordinates.

    Args:
        xi (np.ndarray): x-coordinates at which to evaluate the basis functions.
        eta (np.ndarray): y-coordinates at which to evaluate the basis functions.

    Returns:
        np.ndarray: Values of the basis functions.
    """
    num_shape_func_in_1d = int(np.sqrt(self.num_shape_functions))
    test_x = self.test_fcnx(num_shape_func_in_1d, xi)
    test_y = self.test_fcny(num_shape_func_in_1d, eta)
    values = np.zeros((self.num_shape_functions, len(xi)), dtype=np.float64)

    for i in range(num_shape_func_in_1d):
        values[num_shape_func_in_1d * i : num_shape_func_in_1d * (i + 1), :] = (
            test_x[i, :] * test_y
        )

    return values