Legendre
Module: basis_2d_QN_Legendre.py
This module implements a specialized basis function class for 2D Quad elements using Legendre 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 |
---|---|
Basis2DQNLegendre |
Main class implementing 2D basis functions using Legendre polynomials |
Dependencies
- numpy: For numerical computations
- scipy.special: For Jacobi polynomial calculations
- .basis_function_2d: For base class implementation
Key Features
- Implementation of 2D Q1 element basis functions using Legendre 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
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/
Basis2DQNLegendre
Bases: BasisFunction2D
A specialized implementation of two-dimensional basis functions using Legendre 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 Legendre polynomials implemented through Jacobi polynomial representations with parameters (0,0).
The class inherits from BasisFunction2D and implements all required methods for computing function values and derivatives up to second order.
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 |
test_fcnx |
Computes x-component test functions |
test_fcny |
Computes y-component test functions |
dtest_fcn |
Computes first and 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 use Legendre polynomials via Jacobi polynomials with (0,0) parameters
- Special cases handled for n=1,2 in derivative calculations
- All computations maintain double precision (float64)
- Efficient vectorized operations using numpy arrays
Example
Source code in scirex/core/sciml/fe/basis_2d_qn_legendre.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 |
|
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_legendre.py
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_legendre.py
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_legendre.py
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_legendre.py
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_legendre.py
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_legendre.py
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
|
int
|
First parameter of the Jacobi polynomial. |
required |
b
|
int
|
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. |
Source code in scirex/core/sciml/fe/basis_2d_qn_legendre.py
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_legendre.py
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_legendre.py
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. |