alps_nhds Module

Module including the NHDS implementation for bi-Maxwellian/cold-plasma reference cases. The original NHDS code can be found under github.com/danielver02/NHDS.


Used by

  • module~~alps_nhds~~UsedByGraph module~alps_nhds alps_nhds proc~disp disp proc~disp->module~alps_nhds

Functions

private function besselI(n, x)

Calculates the modified Bessel function of argument x and order n.

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: n

Order of the modified Bessel function.

double precision, intent(in) :: x

Argument of the modified Bessel function.

Return Value doubleprecision

private function dispfunct(zeta, kpos)

Calculates the dispersion function based on the complex error function.

Arguments

Type IntentOptional Attributes Name
double complex, intent(in) :: zeta

Argument of dispersion function.

logical, intent(in) :: kpos

Check whether kpar is positive.

Return Value doublecomplex

private function BESSI(N, X)

Function to calculate the first kind modified Bessel function of integer order N for any real X.

Arguments

Type IntentOptional Attributes Name
integer :: N
real(kind=8) :: X

Return Value doubleprecision

private function BESSI0(X)

Auxiliary Bessel functions for N=0, N=1

Arguments

Type IntentOptional Attributes Name
real(kind=8) :: X

Return Value doubleprecision

private function BESSI1(X)

Modified Bessel function of order 1.

Arguments

Type IntentOptional Attributes Name
real(kind=8) :: X

Return Value doubleprecision


Subroutines

public subroutine calc_chi(chi, j, kz, kperp, x)

Subroutine that calculates the susceptibility of species j based on NHDS.

Arguments

Type IntentOptional Attributes Name
double complex, intent(out) :: chi(3,3)

Susceptibility tensor of species j.

integer, intent(in) :: j

Index for species.

double precision, intent(in) :: kz

Normalised parallel wavenumber.

double precision, intent(in) :: kperp

Normalised perpendicular wavenumber.

double complex, intent(in) :: x

Normalised complex frequency.

private subroutine calc_ypsilon(Y, j, n, kz, kperp, x)

Calculates the Y-tensor according to Stix for a bi-Maxwelling, using the NHDS calculation.

Arguments

Type IntentOptional Attributes Name
double complex, intent(out) :: Y(3,3)

Y-tensor as defined by Stix.

integer, intent(in) :: j

Index for species.

integer, intent(in) :: n

Index of sum over Bessel functions.

double precision, intent(in) :: kz

Normalised parallel wavenumber.

double precision, intent(in) :: kperp

Normalised perpendicular wavenumber.

double complex, intent(in) :: x

Normalised complex frequency.

private subroutine calc_chi_cold(chi, j, kz, kperp, x)

Subroutine that calculates the susceptibility of species j based on the cold-plasma dispersion relation based on the paper Verscharen & Chandran, ApJ 764, 88, 2013.

Arguments

Type IntentOptional Attributes Name
double complex, intent(out) :: chi(3,3)

Susceptibility tensor of species j.

integer, intent(in) :: j

Index for species.

double precision, intent(in) :: kz

Normalised parallel wavenumber.

double precision, intent(in) :: kperp

Normalised perpendicular wavenumber.

double complex, intent(in) :: x

Normalised complex frequency.

private subroutine WOFZ(XI, YI, U, V, FLAG)

Given a complex number Z = (XI,YI), this subroutine computes the value of the Faddeeva function W(Z) = exp(-Z2)erfc(-IZ), where erfc is the complex complementary error function and I is the imaginary unit.

Arguments

Type IntentOptional Attributes Name
real :: XI
real :: YI
real :: U
real :: V
logical :: FLAG