split_nonlinear_terms Module

A module to deal with advancing the nonlinear term separately from the linear terms.



Contents


Variables

Type Visibility Attributes Name Initial
logical, private :: initialized = .false.
integer, private :: nfailed_steps = 0

Statistics for the number of steps taken in the integrators

integer, private :: nrhs_calls = 0

Statistics for the number of steps taken in the integrators

integer, private :: nsuccessful_steps = 0

Statistics for the number of steps taken in the integrators

integer, private :: split_method_switch
integer, private, parameter :: split_method_ab3 = 1
integer, private, parameter :: split_method_backwards_euler = 2
integer, private, parameter :: split_method_rk = 3
integer, private, parameter :: split_method_picard = 4
type(rk_scheme_type), private :: the_rk_scheme
logical, private :: show_statistics
logical, public :: advance_nonadiabatic_dfn
logical, public :: strang_split
real, private :: convergence_tolerance
real, private :: time_step_safety_factor
real, private :: relative_tolerance
real, private :: absolute_tolerance
type(split_nonlinear_terms_config_type), private :: split_nonlinear_terms_config

Interfaces

interface

Interface to describe the call signatures for routines to calculate the field from gnew (invert_field_func) and to calculate the source term used in the advance_nonlinear_term methods. Used to allow us to pass in these methods, giving us a way to test the integrators on different problems.

  • private subroutine invert_field_func(g_in, phi, apar, bpar, gf_lo, local_only)

    Arguments

    Type IntentOptional Attributes Name
    complex, intent(in), dimension (:,:,:) :: g_in
    complex, intent(out), dimension (:,:,:) :: phi
    complex, intent(out), dimension (:,:,:) :: apar
    complex, intent(out), dimension (:,:,:) :: bpar
    logical, intent(in), optional :: gf_lo
    logical, intent(in), optional :: local_only

interface

Interface to describe the call signatures for routines to calculate the field from gnew (invert_field_func) and to calculate the source term used in the advance_nonlinear_term methods. Used to allow us to pass in these methods, giving us a way to test the integrators on different problems.

  • private subroutine source_term_func(g_in, g1, phi, apar, bpar, max_vel_local, need_to_adjust, calculate_cfl_limit)

    Arguments

    Type IntentOptional Attributes Name
    complex, intent(in), dimension (:,:,:) :: g_in
    complex, intent(out), dimension (:,:,:) :: g1
    complex, intent(in), dimension (:,:,:) :: phi
    complex, intent(in), dimension (:,:,:) :: apar
    complex, intent(in), dimension (:,:,:) :: bpar
    real, intent(out) :: max_vel_local
    logical, intent(in) :: need_to_adjust
    logical, intent(in), optional :: calculate_cfl_limit

Derived Types

Used to represent the input configuration of split_nonlinear_terms

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 :: skip_smart_defaults = .false.

Do we want to skip the smaart defaults in init?

real, public :: absolute_tolerance = 1.0e-3

The absolute tolerance used in error controlled schemes. Attempt to adjust time step to keep error below this tolerance.

logical, public :: advance_nonadiabatic_dfn = .false.

If true then nonlinear integrator advances the nonadiabatic dfn, h, rather than the modified distribution function, g.

real, public :: convergence_tolerance = 1.0e-1

The tolerance used in convergence checks. Currently only used for halting the fixed-point iteration in the beuler method

real, public :: relative_tolerance = 1.0e-2

The relative tolerance used in error controlled schemes. Attempt to adjust time step to keep error below this tolerance.

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

Choose which specific RK scheme to use in the RK method. Generally if the tolerances are small a higher order scheme will be more effective than low order. At loose tolerances low order schemes are generally more efficient Valid options include:

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

If true then reports the number of right hand side (NL term) calculations, the number of failed internal steps and the number of successful steps at the end of each linear size step.

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

What algorithm do we use to advance the nonlinear term if split_nonlinear is true. Valid options include:

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

If true then indicates that we want to use a strang split approach where we advance the NL operator by a half step, linear by a full step and then NL by another half step. This should be second order accurate whilst the regular split is only first order accurate.

real, public :: time_step_safety_factor = 0.9

Multiplies the new time step calculated from either the error or cfl limits. Smaller values are more conservative but may help avoid repeated violation of the error/cfl limits and hence reduce the number of failed steps.

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 , :: set_smart_defaults => set_smart_defaults_null Subroutine
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_split_nonlinear_terms_config Subroutine
procedure , public :: write => write_split_nonlinear_terms_config Subroutine
procedure , public :: reset => reset_split_nonlinear_terms_config Subroutine
procedure , public :: broadcast => broadcast_split_nonlinear_terms_config Subroutine
procedure , public , nopass :: get_default_name => get_default_name_split_nonlinear_terms_config Function
procedure , public , nopass :: get_default_requires_index => get_default_requires_index_split_nonlinear_terms_config Function

Functions

private pure function get_inverse_error(errors) result(inverse_error)

Returns the normalised inverse error estimate. In other words the error tolerance, rtol*errors(2) + atol, divided by the error estimate (plus a small number to avoid divide by zero).

Read more…

Arguments

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

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_split_nonlinear_terms_config()

Get the module level config instance

Arguments

None

Return Value type(split_nonlinear_terms_config_type)


Subroutines

public subroutine init_split_nonlinear_terms(split_nonlinear_terms_config_in)

Initialises the split_nonlinear_terms module. Primarily just reading input

Arguments

Type IntentOptional Attributes Name
type(split_nonlinear_terms_config_type), intent(in), optional :: split_nonlinear_terms_config_in

private subroutine read_parameters(split_nonlinear_terms_config_in)

Read the split_nonlinear_terms namelist

Arguments

Type IntentOptional Attributes Name
type(split_nonlinear_terms_config_type), intent(in), optional :: split_nonlinear_terms_config_in

public subroutine finish_split_nonlinear_terms()

Reset the module, freeing memory etc.

Arguments

None

private subroutine reset_step_statistics()

Reset the step statistics

Arguments

None

private subroutine report_step_statistics(istep, unit_in)

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: istep
integer, intent(in), optional :: unit_in

public subroutine advance_nonlinear_term(g_state, istep, phi, apar, bpar, dt, fields_func, source_func)

Advances dg/dt = NL(g,chi) from g_state, chi from t -> t+dt by calling specific requested method. On output g_state contains new g at next step, chi is left unchanged

Arguments

Type IntentOptional Attributes Name
complex, intent(inout), dimension (:,:,:) :: g_state
integer, intent(in) :: istep
complex, intent(in), dimension (:,:,:) :: phi
complex, intent(in), dimension (:,:,:) :: apar
complex, intent(in), dimension (:,:,:) :: bpar
real, intent(in) :: dt
procedure(invert_field_func) :: fields_func
procedure(source_term_func) :: source_func

private subroutine advance_nonlinear_term_ab3(g_state, phinew, aparnew, bparnew, dt, fields_func, source_func)

Advances dg/dt = NL(g,chi) from current g_state, chi from t -> t+dt by using AB scheme of orders up to 3rd.

Arguments

Type IntentOptional Attributes Name
complex, intent(inout), dimension (-ntgrid:, :, g_lo%llim_proc:) :: g_state
complex, intent(inout), dimension (-ntgrid:,:,:) :: phinew
complex, intent(inout), dimension (-ntgrid:,:,:) :: aparnew
complex, intent(inout), dimension (-ntgrid:,:,:) :: bparnew
real, intent(in) :: dt
procedure(invert_field_func) :: fields_func
procedure(source_term_func) :: source_func

private subroutine advance_nonlinear_term_beuler(g_state, phinew, aparnew, bparnew, dt, fields_func, source_func)

Advances dg/dt = NL(g,chi) from current g_state, chi from t -> t+dt by using backwards Euler with fixed-point iteration.

Arguments

Type IntentOptional Attributes Name
complex, intent(inout), dimension (-ntgrid:, :, g_lo%llim_proc:) :: g_state
complex, intent(inout), dimension (-ntgrid:,:,:) :: phinew
complex, intent(inout), dimension (-ntgrid:,:,:) :: aparnew
complex, intent(inout), dimension (-ntgrid:,:,:) :: bparnew
real, intent(in) :: dt
procedure(invert_field_func) :: fields_func
procedure(source_term_func) :: source_func

private subroutine advance_nonlinear_term_rk(g_state, phinew, aparnew, bparnew, dt, fields_func, source_func)

Advances dg/dt = NL(g,chi) from current g_state, chi from t -> t+dt by using an Runge-Kutta scheme with embedded error estimate for error control. This method is a wrapper to the true generic RK implementation.

Arguments

Type IntentOptional Attributes Name
complex, intent(inout), dimension (:,:,:) :: g_state
complex, intent(inout), dimension (:,:,:) :: phinew
complex, intent(inout), dimension (:,:,:) :: aparnew
complex, intent(inout), dimension (:,:,:) :: bparnew
real, intent(in) :: dt
procedure(invert_field_func) :: fields_func
procedure(source_term_func) :: source_func

private subroutine advance_nonlinear_term_rk_implementation(g_state, phinew, aparnew, bparnew, dt, fields_func, source_func, scheme)

Advances dg/dt = NL(g,chi) from current g_state, chi from t -> t+dt by using an Runge-Kutta scheme with embedded error estimate for error control

Arguments

Type IntentOptional Attributes Name
complex, intent(inout), dimension (-ntgrid:, :, g_lo%llim_proc:) :: g_state
complex, intent(inout), dimension (-ntgrid:,:,:) :: phinew
complex, intent(inout), dimension (-ntgrid:,:,:) :: aparnew
complex, intent(inout), dimension (-ntgrid:,:,:) :: bparnew
real, intent(in) :: dt
procedure(invert_field_func) :: fields_func
procedure(source_term_func) :: source_func
type(rk_scheme_type), intent(in) :: scheme

private subroutine advance_nonlinear_term_picard(g_state, phinew, aparnew, bparnew, dt, fields_func, source_func)

Advances dg/dt = NL(g,chi) from current g_state, chi from t -> t+dt by using a simple Picard iteration scheme

Arguments

Type IntentOptional Attributes Name
complex, intent(inout), dimension (-ntgrid:, :, g_lo%llim_proc:) :: g_state
complex, intent(inout), dimension (-ntgrid:,:,:) :: phinew
complex, intent(inout), dimension (-ntgrid:,:,:) :: aparnew
complex, intent(inout), dimension (-ntgrid:,:,:) :: bparnew
real, intent(in) :: dt
procedure(invert_field_func) :: fields_func
procedure(source_term_func) :: source_func

private subroutine calculate_fields(fields_func, g_in, phi_out, apar_out, bpar_out)

Calculates the fields consistent with the passed distribution function

Arguments

Type IntentOptional Attributes Name
procedure(invert_field_func) :: fields_func
complex, intent(in), dimension(:, :, :) :: g_in
complex, intent(out), dimension (:,:,:) :: phi_out
complex, intent(out), dimension (:,:,:) :: apar_out
complex, intent(out), dimension (:,:,:) :: bpar_out

private subroutine read_split_nonlinear_terms_config(self)

Reads in the split_nonlinear_terms_knobs namelist and populates the member variables

Arguments

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

private subroutine write_split_nonlinear_terms_config(self, unit)

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

Arguments

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

private subroutine reset_split_nonlinear_terms_config(self)

Resets the config object to the initial empty state

Arguments

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

private subroutine broadcast_split_nonlinear_terms_config(self)

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

Arguments

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

public subroutine set_split_nonlinear_terms_config(split_nonlinear_terms_config_in)

Set the module level config instance

Arguments

Type IntentOptional Attributes Name
type(split_nonlinear_terms_config_type), intent(in), optional :: split_nonlinear_terms_config_in