Create the actual Adams-Bashforth-Moulton integrator: initialize the b coefficients.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(integrator_adams_bashforth_moulton), | intent(inout) | :: | self | Integrator. |
||
character(len=*), | intent(in) | :: | scheme | Selected scheme. |
||
integer(kind=I_P), | intent(in), | optional | :: | iterations | Implicit iterations. |
|
logical, | intent(in), | optional | :: | autoupdate | Enable cyclic autoupdate of previous time steps. |
|
class(integrand_object), | intent(in), | optional | :: | U | Integrand molding prototype. |
|
logical, | intent(in), | optional | :: | stop_on_fail | Stop execution if initialization fail. |
subroutine initialize(self, scheme, iterations, autoupdate, U, stop_on_fail)
!< Create the actual Adams-Bashforth-Moulton integrator: initialize the *b* coefficients.
class(integrator_adams_bashforth_moulton), intent(inout) :: self !< Integrator.
character(*), intent(in) :: scheme !< Selected scheme.
integer(I_P), intent(in), optional :: iterations !< Implicit iterations.
logical, intent(in), optional :: autoupdate !< Enable cyclic autoupdate of previous time steps.
class(integrand_object), intent(in), optional :: U !< Integrand molding prototype.
logical, intent(in), optional :: stop_on_fail !< Stop execution if initialization fail.
character(len=99), allocatable :: schemes_ab(:) !< Adams-Bashforth schemes.
character(len=99), allocatable :: schemes_am(:) !< Adams-Moulton schemes.
integer(I_P) :: scheme_number_ !< Scheme number in the list of supported schemes.
if (self%is_supported(scheme=scheme)) then
call self%destroy
self%description_ = trim(adjustl(scheme))
scheme_number_ = self%scheme_number(scheme=scheme)
schemes_ab = self%predictor%supported_schemes()
schemes_am = self%corrector%supported_schemes()
self%autoupdate = .true. ; if (present(autoupdate)) self%autoupdate = autoupdate
self%iterations = 1 ; if (present(iterations)) self%iterations = iterations
call self%predictor%initialize(scheme=schemes_ab(scheme_number_), U=U, autoupdate=.false.)
call self%corrector%initialize(scheme=schemes_am(scheme_number_), U=U, iterations=self%iterations, autoupdate=.false.)
self%steps = self%predictor%steps_number()
self%registers = self%steps
if (present(U)) call self%allocate_integrand_members(U=U)
else
call self%trigger_error(error=ERROR_UNSUPPORTED_SCHEME, &
error_message='"'//trim(adjustl(scheme))//'" unsupported scheme', &
is_severe=stop_on_fail)
endif
endsubroutine initialize