Oasis3 4.0.2
|
00001 #ifndef key_noIO 00002 !----------------------------------------------------------------------- 00003 ! Domain decomposition and domain update for message-passing codes 00004 ! 00005 ! AUTHOR: V. Balaji (vbalaji@noaa.gov) 00006 ! Princeton University/GFDL 00007 ! 00008 ! MODIFICATIONS: Reiner Vogelsang (reiner@sgi.com) 00009 ! Sophie Valcke: renamed the routine with _oa suffix 00010 ! 00011 ! This program is free software; The author agrees that you can 00012 ! redistribute and/or modify this version of the program under the 00013 ! terms of the Lesser GNU General Public License as published 00014 ! by the Free Software Foundation. 00015 ! 00016 ! This program is distributed in the hope that it will be useful, 00017 ! but WITHOUT ANY WARRANTY; without even the implied warranty of 00018 ! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the 00019 ! Lesser GNU General Public License for more details 00020 ! (http://www.gnu.org/copyleft/lesser.html). 00021 !----------------------------------------------------------------------- 00022 !#include <os.h> 00023 #if defined(_CRAYT3E) || defined(_CRAYT3D) || defined(sgi_mipspro) 00024 #define SGICRAY_MPP 00025 #endif 00026 !shmalloc is used on MPP SGI/Cray systems for shmem 00027 #if defined(use_libSMA) && defined(SGICRAY_MPP) 00028 #define use_shmalloc 00029 #endif 00030 00031 module mpp_domains_mod_oa 00032 use mod_kinds_mpp 00033 !a generalized domain decomposition package for use with mpp_mod_oa 00034 !Balaji (vb@gfdl.gov) 15 March 1999 00035 use mpp_mod_oa 00036 #include <os.h> 00037 implicit none 00038 private 00039 character(len=128), private :: version= 00040 '$Id: mpp_domains_mod_oa.F90 2826 2010-12-10 11:14:21Z valcke $' 00041 character(len=128), private :: tagname= 00042 '$Name$' 00043 character(len=128), private :: version_update_domains2D, version_global_reduce, version_global_sum, version_global_field 00044 00045 !parameters used to define domains: these are passed to the flags argument of mpp_define_domains 00046 ! if data domain is to have global extent, set GLOBAL_DATA_DOMAIN 00047 ! if global domain has periodic boundaries, set CYCLIC_GLOBAL_DOMAIN 00048 ! sum flags together if more than one of the above conditions is to be met. 00049 integer, parameter, private :: GLOBAL=0, CYCLIC=1 00050 integer, parameter, private :: WEST=2, EAST=3, SOUTH=4, NORTH=5 00051 integer, parameter, private :: SEND=1, RECV=2 00052 integer, parameter, public :: GLOBAL_DATA_DOMAIN=2**GLOBAL, CYCLIC_GLOBAL_DOMAIN=2**CYCLIC 00053 !gridtypes 00054 integer, parameter, private :: AGRID=0, BGRID=1, CGRID=2 00055 integer, parameter, public :: BGRID_NE=BGRID+2**NORTH+2**EAST 00056 integer, parameter, public :: BGRID_SW=BGRID+2**SOUTH+2**WEST 00057 integer, parameter, public :: CGRID_NE=CGRID+2**NORTH+2**EAST 00058 integer, parameter, public :: CGRID_SW=CGRID+2**SOUTH+2**WEST 00059 integer, private :: grid_offset_type=AGRID 00060 !folds 00061 integer, parameter, public :: FOLD_WEST_EDGE = 2**WEST, FOLD_EAST_EDGE = 2**EAST 00062 integer, parameter, public :: FOLD_SOUTH_EDGE=2**SOUTH, FOLD_NORTH_EDGE=2**NORTH 00063 !update 00064 integer, parameter, public :: WUPDATE=2**WEST, EUPDATE=2**EAST, SUPDATE=2**SOUTH, NUPDATE=2**NORTH 00065 integer, parameter, public :: XUPDATE=WUPDATE+EUPDATE, YUPDATE=SUPDATE+NUPDATE 00066 !used by mpp_global_sum 00067 integer, parameter, public :: BITWISE_EXACT_SUM=1 00068 00069 type, public :: domain_axis_spec !type used to specify index limits along an axis of a domain 00070 private 00071 integer :: begin, end, size, max_size !start, end of domain axis, size, max size in set 00072 logical :: is_global !TRUE if domain axis extent covers global domain 00073 end type domain_axis_spec 00074 type, public :: domain1D 00075 private 00076 type(domain_axis_spec) :: compute, data, global 00077 logical :: cyclic 00078 type(domain1D), pointer :: list(:) 00079 integer :: pe !PE to which this domain is assigned 00080 integer :: pos !position of this PE within link list, i.e domain%list(pos)%pe = pe 00081 end type domain1D 00082 !domaintypes of higher rank can be constructed from type domain1D 00083 !typically we only need 1 and 2D, but could need higher (e.g 3D LES) 00084 !some elements are repeated below if they are needed once per domain, not once per axis 00085 type, private :: rectangle 00086 integer :: is, ie, js, je 00087 logical :: overlap, folded 00088 end type rectangle 00089 type, public :: domain2D 00090 private 00091 type(domain1D) :: x 00092 type(domain1D) :: y 00093 type(domain2D), pointer :: list(:) 00094 integer :: pe !PE to which this domain is assigned 00095 integer :: pos !position of this PE within link list, i.e domain%list(pos)%pe = pe 00096 integer :: fold, gridtype 00097 logical :: overlap 00098 type(rectangle) :: recv_e, recv_se, recv_s, recv_sw, 00099 recv_w, recv_nw, recv_n, recv_ne 00100 type(rectangle) :: send_e, send_se, send_s, send_sw, 00101 send_w, send_nw, send_n, send_ne 00102 logical :: remote_domains_initialized 00103 type(rectangle) :: recv_e_off, recv_se_off, recv_s_off, recv_sw_off, 00104 recv_w_off, recv_nw_off, recv_n_off, recv_ne_off 00105 type(rectangle) :: send_e_off, send_se_off, send_s_off, send_sw_off, 00106 send_w_off, send_nw_off, send_n_off, send_ne_off 00107 logical :: remote_off_domains_initialized 00108 end type domain2D 00109 type(domain1D), public :: NULL_DOMAIN1D 00110 type(domain2D), public :: NULL_DOMAIN2D 00111 00112 integer, private :: pe 00113 00114 integer, private :: tk 00115 #ifdef DEBUG 00116 logical, private :: verbose=.TRUE., debug=.TRUE., domain_clocks_on=.FALSE. 00117 #else 00118 logical, private :: verbose=.FALSE., debug=.FALSE., domain_clocks_on=.FALSE. 00119 #endif 00120 logical, private :: module_is_initialized=.FALSE. 00121 integer, parameter, public :: MPP_DOMAIN_TIME=MPP_DEBUG+1 00122 integer :: send_clock=0, recv_clock=0, unpk_clock=0, wait_clock=0, pack_clock=0, pack_loop_clock=0 00123 00124 !stack for internal buffers 00125 !allocated differently if use_shmalloc 00126 #ifdef use_shmalloc 00127 real(DOUBLE_KIND), private :: mpp_domains_stack(1) 00128 pointer( ptr_stack, mpp_domains_stack ) 00129 #else 00130 real(DOUBLE_KIND), private, allocatable :: mpp_domains_stack(:) 00131 #endif 00132 integer, private :: mpp_domains_stack_size=0, mpp_domains_stack_hwm=0 00133 00134 !used by mpp_define_domains2D_new to transmit data 00135 integer :: domain_info_buf(16) 00136 #ifdef use_shmalloc 00137 pointer( ptr_info, domain_info_buf ) 00138 #endif 00139 00140 !public interfaces 00141 00142 ! this interface can add halo to a no-halo domain and copy everything else. 00143 00144 interface mpp_copy_domains 00145 module procedure mpp_copy_domains1D 00146 module procedure mpp_copy_domains2D 00147 end interface 00148 00149 interface mpp_define_domains 00150 module procedure mpp_define_domains1D 00151 module procedure mpp_define_domains2D 00152 end interface 00153 00154 interface mpp_update_domains 00155 module procedure mpp_update_domain2D_r8_2d 00156 module procedure mpp_update_domain2D_r8_3d 00157 module procedure mpp_update_domain2D_r8_4d 00158 module procedure mpp_update_domain2D_r8_5d 00159 module procedure mpp_update_domain2D_r8_2dv 00160 module procedure mpp_update_domain2D_r8_3dv 00161 module procedure mpp_update_domain2D_r8_4dv 00162 module procedure mpp_update_domain2D_r8_5dv 00163 module procedure mpp_update_domain2D_c8_2d 00164 module procedure mpp_update_domain2D_c8_3d 00165 module procedure mpp_update_domain2D_c8_4d 00166 module procedure mpp_update_domain2D_c8_5d 00167 #ifndef no_8byte_integers 00168 module procedure mpp_update_domain2D_i8_2d 00169 module procedure mpp_update_domain2D_i8_3d 00170 module procedure mpp_update_domain2D_i8_4d 00171 module procedure mpp_update_domain2D_i8_5d 00172 module procedure mpp_update_domain2D_l8_2d 00173 module procedure mpp_update_domain2D_l8_3d 00174 module procedure mpp_update_domain2D_l8_4d 00175 module procedure mpp_update_domain2D_l8_5d 00176 #endif 00177 #ifndef no_4byte_reals 00178 module procedure mpp_update_domain2D_r4_2d 00179 module procedure mpp_update_domain2D_r4_3d 00180 module procedure mpp_update_domain2D_r4_4d 00181 module procedure mpp_update_domain2D_r4_5d 00182 #endif 00183 #ifndef no_4byte_cmplx 00184 module procedure mpp_update_domain2D_c4_2d 00185 module procedure mpp_update_domain2D_c4_3d 00186 module procedure mpp_update_domain2D_c4_4d 00187 module procedure mpp_update_domain2D_c4_5d 00188 #endif 00189 module procedure mpp_update_domain2D_i4_2d 00190 module procedure mpp_update_domain2D_i4_3d 00191 module procedure mpp_update_domain2D_i4_4d 00192 module procedure mpp_update_domain2D_i4_5d 00193 module procedure mpp_update_domain2D_l4_2d 00194 module procedure mpp_update_domain2D_l4_3d 00195 module procedure mpp_update_domain2D_l4_4d 00196 module procedure mpp_update_domain2D_l4_5d 00197 #ifndef no_4byte_reals 00198 module procedure mpp_update_domain2D_r4_2dv 00199 module procedure mpp_update_domain2D_r4_3dv 00200 module procedure mpp_update_domain2D_r4_4dv 00201 module procedure mpp_update_domain2D_r4_5dv 00202 #endif 00203 00204 ! module procedure mpp_update_domain1D_r8_2d 00205 ! module procedure mpp_update_domain1D_r8_3d 00206 ! module procedure mpp_update_domain1D_r8_4d 00207 ! module procedure mpp_update_domain1D_r8_5d 00208 ! module procedure mpp_update_domain1D_c8_2d 00209 ! module procedure mpp_update_domain1D_c8_3d 00210 ! module procedure mpp_update_domain1D_c8_4d 00211 ! module procedure mpp_update_domain1D_c8_5d 00212 !#ifndef no_8byte_integers 00213 ! module procedure mpp_update_domain1D_i8_2d 00214 ! module procedure mpp_update_domain1D_i8_3d 00215 ! module procedure mpp_update_domain1D_i8_4d 00216 ! module procedure mpp_update_domain1D_i8_5d 00217 ! module procedure mpp_update_domain1D_l8_2d 00218 ! module procedure mpp_update_domain1D_l8_3d 00219 ! module procedure mpp_update_domain1D_l8_4d 00220 ! module procedure mpp_update_domain1D_l8_5d 00221 !#endif 00222 ! module procedure mpp_update_domain1D_r4_2d 00223 ! module procedure mpp_update_domain1D_r4_3d 00224 ! module procedure mpp_update_domain1D_r4_4d 00225 ! module procedure mpp_update_domain1D_r4_5d 00226 ! module procedure mpp_update_domain1D_c4_2d 00227 ! module procedure mpp_update_domain1D_c4_3d 00228 ! module procedure mpp_update_domain1D_c4_4d 00229 ! module procedure mpp_update_domain1D_c4_5d 00230 ! module procedure mpp_update_domain1D_i4_2d 00231 ! module procedure mpp_update_domain1D_i4_3d 00232 ! module procedure mpp_update_domain1D_i4_4d 00233 ! module procedure mpp_update_domain1D_i4_5d 00234 ! module procedure mpp_update_domain1D_l4_2d 00235 ! module procedure mpp_update_domain1D_l4_3d 00236 ! module procedure mpp_update_domain1D_l4_4d 00237 ! module procedure mpp_update_domain1D_l4_5d 00238 end interface 00239 00240 interface mpp_redistribute 00241 module procedure mpp_redistribute_r8_2D 00242 module procedure mpp_redistribute_r8_3D 00243 module procedure mpp_redistribute_r8_4D 00244 module procedure mpp_redistribute_r8_5D 00245 module procedure mpp_redistribute_c8_2D 00246 module procedure mpp_redistribute_c8_3D 00247 module procedure mpp_redistribute_c8_4D 00248 module procedure mpp_redistribute_c8_5D 00249 #ifndef no_8byte_integers 00250 module procedure mpp_redistribute_i8_2D 00251 module procedure mpp_redistribute_i8_3D 00252 module procedure mpp_redistribute_i8_4D 00253 module procedure mpp_redistribute_i8_5D 00254 module procedure mpp_redistribute_l8_2D 00255 module procedure mpp_redistribute_l8_3D 00256 module procedure mpp_redistribute_l8_4D 00257 module procedure mpp_redistribute_l8_5D 00258 #endif 00259 #ifndef no_4byte_reals 00260 module procedure mpp_redistribute_r4_2D 00261 module procedure mpp_redistribute_r4_3D 00262 module procedure mpp_redistribute_r4_4D 00263 module procedure mpp_redistribute_r4_5D 00264 #endif 00265 #ifndef no_4byte_cmplx 00266 module procedure mpp_redistribute_c4_2D 00267 module procedure mpp_redistribute_c4_3D 00268 module procedure mpp_redistribute_c4_4D 00269 module procedure mpp_redistribute_c4_5D 00270 #endif 00271 module procedure mpp_redistribute_i4_2D 00272 module procedure mpp_redistribute_i4_3D 00273 module procedure mpp_redistribute_i4_4D 00274 module procedure mpp_redistribute_i4_5D 00275 module procedure mpp_redistribute_l4_2D 00276 module procedure mpp_redistribute_l4_3D 00277 module procedure mpp_redistribute_l4_4D 00278 module procedure mpp_redistribute_l4_5D 00279 end interface 00280 00281 interface mpp_global_field 00282 module procedure mpp_global_field2D_r8_2d 00283 module procedure mpp_global_field2D_r8_3d 00284 module procedure mpp_global_field2D_r8_4d 00285 module procedure mpp_global_field2D_r8_5d 00286 module procedure mpp_global_field2D_c8_2d 00287 module procedure mpp_global_field2D_c8_3d 00288 module procedure mpp_global_field2D_c8_4d 00289 module procedure mpp_global_field2D_c8_5d 00290 #ifndef no_8byte_integers 00291 module procedure mpp_global_field2D_i8_2d 00292 module procedure mpp_global_field2D_i8_3d 00293 module procedure mpp_global_field2D_i8_4d 00294 module procedure mpp_global_field2D_i8_5d 00295 module procedure mpp_global_field2D_l8_2d 00296 module procedure mpp_global_field2D_l8_3d 00297 module procedure mpp_global_field2D_l8_4d 00298 module procedure mpp_global_field2D_l8_5d 00299 #endif 00300 #ifndef no_4byte_reals 00301 module procedure mpp_global_field2D_r4_2d 00302 module procedure mpp_global_field2D_r4_3d 00303 module procedure mpp_global_field2D_r4_4d 00304 module procedure mpp_global_field2D_r4_5d 00305 #endif 00306 #ifndef no_4byte_cmplx 00307 module procedure mpp_global_field2D_c4_2d 00308 module procedure mpp_global_field2D_c4_3d 00309 module procedure mpp_global_field2D_c4_4d 00310 module procedure mpp_global_field2D_c4_5d 00311 #endif 00312 module procedure mpp_global_field2D_i4_2d 00313 module procedure mpp_global_field2D_i4_3d 00314 module procedure mpp_global_field2D_i4_4d 00315 module procedure mpp_global_field2D_i4_5d 00316 module procedure mpp_global_field2D_l4_2d 00317 module procedure mpp_global_field2D_l4_3d 00318 module procedure mpp_global_field2D_l4_4d 00319 module procedure mpp_global_field2D_l4_5d 00320 module procedure mpp_global_field1D_r8_2d 00321 module procedure mpp_global_field1D_c8_2d 00322 #ifndef no_8byte_integers 00323 module procedure mpp_global_field1D_i8_2d 00324 module procedure mpp_global_field1D_l8_2d 00325 #endif 00326 #ifndef no_4byte_reals 00327 module procedure mpp_global_field1D_r4_2d 00328 #endif 00329 #ifndef no_4byte_cmplx 00330 module procedure mpp_global_field1D_c4_2d 00331 #endif 00332 module procedure mpp_global_field1D_i4_2d 00333 module procedure mpp_global_field1D_l4_2d 00334 end interface 00335 00336 interface mpp_global_max 00337 module procedure mpp_global_max_r8_2d 00338 module procedure mpp_global_max_r8_3d 00339 module procedure mpp_global_max_r8_4d 00340 module procedure mpp_global_max_r8_5d 00341 #ifndef no_4byte_reals 00342 module procedure mpp_global_max_r4_2d 00343 module procedure mpp_global_max_r4_3d 00344 module procedure mpp_global_max_r4_4d 00345 module procedure mpp_global_max_r4_5d 00346 #endif 00347 #ifndef no_8byte_integers 00348 module procedure mpp_global_max_i8_2d 00349 module procedure mpp_global_max_i8_3d 00350 module procedure mpp_global_max_i8_4d 00351 module procedure mpp_global_max_i8_5d 00352 #endif 00353 module procedure mpp_global_max_i4_2d 00354 module procedure mpp_global_max_i4_3d 00355 module procedure mpp_global_max_i4_4d 00356 module procedure mpp_global_max_i4_5d 00357 end interface 00358 00359 interface mpp_global_min 00360 module procedure mpp_global_min_r8_2d 00361 module procedure mpp_global_min_r8_3d 00362 module procedure mpp_global_min_r8_4d 00363 module procedure mpp_global_min_r8_5d 00364 #ifndef no_4byte_reals 00365 module procedure mpp_global_min_r4_2d 00366 module procedure mpp_global_min_r4_3d 00367 module procedure mpp_global_min_r4_4d 00368 module procedure mpp_global_min_r4_5d 00369 #endif 00370 #ifndef no_8byte_integers 00371 module procedure mpp_global_min_i8_2d 00372 module procedure mpp_global_min_i8_3d 00373 module procedure mpp_global_min_i8_4d 00374 module procedure mpp_global_min_i8_5d 00375 #endif 00376 module procedure mpp_global_min_i4_2d 00377 module procedure mpp_global_min_i4_3d 00378 module procedure mpp_global_min_i4_4d 00379 module procedure mpp_global_min_i4_5d 00380 end interface 00381 00382 interface mpp_global_sum 00383 module procedure mpp_global_sum_r8_2d 00384 module procedure mpp_global_sum_r8_3d 00385 module procedure mpp_global_sum_r8_4d 00386 module procedure mpp_global_sum_r8_5d 00387 #ifndef no_4byte_reals 00388 module procedure mpp_global_sum_r4_2d 00389 module procedure mpp_global_sum_r4_3d 00390 module procedure mpp_global_sum_r4_4d 00391 module procedure mpp_global_sum_r4_5d 00392 #endif 00393 module procedure mpp_global_sum_c8_2d 00394 module procedure mpp_global_sum_c8_3d 00395 module procedure mpp_global_sum_c8_4d 00396 module procedure mpp_global_sum_c8_5d 00397 #ifndef no_4byte_cmplx 00398 module procedure mpp_global_sum_c4_2d 00399 module procedure mpp_global_sum_c4_3d 00400 module procedure mpp_global_sum_c4_4d 00401 module procedure mpp_global_sum_c4_5d 00402 #endif 00403 #ifndef no_8byte_integers 00404 module procedure mpp_global_sum_i8_2d 00405 module procedure mpp_global_sum_i8_3d 00406 module procedure mpp_global_sum_i8_4d 00407 module procedure mpp_global_sum_i8_5d 00408 #endif 00409 module procedure mpp_global_sum_i4_2d 00410 module procedure mpp_global_sum_i4_3d 00411 module procedure mpp_global_sum_i4_4d 00412 module procedure mpp_global_sum_i4_5d 00413 end interface 00414 00415 interface operator(.EQ.) 00416 module procedure mpp_domain1D_eq 00417 module procedure mpp_domain2D_eq 00418 end interface 00419 00420 interface operator(.NE.) 00421 module procedure mpp_domain1D_ne 00422 module procedure mpp_domain2D_ne 00423 end interface 00424 00425 interface mpp_get_compute_domain 00426 module procedure mpp_get_compute_domain1D 00427 module procedure mpp_get_compute_domain2D 00428 end interface 00429 00430 interface mpp_get_compute_domains 00431 module procedure mpp_get_compute_domains1D 00432 module procedure mpp_get_compute_domains2D 00433 end interface 00434 00435 interface mpp_get_data_domain 00436 module procedure mpp_get_data_domain1D 00437 module procedure mpp_get_data_domain2D 00438 end interface 00439 00440 interface mpp_get_global_domain 00441 module procedure mpp_get_global_domain1D 00442 module procedure mpp_get_global_domain2D 00443 end interface 00444 00445 interface mpp_define_layout 00446 module procedure mpp_define_layout2D 00447 end interface 00448 00449 interface mpp_get_pelist 00450 module procedure mpp_get_pelist1D 00451 module procedure mpp_get_pelist2D 00452 end interface 00453 00454 interface mpp_get_layout 00455 module procedure mpp_get_layout1D 00456 module procedure mpp_get_layout2D 00457 end interface 00458 00459 public :: mpp_broadcast_domain, mpp_define_layout, mpp_define_domains, mpp_domains_init, mpp_domains_set_stack_size, & 00460 mpp_domains_exit, mpp_get_compute_domain, mpp_get_compute_domains, mpp_get_data_domain, mpp_get_global_domain, & 00461 mpp_get_domain_components, mpp_get_layout, mpp_get_pelist, mpp_redistribute, mpp_update_domains, & 00462 mpp_global_field, mpp_global_max, mpp_global_min, mpp_global_sum, operator(.EQ.), operator(.NE.), & 00463 mpp_copy_domains 00464 00465 contains 00466 00467 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00468 ! ! 00469 ! MPP_DOMAINS: initialization and termination ! 00470 ! ! 00471 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00472 00473 subroutine mpp_domains_init(flags) 00474 !initialize domain decomp package 00475 integer, intent(in), optional :: flags 00476 integer :: l=0 00477 00478 if( module_is_initialized )return 00479 call mpp_init(flags) !this is a no-op if already initialized 00480 pe = mpp_pe() 00481 module_is_initialized = .TRUE. 00482 if( pe.EQ.mpp_root_pe() )then 00483 write( nulprt,'(/a)' )'MPP_DOMAINS module '//trim(version) & 00484 //trim(tagname) 00485 ! write( nulprt,'(a)' )trim(version_update_domains2D) 00486 end if 00487 00488 if( PRESENT(flags) )then 00489 debug = flags.EQ.MPP_DEBUG 00490 verbose = flags.EQ.MPP_VERBOSE .OR. debug 00491 domain_clocks_on = flags.EQ.MPP_DOMAIN_TIME 00492 end if 00493 00494 call mpp_domains_set_stack_size(983040) !default, pretty arbitrary 00495 #ifdef use_shmalloc 00496 call mpp_malloc( ptr_info, 16, l ) 00497 #endif 00498 00499 !NULL_DOMAIN is a domaintype that can be used to initialize to undef 00500 NULL_DOMAIN1D%global%begin = -1; NULL_DOMAIN1D%global%end = -1; NULL_DOMAIN1D%global%size = 0 00501 NULL_DOMAIN1D%data%begin = -1; NULL_DOMAIN1D%data%end = -1; NULL_DOMAIN1D%data%size = 0 00502 NULL_DOMAIN1D%compute%begin = -1; NULL_DOMAIN1D%compute%end = -1; NULL_DOMAIN1D%compute%size = 0 00503 NULL_DOMAIN1D%pe = ANY_PE 00504 NULL_DOMAIN2D%x = NULL_DOMAIN1D 00505 NULL_DOMAIN2D%y = NULL_DOMAIN1D 00506 NULL_DOMAIN2D%pe = ANY_PE 00507 00508 if( domain_clocks_on )then 00509 pack_clock = mpp_clock_id( 'Halo pack' ) 00510 pack_loop_clock = mpp_clock_id( 'Halo pack loop' ) 00511 send_clock = mpp_clock_id( 'Halo send' ) 00512 recv_clock = mpp_clock_id( 'Halo recv' ) 00513 unpk_clock = mpp_clock_id( 'Halo unpk' ) 00514 wait_clock = mpp_clock_id( 'Halo wait' ) 00515 end if 00516 return 00517 end subroutine mpp_domains_init 00518 00519 subroutine mpp_domains_set_stack_size(n) 00520 !set the mpp_domains_stack variable to be at least n LONG words long 00521 integer, intent(in) :: n 00522 character(len=8) :: text 00523 00524 if( n.LE.mpp_domains_stack_size )return 00525 #ifdef use_shmalloc 00526 call mpp_malloc( ptr_stack, n, mpp_domains_stack_size ) 00527 #else 00528 if( allocated(mpp_domains_stack) )deallocate(mpp_domains_stack) 00529 allocate( mpp_domains_stack(n) ) 00530 mpp_domains_stack_size = n 00531 #endif 00532 write( text,'(i8)' )n 00533 if( pe.EQ.mpp_root_pe() )call mpp_error( NOTE, 'MPP_DOMAINS_SET_STACK_SIZE: stack size set to '//text//'.' ) 00534 00535 return 00536 end subroutine mpp_domains_set_stack_size 00537 00538 subroutine mpp_domains_exit() 00539 !currently does not have much to do, but provides the possibility of re-initialization 00540 if( .NOT.module_is_initialized )return 00541 call mpp_max(mpp_domains_stack_hwm) 00542 if( pe.EQ.mpp_root_pe() )write( nulprt,* )'MPP_DOMAINS_STACK high water mark=', mpp_domains_stack_hwm 00543 module_is_initialized = .FALSE. 00544 return 00545 end subroutine mpp_domains_exit 00546 00547 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00548 ! ! 00549 ! MPP_DOMAINS: overloaded operators (==, /=) ! 00550 ! ! 00551 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00552 00553 function mpp_domain1D_eq( a, b ) 00554 logical :: mpp_domain1D_eq 00555 type(domain1D), intent(in) :: a, b 00556 00557 mpp_domain1D_eq = ( a%compute%begin.EQ.b%compute%begin .AND. & 00558 a%compute%end .EQ.b%compute%end .AND. & 00559 a%data%begin .EQ.b%data%begin .AND. & 00560 a%data%end .EQ.b%data%end .AND. & 00561 a%global%begin .EQ.b%global%begin .AND. & 00562 a%global%end .EQ.b%global%end ) 00563 !compare pelists 00564 ! if( mpp_domain1D_eq )mpp_domain1D_eq = ASSOCIATED(a%list) .AND. ASSOCIATED(b%list) 00565 ! if( mpp_domain1D_eq )mpp_domain1D_eq = size(a%list).EQ.size(b%list) 00566 ! if( mpp_domain1D_eq )mpp_domain1D_eq = ALL(a%list%pe.EQ.b%list%pe) 00567 00568 return 00569 end function mpp_domain1D_eq 00570 00571 function mpp_domain1D_ne( a, b ) 00572 logical :: mpp_domain1D_ne 00573 type(domain1D), intent(in) :: a, b 00574 00575 mpp_domain1D_ne = .NOT. ( a.EQ.b ) 00576 return 00577 end function mpp_domain1D_ne 00578 00579 function mpp_domain2D_eq( a, b ) 00580 logical :: mpp_domain2D_eq 00581 type(domain2D), intent(in) :: a, b 00582 00583 mpp_domain2D_eq = a%x.EQ.b%x .AND. a%y.EQ.b%y 00584 if( mpp_domain2D_eq .AND. ((a%pe.EQ.ANY_PE).OR.(b%pe.EQ.ANY_PE)) )return !NULL_DOMAIN2D 00585 !compare pelists 00586 if( mpp_domain2D_eq )mpp_domain2D_eq = ASSOCIATED(a%list) .AND. ASSOCIATED(b%list) 00587 if( mpp_domain2D_eq )mpp_domain2D_eq = size(a%list).EQ.size(b%list) 00588 if( mpp_domain2D_eq )mpp_domain2D_eq = ALL(a%list%pe.EQ.b%list%pe) 00589 return 00590 end function mpp_domain2D_eq 00591 00592 function mpp_domain2D_ne( a, b ) 00593 logical :: mpp_domain2D_ne 00594 type(domain2D), intent(in) :: a, b 00595 00596 mpp_domain2D_ne = .NOT. ( a.EQ.b ) 00597 return 00598 end function mpp_domain2D_ne 00599 00600 00601 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00602 ! ! 00603 ! MPP_COPY_DOMAINS: Copy domain and add halo ! 00604 ! ! 00605 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00606 00607 subroutine mpp_copy_domains1D(domain_in, domain_out, halo) 00608 ! ARGUMENTS 00609 ! domain_in is the input domain 00610 ! domain_out is the returned domain 00611 ! halo (optional) defines halo width (currently the same on both sides) 00612 00613 type(domain1D), intent(in) :: domain_in 00614 type(domain1D), intent(out) :: domain_out 00615 integer, intent(in), optional :: halo 00616 00617 integer, dimension(:), allocatable :: pelist, extent 00618 integer :: ndivs, global_indices(2) !(/ isg, ieg /) 00619 integer :: isg, ieg, halosz 00620 00621 00622 ! get the global indices of the input domain 00623 call mpp_get_global_domain(domain_in, isg, ieg ) 00624 global_indices(1) = isg; global_indices(2) = ieg 00625 00626 ! get the number of divisions of domain 00627 call mpp_get_layout(domain_in, ndivs) 00628 00629 allocate(pelist(0:ndivs-1),extent(0:ndivs-1)) 00630 00631 ! get the pe list of the input domain 00632 call mpp_get_pelist( domain_in, pelist) 00633 00634 ! get the halo 00635 halosz = 0 00636 if(present(halo)) halosz = halo 00637 00638 ! get the extent 00639 call mpp_get_compute_domains(Domain_in, size = extent(0:ndivs-1)) 00640 00641 call mpp_define_domains( global_indices, ndivs, domain_out, pelist = pelist, & 00642 halo = halosz, extent = extent ) 00643 00644 end subroutine mpp_copy_domains1D 00645 00646 00647 !---------------------------------------------------------------------------------- 00648 00649 00650 subroutine mpp_copy_domains2D(domain_in, domain_out, xhalo, yhalo) 00651 ! ARGUMENTS 00652 ! domain_in is the input domain 00653 ! domain_out is the returned domain 00654 ! xhalo, yhalo (optional) defines halo width. 00655 00656 type(domain2D), intent(in) :: domain_in 00657 type(domain2D), intent(out) :: domain_out 00658 integer, intent(in), optional :: xhalo, yhalo 00659 00660 integer, dimension(:), allocatable :: pelist, xextent, yextent, xbegin, xend, 00661 xsize, ybegin, yend, ysize 00662 integer :: ndivx, ndivy, ndivs, npes, global_indices(4), layout(2) 00663 integer :: isg, ieg, jsg, jeg, xhalosz, yhalosz, i, j, m, n, pe 00664 00665 ! get the global indices of the input domain 00666 call mpp_get_global_domain(domain_in, isg, ieg, jsg, jeg ) 00667 global_indices(1) = isg; global_indices(2) = ieg 00668 global_indices(3) = jsg; global_indices(4) = jeg 00669 00670 ! get the number of divisions of domain 00671 call mpp_get_layout(domain_in, layout) 00672 ndivx = layout(1); ndivy = layout(2) 00673 ndivs = ndivx * ndivy 00674 npes = mpp_npes() 00675 00676 if( ndivs.NE.npes )call mpp_error( 'mpp_domains_mod_oa', & 00677 'mpp_copy_domains: number of PEs is not consistent with the layout', FATAL ) 00678 allocate(pelist(0:npes-1), xsize(0:npes-1), ysize(0:npes-1), & 00679 xbegin(0:npes-1), xend(0:npes-1), & 00680 ybegin(0:npes-1), yend(0:npes-1), & 00681 xextent(0:ndivx), yextent(0:ndivy)) 00682 00683 ! get the pe list of the input domain 00684 call mpp_get_pelist2D( domain_in, pelist) 00685 00686 ! get the halo 00687 xhalosz = 0; yhalosz = 0 00688 if(present(xhalo)) xhalosz = xhalo 00689 if(present(yhalo)) yhalosz = yhalo 00690 00691 ! get the extent 00692 call mpp_get_compute_domains(Domain_in, xbegin, xend, xsize, & 00693 ybegin, yend, ysize ) 00694 00695 ! compute the exact x-axis decomposition 00696 i = isg 00697 m = 0 00698 xextent = 0 00699 do pe = 0, npes-1 00700 if ( xbegin(pe) == i ) then 00701 m = m+1 00702 xextent (m) = xsize (pe) 00703 i = i + xsize(pe) 00704 if ( m == size(xextent) .or. i > ieg ) exit 00705 endif 00706 enddo 00707 00708 ! compute the exact y-axis decomposition 00709 j = jsg 00710 n = 0 00711 yextent = 0 00712 do pe = 0, npes-1 00713 if ( ybegin(pe) == j ) then 00714 n = n+1 00715 yextent (n) = ysize (pe) 00716 j = j + ysize(pe) 00717 if ( n == size(yextent) .or. j > jeg ) exit 00718 endif 00719 enddo 00720 00721 call mpp_define_domains( global_indices, layout, domain_out, pelist = pelist, & 00722 xhalo = xhalosz, yhalo = yhalosz, xextent = xextent, & 00723 yextent = yextent) 00724 00725 end subroutine mpp_copy_domains2D 00726 00727 00728 00729 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00730 ! ! 00731 ! MPP_DEFINE_DOMAINS: define layout and decomposition ! 00732 ! ! 00733 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00734 00735 subroutine mpp_define_domains1D( global_indices, ndivs, domain, pelist, & 00736 flags, halo, extent, maskmap,offset) 00737 !routine to divide global array indices among domains, and assign domains to PEs 00738 !domain is of type domain1D 00739 !ARGUMENTS: 00740 ! global_indices(2)=(isg,ieg) gives the extent of global domain 00741 ! ndivs is number of divisions of domain: even divisions unless extent is present. 00742 ! domain is the returned domain1D 00743 ! pelist (optional) list of PEs to which domains are to be assigned (default 0...npes-1) 00744 ! size of pelist must correspond to number of mask=.TRUE. divisions 00745 ! flags define whether compute and data domains are global (undecomposed) and whether global domain has periodic boundaries 00746 ! halo (optional) defines halo width (currently the same on both sides) 00747 ! extent (optional) array defines width of each division (used for non-uniform domain decomp, for e.g load-balancing) 00748 ! maskmap (optional) a division whose maskmap=.FALSE. is not assigned to any domain 00749 ! By default we assume decomposition of compute and data domains, non-periodic boundaries, no halo, as close to uniform extents 00750 ! as the input parameters permit 00751 integer, intent(in) :: global_indices(2) !(/ isg, ieg /) 00752 integer, intent(in) :: ndivs 00753 type(domain1D), intent(inout) :: domain !declared inout so that existing links, if any, can be nullified 00754 integer, intent(in), optional :: pelist(0:) 00755 integer, intent(in), optional :: flags, halo 00756 integer, intent(in), optional :: extent(0:) 00757 logical, intent(in), optional :: maskmap(0:) 00758 integer, intent(in), optional :: offset(0:) 00759 00760 logical :: compute_domain_is_global, data_domain_is_global 00761 integer :: ndiv, n, isg, ieg, is, ie, i, off, pos, hs, he 00762 integer, allocatable :: pes(:) 00763 logical, allocatable :: mask(:) 00764 !rv 00765 logical ::blocks 00766 integer :: lastie,il_compute_max_size 00767 !rv 00768 integer :: halosz 00769 !used by symmetry algorithm 00770 integer :: imax, ndmax, ndmirror 00771 logical :: symmetrize 00772 !statement functions 00773 logical :: even, odd 00774 even(n) = (mod(n,2).EQ.0) 00775 odd (n) = (mod(n,2).EQ.1) 00776 00777 if( .NOT.module_is_initialized )call mpp_error( FATAL, 'MPP_DEFINE_DOMAINS1D: You must first call mpp_domains_init.' ) 00778 !get global indices 00779 !RV 00780 is=0 00781 ie=0 00782 !RV 00783 isg = global_indices(1) 00784 ieg = global_indices(2) 00785 if( ndivs.GT.ieg-isg+1 )call mpp_error( WARNING, 'MPP_DEFINE_DOMAINS1D: more divisions requested than rows available.' ) 00786 !get the list of PEs on which to assign domains; if pelist is absent use 0..npes-1 00787 if( PRESENT(pelist) )then 00788 if( .NOT.any(pelist.EQ.pe) )then 00789 write( stderr(),* )'pe=', pe, ' pelist=', pelist 00790 call mpp_error( FATAL, 'MPP_DEFINE_DOMAINS1D: pe must be in pelist.' ) 00791 end if 00792 allocate( pes(0:size(pelist)-1) ) 00793 pes(:) = pelist(:) 00794 else 00795 allocate( pes(0:mpp_npes()-1) ) 00796 pes(:) = (/ (i,i=0,mpp_npes()-1) /) 00797 end if 00798 00799 !rv Block distribution ? 00800 blocks=.false. 00801 if(present(offset)) blocks=.true. 00802 if(blocks) then 00803 if((.not.present(maskmap)).and.(.not.present(extent))) & 00804 call mpp_error(FATAL, 'MPP_DEFINE_DOMAINS1D: you have to define maskmap and extent!') 00805 if(size(offset).ne.ndivs) & 00806 call mpp_error( FATAL, 'MPP_DEFINE_DOMAINS1D: offset array size must equal number of domain divisions.' ) 00807 endif 00808 !rv 00809 !get number of real domains: 1 mask domain per PE in pes 00810 allocate( mask(0:ndivs-1) ) 00811 mask = .TRUE. !default mask 00812 if( PRESENT(maskmap) )then 00813 if( size(maskmap).NE.ndivs ) & 00814 call mpp_error( FATAL, 'MPP_DEFINE_DOMAINS1D: maskmap array size must equal number of domain divisions.' ) 00815 mask(:) = maskmap(:) 00816 end if 00817 if( count(mask).NE.size(pes) ) & 00818 call mpp_error( FATAL, 'MPP_DEFINE_DOMAINS1D: number of TRUEs in maskmap array must match PE count.' ) 00819 00820 if( PRESENT(extent) )then 00821 if( size(extent).NE.ndivs ) & 00822 call mpp_error( FATAL, 'MPP_DEFINE_DOMAINS1D: extent array size must equal number of domain divisions.' ) 00823 end if 00824 !get halosize 00825 halosz = 0 00826 if( PRESENT(halo) )halosz = halo 00827 00828 !get flags 00829 compute_domain_is_global = .FALSE. 00830 data_domain_is_global = .FALSE. 00831 domain%cyclic = .FALSE. 00832 if( PRESENT(flags) )then 00833 !NEW: obsolete flag global_compute_domain, since ndivs is non-optional and you cannot have global compute and ndivs.NE.1 00834 compute_domain_is_global = ndivs.EQ.1 00835 !if compute domain is global, data domain must also be 00836 data_domain_is_global = BTEST(flags,GLOBAL) .OR. compute_domain_is_global 00837 domain%cyclic = BTEST(flags,CYCLIC) .AND. halosz.NE.0 00838 end if 00839 00840 !set up links list 00841 allocate( domain%list(0:ndivs-1) ) 00842 00843 !set global domain 00844 domain%list(:)%global%begin = isg 00845 domain%list(:)%global%end = ieg 00846 domain%list(:)%global%size = ieg-isg+1 00847 domain%list(:)%global%max_size = ieg-isg+1 00848 domain%list(:)%global%is_global = .TRUE. !always 00849 00850 !get compute domain 00851 if( compute_domain_is_global )then 00852 domain%list(:)%compute%begin = isg 00853 domain%list(:)%compute%end = ieg 00854 domain%list(:)%compute%is_global = .TRUE. 00855 domain%list(:)%pe = pes(:) 00856 domain%pos = 0 00857 else 00858 domain%list(:)%compute%is_global = .FALSE. 00859 is = isg 00860 lastie=-1 00861 n = 0 00862 !RV 00863 il_compute_max_size=0 00864 !RV 00865 do ndiv=0,ndivs-1 00866 if( PRESENT(extent) )then 00867 !rv build the grid by individual block sizes 00868 if(blocks) then 00869 if(mask(ndiv)) then 00870 is=isg+offset(ndiv) 00871 ie = is + extent(ndiv) - 1 00872 ! print*,'MPP_DEFINE_DOMAINS: ',pe,':',ndiv,blocks,is,ie,ieg,offset 00873 endif 00874 ! if( ndiv.EQ.ndivs-1 .AND. ie.NE.ieg ) & 00875 ! call mpp_error( FATAL, 'MPP_DEFINE_DOMAINS: extent array limits do not match global domain.' ) 00876 else 00877 ie = is + extent(ndiv) - 1 00878 if( ndiv.EQ.ndivs-1 .AND. ie.NE.ieg ) & 00879 call mpp_error( FATAL, 'MPP_DEFINE_DOMAINS: extent array limits do not match global domain.' ) 00880 endif 00881 else 00882 !modified for mirror-symmetry 00883 !original line 00884 ! ie = is + CEILING( float(ieg-is+1)/(ndivs-ndiv) ) - 1 00885 00886 !problem of dividing nx points into n domains maintaining symmetry 00887 !i.e nx=18 n=4 4554 and 5445 are solutions but 4455 is not. 00888 !this will always work for nx even n even or odd 00889 !this will always work for nx odd, n odd 00890 !this will never work for nx odd, n even: for this case we supersede the mirror calculation 00891 ! symmetrize = .NOT. ( mod(ndivs,2).EQ.0 .AND. mod(ieg-isg+1,2).EQ.1 ) 00892 !nx even n odd fails if n>nx/2 00893 symmetrize = ( even(ndivs) .AND. even(ieg-isg+1) ) .OR. & 00894 ( odd(ndivs) .AND. odd(ieg-isg+1) ) .OR. & 00895 ( odd(ndivs) .AND. even(ieg-isg+1) .AND. ndivs.LT.(ieg-isg+1)/2 ) 00896 00897 !mirror domains are stored in the list and retrieved if required. 00898 if( ndiv.EQ.0 )then 00899 !initialize max points and max domains 00900 imax = ieg 00901 ndmax = ndivs 00902 end if 00903 !do bottom half of decomposition, going over the midpoint for odd ndivs 00904 if( ndiv.LT.(ndivs-1)/2+1 )then 00905 !domain is sized by dividing remaining points by remaining domains 00906 ie = is + CEILING( REAL(imax-is+1)/(ndmax-ndiv) ) - 1 00907 ndmirror = (ndivs-1) - ndiv !mirror domain 00908 if( ndmirror.GT.ndiv .AND. symmetrize )then !only for domains over the midpoint 00909 !mirror extents, the max(,) is to eliminate overlaps 00910 domain%list(ndmirror)%compute%begin = max( isg+ieg-ie, ie+1 ) 00911 domain%list(ndmirror)%compute%end = max( isg+ieg-is, ie+1 ) 00912 imax = domain%list(ndmirror)%compute%begin - 1 00913 ndmax = ndmax - 1 00914 end if 00915 else 00916 if( symmetrize )then 00917 !do top half of decomposition by retrieving saved values 00918 is = domain%list(ndiv)%compute%begin 00919 ie = domain%list(ndiv)%compute%end 00920 else 00921 ie = is + CEILING( REAL(imax-is+1)/(ndmax-ndiv) ) - 1 00922 end if 00923 end if 00924 end if 00925 if(blocks) then 00926 if(mask(ndiv))then 00927 if( ie.LT.is )call mpp_error( FATAL, 'MPP_DEFINE_DOMAINS: domain extents must be positive definite.' ) 00928 domain%list(ndiv)%compute%begin = is 00929 domain%list(ndiv)%compute%end = ie 00930 ! if(lastie.gt.0.and.is.ne.lastie+1) & 00931 ! call mpp_error( FATAL, 'MPP_DEFINE_DOMAINS: block extents do not span space completely.' ) 00932 lastie=ie 00933 domain%list(ndiv)%pe = pes(n) 00934 if( pe.EQ.pes(n) )domain%pos = ndiv 00935 n = n + 1 00936 is = ie + 1 00937 else 00938 !rv for Arnaud 00939 domain%list(ndiv)%compute%begin = 0 00940 domain%list(ndiv)%compute%end = 0 00941 endif 00942 if(ndiv.eq.ndivs-1.and.lastie.gt.ieg) & 00943 call mpp_error( FATAL, 'MPP_DEFINE_DOMAINS: block extents exceed max. space.' ) 00944 else 00945 if( ie.LT.is )call mpp_error( FATAL, 'MPP_DEFINE_DOMAINS: domain extents must be positive definite.' ) 00946 domain%list(ndiv)%compute%begin = is 00947 domain%list(ndiv)%compute%end = ie 00948 !rr 2005-07-20, Luis Kornblueh, MPIMet 00949 if( ndiv.GT.0 ) THEN 00950 if( is.NE.domain%list(ndiv-1)%compute%end+1 ) & 00951 call mpp_error( FATAL, 'MPP_DEFINE_DOMAINS: domain extents do not span space completely.' ) 00952 end if 00953 if( ndiv.EQ.ndivs-1 .AND. domain%list(ndiv)%compute%end.NE.ieg ) & 00954 call mpp_error( FATAL, 'MPP_DEFINE_DOMAINS: domain extents do not span space completely.' ) 00955 if( mask(ndiv) )then 00956 domain%list(ndiv)%pe = pes(n) 00957 if( pe.EQ.pes(n) )domain%pos = ndiv 00958 n = n + 1 00959 end if 00960 is = ie + 1 00961 endif 00962 end do 00963 end if 00964 00965 domain%list(:)%compute%size = domain%list(:)%compute%end - domain%list(:)%compute%begin + 1 00966 !get data domain 00967 !data domain is at least equal to compute domain 00968 domain%list(:)%data%begin = domain%list(:)%compute%begin 00969 domain%list(:)%data%end = domain%list(:)%compute%end 00970 domain%list(:)%data%is_global = .FALSE. 00971 !apply global flags 00972 if( data_domain_is_global )then 00973 domain%list(:)%data%begin = isg 00974 domain%list(:)%data%end = ieg 00975 domain%list(:)%data%is_global = .TRUE. 00976 end if 00977 !apply margins 00978 domain%list(:)%data%begin = domain%list(:)%data%begin - halosz 00979 domain%list(:)%data%end = domain%list(:)%data%end + halosz 00980 domain%list(:)%data%size = domain%list(:)%data%end - domain%list(:)%data%begin + 1 00981 00982 ! domain = domain%list(pos) !load domain from domain%list(pos) 00983 domain%compute = domain%list(domain%pos)%compute 00984 domain%data = domain%list(domain%pos)%data 00985 domain%global = domain%list(domain%pos)%global 00986 domain%compute%max_size = MAXVAL( domain%list(:)%compute%size ) 00987 !RV 00988 if(blocks) then 00989 call mpp_max(domain%compute%max_size) 00990 endif 00991 !RV 00992 domain%data%max_size = MAXVAL( domain%list(:)%data%size ) 00993 domain%global%max_size = domain%global%size 00994 00995 !PV786667: the deallocate stmts can be removed when fixed (7.3.1.3m) 00996 deallocate( pes, mask ) 00997 return 00998 00999 contains 01000 01001 function if_overlap( hs, he, cs, ce, os, oe ) 01002 !function to look if a portion of the halo [hs,he] lies with in the compute region [cs,ce] 01003 !if yes, if_overlap returns true, and the overlap region is returned in [os,oe] 01004 logical :: if_overlap 01005 integer, intent(in) :: hs, he, cs, ce 01006 integer, intent(out) :: os, oe 01007 os = max(hs,cs) 01008 oe = min(he,ce) 01009 if( debug )write( stderr(),'(a,7i4)' ) & 01010 'MPP_DEFINE_DOMAINS1D: pe, hs, he, cs, ce, os, oe=', pe, hs, he, cs, ce, os, oe 01011 if_overlap = oe.GE.os 01012 return 01013 end function if_overlap 01014 01015 end subroutine mpp_define_domains1D 01016 01017 subroutine mpp_define_domains2D( global_indices, layout, domain, pelist, & 01018 xflags, yflags, & 01019 xhalo, yhalo, xextent, yextent, maskmap,offsetx,offsety, name ) 01020 !define 2D data and computational domain on global rectilinear cartesian domain (isg:ieg,jsg:jeg) and assign them to PEs 01021 integer, intent(in) :: global_indices(4) !(/ isg, ieg, jsg, jeg /) 01022 integer, intent(in) :: layout(2) 01023 type(domain2D), intent(inout) :: domain 01024 integer, intent(in), optional :: pelist(0:) 01025 integer, intent(in), optional :: xflags, yflags, xhalo, yhalo 01026 integer, intent(in), optional :: xextent(0:), yextent(0:) 01027 logical, intent(in), optional :: maskmap(0:,0:) 01028 integer, intent(in), optional :: offsetx(0:), offsety(0:) 01029 character(len=*), intent(in), optional :: name 01030 integer :: i, j, m, n 01031 integer :: ipos, jpos, pos 01032 integer :: ndivx, ndivy, isg, ieg, jsg, jeg, isd, ied, jsd, jed 01033 01034 logical, allocatable :: mask(:,:) 01035 !rv 01036 logical :: blocks 01037 !rv 01038 integer, allocatable :: pes(:), pearray(:,:) 01039 character(len=8) :: text 01040 01041 if( .NOT.module_is_initialized )call mpp_error( FATAL, 'MPP_DEFINE_DOMAINS2D: You must first call mpp_domains_init.' ) 01042 ndivx = layout(1); ndivy = layout(2) 01043 isg = global_indices(1); ieg = global_indices(2); jsg = global_indices(3); jeg = global_indices(4) 01044 01045 if( PRESENT(pelist) )then 01046 allocate( pes(0:size(pelist)-1) ) 01047 pes = pelist 01048 else 01049 allocate( pes(0:mpp_npes()-1) ) 01050 call mpp_get_current_pelist(pes) 01051 ! pes = (/ (i,i=0,mpp_npes()-1) /) 01052 end if 01053 !rv Block distribution ? 01054 blocks=.false. 01055 if(present(offsetx).and.present(offsety)) blocks=.true. 01056 if(present(offsetx).and.(.not.present(offsety))) then 01057 call mpp_error(FATAL,'MPP_DEFINE_DOMAINS1D: you have to define offsetx AND offsety!') 01058 endif 01059 if((.not.present(offsetx)).and.present(offsety)) then 01060 call mpp_error(FATAL,'MPP_DEFINE_DOMAINS1D: you have to define offsetx AND offsety!') 01061 endif 01062 if(blocks) then 01063 if((.not.present(maskmap)).and.(.not.present(xextent))) & 01064 call mpp_error(FATAL, 'MPP_DEFINE_DOMAINS1D: you have to define maskmap and xextent!') 01065 if((.not.present(maskmap)).and.(.not.present(yextent))) & 01066 call mpp_error(FATAL, 'MPP_DEFINE_DOMAINS1D: you have to define maskmap and yextent!') 01067 if(size(offsetx).ne.ndivx) & 01068 call mpp_error(FATAL,'MPP_DEFINE_DOMAINS2D: offsetx array does not match layout!') 01069 if(size(offsety).ne.ndivy) & 01070 call mpp_error(FATAL,'MPP_DEFINE_DOMAINS2D: offsety array does not match layout!') 01071 endif 01072 !rv 01073 01074 allocate( mask(0:ndivx-1,0:ndivy-1) ) 01075 mask = .TRUE. 01076 if( PRESENT(maskmap) )then 01077 if( size(maskmap,1).NE.ndivx .OR. size(maskmap,2).NE.ndivy ) & 01078 call mpp_error( FATAL, 'MPP_DEFINE_DOMAINS2D: maskmap array does not match layout.' ) 01079 mask(:,:) = maskmap(:,:) 01080 end if 01081 !number of unmask domains in layout must equal number of PEs assigned 01082 n = count(mask) 01083 if( n.NE.size(pes) )then 01084 write( text,'(i8)' )n 01085 call mpp_error( FATAL, 'MPP_DEFINE_DOMAINS2D: incorrect number of PEs assigned for this layout and maskmap. Use ' & 01086 //text//' PEs for this domain decomposition.' ) 01087 end if 01088 01089 !place on PE array; need flag to assign them to j first and then i 01090 allocate( pearray(0:ndivx-1,0:ndivy-1) ) 01091 pearray(:,:) = NULL_PE 01092 ipos = NULL_PE; jpos = NULL_PE; pos = NULL_PE 01093 n = 0 01094 do j = 0,ndivy-1 01095 do i = 0,ndivx-1 01096 if( mask(i,j) )then 01097 pearray(i,j) = pes(n) 01098 if( pes(n).EQ.pe )then 01099 pos = n 01100 ipos = i 01101 jpos = j 01102 end if 01103 n = n + 1 01104 end if 01105 end do 01106 end do 01107 if( ipos.EQ.NULL_PE .OR. jpos.EQ.NULL_PE .or. pos.EQ.NULL_PE ) & 01108 call mpp_error( FATAL, 'MPP_DEFINE_DOMAINS2D: pelist must include this PE.' ) 01109 if( debug )write( stderr(), * )'pe, ipos, jpos=', pe, ipos, jpos, ' pearray(:,jpos)=', pearray(:,jpos), & 01110 ' pearray(ipos,:)=', pearray(ipos,:) 01111 01112 if(blocks) then 01113 !do domain decomposition using 1D versions in X and Y 01114 call mpp_define_domains( global_indices(1:2), ndivx, domain%x, & 01115 pack(pearray(:,jpos),mask(:,jpos)), xflags, xhalo, xextent, mask(:,jpos),offset=offsetx ) 01116 call mpp_define_domains( global_indices(3:4), ndivy, domain%y, & 01117 pack(pearray(ipos,:),mask(ipos,:)), yflags, yhalo, yextent, mask(ipos,:),offset=offsety ) 01118 else 01119 !do domain decomposition using 1D versions in X and Y 01120 call mpp_define_domains( global_indices(1:2), ndivx, domain%x, & 01121 pack(pearray(:,jpos),mask(:,jpos)), xflags, xhalo, xextent, mask(:,jpos)) 01122 call mpp_define_domains( global_indices(3:4), ndivy, domain%y, & 01123 pack(pearray(ipos,:),mask(ipos,:)), yflags, yhalo, yextent, mask(ipos,:)) 01124 endif 01125 if( domain%x%list(domain%x%pos)%pe.NE.domain%y%list(domain%y%pos)%pe ) & 01126 call mpp_error( FATAL, 'MPP_DEFINE_DOMAINS2D: domain%x%list(ipos)%pe.NE.domain%y%list(jpos)%pe.' ) 01127 domain%pos = pos 01128 domain%pe = pe 01129 01130 !set up fold 01131 domain%fold = 0 01132 if( PRESENT(xflags) )then 01133 if( BTEST(xflags,WEST) )domain%fold = domain%fold + FOLD_WEST_EDGE 01134 if( BTEST(xflags,EAST) )domain%fold = domain%fold + FOLD_EAST_EDGE 01135 end if 01136 if( PRESENT(yflags) )then 01137 if( BTEST(yflags,SOUTH) )domain%fold = domain%fold + FOLD_SOUTH_EDGE 01138 if( BTEST(yflags,NORTH) )domain%fold = domain%fold + FOLD_NORTH_EDGE 01139 end if 01140 01141 if( BTEST(domain%fold,SOUTH) .OR. BTEST(domain%fold,NORTH) )then 01142 if( domain%y%cyclic )call mpp_error( FATAL, 'MPP_DEFINE_DOMAINS: an axis cannot be both folded and cyclic.' ) 01143 if( modulo(domain%x%global%size,2).NE.0 ) & 01144 call mpp_error( FATAL, 'MPP_DEFINE_DOMAINS: number of points in X must be even when there is a fold in Y.' ) 01145 !check if folded domain boundaries line up in X: compute domains lining up is a sufficient condition for symmetry 01146 n = ndivx - 1 01147 do i = 0,n/2 01148 if( domain%x%list(i)%compute%size.NE.domain%x%list(n-i)%compute%size ) & 01149 call mpp_error( FATAL, 'MPP_DEFINE_DOMAINS: Folded domain boundaries must line up (mirror-symmetric extents).' ) 01150 end do 01151 end if 01152 if( BTEST(domain%fold,WEST) .OR. BTEST(domain%fold,EAST) )then 01153 if( domain%x%cyclic )call mpp_error( FATAL, 'MPP_DEFINE_DOMAINS: an axis cannot be both folded and cyclic.' ) 01154 if( modulo(domain%y%global%size,2).NE.0 ) & 01155 call mpp_error( FATAL, 'MPP_DEFINE_DOMAINS: number of points in Y must be even when there is a fold in X.' ) 01156 !check if folded domain boundaries line up in Y: compute domains lining up is a sufficient condition for symmetry 01157 n = ndivy - 1 01158 do i = 0,n/2 01159 if( domain%y%list(i)%compute%size.NE.domain%y%list(n-i)%compute%size ) & 01160 call mpp_error( FATAL, 'MPP_DEFINE_DOMAINS: Folded domain boundaries must line up (mirror-symmetric extents).' ) 01161 end do 01162 end if 01163 01164 !set up domain%list 01165 if( debug )write( stderr(),'(a,9i4)' )'pe, domain=', pe, domain_info_buf(1:8) 01166 if( pe.EQ.mpp_root_pe() .AND. PRESENT(name) )then 01167 write( nulprt, '(/a,i3,a,i3)' )trim(name)//' domain decomposition: ', ndivx, ' X', ndivy 01168 write( nulprt, '(3x,a)' )'pe, is, ie, js, je, isd, ied, jsd, jed' 01169 end if 01170 call mpp_sync() 01171 call mpp_get_compute_domain( domain, domain_info_buf(1), domain_info_buf(2), domain_info_buf(3), domain_info_buf(4) ) 01172 call mpp_get_data_domain ( domain, domain_info_buf(5), domain_info_buf(6), domain_info_buf(7), domain_info_buf(8) ) 01173 n = size(pes) 01174 allocate( domain%list(0:n-1) ) !this is only used for storage of remote compute and data domain info 01175 do i = 0,n-1 01176 m = mod(pos+i,n) 01177 domain%list(m)%pe = pes(m) 01178 call mpp_transmit( domain_info_buf(1:8), 8, pes(mod(pos+n-i,n)), domain_info_buf(9:16), 8, pes(m) ) 01179 domain%list(m)%x%compute%begin = domain_info_buf(9) 01180 domain%list(m)%x%compute%end = domain_info_buf(10) 01181 domain%list(m)%y%compute%begin = domain_info_buf(11) 01182 domain%list(m)%y%compute%end = domain_info_buf(12) 01183 domain%list(m)%x%data%begin = domain_info_buf(13) 01184 domain%list(m)%x%data%end = domain_info_buf(14) 01185 domain%list(m)%y%data%begin = domain_info_buf(15) 01186 domain%list(m)%y%data%end = domain_info_buf(16) 01187 if( pe.EQ.mpp_root_pe() .AND. PRESENT(name) )write( nulprt, '(2x,i3,1x,4i5,3x,4i5)' )pes(m), domain_info_buf(9:) 01188 end do 01189 call mpp_sync_self(pes) 01190 domain%list(:)%x%compute%size = domain%list(:)%x%compute%end - domain%list(:)%x%compute%begin + 1 01191 domain%list(:)%y%compute%size = domain%list(:)%y%compute%end - domain%list(:)%y%compute%begin + 1 01192 domain%list(:)%x%data%size = domain%list(:)%x%data%end - domain%list(:)%x%data%begin + 1 01193 domain%list(:)%y%data%size = domain%list(:)%y%data%end - domain%list(:)%y%data%begin + 1 01194 01195 domain%remote_domains_initialized = .FALSE. 01196 domain%remote_off_domains_initialized = .FALSE. 01197 call compute_overlaps(domain) 01198 !PV786667: the deallocate stmts can be removed when fixed (7.3.1.3m) 01199 deallocate( pes, mask, pearray ) 01200 01201 return 01202 end subroutine mpp_define_domains2D 01203 01204 subroutine mpp_broadcast_domain( domain ) 01205 !broadcast domain (useful only outside the context of its own pelist) 01206 type(domain2D), intent(inout) :: domain 01207 integer, allocatable :: pes(:) 01208 logical :: native !true if I'm on the pelist of this domain 01209 integer :: listsize, listpos 01210 integer :: n 01211 integer, dimension(5) :: msg, info !pe and compute domain of each item in list 01212 01213 if( .NOT.module_is_initialized )call mpp_error( FATAL, 'MPP_BROADCAST_DOMAIN: You must first call mpp_domains_init.' ) 01214 01215 !get the current pelist 01216 allocate( pes(0:mpp_npes()-1) ) 01217 call mpp_get_current_pelist(pes) 01218 01219 !am I part of this domain? 01220 native = ASSOCIATED(domain%list) 01221 01222 !set local list size 01223 if( native )then 01224 listsize = size(domain%list) 01225 else 01226 listsize = 0 01227 end if 01228 call mpp_max(listsize) 01229 01230 if( .NOT.native )then 01231 !initialize domain%list and set null values in message 01232 allocate( domain%list(0:listsize-1) ) 01233 domain%pe = NULL_PE 01234 domain%pos = -1 01235 domain%x%compute%begin = HUGE(1) 01236 domain%x%compute%end = -HUGE(1) 01237 domain%y%compute%begin = HUGE(1) 01238 domain%y%compute%end = -HUGE(1) 01239 end if 01240 !initialize values in info 01241 info(1) = domain%pe 01242 call mpp_get_compute_domain( domain, info(2), info(3), info(4), info(5) ) 01243 01244 !broadcast your info across current pelist and unpack if needed 01245 listpos = 0 01246 do n = 0,mpp_npes()-1 01247 msg = info 01248 if( pe.EQ.pes(n) .AND. debug )write( stderr(),* )'PE ', pe, 'broadcasting msg ', msg 01249 call mpp_broadcast( msg, 5, pes(n) ) 01250 !no need to unpack message if native 01251 !no need to unpack message from non-native PE 01252 if( .NOT.native .AND. msg(1).NE.NULL_PE )then 01253 domain%list(listpos)%pe = msg(1) 01254 domain%list(listpos)%x%compute%begin = msg(2) 01255 domain%list(listpos)%x%compute%end = msg(3) 01256 domain%list(listpos)%y%compute%begin = msg(4) 01257 domain%list(listpos)%y%compute%end = msg(5) 01258 listpos = listpos + 1 01259 if( debug )write( stderr(),* )'PE ', pe, 'received domain from PE ', msg(1), 'is,ie,js,je=', msg(2:5) 01260 end if 01261 end do 01262 end subroutine mpp_broadcast_domain 01263 01264 subroutine compute_overlaps( domain ) 01265 !computes remote domain overlaps 01266 !assumes only one in each direction 01267 type(domain2D), intent(inout) :: domain 01268 integer :: i, j, k, m, n 01269 integer :: is, ie, js, je, isc, iec, jsc, jec, isd, ied, jsd, jed, isg, ieg, jsg, jeg, ioff, joff 01270 integer :: list 01271 01272 if( grid_offset_type.EQ.AGRID .AND. domain%remote_domains_initialized )return 01273 if( grid_offset_type.NE.AGRID .AND. domain%remote_off_domains_initialized )return 01274 domain%gridtype = grid_offset_type 01275 n = size(domain%list) 01276 !send 01277 call mpp_get_compute_domain( domain, isc, iec, jsc, jec ) 01278 call mpp_get_global_domain ( domain, isg, ieg, jsg, jeg, xsize=ioff, ysize=joff ) !cyclic offsets 01279 domain%list(:)%overlap = .FALSE. 01280 do list = 0,n-1 01281 m = mod( domain%pos+list, n ) 01282 !to_pe's eastern halo 01283 is = domain%list(m)%x%compute%end+1; ie = domain%list(m)%x%data%end 01284 js = domain%list(m)%y%compute%begin; je = domain%list(m)%y%compute%end 01285 if( is.GT.ieg )then 01286 if( domain%x%cyclic )then !try cyclic offset 01287 is = is-ioff; ie = ie-ioff 01288 else if( BTEST(domain%fold,EAST) )then 01289 i=is; is = 2*ieg-ie+1; ie = 2*ieg-i+1 01290 j=js; js = jsg+jeg-je; je = jsg+jeg-j 01291 if( BTEST(grid_offset_type,EAST) )then 01292 is = is - 1; ie = ie - 1 01293 end if 01294 end if 01295 end if 01296 is = max(is,isc); ie = min(ie,iec) 01297 js = max(js,jsc); je = min(je,jec) 01298 if( ie.GE.is .AND. je.GE.js )then 01299 domain%list(m)%overlap = .TRUE. 01300 if( grid_offset_type.NE.AGRID )then 01301 domain%list(m)%send_w_off%overlap = .TRUE. 01302 domain%list(m)%send_w_off%is = is 01303 domain%list(m)%send_w_off%ie = ie 01304 domain%list(m)%send_w_off%js = js 01305 domain%list(m)%send_w_off%je = je 01306 else 01307 domain%list(m)%send_w%overlap = .TRUE. 01308 domain%list(m)%send_w%is = is 01309 domain%list(m)%send_w%ie = ie 01310 domain%list(m)%send_w%js = js 01311 domain%list(m)%send_w%je = je 01312 end if 01313 else 01314 domain%list(m)%send_w%overlap = .FALSE. 01315 domain%list(m)%send_w_off%overlap = .FALSE. 01316 end if 01317 !to_pe's SE halo 01318 is = domain%list(m)%x%compute%end+1; ie = domain%list(m)%x%data%end 01319 js = domain%list(m)%y%data%begin; je = domain%list(m)%y%compute%begin-1 01320 if( is.GT.ieg )then 01321 if( domain%x%cyclic )then !try cyclic offset 01322 is = is-ioff; ie = ie-ioff 01323 else if( BTEST(domain%fold,EAST) )then 01324 i=is; is = 2*ieg-ie+1; ie = 2*ieg-i+1 01325 j=js; js = jsg+jeg-je; je = jsg+jeg-j 01326 if( BTEST(grid_offset_type,EAST) )then 01327 is = is - 1; ie = ie - 1 01328 end if 01329 end if 01330 end if 01331 if( jsg.GT.je )then 01332 if( domain%y%cyclic )then !try cyclic offset 01333 js = js+joff; je = je+joff 01334 else if( BTEST(domain%fold,SOUTH) )then 01335 i=is; is = isg+ieg-ie; ie = isg+ieg-i 01336 j=js; js = 2*jsg-je-1; je = 2*jsg-j-1 01337 if( BTEST(grid_offset_type,SOUTH) )then 01338 js = js + 1; je = je + 1 01339 end if 01340 end if 01341 end if 01342 is = max(is,isc); ie = min(ie,iec) 01343 js = max(js,jsc); je = min(je,jec) 01344 if( ie.GE.is .AND. je.GE.js )then 01345 domain%list(m)%overlap = .TRUE. 01346 if( grid_offset_type.NE.AGRID )then 01347 domain%list(m)%send_nw_off%overlap = .TRUE. 01348 domain%list(m)%send_nw_off%is = is 01349 domain%list(m)%send_nw_off%ie = ie 01350 domain%list(m)%send_nw_off%js = js 01351 domain%list(m)%send_nw_off%je = je 01352 else 01353 domain%list(m)%send_nw%overlap = .TRUE. 01354 domain%list(m)%send_nw%is = is 01355 domain%list(m)%send_nw%ie = ie 01356 domain%list(m)%send_nw%js = js 01357 domain%list(m)%send_nw%je = je 01358 end if 01359 else 01360 domain%list(m)%send_nw%overlap = .FALSE. 01361 domain%list(m)%send_nw_off%overlap = .FALSE. 01362 end if 01363 !to_pe's southern halo 01364 is = domain%list(m)%x%compute%begin; ie = domain%list(m)%x%compute%end 01365 js = domain%list(m)%y%data%begin; je = domain%list(m)%y%compute%begin-1 01366 if( jsg.GT.je )then 01367 if( domain%y%cyclic )then !try cyclic offset 01368 js = js+joff; je = je+joff 01369 else if( BTEST(domain%fold,SOUTH) )then 01370 i=is; is = isg+ieg-ie; ie = isg+ieg-i 01371 j=js; js = 2*jsg-je-1; je = 2*jsg-j-1 01372 if( BTEST(grid_offset_type,SOUTH) )then 01373 js = js + 1; je = je + 1 01374 end if 01375 end if 01376 end if 01377 is = max(is,isc); ie = min(ie,iec) 01378 js = max(js,jsc); je = min(je,jec) 01379 if( ie.GE.is .AND. je.GE.js )then 01380 domain%list(m)%overlap = .TRUE. 01381 if( grid_offset_type.NE.AGRID )then 01382 domain%list(m)%send_n_off%overlap = .TRUE. 01383 domain%list(m)%send_n_off%is = is 01384 domain%list(m)%send_n_off%ie = ie 01385 domain%list(m)%send_n_off%js = js 01386 domain%list(m)%send_n_off%je = je 01387 else 01388 domain%list(m)%send_n%overlap = .TRUE. 01389 domain%list(m)%send_n%is = is 01390 domain%list(m)%send_n%ie = ie 01391 domain%list(m)%send_n%js = js 01392 domain%list(m)%send_n%je = je 01393 end if 01394 else 01395 domain%list(m)%send_n%overlap = .FALSE. 01396 domain%list(m)%send_n_off%overlap = .FALSE. 01397 end if 01398 !to_pe's SW halo 01399 is = domain%list(m)%x%data%begin; ie = domain%list(m)%x%compute%begin-1 01400 js = domain%list(m)%y%data%begin; je = domain%list(m)%y%compute%begin-1 01401 if( isg.GT.ie )then 01402 if( domain%x%cyclic )then !try cyclic offset 01403 is = is+ioff; ie = ie+ioff 01404 else if( BTEST(domain%fold,WEST) )then 01405 i=is; is = 2*isg-ie-1; ie = 2*isg-i-1 01406 j=js; js = jsg+jeg-je; je = jsg+jeg-j 01407 if( BTEST(grid_offset_type,WEST) )then 01408 is = is + 1; ie = ie + 1 01409 end if 01410 end if 01411 end if 01412 if( jsg.GT.je )then 01413 if( domain%y%cyclic )then !try cyclic offset 01414 js = js+joff; je = je+joff 01415 else if( BTEST(domain%fold,SOUTH) )then 01416 i=is; is = isg+ieg-ie; ie = isg+ieg-i 01417 j=js; js = 2*jsg-je-1; je = 2*jsg-j-1 01418 if( BTEST(grid_offset_type,SOUTH) )then 01419 js = js + 1; je = je + 1 01420 end if 01421 end if 01422 end if 01423 is = max(is,isc); ie = min(ie,iec) 01424 js = max(js,jsc); je = min(je,jec) 01425 if( ie.GE.is .AND. je.GE.js )then 01426 domain%list(m)%overlap = .TRUE. 01427 if( grid_offset_type.NE.AGRID )then 01428 domain%list(m)%send_ne_off%overlap = .TRUE. 01429 domain%list(m)%send_ne_off%is = is 01430 domain%list(m)%send_ne_off%ie = ie 01431 domain%list(m)%send_ne_off%js = js 01432 domain%list(m)%send_ne_off%je = je 01433 else 01434 domain%list(m)%send_ne%overlap = .TRUE. 01435 domain%list(m)%send_ne%is = is 01436 domain%list(m)%send_ne%ie = ie 01437 domain%list(m)%send_ne%js = js 01438 domain%list(m)%send_ne%je = je 01439 end if 01440 else 01441 domain%list(m)%send_ne%overlap = .FALSE. 01442 domain%list(m)%send_ne_off%overlap = .FALSE. 01443 end if 01444 !to_pe's western halo 01445 is = domain%list(m)%x%data%begin; ie = domain%list(m)%x%compute%begin-1 01446 js = domain%list(m)%y%compute%begin; je = domain%list(m)%y%compute%end 01447 if( isg.GT.ie )then 01448 if( domain%x%cyclic )then !try cyclic offset 01449 is = is+ioff; ie = ie+ioff 01450 else if( BTEST(domain%fold,WEST) )then 01451 i=is; is = 2*isg-ie-1; ie = 2*isg-i-1 01452 j=js; js = jsg+jeg-je; je = jsg+jeg-j 01453 if( BTEST(grid_offset_type,WEST) )then 01454 is = is + 1; ie = ie + 1 01455 end if 01456 end if 01457 end if 01458 is = max(is,isc); ie = min(ie,iec) 01459 js = max(js,jsc); je = min(je,jec) 01460 if( ie.GE.is .AND. je.GE.js )then 01461 domain%list(m)%overlap = .TRUE. 01462 if( grid_offset_type.NE.AGRID )then 01463 domain%list(m)%send_e_off%overlap = .TRUE. 01464 domain%list(m)%send_e_off%is = is 01465 domain%list(m)%send_e_off%ie = ie 01466 domain%list(m)%send_e_off%js = js 01467 domain%list(m)%send_e_off%je = je 01468 else 01469 domain%list(m)%send_e%overlap = .TRUE. 01470 domain%list(m)%send_e%is = is 01471 domain%list(m)%send_e%ie = ie 01472 domain%list(m)%send_e%js = js 01473 domain%list(m)%send_e%je = je 01474 end if 01475 else 01476 domain%list(m)%send_e%overlap = .FALSE. 01477 domain%list(m)%send_e_off%overlap = .FALSE. 01478 end if 01479 !to_pe's NW halo 01480 is = domain%list(m)%x%data%begin; ie = domain%list(m)%x%compute%begin-1 01481 js = domain%list(m)%y%compute%end+1; je = domain%list(m)%y%data%end 01482 if( isg.GT.ie )then 01483 if( domain%x%cyclic )then !try cyclic offset 01484 is = is+ioff; ie = ie+ioff 01485 else if( BTEST(domain%fold,WEST) )then 01486 i=is; is = 2*isg-ie-1; ie = 2*isg-i-1 01487 j=js; js = jsg+jeg-je; je = jsg+jeg-j 01488 if( BTEST(grid_offset_type,WEST) )then 01489 is = is + 1; ie = ie + 1 01490 end if 01491 end if 01492 end if 01493 if( js.GT.jeg )then 01494 if( domain%y%cyclic )then !try cyclic offset 01495 js = js-joff; je = je-joff 01496 else if( BTEST(domain%fold,NORTH) )then 01497 i=is; is = isg+ieg-ie; ie = isg+ieg-i 01498 j=js; js = 2*jeg-je+1; je = 2*jeg-j+1 01499 if( BTEST(grid_offset_type,NORTH) )then 01500 js = js - 1; je = je - 1 01501 end if 01502 end if 01503 end if 01504 is = max(is,isc); ie = min(ie,iec) 01505 js = max(js,jsc); je = min(je,jec) 01506 if( ie.GE.is .AND. je.GE.js )then 01507 domain%list(m)%overlap = .TRUE. 01508 if( grid_offset_type.NE.AGRID )then 01509 domain%list(m)%send_se_off%overlap = .TRUE. 01510 domain%list(m)%send_se_off%is = is 01511 domain%list(m)%send_se_off%ie = ie 01512 domain%list(m)%send_se_off%js = js 01513 domain%list(m)%send_se_off%je = je 01514 else 01515 domain%list(m)%send_se%overlap = .TRUE. 01516 domain%list(m)%send_se%is = is 01517 domain%list(m)%send_se%ie = ie 01518 domain%list(m)%send_se%js = js 01519 domain%list(m)%send_se%je = je 01520 end if 01521 else 01522 domain%list(m)%send_se%overlap = .FALSE. 01523 domain%list(m)%send_se_off%overlap = .FALSE. 01524 end if 01525 !to_pe's northern halo 01526 is = domain%list(m)%x%compute%begin; ie = domain%list(m)%x%compute%end 01527 js = domain%list(m)%y%compute%end+1; je = domain%list(m)%y%data%end 01528 if( js.GT.jeg )then 01529 if( domain%y%cyclic )then !try cyclic offset 01530 js = js-joff; je = je-joff 01531 else if( BTEST(domain%fold,NORTH) )then 01532 i=is; is = isg+ieg-ie; ie = isg+ieg-i 01533 j=js; js = 2*jeg-je+1; je = 2*jeg-j+1 01534 if( BTEST(grid_offset_type,NORTH) )then 01535 js = js - 1; je = je - 1 01536 end if 01537 end if 01538 end if 01539 is = max(is,isc); ie = min(ie,iec) 01540 js = max(js,jsc); je = min(je,jec) 01541 if( ie.GE.is .AND. je.GE.js )then 01542 domain%list(m)%overlap = .TRUE. 01543 if( grid_offset_type.NE.AGRID )then 01544 domain%list(m)%send_s_off%overlap = .TRUE. 01545 domain%list(m)%send_s_off%is = is 01546 domain%list(m)%send_s_off%ie = ie 01547 domain%list(m)%send_s_off%js = js 01548 domain%list(m)%send_s_off%je = je 01549 else 01550 domain%list(m)%send_s%overlap = .TRUE. 01551 domain%list(m)%send_s%is = is 01552 domain%list(m)%send_s%ie = ie 01553 domain%list(m)%send_s%js = js 01554 domain%list(m)%send_s%je = je 01555 end if 01556 else 01557 domain%list(m)%send_s%overlap = .FALSE. 01558 domain%list(m)%send_s_off%overlap = .FALSE. 01559 end if 01560 !to_pe's NE halo 01561 is = domain%list(m)%x%compute%end+1; ie = domain%list(m)%x%data%end 01562 js = domain%list(m)%y%compute%end+1; je = domain%list(m)%y%data%end 01563 if( is.GT.ieg )then 01564 if( domain%x%cyclic )then !try cyclic offset 01565 is = is-ioff; ie = ie-ioff 01566 else if( BTEST(domain%fold,EAST) )then 01567 i=is; is = 2*ieg-ie+1; ie = 2*ieg-i+1 01568 j=js; js = jsg+jeg-je; je = jsg+jeg-j 01569 end if 01570 end if 01571 if( js.GT.jeg )then 01572 if( domain%y%cyclic )then !try cyclic offset 01573 js = js-joff; je = je-joff 01574 else if( BTEST(domain%fold,NORTH) )then 01575 i=is; is = isg+ieg-ie; ie = isg+ieg-i 01576 j=js; js = 2*jeg-je+1; je = 2*jeg-j+1 01577 if( BTEST(grid_offset_type,NORTH) )then 01578 js = js - 1; je = je - 1 01579 end if 01580 end if 01581 end if 01582 is = max(is,isc); ie = min(ie,iec) 01583 js = max(js,jsc); je = min(je,jec) 01584 if( ie.GE.is .AND. je.GE.js )then 01585 domain%list(m)%overlap = .TRUE. 01586 if( grid_offset_type.NE.AGRID )then 01587 domain%list(m)%send_sw_off%overlap = .TRUE. 01588 domain%list(m)%send_sw_off%is = is 01589 domain%list(m)%send_sw_off%ie = ie 01590 domain%list(m)%send_sw_off%js = js 01591 domain%list(m)%send_sw_off%je = je 01592 else 01593 domain%list(m)%send_sw%overlap = .TRUE. 01594 domain%list(m)%send_sw%is = is 01595 domain%list(m)%send_sw%ie = ie 01596 domain%list(m)%send_sw%js = js 01597 domain%list(m)%send_sw%je = je 01598 end if 01599 else 01600 domain%list(m)%send_sw%overlap = .FALSE. 01601 domain%list(m)%send_sw_off%overlap = .FALSE. 01602 end if 01603 end do 01604 01605 !recv 01606 do list = 0,n-1 01607 m = mod( domain%pos+n-list, n ) 01608 call mpp_get_compute_domain( domain%list(m), isc, iec, jsc, jec ) 01609 !recv_e 01610 isd = domain%x%compute%end+1; ied = domain%x%data%end 01611 jsd = domain%y%compute%begin; jed = domain%y%compute%end 01612 is=isc; ie=iec; js=jsc; je=jec 01613 domain%list(m)%recv_e%folded = .FALSE. 01614 if( isd.GT.ieg )then 01615 if( domain%x%cyclic )then !try cyclic offset 01616 is = is+ioff; ie = ie+ioff 01617 else if( BTEST(domain%fold,EAST) )then 01618 domain%list(m)%recv_e%folded = .TRUE. 01619 i=is; is = 2*ieg-ie+1; ie = 2*ieg-i+1 01620 j=js; js = jsg+jeg-je; je = jsg+jeg-j 01621 if( BTEST(grid_offset_type,EAST) )then 01622 is = is - 1; ie = ie - 1 01623 end if 01624 end if 01625 end if 01626 is = max(isd,is); ie = min(ied,ie) 01627 js = max(jsd,js); je = min(jed,je) 01628 if( ie.GE.is .AND. je.GE.js )then 01629 domain%list(m)%overlap = .TRUE. 01630 if( grid_offset_type.NE.AGRID )then 01631 domain%list(m)%recv_e_off%overlap = .TRUE. 01632 domain%list(m)%recv_e_off%is = is 01633 domain%list(m)%recv_e_off%ie = ie 01634 domain%list(m)%recv_e_off%js = js 01635 domain%list(m)%recv_e_off%je = je 01636 else 01637 domain%list(m)%recv_e%overlap = .TRUE. 01638 domain%list(m)%recv_e%is = is 01639 domain%list(m)%recv_e%ie = ie 01640 domain%list(m)%recv_e%js = js 01641 domain%list(m)%recv_e%je = je 01642 endif 01643 else 01644 domain%list(m)%recv_e%overlap = .FALSE. 01645 domain%list(m)%recv_e_off%overlap = .FALSE. 01646 end if 01647 !recv_se 01648 isd = domain%x%compute%end+1; ied = domain%x%data%end 01649 jsd = domain%y%data%begin; jed = domain%y%compute%begin-1 01650 is=isc; ie=iec; js=jsc; je=jec 01651 domain%list(m)%recv_se%folded = .FALSE. 01652 if( jed.LT.jsg )then 01653 if( domain%y%cyclic )then !try cyclic offset 01654 js = js-joff; je = je-joff 01655 else if( BTEST(domain%fold,SOUTH) )then 01656 i=is; is = isg+ieg-ie; ie = isg+ieg-i 01657 j=js; js = 2*jsg-je-1; je = 2*jsg-j-1 01658 domain%list(m)%recv_se%folded = .TRUE. 01659 if( BTEST(grid_offset_type,SOUTH) )then 01660 js = js + 1; je = je + 1 01661 end if 01662 end if 01663 end if 01664 if( isd.GT.ieg )then 01665 if( domain%x%cyclic )then !try cyclic offset 01666 is = is+ioff; ie = ie+ioff 01667 else if( BTEST(domain%fold,EAST) )then 01668 i=is; is = 2*ieg-ie+1; ie = 2*ieg-i+1 01669 j=js; js = jsg+jeg-je; je = jsg+jeg-j 01670 domain%list(m)%recv_se%folded = .TRUE. 01671 if( BTEST(grid_offset_type,EAST) )then 01672 is = is - 1; ie = ie - 1 01673 end if 01674 end if 01675 end if 01676 is = max(isd,is); ie = min(ied,ie) 01677 js = max(jsd,js); je = min(jed,je) 01678 if( ie.GE.is .AND. je.GE.js )then 01679 domain%list(m)%overlap = .TRUE. 01680 if( grid_offset_type.NE.AGRID )then 01681 domain%list(m)%recv_se_off%overlap = .TRUE. 01682 domain%list(m)%recv_se_off%is = is 01683 domain%list(m)%recv_se_off%ie = ie 01684 domain%list(m)%recv_se_off%js = js 01685 domain%list(m)%recv_se_off%je = je 01686 else 01687 domain%list(m)%recv_se%overlap = .TRUE. 01688 domain%list(m)%recv_se%is = is 01689 domain%list(m)%recv_se%ie = ie 01690 domain%list(m)%recv_se%js = js 01691 domain%list(m)%recv_se%je = je 01692 endif 01693 else 01694 domain%list(m)%recv_se%overlap = .FALSE. 01695 domain%list(m)%recv_se_off%overlap = .FALSE. 01696 end if 01697 !recv_s 01698 isd = domain%x%compute%begin; ied = domain%x%compute%end 01699 jsd = domain%y%data%begin; jed = domain%y%compute%begin-1 01700 is=isc; ie=iec; js=jsc; je=jec 01701 domain%list(m)%recv_s%folded = .FALSE. 01702 if( jed.LT.jsg )then 01703 if( domain%y%cyclic )then !try cyclic offset 01704 js = js-joff; je = je-joff 01705 else if( BTEST(domain%fold,SOUTH) )then 01706 i=is; is = isg+ieg-ie; ie = isg+ieg-i 01707 j=js; js = 2*jsg-je-1; je = 2*jsg-j-1 01708 domain%list(m)%recv_s%folded = .TRUE. 01709 if( BTEST(grid_offset_type,SOUTH) )then 01710 js = js + 1; je = je + 1 01711 end if 01712 end if 01713 end if 01714 is = max(isd,is); ie = min(ied,ie) 01715 js = max(jsd,js); je = min(jed,je) 01716 if( ie.GE.is .AND. je.GE.js )then 01717 domain%list(m)%overlap = .TRUE. 01718 if( grid_offset_type.NE.AGRID )then 01719 domain%list(m)%recv_s_off%overlap = .TRUE. 01720 domain%list(m)%recv_s_off%is = is 01721 domain%list(m)%recv_s_off%ie = ie 01722 domain%list(m)%recv_s_off%js = js 01723 domain%list(m)%recv_s_off%je = je 01724 else 01725 domain%list(m)%recv_s%overlap = .TRUE. 01726 domain%list(m)%recv_s%is = is 01727 domain%list(m)%recv_s%ie = ie 01728 domain%list(m)%recv_s%js = js 01729 domain%list(m)%recv_s%je = je 01730 endif 01731 else 01732 domain%list(m)%recv_s%overlap = .FALSE. 01733 domain%list(m)%recv_s_off%overlap = .FALSE. 01734 01735 end if 01736 !recv_sw 01737 isd = domain%x%data%begin; ied = domain%x%compute%begin-1 01738 jsd = domain%y%data%begin; jed = domain%y%compute%begin-1 01739 is=isc; ie=iec; js=jsc; je=jec 01740 domain%list(m)%recv_sw%folded = .FALSE. 01741 if( jed.LT.jsg )then 01742 if( domain%y%cyclic )then !try cyclic offset 01743 js = js-joff; je = je-joff 01744 else if( BTEST(domain%fold,SOUTH) )then 01745 i=is; is = isg+ieg-ie; ie = isg+ieg-i 01746 j=js; js = 2*jsg-je-1; je = 2*jsg-j-1 01747 domain%list(m)%recv_sw%folded = .TRUE. 01748 if( BTEST(grid_offset_type,SOUTH) )then 01749 js = js + 1; je = je + 1 01750 end if 01751 end if 01752 end if 01753 if( ied.LT.isg )then 01754 if( domain%x%cyclic )then !try cyclic offset 01755 is = is-ioff; ie = ie-ioff 01756 else if( BTEST(domain%fold,WEST) )then 01757 i=is; is = 2*isg-ie-1; ie = 2*isg-i-1 01758 j=js; js = jsg+jeg-je; je = jsg+jeg-j 01759 domain%list(m)%recv_sw%folded = .TRUE. 01760 if( BTEST(grid_offset_type,WEST) )then 01761 is = is + 1; ie = ie + 1 01762 end if 01763 end if 01764 end if 01765 is = max(isd,is); ie = min(ied,ie) 01766 js = max(jsd,js); je = min(jed,je) 01767 if( ie.GE.is .AND. je.GE.js )then 01768 domain%list(m)%overlap = .TRUE. 01769 if( grid_offset_type.NE.AGRID )then 01770 domain%list(m)%recv_sw_off%overlap = .TRUE. 01771 domain%list(m)%recv_sw_off%is = is 01772 domain%list(m)%recv_sw_off%ie = ie 01773 domain%list(m)%recv_sw_off%js = js 01774 domain%list(m)%recv_sw_off%je = je 01775 else 01776 domain%list(m)%recv_sw%overlap = .TRUE. 01777 domain%list(m)%recv_sw%is = is 01778 domain%list(m)%recv_sw%ie = ie 01779 domain%list(m)%recv_sw%js = js 01780 domain%list(m)%recv_sw%je = je 01781 endif 01782 else 01783 domain%list(m)%recv_sw%overlap = .FALSE. 01784 domain%list(m)%recv_sw_off%overlap = .FALSE. 01785 end if 01786 !recv_w 01787 isd = domain%x%data%begin; ied = domain%x%compute%begin-1 01788 jsd = domain%y%compute%begin; jed = domain%y%compute%end 01789 is=isc; ie=iec; js=jsc; je=jec 01790 domain%list(m)%recv_w%folded = .FALSE. 01791 if( ied.LT.isg )then 01792 if( domain%x%cyclic )then !try cyclic offset 01793 is = is-ioff; ie = ie-ioff 01794 else if( BTEST(domain%fold,WEST) )then 01795 i=is; is = 2*isg-ie-1; ie = 2*isg-i-1 01796 j=js; js = jsg+jeg-je; je = jsg+jeg-j 01797 domain%list(m)%recv_w%folded = .TRUE. 01798 if( BTEST(grid_offset_type,WEST) )then 01799 is = is + 1; ie = ie + 1 01800 end if 01801 end if 01802 end if 01803 is = max(isd,is); ie = min(ied,ie) 01804 js = max(jsd,js); je = min(jed,je) 01805 if( ie.GE.is .AND. je.GE.js )then 01806 domain%list(m)%overlap = .TRUE. 01807 if( grid_offset_type.NE.AGRID )then 01808 domain%list(m)%recv_w_off%overlap = .TRUE. 01809 domain%list(m)%recv_w_off%is = is 01810 domain%list(m)%recv_w_off%ie = ie 01811 domain%list(m)%recv_w_off%js = js 01812 domain%list(m)%recv_w_off%je = je 01813 else 01814 domain%list(m)%recv_w%overlap = .TRUE. 01815 domain%list(m)%recv_w%is = is 01816 domain%list(m)%recv_w%ie = ie 01817 domain%list(m)%recv_w%js = js 01818 domain%list(m)%recv_w%je = je 01819 endif 01820 else 01821 domain%list(m)%recv_w%overlap = .FALSE. 01822 domain%list(m)%recv_w_off%overlap = .FALSE. 01823 end if 01824 !recv_nw 01825 isd = domain%x%data%begin; ied = domain%x%compute%begin-1 01826 jsd = domain%y%compute%end+1; jed = domain%y%data%end 01827 is=isc; ie=iec; js=jsc; je=jec 01828 domain%list(m)%recv_nw%folded = .FALSE. 01829 if( jsd.GT.jeg )then 01830 if( domain%y%cyclic )then !try cyclic offset 01831 js = js+joff; je = je+joff 01832 else if( BTEST(domain%fold,NORTH) )then 01833 i=is; is = isg+ieg-ie; ie = isg+ieg-i 01834 j=js; js = 2*jeg-je+1; je = 2*jeg-j+1 01835 domain%list(m)%recv_nw%folded = .TRUE. 01836 if( BTEST(grid_offset_type,NORTH) )then 01837 js = js - 1; je = je - 1 01838 end if 01839 end if 01840 end if 01841 if( ied.LT.isg )then 01842 if( domain%x%cyclic )then !try cyclic offset 01843 is = is-ioff; ie = ie-ioff 01844 else if( BTEST(domain%fold,WEST) )then 01845 i=is; is = 2*isg-ie-1; ie = 2*isg-i-1 01846 j=js; js = jsg+jeg-je; je = jsg+jeg-j 01847 domain%list(m)%recv_nw%folded = .TRUE. 01848 if( BTEST(grid_offset_type,WEST) )then 01849 is = is + 1; ie = ie + 1 01850 end if 01851 end if 01852 end if 01853 is = max(isd,is); ie = min(ied,ie) 01854 js = max(jsd,js); je = min(jed,je) 01855 if( ie.GE.is .AND. je.GE.js )then 01856 domain%list(m)%overlap = .TRUE. 01857 if( grid_offset_type.NE.AGRID )then 01858 domain%list(m)%recv_nw_off%overlap = .TRUE. 01859 domain%list(m)%recv_nw_off%is = is 01860 domain%list(m)%recv_nw_off%ie = ie 01861 domain%list(m)%recv_nw_off%js = js 01862 domain%list(m)%recv_nw_off%je = je 01863 else 01864 domain%list(m)%recv_nw%overlap = .TRUE. 01865 domain%list(m)%recv_nw%is = is 01866 domain%list(m)%recv_nw%ie = ie 01867 domain%list(m)%recv_nw%js = js 01868 domain%list(m)%recv_nw%je = je 01869 endif 01870 else 01871 domain%list(m)%recv_nw%overlap = .FALSE. 01872 domain%list(m)%recv_nw_off%overlap = .FALSE. 01873 end if 01874 !recv_n 01875 isd = domain%x%compute%begin; ied = domain%x%compute%end 01876 jsd = domain%y%compute%end+1; jed = domain%y%data%end 01877 is=isc; ie=iec; js=jsc; je=jec 01878 domain%list(m)%recv_n%folded = .FALSE. 01879 if( jsd.GT.jeg )then 01880 if( domain%y%cyclic )then !try cyclic offset 01881 js = js+joff; je = je+joff 01882 else if( BTEST(domain%fold,NORTH) )then 01883 i=is; is = isg+ieg-ie; ie = isg+ieg-i 01884 j=js; js = 2*jeg-je+1; je = 2*jeg-j+1 01885 domain%list(m)%recv_n%folded = .TRUE. 01886 if( BTEST(grid_offset_type,NORTH) )then 01887 js = js - 1; je = je - 1 01888 end if 01889 end if 01890 end if 01891 is = max(isd,is); ie = min(ied,ie) 01892 js = max(jsd,js); je = min(jed,je) 01893 if( ie.GE.is .AND. je.GE.js )then 01894 domain%list(m)%overlap = .TRUE. 01895 if( grid_offset_type.NE.AGRID )then 01896 domain%list(m)%recv_n_off%overlap = .TRUE. 01897 domain%list(m)%recv_n_off%is = is 01898 domain%list(m)%recv_n_off%ie = ie 01899 domain%list(m)%recv_n_off%js = js 01900 domain%list(m)%recv_n_off%je = je 01901 else 01902 domain%list(m)%recv_n%overlap = .TRUE. 01903 domain%list(m)%recv_n%is = is 01904 domain%list(m)%recv_n%ie = ie 01905 domain%list(m)%recv_n%js = js 01906 domain%list(m)%recv_n%je = je 01907 end if 01908 else 01909 domain%list(m)%recv_n%overlap = .FALSE. 01910 domain%list(m)%recv_n_off%overlap = .FALSE. 01911 end if 01912 !recv_ne 01913 isd = domain%x%compute%end+1; ied = domain%x%data%end 01914 jsd = domain%y%compute%end+1; jed = domain%y%data%end 01915 is=isc; ie=iec; js=jsc; je=jec 01916 domain%list(m)%recv_ne%folded = .FALSE. 01917 if( jsd.GT.jeg )then 01918 if( domain%y%cyclic )then !try cyclic offset 01919 js = js+joff; je = je+joff 01920 else if( BTEST(domain%fold,NORTH) )then 01921 i=is; is = isg+ieg-ie; ie = isg+ieg-i 01922 j=js; js = 2*jeg-je+1; je = 2*jeg-j+1 01923 domain%list(m)%recv_ne%folded = .TRUE. 01924 if( BTEST(grid_offset_type,NORTH) )then 01925 js = js - 1; je = je - 1 01926 end if 01927 end if 01928 end if 01929 if( isd.GT.ieg )then 01930 if( domain%x%cyclic )then !try cyclic offset 01931 is = is+ioff; ie = ie+ioff 01932 else if( BTEST(domain%fold,EAST) )then 01933 i=is; is = 2*ieg-ie+1; ie = 2*ieg-i+1 01934 j=js; js = jsg+jeg-je; je = jsg+jeg-j 01935 domain%list(m)%recv_ne%folded = .TRUE. 01936 if( BTEST(grid_offset_type,EAST) )then 01937 is = is - 1; ie = ie - 1 01938 end if 01939 end if 01940 end if 01941 is = max(isd,is); ie = min(ied,ie) 01942 js = max(jsd,js); je = min(jed,je) 01943 if( ie.GE.is .AND. je.GE.js )then 01944 domain%list(m)%overlap = .TRUE. 01945 if( grid_offset_type.NE.AGRID )then 01946 domain%list(m)%recv_ne_off%overlap = .TRUE. 01947 domain%list(m)%recv_ne_off%is = is 01948 domain%list(m)%recv_ne_off%ie = ie 01949 domain%list(m)%recv_ne_off%js = js 01950 domain%list(m)%recv_ne_off%je = je 01951 else 01952 domain%list(m)%recv_ne%overlap = .TRUE. 01953 domain%list(m)%recv_ne%is = is 01954 domain%list(m)%recv_ne%ie = ie 01955 domain%list(m)%recv_ne%js = js 01956 domain%list(m)%recv_ne%je = je 01957 end if 01958 else 01959 domain%list(m)%recv_ne%overlap = .FALSE. 01960 domain%list(m)%recv_ne_off%overlap = .FALSE. 01961 end if 01962 end do 01963 if( grid_offset_type.EQ.AGRID )domain%remote_domains_initialized = .TRUE. 01964 if( grid_offset_type.NE.AGRID )domain%remote_off_domains_initialized = .TRUE. 01965 return 01966 end subroutine compute_overlaps 01967 01968 subroutine mpp_define_layout2D( global_indices, ndivs, layout ) 01969 integer, intent(in) :: global_indices(4) !(/ isg, ieg, jsg, jeg /) 01970 integer, intent(in) :: ndivs !number of divisions to divide global domain 01971 integer, intent(out) :: layout(2) 01972 01973 integer :: isg, ieg, jsg, jeg, isz, jsz, idiv, jdiv 01974 01975 isg = global_indices(1) 01976 ieg = global_indices(2) 01977 jsg = global_indices(3) 01978 jeg = global_indices(4) 01979 01980 isz = ieg - isg + 1 01981 jsz = jeg - jsg + 1 01982 !first try to divide ndivs in the domain aspect ratio: if imperfect aspect, reduce idiv till it divides ndivs 01983 idiv = nint( sqrt(float(ndivs*isz)/jsz) ) 01984 idiv = max(idiv,1) !for isz=1 line above can give 0 01985 do while( mod(ndivs,idiv).NE.0 ) 01986 idiv = idiv - 1 01987 end do !will terminate at idiv=1 if not before 01988 jdiv = ndivs/idiv 01989 01990 layout = (/ idiv, jdiv /) 01991 return 01992 end subroutine mpp_define_layout2D 01993 01994 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 01995 ! ! 01996 ! MPP_GET and SET routiness: retrieve various components of domains ! 01997 ! ! 01998 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 01999 02000 subroutine mpp_get_compute_domain1D( domain, begin, end, size, max_size, is_global ) 02001 type(domain1D), intent(in) :: domain 02002 integer, intent(out), optional :: begin, end, size, max_size 02003 logical, intent(out), optional :: is_global 02004 02005 if( PRESENT(begin) )begin = domain%compute%begin 02006 if( PRESENT(end) )end = domain%compute%end 02007 if( PRESENT(size) )size = domain%compute%size 02008 if( PRESENT(max_size) )max_size = domain%compute%max_size 02009 if( PRESENT(is_global) )is_global = domain%compute%is_global 02010 return 02011 end subroutine mpp_get_compute_domain1D 02012 02013 subroutine mpp_get_data_domain1D( domain, begin, end, size, max_size, is_global ) 02014 type(domain1D), intent(in) :: domain 02015 integer, intent(out), optional :: begin, end, size, max_size 02016 logical, intent(out), optional :: is_global 02017 02018 if( PRESENT(begin) )begin = domain%data%begin 02019 if( PRESENT(end) )end = domain%data%end 02020 if( PRESENT(size) )size = domain%data%size 02021 if( PRESENT(max_size) )max_size = domain%data%max_size 02022 if( PRESENT(is_global) )is_global = domain%data%is_global 02023 return 02024 end subroutine mpp_get_data_domain1D 02025 02026 subroutine mpp_get_global_domain1D( domain, begin, end, size, max_size ) 02027 type(domain1D), intent(in) :: domain 02028 integer, intent(out), optional :: begin, end, size, max_size 02029 02030 if( PRESENT(begin) )begin = domain%global%begin 02031 if( PRESENT(end) )end = domain%global%end 02032 if( PRESENT(size) )size = domain%global%size 02033 if( PRESENT(max_size) )max_size = domain%global%max_size 02034 return 02035 end subroutine mpp_get_global_domain1D 02036 02037 subroutine mpp_get_compute_domain2D( domain, xbegin, xend, ybegin, yend, xsize, xmax_size, ysize, ymax_size, & 02038 x_is_global, y_is_global ) 02039 type(domain2D), intent(in) :: domain 02040 integer, intent(out), optional :: xbegin, xend, ybegin, yend, xsize, xmax_size, ysize, ymax_size 02041 logical, intent(out), optional :: x_is_global, y_is_global 02042 call mpp_get_compute_domain( domain%x, xbegin, xend, xsize, xmax_size, x_is_global ) 02043 call mpp_get_compute_domain( domain%y, ybegin, yend, ysize, ymax_size, y_is_global ) 02044 return 02045 end subroutine mpp_get_compute_domain2D 02046 02047 subroutine mpp_get_data_domain2D( domain, xbegin, xend, ybegin, yend, xsize, xmax_size, ysize, ymax_size, & 02048 x_is_global, y_is_global ) 02049 type(domain2D), intent(in) :: domain 02050 integer, intent(out), optional :: xbegin, xend, ybegin, yend, xsize, xmax_size, ysize, ymax_size 02051 logical, intent(out), optional :: x_is_global, y_is_global 02052 call mpp_get_data_domain( domain%x, xbegin, xend, xsize, xmax_size, x_is_global ) 02053 call mpp_get_data_domain( domain%y, ybegin, yend, ysize, ymax_size, y_is_global ) 02054 return 02055 end subroutine mpp_get_data_domain2D 02056 02057 subroutine mpp_get_global_domain2D( domain, xbegin, xend, ybegin, yend, xsize, xmax_size, ysize, ymax_size ) 02058 type(domain2D), intent(in) :: domain 02059 integer, intent(out), optional :: xbegin, xend, ybegin, yend, xsize, xmax_size, ysize, ymax_size 02060 call mpp_get_global_domain( domain%x, xbegin, xend, xsize, xmax_size ) 02061 call mpp_get_global_domain( domain%y, ybegin, yend, ysize, ymax_size ) 02062 return 02063 end subroutine mpp_get_global_domain2D 02064 02065 subroutine mpp_get_domain_components( domain, x, y ) 02066 type(domain2D), intent(in) :: domain 02067 type(domain1D), intent(out), optional :: x, y 02068 if( PRESENT(x) )x = domain%x 02069 if( PRESENT(y) )y = domain%y 02070 return 02071 end subroutine mpp_get_domain_components 02072 02073 subroutine mpp_get_compute_domains1D( domain, begin, end, size ) 02074 type(domain1D), intent(in) :: domain 02075 integer, intent(out), optional, dimension(:) :: begin, end, size 02076 02077 if( .NOT.module_is_initialized ) & 02078 call mpp_error( FATAL, 'MPP_GET_COMPUTE_DOMAINS: must first call mpp_domains_init.' ) 02079 !we use shape instead of size for error checks because size is used as an argument 02080 if( PRESENT(begin) )then 02081 if( any(shape(begin).NE.shape(domain%list)) ) & 02082 call mpp_error( FATAL, 'MPP_GET_COMPUTE_DOMAINS: begin array size does not match domain.' ) 02083 begin(:) = domain%list(:)%compute%begin 02084 end if 02085 if( PRESENT(end) )then 02086 if( any(shape(end).NE.shape(domain%list)) ) & 02087 call mpp_error( FATAL, 'MPP_GET_COMPUTE_DOMAINS: end array size does not match domain.' ) 02088 end(:) = domain%list(:)%compute%end 02089 end if 02090 if( PRESENT(size) )then 02091 if( any(shape(size).NE.shape(domain%list)) ) & 02092 call mpp_error( FATAL, 'MPP_GET_COMPUTE_DOMAINS: size array size does not match domain.' ) 02093 size(:) = domain%list(:)%compute%size 02094 end if 02095 return 02096 end subroutine mpp_get_compute_domains1D 02097 02098 subroutine mpp_get_compute_domains2D( domain, xbegin, xend, xsize, ybegin, yend, ysize ) 02099 type(domain2D), intent(in) :: domain 02100 integer, intent(out), optional, dimension(:) :: xbegin, xend, xsize, ybegin, yend, ysize 02101 02102 if( .NOT.module_is_initialized ) & 02103 call mpp_error( FATAL, 'MPP_GET_COMPUTE_DOMAINS: must first call mpp_domains_init.' ) 02104 02105 if( PRESENT(xbegin) )then 02106 if( size(xbegin).NE.size(domain%list) ) & 02107 call mpp_error( FATAL, 'MPP_GET_COMPUTE_DOMAINS: xbegin array size does not match domain.' ) 02108 xbegin(:) = domain%list(:)%x%compute%begin 02109 end if 02110 if( PRESENT(xend) )then 02111 if( size(xend).NE.size(domain%list) ) & 02112 call mpp_error( FATAL, 'MPP_GET_COMPUTE_DOMAINS: xend array size does not match domain.' ) 02113 xend(:) = domain%list(:)%x%compute%end 02114 end if 02115 if( PRESENT(xsize) )then 02116 if( size(xsize).NE.size(domain%list) ) & 02117 call mpp_error( FATAL, 'MPP_GET_COMPUTE_DOMAINS: xsize array size does not match domain.' ) 02118 xsize(:) = domain%list(:)%x%compute%size 02119 end if 02120 if( PRESENT(ybegin) )then 02121 if( size(ybegin).NE.size(domain%list) ) & 02122 call mpp_error( FATAL, 'MPP_GET_COMPUTE_DOMAINS: ybegin array size does not match domain.' ) 02123 ybegin(:) = domain%list(:)%y%compute%begin 02124 end if 02125 if( PRESENT(yend) )then 02126 if( size(yend).NE.size(domain%list) ) & 02127 call mpp_error( FATAL, 'MPP_GET_COMPUTE_DOMAINS: yend array size does not match domain.' ) 02128 yend(:) = domain%list(:)%y%compute%end 02129 end if 02130 if( PRESENT(ysize) )then 02131 if( size(ysize).NE.size(domain%list) ) & 02132 call mpp_error( FATAL, 'MPP_GET_COMPUTE_DOMAINS: ysize array size does not match domain.' ) 02133 ysize(:) = domain%list(:)%y%compute%size 02134 end if 02135 return 02136 end subroutine mpp_get_compute_domains2D 02137 02138 subroutine mpp_get_pelist1D( domain, pelist, pos ) 02139 type(domain1D), intent(in) :: domain 02140 integer, intent(out) :: pelist(:) 02141 integer, intent(out), optional :: pos 02142 integer :: ndivs 02143 02144 if( .NOT.module_is_initialized ) & 02145 call mpp_error( FATAL, 'MPP_GET_PELIST: must first call mpp_domains_init.' ) 02146 ndivs = size(domain%list) 02147 02148 if( size(pelist).NE.ndivs ) & 02149 call mpp_error( FATAL, 'MPP_GET_PELIST: pelist array size does not match domain.' ) 02150 02151 pelist(:) = domain%list(0:ndivs-1)%pe 02152 if( PRESENT(pos) )pos = domain%pos 02153 return 02154 end subroutine mpp_get_pelist1D 02155 02156 subroutine mpp_get_pelist2D( domain, pelist, pos ) 02157 type(domain2D), intent(in) :: domain 02158 integer, intent(out) :: pelist(:) 02159 integer, intent(out), optional :: pos 02160 02161 if( .NOT.module_is_initialized ) & 02162 call mpp_error( FATAL, 'MPP_GET_PELIST: must first call mpp_domains_init.' ) 02163 if( size(pelist).NE.size(domain%list) ) & 02164 call mpp_error( FATAL, 'MPP_GET_PELIST: pelist array size does not match domain.' ) 02165 02166 pelist(:) = domain%list(:)%pe 02167 if( PRESENT(pos) )pos = domain%pos 02168 return 02169 end subroutine mpp_get_pelist2D 02170 02171 subroutine mpp_get_layout1D( domain, layout ) 02172 type(domain1D), intent(in) :: domain 02173 integer, intent(out) :: layout 02174 02175 if( .NOT.module_is_initialized ) & 02176 call mpp_error( FATAL, 'MPP_GET_LAYOUT: must first call mpp_domains_init.' ) 02177 02178 layout = size(domain%list) 02179 return 02180 end subroutine mpp_get_layout1D 02181 02182 subroutine mpp_get_layout2D( domain, layout ) 02183 type(domain2D), intent(in) :: domain 02184 integer, intent(out) :: layout(2) 02185 02186 if( .NOT.module_is_initialized ) & 02187 call mpp_error( FATAL, 'MPP_GET_LAYOUT: must first call mpp_domains_init.' ) 02188 02189 layout(1) = size(domain%x%list) 02190 layout(2) = size(domain%y%list) 02191 return 02192 end subroutine mpp_get_layout2D 02193 02194 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 02195 ! ! 02196 ! MPP_UPDATE_DOMAINS: fill halos for 2D decomposition ! 02197 ! ! 02198 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 02199 02200 #define VECTOR_FIELD_ 02201 #define MPP_TYPE_ real(DOUBLE_KIND) 02202 #define MPP_UPDATE_DOMAINS_2D_ mpp_update_domain2D_r8_2D 02203 #define MPP_UPDATE_DOMAINS_3D_ mpp_update_domain2D_r8_3D 02204 #define MPP_UPDATE_DOMAINS_4D_ mpp_update_domain2D_r8_4D 02205 #define MPP_UPDATE_DOMAINS_5D_ mpp_update_domain2D_r8_5D 02206 02207 #ifdef VECTOR_FIELD_ 02208 #define MPP_UPDATE_DOMAINS_2D_V_ mpp_update_domain2D_r8_2Dv 02209 #define MPP_UPDATE_DOMAINS_3D_V_ mpp_update_domain2D_r8_3Dv 02210 #define MPP_UPDATE_DOMAINS_4D_V_ mpp_update_domain2D_r8_4Dv 02211 #define MPP_UPDATE_DOMAINS_5D_V_ mpp_update_domain2D_r8_5Dv 02212 #endif 02213 #define MPP_REDISTRIBUTE_2D_ mpp_redistribute_r8_2D 02214 #define MPP_REDISTRIBUTE_3D_ mpp_redistribute_r8_3D 02215 #define MPP_REDISTRIBUTE_4D_ mpp_redistribute_r8_4D 02216 #define MPP_REDISTRIBUTE_5D_ mpp_redistribute_r8_5D 02217 #include <mpp_update_domains2D.h> 02218 #undef VECTOR_FIELD_ 02219 02220 #define MPP_TYPE_ complex(DOUBLE_KIND) 02221 #define MPP_UPDATE_DOMAINS_2D_ mpp_update_domain2D_c8_2D 02222 #define MPP_UPDATE_DOMAINS_3D_ mpp_update_domain2D_c8_3D 02223 #define MPP_UPDATE_DOMAINS_4D_ mpp_update_domain2D_c8_4D 02224 #define MPP_UPDATE_DOMAINS_5D_ mpp_update_domain2D_c8_5D 02225 #define MPP_REDISTRIBUTE_2D_ mpp_redistribute_c8_2D 02226 #define MPP_REDISTRIBUTE_3D_ mpp_redistribute_c8_3D 02227 #define MPP_REDISTRIBUTE_4D_ mpp_redistribute_c8_4D 02228 #define MPP_REDISTRIBUTE_5D_ mpp_redistribute_c8_5D 02229 #include <mpp_update_domains2D.h> 02230 02231 #ifndef no_8byte_integers 02232 #define MPP_TYPE_ integer(LONG_KIND) 02233 #define MPP_UPDATE_DOMAINS_2D_ mpp_update_domain2D_i8_2D 02234 #define MPP_UPDATE_DOMAINS_3D_ mpp_update_domain2D_i8_3D 02235 #define MPP_UPDATE_DOMAINS_4D_ mpp_update_domain2D_i8_4D 02236 #define MPP_UPDATE_DOMAINS_5D_ mpp_update_domain2D_i8_5D 02237 #define MPP_REDISTRIBUTE_2D_ mpp_redistribute_i8_2D 02238 #define MPP_REDISTRIBUTE_3D_ mpp_redistribute_i8_3D 02239 #define MPP_REDISTRIBUTE_4D_ mpp_redistribute_i8_4D 02240 #define MPP_REDISTRIBUTE_5D_ mpp_redistribute_i8_5D 02241 #include <mpp_update_domains2D.h> 02242 02243 #define MPP_TYPE_ logical(LONG_KIND) 02244 #define MPP_UPDATE_DOMAINS_2D_ mpp_update_domain2D_l8_2D 02245 #define MPP_UPDATE_DOMAINS_3D_ mpp_update_domain2D_l8_3D 02246 #define MPP_UPDATE_DOMAINS_4D_ mpp_update_domain2D_l8_4D 02247 #define MPP_UPDATE_DOMAINS_5D_ mpp_update_domain2D_l8_5D 02248 #define MPP_REDISTRIBUTE_2D_ mpp_redistribute_l8_2D 02249 #define MPP_REDISTRIBUTE_3D_ mpp_redistribute_l8_3D 02250 #define MPP_REDISTRIBUTE_4D_ mpp_redistribute_l8_4D 02251 #define MPP_REDISTRIBUTE_5D_ mpp_redistribute_l8_5D 02252 #include <mpp_update_domains2D.h> 02253 #endif 02254 02255 #ifndef no_4byte_reals 02256 #define VECTOR_FIELD_ 02257 #define MPP_TYPE_ real(FLOAT_KIND) 02258 #define MPP_UPDATE_DOMAINS_2D_ mpp_update_domain2D_r4_2D 02259 #define MPP_UPDATE_DOMAINS_3D_ mpp_update_domain2D_r4_3D 02260 #define MPP_UPDATE_DOMAINS_4D_ mpp_update_domain2D_r4_4D 02261 #define MPP_UPDATE_DOMAINS_5D_ mpp_update_domain2D_r4_5D 02262 #ifdef VECTOR_FIELD_ 02263 #define MPP_UPDATE_DOMAINS_2D_V_ mpp_update_domain2D_r4_2Dv 02264 #define MPP_UPDATE_DOMAINS_3D_V_ mpp_update_domain2D_r4_3Dv 02265 #define MPP_UPDATE_DOMAINS_4D_V_ mpp_update_domain2D_r4_4Dv 02266 #define MPP_UPDATE_DOMAINS_5D_V_ mpp_update_domain2D_r4_5Dv 02267 #endif 02268 #define MPP_REDISTRIBUTE_2D_ mpp_redistribute_r4_2D 02269 #define MPP_REDISTRIBUTE_3D_ mpp_redistribute_r4_3D 02270 #define MPP_REDISTRIBUTE_4D_ mpp_redistribute_r4_4D 02271 #define MPP_REDISTRIBUTE_5D_ mpp_redistribute_r4_5D 02272 #include <mpp_update_domains2D.h> 02273 #undef VECTOR_FIELD_ 02274 #endif 02275 02276 #ifndef no_4byte_cmplx 02277 #define MPP_TYPE_ complex(FLOAT_KIND) 02278 #define MPP_UPDATE_DOMAINS_2D_ mpp_update_domain2D_c4_2D 02279 #define MPP_UPDATE_DOMAINS_3D_ mpp_update_domain2D_c4_3D 02280 #define MPP_UPDATE_DOMAINS_4D_ mpp_update_domain2D_c4_4D 02281 #define MPP_UPDATE_DOMAINS_5D_ mpp_update_domain2D_c4_5D 02282 #define MPP_REDISTRIBUTE_2D_ mpp_redistribute_c4_2D 02283 #define MPP_REDISTRIBUTE_3D_ mpp_redistribute_c4_3D 02284 #define MPP_REDISTRIBUTE_4D_ mpp_redistribute_c4_4D 02285 #define MPP_REDISTRIBUTE_5D_ mpp_redistribute_c4_5D 02286 #include <mpp_update_domains2D.h> 02287 #endif 02288 02289 #define MPP_TYPE_ integer(INT_KIND) 02290 #define MPP_UPDATE_DOMAINS_2D_ mpp_update_domain2D_i4_2D 02291 #define MPP_UPDATE_DOMAINS_3D_ mpp_update_domain2D_i4_3D 02292 #define MPP_UPDATE_DOMAINS_4D_ mpp_update_domain2D_i4_4D 02293 #define MPP_UPDATE_DOMAINS_5D_ mpp_update_domain2D_i4_5D 02294 #define MPP_REDISTRIBUTE_2D_ mpp_redistribute_i4_2D 02295 #define MPP_REDISTRIBUTE_3D_ mpp_redistribute_i4_3D 02296 #define MPP_REDISTRIBUTE_4D_ mpp_redistribute_i4_4D 02297 #define MPP_REDISTRIBUTE_5D_ mpp_redistribute_i4_5D 02298 #include <mpp_update_domains2D.h> 02299 02300 #define MPP_TYPE_ logical(INT_KIND) 02301 #define MPP_UPDATE_DOMAINS_2D_ mpp_update_domain2D_l4_2D 02302 #define MPP_UPDATE_DOMAINS_3D_ mpp_update_domain2D_l4_3D 02303 #define MPP_UPDATE_DOMAINS_4D_ mpp_update_domain2D_l4_4D 02304 #define MPP_UPDATE_DOMAINS_5D_ mpp_update_domain2D_l4_5D 02305 #define MPP_REDISTRIBUTE_2D_ mpp_redistribute_l4_2D 02306 #define MPP_REDISTRIBUTE_3D_ mpp_redistribute_l4_3D 02307 #define MPP_REDISTRIBUTE_4D_ mpp_redistribute_l4_4D 02308 #define MPP_REDISTRIBUTE_5D_ mpp_redistribute_l4_5D 02309 #include <mpp_update_domains2D.h> 02310 02311 !#define MPP_UPDATE_DOMAINS_2D_ mpp_update_domain1D_r8_2D 02312 !#define MPP_UPDATE_DOMAINS_3D_ mpp_update_domain1D_r8_3D 02313 !#define MPP_UPDATE_DOMAINS_4D_ mpp_update_domain1D_r8_4D 02314 !#define MPP_UPDATE_DOMAINS_5D_ mpp_update_domain1D_r8_5D 02315 !#define MPP_TYPE_ real(DOUBLE_KIND) 02316 !#include <mpp_update_domains1D.h> 02317 ! 02318 !#define MPP_UPDATE_DOMAINS_2D_ mpp_update_domain1D_c8_2D 02319 !#define MPP_UPDATE_DOMAINS_3D_ mpp_update_domain1D_c8_3D 02320 !#define MPP_UPDATE_DOMAINS_4D_ mpp_update_domain1D_c8_4D 02321 !#define MPP_UPDATE_DOMAINS_5D_ mpp_update_domain1D_c8_5D 02322 !#define MPP_TYPE_ complex(DOUBLE_KIND) 02323 !#include <mpp_update_domains1D.h> 02324 ! 02325 !#ifndef no_8byte_integers 02326 !#define MPP_UPDATE_DOMAINS_2D_ mpp_update_domain1D_i8_2D 02327 !#define MPP_UPDATE_DOMAINS_3D_ mpp_update_domain1D_i8_3D 02328 !#define MPP_UPDATE_DOMAINS_4D_ mpp_update_domain1D_i8_4D 02329 !#define MPP_UPDATE_DOMAINS_5D_ mpp_update_domain1D_i8_5D 02330 !#define MPP_TYPE_ integer(LONG_KIND) 02331 !#include <mpp_update_domains1D.h> 02332 ! 02333 !#define MPP_UPDATE_DOMAINS_2D_ mpp_update_domain1D_l8_2D 02334 !#define MPP_UPDATE_DOMAINS_3D_ mpp_update_domain1D_l8_3D 02335 !#define MPP_UPDATE_DOMAINS_4D_ mpp_update_domain1D_l8_4D 02336 !#define MPP_UPDATE_DOMAINS_5D_ mpp_update_domain1D_l8_5D 02337 !#define MPP_TYPE_ logical(LONG_KIND) 02338 !#include <mpp_update_domains1D.h> 02339 !#endif 02340 ! 02341 !#define MPP_UPDATE_DOMAINS_2D_ mpp_update_domain1D_r4_2D 02342 !#define MPP_UPDATE_DOMAINS_3D_ mpp_update_domain1D_r4_3D 02343 !#define MPP_UPDATE_DOMAINS_4D_ mpp_update_domain1D_r4_4D 02344 !#define MPP_UPDATE_DOMAINS_5D_ mpp_update_domain1D_r4_5D 02345 !#define MPP_TYPE_ real(FLOAT_KIND) 02346 !#include <mpp_update_domains1D.h> 02347 ! 02348 !#define MPP_UPDATE_DOMAINS_2D_ mpp_update_domain1D_c4_2D 02349 !#define MPP_UPDATE_DOMAINS_3D_ mpp_update_domain1D_c4_3D 02350 !#define MPP_UPDATE_DOMAINS_4D_ mpp_update_domain1D_c4_4D 02351 !#define MPP_UPDATE_DOMAINS_5D_ mpp_update_domain1D_c4_5D 02352 !#define MPP_TYPE_ complex(FLOAT_KIND) 02353 !#include <mpp_update_domains1D.h> 02354 ! 02355 !#define MPP_UPDATE_DOMAINS_2D_ mpp_update_domain1D_i4_2D 02356 !#define MPP_UPDATE_DOMAINS_3D_ mpp_update_domain1D_i4_3D 02357 !#define MPP_UPDATE_DOMAINS_4D_ mpp_update_domain1D_i4_4D 02358 !#define MPP_UPDATE_DOMAINS_5D_ mpp_update_domain1D_i4_5D 02359 !#define MPP_TYPE_ integer(INT_KIND) 02360 !#include <mpp_update_domains1D.h> 02361 ! 02362 !#define MPP_UPDATE_DOMAINS_2D_ mpp_update_domain1D_l4_2D 02363 !#define MPP_UPDATE_DOMAINS_3D_ mpp_update_domain1D_l4_3D 02364 !#define MPP_UPDATE_DOMAINS_4D_ mpp_update_domain1D_l4_4D 02365 !#define MPP_UPDATE_DOMAINS_5D_ mpp_update_domain1D_l4_5D 02366 !#define MPP_TYPE_ logical(INT_KIND) 02367 !#include <mpp_update_domains1D.h> 02368 02369 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 02370 ! ! 02371 ! MPP_GLOBAL_REDUCE: get global max/min of field ! 02372 ! ! 02373 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 02374 02375 #define MPP_GLOBAL_REDUCE_2D_ mpp_global_max_r8_2d 02376 #define MPP_GLOBAL_REDUCE_3D_ mpp_global_max_r8_3d 02377 #define MPP_GLOBAL_REDUCE_4D_ mpp_global_max_r8_4d 02378 #define MPP_GLOBAL_REDUCE_5D_ mpp_global_max_r8_5d 02379 #define MPP_TYPE_ real(DOUBLE_KIND) 02380 #define REDUCE_VAL_ maxval 02381 #define REDUCE_LOC_ maxloc 02382 #define MPP_REDUCE_ mpp_max 02383 #include <mpp_global_reduce.h> 02384 02385 #define MPP_GLOBAL_REDUCE_2D_ mpp_global_min_r8_2d 02386 #define MPP_GLOBAL_REDUCE_3D_ mpp_global_min_r8_3d 02387 #define MPP_GLOBAL_REDUCE_4D_ mpp_global_min_r8_4d 02388 #define MPP_GLOBAL_REDUCE_5D_ mpp_global_min_r8_5d 02389 #define MPP_TYPE_ real(DOUBLE_KIND) 02390 #define REDUCE_VAL_ minval 02391 #define REDUCE_LOC_ minloc 02392 #define MPP_REDUCE_ mpp_min 02393 #include <mpp_global_reduce.h> 02394 02395 #ifndef no_4byte_reals 02396 #define MPP_GLOBAL_REDUCE_2D_ mpp_global_max_r4_2d 02397 #define MPP_GLOBAL_REDUCE_3D_ mpp_global_max_r4_3d 02398 #define MPP_GLOBAL_REDUCE_4D_ mpp_global_max_r4_4d 02399 #define MPP_GLOBAL_REDUCE_5D_ mpp_global_max_r4_5d 02400 #define MPP_TYPE_ real(FLOAT_KIND) 02401 #define REDUCE_VAL_ maxval 02402 #define REDUCE_LOC_ maxloc 02403 #define MPP_REDUCE_ mpp_max 02404 #include <mpp_global_reduce.h> 02405 02406 #define MPP_GLOBAL_REDUCE_2D_ mpp_global_min_r4_2d 02407 #define MPP_GLOBAL_REDUCE_3D_ mpp_global_min_r4_3d 02408 #define MPP_GLOBAL_REDUCE_4D_ mpp_global_min_r4_4d 02409 #define MPP_GLOBAL_REDUCE_5D_ mpp_global_min_r4_5d 02410 #define MPP_TYPE_ real(FLOAT_KIND) 02411 #define REDUCE_VAL_ minval 02412 #define REDUCE_LOC_ minloc 02413 #define MPP_REDUCE_ mpp_min 02414 #include <mpp_global_reduce.h> 02415 #endif 02416 02417 #ifndef no_8byte_integers 02418 #define MPP_GLOBAL_REDUCE_2D_ mpp_global_max_i8_2d 02419 #define MPP_GLOBAL_REDUCE_3D_ mpp_global_max_i8_3d 02420 #define MPP_GLOBAL_REDUCE_4D_ mpp_global_max_i8_4d 02421 #define MPP_GLOBAL_REDUCE_5D_ mpp_global_max_i8_5d 02422 #define MPP_TYPE_ integer(LONG_KIND) 02423 #define REDUCE_VAL_ maxval 02424 #define REDUCE_LOC_ maxloc 02425 #define MPP_REDUCE_ mpp_max 02426 #include <mpp_global_reduce.h> 02427 02428 #define MPP_GLOBAL_REDUCE_2D_ mpp_global_min_i8_2d 02429 #define MPP_GLOBAL_REDUCE_3D_ mpp_global_min_i8_3d 02430 #define MPP_GLOBAL_REDUCE_4D_ mpp_global_min_i8_4d 02431 #define MPP_GLOBAL_REDUCE_5D_ mpp_global_min_i8_5d 02432 #define MPP_TYPE_ integer(LONG_KIND) 02433 #define REDUCE_VAL_ minval 02434 #define REDUCE_LOC_ minloc 02435 #define MPP_REDUCE_ mpp_min 02436 #include <mpp_global_reduce.h> 02437 #endif 02438 02439 #define MPP_GLOBAL_REDUCE_2D_ mpp_global_max_i4_2d 02440 #define MPP_GLOBAL_REDUCE_3D_ mpp_global_max_i4_3d 02441 #define MPP_GLOBAL_REDUCE_4D_ mpp_global_max_i4_4d 02442 #define MPP_GLOBAL_REDUCE_5D_ mpp_global_max_i4_5d 02443 #define MPP_TYPE_ integer(INT_KIND) 02444 #define REDUCE_VAL_ maxval 02445 #define REDUCE_LOC_ maxloc 02446 #define MPP_REDUCE_ mpp_max 02447 #include <mpp_global_reduce.h> 02448 02449 #define MPP_GLOBAL_REDUCE_2D_ mpp_global_min_i4_2d 02450 #define MPP_GLOBAL_REDUCE_3D_ mpp_global_min_i4_3d 02451 #define MPP_GLOBAL_REDUCE_4D_ mpp_global_min_i4_4d 02452 #define MPP_GLOBAL_REDUCE_5D_ mpp_global_min_i4_5d 02453 #define MPP_TYPE_ integer(INT_KIND) 02454 #define REDUCE_VAL_ minval 02455 #define REDUCE_LOC_ minloc 02456 #define MPP_REDUCE_ mpp_min 02457 #include <mpp_global_reduce.h> 02458 02459 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 02460 ! ! 02461 ! MPP_GLOBAL_SUM: global sum of field ! 02462 ! ! 02463 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 02464 02465 #define MPP_GLOBAL_SUM_ mpp_global_sum_r8_2d 02466 #define MPP_EXTRA_INDICES_ 02467 #define MPP_TYPE_ real(DOUBLE_KIND) 02468 #include <mpp_global_sum.h> 02469 02470 #define MPP_GLOBAL_SUM_ mpp_global_sum_r8_3d 02471 #define MPP_EXTRA_INDICES_ ,: 02472 #define MPP_TYPE_ real(DOUBLE_KIND) 02473 #include <mpp_global_sum.h> 02474 02475 #define MPP_GLOBAL_SUM_ mpp_global_sum_r8_4d 02476 #define MPP_EXTRA_INDICES_ ,:,: 02477 #define MPP_TYPE_ real(DOUBLE_KIND) 02478 #include <mpp_global_sum.h> 02479 02480 #define MPP_GLOBAL_SUM_ mpp_global_sum_r8_5d 02481 #define MPP_EXTRA_INDICES_ ,:,:,: 02482 #define MPP_TYPE_ real(DOUBLE_KIND) 02483 #include <mpp_global_sum.h> 02484 02485 #ifndef no_4byte_reals 02486 #define MPP_GLOBAL_SUM_ mpp_global_sum_r4_2d 02487 #define MPP_EXTRA_INDICES_ 02488 #define MPP_TYPE_ real(FLOAT_KIND) 02489 #include <mpp_global_sum.h> 02490 02491 #define MPP_GLOBAL_SUM_ mpp_global_sum_r4_3d 02492 #define MPP_EXTRA_INDICES_ ,: 02493 #define MPP_TYPE_ real(FLOAT_KIND) 02494 #include <mpp_global_sum.h> 02495 02496 #define MPP_GLOBAL_SUM_ mpp_global_sum_r4_4d 02497 #define MPP_EXTRA_INDICES_ ,:,: 02498 #define MPP_TYPE_ real(FLOAT_KIND) 02499 #include <mpp_global_sum.h> 02500 02501 #define MPP_GLOBAL_SUM_ mpp_global_sum_r4_5d 02502 #define MPP_EXTRA_INDICES_ ,:,:,: 02503 #define MPP_TYPE_ real(FLOAT_KIND) 02504 #include <mpp_global_sum.h> 02505 #endif 02506 02507 #define MPP_GLOBAL_SUM_ mpp_global_sum_c8_2d 02508 #define MPP_EXTRA_INDICES_ 02509 #define MPP_TYPE_ complex(DOUBLE_KIND) 02510 #include <mpp_global_sum.h> 02511 02512 #define MPP_GLOBAL_SUM_ mpp_global_sum_c8_3d 02513 #define MPP_EXTRA_INDICES_ ,: 02514 #define MPP_TYPE_ complex(DOUBLE_KIND) 02515 #include <mpp_global_sum.h> 02516 02517 #define MPP_GLOBAL_SUM_ mpp_global_sum_c8_4d 02518 #define MPP_EXTRA_INDICES_ ,:,: 02519 #define MPP_TYPE_ complex(DOUBLE_KIND) 02520 #include <mpp_global_sum.h> 02521 02522 #define MPP_GLOBAL_SUM_ mpp_global_sum_c8_5d 02523 #define MPP_EXTRA_INDICES_ ,:,:,: 02524 #define MPP_TYPE_ complex(DOUBLE_KIND) 02525 #include <mpp_global_sum.h> 02526 02527 #ifndef no_4byte_cmplx 02528 #define MPP_GLOBAL_SUM_ mpp_global_sum_c4_2d 02529 #define MPP_EXTRA_INDICES_ 02530 #define MPP_TYPE_ complex(FLOAT_KIND) 02531 #include <mpp_global_sum.h> 02532 02533 #define MPP_GLOBAL_SUM_ mpp_global_sum_c4_3d 02534 #define MPP_EXTRA_INDICES_ ,: 02535 #define MPP_TYPE_ complex(FLOAT_KIND) 02536 #include <mpp_global_sum.h> 02537 02538 #define MPP_GLOBAL_SUM_ mpp_global_sum_c4_4d 02539 #define MPP_EXTRA_INDICES_ ,:,: 02540 #define MPP_TYPE_ complex(FLOAT_KIND) 02541 #include <mpp_global_sum.h> 02542 02543 #define MPP_GLOBAL_SUM_ mpp_global_sum_c4_5d 02544 #define MPP_EXTRA_INDICES_ ,:,:,: 02545 #define MPP_TYPE_ complex(FLOAT_KIND) 02546 #include <mpp_global_sum.h> 02547 #endif 02548 02549 #ifndef no_8byte_integers 02550 #define MPP_GLOBAL_SUM_ mpp_global_sum_i8_2d 02551 #define MPP_EXTRA_INDICES_ 02552 #define MPP_TYPE_ integer(LONG_KIND) 02553 #include <mpp_global_sum.h> 02554 02555 #define MPP_GLOBAL_SUM_ mpp_global_sum_i8_3d 02556 #define MPP_EXTRA_INDICES_ ,: 02557 #define MPP_TYPE_ integer(LONG_KIND) 02558 #include <mpp_global_sum.h> 02559 02560 #define MPP_GLOBAL_SUM_ mpp_global_sum_i8_4d 02561 #define MPP_EXTRA_INDICES_ ,:,: 02562 #define MPP_TYPE_ integer(LONG_KIND) 02563 #include <mpp_global_sum.h> 02564 02565 #define MPP_GLOBAL_SUM_ mpp_global_sum_i8_5d 02566 #define MPP_EXTRA_INDICES_ ,:,:,: 02567 #define MPP_TYPE_ integer(LONG_KIND) 02568 #include <mpp_global_sum.h> 02569 #endif 02570 02571 #define MPP_GLOBAL_SUM_ mpp_global_sum_i4_2d 02572 #define MPP_EXTRA_INDICES_ 02573 #define MPP_TYPE_ integer(INT_KIND) 02574 #include <mpp_global_sum.h> 02575 02576 #define MPP_GLOBAL_SUM_ mpp_global_sum_i4_3d 02577 #define MPP_EXTRA_INDICES_ ,: 02578 #define MPP_TYPE_ integer(INT_KIND) 02579 #include <mpp_global_sum.h> 02580 02581 #define MPP_GLOBAL_SUM_ mpp_global_sum_i4_4d 02582 #define MPP_EXTRA_INDICES_ ,:,: 02583 #define MPP_TYPE_ integer(INT_KIND) 02584 #include <mpp_global_sum.h> 02585 02586 #define MPP_GLOBAL_SUM_ mpp_global_sum_i4_5d 02587 #define MPP_EXTRA_INDICES_ ,:,:,: 02588 #define MPP_TYPE_ integer(INT_KIND) 02589 #include <mpp_global_sum.h> 02590 02591 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 02592 ! ! 02593 ! MPP_GLOBAL_FIELD: get global field from domain field ! 02594 ! ! 02595 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 02596 02597 #define MPP_GLOBAL_FIELD_2D_ mpp_global_field2D_r8_2d 02598 #define MPP_GLOBAL_FIELD_3D_ mpp_global_field2D_r8_3d 02599 #define MPP_GLOBAL_FIELD_4D_ mpp_global_field2D_r8_4d 02600 #define MPP_GLOBAL_FIELD_5D_ mpp_global_field2D_r8_5d 02601 #define MPP_GLOBAL1D_FIELD_2D_ mpp_global_field1D_r8_2d 02602 #define MPP_TYPE_ real(DOUBLE_KIND) 02603 #include <mpp_global_field.h> 02604 02605 #define MPP_GLOBAL_FIELD_2D_ mpp_global_field2D_c8_2d 02606 #define MPP_GLOBAL_FIELD_3D_ mpp_global_field2D_c8_3d 02607 #define MPP_GLOBAL_FIELD_4D_ mpp_global_field2D_c8_4d 02608 #define MPP_GLOBAL_FIELD_5D_ mpp_global_field2D_c8_5d 02609 #define MPP_GLOBAL1D_FIELD_2D_ mpp_global_field1D_c8_2d 02610 #define MPP_TYPE_ complex(DOUBLE_KIND) 02611 #include <mpp_global_field.h> 02612 02613 #ifndef no_8byte_integers 02614 #define MPP_GLOBAL_FIELD_2D_ mpp_global_field2D_i8_2d 02615 #define MPP_GLOBAL_FIELD_3D_ mpp_global_field2D_i8_3d 02616 #define MPP_GLOBAL_FIELD_4D_ mpp_global_field2D_i8_4d 02617 #define MPP_GLOBAL_FIELD_5D_ mpp_global_field2D_i8_5d 02618 #define MPP_GLOBAL1D_FIELD_2D_ mpp_global_field1D_i8_2d 02619 #define MPP_TYPE_ integer(LONG_KIND) 02620 #include <mpp_global_field.h> 02621 02622 #define MPP_GLOBAL_FIELD_2D_ mpp_global_field2D_l8_2d 02623 #define MPP_GLOBAL_FIELD_3D_ mpp_global_field2D_l8_3d 02624 #define MPP_GLOBAL_FIELD_4D_ mpp_global_field2D_l8_4d 02625 #define MPP_GLOBAL_FIELD_5D_ mpp_global_field2D_l8_5d 02626 #define MPP_GLOBAL1D_FIELD_2D_ mpp_global_field1D_l8_2d 02627 #define MPP_TYPE_ logical(LONG_KIND) 02628 #include <mpp_global_field.h> 02629 #endif 02630 02631 #ifndef no_4byte_reals 02632 #define MPP_GLOBAL_FIELD_2D_ mpp_global_field2D_r4_2d 02633 #define MPP_GLOBAL_FIELD_3D_ mpp_global_field2D_r4_3d 02634 #define MPP_GLOBAL_FIELD_4D_ mpp_global_field2D_r4_4d 02635 #define MPP_GLOBAL_FIELD_5D_ mpp_global_field2D_r4_5d 02636 #define MPP_GLOBAL1D_FIELD_2D_ mpp_global_field1D_r4_2d 02637 #define MPP_TYPE_ real(FLOAT_KIND) 02638 #include <mpp_global_field.h> 02639 #endif 02640 02641 #ifndef no_4byte_cmplx 02642 #define MPP_GLOBAL_FIELD_2D_ mpp_global_field2D_c4_2d 02643 #define MPP_GLOBAL_FIELD_3D_ mpp_global_field2D_c4_3d 02644 #define MPP_GLOBAL_FIELD_4D_ mpp_global_field2D_c4_4d 02645 #define MPP_GLOBAL_FIELD_5D_ mpp_global_field2D_c4_5d 02646 #define MPP_GLOBAL1D_FIELD_2D_ mpp_global_field1D_c4_2d 02647 #define MPP_TYPE_ complex(FLOAT_KIND) 02648 #include <mpp_global_field.h> 02649 #endif 02650 02651 #define MPP_GLOBAL_FIELD_2D_ mpp_global_field2D_i4_2d 02652 #define MPP_GLOBAL_FIELD_3D_ mpp_global_field2D_i4_3d 02653 #define MPP_GLOBAL_FIELD_4D_ mpp_global_field2D_i4_4d 02654 #define MPP_GLOBAL_FIELD_5D_ mpp_global_field2D_i4_5d 02655 #define MPP_GLOBAL1D_FIELD_2D_ mpp_global_field1D_i4_2d 02656 #define MPP_TYPE_ integer(INT_KIND) 02657 #include <mpp_global_field.h> 02658 02659 #define MPP_GLOBAL_FIELD_2D_ mpp_global_field2D_l4_2d 02660 #define MPP_GLOBAL_FIELD_3D_ mpp_global_field2D_l4_3d 02661 #define MPP_GLOBAL_FIELD_4D_ mpp_global_field2D_l4_4d 02662 #define MPP_GLOBAL_FIELD_5D_ mpp_global_field2D_l4_5d 02663 #define MPP_GLOBAL1D_FIELD_2D_ mpp_global_field1D_l4_2d 02664 #define MPP_TYPE_ logical(INT_KIND) 02665 #include <mpp_global_field.h> 02666 02667 end module mpp_domains_mod_oa 02668 02669 #endif 02670