This routine builds a description of the faces on each level. The resulting array of tem_face_type ranges from the minimal level of the mesh to the maximal level of the mesh.
\todo Remove all the unnecessary faces from the face descriptor
\todo JZ: this is very memory consuming since we store them in different arrays several times.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
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 |
subroutine tem_build_face_info( tree, boundary, commPattern, proc, faces, &
& nEligibleChildren )
! --------------------------------------------------------------------------
!> Tree representation of the mesh to buid the faces for.
type(treelmesh_type), intent(inout) :: tree
!> The boundaries of your simulation domain.
type(tem_bc_prop_type), intent(in) :: boundary
!> The communication pattern you use for the buffer.
type(tem_commpattern_type), intent(in) :: commPattern
!> Process description to use.
type(tem_comm_env_type), intent(in) :: proc
!> The number of eligible children for the vertical face dependency
integer, intent(in) :: nEligibleChildren
!> The created face descriptors (one for each level).
type(tem_face_type), allocatable, intent(out) :: faces(:)
! --------------------------------------------------------------------------
! 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) :: levelDesc( 1:3, tree%global%minLevel &
& :tree%global%maxLevel )
type(tem_levelDesc_type), allocatable :: levelDescX(:)
type(tem_levelDesc_type), allocatable :: levelDescY(:)
type(tem_levelDesc_type), allocatable :: levelDescZ(:)
integer :: iLevel
! --------------------------------------------------------------------------
allocate(faces(tree%global%minLevel:tree%global%maxLevel))
! create the level desriptor in a dim-by-dim manner.
write(logUnit(3),*)'Building dim-by-dim level descriptors for faces ...'
call tem_dimByDim_construction( tree, boundary, commPattern, proc, &
& levelDescX, levelDescY, levelDescZ )
do iLevel = tree%global%minLevel,tree%global%maxLevel
levelDesc(1,iLevel) = levelDescX(iLevel)
levelDesc(2,iLevel) = levelDescY(iLevel)
levelDesc(3,iLevel) = levelDescZ(iLevel)
end do
deallocate(levelDescX, levelDescY, levelDescZ)
! First we collect all the faces. A face is always identified by its left
! element tree id.
! We collect all the faces, even if only one element next to the face
! exists.
! Therefore, we run over fluids, ghostFromCoarser, ghostFromFiner and halo
! elements.
! We attach two-sided properties to the faces, too (depending on the element
! left and right to the face).
call tem_collect_faces( levelDesc, tree%global%minLevel, &
& tree%global%maxLevel, faces )
! Now, we have to check all the faces which have the communication property
! again.
! We have to check if the halo element next to it is refined/from coarser on
! the remote rank. So, we extend the faces' properties.
call tem_extend_remotePrp( levelDesc, tree%global%minLevel, &
& tree%global%maxLevel, faces )
! Now, we build the vertical dependency of the faces.
call tem_faceDep_vertical( tree%global%minLevel, tree%global%maxLevel, &
& faces, nEligibleChildren )
!> \todo Remove all the unnecessary faces from the face descriptor
! Build the buffers for sent/received faces.
call tem_build_faceBuffers( tree%global%minLevel, tree%global%maxLevel, &
& faces, levelDesc )
! Build the iterators for the faces. These ones will be used by the solvers.
call tem_build_faceLists(tree%global%minLevel, tree%global%maxLevel, &
& nEligibleChildren, faces )
!> \todo JZ: this is very memory consuming since we store them in different
!! arrays several times.
! Buffer the dimension by dimension descriptor
do iLevel = tree%global%minLevel,tree%global%maxLevel
faces(iLevel)%dimByDimDesc(1) = levelDesc(1,iLevel)
faces(iLevel)%dimByDimDesc(2) = levelDesc(2,iLevel)
faces(iLevel)%dimByDimDesc(3) = levelDesc(3,iLevel)
end do
end subroutine tem_build_face_info