2 ! WRF:MODEL_LAYER:PHYSICS
4 MODULE module_diffusion_em
6 USE module_bc, only: set_physical_bc3d, set_physical_bc2d ! XZ
7 USE module_state_description, only: p_m23, p_m13, p_m22, p_m33, p_r23, p_r13, p_r12, p_m12, p_m11, &
8 P_QNS, P_QNR, P_QNG, P_QT, P_QNH, P_QVOLG ! XZ
9 USE module_big_step_utilities_em, only: grid_config_rec_type, param_first_scalar, p_qv, p_qi, p_qc
10 USE module_model_constants
14 !=======================================================================
15 !=======================================================================
17 SUBROUTINE cal_deform_and_div( config_flags, u, v, w, div, &
18 defor11, defor22, defor33, &
19 defor12, defor13, defor23, &
21 u_base, v_base, msfux, msfuy, &
22 msfvx, msfvy, msftx, msfty, &
23 rdx, rdy, dn, dnw, rdz, rdzw, &
24 fnm, fnp, cf1, cf2, cf3, zx, zy, &
25 ids, ide, jds, jde, kds, kde, &
26 ims, ime, jms, jme, kms, kme, &
27 its, ite, jts, jte, kts, kte )
29 ! History: Sep 2003 Changes by Jason Knievel and George Bryan, NCAR
30 ! Oct 2001 Converted to mass core by Bill Skamarock, NCAR
33 ! Purpose: This routine calculates deformation and 3-d divergence.
35 ! References: Klemp and Wilhelmson (JAS 1978)
36 ! Chen and Dudhia (NCAR WRF physics report 2000)
38 !-----------------------------------------------------------------------
40 ! Equations 13a-f, Chen and Dudhia 2000, Appendix A:
41 ! Eqn 13a: D11=defor11= 2m^2 * (partial du^/dX + partial dpsi/dx * partial du^/dpsi)
42 ! Eqn 13b: D22=defor22= 2m^2 * (partial dv^/dY + partial dpsi/dy * partial dv^/dpsi)
43 ! Eqn 13c: D33=defor33= 2 * partial dw/dz [SIMPLER FORM]
44 ! Eqn 13d: D12=defor12= m^2 * (partial dv^/dX + partial du^/dY +
45 ! partial dpsi/dx * partial dv^/dpsi +
46 ! partial dpsi/dy * partial du^/dpsi)
47 ! Eqn 13e: D13=defor13= m^2 * (partial dw^/dX + partial dpsi/dx * partial dw^/dpsi)
48 ! + partial du/dz [SIMPLER FORM]
49 ! Eqn 13f: D23=defor23= m^2 * (partial dw^/dY + partial dpsi/dy * partial dw^/dpsi)
50 ! + partial dv/dz [SIMPLER FORM]
51 !-----------------------------------------------------------------------
56 TYPE( grid_config_rec_type ), INTENT( IN ) &
59 INTEGER, INTENT( IN ) &
60 :: ids, ide, jds, jde, kds, kde, &
61 ims, ime, jms, jme, kms, kme, &
62 its, ite, jts, jte, kts, kte
65 :: rdx, rdy, cf1, cf2, cf3
67 REAL, DIMENSION( kms:kme ), INTENT( IN ) &
68 :: fnm, fnp, dn, dnw, u_base, v_base
70 REAL, DIMENSION( ims:ime , jms:jme ), INTENT( IN ) &
71 :: msfux, msfuy, msfvx, msfvy, msftx, msfty
73 REAL, DIMENSION( ims:ime, kms:kme, jms:jme ), INTENT( IN ) &
74 :: u, v, w, zx, zy, rdz, rdzw
76 REAL, DIMENSION( ims:ime, kms:kme, jms:jme ), INTENT( INOUT ) &
77 :: defor11, defor22, defor33, defor12, defor13, defor23, div
79 INTEGER, INTENT( IN ) :: n_nba_rij
81 REAL , DIMENSION(ims:ime, kms:kme, jms:jme, n_nba_rij), INTENT(INOUT) &
88 :: i, j, k, ktf, ktes1, ktes2, i_start, i_end, j_start, j_end
91 :: tmp, tmpzx, tmpzy, tmpzeta_z, cft1, cft2
93 REAL, DIMENSION( its:ite, jts:jte ) &
94 :: mm, zzavg, zeta_zd12
96 REAL, DIMENSION( its-2:ite+2, kts:kte, jts-2:jte+2 ) &
100 !-----------------------------------------------------------------------
102 ! Comments 10-MAR-2005
103 ! Treat all differentials as 'div-style' [or 'curl-style'],
104 ! i.e., du/dX becomes (in map coordinate space) mx*my * d(u/my)/dx,
105 ! NB - all equations referred to here are from Chen and Dudhia 2002, from the
106 ! WRF physics documents web pages:
107 ! http://www.mmm.ucar.edu/wrf/users/docs/wrf-doc-physics.pdf
109 !=======================================================================
110 ! In the following section, calculate 3-d divergence and the first three
111 ! (defor11, defor22, defor33) of six deformation terms.
116 cft2 = - 0.5 * dnw(ktes1) / dn(ktes1)
119 ktf = MIN( kte, kde-1 )
122 i_end = MIN( ite, ide-1 )
124 j_end = MIN( jte, jde-1 )
126 ! Square the map scale factor.
128 DO j = j_start, j_end
129 DO i = i_start, i_end
130 mm(i,j) = msftx(i,j) * msfty(i,j)
134 !-----------------------------------------------------------------------
137 ! Apply a coordinate transformation to zonal velocity, u.
139 DO j = j_start, j_end
141 DO i = i_start, i_end+1
142 hat(i,k,j) = u(i,k,j) / msfuy(i,j)
147 ! Average in x and z.
152 hatavg(i,k,j) = 0.5 * &
153 ( fnm(k) * ( hat(i,k ,j) + hat(i+1, k,j) ) + &
154 fnp(k) * ( hat(i,k-1,j) + hat(i+1,k-1,j) ) )
159 ! Extrapolate to top and bottom of domain (to w levels).
161 DO j = j_start, j_end
162 DO i = i_start, i_end
163 hatavg(i,1,j) = 0.5 * ( &
164 cf1 * hat(i ,1,j) + &
165 cf2 * hat(i ,2,j) + &
166 cf3 * hat(i ,3,j) + &
167 cf1 * hat(i+1,1,j) + &
168 cf2 * hat(i+1,2,j) + &
170 hatavg(i,kte,j) = 0.5 * ( &
171 cft1 * ( hat(i,ktes1,j) + hat(i+1,ktes1,j) ) + &
172 cft2 * ( hat(i,ktes2,j) + hat(i+1,ktes2,j) ) )
177 ! Eqn 13a: D11=defor11= 2m^2 * (partial du^/dX + partial dpsi/dx * partial du^/dpsi)
178 ! Below, D11 is set = 2*tmp1
179 ! => tmp1 = m^2 * (partial du^/dX + partial dpsi/dx * partial du^/dpsi)
180 ! tmpzx = averaged value of dpsi/dx (=zx)
182 DO j = j_start, j_end
184 DO i = i_start, i_end
186 zx(i,k ,j) + zx(i+1,k ,j) + &
187 zx(i,k+1,j) + zx(i+1,k+1,j) )
188 tmp1(i,k,j) = ( hatavg(i,k+1,j) - hatavg(i,k,j) ) *tmpzx * rdzw(i,k,j)
189 ! tmp1 to here = partial dpsi/dx * partial du^/dpsi:
194 DO j = j_start, j_end
196 DO i = i_start, i_end
197 tmp1(i,k,j) = mm(i,j) * ( rdx * ( hat(i+1,k,j) - hat(i,k,j) ) - &
203 ! End calculation of du/dx.
204 !-----------------------------------------------------------------------
206 !-----------------------------------------------------------------------
207 ! Calculate defor11 (2*du/dx).
209 ! Eqn 13a: D11=defor11= 2 m^2 * (partial du^/dX + partial dpsi/dx * partial du^/dpsi)
212 DO j = j_start, j_end
214 DO i = i_start, i_end
215 defor11(i,k,j) = 2.0 * tmp1(i,k,j)
220 ! End calculation of defor11.
221 !-----------------------------------------------------------------------
223 !-----------------------------------------------------------------------
224 ! Calculate zonal divergence (du/dx) and add it to the divergence array.
226 DO j = j_start, j_end
228 DO i = i_start, i_end
229 div(i,k,j) = tmp1(i,k,j)
234 ! End calculation of zonal divergence.
235 !-----------------------------------------------------------------------
237 !-----------------------------------------------------------------------
240 ! Apply a coordinate transformation to meridional velocity, v.
242 DO j = j_start, j_end+1
244 DO i = i_start, i_end
245 ! Because msfvx at the poles will be undefined (1./0.), we will have
246 ! trouble. But we are OK since v at the poles is 0., and that takes
247 ! precedence in this case.
248 IF ((config_flags%polar) .AND. ((j == jds) .OR. (j == jde))) THEN
251 hat(i,k,j) = v(i,k,j) / msfvx(i,j)
257 ! Account for the slope in y of eta surfaces.
262 hatavg(i,k,j) = 0.5 * ( &
263 fnm(k) * ( hat(i,k ,j) + hat(i,k ,j+1) ) + &
264 fnp(k) * ( hat(i,k-1,j) + hat(i,k-1,j+1) ) )
269 ! Extrapolate to top and bottom of domain (to w levels).
271 DO j = j_start, j_end
272 DO i = i_start, i_end
273 hatavg(i,1,j) = 0.5 * ( &
274 cf1 * hat(i,1,j ) + &
275 cf2 * hat(i,2,j ) + &
276 cf3 * hat(i,3,j ) + &
277 cf1 * hat(i,1,j+1) + &
278 cf2 * hat(i,2,j+1) + &
280 hatavg(i,kte,j) = 0.5 * ( &
281 cft1 * ( hat(i,ktes1,j) + hat(i,ktes1,j+1) ) + &
282 cft2 * ( hat(i,ktes2,j) + hat(i,ktes2,j+1) ) )
287 ! Eqn 13b: D22=defor22= 2m^2 * (partial dv^/dY + partial dpsi/dy * partial dv^/dpsi)
288 ! Below, D22 is set = 2*tmp1
289 ! => tmp1 = m^2 * (partial dv^/dY + partial dpsi/dy * partial dv^/dpsi)
290 ! tmpzy = averaged value of dpsi/dy (=zy)
292 DO j = j_start, j_end
294 DO i = i_start, i_end
296 zy(i,k ,j) + zy(i,k ,j+1) + &
297 zy(i,k+1,j) + zy(i,k+1,j+1) )
298 tmp1(i,k,j) = ( hatavg(i,k+1,j) - hatavg(i,k,j) ) * tmpzy * rdzw(i,k,j)
299 ! tmp1 to here = partial dpsi/dy * partial dv^/dpsi:
304 DO j = j_start, j_end
306 DO i = i_start, i_end
307 tmp1(i,k,j) = mm(i,j) * ( &
308 rdy * ( hat(i,k,j+1) - hat(i,k,j) ) - tmp1(i,k,j) )
313 ! End calculation of dv/dy.
314 !-----------------------------------------------------------------------
316 !-----------------------------------------------------------------------
317 ! Calculate defor22 (2*dv/dy).
319 ! Eqn 13b: D22=defor22= 2 m^2 * (partial dv^/dY + partial dpsi/dy * partial dv^/dpsi)
322 DO j = j_start, j_end
324 DO i = i_start, i_end
325 defor22(i,k,j) = 2.0 * tmp1(i,k,j)
330 ! End calculation of defor22.
331 !-----------------------------------------------------------------------
333 !-----------------------------------------------------------------------
334 ! Calculate meridional divergence (dv/dy) and add it to the divergence
337 DO j = j_start, j_end
339 DO i = i_start, i_end
340 div(i,k,j) = div(i,k,j) + tmp1(i,k,j)
345 ! End calculation of meridional divergence.
346 !-----------------------------------------------------------------------
348 !-----------------------------------------------------------------------
350 ! Eqn 13c: D33=defor33= 2 * partial dw/dz
351 ! Below, D33 is set = 2*tmp1
352 ! => tmp1 = partial dw/dz
356 DO j = j_start, j_end
358 DO i = i_start, i_end
359 tmp1(i,k,j) = ( w(i,k+1,j) - w(i,k,j) ) * rdzw(i,k,j)
364 ! End calculation of dw/dz.
365 !-----------------------------------------------------------------------
367 !-----------------------------------------------------------------------
368 ! Calculate defor33 (2*dw/dz).
370 DO j = j_start, j_end
372 DO i = i_start, i_end
373 defor33(i,k,j) = 2.0 * tmp1(i,k,j)
378 ! End calculation of defor33.
379 !-----------------------------------------------------------------------
381 !-----------------------------------------------------------------------
382 ! Calculate vertical divergence (dw/dz) and add it to the divergence
385 DO j = j_start, j_end
387 DO i = i_start, i_end
388 div(i,k,j) = div(i,k,j) + tmp1(i,k,j)
393 ! End calculation of vertical divergence.
394 !-----------------------------------------------------------------------
396 ! Three-dimensional divergence is now finished and values are in array
397 ! "div." Also, the first three (defor11, defor22, defor33) of six
398 ! deformation terms are now calculated at pressure points.
399 !=======================================================================
401 ! Comments 10-MAR-2005
402 ! Treat all differentials as 'div-style' [or 'curl-style'],
403 ! i.e., du/dY becomes (in map coordinate space) mx*my * d(u/mx)/dy,
404 ! dv/dX becomes (in map coordinate space) mx*my * d(v/my)/dx,
405 ! (see e.g. Haltiner and Williams p. 441)
407 !=======================================================================
408 ! Calculate the final three deformations (defor12, defor13, defor23) at
416 IF ( config_flags%open_xs .OR. config_flags%specified .OR. &
417 config_flags%nested) i_start = MAX( ids+1, its )
418 IF ( config_flags%open_xe .OR. config_flags%specified .OR. &
419 config_flags%nested) i_end = MIN( ide-1, ite )
420 IF ( config_flags%open_ys .OR. config_flags%specified .OR. &
421 config_flags%nested) j_start = MAX( jds+1, jts )
422 IF ( config_flags%open_ye .OR. config_flags%specified .OR. &
423 config_flags%nested) j_end = MIN( jde-1, jte )
424 IF ( config_flags%periodic_x ) i_start = its
425 IF ( config_flags%periodic_x ) i_end = ite
428 !-----------------------------------------------------------------------
431 ! First, calculate an average mapscale factor.
434 ! du/dy => need u map scale factor in x (which is defined at u points)
435 ! averaged over j and j-1
436 ! dv/dx => need v map scale factor in y (which is defined at v points)
437 ! averaged over i and i-1
439 DO j = j_start, j_end
440 DO i = i_start, i_end
441 mm(i,j) = 0.25 * ( msfux(i,j-1) + msfux(i,j) ) * ( msfvy(i-1,j) + msfvy(i,j) )
445 ! Apply a coordinate transformation to zonal velocity, u.
447 DO j =j_start-1, j_end
450 ! Fixes to set_physical_bc2/3d for polar boundary conditions
451 ! remove issues with loop over j
452 hat(i,k,j) = u(i,k,j) / msfux(i,j)
457 ! Average in y and z.
462 hatavg(i,k,j) = 0.5 * ( &
463 fnm(k) * ( hat(i,k ,j-1) + hat(i,k ,j) ) + &
464 fnp(k) * ( hat(i,k-1,j-1) + hat(i,k-1,j) ) )
469 ! Extrapolate to top and bottom of domain (to w levels).
471 DO j = j_start, j_end
472 DO i = i_start, i_end
473 hatavg(i,1,j) = 0.5 * ( &
474 cf1 * hat(i,1,j-1) + &
475 cf2 * hat(i,2,j-1) + &
476 cf3 * hat(i,3,j-1) + &
477 cf1 * hat(i,1,j ) + &
478 cf2 * hat(i,2,j ) + &
480 hatavg(i,kte,j) = 0.5 * ( &
481 cft1 * ( hat(i,ktes1,j-1) + hat(i,ktes1,j) ) + &
482 cft2 * ( hat(i,ktes2,j-1) + hat(i,ktes2,j) ) )
486 ! tmpzy = averaged value of dpsi/dy (=zy) on vorticity grid
487 ! tmp1 = partial dpsi/dy * partial du^/dpsi
488 DO j = j_start, j_end
490 DO i = i_start, i_end
492 zy(i-1,k ,j) + zy(i,k ,j) + &
493 zy(i-1,k+1,j) + zy(i,k+1,j) )
494 tmp1(i,k,j) = ( hatavg(i,k+1,j) - hatavg(i,k,j) ) * &
495 0.25 * tmpzy * ( rdzw(i,k,j) + rdzw(i-1,k,j) + &
496 rdzw(i-1,k,j-1) + rdzw(i,k,j-1) )
501 ! End calculation of du/dy.
502 !----------------------------------------------------------------------
504 !-----------------------------------------------------------------------
505 ! Add the first term to defor12 (du/dy+dv/dx) at vorticity points.
508 ! Eqn 13d: D12=defor12= m^2 * (partial dv^/dX + partial du^/dY +
509 ! partial dpsi/dx * partial dv^/dpsi +
510 ! partial dpsi/dy * partial du^/dpsi)
511 ! Here deal with m^2 * (partial du^/dY + partial dpsi/dy * partial du^/dpsi)
512 ! Still need to add v^ terms:
513 ! m^2 * (partial dv^/dX + partial dpsi/dx * partial dv^/dpsi)
515 DO j = j_start, j_end
517 DO i = i_start, i_end
518 defor12(i,k,j) = mm(i,j) * ( &
519 rdy * ( hat(i,k,j) - hat(i,k,j-1) ) - tmp1(i,k,j) )
524 ! End addition of the first term to defor12.
525 !-----------------------------------------------------------------------
527 !-----------------------------------------------------------------------
530 ! Apply a coordinate transformation to meridional velocity, v.
532 DO j = j_start, j_end
534 DO i = i_start-1, i_end
535 hat(i,k,j) = v(i,k,j) / msfvy(i,j)
540 ! Account for the slope in x of eta surfaces.
542 DO j = j_start, j_end
544 DO i = i_start, i_end
545 hatavg(i,k,j) = 0.5 * ( &
546 fnm(k) * ( hat(i-1,k ,j) + hat(i,k ,j) ) + &
547 fnp(k) * ( hat(i-1,k-1,j) + hat(i,k-1,j) ) )
552 ! Extrapolate to top and bottom of domain (to w levels).
554 DO j = j_start, j_end
555 DO i = i_start, i_end
556 hatavg(i,1,j) = 0.5 * ( &
557 cf1 * hat(i-1,1,j) + &
558 cf2 * hat(i-1,2,j) + &
559 cf3 * hat(i-1,3,j) + &
560 cf1 * hat(i ,1,j) + &
561 cf2 * hat(i ,2,j) + &
563 hatavg(i,kte,j) = 0.5 * ( &
564 cft1 * ( hat(i,ktes1,j) + hat(i-1,ktes1,j) ) + &
565 cft2 * ( hat(i,ktes2,j) + hat(i-1,ktes2,j) ) )
569 ! Fixes to set_physical_bc2/3d have made any check for polar B.C.'s
570 ! unnecessary in this place. zx, rdzw, and hatavg are all defined
571 ! in places they need to be and the values at the poles are replications
572 ! of the values one grid point in, so the averaging over j and j-1 works
573 ! to act as just using the value at j or j-1 (with out extra code).
575 ! tmpzx = averaged value of dpsi/dx (=zx) on vorticity grid
576 ! tmp1 = partial dpsi/dx * partial dv^/dpsi
577 DO j = j_start, j_end
579 DO i = i_start, i_end
581 zx(i,k ,j-1) + zx(i,k ,j) + &
582 zx(i,k+1,j-1) + zx(i,k+1,j) )
583 tmp1(i,k,j) = ( hatavg(i,k+1,j) - hatavg(i,k,j) ) * &
584 0.25 * tmpzx * ( rdzw(i,k,j) + rdzw(i,k,j-1) + &
585 rdzw(i-1,k,j-1) + rdzw(i-1,k,j) )
590 ! End calculation of dv/dx.
591 !-----------------------------------------------------------------------
593 !-----------------------------------------------------------------------
594 ! Add the second term to defor12 (du/dy+dv/dx) at vorticity points.
597 ! Eqn 13d: D12=defor12= m^2 * (partial dv^/dX + partial du^/dY +
598 ! partial dpsi/dx * partial dv^/dpsi +
599 ! partial dpsi/dy * partial du^/dpsi)
600 ! Here adding v^ terms:
601 ! m^2 * (partial dv^/dX + partial dpsi/dx * partial dv^/dpsi)
603 IF ( config_flags%sfs_opt .GT. 0 ) THEN ! NBA--
605 !JDM____________________________________________________________________
607 ! s12 = du/dy + dv/dx
608 ! = (du/dy - dz/dy*du/dz) + (dv/dx - dz/dx*dv/dz)
609 ! ______defor12______ ___tmp1___
611 ! r12 = du/dy - dv/dx
612 ! = (du/dy - dz/dy*du/dz) - (dv/dx - dz/dx*dv/dz)
613 ! ______defor12______ ___tmp1___
614 !_______________________________________________________________________
617 DO j = j_start, j_end
619 DO i = i_start, i_end
621 nba_rij(i,k,j,P_r12) = defor12(i,k,j) - &
623 rdx * ( hat(i,k,j) - hat(i-1,k,j) ) - tmp1(i,k,j) )
625 defor12(i,k,j) = defor12(i,k,j) + &
627 rdx * ( hat(i,k,j) - hat(i-1,k,j) ) - tmp1(i,k,j) )
632 ! End addition of the second term to defor12.
633 !-----------------------------------------------------------------------
635 !-----------------------------------------------------------------------
636 ! Update the boundary for defor12 (might need to change later).
638 IF ( .NOT. config_flags%periodic_x .AND. i_start .EQ. ids+1 ) THEN
641 defor12(ids,k,j) = defor12(ids+1,k,j)
642 nba_rij(ids,k,j,P_r12) = nba_rij(ids+1,k,j,P_r12)
647 IF ( .NOT. config_flags%periodic_y .AND. j_start .EQ. jds+1) THEN
650 defor12(i,k,jds) = defor12(i,k,jds+1)
651 nba_rij(i,k,jds,P_r12) = nba_rij(i,k,jds+1,P_r12)
656 IF ( .NOT. config_flags%periodic_x .AND. i_end .EQ. ide-1) THEN
659 defor12(ide,k,j) = defor12(ide-1,k,j)
660 nba_rij(ide,k,j,P_r12) = nba_rij(ide-1,k,j,P_r12)
665 IF ( .NOT. config_flags%periodic_y .AND. j_end .EQ. jde-1) THEN
668 defor12(i,k,jde) = defor12(i,k,jde-1)
669 nba_rij(i,k,jde,P_r12) = nba_rij(i,k,jde-1,P_r12)
674 ELSE ! NOT NBA--------------------------------------------------------
676 DO j = j_start, j_end
678 DO i = i_start, i_end
679 defor12(i,k,j) = defor12(i,k,j) + &
681 rdx * ( hat(i,k,j) - hat(i-1,k,j) ) - tmp1(i,k,j) )
686 ! End addition of the second term to defor12.
687 !-----------------------------------------------------------------------
689 !-----------------------------------------------------------------------
690 ! Update the boundary for defor12 (might need to change later).
692 IF ( .NOT. config_flags%periodic_x .AND. i_start .EQ. ids+1 ) THEN
695 defor12(ids,k,j) = defor12(ids+1,k,j)
700 IF ( .NOT. config_flags%periodic_y .AND. j_start .EQ. jds+1) THEN
703 defor12(i,k,jds) = defor12(i,k,jds+1)
708 IF ( .NOT. config_flags%periodic_x .AND. i_end .EQ. ide-1) THEN
711 defor12(ide,k,j) = defor12(ide-1,k,j)
716 IF ( .NOT. config_flags%periodic_y .AND. j_end .EQ. jde-1) THEN
719 defor12(i,k,jde) = defor12(i,k,jde-1)
724 ENDIF ! NBA-----------------------------------------------------------
726 ! End update of boundary for defor12.
727 !-----------------------------------------------------------------------
730 ! Further deformation terms not needed for 2-dimensional Smagorinsky diffusion,
731 ! so those terms have not been dealt with yet.
732 ! A "y" has simply been added to all map scale factors to allow the model to
733 ! compile without errors.
735 !-----------------------------------------------------------------------
739 i_end = MIN( ite, ide-1 )
741 j_end = MIN( jte, jde-1 )
743 IF ( config_flags%open_xs .OR. config_flags%specified .OR. &
744 config_flags%nested) i_start = MAX( ids+1, its )
745 IF ( config_flags%open_ys .OR. config_flags%specified .OR. &
746 config_flags%nested) j_start = MAX( jds+1, jts )
748 IF ( config_flags%periodic_x ) i_start = its
749 IF ( config_flags%periodic_x ) i_end = MIN( ite, ide )
750 IF ( config_flags%periodic_y ) j_end = MIN( jte, jde )
752 ! Square the mapscale factor.
756 mm(i,j) = msfux(i,j) * msfuy(i,j)
760 ! Apply a coordinate transformation to vertical velocity, w. This is for both
761 ! defor13 and defor23.
763 DO j = j_start, j_end
765 DO i = i_start, i_end
766 hat(i,k,j) = w(i,k,j) / msfty(i,j)
772 DO j = j_start, MIN( jte, jde-1 )
774 hat(i,k,j) = w(i,k,j) / msfty(i,j)
780 DO i = i_start, MIN( ite, ide-1 )
781 hat(i,k,j) = w(i,k,j) / msfty(i,j)
785 ! QUESTION: What is this for?
787 DO j = j_start, j_end
789 DO i = i_start, i_end
790 hatavg(i,k,j) = 0.25 * ( &
801 DO j = j_start, j_end
803 DO i = i_start, i_end
804 tmp1(i,k,j) = ( hatavg(i,k,j) - hatavg(i,k-1,j) ) * zx(i,k,j) * &
805 0.5 * ( rdz(i,k,j) + rdz(i-1,k,j) )
810 ! End calculation of dw/dx.
811 !-----------------------------------------------------------------------
813 !-----------------------------------------------------------------------
814 ! Add the first term (dw/dx) to defor13 (dw/dx+du/dz) at vorticity
817 DO j = j_start, j_end
819 DO i = i_start, i_end
820 defor13(i,k,j) = mm(i,j) * ( &
821 rdx * ( hat(i,k,j) - hat(i-1,k,j) ) - tmp1(i,k,j) )
826 DO j = j_start, j_end
827 DO i = i_start, i_end
828 defor13(i,kts,j ) = 0.0
829 defor13(i,ktf+1,j) = 0.0
833 ! End addition of the first term to defor13.
834 !-----------------------------------------------------------------------
836 !-----------------------------------------------------------------------
839 IF ( config_flags%mix_full_fields ) THEN
841 DO j = j_start, j_end
843 DO i = i_start, i_end
844 tmp1(i,k,j) = ( u(i,k,j) - u(i,k-1,j) ) * &
845 0.5 * ( rdz(i,k,j) + rdz(i-1,k,j) )
852 DO j = j_start, j_end
854 DO i = i_start, i_end
855 tmp1(i,k,j) = ( u(i,k,j) - u_base(k) - u(i,k-1,j) + u_base(k-1) ) * &
856 0.5 * ( rdz(i,k,j) + rdz(i-1,k,j) )
863 !-----------------------------------------------------------------------
864 ! Add the second term (du/dz) to defor13 (dw/dx+du/dz) at vorticity
868 IF ( config_flags%sfs_opt .GT. 0 ) THEN ! NBA--
870 !JDM____________________________________________________________________
872 ! s13 = du/dz + dw/dx
873 ! = du/dz + (dw/dx - dz/dx*dw/dz)
874 ! = tmp1 + ______defor13______
876 ! r13 = du/dz - dw/dx
877 ! = du/dz - (dw/dx - dz/dx*dw/dz)
878 ! = tmp1 - ______defor13______
879 !_______________________________________________________________________
881 DO j = j_start, j_end
883 DO i = i_start, i_end
884 nba_rij(i,k,j,P_r13) = tmp1(i,k,j) - defor13(i,k,j)
885 defor13(i,k,j) = defor13(i,k,j) + tmp1(i,k,j)
890 DO j = j_start, j_end !change for different surface B. C.
891 DO i = i_start, i_end
892 nba_rij(i,kts ,j,P_r13) = 0.0
893 nba_rij(i,ktf+1,j,P_r13) = 0.0
897 ELSE ! NOT NBA--------------------------------------------------------
899 DO j = j_start, j_end
901 DO i = i_start, i_end
902 defor13(i,k,j) = defor13(i,k,j) + tmp1(i,k,j)
907 ENDIF ! NBA-----------------------------------------------------------
909 ! End addition of the second term to defor13.
910 !-----------------------------------------------------------------------
912 !-----------------------------------------------------------------------
916 i_end = MIN( ite, ide-1 )
918 j_end = MIN( jte, jde-1 )
920 IF ( config_flags%open_xs .OR. config_flags%specified .OR. &
921 config_flags%nested) i_start = MAX( ids+1, its )
922 IF ( config_flags%open_ys .OR. config_flags%specified .OR. &
923 config_flags%nested) j_start = MAX( jds+1, jts )
924 IF ( config_flags%periodic_y ) j_end = MIN( jte, jde )
925 IF ( config_flags%periodic_x ) i_start = its
926 IF ( config_flags%periodic_x ) i_end = MIN( ite, ide-1 )
928 ! Square mapscale factor.
932 mm(i,j) = msfvx(i,j) * msfvy(i,j)
936 ! Apply a coordinate transformation to vertical velocity, w. Added by CW 7/19/07
938 DO j = j_start, j_end
940 DO i = i_start, i_end
941 hat(i,k,j) = w(i,k,j) / msftx(i,j)
947 DO j = j_start, MIN( jte, jde-1 )
949 hat(i,k,j) = w(i,k,j) / msftx(i,j)
955 DO i = i_start, MIN( ite, ide-1 )
956 hat(i,k,j) = w(i,k,j) / msftx(i,j)
960 ! QUESTION: What is this for?
962 DO j = j_start, j_end
964 DO i = i_start, i_end
965 hatavg(i,k,j) = 0.25 * ( &
974 ! Calculate dw/dy and store in tmp1.
976 DO j = j_start, j_end
978 DO i = i_start, i_end
979 tmp1(i,k,j) = ( hatavg(i,k,j) - hatavg(i,k-1,j) ) * zy(i,k,j) * &
980 0.5 * ( rdz(i,k,j) + rdz(i,k,j-1) )
985 ! End calculation of dw/dy.
986 !-----------------------------------------------------------------------
988 !-----------------------------------------------------------------------
989 ! Add the first term (dw/dy) to defor23 (dw/dy+dv/dz) at vorticity
992 DO j = j_start, j_end
994 DO i = i_start, i_end
995 defor23(i,k,j) = mm(i,j) * ( &
996 rdy * ( hat(i,k,j) - hat(i,k,j-1) ) - tmp1(i,k,j) )
1001 DO j = j_start, j_end
1002 DO i = i_start, i_end
1003 defor23(i,kts,j ) = 0.0
1004 defor23(i,ktf+1,j) = 0.0
1008 ! End addition of the first term to defor23.
1009 !-----------------------------------------------------------------------
1011 !-----------------------------------------------------------------------
1014 IF ( config_flags%mix_full_fields ) THEN
1016 DO j = j_start, j_end
1018 DO i = i_start, i_end
1019 tmp1(i,k,j) = ( v(i,k,j) - v(i,k-1,j) ) * &
1020 0.5 * ( rdz(i,k,j) + rdz(i,k,j-1) )
1027 DO j = j_start, j_end
1029 DO i = i_start, i_end
1030 tmp1(i,k,j) = ( v(i,k,j) - v_base(k) - v(i,k-1,j) + v_base(k-1) ) * &
1031 0.5 * ( rdz(i,k,j) + rdz(i,k,j-1) )
1038 ! End calculation of dv/dz.
1039 !-----------------------------------------------------------------------
1041 !-----------------------------------------------------------------------
1042 ! Add the second term (dv/dz) to defor23 (dw/dy+dv/dz) at vorticity
1045 IF ( config_flags%sfs_opt .GT. 0 ) THEN ! NBA--
1047 !JDM___________________________________________________________________
1049 ! s23 = dv/dz + dw/dy
1050 ! = dv/dz + (dw/dy - dz/dy*dw/dz)
1051 ! tmp1 + ______defor23______
1053 ! r23 = dv/dz - dw/dy
1054 ! = dv/dz - (dw/dy - dz/dy*dw/dz)
1055 ! = tmp1 - ______defor23______
1057 ! Add tmp1 to defor23.
1059 DO j = j_start, j_end
1061 DO i = i_start, i_end
1062 nba_rij(i,k,j,P_r23) = tmp1(i,k,j) - defor23(i,k,j)
1063 defor23(i,k,j) = defor23(i,k,j) + tmp1(i,k,j)
1068 DO j = j_start, j_end
1069 DO i = i_start, i_end
1070 nba_rij(i,kts ,j,P_r23) = 0.0
1071 nba_rij(i,ktf+1,j,P_r23) = 0.0
1075 ! End addition of the second term to defor23.
1076 !-----------------------------------------------------------------------
1078 !-----------------------------------------------------------------------
1079 ! Update the boundary for defor13 and defor23 (might need to change
1082 IF ( .NOT. config_flags%periodic_x .AND. i_start .EQ. ids+1) THEN
1085 defor13(ids,k,j) = defor13(ids+1,k,j)
1086 defor23(ids,k,j) = defor23(ids+1,k,j)
1087 nba_rij(ids,k,j,P_r13) = nba_rij(ids+1,k,j,P_r13)
1088 nba_rij(ids,k,j,P_r23) = nba_rij(ids+1,k,j,P_r23)
1093 IF ( .NOT. config_flags%periodic_y .AND. j_start .EQ. jds+1) THEN
1096 defor13(i,k,jds) = defor13(i,k,jds+1)
1097 defor23(i,k,jds) = defor23(i,k,jds+1)
1098 nba_rij(i,k,jds,P_r13) = nba_rij(i,k,jds+1,P_r13)
1099 nba_rij(i,k,jds,P_r23) = nba_rij(i,k,jds+1,P_r23)
1104 IF ( .NOT. config_flags%periodic_x .AND. i_end .EQ. ide-1) THEN
1107 defor13(ide,k,j) = defor13(ide-1,k,j)
1108 defor23(ide,k,j) = defor23(ide-1,k,j)
1109 nba_rij(ide,k,j,P_r13) = nba_rij(ide-1,k,j,P_r13)
1110 nba_rij(ide,k,j,P_r23) = nba_rij(ide-1,k,j,P_r23)
1115 IF ( .NOT. config_flags%periodic_y .AND. j_end .EQ. jde-1) THEN
1118 defor13(i,k,jde) = defor13(i,k,jde-1)
1119 defor23(i,k,jde) = defor23(i,k,jde-1)
1120 nba_rij(i,k,jde,P_r13) = nba_rij(i,k,jde-1,P_r13)
1121 nba_rij(i,k,jde,P_r23) = nba_rij(i,k,jde-1,P_r23)
1126 ELSE ! NOT NBA--------------------------------------------------------
1128 ! Add tmp1 to defor23.
1130 DO j = j_start, j_end
1132 DO i = i_start, i_end
1133 defor23(i,k,j) = defor23(i,k,j) + tmp1(i,k,j)
1138 ! End addition of the second term to defor23.
1139 !-----------------------------------------------------------------------
1141 !-----------------------------------------------------------------------
1142 ! Update the boundary for defor13 and defor23 (might need to change
1145 IF ( .NOT. config_flags%periodic_x .AND. i_start .EQ. ids+1) THEN
1148 defor13(ids,k,j) = defor13(ids+1,k,j)
1149 defor23(ids,k,j) = defor23(ids+1,k,j)
1154 IF ( .NOT. config_flags%periodic_y .AND. j_start .EQ. jds+1) THEN
1157 defor13(i,k,jds) = defor13(i,k,jds+1)
1158 defor23(i,k,jds) = defor23(i,k,jds+1)
1163 IF ( .NOT. config_flags%periodic_x .AND. i_end .EQ. ide-1) THEN
1166 defor13(ide,k,j) = defor13(ide-1,k,j)
1167 defor23(ide,k,j) = defor23(ide-1,k,j)
1172 IF ( .NOT. config_flags%periodic_y .AND. j_end .EQ. jde-1) THEN
1175 defor13(i,k,jde) = defor13(i,k,jde-1)
1176 defor23(i,k,jde) = defor23(i,k,jde-1)
1181 ENDIF ! NBA-----------------------------------------------------------
1183 ! End update of boundary for defor13 and defor23.
1184 !-----------------------------------------------------------------------
1186 ! The second three (defor12, defor13, defor23) of six deformation terms
1187 ! are now calculated at vorticity points.
1188 !=======================================================================
1190 END SUBROUTINE cal_deform_and_div
1192 !=======================================================================
1193 !=======================================================================
1195 SUBROUTINE calculate_km_kh( config_flags, dt, &
1196 dampcoef, zdamp, damp_opt, &
1197 xkmh, xkmv, xkhh, xkhv, &
1198 BN2, khdif, kvdif, div, &
1199 defor11, defor22, defor33, &
1200 defor12, defor13, defor23, &
1201 tke, p8w, t8w, theta, t, p, moist, &
1202 dn, dnw, dx, dy, rdz, rdzw, isotropic, &
1203 n_moist, cf1, cf2, cf3, warm_rain, &
1207 hpbl,dlk,xkmv_meso, & !XZ
1208 ids, ide, jds, jde, kds, kde, &
1209 ims, ime, jms, jme, kms, kme, &
1210 its, ite, jts, jte, kts, kte )
1212 ! History: Sep 2003 Changes by George Bryan and Jason Knievel, NCAR
1213 ! Oct 2001 Converted to mass core by Bill Skamarock, NCAR
1216 ! Purpose: This routine calculates exchange coefficients for the TKE
1219 ! References: Klemp and Wilhelmson (JAS 1978)
1220 ! Deardorff (B-L Meteor 1980)
1221 ! Chen and Dudhia (NCAR WRF physics report 2000)
1223 !-----------------------------------------------------------------------
1224 ! Begin declarations.
1228 TYPE( grid_config_rec_type ), INTENT( IN ) &
1231 INTEGER, INTENT( IN ) &
1232 :: n_moist, damp_opt, isotropic, &
1233 ids, ide, jds, jde, kds, kde, &
1234 ims, ime, jms, jme, kms, kme, &
1235 its, ite, jts, jte, kts, kte
1237 LOGICAL, INTENT( IN ) &
1240 REAL, INTENT( IN ) &
1241 :: dx, dy, zdamp, dt, dampcoef, cf1, cf2, cf3, khdif, kvdif
1243 REAL, DIMENSION( kms:kme ), INTENT( IN ) &
1246 REAL, DIMENSION( ims:ime, kms:kme, jms:jme, n_moist ), INTENT( INOUT ) &
1249 REAL, DIMENSION( ims:ime, kms:kme, jms:jme ), INTENT( INOUT ) &
1250 :: xkmv, xkmh, xkhv, xkhh, BN2
1252 REAL, DIMENSION( ims:ime , kms:kme, jms:jme ), INTENT( IN ) &
1253 :: defor11, defor22, defor33, defor12, defor13, defor23, &
1254 div, rdz, rdzw, p8w, t8w, theta, t, p, zx, zy
1256 REAL, DIMENSION( ims:ime, kms:kme, jms:jme ), INTENT( INOUT ) &
1259 REAL, INTENT( IN ) &
1262 REAL, DIMENSION( ims:ime, jms:jme ), INTENT( IN ) &
1265 REAL, DIMENSION( ims:ime, jms:jme ), INTENT( INOUT ) &
1268 REAL, DIMENSION( ims:ime, kms:kme, jms:jme ), INTENT( IN ) &
1271 REAL, DIMENSION( ims:ime, kms:kme, jms:jme ), INTENT( INOUT ) &
1278 :: i_start, i_end, j_start, j_end, ktf, i, j, k
1281 !-----------------------------------------------------------------------
1283 ktf = MIN( kte, kde-1 )
1285 i_end = MIN( ite, ide-1 )
1287 j_end = MIN( jte, jde-1 )
1289 CALL calculate_N2( config_flags, BN2, moist, &
1290 theta, t, p, p8w, t8w, &
1291 dnw, dn, rdz, rdzw, &
1292 n_moist, cf1, cf2, cf3, warm_rain, &
1293 ids, ide, jds, jde, kds, kde, &
1294 ims, ime, jms, jme, kms, kme, &
1295 its, ite, jts, jte, kts, kte )
1297 ! Select a scheme for calculating diffusion coefficients.
1299 km_coef: SELECT CASE( config_flags%km_opt )
1302 CALL isotropic_km( config_flags, xkmh, xkmv, &
1303 xkhh, xkhv, khdif, kvdif, &
1304 ids, ide, jds, jde, kds, kde, &
1305 ims, ime, jms, jme, kms, kme, &
1306 its, ite, jts, jte, kts, kte )
1308 CALL tke_km( config_flags, xkmh, xkmv, &
1309 xkhh, xkhv, BN2, tke, p8w, t8w, theta, &
1310 rdz, rdzw, dx, dy, dt, isotropic, &
1311 mix_upper_bound, msftx, msfty, &
1312 hpbl,dlk,xkmv_meso, & ! XZ
1313 defor11,defor22,defor12,zx,zy, & ! XZ
1314 ids, ide, jds, jde, kds, kde, &
1315 ims, ime, jms, jme, kms, kme, &
1316 its, ite, jts, jte, kts, kte )
1318 CALL smag_km( config_flags, xkmh, xkmv, &
1319 xkhh, xkhv, BN2, div, &
1320 defor11, defor22, defor33, &
1321 defor12, defor13, defor23, &
1322 rdzw, dx, dy, dt, isotropic, &
1323 mix_upper_bound, msftx, msfty, &
1324 ids, ide, jds, jde, kds, kde, &
1325 ims, ime, jms, jme, kms, kme, &
1326 its, ite, jts, jte, kts, kte )
1328 CALL smag2d_km( config_flags, xkmh, xkmv, &
1329 xkhh, xkhv, defor11, defor22, defor12, &
1330 rdzw, dx, dy, msftx, msfty, &
1332 ids, ide, jds, jde, kds, kde, &
1333 ims, ime, jms, jme, kms, kme, &
1334 its, ite, jts, jte, kts, kte )
1336 CALL wrf_error_fatal( 'Please choose diffusion coefficient scheme' )
1340 IF ( damp_opt .eq. 1 ) THEN
1341 CALL cal_dampkm( config_flags, xkmh, xkhh, xkmv, xkhv, &
1342 dx, dy, dt, dampcoef, rdz, rdzw, zdamp, &
1344 ids, ide, jds, jde, kds, kde, &
1345 ims, ime, jms, jme, kms, kme, &
1346 its, ite, jts, jte, kts, kte )
1349 END SUBROUTINE calculate_km_kh
1351 !=======================================================================
1353 SUBROUTINE cal_dampkm( config_flags,xkmh,xkhh,xkmv,xkhv, &
1354 dx,dy,dt,dampcoef, &
1357 ids,ide, jds,jde, kds,kde, &
1358 ims,ime, jms,jme, kms,kme, &
1359 its,ite, jts,jte, kts,kte )
1361 !-----------------------------------------------------------------------
1362 ! Begin declarations.
1366 TYPE(grid_config_rec_type) , INTENT(IN ) :: config_flags
1368 INTEGER , INTENT(IN ) :: ids, ide, jds, jde, kds, kde, &
1369 ims, ime, jms, jme, kms, kme, &
1370 its, ite, jts, jte, kts, kte
1372 REAL , INTENT(IN ) :: zdamp,dx,dy,dt,dampcoef
1375 REAL , DIMENSION( ims:ime, kms:kme, jms:jme), INTENT(INOUT) :: xkmh , &
1380 REAL , DIMENSION( ims:ime, kms:kme, jms:jme), INTENT(IN ) :: rdz, &
1383 REAL , DIMENSION( ims:ime, jms:jme), INTENT(IN ) :: msftx, &
1387 INTEGER :: i_start, i_end, j_start, j_end, ktf, ktfm1, i, j, k
1388 REAL :: kmmax,kmmvmax,degrad90,dz,tmp
1390 REAL , DIMENSION( its:ite ) :: deltaz
1391 REAL , DIMENSION( its:ite, kts:kte, jts:jte) :: dampk,dampkv
1394 !-----------------------------------------------------------------------
1396 ktf = min(kte,kde-1)
1400 i_end = MIN(ite,ide-1)
1402 j_end = MIN(jte,jde-1)
1404 ! keep upper damping diffusion away from relaxation zones at boundaries if used
1405 IF(config_flags%specified .OR. config_flags%nested)THEN
1406 i_start = MAX(i_start,ids+config_flags%spec_bdy_width-1)
1407 i_end = MIN(i_end,ide-config_flags%spec_bdy_width)
1408 j_start = MAX(j_start,jds+config_flags%spec_bdy_width-1)
1409 j_end = MIN(j_end,jde-config_flags%spec_bdy_width)
1414 DO j = j_start, j_end
1417 DO i = i_start, i_end
1418 ! Unmodified dx used above may produce very large diffusivities
1419 ! when msftx is very large. And the above formula ignores the fact
1420 ! that dy may now be different from dx as well. Let's fix that by
1421 ! defining a "ds" as the minimum of the "real-space" (physical
1422 ! distance) values of dx and dy, and then using that smallest value
1423 ! to calculate a point-by-point kmmax
1424 ds = MIN(dx/msftx(i,j),dy/msfty(i,j))
1427 ! deltaz(i)=0.5*dnw(k)/zeta_z(i,j)
1428 ! dz=dnw(k)/zeta_z(i,j)
1433 tmp=min(deltaz(i)/zdamp,1.)
1434 dampk(i,k,j)=cos(degrad90*tmp)*cos(degrad90*tmp)*kmmax*dampcoef
1435 dampkv(i,k,j)=cos(degrad90*tmp)*cos(degrad90*tmp)*kmmvmax*dampcoef
1436 ! set upper limit on vertical K (based on horizontal K)
1437 dampkv(i,k,j)=min(dampkv(i,k,j),dampk(i,k,j))
1442 DO i = i_start, i_end
1443 ! Unmodified dx used above may produce very large diffusivities
1444 ! when msftx is very large. And the above formula ignores the fact
1445 ! that dy may now be different from dx as well. Let's fix that by
1446 ! defining a "ds" as the minimum of the "real-space" (physical
1447 ! distance) values of dx and dy, and then using that smallest value
1448 ! to calculate a point-by-point kmmax
1449 ds = MIN(dx/msftx(i,j),dy/msfty(i,j))
1452 ! deltaz(i)=deltaz(i)+dn(k)/zeta_z(i,j)
1453 ! dz=dnw(k)/zeta_z(i,j)
1455 deltaz(i) = deltaz(i) + dz
1459 tmp=min(deltaz(i)/zdamp,1.)
1460 dampk(i,k,j)=cos(degrad90*tmp)*cos(degrad90*tmp)*kmmax*dampcoef
1461 dampkv(i,k,j)=cos(degrad90*tmp)*cos(degrad90*tmp)*kmmvmax*dampcoef
1462 ! set upper limit on vertical K (based on horizontal K)
1463 dampkv(i,k,j)=min(dampkv(i,k,j),dampk(i,k,j))
1469 DO j = j_start, j_end
1471 DO i = i_start, i_end
1472 xkmh(i,k,j)=max(xkmh(i,k,j),dampk(i,k,j))
1473 xkhh(i,k,j)=max(xkhh(i,k,j),dampk(i,k,j))
1474 xkmv(i,k,j)=max(xkmv(i,k,j),dampkv(i,k,j))
1475 xkhv(i,k,j)=max(xkhv(i,k,j),dampkv(i,k,j))
1480 END SUBROUTINE cal_dampkm
1482 !=======================================================================
1483 !=======================================================================
1485 SUBROUTINE calculate_N2( config_flags, BN2, moist, &
1486 theta, t, p, p8w, t8w, &
1487 dnw, dn, rdz, rdzw, &
1488 n_moist, cf1, cf2, cf3, warm_rain, &
1489 ids, ide, jds, jde, kds, kde, &
1490 ims, ime, jms, jme, kms, kme, &
1491 its, ite, jts, jte, kts, kte )
1493 !-----------------------------------------------------------------------
1494 ! Begin declarations.
1498 TYPE( grid_config_rec_type ), INTENT( IN ) &
1501 INTEGER, INTENT( IN ) &
1503 ids, ide, jds, jde, kds, kde, &
1504 ims, ime, jms, jme, kms, kme, &
1505 its, ite, jts, jte, kts, kte
1507 LOGICAL, INTENT( IN ) &
1510 REAL, INTENT( IN ) &
1513 REAL, DIMENSION( ims:ime, kms:kme, jms:jme ), INTENT( INOUT ) &
1516 REAL, DIMENSION( ims:ime, kms:kme, jms:jme ), INTENT( IN ) &
1517 :: rdz, rdzw, theta, t, p, p8w, t8w
1519 REAL, DIMENSION( kms:kme ), INTENT( IN ) &
1522 REAL, DIMENSION( ims:ime, kms:kme, jms:jme, n_moist), INTENT( INOUT ) &
1528 :: i, j, k, ktf, ispe, ktes1, ktes2, &
1529 i_start, i_end, j_start, j_end
1532 :: coefa, thetaep1, thetaem1, qc_cr, es, tc, qlpqi, qsw, qsi, &
1533 tmpdz, xlvqv, thetaesfc, thetasfc, qvtop, qvsfc, thetatop, thetaetop
1535 REAL, DIMENSION( its:ite, jts:jte ) &
1538 REAL, DIMENSION( its:ite, kts:kte, jts:jte ) &
1542 !-----------------------------------------------------------------------
1544 qc_cr = 0.00001 ! in Kg/Kg
1546 ktf = MIN( kte, kde-1 )
1551 i_end = MIN( ite, ide-1 )
1553 j_end = MIN( jte, jde-1 )
1555 IF ( config_flags%open_xs .OR. config_flags%specified .OR. &
1556 config_flags%nested) i_start = MAX( ids+1, its )
1557 IF ( config_flags%open_xe .OR. config_flags%specified .OR. &
1558 config_flags%nested) i_end = MIN( ide-2, ite )
1559 IF ( config_flags%open_ys .OR. config_flags%specified .OR. &
1560 config_flags%nested) j_start = MAX( jds+1, jts )
1561 IF ( config_flags%open_ye .OR. config_flags%specified .OR. &
1562 config_flags%nested) j_end = MIN( jde-2 ,jte )
1563 IF ( config_flags%periodic_x ) i_start = its
1564 IF ( config_flags%periodic_x ) i_end = MIN( ite, ide-1 )
1566 IF ( P_QC .GT. PARAM_FIRST_SCALAR) THEN
1567 DO j = j_start, j_end
1569 DO i = i_start, i_end
1570 qctmp(i,k,j) = moist(i,k,j,P_QC)
1575 DO j = j_start, j_end
1577 DO i = i_start, i_end
1599 DO ispe = PARAM_FIRST_SCALAR, n_moist
1600 IF ( ispe .EQ. P_QV .OR. ispe .EQ. P_QC .OR. ispe .EQ. P_QI) THEN
1601 DO j = j_start, j_end
1603 DO i = i_start, i_end
1604 tmp1(i,k,j) = tmp1(i,k,j) + moist(i,k,j,ispe)
1609 DO j = j_start, j_end
1610 DO i = i_start, i_end
1611 tmp1sfc(i,j) = tmp1sfc(i,j) + &
1612 cf1 * moist(i,1,j,ispe) + &
1613 cf2 * moist(i,2,j,ispe) + &
1614 cf3 * moist(i,3,j,ispe)
1615 tmp1top(i,j) = tmp1top(i,j) + &
1616 moist(i,ktes1,j,ispe) + &
1617 ( moist(i,ktes1,j,ispe) - moist(i,ktes2,j,ispe) ) * &
1618 0.5 * dnw(ktes1) / dn(ktes1)
1624 ! Calculate saturation mixing ratio.
1626 DO j = j_start, j_end
1628 DO i = i_start, i_end
1629 tc = t(i,k,j) - SVPT0
1630 es = 1000.0 * SVP1 * EXP( SVP2 * tc / ( t(i,k,j) - SVP3 ) )
1631 qvs(i,k,j) = EP_2 * es / ( p(i,k,j) - es )
1636 DO j = j_start, j_end
1638 DO i = i_start, i_end
1639 tmpdz = 1.0 / rdz(i,k,j) + 1.0 / rdz(i,k+1,j)
1640 IF ( moist(i,k,j,P_QV) .GE. qvs(i,k,j) .OR. qctmp(i,k,j) .GE. qc_cr) THEN
1641 xlvqv = XLV * moist(i,k,j,P_QV)
1642 coefa = ( 1.0 + xlvqv / R_d / t(i,k,j) ) / &
1643 ( 1.0 + XLV * xlvqv / Cp / R_v / t(i,k,j) / t(i,k,j) ) / &
1645 thetaep1 = theta(i,k+1,j) * &
1646 ( 1.0 + XLV * qvs(i,k+1,j) / Cp / t(i,k+1,j) )
1647 thetaem1 = theta(i,k-1,j) * &
1648 ( 1.0 + XLV * qvs(i,k-1,j) / Cp / t(i,k-1,j) )
1649 BN2(i,k,j) = g * ( coefa * ( thetaep1 - thetaem1 ) / tmpdz - &
1650 ( tmp1(i,k+1,j) - tmp1(i,k-1,j) ) / tmpdz )
1652 BN2(i,k,j) = g * ( (theta(i,k+1,j) - theta(i,k-1,j) ) / &
1653 theta(i,k,j) / tmpdz + &
1654 1.61 * ( moist(i,k+1,j,P_QV) - moist(i,k-1,j,P_QV) ) / &
1656 ( tmp1(i,k+1,j) - tmp1(i,k-1,j) ) / tmpdz )
1663 DO j = j_start, j_end
1664 DO i = i_start, i_end
1665 tmpdz = 1.0 / rdz(i,k+1,j) + 0.5 / rdzw(i,k,j)
1666 thetasfc = T8w(i,kts,j) / ( p8w(i,k,j) / p1000mb )**( R_d / Cp )
1667 IF ( moist(i,k,j,P_QV) .GE. qvs(i,k,j) .OR. qctmp(i,k,j) .GE. qc_cr) THEN
1668 qvsfc = cf1 * qvs(i,1,j) + &
1669 cf2 * qvs(i,2,j) + &
1671 xlvqv = XLV * moist(i,k,j,P_QV)
1672 coefa = ( 1.0 + xlvqv / R_d / t(i,k,j) ) / &
1673 ( 1.0 + XLV * xlvqv / Cp / R_v / t(i,k,j) / t(i,k,j) ) / &
1675 thetaep1 = theta(i,k+1,j) * &
1676 ( 1.0 + XLV * qvs(i,k+1,j) / Cp / t(i,k+1,j) )
1677 thetaesfc = thetasfc * &
1678 ( 1.0 + XLV * qvsfc / Cp / t8w(i,kts,j) )
1679 BN2(i,k,j) = g * ( coefa * ( thetaep1 - thetaesfc ) / tmpdz - &
1680 ( tmp1(i,k+1,j) - tmp1sfc(i,j) ) / tmpdz )
1682 qvsfc = cf1 * moist(i,1,j,P_QV) + &
1683 cf2 * moist(i,2,j,P_QV) + &
1684 cf3 * moist(i,3,j,P_QV)
1685 ! BN2(i,k,j) = g * ( ( theta(i,k+1,j) - thetasfc ) / &
1686 ! theta(i,k,j) / tmpdz + &
1687 ! 1.61 * ( moist(i,k+1,j,P_QV) - qvsfc ) / &
1689 ! ( tmp1(i,k+1,j) - tmp1sfc(i,j) ) / tmpdz )
1690 !...... MARTA: change in computation of BN2 at the surface, WCS 040331
1692 tmpdz= 1./rdzw(i,k,j) ! controlare come calcola rdzw
1693 BN2(i,k,j) = g * ( ( theta(i,k+1,j) - theta(i,k,j)) / &
1694 theta(i,k,j) / tmpdz + &
1695 1.61 * ( moist(i,k+1,j,P_QV) - qvsfc ) / &
1697 ( tmp1(i,k+1,j) - tmp1sfc(i,j) ) / tmpdz )
1698 ! end of MARTA/WCS change
1705 !...... MARTA: change in computation of BN2 at the top, WCS 040331
1706 DO j = j_start, j_end
1707 DO i = i_start, i_end
1708 BN2(i,ktf,j)=BN2(i,ktf-1,j)
1711 ! end of MARTA/WCS change
1713 END SUBROUTINE calculate_N2
1715 !=======================================================================
1716 !=======================================================================
1718 SUBROUTINE isotropic_km( config_flags, &
1719 xkmh,xkmv,xkhh,xkhv,khdif,kvdif, &
1720 ids,ide, jds,jde, kds,kde, &
1721 ims,ime, jms,jme, kms,kme, &
1722 its,ite, jts,jte, kts,kte )
1724 !-----------------------------------------------------------------------
1725 ! Begin declarations.
1729 TYPE(grid_config_rec_type) , INTENT(IN ) :: config_flags
1731 INTEGER , INTENT(IN ) :: ids, ide, jds, jde, kds, kde, &
1732 ims, ime, jms, jme, kms, kme, &
1733 its, ite, jts, jte, kts, kte
1735 REAL , INTENT(IN ) :: khdif,kvdif
1737 REAL, DIMENSION( ims:ime , kms:kme , jms:jme ) , INTENT(INOUT) :: xkmh, &
1743 INTEGER :: i_start, i_end, j_start, j_end, ktf, i, j, k
1744 REAL :: khdif3,kvdif3
1747 !-----------------------------------------------------------------------
1752 i_end = MIN(ite,ide-1)
1754 j_end = MIN(jte,jde-1)
1758 khdif3=khdif/prandtl
1759 kvdif3=kvdif/prandtl
1761 DO j = j_start, j_end
1763 DO i = i_start, i_end
1772 END SUBROUTINE isotropic_km
1774 !=======================================================================
1775 !=======================================================================
1777 SUBROUTINE smag_km( config_flags,xkmh,xkmv,xkhh,xkhv,BN2, &
1778 div,defor11,defor22,defor33,defor12, &
1780 rdzw,dx,dy,dt,isotropic, &
1781 mix_upper_bound, msftx, msfty, &
1782 ids,ide, jds,jde, kds,kde, &
1783 ims,ime, jms,jme, kms,kme, &
1784 its,ite, jts,jte, kts,kte )
1786 !-----------------------------------------------------------------------
1787 ! Begin declarations.
1791 TYPE(grid_config_rec_type) , INTENT(IN ) :: config_flags
1793 INTEGER , INTENT(IN ) :: ids, ide, jds, jde, kds, kde, &
1794 ims, ime, jms, jme, kms, kme, &
1795 its, ite, jts, jte, kts, kte
1797 INTEGER , INTENT(IN ) :: isotropic
1798 REAL , INTENT(IN ) :: dx, dy, dt, mix_upper_bound
1801 REAL, DIMENSION( ims:ime , kms:kme , jms:jme ) , INTENT(IN ) :: BN2, &
1804 REAL, DIMENSION( ims:ime , kms:kme , jms:jme ) , INTENT(INOUT) :: xkmh, &
1809 REAL , DIMENSION( ims:ime , kms:kme, jms:jme ), INTENT(IN ) :: &
1817 REAL , DIMENSION( ims:ime, jms:jme), INTENT(IN ) :: msftx, &
1821 INTEGER :: i_start, i_end, j_start, j_end, ktf, i, j, k
1822 REAL :: deltas, tmp, pr, mlen_h, mlen_v, c_s
1824 REAL, DIMENSION( its:ite , kts:kte , jts:jte ) :: def2
1827 !-----------------------------------------------------------------------
1829 ktf = min(kte,kde-1)
1832 i_end = MIN(ite,ide-1)
1834 j_end = MIN(jte,jde-1)
1836 IF ( config_flags%open_xs .or. config_flags%specified .or. &
1837 config_flags%nested) i_start = MAX(ids+1,its)
1838 IF ( config_flags%open_xe .or. config_flags%specified .or. &
1839 config_flags%nested) i_end = MIN(ide-2,ite)
1840 IF ( config_flags%open_ys .or. config_flags%specified .or. &
1841 config_flags%nested) j_start = MAX(jds+1,jts)
1842 IF ( config_flags%open_ye .or. config_flags%specified .or. &
1843 config_flags%nested) j_end = MIN(jde-2,jte)
1844 IF ( config_flags%periodic_x ) i_start = its
1845 IF ( config_flags%periodic_x ) i_end = MIN( ite, ide-1 )
1848 c_s = config_flags%c_s
1853 def2(i,k,j)=0.5*(defor11(i,k,j)*defor11(i,k,j) + &
1854 defor22(i,k,j)*defor22(i,k,j) + &
1855 defor33(i,k,j)*defor33(i,k,j))
1863 tmp=0.25*(defor12(i ,k,j)+defor12(i ,k,j+1)+ &
1864 defor12(i+1,k,j)+defor12(i+1,k,j+1))
1865 def2(i,k,j)=def2(i,k,j)+tmp*tmp
1873 tmp=0.25*(defor13(i ,k+1,j)+defor13(i ,k,j)+ &
1874 defor13(i+1,k+1,j)+defor13(i+1,k,j))
1875 def2(i,k,j)=def2(i,k,j)+tmp*tmp
1883 tmp=0.25*(defor23(i,k+1,j )+defor23(i,k,j )+ &
1884 defor23(i,k+1,j+1)+defor23(i,k,j+1))
1885 def2(i,k,j)=def2(i,k,j)+tmp*tmp
1890 IF (isotropic .EQ. 0) THEN
1891 DO j = j_start, j_end
1893 DO i = i_start, i_end
1894 mlen_h=sqrt(dx/msftx(i,j) * dy/msfty(i,j))
1895 mlen_v= 1./rdzw(i,k,j)
1896 tmp=max(0.,def2(i,k,j)-BN2(i,k,j)/pr)
1898 xkmh(i,k,j)=max(c_s*c_s*mlen_h*mlen_h*tmp, 1.0E-6*mlen_h*mlen_h )
1899 xkmh(i,k,j)=min(xkmh(i,k,j), mix_upper_bound * mlen_h * mlen_h / dt )
1900 xkmv(i,k,j)=max(c_s*c_s*mlen_v*mlen_v*tmp, 1.0E-6*mlen_v*mlen_v )
1901 xkmv(i,k,j)=min(xkmv(i,k,j), mix_upper_bound * mlen_v * mlen_v / dt )
1902 xkhh(i,k,j)=xkmh(i,k,j)/pr
1903 xkhh(i,k,j)=min(xkhh(i,k,j), mix_upper_bound * mlen_h * mlen_h / dt )
1904 xkhv(i,k,j)=xkmv(i,k,j)/pr
1905 xkhv(i,k,j)=min(xkhv(i,k,j), mix_upper_bound * mlen_v * mlen_v / dt )
1910 DO j = j_start, j_end
1912 DO i = i_start, i_end
1913 deltas=(dx/msftx(i,j) * dy/msfty(i,j)/rdzw(i,k,j))**0.33333333
1914 tmp=max(0.,def2(i,k,j)-BN2(i,k,j)/pr)
1916 xkmh(i,k,j)=max(c_s*c_s*deltas*deltas*tmp, 1.0E-6*deltas*deltas )
1917 xkmh(i,k,j)=min(xkmh(i,k,j), mix_upper_bound * dx/msftx(i,j) * dy/msfty(i,j) / dt )
1918 xkmv(i,k,j)=xkmh(i,k,j)
1919 xkmv(i,k,j)=min(xkmv(i,k,j), mix_upper_bound / rdzw(i,k,j) / rdzw(i,k,j) / dt )
1920 xkhh(i,k,j)=xkmh(i,k,j)/pr
1921 xkhh(i,k,j)=min(xkhh(i,k,j), mix_upper_bound * dx/msftx(i,j) * dy/msfty(i,j) / dt )
1922 xkhv(i,k,j)=xkmv(i,k,j)/pr
1923 xkhv(i,k,j)=min(xkhv(i,k,j), mix_upper_bound / rdzw(i,k,j) / rdzw(i,k,j) / dt )
1929 END SUBROUTINE smag_km
1931 !=======================================================================
1932 !=======================================================================
1934 SUBROUTINE smag2d_km( config_flags,xkmh,xkmv,xkhh,xkhv, &
1935 defor11,defor22,defor12, &
1936 rdzw,dx,dy,msftx, msfty,zx,zy, &
1937 ids,ide, jds,jde, kds,kde, &
1938 ims,ime, jms,jme, kms,kme, &
1939 its,ite, jts,jte, kts,kte )
1941 !-----------------------------------------------------------------------
1942 ! Begin declarations.
1946 TYPE(grid_config_rec_type) , INTENT(IN ) :: config_flags
1948 INTEGER , INTENT(IN ) :: ids, ide, jds, jde, kds, kde, &
1949 ims, ime, jms, jme, kms, kme, &
1950 its, ite, jts, jte, kts, kte
1952 REAL , INTENT(IN ) :: dx, dy
1955 REAL, DIMENSION( ims:ime , kms:kme , jms:jme ) , INTENT(IN ) :: rdzw,zx,zy
1957 REAL, DIMENSION( ims:ime , kms:kme , jms:jme ) , INTENT(INOUT) :: xkmh, &
1962 REAL , DIMENSION( ims:ime , kms:kme, jms:jme ), INTENT(IN ) :: &
1967 REAL , DIMENSION( ims:ime, jms:jme), INTENT(IN ) :: msftx, &
1971 INTEGER :: i_start, i_end, j_start, j_end, ktf, i, j, k
1972 REAL :: deltas, tmp, pr, mlen_h, c_s
1973 REAL :: dxm, dym, tmpzx, tmpzy, alpha, def_limit
1975 REAL, DIMENSION( its:ite , kts:kte , jts:jte ) :: def2
1978 !-----------------------------------------------------------------------
1980 ktf = min(kte,kde-1)
1983 i_end = MIN(ite,ide-1)
1985 j_end = MIN(jte,jde-1)
1987 IF ( config_flags%open_xs .or. config_flags%specified .or. &
1988 config_flags%nested) i_start = MAX(ids+1,its)
1989 IF ( config_flags%open_xe .or. config_flags%specified .or. &
1990 config_flags%nested) i_end = MIN(ide-2,ite)
1991 IF ( config_flags%open_ys .or. config_flags%specified .or. &
1992 config_flags%nested) j_start = MAX(jds+1,jts)
1993 IF ( config_flags%open_ye .or. config_flags%specified .or. &
1994 config_flags%nested) j_end = MIN(jde-2,jte)
1995 IF ( config_flags%periodic_x ) i_start = its
1996 IF ( config_flags%periodic_x ) i_end = MIN( ite, ide-1 )
1999 c_s = config_flags%c_s
2004 def2(i,k,j)=0.25*((defor11(i,k,j)-defor22(i,k,j))*(defor11(i,k,j)-defor22(i,k,j)))
2005 tmp=0.25*(defor12(i ,k,j)+defor12(i ,k,j+1)+ &
2006 defor12(i+1,k,j)+defor12(i+1,k,j+1))
2007 def2(i,k,j)=def2(i,k,j)+tmp*tmp
2012 DO j = j_start, j_end
2014 DO i = i_start, i_end
2015 mlen_h=sqrt(dx/msftx(i,j) * dy/msfty(i,j))
2016 tmp=sqrt(def2(i,k,j))
2017 ! xkmh(i,k,j)=max(c_s*c_s*mlen_h*mlen_h*tmp, 1.0E-6*mlen_h*mlen_h )
2018 xkmh(i,k,j)=c_s*c_s*mlen_h*mlen_h*tmp
2019 xkmh(i,k,j)=min(xkmh(i,k,j), 10.*mlen_h )
2020 xkmv(i,k,j)=xkmh(i,k,j) ! v4.2 and later, this is used for hor. diff. of w
2021 xkhh(i,k,j)=xkmh(i,k,j)/pr
2023 IF(config_flags%diff_opt .EQ. 2)THEN
2024 ! jd: reduce diffusion coefficient by slope factor (modified by JB August 2014)
2027 tmpzx = (0.25*( abs(zx(i,k,j))+ abs(zx(i+1,k,j )) + abs(zx(i,k+1,j))+ abs(zx(i+1,k+1,j )))*rdzw(i,k,j)*dxm)
2028 tmpzy = (0.25*( abs(zy(i,k,j))+ abs(zy(i ,k,j+1)) + abs(zy(i,k+1,j))+ abs(zy(i ,k+1,j+1)))*rdzw(i,k,j)*dym)
2029 alpha = max(sqrt(tmpzx*tmpzx+tmpzy*tmpzy),1.0)
2030 ! If deformation is large, further reduce the diffusion coefficient
2031 def_limit = max(10./mlen_h,1.e-3)
2032 if ( tmp .gt. def_limit ) then
2033 xkmh(i,k,j)=xkmh(i,k,j)/(alpha*alpha)
2035 xkmh(i,k,j)=xkmh(i,k,j)/(alpha)
2037 xkhh(i,k,j)=xkmh(i,k,j)/pr
2038 xkmv(i,k,j)=xkmh(i,k,j) ! v4.2 and later, this is used for hor. diff. of w
2044 END SUBROUTINE smag2d_km
2046 !=======================================================================
2047 !=======================================================================
2049 SUBROUTINE tke_km( config_flags, xkmh, xkmv, xkhh, xkhv, &
2050 bn2, tke, p8w, t8w, theta, &
2051 rdz, rdzw, dx,dy, dt, isotropic, &
2052 mix_upper_bound, msftx, msfty, &
2053 hpbl,dlk,xkmv_meso, & ! XZ
2054 defor11,defor22,defor12,zx,zy, & ! XZ
2055 ids, ide, jds, jde, kds, kde, &
2056 ims, ime, jms, jme, kms, kme, &
2057 its, ite, jts, jte, kts, kte )
2059 ! History: Sep 2003 Changes by Jason Knievel and George Bryan, NCAR
2060 ! Oct 2001 Converted to mass core by Bill Skamarock, NCAR
2063 ! Purpose: This routine calculates the exchange coefficients for the
2064 ! TKE turbulence parameterization.
2066 ! References: Klemp and Wilhelmson (JAS 1978)
2067 ! Chen and Dudhia (NCAR WRF physics report 2000)
2069 !-----------------------------------------------------------------------
2070 ! Begin declarations.
2074 TYPE( grid_config_rec_type ), INTENT( IN ) &
2077 INTEGER, INTENT( IN ) &
2078 :: ids, ide, jds, jde, kds, kde, &
2079 ims, ime, jms, jme, kms, kme, &
2080 its, ite, jts, jte, kts, kte
2082 INTEGER, INTENT( IN ) :: isotropic
2083 REAL, INTENT( IN ) &
2086 REAL, DIMENSION( ims:ime, kms:kme, jms:jme ), INTENT( IN ) &
2087 :: tke, p8w, t8w, theta, rdz, rdzw, bn2
2089 REAL, DIMENSION( ims:ime, kms:kme, jms:jme ), INTENT( INOUT ) &
2090 :: xkmh, xkmv, xkhh, xkhv
2092 REAL, INTENT( IN ) &
2095 REAL , DIMENSION( ims:ime, jms:jme), INTENT(IN ) :: msftx, &
2098 REAL, DIMENSION( ims:ime, kms:kme, jms:jme ), INTENT( IN ) &
2101 REAL, DIMENSION( ims:ime, kms:kme, jms:jme ), INTENT( IN ) &
2102 :: defor11, defor22, defor12
2104 REAL, DIMENSION( ims:ime, kms:kme, jms:jme ), INTENT( INOUT ) &
2107 REAL, DIMENSION( ims:ime, jms:jme), INTENT( INOUT ) &
2110 REAL, DIMENSION( ims:ime, kms:kme, jms:jme ), INTENT( IN ) &
2115 REAL, DIMENSION( its:ite, kts:kte, jts:jte ) &
2118 REAL, DIMENSION( its:ite, kts:kte, jts:jte ) &
2122 :: deltas, tmp, mlen_s, mlen_h, mlen_v, tmpdz, &
2123 thetasfc, thetatop, minkx, pr_inv, pr_inv_h, pr_inv_v, c_k
2126 :: i_start, i_end, j_start, j_end, ktf, i, j, k
2128 REAL, PARAMETER :: tke_seed_value = 1.e-06
2130 REAL, PARAMETER :: epsilon = 1.e-10
2132 REAL :: pth1, delxy, pu1, xkmv_les, pr_inv_v_meso, pr_inv_v_les, pr
2133 REAL :: dxm, dym, tmpzx, tmpzy, alpha, def_limit, c_s, xkmh_t, xkhh_t
2134 REAL, DIMENSION( its:ite , kts:kte , jts:jte ) :: def2
2135 REAL, DIMENSION( its:ite , kts:kte , jts:jte ) :: xkmh_s, xkhh_s
2136 REAL, PARAMETER :: xkmvo = 0.0, xkhvo = 0.0 !background diffusivity
2139 !-----------------------------------------------------------------------
2141 ktf = MIN( kte, kde-1 )
2143 i_end = MIN( ite, ide-1 )
2145 j_end = MIN( jte, jde-1 )
2147 IF ( config_flags%open_xs .OR. config_flags%specified .OR. &
2148 config_flags%nested) i_start = MAX( ids+1, its )
2149 IF ( config_flags%open_xe .OR. config_flags%specified .OR. &
2150 config_flags%nested) i_end = MIN( ide-2, ite )
2151 IF ( config_flags%open_ys .OR. config_flags%specified .OR. &
2152 config_flags%nested) j_start = MAX( jds+1, jts )
2153 IF ( config_flags%open_ye .OR. config_flags%specified .OR. &
2154 config_flags%nested) j_end = MIN( jde-2, jte)
2155 IF ( config_flags%periodic_x ) i_start = its
2156 IF ( config_flags%periodic_x ) i_end = MIN( ite, ide-1 )
2158 ! in the absence of surface drag or a surface heat flux, there
2159 ! is no way to generate tke without pre-existing tke. Use
2160 ! tke_seed if the drag and flux are off.
2162 c_k = config_flags%c_k
2164 IF (config_flags%isfflx .eq. 0) THEN
2165 IF ((config_flags%diff_opt .eq. 2) .and. (config_flags%bl_pbl_physics .eq. 0)) THEN
2166 IF( (config_flags%tke_drag_coefficient .lt. epsilon) .and. &
2167 (config_flags%tke_heat_flux .lt. epsilon) ) THEN
2168 tke_seed = tke_seed_value
2171 !tke_drag_coefficient and tke_heat_flux are irrelevant here
2172 tke_seed = tke_seed_value
2176 DO j = j_start, j_end
2178 DO i = i_start, i_end
2179 tmpdz = 1.0 / rdz(i,k+1,j) + 1.0 / rdz(i,k,j)
2180 dthrdn(i,k,j) = ( theta(i,k+1,j) - theta(i,k-1,j) ) / tmpdz
2186 DO j = j_start, j_end
2187 DO i = i_start, i_end
2188 tmpdz = 1.0 / rdzw(i,k+1,j) + 1.0 / rdzw(i,k,j)
2189 thetasfc = T8w(i,kts,j) / ( p8w(i,k,j) / p1000mb )**( R_d / Cp )
2190 dthrdn(i,k,j) = ( theta(i,k+1,j) - thetasfc ) / tmpdz
2195 DO j = j_start, j_end
2196 DO i = i_start, i_end
2197 tmpdz = 1.0 / rdz(i,k,j) + 0.5 / rdzw(i,k,j)
2198 thetatop = T8w(i,kde,j) / ( p8w(i,kde,j) / p1000mb )**( R_d / Cp )
2199 dthrdn(i,k,j) = ( thetatop - theta(i,k-1,j) ) / tmpdz
2203 IF ( config_flags%km_opt .EQ. 2 .and. isotropic .EQ. 0 ) THEN !XZ
2204 DO j = j_start, j_end
2206 DO i = i_start, i_end
2207 mlen_h = SQRT( dx/msftx(i,j) * dy/msfty(i,j) )
2208 tmp = SQRT( MAX( tke(i,k,j), tke_seed ) )
2209 deltas = 1.0 / rdzw(i,k,j)
2211 IF ( dthrdn(i,k,j) .GT. 0.) THEN
2212 mlen_s = 0.76 * tmp / ( ABS( g / theta(i,k,j) * dthrdn(i,k,j) ) )**0.5
2213 mlen_v = MIN( mlen_v, mlen_s )
2215 xkmh(i,k,j) = MAX( c_k * tmp * mlen_h, 1.0E-6 * mlen_h * mlen_h )
2216 xkmh(i,k,j) = MIN( xkmh(i,k,j), mix_upper_bound * mlen_h *mlen_h / dt )
2217 xkmv(i,k,j) = MAX( c_k * tmp * mlen_v, 1.0E-6 * deltas * deltas )
2218 xkmv(i,k,j) = MIN( xkmv(i,k,j), mix_upper_bound * deltas *deltas / dt )
2219 pr_inv_h = 1./prandtl
2220 pr_inv_v = 1.0 + 2.0 * mlen_v / deltas
2221 xkhh(i,k,j) = xkmh(i,k,j) * pr_inv_h
2222 xkhv(i,k,j) = xkmv(i,k,j) * pr_inv_v
2228 IF ( config_flags%km_opt .eq.2.and.isotropic .NE. 0 ) THEN ! XZ
2229 CALL calc_l_scale( config_flags, tke, BN2, l_scale, &
2230 i_start, i_end, ktf, j_start, j_end, &
2231 dx, dy, rdzw, msftx, msfty, &
2232 ids, ide, jds, jde, kds, kde, &
2233 ims, ime, jms, jme, kms, kme, &
2234 its, ite, jts, jte, kts, kte )
2235 DO j = j_start, j_end
2237 DO i = i_start, i_end
2238 tmp = SQRT( MAX( tke(i,k,j), tke_seed ) )
2239 deltas = ( dx/msftx(i,j) * dy/msfty(i,j) / rdzw(i,k,j) )**0.33333333
2240 xkmh(i,k,j) = c_k * tmp * l_scale(i,k,j)
2241 xkmh(i,k,j) = MIN( mix_upper_bound * dx/msftx(i,j) * dy/msfty(i,j) / dt, xkmh(i,k,j) )
2242 xkmv(i,k,j) = c_k * tmp * l_scale(i,k,j)
2243 xkmv(i,k,j) = MIN( mix_upper_bound / rdzw(i,k,j) / rdzw(i,k,j) / dt , xkmv(i,k,j) )
2244 pr_inv = 1.0 + 2.0 * l_scale(i,k,j) / deltas
2245 xkhh(i,k,j) = MIN( mix_upper_bound * dx/msftx(i,j) * dy/msfty(i,j) / dt, xkmh(i,k,j) * pr_inv )
2246 xkhv(i,k,j) = MIN( mix_upper_bound / rdzw(i,k,j) / rdzw(i,k,j) / dt, xkmv(i,k,j) * pr_inv )
2253 IF(config_flags%km_opt .eq.5 ) THEN ! 3DTKE diffusivity
2255 ! smag_2d horizontal diffusivity
2257 c_s = config_flags%c_s
2262 def2(i,k,j) = 0.25*((defor11(i,k,j)-defor22(i,k,j))*(defor11(i,k,j)-defor22(i,k,j)))
2263 tmp = 0.25*(defor12(i ,k,j)+defor12(i ,k,j+1)+ &
2264 defor12(i+1,k,j)+defor12(i+1,k,j+1))
2265 def2(i,k,j)= def2(i,k,j)+tmp*tmp
2270 DO j = j_start, j_end
2272 DO i = i_start, i_end
2273 mlen_h = sqrt(dx/msftx(i,j) * dy/msfty(i,j))
2274 tmp = sqrt(def2(i,k,j))
2275 xkmh_s(i,k,j) = c_s*c_s*mlen_h*mlen_h*tmp
2276 xkmh_s(i,k,j) = min(xkmh_s(i,k,j), 10.*mlen_h )
2277 xkhh_s(i,k,j) = xkmh_s(i,k,j)/pr
2278 IF(config_flags%diff_opt .EQ. 2)THEN
2281 tmpzx = (0.25*( abs(zx(i,k,j))+ abs(zx(i+1,k,j )) + abs(zx(i,k+1,j))+ abs(zx(i+1,k+1,j )))*rdzw(i,k,j)*dxm)
2282 tmpzy = (0.25*( abs(zy(i,k,j))+ abs(zy(i ,k,j+1)) + abs(zy(i,k+1,j))+ abs(zy(i ,k+1,j+1)))*rdzw(i,k,j)*dym)
2283 alpha = max(sqrt(tmpzx*tmpzx+tmpzy*tmpzy),1.0)
2284 def_limit = max(10./mlen_h,1.e-3)
2285 IF ( tmp .GT. def_limit ) THEN
2286 xkmh_s(i,k,j)=xkmh_s(i,k,j)/(alpha*alpha)
2288 xkmh_s(i,k,j)=xkmh_s(i,k,j)/(alpha)
2290 xkhh_s(i,k,j)=xkmh_s(i,k,j)/pr
2296 !! difusivity based on TKE
2297 DO j = j_start, j_end
2299 DO i = i_start, i_end
2300 mlen_h = SQRT( dx/msftx(i,j) * dy/msfty(i,j) ) ! horizontal mixing length
2301 tmp = SQRT( MAX( tke(i,k,j), tke_seed ) )
2302 deltas = 1.0 / rdzw(i,k,j)
2304 IF ( dthrdn(i,k,j) .GT. 0.) THEN
2305 mlen_s = 0.76 * tmp / ( ABS( g / theta(i,k,j) * dthrdn(i,k,j)))**0.5
2306 mlen_v = MIN( mlen_v, mlen_s )
2308 ! length scales multiplied by partial function for grid-size dependency
2309 delxy = SQRT(dx/msftx(i,j)*dy/msfty(i,j))
2310 pth1= pthl(delxy,hpbl(i,j))
2311 pu1 = pu (delxy,hpbl(i,j))
2313 xkmh_t = c_k * tmp * mlen_h
2314 xkmv_meso(i,k,j) = 0.4 * tmp * dlk(i,k,j) ! diffusivity using mesoscale length scale
2315 xkmv_meso(i,k,j) = xkmv_meso(i,k,j) + xkmvo
2316 xkmv_les = c_k * tmp * mlen_v ! diffusivity using LES length scale
2317 xkmv(i,k,j) = ( 1.0 - pu1 ) * xkmv_les + pu1 * xkmv_meso(i,k,j)
2318 xkmv(i,k,j) = MIN(xkmv(i,k,j), 1000.)
2320 pr_inv_h = 1./prandtl
2321 pr_inv_v_les = 1.0 + 2.0 * mlen_v / deltas
2323 xkhh_t = xkmh_t * pr_inv_h
2324 xkhv(i,k,j) = xkmv(i,k,j) * (1.0-pth1) * pr_inv_v_les &
2325 + pth1 * (pr_inv_v_meso * xkmv(i,k,j) + xkhvo)
2326 xkhv(i,k,j) = MIN(xkhv(i,k,j), 1000.)
2327 ! blending of horizontal K(deformation) and K(TKE)
2328 xkmh(i,k,j) = pth1 * xkmh_s(i,k,j) + ( 1.0 - pth1 ) * xkmh_t
2329 xkhh(i,k,j) = pth1 * xkhh_s(i,k,j) + ( 1.0 - pth1 ) * xkhh_t
2336 END SUBROUTINE tke_km
2338 !=======================================================================
2339 !=======================================================================
2341 SUBROUTINE calc_l_scale( config_flags, tke, BN2, l_scale, &
2342 i_start, i_end, ktf, j_start, j_end, &
2343 dx, dy, rdzw, msftx, msfty, &
2344 ids, ide, jds, jde, kds, kde, &
2345 ims, ime, jms, jme, kms, kme, &
2346 its, ite, jts, jte, kts, kte )
2348 ! History: Sep 2003 Written by Bryan and Knievel, NCAR
2350 ! Purpose: This routine calculates the length scale, based on stability,
2351 ! for TKE parameterization of subgrid-scale turbulence.
2353 !-----------------------------------------------------------------------
2354 ! Begin declarations.
2358 TYPE( grid_config_rec_type ), INTENT( IN ) &
2361 INTEGER, INTENT( IN ) &
2362 :: i_start, i_end, ktf, j_start, j_end, &
2363 ids, ide, jds, jde, kds, kde, &
2364 ims, ime, jms, jme, kms, kme, &
2365 its, ite, jts, jte, kts, kte
2367 REAL, DIMENSION( ims:ime, kms:kme, jms:jme ), INTENT( IN ) &
2370 REAL, INTENT( IN ) &
2373 REAL, DIMENSION( its:ite, kts:kte, jts:jte ), INTENT( OUT ) &
2376 REAL , DIMENSION( ims:ime, jms:jme), INTENT(IN ) :: msftx, &
2387 !-----------------------------------------------------------------------
2389 DO j = j_start, j_end
2391 DO i = i_start, i_end
2392 deltas = ( dx/msftx(i,j) * dy/msfty(i,j) / rdzw(i,k,j) )**0.33333333
2393 l_scale(i,k,j) = deltas
2395 IF ( BN2(i,k,j) .gt. 1.0e-6 ) THEN
2396 tmp = SQRT( MAX( tke(i,k,j), 1.0e-6 ) )
2397 l_scale(i,k,j) = 0.76 * tmp / SQRT( BN2(i,k,j) )
2398 l_scale(i,k,j) = MIN( l_scale(i,k,j), deltas)
2399 l_scale(i,k,j) = MAX( l_scale(i,k,j), 0.001 * deltas )
2406 END SUBROUTINE calc_l_scale
2408 !=======================================================================
2409 !=======================================================================
2411 SUBROUTINE meso_length_scale(config_flags,dx,dy,rdzw,rdz,tke, &
2412 p8w,t8w,theta,dlk,hpbl,elmin, &
2413 rmol,rho,hfx,qfx,moist,n_moist, &
2414 ids, ide, jds, jde, kds, kde, &
2415 ims, ime, jms, jme, kms, kme, &
2416 its, ite, jts, jte, kts, kte )
2418 ! The mesoscale length scale based on Nakanishi and Niino (2009)
2419 ! and modified by X. Zhang
2423 TYPE( grid_config_rec_type ), INTENT( IN ) &
2426 INTEGER, INTENT( IN ) &
2427 :: ids, ide, jds, jde, kds, kde, &
2428 ims, ime, jms, jme, kms, kme, &
2429 its, ite, jts, jte, kts, kte
2431 INTEGER, INTENT( IN ) &
2434 REAL, INTENT( IN ) &
2437 REAL, DIMENSION( ims:ime, kms:kme, jms:jme, n_moist ), INTENT( INOUT ) &
2440 REAL, DIMENSION( ims:ime, kms:kme, jms:jme ), INTENT( IN ) &
2441 :: tke, rdzw, rdz, p8w, t8w
2443 REAL, DIMENSION( ims:ime, kms:kme, jms:jme ), INTENT( IN ) &
2446 REAL, DIMENSION( ims:ime, kms:kme, jms:jme ), INTENT( INOUT ) &
2449 REAL, DIMENSION( ims:ime, jms:jme ), INTENT( INOUT ) &
2450 :: rmol,hfx,qfx,hpbl
2453 REAL, DIMENSION( its:ite, kts:kte, jts:jte ) :: zfull,za
2454 REAL, DIMENSION( its:ite, kts:kte, jts:jte ) :: elb,qtke,els,elf
2455 REAL, DIMENSION( its:ite, jts:jte ) :: sflux,elt,vsc
2456 REAL, DIMENSION( its:ite, kts:kte, jts:jte ) :: dthrdn
2457 INTEGER :: izz,i,k,j
2458 INTEGER :: i_start, i_end, ktf, j_start, j_end
2459 REAL,PARAMETER :: alp1 = 0.8, alp2 = 1.0, alp3 = 1.0, alp4 = 100.0 !
2460 REAL :: cpm,qcv,N2,tmpdz,thetasfc,thetatop,heat_flux,gtr,qdz,coe
2461 REAL :: zi2,h1,h2,wt
2463 REAL, PARAMETER :: minzi = 300. !min mixed-layer height
2464 REAL, PARAMETER :: maxdz = 750. !max (half) transition layer depth
2465 !=0.3*2500 m PBLH, so the transition
2466 !layer stops growing for PBLHs > 2.5 km.
2467 REAL, PARAMETER :: mindz = 300. !min (half) transition layer depth
2470 ktf = MIN( kte, kde-1 )
2472 i_end = MIN( ite, ide-1 )
2474 j_end = MIN( jte, jde-1 )
2476 IF ( config_flags%open_xs .OR. config_flags%specified .OR. &
2477 config_flags%nested) i_start = MAX( ids+1, its )
2478 IF ( config_flags%open_xe .OR. config_flags%specified .OR. &
2479 config_flags%nested) i_end = MIN( ide-2, ite )
2480 IF ( config_flags%open_ys .OR. config_flags%specified .OR. &
2481 config_flags%nested) j_start = MAX( jds+1, jts )
2482 IF ( config_flags%open_ye .OR. config_flags%specified .OR. &
2483 config_flags%nested) j_end = MIN( jde-2, jte)
2484 IF ( config_flags%periodic_x ) i_start = its
2485 IF ( config_flags%periodic_x ) i_end = MIN( ite, ide-1 )
2488 DO j = j_start, j_end
2489 DO i = i_start, i_end
2494 DO j = j_start, j_end
2496 DO i = i_start, i_end
2497 zfull(i,k+1,j) = 1.0/rdzw(i,k,j) + zfull(i,k,j)
2502 DO j = j_start, j_end
2504 DO i = i_start, i_end
2505 za(i,k,j) = (zfull(i,k,j) + zfull(i,k+1,j))/2.0
2510 DO j = j_start, j_end
2512 DO i = i_start, i_end
2518 DO j = j_start, j_end
2520 DO i = i_start, i_end
2521 tmpdz = 1.0 / rdz(i,k+1,j) + 1.0 / rdz(i,k,j)
2522 dthrdn(i,k,j) = ( theta(i,k+1,j) - theta(i,k-1,j) ) / tmpdz
2528 DO j = j_start, j_end
2529 DO i = i_start, i_end
2530 tmpdz = 1.0 / rdzw(i,k+1,j) + 1.0 / rdzw(i,k,j)
2531 thetasfc = T8w(i,kts,j) / ( p8w(i,k,j) / p1000mb )**( R_d / Cp )
2532 dthrdn(i,k,j) = ( theta(i,k+1,j) - thetasfc ) / tmpdz
2537 DO j = j_start, j_end
2538 DO i = i_start, i_end
2539 tmpdz = 1.0 / rdz(i,k,j) + 0.5 / rdzw(i,k,j)
2540 thetatop = T8w(i,kde,j) / ( p8w(i,kde,j) / p1000mb )**( R_d / Cp )
2541 dthrdn(i,k,j) = ( thetatop - theta(i,k-1,j) ) / tmpdz
2545 DO j = j_start, j_end
2547 DO i = i_start, i_end
2548 qtke(i,k,j) = sqrt(2.0*tke(i,k,j))
2553 DO j = j_start, j_end
2554 DO i = i_start, i_end
2560 DO j = j_start, j_end
2562 DO i = i_start, i_end
2563 qdz = qtke(i,k,j)*(1.0/rdzw(i,k,j))
2564 elt(i,j) = elt(i,j) + qdz*za(i,k,j)
2565 vsc(i,j) = vsc(i,j) + qdz
2570 DO j = j_start, j_end
2571 DO i = i_start, i_end
2572 elt(i,j) = alp1*elt(i,j)/vsc(i,j)
2576 hflux: SELECT CASE( config_flags%isfflx )
2577 CASE (0,2) ! with fixed surface heat flux given in the namelist
2578 heat_flux = config_flags%tke_heat_flux ! constant heat flux value
2579 DO j = j_start, j_end
2580 DO i = i_start, i_end
2581 cpm = cp * (1. + 0.8 * moist(i,kts,j,P_QV))
2582 hfx(i,j)= heat_flux*cpm*rho(i,1,j)
2586 CASE (1) ! use surface heat flux computed from surface routine
2587 DO j = j_start, j_end
2588 DO i = i_start, i_end
2589 cpm = cp * (1. + 0.8 * moist(i,kts,j,P_QV))
2590 heat_flux = hfx(i,j)/cpm/rho(i,1,j)
2595 CALL wrf_error_fatal( 'isfflx value invalid for diff_opt=2' )
2598 DO j = j_start,j_end
2599 DO i = i_start,i_end
2600 cpm = cp * (1. + 0.8*moist(i,1,j,P_QV))
2601 sflux(i,j) = (hfx(i,j)/cpm)/rho(i,1,j)
2605 !-------Length scale limited by bouyancy effect-----
2607 DO j = j_start, j_end
2609 DO i = i_start, i_end
2610 IF( dthrdn(i,k,j).GT.0.0 ) THEN
2611 N2 = gtr*dthrdn(i,k,j)
2612 qcv = (gtr*MAX(sflux(i,j),0.0)*elt(i,j))**(1.0/3.0)
2613 elb(i,k,j) = qtke(i,k,j)/sqrt(N2)*(alp2 + alp3*sqrt(qcv/(elt(i,j)*sqrt(N2))))
2614 elf(i,k,j) = alp2*qtke(i,k,j)/sqrt(N2)
2617 elf(i,k,j) = elb(i,k,j)
2623 !------Length scale in the surface layer-------
2624 DO j = j_start, j_end
2626 DO i = i_start,i_end
2627 IF (rmol(i,j) .GT. 0.0) THEN
2628 els(i,k,j) = karman*za(i,k,j)/(1.0+2.7*min(za(i,k,j)*rmol(i,j),1.0))
2630 coe = (1.0 - alp4*za(i,k,j)*rmol(i,j))**0.2
2631 els(i,k,j) = 1.0*karman*za(i,k,j)*coe
2637 !----Harmonic averaging of mixing length scales-----
2638 DO j = j_start, j_end
2640 DO i = i_start, i_end
2641 dlk(i,k,j) = MIN(elb(i,k,j)/(elb(i,k,j)/elt(i,j)+elb(i,k,j)/els(i,k,j)+1.0),elf(i,k,j))
2646 !add blending to use BouLac mixing length in free atmos;
2647 !defined relative to the PBLH (zi) + transition layer (h1)
2648 DO j = j_start, j_end
2650 DO i = i_start, i_end
2651 zi2=MAX(hpbl(i,j),minzi)
2652 h1=MAX(0.3*zi2,mindz)
2653 h1=MIN(h1,maxdz) ! 1/2 transition layer depth
2654 h2=h1/2.0 ! 1/4 transition layer depth
2656 wt=.5*TANH((za(i,k,j) - (zi2+h1))/h2) + .5
2657 dlk(i,k,j) = dlk(i,k,j)*(1.-wt) + 0.4*elmin(i,k,j)*wt
2662 END SUBROUTINE meso_length_scale
2663 !=======================================================================
2664 !=======================================================================
2665 SUBROUTINE free_atmos_length(config_flags,dx,dy,rdzw,rdz,tke,theta,elmin, &
2666 hfx,qfx,moist,n_moist, &
2667 ids, ide, jds, jde, kds, kde, &
2668 ims, ime, jms, jme, kms, kme, &
2669 its, ite, jts, jte, kts, kte )
2671 ! NOTE: This subroutine is based on the BouLac PBL to calculate mixing length
2672 ! in the free atmosphere and modified for integration into the 3DTKE scheme.
2673 ! Modified by X.Zhang
2675 ! dlu = the distance a parcel can be lifted upwards give a finite
2677 ! dld = the distance a parcel can be displaced downwards given a
2678 ! finite amount of TKE.
2680 ! compute the length scales up and down
2684 TYPE( grid_config_rec_type ), INTENT( IN ) &
2687 INTEGER, INTENT( IN ) &
2688 :: ids, ide, jds, jde, kds, kde, &
2689 ims, ime, jms, jme, kms, kme, &
2690 its, ite, jts, jte, kts, kte
2692 INTEGER, INTENT( IN ) &
2695 REAL, INTENT( IN ) &
2698 REAL, DIMENSION( ims:ime, kms:kme, jms:jme, n_moist ), INTENT( INOUT ) &
2701 REAL, DIMENSION( ims:ime, kms:kme, jms:jme ), INTENT( IN ) &
2704 REAL, DIMENSION( ims:ime, kms:kme, jms:jme ), INTENT( IN ) &
2707 REAL, DIMENSION( ims:ime, kms:kme, jms:jme ), INTENT( INOUT ) &
2710 REAL, DIMENSION( ims:ime, jms:jme ), INTENT( INOUT ) &
2714 REAL, DIMENSION( its:ite, kts:kte, jts:jte ) :: zfull,za
2716 REAL, DIMENSION( its:ite, kts:kte, jts:jte ) :: dlg,dlu,dld
2718 INTEGER :: izz, found, i, k, j
2720 REAL :: dzt,zup,beta,zup_inf,bbb,tl,zdo,zdo_sup,zzz
2722 INTEGER :: i_start, i_end, ktf, j_start, j_end
2724 ktf = MIN( kte, kde-1 )
2726 i_end = MIN( ite, ide-1 )
2728 j_end = MIN( jte, jde-1 )
2730 IF ( config_flags%open_xs .OR. config_flags%specified .OR. &
2731 config_flags%nested) i_start = MAX( ids+1, its )
2732 IF ( config_flags%open_xe .OR. config_flags%specified .OR. &
2733 config_flags%nested) i_end = MIN( ide-2, ite )
2734 IF ( config_flags%open_ys .OR. config_flags%specified .OR. &
2735 config_flags%nested) j_start = MAX( jds+1, jts )
2736 IF ( config_flags%open_ye .OR. config_flags%specified .OR. &
2737 config_flags%nested) j_end = MIN( jde-2, jte)
2738 IF ( config_flags%periodic_x ) i_start = its
2739 IF ( config_flags%periodic_x ) i_end = MIN( ite, ide-1 )
2741 DO j = j_start, j_end
2742 DO i = i_start, i_end
2747 DO j = j_start, j_end
2749 DO i = i_start, i_end
2750 zfull(i,k+1,j) = 1.0/rdzw(i,k,j) + zfull(i,k,j)
2755 DO j = j_start, j_end
2757 DO i = i_start, i_end
2758 za(i,k,j) = (zfull(i,k,j) + zfull(i,k+1,j))/2.0
2763 DO j = j_start, j_end
2765 DO i = i_start, i_end
2771 DO j = j_start, j_end
2773 DO i = i_start, i_end
2775 dlu(i,k,j) = zfull(i,ktf+1,j) - zfull(i,k,j) - 1.0/rdzw(i,k,j)/2.
2778 beta=g/300. !th(i,1,j) !Buoyancy coefficient
2783 DO WHILE (found.EQ.0)
2784 IF (izz .LT. ktf) THEN
2785 dzt = (1.0/rdzw(i,izz+1,j) + 1.0/rdzw(i,izz,j))/2.
2786 zup = zup - beta * theta(i,k,j) * dzt
2787 zup = zup + beta * (theta(i,izz+1,j) + theta(i,izz,j)) * dzt/2.
2789 IF (tke(i,k,j).LT.zup.and.tke(i,k,j).GE.zup_inf) THEN
2790 bbb=(theta(i,izz+1,j)-theta(i,izz,j))/dzt
2792 tl = (-beta * (theta(i,izz,j) - theta(i,k,j)) &
2793 + sqrt(max(0.,(beta*(theta(i,izz,j)-theta(i,k,j)))**2. &
2794 + 2.*bbb*beta*(tke(i,k,j) - zup_inf))))/bbb/beta
2796 IF (theta(i,izz,j) .NE. theta(i,k,j)) THEN
2797 tl = (tke(i,k,j) - zup_inf)/(beta*(theta(i,izz,j) - theta(i,k,j)))
2802 dlu(i,k,j) = zzz - dzt + tl
2815 dld(i,k,j) = zfull(i,k,j) + 1.0/rdzw(i,k,j)/2.
2817 IF (k .GT. kts) THEN
2820 DO WHILE (found.EQ.0)
2821 IF (izz .GT. kts) THEN
2822 dzt = (1.0/rdzw(i,izz-1,j) + 1.0/rdzw(i,izz,j))/2.
2823 zdo = zdo + beta*theta(i,k,j)*dzt
2824 zdo = zdo - beta*(theta(i,izz-1,j) + theta(i,izz,j))*dzt/2.
2826 IF(tke(i,k,j) .LT. zdo .and. tke(i,k,j).GE.zdo_sup) THEN
2827 bbb = (theta(i,izz,j) - theta(i,izz-1,j))/dzt
2829 tl = (beta*(theta(i,izz,j) - theta(i,k,j)) &
2830 + sqrt(max(0.,(beta*(theta(i,izz,j) - theta(i,k,j)))**2. &
2831 + 2.*bbb*beta*(tke(i,k,j) - zdo_sup))))/bbb/beta
2833 IF(theta(i,izz,j).NE.theta(i,k,j)) THEN
2834 tl = (tke(i,k,j) - zdo_sup)/(beta*(theta(i,izz,j) - theta(i,k,j)))
2840 dld(i,k,j) = zzz - dzt + tl
2851 dlg(i,k,j) = (zfull(i,k,j)+zfull(i,k+1,j))/2.0
2852 dld(i,k,j) = min(dld(i,k,j),dlg(i,k,j))
2854 elmin(i,k,j) = min(dlu(i,k,j),dld(i,k,j))
2855 elmin(i,k,j) = elmin(i,k,j)/(1. + (elmin(i,k,j)/2000.))
2860 END SUBROUTINE free_atmos_length
2861 !=======================================================================
2862 !=======================================================================
2864 SUBROUTINE horizontal_diffusion_2 ( rt_tendf, ru_tendf, rv_tendf, rw_tendf, &
2866 moist_tendf, n_moist, &
2867 chem_tendf, n_chem, &
2868 scalar_tendf, n_scalar, &
2869 tracer_tendf, n_tracer, &
2870 thp, theta, tke, config_flags, &
2871 defor11, defor22, defor12, &
2873 nba_mij, n_nba_mij, &
2875 moist, chem, scalar,tracer, &
2876 msfux, msfuy, msfvx, msfvy, &
2877 msftx, msfty, xkmh, xkmv, xkhh, km_opt, &
2878 rdx, rdy, rdz, rdzw, fnm, fnp, &
2879 cf1, cf2, cf3, zx, zy, dn, dnw, rho, &
2880 ids, ide, jds, jde, kds, kde, &
2881 ims, ime, jms, jme, kms, kme, &
2882 its, ite, jts, jte, kts, kte )
2884 !-----------------------------------------------------------------------
2885 ! Begin declarations.
2889 TYPE(grid_config_rec_type), INTENT(IN ) :: config_flags
2891 INTEGER , INTENT(IN ) :: ids, ide, jds, jde, kds, kde, &
2892 ims, ime, jms, jme, kms, kme, &
2893 its, ite, jts, jte, kts, kte
2895 INTEGER , INTENT(IN ) :: n_moist,n_chem,n_scalar,n_tracer,km_opt
2897 REAL , INTENT(IN ) :: cf1, cf2, cf3
2899 REAL , DIMENSION( kms:kme ) , INTENT(IN ) :: fnm
2900 REAL , DIMENSION( kms:kme ) , INTENT(IN ) :: fnp
2901 REAL , DIMENSION( kms:kme ) , INTENT(IN ) :: dnw
2902 REAL , DIMENSION( kms:kme ) , INTENT(IN ) :: dn
2904 REAL , DIMENSION( ims:ime, jms:jme) , INTENT(IN ) :: msfux, &
2911 REAL , DIMENSION( ims:ime, kms:kme, jms:jme), INTENT(INOUT) ::rt_tendf,&
2917 REAL , DIMENSION( ims:ime, kms:kme, jms:jme, n_moist), &
2918 INTENT(INOUT) :: moist_tendf
2920 REAL , DIMENSION( ims:ime, kms:kme, jms:jme, n_chem), &
2921 INTENT(INOUT) :: chem_tendf
2923 REAL , DIMENSION( ims:ime, kms:kme, jms:jme, n_scalar), &
2924 INTENT(INOUT) :: scalar_tendf
2926 REAL , DIMENSION( ims:ime, kms:kme, jms:jme, n_tracer), &
2927 INTENT(INOUT) :: tracer_tendf
2929 REAL , DIMENSION( ims:ime, kms:kme, jms:jme, n_moist), &
2930 INTENT(IN ) :: moist
2932 REAL , DIMENSION( ims:ime, kms:kme, jms:jme, n_chem), &
2935 REAL , DIMENSION( ims:ime, kms:kme, jms:jme, n_scalar) , &
2936 INTENT(IN ) :: scalar
2938 REAL , DIMENSION( ims:ime, kms:kme, jms:jme, n_tracer) , &
2939 INTENT(IN ) :: tracer
2941 REAL , DIMENSION( ims:ime, kms:kme, jms:jme), INTENT(IN ) ::defor11, &
2960 REAL , INTENT(IN ) :: rdx, &
2962 INTEGER, INTENT( IN ) :: n_nba_mij
2964 REAL , DIMENSION(ims:ime, kms:kme, jms:jme, n_nba_mij), INTENT(INOUT) &
2969 INTEGER :: im, ic, is
2971 ! REAL , DIMENSION(its-1:ite+1, kts:kte, jts-1:jte+1) :: xkhh
2974 !-----------------------------------------------------------------------
2976 !-----------------------------------------------------------------------
2977 ! Call diffusion subroutines.
2979 CALL horizontal_diffusion_u_2( ru_tendf, config_flags, &
2980 defor11, defor12, div, &
2981 nba_mij, n_nba_mij, &
2983 msfux, msfuy, xkmh, rdx, rdy, fnm, fnp, &
2984 dnw, zx, zy, rdzw, rho, &
2985 ids, ide, jds, jde, kds, kde, &
2986 ims, ime, jms, jme, kms, kme, &
2987 its, ite, jts, jte, kts, kte )
2989 CALL horizontal_diffusion_v_2( rv_tendf, config_flags, &
2990 defor12, defor22, div, &
2991 nba_mij, n_nba_mij, &
2993 msfvx, msfvy, xkmh, rdx, rdy, fnm, fnp, &
2994 dnw, zx, zy, rdzw, rho, &
2995 ids, ide, jds, jde, kds, kde, &
2996 ims, ime, jms, jme, kms, kme, &
2997 its, ite, jts, jte, kts, kte )
2999 CALL horizontal_diffusion_w_2( rw_tendf, config_flags, &
3000 defor13, defor23, div, &
3001 nba_mij, n_nba_mij, &
3003 msftx, msfty, xkmv, rdx, rdy, fnm, fnp, &
3004 dn, zx, zy, rdz, rho, &
3005 ids, ide, jds, jde, kds, kde, &
3006 ims, ime, jms, jme, kms, kme, &
3007 its, ite, jts, jte, kts, kte )
3009 CALL horizontal_diffusion_s ( rt_tendf, config_flags, &
3011 msftx, msfty, msfux, msfuy, &
3012 msfvx, msfvy, xkhh, rdx, rdy, &
3013 fnm, fnp, cf1, cf2, cf3, &
3014 zx, zy, rdz, rdzw, dnw, dn, rho, &
3016 ids, ide, jds, jde, kds, kde, &
3017 ims, ime, jms, jme, kms, kme, &
3018 its, ite, jts, jte, kts, kte )
3020 IF ( ((km_opt .eq. 2) .and. (.not.config_flags%tke_mix2_off)) & !XZ
3021 .or. (km_opt .eq. 5) ) & !XZ
3022 CALL horizontal_diffusion_s ( tke_tendf(ims,kms,jms), &
3025 msftx, msfty, msfux, msfuy, &
3026 msfvx, msfvy, xkhh, rdx, rdy, &
3027 fnm, fnp, cf1, cf2, cf3, &
3028 zx, zy, rdz, rdzw, dnw, dn, rho, &
3030 ids, ide, jds, jde, kds, kde, &
3031 ims, ime, jms, jme, kms, kme, &
3032 its, ite, jts, jte, kts, kte )
3034 IF ((n_moist .ge. PARAM_FIRST_SCALAR) .and. (.not.config_flags%moist_mix2_off)) THEN
3036 moist_loop: do im = PARAM_FIRST_SCALAR, n_moist
3038 CALL horizontal_diffusion_s( moist_tendf(ims,kms,jms,im), &
3040 moist(ims,kms,jms,im), &
3041 msftx, msfty, msfux, msfuy, &
3042 msfvx, msfvy, xkhh, rdx, rdy, &
3043 fnm, fnp, cf1, cf2, cf3, &
3044 zx, zy, rdz, rdzw, dnw, dn, rho, &
3046 ids, ide, jds, jde, kds, kde, &
3047 ims, ime, jms, jme, kms, kme, &
3048 its, ite, jts, jte, kts, kte )
3054 IF ((n_chem .ge. PARAM_FIRST_SCALAR) .and. (.not.config_flags%chem_mix2_off)) THEN
3056 chem_loop: do ic = PARAM_FIRST_SCALAR, n_chem
3058 CALL horizontal_diffusion_s( chem_tendf(ims,kms,jms,ic), &
3060 chem(ims,kms,jms,ic), &
3061 msftx, msfty, msfux, msfuy, &
3062 msfvx, msfvy, xkhh, rdx, rdy, &
3063 fnm, fnp, cf1, cf2, cf3, &
3064 zx, zy, rdz, rdzw, dnw, dn, rho,&
3066 ids, ide, jds, jde, kds, kde, &
3067 ims, ime, jms, jme, kms, kme, &
3068 its, ite, jts, jte, kts, kte )
3074 IF ((n_tracer .ge. PARAM_FIRST_SCALAR) .and. (.not.config_flags%tracer_mix2_off)) THEN
3076 tracer_loop: do ic = PARAM_FIRST_SCALAR, n_tracer
3078 CALL horizontal_diffusion_s( tracer_tendf(ims,kms,jms,ic), &
3080 tracer(ims,kms,jms,ic), &
3081 msftx, msfty, msfux, msfuy, &
3082 msfvx, msfvy, xkhh, rdx, rdy, &
3083 fnm, fnp, cf1, cf2, cf3, &
3084 zx, zy, rdz, rdzw, dnw, dn, rho, &
3086 ids, ide, jds, jde, kds, kde, &
3087 ims, ime, jms, jme, kms, kme, &
3088 its, ite, jts, jte, kts, kte )
3093 IF ((n_scalar .ge. PARAM_FIRST_SCALAR) .and. (.not.config_flags%scalar_mix2_off)) THEN
3095 scalar_loop: do is = PARAM_FIRST_SCALAR, n_scalar
3097 CALL horizontal_diffusion_s( scalar_tendf(ims,kms,jms,is), &
3099 scalar(ims,kms,jms,is), &
3100 msftx, msfty, msfux, msfuy, &
3101 msfvx, msfvy, xkhh, rdx, rdy, &
3102 fnm, fnp, cf1, cf2, cf3, &
3103 zx, zy, rdz, rdzw, dnw, dn, rho, &
3105 ids, ide, jds, jde, kds, kde, &
3106 ims, ime, jms, jme, kms, kme, &
3107 its, ite, jts, jte, kts, kte )
3113 END SUBROUTINE horizontal_diffusion_2
3115 !=======================================================================
3116 !=======================================================================
3118 SUBROUTINE horizontal_diffusion_u_2( tendency, config_flags, &
3119 defor11, defor12, div, &
3120 nba_mij, n_nba_mij, &
3123 xkmh, rdx, rdy, fnm, fnp, &
3124 dnw, zx, zy, rdzw, rho, &
3125 ids, ide, jds, jde, kds, kde, &
3126 ims, ime, jms, jme, kms, kme, &
3127 its, ite, jts, jte, kts, kte )
3129 !-----------------------------------------------------------------------
3130 ! Begin declarations.
3134 TYPE(grid_config_rec_type), INTENT(IN ) :: config_flags
3136 INTEGER , INTENT(IN ) :: ids, ide, jds, jde, kds, kde, &
3137 ims, ime, jms, jme, kms, kme, &
3138 its, ite, jts, jte, kts, kte
3140 REAL , DIMENSION( kms:kme ) , INTENT(IN ) :: fnm
3141 REAL , DIMENSION( kms:kme ) , INTENT(IN ) :: fnp
3142 REAL , DIMENSION( kms:kme ) , INTENT(IN ) :: dnw
3144 REAL , DIMENSION( ims:ime, jms:jme) , INTENT(IN ) :: msfux, &
3147 REAL , DIMENSION( ims:ime, kms:kme, jms:jme), INTENT(INOUT) ::tendency
3149 REAL , DIMENSION( ims:ime, kms:kme, jms:jme), INTENT(IN ) :: rdzw, &
3152 REAL , DIMENSION( ims:ime, kms:kme, jms:jme), INTENT(IN ) ::defor11, &
3160 INTEGER, INTENT( IN ) :: n_nba_mij
3162 REAL , DIMENSION(ims:ime, kms:kme, jms:jme, n_nba_mij), INTENT(INOUT) &
3165 REAL , INTENT(IN ) :: rdx, &
3169 INTEGER :: i, j, k, ktf
3171 INTEGER :: i_start, i_end, j_start, j_end
3172 INTEGER :: is_ext,ie_ext,js_ext,je_ext
3174 REAL , DIMENSION( its-1:ite+1, kts:kte, jts-1:jte+1) :: titau1avg, &
3186 REAL :: mrdx, mrdy, rcoup
3188 REAL :: tmpzy, tmpzeta_z
3192 REAL :: term1, term2, term3
3195 !-----------------------------------------------------------------------
3199 !-----------------------------------------------------------------------
3200 ! u : p (.), u(|), w(-)
3204 ! p | . | . | . | k+1 | . | . | . | k+1
3206 ! w - 13 - - k+1 13 k+1
3208 ! p | 11 O 11 | . | k | 12 O 12 | . | k
3212 ! p | . | . | . | k-1 | . | . | . | k-1
3214 ! i-1 i i i+1 j-1 j j j+1 j+1
3220 j_end = MIN(jte,jde-1)
3222 IF ( config_flags%open_xs .or. config_flags%specified .or. &
3223 config_flags%nested) i_start = MAX(ids+1,its)
3224 IF ( config_flags%open_xe .or. config_flags%specified .or. &
3225 config_flags%nested) i_end = MIN(ide-1,ite)
3226 IF ( config_flags%open_ys .or. config_flags%specified .or. &
3227 config_flags%nested) j_start = MAX(jds+1,jts)
3228 IF ( config_flags%open_ye .or. config_flags%specified .or. &
3229 config_flags%nested) j_end = MIN(jde-2,jte)
3230 IF ( config_flags%periodic_x ) i_start = its
3231 IF ( config_flags%periodic_x ) i_end = ite
3238 CALL cal_titau_11_22_33( config_flags, titau1, &
3239 tke, xkmh, defor11, &
3240 nba_mij(ims,kms,jms,P_m11), rho, &
3241 is_ext, ie_ext, js_ext, je_ext, &
3242 ids, ide, jds, jde, kds, kde, &
3243 ims, ime, jms, jme, kms, kme, &
3244 its, ite, jts, jte, kts, kte )
3251 CALL cal_titau_12_21( config_flags, titau2, &
3253 nba_mij(ims,kms,jms,P_m12), rho, &
3254 is_ext, ie_ext, js_ext, je_ext, &
3255 ids, ide, jds, jde, kds, kde, &
3256 ims, ime, jms, jme, kms, kme, &
3257 its, ite, jts, jte, kts, kte )
3259 ! titau1avg = titau11avg
3260 ! titau2avg = titau12avg
3262 DO j = j_start, j_end
3264 DO i = i_start, i_end
3265 titau1avg(i,k,j)=0.5*(fnm(k)*(titau1(i-1,k ,j)+titau1(i,k ,j))+ &
3266 fnp(k)*(titau1(i-1,k-1,j)+titau1(i,k-1,j)))
3267 titau2avg(i,k,j)=0.5*(fnm(k)*(titau2(i,k ,j+1)+titau2(i,k ,j))+ &
3268 fnp(k)*(titau2(i,k-1,j+1)+titau2(i,k-1,j)))
3269 ! tmpzy = 0.25*( zy(i-1,k,j )+zy(i,k,j )+ &
3270 ! zy(i-1,k,j+1)+zy(i,k,j+1) )
3271 ! tmpzeta_z = 0.5*(zeta_z(i,j)+zeta_z(i-1,j))
3272 ! titau1avg(i,k,j)=titau1avg(i,k,j)*zx(i,k,j)*tmpzeta_z
3273 ! titau2avg(i,k,j)=titau2avg(i,k,j)*tmpzy *tmpzeta_z
3275 zx_at_u(i, k, j) = 0.5 * (zx(i, k, j) + zx(i, k + 1 , j))
3276 zy_at_u(i, k, j) = 0.125 * (zy(i - 1, k, j) + zy(i, k, j) + &
3277 zy(i - 1, k, j + 1) + zy(i, k, j + 1) + zy(i - 1, k + 1, j) + &
3278 zy(i, k + 1, j) + zy(i - 1, k + 1, j + 1) + zy(i, k + 1, j + 1))
3280 ! titau1avg(i,k,j)=titau1avg(i,k,j)*zx(i,k,j)
3281 ! titau2avg(i,k,j)=titau2avg(i,k,j)*tmpzy
3287 DO j = j_start, j_end
3288 DO i = i_start, i_end
3289 titau1avg(i,kts,j)=0.
3290 titau1avg(i,ktf+1,j)=0.
3291 titau2avg(i,kts,j)=0.
3292 titau2avg(i,ktf+1,j)=0.
3293 zx_at_u(i, kts, j) = 0.5 * (zx(i, kts, j) + zx(i, kts + 1 , j))
3294 zy_at_u(i, kts, j) = 0.125 * (zy(i - 1, kts, j) + zy(i, kts, j) + &
3295 zy(i - 1, kts, j + 1) + zy(i, kts, j + 1) + zy(i - 1, kts + 1, j) + &
3296 zy(i, kts + 1, j) + zy(i - 1, kts + 1, j + 1) + zy(i, kts + 1, j + 1))
3300 DO j = j_start, j_end
3302 DO i = i_start, i_end
3307 tmpdz = (1./rdzw(i,k,j)+1./rdzw(i-1,k,j))/2.
3308 tendency(i,k,j)=tendency(i,k,j) + g*tmpdz/dnw(k) * &
3309 (mrdx*(titau1(i,k,j ) - titau1(i-1,k,j)) + &
3310 mrdy*(titau2(i,k,j+1) - titau2(i ,k,j)) - &
3311 msfux(i,j)*zx_at_u(i,k,j)*(titau1avg(i,k+1,j)-titau1avg(i,k,j)) / tmpdz - &
3312 msfuy(i,j)*zy_at_u(i,k,j)*(titau2avg(i,k+1,j)-titau2avg(i,k,j)) / tmpdz &
3318 END SUBROUTINE horizontal_diffusion_u_2
3320 !=======================================================================
3321 !=======================================================================
3323 SUBROUTINE horizontal_diffusion_v_2( tendency, config_flags, &
3324 defor12, defor22, div, &
3325 nba_mij, n_nba_mij, &
3328 xkmh, rdx, rdy, fnm, fnp, &
3329 dnw, zx, zy, rdzw, rho, &
3330 ids, ide, jds, jde, kds, kde, &
3331 ims, ime, jms, jme, kms, kme, &
3332 its, ite, jts, jte, kts, kte )
3334 !-----------------------------------------------------------------------
3335 ! Begin declarations.
3339 TYPE(grid_config_rec_type), INTENT(IN ) :: config_flags
3341 INTEGER , INTENT(IN ) :: ids, ide, jds, jde, kds, kde, &
3342 ims, ime, jms, jme, kms, kme, &
3343 its, ite, jts, jte, kts, kte
3345 REAL , DIMENSION( kms:kme ) , INTENT(IN ) :: fnm
3346 REAL , DIMENSION( kms:kme ) , INTENT(IN ) :: fnp
3347 REAL , DIMENSION( kms:kme ) , INTENT(IN ) :: dnw
3349 REAL , DIMENSION( ims:ime, jms:jme) , INTENT(IN ) :: msfvx, &
3352 REAL , DIMENSION( ims:ime, kms:kme, jms:jme), INTENT(INOUT) :: tendency
3354 REAL , DIMENSION( ims:ime, kms:kme, jms:jme), INTENT(IN ) ::defor12, &
3364 INTEGER, INTENT( IN ) :: n_nba_mij
3366 REAL , DIMENSION(ims:ime, kms:kme, jms:jme, n_nba_mij), INTENT(INOUT) &
3369 REAL , INTENT(IN ) :: rdx, &
3374 INTEGER :: i, j, k, ktf
3376 INTEGER :: i_start, i_end, j_start, j_end
3377 INTEGER :: is_ext,ie_ext,js_ext,je_ext
3379 REAL , DIMENSION( its-1:ite+1, kts:kte, jts-1:jte+1) :: titau1avg, &
3391 REAL :: mrdx, mrdy, rcoup
3393 REAL :: tmpzx, tmpzeta_z
3396 !-----------------------------------------------------------------------
3400 !-----------------------------------------------------------------------
3401 ! v : p (.), v(+), w(-)
3405 ! p + . + . + . + k+1 + . + . + . + k+1
3407 ! w - 23 - - k+1 23 k+1
3409 ! p + 22 O 22 + . + k + 21 O 21 + . + k
3413 ! p + . + . + . + k-1 + . + . + . + k-1
3415 ! j-1 j j j+1 i-1 i i i+1 i+1
3419 i_end = MIN(ite,ide-1)
3423 IF ( config_flags%open_xs .or. config_flags%specified .or. &
3424 config_flags%nested) i_start = MAX(ids+1,its)
3425 IF ( config_flags%open_xe .or. config_flags%specified .or. &
3426 config_flags%nested) i_end = MIN(ide-2,ite)
3427 IF ( config_flags%open_ys .or. config_flags%specified .or. &
3428 config_flags%nested) j_start = MAX(jds+1,jts)
3429 IF ( config_flags%open_ye .or. config_flags%specified .or. &
3430 config_flags%nested) j_end = MIN(jde-1,jte)
3431 IF ( config_flags%periodic_x ) i_start = its
3432 IF ( config_flags%periodic_x ) i_end = MIN(ite,ide-1)
3439 CALL cal_titau_12_21( config_flags, titau1, &
3441 nba_mij(ims,kms,jms,P_m12),rho, &
3442 is_ext,ie_ext,js_ext,je_ext, &
3443 ids, ide, jds, jde, kds, kde, &
3444 ims, ime, jms, jme, kms, kme, &
3445 its, ite, jts, jte, kts, kte )
3452 CALL cal_titau_11_22_33( config_flags, titau2, &
3453 tke, xkmh, defor22, &
3454 nba_mij(ims,kms,jms,P_m22),rho, &
3455 is_ext, ie_ext, js_ext, je_ext, &
3456 ids, ide, jds, jde, kds, kde, &
3457 ims, ime, jms, jme, kms, kme, &
3458 its, ite, jts, jte, kts, kte )
3460 DO j = j_start, j_end
3462 DO i = i_start, i_end
3463 titau1avg(i,k,j)=0.5*(fnm(k)*(titau1(i+1,k ,j)+titau1(i,k ,j))+ &
3464 fnp(k)*(titau1(i+1,k-1,j)+titau1(i,k-1,j)))
3465 titau2avg(i,k,j)=0.5*(fnm(k)*(titau2(i,k ,j-1)+titau2(i,k ,j))+ &
3466 fnp(k)*(titau2(i,k-1,j-1)+titau2(i,k-1,j)))
3468 ! tmpzx = 0.25*( zx(i,k,j )+zx(i+1,k,j )+ &
3469 ! zx(i,k,j-1)+zx(i+1,k,j-1) )
3470 zx_at_v(i, k, j) = 0.125 * (zx(i, k, j) + zx(i + 1, k, j) + &
3471 zx(i, k, j - 1) + zx(i + 1, k, j - 1) + zx(i, k + 1, j) + &
3472 zx(i + 1, k + 1, j) + zx(i, k + 1, j - 1) + zx(i + 1, k + 1, j - 1))
3473 zy_at_v(i, k, j) = 0.5 * (zy(i, k, j) + zy(i, k + 1 , j))
3475 ! titau1avg(i,k,j)=titau1avg(i,k,j)*tmpzx
3476 ! titau2avg(i,k,j)=titau2avg(i,k,j)*zy(i,k,j)
3483 DO j = j_start, j_end
3484 DO i = i_start, i_end
3485 titau1avg(i,kts,j)=0.
3486 titau1avg(i,ktf+1,j)=0.
3487 titau2avg(i,kts,j)=0.
3488 titau2avg(i,ktf+1,j)=0.
3489 zx_at_v(i, kts, j) = 0.125 * (zx(i, kts, j) + zx(i + 1, kts, j) + &
3490 zx(i, kts, j - 1) + zx(i + 1, kts, j - 1) + zx(i, kts + 1, j) + &
3491 zx(i + 1, kts + 1, j) + zx(i, kts + 1, j - 1) + zx(i + 1, kts + 1, j - 1))
3492 zy_at_v(i, kts, j) = 0.5 * (zy(i, kts, j) + zy(i, kts + 1 , j))
3496 DO j = j_start, j_end
3498 DO i = i_start, i_end
3502 tmpdz = (1./rdzw(i,k,j)+1./rdzw(i,k,j-1))/2.
3503 tendency(i,k,j)=tendency(i,k,j) + g*tmpdz/dnw(k) * &
3504 (mrdy*(titau2(i,k,j ) - titau2(i,k,j-1)) + &
3505 mrdx*(titau1(i+1,k,j) - titau1(i ,k,j)) - &
3506 msfvx(i,j)*zx_at_v(i,k,j)*(titau1avg(i,k+1,j)-titau1avg(i,k,j)) / tmpdz - &
3507 msfvy(i,j)*zy_at_v(i,k,j)*(titau2avg(i,k+1,j)-titau2avg(i,k,j)) / tmpdz &
3514 END SUBROUTINE horizontal_diffusion_v_2
3516 !=======================================================================
3517 !=======================================================================
3519 SUBROUTINE horizontal_diffusion_w_2( tendency, config_flags, &
3520 defor13, defor23, div, &
3521 nba_mij, n_nba_mij, &
3524 xkmv, rdx, rdy, fnm, fnp, &
3525 dn, zx, zy, rdz, rho, &
3526 ids, ide, jds, jde, kds, kde, &
3527 ims, ime, jms, jme, kms, kme, &
3528 its, ite, jts, jte, kts, kte )
3530 !-----------------------------------------------------------------------
3531 ! Begin declarations.
3535 TYPE(grid_config_rec_type), INTENT(IN ) :: config_flags
3537 INTEGER , INTENT(IN ) :: ids, ide, jds, jde, kds, kde, &
3538 ims, ime, jms, jme, kms, kme, &
3539 its, ite, jts, jte, kts, kte
3541 REAL , DIMENSION( kms:kme ) , INTENT(IN ) :: fnm
3542 REAL , DIMENSION( kms:kme ) , INTENT(IN ) :: fnp
3543 REAL , DIMENSION( kms:kme ) , INTENT(IN ) :: dn
3545 REAL , DIMENSION( ims:ime, jms:jme) , INTENT(IN ) :: msftx, &
3548 REAL , DIMENSION( ims:ime, kms:kme, jms:jme), INTENT(INOUT) :: tendency
3550 REAL , DIMENSION( ims:ime, kms:kme, jms:jme), INTENT(IN ) ::defor13, &
3560 INTEGER, INTENT( IN ) :: n_nba_mij
3562 REAL , DIMENSION(ims:ime, kms:kme, jms:jme, n_nba_mij), INTENT(INOUT) &
3565 REAL , INTENT(IN ) :: rdx, &
3570 INTEGER :: i, j, k, ktf
3572 INTEGER :: i_start, i_end, j_start, j_end
3573 INTEGER :: is_ext,ie_ext,js_ext,je_ext
3575 REAL , DIMENSION( its-1:ite+1, kts:kte, jts-1:jte+1) :: titau1avg, &
3587 REAL :: mrdx, mrdy, rcoup
3589 REAL :: tmpzx, tmpzy, tmpzeta_z
3592 !-----------------------------------------------------------------------
3596 !-----------------------------------------------------------------------
3597 ! w : p (.), u(|), v(+), w(-)
3601 ! w - - - k+1 w - - - k+1
3603 ! p . | 33 | . k p . + 33 + . k
3605 ! w - 31 O 31 - k w - 32 O 32 - k
3607 ! p . | 33 | . k-1 p . | 33 | . k-1
3609 ! w - - - k-1 w - - - k-1
3611 ! i-1 i i i+1 j-1 j j j+1
3614 i_end = MIN(ite,ide-1)
3616 j_end = MIN(jte,jde-1)
3618 IF ( config_flags%open_xs .or. config_flags%specified .or. &
3619 config_flags%nested) i_start = MAX(ids+1,its)
3620 IF ( config_flags%open_xe .or. config_flags%specified .or. &
3621 config_flags%nested) i_end = MIN(ide-2,ite)
3622 IF ( config_flags%open_ys .or. config_flags%specified .or. &
3623 config_flags%nested) j_start = MAX(jds+1,jts)
3624 IF ( config_flags%open_ye .or. config_flags%specified .or. &
3625 config_flags%nested) j_end = MIN(jde-2,jte)
3626 IF ( config_flags%periodic_x ) i_start = its
3627 IF ( config_flags%periodic_x ) i_end = MIN(ite,ide-1)
3634 CALL cal_titau_13_31( config_flags, titau1, defor13, &
3635 nba_mij(ims,kms,jms,P_m13), &
3636 xkmv, fnm, fnp, rho, &
3637 is_ext, ie_ext, js_ext, je_ext, &
3638 ids, ide, jds, jde, kds, kde, &
3639 ims, ime, jms, jme, kms, kme, &
3640 its, ite, jts, jte, kts, kte )
3647 CALL cal_titau_23_32( config_flags, titau2, defor23, &
3648 nba_mij(ims,kms,jms,P_m23), &
3649 xkmv, fnm, fnp, rho, &
3650 is_ext, ie_ext, js_ext, je_ext, &
3651 ids, ide, jds, jde, kds, kde, &
3652 ims, ime, jms, jme, kms, kme, &
3653 its, ite, jts, jte, kts, kte )
3655 ! titau1avg = titau31avg * zx * zeta_z = titau13avg * zx * zeta_z
3656 ! titau2avg = titau32avg * zy * zeta_z = titau23avg * zy * zeta_z
3658 DO j = j_start, j_end
3660 DO i = i_start, i_end
3661 titau1avg(i,k,j)=0.25*(titau1(i+1,k+1,j)+titau1(i,k+1,j)+ &
3662 titau1(i+1,k ,j)+titau1(i,k ,j))
3663 titau2avg(i,k,j)=0.25*(titau2(i,k+1,j+1)+titau2(i,k+1,j)+ &
3664 titau2(i,k ,j+1)+titau2(i,k ,j))
3665 zx_at_w(i, k, j) = 0.5 * (zx(i, k, j) + zx(i + 1, k, j))
3666 zy_at_w(i, k, j) = 0.5 * (zy(i, k, j) + zy(i, k, j + 1))
3668 ! tmpzx =0.25*( zx(i,k ,j)+zx(i+1,k ,j)+ &
3669 ! zx(i,k+1,j)+zx(i+1,k+1,j) )
3670 ! tmpzy =0.25*( zy(i,k ,j)+zy(i,k ,j+1)+ &
3671 ! zy(i,k+1,j)+zy(i,k+1,j+1) )
3673 ! titau1avg(i,k,j)=titau1avg(i,k,j)*tmpzx
3674 ! titau2avg(i,k,j)=titau2avg(i,k,j)*tmpzy
3675 ! titau1avg(i,k,j)=titau1avg(i,k,j)*tmpzx*zeta_z(i,j)
3676 ! titau2avg(i,k,j)=titau2avg(i,k,j)*tmpzy*zeta_z(i,j)
3681 DO j = j_start, j_end
3682 DO i = i_start, i_end
3683 titau1avg(i,ktf+1,j)=0.
3684 titau2avg(i,ktf+1,j)=0.
3688 DO j = j_start, j_end
3690 DO i = i_start, i_end
3695 tendency(i,k,j)=tendency(i,k,j) + g/(dn(k)*rdz(i,k,j)) * &
3696 (mrdx*(titau1(i+1,k,j)-titau1(i,k,j))+ &
3697 mrdy*(titau2(i,k,j+1)-titau2(i,k,j))- &
3698 msfty(i,j)*rdz(i,k,j)*(zx_at_w(i,k,j)*(titau1avg(i,k,j)-titau1avg(i,k-1,j))+ &
3699 zy_at_w(i,k,j)*(titau2avg(i,k,j)-titau2avg(i,k-1,j)) &
3706 END SUBROUTINE horizontal_diffusion_w_2
3708 !=======================================================================
3709 !=======================================================================
3711 SUBROUTINE horizontal_diffusion_s (tendency, config_flags, &
3713 msftx, msfty, msfux, msfuy, &
3714 msfvx, msfvy, xkhh, rdx, rdy, &
3715 fnm, fnp, cf1, cf2, cf3, &
3716 zx, zy, rdz, rdzw, dnw, dn, rho, &
3718 ids, ide, jds, jde, kds, kde, &
3719 ims, ime, jms, jme, kms, kme, &
3720 its, ite, jts, jte, kts, kte )
3722 !-----------------------------------------------------------------------
3723 ! Begin declarations.
3727 TYPE(grid_config_rec_type), INTENT(IN ) :: config_flags
3729 INTEGER , INTENT(IN ) :: ids, ide, jds, jde, kds, kde, &
3730 ims, ime, jms, jme, kms, kme, &
3731 its, ite, jts, jte, kts, kte
3733 LOGICAL, INTENT(IN ) :: doing_tke
3735 REAL , INTENT(IN ) :: cf1, cf2, cf3
3737 REAL , DIMENSION( kms:kme ) , INTENT(IN ) :: fnm
3738 REAL , DIMENSION( kms:kme ) , INTENT(IN ) :: fnp
3739 REAL , DIMENSION( kms:kme ) , INTENT(IN ) :: dn
3740 REAL , DIMENSION( kms:kme ) , INTENT(IN ) :: dnw
3742 REAL , DIMENSION( ims:ime, jms:jme) , INTENT(IN ) :: msfux
3743 REAL , DIMENSION( ims:ime, jms:jme) , INTENT(IN ) :: msfuy
3744 REAL , DIMENSION( ims:ime, jms:jme) , INTENT(IN ) :: msfvx
3745 REAL , DIMENSION( ims:ime, jms:jme) , INTENT(IN ) :: msfvy
3746 REAL , DIMENSION( ims:ime, jms:jme) , INTENT(IN ) :: msftx
3747 REAL , DIMENSION( ims:ime, jms:jme) , INTENT(IN ) :: msfty
3749 REAL , DIMENSION( ims:ime, kms:kme, jms:jme), INTENT(INOUT) :: tendency
3751 REAL , DIMENSION( ims:ime, kms:kme, jms:jme), INTENT(IN ) :: &
3757 REAL , DIMENSION( ims:ime, kms:kme, jms:jme), INTENT(IN ) :: var, &
3761 REAL , INTENT(IN ) :: rdx, &
3766 INTEGER :: i, j, k, ktf
3768 INTEGER :: i_start, i_end, j_start, j_end
3770 REAL , DIMENSION( its-1:ite+1, kts:kte, jts-1:jte+1) :: H1avg, &
3778 REAL , DIMENSION( its:ite, kts:kte, jts:jte) :: tmptendf
3780 REAL :: mrdx, mrdy, rcoup
3781 REAL :: tmpzx, tmpzy, tmpzeta_z, rdzu, rdzv
3782 INTEGER :: ktes1,ktes2
3785 !-----------------------------------------------------------------------
3789 !-----------------------------------------------------------------------
3790 ! scalars: t (.), u(|), v(+), w(-)
3794 ! w - 3 - k+1 w - 3 - k+1
3796 ! t . 1 O 1 . k t . 2 O 2 . k
3798 ! w - 3 - k w - 3 - k
3800 ! t . | . | . k-1 t . + . + . k-1
3802 ! w - - - k-1 w - - - k-1
3804 ! t i-1 i i i+1 j-1 j j j+1
3811 i_end = MIN(ite,ide-1)
3813 j_end = MIN(jte,jde-1)
3815 IF ( config_flags%open_xs .or. config_flags%specified .or. &
3816 config_flags%nested) i_start = MAX(ids+1,its)
3817 IF ( config_flags%open_xe .or. config_flags%specified .or. &
3818 config_flags%nested) i_end = MIN(ide-2,ite)
3819 IF ( config_flags%open_ys .or. config_flags%specified .or. &
3820 config_flags%nested) j_start = MAX(jds+1,jts)
3821 IF ( config_flags%open_ye .or. config_flags%specified .or. &
3822 config_flags%nested) j_end = MIN(jde-2,jte)
3823 IF ( config_flags%periodic_x ) i_start = its
3824 IF ( config_flags%periodic_x ) i_end = MIN(ite,ide-1)
3826 ! diffusion of the TKE needs mutiple 2
3828 IF ( doing_tke ) THEN
3829 DO j = j_start, j_end
3831 DO i = i_start, i_end
3832 tmptendf(i,k,j)=tendency(i,k,j)
3838 ! H1 = partial var over partial x
3840 DO j = j_start, j_end
3842 DO i = i_start, i_end + 1
3843 xkxavg(i,k,j)=0.5*(xkhh(i-1,k,j)+xkhh(i,k,j))*0.5*(rho(i-1,k,j)+rho(i,k,j))
3848 DO j = j_start, j_end
3850 DO i = i_start, i_end + 1
3851 H1avg(i,k,j)=0.5*(fnm(k)*(var(i-1,k ,j)+var(i,k ,j))+ &
3852 fnp(k)*(var(i-1,k-1,j)+var(i,k-1,j)))
3857 DO j = j_start, j_end
3858 DO i = i_start, i_end + 1
3859 H1avg(i,kts ,j)=0.5*(cf1*var(i ,1,j)+cf2*var(i ,2,j)+ &
3860 cf3*var(i ,3,j)+cf1*var(i-1,1,j)+ &
3861 cf2*var(i-1,2,j)+cf3*var(i-1,3,j))
3862 H1avg(i,ktf+1,j)=0.5*(var(i,ktes1,j)+(var(i,ktes1,j)- &
3863 var(i,ktes2,j))*0.5*dnw(ktes1)/dn(ktes1)+ &
3864 var(i-1,ktes1,j)+(var(i-1,ktes1,j)- &
3865 var(i-1,ktes2,j))*0.5*dnw(ktes1)/dn(ktes1))
3869 DO j = j_start, j_end
3871 DO i = i_start, i_end + 1
3873 tmpzx = 0.5*( zx(i,k,j)+ zx(i,k+1,j))
3874 rdzu = 2./(1./rdzw(i,k,j) + 1./rdzw(i-1,k,j))
3875 H1(i,k,j)=-msfux(i,j)*xkxavg(i,k,j)*( &
3876 rdx*(var(i,k,j)-var(i-1,k,j)) - tmpzx* &
3877 (H1avg(i,k+1,j)-H1avg(i,k,j))*rdzu )
3879 ! tmpzeta_z = 0.5*(zeta_z(i,j)+zeta_z(i-1,j))
3880 ! H1(i,k,j)=-msfuy(i,j)*xkxavg(i,k,j)*( &
3881 ! rdx*(var(i,k,j)-var(i-1,k,j)) - tmpzx*tmpzeta_z* &
3882 ! (H1avg(i,k+1,j)-H1avg(i,k,j))/dnw(k))
3887 ! H2 = partial var over partial y
3889 DO j = j_start, j_end + 1
3891 DO i = i_start, i_end
3892 xkxavg(i,k,j)=0.5*(xkhh(i,k,j-1)+xkhh(i,k,j))*0.5*(rho(i,k,j-1)+rho(i,k,j))
3897 DO j = j_start, j_end + 1
3899 DO i = i_start, i_end
3900 H2avg(i,k,j)=0.5*(fnm(k)*(var(i,k ,j-1)+var(i,k ,j))+ &
3901 fnp(k)*(var(i,k-1,j-1)+var(i,k-1,j)))
3906 DO j = j_start, j_end + 1
3907 DO i = i_start, i_end
3908 H2avg(i,kts ,j)=0.5*(cf1*var(i,1,j )+cf2*var(i ,2,j)+ &
3909 cf3*var(i,3,j )+cf1*var(i,1,j-1)+ &
3910 cf2*var(i,2,j-1)+cf3*var(i,3,j-1))
3911 H2avg(i,ktf+1,j)=0.5*(var(i,ktes1,j)+(var(i,ktes1,j)- &
3912 var(i,ktes2,j))*0.5*dnw(ktes1)/dn(ktes1)+ &
3913 var(i,ktes1,j-1)+(var(i,ktes1,j-1)- &
3914 var(i,ktes2,j-1))*0.5*dnw(ktes1)/dn(ktes1))
3918 DO j = j_start, j_end + 1
3920 DO i = i_start, i_end
3922 tmpzy = 0.5*( zy(i,k,j)+ zy(i,k+1,j))
3923 rdzv = 2./(1./rdzw(i,k,j) + 1./rdzw(i,k,j-1))
3924 H2(i,k,j)=-msfvy(i,j)*xkxavg(i,k,j)*( &
3925 rdy*(var(i,k,j)-var(i,k,j-1)) - tmpzy* &
3926 (H2avg(i ,k+1,j)-H2avg(i,k,j))*rdzv)
3928 ! tmpzeta_z = 0.5*(zeta_z(i,j)+zeta_z(i,j-1))
3929 ! H2(i,k,j)=-msfvy(i,j)*xkxavg(i,k,j)*( &
3930 ! rdy*(var(i,k,j)-var(i,k,j-1)) - tmpzy*tmpzeta_z* &
3931 ! (H2avg(i ,k+1,j)-H2avg(i,k,j))/dnw(k))
3936 DO j = j_start, j_end
3938 DO i = i_start, i_end
3939 H1avg(i,k,j)=0.5*(fnm(k)*(H1(i+1,k ,j)+H1(i,k ,j))+ &
3940 fnp(k)*(H1(i+1,k-1,j)+H1(i,k-1,j)))
3941 H2avg(i,k,j)=0.5*(fnm(k)*(H2(i,k ,j+1)+H2(i,k ,j))+ &
3942 fnp(k)*(H2(i,k-1,j+1)+H2(i,k-1,j)))
3944 ! zxavg(i,k,j)=fnm(k)*zx(i,k,j)+fnp(k)*zx(i,k-1,j)
3945 ! zyavg(i,k,j)=fnm(k)*zy(i,k,j)+fnp(k)*zy(i,k-1,j)
3947 ! H1avg(i,k,j)=zx*H1avg*zeta_z
3948 ! H2avg(i,k,j)=zy*H2avg*zeta_z
3950 zx_at_m(i, k, j) = 0.25*(zx(i,k,j)+ zx(i+1,k,j) + zx(i,k+1,j)+ zx(i+1,k+1,j))
3951 zy_at_m(i, k, j) = 0.25*(zy(i,k,j)+ zy(i,k,j+1) + zy(i,k+1,j)+ zy(i,k+1,j+1))
3953 ! H1avg(i,k,j)=H1avg(i,k,j)*tmpzx*zeta_z(i,j)
3954 ! H2avg(i,k,j)=H2avg(i,k,j)*tmpzy*zeta_z(i,j)
3959 DO j = j_start, j_end
3960 DO i = i_start, i_end
3965 zx_at_m(i, kts, j) = 0.25*(zx(i,kts,j)+ zx(i+1,kts,j) + zx(i,kts+1,j)+ zx(i+1,kts+1,j))
3966 zy_at_m(i, kts, j) = 0.25*(zy(i,kts,j)+ zy(i,kts,j+1) + zy(i,kts+1,j)+ zy(i,kts+1,j+1))
3971 DO j = j_start, j_end
3973 DO i = i_start, i_end
3978 tendency(i,k,j)=tendency(i,k,j) + g/(dnw(k)*rdzw(i,k,j)) * &
3979 (mrdx*(H1(i+1,k,j)-H1(i ,k,j)) + &
3980 mrdy*(H2(i,k,j+1)-H2(i,k,j )) - &
3981 msftx(i,j)*zx_at_m(i, k, j)*(H1avg(i,k+1,j)-H1avg(i,k,j))*rdzw(i,k,j) - &
3982 msfty(i,j)*zy_at_m(i, k, j)*(H2avg(i,k+1,j)-H2avg(i,k,j))*rdzw(i,k,j) &
3988 IF ( doing_tke ) THEN
3989 DO j = j_start, j_end
3991 DO i = i_start, i_end
3992 tendency(i,k,j)=tmptendf(i,k,j)+2.* &
3993 (tendency(i,k,j)-tmptendf(i,k,j))
3999 END SUBROUTINE horizontal_diffusion_s
4001 !=======================================================================
4002 !=======================================================================
4004 SUBROUTINE vertical_diffusion_2 ( ru_tendf, rv_tendf, rw_tendf, rt_tendf, &
4005 tke_tendf, moist_tendf, n_moist, &
4006 chem_tendf, n_chem, &
4007 scalar_tendf, n_scalar, &
4008 tracer_tendf, n_tracer, &
4010 thp,u_base,v_base,t_base,qv_base,tke, &
4012 config_flags,defor13,defor23,defor33, &
4013 nba_mij, n_nba_mij, &
4015 moist,chem,scalar,tracer, &
4016 xkmv,xkhv,xkmh,km_opt, & ! xkmh added
4017 fnm, fnp, dn, dnw, rdz, rdzw, &
4018 hfx, qfx, ust, rho, &
4019 ids, ide, jds, jde, kds, kde, &
4020 ims, ime, jms, jme, kms, kme, &
4021 its, ite, jts, jte, kts, kte )
4023 !-----------------------------------------------------------------------
4024 ! Begin declarations.
4028 TYPE(grid_config_rec_type), INTENT(IN ) :: config_flags
4030 INTEGER , INTENT(IN ) :: ids, ide, jds, jde, kds, kde, &
4031 ims, ime, jms, jme, kms, kme, &
4032 its, ite, jts, jte, kts, kte
4034 INTEGER , INTENT(IN ) :: n_moist,n_chem,n_scalar,n_tracer,km_opt
4036 REAL , DIMENSION( kms:kme ) , INTENT(IN ) :: fnm
4037 REAL , DIMENSION( kms:kme ) , INTENT(IN ) :: fnp
4038 REAL , DIMENSION( kms:kme ) , INTENT(IN ) :: dnw
4039 REAL , DIMENSION( kms:kme ) , INTENT(IN ) :: dn
4041 REAL , DIMENSION( kms:kme ) , INTENT(IN ) :: qv_base
4042 REAL , DIMENSION( kms:kme ) , INTENT(IN ) :: u_base
4043 REAL , DIMENSION( kms:kme ) , INTENT(IN ) :: v_base
4044 REAL , DIMENSION( kms:kme ) , INTENT(IN ) :: t_base
4046 REAL , DIMENSION( ims:ime, kms:kme, jms:jme), INTENT(INOUT) ::ru_tendf,&
4052 REAL , DIMENSION( ims:ime, kms:kme, jms:jme, n_moist), &
4053 INTENT(INOUT) :: moist_tendf
4055 REAL , DIMENSION( ims:ime, kms:kme, jms:jme, n_chem), &
4056 INTENT(INOUT) :: chem_tendf
4058 REAL , DIMENSION( ims:ime, kms:kme, jms:jme, n_scalar) , &
4059 INTENT(INOUT) :: scalar_tendf
4060 REAL , DIMENSION( ims:ime, kms:kme, jms:jme, n_tracer) , &
4061 INTENT(INOUT) :: tracer_tendf
4063 REAL , DIMENSION( ims:ime, kms:kme, jms:jme, n_moist), &
4064 INTENT(INOUT) :: moist
4066 REAL , DIMENSION( ims:ime, kms:kme, jms:jme, n_chem), &
4067 INTENT(INOUT) :: chem
4069 REAL , DIMENSION( ims:ime, kms:kme, jms:jme, n_scalar) , &
4070 INTENT(IN ) :: scalar
4071 REAL , DIMENSION( ims:ime, kms:kme, jms:jme, n_tracer) , &
4072 INTENT(IN ) :: tracer
4074 REAL , DIMENSION( ims:ime, kms:kme, jms:jme), INTENT(IN ) ::defor13, &
4088 INTEGER, INTENT( IN ) :: n_nba_mij
4090 REAL , DIMENSION(ims:ime, kms:kme, jms:jme, n_nba_mij), INTENT(INOUT) &
4093 REAL , DIMENSION( ims:ime, kms:kme, jms:jme), INTENT(IN ) :: rho
4094 REAL , DIMENSION( ims:ime, jms:jme), INTENT(INOUT) :: hfx, &
4096 REAL , DIMENSION( ims:ime, jms:jme), INTENT(IN ) :: ust
4097 REAL , DIMENSION( ims:ime, kms:kme, jms:jme), INTENT(INOUT) :: thp
4101 REAL , DIMENSION( ims:ime, kms:kme, jms:jme) :: var_mix
4103 INTEGER :: im, i,j,k
4104 INTEGER :: i_start, i_end, j_start, j_end
4106 ! REAL , DIMENSION( its:ite, kts:kte, jts:jte) :: xkhv
4108 !***************************************************************************
4109 !***************************************************************************
4110 !MODIFICA VARIABILI PER I FLUSSI
4112 REAL :: V0_u,V0_v,tao_xz,tao_yz,ustar,cd0
4113 REAL :: xsfc,psi1,vk2,zrough,lnz
4114 REAL :: heat_flux, moist_flux, heat_flux0
4117 !FINE MODIFICA VARIABILI PER I FLUSSI
4118 !***************************************************************************
4122 !-----------------------------------------------------------------------
4125 i_end = MIN(ite,ide-1)
4127 j_end = MIN(jte,jde-1)
4129 !-----------------------------------------------------------------------
4131 CALL vertical_diffusion_u_2( ru_tendf, config_flags, &
4133 nba_mij, n_nba_mij, &
4134 dnw, rdzw, fnm, fnp, rho, &
4135 ids, ide, jds, jde, kds, kde, &
4136 ims, ime, jms, jme, kms, kme, &
4137 its, ite, jts, jte, kts, kte )
4140 CALL vertical_diffusion_v_2( rv_tendf, config_flags, &
4142 nba_mij, n_nba_mij, &
4143 dnw, rdzw, fnm, fnp, rho, &
4144 ids, ide, jds, jde, kds, kde, &
4145 ims, ime, jms, jme, kms, kme, &
4146 its, ite, jts, jte, kts, kte )
4148 CALL vertical_diffusion_w_2( rw_tendf, config_flags, &
4149 defor33, tke(ims,kms,jms), &
4150 nba_mij, n_nba_mij, &
4152 dn, rdz, fnm, fnp, rho, &
4153 ids, ide, jds, jde, kds, kde, &
4154 ims, ime, jms, jme, kms, kme, &
4155 its, ite, jts, jte, kts, kte )
4157 !*****************************************
4158 !*****************************************
4159 ! MODIFICA al flusso di momento alla parete
4161 vflux: SELECT CASE( config_flags%isfflx )
4162 CASE (0) ! Assume cd a constant, specified in namelist
4163 cd0 = config_flags%tke_drag_coefficient ! constant drag coefficient
4164 ! set in namelist.input
4166 !calcolo del modulo della velocita
4167 DO j = j_start, j_end
4171 V0_u= sqrt((u_2(i,kts,j)**2) + &
4172 (((v_2(i ,kts,j )+ &
4175 v_2(i-1,kts,j+1))/4)**2))+epsilon
4177 tao_xz=cd0*V0_u*u_2(i,kts,j)
4178 ru_tendf(i,kts,j)=ru_tendf(i,kts,j) + g*tao_xz*0.5*(rho(i,kts,j)+rho(i-1,kts,j))/dnw(kts)
4179 IF ( (config_flags%m_opt .EQ. 1) .OR. (config_flags%sfs_opt .GT. 0) ) THEN
4180 nba_mij(i,kts,j,P_m13) = -tao_xz
4186 DO i = i_start, i_end
4189 V0_v= sqrt((v_2(i,kts,j)**2) + &
4190 (((u_2(i ,kts,j )+ &
4193 u_2(i+1,kts,j-1))/4)**2))+epsilon
4195 tao_yz=cd0*V0_v*v_2(i,kts,j)
4196 rv_tendf(i,kts,j)=rv_tendf(i,kts,j) + g*tao_yz*0.5*(rho(i,kts,j)+rho(i,kts,j-1))/dnw(kts)
4197 IF ( (config_flags%m_opt .EQ. 1) .OR. (config_flags%sfs_opt .GT. 0) ) THEN
4198 nba_mij(i,kts,j,P_m23) = -tao_yz
4203 CASE (1,2) ! ustar computed from surface routine
4204 DO j = j_start, j_end
4208 V0_u= sqrt((u_2(i,kts,j)**2) + &
4209 (((v_2(i ,kts,j )+ &
4212 v_2(i-1,kts,j+1))/4)**2))+epsilon
4213 ustar=0.5*(ust(i,j)+ust(i-1,j))
4215 tao_xz=ustar*ustar*u_2(i,kts,j)/V0_u
4216 ru_tendf(i,kts,j)=ru_tendf(i,kts,j) + g*tao_xz*0.5*(rho(i,kts,j)+rho(i-1,kts,j))/dnw(kts)
4217 IF ( (config_flags%m_opt .EQ. 1) .OR. (config_flags%sfs_opt .GT. 0) ) THEN
4218 nba_mij(i,kts,j,P_m13) = -tao_xz
4224 DO i = i_start, i_end
4227 V0_v= sqrt((v_2(i,kts,j)**2) + &
4228 (((u_2(i ,kts,j )+ &
4231 u_2(i+1,kts,j-1))/4)**2))+epsilon
4232 ustar=0.5*(ust(i,j)+ust(i,j-1))
4234 tao_yz=ustar*ustar*v_2(i,kts,j)/V0_v
4235 rv_tendf(i,kts,j)=rv_tendf(i,kts,j) + g*tao_yz*0.5*(rho(i,kts,j)+rho(i,kts,j-1))/dnw(kts)
4236 IF ( (config_flags%m_opt .EQ. 1) .OR. (config_flags%sfs_opt .GT. 0) ) THEN
4237 nba_mij(i,kts,j,P_m23) = -tao_yz
4243 CALL wrf_error_fatal( 'isfflx value invalid for diff_opt=2' )
4247 ! FINE MODIFICA al flusso di momento alla parete
4248 !*****************************************
4249 !*****************************************
4251 IF ( config_flags%mix_full_fields ) THEN
4253 DO j=jts,min(jte,jde-1)
4255 DO i=its,min(ite,ide-1)
4256 var_mix(i,k,j) = thp(i,k,j)
4263 DO j=jts,min(jte,jde-1)
4265 DO i=its,min(ite,ide-1)
4266 var_mix(i,k,j) = thp(i,k,j) - t_base(k)
4273 CALL vertical_diffusion_s( rt_tendf, config_flags, var_mix, xkhv, &
4274 dn, dnw, rdz, rdzw, fnm, fnp, rho, &
4276 ids, ide, jds, jde, kds, kde, &
4277 ims, ime, jms, jme, kms, kme, &
4278 its, ite, jts, jte, kts, kte )
4281 !*****************************************
4282 !*****************************************
4283 !MODIFICA al flusso di calore
4286 hflux: SELECT CASE( config_flags%isfflx )
4287 CASE (0,2) ! with fixed surface heat flux given in the namelist
4288 heat_flux = config_flags%tke_heat_flux ! constant heat flux value
4289 ! set in namelist.input
4290 DO j = j_start, j_end
4291 DO i = i_start, i_end
4292 cpm = cp * (1. + 0.8 * moist(i,kts,j,P_QV))
4293 hfx(i,j)=heat_flux*cpm*rho(i,kts,j) ! provided for output only
4294 if(config_flags%use_theta_m == 1)THEN
4295 rt_tendf(i,kts,j)=rt_tendf(i,kts,j) &
4296 -g*heat_flux*(1.+rvovrd*moist(i,kts,j,P_QV))*rho(i,kts,j)/dnw(kts) &
4297 -g*1.61*theta(i,kts,j)*qfx(i,j)/dnw(kts)
4299 rt_tendf(i,kts,j)=rt_tendf(i,kts,j) &
4300 -g*heat_flux*rho(i,kts,j)/dnw(kts)
4305 CASE (1) ! use surface heat flux computed from surface routine
4306 DO j = j_start, j_end
4307 DO i = i_start, i_end
4309 cpm = cp * (1. + 0.8 * moist(i,kts,j,P_QV))
4310 heat_flux = hfx(i,j)/cpm
4311 if(config_flags%use_theta_m == 1)THEN
4312 rt_tendf(i,kts,j)=rt_tendf(i,kts,j) &
4313 -g*heat_flux*(1.+rvovrd*moist(i,kts,j,P_QV))/dnw(kts) &
4314 -g*1.61*theta(i,kts,j)*qfx(i,j)/dnw(kts)
4316 rt_tendf(i,kts,j)=rt_tendf(i,kts,j) &
4317 -g*heat_flux/dnw(kts)
4324 CALL wrf_error_fatal( 'isfflx value invalid for diff_opt=2' )
4328 ! FINE MODIFICA al flusso di calore
4329 !*****************************************
4330 !*****************************************
4332 If (km_opt .eq. 2) then
4333 CALL vertical_diffusion_s( tke_tendf(ims,kms,jms), &
4334 config_flags, tke(ims,kms,jms), &
4336 dn, dnw, rdz, rdzw, fnm, fnp, rho, &
4338 ids, ide, jds, jde, kds, kde, &
4339 ims, ime, jms, jme, kms, kme, &
4340 its, ite, jts, jte, kts, kte )
4343 IF (n_moist .ge. PARAM_FIRST_SCALAR) THEN
4345 moist_loop: do im = PARAM_FIRST_SCALAR, n_moist
4347 IF ( (.not. config_flags%mix_full_fields) .and. (im == P_QV) ) THEN
4349 DO j=jts,min(jte,jde-1)
4351 DO i=its,min(ite,ide-1)
4352 var_mix(i,k,j) = moist(i,k,j,im) - qv_base(k)
4359 DO j=jts,min(jte,jde-1)
4361 DO i=its,min(ite,ide-1)
4362 var_mix(i,k,j) = moist(i,k,j,im)
4370 CALL vertical_diffusion_s( moist_tendf(ims,kms,jms,im), &
4371 config_flags, var_mix, &
4373 dn, dnw, rdz, rdzw, fnm, fnp, rho, &
4375 ids, ide, jds, jde, kds, kde, &
4376 ims, ime, jms, jme, kms, kme, &
4377 its, ite, jts, jte, kts, kte )
4379 !*****************************************
4380 !*****************************************
4381 !MODIFICATIONS for water vapor flux
4384 qflux: SELECT CASE( config_flags%isfflx )
4387 CASE (1,2) ! with surface moisture flux
4388 IF ( im == P_QV ) THEN
4389 DO j = j_start, j_end
4390 DO i = i_start, i_end
4391 moist_flux = qfx(i,j)
4392 moist_tendf(i,kts,j,im)=moist_tendf(i,kts,j,im) &
4393 -g*moist_flux/dnw(kts)
4399 CALL wrf_error_fatal( 'isfflx value invalid for diff_opt=2' )
4402 ! END MODIFICATIONS for water vapor flux
4403 !*****************************************
4404 !*****************************************
4409 IF (n_chem .ge. PARAM_FIRST_SCALAR) THEN
4411 chem_loop: do im = PARAM_FIRST_SCALAR, n_chem
4413 CALL vertical_diffusion_s( chem_tendf(ims,kms,jms,im), &
4414 config_flags, chem(ims,kms,jms,im), &
4416 dn, dnw, rdz, rdzw, fnm, fnp, rho, &
4418 ids, ide, jds, jde, kds, kde, &
4419 ims, ime, jms, jme, kms, kme, &
4420 its, ite, jts, jte, kts, kte )
4425 IF (n_tracer .ge. PARAM_FIRST_SCALAR) THEN
4427 tracer_loop: do im = PARAM_FIRST_SCALAR, n_tracer
4429 CALL vertical_diffusion_s( tracer_tendf(ims,kms,jms,im), &
4430 config_flags, tracer(ims,kms,jms,im), &
4432 dn, dnw, rdz, rdzw, fnm, fnp, rho, &
4434 ids, ide, jds, jde, kds, kde, &
4435 ims, ime, jms, jme, kms, kme, &
4436 its, ite, jts, jte, kts, kte )
4442 IF (n_scalar .ge. PARAM_FIRST_SCALAR) THEN
4444 scalar_loop: do im = PARAM_FIRST_SCALAR, n_scalar
4446 CALL vertical_diffusion_s( scalar_tendf(ims,kms,jms,im), &
4447 config_flags, scalar(ims,kms,jms,im), &
4449 dn, dnw, rdz, rdzw, fnm, fnp, rho, &
4451 ids, ide, jds, jde, kds, kde, &
4452 ims, ime, jms, jme, kms, kme, &
4453 its, ite, jts, jte, kts, kte )
4458 END SUBROUTINE vertical_diffusion_2
4460 !=======================================================================
4461 !=======================================================================
4463 SUBROUTINE vertical_diffusion_u_2( tendency, config_flags, &
4465 nba_mij, n_nba_mij, &
4466 dnw, rdzw, fnm, fnp, rho, &
4467 ids, ide, jds, jde, kds, kde, &
4468 ims, ime, jms, jme, kms, kme, &
4469 its, ite, jts, jte, kts, kte )
4471 !-----------------------------------------------------------------------
4472 ! Begin declarations.
4476 TYPE(grid_config_rec_type), INTENT(IN ) :: config_flags
4478 INTEGER , INTENT(IN ) :: ids, ide, jds, jde, kds, kde, &
4479 ims, ime, jms, jme, kms, kme, &
4480 its, ite, jts, jte, kts, kte
4482 REAL , DIMENSION( kms:kme ) , INTENT(IN ) :: fnm
4483 REAL , DIMENSION( kms:kme ) , INTENT(IN ) :: fnp
4484 REAL , DIMENSION( kms:kme ) , INTENT(IN ) :: dnw
4485 ! REAL , DIMENSION( ims:ime , jms:jme ) , INTENT(IN ) :: zeta_z
4487 REAL , DIMENSION( ims:ime, kms:kme, jms:jme), INTENT(INOUT) ::tendency
4489 REAL , DIMENSION( ims:ime , kms:kme, jms:jme ) , &
4490 INTENT(IN ) ::defor13, &
4495 INTEGER, INTENT( IN ) :: n_nba_mij
4497 REAL , DIMENSION(ims:ime, kms:kme, jms:jme, n_nba_mij), INTENT(INOUT) &
4502 INTEGER :: i, j, k, ktf
4504 INTEGER :: i_start, i_end, j_start, j_end
4505 INTEGER :: is_ext,ie_ext,js_ext,je_ext
4507 REAL , DIMENSION( its-1:ite+1, kts:kte, jts-1:jte+1) :: titau3
4509 REAL , DIMENSION( its:ite, jts:jte) :: zzavg
4514 !-----------------------------------------------------------------------
4521 j_end = MIN(jte,jde-1)
4523 IF ( config_flags%open_xs .or. config_flags%specified .or. &
4524 config_flags%nested) i_start = MAX(ids+1,its)
4525 IF ( config_flags%open_xe .or. config_flags%specified .or. &
4526 config_flags%nested) i_end = MIN(ide-1,ite)
4527 IF ( config_flags%open_ys .or. config_flags%specified .or. &
4528 config_flags%nested) j_start = MAX(jds+1,jts)
4529 IF ( config_flags%open_ye .or. config_flags%specified .or. &
4530 config_flags%nested) j_end = MIN(jde-2,jte)
4531 IF ( config_flags%periodic_x ) i_start = its
4532 IF ( config_flags%periodic_x ) i_end = ite
4539 CALL cal_titau_13_31( config_flags, titau3, defor13, &
4540 nba_mij(ims,kms,jms,P_m13), &
4541 xkmv, fnm, fnp, rho, &
4542 is_ext, ie_ext, js_ext, je_ext, &
4543 ids, ide, jds, jde, kds, kde, &
4544 ims, ime, jms, jme, kms, kme, &
4545 its, ite, jts, jte, kts, kte )
4547 DO j = j_start, j_end
4549 DO i = i_start, i_end
4552 tendency(i,k,j)=tendency(i,k,j)-rdzu*(titau3(i,k+1,j)-titau3(i,k,j))
4559 ! we will pick up the surface drag (titau3(i,kts,j)) later
4561 DO j = j_start, j_end
4563 DO i = i_start, i_end
4566 tendency(i,k,j)=tendency(i,k,j)-rdzu*(titau3(i,k+1,j))
4571 END SUBROUTINE vertical_diffusion_u_2
4573 !=======================================================================
4574 !=======================================================================
4576 SUBROUTINE vertical_diffusion_v_2( tendency, config_flags, &
4578 nba_mij, n_nba_mij, &
4579 dnw, rdzw, fnm, fnp, rho, &
4580 ids, ide, jds, jde, kds, kde, &
4581 ims, ime, jms, jme, kms, kme, &
4582 its, ite, jts, jte, kts, kte )
4583 !-----------------------------------------------------------------------
4584 ! Begin declarations.
4588 TYPE(grid_config_rec_type), INTENT(IN ) :: config_flags
4590 INTEGER , INTENT(IN ) :: ids, ide, jds, jde, kds, kde, &
4591 ims, ime, jms, jme, kms, kme, &
4592 its, ite, jts, jte, kts, kte
4593 REAL , DIMENSION( kms:kme ) , INTENT(IN ) :: fnm
4594 REAL , DIMENSION( kms:kme ) , INTENT(IN ) :: fnp
4595 REAL , DIMENSION( kms:kme ) , INTENT(IN ) :: dnw
4596 ! REAL , DIMENSION( ims:ime , jms:jme ) , INTENT(IN ) :: zeta_z
4598 REAL , DIMENSION( ims:ime, kms:kme, jms:jme), INTENT(INOUT) ::tendency
4600 REAL , DIMENSION( ims:ime , kms:kme, jms:jme ) , &
4601 INTENT(IN ) ::defor23, &
4606 INTEGER, INTENT( IN ) :: n_nba_mij
4608 REAL , DIMENSION(ims:ime, kms:kme, jms:jme, n_nba_mij), INTENT(INOUT) &
4613 INTEGER :: i, j, k, ktf
4615 INTEGER :: i_start, i_end, j_start, j_end
4616 INTEGER :: is_ext,ie_ext,js_ext,je_ext
4618 REAL , DIMENSION( its-1:ite+1, kts:kte, jts-1:jte+1) :: titau3
4620 REAL , DIMENSION( its:ite, jts:jte) :: zzavg
4625 !-----------------------------------------------------------------------
4630 i_end = MIN(ite,ide-1)
4634 IF ( config_flags%open_xs .or. config_flags%specified .or. &
4635 config_flags%nested) i_start = MAX(ids+1,its)
4636 IF ( config_flags%open_xe .or. config_flags%specified .or. &
4637 config_flags%nested) i_end = MIN(ide-2,ite)
4638 IF ( config_flags%open_ys .or. config_flags%specified .or. &
4639 config_flags%nested) j_start = MAX(jds+1,jts)
4640 IF ( config_flags%open_ye .or. config_flags%specified .or. &
4641 config_flags%nested) j_end = MIN(jde-1,jte)
4642 IF ( config_flags%periodic_x ) i_start = its
4643 IF ( config_flags%periodic_x ) i_end = MIN(ite,ide-1)
4650 CALL cal_titau_23_32( config_flags, titau3, defor23, &
4651 nba_mij(ims,kms,jms,P_m23), &
4652 xkmv, fnm, fnp, rho, &
4653 is_ext, ie_ext, js_ext, je_ext, &
4654 ids, ide, jds, jde, kds, kde, &
4655 ims, ime, jms, jme, kms, kme, &
4656 its, ite, jts, jte, kts, kte )
4658 DO j = j_start, j_end
4660 DO i = i_start, i_end
4663 tendency(i,k,j)=tendency(i,k,j)-rdzv*(titau3(i,k+1,j)-titau3(i,k,j))
4670 ! we will pick up the surface drag (titau3(i,kts,j)) later
4672 DO j = j_start, j_end
4674 DO i = i_start, i_end
4677 tendency(i,k,j)=tendency(i,k,j)-rdzv*(titau3(i,k+1,j))
4683 END SUBROUTINE vertical_diffusion_v_2
4685 !=======================================================================
4686 !=======================================================================
4688 SUBROUTINE vertical_diffusion_w_2(tendency, config_flags, &
4690 nba_mij, n_nba_mij, &
4692 dn, rdz, fnm, fnp, rho, &
4693 ids, ide, jds, jde, kds, kde, &
4694 ims, ime, jms, jme, kms, kme, &
4695 its, ite, jts, jte, kts, kte )
4697 !-----------------------------------------------------------------------
4698 ! Begin declarations.
4702 TYPE(grid_config_rec_type), INTENT(IN ) :: config_flags
4704 INTEGER , INTENT(IN ) :: ids, ide, jds, jde, kds, kde, &
4705 ims, ime, jms, jme, kms, kme, &
4706 its, ite, jts, jte, kts, kte
4708 REAL , DIMENSION( kms:kme ) , INTENT(IN ) :: dn, fnm, fnp
4710 REAL , DIMENSION( ims:ime, kms:kme, jms:jme), INTENT(INOUT) ::tendency
4712 REAL , DIMENSION( ims:ime , kms:kme, jms:jme ) , &
4713 INTENT(IN ) ::defor33, &
4720 INTEGER, INTENT( IN ) :: n_nba_mij
4722 REAL , DIMENSION(ims:ime, kms:kme, jms:jme, n_nba_mij), INTENT(INOUT) &
4727 INTEGER :: i, j, k, ktf
4729 INTEGER :: i_start, i_end, j_start, j_end
4730 INTEGER :: is_ext,ie_ext,js_ext,je_ext
4732 REAL , DIMENSION( its-1:ite+1, kts:kte, jts-1:jte+1) :: titau3
4735 !-----------------------------------------------------------------------
4740 i_end = MIN(ite,ide-1)
4742 j_end = MIN(jte,jde-1)
4744 IF ( config_flags%open_xs .or. config_flags%specified .or. &
4745 config_flags%nested) i_start = MAX(ids+1,its)
4746 IF ( config_flags%open_xe .or. config_flags%specified .or. &
4747 config_flags%nested) i_end = MIN(ide-2,ite)
4748 IF ( config_flags%open_ys .or. config_flags%specified .or. &
4749 config_flags%nested) j_start = MAX(jds+1,jts)
4750 IF ( config_flags%open_ye .or. config_flags%specified .or. &
4751 config_flags%nested) j_end = MIN(jde-2,jte)
4752 IF ( config_flags%periodic_x ) i_start = its
4753 IF ( config_flags%periodic_x ) i_end = MIN(ite,ide-1)
4760 CALL cal_titau_11_22_33( config_flags, titau3, &
4761 tke, xkmh, defor33, &
4762 nba_mij(ims,kms,jms,P_m33), rho, &
4763 is_ext, ie_ext, js_ext, je_ext, &
4764 ids, ide, jds, jde, kds, kde, &
4765 ims, ime, jms, jme, kms, kme, &
4766 its, ite, jts, jte, kts, kte )
4768 ! DO j = j_start, j_end
4770 ! DO i = i_start, i_end
4771 ! titau3(i,k,j)=titau3(i,k,j)*zeta_z(i,j)
4776 DO j = j_start, j_end
4778 DO i = i_start, i_end
4779 tendency(i,k,j)=tendency(i,k,j)+ g*(titau3(i,k,j)-titau3(i,k-1,j))/dn(k)
4784 END SUBROUTINE vertical_diffusion_w_2
4786 !=======================================================================
4787 !=======================================================================
4789 SUBROUTINE vertical_diffusion_s( tendency, config_flags, var, xkhv, &
4790 dn, dnw, rdz, rdzw, fnm, fnp, rho, &
4792 ids, ide, jds, jde, kds, kde, &
4793 ims, ime, jms, jme, kms, kme, &
4794 its, ite, jts, jte, kts, kte )
4796 !-----------------------------------------------------------------------
4797 ! Begin declarations.
4801 TYPE(grid_config_rec_type), INTENT(IN ) :: config_flags
4803 INTEGER , INTENT(IN ) :: ids, ide, jds, jde, kds, kde, &
4804 ims, ime, jms, jme, kms, kme, &
4805 its, ite, jts, jte, kts, kte
4807 LOGICAL, INTENT(IN ) :: doing_tke
4809 REAL , DIMENSION( kms:kme ) , INTENT(IN ) :: fnm
4810 REAL , DIMENSION( kms:kme ) , INTENT(IN ) :: fnp
4811 REAL , DIMENSION( kms:kme ) , INTENT(IN ) :: dn
4812 REAL , DIMENSION( kms:kme ) , INTENT(IN ) :: dnw
4814 REAL , DIMENSION( ims:ime, kms:kme, jms:jme), INTENT(INOUT) ::tendency
4816 REAL , DIMENSION( ims:ime , kms:kme, jms:jme ) , INTENT(IN) :: xkhv
4818 REAL , DIMENSION( ims:ime , kms:kme, jms:jme ) , &
4819 INTENT(IN ) :: var, &
4825 INTEGER :: i, j, k, ktf
4827 INTEGER :: i_start, i_end, j_start, j_end
4829 REAL , DIMENSION( its:ite, kts:kte, jts:jte) :: H3, &
4833 REAL , DIMENSION( its:ite, kts:kte, jts:jte) :: tmptendf
4836 !-----------------------------------------------------------------------
4841 i_end = MIN(ite,ide-1)
4843 j_end = MIN(jte,jde-1)
4845 IF ( config_flags%open_xs .or. config_flags%specified .or. &
4846 config_flags%nested) i_start = MAX(ids+1,its)
4847 IF ( config_flags%open_xe .or. config_flags%specified .or. &
4848 config_flags%nested) i_end = MIN(ide-2,ite)
4849 IF ( config_flags%open_ys .or. config_flags%specified .or. &
4850 config_flags%nested) j_start = MAX(jds+1,jts)
4851 IF ( config_flags%open_ye .or. config_flags%specified .or. &
4852 config_flags%nested) j_end = MIN(jde-2,jte)
4853 IF ( config_flags%periodic_x ) i_start = its
4854 IF ( config_flags%periodic_x ) i_end = MIN(ite,ide-1)
4857 DO j = j_start, j_end
4859 DO i = i_start, i_end
4860 tmptendf(i,k,j)=tendency(i,k,j)
4870 DO j = j_start, j_end
4872 DO i = i_start, i_end
4873 xkxavg(i,k,j)=fnm(k)*xkhv(i,k,j)+fnp(k)*xkhv(i,k-1,j)
4874 xkxavg(i,k,j)=xkxavg(i,k,j)*(fnm(k)*rho(i,k,j)+fnp(k)*rho(i,k-1,j))
4875 H3(i,k,j)=-xkxavg(i,k,j)*(var(i,k,j)-var(i,k-1,j))*rdz(i,k,j)
4880 DO j = j_start, j_end
4881 DO i = i_start, i_end
4887 DO j = j_start, j_end
4889 DO i = i_start, i_end
4890 tendency(i,k,j)=tendency(i,k,j) &
4891 + g * (H3(i,k+1,j)-H3(i,k,j))/dnw(k)
4897 DO j = j_start, j_end
4899 DO i = i_start, i_end
4900 tendency(i,k,j)=tmptendf(i,k,j)+2.* &
4901 (tendency(i,k,j)-tmptendf(i,k,j))
4907 END SUBROUTINE vertical_diffusion_s
4909 !=======================================================================
4910 !=======================================================================
4911 SUBROUTINE nonlocal_flux (config_flags,nlflux,gamu,gamv,hpbl,kpbl, &
4912 dx,dy,dt,ust,hfx,qfx,br,ep1,ep2, &
4913 karman,u_phy,v_phy,th_phy,rho,moist,n_moist, &
4916 ids, ide, jds, jde, kds, kde, &
4917 ims, ime, jms, jme, kms, kme, &
4918 its, ite, jts, jte, kts, kte )
4919 !-----------------------------------------------------------------------
4920 !This subroutine prescibes the nonlocal heat flux profile based on LES analysis
4921 !And compute the nonlocal momentum gamma term
4923 ! Shin and Hong 2015, MWR
4924 ! Xu Zhang et al. 2018, MWR
4925 !-----------------------------------------------------------------------
4926 ! Begin declarations.
4930 TYPE(grid_config_rec_type), INTENT(IN ) :: config_flags
4932 INTEGER , INTENT(IN ) :: ids, ide, jds, jde, kds, kde, &
4933 ims, ime, jms, jme, kms, kme, &
4934 its, ite, jts, jte, kts, kte
4936 INTEGER , INTENT(IN ) :: n_moist
4938 REAL, INTENT(IN ) :: ep1,ep2,karman
4940 REAL, INTENT(IN ) :: dx,dy,dt
4942 REAL, DIMENSION( ims:ime, kms:kme, jms:jme ), INTENT(IN) :: &
4945 REAL, DIMENSION( ims:ime, kms:kme, jms:jme ), INTENT(IN) :: rdzw
4947 REAL, DIMENSION( ims:ime, kms:kme, jms:jme, n_moist ), INTENT( INOUT ) :: &
4950 REAL, DIMENSION( ims:ime, kms:kme, jms:jme ), INTENT(IN) :: rho
4952 REAL, DIMENSION( ims:ime, jms:jme ), INTENT(INOUT) :: hfx,qfx,br,ust,hpbl
4954 REAL, DIMENSION( ims:ime, jms:jme ), INTENT(IN ) :: u10,v10,wspd
4956 REAL, DIMENSION( ims:ime, jms:jme ), INTENT(IN ) :: msftx,msfty
4958 REAL, DIMENSION( ims:ime, kms:kme, jms:jme ), INTENT(OUT) :: nlflux
4960 REAL, DIMENSION( ims:ime, jms:jme ), INTENT(OUT ) :: gamu,gamv
4962 INTEGER, DIMENSION( ims:ime, jms:jme ), INTENT(INOUT ) :: kpbl
4965 REAL, DIMENSION( its:ite, kts:kte, jts:jte ) :: zq
4966 REAL, DIMENSION( its:ite, kts:kte-1, jts:jte ) :: za,thv
4967 REAL, DIMENSION( its:ite, kts:kte, jts:jte ) :: zfacmf,entfacmf
4968 REAL, DIMENSION( its:ite, jts:jte ) :: govrth,sflux,wstar3,wstar,rigs
4969 LOGICAL, DIMENSION( its:ite,jts:jte ) :: pblflg, sfcflg, stable
4970 REAL, DIMENSION( its:ite, jts:jte ) :: deltaoh,we,enlfrac2
4971 REAL, DIMENSION( its:ite, jts:jte ) :: hfxpbl,bfxpbl
4972 REAL, DIMENSION( its:ite, jts:jte ) :: dthv,wm2
4973 REAL, DIMENSION( its:ite, jts:jte ) :: wscale,thermal
4974 REAL :: tvcon,delb,cpm,wm3,bfx0,ust3
4975 REAL :: dtheta,du,dv,thsfc,brint
4976 INTEGER :: i,j,k,i_start,i_end,j_start,j_end,ktf
4977 REAL :: mlfrac,ezfrac,sfcfracn,sflux0,snlflux0
4978 REAL :: amf1,amf2,amf3,bmf2,bmf3,pth1,delxy,pu1
4980 REAL,PARAMETER :: phifac = 8.,sfcfrac = 0.1,d1 = 0.02, d2 = 0.05, zfmin = 1.e-8
4981 REAL,PARAMETER :: h1 = 0.33333335, h2 = 0.6666667, tmin = 1.e-2
4982 ! tunable parameters
4983 REAL,PARAMETER :: mltop = 1.0, sfcfracn1 = 0.075
4984 REAL,PARAMETER :: nlfrac = 0.68 ! the ratio of nonlocal heat flux to total heat flux
4985 REAL,PARAMETER :: enlfrac = -0.15
4986 REAL,PARAMETER :: a11 = 1.0, a12 = -1.15
4987 REAL,PARAMETER :: ezfac = 1.5
4988 REAL,PARAMETER :: cpent = -0.4, rigsmax = 100., rimin = -100.
4989 REAL,PARAMETER :: entfmin = 1.0, entfmax = 5.0, sm = 10.9
4994 i_end = MIN(ite,ide-1)
4996 j_end = MIN(jte,jde-1)
4998 IF ( config_flags%open_xs .or. config_flags%specified .or. &
4999 config_flags%nested) i_start = MAX(ids+1,its)
5000 IF ( config_flags%open_xe .or. config_flags%specified .or. &
5001 config_flags%nested) i_end = MIN(ide-2,ite)
5002 IF ( config_flags%open_ys .or. config_flags%specified .or. &
5003 config_flags%nested) j_start = MAX(jds+1,jts)
5004 IF ( config_flags%open_ye .or. config_flags%specified .or. &
5005 config_flags%nested) j_end = MIN(jde-2,jte)
5006 IF ( config_flags%periodic_x ) i_start = its
5007 IF ( config_flags%periodic_x ) i_end = MIN(ite,ide-1)
5009 DO j = j_start, j_end
5010 DO i = i_start, i_end
5015 DO j = j_start, j_end
5017 DO i = i_start, i_end
5018 zq(i,k+1,j) = 1.0/rdzw(i,k,j) + zq(i,k,j)
5023 DO j = j_start,j_end
5025 DO i = i_start,i_end
5026 za(i,k,j) = 0.5*(zq(i,k,j) + zq(i,k+1,j))
5031 DO j = j_start, j_end
5032 DO i = i_start, i_end
5039 DO j = j_start, j_end
5041 DO i = i_start, i_end
5047 DO j = j_start, j_end
5049 DO i = i_start, i_end
5055 DO j = j_start, j_end
5057 DO i = i_start, i_end
5058 tvcon = 1. + ep1*moist(i,k,j,P_QV)
5059 thv(i,k,j) = th_phy(i,k,j)*tvcon
5064 DO j = j_start,j_end
5065 DO i = i_start,i_end
5066 govrth(i,j) = g/th_phy(i,1,j)
5070 hflux: SELECT CASE( config_flags%isfflx )
5071 CASE (0,2) ! with fixed surface heat flux given in the namelist
5072 heat_flux = config_flags%tke_heat_flux ! constant heat flux value
5073 DO j = j_start, j_end
5074 DO i = i_start, i_end
5075 cpm = cp * (1. + 0.8 * moist(i,kts,j,P_QV))
5076 hfx(i,j)= heat_flux*cpm*rho(i,1,j)
5080 CASE (1) ! use surface heat flux computed from surface routine
5081 DO j = j_start, j_end
5082 DO i = i_start, i_end
5083 cpm = cp * (1. + 0.8 * moist(i,kts,j,P_QV))
5084 heat_flux = hfx(i,j)/cpm/rho(i,1,j)
5089 CALL wrf_error_fatal( 'isfflx value invalid for diff_opt=2' )
5092 DO j = j_start, j_end
5093 DO i = i_start, i_end
5095 hpbl(i,j) = zq(i,1,j)
5096 pblflg(i,j) = .true.
5097 sfcflg(i,j) = .true.
5099 cpm = cp * (1. + 0.8*moist(i,1,j,P_QV))
5100 sflux(i,j) = (hfx(i,j)/cpm)/rho(i,1,j)
5102 IF(br(i,j).GT.0.0) sfcflg(i,j) = .false.
5106 ! get pbl height begin based on theta method
5107 DO j = j_start, j_end
5108 DO i = i_start, i_end
5109 thsfc = thv(i,kts,j) + 0.5
5111 IF(thv(i,k,j).GT.thsfc) THEN
5112 hpbl(i,j)=za(i,k-1,j)+(thsfc-thv(i,k-1,j))/(max(0.01,thv(i,k,j)-thv(i,k-1,j)))*(za(i,k,j)-za(i,k-1,j))
5121 DO j = j_start, j_end
5122 DO i = i_start, i_end
5123 IF(hpbl(i,j).LT.zq(i,2,j)) kpbl(i,j) = 1
5124 IF(kpbl(i,j).LE.1) pblflg(i,j) = .false.
5128 DO j = j_start, j_end
5129 DO i = i_start, i_end
5131 bfx0 = max(sflux(i,j),0.)
5132 wstar3(i,j) = (govrth(i,j)*bfx0*hpbl(i,j))
5133 wstar(i,j) = (wstar3(i,j))**h1
5139 wscale(i,j) = (ust3 + phifac*karman*wstar3(i,j)*0.5)**h1
5140 wscale(i,j) = MIN(wscale(i,j), ust(i,j)*16.0)
5141 wscale(i,j) = MAX(wscale(i,j), ust(i,j)/5.0 )
5145 DO j = j_start, j_end
5146 DO i = i_start, i_end
5147 delxy = sqrt(dx/msftx(i,j)*dy/msfty(i,j))
5148 pu1=pu(delxy,hpbl(i,j))
5149 IF((sfcflg(i,j).and.sflux(i,j).GT.0.0) .and. (hpbl(i,j) .GT. 0)) THEN
5150 !nonlocal momentum flux based on Brown and Grant (1997)
5151 brint = -sm*ust(i,j)*ust(i,j)*wstar3(i,j)/(hpbl(i,j)*wscale(i,j)**4)
5152 gamu(i,j) = pu1 * brint*u_phy(i,1,j)/wspd(i,j)
5153 gamv(i,j) = pu1 * brint*v_phy(i,1,j)/wspd(i,j)
5155 pblflg(i,j) = .false.
5160 DO j = j_start,j_end
5161 DO i = i_start,i_end
5164 wm3 = wstar3(i,j) + 5. * ust(i,j)**3
5166 bfxpbl(i,j) = -0.15*thv(i,1,j)/g*wm3/hpbl(i,j)
5167 dthv(i,j) = max(thv(i,k+1,j) - thv(i,k,j),tmin)
5168 dtheta = max(th_phy(i,k+1,j) - th_phy(i,k,j),tmin)
5170 du = u_phy(i,k+1,j)-u_phy(i,k,j)
5171 dv = v_phy(i,k+1,j)-v_phy(i,k,j)
5173 we(i,j) = max(bfxpbl(i,j)/dthv(i,j),-sqrt(wm2(i,j)))
5174 hfxpbl(i,j) = we(i,j)*dtheta
5175 delb = govrth(i,j)*dthv(i,j)
5177 deltaoh(i,j) = d1*hpbl(i,j) + d2*wm2(i,j)/delb
5178 deltaoh(i,j) = max(ezfac*deltaoh(i,j),hpbl(i,j)-za(i,k,j)-1.)
5179 deltaoh(i,j) = min(deltaoh(i,j), hpbl(i,j))
5181 if ((du .ne. 0) .or. (dv .ne. 0)) then
5182 rigs(i,j) = govrth(i,j)*dthv(i,j)*deltaoh(i,j)/(du**2.+dv**2.)
5183 rigs(i,j) = max(min(rigs(i,j), rigsmax),rimin)
5187 enlfrac2(i,j) = max(min(wm3/wstar3(i,j)/(1.+cpent/rigs(i,j)),entfmax),entfmin)
5188 enlfrac2(i,j) = enlfrac2(i,j)*enlfrac
5193 DO j = j_start, j_end
5195 DO i = i_start, i_end
5197 entfacmf(i,k,j) = sqrt(((zq(i,k+1,j) - hpbl(i,j))/deltaoh(i,j))**2.)
5203 DO j = j_start, j_end
5204 DO i = i_start, i_end
5205 IF (pblflg(i,j)) THEN
5206 deltaoh(i,j) = deltaoh(i,j)/hpbl(i,j)
5207 delxy = sqrt(dx/msftx(i,j)*dy/msfty(i,j))
5208 mlfrac = mltop-deltaoh(i,j)
5209 ezfrac = mltop+deltaoh(i,j)
5210 zfacmf(i,1,j) = min(max((zq(i,2,j)/hpbl(i,j)),zfmin),1.)
5211 sfcfracn = max(sfcfracn1,zfacmf(i,1,j))
5213 sflux0 = (a11+a12*sfcfracn)*sflux(i,j)
5214 snlflux0 = nlfrac*sflux0
5215 amf1 = snlflux0/sfcfracn
5216 amf2 = -snlflux0/(mlfrac-sfcfracn)
5218 amf3 = snlflux0*enlfrac2(i,j)/deltaoh(i,j)
5220 hfxpbl(i,j) = amf3+bmf3
5221 pth1 = pthnl(delxy,hpbl(i,j))
5222 hfxpbl(i,j) = hfxpbl(i,j)*pth1
5225 zfacmf(i,k,j) = max((zq(i,k+1,j)/hpbl(i,j)),zfmin)
5226 IF(k.LT.kpbl(i,j)) THEN
5227 IF(zfacmf(i,k,j).LE.sfcfracn) THEN
5228 nlflux(i,k,j) = amf1*zfacmf(i,k,j)
5229 ELSE IF (zfacmf(i,k,j).LE.mlfrac) THEN
5230 nlflux(i,k,j) = amf2*zfacmf(i,k,j)+bmf2
5232 nlflux(i,k,j) = nlflux(i,k,j) + hfxpbl(i,j)*exp(-entfacmf(i,k,j))
5233 nlflux(i,k,j) = nlflux(i,k,j)*pth1
5239 END SUBROUTINE nonlocal_flux
5241 !==============================================================================
5242 !==============================================================================
5244 ! partial function for nonlocal heat flux
5248 REAL,PARAMETER :: pmin = 0.0,pmax = 1.0
5249 REAL,PARAMETER :: a1 = 1.000, a2 = 0.936, a3 = -1.110, &
5250 a4 = 1.000, a5 = 0.312, a6 = 0.329, a7 = 0.243
5251 REAL,PARAMETER :: b1 = 2.0, b2 = 0.875
5252 real :: d,h,doh,num,den
5256 num = a1*(doh)**b1 + a2*(doh)**b2+a3
5257 den = a4*(doh)**b1 + a5*(doh)**b2+a6
5258 pthnl = a7*num/den + (1. - a7)
5263 pthnl = max(pthnl,pmin)
5264 pthnl = min(pthnl,pmax)
5266 IF(d.LE.100.) THEN ! assume dx<=100m as LES
5273 ! partial function for local heat flux
5277 REAL,PARAMETER :: pmin = 0.0,pmax = 1.0
5278 REAL,PARAMETER :: a1 = 1.000, a2 = 0.870, a3 = -0.913, &
5279 a4 = 1.000, a5 = 0.153, a6 = 0.278, a7 = 0.280
5280 REAL,PARAMETER :: b1 = 2.0, b2 = 0.5
5281 REAL :: d,h,doh,num,den
5285 num = a1*(doh)**b1 + a2*(doh)**b2+a3
5286 den = a4*(doh)**b1 + a5*(doh)**b2+a6
5287 pthl = a7*num/den+(1. - a7)
5291 pthl = max(pthl,pmin)
5292 pthl = min(pthl,pmax)
5294 IF(d.LE.100.) THEN ! assume dx<=100m as LES
5301 ! partial function for momentum flux
5305 REAL,PARAMETER :: pmin = 0.0,pmax = 1.0
5306 REAL,PARAMETER :: a1 = 1.0, a2 = 0.070, a3 = 1.0, a4 = 0.142, a5 = 0.071
5307 REAL,PARAMETER :: b1 = 2.0, b2 = 0.6666667
5308 REAL :: d,h,doh,num,den
5312 num = a1*(doh)**b1 + a2*(doh)**b2
5313 den = a3*(doh)**b1 + a4*(doh)**b2+a5
5321 IF(d.LE.100.) THEN ! assume dx<=100 as LES
5328 !=======================================================================
5329 !=======================================================================
5331 SUBROUTINE cal_titau_11_22_33( config_flags, titau, &
5334 is_ext, ie_ext, js_ext, je_ext, &
5335 ids, ide, jds, jde, kds, kde, &
5336 ims, ime, jms, jme, kms, kme, &
5337 its, ite, jts, jte, kts, kte )
5339 ! History: Sep 2003 Changes by George Bryan and Jason Knievel, NCAR
5340 ! Oct 2001 Converted to mass core by Bill Skamarock, NCAR
5341 ! Aug 2000 Original code by Shu-Hua Chen, UC-Davis
5343 ! Purpose: This routine calculates stress terms (taus) for use in
5344 ! the calculation of production of TKE by sheared wind
5346 ! References: Klemp and Wilhelmson (JAS 1978)
5347 ! Deardorff (B-L Meteor 1980)
5348 ! Chen and Dudhia (NCAR WRF physics report 2000)
5352 !-----------------------------------------------------------------------
5353 ! Begin declarations.
5357 TYPE( grid_config_rec_type ), INTENT( IN ) &
5360 INTEGER, INTENT( IN ) &
5361 :: ids, ide, jds, jde, kds, kde, &
5362 ims, ime, jms, jme, kms, kme, &
5363 its, ite, jts, jte, kts, kte
5365 INTEGER, INTENT( IN ) &
5366 :: is_ext, ie_ext, js_ext, je_ext
5368 REAL, DIMENSION( its-1:ite+1, kts:kte, jts-1:jte+1 ), INTENT( INOUT ) &
5371 REAL, DIMENSION( ims:ime, kms:kme, jms:jme ), INTENT( IN ) &
5372 :: defor, xkx, tke, rho
5374 REAL, DIMENSION( ims:ime, kms:kme, jms:jme ), INTENT( INOUT ) &
5380 :: i, j, k, ktf, i_start, i_end, j_start, j_end
5383 !-----------------------------------------------------------------------
5385 ktf = MIN( kte, kde-1 )
5392 IF ( config_flags%open_xs .OR. config_flags%specified .OR. &
5393 config_flags%nested) i_start = MAX( ids+1, its )
5394 IF ( config_flags%open_xe .OR. config_flags%specified .OR. &
5395 config_flags%nested) i_end = MIN( ide-1, ite )
5396 IF ( config_flags%open_ys .OR. config_flags%specified .OR. &
5397 config_flags%nested) j_start = MAX( jds+1, jts )
5398 IF ( config_flags%open_ye .OR. config_flags%specified .OR. &
5399 config_flags%nested) j_end = MIN( jde-1, jte )
5400 IF ( config_flags%periodic_x ) i_start = its
5401 IF ( config_flags%periodic_x ) i_end = ite
5403 i_start = i_start - is_ext
5404 i_end = i_end + ie_ext
5405 j_start = j_start - js_ext
5406 j_end = j_end + je_ext
5408 IF ( config_flags%sfs_opt .GT. 0 ) THEN ! USE NBA MODEL SFS STRESSES
5410 DO j = j_start, j_end
5412 DO i = i_start, i_end
5414 titau(i,k,j) = rho(i,k,j) * mtau(i,k,j)
5422 IF ( config_flags%m_opt .EQ. 1 ) THEN ! ASSIGN STRESS TO MTAU FOR OUTPUT
5424 DO j = j_start, j_end
5426 DO i = i_start, i_end
5428 titau(i,k,j) = - rho(i,k,j) * xkx(i,k,j) * defor(i,k,j)
5429 mtau(i,k,j) = - xkx(i,k,j) * defor(i,k,j)
5435 ELSE !NO STRESS OUTPUT
5437 DO j = j_start, j_end
5439 DO i = i_start, i_end
5441 titau(i,k,j) = - rho(i,k,j) * xkx(i,k,j) * defor(i,k,j)
5451 END SUBROUTINE cal_titau_11_22_33
5453 !=======================================================================
5454 !=======================================================================
5456 SUBROUTINE cal_titau_12_21( config_flags, titau, &
5459 is_ext, ie_ext, js_ext, je_ext, &
5460 ids, ide, jds, jde, kds, kde, &
5461 ims, ime, jms, jme, kms, kme, &
5462 its, ite, jts, jte, kts, kte )
5464 ! History: Sep 2003 Modifications by George Bryan and Jason Knievel, NCAR
5465 ! Oct 2001 Converted to mass core by Bill Skamarock, NCAR
5466 ! Aug 2000 Original code by Shu-Hua Chen, UC-Davis
5468 ! Pusrpose This routine calculates the stress terms (taus) for use in
5469 ! the calculation of production of TKE by sheared wind
5471 ! References: Klemp and Wilhelmson (JAS 1978)
5472 ! Deardorff (B-L Meteor 1980)
5473 ! Chen and Dudhia (NCAR WRF physics report 2000)
5477 !-----------------------------------------------------------------------
5478 ! Begin declarations.
5482 TYPE( grid_config_rec_type), INTENT( IN ) &
5485 INTEGER, INTENT( IN ) &
5486 :: ids, ide, jds, jde, kds, kde, &
5487 ims, ime, jms, jme, kms, kme, &
5488 its, ite, jts, jte, kts, kte
5490 INTEGER, INTENT( IN ) &
5491 :: is_ext, ie_ext, js_ext, je_ext
5493 REAL, DIMENSION( its-1:ite+1, kts:kte, jts-1:jte+1 ), INTENT( INOUT ) &
5496 REAL, DIMENSION( ims:ime, kms:kme, jms:jme ), INTENT( IN ) &
5499 REAL, DIMENSION( ims:ime, kms:kme, jms:jme ), INTENT( INOUT ) &
5505 :: i, j, k, ktf, i_start, i_end, j_start, j_end
5507 REAL, DIMENSION( its-1:ite+1, kts:kte, jts-1:jte+1 ) &
5511 !-----------------------------------------------------------------------
5513 ktf = MIN( kte, kde-1 )
5515 ! Needs one more point in the x and y directions.
5522 IF ( config_flags%open_xs .OR. config_flags%specified .OR. &
5523 config_flags%nested ) i_start = MAX( ids+1, its )
5524 IF ( config_flags%open_xe .OR. config_flags%specified .OR. &
5525 config_flags%nested ) i_end = MIN( ide-1, ite )
5526 IF ( config_flags%open_ys .OR. config_flags%specified .OR. &
5527 config_flags%nested ) j_start = MAX( jds+1, jts )
5528 IF ( config_flags%open_ye .OR. config_flags%specified .OR. &
5529 config_flags%nested ) j_end = MIN( jde-1, jte )
5530 IF ( config_flags%periodic_x ) i_start = its
5531 IF ( config_flags%periodic_x ) i_end = ite
5533 i_start = i_start - is_ext
5534 i_end = i_end + ie_ext
5535 j_start = j_start - js_ext
5536 j_end = j_end + je_ext
5538 DO j = j_start, j_end
5540 DO i = i_start, i_end
5541 rhoavg(i,k,j) = 0.25 * ( rho(i-1,k,j ) + rho(i,k,j ) + &
5542 rho(i-1,k,j-1) + rho(i,k,j-1) )
5543 xkxavg(i,k,j) = rhoavg(i,k,j) * &
5544 0.25 * ( xkx(i-1,k,j ) + xkx(i,k,j ) + &
5545 xkx(i-1,k,j-1) + xkx(i,k,j-1) )
5550 ! titau12 or titau21
5552 IF ( config_flags%sfs_opt .GT. 0 ) THEN ! USE NBA MODEL SFS STRESSES
5554 DO j = j_start, j_end
5556 DO i = i_start, i_end
5558 titau(i,k,j) = rhoavg(i,k,j) * mtau(i,k,j)
5566 IF ( config_flags%m_opt .EQ. 1 ) THEN ! ASSIGN STRESS TO MTAU FOR OUTPUT
5568 DO j = j_start, j_end
5570 DO i = i_start, i_end
5571 titau(i,k,j) = - xkxavg(i,k,j) * defor(i,k,j)
5572 mtau(i,k,j) = - xkxavg(i,k,j) * defor(i,k,j) / rhoavg(i,k,j)
5578 ELSE ! NO STRESS OUTPUT
5580 DO j = j_start, j_end
5582 DO i = i_start, i_end
5584 titau(i,k,j) = - xkxavg(i,k,j) * defor(i,k,j)
5594 END SUBROUTINE cal_titau_12_21
5596 !=======================================================================
5598 SUBROUTINE cal_titau_13_31( config_flags, titau, &
5601 xkx, fnm, fnp, rho, &
5602 is_ext, ie_ext, js_ext, je_ext, &
5603 ids, ide, jds, jde, kds, kde, &
5604 ims, ime, jms, jme, kms, kme, &
5605 its, ite, jts, jte, kts, kte )
5607 ! History: Sep 2003 Modifications by George Bryan and Jason Knievel, NCAR
5608 ! Oct 2001 Converted to mass core by Bill Skamarock, NCAR
5609 ! Aug 2000 Original code by Shu-Hua Chen, UC-Davis
5611 ! Purpose: This routine calculates the stress terms (taus) for use in
5612 ! the calculation of production of TKE by sheared wind
5614 ! References: Klemp and Wilhelmson (JAS 1978)
5615 ! Deardorff (B-L Meteor 1980)
5616 ! Chen and Dudhia (NCAR WRF physics report 2000)
5620 !-----------------------------------------------------------------------
5621 ! Begin declarations.
5625 TYPE( grid_config_rec_type), INTENT( IN ) &
5628 INTEGER, INTENT( IN ) &
5629 :: ids, ide, jds, jde, kds, kde, &
5630 ims, ime, jms, jme, kms, kme, &
5631 its, ite, jts, jte, kts, kte
5633 INTEGER, INTENT( IN ) &
5634 :: is_ext, ie_ext, js_ext, je_ext
5636 REAL, DIMENSION( kms:kme ), INTENT( IN ) &
5639 REAL, DIMENSION( its-1:ite+1, kts:kte, jts-1:jte+1 ), INTENT( INOUT ) &
5642 REAL, DIMENSION( ims:ime, kms:kme, jms:jme), INTENT( IN ) &
5645 REAL, DIMENSION( ims:ime, kms:kme, jms:jme ), INTENT( INOUT ) &
5651 :: i, j, k, ktf, i_start, i_end, j_start, j_end
5653 REAL, DIMENSION( its-1:ite+1, kts:kte, jts-1:jte+1 ) &
5657 !-----------------------------------------------------------------------
5659 ktf = MIN( kte, kde-1 )
5661 ! Find ide-1 and jde-1 for averaging to p point.
5666 j_end = MIN( jte, jde-1 )
5668 IF ( config_flags%open_xs .OR. config_flags%specified .OR. &
5669 config_flags%nested) i_start = MAX( ids+1, its )
5670 IF ( config_flags%open_xe .OR. config_flags%specified .OR. &
5671 config_flags%nested) i_end = MIN( ide-1, ite )
5672 IF ( config_flags%open_ys .OR. config_flags%specified .OR. &
5673 config_flags%nested) j_start = MAX( jds+1, jts )
5674 IF ( config_flags%open_ye .OR. config_flags%specified .OR. &
5675 config_flags%nested) j_end = MIN( jde-2, jte )
5676 IF ( config_flags%periodic_x ) i_start = its
5677 IF ( config_flags%periodic_x ) i_end = ite
5679 i_start = i_start - is_ext
5680 i_end = i_end + ie_ext
5681 j_start = j_start - js_ext
5682 j_end = j_end + je_ext
5684 DO j = j_start, j_end
5686 DO i = i_start, i_end
5687 rhoavg(i,k,j) = 0.5 * ( fnm(k) * ( rho(i-1,k ,j) + rho(i,k ,j) ) + &
5688 fnp(k) * ( rho(i-1,k-1,j) + rho(i,k-1,j) ) )
5689 xkxavg(i,k,j) = rhoavg(i,k,j) * &
5690 0.5 * ( fnm(k) * ( xkx(i,k ,j) + xkx(i-1,k ,j) ) + &
5691 fnp(k) * ( xkx(i,k-1,j) + xkx(i-1,k-1,j) ) )
5696 IF ( config_flags%sfs_opt .GT. 0 ) THEN ! USE NBA MODEL SFS STRESSES
5698 DO j = j_start, j_end
5700 DO i = i_start, i_end
5701 titau(i,k,j) = rhoavg(i,k,j) * mtau(i,k,j)
5708 IF ( config_flags%m_opt .EQ. 1 ) THEN ! ASSIGN STRESS TO MTAU FOR OUTPUT
5710 DO j = j_start, j_end
5712 DO i = i_start, i_end
5714 titau(i,k,j) = - xkxavg(i,k,j) * defor(i,k,j)
5715 mtau(i,k,j) = - xkxavg(i,k,j) * defor(i,k,j) / rhoavg(i,k,j)
5721 ELSE ! NO STRESS OUTPUT
5723 DO j = j_start, j_end
5725 DO i = i_start, i_end
5727 titau(i,k,j) = - xkxavg(i,k,j) * defor(i,k,j)
5737 DO j = j_start, j_end
5738 DO i = i_start, i_end
5739 titau(i,kts ,j) = 0.0
5740 titau(i,ktf+1,j) = 0.0
5744 END SUBROUTINE cal_titau_13_31
5746 !=======================================================================
5747 !=======================================================================
5749 SUBROUTINE cal_titau_23_32( config_flags, titau, defor, &
5751 xkx, fnm, fnp, rho, &
5752 is_ext, ie_ext, js_ext, je_ext, &
5753 ids, ide, jds, jde, kds, kde, &
5754 ims, ime, jms, jme, kms, kme, &
5755 its, ite, jts, jte, kts, kte )
5757 ! History: Sep 2003 Changes by George Bryan and Jason Knievel, NCAR
5758 ! Oct 2001 Converted to mass core by Bill Skamarock, NCAR
5759 ! Aug 2000 Original code by Shu-Hua Chen, UC-Davis
5761 ! Purpose: This routine calculates stress terms (taus) for use in
5762 ! the calculation of production of TKE by sheared wind
5764 ! References: Klemp and Wilhelmson (JAS 1978)
5765 ! Deardorff (B-L Meteor 1980)
5766 ! Chen and Dudhia (NCAR WRF physics report 2000)
5770 !-----------------------------------------------------------------------
5771 ! Begin declarations.
5775 TYPE( grid_config_rec_type ), INTENT( IN ) &
5778 INTEGER, INTENT( IN ) &
5779 :: ids, ide, jds, jde, kds, kde, &
5780 ims, ime, jms, jme, kms, kme, &
5781 its, ite, jts, jte, kts, kte
5783 INTEGER, INTENT( IN ) &
5784 :: is_ext,ie_ext,js_ext,je_ext
5786 REAL, DIMENSION( kms:kme ), INTENT( IN ) &
5789 REAL, DIMENSION( its-1:ite+1, kts:kte, jts-1:jte+1 ), INTENT( INOUT ) &
5792 REAL, DIMENSION( ims:ime, kms:kme, jms:jme ), INTENT( IN ) &
5795 REAL, DIMENSION( ims:ime, kms:kme, jms:jme ), INTENT( INOUT ) &
5801 :: i, j, k, ktf, i_start, i_end, j_start, j_end
5803 REAL, DIMENSION( its-1:ite+1, kts:kte, jts-1:jte+1 ) &
5807 !-----------------------------------------------------------------------
5809 ktf = MIN( kte, kde-1 )
5811 ! Find ide-1 and jde-1 for averaging to p point.
5814 i_end = MIN( ite, ide-1 )
5818 IF ( config_flags%open_xs .OR. config_flags%specified .OR. &
5819 config_flags%nested) i_start = MAX( ids+1, its )
5820 IF ( config_flags%open_xe .OR. config_flags%specified .OR. &
5821 config_flags%nested) i_end = MIN( ide-2, ite )
5822 IF ( config_flags%open_ys .OR. config_flags%specified .OR. &
5823 config_flags%nested) j_start = MAX( jds+1, jts )
5824 IF ( config_flags%open_ye .OR. config_flags%specified .OR. &
5825 config_flags%nested) j_end = MIN( jde-1, jte )
5826 IF ( config_flags%periodic_x ) i_start = its
5827 IF ( config_flags%periodic_x ) i_end = MIN( ite, ide-1 )
5829 i_start = i_start - is_ext
5830 i_end = i_end + ie_ext
5831 j_start = j_start - js_ext
5832 j_end = j_end + je_ext
5834 DO j = j_start, j_end
5836 DO i = i_start, i_end
5837 rhoavg(i,k,j) = 0.5 * ( fnm(k) * ( rho(i,k ,j) + rho(i,k ,j-1) ) + &
5838 fnp(k) * ( rho(i,k-1,j) + rho(i,k-1,j-1) ) )
5839 xkxavg(i,k,j) = rhoavg(i,k,j) * &
5840 0.5 * ( fnm(k) * ( xkx(i,k ,j) + xkx(i,k ,j-1) ) + &
5841 fnp(k) * ( xkx(i,k-1,j) + xkx(i,k-1,j-1) ) )
5846 IF ( config_flags%sfs_opt .GT. 0 ) THEN ! USE NBA MODEL SFS STRESSES
5848 DO j = j_start, j_end
5850 DO i = i_start, i_end
5852 titau(i,k,j) = rhoavg(i,k,j) * mtau(i,k,j)
5860 IF ( config_flags%m_opt .EQ. 1 ) THEN ! ASSIGN STRESS TO MTAU FOR OUTPUT
5862 DO j = j_start, j_end
5864 DO i = i_start, i_end
5866 titau(i,k,j) = - xkxavg(i,k,j) * defor(i,k,j)
5867 mtau(i,k,j) = - xkxavg(i,k,j) * defor(i,k,j) / rhoavg(i,k,j)
5873 ELSE ! NO STRESS OUTPUT
5875 DO j = j_start, j_end
5877 DO i = i_start, i_end
5879 titau(i,k,j) = - xkxavg(i,k,j) * defor(i,k,j)
5889 DO j = j_start, j_end
5890 DO i = i_start, i_end
5891 titau(i,kts ,j) = 0.0
5892 titau(i,ktf+1,j) = 0.0
5896 END SUBROUTINE cal_titau_23_32
5898 !=======================================================================
5899 !=======================================================================
5901 SUBROUTINE phy_bc ( config_flags,div,defor11,defor22,defor33, &
5902 defor12,defor13,defor23,xkmh,xkmv,xkhh,xkhv,tke,rho, &
5906 gamu, gamv, xkmv_meso, & ! XZ
5907 ids, ide, jds, jde, kds, kde, &
5908 ims, ime, jms, jme, kms, kme, &
5909 ips, ipe, jps, jpe, kps, kpe, &
5910 its, ite, jts, jte, kts, kte )
5912 !------------------------------------------------------------------------------
5913 ! Begin declarations.
5917 TYPE(grid_config_rec_type), INTENT(IN ) :: config_flags
5919 INTEGER , INTENT(IN ) :: ids, ide, jds, jde, kds, kde, &
5920 ims, ime, jms, jme, kms, kme, &
5921 ips, ipe, jps, jpe, kps, kpe, &
5922 its, ite, jts, jte, kts, kte
5924 REAL , DIMENSION( ims:ime, kms:kme, jms:jme), INTENT(INOUT) ::RUBLTEN, &
5945 REAL , DIMENSION( ims:ime, jms:jme), INTENT(INOUT) :: gamu, gamv
5947 !-----------------------------------------------------------------------
5949 IF(config_flags%bl_pbl_physics .GT. 0) THEN
5951 CALL set_physical_bc3d( RUBLTEN , 't', config_flags, &
5952 ids, ide, jds, jde, kds, kde, &
5953 ims, ime, jms, jme, kms, kme, &
5954 ips, ipe, jps, jpe, kps, kpe, &
5955 its, ite, jts, jte, kts, kte )
5957 CALL set_physical_bc3d( RVBLTEN , 't', config_flags, &
5958 ids, ide, jds, jde, kds, kde, &
5959 ims, ime, jms, jme, kms, kme, &
5960 ips, ipe, jps, jpe, kps, kpe, &
5961 its, ite, jts, jte, kts, kte )
5965 IF(config_flags%cu_physics .GT. 0) THEN
5967 CALL set_physical_bc3d( RUCUTEN , 't', config_flags, &
5968 ids, ide, jds, jde, kds, kde, &
5969 ims, ime, jms, jme, kms, kme, &
5970 ips, ipe, jps, jpe, kps, kpe, &
5971 its, ite, jts, jte, kts, kte )
5973 CALL set_physical_bc3d( RVCUTEN , 't', config_flags, &
5974 ids, ide, jds, jde, kds, kde, &
5975 ims, ime, jms, jme, kms, kme, &
5976 ips, ipe, jps, jpe, kps, kpe, &
5977 its, ite, jts, jte, kts, kte )
5981 IF(config_flags%shcu_physics .GT. 0) THEN
5983 CALL set_physical_bc3d( RUSHTEN , 't', config_flags, &
5984 ids, ide, jds, jde, kds, kde, &
5985 ims, ime, jms, jme, kms, kme, &
5986 ips, ipe, jps, jpe, kps, kpe, &
5987 its, ite, jts, jte, kts, kte )
5989 CALL set_physical_bc3d( RVSHTEN , 't', config_flags, &
5990 ids, ide, jds, jde, kds, kde, &
5991 ims, ime, jms, jme, kms, kme, &
5992 ips, ipe, jps, jpe, kps, kpe, &
5993 its, ite, jts, jte, kts, kte )
5997 ! move out of the conditional, below; horiz coeffs needed for
5998 ! all diff_opt cases. JM
6000 CALL set_physical_bc3d( xkmh , 't', config_flags, &
6001 ids, ide, jds, jde, kds, kde, &
6002 ims, ime, jms, jme, kms, kme, &
6003 ips, ipe, jps, jpe, kps, kpe, &
6004 its, ite, jts, jte, kts, kte )
6006 CALL set_physical_bc3d( xkhh , 't', config_flags, &
6007 ids, ide, jds, jde, kds, kde, &
6008 ims, ime, jms, jme, kms, kme, &
6009 ips, ipe, jps, jpe, kps, kpe, &
6010 its, ite, jts, jte, kts, kte )
6012 IF(config_flags%diff_opt .eq. 2) THEN
6014 CALL set_physical_bc3d( xkmv , 't', config_flags, &
6015 ids, ide, jds, jde, kds, kde, &
6016 ims, ime, jms, jme, kms, kme, &
6017 ips, ipe, jps, jpe, kps, kpe, &
6018 its, ite, jts, jte, kts, kte )
6020 CALL set_physical_bc3d( xkhv , 't', config_flags, &
6021 ids, ide, jds, jde, kds, kde, &
6022 ims, ime, jms, jme, kms, kme, &
6023 ips, ipe, jps, jpe, kps, kpe, &
6024 its, ite, jts, jte, kts, kte )
6026 CALL set_physical_bc3d( div , 't', config_flags, &
6027 ids, ide, jds, jde, kds, kde, &
6028 ims, ime, jms, jme, kms, kme, &
6029 ips, ipe, jps, jpe, kps, kpe, &
6030 its, ite, jts, jte, kts, kte )
6032 CALL set_physical_bc3d( defor11 , 't', config_flags, &
6033 ids, ide, jds, jde, kds, kde, &
6034 ims, ime, jms, jme, kms, kme, &
6035 ips, ipe, jps, jpe, kps, kpe, &
6036 its, ite, jts, jte, kts, kte )
6038 CALL set_physical_bc3d( defor22 , 't', config_flags, &
6039 ids, ide, jds, jde, kds, kde, &
6040 ims, ime, jms, jme, kms, kme, &
6041 ips, ipe, jps, jpe, kps, kpe, &
6042 its, ite, jts, jte, kts, kte )
6044 CALL set_physical_bc3d( defor33 , 't', config_flags, &
6045 ids, ide, jds, jde, kds, kde, &
6046 ims, ime, jms, jme, kms, kme, &
6047 ips, ipe, jps, jpe, kps, kpe, &
6048 its, ite, jts, jte, kts, kte )
6050 CALL set_physical_bc3d( defor12 , 'd', config_flags, &
6051 ids, ide, jds, jde, kds, kde, &
6052 ims, ime, jms, jme, kms, kme, &
6053 ips, ipe, jps, jpe, kps, kpe, &
6054 its, ite, jts, jte, kts, kte )
6056 CALL set_physical_bc3d( defor13 , 'e', config_flags, &
6057 ids, ide, jds, jde, kds, kde, &
6058 ims, ime, jms, jme, kms, kme, &
6059 ips, ipe, jps, jpe, kps, kpe, &
6060 its, ite, jts, jte, kts, kte )
6062 CALL set_physical_bc3d( defor23 , 'f', config_flags, &
6063 ids, ide, jds, jde, kds, kde, &
6064 ims, ime, jms, jme, kms, kme, &
6065 ips, ipe, jps, jpe, kps, kpe, &
6066 its, ite, jts, jte, kts, kte )
6068 CALL set_physical_bc3d( rho , 't', config_flags, &
6069 ids, ide, jds, jde, kds, kde, &
6070 ims, ime, jms, jme, kms, kme, &
6071 ips, ipe, jps, jpe, kps, kpe, &
6072 its, ite, jts, jte, kts, kte )
6076 IF(config_flags%diff_opt .eq. 2 .and. config_flags%km_opt .eq. 5) THEN
6077 CALL set_physical_bc2d( gamu, 't', config_flags, &
6078 ids, ide, jds, jde, &
6079 ims, ime, jms, jme, &
6080 ips, ipe, jps, jpe, &
6081 its, ite, jts, jte )
6082 CALL set_physical_bc2d( gamv, 't', config_flags, &
6083 ids, ide, jds, jde, &
6084 ims, ime, jms, jme, &
6085 ips, ipe, jps, jpe, &
6086 its, ite, jts, jte )
6087 CALL set_physical_bc3d( xkmv_meso,'t', config_flags, &
6088 ids, ide, jds, jde, kds, kde, &
6089 ims, ime, jms, jme, kms, kme, &
6090 ips, ipe, jps, jpe, kps, kpe, &
6091 its, ite, jts, jte, kts, kte )
6094 END SUBROUTINE phy_bc
6096 !=======================================================================
6097 !=======================================================================
6099 SUBROUTINE tke_rhs( tendency, BN2, config_flags, &
6100 defor11, defor22, defor33, &
6101 defor12, defor13, defor23, &
6102 u, v, w, div, tke, mu, c1, c2, &
6103 theta, p, p8w, t8w, z, fnm, fnp, &
6104 cf1, cf2, cf3, msftx, msfty, &
6106 rdx, rdy, dx, dy, dt, zx, zy, &
6107 rdz, rdzw, dn, dnw, isotropic, &
6108 hfx, qfx, qv, ust, rho, &
6109 l_diss, nlflux, hpbl, dlk, & !XZ
6110 ids, ide, jds, jde, kds, kde, &
6111 ims, ime, jms, jme, kms, kme, &
6112 its, ite, jts, jte, kts, kte )
6114 !-----------------------------------------------------------------------
6115 ! Begin declarations.
6119 TYPE( grid_config_rec_type ), INTENT( IN ) &
6122 INTEGER, INTENT( IN ) &
6123 :: ids, ide, jds, jde, kds, kde, &
6124 ims, ime, jms, jme, kms, kme, &
6125 its, ite, jts, jte, kts, kte
6127 INTEGER, INTENT( IN ) :: isotropic
6128 REAL, INTENT( IN ) &
6129 :: cf1, cf2, cf3, dt, rdx, rdy, dx, dy
6131 REAL, DIMENSION( kms:kme ), INTENT( IN ) &
6132 :: fnm, fnp, dnw, dn
6134 REAL, DIMENSION( ims:ime, jms:jme ), INTENT( IN ) &
6137 REAL, DIMENSION( ims:ime, kms:kme, jms:jme ), INTENT( INOUT ) &
6140 REAL, DIMENSION( ims:ime, kms:kme, jms:jme ), INTENT( IN ) &
6141 :: defor11, defor22, defor33, defor12, defor13, defor23, &
6142 div, BN2, tke, xkmh, xkmv, xkhv, zx, zy, u, v, w, theta, &
6143 p, p8w, t8w, z, rdz, rdzw
6145 REAL, DIMENSION( ims:ime, jms:jme ), INTENT( IN ) &
6147 REAL, DIMENSION( kms:kme) , INTENT( IN ) &
6150 REAL, DIMENSION ( ims:ime, jms:jme ), INTENT( IN ) &
6152 REAL, DIMENSION ( ims:ime, kms:kme, jms:jme ), INTENT ( IN ) &
6155 REAL, DIMENSION( ims:ime, kms:kme, jms:jme ), INTENT( OUT ) &
6158 REAL, DIMENSION( ims:ime, kms:kme, jms:jme ), INTENT ( IN ) &
6161 REAL, DIMENSION( ims:ime, jms:jme ), INTENT(INOUT) &
6167 :: i, j, k, ktf, i_start, i_end, j_start, j_end
6170 !-----------------------------------------------------------------------
6172 CALL tke_shear( tendency, config_flags, &
6173 defor11, defor22, defor33, &
6174 defor12, defor13, defor23, &
6175 u, v, w, tke, ust, mu, &
6177 cf1, cf2, cf3, msftx, msfty, &
6179 rdx, rdy, zx, zy, rdz, rdzw, dnw, dn, &
6180 ids, ide, jds, jde, kds, kde, &
6181 ims, ime, jms, jme, kms, kme, &
6182 its, ite, jts, jte, kts, kte )
6184 CALL tke_buoyancy( tendency, config_flags, mu, &
6186 tke, xkhv, BN2, theta, dt, &
6187 hfx, qfx, qv, rho, &
6189 ids, ide, jds, jde, kds, kde, &
6190 ims, ime, jms, jme, kms, kme, &
6191 its, ite, jts, jte, kts, kte )
6193 CALL tke_dissip( tendency, config_flags, mu, c1, c2, &
6194 tke, bn2, theta, p8w, t8w, z, &
6195 dx, dy,rdz, rdzw, isotropic, &
6197 hpbl, dlk, l_diss, & ! XZ
6198 ids, ide, jds, jde, kds, kde, &
6199 ims, ime, jms, jme, kms, kme, &
6200 its, ite, jts, jte, kts, kte )
6202 ! Set a lower limit on TKE.
6204 ktf = MIN( kte, kde-1 )
6206 i_end = MIN( ite, ide-1 )
6208 j_end = MIN( jte, jde-1 )
6210 IF ( config_flags%open_xs .or. config_flags%specified .or. &
6211 config_flags%nested) i_start = MAX(ids+1,its)
6212 IF ( config_flags%open_xe .or. config_flags%specified .or. &
6213 config_flags%nested) i_end = MIN(ide-2,ite)
6214 IF ( config_flags%open_ys .or. config_flags%specified .or. &
6215 config_flags%nested) j_start = MAX(jds+1,jts)
6216 IF ( config_flags%open_ye .or. config_flags%specified .or. &
6217 config_flags%nested) j_end = MIN(jde-2,jte)
6218 IF ( config_flags%periodic_x ) i_start = its
6219 IF ( config_flags%periodic_x ) i_end = MIN( ite, ide-1 )
6221 DO j = j_start, j_end
6223 DO i = i_start, i_end
6224 tendency(i,k,j) = max( tendency(i,k,j), -(c1(k)*mu(i,j)+c2(k)) * max( 0.0 , tke(i,k,j) ) / dt )
6229 END SUBROUTINE tke_rhs
6231 !=======================================================================
6232 !=======================================================================
6234 SUBROUTINE tke_buoyancy( tendency, config_flags, mu, &
6236 tke, xkhv, BN2, theta, dt, &
6237 hfx, qfx, qv, rho, &
6239 ids, ide, jds, jde, kds, kde, &
6240 ims, ime, jms, jme, kms, kme, &
6241 its, ite, jts, jte, kts, kte )
6243 !-----------------------------------------------------------------------
6244 ! Begin declarations.
6248 TYPE( grid_config_rec_type ), INTENT( IN ) &
6251 INTEGER, INTENT( IN ) &
6252 :: ids, ide, jds, jde, kds, kde, &
6253 ims, ime, jms, jme, kms, kme, &
6254 its, ite, jts, jte, kts, kte
6256 REAL, INTENT( IN ) &
6259 REAL, DIMENSION( ims:ime, kms:kme, jms:jme ), INTENT( INOUT ) &
6262 REAL, DIMENSION( ims:ime, kms:kme, jms:jme ), INTENT( IN ) &
6263 :: xkhv, tke, BN2, theta
6265 REAL, DIMENSION( ims:ime, jms:jme ), INTENT( IN ) &
6267 REAL, DIMENSION( kms:kme) , INTENT( IN ) &
6270 REAL, DIMENSION( ims:ime, kms:kme, jms:jme ), INTENT ( IN ) &
6273 REAL, DIMENSION(ims:ime, jms:jme ), INTENT ( IN ) :: hfx, qfx
6275 REAL, DIMENSION( ims:ime, kms:kme, jms:jme ), INTENT ( IN ) &
6284 :: i_start, i_end, j_start, j_end
6286 REAL :: heat_flux, heat_flux0
6291 !-----------------------------------------------------------------------
6293 !-----------------------------------------------------------------------
6294 ! Add to the TKE tendency the term that accounts for production of TKE
6295 ! due to buoyant motions.
6297 ktf = MIN( kte, kde-1 )
6299 i_end = MIN( ite, ide-1 )
6301 j_end = MIN( jte, jde-1 )
6303 IF ( config_flags%open_xs .OR. config_flags%specified .OR. &
6304 config_flags%nested ) i_start = MAX( ids+1, its )
6305 IF ( config_flags%open_xe .OR. config_flags%specified .OR. &
6306 config_flags%nested ) i_end = MIN( ide-2, ite )
6307 IF ( config_flags%open_ys .OR. config_flags%specified .OR. &
6308 config_flags%nested ) j_start = MAX( jds+1, jts )
6309 IF ( config_flags%open_ye .OR. config_flags%specified .OR. &
6310 config_flags%nested ) j_end = MIN( jde-2, jte )
6311 IF ( config_flags%periodic_x ) i_start = its
6312 IF ( config_flags%periodic_x ) i_end = MIN( ite, ide-1 )
6315 IF(config_flags%km_opt.eq.5) THEN
6316 DO j = j_start, j_end
6318 DO i = i_start, i_end
6319 tendency(i,k,j) = tendency(i,k,j) - (c1(k)*mu(i,j)+c2(k)) * (xkhv(i,k,j) * BN2(i,k,j) &
6320 - 0.5*(nlflux(i,k+1,j)+nlflux(i,k,j))*g/theta(i,k,j))
6325 DO j = j_start, j_end
6327 DO i = i_start, i_end
6328 tendency(i,k,j) = tendency(i,k,j) - (c1(k)*mu(i,j)+c2(k)) * xkhv(i,k,j) * BN2(i,k,j)
6334 ! MARTA: change in the computation of the tke's tendency at the surface.
6335 ! the buoyancy flux is the average of the surface heat flux (0.06) and the
6336 ! flux at the first w level
6340 hflux: SELECT CASE( config_flags%isfflx )
6341 CASE (0,2) ! with fixed surface heat flux given in the namelist
6342 heat_flux0 = config_flags%tke_heat_flux ! constant heat flux value
6343 ! set in namelist.input
6347 DO j = j_start, j_end
6348 DO i = i_start, i_end
6349 heat_flux = heat_flux0
6350 tendency(i,k,j)= tendency(i,k,j) - &
6351 (c1(k)*mu(i,j)+c2(k))*((xkhv(i,k,j)*BN2(i,k,j))- (g/theta(i,k,j))*heat_flux)/2.
6356 CASE (1) ! use surface heat flux computed from surface routine
6359 DO j = j_start, j_end
6360 DO i = i_start, i_end
6361 cpm = cp * (1. + 0.8*qv(i,k,j))
6362 heat_flux = (hfx(i,j)/cpm)/rho(i,k,j)
6364 tendency(i,k,j)= tendency(i,k,j) - &
6365 (c1(k)*mu(i,j)+c2(k))*((xkhv(i,k,j)*BN2(i,k,j))- (g/theta(i,k,j))*heat_flux)/2.
6370 CALL wrf_error_fatal( 'isfflx value invalid for diff_opt=2' )
6373 ! end of MARTA/WCS change
6375 ! The tendency array now includes production of TKE from buoyant
6377 !-----------------------------------------------------------------------
6379 END SUBROUTINE tke_buoyancy
6381 !=======================================================================
6382 !=======================================================================
6384 SUBROUTINE tke_dissip( tendency, config_flags, &
6386 tke, bn2, theta, p8w, t8w, z, &
6387 dx, dy, rdz, rdzw, isotropic, &
6389 hpbl, dlk, l_diss, & !XZ
6390 ids, ide, jds, jde, kds, kde, &
6391 ims, ime, jms, jme, kms, kme, &
6392 its, ite, jts, jte, kts, kte )
6394 ! History: Sep 2003 Changes by George Bryan and Jason Knievel, NCAR
6395 ! Oct 2001 Converted to mass core by Bill Skamarock, NCAR
6396 ! Aug 2000 Original code by Shu-Hua Chen, UC-Davis
6398 ! Purpose: This routine calculates dissipation of turbulent kinetic
6401 ! References: Klemp and Wilhelmson (JAS 1978)
6402 ! Deardorff (B-L Meteor 1980)
6403 ! Chen and Dudhia (NCAR WRF physics report 2000)
6405 !-----------------------------------------------------------------------
6406 ! Begin declarations.
6410 TYPE( grid_config_rec_type ), INTENT( IN ) &
6413 INTEGER, INTENT( IN ) &
6414 :: ids, ide, jds, jde, kds, kde, &
6415 ims, ime, jms, jme, kms, kme, &
6416 its, ite, jts, jte, kts, kte
6418 INTEGER, INTENT( IN ) :: isotropic
6419 REAL, INTENT( IN ) &
6422 REAL, DIMENSION( ims:ime, kms:kme, jms:jme ), INTENT( INOUT ) &
6425 REAL, DIMENSION( ims:ime, kms:kme, jms:jme ), INTENT( IN ) &
6426 :: tke, bn2, theta, p8w, t8w, z, rdz, rdzw
6428 REAL, DIMENSION( ims:ime, jms:jme ), INTENT( IN ) &
6430 REAL, DIMENSION( kms:kme) , INTENT( IN ) &
6433 REAL, DIMENSION( ims:ime, jms:jme ), INTENT( IN ) &
6436 REAL, DIMENSION( ims:ime, jms:jme ), INTENT(INOUT) &
6439 REAL, DIMENSION( ims:ime, kms:kme, jms:jme ), INTENT(IN) &
6442 REAL, DIMENSION( ims:ime, kms:kme, jms:jme ), INTENT(OUT) &
6447 REAL, DIMENSION( its:ite, kts:kte, jts:jte ) &
6450 REAL, DIMENSION( its:ite, kts:kte, jts:jte ) &
6453 REAL, DIMENSION( its:ite ) &
6457 :: i, j, k, ktf, i_start, i_end, j_start, j_end
6460 :: disp_len, deltas, coefc, tmpdz, len_s, thetasfc, &
6461 thetatop, len_0, tketmp, tmp, ce1, ce2, c_k
6464 :: delxy,coefc_m,pu1
6467 !-----------------------------------------------------------------------
6468 c_k = config_flags%c_k
6470 ce1 = ( c_k / 0.10 ) * 0.19
6471 ce2 = max( 0.0 , 0.93 - ce1 )
6473 ktf = MIN( kte, kde-1 )
6475 i_end = MIN(ite,ide-1)
6477 j_end = MIN(jte,jde-1)
6479 IF ( config_flags%open_xs .OR. config_flags%specified .OR. &
6480 config_flags%nested) i_start = MAX( ids+1, its )
6481 IF ( config_flags%open_xe .OR. config_flags%specified .OR. &
6482 config_flags%nested) i_end = MIN( ide-2, ite )
6483 IF ( config_flags%open_ys .OR. config_flags%specified .OR. &
6484 config_flags%nested) j_start = MAX( jds+1, jts )
6485 IF ( config_flags%open_ye .OR. config_flags%specified .OR. &
6486 config_flags%nested) j_end = MIN( jde-2, jte )
6487 IF ( config_flags%periodic_x ) i_start = its
6488 IF ( config_flags%periodic_x ) i_end = MIN( ite, ide-1 )
6490 CALL calc_l_scale( config_flags, tke, BN2, l_scale, &
6491 i_start, i_end, ktf, j_start, j_end, &
6492 dx, dy, rdzw, msftx, msfty, &
6493 ids, ide, jds, jde, kds, kde, &
6494 ims, ime, jms, jme, kms, kme, &
6495 its, ite, jts, jte, kts, kte )
6496 DO j = j_start, j_end
6498 DO i = i_start, i_end
6499 deltas = ( dx/msftx(i,j) * dy/msfty(i,j) / rdzw(i,k,j) )**0.33333333
6500 tketmp = MAX( tke(i,k,j), 1.0e-6 )
6502 ! Apply Deardorff's (1980) "wall effect" at the bottom of the domain.
6503 ! For LES with fine grid, no need for this wall effect!
6505 IF ( k .eq. kts .or. k .eq. ktf ) then
6508 coefc = ce1 + ce2 * l_scale(i,k,j) / deltas
6512 IF(config_flags%km_opt.eq.5) THEN
6513 delxy = SQRT(dx/msftx(i,j)*dy/msfty(i,j))
6514 pu1 = pu(delxy,hpbl(i,j))
6516 l_diss(i,k,j) = (1.0-pu1)*l_scale(i,k,j)/coefc + pu1*dlk(i,k,j)/coefc_m
6518 tendency(i,k,j) = tendency(i,k,j) - &
6519 (c1(k)*mu(i,j)+c2(k)) * coefc * tketmp**1.5 / l_scale(i,k,j)
6524 END SUBROUTINE tke_dissip
6526 !=======================================================================
6527 !=======================================================================
6529 SUBROUTINE tke_shear( tendency, config_flags, &
6530 defor11, defor22, defor33, &
6531 defor12, defor13, defor23, &
6532 u, v, w, tke, ust, mu, c1, c2, &
6534 cf1, cf2, cf3, msftx, msfty, &
6536 rdx, rdy, zx, zy, rdz, rdzw, dn, dnw, &
6537 ids, ide, jds, jde, kds, kde, &
6538 ims, ime, jms, jme, kms, kme, &
6539 its, ite, jts, jte, kts, kte )
6541 ! History: Sep 2003 Rewritten by George Bryan and Jason Knievel,
6543 ! Oct 2001 Converted to mass core by Bill Skamarock, NCAR
6544 ! Aug 2000 Original code by Shu-Hua Chen, UC-Davis
6546 ! Purpose: This routine calculates the production of turbulent
6547 ! kinetic energy by stresses due to sheared wind.
6549 ! References: Klemp and Wilhelmson (JAS 1978)
6550 ! Deardorff (B-L Meteor 1980)
6551 ! Chen and Dudhia (NCAR WRF physics report 2000)
6555 ! avg temporary working array
6559 ! defor11 deformation term ( du/dx + du/dx )
6560 ! defor12 deformation term ( dv/dx + du/dy ); same as defor21
6561 ! defor13 deformation term ( dw/dx + du/dz ); same as defor31
6562 ! defor22 deformation term ( dv/dy + dv/dy )
6563 ! defor23 deformation term ( dw/dy + dv/dz ); same as defor32
6564 ! defor33 deformation term ( dw/dz + dw/dz )
6565 ! div 3-d divergence
6575 ! titau tau (stress tensor) with a tilde, indicating division by
6576 ! a map-scale factor and the fraction of the total modeled
6577 ! atmosphere beneath a given altitude (titau = tau/m/zeta)
6578 ! tke turbulent kinetic energy
6580 !-----------------------------------------------------------------------
6581 ! Begin declarations.
6585 TYPE( grid_config_rec_type ), INTENT( IN ) &
6588 INTEGER, INTENT( IN ) &
6589 :: ids, ide, jds, jde, kds, kde, &
6590 ims, ime, jms, jme, kms, kme, &
6591 its, ite, jts, jte, kts, kte
6593 REAL, INTENT( IN ) &
6594 :: cf1, cf2, cf3, rdx, rdy
6596 REAL, DIMENSION( kms:kme ), INTENT( IN ) &
6597 :: fnm, fnp, dn, dnw
6599 REAL, DIMENSION( ims:ime, jms:jme ), INTENT( IN ) &
6602 REAL, DIMENSION( ims:ime, kms:kme, jms:jme ), INTENT( INOUT ) &
6605 REAL, DIMENSION( ims:ime, kms:kme, jms:jme ), INTENT( IN ) &
6606 :: defor11, defor22, defor33, defor12, defor13, defor23, &
6607 tke, xkmh, xkmv, zx, zy, u, v, w, rdz, rdzw
6609 REAL, DIMENSION( ims:ime, jms:jme ), INTENT( IN ) &
6611 REAL, DIMENSION( kms:kme) , INTENT( IN ) &
6614 REAL, DIMENSION( ims:ime, jms:jme ), INTENT( IN ) &
6620 :: i, j, k, ktf, ktes1, ktes2, &
6621 i_start, i_end, j_start, j_end, &
6622 is_ext, ie_ext, js_ext, je_ext
6627 REAL, DIMENSION( its-1:ite+1, kts:kte, jts-1:jte+1 ) &
6630 REAL, DIMENSION( its:ite, kts:kte, jts:jte ) &
6631 :: titau12, tmp1, zxavg, zyavg
6633 REAL :: absU, cd0, Cd
6636 !-----------------------------------------------------------------------
6638 ktf = MIN( kte, kde-1 )
6643 i_end = MIN( ite, ide-1 )
6645 j_end = MIN( jte, jde-1 )
6647 IF ( config_flags%open_xs .OR. config_flags%specified .OR. &
6648 config_flags%nested ) i_start = MAX( ids+1, its )
6649 IF ( config_flags%open_xe .OR. config_flags%specified .OR. &
6650 config_flags%nested ) i_end = MIN( ide-2, ite )
6651 IF ( config_flags%open_ys .OR. config_flags%specified .OR. &
6652 config_flags%nested ) j_start = MAX( jds+1, jts )
6653 IF ( config_flags%open_ye .OR. config_flags%specified .OR. &
6654 config_flags%nested ) j_end = MIN( jde-2, jte )
6655 IF ( config_flags%periodic_x ) i_start = its
6656 IF ( config_flags%periodic_x ) i_end = MIN( ite, ide-1 )
6658 DO j = j_start, j_end
6660 DO i = i_start, i_end
6661 zxavg(i,k,j) = 0.25 * ( zx(i,k ,j) + zx(i+1,k ,j) + &
6662 zx(i,k+1,j) + zx(i+1,k+1,j) )
6663 zyavg(i,k,j) = 0.25 * ( zy(i,k ,j) + zy(i,k ,j+1) + &
6664 zy(i,k+1,j) + zy(i,k+1,j+1) )
6669 ! Begin calculating production of turbulence due to shear. The approach
6670 ! is to add together contributions from six terms, each of which is the
6671 ! square of a deformation that is then multiplied by an exchange
6672 ! coefficiant. The same exchange coefficient is assumed for horizontal
6673 ! and vertical coefficients for some of the terms (the vertical value is
6678 DO j = j_start, j_end
6680 DO i = i_start, i_end
6681 tendency(i,k,j) = tendency(i,k,j) + 0.5 * &
6682 (c1(k)*mu(i,j)+c2(k)) * xkmh(i,k,j) * ( ( defor11(i,k,j) )**2 )
6689 DO j = j_start, j_end
6691 DO i = i_start, i_end
6692 tendency(i,k,j) = tendency(i,k,j) + 0.5 * &
6693 (c1(k)*mu(i,j)+c2(k)) * xkmh(i,k,j) * ( ( defor22(i,k,j) )**2 )
6700 DO j = j_start, j_end
6702 DO i = i_start, i_end
6703 tendency(i,k,j) = tendency(i,k,j) + 0.5 * &
6704 (c1(k)*mu(i,j)+c2(k)) * xkmv(i,k,j) * ( ( defor33(i,k,j) )**2 )
6711 DO j = j_start, j_end
6713 DO i = i_start, i_end
6714 avg(i,k,j) = 0.25 * &
6715 ( ( defor12(i ,k,j)**2 ) + ( defor12(i ,k,j+1)**2 ) + &
6716 ( defor12(i+1,k,j)**2 ) + ( defor12(i+1,k,j+1)**2 ) )
6721 DO j = j_start, j_end
6723 DO i = i_start, i_end
6724 tendency(i,k,j) = tendency(i,k,j) + (c1(k)*mu(i,j)+c2(k)) * xkmh(i,k,j) * avg(i,k,j)
6731 DO j = j_start, j_end
6733 DO i = i_start, i_end+1
6734 tmp2(i,k,j) = defor13(i,k,j)
6739 DO j = j_start, j_end
6740 DO i = i_start, i_end+1
6741 tmp2(i,kts ,j) = 0.0
6742 tmp2(i,ktf+1,j) = 0.0
6746 DO j = j_start, j_end
6748 DO i = i_start, i_end
6749 avg(i,k,j) = 0.25 * &
6750 ( ( tmp2(i ,k+1,j)**2 ) + ( tmp2(i ,k,j)**2 ) + &
6751 ( tmp2(i+1,k+1,j)**2 ) + ( tmp2(i+1,k,j)**2 ) )
6756 DO j = j_start, j_end
6758 DO i = i_start, i_end
6759 tendency(i,k,j) = tendency(i,k,j) + (c1(k)*mu(i,j)+c2(k)) * xkmv(i,k,j) * avg(i,k,j)
6764 !MARTA: add the drag at the surface; WCS 040331
6767 uflux: SELECT CASE( config_flags%isfflx )
6768 CASE (0) ! Assume cd a constant, specified in namelist
6770 cd0 = config_flags%tke_drag_coefficient ! drag coefficient set
6772 DO j = j_start, j_end
6773 DO i = i_start, i_end
6775 absU=0.5*sqrt((u(i,k,j)+u(i+1,k,j))**2+(v(i,k,j)+v(i,k,j+1))**2)
6777 tendency(i,k,j) = tendency(i,k,j) + &
6778 (c1(k)*mu(i,j)+c2(k))*( (u(i,k,j)+u(i+1,k,j))*0.5* &
6779 Cd*absU*(defor13(i,kts+1,j)+defor13(i+1,kts+1,j))*0.5 )
6784 CASE (1,2) ! ustar computed from surface routine
6786 DO j = j_start, j_end
6787 DO i = i_start, i_end
6789 absU=0.5*sqrt((u(i,k,j)+u(i+1,k,j))**2+(v(i,k,j)+v(i,k,j+1))**2)+epsilon
6790 Cd = (ust(i,j)**2)/(absU**2)
6791 tendency(i,k,j) = tendency(i,k,j) + &
6792 (c1(k)*mu(i,j)+c2(k))*( (u(i,k,j)+u(i+1,k,j))*0.5* &
6793 Cd*absU*(defor13(i,kts+1,j)+defor13(i+1,kts+1,j))*0.5 )
6799 CALL wrf_error_fatal( 'isfflx value invalid for diff_opt=2' )
6801 ! end of MARTA/WCS change
6805 DO j = j_start, j_end+1
6807 DO i = i_start, i_end
6808 tmp2(i,k,j) = defor23(i,k,j)
6813 DO j = j_start, j_end+1
6814 DO i = i_start, i_end
6815 tmp2(i,kts, j) = 0.0
6816 tmp2(i,ktf+1,j) = 0.0
6820 DO j = j_start, j_end
6822 DO i = i_start, i_end
6823 avg(i,k,j) = 0.25 * &
6824 ( ( tmp2(i,k+1,j )**2 ) + ( tmp2(i,k,j )**2) + &
6825 ( tmp2(i,k+1,j+1)**2 ) + ( tmp2(i,k,j+1)**2) )
6830 DO j = j_start, j_end
6832 DO i = i_start, i_end
6833 tendency(i,k,j) = tendency(i,k,j) + (c1(k)*mu(i,j)+c2(k)) * xkmv(i,k,j) * avg(i,k,j)
6838 !MARTA: add the drag at the surface; WCS 040331
6841 vflux: SELECT CASE( config_flags%isfflx )
6842 CASE (0) ! Assume cd a constant, specified in namelist
6844 cd0 = config_flags%tke_drag_coefficient ! constant drag coefficient
6845 ! set in namelist.input
6846 DO j = j_start, j_end
6847 DO i = i_start, i_end
6849 absU=0.5*sqrt((u(i,k,j)+u(i+1,k,j))**2+(v(i,k,j)+v(i,k,j+1))**2)
6851 tendency(i,k,j) = tendency(i,k,j) + &
6852 (c1(k)*mu(i,j)+c2(k))*( (v(i,k,j)+v(i,k,j+1))*0.5* &
6853 Cd*absU*(defor23(i,kts+1,j)+defor23(i,kts+1,j+1))*0.5 )
6858 CASE (1,2) ! ustar computed from surface routine
6860 DO j = j_start, j_end
6861 DO i = i_start, i_end
6863 absU=0.5*sqrt((u(i,k,j)+u(i+1,k,j))**2+(v(i,k,j)+v(i,k,j+1))**2)+epsilon
6864 Cd = (ust(i,j)**2)/(absU**2)
6865 tendency(i,k,j) = tendency(i,k,j) + &
6866 (c1(k)*mu(i,j)+c2(k))*( (v(i,k,j)+v(i,k,j+1))*0.5* &
6867 Cd*absU*(defor23(i,kts+1,j)+defor23(i,kts+1,j+1))*0.5 )
6873 CALL wrf_error_fatal( 'isfflx value invalid for diff_opt=2' )
6875 ! end of MARTA/WCS change
6877 END SUBROUTINE tke_shear
6879 !=======================================================================
6880 !=======================================================================
6882 SUBROUTINE compute_diff_metrics( config_flags, ph, phb, z, rdz, rdzw, &
6884 ids, ide, jds, jde, kds, kde, &
6885 ims, ime, jms, jme, kms, kme, &
6886 its, ite, jts, jte, kts, kte )
6888 !-----------------------------------------------------------------------
6889 ! Begin declarations.
6893 TYPE( grid_config_rec_type ), INTENT( IN ) &
6896 INTEGER, INTENT( IN ) &
6897 :: ids, ide, jds, jde, kds, kde, &
6898 ims, ime, jms, jme, kms, kme, &
6899 its, ite, jts, jte, kts, kte
6901 REAL, DIMENSION( ims:ime, kms:kme, jms:jme ), INTENT( IN ) &
6904 REAL, DIMENSION( ims:ime, kms:kme, jms:jme ), INTENT( OUT ) &
6905 :: rdz, rdzw, zx, zy, z
6908 REAL, INTENT( IN ) &
6913 REAL, DIMENSION( its-1:ite, kts:kte, jts-1:jte ) &
6917 :: i, j, k, i_start, i_end, j_start, j_end, ktf
6920 !-----------------------------------------------------------------------
6922 ktf = MIN( kte, kde-1 )
6924 ! Bug fix, WCS, 22 april 2002.
6926 ! We need rdzw in halo for average to u and v points.
6931 ! Begin with dz computations.
6933 DO j = j_start, j_end
6935 IF ( ( j_start >= jts ) .AND. ( j_end <= MIN( jte, jde-1 ) ) ) THEN
6940 i_end = MIN( ite, ide-1 )
6943 ! Compute z at w points for rdz and rdzw computations. We'll switch z
6944 ! to z at p points before returning
6948 ! Bug fix, WCS, 22 april 2002
6950 DO i = i_start, i_end
6951 z_at_w(i,k,j) = ( ph(i,k,j) + phb(i,k,j) ) / g
6956 DO i = i_start, i_end
6957 rdzw(i,k,j) = 1.0 / ( z_at_w(i,k+1,j) - z_at_w(i,k,j) )
6962 DO i = i_start, i_end
6963 rdz(i,k,j) = 2.0 / ( z_at_w(i,k+1,j) - z_at_w(i,k-1,j) )
6967 ! Bug fix, WCS, 22 april 2002; added the following code
6969 DO i = i_start, i_end
6970 rdz(i,1,j) = 2./(z_at_w(i,2,j)-z_at_w(i,1,j))
6977 ! Now compute zx and zy; we'll assume that the halo for ph and phb is
6981 i_end = MIN( ite, ide-1 )
6983 j_end = MIN( jte, jde-1 )
6985 DO j = j_start, j_end
6987 DO i = MAX( ids+1, its ), i_end
6988 zx(i,k,j) = rdx * ( phb(i,k,j) - phb(i-1,k,j) ) / g
6993 DO j = j_start, j_end
6995 DO i = MAX( ids+1, its ), i_end
6996 zx(i,k,j) = zx(i,k,j) + rdx * ( ph(i,k,j) - ph(i-1,k,j) ) / g
7001 DO j = MAX( jds+1, jts ), j_end
7003 DO i = i_start, i_end
7004 zy(i,k,j) = rdy * ( phb(i,k,j) - phb(i,k,j-1) ) / g
7009 DO j = MAX( jds+1, jts ), j_end
7011 DO i = i_start, i_end
7012 zy(i,k,j) = zy(i,k,j) + rdy * ( ph(i,k,j) - ph(i,k,j-1) ) / g
7017 ! Some b.c. on zx and zy.
7019 IF ( .NOT. config_flags%periodic_x ) THEN
7021 IF ( ite == ide ) THEN
7022 DO j = j_start, j_end
7029 IF ( its == ids ) THEN
7030 DO j = j_start, j_end
7039 IF ( ite == ide ) THEN
7042 zx(ide,k,j) = rdx * ( phb(ide,k,j) - phb(ide-1,k,j) ) / g
7046 DO j = j_start, j_end
7048 zx(ide,k,j) = zx(ide,k,j) + rdx * ( ph(ide,k,j) - ph(ide-1,k,j) ) / g
7053 IF ( its == ids ) THEN
7054 DO j = j_start, j_end
7056 zx(ids,k,j) = rdx * ( phb(ids,k,j) - phb(ids-1,k,j) ) / g
7062 zx(ids,k,j) = zx(ids,k,j) + rdx * ( ph(ids,k,j) - ph(ids-1,k,j) ) / g
7069 IF ( .NOT. config_flags%periodic_y ) THEN
7071 IF ( jte == jde ) THEN
7073 DO i =i_start, i_end
7079 IF ( jts == jds ) THEN
7081 DO i =i_start, i_end
7089 IF ( jte == jde ) THEN
7091 DO i =i_start, i_end
7092 zy(i,k,jde) = rdy * ( phb(i,k,jde) - phb(i,k,jde-1) ) / g
7097 DO i =i_start, i_end
7098 zy(i,k,jde) = zy(i,k,jde) + rdy * ( ph(i,k,jde) - ph(i,k,jde-1) ) / g
7103 IF ( jts == jds ) THEN
7105 DO i =i_start, i_end
7106 zy(i,k,jds) = rdy * ( phb(i,k,jds) - phb(i,k,jds-1) ) / g
7111 DO i =i_start, i_end
7112 zy(i,k,jds) = zy(i,k,jds) + rdy * ( ph(i,k,jds) - ph(i,k,jds-1) ) / g
7119 ! Calculate z at p points.
7121 DO j = j_start, j_end
7123 DO i = i_start, i_end
7125 ( ph(i,k,j) + phb(i,k,j) + ph(i,k+1,j) + phb(i,k+1,j) ) / g
7130 END SUBROUTINE compute_diff_metrics
7132 SUBROUTINE cal_helicity ( config_flags, u, v, w, uh, up_heli_max,&
7137 rdx, rdy, dn, dnw, rdz, rdzw, &
7138 fnm, fnp, cf1, cf2, cf3, zx, zy, &
7139 ids, ide, jds, jde, kds, kde, &
7140 ims, ime, jms, jme, kms, kme, &
7141 its, ite, jts, jte, kts, kte )
7143 ! History: Apr 2007 Coded by Scott Dembek, USRA, basically using
7144 ! code segments stolen from the deformation
7145 ! and divergence subroutine.
7148 ! Purpose: This routine calculates updraft helicity.
7149 ! w ( dv/dx - du/dy )
7153 !-----------------------------------------------------------------------
7154 ! Begin declarations.
7158 TYPE( grid_config_rec_type ), INTENT( IN ) &
7161 INTEGER, INTENT( IN ) &
7162 :: ids, ide, jds, jde, kds, kde, &
7163 ims, ime, jms, jme, kms, kme, &
7164 its, ite, jts, jte, kts, kte
7166 REAL, INTENT( IN ) &
7167 :: rdx, rdy, cf1, cf2, cf3
7169 REAL, DIMENSION( kms:kme ), INTENT( IN ) &
7170 :: fnm, fnp, dn, dnw
7172 REAL, DIMENSION( ims:ime , jms:jme ), INTENT( IN ) &
7173 :: msfux, msfuy, msfvx, msfvy, ht
7175 REAL, DIMENSION( ims:ime, kms:kme, jms:jme ), INTENT( IN ) &
7176 :: u, v, w, ph, phb, zx, zy, rdz, rdzw
7178 REAL, DIMENSION( ims:ime, jms:jme ), INTENT( INOUT ) &
7184 :: i, j, k, ktf, ktes1, ktes2, i_start, i_end, j_start, j_end
7187 :: tmp, tmpzx, tmpzy, tmpzeta_z, cft1, cft2
7192 REAL, DIMENSION( its-3:ite+2, jts-3:jte+2 ) :: mm
7194 REAL, DIMENSION( its-3:ite+2, kts:kte, jts-3:jte+2 ) &
7195 :: tmp1, hat, hatavg
7197 REAL, DIMENSION( its-3:ite+2, kts:kte, jts-3:jte+2 ) &
7200 LOGICAL, DIMENSION( its-3:ite+2, jts-3:jte+2 ) &
7204 !-----------------------------------------------------------------------
7206 !=======================================================================
7207 ! In the following section, calculate the vertical component of the
7208 ! relative vorticity as the first two terms of the updraft helicity
7209 ! w ( dv/dx - du/dy ).
7214 cft2 = - 0.5 * dnw(ktes1) / dn(ktes1)
7217 ktf = MIN( kte, kde-1 )
7219 !=======================================================================
7226 IF ( config_flags%open_xs .OR. config_flags%specified .OR. &
7227 config_flags%nested) i_start = MAX( ids+1, its-2 )
7228 IF ( config_flags%open_xe .OR. config_flags%specified .OR. &
7229 config_flags%nested) i_end = MIN( ide-1, ite+2 )
7230 IF ( config_flags%open_ys .OR. config_flags%specified .OR. &
7231 config_flags%nested) j_start = MAX( jds+1, jts-2 )
7232 IF ( config_flags%open_ye .OR. config_flags%specified .OR. &
7233 config_flags%nested) j_end = MIN( jde-1, jte+2 )
7234 IF ( config_flags%periodic_x ) i_start = its
7235 IF ( config_flags%periodic_x ) i_end = ite
7238 !-----------------------------------------------------------------------
7241 ! First, calculate an average mapscale factor.
7243 ! Comments 10-MAR-05
7244 ! du/dy => need u map scale factor in x (which is defined at u points)
7245 ! averaged over j and j-1
7246 ! dv/dx => need v map scale factor in y (which is defined at v points)
7247 ! averaged over i and i-1
7249 DO j = j_start, j_end
7250 DO i = i_start, i_end
7251 mm(i,j) = 0.25 * ( msfux(i,j-1) + msfux(i,j) ) * ( msfvy(i-1,j) + msfvy(i,j) )
7255 ! Apply a coordinate transformation to meridional velocity, v.
7257 DO j = j_start, j_end
7259 DO i = i_start-1, i_end
7260 hat(i,k,j) = v(i,k,j) / msfvy(i,j)
7265 ! Account for the slope in x of eta surfaces.
7267 DO j = j_start, j_end
7269 DO i = i_start, i_end
7270 hatavg(i,k,j) = 0.5 * ( &
7271 fnm(k) * ( hat(i-1,k ,j) + hat(i,k ,j) ) + &
7272 fnp(k) * ( hat(i-1,k-1,j) + hat(i,k-1,j) ) )
7277 ! Extrapolate to top and bottom of domain (to w levels).
7279 DO j = j_start, j_end
7280 DO i = i_start, i_end
7281 hatavg(i,1,j) = 0.5 * ( &
7282 cf1 * hat(i-1,1,j) + &
7283 cf2 * hat(i-1,2,j) + &
7284 cf3 * hat(i-1,3,j) + &
7285 cf1 * hat(i ,1,j) + &
7286 cf2 * hat(i ,2,j) + &
7288 hatavg(i,kte,j) = 0.5 * ( &
7289 cft1 * ( hat(i,ktes1,j) + hat(i-1,ktes1,j) ) + &
7290 cft2 * ( hat(i,ktes2,j) + hat(i-1,ktes2,j) ) )
7294 ! Fixes to set_physical_bc2/3d have made any check for polar B.C.'s
7295 ! unnecessary in this place. zx, rdzw, and hatavg are all defined
7296 ! in places they need to be and the values at the poles are replications
7297 ! of the values one grid point in, so the averaging over j and j-1 works
7298 ! to act as just using the value at j or j-1 (with out extra code).
7300 ! tmpzx = averaged value of dpsi/dx (=zx) on vorticity grid
7301 ! tmp1 = partial dpsi/dx * partial dv^/dpsi
7302 DO j = j_start, j_end
7304 DO i = i_start, i_end
7306 zx(i,k ,j-1) + zx(i,k ,j) + &
7307 zx(i,k+1,j-1) + zx(i,k+1,j) )
7308 tmp1(i,k,j) = ( hatavg(i,k+1,j) - hatavg(i,k,j) ) * &
7309 0.25 * tmpzx * ( rdzw(i,k,j) + rdzw(i,k,j-1) + &
7310 rdzw(i-1,k,j-1) + rdzw(i-1,k,j) )
7315 ! End calculation of dv/dx.
7316 !-----------------------------------------------------------------------
7318 !-----------------------------------------------------------------------
7319 ! Calculate the first parenthetical term of
7320 ! updraft helicity = w ( dv/dx - du/dy ) at vorticity points.
7322 DO j = j_start, j_end
7324 DO i = i_start, i_end
7325 rvort(i,k,j) = mm(i,j) * ( &
7326 rdx * ( hat(i,k,j) - hat(i-1,k,j) ) - tmp1(i,k,j) )
7331 ! End calculation of the first parenthetical term of updraft helicity.
7332 !-----------------------------------------------------------------------
7334 !-----------------------------------------------------------------------
7337 ! Apply a coordinate transformation to zonal velocity, u.
7339 DO j =j_start-1, j_end
7341 DO i =i_start, i_end
7342 ! Fixes to set_physical_bc2/3d for polar boundary conditions
7343 ! remove issues with loop over j
7344 hat(i,k,j) = u(i,k,j) / msfux(i,j)
7349 ! Average in y and z.
7354 hatavg(i,k,j) = 0.5 * ( &
7355 fnm(k) * ( hat(i,k ,j-1) + hat(i,k ,j) ) + &
7356 fnp(k) * ( hat(i,k-1,j-1) + hat(i,k-1,j) ) )
7361 ! Extrapolate to top and bottom of domain (to w levels).
7363 DO j = j_start, j_end
7364 DO i = i_start, i_end
7365 hatavg(i,1,j) = 0.5 * ( &
7366 cf1 * hat(i,1,j-1) + &
7367 cf2 * hat(i,2,j-1) + &
7368 cf3 * hat(i,3,j-1) + &
7369 cf1 * hat(i,1,j ) + &
7370 cf2 * hat(i,2,j ) + &
7372 hatavg(i,kte,j) = 0.5 * ( &
7373 cft1 * ( hat(i,ktes1,j-1) + hat(i,ktes1,j) ) + &
7374 cft2 * ( hat(i,ktes2,j-1) + hat(i,ktes2,j) ) )
7378 ! tmpzy = averaged value of dpsi/dy (=zy) on vorticity grid
7379 ! tmp1 = partial dpsi/dy * partial du^/dpsi
7380 DO j = j_start, j_end
7382 DO i = i_start, i_end
7384 zy(i-1,k ,j) + zy(i,k ,j) + &
7385 zy(i-1,k+1,j) + zy(i,k+1,j) )
7386 tmp1(i,k,j) = ( hatavg(i,k+1,j) - hatavg(i,k,j) ) * &
7387 0.25 * tmpzy * ( rdzw(i,k,j) + rdzw(i-1,k,j) + &
7388 rdzw(i-1,k,j-1) + rdzw(i,k,j-1) )
7393 ! End calculation of du/dy.
7394 !----------------------------------------------------------------------
7396 !-----------------------------------------------------------------------
7397 ! Add (subtract, actually) the second parenthetical term to
7398 ! updraft helicity = w ( dv/dx - du/dy ) at vorticity points.
7400 DO j = j_start, j_end
7402 DO i = i_start, i_end
7403 rvort(i,k,j) = rvort(i,k,j) - &
7405 rdy * ( hat(i,k,j) - hat(i,k,j-1) ) - tmp1(i,k,j) )
7410 ! End addition of the second parenthetical term to updraft helicity.
7411 !-----------------------------------------------------------------------
7413 !-----------------------------------------------------------------------
7414 ! Update the boundary for the vorticity.
7416 IF ( .NOT. config_flags%periodic_x .AND. i_start .EQ. ids+1 ) THEN
7419 rvort(ids,k,j) = rvort(ids+1,k,j)
7424 IF ( .NOT. config_flags%periodic_y .AND. j_start .EQ. jds+1) THEN
7427 rvort(i,k,jds) = rvort(i,k,jds+1)
7432 IF ( .NOT. config_flags%periodic_x .AND. i_end .EQ. ide-1) THEN
7435 rvort(ide,k,j) = rvort(ide-1,k,j)
7440 IF ( .NOT. config_flags%periodic_y .AND. j_end .EQ. jde-1) THEN
7443 rvort(i,k,jde) = rvort(i,k,jde-1)
7448 ! End update of boundary for the vorticity.
7449 !-----------------------------------------------------------------------
7451 !-----------------------------------------------------------------------
7452 ! Calculate an 8-point average of w used to complete the
7453 ! calculation of updraft helicity = w ( dv/dx - du/dy ).
7455 DO j = j_start, j_end
7457 DO i = i_start, i_end
7458 wavg(i,k,j) = 0.125 * ( &
7459 w(i,k ,j ) + w(i-1,k ,j ) + &
7460 w(i,k ,j-1) + w(i-1,k ,j-1) + &
7461 w(i,k+1,j ) + w(i-1,k+1,j ) + &
7462 w(i,k+1,j-1) + w(i-1,k+1,j-1) )
7467 ! End addition of the average w calculation for updraft helicity.
7468 !-----------------------------------------------------------------------
7470 !-----------------------------------------------------------------------
7471 ! Complete the updraft helicity calculation, only including columns with
7472 ! positive vertical velocity between the 2000 m and 5000 m levels.
7474 DO j = j_start, j_end
7475 DO i = i_start, i_end
7476 use_column(i,j) = .true.
7481 ! We want zl and zu to be the height above ground, so we subtract
7482 ! the terrain elevation from the geopotential/g.
7484 DO j = j_start, j_end
7486 DO i = i_start, i_end
7488 (( ph(i ,k ,j ) + phb(i ,k ,j ) ) / g - ht(i ,j ) ) + &
7489 (( ph(i-1,k ,j ) + phb(i-1,k ,j ) ) / g - ht(i-1,j ) ) + &
7490 (( ph(i ,k ,j-1) + phb(i ,k ,j-1) ) / g - ht(i ,j-1) ) + &
7491 (( ph(i-1,k ,j-1) + phb(i-1,k ,j-1) ) / g - ht(i-1,j-1) ) ) )
7494 (( ph(i ,k+1,j ) + phb(i ,k+1,j ) ) / g - ht(i ,j ) ) + &
7495 (( ph(i-1,k+1,j ) + phb(i-1,k+1,j ) ) / g - ht(i-1,j ) ) + &
7496 (( ph(i ,k+1,j-1) + phb(i ,k+1,j-1) ) / g - ht(i ,j-1) ) + &
7497 (( ph(i-1,k+1,j-1) + phb(i-1,k+1,j-1) ) / g - ht(i-1,j-1) ) ) )
7499 IF ( zl .GE. 2000. .AND. zu .LE. 5000. ) THEN
7500 IF ( wavg(i,k,j) .GT. 0. .AND. wavg(i,k+1,j) .GT. 0. ) THEN
7501 uh(i,j) = uh(i,j) + ( ( wavg(i,k,j) * rvort(i,k,j) + &
7502 wavg(i,k+1,j) * rvort(i,k+1,j) ) * 0.5 ) &
7505 use_column(i,j) = .false.
7515 DO j = MAX(jds+1,jts),MIN(jde-2,jte)
7516 DO i = MAX(ids+1,its),MIN(ide-2,ite)
7517 uh_smth = 0.25 * uh(i ,j ) + &
7518 0.125 * ( uh(i+1,j ) + uh(i-1,j ) + &
7519 uh(i ,j+1) + uh(i ,j-1) ) + &
7520 0.0625 * ( uh(i+1,j+1) + uh(i+1,j-1) + &
7521 uh(i-1,j+1) + uh(i-1,j-1) )
7523 IF ( use_column(i,j) ) THEN
7524 IF ( uh_smth .GT. up_heli_max(i,j) ) THEN
7525 up_heli_max(i,j) = uh_smth
7529 ! IF ( up_heli_max(i,j) .GT. 1.0e+3 ) THEN
7530 ! print *,' i,j,up_heli_max = ', i,j,up_heli_max(i,j)
7535 ! End addition of the average w calculation for updraft helicity.
7536 !-----------------------------------------------------------------------
7538 !-----------------------------------------------------------------------
7539 ! Update the boundary for updraft helicity (might need to change later).
7541 IF ( .NOT. config_flags%periodic_x .AND. i_start .EQ. ids+1 ) THEN
7544 ! rvort(ids,k,j) = rvort(ids+1,k,j)
7545 up_heli_max(ids,j) = up_heli_max(ids+1,j)
7550 IF ( .NOT. config_flags%periodic_y .AND. j_start .EQ. jds+1) THEN
7553 ! rvort(i,k,jds) = rvort(i,k,jds+1)
7554 up_heli_max(i,jds) = up_heli_max(i,jds+1)
7559 IF ( .NOT. config_flags%periodic_x .AND. i_end .EQ. ide-1) THEN
7562 ! rvort(ide,k,j) = rvort(ide-1,k,j)
7563 up_heli_max(ide,j) = up_heli_max(ide-1,j)
7568 IF ( .NOT. config_flags%periodic_y .AND. j_end .EQ. jde-1) THEN
7571 ! rvort(i,k,jde) = rvort(i,k,jde-1)
7572 up_heli_max(i,jde) = up_heli_max(i,jde-1)
7577 ! End update of boundary for updraft helicity.
7579 END SUBROUTINE cal_helicity
7580 !=======================================================================
7581 !=======================================================================
7583 SUBROUTINE tridiag(n,a,b,c,d)
7584 !! to solve system of linear eqs on tridiagonal matrix n times n
7585 !! after Peaceman and Rachford, 1955
7586 !! a,b,c,d - are vectors of order n
7587 !! a,b,c - are coefficients on the LHS
7588 !! d - is initially RHS on the output becomes a solution vector
7589 !-------------------------------------------------------------------
7591 INTEGER, INTENT(IN) :: n
7592 REAL, DIMENSION(n), INTENT(IN) :: a,b
7593 REAL, DIMENSION(n), INTENT(INOUT) :: c,d
7597 REAL, DIMENSION(n) :: q
7604 p = 1./(b(i) + a(i)*q(i-1))
7606 d(i) = (d(i) - a(i)*d(i-1))*p
7610 d(i) = d(i) + q(i)*d(i+1)
7613 END SUBROUTINE tridiag
7614 !=======================================================================
7615 !=======================================================================
7616 SUBROUTINE vertical_diffusion_implicit(ru_tendf, rv_tendf, rw_tendf, rt_tendf,&
7617 tke_tendf, moist_tendf, n_moist, &
7618 chem_tendf, n_chem, &
7619 scalar_tendf, n_scalar, &
7620 tracer_tendf, n_tracer, &
7621 u_2, v_2, w_2, dt, &
7622 thp,u_base,v_base,t_base,qv_base,mu,tke, &
7623 theta, config_flags, &
7624 moist,chem,scalar,tracer, &
7625 xkmv,xkhv,xkmh,km_opt, &
7626 fnm, fnp, dn, dnw, rdz, rdzw, &
7627 hfx, qfx, ust, rho, &
7628 nlflux,gamu,gamv,xkmv_meso,l_diss, &
7629 c1h, c2h, c1f, c2f, &
7630 ids, ide, jds, jde, kds, kde, &
7631 ims, ime, jms, jme, kms, kme, &
7632 its, ite, jts, jte, kts, kte )
7633 !To solve the vertical diffusion equations for u,v,w,th,moist,
7634 !chem,scalar,tracer using implicit method.
7635 !The vertical TKE diffusion and dissipation are implicitly treated.
7636 !-----------------------------------------------------------------------
7637 ! Begin declarations.
7641 TYPE(grid_config_rec_type), INTENT(IN ) :: config_flags
7643 INTEGER , INTENT(IN ) :: ids, ide, jds, jde, kds, kde, &
7644 ims, ime, jms, jme, kms, kme, &
7645 its, ite, jts, jte, kts, kte
7647 INTEGER , INTENT(IN ) :: n_moist,n_chem,n_scalar,n_tracer,km_opt
7648 REAL, INTENT(IN ) :: dt
7650 REAL , DIMENSION( kms:kme ) , INTENT(IN ) :: fnm
7652 REAL , DIMENSION( kms:kme ) , INTENT(IN ) :: fnp
7654 REAL , DIMENSION( kms:kme ) , INTENT(IN ) :: dnw
7656 REAL , DIMENSION( kms:kme ) , INTENT(IN ) :: dn
7658 REAL , DIMENSION( ims:ime , jms:jme ) , INTENT(IN ) :: mu
7660 REAL , DIMENSION( kms:kme ) , INTENT(IN ) :: qv_base
7662 REAL , DIMENSION( kms:kme ) , INTENT(IN ) :: u_base
7664 REAL , DIMENSION( kms:kme ) , INTENT(IN ) :: v_base
7666 REAL , DIMENSION( kms:kme ) , INTENT(IN ) :: t_base
7668 REAL , DIMENSION( kms:kme) , INTENT(IN ) :: &
7671 REAL , DIMENSION( ims:ime, kms:kme, jms:jme), INTENT(INOUT) ::ru_tendf,&
7677 REAL , DIMENSION( ims:ime, kms:kme, jms:jme, n_moist), &
7678 INTENT(INOUT) :: moist_tendf
7680 REAL , DIMENSION( ims:ime, kms:kme, jms:jme, n_chem), &
7681 INTENT(INOUT) :: chem_tendf
7683 REAL , DIMENSION( ims:ime, kms:kme, jms:jme, n_scalar) , &
7684 INTENT(INOUT) :: scalar_tendf
7686 REAL , DIMENSION( ims:ime, kms:kme, jms:jme, n_tracer) , &
7687 INTENT(INOUT) :: tracer_tendf
7689 REAL , DIMENSION( ims:ime, kms:kme, jms:jme, n_moist), &
7690 INTENT(INOUT) :: moist
7692 REAL , DIMENSION( ims:ime, kms:kme, jms:jme, n_chem), &
7693 INTENT(INOUT) :: chem
7695 REAL , DIMENSION( ims:ime, kms:kme, jms:jme, n_scalar) , &
7696 INTENT(IN ) :: scalar
7698 REAL , DIMENSION( ims:ime, kms:kme, jms:jme, n_tracer) , &
7699 INTENT(IN ) :: tracer
7701 REAL , DIMENSION( ims:ime, kms:kme, jms:jme), &
7702 INTENT(IN ) :: xkmv, &
7713 REAL , DIMENSION( ims:ime, kms:kme, jms:jme), INTENT(IN ) :: rho
7715 REAL , DIMENSION( ims:ime, jms:jme), INTENT(INOUT) :: hfx, &
7718 REAL , DIMENSION( ims:ime, jms:jme), INTENT(IN ) :: ust
7720 REAL , DIMENSION( ims:ime, kms:kme, jms:jme), INTENT(INOUT) :: thp
7722 REAL , DIMENSION( ims:ime, kms:kme, jms:jme ), INTENT(INOUT) :: nlflux
7724 REAL , DIMENSION( ims:ime, jms:jme), INTENT(INOUT) :: gamu,&
7727 REAL, DIMENSION( ims:ime, kms:kme, jms:jme ), INTENT(INOUT) :: xkmv_meso
7729 REAL, DIMENSION( ims:ime, kms:kme, jms:jme ), INTENT(INOUT) :: l_diss
7731 REAL , DIMENSION( ims:ime, kms:kme, jms:jme) :: var_mix
7733 REAL, DIMENSION( its-1:ite+1, kts:kte, jts-1:jte+1 ):: xkxavg_m, &
7735 REAL, DIMENSION( its:ite, kts:kte, jts:jte) :: xkxavg, &
7739 REAL, DIMENSION( its-1:ite+1, jts-1:jte+1 ) :: gamvavg, gamuavg
7740 REAL, DIMENSION( its-1:ite+1, jts-1:jte+1 ) :: tao_xz, tao_yz
7741 REAL, DIMENSION( its-1:ite+1, kts:kte, jts-1:jte+1 ):: muavg_u,muavg_v
7742 REAL, DIMENSION( its-1:ite+1, kts:kte, jts-1:jte+1 ):: rdz_u, rdz_v
7743 REAL, DIMENSION(kts:kte) :: a, b, c, d
7744 REAL, DIMENSION(kts:kte-1) :: a1, b1, c1, d1
7746 INTEGER :: im, i,j,k,ktf,nz
7747 INTEGER :: i_start, i_end, j_start, j_end
7749 REAL :: V0_u,V0_v,ustar,beta
7750 REAL :: heat_flux, moist_flux, qfac
7751 REAL :: cpm,rdz_w,rhoavg_u,rhoavg_v,rdz_z
7753 !-----------------------------------------------------------------------
7754 !!============================================
7756 !!============================================
7762 j_end = MIN(jte,jde-1)
7764 IF ( config_flags%open_xs .or. config_flags%specified .or. &
7765 config_flags%nested) i_start = MAX(ids+1,its)
7766 IF ( config_flags%open_xe .or. config_flags%specified .or. &
7767 config_flags%nested) i_end = MIN(ide-1,ite)
7768 IF ( config_flags%open_ys .or. config_flags%specified .or. &
7769 config_flags%nested) j_start = MAX(jds+1,jts)
7770 IF ( config_flags%open_ye .or. config_flags%specified .or. &
7771 config_flags%nested) j_end = MIN(jde-2,jte)
7772 IF ( config_flags%periodic_x ) i_start = its
7773 IF ( config_flags%periodic_x ) i_end = ite
7775 DO j = j_start, j_end
7777 DO i = i_start, i_end
7778 rhoavg_u = 0.5 * ( fnm(k) * ( rho(i,k ,j) + rho(i-1,k,j) ) + &
7779 fnp(k) * ( rho(i,k-1,j) + rho(i-1,k-1,j) ) )
7780 xkxavg_m(i,k,j) = 0.5 * rhoavg_u * &
7781 ( fnm(k) * ( xkmv_meso(i,k ,j) + xkmv_meso(i-1,k,j) ) + &
7782 fnp(k) * ( xkmv_meso(i,k-1,j) + xkmv_meso(i-1,k-1,j) ) )
7783 xkxavg_s(i,k,j) = 0.5 * rhoavg_u * &
7784 ( fnm(k) * ( xkmv(i,k ,j) + xkmv(i-1,k,j) ) + &
7785 fnp(k) * ( xkmv(i,k-1,j) + xkmv(i-1,k-1,j) ) )
7790 DO j = j_start, j_end
7791 DO i = i_start, i_end
7792 xkxavg_m(i,kts, j) = 0.0
7793 xkxavg_m(i,ktf+1,j) = 0.0
7794 xkxavg_s(i,kts, j) = 0.0
7795 xkxavg_s(i,ktf+1,j) = 0.0
7799 DO j = j_start, j_end
7800 DO i = i_start, i_end
7801 gamuavg(i,j) = 0.5 * ( gamu(i,j) + gamu(i-1,j) )
7805 DO j = j_start, j_end
7807 DO i = i_start, i_end
7808 rdz_u(i,k,j) = 2./( 1./rdz(i,k,j) + 1./rdz(i-1,k,j))
7809 muavg_u(i,k,j) = 0.5 * ( (c1h(k)*mu(i,j)+c2h(k)) + (c1h(k)*mu(i-1,j)+c2h(k)) )
7814 DO j = j_start, j_end
7817 V0_u = sqrt((u_2(i,kts,j)**2)+ &
7818 (((v_2(i ,kts,j )+ &
7821 v_2(i-1,kts,j+1))/4)**2))+epsilon
7822 ustar = 0.5*(ust(i,j)+ust(i-1,j))
7823 tao_xz(i,j) = ustar*ustar*(rho(i,kts,j)+rho(i-1,kts,j))/(2.*V0_u)
7829 DO j = j_start, j_end
7830 DO i = i_start, i_end
7832 rdz_w = -g/(dnw(k))/muavg_u(i,k,j)
7834 b(k) = 1.+rdz_w*(rdz_u(i,k+1,j)*xkxavg_s(i,k+1,j)+tao_xz(i,j))*dt
7835 c(k) = -rdz_w*rdz_u(i,k+1,j)*xkxavg_s(i,k+1,j)*dt
7839 rdz_w = -g/(dnw(k))/muavg_u(i,k,j)
7840 a(k) = -rdz_w*rdz_u(i,k,j)*xkxavg_s(i,k,j)*dt
7841 b(k) = 1 + rdz_w*(rdz_u(i,k+1,j)*xkxavg_s(i,k+1,j)+rdz_u(i,k,j)*xkxavg_s(i,k,j))*dt
7842 c(k) = -rdz_w*rdz_u(i,k+1,j)*xkxavg_s(i,k+1,j)*dt
7843 d(k) = -rdz_w*(xkxavg_m(i,k+1,j)-xkxavg_m(i,k,j))*dt*gamuavg(i,j) + u_2(i,k,j)
7849 d(ktf) = u_2(i,ktf,j)
7851 CALL tridiag(nz,a,b,c,d)
7854 ru_tendf(i,k,j) = ru_tendf(i,k,j) + muavg_u(i,k,j) * (d(k) - u_2(i,k,j))/dt
7860 !!============================================
7862 !!============================================
7866 i_end = MIN(ite,ide-1)
7870 IF ( config_flags%open_xs .or. config_flags%specified .or. &
7871 config_flags%nested) i_start = MAX(ids+1,its)
7872 IF ( config_flags%open_xe .or. config_flags%specified .or. &
7873 config_flags%nested) i_end = MIN(ide-2,ite)
7874 IF ( config_flags%open_ys .or. config_flags%specified .or. &
7875 config_flags%nested) j_start = MAX(jds+1,jts)
7876 IF ( config_flags%open_ye .or. config_flags%specified .or. &
7877 config_flags%nested) j_end = MIN(jde-1,jte)
7878 IF ( config_flags%periodic_x ) i_start = its
7879 IF ( config_flags%periodic_x ) i_end = MIN(ite,ide-1)
7881 DO j = j_start, j_end
7883 DO i = i_start, i_end
7884 rhoavg_v = 0.5 * ( fnm(k) * ( rho(i,k ,j) + rho(i,k,j-1) ) + &
7885 fnp(k) * ( rho(i,k-1,j) + rho(i,k-1,j-1) ) )
7886 xkxavg_m(i,k,j) = 0.5 * rhoavg_v * &
7887 ( fnm(k) * ( xkmv_meso(i,k ,j) + xkmv_meso(i,k,j-1) ) + &
7888 fnp(k) * ( xkmv_meso(i,k-1,j) + xkmv_meso(i,k-1,j-1) ) )
7889 xkxavg_s(i,k,j) = 0.5 * rhoavg_v * &
7890 ( fnm(k) * ( xkmv(i,k ,j) + xkmv(i,k,j-1) ) + &
7891 fnp(k) * ( xkmv(i,k-1,j) + xkmv(i,k-1,j-1) ) )
7896 DO j = j_start, j_end
7897 DO i = i_start, i_end
7898 xkxavg_m(i,kts, j) = 0.0
7899 xkxavg_m(i,ktf+1,j) = 0.0
7900 xkxavg_s(i,kts, j) = 0.0
7901 xkxavg_s(i,ktf+1,j) = 0.0
7905 DO j = j_start, j_end
7906 DO i = i_start, i_end
7907 gamvavg(i,j) = 0.5 * ( gamv(i,j) + gamv(i,j-1) )
7911 DO j = j_start, j_end
7913 DO i = i_start, i_end
7914 muavg_v(i,k,j) = 0.5 * ( (c1h(k)*mu(i,j)+c2h(k)) + (c1h(k)*mu(i,j-1)+c2h(k)) )
7919 DO j = j_start, j_end
7921 DO i = i_start, i_end
7922 rdz_v(i,k,j) = 2./( 1./rdz(i,k,j) + 1./rdz(i,k,j-1))
7928 DO i = i_start, i_end
7930 V0_v = sqrt((v_2(i,kts,j)**2) + &
7931 (((u_2(i ,kts,j )+ &
7934 u_2(i+1,kts,j-1))/4)**2))+epsilon
7935 ustar = 0.5*(ust(i,j)+ust(i,j-1))
7936 tao_yz(i,j) = ustar*ustar*(rho(i,kts,j)+rho(i,kts,j-1))/(2.*V0_v)
7942 DO j = j_start, j_end
7943 DO i = i_start, i_end
7945 rdz_w = -g/(dnw(k))/muavg_v(i,k,j)
7947 b(k) = 1. + rdz_w*(rdz_v(i,k+1,j)*xkxavg_s(i,k+1,j)+tao_yz(i,j))*dt
7948 c(k) = -rdz_w* rdz_v(i,k+1,j)*xkxavg_s(i,k+1,j)*dt
7952 rdz_w = -g/(dnw(k))/muavg_v(i,k,j)
7953 a(k) = -rdz_w*rdz_v(i,k,j)*xkxavg_s(i,k,j)*dt
7954 b(k) = 1. + rdz_w*(rdz_v(i,k+1,j)*xkxavg_s(i,k+1,j)+rdz_v(i,k,j)*xkxavg_s(i,k,j))*dt
7955 c(k) = -rdz_w*rdz_v(i,k+1,j)*xkxavg_s(i,k+1,j)*dt
7956 d(k) = -rdz_w*(xkxavg_m(i,k+1,j)-xkxavg_m(i,k,j))*dt*gamvavg(i,j) + v_2(i,k,j)
7962 d(ktf) = v_2(i,ktf,j)
7964 CALL tridiag(nz,a,b,c,d)
7967 rv_tendf(i,k,j) = rv_tendf(i,k,j) + muavg_v(i,k,j) * (d(k) - v_2(i,k,j))/dt
7973 !!============================================
7975 !!============================================
7979 i_end = MIN(ite,ide-1)
7981 j_end = MIN(jte,jde-1)
7983 IF ( config_flags%open_xs .or. config_flags%specified .or. &
7984 config_flags%nested) i_start = MAX(ids+1,its)
7985 IF ( config_flags%open_xe .or. config_flags%specified .or. &
7986 config_flags%nested) i_end = MIN(ide-2,ite)
7987 IF ( config_flags%open_ys .or. config_flags%specified .or. &
7988 config_flags%nested) j_start = MAX(jds+1,jts)
7989 IF ( config_flags%open_ye .or. config_flags%specified .or. &
7990 config_flags%nested) j_end = MIN(jde-2,jte)
7991 IF ( config_flags%periodic_x ) i_start = its
7992 IF ( config_flags%periodic_x ) i_end = MIN(ite,ide-1)
7996 DO j = j_start, j_end
7998 DO i = i_start, i_end
7999 xkxavg_w(i,k,j) = xkmh(i,k,j)*rho(i,k,j)
8004 DO j = j_start, j_end
8005 DO i = i_start, i_end
8008 rdz_z = -g/dn(k)/(c1f(k)*mu(i,j) + c2f(k))
8009 a1(k-1) = - rdz_z* rdzw(i,k-1,j)*xkxavg_w(i,k-1,j)*dt
8010 b1(k-1) = 1. + rdz_z*(rdzw(i,k ,j)*xkxavg_w(i,k,j)+rdzw(i,k-1,j)*xkxavg_w(i,k-1,j))*dt
8011 c1(k-1) = - rdz_z* rdzw(i,k ,j)*xkxavg_w(i,k,j)*dt
8012 d1(k-1) = w_2(i,k,j)
8015 CALL tridiag(nz,a1,b1,c1,d1)
8018 rw_tendf(i,k,j) = rw_tendf(i,k,j) + (c1f(k)*mu(i,j) + c2f(k)) * (d1(k-1)-w_2(i,k,j))/dt
8024 !!============================================
8026 !!============================================
8030 i_end = MIN(ite,ide-1)
8032 j_end = MIN(jte,jde-1)
8034 IF ( config_flags%open_xs .or. config_flags%specified .or. &
8035 config_flags%nested) i_start = MAX(ids+1,its)
8036 IF ( config_flags%open_xe .or. config_flags%specified .or. &
8037 config_flags%nested) i_end = MIN(ide-2,ite)
8038 IF ( config_flags%open_ys .or. config_flags%specified .or. &
8039 config_flags%nested) j_start = MAX(jds+1,jts)
8040 IF ( config_flags%open_ye .or. config_flags%specified .or. &
8041 config_flags%nested) j_end = MIN(jde-2,jte)
8042 IF ( config_flags%periodic_x ) i_start = its
8043 IF ( config_flags%periodic_x ) i_end = MIN(ite,ide-1)
8045 IF ( config_flags%mix_full_fields ) THEN
8046 DO j = jts,min(jte,jde-1)
8048 DO i = its,min(ite,ide-1)
8049 var_mix(i,k,j) = thp(i,k,j)
8054 DO j = jts,min(jte,jde-1)
8056 DO i = its,min(ite,ide-1)
8057 var_mix(i,k,j) = thp(i,k,j) - t_base(k)
8063 DO j = j_start, j_end
8065 DO i = i_start, i_end
8066 xkxavg(i,k,j) = (fnm(k)*xkhv(i,k,j)+fnp(k)*xkhv(i,k-1,j)) * &
8067 (fnm(k)*rho (i,k,j)+fnp(k)*rho (i,k-1,j))
8072 DO j = j_start, j_end
8074 DO i = i_start, i_end
8075 nlflux_rho(i,k,j) = nlflux(i,k,j) * &
8076 (fnm(k)*rho (i,k,j)+fnp(k)*rho (i,k-1,j))
8081 DO j = j_start, j_end
8082 DO i = i_start, i_end
8083 xkxavg(i,kts ,j) = 0.0
8084 xkxavg(i,ktf+1,j) = 0.0
8090 hflux: SELECT CASE( config_flags%isfflx )
8091 CASE (0,2) ! with fixed surface heat flux given in the namelist
8092 heat_flux = config_flags%tke_heat_flux ! constant heat flux value
8093 ! set in namelist.input
8094 DO j = j_start, j_end
8095 DO i = i_start, i_end
8096 cpm = cp * (1. + 0.8 * moist(i,kts,j,P_QV))
8097 hfx(i,j)=heat_flux*cpm*rho(i,1,j)
8101 CASE (1) ! use surface heat flux computed from surface routine
8104 CALL wrf_error_fatal( 'isfflx value invalid for diff_opt=2' )
8107 DO j = j_start, j_end
8108 DO i = i_start, i_end
8109 cpm = cp * (1. + 0.8 * moist(i,kts,j,P_QV))
8111 rdz_w = -g/dnw(k)/(c1h(k)*mu(i,j) + c2h(k))
8113 b(k) = 1. + rdz_w*(rdz(i,k+1,j)*xkxavg(i,k+1,j))*dt
8114 c(k) = - rdz_w* rdz(i,k+1,j)*xkxavg(i,k+1,j) *dt
8116 IF(config_flags%use_theta_m == 0)THEN
8117 d(k) = var_mix(i,kts,j) &
8118 - dt*rdz_w*nlflux_rho(i,k+1,j) &
8119 + dt*rdz_w*hfx(i,j)/cpm
8121 qfac = 1.+rvovrd*moist(i,kts,j,P_QV)
8122 d(k) = var_mix(i,kts,j) &
8123 - dt*rdz_w*nlflux_rho(i,k+1,j)*qfac &
8124 + dt*rdz_w*(qfac*hfx(i,j)/cpm+1.61*theta(i,kts,j)*qfx(i,j))
8127 IF(config_flags%use_theta_m == 1)THEN
8128 qfac = 1.+rvovrd*moist(i,kts,j,P_QV)
8133 rdz_w = -g/dnw(k)/(c1h(k)*mu(i,j) + c2h(k))
8134 a(k) = -rdz_w*rdz(i,k,j)*xkxavg(i,k,j)*dt
8135 b(k) = 1.+rdz_w*(rdz(i,k+1,j)*xkxavg(i,k+1,j)+rdz(i,k,j)*xkxavg(i,k,j))*dt
8136 c(k) = -rdz_w*rdz(i,k+1,j)*xkxavg(i,k+1,j)*dt
8137 d(k) = -rdz_w*(nlflux_rho(i,k+1,j)-nlflux_rho(i,k,j))*qfac*dt + var_mix(i,k,j)
8143 d(ktf) = var_mix(i,ktf,j)
8145 CALL tridiag(nz,a,b,c,d)
8148 rt_tendf(i,k,j) = rt_tendf(i,k,j) + (c1h(k)*mu(i,j) + c2h(k)) * (d(k) - var_mix(i,k,j))/dt
8154 !!============================================
8155 !! tke(including implicit treatments of vertical diffusion and disspation)
8156 !!============================================
8157 DO j = j_start, j_end
8159 DO i = i_start, i_end
8160 xkxavg(i,k,j) = ( fnm(k) * xkhv(i,k,j) + fnp(k) * xkhv(i,k-1,j) ) &
8161 *( fnm(k) * rho (i,k,j) + fnp(k) * rho (i,k-1,j) )
8166 DO j = j_start, j_end
8167 DO i = i_start, i_end
8168 xkxavg(i,kts, j) = 0.0
8169 xkxavg(i,ktf+1,j) = 0.0
8175 DO j = j_start, j_end
8176 DO i = i_start, i_end
8178 rdz_w = - g/dnw(k)/(c1h(k)*mu(i,j)+c2h(k))
8179 IF (l_diss(i,k,j) .ne. 0) THEN
8180 beta = 1.5*sqrt(tke(i,k,j))/l_diss(i,k,j)
8184 a(k) = - 2.0*xkxavg(i,k,j)*dt*rdz_w*rdz(i,k,j)
8185 b(k) = 1.0 + 2.0*dt*rdz_w*(rdz(i,k,j)*xkxavg(i,k,j) &
8186 + rdz(i,k+1,j)*xkxavg(i,k+1,j)) &
8188 c(k) = - 2.0*xkxavg(i,k+1,j)*dt*rdz(i,k+1,j)*rdz_w
8191 IF (l_diss(i,k,j) .ne. 0) THEN
8192 d(k) = d(k) + 0.5*dt*tke(i,k,j)**1.5/l_diss(i,k,j)
8201 CALL tridiag(nz,a,b,c,d)
8204 tke_tendf(i,k,j) = tke_tendf(i,k,j) + (c1h(k)*mu(i,j)+c2h(k)) * (d(k)-tke(i,k,j))/dt
8210 !!============================================
8212 !!============================================
8213 DO j = j_start, j_end
8215 DO i = i_start, i_end
8216 xkxavg(i,k,j) = ( fnm(k) * xkhv(i,k,j) + fnp(k) * xkhv(i,k-1,j) ) &
8217 *( fnm(k) * rho (i,k,j) + fnp(k) * rho (i,k-1,j) )
8222 DO j = j_start, j_end
8223 DO i = i_start, i_end
8224 xkxavg(i,kts, j) = 0.0
8225 xkxavg(i,ktf+1,j) = 0.0
8229 IF (n_moist .ge. PARAM_FIRST_SCALAR) THEN
8230 moist_loop: do im = PARAM_FIRST_SCALAR, n_moist
8231 IF ( (.not. config_flags%mix_full_fields) .and. (im == P_QV) ) THEN
8232 DO j = jts,min(jte,jde-1)
8234 DO i = its,min(ite,ide-1)
8235 var_mix(i,k,j) = moist(i,k,j,im) - qv_base(k)
8240 DO j = jts,min(jte,jde-1)
8242 DO i = its,min(ite,ide-1)
8243 var_mix(i,k,j) = moist(i,k,j,im)
8251 IF ( im == P_QV ) THEN
8252 DO j = j_start, j_end
8253 DO i = i_start, i_end
8255 rdz_w = -g/dnw(k)/(c1h(k)*mu(i,j)+c2h(k))
8257 b(k) = 1.+rdz_w*(rdz(i,k+1,j)*xkxavg(i,k+1,j))*dt
8258 c(k) = -rdz_w* rdz(i,k+1,j)*xkxavg(i,k+1,j) *dt
8259 ! d(k) = var_mix(i,k,j) + dt*rdz_w*qfx(i,j)/(1.+moist(i,k,j,P_QV)) !??
8260 d(k) = var_mix(i,k,j) + dt*rdz_w*qfx(i,j)
8263 rdz_w = -g/dnw(k)/(c1h(k)*mu(i,j)+c2h(k))
8264 a(k) = -rdz_w*rdz(i,k,j)*xkxavg(i,k,j)*dt
8265 b(k) = 1.+rdz_w*(rdz(i,k+1,j)*xkxavg(i,k+1,j)+rdz(i,k,j)*xkxavg(i,k,j))*dt
8266 c(k) = -rdz_w*rdz(i,k+1,j)*xkxavg(i,k+1,j)*dt
8267 d(k) = var_mix(i,k,j)
8273 d(ktf) = var_mix(i,ktf,j)
8275 CALL tridiag(nz,a,b,c,d)
8278 moist_tendf(i,k,j,im) = moist_tendf(i,k,j,im) + (c1h(k)*mu(i,j)+c2h(k)) * (d(k)-var_mix(i,k,j))/dt
8285 IF ( im /= P_QV ) THEN
8286 DO j = j_start, j_end
8287 DO i = i_start, i_end
8289 rdz_w = -g/dnw(k)/(c1h(k)*mu(i,j)+c2h(k))
8291 b(k) = 1.+rdz_w*(rdz(i,k+1,j)*xkxavg(i,k+1,j))*dt
8292 c(k) = -rdz_w* rdz(i,k+1,j)*xkxavg(i,k+1,j)*dt
8293 d(k) = var_mix(i,k,j)
8296 rdz_w = -g/dnw(k)/(c1h(k)*mu(i,j)+c2h(k))
8297 a(k) = -rdz_w*rdz(i,k,j)*xkxavg(i,k,j)*dt
8298 b(k) = 1.+rdz_w*(rdz(i,k+1,j)*xkxavg(i,k+1,j)+rdz(i,k,j)*xkxavg(i,k,j))*dt
8299 c(k) = -rdz_w*rdz(i,k+1,j)*xkxavg(i,k+1,j)*dt
8300 d(k) = var_mix(i,k,j)
8306 d(ktf) = var_mix(i,ktf,j)
8308 CALL tridiag(nz,a,b,c,d)
8311 moist_tendf(i,k,j,im) = moist_tendf(i,k,j,im) + (c1h(k)*mu(i,j)+c2h(k)) * (d(k)-var_mix(i,k,j))/dt
8320 !!============================================
8322 !!============================================
8323 DO j = j_start, j_end
8325 DO i = i_start, i_end
8326 xkxavg(i,k,j) = ( fnm(k) * xkhv(i,k,j) + fnp(k) * xkhv(i,k-1,j) ) &
8327 *( fnm(k) * rho (i,k,j) + fnp(k) * rho (i,k-1,j) )
8332 DO j = j_start, j_end
8333 DO i = i_start, i_end
8334 xkxavg(i,kts ,j) = 0.0
8335 xkxavg(i,ktf+1,j) = 0.0
8339 IF (n_scalar .ge. PARAM_FIRST_SCALAR) THEN
8340 scalar_loop: do im = PARAM_FIRST_SCALAR, n_scalar
8341 ! need to avoid mixing scalars associated with precipitating species (e.g. Nr)
8342 IF(im.NE.P_QNS .AND. im.NE.P_QNR .AND. im.NE.P_QNG .AND. im.NE.P_QNH .AND. &
8343 im.NE.P_QT .AND. im.NE.P_QVOLG) THEN
8344 DO j = jts, min(jte,jde-1)
8346 DO i = its, min(ite,ide-1)
8347 var_mix(i,k,j) = scalar(i,k,j,im)
8352 DO j = j_start, j_end
8353 DO i = i_start, i_end
8355 rdz_w = -g/dnw(k)/(c1h(k)*mu(i,j)+c2h(k))
8357 b(k) = 1.+rdz_w*(rdz(i,k+1,j)*xkxavg(i,k+1,j))*dt
8358 c(k) = -rdz_w* rdz(i,k+1,j)*xkxavg(i,k+1,j)*dt
8359 d(k) = var_mix(i,k,j)
8362 rdz_w = -g/dnw(k)/(c1h(k)*mu(i,j)+c2h(k))
8363 a(k) = -rdz_w*rdz(i,k,j)*xkxavg(i,k,j)*dt
8364 b(k) = 1.+rdz_w*(rdz(i,k+1,j)*xkxavg(i,k+1,j)+rdz(i,k,j)*xkxavg(i,k,j))*dt
8365 c(k) = -rdz_w*rdz(i,k+1,j)*xkxavg(i,k+1,j)*dt
8366 d(k) = var_mix(i,k,j)
8372 d(ktf) = var_mix(i,ktf,j)
8374 CALL tridiag(nz,a,b,c,d)
8377 scalar_tendf(i,k,j,im) = scalar_tendf(i,k,j,im) + (c1h(k)*mu(i,j)+c2h(k)) * (d(k)-var_mix(i,k,j))/dt
8385 !!============================================
8387 !!============================================
8388 IF (n_tracer .ge. PARAM_FIRST_SCALAR) THEN
8389 tracer_loop: do im = PARAM_FIRST_SCALAR, n_tracer
8390 DO j = jts, min(jte,jde-1)
8392 DO i = its, min(ite,ide-1)
8393 var_mix(i,k,j) = tracer(i,k,j,im)
8398 DO j = j_start, j_end
8399 DO i = i_start, i_end
8401 rdz_w = -g/dnw(k)/(c1h(k)*mu(i,j)+c2h(k))
8403 b(k) = 1.+rdz_w*(rdz(i,k+1,j)*xkxavg(i,k+1,j))*dt
8404 c(k) = -rdz_w* rdz(i,k+1,j)*xkxavg(i,k+1,j)*dt
8405 d(k) = var_mix(i,kts,j)
8408 rdz_w = -g/dnw(k)/(c1h(k)*mu(i,j)+c2h(k))
8409 a(k) = -rdz_w*rdz(i,k,j)*xkxavg(i,k,j)*dt
8410 b(k) = 1.+rdz_w*(rdz(i,k+1,j)*xkxavg(i,k+1,j)+rdz(i,k,j)*xkxavg(i,k,j))*dt
8411 c(k) = -rdz_w*rdz(i,k+1,j)*xkxavg(i,k+1,j)*dt
8412 d(k) = var_mix(i,k,j)
8418 d(ktf) = var_mix(i,ktf,j)
8420 CALL tridiag(nz,a,b,c,d)
8423 tracer_tendf(i,k,j,im) = tracer_tendf(i,k,j,im) + (c1h(k)*mu(i,j)+c2h(k)) * (d(k) - var_mix(i,k,j))/dt
8430 !!============================================
8432 !!============================================
8433 IF (n_chem .ge. PARAM_FIRST_SCALAR) THEN
8434 chem_loop: do im = PARAM_FIRST_SCALAR, n_chem
8435 DO j = jts,min(jte,jde-1)
8437 DO i = its,min(ite,ide-1)
8438 var_mix(i,k,j) = chem(i,k,j,im)
8443 DO j = j_start, j_end
8444 DO i = i_start, i_end
8446 rdz_w = -g/dnw(k)/(c1h(k)*mu(i,j)+c2h(k))
8448 b(k) = 1.+rdz_w*(rdz(i,k+1,j)*xkxavg(i,k+1,j))*dt
8449 c(k) = -rdz_w* rdz(i,k+1,j)*xkxavg(i,k+1,j)*dt
8450 d(k) = var_mix(i,kts,j)
8453 rdz_w = -g/dnw(k)/(c1h(k)*mu(i,j)+c2h(k))
8454 a(k) = -rdz_w*rdz(i,k,j)*xkxavg(i,k,j)*dt
8455 b(k) = 1.+rdz_w*(rdz(i,k+1,j)*xkxavg(i,k+1,j)+rdz(i,k,j)*xkxavg(i,k,j))*dt
8456 c(k) = -rdz_w*rdz(i,k+1,j)*xkxavg(i,k+1,j)*dt
8457 d(k) = var_mix(i,k,j)
8463 d(ktf) = var_mix(i,ktf,j)
8465 CALL tridiag(nz,a,b,c,d)
8468 chem_tendf(i,k,j,im) = chem_tendf(i,k,j,im) + (c1h(k)*mu(i,j)+c2h(k)) * (d(k) - var_mix(i,k,j))/dt
8475 END SUBROUTINE vertical_diffusion_implicit
8477 !=======================================================================
8478 !=======================================================================
8479 END MODULE module_diffusion_em
8481 !=======================================================================
8482 !=======================================================================