Report Typos and Errors    
Semi-Lagrangian Library
Modular library for kinetic and gyrokinetic simulations of plasmas in fusion energy devices.
sll_m_ode_integrator_base.F90
Go to the documentation of this file.
1 !------------------------------------------------------------------------------
2 ! SELALIB
3 !------------------------------------------------------------------------------
4 ! MODULE: sll_m_ode_integrator_base
5 !
6 ! DESCRIPTION:
12 !------------------------------------------------------------------------------
14 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
15 
16  use sll_m_working_precision, only: &
17  f64
18 
19  use sll_m_vector_space_base, only: &
21 
22  implicit none
23 
24  public :: &
25  sll_c_ode, &
27 
28  private
29 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
30 
32  integer, parameter :: wp = f64
33 
34  !----------------------------------------------------------------------------
35  ! Abstract type: sll_c_ode
36  !----------------------------------------------------------------------------
45  type, abstract :: sll_c_ode
46  contains
47  procedure(i_rhs), deferred :: rhs
48 ! procedure :: solve !< solve an implicit problem
49  end type sll_c_ode
50 
51  !----------------------------------------------------------------------------
52  ! Abstract type: sll_c_ode_integrator
53  !----------------------------------------------------------------------------
62  type, abstract :: sll_c_ode_integrator
63 
65  class(sll_c_ode), pointer :: ode => null()
66 
68  class(sll_c_vector_space), allocatable :: work(:)
69 
70  contains
71  procedure(i_init), deferred :: init
72  procedure(i_step), deferred :: step
73  procedure(i_clean), deferred :: clean
74 
75  end type sll_c_ode_integrator
76 
77  !----------------------------------------------------------------------------
78  ! Abstract interface: i_rhs
79  !----------------------------------------------------------------------------
86  abstract interface
87  subroutine i_rhs(self, t, y, ydot)
89  class(sll_c_ode), intent(inout) :: self
90  real(wp), intent(in) :: t
91  class(sll_c_vector_space), intent(in) :: y
92  class(sll_c_vector_space), intent(inout) :: ydot
93  end subroutine i_rhs
94  end interface
95 
96  !----------------------------------------------------------------------------
97  ! Abstract interface: i_init
98  !----------------------------------------------------------------------------
110  abstract interface
111  subroutine i_init(self, ode, t0, y0)
113  class(sll_c_ode_integrator), intent(out) :: self
114  class(sll_c_ode), pointer, intent(in) :: ode
115  real(wp), intent(in) :: t0
116  class(sll_c_vector_space), intent(inout) :: y0
117  end subroutine i_init
118  end interface
119 
120  !----------------------------------------------------------------------------
121  ! Abstract interface: i_step
122  !----------------------------------------------------------------------------
130  abstract interface
131  subroutine i_step(self, t, y, h, ynew)
133  class(sll_c_ode_integrator), intent(inout) :: self
134  real(wp), intent(in) :: t
135  class(sll_c_vector_space), intent(in) :: y
136  real(wp), intent(in) :: h
137  class(sll_c_vector_space), intent(inout) :: ynew
138  end subroutine i_step
139  end interface
140 
141  !----------------------------------------------------------------------------
142  ! Abstract interface: i_clean
143  !----------------------------------------------------------------------------
147  abstract interface
148  subroutine i_clean(self)
149  import sll_c_ode_integrator
150  class(sll_c_ode_integrator), intent(inout) :: self
151  end subroutine i_clean
152  end interface
153 
154 !==============================================================================
155 !contains
156 !==============================================================================
157 
158 ! subroutine solve( self, t, y, ynew, b, sigma )
159 ! class( sll_c_ode ), intent( inout ) :: self
160 ! real(wp) , intent( in ) :: t
161 ! end subroutine solve
162 
163 end module sll_m_ode_integrator_base
Compute the time derivative of the state vector.
Advance the solution by one time step.
Abstract types for: 1) generic ODE system, and 2) ODE integrator.
integer, parameter wp
Working precision.
Abstract type implementing a generic vector space.
Module to select the kind parameter.
integer, parameter, public f64
f64 is the kind type for 64-bit reals (double precision)
Abstract base class for all vector spaces.
    Report Typos and Errors