direction. It attaches also two-sided properties to the face.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
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. |
subroutine tem_get_faces(levelDesc, direction, faces) ! -------------------------------------------------------------------------- !> Level descriptor of the level of the mesh you want to collect the !! faces for. type(tem_levelDesc_type), intent(in) :: levelDesc !> The spatial direction to collect the faces for: !! 1 -> x direction !! 2 -> y direction !! 3 -> z direction integer, intent(in) :: direction !> Description of the faces on the current level. type(tem_face_descriptor_type), intent(inout) :: faces ! -------------------------------------------------------------------------- integer :: iElem integer :: neighIndex integer :: neighPos integer :: elemPrp, neighPrp integer(kind=long_k) :: elemId, neighId ! -------------------------------------------------------------------------- ! Initialize the facelist (with a dynamic array for the faceIDs) call tem_init_faceList( faces%faceList ) ! Loop over all the elements (includes fluid, ghost and halos) have a look ! to the left and to the right face. ! Try to add both faces of the element to the list of all faces. elemLoop: do iElem = 1, levelDesc%nElems ! Get the current element property elemPrp = tem_get_elemPrp(levelDesc, iElem) elemId = levelDesc%total(iElem) ! Look for the left face of the current element. ! We use the compute stencil to determine the position of the neighboring ! left element. neighIndex = tem_left call tem_get_faceNeigh( levelDesc, iElem, direction, neighIndex, & & neighId, neighPos ) neighPrp = tem_get_elemPrp(levelDesc, neighPos) ! For the left face, the current element will be on its right side. call tem_addFace( leftElemId = neighId, & & leftElemPos = neighPos, & & leftElemPrp = neighPrp, & & rightElemId = elemId, & & rightElemPos = iElem, & & rightElemPrp = elemPrp, & & faces = faces ) ! Look for the right face of the current element. ! Get the right neighbor of the element. neighIndex = tem_right call tem_get_faceNeigh( levelDesc, iElem, direction, neighIndex, & & neighId, neighPos ) neighPrp = tem_get_elemPrp(levelDesc, neighPos) ! For the right face, the current element will be on its left side. call tem_addFace( leftElemId = elemId, & & leftElemPos = iElem, & & leftElemPrp = elemPrp, & & rightElemId = neighId, & & rightElemPos = neighPos, & & rightElemPrp = neighPrp, & & faces = faces ) end do elemLoop end subroutine tem_get_faces