initibmnorm Subroutine

public subroutine initibmnorm(fname, solid_info)

Uses

  • proc~~initibmnorm~~UsesGraph proc~initibmnorm initibmnorm decomp_2d decomp_2d proc~initibmnorm->decomp_2d module~modglobal modglobal proc~initibmnorm->module~modglobal module~modmpi modmpi proc~initibmnorm->module~modmpi mpi mpi module~modmpi->mpi

Arguments

Type IntentOptional Attributes Name
character(len=11), intent(in) :: fname
type(solid_info_type), intent(inout) :: solid_info

Calls

proc~~initibmnorm~~CallsGraph proc~initibmnorm initibmnorm mpi_bcast mpi_bcast proc~initibmnorm->mpi_bcast zend zend proc~initibmnorm->zend zstart zstart proc~initibmnorm->zstart

Called by

proc~~initibmnorm~~CalledByGraph proc~initibmnorm initibmnorm proc~initibm initibm proc~initibm->proc~initibmnorm program~dalesurban DALESURBAN program~dalesurban->proc~initibm

Source Code

   subroutine initibmnorm(fname, solid_info)
     use modglobal, only : ifinput
     use modmpi,    only : myid, comm3d, mpierr
     use decomp_2d, only : zstart, zend

     character(11), intent(in) :: fname

     type(solid_info_type), intent(inout) :: solid_info

     logical :: lsolptsrank(solid_info%nsolpts)
     integer n, m

     character(80) chmess

     allocate(solid_info%solpts(solid_info%nsolpts,3))

     ! read u points
     if (myid == 0) then
       open (ifinput, file=fname)
       read (ifinput, '(a80)') chmess
       do n = 1, solid_info%nsolpts
         read (ifinput, *) solid_info%solpts(n,1), solid_info%solpts(n,2), solid_info%solpts(n,3)
       end do
       close (ifinput)
     end if

     call MPI_BCAST(solid_info%solpts, solid_info%nsolpts*3, MPI_INTEGER, 0, comm3d, mpierr)

     ! Determine whether points are on this rank
     solid_info%nsolptsrank = 0
     do n = 1, solid_info%nsolpts
       if ((solid_info%solpts(n,1) >= zstart(1) .and. solid_info%solpts(n,1) <= zend(1)) .and. &
           (solid_info%solpts(n,2) >= zstart(2) .and. solid_info%solpts(n,2) <= zend(2))) then
          lsolptsrank(n) = .true.
          solid_info%nsolptsrank = solid_info%nsolptsrank + 1
        else
          lsolptsrank(n) = .false.
       end if
     end do

     ! Store points on current rank - only loop through these points
     allocate(solid_info%solptsrank(solid_info%nsolptsrank))
     allocate(solid_info%solpts_loc(solid_info%nsolptsrank,3))
     m = 0
     do n = 1, solid_info%nsolpts
       if (lsolptsrank(n)) then
          m = m + 1
          solid_info%solptsrank(m) = n
          solid_info%solpts_loc(m,:) = (/solid_info%solpts(n,1),solid_info%solpts(n,2),solid_info%solpts(n,3)/)
       end if
     end do

     !write(*,*) "rank ", myid, " has ", solid_info%nsolptsrank, " solid points from ", fname

     deallocate(solid_info%solpts)

   end subroutine initibmnorm