FV3 Bundle
boundary_tlm.F90
Go to the documentation of this file.
1 !***********************************************************************
2 !* GNU General Public License *
3 !* This file is a part of fvGFS. *
4 !* *
5 !* fvGFS is free software; you can redistribute it and/or modify it *
6 !* and are expected to follow the terms of the GNU General Public *
7 !* License as published by the Free Software Foundation; either *
8 !* version 2 of the License, or (at your option) any later version. *
9 !* *
10 !* fvGFS is distributed in the hope that it will be useful, but *
11 !* WITHOUT ANY WARRANTY; without even the implied warranty of *
12 !* MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU *
13 !* General Public License for more details. *
14 !* *
15 !* For the full text of the GNU General Public License, *
16 !* write to: Free Software Foundation, Inc., *
17 !* 675 Mass Ave, Cambridge, MA 02139, USA. *
18 !* or see: http://www.gnu.org/licenses/gpl.html *
19 !***********************************************************************
21 
22  use fv_mp_nlm_mod, only: ng, isc,jsc,iec,jec, isd,jsd,ied,jed, is,js,ie,je, is_master
23  use constants_mod, only: grav
24 
26  use mpp_domains_mod, only: center, corner, north, east
28  use mpp_mod, only: mpp_error, fatal, mpp_sum, mpp_sync, mpp_npes, mpp_broadcast, warning, mpp_pe
29 
30  use fv_mp_nlm_mod, only: mp_bcst
32  use mpp_mod, only: mpp_send, mpp_recv
34  use mpp_domains_mod, only : nest_domain_type, west, south
35  use mpp_domains_mod, only : mpp_get_c2f_index, mpp_update_nest_fine
36  use mpp_domains_mod, only : mpp_get_f2c_index, mpp_update_nest_coarse
37  !use mpp_domains_mod, only : mpp_get_domain_shift
38 
39  implicit none
40  public extrapolation_bc
44 
46 
47  interface nested_grid_bc
48  module procedure nested_grid_bc_2d
49  module procedure nested_grid_bc_mpp
50  module procedure nested_grid_bc_mpp_send
51  module procedure nested_grid_bc_2d_mpp
52  module procedure nested_grid_bc_3d
53  end interface
54 
55 
56  interface fill_nested_grid
57  module procedure fill_nested_grid_2d
58  module procedure fill_nested_grid_3d
59  end interface
60 
62  module procedure update_coarse_grid_mpp
63  module procedure update_coarse_grid_mpp_2d
64  end interface
65 
66 CONTAINS
67 !Linear extrapolation into halo region
68 !Not to be confused with extrapolated-in-time nested BCs
69  SUBROUTINE extrapolation_bc(q, istag, jstag, npx, npy, bd, pd_in, &
70 & debug_in)
71  IMPLICIT NONE
72  TYPE(fv_grid_bounds_type), INTENT(IN) :: bd
73  INTEGER, INTENT(IN) :: istag, jstag, npx, npy
74  REAL, DIMENSION(bd%isd:bd%ied+istag, bd%jsd:bd%jed+jstag), intent(&
75 & INOUT) :: q
76  LOGICAL, INTENT(IN), OPTIONAL :: pd_in, debug_in
77  INTEGER :: i, j, istart, iend, jstart, jend
78  LOGICAL :: pd, debug
79  INTEGER :: is, ie, js, je
80  INTEGER :: isd, ied, jsd, jed
81  INTRINSIC max
82  INTRINSIC min
83  INTRINSIC PRESENT
84  INTRINSIC real
85  is = bd%is
86  ie = bd%ie
87  js = bd%js
88  je = bd%je
89  isd = bd%isd
90  ied = bd%ied
91  jsd = bd%jsd
92  jed = bd%jed
93  IF (isd .LT. 1) THEN
94  istart = 1
95  ELSE
96  istart = isd
97  END IF
98  IF (ied .GT. npx - 1) THEN
99  iend = npx - 1
100  ELSE
101  iend = ied
102  END IF
103  IF (jsd .LT. 1) THEN
104  jstart = 1
105  ELSE
106  jstart = jsd
107  END IF
108  IF (jed .GT. npy - 1) THEN
109  jend = npy - 1
110  ELSE
111  jend = jed
112  END IF
113 !Positive-definite extrapolation: shift from linear extrapolation to zero-gradient when the extrapolated value turns negative.
114  IF (PRESENT(pd_in)) THEN
115  pd = pd_in
116  ELSE
117  pd = .false.
118  END IF
119  IF (PRESENT(debug_in)) THEN
120  debug = debug_in
121  ELSE
122  debug = .false.
123  END IF
124  IF (is .EQ. 1) THEN
125  IF (pd) THEN
126  DO j=jstart,jend+jstag
127  DO i=0,isd,-1
128  IF (REAL(i) .LE. 1. - q(1, j)/(q(2, j)-q(1, j)+1.e-12) .AND.&
129 & q(1, j) .LT. q(2, j)) THEN
130  q(i, j) = q(i+1, j)
131  ELSE
132  q(i, j) = REAL(2-i)*q(1, j) - REAL(1-i)*q(2, j)
133  END IF
134  END DO
135  END DO
136  ELSE
137  DO j=jstart,jend+jstag
138  DO i=0,isd,-1
139  q(i, j) = REAL(2-i)*q(1, j) - REAL(1-i)*q(2, j)
140  END DO
141  END DO
142  END IF
143  END IF
144  IF (js .EQ. 1) THEN
145  IF (pd) THEN
146  DO j=0,jsd,-1
147  DO i=istart,iend+istag
148  IF (REAL(j) .LE. 1. - q(i, 1)/(q(i, 2)-q(i, 1)+1.e-12) .AND.&
149 & q(i, 1) .LT. q(i, 2)) THEN
150  q(i, j) = q(i, j+1)
151  ELSE
152  q(i, j) = REAL(2-j)*q(i, 1) - REAL(1-j)*q(i, 2)
153  END IF
154  END DO
155  END DO
156  ELSE
157  DO j=0,jsd,-1
158  DO i=istart,iend+istag
159  q(i, j) = REAL(2-j)*q(i, 1) - REAL(1-j)*q(i, 2)
160  END DO
161  END DO
162  END IF
163  END IF
164  IF (ie .EQ. npx - 1) THEN
165  IF (pd) THEN
166  DO j=jstart,jend+jstag
167  DO i=ie+1+istag,ied+istag
168  IF (REAL(i) .GE. ie + istag + q(ie+istag, j)/(q(ie+istag-1, &
169 & j)-q(ie+istag, j)+1.e-12) .AND. q(ie+istag, j) .LT. q(ie&
170 & +istag-1, j)) THEN
171  q(i, j) = q(i-1, j)
172  ELSE
173  q(i, j) = REAL(i-(ie+istag-1))*q(ie+istag, j) + REAL(ie+&
174 & istag-i)*q(ie+istag-1, j)
175  end if
176  END DO
177  END DO
178  ELSE
179  DO j=jstart,jend+jstag
180  DO i=ie+1+istag,ied+istag
181  q(i, j) = REAL(i-(ie+istag-1))*q(ie+istag, j) + REAL(ie+&
182 & istag-i)*q(ie+istag-1, j)
183  end do
184  END DO
185  END IF
186  END IF
187  IF (je .EQ. npy - 1) THEN
188  IF (pd) THEN
189  DO j=je+1+jstag,jed+jstag
190  DO i=istart,iend+istag
191  IF (REAL(j) .GE. je + jstag + q(i, je+jstag)/(q(i, je+jstag-&
192 & 1)-q(i, je+jstag)+1.e-12) .AND. q(i, je+jstag-1) .GT. q(&
193 & i, je+jstag)) THEN
194  q(i, j) = q(i, j-1)
195  ELSE
196  q(i, j) = REAL(j-(je+jstag-1))*q(i, je+jstag) + REAL(je+&
197 & jstag-j)*q(i, je+jstag-1)
198  end if
199  END DO
200  END DO
201  ELSE
202  DO j=je+1+jstag,jed+jstag
203  DO i=istart,iend+istag
204  q(i, j) = REAL(j-(je+jstag-1))*q(i, je+jstag) + REAL(je+&
205 & jstag-j)*q(i, je+jstag-1)
206  end do
207  END DO
208  END IF
209  END IF
210 !CORNERS: Average of extrapolations
211  IF (is .EQ. 1 .AND. js .EQ. 1) THEN
212  IF (pd) THEN
213  DO j=0,jsd,-1
214  DO i=0,isd,-1
215  IF (REAL(i) .LE. 1. - q(1, j)/(q(2, j)-q(1, j)+1.e-12) .AND.&
216 & q(2, j) .GT. q(1, j)) THEN
217  q(i, j) = 0.5*q(i+1, j)
218  ELSE
219  q(i, j) = 0.5*(REAL(2-i)*q(1, j)-REAL(1-i)*q(2, j))
220  END IF
221  IF (REAL(j) .LE. 1. - q(i, 1)/(q(i, 2)-q(i, 1)+1.e-12) .AND.&
222 & q(i, 2) .GT. q(i, 1)) THEN
223  q(i, j) = q(i, j) + 0.5*q(i, j+1)
224  ELSE
225  q(i, j) = q(i, j) + 0.5*(REAL(2-j)*q(i, 1)-REAL(1-j)*q(i, &
226 & 2))
227  end if
228  END DO
229  END DO
230  ELSE
231  DO j=jsd,0
232  DO i=isd,0
233  q(i, j) = 0.5*(REAL(2-i)*q(1, j)-REAL(1-i)*q(2, j)) + 0.5*(&
234 & REAL(2-j)*q(i, 1)-REAL(1-j)*q(i, 2))
235  end do
236  END DO
237  END IF
238  END IF
239  IF (is .EQ. 1 .AND. je .EQ. npy - 1) THEN
240  IF (pd) THEN
241  DO j=je+1+jstag,jed+jstag
242  DO i=0,isd,-1
243  IF (REAL(i) .LE. 1. - q(1, j)/(q(2, j)-q(1, j)+1.e-12) .AND.&
244 & q(2, j) .GT. q(1, j)) THEN
245  q(i, j) = 0.5*q(i+1, j)
246  ELSE
247  q(i, j) = 0.5*(REAL(2-i)*q(1, j)-REAL(1-i)*q(2, j))
248  END IF
249 !'Unary plus' removed to appease IBM compiler
250 !if (real(j) >= je+jstag - q(i,je+jstag)/(q(i,je+jstag-1)-q(i,je+jstag)+1.e-12) .and. &
251  IF (REAL(j) .GE. je + jstag - q(i, je+jstag)/(q(i, je+jstag-&
252 & 1)-q(i, je+jstag)+1.e-12) .AND. q(i, je+jstag-1) .GT. q(&
253 & i, je+jstag)) THEN
254  q(i, j) = q(i, j) + 0.5*q(i, j-1)
255  ELSE
256  q(i, j) = q(i, j) + 0.5*(REAL(j-(je+jstag-1))*q(i, je+&
257 & jstag)+REAL(je+jstag-j)*q(i, je+jstag-1))
258  end if
259  END DO
260  END DO
261  ELSE
262  DO j=je+1+jstag,jed+jstag
263  DO i=isd,0
264  q(i, j) = 0.5*(REAL(2-i)*q(1, j)-REAL(1-i)*q(2, j)) + 0.5*(&
265 & REAL(j-(je+jstag-1))*q(i, je+jstag)+REAL(je+jstag-j)*q(i, &
266 & je+jstag-1))
267  end do
268  END DO
269  END IF
270  END IF
271  IF (ie .EQ. npx - 1 .AND. je .EQ. npy - 1) THEN
272  IF (pd) THEN
273  DO j=je+1+jstag,jed+jstag
274  DO i=ie+1+istag,ied+istag
275  IF (REAL(i) .GE. ie + istag + q(ie+istag, j)/(q(ie+istag-1, &
276 & j)-q(ie+istag, j)+1.e-12) .AND. q(ie+istag-1, j) .GT. q(&
277 & ie+istag, j)) THEN
278  q(i, j) = 0.5*q(i-1, j)
279  ELSE
280  q(i, j) = 0.5*(REAL(i-(ie+istag-1))*q(ie+istag, j)+REAL(ie&
281 & +istag-i)*q(ie+istag-1, j))
282  end if
283  IF (REAL(j) .GE. je + jstag + q(i, je+jstag)/(q(i, je+jstag-&
284 & 1)-q(i, je+jstag)+1.e-12) .AND. q(i, je+jstag-1) .GT. q(&
285 & i, je+jstag)) THEN
286  q(i, j) = q(i, j) + 0.5*q(i, j-1)
287  ELSE
288  q(i, j) = q(i, j) + 0.5*(REAL(j-(je+jstag-1))*q(i, je+&
289 & jstag)+REAL(je+jstag-j)*q(i, je+jstag-1))
290  end if
291  END DO
292  END DO
293  ELSE
294  DO j=je+1+jstag,jed+jstag
295  DO i=ie+1+istag,ied+istag
296  q(i, j) = 0.5*(REAL(i-(ie+istag-1))*q(ie+istag, j)+REAL(ie+&
297 & istag-i)*q(ie+istag-1, j)) + 0.5*(REAL(j-(je+jstag-1))*q(i&
298 & , je+jstag)+REAL(je+jstag-j)*q(i, je+jstag-1))
299  end do
300  END DO
301  END IF
302  END IF
303  IF (ie .EQ. npx - 1 .AND. js .EQ. 1) THEN
304  IF (pd) THEN
305  DO j=0,jsd,-1
306  DO i=ie+1+istag,ied+istag
307  IF (REAL(i) .GE. ie + istag + q(ie+istag, j)/(q(ie+istag-1, &
308 & j)-q(ie+istag, j)+1.e-12) .AND. q(ie+istag-1, j) .GT. q(&
309 & ie+istag, j)) THEN
310  q(i, j) = 0.5*q(i-1, j)
311  ELSE
312  q(i, j) = 0.5*(REAL(i-(ie+istag-1))*q(ie+istag, j)+REAL(ie&
313 & +istag-i)*q(ie+istag-1, j))
314  end if
315  IF (REAL(j) .LE. 1. - q(i, 1)/(q(i, 2)-q(i, 1)+1.e-12) .AND.&
316 & q(i, 2) .GT. q(i, 1)) THEN
317  q(i, j) = q(i, j) + 0.5*q(i, j+1)
318  ELSE
319  q(i, j) = q(i, j) + 0.5*(REAL(2-j)*q(i, 1)-REAL(1-j)*q(i, &
320 & 2))
321  end if
322  END DO
323  END DO
324  ELSE
325  DO j=jsd,0
326  DO i=ie+1+istag,ied+istag
327  q(i, j) = 0.5*(REAL(i-(ie+istag-1))*q(ie+istag, j)+REAL(ie+&
328 & istag-i)*q(ie+istag-1, j)) + 0.5*(REAL(2-j)*q(i, 1)-REAL(1&
329 & -j)*q(i, 2))
330  end do
331  END DO
332  END IF
333  END IF
334  END SUBROUTINE extrapolation_bc
335  SUBROUTINE fill_nested_grid_2d(var_nest, var_coarse, ind, wt, istag, &
336 & jstag, isg, ieg, jsg, jeg, bd, istart_in, iend_in, jstart_in, &
337 & jend_in)
338  IMPLICIT NONE
339  TYPE(FV_GRID_BOUNDS_TYPE), INTENT(IN) :: bd
340  INTEGER, INTENT(IN) :: istag, jstag, isg, ieg, jsg, jeg
341  REAL, DIMENSION(bd%isd:bd%ied+istag, bd%jsd:bd%jed+jstag), INTENT(&
342 & INOUT) :: var_nest
343  REAL, DIMENSION(isg:ieg+istag, jsg:jeg+jstag), INTENT(IN) :: &
344 & var_coarse
345  INTEGER, DIMENSION(bd%isd:bd%ied+istag, bd%jsd:bd%jed+jstag, 2), &
346 & INTENT(IN) :: ind
347  REAL, DIMENSION(bd%isd:bd%ied+istag, bd%jsd:bd%jed+jstag, 4), &
348 & INTENT(IN) :: wt
349  INTEGER, INTENT(IN), OPTIONAL :: istart_in, iend_in, jstart_in, &
350 & jend_in
351  INTEGER :: i, j, ic, jc
352  INTEGER :: istart, iend, jstart, jend
353  INTEGER :: is, ie, js, je
354  INTEGER :: isd, ied, jsd, jed
355  INTRINSIC PRESENT
356  is = bd%is
357  ie = bd%ie
358  js = bd%js
359  je = bd%je
360  isd = bd%isd
361  ied = bd%ied
362  jsd = bd%jsd
363  jed = bd%jed
364  IF (PRESENT(istart_in)) THEN
365  istart = istart_in
366  ELSE
367  istart = isd
368  END IF
369  IF (PRESENT(iend_in)) THEN
370  iend = iend_in + istag
371  ELSE
372  iend = ied + istag
373  END IF
374  IF (PRESENT(jstart_in)) THEN
375  jstart = jstart_in
376  ELSE
377  jstart = jsd
378  END IF
379  IF (PRESENT(jend_in)) THEN
380  jend = jend_in + jstag
381  ELSE
382  jend = jed + jstag
383  END IF
384  DO j=jstart,jend
385  DO i=istart,iend
386  ic = ind(i, j, 1)
387  jc = ind(i, j, 2)
388  var_nest(i, j) = wt(i, j, 1)*var_coarse(ic, jc) + wt(i, j, 2)*&
389 & var_coarse(ic, jc+1) + wt(i, j, 3)*var_coarse(ic+1, jc+1) + wt&
390 & (i, j, 4)*var_coarse(ic+1, jc)
391  END DO
392  END DO
393  END SUBROUTINE fill_nested_grid_2d
394  SUBROUTINE fill_nested_grid_3d(var_nest, var_coarse, ind, wt, istag, &
395 & jstag, isg, ieg, jsg, jeg, npz, bd, istart_in, iend_in, jstart_in, &
396 & jend_in)
397  IMPLICIT NONE
398  TYPE(FV_GRID_BOUNDS_TYPE), INTENT(IN) :: bd
399  INTEGER, INTENT(IN) :: istag, jstag, isg, ieg, jsg, jeg, npz
400  REAL, DIMENSION(bd%isd:bd%ied+istag, bd%jsd:bd%jed+jstag, npz), &
401 & INTENT(INOUT) :: var_nest
402  REAL, DIMENSION(isg:ieg+istag, jsg:jeg+jstag, npz), INTENT(IN) :: &
403 & var_coarse
404  INTEGER, DIMENSION(bd%isd:bd%ied+istag, bd%jsd:bd%jed+jstag, 2), &
405 & INTENT(IN) :: ind
406  REAL, DIMENSION(bd%isd:bd%ied+istag, bd%jsd:bd%jed+jstag, 4), &
407 & INTENT(IN) :: wt
408  INTEGER, INTENT(IN), OPTIONAL :: istart_in, iend_in, jstart_in, &
409 & jend_in
410  INTEGER :: i, j, ic, jc, k
411  INTEGER :: istart, iend, jstart, jend
412  INTEGER :: is, ie, js, je
413  INTEGER :: isd, ied, jsd, jed
414  INTRINSIC PRESENT
415  is = bd%is
416  ie = bd%ie
417  js = bd%js
418  je = bd%je
419  isd = bd%isd
420  ied = bd%ied
421  jsd = bd%jsd
422  jed = bd%jed
423  IF (PRESENT(istart_in)) THEN
424  istart = istart_in
425  ELSE
426  istart = isd
427  END IF
428  IF (PRESENT(iend_in)) THEN
429  iend = iend_in + istag
430  ELSE
431  iend = ied + istag
432  END IF
433  IF (PRESENT(jstart_in)) THEN
434  jstart = jstart_in
435  ELSE
436  jstart = jsd
437  END IF
438  IF (PRESENT(jend_in)) THEN
439  jend = jend_in + jstag
440  ELSE
441  jend = jed + jstag
442  END IF
443  DO k=1,npz
444  DO j=jstart,jend
445  DO i=istart,iend
446  ic = ind(i, j, 1)
447  jc = ind(i, j, 2)
448  var_nest(i, j, k) = wt(i, j, 1)*var_coarse(ic, jc, k) + wt(i, &
449 & j, 2)*var_coarse(ic, jc+1, k) + wt(i, j, 3)*var_coarse(ic+1&
450 & , jc+1, k) + wt(i, j, 4)*var_coarse(ic+1, jc, k)
451  END DO
452  END DO
453  END DO
454  END SUBROUTINE fill_nested_grid_3d
455  SUBROUTINE nested_grid_bc_mpp(var_nest, var_coarse, nest_domain, ind, &
456 & wt, istag, jstag, npx, npy, npz, bd, isg, ieg, jsg, jeg, nstep_in, &
457 & nsplit_in, proc_in)
458  IMPLICIT NONE
459  TYPE(FV_GRID_BOUNDS_TYPE), INTENT(IN) :: bd
460  INTEGER, INTENT(IN) :: istag, jstag, npx, npy, npz, isg, ieg, jsg, &
461 & jeg
462  REAL, DIMENSION(bd%isd:bd%ied+istag, bd%jsd:bd%jed+jstag, npz), &
463 & INTENT(INOUT) :: var_nest
464  REAL, DIMENSION(isg:ieg+istag, jsg:jeg+jstag, npz), INTENT(IN) :: &
465 & var_coarse
466  TYPE(NEST_DOMAIN_TYPE), INTENT(INOUT) :: nest_domain
467  INTEGER, DIMENSION(bd%isd:bd%ied+istag, bd%jsd:bd%jed+jstag, 2), &
468 & INTENT(IN) :: ind
469  REAL, DIMENSION(bd%isd:bd%ied+istag, bd%jsd:bd%jed+jstag, 4), &
470 & INTENT(IN) :: wt
471  INTEGER, INTENT(IN), OPTIONAL :: nstep_in, nsplit_in
472  LOGICAL, INTENT(IN), OPTIONAL :: proc_in
473  INTEGER :: isw_f, iew_f, jsw_f, jew_f, isw_c, iew_c, jsw_c, jew_c
474  INTEGER :: ise_f, iee_f, jse_f, jee_f, ise_c, iee_c, jse_c, jee_c
475  INTEGER :: iss_f, ies_f, jss_f, jes_f, iss_c, ies_c, jss_c, jes_c
476  INTEGER :: isn_f, ien_f, jsn_f, jen_f, isn_c, ien_c, jsn_c, jen_c
477  REAL, ALLOCATABLE :: wbuffer(:, :, :)
478  REAL, ALLOCATABLE :: ebuffer(:, :, :)
479  REAL, ALLOCATABLE :: sbuffer(:, :, :)
480  REAL, ALLOCATABLE :: nbuffer(:, :, :)
481  INTEGER :: i, j, ic, jc, istart, iend, k
482  INTEGER :: position
483  LOGICAL :: process
484  INTEGER :: is, ie, js, je
485  INTEGER :: isd, ied, jsd, jed
486  INTRINSIC PRESENT
487  is = bd%is
488  ie = bd%ie
489  js = bd%js
490  je = bd%je
491  isd = bd%isd
492  ied = bd%ied
493  jsd = bd%jsd
494  jed = bd%jed
495  IF (PRESENT(proc_in)) THEN
496  process = proc_in
497  ELSE
498  process = .true.
499  END IF
500  IF (istag .EQ. 1 .AND. jstag .EQ. 1) THEN
501  position = corner
502  ELSE IF (istag .EQ. 0 .AND. jstag .EQ. 1) THEN
503  position = north
504  ELSE IF (istag .EQ. 1 .AND. jstag .EQ. 0) THEN
505  position = east
506  ELSE
507  position = center
508  END IF
509  CALL mpp_get_c2f_index(nest_domain, isw_f, iew_f, jsw_f, jew_f, &
510 & isw_c, iew_c, jsw_c, jew_c, west, position)
511  CALL mpp_get_c2f_index(nest_domain, ise_f, iee_f, jse_f, jee_f, &
512 & ise_c, iee_c, jse_c, jee_c, east, position)
513  CALL mpp_get_c2f_index(nest_domain, iss_f, ies_f, jss_f, jes_f, &
514 & iss_c, ies_c, jss_c, jes_c, south, position)
515  CALL mpp_get_c2f_index(nest_domain, isn_f, ien_f, jsn_f, jen_f, &
516 & isn_c, ien_c, jsn_c, jen_c, north, position)
517  IF (iew_c .GE. isw_c .AND. jew_c .GE. jsw_c) THEN
518  ALLOCATE(wbuffer(isw_c:iew_c, jsw_c:jew_c, npz))
519  ELSE
520  ALLOCATE(wbuffer(1, 1, 1))
521  END IF
522  wbuffer = 0.0
523  IF (iee_c .GE. ise_c .AND. jee_c .GE. jse_c) THEN
524  ALLOCATE(ebuffer(ise_c:iee_c, jse_c:jee_c, npz))
525  ELSE
526  ALLOCATE(ebuffer(1, 1, 1))
527  END IF
528  ebuffer = 0.0
529  IF (ies_c .GE. iss_c .AND. jes_c .GE. jss_c) THEN
530  ALLOCATE(sbuffer(iss_c:ies_c, jss_c:jes_c, npz))
531  ELSE
532  ALLOCATE(sbuffer(1, 1, 1))
533  END IF
534  sbuffer = 0.0
535  IF (ien_c .GE. isn_c .AND. jen_c .GE. jsn_c) THEN
536  ALLOCATE(nbuffer(isn_c:ien_c, jsn_c:jen_c, npz))
537  ELSE
538  ALLOCATE(nbuffer(1, 1, 1))
539  END IF
540  nbuffer = 0.0
541  CALL timing_on('COMM_TOTAL')
542  CALL mpp_update_nest_fine(var_coarse, nest_domain, wbuffer, sbuffer&
543 & , ebuffer, nbuffer, position)
544  CALL timing_off('COMM_TOTAL')
545 !process
546  IF (process) THEN
547  IF (is .EQ. 1) THEN
548  DO k=1,npz
549  DO j=jsd,jed+jstag
550  DO i=isd,0
551  ic = ind(i, j, 1)
552  jc = ind(i, j, 2)
553  var_nest(i, j, k) = wt(i, j, 1)*wbuffer(ic, jc, k) + wt(i&
554 & , j, 2)*wbuffer(ic, jc+1, k) + wt(i, j, 3)*wbuffer(ic+1&
555 & , jc+1, k) + wt(i, j, 4)*wbuffer(ic+1, jc, k)
556  END DO
557  END DO
558  END DO
559  END IF
560  IF (js .EQ. 1) THEN
561  IF (is .EQ. 1) THEN
562  istart = is
563  ELSE
564  istart = isd
565  END IF
566  IF (ie .EQ. npx - 1) THEN
567  iend = ie
568  ELSE
569  iend = ied
570  END IF
571  DO k=1,npz
572  DO j=jsd,0
573  DO i=istart,iend+istag
574  ic = ind(i, j, 1)
575  jc = ind(i, j, 2)
576  var_nest(i, j, k) = wt(i, j, 1)*sbuffer(ic, jc, k) + wt(i&
577 & , j, 2)*sbuffer(ic, jc+1, k) + wt(i, j, 3)*sbuffer(ic+1&
578 & , jc+1, k) + wt(i, j, 4)*sbuffer(ic+1, jc, k)
579  END DO
580  END DO
581  END DO
582  END IF
583  IF (ie .EQ. npx - 1) THEN
584  DO k=1,npz
585  DO j=jsd,jed+jstag
586  DO i=npx+istag,ied+istag
587  ic = ind(i, j, 1)
588  jc = ind(i, j, 2)
589  var_nest(i, j, k) = wt(i, j, 1)*ebuffer(ic, jc, k) + wt(i&
590 & , j, 2)*ebuffer(ic, jc+1, k) + wt(i, j, 3)*ebuffer(ic+1&
591 & , jc+1, k) + wt(i, j, 4)*ebuffer(ic+1, jc, k)
592  END DO
593  END DO
594  END DO
595  END IF
596  IF (je .EQ. npy - 1) THEN
597  IF (is .EQ. 1) THEN
598  istart = is
599  ELSE
600  istart = isd
601  END IF
602  IF (ie .EQ. npx - 1) THEN
603  iend = ie
604  ELSE
605  iend = ied
606  END IF
607  DO k=1,npz
608  DO j=npy+jstag,jed+jstag
609  DO i=istart,iend+istag
610  ic = ind(i, j, 1)
611  jc = ind(i, j, 2)
612  var_nest(i, j, k) = wt(i, j, 1)*nbuffer(ic, jc, k) + wt(i&
613 & , j, 2)*nbuffer(ic, jc+1, k) + wt(i, j, 3)*nbuffer(ic+1&
614 & , jc+1, k) + wt(i, j, 4)*nbuffer(ic+1, jc, k)
615  END DO
616  END DO
617  END DO
618  END IF
619  END IF
620  DEALLOCATE(wbuffer)
621  DEALLOCATE(ebuffer)
622  DEALLOCATE(sbuffer)
623  DEALLOCATE(nbuffer)
624  END SUBROUTINE nested_grid_bc_mpp
625  SUBROUTINE nested_grid_bc_mpp_send(var_coarse, nest_domain, istag, &
626 & jstag)
627  IMPLICIT NONE
628  REAL, DIMENSION(:, :, :), INTENT(IN) :: var_coarse
629  TYPE(NEST_DOMAIN_TYPE), INTENT(INOUT) :: nest_domain
630  INTEGER, INTENT(IN) :: istag, jstag
631  REAL, ALLOCATABLE :: wbuffer(:, :, :)
632  REAL, ALLOCATABLE :: ebuffer(:, :, :)
633  REAL, ALLOCATABLE :: sbuffer(:, :, :)
634  REAL, ALLOCATABLE :: nbuffer(:, :, :)
635  INTEGER :: i, j, ic, jc, istart, iend, k
636  INTEGER :: position
637  IF (istag .EQ. 1 .AND. jstag .EQ. 1) THEN
638  position = corner
639  ELSE IF (istag .EQ. 0 .AND. jstag .EQ. 1) THEN
640  position = north
641  ELSE IF (istag .EQ. 1 .AND. jstag .EQ. 0) THEN
642  position = east
643  ELSE
644  position = center
645  END IF
646  ALLOCATE(wbuffer(1, 1, 1))
647  ALLOCATE(ebuffer(1, 1, 1))
648  ALLOCATE(sbuffer(1, 1, 1))
649  ALLOCATE(nbuffer(1, 1, 1))
650  CALL timing_on('COMM_TOTAL')
651  CALL mpp_update_nest_fine(var_coarse, nest_domain, wbuffer, sbuffer&
652 & , ebuffer, nbuffer, position)
653  CALL timing_off('COMM_TOTAL')
654  DEALLOCATE(wbuffer)
655  DEALLOCATE(ebuffer)
656  DEALLOCATE(sbuffer)
657  DEALLOCATE(nbuffer)
658  END SUBROUTINE nested_grid_bc_mpp_send
659  SUBROUTINE nested_grid_bc_2d_mpp(var_nest, var_coarse, nest_domain, &
660 & ind, wt, istag, jstag, npx, npy, bd, isg, ieg, jsg, jeg, nstep_in, &
661 & nsplit_in, proc_in)
662  IMPLICIT NONE
663  TYPE(FV_GRID_BOUNDS_TYPE), INTENT(IN) :: bd
664  INTEGER, INTENT(IN) :: istag, jstag, npx, npy, isg, ieg, jsg, jeg
665  REAL, DIMENSION(bd%isd:bd%ied+istag, bd%jsd:bd%jed+jstag), INTENT(&
666 & INOUT) :: var_nest
667  REAL, DIMENSION(isg:ieg+istag, jsg:jeg+jstag), INTENT(IN) :: &
668 & var_coarse
669  TYPE(NEST_DOMAIN_TYPE), INTENT(INOUT) :: nest_domain
670  INTEGER, DIMENSION(bd%isd:bd%ied+istag, bd%jsd:bd%jed+jstag, 2), &
671 & INTENT(IN) :: ind
672  REAL, DIMENSION(bd%isd:bd%ied+istag, bd%jsd:bd%jed+jstag, 4), &
673 & INTENT(IN) :: wt
674  INTEGER, INTENT(IN), OPTIONAL :: nstep_in, nsplit_in
675  LOGICAL, INTENT(IN), OPTIONAL :: proc_in
676  INTEGER :: isw_f, iew_f, jsw_f, jew_f, isw_c, iew_c, jsw_c, jew_c
677  INTEGER :: ise_f, iee_f, jse_f, jee_f, ise_c, iee_c, jse_c, jee_c
678  INTEGER :: iss_f, ies_f, jss_f, jes_f, iss_c, ies_c, jss_c, jes_c
679  INTEGER :: isn_f, ien_f, jsn_f, jen_f, isn_c, ien_c, jsn_c, jen_c
680  REAL, ALLOCATABLE :: wbuffer(:, :)
681  REAL, ALLOCATABLE :: ebuffer(:, :)
682  REAL, ALLOCATABLE :: sbuffer(:, :)
683  REAL, ALLOCATABLE :: nbuffer(:, :)
684  INTEGER :: i, j, ic, jc, istart, iend, k
685  INTEGER :: position
686  LOGICAL :: process
687  INTEGER :: is, ie, js, je
688  INTEGER :: isd, ied, jsd, jed
689  INTRINSIC PRESENT
690  is = bd%is
691  ie = bd%ie
692  js = bd%js
693  je = bd%je
694  isd = bd%isd
695  ied = bd%ied
696  jsd = bd%jsd
697  jed = bd%jed
698  IF (PRESENT(proc_in)) THEN
699  process = proc_in
700  ELSE
701  process = .true.
702  END IF
703  IF (istag .EQ. 1 .AND. jstag .EQ. 1) THEN
704  position = corner
705  ELSE IF (istag .EQ. 0 .AND. jstag .EQ. 1) THEN
706  position = north
707  ELSE IF (istag .EQ. 1 .AND. jstag .EQ. 0) THEN
708  position = east
709  ELSE
710  position = center
711  END IF
712  CALL mpp_get_c2f_index(nest_domain, isw_f, iew_f, jsw_f, jew_f, &
713 & isw_c, iew_c, jsw_c, jew_c, west, position)
714  CALL mpp_get_c2f_index(nest_domain, ise_f, iee_f, jse_f, jee_f, &
715 & ise_c, iee_c, jse_c, jee_c, east, position)
716  CALL mpp_get_c2f_index(nest_domain, iss_f, ies_f, jss_f, jes_f, &
717 & iss_c, ies_c, jss_c, jes_c, south, position)
718  CALL mpp_get_c2f_index(nest_domain, isn_f, ien_f, jsn_f, jen_f, &
719 & isn_c, ien_c, jsn_c, jen_c, north, position)
720  IF (iew_c .GE. isw_c .AND. jew_c .GE. jsw_c) THEN
721  ALLOCATE(wbuffer(isw_c:iew_c, jsw_c:jew_c))
722  ELSE
723  ALLOCATE(wbuffer(1, 1))
724  END IF
725  wbuffer = 0.0
726  IF (iee_c .GE. ise_c .AND. jee_c .GE. jse_c) THEN
727  ALLOCATE(ebuffer(ise_c:iee_c, jse_c:jee_c))
728  ELSE
729  ALLOCATE(ebuffer(1, 1))
730  END IF
731  ebuffer = 0.0
732  IF (ies_c .GE. iss_c .AND. jes_c .GE. jss_c) THEN
733  ALLOCATE(sbuffer(iss_c:ies_c, jss_c:jes_c))
734  ELSE
735  ALLOCATE(sbuffer(1, 1))
736  END IF
737  sbuffer = 0.0
738  IF (ien_c .GE. isn_c .AND. jen_c .GE. jsn_c) THEN
739  ALLOCATE(nbuffer(isn_c:ien_c, jsn_c:jen_c))
740  ELSE
741  ALLOCATE(nbuffer(1, 1))
742  END IF
743  nbuffer = 0.0
744  CALL timing_on('COMM_TOTAL')
745  CALL mpp_update_nest_fine(var_coarse, nest_domain, wbuffer, sbuffer&
746 & , ebuffer, nbuffer, position)
747  CALL timing_off('COMM_TOTAL')
748 !process
749  IF (process) THEN
750  IF (is .EQ. 1) THEN
751  DO j=jsd,jed+jstag
752  DO i=isd,0
753  ic = ind(i, j, 1)
754  jc = ind(i, j, 2)
755  var_nest(i, j) = wt(i, j, 1)*wbuffer(ic, jc) + wt(i, j, 2)*&
756 & wbuffer(ic, jc+1) + wt(i, j, 3)*wbuffer(ic+1, jc+1) + wt(i&
757 & , j, 4)*wbuffer(ic+1, jc)
758  END DO
759  END DO
760  END IF
761  IF (js .EQ. 1) THEN
762  IF (is .EQ. 1) THEN
763  istart = is
764  ELSE
765  istart = isd
766  END IF
767  IF (ie .EQ. npx - 1) THEN
768  iend = ie
769  ELSE
770  iend = ied
771  END IF
772  DO j=jsd,0
773  DO i=istart,iend+istag
774  ic = ind(i, j, 1)
775  jc = ind(i, j, 2)
776  var_nest(i, j) = wt(i, j, 1)*sbuffer(ic, jc) + wt(i, j, 2)*&
777 & sbuffer(ic, jc+1) + wt(i, j, 3)*sbuffer(ic+1, jc+1) + wt(i&
778 & , j, 4)*sbuffer(ic+1, jc)
779  END DO
780  END DO
781  END IF
782  IF (ie .EQ. npx - 1) THEN
783  DO j=jsd,jed+jstag
784  DO i=npx+istag,ied+istag
785  ic = ind(i, j, 1)
786  jc = ind(i, j, 2)
787  var_nest(i, j) = wt(i, j, 1)*ebuffer(ic, jc) + wt(i, j, 2)*&
788 & ebuffer(ic, jc+1) + wt(i, j, 3)*ebuffer(ic+1, jc+1) + wt(i&
789 & , j, 4)*ebuffer(ic+1, jc)
790  END DO
791  END DO
792  END IF
793  IF (je .EQ. npy - 1) THEN
794  IF (is .EQ. 1) THEN
795  istart = is
796  ELSE
797  istart = isd
798  END IF
799  IF (ie .EQ. npx - 1) THEN
800  iend = ie
801  ELSE
802  iend = ied
803  END IF
804  DO j=npy+jstag,jed+jstag
805  DO i=istart,iend+istag
806  ic = ind(i, j, 1)
807  jc = ind(i, j, 2)
808  var_nest(i, j) = wt(i, j, 1)*nbuffer(ic, jc) + wt(i, j, 2)*&
809 & nbuffer(ic, jc+1) + wt(i, j, 3)*nbuffer(ic+1, jc+1) + wt(i&
810 & , j, 4)*nbuffer(ic+1, jc)
811  END DO
812  END DO
813  END IF
814  END IF
815  DEALLOCATE(wbuffer)
816  DEALLOCATE(ebuffer)
817  DEALLOCATE(sbuffer)
818  DEALLOCATE(nbuffer)
819  END SUBROUTINE nested_grid_bc_2d_mpp
820  SUBROUTINE nested_grid_bc_2d(var_nest, var_coarse, ind, wt, istag, &
821 & jstag, npx, npy, bd, isg, ieg, jsg, jeg, nstep_in, nsplit_in)
822  IMPLICIT NONE
823  TYPE(FV_GRID_BOUNDS_TYPE), INTENT(IN) :: bd
824  INTEGER, INTENT(IN) :: istag, jstag, npx, npy, isg, ieg, jsg, jeg
825  REAL, DIMENSION(bd%isd:bd%ied+istag, bd%jsd:bd%jed+jstag), INTENT(&
826 & INOUT) :: var_nest
827  REAL, DIMENSION(isg:ieg+istag, jsg:jeg+jstag), INTENT(IN) :: &
828 & var_coarse
829  INTEGER, DIMENSION(bd%isd:bd%ied+istag, bd%jsd:bd%jed+jstag, 2), &
830 & INTENT(IN) :: ind
831  REAL, DIMENSION(bd%isd:bd%ied+istag, bd%jsd:bd%jed+jstag, 4), &
832 & INTENT(IN) :: wt
833  INTEGER, INTENT(IN), OPTIONAL :: nstep_in, nsplit_in
834  INTEGER :: nstep, nsplit
835  INTEGER :: i, j, ic, jc, istart, iend
836  INTEGER :: is, ie, js, je
837  INTEGER :: isd, ied, jsd, jed
838  INTRINSIC PRESENT
839  is = bd%is
840  ie = bd%ie
841  js = bd%js
842  je = bd%je
843  isd = bd%isd
844  ied = bd%ied
845  jsd = bd%jsd
846  jed = bd%jed
847  IF ((.NOT.PRESENT(nstep_in)) .OR. (.NOT.PRESENT(nsplit_in))) THEN
848  nstep = 1
849  nsplit = 2
850  ELSE
851  nstep = nstep_in
852  nsplit = nsplit_in
853  END IF
854  IF (is .EQ. 1) THEN
855  DO j=jsd,jed+jstag
856  DO i=isd,0
857  ic = ind(i, j, 1)
858  jc = ind(i, j, 2)
859  var_nest(i, j) = wt(i, j, 1)*var_coarse(ic, jc) + wt(i, j, 2)*&
860 & var_coarse(ic, jc+1) + wt(i, j, 3)*var_coarse(ic+1, jc+1) + &
861 & wt(i, j, 4)*var_coarse(ic+1, jc)
862  END DO
863  END DO
864  END IF
865  IF (js .EQ. 1) THEN
866  IF (is .EQ. 1) THEN
867  istart = is
868  ELSE
869  istart = isd
870  END IF
871  IF (ie .EQ. npx - 1) THEN
872  iend = ie
873  ELSE
874  iend = ied
875  END IF
876  DO j=jsd,0
877  DO i=istart,iend+istag
878  ic = ind(i, j, 1)
879  jc = ind(i, j, 2)
880  var_nest(i, j) = wt(i, j, 1)*var_coarse(ic, jc) + wt(i, j, 2)*&
881 & var_coarse(ic, jc+1) + wt(i, j, 3)*var_coarse(ic+1, jc+1) + &
882 & wt(i, j, 4)*var_coarse(ic+1, jc)
883  END DO
884  END DO
885  END IF
886  IF (ie .EQ. npx - 1) THEN
887  DO j=jsd,jed+jstag
888  DO i=npx+istag,ied+istag
889  ic = ind(i, j, 1)
890  jc = ind(i, j, 2)
891  var_nest(i, j) = wt(i, j, 1)*var_coarse(ic, jc) + wt(i, j, 2)*&
892 & var_coarse(ic, jc+1) + wt(i, j, 3)*var_coarse(ic+1, jc+1) + &
893 & wt(i, j, 4)*var_coarse(ic+1, jc)
894  END DO
895  END DO
896  END IF
897  IF (je .EQ. npy - 1) THEN
898  IF (is .EQ. 1) THEN
899  istart = is
900  ELSE
901  istart = isd
902  END IF
903  IF (ie .EQ. npx - 1) THEN
904  iend = ie
905  ELSE
906  iend = ied
907  END IF
908  DO j=npy+jstag,jed+jstag
909  DO i=istart,iend+istag
910  ic = ind(i, j, 1)
911  jc = ind(i, j, 2)
912  var_nest(i, j) = wt(i, j, 1)*var_coarse(ic, jc) + wt(i, j, 2)*&
913 & var_coarse(ic, jc+1) + wt(i, j, 3)*var_coarse(ic+1, jc+1) + &
914 & wt(i, j, 4)*var_coarse(ic+1, jc)
915  END DO
916  END DO
917  END IF
918  END SUBROUTINE nested_grid_bc_2d
919  SUBROUTINE nested_grid_bc_3d(var_nest, var_coarse, ind, wt, istag, &
920 & jstag, npx, npy, npz, bd, isg, ieg, jsg, jeg, nstep_in, nsplit_in)
921  IMPLICIT NONE
922  TYPE(FV_GRID_BOUNDS_TYPE), INTENT(IN) :: bd
923  INTEGER, INTENT(IN) :: istag, jstag, npx, npy, isg, ieg, jsg, jeg, &
924 & npz
925  REAL, DIMENSION(bd%isd:bd%ied+istag, bd%jsd:bd%jed+jstag, npz), &
926 & INTENT(INOUT) :: var_nest
927  REAL, DIMENSION(isg:ieg+istag, jsg:jeg+jstag, npz), INTENT(IN) :: &
928 & var_coarse
929  INTEGER, DIMENSION(bd%isd:bd%ied+istag, bd%jsd:bd%jed+jstag, 2), &
930 & INTENT(IN) :: ind
931  REAL, DIMENSION(bd%isd:bd%ied+istag, bd%jsd:bd%jed+jstag, 4), &
932 & INTENT(IN) :: wt
933  INTEGER, INTENT(IN), OPTIONAL :: nstep_in, nsplit_in
934  INTEGER :: nstep, nsplit
935  INTEGER :: i, j, ic, jc, istart, iend, k
936  INTEGER :: is, ie, js, je
937  INTEGER :: isd, ied, jsd, jed
938  INTRINSIC PRESENT
939  is = bd%is
940  ie = bd%ie
941  js = bd%js
942  je = bd%je
943  isd = bd%isd
944  ied = bd%ied
945  jsd = bd%jsd
946  jed = bd%jed
947  IF ((.NOT.PRESENT(nstep_in)) .OR. (.NOT.PRESENT(nsplit_in))) THEN
948  nstep = 1
949  nsplit = 2
950  ELSE
951  nstep = nstep_in
952  nsplit = nsplit_in
953  END IF
954  IF (is .EQ. 1) THEN
955  DO k=1,npz
956  DO j=jsd,jed+jstag
957  DO i=isd,0
958  ic = ind(i, j, 1)
959  jc = ind(i, j, 2)
960  var_nest(i, j, k) = wt(i, j, 1)*var_coarse(ic, jc, k) + wt(i&
961 & , j, 2)*var_coarse(ic, jc+1, k) + wt(i, j, 3)*var_coarse(&
962 & ic+1, jc+1, k) + wt(i, j, 4)*var_coarse(ic+1, jc, k)
963  END DO
964  END DO
965  END DO
966  END IF
967  IF (js .EQ. 1) THEN
968  IF (is .EQ. 1) THEN
969  istart = is
970  ELSE
971  istart = isd
972  END IF
973  IF (ie .EQ. npx - 1) THEN
974  iend = ie
975  ELSE
976  iend = ied
977  END IF
978  DO k=1,npz
979  DO j=jsd,0
980  DO i=istart,iend+istag
981  ic = ind(i, j, 1)
982  jc = ind(i, j, 2)
983  var_nest(i, j, k) = wt(i, j, 1)*var_coarse(ic, jc, k) + wt(i&
984 & , j, 2)*var_coarse(ic, jc+1, k) + wt(i, j, 3)*var_coarse(&
985 & ic+1, jc+1, k) + wt(i, j, 4)*var_coarse(ic+1, jc, k)
986  END DO
987  END DO
988  END DO
989  END IF
990  IF (ie .EQ. npx - 1) THEN
991  DO k=1,npz
992  DO j=jsd,jed+jstag
993  DO i=npx+istag,ied+istag
994  ic = ind(i, j, 1)
995  jc = ind(i, j, 2)
996  var_nest(i, j, k) = wt(i, j, 1)*var_coarse(ic, jc, k) + wt(i&
997 & , j, 2)*var_coarse(ic, jc+1, k) + wt(i, j, 3)*var_coarse(&
998 & ic+1, jc+1, k) + wt(i, j, 4)*var_coarse(ic+1, jc, k)
999  END DO
1000  END DO
1001  END DO
1002  END IF
1003  IF (je .EQ. npy - 1) THEN
1004  IF (is .EQ. 1) THEN
1005  istart = is
1006  ELSE
1007  istart = isd
1008  END IF
1009  IF (ie .EQ. npx - 1) THEN
1010  iend = ie
1011  ELSE
1012  iend = ied
1013  END IF
1014  DO k=1,npz
1015  DO j=npy+jstag,jed+jstag
1016  DO i=istart,iend+istag
1017  ic = ind(i, j, 1)
1018  jc = ind(i, j, 2)
1019  var_nest(i, j, k) = wt(i, j, 1)*var_coarse(ic, jc, k) + wt(i&
1020 & , j, 2)*var_coarse(ic, jc+1, k) + wt(i, j, 3)*var_coarse(&
1021 & ic+1, jc+1, k) + wt(i, j, 4)*var_coarse(ic+1, jc, k)
1022  END DO
1023  END DO
1024  END DO
1025  END IF
1026  END SUBROUTINE nested_grid_bc_3d
1027  SUBROUTINE nested_grid_bc_send(var_coarse, nest_domain, istag, jstag)
1028  IMPLICIT NONE
1029  REAL, DIMENSION(:, :, :), INTENT(IN) :: var_coarse
1030  TYPE(nest_domain_type), INTENT(INOUT) :: nest_domain
1031  INTEGER, INTENT(IN) :: istag, jstag
1032  INTEGER :: position
1033  REAL :: wbuffer(1, 1, 1)
1034  REAL :: ebuffer(1, 1, 1)
1035  REAL :: sbuffer(1, 1, 1)
1036  REAL :: nbuffer(1, 1, 1)
1037  IF (istag .EQ. 1 .AND. jstag .EQ. 1) THEN
1038  position = corner
1039  ELSE IF (istag .EQ. 0 .AND. jstag .EQ. 1) THEN
1040  position = north
1041  ELSE IF (istag .EQ. 1 .AND. jstag .EQ. 0) THEN
1042  position = east
1043  ELSE
1044  position = center
1045  END IF
1046  CALL timing_on('COMM_TOTAL')
1047  CALL mpp_update_nest_fine(var_coarse, nest_domain, wbuffer, sbuffer&
1048 & , ebuffer, nbuffer, position)
1049  CALL timing_off('COMM_TOTAL')
1050  END SUBROUTINE nested_grid_bc_send
1051  SUBROUTINE nested_grid_bc_recv(nest_domain, istag, jstag, npz, bd, &
1052 & nest_bc_buffers)
1053  IMPLICIT NONE
1054  TYPE(fv_grid_bounds_type), INTENT(IN) :: bd
1055  TYPE(nest_domain_type), INTENT(INOUT) :: nest_domain
1056  INTEGER, INTENT(IN) :: istag, jstag, npz
1057  TYPE(fv_nest_bc_type_3d), INTENT(INOUT), TARGET :: nest_bc_buffers
1058  REAL, DIMENSION(bd%isd:bd%ied+istag, bd%jsd:bd%jed+jstag, npz) :: &
1059 & var_coarse_dummy
1060  INTEGER :: position
1061  INTEGER :: isw_f, iew_f, jsw_f, jew_f, isw_c, iew_c, jsw_c, jew_c
1062  INTEGER :: ise_f, iee_f, jse_f, jee_f, ise_c, iee_c, jse_c, jee_c
1063  INTEGER :: iss_f, ies_f, jss_f, jes_f, iss_c, ies_c, jss_c, jes_c
1064  INTEGER :: isn_f, ien_f, jsn_f, jen_f, isn_c, ien_c, jsn_c, jen_c
1065  INTEGER :: i, j, k
1066  INTRINSIC ALLOCATED
1067  var_coarse_dummy = 0.0
1068  IF (istag .EQ. 1 .AND. jstag .EQ. 1) THEN
1069  position = corner
1070  ELSE IF (istag .EQ. 0 .AND. jstag .EQ. 1) THEN
1071  position = north
1072  ELSE IF (istag .EQ. 1 .AND. jstag .EQ. 0) THEN
1073  position = east
1074  ELSE
1075  position = center
1076  END IF
1077  IF (.NOT.ALLOCATED(nest_bc_buffers%west_t1)) THEN
1078  CALL mpp_get_c2f_index(nest_domain, isw_f, iew_f, jsw_f, jew_f, &
1079 & isw_c, iew_c, jsw_c, jew_c, west, position)
1080  CALL mpp_get_c2f_index(nest_domain, ise_f, iee_f, jse_f, jee_f, &
1081 & ise_c, iee_c, jse_c, jee_c, east, position)
1082  CALL mpp_get_c2f_index(nest_domain, iss_f, ies_f, jss_f, jes_f, &
1083 & iss_c, ies_c, jss_c, jes_c, south, position)
1084  CALL mpp_get_c2f_index(nest_domain, isn_f, ien_f, jsn_f, jen_f, &
1085 & isn_c, ien_c, jsn_c, jen_c, north, position)
1086  IF (iew_c .GE. isw_c .AND. jew_c .GE. jsw_c) THEN
1087  IF (.NOT.ALLOCATED(nest_bc_buffers%west_t1)) THEN
1088  ALLOCATE(nest_bc_buffers%west_t1(isw_c:iew_c, jsw_c:jew_c, npz&
1089 & ))
1090  END IF
1091 !compatible with first touch principle
1092  DO k=1,npz
1093  DO j=jsw_c,jew_c
1094  DO i=isw_c,iew_c
1095  nest_bc_buffers%west_t1(i, j, k) = 0.
1096  END DO
1097  END DO
1098  END DO
1099  ELSE
1100  ALLOCATE(nest_bc_buffers%west_t1(1, 1, 1))
1101  nest_bc_buffers%west_t1(1, 1, 1) = 0.
1102  END IF
1103  IF (iee_c .GE. ise_c .AND. jee_c .GE. jse_c) THEN
1104  IF (.NOT.ALLOCATED(nest_bc_buffers%east_t1)) THEN
1105  ALLOCATE(nest_bc_buffers%east_t1(ise_c:iee_c, jse_c:jee_c, npz&
1106 & ))
1107  END IF
1108  DO k=1,npz
1109  DO j=jse_c,jee_c
1110  DO i=ise_c,iee_c
1111  nest_bc_buffers%east_t1(i, j, k) = 0.
1112  END DO
1113  END DO
1114  END DO
1115  ELSE
1116  ALLOCATE(nest_bc_buffers%east_t1(1, 1, 1))
1117  nest_bc_buffers%east_t1(1, 1, 1) = 0.
1118  END IF
1119  IF (ies_c .GE. iss_c .AND. jes_c .GE. jss_c) THEN
1120  IF (.NOT.ALLOCATED(nest_bc_buffers%south_t1)) THEN
1121  ALLOCATE(nest_bc_buffers%south_t1(iss_c:ies_c, jss_c:jes_c, &
1122 & npz))
1123  END IF
1124  DO k=1,npz
1125  DO j=jss_c,jes_c
1126  DO i=iss_c,ies_c
1127  nest_bc_buffers%south_t1(i, j, k) = 0.
1128  END DO
1129  END DO
1130  END DO
1131  ELSE
1132  ALLOCATE(nest_bc_buffers%south_t1(1, 1, 1))
1133  nest_bc_buffers%south_t1(1, 1, 1) = 0.
1134  END IF
1135  IF (ien_c .GE. isn_c .AND. jen_c .GE. jsn_c) THEN
1136  IF (.NOT.ALLOCATED(nest_bc_buffers%north_t1)) THEN
1137  ALLOCATE(nest_bc_buffers%north_t1(isn_c:ien_c, jsn_c:jen_c, &
1138 & npz))
1139  END IF
1140  DO k=1,npz
1141  DO j=jsn_c,jen_c
1142  DO i=isn_c,ien_c
1143  nest_bc_buffers%north_t1(i, j, k) = 0.
1144  END DO
1145  END DO
1146  END DO
1147  ELSE
1148  ALLOCATE(nest_bc_buffers%north_t1(1, 1, 1))
1149  nest_bc_buffers%north_t1(1, 1, 1) = 0
1150  END IF
1151  END IF
1152  CALL timing_on('COMM_TOTAL')
1153  CALL mpp_update_nest_fine(var_coarse_dummy, nest_domain, &
1154 & nest_bc_buffers%west_t1, nest_bc_buffers%&
1155 & south_t1, nest_bc_buffers%east_t1, &
1156 & nest_bc_buffers%north_t1, position)
1157  CALL timing_off('COMM_TOTAL')
1158  END SUBROUTINE nested_grid_bc_recv
1159  SUBROUTINE nested_grid_bc_save_proc(nest_domain, ind, wt, istag, jstag&
1160 & , npx, npy, npz, bd, nest_bc, nest_bc_buffers, pd_in)
1161  IMPLICIT NONE
1162  TYPE(fv_grid_bounds_type), INTENT(IN) :: bd
1163  TYPE(nest_domain_type), INTENT(INOUT) :: nest_domain
1164  INTEGER, INTENT(IN) :: istag, jstag, npx, npy, npz
1165  INTEGER, DIMENSION(bd%isd:bd%ied+istag, bd%jsd:bd%jed+jstag, 2), &
1166 & INTENT(IN) :: ind
1167  REAL, DIMENSION(bd%isd:bd%ied+istag, bd%jsd:bd%jed+jstag, 4), &
1168 & INTENT(IN) :: wt
1169  LOGICAL, INTENT(IN), OPTIONAL :: pd_in
1170 !!NOTE: if declaring an ALLOCATABLE array with intent(OUT), the resulting dummy array
1171 !! will NOT be allocated! This goes for allocatable members of derived types as well.
1172  TYPE(fv_nest_bc_type_3d), INTENT(INOUT), TARGET :: nest_bc, &
1173 & nest_bc_buffers
1174  REAL, DIMENSION(bd%isd:bd%ied+istag, bd%jsd:bd%jed+jstag, npz) :: &
1175 & var_coarse_dummy
1176  REAL, DIMENSION(:, :, :), POINTER :: var_east, var_west, var_south, &
1177 & var_north
1178  REAL, DIMENSION(:, :, :), POINTER :: buf_east, buf_west, buf_south, &
1179 & buf_north
1180  INTEGER :: position
1181  INTEGER :: i, j, k, ic, jc, istart, iend
1182  LOGICAL :: process
1183  LOGICAL, SAVE :: pd=.false.
1184  INTEGER :: is, ie, js, je
1185  INTEGER :: isd, ied, jsd, jed
1186  INTRINSIC PRESENT
1187  INTRINSIC max
1188  is = bd%is
1189  ie = bd%ie
1190  js = bd%js
1191  je = bd%je
1192  isd = bd%isd
1193  ied = bd%ied
1194  jsd = bd%jsd
1195  jed = bd%jed
1196  IF (PRESENT(pd_in)) THEN
1197  pd = pd_in
1198  ELSE
1199  pd = .false.
1200  END IF
1201  var_east => nest_bc%east_t1
1202  var_west => nest_bc%west_t1
1203  var_north => nest_bc%north_t1
1204  var_south => nest_bc%south_t1
1205  buf_east => nest_bc_buffers%east_t1
1206  buf_west => nest_bc_buffers%west_t1
1207  buf_north => nest_bc_buffers%north_t1
1208  buf_south => nest_bc_buffers%south_t1
1209 ! ?buffer has uninterpolated coarse-grid data; need to perform interpolation ourselves
1210 !To do this more securely, instead of using is/etc we could use the fine-grid indices defined above
1211  IF (is .EQ. 1) THEN
1212 !$NO-MP parallel do default(none) shared(npz,isd,ied,jsd,jed,jstag,ind,var_west,wt,buf_west) private(ic,jc)
1213  DO k=1,npz
1214  DO j=jsd,jed+jstag
1215  DO i=isd,0
1216  ic = ind(i, j, 1)
1217  jc = ind(i, j, 2)
1218  var_west(i, j, k) = wt(i, j, 1)*buf_west(ic, jc, k) + wt(i, &
1219 & j, 2)*buf_west(ic, jc+1, k) + wt(i, j, 3)*buf_west(ic+1, &
1220 & jc+1, k) + wt(i, j, 4)*buf_west(ic+1, jc, k)
1221  END DO
1222  END DO
1223  END DO
1224  IF (pd) THEN
1225 !$NO-MP parallel do default(none) shared(npz,jsd,jed,jstag,isd,var_west,nest_BC)
1226  DO k=1,npz
1227  DO j=jsd,jed+jstag
1228  DO i=isd,0
1229  IF (var_west(i, j, k) .LT. 0.5*nest_bc%west_t0(i, j, k)) &
1230 & THEN
1231  var_west(i, j, k) = 0.5*nest_bc%west_t0(i, j, k)
1232  ELSE
1233  var_west(i, j, k) = var_west(i, j, k)
1234  END IF
1235  END DO
1236  END DO
1237  END DO
1238  END IF
1239  END IF
1240  IF (js .EQ. 1) THEN
1241  IF (is .EQ. 1) THEN
1242  istart = is
1243  ELSE
1244  istart = isd
1245  END IF
1246  IF (ie .EQ. npx - 1) THEN
1247  iend = ie
1248  ELSE
1249  iend = ied
1250  END IF
1251 !$NO-MP parallel do default(none) shared(npz,istart,iend,jsd,jed,istag,ind,var_south,wt,buf_south) private(ic,jc)
1252  DO k=1,npz
1253  DO j=jsd,0
1254  DO i=istart,iend+istag
1255  ic = ind(i, j, 1)
1256  jc = ind(i, j, 2)
1257  var_south(i, j, k) = wt(i, j, 1)*buf_south(ic, jc, k) + wt(i&
1258 & , j, 2)*buf_south(ic, jc+1, k) + wt(i, j, 3)*buf_south(ic+&
1259 & 1, jc+1, k) + wt(i, j, 4)*buf_south(ic+1, jc, k)
1260  END DO
1261  END DO
1262  END DO
1263  IF (pd) THEN
1264 !$NO-MP parallel do default(none) shared(npz,jsd,jed,istart,iend,istag,var_south,nest_BC)
1265  DO k=1,npz
1266  DO j=jsd,0
1267  DO i=istart,iend+istag
1268  IF (var_south(i, j, k) .LT. 0.5*nest_bc%south_t0(i, j, k)&
1269 & ) THEN
1270  var_south(i, j, k) = 0.5*nest_bc%south_t0(i, j, k)
1271  ELSE
1272  var_south(i, j, k) = var_south(i, j, k)
1273  END IF
1274  END DO
1275  END DO
1276  END DO
1277  END IF
1278  END IF
1279  IF (ie .EQ. npx - 1) THEN
1280 !$NO-MP parallel do default(none) shared(npx,npz,isd,ied,jsd,jed,istag,jstag,ind,var_east,wt,buf_east) private(ic,jc)
1281  DO k=1,npz
1282  DO j=jsd,jed+jstag
1283  DO i=npx+istag,ied+istag
1284  ic = ind(i, j, 1)
1285  jc = ind(i, j, 2)
1286  var_east(i, j, k) = wt(i, j, 1)*buf_east(ic, jc, k) + wt(i, &
1287 & j, 2)*buf_east(ic, jc+1, k) + wt(i, j, 3)*buf_east(ic+1, &
1288 & jc+1, k) + wt(i, j, 4)*buf_east(ic+1, jc, k)
1289  END DO
1290  END DO
1291  END DO
1292  IF (pd) THEN
1293 !$NO-MP parallel do default(none) shared(npx,npz,jsd,jed,istag,jstag,ied,var_east,nest_BC)
1294  DO k=1,npz
1295  DO j=jsd,jed+jstag
1296  DO i=npx+istag,ied+istag
1297  IF (var_east(i, j, k) .LT. 0.5*nest_bc%east_t0(i, j, k)) &
1298 & THEN
1299  var_east(i, j, k) = 0.5*nest_bc%east_t0(i, j, k)
1300  ELSE
1301  var_east(i, j, k) = var_east(i, j, k)
1302  END IF
1303  END DO
1304  END DO
1305  END DO
1306  END IF
1307  END IF
1308  IF (je .EQ. npy - 1) THEN
1309  IF (is .EQ. 1) THEN
1310  istart = is
1311  ELSE
1312  istart = isd
1313  END IF
1314  IF (ie .EQ. npx - 1) THEN
1315  iend = ie
1316  ELSE
1317  iend = ied
1318  END IF
1319 !$NO-MP parallel do default(none) shared(npy,npz,istart,iend,jsd,jed,istag,jstag,ind,var_north,wt,buf_north) private(ic,jc)
1320  DO k=1,npz
1321  DO j=npy+jstag,jed+jstag
1322  DO i=istart,iend+istag
1323  ic = ind(i, j, 1)
1324  jc = ind(i, j, 2)
1325  var_north(i, j, k) = wt(i, j, 1)*buf_north(ic, jc, k) + wt(i&
1326 & , j, 2)*buf_north(ic, jc+1, k) + wt(i, j, 3)*buf_north(ic+&
1327 & 1, jc+1, k) + wt(i, j, 4)*buf_north(ic+1, jc, k)
1328  END DO
1329  END DO
1330  END DO
1331  IF (pd) THEN
1332 !$NO-MP parallel do default(none) shared(npy,npz,jsd,jed,istart,iend,istag,jstag,ied,var_north,nest_BC)
1333  DO k=1,npz
1334  DO j=npy+jstag,jed+jstag
1335  DO i=istart,iend+istag
1336  IF (var_north(i, j, k) .LT. 0.5*nest_bc%north_t0(i, j, k)&
1337 & ) THEN
1338  var_north(i, j, k) = 0.5*nest_bc%north_t0(i, j, k)
1339  ELSE
1340  var_north(i, j, k) = var_north(i, j, k)
1341  END IF
1342  END DO
1343  END DO
1344  END DO
1345  END IF
1346  END IF
1347  END SUBROUTINE nested_grid_bc_save_proc
1348 ! Differentiation of nested_grid_bc_apply_intt in forward (tangent) mode:
1349 ! variations of useful results: var_nest
1350 ! with respect to varying inputs: var_nest
1351 ! A NOTE ON BCTYPE: currently only an interpolation BC is implemented,
1352 ! bctype >= 2 currently correspond
1353 ! to a flux BC on the tracers ONLY, which is implemented in fv_tracer.
1354  SUBROUTINE nested_grid_bc_apply_intt_tlm(var_nest, var_nest_tl, istag&
1355 & , jstag, npx, npy, npz, bd, step, split, bc, bctype)
1356  IMPLICIT NONE
1357  TYPE(fv_grid_bounds_type), INTENT(IN) :: bd
1358  INTEGER, INTENT(IN) :: istag, jstag, npx, npy, npz
1359  REAL, DIMENSION(bd%isd:bd%ied+istag, bd%jsd:bd%jed+jstag, npz), &
1360 & INTENT(INOUT) :: var_nest
1361  REAL, DIMENSION(bd%isd:bd%ied+istag, bd%jsd:bd%jed+jstag, npz), &
1362 & INTENT(INOUT) :: var_nest_tl
1363  REAL, INTENT(IN) :: split, step
1364  INTEGER, INTENT(IN) :: bctype
1365  TYPE(fv_nest_bc_type_3d), INTENT(IN), TARGET :: bc
1366  REAL, DIMENSION(:, :, :), POINTER :: var_t0, var_t1
1367  INTEGER :: i, j, istart, iend, k
1368  REAL :: denom
1369  LOGICAL, SAVE :: printdiag=.true.
1370  INTEGER :: is, ie, js, je
1371  INTEGER :: isd, ied, jsd, jed
1372  is = bd%is
1373  ie = bd%ie
1374  js = bd%js
1375  je = bd%je
1376  isd = bd%isd
1377  ied = bd%ied
1378  jsd = bd%jsd
1379  jed = bd%jed
1380  denom = 1./split
1381  IF (is .EQ. 1) THEN
1382  var_t0 => bc%west_t0
1383  var_t1 => bc%west_t1
1384  DO k=1,npz
1385  DO j=jsd,jed+jstag
1386  DO i=isd,0
1387  var_nest_tl(i, j, k) = 0.0
1388  var_nest(i, j, k) = (var_t0(i, j, k)*(split-step)+step*&
1389 & var_t1(i, j, k))*denom
1390  END DO
1391  END DO
1392  END DO
1393  END IF
1394  IF (js .EQ. 1) THEN
1395  IF (is .EQ. 1) THEN
1396  istart = is
1397  ELSE
1398  istart = isd
1399  END IF
1400  IF (ie .EQ. npx - 1) THEN
1401  iend = ie
1402  ELSE
1403  iend = ied
1404  END IF
1405  var_t0 => bc%south_t0
1406  var_t1 => bc%south_t1
1407  DO k=1,npz
1408  DO j=jsd,0
1409  DO i=istart,iend+istag
1410  var_nest_tl(i, j, k) = 0.0
1411  var_nest(i, j, k) = (var_t0(i, j, k)*(split-step)+step*&
1412 & var_t1(i, j, k))*denom
1413  END DO
1414  END DO
1415  END DO
1416  END IF
1417  IF (ie .EQ. npx - 1) THEN
1418  var_t0 => bc%east_t0
1419  var_t1 => bc%east_t1
1420  DO k=1,npz
1421  DO j=jsd,jed+jstag
1422  DO i=npx+istag,ied+istag
1423  var_nest_tl(i, j, k) = 0.0
1424  var_nest(i, j, k) = (var_t0(i, j, k)*(split-step)+step*&
1425 & var_t1(i, j, k))*denom
1426  END DO
1427  END DO
1428  END DO
1429  END IF
1430  IF (je .EQ. npy - 1) THEN
1431  IF (is .EQ. 1) THEN
1432  istart = is
1433  ELSE
1434  istart = isd
1435  END IF
1436  IF (ie .EQ. npx - 1) THEN
1437  iend = ie
1438  ELSE
1439  iend = ied
1440  END IF
1441  var_t0 => bc%north_t0
1442  var_t1 => bc%north_t1
1443  DO k=1,npz
1444  DO j=npy+jstag,jed+jstag
1445  DO i=istart,iend+istag
1446  var_nest_tl(i, j, k) = 0.0
1447  var_nest(i, j, k) = (var_t0(i, j, k)*(split-step)+step*&
1448 & var_t1(i, j, k))*denom
1449  END DO
1450  END DO
1451  END DO
1452  END IF
1453  END SUBROUTINE nested_grid_bc_apply_intt_tlm
1454 ! A NOTE ON BCTYPE: currently only an interpolation BC is implemented,
1455 ! bctype >= 2 currently correspond
1456 ! to a flux BC on the tracers ONLY, which is implemented in fv_tracer.
1457  SUBROUTINE nested_grid_bc_apply_intt(var_nest, istag, jstag, npx, npy&
1458 & , npz, bd, step, split, bc, bctype)
1459  IMPLICIT NONE
1460  TYPE(fv_grid_bounds_type), INTENT(IN) :: bd
1461  INTEGER, INTENT(IN) :: istag, jstag, npx, npy, npz
1462  REAL, DIMENSION(bd%isd:bd%ied+istag, bd%jsd:bd%jed+jstag, npz), &
1463 & INTENT(INOUT) :: var_nest
1464  REAL, INTENT(IN) :: split, step
1465  INTEGER, INTENT(IN) :: bctype
1466  TYPE(fv_nest_bc_type_3d), INTENT(IN), TARGET :: bc
1467  REAL, DIMENSION(:, :, :), POINTER :: var_t0, var_t1
1468  INTEGER :: i, j, istart, iend, k
1469  REAL :: denom
1470  LOGICAL, SAVE :: printdiag=.true.
1471  INTEGER :: is, ie, js, je
1472  INTEGER :: isd, ied, jsd, jed
1473  is = bd%is
1474  ie = bd%ie
1475  js = bd%js
1476  je = bd%je
1477  isd = bd%isd
1478  ied = bd%ied
1479  jsd = bd%jsd
1480  jed = bd%jed
1481  denom = 1./split
1482  IF (is .EQ. 1) THEN
1483  var_t0 => bc%west_t0
1484  var_t1 => bc%west_t1
1485  DO k=1,npz
1486  DO j=jsd,jed+jstag
1487  DO i=isd,0
1488  var_nest(i, j, k) = (var_t0(i, j, k)*(split-step)+step*&
1489 & var_t1(i, j, k))*denom
1490  END DO
1491  END DO
1492  END DO
1493  END IF
1494  IF (js .EQ. 1) THEN
1495  IF (is .EQ. 1) THEN
1496  istart = is
1497  ELSE
1498  istart = isd
1499  END IF
1500  IF (ie .EQ. npx - 1) THEN
1501  iend = ie
1502  ELSE
1503  iend = ied
1504  END IF
1505  var_t0 => bc%south_t0
1506  var_t1 => bc%south_t1
1507  DO k=1,npz
1508  DO j=jsd,0
1509  DO i=istart,iend+istag
1510  var_nest(i, j, k) = (var_t0(i, j, k)*(split-step)+step*&
1511 & var_t1(i, j, k))*denom
1512  END DO
1513  END DO
1514  END DO
1515  END IF
1516  IF (ie .EQ. npx - 1) THEN
1517  var_t0 => bc%east_t0
1518  var_t1 => bc%east_t1
1519  DO k=1,npz
1520  DO j=jsd,jed+jstag
1521  DO i=npx+istag,ied+istag
1522  var_nest(i, j, k) = (var_t0(i, j, k)*(split-step)+step*&
1523 & var_t1(i, j, k))*denom
1524  END DO
1525  END DO
1526  END DO
1527  END IF
1528  IF (je .EQ. npy - 1) THEN
1529  IF (is .EQ. 1) THEN
1530  istart = is
1531  ELSE
1532  istart = isd
1533  END IF
1534  IF (ie .EQ. npx - 1) THEN
1535  iend = ie
1536  ELSE
1537  iend = ied
1538  END IF
1539  var_t0 => bc%north_t0
1540  var_t1 => bc%north_t1
1541  DO k=1,npz
1542  DO j=npy+jstag,jed+jstag
1543  DO i=istart,iend+istag
1544  var_nest(i, j, k) = (var_t0(i, j, k)*(split-step)+step*&
1545 & var_t1(i, j, k))*denom
1546  END DO
1547  END DO
1548  END DO
1549  END IF
1550  END SUBROUTINE nested_grid_bc_apply_intt
1551  SUBROUTINE update_coarse_grid_mpp_2d(var_coarse, var_nest, nest_domain&
1552 & , ind_update, dx, dy, area, isd_p, ied_p, jsd_p, jed_p, is_n, ie_n, &
1553 & js_n, je_n, isu, ieu, jsu, jeu, npx, npy, istag, jstag, r, &
1554 & nestupdate, upoff, nsponge, parent_proc, child_proc, parent_grid)
1555  IMPLICIT NONE
1556  INTEGER, INTENT(IN) :: isd_p, ied_p, jsd_p, jed_p, is_n, ie_n, js_n&
1557 & , je_n
1558  INTEGER, INTENT(IN) :: isu, ieu, jsu, jeu
1559  INTEGER, INTENT(IN) :: istag, jstag, r, nestupdate, upoff, nsponge
1560  INTEGER, INTENT(IN) :: ind_update(isd_p:ied_p+1, jsd_p:jed_p+1, 2)
1561  INTEGER, INTENT(IN) :: npx, npy
1562  REAL, INTENT(IN) :: var_nest(is_n:ie_n+istag, js_n:je_n+jstag)
1563  REAL, INTENT(INOUT) :: var_coarse(isd_p:ied_p+istag, jsd_p:jed_p+&
1564 & jstag)
1565  REAL, INTENT(IN) :: dx(isd:ied, jsd:jed+1)
1566  REAL, INTENT(IN) :: dy(isd:ied+1, jsd:jed)
1567  REAL, INTENT(IN) :: area(isd:ied, jsd:jed)
1568  LOGICAL, INTENT(IN) :: parent_proc, child_proc
1569  TYPE(FV_ATMOS_TYPE), INTENT(INOUT) :: parent_grid
1570  TYPE(NEST_DOMAIN_TYPE), INTENT(INOUT) :: nest_domain
1571  REAL :: var_nest_3d(is_n:ie_n+istag, js_n:je_n+jstag, 1)
1572  REAL :: var_coarse_3d(isd_p:ied_p+istag, jsd_p:jed_p+jstag, 1)
1573  INTRINSIC SIZE
1574  IF (child_proc .AND. SIZE(var_nest) .GT. 1) var_nest_3d(is_n:ie_n+&
1575 & istag, js_n:je_n+jstag, 1) = var_nest(is_n:ie_n+istag, js_n:je_n+&
1576 & jstag)
1577  IF (parent_proc .AND. SIZE(var_coarse) .GT. 1) var_coarse_3d(isd_p:&
1578 & ied_p+istag, jsd_p:jed_p, 1) = var_coarse(isd_p:ied_p+istag, jsd_p&
1579 & :jed_p+jstag)
1580  CALL update_coarse_grid_mpp(var_coarse_3d, var_nest_3d, nest_domain&
1581 & , ind_update, dx, dy, area, isd_p, ied_p, &
1582 & jsd_p, jed_p, is_n, ie_n, js_n, je_n, isu, ieu&
1583 & , jsu, jeu, npx, npy, 1, istag, jstag, r, &
1584 & nestupdate, upoff, nsponge, parent_proc, &
1585 & child_proc, parent_grid)
1586  IF (SIZE(var_coarse) .GT. 1 .AND. parent_proc) var_coarse(isd_p:&
1587 & ied_p+istag, jsd_p:jed_p+jstag) = var_coarse_3d(isd_p:ied_p+istag&
1588 & , jsd_p:jed_p, 1)
1589  END SUBROUTINE update_coarse_grid_mpp_2d
1590  SUBROUTINE update_coarse_grid_mpp(var_coarse, var_nest, nest_domain, &
1591 & ind_update, dx, dy, area, isd_p, ied_p, jsd_p, jed_p, is_n, ie_n, &
1592 & js_n, je_n, isu, ieu, jsu, jeu, npx, npy, npz, istag, jstag, r, &
1593 & nestupdate, upoff, nsponge, parent_proc, child_proc, parent_grid)
1594  IMPLICIT NONE
1595 !This routine assumes the coarse and nested grids are properly
1596 ! aligned, and that in particular for odd refinement ratios all
1597 ! coarse-grid points coincide with nested-grid points
1598  INTEGER, INTENT(IN) :: isd_p, ied_p, jsd_p, jed_p, is_n, ie_n, js_n&
1599 & , je_n
1600  INTEGER, INTENT(IN) :: isu, ieu, jsu, jeu
1601  INTEGER, INTENT(IN) :: istag, jstag, npx, npy, npz, r, nestupdate, &
1602 & upoff, nsponge
1603  INTEGER, INTENT(IN) :: ind_update(isd_p:ied_p+1, jsd_p:jed_p+1, 2)
1604  REAL, INTENT(IN) :: var_nest(is_n:ie_n+istag, js_n:je_n+jstag, npz)
1605  REAL, INTENT(INOUT) :: var_coarse(isd_p:ied_p+istag, jsd_p:jed_p+&
1606 & jstag, npz)
1607  REAL, INTENT(IN) :: area(isd:ied, jsd:jed)
1608  REAL, INTENT(IN) :: dx(isd:ied, jsd:jed+1)
1609  REAL, INTENT(IN) :: dy(isd:ied+1, jsd:jed)
1610  LOGICAL, INTENT(IN) :: parent_proc, child_proc
1611  TYPE(FV_ATMOS_TYPE), INTENT(INOUT) :: parent_grid
1612  TYPE(NEST_DOMAIN_TYPE), INTENT(INOUT) :: nest_domain
1613  INTEGER :: in, jn, ini, jnj, s, qr
1614  INTEGER :: is_c, ie_c, js_c, je_c, is_f, ie_f, js_f, je_f
1615  INTEGER :: istart, istop, jstart, jstop, ishift, jshift, j, i, k
1616  REAL :: val
1617  REAL, DIMENSION(:, :, :), ALLOCATABLE :: nest_dat
1618  REAL :: var_nest_send(is_n:ie_n+istag, js_n:je_n+jstag, npz)
1619  INTEGER :: position
1620  IF (istag .EQ. 1 .AND. jstag .EQ. 1) THEN
1621  position = corner
1622  ELSE IF (istag .EQ. 0 .AND. jstag .EQ. 1) THEN
1623  position = north
1624  ELSE IF (istag .EQ. 1 .AND. jstag .EQ. 0) THEN
1625  position = east
1626  ELSE
1627  position = center
1628  END IF
1629  CALL mpp_get_f2c_index(nest_domain, is_c, ie_c, js_c, je_c, is_f, &
1630 & ie_f, js_f, je_f, position)
1631  IF (ie_f .GT. is_f .AND. je_f .GT. js_f) THEN
1632  ALLOCATE(nest_dat(is_f:ie_f, js_f:je_f, npz))
1633  ELSE
1634  ALLOCATE(nest_dat(1, 1, 1))
1635  END IF
1636  nest_dat = -600.0
1637  IF (child_proc) THEN
1638 !! IF an area average (for istag == jstag == 0) or a linear average then multiply in the areas before sending data
1639  IF (istag .EQ. 0 .AND. jstag .EQ. 0) THEN
1640  SELECT CASE (nestupdate)
1641  CASE (1, 2, 6, 7, 8)
1642 !$NO-MP parallel do default(none) shared(npz,js_n,je_n,is_n,ie_n,var_nest_send,var_nest,area)
1643  DO k=1,npz
1644  DO j=js_n,je_n
1645  DO i=is_n,ie_n
1646  var_nest_send(i, j, k) = var_nest(i, j, k)*area(i, j)
1647  END DO
1648  END DO
1649  END DO
1650  END SELECT
1651  ELSE IF (istag .EQ. 0 .AND. jstag .GT. 0) THEN
1652  SELECT CASE (nestupdate)
1653  CASE (1, 6, 7, 8)
1654 !$NO-MP parallel do default(none) shared(npz,js_n,je_n,is_n,ie_n,var_nest_send,var_nest,dx)
1655  DO k=1,npz
1656  DO j=js_n,je_n+1
1657  DO i=is_n,ie_n
1658  var_nest_send(i, j, k) = var_nest(i, j, k)*dx(i, j)
1659  END DO
1660  END DO
1661  END DO
1662  CASE DEFAULT
1663  CALL mpp_error(fatal, 'nestupdate type not implemented')
1664  END SELECT
1665  ELSE IF (istag .GT. 0 .AND. jstag .EQ. 0) THEN
1666  SELECT CASE (nestupdate)
1667  CASE (1, 6, 7, 8)
1668 !averaging update; in-line average for face-averaged values instead of areal average
1669 !$NO-MP parallel do default(none) shared(npz,js_n,je_n,is_n,ie_n,var_nest_send,var_nest,dy)
1670  DO k=1,npz
1671  DO j=js_n,je_n
1672  DO i=is_n,ie_n+1
1673  var_nest_send(i, j, k) = var_nest(i, j, k)*dy(i, j)
1674  END DO
1675  END DO
1676  END DO
1677  CASE DEFAULT
1678  CALL mpp_error(fatal, 'nestupdate type not implemented')
1679  END SELECT
1680  ELSE
1681  CALL mpp_error(fatal, &
1682 & 'Cannot have both nonzero istag and jstag.')
1683  END IF
1684  END IF
1685  CALL timing_on('COMM_TOTAL')
1686  CALL mpp_update_nest_coarse(var_nest_send, nest_domain, nest_dat, &
1687 & position=position)
1688  CALL timing_off('COMM_TOTAL')
1689 !rounds down (since r > 0)
1690  s = r/2
1691  qr = r*upoff + nsponge - s
1692  IF (parent_proc .AND. (.NOT.(ieu .LT. isu .OR. jeu .LT. jsu))) THEN
1693  IF (istag .EQ. 0 .AND. jstag .EQ. 0) THEN
1694  SELECT CASE (nestupdate)
1695  CASE (1, 2, 6, 7, 8)
1696 ! 1 = Conserving update on all variables; 2 = conserving update for cell-centered values; 6 = conserving remap-update
1697 !$NO-MP parallel do default(none) shared(npz,jsu,jeu,isu,ieu,ind_update,nest_dat,parent_grid,var_coarse,r) &
1698 !$NO-MP private(in,jn,val)
1699  DO k=1,npz
1700  DO j=jsu,jeu
1701  DO i=isu,ieu
1702  in = ind_update(i, j, 1)
1703  jn = ind_update(i, j, 2)
1704 !!$ if (in < max(1+qr,is_f) .or. in > min(npx-1-qr-r+1,ie_f) .or. &
1705 !!$ jn < max(1+qr,js_f) .or. jn > min(npy-1-qr-r+1,je_f)) then
1706 !!$ write(mpp_pe()+3000,'(A, 14I6)') 'SKIP: ', i, j, in, jn, 1+qr, is_f, ie_f, js_f, je_f, npy-1-qr-r+1, isu, ieu,
1707 !jsu, jeu
1708 !!$ cycle
1709 !!$ endif
1710  val = 0.
1711  DO jnj=jn,jn+r-1
1712  DO ini=in,in+r-1
1713  val = val + nest_dat(ini, jnj, k)
1714  END DO
1715  END DO
1716 !var_coarse(i,j,k) = val/r**2.
1717 !!! CLEANUP: Couldn't rarea and rdx and rdy be built into the weight arrays?
1718 !!! Two-way updates do not yet have weights, tho
1719  var_coarse(i, j, k) = val*parent_grid%gridstruct%rarea(i&
1720 & , j)
1721  END DO
1722  END DO
1723  END DO
1724  CASE DEFAULT
1725  CALL mpp_error(fatal, 'nestupdate type not implemented')
1726  END SELECT
1727  ELSE IF (istag .EQ. 0 .AND. jstag .GT. 0) THEN
1728  SELECT CASE (nestupdate)
1729  CASE (1, 6, 7, 8)
1730 !$NO-MP parallel do default(none) shared(npz,jsu,jeu,isu,ieu,ind_update,nest_dat,parent_grid,var_coarse,r) &
1731 !$NO-MP private(in,jn,val)
1732  DO k=1,npz
1733  DO j=jsu,jeu+1
1734  DO i=isu,ieu
1735  in = ind_update(i, j, 1)
1736  jn = ind_update(i, j, 2)
1737 !!$ if (in < max(1+qr,is_f) .or. in > min(npx-1-qr-r+1,ie_f) .or. &
1738 !!$ jn < max(1+qr+s,js_f) .or. jn > min(npy-1-qr-s+1,je_f)) then
1739 !!$ write(mpp_pe()+3000,'(A, 14I)') 'SKIP u: ', i, j, in, jn, 1+qr, is_f, ie_f, js_f, je_f, npy-1-qr-s+1, isu, ieu,
1740 ! jsu, jeu
1741 !!$ cycle
1742 !!$ endif
1743  val = 0.
1744  DO ini=in,in+r-1
1745  val = val + nest_dat(ini, jn, k)
1746  END DO
1747 ! var_coarse(i,j,k) = val/r
1748  var_coarse(i, j, k) = val*parent_grid%gridstruct%rdx(i, &
1749 & j)
1750  END DO
1751  END DO
1752  END DO
1753  CASE DEFAULT
1754  CALL mpp_error(fatal, 'nestupdate type not implemented')
1755  END SELECT
1756  ELSE IF (istag .GT. 0 .AND. jstag .EQ. 0) THEN
1757  SELECT CASE (nestupdate)
1758  CASE (1, 6, 7, 8)
1759 !averaging update; in-line average for face-averaged values instead of areal average
1760 !$NO-MP parallel do default(none) shared(npz,jsu,jeu,isu,ieu,ind_update,nest_dat,parent_grid,var_coarse,r) &
1761 !$NO-MP private(in,jn,val)
1762  DO k=1,npz
1763  DO j=jsu,jeu
1764  DO i=isu,ieu+1
1765  in = ind_update(i, j, 1)
1766  jn = ind_update(i, j, 2)
1767 !!$ if (in < max(1+qr+s,is_f) .or. in > min(npx-1-qr-s+1,ie_f) .or. &
1768 !!$ jn < max(1+qr,js_f) .or. jn > min(npy-1-qr-r+1,je_f)) then
1769 !!$ write(mpp_pe()+3000,'(A, 14I6)') 'SKIP v: ', i, j, in, jn, 1+qr, is_f, ie_f, js_f, je_f, npx-1-qr-s+1, isu, ieu
1770 !, jsu, jeu
1771 !!$ cycle
1772 !!$ endif
1773  val = 0.
1774  DO jnj=jn,jn+r-1
1775  val = val + nest_dat(in, jnj, k)
1776  END DO
1777 ! var_coarse(i,j,k) = val/r
1778  var_coarse(i, j, k) = val*parent_grid%gridstruct%rdy(i, &
1779 & j)
1780  END DO
1781  END DO
1782  END DO
1783  CASE DEFAULT
1784  CALL mpp_error(fatal, 'nestupdate type not implemented')
1785  END SELECT
1786  END IF
1787  END IF
1788  DEALLOCATE(nest_dat)
1789  END SUBROUTINE update_coarse_grid_mpp
1790 
1791 end module boundary_tlm_mod
subroutine, public nested_grid_bc_apply_intt(var_nest, istag, jstag, npx, npy, npz, bd, step, split, bc, bctype)
subroutine nested_grid_bc_2d(var_nest, var_coarse, ind, wt, istag, jstag, npx, npy, bd, isg, ieg, jsg, jeg, nstep_in, nsplit_in)
subroutine, public nested_grid_bc_send(var_coarse, nest_domain, istag, jstag)
Definition: mpp.F90:39
l_size ! loop over number of fields ke do je do ie to je n if(.NOT. d_comm%R_do_buf(list)) cycle from_pe
subroutine update_coarse_grid_mpp_2d(var_coarse, var_nest, nest_domain, ind_update, dx, dy, area, isd_p, ied_p, jsd_p, jed_p, is_n, ie_n, js_n, je_n, isu, ieu, jsu, jeu, npx, npy, istag, jstag, r, nestupdate, upoff, nsponge, parent_proc, child_proc, parent_grid)
integer, parameter, public ng
subroutine timing_on(blk_name)
subroutine nested_grid_bc_mpp_send(var_coarse, nest_domain, istag, jstag)
subroutine fill_nested_grid_2d(var_nest, var_coarse, ind, wt, istag, jstag, isg, ieg, jsg, jeg, bd, istart_in, iend_in, jstart_in, jend_in)
subroutine update_coarse_grid_mpp(var_coarse, var_nest, nest_domain, ind_update, dx, dy, area, isd_p, ied_p, jsd_p, jed_p, is_n, ie_n, js_n, je_n, isu, ieu, jsu, jeu, npx, npy, npz, istag, jstag, r, nestupdate, upoff, nsponge, parent_proc, child_proc, parent_grid)
subroutine fill_nested_grid_3d(var_nest, var_coarse, ind, wt, istag, jstag, isg, ieg, jsg, jeg, npz, bd, istart_in, iend_in, jstart_in, jend_in)
subroutine nested_grid_bc_3d(var_nest, var_coarse, ind, wt, istag, jstag, npx, npy, npz, bd, isg, ieg, jsg, jeg, nstep_in, nsplit_in)
real, parameter, public grav
Acceleration due to gravity [m/s^2].
Definition: constants.F90:76
subroutine, public nested_grid_bc_save_proc(nest_domain, ind, wt, istag, jstag, npx, npy, npz, bd, nest_bc, nest_bc_buffers, pd_in)
#define max(a, b)
Definition: mosaic_util.h:33
subroutine, public nested_grid_bc_apply_intt_tlm(var_nest, var_nest_tl, istag, jstag, npx, npy, npz, bd, step, split, bc, bctype)
subroutine, public extrapolation_bc(q, istag, jstag, npx, npy, bd, pd_in, debug_in)
subroutine nested_grid_bc_mpp(var_nest, var_coarse, nest_domain, ind, wt, istag, jstag, npx, npy, npz, bd, isg, ieg, jsg, jeg, nstep_in, nsplit_in, proc_in)
#define min(a, b)
Definition: mosaic_util.h:32
subroutine, public nested_grid_bc_recv(nest_domain, istag, jstag, npz, bd, nest_bc_buffers)
subroutine nested_grid_bc_2d_mpp(var_nest, var_coarse, nest_domain, ind, wt, istag, jstag, npx, npy, bd, isg, ieg, jsg, jeg, nstep_in, nsplit_in, proc_in)
subroutine timing_off(blk_name)