initialize_rk_class Subroutine

private subroutine initialize_rk_class(me, n, f, report, g, stop_on_errors, max_number_of_steps, report_rate, solver)

Initialize the rk_class.

Type Bound

rk_class

Arguments

Type IntentOptional Attributes Name
class(rk_class), intent(inout) :: me
integer, intent(in) :: n

number of variables

procedure(deriv_func) :: f

derivative function

procedure(report_func), optional :: report

for reporting the steps

procedure(event_func), optional :: g

for stopping at an event

logical, intent(in), optional :: stop_on_errors

stop the program for any errors (default is False)

integer, intent(in), optional :: max_number_of_steps

max number of steps allowed

integer, intent(in), optional :: report_rate

how often to call report function: 0 : no reporting (same as not associating report), 1 : report every point, 2 : report every other point, etc. The first and last point are always reported.

class(root_method), intent(in), optional :: solver

the root-finding method to use for even finding. if not present, then brent_solver is used.


Calls

proc~~initialize_rk_class~~CallsGraph proc~initialize_rk_class rklib_module::rk_class%initialize_rk_class proc~destroy rklib_module::rk_class%destroy proc~initialize_rk_class->proc~destroy properties properties proc~initialize_rk_class->properties

Source Code

    subroutine initialize_rk_class(me,n,f,report,g,stop_on_errors,&
                                   max_number_of_steps,report_rate,&
                                   solver)

    implicit none

    class(rk_class),intent(inout)   :: me
    integer,intent(in)              :: n       !! number of variables
    procedure(deriv_func)           :: f       !! derivative function
    procedure(event_func),optional  :: g       !! for stopping at an event
    procedure(report_func),optional :: report  !! for reporting the steps
    logical,intent(in),optional     :: stop_on_errors !! stop the program for
                                                      !! any errors (default is False)
    integer,intent(in),optional     :: max_number_of_steps !! max number of steps allowed
    integer,intent(in),optional     :: report_rate !! how often to call report function:
                                                   !! `0` : no reporting (same as not associating `report`),
                                                   !! `1` : report every point,
                                                   !! `2` : report every other point, etc.
                                                   !! The first and last point are always reported.
    class(root_method),intent(in),optional :: solver !! the root-finding method to use for even finding.
                                                     !! if not present, then `brent_solver` is used.

    type(rklib_properties) :: props !! to get the method properties

    call me%destroy()

    me%n = n
    me%f => f
    if (present(report)) me%report => report
    if (present(g))      me%g      => g
    if (present(stop_on_errors)) me%stop_on_errors = stop_on_errors
    if (present(max_number_of_steps)) me%max_number_of_steps = abs(max_number_of_steps)
    if (present(report_rate)) me%report_rate = abs(report_rate)
    if (present(solver)) me%solver = solver

    ! allocate the registers:
    props = me%properties()
    if (allocated(me%funcs)) deallocate(me%funcs)
    allocate(me%funcs(n, props%number_of_registers))
    me%funcs = zero

    ! reset internal variables:
    me%num_steps = 0
    me%stopped = .false.

    end subroutine initialize_rk_class