inittimedep Subroutine

public subroutine inittimedep()

Uses

  • proc~~inittimedep~~UsesGraph proc~inittimedep inittimedep module~modfields modfields proc~inittimedep->module~modfields module~modglobal modglobal proc~inittimedep->module~modglobal module~modibmdata modibmdata proc~inittimedep->module~modibmdata module~modmpi modmpi proc~inittimedep->module~modmpi decomp_2d decomp_2d module~modfields->decomp_2d mpi mpi module~modmpi->mpi

Arguments

None

Calls

proc~~inittimedep~~CallsGraph proc~inittimedep inittimedep mpi_bcast mpi_bcast proc~inittimedep->mpi_bcast proc~timedep timedep proc~inittimedep->proc~timedep proc~timedeplw timedeplw proc~timedep->proc~timedeplw proc~timedepnudge timedepnudge proc~timedep->proc~timedepnudge proc~timedepsurf timedepsurf proc~timedep->proc~timedepsurf proc~timedepsw timedepsw proc~timedep->proc~timedepsw

Called by

proc~~inittimedep~~CalledByGraph proc~inittimedep inittimedep program~dalesurban DALESURBAN program~dalesurban->proc~inittimedep

Source Code

  subroutine inittimedep
    use modmpi,    only :myid,my_real,mpi_logical,mpierr,comm3d
    use modglobal, only :cexpnr,kb,ke,kh,kmax,ifinput,runtime,zf,skyLW,nfcts
    use modibmdata, only : bctfxm, bctfxp, bctfym, bctfyp, bctfz!, bctfzf
    use modfields, only: thlprof
    !use initfac, only : netsw !Should probably be moved to somewhere else

    implicit none

    character (80):: chmess
    character (1) :: chmess1
    integer :: k,t,n,ierr
    real :: dummyr
    real, allocatable, dimension (:) :: height

    ltimedep = (ltimedepsurf .or. ltimedepnudge) .or. (ltimedeplw .or. ltimedepsw)

    if (.not. ltimedep) return

    if (ltimedepsurf) then
      allocate(timeflux (1:ntimedepsurf))
      allocate(bctfxmt (1:ntimedepsurf))
      allocate(bctfxpt (1:ntimedepsurf))
      allocate(bctfymt (1:ntimedepsurf))
      allocate(bctfypt (1:ntimedepsurf))
      allocate(bctfzt (1:ntimedepsurf))
      !allocate(bctfzft (1:ntimedepsurf))

      timeflux = 0.
      bctfxmt = bctfxm
      bctfxpt = bctfxp
      bctfymt = bctfym
      bctfypt = bctfyp
      bctfzt = bctfz
      !bctfzft = bctfzf

      if (myid==0) then
        open(ifinput,file='timedepsurf.inp.'//cexpnr)
        read(ifinput,'(a80)') chmess
        !write(6,*) chmess
        read(ifinput,'(a80)') chmess
        !write(6,*) chmess

        !--- load fluxes---
        !t    = 1
        ierr = 0
        do t = 1,ntimedepsurf
          read(ifinput,*, iostat = ierr) timeflux(t), bctfxmt(t), bctfxpt(t), bctfymt(t), bctfypt(t), bctfzt(t)!, bctfzft(t)
          !write(*,*) t, timeflux(t), bctfxmt(t), bctfxpt(t), bctfymt(t), bctfypt(t), bctfzt(t)!, bctfzft(t)
          !if (ierr < 0) then
            !stop 'STOP: No time dependend data for end of run (surface fluxes)'
          !end if
        end do
        !if(timeflux(1)>runtime) then
         !write(6,*) 'Time dependent surface variables do not change before end of simulation'
         !ltimedepsurf=.false.
        !endif
        ! flush to the end of fluxlist
        !do while (ierr ==0)
          !read (ifinput,*,iostat=ierr) dummyr
        !end do
        !backspace (ifinput)
        !close(ifinput)

      end if !myid==0

      call MPI_BCAST(timeflux, ntimedepsurf, MY_REAL, 0, comm3d, mpierr)
      call MPI_BCAST(bctfxmt , ntimedepsurf, MY_REAL, 0, comm3d, mpierr)
      call MPI_BCAST(bctfxpt , ntimedepsurf, MY_REAL, 0, comm3d,mpierr)
      call MPI_BCAST(bctfymt , ntimedepsurf, MY_REAL, 0, comm3d,mpierr)
      call MPI_BCAST(bctfypt , ntimedepsurf, MY_REAL, 0, comm3d,mpierr)
      call MPI_BCAST(bctfzt  , ntimedepsurf, MY_REAL, 0, comm3d,mpierr)
      !call MPI_BCAST(bctfzft , ntimedepsurf, MY_REAL, 0, comm3d,mpierr)

    end if

    if (ltimedepnudge) then

      allocate(timenudge (1:ntimedepnudge))
      allocate(height    (kb:ke+kh))
      allocate(uproft    (kb:ke+kh,ntimedepnudge))
      allocate(vproft    (kb:ke+kh,ntimedepnudge))
      allocate(thlproft  (kb:ke+kh,ntimedepnudge))
      allocate(qtproft   (kb:ke+kh,ntimedepnudge))

      timenudge = 0
      thlproft = 0
      qtproft = 0
      uproft = 0
      vproft = 0

      if (myid == 0) then
        !---load nudging profiles----
        open(ifinput,file='timedepnudge.inp.'//cexpnr)
        read(ifinput,'(a80)') chmess
        !write(6,*) chmess

        !t = 0
        do t = 1,ntimedepnudge
          !t = t + 1
          chmess1 = "#"
          ierr = 1 ! not zero
          do while (.not.(chmess1 == "#" .and. ierr == 0)) !search for the next line consisting of "# time", from there onwards the profiles will be read
            read(ifinput,*,iostat=ierr) chmess1, timenudge(t)
            !if (ierr < 0) then
              !stop 'STOP: No time dependend data for end of run'
            !end if
          end do

          !write (*,*) 'timenudge = ',timenudge(t)
          !write(*,*) 'Nudging profiles'
          do k=kb,ke
            read (ifinput,*) &
            height   (k)  , &
            thlproft (k,t), &
            qtproft  (k,t), &
            uproft   (k,t), &
            vproft   (k,t)

            !write(*,*) height(k), thlproft (k,t), qtproft(k,t), uproft(k,t), vproft(k,t)
          end do
        end do
      end if !myid == 0

      call MPI_BCAST(timenudge, ntimedepnudge, MY_REAL, 0, comm3d, mpierr)
      call MPI_BCAST(thlproft, (kmax+1)*ntimedepnudge,MY_REAL,0,comm3d,mpierr)
      call MPI_BCAST(qtproft,  (kmax+1)*ntimedepnudge,MY_REAL,0,comm3d,mpierr)
      call MPI_BCAST(uproft,   (kmax+1)*ntimedepnudge,MY_REAL,0,comm3d,mpierr)
      call MPI_BCAST(vproft,   (kmax+1)*ntimedepnudge,MY_REAL,0,comm3d,mpierr)

      deallocate(height)

    end if

    if (ltimedeplw) then
      allocate(timelw (1:ntimedeplw))
      allocate(skyLWt (1:ntimedeplw))

      ! Read longwave
      timelw = 0.
      skyLWt = skyLW

      if (myid==0) then
       open(ifinput,file='timedeplw.inp.'//cexpnr)
       read(ifinput,'(a80)') chmess
       !write(6,*) chmess
       read(ifinput,'(a80)') chmess
       !write(6,*) chmess

       !--- load fluxes---
       !t    = 1
       ierr = 0
       do t = 1,ntimedeplw
          read(ifinput,*, iostat = ierr) timelw(t), skyLWt(t)
          !write(*,*) t, timelw(t), skyLWt(t)
          !if (ierr < 0) then
            !stop 'STOP: No time dependend data for end of run (surface fluxes)'
          !end if
       end do
       !if(timelw(1)>runtime) then
         !write(6,*) 'Time dependent surface variables do not change before end of simulation'
         !ltimedeplw=.false.
       !endif
       ! flush to the end of fluxlist
       !do while (ierr ==0)
          !read (ifinput,*,iostat=ierr) dummyr
       !end do
       !backspace (ifinput)
       !close(ifinput)

      end if ! myid = 0

      call MPI_BCAST(timelw, ntimedeplw, MY_REAL, 0, comm3d, mpierr)
      call MPI_BCAST(skyLWt, ntimedeplw, MY_REAL, 0, comm3d, mpierr)

    end if !ltimedeplw

    if (ltimedepsw) then
      allocate(timesw(1:ntimedepsw))
      allocate(netswt(1:nfcts, 1:ntimedepsw))

      timesw = 0.
      netswt = 0.

      if (myid == 0) then
        ! Read shortwave
        open (ifinput,file='timedepsw.inp.'//cexpnr)
        read (ifinput,'(a80)') chmess ! first line is a description of the file
        read (ifinput, *) (timesw(t), t=1,ntimedepsw) ! second line is the times

        do n = 1,nfcts
         read (ifinput, *) (netswt(n,t), t=1,ntimedepsw)
        end do

        !write(*,*) "read timedepsw"

      end if !myid==0

      call MPI_BCAST(timesw, ntimedepsw, MY_REAL, 0, comm3d, mpierr)
      call MPI_BCAST(netswt, ntimedepsw*nfcts, MY_REAL, 0, comm3d, mpierr)

    end if ! timedepsw

    call MPI_BCAST(ltimedepsurf ,1,MPI_LOGICAL,0,comm3d,mpierr)
    call MPI_BCAST(ltimedepnudge,1,MPI_LOGICAL,0,comm3d,mpierr)
    call MPI_BCAST(ltimedeplw,1,MPI_LOGICAL,0,comm3d,mpierr)
    call MPI_BCAST(ltimedepsw,1,MPI_LOGICAL,0,comm3d,mpierr)

    call timedep

  end subroutine inittimedep