63 use constants_mod
, only: grav
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
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
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
139 integer :: i,j, istart, iend, jstart, jend
142 integer :: is, ie, js, je
143 integer :: isd, ied, jsd, jed
155 iend = min(ied,npx-1)
157 jend = min(jed,npy-1)
160 if (
present(pd_in))
then 166 if (
present(debug_in))
then 176 do j = jstart,jend+jstag
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 182 q(i,j) =
real(2-i)*q(1,j) -
real(1-i)*q(2,j)
190 do j = jstart,jend+jstag
193 q(i,j) =
real(2-i)*q(1,j) -
real(1-i)*q(2,j)
207 do i = istart,iend+istag
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 212 q(i,j) =
real(2-j)*q(i,1) -
real(1-j)*q(i,2)
221 do i = istart,iend+istag
223 q(i,j) =
real(2-j)*q(i,1) -
real(1-j)*q(i,2)
232 if (ie == npx - 1)
then 236 do j=jstart,jend+jstag
237 do i=ie+1+istag,ied+istag
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 243 q(i,j) =
real(i - (ie+istag-1))*q(ie+istag,j) +
real((ie+istag) - i)*q(ie+istag-1,j)
251 do j=jstart,jend+jstag
252 do i=ie+1+istag,ied+istag
254 q(i,j) =
real(i - (ie+istag-1))*q(ie+istag,j) +
real((ie+istag) - i)*q(ie+istag-1,j)
263 if (je == npy - 1)
then 267 do j=je+1+jstag,jed+jstag
268 do i=istart,iend+istag
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 274 q(i,j) =
real(j - (je+jstag-1))*q(i,je+jstag) +
real((je+jstag) - j)*q(i,je+jstag-1)
282 do j=je+1+jstag,jed+jstag
283 do i=istart,iend+istag
285 q(i,j) =
real(j - (je+jstag-1))*q(i,je+jstag) +
real((je+jstag) - j)*q(i,je+jstag-1)
297 if (is == 1 .and. js == 1)
then 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)
307 q(i,j) = 0.5*(
real(2-i)*q(1,j) -
real(1-i)*q(2,j) )
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)
314 q(i,j) = q(i,j) + 0.5*(
real(2-j)*q(i,1) -
real(1-j)*q(i,2))
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) )
335 if (is == 1 .and. je == npy-1)
then 339 do j=je+1+jstag,jed+jstag
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)
345 q(i,j) = 0.5*(
real(2-i)*q(1,j) -
real(1-i)*q(2,j) )
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)
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) )
362 do j=je+1+jstag,jed+jstag
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) )
375 if (ie == npx-1 .and. je == npy-1)
then 379 do j=je+1+jstag,jed+jstag
380 do i=ie+1+istag,ied+istag
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)
387 q(i,j) = 0.5*(
real(i - (ie+istag-1))*q(ie+istag,j) +
real((ie+istag) - i)*q(ie+istag-1,j))
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)
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) )
402 do j=je+1+jstag,jed+jstag
403 do i=ie+1+istag,ied+istag
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) )
415 if (ie == npx-1 .and. js == 1)
then 420 do i=ie+1+istag,ied+istag
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)
427 q(i,j) = 0.5*(
real(i - (ie+istag-1))*q(ie+istag,j) +
real((ie+istag) - i)*q(ie+istag-1,j))
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)
434 q(i,j) = q(i,j) + 0.5*(
real(2-j)*q(i,1) -
real(1-j)*q(i,2))
444 do i=ie+1+istag,ied+istag
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) )
460 isg, ieg, jsg, jeg, bd, istart_in, iend_in, jstart_in, jend_in)
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
470 integer :: i,j, ic, jc
471 integer :: istart, iend, jstart, jend
473 integer :: is, ie, js, je
474 integer :: isd, ied, jsd, jed
485 if (
present(istart_in))
then 490 if (
present(iend_in))
then 496 if (
present(jstart_in))
then 501 if (
present(jend_in))
then 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)
525 isg, ieg, jsg, jeg, npz, bd, istart_in, iend_in, jstart_in, jend_in)
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
535 integer :: i,j, ic, jc, k
536 integer :: istart, iend, jstart, jend
538 integer :: is, ie, js, je
539 integer :: isd, ied, jsd, jed
550 if (
present(istart_in))
then 555 if (
present(iend_in))
then 561 if (
present(jstart_in))
then 566 if (
present(jend_in))
then 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)
623 npx, npy, npz, bd, isg, ieg, jsg, jeg, nest_level, nstep_in, nsplit_in, proc_in)
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
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(:,:,:)
645 integer :: i,j, ic, jc, istart, iend, k
650 integer :: is, ie, js, je
651 integer :: isd, ied, jsd, jed
662 if (
PRESENT(proc_in))
then 668 if (istag == 1 .and. jstag == 1)
then 670 else if (istag == 0 .and. jstag == 1)
then 672 else if (istag == 1 .and. jstag == 0)
then 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)
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))
690 allocate(wbuffer(1,1,1))
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))
697 allocate(ebuffer(1,1,1))
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))
704 allocate(sbuffer(1,1,1))
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))
711 allocate(nbuffer(1,1,1))
717 call mpp_update_nest_fine(var_coarse, nest_domain, wbuffer, sbuffer, ebuffer, nbuffer, &
718 nest_level=nest_level, position=position)
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)
751 if (ie == npx-1)
then 760 do i=istart,iend+istag
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)
777 if (ie == npx-1)
then 781 do i=npx+istag,ied+istag
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)
797 if (je == npy-1)
then 805 if (ie == npx-1)
then 813 do j=npy+jstag,jed+jstag
814 do i=istart,iend+istag
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)
832 deallocate(wbuffer, ebuffer, sbuffer, nbuffer)
837 integer,
intent(OUT) :: position_x, position_y
838 integer,
optional,
intent(IN) :: gridtype
840 integer :: grid_offset_type
842 grid_offset_type = agrid
843 if(
present(gridtype)) grid_offset_type = gridtype
845 select case(grid_offset_type)
859 call mpp_error(fatal,
"get_vector_position: invalid value of gridtype")
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
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)
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))
886 allocate(wbuffer(1,1,1))
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))
893 allocate(ebuffer(1,1,1))
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))
900 allocate(sbuffer(1,1,1))
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))
907 allocate(nbuffer(1,1,1))
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, &
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
934 real,
allocatable :: wbufferx(:,:,:), wbuffery(:,:,:)
935 real,
allocatable :: ebufferx(:,:,:), ebuffery(:,:,:)
936 real,
allocatable :: sbufferx(:,:,:), sbuffery(:,:,:)
937 real,
allocatable :: nbufferx(:,:,:), nbuffery(:,:,:)
939 integer :: i,j, ic, jc, istart, iend, k
941 integer :: position_x, position_y
944 integer :: is, ie, js, je
945 integer :: isd, ied, jsd, jed
956 if (
PRESENT(proc_in))
then 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)
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)
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)
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)
1016 if (ie == npx-1)
then 1025 do i=istart,iend+istag_u
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)
1039 do i=istart,iend+istag_v
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)
1056 if (ie == npx-1)
then 1059 do j=jsd,jed+jstag_u
1060 do i=npx+istag_u,ied+istag_u
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)
1073 do j=jsd,jed+jstag_v
1074 do i=npx+istag_v,ied+istag_v
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)
1090 if (je == npy-1)
then 1098 if (ie == npx-1)
then 1106 do j=npy+jstag_u,jed+jstag_u
1107 do i=istart,iend+istag_u
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)
1120 do j=npy+jstag_v,jed+jstag_v
1121 do i=istart,iend+istag_v
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)
1139 deallocate(wbufferx, ebufferx, sbufferx, nbufferx)
1140 deallocate(wbuffery, ebuffery, sbuffery, nbuffery)
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
1152 real,
allocatable :: wbuffer(:,:,:)
1153 real,
allocatable :: ebuffer(:,:,:)
1154 real,
allocatable :: sbuffer(:,:,:)
1155 real,
allocatable :: nbuffer(:,:,:)
1157 integer :: i,j, ic, jc, istart, iend, k
1162 if (istag == 1 .and. jstag == 1)
then 1164 else if (istag == 0 .and. jstag == 1)
then 1166 else if (istag == 1 .and. jstag == 0)
then 1173 allocate(wbuffer(1,1,1))
1175 allocate(ebuffer(1,1,1))
1177 allocate(sbuffer(1,1,1))
1179 allocate(nbuffer(1,1,1))
1183 call mpp_update_nest_fine(var_coarse, nest_domain, wbuffer, sbuffer, ebuffer, nbuffer, nest_level=nest_level, position=position)
1187 deallocate(wbuffer, ebuffer, sbuffer, nbuffer)
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
1198 real,
allocatable :: wbuffer(:,:)
1199 real,
allocatable :: ebuffer(:,:)
1200 real,
allocatable :: sbuffer(:,:)
1201 real,
allocatable :: nbuffer(:,:)
1203 integer :: i,j, ic, jc, istart, iend, k
1208 if (istag == 1 .and. jstag == 1)
then 1210 else if (istag == 0 .and. jstag == 1)
then 1212 else if (istag == 1 .and. jstag == 0)
then 1219 allocate(wbuffer(1,1))
1221 allocate(ebuffer(1,1))
1223 allocate(sbuffer(1,1))
1225 allocate(nbuffer(1,1))
1229 call mpp_update_nest_fine(var_coarse, nest_domain, wbuffer, sbuffer, ebuffer, nbuffer, nest_level=nest_level, position=position)
1233 deallocate(wbuffer, ebuffer, sbuffer, nbuffer)
1238 npx, npy, bd, isg, ieg, jsg, jeg, nest_level, nstep_in, nsplit_in, proc_in)
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
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(:,:)
1260 integer :: i,j, ic, jc, istart, iend, k
1266 integer :: is, ie, js, je
1267 integer :: isd, ied, jsd, jed
1278 if (
PRESENT(proc_in))
then 1284 if (
PRESENT(nest_level))
then 1288 if (istag == 1 .and. jstag == 1)
then 1290 else if (istag == 0 .and. jstag == 1)
then 1292 else if (istag == 1 .and. jstag == 0)
then 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)
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))
1310 allocate(wbuffer(1,1))
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))
1317 allocate(ebuffer(1,1))
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))
1324 allocate(sbuffer(1,1))
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))
1331 allocate(nbuffer(1,1))
1336 call mpp_update_nest_fine(var_coarse, nest_domain, wbuffer, sbuffer, ebuffer, nbuffer, nest_level=nl, position=position)
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)
1366 if (ie == npx-1)
then 1373 do i=istart,iend+istag
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)
1389 if (ie == npx-1)
then 1391 do i=npx+istag,ied+istag
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)
1406 if (je == npy-1)
then 1414 if (ie == npx-1)
then 1420 do j=npy+jstag,jed+jstag
1421 do i=istart,iend+istag
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)
1438 deallocate(wbuffer, ebuffer, sbuffer, nbuffer)
1443 npx, npy, bd, isg, ieg, jsg, jeg, nstep_in, nsplit_in)
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
1453 integer :: nstep, nsplit
1455 integer :: i,j, ic, jc, istart, iend
1457 integer :: is, ie, js, je
1458 integer :: isd, ied, jsd, jed
1469 if ( .not.
present(nstep_in) .or. .not.
present(nsplit_in) )
then 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)
1502 if (ie == npx-1)
then 1509 do i=istart,iend+istag
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)
1525 if (ie == npx-1)
then 1527 do i=npx+istag,ied+istag
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)
1542 if (je == npy-1)
then 1550 if (ie == npx-1)
then 1557 do j=npy+jstag,jed+jstag
1558 do i=istart,iend+istag
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)
1578 npx, npy, npz, bd, isg, ieg, jsg, jeg, nstep_in, nsplit_in)
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
1588 integer :: nstep, nsplit
1590 integer :: i,j, ic, jc, istart, iend, k
1592 integer :: is, ie, js, je
1593 integer :: isd, ied, jsd, jed
1604 if ( .not.
present(nstep_in) .or. .not.
present(nsplit_in) )
then 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)
1640 if (ie == npx-1)
then 1649 do i=istart,iend+istag
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)
1666 if (ie == npx-1)
then 1670 do i=npx+istag,ied+istag
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)
1686 if (je == npy-1)
then 1694 if (ie == npx-1)
then 1702 do j=npy+jstag,jed+jstag
1703 do i=istart,iend+istag
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)
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
1732 real :: wbuffer(1,1,1)
1733 real :: ebuffer(1,1,1)
1734 real :: sbuffer(1,1,1)
1735 real :: nbuffer(1,1,1)
1738 if (istag == 1 .and. jstag == 1)
then 1740 else if (istag == 0 .and. jstag == 1)
then 1742 else if (istag == 1 .and. jstag == 0)
then 1749 call mpp_update_nest_fine(var_coarse, nest_domain, wbuffer, sbuffer, ebuffer, nbuffer, nest_level=nest_level, position=position)
1755 bd, nest_BC_buffers, nest_level)
1758 type(nest_domain_type),
intent(INOUT) :: nest_domain
1759 integer,
intent(IN) :: istag, jstag, npz
1760 integer,
intent(IN) :: nest_level
1764 real,
dimension(bd%isd:bd%ied+istag,bd%jsd:bd%jed+jstag,npz) :: var_coarse_dummy
1775 if (istag == 1 .and. jstag == 1)
then 1777 else if (istag == 0 .and. jstag == 1)
then 1779 else if (istag == 1 .and. jstag == 0)
then 1785 if (.not.
allocated(nest_bc_buffers%west_t1) )
then 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)
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
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)
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)
1816 subroutine init_nest_bc_type(nest_domain, nest_BC_buffers, npz, nest_level, position)
1817 type(nest_domain_type),
intent(INOUT) :: nest_domain
1819 integer,
intent(IN) :: npz, position, nest_level
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
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)
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))
1843 nest_bc_buffers%west_t1(i,j,k) = 1.e24
1848 allocate(nest_bc_buffers%west_t1(1,1,1))
1849 nest_bc_buffers%west_t1(1,1,1) = 1.e24
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))
1858 nest_bc_buffers%east_t1(i,j,k) = 1.e24
1863 allocate(nest_bc_buffers%east_t1(1,1,1))
1864 nest_bc_buffers%east_t1(1,1,1) = 1.e24
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))
1873 nest_bc_buffers%south_t1(i,j,k) = 1.e24
1878 allocate(nest_bc_buffers%south_t1(1,1,1))
1879 nest_bc_buffers%south_t1(1,1,1) = 1.e24
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))
1888 nest_bc_buffers%north_t1(i,j,k) = 1.e24
1893 allocate(nest_bc_buffers%north_t1(1,1,1))
1894 nest_bc_buffers%north_t1(1,1,1) = 1.e24
1900 subroutine nested_grid_bc_recv_vector(nest_domain, npz, bd, nest_BC_u_buffers, nest_BC_v_buffers, nest_level, flags, gridtype)
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
1909 real,
dimension(1,1,npz) :: u_coarse_dummy, v_coarse_dummy
1912 integer :: position_x, position_y
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)
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)
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)
1935 npx, npy, npz, bd, nest_BC, nest_BC_buffers, pd_in)
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
1948 real,
dimension(bd%isd:bd%ied+istag,bd%jsd:bd%jed+jstag,npz) :: var_coarse_dummy
1950 real,
dimension(:,:,:),
pointer :: var_east, var_west, var_south, var_north
1951 real,
dimension(:,:,:),
pointer :: buf_east, buf_west, buf_south, buf_north
1956 integer :: i,j, k, ic, jc, istart, iend
1957 logical :: process, pd = .false.
1959 integer :: is, ie, js, je
1960 integer :: isd, ied, jsd, jed
1972 if (
present(pd_in))
then 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
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
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)
2017 var_west(i,j,k) = max(var_west(i,j,k), 0.5*nest_bc%west_t0(i,j,k))
2033 if (ie == npx-1)
then 2042 do i=istart,iend+istag
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)
2062 do i=istart,iend+istag
2064 var_south(i,j,k) = max(var_south(i,j,k), 0.5*nest_bc%south_t0(i,j,k))
2074 if (ie == npx-1 )
then 2079 do i=npx+istag,ied+istag
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)
2099 do i=npx+istag,ied+istag
2101 var_east(i,j,k) = max(var_east(i,j,k), 0.5*nest_bc%east_t0(i,j,k))
2110 if (je == npy-1 )
then 2118 if (ie == npx-1)
then 2126 do j=npy+jstag,jed+jstag
2127 do i=istart,iend+istag
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)
2146 do j=npy+jstag,jed+jstag
2147 do i=istart,iend+istag
2149 var_north(i,j,k) = max(var_north(i,j,k), 0.5*nest_bc%north_t0(i,j,k))
2169 npx, npy, npz, bd, step, split, &
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
2179 real,
pointer,
dimension(:,:,:) :: var_t0, var_t1
2181 integer :: i,j, istart, iend, k
2184 logical,
save :: printdiag = .true.
2186 integer :: is, ie, js, je
2187 integer :: isd, ied, jsd, jed
2200 var_t0 => bc%west_t0
2201 var_t1 => bc%west_t1
2206 var_nest(i,j,k) = (var_t0(i,j,k)*(split-step) + step*var_t1(i,j,k))*denom
2221 if (ie == npx-1)
then 2227 var_t0 => bc%south_t0
2228 var_t1 => bc%south_t1
2232 do i=istart,iend+istag
2234 var_nest(i,j,k) = (var_t0(i,j,k)*(split-step) + step*var_t1(i,j,k))*denom
2241 if (ie == npx-1 )
then 2242 var_t0 => bc%east_t0
2243 var_t1 => bc%east_t1
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
2256 if (je == npy-1 )
then 2264 if (ie == npx-1)
then 2270 var_t0 => bc%north_t0
2271 var_t1 => bc%north_t1
2274 do j=npy+jstag,jed+jstag
2275 do i=istart,iend+istag
2277 var_nest(i,j,k) = (var_t0(i,j,k)*(split-step) + step*var_t1(i,j,k))*denom
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)
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
2304 type(nest_domain_type),
intent(INOUT) :: nest_domain
2305 integer,
intent(IN) :: nest_level
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)
2314 ptr_nest = loc(var_nest)
2315 ptr_coarse = loc(var_coarse)
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 )
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)
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
2347 type(nest_domain_type),
intent(INOUT) :: nest_domain
2348 integer,
intent(IN) :: nest_level
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
2354 real,
allocatable,
dimension(:,:,:) :: coarse_dat_send
2355 real,
allocatable :: coarse_dat_recv(:,:,:)
2358 if (istag == 1 .and. jstag == 1)
then 2360 else if (istag == 0 .and. jstag == 1)
then 2362 else if (istag == 1 .and. jstag == 0)
then 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.
2375 allocate(coarse_dat_recv(isd_p:ied_p+istag, jsd_p:jed_p+jstag, npz))
2377 if (child_proc)
then 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)
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)
2387 if (
allocated(coarse_dat_send))
then 2388 deallocate(coarse_dat_send)
2394 qr = r*upoff + nsponge - s
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)
2401 if (
allocated(coarse_dat_recv))
deallocate(coarse_dat_recv)
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)
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
2422 if (istag == 0 .and. jstag == 0)
then 2423 select case (nestupdate)
2436 val = val + var_nest(ini,jnj,k)*area(ini,jnj)
2439 coarse_dat_send(i,j,k) = val
2448 else if (istag == 0 .and. jstag > 0)
then 2450 select case (nestupdate)
2462 val = val + var_nest(ini,jn,k)*dx(ini,jn)
2464 coarse_dat_send(i,j,k) = val
2474 call mpp_error(fatal,
'nestupdate type not implemented')
2478 else if (istag > 0 .and. jstag == 0)
then 2479 select case (nestupdate)
2492 val = val + var_nest(in,jnj,k)*dy(in,jnj)
2494 coarse_dat_send(i,j,k) = val
2504 call mpp_error(fatal,
'nestupdate type not implemented')
2510 call mpp_error(fatal,
"Cannot have both nonzero istag and jstag.")
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)
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)
2534 if (istag == 0 .and. jstag == 0)
then 2536 select case (nestupdate)
2543 var_coarse(i,j,k) = coarse_dat_recv(i,j,k)*parent_grid%gridstruct%rarea(i,j)
2551 call mpp_error(fatal,
'nestupdate type not implemented')
2556 else if (istag == 0 .and. jstag > 0)
then 2559 select case (nestupdate)
2566 var_coarse(i,j,k) = coarse_dat_recv(i,j,k)*parent_grid%gridstruct%rdx(i,j)
2573 call mpp_error(fatal,
'nestupdate type not implemented')
2577 else if (istag > 0 .and. jstag == 0)
then 2579 select case (nestupdate)
2586 var_coarse(i,j,k) = coarse_dat_recv(i,j,k)*parent_grid%gridstruct%rdy(i,j)
2593 call mpp_error(fatal,
'nestupdate type not implemented')
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)
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
2626 integer,
intent(IN) :: nest_level
2627 type(nest_domain_type),
intent(INOUT) :: nest_domain
2628 integer,
optional,
intent(IN) :: flags, gridtype
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
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
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.
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))
2654 if (child_proc)
then 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)
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)
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)
2667 if (
allocated(coarse_dat_send_u))
deallocate(coarse_dat_send_u)
2668 if (
allocated(coarse_dat_send_v))
deallocate(coarse_dat_send_v)
2673 qr = r*upoff + nsponge - s
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)
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)
2684 if (
allocated(coarse_dat_recv_u))
deallocate(coarse_dat_recv_u)
2685 if (
allocated(coarse_dat_recv_v))
deallocate(coarse_dat_recv_v)
The module 'fv_mp_mod' is a single program multiple data (SPMD) parallel decompostion/communication m...
subroutine timing_off(blk_name)
The subroutine 'timing_off' stops a timer.
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)
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)
The interface'update_coarse_grid_mpp'contains subroutines that fetch data from the nested grid and in...
subroutine get_vector_position(position_x, position_y, gridtype)
subroutine nested_grid_bc_send_vector(u_coarse, v_coarse, nest_domain, nest_level, flags, gridtype)
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)
subroutine init_buffer(nest_domain, wbuffer, sbuffer, ebuffer, nbuffer, npz, nest_level, position)
The module 'fv_timing' contains FV3 timers.
The module 'boundary' contains utility routines for grid nesting and boundary conditions.
subroutine init_nest_bc_type(nest_domain, nest_BC_buffers, npz, nest_level, position)
The module 'fv_arrays' contains the 'fv_atmos_type' and associated datatypes.
subroutine nested_grid_bc_2d(var_nest, var_coarse, ind, wt, istag, jstag, npx, npy, bd, isg, ieg, jsg, jeg, nstep_in, nsplit_in)
interface 'nested_grid_BC' includes subroutines 'nested_grid_BC_2d' and 'nested_grid_BC_3d' that fetc...
subroutine nested_grid_bc_recv_vector(nest_domain, npz, bd, nest_BC_u_buffers, nest_BC_v_buffers, nest_level, flags, gridtype)
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)
subroutine nested_grid_bc_send_scalar(var_coarse, nest_domain, istag, jstag, nest_level)
The subroutine 'nested_grid_BC_send' sends coarse-grid data to create boundary conditions.
The interface 'fill_nested_grid' includes subroutines 'fill_nested_grid_2d' and 'fill_nested_grid_3d'...
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)
subroutine timing_on(blk_name)
The subroutine 'timing_on' starts a timer.
subroutine nested_grid_bc_mpp_send_3d(var_coarse, nest_domain, istag, jstag, nest_level)
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)
subroutine nested_grid_bc_mpp_send_2d(var_coarse, nest_domain, istag, jstag, nest_level)
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)
subroutine nested_grid_bc_recv_scalar(nest_domain, istag, jstag, npz, bd, nest_BC_buffers, nest_level)
subroutine, public nested_grid_bc_apply_intt(var_nest, istag, jstag, npx, npy, npz, bd, step, split, BC, bctype)
The subroutine 'nested_grid_BC_apply_intT' performs linear interpolation or extrapolation in time for...
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 'nested_grid_BC_save_proc' saves data received by 'nested_grid_BC_recv' into the datat...
subroutine, public extrapolation_bc(q, istag, jstag, npx, npy, bd, pd_in, debug_in)
The subroutine 'extrapolation_BC' performs linear extrapolation into the halo region.
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)
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)
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)
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)