Skip to content

Chebyshev

Module: basis_2d_qn_chebyshev_2.py

This module implements a specialized basis function class for 2D finite elements using Chebyshev 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:

Name Description
Basis2DQNChebyshev2

Main class implementing 2D basis functions using Chebyshev polynomials

Dependencies
  • numpy: For numerical computations and array operations
  • scipy.special: For Jacobi polynomial calculations and evaluations
  • .basis_function_2d: For base class BasisFunction2D implementation
Key Features
  • Implementation of 2D element basis functions using Chebyshev 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 through initialization
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/

Basis2DQNChebyshev2

Bases: BasisFunction2D

A specialized implementation of two-dimensional basis functions using Chebyshev polynomials for Q1 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 Chebyshev polynomials through Jacobi polynomial representations.

The class inherits from BasisFunction2D and implements all required methods for computing function values and derivatives. The implementation follows the methodology described in hp-VPINNs research by Ehsan Kharazmi et al.

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
value

Computes values of all basis functions at given points

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 normalized Jacobi polynomials
  • Special cases are handled for first few polynomial degrees in derivatives
  • All computations maintain double precision (float64)
  • Efficient vectorized operations using numpy arrays
Example
basis = Basis2DQNChebyshev2(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)
Notes
  • num_shape_functions must be a perfect square
  • All coordinate inputs (xi, eta) should be in the range [-1, 1]
  • Implementation optimized for vectorized operations on multiple points
  • Based on hp-VPINNs methodology: https://github.com/ehsankharazmi/hp-VPINNs/
References

Kharazmi, E., et al. "hp-VPINNs: Variational Physics-Informed Neural Networks With Domain Decomposition"

Source code in scirex/core/sciml/fe/basis_2d_qn_chebyshev_2.py
 64
 65
 66
 67
 68
 69
 70
 71
 72
 73
 74
 75
 76
 77
 78
 79
 80
 81
 82
 83
 84
 85
 86
 87
 88
 89
 90
 91
 92
 93
 94
 95
 96
 97
 98
 99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
class Basis2DQNChebyshev2(BasisFunction2D):
    """A specialized implementation of two-dimensional basis functions using Chebyshev polynomials for Q1 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 Chebyshev polynomials
    through Jacobi polynomial representations.

    The class inherits from BasisFunction2D and implements all required methods for computing
    function values and derivatives. The implementation follows the methodology described in
    hp-VPINNs research by Ehsan Kharazmi et al.

    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:
        value(xi, eta): Computes values of all basis functions at given points
        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 normalized Jacobi polynomials
        - Special cases are handled for first few polynomial degrees in derivatives
        - All computations maintain double precision (float64)
        - Efficient vectorized operations using numpy arrays

    Example:
        ```python
        basis = Basis2DQNChebyshev2(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)
        ```

    Notes:
        - num_shape_functions must be a perfect square
        - All coordinate inputs (xi, eta) should be in the range [-1, 1]
        - Implementation optimized for vectorized operations on multiple points
        - Based on hp-VPINNs methodology: https://github.com/ehsankharazmi/hp-VPINNs/

    References:
        Kharazmi, E., et al. "hp-VPINNs: Variational Physics-Informed Neural Networks
        With Domain Decomposition"
    """

    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:
        """Evaluates Jacobi polynomial at specified points.

        Computes values of nth degree Jacobi polynomial with parameters (a,b)
        at given points x.

        Args:
            n: Degree of Jacobi polynomial. Must be non-negative integer.
            a: First parameter of Jacobi polynomial
            b: Second parameter of Jacobi polynomial
            x: Points at which to evaluate polynomial
                Shape: (n_points,)

        Returns:
            np.ndarray: Values of Jacobi polynomial at input points
                Shape: Same as input x

        Notes:
            Wrapper around scipy.special.jacobi that ensures float64 precision
            and proper array handling.
        """
        x = np.array(x, dtype=np.float64)
        return jacobi(n, a, b)(x)

    ## Helper Function
    def test_fcnx(self, n_test: int, x: np.ndarray) -> np.ndarray:
        """Computes x-component test functions.

        Evaluates the x-direction test functions constructed as differences
        of normalized Jacobi polynomials.

        Args:
            n_test: Number of test functions to compute
            x: Points at which to evaluate functions
                Shape: (n_points,)

        Returns:
            np.ndarray: Values of test functions at input points
                Shape: (n_test, n_points)

        Notes:
            Test functions are constructed as differences of normalized Jacobi
            polynomials following hp-VPINNs methodology.
        """
        test_total = []
        for n in range(1, n_test + 1):
            test = self.jacobi_wrapper(n + 1, -1 / 2, -1 / 2, x) / self.jacobi_wrapper(
                n + 1, -1 / 2, -1 / 2, 1
            ) - self.jacobi_wrapper(n - 1, -1 / 2, -1 / 2, x) / self.jacobi_wrapper(
                n - 1, -1 / 2, -1 / 2, 1
            )
            test_total.append(test)
        return np.asarray(test_total, np.float64)

    def test_fcny(self, n_test: int, y: np.ndarray) -> np.ndarray:
        """Computes y-component test functions.

        Evaluates the y-direction test functions constructed as differences
        of normalized Jacobi polynomials.

        Args:
            n_test: Number of test functions to compute
            y: Points at which to evaluate functions
                Shape: (n_points,)

        Returns:
            np.ndarray: Values of test functions at input points
                Shape: (n_test, n_points)

        Notes:
            Test functions are constructed as differences of normalized Jacobi
            polynomials following hp-VPINNs methodology.
        """
        test_total = []
        for n in range(1, n_test + 1):
            test = self.jacobi_wrapper(n + 1, -1 / 2, -1 / 2, y) / self.jacobi_wrapper(
                n + 1, -1 / 2, -1 / 2, 1
            ) - self.jacobi_wrapper(n - 1, -1 / 2, -1 / 2, y) / self.jacobi_wrapper(
                n - 1, -1 / 2, -1 / 2, 1
            )
            test_total.append(test)
        return np.asarray(test_total, np.float64)

    def dtest_fcn(self, n_test: int, x: np.ndarray) -> np.ndarray:
        """Computes first and second derivatives of test functions.

        Calculates derivatives of test functions constructed from Jacobi
        polynomials, handling special cases for n=1,2 separately.

        Args:
            n_test: Number of test functions
            x: Points at which to evaluate derivatives
                Shape: (n_points,)

        Returns:
            tuple(np.ndarray, np.ndarray): First and second derivatives
                First element: First derivatives, shape (n_test, n_points)
                Second element: Second derivatives, shape (n_test, n_points)

        Notes:
            Special cases for n=1,2 ensure proper derivative calculations
            following hp-VPINNs methodology.
        """
        d1test_total = []
        d2test_total = []
        for n in range(1, n_test + 1):
            if n == 1:
                d1test = (
                    ((n + 1) / 2)
                    * self.jacobi_wrapper(n, 1 / 2, 1 / 2, x)
                    / self.jacobi_wrapper(n + 1, -1 / 2, -1 / 2, 1)
                )
                d2test = (
                    ((n + 2) * (n + 1) / (2 * 2))
                    * self.jacobi_wrapper(n - 1, 3 / 2, 3 / 2, x)
                    / self.jacobi_wrapper(n + 1, -1 / 2, -1 / 2, 1)
                )
                d1test_total.append(d1test)
                d2test_total.append(d2test)
            elif n == 2:
                d1test = ((n + 1) / 2) * self.jacobi_wrapper(
                    n, 1 / 2, 1 / 2, x
                ) / self.jacobi_wrapper(n + 1, -1 / 2, -1 / 2, 1) - (
                    (n - 1) / 2
                ) * self.jacobi_wrapper(
                    n - 2, 1 / 2, 1 / 2, x
                ) / self.jacobi_wrapper(
                    n - 1, -1 / 2, -1 / 2, 1
                )
                d2test = (
                    ((n + 2) * (n + 1) / (2 * 2))
                    * self.jacobi_wrapper(n - 1, 3 / 2, 3 / 2, x)
                    / self.jacobi_wrapper(n + 1, -1 / 2, -1 / 2, 1)
                )
                d1test_total.append(d1test)
                d2test_total.append(d2test)
            else:
                d1test = ((n + 1) / 2) * self.jacobi_wrapper(
                    n, 1 / 2, 1 / 2, x
                ) / self.jacobi_wrapper(n + 1, -1 / 2, -1 / 2, 1) - (
                    (n - 1) / 2
                ) * self.jacobi_wrapper(
                    n - 2, 1 / 2, 1 / 2, x
                ) / self.jacobi_wrapper(
                    n - 1, -1 / 2, -1 / 2, 1
                )
                d2test = ((n + 2) * (n + 1) / (2 * 2)) * self.jacobi_wrapper(
                    n - 1, 3 / 2, 3 / 2, x
                ) / self.jacobi_wrapper(n + 1, -1 / 2, -1 / 2, 1) - (
                    (n) * (n - 1) / (2 * 2)
                ) * self.jacobi_wrapper(
                    n - 3, 3 / 2, 3 / 2, x
                ) / self.jacobi_wrapper(
                    n - 1, -1 / 2, -1 / 2, 1
                )
                d1test_total.append(d1test)
                d2test_total.append(d2test)
        return np.asarray(d1test_total), np.asarray(d2test_total)

    def value(self, xi: np.ndarray, eta: np.ndarray) -> np.ndarray:
        """Evaluates basis functions at given coordinates.

        Computes values of all basis functions at specified (xi,eta) points
        using tensor product of 1D test functions.

        Args:
            xi: x-coordinates at which to evaluate functions
                Shape: (n_points,)
            eta: y-coordinates at which to evaluate functions
                Shape: (n_points,)

        Returns:
            np.ndarray: Values of all basis functions
                Shape: (num_shape_functions, n_points)

        Notes:
            Basis functions are constructed as products of 1D test functions
            in x and y directions.
        """
        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:
        """Computes x-derivatives of basis functions.

        Evaluates partial derivatives with respect to x of all basis
        functions at given coordinates.

        Args:
            xi: x-coordinates at which to evaluate derivatives
                Shape: (n_points,)
            eta: y-coordinates at which to evaluate derivatives
                Shape: (n_points,)

        Returns:
            np.ndarray: Values of x-derivatives
                Shape: (num_shape_functions, n_points)

        Notes:
            Uses product rule with x-derivatives of test functions in
            x-direction and values in y-direction.
        """
        num_shape_func_in_1d = int(np.sqrt(self.num_shape_functions))
        grad_test_x = self.dtest_fcn(num_shape_func_in_1d, xi)[0]
        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:
        """Computes y-derivatives of basis functions.

        Evaluates partial derivatives with respect to y of all basis
        functions at given coordinates.

        Args:
            xi: x-coordinates at which to evaluate derivatives
                Shape: (n_points,)
            eta: y-coordinates at which to evaluate derivatives
                Shape: (n_points,)

        Returns:
            np.ndarray: Values of y-derivatives
                Shape: (num_shape_functions, n_points)

        Notes:
            Uses product rule with values in x-direction and y-derivatives
            of test functions in y-direction.
        """
        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)[0]
        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:
        """Computes second x-derivatives of basis functions.

        Evaluates second partial derivatives with respect to x of all basis
        functions at given coordinates.

        Args:
            xi: x-coordinates at which to evaluate derivatives
                Shape: (n_points,)
            eta: y-coordinates at which to evaluate derivatives
                Shape: (n_points,)

        Returns:
            np.ndarray: Values of second x-derivatives
                Shape: (num_shape_functions, n_points)

        Notes:
            Uses product rule with second x-derivatives of test functions in
            x-direction and values in y-direction.
        """
        num_shape_func_in_1d = int(np.sqrt(self.num_shape_functions))
        grad_grad_x = self.dtest_fcn(num_shape_func_in_1d, xi)[1]
        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:
        """Computes second x-derivatives of basis functions.

        Evaluates second partial derivatives with respect to x of all basis
        functions at given coordinates.

        Args:
            xi: x-coordinates at which to evaluate derivatives
                Shape: (n_points,)
            eta: y-coordinates at which to evaluate derivatives
                Shape: (n_points,)

        Returns:
            np.ndarray: Values of second x-derivatives
                Shape: (num_shape_functions, n_points)

        Notes:
            Uses product rule with second x-derivatives of test functions in
            x-direction and y derivative values in y-direction.
        """
        num_shape_func_in_1d = int(np.sqrt(self.num_shape_functions))
        grad_test_x = self.dtest_fcn(num_shape_func_in_1d, xi)[0]
        grad_test_y = self.dtest_fcn(num_shape_func_in_1d, eta)[0]
        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:
        """Computes second x-derivatives of basis functions.

        Evaluates second partial derivatives with respect to x of all basis
        functions at given coordinates.

        Args:
            xi: x-coordinates at which to evaluate derivatives
                Shape: (n_points,)
            eta: y-coordinates at which to evaluate derivatives
                Shape: (n_points,)

        Returns:
            np.ndarray: Values of second x-derivatives
                Shape: (num_shape_functions, n_points)

        Notes:
            Uses product rule with second y-derivatives of test functions in
            x-direction and values in y-direction.
        """
        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.dtest_fcn(num_shape_func_in_1d, eta)[1]
        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

dtest_fcn(n_test, x)

Computes first and second derivatives of test functions.

Calculates derivatives of test functions constructed from Jacobi polynomials, handling special cases for n=1,2 separately.

Parameters:

Name Type Description Default
n_test int

Number of test functions

required
x ndarray

Points at which to evaluate derivatives Shape: (n_points,)

required

Returns:

Name Type Description
tuple (ndarray, ndarray)

First and second derivatives First element: First derivatives, shape (n_test, n_points) Second element: Second derivatives, shape (n_test, n_points)

Notes

Special cases for n=1,2 ensure proper derivative calculations following hp-VPINNs methodology.

Source code in scirex/core/sciml/fe/basis_2d_qn_chebyshev_2.py
def dtest_fcn(self, n_test: int, x: np.ndarray) -> np.ndarray:
    """Computes first and second derivatives of test functions.

    Calculates derivatives of test functions constructed from Jacobi
    polynomials, handling special cases for n=1,2 separately.

    Args:
        n_test: Number of test functions
        x: Points at which to evaluate derivatives
            Shape: (n_points,)

    Returns:
        tuple(np.ndarray, np.ndarray): First and second derivatives
            First element: First derivatives, shape (n_test, n_points)
            Second element: Second derivatives, shape (n_test, n_points)

    Notes:
        Special cases for n=1,2 ensure proper derivative calculations
        following hp-VPINNs methodology.
    """
    d1test_total = []
    d2test_total = []
    for n in range(1, n_test + 1):
        if n == 1:
            d1test = (
                ((n + 1) / 2)
                * self.jacobi_wrapper(n, 1 / 2, 1 / 2, x)
                / self.jacobi_wrapper(n + 1, -1 / 2, -1 / 2, 1)
            )
            d2test = (
                ((n + 2) * (n + 1) / (2 * 2))
                * self.jacobi_wrapper(n - 1, 3 / 2, 3 / 2, x)
                / self.jacobi_wrapper(n + 1, -1 / 2, -1 / 2, 1)
            )
            d1test_total.append(d1test)
            d2test_total.append(d2test)
        elif n == 2:
            d1test = ((n + 1) / 2) * self.jacobi_wrapper(
                n, 1 / 2, 1 / 2, x
            ) / self.jacobi_wrapper(n + 1, -1 / 2, -1 / 2, 1) - (
                (n - 1) / 2
            ) * self.jacobi_wrapper(
                n - 2, 1 / 2, 1 / 2, x
            ) / self.jacobi_wrapper(
                n - 1, -1 / 2, -1 / 2, 1
            )
            d2test = (
                ((n + 2) * (n + 1) / (2 * 2))
                * self.jacobi_wrapper(n - 1, 3 / 2, 3 / 2, x)
                / self.jacobi_wrapper(n + 1, -1 / 2, -1 / 2, 1)
            )
            d1test_total.append(d1test)
            d2test_total.append(d2test)
        else:
            d1test = ((n + 1) / 2) * self.jacobi_wrapper(
                n, 1 / 2, 1 / 2, x
            ) / self.jacobi_wrapper(n + 1, -1 / 2, -1 / 2, 1) - (
                (n - 1) / 2
            ) * self.jacobi_wrapper(
                n - 2, 1 / 2, 1 / 2, x
            ) / self.jacobi_wrapper(
                n - 1, -1 / 2, -1 / 2, 1
            )
            d2test = ((n + 2) * (n + 1) / (2 * 2)) * self.jacobi_wrapper(
                n - 1, 3 / 2, 3 / 2, x
            ) / self.jacobi_wrapper(n + 1, -1 / 2, -1 / 2, 1) - (
                (n) * (n - 1) / (2 * 2)
            ) * self.jacobi_wrapper(
                n - 3, 3 / 2, 3 / 2, x
            ) / self.jacobi_wrapper(
                n - 1, -1 / 2, -1 / 2, 1
            )
            d1test_total.append(d1test)
            d2test_total.append(d2test)
    return np.asarray(d1test_total), np.asarray(d2test_total)

gradx(xi, eta)

Computes x-derivatives of basis functions.

Evaluates partial derivatives with respect to x of all basis functions at given coordinates.

Parameters:

Name Type Description Default
xi ndarray

x-coordinates at which to evaluate derivatives Shape: (n_points,)

required
eta ndarray

y-coordinates at which to evaluate derivatives Shape: (n_points,)

required

Returns:

Type Description
ndarray

np.ndarray: Values of x-derivatives Shape: (num_shape_functions, n_points)

Notes

Uses product rule with x-derivatives of test functions in x-direction and values in y-direction.

Source code in scirex/core/sciml/fe/basis_2d_qn_chebyshev_2.py
def gradx(self, xi: np.ndarray, eta: np.ndarray) -> np.ndarray:
    """Computes x-derivatives of basis functions.

    Evaluates partial derivatives with respect to x of all basis
    functions at given coordinates.

    Args:
        xi: x-coordinates at which to evaluate derivatives
            Shape: (n_points,)
        eta: y-coordinates at which to evaluate derivatives
            Shape: (n_points,)

    Returns:
        np.ndarray: Values of x-derivatives
            Shape: (num_shape_functions, n_points)

    Notes:
        Uses product rule with x-derivatives of test functions in
        x-direction and values in y-direction.
    """
    num_shape_func_in_1d = int(np.sqrt(self.num_shape_functions))
    grad_test_x = self.dtest_fcn(num_shape_func_in_1d, xi)[0]
    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)

Computes second x-derivatives of basis functions.

Evaluates second partial derivatives with respect to x of all basis functions at given coordinates.

Parameters:

Name Type Description Default
xi ndarray

x-coordinates at which to evaluate derivatives Shape: (n_points,)

required
eta ndarray

y-coordinates at which to evaluate derivatives Shape: (n_points,)

required

Returns:

Type Description
ndarray

np.ndarray: Values of second x-derivatives Shape: (num_shape_functions, n_points)

Notes

Uses product rule with second x-derivatives of test functions in x-direction and values in y-direction.

Source code in scirex/core/sciml/fe/basis_2d_qn_chebyshev_2.py
def gradxx(self, xi: np.ndarray, eta: np.ndarray) -> np.ndarray:
    """Computes second x-derivatives of basis functions.

    Evaluates second partial derivatives with respect to x of all basis
    functions at given coordinates.

    Args:
        xi: x-coordinates at which to evaluate derivatives
            Shape: (n_points,)
        eta: y-coordinates at which to evaluate derivatives
            Shape: (n_points,)

    Returns:
        np.ndarray: Values of second x-derivatives
            Shape: (num_shape_functions, n_points)

    Notes:
        Uses product rule with second x-derivatives of test functions in
        x-direction and values in y-direction.
    """
    num_shape_func_in_1d = int(np.sqrt(self.num_shape_functions))
    grad_grad_x = self.dtest_fcn(num_shape_func_in_1d, xi)[1]
    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)

Computes second x-derivatives of basis functions.

Evaluates second partial derivatives with respect to x of all basis functions at given coordinates.

Parameters:

Name Type Description Default
xi ndarray

x-coordinates at which to evaluate derivatives Shape: (n_points,)

required
eta ndarray

y-coordinates at which to evaluate derivatives Shape: (n_points,)

required

Returns:

Type Description
ndarray

np.ndarray: Values of second x-derivatives Shape: (num_shape_functions, n_points)

Notes

Uses product rule with second x-derivatives of test functions in x-direction and y derivative values in y-direction.

Source code in scirex/core/sciml/fe/basis_2d_qn_chebyshev_2.py
def gradxy(self, xi: np.ndarray, eta: np.ndarray) -> np.ndarray:
    """Computes second x-derivatives of basis functions.

    Evaluates second partial derivatives with respect to x of all basis
    functions at given coordinates.

    Args:
        xi: x-coordinates at which to evaluate derivatives
            Shape: (n_points,)
        eta: y-coordinates at which to evaluate derivatives
            Shape: (n_points,)

    Returns:
        np.ndarray: Values of second x-derivatives
            Shape: (num_shape_functions, n_points)

    Notes:
        Uses product rule with second x-derivatives of test functions in
        x-direction and y derivative values in y-direction.
    """
    num_shape_func_in_1d = int(np.sqrt(self.num_shape_functions))
    grad_test_x = self.dtest_fcn(num_shape_func_in_1d, xi)[0]
    grad_test_y = self.dtest_fcn(num_shape_func_in_1d, eta)[0]
    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)

Computes y-derivatives of basis functions.

Evaluates partial derivatives with respect to y of all basis functions at given coordinates.

Parameters:

Name Type Description Default
xi ndarray

x-coordinates at which to evaluate derivatives Shape: (n_points,)

required
eta ndarray

y-coordinates at which to evaluate derivatives Shape: (n_points,)

required

Returns:

Type Description
ndarray

np.ndarray: Values of y-derivatives Shape: (num_shape_functions, n_points)

Notes

Uses product rule with values in x-direction and y-derivatives of test functions in y-direction.

Source code in scirex/core/sciml/fe/basis_2d_qn_chebyshev_2.py
def grady(self, xi: np.ndarray, eta: np.ndarray) -> np.ndarray:
    """Computes y-derivatives of basis functions.

    Evaluates partial derivatives with respect to y of all basis
    functions at given coordinates.

    Args:
        xi: x-coordinates at which to evaluate derivatives
            Shape: (n_points,)
        eta: y-coordinates at which to evaluate derivatives
            Shape: (n_points,)

    Returns:
        np.ndarray: Values of y-derivatives
            Shape: (num_shape_functions, n_points)

    Notes:
        Uses product rule with values in x-direction and y-derivatives
        of test functions in y-direction.
    """
    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)[0]
    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)

Computes second x-derivatives of basis functions.

Evaluates second partial derivatives with respect to x of all basis functions at given coordinates.

Parameters:

Name Type Description Default
xi ndarray

x-coordinates at which to evaluate derivatives Shape: (n_points,)

required
eta ndarray

y-coordinates at which to evaluate derivatives Shape: (n_points,)

required

Returns:

Type Description
ndarray

np.ndarray: Values of second x-derivatives Shape: (num_shape_functions, n_points)

Notes

Uses product rule with second y-derivatives of test functions in x-direction and values in y-direction.

Source code in scirex/core/sciml/fe/basis_2d_qn_chebyshev_2.py
def gradyy(self, xi: np.ndarray, eta: np.ndarray) -> np.ndarray:
    """Computes second x-derivatives of basis functions.

    Evaluates second partial derivatives with respect to x of all basis
    functions at given coordinates.

    Args:
        xi: x-coordinates at which to evaluate derivatives
            Shape: (n_points,)
        eta: y-coordinates at which to evaluate derivatives
            Shape: (n_points,)

    Returns:
        np.ndarray: Values of second x-derivatives
            Shape: (num_shape_functions, n_points)

    Notes:
        Uses product rule with second y-derivatives of test functions in
        x-direction and values in y-direction.
    """
    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.dtest_fcn(num_shape_func_in_1d, eta)[1]
    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)

Evaluates Jacobi polynomial at specified points.

Computes values of nth degree Jacobi polynomial with parameters (a,b) at given points x.

Parameters:

Name Type Description Default
n int

Degree of Jacobi polynomial. Must be non-negative integer.

required
a int

First parameter of Jacobi polynomial

required
b int

Second parameter of Jacobi polynomial

required
x ndarray

Points at which to evaluate polynomial Shape: (n_points,)

required

Returns:

Type Description
ndarray

np.ndarray: Values of Jacobi polynomial at input points Shape: Same as input x

Notes

Wrapper around scipy.special.jacobi that ensures float64 precision and proper array handling.

Source code in scirex/core/sciml/fe/basis_2d_qn_chebyshev_2.py
def jacobi_wrapper(self, n: int, a: int, b: int, x: np.ndarray) -> np.ndarray:
    """Evaluates Jacobi polynomial at specified points.

    Computes values of nth degree Jacobi polynomial with parameters (a,b)
    at given points x.

    Args:
        n: Degree of Jacobi polynomial. Must be non-negative integer.
        a: First parameter of Jacobi polynomial
        b: Second parameter of Jacobi polynomial
        x: Points at which to evaluate polynomial
            Shape: (n_points,)

    Returns:
        np.ndarray: Values of Jacobi polynomial at input points
            Shape: Same as input x

    Notes:
        Wrapper around scipy.special.jacobi that ensures float64 precision
        and proper array handling.
    """
    x = np.array(x, dtype=np.float64)
    return jacobi(n, a, b)(x)

test_fcnx(n_test, x)

Computes x-component test functions.

Evaluates the x-direction test functions constructed as differences of normalized Jacobi polynomials.

Parameters:

Name Type Description Default
n_test int

Number of test functions to compute

required
x ndarray

Points at which to evaluate functions Shape: (n_points,)

required

Returns:

Type Description
ndarray

np.ndarray: Values of test functions at input points Shape: (n_test, n_points)

Notes

Test functions are constructed as differences of normalized Jacobi polynomials following hp-VPINNs methodology.

Source code in scirex/core/sciml/fe/basis_2d_qn_chebyshev_2.py
def test_fcnx(self, n_test: int, x: np.ndarray) -> np.ndarray:
    """Computes x-component test functions.

    Evaluates the x-direction test functions constructed as differences
    of normalized Jacobi polynomials.

    Args:
        n_test: Number of test functions to compute
        x: Points at which to evaluate functions
            Shape: (n_points,)

    Returns:
        np.ndarray: Values of test functions at input points
            Shape: (n_test, n_points)

    Notes:
        Test functions are constructed as differences of normalized Jacobi
        polynomials following hp-VPINNs methodology.
    """
    test_total = []
    for n in range(1, n_test + 1):
        test = self.jacobi_wrapper(n + 1, -1 / 2, -1 / 2, x) / self.jacobi_wrapper(
            n + 1, -1 / 2, -1 / 2, 1
        ) - self.jacobi_wrapper(n - 1, -1 / 2, -1 / 2, x) / self.jacobi_wrapper(
            n - 1, -1 / 2, -1 / 2, 1
        )
        test_total.append(test)
    return np.asarray(test_total, np.float64)

test_fcny(n_test, y)

Computes y-component test functions.

Evaluates the y-direction test functions constructed as differences of normalized Jacobi polynomials.

Parameters:

Name Type Description Default
n_test int

Number of test functions to compute

required
y ndarray

Points at which to evaluate functions Shape: (n_points,)

required

Returns:

Type Description
ndarray

np.ndarray: Values of test functions at input points Shape: (n_test, n_points)

Notes

Test functions are constructed as differences of normalized Jacobi polynomials following hp-VPINNs methodology.

Source code in scirex/core/sciml/fe/basis_2d_qn_chebyshev_2.py
def test_fcny(self, n_test: int, y: np.ndarray) -> np.ndarray:
    """Computes y-component test functions.

    Evaluates the y-direction test functions constructed as differences
    of normalized Jacobi polynomials.

    Args:
        n_test: Number of test functions to compute
        y: Points at which to evaluate functions
            Shape: (n_points,)

    Returns:
        np.ndarray: Values of test functions at input points
            Shape: (n_test, n_points)

    Notes:
        Test functions are constructed as differences of normalized Jacobi
        polynomials following hp-VPINNs methodology.
    """
    test_total = []
    for n in range(1, n_test + 1):
        test = self.jacobi_wrapper(n + 1, -1 / 2, -1 / 2, y) / self.jacobi_wrapper(
            n + 1, -1 / 2, -1 / 2, 1
        ) - self.jacobi_wrapper(n - 1, -1 / 2, -1 / 2, y) / self.jacobi_wrapper(
            n - 1, -1 / 2, -1 / 2, 1
        )
        test_total.append(test)
    return np.asarray(test_total, np.float64)

value(xi, eta)

Evaluates basis functions at given coordinates.

Computes values of all basis functions at specified (xi,eta) points using tensor product of 1D test functions.

Parameters:

Name Type Description Default
xi ndarray

x-coordinates at which to evaluate functions Shape: (n_points,)

required
eta ndarray

y-coordinates at which to evaluate functions Shape: (n_points,)

required

Returns:

Type Description
ndarray

np.ndarray: Values of all basis functions Shape: (num_shape_functions, n_points)

Notes

Basis functions are constructed as products of 1D test functions in x and y directions.

Source code in scirex/core/sciml/fe/basis_2d_qn_chebyshev_2.py
def value(self, xi: np.ndarray, eta: np.ndarray) -> np.ndarray:
    """Evaluates basis functions at given coordinates.

    Computes values of all basis functions at specified (xi,eta) points
    using tensor product of 1D test functions.

    Args:
        xi: x-coordinates at which to evaluate functions
            Shape: (n_points,)
        eta: y-coordinates at which to evaluate functions
            Shape: (n_points,)

    Returns:
        np.ndarray: Values of all basis functions
            Shape: (num_shape_functions, n_points)

    Notes:
        Basis functions are constructed as products of 1D test functions
        in x and y directions.
    """
    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