This module contains the key numerical functions of ALPS.
This function returns the determinant of the dispersion tensor for a given frequency om.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
double complex, | intent(in) | :: | om |
Complex wave frequency . |
This function returns the full integral expression according to Eq. (2.9) in the code paper.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
double complex, | intent(in) | :: | om |
Complex wave frequency . |
||
integer, | intent(in) | :: | nn |
Order of the Bessel function. |
||
integer, | intent(in) | :: | mode |
Index of the entries in the T-tensor of Eq. (2.10). |
||
logical, | intent(in) | :: | found_res |
Check whether a resonance is found. |
This function performs the integral in Eq. (2.9) of the code paper, but without accounting for the Landau contour integral. It is called by full_integrate.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
double complex, | intent(in) | :: | om |
Complex wave frequency . |
||
integer, | intent(in) | :: | nn |
Order of the Bessel function. |
||
integer, | intent(in) | :: | mode |
Index of the entries in the T-tensor of Eq. (2.10). |
||
integer | :: | iparmin |
Minimum limit index of parallel momentum for integration. |
|||
integer | :: | iparmax |
Maximum limit index of parallel momentum for integration. |
This function performs the integration near resonances as described in Section 3.1 of the code paper. It is only called if resonances are present in or near the integration domain.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
double complex, | intent(in) | :: | om |
Complex wave frequency . |
||
integer, | intent(in) | :: | nn |
Order of the Bessel function. |
||
integer, | intent(in) | :: | mode |
Index of the entries in the T-tensor of Eq. (2.10). |
This function returns the function from Eq. (3.2) of the code paper.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
double precision, | intent(in) | :: | ppar_real |
Real part of the momentum at which is evaluated. |
||
integer, | intent(in) | :: | iperp |
Index of the perpendicular momentum. |
||
double complex, | intent(in) | :: | om |
Complex wave frequency . |
||
integer, | intent(in) | :: | nn |
Order of the Bessel function. |
||
integer, | intent(in) | :: | mode |
Index of the entries in the T-tensor of Eq. (2.10). |
This function evaluates the Landau contour according to Eqs. (3.8) and (3.9) of the code paper.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
double complex, | intent(in) | :: | om |
Complex wave frequency . |
||
integer, | intent(in) | :: | nn |
Order of the Bessel function. |
||
integer, | intent(in) | :: | mode |
Index of the entries in the T-tensor of Eq. (2.10). |
This function returns the ee term in Eq. (2.9).
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
double complex, | intent(in) | :: | om |
Complex wave frequency . |
This function evaluates the term proportional to in Eq. (2.9) of the code paper.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
double complex, | intent(in) | :: | om |
Complex wave frequency . |
||
integer, | intent(in) | :: | nn |
Order of the Bessel function. |
||
integer, | intent(in) | :: | iperp |
Index to loop over perpendicular momentum. |
||
integer, | intent(in) | :: | ipar |
Index to loop over parallel momentum. |
This function returns the T-tensor according to Eq. (2.10) of the code paper.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
integer, | intent(in) | :: | nn |
Order of the Bessel function. |
||
integer, | intent(in) | :: | iperp |
Index to loop over perpendicular momentum. |
||
integer, | intent(in) | :: | ipar |
Index to loop over parallel momentum. |
||
integer, | intent(in) | :: | mode |
Index of the entries in the T-tensor of Eq. (2.10). |
This function returns the T-tensor according to Eq. (2.10) of the code paper for the case in which it is evaluated at the complex resonance momentum.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
integer, | intent(in) | :: | nn |
Order of the Bessel function. |
||
integer, | intent(in) | :: | iperp |
Index to loop over perpendicular momentum. |
||
double complex, | intent(in) | :: | p_res |
Complex resonance momentum. |
||
integer, | intent(in) | :: | mode |
Index of the entries in the T-tensor of Eq. (2.10). |
This subroutine calculates the perpendicular and parallel derivatives of the background velocity distribution function f0.
This subroutine determines whether any kinetic resonances are located in the integration domain.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
double complex, | intent(in) | :: | om |
Complex wave frequency . |
||
integer, | intent(in) | :: | nn |
Order of Bessel function. |
||
logical, | intent(out) | :: | found_res_plus |
Check whether a resonance is found at positive n. |
||
logical, | intent(out) | :: | found_res_minus |
Check whether a resonance is found at negative n. |
This subroutine applies the secant method to find the roots of the dispersion tensor.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
double complex, | intent(inout) | :: | om |
Complex wave frequency . |
||
integer, | intent(in) | :: | in |
Root number |
This subroutine scans solutions along a single prescribed path in wavevector space. KGK: This line causes the solution to (occasionally) smoothly transition to unphysical values. Suppressing until we understand the error.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
integer, | intent(in) | :: | ik |
Index of scan number. |
This subroutine calculates the relative electric and magnetic field amplitudes, the relative fluctuations in the density and velocity of all species, and the heating rates of the given solution. It is based on the calc_eigen routine by Greg Howes and Kris Klein.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
double complex, | intent(in) | :: | omega |
Complex wave frequency . |
||
double complex, | intent(out), | dimension(1:3) | :: | electric |
Relative electric field amplitude (eigenfunction). |
|
double complex, | intent(out), | dimension(1:3) | :: | magnetic |
Relative magnetic field amplitude (eigenfunction). |
|
double complex, | intent(out), | dimension(1:3,1:nspec) | :: | vmean |
Relative velocity-fluctuation amplitude (eigenfunction). |
|
double complex, | intent(out), | dimension(1:nspec) | :: | ds |
Relative density-fluctuation amplitude (eigenfunction). |
|
double precision, | intent(out), | dimension(1:nspec) | :: | Ps |
Relative heating rate of a given species. |
|
logical, | intent(in) | :: | eigen_L |
Check whether eigenfunction calculation is requested. |
||
logical, | intent(in) | :: | heat_L |
Check whether eigenfunction calculation is requested. |
This subroutine scans along a prescribed plane in wavevector space to map out in this space. It is required that n_scan=2. KGK: This line causes the solution to (occasionally) smoothly transition to unphysical values. Suppressing until we understand the error.
This subroutine calculates the map of the determinant of the dispersion tensor in complex frequency space. check
This subroutine refines the guess at the starting point of the search for solutions to the dispersion relation when scanning. It is also used by map_search to identify the roots on the map.
This subroutine identifies the minima of the coarse map grid. It is called by map_search. The code is based on a routine by Greg Howes, 2006.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
double precision, | intent(in), | dimension(:,:), pointer | :: | val |
Array of determinant of the dispersion tensor. |
|
integer, | intent(in) | :: | numroots |
Number of roots. |
||
integer, | intent(out), | dimension(1:2,1:numroots) | :: | iroots |
Indices of roots. |
|
integer, | intent(out) | :: | nroots |
Number of roots found. |
This subroutine determines the maximum required order of the Bessel functions in Eq. (2.9) of the code paper.
This subroutine defines the tasks for the individual processes. It uses the number of species and the required orders of the Bessel functions to define the splitting across the MPI processes.
This subroutine determines the array of Bessel functions that is used in the T-tensor of Eq. (2.10) of the code paper.