readrestartfiles Subroutine

public subroutine readrestartfiles()

Uses

  • proc~~readrestartfiles~~UsesGraph proc~readrestartfiles modstartup::readrestartfiles module~modfields modfields proc~readrestartfiles->module~modfields module~modglobal modglobal proc~readrestartfiles->module~modglobal module~modinlet modinlet proc~readrestartfiles->module~modinlet module~modinletdata modinletdata proc~readrestartfiles->module~modinletdata module~modmpi modmpi proc~readrestartfiles->module~modmpi module~modsubgriddata modsubgriddata proc~readrestartfiles->module~modsubgriddata module~modsurfdata modsurfdata proc~readrestartfiles->module~modsurfdata module~modinlet->module~modinletdata mpi mpi module~modmpi->mpi

Arguments

None

Calls

proc~~readrestartfiles~~CallsGraph proc~readrestartfiles modstartup::readrestartfiles proc~zinterpolate1d modinlet::zinterpolate1d proc~readrestartfiles->proc~zinterpolate1d proc~zinterpolate2d modinlet::zinterpolate2d proc~readrestartfiles->proc~zinterpolate2d proc~zinterpolatet1d modinlet::zinterpolatet1d proc~readrestartfiles->proc~zinterpolatet1d proc~zinterpolatew1d modinlet::zinterpolatew1d proc~readrestartfiles->proc~zinterpolatew1d

Called by

proc~~readrestartfiles~~CalledByGraph proc~readrestartfiles modstartup::readrestartfiles proc~readinitfiles modstartup::readinitfiles proc~readinitfiles->proc~readrestartfiles proc~startup modstartup::startup proc~startup->proc~readinitfiles program~dalesurban DALESURBAN program~dalesurban->proc~startup

Contents

Source Code


Source Code

   subroutine readrestartfiles

      use modsurfdata, only:ustar, thlflux, qtflux, svflux, dudz, dvdz, dthldz, dqtdz, ps, thls, qts, thvs, oblav, &
         wtsurf
      use modfields, only:u0, v0, w0, thl0, qt0, ql0, ql0h, qtav, qlav, e120, dthvdz, presf, presh, sv0, mindist, wall, &
         uav, vav, wav, uuav, vvav, wwav, uvav, uwav, vwav, svav, thlav, thl2av, sv2av, pres0, svm, &
         svprof, viscratioav, thluav, thlvav, thlwav, svuav, svvav, svwav, presav, &
         uusgsav, vvsgsav, wwsgsav, uwsgsav, thlusgsav, thlwsgsav, svusgsav, svwsgsav, tkesgsav, &
         strain2av, nusgsav
      use modglobal, only:ib, ie, ih, jb, je, jh, kb, ke, kh, dtheta, dqt, dsv, startfile, timee, totavtime, runavtime, &
         iexpnr, ntimee, rk3step, ifinput, nsv, runtime, dt, cexpnr, lreadmean, lreadminl, &
         totinletav, lreadscal, ltempeq, dzf, numol, prandtlmoli
      use modmpi, only:cmyid, myid
      use modsubgriddata, only:ekm
      use modinlet, only:zinterpolate1d, zinterpolatet1d, zinterpolatew1d, zinterpolate2d
      use modinletdata, only:Uinl, Urec, Wrec, Utav, Tinl, Trec, linuf, linuh, &
         kbin, kein, lzinzsim, utaui, Ttav, ttaui

      real, dimension(ib:ie, jb:je, kb:ke)  ::  dummy3d
      real, dimension(ib:ie, kbin:kein)    ::  Utavin
      real, dimension(ib:ie, kbin:kein)    ::  Ttavin
      real, dimension(kbin:kein)          ::  Uinlin
      real, dimension(kbin:kein)          ::  Urecin
      real, dimension(kbin:kein)          ::  Tinlin
      real, dimension(kbin:kein)          ::  Trecin
      real, dimension(kbin:kein + 1)        ::  Wrecin
      character(50) :: name, name2, name4
      real dummy
      integer i, j, k, n
      !********************************************************************

      !    1.0 Read initfiles
      !-----------------------------------------------------------------
      name = startfile
      name(5:5) = 'd'
      name(15:17) = cmyid
      write (6, *) 'loading ', name
      open (unit=ifinput, file=name, form='unformatted', status='old')

      read (ifinput) (((mindist(i, j, k), i=ib, ie), j=jb, je), k=kb, ke)
      read (ifinput) ((((wall(i, j, k, n), i=ib, ie), j=jb, je), k=kb, ke), n=1, 5)
      read (ifinput) (((u0(i, j, k), i=ib - ih, ie + ih), j=jb - jh, je + jh), k=kb, ke + kh)
      read (ifinput) (((v0(i, j, k), i=ib - ih, ie + ih), j=jb - jh, je + jh), k=kb, ke + kh)
      read (ifinput) (((w0(i, j, k), i=ib - ih, ie + ih), j=jb - jh, je + jh), k=kb, ke + kh)
      read (ifinput) (((pres0(i, j, k), i=ib - ih, ie + ih), j=jb - jh, je + jh), k=kb, ke + kh)
      read (ifinput) (((thl0(i, j, k), i=ib - ih, ie + ih), j=jb - jh, je + jh), k=kb, ke + kh)
      read (ifinput) (((e120(i, j, k), i=ib - ih, ie + ih), j=jb - jh, je + jh), k=kb, ke + kh)
      read (ifinput) (((ekm(i, j, k), i=ib - ih, ie + ih), j=jb - jh, je + jh), k=kb, ke + kh)
      read (ifinput) (((qt0(i, j, k), i=ib - ih, ie + ih), j=jb - jh, je + jh), k=kb, ke + kh)
      read (ifinput) (((ql0(i, j, k), i=ib - ih, ie + ih), j=jb - jh, je + jh), k=kb, ke + kh)
      read (ifinput) (((ql0h(i, j, k), i=ib - ih, ie + ih), j=jb - jh, je + jh), k=kb, ke + kh)
      read (ifinput) timee, dt
      close (ifinput)
      write (6, *) 'finished loading ', name

      if ((nsv > 0) .and. (lreadscal)) then
         name(5:5) = 's'
         write (6, *) 'loading ', name
         open (unit=ifinput, file=name, form='unformatted')
         read (ifinput) ((((sv0(i, j, k, n), i=ib - ih, ie + ih), j=jb - jh, je + jh), k=kb, ke + kh), n=1, nsv)
         read (ifinput) timee
         close (ifinput)
         write (6, *) 'finished loading ', name
      elseif ((nsv > 0) .and. (.not. lreadscal)) then
         sv0 = 0.
         svprof = 0.
      end if

      ! read mean variables if asked for by lreadmean
      name2 = 'means   .'
      name2(6:8) = cmyid
      name2(10:12) = cexpnr
      if (lreadmean) then
         write (6, *) 'Reading meansXXX.XXX, proc = ', myid
         open (unit=ifinput, file=name2, form='unformatted')
         read (ifinput) totavtime, nsv
         read (ifinput) (((uav(i, j, k), i=ib - ih, ie + ih), j=jb - jh, je + jh), k=kb - kh, ke + kh)
         read (ifinput) (((vav(i, j, k), i=ib - ih, ie + ih), j=jb - jh, je + jh), k=kb - kh, ke + kh)
         read (ifinput) (((wav(i, j, k), i=ib - ih, ie + ih), j=jb - jh, je + jh), k=kb - kh, ke + kh)
         read (ifinput) (((thlav(i, j, k), i=ib - ih, ie + ih), j=jb - jh, je + jh), k=kb - kh, ke + kh)
         read (ifinput) (((qtav(i, j, k), i=ib - ih, ie + ih), j=jb - jh, je + jh), k=kb - kh, ke + kh)
         read (ifinput) (((qlav(i, j, k), i=ib - ih, ie + ih), j=jb - jh, je + jh), k=kb - kh, ke + kh)
         read (ifinput) (((presav(i, j, k), i=ib - ih, ie + ih), j=jb - jh, je + jh), k=kb - kh, ke + kh)
         read (ifinput) ((((svav(i, j, k, n), i=ib - ih, ie + ih), j=jb - jh, je + jh), k=kb - kh, ke + kh), n=1, nsv)
         read (ifinput) (((viscratioav(i, j, k), i=ib - ih, ie + ih), j=jb - jh, je + jh), k=kb - kh, ke + kh)
         read (ifinput) (((uuav(i, j, k), i=ib - ih, ie + ih), j=jb - jh, je + jh), k=kb - kh, ke + kh)
         read (ifinput) (((vvav(i, j, k), i=ib - ih, ie + ih), j=jb - jh, je + jh), k=kb - kh, ke + kh)
         read (ifinput) (((wwav(i, j, k), i=ib - ih, ie + ih), j=jb - jh, je + jh), k=kb - kh, ke + kh)
         read (ifinput) (((thl2av(i, j, k), i=ib - ih, ie + ih), j=jb - jh, je + jh), k=kb - kh, ke + kh)
         read (ifinput) ((((sv2av(i, j, k, n), i=ib - ih, ie + ih), j=jb - jh, je + jh), k=kb - kh, ke + kh), n=1, nsv)
         read (ifinput) (((uvav(i, j, k), i=ib, ie + ih), j=jb, je + jh), k=kb, ke)
         read (ifinput) (((uwav(i, j, k), i=ib, ie + ih), j=jb, je), k=kb, ke + kh)
         read (ifinput) (((vwav(i, j, k), i=ib, ie), j=jb, je + jh), k=kb, ke + kh)
         read (ifinput) (((thluav(i, j, k), i=ib, ie), j=jb, je), k=kb, ke)
         read (ifinput) (((thlvav(i, j, k), i=ib, ie), j=jb, je + jh), k=kb, ke)
         read (ifinput) (((thlwav(i, j, k), i=ib, ie), j=jb, je), k=kb, ke + kh)
         read (ifinput) ((((svuav(i, j, k, n), i=ib, ie), j=jb, je), k=kb, ke), n=1, nsv)
         read (ifinput) ((((svvav(i, j, k, n), i=ib, ie), j=jb, je + jh), k=kb, ke), n=1, nsv)
         read (ifinput) ((((svwav(i, j, k, n), i=ib, ie), j=jb, je), k=kb, ke + kh), n=1, nsv)
         close (ifinput)
         write (6, *) 'Total averaging time so far: ', totavtime

         ! read <x'y'>_SGS to file.
         name2 = 'SGS__   .'
         name2(6:8) = cmyid
         name2(10:12) = cexpnr
         open (unit=ifinput, file=name2, form='unformatted')
         read (ifinput) dummy, dummy
         read (ifinput) (((uusgsav(i, j, k), i=ib - ih, ie + ih), j=jb - jh, je + jh), k=kb - kh, ke + kh)
         read (ifinput) (((vvsgsav(i, j, k), i=ib - ih, ie + ih), j=jb - jh, je + jh), k=kb - kh, ke + kh)
         read (ifinput) (((wwsgsav(i, j, k), i=ib - ih, ie + ih), j=jb - jh, je + jh), k=kb - kh, ke + kh)
         read (ifinput) (((uwsgsav(i, j, k), i=ib, ie + ih), j=jb, je), k=kb, ke + kh)
         read (ifinput) (((dummy3d(i, j, k), i=ib, ie), j=jb, je), k=kb, ke) ! this is dissresav, which will be computed using other mean quantities
         read (ifinput) (((tkesgsav(i, j, k), i=ib, ie), j=jb, je), k=kb, ke)
         read (ifinput) (((dummy3d(i, j, k), i=ib, ie), j=jb, je), k=kb, ke) ! this is disssgsav, which will be computed using other mean quantities
         read (ifinput) (((strain2av(i, j, k), i=ib, ie), j=jb, je), k=kb, ke) ! <SijSij> (NOT <Sij><Sij> !!) (average over time)
         read (ifinput) (((nusgsav(i, j, k), i=ib, ie), j=jb, je), k=kb, ke) ! <nu_sgs> (average over time)
         read (ifinput) (((thlusgsav(i, j, k), i=ib, ie + ih), j=jb, je), k=kb, ke)
         read (ifinput) (((thlwsgsav(i, j, k), i=ib, ie), j=jb, je), k=kb, ke + kh)
         read (ifinput) ((((svusgsav(i, j, k, n), i=ib, ie + ih), j=jb, je), k=kb, ke), n=1, nsv)
         read (ifinput) ((((svwsgsav(i, j, k, n), i=ib, ie), j=jb, je), k=kb, ke + kh), n=1, nsv)
         close (ifinput)

      end if

      ! read mean profiles for inlet generator
      if (lreadminl) then
         if (.not. lzinzsim) then
            name4 = 'meaninlet.   '
            name4(11:13) = cexpnr
            open (unit=ifinput, file=name4, form='unformatted')
            read (ifinput) totinletav ! interval of time-average
            read (ifinput) (Uinlin(k), k=kbin, kein)
            read (ifinput) (Urecin(k), k=kbin, kein)
            read (ifinput) (Wrecin(k), k=kbin, kein + 1)
            read (ifinput) ((Utavin(i, k), i=ib, ie), k=kbin, kein)
            close (ifinput)

            call zinterpolate1d(Uinlin, Uinl) ! interpolate inlet profile to zgrid
            call zinterpolate1d(Urecin, Urec)
            call zinterpolatew1d(Wrecin, Wrec)
            call zinterpolate2d(Utavin, Utav)

            if (ltempeq) then
               name4 = 'tempinlet.   '
               name4(11:13) = cexpnr
               open (unit=ifinput, file=name4, form='unformatted')
               read (ifinput) totinletav ! interval of time-average
               read (ifinput) (Tinlin(k), k=kbin, kein)
               read (ifinput) (Trecin(k), k=kbin, kein)
               read (ifinput) ((Ttavin(i, k), i=ib, ie), k=kbin, kein)
               close (ifinput)

               call zinterpolatet1d(Tinlin, Tinl)
               call zinterpolatet1d(Trecin, Trec)
               call zinterpolate2d(Ttavin, Ttav)
            end if ! ltempeq

         else !lzinzsim=.true. -> inlet grid equals sim grid
            name4 = 'meaninlet.   '
            name4(11:13) = cexpnr
            open (unit=ifinput, file=name4, form='unformatted')
            read (ifinput) totinletav ! interval of time-average
            read (ifinput) (Uinl(k), k=kb, ke)
            read (ifinput) (Urec(k), k=kb, ke)
            read (ifinput) (Wrec(k), k=kb, ke + 1)
            read (ifinput) ((Utav(i, k), i=ib, ie), k=kb, ke)

            close (ifinput)

            if (ltempeq) then
               name4 = 'tempinlet.   '
               name4(11:13) = cexpnr
               open (unit=ifinput, file=name4, form='unformatted')
               read (ifinput) totinletav ! interval of time-average
               read (ifinput) (Tinl(k), k=kb, ke)
               read (ifinput) (Trec(k), k=kb, ke)
               read (ifinput) ((Ttav(i, k), i=ib, ie), k=kb, ke)

               close (ifinput)
            end if ! ltempeq
         end if ! lzinzsim

         utaui = sqrt(abs(2*numol*Uinl(kb)/dzf(kb))) ! average streamwise friction at inlet (need for first time step)
         if (ltempeq) then
            ttaui = numol*prandtlmoli*2.*(Tinl(kb) - thls)/(dzf(kb)*utaui)
         end if
      end if !(lreadminl)

   end subroutine readrestartfiles