multistep_m.f90 Source File


Contents

Source Code


Source Code

module multistep_m
    !! Time-integration methods based on multistep methods
    !!    
    !! multistep_storage_t:
    !! - takes care of storage and shift of data at previous timestep
    !!
    !! karniadakis_t
    !! - Explicit part of time-stepping according to Karniadakis scheme
    !! - Implicit part not treated within this module but must be handled explicitly by yourself
    !! - Automatic order reduction for initial time-stepd
    !! - Generic implementation, i.e. independent of mesh, equilibria, variable, etc...
    use precision_grillix_m, only : GP
    use screen_io_m, only :  get_stdout
    use parallelisation_setup_m, only : is_rank_info_writer
    use error_handling_grillix_m, only : handle_error, error_info_t
    use status_codes_grillix_m, only : GRILLIX_ERR_OTHER
    implicit none

    type, public :: multistep_storage_t
        !! Datatype for storage
        integer, private :: ndim
        !! Dimension of storage values
        integer, private :: nstorage
        !! Number of storage points
        integer, private :: nfields
        !! Number of fields to be stored
        real(GP), public, allocatable, dimension(:,:,:) :: vstore
        !! Values of storage (ndim, nstorage, nfields)
    contains
        procedure, public :: init_storage
        procedure, public :: get_ndim
        procedure, public :: get_nstorage
        procedure, public :: get_nfields
        procedure, public :: shift_storage
        final :: destructor_storage                        
    end type 
    
    interface

        pure integer module function get_ndim(self)
            !! Returns ndim
            class(multistep_storage_t), intent(in) :: self
            !! Instance of the type                   
        end function

        pure integer module function get_nstorage(self)
            !! Returns nstorage
            class(multistep_storage_t), intent(in) :: self
            !! Instance of the type                   
        end function

        pure integer module function get_nfields(self)
            !! Returns nfields
            class(multistep_storage_t), intent(in) :: self
            !! Instance of the type                   
        end function

        module subroutine init_storage(self, ndim, nstorage, nfields, vstore_firsts)
            !! Initialises a Karniadakis scheme
            class(multistep_storage_t), intent(inout) :: self
            !! Instance of the type 
            integer, intent(in) :: ndim
            !! Dimension of storage values
            integer, intent(in) :: nstorage
            !! Number of storage points 
            integer, intent(in) :: nfields
            !! Number of fields to be stored
            real(GP), dimension(ndim, nstorage, nfields), intent(in), optional :: vstore_firsts
            !! Initial values for storage (default = 0)
        end subroutine

        module subroutine shift_storage(self, vals_new)
            !! Shifts storage and vals_new will be placed at vsore(:,1,:) 
            class(multistep_storage_t), intent(inout) :: self
            !! Instance of the type 
             real(GP), dimension(self%ndim, self%nfields), intent(in) :: vals_new
            !! New storage value
        end subroutine    

        module subroutine destructor_storage(self)
            !! Frees memory associated with multistep_storage
            type(multistep_storage_t), intent(inout) :: self
            !! Instance of the type
        end subroutine

    end interface



    type, public, extends(multistep_storage_t) :: karniadakis_t
        !! Datatype for karniadakis
        integer, private :: ncount
        !! Counter for Karniadakis integrator
        integer, private :: order
        !! Order of time-stepping scheme
        integer, private :: current_order
        !! Current order of timestep (reduced for initial time-steps)     
        real(GP), private :: dtau
        !! Time step size
    contains
        procedure, public :: init => init_karniadakis
        procedure, public :: get_ncount
        procedure, public :: get_order
        procedure, public :: get_current_order
        procedure, public :: get_dtau
        procedure, public :: get_dtau_bdf
        procedure, public :: advance
        procedure, public :: display => display_karniadakis
        final :: destructor_karniadakis                        
    end type 
    
    interface
    
        pure integer module function get_ncount(self)
            !! Returns ncount 
            class(karniadakis_t), intent(in) :: self
            !! Instance of the type                   
        end function

        pure integer module function get_order(self)
            !! Returns order that eventually approached (after initial time-steps) 
            class(karniadakis_t), intent(in) :: self
            !! Instance of the type                   
        end function

        pure integer module function get_current_order(self)
            !! Returns order used for last time-step
            class(karniadakis_t), intent(in) :: self
            !! Instance of the type                   
        end function

        real(GP) module function get_dtau(self)
            !! Returns dtau
            class(karniadakis_t), intent(in) :: self
            !! Instance of the type                   
        end function

        real(GP) module function get_dtau_bdf(self)
            !! Returns dtau_bdf
            class(karniadakis_t), intent(in) :: self
            !! Instance of the type                   
        end function

        module subroutine init_karniadakis(self, ndim, order, dtau, y_firsts, dy_firsts)
            !! Initialises a Karniadakis scheme
            class(karniadakis_t), intent(inout) :: self
            !! Instance of the type 
            integer, intent(in) :: ndim
            !! Dimension of problem
            integer, intent(in) :: order
            !! Order of time-stepping scheme  
            real(GP), intent(in) :: dtau  
            !! Time step size   
            real(GP), dimension(ndim, order-1), intent(in), optional :: y_firsts
            !! Values for initial time steps in descending order 
            real(GP), dimension(ndim, order-1), intent(in), optional :: dy_firsts  
            !! Right hand side for initial time-steps in descending order   
        end subroutine

        module subroutine advance(self, dy, y, yn, dbgout)
            !! Advances variable from t -> t+1
            class(karniadakis_t), intent(inout) :: self
            !! Instance of the type 
            real(GP), dimension(self%ndim), intent(in) :: dy
            !! Right hand side of differential equation at time-step t
            real(GP), dimension(self%ndim), intent(in) :: y
            !! Variable at time point t, on output advanced to t+1
            real(GP), dimension(self%ndim), intent(out) :: yn
            !! Variable advanced to time point t+1
            integer, intent(in), optional :: dbgout 
            !! debug output level, default: 0
        end subroutine    

        module subroutine display_karniadakis(self)
            !! Displays basic information of karniadakis
            class(karniadakis_t), intent(in) :: self
            !! Instance of the type
        end subroutine     

        module subroutine destructor_karniadakis(self)
            !! Frees memory associated with karniadakis
            type(karniadakis_t), intent(inout) :: self
            !! Instance of the type
        end subroutine

    end interface

end module