Create the actual RK integrator: initialize the Butcher' table coefficients.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(integrator_runge_kutta_ssp), | intent(inout) | :: | self | Integrator. |
||
character(len=*), | intent(in) | :: | scheme | Selected scheme. |
||
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, U, stop_on_fail)
!< Create the actual RK integrator: initialize the Butcher' table coefficients.
class(integrator_runge_kutta_ssp), intent(inout) :: self !< Integrator.
character(*), intent(in) :: scheme !< Selected scheme.
class(integrand_object), intent(in), optional :: U !< Integrand molding prototype.
logical, intent(in), optional :: stop_on_fail !< Stop execution if initialization fail.
if (self%is_supported(scheme=scheme)) then
call self%destroy
self%description_ = trim(adjustl(scheme))
select case(trim(adjustl(scheme)))
case('runge_kutta_ssp_stages_1_order_1')
self%stages = 1
allocate(self%beta(1:self%stages )) ; self%beta = 0._R_P
allocate(self%alph(1:self%stages, 1:self%stages)) ; self%alph = 0._R_P
allocate(self%gamm( 1:self%stages)) ; self%gamm = 0._R_P
self%beta(1) = 1._R_P
case('runge_kutta_ssp_stages_2_order_2')
self%stages = 2
allocate(self%beta(1:self%stages )) ; self%beta = 0._R_P
allocate(self%alph(1:self%stages, 1:self%stages)) ; self%alph = 0._R_P
allocate(self%gamm( 1:self%stages)) ; self%gamm = 0._R_P
self%beta(1) = 0.5_R_P
self%beta(2) = 0.5_R_P
self%alph(2, 1) = 1._R_P
self%gamm(2) = 1._R_P
case('runge_kutta_ssp_stages_3_order_3')
self%stages = 3
allocate(self%beta(1:self%stages )) ; self%beta = 0._R_P
allocate(self%alph(1:self%stages, 1:self%stages)) ; self%alph = 0._R_P
allocate(self%gamm( 1:self%stages)) ; self%gamm = 0._R_P
self%beta(1) = 1._R_P/6._R_P
self%beta(2) = 1._R_P/6._R_P
self%beta(3) = 2._R_P/3._R_P
self%alph(2, 1) = 1._R_P
self%alph(3, 1) = 0.25_R_P ; self%alph(3, 2) = 0.25_R_P
self%gamm(2) = 1._R_P
self%gamm(3) = 0.5_R_P
case('runge_kutta_ssp_stages_5_order_4')
self%stages = 5
allocate(self%beta(1:self%stages )) ; self%beta = 0._R_P
allocate(self%alph(1:self%stages, 1:self%stages)) ; self%alph = 0._R_P
allocate(self%gamm( 1:self%stages)) ; self%gamm = 0._R_P
self%beta(1) = 0.14681187618661_R_P
self%beta(2) = 0.24848290924556_R_P
self%beta(3) = 0.10425883036650_R_P
self%beta(4) = 0.27443890091960_R_P
self%beta(5) = 0.22600748319395_R_P
self%alph(2, 1)=0.39175222700392_R_P
self%alph(3, 1)=0.21766909633821_R_P;self%alph(3, 2)=0.36841059262959_R_P
self%alph(4, 1)=0.08269208670950_R_P;self%alph(4, 2)=0.13995850206999_R_P;self%alph(4, 3)=0.25189177424738_R_P
self%alph(5, 1)=0.06796628370320_R_P;self%alph(5, 2)=0.11503469844438_R_P;self%alph(5, 3)=0.20703489864929_R_P
self%alph(5, 4)=0.54497475021237_R_P
self%gamm(2) = 0.39175222700392_R_P
self%gamm(3) = 0.58607968896780_R_P
self%gamm(4) = 0.47454236302687_R_P
self%gamm(5) = 0.93501063100924_R_P
endselect
self%registers = self%stages
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