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