73 use fms_mod
, only : fms_init, fms_end
74 use mpp_mod
, only : fatal, mpp_debug, note, mpp_clock_sync,mpp_clock_detailed, warning
75 use mpp_mod
, only : mpp_pe, mpp_npes, mpp_node, mpp_root_pe, mpp_error, mpp_set_warn_level
76 use mpp_mod
, only : mpp_declare_pelist, mpp_set_current_pelist, mpp_sync
77 use mpp_mod
, only : mpp_clock_begin, mpp_clock_end, mpp_clock_id
78 use mpp_mod
, only : mpp_chksum, stdout, stderr, mpp_broadcast
79 use mpp_mod
, only : mpp_send, mpp_recv, mpp_sync_self, event_recv, mpp_gather
80 use mpp_domains_mod
, only : global_data_domain, bitwise_exact_sum, bgrid_ne, fold_north_edge, cgrid_ne
81 use mpp_domains_mod
, only : mpp_domain_time, cyclic_global_domain, nupdate,eupdate, xupdate, yupdate, scalar_pair
82 use mpp_domains_mod
, only : domain1d, domain2d, domaincommunicator2d, mpp_get_ntile_count
83 use mpp_domains_mod
, only : mpp_get_compute_domain, mpp_get_data_domain
84 use mpp_domains_mod
, only : mpp_global_field, mpp_global_sum, mpp_global_max, mpp_global_min
85 use mpp_domains_mod
, only : mpp_domains_init, mpp_domains_exit, mpp_broadcast_domain
86 use mpp_domains_mod
, only : mpp_check_field, mpp_define_layout
87 use mpp_domains_mod
, only : mpp_get_neighbor_pe, mpp_define_mosaic, mpp_define_io_domain
88 use mpp_domains_mod
, only : north, north_east, east, south_east
89 use mpp_domains_mod
, only : south, south_west, west, north_west
90 use mpp_domains_mod
, only : mpp_start_group_update, mpp_complete_group_update
91 use mpp_domains_mod
, only : mpp_group_update_initialized, mpp_do_group_update
92 use mpp_domains_mod
, only : mpp_create_group_update,mpp_reset_group_update_field
93 use mpp_domains_mod
, only : group_halo_update_type => mpp_group_update_type
94 use mpp_parameter_mod
, only : wupdate, eupdate, supdate, nupdate, xupdate, yupdate
96 use fms_io_mod
, only: set_domain
97 use mpp_mod
, only : mpp_get_current_pelist
98 use mpp_domains_mod
, only : mpp_define_domains
99 use mpp_domains_mod
, only : mpp_define_nest_domains, nest_domain_type
100 use mpp_domains_mod
, only : mpp_get_c2f_index, mpp_update_nest_fine
101 use mpp_domains_mod
, only : mpp_get_f2c_index, mpp_update_nest_coarse
102 use mpp_domains_mod
, only : mpp_get_domain_shift
103 use ensemble_manager_mod
, only : get_ensemble_id
108 integer,
parameter::
ng = 3
111 integer,
parameter :: xdir=1
112 integer,
parameter :: ydir=2
113 integer :: commglobal, ierror, npes
118 integer,
allocatable,
dimension(:) :: npes_tile, tile1, tile2
119 integer,
allocatable,
dimension(:) :: istart1, iend1, jstart1, jend1
120 integer,
allocatable,
dimension(:) :: istart2, iend2, jstart2, jend2
121 integer,
allocatable,
dimension(:,:) :: layout2d, global_indices
126 type(nest_domain_type),
allocatable,
dimension(:) :: nest_domain
127 integer :: this_pe_grid = 0
128 integer,
EXTERNAL :: omp_get_thread_num, omp_get_num_threads
130 integer :: npes_this_grid
135 integer :: is, ie, js, je
136 integer :: isd, ied, jsd, jed
137 integer :: isc, iec, jsc, jec
142 public mp_start, mp_assign_gid, mp_barrier, mp_stop
143 public domain_decomp, mp_bcst, mp_reduce_max, mp_reduce_sum, mp_gather
145 public fill_corners, xdir, ydir
146 public switch_current_domain, switch_current_atm, broadcast_domains
147 public is_master, setup_master
151 public is, ie, js, je, isd, ied, jsd, jed, isc, iec, jsc, jec,
ng 152 public start_group_halo_update, complete_group_halo_update
153 public group_halo_update_type
155 interface start_group_halo_update
156 module procedure start_var_group_update_2d
157 module procedure start_var_group_update_3d
158 module procedure start_var_group_update_4d
159 module procedure start_vector_group_update_2d
160 module procedure start_vector_group_update_3d
161 end interface start_group_halo_update
163 INTERFACE fill_corners
164 MODULE PROCEDURE fill_corners_2d_r4
165 MODULE PROCEDURE fill_corners_2d_r8
166 MODULE PROCEDURE fill_corners_xy_2d_r4
167 MODULE PROCEDURE fill_corners_xy_2d_r8
168 MODULE PROCEDURE fill_corners_xy_3d_r4
169 MODULE PROCEDURE fill_corners_xy_3d_r8
172 INTERFACE fill_corners_agrid
173 MODULE PROCEDURE fill_corners_agrid_r4
174 MODULE PROCEDURE fill_corners_agrid_r8
177 INTERFACE fill_corners_cgrid
178 MODULE PROCEDURE fill_corners_cgrid_r4
179 MODULE PROCEDURE fill_corners_cgrid_r8
182 INTERFACE fill_corners_dgrid
183 MODULE PROCEDURE fill_corners_dgrid_r4
184 MODULE PROCEDURE fill_corners_dgrid_r8
190 MODULE PROCEDURE mp_bcst_i
191 MODULE PROCEDURE mp_bcst_r4
192 MODULE PROCEDURE mp_bcst_r8
193 MODULE PROCEDURE mp_bcst_1d_r4
194 MODULE PROCEDURE mp_bcst_1d_r8
195 MODULE PROCEDURE mp_bcst_2d_r4
196 MODULE PROCEDURE mp_bcst_2d_r8
197 MODULE PROCEDURE mp_bcst_3d_r4
198 MODULE PROCEDURE mp_bcst_3d_r8
199 MODULE PROCEDURE mp_bcst_4d_r4
200 MODULE PROCEDURE mp_bcst_4d_r8
201 MODULE PROCEDURE mp_bcst_1d_i
202 MODULE PROCEDURE mp_bcst_2d_i
203 MODULE PROCEDURE mp_bcst_3d_i
204 MODULE PROCEDURE mp_bcst_4d_i
210 INTERFACE mp_reduce_min
211 MODULE PROCEDURE mp_reduce_min_r4
212 MODULE PROCEDURE mp_reduce_min_r8
218 INTERFACE mp_reduce_max
219 MODULE PROCEDURE mp_reduce_max_r4_1d
220 MODULE PROCEDURE mp_reduce_max_r4
221 MODULE PROCEDURE mp_reduce_max_r8_1d
222 MODULE PROCEDURE mp_reduce_max_r8
223 MODULE PROCEDURE mp_reduce_max_i
230 INTERFACE mp_reduce_sum
231 MODULE PROCEDURE mp_reduce_sum_r4
232 MODULE PROCEDURE mp_reduce_sum_r4_1d
233 MODULE PROCEDURE mp_reduce_sum_r4_1darr
234 MODULE PROCEDURE mp_reduce_sum_r4_2darr
235 MODULE PROCEDURE mp_reduce_sum_r8
236 MODULE PROCEDURE mp_reduce_sum_r8_1d
237 MODULE PROCEDURE mp_reduce_sum_r8_1darr
238 MODULE PROCEDURE mp_reduce_sum_r8_2darr
245 MODULE PROCEDURE mp_gather_4d_r4
246 MODULE PROCEDURE mp_gather_3d_r4
247 MODULE PROCEDURE mp_gather_3d_r8
250 integer :: halo_update_type = 1
254 subroutine mp_assign_gid
259 end subroutine mp_assign_gid
264 subroutine mp_start(commID, halo_update_type_in)
265 integer,
intent(in),
optional :: commid
266 integer,
intent(in),
optional :: halo_update_type_in
272 commglobal = mpi_comm_world
273 if(
PRESENT(commid) )
then 276 halo_update_type = halo_update_type_in
285 if ( mpp_pe()==mpp_root_pe() )
then 288 write(unit,*)
'Starting PEs : ', mpp_npes()
289 write(unit,*)
'Starting Threads : ', numthreads
294 if (mpp_npes() > 1)
call mpi_barrier(commglobal, ierror)
296 end subroutine mp_start
301 logical function is_master()
305 end function is_master
307 subroutine setup_master(pelist_local)
309 integer,
intent(IN) :: pelist_local(:)
311 if (any(
gid == pelist_local))
then 318 end subroutine setup_master
323 subroutine mp_barrier()
325 call mpi_barrier(commglobal, ierror)
327 end subroutine mp_barrier
337 call mpi_barrier(commglobal, ierror)
342 end subroutine mp_stop
350 subroutine domain_decomp(npx,npy,nregions,grid_type,nested,Atm,layout,io_layout)
352 integer,
intent(IN) :: npx,npy,grid_type
353 integer,
intent(INOUT) :: nregions
354 logical,
intent(IN):: nested
356 integer,
intent(INOUT) :: layout(2), io_layout(2)
358 integer,
allocatable :: pe_start(:), pe_end(:)
360 integer :: nx,ny,n,num_alloc
361 character(len=32) ::
type =
"unknown" 362 logical :: is_symmetry
363 logical :: debug=.false.
364 integer,
allocatable :: tile_id(:)
367 integer :: npes_x, npes_y
369 integer,
pointer :: pelist(:), grid_number, num_contact, npes_per_tile
370 logical,
pointer :: square_domain
371 type(domain2d),
pointer :: domain, domain_for_coupler
378 grid_number => atm%grid_number
379 num_contact => atm%num_contact
381 domain_for_coupler => atm%domain_for_coupler
382 npes_per_tile => atm%npes_per_tile
388 call mpp_domains_init(mpp_domain_time)
390 select case(nregions)
393 select case (grid_type)
396 type =
"Cubed-sphere nested grid" 398 type =
"Cubed-sphere, single face" 402 npes_per_tile = npes_x*npes_y
404 call mpp_define_layout( (/1,npx-1,1,npy-1/), npes_per_tile, layout )
406 if ( npes_x == 0 )
then 409 if ( npes_y == 0 )
then 413 if ( npes_x==npes_y .and. (npx-1)==((npx-1)/npes_x)*npes_x ) atm%gridstruct%square_domain = .true.
415 if ( (npx/npes_x <
ng) .or. (npy/npes_y <
ng) )
then 416 write(*,310) npes_x, npes_y, npx/npes_x, npy/npes_y
421 layout = (/npes_x,npes_y/)
423 type=
"Lat-Lon: cyclic" 426 if( mod(npes,nregions) .NE. 0 )
then 427 call mpp_error(note,
'TEST_MPP_DOMAINS: for Cyclic mosaic, npes should be multiple of nregions. ' // &
428 'No test is done for Cyclic mosaic. ' )
431 npes_per_tile = npes/nregions
432 call mpp_define_layout( (/1,npx-1,1,npy-1/), npes_per_tile, layout )
433 layout = (/1,npes_per_tile/)
435 type=
"Cartesian: double periodic" 438 npes_per_tile = npes/nregions
439 if(npes_x*npes_y == npes_per_tile)
then 440 layout = (/npes_x,npes_y/)
442 call mpp_define_layout( (/1,npx-1,1,npy-1/), npes_per_tile, layout )
445 type=
"Lat-Lon: patch" 448 npes_per_tile = npes/nregions
449 call mpp_define_layout( (/1,npx-1,1,npy-1/), npes_per_tile, layout )
451 type=
"Lat-Lon: strip" 454 npes_per_tile = npes/nregions
455 call mpp_define_layout( (/1,npx-1,1,npy-1/), npes_per_tile, layout )
457 type=
"Cartesian: channel" 460 npes_per_tile = npes/nregions
461 call mpp_define_layout( (/1,npx-1,1,npy-1/), npes_per_tile, layout )
465 type=
"Cubic: cubed-sphere" 467 call mpp_error(fatal,
'For a nested grid with grid_type < 3 nregions_domain must equal 1.')
472 npes_per_tile = npes_x*npes_y
473 call mpp_define_layout( (/1,npx-1,1,npy-1/), npes_per_tile, layout )
475 if ( npes_x == 0 )
then 478 if ( npes_y == 0 )
then 482 if ( npes_x==npes_y .and. (npx-1)==((npx-1)/npes_x)*npes_x ) atm%gridstruct%square_domain = .true.
484 if ( (npx/npes_x <
ng) .or. (npy/npes_y <
ng) )
then 485 write(*,310) npes_x, npes_y, npx/npes_x, npy/npes_y
486 310
format(
'Invalid layout, NPES_X:',i4.4,
'NPES_Y:',i4.4,
'ncells_X:',i4.4,
'ncells_Y:',i4.4)
491 layout = (/npes_x,npes_y/)
493 call mpp_error(fatal,
'domain_decomp: no such test: '//type)
496 allocate(layout2d(2,nregions), global_indices(4,nregions), npes_tile(nregions) )
497 allocate(pe_start(nregions),pe_end(nregions))
498 npes_tile = npes_per_tile
500 global_indices(:,n) = (/1,npx-1,1,npy-1/)
501 layout2d(:,n) = layout
502 pe_start(n) = pelist(1) + (n-1)*layout(1)*layout(2)
503 pe_end(n) = pe_start(n) + layout(1)*layout(2) -1
505 num_alloc=max(1,num_contact)
506 allocate(tile1(num_alloc), tile2(num_alloc) )
507 allocate(istart1(num_alloc), iend1(num_alloc), jstart1(num_alloc), jend1(num_alloc) )
508 allocate(istart2(num_alloc), iend2(num_alloc), jstart2(num_alloc), jend2(num_alloc) )
511 select case(nregions)
514 select case (grid_type)
519 tile1(1) = 1; tile2(1) = 2
520 istart1(1) = nx; iend1(1) = nx; jstart1(1) = 1; jend1(1) = ny
521 istart2(1) = 1; iend2(1) = 1; jstart2(1) = 1; jend2(1) = ny
523 tile1(2) = 1; tile2(2) = 3
524 istart1(2) = 1; iend1(2) = nx; jstart1(2) = 1; jend1(2) = 1
525 istart2(2) = 1; iend2(2) = nx; jstart2(2) = ny; jend2(2) = ny
527 tile1(3) = 1; tile2(3) = 2
528 istart1(3) = 1; iend1(3) = 1; jstart1(3) = 1; jend1(3) = ny
529 istart2(3) = nx; iend2(3) = nx; jstart2(3) = 1; jend2(3) = ny
531 tile1(4) = 1; tile2(4) = 3
532 istart1(4) = 1; iend1(4) = nx; jstart1(4) = ny; jend1(4) = ny
533 istart2(4) = 1; iend2(4) = nx; jstart2(4) = 1; jend2(4) = 1
535 tile1(5) = 2; tile2(5) = 4
536 istart1(5) = 1; iend1(5) = nx; jstart1(5) = 1; jend1(5) = 1
537 istart2(5) = 1; iend2(5) = nx; jstart2(5) = ny; jend2(5) = ny
539 tile1(6) = 2; tile2(6) = 4
540 istart1(6) = 1; iend1(6) = nx; jstart1(6) = ny; jend1(6) = ny
541 istart2(6) = 1; iend2(6) = nx; jstart2(6) = 1; jend2(6) = 1
543 tile1(7) = 3; tile2(7) = 4
544 istart1(7) = nx; iend1(7) = nx; jstart1(7) = 1; jend1(7) = ny
545 istart2(7) = 1; iend2(7) = 1; jstart2(7) = 1; jend2(7) = ny
547 tile1(8) = 3; tile2(8) = 4
548 istart1(8) = 1; iend1(8) = 1; jstart1(8) = 1; jend1(8) = ny
549 istart2(8) = nx; iend2(8) = nx; jstart2(8) = 1; jend2(8) = ny
550 is_symmetry = .false.
553 tile1(1) = 1; tile2(1) = 1
554 istart1(1) = nx; iend1(1) = nx; jstart1(1) = 1; jend1(1) = ny
555 istart2(1) = 1; iend2(1) = 1; jstart2(1) = 1; jend2(1) = ny
557 tile1(2) = 1; tile2(2) = 1
558 istart1(2) = 1; iend1(2) = nx; jstart1(2) = 1; jend1(2) = 1
559 istart2(2) = 1; iend2(2) = nx; jstart2(2) = ny; jend2(2) = ny
564 tile1(1) = 1; tile2(1) = 1
565 istart1(1) = nx; iend1(1) = nx; jstart1(1) = 1; jend1(1) = ny
566 istart2(1) = 1; iend2(1) = 1; jstart2(1) = 1; jend2(1) = ny
569 tile1(1) = 1; tile2(1) = 1
570 istart1(1) = nx; iend1(1) = nx; jstart1(1) = 1; jend1(1) = ny
571 istart2(1) = 1; iend2(1) = 1; jstart2(1) = 1; jend2(1) = ny
576 tile1(1) = 1; tile2(1) = 2
577 istart1(1) = nx; iend1(1) = nx; jstart1(1) = 1; jend1(1) = ny
578 istart2(1) = 1; iend2(1) = 1; jstart2(1) = 1; jend2(1) = ny
580 tile1(2) = 1; tile2(2) = 3
581 istart1(2) = 1; iend1(2) = nx; jstart1(2) = ny; jend1(2) = ny
582 istart2(2) = 1; iend2(2) = 1; jstart2(2) = ny; jend2(2) = 1
584 tile1(3) = 1; tile2(3) = 5
585 istart1(3) = 1; iend1(3) = 1; jstart1(3) = 1; jend1(3) = ny
586 istart2(3) = nx; iend2(3) = 1; jstart2(3) = ny; jend2(3) = ny
588 tile1(4) = 1; tile2(4) = 6
589 istart1(4) = 1; iend1(4) = nx; jstart1(4) = 1; jend1(4) = 1
590 istart2(4) = 1; iend2(4) = nx; jstart2(4) = ny; jend2(4) = ny
592 tile1(5) = 2; tile2(5) = 3
593 istart1(5) = 1; iend1(5) = nx; jstart1(5) = ny; jend1(5) = ny
594 istart2(5) = 1; iend2(5) = nx; jstart2(5) = 1; jend2(5) = 1
596 tile1(6) = 2; tile2(6) = 4
597 istart1(6) = nx; iend1(6) = nx; jstart1(6) = 1; jend1(6) = ny
598 istart2(6) = nx; iend2(6) = 1; jstart2(6) = 1; jend2(6) = 1
600 tile1(7) = 2; tile2(7) = 6
601 istart1(7) = 1; iend1(7) = nx; jstart1(7) = 1; jend1(7) = 1
602 istart2(7) = nx; iend2(7) = nx; jstart2(7) = ny; jend2(7) = 1
604 tile1(8) = 3; tile2(8) = 4
605 istart1(8) = nx; iend1(8) = nx; jstart1(8) = 1; jend1(8) = ny
606 istart2(8) = 1; iend2(8) = 1; jstart2(8) = 1; jend2(8) = ny
608 tile1(9) = 3; tile2(9) = 5
609 istart1(9) = 1; iend1(9) = nx; jstart1(9) = ny; jend1(9) = ny
610 istart2(9) = 1; iend2(9) = 1; jstart2(9) = ny; jend2(9) = 1
612 tile1(10) = 4; tile2(10) = 5
613 istart1(10) = 1; iend1(10) = nx; jstart1(10) = ny; jend1(10) = ny
614 istart2(10) = 1; iend2(10) = nx; jstart2(10) = 1; jend2(10) = 1
616 tile1(11) = 4; tile2(11) = 6
617 istart1(11) = nx; iend1(11) = nx; jstart1(11) = 1; jend1(11) = ny
618 istart2(11) = nx; iend2(11) = 1; jstart2(11) = 1; jend2(11) = 1
620 tile1(12) = 5; tile2(12) = 6
621 istart1(12) = nx; iend1(12) = nx; jstart1(12) = 1; jend1(12) = ny
622 istart2(12) = 1; iend2(12) = 1; jstart2(12) = 1; jend2(12) = ny
625 if ( any(pelist ==
gid) )
then 626 allocate(tile_id(nregions))
628 if( nregions .NE. 1 )
then 629 call mpp_error(fatal,
'domain_decomp: nregions should be 1 for nested region, contact developer')
638 call mpp_define_mosaic(global_indices, layout2d, domain, nregions, num_contact, tile1, tile2, &
639 istart1, iend1, jstart1, jend1, istart2, iend2, jstart2, jend2, &
640 pe_start=pe_start, pe_end=pe_end, symmetry=is_symmetry, &
641 shalo =
ng, nhalo =
ng, whalo =
ng, ehalo =
ng, tile_id=tile_id, name = type)
642 call mpp_define_mosaic(global_indices, layout2d, domain_for_coupler, nregions, num_contact, tile1, tile2, &
643 istart1, iend1, jstart1, jend1, istart2, iend2, jstart2, jend2, &
644 pe_start=pe_start, pe_end=pe_end, symmetry=is_symmetry, &
645 shalo = 1, nhalo = 1, whalo = 1, ehalo = 1, tile_id=tile_id, name = type)
647 call mpp_define_io_domain(domain, io_layout)
648 call mpp_define_io_domain(domain_for_coupler, io_layout)
652 deallocate(pe_start,pe_end)
653 deallocate(layout2d, global_indices, npes_tile)
654 deallocate(tile1, tile2)
655 deallocate(istart1, iend1, jstart1, jend1)
656 deallocate(istart2, iend2, jstart2, jend2)
659 atm%tile = (
gid-pelist(1))/npes_per_tile+1
660 if (any(pelist ==
gid))
then 661 npes_this_grid = npes_per_tile*nregions
663 call mpp_get_compute_domain( domain, is, ie, js, je )
664 call mpp_get_data_domain ( domain, isd, ied, jsd, jed )
681 if (debug .and. nregions==1)
then 683 write(*,200) tile, is, ie, js, je
687 200
format(i4.4,
' ', i4.4,
' ', i4.4,
' ', i4.4,
' ', i4.4,
' ')
707 end subroutine domain_decomp
712 subroutine start_var_group_update_2d(group, array, domain, flags, position, whalo, ehalo, shalo, nhalo, complete)
713 type(group_halo_update_type),
intent(inout) :: group
714 real,
dimension(:,:),
intent(inout) :: array
715 type(domain2d),
intent(inout) :: domain
716 integer,
optional,
intent(in) :: flags
717 integer,
optional,
intent(in) :: position
718 integer,
optional,
intent(in) :: whalo, ehalo, shalo, nhalo
719 logical,
optional,
intent(in) :: complete
722 logical :: is_complete
737 if (mpp_group_update_initialized(group))
then 738 call mpp_reset_group_update_field(group,array)
740 call mpp_create_group_update(group, array, domain, flags=flags, position=position, &
741 whalo=whalo, ehalo=ehalo, shalo=shalo, nhalo=nhalo)
745 if(
present(complete)) is_complete = complete
746 if(is_complete .and. halo_update_type == 1)
then 747 call mpp_start_group_update(group, domain, d_type)
750 end subroutine start_var_group_update_2d
753 subroutine start_var_group_update_3d(group, array, domain, flags, position, whalo, ehalo, shalo, nhalo, complete)
754 type(group_halo_update_type),
intent(inout) :: group
755 real,
dimension(:,:,:),
intent(inout) :: array
756 type(domain2d),
intent(inout) :: domain
757 integer,
optional,
intent(in) :: flags
758 integer,
optional,
intent(in) :: position
759 integer,
optional,
intent(in) :: whalo, ehalo, shalo, nhalo
760 logical,
optional,
intent(in) :: complete
763 logical :: is_complete
779 if (mpp_group_update_initialized(group))
then 780 call mpp_reset_group_update_field(group,array)
782 call mpp_create_group_update(group, array, domain, flags=flags, position=position, &
783 whalo=whalo, ehalo=ehalo, shalo=shalo, nhalo=nhalo)
787 if(
present(complete)) is_complete = complete
788 if(is_complete .and. halo_update_type == 1 )
then 789 call mpp_start_group_update(group, domain, d_type)
792 end subroutine start_var_group_update_3d
794 subroutine start_var_group_update_4d(group, array, domain, flags, position, whalo, ehalo, shalo, nhalo, complete)
795 type(group_halo_update_type),
intent(inout) :: group
796 real,
dimension(:,:,:,:),
intent(inout) :: array
797 type(domain2d),
intent(inout) :: domain
798 integer,
optional,
intent(in) :: flags
799 integer,
optional,
intent(in) :: position
801 integer,
optional,
intent(in) :: whalo, ehalo, shalo, nhalo
802 logical,
optional,
intent(in) :: complete
805 logical :: is_complete
823 if (mpp_group_update_initialized(group))
then 824 call mpp_reset_group_update_field(group,array)
826 call mpp_create_group_update(group, array, domain, flags=flags, position=position, &
827 whalo=whalo, ehalo=ehalo, shalo=shalo, nhalo=nhalo)
831 if(
present(complete)) is_complete = complete
832 if(is_complete .and. halo_update_type == 1 )
then 833 call mpp_start_group_update(group, domain, d_type)
836 end subroutine start_var_group_update_4d
840 subroutine start_vector_group_update_2d(group, u_cmpt, v_cmpt, domain, flags, gridtype, whalo, ehalo, shalo, nhalo, complete)
841 type(group_halo_update_type),
intent(inout) :: group
842 real,
dimension(:,:),
intent(inout) :: u_cmpt, v_cmpt
845 type(domain2d),
intent(inout) :: domain
846 integer,
optional,
intent(in) :: flags
847 integer,
optional,
intent(in) :: gridtype
850 integer,
optional,
intent(in) :: whalo, ehalo, shalo, nhalo
851 logical,
optional,
intent(in) :: complete
854 logical :: is_complete
874 if (mpp_group_update_initialized(group))
then 875 call mpp_reset_group_update_field(group,u_cmpt, v_cmpt)
877 call mpp_create_group_update(group, u_cmpt, v_cmpt, domain, &
878 flags=flags, gridtype=gridtype, &
879 whalo=whalo, ehalo=ehalo, shalo=shalo, nhalo=nhalo)
883 if(
present(complete)) is_complete = complete
884 if(is_complete .and. halo_update_type == 1 )
then 885 call mpp_start_group_update(group, domain, d_type)
888 end subroutine start_vector_group_update_2d
890 subroutine start_vector_group_update_3d(group, u_cmpt, v_cmpt, domain, flags, gridtype, whalo, ehalo, shalo, nhalo, complete)
891 type(group_halo_update_type),
intent(inout) :: group
892 real,
dimension(:,:,:),
intent(inout) :: u_cmpt, v_cmpt
895 type(domain2d),
intent(inout) :: domain
896 integer,
optional,
intent(in) :: flags
897 integer,
optional,
intent(in) :: gridtype
900 integer,
optional,
intent(in) :: whalo, ehalo, shalo, nhalo
901 logical,
optional,
intent(in) :: complete
904 logical :: is_complete
924 if (mpp_group_update_initialized(group))
then 925 call mpp_reset_group_update_field(group,u_cmpt, v_cmpt)
927 call mpp_create_group_update(group, u_cmpt, v_cmpt, domain, &
928 flags=flags, gridtype=gridtype, &
929 whalo=whalo, ehalo=ehalo, shalo=shalo, nhalo=nhalo)
933 if(
present(complete)) is_complete = complete
934 if(is_complete .and. halo_update_type == 1)
then 935 call mpp_start_group_update(group, domain, d_type)
938 end subroutine start_vector_group_update_3d
941 subroutine complete_group_halo_update(group, domain)
942 type(group_halo_update_type),
intent(inout) :: group
943 type(domain2d),
intent(inout) :: domain
950 if( halo_update_type == 1 )
then 951 call mpp_complete_group_update(group, domain, d_type)
953 call mpp_do_group_update(group, domain, d_type)
956 end subroutine complete_group_halo_update
961 subroutine broadcast_domains(Atm)
965 integer :: n, i1, i2, j1, j2, i
966 integer :: ens_root_pe, ensemble_id
973 ensemble_id = get_ensemble_id()
974 ens_root_pe = (ensemble_id-1)*npes
977 call mpp_set_current_pelist((/ (i,i=ens_root_pe,npes-1+ens_root_pe) /))
979 call mpp_broadcast_domain(atm(n)%domain)
980 call mpp_broadcast_domain(atm(n)%domain_for_coupler)
983 end subroutine broadcast_domains
985 subroutine switch_current_domain(new_domain,new_domain_for_coupler)
987 type(domain2d),
intent(in),
target :: new_domain, new_domain_for_coupler
988 logical,
parameter :: debug = .false.
993 call mpp_get_compute_domain( new_domain, is, ie, js, je )
996 call mpp_get_data_domain ( new_domain, isd, ied, jsd, jed )
1002 call set_domain(new_domain)
1005 end subroutine switch_current_domain
1007 subroutine switch_current_atm(new_Atm, switch_domain)
1010 logical,
intent(IN),
optional :: switch_domain
1011 logical,
parameter :: debug = .false.
1014 if (debug .AND. (
gid==
masterproc)) print*,
'SWITCHING ATM STRUCTURES', new_atm%grid_number
1015 if (
present(switch_domain))
then 1020 if (swd)
call switch_current_domain(new_atm%domain, new_atm%domain_for_coupler)
1025 end subroutine switch_current_atm
1030 subroutine fill_corners_2d_r4(q, npx, npy, FILL, AGRID, BGRID)
1031 real(kind=4),
DIMENSION(isd:,jsd:),
intent(INOUT):: q
1032 integer,
intent(IN):: npx,npy
1033 integer,
intent(IN):: fill
1034 logical,
OPTIONAL,
intent(IN) :: agrid, bgrid
1037 if (
present(bgrid))
then 1043 if ((is== 1) .and. (js== 1)) q(1-i ,1-j ) = q(1-j ,i+1 )
1044 if ((is== 1) .and. (je==npy-1)) q(1-i ,npy+j) = q(1-j ,npy-i )
1045 if ((ie==npx-1) .and. (js== 1)) q(npx+i,1-j ) = q(npx+j,i+1 )
1046 if ((ie==npx-1) .and. (je==npy-1)) q(npx+i,npy+j) = q(npx+j,npy-i )
1052 if ((is== 1) .and. (js== 1)) q(1-j ,1-i ) = q(i+1 ,1-j )
1053 if ((is== 1) .and. (je==npy-1)) q(1-j ,npy+i) = q(i+1 ,npy+j )
1054 if ((ie==npx-1) .and. (js== 1)) q(npx+j,1-i ) = q(npx-i,1-j )
1055 if ((ie==npx-1) .and. (je==npy-1)) q(npx+j,npy+i) = q(npx-i,npy+j )
1061 if ((is== 1) .and. (js== 1)) q(1-i ,1-j ) = q(1-j ,i+1 )
1062 if ((is== 1) .and. (je==npy-1)) q(1-i ,npy+j) = q(1-j ,npy-i )
1063 if ((ie==npx-1) .and. (js== 1)) q(npx+i,1-j ) = q(npx+j,i+1 )
1064 if ((ie==npx-1) .and. (je==npy-1)) q(npx+i,npy+j) = q(npx+j,npy-i )
1069 elseif (
present(agrid))
then 1075 if ((is== 1) .and. (js== 1)) q(1-i ,1-j ) = q(1-j ,i )
1076 if ((is== 1) .and. (je==npy-1)) q(1-i ,npy-1+j) = q(1-j ,npy-1-i+1)
1077 if ((ie==npx-1) .and. (js== 1)) q(npx-1+i,1-j ) = q(npx-1+j,i )
1078 if ((ie==npx-1) .and. (je==npy-1)) q(npx-1+i,npy-1+j) = q(npx-1+j,npy-1-i+1)
1084 if ((is== 1) .and. (js== 1)) q(1-j ,1-i ) = q(i ,1-j )
1085 if ((is== 1) .and. (je==npy-1)) q(1-j ,npy-1+i) = q(i ,npy-1+j)
1086 if ((ie==npx-1) .and. (js== 1)) q(npx-1+j,1-i ) = q(npx-1-i+1,1-j )
1087 if ((ie==npx-1) .and. (je==npy-1)) q(npx-1+j,npy-1+i) = q(npx-1-i+1,npy-1+j)
1093 if ((is== 1) .and. (js== 1)) q(1-j ,1-i ) = q(i ,1-j )
1094 if ((is== 1) .and. (je==npy-1)) q(1-j ,npy-1+i) = q(i ,npy-1+j)
1095 if ((ie==npx-1) .and. (js== 1)) q(npx-1+j,1-i ) = q(npx-1-i+1,1-j )
1096 if ((ie==npx-1) .and. (je==npy-1)) q(npx-1+j,npy-1+i) = q(npx-1-i+1,npy-1+j)
1103 end subroutine fill_corners_2d_r4
1110 subroutine fill_corners_2d_r8(q, npx, npy, FILL, AGRID, BGRID)
1111 real(kind=8),
DIMENSION(isd:,jsd:),
intent(INOUT):: q
1112 integer,
intent(IN):: npx,npy
1113 integer,
intent(IN):: fill
1114 logical,
OPTIONAL,
intent(IN) :: agrid, bgrid
1117 if (
present(bgrid))
then 1123 if ((is== 1) .and. (js== 1)) q(1-i ,1-j ) = q(1-j ,i+1 )
1124 if ((is== 1) .and. (je==npy-1)) q(1-i ,npy+j) = q(1-j ,npy-i )
1125 if ((ie==npx-1) .and. (js== 1)) q(npx+i,1-j ) = q(npx+j,i+1 )
1126 if ((ie==npx-1) .and. (je==npy-1)) q(npx+i,npy+j) = q(npx+j,npy-i )
1132 if ((is== 1) .and. (js== 1)) q(1-j ,1-i ) = q(i+1 ,1-j )
1133 if ((is== 1) .and. (je==npy-1)) q(1-j ,npy+i) = q(i+1 ,npy+j )
1134 if ((ie==npx-1) .and. (js== 1)) q(npx+j,1-i ) = q(npx-i,1-j )
1135 if ((ie==npx-1) .and. (je==npy-1)) q(npx+j,npy+i) = q(npx-i,npy+j )
1141 if ((is== 1) .and. (js== 1)) q(1-i ,1-j ) = q(1-j ,i+1 )
1142 if ((is== 1) .and. (je==npy-1)) q(1-i ,npy+j) = q(1-j ,npy-i )
1143 if ((ie==npx-1) .and. (js== 1)) q(npx+i,1-j ) = q(npx+j,i+1 )
1144 if ((ie==npx-1) .and. (je==npy-1)) q(npx+i,npy+j) = q(npx+j,npy-i )
1149 elseif (
present(agrid))
then 1155 if ((is== 1) .and. (js== 1)) q(1-i ,1-j ) = q(1-j ,i )
1156 if ((is== 1) .and. (je==npy-1)) q(1-i ,npy-1+j) = q(1-j ,npy-1-i+1)
1157 if ((ie==npx-1) .and. (js== 1)) q(npx-1+i,1-j ) = q(npx-1+j,i )
1158 if ((ie==npx-1) .and. (je==npy-1)) q(npx-1+i,npy-1+j) = q(npx-1+j,npy-1-i+1)
1164 if ((is== 1) .and. (js== 1)) q(1-j ,1-i ) = q(i ,1-j )
1165 if ((is== 1) .and. (je==npy-1)) q(1-j ,npy-1+i) = q(i ,npy-1+j)
1166 if ((ie==npx-1) .and. (js== 1)) q(npx-1+j,1-i ) = q(npx-1-i+1,1-j )
1167 if ((ie==npx-1) .and. (je==npy-1)) q(npx-1+j,npy-1+i) = q(npx-1-i+1,npy-1+j)
1173 if ((is== 1) .and. (js== 1)) q(1-j ,1-i ) = q(i ,1-j )
1174 if ((is== 1) .and. (je==npy-1)) q(1-j ,npy-1+i) = q(i ,npy-1+j)
1175 if ((ie==npx-1) .and. (js== 1)) q(npx-1+j,1-i ) = q(npx-1-i+1,1-j )
1176 if ((ie==npx-1) .and. (je==npy-1)) q(npx-1+j,npy-1+i) = q(npx-1-i+1,npy-1+j)
1183 end subroutine fill_corners_2d_r8
1191 subroutine fill_corners_xy_2d_r8(x, y, npx, npy, DGRID, AGRID, CGRID, VECTOR)
1192 real(kind=8),
DIMENSION(isd:,jsd:),
intent(INOUT):: x
1193 real(kind=8),
DIMENSION(isd:,jsd:),
intent(INOUT):: y
1194 integer,
intent(IN):: npx,npy
1195 logical,
OPTIONAL,
intent(IN) :: dgrid, agrid, cgrid, vector
1198 real(kind=8) :: mysign
1201 if (
present(vector))
then 1202 if (vector) mysign = -1.0
1205 if (
present(dgrid))
then 1206 call fill_corners_dgrid(x, y, npx, npy, mysign)
1207 elseif (
present(cgrid))
then 1208 call fill_corners_cgrid(x, y, npx, npy, mysign)
1209 elseif (
present(agrid))
then 1210 call fill_corners_agrid(x, y, npx, npy, mysign)
1212 call fill_corners_agrid(x, y, npx, npy, mysign)
1215 end subroutine fill_corners_xy_2d_r8
1223 subroutine fill_corners_xy_2d_r4(x, y, npx, npy, DGRID, AGRID, CGRID, VECTOR)
1224 real(kind=4),
DIMENSION(isd:,jsd:),
intent(INOUT):: x
1225 real(kind=4),
DIMENSION(isd:,jsd:),
intent(INOUT):: y
1226 integer,
intent(IN):: npx,npy
1227 logical,
OPTIONAL,
intent(IN) :: dgrid, agrid, cgrid, vector
1230 real(kind=4) :: mysign
1233 if (
present(vector))
then 1234 if (vector) mysign = -1.0
1237 if (
present(dgrid))
then 1238 call fill_corners_dgrid(x, y, npx, npy, mysign)
1239 elseif (
present(cgrid))
then 1240 call fill_corners_cgrid(x, y, npx, npy, mysign)
1241 elseif (
present(agrid))
then 1242 call fill_corners_agrid(x, y, npx, npy, mysign)
1244 call fill_corners_agrid(x, y, npx, npy, mysign)
1247 end subroutine fill_corners_xy_2d_r4
1255 subroutine fill_corners_xy_3d_r8(x, y, npx, npy, npz, DGRID, AGRID, CGRID, VECTOR)
1256 real(kind=8),
DIMENSION(isd:,jsd:,:),
intent(INOUT):: x
1257 real(kind=8),
DIMENSION(isd:,jsd:,:),
intent(INOUT):: y
1258 integer,
intent(IN):: npx,npy,npz
1259 logical,
OPTIONAL,
intent(IN) :: dgrid, agrid, cgrid, vector
1262 real(kind=8) :: mysign
1265 if (
present(vector))
then 1266 if (vector) mysign = -1.0
1269 if (
present(dgrid))
then 1271 call fill_corners_dgrid(x(:,:,k), y(:,:,k), npx, npy, mysign)
1273 elseif (
present(cgrid))
then 1275 call fill_corners_cgrid(x(:,:,k), y(:,:,k), npx, npy, mysign)
1277 elseif (
present(agrid))
then 1279 call fill_corners_agrid(x(:,:,k), y(:,:,k), npx, npy, mysign)
1283 call fill_corners_agrid(x(:,:,k), y(:,:,k), npx, npy, mysign)
1287 end subroutine fill_corners_xy_3d_r8
1295 subroutine fill_corners_xy_3d_r4(x, y, npx, npy, npz, DGRID, AGRID, CGRID, VECTOR)
1296 real(kind=4),
DIMENSION(isd:,jsd:,:),
intent(INOUT):: x
1297 real(kind=4),
DIMENSION(isd:,jsd:,:),
intent(INOUT):: y
1298 integer,
intent(IN):: npx,npy,npz
1299 logical,
OPTIONAL,
intent(IN) :: dgrid, agrid, cgrid, vector
1302 real(kind=4) :: mysign
1305 if (
present(vector))
then 1306 if (vector) mysign = -1.0
1309 if (
present(dgrid))
then 1311 call fill_corners_dgrid(x(:,:,k), y(:,:,k), npx, npy, mysign)
1313 elseif (
present(cgrid))
then 1315 call fill_corners_cgrid(x(:,:,k), y(:,:,k), npx, npy, mysign)
1317 elseif (
present(agrid))
then 1319 call fill_corners_agrid(x(:,:,k), y(:,:,k), npx, npy, mysign)
1323 call fill_corners_agrid(x(:,:,k), y(:,:,k), npx, npy, mysign)
1327 end subroutine fill_corners_xy_3d_r4
1335 subroutine fill_corners_dgrid_r8(x, y, npx, npy, mySign)
1336 real(kind=8),
DIMENSION(isd:,jsd:),
intent(INOUT):: x
1337 real(kind=8),
DIMENSION(isd:,jsd:),
intent(INOUT):: y
1338 integer,
intent(IN):: npx,npy
1339 real(kind=8),
intent(IN) :: mysign
1348 if ((is == 1) .and. (js == 1)) x(1-i ,1-j ) = mysign*y(1-j ,i )
1349 if ((is == 1) .and. (je+1==npy)) x(1-i ,npy+j) = y(1-j ,npy-i)
1350 if ((ie+1==npx) .and. (js == 1)) x(npx-1+i,1-j ) = y(npx+j,i )
1351 if ((ie+1==npx) .and. (je+1==npy)) x(npx-1+i,npy+j) = mysign*y(npx+j,npy-i)
1360 if ((is == 1) .and. (js == 1)) y(1-i ,1-j ) = mysign*x(j ,1-i )
1361 if ((is == 1) .and. (je+1==npy)) y(1-i ,npy-1+j) = x(j ,npy+i)
1362 if ((ie+1==npx) .and. (js == 1)) y(npx+i ,1-j ) = x(npx-j ,1-i )
1363 if ((ie+1==npx) .and. (je+1==npy)) y(npx+i ,npy-1+j) = mysign*x(npx-j ,npy+i)
1367 end subroutine fill_corners_dgrid_r8
1375 subroutine fill_corners_dgrid_r4(x, y, npx, npy, mySign)
1376 real(kind=4),
DIMENSION(isd:,jsd:),
intent(INOUT):: x
1377 real(kind=4),
DIMENSION(isd:,jsd:),
intent(INOUT):: y
1378 integer,
intent(IN):: npx,npy
1379 real(kind=4),
intent(IN) :: mysign
1388 if ((is == 1) .and. (js == 1)) x(1-i ,1-j ) = mysign*y(1-j ,i )
1389 if ((is == 1) .and. (je+1==npy)) x(1-i ,npy+j) = y(1-j ,npy-i)
1390 if ((ie+1==npx) .and. (js == 1)) x(npx-1+i,1-j ) = y(npx+j,i )
1391 if ((ie+1==npx) .and. (je+1==npy)) x(npx-1+i,npy+j) = mysign*y(npx+j,npy-i)
1400 if ((is == 1) .and. (js == 1)) y(1-i ,1-j ) = mysign*x(j ,1-i )
1401 if ((is == 1) .and. (je+1==npy)) y(1-i ,npy-1+j) = x(j ,npy+i)
1402 if ((ie+1==npx) .and. (js == 1)) y(npx+i ,1-j ) = x(npx-j ,1-i )
1403 if ((ie+1==npx) .and. (je+1==npy)) y(npx+i ,npy-1+j) = mysign*x(npx-j ,npy+i)
1407 end subroutine fill_corners_dgrid_r4
1415 subroutine fill_corners_cgrid_r4(x, y, npx, npy, mySign)
1416 real(kind=4),
DIMENSION(isd:,jsd:),
intent(INOUT):: x
1417 real(kind=4),
DIMENSION(isd:,jsd:),
intent(INOUT):: y
1418 integer,
intent(IN):: npx,npy
1419 real(kind=4),
intent(IN) :: mysign
1424 if ((is == 1) .and. (js == 1)) x(1-i ,1-j ) = y(j ,1-i )
1425 if ((is == 1) .and. (je+1==npy)) x(1-i ,npy-1+j) = mysign*y(j ,npy+i)
1426 if ((ie+1==npx) .and. (js == 1)) x(npx+i ,1-j ) = mysign*y(npx-j ,1-i )
1427 if ((ie+1==npx) .and. (je+1==npy)) x(npx+i ,npy-1+j) = y(npx-j ,npy+i)
1432 if ((is == 1) .and. (js == 1)) y(1-i ,1-j ) = x(1-j ,i )
1433 if ((is == 1) .and. (je+1==npy)) y(1-i ,npy+j) = mysign*x(1-j ,npy-i)
1434 if ((ie+1==npx) .and. (js == 1)) y(npx-1+i,1-j ) = mysign*x(npx+j,i )
1435 if ((ie+1==npx) .and. (je+1==npy)) y(npx-1+i,npy+j) = x(npx+j,npy-i)
1439 end subroutine fill_corners_cgrid_r4
1447 subroutine fill_corners_cgrid_r8(x, y, npx, npy, mySign)
1448 real(kind=8),
DIMENSION(isd:,jsd:),
intent(INOUT):: x
1449 real(kind=8),
DIMENSION(isd:,jsd:),
intent(INOUT):: y
1450 integer,
intent(IN):: npx,npy
1451 real(kind=8),
intent(IN) :: mysign
1456 if ((is == 1) .and. (js == 1)) x(1-i ,1-j ) = y(j ,1-i )
1457 if ((is == 1) .and. (je+1==npy)) x(1-i ,npy-1+j) = mysign*y(j ,npy+i)
1458 if ((ie+1==npx) .and. (js == 1)) x(npx+i ,1-j ) = mysign*y(npx-j ,1-i )
1459 if ((ie+1==npx) .and. (je+1==npy)) x(npx+i ,npy-1+j) = y(npx-j ,npy+i)
1464 if ((is == 1) .and. (js == 1)) y(1-i ,1-j ) = x(1-j ,i )
1465 if ((is == 1) .and. (je+1==npy)) y(1-i ,npy+j) = mysign*x(1-j ,npy-i)
1466 if ((ie+1==npx) .and. (js == 1)) y(npx-1+i,1-j ) = mysign*x(npx+j,i )
1467 if ((ie+1==npx) .and. (je+1==npy)) y(npx-1+i,npy+j) = x(npx+j,npy-i)
1471 end subroutine fill_corners_cgrid_r8
1479 subroutine fill_corners_agrid_r4(x, y, npx, npy, mySign)
1480 real(kind=4),
DIMENSION(isd:,jsd:),
intent(INOUT):: x
1481 real(kind=4),
DIMENSION(isd:,jsd:),
intent(INOUT):: y
1482 integer,
intent(IN):: npx,npy
1483 real(kind=4),
intent(IN) :: mysign
1488 if ((is== 1) .and. (js== 1)) x(1-i ,1-j ) = mysign*y(1-j ,i )
1489 if ((is== 1) .and. (je==npy-1)) x(1-i ,npy-1+j) = y(1-j ,npy-1-i+1)
1490 if ((ie==npx-1) .and. (js== 1)) x(npx-1+i,1-j ) = y(npx-1+j,i )
1491 if ((ie==npx-1) .and. (je==npy-1)) x(npx-1+i,npy-1+j) = mysign*y(npx-1+j,npy-1-i+1)
1496 if ((is== 1) .and. (js== 1)) y(1-j ,1-i ) = mysign*x(i ,1-j )
1497 if ((is== 1) .and. (je==npy-1)) y(1-j ,npy-1+i) = x(i ,npy-1+j)
1498 if ((ie==npx-1) .and. (js== 1)) y(npx-1+j,1-i ) = x(npx-1-i+1,1-j )
1499 if ((ie==npx-1) .and. (je==npy-1)) y(npx-1+j,npy-1+i) = mysign*x(npx-1-i+1,npy-1+j)
1503 end subroutine fill_corners_agrid_r4
1511 subroutine fill_corners_agrid_r8(x, y, npx, npy, mySign)
1512 real(kind=8),
DIMENSION(isd:,jsd:),
intent(INOUT):: x
1513 real(kind=8),
DIMENSION(isd:,jsd:),
intent(INOUT):: y
1514 integer,
intent(IN):: npx,npy
1515 real(kind=8),
intent(IN) :: mysign
1520 if ((is== 1) .and. (js== 1)) x(1-i ,1-j ) = mysign*y(1-j ,i )
1521 if ((is== 1) .and. (je==npy-1)) x(1-i ,npy-1+j) = y(1-j ,npy-1-i+1)
1522 if ((ie==npx-1) .and. (js== 1)) x(npx-1+i,1-j ) = y(npx-1+j,i )
1523 if ((ie==npx-1) .and. (je==npy-1)) x(npx-1+i,npy-1+j) = mysign*y(npx-1+j,npy-1-i+1)
1528 if ((is== 1) .and. (js== 1)) y(1-j ,1-i ) = mysign*x(i ,1-j )
1529 if ((is== 1) .and. (je==npy-1)) y(1-j ,npy-1+i) = x(i ,npy-1+j)
1530 if ((ie==npx-1) .and. (js== 1)) y(npx-1+j,1-i ) = x(npx-1-i+1,1-j )
1531 if ((ie==npx-1) .and. (je==npy-1)) y(npx-1+j,npy-1+i) = mysign*x(npx-1-i+1,npy-1+j)
1535 end subroutine fill_corners_agrid_r8
1920 subroutine mp_gather_4d_r4(q, i1,i2, j1,j2, idim, jdim, kdim, ldim)
1921 integer,
intent(IN) :: i1,i2, j1,j2
1922 integer,
intent(IN) :: idim, jdim, kdim, ldim
1923 real(kind=4),
intent(INOUT):: q(idim,jdim,kdim,ldim)
1924 integer :: i,j,k,l,n,icnt
1925 integer :: lsize, lsize_buf(1)
1927 integer :: lsizes(npes_this_grid), ldispl(npes_this_grid), cnts(npes_this_grid)
1928 integer :: ldims(5), gdims(5*npes_this_grid)
1929 real(kind=4),
allocatable,
dimension(:) :: larr, garr
1936 do l=1,npes_this_grid
1940 call mpp_gather(ldims, gdims)
1943 lsize = ( (i2 - i1 + 1) * (j2 - j1 + 1) ) * kdim
1944 do l=1,npes_this_grid
1949 lsize_buf(1) = lsize
1950 call mpp_gather(lsize_buf, lsizes)
1953 allocate ( larr(lsize) )
1958 larr(icnt) = q(i,j,k,tile)
1965 call mpp_broadcast(lsizes, npes_this_grid,
masterproc)
1967 do l=2,npes_this_grid
1969 ldispl(l) = ldispl(l-1) + lsizes(l-1)
1970 gsize = gsize + lsizes(l)
1972 allocate ( garr(gsize) )
1974 call mpp_gather(larr, lsize, garr, lsizes)
1978 do n=2,npes_this_grid
1980 do l=gdims( (n-1)*5 + 5 ), gdims( (n-1)*5 + 5 )
1982 do j=gdims( (n-1)*5 + 3 ), gdims( (n-1)*5 + 4 )
1983 do i=gdims( (n-1)*5 + 1 ), gdims( (n-1)*5 + 2 )
1984 q(i,j,k,l) = garr(ldispl(n)+icnt)
1995 end subroutine mp_gather_4d_r4
2006 subroutine mp_gather_3d_r4(q, i1,i2, j1,j2, idim, jdim, ldim)
2007 integer,
intent(IN) :: i1,i2, j1,j2
2008 integer,
intent(IN) :: idim, jdim, ldim
2009 real(kind=4),
intent(INOUT):: q(idim,jdim,ldim)
2010 integer :: i,j,l,n,icnt
2011 integer :: lsize, lsize_buf(1)
2013 integer :: lsizes(npes_this_grid), ldispl(npes_this_grid), cnts(npes_this_grid)
2014 integer :: ldims(5), gdims(5*npes_this_grid)
2015 real(kind=4),
allocatable,
dimension(:) :: larr, garr
2022 do l=1,npes_this_grid
2027 call mpp_gather(ldims, gdims)
2029 lsize = ( (i2 - i1 + 1) * (j2 - j1 + 1) )
2030 do l=1,npes_this_grid
2035 lsize_buf(1) = lsize
2036 call mpp_gather(lsize_buf, lsizes)
2039 allocate ( larr(lsize) )
2043 larr(icnt) = q(i,j,tile)
2049 call mpp_broadcast(lsizes, npes_this_grid,
masterproc)
2051 do l=2,npes_this_grid
2053 ldispl(l) = ldispl(l-1) + lsizes(l-1)
2054 gsize = gsize + lsizes(l)
2056 allocate ( garr(gsize) )
2057 call mpp_gather(larr, lsize, garr, lsizes)
2060 do n=2,npes_this_grid
2062 do l=gdims( (n-1)*5 + 5 ), gdims( (n-1)*5 + 5 )
2063 do j=gdims( (n-1)*5 + 3 ), gdims( (n-1)*5 + 4 )
2064 do i=gdims( (n-1)*5 + 1 ), gdims( (n-1)*5 + 2 )
2065 q(i,j,l) = garr(ldispl(n)+icnt)
2075 end subroutine mp_gather_3d_r4
2085 subroutine mp_gather_3d_r8(q, i1,i2, j1,j2, idim, jdim, ldim)
2086 integer,
intent(IN) :: i1,i2, j1,j2
2087 integer,
intent(IN) :: idim, jdim, ldim
2088 real(kind=8),
intent(INOUT):: q(idim,jdim,ldim)
2089 integer :: i,j,l,n,icnt
2090 integer :: lsize, lsize_buf(1)
2092 integer :: lsizes(npes_this_grid), ldispl(npes_this_grid), cnts(npes_this_grid)
2093 integer :: ldims(5), gdims(5*npes_this_grid)
2094 real(kind=8),
allocatable,
dimension(:) :: larr, garr
2101 do l=1,npes_this_grid
2106 call mpp_gather(ldims, gdims)
2107 lsize = ( (i2 - i1 + 1) * (j2 - j1 + 1) )
2108 do l=1,npes_this_grid
2115 lsize_buf(1) = lsize
2116 call mpp_gather(lsize_buf, lsizes)
2118 allocate ( larr(lsize) )
2122 larr(icnt) = q(i,j,tile)
2127 call mpp_broadcast(lsizes, npes_this_grid,
masterproc)
2130 do l=2,npes_this_grid
2132 ldispl(l) = ldispl(l-1) + lsizes(l-1)
2133 gsize = gsize + lsizes(l)
2136 allocate ( garr(gsize) )
2137 call mpp_gather(larr, lsize, garr, lsizes)
2140 do n=2,npes_this_grid
2142 do l=gdims( (n-1)*5 + 5 ), gdims( (n-1)*5 + 5 )
2143 do j=gdims( (n-1)*5 + 3 ), gdims( (n-1)*5 + 4 )
2144 do i=gdims( (n-1)*5 + 1 ), gdims( (n-1)*5 + 2 )
2145 q(i,j,l) = garr(ldispl(n)+icnt)
2155 end subroutine mp_gather_3d_r8
2165 subroutine mp_bcst_i(q)
2166 integer,
intent(INOUT) :: q
2168 call mpi_bcast(q, 1, mpi_integer,
masterproc, commglobal, ierror)
2170 end subroutine mp_bcst_i
2180 subroutine mp_bcst_r4(q)
2181 real(kind=4),
intent(INOUT) :: q
2183 call mpi_bcast(q, 1, mpi_real,
masterproc, commglobal, ierror)
2185 end subroutine mp_bcst_r4
2195 subroutine mp_bcst_r8(q)
2196 real(kind=8),
intent(INOUT) :: q
2198 call mpi_bcast(q, 1, mpi_double_precision,
masterproc, commglobal, ierror)
2200 end subroutine mp_bcst_r8
2210 subroutine mp_bcst_1d_r4(q, idim)
2211 integer,
intent(IN) :: idim
2212 real(kind=4),
intent(INOUT) :: q(idim)
2214 call mpi_bcast(q, idim, mpi_real,
masterproc, commglobal, ierror)
2216 end subroutine mp_bcst_1d_r4
2226 subroutine mp_bcst_1d_r8(q, idim)
2227 integer,
intent(IN) :: idim
2228 real(kind=8),
intent(INOUT) :: q(idim)
2230 call mpi_bcast(q, idim, mpi_double_precision,
masterproc, commglobal, ierror)
2232 end subroutine mp_bcst_1d_r8
2242 subroutine mp_bcst_2d_r4(q, idim, jdim)
2243 integer,
intent(IN) :: idim, jdim
2244 real(kind=4),
intent(INOUT) :: q(idim,jdim)
2246 call mpi_bcast(q, idim*jdim, mpi_real,
masterproc, commglobal, ierror)
2248 end subroutine mp_bcst_2d_r4
2258 subroutine mp_bcst_2d_r8(q, idim, jdim)
2259 integer,
intent(IN) :: idim, jdim
2260 real(kind=8),
intent(INOUT) :: q(idim,jdim)
2262 call mpi_bcast(q, idim*jdim, mpi_double_precision,
masterproc, commglobal, ierror)
2264 end subroutine mp_bcst_2d_r8
2274 subroutine mp_bcst_3d_r4(q, idim, jdim, kdim)
2275 integer,
intent(IN) :: idim, jdim, kdim
2276 real(kind=4),
intent(INOUT) :: q(idim,jdim,kdim)
2278 call mpi_bcast(q, idim*jdim*kdim, mpi_real,
masterproc, commglobal, ierror)
2280 end subroutine mp_bcst_3d_r4
2290 subroutine mp_bcst_3d_r8(q, idim, jdim, kdim)
2291 integer,
intent(IN) :: idim, jdim, kdim
2292 real(kind=8),
intent(INOUT) :: q(idim,jdim,kdim)
2294 call mpi_bcast(q, idim*jdim*kdim, mpi_double_precision,
masterproc, commglobal, ierror)
2296 end subroutine mp_bcst_3d_r8
2306 subroutine mp_bcst_4d_r4(q, idim, jdim, kdim, ldim)
2307 integer,
intent(IN) :: idim, jdim, kdim, ldim
2308 real(kind=4),
intent(INOUT) :: q(idim,jdim,kdim,ldim)
2310 call mpi_bcast(q, idim*jdim*kdim*ldim, mpi_real,
masterproc, commglobal, ierror)
2312 end subroutine mp_bcst_4d_r4
2322 subroutine mp_bcst_4d_r8(q, idim, jdim, kdim, ldim)
2323 integer,
intent(IN) :: idim, jdim, kdim, ldim
2324 real(kind=8),
intent(INOUT) :: q(idim,jdim,kdim,ldim)
2326 call mpi_bcast(q, idim*jdim*kdim*ldim, mpi_double_precision,
masterproc, commglobal, ierror)
2328 end subroutine mp_bcst_4d_r8
2338 subroutine mp_bcst_3d_i(q, idim, jdim, kdim)
2339 integer,
intent(IN) :: idim, jdim, kdim
2340 integer,
intent(INOUT) :: q(idim,jdim,kdim)
2342 call mpi_bcast(q, idim*jdim*kdim, mpi_integer,
masterproc, commglobal, ierror)
2344 end subroutine mp_bcst_3d_i
2354 subroutine mp_bcst_1d_i(q, idim)
2355 integer,
intent(IN) :: idim
2356 integer,
intent(INOUT) :: q(idim)
2358 call mpi_bcast(q, idim, mpi_integer,
masterproc, commglobal, ierror)
2360 end subroutine mp_bcst_1d_i
2370 subroutine mp_bcst_2d_i(q, idim, jdim)
2371 integer,
intent(IN) :: idim, jdim
2372 integer,
intent(INOUT) :: q(idim,jdim)
2374 call mpi_bcast(q, idim*jdim, mpi_integer,
masterproc, commglobal, ierror)
2376 end subroutine mp_bcst_2d_i
2385 subroutine mp_bcst_4d_i(q, idim, jdim, kdim, ldim)
2386 integer,
intent(IN) :: idim, jdim, kdim, ldim
2387 integer,
intent(INOUT) :: q(idim,jdim,kdim,ldim)
2389 call mpi_bcast(q, idim*jdim*kdim*ldim, mpi_integer,
masterproc, commglobal, ierror)
2391 end subroutine mp_bcst_4d_i
2402 subroutine mp_reduce_max_r4_1d(mymax,npts)
2403 integer,
intent(IN) :: npts
2404 real(kind=4),
intent(INOUT) :: mymax(npts)
2406 real(kind=4) :: gmax(npts)
2408 call mpi_allreduce( mymax, gmax, npts, mpi_real, mpi_max, &
2409 commglobal, ierror )
2413 end subroutine mp_reduce_max_r4_1d
2424 subroutine mp_reduce_max_r8_1d(mymax,npts)
2425 integer,
intent(IN) :: npts
2426 real(kind=8),
intent(INOUT) :: mymax(npts)
2428 real(kind=8) :: gmax(npts)
2430 call mpi_allreduce( mymax, gmax, npts, mpi_double_precision, mpi_max, &
2431 commglobal, ierror )
2435 end subroutine mp_reduce_max_r8_1d
2446 subroutine mp_reduce_max_r4(mymax)
2447 real(kind=4),
intent(INOUT) :: mymax
2449 real(kind=4) :: gmax
2451 call mpi_allreduce( mymax, gmax, 1, mpi_real, mpi_max, &
2452 commglobal, ierror )
2456 end subroutine mp_reduce_max_r4
2463 subroutine mp_reduce_max_r8(mymax)
2464 real(kind=8),
intent(INOUT) :: mymax
2466 real(kind=8) :: gmax
2468 call mpi_allreduce( mymax, gmax, 1, mpi_double_precision, mpi_max, &
2469 commglobal, ierror )
2473 end subroutine mp_reduce_max_r8
2475 subroutine mp_reduce_min_r4(mymin)
2476 real(kind=4),
intent(INOUT) :: mymin
2478 real(kind=4) :: gmin
2480 call mpi_allreduce( mymin, gmin, 1, mpi_real, mpi_min, &
2481 commglobal, ierror )
2485 end subroutine mp_reduce_min_r4
2487 subroutine mp_reduce_min_r8(mymin)
2488 real(kind=8),
intent(INOUT) :: mymin
2490 real(kind=8) :: gmin
2492 call mpi_allreduce( mymin, gmin, 1, mpi_double_precision, mpi_min, &
2493 commglobal, ierror )
2497 end subroutine mp_reduce_min_r8
2507 subroutine mp_reduce_max_i(mymax)
2508 integer,
intent(INOUT) :: mymax
2512 call mpi_allreduce( mymax, gmax, 1, mpi_integer, mpi_max, &
2513 commglobal, ierror )
2517 end subroutine mp_reduce_max_i
2527 subroutine mp_reduce_sum_r4(mysum)
2528 real(kind=4),
intent(INOUT) :: mysum
2530 real(kind=4) :: gsum
2532 call mpi_allreduce( mysum, gsum, 1, mpi_real, mpi_sum, &
2533 commglobal, ierror )
2537 end subroutine mp_reduce_sum_r4
2547 subroutine mp_reduce_sum_r8(mysum)
2548 real(kind=8),
intent(INOUT) :: mysum
2550 real(kind=8) :: gsum
2552 call mpi_allreduce( mysum, gsum, 1, mpi_double_precision, mpi_sum, &
2553 commglobal, ierror )
2557 end subroutine mp_reduce_sum_r8
2569 subroutine mp_reduce_sum_r4_1darr(mysum, npts)
2570 integer,
intent(in) :: npts
2571 real(kind=4),
intent(inout) :: mysum(npts)
2572 real(kind=4) :: gsum(npts)
2575 call mpi_allreduce( mysum, gsum, npts, mpi_real, mpi_sum, &
2576 commglobal, ierror )
2580 end subroutine mp_reduce_sum_r4_1darr
2592 subroutine mp_reduce_sum_r4_2darr(mysum, npts1,npts2)
2593 integer,
intent(in) :: npts1,npts2
2594 real(kind=4),
intent(inout) :: mysum(npts1,npts2)
2595 real(kind=4) :: gsum(npts1,npts2)
2598 call mpi_allreduce( mysum, gsum, npts1*npts2, mpi_real, mpi_sum, &
2599 commglobal, ierror )
2603 end subroutine mp_reduce_sum_r4_2darr
2614 subroutine mp_reduce_sum_r4_1d(mysum, sum1d, npts)
2615 integer,
intent(in) :: npts
2616 real(kind=4),
intent(in) :: sum1d(npts)
2617 real(kind=4),
intent(INOUT) :: mysum
2619 real(kind=4) :: gsum
2624 mysum = mysum + sum1d(i)
2627 call mpi_allreduce( mysum, gsum, 1, mpi_double_precision, mpi_sum, &
2628 commglobal, ierror )
2632 end subroutine mp_reduce_sum_r4_1d
2642 subroutine mp_reduce_sum_r8_1d(mysum, sum1d, npts)
2643 integer,
intent(in) :: npts
2644 real(kind=8),
intent(in) :: sum1d(npts)
2645 real(kind=8),
intent(INOUT) :: mysum
2647 real(kind=8) :: gsum
2652 mysum = mysum + sum1d(i)
2655 call mpi_allreduce( mysum, gsum, 1, mpi_double_precision, mpi_sum, &
2656 commglobal, ierror )
2660 end subroutine mp_reduce_sum_r8_1d
2669 subroutine mp_reduce_sum_r8_1darr(mysum, npts)
2670 integer,
intent(in) :: npts
2671 real(kind=8),
intent(inout) :: mysum(npts)
2672 real(kind=8) :: gsum(npts)
2676 call mpi_allreduce( mysum, gsum, npts, mpi_double_precision, &
2678 commglobal, ierror )
2682 end subroutine mp_reduce_sum_r8_1darr
2693 subroutine mp_reduce_sum_r8_2darr(mysum, npts1,npts2)
2694 integer,
intent(in) :: npts1,npts2
2695 real(kind=8),
intent(inout) :: mysum(npts1,npts2)
2696 real(kind=8) :: gsum(npts1,npts2)
2700 call mpi_allreduce( mysum, gsum, npts1*npts2, &
2701 mpi_double_precision, mpi_sum, &
2702 commglobal, ierror )
2706 end subroutine mp_reduce_sum_r8_2darr
2716 integer,
parameter::
ng = 3
The module 'fv_mp_mod' is a single program multiple data (SPMD) parallel decompostion/communication m...
The module 'fv_arrays' contains the 'fv_atmos_type' and associated datatypes.
integer, parameter, public ng
integer, public masterproc