m_base_backend Module



Abstract Interfaces

abstract interface

  • public subroutine transeq_ders(self, du, dv, dw, u, v, w, dirps)

    transeq equation obtains the derivatives direction by direction, and the exact algorithm used to obtain these derivatives are decided at runtime. Backend implementations are responsible from directing calls to transeq_ders into the correct algorithm.

    Arguments

    Type IntentOptional Attributes Name
    class(base_backend_t) :: self
    class(field_t), intent(inout) :: du
    class(field_t), intent(inout) :: dv
    class(field_t), intent(inout) :: dw
    class(field_t), intent(in) :: u
    class(field_t), intent(in) :: v
    class(field_t), intent(in) :: w
    type(dirps_t), intent(in) :: dirps

abstract interface

  • public subroutine tds_solve(self, du, u, tdsops)

    transeq equation obtains the derivatives direction by direction, and the exact algorithm used to obtain these derivatives are decided at runtime. Backend implementations are responsible from directing calls to tds_solve to the correct algorithm.

    Arguments

    Type IntentOptional Attributes Name
    class(base_backend_t) :: self
    class(field_t), intent(inout) :: du
    class(field_t), intent(in) :: u
    class(tdsops_t), intent(in) :: tdsops

abstract interface

  • public subroutine reorder(self, u_, u, direction)

    reorder subroutines are straightforward, they rearrange data into our specialist data structure so that regardless of the direction tridiagonal systems are solved efficiently and fast.

    Arguments

    Type IntentOptional Attributes Name
    class(base_backend_t) :: self
    class(field_t), intent(inout) :: u_
    class(field_t), intent(in) :: u
    integer, intent(in) :: direction

abstract interface

  • public subroutine sum_intox(self, u, u_)

    sum9into3 subroutine combines all the directional velocity derivatives into the corresponding x directional fields.

    Arguments

    Type IntentOptional Attributes Name
    class(base_backend_t) :: self
    class(field_t), intent(inout) :: u
    class(field_t), intent(in) :: u_

abstract interface

  • public subroutine vecadd(self, a, x, b, y)

    adds two vectors together: y = ax + by

    Arguments

    Type IntentOptional Attributes Name
    class(base_backend_t) :: self
    real(kind=dp), intent(in) :: a
    class(field_t), intent(in) :: x
    real(kind=dp), intent(in) :: b
    class(field_t), intent(inout) :: y

abstract interface

  • public function scalar_product(self, x, y) result(s)

    Calculates the scalar product of two input fields

    Arguments

    Type IntentOptional Attributes Name
    class(base_backend_t) :: self
    class(field_t), intent(in) :: x
    class(field_t), intent(in) :: y

    Return Value real(kind=dp)

abstract interface

  • public subroutine copy_data_to_f(self, f, data)

    Copy the specialist data structure from device or host back to a regular 3D data array in host memory.

    Arguments

    Type IntentOptional Attributes Name
    class(base_backend_t), intent(inout) :: self
    class(field_t), intent(inout) :: f
    real(kind=dp), intent(in), dimension(:, :, :) :: data

abstract interface

  • public subroutine copy_f_to_data(self, data, f)

    Copy a regular 3D array in host memory into the specialist data structure field that lives on device or host

    Arguments

    Type IntentOptional Attributes Name
    class(base_backend_t), intent(inout) :: self
    real(kind=dp), intent(out), dimension(:, :, :) :: data
    class(field_t), intent(in) :: f

abstract interface

  • public subroutine alloc_tdsops(self, tdsops, dir, operation, scheme, n_halo, from_to, bc_start, bc_end, sym, c_nu, nu0_nu)

    Arguments

    Type IntentOptional Attributes Name
    class(base_backend_t) :: self
    class(tdsops_t), intent(inout), allocatable :: tdsops
    integer, intent(in) :: dir
    character(len=*), intent(in) :: operation
    character(len=*), intent(in) :: scheme
    integer, intent(in), optional :: n_halo
    character(len=*), intent(in), optional :: from_to
    character(len=*), intent(in), optional :: bc_start
    character(len=*), intent(in), optional :: bc_end
    logical, intent(in), optional :: sym
    real(kind=dp), intent(in), optional :: c_nu
    real(kind=dp), intent(in), optional :: nu0_nu

abstract interface

  • public subroutine init_poisson_fft(self, mesh, xdirps, ydirps, zdirps)

    Arguments

    Type IntentOptional Attributes Name
    class(base_backend_t) :: self
    class(mesh_t), intent(in) :: mesh
    type(dirps_t), intent(in) :: xdirps
    type(dirps_t), intent(in) :: ydirps
    type(dirps_t), intent(in) :: zdirps

Derived Types

type, public, abstract ::  base_backend_t

base_backend class defines all the abstract operations that the solver class requires.

Read more…

Components

Type Visibility Attributes Name Initial
real(kind=dp), public :: nu
class(mesh_t), public, pointer :: mesh
class(allocator_t), public, pointer :: allocator
class(poisson_fft_t), public, pointer :: poisson_fft

Type-Bound Procedures

procedure(transeq_ders), public, deferred :: transeq_x
procedure(transeq_ders), public, deferred :: transeq_y
procedure(transeq_ders), public, deferred :: transeq_z
procedure(tds_solve), public, deferred :: tds_solve
procedure(reorder), public, deferred :: reorder
procedure(sum_intox), public, deferred :: sum_yintox
procedure(sum_intox), public, deferred :: sum_zintox
procedure(vecadd), public, deferred :: vecadd
procedure(scalar_product), public, deferred :: scalar_product
procedure(copy_data_to_f), public, deferred :: copy_data_to_f
procedure(copy_f_to_data), public, deferred :: copy_f_to_data
procedure(alloc_tdsops), public, deferred :: alloc_tdsops
procedure(init_poisson_fft), public, deferred :: init_poisson_fft
procedure, public :: base_init
procedure, public :: get_field_data
procedure, public :: set_field_data

Subroutines

public subroutine base_init(self)

Arguments

Type IntentOptional Attributes Name
class(base_backend_t) :: self

public subroutine get_field_data(self, data, f, dir)

Extract data from field f optionally reordering into dir orientation. To output in same orientation as f, use call ...%get_field_data(data, f, f%dir)

Arguments

Type IntentOptional Attributes Name
class(base_backend_t) :: self
real(kind=dp), intent(out), dimension(:, :, :) :: data

Output array

class(field_t), intent(in) :: f

Field

integer, intent(in), optional :: dir

Desired orientation of output array (defaults to Cartesian)

public subroutine set_field_data(self, f, data, dir)

Arguments

Type IntentOptional Attributes Name
class(base_backend_t) :: self
class(field_t), intent(inout) :: f

Field

real(kind=dp), intent(in), dimension(:, :, :) :: data

Input array

integer, intent(in), optional :: dir

Orientation of input array (defaults to Cartesian)