Spherical Derivative Operators¶
SphericalDerivative1D
¶
Bases: Module
1D spectral derivative operator on Gauss-Legendre latitude grid.
Supports colatitude (theta) gradient and Laplacian using Legendre polynomial spectral coefficients.
Mathematical Formulation:¶
Expand u in normalised Legendre polynomials: u(theta) = sum_l c_l * P_l_norm(cos(theta))
Colatitude gradient (using recurrence P_l'(cos(theta)) = -P_l^1(cos(theta))): du/d_theta = -sum_l c_l * sqrt(l*(l+1)) * P_l^1_norm(cos(theta))
Laplacian on the unit sphere (eigenvalue relation): nabla^2_sphere u = -sum_l c_l * l*(l+1) / R^2 * P_l_norm(cos(theta))
Attributes:¶
grid : SphericalGrid1D The 1D Gauss-Legendre grid.
Source code in spectraldiffx/_src/spherical/operators.py
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 | |
Functions¶
to_spectral(u)
¶
Forward Discrete Legendre Transform.
c_l = sum_j w_j * P_l_norm(cos(theta_j)) * u(theta_j)
Parameters:¶
u : Float[Array, "N"] Physical field at GL nodes.
Returns:¶
c : Float[Array, "N"] Legendre spectral coefficients.
Source code in spectraldiffx/_src/spherical/operators.py
from_spectral(c)
¶
Inverse Discrete Legendre Transform.
u(theta_j) ≈ sum_l c_l * P_l_norm(cos(theta_j))
Parameters:¶
c : Float[Array, "N"] Legendre spectral coefficients.
Returns:¶
u : Float[Array, "N"] Physical field at GL nodes.
Source code in spectraldiffx/_src/spherical/operators.py
gradient(u, spectral=False)
¶
Colatitude gradient: du/d_theta.
Algorithm
- Forward DLT: c_l = P_norm_matrix @ (w * u)
- Multiply: c_grad_l = -sqrt(l*(l+1)) * c_l
- Reconstruct: du/d_theta = P1_norm_matrix.T @ c_grad
Parameters:¶
u : Float[Array, "N"] Physical field or spectral coefficients (if spectral=True). spectral : bool If True, u contains Legendre coefficients.
Returns:¶
du_dtheta : Float[Array, "N"] Colatitude derivative du/d_theta at GL nodes.
Source code in spectraldiffx/_src/spherical/operators.py
laplacian(u, spectral=False)
¶
Spherical Laplacian: (1/sin(theta)) * d/d_theta [sin(theta) * du/d_theta] (zonal, m=0 case): nabla^2_sphere u = -l*(l+1)/R^2 * u in spectral space.
Parameters:¶
u : Float[Array, "N"] Physical field or spectral coefficients. spectral : bool If True, u contains Legendre coefficients.
Returns:¶
lap_u : Float[Array, "N"] Laplacian at GL nodes.
Source code in spectraldiffx/_src/spherical/operators.py
__call__(u, order=1, spectral=False)
¶
Apply derivative operator.
Parameters:¶
u : Float[Array, "N"] Physical field or spectral coefficients. order : int 1 = gradient (du/d_theta), 2 = Laplacian. spectral : bool If True, u contains Legendre coefficients.
Returns:¶
Array [N]
Source code in spectraldiffx/_src/spherical/operators.py
SphericalDerivative2D
¶
Bases: Module
2D pseudo-spectral derivative operators on the sphere.
Computes physical-space differential operators (gradient, divergence, curl, Laplacian) for a field u(theta, phi) on the sphere using the full Spherical Harmonic Transform.
Mathematical Formulation:¶
For a field u(theta, phi) on a sphere of radius R = Ly / pi:
Gradient (covariant components):
grad_theta u = (1/R) * du/d_theta
grad_phi u = 1 / (R * sin(theta)) * du/d_phi
Divergence of vector (V_theta, V_phi) (in physical coordinates):
div V = 1/(R*sin(theta)) * [d(V_theta*sin(theta))/d_theta + d(V_phi)/d_phi]
Scalar curl (vertical component of curl):
curl V = 1/(R*sin(theta)) * [d(V_theta)/d_phi - d(V_phi*sin(theta))/d_theta]
Laplacian (via eigenvalue in spectral space):
nabla^2 u = sum_{l,m} -l*(l+1)/R^2 * u_hat(l,m) * Y_l^m
Attributes:¶
grid : SphericalGrid2D The full lat-lon grid. deriv_theta : SphericalDerivative1D 1D operator reused for latitude derivatives.
Source code in spectraldiffx/_src/spherical/operators.py
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 | |
Functions¶
gradient(u, spectral=False)
¶
Compute the gradient of u on the sphere.
grad_theta u = (1/R) * du/d_theta
grad_phi u = 1 / (R * sin(theta)) * du/d_phi
Parameters:¶
u : Float[Array, "Ny Nx"] Physical field. spectral : bool If True, u is already in spectral (SHT) space.
Returns:¶
(grad_theta, grad_phi) : tuple of Float[Array, "Ny Nx"] Covariant gradient components.
Source code in spectraldiffx/_src/spherical/operators.py
divergence(v_theta, v_phi, spectral=False)
¶
Divergence of a vector field (V_theta, V_phi) on the sphere.
div V = 1/(R*sin(theta)) * [d(V_theta*sin(theta))/d_theta + dV_phi/d_phi]
Parameters:¶
v_theta : Float[Array, "Ny Nx"] Colatitude component (physical space). v_phi : Float[Array, "Ny Nx"] Longitude component (physical space). spectral : bool Unused (reserved for API consistency; inputs expected in physical space).
Returns:¶
div : Float[Array, "Ny Nx"] Scalar divergence field.
Source code in spectraldiffx/_src/spherical/operators.py
curl(v_theta, v_phi, spectral=False)
¶
Vertical curl (scalar vorticity) of a 2D vector field on the sphere.
curl V = 1/(R*sin(theta)) * [dV_theta/d_phi - d(V_phi*sin(theta))/d_theta]
Parameters:¶
v_theta : Float[Array, "Ny Nx"] Colatitude component. v_phi : Float[Array, "Ny Nx"] Longitude component. spectral : bool Unused (reserved for API consistency).
Returns:¶
zeta : Float[Array, "Ny Nx"] Scalar vorticity field.
Source code in spectraldiffx/_src/spherical/operators.py
laplacian(u, spectral=False)
¶
Scalar Laplace-Beltrami operator on the sphere.
In spectral space
nabla^2 u_hat(l, m) = -l*(l+1)/R^2 * u_hat(l, m)
Parameters:¶
u : Float[Array, "Ny Nx"] Physical field or spectral coefficients. spectral : bool If True, u is already in spectral space.
Returns:¶
lap_u : Float[Array, "Ny Nx"] Laplacian in physical space.
Source code in spectraldiffx/_src/spherical/operators.py
advection_scalar(v_theta, v_phi, q)
¶
Pseudo-spectral scalar advection: (V · nabla) q.
Computes the gradient of q spectrally, then multiplies by V in physical space: adv = V_theta * grad_theta(q) + V_phi * grad_phi(q)
Parameters:¶
v_theta : Float[Array, "Ny Nx"] Colatitude velocity component. v_phi : Float[Array, "Ny Nx"] Longitude velocity component. q : Float[Array, "Ny Nx"] Scalar field to advect.
Returns:¶
adv : Float[Array, "Ny Nx"] Advection term at each grid point.