mus_particle_aux_module Module

Module containing some auxiliary routines for LBM-DEM simulations for particles in a flow.



Functions

public function coordLocalOnMyRank(coord, scheme, myRank)

Routine which checks whether integer coordinate is local on the current rank

Arguments

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

integer coordinate (x,y,z,level)

type(mus_scheme_type), intent(in) :: scheme

Scheme for access to level descriptor

integer, intent(in) :: myRank

This process rank

Return Value logical

public function positionLocalOnMyRank(pos, geometry, scheme, myRank)

Routine which checks whether a spatial position is local on the current rank

Arguments

Type IntentOptional Attributes Name
real(kind=rk), intent(in) :: pos(3)

Cartesian position (x, y, z)

type(mus_geom_type), intent(in) :: geometry

Geometry to look in

type(mus_scheme_type), intent(in) :: scheme

Scheme for access to level descriptor

integer, intent(in) :: myRank

This process rank

Return Value logical

public function coordExistsOnMyRank(coord, scheme, myRank)

Checks if integer coordinate is in the total list of this rank

Arguments

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

integer coordinate (x,y,z,level)

type(mus_scheme_type), intent(in) :: scheme

Scheme for access to level descriptor

integer, intent(in) :: myRank

This process rank

Return Value logical

public pure function getBaryOfCoord(coord, origin, dx) result(bary)

Convenience function to get barycenter of integer coordinate

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: coord(3)
real(kind=rk), intent(in) :: origin(3)
real(kind=rk), intent(in) :: dx

Return Value real(kind=rk), (3)

public function getPosOfCoord(coord, scheme) result(ldPos)

Get the position in the levelDesc%total list of an element with integer coordinate coord. This routine returns ldPos <= 0 if the element could not be found on this process.

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: coord(4)
type(mus_scheme_type), intent(in) :: scheme

Return Value integer

public function followNeighPath(neighPath, startPos, neigh, nElems_fluid, zeroDir) result(endPos)

followNeighPath can be used to determine the position of a neighboring lattice site in the total list. For this we follow one or more entries of the levelDesc%neigh array. Since this array is only defined for local fluid elements, startPos must point to a local fluid element.

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: neighPath(:)

Array containing the sequence of indices of directions that must be followed to get to the direction described by neighDir. This depends on the stencil being used for the simulation. For example in the d3q19 stencil, to reach the neighbor located with displacement (1,1,1) from the element at startPos, we follow first (1,0,0) and then (0,1,1). For d3q27 neighPath will only be one entry (with value 26) corresponding to the index (1,1,1).

integer, intent(in) :: startPos

Position of the starting element in the total list This element must be a local fluid!

type(tem_levelNeighbor_type), intent(in) :: neigh

Neighbor list

integer, intent(in) :: nElems_fluid

Number of local fluid elements on this process We need this to check if intermediate values of endPos are still local fluids

integer, intent(in) :: zeroDir

Index of the "zero direction" of the stencil i.e. stencil%cxdir(1:3,zeroDir) = [0,0,0]

Return Value integer

consecutively neighPath(1), neighPath(2) etc.

public function findPropIndex(prp_bitpos, glob) result(propIndex)

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: prp_bitpos

Property bit of property to find index for

type(tem_global_type), intent(in) :: glob

tree%global mesh information

Return Value integer

private function treeIDlocalOnMyRank(TreeID, scheme, lev)

Arguments

Type IntentOptional Attributes Name
integer(kind=long_k), intent(in) :: TreeID

TreeID to check

type(mus_scheme_type), intent(in) :: scheme

Scheme for access to level descriptor

integer, intent(in) :: lev

Level

Return Value logical

private function has_remote_neighbors(iElem, scheme, lev)

Function which checks if element has remote neighbors

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: iElem

Index of this element in the total list

type(mus_scheme_type), intent(in) :: scheme

Scheme for access to levelDescriptor

integer, intent(in) :: lev

Level

Return Value logical


Subroutines

public subroutine findPartitionOfTreeID(sTreeID, geometry, procs, nProcs, outProc, procIndex)

Routine to find the process a certain treeID is on. Looks in the ranks in the array of argument procs Returns the INDEX of the found proc in this array Example: procs = [4 8 7] and we find our treeID on rank 8 will return procIndex = 2

Arguments

Type IntentOptional Attributes Name
integer(kind=long_k), intent(in) :: sTreeID

TreeID to look for

type(mus_geom_type), intent(in) :: geometry

Geometry to look in

integer, intent(in), optional :: procs(:)

Only look on these ranks (look in all ranks if omitted)

integer, intent(in) :: nProcs

Number of elements in procs

integer, intent(out) :: outProc

Process which the treeID is on

integer, intent(out) :: procIndex

Index of that process

public subroutine cross_product(a, b, res)

computes cross-product a x b = res

Arguments

Type IntentOptional Attributes Name
real(kind=rk), intent(in) :: a(3)
real(kind=rk), intent(in) :: b(3)
real(kind=rk), intent(out) :: res(3)

public subroutine getCartesianStencilIndices(cxDir, cartIndices)

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: cxDir(:,:)

Integer stencil direction

integer, intent(out) :: cartIndices(6)

Output indices of the Cartesian directions [-x +x -y +y -z +z] in cxDir

public subroutine getPathToVec(vec, cxDir, path)

getPathToVec builds a series of 3 indices which can be followed using the levelDesc(lev)%neigh array to get the position of a neighboring element in the total list. Usage: given an element at coordinate [0,0,0] with position startPos in the total list, to reach the element at [1,1,1] we do: call getPathToVec([1,1,1], cxDir, path) pos = startPos do k = 1,3 pos = levelDesc(lev)%neigh%nghElems(path(k),pos) end do At the end of the loop pos will point to the requested element

Arguments

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

Vector we'd like to build the path to

integer, intent(in) :: cxDir(:,:)

Integer stencil direction

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

Output path of iDirs

public subroutine getProcessBoundingBox(scheme, geometry, boundingBox)

Routine to find the bounding box of the part of the spatial domain on this process

Arguments

Type IntentOptional Attributes Name
type(mus_scheme_type), intent(in) :: scheme

Scheme for access to level descriptor

type(mus_geom_type), intent(in) :: geometry

Geometry information to determine TreeIDs of elements 'covered' by particle

real(kind=rk), intent(out) :: boundingBox(6)

public subroutine mus_particles_global_errorcheck(flag, comm)

Routine to aggregate local error flags and use MPI allreduce to set a global flag if any of the local flags is true

Arguments

Type IntentOptional Attributes Name
logical, intent(inout) :: flag

Error flag: at input this will be the local error code of this rank only which is set to TRUE if there is an error and FALSE if not. Upon output this will be set to flag = TRUE if ANY of the local error codes (on all processes in comm) were TRUE.

integer, intent(in) :: comm

MPI communicator