block_obj Derived Type

type, public :: block_obj

A block object


Inherits

type~~block_obj~~InheritsGraph type~block_obj block_obj MPI_Datatype MPI_Datatype type~block_obj->MPI_Datatype gc_slab_r, gc_slab_i type~axis_obj axis_obj type~block_obj->type~axis_obj axis, axis_partition type~parallel_obj parallel_obj type~block_obj->type~parallel_obj parallel type~parallel_obj->MPI_Datatype REAL_SP, REAL_DP, REAL_WP, COMPLEX_SP, COMPLEX_DP, COMPLEX_WP, INTEGER, INT8, LOGICAL MPI_Info MPI_Info type~parallel_obj->MPI_Info mpi_info type~communicators communicators type~parallel_obj->type~communicators comm type~patch patch type~parallel_obj->type~patch rank MPI_Comm MPI_Comm type~communicators->MPI_Comm w, g

Inherited by

type~~block_obj~~InheritedByGraph type~block_obj block_obj type~bc_set bc_set type~bc_set->type~block_obj block type~region_obj region_obj type~bc_set->type~region_obj region type~case_obj case_obj type~case_obj->type~block_obj block type~cdifs_obj cdifs_obj type~cdifs_obj->type~block_obj block type~cdifs_obj->type~bc_set bcs type~eulerian_set eulerian_set type~cdifs_obj->type~eulerian_set fields type~hypre_obj hypre_obj type~cdifs_obj->type~hypre_obj hypre type~op_obj op_obj type~cdifs_obj->type~op_obj op type~cdifs_case_obj cdifs_case_obj type~cdifs_obj->type~cdifs_case_obj case type~eulerian_obj_i eulerian_obj_i type~cdifs_obj->type~eulerian_obj_i maskV type~eulerian_obj_r eulerian_obj_r type~cdifs_obj->type~eulerian_obj_r V, P, dP, SA, ibVF, ibF, ibN, Vold, resV, rhs, divu, Vm type~marker_set marker_set type~cdifs_obj->type~marker_set IB type~respart_set ResPart_set type~cdifs_obj->type~respart_set RP type~eulerian_obj_base eulerian_obj_base type~eulerian_obj_base->type~block_obj block type~eulerian_set->type~block_obj block type~eulerian_ptr eulerian_ptr type~eulerian_set->type~eulerian_ptr field type~grans_obj grans_obj type~grans_obj->type~block_obj block type~grans_obj->type~eulerian_set fields type~grans_obj->type~hypre_obj hypre type~grans_obj->type~eulerian_obj_r ibVF, PVF, SA, Fp type~grans_case_obj grans_case_obj type~grans_obj->type~grans_case_obj case type~grans_obj->type~marker_set IB type~particle_set particle_set type~grans_obj->type~particle_set PP type~grans_obj->type~respart_set RP type~hypre_obj->type~block_obj block type~hypre_obj->type~eulerian_obj_i irow type~lagrangian_set lagrangian_set type~lagrangian_set->type~block_obj cblock, block type~ngadata_obj ngadata_obj type~ngadata_obj->type~block_obj block type~op_obj->type~block_obj block type~op_obj->type~eulerian_obj_i mask type~region_obj->type~block_obj region type~solid_set solid_set type~solid_set->type~block_obj block type~solid_obj solid_obj type~solid_set->type~solid_obj p type~cdifs_case_obj->type~case_obj type~eulerian_obj_i->type~eulerian_obj_base type~eulerian_obj_r->type~eulerian_obj_base type~eulerian_ptr->type~eulerian_obj_base p type~grans_case_obj->type~case_obj type~marker_set->type~lagrangian_set type~marker_set->type~op_obj op type~particle_set->type~lagrangian_set type~respart_set->type~lagrangian_set type~respart_set->type~op_obj op type~respart_set->type~marker_set ib type~solid_obj->type~marker_set

Components

Type Visibility Attributes Name Initial
type(axis_obj), public :: axis(3)

Axes in x1, x2, and x3 directions

type(axis_obj), public :: axis_partition(3)

Axes of the main partition, (spacing equal to 1 sub-block)

real(kind=wp), public :: dx(3)

Minimum mesh spacing in each direction

real(kind=wp), public, pointer, contiguous :: dxm(:) => null()

Convenience pointer to x-cell spacings

real(kind=wp), public, pointer, contiguous :: dym(:) => null()

Convenience pointer to y-cell spacings

real(kind=wp), public, pointer, contiguous :: dzm(:) => null()

Convenience pointer to z-cell spacings

type(MPI_Datatype), public :: gc_slab_i(3) = MPI_DATATYPE_NULL

MPI derived type for integer gc

type(MPI_Datatype), public :: gc_slab_r(3) = MPI_DATATYPE_NULL

MPI derived type for wp real gc

integer, public :: hi(3)

Array upper bound

logical, public :: is_partitioned = .false.

Flag for parallel partitioning

integer, public :: lo(3)

Array lower bound

integer, public :: ngc = 2

Number of ghostcells

type(parallel_obj), public, pointer :: parallel => null()

Associated parallel structure

logical, public :: periods(3) = .false.

Periodicity in each direction

real(kind=wp), public :: pmax(3)

Min and max locations across all blocks

real(kind=wp), public :: pmin(3)

Min and max locations across all blocks

real(kind=wp), public, pointer, contiguous :: x(:) => null()

Convenience pointer to x-nodal points

real(kind=wp), public :: xhi(3)

Coordinate of the top right corner

real(kind=wp), public :: xlo(3)

Coordinate of the bottom left corner

real(kind=wp), public, pointer, contiguous :: xm(:) => null()

Convenience pointer to x-mid points

real(kind=wp), public :: xmax

Min and max locations across all blocks

real(kind=wp), public :: xmin

Min and max locations across all blocks

real(kind=wp), public, pointer, contiguous :: y(:) => null()

Convenience pointer to y-nodal points

real(kind=wp), public, pointer, contiguous :: ym(:) => null()

Convenience pointer to y-mid points

real(kind=wp), public :: ymax

Min and max locations across all blocks

real(kind=wp), public :: ymin

Min and max locations across all blocks

real(kind=wp), public, pointer, contiguous :: z(:) => null()

Convenience pointer to z-nodal points

real(kind=wp), public, pointer, contiguous :: zm(:) => null()

Convenience pointer to z-mid points

real(kind=wp), public :: zmax

Min and max locations across all blocks

real(kind=wp), public :: zmin

Min and max locations across all blocks


Type-Bound Procedures

procedure, public :: Finalize => block_obj_Final

  • private impure subroutine block_obj_Final(this)

    Finalize the block object

    Arguments

    Type IntentOptional Attributes Name
    class(block_obj), intent(inout) :: this

    A block object

procedure, public :: Info => block_obj_Info

  • private impure subroutine block_obj_Info(this)

    Print to stdout information about this block

    Arguments

    Type IntentOptional Attributes Name
    class(block_obj), intent(in) :: this

    A block object

generic, public :: Initialize => block_obj_Init, block_obj_Init2

  • private impure subroutine block_obj_Init(this, ngc, parallel)

    Initialize block object

    Arguments

    Type IntentOptional Attributes Name
    class(block_obj), intent(inout) :: this

    A block object

    integer, intent(in) :: ngc

    Number of ghostcells

    type(parallel_obj), intent(in), target :: parallel

    parallel structure from main program

  • private impure subroutine block_obj_Init2(this, xlo, xhi, lo, hi, ngc, parallel)

    [DEPRECATED] Initialize block object

    Arguments

    Type IntentOptional Attributes Name
    class(block_obj), intent(inout) :: this

    A block object

    real(kind=wp), intent(in) :: xlo(3)

    Coordinates of the bottom left corner

    real(kind=wp), intent(in) :: xhi(3)

    Coordinates of the top right corner

    integer, intent(in) :: lo(3)

    Array lower bound

    integer, intent(in) :: hi(3)

    Array upper bound

    integer, intent(in) :: ngc

    Number of ghostcells

    type(parallel_obj), intent(in), target :: parallel

    parallel structure from main program

procedure, public :: Locate => block_obj_Locate

  • private impure function block_obj_Locate(this, p) result(rank)

    Return block ID and rank of the block where the point is located using a binary search alogirthm. Note that this function assumes that the point is within the domain, i.e., (pmin <= p <= pmax) and that any treatment for periodicity has been previously applied.

    Arguments

    Type IntentOptional Attributes Name
    class(block_obj), intent(in) :: this

    A block object

    real(kind=wp), intent(in) :: p(3)

    Position to locate

    Return Value integer

    MPI rank

procedure, public :: Partition => block_obj_Partition

  • private impure subroutine block_obj_Partition(this, Nb)

    Partition a parent block into sub-blocks based on a given decomposition Nb(3) for parallel simulations. This will also define the partition axes, and update local bounds/extents, and global bounds/extents.

    Arguments

    Type IntentOptional Attributes Name
    class(block_obj), intent(inout) :: this

    A block object

    integer, intent(in) :: Nb(3)

    Number of blocks in each direction

procedure, public :: Read => block_obj_Read

  • private impure subroutine block_obj_Read(this, name)

    Read block data using HDF5

    Arguments

    Type IntentOptional Attributes Name
    class(block_obj), intent(inout) :: this

    A block object

    character(len=*), intent(in) :: name

    Name of file to write

procedure, public :: SetConveniencePointers => block_obj_SetConveniencePointers

procedure, public :: SetPeriodicity => block_obj_SetPeriodicity

  • private pure subroutine block_obj_SetPeriodicity(this, periods)

    Set block periodicity in each direction

    Arguments

    Type IntentOptional Attributes Name
    class(block_obj), intent(inout) :: this

    A block object

    logical, intent(in) :: periods(3)

    Periodicity

procedure, public :: SetupMPITypes => block_obj_SetupMPITypes

  • private impure subroutine block_obj_SetupMPITypes(this)

    Define MPI derived type for communicating ghostcells

    Arguments

    Type IntentOptional Attributes Name
    class(block_obj), intent(inout) :: this

    A block object

procedure, public :: SetupUniformGrid => block_obj_SetupUniformGrid

  • private impure subroutine block_obj_SetupUniformGrid(this, xlo, xhi, lo, hi)

    Initialize a uniform grid on this block

    Arguments

    Type IntentOptional Attributes Name
    class(block_obj), intent(inout) :: this

    A block object

    real(kind=wp), intent(in) :: xlo(3)

    Coordinates of the bottom left corner

    real(kind=wp), intent(in) :: xhi(3)

    Coordinates of the top right corner

    integer, intent(in) :: lo(3)

    Array lower bound

    integer, intent(in) :: hi(3)

    Array upper bound

procedure, public, nopass :: SubDivideBlock => block_obj_SubDivideBlock

  • private pure subroutine block_obj_SubDivideBlock(coord, Nb, minlo, maxhi, sublo, subhi)

    Computes the bounds of the sub-block form the bounds of the parent block. Each sub-block gets about the same number of grid points in each direction.

    Arguments

    Type IntentOptional Attributes Name
    integer, intent(in) :: coord(3)

    Coordinates of MPI rank on MPI grid

    integer, intent(in) :: Nb(3)

    Number of blocks in each direction

    integer, intent(in) :: minlo(3)

    Global lo bounds

    integer, intent(in) :: maxhi(3)

    Global hi bounds

    integer, intent(out) :: sublo(3)

    lo bounds of sub-block

    integer, intent(out) :: subhi(3)

    hi bounds of sub-block

procedure, public :: UpdateExtents => block_obj_UpdateExtents

  • private impure subroutine block_obj_UpdateExtents(this)

    Updates the dimensional extents of the block

    Arguments

    Type IntentOptional Attributes Name
    class(block_obj), intent(inout) :: this

    A block object

generic, public :: UpdateGridGhostCells => block_obj_UpdateGridGhostCells, block_obj_UpdateGridGhostCells2

  • private impure subroutine block_obj_UpdateGridGhostCells(this)

    Updates the ghostcell values of local grid owned by the current MPI rank. Note that each MPI rank stores only its portion of the grid, thus needs to have proper ghostcell values.

    Arguments

    Type IntentOptional Attributes Name
    class(block_obj), intent(inout) :: this

    A block object

  • private subroutine block_obj_UpdateGridGhostCells2(this, axis, idir)

    [DEPRECATED] Updates the ghostcell values of local grid owned by the current MPI rank. Note that each MPI rank stores only its portion of the grid, thus needs to have proper ghostcell values. Although SetupUniformGrid fills the ghostcells of x/y/z, it does it assuming fixed grid spacing, which may not be the correct if a non-uniform grid is used.

    Arguments

    Type IntentOptional Attributes Name
    class(block_obj), intent(inout) :: this

    A block object

    real(kind=wp), intent(inout) :: axis(:)
    integer, intent(in) :: idir

procedure, public :: UpdateMidPoints => block_obj_UpdateMidPoints

  • private pure subroutine block_obj_UpdateMidPoints(this)

    Update mid points

    Arguments

    Type IntentOptional Attributes Name
    class(block_obj), intent(inout) :: this

    A block object

procedure, public :: UpdateSpacing => block_obj_UpdateSpacing

  • private impure subroutine block_obj_UpdateSpacing(this)

    Update grid spacing arrays

    Arguments

    Type IntentOptional Attributes Name
    class(block_obj), intent(inout) :: this

    A block object

procedure, public :: Write => block_obj_Write

  • private impure subroutine block_obj_Write(this, name)

    Write block data using HDF5

    Arguments

    Type IntentOptional Attributes Name
    class(block_obj), intent(in) :: this

    A block object

    character(len=*), intent(in) :: name

    Name of file to write

procedure, private :: block_obj_Init

  • private impure subroutine block_obj_Init(this, ngc, parallel)

    Initialize block object

    Arguments

    Type IntentOptional Attributes Name
    class(block_obj), intent(inout) :: this

    A block object

    integer, intent(in) :: ngc

    Number of ghostcells

    type(parallel_obj), intent(in), target :: parallel

    parallel structure from main program

procedure, private :: block_obj_Init2

  • private impure subroutine block_obj_Init2(this, xlo, xhi, lo, hi, ngc, parallel)

    [DEPRECATED] Initialize block object

    Arguments

    Type IntentOptional Attributes Name
    class(block_obj), intent(inout) :: this

    A block object

    real(kind=wp), intent(in) :: xlo(3)

    Coordinates of the bottom left corner

    real(kind=wp), intent(in) :: xhi(3)

    Coordinates of the top right corner

    integer, intent(in) :: lo(3)

    Array lower bound

    integer, intent(in) :: hi(3)

    Array upper bound

    integer, intent(in) :: ngc

    Number of ghostcells

    type(parallel_obj), intent(in), target :: parallel

    parallel structure from main program

procedure, private :: block_obj_UpdateGridGhostCells

  • private impure subroutine block_obj_UpdateGridGhostCells(this)

    Updates the ghostcell values of local grid owned by the current MPI rank. Note that each MPI rank stores only its portion of the grid, thus needs to have proper ghostcell values.

    Arguments

    Type IntentOptional Attributes Name
    class(block_obj), intent(inout) :: this

    A block object

procedure, private :: block_obj_UpdateGridGhostCells2

  • private subroutine block_obj_UpdateGridGhostCells2(this, axis, idir)

    [DEPRECATED] Updates the ghostcell values of local grid owned by the current MPI rank. Note that each MPI rank stores only its portion of the grid, thus needs to have proper ghostcell values. Although SetupUniformGrid fills the ghostcells of x/y/z, it does it assuming fixed grid spacing, which may not be the correct if a non-uniform grid is used.

    Arguments

    Type IntentOptional Attributes Name
    class(block_obj), intent(inout) :: this

    A block object

    real(kind=wp), intent(inout) :: axis(:)
    integer, intent(in) :: idir