1 !WRF:MEDIATION_LAYER:INTERPOLATIONFUNCTION
9 ! Note, NMM-specific routines moved to end. 20080612. JM
11 SUBROUTINE interp_fcn ( cfld, & ! CD field
12 cids, cide, ckds, ckde, cjds, cjde, &
13 cims, cime, ckms, ckme, cjms, cjme, &
14 cits, cite, ckts, ckte, cjts, cjte, &
16 nids, nide, nkds, nkde, njds, njde, &
17 nims, nime, nkms, nkme, njms, njme, &
18 nits, nite, nkts, nkte, njts, njte, &
19 shw, & ! stencil half width for interp
20 imask, & ! interpolation mask
21 xstag, ystag, & ! staggering of field
22 ipos, jpos, & ! Position of lower left of nest in CD
23 nri, nrj ) ! nest ratios
29 INTEGER, INTENT(IN) :: cids, cide, ckds, ckde, cjds, cjde, &
30 cims, cime, ckms, ckme, cjms, cjme, &
31 cits, cite, ckts, ckte, cjts, cjte, &
32 nids, nide, nkds, nkde, njds, njde, &
33 nims, nime, nkms, nkme, njms, njme, &
34 nits, nite, nkts, nkte, njts, njte, &
38 LOGICAL, INTENT(IN) :: xstag, ystag
40 REAL, DIMENSION ( cims:cime, ckms:ckme, cjms:cjme ) :: cfld
41 REAL, DIMENSION ( nims:nime, nkms:nkme, njms:njme ) :: nfld
42 INTEGER, DIMENSION ( nims:nime, njms:njme ) :: imask
48 INTEGER ci, cj, ck, ni, nj, nk, ip, jp, ioff, joff, nioff, njoff
53 REAL psca(cims:cime,cjms:cjme,nri*nrj)
54 LOGICAL icmask( cims:cime, cjms:cjme )
59 ! Iterate over the ND tile and compute the values
80 !$OMP PRIVATE ( i,j,k,ni,nj,ci,cj,ip,jp,nk,ck,nf,icmask,psca )
85 nj = (j-jpos) * nrj + ( nrjo2 + 1 ) ! j point on nest
87 ni = (i-ipos) * nri + ( nrio2 + 1 ) ! i point on nest
88 if ( ni .ge. nits-nioff-nrio2 .and. &
89 ni .le. nite+nioff+nrio2 .and. &
90 nj .ge. njts-njoff-nrjo2 .and. &
91 nj .le. njte+njoff+nrjo2 ) then
92 ! if ( imask(ni,nj) .eq. 1 .or. imask(ni-nioff,nj-njoff) .eq. 1 ) then
93 ! icmask( i, j ) = .TRUE.
95 if ( ni.ge.nims.and.ni.le.nime.and.nj.ge.njms.and.nj.le.njme) then
96 if ( imask(ni,nj) .eq. 1 ) then
97 icmask( i, j ) = .TRUE.
100 if ( ni-nioff.ge.nims.and.ni.le.nime.and.nj-njoff.ge.njms.and.nj.le.njme) then
101 if (ni .ge. nits-nioff .and. nj .ge. njts-njoff ) then
102 if ( imask(ni-nioff,nj-njoff) .eq. 1) then
103 icmask( i, j ) = .TRUE.
108 psca(i,j,nf) = cfld(i,k,j)
113 ! tile dims in this call to sint are 1-over to account for the fact
114 ! that the number of cells on the nest local subdomain is not
115 ! necessarily a multiple of the nest ratio in a given dim.
116 ! this could be a little less ham-handed.
121 cims, cime, cjms, cjme, icmask, &
122 cits-1, cite+1, cjts-1, cjte+1, nrj*nri, xstag, ystag )
124 !call end_timing( ' sint ' )
126 DO nj = njts, njte+joff
127 cj = jpos + (nj-1) / nrj ! j coord of CD point
128 jp = mod ( nj-1 , nrj ) ! coord of ND w/i CD point
131 DO ni = nits, nite+ioff
132 ci = ipos + (ni-1) / nri ! i coord of CD point
133 ip = mod ( ni-1 , nri ) ! coord of ND w/i CD point
134 if ( imask ( ni, nj ) .eq. 1 .or. imask ( ni-ioff, nj-joff ) .eq. 1 ) then
135 nfld( ni-ioff, nk, nj-joff ) = psca( ci , cj, ip+1 + (jp)*nri )
140 !$OMP END PARALLEL DO
144 !write(0,'(") cims:cime, ckms:ckme, cjms:cjme ",6i4)')cims,cime, ckms,ckme, cjms,cjme
145 !write(0,'(") nims:nime, nkms:nkme, njms:njme ",6i4)')nims,nime, nkms,nkme, njms,njme
146 !write(0,'(") cits:cite, ckts:ckte, cjts:cjte ",6i4)')cits,cite, ckts,ckte, cjts,cjte
147 !write(0,'(") nits:nite, nkts:nkte, njts:njte ",6i4)')nits,nite, nkts,nkte, njts,njte
150 cj = jpos + (nj-1) / nrj ! j coord of CD point
151 jp = mod ( nj , nrj ) ! coord of ND w/i CD point
155 ci = ipos + (ni-1) / nri ! j coord of CD point
156 ip = mod ( ni , nri ) ! coord of ND w/i CD point
157 ! This is a trivial implementation of the interp_fcn; just copies
158 ! the values from the CD into the ND
159 if ( imask ( ni, nj ) .eq. 1 ) then
160 nfld( ni, nk, nj ) = cfld( ci , ck , cj )
169 END SUBROUTINE interp_fcn
171 !==================================
172 ! this is the default function used in feedback.
174 SUBROUTINE copy_fcn ( cfld, & ! CD field
175 cids, cide, ckds, ckde, cjds, cjde, &
176 cims, cime, ckms, ckme, cjms, cjme, &
177 cits, cite, ckts, ckte, cjts, cjte, &
179 nids, nide, nkds, nkde, njds, njde, &
180 nims, nime, nkms, nkme, njms, njme, &
181 nits, nite, nkts, nkte, njts, njte, &
182 shw, & ! stencil half width for interp
183 imask, & ! interpolation mask
184 xstag, ystag, & ! staggering of field
185 ipos, jpos, & ! Position of lower left of nest in CD
186 nri, nrj ) ! nest ratios
191 INTEGER, INTENT(IN) :: cids, cide, ckds, ckde, cjds, cjde, &
192 cims, cime, ckms, ckme, cjms, cjme, &
193 cits, cite, ckts, ckte, cjts, cjte, &
194 nids, nide, nkds, nkde, njds, njde, &
195 nims, nime, nkms, nkme, njms, njme, &
196 nits, nite, nkts, nkte, njts, njte, &
200 LOGICAL, INTENT(IN) :: xstag, ystag
202 REAL, DIMENSION ( cims:cime, ckms:ckme, cjms:cjme ), INTENT(OUT) :: cfld
203 REAL, DIMENSION ( nims:nime, nkms:nkme, njms:njme ),INTENT(IN) :: nfld
204 INTEGER, DIMENSION ( nims:nime, njms:njme ),INTENT(IN) :: imask
208 INTEGER ci, cj, ck, ni, nj, nk, ip, jp, ioff, joff, ioffa, joffa
209 INTEGER :: icmin,icmax,jcmin,jcmax
210 INTEGER :: istag,jstag, ipoints,jpoints,ijpoints
211 INTEGER , PARAMETER :: passes = 2
214 ! Loop over the coarse grid in the area of the fine mesh. Do not
215 ! process the coarse grid values that are along the lateral BC
216 ! provided to the fine grid. Since that is in the specified zone
217 ! for the fine grid, it should not be used in any feedback to the
218 ! coarse grid as it should not have changed.
220 ! Due to peculiarities of staggering, it is simpler to handle the feedback
221 ! for the staggerings based upon whether it is a even ratio (2::1, 4::1, etc.) or
222 ! an odd staggering ratio (3::1, 5::1, etc.).
224 ! Though there are separate grid ratios for the i and j directions, this code
225 ! is not general enough to handle aspect ratios .NE. 1 for the fine grid cell.
227 ! These are local integer increments in the looping. Basically, istag=1 means
228 ! that we will assume one less point in the i direction. Note that ci and cj
229 ! have a maximum value that is decreased by istag and jstag, respectively.
231 ! Horizontal momentum feedback is along the face, not within the cell. For a
232 ! 3::1 ratio, temperature would use 9 pts for feedback, while u and v use
233 ! only 3 points for feedback from the nest to the parent.
235 CALL nl_get_spec_zone( 1 , spec_zone )
236 istag = 1 ; jstag = 1
237 IF ( xstag ) istag = 0
238 IF ( ystag ) jstag = 0
240 IF( MOD(nrj,2) .NE. 0) THEN ! odd refinement ratio
242 IF ( ( .NOT. xstag ) .AND. ( .NOT. ystag ) ) THEN
243 DO cj = MAX(jpos+spec_zone,cjts),MIN(jpos+(njde-njds)/nrj-jstag-spec_zone,cjte)
244 nj = (cj-jpos)*nrj + jstag + 1
247 DO ci = MAX(ipos+spec_zone,cits),MIN(ipos+(nide-nids)/nri-istag-spec_zone,cite)
248 ni = (ci-ipos)*nri + istag + 1
249 cfld( ci, ck, cj ) = 0.
250 DO ijpoints = 1 , nri * nrj
251 ipoints = MOD((ijpoints-1),nri) + 1 - nri/2 - 1
252 jpoints = (ijpoints-1)/nri + 1 - nrj/2 - 1
253 cfld( ci, ck, cj ) = cfld( ci, ck, cj ) + &
254 1./REAL(nri*nrj) * nfld( ni+ipoints , nk , nj+jpoints )
256 ! cfld( ci, ck, cj ) = 1./9. * &
257 ! ( nfld( ni-1, nk , nj-1) + &
258 ! nfld( ni , nk , nj-1) + &
259 ! nfld( ni+1, nk , nj-1) + &
260 ! nfld( ni-1, nk , nj ) + &
261 ! nfld( ni , nk , nj ) + &
262 ! nfld( ni+1, nk , nj ) + &
263 ! nfld( ni-1, nk , nj+1) + &
264 ! nfld( ni , nk , nj+1) + &
265 ! nfld( ni+1, nk , nj+1) )
270 ELSE IF ( ( xstag ) .AND. ( .NOT. ystag ) ) THEN
271 DO cj = MAX(jpos+spec_zone,cjts),MIN(jpos+(njde-njds)/nrj-jstag-spec_zone,cjte)
272 nj = (cj-jpos)*nrj + jstag + 1
275 DO ci = MAX(ipos+spec_zone,cits),MIN(ipos+(nide-nids)/nri-istag-spec_zone,cite)
276 ni = (ci-ipos)*nri + istag + 1
277 cfld( ci, ck, cj ) = 0.
278 DO ijpoints = (nri+1)/2 , (nri+1)/2 + nri*(nri-1) , nri
279 ipoints = MOD((ijpoints-1),nri) + 1 - nri/2 - 1
280 jpoints = (ijpoints-1)/nri + 1 - nrj/2 - 1
281 cfld( ci, ck, cj ) = cfld( ci, ck, cj ) + &
282 1./REAL(nri ) * nfld( ni+ipoints , nk , nj+jpoints )
284 ! cfld( ci, ck, cj ) = 1./3. * &
285 ! ( nfld( ni , nk , nj-1) + &
286 ! nfld( ni , nk , nj ) + &
287 ! nfld( ni , nk , nj+1) )
292 ELSE IF ( ( .NOT. xstag ) .AND. ( ystag ) ) THEN
293 DO cj = MAX(jpos+spec_zone,cjts),MIN(jpos+(njde-njds)/nrj-jstag-spec_zone,cjte)
294 nj = (cj-jpos)*nrj + jstag + 1
297 DO ci = MAX(ipos+spec_zone,cits),MIN(ipos+(nide-nids)/nri-istag-spec_zone,cite)
298 ni = (ci-ipos)*nri + istag + 1
299 cfld( ci, ck, cj ) = 0.
300 DO ijpoints = ( nrj*nrj +1 )/2 - nrj/2 , ( nrj*nrj +1 )/2 - nrj/2 + nrj-1
301 ipoints = MOD((ijpoints-1),nri) + 1 - nri/2 - 1
302 jpoints = (ijpoints-1)/nri + 1 - nrj/2 - 1
303 cfld( ci, ck, cj ) = cfld( ci, ck, cj ) + &
304 1./REAL( nrj) * nfld( ni+ipoints , nk , nj+jpoints )
306 ! cfld( ci, ck, cj ) = 1./3. * &
307 ! ( nfld( ni-1, nk , nj ) + &
308 ! nfld( ni , nk , nj ) + &
309 ! nfld( ni+1, nk , nj ) )
316 ! Even refinement ratio
318 ELSE IF ( MOD(nrj,2) .EQ. 0) THEN
319 IF ( ( .NOT. xstag ) .AND. ( .NOT. ystag ) ) THEN
321 ! This is a simple schematic of the feedback indexing used in the even
322 ! ratio nests. For simplicity, a 2::1 ratio is depicted. Only the
323 ! mass variable staggering is shown.
325 ! the boxes with a "T" and four small "t" represents a coarse grid (CG)
326 ! cell, that is composed of four (2::1 ratio) fine grid (FG) cells.
328 ! Shown below is the area of the CG that is in the area of the FG. The
329 ! first grid point of the depicted CG is the starting location of the nest
330 ! in the parent domain (ipos,jpos - i_parent_start and j_parent_start from
333 ! For each of the CG points, the feedback loop is over each of the FG points
334 ! within the CG cell. For a 2::1 ratio, there are four total points (this is
335 ! the ijpoints loop). The feedback value to the CG is the arithmetic mean of
336 ! all of the FG values within each CG cell.
338 ! |-------------||-------------| |-------------||-------------|
339 ! | t t || t t | | t t || t t |
340 ! jpos+ | || | | || |
341 ! (njde-njds)- | T || T | | T || T |
342 ! jstag | || | | || |
343 ! | t t || t t | | t t || t t |
344 ! |-------------||-------------| |-------------||-------------|
345 ! |-------------||-------------| |-------------||-------------|
346 ! | t t || t t | | t t || t t |
348 ! | T || T | | T || T |
350 ! | t t || t t | | t t || t t |
351 ! |-------------||-------------| |-------------||-------------|
359 ! |-------------||-------------| |-------------||-------------|
360 ! jpoints = 1 | t t || t t | | t t || t t |
362 ! | T || T | | T || T |
364 ! jpoints = 0, | t t || t t | | t t || t t |
365 ! nj=3 |-------------||-------------| |-------------||-------------|
366 ! |-------------||-------------| |-------------||-------------|
367 ! jpoints = 1 | t t || t t | | t t || t t |
369 ! jpos | T || T | ... | T || T |
371 ! jpoints = 0, | t t || t t | ... | t t || t t |
372 ! nj=1 |-------------||-------------| |-------------||-------------|
377 ! ni = 1 3 (nide-nids)/nri
378 ! ipoints= 0 1 0 1 -istag
381 ! For performance benefits, users can comment out the inner most loop (and cfld=0) and
382 ! hardcode the loop feedback. For example, it is set up to run a 2::1 ratio
383 ! if uncommented. This lacks generality, but is likely to gain timing benefits
384 ! with compilers unable to unroll inner loops that do not have parameterized sizes.
386 ! The extra +1 ---------/ and the extra -1 ----\ (both for ci and cj)
387 ! / \ keeps the feedback out of the
388 ! / \ outer row/col, since that CG data
389 ! / \ specified the nest boundary originally
392 ! / \ a sentence to not end a line
393 ! / \ with a stupid backslash
394 DO cj = MAX(jpos+spec_zone,cjts),MIN(jpos+(njde-njds)/nrj-jstag-spec_zone,cjte)
395 nj = (cj-jpos)*nrj + jstag
398 DO ci = MAX(ipos+spec_zone,cits),MIN(ipos+(nide-nids)/nri-istag-spec_zone,cite)
399 ni = (ci-ipos)*nri + istag
400 cfld( ci, ck, cj ) = 0.
401 DO ijpoints = 1 , nri * nrj
402 ipoints = MOD((ijpoints-1),nri)
403 jpoints = (ijpoints-1)/nri
404 cfld( ci, ck, cj ) = cfld( ci, ck, cj ) + &
405 1./REAL(nri*nrj) * nfld( ni+ipoints , nk , nj+jpoints )
407 ! cfld( ci, ck, cj ) = 1./4. * &
408 ! ( nfld( ni , nk , nj ) + &
409 ! nfld( ni+1, nk , nj ) + &
410 ! nfld( ni , nk , nj+1) + &
411 ! nfld( ni+1, nk , nj+1) )
418 ELSE IF ( ( xstag ) .AND. ( .NOT. ystag ) ) THEN
446 DO cj = MAX(jpos+spec_zone,cjts),MIN(jpos+(njde-njds)/nrj-jstag-spec_zone,cjte)
447 nj = (cj-jpos)*nrj + 1
450 DO ci = MAX(ipos+spec_zone,cits),MIN(ipos+(nide-nids)/nri-istag-spec_zone,cite)
451 ni = (ci-ipos)*nri + 1
452 cfld( ci, ck, cj ) = 0.
453 DO ijpoints = 1 , nri*nrj , nri
454 ipoints = MOD((ijpoints-1),nri)
455 jpoints = (ijpoints-1)/nri
456 cfld( ci, ck, cj ) = cfld( ci, ck, cj ) + &
457 1./REAL(nri ) * nfld( ni+ipoints , nk , nj+jpoints )
459 ! cfld( ci, ck, cj ) = 1./2. * &
460 ! ( nfld( ni , nk , nj ) + &
461 ! nfld( ni , nk , nj+1) )
468 ELSE IF ( ( .NOT. xstag ) .AND. ( ystag ) ) THEN
469 DO cj = MAX(jpos+spec_zone,cjts),MIN(jpos+(njde-njds)/nrj-jstag-spec_zone,cjte)
470 nj = (cj-jpos)*nrj + 1
473 DO ci = MAX(ipos+spec_zone,cits),MIN(ipos+(nide-nids)/nri-istag-spec_zone,cite)
474 ni = (ci-ipos)*nri + 1
475 cfld( ci, ck, cj ) = 0.
476 DO ijpoints = 1 , nri
477 ipoints = MOD((ijpoints-1),nri)
478 jpoints = (ijpoints-1)/nri
479 cfld( ci, ck, cj ) = cfld( ci, ck, cj ) + &
480 1./REAL(nri ) * nfld( ni+ipoints , nk , nj+jpoints )
482 ! cfld( ci, ck, cj ) = 1./2. * &
483 ! ( nfld( ni , nk , nj ) + &
484 ! nfld( ni+1, nk , nj ) )
493 END SUBROUTINE copy_fcn
495 !==================================
496 ! this is the 1pt function used in feedback.
498 SUBROUTINE copy_fcnm ( cfld, & ! CD field
499 cids, cide, ckds, ckde, cjds, cjde, &
500 cims, cime, ckms, ckme, cjms, cjme, &
501 cits, cite, ckts, ckte, cjts, cjte, &
503 nids, nide, nkds, nkde, njds, njde, &
504 nims, nime, nkms, nkme, njms, njme, &
505 nits, nite, nkts, nkte, njts, njte, &
506 shw, & ! stencil half width for interp
507 imask, & ! interpolation mask
508 xstag, ystag, & ! staggering of field
509 ipos, jpos, & ! Position of lower left of nest in CD
510 nri, nrj ) ! nest ratios
516 INTEGER, INTENT(IN) :: cids, cide, ckds, ckde, cjds, cjde, &
517 cims, cime, ckms, ckme, cjms, cjme, &
518 cits, cite, ckts, ckte, cjts, cjte, &
519 nids, nide, nkds, nkde, njds, njde, &
520 nims, nime, nkms, nkme, njms, njme, &
521 nits, nite, nkts, nkte, njts, njte, &
525 LOGICAL, INTENT(IN) :: xstag, ystag
527 REAL, DIMENSION ( cims:cime, ckms:ckme, cjms:cjme ), INTENT(OUT) :: cfld
528 REAL, DIMENSION ( nims:nime, nkms:nkme, njms:njme ), INTENT(IN) :: nfld
529 INTEGER, DIMENSION ( nims:nime, njms:njme ), INTENT(IN) :: imask
533 INTEGER ci, cj, ck, ni, nj, nk, ip, jp, ioff, joff, ioffa, joffa
534 INTEGER :: icmin,icmax,jcmin,jcmax
535 INTEGER :: istag,jstag, ipoints,jpoints,ijpoints
536 INTEGER , PARAMETER :: passes = 2
539 CALL nl_get_spec_zone( 1, spec_zone )
540 istag = 1 ; jstag = 1
541 IF ( xstag ) istag = 0
542 IF ( ystag ) jstag = 0
544 IF( MOD(nrj,2) .NE. 0) THEN ! odd refinement ratio
546 DO cj = MAX(jpos+spec_zone,cjts),MIN(jpos+(njde-njds)/nrj-jstag-spec_zone,cjte)
547 nj = (cj-jpos)*nrj + jstag + 1
550 DO ci = MAX(ipos+spec_zone,cits),MIN(ipos+(nide-nids)/nri-istag-spec_zone,cite)
551 ni = (ci-ipos)*nri + istag + 1
552 cfld( ci, ck, cj ) = nfld( ni , nk , nj )
557 ELSE ! even refinement ratio, pick nearest neighbor on SW corner
558 DO cj = MAX(jpos+spec_zone,cjts),MIN(jpos+(njde-njds)/nrj-jstag-spec_zone,cjte)
559 nj = (cj-jpos)*nrj + 1
562 DO ci = MAX(ipos+spec_zone,cits),MIN(ipos+(nide-nids)/nri-istag-spec_zone,cite)
563 ni = (ci-ipos)*nri + 1
566 cfld( ci, ck, cj ) = nfld( ni+ipoints , nk , nj+jpoints )
575 END SUBROUTINE copy_fcnm
577 !==================================
578 ! this is the 1pt function used in feedback for integers
580 SUBROUTINE copy_fcni ( cfld, & ! CD field
581 cids, cide, ckds, ckde, cjds, cjde, &
582 cims, cime, ckms, ckme, cjms, cjme, &
583 cits, cite, ckts, ckte, cjts, cjte, &
585 nids, nide, nkds, nkde, njds, njde, &
586 nims, nime, nkms, nkme, njms, njme, &
587 nits, nite, nkts, nkte, njts, njte, &
588 shw, & ! stencil half width for interp
589 imask, & ! interpolation mask
590 xstag, ystag, & ! staggering of field
591 ipos, jpos, & ! Position of lower left of nest in CD
592 nri, nrj ) ! nest ratios
598 INTEGER, INTENT(IN) :: cids, cide, ckds, ckde, cjds, cjde, &
599 cims, cime, ckms, ckme, cjms, cjme, &
600 cits, cite, ckts, ckte, cjts, cjte, &
601 nids, nide, nkds, nkde, njds, njde, &
602 nims, nime, nkms, nkme, njms, njme, &
603 nits, nite, nkts, nkte, njts, njte, &
607 LOGICAL, INTENT(IN) :: xstag, ystag
609 INTEGER, DIMENSION ( cims:cime, ckms:ckme, cjms:cjme ), INTENT(OUT) :: cfld
610 INTEGER, DIMENSION ( nims:nime, nkms:nkme, njms:njme ), INTENT(IN) :: nfld
611 INTEGER, DIMENSION ( nims:nime, njms:njme ), INTENT(IN) :: imask
615 INTEGER ci, cj, ck, ni, nj, nk, ip, jp, ioff, joff, ioffa, joffa
616 INTEGER :: icmin,icmax,jcmin,jcmax
617 INTEGER :: istag,jstag, ipoints,jpoints,ijpoints
618 INTEGER , PARAMETER :: passes = 2
621 CALL nl_get_spec_zone( 1, spec_zone )
622 istag = 1 ; jstag = 1
623 IF ( xstag ) istag = 0
624 IF ( ystag ) jstag = 0
626 IF( MOD(nrj,2) .NE. 0) THEN ! odd refinement ratio
628 DO cj = MAX(jpos+spec_zone,cjts),MIN(jpos+(njde-njds)/nrj-jstag-spec_zone,cjte)
629 nj = (cj-jpos)*nrj + jstag + 1
632 DO ci = MAX(ipos+spec_zone,cits),MIN(ipos+(nide-nids)/nri-istag-spec_zone,cite)
633 ni = (ci-ipos)*nri + istag + 1
634 cfld( ci, ck, cj ) = nfld( ni , nk , nj )
639 ELSE ! even refinement ratio
640 DO cj = MAX(jpos+spec_zone,cjts),MIN(jpos+(njde-njds)/nrj-jstag-spec_zone,cjte)
641 nj = (cj-jpos)*nrj + 1
644 DO ci = MAX(ipos+spec_zone,cits),MIN(ipos+(nide-nids)/nri-istag-spec_zone,cite)
645 ni = (ci-ipos)*nri + 1
648 cfld( ci, ck, cj ) = nfld( ni+ipoints , nk , nj+jpoints )
657 END SUBROUTINE copy_fcni
659 !==================================
661 SUBROUTINE bdy_interp ( cfld, & ! CD field
662 cids, cide, ckds, ckde, cjds, cjde, &
663 cims, cime, ckms, ckme, cjms, cjme, &
664 cits, cite, ckts, ckte, cjts, cjte, &
666 nids, nide, nkds, nkde, njds, njde, &
667 nims, nime, nkms, nkme, njms, njme, &
668 nits, nite, nkts, nkte, njts, njte, &
669 shw, & ! stencil half width
670 imask, & ! interpolation mask
671 xstag, ystag, & ! staggering of field
672 ipos, jpos, & ! Position of lower left of nest in CD
673 nri, nrj, & ! nest ratios
678 cbdy_txs, nbdy_txs, &
679 cbdy_txe, nbdy_txe, &
680 cbdy_tys, nbdy_tys, &
681 cbdy_tye, nbdy_tye, &
687 INTEGER, INTENT(IN) :: cids, cide, ckds, ckde, cjds, cjde, &
688 cims, cime, ckms, ckme, cjms, cjme, &
689 cits, cite, ckts, ckte, cjts, cjte, &
690 nids, nide, nkds, nkde, njds, njde, &
691 nims, nime, nkms, nkme, njms, njme, &
692 nits, nite, nkts, nkte, njts, njte, &
697 LOGICAL, INTENT(IN) :: xstag, ystag
699 REAL, DIMENSION ( cims:cime, ckms:ckme, cjms:cjme ) :: cfld
700 REAL, DIMENSION ( nims:nime, nkms:nkme, njms:njme ) :: nfld
701 INTEGER, DIMENSION ( nims:nime, njms:njme ) :: imask
702 REAL, DIMENSION( * ), INTENT(INOUT) :: cbdy_xs, cbdy_txs, nbdy_xs, nbdy_txs
703 REAL, DIMENSION( * ), INTENT(INOUT) :: cbdy_xe, cbdy_txe, nbdy_xe, nbdy_txe
704 REAL, DIMENSION( * ), INTENT(INOUT) :: cbdy_ys, cbdy_tys, nbdy_ys, nbdy_tys
705 REAL, DIMENSION( * ), INTENT(INOUT) :: cbdy_ye, cbdy_tye, nbdy_ye, nbdy_tye
710 INTEGER nijds, nijde, spec_bdy_width
712 nijds = min(nids, njds)
713 nijde = max(nide, njde)
714 CALL nl_get_spec_bdy_width( 1, spec_bdy_width )
716 CALL bdy_interp1( cfld, & ! CD field
717 cids, cide, ckds, ckde, cjds, cjde, &
718 cims, cime, ckms, ckme, cjms, cjme, &
719 cits, cite, ckts, ckte, cjts, cjte, &
721 nijds, nijde , spec_bdy_width , &
722 nids, nide, nkds, nkde, njds, njde, &
723 nims, nime, nkms, nkme, njms, njme, &
724 nits, nite, nkts, nkte, njts, njte, &
726 xstag, ystag, & ! staggering of field
727 ipos, jpos, & ! Position of lower left of nest in CD
733 cbdy_txs, nbdy_txs, &
734 cbdy_txe, nbdy_txe, &
735 cbdy_tys, nbdy_tys, &
736 cbdy_tye, nbdy_tye, &
742 END SUBROUTINE bdy_interp
744 SUBROUTINE bdy_interp1( cfld, & ! CD field
745 cids, cide, ckds, ckde, cjds, cjde, &
746 cims, cime, ckms, ckme, cjms, cjme, &
747 cits, cite, ckts, ckte, cjts, cjte, &
749 nijds, nijde, spec_bdy_width , &
750 nids, nide, nkds, nkde, njds, njde, &
751 nims, nime, nkms, nkme, njms, njme, &
752 nits, nite, nkts, nkte, njts, njte, &
754 imask, & ! interpolation mask
755 xstag, ystag, & ! staggering of field
756 ipos, jpos, & ! Position of lower left of nest in CD
770 use module_state_description
773 INTEGER, INTENT(IN) :: cids, cide, ckds, ckde, cjds, cjde, &
774 cims, cime, ckms, ckme, cjms, cjme, &
775 cits, cite, ckts, ckte, cjts, cjte, &
776 nids, nide, nkds, nkde, njds, njde, &
777 nims, nime, nkms, nkme, njms, njme, &
778 nits, nite, nkts, nkte, njts, njte, &
782 INTEGER, INTENT(IN) :: nijds, nijde, spec_bdy_width
783 LOGICAL, INTENT(IN) :: xstag, ystag
785 REAL, DIMENSION ( cims:cime, ckms:ckme, cjms:cjme ), INTENT(INOUT) :: cfld
786 REAL, DIMENSION ( nims:nime, nkms:nkme, njms:njme ), INTENT(INOUT) :: nfld
787 INTEGER, DIMENSION ( nims:nime, njms:njme ) :: imask
788 REAL, DIMENSION ( * ), INTENT(INOUT) :: cbdy_xs, cbdy_txs ! not used
789 REAL, DIMENSION ( * ), INTENT(INOUT) :: cbdy_xe, cbdy_txe ! not used
790 REAL, DIMENSION ( * ), INTENT(INOUT) :: cbdy_ys, cbdy_tys ! not used
791 REAL, DIMENSION ( * ), INTENT(INOUT) :: cbdy_ye, cbdy_tye ! not used
793 REAL, DIMENSION ( njms:njme, nkms:nkme, spec_bdy_width ), INTENT(INOUT) :: bdy_xs, bdy_txs
794 REAL, DIMENSION ( njms:njme, nkms:nkme, spec_bdy_width ), INTENT(INOUT) :: bdy_xe, bdy_txe
795 REAL, DIMENSION ( nims:nime, nkms:nkme, spec_bdy_width ), INTENT(INOUT) :: bdy_ys, bdy_tys
796 REAL, DIMENSION ( nims:nime, nkms:nkme, spec_bdy_width ), INTENT(INOUT) :: bdy_ye, bdy_tye
801 INTEGER ci, cj, ck, ni, nj, nk, ni1, nj1, nk1, ip, jp, ioff, joff
806 REAL psca1(cims:cime,cjms:cjme,nri*nrj)
807 REAL psca(cims:cime,cjms:cjme,nri*nrj)
808 LOGICAL icmask( cims:cime, cjms:cjme )
818 ! statement functions for converting a nest index to coarse
819 n2ci(n) = (n+ipos*nri-1)/nri
820 n2cj(n) = (n+jpos*nrj-1)/nrj
834 ! Iterate over the ND tile and compute the values
838 CALL nl_get_spec_zone( 1, spec_zone )
839 CALL nl_get_relax_zone( 1, relax_zone )
840 sz = MIN(MAX( spec_zone, relax_zone + 1 ),spec_bdy_width)
845 !$OMP PRIVATE ( i,j,k,ni,nj,ni1,nj1,ci,cj,ip,jp,nk,ck,nf,icmask,psca,psca1 )
850 nj = (j-jpos) * nrj + ( nrj / 2 + 1 ) ! j point on nest
852 ni = (i-ipos) * nri + ( nri / 2 + 1 ) ! i point on nest
853 psca1(i,j,nf) = cfld(i,k,j)
857 ! hopefully less ham handed but still correct and more efficient
858 ! sintb ignores icmask so it does not matter that icmask is not set
861 IF ( njts .ge. njds .and. njts .le. njds + sz + joff ) THEN
862 CALL sintb( psca1, psca, &
863 cims, cime, cjms, cjme, icmask, &
864 n2ci(nits)-1, n2ci(nite)+1, n2cj(MAX(njts,njds)), n2cj(MIN(njte,njds+sz+joff)), nrj*nri, xstag, ystag )
867 IF ( njte .le. njde .and. njte .ge. njde - sz - joff ) THEN
868 CALL sintb( psca1, psca, &
869 cims, cime, cjms, cjme, icmask, &
870 n2ci(nits)-1, n2ci(nite)+1, n2cj(MAX(njts,njde-sz-joff)), n2cj(MIN(njte,njde-1+joff)), nrj*nri, xstag, ystag )
873 IF ( nits .ge. nids .and. nits .le. nids + sz + ioff ) THEN
874 CALL sintb( psca1, psca, &
875 cims, cime, cjms, cjme, icmask, &
876 n2ci(MAX(nits,nids)), n2ci(MIN(nite,nids+sz+ioff)), n2cj(njts)-1, n2cj(njte)+1, nrj*nri, xstag, ystag )
879 IF ( nite .le. nide .and. nite .ge. nide - sz - ioff ) THEN
880 CALL sintb( psca1, psca, &
881 cims, cime, cjms, cjme, icmask, &
882 n2ci(MAX(nits,nide-sz-ioff)), n2ci(MIN(nite,nide-1+ioff)), n2cj(njts)-1, n2cj(njte)+1, nrj*nri, xstag, ystag )
885 DO nj1 = MAX(njds,njts-1), MIN(njde+joff,njte+joff+1)
886 cj = jpos + (nj1-1) / nrj ! j coord of CD point
887 jp = mod ( nj1-1 , nrj ) ! coord of ND w/i CD point
890 DO ni1 = MAX(nids,nits-1), MIN(nide+ioff,nite+ioff+1)
891 ci = ipos + (ni1-1) / nri ! j coord of CD point
892 ip = mod ( ni1-1 , nri ) ! coord of ND w/i CD point
897 IF ( ( ni.LT.nids) .OR. (nj.LT.njds) ) THEN
901 !bdy contains the value at t-dt. psca contains the value at t
902 !compute dv/dt and store in bdy_t
903 !afterwards store the new value of v at t into bdy
905 IF ( ni .ge. nids .and. ni .lt. nids + sz ) THEN
906 bdy_txs( nj,k,ni ) = rdt*(psca(ci+shw,cj+shw,ip+1+(jp)*nri)-nfld(ni,k,nj))
907 bdy_xs( nj,k,ni ) = psca(ci+shw,cj+shw,ip+1+(jp)*nri )
911 IF ( nj .ge. njds .and. nj .lt. njds + sz ) THEN
912 bdy_tys( ni,k,nj ) = rdt*(psca(ci+shw,cj+shw,ip+1+(jp)*nri)-nfld(ni,k,nj))
913 bdy_ys( ni,k,nj ) = psca(ci+shw,cj+shw,ip+1+(jp)*nri )
918 IF ( ni .ge. nide - sz + 1 .AND. ni .le. nide ) THEN
919 bdy_txe( nj,k,nide-ni+1 ) = rdt*(psca(ci+shw,cj+shw,ip+1+(jp)*nri)-nfld(ni,k,nj))
920 bdy_xe( nj,k,nide-ni+1 ) = psca(ci+shw,cj+shw,ip+1+(jp)*nri )
923 IF ( ni .ge. nide - sz .AND. ni .le. nide-1 ) THEN
924 bdy_txe( nj,k,nide-ni ) = rdt*(psca(ci+shw,cj+shw,ip+1+(jp)*nri)-nfld(ni,k,nj))
925 bdy_xe( nj,k,nide-ni ) = psca(ci+shw,cj+shw,ip+1+(jp)*nri )
931 IF ( nj .ge. njde - sz + 1 .AND. nj .le. njde ) THEN
932 bdy_tye( ni,k,njde-nj+1 ) = rdt*(psca(ci+shw,cj+shw,ip+1+(jp)*nri)-nfld(ni,k,nj))
933 bdy_ye( ni,k,njde-nj+1 ) = psca(ci+shw,cj+shw,ip+1+(jp)*nri )
936 IF ( nj .ge. njde - sz .AND. nj .le. njde-1 ) THEN
937 bdy_tye(ni,k,njde-nj ) = rdt*(psca(ci+shw,cj+shw,ip+1+(jp)*nri)-nfld(ni,k,nj))
938 bdy_ye( ni,k,njde-nj ) = psca(ci+shw,cj+shw,ip+1+(jp)*nri )
945 !$OMP END PARALLEL DO
950 END SUBROUTINE bdy_interp1
954 SUBROUTINE interp_fcni( cfld, & ! CD field
955 cids, cide, ckds, ckde, cjds, cjde, &
956 cims, cime, ckms, ckme, cjms, cjme, &
957 cits, cite, ckts, ckte, cjts, cjte, &
959 nids, nide, nkds, nkde, njds, njde, &
960 nims, nime, nkms, nkme, njms, njme, &
961 nits, nite, nkts, nkte, njts, njte, &
962 shw, & ! stencil half width
963 imask, & ! interpolation mask
964 xstag, ystag, & ! staggering of field
965 ipos, jpos, & ! Position of lower left of nest in CD
966 nri, nrj ) ! nest ratios
971 INTEGER, INTENT(IN) :: cids, cide, ckds, ckde, cjds, cjde, &
972 cims, cime, ckms, ckme, cjms, cjme, &
973 cits, cite, ckts, ckte, cjts, cjte, &
974 nids, nide, nkds, nkde, njds, njde, &
975 nims, nime, nkms, nkme, njms, njme, &
976 nits, nite, nkts, nkte, njts, njte, &
980 LOGICAL, INTENT(IN) :: xstag, ystag
982 INTEGER, DIMENSION ( cims:cime, ckms:ckme, cjms:cjme ) :: cfld
983 INTEGER, DIMENSION ( nims:nime, nkms:nkme, njms:njme ) :: nfld
984 INTEGER, DIMENSION ( nims:nime, njms:njme ) :: imask
988 INTEGER ci, cj, ck, ni, nj, nk, ip, jp
990 ! Iterate over the ND tile and compute the values
993 !write(0,'("cits:cite, ckts:ckte, cjts:cjte ",6i4)')cits,cite, ckts,ckte, cjts,cjte
994 !write(0,'("nits:nite, nkts:nkte, njts:njte ",6i4)')nits,nite, nkts,nkte, njts,njte
997 cj = jpos + (nj-1) / nrj ! j coord of CD point
998 jp = mod ( nj , nrj ) ! coord of ND w/i CD point
1002 ci = ipos + (ni-1) / nri ! j coord of CD point
1003 ip = mod ( ni , nri ) ! coord of ND w/i CD point
1004 ! This is a trivial implementation of the interp_fcn; just copies
1005 ! the values from the CD into the ND
1006 nfld( ni, nk, nj ) = cfld( ci , ck , cj )
1013 END SUBROUTINE interp_fcni
1015 SUBROUTINE interp_fcnm( cfld, & ! CD field
1016 cids, cide, ckds, ckde, cjds, cjde, &
1017 cims, cime, ckms, ckme, cjms, cjme, &
1018 cits, cite, ckts, ckte, cjts, cjte, &
1020 nids, nide, nkds, nkde, njds, njde, &
1021 nims, nime, nkms, nkme, njms, njme, &
1022 nits, nite, nkts, nkte, njts, njte, &
1023 shw, & ! stencil half width
1024 imask, & ! interpolation mask
1025 xstag, ystag, & ! staggering of field
1026 ipos, jpos, & ! Position of lower left of nest in CD
1027 nri, nrj ) ! nest ratios
1028 USE module_configure
1032 INTEGER, INTENT(IN) :: cids, cide, ckds, ckde, cjds, cjde, &
1033 cims, cime, ckms, ckme, cjms, cjme, &
1034 cits, cite, ckts, ckte, cjts, cjte, &
1035 nids, nide, nkds, nkde, njds, njde, &
1036 nims, nime, nkms, nkme, njms, njme, &
1037 nits, nite, nkts, nkte, njts, njte, &
1041 LOGICAL, INTENT(IN) :: xstag, ystag
1043 REAL, DIMENSION ( cims:cime, ckms:ckme, cjms:cjme ) :: cfld
1044 REAL, DIMENSION ( nims:nime, nkms:nkme, njms:njme ) :: nfld
1045 INTEGER, DIMENSION ( nims:nime, njms:njme ) :: imask
1049 INTEGER ci, cj, ck, ni, nj, nk, ip, jp
1051 ! Iterate over the ND tile and compute the values
1054 !write(0,'("mask cits:cite, ckts:ckte, cjts:cjte ",6i4)')cits,cite, ckts,ckte, cjts,cjte
1055 !write(0,'("mask nits:nite, nkts:nkte, njts:njte ",6i4)')nits,nite, nkts,nkte, njts,njte
1058 cj = jpos + (nj-1) / nrj ! j coord of CD point
1059 jp = mod ( nj , nrj ) ! coord of ND w/i CD point
1063 ci = ipos + (ni-1) / nri ! j coord of CD point
1064 ip = mod ( ni , nri ) ! coord of ND w/i CD point
1065 ! This is a trivial implementation of the interp_fcn; just copies
1066 ! the values from the CD into the ND
1067 nfld( ni, nk, nj ) = cfld( ci , ck , cj )
1074 END SUBROUTINE interp_fcnm
1076 SUBROUTINE interp_mask_land_field ( enable, & ! says whether to allow interpolation or just the bcasts
1078 cids, cide, ckds, ckde, cjds, cjde, &
1079 cims, cime, ckms, ckme, cjms, cjme, &
1080 cits, cite, ckts, ckte, cjts, cjte, &
1082 nids, nide, nkds, nkde, njds, njde, &
1083 nims, nime, nkms, nkme, njms, njme, &
1084 nits, nite, nkts, nkte, njts, njte, &
1085 shw, & ! stencil half width
1086 imask, & ! interpolation mask
1087 xstag, ystag, & ! staggering of field
1088 ipos, jpos, & ! Position of lower left of nest in CD
1089 nri, nrj, & ! nest ratios
1092 USE module_configure
1093 USE module_wrf_error
1098 LOGICAL, INTENT(IN) :: enable
1099 INTEGER, INTENT(IN) :: cids, cide, ckds, ckde, cjds, cjde, &
1100 cims, cime, ckms, ckme, cjms, cjme, &
1101 cits, cite, ckts, ckte, cjts, cjte, &
1102 nids, nide, nkds, nkde, njds, njde, &
1103 nims, nime, nkms, nkme, njms, njme, &
1104 nits, nite, nkts, nkte, njts, njte, &
1108 LOGICAL, INTENT(IN) :: xstag, ystag
1110 REAL, DIMENSION ( cims:cime, ckms:ckme, cjms:cjme ) :: cfld
1111 REAL, DIMENSION ( nims:nime, nkms:nkme, njms:njme ) :: nfld
1112 INTEGER, DIMENSION ( nims:nime, njms:njme ) :: imask
1114 REAL, DIMENSION ( cims:cime, cjms:cjme ) :: clu
1115 REAL, DIMENSION ( nims:nime, njms:njme ) :: nlu
1119 INTEGER ci, cj, ck, ni, nj, nk, ip, jp
1120 INTEGER :: icount , ii , jj , ist , ien , jst , jen , iswater
1121 REAL :: avg , sum , dx , dy
1122 INTEGER , PARAMETER :: max_search = 5
1123 CHARACTER*120 message
1125 ! Find out what the water value is.
1127 CALL nl_get_iswater(1,iswater)
1129 ! Right now, only mass point locations permitted.
1131 IF ( ( .NOT. xstag) .AND. ( .NOT. ystag ) ) THEN
1133 ! Loop over each i,k,j in the nested domain.
1138 IF ( MOD ( nrj , 2 ) .EQ. 0 ) THEN
1139 cj = ( nj + (nrj/2)-1 ) / nrj + jpos -1 ! first coarse position equal to or below nest point
1141 cj = ( nj + (nrj-1)/2 ) / nrj + jpos -1 ! first coarse position equal to or below nest point
1146 IF ( MOD ( nri , 2 ) .EQ. 0 ) THEN
1147 ci = ( ni + (nri/2)-1 ) / nri + ipos -1 ! first coarse position equal to or to the left of nest point
1149 ci = ( ni + (nri-1)/2 ) / nri + ipos -1 ! first coarse position equal to or to the left of nest point
1156 ! (ci,cj+1) (ci+1,cj+1)
1170 ! For odd ratios, at (nri+1)/2, we are on the coarse grid point, so dx = 0
1172 IF ( MOD ( nri , 2 ) .EQ. 0 ) THEN
1173 dx = ( REAL ( MOD ( ni+(nri-1)/2 , nri ) ) + 0.5 ) / REAL ( nri )
1175 dx = REAL ( MOD ( ni+(nri-1)/2 , nri ) ) / REAL ( nri )
1177 IF ( MOD ( nrj , 2 ) .EQ. 0 ) THEN
1178 dy = ( REAL ( MOD ( nj+(nrj-1)/2 , nrj ) ) + 0.5 ) / REAL ( nrj )
1180 dy = REAL ( MOD ( nj+(nrj-1)/2 , nrj ) ) / REAL ( nrj )
1183 ! This is a "land only" field. If this is a water point, no operations required.
1185 IF ( ( NINT(nlu(ni ,nj )) .EQ. iswater ) ) THEN
1187 ! nfld(ni,nk,nj) = 1.e20
1190 ! If this is a nested land point, and the surrounding coarse values are all land points,
1191 ! then this is a simple 4-pt interpolation.
1193 ELSE IF ( ( NINT(nlu(ni ,nj )) .NE. iswater ) .AND. &
1194 ( NINT(clu(ci ,cj )) .NE. iswater ) .AND. &
1195 ( NINT(clu(ci+1,cj )) .NE. iswater ) .AND. &
1196 ( NINT(clu(ci ,cj+1)) .NE. iswater ) .AND. &
1197 ( NINT(clu(ci+1,cj+1)) .NE. iswater ) ) THEN
1198 nfld(ni,nk,nj) = ( 1. - dx ) * ( ( 1. - dy ) * cfld(ci ,ck,cj ) + &
1199 dy * cfld(ci ,ck,cj+1) ) + &
1200 dx * ( ( 1. - dy ) * cfld(ci+1,ck,cj ) + &
1201 dy * cfld(ci+1,ck,cj+1) )
1203 ! If this is a nested land point and there are NO coarse land values surrounding,
1204 ! we temporarily punt.
1206 ELSE IF ( ( NINT(nlu(ni ,nj )) .NE. iswater ) .AND. &
1207 ( NINT(clu(ci ,cj )) .EQ. iswater ) .AND. &
1208 ( NINT(clu(ci+1,cj )) .EQ. iswater ) .AND. &
1209 ( NINT(clu(ci ,cj+1)) .EQ. iswater ) .AND. &
1210 ( NINT(clu(ci+1,cj+1)) .EQ. iswater ) ) THEN
1211 ! nfld(ni,nk,nj) = -1.e20
1214 ! If there are some water points and some land points, take an average.
1216 ELSE IF ( NINT(nlu(ni ,nj )) .NE. iswater ) THEN
1219 IF ( NINT(clu(ci ,cj )) .NE. iswater ) THEN
1221 sum = sum + cfld(ci ,ck,cj )
1223 IF ( NINT(clu(ci+1,cj )) .NE. iswater ) THEN
1225 sum = sum + cfld(ci+1,ck,cj )
1227 IF ( NINT(clu(ci ,cj+1)) .NE. iswater ) THEN
1229 sum = sum + cfld(ci ,ck,cj+1)
1231 IF ( NINT(clu(ci+1,cj+1)) .NE. iswater ) THEN
1233 sum = sum + cfld(ci+1,ck,cj+1)
1235 nfld(ni,nk,nj) = sum / REAL ( icount )
1242 ! Get an average of the whole domain for problem locations.
1249 IF ( ( nfld(ni,nk,nj) .GT. -1.e19 ) .AND. ( nfld(ni,nk,nj) .LT. 1.e19 ) ) THEN
1251 sum = sum + nfld(ni,nk,nj)
1260 CALL wrf_dm_bcast_real( sum, 1 )
1261 CALL wrf_dm_bcast_integer( icount, 1 )
1263 IF ( icount .GT. 0 ) THEN
1264 avg = sum / REAL ( icount )
1266 ! OK, if there were any of those island situations, we try to search a bit broader
1267 ! of an area in the coarse grid.
1272 IF ( nfld(ni,nk,nj) .LT. -1.e19 ) THEN
1273 IF ( MOD ( nrj , 2 ) .EQ. 0 ) THEN
1274 cj = ( nj + (nrj/2)-1 ) / nrj + jpos -1 ! first coarse position equal to or below nest point
1276 cj = ( nj + (nrj-1)/2 ) / nrj + jpos -1 ! first coarse position equal to or below nest point
1278 IF ( MOD ( nri , 2 ) .EQ. 0 ) THEN
1279 ci = ( ni + (nri/2)-1 ) / nri + ipos -1 ! first coarse position equal to or to the left of nest point
1281 ci = ( ni + (nri-1)/2 ) / nri + ipos -1 ! first coarse position equal to or to the left of nest point
1283 ist = MAX (ci-max_search,cits)
1284 ien = MIN (ci+max_search,cite,cide-1)
1285 jst = MAX (cj-max_search,cjts)
1286 jen = MIN (cj+max_search,cjte,cjde-1)
1291 IF ( NINT(clu(ii,jj)) .NE. iswater ) THEN
1293 sum = sum + cfld(ii,nk,jj)
1297 IF ( icount .GT. 0 ) THEN
1298 nfld(ni,nk,nj) = sum / REAL ( icount )
1300 ! CALL wrf_error_fatal ( "horizontal interp error - island" )
1301 write(message,*) 'horizontal interp error - island, using average ', avg
1302 CALL wrf_message ( message )
1303 nfld(ni,nk,nj) = avg
1312 CALL wrf_error_fatal ( "only unstaggered fields right now" )
1315 END SUBROUTINE interp_mask_land_field
1317 SUBROUTINE interp_mask_water_field ( enable, & ! says whether to allow interpolation or just the bcasts
1319 cids, cide, ckds, ckde, cjds, cjde, &
1320 cims, cime, ckms, ckme, cjms, cjme, &
1321 cits, cite, ckts, ckte, cjts, cjte, &
1323 nids, nide, nkds, nkde, njds, njde, &
1324 nims, nime, nkms, nkme, njms, njme, &
1325 nits, nite, nkts, nkte, njts, njte, &
1326 shw, & ! stencil half width
1327 imask, & ! interpolation mask
1328 xstag, ystag, & ! staggering of field
1329 ipos, jpos, & ! Position of lower left of nest in CD
1330 nri, nrj, & ! nest ratios
1333 USE module_configure
1334 USE module_wrf_error
1339 LOGICAL, INTENT(IN) :: enable
1340 INTEGER, INTENT(IN) :: cids, cide, ckds, ckde, cjds, cjde, &
1341 cims, cime, ckms, ckme, cjms, cjme, &
1342 cits, cite, ckts, ckte, cjts, cjte, &
1343 nids, nide, nkds, nkde, njds, njde, &
1344 nims, nime, nkms, nkme, njms, njme, &
1345 nits, nite, nkts, nkte, njts, njte, &
1349 LOGICAL, INTENT(IN) :: xstag, ystag
1351 REAL, DIMENSION ( cims:cime, ckms:ckme, cjms:cjme ) :: cfld
1352 REAL, DIMENSION ( nims:nime, nkms:nkme, njms:njme ) :: nfld
1353 INTEGER, DIMENSION ( nims:nime, njms:njme ) :: imask
1355 REAL, DIMENSION ( cims:cime, cjms:cjme ) :: clu
1356 REAL, DIMENSION ( nims:nime, njms:njme ) :: nlu
1360 INTEGER ci, cj, ck, ni, nj, nk, ip, jp
1361 INTEGER :: icount , ii , jj , ist , ien , jst , jen , iswater
1362 REAL :: avg , sum , dx , dy
1363 INTEGER , PARAMETER :: max_search = 5
1365 ! Find out what the water value is.
1367 CALL nl_get_iswater(1,iswater)
1369 ! Right now, only mass point locations permitted.
1371 IF ( ( .NOT. xstag) .AND. ( .NOT. ystag ) ) THEN
1374 ! Loop over each i,k,j in the nested domain.
1377 IF ( MOD ( nrj , 2 ) .EQ. 0 ) THEN
1378 cj = ( nj + (nrj/2)-1 ) / nrj + jpos -1 ! first coarse position equal to or below nest point
1380 cj = ( nj + (nrj-1)/2 ) / nrj + jpos -1 ! first coarse position equal to or below nest point
1385 IF ( MOD ( nri , 2 ) .EQ. 0 ) THEN
1386 ci = ( ni + (nri/2)-1 ) / nri + ipos -1 ! first coarse position equal to or to the left of nest point
1388 ci = ( ni + (nri-1)/2 ) / nri + ipos -1 ! first coarse position equal to or to the left of nest point
1395 ! (ci,cj+1) (ci+1,cj+1)
1409 ! At ni=2, we are on the coarse grid point, so dx = 0
1411 IF ( MOD ( nri , 2 ) .EQ. 0 ) THEN
1412 dx = ( REAL ( MOD ( ni+(nri-1)/2 , nri ) ) + 0.5 ) / REAL ( nri )
1414 dx = REAL ( MOD ( ni+(nri-1)/2 , nri ) ) / REAL ( nri )
1416 IF ( MOD ( nrj , 2 ) .EQ. 0 ) THEN
1417 dy = ( REAL ( MOD ( nj+(nrj-1)/2 , nrj ) ) + 0.5 ) / REAL ( nrj )
1419 dy = REAL ( MOD ( nj+(nrj-1)/2 , nrj ) ) / REAL ( nrj )
1422 ! This is a "water only" field. If this is a land point, no operations required.
1424 IF ( ( NINT(nlu(ni ,nj )) .NE. iswater ) ) THEN
1426 ! nfld(ni,nk,nj) = 1.e20
1429 ! If this is a nested water point, and the surrounding coarse values are all water points,
1430 ! then this is a simple 4-pt interpolation.
1432 ELSE IF ( ( NINT(nlu(ni ,nj )) .EQ. iswater ) .AND. &
1433 ( NINT(clu(ci ,cj )) .EQ. iswater ) .AND. &
1434 ( NINT(clu(ci+1,cj )) .EQ. iswater ) .AND. &
1435 ( NINT(clu(ci ,cj+1)) .EQ. iswater ) .AND. &
1436 ( NINT(clu(ci+1,cj+1)) .EQ. iswater ) ) THEN
1437 nfld(ni,nk,nj) = ( 1. - dx ) * ( ( 1. - dy ) * cfld(ci ,ck,cj ) + &
1438 dy * cfld(ci ,ck,cj+1) ) + &
1439 dx * ( ( 1. - dy ) * cfld(ci+1,ck,cj ) + &
1440 dy * cfld(ci+1,ck,cj+1) )
1442 ! If this is a nested water point and there are NO coarse water values surrounding,
1443 ! we temporarily punt.
1445 ELSE IF ( ( NINT(nlu(ni ,nj )) .EQ. iswater ) .AND. &
1446 ( NINT(clu(ci ,cj )) .NE. iswater ) .AND. &
1447 ( NINT(clu(ci+1,cj )) .NE. iswater ) .AND. &
1448 ( NINT(clu(ci ,cj+1)) .NE. iswater ) .AND. &
1449 ( NINT(clu(ci+1,cj+1)) .NE. iswater ) ) THEN
1450 ! nfld(ni,nk,nj) = -1.e20
1453 ! If there are some land points and some water points, take an average.
1455 ELSE IF ( NINT(nlu(ni ,nj )) .EQ. iswater ) THEN
1458 IF ( NINT(clu(ci ,cj )) .EQ. iswater ) THEN
1460 sum = sum + cfld(ci ,ck,cj )
1462 IF ( NINT(clu(ci+1,cj )) .EQ. iswater ) THEN
1464 sum = sum + cfld(ci+1,ck,cj )
1466 IF ( NINT(clu(ci ,cj+1)) .EQ. iswater ) THEN
1468 sum = sum + cfld(ci ,ck,cj+1)
1470 IF ( NINT(clu(ci+1,cj+1)) .EQ. iswater ) THEN
1472 sum = sum + cfld(ci+1,ck,cj+1)
1474 nfld(ni,nk,nj) = sum / REAL ( icount )
1480 ! Get an average of the whole domain for problem locations.
1487 IF ( ( nfld(ni,nk,nj) .GT. -1.e19 ) .AND. ( nfld(ni,nk,nj) .LT. 1.e19 ) ) THEN
1489 sum = sum + nfld(ni,nk,nj)
1498 CALL wrf_dm_bcast_real( sum, 1 )
1499 CALL wrf_dm_bcast_integer( icount, 1 )
1501 IF ( icount .NE. 0 ) THEN
1502 avg = sum / REAL ( icount )
1505 ! OK, if there were any of those lake situations, we try to search a bit broader
1506 ! of an area in the coarse grid.
1511 IF ( nfld(ni,nk,nj) .LT. -1.e19 ) THEN
1512 IF ( MOD ( nrj , 2 ) .EQ. 0 ) THEN
1513 cj = ( nj + (nrj/2)-1 ) / nrj + jpos -1 ! first coarse position equal to or below nest point
1515 cj = ( nj + (nrj-1)/2 ) / nrj + jpos -1 ! first coarse position equal to or below nest point
1517 IF ( MOD ( nri , 2 ) .EQ. 0 ) THEN
1518 ci = ( ni + (nri/2)-1 ) / nri + ipos -1 ! first coarse position equal to or to the left of nest point
1520 ci = ( ni + (nri-1)/2 ) / nri + ipos -1 ! first coarse position equal to or to the left of nest point
1522 ist = MAX (ci-max_search,cits)
1523 ien = MIN (ci+max_search,cite,cide-1)
1524 jst = MAX (cj-max_search,cjts)
1525 jen = MIN (cj+max_search,cjte,cjde-1)
1530 IF ( NINT(clu(ii,jj)) .EQ. iswater ) THEN
1532 sum = sum + cfld(ii,nk,jj)
1536 IF ( icount .GT. 0 ) THEN
1537 nfld(ni,nk,nj) = sum / REAL ( icount )
1539 ! CALL wrf_error_fatal ( "horizontal interp error - lake" )
1540 print *,'horizontal interp error - lake, using average ',avg
1541 nfld(ni,nk,nj) = avg
1550 CALL wrf_error_fatal ( "only unstaggered fields right now" )
1553 END SUBROUTINE interp_mask_water_field
1558 SUBROUTINE smoother ( cfld , &
1559 cids, cide, ckds, ckde, cjds, cjde, &
1560 cims, cime, ckms, ckme, cjms, cjme, &
1561 cits, cite, ckts, ckte, cjts, cjte, &
1562 nids, nide, nkds, nkde, njds, njde, &
1563 nims, nime, nkms, nkme, njms, njme, &
1564 nits, nite, nkts, nkte, njts, njte, &
1565 xstag, ystag, & ! staggering of field
1566 ipos, jpos, & ! Position of lower left of nest in
1570 USE module_configure
1573 INTEGER, INTENT(IN) :: cids, cide, ckds, ckde, cjds, cjde, &
1574 cims, cime, ckms, ckme, cjms, cjme, &
1575 cits, cite, ckts, ckte, cjts, cjte, &
1576 nids, nide, nkds, nkde, njds, njde, &
1577 nims, nime, nkms, nkme, njms, njme, &
1578 nits, nite, nkts, nkte, njts, njte, &
1581 LOGICAL, INTENT(IN) :: xstag, ystag
1582 INTEGER :: smooth_option, feedback , spec_zone
1584 REAL, DIMENSION ( cims:cime, ckms:ckme, cjms:cjme ) :: cfld
1586 ! If there is no feedback, there can be no smoothing.
1588 CALL nl_get_feedback ( 1, feedback )
1589 IF ( feedback == 0 ) RETURN
1590 CALL nl_get_spec_zone ( 1, spec_zone )
1592 ! These are the 2d smoothers used on the fedback data. These filters
1593 ! are run on the coarse grid data (after the nested info has been
1594 ! fedback). Only the area of the nest in the coarse grid is filtered.
1596 CALL nl_get_smooth_option ( 1, smooth_option )
1598 IF ( smooth_option == 0 ) THEN
1600 ELSE IF ( smooth_option == 1 ) THEN
1601 CALL sm121 ( cfld , &
1602 cids, cide, ckds, ckde, cjds, cjde, &
1603 cims, cime, ckms, ckme, cjms, cjme, &
1604 cits, cite, ckts, ckte, cjts, cjte, &
1605 xstag, ystag, & ! staggering of field
1606 nids, nide, nkds, nkde, njds, njde, &
1607 nims, nime, nkms, nkme, njms, njme, &
1608 nits, nite, nkts, nkte, njts, njte, &
1610 ipos, jpos & ! Position of lower left of nest in
1612 ELSE IF ( smooth_option == 2 ) THEN
1613 CALL smdsm ( cfld , &
1614 cids, cide, ckds, ckde, cjds, cjde, &
1615 cims, cime, ckms, ckme, cjms, cjme, &
1616 cits, cite, ckts, ckte, cjts, cjte, &
1617 xstag, ystag, & ! staggering of field
1618 nids, nide, nkds, nkde, njds, njde, &
1619 nims, nime, nkms, nkme, njms, njme, &
1620 nits, nite, nkts, nkte, njts, njte, &
1622 ipos, jpos & ! Position of lower left of nest in
1626 END SUBROUTINE smoother
1628 SUBROUTINE sm121 ( cfld , &
1629 cids, cide, ckds, ckde, cjds, cjde, &
1630 cims, cime, ckms, ckme, cjms, cjme, &
1631 cits, cite, ckts, ckte, cjts, cjte, &
1632 xstag, ystag, & ! staggering of field
1633 nids, nide, nkds, nkde, njds, njde, &
1634 nims, nime, nkms, nkme, njms, njme, &
1635 nits, nite, nkts, nkte, njts, njte, &
1637 ipos, jpos & ! Position of lower left of nest in
1640 USE module_configure
1643 INTEGER, INTENT(IN) :: cids, cide, ckds, ckde, cjds, cjde, &
1644 cims, cime, ckms, ckme, cjms, cjme, &
1645 cits, cite, ckts, ckte, cjts, cjte, &
1646 nids, nide, nkds, nkde, njds, njde, &
1647 nims, nime, nkms, nkme, njms, njme, &
1648 nits, nite, nkts, nkte, njts, njte, &
1651 LOGICAL, INTENT(IN) :: xstag, ystag
1653 REAL, DIMENSION ( cims:cime, ckms:ckme, cjms:cjme ) :: cfld
1654 REAL, DIMENSION ( cims:cime, cjms:cjme ) :: cfldnew
1656 INTEGER :: i , j , k , loop
1657 INTEGER :: istag,jstag
1659 INTEGER, PARAMETER :: smooth_passes = 1 ! More passes requires a larger stencil (currently 48 pt)
1661 istag = 1 ; jstag = 1
1662 IF ( xstag ) istag = 0
1663 IF ( ystag ) jstag = 0
1665 ! Simple 1-2-1 smoother.
1667 smoothing_passes : DO loop = 1 , smooth_passes
1671 ! Initialize dummy cfldnew
1673 DO i = MAX(ipos,cits-3) , MIN(ipos+(nide-nids)/nri,cite+3)
1674 DO j = MAX(jpos,cjts-3) , MIN(jpos+(njde-njds)/nrj,cjte+3)
1675 cfldnew(i,j) = cfld(i,k,j)
1679 ! 1-2-1 smoothing in the j direction first,
1681 DO i = MAX(ipos+2,cits-2) , MIN(ipos+(nide-nids)/nri-2-istag,cite+2)
1682 DO j = MAX(jpos+2,cjts-2) , MIN(jpos+(njde-njds)/nrj-2-jstag,cjte+2)
1683 cfldnew(i,j) = 0.25 * ( cfld(i,k,j+1) + 2.*cfld(i,k,j) + cfld(i,k,j-1) )
1687 ! then 1-2-1 smoothing in the i direction last
1689 DO j = MAX(jpos+2,cjts-2) , MIN(jpos+(njde-njds)/nrj-2-jstag,cjte+2)
1690 DO i = MAX(ipos+2,cits-2) , MIN(ipos+(nide-nids)/nri-2-istag,cite+2)
1691 cfld(i,k,j) = 0.25 * ( cfldnew(i+1,j) + 2.*cfldnew(i,j) + cfldnew(i-1,j) )
1697 END DO smoothing_passes
1699 END SUBROUTINE sm121
1701 SUBROUTINE smdsm ( cfld , &
1702 cids, cide, ckds, ckde, cjds, cjde, &
1703 cims, cime, ckms, ckme, cjms, cjme, &
1704 cits, cite, ckts, ckte, cjts, cjte, &
1705 xstag, ystag, & ! staggering of field
1706 nids, nide, nkds, nkde, njds, njde, &
1707 nims, nime, nkms, nkme, njms, njme, &
1708 nits, nite, nkts, nkte, njts, njte, &
1710 ipos, jpos & ! Position of lower left of nest in
1713 USE module_configure
1716 INTEGER, INTENT(IN) :: cids, cide, ckds, ckde, cjds, cjde, &
1717 cims, cime, ckms, ckme, cjms, cjme, &
1718 cits, cite, ckts, ckte, cjts, cjte, &
1719 nids, nide, nkds, nkde, njds, njde, &
1720 nims, nime, nkms, nkme, njms, njme, &
1721 nits, nite, nkts, nkte, njts, njte, &
1724 LOGICAL, INTENT(IN) :: xstag, ystag
1726 REAL, DIMENSION ( cims:cime, ckms:ckme, cjms:cjme ) :: cfld
1727 REAL, DIMENSION ( cims:cime, cjms:cjme ) :: cfldnew
1729 REAL , DIMENSION ( 2 ) :: xnu
1730 INTEGER :: i , j , k , loop , n
1731 INTEGER :: istag,jstag
1733 INTEGER, PARAMETER :: smooth_passes = 1 ! More passes requires a larger stencil (currently 48 pt)
1735 xnu = (/ 0.50 , -0.52 /)
1737 istag = 1 ; jstag = 1
1738 IF ( xstag ) istag = 0
1739 IF ( ystag ) jstag = 0
1741 ! The odd number passes of this are the "smoother", the even
1742 ! number passes are the "de-smoother" (note the different signs on xnu).
1744 smoothing_passes : DO loop = 1 , smooth_passes * 2
1746 n = 2 - MOD ( loop , 2 )
1750 DO i = MAX(ipos+2,cits-2) , MIN(ipos+(nide-nids)/nri-2-istag,cite+2)
1751 DO j = MAX(jpos+2,cjts-2) , MIN(jpos+(njde-njds)/nrj-2-jstag,cjte+2)
1752 cfldnew(i,j) = cfld(i,k,j) + xnu(n) * ((cfld(i,k,j+1) + cfld(i,k,j-1)) * 0.5-cfld(i,k,j))
1756 DO i = MAX(ipos+2,cits-2) , MIN(ipos+(nide-nids)/nri-2-istag,cite+2)
1757 DO j = MAX(jpos+2,cjts-2) , MIN(jpos+(njde-njds)/nrj-2-jstag,cjte+2)
1758 cfld(i,k,j) = cfldnew(i,j)
1762 DO j = MAX(jpos+2,cjts-2) , MIN(jpos+(njde-njds)/nrj-2-jstag,cjte+2)
1763 DO i = MAX(ipos+2,cits-2) , MIN(ipos+(nide-nids)/nri-2-istag,cite+2)
1764 cfldnew(i,j) = cfld(i,k,j) + xnu(n) * ((cfld(i+1,k,j) + cfld(i-1,k,j)) * 0.5-cfld(i,k,j))
1768 DO j = MAX(jpos+2,cjts-2) , MIN(jpos+(njde-njds)/nrj-2-jstag,cjte+2)
1769 DO i = MAX(ipos+2,cits-2) , MIN(ipos+(nide-nids)/nri-2-istag,cite+2)
1770 cfld(i,k,j) = cfldnew(i,j)
1776 END DO smoothing_passes
1778 END SUBROUTINE smdsm
1780 !==================================
1781 ! this is used to modify a field over the nest so we can see where the nest is
1783 SUBROUTINE mark_domain ( cfld, & ! CD field
1784 cids, cide, ckds, ckde, cjds, cjde, &
1785 cims, cime, ckms, ckme, cjms, cjme, &
1786 cits, cite, ckts, ckte, cjts, cjte, &
1788 nids, nide, nkds, nkde, njds, njde, &
1789 nims, nime, nkms, nkme, njms, njme, &
1790 nits, nite, nkts, nkte, njts, njte, &
1791 shw, & ! stencil half width for interp
1792 imask, & ! interpolation mask
1793 xstag, ystag, & ! staggering of field
1794 ipos, jpos, & ! Position of lower left of nest in CD
1795 nri, nrj ) ! nest ratios
1796 USE module_configure
1797 USE module_wrf_error
1801 INTEGER, INTENT(IN) :: cids, cide, ckds, ckde, cjds, cjde, &
1802 cims, cime, ckms, ckme, cjms, cjme, &
1803 cits, cite, ckts, ckte, cjts, cjte, &
1804 nids, nide, nkds, nkde, njds, njde, &
1805 nims, nime, nkms, nkme, njms, njme, &
1806 nits, nite, nkts, nkte, njts, njte, &
1810 LOGICAL, INTENT(IN) :: xstag, ystag
1812 REAL, DIMENSION ( cims:cime, ckms:ckme, cjms:cjme ), INTENT(OUT) :: cfld
1813 REAL, DIMENSION ( nims:nime, nkms:nkme, njms:njme ), INTENT(IN) :: nfld
1814 INTEGER, DIMENSION ( nims:nime, njms:njme ), INTENT(IN) :: imask
1818 INTEGER ci, cj, ck, ni, nj, nk, ip, jp, ioff, joff, ioffa, joffa
1819 INTEGER :: icmin,icmax,jcmin,jcmax
1820 INTEGER :: istag,jstag, ipoints,jpoints,ijpoints
1822 istag = 1 ; jstag = 1
1823 IF ( xstag ) istag = 0
1824 IF ( ystag ) jstag = 0
1826 DO cj = MAX(jpos+1,cjts),MIN(jpos+(njde-njds)/nrj-jstag-1,cjte)
1827 nj = (cj-jpos)*nrj + jstag + 1
1830 DO ci = MAX(ipos+1,cits),MIN(ipos+(nide-nids)/nri-istag-1,cite)
1831 ni = (ci-ipos)*nri + istag + 1
1832 cfld( ci, ck, cj ) = 9021000. !magic number: Beverly Hills * 100.
1837 END SUBROUTINE mark_domain
1839 #if ( NMM_CORE == 1 )
1840 !=======================================================================================
1841 ! E grid interpolation for mass with addition of terrain adjustments. First routine
1842 ! pertains to initial conditions and the next one corresponds to boundary conditions
1843 ! This is gopal's doing
1844 !=======================================================================================
1846 SUBROUTINE interp_mass_nmm (cfld, & ! CD field
1847 cids, cide, ckds, ckde, cjds, cjde, &
1848 cims, cime, ckms, ckme, cjms, cjme, &
1849 cits, cite, ckts, ckte, cjts, cjte, &
1851 nids, nide, nkds, nkde, njds, njde, &
1852 nims, nime, nkms, nkme, njms, njme, &
1853 nits, nite, nkts, nkte, njts, njte, &
1854 shw, & ! stencil half width for interp
1855 imask, & ! interpolation mask
1856 xstag, ystag, & ! staggering of field
1857 ipos, jpos, & ! Position of lower left of nest in CD
1858 nri, nrj, & ! nest ratios
1859 CII, IIH, CJJ, JJH, CBWGT1, HBWGT1, & ! south-western grid locs and weights
1860 CBWGT2, HBWGT2, CBWGT3, HBWGT3, & ! note that "C"ourse grid ones are
1861 CBWGT4, HBWGT4, & ! dummys for weights
1862 CZ3d, Z3d, & ! Z3d interpolated from CZ3d
1863 CFIS,FIS, & ! CFIS dummy on fine domain
1864 CSM,SM, & ! CSM is dummy
1870 USE MODULE_MODEL_CONSTANTS
1874 LOGICAL,INTENT(IN) :: xstag, ystag
1875 INTEGER,INTENT(IN) :: ckzmax,kzmax
1876 INTEGER,INTENT(IN) :: cids, cide, ckds, ckde, cjds, cjde, &
1877 cims, cime, ckms, ckme, cjms, cjme, &
1878 cits, cite, ckts, ckte, cjts, cjte, &
1879 nids, nide, nkds, nkde, njds, njde, &
1880 nims, nime, nkms, nkme, njms, njme, &
1881 nits, nite, nkts, nkte, njts, njte, &
1882 shw,ipos,jpos,nri,nrj
1884 INTEGER,DIMENSION(nims:nime,njms:njme), INTENT(IN) :: IMASK
1888 INTEGER,DIMENSION(cims:cime,cjms:cjme), INTENT(IN) :: CII,CJJ ! dummy
1889 REAL,DIMENSION(cims:cime,cjms:cjme), INTENT(IN) :: CBWGT1,CBWGT2,CBWGT3
1890 REAL,DIMENSION(cims:cime,cjms:cjme), INTENT(IN) :: CBWGT4,CFIS,CSM
1891 REAL,DIMENSION(cims:cime,cjms:cjme,ckms:ckme), INTENT(IN) :: CFLD
1892 REAL,DIMENSION(cims:cime,cjms:cjme,1:KZMAX), INTENT(IN) :: CZ3d
1893 REAL,DIMENSION(1:KZMAX), INTENT(IN) :: CPSTD
1894 REAL,INTENT(IN) :: CPDTOP,CPTOP
1898 INTEGER,DIMENSION(nims:nime,njms:njme), INTENT(IN) :: IIH,JJH
1899 REAL,DIMENSION(nims:nime,njms:njme), INTENT(IN) :: HBWGT1,HBWGT2,HBWGT3
1900 REAL,DIMENSION(nims:nime,njms:njme), INTENT(IN) :: HBWGT4
1901 REAL,DIMENSION(nims:nime,njms:njme), INTENT(IN) :: FIS,SM
1902 REAL,DIMENSION(nims:nime,njms:njme,nkms:nkme), INTENT(INOUT) :: NFLD
1903 REAL,DIMENSION(1:KZMAX), INTENT(IN) :: PSTD
1904 REAL,DIMENSION(nims:nime,njms:njme,1:KZMAX), INTENT(OUT) :: Z3d
1905 REAL,INTENT(IN) :: PDTOP,PTOP
1909 INTEGER,PARAMETER :: JTB=134
1910 REAL, PARAMETER :: LAPSR=6.5E-3,GI=1./G, D608=0.608
1911 REAL, PARAMETER :: COEF3=R_D*GI*LAPSR
1912 INTEGER :: I,J,K,IDUM
1913 REAL :: dlnpdz,tvout,pmo
1914 REAL,DIMENSION(nims:nime,njms:njme) :: ZS,DUM2d
1915 REAL,DIMENSION(JTB) :: PIN,ZIN,Y2,PIO,ZOUT,DUM1,DUM2
1916 !-----------------------------------------------------------------------------------------------------
1918 !*** CHECK DOMAIN BOUNDS BEFORE INTERPOLATION
1920 DO J=NJTS,MIN(NJTE,NJDE-1)
1921 DO I=NITS,MIN(NITE,NIDE-1)
1922 IF(IIH(i,j).LT.(CIDS-shw) .OR. IIH(i,j).GT.(CIDE+shw)) &
1923 CALL wrf_error_fatal ('mass points:check domain bounds along x' )
1924 IF(JJH(i,j).LT.(CJDS-shw) .OR. JJH(i,j).GT.(CJDE+shw)) &
1925 CALL wrf_error_fatal ('mass points:check domain bounds along y' )
1929 IF(KZMAX .GT. (JTB-10)) &
1930 CALL wrf_error_fatal ('mass points: increase JTB in interp_mass_nmm')
1932 ! WRITE(21,*)'------------- MED NEST INITIAL 1 ----------------'
1933 ! DO J=NJTS,MIN(NJTE,NJDE-1)
1934 ! DO I=NITS,MIN(NITE,NIDE-1)
1935 ! WRITE(21,*)I,J,IMASK(I,J),NFLD(I,1,J)
1941 !*** DEFINE LOCAL TOPOGRAPHY ON THE BASIS OF FIS. ALSO CHECK IF SM IS LAND (SM=0) OVER TOPO
1942 !*** YOU DON'T WANT MOUNTAINS INSIDE WATER BODIES!
1945 DO J=NJTS,MIN(NJTE,NJDE-1)
1946 DO I=NITS,MIN(NITE,NIDE-1)
1952 !*** Interpolate GPMs DERIVED FROM STANDARD ATMOSPHERIC LAPSE RATE FROM THE PARENT TO
1953 !*** THE NESTED DOMAIN
1955 !*** INDEX CONVENTIONS
1970 DO K=NKTS,KZMAX ! Please note that we are still in isobaric surfaces
1971 DO J=NJTS,MIN(NJTE,NJDE-1)
1972 DO I=NITS,MIN(NITE,NIDE-1)
1974 IF(MOD(JJH(I,J),2) .NE. 0)THEN ! 1,3,5,7
1975 Z3d(I,J,K) = HBWGT1(I,J)*CZ3d(IIH(I,J), JJH(I,J) ,K) &
1976 + HBWGT2(I,J)*CZ3d(IIH(I,J)+1,JJH(I,J) ,K) &
1977 + HBWGT3(I,J)*CZ3d(IIH(I,J), JJH(I,J)-1,K) &
1978 + HBWGT4(I,J)*CZ3d(IIH(I,J), JJH(I,J)+1,K)
1980 Z3d(I,J,K) = HBWGT1(I,J)*CZ3d(IIH(I,J), JJH(I,J) ,K) &
1981 + HBWGT2(I,J)*CZ3d(IIH(I,J)+1,JJH(I,J) ,K) &
1982 + HBWGT3(I,J)*CZ3d(IIH(I,J)+1,JJH(I,J)-1,K) &
1983 + HBWGT4(I,J)*CZ3d(IIH(I,J)+1,JJH(I,J)+1,K)
1991 ! RECONSTRUCT PDs ON THE BASIS OF TOPOGRAPHY AND THE INTERPOLATED HEIGHTS
1993 DO J=NJTS,MIN(NJTE,NJDE-1)
1994 DO I=NITS,MIN(NITE,NIDE-1)
1996 IF (ZS(I,J) .LT. Z3d(I,J,1)) THEN
1997 dlnpdz = (log(PSTD(1))-log(PSTD(2)) )/(Z3d(i,j,1)-Z3d(i,j,2))
1998 dum2d(i,j) = exp(log(PSTD(1)) + dlnpdz*(ZS(I,J) - Z3d(i,j,1)))
1999 dum2d(i,j) = dum2d(i,j) - PDTOP -PTOP
2000 ELSE ! target level bounded by input levels
2001 DO K =NKTS,KZMAX-1 ! still in the isobaric surfaces
2002 IF(ZS(I,J) .GE. Z3d(I,J,K) .AND. ZS(I,J) .LT. Z3d(I,J,K+1))THEN
2003 dlnpdz = (log(PSTD(K))-log(PSTD(K+1)) ) /(Z3d(I,J,K)-Z3d(I,J,K+1))
2004 dum2d(i,j) = exp(log(PSTD(K)) + dlnpdz*(ZS(I,J)- Z3d(I,J,K)))
2005 dum2d(i,j) = dum2d(i,j) - PDTOP -PTOP
2009 IF(ZS(I,J) .GE. Z3d(I,J,KZMAX))THEN
2010 WRITE(0,*)'I=',I,'J=',J,'K=',KZMAX,'TERRAIN HEIGHT',ZS(I,J),'Z3d',Z3d(I,J,KZMAX)
2011 CALL wrf_error_fatal3 ( "interp_fcn.b" , 176 , "MOUNTAIN TOO HIGH TO FIT THE MODEL DEPTH")
2017 DO K=NKDS,NKDE ! NKTE is 1, nevertheless let us pretend religious
2018 DO J=NJTS,MIN(NJTE,NJDE-1)
2019 DO I=NITS,MIN(NITE,NIDE-1)
2020 IF(IMASK(I,J) .NE. 1)THEN
2021 NFLD(I,J,K)= dum2d(i,j) ! PD defined in the nested domain
2028 END SUBROUTINE interp_mass_nmm
2030 !--------------------------------------------------------------------------------------
2032 SUBROUTINE nmm_bdymass_hinterp ( cfld, & ! CD field
2033 cids, cide, ckds, ckde, cjds, cjde, &
2034 cims, cime, ckms, ckme, cjms, cjme, &
2035 cits, cite, ckts, ckte, cjts, cjte, &
2037 nids, nide, nkds, nkde, njds, njde, &
2038 nims, nime, nkms, nkme, njms, njme, &
2039 nits, nite, nkts, nkte, njts, njte, &
2040 shw, & ! stencil half width
2041 imask, & ! interpolation mask
2042 xstag, ystag, & ! staggering of field
2043 ipos, jpos, & ! Position of lower left of nest in CD
2044 nri, nrj, & ! nest ratios
2053 CTEMP_B,NTEMP_B, & ! These temp arrays should be removed
2054 CTEMP_BT,NTEMP_BT, & ! later on
2055 CII, IIH, CJJ, JJH, CBWGT1, HBWGT1, & ! south-western grid locs and weights
2056 CBWGT2, HBWGT2, CBWGT3, HBWGT3, & ! note that "C"ourse grid ones are
2057 CBWGT4, HBWGT4, & ! dummys
2058 CZ3d, Z3d, & ! Z3d dummy on nested domain
2059 CFIS,FIS, & ! CFIS dummy on fine domain
2060 CSM,SM, & ! CSM is dummy
2067 USE MODULE_MODEL_CONSTANTS
2068 USE module_configure
2069 USE module_wrf_error
2074 INTEGER, INTENT(IN) :: ckzmax,kzmax
2075 INTEGER, INTENT(IN) :: cids, cide, ckds, ckde, cjds, cjde, &
2076 cims, cime, ckms, ckme, cjms, cjme, &
2077 cits, cite, ckts, ckte, cjts, cjte, &
2078 nids, nide, nkds, nkde, njds, njde, &
2079 nims, nime, nkms, nkme, njms, njme, &
2080 nits, nite, nkts, nkte, njts, njte, &
2086 REAL, DIMENSION ( cims:cime, cjms:cjme, ckms:ckme ), INTENT(OUT) :: ctemp_b,ctemp_bt
2087 REAL, DIMENSION ( nims:nime, njms:njme, nkms:nkme ), INTENT(OUT) :: ntemp_b,ntemp_bt
2088 LOGICAL, INTENT(IN) :: xstag, ystag
2089 REAL, DIMENSION( * ), INTENT(INOUT) :: c_bxs,n_bxs,c_bxe,n_bxe,c_bys,n_bys,c_bye,n_bye
2090 REAL, DIMENSION( * ), INTENT(INOUT) :: c_btxs,n_btxs,c_btxe,n_btxe,c_btys,n_btys,c_btye,n_btye
2094 INTEGER,DIMENSION(nims:nime,njms:njme), INTENT(IN) :: IMASK
2095 INTEGER,DIMENSION(cims:cime,cjms:cjme), INTENT(IN) :: CII,CJJ ! dummy
2096 REAL,DIMENSION(cims:cime,cjms:cjme), INTENT(IN) :: CBWGT1,CBWGT2,CBWGT3
2097 REAL,DIMENSION(cims:cime,cjms:cjme), INTENT(IN) :: CBWGT4,CFIS,CSM
2098 REAL,DIMENSION(cims:cime,cjms:cjme,ckms:ckme), INTENT(IN) :: CFLD
2099 REAL,DIMENSION(cims:cime,cjms:cjme,1:KZMAX), INTENT(IN) :: CZ3d
2100 REAL,DIMENSION(1:KZMAX), INTENT(IN) :: CPSTD
2101 REAL,INTENT(IN) :: CPDTOP,CPTOP
2105 INTEGER,DIMENSION(nims:nime,njms:njme), INTENT(IN) :: IIH,JJH
2106 REAL,DIMENSION(nims:nime,njms:njme), INTENT(IN) :: HBWGT1,HBWGT2,HBWGT3
2107 REAL,DIMENSION(nims:nime,njms:njme), INTENT(IN) :: HBWGT4
2108 REAL,DIMENSION(nims:nime,njms:njme), INTENT(IN) :: FIS,SM
2109 REAL,DIMENSION(nims:nime,njms:njme,nkms:nkme), INTENT(INOUT) :: NFLD
2110 REAL,DIMENSION(1:KZMAX), INTENT(IN) :: PSTD
2111 REAL,DIMENSION(nims:nime,njms:njme,1:KZMAX), INTENT(OUT) :: Z3d
2112 REAL,INTENT(IN) :: PDTOP,PTOP
2116 INTEGER :: nijds, nijde, spec_bdy_width,i,j,k
2117 REAL :: dlnpdz,dum2d
2118 REAL,DIMENSION(nims:nime,njms:njme) :: zs
2120 INTEGER,PARAMETER :: JTB=134
2122 REAL, DIMENSION (nims:nime,njms:njme) :: CWK1,CWK2,CWK3,CWK4
2124 nijds = min(nids, njds)
2125 nijde = max(nide, njde)
2126 CALL nl_get_spec_bdy_width( 1, spec_bdy_width )
2129 !*** DEFINE LOCAL TOPOGRAPHY ON THE BASIS OF FIS. ASLO CHECK IF SM IS LAND (SM=0) OVER TOPO
2130 !*** YOU DON'T WANT MOUNTAINS INSIDE WATER BODIES!
2132 DO J=NJTS,MIN(NJTE,NJDE-1)
2133 DO I=NITS,MIN(NITE,NIDE-1)
2140 NMM_XS: IF(NITS .EQ. NIDS)THEN
2141 ! WRITE(0,*)'ENTERING X1 START BOUNDARY AT MASS POINTS',NJTS,MIN(NJTE,NJDE-1)
2145 DO J = NJTS,MIN(NJTE,NJDE-1)
2146 IF(MOD(JJH(I,J),2) .NE. 0)THEN ! 1,3,5,7 of the parent domain
2147 Z3d(I,J,K) = HBWGT1(I,J)*CZ3d(IIH(I,J), JJH(I,J) ,K) &
2148 + HBWGT2(I,J)*CZ3d(IIH(I,J)+1,JJH(I,J) ,K) &
2149 + HBWGT3(I,J)*CZ3d(IIH(I,J), JJH(I,J)-1,K) &
2150 + HBWGT4(I,J)*CZ3d(IIH(I,J), JJH(I,J)+1,K)
2152 Z3d(I,J,K) = HBWGT1(I,J)*CZ3d(IIH(I,J), JJH(I,J) ,K) &
2153 + HBWGT2(I,J)*CZ3d(IIH(I,J)+1,JJH(I,J) ,K) &
2154 + HBWGT3(I,J)*CZ3d(IIH(I,J)+1,JJH(I,J)-1,K) &
2155 + HBWGT4(I,J)*CZ3d(IIH(I,J)+1,JJH(I,J)+1,K)
2160 DO J = NJTS,MIN(NJTE,NJDE-1)
2161 IF(MOD(J,2) .NE. 0)THEN
2162 IF (ZS(I,J) .LT. Z3d(I,J,2)) THEN ! level 2 has to be changed
2163 dlnpdz = (log(PSTD(1))-log(PSTD(2)) )/(Z3d(I,J,1)-Z3d(I,J,2))
2164 dum2d = exp(log(PSTD(1)) + dlnpdz*(ZS(I,J) - Z3d(I,J,1)))
2165 CWK1(I,J) = dum2d -PDTOP -PTOP
2166 ELSE ! target level bounded by input levels
2168 IF(ZS(I,J) .GE. Z3d(I,J,K) .AND. ZS(I,J) .LT. Z3d(I,J,K+1))THEN
2169 dlnpdz = (log(PSTD(K))-log(PSTD(K+1)) ) /(Z3d(I,J,K)-Z3d(I,J,K+1))
2170 dum2d = exp(log(PSTD(K)) + dlnpdz*(ZS(I,J)- Z3d(I,J,K)))
2171 CWK1(I,J) = dum2d -PDTOP -PTOP
2175 IF(ZS(I,J) .GE. Z3d(I,J,KZMAX))THEN
2176 WRITE(0,*)'I=',I,'J=',J,'K=',K,'TERRAIN HEIGHT',ZS(I,J),'Z3d',Z3d(I,J,KZMAX)
2177 CALL wrf_error_fatal("BC:MOUNTAIN TOO HIGH TO FIT THE MODEL DEPTH")
2184 DO J = NJTS,MIN(NJTE,NJDE-1)
2186 ntemp_b(i,j,k) = CWK1(I,J)
2187 ntemp_bt(i,j,k) = 0.0
2194 NMM_XE: IF(NITE-1 .EQ. NIDE-1)THEN
2195 ! WRITE(0,*)'ENTERING X END BOUNDARY AT MASS POINTS',NJTS,MIN(NJTE,NJDE-1)
2200 DO J=NJTS,MIN(NJTE,NJDE-1)
2201 IF(MOD(JJH(I,J),2) .NE. 0)THEN ! 1,3,5,7
2202 Z3d(I,J,K) = HBWGT1(I,J)*CZ3d(IIH(I,J), JJH(I,J) ,K) &
2203 + HBWGT2(I,J)*CZ3d(IIH(I,J)+1,JJH(I,J) ,K) &
2204 + HBWGT3(I,J)*CZ3d(IIH(I,J), JJH(I,J)-1,K) &
2205 + HBWGT4(I,J)*CZ3d(IIH(I,J), JJH(I,J)+1,K)
2207 Z3d(I,J,K) = HBWGT1(I,J)*CZ3d(IIH(I,J), JJH(I,J) ,K) &
2208 + HBWGT2(I,J)*CZ3d(IIH(I,J)+1,JJH(I,J) ,K) &
2209 + HBWGT3(I,J)*CZ3d(IIH(I,J)+1,JJH(I,J)-1,K) &
2210 + HBWGT4(I,J)*CZ3d(IIH(I,J)+1,JJH(I,J)+1,K)
2215 DO J = NJTS,MIN(NJTE,NJDE-1)
2216 IF(MOD(J,2) .NE.0)THEN ! 1,3,5,7 of nested domain
2217 IF (ZS(I,J) .LT. Z3d(I,J,2)) THEN ! level 2 has to be changed
2218 dlnpdz = (log(PSTD(1))-log(PSTD(2)) )/(Z3d(I,J,1)-Z3d(I,J,2))
2219 dum2d = exp(log(PSTD(1)) + dlnpdz*(ZS(I,J) - Z3d(I,J,1)))
2220 CWK2(I,J) = dum2d -PDTOP -PTOP
2221 ELSE ! target level bounded by input levels
2223 IF(ZS(I,J) .GE. Z3d(I,J,K) .AND. ZS(I,J) .LT. Z3d(I,J,K+1))THEN
2224 dlnpdz = (log(PSTD(K))-log(PSTD(K+1)) ) /(Z3d(I,J,K)-Z3d(I,J,K+1))
2225 dum2d = exp(log(PSTD(K)) + dlnpdz*(ZS(I,J)- Z3d(I,J,K)))
2226 CWK2(I,J) = dum2d -PDTOP -PTOP
2230 IF(ZS(I,J) .GE. Z3d(I,J,KZMAX))THEN
2231 WRITE(0,*)'I=',I,'J=',J,'K=',K,'TERRAIN HEIGHT',ZS(I,J),'Z3d',Z3d(I,J,KZMAX)
2232 CALL wrf_error_fatal("BC:MOUNTAIN TOO HIGH TO FIT THE MODEL DEPTH")
2239 DO J = NJTS,MIN(NJTE,NJDE-1)
2241 ntemp_b(i,j,k) = CWK2(I,J)
2242 ntemp_bt(i,j,k) = 0.0
2249 NMM_YS: IF(NJTS .EQ. NJDS)THEN
2250 ! WRITE(20,*)'ENTERING Y START BOUNDARY AT MASS POINTS',NITS,MIN(NITE,NIDE-1)
2253 DO I = NITS,MIN(NITE,NIDE-1)
2254 IF(MOD(JJH(I,J),2) .NE. 0)THEN ! 1,3,5,7
2255 Z3d(I,J,K) = HBWGT1(I,J)*CZ3d(IIH(I,J), JJH(I,J) ,K) &
2256 + HBWGT2(I,J)*CZ3d(IIH(I,J)+1,JJH(I,J) ,K) &
2257 + HBWGT3(I,J)*CZ3d(IIH(I,J), JJH(I,J)-1,K) &
2258 + HBWGT4(I,J)*CZ3d(IIH(I,J), JJH(I,J)+1,K)
2260 Z3d(I,J,K) = HBWGT1(I,J)*CZ3d(IIH(I,J), JJH(I,J) ,K) &
2261 + HBWGT2(I,J)*CZ3d(IIH(I,J)+1,JJH(I,J) ,K) &
2262 + HBWGT3(I,J)*CZ3d(IIH(I,J)+1,JJH(I,J)-1,K) &
2263 + HBWGT4(I,J)*CZ3d(IIH(I,J)+1,JJH(I,J)+1,K)
2268 DO I = NITS,MIN(NITE,NIDE-1)
2269 IF (ZS(I,J) .LT. Z3d(I,J,2)) THEN ! level 2 has to be changed
2270 dlnpdz = (log(PSTD(1))-log(PSTD(2)) )/(Z3d(I,J,1)-Z3d(I,J,2))
2271 dum2d = exp(log(PSTD(1)) + dlnpdz*(ZS(I,J) - Z3d(I,J,1)))
2272 CWK3(I,J) = dum2d -PDTOP -PTOP
2273 ELSE ! target level bounded by input levels
2275 IF(ZS(I,J) .GE. Z3d(I,J,K) .AND. ZS(I,J) .LT. Z3d(I,J,K+1))THEN
2276 dlnpdz = (log(PSTD(K))-log(PSTD(K+1)) ) /(Z3d(I,J,K)-Z3d(I,J,K+1))
2277 dum2d = exp(log(PSTD(K)) + dlnpdz*(ZS(I,J)- Z3d(I,J,K)))
2278 CWK3(I,J) = dum2d -PDTOP -PTOP
2282 IF(ZS(I,J) .GE. Z3d(I,J,KZMAX))THEN
2283 WRITE(0,*)'I=',I,'J=',J,'K=',K,'TERRAIN HEIGHT',ZS(I,J),'Z3d',Z3d(I,J,KZMAX)
2284 CALL wrf_error_fatal("BC:MOUNTAIN TOO HIGH TO FIT THE MODEL DEPTH")
2289 DO I = NITS,MIN(NITE,NIDE-1)
2290 ntemp_b(i,j,k) = CWK3(I,J)
2291 ntemp_bt(i,j,k) = 0.0
2298 NMM_YE: IF(NJTE-1 .EQ. NJDE-1)THEN
2299 ! WRITE(20,*)'ENTERING Y END BOUNDARY AT MASS POINTS',NITS,MIN(NITE,NIDE-1)
2303 DO I = NITS,MIN(NITE,NIDE-1)
2304 IF(MOD(JJH(I,J),2) .NE. 0)THEN ! 1,3,5,7
2305 Z3d(I,J,K) = HBWGT1(I,J)*CZ3d(IIH(I,J), JJH(I,J) ,K) &
2306 + HBWGT2(I,J)*CZ3d(IIH(I,J)+1,JJH(I,J) ,K) &
2307 + HBWGT3(I,J)*CZ3d(IIH(I,J), JJH(I,J)-1,K) &
2308 + HBWGT4(I,J)*CZ3d(IIH(I,J), JJH(I,J)+1,K)
2310 Z3d(I,J,K) = HBWGT1(I,J)*CZ3d(IIH(I,J), JJH(I,J) ,K) &
2311 + HBWGT2(I,J)*CZ3d(IIH(I,J)+1,JJH(I,J) ,K) &
2312 + HBWGT3(I,J)*CZ3d(IIH(I,J)+1,JJH(I,J)-1,K) &
2313 + HBWGT4(I,J)*CZ3d(IIH(I,J)+1,JJH(I,J)+1,K)
2318 DO I = NITS,MIN(NITE,NIDE-1)
2319 IF (ZS(I,J) .LT. Z3d(I,J,2)) THEN ! level 2 has to be changed
2320 dlnpdz = (log(PSTD(1))-log(PSTD(2)) )/(Z3d(I,J,1)-Z3d(I,J,2))
2321 dum2d = exp(log(PSTD(1)) + dlnpdz*(ZS(I,J) - Z3d(I,J,1)))
2322 CWK4(I,J) = dum2d -PDTOP -PTOP
2323 ELSE ! target level bounded by input levels
2325 IF(ZS(I,J) .GE. Z3d(I,J,K) .AND. ZS(I,J) .LT. Z3d(I,J,K+1))THEN
2326 dlnpdz = (log(PSTD(K))-log(PSTD(K+1)) ) /(Z3d(I,J,K)-Z3d(I,J,K+1))
2327 dum2d = exp(log(PSTD(K)) + dlnpdz*(ZS(I,J)- Z3d(I,J,K)))
2328 CWK4(I,J) = dum2d -PDTOP -PTOP
2332 IF(ZS(I,J) .GE. Z3d(I,J,KZMAX))THEN
2333 WRITE(0,*)'I=',I,'J=',J,'K=',K,'TERRAIN HEIGHT',ZS(I,J),'Z3d',Z3d(I,J,KZMAX)
2334 CALL wrf_error_fatal("BC:MOUNTAIN TOO HIGH TO FIT THE MODEL DEPTH")
2339 DO I = NITS,MIN(NITE,NIDE-1)
2340 ntemp_b(i,j,k) = CWK4(I,J)
2341 ntemp_bt(i,j,k) = 0.0
2348 END SUBROUTINE nmm_bdymass_hinterp
2350 !=======================================================================================
2352 ! ADDED FOR INCLUDING MOISTURE AND THERMODYNAMIC ENERGY BALANCE
2354 !=======================================================================================
2356 SUBROUTINE interp_scalar_nmm (cfld, & ! CD field
2357 cids,cide,ckds,ckde,cjds,cjde, &
2358 cims,cime,ckms,ckme,cjms,cjme, &
2359 cits,cite,ckts,ckte,cjts,cjte, &
2361 nids,nide,nkds,nkde,njds,njde, &
2362 nims,nime,nkms,nkme,njms,njme, &
2363 nits,nite,nkts,nkte,njts,njte, &
2364 shw, & ! stencil half width for interp
2365 imask, & ! interpolation mask
2366 xstag,ystag, & ! staggering of field
2367 ipos,jpos, & ! Position of lower left of nest in CD
2368 nri,nrj, & ! nest ratios
2369 CII, IIH, CJJ, JJH, CBWGT1, HBWGT1, & ! south-western grid locs and weights
2370 CBWGT2, HBWGT2, CBWGT3, HBWGT3, & ! note that "C"ourse grid ones are
2371 CBWGT4, HBWGT4, & ! dummys for weights
2377 CETA1,ETA1,CETA2,ETA2 )
2379 USE MODULE_MODEL_CONSTANTS
2383 LOGICAL,INTENT(IN) :: xstag, ystag
2384 INTEGER,INTENT(IN) :: cids, cide, ckds, ckde, cjds, cjde, &
2385 cims, cime, ckms, ckme, cjms, cjme, &
2386 cits, cite, ckts, ckte, cjts, cjte, &
2387 nids, nide, nkds, nkde, njds, njde, &
2388 nims, nime, nkms, nkme, njms, njme, &
2389 nits, nite, nkts, nkte, njts, njte, &
2390 shw,ipos,jpos,nri,nrj
2392 INTEGER,DIMENSION(nims:nime,njms:njme), INTENT(IN) :: IMASK
2396 INTEGER,DIMENSION(cims:cime,cjms:cjme), INTENT(IN) :: CII,CJJ ! dummy
2397 REAL,DIMENSION(cims:cime,cjms:cjme), INTENT(IN) :: CBWGT1,CBWGT2
2398 REAL,DIMENSION(cims:cime,cjms:cjme), INTENT(IN) :: CBWGT3,CBWGT4
2400 REAL,DIMENSION(cims:cime,cjms:cjme,ckms:ckme), INTENT(IN) :: CFLD
2401 REAL,DIMENSION(cims:cime,cjms:cjme,ckms:ckme), INTENT(IN) :: CC3d ! scalar input on constant pressure levels
2402 REAL,DIMENSION(ckms:ckme), INTENT(IN) :: CPSTD
2403 REAL,DIMENSION(cims:cime,cjms:cjme), INTENT(IN) :: CPD
2404 REAL,DIMENSION(ckms:ckme), INTENT(IN) :: CETA1,CETA2
2405 REAL, INTENT(IN) :: CPDTOP,CPTOP
2409 INTEGER,DIMENSION(nims:nime,njms:njme), INTENT(IN) :: IIH,JJH
2410 REAL,DIMENSION(nims:nime,njms:njme), INTENT(IN) :: HBWGT1,HBWGT2
2411 REAL,DIMENSION(nims:nime,njms:njme), INTENT(IN) :: HBWGT3,HBWGT4
2413 REAL,DIMENSION(nims:nime,njms:njme,nkms:nkme), INTENT(OUT):: NFLD ! This is scalar on hybrid levels
2414 REAL,DIMENSION(nims:nime,njms:njme,nkms:nkme), INTENT(OUT):: C3d ! Scalar on constant pressure levels
2415 REAL,DIMENSION(nkms:nkme), INTENT(IN) :: PSTD
2416 REAL,DIMENSION(nims:nime,njms:njme), INTENT(IN) :: PD
2417 REAL,DIMENSION(nkms:nkme), INTENT(IN) :: ETA1,ETA2
2418 REAL,INTENT(IN) :: PDTOP,PTOP
2422 INTEGER,PARAMETER :: JTB=134
2424 REAL,DIMENSION(JTB) :: PIN,CIN,Y2,PIO,PTMP,COUT,DUM1,DUM2
2426 !-----------------------------------------------------------------------------------------------------
2429 ! *** CHECK VERTICAL BOUNDS BEFORE USING SPLINE OR LINEAR INTERPOLATION
2431 IF(nkme .GT. (JTB-10) .OR. NKDE .GT. (JTB-10)) &
2432 CALL wrf_error_fatal ('mass points: increase JTB in interp_mass_nmm')
2435 ! FIRST, HORIZONTALLY INTERPOLATE MOISTURE NOW AVAILABLE ON CONSTANT PRESSURE SURFACE (LEVELS) FROM THE
2436 ! PARENT TO THE NESTED DOMAIN
2438 !*** INDEX CONVENTIONS
2453 DO K=NKDS,NKDE-1 ! Please note that we are still in isobaric surfaces
2454 DO J=NJTS,MIN(NJTE,NJDE-1)
2455 DO I=NITS,MIN(NITE,NIDE-1)
2456 IF(IMASK(I,J) .NE. 1)THEN
2457 IF(MOD(JJH(I,J),2) .NE. 0)THEN ! 1,3,5,7
2458 C3d(I,J,K) = HBWGT1(I,J)*CC3d(IIH(I,J), JJH(I,J) ,K) &
2459 + HBWGT2(I,J)*CC3d(IIH(I,J)+1,JJH(I,J) ,K) &
2460 + HBWGT3(I,J)*CC3d(IIH(I,J), JJH(I,J)-1,K) &
2461 + HBWGT4(I,J)*CC3d(IIH(I,J), JJH(I,J)+1,K)
2464 C3d(I,J,K) = HBWGT1(I,J)*CC3d(IIH(I,J), JJH(I,J) ,K) &
2465 + HBWGT2(I,J)*CC3d(IIH(I,J)+1,JJH(I,J) ,K) &
2466 + HBWGT3(I,J)*CC3d(IIH(I,J)+1,JJH(I,J)-1,K) &
2467 + HBWGT4(I,J)*CC3d(IIH(I,J)+1,JJH(I,J)+1,K)
2476 ! RECOVER THE SCALARS FROM CONSTANT PRESSURE SURFACES (LEVELS) ON TO HYBRID SURFACES
2478 DO J=NJTS,MIN(NJTE,NJDE-1)
2479 DO I=NITS,MIN(NITE,NIDE-1)
2480 IF(IMASK(I,J) .NE. 1)THEN
2482 ! clean local array before use of spline or linear interpolation
2484 CIN=0.;PIN=0.;Y2=0;PIO=0.;PTMP=0.;COUT=0.;DUM1=0.;DUM2=0.
2486 DO K=NKDS+1,NKDE ! inputs at standard levels
2487 PIN(K-1) = EXP((ALOG(PSTD(NKDE-K+1))+ALOG(PSTD(NKDE-K+2)))*0.5)
2488 CIN(K-1) = C3d(I,J,NKDE-K+1)
2494 DO K=NKDS,NKDE ! target points in model interface levels (pint)
2495 PTMP(K) = ETA1(K)*PDTOP + ETA2(K)*PD(I,J) + PTOP
2498 DO K=NKDS,NKDE-1 ! target points in model levels
2499 PIO(K) = EXP((ALOG(PTMP(K))+ALOG(PTMP(K+1)))*0.5)
2502 IF(PTMP(1) .GE. PSTD(1))THEN ! if lower boundary is higher than PMSL(1) re-set lower boundary
2503 PIN(NKDE-1) = PIO(1) ! be consistent with target. This may not happen at all
2504 WRITE(0,*)'WARNING: NESTED DOMAIN PRESSURE AT LOWEST LEVEL HIGHER THAN PSTD'
2505 WRITE(0,*)'I,J,PIO(1),PSTD(1)',I,J,PIO(1),PSTD(1)
2508 CALL SPLINE2(I,J,JTB,NKDE-1,PIN,CIN,Y2,NKDE-1,PIO,COUT,DUM1,DUM2) ! interpolate
2511 NFLD(I,J,K)= COUT(K) ! scalar in the nested domain
2518 END SUBROUTINE interp_scalar_nmm
2520 !===========================================================================================
2522 SUBROUTINE nmm_bdy_scalar (cfld, & ! CD field
2523 cids,cide,ckds,ckde,cjds,cjde, &
2524 cims,cime,ckms,ckme,cjms,cjme, &
2525 cits,cite,ckts,ckte,cjts,cjte, &
2527 nids,nide,nkds,nkde,njds,njde, &
2528 nims,nime,nkms,nkme,njms,njme, &
2529 nits,nite,nkts,nkte,njts,njte, &
2530 shw, & ! stencil half width for interp
2531 imask, & ! interpolation mask
2532 xstag,ystag, & ! staggering of field
2533 ipos,jpos, & ! Position of lower left of nest in CD
2534 nri,nrj, & ! nest ratios
2544 CTEMP_B,NTEMP_B, & ! to be removed
2545 CTEMP_BT,NTEMP_BT, &
2546 CII, IIH, CJJ, JJH, CBWGT1, HBWGT1, & ! south-western grid locs and weights
2547 CBWGT2, HBWGT2, CBWGT3, HBWGT3, & ! note that "C"ourse grid ones are
2548 CBWGT4, HBWGT4, & ! dummys for weights
2554 CETA1,ETA1,CETA2,ETA2 )
2555 USE MODULE_MODEL_CONSTANTS
2559 LOGICAL,INTENT(IN) :: xstag, ystag
2560 REAL, INTENT(INOUT) :: cdt, ndt
2561 INTEGER,INTENT(IN) :: cids, cide, ckds, ckde, cjds, cjde, &
2562 cims, cime, ckms, ckme, cjms, cjme, &
2563 cits, cite, ckts, ckte, cjts, cjte, &
2564 nids, nide, nkds, nkde, njds, njde, &
2565 nims, nime, nkms, nkme, njms, njme, &
2566 nits, nite, nkts, nkte, njts, njte, &
2567 shw,ipos,jpos,nri,nrj
2568 REAL, DIMENSION ( cims:cime, cjms:cjme, ckms:ckme ), INTENT(OUT) :: ctemp_b,ctemp_bt
2569 REAL, DIMENSION ( nims:nime, njms:njme, nkms:nkme ), INTENT(OUT) :: ntemp_b,ntemp_bt
2570 REAL, DIMENSION( * ), INTENT(INOUT) :: c_bxs,n_bxs,c_bxe,n_bxe,c_bys,n_bys,c_bye,n_bye
2571 REAL, DIMENSION( * ), INTENT(INOUT) :: c_btxs,n_btxs,c_btxe,n_btxe,c_btys,n_btys,c_btye,n_btye
2574 INTEGER,DIMENSION(nims:nime,njms:njme), INTENT(IN) :: IMASK
2578 INTEGER,DIMENSION(cims:cime,cjms:cjme), INTENT(IN) :: CII,CJJ ! dummy
2579 REAL,DIMENSION(cims:cime,cjms:cjme), INTENT(IN) :: CBWGT1,CBWGT2
2580 REAL,DIMENSION(cims:cime,cjms:cjme), INTENT(IN) :: CBWGT3,CBWGT4
2581 REAL,DIMENSION(cims:cime,cjms:cjme,ckms:ckme), INTENT(IN) :: CFLD
2582 REAL,DIMENSION(cims:cime,cjms:cjme,ckms:ckme), INTENT(IN) :: CC3d ! scalar input on constant pressure levels
2583 REAL,DIMENSION(ckms:ckme), INTENT(IN) :: CPSTD
2584 REAL,DIMENSION(cims:cime,cjms:cjme), INTENT(IN) :: CPD
2585 REAL,DIMENSION(ckms:ckme), INTENT(IN) :: CETA1,CETA2
2586 REAL, INTENT(IN) :: CPDTOP,CPTOP
2590 INTEGER,DIMENSION(nims:nime,njms:njme), INTENT(IN) :: IIH,JJH
2591 REAL,DIMENSION(nims:nime,njms:njme), INTENT(IN) :: HBWGT1,HBWGT2
2592 REAL,DIMENSION(nims:nime,njms:njme), INTENT(IN) :: HBWGT3,HBWGT4
2593 REAL,DIMENSION(nims:nime,njms:njme,nkms:nkme), INTENT(OUT):: NFLD
2594 REAL,DIMENSION(nims:nime,njms:njme,nkms:nkme), INTENT(OUT):: C3d !Scalar on constant pressure levels
2595 REAL,DIMENSION(nkms:nkme), INTENT(IN) :: PSTD
2596 REAL,DIMENSION(nims:nime,njms:njme), INTENT(IN) :: PD
2597 REAL,DIMENSION(nkms:nkme), INTENT(IN) :: ETA1,ETA2
2598 REAL,INTENT(IN) :: PDTOP,PTOP
2602 INTEGER,PARAMETER :: JTB=134
2603 INTEGER :: I,J,K,II,JJ
2604 REAL,DIMENSION(JTB) :: PIN,CIN,Y2,PIO,PTMP,COUT,DUM1,DUM2
2605 REAL, DIMENSION (nims:nime,njms:njme,nkms:nkme) :: CWK1,CWK2,CWK3,CWK4
2606 !-----------------------------------------------------------------------------------------------------
2609 ! *** CHECK VERTICAL BOUNDS BEFORE USING SPLINE INTERPOLATION
2611 IF(nkme .GT. (JTB-10) .OR. NKDE .GT. (JTB-10)) &
2612 CALL wrf_error_fatal ('mass points: increase JTB in interp_mass_nmm')
2616 NMM_XS: IF(NITS .EQ. NIDS)THEN
2617 ! WRITE(0,*)'ENTERING X1 START BOUNDARY AT T POINTS',NJTS,MIN(NJTE,NJDE-1)
2619 DO K=NKDS,NKDE-1 ! Please note that we are still in isobaric surfaces
2620 DO J = NJTS,MIN(NJTE,NJDE-1)
2621 IF(MOD(JJH(I,J),2) .NE. 0)THEN ! 1,3,5,7 of the parent domain
2622 C3d(I,J,K) = HBWGT1(I,J)*CC3d(IIH(I,J), JJH(I,J) ,K) &
2623 + HBWGT2(I,J)*CC3d(IIH(I,J)+1,JJH(I,J) ,K) &
2624 + HBWGT3(I,J)*CC3d(IIH(I,J), JJH(I,J)-1,K) &
2625 + HBWGT4(I,J)*CC3d(IIH(I,J), JJH(I,J)+1,K)
2627 C3d(I,J,K) = HBWGT1(I,J)*CC3d(IIH(I,J), JJH(I,J) ,K) &
2628 + HBWGT2(I,J)*CC3d(IIH(I,J)+1,JJH(I,J) ,K) &
2629 + HBWGT3(I,J)*CC3d(IIH(I,J)+1,JJH(I,J)-1,K) &
2630 + HBWGT4(I,J)*CC3d(IIH(I,J)+1,JJH(I,J)+1,K)
2635 DO J=NJTS,MIN(NJTE,NJDE-1)
2636 IF(MOD(J,2) .NE. 0)THEN
2637 CIN=0.;PIN=0.;Y2=0;PIO=0.;PTMP=0.;COUT=0.;DUM1=0.;DUM2=0. ! clean up local array
2638 DO K=NKDS+1,NKDE ! inputs at standard levels
2639 PIN(K-1) = EXP((ALOG(PSTD(NKDE-K+1))+ALOG(PSTD(NKDE-K+2)))*0.5)
2640 CIN(K-1) = C3d(I,J,NKDE-K+1)
2644 DO K=NKDS,NKDE ! target points in model interface levels (pint)
2645 PTMP(K) = ETA1(K)*PDTOP + ETA2(K)*PD(I,J) + PTOP
2647 DO K=NKDS,NKDE-1 ! target points in model levels
2648 PIO(K) = EXP((ALOG(PTMP(K))+ALOG(PTMP(K+1)))*0.5)
2650 IF(PTMP(1) .GE. PSTD(1))THEN ! if lower boundary is higher than PMSL(1) re-set lower boundary
2651 PIN(NKDE-1) = PIO(1) ! be consistent with target. This may not happen at all
2652 WRITE(0,*)'WARNING: NESTED DOMAIN PRESSURE AT LOWEST LEVEL HIGHER THAN PSTD'
2653 WRITE(0,*)'I,J,PIO(1),PSTD(1)',I,J,PIO(1),PSTD(1)
2656 CALL SPLINE2(I,J,JTB,NKDE-1,PIN,CIN,Y2,NKDE-1,PIO,COUT,DUM1,DUM2) ! interpolate
2659 CWK1(I,J,K)= COUT(K) ! scalar in the nested domain
2668 DO J = NJTS,MIN(NJTE,NJDE-1)
2670 ntemp_b(i,j,k) = CWK1(I,J,K)
2671 ntemp_bt(i,j,k) = 0.0
2680 NMM_XE: IF(NITE-1 .EQ. NIDE-1)THEN
2681 ! WRITE(0,*)'ENTERING X END BOUNDARY AT T POINTS',NJTS,MIN(NJTE,NJDE-1)
2683 DO K=NKDS,NKDE-1 ! Please note that we are still in isobaric surfaces
2684 DO J = NJTS,MIN(NJTE,NJDE-1)
2685 IF(MOD(JJH(I,J),2) .NE. 0)THEN ! 1,3,5,7 of the parent domain
2686 C3d(I,J,K) = HBWGT1(I,J)*CC3d(IIH(I,J), JJH(I,J) ,K) &
2687 + HBWGT2(I,J)*CC3d(IIH(I,J)+1,JJH(I,J) ,K) &
2688 + HBWGT3(I,J)*CC3d(IIH(I,J), JJH(I,J)-1,K) &
2689 + HBWGT4(I,J)*CC3d(IIH(I,J), JJH(I,J)+1,K)
2691 C3d(I,J,K) = HBWGT1(I,J)*CC3d(IIH(I,J), JJH(I,J) ,K) &
2692 + HBWGT2(I,J)*CC3d(IIH(I,J)+1,JJH(I,J) ,K) &
2693 + HBWGT3(I,J)*CC3d(IIH(I,J)+1,JJH(I,J)-1,K) &
2694 + HBWGT4(I,J)*CC3d(IIH(I,J)+1,JJH(I,J)+1,K)
2699 DO J=NJTS,MIN(NJTE,NJDE-1)
2700 IF(MOD(J,2) .NE. 0)THEN
2701 CIN=0.;PIN=0.;Y2=0;PIO=0.;PTMP=0.;COUT=0.;DUM1=0.;DUM2=0. ! clean up local array
2702 DO K=NKDS+1,NKDE ! inputs at standard levels
2703 PIN(K-1) = EXP((ALOG(PSTD(NKDE-K+1))+ALOG(PSTD(NKDE-K+2)))*0.5)
2704 CIN(K-1) = C3d(I,J,NKDE-K+1)
2708 DO K=NKDS,NKDE ! target points in model interface levels (pint)
2709 PTMP(K) = ETA1(K)*PDTOP + ETA2(K)*PD(I,J) + PTOP
2711 DO K=NKDS,NKDE-1 ! target points in model levels
2712 PIO(K) = EXP((ALOG(PTMP(K))+ALOG(PTMP(K+1)))*0.5)
2714 IF(PTMP(1) .GE. PSTD(1))THEN ! if lower boundary is higher than PMSL(1) re-set lower boundary
2715 PIN(NKDE-1) = PIO(1) ! be consistent with target. This may not happen at all
2716 WRITE(0,*)'WARNING: NESTED DOMAIN PRESSURE AT LOWEST LEVEL HIGHER THAN PSTD'
2717 WRITE(0,*)'I,J,PIO(1),PSTD(1)',I,J,PIO(1),PSTD(1)
2720 CALL SPLINE2(I,J,JTB,NKDE-1,PIN,CIN,Y2,NKDE-1,PIO,COUT,DUM1,DUM2) ! interpolate
2723 CWK2(I,J,K)= COUT(K) ! scalar in the nested domain
2732 DO J = NJTS,MIN(NJTE,NJDE-1)
2733 DO K = NKDS,MIN(NKTE,NKDE-1)
2734 ntemp_b(i,j,k) = CWK2(I,J,K)
2735 ntemp_bt(i,j,k) = 0.0
2743 NMM_YS: IF(NJTS .EQ. NJDS)THEN
2744 ! WRITE(0,*)'ENTERING Y START BOUNDARY AT T POINTS',NITS,MIN(NITE,NIDE-1)
2747 DO I = NITS,MIN(NITE,NIDE-1)
2748 IF(MOD(JJH(I,J),2) .NE. 0)THEN ! 1,3,5,7 of the parent domain
2749 C3d(I,J,K) = HBWGT1(I,J)*CC3d(IIH(I,J), JJH(I,J) ,K) &
2750 + HBWGT2(I,J)*CC3d(IIH(I,J)+1,JJH(I,J) ,K) &
2751 + HBWGT3(I,J)*CC3d(IIH(I,J), JJH(I,J)-1,K) &
2752 + HBWGT4(I,J)*CC3d(IIH(I,J), JJH(I,J)+1,K)
2754 C3d(I,J,K) = HBWGT1(I,J)*CC3d(IIH(I,J), JJH(I,J) ,K) &
2755 + HBWGT2(I,J)*CC3d(IIH(I,J)+1,JJH(I,J) ,K) &
2756 + HBWGT3(I,J)*CC3d(IIH(I,J)+1,JJH(I,J)-1,K) &
2757 + HBWGT4(I,J)*CC3d(IIH(I,J)+1,JJH(I,J)+1,K)
2763 DO I=NITS,MIN(NITE,NIDE-1)
2764 CIN=0.;PIN=0.;Y2=0;PIO=0.;PTMP=0.;COUT=0.;DUM1=0.;DUM2=0. ! clean up local array
2765 DO K=NKDS+1,NKDE ! inputs at standard levels
2766 PIN(K-1) = EXP((ALOG(PSTD(NKDE-K+1))+ALOG(PSTD(NKDE-K+2)))*0.5)
2767 CIN(K-1) = C3d(I,J,NKDE-K+1)
2771 DO K=NKDS,NKDE ! target points in model interface levels (pint)
2772 PTMP(K) = ETA1(K)*PDTOP + ETA2(K)*PD(I,J) + PTOP
2774 DO K=NKDS,NKDE-1 ! target points in model levels
2775 PIO(K) = EXP((ALOG(PTMP(K))+ALOG(PTMP(K+1)))*0.5)
2777 IF(PTMP(1) .GE. PSTD(1))THEN ! if lower boundary is higher than PMSL(1) re-set lower boundary
2778 PIN(NKDE-1) = PIO(1) ! be consistent with target. This may not happen at all
2779 WRITE(0,*)'WARNING: NESTED DOMAIN PRESSURE AT LOWEST LEVEL HIGHER THAN PSTD'
2780 WRITE(0,*)'I,J,PIO(1),PSTD(1)',I,J,PIO(1),PSTD(1)
2783 CALL SPLINE2(I,J,JTB,NKDE-1,PIN,CIN,Y2,NKDE-1,PIO,COUT,DUM1,DUM2) ! interpolate
2786 CWK3(I,J,K)= COUT(K) ! scalar in the nested domain
2791 DO I = NITS,MIN(NITE,NIDE-1)
2792 ntemp_b(i,J,K) = CWK3(I,J,K)
2793 ntemp_bt(i,J,K) = 0.0
2801 NMM_YE: IF(NJTE-1 .EQ. NJDE-1)THEN
2802 ! WRITE(0,*)'ENTERING Y END BOUNDARY AT T POINTS',NITS,MIN(NITE,NIDE-1)
2805 DO I = NITS,MIN(NITE,NIDE-1)
2806 IF(MOD(JJH(I,J),2) .NE. 0)THEN ! 1,3,5,7 of the parent domain
2807 C3d(I,J,K) = HBWGT1(I,J)*CC3d(IIH(I,J), JJH(I,J) ,K) &
2808 + HBWGT2(I,J)*CC3d(IIH(I,J)+1,JJH(I,J) ,K) &
2809 + HBWGT3(I,J)*CC3d(IIH(I,J), JJH(I,J)-1,K) &
2810 + HBWGT4(I,J)*CC3d(IIH(I,J), JJH(I,J)+1,K)
2812 C3d(I,J,K) = HBWGT1(I,J)*CC3d(IIH(I,J), JJH(I,J) ,K) &
2813 + HBWGT2(I,J)*CC3d(IIH(I,J)+1,JJH(I,J) ,K) &
2814 + HBWGT3(I,J)*CC3d(IIH(I,J)+1,JJH(I,J)-1,K) &
2815 + HBWGT4(I,J)*CC3d(IIH(I,J)+1,JJH(I,J)+1,K)
2821 DO I=NITS,MIN(NITE,NIDE-1)
2822 CIN=0.;PIN=0.;Y2=0;PIO=0.;PTMP=0.;COUT=0.;DUM1=0.;DUM2=0. ! clean up local array
2823 DO K=NKDS+1,NKDE ! inputs at standard levels
2824 PIN(K-1) = EXP((ALOG(PSTD(NKDE-K+1))+ALOG(PSTD(NKDE-K+2)))*0.5)
2825 CIN(K-1) = C3d(I,J,NKDE-K+1)
2829 DO K=NKDS,NKDE ! target points in model interface levels (pint)
2830 PTMP(K) = ETA1(K)*PDTOP + ETA2(K)*PD(I,J) + PTOP
2832 DO K=NKDS,NKDE-1 ! target points in model levels
2833 PIO(K) = EXP((ALOG(PTMP(K))+ALOG(PTMP(K+1)))*0.5)
2835 IF(PTMP(1) .GE. PSTD(1))THEN ! if lower boundary is higher than PMSL(1) re-set lower boundary
2836 PIN(NKDE-1) = PIO(1) ! be consistent with target. This may not happen at all
2837 WRITE(0,*)'WARNING: NESTED DOMAIN PRESSURE AT LOWEST LEVEL HIGHER THAN PSTD'
2838 WRITE(0,*)'I,J,PIO(1),PSTD(1)',I,J,PIO(1),PSTD(1)
2841 CALL SPLINE2(I,J,JTB,NKDE-1,PIN,CIN,Y2,NKDE-1,PIO,COUT,DUM1,DUM2) ! interpolate
2844 CWK4(I,J,K)= COUT(K) ! scalar in the nested domain
2849 DO I = NITS,MIN(NITE,NIDE-1)
2850 ntemp_b(i,J,K) = CWK4(I,J,K)
2851 ntemp_bt(i,J,K) = 0.0
2857 END SUBROUTINE nmm_bdy_scalar
2860 !=======================================================================================
2861 SUBROUTINE SPLINE2(I,J,JTBX,NOLD,XOLD,YOLD,Y2,NNEW,XNEW,YNEW,P,Q)
2863 ! ******************************************************************
2865 ! * THIS IS A ONE-DIMENSIONAL CUBIC SPLINE FITTING ROUTINE *
2866 ! * PROGRAMED FOR A SMALL SCALAR MACHINE. *
2868 ! * PROGRAMER Z. JANJIC *
2870 ! * NOLD - NUMBER OF GIVEN VALUES OF THE FUNCTION. MUST BE GE 3. *
2871 ! * XOLD - LOCATIONS OF THE POINTS AT WHICH THE VALUES OF THE *
2872 ! * FUNCTION ARE GIVEN. MUST BE IN ASCENDING ORDER. *
2873 ! * YOLD - THE GIVEN VALUES OF THE FUNCTION AT THE POINTS XOLD. *
2874 ! * Y2 - THE SECOND DERIVATIVES AT THE POINTS XOLD. IF NATURAL *
2875 ! * SPLINE IS FITTED Y2(1)=0. AND Y2(NOLD)=0. MUST BE *
2877 ! * NNEW - NUMBER OF VALUES OF THE FUNCTION TO BE CALCULATED. *
2878 ! * XNEW - LOCATIONS OF THE POINTS AT WHICH THE VALUES OF THE *
2879 ! * FUNCTION ARE CALCULATED. XNEW(K) MUST BE GE XOLD(1) *
2880 ! * AND LE XOLD(NOLD). *
2881 ! * YNEW - THE VALUES OF THE FUNCTION TO BE CALCULATED. *
2882 ! * P, Q - AUXILIARY VECTORS OF THE LENGTH NOLD-2. *
2884 ! ******************************************************************
2885 !---------------------------------------------------------------------
2887 !---------------------------------------------------------------------
2888 INTEGER,INTENT(IN) :: I,J,JTBX,NNEW,NOLD
2889 REAL,DIMENSION(JTBX),INTENT(IN) :: XNEW,XOLD,YOLD
2890 REAL,DIMENSION(JTBX),INTENT(INOUT) :: P,Q,Y2
2891 REAL,DIMENSION(JTBX),INTENT(OUT) :: YNEW
2893 INTEGER :: II,JJ,K,K1,K2,KOLD,NOLDM1
2894 REAL :: AK,BK,CK,DEN,DX,DXC,DXL,DXR,DYDXL,DYDXR &
2895 ,RDX,RTDXC,X,XK,XSQ,Y2K,Y2KP1
2896 !---------------------------------------------------------------------
2902 IF(I.eq.II.and.J.eq.JJ)THEN
2903 WRITE(0,*)'DEBUG in SPLINE2: I,J',I,J
2904 WRITE(0,*)'DEBUG in SPLINE2:HSO= ',xnew(1:nold)
2906 WRITE(0,*)'DEBUG in SPLINE2:L,ZETAI,PINTI= ' &
2915 DYDXL=(YOLD(2)-YOLD(1))/DXL
2916 DYDXR=(YOLD(3)-YOLD(2))/DXR
2919 P(1)= RTDXC*(6.*(DYDXR-DYDXL)-DXL*Y2(1))
2922 IF(NOLD.EQ.3)GO TO 150
2923 !---------------------------------------------------------------------
2928 DXR=XOLD(K+1)-XOLD(K)
2929 DYDXR=(YOLD(K+1)-YOLD(K))/DXR
2931 DEN=1./(DXL*Q(K-2)+DXC+DXC)
2933 P(K-1)= DEN*(6.*(DYDXR-DYDXL)-DXL*P(K-2))
2937 IF(K.LT.NOLD)GO TO 100
2938 !-----------------------------------------------------------------------
2941 200 Y2(K)=P(K-1)+Q(K-1)*Y2(K+1)
2945 !-----------------------------------------------------------------------
2952 IF(XOLD(K2).GT.XK)THEN
2962 450 IF(K1.EQ.1)GO TO 500
2963 IF(K.EQ.KOLD)GO TO 550
2969 DX=XOLD(K+1)-XOLD(K)
2972 AK=.1666667*RDX*(Y2KP1-Y2K)
2974 CK=RDX*(YOLD(K+1)-YOLD(K))-.1666667*DX*(Y2KP1+Y2K+Y2K)
2979 YNEW(K1)=AK*XSQ*X+BK*XSQ+CK*X+YOLD(K)
2983 IF(I.eq.II.and.J.eq.JJ)THEN
2984 WRITE(0,*) 'DEBUG:: k1,xnew(k1),ynew(k1): ', K1,XNEW(k1),YNEW(k1)
2989 IF(K1.LE.NNEW)GO TO 300
2993 END SUBROUTINE SPLINE2
2995 !=======================================================================================
2996 ! E grid interpolation for H and V points
2997 !=======================================================================================
2999 SUBROUTINE interp_h_nmm (cfld, & ! CD field
3000 cids, cide, ckds, ckde, cjds, cjde, &
3001 cims, cime, ckms, ckme, cjms, cjme, &
3002 cits, cite, ckts, ckte, cjts, cjte, &
3004 nids, nide, nkds, nkde, njds, njde, &
3005 nims, nime, nkms, nkme, njms, njme, &
3006 nits, nite, nkts, nkte, njts, njte, &
3007 shw, & ! stencil half width for interp
3008 imask, & ! interpolation mask
3009 xstag, ystag, & ! staggering of field
3010 ipos, jpos, & ! Position of lower left of nest in CD
3011 nri, nrj, & ! nest ratios
3012 CII, IIH, CJJ, JJH, CBWGT1, HBWGT1, & ! south-western grid locs and weights
3013 CBWGT2, HBWGT2, CBWGT3, HBWGT3, & ! note that "C"ourse grid ones are
3014 CBWGT4, HBWGT4 ) ! dummys for weights
3018 INTEGER, INTENT(IN) :: cids, cide, ckds, ckde, cjds, cjde, &
3019 cims, cime, ckms, ckme, cjms, cjme, &
3020 cits, cite, ckts, ckte, cjts, cjte, &
3021 nids, nide, nkds, nkde, njds, njde, &
3022 nims, nime, nkms, nkme, njms, njme, &
3023 nits, nite, nkts, nkte, njts, njte, &
3027 LOGICAL, INTENT(IN) :: xstag, ystag
3029 REAL, DIMENSION ( cims:cime, cjms:cjme, ckms:ckme ) :: cfld
3030 REAL, DIMENSION ( nims:nime, njms:njme, nkms:nkme ) :: nfld
3031 REAL, DIMENSION ( cims:cime, cjms:cjme ), INTENT(IN) :: CBWGT1,CBWGT2,CBWGT3,CBWGT4 ! dummy
3032 REAL, DIMENSION ( nims:nime, njms:njme ), INTENT(IN) :: HBWGT1,HBWGT2,HBWGT3,HBWGT4
3033 INTEGER, DIMENSION ( cims:cime, cjms:cjme ), INTENT(IN) :: CII,CJJ ! dummy
3034 INTEGER, DIMENSION ( nims:nime, njms:njme ), INTENT(IN) :: IIH,JJH
3035 INTEGER, DIMENSION ( nims:nime, njms:njme ) :: imask
3040 !*** CHECK DOMAIN BOUNDS BEFORE INTERPOLATION
3042 DO J=NJTS,MIN(NJTE,NJDE-1)
3043 DO I=NITS,MIN(NITE,NIDE-1)
3044 IF(IIH(i,j).LT.(CIDS-shw) .OR. IIH(i,j).GT.(CIDE+shw)) &
3045 CALL wrf_error_fatal ('hpoints:check domain bounds along x' )
3046 IF(JJH(i,j).LT.(CJDS-shw) .OR. JJH(i,j).GT.(CJDE+shw)) &
3047 CALL wrf_error_fatal ('hpoints:check domain bounds along y' )
3051 !*** INDEX CONVENTIONS
3066 DO J=NJTS,MIN(NJTE,NJDE-1)
3067 DO I=NITS,MIN(NITE,NIDE-1)
3068 IF(IMASK(I,J) .NE. 1)THEN
3070 IF(MOD(JJH(I,J),2) .NE. 0)THEN ! 1,3,5,7
3071 NFLD(I,J,K) = HBWGT1(I,J)*CFLD(IIH(I,J), JJH(I,J) ,K) &
3072 + HBWGT2(I,J)*CFLD(IIH(I,J)+1,JJH(I,J) ,K) &
3073 + HBWGT3(I,J)*CFLD(IIH(I,J), JJH(I,J)-1,K) &
3074 + HBWGT4(I,J)*CFLD(IIH(I,J), JJH(I,J)+1,K)
3076 NFLD(I,J,K) = HBWGT1(I,J)*CFLD(IIH(I,J), JJH(I,J) ,K) &
3077 + HBWGT2(I,J)*CFLD(IIH(I,J)+1,JJH(I,J) ,K) &
3078 + HBWGT3(I,J)*CFLD(IIH(I,J)+1,JJH(I,J)-1,K) &
3079 + HBWGT4(I,J)*CFLD(IIH(I,J)+1,JJH(I,J)+1,K)
3087 END SUBROUTINE interp_h_nmm
3089 SUBROUTINE interp_v_nmm (cfld, & ! CD field
3090 cids, cide, ckds, ckde, cjds, cjde, &
3091 cims, cime, ckms, ckme, cjms, cjme, &
3092 cits, cite, ckts, ckte, cjts, cjte, &
3094 nids, nide, nkds, nkde, njds, njde, &
3095 nims, nime, nkms, nkme, njms, njme, &
3096 nits, nite, nkts, nkte, njts, njte, &
3097 shw, & ! stencil half width for interp
3098 imask, & ! interpolation mask
3099 xstag, ystag, & ! staggering of field
3100 ipos, jpos, & ! Position of lower left of nest in CD
3101 nri, nrj, & ! nest ratios
3102 CII, IIV, CJJ, JJV, CBWGT1, VBWGT1, & ! south-western grid locs and weights
3103 CBWGT2, VBWGT2, CBWGT3, VBWGT3, & ! note that "C"ourse grid ones are
3104 CBWGT4, VBWGT4 ) ! dummys
3108 INTEGER, INTENT(IN) :: cids, cide, ckds, ckde, cjds, cjde, &
3109 cims, cime, ckms, ckme, cjms, cjme, &
3110 cits, cite, ckts, ckte, cjts, cjte, &
3111 nids, nide, nkds, nkde, njds, njde, &
3112 nims, nime, nkms, nkme, njms, njme, &
3113 nits, nite, nkts, nkte, njts, njte, &
3117 LOGICAL, INTENT(IN) :: xstag, ystag
3119 REAL, DIMENSION ( cims:cime, cjms:cjme, ckms:ckme ) :: cfld
3120 REAL, DIMENSION ( nims:nime, njms:njme, nkms:nkme ) :: nfld
3121 REAL, DIMENSION ( cims:cime, cjms:cjme ), INTENT(IN) :: CBWGT1,CBWGT2,CBWGT3,CBWGT4 ! dummy
3122 REAL, DIMENSION ( nims:nime, njms:njme ), INTENT(IN) :: VBWGT1,VBWGT2,VBWGT3,VBWGT4
3123 INTEGER, DIMENSION ( cims:cime, cjms:cjme ), INTENT(IN) :: CII,CJJ ! dummy
3124 INTEGER, DIMENSION ( nims:nime, njms:njme ), INTENT(IN) :: IIV,JJV
3125 INTEGER, DIMENSION ( nims:nime, njms:njme ) :: imask
3132 !*** CHECK DOMAIN BOUNDS BEFORE INTERPOLATION
3134 DO J=NJTS,MIN(NJTE,NJDE-1)
3135 DO I=NITS,MIN(NITE,NIDE-1)
3136 IF(IIV(i,j).LT.(CIDS-shw) .OR. IIV(i,j).GT.(CIDE+shw)) &
3137 CALL wrf_error_fatal ('vpoints:check domain bounds along x' )
3138 IF(JJV(i,j).LT.(CJDS-shw) .OR. JJV(i,j).GT.(CJDE+shw)) &
3139 CALL wrf_error_fatal ('vpoints:check domain bounds along y' )
3143 !*** INDEX CONVENTIONS
3158 DO J=NJTS,MIN(NJTE,NJDE-1)
3159 DO I=NITS,MIN(NITE,NIDE-1)
3160 IF(IMASK(I,J) .NE. 1)THEN
3162 IF(MOD(JJV(I,J),2) .NE. 0)THEN ! 1,3,5,7
3163 NFLD(I,J,K) = VBWGT1(I,J)*CFLD(IIV(I,J), JJV(I,J) ,K) &
3164 + VBWGT2(I,J)*CFLD(IIV(I,J)+1,JJV(I,J) ,K) &
3165 + VBWGT3(I,J)*CFLD(IIV(I,J)+1,JJV(I,J)-1,K) &
3166 + VBWGT4(I,J)*CFLD(IIV(I,J)+1,JJV(I,J)+1,K)
3168 NFLD(I,J,K) = VBWGT1(I,J)*CFLD(IIV(I,J), JJV(I,J) ,K) &
3169 + VBWGT2(I,J)*CFLD(IIV(I,J)+1,JJV(I,J) ,K) &
3170 + VBWGT3(I,J)*CFLD(IIV(I,J), JJV(I,J)-1,K) &
3171 + VBWGT4(I,J)*CFLD(IIV(I,J), JJV(I,J)+1,K)
3179 END SUBROUTINE interp_v_nmm
3181 !=======================================================================================
3182 ! E grid nearest neighbour interpolation for H points.
3183 ! This routine assumes cfld and nfld are in IJK
3184 !=======================================================================================
3186 SUBROUTINE interp_hnear_nmm (cfld, & ! CD field
3187 cids, cide, ckds, ckde, cjds, cjde, &
3188 cims, cime, ckms, ckme, cjms, cjme, &
3189 cits, cite, ckts, ckte, cjts, cjte, &
3191 nids, nide, nkds, nkde, njds, njde, &
3192 nims, nime, nkms, nkme, njms, njme, &
3193 nits, nite, nkts, nkte, njts, njte, &
3194 shw, & ! stencil half width for interp
3195 imask, & ! interpolation mask
3196 xstag, ystag, & ! staggering of field
3197 ipos, jpos, & ! Position of lower left of nest in CD
3198 nri, nrj, & ! nest ratios
3199 CII, IIH, CJJ, JJH, CBWGT1, HBWGT1, & ! south-western grid locs and weights
3200 CBWGT2, HBWGT2, CBWGT3, HBWGT3, & ! note that "C"ourse grid ones are
3201 CBWGT4, HBWGT4 ) ! just dummys
3205 INTEGER, INTENT(IN) :: cids, cide, ckds, ckde, cjds, cjde, &
3206 cims, cime, ckms, ckme, cjms, cjme, &
3207 cits, cite, ckts, ckte, cjts, cjte, &
3208 nids, nide, nkds, nkde, njds, njde, &
3209 nims, nime, nkms, nkme, njms, njme, &
3210 nits, nite, nkts, nkte, njts, njte, &
3214 LOGICAL, INTENT(IN) :: xstag, ystag
3216 REAL, DIMENSION ( cims:cime, cjms:cjme, ckms:ckme ) :: cfld
3217 REAL, DIMENSION ( nims:nime, njms:njme, nkms:nkme ) :: nfld
3218 REAL, DIMENSION ( cims:cime, cjms:cjme ), INTENT(IN) :: CBWGT1,CBWGT2,CBWGT3,CBWGT4 ! dummy
3219 REAL, DIMENSION ( nims:nime, njms:njme ), INTENT(IN) :: HBWGT1,HBWGT2,HBWGT3,HBWGT4
3220 INTEGER, DIMENSION ( cims:cime, cjms:cjme ), INTENT(IN) :: CII,CJJ ! dummy
3221 INTEGER, DIMENSION ( nims:nime, njms:njme ), INTENT(IN) :: IIH,JJH
3222 INTEGER, DIMENSION ( nims:nime, njms:njme ) :: imask
3229 REAL, DIMENSION (4, nims:nime, njms:njme ) :: NBWGT
3232 !*** INDEX CONVENTIONS
3240 !*** NBWGT1=1 NBWGT2=0
3246 DO J=NJTS,MIN(NJTE,NJDE-1)
3247 DO I=NITS,MIN(NITE,NIDE-1)
3248 IF(IMASK(I,J) .NE. 1)THEN
3249 NBWGT(1,I,J)=HBWGT1(I,J)
3250 NBWGT(2,I,J)=HBWGT2(I,J)
3251 NBWGT(3,I,J)=HBWGT3(I,J)
3252 NBWGT(4,I,J)=HBWGT4(I,J)
3257 DO J=NJTS,MIN(NJTE,NJDE-1)
3258 DO I=NITS,MIN(NITE,NIDE-1)
3259 IF(IMASK(I,J) .NE. 1)THEN
3263 AMAXVAL=amax1(NBWGT(N,I,J),AMAXVAL)
3269 IF(AMAXVAL .EQ. NBWGT(N,I,J) .AND. FLIP)THEN
3275 SUM=SUM+NBWGT(N,I,J)
3276 IF(SUM .GT. 1.0)CALL wrf_error_fatal ( "horizontal interp error - interp_hnear_nmm" )
3284 DO J=NJTS,MIN(NJTE,NJDE-1)
3285 DO I=NITS,MIN(NITE,NIDE-1)
3286 IF(IMASK(I,J) .NE. 1)THEN
3287 IF(MOD(JJH(I,J),2) .NE. 0)THEN ! 1,3,5,7
3288 NFLD(I,J,K) = NBWGT(1,I,J)*CFLD(IIH(I,J), JJH(I,J) ,K) &
3289 + NBWGT(2,I,J)*CFLD(IIH(I,J)+1,JJH(I,J) ,K) &
3290 + NBWGT(3,I,J)*CFLD(IIH(I,J), JJH(I,J)-1,K) &
3291 + NBWGT(4,I,J)*CFLD(IIH(I,J), JJH(I,J)+1,K)
3293 NFLD(I,J,K) = NBWGT(1,I,J)*CFLD(IIH(I,J), JJH(I,J) ,K) &
3294 + NBWGT(2,I,J)*CFLD(IIH(I,J)+1,JJH(I,J) ,K) &
3295 + NBWGT(3,I,J)*CFLD(IIH(I,J)+1,JJH(I,J)-1,K) &
3296 + NBWGT(4,I,J)*CFLD(IIH(I,J)+1,JJH(I,J)+1,K)
3303 END SUBROUTINE interp_hnear_nmm
3305 !=======================================================================================
3306 ! E grid nearest neighbour interpolation for H points.
3307 ! This routine assumes cfld and nfld are in IKJ or ILJ
3308 !=======================================================================================
3310 SUBROUTINE interp_hnear_ikj_nmm (cfld, & ! CD field
3311 cids, cide, ckds, ckde, cjds, cjde, &
3312 cims, cime, ckms, ckme, cjms, cjme, &
3313 cits, cite, ckts, ckte, cjts, cjte, &
3315 nids, nide, nkds, nkde, njds, njde, &
3316 nims, nime, nkms, nkme, njms, njme, &
3317 nits, nite, nkts, nkte, njts, njte, &
3318 shw, & ! stencil half width for interp
3319 imask, & ! interpolation mask
3320 xstag, ystag, & ! staggering of field
3321 ipos, jpos, & ! Position of lower left of nest in CD
3322 nri, nrj, & ! nest ratios
3323 CII, IIH, CJJ, JJH, CBWGT1, HBWGT1, & ! south-western grid locs and weights
3324 CBWGT2, HBWGT2, CBWGT3, HBWGT3, & ! note that "C"ourse grid ones are
3325 CBWGT4, HBWGT4 ) ! just dummys
3329 INTEGER, INTENT(IN) :: cids, cide, ckds, ckde, cjds, cjde, &
3330 cims, cime, ckms, ckme, cjms, cjme, &
3331 cits, cite, ckts, ckte, cjts, cjte, &
3332 nids, nide, nkds, nkde, njds, njde, &
3333 nims, nime, nkms, nkme, njms, njme, &
3334 nits, nite, nkts, nkte, njts, njte, &
3338 LOGICAL, INTENT(IN) :: xstag, ystag
3340 REAL, DIMENSION ( cims:cime, ckms:ckme, cjms:cjme ) :: cfld
3341 REAL, DIMENSION ( nims:nime, nkms:nkme, njms:njme ) :: nfld
3342 REAL, DIMENSION ( cims:cime, cjms:cjme ), INTENT(IN) :: CBWGT1,CBWGT2,CBWGT3,CBWGT4 ! dummy
3343 REAL, DIMENSION ( nims:nime, njms:njme ), INTENT(IN) :: HBWGT1,HBWGT2,HBWGT3,HBWGT4
3344 INTEGER, DIMENSION ( cims:cime, cjms:cjme ), INTENT(IN) :: CII,CJJ ! dummy
3345 INTEGER, DIMENSION ( nims:nime, njms:njme ), INTENT(IN) :: IIH,JJH
3346 INTEGER, DIMENSION ( nims:nime, njms:njme ) :: imask
3353 REAL, DIMENSION (4, nims:nime, njms:njme ) :: NBWGT
3356 !*** INDEX CONVENTIONS
3364 !*** NBWGT1=1 NBWGT2=0
3370 DO J=NJTS,MIN(NJTE,NJDE-1)
3371 DO I=NITS,MIN(NITE,NIDE-1)
3372 IF(IMASK(I,J) .NE. 1)THEN
3373 NBWGT(1,I,J)=HBWGT1(I,J)
3374 NBWGT(2,I,J)=HBWGT2(I,J)
3375 NBWGT(3,I,J)=HBWGT3(I,J)
3376 NBWGT(4,I,J)=HBWGT4(I,J)
3381 DO J=NJTS,MIN(NJTE,NJDE-1)
3382 DO I=NITS,MIN(NITE,NIDE-1)
3383 IF(IMASK(I,J) .NE. 1)THEN
3387 AMAXVAL=amax1(NBWGT(N,I,J),AMAXVAL)
3393 IF(AMAXVAL .EQ. NBWGT(N,I,J) .AND. FLIP)THEN
3399 SUM=SUM+NBWGT(N,I,J)
3400 IF(SUM .GT. 1.0)CALL wrf_error_fatal ( "horizontal interp error - interp_hnear_nmm" )
3407 DO J=NJTS,MIN(NJTE,NJDE-1)
3409 DO I=NITS,MIN(NITE,NIDE-1)
3410 IF(IMASK(I,J) .NE. 1)THEN
3411 IF(MOD(JJH(I,J),2) .NE. 0)THEN ! 1,3,5,7
3412 NFLD(I,K,J) = NBWGT(1,I,J)*CFLD(IIH(I,J), K,JJH(I,J) ) &
3413 + NBWGT(2,I,J)*CFLD(IIH(I,J)+1,K,JJH(I,J) ) &
3414 + NBWGT(3,I,J)*CFLD(IIH(I,J), K,JJH(I,J)-1) &
3415 + NBWGT(4,I,J)*CFLD(IIH(I,J), K,JJH(I,J)+1)
3417 NFLD(I,K,J) = NBWGT(1,I,J)*CFLD(IIH(I,J), K,JJH(I,J) ) &
3418 + NBWGT(2,I,J)*CFLD(IIH(I,J)+1,K,JJH(I,J) ) &
3419 + NBWGT(3,I,J)*CFLD(IIH(I,J)+1,K,JJH(I,J)-1) &
3420 + NBWGT(4,I,J)*CFLD(IIH(I,J)+1,K,JJH(I,J)+1)
3427 END SUBROUTINE interp_hnear_ikj_nmm
3429 !=======================================================================================
3430 ! E grid nearest neighbour interpolation for integer H points
3431 !=======================================================================================
3433 SUBROUTINE interp_int_hnear_nmm (cfld, & ! CD field; integers
3434 cids, cide, ckds, ckde, cjds, cjde, &
3435 cims, cime, ckms, ckme, cjms, cjme, &
3436 cits, cite, ckts, ckte, cjts, cjte, &
3437 nfld, & ! ND field; integers
3438 nids, nide, nkds, nkde, njds, njde, &
3439 nims, nime, nkms, nkme, njms, njme, &
3440 nits, nite, nkts, nkte, njts, njte, &
3441 shw, & ! stencil half width for interp
3442 imask, & ! interpolation mask
3443 xstag, ystag, & ! staggering of field
3444 ipos, jpos, & ! lower left of nest in CD
3445 nri, nrj, & ! nest ratios
3446 CII, IIH, CJJ, JJH, CBWGT1, HBWGT1, & ! s-w grid locs and weights
3447 CBWGT2, HBWGT2, CBWGT3, HBWGT3, & ! note that "C"ourse grid ones are
3448 CBWGT4, HBWGT4 ) ! just dummys
3452 INTEGER, INTENT(IN) :: cids, cide, ckds, ckde, cjds, cjde, &
3453 cims, cime, ckms, ckme, cjms, cjme, &
3454 cits, cite, ckts, ckte, cjts, cjte, &
3455 nids, nide, nkds, nkde, njds, njde, &
3456 nims, nime, nkms, nkme, njms, njme, &
3457 nits, nite, nkts, nkte, njts, njte, &
3461 LOGICAL, INTENT(IN) :: xstag, ystag
3463 INTEGER, DIMENSION ( cims:cime, cjms:cjme, ckms:ckme ) :: cfld
3464 INTEGER, DIMENSION ( nims:nime, njms:njme, nkms:nkme ) :: nfld
3465 REAL, DIMENSION ( cims:cime, cjms:cjme ), INTENT(IN) :: CBWGT1,CBWGT2,CBWGT3,CBWGT4 ! dummy
3466 REAL, DIMENSION ( nims:nime, njms:njme ), INTENT(IN) :: HBWGT1,HBWGT2,HBWGT3,HBWGT4
3467 INTEGER, DIMENSION ( cims:cime, cjms:cjme ), INTENT(IN) :: CII,CJJ ! dummy
3468 INTEGER, DIMENSION ( nims:nime, njms:njme ), INTENT(IN) :: IIH,JJH
3469 INTEGER, DIMENSION ( nims:nime, njms:njme ) :: imask
3476 REAL, DIMENSION (4, nims:nime, njms:njme ) :: NBWGT
3479 !*** INDEX CONVENTIONS
3487 !*** NBWGT1=1 NBWGT2=0
3493 DO J=NJTS,MIN(NJTE,NJDE-1)
3494 DO I=NITS,MIN(NITE,NIDE-1)
3495 IF(IMASK(I,J) .NE. 1)THEN
3496 NBWGT(1,I,J)=HBWGT1(I,J)
3497 NBWGT(2,I,J)=HBWGT2(I,J)
3498 NBWGT(3,I,J)=HBWGT3(I,J)
3499 NBWGT(4,I,J)=HBWGT4(I,J)
3504 DO J=NJTS,MIN(NJTE,NJDE-1)
3505 DO I=NITS,MIN(NITE,NIDE-1)
3506 IF(IMASK(I,J) .NE. 1)THEN
3510 AMAXVAL=amax1(NBWGT(N,I,J),AMAXVAL)
3516 IF(AMAXVAL .EQ. NBWGT(N,I,J) .AND. FLIP)THEN
3522 SUM=SUM+NBWGT(N,I,J)
3523 IF(SUM .GT. 1.0)CALL wrf_error_fatal ( "horizontal interp error - interp_hnear_nmm" )
3530 DO J=NJTS,MIN(NJTE,NJDE-1)
3532 DO I=NITS,MIN(NITE,NIDE-1)
3533 IF(IMASK(I,J) .NE. 1)THEN
3534 IF(MOD(JJH(I,J),2) .NE. 0)THEN ! 1,3,5,7
3535 NFLD(I,J,K) = NBWGT(1,I,J)*CFLD(IIH(I,J), JJH(I,J) ,K) &
3536 + NBWGT(2,I,J)*CFLD(IIH(I,J)+1,JJH(I,J) ,K) &
3537 + NBWGT(3,I,J)*CFLD(IIH(I,J), JJH(I,J)-1,K) &
3538 + NBWGT(4,I,J)*CFLD(IIH(I,J), JJH(I,J)+1,K)
3540 NFLD(I,J,K) = NBWGT(1,I,J)*CFLD(IIH(I,J), JJH(I,J) ,K) &
3541 + NBWGT(2,I,J)*CFLD(IIH(I,J)+1,JJH(I,J) ,K) &
3542 + NBWGT(3,I,J)*CFLD(IIH(I,J)+1,JJH(I,J)-1,K) &
3543 + NBWGT(4,I,J)*CFLD(IIH(I,J)+1,JJH(I,J)+1,K)
3550 END SUBROUTINE interp_int_hnear_nmm
3552 !--------------------------------------------------------------------------------------
3554 SUBROUTINE nmm_bdy_hinterp (cfld, & ! CD field
3555 cids, cide, ckds, ckde, cjds, cjde, &
3556 cims, cime, ckms, ckme, cjms, cjme, &
3557 cits, cite, ckts, ckte, cjts, cjte, &
3559 nids, nide, nkds, nkde, njds, njde, &
3560 nims, nime, nkms, nkme, njms, njme, &
3561 nits, nite, nkts, nkte, njts, njte, &
3562 shw, & ! stencil half width
3563 imask, & ! interpolation mask
3564 xstag, ystag, & ! staggering of field
3565 ipos, jpos, & ! Position of lower left of nest in CD
3566 nri, nrj, & ! nest ratios
3575 CTEMP_B,NTEMP_B, & ! These temp arrays should be removed
3576 CTEMP_BT,NTEMP_BT, & ! later on
3577 CII, IIH, CJJ, JJH, CBWGT1, HBWGT1, & ! south-western grid locs and weights
3578 CBWGT2, HBWGT2, CBWGT3, HBWGT3, & ! note that "C"ourse grid ones are
3579 CBWGT4, HBWGT4 ) ! dummys
3581 ! use module_state_description
3582 USE module_configure
3583 USE module_wrf_error
3588 INTEGER, INTENT(IN) :: cids, cide, ckds, ckde, cjds, cjde, &
3589 cims, cime, ckms, ckme, cjms, cjme, &
3590 cits, cite, ckts, ckte, cjts, cjte, &
3591 nids, nide, nkds, nkde, njds, njde, &
3592 nims, nime, nkms, nkme, njms, njme, &
3593 nits, nite, nkts, nkte, njts, njte, &
3598 LOGICAL, INTENT(IN) :: xstag, ystag
3600 REAL, DIMENSION ( cims:cime, cjms:cjme, ckms:ckme ) :: cfld
3601 REAL, DIMENSION ( nims:nime, njms:njme, nkms:nkme ) :: nfld
3603 REAL, DIMENSION ( cims:cime, cjms:cjme, ckms:ckme ) :: ctemp_b,ctemp_bt
3604 REAL, DIMENSION ( nims:nime, njms:njme, nkms:nkme ) :: ntemp_b,ntemp_bt
3606 INTEGER, DIMENSION ( nims:nime, njms:njme ) :: imask
3607 REAL, DIMENSION( * ), INTENT(INOUT) :: c_bxs,n_bxs,c_bxe,n_bxe,c_bys,n_bys,c_bye,n_bye
3608 REAL, DIMENSION( * ), INTENT(INOUT) :: c_btxs,n_btxs,c_btxe,n_btxe,c_btys,n_btys,c_btye,n_btye
3609 REAL, DIMENSION ( cims:cime, cjms:cjme ), INTENT(IN) :: CBWGT1,CBWGT2,CBWGT3,CBWGT4 ! dummy
3610 REAL, DIMENSION ( nims:nime, njms:njme ), INTENT(IN) :: HBWGT1,HBWGT2,HBWGT3,HBWGT4
3611 INTEGER, DIMENSION ( cims:cime, cjms:cjme ), INTENT(IN) :: CII,CJJ ! dummy
3612 INTEGER, DIMENSION ( nims:nime, njms:njme ), INTENT(IN) :: IIH,JJH
3617 REAL, DIMENSION ( nims:nime, njms:njme, nkms:nkme ) :: cwk1,cwk2,cwk3,cwk4
3621 NMM_XS: IF(NITS .EQ. NIDS)THEN
3622 ! WRITE(0,*)'ENTERING X1 START BOUNDARY AT MASS POINTS',NJTS,MIN(NJTE,NJDE-1)
3625 DO J = NJTS,MIN(NJTE,NJDE-1)
3626 IF(MOD(J,2) .NE.0)THEN ! 1,3,5,7 of nested domain
3627 IF(MOD(JJH(I,J),2) .NE. 0)THEN ! 1,3,5,7 of the parent domain
3628 CWK1(I,J,K) = HBWGT1(I,J)*CFLD(IIH(I,J), JJH(I,J) ,K) &
3629 + HBWGT2(I,J)*CFLD(IIH(I,J)+1,JJH(I,J) ,K) &
3630 + HBWGT3(I,J)*CFLD(IIH(I,J), JJH(I,J)-1,K) &
3631 + HBWGT4(I,J)*CFLD(IIH(I,J), JJH(I,J)+1,K)
3635 CWK1(I,J,K) = HBWGT1(I,J)*CFLD(IIH(I,J), JJH(I,J) ,K) &
3636 + HBWGT2(I,J)*CFLD(IIH(I,J)+1,JJH(I,J) ,K) &
3637 + HBWGT3(I,J)*CFLD(IIH(I,J)+1,JJH(I,J)-1,K) &
3638 + HBWGT4(I,J)*CFLD(IIH(I,J)+1,JJH(I,J)+1,K)
3641 CWK1(I,J,K) = 0.0 ! even rows at mass points of the nested domain
3643 ntemp_b(i,J,K) = CWK1(I,J,K)
3644 ntemp_bt(i,J,K) = 0.0
3651 NMM_XE: IF(NITE-1 .EQ. NIDE-1)THEN
3652 ! WRITE(0,*)'ENTERING X END BOUNDARY AT MASS POINTS',NJTS,MIN(NJTE,NJDE-1)
3655 DO J = NJTS,MIN(NJTE,NJDE-1)
3656 IF(MOD(J,2) .NE.0)THEN ! 1,3,5,7 of the nested domain
3657 IF(MOD(JJH(I,J),2) .NE. 0)THEN ! 1,3,5,7 of the parent domain
3658 CWK2(I,J,K) = HBWGT1(I,J)*CFLD(IIH(I,J), JJH(I,J) ,K) &
3659 + HBWGT2(I,J)*CFLD(IIH(I,J)+1,JJH(I,J) ,K) &
3660 + HBWGT3(I,J)*CFLD(IIH(I,J), JJH(I,J)-1,K) &
3661 + HBWGT4(I,J)*CFLD(IIH(I,J), JJH(I,J)+1,K)
3663 CWK2(I,J,K) = HBWGT1(I,J)*CFLD(IIH(I,J), JJH(I,J) ,K) &
3664 + HBWGT2(I,J)*CFLD(IIH(I,J)+1,JJH(I,J) ,K) &
3665 + HBWGT3(I,J)*CFLD(IIH(I,J)+1,JJH(I,J)-1,K) &
3666 + HBWGT4(I,J)*CFLD(IIH(I,J)+1,JJH(I,J)+1,K)
3669 CWK2(I,J,K) = 0.0 ! even rows at mass points
3671 ntemp_b(i,J,K) = CWK2(I,J,K)
3672 ntemp_bt(i,J,K) = 0.0
3679 NMM_YS: IF(NJTS .EQ. NJDS)THEN
3680 ! WRITE(0,*)'ENTERING Y START BOUNDARY AT MASS POINTS',NITS,MIN(NITE,NIDE-1)
3683 DO I = NITS,MIN(NITE,NIDE-1)
3684 IF(MOD(JJH(I,J),2) .NE. 0)THEN ! 1,3,5,7
3685 CWK3(I,J,K) = HBWGT1(I,J)*CFLD(IIH(I,J), JJH(I,J) ,K) &
3686 + HBWGT2(I,J)*CFLD(IIH(I,J)+1,JJH(I,J) ,K) &
3687 + HBWGT3(I,J)*CFLD(IIH(I,J), JJH(I,J)-1,K) &
3688 + HBWGT4(I,J)*CFLD(IIH(I,J), JJH(I,J)+1,K)
3690 CWK3(I,J,K) = HBWGT1(I,J)*CFLD(IIH(I,J), JJH(I,J) ,K) &
3691 + HBWGT2(I,J)*CFLD(IIH(I,J)+1,JJH(I,J) ,K) &
3692 + HBWGT3(I,J)*CFLD(IIH(I,J)+1,JJH(I,J)-1,K) &
3693 + HBWGT4(I,J)*CFLD(IIH(I,J)+1,JJH(I,J)+1,K)
3695 ntemp_b(i,J,K) = CWK3(I,J,K)
3696 ntemp_bt(i,J,K) = 0.0
3703 NMM_YE: IF(NJTE-1 .EQ. NJDE-1)THEN
3704 ! WRITE(0,*)'ENTERING Y END BOUNDARY AT MASS POINTS',NITS,MIN(NITE,NIDE-1)
3707 DO I = NITS,MIN(NITE,NIDE-1)
3708 IF(MOD(JJH(I,J),2) .NE. 0)THEN ! 1,3,5,7
3709 CWK4(I,J,K) = HBWGT1(I,J)*CFLD(IIH(I,J), JJH(I,J) ,K) &
3710 + HBWGT2(I,J)*CFLD(IIH(I,J)+1,JJH(I,J) ,K) &
3711 + HBWGT3(I,J)*CFLD(IIH(I,J), JJH(I,J)-1,K) &
3712 + HBWGT4(I,J)*CFLD(IIH(I,J), JJH(I,J)+1,K)
3714 CWK4(I,J,K) = HBWGT1(I,J)*CFLD(IIH(I,J), JJH(I,J) ,K) &
3715 + HBWGT2(I,J)*CFLD(IIH(I,J)+1,JJH(I,J) ,K) &
3716 + HBWGT3(I,J)*CFLD(IIH(I,J)+1,JJH(I,J)-1,K) &
3717 + HBWGT4(I,J)*CFLD(IIH(I,J)+1,JJH(I,J)+1,K)
3720 ntemp_b(i,J,K) = CWK4(I,J,K)
3721 ntemp_bt(i,J,K) = 0.0
3728 END SUBROUTINE nmm_bdy_hinterp
3730 !--------------------------------------------------------------------------------------
3732 SUBROUTINE nmm_bdy_vinterp ( cfld, & ! CD field
3733 cids, cide, ckds, ckde, cjds, cjde, &
3734 cims, cime, ckms, ckme, cjms, cjme, &
3735 cits, cite, ckts, ckte, cjts, cjte, &
3737 nids, nide, nkds, nkde, njds, njde, &
3738 nims, nime, nkms, nkme, njms, njme, &
3739 nits, nite, nkts, nkte, njts, njte, &
3740 shw, & ! stencil half width
3741 imask, & ! interpolation mask
3742 xstag, ystag, & ! staggering of field
3743 ipos, jpos, & ! Position of lower left of nest in CD
3744 nri, nrj, & ! nest ratios
3753 CTEMP_B,NTEMP_B, & ! These temp arrays should be removed
3754 CTEMP_BT,NTEMP_BT, & ! later on
3755 CII, IIV, CJJ, JJV, CBWGT1, VBWGT1, & ! south-western grid locs and weights
3756 CBWGT2, VBWGT2, CBWGT3, VBWGT3, & ! note that "C"ourse grid ones are
3757 CBWGT4, VBWGT4 ) ! dummys
3759 ! use module_state_description
3760 USE module_configure
3761 USE module_wrf_error
3766 INTEGER, INTENT(IN) :: cids, cide, ckds, ckde, cjds, cjde, &
3767 cims, cime, ckms, ckme, cjms, cjme, &
3768 cits, cite, ckts, ckte, cjts, cjte, &
3769 nids, nide, nkds, nkde, njds, njde, &
3770 nims, nime, nkms, nkme, njms, njme, &
3771 nits, nite, nkts, nkte, njts, njte, &
3776 LOGICAL, INTENT(IN) :: xstag, ystag
3778 REAL, DIMENSION ( cims:cime, cjms:cjme, ckms:ckme ) :: cfld
3779 REAL, DIMENSION ( nims:nime, njms:njme, nkms:nkme ) :: nfld
3781 REAL, DIMENSION ( cims:cime, cjms:cjme, ckms:ckme ) :: ctemp_b,ctemp_bt
3782 REAL, DIMENSION ( nims:nime, njms:njme, nkms:nkme ) :: ntemp_b,ntemp_bt
3784 INTEGER, DIMENSION ( nims:nime, njms:njme ) :: imask
3785 REAL, DIMENSION( * ), INTENT(INOUT) :: c_bxs,n_bxs,c_bxe,n_bxe,c_bys,n_bys,c_bye,n_bye
3786 REAL, DIMENSION( * ), INTENT(INOUT) :: c_btxs,n_btxs,c_btxe,n_btxe,c_btys,n_btys,c_btye,n_btye
3787 REAL, DIMENSION ( cims:cime, cjms:cjme ), INTENT(IN) :: CBWGT1,CBWGT2,CBWGT3,CBWGT4 ! dummy
3788 REAL, DIMENSION ( nims:nime, njms:njme ), INTENT(IN) :: VBWGT1,VBWGT2,VBWGT3,VBWGT4
3789 INTEGER, DIMENSION ( cims:cime, cjms:cjme ), INTENT(IN) :: CII,CJJ ! dummy
3790 INTEGER, DIMENSION ( nims:nime, njms:njme ), INTENT(IN) :: IIV,JJV
3795 REAL, DIMENSION ( nims:nime, njms:njme, nkms:nkme ) :: cwk1,cwk2,cwk3,cwk4
3799 NMM_XS: IF(NITS .EQ. NIDS)THEN
3800 ! WRITE(0,*)'ENTERING X START BOUNDARY AT VELOCITY POINTS',NITS,NIDS,NJTS,MIN(NJTE,NJDE-1)
3803 DO J = NJTS,MIN(NJTE,NJDE-1)
3804 IF(MOD(J,2) .EQ.0)THEN ! 1,3,5,7 of nested domain
3805 IF(MOD(JJV(I,J),2) .NE. 0)THEN ! 1,3,5,7 of the parent domain
3806 CWK1(I,J,K) = VBWGT1(I,J)*CFLD(IIV(I,J), JJV(I,J) ,K) &
3807 + VBWGT2(I,J)*CFLD(IIV(I,J)+1,JJV(I,J) ,K) &
3808 + VBWGT3(I,J)*CFLD(IIV(I,J)+1,JJV(I,J)-1,K) &
3809 + VBWGT4(I,J)*CFLD(IIV(I,J)+1,JJV(I,J)+1,K)
3811 CWK1(I,J,K) = VBWGT1(I,J)*CFLD(IIV(I,J), JJV(I,J) ,K) &
3812 + VBWGT2(I,J)*CFLD(IIV(I,J)+1,JJV(I,J) ,K) &
3813 + VBWGT3(I,J)*CFLD(IIV(I,J), JJV(I,J)-1,K) &
3814 + VBWGT4(I,J)*CFLD(IIV(I,J), JJV(I,J)+1,K)
3817 CWK1(I,J,K) = 0.0 ! odd rows along J, at mass points have zero velocity
3819 ntemp_b(i,J,K) = CWK1(I,J,K)
3820 ntemp_bt(i,J,K) = 0.0
3827 NMM_XE: IF(NITE-1 .EQ. NIDE-1)THEN
3828 ! WRITE(0,*)'ENTERING X END BOUNDARY AT VELOCITY POINTS',NITE-1,NIDE-1,NJTS,MIN(NJTE,NJDE-1)
3831 DO J = NJTS,MIN(NJTE,NJDE-1)
3832 IF(MOD(J,2) .EQ.0)THEN ! 1,3,5,7 of the nested domain
3833 IF(MOD(JJV(I,J),2) .NE. 0)THEN ! 1,3,5,7 of the parent domain
3834 CWK2(I,J,K) = VBWGT1(I,J)*CFLD(IIV(I,J), JJV(I,J) ,K) &
3835 + VBWGT2(I,J)*CFLD(IIV(I,J)+1,JJV(I,J) ,K) &
3836 + VBWGT3(I,J)*CFLD(IIV(I,J)+1,JJV(I,J)-1,K) &
3837 + VBWGT4(I,J)*CFLD(IIV(I,J)+1,JJV(I,J)+1,K)
3839 CWK2(I,J,K) = VBWGT1(I,J)*CFLD(IIV(I,J), JJV(I,J) ,K) &
3840 + VBWGT2(I,J)*CFLD(IIV(I,J)+1,JJV(I,J) ,K) &
3841 + VBWGT3(I,J)*CFLD(IIV(I,J), JJV(I,J)-1,K) &
3842 + VBWGT4(I,J)*CFLD(IIV(I,J), JJV(I,J)+1,K)
3845 CWK2(I,J,K) = 0.0 ! odd rows at mass points
3847 ntemp_b(i,J,K) = CWK2(I,J,K)
3848 ntemp_bt(i,J,K) = 0.0
3855 NMM_YS: IF(NJTS .EQ. NJDS)THEN
3856 ! WRITE(0,*)'ENTERING Y START BOUNDARY AT VELOCITY POINTS',NJTS,NJDS,NITS,MIN(NITE,NIDE-1)
3859 DO I = NITS,MIN(NITE,NIDE-2) ! NIDE-1 SHOULD NOT MATTER IF WE FILL UP PHANTOM CELL
3860 IF(MOD(JJV(I,J),2) .NE. 0)THEN ! 1,3,5,7
3861 CWK3(I,J,K) = VBWGT1(I,J)*CFLD(IIV(I,J), JJV(I,J) ,K) &
3862 + VBWGT2(I,J)*CFLD(IIV(I,J)+1,JJV(I,J) ,K) &
3863 + VBWGT3(I,J)*CFLD(IIV(I,J)+1,JJV(I,J)-1,K) &
3864 + VBWGT4(I,J)*CFLD(IIV(I,J)+1,JJV(I,J)+1,K)
3866 CWK3(I,J,K) = VBWGT1(I,J)*CFLD(IIV(I,J), JJV(I,J) ,K) &
3867 + VBWGT2(I,J)*CFLD(IIV(I,J)+1,JJV(I,J) ,K) &
3868 + VBWGT3(I,J)*CFLD(IIV(I,J), JJV(I,J)-1,K) &
3869 + VBWGT4(I,J)*CFLD(IIV(I,J), JJV(I,J)+1,K)
3871 ntemp_b(i,J,K) = CWK3(I,J,K)
3872 ntemp_bt(i,J,K) = 0.0
3879 NMM_YE: IF(NJTE-1 .EQ. NJDE-1)THEN
3880 ! WRITE(0,*)'ENTERING Y END BOUNDARY AT VELOCITY POINTS',NJTE-1,NJDE-1,NITS,MIN(NITE,NIDE-1)
3883 DO I = NITS,MIN(NITE,NIDE-2) ! NIDE-1 SHOULD NOT MATTER IF WE FILL UP PHANTOM CELL
3884 IF(MOD(JJV(I,J),2) .NE. 0)THEN ! 1,3,5,7
3885 CWK4(I,J,K) = VBWGT1(I,J)*CFLD(IIV(I,J), JJV(I,J) ,K) &
3886 + VBWGT2(I,J)*CFLD(IIV(I,J)+1,JJV(I,J) ,K) &
3887 + VBWGT3(I,J)*CFLD(IIV(I,J)+1,JJV(I,J)-1,K) &
3888 + VBWGT4(I,J)*CFLD(IIV(I,J)+1,JJV(I,J)+1,K)
3890 CWK4(I,J,K) = VBWGT1(I,J)*CFLD(IIV(I,J), JJV(I,J) ,K) &
3891 + VBWGT2(I,J)*CFLD(IIV(I,J)+1,JJV(I,J) ,K) &
3892 + VBWGT3(I,J)*CFLD(IIV(I,J), JJV(I,J)-1,K) &
3893 + VBWGT4(I,J)*CFLD(IIV(I,J), JJV(I,J)+1,K)
3895 ntemp_b(i,J,K) = CWK4(I,J,K)
3896 ntemp_bt(i,J,K) = 0.0
3903 END SUBROUTINE nmm_bdy_vinterp
3906 !=======================================================================================
3907 ! E grid interpolation: simple copy from parent to mother domain
3908 !=======================================================================================
3911 SUBROUTINE nmm_copy ( cfld, & ! CD field
3912 cids, cide, ckds, ckde, cjds, cjde, &
3913 cims, cime, ckms, ckme, cjms, cjme, &
3914 cits, cite, ckts, ckte, cjts, cjte, &
3916 nids, nide, nkds, nkde, njds, njde, &
3917 nims, nime, nkms, nkme, njms, njme, &
3918 nits, nite, nkts, nkte, njts, njte, &
3919 shw, & ! stencil half width
3920 imask, & ! interpolation mask
3921 xstag, ystag, & ! staggering of field
3922 ipos, jpos, & ! Position of lower left of nest in CD
3923 nri, nrj, & ! nest ratios
3924 CII, IIH, CJJ, JJH )
3929 LOGICAL, INTENT(IN) :: xstag, ystag
3930 INTEGER, INTENT(IN) :: cids, cide, ckds, ckde, cjds, cjde, &
3931 cims, cime, ckms, ckme, cjms, cjme, &
3932 cits, cite, ckts, ckte, cjts, cjte, &
3933 nids, nide, nkds, nkde, njds, njde, &
3934 nims, nime, nkms, nkme, njms, njme, &
3935 nits, nite, nkts, nkte, njts, njte, &
3939 REAL, DIMENSION ( cims:cime, cjms:cjme, ckms:ckme ), INTENT(IN) :: cfld
3940 REAL, DIMENSION ( nims:nime, njms:njme, nkms:nkme ), INTENT(INOUT) :: nfld
3941 INTEGER, DIMENSION ( nims:nime, njms:njme ), INTENT(IN) :: imask
3942 INTEGER, DIMENSION ( cims:cime, cjms:cjme ), INTENT(IN) :: CII,CJJ ! dummy
3943 INTEGER, DIMENSION ( nims:nime, njms:njme ), INTENT(IN) :: IIH,JJH
3948 DO J=NJTS,MIN(NJTE,NJDE-1)
3950 DO I=NITS,MIN(NITE,NIDE-1)
3951 NFLD(I,J,K) = CFLD(IIH(I,J),JJH(I,J),K)
3958 END SUBROUTINE nmm_copy
3960 !=======================================================================================
3961 ! E grid test for mass point coincidence
3962 !=======================================================================================
3964 SUBROUTINE test_nmm (cfld, & ! CD field
3965 cids, cide, ckds, ckde, cjds, cjde, &
3966 cims, cime, ckms, ckme, cjms, cjme, &
3967 cits, cite, ckts, ckte, cjts, cjte, &
3969 nids, nide, nkds, nkde, njds, njde, &
3970 nims, nime, nkms, nkme, njms, njme, &
3971 nits, nite, nkts, nkte, njts, njte, &
3972 shw, & ! stencil half width for interp
3973 imask, & ! interpolation mask
3974 xstag, ystag, & ! staggering of field
3975 ipos, jpos, & ! Position of lower left of nest in CD
3976 nri, nrj, & ! nest ratios
3977 CII, IIH, CJJ, JJH, CBWGT1, HBWGT1, & ! south-western grid locs and weights
3978 CBWGT2, HBWGT2, CBWGT3, HBWGT3, & ! note that "C"ourse grid ones are
3979 CBWGT4, HBWGT4 ) ! dummys for weights
3983 INTEGER, INTENT(IN) :: cids, cide, ckds, ckde, cjds, cjde, &
3984 cims, cime, ckms, ckme, cjms, cjme, &
3985 cits, cite, ckts, ckte, cjts, cjte, &
3986 nids, nide, nkds, nkde, njds, njde, &
3987 nims, nime, nkms, nkme, njms, njme, &
3988 nits, nite, nkts, nkte, njts, njte, &
3992 LOGICAL, INTENT(IN) :: xstag, ystag
3994 REAL, DIMENSION ( cims:cime, cjms:cjme, ckms:ckme ) :: cfld
3995 REAL, DIMENSION ( nims:nime, njms:njme, nkms:nkme ) :: nfld
3996 REAL, DIMENSION ( cims:cime, cjms:cjme ), INTENT(IN) :: CBWGT1,CBWGT2,CBWGT3,CBWGT4 ! dummy
3997 REAL, DIMENSION ( nims:nime, njms:njme ), INTENT(IN) :: HBWGT1,HBWGT2,HBWGT3,HBWGT4
3998 INTEGER, DIMENSION ( cims:cime, cjms:cjme ), INTENT(IN) :: CII,CJJ ! dummy
3999 INTEGER, DIMENSION ( nims:nime, njms:njme ), INTENT(IN) :: IIH,JJH
4000 INTEGER, DIMENSION ( nims:nime, njms:njme ) :: imask
4004 REAL,PARAMETER :: error=0.0001,error1=1.0
4007 !*** CHECK DOMAIN BOUNDS BEFORE INTERPOLATION
4009 DO J=NJTS,MIN(NJTE,NJDE-1)
4010 DO I=NITS,MIN(NITE,NIDE-1)
4011 IF(IIH(i,j).LT.(CIDS-shw) .OR. IIH(i,j).GT.(CIDE+shw)) &
4012 CALL wrf_error_fatal ('hpoints:check domain bounds along x' )
4013 IF(JJH(i,j).LT.(CJDS-shw) .OR. JJH(i,j).GT.(CJDE+shw)) &
4014 CALL wrf_error_fatal ('hpoints:check domain bounds along y' )
4019 !*** INDEX CONVENTIONS
4034 ! WRITE(0,*)NITS,MIN(NITE,NIDE-1),CITS,CITE
4035 DO J=NJTS,MIN(NJTE,NJDE-1)
4037 DO I=NITS,MIN(NITE,NIDE-1)
4038 IF(ABS(1.0-HBWGT1(I,J)) .LE. ERROR)THEN
4039 DIFF=ABS(NFLD(I,J,K)-CFLD(IIH(I,J),JJH(I,J),K))
4040 IF(DIFF .GT. ERROR)THEN
4041 CALL wrf_debug(1,"dyn_nmm: NON-COINCIDENT, NESTED MASS POINT")
4042 WRITE(0,*)I,IIH(I,J),J,JJH(I,J),HBWGT1(I,J),NFLD(I,J,K),CFLD(IIH(I,J),JJH(I,J),K),DIFF
4044 IF(DIFF .GT. ERROR1)THEN
4045 WRITE(0,*)I,IIH(I,J),J,JJH(I,J),HBWGT1(I,J),NFLD(I,J,K),CFLD(IIH(I,J),JJH(I,J),K),DIFF
4046 CALL wrf_error_fatal ('dyn_nmm: NON-COINCIDENT, NESTED MASS POINT')
4053 END SUBROUTINE test_nmm
4055 !==================================
4056 ! this is the default function used in nmm feedback at mass points.
4058 SUBROUTINE nmm_feedback ( cfld, & ! CD field
4059 cids, cide, ckds, ckde, cjds, cjde, &
4060 cims, cime, ckms, ckme, cjms, cjme, &
4061 cits, cite, ckts, ckte, cjts, cjte, &
4063 nids, nide, nkds, nkde, njds, njde, &
4064 nims, nime, nkms, nkme, njms, njme, &
4065 nits, nite, nkts, nkte, njts, njte, &
4066 shw, & ! stencil half width for interp
4067 imask, & ! interpolation mask
4068 xstag, ystag, & ! staggering of field
4069 ipos, jpos, & ! Position of lower left of nest in CD
4070 nri, nrj, & ! nest ratios
4071 CII, IIH, CJJ, JJH, &
4072 CBWGT1, HBWGT1, CBWGT2, HBWGT2, &
4073 CBWGT3, HBWGT3, CBWGT4, HBWGT4 )
4074 USE module_configure
4078 INTEGER, INTENT(IN) :: cids, cide, ckds, ckde, cjds, cjde, &
4079 cims, cime, ckms, ckme, cjms, cjme, &
4080 cits, cite, ckts, ckte, cjts, cjte, &
4081 nids, nide, nkds, nkde, njds, njde, &
4082 nims, nime, nkms, nkme, njms, njme, &
4083 nits, nite, nkts, nkte, njts, njte, &
4087 INTEGER,DIMENSION(cims:cime,cjms:cjme), INTENT(IN) :: CII,CJJ ! dummy
4088 INTEGER,DIMENSION(nims:nime,njms:njme), INTENT(IN) :: IIH,JJH
4089 REAL,DIMENSION(cims:cime,cjms:cjme), INTENT(IN) :: CBWGT1,CBWGT2,CBWGT3,CBWGT4
4090 REAL,DIMENSION(nims:nime,njms:njme), INTENT(IN) :: HBWGT1,HBWGT2,HBWGT3,HBWGT4
4091 LOGICAL, INTENT(IN) :: xstag, ystag
4093 REAL, DIMENSION ( cims:cime, cjms:cjme, ckms:ckme ), INTENT(OUT) :: cfld
4094 REAL, DIMENSION ( nims:nime, njms:njme, nkms:nkme ), INTENT(IN) :: nfld
4095 INTEGER, DIMENSION ( nims:nime, njms:njme ),INTENT(IN) :: imask
4099 INTEGER ci, cj, ck, ni, nj, nk, ip, jp, ioff, joff, ioffa, joffa
4100 INTEGER :: icmin,icmax,jcmin,jcmax
4101 INTEGER :: is, ipoints,jpoints,ijpoints
4102 INTEGER , PARAMETER :: passes = 2
4105 !=====================================================================================
4108 IF(nri .ne. 3 .OR. nrj .ne. 3) &
4109 CALL wrf_error_fatal ('Feedback works for only 1:3 ratios, currently. Modify the namelist' )
4111 ! WRITE(0,*)'SIMPLE FEED BACK IS SWITCHED ON FOR MASS'
4117 DO cj = MAX(jpos+1,cjts),MIN(jpos+(njde-njds)/nrj-1,cjte) ! exclude top and bottom BCs
4118 nj = (cj-jpos)*nrj + 1
4119 if(mod(cj,2) .eq. 0)THEN
4120 is=0 ! even rows for mass points (2,4,6,8)
4122 is=1 ! odd rows for mass points (1,3,5,7)
4124 DO ci = MAX(ipos+is,cits),MIN(ipos+(nide-nids)/nri-1,cite) ! excludes LBCs
4125 ni = (ci-ipos)*nri + 2 -is
4126 IF(IS==0)THEN ! (2,4,6,8)
4127 ! AVGH = NFLD(NI,NJ+1,NK) + NFLD(NI,NJ-1,NK) + NFLD(NI+1,NJ+1,NK)+ NFLD(NI+1,NJ-1,NK) &
4128 ! + NFLD(NI+1,NJ,NK) + NFLD(NI-1,NJ,NK) + NFLD(NI,NJ+2,NK) + NFLD(NI,NJ-2,NK) &
4129 ! + NFLD(NI+1,NJ+2,NK)+ NFLD(NI-1,NJ+2,NK)+ NFLD(NI+1,NJ-2,NK)+ NFLD(NI-1,NJ-2,NK)
4131 AVGH = NFLD(NI,NJ+2,NK) &
4132 + NFLD(NI ,NJ+1,NK) + NFLD(NI+1,NJ+1,NK) &
4133 + NFLD(NI-1,NJ ,NK) + NFLD(NI,NJ ,NK) + NFLD(NI+1,NJ ,NK) &
4134 + NFLD(NI ,NJ-1,NK) + NFLD(NI+1,NJ-1,NK) &
4138 ! AVGH = NFLD(NI,NJ+1,NK) + NFLD(NI,NJ-1,NK) + NFLD(NI-1,NJ+1,NK)+ NFLD(NI-1,NJ-1,NK) &
4139 ! + NFLD(NI+1,NJ,NK) + NFLD(NI-1,NJ,NK) + NFLD(NI,NJ+2,NK) + NFLD(NI,NJ-2,NK) &
4140 ! + NFLD(NI+1,NJ+2,NK)+ NFLD(NI-1,NJ+2,NK)+ NFLD(NI+1,NJ-2,NK)+ NFLD(NI-1,NJ-2,NK)
4142 AVGH = NFLD(NI,NJ+2,NK) &
4143 + NFLD(NI-1,NJ+1,NK) + NFLD(NI,NJ+1,NK) &
4144 + NFLD(NI-1,NJ ,NK) + NFLD(NI,NJ ,NK) + NFLD(NI+1,NJ ,NK) &
4145 + NFLD(NI-1,NJ-1,NK) + NFLD(NI,NJ-1,NK) &
4149 !dusan CFLD(CI,CK,CJ) = 0.5*CFLD(CI,CK,CJ) + 0.5*(NFLD(NI,NK,NJ)+AVGH)/13.0
4150 ! CFLD(CI,CJ,CK) = (NFLD(NI,NJ,NK)+AVGH)/13.0
4151 CFLD(CI,CJ,CK) = AVGH/9.0
4156 END SUBROUTINE nmm_feedback
4158 !===========================================================================================
4160 SUBROUTINE nmm_vfeedback ( cfld, & ! CD field
4161 cids, cide, ckds, ckde, cjds, cjde, &
4162 cims, cime, ckms, ckme, cjms, cjme, &
4163 cits, cite, ckts, ckte, cjts, cjte, &
4165 nids, nide, nkds, nkde, njds, njde, &
4166 nims, nime, nkms, nkme, njms, njme, &
4167 nits, nite, nkts, nkte, njts, njte, &
4168 shw, & ! stencil half width for interp
4169 imask, & ! interpolation mask
4170 xstag, ystag, & ! staggering of field
4171 ipos, jpos, & ! Position of lower left of nest in CD
4172 nri, nrj, & ! nest ratios
4173 CII, IIV, CJJ, JJV, &
4174 CBWGT1, VBWGT1, CBWGT2, VBWGT2, &
4175 CBWGT3, VBWGT3, CBWGT4, VBWGT4 )
4176 USE module_configure
4180 INTEGER, INTENT(IN) :: cids, cide, ckds, ckde, cjds, cjde, &
4181 cims, cime, ckms, ckme, cjms, cjme, &
4182 cits, cite, ckts, ckte, cjts, cjte, &
4183 nids, nide, nkds, nkde, njds, njde, &
4184 nims, nime, nkms, nkme, njms, njme, &
4185 nits, nite, nkts, nkte, njts, njte, &
4189 INTEGER,DIMENSION(cims:cime,cjms:cjme), INTENT(IN) :: CII,CJJ ! dummy
4190 INTEGER,DIMENSION(nims:nime,njms:njme), INTENT(IN) :: IIV,JJV
4191 REAL,DIMENSION(cims:cime,cjms:cjme), INTENT(IN) :: CBWGT1,CBWGT2,CBWGT3,CBWGT4
4192 REAL,DIMENSION(nims:nime,njms:njme), INTENT(IN) :: VBWGT1,VBWGT2,VBWGT3,VBWGT4
4193 LOGICAL, INTENT(IN) :: xstag, ystag
4195 REAL, DIMENSION ( cims:cime, cjms:cjme, ckms:ckme ), INTENT(OUT) :: cfld
4196 REAL, DIMENSION ( nims:nime, njms:njme, nkms:nkme ), INTENT(IN) :: nfld
4197 INTEGER, DIMENSION ( nims:nime, njms:njme ),INTENT(IN) :: imask
4201 INTEGER ci, cj, ck, ni, nj, nk, ip, jp, ioff, joff, ioffa, joffa
4202 INTEGER :: icmin,icmax,jcmin,jcmax
4203 INTEGER :: is, ipoints,jpoints,ijpoints
4204 INTEGER , PARAMETER :: passes = 2
4207 !=====================================================================================
4210 IF(nri .ne. 3 .OR. nrj .ne. 3) &
4211 CALL wrf_error_fatal ('Feedback works for only 1:3 ratios, currently. Modify the namelist')
4213 ! WRITE(0,*)'SIMPLE FEED BACK IS SWITCHED ON FOR VELOCITY'
4219 DO cj = MAX(jpos+1,cjts),MIN(jpos+(njde-njds)/nrj-1,cjte) ! exclude top and bottom BCs
4220 nj = (cj-jpos)*nrj + 1
4221 if(mod(cj,2) .eq. 0)THEN
4222 is=1 ! even rows for velocity points (2,4,6,8)
4224 is=0 ! odd rows for velocity points (1,3,5,7)
4226 DO ci = MAX(ipos+is,cits),MIN(ipos+(nide-nids)/nri-1,cite) ! excludes LBCs
4227 ni = (ci-ipos)*nri + 2 -is
4228 IF(IS==0)THEN ! (1,3,5,7)
4229 ! AVGV = NFLD(NI,NK,NJ+1) + NFLD(NI,NK,NJ-1) + NFLD(NI+1,NK,NJ+1)+ NFLD(NI+1,NK,NJ-1) &
4230 ! + NFLD(NI+1,NK,NJ) + NFLD(NI-1,NK,NJ) + NFLD(NI,NK,NJ+2) + NFLD(NI,NK,NJ-2) &
4231 ! + NFLD(NI+1,NK,NJ+2)+ NFLD(NI-1,NK,NJ+2)+ NFLD(NI+1,NK,NJ-2)+ NFLD(NI-1,NK,NJ-2)
4233 AVGV = NFLD(NI,NJ+2,NK) &
4234 + NFLD(NI ,NJ+1,NK) + NFLD(NI+1,NJ+1,NK) &
4235 + NFLD(NI-1,NJ ,NK) + NFLD(NI,NJ ,NK) + NFLD(NI+1,NJ ,NK) &
4236 + NFLD(NI ,NJ-1,NK) + NFLD(NI+1,NJ-1,NK) &
4240 ! AVGV = NFLD(NI,NK,NJ+1) + NFLD(NI,NK,NJ-1) + NFLD(NI-1,NK,NJ+1)+ NFLD(NI-1,NK,NJ-1) &
4241 ! + NFLD(NI+1,NK,NJ) + NFLD(NI-1,NK,NJ) + NFLD(NI,NK,NJ+2) + NFLD(NI,NK,NJ-2) &
4242 ! + NFLD(NI+1,NK,NJ+2)+ NFLD(NI-1,NK,NJ+2)+ NFLD(NI+1,NK,NJ-2)+ NFLD(NI-1,NK,NJ-2)
4244 AVGV = NFLD(NI,NJ+2,NK) &
4245 + NFLD(NI-1,NJ+1,NK) + NFLD(NI,NJ+1,NK) &
4246 + NFLD(NI-1,NJ ,NK) + NFLD(NI,NJ ,NK) + NFLD(NI+1,NJ ,NK) &
4247 + NFLD(NI-1,NJ-1,NK) + NFLD(NI,NJ-1,NK) &
4251 !dusan CFLD(CI,CK,CJ) = 0.5*CFLD(CI,CK,CJ) + 0.5*(NFLD(NI,NK,NJ)+AVGV)/13.0
4252 ! CFLD(CI,CK,CJ) = (NFLD(NI,NK,NJ)+AVGV)/13.0
4253 CFLD(CI,CJ,CK) = AVGV/9.0
4258 END SUBROUTINE nmm_vfeedback
4261 SUBROUTINE nmm_smoother ( cfld , &
4262 cids, cide, ckds, ckde, cjds, cjde, &
4263 cims, cime, ckms, ckme, cjms, cjme, &
4264 cits, cite, ckts, ckte, cjts, cjte, &
4265 nids, nide, nkds, nkde, njds, njde, &
4266 nims, nime, nkms, nkme, njms, njme, &
4267 nits, nite, nkts, nkte, njts, njte, &
4273 USE module_configure
4276 INTEGER, INTENT(IN) :: cids, cide, ckds, ckde, cjds, cjde, &
4277 cims, cime, ckms, ckme, cjms, cjme, &
4278 cits, cite, ckts, ckte, cjts, cjte, &
4279 nids, nide, nkds, nkde, njds, njde, &
4280 nims, nime, nkms, nkme, njms, njme, &
4281 nits, nite, nkts, nkte, njts, njte, &
4284 REAL, DIMENSION ( cims:cime, ckms:ckme, cjms:cjme ), INTENT(INOUT) :: cfld
4285 LOGICAL, INTENT(IN) :: xstag, ystag
4291 INTEGER, PARAMETER :: smooth_passes = 5
4293 REAL, DIMENSION ( cims:cime, ckms:ckme, cjms:cjme ) :: cfldnew
4294 INTEGER :: ci, cj, ck
4295 INTEGER :: is, npass
4299 ! If there is no feedback, there can be no smoothing.
4301 CALL nl_get_feedback ( 1, feedback )
4302 IF ( feedback == 0 ) RETURN
4304 WRITE(0,*)'SIMPLE SMOOTHER IS SWITCHED ON FOR HEIGHT'
4306 DO npass = 1, smooth_passes
4308 DO cj = MAX(jpos+1,cjts),MIN(jpos+(njde-njds)/nrj-1,cjte) ! exclude top and bottom BCs
4309 if(mod(cj,2) .eq. 0)THEN
4310 is=0 ! even rows for mass points (2,4,6,8)
4312 is=1 ! odd rows for mass points (1,3,5,7)
4315 DO ci = MAX(ipos+is,cits),MIN(ipos+(nide-nids)/nri-1,cite) ! excludes LBCs
4316 IF(IS==0)THEN ! (2,4,6,8)
4317 AVGH = CFLD(CI,CK,CJ+1) + CFLD(CI,CK,CJ-1) + CFLD(CI+1,CK,CJ+1) + CFLD(CI+1,CK,CJ-1)
4319 AVGH = CFLD(CI,CK,CJ+1) + CFLD(CI,CK,CJ-1) + CFLD(CI-1,CK,CJ+1) + CFLD(CI-1,CK,CJ-1)
4321 CFLDNEW(CI,CK,CJ) = (AVGH + 4*CFLD(CI,CK,CJ)) / 8.0
4326 DO cj = MAX(jpos+1,cjts),MIN(jpos+(njde-njds)/nrj-1,cjte) ! exclude top and bottom BCs
4327 if(mod(cj,2) .eq. 0)THEN
4328 is=0 ! even rows for mass points (2,4,6,8)
4330 is=1 ! odd rows for mass points (1,3,5,7)
4333 DO ci = MAX(ipos+is,cits),MIN(ipos+(nide-nids)/nri-1,cite) ! excludes LBCs
4334 CFLD(CI,CK,CJ) = CFLDNEW(CI,CK,CJ)
4341 END SUBROUTINE nmm_smoother
4344 SUBROUTINE nmm_vsmoother ( cfld , &
4345 cids, cide, ckds, ckde, cjds, cjde, &
4346 cims, cime, ckms, ckme, cjms, cjme, &
4347 cits, cite, ckts, ckte, cjts, cjte, &
4348 nids, nide, nkds, nkde, njds, njde, &
4349 nims, nime, nkms, nkme, njms, njme, &
4350 nits, nite, nkts, nkte, njts, njte, &
4356 USE module_configure
4359 INTEGER, INTENT(IN) :: cids, cide, ckds, ckde, cjds, cjde, &
4360 cims, cime, ckms, ckme, cjms, cjme, &
4361 cits, cite, ckts, ckte, cjts, cjte, &
4362 nids, nide, nkds, nkde, njds, njde, &
4363 nims, nime, nkms, nkme, njms, njme, &
4364 nits, nite, nkts, nkte, njts, njte, &
4367 REAL, DIMENSION ( cims:cime, ckms:ckme, cjms:cjme ), INTENT(INOUT) :: cfld
4368 LOGICAL, INTENT(IN) :: xstag, ystag
4374 INTEGER, PARAMETER :: smooth_passes = 5
4376 REAL, DIMENSION ( cims:cime, ckms:ckme, cjms:cjme ) :: cfldnew
4377 INTEGER :: ci, cj, ck
4378 INTEGER :: is, npass
4382 ! If there is no feedback, there can be no smoothing.
4384 CALL nl_get_feedback ( 1, feedback )
4385 IF ( feedback == 0 ) RETURN
4387 WRITE(0,*)'SIMPLE SMOOTHER IS SWITCHED ON FOR VELOCITY'
4389 DO npass = 1, smooth_passes
4391 DO cj = MAX(jpos+1,cjts),MIN(jpos+(njde-njds)/nrj-1,cjte) ! exclude top and bottom BCs
4392 if(mod(cj,2) .eq. 0)THEN
4393 is=1 ! even rows for mass points (2,4,6,8)
4395 is=0 ! odd rows for mass points (1,3,5,7)
4398 DO ci = MAX(ipos+is,cits),MIN(ipos+(nide-nids)/nri-1,cite) ! excludes LBCs
4399 IF(IS==0)THEN ! (2,4,6,8)
4400 AVGV = CFLD(CI,CK,CJ+1) + CFLD(CI,CK,CJ-1) + CFLD(CI+1,CK,CJ+1) + CFLD(CI+1,CK,CJ-1)
4402 AVGV = CFLD(CI,CK,CJ+1) + CFLD(CI,CK,CJ-1) + CFLD(CI-1,CK,CJ+1) + CFLD(CI-1,CK,CJ-1)
4404 CFLDNEW(CI,CK,CJ) = (AVGV + 4*CFLD(CI,CK,CJ)) / 8.0
4409 DO cj = MAX(jpos+1,cjts),MIN(jpos+(njde-njds)/nrj-1,cjte) ! exclude top and bottom BCs
4410 if(mod(cj,2) .eq. 0)THEN
4411 is=1 ! even rows for mass points (2,4,6,8)
4413 is=0 ! odd rows for mass points (1,3,5,7)
4416 DO ci = MAX(ipos+is,cits),MIN(ipos+(nide-nids)/nri-1,cite) ! excludes LBCs
4417 CFLD(CI,CK,CJ) = CFLDNEW(CI,CK,CJ)
4424 END SUBROUTINE nmm_vsmoother
4425 !======================================================================================
4426 ! End of gopal's doing
4427 !======================================================================================