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 ( imask(ni,nj) .eq. 1 ) then
96 icmask( i, j ) = .TRUE.
98 if (ni .ge. nits-nioff .and. nj .ge. njts-njoff ) then
99 if ( imask(ni-nioff,nj-njoff) .eq. 1) then
100 icmask( i, j ) = .TRUE.
104 psca(i,j,nf) = cfld(i,k,j)
109 ! tile dims in this call to sint are 1-over to account for the fact
110 ! that the number of cells on the nest local subdomain is not
111 ! necessarily a multiple of the nest ratio in a given dim.
112 ! this could be a little less ham-handed.
117 cims, cime, cjms, cjme, icmask, &
118 cits-1, cite+1, cjts-1, cjte+1, nrj*nri, xstag, ystag )
120 !call end_timing( ' sint ' )
122 DO nj = njts, njte+joff
123 cj = jpos + (nj-1) / nrj ! j coord of CD point
124 jp = mod ( nj-1 , nrj ) ! coord of ND w/i CD point
127 DO ni = nits, nite+ioff
128 ci = ipos + (ni-1) / nri ! i coord of CD point
129 ip = mod ( ni-1 , nri ) ! coord of ND w/i CD point
130 if ( imask ( ni, nj ) .eq. 1 .or. imask ( ni-ioff, nj-joff ) .eq. 1 ) then
131 nfld( ni-ioff, nk, nj-joff ) = psca( ci , cj, ip+1 + (jp)*nri )
136 !$OMP END PARALLEL DO
140 !write(0,'(") cims:cime, ckms:ckme, cjms:cjme ",6i4)')cims,cime, ckms,ckme, cjms,cjme
141 !write(0,'(") nims:nime, nkms:nkme, njms:njme ",6i4)')nims,nime, nkms,nkme, njms,njme
142 !write(0,'(") cits:cite, ckts:ckte, cjts:cjte ",6i4)')cits,cite, ckts,ckte, cjts,cjte
143 !write(0,'(") nits:nite, nkts:nkte, njts:njte ",6i4)')nits,nite, nkts,nkte, njts,njte
146 cj = jpos + (nj-1) / nrj ! j coord of CD point
147 jp = mod ( nj , nrj ) ! coord of ND w/i CD point
151 ci = ipos + (ni-1) / nri ! j coord of CD point
152 ip = mod ( ni , nri ) ! coord of ND w/i CD point
153 ! This is a trivial implementation of the interp_fcn; just copies
154 ! the values from the CD into the ND
155 if ( imask ( ni, nj ) .eq. 1 ) then
156 nfld( ni, nk, nj ) = cfld( ci , ck , cj )
165 END SUBROUTINE interp_fcn
167 !==================================
168 ! this is the default function used in feedback.
170 SUBROUTINE copy_fcn ( cfld, & ! CD field
171 cids, cide, ckds, ckde, cjds, cjde, &
172 cims, cime, ckms, ckme, cjms, cjme, &
173 cits, cite, ckts, ckte, cjts, cjte, &
175 nids, nide, nkds, nkde, njds, njde, &
176 nims, nime, nkms, nkme, njms, njme, &
177 nits, nite, nkts, nkte, njts, njte, &
178 shw, & ! stencil half width for interp
179 imask, & ! interpolation mask
180 xstag, ystag, & ! staggering of field
181 ipos, jpos, & ! Position of lower left of nest in CD
182 nri, nrj ) ! nest ratios
187 INTEGER, INTENT(IN) :: cids, cide, ckds, ckde, cjds, cjde, &
188 cims, cime, ckms, ckme, cjms, cjme, &
189 cits, cite, ckts, ckte, cjts, cjte, &
190 nids, nide, nkds, nkde, njds, njde, &
191 nims, nime, nkms, nkme, njms, njme, &
192 nits, nite, nkts, nkte, njts, njte, &
196 LOGICAL, INTENT(IN) :: xstag, ystag
198 REAL, DIMENSION ( cims:cime, ckms:ckme, cjms:cjme ), INTENT(OUT) :: cfld
199 REAL, DIMENSION ( nims:nime, nkms:nkme, njms:njme ),INTENT(IN) :: nfld
200 INTEGER, DIMENSION ( nims:nime, njms:njme ),INTENT(IN) :: imask
204 INTEGER ci, cj, ck, ni, nj, nk, ip, jp, ioff, joff, ioffa, joffa
205 INTEGER :: icmin,icmax,jcmin,jcmax
206 INTEGER :: istag,jstag, ipoints,jpoints,ijpoints
207 INTEGER , PARAMETER :: passes = 2
210 ! Loop over the coarse grid in the area of the fine mesh. Do not
211 ! process the coarse grid values that are along the lateral BC
212 ! provided to the fine grid. Since that is in the specified zone
213 ! for the fine grid, it should not be used in any feedback to the
214 ! coarse grid as it should not have changed.
216 ! Due to peculiarities of staggering, it is simpler to handle the feedback
217 ! for the staggerings based upon whether it is a even ratio (2::1, 4::1, etc.) or
218 ! an odd staggering ratio (3::1, 5::1, etc.).
220 ! Though there are separate grid ratios for the i and j directions, this code
221 ! is not general enough to handle aspect ratios .NE. 1 for the fine grid cell.
223 ! These are local integer increments in the looping. Basically, istag=1 means
224 ! that we will assume one less point in the i direction. Note that ci and cj
225 ! have a maximum value that is decreased by istag and jstag, respectively.
227 ! Horizontal momentum feedback is along the face, not within the cell. For a
228 ! 3::1 ratio, temperature would use 9 pts for feedback, while u and v use
229 ! only 3 points for feedback from the nest to the parent.
231 CALL nl_get_spec_zone( 1 , spec_zone )
232 istag = 1 ; jstag = 1
233 IF ( xstag ) istag = 0
234 IF ( ystag ) jstag = 0
236 IF( MOD(nrj,2) .NE. 0) THEN ! odd refinement ratio
238 IF ( ( .NOT. xstag ) .AND. ( .NOT. ystag ) ) THEN
239 DO cj = MAX(jpos+spec_zone,cjts),MIN(jpos+(njde-njds)/nrj-jstag-spec_zone,cjte)
240 nj = (cj-jpos)*nrj + jstag + 1
243 DO ci = MAX(ipos+spec_zone,cits),MIN(ipos+(nide-nids)/nri-istag-spec_zone,cite)
244 ni = (ci-ipos)*nri + istag + 1
245 cfld( ci, ck, cj ) = 0.
246 DO ijpoints = 1 , nri * nrj
247 ipoints = MOD((ijpoints-1),nri) + 1 - nri/2 - 1
248 jpoints = (ijpoints-1)/nri + 1 - nrj/2 - 1
249 cfld( ci, ck, cj ) = cfld( ci, ck, cj ) + &
250 1./REAL(nri*nrj) * nfld( ni+ipoints , nk , nj+jpoints )
252 ! cfld( ci, ck, cj ) = 1./9. * &
253 ! ( nfld( ni-1, nk , nj-1) + &
254 ! nfld( ni , nk , nj-1) + &
255 ! nfld( ni+1, nk , nj-1) + &
256 ! nfld( ni-1, nk , nj ) + &
257 ! nfld( ni , nk , nj ) + &
258 ! nfld( ni+1, nk , nj ) + &
259 ! nfld( ni-1, nk , nj+1) + &
260 ! nfld( ni , nk , nj+1) + &
261 ! nfld( ni+1, nk , nj+1) )
266 ELSE IF ( ( xstag ) .AND. ( .NOT. ystag ) ) THEN
267 DO cj = MAX(jpos+spec_zone,cjts),MIN(jpos+(njde-njds)/nrj-jstag-spec_zone,cjte)
268 nj = (cj-jpos)*nrj + jstag + 1
271 DO ci = MAX(ipos+spec_zone,cits),MIN(ipos+(nide-nids)/nri-istag-spec_zone,cite)
272 ni = (ci-ipos)*nri + istag + 1
273 cfld( ci, ck, cj ) = 0.
274 DO ijpoints = (nri+1)/2 , (nri+1)/2 + nri*(nri-1) , nri
275 ipoints = MOD((ijpoints-1),nri) + 1 - nri/2 - 1
276 jpoints = (ijpoints-1)/nri + 1 - nrj/2 - 1
277 cfld( ci, ck, cj ) = cfld( ci, ck, cj ) + &
278 1./REAL(nri ) * nfld( ni+ipoints , nk , nj+jpoints )
280 ! cfld( ci, ck, cj ) = 1./3. * &
281 ! ( nfld( ni , nk , nj-1) + &
282 ! nfld( ni , nk , nj ) + &
283 ! nfld( ni , nk , nj+1) )
288 ELSE IF ( ( .NOT. xstag ) .AND. ( ystag ) ) THEN
289 DO cj = MAX(jpos+spec_zone,cjts),MIN(jpos+(njde-njds)/nrj-jstag-spec_zone,cjte)
290 nj = (cj-jpos)*nrj + jstag + 1
293 DO ci = MAX(ipos+spec_zone,cits),MIN(ipos+(nide-nids)/nri-istag-spec_zone,cite)
294 ni = (ci-ipos)*nri + istag + 1
295 cfld( ci, ck, cj ) = 0.
296 DO ijpoints = ( nrj*nrj +1 )/2 - nrj/2 , ( nrj*nrj +1 )/2 - nrj/2 + nrj-1
297 ipoints = MOD((ijpoints-1),nri) + 1 - nri/2 - 1
298 jpoints = (ijpoints-1)/nri + 1 - nrj/2 - 1
299 cfld( ci, ck, cj ) = cfld( ci, ck, cj ) + &
300 1./REAL( nrj) * nfld( ni+ipoints , nk , nj+jpoints )
302 ! cfld( ci, ck, cj ) = 1./3. * &
303 ! ( nfld( ni-1, nk , nj ) + &
304 ! nfld( ni , nk , nj ) + &
305 ! nfld( ni+1, nk , nj ) )
312 ! Even refinement ratio
314 ELSE IF ( MOD(nrj,2) .EQ. 0) THEN
315 IF ( ( .NOT. xstag ) .AND. ( .NOT. ystag ) ) THEN
317 ! This is a simple schematic of the feedback indexing used in the even
318 ! ratio nests. For simplicity, a 2::1 ratio is depicted. Only the
319 ! mass variable staggering is shown.
321 ! the boxes with a "T" and four small "t" represents a coarse grid (CG)
322 ! cell, that is composed of four (2::1 ratio) fine grid (FG) cells.
324 ! Shown below is the area of the CG that is in the area of the FG. The
325 ! first grid point of the depicted CG is the starting location of the nest
326 ! in the parent domain (ipos,jpos - i_parent_start and j_parent_start from
329 ! For each of the CG points, the feedback loop is over each of the FG points
330 ! within the CG cell. For a 2::1 ratio, there are four total points (this is
331 ! the ijpoints loop). The feedback value to the CG is the arithmetic mean of
332 ! all of the FG values within each CG cell.
334 ! |-------------||-------------| |-------------||-------------|
335 ! | t t || t t | | t t || t t |
336 ! jpos+ | || | | || |
337 ! (njde-njds)- | T || T | | T || T |
338 ! jstag | || | | || |
339 ! | t t || t t | | t t || t t |
340 ! |-------------||-------------| |-------------||-------------|
341 ! |-------------||-------------| |-------------||-------------|
342 ! | t t || t t | | t t || t t |
344 ! | T || T | | T || T |
346 ! | t t || t t | | t t || t t |
347 ! |-------------||-------------| |-------------||-------------|
355 ! |-------------||-------------| |-------------||-------------|
356 ! jpoints = 1 | t t || t t | | t t || t t |
358 ! | T || T | | T || T |
360 ! jpoints = 0, | t t || t t | | t t || t t |
361 ! nj=3 |-------------||-------------| |-------------||-------------|
362 ! |-------------||-------------| |-------------||-------------|
363 ! jpoints = 1 | t t || t t | | t t || t t |
365 ! jpos | T || T | ... | T || T |
367 ! jpoints = 0, | t t || t t | ... | t t || t t |
368 ! nj=1 |-------------||-------------| |-------------||-------------|
373 ! ni = 1 3 (nide-nids)/nri
374 ! ipoints= 0 1 0 1 -istag
377 ! For performance benefits, users can comment out the inner most loop (and cfld=0) and
378 ! hardcode the loop feedback. For example, it is set up to run a 2::1 ratio
379 ! if uncommented. This lacks generality, but is likely to gain timing benefits
380 ! with compilers unable to unroll inner loops that do not have parameterized sizes.
382 ! The extra +1 ---------/ and the extra -1 ----\ (both for ci and cj)
383 ! / \ keeps the feedback out of the
384 ! / \ outer row/col, since that CG data
385 ! / \ specified the nest boundary originally
388 ! / \ a sentence to not end a line
389 ! / \ with a stupid backslash
390 DO cj = MAX(jpos+spec_zone,cjts),MIN(jpos+(njde-njds)/nrj-jstag-spec_zone,cjte)
391 nj = (cj-jpos)*nrj + jstag
394 DO ci = MAX(ipos+spec_zone,cits),MIN(ipos+(nide-nids)/nri-istag-spec_zone,cite)
395 ni = (ci-ipos)*nri + istag
396 cfld( ci, ck, cj ) = 0.
397 DO ijpoints = 1 , nri * nrj
398 ipoints = MOD((ijpoints-1),nri)
399 jpoints = (ijpoints-1)/nri
400 cfld( ci, ck, cj ) = cfld( ci, ck, cj ) + &
401 1./REAL(nri*nrj) * nfld( ni+ipoints , nk , nj+jpoints )
403 ! cfld( ci, ck, cj ) = 1./4. * &
404 ! ( nfld( ni , nk , nj ) + &
405 ! nfld( ni+1, nk , nj ) + &
406 ! nfld( ni , nk , nj+1) + &
407 ! nfld( ni+1, nk , nj+1) )
414 ELSE IF ( ( xstag ) .AND. ( .NOT. ystag ) ) THEN
442 DO cj = MAX(jpos+spec_zone,cjts),MIN(jpos+(njde-njds)/nrj-jstag-spec_zone,cjte)
443 nj = (cj-jpos)*nrj + 1
446 DO ci = MAX(ipos+spec_zone,cits),MIN(ipos+(nide-nids)/nri-istag-spec_zone,cite)
447 ni = (ci-ipos)*nri + 1
448 cfld( ci, ck, cj ) = 0.
449 DO ijpoints = 1 , nri*nrj , nri
450 ipoints = MOD((ijpoints-1),nri)
451 jpoints = (ijpoints-1)/nri
452 cfld( ci, ck, cj ) = cfld( ci, ck, cj ) + &
453 1./REAL(nri ) * nfld( ni+ipoints , nk , nj+jpoints )
455 ! cfld( ci, ck, cj ) = 1./2. * &
456 ! ( nfld( ni , nk , nj ) + &
457 ! nfld( ni , nk , nj+1) )
464 ELSE IF ( ( .NOT. xstag ) .AND. ( ystag ) ) THEN
465 DO cj = MAX(jpos+spec_zone,cjts),MIN(jpos+(njde-njds)/nrj-jstag-spec_zone,cjte)
466 nj = (cj-jpos)*nrj + 1
469 DO ci = MAX(ipos+spec_zone,cits),MIN(ipos+(nide-nids)/nri-istag-spec_zone,cite)
470 ni = (ci-ipos)*nri + 1
471 cfld( ci, ck, cj ) = 0.
472 DO ijpoints = 1 , nri
473 ipoints = MOD((ijpoints-1),nri)
474 jpoints = (ijpoints-1)/nri
475 cfld( ci, ck, cj ) = cfld( ci, ck, cj ) + &
476 1./REAL(nri ) * nfld( ni+ipoints , nk , nj+jpoints )
478 ! cfld( ci, ck, cj ) = 1./2. * &
479 ! ( nfld( ni , nk , nj ) + &
480 ! nfld( ni+1, nk , nj ) )
489 END SUBROUTINE copy_fcn
491 !==================================
492 ! this is the 1pt function used in feedback.
494 SUBROUTINE copy_fcnm ( cfld, & ! CD field
495 cids, cide, ckds, ckde, cjds, cjde, &
496 cims, cime, ckms, ckme, cjms, cjme, &
497 cits, cite, ckts, ckte, cjts, cjte, &
499 nids, nide, nkds, nkde, njds, njde, &
500 nims, nime, nkms, nkme, njms, njme, &
501 nits, nite, nkts, nkte, njts, njte, &
502 shw, & ! stencil half width for interp
503 imask, & ! interpolation mask
504 xstag, ystag, & ! staggering of field
505 ipos, jpos, & ! Position of lower left of nest in CD
506 nri, nrj ) ! nest ratios
512 INTEGER, INTENT(IN) :: cids, cide, ckds, ckde, cjds, cjde, &
513 cims, cime, ckms, ckme, cjms, cjme, &
514 cits, cite, ckts, ckte, cjts, cjte, &
515 nids, nide, nkds, nkde, njds, njde, &
516 nims, nime, nkms, nkme, njms, njme, &
517 nits, nite, nkts, nkte, njts, njte, &
521 LOGICAL, INTENT(IN) :: xstag, ystag
523 REAL, DIMENSION ( cims:cime, ckms:ckme, cjms:cjme ), INTENT(OUT) :: cfld
524 REAL, DIMENSION ( nims:nime, nkms:nkme, njms:njme ), INTENT(IN) :: nfld
525 INTEGER, DIMENSION ( nims:nime, njms:njme ), INTENT(IN) :: imask
529 INTEGER ci, cj, ck, ni, nj, nk, ip, jp, ioff, joff, ioffa, joffa
530 INTEGER :: icmin,icmax,jcmin,jcmax
531 INTEGER :: istag,jstag, ipoints,jpoints,ijpoints
532 INTEGER , PARAMETER :: passes = 2
535 CALL nl_get_spec_zone( 1, spec_zone )
536 istag = 1 ; jstag = 1
537 IF ( xstag ) istag = 0
538 IF ( ystag ) jstag = 0
540 IF( MOD(nrj,2) .NE. 0) THEN ! odd refinement ratio
542 DO cj = MAX(jpos+spec_zone,cjts),MIN(jpos+(njde-njds)/nrj-jstag-spec_zone,cjte)
543 nj = (cj-jpos)*nrj + jstag + 1
546 DO ci = MAX(ipos+spec_zone,cits),MIN(ipos+(nide-nids)/nri-istag-spec_zone,cite)
547 ni = (ci-ipos)*nri + istag + 1
548 cfld( ci, ck, cj ) = nfld( ni , nk , nj )
553 ELSE ! even refinement ratio, pick nearest neighbor on SW corner
554 DO cj = MAX(jpos+spec_zone,cjts),MIN(jpos+(njde-njds)/nrj-jstag-spec_zone,cjte)
555 nj = (cj-jpos)*nrj + 1
558 DO ci = MAX(ipos+spec_zone,cits),MIN(ipos+(nide-nids)/nri-istag-spec_zone,cite)
559 ni = (ci-ipos)*nri + 1
562 cfld( ci, ck, cj ) = nfld( ni+ipoints , nk , nj+jpoints )
571 END SUBROUTINE copy_fcnm
573 !==================================
574 ! this is the 1pt function used in feedback for integers
576 SUBROUTINE copy_fcni ( cfld, & ! CD field
577 cids, cide, ckds, ckde, cjds, cjde, &
578 cims, cime, ckms, ckme, cjms, cjme, &
579 cits, cite, ckts, ckte, cjts, cjte, &
581 nids, nide, nkds, nkde, njds, njde, &
582 nims, nime, nkms, nkme, njms, njme, &
583 nits, nite, nkts, nkte, njts, njte, &
584 shw, & ! stencil half width for interp
585 imask, & ! interpolation mask
586 xstag, ystag, & ! staggering of field
587 ipos, jpos, & ! Position of lower left of nest in CD
588 nri, nrj ) ! nest ratios
594 INTEGER, INTENT(IN) :: cids, cide, ckds, ckde, cjds, cjde, &
595 cims, cime, ckms, ckme, cjms, cjme, &
596 cits, cite, ckts, ckte, cjts, cjte, &
597 nids, nide, nkds, nkde, njds, njde, &
598 nims, nime, nkms, nkme, njms, njme, &
599 nits, nite, nkts, nkte, njts, njte, &
603 LOGICAL, INTENT(IN) :: xstag, ystag
605 INTEGER, DIMENSION ( cims:cime, ckms:ckme, cjms:cjme ), INTENT(OUT) :: cfld
606 INTEGER, DIMENSION ( nims:nime, nkms:nkme, njms:njme ), INTENT(IN) :: nfld
607 INTEGER, DIMENSION ( nims:nime, njms:njme ), INTENT(IN) :: imask
611 INTEGER ci, cj, ck, ni, nj, nk, ip, jp, ioff, joff, ioffa, joffa
612 INTEGER :: icmin,icmax,jcmin,jcmax
613 INTEGER :: istag,jstag, ipoints,jpoints,ijpoints
614 INTEGER , PARAMETER :: passes = 2
617 CALL nl_get_spec_zone( 1, spec_zone )
618 istag = 1 ; jstag = 1
619 IF ( xstag ) istag = 0
620 IF ( ystag ) jstag = 0
622 IF( MOD(nrj,2) .NE. 0) THEN ! odd refinement ratio
624 DO cj = MAX(jpos+spec_zone,cjts),MIN(jpos+(njde-njds)/nrj-jstag-spec_zone,cjte)
625 nj = (cj-jpos)*nrj + jstag + 1
628 DO ci = MAX(ipos+spec_zone,cits),MIN(ipos+(nide-nids)/nri-istag-spec_zone,cite)
629 ni = (ci-ipos)*nri + istag + 1
630 cfld( ci, ck, cj ) = nfld( ni , nk , nj )
635 ELSE ! even refinement ratio
636 DO cj = MAX(jpos+spec_zone,cjts),MIN(jpos+(njde-njds)/nrj-jstag-spec_zone,cjte)
637 nj = (cj-jpos)*nrj + 1
640 DO ci = MAX(ipos+spec_zone,cits),MIN(ipos+(nide-nids)/nri-istag-spec_zone,cite)
641 ni = (ci-ipos)*nri + 1
644 cfld( ci, ck, cj ) = nfld( ni+ipoints , nk , nj+jpoints )
653 END SUBROUTINE copy_fcni
655 !==================================
657 SUBROUTINE bdy_interp ( cfld, & ! CD field
658 cids, cide, ckds, ckde, cjds, cjde, &
659 cims, cime, ckms, ckme, cjms, cjme, &
660 cits, cite, ckts, ckte, cjts, cjte, &
662 nids, nide, nkds, nkde, njds, njde, &
663 nims, nime, nkms, nkme, njms, njme, &
664 nits, nite, nkts, nkte, njts, njte, &
665 shw, & ! stencil half width
666 imask, & ! interpolation mask
667 xstag, ystag, & ! staggering of field
668 ipos, jpos, & ! Position of lower left of nest in CD
669 nri, nrj, & ! nest ratios
674 cbdy_txs, nbdy_txs, &
675 cbdy_txe, nbdy_txe, &
676 cbdy_tys, nbdy_tys, &
677 cbdy_tye, nbdy_tye, &
683 INTEGER, INTENT(IN) :: cids, cide, ckds, ckde, cjds, cjde, &
684 cims, cime, ckms, ckme, cjms, cjme, &
685 cits, cite, ckts, ckte, cjts, cjte, &
686 nids, nide, nkds, nkde, njds, njde, &
687 nims, nime, nkms, nkme, njms, njme, &
688 nits, nite, nkts, nkte, njts, njte, &
693 LOGICAL, INTENT(IN) :: xstag, ystag
695 REAL, DIMENSION ( cims:cime, ckms:ckme, cjms:cjme ) :: cfld
696 REAL, DIMENSION ( nims:nime, nkms:nkme, njms:njme ) :: nfld
697 INTEGER, DIMENSION ( nims:nime, njms:njme ) :: imask
698 REAL, DIMENSION( * ), INTENT(INOUT) :: cbdy_xs, cbdy_txs, nbdy_xs, nbdy_txs
699 REAL, DIMENSION( * ), INTENT(INOUT) :: cbdy_xe, cbdy_txe, nbdy_xe, nbdy_txe
700 REAL, DIMENSION( * ), INTENT(INOUT) :: cbdy_ys, cbdy_tys, nbdy_ys, nbdy_tys
701 REAL, DIMENSION( * ), INTENT(INOUT) :: cbdy_ye, cbdy_tye, nbdy_ye, nbdy_tye
706 INTEGER nijds, nijde, spec_bdy_width
708 nijds = min(nids, njds)
709 nijde = max(nide, njde)
710 CALL nl_get_spec_bdy_width( 1, spec_bdy_width )
712 CALL bdy_interp1( cfld, & ! CD field
713 cids, cide, ckds, ckde, cjds, cjde, &
714 cims, cime, ckms, ckme, cjms, cjme, &
715 cits, cite, ckts, ckte, cjts, cjte, &
717 nijds, nijde , spec_bdy_width , &
718 nids, nide, nkds, nkde, njds, njde, &
719 nims, nime, nkms, nkme, njms, njme, &
720 nits, nite, nkts, nkte, njts, njte, &
722 xstag, ystag, & ! staggering of field
723 ipos, jpos, & ! Position of lower left of nest in CD
729 cbdy_txs, nbdy_txs, &
730 cbdy_txe, nbdy_txe, &
731 cbdy_tys, nbdy_tys, &
732 cbdy_tye, nbdy_tye, &
738 END SUBROUTINE bdy_interp
740 SUBROUTINE bdy_interp1( cfld, & ! CD field
741 cids, cide, ckds, ckde, cjds, cjde, &
742 cims, cime, ckms, ckme, cjms, cjme, &
743 cits, cite, ckts, ckte, cjts, cjte, &
745 nijds, nijde, spec_bdy_width , &
746 nids, nide, nkds, nkde, njds, njde, &
747 nims, nime, nkms, nkme, njms, njme, &
748 nits, nite, nkts, nkte, njts, njte, &
750 imask, & ! interpolation mask
751 xstag, ystag, & ! staggering of field
752 ipos, jpos, & ! Position of lower left of nest in CD
766 use module_state_description
769 INTEGER, INTENT(IN) :: cids, cide, ckds, ckde, cjds, cjde, &
770 cims, cime, ckms, ckme, cjms, cjme, &
771 cits, cite, ckts, ckte, cjts, cjte, &
772 nids, nide, nkds, nkde, njds, njde, &
773 nims, nime, nkms, nkme, njms, njme, &
774 nits, nite, nkts, nkte, njts, njte, &
778 INTEGER, INTENT(IN) :: nijds, nijde, spec_bdy_width
779 LOGICAL, INTENT(IN) :: xstag, ystag
781 REAL, DIMENSION ( cims:cime, ckms:ckme, cjms:cjme ), INTENT(INOUT) :: cfld
782 REAL, DIMENSION ( nims:nime, nkms:nkme, njms:njme ), INTENT(INOUT) :: nfld
783 INTEGER, DIMENSION ( nims:nime, njms:njme ) :: imask
784 REAL, DIMENSION ( * ), INTENT(INOUT) :: cbdy_xs, cbdy_txs ! not used
785 REAL, DIMENSION ( * ), INTENT(INOUT) :: cbdy_xe, cbdy_txe ! not used
786 REAL, DIMENSION ( * ), INTENT(INOUT) :: cbdy_ys, cbdy_tys ! not used
787 REAL, DIMENSION ( * ), INTENT(INOUT) :: cbdy_ye, cbdy_tye ! not used
789 REAL, DIMENSION ( njms:njme, nkms:nkme, spec_bdy_width ), INTENT(INOUT) :: bdy_xs, bdy_txs
790 REAL, DIMENSION ( njms:njme, nkms:nkme, spec_bdy_width ), INTENT(INOUT) :: bdy_xe, bdy_txe
791 REAL, DIMENSION ( nims:nime, nkms:nkme, spec_bdy_width ), INTENT(INOUT) :: bdy_ys, bdy_tys
792 REAL, DIMENSION ( nims:nime, nkms:nkme, spec_bdy_width ), INTENT(INOUT) :: bdy_ye, bdy_tye
797 INTEGER ci, cj, ck, ni, nj, nk, ni1, nj1, nk1, ip, jp, ioff, joff
802 REAL psca1(cims:cime,cjms:cjme,nri*nrj)
803 REAL psca(cims:cime,cjms:cjme,nri*nrj)
804 LOGICAL icmask( cims:cime, cjms:cjme )
814 ! statement functions for converting a nest index to coarse
815 n2ci(n) = (n+ipos*nri-1)/nri
816 n2cj(n) = (n+jpos*nrj-1)/nrj
830 ! Iterate over the ND tile and compute the values
834 CALL nl_get_spec_zone( 1, spec_zone )
835 CALL nl_get_relax_zone( 1, relax_zone )
836 sz = MIN(MAX( spec_zone, relax_zone + 1 ),spec_bdy_width)
841 !$OMP PRIVATE ( i,j,k,ni,nj,ni1,nj1,ci,cj,ip,jp,nk,ck,nf,icmask,psca,psca1 )
846 nj = (j-jpos) * nrj + ( nrj / 2 + 1 ) ! j point on nest
848 ni = (i-ipos) * nri + ( nri / 2 + 1 ) ! i point on nest
849 psca1(i,j,nf) = cfld(i,k,j)
853 ! hopefully less ham handed but still correct and more efficient
854 ! sintb ignores icmask so it does not matter that icmask is not set
857 IF ( njts .ge. njds .and. njts .le. njds + sz + joff ) THEN
858 CALL sintb( psca1, psca, &
859 cims, cime, cjms, cjme, icmask, &
860 n2ci(nits)-1, n2ci(nite)+1, n2cj(MAX(njts,njds)), n2cj(MIN(njte,njds+sz+joff)), nrj*nri, xstag, ystag )
863 IF ( njte .le. njde .and. njte .ge. njde - sz - joff ) THEN
864 CALL sintb( psca1, psca, &
865 cims, cime, cjms, cjme, icmask, &
866 n2ci(nits)-1, n2ci(nite)+1, n2cj(MAX(njts,njde-sz-joff)), n2cj(MIN(njte,njde-1+joff)), nrj*nri, xstag, ystag )
869 IF ( nits .ge. nids .and. nits .le. nids + sz + ioff ) THEN
870 CALL sintb( psca1, psca, &
871 cims, cime, cjms, cjme, icmask, &
872 n2ci(MAX(nits,nids)), n2ci(MIN(nite,nids+sz+ioff)), n2cj(njts)-1, n2cj(njte)+1, nrj*nri, xstag, ystag )
875 IF ( nite .le. nide .and. nite .ge. nide - sz - ioff ) THEN
876 CALL sintb( psca1, psca, &
877 cims, cime, cjms, cjme, icmask, &
878 n2ci(MAX(nits,nide-sz-ioff)), n2ci(MIN(nite,nide-1+ioff)), n2cj(njts)-1, n2cj(njte)+1, nrj*nri, xstag, ystag )
881 DO nj1 = MAX(njds,njts-1), MIN(njde+joff,njte+joff+1)
882 cj = jpos + (nj1-1) / nrj ! j coord of CD point
883 jp = mod ( nj1-1 , nrj ) ! coord of ND w/i CD point
886 DO ni1 = MAX(nids,nits-1), MIN(nide+ioff,nite+ioff+1)
887 ci = ipos + (ni1-1) / nri ! j coord of CD point
888 ip = mod ( ni1-1 , nri ) ! coord of ND w/i CD point
893 IF ( ( ni.LT.nids) .OR. (nj.LT.njds) ) THEN
897 !bdy contains the value at t-dt. psca contains the value at t
898 !compute dv/dt and store in bdy_t
899 !afterwards store the new value of v at t into bdy
901 IF ( ni .ge. nids .and. ni .lt. nids + sz ) THEN
902 bdy_txs( nj,k,ni ) = rdt*(psca(ci+shw,cj+shw,ip+1+(jp)*nri)-nfld(ni,k,nj))
903 bdy_xs( nj,k,ni ) = psca(ci+shw,cj+shw,ip+1+(jp)*nri )
907 IF ( nj .ge. njds .and. nj .lt. njds + sz ) THEN
908 bdy_tys( ni,k,nj ) = rdt*(psca(ci+shw,cj+shw,ip+1+(jp)*nri)-nfld(ni,k,nj))
909 bdy_ys( ni,k,nj ) = psca(ci+shw,cj+shw,ip+1+(jp)*nri )
914 IF ( ni .ge. nide - sz + 1 .AND. ni .le. nide ) THEN
915 bdy_txe( nj,k,nide-ni+1 ) = rdt*(psca(ci+shw,cj+shw,ip+1+(jp)*nri)-nfld(ni,k,nj))
916 bdy_xe( nj,k,nide-ni+1 ) = psca(ci+shw,cj+shw,ip+1+(jp)*nri )
919 IF ( ni .ge. nide - sz .AND. ni .le. nide-1 ) THEN
920 bdy_txe( nj,k,nide-ni ) = rdt*(psca(ci+shw,cj+shw,ip+1+(jp)*nri)-nfld(ni,k,nj))
921 bdy_xe( nj,k,nide-ni ) = psca(ci+shw,cj+shw,ip+1+(jp)*nri )
927 IF ( nj .ge. njde - sz + 1 .AND. nj .le. njde ) THEN
928 bdy_tye( ni,k,njde-nj+1 ) = rdt*(psca(ci+shw,cj+shw,ip+1+(jp)*nri)-nfld(ni,k,nj))
929 bdy_ye( ni,k,njde-nj+1 ) = psca(ci+shw,cj+shw,ip+1+(jp)*nri )
932 IF ( nj .ge. njde - sz .AND. nj .le. njde-1 ) THEN
933 bdy_tye(ni,k,njde-nj ) = rdt*(psca(ci+shw,cj+shw,ip+1+(jp)*nri)-nfld(ni,k,nj))
934 bdy_ye( ni,k,njde-nj ) = psca(ci+shw,cj+shw,ip+1+(jp)*nri )
941 !$OMP END PARALLEL DO
946 END SUBROUTINE bdy_interp1
950 SUBROUTINE interp_fcni( cfld, & ! CD field
951 cids, cide, ckds, ckde, cjds, cjde, &
952 cims, cime, ckms, ckme, cjms, cjme, &
953 cits, cite, ckts, ckte, cjts, cjte, &
955 nids, nide, nkds, nkde, njds, njde, &
956 nims, nime, nkms, nkme, njms, njme, &
957 nits, nite, nkts, nkte, njts, njte, &
958 shw, & ! stencil half width
959 imask, & ! interpolation mask
960 xstag, ystag, & ! staggering of field
961 ipos, jpos, & ! Position of lower left of nest in CD
962 nri, nrj ) ! nest ratios
967 INTEGER, INTENT(IN) :: cids, cide, ckds, ckde, cjds, cjde, &
968 cims, cime, ckms, ckme, cjms, cjme, &
969 cits, cite, ckts, ckte, cjts, cjte, &
970 nids, nide, nkds, nkde, njds, njde, &
971 nims, nime, nkms, nkme, njms, njme, &
972 nits, nite, nkts, nkte, njts, njte, &
976 LOGICAL, INTENT(IN) :: xstag, ystag
978 INTEGER, DIMENSION ( cims:cime, ckms:ckme, cjms:cjme ) :: cfld
979 INTEGER, DIMENSION ( nims:nime, nkms:nkme, njms:njme ) :: nfld
980 INTEGER, DIMENSION ( nims:nime, njms:njme ) :: imask
984 INTEGER ci, cj, ck, ni, nj, nk, ip, jp
986 ! Iterate over the ND tile and compute the values
989 !write(0,'("cits:cite, ckts:ckte, cjts:cjte ",6i4)')cits,cite, ckts,ckte, cjts,cjte
990 !write(0,'("nits:nite, nkts:nkte, njts:njte ",6i4)')nits,nite, nkts,nkte, njts,njte
993 cj = jpos + (nj-1) / nrj ! j coord of CD point
994 jp = mod ( nj , nrj ) ! coord of ND w/i CD point
998 ci = ipos + (ni-1) / nri ! j coord of CD point
999 ip = mod ( ni , nri ) ! coord of ND w/i CD point
1000 ! This is a trivial implementation of the interp_fcn; just copies
1001 ! the values from the CD into the ND
1002 nfld( ni, nk, nj ) = cfld( ci , ck , cj )
1009 END SUBROUTINE interp_fcni
1011 SUBROUTINE interp_fcnm( cfld, & ! CD field
1012 cids, cide, ckds, ckde, cjds, cjde, &
1013 cims, cime, ckms, ckme, cjms, cjme, &
1014 cits, cite, ckts, ckte, cjts, cjte, &
1016 nids, nide, nkds, nkde, njds, njde, &
1017 nims, nime, nkms, nkme, njms, njme, &
1018 nits, nite, nkts, nkte, njts, njte, &
1019 shw, & ! stencil half width
1020 imask, & ! interpolation mask
1021 xstag, ystag, & ! staggering of field
1022 ipos, jpos, & ! Position of lower left of nest in CD
1023 nri, nrj ) ! nest ratios
1024 USE module_configure
1028 INTEGER, INTENT(IN) :: cids, cide, ckds, ckde, cjds, cjde, &
1029 cims, cime, ckms, ckme, cjms, cjme, &
1030 cits, cite, ckts, ckte, cjts, cjte, &
1031 nids, nide, nkds, nkde, njds, njde, &
1032 nims, nime, nkms, nkme, njms, njme, &
1033 nits, nite, nkts, nkte, njts, njte, &
1037 LOGICAL, INTENT(IN) :: xstag, ystag
1039 REAL, DIMENSION ( cims:cime, ckms:ckme, cjms:cjme ) :: cfld
1040 REAL, DIMENSION ( nims:nime, nkms:nkme, njms:njme ) :: nfld
1041 INTEGER, DIMENSION ( nims:nime, njms:njme ) :: imask
1045 INTEGER ci, cj, ck, ni, nj, nk, ip, jp
1047 ! Iterate over the ND tile and compute the values
1050 !write(0,'("mask cits:cite, ckts:ckte, cjts:cjte ",6i4)')cits,cite, ckts,ckte, cjts,cjte
1051 !write(0,'("mask nits:nite, nkts:nkte, njts:njte ",6i4)')nits,nite, nkts,nkte, njts,njte
1054 cj = jpos + (nj-1) / nrj ! j coord of CD point
1055 jp = mod ( nj , nrj ) ! coord of ND w/i CD point
1059 ci = ipos + (ni-1) / nri ! j coord of CD point
1060 ip = mod ( ni , nri ) ! coord of ND w/i CD point
1061 ! This is a trivial implementation of the interp_fcn; just copies
1062 ! the values from the CD into the ND
1063 nfld( ni, nk, nj ) = cfld( ci , ck , cj )
1070 END SUBROUTINE interp_fcnm
1072 SUBROUTINE interp_mask_land_field ( enable, & ! says whether to allow interpolation or just the bcasts
1074 cids, cide, ckds, ckde, cjds, cjde, &
1075 cims, cime, ckms, ckme, cjms, cjme, &
1076 cits, cite, ckts, ckte, cjts, cjte, &
1078 nids, nide, nkds, nkde, njds, njde, &
1079 nims, nime, nkms, nkme, njms, njme, &
1080 nits, nite, nkts, nkte, njts, njte, &
1081 shw, & ! stencil half width
1082 imask, & ! interpolation mask
1083 xstag, ystag, & ! staggering of field
1084 ipos, jpos, & ! Position of lower left of nest in CD
1085 nri, nrj, & ! nest ratios
1088 USE module_configure
1089 USE module_wrf_error
1094 LOGICAL, INTENT(IN) :: enable
1095 INTEGER, INTENT(IN) :: cids, cide, ckds, ckde, cjds, cjde, &
1096 cims, cime, ckms, ckme, cjms, cjme, &
1097 cits, cite, ckts, ckte, cjts, cjte, &
1098 nids, nide, nkds, nkde, njds, njde, &
1099 nims, nime, nkms, nkme, njms, njme, &
1100 nits, nite, nkts, nkte, njts, njte, &
1104 LOGICAL, INTENT(IN) :: xstag, ystag
1106 REAL, DIMENSION ( cims:cime, ckms:ckme, cjms:cjme ) :: cfld
1107 REAL, DIMENSION ( nims:nime, nkms:nkme, njms:njme ) :: nfld
1108 INTEGER, DIMENSION ( nims:nime, njms:njme ) :: imask
1110 REAL, DIMENSION ( cims:cime, cjms:cjme ) :: clu
1111 REAL, DIMENSION ( nims:nime, njms:njme ) :: nlu
1115 INTEGER ci, cj, ck, ni, nj, nk, ip, jp
1116 INTEGER :: icount , ii , jj , ist , ien , jst , jen , iswater
1117 REAL :: avg , sum , dx , dy
1118 INTEGER , PARAMETER :: max_search = 5
1119 CHARACTER*120 message
1121 ! Find out what the water value is.
1123 CALL nl_get_iswater(1,iswater)
1125 ! Right now, only mass point locations permitted.
1127 IF ( ( .NOT. xstag) .AND. ( .NOT. ystag ) ) THEN
1129 ! Loop over each i,k,j in the nested domain.
1134 IF ( MOD ( nrj , 2 ) .EQ. 0 ) THEN
1135 cj = ( nj + (nrj/2)-1 ) / nrj + jpos -1 ! first coarse position equal to or below nest point
1137 cj = ( nj + (nrj-1)/2 ) / nrj + jpos -1 ! first coarse position equal to or below nest point
1142 IF ( MOD ( nri , 2 ) .EQ. 0 ) THEN
1143 ci = ( ni + (nri/2)-1 ) / nri + ipos -1 ! first coarse position equal to or to the left of nest point
1145 ci = ( ni + (nri-1)/2 ) / nri + ipos -1 ! first coarse position equal to or to the left of nest point
1152 ! (ci,cj+1) (ci+1,cj+1)
1166 ! For odd ratios, at (nri+1)/2, we are on the coarse grid point, so dx = 0
1168 IF ( MOD ( nri , 2 ) .EQ. 0 ) THEN
1169 dx = ( REAL ( MOD ( ni+(nri-1)/2 , nri ) ) + 0.5 ) / REAL ( nri )
1171 dx = REAL ( MOD ( ni+(nri-1)/2 , nri ) ) / REAL ( nri )
1173 IF ( MOD ( nrj , 2 ) .EQ. 0 ) THEN
1174 dy = ( REAL ( MOD ( nj+(nrj-1)/2 , nrj ) ) + 0.5 ) / REAL ( nrj )
1176 dy = REAL ( MOD ( nj+(nrj-1)/2 , nrj ) ) / REAL ( nrj )
1179 ! This is a "land only" field. If this is a water point, no operations required.
1181 IF ( ( NINT(nlu(ni ,nj )) .EQ. iswater ) ) THEN
1183 ! nfld(ni,nk,nj) = 1.e20
1186 ! If this is a nested land point, and the surrounding coarse values are all land points,
1187 ! then this is a simple 4-pt interpolation.
1189 ELSE IF ( ( NINT(nlu(ni ,nj )) .NE. iswater ) .AND. &
1190 ( NINT(clu(ci ,cj )) .NE. iswater ) .AND. &
1191 ( NINT(clu(ci+1,cj )) .NE. iswater ) .AND. &
1192 ( NINT(clu(ci ,cj+1)) .NE. iswater ) .AND. &
1193 ( NINT(clu(ci+1,cj+1)) .NE. iswater ) ) THEN
1194 nfld(ni,nk,nj) = ( 1. - dx ) * ( ( 1. - dy ) * cfld(ci ,ck,cj ) + &
1195 dy * cfld(ci ,ck,cj+1) ) + &
1196 dx * ( ( 1. - dy ) * cfld(ci+1,ck,cj ) + &
1197 dy * cfld(ci+1,ck,cj+1) )
1199 ! If this is a nested land point and there are NO coarse land values surrounding,
1200 ! we temporarily punt.
1202 ELSE IF ( ( NINT(nlu(ni ,nj )) .NE. iswater ) .AND. &
1203 ( NINT(clu(ci ,cj )) .EQ. iswater ) .AND. &
1204 ( NINT(clu(ci+1,cj )) .EQ. iswater ) .AND. &
1205 ( NINT(clu(ci ,cj+1)) .EQ. iswater ) .AND. &
1206 ( NINT(clu(ci+1,cj+1)) .EQ. iswater ) ) THEN
1207 ! nfld(ni,nk,nj) = -1.e20
1210 ! If there are some water points and some land points, take an average.
1212 ELSE IF ( NINT(nlu(ni ,nj )) .NE. iswater ) THEN
1215 IF ( NINT(clu(ci ,cj )) .NE. iswater ) THEN
1217 sum = sum + cfld(ci ,ck,cj )
1219 IF ( NINT(clu(ci+1,cj )) .NE. iswater ) THEN
1221 sum = sum + cfld(ci+1,ck,cj )
1223 IF ( NINT(clu(ci ,cj+1)) .NE. iswater ) THEN
1225 sum = sum + cfld(ci ,ck,cj+1)
1227 IF ( NINT(clu(ci+1,cj+1)) .NE. iswater ) THEN
1229 sum = sum + cfld(ci+1,ck,cj+1)
1231 nfld(ni,nk,nj) = sum / REAL ( icount )
1238 ! Get an average of the whole domain for problem locations.
1245 IF ( ( nfld(ni,nk,nj) .GT. -1.e19 ) .AND. ( nfld(ni,nk,nj) .LT. 1.e19 ) ) THEN
1247 sum = sum + nfld(ni,nk,nj)
1256 CALL wrf_dm_bcast_real( sum, 1 )
1257 CALL wrf_dm_bcast_integer( icount, 1 )
1259 IF ( icount .GT. 0 ) THEN
1260 avg = sum / REAL ( icount )
1262 ! OK, if there were any of those island situations, we try to search a bit broader
1263 ! of an area in the coarse grid.
1268 IF ( nfld(ni,nk,nj) .LT. -1.e19 ) THEN
1269 IF ( MOD ( nrj , 2 ) .EQ. 0 ) THEN
1270 cj = ( nj + (nrj/2)-1 ) / nrj + jpos -1 ! first coarse position equal to or below nest point
1272 cj = ( nj + (nrj-1)/2 ) / nrj + jpos -1 ! first coarse position equal to or below nest point
1274 IF ( MOD ( nri , 2 ) .EQ. 0 ) THEN
1275 ci = ( ni + (nri/2)-1 ) / nri + ipos -1 ! first coarse position equal to or to the left of nest point
1277 ci = ( ni + (nri-1)/2 ) / nri + ipos -1 ! first coarse position equal to or to the left of nest point
1279 ist = MAX (ci-max_search,cits)
1280 ien = MIN (ci+max_search,cite,cide-1)
1281 jst = MAX (cj-max_search,cjts)
1282 jen = MIN (cj+max_search,cjte,cjde-1)
1287 IF ( NINT(clu(ii,jj)) .NE. iswater ) THEN
1289 sum = sum + cfld(ii,nk,jj)
1293 IF ( icount .GT. 0 ) THEN
1294 nfld(ni,nk,nj) = sum / REAL ( icount )
1296 ! CALL wrf_error_fatal ( "horizontal interp error - island" )
1297 write(message,*) 'horizontal interp error - island, using average ', avg
1298 CALL wrf_message ( message )
1299 nfld(ni,nk,nj) = avg
1308 CALL wrf_error_fatal ( "only unstaggered fields right now" )
1311 END SUBROUTINE interp_mask_land_field
1313 SUBROUTINE interp_mask_water_field ( enable, & ! says whether to allow interpolation or just the bcasts
1315 cids, cide, ckds, ckde, cjds, cjde, &
1316 cims, cime, ckms, ckme, cjms, cjme, &
1317 cits, cite, ckts, ckte, cjts, cjte, &
1319 nids, nide, nkds, nkde, njds, njde, &
1320 nims, nime, nkms, nkme, njms, njme, &
1321 nits, nite, nkts, nkte, njts, njte, &
1322 shw, & ! stencil half width
1323 imask, & ! interpolation mask
1324 xstag, ystag, & ! staggering of field
1325 ipos, jpos, & ! Position of lower left of nest in CD
1326 nri, nrj, & ! nest ratios
1329 USE module_configure
1330 USE module_wrf_error
1335 LOGICAL, INTENT(IN) :: enable
1336 INTEGER, INTENT(IN) :: cids, cide, ckds, ckde, cjds, cjde, &
1337 cims, cime, ckms, ckme, cjms, cjme, &
1338 cits, cite, ckts, ckte, cjts, cjte, &
1339 nids, nide, nkds, nkde, njds, njde, &
1340 nims, nime, nkms, nkme, njms, njme, &
1341 nits, nite, nkts, nkte, njts, njte, &
1345 LOGICAL, INTENT(IN) :: xstag, ystag
1347 REAL, DIMENSION ( cims:cime, ckms:ckme, cjms:cjme ) :: cfld
1348 REAL, DIMENSION ( nims:nime, nkms:nkme, njms:njme ) :: nfld
1349 INTEGER, DIMENSION ( nims:nime, njms:njme ) :: imask
1351 REAL, DIMENSION ( cims:cime, cjms:cjme ) :: clu
1352 REAL, DIMENSION ( nims:nime, njms:njme ) :: nlu
1356 INTEGER ci, cj, ck, ni, nj, nk, ip, jp
1357 INTEGER :: icount , ii , jj , ist , ien , jst , jen , iswater
1358 REAL :: avg , sum , dx , dy
1359 INTEGER , PARAMETER :: max_search = 5
1361 ! Find out what the water value is.
1363 CALL nl_get_iswater(1,iswater)
1365 ! Right now, only mass point locations permitted.
1367 IF ( ( .NOT. xstag) .AND. ( .NOT. ystag ) ) THEN
1370 ! Loop over each i,k,j in the nested domain.
1373 IF ( MOD ( nrj , 2 ) .EQ. 0 ) THEN
1374 cj = ( nj + (nrj/2)-1 ) / nrj + jpos -1 ! first coarse position equal to or below nest point
1376 cj = ( nj + (nrj-1)/2 ) / nrj + jpos -1 ! first coarse position equal to or below nest point
1381 IF ( MOD ( nri , 2 ) .EQ. 0 ) THEN
1382 ci = ( ni + (nri/2)-1 ) / nri + ipos -1 ! first coarse position equal to or to the left of nest point
1384 ci = ( ni + (nri-1)/2 ) / nri + ipos -1 ! first coarse position equal to or to the left of nest point
1391 ! (ci,cj+1) (ci+1,cj+1)
1405 ! At ni=2, we are on the coarse grid point, so dx = 0
1407 IF ( MOD ( nri , 2 ) .EQ. 0 ) THEN
1408 dx = ( REAL ( MOD ( ni+(nri-1)/2 , nri ) ) + 0.5 ) / REAL ( nri )
1410 dx = REAL ( MOD ( ni+(nri-1)/2 , nri ) ) / REAL ( nri )
1412 IF ( MOD ( nrj , 2 ) .EQ. 0 ) THEN
1413 dy = ( REAL ( MOD ( nj+(nrj-1)/2 , nrj ) ) + 0.5 ) / REAL ( nrj )
1415 dy = REAL ( MOD ( nj+(nrj-1)/2 , nrj ) ) / REAL ( nrj )
1418 ! This is a "water only" field. If this is a land point, no operations required.
1420 IF ( ( NINT(nlu(ni ,nj )) .NE. iswater ) ) THEN
1422 ! nfld(ni,nk,nj) = 1.e20
1425 ! If this is a nested water point, and the surrounding coarse values are all water points,
1426 ! then this is a simple 4-pt interpolation.
1428 ELSE IF ( ( NINT(nlu(ni ,nj )) .EQ. iswater ) .AND. &
1429 ( NINT(clu(ci ,cj )) .EQ. iswater ) .AND. &
1430 ( NINT(clu(ci+1,cj )) .EQ. iswater ) .AND. &
1431 ( NINT(clu(ci ,cj+1)) .EQ. iswater ) .AND. &
1432 ( NINT(clu(ci+1,cj+1)) .EQ. iswater ) ) THEN
1433 nfld(ni,nk,nj) = ( 1. - dx ) * ( ( 1. - dy ) * cfld(ci ,ck,cj ) + &
1434 dy * cfld(ci ,ck,cj+1) ) + &
1435 dx * ( ( 1. - dy ) * cfld(ci+1,ck,cj ) + &
1436 dy * cfld(ci+1,ck,cj+1) )
1438 ! If this is a nested water point and there are NO coarse water values surrounding,
1439 ! we temporarily punt.
1441 ELSE IF ( ( NINT(nlu(ni ,nj )) .EQ. iswater ) .AND. &
1442 ( NINT(clu(ci ,cj )) .NE. iswater ) .AND. &
1443 ( NINT(clu(ci+1,cj )) .NE. iswater ) .AND. &
1444 ( NINT(clu(ci ,cj+1)) .NE. iswater ) .AND. &
1445 ( NINT(clu(ci+1,cj+1)) .NE. iswater ) ) THEN
1446 ! nfld(ni,nk,nj) = -1.e20
1449 ! If there are some land points and some water points, take an average.
1451 ELSE IF ( NINT(nlu(ni ,nj )) .EQ. iswater ) THEN
1454 IF ( NINT(clu(ci ,cj )) .EQ. iswater ) THEN
1456 sum = sum + cfld(ci ,ck,cj )
1458 IF ( NINT(clu(ci+1,cj )) .EQ. iswater ) THEN
1460 sum = sum + cfld(ci+1,ck,cj )
1462 IF ( NINT(clu(ci ,cj+1)) .EQ. iswater ) THEN
1464 sum = sum + cfld(ci ,ck,cj+1)
1466 IF ( NINT(clu(ci+1,cj+1)) .EQ. iswater ) THEN
1468 sum = sum + cfld(ci+1,ck,cj+1)
1470 nfld(ni,nk,nj) = sum / REAL ( icount )
1476 ! Get an average of the whole domain for problem locations.
1483 IF ( ( nfld(ni,nk,nj) .GT. -1.e19 ) .AND. ( nfld(ni,nk,nj) .LT. 1.e19 ) ) THEN
1485 sum = sum + nfld(ni,nk,nj)
1494 CALL wrf_dm_bcast_real( sum, 1 )
1495 CALL wrf_dm_bcast_integer( icount, 1 )
1497 IF ( icount .NE. 0 ) THEN
1498 avg = sum / REAL ( icount )
1501 ! OK, if there were any of those lake situations, we try to search a bit broader
1502 ! of an area in the coarse grid.
1507 IF ( nfld(ni,nk,nj) .LT. -1.e19 ) THEN
1508 IF ( MOD ( nrj , 2 ) .EQ. 0 ) THEN
1509 cj = ( nj + (nrj/2)-1 ) / nrj + jpos -1 ! first coarse position equal to or below nest point
1511 cj = ( nj + (nrj-1)/2 ) / nrj + jpos -1 ! first coarse position equal to or below nest point
1513 IF ( MOD ( nri , 2 ) .EQ. 0 ) THEN
1514 ci = ( ni + (nri/2)-1 ) / nri + ipos -1 ! first coarse position equal to or to the left of nest point
1516 ci = ( ni + (nri-1)/2 ) / nri + ipos -1 ! first coarse position equal to or to the left of nest point
1518 ist = MAX (ci-max_search,cits)
1519 ien = MIN (ci+max_search,cite,cide-1)
1520 jst = MAX (cj-max_search,cjts)
1521 jen = MIN (cj+max_search,cjte,cjde-1)
1526 IF ( NINT(clu(ii,jj)) .EQ. iswater ) THEN
1528 sum = sum + cfld(ii,nk,jj)
1532 IF ( icount .GT. 0 ) THEN
1533 nfld(ni,nk,nj) = sum / REAL ( icount )
1535 ! CALL wrf_error_fatal ( "horizontal interp error - lake" )
1536 print *,'horizontal interp error - lake, using average ',avg
1537 nfld(ni,nk,nj) = avg
1546 CALL wrf_error_fatal ( "only unstaggered fields right now" )
1549 END SUBROUTINE interp_mask_water_field
1554 SUBROUTINE smoother ( cfld , &
1555 cids, cide, ckds, ckde, cjds, cjde, &
1556 cims, cime, ckms, ckme, cjms, cjme, &
1557 cits, cite, ckts, ckte, cjts, cjte, &
1558 nids, nide, nkds, nkde, njds, njde, &
1559 nims, nime, nkms, nkme, njms, njme, &
1560 nits, nite, nkts, nkte, njts, njte, &
1561 xstag, ystag, & ! staggering of field
1562 ipos, jpos, & ! Position of lower left of nest in
1566 USE module_configure
1569 INTEGER, INTENT(IN) :: cids, cide, ckds, ckde, cjds, cjde, &
1570 cims, cime, ckms, ckme, cjms, cjme, &
1571 cits, cite, ckts, ckte, cjts, cjte, &
1572 nids, nide, nkds, nkde, njds, njde, &
1573 nims, nime, nkms, nkme, njms, njme, &
1574 nits, nite, nkts, nkte, njts, njte, &
1577 LOGICAL, INTENT(IN) :: xstag, ystag
1578 INTEGER :: smooth_option, feedback , spec_zone
1580 REAL, DIMENSION ( cims:cime, ckms:ckme, cjms:cjme ) :: cfld
1582 ! If there is no feedback, there can be no smoothing.
1584 CALL nl_get_feedback ( 1, feedback )
1585 IF ( feedback == 0 ) RETURN
1586 CALL nl_get_spec_zone ( 1, spec_zone )
1588 ! These are the 2d smoothers used on the fedback data. These filters
1589 ! are run on the coarse grid data (after the nested info has been
1590 ! fedback). Only the area of the nest in the coarse grid is filtered.
1592 CALL nl_get_smooth_option ( 1, smooth_option )
1594 IF ( smooth_option == 0 ) THEN
1596 ELSE IF ( smooth_option == 1 ) THEN
1597 CALL sm121 ( cfld , &
1598 cids, cide, ckds, ckde, cjds, cjde, &
1599 cims, cime, ckms, ckme, cjms, cjme, &
1600 cits, cite, ckts, ckte, cjts, cjte, &
1601 xstag, ystag, & ! staggering of field
1602 nids, nide, nkds, nkde, njds, njde, &
1603 nims, nime, nkms, nkme, njms, njme, &
1604 nits, nite, nkts, nkte, njts, njte, &
1606 ipos, jpos & ! Position of lower left of nest in
1608 ELSE IF ( smooth_option == 2 ) THEN
1609 CALL smdsm ( cfld , &
1610 cids, cide, ckds, ckde, cjds, cjde, &
1611 cims, cime, ckms, ckme, cjms, cjme, &
1612 cits, cite, ckts, ckte, cjts, cjte, &
1613 xstag, ystag, & ! staggering of field
1614 nids, nide, nkds, nkde, njds, njde, &
1615 nims, nime, nkms, nkme, njms, njme, &
1616 nits, nite, nkts, nkte, njts, njte, &
1618 ipos, jpos & ! Position of lower left of nest in
1622 END SUBROUTINE smoother
1624 SUBROUTINE sm121 ( cfld , &
1625 cids, cide, ckds, ckde, cjds, cjde, &
1626 cims, cime, ckms, ckme, cjms, cjme, &
1627 cits, cite, ckts, ckte, cjts, cjte, &
1628 xstag, ystag, & ! staggering of field
1629 nids, nide, nkds, nkde, njds, njde, &
1630 nims, nime, nkms, nkme, njms, njme, &
1631 nits, nite, nkts, nkte, njts, njte, &
1633 ipos, jpos & ! Position of lower left of nest in
1636 USE module_configure
1639 INTEGER, INTENT(IN) :: cids, cide, ckds, ckde, cjds, cjde, &
1640 cims, cime, ckms, ckme, cjms, cjme, &
1641 cits, cite, ckts, ckte, cjts, cjte, &
1642 nids, nide, nkds, nkde, njds, njde, &
1643 nims, nime, nkms, nkme, njms, njme, &
1644 nits, nite, nkts, nkte, njts, njte, &
1647 LOGICAL, INTENT(IN) :: xstag, ystag
1649 REAL, DIMENSION ( cims:cime, ckms:ckme, cjms:cjme ) :: cfld
1650 REAL, DIMENSION ( cims:cime, cjms:cjme ) :: cfldnew
1652 INTEGER :: i , j , k , loop
1653 INTEGER :: istag,jstag
1655 INTEGER, PARAMETER :: smooth_passes = 1 ! More passes requires a larger stencil (currently 48 pt)
1657 istag = 1 ; jstag = 1
1658 IF ( xstag ) istag = 0
1659 IF ( ystag ) jstag = 0
1661 ! Simple 1-2-1 smoother.
1663 smoothing_passes : DO loop = 1 , smooth_passes
1667 ! Initialize dummy cfldnew
1669 DO i = MAX(ipos,cits-3) , MIN(ipos+(nide-nids)/nri,cite+3)
1670 DO j = MAX(jpos,cjts-3) , MIN(jpos+(njde-njds)/nrj,cjte+3)
1671 cfldnew(i,j) = cfld(i,k,j)
1675 ! 1-2-1 smoothing in the j direction first,
1677 DO i = MAX(ipos+2,cits-2) , MIN(ipos+(nide-nids)/nri-2-istag,cite+2)
1678 DO j = MAX(jpos+2,cjts-2) , MIN(jpos+(njde-njds)/nrj-2-jstag,cjte+2)
1679 cfldnew(i,j) = 0.25 * ( cfld(i,k,j+1) + 2.*cfld(i,k,j) + cfld(i,k,j-1) )
1683 ! then 1-2-1 smoothing in the i direction last
1685 DO j = MAX(jpos+2,cjts-2) , MIN(jpos+(njde-njds)/nrj-2-jstag,cjte+2)
1686 DO i = MAX(ipos+2,cits-2) , MIN(ipos+(nide-nids)/nri-2-istag,cite+2)
1687 cfld(i,k,j) = 0.25 * ( cfldnew(i+1,j) + 2.*cfldnew(i,j) + cfldnew(i-1,j) )
1693 END DO smoothing_passes
1695 END SUBROUTINE sm121
1697 SUBROUTINE smdsm ( cfld , &
1698 cids, cide, ckds, ckde, cjds, cjde, &
1699 cims, cime, ckms, ckme, cjms, cjme, &
1700 cits, cite, ckts, ckte, cjts, cjte, &
1701 xstag, ystag, & ! staggering of field
1702 nids, nide, nkds, nkde, njds, njde, &
1703 nims, nime, nkms, nkme, njms, njme, &
1704 nits, nite, nkts, nkte, njts, njte, &
1706 ipos, jpos & ! Position of lower left of nest in
1709 USE module_configure
1712 INTEGER, INTENT(IN) :: cids, cide, ckds, ckde, cjds, cjde, &
1713 cims, cime, ckms, ckme, cjms, cjme, &
1714 cits, cite, ckts, ckte, cjts, cjte, &
1715 nids, nide, nkds, nkde, njds, njde, &
1716 nims, nime, nkms, nkme, njms, njme, &
1717 nits, nite, nkts, nkte, njts, njte, &
1720 LOGICAL, INTENT(IN) :: xstag, ystag
1722 REAL, DIMENSION ( cims:cime, ckms:ckme, cjms:cjme ) :: cfld
1723 REAL, DIMENSION ( cims:cime, cjms:cjme ) :: cfldnew
1725 REAL , DIMENSION ( 2 ) :: xnu
1726 INTEGER :: i , j , k , loop , n
1727 INTEGER :: istag,jstag
1729 INTEGER, PARAMETER :: smooth_passes = 1 ! More passes requires a larger stencil (currently 48 pt)
1731 xnu = (/ 0.50 , -0.52 /)
1733 istag = 1 ; jstag = 1
1734 IF ( xstag ) istag = 0
1735 IF ( ystag ) jstag = 0
1737 ! The odd number passes of this are the "smoother", the even
1738 ! number passes are the "de-smoother" (note the different signs on xnu).
1740 smoothing_passes : DO loop = 1 , smooth_passes * 2
1742 n = 2 - MOD ( loop , 2 )
1746 DO i = MAX(ipos+2,cits-2) , MIN(ipos+(nide-nids)/nri-2-istag,cite+2)
1747 DO j = MAX(jpos+2,cjts-2) , MIN(jpos+(njde-njds)/nrj-2-jstag,cjte+2)
1748 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))
1752 DO i = MAX(ipos+2,cits-2) , MIN(ipos+(nide-nids)/nri-2-istag,cite+2)
1753 DO j = MAX(jpos+2,cjts-2) , MIN(jpos+(njde-njds)/nrj-2-jstag,cjte+2)
1754 cfld(i,k,j) = cfldnew(i,j)
1758 DO j = MAX(jpos+2,cjts-2) , MIN(jpos+(njde-njds)/nrj-2-jstag,cjte+2)
1759 DO i = MAX(ipos+2,cits-2) , MIN(ipos+(nide-nids)/nri-2-istag,cite+2)
1760 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))
1764 DO j = MAX(jpos+2,cjts-2) , MIN(jpos+(njde-njds)/nrj-2-jstag,cjte+2)
1765 DO i = MAX(ipos+2,cits-2) , MIN(ipos+(nide-nids)/nri-2-istag,cite+2)
1766 cfld(i,k,j) = cfldnew(i,j)
1772 END DO smoothing_passes
1774 END SUBROUTINE smdsm
1776 !==================================
1777 ! this is used to modify a field over the nest so we can see where the nest is
1779 SUBROUTINE mark_domain ( cfld, & ! CD field
1780 cids, cide, ckds, ckde, cjds, cjde, &
1781 cims, cime, ckms, ckme, cjms, cjme, &
1782 cits, cite, ckts, ckte, cjts, cjte, &
1784 nids, nide, nkds, nkde, njds, njde, &
1785 nims, nime, nkms, nkme, njms, njme, &
1786 nits, nite, nkts, nkte, njts, njte, &
1787 shw, & ! stencil half width for interp
1788 imask, & ! interpolation mask
1789 xstag, ystag, & ! staggering of field
1790 ipos, jpos, & ! Position of lower left of nest in CD
1791 nri, nrj ) ! nest ratios
1792 USE module_configure
1793 USE module_wrf_error
1797 INTEGER, INTENT(IN) :: cids, cide, ckds, ckde, cjds, cjde, &
1798 cims, cime, ckms, ckme, cjms, cjme, &
1799 cits, cite, ckts, ckte, cjts, cjte, &
1800 nids, nide, nkds, nkde, njds, njde, &
1801 nims, nime, nkms, nkme, njms, njme, &
1802 nits, nite, nkts, nkte, njts, njte, &
1806 LOGICAL, INTENT(IN) :: xstag, ystag
1808 REAL, DIMENSION ( cims:cime, ckms:ckme, cjms:cjme ), INTENT(OUT) :: cfld
1809 REAL, DIMENSION ( nims:nime, nkms:nkme, njms:njme ), INTENT(IN) :: nfld
1810 INTEGER, DIMENSION ( nims:nime, njms:njme ), INTENT(IN) :: imask
1814 INTEGER ci, cj, ck, ni, nj, nk, ip, jp, ioff, joff, ioffa, joffa
1815 INTEGER :: icmin,icmax,jcmin,jcmax
1816 INTEGER :: istag,jstag, ipoints,jpoints,ijpoints
1818 istag = 1 ; jstag = 1
1819 IF ( xstag ) istag = 0
1820 IF ( ystag ) jstag = 0
1822 DO cj = MAX(jpos+1,cjts),MIN(jpos+(njde-njds)/nrj-jstag-1,cjte)
1823 nj = (cj-jpos)*nrj + jstag + 1
1826 DO ci = MAX(ipos+1,cits),MIN(ipos+(nide-nids)/nri-istag-1,cite)
1827 ni = (ci-ipos)*nri + istag + 1
1828 cfld( ci, ck, cj ) = 9021000. !magic number: Beverly Hills * 100.
1833 END SUBROUTINE mark_domain
1835 #if ( NMM_CORE == 1 )
1836 !=======================================================================================
1837 ! E grid interpolation for mass with addition of terrain adjustments. First routine
1838 ! pertains to initial conditions and the next one corresponds to boundary conditions
1839 ! This is gopal's doing
1840 !=======================================================================================
1842 SUBROUTINE interp_mass_nmm (cfld, & ! CD field
1843 cids, cide, ckds, ckde, cjds, cjde, &
1844 cims, cime, ckms, ckme, cjms, cjme, &
1845 cits, cite, ckts, ckte, cjts, cjte, &
1847 nids, nide, nkds, nkde, njds, njde, &
1848 nims, nime, nkms, nkme, njms, njme, &
1849 nits, nite, nkts, nkte, njts, njte, &
1850 shw, & ! stencil half width for interp
1851 imask, & ! interpolation mask
1852 xstag, ystag, & ! staggering of field
1853 ipos, jpos, & ! Position of lower left of nest in CD
1854 nri, nrj, & ! nest ratios
1855 CII, IIH, CJJ, JJH, CBWGT1, HBWGT1, & ! south-western grid locs and weights
1856 CBWGT2, HBWGT2, CBWGT3, HBWGT3, & ! note that "C"ourse grid ones are
1857 CBWGT4, HBWGT4, & ! dummys for weights
1858 CZ3d, Z3d, & ! Z3d interpolated from CZ3d
1859 CFIS,FIS, & ! CFIS dummy on fine domain
1860 CSM,SM, & ! CSM is dummy
1866 USE MODULE_MODEL_CONSTANTS
1870 LOGICAL,INTENT(IN) :: xstag, ystag
1871 INTEGER,INTENT(IN) :: ckzmax,kzmax
1872 INTEGER,INTENT(IN) :: cids, cide, ckds, ckde, cjds, cjde, &
1873 cims, cime, ckms, ckme, cjms, cjme, &
1874 cits, cite, ckts, ckte, cjts, cjte, &
1875 nids, nide, nkds, nkde, njds, njde, &
1876 nims, nime, nkms, nkme, njms, njme, &
1877 nits, nite, nkts, nkte, njts, njte, &
1878 shw,ipos,jpos,nri,nrj
1880 INTEGER,DIMENSION(nims:nime,njms:njme), INTENT(IN) :: IMASK
1884 INTEGER,DIMENSION(cims:cime,cjms:cjme), INTENT(IN) :: CII,CJJ ! dummy
1885 REAL,DIMENSION(cims:cime,cjms:cjme), INTENT(IN) :: CBWGT1,CBWGT2,CBWGT3
1886 REAL,DIMENSION(cims:cime,cjms:cjme), INTENT(IN) :: CBWGT4,CFIS,CSM
1887 REAL,DIMENSION(cims:cime,cjms:cjme,ckms:ckme), INTENT(IN) :: CFLD
1888 REAL,DIMENSION(cims:cime,cjms:cjme,1:KZMAX), INTENT(IN) :: CZ3d
1889 REAL,DIMENSION(1:KZMAX), INTENT(IN) :: CPSTD
1890 REAL,INTENT(IN) :: CPDTOP,CPTOP
1894 INTEGER,DIMENSION(nims:nime,njms:njme), INTENT(IN) :: IIH,JJH
1895 REAL,DIMENSION(nims:nime,njms:njme), INTENT(IN) :: HBWGT1,HBWGT2,HBWGT3
1896 REAL,DIMENSION(nims:nime,njms:njme), INTENT(IN) :: HBWGT4
1897 REAL,DIMENSION(nims:nime,njms:njme), INTENT(IN) :: FIS,SM
1898 REAL,DIMENSION(nims:nime,njms:njme,nkms:nkme), INTENT(INOUT) :: NFLD
1899 REAL,DIMENSION(1:KZMAX), INTENT(IN) :: PSTD
1900 REAL,DIMENSION(nims:nime,njms:njme,1:KZMAX), INTENT(OUT) :: Z3d
1901 REAL,INTENT(IN) :: PDTOP,PTOP
1905 INTEGER,PARAMETER :: JTB=134
1906 REAL, PARAMETER :: LAPSR=6.5E-3,GI=1./G, D608=0.608
1907 REAL, PARAMETER :: COEF3=R_D*GI*LAPSR
1908 INTEGER :: I,J,K,IDUM
1909 REAL :: dlnpdz,tvout,pmo
1910 REAL,DIMENSION(nims:nime,njms:njme) :: ZS,DUM2d
1911 REAL,DIMENSION(JTB) :: PIN,ZIN,Y2,PIO,ZOUT,DUM1,DUM2
1912 !-----------------------------------------------------------------------------------------------------
1914 !*** CHECK DOMAIN BOUNDS BEFORE INTERPOLATION
1916 DO J=NJTS,MIN(NJTE,NJDE-1)
1917 DO I=NITS,MIN(NITE,NIDE-1)
1918 IF(IIH(i,j).LT.(CIDS-shw) .OR. IIH(i,j).GT.(CIDE+shw)) &
1919 CALL wrf_error_fatal ('mass points:check domain bounds along x' )
1920 IF(JJH(i,j).LT.(CJDS-shw) .OR. JJH(i,j).GT.(CJDE+shw)) &
1921 CALL wrf_error_fatal ('mass points:check domain bounds along y' )
1925 IF(KZMAX .GT. (JTB-10)) &
1926 CALL wrf_error_fatal ('mass points: increase JTB in interp_mass_nmm')
1928 ! WRITE(21,*)'------------- MED NEST INITIAL 1 ----------------'
1929 ! DO J=NJTS,MIN(NJTE,NJDE-1)
1930 ! DO I=NITS,MIN(NITE,NIDE-1)
1931 ! WRITE(21,*)I,J,IMASK(I,J),NFLD(I,1,J)
1937 !*** DEFINE LOCAL TOPOGRAPHY ON THE BASIS OF FIS. ALSO CHECK IF SM IS LAND (SM=0) OVER TOPO
1938 !*** YOU DON'T WANT MOUNTAINS INSIDE WATER BODIES!
1941 DO J=NJTS,MIN(NJTE,NJDE-1)
1942 DO I=NITS,MIN(NITE,NIDE-1)
1948 !*** Interpolate GPMs DERIVED FROM STANDARD ATMOSPHERIC LAPSE RATE FROM THE PARENT TO
1949 !*** THE NESTED DOMAIN
1951 !*** INDEX CONVENTIONS
1966 DO K=NKTS,KZMAX ! Please note that we are still in isobaric surfaces
1967 DO J=NJTS,MIN(NJTE,NJDE-1)
1968 DO I=NITS,MIN(NITE,NIDE-1)
1970 IF(MOD(JJH(I,J),2) .NE. 0)THEN ! 1,3,5,7
1971 Z3d(I,J,K) = HBWGT1(I,J)*CZ3d(IIH(I,J), JJH(I,J) ,K) &
1972 + HBWGT2(I,J)*CZ3d(IIH(I,J)+1,JJH(I,J) ,K) &
1973 + HBWGT3(I,J)*CZ3d(IIH(I,J), JJH(I,J)-1,K) &
1974 + HBWGT4(I,J)*CZ3d(IIH(I,J), JJH(I,J)+1,K)
1976 Z3d(I,J,K) = HBWGT1(I,J)*CZ3d(IIH(I,J), JJH(I,J) ,K) &
1977 + HBWGT2(I,J)*CZ3d(IIH(I,J)+1,JJH(I,J) ,K) &
1978 + HBWGT3(I,J)*CZ3d(IIH(I,J)+1,JJH(I,J)-1,K) &
1979 + HBWGT4(I,J)*CZ3d(IIH(I,J)+1,JJH(I,J)+1,K)
1987 ! RECONSTRUCT PDs ON THE BASIS OF TOPOGRAPHY AND THE INTERPOLATED HEIGHTS
1989 DO J=NJTS,MIN(NJTE,NJDE-1)
1990 DO I=NITS,MIN(NITE,NIDE-1)
1992 IF (ZS(I,J) .LT. Z3d(I,J,1)) THEN
1993 dlnpdz = (log(PSTD(1))-log(PSTD(2)) )/(Z3d(i,j,1)-Z3d(i,j,2))
1994 dum2d(i,j) = exp(log(PSTD(1)) + dlnpdz*(ZS(I,J) - Z3d(i,j,1)))
1995 dum2d(i,j) = dum2d(i,j) - PDTOP -PTOP
1996 ELSE ! target level bounded by input levels
1997 DO K =NKTS,KZMAX-1 ! still in the isobaric surfaces
1998 IF(ZS(I,J) .GE. Z3d(I,J,K) .AND. ZS(I,J) .LT. Z3d(I,J,K+1))THEN
1999 dlnpdz = (log(PSTD(K))-log(PSTD(K+1)) ) /(Z3d(I,J,K)-Z3d(I,J,K+1))
2000 dum2d(i,j) = exp(log(PSTD(K)) + dlnpdz*(ZS(I,J)- Z3d(I,J,K)))
2001 dum2d(i,j) = dum2d(i,j) - PDTOP -PTOP
2005 IF(ZS(I,J) .GE. Z3d(I,J,KZMAX))THEN
2006 WRITE(0,*)'I=',I,'J=',J,'K=',KZMAX,'TERRAIN HEIGHT',ZS(I,J),'Z3d',Z3d(I,J,KZMAX)
2007 CALL wrf_error_fatal3 ( "interp_fcn.b" , 176 , "MOUNTAIN TOO HIGH TO FIT THE MODEL DEPTH")
2013 DO K=NKDS,NKDE ! NKTE is 1, nevertheless let us pretend religious
2014 DO J=NJTS,MIN(NJTE,NJDE-1)
2015 DO I=NITS,MIN(NITE,NIDE-1)
2016 IF(IMASK(I,J) .NE. 1)THEN
2017 NFLD(I,J,K)= dum2d(i,j) ! PD defined in the nested domain
2024 END SUBROUTINE interp_mass_nmm
2026 !--------------------------------------------------------------------------------------
2028 SUBROUTINE nmm_bdymass_hinterp ( cfld, & ! CD field
2029 cids, cide, ckds, ckde, cjds, cjde, &
2030 cims, cime, ckms, ckme, cjms, cjme, &
2031 cits, cite, ckts, ckte, cjts, cjte, &
2033 nids, nide, nkds, nkde, njds, njde, &
2034 nims, nime, nkms, nkme, njms, njme, &
2035 nits, nite, nkts, nkte, njts, njte, &
2036 shw, & ! stencil half width
2037 imask, & ! interpolation mask
2038 xstag, ystag, & ! staggering of field
2039 ipos, jpos, & ! Position of lower left of nest in CD
2040 nri, nrj, & ! nest ratios
2049 CTEMP_B,NTEMP_B, & ! These temp arrays should be removed
2050 CTEMP_BT,NTEMP_BT, & ! later on
2051 CII, IIH, CJJ, JJH, CBWGT1, HBWGT1, & ! south-western grid locs and weights
2052 CBWGT2, HBWGT2, CBWGT3, HBWGT3, & ! note that "C"ourse grid ones are
2053 CBWGT4, HBWGT4, & ! dummys
2054 CZ3d, Z3d, & ! Z3d dummy on nested domain
2055 CFIS,FIS, & ! CFIS dummy on fine domain
2056 CSM,SM, & ! CSM is dummy
2063 USE MODULE_MODEL_CONSTANTS
2064 USE module_configure
2065 USE module_wrf_error
2070 INTEGER, INTENT(IN) :: ckzmax,kzmax
2071 INTEGER, INTENT(IN) :: cids, cide, ckds, ckde, cjds, cjde, &
2072 cims, cime, ckms, ckme, cjms, cjme, &
2073 cits, cite, ckts, ckte, cjts, cjte, &
2074 nids, nide, nkds, nkde, njds, njde, &
2075 nims, nime, nkms, nkme, njms, njme, &
2076 nits, nite, nkts, nkte, njts, njte, &
2082 REAL, DIMENSION ( cims:cime, cjms:cjme, ckms:ckme ), INTENT(OUT) :: ctemp_b,ctemp_bt
2083 REAL, DIMENSION ( nims:nime, njms:njme, nkms:nkme ), INTENT(OUT) :: ntemp_b,ntemp_bt
2084 LOGICAL, INTENT(IN) :: xstag, ystag
2085 REAL, DIMENSION( * ), INTENT(INOUT) :: c_bxs,n_bxs,c_bxe,n_bxe,c_bys,n_bys,c_bye,n_bye
2086 REAL, DIMENSION( * ), INTENT(INOUT) :: c_btxs,n_btxs,c_btxe,n_btxe,c_btys,n_btys,c_btye,n_btye
2090 INTEGER,DIMENSION(nims:nime,njms:njme), INTENT(IN) :: IMASK
2091 INTEGER,DIMENSION(cims:cime,cjms:cjme), INTENT(IN) :: CII,CJJ ! dummy
2092 REAL,DIMENSION(cims:cime,cjms:cjme), INTENT(IN) :: CBWGT1,CBWGT2,CBWGT3
2093 REAL,DIMENSION(cims:cime,cjms:cjme), INTENT(IN) :: CBWGT4,CFIS,CSM
2094 REAL,DIMENSION(cims:cime,cjms:cjme,ckms:ckme), INTENT(IN) :: CFLD
2095 REAL,DIMENSION(cims:cime,cjms:cjme,1:KZMAX), INTENT(IN) :: CZ3d
2096 REAL,DIMENSION(1:KZMAX), INTENT(IN) :: CPSTD
2097 REAL,INTENT(IN) :: CPDTOP,CPTOP
2101 INTEGER,DIMENSION(nims:nime,njms:njme), INTENT(IN) :: IIH,JJH
2102 REAL,DIMENSION(nims:nime,njms:njme), INTENT(IN) :: HBWGT1,HBWGT2,HBWGT3
2103 REAL,DIMENSION(nims:nime,njms:njme), INTENT(IN) :: HBWGT4
2104 REAL,DIMENSION(nims:nime,njms:njme), INTENT(IN) :: FIS,SM
2105 REAL,DIMENSION(nims:nime,njms:njme,nkms:nkme), INTENT(INOUT) :: NFLD
2106 REAL,DIMENSION(1:KZMAX), INTENT(IN) :: PSTD
2107 REAL,DIMENSION(nims:nime,njms:njme,1:KZMAX), INTENT(OUT) :: Z3d
2108 REAL,INTENT(IN) :: PDTOP,PTOP
2112 INTEGER :: nijds, nijde, spec_bdy_width,i,j,k
2113 REAL :: dlnpdz,dum2d
2114 REAL,DIMENSION(nims:nime,njms:njme) :: zs
2116 INTEGER,PARAMETER :: JTB=134
2118 REAL, DIMENSION (nims:nime,njms:njme) :: CWK1,CWK2,CWK3,CWK4
2120 nijds = min(nids, njds)
2121 nijde = max(nide, njde)
2122 CALL nl_get_spec_bdy_width( 1, spec_bdy_width )
2125 !*** DEFINE LOCAL TOPOGRAPHY ON THE BASIS OF FIS. ASLO CHECK IF SM IS LAND (SM=0) OVER TOPO
2126 !*** YOU DON'T WANT MOUNTAINS INSIDE WATER BODIES!
2128 DO J=NJTS,MIN(NJTE,NJDE-1)
2129 DO I=NITS,MIN(NITE,NIDE-1)
2136 NMM_XS: IF(NITS .EQ. NIDS)THEN
2137 ! WRITE(0,*)'ENTERING X1 START BOUNDARY AT MASS POINTS',NJTS,MIN(NJTE,NJDE-1)
2141 DO J = NJTS,MIN(NJTE,NJDE-1)
2142 IF(MOD(JJH(I,J),2) .NE. 0)THEN ! 1,3,5,7 of the parent domain
2143 Z3d(I,J,K) = HBWGT1(I,J)*CZ3d(IIH(I,J), JJH(I,J) ,K) &
2144 + HBWGT2(I,J)*CZ3d(IIH(I,J)+1,JJH(I,J) ,K) &
2145 + HBWGT3(I,J)*CZ3d(IIH(I,J), JJH(I,J)-1,K) &
2146 + HBWGT4(I,J)*CZ3d(IIH(I,J), JJH(I,J)+1,K)
2148 Z3d(I,J,K) = HBWGT1(I,J)*CZ3d(IIH(I,J), JJH(I,J) ,K) &
2149 + HBWGT2(I,J)*CZ3d(IIH(I,J)+1,JJH(I,J) ,K) &
2150 + HBWGT3(I,J)*CZ3d(IIH(I,J)+1,JJH(I,J)-1,K) &
2151 + HBWGT4(I,J)*CZ3d(IIH(I,J)+1,JJH(I,J)+1,K)
2156 DO J = NJTS,MIN(NJTE,NJDE-1)
2157 IF(MOD(J,2) .NE. 0)THEN
2158 IF (ZS(I,J) .LT. Z3d(I,J,2)) THEN ! level 2 has to be changed
2159 dlnpdz = (log(PSTD(1))-log(PSTD(2)) )/(Z3d(I,J,1)-Z3d(I,J,2))
2160 dum2d = exp(log(PSTD(1)) + dlnpdz*(ZS(I,J) - Z3d(I,J,1)))
2161 CWK1(I,J) = dum2d -PDTOP -PTOP
2162 ELSE ! target level bounded by input levels
2164 IF(ZS(I,J) .GE. Z3d(I,J,K) .AND. ZS(I,J) .LT. Z3d(I,J,K+1))THEN
2165 dlnpdz = (log(PSTD(K))-log(PSTD(K+1)) ) /(Z3d(I,J,K)-Z3d(I,J,K+1))
2166 dum2d = exp(log(PSTD(K)) + dlnpdz*(ZS(I,J)- Z3d(I,J,K)))
2167 CWK1(I,J) = dum2d -PDTOP -PTOP
2171 IF(ZS(I,J) .GE. Z3d(I,J,KZMAX))THEN
2172 WRITE(0,*)'I=',I,'J=',J,'K=',K,'TERRAIN HEIGHT',ZS(I,J),'Z3d',Z3d(I,J,KZMAX)
2173 CALL wrf_error_fatal("BC:MOUNTAIN TOO HIGH TO FIT THE MODEL DEPTH")
2180 DO J = NJTS,MIN(NJTE,NJDE-1)
2182 ntemp_b(i,j,k) = CWK1(I,J)
2183 ntemp_bt(i,j,k) = 0.0
2190 NMM_XE: IF(NITE-1 .EQ. NIDE-1)THEN
2191 ! WRITE(0,*)'ENTERING X END BOUNDARY AT MASS POINTS',NJTS,MIN(NJTE,NJDE-1)
2196 DO J=NJTS,MIN(NJTE,NJDE-1)
2197 IF(MOD(JJH(I,J),2) .NE. 0)THEN ! 1,3,5,7
2198 Z3d(I,J,K) = HBWGT1(I,J)*CZ3d(IIH(I,J), JJH(I,J) ,K) &
2199 + HBWGT2(I,J)*CZ3d(IIH(I,J)+1,JJH(I,J) ,K) &
2200 + HBWGT3(I,J)*CZ3d(IIH(I,J), JJH(I,J)-1,K) &
2201 + HBWGT4(I,J)*CZ3d(IIH(I,J), JJH(I,J)+1,K)
2203 Z3d(I,J,K) = HBWGT1(I,J)*CZ3d(IIH(I,J), JJH(I,J) ,K) &
2204 + HBWGT2(I,J)*CZ3d(IIH(I,J)+1,JJH(I,J) ,K) &
2205 + HBWGT3(I,J)*CZ3d(IIH(I,J)+1,JJH(I,J)-1,K) &
2206 + HBWGT4(I,J)*CZ3d(IIH(I,J)+1,JJH(I,J)+1,K)
2211 DO J = NJTS,MIN(NJTE,NJDE-1)
2212 IF(MOD(J,2) .NE.0)THEN ! 1,3,5,7 of nested domain
2213 IF (ZS(I,J) .LT. Z3d(I,J,2)) THEN ! level 2 has to be changed
2214 dlnpdz = (log(PSTD(1))-log(PSTD(2)) )/(Z3d(I,J,1)-Z3d(I,J,2))
2215 dum2d = exp(log(PSTD(1)) + dlnpdz*(ZS(I,J) - Z3d(I,J,1)))
2216 CWK2(I,J) = dum2d -PDTOP -PTOP
2217 ELSE ! target level bounded by input levels
2219 IF(ZS(I,J) .GE. Z3d(I,J,K) .AND. ZS(I,J) .LT. Z3d(I,J,K+1))THEN
2220 dlnpdz = (log(PSTD(K))-log(PSTD(K+1)) ) /(Z3d(I,J,K)-Z3d(I,J,K+1))
2221 dum2d = exp(log(PSTD(K)) + dlnpdz*(ZS(I,J)- Z3d(I,J,K)))
2222 CWK2(I,J) = dum2d -PDTOP -PTOP
2226 IF(ZS(I,J) .GE. Z3d(I,J,KZMAX))THEN
2227 WRITE(0,*)'I=',I,'J=',J,'K=',K,'TERRAIN HEIGHT',ZS(I,J),'Z3d',Z3d(I,J,KZMAX)
2228 CALL wrf_error_fatal("BC:MOUNTAIN TOO HIGH TO FIT THE MODEL DEPTH")
2235 DO J = NJTS,MIN(NJTE,NJDE-1)
2237 ntemp_b(i,j,k) = CWK2(I,J)
2238 ntemp_bt(i,j,k) = 0.0
2245 NMM_YS: IF(NJTS .EQ. NJDS)THEN
2246 ! WRITE(20,*)'ENTERING Y START BOUNDARY AT MASS POINTS',NITS,MIN(NITE,NIDE-1)
2249 DO I = NITS,MIN(NITE,NIDE-1)
2250 IF(MOD(JJH(I,J),2) .NE. 0)THEN ! 1,3,5,7
2251 Z3d(I,J,K) = HBWGT1(I,J)*CZ3d(IIH(I,J), JJH(I,J) ,K) &
2252 + HBWGT2(I,J)*CZ3d(IIH(I,J)+1,JJH(I,J) ,K) &
2253 + HBWGT3(I,J)*CZ3d(IIH(I,J), JJH(I,J)-1,K) &
2254 + HBWGT4(I,J)*CZ3d(IIH(I,J), JJH(I,J)+1,K)
2256 Z3d(I,J,K) = HBWGT1(I,J)*CZ3d(IIH(I,J), JJH(I,J) ,K) &
2257 + HBWGT2(I,J)*CZ3d(IIH(I,J)+1,JJH(I,J) ,K) &
2258 + HBWGT3(I,J)*CZ3d(IIH(I,J)+1,JJH(I,J)-1,K) &
2259 + HBWGT4(I,J)*CZ3d(IIH(I,J)+1,JJH(I,J)+1,K)
2264 DO I = NITS,MIN(NITE,NIDE-1)
2265 IF (ZS(I,J) .LT. Z3d(I,J,2)) THEN ! level 2 has to be changed
2266 dlnpdz = (log(PSTD(1))-log(PSTD(2)) )/(Z3d(I,J,1)-Z3d(I,J,2))
2267 dum2d = exp(log(PSTD(1)) + dlnpdz*(ZS(I,J) - Z3d(I,J,1)))
2268 CWK3(I,J) = dum2d -PDTOP -PTOP
2269 ELSE ! target level bounded by input levels
2271 IF(ZS(I,J) .GE. Z3d(I,J,K) .AND. ZS(I,J) .LT. Z3d(I,J,K+1))THEN
2272 dlnpdz = (log(PSTD(K))-log(PSTD(K+1)) ) /(Z3d(I,J,K)-Z3d(I,J,K+1))
2273 dum2d = exp(log(PSTD(K)) + dlnpdz*(ZS(I,J)- Z3d(I,J,K)))
2274 CWK3(I,J) = dum2d -PDTOP -PTOP
2278 IF(ZS(I,J) .GE. Z3d(I,J,KZMAX))THEN
2279 WRITE(0,*)'I=',I,'J=',J,'K=',K,'TERRAIN HEIGHT',ZS(I,J),'Z3d',Z3d(I,J,KZMAX)
2280 CALL wrf_error_fatal("BC:MOUNTAIN TOO HIGH TO FIT THE MODEL DEPTH")
2285 DO I = NITS,MIN(NITE,NIDE-1)
2286 ntemp_b(i,j,k) = CWK3(I,J)
2287 ntemp_bt(i,j,k) = 0.0
2294 NMM_YE: IF(NJTE-1 .EQ. NJDE-1)THEN
2295 ! WRITE(20,*)'ENTERING Y END BOUNDARY AT MASS POINTS',NITS,MIN(NITE,NIDE-1)
2299 DO I = NITS,MIN(NITE,NIDE-1)
2300 IF(MOD(JJH(I,J),2) .NE. 0)THEN ! 1,3,5,7
2301 Z3d(I,J,K) = HBWGT1(I,J)*CZ3d(IIH(I,J), JJH(I,J) ,K) &
2302 + HBWGT2(I,J)*CZ3d(IIH(I,J)+1,JJH(I,J) ,K) &
2303 + HBWGT3(I,J)*CZ3d(IIH(I,J), JJH(I,J)-1,K) &
2304 + HBWGT4(I,J)*CZ3d(IIH(I,J), JJH(I,J)+1,K)
2306 Z3d(I,J,K) = HBWGT1(I,J)*CZ3d(IIH(I,J), JJH(I,J) ,K) &
2307 + HBWGT2(I,J)*CZ3d(IIH(I,J)+1,JJH(I,J) ,K) &
2308 + HBWGT3(I,J)*CZ3d(IIH(I,J)+1,JJH(I,J)-1,K) &
2309 + HBWGT4(I,J)*CZ3d(IIH(I,J)+1,JJH(I,J)+1,K)
2314 DO I = NITS,MIN(NITE,NIDE-1)
2315 IF (ZS(I,J) .LT. Z3d(I,J,2)) THEN ! level 2 has to be changed
2316 dlnpdz = (log(PSTD(1))-log(PSTD(2)) )/(Z3d(I,J,1)-Z3d(I,J,2))
2317 dum2d = exp(log(PSTD(1)) + dlnpdz*(ZS(I,J) - Z3d(I,J,1)))
2318 CWK4(I,J) = dum2d -PDTOP -PTOP
2319 ELSE ! target level bounded by input levels
2321 IF(ZS(I,J) .GE. Z3d(I,J,K) .AND. ZS(I,J) .LT. Z3d(I,J,K+1))THEN
2322 dlnpdz = (log(PSTD(K))-log(PSTD(K+1)) ) /(Z3d(I,J,K)-Z3d(I,J,K+1))
2323 dum2d = exp(log(PSTD(K)) + dlnpdz*(ZS(I,J)- Z3d(I,J,K)))
2324 CWK4(I,J) = dum2d -PDTOP -PTOP
2328 IF(ZS(I,J) .GE. Z3d(I,J,KZMAX))THEN
2329 WRITE(0,*)'I=',I,'J=',J,'K=',K,'TERRAIN HEIGHT',ZS(I,J),'Z3d',Z3d(I,J,KZMAX)
2330 CALL wrf_error_fatal("BC:MOUNTAIN TOO HIGH TO FIT THE MODEL DEPTH")
2335 DO I = NITS,MIN(NITE,NIDE-1)
2336 ntemp_b(i,j,k) = CWK4(I,J)
2337 ntemp_bt(i,j,k) = 0.0
2344 END SUBROUTINE nmm_bdymass_hinterp
2346 !=======================================================================================
2348 ! ADDED FOR INCLUDING MOISTURE AND THERMODYNAMIC ENERGY BALANCE
2350 !=======================================================================================
2352 SUBROUTINE interp_scalar_nmm (cfld, & ! CD field
2353 cids,cide,ckds,ckde,cjds,cjde, &
2354 cims,cime,ckms,ckme,cjms,cjme, &
2355 cits,cite,ckts,ckte,cjts,cjte, &
2357 nids,nide,nkds,nkde,njds,njde, &
2358 nims,nime,nkms,nkme,njms,njme, &
2359 nits,nite,nkts,nkte,njts,njte, &
2360 shw, & ! stencil half width for interp
2361 imask, & ! interpolation mask
2362 xstag,ystag, & ! staggering of field
2363 ipos,jpos, & ! Position of lower left of nest in CD
2364 nri,nrj, & ! nest ratios
2365 CII, IIH, CJJ, JJH, CBWGT1, HBWGT1, & ! south-western grid locs and weights
2366 CBWGT2, HBWGT2, CBWGT3, HBWGT3, & ! note that "C"ourse grid ones are
2367 CBWGT4, HBWGT4, & ! dummys for weights
2373 CETA1,ETA1,CETA2,ETA2 )
2375 USE MODULE_MODEL_CONSTANTS
2379 LOGICAL,INTENT(IN) :: xstag, ystag
2380 INTEGER,INTENT(IN) :: cids, cide, ckds, ckde, cjds, cjde, &
2381 cims, cime, ckms, ckme, cjms, cjme, &
2382 cits, cite, ckts, ckte, cjts, cjte, &
2383 nids, nide, nkds, nkde, njds, njde, &
2384 nims, nime, nkms, nkme, njms, njme, &
2385 nits, nite, nkts, nkte, njts, njte, &
2386 shw,ipos,jpos,nri,nrj
2388 INTEGER,DIMENSION(nims:nime,njms:njme), INTENT(IN) :: IMASK
2392 INTEGER,DIMENSION(cims:cime,cjms:cjme), INTENT(IN) :: CII,CJJ ! dummy
2393 REAL,DIMENSION(cims:cime,cjms:cjme), INTENT(IN) :: CBWGT1,CBWGT2
2394 REAL,DIMENSION(cims:cime,cjms:cjme), INTENT(IN) :: CBWGT3,CBWGT4
2396 REAL,DIMENSION(cims:cime,cjms:cjme,ckms:ckme), INTENT(IN) :: CFLD
2397 REAL,DIMENSION(cims:cime,cjms:cjme,ckms:ckme), INTENT(IN) :: CC3d ! scalar input on constant pressure levels
2398 REAL,DIMENSION(ckms:ckme), INTENT(IN) :: CPSTD
2399 REAL,DIMENSION(cims:cime,cjms:cjme), INTENT(IN) :: CPD
2400 REAL,DIMENSION(ckms:ckme), INTENT(IN) :: CETA1,CETA2
2401 REAL, INTENT(IN) :: CPDTOP,CPTOP
2405 INTEGER,DIMENSION(nims:nime,njms:njme), INTENT(IN) :: IIH,JJH
2406 REAL,DIMENSION(nims:nime,njms:njme), INTENT(IN) :: HBWGT1,HBWGT2
2407 REAL,DIMENSION(nims:nime,njms:njme), INTENT(IN) :: HBWGT3,HBWGT4
2409 REAL,DIMENSION(nims:nime,njms:njme,nkms:nkme), INTENT(OUT):: NFLD ! This is scalar on hybrid levels
2410 REAL,DIMENSION(nims:nime,njms:njme,nkms:nkme), INTENT(OUT):: C3d ! Scalar on constant pressure levels
2411 REAL,DIMENSION(nkms:nkme), INTENT(IN) :: PSTD
2412 REAL,DIMENSION(nims:nime,njms:njme), INTENT(IN) :: PD
2413 REAL,DIMENSION(nkms:nkme), INTENT(IN) :: ETA1,ETA2
2414 REAL,INTENT(IN) :: PDTOP,PTOP
2418 INTEGER,PARAMETER :: JTB=134
2420 REAL,DIMENSION(JTB) :: PIN,CIN,Y2,PIO,PTMP,COUT,DUM1,DUM2
2422 !-----------------------------------------------------------------------------------------------------
2425 ! *** CHECK VERTICAL BOUNDS BEFORE USING SPLINE OR LINEAR INTERPOLATION
2427 IF(nkme .GT. (JTB-10) .OR. NKDE .GT. (JTB-10)) &
2428 CALL wrf_error_fatal ('mass points: increase JTB in interp_mass_nmm')
2431 ! FIRST, HORIZONTALLY INTERPOLATE MOISTURE NOW AVAILABLE ON CONSTANT PRESSURE SURFACE (LEVELS) FROM THE
2432 ! PARENT TO THE NESTED DOMAIN
2434 !*** INDEX CONVENTIONS
2449 DO K=NKDS,NKDE-1 ! Please note that we are still in isobaric surfaces
2450 DO J=NJTS,MIN(NJTE,NJDE-1)
2451 DO I=NITS,MIN(NITE,NIDE-1)
2452 IF(IMASK(I,J) .NE. 1)THEN
2453 IF(MOD(JJH(I,J),2) .NE. 0)THEN ! 1,3,5,7
2454 C3d(I,J,K) = HBWGT1(I,J)*CC3d(IIH(I,J), JJH(I,J) ,K) &
2455 + HBWGT2(I,J)*CC3d(IIH(I,J)+1,JJH(I,J) ,K) &
2456 + HBWGT3(I,J)*CC3d(IIH(I,J), JJH(I,J)-1,K) &
2457 + HBWGT4(I,J)*CC3d(IIH(I,J), JJH(I,J)+1,K)
2460 C3d(I,J,K) = HBWGT1(I,J)*CC3d(IIH(I,J), JJH(I,J) ,K) &
2461 + HBWGT2(I,J)*CC3d(IIH(I,J)+1,JJH(I,J) ,K) &
2462 + HBWGT3(I,J)*CC3d(IIH(I,J)+1,JJH(I,J)-1,K) &
2463 + HBWGT4(I,J)*CC3d(IIH(I,J)+1,JJH(I,J)+1,K)
2472 ! RECOVER THE SCALARS FROM CONSTANT PRESSURE SURFACES (LEVELS) ON TO HYBRID SURFACES
2474 DO J=NJTS,MIN(NJTE,NJDE-1)
2475 DO I=NITS,MIN(NITE,NIDE-1)
2476 IF(IMASK(I,J) .NE. 1)THEN
2478 ! clean local array before use of spline or linear interpolation
2480 CIN=0.;PIN=0.;Y2=0;PIO=0.;PTMP=0.;COUT=0.;DUM1=0.;DUM2=0.
2482 DO K=NKDS+1,NKDE ! inputs at standard levels
2483 PIN(K-1) = EXP((ALOG(PSTD(NKDE-K+1))+ALOG(PSTD(NKDE-K+2)))*0.5)
2484 CIN(K-1) = C3d(I,J,NKDE-K+1)
2490 DO K=NKDS,NKDE ! target points in model interface levels (pint)
2491 PTMP(K) = ETA1(K)*PDTOP + ETA2(K)*PD(I,J) + PTOP
2494 DO K=NKDS,NKDE-1 ! target points in model levels
2495 PIO(K) = EXP((ALOG(PTMP(K))+ALOG(PTMP(K+1)))*0.5)
2498 IF(PTMP(1) .GE. PSTD(1))THEN ! if lower boundary is higher than PMSL(1) re-set lower boundary
2499 PIN(NKDE-1) = PIO(1) ! be consistent with target. This may not happen at all
2500 WRITE(0,*)'WARNING: NESTED DOMAIN PRESSURE AT LOWEST LEVEL HIGHER THAN PSTD'
2501 WRITE(0,*)'I,J,PIO(1),PSTD(1)',I,J,PIO(1),PSTD(1)
2504 CALL SPLINE2(I,J,JTB,NKDE-1,PIN,CIN,Y2,NKDE-1,PIO,COUT,DUM1,DUM2) ! interpolate
2507 NFLD(I,J,K)= COUT(K) ! scalar in the nested domain
2514 END SUBROUTINE interp_scalar_nmm
2516 !===========================================================================================
2518 SUBROUTINE nmm_bdy_scalar (cfld, & ! CD field
2519 cids,cide,ckds,ckde,cjds,cjde, &
2520 cims,cime,ckms,ckme,cjms,cjme, &
2521 cits,cite,ckts,ckte,cjts,cjte, &
2523 nids,nide,nkds,nkde,njds,njde, &
2524 nims,nime,nkms,nkme,njms,njme, &
2525 nits,nite,nkts,nkte,njts,njte, &
2526 shw, & ! stencil half width for interp
2527 imask, & ! interpolation mask
2528 xstag,ystag, & ! staggering of field
2529 ipos,jpos, & ! Position of lower left of nest in CD
2530 nri,nrj, & ! nest ratios
2540 CTEMP_B,NTEMP_B, & ! to be removed
2541 CTEMP_BT,NTEMP_BT, &
2542 CII, IIH, CJJ, JJH, CBWGT1, HBWGT1, & ! south-western grid locs and weights
2543 CBWGT2, HBWGT2, CBWGT3, HBWGT3, & ! note that "C"ourse grid ones are
2544 CBWGT4, HBWGT4, & ! dummys for weights
2550 CETA1,ETA1,CETA2,ETA2 )
2551 USE MODULE_MODEL_CONSTANTS
2555 LOGICAL,INTENT(IN) :: xstag, ystag
2556 REAL, INTENT(INOUT) :: cdt, ndt
2557 INTEGER,INTENT(IN) :: cids, cide, ckds, ckde, cjds, cjde, &
2558 cims, cime, ckms, ckme, cjms, cjme, &
2559 cits, cite, ckts, ckte, cjts, cjte, &
2560 nids, nide, nkds, nkde, njds, njde, &
2561 nims, nime, nkms, nkme, njms, njme, &
2562 nits, nite, nkts, nkte, njts, njte, &
2563 shw,ipos,jpos,nri,nrj
2564 REAL, DIMENSION ( cims:cime, cjms:cjme, ckms:ckme ), INTENT(OUT) :: ctemp_b,ctemp_bt
2565 REAL, DIMENSION ( nims:nime, njms:njme, nkms:nkme ), INTENT(OUT) :: ntemp_b,ntemp_bt
2566 REAL, DIMENSION( * ), INTENT(INOUT) :: c_bxs,n_bxs,c_bxe,n_bxe,c_bys,n_bys,c_bye,n_bye
2567 REAL, DIMENSION( * ), INTENT(INOUT) :: c_btxs,n_btxs,c_btxe,n_btxe,c_btys,n_btys,c_btye,n_btye
2570 INTEGER,DIMENSION(nims:nime,njms:njme), INTENT(IN) :: IMASK
2574 INTEGER,DIMENSION(cims:cime,cjms:cjme), INTENT(IN) :: CII,CJJ ! dummy
2575 REAL,DIMENSION(cims:cime,cjms:cjme), INTENT(IN) :: CBWGT1,CBWGT2
2576 REAL,DIMENSION(cims:cime,cjms:cjme), INTENT(IN) :: CBWGT3,CBWGT4
2577 REAL,DIMENSION(cims:cime,cjms:cjme,ckms:ckme), INTENT(IN) :: CFLD
2578 REAL,DIMENSION(cims:cime,cjms:cjme,ckms:ckme), INTENT(IN) :: CC3d ! scalar input on constant pressure levels
2579 REAL,DIMENSION(ckms:ckme), INTENT(IN) :: CPSTD
2580 REAL,DIMENSION(cims:cime,cjms:cjme), INTENT(IN) :: CPD
2581 REAL,DIMENSION(ckms:ckme), INTENT(IN) :: CETA1,CETA2
2582 REAL, INTENT(IN) :: CPDTOP,CPTOP
2586 INTEGER,DIMENSION(nims:nime,njms:njme), INTENT(IN) :: IIH,JJH
2587 REAL,DIMENSION(nims:nime,njms:njme), INTENT(IN) :: HBWGT1,HBWGT2
2588 REAL,DIMENSION(nims:nime,njms:njme), INTENT(IN) :: HBWGT3,HBWGT4
2589 REAL,DIMENSION(nims:nime,njms:njme,nkms:nkme), INTENT(OUT):: NFLD
2590 REAL,DIMENSION(nims:nime,njms:njme,nkms:nkme), INTENT(OUT):: C3d !Scalar on constant pressure levels
2591 REAL,DIMENSION(nkms:nkme), INTENT(IN) :: PSTD
2592 REAL,DIMENSION(nims:nime,njms:njme), INTENT(IN) :: PD
2593 REAL,DIMENSION(nkms:nkme), INTENT(IN) :: ETA1,ETA2
2594 REAL,INTENT(IN) :: PDTOP,PTOP
2598 INTEGER,PARAMETER :: JTB=134
2599 INTEGER :: I,J,K,II,JJ
2600 REAL,DIMENSION(JTB) :: PIN,CIN,Y2,PIO,PTMP,COUT,DUM1,DUM2
2601 REAL, DIMENSION (nims:nime,njms:njme,nkms:nkme) :: CWK1,CWK2,CWK3,CWK4
2602 !-----------------------------------------------------------------------------------------------------
2605 ! *** CHECK VERTICAL BOUNDS BEFORE USING SPLINE INTERPOLATION
2607 IF(nkme .GT. (JTB-10) .OR. NKDE .GT. (JTB-10)) &
2608 CALL wrf_error_fatal ('mass points: increase JTB in interp_mass_nmm')
2612 NMM_XS: IF(NITS .EQ. NIDS)THEN
2613 ! WRITE(0,*)'ENTERING X1 START BOUNDARY AT T POINTS',NJTS,MIN(NJTE,NJDE-1)
2615 DO K=NKDS,NKDE-1 ! Please note that we are still in isobaric surfaces
2616 DO J = NJTS,MIN(NJTE,NJDE-1)
2617 IF(MOD(JJH(I,J),2) .NE. 0)THEN ! 1,3,5,7 of the parent domain
2618 C3d(I,J,K) = HBWGT1(I,J)*CC3d(IIH(I,J), JJH(I,J) ,K) &
2619 + HBWGT2(I,J)*CC3d(IIH(I,J)+1,JJH(I,J) ,K) &
2620 + HBWGT3(I,J)*CC3d(IIH(I,J), JJH(I,J)-1,K) &
2621 + HBWGT4(I,J)*CC3d(IIH(I,J), JJH(I,J)+1,K)
2623 C3d(I,J,K) = HBWGT1(I,J)*CC3d(IIH(I,J), JJH(I,J) ,K) &
2624 + HBWGT2(I,J)*CC3d(IIH(I,J)+1,JJH(I,J) ,K) &
2625 + HBWGT3(I,J)*CC3d(IIH(I,J)+1,JJH(I,J)-1,K) &
2626 + HBWGT4(I,J)*CC3d(IIH(I,J)+1,JJH(I,J)+1,K)
2631 DO J=NJTS,MIN(NJTE,NJDE-1)
2632 IF(MOD(J,2) .NE. 0)THEN
2633 CIN=0.;PIN=0.;Y2=0;PIO=0.;PTMP=0.;COUT=0.;DUM1=0.;DUM2=0. ! clean up local array
2634 DO K=NKDS+1,NKDE ! inputs at standard levels
2635 PIN(K-1) = EXP((ALOG(PSTD(NKDE-K+1))+ALOG(PSTD(NKDE-K+2)))*0.5)
2636 CIN(K-1) = C3d(I,J,NKDE-K+1)
2640 DO K=NKDS,NKDE ! target points in model interface levels (pint)
2641 PTMP(K) = ETA1(K)*PDTOP + ETA2(K)*PD(I,J) + PTOP
2643 DO K=NKDS,NKDE-1 ! target points in model levels
2644 PIO(K) = EXP((ALOG(PTMP(K))+ALOG(PTMP(K+1)))*0.5)
2646 IF(PTMP(1) .GE. PSTD(1))THEN ! if lower boundary is higher than PMSL(1) re-set lower boundary
2647 PIN(NKDE-1) = PIO(1) ! be consistent with target. This may not happen at all
2648 WRITE(0,*)'WARNING: NESTED DOMAIN PRESSURE AT LOWEST LEVEL HIGHER THAN PSTD'
2649 WRITE(0,*)'I,J,PIO(1),PSTD(1)',I,J,PIO(1),PSTD(1)
2652 CALL SPLINE2(I,J,JTB,NKDE-1,PIN,CIN,Y2,NKDE-1,PIO,COUT,DUM1,DUM2) ! interpolate
2655 CWK1(I,J,K)= COUT(K) ! scalar in the nested domain
2664 DO J = NJTS,MIN(NJTE,NJDE-1)
2666 ntemp_b(i,j,k) = CWK1(I,J,K)
2667 ntemp_bt(i,j,k) = 0.0
2676 NMM_XE: IF(NITE-1 .EQ. NIDE-1)THEN
2677 ! WRITE(0,*)'ENTERING X END BOUNDARY AT T POINTS',NJTS,MIN(NJTE,NJDE-1)
2679 DO K=NKDS,NKDE-1 ! Please note that we are still in isobaric surfaces
2680 DO J = NJTS,MIN(NJTE,NJDE-1)
2681 IF(MOD(JJH(I,J),2) .NE. 0)THEN ! 1,3,5,7 of the parent domain
2682 C3d(I,J,K) = HBWGT1(I,J)*CC3d(IIH(I,J), JJH(I,J) ,K) &
2683 + HBWGT2(I,J)*CC3d(IIH(I,J)+1,JJH(I,J) ,K) &
2684 + HBWGT3(I,J)*CC3d(IIH(I,J), JJH(I,J)-1,K) &
2685 + HBWGT4(I,J)*CC3d(IIH(I,J), JJH(I,J)+1,K)
2687 C3d(I,J,K) = HBWGT1(I,J)*CC3d(IIH(I,J), JJH(I,J) ,K) &
2688 + HBWGT2(I,J)*CC3d(IIH(I,J)+1,JJH(I,J) ,K) &
2689 + HBWGT3(I,J)*CC3d(IIH(I,J)+1,JJH(I,J)-1,K) &
2690 + HBWGT4(I,J)*CC3d(IIH(I,J)+1,JJH(I,J)+1,K)
2695 DO J=NJTS,MIN(NJTE,NJDE-1)
2696 IF(MOD(J,2) .NE. 0)THEN
2697 CIN=0.;PIN=0.;Y2=0;PIO=0.;PTMP=0.;COUT=0.;DUM1=0.;DUM2=0. ! clean up local array
2698 DO K=NKDS+1,NKDE ! inputs at standard levels
2699 PIN(K-1) = EXP((ALOG(PSTD(NKDE-K+1))+ALOG(PSTD(NKDE-K+2)))*0.5)
2700 CIN(K-1) = C3d(I,J,NKDE-K+1)
2704 DO K=NKDS,NKDE ! target points in model interface levels (pint)
2705 PTMP(K) = ETA1(K)*PDTOP + ETA2(K)*PD(I,J) + PTOP
2707 DO K=NKDS,NKDE-1 ! target points in model levels
2708 PIO(K) = EXP((ALOG(PTMP(K))+ALOG(PTMP(K+1)))*0.5)
2710 IF(PTMP(1) .GE. PSTD(1))THEN ! if lower boundary is higher than PMSL(1) re-set lower boundary
2711 PIN(NKDE-1) = PIO(1) ! be consistent with target. This may not happen at all
2712 WRITE(0,*)'WARNING: NESTED DOMAIN PRESSURE AT LOWEST LEVEL HIGHER THAN PSTD'
2713 WRITE(0,*)'I,J,PIO(1),PSTD(1)',I,J,PIO(1),PSTD(1)
2716 CALL SPLINE2(I,J,JTB,NKDE-1,PIN,CIN,Y2,NKDE-1,PIO,COUT,DUM1,DUM2) ! interpolate
2719 CWK2(I,J,K)= COUT(K) ! scalar in the nested domain
2728 DO J = NJTS,MIN(NJTE,NJDE-1)
2729 DO K = NKDS,MIN(NKTE,NKDE-1)
2730 ntemp_b(i,j,k) = CWK2(I,J,K)
2731 ntemp_bt(i,j,k) = 0.0
2739 NMM_YS: IF(NJTS .EQ. NJDS)THEN
2740 ! WRITE(0,*)'ENTERING Y START BOUNDARY AT T POINTS',NITS,MIN(NITE,NIDE-1)
2743 DO I = NITS,MIN(NITE,NIDE-1)
2744 IF(MOD(JJH(I,J),2) .NE. 0)THEN ! 1,3,5,7 of the parent domain
2745 C3d(I,J,K) = HBWGT1(I,J)*CC3d(IIH(I,J), JJH(I,J) ,K) &
2746 + HBWGT2(I,J)*CC3d(IIH(I,J)+1,JJH(I,J) ,K) &
2747 + HBWGT3(I,J)*CC3d(IIH(I,J), JJH(I,J)-1,K) &
2748 + HBWGT4(I,J)*CC3d(IIH(I,J), JJH(I,J)+1,K)
2750 C3d(I,J,K) = HBWGT1(I,J)*CC3d(IIH(I,J), JJH(I,J) ,K) &
2751 + HBWGT2(I,J)*CC3d(IIH(I,J)+1,JJH(I,J) ,K) &
2752 + HBWGT3(I,J)*CC3d(IIH(I,J)+1,JJH(I,J)-1,K) &
2753 + HBWGT4(I,J)*CC3d(IIH(I,J)+1,JJH(I,J)+1,K)
2759 DO I=NITS,MIN(NITE,NIDE-1)
2760 CIN=0.;PIN=0.;Y2=0;PIO=0.;PTMP=0.;COUT=0.;DUM1=0.;DUM2=0. ! clean up local array
2761 DO K=NKDS+1,NKDE ! inputs at standard levels
2762 PIN(K-1) = EXP((ALOG(PSTD(NKDE-K+1))+ALOG(PSTD(NKDE-K+2)))*0.5)
2763 CIN(K-1) = C3d(I,J,NKDE-K+1)
2767 DO K=NKDS,NKDE ! target points in model interface levels (pint)
2768 PTMP(K) = ETA1(K)*PDTOP + ETA2(K)*PD(I,J) + PTOP
2770 DO K=NKDS,NKDE-1 ! target points in model levels
2771 PIO(K) = EXP((ALOG(PTMP(K))+ALOG(PTMP(K+1)))*0.5)
2773 IF(PTMP(1) .GE. PSTD(1))THEN ! if lower boundary is higher than PMSL(1) re-set lower boundary
2774 PIN(NKDE-1) = PIO(1) ! be consistent with target. This may not happen at all
2775 WRITE(0,*)'WARNING: NESTED DOMAIN PRESSURE AT LOWEST LEVEL HIGHER THAN PSTD'
2776 WRITE(0,*)'I,J,PIO(1),PSTD(1)',I,J,PIO(1),PSTD(1)
2779 CALL SPLINE2(I,J,JTB,NKDE-1,PIN,CIN,Y2,NKDE-1,PIO,COUT,DUM1,DUM2) ! interpolate
2782 CWK3(I,J,K)= COUT(K) ! scalar in the nested domain
2787 DO I = NITS,MIN(NITE,NIDE-1)
2788 ntemp_b(i,J,K) = CWK3(I,J,K)
2789 ntemp_bt(i,J,K) = 0.0
2797 NMM_YE: IF(NJTE-1 .EQ. NJDE-1)THEN
2798 ! WRITE(0,*)'ENTERING Y END BOUNDARY AT T POINTS',NITS,MIN(NITE,NIDE-1)
2801 DO I = NITS,MIN(NITE,NIDE-1)
2802 IF(MOD(JJH(I,J),2) .NE. 0)THEN ! 1,3,5,7 of the parent domain
2803 C3d(I,J,K) = HBWGT1(I,J)*CC3d(IIH(I,J), JJH(I,J) ,K) &
2804 + HBWGT2(I,J)*CC3d(IIH(I,J)+1,JJH(I,J) ,K) &
2805 + HBWGT3(I,J)*CC3d(IIH(I,J), JJH(I,J)-1,K) &
2806 + HBWGT4(I,J)*CC3d(IIH(I,J), JJH(I,J)+1,K)
2808 C3d(I,J,K) = HBWGT1(I,J)*CC3d(IIH(I,J), JJH(I,J) ,K) &
2809 + HBWGT2(I,J)*CC3d(IIH(I,J)+1,JJH(I,J) ,K) &
2810 + HBWGT3(I,J)*CC3d(IIH(I,J)+1,JJH(I,J)-1,K) &
2811 + HBWGT4(I,J)*CC3d(IIH(I,J)+1,JJH(I,J)+1,K)
2817 DO I=NITS,MIN(NITE,NIDE-1)
2818 CIN=0.;PIN=0.;Y2=0;PIO=0.;PTMP=0.;COUT=0.;DUM1=0.;DUM2=0. ! clean up local array
2819 DO K=NKDS+1,NKDE ! inputs at standard levels
2820 PIN(K-1) = EXP((ALOG(PSTD(NKDE-K+1))+ALOG(PSTD(NKDE-K+2)))*0.5)
2821 CIN(K-1) = C3d(I,J,NKDE-K+1)
2825 DO K=NKDS,NKDE ! target points in model interface levels (pint)
2826 PTMP(K) = ETA1(K)*PDTOP + ETA2(K)*PD(I,J) + PTOP
2828 DO K=NKDS,NKDE-1 ! target points in model levels
2829 PIO(K) = EXP((ALOG(PTMP(K))+ALOG(PTMP(K+1)))*0.5)
2831 IF(PTMP(1) .GE. PSTD(1))THEN ! if lower boundary is higher than PMSL(1) re-set lower boundary
2832 PIN(NKDE-1) = PIO(1) ! be consistent with target. This may not happen at all
2833 WRITE(0,*)'WARNING: NESTED DOMAIN PRESSURE AT LOWEST LEVEL HIGHER THAN PSTD'
2834 WRITE(0,*)'I,J,PIO(1),PSTD(1)',I,J,PIO(1),PSTD(1)
2837 CALL SPLINE2(I,J,JTB,NKDE-1,PIN,CIN,Y2,NKDE-1,PIO,COUT,DUM1,DUM2) ! interpolate
2840 CWK4(I,J,K)= COUT(K) ! scalar in the nested domain
2845 DO I = NITS,MIN(NITE,NIDE-1)
2846 ntemp_b(i,J,K) = CWK4(I,J,K)
2847 ntemp_bt(i,J,K) = 0.0
2853 END SUBROUTINE nmm_bdy_scalar
2856 !=======================================================================================
2857 SUBROUTINE SPLINE2(I,J,JTBX,NOLD,XOLD,YOLD,Y2,NNEW,XNEW,YNEW,P,Q)
2859 ! ******************************************************************
2861 ! * THIS IS A ONE-DIMENSIONAL CUBIC SPLINE FITTING ROUTINE *
2862 ! * PROGRAMED FOR A SMALL SCALAR MACHINE. *
2864 ! * PROGRAMER Z. JANJIC *
2866 ! * NOLD - NUMBER OF GIVEN VALUES OF THE FUNCTION. MUST BE GE 3. *
2867 ! * XOLD - LOCATIONS OF THE POINTS AT WHICH THE VALUES OF THE *
2868 ! * FUNCTION ARE GIVEN. MUST BE IN ASCENDING ORDER. *
2869 ! * YOLD - THE GIVEN VALUES OF THE FUNCTION AT THE POINTS XOLD. *
2870 ! * Y2 - THE SECOND DERIVATIVES AT THE POINTS XOLD. IF NATURAL *
2871 ! * SPLINE IS FITTED Y2(1)=0. AND Y2(NOLD)=0. MUST BE *
2873 ! * NNEW - NUMBER OF VALUES OF THE FUNCTION TO BE CALCULATED. *
2874 ! * XNEW - LOCATIONS OF THE POINTS AT WHICH THE VALUES OF THE *
2875 ! * FUNCTION ARE CALCULATED. XNEW(K) MUST BE GE XOLD(1) *
2876 ! * AND LE XOLD(NOLD). *
2877 ! * YNEW - THE VALUES OF THE FUNCTION TO BE CALCULATED. *
2878 ! * P, Q - AUXILIARY VECTORS OF THE LENGTH NOLD-2. *
2880 ! ******************************************************************
2881 !---------------------------------------------------------------------
2883 !---------------------------------------------------------------------
2884 INTEGER,INTENT(IN) :: I,J,JTBX,NNEW,NOLD
2885 REAL,DIMENSION(JTBX),INTENT(IN) :: XNEW,XOLD,YOLD
2886 REAL,DIMENSION(JTBX),INTENT(INOUT) :: P,Q,Y2
2887 REAL,DIMENSION(JTBX),INTENT(OUT) :: YNEW
2889 INTEGER :: II,JJ,K,K1,K2,KOLD,NOLDM1
2890 REAL :: AK,BK,CK,DEN,DX,DXC,DXL,DXR,DYDXL,DYDXR &
2891 ,RDX,RTDXC,X,XK,XSQ,Y2K,Y2KP1
2892 !---------------------------------------------------------------------
2898 IF(I.eq.II.and.J.eq.JJ)THEN
2899 WRITE(0,*)'DEBUG in SPLINE2: I,J',I,J
2900 WRITE(0,*)'DEBUG in SPLINE2:HSO= ',xnew(1:nold)
2902 WRITE(0,*)'DEBUG in SPLINE2:L,ZETAI,PINTI= ' &
2911 DYDXL=(YOLD(2)-YOLD(1))/DXL
2912 DYDXR=(YOLD(3)-YOLD(2))/DXR
2915 P(1)= RTDXC*(6.*(DYDXR-DYDXL)-DXL*Y2(1))
2918 IF(NOLD.EQ.3)GO TO 150
2919 !---------------------------------------------------------------------
2924 DXR=XOLD(K+1)-XOLD(K)
2925 DYDXR=(YOLD(K+1)-YOLD(K))/DXR
2927 DEN=1./(DXL*Q(K-2)+DXC+DXC)
2929 P(K-1)= DEN*(6.*(DYDXR-DYDXL)-DXL*P(K-2))
2933 IF(K.LT.NOLD)GO TO 100
2934 !-----------------------------------------------------------------------
2937 200 Y2(K)=P(K-1)+Q(K-1)*Y2(K+1)
2941 !-----------------------------------------------------------------------
2948 IF(XOLD(K2).GT.XK)THEN
2958 450 IF(K1.EQ.1)GO TO 500
2959 IF(K.EQ.KOLD)GO TO 550
2965 DX=XOLD(K+1)-XOLD(K)
2968 AK=.1666667*RDX*(Y2KP1-Y2K)
2970 CK=RDX*(YOLD(K+1)-YOLD(K))-.1666667*DX*(Y2KP1+Y2K+Y2K)
2975 YNEW(K1)=AK*XSQ*X+BK*XSQ+CK*X+YOLD(K)
2979 IF(I.eq.II.and.J.eq.JJ)THEN
2980 WRITE(0,*) 'DEBUG:: k1,xnew(k1),ynew(k1): ', K1,XNEW(k1),YNEW(k1)
2985 IF(K1.LE.NNEW)GO TO 300
2989 END SUBROUTINE SPLINE2
2991 !=======================================================================================
2992 ! E grid interpolation for H and V points
2993 !=======================================================================================
2995 SUBROUTINE interp_h_nmm (cfld, & ! CD field
2996 cids, cide, ckds, ckde, cjds, cjde, &
2997 cims, cime, ckms, ckme, cjms, cjme, &
2998 cits, cite, ckts, ckte, cjts, cjte, &
3000 nids, nide, nkds, nkde, njds, njde, &
3001 nims, nime, nkms, nkme, njms, njme, &
3002 nits, nite, nkts, nkte, njts, njte, &
3003 shw, & ! stencil half width for interp
3004 imask, & ! interpolation mask
3005 xstag, ystag, & ! staggering of field
3006 ipos, jpos, & ! Position of lower left of nest in CD
3007 nri, nrj, & ! nest ratios
3008 CII, IIH, CJJ, JJH, CBWGT1, HBWGT1, & ! south-western grid locs and weights
3009 CBWGT2, HBWGT2, CBWGT3, HBWGT3, & ! note that "C"ourse grid ones are
3010 CBWGT4, HBWGT4 ) ! dummys for weights
3014 INTEGER, INTENT(IN) :: cids, cide, ckds, ckde, cjds, cjde, &
3015 cims, cime, ckms, ckme, cjms, cjme, &
3016 cits, cite, ckts, ckte, cjts, cjte, &
3017 nids, nide, nkds, nkde, njds, njde, &
3018 nims, nime, nkms, nkme, njms, njme, &
3019 nits, nite, nkts, nkte, njts, njte, &
3023 LOGICAL, INTENT(IN) :: xstag, ystag
3025 REAL, DIMENSION ( cims:cime, cjms:cjme, ckms:ckme ) :: cfld
3026 REAL, DIMENSION ( nims:nime, njms:njme, nkms:nkme ) :: nfld
3027 REAL, DIMENSION ( cims:cime, cjms:cjme ), INTENT(IN) :: CBWGT1,CBWGT2,CBWGT3,CBWGT4 ! dummy
3028 REAL, DIMENSION ( nims:nime, njms:njme ), INTENT(IN) :: HBWGT1,HBWGT2,HBWGT3,HBWGT4
3029 INTEGER, DIMENSION ( cims:cime, cjms:cjme ), INTENT(IN) :: CII,CJJ ! dummy
3030 INTEGER, DIMENSION ( nims:nime, njms:njme ), INTENT(IN) :: IIH,JJH
3031 INTEGER, DIMENSION ( nims:nime, njms:njme ) :: imask
3036 !*** CHECK DOMAIN BOUNDS BEFORE INTERPOLATION
3038 DO J=NJTS,MIN(NJTE,NJDE-1)
3039 DO I=NITS,MIN(NITE,NIDE-1)
3040 IF(IIH(i,j).LT.(CIDS-shw) .OR. IIH(i,j).GT.(CIDE+shw)) &
3041 CALL wrf_error_fatal ('hpoints:check domain bounds along x' )
3042 IF(JJH(i,j).LT.(CJDS-shw) .OR. JJH(i,j).GT.(CJDE+shw)) &
3043 CALL wrf_error_fatal ('hpoints:check domain bounds along y' )
3047 !*** INDEX CONVENTIONS
3062 DO J=NJTS,MIN(NJTE,NJDE-1)
3063 DO I=NITS,MIN(NITE,NIDE-1)
3064 IF(IMASK(I,J) .NE. 1)THEN
3066 IF(MOD(JJH(I,J),2) .NE. 0)THEN ! 1,3,5,7
3067 NFLD(I,J,K) = HBWGT1(I,J)*CFLD(IIH(I,J), JJH(I,J) ,K) &
3068 + HBWGT2(I,J)*CFLD(IIH(I,J)+1,JJH(I,J) ,K) &
3069 + HBWGT3(I,J)*CFLD(IIH(I,J), JJH(I,J)-1,K) &
3070 + HBWGT4(I,J)*CFLD(IIH(I,J), JJH(I,J)+1,K)
3072 NFLD(I,J,K) = HBWGT1(I,J)*CFLD(IIH(I,J), JJH(I,J) ,K) &
3073 + HBWGT2(I,J)*CFLD(IIH(I,J)+1,JJH(I,J) ,K) &
3074 + HBWGT3(I,J)*CFLD(IIH(I,J)+1,JJH(I,J)-1,K) &
3075 + HBWGT4(I,J)*CFLD(IIH(I,J)+1,JJH(I,J)+1,K)
3083 END SUBROUTINE interp_h_nmm
3085 SUBROUTINE interp_v_nmm (cfld, & ! CD field
3086 cids, cide, ckds, ckde, cjds, cjde, &
3087 cims, cime, ckms, ckme, cjms, cjme, &
3088 cits, cite, ckts, ckte, cjts, cjte, &
3090 nids, nide, nkds, nkde, njds, njde, &
3091 nims, nime, nkms, nkme, njms, njme, &
3092 nits, nite, nkts, nkte, njts, njte, &
3093 shw, & ! stencil half width for interp
3094 imask, & ! interpolation mask
3095 xstag, ystag, & ! staggering of field
3096 ipos, jpos, & ! Position of lower left of nest in CD
3097 nri, nrj, & ! nest ratios
3098 CII, IIV, CJJ, JJV, CBWGT1, VBWGT1, & ! south-western grid locs and weights
3099 CBWGT2, VBWGT2, CBWGT3, VBWGT3, & ! note that "C"ourse grid ones are
3100 CBWGT4, VBWGT4 ) ! dummys
3104 INTEGER, INTENT(IN) :: cids, cide, ckds, ckde, cjds, cjde, &
3105 cims, cime, ckms, ckme, cjms, cjme, &
3106 cits, cite, ckts, ckte, cjts, cjte, &
3107 nids, nide, nkds, nkde, njds, njde, &
3108 nims, nime, nkms, nkme, njms, njme, &
3109 nits, nite, nkts, nkte, njts, njte, &
3113 LOGICAL, INTENT(IN) :: xstag, ystag
3115 REAL, DIMENSION ( cims:cime, cjms:cjme, ckms:ckme ) :: cfld
3116 REAL, DIMENSION ( nims:nime, njms:njme, nkms:nkme ) :: nfld
3117 REAL, DIMENSION ( cims:cime, cjms:cjme ), INTENT(IN) :: CBWGT1,CBWGT2,CBWGT3,CBWGT4 ! dummy
3118 REAL, DIMENSION ( nims:nime, njms:njme ), INTENT(IN) :: VBWGT1,VBWGT2,VBWGT3,VBWGT4
3119 INTEGER, DIMENSION ( cims:cime, cjms:cjme ), INTENT(IN) :: CII,CJJ ! dummy
3120 INTEGER, DIMENSION ( nims:nime, njms:njme ), INTENT(IN) :: IIV,JJV
3121 INTEGER, DIMENSION ( nims:nime, njms:njme ) :: imask
3128 !*** CHECK DOMAIN BOUNDS BEFORE INTERPOLATION
3130 DO J=NJTS,MIN(NJTE,NJDE-1)
3131 DO I=NITS,MIN(NITE,NIDE-1)
3132 IF(IIV(i,j).LT.(CIDS-shw) .OR. IIV(i,j).GT.(CIDE+shw)) &
3133 CALL wrf_error_fatal ('vpoints:check domain bounds along x' )
3134 IF(JJV(i,j).LT.(CJDS-shw) .OR. JJV(i,j).GT.(CJDE+shw)) &
3135 CALL wrf_error_fatal ('vpoints:check domain bounds along y' )
3139 !*** INDEX CONVENTIONS
3154 DO J=NJTS,MIN(NJTE,NJDE-1)
3155 DO I=NITS,MIN(NITE,NIDE-1)
3156 IF(IMASK(I,J) .NE. 1)THEN
3158 IF(MOD(JJV(I,J),2) .NE. 0)THEN ! 1,3,5,7
3159 NFLD(I,J,K) = VBWGT1(I,J)*CFLD(IIV(I,J), JJV(I,J) ,K) &
3160 + VBWGT2(I,J)*CFLD(IIV(I,J)+1,JJV(I,J) ,K) &
3161 + VBWGT3(I,J)*CFLD(IIV(I,J)+1,JJV(I,J)-1,K) &
3162 + VBWGT4(I,J)*CFLD(IIV(I,J)+1,JJV(I,J)+1,K)
3164 NFLD(I,J,K) = VBWGT1(I,J)*CFLD(IIV(I,J), JJV(I,J) ,K) &
3165 + VBWGT2(I,J)*CFLD(IIV(I,J)+1,JJV(I,J) ,K) &
3166 + VBWGT3(I,J)*CFLD(IIV(I,J), JJV(I,J)-1,K) &
3167 + VBWGT4(I,J)*CFLD(IIV(I,J), JJV(I,J)+1,K)
3175 END SUBROUTINE interp_v_nmm
3177 !=======================================================================================
3178 ! E grid nearest neighbour interpolation for H points.
3179 ! This routine assumes cfld and nfld are in IJK
3180 !=======================================================================================
3182 SUBROUTINE interp_hnear_nmm (cfld, & ! CD field
3183 cids, cide, ckds, ckde, cjds, cjde, &
3184 cims, cime, ckms, ckme, cjms, cjme, &
3185 cits, cite, ckts, ckte, cjts, cjte, &
3187 nids, nide, nkds, nkde, njds, njde, &
3188 nims, nime, nkms, nkme, njms, njme, &
3189 nits, nite, nkts, nkte, njts, njte, &
3190 shw, & ! stencil half width for interp
3191 imask, & ! interpolation mask
3192 xstag, ystag, & ! staggering of field
3193 ipos, jpos, & ! Position of lower left of nest in CD
3194 nri, nrj, & ! nest ratios
3195 CII, IIH, CJJ, JJH, CBWGT1, HBWGT1, & ! south-western grid locs and weights
3196 CBWGT2, HBWGT2, CBWGT3, HBWGT3, & ! note that "C"ourse grid ones are
3197 CBWGT4, HBWGT4 ) ! just dummys
3201 INTEGER, INTENT(IN) :: cids, cide, ckds, ckde, cjds, cjde, &
3202 cims, cime, ckms, ckme, cjms, cjme, &
3203 cits, cite, ckts, ckte, cjts, cjte, &
3204 nids, nide, nkds, nkde, njds, njde, &
3205 nims, nime, nkms, nkme, njms, njme, &
3206 nits, nite, nkts, nkte, njts, njte, &
3210 LOGICAL, INTENT(IN) :: xstag, ystag
3212 REAL, DIMENSION ( cims:cime, cjms:cjme, ckms:ckme ) :: cfld
3213 REAL, DIMENSION ( nims:nime, njms:njme, nkms:nkme ) :: nfld
3214 REAL, DIMENSION ( cims:cime, cjms:cjme ), INTENT(IN) :: CBWGT1,CBWGT2,CBWGT3,CBWGT4 ! dummy
3215 REAL, DIMENSION ( nims:nime, njms:njme ), INTENT(IN) :: HBWGT1,HBWGT2,HBWGT3,HBWGT4
3216 INTEGER, DIMENSION ( cims:cime, cjms:cjme ), INTENT(IN) :: CII,CJJ ! dummy
3217 INTEGER, DIMENSION ( nims:nime, njms:njme ), INTENT(IN) :: IIH,JJH
3218 INTEGER, DIMENSION ( nims:nime, njms:njme ) :: imask
3225 REAL, DIMENSION (4, nims:nime, njms:njme ) :: NBWGT
3228 !*** INDEX CONVENTIONS
3236 !*** NBWGT1=1 NBWGT2=0
3242 DO J=NJTS,MIN(NJTE,NJDE-1)
3243 DO I=NITS,MIN(NITE,NIDE-1)
3244 IF(IMASK(I,J) .NE. 1)THEN
3245 NBWGT(1,I,J)=HBWGT1(I,J)
3246 NBWGT(2,I,J)=HBWGT2(I,J)
3247 NBWGT(3,I,J)=HBWGT3(I,J)
3248 NBWGT(4,I,J)=HBWGT4(I,J)
3253 DO J=NJTS,MIN(NJTE,NJDE-1)
3254 DO I=NITS,MIN(NITE,NIDE-1)
3255 IF(IMASK(I,J) .NE. 1)THEN
3259 AMAXVAL=amax1(NBWGT(N,I,J),AMAXVAL)
3265 IF(AMAXVAL .EQ. NBWGT(N,I,J) .AND. FLIP)THEN
3271 SUM=SUM+NBWGT(N,I,J)
3272 IF(SUM .GT. 1.0)CALL wrf_error_fatal ( "horizontal interp error - interp_hnear_nmm" )
3280 DO J=NJTS,MIN(NJTE,NJDE-1)
3281 DO I=NITS,MIN(NITE,NIDE-1)
3282 IF(IMASK(I,J) .NE. 1)THEN
3283 IF(MOD(JJH(I,J),2) .NE. 0)THEN ! 1,3,5,7
3284 NFLD(I,J,K) = NBWGT(1,I,J)*CFLD(IIH(I,J), JJH(I,J) ,K) &
3285 + NBWGT(2,I,J)*CFLD(IIH(I,J)+1,JJH(I,J) ,K) &
3286 + NBWGT(3,I,J)*CFLD(IIH(I,J), JJH(I,J)-1,K) &
3287 + NBWGT(4,I,J)*CFLD(IIH(I,J), JJH(I,J)+1,K)
3289 NFLD(I,J,K) = NBWGT(1,I,J)*CFLD(IIH(I,J), JJH(I,J) ,K) &
3290 + NBWGT(2,I,J)*CFLD(IIH(I,J)+1,JJH(I,J) ,K) &
3291 + NBWGT(3,I,J)*CFLD(IIH(I,J)+1,JJH(I,J)-1,K) &
3292 + NBWGT(4,I,J)*CFLD(IIH(I,J)+1,JJH(I,J)+1,K)
3299 END SUBROUTINE interp_hnear_nmm
3301 !=======================================================================================
3302 ! E grid nearest neighbour interpolation for H points.
3303 ! This routine assumes cfld and nfld are in IKJ or ILJ
3304 !=======================================================================================
3306 SUBROUTINE interp_hnear_ikj_nmm (cfld, & ! CD field
3307 cids, cide, ckds, ckde, cjds, cjde, &
3308 cims, cime, ckms, ckme, cjms, cjme, &
3309 cits, cite, ckts, ckte, cjts, cjte, &
3311 nids, nide, nkds, nkde, njds, njde, &
3312 nims, nime, nkms, nkme, njms, njme, &
3313 nits, nite, nkts, nkte, njts, njte, &
3314 shw, & ! stencil half width for interp
3315 imask, & ! interpolation mask
3316 xstag, ystag, & ! staggering of field
3317 ipos, jpos, & ! Position of lower left of nest in CD
3318 nri, nrj, & ! nest ratios
3319 CII, IIH, CJJ, JJH, CBWGT1, HBWGT1, & ! south-western grid locs and weights
3320 CBWGT2, HBWGT2, CBWGT3, HBWGT3, & ! note that "C"ourse grid ones are
3321 CBWGT4, HBWGT4 ) ! just dummys
3325 INTEGER, INTENT(IN) :: cids, cide, ckds, ckde, cjds, cjde, &
3326 cims, cime, ckms, ckme, cjms, cjme, &
3327 cits, cite, ckts, ckte, cjts, cjte, &
3328 nids, nide, nkds, nkde, njds, njde, &
3329 nims, nime, nkms, nkme, njms, njme, &
3330 nits, nite, nkts, nkte, njts, njte, &
3334 LOGICAL, INTENT(IN) :: xstag, ystag
3336 REAL, DIMENSION ( cims:cime, ckms:ckme, cjms:cjme ) :: cfld
3337 REAL, DIMENSION ( nims:nime, nkms:nkme, njms:njme ) :: nfld
3338 REAL, DIMENSION ( cims:cime, cjms:cjme ), INTENT(IN) :: CBWGT1,CBWGT2,CBWGT3,CBWGT4 ! dummy
3339 REAL, DIMENSION ( nims:nime, njms:njme ), INTENT(IN) :: HBWGT1,HBWGT2,HBWGT3,HBWGT4
3340 INTEGER, DIMENSION ( cims:cime, cjms:cjme ), INTENT(IN) :: CII,CJJ ! dummy
3341 INTEGER, DIMENSION ( nims:nime, njms:njme ), INTENT(IN) :: IIH,JJH
3342 INTEGER, DIMENSION ( nims:nime, njms:njme ) :: imask
3349 REAL, DIMENSION (4, nims:nime, njms:njme ) :: NBWGT
3352 !*** INDEX CONVENTIONS
3360 !*** NBWGT1=1 NBWGT2=0
3366 DO J=NJTS,MIN(NJTE,NJDE-1)
3367 DO I=NITS,MIN(NITE,NIDE-1)
3368 IF(IMASK(I,J) .NE. 1)THEN
3369 NBWGT(1,I,J)=HBWGT1(I,J)
3370 NBWGT(2,I,J)=HBWGT2(I,J)
3371 NBWGT(3,I,J)=HBWGT3(I,J)
3372 NBWGT(4,I,J)=HBWGT4(I,J)
3377 DO J=NJTS,MIN(NJTE,NJDE-1)
3378 DO I=NITS,MIN(NITE,NIDE-1)
3379 IF(IMASK(I,J) .NE. 1)THEN
3383 AMAXVAL=amax1(NBWGT(N,I,J),AMAXVAL)
3389 IF(AMAXVAL .EQ. NBWGT(N,I,J) .AND. FLIP)THEN
3395 SUM=SUM+NBWGT(N,I,J)
3396 IF(SUM .GT. 1.0)CALL wrf_error_fatal ( "horizontal interp error - interp_hnear_nmm" )
3403 DO J=NJTS,MIN(NJTE,NJDE-1)
3405 DO I=NITS,MIN(NITE,NIDE-1)
3406 IF(IMASK(I,J) .NE. 1)THEN
3407 IF(MOD(JJH(I,J),2) .NE. 0)THEN ! 1,3,5,7
3408 NFLD(I,K,J) = NBWGT(1,I,J)*CFLD(IIH(I,J), K,JJH(I,J) ) &
3409 + NBWGT(2,I,J)*CFLD(IIH(I,J)+1,K,JJH(I,J) ) &
3410 + NBWGT(3,I,J)*CFLD(IIH(I,J), K,JJH(I,J)-1) &
3411 + NBWGT(4,I,J)*CFLD(IIH(I,J), K,JJH(I,J)+1)
3413 NFLD(I,K,J) = NBWGT(1,I,J)*CFLD(IIH(I,J), K,JJH(I,J) ) &
3414 + NBWGT(2,I,J)*CFLD(IIH(I,J)+1,K,JJH(I,J) ) &
3415 + NBWGT(3,I,J)*CFLD(IIH(I,J)+1,K,JJH(I,J)-1) &
3416 + NBWGT(4,I,J)*CFLD(IIH(I,J)+1,K,JJH(I,J)+1)
3423 END SUBROUTINE interp_hnear_ikj_nmm
3425 !=======================================================================================
3426 ! E grid nearest neighbour interpolation for integer H points
3427 !=======================================================================================
3429 SUBROUTINE interp_int_hnear_nmm (cfld, & ! CD field; integers
3430 cids, cide, ckds, ckde, cjds, cjde, &
3431 cims, cime, ckms, ckme, cjms, cjme, &
3432 cits, cite, ckts, ckte, cjts, cjte, &
3433 nfld, & ! ND field; integers
3434 nids, nide, nkds, nkde, njds, njde, &
3435 nims, nime, nkms, nkme, njms, njme, &
3436 nits, nite, nkts, nkte, njts, njte, &
3437 shw, & ! stencil half width for interp
3438 imask, & ! interpolation mask
3439 xstag, ystag, & ! staggering of field
3440 ipos, jpos, & ! lower left of nest in CD
3441 nri, nrj, & ! nest ratios
3442 CII, IIH, CJJ, JJH, CBWGT1, HBWGT1, & ! s-w grid locs and weights
3443 CBWGT2, HBWGT2, CBWGT3, HBWGT3, & ! note that "C"ourse grid ones are
3444 CBWGT4, HBWGT4 ) ! just dummys
3448 INTEGER, INTENT(IN) :: cids, cide, ckds, ckde, cjds, cjde, &
3449 cims, cime, ckms, ckme, cjms, cjme, &
3450 cits, cite, ckts, ckte, cjts, cjte, &
3451 nids, nide, nkds, nkde, njds, njde, &
3452 nims, nime, nkms, nkme, njms, njme, &
3453 nits, nite, nkts, nkte, njts, njte, &
3457 LOGICAL, INTENT(IN) :: xstag, ystag
3459 INTEGER, DIMENSION ( cims:cime, cjms:cjme, ckms:ckme ) :: cfld
3460 INTEGER, DIMENSION ( nims:nime, njms:njme, nkms:nkme ) :: nfld
3461 REAL, DIMENSION ( cims:cime, cjms:cjme ), INTENT(IN) :: CBWGT1,CBWGT2,CBWGT3,CBWGT4 ! dummy
3462 REAL, DIMENSION ( nims:nime, njms:njme ), INTENT(IN) :: HBWGT1,HBWGT2,HBWGT3,HBWGT4
3463 INTEGER, DIMENSION ( cims:cime, cjms:cjme ), INTENT(IN) :: CII,CJJ ! dummy
3464 INTEGER, DIMENSION ( nims:nime, njms:njme ), INTENT(IN) :: IIH,JJH
3465 INTEGER, DIMENSION ( nims:nime, njms:njme ) :: imask
3472 REAL, DIMENSION (4, nims:nime, njms:njme ) :: NBWGT
3475 !*** INDEX CONVENTIONS
3483 !*** NBWGT1=1 NBWGT2=0
3489 DO J=NJTS,MIN(NJTE,NJDE-1)
3490 DO I=NITS,MIN(NITE,NIDE-1)
3491 IF(IMASK(I,J) .NE. 1)THEN
3492 NBWGT(1,I,J)=HBWGT1(I,J)
3493 NBWGT(2,I,J)=HBWGT2(I,J)
3494 NBWGT(3,I,J)=HBWGT3(I,J)
3495 NBWGT(4,I,J)=HBWGT4(I,J)
3500 DO J=NJTS,MIN(NJTE,NJDE-1)
3501 DO I=NITS,MIN(NITE,NIDE-1)
3502 IF(IMASK(I,J) .NE. 1)THEN
3506 AMAXVAL=amax1(NBWGT(N,I,J),AMAXVAL)
3512 IF(AMAXVAL .EQ. NBWGT(N,I,J) .AND. FLIP)THEN
3518 SUM=SUM+NBWGT(N,I,J)
3519 IF(SUM .GT. 1.0)CALL wrf_error_fatal ( "horizontal interp error - interp_hnear_nmm" )
3526 DO J=NJTS,MIN(NJTE,NJDE-1)
3528 DO I=NITS,MIN(NITE,NIDE-1)
3529 IF(IMASK(I,J) .NE. 1)THEN
3530 IF(MOD(JJH(I,J),2) .NE. 0)THEN ! 1,3,5,7
3531 NFLD(I,J,K) = NBWGT(1,I,J)*CFLD(IIH(I,J), JJH(I,J) ,K) &
3532 + NBWGT(2,I,J)*CFLD(IIH(I,J)+1,JJH(I,J) ,K) &
3533 + NBWGT(3,I,J)*CFLD(IIH(I,J), JJH(I,J)-1,K) &
3534 + NBWGT(4,I,J)*CFLD(IIH(I,J), JJH(I,J)+1,K)
3536 NFLD(I,J,K) = NBWGT(1,I,J)*CFLD(IIH(I,J), JJH(I,J) ,K) &
3537 + NBWGT(2,I,J)*CFLD(IIH(I,J)+1,JJH(I,J) ,K) &
3538 + NBWGT(3,I,J)*CFLD(IIH(I,J)+1,JJH(I,J)-1,K) &
3539 + NBWGT(4,I,J)*CFLD(IIH(I,J)+1,JJH(I,J)+1,K)
3546 END SUBROUTINE interp_int_hnear_nmm
3548 !--------------------------------------------------------------------------------------
3550 SUBROUTINE nmm_bdy_hinterp (cfld, & ! CD field
3551 cids, cide, ckds, ckde, cjds, cjde, &
3552 cims, cime, ckms, ckme, cjms, cjme, &
3553 cits, cite, ckts, ckte, cjts, cjte, &
3555 nids, nide, nkds, nkde, njds, njde, &
3556 nims, nime, nkms, nkme, njms, njme, &
3557 nits, nite, nkts, nkte, njts, njte, &
3558 shw, & ! stencil half width
3559 imask, & ! interpolation mask
3560 xstag, ystag, & ! staggering of field
3561 ipos, jpos, & ! Position of lower left of nest in CD
3562 nri, nrj, & ! nest ratios
3571 CTEMP_B,NTEMP_B, & ! These temp arrays should be removed
3572 CTEMP_BT,NTEMP_BT, & ! later on
3573 CII, IIH, CJJ, JJH, CBWGT1, HBWGT1, & ! south-western grid locs and weights
3574 CBWGT2, HBWGT2, CBWGT3, HBWGT3, & ! note that "C"ourse grid ones are
3575 CBWGT4, HBWGT4 ) ! dummys
3577 ! use module_state_description
3578 USE module_configure
3579 USE module_wrf_error
3584 INTEGER, INTENT(IN) :: cids, cide, ckds, ckde, cjds, cjde, &
3585 cims, cime, ckms, ckme, cjms, cjme, &
3586 cits, cite, ckts, ckte, cjts, cjte, &
3587 nids, nide, nkds, nkde, njds, njde, &
3588 nims, nime, nkms, nkme, njms, njme, &
3589 nits, nite, nkts, nkte, njts, njte, &
3594 LOGICAL, INTENT(IN) :: xstag, ystag
3596 REAL, DIMENSION ( cims:cime, cjms:cjme, ckms:ckme ) :: cfld
3597 REAL, DIMENSION ( nims:nime, njms:njme, nkms:nkme ) :: nfld
3599 REAL, DIMENSION ( cims:cime, cjms:cjme, ckms:ckme ) :: ctemp_b,ctemp_bt
3600 REAL, DIMENSION ( nims:nime, njms:njme, nkms:nkme ) :: ntemp_b,ntemp_bt
3602 INTEGER, DIMENSION ( nims:nime, njms:njme ) :: imask
3603 REAL, DIMENSION( * ), INTENT(INOUT) :: c_bxs,n_bxs,c_bxe,n_bxe,c_bys,n_bys,c_bye,n_bye
3604 REAL, DIMENSION( * ), INTENT(INOUT) :: c_btxs,n_btxs,c_btxe,n_btxe,c_btys,n_btys,c_btye,n_btye
3605 REAL, DIMENSION ( cims:cime, cjms:cjme ), INTENT(IN) :: CBWGT1,CBWGT2,CBWGT3,CBWGT4 ! dummy
3606 REAL, DIMENSION ( nims:nime, njms:njme ), INTENT(IN) :: HBWGT1,HBWGT2,HBWGT3,HBWGT4
3607 INTEGER, DIMENSION ( cims:cime, cjms:cjme ), INTENT(IN) :: CII,CJJ ! dummy
3608 INTEGER, DIMENSION ( nims:nime, njms:njme ), INTENT(IN) :: IIH,JJH
3613 REAL, DIMENSION ( nims:nime, njms:njme, nkms:nkme ) :: cwk1,cwk2,cwk3,cwk4
3617 NMM_XS: IF(NITS .EQ. NIDS)THEN
3618 ! WRITE(0,*)'ENTERING X1 START BOUNDARY AT MASS POINTS',NJTS,MIN(NJTE,NJDE-1)
3621 DO J = NJTS,MIN(NJTE,NJDE-1)
3622 IF(MOD(J,2) .NE.0)THEN ! 1,3,5,7 of nested domain
3623 IF(MOD(JJH(I,J),2) .NE. 0)THEN ! 1,3,5,7 of the parent domain
3624 CWK1(I,J,K) = HBWGT1(I,J)*CFLD(IIH(I,J), JJH(I,J) ,K) &
3625 + HBWGT2(I,J)*CFLD(IIH(I,J)+1,JJH(I,J) ,K) &
3626 + HBWGT3(I,J)*CFLD(IIH(I,J), JJH(I,J)-1,K) &
3627 + HBWGT4(I,J)*CFLD(IIH(I,J), JJH(I,J)+1,K)
3631 CWK1(I,J,K) = HBWGT1(I,J)*CFLD(IIH(I,J), JJH(I,J) ,K) &
3632 + HBWGT2(I,J)*CFLD(IIH(I,J)+1,JJH(I,J) ,K) &
3633 + HBWGT3(I,J)*CFLD(IIH(I,J)+1,JJH(I,J)-1,K) &
3634 + HBWGT4(I,J)*CFLD(IIH(I,J)+1,JJH(I,J)+1,K)
3637 CWK1(I,J,K) = 0.0 ! even rows at mass points of the nested domain
3639 ntemp_b(i,J,K) = CWK1(I,J,K)
3640 ntemp_bt(i,J,K) = 0.0
3647 NMM_XE: IF(NITE-1 .EQ. NIDE-1)THEN
3648 ! WRITE(0,*)'ENTERING X END BOUNDARY AT MASS POINTS',NJTS,MIN(NJTE,NJDE-1)
3651 DO J = NJTS,MIN(NJTE,NJDE-1)
3652 IF(MOD(J,2) .NE.0)THEN ! 1,3,5,7 of the nested domain
3653 IF(MOD(JJH(I,J),2) .NE. 0)THEN ! 1,3,5,7 of the parent domain
3654 CWK2(I,J,K) = HBWGT1(I,J)*CFLD(IIH(I,J), JJH(I,J) ,K) &
3655 + HBWGT2(I,J)*CFLD(IIH(I,J)+1,JJH(I,J) ,K) &
3656 + HBWGT3(I,J)*CFLD(IIH(I,J), JJH(I,J)-1,K) &
3657 + HBWGT4(I,J)*CFLD(IIH(I,J), JJH(I,J)+1,K)
3659 CWK2(I,J,K) = HBWGT1(I,J)*CFLD(IIH(I,J), JJH(I,J) ,K) &
3660 + HBWGT2(I,J)*CFLD(IIH(I,J)+1,JJH(I,J) ,K) &
3661 + HBWGT3(I,J)*CFLD(IIH(I,J)+1,JJH(I,J)-1,K) &
3662 + HBWGT4(I,J)*CFLD(IIH(I,J)+1,JJH(I,J)+1,K)
3665 CWK2(I,J,K) = 0.0 ! even rows at mass points
3667 ntemp_b(i,J,K) = CWK2(I,J,K)
3668 ntemp_bt(i,J,K) = 0.0
3675 NMM_YS: IF(NJTS .EQ. NJDS)THEN
3676 ! WRITE(0,*)'ENTERING Y START BOUNDARY AT MASS POINTS',NITS,MIN(NITE,NIDE-1)
3679 DO I = NITS,MIN(NITE,NIDE-1)
3680 IF(MOD(JJH(I,J),2) .NE. 0)THEN ! 1,3,5,7
3681 CWK3(I,J,K) = HBWGT1(I,J)*CFLD(IIH(I,J), JJH(I,J) ,K) &
3682 + HBWGT2(I,J)*CFLD(IIH(I,J)+1,JJH(I,J) ,K) &
3683 + HBWGT3(I,J)*CFLD(IIH(I,J), JJH(I,J)-1,K) &
3684 + HBWGT4(I,J)*CFLD(IIH(I,J), JJH(I,J)+1,K)
3686 CWK3(I,J,K) = HBWGT1(I,J)*CFLD(IIH(I,J), JJH(I,J) ,K) &
3687 + HBWGT2(I,J)*CFLD(IIH(I,J)+1,JJH(I,J) ,K) &
3688 + HBWGT3(I,J)*CFLD(IIH(I,J)+1,JJH(I,J)-1,K) &
3689 + HBWGT4(I,J)*CFLD(IIH(I,J)+1,JJH(I,J)+1,K)
3691 ntemp_b(i,J,K) = CWK3(I,J,K)
3692 ntemp_bt(i,J,K) = 0.0
3699 NMM_YE: IF(NJTE-1 .EQ. NJDE-1)THEN
3700 ! WRITE(0,*)'ENTERING Y END BOUNDARY AT MASS POINTS',NITS,MIN(NITE,NIDE-1)
3703 DO I = NITS,MIN(NITE,NIDE-1)
3704 IF(MOD(JJH(I,J),2) .NE. 0)THEN ! 1,3,5,7
3705 CWK4(I,J,K) = HBWGT1(I,J)*CFLD(IIH(I,J), JJH(I,J) ,K) &
3706 + HBWGT2(I,J)*CFLD(IIH(I,J)+1,JJH(I,J) ,K) &
3707 + HBWGT3(I,J)*CFLD(IIH(I,J), JJH(I,J)-1,K) &
3708 + HBWGT4(I,J)*CFLD(IIH(I,J), JJH(I,J)+1,K)
3710 CWK4(I,J,K) = HBWGT1(I,J)*CFLD(IIH(I,J), JJH(I,J) ,K) &
3711 + HBWGT2(I,J)*CFLD(IIH(I,J)+1,JJH(I,J) ,K) &
3712 + HBWGT3(I,J)*CFLD(IIH(I,J)+1,JJH(I,J)-1,K) &
3713 + HBWGT4(I,J)*CFLD(IIH(I,J)+1,JJH(I,J)+1,K)
3716 ntemp_b(i,J,K) = CWK4(I,J,K)
3717 ntemp_bt(i,J,K) = 0.0
3724 END SUBROUTINE nmm_bdy_hinterp
3726 !--------------------------------------------------------------------------------------
3728 SUBROUTINE nmm_bdy_vinterp ( cfld, & ! CD field
3729 cids, cide, ckds, ckde, cjds, cjde, &
3730 cims, cime, ckms, ckme, cjms, cjme, &
3731 cits, cite, ckts, ckte, cjts, cjte, &
3733 nids, nide, nkds, nkde, njds, njde, &
3734 nims, nime, nkms, nkme, njms, njme, &
3735 nits, nite, nkts, nkte, njts, njte, &
3736 shw, & ! stencil half width
3737 imask, & ! interpolation mask
3738 xstag, ystag, & ! staggering of field
3739 ipos, jpos, & ! Position of lower left of nest in CD
3740 nri, nrj, & ! nest ratios
3749 CTEMP_B,NTEMP_B, & ! These temp arrays should be removed
3750 CTEMP_BT,NTEMP_BT, & ! later on
3751 CII, IIV, CJJ, JJV, CBWGT1, VBWGT1, & ! south-western grid locs and weights
3752 CBWGT2, VBWGT2, CBWGT3, VBWGT3, & ! note that "C"ourse grid ones are
3753 CBWGT4, VBWGT4 ) ! dummys
3755 ! use module_state_description
3756 USE module_configure
3757 USE module_wrf_error
3762 INTEGER, INTENT(IN) :: cids, cide, ckds, ckde, cjds, cjde, &
3763 cims, cime, ckms, ckme, cjms, cjme, &
3764 cits, cite, ckts, ckte, cjts, cjte, &
3765 nids, nide, nkds, nkde, njds, njde, &
3766 nims, nime, nkms, nkme, njms, njme, &
3767 nits, nite, nkts, nkte, njts, njte, &
3772 LOGICAL, INTENT(IN) :: xstag, ystag
3774 REAL, DIMENSION ( cims:cime, cjms:cjme, ckms:ckme ) :: cfld
3775 REAL, DIMENSION ( nims:nime, njms:njme, nkms:nkme ) :: nfld
3777 REAL, DIMENSION ( cims:cime, cjms:cjme, ckms:ckme ) :: ctemp_b,ctemp_bt
3778 REAL, DIMENSION ( nims:nime, njms:njme, nkms:nkme ) :: ntemp_b,ntemp_bt
3780 INTEGER, DIMENSION ( nims:nime, njms:njme ) :: imask
3781 REAL, DIMENSION( * ), INTENT(INOUT) :: c_bxs,n_bxs,c_bxe,n_bxe,c_bys,n_bys,c_bye,n_bye
3782 REAL, DIMENSION( * ), INTENT(INOUT) :: c_btxs,n_btxs,c_btxe,n_btxe,c_btys,n_btys,c_btye,n_btye
3783 REAL, DIMENSION ( cims:cime, cjms:cjme ), INTENT(IN) :: CBWGT1,CBWGT2,CBWGT3,CBWGT4 ! dummy
3784 REAL, DIMENSION ( nims:nime, njms:njme ), INTENT(IN) :: VBWGT1,VBWGT2,VBWGT3,VBWGT4
3785 INTEGER, DIMENSION ( cims:cime, cjms:cjme ), INTENT(IN) :: CII,CJJ ! dummy
3786 INTEGER, DIMENSION ( nims:nime, njms:njme ), INTENT(IN) :: IIV,JJV
3791 REAL, DIMENSION ( nims:nime, njms:njme, nkms:nkme ) :: cwk1,cwk2,cwk3,cwk4
3795 NMM_XS: IF(NITS .EQ. NIDS)THEN
3796 ! WRITE(0,*)'ENTERING X START BOUNDARY AT VELOCITY POINTS',NITS,NIDS,NJTS,MIN(NJTE,NJDE-1)
3799 DO J = NJTS,MIN(NJTE,NJDE-1)
3800 IF(MOD(J,2) .EQ.0)THEN ! 1,3,5,7 of nested domain
3801 IF(MOD(JJV(I,J),2) .NE. 0)THEN ! 1,3,5,7 of the parent domain
3802 CWK1(I,J,K) = VBWGT1(I,J)*CFLD(IIV(I,J), JJV(I,J) ,K) &
3803 + VBWGT2(I,J)*CFLD(IIV(I,J)+1,JJV(I,J) ,K) &
3804 + VBWGT3(I,J)*CFLD(IIV(I,J)+1,JJV(I,J)-1,K) &
3805 + VBWGT4(I,J)*CFLD(IIV(I,J)+1,JJV(I,J)+1,K)
3807 CWK1(I,J,K) = VBWGT1(I,J)*CFLD(IIV(I,J), JJV(I,J) ,K) &
3808 + VBWGT2(I,J)*CFLD(IIV(I,J)+1,JJV(I,J) ,K) &
3809 + VBWGT3(I,J)*CFLD(IIV(I,J), JJV(I,J)-1,K) &
3810 + VBWGT4(I,J)*CFLD(IIV(I,J), JJV(I,J)+1,K)
3813 CWK1(I,J,K) = 0.0 ! odd rows along J, at mass points have zero velocity
3815 ntemp_b(i,J,K) = CWK1(I,J,K)
3816 ntemp_bt(i,J,K) = 0.0
3823 NMM_XE: IF(NITE-1 .EQ. NIDE-1)THEN
3824 ! WRITE(0,*)'ENTERING X END BOUNDARY AT VELOCITY POINTS',NITE-1,NIDE-1,NJTS,MIN(NJTE,NJDE-1)
3827 DO J = NJTS,MIN(NJTE,NJDE-1)
3828 IF(MOD(J,2) .EQ.0)THEN ! 1,3,5,7 of the nested domain
3829 IF(MOD(JJV(I,J),2) .NE. 0)THEN ! 1,3,5,7 of the parent domain
3830 CWK2(I,J,K) = VBWGT1(I,J)*CFLD(IIV(I,J), JJV(I,J) ,K) &
3831 + VBWGT2(I,J)*CFLD(IIV(I,J)+1,JJV(I,J) ,K) &
3832 + VBWGT3(I,J)*CFLD(IIV(I,J)+1,JJV(I,J)-1,K) &
3833 + VBWGT4(I,J)*CFLD(IIV(I,J)+1,JJV(I,J)+1,K)
3835 CWK2(I,J,K) = VBWGT1(I,J)*CFLD(IIV(I,J), JJV(I,J) ,K) &
3836 + VBWGT2(I,J)*CFLD(IIV(I,J)+1,JJV(I,J) ,K) &
3837 + VBWGT3(I,J)*CFLD(IIV(I,J), JJV(I,J)-1,K) &
3838 + VBWGT4(I,J)*CFLD(IIV(I,J), JJV(I,J)+1,K)
3841 CWK2(I,J,K) = 0.0 ! odd rows at mass points
3843 ntemp_b(i,J,K) = CWK2(I,J,K)
3844 ntemp_bt(i,J,K) = 0.0
3851 NMM_YS: IF(NJTS .EQ. NJDS)THEN
3852 ! WRITE(0,*)'ENTERING Y START BOUNDARY AT VELOCITY POINTS',NJTS,NJDS,NITS,MIN(NITE,NIDE-1)
3855 DO I = NITS,MIN(NITE,NIDE-2) ! NIDE-1 SHOULD NOT MATTER IF WE FILL UP PHANTOM CELL
3856 IF(MOD(JJV(I,J),2) .NE. 0)THEN ! 1,3,5,7
3857 CWK3(I,J,K) = VBWGT1(I,J)*CFLD(IIV(I,J), JJV(I,J) ,K) &
3858 + VBWGT2(I,J)*CFLD(IIV(I,J)+1,JJV(I,J) ,K) &
3859 + VBWGT3(I,J)*CFLD(IIV(I,J)+1,JJV(I,J)-1,K) &
3860 + VBWGT4(I,J)*CFLD(IIV(I,J)+1,JJV(I,J)+1,K)
3862 CWK3(I,J,K) = VBWGT1(I,J)*CFLD(IIV(I,J), JJV(I,J) ,K) &
3863 + VBWGT2(I,J)*CFLD(IIV(I,J)+1,JJV(I,J) ,K) &
3864 + VBWGT3(I,J)*CFLD(IIV(I,J), JJV(I,J)-1,K) &
3865 + VBWGT4(I,J)*CFLD(IIV(I,J), JJV(I,J)+1,K)
3867 ntemp_b(i,J,K) = CWK3(I,J,K)
3868 ntemp_bt(i,J,K) = 0.0
3875 NMM_YE: IF(NJTE-1 .EQ. NJDE-1)THEN
3876 ! WRITE(0,*)'ENTERING Y END BOUNDARY AT VELOCITY POINTS',NJTE-1,NJDE-1,NITS,MIN(NITE,NIDE-1)
3879 DO I = NITS,MIN(NITE,NIDE-2) ! NIDE-1 SHOULD NOT MATTER IF WE FILL UP PHANTOM CELL
3880 IF(MOD(JJV(I,J),2) .NE. 0)THEN ! 1,3,5,7
3881 CWK4(I,J,K) = VBWGT1(I,J)*CFLD(IIV(I,J), JJV(I,J) ,K) &
3882 + VBWGT2(I,J)*CFLD(IIV(I,J)+1,JJV(I,J) ,K) &
3883 + VBWGT3(I,J)*CFLD(IIV(I,J)+1,JJV(I,J)-1,K) &
3884 + VBWGT4(I,J)*CFLD(IIV(I,J)+1,JJV(I,J)+1,K)
3886 CWK4(I,J,K) = VBWGT1(I,J)*CFLD(IIV(I,J), JJV(I,J) ,K) &
3887 + VBWGT2(I,J)*CFLD(IIV(I,J)+1,JJV(I,J) ,K) &
3888 + VBWGT3(I,J)*CFLD(IIV(I,J), JJV(I,J)-1,K) &
3889 + VBWGT4(I,J)*CFLD(IIV(I,J), JJV(I,J)+1,K)
3891 ntemp_b(i,J,K) = CWK4(I,J,K)
3892 ntemp_bt(i,J,K) = 0.0
3899 END SUBROUTINE nmm_bdy_vinterp
3902 !=======================================================================================
3903 ! E grid interpolation: simple copy from parent to mother domain
3904 !=======================================================================================
3907 SUBROUTINE nmm_copy ( cfld, & ! CD field
3908 cids, cide, ckds, ckde, cjds, cjde, &
3909 cims, cime, ckms, ckme, cjms, cjme, &
3910 cits, cite, ckts, ckte, cjts, cjte, &
3912 nids, nide, nkds, nkde, njds, njde, &
3913 nims, nime, nkms, nkme, njms, njme, &
3914 nits, nite, nkts, nkte, njts, njte, &
3915 shw, & ! stencil half width
3916 imask, & ! interpolation mask
3917 xstag, ystag, & ! staggering of field
3918 ipos, jpos, & ! Position of lower left of nest in CD
3919 nri, nrj, & ! nest ratios
3920 CII, IIH, CJJ, JJH )
3925 LOGICAL, INTENT(IN) :: xstag, ystag
3926 INTEGER, INTENT(IN) :: cids, cide, ckds, ckde, cjds, cjde, &
3927 cims, cime, ckms, ckme, cjms, cjme, &
3928 cits, cite, ckts, ckte, cjts, cjte, &
3929 nids, nide, nkds, nkde, njds, njde, &
3930 nims, nime, nkms, nkme, njms, njme, &
3931 nits, nite, nkts, nkte, njts, njte, &
3935 REAL, DIMENSION ( cims:cime, cjms:cjme, ckms:ckme ), INTENT(IN) :: cfld
3936 REAL, DIMENSION ( nims:nime, njms:njme, nkms:nkme ), INTENT(INOUT) :: nfld
3937 INTEGER, DIMENSION ( nims:nime, njms:njme ), INTENT(IN) :: imask
3938 INTEGER, DIMENSION ( cims:cime, cjms:cjme ), INTENT(IN) :: CII,CJJ ! dummy
3939 INTEGER, DIMENSION ( nims:nime, njms:njme ), INTENT(IN) :: IIH,JJH
3944 DO J=NJTS,MIN(NJTE,NJDE-1)
3946 DO I=NITS,MIN(NITE,NIDE-1)
3947 NFLD(I,J,K) = CFLD(IIH(I,J),JJH(I,J),K)
3954 END SUBROUTINE nmm_copy
3956 !=======================================================================================
3957 ! E grid test for mass point coincidence
3958 !=======================================================================================
3960 SUBROUTINE test_nmm (cfld, & ! CD field
3961 cids, cide, ckds, ckde, cjds, cjde, &
3962 cims, cime, ckms, ckme, cjms, cjme, &
3963 cits, cite, ckts, ckte, cjts, cjte, &
3965 nids, nide, nkds, nkde, njds, njde, &
3966 nims, nime, nkms, nkme, njms, njme, &
3967 nits, nite, nkts, nkte, njts, njte, &
3968 shw, & ! stencil half width for interp
3969 imask, & ! interpolation mask
3970 xstag, ystag, & ! staggering of field
3971 ipos, jpos, & ! Position of lower left of nest in CD
3972 nri, nrj, & ! nest ratios
3973 CII, IIH, CJJ, JJH, CBWGT1, HBWGT1, & ! south-western grid locs and weights
3974 CBWGT2, HBWGT2, CBWGT3, HBWGT3, & ! note that "C"ourse grid ones are
3975 CBWGT4, HBWGT4 ) ! dummys for weights
3979 INTEGER, INTENT(IN) :: cids, cide, ckds, ckde, cjds, cjde, &
3980 cims, cime, ckms, ckme, cjms, cjme, &
3981 cits, cite, ckts, ckte, cjts, cjte, &
3982 nids, nide, nkds, nkde, njds, njde, &
3983 nims, nime, nkms, nkme, njms, njme, &
3984 nits, nite, nkts, nkte, njts, njte, &
3988 LOGICAL, INTENT(IN) :: xstag, ystag
3990 REAL, DIMENSION ( cims:cime, cjms:cjme, ckms:ckme ) :: cfld
3991 REAL, DIMENSION ( nims:nime, njms:njme, nkms:nkme ) :: nfld
3992 REAL, DIMENSION ( cims:cime, cjms:cjme ), INTENT(IN) :: CBWGT1,CBWGT2,CBWGT3,CBWGT4 ! dummy
3993 REAL, DIMENSION ( nims:nime, njms:njme ), INTENT(IN) :: HBWGT1,HBWGT2,HBWGT3,HBWGT4
3994 INTEGER, DIMENSION ( cims:cime, cjms:cjme ), INTENT(IN) :: CII,CJJ ! dummy
3995 INTEGER, DIMENSION ( nims:nime, njms:njme ), INTENT(IN) :: IIH,JJH
3996 INTEGER, DIMENSION ( nims:nime, njms:njme ) :: imask
4000 REAL,PARAMETER :: error=0.0001,error1=1.0
4003 !*** CHECK DOMAIN BOUNDS BEFORE INTERPOLATION
4005 DO J=NJTS,MIN(NJTE,NJDE-1)
4006 DO I=NITS,MIN(NITE,NIDE-1)
4007 IF(IIH(i,j).LT.(CIDS-shw) .OR. IIH(i,j).GT.(CIDE+shw)) &
4008 CALL wrf_error_fatal ('hpoints:check domain bounds along x' )
4009 IF(JJH(i,j).LT.(CJDS-shw) .OR. JJH(i,j).GT.(CJDE+shw)) &
4010 CALL wrf_error_fatal ('hpoints:check domain bounds along y' )
4015 !*** INDEX CONVENTIONS
4030 ! WRITE(0,*)NITS,MIN(NITE,NIDE-1),CITS,CITE
4031 DO J=NJTS,MIN(NJTE,NJDE-1)
4033 DO I=NITS,MIN(NITE,NIDE-1)
4034 IF(ABS(1.0-HBWGT1(I,J)) .LE. ERROR)THEN
4035 DIFF=ABS(NFLD(I,J,K)-CFLD(IIH(I,J),JJH(I,J),K))
4036 IF(DIFF .GT. ERROR)THEN
4037 CALL wrf_debug(1,"dyn_nmm: NON-COINCIDENT, NESTED MASS POINT")
4038 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
4040 IF(DIFF .GT. ERROR1)THEN
4041 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
4042 CALL wrf_error_fatal ('dyn_nmm: NON-COINCIDENT, NESTED MASS POINT')
4049 END SUBROUTINE test_nmm
4051 !==================================
4052 ! this is the default function used in nmm feedback at mass points.
4054 SUBROUTINE nmm_feedback ( cfld, & ! CD field
4055 cids, cide, ckds, ckde, cjds, cjde, &
4056 cims, cime, ckms, ckme, cjms, cjme, &
4057 cits, cite, ckts, ckte, cjts, cjte, &
4059 nids, nide, nkds, nkde, njds, njde, &
4060 nims, nime, nkms, nkme, njms, njme, &
4061 nits, nite, nkts, nkte, njts, njte, &
4062 shw, & ! stencil half width for interp
4063 imask, & ! interpolation mask
4064 xstag, ystag, & ! staggering of field
4065 ipos, jpos, & ! Position of lower left of nest in CD
4066 nri, nrj, & ! nest ratios
4067 CII, IIH, CJJ, JJH, &
4068 CBWGT1, HBWGT1, CBWGT2, HBWGT2, &
4069 CBWGT3, HBWGT3, CBWGT4, HBWGT4 )
4070 USE module_configure
4074 INTEGER, INTENT(IN) :: cids, cide, ckds, ckde, cjds, cjde, &
4075 cims, cime, ckms, ckme, cjms, cjme, &
4076 cits, cite, ckts, ckte, cjts, cjte, &
4077 nids, nide, nkds, nkde, njds, njde, &
4078 nims, nime, nkms, nkme, njms, njme, &
4079 nits, nite, nkts, nkte, njts, njte, &
4083 INTEGER,DIMENSION(cims:cime,cjms:cjme), INTENT(IN) :: CII,CJJ ! dummy
4084 INTEGER,DIMENSION(nims:nime,njms:njme), INTENT(IN) :: IIH,JJH
4085 REAL,DIMENSION(cims:cime,cjms:cjme), INTENT(IN) :: CBWGT1,CBWGT2,CBWGT3,CBWGT4
4086 REAL,DIMENSION(nims:nime,njms:njme), INTENT(IN) :: HBWGT1,HBWGT2,HBWGT3,HBWGT4
4087 LOGICAL, INTENT(IN) :: xstag, ystag
4089 REAL, DIMENSION ( cims:cime, cjms:cjme, ckms:ckme ), INTENT(OUT) :: cfld
4090 REAL, DIMENSION ( nims:nime, njms:njme, nkms:nkme ), INTENT(IN) :: nfld
4091 INTEGER, DIMENSION ( nims:nime, njms:njme ),INTENT(IN) :: imask
4095 INTEGER ci, cj, ck, ni, nj, nk, ip, jp, ioff, joff, ioffa, joffa
4096 INTEGER :: icmin,icmax,jcmin,jcmax
4097 INTEGER :: is, ipoints,jpoints,ijpoints
4098 INTEGER , PARAMETER :: passes = 2
4101 !=====================================================================================
4104 IF(nri .ne. 3 .OR. nrj .ne. 3) &
4105 CALL wrf_error_fatal ('Feedback works for only 1:3 ratios, currently. Modify the namelist' )
4107 ! WRITE(0,*)'SIMPLE FEED BACK IS SWITCHED ON FOR MASS'
4113 DO cj = MAX(jpos+1,cjts),MIN(jpos+(njde-njds)/nrj-1,cjte) ! exclude top and bottom BCs
4114 nj = (cj-jpos)*nrj + 1
4115 if(mod(cj,2) .eq. 0)THEN
4116 is=0 ! even rows for mass points (2,4,6,8)
4118 is=1 ! odd rows for mass points (1,3,5,7)
4120 DO ci = MAX(ipos+is,cits),MIN(ipos+(nide-nids)/nri-1,cite) ! excludes LBCs
4121 ni = (ci-ipos)*nri + 2 -is
4122 IF(IS==0)THEN ! (2,4,6,8)
4123 ! AVGH = NFLD(NI,NJ+1,NK) + NFLD(NI,NJ-1,NK) + NFLD(NI+1,NJ+1,NK)+ NFLD(NI+1,NJ-1,NK) &
4124 ! + NFLD(NI+1,NJ,NK) + NFLD(NI-1,NJ,NK) + NFLD(NI,NJ+2,NK) + NFLD(NI,NJ-2,NK) &
4125 ! + NFLD(NI+1,NJ+2,NK)+ NFLD(NI-1,NJ+2,NK)+ NFLD(NI+1,NJ-2,NK)+ NFLD(NI-1,NJ-2,NK)
4127 AVGH = NFLD(NI,NJ+2,NK) &
4128 + NFLD(NI ,NJ+1,NK) + NFLD(NI+1,NJ+1,NK) &
4129 + NFLD(NI-1,NJ ,NK) + NFLD(NI,NJ ,NK) + NFLD(NI+1,NJ ,NK) &
4130 + NFLD(NI ,NJ-1,NK) + NFLD(NI+1,NJ-1,NK) &
4134 ! AVGH = NFLD(NI,NJ+1,NK) + NFLD(NI,NJ-1,NK) + NFLD(NI-1,NJ+1,NK)+ NFLD(NI-1,NJ-1,NK) &
4135 ! + NFLD(NI+1,NJ,NK) + NFLD(NI-1,NJ,NK) + NFLD(NI,NJ+2,NK) + NFLD(NI,NJ-2,NK) &
4136 ! + NFLD(NI+1,NJ+2,NK)+ NFLD(NI-1,NJ+2,NK)+ NFLD(NI+1,NJ-2,NK)+ NFLD(NI-1,NJ-2,NK)
4138 AVGH = NFLD(NI,NJ+2,NK) &
4139 + NFLD(NI-1,NJ+1,NK) + NFLD(NI,NJ+1,NK) &
4140 + NFLD(NI-1,NJ ,NK) + NFLD(NI,NJ ,NK) + NFLD(NI+1,NJ ,NK) &
4141 + NFLD(NI-1,NJ-1,NK) + NFLD(NI,NJ-1,NK) &
4145 !dusan CFLD(CI,CK,CJ) = 0.5*CFLD(CI,CK,CJ) + 0.5*(NFLD(NI,NK,NJ)+AVGH)/13.0
4146 ! CFLD(CI,CJ,CK) = (NFLD(NI,NJ,NK)+AVGH)/13.0
4147 CFLD(CI,CJ,CK) = AVGH/9.0
4152 END SUBROUTINE nmm_feedback
4154 !===========================================================================================
4156 SUBROUTINE nmm_vfeedback ( cfld, & ! CD field
4157 cids, cide, ckds, ckde, cjds, cjde, &
4158 cims, cime, ckms, ckme, cjms, cjme, &
4159 cits, cite, ckts, ckte, cjts, cjte, &
4161 nids, nide, nkds, nkde, njds, njde, &
4162 nims, nime, nkms, nkme, njms, njme, &
4163 nits, nite, nkts, nkte, njts, njte, &
4164 shw, & ! stencil half width for interp
4165 imask, & ! interpolation mask
4166 xstag, ystag, & ! staggering of field
4167 ipos, jpos, & ! Position of lower left of nest in CD
4168 nri, nrj, & ! nest ratios
4169 CII, IIV, CJJ, JJV, &
4170 CBWGT1, VBWGT1, CBWGT2, VBWGT2, &
4171 CBWGT3, VBWGT3, CBWGT4, VBWGT4 )
4172 USE module_configure
4176 INTEGER, INTENT(IN) :: cids, cide, ckds, ckde, cjds, cjde, &
4177 cims, cime, ckms, ckme, cjms, cjme, &
4178 cits, cite, ckts, ckte, cjts, cjte, &
4179 nids, nide, nkds, nkde, njds, njde, &
4180 nims, nime, nkms, nkme, njms, njme, &
4181 nits, nite, nkts, nkte, njts, njte, &
4185 INTEGER,DIMENSION(cims:cime,cjms:cjme), INTENT(IN) :: CII,CJJ ! dummy
4186 INTEGER,DIMENSION(nims:nime,njms:njme), INTENT(IN) :: IIV,JJV
4187 REAL,DIMENSION(cims:cime,cjms:cjme), INTENT(IN) :: CBWGT1,CBWGT2,CBWGT3,CBWGT4
4188 REAL,DIMENSION(nims:nime,njms:njme), INTENT(IN) :: VBWGT1,VBWGT2,VBWGT3,VBWGT4
4189 LOGICAL, INTENT(IN) :: xstag, ystag
4191 REAL, DIMENSION ( cims:cime, cjms:cjme, ckms:ckme ), INTENT(OUT) :: cfld
4192 REAL, DIMENSION ( nims:nime, njms:njme, nkms:nkme ), INTENT(IN) :: nfld
4193 INTEGER, DIMENSION ( nims:nime, njms:njme ),INTENT(IN) :: imask
4197 INTEGER ci, cj, ck, ni, nj, nk, ip, jp, ioff, joff, ioffa, joffa
4198 INTEGER :: icmin,icmax,jcmin,jcmax
4199 INTEGER :: is, ipoints,jpoints,ijpoints
4200 INTEGER , PARAMETER :: passes = 2
4203 !=====================================================================================
4206 IF(nri .ne. 3 .OR. nrj .ne. 3) &
4207 CALL wrf_error_fatal ('Feedback works for only 1:3 ratios, currently. Modify the namelist')
4209 ! WRITE(0,*)'SIMPLE FEED BACK IS SWITCHED ON FOR VELOCITY'
4215 DO cj = MAX(jpos+1,cjts),MIN(jpos+(njde-njds)/nrj-1,cjte) ! exclude top and bottom BCs
4216 nj = (cj-jpos)*nrj + 1
4217 if(mod(cj,2) .eq. 0)THEN
4218 is=1 ! even rows for velocity points (2,4,6,8)
4220 is=0 ! odd rows for velocity points (1,3,5,7)
4222 DO ci = MAX(ipos+is,cits),MIN(ipos+(nide-nids)/nri-1,cite) ! excludes LBCs
4223 ni = (ci-ipos)*nri + 2 -is
4224 IF(IS==0)THEN ! (1,3,5,7)
4225 ! AVGV = NFLD(NI,NK,NJ+1) + NFLD(NI,NK,NJ-1) + NFLD(NI+1,NK,NJ+1)+ NFLD(NI+1,NK,NJ-1) &
4226 ! + NFLD(NI+1,NK,NJ) + NFLD(NI-1,NK,NJ) + NFLD(NI,NK,NJ+2) + NFLD(NI,NK,NJ-2) &
4227 ! + NFLD(NI+1,NK,NJ+2)+ NFLD(NI-1,NK,NJ+2)+ NFLD(NI+1,NK,NJ-2)+ NFLD(NI-1,NK,NJ-2)
4229 AVGV = NFLD(NI,NJ+2,NK) &
4230 + NFLD(NI ,NJ+1,NK) + NFLD(NI+1,NJ+1,NK) &
4231 + NFLD(NI-1,NJ ,NK) + NFLD(NI,NJ ,NK) + NFLD(NI+1,NJ ,NK) &
4232 + NFLD(NI ,NJ-1,NK) + NFLD(NI+1,NJ-1,NK) &
4236 ! AVGV = NFLD(NI,NK,NJ+1) + NFLD(NI,NK,NJ-1) + NFLD(NI-1,NK,NJ+1)+ NFLD(NI-1,NK,NJ-1) &
4237 ! + NFLD(NI+1,NK,NJ) + NFLD(NI-1,NK,NJ) + NFLD(NI,NK,NJ+2) + NFLD(NI,NK,NJ-2) &
4238 ! + NFLD(NI+1,NK,NJ+2)+ NFLD(NI-1,NK,NJ+2)+ NFLD(NI+1,NK,NJ-2)+ NFLD(NI-1,NK,NJ-2)
4240 AVGV = NFLD(NI,NJ+2,NK) &
4241 + NFLD(NI-1,NJ+1,NK) + NFLD(NI,NJ+1,NK) &
4242 + NFLD(NI-1,NJ ,NK) + NFLD(NI,NJ ,NK) + NFLD(NI+1,NJ ,NK) &
4243 + NFLD(NI-1,NJ-1,NK) + NFLD(NI,NJ-1,NK) &
4247 !dusan CFLD(CI,CK,CJ) = 0.5*CFLD(CI,CK,CJ) + 0.5*(NFLD(NI,NK,NJ)+AVGV)/13.0
4248 ! CFLD(CI,CK,CJ) = (NFLD(NI,NK,NJ)+AVGV)/13.0
4249 CFLD(CI,CJ,CK) = AVGV/9.0
4254 END SUBROUTINE nmm_vfeedback
4257 SUBROUTINE nmm_smoother ( cfld , &
4258 cids, cide, ckds, ckde, cjds, cjde, &
4259 cims, cime, ckms, ckme, cjms, cjme, &
4260 cits, cite, ckts, ckte, cjts, cjte, &
4261 nids, nide, nkds, nkde, njds, njde, &
4262 nims, nime, nkms, nkme, njms, njme, &
4263 nits, nite, nkts, nkte, njts, njte, &
4269 USE module_configure
4272 INTEGER, INTENT(IN) :: cids, cide, ckds, ckde, cjds, cjde, &
4273 cims, cime, ckms, ckme, cjms, cjme, &
4274 cits, cite, ckts, ckte, cjts, cjte, &
4275 nids, nide, nkds, nkde, njds, njde, &
4276 nims, nime, nkms, nkme, njms, njme, &
4277 nits, nite, nkts, nkte, njts, njte, &
4280 REAL, DIMENSION ( cims:cime, ckms:ckme, cjms:cjme ), INTENT(INOUT) :: cfld
4281 LOGICAL, INTENT(IN) :: xstag, ystag
4287 INTEGER, PARAMETER :: smooth_passes = 5
4289 REAL, DIMENSION ( cims:cime, ckms:ckme, cjms:cjme ) :: cfldnew
4290 INTEGER :: ci, cj, ck
4291 INTEGER :: is, npass
4295 ! If there is no feedback, there can be no smoothing.
4297 CALL nl_get_feedback ( 1, feedback )
4298 IF ( feedback == 0 ) RETURN
4300 WRITE(0,*)'SIMPLE SMOOTHER IS SWITCHED ON FOR HEIGHT'
4302 DO npass = 1, smooth_passes
4304 DO cj = MAX(jpos+1,cjts),MIN(jpos+(njde-njds)/nrj-1,cjte) ! exclude top and bottom BCs
4305 if(mod(cj,2) .eq. 0)THEN
4306 is=0 ! even rows for mass points (2,4,6,8)
4308 is=1 ! odd rows for mass points (1,3,5,7)
4311 DO ci = MAX(ipos+is,cits),MIN(ipos+(nide-nids)/nri-1,cite) ! excludes LBCs
4312 IF(IS==0)THEN ! (2,4,6,8)
4313 AVGH = CFLD(CI,CK,CJ+1) + CFLD(CI,CK,CJ-1) + CFLD(CI+1,CK,CJ+1) + CFLD(CI+1,CK,CJ-1)
4315 AVGH = CFLD(CI,CK,CJ+1) + CFLD(CI,CK,CJ-1) + CFLD(CI-1,CK,CJ+1) + CFLD(CI-1,CK,CJ-1)
4317 CFLDNEW(CI,CK,CJ) = (AVGH + 4*CFLD(CI,CK,CJ)) / 8.0
4322 DO cj = MAX(jpos+1,cjts),MIN(jpos+(njde-njds)/nrj-1,cjte) ! exclude top and bottom BCs
4323 if(mod(cj,2) .eq. 0)THEN
4324 is=0 ! even rows for mass points (2,4,6,8)
4326 is=1 ! odd rows for mass points (1,3,5,7)
4329 DO ci = MAX(ipos+is,cits),MIN(ipos+(nide-nids)/nri-1,cite) ! excludes LBCs
4330 CFLD(CI,CK,CJ) = CFLDNEW(CI,CK,CJ)
4337 END SUBROUTINE nmm_smoother
4340 SUBROUTINE nmm_vsmoother ( cfld , &
4341 cids, cide, ckds, ckde, cjds, cjde, &
4342 cims, cime, ckms, ckme, cjms, cjme, &
4343 cits, cite, ckts, ckte, cjts, cjte, &
4344 nids, nide, nkds, nkde, njds, njde, &
4345 nims, nime, nkms, nkme, njms, njme, &
4346 nits, nite, nkts, nkte, njts, njte, &
4352 USE module_configure
4355 INTEGER, INTENT(IN) :: cids, cide, ckds, ckde, cjds, cjde, &
4356 cims, cime, ckms, ckme, cjms, cjme, &
4357 cits, cite, ckts, ckte, cjts, cjte, &
4358 nids, nide, nkds, nkde, njds, njde, &
4359 nims, nime, nkms, nkme, njms, njme, &
4360 nits, nite, nkts, nkte, njts, njte, &
4363 REAL, DIMENSION ( cims:cime, ckms:ckme, cjms:cjme ), INTENT(INOUT) :: cfld
4364 LOGICAL, INTENT(IN) :: xstag, ystag
4370 INTEGER, PARAMETER :: smooth_passes = 5
4372 REAL, DIMENSION ( cims:cime, ckms:ckme, cjms:cjme ) :: cfldnew
4373 INTEGER :: ci, cj, ck
4374 INTEGER :: is, npass
4378 ! If there is no feedback, there can be no smoothing.
4380 CALL nl_get_feedback ( 1, feedback )
4381 IF ( feedback == 0 ) RETURN
4383 WRITE(0,*)'SIMPLE SMOOTHER IS SWITCHED ON FOR VELOCITY'
4385 DO npass = 1, smooth_passes
4387 DO cj = MAX(jpos+1,cjts),MIN(jpos+(njde-njds)/nrj-1,cjte) ! exclude top and bottom BCs
4388 if(mod(cj,2) .eq. 0)THEN
4389 is=1 ! even rows for mass points (2,4,6,8)
4391 is=0 ! odd rows for mass points (1,3,5,7)
4394 DO ci = MAX(ipos+is,cits),MIN(ipos+(nide-nids)/nri-1,cite) ! excludes LBCs
4395 IF(IS==0)THEN ! (2,4,6,8)
4396 AVGV = CFLD(CI,CK,CJ+1) + CFLD(CI,CK,CJ-1) + CFLD(CI+1,CK,CJ+1) + CFLD(CI+1,CK,CJ-1)
4398 AVGV = CFLD(CI,CK,CJ+1) + CFLD(CI,CK,CJ-1) + CFLD(CI-1,CK,CJ+1) + CFLD(CI-1,CK,CJ-1)
4400 CFLDNEW(CI,CK,CJ) = (AVGV + 4*CFLD(CI,CK,CJ)) / 8.0
4405 DO cj = MAX(jpos+1,cjts),MIN(jpos+(njde-njds)/nrj-1,cjte) ! exclude top and bottom BCs
4406 if(mod(cj,2) .eq. 0)THEN
4407 is=1 ! even rows for mass points (2,4,6,8)
4409 is=0 ! odd rows for mass points (1,3,5,7)
4412 DO ci = MAX(ipos+is,cits),MIN(ipos+(nide-nids)/nri-1,cite) ! excludes LBCs
4413 CFLD(CI,CK,CJ) = CFLDNEW(CI,CK,CJ)
4420 END SUBROUTINE nmm_vsmoother
4421 !======================================================================================
4422 ! End of gopal's doing
4423 !======================================================================================