collisions Module

Routines for implementing the model collision operator defined by Barnes, Abel et al. 2009 or on arxiv. The collision operator causes physically motivated smoothing of structure in velocity space which is necessary to prevent buildup of structure at fine scales in velocity space, while conserving energy and momentum.



Contents


Variables

Type Visibility Attributes Name Initial
logical, public :: use_le_layout
logical, private :: const_v
logical, private :: conserve_moments
logical, private :: conservative
logical, private :: resistivity
integer, public :: collision_model_switch
integer, private :: lorentz_switch
integer, private :: ediff_switch
logical, public :: adjust
logical, public :: heating
logical, public :: hyper_colls
logical, private :: ei_coll_only
logical, private :: test
logical, private :: special_wfb_lorentz
logical, private :: vpar_zero_mean
logical, private :: conserve_forbid_zero
integer, public, parameter :: collision_model_lorentz = 1
integer, public, parameter :: collision_model_none = 3
integer, public, parameter :: collision_model_lorentz_test = 5
integer, public, parameter :: collision_model_full = 6
integer, public, parameter :: collision_model_ediffuse = 7
integer, private, parameter :: lorentz_scheme_default = 1
integer, private, parameter :: lorentz_scheme_old = 2
integer, private, parameter :: ediff_scheme_default = 1
integer, private, parameter :: ediff_scheme_old = 2
real, public, dimension (2) :: vnmult = -1.0
integer, public :: ncheck
logical, public :: vary_vnew
real, private :: vnfac
real, private :: vnslow
real, private :: etol
real, private :: ewindow
real, private :: etola
real, private :: ewindowa
integer, private :: timesteps_between_collisions
logical, private :: force_collisions
logical, private :: has_lorentz
logical, private :: has_diffuse
real, private, dimension (2) :: vnm_init = 1.0
real, public, dimension (:,:,:), allocatable :: dtot
real, public, dimension (:,:), allocatable :: fdf
real, public, dimension (:,:), allocatable :: fdb
real, private, dimension (:,:,:), allocatable :: vnew
real, private, dimension (:,:,:), allocatable :: vnew_s
real, private, dimension (:,:,:), allocatable :: vnew_D
real, private, dimension (:,:,:), allocatable :: vnew_E
real, private, dimension (:,:,:), allocatable :: delvnew
real, private, dimension (:,:), allocatable :: vpdiffle
real, private, dimension (:,:,:), allocatable :: vpdiff
real, private, dimension (:,:,:,:), allocatable :: vnewh
real, private, dimension(:,:,:), allocatable :: s0
real, private, dimension(:,:,:), allocatable :: w0
real, private, dimension(:,:,:), allocatable :: z0
real, private, dimension (:,:,:), allocatable :: s0le
real, private, dimension (:,:,:), allocatable :: w0le
real, private, dimension (:,:,:), allocatable :: z0le
real, private, dimension (:,:,:), allocatable :: aj0le
real, private, dimension (:,:,:), allocatable :: vperp_aj1le
real, private, dimension (:,:,:), allocatable :: vpa_aj0_le
real, private, dimension(:,:,:), allocatable :: bs0
real, private, dimension(:,:,:), allocatable :: bw0
real, private, dimension(:,:,:), allocatable :: bz0
real, private, dimension (:,:,:), allocatable :: bs0le
real, private, dimension (:,:,:), allocatable :: bw0le
real, private, dimension (:,:,:), allocatable :: bz0le
real, private :: cfac
integer, private :: nxi_lim

Sets the upper xi index which we consider in loops

real, private, dimension (:,:,:), allocatable :: c1le
real, private, dimension (:,:,:), allocatable :: betaale
real, private, dimension (:,:,:), allocatable :: qle
real, private, dimension (:,:,:), allocatable :: d1le
real, private, dimension (:,:,:), allocatable :: h1le
real, private, dimension (:,:,:), allocatable :: ec1le
real, private, dimension (:,:,:), allocatable :: ebetaale
real, private, dimension (:,:,:), allocatable :: eqle
real, private, dimension (:,:), allocatable :: c1
real, private, dimension (:,:), allocatable :: betaa
real, private, dimension (:,:), allocatable :: ql
real, private, dimension (:,:), allocatable :: d1
real, private, dimension (:,:), allocatable :: h1
real, private, dimension (:,:), allocatable :: ec1
real, private, dimension (:,:), allocatable :: ebetaa
real, private, dimension (:,:), allocatable :: eql
real, private, dimension (:, :), allocatable :: pitch_weights
logical, private :: drag = .false.
logical, public :: colls = .true.
logical, public :: split_collisions
logical, private :: hypermult
logical, private :: initialized = .false.
logical, private :: exist
type(collisions_config_type), private :: collisions_config

Interfaces

public interface solfp1

  • private subroutine solfp1_le_layout(gle, diagnostics)

    FIXME : Add documentation

    Arguments

    Type IntentOptional Attributes Name
    complex, intent(inout), dimension (:,:,le_lo%llim_proc:) :: gle
    integer, intent(in), optional :: diagnostics
  • private subroutine solfp1_standard_layout(g, g1, gc1, gc2, diagnostics, gtoc, ctog)

    FIXME : Add documentation

    Arguments

    Type IntentOptional Attributes Name
    complex, intent(inout), dimension (-ntgrid:,:,g_lo%llim_proc:) :: g
    complex, intent(inout), dimension (-ntgrid:,:,g_lo%llim_proc:) :: g1
    complex, intent(inout), dimension (-ntgrid:,:,g_lo%llim_proc:) :: gc1
    complex, intent(inout), dimension (-ntgrid:,:,g_lo%llim_proc:) :: gc2
    integer, intent(in), optional :: diagnostics
    logical, intent(in), optional :: gtoc
    logical, intent(in), optional :: ctog

private interface solfp_lorentz

  • private subroutine solfp_lorentz_le_layout(gle, diagnostics)

    FIXME : Add documentation

    Arguments

    Type IntentOptional Attributes Name
    complex, intent(inout), dimension (:,:,le_lo%llim_proc:) :: gle
    integer, intent(in), optional :: diagnostics
  • private subroutine solfp_lorentz_standard_layout(g, gc, gh, diagnostics)

    FIXME : Add documentation

    Arguments

    Type IntentOptional Attributes Name
    complex, intent(inout), dimension (-ntgrid:,:,g_lo%llim_proc:) :: g
    complex, intent(inout), dimension (-ntgrid:,:,g_lo%llim_proc:) :: gc
    complex, intent(inout), dimension (-ntgrid:,:,g_lo%llim_proc:) :: gh
    integer, intent(in), optional :: diagnostics

private interface conserve_lorentz

  • private subroutine conserve_lorentz_le_layout(gle)

    FIXME : Add documentation !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! First get v0y0 and y1 = y0 - v0y0 * z0 / (1 + v0z0) v0 = vpa J0 f0, y0 = gle $OMP PARALLEL DO DEFAULT(none) & $OMP PRIVATE(ile, it, ik, ixi, ie, is, ig, v0y0) & $OMP SHARED(le_lo, kwork_filter, ieqzip, negrid, nxi_lim, & $OMP z0le, w, wxi, vpa_aj0_le, gle) & $OMP SCHEDULE(static) !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! Now get v1y1 and y2 = y1 - v1y1 * s1 / (1 + v1s1) v1 = nud vpa J0 f0, y1 = gle

    !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! Now get v2y2 and x = y2 - v2y2 * w2 / (1 + v2w2) $OMP PARALLEL DO DEFAULT(none) & $OMP PRIVATE(ile, it, ik, is, ie, ig, ixi, nud, v2y2) & $OMP SHARED(le_lo, kwork_filter, ieqzip, negrid, nxi_lim, & $OMP vnmult, vnew_D, vperp_aj1le, gle, w, wxi, w0le) & $OMP SCHEDULE(static)

    Arguments

    Type IntentOptional Attributes Name
    complex, intent(inout), dimension (:,:,le_lo%llim_proc:) :: gle
  • private subroutine conserve_lorentz_standard_layout(g, g1)

    FIXME : Add documentation !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! First get v0y0 $OMP PARALLEL DO DEFAULT(none) & $OMP PRIVATE(iglo, it, ik, isgn) & $OMP SHARED(g_lo, kwork_filter, gtmp, vpa, aj0, g) & $OMP SCHEDULE(static) !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! Get y1 = y0 - v0y0 * z0 / (1 + v0z0) $OMP PARALLEL DO DEFAULT(none) & $OMP PRIVATE(iglo, it, ik, is, isgn) & $OMP SHARED(g_lo, kwork_filter, g1, g, ieqzip, v0y0, z0) & $OMP SCHEDULE(static) !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! Now get v1y1 $OMP PARALLEL DO DEFAULT(none) & $OMP PRIVATE(iglo, it, ik, il, ie, is, isgn) & $OMP SHARED(g_lo, kwork_filter, conservative, gtmp, vns, speed, vpdiff, aj0, g1, vpa) & $OMP SCHEDULE(static) !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! Get y2 = y1 - v1y1 * s1 / (1 + v1s1) $OMP PARALLEL DO DEFAULT(none) & $OMP PRIVATE(iglo, it, ik, is, isgn) & $OMP SHARED(g_lo, kwork_filter, g1, ieqzip, v1y1, s0) & $OMP SCHEDULE(static) !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! Now get v2y2 $OMP PARALLEL DO DEFAULT(none) & $OMP PRIVATE(iglo, it, ik, ie, il, is, isgn) & $OMP SHARED(g_lo, kwork_filter, gtmp, vns, energy, al, aj1, g1) & $OMP SCHEDULE(static) !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! Finally get x = y2 - v2y2 * w2 / (1 + v2w2) $OMP PARALLEL DO DEFAULT(none) & $OMP PRIVATE(iglo, it, ik, is, isgn) & $OMP SHARED(g_lo, kwork_filter, g, g1, ieqzip, v2y2, w0) & $OMP SCHEDULE(static)

    Arguments

    Type IntentOptional Attributes Name
    complex, intent(inout), dimension (-ntgrid:,:,g_lo%llim_proc:) :: g
    complex, intent(inout), dimension (-ntgrid:,:,g_lo%llim_proc:) :: g1

private interface conserve_diffuse

  • private subroutine conserve_diffuse_le_layout(gle)

    FIXME : Add documentation !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! First get v0y0 and then y1 = y0 - v0y0 * z0 / (1 + v0z0) $OMP PARALLEL DO DEFAULT(none) & $OMP PRIVATE(ile, it, ik, is, ie, ig, ixi, vnE, v0y0) & $OMP SHARED(le_lo, kwork_filter, negrid, nxi_lim, w, wxi, ieqzip, bz0le, vnmult, vnew_E, aj0le, gle) & $OMP SCHEDULE(static) !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! Now get v1y1 and then y2 = y1 - v1y1 * s1 / (1 + v1s1) $OMP PARALLEL DO DEFAULT(none) & $OMP PRIVATE(ile, it, ik, is, ig, delnu, ie, ixi, v1y1) & $OMP SHARED(le_lo, kwork_filter, negrid, vnmult, delvnew, nxi_lim, & $OMP bs0le, vpa_aj0_le, gle, w, wxi, ieqzip) & $OMP SCHEDULE(static) !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! Now get v2y2 and then get x = y2 - v2y2 * w2 / (1 + v2w2) $OMP PARALLEL DO DEFAULT(none) & $OMP PRIVATE(ile, it, ik, is, ie, ig, ixi, delnu, v2y2) & $OMP SHARED(le_lo, kwork_filter, negrid, nxi_lim, ieqzip, & $OMP bw0le, delvnew, vperp_aj1le, gle, vnmult, w, wxi) & $OMP SCHEDULE(static)

    Arguments

    Type IntentOptional Attributes Name
    complex, intent(inout), dimension (:,:,le_lo%llim_proc:) :: gle
  • private subroutine conserve_diffuse_standard_layout(g, g1)

    FIXME : Add documentation !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! First get v0y0 $OMP PARALLEL DO DEFAULT(none) & $OMP PRIVATE(iglo, it, ik, ie, is, isgn) & $OMP SHARED(g_lo, kwork_filter, gtmp, vnmult, vnew_E, aj0, g) & $OMP SCHEDULE(static) !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! Get y1 = y0 - v0y0 * z0 / (1 + v0z0) $OMP PARALLEL DO DEFAULT(none) & $OMP PRIVATE(iglo, it, ik, is, isgn) & $OMP SHARED(g_lo, kwork_filter, g1, g, ieqzip, v0y0, bz0) & $OMP SCHEDULE(static) !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! Now get v1y1 $OMP PARALLEL DO DEFAULT(none) & $OMP PRIVATE(iglo, it, ik, ie, is, isgn) & $OMP SHARED(g_lo, kwork_filter, gtmp, vns, vpa, aj0, g1) & $OMP SCHEDULE(static) !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! Get y2 = y1 - v1y1 * s1 / (1 + v1s1) $OMP PARALLEL DO DEFAULT(none) & $OMP PRIVATE(iglo, it, ik, is, isgn) & $OMP SHARED(g_lo, kwork_filter, g1, ieqzip, v1y1, bs0) & $OMP SCHEDULE(static) !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! Now get v2y2 $OMP PARALLEL DO DEFAULT(none) & $OMP PRIVATE(iglo, it, ik, ie, il, is, isgn) & $OMP SHARED(g_lo, kwork_filter, gtmp, vns, energy, al, aj1, g1) & $OMP SCHEDULE(static) !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! Finally get x = y2 - v2y2 * w2 / (1 + v2w2) $OMP PARALLEL DO DEFAULT(none) & $OMP PRIVATE(iglo, it, ik, is, isgn) & $OMP SHARED(g_lo, kwork_filter, g, g1, ieqzip, v2y2, bw0) & $OMP SCHEDULE(static)

    Arguments

    Type IntentOptional Attributes Name
    complex, intent(inout), dimension (-ntgrid:,:,g_lo%llim_proc:) :: g
    complex, intent(inout), dimension (-ntgrid:,:,g_lo%llim_proc:) :: g1

Derived Types

type, public, extends(abstract_config_type) ::  collisions_config_type

Used to represent the input configuration of collisions

Components

Type Visibility Attributes Name Initial
logical, public :: exist = .false.

Does the related namelist exist in the target input file?

integer, public :: index = 0

Used to hold the specific index of numbered namelists

logical, public :: skip_read = .false.

Do we want to skip the read step in init?

logical, public :: skip_broadcast = .false.

Do we want to skip the broadcast step in init?

logical, public :: adjust = .true.

If true (default) then transform from the gyro-averaged distribution function evolved by GS2 to the non-Boltzmann part of , (), when applying the collision operator. This is the physically appropriate choice, this parameter is primarily for numerical testing.

real, public :: cfac = 1.0

Factor multipyling the finite Larmor radius terms in the collision operator. This term is essentially just classical diffusion. Set cfac to 0 to turn off this diffusion.

Read more…
character(len=20), public :: collision_model = 'default'

Selects the collision model used in the simulation. Can be one of

Read more…
logical, public :: conservative = .true.

If true (default) then guarantee exact conservation properties.

logical, public :: conserve_forbid_zero = .true.

If true (default) then forces conservation corrections to zero in the forbidden region to avoid introducing unphysical non-zero values for the distribution function in the forbidden region of trapped particles.

Read more…
logical, public :: conserve_moments = .true.

If true (default) then guarantee collision operator conserves momentum and energy.

Read more…
logical, public :: const_v = .false.

If true (not the default) then use the thermal velocity to evaluate the collision frequencies to be used. This results in an energy independent collision frequency being used for all species.

character(len=20), public :: ediff_scheme = 'default'

Controls how the coefficients in the matrix representing the energy diffusion operator are obtained. Can be one of

Read more…
logical, public :: ei_coll_only = .false.

If true (not the default) then force the collision frequency used for all non-electron species to zero and force all electron-electron terms to zero.

real, public :: etol = 2.e-2

Only used in get_verr as a part of the adaptive collisionality algorithm. Sets the maximum relative error allowed, above which the collision frequency must be increased.

Read more…
real, public :: etola = 2.e-2

Only used in get_verr as a part of the adaptive collisionality algorithm. Sets the maximum absolute error allowed, above which the collision frequency must be increased.

Read more…
real, public :: ewindow = 1.e-2

Only used in get_verr as a part of the adaptive collisionality algorithm. Sets an offset to apply to the relative error limit set by etol. This is used to provide hysteresis is the adaptive collisionality algorithm so to avoid adjusting the collision frequency up and down every step (similar to delt_cushion).

Read more…
real, public :: ewindowa = 1.e-2

Only used in get_verr as a part of the adaptive collisionality algorithm. Sets an offset to apply to the absolute error limit set by etola. This is used to provide hysteresis is the adaptive collisionality algorithm so to avoid adjusting the collision frequency up and down every step (similar to the delt_cushion).

Read more…
logical, public :: force_collisions = .false.

Currently we skip the application of the selected collision operator for species with zero collision frequency and disable collisions entirely if no species have non-zero collision frequency. This is generally sensible, but it can be useful to force the use of the collisional code path in some situations such as in code testing. Setting this flag to .true. forces the selected collision model operator to be applied even if the collision frequency is zero.

logical, public :: heating = .false.

If true (not the default) then calculate collisional heating when applying the collion operator. This is purely a diagnostic calculation. It should not change the evolution.

Read more…
logical, public :: hypermult = .false.

If true (not the default) then multiply the hyper collision frequency by the species' collision frequency nu_h. This only impacts the pitch angle scattering operator.

Read more…
character(len=20), public :: lorentz_scheme = 'default'

Controls how the coefficients in the matrix representing the pitch angle scattering operator are obtained. Can be one of

Read more…
integer, public :: ncheck = 100

Used as a part of the adaptive collisionality algorithm. When active we check the velocity space error with get_verr every ncheck time steps. This check can be relatively expensive so it is recommended to avoid small values of ncheck.

Read more…
logical, public :: resistivity = .true.

If true (default) then potentially include the drag term in the pitch angle scattering operator. This is a necessary but not sufficient criteria. For the drag term to be included we also require , more than one simulated species and finite perturbations included in the simulation (i.e. fapar /= 0).

logical, public :: special_wfb_lorentz = .false.

If true (not the default) then use special handling for the wfb particle in the pitch angle scattering operator.

Read more…
logical, public :: split_collisions = .false.

If true (not the default) then remove the collision operator from the usual time advance algorithm. Instead the collision operator is only applied every timesteps_between_collisions timesteps. This can potentially substantially speed up collisional simulations, both in the initialisation and advance phases.

Read more…
logical, public :: test = .false.

If true (not the default) then performs some additional checks of the data redistribution routines associated with transforming being the standard and collisional data decompositions.

integer, public :: timesteps_between_collisions = 1

Should set the number of timesteps between application of the collision operator if split_collisions is true. Currently this is ignored.

Read more…
logical, public :: use_le_layout = .true.

If true (default) then use a data decomposition for collisions that brings both pitch angle and energy local to a processor. This is typically the most efficient option, however for collisional simulations that only use part of the collision operator (either energy diffusion or pitch angle scattering) then it may be beneficial to set this flag to false such that we use a decomposition that only brings either energy or pitch angle local.

logical, public :: vary_vnew = .false.

Set to true (not the default) to activate the adaptive collisionality algorithm.

Read more…
real, public :: vnfac = 1.1

If the collisionality is to be increased as a part of the adaptive collisionality algorithm then increase it by this factor.

real, public :: vnslow = 0.9

If the collisionality is to be decreased as a part of the adaptive collisionality algorithm then decrease it by this factor.

logical, public :: vpar_zero_mean = .true.

Controls how the duplicate vpar = 0 point is handled. When vpar_zero_mean = .true. (the default) the average of g(vpar = 0) for both signs of the parallel velcoity (isgn) is used in the collision operator instead of just g(vpar = 0) at isgn=2. This is seen to suppress a numerical instability when special_wfb_lorentz =.false.. With these defaults pitch angle scattering at is now being handled physically i.e. vpar = 0 at this theta location is no longer being skipped.

Read more…

Type-Bound Procedures

procedure , public , :: is_initialised => is_initialised_generic Function
procedure , public , :: init => init_generic Subroutine
procedure , public , :: setup => setup_generic Subroutine
procedure , public , :: write_namelist_header Subroutine
procedure , public , :: get_name => get_name_generic Function
procedure , public , :: get_requires_index => get_requires_index_generic Function
procedure , public , nopass :: write_namelist_footer Subroutine
generic, public , :: write_key_val => write_key_val_string, write_key_val_real, write_key_val_complex, write_key_val_integer, write_key_val_logical, write_key_val_real_array, write_key_val_complex_array, write_key_val_integer_array
procedure , public :: read => read_collisions_config Subroutine
procedure , public :: write => write_collisions_config Subroutine
procedure , public :: reset => reset_collisions_config Subroutine
procedure , public :: broadcast => broadcast_collisions_config Subroutine
procedure , public , nopass :: get_default_name => get_default_name_collisions_config Function
procedure , public , nopass :: get_default_requires_index => get_default_requires_index_collisions_config Function

Functions

private elemental function safe_sqrt(arg)

A wrapper to sqrt which replaces -ve values with 0.0 to avoid NaNs arising from slight floating point discrepancies. We could consider adding a debug check to abort/warn if the passed value is too negative (i.e. if it looks like an error rather than small round off).

Arguments

Type IntentOptional Attributes Name
real, intent(in) :: arg

Return Value real

private elemental function hsg_func(vel)

FIXME : Add documentation

Arguments

Type IntentOptional Attributes Name
real, intent(in) :: vel

Return Value real

Gets the default name for this namelist

Arguments

None

Return Value character(len=CONFIG_MAX_NAME_LEN)

Gets the default requires index for this namelist

Arguments

None

Return Value logical

public pure function get_collisions_config()

Get the module level config instance

Arguments

None

Return Value type(collisions_config_type)


Subroutines

public subroutine set_vnmult(vnmult_in)

Arguments

Type IntentOptional Attributes Name
real, intent(in), dimension(2) :: vnmult_in

public subroutine set_heating(heating_in)

Arguments

Type IntentOptional Attributes Name
logical, intent(in) :: heating_in

public subroutine check_collisions(report_unit)

FIXME : Add documentation

Arguments

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

public subroutine wnml_collisions(unit)

FIXME : Add documentation

Arguments

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

public subroutine init_collisions(collisions_config_in)

FIXME : Add documentation

Arguments

Type IntentOptional Attributes Name
type(collisions_config_type), intent(in), optional :: collisions_config_in

public subroutine read_parameters(collisions_config_in)

FIXME : Add documentation

Arguments

Type IntentOptional Attributes Name
type(collisions_config_type), intent(in), optional :: collisions_config_in

private subroutine init_arrays()

FIXME : Add documentation

Read more…

Arguments

None

private subroutine init_le_bessel()

Communicate Bessel functions from g_lo to le_lo Currently just aj1

Arguments

None

private subroutine init_lorentz_conserve()

Precompute three quantities needed for momentum and energy conservation: z0, w0, s0 (z0le, w0le, s0le if use_le_layout chosen in the input file defined) !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! Now get z0 (first form) $OMP PARALLEL DO DEFAULT(none) & $OMP PRIVATE(iglo, ik, ie, il, is, isgn) & $OMP SHARED(g_lo, conservative, z0tmp, code_dt, vns, vpdiff, speed, aj0, vpa) & $OMP SCHEDULE(static) !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! Now get v0z0 $OMP PARALLEL DO DEFAULT(none) & $OMP PRIVATE(iglo, isgn) & $OMP SHARED(g_lo, vpa, gtmp, aj0, z0tmp) & $OMP COLLAPSE(2) & $OMP SCHEDULE(static) !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! Redefine z0 = z0 / (1 + v0z0) $OMP PARALLEL DO DEFAULT(none) & $OMP PRIVATE(iglo, ik, it, is, isgn) & $OMP SHARED(g_lo, z0tmp, dtmp) & $OMP SCHEDULE(static) !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! du == int (E nu_s f_0); du = du(z, kx, ky, s) duinv = 1/du

Read more…

Arguments

None

private subroutine init_diffuse_conserve()

Precompute three quantities needed for momentum and energy conservation: bz0, bw0, bs0 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! Now get z0 (first form) $OMP PARALLEL DO DEFAULT(none) & $OMP PRIVATE(iglo, ik, it, ie, il, is, isgn, ig) & $OMP SHARED(g_lo, ntgrid, forbid, bz0tmp, code_dt, vnmult, vnew_E, & $OMP aj0, duinv, conserve_forbid_zero) & $OMP SCHEDULE(static) !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! Now get v0z0 $OMP PARALLEL DO DEFAULT(none) & $OMP PRIVATE(iglo, ik, ie, is, isgn) & $OMP SHARED(g_lo, gtmp, vnmult, vnew_E, aj0, bz0tmp) & $OMP SCHEDULE(static) !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! Redefine z0 = z0 / (1 + v0z0) $OMP PARALLEL DO DEFAULT(none) & $OMP PRIVATE(iglo, ik, it, is, isgn) & $OMP SHARED(g_lo, bz0tmp, dtmp) & $OMP SCHEDULE(static) !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! redefine dq = du (for momentum-conserving terms) du == int (E nu_s f_0); du = du(z, kx, ky, s) duinv = 1/du $OMP PARALLEL DO DEFAULT(none) & $OMP PRIVATE(iglo, ik, ie, is, isgn) & $OMP SHARED(g_lo, gtmp, vns, vpa) & $OMP SCHEDULE(static) !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! Now get s0 (first form) $OMP PARALLEL DO DEFAULT(none) & $OMP PRIVATE(iglo, ik, it, ie, is, isgn) & $OMP SHARED(g_lo, bs0tmp, vns, vpa, aj0, code_dt, duinv) & $OMP SCHEDULE(static) !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! Now get v0s0 $OMP PARALLEL DO DEFAULT(none) & $OMP PRIVATE(iglo, ik, ie, is, isgn) & $OMP SHARED(g_lo, gtmp, vnmult, vnew_E, aj0, bs0tmp) & $OMP SCHEDULE(static) !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! Redefine s0 = s0 - v0s0 * z0 / (1 + v0z0) $OMP PARALLEL DO DEFAULT(none) & $OMP PRIVATE(iglo, ik, it, is, isgn) & $OMP SHARED(g_lo, bs0tmp, dtmp, bz0tmp) & $OMP SCHEDULE(static) !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! Now get v1s0 $OMP PARALLEL DO DEFAULT(none) & $OMP PRIVATE(iglo, ik, ie, is, isgn) & $OMP SHARED(g_lo, gtmp, vns, vpa, aj0, bs0tmp) & $OMP SCHEDULE(static) !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! Redefine s0 = s0 / (1 + v0s0) $OMP PARALLEL DO DEFAULT(none) & $OMP PRIVATE(iglo, ik, it, is, isgn) & $OMP SHARED(g_lo, bs0tmp, dtmp) & $OMP SCHEDULE(static) !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! Now get w0 $OMP PARALLEL DO DEFAULT(none) & $OMP PRIVATE(iglo, ik, it, ie, il, is, isgn, ig) & $OMP SHARED(g_lo, ntgrid, forbid, bw0tmp, vns, energy, al, aj1, code_dt, & $OMP spec, kperp2, duinv, bmag, conserve_forbid_zero) & $OMP SCHEDULE(static) !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! Now get v0w0 $OMP PARALLEL DO DEFAULT(none) & $OMP PRIVATE(iglo, ik, ie, is, isgn) & $OMP SHARED(g_lo, gtmp, vnmult, vnew_E, aj0, bw0tmp) & $OMP SCHEDULE(static) !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! Redefine w0 = w0 - v0w0 * z0 / (1 + v0z0) (this is w1 from MAB notes) $OMP PARALLEL DO DEFAULT(none) & $OMP PRIVATE(iglo, ik, it, is, isgn) & $OMP SHARED(g_lo, bw0tmp, bz0tmp, dtmp) & $OMP SCHEDULE(static) !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! Get v1w1 $OMP PARALLEL DO DEFAULT(none) & $OMP PRIVATE(iglo, ik, ie, is, isgn) & $OMP SHARED(g_lo, gtmp, vns, vpa, aj0, bw0tmp) & $OMP SCHEDULE(static) !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! Redefine w0 = w1 - v1w1 * s1 / (1 + v1s1) (this is w2 from MAB notes) $OMP PARALLEL DO DEFAULT(none) & $OMP PRIVATE(iglo, ik, it, is, isgn) & $OMP SHARED(g_lo, bw0tmp, bs0tmp, dtmp) & $OMP SCHEDULE(static) !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! Get v2w2 $OMP PARALLEL DO DEFAULT(none) & $OMP PRIVATE(iglo, ik, ie, il, is, isgn) & $OMP SHARED(g_lo, gtmp, vns, energy, al, aj1, bw0tmp) & $OMP SCHEDULE(static) !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! Redefine w0 = w2 / (1 + v2w2) $OMP PARALLEL DO DEFAULT(none) & $OMP PRIVATE(iglo, ik, it, is, isgn) & $OMP SHARED(g_lo, bw0tmp, dtmp) & $OMP SCHEDULE(static)

Arguments

None

private subroutine zero_out_passing_hybrid_electrons(arr)

If we have hybrid electrons then we need to remove the contribution to the conservation terms from the adiabatic passing part. This routine does this for us.

Arguments

Type IntentOptional Attributes Name
complex, intent(inout), dimension(:, :, g_lo%llim_proc:) :: arr

private subroutine init_vnew()

FIXME : Add documentation MAB MAB

Arguments

None

public subroutine adjust_vnmult(errest, consider_trapped_error)

Given estimates of the velocity space integration errors adjust the collision frequency scaling factor vnmult.

Arguments

Type IntentOptional Attributes Name
real, intent(in), dimension(5, 2) :: errest
logical, intent(in) :: consider_trapped_error

private subroutine init_ediffuse(vnmult_target)

FIXME : Add documentation

Arguments

Type IntentOptional Attributes Name
real, intent(in), optional :: vnmult_target

private subroutine get_ediffuse_matrix(aa, bb, cc, ig, ik, it, il, is, xe)

FIXME : Add documentation

Arguments

Type IntentOptional Attributes Name
real, intent(out), dimension (:) :: aa
real, intent(out), dimension (:) :: bb
real, intent(out), dimension (:) :: cc
integer, intent(in) :: ig
integer, intent(in) :: ik
integer, intent(in) :: it
integer, intent(in) :: il
integer, intent(in) :: is
real, intent(in), dimension (:) :: xe

private subroutine init_lorentz(vnmult_target)

FIXME : Add documentation

Arguments

Type IntentOptional Attributes Name
real, intent(in), optional :: vnmult_target

private subroutine set_hzero_lorentz_collisions_matrix(aa, bb, cc, je, ik, is)

Special behaviour when h=0 for passing non-zonal electrons

Read more…

Arguments

Type IntentOptional Attributes Name
real, intent(inout), dimension(:) :: aa
real, intent(inout), dimension(:) :: bb
real, intent(inout), dimension(:) :: cc
integer, intent(in) :: je
integer, intent(in) :: ik
integer, intent(in) :: is

private subroutine get_lorentz_matrix(aa, bb, cc, dd, hh, ig, ik, it, ie, is)

FIXME : Add documentation

Arguments

Type IntentOptional Attributes Name
real, intent(out), dimension (:) :: aa
real, intent(out), dimension (:) :: bb
real, intent(out), dimension (:) :: cc
real, intent(out), dimension (:) :: dd
real, intent(out), dimension (:) :: hh
integer, intent(in) :: ig
integer, intent(in) :: ik
integer, intent(in) :: it
integer, intent(in) :: ie
integer, intent(in) :: is

public subroutine init_lorentz_error()

FIXME : Add documentation

Read more…

Arguments

None

private subroutine solfp1_standard_layout(g, g1, gc1, gc2, diagnostics, gtoc, ctog)

FIXME : Add documentation

Arguments

Type IntentOptional Attributes Name
complex, intent(inout), dimension (-ntgrid:,:,g_lo%llim_proc:) :: g
complex, intent(inout), dimension (-ntgrid:,:,g_lo%llim_proc:) :: g1
complex, intent(inout), dimension (-ntgrid:,:,g_lo%llim_proc:) :: gc1
complex, intent(inout), dimension (-ntgrid:,:,g_lo%llim_proc:) :: gc2
integer, intent(in), optional :: diagnostics
logical, intent(in), optional :: gtoc
logical, intent(in), optional :: ctog

private subroutine solfp1_le_layout(gle, diagnostics)

FIXME : Add documentation

Arguments

Type IntentOptional Attributes Name
complex, intent(inout), dimension (:,:,le_lo%llim_proc:) :: gle
integer, intent(in), optional :: diagnostics

private subroutine conserve_lorentz_standard_layout(g, g1)

FIXME : Add documentation !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! First get v0y0 $OMP PARALLEL DO DEFAULT(none) & $OMP PRIVATE(iglo, it, ik, isgn) & $OMP SHARED(g_lo, kwork_filter, gtmp, vpa, aj0, g) & $OMP SCHEDULE(static) !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! Get y1 = y0 - v0y0 * z0 / (1 + v0z0) $OMP PARALLEL DO DEFAULT(none) & $OMP PRIVATE(iglo, it, ik, is, isgn) & $OMP SHARED(g_lo, kwork_filter, g1, g, ieqzip, v0y0, z0) & $OMP SCHEDULE(static) !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! Now get v1y1 $OMP PARALLEL DO DEFAULT(none) & $OMP PRIVATE(iglo, it, ik, il, ie, is, isgn) & $OMP SHARED(g_lo, kwork_filter, conservative, gtmp, vns, speed, vpdiff, aj0, g1, vpa) & $OMP SCHEDULE(static) !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! Get y2 = y1 - v1y1 * s1 / (1 + v1s1) $OMP PARALLEL DO DEFAULT(none) & $OMP PRIVATE(iglo, it, ik, is, isgn) & $OMP SHARED(g_lo, kwork_filter, g1, ieqzip, v1y1, s0) & $OMP SCHEDULE(static) !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! Now get v2y2 $OMP PARALLEL DO DEFAULT(none) & $OMP PRIVATE(iglo, it, ik, ie, il, is, isgn) & $OMP SHARED(g_lo, kwork_filter, gtmp, vns, energy, al, aj1, g1) & $OMP SCHEDULE(static) !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! Finally get x = y2 - v2y2 * w2 / (1 + v2w2) $OMP PARALLEL DO DEFAULT(none) & $OMP PRIVATE(iglo, it, ik, is, isgn) & $OMP SHARED(g_lo, kwork_filter, g, g1, ieqzip, v2y2, w0) & $OMP SCHEDULE(static)

Arguments

Type IntentOptional Attributes Name
complex, intent(inout), dimension (-ntgrid:,:,g_lo%llim_proc:) :: g
complex, intent(inout), dimension (-ntgrid:,:,g_lo%llim_proc:) :: g1

private subroutine conserve_lorentz_le_layout(gle)

FIXME : Add documentation !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! First get v0y0 and y1 = y0 - v0y0 * z0 / (1 + v0z0) v0 = vpa J0 f0, y0 = gle $OMP PARALLEL DO DEFAULT(none) & $OMP PRIVATE(ile, it, ik, ixi, ie, is, ig, v0y0) & $OMP SHARED(le_lo, kwork_filter, ieqzip, negrid, nxi_lim, & $OMP z0le, w, wxi, vpa_aj0_le, gle) & $OMP SCHEDULE(static) !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! Now get v1y1 and y2 = y1 - v1y1 * s1 / (1 + v1s1) v1 = nud vpa J0 f0, y1 = gle

Read more…

Arguments

Type IntentOptional Attributes Name
complex, intent(inout), dimension (:,:,le_lo%llim_proc:) :: gle

private subroutine conserve_diffuse_standard_layout(g, g1)

FIXME : Add documentation !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! First get v0y0 $OMP PARALLEL DO DEFAULT(none) & $OMP PRIVATE(iglo, it, ik, ie, is, isgn) & $OMP SHARED(g_lo, kwork_filter, gtmp, vnmult, vnew_E, aj0, g) & $OMP SCHEDULE(static) !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! Get y1 = y0 - v0y0 * z0 / (1 + v0z0) $OMP PARALLEL DO DEFAULT(none) & $OMP PRIVATE(iglo, it, ik, is, isgn) & $OMP SHARED(g_lo, kwork_filter, g1, g, ieqzip, v0y0, bz0) & $OMP SCHEDULE(static) !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! Now get v1y1 $OMP PARALLEL DO DEFAULT(none) & $OMP PRIVATE(iglo, it, ik, ie, is, isgn) & $OMP SHARED(g_lo, kwork_filter, gtmp, vns, vpa, aj0, g1) & $OMP SCHEDULE(static) !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! Get y2 = y1 - v1y1 * s1 / (1 + v1s1) $OMP PARALLEL DO DEFAULT(none) & $OMP PRIVATE(iglo, it, ik, is, isgn) & $OMP SHARED(g_lo, kwork_filter, g1, ieqzip, v1y1, bs0) & $OMP SCHEDULE(static) !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! Now get v2y2 $OMP PARALLEL DO DEFAULT(none) & $OMP PRIVATE(iglo, it, ik, ie, il, is, isgn) & $OMP SHARED(g_lo, kwork_filter, gtmp, vns, energy, al, aj1, g1) & $OMP SCHEDULE(static) !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! Finally get x = y2 - v2y2 * w2 / (1 + v2w2) $OMP PARALLEL DO DEFAULT(none) & $OMP PRIVATE(iglo, it, ik, is, isgn) & $OMP SHARED(g_lo, kwork_filter, g, g1, ieqzip, v2y2, bw0) & $OMP SCHEDULE(static)

Arguments

Type IntentOptional Attributes Name
complex, intent(inout), dimension (-ntgrid:,:,g_lo%llim_proc:) :: g
complex, intent(inout), dimension (-ntgrid:,:,g_lo%llim_proc:) :: g1

private subroutine conserve_diffuse_le_layout(gle)

FIXME : Add documentation !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! First get v0y0 and then y1 = y0 - v0y0 * z0 / (1 + v0z0) $OMP PARALLEL DO DEFAULT(none) & $OMP PRIVATE(ile, it, ik, is, ie, ig, ixi, vnE, v0y0) & $OMP SHARED(le_lo, kwork_filter, negrid, nxi_lim, w, wxi, ieqzip, bz0le, vnmult, vnew_E, aj0le, gle) & $OMP SCHEDULE(static) !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! Now get v1y1 and then y2 = y1 - v1y1 * s1 / (1 + v1s1) $OMP PARALLEL DO DEFAULT(none) & $OMP PRIVATE(ile, it, ik, is, ig, delnu, ie, ixi, v1y1) & $OMP SHARED(le_lo, kwork_filter, negrid, vnmult, delvnew, nxi_lim, & $OMP bs0le, vpa_aj0_le, gle, w, wxi, ieqzip) & $OMP SCHEDULE(static) !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! Now get v2y2 and then get x = y2 - v2y2 * w2 / (1 + v2w2) $OMP PARALLEL DO DEFAULT(none) & $OMP PRIVATE(ile, it, ik, is, ie, ig, ixi, delnu, v2y2) & $OMP SHARED(le_lo, kwork_filter, negrid, nxi_lim, ieqzip, & $OMP bw0le, delvnew, vperp_aj1le, gle, vnmult, w, wxi) & $OMP SCHEDULE(static)

Arguments

Type IntentOptional Attributes Name
complex, intent(inout), dimension (:,:,le_lo%llim_proc:) :: gle

private subroutine solfp_lorentz_standard_layout(g, gc, gh, diagnostics)

FIXME : Add documentation

Arguments

Type IntentOptional Attributes Name
complex, intent(inout), dimension (-ntgrid:,:,g_lo%llim_proc:) :: g
complex, intent(inout), dimension (-ntgrid:,:,g_lo%llim_proc:) :: gc
complex, intent(inout), dimension (-ntgrid:,:,g_lo%llim_proc:) :: gh
integer, intent(in), optional :: diagnostics

private subroutine solfp_lorentz_le_layout(gle, diagnostics)

FIXME : Add documentation

Arguments

Type IntentOptional Attributes Name
complex, intent(inout), dimension (:,:,le_lo%llim_proc:) :: gle
integer, intent(in), optional :: diagnostics

private subroutine solfp_ediffuse_standard_layout(g)

Energy diffusion subroutine used with energy layout (not le_layout) this is always the case when initializing the conserving terms, otherwise is the case if use_le_layout is no specified in the input file.

Arguments

Type IntentOptional Attributes Name
complex, intent(inout), dimension (-ntgrid:,:,g_lo%llim_proc:) :: g

private subroutine solfp_ediffuse_le_layout(gle)

FIXME : Add documentation

Arguments

Type IntentOptional Attributes Name
complex, intent(inout), dimension (:,:,le_lo%llim_proc:) :: gle

private subroutine init_vpdiff()

FIXME : Add documentation

Arguments

None

public subroutine reset_init()

Forces recalculation of coefficients in collision operator when timestep changes. Currently just calls finish_collisions

Arguments

None

public subroutine finish_collisions()

Forces recalculation of coefficients in collision operator when timestep changes.

Arguments

None

public subroutine set_collisions_config(collisions_config_in)

Set the module level config type Will abort if the module has already been initialised to avoid inconsistencies.

Arguments

Type IntentOptional Attributes Name
type(collisions_config_type), intent(in), optional :: collisions_config_in

private subroutine read_collisions_config(self)

Reads in the collisions_knobs namelist and populates the member variables

Arguments

Type IntentOptional Attributes Name
class(collisions_config_type), intent(inout) :: self

private subroutine write_collisions_config(self, unit)

Writes out a namelist representing the current state of the config object

Arguments

Type IntentOptional Attributes Name
class(collisions_config_type), intent(in) :: self
integer, intent(in), optional :: unit

private subroutine reset_collisions_config(self)

Resets the config object to the initial empty state

Arguments

Type IntentOptional Attributes Name
class(collisions_config_type), intent(inout) :: self

private subroutine broadcast_collisions_config(self)

Broadcasts all config parameters so object is populated identically on all processors

Arguments

Type IntentOptional Attributes Name
class(collisions_config_type), intent(inout) :: self