FV3DYCORE  Version 2.0.0
boundary.F90
Go to the documentation of this file.
1 !***********************************************************************
2 !* GNU Lesser General Public License
3 !*
4 !* This file is part of the FV3 dynamical core.
5 !*
6 !* The FV3 dynamical core is free software: you can redistribute it
7 !* and/or modify it under the terms of the
8 !* GNU Lesser General Public License as published by the
9 !* Free Software Foundation, either version 3 of the License, or
10 !* (at your option) any later version.
11 !*
12 !* The FV3 dynamical core is distributed in the hope that it will be
13 !* useful, but WITHOUT ANYWARRANTY; without even the implied warranty
14 !* of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.
15 !* See the GNU General Public License for more details.
16 !*
17 !* You should have received a copy of the GNU Lesser General Public
18 !* License along with the FV3 dynamical core.
19 !* If not, see <http://www.gnu.org/licenses/>.
20 !***********************************************************************
21 
24 
26 
27 ! Modules Included:
28 ! <table>
29 ! <tr>
30 ! <th>Module Name</th>
31 ! <th>Functions Included</th>
32 ! </tr>
33 ! <tr>
34 ! <td>constants_mod</td>
35 ! <td>grav</td>
36 ! <tr>
37 ! <td>fv_arrays_mod</td>
38 ! <td>fv_atmos_type, fv_nest_BC_type_3D, fv_grid_bounds_type</td>
39 ! </tr>
40 ! <tr>
41 ! <td>fv_mp_mod</td>
42 ! <td>ng, isc,jsc,iec,jec, isd,jsd,ied,jed, is,js,ie,je, is_master, mp_bcst</td>
43 ! </tr>
44 ! <tr>
45 ! <td>fv_timing_mod</td>
46 ! <td>timing_on, timing_off</td>
47 ! </tr>
48 ! <tr>
49 ! <td>mpp_mod/td>
50 ! <td>mpp_error, FATAL, mpp_sum, mpp_sync, mpp_npes, mpp_broadcast, WARNING, mpp_pe,
51 ! mpp_send, mpp_recv</td>
52 ! </tr>
53 ! <tr>
54 ! <td>mpp_domains_mod/td>
55 ! <td>mpp_get_compute_domain, mpp_get_data_domain, mpp_get_global_domain,
56 ! ENTER, CORNER, NORTH, EAST,nest_domain_type, WEST, SOUTH,
57 ! mpp_get_C2F_index, mpp_update_nest_fine,mpp_global_field, mpp_get_pelist
58 ! mpp_get_F2C_index, mpp_update_nest_coarse</td>
59 ! </tr>
60 ! </table>
61 
62  use fv_mp_mod, only: is_master
63  use constants_mod, only: grav
64 
65  use mpp_domains_mod, only: mpp_get_compute_domain, mpp_get_data_domain, mpp_get_global_domain
66  use mpp_domains_mod, only: center, corner, north, east
67  use mpp_domains_mod, only: mpp_global_field, mpp_get_pelist
68  use mpp_domains_mod, only: agrid, bgrid_ne, cgrid_ne, dgrid_ne
69  use mpp_mod, only: mpp_error, fatal, mpp_sum, mpp_sync, mpp_npes, mpp_broadcast, warning, mpp_pe
70 
71  use fv_mp_mod, only: mp_bcst
73  use mpp_mod, only: mpp_send, mpp_recv
75  use mpp_domains_mod, only : nest_domain_type, west, south
76  use mpp_domains_mod, only : mpp_get_c2f_index, mpp_update_nest_fine
77  use mpp_domains_mod, only : mpp_get_f2c_index, mpp_update_nest_coarse
78  !use mpp_domains_mod, only : mpp_get_domain_shift
79 
80  implicit none
81  public extrapolation_bc
85 
89  interface nested_grid_bc
90  module procedure nested_grid_bc_2d
91 ! module procedure nested_grid_BC_mpp_2d
92  module procedure nested_grid_bc_mpp_3d
93  module procedure nested_grid_bc_mpp_send_2d
94  module procedure nested_grid_bc_mpp_send_3d
95  module procedure nested_grid_bc_2d_mpp
96  module procedure nested_grid_bc_3d
97  module procedure nested_grid_bc_mpp_3d_vector
98  end interface
99 
101  module procedure nested_grid_bc_send_scalar
102  module procedure nested_grid_bc_send_vector
103  end interface
104 
106  module procedure nested_grid_bc_recv_scalar
107  module procedure nested_grid_bc_recv_vector
108  end interface
112 
114  module procedure fill_nested_grid_2d
115  module procedure fill_nested_grid_3d
116  end interface
117 
123  module procedure update_coarse_grid_mpp
124  module procedure update_coarse_grid_mpp_2d
125  module procedure update_coarse_grid_mpp_vector
126  end interface
127 
128 contains
129 
131 !Not to be confused with extrapolated-in-time nested BCs
132  subroutine extrapolation_bc(q, istag, jstag, npx, npy, bd, pd_in, debug_in)
134  type(fv_grid_bounds_type), intent(IN) :: bd
135  integer, intent(in) :: istag, jstag, npx, npy
136  real, intent(inout), dimension(bd%isd:bd%ied+istag, bd%jsd:bd%jed+jstag) :: q
137  logical, intent(in), OPTIONAL :: pd_in, debug_in
138 
139  integer :: i,j, istart, iend, jstart, jend
140  logical :: pd, debug
141 
142  integer :: is, ie, js, je
143  integer :: isd, ied, jsd, jed
144 
145  is = bd%is
146  ie = bd%ie
147  js = bd%js
148  je = bd%je
149  isd = bd%isd
150  ied = bd%ied
151  jsd = bd%jsd
152  jed = bd%jed
153 
154  istart = max(isd, 1)
155  iend = min(ied,npx-1)
156  jstart = max(jsd, 1)
157  jend = min(jed,npy-1)
158 
159  !Positive-definite extrapolation: shift from linear extrapolation to zero-gradient when the extrapolated value turns negative.
160  if (present(pd_in)) then
161  pd = pd_in
162  else
163  pd = .false.
164  end if
165 
166  if (present(debug_in)) then
167  debug = debug_in
168  else
169  debug = .false.
170  end if
171 
172  if (is == 1) then
173 
174  if (pd) then
175 
176  do j = jstart,jend+jstag
177  do i = 0,isd,-1
178 
179  if (real(i) <= 1. - q(1,j)/(q(2,j) - q(1,j) + 1.e-12) .and. q(1,j) < q(2,j)) then
180  q(i,j) = q(i+1,j)
181  else
182  q(i,j) = real(2-i)*q(1,j) - real(1-i)*q(2,j)
183  end if
184 
185  end do
186  end do
187 
188  else
189 
190  do j = jstart,jend+jstag
191  do i = 0,isd,-1
192 
193  q(i,j) = real(2-i)*q(1,j) - real(1-i)*q(2,j)
194 
195  end do
196  end do
197 
198  end if
199 
200  end if
201 
202  if (js == 1) then
203 
204  if (pd) then
205 
206  do j = 0,jsd,-1
207  do i = istart,iend+istag
208 
209  if (real(j) <= 1. - q(i,1)/(q(i,2) - q(i,1) + 1.e-12) .and. q(i,1) < q(i,2)) then
210  q(i,j) = q(i,j+1)
211  else
212  q(i,j) = real(2-j)*q(i,1) - real(1-j)*q(i,2)
213  end if
214 
215  end do
216  end do
217 
218  else
219 
220  do j = 0,jsd,-1
221  do i = istart,iend+istag
222 
223  q(i,j) = real(2-j)*q(i,1) - real(1-j)*q(i,2)
224 
225  end do
226  end do
227 
228  end if
229 
230  end if
231 
232  if (ie == npx - 1) then
233 
234  if (pd) then
235 
236  do j=jstart,jend+jstag
237  do i=ie+1+istag,ied+istag
238 
239  if (real(i) >= ie+istag + q(ie+istag,j)/(q(ie+istag-1,j)-q(ie+istag,j)+1.e-12) .and. &
240  q(ie+istag,j) < q(ie+istag-1,j)) then
241  q(i,j) = q(i-1,j)
242  else
243  q(i,j) = real(i - (ie+istag-1))*q(ie+istag,j) + real((ie+istag) - i)*q(ie+istag-1,j)
244  end if
245 
246  end do
247  end do
248 
249  else
250 
251  do j=jstart,jend+jstag
252  do i=ie+1+istag,ied+istag
253 
254  q(i,j) = real(i - (ie+istag-1))*q(ie+istag,j) + real((ie+istag) - i)*q(ie+istag-1,j)
255 
256  end do
257  end do
258 
259  end if
260 
261  end if
262 
263  if (je == npy - 1) then
264 
265  if (pd) then
266 
267  do j=je+1+jstag,jed+jstag
268  do i=istart,iend+istag
269 
270  if (real(j) >= je+jstag + q(i,je+jstag)/(q(i,je+jstag-1)-q(i,je+jstag)+1.e-12) .and. &
271  q(i,je+jstag-1) > q(i,je+jstag)) then
272  q(i,j) = q(i,j-1)
273  else
274  q(i,j) = real(j - (je+jstag-1))*q(i,je+jstag) + real((je+jstag) - j)*q(i,je+jstag-1)
275  end if
276 
277  end do
278  end do
279 
280  else
281 
282  do j=je+1+jstag,jed+jstag
283  do i=istart,iend+istag
284 
285  q(i,j) = real(j - (je+jstag-1))*q(i,je+jstag) + real((je+jstag) - j)*q(i,je+jstag-1)
286 
287  end do
288  end do
289 
290  end if
291 
292  end if
293 
294 
295  !CORNERS: Average of extrapolations
296 
297  if (is == 1 .and. js == 1) then
298 
299  if (pd) then
300 
301  do j=0,jsd,-1
302  do i=0,isd,-1
303 
304  if (real(i) <= 1. - q(1,j)/(q(2,j) - q(1,j) + 1.e-12) .and. q(2,j) > q(1,j)) then
305  q(i,j) = 0.5*q(i+1,j)
306  else
307  q(i,j) = 0.5*( real(2-i)*q(1,j) - real(1-i)*q(2,j) )
308  end if
309 
310  if (real(j) <= 1. - q(i,1)/(q(i,2) - q(i,1) + 1.e-12) .and. q(i,2) > q(i,1)) then
311  q(i,j) = q(i,j) + 0.5*q(i,j+1)
312 
313  else
314  q(i,j) = q(i,j) + 0.5*(real(2-j)*q(i,1) - real(1-j)*q(i,2))
315  end if
316 
317  end do
318  end do
319 
320  else
321 
322  do j=jsd,0
323  do i=isd,0
324 
325  q(i,j) = 0.5*( real(2-i)*q(1,j) - real(1-i)*q(2,j) ) + &
326  0.5*( real(2-j)*q(i,1) - real(1-j)*q(i,2) )
327 
328  end do
329  end do
330 
331  end if
332 
333  end if
334 
335  if (is == 1 .and. je == npy-1) then
336 
337  if (pd) then
338 
339  do j=je+1+jstag,jed+jstag
340  do i=0,isd,-1
341 
342  if (real(i) <= 1. - q(1,j)/(q(2,j) - q(1,j) + 1.e-12) .and. q(2,j) > q(1,j)) then
343  q(i,j) = 0.5*q(i+1,j)
344  else
345  q(i,j) = 0.5*( real(2-i)*q(1,j) - real(1-i)*q(2,j) )
346  end if
347 
348  !'Unary plus' removed to appease IBM compiler
349  !if (real(j) >= je+jstag - q(i,je+jstag)/(q(i,je+jstag-1)-q(i,je+jstag)+1.e-12) .and. &
350  if (real(j) >= je+jstag - q(i,je+jstag)/(q(i,je+jstag-1)-q(i,je+jstag)+1.e-12) .and. &
351  q(i,je+jstag-1) > q(i,je+jstag) ) then
352  q(i,j) = q(i,j) + 0.5*q(i,j-1)
353  else
354  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) )
355  end if
356 
357  end do
358  end do
359 
360  else
361 
362  do j=je+1+jstag,jed+jstag
363  do i=isd,0
364 
365  q(i,j) = 0.5*( real(2-i)*q(1,j) - real(1-i)*q(2,j) ) + &
366  0.5*( real(j - (je+jstag-1))*q(i,je+jstag) + real((je+jstag) - j)*q(i,je+jstag-1) )
367 
368  end do
369  end do
370 
371  end if
372 
373  end if
374 
375  if (ie == npx-1 .and. je == npy-1) then
376 
377  if (pd) then
378 
379  do j=je+1+jstag,jed+jstag
380  do i=ie+1+istag,ied+istag
381 
382 
383  if (real(i) >= ie+istag + q(ie+istag,j)/(q(ie+istag-1,j)-q(ie+istag,j)+1.e-12) .and. &
384  q(ie+istag-1,j) > q(ie+istag,j)) then
385  q(i,j) = 0.5*q(i-1,j)
386  else
387  q(i,j) = 0.5*(real(i - (ie+istag-1))*q(ie+istag,j) + real((ie+istag) - i)*q(ie+istag-1,j))
388  end if
389 
390  if (real(j) >= je+jstag + q(i,je+jstag)/(q(i,je+jstag-1)-q(i,je+jstag)+1.e-12) .and. &
391  q(i,je+jstag-1) > q(i,je+jstag)) then
392  q(i,j) = q(i,j) + 0.5*q(i,j-1)
393  else
394  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) )
395  end if
396 
397  end do
398  end do
399 
400  else
401 
402  do j=je+1+jstag,jed+jstag
403  do i=ie+1+istag,ied+istag
404 
405  q(i,j) = 0.5*( real(i - (ie+istag-1))*q(ie+istag,j) + real((ie+istag) - i)*q(ie+istag-1,j) ) + &
406  0.5*( real(j - (je+jstag-1))*q(i,je+jstag) + real((je+jstag) - j)*q(i,je+jstag-1) )
407 
408  end do
409  end do
410 
411  end if
412 
413  end if
414 
415  if (ie == npx-1 .and. js == 1) then
416 
417  if (pd) then
418 
419  do j=0,jsd,-1
420  do i=ie+1+istag,ied+istag
421 
422 
423  if (real(i) >= ie+istag + q(ie+istag,j)/(q(ie+istag-1,j)-q(ie+istag,j)+1.e-12) .and. &
424  q(ie+istag-1,j) > q(ie+istag,j)) then
425  q(i,j) = 0.5*q(i-1,j)
426  else
427  q(i,j) = 0.5*(real(i - (ie+istag-1))*q(ie+istag,j) + real((ie+istag) - i)*q(ie+istag-1,j))
428  end if
429 
430  if (real(j) <= 1. - q(i,1)/(q(i,2) - q(i,1) + 1.e-12) .and. &
431  q(i,2) > q(i,1)) then
432  q(i,j) = q(i,j) + 0.5*q(i,j+1)
433  else
434  q(i,j) = q(i,j) + 0.5*(real(2-j)*q(i,1) - real(1-j)*q(i,2))
435  end if
436 
437  end do
438  end do
439 
440 
441  else
442 
443  do j=jsd,0
444  do i=ie+1+istag,ied+istag
445 
446  q(i,j) = 0.5*( real(i - (ie+istag-1))*q(ie+istag,j) + real((ie+istag) - i)*q(ie+istag-1,j) ) + &
447  0.5*( real(2-j)*q(i,1) - real(1-j)*q(i,2) )
448 
449  end do
450  end do
451 
452  end if
453 
454  end if
455 
456 
457  end subroutine extrapolation_bc
458 
459  subroutine fill_nested_grid_2d(var_nest, var_coarse, ind, wt, istag, jstag, &
460  isg, ieg, jsg, jeg, bd, istart_in, iend_in, jstart_in, jend_in)
462  type(fv_grid_bounds_type), intent(IN) :: bd
463  real, dimension(bd%isd:bd%ied+istag,bd%jsd:bd%jed+jstag), intent(INOUT) :: var_nest
464  real, dimension(isg:ieg+istag,jsg:jeg+jstag), intent(IN) :: var_coarse
465  integer, dimension(bd%isd:bd%ied+istag,bd%jsd:bd%jed+jstag,2), intent(IN) :: ind
466  real, dimension(bd%isd:bd%ied+istag,bd%jsd:bd%jed+jstag,4), intent(IN) :: wt
467  integer, intent(IN) :: istag, jstag, isg, ieg, jsg, jeg
468  integer, intent(IN), OPTIONAL :: istart_in, iend_in, jstart_in, jend_in
469 
470  integer :: i,j, ic, jc
471  integer :: istart, iend, jstart, jend
472 
473  integer :: is, ie, js, je
474  integer :: isd, ied, jsd, jed
475 
476  is = bd%is
477  ie = bd%ie
478  js = bd%js
479  je = bd%je
480  isd = bd%isd
481  ied = bd%ied
482  jsd = bd%jsd
483  jed = bd%jed
484 
485  if (present(istart_in)) then
486  istart = istart_in
487  else
488  istart = isd
489  end if
490  if (present(iend_in)) then
491  iend = iend_in+istag
492  else
493  iend = ied+istag
494  end if
495 
496  if (present(jstart_in)) then
497  jstart = jstart_in
498  else
499  jstart = jsd
500  end if
501  if (present(jend_in)) then
502  jend = jend_in+jstag
503  else
504  jend = jed+jstag
505  end if
506 
507  do j=jstart,jend
508  do i=istart,iend
509 
510  ic = ind(i,j,1)
511  jc = ind(i,j,2)
512 
513  var_nest(i,j) = &
514  wt(i,j,1)*var_coarse(ic, jc) + &
515  wt(i,j,2)*var_coarse(ic, jc+1) + &
516  wt(i,j,3)*var_coarse(ic+1,jc+1) + &
517  wt(i,j,4)*var_coarse(ic+1,jc)
518 
519  end do
520  end do
521 
522  end subroutine fill_nested_grid_2d
523 
524  subroutine fill_nested_grid_3d(var_nest, var_coarse, ind, wt, istag, jstag, &
525  isg, ieg, jsg, jeg, npz, bd, istart_in, iend_in, jstart_in, jend_in)
527  type(fv_grid_bounds_type), intent(IN) :: bd
528  real, dimension(bd%isd:bd%ied+istag,bd%jsd:bd%jed+jstag,npz), intent(INOUT) :: var_nest
529  real, dimension(isg:ieg+istag,jsg:jeg+jstag,npz), intent(IN) :: var_coarse
530  integer, dimension(bd%isd:bd%ied+istag,bd%jsd:bd%jed+jstag,2), intent(IN) :: ind
531  real, dimension(bd%isd:bd%ied+istag,bd%jsd:bd%jed+jstag,4), intent(IN) :: wt
532  integer, intent(IN) :: istag, jstag, isg, ieg, jsg, jeg, npz
533  integer, intent(IN), OPTIONAL :: istart_in, iend_in, jstart_in, jend_in
534 
535  integer :: i,j, ic, jc, k
536  integer :: istart, iend, jstart, jend
537 
538  integer :: is, ie, js, je
539  integer :: isd, ied, jsd, jed
540 
541  is = bd%is
542  ie = bd%ie
543  js = bd%js
544  je = bd%je
545  isd = bd%isd
546  ied = bd%ied
547  jsd = bd%jsd
548  jed = bd%jed
549 
550  if (present(istart_in)) then
551  istart = istart_in
552  else
553  istart = isd
554  end if
555  if (present(iend_in)) then
556  iend = iend_in+istag
557  else
558  iend = ied+istag
559  end if
560 
561  if (present(jstart_in)) then
562  jstart = jstart_in
563  else
564  jstart = jsd
565  end if
566  if (present(jend_in)) then
567  jend = jend_in+jstag
568  else
569  jend = jed+jstag
570  end if
571 
572  do k=1,npz
573 
574  do j=jstart,jend
575  do i=istart,iend
576 
577  ic = ind(i,j,1)
578  jc = ind(i,j,2)
579 
580  var_nest(i,j,k) = &
581  wt(i,j,1)*var_coarse(ic, jc, k) + &
582  wt(i,j,2)*var_coarse(ic, jc+1,k) + &
583  wt(i,j,3)*var_coarse(ic+1,jc+1,k) + &
584  wt(i,j,4)*var_coarse(ic+1,jc, k)
585 
586  end do
587  end do
588 
589  end do
590 
591  end subroutine fill_nested_grid_3d
592 
593 !!$ subroutine nested_grid_BC_mpp_2d(var_nest, nest_domain, ind, wt, istag, jstag, &
594 !!$ npx, npy, bd, isg, ieg, jsg, jeg, nstep_in, nsplit_in, proc_in)
595 !!$
596 !!$ type(fv_grid_bounds_type), intent(IN) :: bd
597 !!$ real, dimension(bd%isd:bd%ied+istag,bd%jsd:bd%jed+jstag), intent(INOUT) :: var_nest
598 !!$ real, dimension(isg:ieg+istag,jsg:jeg+jstag), intent(IN) :: var_coarse
599 !!$ type(nest_domain_type), intent(INOUT) :: nest_domain
600 !!$ integer, dimension(bd%isd:bd%ied+istag,bd%jsd:bd%jed+jstag,2), intent(IN) :: ind
601 !!$ real, dimension(bd%isd:bd%ied+istag,bd%jsd:bd%jed+jstag,4), intent(IN) :: wt
602 !!$ integer, intent(IN) :: istag, jstag, npx, npy, isg, ieg, jsg, jeg
603 !!$ integer, intent(IN), OPTIONAL :: nstep_in, nsplit_in
604 !!$ logical, intent(IN), OPTIONAL :: proc_in
605 !!$
606 !!$ real, dimension(bd%isd:bd%ied+istag,bd%jsd:bd%jed+jstag,1) :: var_nest_3d
607 !!$
608 !!$ integer :: i,j
609 !!$
610 !!$ do j=bd%jsd,bd%jed+jstag
611 !!$ do i=bd%isd,bd%ied+istag
612 !!$ var_nest_3d(i,j,1) = var_nest(i,j)
613 !!$ enddo
614 !!$ enddo
615 !!$
616 !!$ call nested_grid_BC_mpp_3d(var_nest_3d, nest_domain, ind, wt, istag, jstag, &
617 !!$ npx, npy, 1, bd, isg, ieg, jsg, jeg, nstep_in, nsplit_in, proc_in)
618 !!$
619 !!$
620 !!$ end subroutine nested_grid_BC_mpp_2d
621 
622  subroutine nested_grid_bc_mpp_3d(var_nest, var_coarse, nest_domain, ind, wt, istag, jstag, &
623  npx, npy, npz, bd, isg, ieg, jsg, jeg, nest_level, nstep_in, nsplit_in, proc_in)
625  type(fv_grid_bounds_type), intent(IN) :: bd
626  real, dimension(bd%isd:bd%ied+istag,bd%jsd:bd%jed+jstag,npz), intent(INOUT) :: var_nest
627  real, dimension(isg:ieg+istag,jsg:jeg+jstag,npz), intent(IN) :: var_coarse
628  type(nest_domain_type), intent(INOUT) :: nest_domain
629  integer, dimension(bd%isd:bd%ied+istag,bd%jsd:bd%jed+jstag,2), intent(IN) :: ind
630  real, dimension(bd%isd:bd%ied+istag,bd%jsd:bd%jed+jstag,4), intent(IN) :: wt
631  integer, intent(IN) :: istag, jstag, npx, npy, npz, isg, ieg, jsg, jeg
632  integer, intent(IN) :: nest_level
633  integer, intent(IN), OPTIONAL :: nstep_in, nsplit_in
634  logical, intent(IN), OPTIONAL :: proc_in
635 
636  integer :: isw_f, iew_f, jsw_f, jew_f, isw_c, iew_c, jsw_c, jew_c
637  integer :: ise_f, iee_f, jse_f, jee_f, ise_c, iee_c, jse_c, jee_c
638  integer :: iss_f, ies_f, jss_f, jes_f, iss_c, ies_c, jss_c, jes_c
639  integer :: isn_f, ien_f, jsn_f, jen_f, isn_c, ien_c, jsn_c, jen_c
640  real, allocatable :: wbuffer(:,:,:)
641  real, allocatable :: ebuffer(:,:,:)
642  real, allocatable :: sbuffer(:,:,:)
643  real, allocatable :: nbuffer(:,:,:)
644 
645  integer :: i,j, ic, jc, istart, iend, k
646 
647  integer :: position
648  logical :: process
649 
650  integer :: is, ie, js, je
651  integer :: isd, ied, jsd, jed
652 
653  is = bd%is
654  ie = bd%ie
655  js = bd%js
656  je = bd%je
657  isd = bd%isd
658  ied = bd%ied
659  jsd = bd%jsd
660  jed = bd%jed
661 
662  if (PRESENT(proc_in)) then
663  process = proc_in
664  else
665  process = .true.
666  endif
667 
668  if (istag == 1 .and. jstag == 1) then
669  position = corner
670  else if (istag == 0 .and. jstag == 1) then
671  position = north
672  else if (istag == 1 .and. jstag == 0) then
673  position = east
674  else
675  position = center
676  end if
677 
678  call mpp_get_c2f_index(nest_domain, isw_f, iew_f, jsw_f, jew_f, isw_c, iew_c, jsw_c, jew_c, &
679  west, nest_level=nest_level, position=position)
680  call mpp_get_c2f_index(nest_domain, ise_f, iee_f, jse_f, jee_f, ise_c, iee_c, jse_c, jee_c, &
681  east, nest_level=nest_level, position=position)
682  call mpp_get_c2f_index(nest_domain, iss_f, ies_f, jss_f, jes_f, iss_c, ies_c, jss_c, jes_c, &
683  south, nest_level=nest_level, position=position)
684  call mpp_get_c2f_index(nest_domain, isn_f, ien_f, jsn_f, jen_f, isn_c, ien_c, jsn_c, jen_c, &
685  north, nest_level=nest_level, position=position)
686 
687  if( iew_c .GE. isw_c .AND. jew_c .GE. jsw_c ) then
688  allocate(wbuffer(isw_c:iew_c, jsw_c:jew_c,npz))
689  else
690  allocate(wbuffer(1,1,1))
691  endif
692  wbuffer = 0
693 
694  if( iee_c .GE. ise_c .AND. jee_c .GE. jse_c ) then
695  allocate(ebuffer(ise_c:iee_c, jse_c:jee_c,npz))
696  else
697  allocate(ebuffer(1,1,1))
698  endif
699  ebuffer = 0
700 
701  if( ies_c .GE. iss_c .AND. jes_c .GE. jss_c ) then
702  allocate(sbuffer(iss_c:ies_c, jss_c:jes_c,npz))
703  else
704  allocate(sbuffer(1,1,1))
705  endif
706  sbuffer = 0
707 
708  if( ien_c .GE. isn_c .AND. jen_c .GE. jsn_c ) then
709  allocate(nbuffer(isn_c:ien_c, jsn_c:jen_c,npz))
710  else
711  allocate(nbuffer(1,1,1))
712  endif
713  nbuffer = 0
714 
715 
716  call timing_on ('COMM_TOTAL')
717  call mpp_update_nest_fine(var_coarse, nest_domain, wbuffer, sbuffer, ebuffer, nbuffer, &
718  nest_level=nest_level, position=position)
719  call timing_off('COMM_TOTAL')
720 
721  if (process) then
722 
723  if (is == 1) then
724 !OMP parallel do default(none) shared(npz,jsd,jed,jstag,isd,ind,var_nest,wt,wbuffer) private(ic,jc)
725  do k=1,npz
726  do j=jsd,jed+jstag
727  do i=isd,0
728 
729  ic = ind(i,j,1)
730  jc = ind(i,j,2)
731 
732  var_nest(i,j,k) = &
733  wt(i,j,1)*wbuffer(ic, jc, k) + &
734  wt(i,j,2)*wbuffer(ic, jc+1,k) + &
735  wt(i,j,3)*wbuffer(ic+1,jc+1,k) + &
736  wt(i,j,4)*wbuffer(ic+1,jc, k)
737 
738  end do
739  end do
740  end do
741  end if
742 
743  if (js == 1) then
744 
745  if (is == 1) then
746  istart = is
747  else
748  istart = isd
749  end if
750 
751  if (ie == npx-1) then
752  iend = ie
753  else
754  iend = ied
755  end if
756 
757 !OMP parallel do default(none) shared(npz,jsd,istart,iend,istag,ind,var_nest,wt,sbuffer) private(ic,jc)
758  do k=1,npz
759  do j=jsd,0
760  do i=istart,iend+istag
761 
762  ic = ind(i,j,1)
763  jc = ind(i,j,2)
764 
765  var_nest(i,j,k) = &
766  wt(i,j,1)*sbuffer(ic, jc, k) + &
767  wt(i,j,2)*sbuffer(ic, jc+1,k) + &
768  wt(i,j,3)*sbuffer(ic+1,jc+1,k) + &
769  wt(i,j,4)*sbuffer(ic+1,jc, k)
770 
771  end do
772  end do
773  end do
774  end if
775 
776 
777  if (ie == npx-1) then
778 !OMP parallel do default(none) shared(npz,jsd,jed,jstag,npx,ied,istag,ind,var_nest,wt,ebuffer) private(ic,jc)
779  do k=1,npz
780  do j=jsd,jed+jstag
781  do i=npx+istag,ied+istag
782 
783  ic = ind(i,j,1)
784  jc = ind(i,j,2)
785 
786  var_nest(i,j,k) = &
787  wt(i,j,1)*ebuffer(ic, jc, k) + &
788  wt(i,j,2)*ebuffer(ic, jc+1,k) + &
789  wt(i,j,3)*ebuffer(ic+1,jc+1,k) + &
790  wt(i,j,4)*ebuffer(ic+1,jc, k)
791 
792  end do
793  end do
794  end do
795  end if
796 
797  if (je == npy-1) then
798 
799  if (is == 1) then
800  istart = is
801  else
802  istart = isd
803  end if
804 
805  if (ie == npx-1) then
806  iend = ie
807  else
808  iend = ied
809  end if
810 
811 !OMP parallel do default(none) shared(npz,jstag,npy,jed,istart,iend,istag,ind,var_nest,wt,nbuffer) private(ic,jc)
812  do k=1,npz
813  do j=npy+jstag,jed+jstag
814  do i=istart,iend+istag
815 
816  ic = ind(i,j,1)
817  jc = ind(i,j,2)
818 
819  var_nest(i,j,k) = &
820  wt(i,j,1)*nbuffer(ic, jc, k) + &
821  wt(i,j,2)*nbuffer(ic, jc+1,k) + &
822  wt(i,j,3)*nbuffer(ic+1,jc+1,k) + &
823  wt(i,j,4)*nbuffer(ic+1,jc, k)
824 
825  end do
826  end do
827  end do
828  end if
829 
830  endif !process
831 
832  deallocate(wbuffer, ebuffer, sbuffer, nbuffer)
833 
834  end subroutine nested_grid_bc_mpp_3d
835 
836  subroutine get_vector_position(position_x, position_y, gridtype)
837  integer, intent(OUT) :: position_x, position_y
838  integer, optional, intent(IN) :: gridtype
839 
840  integer :: grid_offset_type
841 
842  grid_offset_type = agrid
843  if(present(gridtype)) grid_offset_type = gridtype
844 
845  select case(grid_offset_type)
846  case (agrid)
847  position_x = center
848  position_y = center
849  case (bgrid_ne)
850  position_x = corner
851  position_y = corner
852  case (cgrid_ne)
853  position_x = east
854  position_y = north
855  case (dgrid_ne)
856  position_y = east
857  position_x = north
858  case default
859  call mpp_error(fatal, "get_vector_position: invalid value of gridtype")
860  end select
861 
862 
863  end subroutine get_vector_position
864 
865  subroutine init_buffer(nest_domain, wbuffer, sbuffer, ebuffer, nbuffer, npz, nest_level, position)
866  type(nest_domain_type), intent(INOUT) :: nest_domain
867  real, allocatable, dimension(:,:,:), intent(OUT) :: wbuffer, sbuffer, ebuffer, nbuffer
868  integer, intent(IN) :: npz, position, nest_level
869  integer :: isw_f, iew_f, jsw_f, jew_f, isw_c, iew_c, jsw_c, jew_c
870  integer :: ise_f, iee_f, jse_f, jee_f, ise_c, iee_c, jse_c, jee_c
871  integer :: iss_f, ies_f, jss_f, jes_f, iss_c, ies_c, jss_c, jes_c
872  integer :: isn_f, ien_f, jsn_f, jen_f, isn_c, ien_c, jsn_c, jen_c
873 
874  call mpp_get_c2f_index(nest_domain, isw_f, iew_f, jsw_f, jew_f, isw_c, iew_c, jsw_c, jew_c, &
875  west, nest_level=nest_level, position=position)
876  call mpp_get_c2f_index(nest_domain, ise_f, iee_f, jse_f, jee_f, ise_c, iee_c, jse_c, jee_c, &
877  east, nest_level=nest_level, position=position)
878  call mpp_get_c2f_index(nest_domain, iss_f, ies_f, jss_f, jes_f, iss_c, ies_c, jss_c, jes_c, &
879  south, nest_level=nest_level, position=position)
880  call mpp_get_c2f_index(nest_domain, isn_f, ien_f, jsn_f, jen_f, isn_c, ien_c, jsn_c, jen_c, &
881  north, nest_level=nest_level, position=position)
882 
883  if( iew_c .GE. isw_c .AND. jew_c .GE. jsw_c ) then
884  allocate(wbuffer(isw_c:iew_c, jsw_c:jew_c,npz))
885  else
886  allocate(wbuffer(1,1,1))
887  endif
888  wbuffer = 0
889 
890  if( iee_c .GE. ise_c .AND. jee_c .GE. jse_c ) then
891  allocate(ebuffer(ise_c:iee_c, jse_c:jee_c,npz))
892  else
893  allocate(ebuffer(1,1,1))
894  endif
895  ebuffer = 0
896 
897  if( ies_c .GE. iss_c .AND. jes_c .GE. jss_c ) then
898  allocate(sbuffer(iss_c:ies_c, jss_c:jes_c,npz))
899  else
900  allocate(sbuffer(1,1,1))
901  endif
902  sbuffer = 0
903 
904  if( ien_c .GE. isn_c .AND. jen_c .GE. jsn_c ) then
905  allocate(nbuffer(isn_c:ien_c, jsn_c:jen_c,npz))
906  else
907  allocate(nbuffer(1,1,1))
908  endif
909  nbuffer = 0
910 
911  end subroutine init_buffer
912 
913 
914  subroutine nested_grid_bc_mpp_3d_vector(u_nest, v_nest, u_coarse, v_coarse, nest_domain, ind_u, ind_v, wt_u, wt_v, &
915  istag_u, jstag_u, istag_v, jstag_v, npx, npy, npz, bd, isg, ieg, jsg, jeg, nest_level, nstep_in, nsplit_in, proc_in, &
916  flags, gridtype)
918  type(fv_grid_bounds_type), intent(IN) :: bd
919  integer, intent(IN) :: istag_u, jstag_u, istag_v, jstag_v, npx, npy, npz, isg, ieg, jsg, jeg
920  real, dimension(bd%isd:bd%ied+istag_u,bd%jsd:bd%jed+jstag_u,npz), intent(INOUT) :: u_nest
921  real, dimension(bd%isd:bd%ied+istag_v,bd%jsd:bd%jed+jstag_v,npz), intent(INOUT) :: v_nest
922  real, dimension(isg:ieg+istag_u,jsg:jeg+jstag_u,npz), intent(IN) :: u_coarse
923  real, dimension(isg:ieg+istag_v,jsg:jeg+jstag_v,npz), intent(IN) :: v_coarse
924  type(nest_domain_type), intent(INOUT) :: nest_domain
925  integer, dimension(bd%isd:bd%ied+istag_u,bd%jsd:bd%jed+jstag_u,2), intent(IN) :: ind_u
926  integer, dimension(bd%isd:bd%ied+istag_v,bd%jsd:bd%jed+jstag_v,2), intent(IN) :: ind_v
927  real, dimension(bd%isd:bd%ied+istag_u,bd%jsd:bd%jed+jstag_u,4), intent(IN) :: wt_u
928  real, dimension(bd%isd:bd%ied+istag_v,bd%jsd:bd%jed+jstag_v,4), intent(IN) :: wt_v
929  integer, intent(IN) :: nest_level
930  integer, intent(IN), OPTIONAL :: nstep_in, nsplit_in
931  logical, intent(IN), OPTIONAL :: proc_in
932  integer, intent(IN), OPTIONAL :: flags, gridtype
933 
934  real, allocatable :: wbufferx(:,:,:), wbuffery(:,:,:)
935  real, allocatable :: ebufferx(:,:,:), ebuffery(:,:,:)
936  real, allocatable :: sbufferx(:,:,:), sbuffery(:,:,:)
937  real, allocatable :: nbufferx(:,:,:), nbuffery(:,:,:)
938 
939  integer :: i,j, ic, jc, istart, iend, k
940 
941  integer :: position_x, position_y
942  logical :: process
943 
944  integer :: is, ie, js, je
945  integer :: isd, ied, jsd, jed
946 
947  is = bd%is
948  ie = bd%ie
949  js = bd%js
950  je = bd%je
951  isd = bd%isd
952  ied = bd%ied
953  jsd = bd%jsd
954  jed = bd%jed
955 
956  if (PRESENT(proc_in)) then
957  process = proc_in
958  else
959  process = .true.
960  endif
961 
962  call get_vector_position(position_x, position_y, gridtype)
963  call init_buffer(nest_domain, wbufferx, sbufferx, ebufferx, nbufferx, npz, nest_level, position_x)
964  call init_buffer(nest_domain, wbuffery, sbuffery, ebuffery, nbuffery, npz, nest_level, position_x)
965 
966  call timing_on ('COMM_TOTAL')
967  call mpp_update_nest_fine(u_coarse, v_coarse, nest_domain, wbufferx, wbuffery, sbufferx, sbuffery, &
968  ebufferx, ebuffery, nbufferx, nbuffery, flags=flags, nest_level=nest_level, gridtype=gridtype)
969  call timing_off('COMM_TOTAL')
970 
971  if (process) then
972 
973  if (is == 1) then
974 !OMP parallel do default(none) shared(npz,jsd,jed,jstag,isd,ind,var_nest,wt,wbuffer) private(ic,jc)
975  do k=1,npz
976  do j=jsd,jed+jstag_u
977  do i=isd,0
978 
979  ic = ind_u(i,j,1)
980  jc = ind_u(i,j,2)
981 
982  u_nest(i,j,k) = &
983  wt_u(i,j,1)*wbufferx(ic, jc, k) + &
984  wt_u(i,j,2)*wbufferx(ic, jc+1,k) + &
985  wt_u(i,j,3)*wbufferx(ic+1,jc+1,k) + &
986  wt_u(i,j,4)*wbufferx(ic+1,jc, k)
987 
988  end do
989  end do
990  do j=jsd,jed+jstag_v
991  do i=isd,0
992 
993  ic = ind_v(i,j,1)
994  jc = ind_v(i,j,2)
995 
996  v_nest(i,j,k) = &
997  wt_v(i,j,1)*wbuffery(ic, jc, k) + &
998  wt_v(i,j,2)*wbuffery(ic, jc+1,k) + &
999  wt_v(i,j,3)*wbuffery(ic+1,jc+1,k) + &
1000  wt_v(i,j,4)*wbuffery(ic+1,jc, k)
1001 
1002  end do
1003  end do
1004  end do
1005 
1006  end if
1007 
1008  if (js == 1) then
1009 
1010  if (is == 1) then
1011  istart = is
1012  else
1013  istart = isd
1014  end if
1015 
1016  if (ie == npx-1) then
1017  iend = ie
1018  else
1019  iend = ied
1020  end if
1021 
1022 !OMP parallel do default(none) shared(npz,jsd,istart,iend,istag,ind,var_nest,wt,sbuffer) private(ic,jc)
1023  do k=1,npz
1024  do j=jsd,0
1025  do i=istart,iend+istag_u
1026 
1027  ic = ind_u(i,j,1)
1028  jc = ind_u(i,j,2)
1029 
1030  u_nest(i,j,k) = &
1031  wt_u(i,j,1)*sbufferx(ic, jc, k) + &
1032  wt_u(i,j,2)*sbufferx(ic, jc+1,k) + &
1033  wt_u(i,j,3)*sbufferx(ic+1,jc+1,k) + &
1034  wt_u(i,j,4)*sbufferx(ic+1,jc, k)
1035 
1036  end do
1037  end do
1038  do j=jsd,0
1039  do i=istart,iend+istag_v
1040 
1041  ic = ind_v(i,j,1)
1042  jc = ind_v(i,j,2)
1043 
1044  v_nest(i,j,k) = &
1045  wt_v(i,j,1)*sbuffery(ic, jc, k) + &
1046  wt_v(i,j,2)*sbuffery(ic, jc+1,k) + &
1047  wt_v(i,j,3)*sbuffery(ic+1,jc+1,k) + &
1048  wt_v(i,j,4)*sbuffery(ic+1,jc, k)
1049 
1050  end do
1051  end do
1052  end do
1053  end if
1054 
1055 
1056  if (ie == npx-1) then
1057 !OMP parallel do default(none) shared(npz,jsd,jed,jstag,npx,ied,istag,ind,var_nest,wt,ebuffer) private(ic,jc)
1058  do k=1,npz
1059  do j=jsd,jed+jstag_u
1060  do i=npx+istag_u,ied+istag_u
1061 
1062  ic = ind_u(i,j,1)
1063  jc = ind_u(i,j,2)
1064 
1065  u_nest(i,j,k) = &
1066  wt_u(i,j,1)*ebufferx(ic, jc, k) + &
1067  wt_u(i,j,2)*ebufferx(ic, jc+1,k) + &
1068  wt_u(i,j,3)*ebufferx(ic+1,jc+1,k) + &
1069  wt_u(i,j,4)*ebufferx(ic+1,jc, k)
1070 
1071  end do
1072  end do
1073  do j=jsd,jed+jstag_v
1074  do i=npx+istag_v,ied+istag_v
1075 
1076  ic = ind_v(i,j,1)
1077  jc = ind_v(i,j,2)
1078 
1079  v_nest(i,j,k) = &
1080  wt_v(i,j,1)*ebuffery(ic, jc, k) + &
1081  wt_v(i,j,2)*ebuffery(ic, jc+1,k) + &
1082  wt_v(i,j,3)*ebuffery(ic+1,jc+1,k) + &
1083  wt_v(i,j,4)*ebuffery(ic+1,jc, k)
1084 
1085  end do
1086  end do
1087  end do
1088  end if
1089 
1090  if (je == npy-1) then
1091 
1092  if (is == 1) then
1093  istart = is
1094  else
1095  istart = isd
1096  end if
1097 
1098  if (ie == npx-1) then
1099  iend = ie
1100  else
1101  iend = ied
1102  end if
1103 
1104 !OMP parallel do default(none) shared(npz,jstag,npy,jed,istart,iend,istag,ind,var_nest,wt,nbuffer) private(ic,jc)
1105  do k=1,npz
1106  do j=npy+jstag_u,jed+jstag_u
1107  do i=istart,iend+istag_u
1108 
1109  ic = ind_u(i,j,1)
1110  jc = ind_u(i,j,2)
1111 
1112  u_nest(i,j,k) = &
1113  wt_u(i,j,1)*nbufferx(ic, jc, k) + &
1114  wt_u(i,j,2)*nbufferx(ic, jc+1,k) + &
1115  wt_u(i,j,3)*nbufferx(ic+1,jc+1,k) + &
1116  wt_u(i,j,4)*nbufferx(ic+1,jc, k)
1117 
1118  end do
1119  end do
1120  do j=npy+jstag_v,jed+jstag_v
1121  do i=istart,iend+istag_v
1122 
1123  ic = ind_v(i,j,1)
1124  jc = ind_v(i,j,2)
1125 
1126  v_nest(i,j,k) = &
1127  wt_v(i,j,1)*nbuffery(ic, jc, k) + &
1128  wt_v(i,j,2)*nbuffery(ic, jc+1,k) + &
1129  wt_v(i,j,3)*nbuffery(ic+1,jc+1,k) + &
1130  wt_v(i,j,4)*nbuffery(ic+1,jc, k)
1131 
1132  end do
1133  end do
1134  end do
1135  end if
1136 
1137  endif !process
1138 
1139  deallocate(wbufferx, ebufferx, sbufferx, nbufferx)
1140  deallocate(wbuffery, ebuffery, sbuffery, nbuffery)
1141 
1142  end subroutine nested_grid_bc_mpp_3d_vector
1143 
1144 
1145  subroutine nested_grid_bc_mpp_send_3d(var_coarse, nest_domain, istag, jstag, nest_level)
1147  real, dimension(:,:,:), intent(IN) :: var_coarse
1148  type(nest_domain_type), intent(INOUT) :: nest_domain
1149  integer, intent(IN) :: istag, jstag
1150  integer, intent(IN) :: nest_level
1151 
1152  real, allocatable :: wbuffer(:,:,:)
1153  real, allocatable :: ebuffer(:,:,:)
1154  real, allocatable :: sbuffer(:,:,:)
1155  real, allocatable :: nbuffer(:,:,:)
1156 
1157  integer :: i,j, ic, jc, istart, iend, k
1158 
1159  integer :: position
1160 
1161 
1162  if (istag == 1 .and. jstag == 1) then
1163  position = corner
1164  else if (istag == 0 .and. jstag == 1) then
1165  position = north
1166  else if (istag == 1 .and. jstag == 0) then
1167  position = east
1168  else
1169  position = center
1170  end if
1171 
1172 
1173  allocate(wbuffer(1,1,1))
1174 
1175  allocate(ebuffer(1,1,1))
1176 
1177  allocate(sbuffer(1,1,1))
1178 
1179  allocate(nbuffer(1,1,1))
1180 
1181 
1182  call timing_on ('COMM_TOTAL')
1183  call mpp_update_nest_fine(var_coarse, nest_domain, wbuffer, sbuffer, ebuffer, nbuffer, nest_level=nest_level, position=position)
1184  call timing_off('COMM_TOTAL')
1185 
1186 
1187  deallocate(wbuffer, ebuffer, sbuffer, nbuffer)
1188 
1189  end subroutine nested_grid_bc_mpp_send_3d
1190 
1191  subroutine nested_grid_bc_mpp_send_2d(var_coarse, nest_domain, istag, jstag, nest_level)
1193  real, dimension(:,:), intent(IN) :: var_coarse
1194  type(nest_domain_type), intent(INOUT) :: nest_domain
1195  integer, intent(IN) :: istag, jstag
1196  integer, intent(IN) :: nest_level
1197 
1198  real, allocatable :: wbuffer(:,:)
1199  real, allocatable :: ebuffer(:,:)
1200  real, allocatable :: sbuffer(:,:)
1201  real, allocatable :: nbuffer(:,:)
1202 
1203  integer :: i,j, ic, jc, istart, iend, k
1204 
1205  integer :: position
1206 
1207 
1208  if (istag == 1 .and. jstag == 1) then
1209  position = corner
1210  else if (istag == 0 .and. jstag == 1) then
1211  position = north
1212  else if (istag == 1 .and. jstag == 0) then
1213  position = east
1214  else
1215  position = center
1216  end if
1217 
1218 
1219  allocate(wbuffer(1,1))
1220 
1221  allocate(ebuffer(1,1))
1222 
1223  allocate(sbuffer(1,1))
1224 
1225  allocate(nbuffer(1,1))
1226 
1227 
1228  call timing_on ('COMM_TOTAL')
1229  call mpp_update_nest_fine(var_coarse, nest_domain, wbuffer, sbuffer, ebuffer, nbuffer, nest_level=nest_level, position=position)
1230  call timing_off('COMM_TOTAL')
1231 
1232 
1233  deallocate(wbuffer, ebuffer, sbuffer, nbuffer)
1234 
1235  end subroutine nested_grid_bc_mpp_send_2d
1236 
1237  subroutine nested_grid_bc_2d_mpp(var_nest, var_coarse, nest_domain, ind, wt, istag, jstag, &
1238  npx, npy, bd, isg, ieg, jsg, jeg, nest_level, nstep_in, nsplit_in, proc_in)
1240  type(fv_grid_bounds_type), intent(IN) :: bd
1241  real, dimension(bd%isd:bd%ied+istag,bd%jsd:bd%jed+jstag), intent(INOUT) :: var_nest
1242  real, dimension(isg:ieg+istag,jsg:jeg+jstag), intent(IN) :: var_coarse
1243  type(nest_domain_type), intent(INOUT) :: nest_domain
1244  integer, dimension(bd%isd:bd%ied+istag,bd%jsd:bd%jed+jstag,2), intent(IN) :: ind
1245  real, dimension(bd%isd:bd%ied+istag,bd%jsd:bd%jed+jstag,4), intent(IN) :: wt
1246  integer, intent(IN) :: istag, jstag, npx, npy, isg, ieg, jsg, jeg
1247  integer, intent(IN), OPTIONAL :: nest_level
1248  integer, intent(IN), OPTIONAL :: nstep_in, nsplit_in
1249  logical, intent(IN), OPTIONAL :: proc_in
1250 
1251  integer :: isw_f, iew_f, jsw_f, jew_f, isw_c, iew_c, jsw_c, jew_c
1252  integer :: ise_f, iee_f, jse_f, jee_f, ise_c, iee_c, jse_c, jee_c
1253  integer :: iss_f, ies_f, jss_f, jes_f, iss_c, ies_c, jss_c, jes_c
1254  integer :: isn_f, ien_f, jsn_f, jen_f, isn_c, ien_c, jsn_c, jen_c
1255  real, allocatable :: wbuffer(:,:)
1256  real, allocatable :: ebuffer(:,:)
1257  real, allocatable :: sbuffer(:,:)
1258  real, allocatable :: nbuffer(:,:)
1259 
1260  integer :: i,j, ic, jc, istart, iend, k
1261  integer :: nl = 1 !nest_level
1262 
1263  integer :: position
1264  logical :: process
1265 
1266  integer :: is, ie, js, je
1267  integer :: isd, ied, jsd, jed
1268 
1269  is = bd%is
1270  ie = bd%ie
1271  js = bd%js
1272  je = bd%je
1273  isd = bd%isd
1274  ied = bd%ied
1275  jsd = bd%jsd
1276  jed = bd%jed
1277 
1278  if (PRESENT(proc_in)) then
1279  process = proc_in
1280  else
1281  process = .true.
1282  endif
1283 
1284  if (PRESENT(nest_level)) then
1285  nl = nest_level
1286  endif
1287 
1288  if (istag == 1 .and. jstag == 1) then
1289  position = corner
1290  else if (istag == 0 .and. jstag == 1) then
1291  position = north
1292  else if (istag == 1 .and. jstag == 0) then
1293  position = east
1294  else
1295  position = center
1296  end if
1297 
1298  call mpp_get_c2f_index(nest_domain, isw_f, iew_f, jsw_f, jew_f, isw_c, iew_c, jsw_c, jew_c, &
1299  west, nest_level=nl, position=position)
1300  call mpp_get_c2f_index(nest_domain, ise_f, iee_f, jse_f, jee_f, ise_c, iee_c, jse_c, jee_c, &
1301  east, nest_level=nl, position=position)
1302  call mpp_get_c2f_index(nest_domain, iss_f, ies_f, jss_f, jes_f, iss_c, ies_c, jss_c, jes_c, &
1303  south, nest_level=nl, position=position)
1304  call mpp_get_c2f_index(nest_domain, isn_f, ien_f, jsn_f, jen_f, isn_c, ien_c, jsn_c, jen_c, &
1305  north, nest_level=nl, position=position)
1306 
1307  if( iew_c .GE. isw_c .AND. jew_c .GE. jsw_c ) then
1308  allocate(wbuffer(isw_c:iew_c, jsw_c:jew_c))
1309  else
1310  allocate(wbuffer(1,1))
1311  endif
1312  wbuffer = 0
1313 
1314  if( iee_c .GE. ise_c .AND. jee_c .GE. jse_c ) then
1315  allocate(ebuffer(ise_c:iee_c, jse_c:jee_c))
1316  else
1317  allocate(ebuffer(1,1))
1318  endif
1319  ebuffer = 0
1320 
1321  if( ies_c .GE. iss_c .AND. jes_c .GE. jss_c ) then
1322  allocate(sbuffer(iss_c:ies_c, jss_c:jes_c))
1323  else
1324  allocate(sbuffer(1,1))
1325  endif
1326  sbuffer = 0
1327 
1328  if( ien_c .GE. isn_c .AND. jen_c .GE. jsn_c ) then
1329  allocate(nbuffer(isn_c:ien_c, jsn_c:jen_c))
1330  else
1331  allocate(nbuffer(1,1))
1332  endif
1333  nbuffer = 0
1334 
1335  call timing_on ('COMM_TOTAL')
1336  call mpp_update_nest_fine(var_coarse, nest_domain, wbuffer, sbuffer, ebuffer, nbuffer, nest_level=nl, position=position)
1337  call timing_off('COMM_TOTAL')
1338 
1339  if (process) then
1340 
1341  if (is == 1) then
1342  do j=jsd,jed+jstag
1343  do i=isd,0
1344 
1345  ic = ind(i,j,1)
1346  jc = ind(i,j,2)
1347 
1348  var_nest(i,j) = &
1349  wt(i,j,1)*wbuffer(ic, jc) + &
1350  wt(i,j,2)*wbuffer(ic, jc+1) + &
1351  wt(i,j,3)*wbuffer(ic+1,jc+1) + &
1352  wt(i,j,4)*wbuffer(ic+1,jc)
1353 
1354  end do
1355  end do
1356  end if
1357 
1358  if (js == 1) then
1359 
1360  if (is == 1) then
1361  istart = is
1362  else
1363  istart = isd
1364  end if
1365 
1366  if (ie == npx-1) then
1367  iend = ie
1368  else
1369  iend = ied
1370  end if
1371 
1372  do j=jsd,0
1373  do i=istart,iend+istag
1374 
1375  ic = ind(i,j,1)
1376  jc = ind(i,j,2)
1377 
1378  var_nest(i,j) = &
1379  wt(i,j,1)*sbuffer(ic, jc) + &
1380  wt(i,j,2)*sbuffer(ic, jc+1) + &
1381  wt(i,j,3)*sbuffer(ic+1,jc+1) + &
1382  wt(i,j,4)*sbuffer(ic+1,jc)
1383 
1384  end do
1385  end do
1386  end if
1387 
1388 
1389  if (ie == npx-1) then
1390  do j=jsd,jed+jstag
1391  do i=npx+istag,ied+istag
1392 
1393  ic = ind(i,j,1)
1394  jc = ind(i,j,2)
1395 
1396  var_nest(i,j) = &
1397  wt(i,j,1)*ebuffer(ic, jc) + &
1398  wt(i,j,2)*ebuffer(ic, jc+1) + &
1399  wt(i,j,3)*ebuffer(ic+1,jc+1) + &
1400  wt(i,j,4)*ebuffer(ic+1,jc)
1401 
1402  end do
1403  end do
1404  end if
1405 
1406  if (je == npy-1) then
1407 
1408  if (is == 1) then
1409  istart = is
1410  else
1411  istart = isd
1412  end if
1413 
1414  if (ie == npx-1) then
1415  iend = ie
1416  else
1417  iend = ied
1418  end if
1419 
1420  do j=npy+jstag,jed+jstag
1421  do i=istart,iend+istag
1422 
1423  ic = ind(i,j,1)
1424  jc = ind(i,j,2)
1425 
1426  var_nest(i,j) = &
1427  wt(i,j,1)*nbuffer(ic, jc) + &
1428  wt(i,j,2)*nbuffer(ic, jc+1) + &
1429  wt(i,j,3)*nbuffer(ic+1,jc+1) + &
1430  wt(i,j,4)*nbuffer(ic+1,jc)
1431 
1432  end do
1433  end do
1434  end if
1435 
1436  endif !process
1437 
1438  deallocate(wbuffer, ebuffer, sbuffer, nbuffer)
1439 
1440  end subroutine nested_grid_bc_2d_mpp
1441 
1442  subroutine nested_grid_bc_2d(var_nest, var_coarse, ind, wt, istag, jstag, &
1443  npx, npy, bd, isg, ieg, jsg, jeg, nstep_in, nsplit_in)
1445  type(fv_grid_bounds_type), intent(IN) :: bd
1446  real, dimension(bd%isd:bd%ied+istag,bd%jsd:bd%jed+jstag), intent(INOUT) :: var_nest
1447  real, dimension(isg:ieg+istag,jsg:jeg+jstag), intent(IN) :: var_coarse
1448  integer, dimension(bd%isd:bd%ied+istag,bd%jsd:bd%jed+jstag,2), intent(IN) :: ind
1449  real, dimension(bd%isd:bd%ied+istag,bd%jsd:bd%jed+jstag,4), intent(IN) :: wt
1450  integer, intent(IN) :: istag, jstag, npx, npy, isg, ieg, jsg, jeg
1451  integer, intent(IN), OPTIONAL :: nstep_in, nsplit_in
1452 
1453  integer :: nstep, nsplit
1454 
1455  integer :: i,j, ic, jc, istart, iend
1456 
1457  integer :: is, ie, js, je
1458  integer :: isd, ied, jsd, jed
1459 
1460  is = bd%is
1461  ie = bd%ie
1462  js = bd%js
1463  je = bd%je
1464  isd = bd%isd
1465  ied = bd%ied
1466  jsd = bd%jsd
1467  jed = bd%jed
1468 
1469  if ( .not. present(nstep_in) .or. .not. present(nsplit_in) ) then
1470  nstep = 1
1471  nsplit = 2
1472  else
1473  nstep = nstep_in
1474  nsplit = nsplit_in
1475  end if
1476 
1477  if (is == 1) then
1478  do j=jsd,jed+jstag
1479  do i=isd,0
1480 
1481  ic = ind(i,j,1)
1482  jc = ind(i,j,2)
1483 
1484  var_nest(i,j) = &
1485  wt(i,j,1)*var_coarse(ic, jc) + &
1486  wt(i,j,2)*var_coarse(ic, jc+1) + &
1487  wt(i,j,3)*var_coarse(ic+1,jc+1) + &
1488  wt(i,j,4)*var_coarse(ic+1,jc)
1489 
1490  end do
1491  end do
1492  end if
1493 
1494  if (js == 1) then
1495 
1496  if (is == 1) then
1497  istart = is
1498  else
1499  istart = isd
1500  end if
1501 
1502  if (ie == npx-1) then
1503  iend = ie
1504  else
1505  iend = ied
1506  end if
1507 
1508  do j=jsd,0
1509  do i=istart,iend+istag
1510 
1511  ic = ind(i,j,1)
1512  jc = ind(i,j,2)
1513 
1514  var_nest(i,j) = &
1515  wt(i,j,1)*var_coarse(ic, jc) + &
1516  wt(i,j,2)*var_coarse(ic, jc+1) + &
1517  wt(i,j,3)*var_coarse(ic+1,jc+1) + &
1518  wt(i,j,4)*var_coarse(ic+1,jc)
1519 
1520  end do
1521  end do
1522  end if
1523 
1524 
1525  if (ie == npx-1) then
1526  do j=jsd,jed+jstag
1527  do i=npx+istag,ied+istag
1528 
1529  ic = ind(i,j,1)
1530  jc = ind(i,j,2)
1531 
1532  var_nest(i,j) = &
1533  wt(i,j,1)*var_coarse(ic, jc) + &
1534  wt(i,j,2)*var_coarse(ic, jc+1) + &
1535  wt(i,j,3)*var_coarse(ic+1,jc+1) + &
1536  wt(i,j,4)*var_coarse(ic+1,jc)
1537 
1538  end do
1539  end do
1540  end if
1541 
1542  if (je == npy-1) then
1543 
1544  if (is == 1) then
1545  istart = is
1546  else
1547  istart = isd
1548  end if
1549 
1550  if (ie == npx-1) then
1551  iend = ie
1552  else
1553  iend = ied
1554  end if
1555 
1556 
1557  do j=npy+jstag,jed+jstag
1558  do i=istart,iend+istag
1559 
1560  ic = ind(i,j,1)
1561  jc = ind(i,j,2)
1562 
1563  var_nest(i,j) = &
1564  wt(i,j,1)*var_coarse(ic, jc) + &
1565  wt(i,j,2)*var_coarse(ic, jc+1) + &
1566  wt(i,j,3)*var_coarse(ic+1,jc+1) + &
1567  wt(i,j,4)*var_coarse(ic+1,jc)
1568 
1569  end do
1570  end do
1571  end if
1572 
1573 
1574 
1575  end subroutine nested_grid_bc_2d
1576 
1577  subroutine nested_grid_bc_3d(var_nest, var_coarse, ind, wt, istag, jstag, &
1578  npx, npy, npz, bd, isg, ieg, jsg, jeg, nstep_in, nsplit_in)
1580  type(fv_grid_bounds_type), intent(IN) :: bd
1581  real, dimension(bd%isd:bd%ied+istag,bd%jsd:bd%jed+jstag,npz), intent(INOUT) :: var_nest
1582  real, dimension(isg:ieg+istag,jsg:jeg+jstag,npz), intent(IN) :: var_coarse
1583  integer, dimension(bd%isd:bd%ied+istag,bd%jsd:bd%jed+jstag,2), intent(IN) :: ind
1584  real, dimension(bd%isd:bd%ied+istag,bd%jsd:bd%jed+jstag,4), intent(IN) :: wt
1585  integer, intent(IN) :: istag, jstag, npx, npy, isg, ieg, jsg, jeg, npz
1586  integer, intent(IN), OPTIONAL :: nstep_in, nsplit_in
1587 
1588  integer :: nstep, nsplit
1589 
1590  integer :: i,j, ic, jc, istart, iend, k
1591 
1592  integer :: is, ie, js, je
1593  integer :: isd, ied, jsd, jed
1594 
1595  is = bd%is
1596  ie = bd%ie
1597  js = bd%js
1598  je = bd%je
1599  isd = bd%isd
1600  ied = bd%ied
1601  jsd = bd%jsd
1602  jed = bd%jed
1603 
1604  if ( .not. present(nstep_in) .or. .not. present(nsplit_in) ) then
1605  nstep = 1
1606  nsplit = 2
1607  else
1608  nstep = nstep_in
1609  nsplit = nsplit_in
1610  end if
1611 
1612  if (is == 1) then
1613 !OMP parallel do default(none) shared(npz,jsd,jed,jstag,isd,ind,var_nest,wt,var_coarse) private(ic,jc)
1614  do k=1,npz
1615  do j=jsd,jed+jstag
1616  do i=isd,0
1617 
1618  ic = ind(i,j,1)
1619  jc = ind(i,j,2)
1620 
1621  var_nest(i,j,k) = &
1622  wt(i,j,1)*var_coarse(ic, jc, k) + &
1623  wt(i,j,2)*var_coarse(ic, jc+1,k) + &
1624  wt(i,j,3)*var_coarse(ic+1,jc+1,k) + &
1625  wt(i,j,4)*var_coarse(ic+1,jc, k)
1626 
1627  end do
1628  end do
1629  end do
1630  end if
1631 
1632  if (js == 1) then
1633 
1634  if (is == 1) then
1635  istart = is
1636  else
1637  istart = isd
1638  end if
1639 
1640  if (ie == npx-1) then
1641  iend = ie
1642  else
1643  iend = ied
1644  end if
1645 
1646 !OMP parallel do default(none) shared(npz,jsd,istart,iend,istag,ind,var_nest,wt,var_coarse) private(ic,jc)
1647  do k=1,npz
1648  do j=jsd,0
1649  do i=istart,iend+istag
1650 
1651  ic = ind(i,j,1)
1652  jc = ind(i,j,2)
1653 
1654  var_nest(i,j,k) = &
1655  wt(i,j,1)*var_coarse(ic, jc, k) + &
1656  wt(i,j,2)*var_coarse(ic, jc+1,k) + &
1657  wt(i,j,3)*var_coarse(ic+1,jc+1,k) + &
1658  wt(i,j,4)*var_coarse(ic+1,jc, k)
1659 
1660  end do
1661  end do
1662  end do
1663  end if
1664 
1665 
1666  if (ie == npx-1) then
1667 !OMP parallel do default(none) shared(npz,jsd,jed,jstag,npx,ied,istag,ind,var_nest,wt,var_coarse) private(ic,jc)
1668  do k=1,npz
1669  do j=jsd,jed+jstag
1670  do i=npx+istag,ied+istag
1671 
1672  ic = ind(i,j,1)
1673  jc = ind(i,j,2)
1674 
1675  var_nest(i,j,k) = &
1676  wt(i,j,1)*var_coarse(ic, jc, k) + &
1677  wt(i,j,2)*var_coarse(ic, jc+1,k) + &
1678  wt(i,j,3)*var_coarse(ic+1,jc+1,k) + &
1679  wt(i,j,4)*var_coarse(ic+1,jc, k)
1680 
1681  end do
1682  end do
1683  end do
1684  end if
1685 
1686  if (je == npy-1) then
1687 
1688  if (is == 1) then
1689  istart = is
1690  else
1691  istart = isd
1692  end if
1693 
1694  if (ie == npx-1) then
1695  iend = ie
1696  else
1697  iend = ied
1698  end if
1699 
1700 !OMP parallel do default(none) shared(npz,npy,jed,jstag,istart,iend,istag,ind,var_nest,wt,var_coarse) private(ic,jc)
1701  do k=1,npz
1702  do j=npy+jstag,jed+jstag
1703  do i=istart,iend+istag
1704 
1705  ic = ind(i,j,1)
1706  jc = ind(i,j,2)
1707 
1708  var_nest(i,j,k) = &
1709  wt(i,j,1)*var_coarse(ic, jc, k) + &
1710  wt(i,j,2)*var_coarse(ic, jc+1,k) + &
1711  wt(i,j,3)*var_coarse(ic+1,jc+1,k) + &
1712  wt(i,j,4)*var_coarse(ic+1,jc, k)
1713 
1714  end do
1715  end do
1716  end do
1717  end if
1718 
1719 
1720 
1721  end subroutine nested_grid_bc_3d
1723  subroutine nested_grid_bc_send_scalar(var_coarse, nest_domain, istag, jstag, nest_level)
1725  real, dimension(:,:,:), intent(IN) :: var_coarse
1726  type(nest_domain_type), intent(INOUT) :: nest_domain
1727  integer, intent(IN) :: istag, jstag
1728  integer, intent(IN) :: nest_level
1729 
1730  integer :: position
1731 
1732  real :: wbuffer(1,1,1)
1733  real :: ebuffer(1,1,1)
1734  real :: sbuffer(1,1,1)
1735  real :: nbuffer(1,1,1)
1736 
1737 
1738  if (istag == 1 .and. jstag == 1) then
1739  position = corner
1740  else if (istag == 0 .and. jstag == 1) then
1741  position = north
1742  else if (istag == 1 .and. jstag == 0) then
1743  position = east
1744  else
1745  position = center
1746  end if
1747 
1748  call timing_on ('COMM_TOTAL')
1749  call mpp_update_nest_fine(var_coarse, nest_domain, wbuffer, sbuffer, ebuffer, nbuffer, nest_level=nest_level, position=position)
1750  call timing_off('COMM_TOTAL')
1751 
1752  end subroutine nested_grid_bc_send_scalar
1753 
1754  subroutine nested_grid_bc_recv_scalar(nest_domain, istag, jstag, npz, &
1755  bd, nest_BC_buffers, nest_level)
1757  type(fv_grid_bounds_type), intent(IN) :: bd
1758  type(nest_domain_type), intent(INOUT) :: nest_domain
1759  integer, intent(IN) :: istag, jstag, npz
1760  integer, intent(IN) :: nest_level
1761 
1762  type(fv_nest_bc_type_3d), intent(INOUT), target :: nest_BC_buffers
1763 
1764  real, dimension(bd%isd:bd%ied+istag,bd%jsd:bd%jed+jstag,npz) :: var_coarse_dummy
1765 
1766  integer :: position
1767 
1768 !!$ integer :: isw_f, iew_f, jsw_f, jew_f, isw_c, iew_c, jsw_c, jew_c
1769 !!$ integer :: ise_f, iee_f, jse_f, jee_f, ise_c, iee_c, jse_c, jee_c
1770 !!$ integer :: iss_f, ies_f, jss_f, jes_f, iss_c, ies_c, jss_c, jes_c
1771 !!$ integer :: isn_f, ien_f, jsn_f, jen_f, isn_c, ien_c, jsn_c, jen_c
1772 
1773  integer :: i,j, k
1774 
1775  if (istag == 1 .and. jstag == 1) then
1776  position = corner
1777  else if (istag == 0 .and. jstag == 1) then
1778  position = north
1779  else if (istag == 1 .and. jstag == 0) then
1780  position = east
1781  else
1782  position = center
1783  end if
1784 
1785  if (.not. allocated(nest_bc_buffers%west_t1) ) then
1786  call init_nest_bc_type(nest_domain, nest_bc_buffers, npz, nest_level, position)
1787  endif
1788 
1789  call timing_on ('COMM_TOTAL')
1790  call mpp_update_nest_fine(var_coarse_dummy, nest_domain, nest_bc_buffers%west_t1, nest_bc_buffers%south_t1, &
1791  nest_bc_buffers%east_t1, nest_bc_buffers%north_t1, nest_level=nest_level, position=position)
1792  call timing_off('COMM_TOTAL')
1793 
1794  end subroutine nested_grid_bc_recv_scalar
1795 
1796  subroutine nested_grid_bc_send_vector(u_coarse, v_coarse, nest_domain, nest_level, flags, gridtype)
1797  real, dimension(:,:,:), intent(IN) :: u_coarse, v_coarse
1798  type(nest_domain_type), intent(INOUT) :: nest_domain
1799  integer, intent(IN) :: nest_level
1800  integer, optional, intent(IN) :: flags, gridtype
1801 
1802  real :: wbufferx(1,1,1), wbuffery(1,1,1)
1803  real :: ebufferx(1,1,1), ebuffery(1,1,1)
1804  real :: sbufferx(1,1,1), sbuffery(1,1,1)
1805  real :: nbufferx(1,1,1), nbuffery(1,1,1)
1806 
1807  integer :: nl = 1
1808 
1809  call timing_on ('COMM_TOTAL')
1810  call mpp_update_nest_fine(u_coarse, v_coarse, nest_domain, wbufferx,wbuffery, sbufferx, sbuffery, &
1811  ebufferx, ebuffery, nbufferx, nbuffery, nest_level=nest_level, flags=flags, gridtype=gridtype)
1812  call timing_off('COMM_TOTAL')
1813 
1814  end subroutine nested_grid_bc_send_vector
1815 
1816  subroutine init_nest_bc_type(nest_domain, nest_BC_buffers, npz, nest_level, position)
1817  type(nest_domain_type), intent(INOUT) :: nest_domain
1818  type(fv_nest_bc_type_3d), intent(INOUT) :: nest_BC_buffers
1819  integer, intent(IN) :: npz, position, nest_level
1820 
1821  integer :: isw_f, iew_f, jsw_f, jew_f, isw_c, iew_c, jsw_c, jew_c
1822  integer :: ise_f, iee_f, jse_f, jee_f, ise_c, iee_c, jse_c, jee_c
1823  integer :: iss_f, ies_f, jss_f, jes_f, iss_c, ies_c, jss_c, jes_c
1824  integer :: isn_f, ien_f, jsn_f, jen_f, isn_c, ien_c, jsn_c, jen_c
1825  integer :: i, j, k
1826 
1827  call mpp_get_c2f_index(nest_domain, isw_f, iew_f, jsw_f, jew_f, isw_c, iew_c, jsw_c, jew_c, &
1828  west, nest_level=nest_level, position=position)
1829  call mpp_get_c2f_index(nest_domain, ise_f, iee_f, jse_f, jee_f, ise_c, iee_c, jse_c, jee_c, &
1830  east, nest_level=nest_level, position=position)
1831  call mpp_get_c2f_index(nest_domain, iss_f, ies_f, jss_f, jes_f, iss_c, ies_c, jss_c, jes_c, &
1832  south, nest_level=nest_level, position=position)
1833  call mpp_get_c2f_index(nest_domain, isn_f, ien_f, jsn_f, jen_f, isn_c, ien_c, jsn_c, jen_c, &
1834  north, nest_level=nest_level, position=position)
1835 
1836  if( iew_c .GE. isw_c .AND. jew_c .GE. jsw_c ) then
1837  If (.not. allocated(nest_bc_buffers%west_t1)) allocate(nest_bc_buffers%west_t1(isw_c:iew_c, jsw_c:jew_c,npz))
1838  !compatible with first touch principle
1839 !OMP parallel do default(none) shared(npz,jsw_c,jew_c,isw_c,iew_c,nest_BC_buffers)
1840  do k=1,npz
1841  do j=jsw_c,jew_c
1842  do i=isw_c,iew_c
1843  nest_bc_buffers%west_t1(i,j,k) = 1.e24
1844  enddo
1845  enddo
1846  enddo
1847  else
1848  allocate(nest_bc_buffers%west_t1(1,1,1))
1849  nest_bc_buffers%west_t1(1,1,1) = 1.e24
1850  endif
1851 
1852  if( iee_c .GE. ise_c .AND. jee_c .GE. jse_c ) then
1853  If (.not. allocated(nest_bc_buffers%east_t1)) allocate(nest_bc_buffers%east_t1(ise_c:iee_c, jse_c:jee_c,npz))
1854 !OMP parallel do default(none) shared(npz,jse_c,jee_c,ise_c,iee_c,nest_BC_buffers)
1855  do k=1,npz
1856  do j=jse_c,jee_c
1857  do i=ise_c,iee_c
1858  nest_bc_buffers%east_t1(i,j,k) = 1.e24
1859  enddo
1860  enddo
1861  enddo
1862  else
1863  allocate(nest_bc_buffers%east_t1(1,1,1))
1864  nest_bc_buffers%east_t1(1,1,1) = 1.e24
1865  endif
1866 
1867  if( ies_c .GE. iss_c .AND. jes_c .GE. jss_c ) then
1868  If (.not. allocated(nest_bc_buffers%south_t1)) allocate(nest_bc_buffers%south_t1(iss_c:ies_c, jss_c:jes_c,npz))
1869 !OMP parallel do default(none) shared(npz,jss_c,jes_c,iss_c,ies_c,nest_BC_buffers)
1870  do k=1,npz
1871  do j=jss_c,jes_c
1872  do i=iss_c,ies_c
1873  nest_bc_buffers%south_t1(i,j,k) = 1.e24
1874  enddo
1875  enddo
1876  enddo
1877  else
1878  allocate(nest_bc_buffers%south_t1(1,1,1))
1879  nest_bc_buffers%south_t1(1,1,1) = 1.e24
1880  endif
1881 
1882  if( ien_c .GE. isn_c .AND. jen_c .GE. jsn_c ) then
1883  If (.not. allocated(nest_bc_buffers%north_t1)) allocate(nest_bc_buffers%north_t1(isn_c:ien_c, jsn_c:jen_c,npz))
1884 !OMP parallel do default(none) shared(npz,jsn_c,jen_c,isn_c,ien_c,nest_BC_buffers)
1885  do k=1,npz
1886  do j=jsn_c,jen_c
1887  do i=isn_c,ien_c
1888  nest_bc_buffers%north_t1(i,j,k) = 1.e24
1889  enddo
1890  enddo
1891  enddo
1892  else
1893  allocate(nest_bc_buffers%north_t1(1,1,1))
1894  nest_bc_buffers%north_t1(1,1,1) = 1.e24
1895  endif
1896 
1897 
1898  end subroutine init_nest_bc_type
1899 
1900  subroutine nested_grid_bc_recv_vector(nest_domain, npz, bd, nest_BC_u_buffers, nest_BC_v_buffers, nest_level, flags, gridtype)
1902  type(fv_grid_bounds_type), intent(IN) :: bd
1903  type(nest_domain_type), intent(INOUT) :: nest_domain
1904  integer, intent(IN) :: npz
1905  type(fv_nest_bc_type_3d), intent(INOUT), target :: nest_BC_u_buffers, nest_BC_v_buffers
1906  integer, intent(IN) :: nest_level
1907  integer, optional, intent(IN) :: flags, gridtype
1908 
1909  real, dimension(1,1,npz) :: u_coarse_dummy, v_coarse_dummy
1910 
1911  integer :: i,j, k
1912  integer :: position_x, position_y
1913 
1914  call get_vector_position(position_x, position_y, gridtype)
1915 
1916  if (.not. allocated(nest_bc_u_buffers%west_t1) ) then
1917  call init_nest_bc_type(nest_domain, nest_bc_u_buffers, npz, nest_level, position_x)
1918  endif
1919  if (.not. allocated(nest_bc_v_buffers%west_t1) ) then
1920  call init_nest_bc_type(nest_domain, nest_bc_v_buffers, npz, nest_level, position_y)
1921  endif
1922 
1923  call timing_on ('COMM_TOTAL')
1924  call mpp_update_nest_fine(u_coarse_dummy, v_coarse_dummy, nest_domain, &
1925  nest_bc_u_buffers%west_t1, nest_bc_v_buffers%west_t1, nest_bc_u_buffers%south_t1, nest_bc_v_buffers%south_t1, &
1926  nest_bc_u_buffers%east_t1, nest_bc_v_buffers%east_t1, nest_bc_u_buffers%north_t1, nest_bc_v_buffers%north_t1, &
1927  nest_level, flags, gridtype)
1928  call timing_off('COMM_TOTAL')
1929 
1930  end subroutine nested_grid_bc_recv_vector
1931 
1934  subroutine nested_grid_bc_save_proc(nest_domain, ind, wt, istag, jstag, &
1935  npx, npy, npz, bd, nest_BC, nest_BC_buffers, pd_in)
1937  type(fv_grid_bounds_type), intent(IN) :: bd
1938  type(nest_domain_type), intent(INOUT) :: nest_domain
1939  integer, dimension(bd%isd:bd%ied+istag,bd%jsd:bd%jed+jstag,2), intent(IN) :: ind
1940  real, dimension(bd%isd:bd%ied+istag,bd%jsd:bd%jed+jstag,4), intent(IN) :: wt
1941  integer, intent(IN) :: istag, jstag, npx, npy, npz
1942  logical, intent(IN), OPTIONAL :: pd_in
1943 
1944  !!NOTE: if declaring an ALLOCATABLE array with intent(OUT), the resulting dummy array
1945  !! will NOT be allocated! This goes for allocatable members of derived types as well.
1946  type(fv_nest_bc_type_3d), intent(INOUT), target :: nest_BC, nest_BC_buffers
1947 
1948  real, dimension(bd%isd:bd%ied+istag,bd%jsd:bd%jed+jstag,npz) :: var_coarse_dummy
1949 
1950  real, dimension(:,:,:), pointer :: var_east, var_west, var_south, var_north
1951  real, dimension(:,:,:), pointer :: buf_east, buf_west, buf_south, buf_north
1952 
1953  integer :: position
1954 
1955 
1956  integer :: i,j, k, ic, jc, istart, iend
1957  logical :: process, pd = .false.
1958 
1959  integer :: is, ie, js, je
1960  integer :: isd, ied, jsd, jed
1961 
1962  is = bd%is
1963  ie = bd%ie
1964  js = bd%js
1965  je = bd%je
1966  isd = bd%isd
1967  ied = bd%ied
1968  jsd = bd%jsd
1969  jed = bd%jed
1970 
1971 
1972  if (present(pd_in)) then
1973  pd = pd_in
1974  else
1975  pd = .false.
1976  endif
1977 
1978 
1979  var_east => nest_bc%east_t1
1980  var_west => nest_bc%west_t1
1981  var_north => nest_bc%north_t1
1982  var_south => nest_bc%south_t1
1983 
1984  buf_east => nest_bc_buffers%east_t1
1985  buf_west => nest_bc_buffers%west_t1
1986  buf_north => nest_bc_buffers%north_t1
1987  buf_south => nest_bc_buffers%south_t1
1988  ! ?buffer has uninterpolated coarse-grid data; need to perform interpolation ourselves
1989  !To do this more securely, instead of using is/etc we could use the fine-grid indices defined above
1990  if (is == 1 ) then
1991 
1992 !$OMP parallel do default(none) shared(npz,isd,ied,jsd,jed,jstag,ind,var_west,wt,buf_west) private(ic,jc)
1993  do k=1,npz
1994  do j=jsd,jed+jstag
1995  do i=isd,0
1996 
1997  ic = ind(i,j,1)
1998  jc = ind(i,j,2)
1999 
2000 
2001  var_west(i,j,k) = &
2002  wt(i,j,1)*buf_west(ic, jc,k) + &
2003  wt(i,j,2)*buf_west(ic, jc+1,k) + &
2004  wt(i,j,3)*buf_west(ic+1,jc+1,k) + &
2005  wt(i,j,4)*buf_west(ic+1,jc,k)
2006 
2007  end do
2008  end do
2009  end do
2010 
2011  if (pd) then
2012 !$OMP parallel do default(none) shared(npz,jsd,jed,jstag,isd,var_west,nest_BC)
2013  do k=1,npz
2014  do j=jsd,jed+jstag
2015  do i=isd,0
2016 
2017  var_west(i,j,k) = max(var_west(i,j,k), 0.5*nest_bc%west_t0(i,j,k))
2018  end do
2019  end do
2020  end do
2021  endif
2022 
2023  end if
2024 
2025  if (js == 1 ) then
2026 
2027  if (is == 1) then
2028  istart = is
2029  else
2030  istart = isd
2031  end if
2032 
2033  if (ie == npx-1) then
2034  iend = ie
2035  else
2036  iend = ied
2037  end if
2038 
2039 !$OMP parallel do default(none) shared(npz,istart,iend,jsd,jed,istag,ind,var_south,wt,buf_south) private(ic,jc)
2040  do k=1,npz
2041  do j=jsd,0
2042  do i=istart,iend+istag
2043 
2044  ic = ind(i,j,1)
2045  jc = ind(i,j,2)
2046 
2047 
2048  var_south(i,j,k) = &
2049  wt(i,j,1)*buf_south(ic, jc,k) + &
2050  wt(i,j,2)*buf_south(ic, jc+1,k) + &
2051  wt(i,j,3)*buf_south(ic+1,jc+1,k) + &
2052  wt(i,j,4)*buf_south(ic+1,jc,k)
2053 
2054  end do
2055  end do
2056  end do
2057 
2058  if (pd) then
2059 !$OMP parallel do default(none) shared(npz,jsd,jed,istart,iend,istag,var_south,nest_BC)
2060  do k=1,npz
2061  do j=jsd,0
2062  do i=istart,iend+istag
2063 
2064  var_south(i,j,k) = max(var_south(i,j,k), 0.5*nest_bc%south_t0(i,j,k))
2065 
2066  end do
2067  end do
2068  end do
2069  endif
2070 
2071  end if
2072 
2073 
2074  if (ie == npx-1 ) then
2075 
2076 !$OMP parallel do default(none) shared(npx,npz,isd,ied,jsd,jed,istag,jstag,ind,var_east,wt,buf_east) private(ic,jc)
2077  do k=1,npz
2078  do j=jsd,jed+jstag
2079  do i=npx+istag,ied+istag
2080 
2081  ic = ind(i,j,1)
2082  jc = ind(i,j,2)
2083 
2084 
2085  var_east(i,j,k) = &
2086  wt(i,j,1)*buf_east(ic, jc,k) + &
2087  wt(i,j,2)*buf_east(ic, jc+1,k) + &
2088  wt(i,j,3)*buf_east(ic+1,jc+1,k) + &
2089  wt(i,j,4)*buf_east(ic+1,jc,k)
2090 
2091  end do
2092  end do
2093  end do
2094 
2095  if (pd) then
2096 !$OMP parallel do default(none) shared(npx,npz,jsd,jed,istag,jstag,ied,var_east,nest_BC)
2097  do k=1,npz
2098  do j=jsd,jed+jstag
2099  do i=npx+istag,ied+istag
2100 
2101  var_east(i,j,k) = max(var_east(i,j,k), 0.5*nest_bc%east_t0(i,j,k))
2102 
2103  end do
2104  end do
2105  end do
2106  endif
2107 
2108  end if
2109 
2110  if (je == npy-1 ) then
2111 
2112  if (is == 1) then
2113  istart = is
2114  else
2115  istart = isd
2116  end if
2117 
2118  if (ie == npx-1) then
2119  iend = ie
2120  else
2121  iend = ied
2122  end if
2123 
2124 !$OMP parallel do default(none) shared(npy,npz,istart,iend,jsd,jed,istag,jstag,ind,var_north,wt,buf_north) private(ic,jc)
2125  do k=1,npz
2126  do j=npy+jstag,jed+jstag
2127  do i=istart,iend+istag
2128 
2129  ic = ind(i,j,1)
2130  jc = ind(i,j,2)
2131 
2132 
2133  var_north(i,j,k) = &
2134  wt(i,j,1)*buf_north(ic, jc,k) + &
2135  wt(i,j,2)*buf_north(ic, jc+1,k) + &
2136  wt(i,j,3)*buf_north(ic+1,jc+1,k) + &
2137  wt(i,j,4)*buf_north(ic+1,jc,k)
2138 
2139  end do
2140  end do
2141  end do
2142 
2143  if (pd) then
2144 !$OMP parallel do default(none) shared(npy,npz,jsd,jed,istart,iend,istag,jstag,ied,var_north,nest_BC)
2145  do k=1,npz
2146  do j=npy+jstag,jed+jstag
2147  do i=istart,iend+istag
2148 
2149  var_north(i,j,k) = max(var_north(i,j,k), 0.5*nest_bc%north_t0(i,j,k))
2150 
2151  end do
2152  end do
2153  end do
2154  endif
2155 
2156  end if
2157 
2158  end subroutine nested_grid_bc_save_proc
2159 
2160 
2161  ! A NOTE ON BCTYPE: currently only an interpolation BC is implemented,
2162  ! bctype >= 2 currently correspond
2163  ! to a flux BC on the tracers ONLY, which is implemented in fv_tracer.
2164 
2168  subroutine nested_grid_bc_apply_intt(var_nest, istag, jstag, &
2169  npx, npy, npz, bd, step, split, &
2170  BC, bctype)
2172  type(fv_grid_bounds_type), intent(IN) :: bd
2173  real, dimension(bd%isd:bd%ied+istag,bd%jsd:bd%jed+jstag, npz), intent(INOUT) :: var_nest
2174  integer, intent(IN) :: istag, jstag, npx, npy, npz
2175  real, intent(IN) :: split, step
2176  integer, intent(IN) :: bctype
2177 
2178  type(fv_nest_bc_type_3d), intent(IN), target :: BC
2179  real, pointer, dimension(:,:,:) :: var_t0, var_t1
2180 
2181  integer :: i,j, istart, iend, k
2182  real :: denom
2183 
2184  logical, save :: printdiag = .true.
2185 
2186  integer :: is, ie, js, je
2187  integer :: isd, ied, jsd, jed
2188 
2189  is = bd%is
2190  ie = bd%ie
2191  js = bd%js
2192  je = bd%je
2193  isd = bd%isd
2194  ied = bd%ied
2195  jsd = bd%jsd
2196  jed = bd%jed
2197 
2198  denom = 1./split
2199  if (is == 1 ) then
2200  var_t0 => bc%west_t0
2201  var_t1 => bc%west_t1
2202 !OMP parallel do default(none) shared(npz,jsd,jed,jstag,isd,var_nest,var_t0,var_t1,split,step,denom)
2203  do k=1,npz
2204  do j=jsd,jed+jstag
2205  do i=isd,0
2206  var_nest(i,j,k) = (var_t0(i,j,k)*(split-step) + step*var_t1(i,j,k))*denom
2207  end do
2208 
2209  end do
2210  end do
2211  end if
2212 
2213  if (js == 1 ) then
2214 
2215  if (is == 1) then
2216  istart = is
2217  else
2218  istart = isd
2219  end if
2220 
2221  if (ie == npx-1) then
2222  iend = ie
2223  else
2224  iend = ied
2225  end if
2226 
2227  var_t0 => bc%south_t0
2228  var_t1 => bc%south_t1
2229 !OMP parallel do default(none) shared(npz,jsd,istart,iend,istag,var_nest,var_t0,var_t1,split,step,denom)
2230  do k=1,npz
2231  do j=jsd,0
2232  do i=istart,iend+istag
2233 
2234  var_nest(i,j,k) = (var_t0(i,j,k)*(split-step) + step*var_t1(i,j,k))*denom
2235  end do
2236  end do
2237  end do
2238  end if
2239 
2240 
2241  if (ie == npx-1 ) then
2242  var_t0 => bc%east_t0
2243  var_t1 => bc%east_t1
2244 !OMP parallel do default(none) shared(npz,jsd,jed,jstag,npx,isd,istag,var_nest,var_t0,var_t1,split,step,denom)
2245  do k=1,npz
2246  do j=jsd,jed+jstag
2247  do i=npx+istag,ied+istag
2248  var_nest(i,j,k) = (var_t0(i,j,k)*(split-step) + step*var_t1(i,j,k))*denom
2249 
2250  end do
2251  end do
2252  end do
2253 
2254  end if
2255 
2256  if (je == npy-1 ) then
2257 
2258  if (is == 1) then
2259  istart = is
2260  else
2261  istart = isd
2262  end if
2263 
2264  if (ie == npx-1) then
2265  iend = ie
2266  else
2267  iend = ied
2268  end if
2269 
2270  var_t0 => bc%north_t0
2271  var_t1 => bc%north_t1
2272 !OMP parallel do default(none) shared(npz,npy,jed,jstag,istart,iend,istag,var_nest,var_t0,var_t1,split,step,denom)
2273  do k=1,npz
2274  do j=npy+jstag,jed+jstag
2275  do i=istart,iend+istag
2276 
2277  var_nest(i,j,k) = (var_t0(i,j,k)*(split-step) + step*var_t1(i,j,k))*denom
2278 
2279  end do
2280  end do
2281  end do
2282 
2283  end if
2284 
2285 
2286  end subroutine nested_grid_bc_apply_intt
2287 
2288  subroutine update_coarse_grid_mpp_2d(var_coarse, var_nest, nest_domain, dx, dy, area, &
2289  bd, isd_p, ied_p, jsd_p, jed_p, is_n, ie_n, js_n, je_n, isu, ieu, jsu, jeu, npx, npy, &
2290  istag, jstag, r, nestupdate, upoff, nsponge, parent_proc, child_proc, parent_grid, nest_level)
2292  type(fv_grid_bounds_type), intent(IN) :: bd
2293  integer, intent(IN) :: isd_p, ied_p, jsd_p, jed_p, is_n, ie_n, js_n, je_n
2294  integer, intent(IN) :: isu, ieu, jsu, jeu
2295  integer, intent(IN) :: istag, jstag, r, nestupdate, upoff, nsponge
2296  integer, intent(IN) :: npx, npy
2297  real, intent(IN), target :: var_nest(is_n:ie_n+istag,js_n:je_n+jstag)
2298  real, intent(INOUT), target :: var_coarse(isd_p:ied_p+istag,jsd_p:jed_p+jstag)
2299  real, intent(IN) :: dx(bd%isd:bd%ied, bd%jsd:bd%jed+1)
2300  real, intent(IN) :: dy(bd%isd:bd%ied+1,bd%jsd:bd%jed)
2301  real, intent(IN) :: area(bd%isd:bd%ied,bd%jsd:bd%jed)
2302  logical, intent(IN) :: parent_proc, child_proc
2303  type(fv_atmos_type), pointer, intent(IN) :: parent_grid
2304  type(nest_domain_type), intent(INOUT) :: nest_domain
2305  integer, intent(IN) :: nest_level
2306 
2307  real :: var_nest_3d(is_n:ie_n+istag,js_n:je_n+jstag,1)
2308  real :: var_coarse_3d(isd_p:ied_p+istag,jsd_p:jed_p+jstag,1)
2309  pointer(ptr_nest, var_nest_3d)
2310  pointer(ptr_coarse, var_coarse_3d)
2311 
2312 ! if (child_proc .and. size(var_nest) > 1) ptr_nest = LOC(var_nest)
2313 ! if (parent_proc .and. size(var_coarse) > 1) ptr_coarse = LOC(var_coarse)
2314  ptr_nest = loc(var_nest)
2315  ptr_coarse = loc(var_coarse)
2316  call update_coarse_grid_mpp(var_coarse_3d, var_nest_3d, &
2317  nest_domain, dx, dy, area, &
2318  bd, isd_p, ied_p, jsd_p, jed_p, is_n, ie_n, js_n, je_n, &
2319  isu, ieu, jsu, jeu, npx, npy, 1, &
2320  istag, jstag, r, nestupdate, upoff, nsponge, &
2321  parent_proc, child_proc, parent_grid, nest_level )
2322 
2323  end subroutine update_coarse_grid_mpp_2d
2324 
2325 
2326  subroutine update_coarse_grid_mpp(var_coarse, var_nest, nest_domain, dx, dy, area, &
2327  bd, isd_p, ied_p, jsd_p, jed_p, is_n, ie_n, js_n, je_n, &
2328  isu, ieu, jsu, jeu, npx, npy, npz, &
2329  istag, jstag, r, nestupdate, upoff, nsponge, &
2330  parent_proc, child_proc, parent_grid, nest_level)
2332  !This routine assumes the coarse and nested grids are properly
2333  ! aligned, and that in particular for odd refinement ratios all
2334  ! coarse-grid cells (faces) coincide with nested-grid cells (faces)
2335 
2336  type(fv_grid_bounds_type), intent(IN) :: bd
2337  integer, intent(IN) :: isd_p, ied_p, jsd_p, jed_p, is_n, ie_n, js_n, je_n
2338  integer, intent(IN) :: isu, ieu, jsu, jeu
2339  integer, intent(IN) :: istag, jstag, npx, npy, npz, r, nestupdate, upoff, nsponge
2340  real, intent(IN) :: var_nest(is_n:ie_n+istag,js_n:je_n+jstag,npz)
2341  real, intent(INOUT) :: var_coarse(isd_p:ied_p+istag,jsd_p:jed_p+jstag,npz)
2342  real, intent(IN) :: area(bd%isd:bd%ied,bd%jsd:bd%jed)
2343  real, intent(IN) :: dx(bd%isd:bd%ied,bd%jsd:bd%jed+1)
2344  real, intent(IN) :: dy(bd%isd:bd%ied+1,bd%jsd:bd%jed)
2345  logical, intent(IN) :: parent_proc, child_proc
2346  type(fv_atmos_type), pointer, intent(IN) :: parent_grid
2347  type(nest_domain_type), intent(INOUT) :: nest_domain
2348  integer, intent(IN) :: nest_level
2349 
2350  integer :: in, jn, ini, jnj, s, qr
2351  integer :: is_c, ie_c, js_c, je_c, is_f, ie_f, js_f, je_f
2352  integer :: istart, istop, jstart, jstop, ishift, jshift, j, i, k
2353  real :: val
2354  real, allocatable, dimension(:,:,:) :: coarse_dat_send
2355  real, allocatable :: coarse_dat_recv(:,:,:)
2356  integer :: position
2357 
2358  if (istag == 1 .and. jstag == 1) then
2359  position = corner
2360  else if (istag == 0 .and. jstag == 1) then
2361  position = north
2362  else if (istag == 1 .and. jstag == 0) then
2363  position = east
2364  else
2365  position = center
2366  end if
2367 
2368  !Note that *_c does not have values on the parent_proc.
2369  !Must use isu, etc. to get bounds of update region on parent.
2370  call mpp_get_f2c_index(nest_domain, is_c, ie_c, js_c, je_c, is_f, ie_f, js_f, je_f, nest_level=nest_level, position=position)
2371  if (child_proc) then
2372  allocate(coarse_dat_send(is_c:ie_c, js_c:je_c,npz))
2373  coarse_dat_send = -1200.
2374  endif
2375  allocate(coarse_dat_recv(isd_p:ied_p+istag, jsd_p:jed_p+jstag, npz))
2376 
2377  if (child_proc) then
2378  call fill_coarse_data_send(coarse_dat_send, var_nest, dx, dy, area, &
2379  bd, is_c, ie_c, js_c, je_c, is_f, js_f, is_n, ie_n, js_n, je_n, &
2380  npx, npy, npz, istag, jstag, r, nestupdate)
2381  endif
2382 
2383  call timing_on('COMM_TOTAL')
2384  call mpp_update_nest_coarse(field_in=coarse_dat_send, nest_domain=nest_domain, field_out=coarse_dat_recv, &
2385  nest_level=nest_level, position=position)
2386 
2387  if (allocated(coarse_dat_send)) then
2388  deallocate(coarse_dat_send)
2389  end if
2390 
2391  call timing_off('COMM_TOTAL')
2392 
2393  s = r/2 !rounds down (since r > 0)
2394  qr = r*upoff + nsponge - s
2395 
2396  if (parent_proc .and. .not. (ieu < isu .or. jeu < jsu)) then
2397  call fill_var_coarse(var_coarse, coarse_dat_recv, isd_p, ied_p, jsd_p, jed_p, &
2398  isu, ieu, jsu, jeu, npx, npy, npz, istag, jstag, nestupdate, parent_grid)
2399  endif
2400 
2401  if (allocated(coarse_dat_recv)) deallocate(coarse_dat_recv)
2402 
2403  end subroutine update_coarse_grid_mpp
2404 
2405  subroutine fill_coarse_data_send(coarse_dat_send, var_nest, dx, dy, area, &
2406  bd, is_c, ie_c, js_c, je_c, is_f, js_f, is_n, ie_n, js_n, je_n, &
2407  npx, npy, npz, istag, jstag, r, nestupdate)
2408  type(fv_grid_bounds_type), intent(IN) :: bd
2409  integer, intent(IN) :: is_c, ie_c, js_c, je_c, is_n, ie_n, js_n, je_n
2410  integer, intent(IN) :: is_f, js_f
2411  integer, intent(IN) :: istag, jstag
2412  integer, intent(IN) :: npx, npy, npz, r, nestupdate
2413  real, intent(INOUT) :: coarse_dat_send(is_c:ie_c,js_c:je_c,npz)
2414  real, intent(IN) :: var_nest(is_n:ie_n+istag,js_n:je_n+jstag,npz)
2415  real, intent(IN) :: area(bd%isd:bd%ied,bd%jsd:bd%jed)
2416  real, intent(IN) :: dx(bd%isd:bd%ied,bd%jsd:bd%jed+1)
2417  real, intent(IN) :: dy(bd%isd:bd%ied+1,bd%jsd:bd%jed)
2418  integer :: in, jn, ini, jnj, k, j, i
2419  real :: val
2420 
2421 
2422  if (istag == 0 .and. jstag == 0) then
2423  select case (nestupdate)
2424  case (1,2,6,7,8)
2425 
2426 !$OMP parallel do default(none) shared(npz,js_c,je_c,is_c,ie_c,js_f,is_f,coarse_dat_send,var_nest,area,r) private(in,jn,val)
2427  do k=1,npz
2428  jn = js_f
2429  do j=js_c,je_c
2430  in = is_f
2431  do i=is_c,ie_c
2432 
2433  val = 0.
2434  do jnj=jn,jn+r-1
2435  do ini=in,in+r-1
2436  val = val + var_nest(ini,jnj,k)*area(ini,jnj)
2437  end do
2438  end do
2439  coarse_dat_send(i,j,k) = val !divide area on coarse grid
2440 
2441  in = in + r
2442  end do
2443  jn = jn + r
2444  end do
2445  end do
2446 
2447  end select
2448  else if (istag == 0 .and. jstag > 0) then
2449 
2450  select case (nestupdate)
2451  case (1,6,7,8)
2452 
2453 !$OMP parallel do default(none) shared(npz,js_c,je_c,is_c,ie_c,js_f,is_f,coarse_dat_send,var_nest,dx,r) private(in,jn,val)
2454  do k=1,npz
2455  jn = js_f
2456  do j=js_c,je_c!+1
2457  in = is_f
2458  do i=is_c,ie_c
2459 
2460  val = 0.
2461  do ini=in,in+r-1
2462  val = val + var_nest(ini,jn,k)*dx(ini,jn)
2463  end do
2464  coarse_dat_send(i,j,k) = val
2465 
2466  in = in + r
2467  end do
2468  jn = jn + r
2469  end do
2470  end do
2471 
2472  case default
2473 
2474  call mpp_error(fatal, 'nestupdate type not implemented')
2475 
2476  end select
2477 
2478  else if (istag > 0 .and. jstag == 0) then
2479  select case (nestupdate)
2480 
2481  case (1,6,7,8) !averaging update; in-line average for face-averaged values instead of areal average
2482 
2483 !$OMP parallel do default(none) shared(npz,js_c,je_c,is_c,ie_c,js_f,is_f,coarse_dat_send,var_nest,dy,r) private(in,jn,val)
2484  do k=1,npz
2485  jn = js_f
2486  do j=js_c,je_c
2487  in = is_f
2488  do i=is_c,ie_c!+1
2489 
2490  val = 0.
2491  do jnj=jn,jn+r-1
2492  val = val + var_nest(in,jnj,k)*dy(in,jnj)
2493  end do
2494  coarse_dat_send(i,j,k) = val
2495 
2496  in = in + r
2497  end do
2498  jn = jn + r
2499  end do
2500  end do
2501 
2502  case default
2503 
2504  call mpp_error(fatal, 'nestupdate type not implemented')
2505 
2506  end select
2507 
2508  else
2509 
2510  call mpp_error(fatal, "Cannot have both nonzero istag and jstag.")
2511 
2512  endif
2513 
2514 
2515  end subroutine fill_coarse_data_send
2516 
2517  subroutine fill_var_coarse(var_coarse, coarse_dat_recv, isd_p, ied_p, jsd_p, jed_p, &
2518  isu, ieu, jsu, jeu, npx, npy, npz, istag, jstag, nestupdate, parent_grid)
2520  !This routine assumes the coarse and nested grids are properly
2521  ! aligned, and that in particular for odd refinement ratios all
2522  ! coarse-grid cells (faces) coincide with nested-grid cells (faces)
2523 
2524  integer, intent(IN) :: isd_p, ied_p, jsd_p, jed_p
2525  integer, intent(IN) :: isu, ieu, jsu, jeu
2526  integer, intent(IN) :: istag, jstag
2527  integer, intent(IN) :: npx, npy, npz, nestupdate
2528  real, intent(INOUT) :: var_coarse(isd_p:ied_p+istag,jsd_p:jed_p+jstag,npz)
2529  real, intent(INOUT) :: coarse_dat_recv(isd_p:ied_p+istag,jsd_p:jed_p+jstag,npz)
2530  type(fv_atmos_type), intent(IN) :: parent_grid
2531 
2532  integer :: i, j, k
2533 
2534  if (istag == 0 .and. jstag == 0) then
2535 
2536  select case (nestupdate)
2537  case (1,2,6,7,8) ! 1 = Conserving update on all variables; 2 = conserving update for cell-centered values; 6 = conserving remap-update
2538 
2539 !$OMP parallel do default(none) shared(npz,jsu,jeu,isu,ieu,coarse_dat_recv,parent_grid,var_coarse)
2540  do k=1,npz
2541  do j=jsu,jeu
2542  do i=isu,ieu
2543  var_coarse(i,j,k) = coarse_dat_recv(i,j,k)*parent_grid%gridstruct%rarea(i,j)
2544  end do
2545  end do
2546  end do
2547 
2548 
2549  case default
2550 
2551  call mpp_error(fatal, 'nestupdate type not implemented')
2552 
2553 
2554  end select
2555 
2556  else if (istag == 0 .and. jstag > 0) then
2557 
2558 
2559  select case (nestupdate)
2560  case (1,6,7,8)
2561 
2562 !$OMP parallel do default(none) shared(npz,jsu,jeu,isu,ieu,coarse_dat_recv,parent_grid,var_coarse)
2563  do k=1,npz
2564  do j=jsu,jeu+1
2565  do i=isu,ieu
2566  var_coarse(i,j,k) = coarse_dat_recv(i,j,k)*parent_grid%gridstruct%rdx(i,j)
2567  end do
2568  end do
2569  end do
2570 
2571  case default
2572 
2573  call mpp_error(fatal, 'nestupdate type not implemented')
2574 
2575  end select
2576 
2577  else if (istag > 0 .and. jstag == 0) then
2578 
2579  select case (nestupdate)
2580  case (1,6,7,8) !averaging update; in-line average for face-averaged values instead of areal average
2581 
2582 !$OMP parallel do default(none) shared(npz,jsu,jeu,isu,ieu,coarse_dat_recv,parent_grid,var_coarse)
2583  do k=1,npz
2584  do j=jsu,jeu
2585  do i=isu,ieu+1
2586  var_coarse(i,j,k) = coarse_dat_recv(i,j,k)*parent_grid%gridstruct%rdy(i,j)
2587  end do
2588  end do
2589  end do
2590 
2591  case default
2592 
2593  call mpp_error(fatal, 'nestupdate type not implemented')
2594 
2595  end select
2596 
2597  end if
2598 
2599 
2600  end subroutine fill_var_coarse
2601 
2602  subroutine update_coarse_grid_mpp_vector(u_coarse, v_coarse, u_nest, v_nest, nest_domain, dx, dy, area, &
2603  bd, isd_p, ied_p, jsd_p, jed_p, is_n, ie_n, js_n, je_n, &
2604  isu, ieu, jsu, jeu, npx, npy, npz, istag_u, jstag_u, istag_v, jstag_v, &
2605  r, nestupdate, upoff, nsponge, &
2606  parent_proc, child_proc, parent_grid, nest_level, flags, gridtype)
2608  !This routine assumes the coarse and nested grids are properly
2609  ! aligned, and that in particular for odd refinement ratios all
2610  ! coarse-grid cells (faces) coincide with nested-grid cells (faces)
2611 
2612  type(fv_grid_bounds_type), intent(IN) :: bd
2613  integer, intent(IN) :: isd_p, ied_p, jsd_p, jed_p, is_n, ie_n, js_n, je_n
2614  integer, intent(IN) :: isu, ieu, jsu, jeu
2615  integer, intent(IN) :: istag_u, jstag_u, istag_v, jstag_v
2616  integer, intent(IN) :: npx, npy, npz, r, nestupdate, upoff, nsponge
2617  real, intent(IN) :: u_nest(is_n:ie_n+istag_u,js_n:je_n+jstag_u,npz)
2618  real, intent(INOUT) :: u_coarse(isd_p:ied_p+istag_u,jsd_p:jed_p+jstag_u,npz)
2619  real, intent(IN) :: v_nest(is_n:ie_n+istag_v,js_n:je_n+jstag_v,npz)
2620  real, intent(INOUT) :: v_coarse(isd_p:ied_p+istag_v,jsd_p:jed_p+jstag_v,npz)
2621  real, intent(IN) :: area(bd%isd:bd%ied,bd%jsd:bd%jed)
2622  real, intent(IN) :: dx(bd%isd:bd%ied,bd%jsd:bd%jed+1)
2623  real, intent(IN) :: dy(bd%isd:bd%ied+1,bd%jsd:bd%jed)
2624  logical, intent(IN) :: parent_proc, child_proc
2625  type(fv_atmos_type), intent(INOUT) :: parent_grid
2626  integer, intent(IN) :: nest_level
2627  type(nest_domain_type), intent(INOUT) :: nest_domain
2628  integer, optional, intent(IN) :: flags, gridtype
2629 
2630  integer :: s, qr
2631  integer :: is_cx, ie_cx, js_cx, je_cx, is_fx, ie_fx, js_fx, je_fx
2632  integer :: is_cy, ie_cy, js_cy, je_cy, is_fy, ie_fy, js_fy, je_fy
2633  integer :: istart, istop, jstart, jstop, ishift, jshift, j, i, k
2634  real :: val
2635  real, allocatable, dimension(:,:,:) :: coarse_dat_send_u, coarse_dat_send_v
2636  real, allocatable :: coarse_dat_recv_u(:,:,:), coarse_dat_recv_v(:,:,:)
2637  integer :: position_x, position_y
2638 
2639  call get_vector_position(position_x, position_y, gridtype)
2640 
2641  call mpp_get_f2c_index(nest_domain, is_cx, ie_cx, js_cx, je_cx, is_fx, ie_fx, js_fx, je_fx, &
2642  nest_level=nest_level, position=position_x)
2643  call mpp_get_f2c_index(nest_domain, is_cy, ie_cy, js_cy, je_cy, is_fy, ie_fy, js_fy, je_fy, &
2644  nest_level=nest_level, position=position_y)
2645  if (child_proc) then
2646  allocate(coarse_dat_send_u(is_cx:ie_cx, js_cx:je_cx,npz))
2647  allocate(coarse_dat_send_v(is_cy:ie_cy, js_cy:je_cy,npz))
2648  coarse_dat_send_u = -1200.
2649  coarse_dat_send_v = -1200.
2650  endif
2651  allocate(coarse_dat_recv_u(isd_p:ied_p+istag_u, jsd_p:jed_p+jstag_u, npz))
2652  allocate(coarse_dat_recv_v(isd_p:ied_p+istag_v, jsd_p:jed_p+jstag_v, npz))
2653 
2654  if (child_proc) then
2655  call fill_coarse_data_send(coarse_dat_send_u, u_nest, dx, dy, area, &
2656  bd, is_cx, ie_cx, js_cx, je_cx, is_fx, js_fx, is_n, ie_n, js_n, je_n, &
2657  npx, npy, npz, istag_u, jstag_u, r, nestupdate)
2658  call fill_coarse_data_send(coarse_dat_send_v, v_nest, dx, dy, area, &
2659  bd, is_cy, ie_cy, js_cy, je_cy, is_fy, js_fy, is_n, ie_n, js_n, je_n, &
2660  npx, npy, npz, istag_v, jstag_v, r, nestupdate)
2661  endif
2662 
2663  call timing_on('COMM_TOTAL')
2664  call mpp_update_nest_coarse(coarse_dat_send_u, coarse_dat_send_v, nest_domain, coarse_dat_recv_u, &
2665  coarse_dat_recv_v, nest_level, flags, gridtype)
2666 
2667  if (allocated(coarse_dat_send_u)) deallocate(coarse_dat_send_u)
2668  if (allocated(coarse_dat_send_v)) deallocate(coarse_dat_send_v)
2669 
2670  call timing_off('COMM_TOTAL')
2671 
2672  s = r/2 !rounds down (since r > 0)
2673  qr = r*upoff + nsponge - s
2674 
2675  if (parent_proc .and. .not. (ieu < isu .or. jeu < jsu)) then
2676  call fill_var_coarse(u_coarse, coarse_dat_recv_u, isd_p, ied_p, jsd_p, jed_p, &
2677  isu, ieu, jsu, jeu, npx, npy, npz, istag_u, jstag_u, nestupdate, parent_grid)
2678  endif
2679  if (parent_proc .and. .not. (ieu < isu .or. jeu < jsu)) then
2680  call fill_var_coarse(v_coarse, coarse_dat_recv_v, isd_p, ied_p, jsd_p, jed_p, &
2681  isu, ieu, jsu, jeu, npx, npy, npz, istag_v, jstag_v, nestupdate, parent_grid)
2682  endif
2683 
2684  if (allocated(coarse_dat_recv_u)) deallocate(coarse_dat_recv_u)
2685  if (allocated(coarse_dat_recv_v)) deallocate(coarse_dat_recv_v)
2686 
2687  end subroutine update_coarse_grid_mpp_vector
2688 
2689 end module boundary_mod
The module &#39;fv_mp_mod&#39; is a single program multiple data (SPMD) parallel decompostion/communication m...
Definition: fv_mp_mod.F90:24
subroutine timing_off(blk_name)
The subroutine &#39;timing_off&#39; stops a timer.
Definition: fv_timing.F90:180
subroutine nested_grid_bc_mpp_3d(var_nest, var_coarse, nest_domain, ind, wt, istag, jstag, npx, npy, npz, bd, isg, ieg, jsg, jeg, nest_level, nstep_in, nsplit_in, proc_in)
Definition: boundary.F90:624
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)
Definition: boundary.F90:526
The interface&#39;update_coarse_grid_mpp&#39;contains subroutines that fetch data from the nested grid and in...
Definition: boundary.F90:122
subroutine get_vector_position(position_x, position_y, gridtype)
Definition: boundary.F90:837
subroutine nested_grid_bc_send_vector(u_coarse, v_coarse, nest_domain, nest_level, flags, gridtype)
Definition: boundary.F90:1797
subroutine update_coarse_grid_mpp_2d(var_coarse, var_nest, nest_domain, dx, dy, area, bd, 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, nest_level)
Definition: boundary.F90:2291
subroutine init_buffer(nest_domain, wbuffer, sbuffer, ebuffer, nbuffer, npz, nest_level, position)
Definition: boundary.F90:866
The module &#39;fv_timing&#39; contains FV3 timers.
Definition: fv_timing.F90:24
The module &#39;boundary&#39; contains utility routines for grid nesting and boundary conditions.
Definition: boundary.F90:25
subroutine init_nest_bc_type(nest_domain, nest_BC_buffers, npz, nest_level, position)
Definition: boundary.F90:1817
The module &#39;fv_arrays&#39; contains the &#39;fv_atmos_type&#39; and associated datatypes.
Definition: fv_arrays.F90:24
subroutine nested_grid_bc_2d(var_nest, var_coarse, ind, wt, istag, jstag, npx, npy, bd, isg, ieg, jsg, jeg, nstep_in, nsplit_in)
Definition: boundary.F90:1444
interface &#39;nested_grid_BC&#39; includes subroutines &#39;nested_grid_BC_2d&#39; and &#39;nested_grid_BC_3d&#39; that fetc...
Definition: boundary.F90:89
subroutine nested_grid_bc_recv_vector(nest_domain, npz, bd, nest_BC_u_buffers, nest_BC_v_buffers, nest_level, flags, gridtype)
Definition: boundary.F90:1901
subroutine nested_grid_bc_2d_mpp(var_nest, var_coarse, nest_domain, ind, wt, istag, jstag, npx, npy, bd, isg, ieg, jsg, jeg, nest_level, nstep_in, nsplit_in, proc_in)
Definition: boundary.F90:1239
subroutine nested_grid_bc_send_scalar(var_coarse, nest_domain, istag, jstag, nest_level)
The subroutine &#39;nested_grid_BC_send&#39; sends coarse-grid data to create boundary conditions.
Definition: boundary.F90:1724
The interface &#39;fill_nested_grid&#39; includes subroutines &#39;fill_nested_grid_2d&#39; and &#39;fill_nested_grid_3d&#39;...
Definition: boundary.F90:113
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: boundary.F90:1579
subroutine timing_on(blk_name)
The subroutine &#39;timing_on&#39; starts a timer.
Definition: fv_timing.F90:116
subroutine nested_grid_bc_mpp_send_3d(var_coarse, nest_domain, istag, jstag, nest_level)
Definition: boundary.F90:1146
subroutine update_coarse_grid_mpp_vector(u_coarse, v_coarse, u_nest, v_nest, nest_domain, dx, dy, area, bd, isd_p, ied_p, jsd_p, jed_p, is_n, ie_n, js_n, je_n, isu, ieu, jsu, jeu, npx, npy, npz, istag_u, jstag_u, istag_v, jstag_v, r, nestupdate, upoff, nsponge, parent_proc, child_proc, parent_grid, nest_level, flags, gridtype)
Definition: boundary.F90:2607
subroutine nested_grid_bc_mpp_send_2d(var_coarse, nest_domain, istag, jstag, nest_level)
Definition: boundary.F90:1192
subroutine nested_grid_bc_mpp_3d_vector(u_nest, v_nest, u_coarse, v_coarse, nest_domain, ind_u, ind_v, wt_u, wt_v, istag_u, jstag_u, istag_v, jstag_v, npx, npy, npz, bd, isg, ieg, jsg, jeg, nest_level, nstep_in, nsplit_in, proc_in, flags, gridtype)
Definition: boundary.F90:917
subroutine nested_grid_bc_recv_scalar(nest_domain, istag, jstag, npz, bd, nest_BC_buffers, nest_level)
Definition: boundary.F90:1756
subroutine, public nested_grid_bc_apply_intt(var_nest, istag, jstag, npx, npy, npz, bd, step, split, BC, bctype)
The subroutine &#39;nested_grid_BC_apply_intT&#39; performs linear interpolation or extrapolation in time for...
Definition: boundary.F90:2171
subroutine, public nested_grid_bc_save_proc(nest_domain, ind, wt, istag, jstag, npx, npy, npz, bd, nest_BC, nest_BC_buffers, pd_in)
The subroutine &#39;nested_grid_BC_save_proc&#39; saves data received by &#39;nested_grid_BC_recv&#39; into the datat...
Definition: boundary.F90:1936
subroutine, public extrapolation_bc(q, istag, jstag, npx, npy, bd, pd_in, debug_in)
The subroutine &#39;extrapolation_BC&#39; performs linear extrapolation into the halo region.
Definition: boundary.F90:133
subroutine fill_coarse_data_send(coarse_dat_send, var_nest, dx, dy, area, bd, is_c, ie_c, js_c, je_c, is_f, js_f, is_n, ie_n, js_n, je_n, npx, npy, npz, istag, jstag, r, nestupdate)
Definition: boundary.F90:2408
subroutine fill_var_coarse(var_coarse, coarse_dat_recv, isd_p, ied_p, jsd_p, jed_p, isu, ieu, jsu, jeu, npx, npy, npz, istag, jstag, nestupdate, parent_grid)
Definition: boundary.F90:2519
subroutine update_coarse_grid_mpp(var_coarse, var_nest, nest_domain, dx, dy, area, bd, 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, nest_level)
Definition: boundary.F90:2331
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)
Definition: boundary.F90:461