For computing direct and inverse Legendre transforms.
Legendre module variables
| Type | Visibility | Attributes | Name | Initial | |||
|---|---|---|---|---|---|---|---|
| logical, | private | :: | mod_legendre_initialized | = | .false. | ||
| real(kind=8), | private, | allocatable, dimension(:, :) | :: | epsi | |||
| real(kind=8), | private, | allocatable, dimension(:, :, :) | :: | cpol | |||
| real(kind=8), | private, | allocatable, dimension(:, :) | :: | repsi | |||
| integer, | private, | allocatable, dimension(:) | :: | nsh2 | |||
| real(kind=8), | private, | allocatable, dimension(:) | :: | wt | |||
| class(ModGeometry_t), | private, | pointer | :: | mod_geometry | => | null() | Spectral module instance |
| logical, | private | :: | mod_geometry_initialized | = | .false. |
| procedure, public :: initialize => ModLegendre_initialize | |
| procedure, public :: legendre => ModLegendre_legendre | |
| procedure, public :: legendre_inv => ModLegendre_legendre_inv | |
| procedure, public :: delete => ModLegendre_delete | |
| procedure, public :: legendre_poly => ModLegendre_legendre_poly |
Computes inverse Legendre transformation. The Legendre polynomials (cpol) and the triangular shape definition (nsh2) needs to be initialized and passed to the function.
| Type | Intent | Optional | Attributes | Name | ||
|---|---|---|---|---|---|---|
| class(ModLegendre_t), | intent(in) | :: | this | |||
| real(kind=p), | intent(in) | :: | input(2*mx,nx) | Input field |
Output field
Computes direct Legendre transformation. The Legendre polynomials (cpol), the triangular shape definition (nsh2), and the gaussian weights (wt) used for the integration of the Legendre transform needs to be initialized and passed to the function.
| Type | Intent | Optional | Attributes | Name | ||
|---|---|---|---|---|---|---|
| class(ModLegendre_t), | intent(in) | :: | this | |||
| real(kind=p), | intent(in) | :: | input(2*mx,il) | Input field |
Output field
Compute Gaussian weights for direct Legendre transform
Weights in gaussian quadrature (sum should equal 1.0)
Compute associated Legendre polynomials at given latitude.
| Type | Intent | Optional | Attributes | Name | ||
|---|---|---|---|---|---|---|
| class(ModLegendre_t), | intent(in) | :: | this | |||
| integer, | intent(in) | :: | j | The latitude index to compute the polynomials at |
The Legendre polynomials
Initializes Legendre transforms and constants used for other subroutines that manipulate spherical harmonics. The Legendre polynomials Epsilon function used for various spectral calculations 1/epsi Used for defining shape of spectral triangle Gaussian weights used for integration in direct Legendre transform
| Type | Intent | Optional | Attributes | Name | ||
|---|---|---|---|---|---|---|
| class(ModLegendre_t), | intent(inout) | :: | this | |||
| class(ModGeometry_t), | intent(in), | target | :: | mod_geometry |
| Type | Intent | Optional | Attributes | Name | ||
|---|---|---|---|---|---|---|
| class(ModLegendre_t), | intent(inout) | :: | this |