time_advance Module



Variables

Type Visibility Attributes Name Initial
logical, private :: time_advance_initialized = .false.
logical, private :: wdriftinit = .false.
logical, private :: wstarinit = .false.
logical, private :: parnlinit = .false.
logical, private :: readinit = .false.
logical, private :: radialinit = .false.
logical, private :: driftimpinit = .false.
real, private, dimension(:, :), allocatable :: par_nl_fac
real, private, dimension(:, :), allocatable :: d_par_nl_fac_dr
real, private, dimension(:, :), allocatable :: par_nl_curv
real, private, dimension(:, :), allocatable :: d_par_nl_curv_dr
real, private, dimension(:), allocatable :: par_nl_driftx
real, private, dimension(:), allocatable :: par_nl_drifty
real, private, dimension(:), allocatable :: d_par_nl_driftx_dr
real, private, dimension(:), allocatable :: d_par_nl_drifty_dr
real, public, dimension(2, 10) :: time_gke = 0.
real, public, dimension(2, 2) :: time_parallel_nl = 0.

Interfaces

private interface get_dgdy

  • private subroutine get_dgdy_2d(g, dgdy)

    compute dg/dy in k-space accepts g(ky,kx)

    Arguments

    Type IntentOptional Attributes Name
    complex, intent(in), dimension(:, :) :: g
    complex, intent(out), dimension(:, :) :: dgdy
  • private subroutine get_dgdy_3d(g, dgdy)

    compute dg/dy in k-space accepts g(ky,kx,z,tube)

    Arguments

    Type IntentOptional Attributes Name
    complex, intent(in), dimension(:, :, -nzgrid:, :) :: g
    complex, intent(out), dimension(:, :, -nzgrid:, :) :: dgdy
  • private subroutine get_dgdy_4d(g, dgdy)

    compute dg/dy in k-space accepts g(ky,kx,z,tube,(vpa,mu,spec))

    Arguments

    Type IntentOptional Attributes Name
    complex, intent(in), dimension(:, :, -nzgrid:, :, vmu_lo%llim_proc:) :: g
    complex, intent(out), dimension(:, :, -nzgrid:, :, vmu_lo%llim_proc:) :: dgdy

private interface get_dgdx

  • private subroutine get_dgdx_2d(g, dgdx)

    compute dg/dx in k-space accepts g(ky,kx)

    Arguments

    Type IntentOptional Attributes Name
    complex, intent(in), dimension(:, :) :: g
    complex, intent(out), dimension(:, :) :: dgdx
  • private subroutine get_dgdx_3d(g, dgdx)

    compute dg/dx in k-space accepts g(ky,kx,z,tube)

    Arguments

    Type IntentOptional Attributes Name
    complex, intent(in), dimension(:, :, -nzgrid:, :) :: g
    complex, intent(out), dimension(:, :, -nzgrid:, :) :: dgdx
  • private subroutine get_dgdx_4d(g, dgdx)

    compute dg/dx in k-space accepts g(ky,kx,z,tube,(vpa,mu,spec))

    Arguments

    Type IntentOptional Attributes Name
    complex, intent(in), dimension(:, :, -nzgrid:, :, vmu_lo%llim_proc:) :: g
    complex, intent(out), dimension(:, :, -nzgrid:, :, vmu_lo%llim_proc:) :: dgdx

public interface checksum

  • private subroutine checksum_field(field, total)

    Arguments

    Type IntentOptional Attributes Name
    complex, intent(in), dimension(:, :, -nzgrid:, :) :: field
    real, intent(out) :: total
  • private subroutine checksum_dist(dist, total, norm)

    Arguments

    Type IntentOptional Attributes Name
    complex, intent(in), dimension(:, :, -nzgrid:, :, vmu_lo%llim_proc:) :: dist
    real, intent(out) :: total
    logical, intent(in), optional :: norm

Subroutines

public subroutine init_time_advance()

read time_advance_knobs namelist from the input file; sets the explicit time advance option, as well as allows for scaling of the x and y components of the magnetic drifts and of the drive term allocate distribution function sized arrays needed, e.g., for Runge-Kutta time advance set up neoclassical corrections to the equilibrium Maxwellian; only calculated/needed when simulating higher order terms in rhostar for intrinsic rotation calculate the term multiplying dg/dvpa in the mirror term and set up either the semi-Lagrange machinery or the tridiagonal matrix to be inverted if solving implicitly calculate the term multiplying dg/dz in the parallel streaming term and set up the tridiagonal matrix to be inverted if solving implicitly allocate and calculate the factors multiplying dg/dx, dg/dy, dphi/dx and dphi/dy in the magnetic drift terms allocate and calculate the factor multiplying dphi/dy in the gradient drive term

Arguments

None

private subroutine init_wdrift()

allocate wdriftx_phi, the factor multiplying dphi/dx in the magnetic drift term allocate wdrifty_phi, the factor multiplying dphi/dy in the magnetic drift term allocate wdriftx_bpar, the factor multiplying dbpar/dx in the magnetic drift term allocate wdrifty_bpar, the factor multiplying dbpar/dy in the magnetic drift term allocate wdriftx_g, the factor multiplying dg/dx in the magnetic drift term allocate wdrifty_g, the factor multiplying dg/dy in the magnetic drift term this is the curvature drift piece of wdrifty with missing factor of vpa vpa factor is missing to avoid singularity when including non-Maxwellian corrections to equilibrium this is the grad-B drift piece of wdrifty if including neoclassical correction to equilibrium Maxwellian, then add in v_E^{nc} . grad y dg/dy coefficient here if maxwwellian_normalization = .true., evolved distribution function is normalised by a Maxwellian otherwise, it is not; a Maxwellian weighting factor must thus be included assign wdrifty_bpar, neoclassical terms not supported if including neoclassical corrections to equilibrium, add in -(Ze/m) * v_curv/vpa . grad y d/dy * dF^{nc}/dvpa term and v_E . grad z dF^{nc}/dz (here get the dphi/dy part of v_E) NB: the below neoclassical correction needs to be divided by an equilibrium Maxwellian if maxwellian_normalization = .true. this is the curvature drift piece of wdriftx with missing factor of vpa vpa factor is missing to avoid singularity when including non-Maxwellian corrections to equilibrium this is the grad-B drift piece of wdriftx if including neoclassical correction to equilibrium Maxwellian, then add in v_E^{nc} . grad x dg/dx coefficient here if maxwellian_normalizatiion = .true., evolved distribution function is normalised by a Maxwellian otherwise, it is not; a Maxwellian weighting factor must thus be included assign wdriftx_bpar, neoclassical terms not supported if including neoclassical corrections to equilibrium, add in (Ze/m) * v_curv/vpa . grad x d/dx * dF^{nc}/dvpa term and v_E . grad z dF^{nc}/dz (here get the dphi/dx part of v_E) and v_E . grad alpha dF^{nc}/dalpha (dphi/dx part of v_E) NB: the below neoclassical correction needs to be divided by an equilibrium Maxwellian if running with maxwellian_normalzation = .true.

Arguments

None

private subroutine init_wstar()

Arguments

None

private subroutine init_parallel_nonlinearity()

Arguments

None

private subroutine init_radial_variation()

Arguments

None

private subroutine allocate_arrays()

Arguments

None

private subroutine init_cfl()

TODO:GA- add correct CFL condition

Arguments

None

private subroutine reset_dt()

Arguments

None

public subroutine advance_stella(istep, stop_stella)

unless running in multibox mode, no need to worry about mb_communicate calls as the subroutine is immediately exited if not in multibox mode. save value of phi & apar for use in diagnostics (to obtain frequency) reverse the order of operations every time step as part of alternating direction operator splitting this is needed to ensure 2nd order accuracy in time

Read more…

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: istep
logical, intent(inout) :: stop_stella

private subroutine advance_explicit(g, restart_time_step, istep)

advance_explicit takes as input the guiding centre distribution function in k-space and updates it to account for all of the terms in the GKE that are advanced explicitly in time

Read more…

Arguments

Type IntentOptional Attributes Name
complex, intent(inout), dimension(:, :, -nzgrid:, :, vmu_lo%llim_proc:) :: g
logical, intent(inout) :: restart_time_step
integer, intent(in) :: istep

private subroutine advance_explicit_euler(g, restart_time_step, istep)

advance_explicit_euler uses forward Euler to advance one time step

Read more…

Arguments

Type IntentOptional Attributes Name
complex, intent(inout), dimension(:, :, -nzgrid:, :, vmu_lo%llim_proc:) :: g
logical, intent(inout) :: restart_time_step
integer, intent(in) :: istep

private subroutine advance_explicit_rk2(g, restart_time_step, istep)

advance_expliciit_rk2 uses strong stability-preserving RK2 to advance one time step

Read more…

Arguments

Type IntentOptional Attributes Name
complex, intent(inout), dimension(:, :, -nzgrid:, :, vmu_lo%llim_proc:) :: g
logical, intent(inout) :: restart_time_step
integer, intent(in) :: istep

private subroutine advance_explicit_rk3(g, restart_time_step, istep)

strong stability-preserving RK3

Read more…

Arguments

Type IntentOptional Attributes Name
complex, intent(inout), dimension(:, :, -nzgrid:, :, vmu_lo%llim_proc:) :: g
logical, intent(inout) :: restart_time_step
integer, intent(in) :: istep

private subroutine advance_explicit_rk4(g, restart_time_step, istep)

standard RK4

Read more…

Arguments

Type IntentOptional Attributes Name
complex, intent(inout), dimension(:, :, -nzgrid:, :, vmu_lo%llim_proc:) :: g
logical, intent(inout) :: restart_time_step
integer, intent(in) :: istep

private subroutine solve_gke(pdf, rhs_ky, restart_time_step, istep)

solve_gke accepts as argument pdf, the guiding centre distribution function in k-space, and returns rhs_ky, the right-hand side of the gyrokinetic equation in k-space; i.e., if dg/dt = r, then rhs_ky = rdt; note that if include_apar = T, then the input pdf is actually gbar = g + (Ze/T)(vpa/c)F0

Read more…

Arguments

Type IntentOptional Attributes Name
complex, intent(inout), dimension(:, :, -nzgrid:, :, vmu_lo%llim_proc:) :: pdf
complex, intent(out), dimension(:, :, -nzgrid:, :, vmu_lo%llim_proc:), target :: rhs_ky
logical, intent(out) :: restart_time_step
integer, intent(in) :: istep

private subroutine advance_hyper_explicit(gin, gout)

Arguments

Type IntentOptional Attributes Name
complex, intent(in), dimension(:, :, -nzgrid:, :, vmu_lo%llim_proc:) :: gin
complex, intent(inout), dimension(:, :, -nzgrid:, :, vmu_lo%llim_proc:) :: gout

private subroutine advance_wstar_explicit(phi, gout)

start timing the time advance due to the driving gradients

Read more…

Arguments

Type IntentOptional Attributes Name
complex, intent(in), dimension(:, :, -nzgrid:, :) :: phi
complex, intent(inout), dimension(:, :, -nzgrid:, :, vmu_lo%llim_proc:) :: gout

private subroutine advance_wdrifty_explicit(g, phi, bpar, gout)

advance_wdrifty_explicit subroutine calculates and adds the y-component of the magnetic drift term to the RHS of the GK equation

Read more…

Arguments

Type IntentOptional Attributes Name
complex, intent(in), dimension(:, :, -nzgrid:, :, vmu_lo%llim_proc:) :: g
complex, intent(in), dimension(:, :, -nzgrid:, :) :: phi
complex, intent(in), dimension(:, :, -nzgrid:, :) :: bpar
complex, intent(inout), dimension(:, :, -nzgrid:, :, vmu_lo%llim_proc:) :: gout

private subroutine advance_wdriftx_explicit(g, phi, bpar, gout)

advance_wdriftx_explicit subroutine calculates and adds the x-component of the magnetic drift term to the RHS of the GK equation

Read more…

Arguments

Type IntentOptional Attributes Name
complex, intent(in), dimension(:, :, -nzgrid:, :, vmu_lo%llim_proc:) :: g
complex, intent(in), dimension(:, :, -nzgrid:, :) :: phi
complex, intent(in), dimension(:, :, -nzgrid:, :) :: bpar
complex, intent(inout), dimension(:, :, -nzgrid:, :, vmu_lo%llim_proc:) :: gout

private subroutine advance_ExB_nonlinearity(g, gout, restart_time_step, istep)

compute phase factor needed when running with equilibrium flow shear compute ikyg FFT to get dg/dy in (y,x) space compute ikx zero out the zonal contribution to d/dx if requested if running with equilibrium flow shear, make adjustment to the term multiplying dg/dy FFT to get d/dx in (y,x) space multiply by the geometric factor appearing in the Poisson bracket; i.e., (dx/dpsidy/dalpha)0.5 compute the contribution to the Poisson bracket from dg/dy*d/dx

Read more…

Arguments

Type IntentOptional Attributes Name
complex, intent(inout), dimension(:, :, -nzgrid:, :, vmu_lo%llim_proc:) :: g
complex, intent(inout), dimension(:, :, -nzgrid:, :, vmu_lo%llim_proc:) :: gout
logical, intent(out) :: restart_time_step
integer, intent(in) :: istep

private subroutine advance_parallel_nonlinearity(g, gout, restart_time_step)

check estimated cfl_dt to see if the time step size needs to be changed

Arguments

Type IntentOptional Attributes Name
complex, intent(in), dimension(:, :, -nzgrid:, :, vmu_lo%llim_proc:) :: g
complex, intent(inout), dimension(:, :, -nzgrid:, :, vmu_lo%llim_proc:) :: gout
logical, intent(out) :: restart_time_step

private subroutine advance_radial_variation(g, gout)

Arguments

Type IntentOptional Attributes Name
complex, intent(in), dimension(:, :, -nzgrid:, :, vmu_lo%llim_proc:) :: g
complex, intent(inout), dimension(:, :, -nzgrid:, :, vmu_lo%llim_proc:) :: gout

private subroutine get_dgdy_2d(g, dgdy)

compute dg/dy in k-space accepts g(ky,kx)

Arguments

Type IntentOptional Attributes Name
complex, intent(in), dimension(:, :) :: g
complex, intent(out), dimension(:, :) :: dgdy

private subroutine get_dgdy_3d(g, dgdy)

compute dg/dy in k-space accepts g(ky,kx,z,tube)

Arguments

Type IntentOptional Attributes Name
complex, intent(in), dimension(:, :, -nzgrid:, :) :: g
complex, intent(out), dimension(:, :, -nzgrid:, :) :: dgdy

private subroutine get_dgdy_4d(g, dgdy)

compute dg/dy in k-space accepts g(ky,kx,z,tube,(vpa,mu,spec))

Arguments

Type IntentOptional Attributes Name
complex, intent(in), dimension(:, :, -nzgrid:, :, vmu_lo%llim_proc:) :: g
complex, intent(out), dimension(:, :, -nzgrid:, :, vmu_lo%llim_proc:) :: dgdy

private subroutine get_dgdx_2d(g, dgdx)

compute dg/dx in k-space accepts g(ky,kx)

Arguments

Type IntentOptional Attributes Name
complex, intent(in), dimension(:, :) :: g
complex, intent(out), dimension(:, :) :: dgdx

private subroutine get_dgdx_3d(g, dgdx)

compute dg/dx in k-space accepts g(ky,kx,z,tube)

Arguments

Type IntentOptional Attributes Name
complex, intent(in), dimension(:, :, -nzgrid:, :) :: g
complex, intent(out), dimension(:, :, -nzgrid:, :) :: dgdx

private subroutine get_dgdx_4d(g, dgdx)

compute dg/dx in k-space accepts g(ky,kx,z,tube,(vpa,mu,spec))

Arguments

Type IntentOptional Attributes Name
complex, intent(in), dimension(:, :, -nzgrid:, :, vmu_lo%llim_proc:) :: g
complex, intent(out), dimension(:, :, -nzgrid:, :, vmu_lo%llim_proc:) :: dgdx

private subroutine add_explicit_term(g, pre_factor, src)

Arguments

Type IntentOptional Attributes Name
complex, intent(in), dimension(:, :, -nzgrid:, :, vmu_lo%llim_proc:) :: g
real, intent(in), dimension(-nzgrid:, vmu_lo%llim_proc:) :: pre_factor
complex, intent(inout), dimension(:, :, -nzgrid:, :, vmu_lo%llim_proc:) :: src

private subroutine add_explicit_term_ffs(g, pre_factor, src)

add vM . grad y d/dy or vM . grad x d/dx (or equivalents with g) or omega_* * d/dy term to RHS of GK equation

Arguments

Type IntentOptional Attributes Name
complex, intent(in), dimension(:, :, -nzgrid:, :, vmu_lo%llim_proc:) :: g
real, intent(in), dimension(:, -nzgrid:, vmu_lo%llim_proc:) :: pre_factor
complex, intent(inout), dimension(:, :, -nzgrid:, :, vmu_lo%llim_proc:) :: src

private subroutine advance_implicit(istep, phi, apar, bpar, g)

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: istep
complex, intent(inout), dimension(:, :, -nzgrid:, :) :: phi
complex, intent(inout), dimension(:, :, -nzgrid:, :) :: apar
complex, intent(inout), dimension(:, :, -nzgrid:, :) :: bpar
complex, intent(inout), dimension(:, :, -nzgrid:, :, vmu_lo%llim_proc:) :: g

private subroutine mb_communicate(g_in)

Arguments

Type IntentOptional Attributes Name
complex, intent(inout), dimension(:, :, -nzgrid:, :, vmu_lo%llim_proc:) :: g_in

private subroutine checksum_field(field, total)

Arguments

Type IntentOptional Attributes Name
complex, intent(in), dimension(:, :, -nzgrid:, :) :: field
real, intent(out) :: total

private subroutine checksum_dist(dist, total, norm)

Arguments

Type IntentOptional Attributes Name
complex, intent(in), dimension(:, :, -nzgrid:, :, vmu_lo%llim_proc:) :: dist
real, intent(out) :: total
logical, intent(in), optional :: norm

public subroutine finish_time_advance()

Arguments

None

private subroutine finish_parallel_nonlinearity()

Arguments

None

private subroutine finish_wdrift()

Arguments

None

private subroutine finish_wstar()

Arguments

None

private subroutine deallocate_arrays()

Arguments

None