ddeabm_with_event_initialize_vec Subroutine

private subroutine ddeabm_with_event_initialize_vec(me, neq, maxnum, df, rtol, atol, ng, g, root_tol, report, bracket, initial_step_mode, initial_step_size)

Initialize ddeabm_with_event_class_vec class, and set the variables that cannot be changed during a problem.

See also

Type Bound

ddeabm_with_event_class_vec

Arguments

Type IntentOptional Attributes Name
class(ddeabm_with_event_class_vec), intent(inout) :: me
integer, intent(in) :: neq

number of equations to be integrated

integer, intent(in) :: maxnum

maximum number of integration steps allowed

procedure(deriv_func) :: df

derivative function

real(kind=wp), intent(in), dimension(:) :: rtol

relative tolerance for integration (see ddeabm)

real(kind=wp), intent(in), dimension(:) :: atol

absolution tolerance for integration (see ddeabm)

integer, intent(in) :: ng

number of event functions in g

procedure(event_func_vec) :: g

Event function . This should be a smooth function than can have values and . When (within the tolerance), then a root has been located and the integration will stop.

real(kind=wp), intent(in), dimension(:) :: root_tol

tolerance for the root finding. This should be sized ng or 1 (in which case the value is used for all elements)

procedure(report_func), optional :: report

reporting function

procedure(bracket_func_vec), optional :: bracket

root bracketing function. if not present, the default is used.

integer, intent(in), optional :: initial_step_mode

how to choose the initial step h:

  1. Use dhstrt.
  2. Use the older (quicker) algorithm.
  3. Use the user-specified value initial_step_size (>0).
real(kind=wp), intent(in), optional :: initial_step_size

for initial_step_mode=3


Calls

proc~~ddeabm_with_event_initialize_vec~~CallsGraph proc~ddeabm_with_event_initialize_vec ddeabm_with_event_class_vec%ddeabm_with_event_initialize_vec proc~ddeabm_initialize ddeabm_class%ddeabm_initialize proc~ddeabm_with_event_initialize_vec->proc~ddeabm_initialize proc~destroy_ddeabm ddeabm_class%destroy_ddeabm proc~ddeabm_initialize->proc~destroy_ddeabm

Source Code

   subroutine ddeabm_with_event_initialize_vec(me, neq, maxnum, df, rtol, atol, &
                                               ng, g, root_tol, report, bracket, &
                                               initial_step_mode, initial_step_size)

      implicit none

      class(ddeabm_with_event_class_vec), intent(inout) :: me
      integer, intent(in)                           :: neq       !! number of equations to be integrated
      integer, intent(in)                           :: maxnum    !! maximum number of integration steps allowed
      procedure(deriv_func)                        :: df        !! derivative function \( dx/dt \)
      real(wp), dimension(:), intent(in)             :: rtol      !! relative tolerance for integration (see [[ddeabm]])
      real(wp), dimension(:), intent(in)             :: atol      !! absolution tolerance for integration (see [[ddeabm]])
      integer, intent(in)                           :: ng        !! number of event functions in `g`
      procedure(event_func_vec)                    :: g         !! Event function \( g(t,x) \). This should be a smooth function
                                                              !! than can have values \( <0 \) and \( \ge 0 \).
                                                              !! When \( g=0 \) (within the tolerance),
                                                              !! then a root has been located and
                                                              !! the integration will stop.
      real(wp), dimension(:), intent(in)             :: root_tol  !! tolerance for the root finding.
                                                              !! This should be sized `ng` or `1` (in which
                                                              !! case the value is used for all elements)
      procedure(report_func), optional              :: report    !! reporting function
      procedure(bracket_func_vec), optional         :: bracket   !! root bracketing function. if not present,
                                                              !! the default is used.
      integer, intent(in), optional         :: initial_step_mode !! how to choose the initial step `h`:
                                                             !!
                                                             !! 1. Use [[dhstrt]].
                                                             !! 2. Use the older (quicker) algorithm.
                                                             !! 3. Use the user-specified value `initial_step_size` (>0).
      real(wp), intent(in), optional        :: initial_step_size !! for `initial_step_mode=3`

      !base class initialization:
      call me%initialize(neq, maxnum, df, rtol, atol, report, initial_step_mode, initial_step_size)

      !event finding variables:
      me%n_g_eqns = ng
      allocate (me%tol(ng))
      if (size(root_tol) == 1) then
         me%tol = root_tol(1) ! use this one for all of them
      elseif (size(root_tol) == ng) then
         me%tol = root_tol
      else
         write (*, *) 'Error in ddeabm_with_event_initialize_vec: '// &
            'incorrect size for root_tol input. Using first element only.'
         me%tol = root_tol(1)
      end if
      me%gfunc => g

      ! bracketing function
      ! (if not specified, use the default):
      if (present(bracket)) then
         me%bracket => bracket
      else
         me%bracket => null()
      end if

   end subroutine ddeabm_with_event_initialize_vec