leapBlock Module

Defines structured blocks. A block contains three axes and may be standalone or as part of an domain partition.

Axes define mid points, nodal points, and spacing as following:

                 xm(lo)               xm(hi)
                   |                   |
                |-----|-----|------|------|
                |     <-dxm->             |
              x(lo)                     x(hi+1)

Additionally, ghostcells can be added to the extremities of the axes to faciliate parallel operations.

A standalone (serial) block contains three axes for each direction, information on the number of grid points in each direction, information on periodicity, and positions of corner points. Many serial blocks can be instantiated as the same time.

                                           pmax
       x2(hi+1) |-----|-----|------|------|
                |     |     |      |      |
                |-----|-----|------|------|
                |     |     |      |      |
                |-----|-----|------|------|
                |     |     |      |      |
                |-----|-----|------|------|
                |     |     |      |      |
                |-----|-----|------|------|
                |     |     |      |      |
         x2(lo) |-----|-----|------|------|
            pmin
              x(lo)                     x(hi+1)

A partitioned block represents a sub-block from a larger global block. It only contains information about its own sub-grid. E.g.:

               2-by-2 domain decomposition
                |-----------|-----------|
                |           |           |
                |  RANK=3   |  RANK=4   |
                |           |           |
                |-----------|-----------|
                |           |           |
                |  RANK=1   |  RANK=2   |
                |           |           |
                |-----------|-----------|


Uses

  • module~~leapblock~~UsesGraph module~leapblock leapBlock iso_fortran_env iso_fortran_env module~leapblock->iso_fortran_env module~leapio_hdf5 leapIO_hdf5 module~leapblock->module~leapio_hdf5 module~leapkinds leapKinds module~leapblock->module~leapkinds module~leapparallel leapParallel module~leapblock->module~leapparallel mpi_f08 mpi_f08 module~leapblock->mpi_f08 module~leapio_hdf5->module~leapkinds module~leapio_hdf5->module~leapparallel hdf5 hdf5 module~leapio_hdf5->hdf5 module~leaputils leapUtils module~leapio_hdf5->module~leaputils module~leapkinds->iso_fortran_env module~leapparallel->iso_fortran_env module~leapparallel->module~leapkinds module~leapparallel->mpi_f08 module~leaputils->module~leapkinds

Used by

  • module~~leapblock~~UsedByGraph module~leapblock leapBlock module~cdifs cdifs module~cdifs->module~leapblock module~collisions collisions module~cdifs->module~collisions module~leapbc leapBC module~cdifs->module~leapbc module~leapdiffop leapDiffOp module~cdifs->module~leapdiffop module~leapeulerian leapEulerian module~cdifs->module~leapeulerian module~leaphypre leapHypre module~cdifs->module~leaphypre module~particles_resolved particles_resolved module~cdifs->module~particles_resolved module~immersed_boundaries immersed_boundaries module~cdifs->module~immersed_boundaries module~collisions->module~leapblock module~particles_point particles_point module~collisions->module~particles_point module~collisions->module~particles_resolved module~collisions->module~immersed_boundaries module~grans grans module~grans->module~leapblock module~grans->module~collisions module~grans->module~leapbc module~grans->module~leapdiffop module~grans->module~leapeulerian module~grans->module~leaphypre module~grans->module~particles_point module~grans->module~particles_resolved module~grans->module~immersed_boundaries module~immersed_boundaries_markers immersed_boundaries_markers module~immersed_boundaries_markers->module~leapblock module~immersed_boundaries_markers->module~leapbc module~immersed_boundaries_markers->module~leapdiffop module~immersed_boundaries_markers->module~leapeulerian module~immersed_boundaries_markers->module~leaphypre module~leaplagrangian leapLagrangian module~immersed_boundaries_markers->module~leaplagrangian module~immersed_boundaries_solids immersed_boundaries_solids module~immersed_boundaries_solids->module~leapblock module~immersed_boundaries_solids->module~immersed_boundaries_markers module~immersed_boundaries_solids->module~leapeulerian module~immersed_boundaries_solids->module~leaplagrangian module~leapbc->module~leapblock module~leapbc->module~leapeulerian module~leapdiffop->module~leapblock module~leapdiffop->module~leapbc module~leapdiffop->module~leapeulerian module~leapeulerian->module~leapblock module~leaphypre->module~leapblock module~leaphypre->module~leapeulerian module~leaplagrangian->module~leapblock module~particles_point->module~leapblock module~particles_point->module~leapdiffop module~particles_point->module~leapeulerian module~particles_point->module~leaplagrangian module~particles_point->module~immersed_boundaries module~particles_resolved->module~leapblock module~particles_resolved->module~leapbc module~particles_resolved->module~leapdiffop module~particles_resolved->module~leapeulerian module~particles_resolved->module~leaplagrangian module~particles_resolved->module~immersed_boundaries proc~collision_obj_setupcollisionblock2 collision_obj%collision_obj_SetupCollisionBlock2 proc~collision_obj_setupcollisionblock2->module~leapblock module~cdifs_advancesolution_smod cdifs_AdvanceSolution_smod module~cdifs_advancesolution_smod->module~cdifs module~cdifs_monitor_smod cdifs_Monitor_smod module~cdifs_monitor_smod->module~cdifs module~cdifs_preparesolver_smod cdifs_PrepareSolver_smod module~cdifs_preparesolver_smod->module~cdifs module~cdifs_writeoutputdata_smod cdifs_WriteOutputData_smod module~cdifs_writeoutputdata_smod->module~cdifs module~cdifs_writerestartdata_smod cdifs_WriteRestartData_smod module~cdifs_writerestartdata_smod->module~cdifs module~grans_advancesolution_smod grans_AdvanceSolution_smod module~grans_advancesolution_smod->module~grans module~grans_module_smod grans_module_smod module~grans_module_smod->module~grans module~grans_preparesolver_smod grans_PrepareSolver_smod module~grans_preparesolver_smod->module~grans module~grans_writeoutputdata_smod grans_WriteOutputData_smod module~grans_writeoutputdata_smod->module~grans module~grans_writerestartdata_smod grans_WriteRestartData_smod module~grans_writerestartdata_smod->module~grans module~immersed_boundaries->module~immersed_boundaries_markers module~immersed_boundaries->module~immersed_boundaries_solids proc~bc_set_buildmask bc_set%bc_set_BuildMask proc~bc_set_buildmask->module~leapeulerian proc~cdifs_obj_preparesolverbodyforce cdifs_obj_PrepareSolverBodyforce proc~cdifs_obj_preparesolverbodyforce->module~leapbc proc~cdifs_obj_preparesolveroperators cdifs_obj_PrepareSolverOperators proc~cdifs_obj_preparesolveroperators->module~leapbc proc~cdifs_obj_preparesolveroperatorsdiv cdifs_obj_PrepareSolverOperatorsDIV proc~cdifs_obj_preparesolveroperatorsdiv->module~leapbc proc~cdifs_obj_preparesolveroperatorspgrad cdifs_obj_PrepareSolverOperatorsPGRAD proc~cdifs_obj_preparesolveroperatorspgrad->module~leapbc proc~cdifs_obj_preparesolveroperatorsvlap cdifs_obj_PrepareSolverOperatorsVLAP proc~cdifs_obj_preparesolveroperatorsvlap->module~leapbc proc~grans_obj_preparesolveroperators grans_obj_PrepareSolverOperators proc~grans_obj_preparesolveroperators->module~leapdiffop proc~marker_set_computesolidvolfrac marker_set%marker_set_ComputeSolidVolFrac proc~marker_set_computesolidvolfrac->module~leapdiffop program~main main program~main->module~cdifs program~main->module~grans

Derived Types

type, public ::  axis_obj

Defines a 1D axis

Components

Type Visibility Attributes Name Initial
real(kind=wp), public, pointer, contiguous :: dxm(:) => null()

Cell spacings

integer, public :: hi

Higher bound

integer, public :: lo

Lower bound

integer, public :: ngc

Number of ghostcells

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

Nodal points

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

Mid points

Type-Bound Procedures

procedure, public :: Finalize => axis_obj_Final
procedure, public :: Initialize => axis_obj_Init

type, public ::  block_obj

A block object

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, contiguous, pointer :: dxm(:) => null()

Convenience pointer to x-cell spacings

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

Convenience pointer to y-cell spacings

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

Convenience pointer to z-cell spacings

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

MPI derived type for integer gc

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

MPI derived type for wp real gc

type(hdf5_obj), public, pointer :: hdf5 => null()

HDF5 object for IO

integer, public :: hi(3)

Array upper bound

logical, public :: is_initialized = .false.

Flag to determine whether this has been initialized

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, contiguous, pointer :: 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, contiguous, pointer :: 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, contiguous, pointer :: y(:) => null()

Convenience pointer to y-nodal points

real(kind=wp), public, contiguous, pointer :: 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, contiguous, pointer :: z(:) => null()

Convenience pointer to z-nodal points

real(kind=wp), public, contiguous, pointer :: 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
procedure, public :: GetOwningCell => block_obj_GetOwningCell
procedure, public :: GetOwningCellWGC => block_obj_GetOwningCellWGC
procedure, public :: Info => block_obj_Info
generic, public :: Initialize => block_obj_Init, block_obj_Init2
procedure, public :: Locate => block_obj_Locate
procedure, public :: Partition => block_obj_Partition
procedure, public :: Read => block_obj_Read
procedure, public :: SetConveniencePointers => block_obj_SetConveniencePointers
procedure, public :: SetPeriodicity => block_obj_SetPeriodicity
procedure, public :: SetupMPITypes => block_obj_SetupMPITypes
procedure, public :: SetupUniformGrid => block_obj_SetupUniformGrid
procedure, public, nopass :: SubDivideBlock => block_obj_SubDivideBlock
procedure, public :: UpdateExtents => block_obj_UpdateExtents
procedure, public :: UpdateGridGhostCells => block_obj_UpdateGridGhostCells
procedure, public :: UpdateMidPoints => block_obj_UpdateMidPoints
procedure, public :: UpdateSpacing => block_obj_UpdateSpacing
procedure, public :: Write => block_obj_Write
procedure, private :: block_obj_Init
procedure, private :: block_obj_Init2

Functions

private pure function block_obj_GetOwningCell(this, p) result(cell)

Returns the coordinates of the cell that contains the position p. This subroutine ignores gost cells. As such, it will always return an internal cell of the domain. If the point is located in a ghost cell, the returned cell will be an internal boundary cell that is closest to that ghost cell.

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, (3)

Cell coordinates

private pure function block_obj_GetOwningCellWGC(this, p) result(cell)

Returns the coordinates of the cell that contains the position p, inclusive of ghost cells.

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, (3)

Cell coordinates

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

Returns 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


Subroutines

private pure subroutine axis_obj_Final(this)

Finalizes object and frees data.

Arguments

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

A axis object

private pure subroutine axis_obj_Init(this, lo, hi, ngc)

Initialize axis.

Arguments

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

A axis object

integer, intent(in) :: lo

Array lower bound

integer, intent(in) :: hi

Array higher bound

integer, intent(in) :: ngc

Number of ghost cells

private impure subroutine block_obj_Final(this)

Finalizes the block object.

Arguments

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

A block object

private impure subroutine block_obj_Info(this)

Prints to stdout information about this block

Arguments

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

A block object

private impure subroutine block_obj_Init(this, ngc, parallel)

Initializes 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 to link with

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

[DEPRECATED] Initializes 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 to link with

private impure subroutine block_obj_Partition(this, Nb)

Partitions 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

private impure subroutine block_obj_Read(this, name)

Reads 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

private pure subroutine block_obj_SetConveniencePointers(this)

Associate convenience pointers.

Arguments

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

A block object

private pure subroutine block_obj_SetPeriodicity(this, periods)

Sets 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

private impure subroutine block_obj_SetupMPITypes(this)

Defines MPI derived type for communicating ghostcells.

Arguments

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

A block object

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

Initializes 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

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

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

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 pure subroutine block_obj_UpdateMidPoints(this)

Updates mid points.

Arguments

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

A block object

private impure subroutine block_obj_UpdateSpacing(this)

Updates grid spacing arrays.

Arguments

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

A block object

private impure subroutine block_obj_Write(this, name)

Writes 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