Merge pull request #22 from wirc-sjsu/develop-w21
[WRF-Fire-merge.git] / phys / module_sf_scmflux.F
blobe060387b1ae2551b35acb47ad3334d382648bf3b
5 module module_sf_scmflux
6 contains
8 !-------------------------------------------------------------------
10    subroutine scmflux(u3d, v3d, t3d, qv3d, p3d, dz8w,                           &
11                      cp, rovcp, xlv, psfc, cpm, xland,                          &
12                      psim, psih, hfx, qfx, lh, tsk, flhc, flqc,                 &
13                      znt, gz1oz0, wspd,                                         &
14                      julian_in, karman, p1000mb,                                &
15                      itimestep,chklowq,                                          &
16                      ids, ide, jds, jde, kds, kde,                              &
17                      ims, ime, jms, jme, kms, kme,                              &
18                      its, ite, jts, jte, kts, kte   )
19 !-------------------------------------------------------------------
20       implicit none
21 !-------------------------------------------------------------------
23    integer, intent(in)   ::                       ids, ide, jds, jde, kds, kde, &
24                                                   ims, ime, jms, jme, kms, kme, &
25                                        its, ite, jts, jte, kts, kte, itimestep        
26 !   
27    real, intent(in)      ::         cp, rovcp, xlv, julian_in, karman, p1000mb
29    real, dimension( ims:ime, kms:kme, jms:jme )                               , &
30             intent(in)   ::                                                u3d, &
31                                                                            v3d, &
32                                                                            t3d, &
33                                                                           qv3d, &
34                                                                            p3d, &
35                                                                           dz8w
36    real, dimension( ims:ime, jms:jme )                                        , &
37             intent(in)   ::                                               psfc, &
38                                                                          xland, &
39                                                                           flhc, &
40                                                                           flqc 
42    real, dimension( ims:ime, jms:jme )                                        , &
43             intent(inout)::                                                cpm, &
44                                                                            znt, &
45                                                                         gz1oz0, &
46                                                                           wspd, &
47                                                                           psim, &
48                                                                           psih, &
49                                                                            hfx, &
50                                                                            qfx, &
51                                                                             lh, &
52                                                                            tsk,&
53                                                                         chklowq           
54 ! local vars
56    integer, parameter    ::                                       n_max = 1200
57    integer               ::                                 i, j, n, nm, nt, m
58    real, parameter       ::                                           den = 1.
59    real                  ::                         julian_s, julian_e, fc_int, &
60                                                             fm, fh, ch, dtdiff
61    real, dimension( 1:n_max ) ::                     fc_qfx, fc_hfx, fc_julian !JP 0 ->1
62    real                       ::                     qfx_interp,hfx_interp ! JP
63    real, dimension( its:ite, jts:jte) ::                                   u2d, &
64                                                                            v2d, &
65                                                                            t2d, &
66                                                                           qv2d, &
67                                                                            p2d, &
68                                                                         dz8w1d, &
69                                                                             za, &
70                                                                            thx, &
71                                                                           thgb
72    logical               ::                                        end_of_file
74 !-----open scmflx_bdy and read the julian_s, julian_e, fc_int
76    open(unit=11, file='scmflx_bdy', form='formatted', status='old')
77    print*,'scmflx_bdy' 
78    read(11,*) julian_s, julian_e, fc_int
80      end_of_file = .false.
81      n=1
82      do while (.not. end_of_file)
83        read(11,*,end=100) fc_hfx(n), fc_qfx(n)
84        fc_julian(n)=julian_s+(n-1)*fc_int/86400.
85        n=n+1
86        go to 110
87  100   end_of_file = .true.  
88  110   continue
89      enddo
90      nt=n-1
91    close(11)
93 !-----linear interpolation of the fluxes for each time step
95    do n=1,nt 
96      if (julian_in.ge.fc_julian(n) .and. julian_in.lt.fc_julian(n+1)) then
97        qfx_interp= fc_qfx(n)                                                      &
98                 +(fc_qfx(n+1)-fc_qfx(n))*(julian_in-fc_julian(n))/(fc_int/86400.)
99        hfx_interp= fc_hfx(n)                                                      &
100                 +(fc_hfx(n+1)-fc_hfx(n))*(julian_in-fc_julian(n))/(fc_int/86400.)
101      endif
102    enddo
104 !-----compute surface moisture and heat fluxes, in the unit of [W m-2]
107 !-----compute skin temperature
109    do j=jts,jte
110      do i=its,ite
111        u2d(i,j)=u3d(i,1,j)
112        v2d(i,j)=v3d(i,1,j)
113        t2d(i,j)=t3d(i,1,j)
114        qv2d(i,j)=qv3d(i,1,j)
115        p2d(i,j)=p3d(i,1,j)
116        dz8w1d(i,j)=dz8w(i,1,j)
117        za(i,j)=0.5*dz8w1d(i,j)
118      enddo
119    enddo 
121    do j=jts, jte
122      do i=its, ite
124 !-----compute surface moisture flux
126        qfx(i,j)=qfx_interp/1000.
127        qfx(i,j)=amax1(qfx(i,j),0.)
128        lh(i,j)=xlv*qfx(i,j)
132 !-----compute surface heat flux
134        cpm(i,j)=cp*(1.+0.8*qv2d(i,j))
135 !       print*,'i j cpm xland qv2d',i,j,cpm(i,j),xland(i,j), qv2d(i,j)
136 !       print*,hfx_interp
137        if(xland(i,j)-1.5 .gt. 0.)then
138          hfx(i,j)=hfx_interp*cpm(i,j)
139        elseif(xland(i,j)-1.5 .lt. 0.)then
140          hfx(i,j)=hfx_interp*cpm(i,j)
141          hfx(i,j)=amax1(hfx(i,j),-250.)
142        endif
143      enddo
144    enddo
146    
147    if (itimestep .eq. 1) then
148      psih=0.0
149      psim=0.0
151    endif
152      chklowq=1.0 !JP
154    
155    do j=jts,jte
156      do i=its,ite
157        gz1oz0(i,j)=alog(za(i,j)/znt(i,j))
158        fh=gz1oz0(i,j)-psih(i,j)
159        fm=gz1oz0(i,j)-psim(i,j)
160        ch=karman**2/fh/fm
161        wspd(i,j)=sqrt(u2d(i,j)**2+v2d(i,j)**2)
162        dtdiff=-hfx(i,j)/den/cpm(i,j)/ch/wspd(i,j)
163        tsk(i,j)=t2d(i,j)-dtdiff
164      enddo
165    enddo
166    
167    end subroutine scmflux
168 !-------------------------------------------------------------------
169 end module module_sf_scmflux