base_backend_t Derived Type

type, public, abstract :: base_backend_t

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

For example, transport equation in solver class evaluates the derivatives in x, y, and z directions, and reorders the input fields as required. Then finally, combines all the directional derivatives to obtain the divergence of U*.

All these high level operations solver class executes are defined here using the abstract interfaces. Every backend implementation extends the present abstact backend class to define the specifics of these operations based on the target architecture.


Inherits

type~~base_backend_t~~InheritsGraph type~base_backend_t base_backend_t type~allocator_t allocator_t type~base_backend_t->type~allocator_t allocator type~mesh_t mesh_t type~base_backend_t->type~mesh_t mesh type~poisson_fft_t poisson_fft_t type~base_backend_t->type~poisson_fft_t poisson_fft type~allocator_t->type~mesh_t mesh type~field_t field_t type~allocator_t->type~field_t first type~geo_t geo_t type~mesh_t->type~geo_t geo type~grid_t grid_t type~mesh_t->type~grid_t grid type~par_t par_t type~mesh_t->type~par_t par type~field_t->type~field_t next

Inherited by

type~~base_backend_t~~InheritedByGraph type~base_backend_t base_backend_t type~cuda_backend_t cuda_backend_t type~cuda_backend_t->type~base_backend_t type~omp_backend_t omp_backend_t type~omp_backend_t->type~base_backend_t type~solver_t solver_t type~solver_t->type~base_backend_t backend type~time_intg_t time_intg_t type~solver_t->type~time_intg_t time_integrator type~vector_calculus_t vector_calculus_t type~solver_t->type~vector_calculus_t vector_calculus type~time_intg_t->type~base_backend_t backend type~vector_calculus_t->type~base_backend_t backend type~base_case_t base_case_t type~base_case_t->type~solver_t solver type~case_channel_t case_channel_t type~case_channel_t->type~base_case_t type~case_generic_t case_generic_t type~case_generic_t->type~base_case_t type~case_tgv_t case_tgv_t type~case_tgv_t->type~base_case_t

Components

Type Visibility Attributes Name Initial
real(kind=dp), public :: nu
type(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

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

    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

procedure(transeq_ders), public, deferred :: transeq_y

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

    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

procedure(transeq_ders), public, deferred :: transeq_z

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

    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

procedure(tds_solve), public, deferred :: tds_solve

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

    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

procedure(reorder), public, deferred :: reorder

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

    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

procedure(sum_intox), public, deferred :: sum_yintox

  • subroutine sum_intox(self, u, u_) Prototype

    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_

procedure(sum_intox), public, deferred :: sum_zintox

  • subroutine sum_intox(self, u, u_) Prototype

    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_

procedure(vecadd), public, deferred :: vecadd

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

    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

procedure(scalar_product), public, deferred :: scalar_product

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

    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)

procedure(field_max_mean), public, deferred :: field_max_mean

  • subroutine field_max_mean(self, max_val, mean_val, f, enforced_data_loc) Prototype

    Obtains maximum and mean values in a field

    Arguments

    Type IntentOptional Attributes Name
    class(base_backend_t) :: self
    real(kind=dp), intent(out) :: max_val
    real(kind=dp), intent(out) :: mean_val
    class(field_t), intent(in) :: f
    integer, intent(in), optional :: enforced_data_loc

procedure(field_ops), public, deferred :: field_scale

  • subroutine field_ops(self, f, a) Prototype

    Scales or shifts a field by a

    Arguments

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

procedure(field_ops), public, deferred :: field_shift

  • subroutine field_ops(self, f, a) Prototype

    Scales or shifts a field by a

    Arguments

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

procedure(field_reduce), public, deferred :: field_volume_integral

  • function field_reduce(self, f) result(s) Prototype

    Reduces field to a scalar, example: volume integral

    Arguments

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

    Return Value real(kind=dp)

procedure(field_set_face), public, deferred :: field_set_face

  • subroutine field_set_face(self, f, c_start, c_end, face) Prototype

    A field is a subdomain with a rectangular cuboid shape. It has 6 faces, and these faces are either a subdomain boundary or a global domain boundary based on the location of the subdomain. This subroutine allows us to set any of these faces to a value, 'c_start' and 'c_end' for faces at opposite sides. 'face' is one of X_FACE, Y_FACE, Z_FACE from common.f90

    Arguments

    Type IntentOptional Attributes Name
    class(base_backend_t) :: self
    class(field_t), intent(inout) :: f
    real(kind=dp), intent(in) :: c_start
    real(kind=dp), intent(in) :: c_end
    integer, intent(in) :: face

procedure(copy_data_to_f), public, deferred :: copy_data_to_f

  • subroutine copy_data_to_f(self, f, data) Prototype

    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

procedure(copy_f_to_data), public, deferred :: copy_f_to_data

  • subroutine copy_f_to_data(self, data, f) Prototype

    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

procedure(alloc_tdsops), public, deferred :: alloc_tdsops

  • subroutine alloc_tdsops(self, tdsops, n_tds, delta, operation, scheme, bc_start, bc_end, stretch, stretch_correct, n_halo, from_to, sym, c_nu, nu0_nu) Prototype

    Arguments

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

procedure(init_poisson_fft), public, deferred :: init_poisson_fft

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

    Arguments

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

procedure, public :: base_init

procedure, public :: get_field_data

  • 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)

procedure, public :: set_field_data

  • 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)