A module to deal with advancing the nonlinear term separately from the linear terms.
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 |
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.
Type | Intent | Optional | 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 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.
Type | Intent | Optional | 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 |
Used to represent the input configuration of split_nonlinear_terms
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: |
|
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: |
|
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. |
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).
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
real, | intent(in), | dimension(:) | :: | errors |
Gets the default name for this namelist
Gets the default requires index for this namelist
Get the module level config instance
Initialises the split_nonlinear_terms module. Primarily just reading input
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
type(split_nonlinear_terms_config_type), | intent(in), | optional | :: | split_nonlinear_terms_config_in |
Read the split_nonlinear_terms namelist
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
type(split_nonlinear_terms_config_type), | intent(in), | optional | :: | split_nonlinear_terms_config_in |
Reset the module, freeing memory etc.
Reset the step statistics
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
integer, | intent(in) | :: | istep | |||
integer, | intent(in), | optional | :: | unit_in |
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
Type | Intent | Optional | 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 |
Advances dg/dt = NL(g,chi) from current g_state, chi from t -> t+dt by using AB scheme of orders up to 3rd.
Type | Intent | Optional | 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 |
Advances dg/dt = NL(g,chi) from current g_state, chi from t -> t+dt by using backwards Euler with fixed-point iteration.
Type | Intent | Optional | 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 |
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.
Type | Intent | Optional | 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 |
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
Type | Intent | Optional | 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 |
Advances dg/dt = NL(g,chi) from current g_state, chi from t -> t+dt by using a simple Picard iteration scheme
Type | Intent | Optional | 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 |
Calculates the fields consistent with the passed distribution function
Type | Intent | Optional | 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 |
Reads in the split_nonlinear_terms_knobs namelist and populates the member variables
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(split_nonlinear_terms_config_type), | intent(inout) | :: | self |
Writes out a namelist representing the current state of the config object
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(split_nonlinear_terms_config_type), | intent(in) | :: | self | |||
integer, | intent(in), | optional | :: | unit |
Resets the config object to the initial empty state
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(split_nonlinear_terms_config_type), | intent(inout) | :: | self |
Broadcasts all config parameters so object is populated identically on all processors
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(split_nonlinear_terms_config_type), | intent(inout) | :: | self |
Set the module level config instance
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
type(split_nonlinear_terms_config_type), | intent(in), | optional | :: | split_nonlinear_terms_config_in |