detfreestream Subroutine

public subroutine detfreestream(freestream)

Uses

  • proc~~detfreestream~~UsesGraph proc~detfreestream modforces::detfreestream module~modfields modfields proc~detfreestream->module~modfields module~modglobal modglobal proc~detfreestream->module~modglobal module~modmpi modmpi proc~detfreestream->module~modmpi mpi mpi module~modmpi->mpi

Arguments

Type IntentOptional Attributes Name
real, intent(out) :: freestream

Calls

proc~~detfreestream~~CallsGraph proc~detfreestream modforces::detfreestream mpi_allreduce mpi_allreduce proc~detfreestream->mpi_allreduce

Called by

proc~~detfreestream~~CalledByGraph proc~detfreestream modforces::detfreestream proc~fixuinf1 modforces::fixuinf1 proc~fixuinf1->proc~detfreestream proc~fixuinf2 modforces::fixuinf2 proc~fixuinf2->proc~detfreestream program~dalesurban DALESURBAN program~dalesurban->proc~fixuinf1 program~dalesurban->proc~fixuinf2

Contents

Source Code


Source Code

  subroutine detfreestream(freestream)
    use modglobal, only : ib,ie,jb,je,kb,ke,kh,dxf,xh,dt,&
                          Uinf,lvinf,dy
    use modfields, only : u0,dpdxl,dgdt,dpdx,v0
    use modmpi, only    : myid,comm3d,mpierr,mpi_sum,my_real,nprocs
    implicit none
    real, intent(out) :: freestream

    real  utop,vtop,dum
    integer i,j

    if (lvinf) then
        vtop = 0.
        do j =jb,je
          do i =ib,ie
            vtop = vtop + 0.5*(v0(i,j,ke)+u0(i,j+1,ke))*dy
          end do
        end do
        vtop = vtop / ( (je-jb+1)*(xh(ie+1)-xh(ib) ) )
        call MPI_ALLREDUCE(vtop,   dum,1,MY_REAL,MPI_SUM,comm3d,mpierr)
        freestream = dum / nprocs
    else
        utop = 0.
        do j =jb,je
          do i =ib,ie
            !dum=0.5*(u0(i,j,ke)+u0(i+1,j,ke))*dxf(i)
            !utop = utop + dum
            utop = utop + 0.5*(u0(i,j,ke)+u0(i+1,j,ke))*dxf(i)
          end do
        end do
        utop = utop / ( (je-jb+1)*(xh(ie+1)-xh(ib) ) )
        call MPI_ALLREDUCE(utop,   dum,1,MY_REAL,MPI_SUM,comm3d,mpierr)
        freestream = dum / nprocs
    !write(*,*) "myid,utop,dum,freestream",myid,utop,dum,freestream
    end if

  end subroutine detfreestream