Merge remote-tracking branch 'origin/release-v4.6.1'
[WRF.git] / phys / module_cu_ntiedtke.F
blob3b56132b667600b8a932389ad0262357d2a3d97f
1 !=================================================================================================================
2  module module_cu_ntiedtke
3  use ccpp_kind_types,only: kind_phys
5  use cu_ntiedtke,only: cu_ntiedtke_run, &
6                        cu_ntiedtke_init
7  use cu_ntiedtke_common
9  implicit none
10  private
11  public:: cu_ntiedtke_driver, &
12           ntiedtkeinit
15  contains
18 !=================================================================================================================
19  subroutine cu_ntiedtke_driver(                          &
20                  dt,itimestep,stepcu                     &
21                 ,raincv,pratec,qfx,hfx                   &
22                 ,u3d,v3d,w,t3d,qv3d,qc3d,qi3d,pi3d,rho3d &
23                 ,qvften,thften                           &
24                 ,dz8w,pcps,p8w,xland,cu_act_flag,dx      &
25                 ,f_qv,f_qc,f_qr,f_qi,f_qs                &
26                 ,grav,xlf,xls,xlv,rd,rv,cp               &
27                 ,rthcuten,rqvcuten,rqccuten,rqicuten     &
28                 ,rucuten,rvcuten                         &
29                 ,ids,ide,jds,jde,kds,kde                 &
30                 ,ims,ime,jms,jme,kms,kme                 &
31                 ,its,ite,jts,jte,kts,kte                 &
32                 ,errmsg,errflg)
33 !=================================================================================================================
34 !-- u3d         3d u-velocity interpolated to theta points (m/s)
35 !-- v3d         3d v-velocity interpolated to theta points (m/s)
36 !-- th3d        3d potential temperature (k)
37 !-- t3d         temperature (k)
38 !-- qv3d        3d water vapor mixing ratio (kg/kg)
39 !-- qc3d        3d cloud mixing ratio (kg/kg)
40 !-- qi3d        3d ice mixing ratio (kg/kg)
41 !-- rho3d       3d air density (kg/m^3)
42 !-- p8w         3d hydrostatic pressure at full levels (pa)
43 !-- pcps        3d hydrostatic pressure at half levels (pa)
44 !-- pi3d        3d exner function (dimensionless)
45 !-- qvften      3d total advective + PBL moisture tendency (kg kg-1 s-1)
46 !-- thften      3d total advective + PBL + radiative temperature tendency (k s-1)
47 !-- rthcuten      theta tendency due to 
48 !                 cumulus scheme precipitation (k/s)
49 !-- rucuten       u wind tendency due to 
50 !                 cumulus scheme precipitation (k/s)
51 !-- rvcuten       v wind tendency due to 
52 !                 cumulus scheme precipitation (k/s)
53 !-- rqvcuten      qv tendency due to 
54 !                 cumulus scheme precipitation (kg/kg/s)
55 !-- rqccuten      qc tendency due to 
56 !                 cumulus scheme precipitation (kg/kg/s)
57 !-- rqicuten      qi tendency due to 
58 !                 cumulus scheme precipitation (kg/kg/s)
59 !-- rainc         accumulated total cumulus scheme precipitation (mm)
60 !-- raincv        cumulus scheme precipitation (mm)
61 !-- pratec        precipitiation rate from cumulus scheme (mm/s)
62 !-- dz8w        dz between full levels (m)
63 !-- qfx         upward moisture flux at the surface (kg/m^2/s)
64 !-- hfx         upward heat flux at the surface (w/m^2) 
65 !-- dt          time step (s)
66 !-- ids         start index for i in domain
67 !-- ide         end index for i in domain
68 !-- jds         start index for j in domain
69 !-- jde         end index for j in domain
70 !-- kds         start index for k in domain
71 !-- kde         end index for k in domain
72 !-- ims         start index for i in memory
73 !-- ime         end index for i in memory
74 !-- jms         start index for j in memory
75 !-- jme         end index for j in memory
76 !-- kms         start index for k in memory
77 !-- kme         end index for k in memory
78 !-- its         start index for i in tile
79 !-- ite         end index for i in tile
80 !-- jts         start index for j in tile
81 !-- jte         end index for j in tile
82 !-- kts         start index for k in tile
83 !-- kte         end index for k in tile
84 !-----------------------------------------------------------------------------------------------------------------
86 !--- input arguments:
87  logical,intent(in),optional:: f_qv,f_qc,f_qr,f_qi,f_qs
89  integer,intent(in):: ids,ide,jds,jde,kds,kde, &
90                       ims,ime,jms,jme,kms,kme, &
91                       its,ite,jts,jte,kts,kte
93  integer,intent(in):: itimestep,stepcu
95  real(kind=kind_phys),intent(in):: cp,grav,rd,rv,xlf,xls,xlv
97  real(kind=kind_phys),intent(in):: dt
99  real(kind=kind_phys),intent(in),dimension(ims:ime,jms:jme):: dx,hfx,qfx,xland
101  real(kind=kind_phys),intent(in),dimension(ims:ime,kms:kme,jms:jme):: &
102     dz8w,   &
103     pcps,   &
104     p8w,    &
105     pi3d,   &
106     qc3d,   &
107     qvften, &
108     thften, &
109     qi3d,   &
110     qv3d,   &
111     rho3d,  &
112     t3d,    &
113     u3d,    &
114     v3d,    &
115     w
117 !--- inout arguments:
118  logical,intent(inout),dimension(ims:ime,jms:jme):: cu_act_flag
120  real(kind=kind_phys),intent(inout),dimension(ims:ime,jms:jme):: raincv, pratec
122  real(kind=kind_phys),intent(inout),dimension(ims:ime,kms:kme,jms:jme),optional:: &
123  rqccuten,  &
124  rqicuten,  &
125  rqvcuten,  &
126  rthcuten,  &
127  rucuten,   &
128  rvcuten
130 !--- output arguments:
131  character(len=*),intent(out):: errmsg
132  integer,intent(out):: errflg
134 !--- local variables and arrays:
135  integer:: i,im,j,k,kx,kx1
136  integer,dimension(its:ite)::  slimsk
138  real(kind=kind_phys):: delt
139  real(kind=kind_phys),dimension(its:ite):: rn
140  real(kind=kind_phys),dimension(its:ite,kts:kte):: prsl,omg,ghtl
141  real(kind=kind_phys),dimension(its:ite,kts:kte):: uf,vf,tf,qvf,qcf,qif
142  real(kind=kind_phys),dimension(its:ite,kts:kte):: qvftenz,thftenz
143  real(kind=kind_phys),dimension(its:ite,kts:kte+1):: prsi,ghti,zi
145  real(kind=kind_phys),dimension(its:ite):: dx_hv,hfx_hv,qfx_hv,xland_hv
146  real(kind=kind_phys),dimension(its:ite,kts:kte):: dz_hv,pi_hv,prsl_hv
147  real(kind=kind_phys),dimension(its:ite,kts:kte):: qv_hv,qc_hv,qi_hv,rho_hv,t_hv,u_hv,v_hv
148  real(kind=kind_phys),dimension(its:ite,kts:kte):: qvften_hv,thften_hv
149  real(kind=kind_phys),dimension(its:ite,kts:kte+1):: prsi_hv,w_hv
151  real(kind=kind_phys),dimension(its:ite):: raincv_hv,pratec_hv
152  real(kind=kind_phys),dimension(its:ite,kts:kte):: rthcuten_hv,rqvcuten_hv,rqccuten_hv,rqicuten_hv, &
153                                                    rucuten_hv,rvcuten_hv
155 !-----------------------------------------------------------------------------------------------------------------
157  errmsg = ' '
158  errflg = 0
160  call cu_ntiedtke_init( &
161                        con_cp  = cp  , con_rd = rd   , con_rv   = rv   , con_xlv = xlv    , &
162                        con_xls = xls , con_xlf = xlf , con_grav = grav , errmsg  = errmsg , &
163                        errflg  = errflg                                                     &
164                       )
166  do j = jts,jte
167     do i = its,ite
168        cu_act_flag(i,j)=.true.
169     enddo
170  enddo
172  do j = jts,jte
174     do i = its,ite
175        dx_hv(i)    = dx(i,j)
176        hfx_hv(i)   = hfx(i,j)
177        qfx_hv(i)   = qfx(i,j)
178        xland_hv(i) = xland(i,j)
179     enddo
181     do k = kts,kte
182        do i = its,ite
183           dz_hv(i,k)   = dz8w(i,k,j)
184           pi_hv(i,k)   = pi3d(i,k,j)
185           prsl_hv(i,k) = pcps(i,k,j)
186           qv_hv(i,k)   = qv3d(i,k,j)
187           qc_hv(i,k)   = qc3d(i,k,j)
188           qi_hv(i,k)   = qi3d(i,k,j)
189           rho_hv(i,k)  = rho3d(i,k,j)
190           t_hv(i,k)    = t3d(i,k,j)
191           u_hv(i,k)    = u3d(i,k,j)
192           v_hv(i,k)    = v3d(i,k,j)
194           qvften_hv(i,k) = qvften(i,k,j)
195           thften_hv(i,k) = thften(i,k,j)
196        enddo
197     enddo
198     do k = kts,kte+1
199        do i = its,ite
200           prsi_hv(i,k) = p8w(i,k,j)
201           w_hv(i,k)    = w(i,k,j)
202        enddo
203     enddo
205     call cu_ntiedtke_pre_run( &
206                      its    = its       , ite     = ite     , kts     = kts     , kte       = kte       , &
207                      im     = im        , kx      = kx      , kx1     = kx1     , itimestep = itimestep , &
208                      stepcu = stepcu    , dt      = dt      , grav    = grav    , xland     = xland_hv  , &
209                      dz     = dz_hv     , pres    = prsl_hv , presi   = prsi_hv , t         = t_hv      , &
210                      rho    = rho_hv    , qv      = qv_hv   , qc      = qc_hv   , qi        = qi_hv     , &
211                      u      = u_hv      , v       = v_hv    , w       = w_hv    , qvften    = qvften_hv , &
212                      thften = thften_hv , qvftenz = qvftenz , thftenz = thftenz , slimsk    = slimsk    , &
213                      delt   = delt      , prsl    = prsl    , ghtl    = ghtl    , tf        = tf        , &
214                      qvf    = qvf       , qcf     = qcf     , qif     = qif     , uf        = uf        , &
215                      vf     = vf        , prsi    = prsi    , ghti    = ghti    , omg       = omg       , &
216                      errmsg = errmsg    , errflg  = errflg                                                &
217                             )
219     call cu_ntiedtke_run( &
220                      pu   = uf     , pv   = vf     , pt     = tf      , pqv    = qvf     , &
221                      pqc  = qcf    , pqi  = qif    , pqvf   = qvftenz , ptf    = thftenz , &
222                      poz  = ghtl   , pzz  = ghti   , pomg   = omg     , pap    = prsl    , &
223                      paph = prsi   , evap = qfx_hv , hfx    = hfx_hv  , zprecc = rn      , &
224                      lndj = slimsk , lq   = im     , km     = kx      , km1    = kx1     , &
225                      dt   = delt   , dx   = dx_hv  , errmsg = errmsg  , errflg = errflg    &
226                         )
228     call cu_ntiedtke_post_run( &
229                      its      = its         , ite      = ite         , kts      = kts         , kte      = kte         , &
230                      stepcu   = stepcu      , dt       = dt          , exner    = pi_hv       , qv       = qv_hv       , &
231                      qc       = qc_hv       , qi       = qi_hv       , t        = t_hv        , u        = u_hv        , &
232                      v        = v_hv        , qvf      = qvf         , qcf      = qcf         , qif      = qif         , &
233                      tf       = tf          , uf       = uf          , vf       = vf          , rn       = rn          , &
234                      raincv   = raincv_hv   , pratec   = pratec_hv   , rthcuten = rthcuten_hv , rqvcuten = rqvcuten_hv , &
235                      rqccuten = rqccuten_hv , rqicuten = rqicuten_hv , rucuten  = rucuten_hv  , rvcuten  = rvcuten_hv  , &
236                      errmsg   = errmsg      , errflg   = errflg                                                          &
237                              )
239     do i = its,ite
240        raincv(i,j) = raincv_hv(i)
241        pratec(i,j) = pratec_hv(i)
242     enddo
244     do k = kts,kte
245        do i = its,ite
246           rucuten(i,k,j)  = rucuten_hv(i,k)
247           rvcuten(i,k,j)  = rvcuten_hv(i,k)
248           rthcuten(i,k,j) = rthcuten_hv(i,k)
249           rqvcuten(i,k,j) = rqvcuten_hv(i,k)
250        enddo
251     enddo
253     if(present(rqccuten))then
254        if(f_qc) then
255           do k = kts,kte
256              do i = its,ite
257                 rqccuten(i,k,j) = rqccuten_hv(i,k)
258              enddo
259           enddo
260        endif
261     endif
263     if(present(rqicuten))then
264        if(f_qi) then
265           do k = kts,kte
266              do i = its,ite
267                 rqicuten(i,k,j) = rqicuten_hv(i,k)
268              enddo
269           enddo
270        endif
271     endif
273  enddo
275  end subroutine cu_ntiedtke_driver
277 !=================================================================================================================
278  subroutine ntiedtkeinit(rthcuten,rqvcuten,rqccuten,rqicuten, &
279                          rucuten,rvcuten,rthften,rqvften,     &
280                          restart,p_qc,p_qi,p_first_scalar,    &
281                          allowed_to_read,                     &
282                          ids, ide, jds, jde, kds, kde,        &
283                          ims, ime, jms, jme, kms, kme,        &
284                          its, ite, jts, jte, kts, kte)
285 !=================================================================================================================
287 !--- input arguments:
288  logical,intent(in):: allowed_to_read,restart
290  integer,intent(in):: ids, ide, jds, jde, kds, kde, &
291                       ims, ime, jms, jme, kms, kme, &
292                       its, ite, jts, jte, kts, kte
293  integer,intent(in):: p_first_scalar,p_qi,p_qc
295 !--- output arguments:
296  real(kind=kind_phys),intent(out),dimension(ims:ime,kms:kme,jms:jme )::  &
297     rthcuten,rqvcuten,rqccuten,rqicuten,rucuten,rvcuten,rthften,rqvften
299 !--- local variables and arrays:
300  integer:: i,j,k,itf,jtf,ktf
302 !-----------------------------------------------------------------------------------------------------------------
304  jtf = min0(jte,jde-1)
305  ktf = min0(kte,kde-1)
306  itf = min0(ite,ide-1)
308  if(.not.restart)then
309     do j = jts,jtf
310        do k = kts,ktf
311           do i = its,itf
312              rthcuten(i,k,j) = 0.
313              rqvcuten(i,k,j) = 0.
314              rucuten(i,k,j)  = 0.
315              rvcuten(i,k,j)  = 0.
316           enddo
317        enddo
318     enddo
320     do j = jts,jtf
321        do k = kts,ktf
322           do i = its,itf
323              rthften(i,k,j)=0.
324              rqvften(i,k,j)=0.
325           enddo
326        enddo
327     enddo
329     if(p_qc .ge. p_first_scalar) then
330        do j = jts,jtf
331           do k = kts,ktf
332              do i = its,itf
333                 rqccuten(i,k,j)=0.
334              enddo
335           enddo
336        enddo
337     endif
339     if(p_qi .ge. p_first_scalar) then
340        do j = jts,jtf
341           do k = kts,ktf
342              do i = its,itf
343                 rqicuten(i,k,j)=0.
344              enddo
345           enddo
346        enddo
347     endif
348  endif
350  end subroutine ntiedtkeinit
352 !=================================================================================================================
353  subroutine cu_ntiedtke_pre_run(its,ite,kts,kte,im,kx,kx1,itimestep,stepcu,dt,grav,xland,dz,pres,presi,   &
354                                 t,rho,qv,qc,qi,u,v,w,qvften,thften,qvftenz,thftenz,slimsk,delt,prsl,ghtl, &
355                                 tf,qvf,qcf,qif,uf,vf,prsi,ghti,omg,errmsg,errflg)
356 !=================================================================================================================
358 !--- input arguments:
359  integer,intent(in):: its,ite,kts,kte
360  integer,intent(in):: itimestep
361  integer,intent(in):: stepcu
363  real(kind=kind_phys),intent(in):: dt,grav
364  real(kind=kind_phys),intent(in),dimension(its:ite):: xland
365  real(kind=kind_phys),intent(in),dimension(its:ite,kts:kte):: dz,pres,t,rho,qv,qc,qi,u,v
366  real(kind=kind_phys),intent(inout),dimension(its:ite,kts:kte):: qvften,thften
367  real(kind=kind_phys),intent(in),dimension(its:ite,kts:kte+1):: presi,w
369 !--- inout arguments:
370  integer,intent(inout):: im,kx,kx1
371  integer,intent(inout),dimension(its:ite):: slimsk
373  real(kind=kind_phys),intent(inout):: delt
374  real(kind=kind_phys),intent(inout),dimension(its:ite,kts:kte):: tf,qvf,qcf,qif,uf,vf
375  real(kind=kind_phys),intent(inout),dimension(its:ite,kts:kte):: ghtl,omg,prsl
376  real(kind=kind_phys),intent(inout),dimension(its:ite,kts:kte):: qvftenz,thftenz
377  real(kind=kind_phys),intent(inout),dimension(its:ite,kts:kte+1):: ghti,prsi
379 !--- output arguments:
380  character(len=*),intent(out):: errmsg
381  integer,intent(out):: errflg
383 !--- local variables and arrays:
384  integer:: i,k,pp,zz
386  real(kind=kind_phys),dimension(its:ite,kts:kte):: zl,dot
387  real(kind=kind_phys),dimension(its:ite,kts:kte+1):: zi
389 !-----------------------------------------------------------------------------------------------------------------
391  im  = ite-its+1
392  kx  = kte-kts+1
393  kx1 = kx+1
395  delt  = dt*stepcu
397  do i = its,ite
398     slimsk(i) = (abs(xland(i)-2.))
399  enddo
401  k = kts
402  do i = its,ite
403     zi(i,k) = 0.
404  enddo
405  do k = kts,kte
406     do i = its,ite
407        zi(i,k+1) = zi(i,k)+dz(i,k)
408     enddo
409  enddo
410  do k = kts,kte
411     do i = its,ite
412        zl(i,k)  = 0.5*(zi(i,k)+zi(i,k+1))
413        dot(i,k) = -0.5*grav*rho(i,k)*(w(i,k)+w(i,k+1))
414     enddo
415  enddo
417  pp = 0
418  do k = kts,kte+1
419     zz = kte + 1 - pp
420     do i = its,ite
421        ghti(i,zz) = zi(i,k)
422        prsi(i,zz) = presi(i,k)
423     enddo
424     pp = pp + 1
425  enddo
426  pp = 0
427  do k = kts,kte
428     zz = kte-pp
429     do i = its,ite
430        ghtl(i,zz) = zl(i,k)
431        omg(i,zz)  = dot(i,k)
432        prsl(i,zz) = pres(i,k)
433     enddo
434     pp = pp + 1
435  enddo
437  pp = 0
438  do k = kts,kte
439     zz = kte-pp
440     do i = its,ite
441        tf(i,zz)  = t(i,k)
442        qvf(i,zz) = qv(i,k)
443        qcf(i,zz) = qc(i,k)
444        qif(i,zz) = qi(i,k)
445        uf(i,zz)  = u(i,k)
446        vf(i,zz)  = v(i,k)
447     enddo
448     pp = pp + 1
449  enddo
451  if(itimestep == 1) then
452     do k = kts,kte
453        do i = its,ite
454           qvftenz(i,k) = 0.
455           thftenz(i,k) = 0.
456        enddo
457     enddo
458  else
459     pp = 0
460     do k = kts,kte
461        zz = kte-pp
462        do i = its,ite
463           qvftenz(i,zz) = qvften(i,k)
464           thftenz(i,zz) = thften(i,k)
465        enddo
466        pp = pp + 1
467     enddo
468  endif
470  errmsg = 'cu_ntiedtke_pre_run OK'
471  errflg = 0
473  end subroutine cu_ntiedtke_pre_run
475 !=================================================================================================================
476  subroutine cu_ntiedtke_post_run(its,ite,kts,kte,stepcu,dt,exner,qv,qc,qi,t,u,v,qvf,qcf,qif,tf,uf,vf,rn,raincv, &
477                                  pratec,rthcuten,rqvcuten,rqccuten,rqicuten,rucuten,rvcuten,errmsg,errflg)
478 !=================================================================================================================
480 !--- input arguments:
481  integer,intent(in):: its,ite,kts,kte
482  integer,intent(in):: stepcu
484  real(kind=kind_phys),intent(in):: dt
485  real(kind=kind_phys),intent(in),dimension(its:ite):: rn
486  real(kind=kind_phys),intent(in),dimension(its:ite,kts:kte):: exner,qv,qc,qi,t,u,v,qvf,qcf,qif,tf,uf,vf
488 !--- inout arguments:
489  real(kind=kind_phys),intent(inout),dimension(its:ite):: raincv,pratec
490  real(kind=kind_phys),intent(inout),dimension(its:ite,kts:kte):: rqvcuten,rqccuten,rqicuten
491  real(kind=kind_phys),intent(inout),dimension(its:ite,kts:kte):: rthcuten,rucuten,rvcuten
493 !--- output arguments:
494  character(len=*),intent(out):: errmsg
495  integer,intent(out):: errflg
497 !--- local variables and arrays:
498  integer:: i,k,pp,zz
500  real(kind=kind_phys):: delt,rdelt
502 !-----------------------------------------------------------------------------------------------------------------
504  delt  = dt*stepcu
505  rdelt = 1./delt
507  do i = its,ite
508     raincv(i) = rn(i)/stepcu
509     pratec(i) = rn(i)/(stepcu*dt)
510  enddo
512  pp = 0
513  do k = kts,kte
514     zz = kte - pp
515     do i = its,ite
516        rthcuten(i,k) = (tf(i,zz)-t(i,k))/exner(i,k)*rdelt
517        rqvcuten(i,k) = (qvf(i,zz)-qv(i,k))*rdelt
518        rqccuten(i,k) = (qcf(i,zz)-qc(i,k))*rdelt
519        rqicuten(i,k) = (qif(i,zz)-qi(i,k))*rdelt
520        rucuten(i,k)  = (uf(i,zz)-u(i,k))*rdelt
521        rvcuten(i,k)  = (vf(i,zz)-v(i,k))*rdelt
522     enddo
523     pp = pp + 1
524  enddo
526  errmsg = 'cu_ntiedtke_timestep_final OK'
527  errflg = 0
529  end subroutine cu_ntiedtke_post_run
531 !=================================================================================================================
532  end module module_cu_ntiedtke
533 !=================================================================================================================