Report the actually existing elements, which were requested as halos from remote
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
type(tem_levelDesc_type), | intent(inout) | :: | levelDesc |
the level descriptor of specific level |
||
type(tem_communication_type), | intent(inout) | :: | sendbuffer |
send and receive communication buffer type |
||
type(tem_communication_type), | intent(inout) | :: | recvbuffer |
send and receive communication buffer type |
||
type(tem_comm_env_type), | intent(in) | :: | proc |
Process description to use. |
||
type(tem_commPattern_type) | :: | commPattern |
communication pattern |
|||
type(tem_stencilHeader_type) | :: | computeStencil(:) |
array of all stencils used in the simulation |
subroutine redefine_halos( levelDesc, sendbuffer, recvbuffer, proc, & & commPattern, computeStencil ) ! --------------------------------------------------------------------------- !> the level descriptor of specific level type(tem_levelDesc_type), intent(inout) :: levelDesc !> send and receive communication buffer type type(tem_communication_type), intent(inout) :: sendbuffer, recvbuffer !> Process description to use. type(tem_comm_env_type), intent(in) :: proc !> communication pattern type(tem_commPattern_type) :: commPattern !> array of all stencils used in the simulation type(tem_stencilHeader_type) :: computeStencil(:) ! --------------------------------------------------------------------------- integer :: iProc, tPos, iElem, nElems integer :: elemCount integer(kind=long_k), allocatable :: recvBuf(:) integer :: message_flag logical :: wasAdded ! --------------------------------------------------------------------------- message_flag = 3 allocate( recvBuf(sum(recvbuffer%nElemsProc)) ) elemCount = 0 do iProc = 1, recvbuffer%nProcs nElems = recvbuffer%nElemsProc(iProc) call init( me = recvBuffer%elemPos( iProc ), & & length = nElems ) do iElem=1,nElems elemCount = elemCount + 1 call append( me = recvbuffer%elemPos( iProc ), & & val = elemCount ) end do call commPattern%initBuf_long( me = recvBuffer%buf_long( iProc ), & & pos = recvBuffer%elemPos( iProc )%val, & & nVals = recvBuffer%nElemsProc( iProc )) end do if (tem_logging_isActive(main_debug%logger, 6)) then write(dbgUnit(0),*) ' SEND BUFFER ' call tem_comm_dumpType( me = sendbuffer, & & nUnit = dbgUnit(6) ) write(dbgUnit(0),*) ' ' write(dbgUnit(0),*) ' RECV BUFFER ' call tem_comm_dumpType( me = recvbuffer, & & nUnit = dbgUnit(6) ) end if ! --------------------------------------------------------------------- ! initialize send buffer do iProc = 1, sendbuffer%nProcs call commPattern%initBuf_long( me = sendBuffer%buf_long( iProc ), & & pos = sendBuffer%elemPos( iProc )%val, & & nVals = sendBuffer%nElemsProc( iProc )) end do call commPattern%exchange_long( send = sendbuffer, & & recv = recvbuffer, & & state = recvBuf, & & message_flag = message_flag, & & comm = proc%comm, & & send_state = levelDesc%elem%tID%val ) ! Now I received all the halos, which are really existing ! and required on the remote processes ! --------------------------------------------------------------------- ! --------------------------------------------------------------------- ! Add all received elements ! Set the eType accordingly or leave as nonExisting if it was not received elemCount = 0 do iProc = 1, recvbuffer%nProcs nElems = recvBuffer%nElemsProc( iProc ) call commPattern%finBuf_long(me = recvbuffer%buf_long(iProc)) call init( me = recvBuffer%elemPos( iProc ), & & length = nElems ) do iElem = 1, nElems elemCount = elemCount + 1 call append( me = levelDesc%elem, & & tID = recvBuf( elemCount ), & & sourceProc = recvbuffer%proc( iProc )+1, & & eType = eT_halo, & & pos = tPos, & & wasAdded = wasAdded ) if ( .not. wasAdded ) then ! If it was not added, it was there before. ! Then simply set the element type to halo (was initialized with ! eT_nonExistingElem above ) if ( levelDesc%elem%eType%val( tPos ) == eT_nonExisting ) then call changeType( levelDesc%elem, tPos, eT_halo ) end if end if call append( me = recvbuffer%elemPos( iProc ), & & val = tPos ) end do call commPattern%initBuf_long( me = recvbuffer%buf_long( iProc ), & & pos = recvbuffer%elemPos( iProc )%val, & & nVals = recvBuffer%nElemsProc( iProc )) end do deallocate(recvBuf) ! --------------------------------------------------------------------- ! communicate the properties of the halos call commPattern%exchange_long( & & send = sendbuffer, & & recv = recvbuffer, & & state = levelDesc%elem%property%val(:), & & message_flag = 5, & & comm = proc%comm ) ! communicate the stencil neighbors ! the stencil which includes the boundary information if ( recvbuffer%nProcs > 0 .or. sendbuffer%nProcs > 0 ) then ! communicate only fluid stencil which is defined for all ! elements in elem list including halo elements. ! In our case, fluid stencil is 1st Stencil so setting iStencil=1 call tem_stencil_communicate( send = sendbuffer, & & recv = recvbuffer, & & elem = levelDesc%elem, & & proc = proc, & & commPattern = commPattern, & & computeStencil = computeStencil(1), & & iStencil = 1 ) end if ! Deallocate long type buffer do iProc = 1, sendBuffer%nProcs call commPattern%finbuf_long( sendBuffer%buf_long(iProc) ) end do do iProc = 1, recvBuffer%nProcs call commPattern%finbuf_long( recvBuffer%buf_long(iProc) ) end do write(logUnit(5),*) ' Done communicating stencil buffer' ! --------------------------------------------------------------------- end subroutine redefine_halos