tem_face_module Module

This module provides information about the existence and properties of faces in an Octree based mesh. It is able to deliver info about the horizontal (i.e. on the same refinement level) and vertical (i.e. between the different levels) dependencies.

If you want to make use of face descriptions in your code, you just have to call tem_build_face_info to buid up the face descriptions in your solver. You have to make use of the created iterators and face descriptions. The faces structure created by tem_build_face_info is an array of tem_face_type with the face information for each level. On each level there are three tem_face_descriptor_type for the three spatial dimensions faces(iLevel)%faces(1:3). Each of these descriptors holds the complete list of faces on that level for this direction in the faceList.

In general, the construction of the face description follows the following procedure:

  1. Generate level descriptor in a dimension-by-dimension approach.
  2. Collect all the faces and attach properties to the faces.
  3. Extend the collected properties of the remote elements by from finer or from coarser properties (to indicate if a halo element is from coarser or from finer).
  4. Create vertical dependencies between the faces.
  5. Create list of faces which are relevant for computation and interpolation.
  6. Create communication buffers for send and receive of face information.

Each face has two neighboring elements. To identify each face uniquely, we use the tree ID of its left element.

       ---------------------------------------
      /                  /                  /
     /                  /                  /
    /   left elem.     /   right elem.    /
   /                  /                  /
  /                  /                  /
  --------------------------------------

This module provides a basic description of the faces per spatial direction. Building this description is achieved in several steps:

  1. Collect all the faces (each level independently).
  2. Determine the properties of the faces. a. Check if halo elements are refined on the remote partition.
  3. Build the vertical dependencies between the faces of the different levels.
  4. Build buffers and lists of faces from the previous steps.

To define a correct list of faces for the solvers, we attach left and right properties to each face (in part 1 of the upper algorithm). We differentiate between a left and right property of each face. The left one is determined by it's left element and the right one by it's right element respectively.

              tem_left / / / tem_right
       -----------------------------------------------------
      /                      / / /                        /
     /                      / / /                        /
    /   left elem. --->    / / /    <--- right elem.    /
   /                      / / /                        /
  /                      / / /                        /
  ----------------------------------------------------
                       / / /

We use one or more of the following properties to the left and right properties of the faces:

a. tem_fluidFace_prp b. tem_fromFinerFace_prp c. tem_fromCoarserFace_prp d. tem_remoteFace_prp

These properties are the basis for building the correct list of faces in step 4 of the above algorithm.

See [[Face Details]] for examples on face configurations.


Uses

Used by

  • module~~tem_face_module~~UsedByGraph module~tem_face_module tem_face_module program~tem_face_test~3 tem_face_test program~tem_face_test~3->module~tem_face_module program~tem_face_test tem_face_test program~tem_face_test->module~tem_face_module program~tem_face_test~2 tem_face_test program~tem_face_test~2->module~tem_face_module

Contents


Functions

private pure function tem_isRecvFace(facePos, faces) result(isRecvFace)

Subroutine to check if the current rank will receive information

Read more…

Arguments

TypeIntentOptionalAttributesName
integer, intent(in) :: facePos

The position of the face in faces to be checked.

type(tem_face_descriptor_type), intent(in) :: faces

The description of the faces.

Return Value integer

Returns 0 if no information has to be received from remote ranks before a compute step can be done. In case that information has to be received it will return a non-zero positive value. Furhtermore, it determines whether information for the left face value (left adjacent element) or the right face value (right adjacent element) has to be received. In the first case it will return tem_left and in the second case it will return tem_right. If something went wrong, e.g. an unknown combination occurs, it will return -1.

private pure function tem_isSendFace(facePos, faces) result(isSendFace)

Subroutine to check if the current rank will send information

Read more…

Arguments

TypeIntentOptionalAttributesName
integer, intent(in) :: facePos

The position of the face in faces to be checked.

type(tem_face_descriptor_type), intent(in) :: faces

The description of the faces.

Return Value integer

Returns 0 if no information has to be received from remote ranks before a compute step can be done. In case that information has to be received it will return a non-zero positive value. Furhtermore, it determines whether information for the left face value (left adjacent element) or the right face value (right adjacent element) has to be received. In the first case it will return tem_left and in the second case it will return tem_right. If something went wrong, e.g. an unknown combination occurs, it will return -1.

private function tem_isFromFinerFace(facePos, faces) result(isFromFiner)

Function to decide if a certain face is constructed from finer

Read more…

Arguments

TypeIntentOptionalAttributesName
integer, intent(in) :: facePos
type(tem_face_descriptor_type), intent(in) :: faces

Return Value integer

private function tem_reqDep_down(leftPrp, rightPrp) result(reqDep)

Function to check if a face with given left and right property requires a downward (coarse->fine) dependency.

Arguments

TypeIntentOptionalAttributesName
integer, intent(in) :: leftPrp

Left property of the face.

integer, intent(in) :: rightPrp

Right property of the face.

Return Value logical

Does the face requires a vertical dependency.

private function tem_reqDep_up(leftPrp, rightPrp) result(reqDep)

Function to check if a face with given left and right property requires a upward (fine->coarse) dependency.

Arguments

TypeIntentOptionalAttributesName
integer, intent(in) :: leftPrp

Left property of the face.

integer, intent(in) :: rightPrp

Right property of the face.

Return Value logical

Does the face requires a vertical dependency.

private function tem_melt_facePrp(firstPrp, secondPrp) result(meltedPrp)

Function to melt two properties together. The resulting property holds the union of firstPrp and secondPrp.

Arguments

TypeIntentOptionalAttributesName
integer, intent(in) :: firstPrp

The first property to be melted.

integer, intent(in) :: secondPrp

The second property to be melted.

Return Value integer

The resulting property (the union of the first and second property)

private function tem_get_elemPrp(levelDesc, elemPos) result(elemPrp)

Returns the property of a given element in the level descriptor

Read more…

Arguments

TypeIntentOptionalAttributesName
type(tem_levelDesc_type), intent(in) :: levelDesc

The level descriptor that contains the investigated element.

integer, intent(in) :: elemPos

The position of the element in the total list of the level descriptor.

Return Value integer

The element property of the element.


Subroutines

public subroutine tem_build_face_info(tree, boundary, commPattern, proc, faces, nEligibleChildren)

Routine to get face information for your whole mesh.

Read more…

Arguments

TypeIntentOptionalAttributesName
type(treelmesh_type), intent(inout) :: tree

Tree representation of the mesh to buid the faces for.

type(tem_BC_prop_type), intent(in) :: boundary

The boundaries of your simulation domain.

type(tem_commPattern_type), intent(in) :: commPattern

The communication pattern you use for the buffer.

type(tem_comm_env_type), intent(in) :: proc

Process description to use.

type(tem_face_type), intent(out), allocatable:: faces(:)

The created face descriptors (one for each level).

integer, intent(in) :: nEligibleChildren

The number of eligible children for the vertical face dependency

private subroutine tem_build_faceBuffers(minLevel, maxLevel, faces, levelDesc)

Arguments

TypeIntentOptionalAttributesName
integer, intent(in) :: minLevel

The minimum refinement level of your mesh.

integer, intent(in) :: maxLevel

The maximum refinement level of your mesh.

type(tem_face_type), intent(inout) :: faces(minLevel:maxLevel)

The communication pattern you want use for the buffer. The created face descriptor.

type(tem_levelDesc_type), intent(in) :: levelDesc(1:3,minLevel:maxLevel)

Dimension-by-dimension level descriptors

private subroutine tem_build_faceRecvBuffers(levelDesc, faces, buf)

Subroutine to build communication buffers for the faces the

Read more…

Arguments

TypeIntentOptionalAttributesName
type(tem_levelDesc_type), intent(in) :: levelDesc

Dimension-by-dimension level descriptor for the current level and direction.

type(tem_face_descriptor_type), intent(inout) :: faces

The communication pattern you want use for the buffer. The created face descriptor.

type(tem_communication_type), intent(out) :: buf(2)

The created receive buffer. Size is two, due to left and right limes of the face values. To access it use tem_left and tem_right.

private subroutine tem_build_faceSendBuffers(levelDesc, faces, buf)

Subroutine to build communication buffers for the faces the

Read more…

Arguments

TypeIntentOptionalAttributesName
type(tem_levelDesc_type), intent(in) :: levelDesc

Dimension-by-dimension level descriptor for the current level and direction.

type(tem_face_descriptor_type), intent(inout) :: faces

The communication pattern you want use for the buffer. The created face descriptor.

type(tem_communication_type), intent(out) :: buf(2)

The created send buffer. Size is two, due to left and right limes of the face values. To access it use tem_left and tem_right.

private subroutine tem_dimByDim_construction(tree, boundary, commPattern, proc, levelDescX, levelDescY, levelDescZ)

Creates dimension by dimension level descriptors.

Arguments

TypeIntentOptionalAttributesName
type(treelmesh_type), intent(inout) :: tree

Tree representation of the mesh.

type(tem_BC_prop_type), intent(in) :: boundary

The boundaries of your simulation domain

type(tem_commPattern_type), intent(in) :: commPattern

The communication pattern you use for the buffer.

type(tem_comm_env_type), intent(in) :: proc

Process description to use.

type(tem_levelDesc_type), intent(out), allocatable:: levelDescX(:)

Level descriptor for each spatial direction and each level of your mesh. The level descriptor have to be constructed with the dimension by dimension stencils (+1, 0, -1) for each spatial direction.

type(tem_levelDesc_type), intent(out), allocatable:: levelDescY(:)

Level descriptor for each spatial direction and each level of your mesh. The level descriptor have to be constructed with the dimension by dimension stencils (+1, 0, -1) for each spatial direction.

type(tem_levelDesc_type), intent(out), allocatable:: levelDescZ(:)

Level descriptor for each spatial direction and each level of your mesh. The level descriptor have to be constructed with the dimension by dimension stencils (+1, 0, -1) for each spatial direction.

private subroutine tem_create_levelDesc(tree, stencil, boundary, commPattern, levelDesc, proc)

Create a face level descriptor

Arguments

TypeIntentOptionalAttributesName
type(treelmesh_type), intent(inout) :: tree

Tree representation of the mesh.

type(tem_stencilHeader_type), intent(inout) :: stencil(1)

The stencil to create the level descriptor for.

type(tem_BC_prop_type), intent(in) :: boundary

The boundaries of your simulation domain

type(tem_commPattern_type), intent(in) :: commPattern

The communication pattern you use for the buffer.

type(tem_levelDesc_type), intent(out), allocatable:: levelDesc(:)

The created level descriptor.

type(tem_comm_env_type), intent(in) :: proc

Process description to use.

private subroutine tem_define_dimStencil(stencil, nElems, direction)

Creates the dimension-by-dimension stencils for a given spatial

Read more…

Arguments

TypeIntentOptionalAttributesName
type(tem_stencilHeader_type), intent(inout) :: stencil

The stencil layout to set.

integer, intent(in) :: nElems

The number of elements which share this stencil

integer, intent(in) :: direction

The spatial direction: 1 -> x direction 2 -> y direction 3 -> z direction

private subroutine tem_build_faceLists(minLevel, maxLevel, nEligibleChildren, faces)

Routine to generate the separated lists of faces for

Read more…

Arguments

TypeIntentOptionalAttributesName
integer, intent(in) :: minLevel

Minimum level of your mesh.

integer, intent(in) :: maxLevel

Maximum level of your mesh.

integer, intent(in) :: nEligibleChildren

The number of eligible children for the vertical face dependency

type(tem_face_type), intent(inout) :: faces(minLevel:maxLevel)

The created face descriptor.

private subroutine tem_build_fromFinerList(minLevel, maxLevel, nEligibleChildren, faces)

Extracts all the faces which will be computed by current rank.

Arguments

TypeIntentOptionalAttributesName
integer, intent(in) :: minLevel

The min refinement level of the mesh.

integer, intent(in) :: maxLevel

The max refinement level of the mesh.

integer, intent(in) :: nEligibleChildren

The number of eligible children for the vertical face dependency

type(tem_face_type), intent(inout) :: faces(minLevel:maxLevel)

The created face descriptor.

private subroutine tem_build_computeList(faces, nEligibleChildren)

Extracts all the faces which will be computed by current rank.

Arguments

TypeIntentOptionalAttributesName
type(tem_face_type), intent(inout) :: faces

The created face descriptor.

integer, intent(in) :: nEligibleChildren

The number of eligible children for the vertical face dependency

private subroutine tem_faceDep_vertical(minLevel, maxLevel, faces, nEligibleChildren)

Arguments

TypeIntentOptionalAttributesName
integer, intent(in) :: minLevel

Minimum level of your mesh.

integer, intent(in) :: maxLevel

Maximum level of your mesh.

type(tem_face_type), intent(inout) :: faces(minLevel:maxLevel)

Face descriptor where the faces will be appended to.

integer, intent(in) :: nEligibleChildren

The number of eligible children for the vertical face dependency

private subroutine tem_init_faceDep(faces, nEligibleChildren)

Subroutine to initialize the container for the vertical

Read more…

Arguments

TypeIntentOptionalAttributesName
type(tem_face_descriptor_type), intent(inout) :: faces

The face description in which you want to initialize the vertical face dependency container.

integer, intent(in) :: nEligibleChildren

The number of eligible children for the vertical face dependency

private subroutine tem_faceDep_verticalDown(coarseFaces, fineFaces, dir, nEligibleChildren)

Arguments

TypeIntentOptionalAttributesName
type(tem_face_descriptor_type), intent(inout) :: coarseFaces

Face description on the coarse level. The dependencies to the finer level will be appended to this face descriptor.

type(tem_face_descriptor_type), intent(in) :: fineFaces

Face description on the finer level.

integer, intent(in) :: dir

The spatial direction of the face to add the downward dependencies for. 1 --> x direction 2 --> y direction 3 --> z direction

integer, intent(in) :: nEligibleChildren

The number of eligible children for the vertical face dependency

private subroutine tem_faceDep_verticalUp(fineFaces, coarseFaces)

Arguments

TypeIntentOptionalAttributesName
type(tem_face_descriptor_type), intent(inout) :: fineFaces

Face description on the fine level. The dependencies to the coarser level will be appended to this face descriptor.

type(tem_face_descriptor_type), intent(in) :: coarseFaces

Face description on the coarse level.

private subroutine tem_getFace_prp(faces, facePos, leftPrp, rightPrp)

Subroutine returns left and right propery of a face.

Arguments

TypeIntentOptionalAttributesName
type(tem_face_descriptor_type), intent(in) :: faces

The face descriptor the face is located in.

integer, intent(in) :: facePos

The position of the face in the face descriptor.

integer, intent(out) :: leftPrp

The left property of the face.

integer, intent(out) :: rightPrp

The right property of the face.

private subroutine tem_addDep_down(coarseFacePos, coarseFaces, fineFaces, dir, nEligibleChildren)

Adds an downward face dependency (coarse->fine) to the coarse face

Read more…

Arguments

TypeIntentOptionalAttributesName
integer, intent(in) :: coarseFacePos

The position of the coarse face in coarseFaces you want to add the child dependencies for.

type(tem_face_descriptor_type), intent(inout) :: coarseFaces

The description of the faces on the coarse level. The dependency will be added to this face descriptor.

type(tem_face_descriptor_type), intent(in) :: fineFaces

The descriptor of the faces on the fine level.

integer, intent(in) :: dir

The spatial direction of the face to add the downward dependencies for. 1 --> x direction 2 --> y direction 3 --> z direction

integer, intent(in) :: nEligibleChildren

The number of eligible children for the vertical face dependency

private subroutine tem_addDep_up(fineFacePos, fineFaces, coarseFaces)

Adds an upward face dependency (fine->coarse) to the fine face

Read more…

Arguments

TypeIntentOptionalAttributesName
integer, intent(in) :: fineFacePos

The position of the face in the fineFaces descriptor you want to add.

type(tem_face_descriptor_type), intent(inout) :: fineFaces

The description of the faces on the fine level. The dependency will be added to this face descriptor.

type(tem_face_descriptor_type), intent(in) :: coarseFaces

The descriptor of the faces on the coarse level.

private subroutine tem_collect_faces(levelDesc, minLevel, maxLevel, faces)

Collects all the faces which have at least one neighboring

Read more…

Arguments

TypeIntentOptionalAttributesName
type(tem_levelDesc_type), intent(in) :: levelDesc(1:3,minLevel:maxLevel)

Level descriptor for each spatial direction and each level of your mesh.

integer, intent(in) :: minLevel

Minimum level of your mesh.

integer, intent(in) :: maxLevel

Maximum level of your mesh.

type(tem_face_type), intent(inout) :: faces(minLevel:maxLevel)

Face descriptor where the faces will be appended to.

private subroutine tem_init_faceList(faceList)

Subroutine to initialize the dyanmic content of a face list.

Arguments

TypeIntentOptionalAttributesName
type(tem_faceList_type), intent(inout) :: faceList

The face list to initialize.

private subroutine tem_get_faces(levelDesc, direction, faces)

Collects all the faces on a given level of the mesh in a given

Read more…

Arguments

TypeIntentOptionalAttributesName
type(tem_levelDesc_type), intent(in) :: levelDesc

Level descriptor of the level of the mesh you want to collect the faces for.

integer, intent(in) :: direction

The spatial direction to collect the faces for: 1 -> x direction 2 -> y direction 3 -> z direction

type(tem_face_descriptor_type), intent(inout) :: faces

Description of the faces on the current level.

private subroutine tem_get_faceNeigh(levelDesc, elemPos, dir, leftOrRight, neighId, neighPos)

Function to get the face neighbor of a certain element in the level descriptor. Even non-existing face neighbors can be handled by this routine.

Arguments

TypeIntentOptionalAttributesName
type(tem_levelDesc_type), intent(in) :: levelDesc

The level descriptor the element is located in.

integer, intent(in) :: elemPos

The element position in the level descriptor.

integer, intent(in) :: dir

The spatial direction of the face.

  1. -> x direction
  2. -> y direction
  3. -> z direction
integer, intent(in) :: leftOrRight

Find left or right face neighbor of the element. Use \ref tem_face_module::tem_left or \ref tem_face_module::tem_right.

integer(kind=long_k), intent(out) :: neighId

The treeid of the face neighbor element.

integer, intent(out) :: neighPos

Position of the neighbor in the total list, or 0 if element does not exist.

private subroutine tem_addFace(leftElemId, leftElemPos, leftElemPrp, rightElemId, rightElemPos, rightElemPrp, faces)

Adds a new face to the face description.

Read more…

Arguments

TypeIntentOptionalAttributesName
integer(kind=long_k), intent(in) :: leftElemId

Element id of the left element

integer, intent(in) :: leftElemPos

Position of the left element in the level descriptor's total list.

integer, intent(in) :: leftElemPrp

Properties of the left element.

integer(kind=long_k), intent(in) :: rightElemId

Element id of the right element

integer, intent(in) :: rightElemPos

Position of the right element in the level desriptor's total list.

integer, intent(in) :: rightElemPrp

Properties of the right element

type(tem_face_descriptor_type), intent(inout) :: faces

The face description the new face will be added to. If the face already exists in this face description. The existing face's property will be overwritten by the new ones.

private subroutine tem_extend_remotePrp(levelDesc, minLevel, maxLevel, faces)

Extends the properties of the faces which are marked as

Read more…

Arguments

TypeIntentOptionalAttributesName
type(tem_levelDesc_type), intent(in) :: levelDesc(1:3,minLevel:maxLevel)

Level descriptor for each level of your mesh.

integer, intent(in) :: minLevel

Minimum level of your mesh.

integer, intent(in) :: maxLevel

Maximum level of your mesh.

type(tem_face_type), intent(inout) :: faces(minLevel:maxLevel)

The face descriptor to be corrected.

private subroutine tem_extend_commFromFinerPrp(levelDesc, direction, faces)

Extend communication property for faces by the finer property if

Read more…

Arguments

TypeIntentOptionalAttributesName
type(tem_levelDesc_type), intent(in) :: levelDesc

Level descriptor of the level of the mesh you want to collect the faces for.

integer, intent(in) :: direction

The face direction.

type(tem_face_descriptor_type), intent(inout) :: faces

Description of the faces on the current level.

private subroutine tem_extend_commFromCoarserPrp(levelDesc, direction, faces)

Extend communication property for faces by the coarser property if

Read more…

Arguments

TypeIntentOptionalAttributesName
type(tem_levelDesc_type), intent(in) :: levelDesc

Level descriptor of the level of the mesh you want to collect the faces for.

integer, intent(in) :: direction

The face direction.

type(tem_face_descriptor_type), intent(inout) :: faces

Description of the faces on the current level.

private subroutine tem_appendFace_prp(faceId, prp, faces, prp_dir)

Attaches another property to a given face (from left or right). If the face does not exist this routine will not do anything.

Arguments

TypeIntentOptionalAttributesName
integer(kind=long_k), intent(in) :: faceId

The face identifier to append the face property to.

integer, intent(in) :: prp

The property to attach

type(tem_face_descriptor_type), intent(inout) :: faces

The face desriptor to update.

integer, intent(in) :: prp_dir

Attach the property to the left or right property of the face.