V3FIT
Data Types | Functions/Subroutines | Variables
fourier Module Reference

Module contains subroutines for computing FFT on parallel radial intervals. Converts quantities between real and fourier space. Note fixarray must be called once before call any of the Fourier subroutines to calculate all necessary sine and cosine factors on a fixed mesh of collocation angles. More...

Data Types

interface  check
 Interface for checking results of the unit tests. More...
 
interface  fourier_class
 Base class containing fourier memory. More...
 

Functions/Subroutines

class(fourier_class) function, pointer fourier_construct (mpol, ntor, ntheta, nzeta, nfp, sym)
 Construct a fourier_class object.. More...
 
subroutine fourier_destruct (this)
 Deconstruct a fourier_class object. More...
 
subroutine fourier_tomnsp_3d (this, xuv, xmn, parity)
 Convert a quantity from real to fourier space in parallel. More...
 
subroutine fourier_tomnsp_3d_pest (this, xuv, xmn, lmns, lmnc,
 Convert a quantity from real to fourier space in pest coordinates. More...
 
pure subroutine fourier_tomnsp_2d (this, xuv, xmn, parity)
 Fourier transform a 2D quantity to fourier space. More...
 
pure subroutine fourier_tomnsp_2d_pest (this, xuv, xmn, lambda, dup
 Fourier transform a 2D quantity to fourier space. More...
 
pure subroutine fourier_tomnsp_2d_u (this, xuv)
 Sum over the poloidal part. More...
 
pure subroutine fourier_tomnsp_2d_u_pest (this, xuv, lambda, dupdu,
 Sum over the poloidal part. More...
 
pure subroutine fourier_tomnsp_2d_v (this, xmn, parity)
 Sum over toroidal angle. More...
 
subroutine fourier_toijsp_3d (this, xmn, xuv, dflag, parity)
 Convert a quantity from fourier to real space in parallel. More...
 
pure subroutine fourier_toijsp_2d (this, xmn, xuv, dflag, parity)
 Fourier transform a 2D quantity to real space. More...
 
subroutine fourier_get_origin (this, xuv, mode, asym)
 Computes the origin value of a half-mesh quantity. More...
 
subroutine fourier_round_trig (cosi, sini)
 Round trig values to whole number for near numbers. More...
 
logical function test_fourier ()
 Test fourier routines. More...
 
logical function, private check_mn (expected, received, m, n, name)
 Check test values. More...
 
logical function, private check_all (expected, received, mpol, ntor, name)
 Check all test values. More...
 

Variables

integer, parameter f_cos = 0
 Cosine parity.
 
integer, parameter f_sin = 1
 Sine parity.
 
integer, parameter f_none = 0
 Do not sum fouier real space transformed quantity. This is used when a stellarator symmetric parity is being converted to real space.
 
integer, parameter f_ds = 1
 Radial derivative flag.
 
integer, parameter f_du = 2
 Poloidal derivative flag.
 
integer, parameter f_dv = 4
 Toroidal derivative flag.
 
integer, parameter f_dudv = IOR(f_du, f_dv)
 Combined toroidal and poloidal derivatives.
 
integer, parameter f_sum = 8
 Sum fouier real space transformed quantity. This is used when a non-stellarator symmetric parity is being converted to real space.
 
integer, parameter b_ds = 0
 Bit position of the f_ds flag.
 
integer, parameter b_du = 1
 Bit position of the f_du flag.
 
integer, parameter b_dv = 2
 Bit position of the f_dv flag.
 
integer, parameter b_sum = 3
 Bit position of the f_sum flag.
 
integer, parameter m0 = 0
 m = 0 mode.
 
integer, parameter m1 = 1
 m = 1 mode.
 
integer, parameter m2 = 2
 m = 2 mode.
 
integer, parameter n0 = 0
 n = 0 mode.
 
integer, parameter n1 = 1
 n = 1 mode.
 

Detailed Description

Module contains subroutines for computing FFT on parallel radial intervals. Converts quantities between real and fourier space. Note fixarray must be called once before call any of the Fourier subroutines to calculate all necessary sine and cosine factors on a fixed mesh of collocation angles.

Function/Subroutine Documentation

◆ check_all()

logical function, private fourier::check_all ( real (rprec), dimension(:,:), intent(in)  expected,
real (rprec), dimension(:,:), intent(in)  received,
integer, intent(in)  mpol,
integer, intent(in)  ntor,
character (len=*), intent(in)  name 
)
private

Check all test values.

Parameters
[in]expectedExpected value for the test.
[in]receivedRecieved value for the test.
[in]nameName of the test.

Definition at line 1203 of file fourier.f90.

◆ check_mn()

logical function, private fourier::check_mn ( real (rprec), intent(in)  expected,
real (rprec), intent(in)  received,
integer, intent(in)  m,
integer, intent(in)  n,
character (len=*), intent(in)  name 
)
private

Check test values.

Parameters
[in]expectedExpected value for the test.
[in]receivedRecieved value for the test.
[in]mPoloidal mode to check.
[in]nToroidal mode to check.
[in]nameName of the test.

Definition at line 1170 of file fourier.f90.

◆ fourier_construct()

class (fourier_class) function, pointer fourier::fourier_construct ( integer, intent(in)  mpol,
integer, intent(in)  ntor,
integer, intent(in)  ntheta,
integer, intent(in)  nzeta,
integer, intent(in)  nfp,
logical, intent(in)  sym 
)

Construct a fourier_class object..

Allocates memory and initializes a fourier_class object. Computes the othronorm and cosine and sine buffers.

This subroutine computes the cosine-sine factors that will be needed when moving between Fourier and real space. All normalizations are contained in the poloidal quantities used in the Fourier to real space transformation:

  • sinmui
  • cosmui

Fourier representations are assumed to have stellarator symmetry: *# cosine (iparity=0): C(u, v) = sum_u sum_v (C_mn COS(mu + n*nfp*v)) *# sine (iparity=1): S(u, v) = sum_u sum_v (S_mn SIN(mu + n*nfp*v))

The number of collocation points have been set initalially equal to the number of modes: theta_j = j*pi/M, j = 0,...,M Where M = mpol - 1 zeta_k = k*2*pi/(2N + 1) k = 0,...,2N Where N = ntor.

Parameters
[in]mpolNumber of poloidal modes.
[in]ntorNumber of toroidal modes.
[in]nthetaNumber of poloidal real grid points.
[in]nzetaNumber of toroidal real grid points.
[in]nfpNumber of field periods.
[in]symSymmetry flag.

Definition at line 180 of file fourier.f90.

◆ fourier_destruct()

subroutine fourier::fourier_destruct ( type (fourier_class), intent(inout)  this)

Deconstruct a fourier_class object.

Deallocates memory and uninitializes a fourier_class object.

Parameters
[in,out]thisA fourier_class instance.

Definition at line 300 of file fourier.f90.

◆ fourier_get_origin()

subroutine fourier::fourier_get_origin ( class (fourier_class), intent(inout)  this,
real(dp), dimension(:,:,:), intent(inout), allocatable  xuv,
integer, intent(in)  mode,
logical, intent(in)  asym 
)

Computes the origin value of a half-mesh quantity.

The origin quantity is found from the m mode of the js = 2 value.

Parameters
[in,out]thisA fourier_class instance.
[in,out]xuvReal space quantity.
[in]modePoloidal mode number to be retained.
[in]asymAdd stellarator asymmetric terms.

Definition at line 949 of file fourier.f90.

◆ fourier_round_trig()

subroutine fourier::fourier_round_trig ( real (dp), intent(inout)  cosi,
real (dp), intent(inout)  sini 
)

Round trig values to whole number for near numbers.

Parameters
[in,out]cosiCosine value to round.
[in,out]siniSine value to round.

Definition at line 1011 of file fourier.f90.

◆ fourier_toijsp_2d()

pure subroutine fourier::fourier_toijsp_2d ( class (fourier_class), intent(inout)  this,
real (dp), dimension(:,:), intent(in)  xmn,
real (dp), dimension(:,:), intent(inout)  xuv,
integer, intent(in)  dflag,
integer, intent(in)  parity 
)

Fourier transform a 2D quantity to real space.

Works by reshaping to a 3D array then passing to the normal 3D version.

Parameters
[in]thisA fourier_class instance.
[in]xmnFourier space quantity.
[in,out]xuvReal space quantity.
[in]dflagDerivative and sum control flags.
[in]parityFourier parity flag.

Definition at line 826 of file fourier.f90.

◆ fourier_toijsp_3d()

subroutine fourier::fourier_toijsp_3d ( class (fourier_class), intent(inout)  this,
real (dp), dimension(:,:,:), intent(in)  xmn,
real (dp), dimension(:,:,:), intent(inout)  xuv,
integer, intent(in)  dflag,
integer, intent(in)  parity 
)

Convert a quantity from fourier to real space in parallel.

This subroutine moves a quantity to real space by summing over its Fourier harmonics. If the f_du flag is set take the poloidal derivative. If the f_dv flag is set take the toroidal derivative. If the f_sum flag is set add the real quanity to the previous value.

Parameters
[in]thisA fourier_class instance.
[in]xmnFourier space quantity.
[in,out]xuvReal space quantity.
[in]dflagDerivative and sum control flags.
[in]parityFourier parity flag.

Definition at line 783 of file fourier.f90.

◆ fourier_tomnsp_2d()

pure subroutine fourier::fourier_tomnsp_2d ( class (fourier_class), intent(inout)  this,
real (dp), dimension(:,:), intent(in)  xuv,
real (dp), dimension(:,:), intent(out)  xmn,
integer, intent(in)  parity 
)

Fourier transform a 2D quantity to fourier space.

Works by reshaping to a 3D array then passing to the normal 3D version. Have the 3D function call the 2D function for each slice.

Parameters
[in,out]thisA fourier_class instance.
[in]xuvFourier space quantity.
[out]xmnReal space quantity.
[in]parityFourier parity flag.

Definition at line 523 of file fourier.f90.

◆ fourier_tomnsp_2d_pest()

pure subroutine fourier::fourier_tomnsp_2d_pest ( class (fourier_class), intent(inout)  this,
real (dp), dimension(:,:), intent(in)  xuv,
real (dp), dimension(:,:), intent(out)  xmn,
real (dp), dimension(:,:), intent(in)  lambda,
  dup 
)

Fourier transform a 2D quantity to fourier space.

Works by reshaping to a 3D array then passing to the normal 3D version. Have the 3D function call the 2D function for each slice.

Parameters
[in,out]thisA fourier_class instance.
[in]xuvFourier space quantity.
[out]xmnReal space quantity.
[in]lambdaStellarator symmetric lambda.
[in]dupduStellarator asymmetric lambda.
[in]parityFourier parity flag.
[in]asymAdd stellarator asymmetric terms.

Definition at line 555 of file fourier.f90.

◆ fourier_tomnsp_2d_u()

pure subroutine fourier::fourier_tomnsp_2d_u ( class (fourier_class), intent(inout)  this,
real (dp), dimension(:,:), intent(in)  xuv 
)

Sum over the poloidal part.

Pest transforms use the same sum over v but different sums over u. To encourage code reuse break the u and v sums into separate methods. This version uses the regular poloidal sum.

Parameters
[in,out]thisA fourier_class instance.
[in]xuvFourier space quantity.

Definition at line 587 of file fourier.f90.

◆ fourier_tomnsp_2d_u_pest()

pure subroutine fourier::fourier_tomnsp_2d_u_pest ( class (fourier_class), intent(inout)  this,
real (dp), dimension(:,:), intent(in)  xuv,
real (dp), dimension(:,:), intent(in)  lambda,
real (dp), dimension(:,:), intent(in)  dupdu 
)

Sum over the poloidal part.

Pest transforms use the same sum over v but different sums over u. To encourage code reuse break the u and v sums into separate methods. This version uses the lambda poloidal sum.

Parameters
[in,out]thisA fourier_class instance.
[in]xuvFourier space quantity.
[in]lmnsStellarator symmetric lambda.
[in]lmncStellarator asymmetric lambda.
[in]asymAdd stellarator asymmetric terms.

Definition at line 625 of file fourier.f90.

◆ fourier_tomnsp_2d_v()

pure subroutine fourier::fourier_tomnsp_2d_v ( class (fourier_class), intent(inout)  this,
real (dp), dimension(:,:), intent(out)  xmn,
integer, intent(in)  parity 
)

Sum over toroidal angle.

Pest transforms use the same sum over v but different sums over u. To encourage code reuse break the u and v sums into separate methods. This impliments the common toroidal sum.

Parameters
[in]thisA fourier_class instance.
[in]xuvFourier space quantity.
[out]xmnReal space quantity.
[in]parityFourier parity flag.

Definition at line 707 of file fourier.f90.

◆ fourier_tomnsp_3d()

subroutine fourier::fourier_tomnsp_3d ( class (fourier_class), intent(inout)  this,
real (dp), dimension(:,:,:), intent(in)  xuv,
real (dp), dimension(:,:,:), intent(out)  xmn,
integer, intent(in)  parity 
)

Convert a quantity from real to fourier space in parallel.

This subroutine moves a 3D quantity to fourier space by summing over its real represenetation for each surface.

See also
fourier_tomnsp_2d
Parameters
[in,out]thisA fourier_class instance.
[in]xuvReal space quantity.
[out]xmnFourier space quantity.
[in]parityFourier parity flag.

Definition at line 399 of file fourier.f90.

◆ fourier_tomnsp_3d_pest()

subroutine fourier::fourier_tomnsp_3d_pest ( class (fourier_class), intent(inout)  this,
real (dp), dimension(:,:,:), intent(in)  xuv,
real (dp), dimension(:,:,:), intent(out)  xmn,
real (dp), dimension(:,:,:), intent(in)  lmns,
real (dp), dimension(:,:,:), intent(in)  lmnc 
)

Convert a quantity from real to fourier space in pest coordinates.

This subroutine moves a 3D quantity to fourier space by summing over its real represenetation for each surface.

See also
fourier_tomnsp_2d_pest.
Parameters
[in,out]thisA fourier_class instance.
[in]xuvReal space quantity.
[out]xmnFourier space quantity.
[in]lmnsStellarator symmetric lambda.
[in]lmncStellarator asymmetric lambda.
[in]parityFourier parity flag.
[in]asymAdd stellarator asymmetric terms.

Definition at line 444 of file fourier.f90.

◆ test_fourier()

logical function fourier::test_fourier

Test fourier routines.

This tests the routines by fourier transforming form flux to real space then back again.

Definition at line 1050 of file fourier.f90.