rklib_example Program

Uses

  • program~~rklib_example~~UsesGraph program~rklib_example rklib_example iso_fortran_env iso_fortran_env program~rklib_example->iso_fortran_env module~rklib_module rklib_module program~rklib_example->module~rklib_module module~rklib_module->iso_fortran_env root_module root_module module~rklib_module->root_module

Calls

program~~rklib_example~~CallsGraph program~rklib_example rklib_example proc~initialize_variable_step rklib_module::rk_variable_step_class%initialize_variable_step program~rklib_example->proc~initialize_variable_step proc~integrate_variable_step rklib_module::rk_variable_step_class%integrate_variable_step program~rklib_example->proc~integrate_variable_step proc~rk_class_status rklib_module::rk_class%rk_class_status program~rklib_example->proc~rk_class_status init init proc~initialize_variable_step->init proc~raise_exception rklib_module::rk_class%raise_exception proc~initialize_variable_step->proc~raise_exception proc~begin_integration_rk_variable_step_class rklib_module::rk_variable_step_class%begin_integration_rk_variable_step_class proc~integrate_variable_step->proc~begin_integration_rk_variable_step_class proc~compute_initial_step rklib_module::rk_variable_step_class%compute_initial_step proc~integrate_variable_step->proc~compute_initial_step proc~compute_stepsize rklib_module::stepsize_class%compute_stepsize proc~integrate_variable_step->proc~compute_stepsize proc~export_point rklib_module::rk_class%export_point proc~integrate_variable_step->proc~export_point proc~order rklib_module::rk_variable_step_class%order proc~integrate_variable_step->proc~order proc~integrate_variable_step->proc~raise_exception step step proc~integrate_variable_step->step begin begin proc~begin_integration_rk_variable_step_class->begin proc~destroy_fsal_cache rklib_module::rk_variable_step_fsal_class%destroy_fsal_cache proc~begin_integration_rk_variable_step_class->proc~destroy_fsal_cache f f proc~compute_initial_step->f proc~hinit rklib_module::rk_variable_step_class%hinit proc~compute_initial_step->proc~hinit proc~hstart rklib_module::rk_variable_step_class%hstart proc~compute_initial_step->proc~hstart properties properties proc~order->properties proc~hinit->proc~order proc~hinit->f proc~hstart->proc~order proc~hstart->f

Variables

Type Attributes Name Initial
integer, parameter :: n = 2

dimension of the system

real(kind=wp), parameter :: tol = 1.0e-12_wp

integration tolerance

real(kind=wp), parameter :: t0 = 0.0_wp

initial t value

real(kind=wp), parameter :: dt = 1.0_wp

initial step size

real(kind=wp), parameter :: tf = 100.0_wp

endpoint of integration

real(kind=wp), parameter, dimension(n) :: x0 = [0.0_wp, 0.1_wp]

initial x value

real(kind=wp), dimension(n) :: xf

final x value

type(rktp86_class) :: prop
character(len=:), allocatable :: message

Subroutines

subroutine fvpol(me, t, x, f)

Right-hand side of van der Pol equation

Arguments

Type IntentOptional Attributes Name
class(rk_class), intent(inout) :: me
real(kind=wp), intent(in) :: t
real(kind=wp), intent(in), dimension(:) :: x
real(kind=wp), intent(out), dimension(:) :: f

Source Code

program rklib_example

  use rklib_module, wp => rk_module_rk
  use iso_fortran_env, only: output_unit

  implicit none

  integer,parameter :: n = 2 !! dimension of the system
  real(wp),parameter :: tol = 1.0e-12_wp !! integration tolerance
  real(wp),parameter :: t0 = 0.0_wp !! initial t value
  real(wp),parameter :: dt = 1.0_wp !! initial step size
  real(wp),parameter :: tf = 100.0_wp !! endpoint of integration
  real(wp),dimension(n),parameter :: x0 = [0.0_wp,0.1_wp] !! initial x value

  real(wp),dimension(n) :: xf !! final x value
  type(rktp86_class) :: prop
  character(len=:),allocatable :: message

  call prop%initialize(n=n,f=fvpol,rtol=[tol],atol=[tol])
  call prop%integrate(t0,x0,dt,tf,xf)
  call prop%status(message=message)

  write (output_unit,'(A)') message
  write (output_unit,'(A,F7.2/,A,2E18.10)') &
              'tf =',tf ,'xf =',xf(1),xf(2)

contains

  subroutine fvpol(me,t,x,f)
    !! Right-hand side of van der Pol equation

    class(rk_class),intent(inout)     :: me
    real(wp),intent(in)               :: t
    real(wp),dimension(:),intent(in)  :: x
    real(wp),dimension(:),intent(out) :: f

    f(1) = x(2)
    f(2) = 0.2_wp*(1.0_wp-x(1)**2)*x(2) - x(1)

  end subroutine fvpol

end program rklib_example