Set the optional inputs for dop853.
In the original code, these were part of the work
and iwork
arrays.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(dop853_class), | intent(inout) | :: | me | |||
integer, | intent(in) | :: | n |
the dimension of the system (size of and vectors) |
||
procedure(deriv_func) | :: | fcn |
subroutine computing the value of |
|||
procedure(solout_func), | optional | :: | solout |
subroutine providing the
numerical solution during integration.
if |
||
integer, | intent(in), | optional | :: | iprint |
switch for printing error messages
if |
|
integer, | intent(in), | optional | :: | nstiff |
test for stiffness is activated after step number
|
|
integer, | intent(in), | optional | :: | nmax |
the maximal number of allowed steps. |
|
real(kind=wp), | intent(in), | optional | :: | hinitial |
initial step size, for |
|
real(kind=wp), | intent(in), | optional | :: | hmax |
maximal step size, defaults to |
|
real(kind=wp), | intent(in), | optional | :: | safe |
safety factor in step size prediction |
|
real(kind=wp), | intent(in), | optional | :: | fac1 |
parameter for step size selection.
the new step size is chosen subject to the restriction
|
|
real(kind=wp), | intent(in), | optional | :: | fac2 |
parameter for step size selection.
the new step size is chosen subject to the restriction
|
|
real(kind=wp), | intent(in), | optional | :: | beta |
is the |
|
integer, | intent(in), | optional, | dimension(:) | :: | icomp |
the components for which dense output is required (size from 0 to |
logical, | intent(out) | :: | status_ok |
will be false for invalid inputs. |
subroutine set_parameters(me,n,fcn,solout,iprint,nstiff,nmax,hinitial,&
hmax,safe,fac1,fac2,beta,icomp,status_ok)
implicit none
class(dop853_class),intent(inout) :: me
integer,intent(in) :: n !! the dimension of the system (size of \(y\) and \(y'\) vectors)
procedure(deriv_func) :: fcn !! subroutine computing the value of \( y' = f(x,y) \)
procedure(solout_func),optional :: solout !! subroutine providing the
!! numerical solution during integration.
!! if `iout>=1`, it is called during integration.
!! supply a dummy subroutine if `iout=0`.
integer,intent(in),optional :: iprint !! switch for printing error messages
!! if `iprint==0` no messages are being printed
!! if `iprint/=0` messages are printed with
!! `write (iprint,*)` ...
integer,intent(in),optional :: nstiff !! test for stiffness is activated after step number
!! `j*nstiff` (`j` integer), provided `nstiff>0`.
!! for negative `nstiff` the stiffness test is
!! never activated.
integer,intent(in),optional :: nmax !! the maximal number of allowed steps.
real(wp),intent(in),optional :: hinitial !! initial step size, for `hinitial=0` an initial guess
!! is computed with help of the function [[hinit]].
real(wp),intent(in),optional :: hmax !! maximal step size, defaults to `xend-x` if `hmax=0`.
real(wp),intent(in),optional :: safe !! safety factor in step size prediction
real(wp),intent(in),optional :: fac1 !! parameter for step size selection.
!! the new step size is chosen subject to the restriction
!! `fac1 <= hnew/hold <= fac2`
real(wp),intent(in),optional :: fac2 !! parameter for step size selection.
!! the new step size is chosen subject to the restriction
!! `fac1 <= hnew/hold <= fac2`
real(wp),intent(in),optional :: beta !! is the `beta` for stabilized step size control
!! (see section iv.2). positive values of `beta` ( <= 0.04 )
!! make the step size control more stable.
integer,dimension(:),intent(in),optional :: icomp !! the components for which dense output is required (size from 0 to `n`).
logical,intent(out) :: status_ok !! will be false for invalid inputs.
call me%destroy()
status_ok = .true.
!required inputs:
me%n = n
me%fcn => fcn
!optional inputs:
if (present(solout)) me%solout => solout
if (present(iprint)) me%iprint = iprint
if (present(nstiff)) me%nstiff = nstiff
if (present(hinitial)) me%hinitial = hinitial
if (present(hmax)) me%hmax = hmax
if (present(fac1)) me%fac1 = fac1
if (present(fac2)) me%fac2 = fac2
if (present(nmax)) then
if ( nmax<=0 ) then
if ( me%iprint/=0 ) &
write (me%iprint,*) ' wrong input nmax=', nmax
status_ok = .false.
else
me%nmax = nmax
end if
end if
if (present(safe)) then
if ( safe>=1.0_wp .or. safe<=1.0e-4_wp ) then
if ( me%iprint/=0 ) &
write (me%iprint,*) ' curious input for safety factor safe:', &
safe
status_ok = .false.
else
me%safe = safe
end if
end if
if (present(beta)) then
if ( beta<=0.0_wp ) then
me%beta = 0.0_wp
else
if ( beta>0.2_wp ) then
if ( me%iprint/=0 ) write (me%iprint,*) &
' curious input for beta: ', beta
status_ok = .false.
else
me%beta = beta
end if
end if
end if
if (present(icomp)) then
me%nrdens = size(icomp)
!check validity of icomp array:
if (size(icomp)<=me%n .and. all(icomp>0 .and. icomp<=me%n)) then
allocate(me%icomp(me%nrdens)); me%icomp = icomp
allocate(me%cont(8*me%nrdens)); me%cont = 0.0_wp
else
if ( me%iprint/=0 ) write (me%iprint,*) &
' invalid icomp array: ',icomp
status_ok = .false.
end if
end if
end subroutine set_parameters