[1682] | 1 | !> @file init_pegrid.f90 |
---|
[2000] | 2 | !------------------------------------------------------------------------------! |
---|
[2696] | 3 | ! This file is part of the PALM model system. |
---|
[1036] | 4 | ! |
---|
[2000] | 5 | ! PALM is free software: you can redistribute it and/or modify it under the |
---|
| 6 | ! terms of the GNU General Public License as published by the Free Software |
---|
| 7 | ! Foundation, either version 3 of the License, or (at your option) any later |
---|
| 8 | ! version. |
---|
[1036] | 9 | ! |
---|
| 10 | ! PALM is distributed in the hope that it will be useful, but WITHOUT ANY |
---|
| 11 | ! WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR |
---|
| 12 | ! A PARTICULAR PURPOSE. See the GNU General Public License for more details. |
---|
| 13 | ! |
---|
| 14 | ! You should have received a copy of the GNU General Public License along with |
---|
| 15 | ! PALM. If not, see <http://www.gnu.org/licenses/>. |
---|
| 16 | ! |
---|
[4360] | 17 | ! Copyright 1997-2020 Leibniz Universitaet Hannover |
---|
[2000] | 18 | !------------------------------------------------------------------------------! |
---|
[1036] | 19 | ! |
---|
[254] | 20 | ! Current revisions: |
---|
[1322] | 21 | ! ------------------ |
---|
[2201] | 22 | ! |
---|
[3553] | 23 | ! |
---|
[1321] | 24 | ! Former revisions: |
---|
| 25 | ! ----------------- |
---|
| 26 | ! $Id: init_pegrid.f90 4564 2020-06-12 14:03:36Z eckhard $ |
---|
[4564] | 27 | ! Vertical nesting method of Huq et al. (2019) removed |
---|
| 28 | ! |
---|
| 29 | ! 4461 2020-03-12 16:51:59Z raasch |
---|
[4461] | 30 | ! communicator configurations for four virtual pe grids defined |
---|
| 31 | ! |
---|
| 32 | ! 4444 2020-03-05 15:59:50Z raasch |
---|
[4444] | 33 | ! bugfix: cpp-directives for serial mode added |
---|
| 34 | ! |
---|
| 35 | ! 4360 2020-01-07 11:25:50Z suehring |
---|
[4301] | 36 | ! changed message PA0467 |
---|
| 37 | ! |
---|
| 38 | ! 4264 2019-10-15 16:00:23Z scharf |
---|
[4264] | 39 | ! corrected error message string |
---|
| 40 | ! |
---|
| 41 | ! 4241 2019-09-27 06:32:47Z raasch |
---|
[4241] | 42 | ! Check added to ensure that subdomain grid has at least the size as given by the number |
---|
| 43 | ! of ghost points |
---|
| 44 | ! |
---|
| 45 | ! 4182 2019-08-22 15:20:23Z scharf |
---|
[4182] | 46 | ! Corrected "Former revisions" section |
---|
| 47 | ! |
---|
| 48 | ! 4045 2019-06-21 10:58:47Z raasch |
---|
[4045] | 49 | ! bugfix: kind attribute added to nint function to allow for large integers which may appear in |
---|
| 50 | ! case of default recycling width and small grid spacings |
---|
| 51 | ! |
---|
| 52 | ! 3999 2019-05-23 16:09:37Z suehring |
---|
[3999] | 53 | ! Spend 3 ghost points also in case of pw-scheme when nesting is applied |
---|
| 54 | ! |
---|
| 55 | ! 3897 2019-04-15 11:51:14Z suehring |
---|
[3897] | 56 | ! Minor revision of multigrid check; give warning instead of an abort. |
---|
| 57 | ! |
---|
| 58 | ! 3890 2019-04-12 15:59:20Z suehring |
---|
[3890] | 59 | ! Check if grid coarsening is possible on subdomain, in order to avoid that |
---|
| 60 | ! multigrid approach effectively reduces to a Gauss-Seidel scheme. |
---|
| 61 | ! |
---|
| 62 | ! 3885 2019-04-11 11:29:34Z kanani |
---|
[3885] | 63 | ! Changes related to global restructuring of location messages and introduction |
---|
| 64 | ! of additional debug messages |
---|
| 65 | ! |
---|
| 66 | ! 3884 2019-04-10 13:31:55Z Giersch |
---|
[3884] | 67 | ! id_recycling is only calculated in case of tubulent inflow |
---|
| 68 | ! |
---|
| 69 | ! 3761 2019-02-25 15:31:42Z raasch |
---|
[3761] | 70 | ! unused variable removed |
---|
| 71 | ! |
---|
| 72 | ! 3655 2019-01-07 16:51:22Z knoop |
---|
[3553] | 73 | ! variables documented |
---|
[2259] | 74 | ! |
---|
[4182] | 75 | ! Revision 1.1 1997/07/24 11:15:09 raasch |
---|
| 76 | ! Initial revision |
---|
| 77 | ! |
---|
| 78 | ! |
---|
[1] | 79 | ! Description: |
---|
| 80 | ! ------------ |
---|
[1682] | 81 | !> Determination of the virtual processor topology (if not prescribed by the |
---|
| 82 | !> user)and computation of the grid point number and array bounds of the local |
---|
| 83 | !> domains. |
---|
[2696] | 84 | !> @todo: remove MPI-data types for 2D exchange on coarse multigrid level (not |
---|
| 85 | !> used any more) |
---|
[1] | 86 | !------------------------------------------------------------------------------! |
---|
[1682] | 87 | SUBROUTINE init_pegrid |
---|
| 88 | |
---|
[1] | 89 | |
---|
[1320] | 90 | USE control_parameters, & |
---|
[3182] | 91 | ONLY: bc_dirichlet_l, bc_dirichlet_n, bc_dirichlet_r, bc_dirichlet_s, & |
---|
| 92 | bc_lr, bc_ns, bc_radiation_l, bc_radiation_n, bc_radiation_r, & |
---|
[4444] | 93 | bc_radiation_s, & |
---|
[3241] | 94 | grid_level, grid_level_count, maximum_grid_level, & |
---|
[4444] | 95 | message_string, mg_switch_to_pe0_level, & |
---|
| 96 | psolver |
---|
[1] | 97 | |
---|
[4444] | 98 | |
---|
| 99 | #if defined( __parallel ) |
---|
| 100 | USE control_parameters, & |
---|
| 101 | ONLY: coupling_mode, coupling_topology, gathered_size, momentum_advec, & |
---|
| 102 | outflow_source_plane, recycling_width, scalar_advec, subdomain_size, & |
---|
| 103 | turbulent_inflow, turbulent_outflow, y_shift |
---|
| 104 | |
---|
[1320] | 105 | USE grid_variables, & |
---|
| 106 | ONLY: dx |
---|
[4444] | 107 | #endif |
---|
[1320] | 108 | |
---|
| 109 | USE indices, & |
---|
[4444] | 110 | ONLY: nnx, nny, nnz, nx, nxl, nxl_mg, & |
---|
| 111 | nxlu, nxr, nxr_mg, ny, nyn, nyn_mg, nys, nys_mg, & |
---|
[1320] | 112 | nysv, nz, nzb, nzt, nzt_mg, wall_flags_1, wall_flags_2, & |
---|
| 113 | wall_flags_3, wall_flags_4, wall_flags_5, wall_flags_6, & |
---|
| 114 | wall_flags_7, wall_flags_8, wall_flags_9, wall_flags_10 |
---|
[1] | 115 | |
---|
[4444] | 116 | #if defined( __parallel ) |
---|
| 117 | USE indices, & |
---|
| 118 | ONLY: mg_loc_ind, nbgp, nx_a, nx_o, ny_a, ny_o |
---|
| 119 | #endif |
---|
| 120 | |
---|
[1320] | 121 | USE kinds |
---|
| 122 | |
---|
| 123 | USE pegrid |
---|
[3999] | 124 | |
---|
[4444] | 125 | #if defined( __parallel ) |
---|
[3999] | 126 | USE pmc_interface, & |
---|
| 127 | ONLY: nested_run |
---|
[4444] | 128 | |
---|
[1833] | 129 | USE spectra_mod, & |
---|
[4444] | 130 | ONLY: calculate_spectra |
---|
[1833] | 131 | |
---|
[2259] | 132 | USE synthetic_turbulence_generator_mod, & |
---|
[2938] | 133 | ONLY: id_stg_left, id_stg_north, id_stg_right, id_stg_south, & |
---|
| 134 | use_syn_turb_gen |
---|
[4444] | 135 | #endif |
---|
[2259] | 136 | |
---|
[1320] | 137 | USE transpose_indices, & |
---|
[4444] | 138 | ONLY: nxl_y, nxl_z, nxr_y, nxr_z, nyn_x, nyn_z, nys_x,& |
---|
| 139 | nys_z, nzb_x, nzb_y, nzt_x, nzt_y |
---|
[667] | 140 | |
---|
[4444] | 141 | #if defined( __parallel ) |
---|
| 142 | USE transpose_indices, & |
---|
| 143 | ONLY: nxl_yd, nxr_yd, nzb_yd, nzt_yd |
---|
| 144 | #endif |
---|
[2365] | 145 | |
---|
[1] | 146 | IMPLICIT NONE |
---|
| 147 | |
---|
[3552] | 148 | INTEGER(iwp) :: i !< running index over number of processors or number of multigrid level |
---|
[4444] | 149 | #if defined( __parallel ) |
---|
[3552] | 150 | INTEGER(iwp) :: id_inflow_l !< ID indicating processors located at the left inflow boundary |
---|
[2050] | 151 | INTEGER(iwp) :: id_outflow_l !< local value of id_outflow |
---|
| 152 | INTEGER(iwp) :: id_outflow_source_l !< local value of id_outflow_source |
---|
[3552] | 153 | INTEGER(iwp) :: id_recycling_l !< ID indicating processors located at the recycling plane |
---|
[2938] | 154 | INTEGER(iwp) :: id_stg_left_l !< left lateral boundary local core id in case of turbulence generator |
---|
| 155 | INTEGER(iwp) :: id_stg_north_l !< north lateral boundary local core id in case of turbulence generator |
---|
| 156 | INTEGER(iwp) :: id_stg_right_l !< right lateral boundary local core id in case of turbulence generator |
---|
| 157 | INTEGER(iwp) :: id_stg_south_l !< south lateral boundary local core id in case of turbulence generator |
---|
[3552] | 158 | INTEGER(iwp) :: ind(5) !< array containing the subdomain bounds |
---|
[4444] | 159 | #endif |
---|
[3552] | 160 | INTEGER(iwp) :: j !< running index, used for various loops |
---|
| 161 | INTEGER(iwp) :: k !< number of vertical grid points in different multigrid level |
---|
| 162 | INTEGER(iwp) :: maximum_grid_level_l !< maximum number of grid level without switching to PE 0 |
---|
| 163 | INTEGER(iwp) :: mg_levels_x !< maximum number of grid level allowed along x-direction |
---|
| 164 | INTEGER(iwp) :: mg_levels_y !< maximum number of grid level allowed along y-direction |
---|
| 165 | INTEGER(iwp) :: mg_levels_z !< maximum number of grid level allowed along z-direction |
---|
| 166 | INTEGER(iwp) :: mg_switch_to_pe0_level_l !< maximum number of grid level with switching to PE 0 |
---|
[4444] | 167 | #if defined( __parallel ) |
---|
[3552] | 168 | INTEGER(iwp) :: nnx_y !< quotient of number of grid points along x-direction and number of PEs used along y-direction |
---|
| 169 | INTEGER(iwp) :: nny_x !< quotient of number of grid points along y-direction and number of PEs used along x-direction |
---|
| 170 | INTEGER(iwp) :: nny_z !< quotient of number of grid points along y-direction and number of PEs used along x-direction |
---|
| 171 | INTEGER(iwp) :: nnz_x !< quotient of number of grid points along z-direction and number of PEs used along x-direction |
---|
| 172 | INTEGER(iwp) :: nnz_y !< quotient of number of grid points along z-direction and number of PEs used along x-direction |
---|
| 173 | INTEGER(iwp) :: numproc_sqr !< square root of the number of processors |
---|
[4444] | 174 | #endif |
---|
[3552] | 175 | INTEGER(iwp) :: nxl_l !< lower index bound along x-direction on subdomain and different multigrid level |
---|
| 176 | INTEGER(iwp) :: nxr_l !< upper index bound along x-direction on subdomain and different multigrid level |
---|
| 177 | INTEGER(iwp) :: nyn_l !< lower index bound along y-direction on subdomain and different multigrid level |
---|
| 178 | INTEGER(iwp) :: nys_l !< upper index bound along y-direction on subdomain and different multigrid level |
---|
[4444] | 179 | #if defined( __parallel ) |
---|
[3552] | 180 | INTEGER(iwp) :: nzb_l !< lower index bound along z-direction on subdomain and different multigrid level |
---|
[4444] | 181 | #endif |
---|
[3552] | 182 | INTEGER(iwp) :: nzt_l !< upper index bound along z-direction on subdomain and different multigrid level |
---|
| 183 | !$ INTEGER(iwp) :: omp_get_num_threads !< number of OpenMP threads |
---|
[1] | 184 | |
---|
[4444] | 185 | #if defined( __parallel ) |
---|
[3552] | 186 | INTEGER(iwp), DIMENSION(:), ALLOCATABLE :: ind_all !< dummy array containing index bounds on subdomain, used for gathering |
---|
| 187 | INTEGER(iwp), DIMENSION(:), ALLOCATABLE :: nxlf !< lower index bound allong x-direction for every PE |
---|
| 188 | INTEGER(iwp), DIMENSION(:), ALLOCATABLE :: nxrf !< upper index bound allong x-direction for every PE |
---|
| 189 | INTEGER(iwp), DIMENSION(:), ALLOCATABLE :: nynf !< lower index bound allong y-direction for every PE |
---|
| 190 | INTEGER(iwp), DIMENSION(:), ALLOCATABLE :: nysf !< lower index bound allong y-direction for every PE |
---|
[1] | 191 | |
---|
[3552] | 192 | INTEGER(iwp), DIMENSION(2) :: pdims_remote !< number of PEs used for coupled model (only in atmospher-ocean coupling) |
---|
[2372] | 193 | INTEGER(iwp) :: lcoord(2) !< PE coordinates of left neighbor along x and y |
---|
| 194 | INTEGER(iwp) :: rcoord(2) !< PE coordinates of right neighbor along x and y |
---|
[4444] | 195 | #endif |
---|
[667] | 196 | |
---|
[1] | 197 | ! |
---|
| 198 | !-- Get the number of OpenMP threads |
---|
| 199 | !$OMP PARALLEL |
---|
| 200 | !$ threads_per_task = omp_get_num_threads() |
---|
| 201 | !$OMP END PARALLEL |
---|
| 202 | |
---|
| 203 | |
---|
| 204 | #if defined( __parallel ) |
---|
[667] | 205 | |
---|
[3885] | 206 | CALL location_message( 'creating virtual PE grids + MPI derived data types', 'start' ) |
---|
[1764] | 207 | |
---|
[1] | 208 | ! |
---|
[2177] | 209 | !-- Determine the processor topology or check it, if prescribed by the user |
---|
[1779] | 210 | IF ( npex == -1 .AND. npey == -1 ) THEN |
---|
[1] | 211 | |
---|
| 212 | ! |
---|
[2177] | 213 | !-- Automatic determination of the topology |
---|
[1779] | 214 | numproc_sqr = SQRT( REAL( numprocs, KIND=wp ) ) |
---|
| 215 | pdims(1) = MAX( numproc_sqr , 1 ) |
---|
[2180] | 216 | DO WHILE ( MOD( numprocs , pdims(1) ) /= 0 ) |
---|
[1779] | 217 | pdims(1) = pdims(1) - 1 |
---|
| 218 | ENDDO |
---|
[2180] | 219 | pdims(2) = numprocs / pdims(1) |
---|
[1] | 220 | |
---|
[1779] | 221 | ELSEIF ( npex /= -1 .AND. npey /= -1 ) THEN |
---|
[1] | 222 | |
---|
| 223 | ! |
---|
[1779] | 224 | !-- Prescribed by user. Number of processors on the prescribed topology |
---|
| 225 | !-- must be equal to the number of PEs available to the job |
---|
| 226 | IF ( ( npex * npey ) /= numprocs ) THEN |
---|
[3045] | 227 | WRITE( message_string, * ) 'number of PEs of the prescribed ', & |
---|
[3046] | 228 | 'topology (', npex*npey,') does not match & the number of ', & |
---|
[1779] | 229 | 'PEs available to the job (', numprocs, ')' |
---|
| 230 | CALL message( 'init_pegrid', 'PA0221', 1, 2, 0, 6, 0 ) |
---|
| 231 | ENDIF |
---|
| 232 | pdims(1) = npex |
---|
| 233 | pdims(2) = npey |
---|
[1] | 234 | |
---|
[1779] | 235 | ELSE |
---|
[1] | 236 | ! |
---|
[1779] | 237 | !-- If the processor topology is prescribed by the user, the number of |
---|
| 238 | !-- PEs must be given in both directions |
---|
[3045] | 239 | message_string = 'if the processor topology is prescribed by th' // & |
---|
[3046] | 240 | 'e user & both values of "npex" and "npey" must be given' // & |
---|
[1779] | 241 | ' in the &NAMELIST-parameter file' |
---|
| 242 | CALL message( 'init_pegrid', 'PA0222', 1, 2, 0, 6, 0 ) |
---|
[1] | 243 | |
---|
| 244 | ENDIF |
---|
| 245 | |
---|
| 246 | ! |
---|
[4461] | 247 | !-- Create four default MPI communicators for the 2d virtual PE grid. One of them will be used |
---|
| 248 | !-- as the main communicator for this run, while others might be used for specific quantities like |
---|
| 249 | !-- aerosol, chemical species, or passive scalars), if their horizontal boundary conditions shall |
---|
| 250 | !-- be different from those of the other quantities (e.g. non-cyclic conditions for aerosols, and |
---|
| 251 | !-- cyclic conditions for all others). |
---|
| 252 | DO i = 1, 4 |
---|
| 253 | |
---|
| 254 | IF ( i == 1 ) cyclic = (/ .TRUE., .TRUE. /) ! cyclic along x and y |
---|
| 255 | IF ( i == 2 ) cyclic = (/ .TRUE., .FALSE. /) ! cyclic along x |
---|
| 256 | IF ( i == 3 ) cyclic = (/ .FALSE., .TRUE. /) ! cyllic along y |
---|
| 257 | IF ( i == 4 ) cyclic = (/ .FALSE., .FALSE. /) ! non-cyclic |
---|
| 258 | |
---|
| 259 | CALL MPI_CART_CREATE( comm_palm, ndim, pdims, cyclic, reorder, & |
---|
| 260 | communicator_configurations(i)%mpi_communicator, ierr ) |
---|
| 261 | |
---|
| 262 | CALL MPI_CART_SHIFT( communicator_configurations(i)%mpi_communicator, 0, 1, & |
---|
| 263 | communicator_configurations(i)%pleft, & |
---|
| 264 | communicator_configurations(i)%pright, ierr ) |
---|
| 265 | |
---|
| 266 | CALL MPI_CART_SHIFT( communicator_configurations(i)%mpi_communicator, 1, 1, & |
---|
| 267 | communicator_configurations(i)%psouth, & |
---|
| 268 | communicator_configurations(i)%pnorth, ierr ) |
---|
| 269 | |
---|
| 270 | ENDDO |
---|
| 271 | |
---|
| 272 | ! |
---|
[1] | 273 | !-- If necessary, set horizontal boundary conditions to non-cyclic |
---|
[722] | 274 | IF ( bc_lr /= 'cyclic' ) cyclic(1) = .FALSE. |
---|
| 275 | IF ( bc_ns /= 'cyclic' ) cyclic(2) = .FALSE. |
---|
[1] | 276 | |
---|
[807] | 277 | |
---|
[1] | 278 | ! |
---|
[4461] | 279 | !-- Set the main communicator (virtual pe grid) for this run |
---|
| 280 | IF ( bc_lr == 'cyclic' .AND. bc_ns == 'cyclic' ) i = 1 |
---|
| 281 | IF ( bc_lr == 'cyclic' .AND. bc_ns /= 'cyclic' ) i = 2 |
---|
| 282 | IF ( bc_lr /= 'cyclic' .AND. bc_ns == 'cyclic' ) i = 3 |
---|
| 283 | IF ( bc_lr /= 'cyclic' .AND. bc_ns /= 'cyclic' ) i = 4 |
---|
| 284 | |
---|
| 285 | comm2d = communicator_configurations(i)%mpi_communicator |
---|
| 286 | pleft = communicator_configurations(i)%pleft |
---|
| 287 | pright = communicator_configurations(i)%pright |
---|
| 288 | psouth = communicator_configurations(i)%psouth |
---|
| 289 | pnorth = communicator_configurations(i)%pnorth |
---|
| 290 | |
---|
| 291 | ! |
---|
| 292 | !-- Set rank and coordinates of the main communicator |
---|
[1] | 293 | CALL MPI_COMM_RANK( comm2d, myid, ierr ) |
---|
[1468] | 294 | WRITE (myid_char,'(''_'',I6.6)') myid |
---|
[1] | 295 | |
---|
| 296 | CALL MPI_CART_COORDS( comm2d, myid, ndim, pcoord, ierr ) |
---|
[4461] | 297 | |
---|
[1] | 298 | ! |
---|
[2372] | 299 | !-- In case of cyclic boundary conditions, a y-shift at the boundaries in |
---|
| 300 | !-- x-direction can be introduced via parameter y_shift. The shift is done |
---|
| 301 | !-- by modifying the processor grid in such a way that processors located |
---|
| 302 | !-- at the x-boundary communicate across it to processors with y-coordinate |
---|
| 303 | !-- shifted by y_shift relative to their own. This feature can not be used |
---|
| 304 | !-- in combination with an fft pressure solver. It has been implemented to |
---|
| 305 | !-- counter the effect of streak structures in case of cyclic boundary |
---|
| 306 | !-- conditions. For a description of these see Munters |
---|
| 307 | !-- (2016; dx.doi.org/10.1063/1.4941912) |
---|
| 308 | !-- |
---|
| 309 | !-- Get coordinates of left and right neighbor on PE grid |
---|
| 310 | IF ( y_shift /= 0 ) THEN |
---|
[4301] | 311 | IF ( bc_lr == 'cyclic' ) THEN |
---|
| 312 | IF ( TRIM( psolver ) /= 'multigrid' .AND. & |
---|
| 313 | TRIM( psolver ) /= 'multigrid_noopt') & |
---|
| 314 | THEN |
---|
| 315 | message_string = 'y_shift /= 0 requires a multigrid pressure solver ' |
---|
| 316 | CALL message( 'check_parameters', 'PA0468', 1, 2, 0, 6, 0 ) |
---|
| 317 | ENDIF |
---|
[2372] | 318 | |
---|
[4301] | 319 | CALL MPI_CART_COORDS( comm2d, pright, ndim, rcoord, ierr ) |
---|
| 320 | CALL MPI_CART_COORDS( comm2d, pleft, ndim, lcoord, ierr ) |
---|
[2372] | 321 | |
---|
[4301] | 322 | ! |
---|
| 323 | !-- If the x(y)-coordinate of the right (left) neighbor is smaller (greater) |
---|
| 324 | !-- than that of the calling process, then the calling process is located on |
---|
| 325 | !-- the right (left) boundary of the processor grid. In that case, |
---|
| 326 | !-- the y-coordinate of that neighbor is increased (decreased) by y_shift. |
---|
| 327 | !-- The rank of the process with that coordinate is then inquired and the |
---|
| 328 | !-- neighbor rank for MPI_SENDRECV, pright (pleft) is set to it. |
---|
| 329 | !-- In this way, the calling process receives a new right (left) neighbor |
---|
| 330 | !-- for all future MPI_SENDRECV calls. That neighbor has a y-coordinate |
---|
| 331 | !-- of y+(-)y_shift, where y is the original right (left) neighbor's |
---|
| 332 | !-- y-coordinate. The modulo-operation ensures that if the neighbor's |
---|
| 333 | !-- y-coordinate exceeds the grid-boundary, it will be relocated to |
---|
| 334 | !-- the opposite part of the grid cyclicly. |
---|
| 335 | IF ( rcoord(1) < pcoord(1) ) THEN |
---|
| 336 | rcoord(2) = MODULO( rcoord(2) + y_shift, pdims(2) ) |
---|
| 337 | CALL MPI_CART_RANK( comm2d, rcoord, pright, ierr ) |
---|
| 338 | ENDIF |
---|
[2372] | 339 | |
---|
[4301] | 340 | IF ( lcoord(1) > pcoord(1) ) THEN |
---|
| 341 | lcoord(2) = MODULO( lcoord(2) - y_shift, pdims(2) ) |
---|
| 342 | CALL MPI_CART_RANK( comm2d, lcoord, pleft, ierr ) |
---|
| 343 | ENDIF |
---|
| 344 | |
---|
| 345 | ELSE |
---|
[2372] | 346 | ! |
---|
[4301] | 347 | !-- y-shift for non-cyclic boundary conditions is only implemented |
---|
| 348 | !-- for the turbulence recycling method in inflow_turbulence.f90 |
---|
| 349 | IF ( .NOT. turbulent_inflow ) THEN |
---|
| 350 | message_string = 'y_shift /= 0 is only allowed for cyclic ' // & |
---|
| 351 | 'boundary conditions in both directions ' // & |
---|
| 352 | 'or with turbulent_inflow == .TRUE.' |
---|
| 353 | CALL message( 'check_parameters', 'PA0467', 1, 2, 0, 6, 0 ) |
---|
| 354 | ENDIF |
---|
[2372] | 355 | ENDIF |
---|
| 356 | ENDIF |
---|
[4461] | 357 | |
---|
[2372] | 358 | ! |
---|
[1] | 359 | !-- Determine sub-topologies for transpositions |
---|
| 360 | !-- Transposition from z to x: |
---|
| 361 | remain_dims(1) = .TRUE. |
---|
| 362 | remain_dims(2) = .FALSE. |
---|
| 363 | CALL MPI_CART_SUB( comm2d, remain_dims, comm1dx, ierr ) |
---|
| 364 | CALL MPI_COMM_RANK( comm1dx, myidx, ierr ) |
---|
| 365 | ! |
---|
| 366 | !-- Transposition from x to y |
---|
| 367 | remain_dims(1) = .FALSE. |
---|
| 368 | remain_dims(2) = .TRUE. |
---|
| 369 | CALL MPI_CART_SUB( comm2d, remain_dims, comm1dy, ierr ) |
---|
| 370 | CALL MPI_COMM_RANK( comm1dy, myidy, ierr ) |
---|
| 371 | |
---|
| 372 | |
---|
| 373 | ! |
---|
[1003] | 374 | !-- Calculate array bounds along x-direction for every PE. |
---|
[3045] | 375 | ALLOCATE( nxlf(0:pdims(1)-1), nxrf(0:pdims(1)-1), nynf(0:pdims(2)-1), & |
---|
[1003] | 376 | nysf(0:pdims(2)-1) ) |
---|
[1] | 377 | |
---|
[1003] | 378 | IF ( MOD( nx+1 , pdims(1) ) /= 0 ) THEN |
---|
[3045] | 379 | WRITE( message_string, * ) 'x-direction: gridpoint number (',nx+1,') ', & |
---|
[4264] | 380 | 'is not an& integral multiple of the number',& |
---|
| 381 | ' of processors (',pdims(1),')' |
---|
[254] | 382 | CALL message( 'init_pegrid', 'PA0225', 1, 2, 0, 6, 0 ) |
---|
[1] | 383 | ELSE |
---|
[1003] | 384 | nnx = ( nx + 1 ) / pdims(1) |
---|
[4264] | 385 | ENDIF |
---|
[1] | 386 | |
---|
| 387 | ! |
---|
| 388 | !-- Left and right array bounds, number of gridpoints |
---|
| 389 | DO i = 0, pdims(1)-1 |
---|
| 390 | nxlf(i) = i * nnx |
---|
| 391 | nxrf(i) = ( i + 1 ) * nnx - 1 |
---|
| 392 | ENDDO |
---|
| 393 | |
---|
| 394 | ! |
---|
| 395 | !-- Calculate array bounds in y-direction for every PE. |
---|
[1003] | 396 | IF ( MOD( ny+1 , pdims(2) ) /= 0 ) THEN |
---|
[274] | 397 | WRITE( message_string, * ) 'y-direction: gridpoint number (',ny+1,') ', & |
---|
[4264] | 398 | 'is not an& integral multiple of the number',& |
---|
| 399 | ' of processors (',pdims(2),')' |
---|
[254] | 400 | CALL message( 'init_pegrid', 'PA0227', 1, 2, 0, 6, 0 ) |
---|
[1] | 401 | ELSE |
---|
[1003] | 402 | nny = ( ny + 1 ) / pdims(2) |
---|
[4264] | 403 | ENDIF |
---|
[1] | 404 | |
---|
| 405 | ! |
---|
| 406 | !-- South and north array bounds |
---|
| 407 | DO j = 0, pdims(2)-1 |
---|
| 408 | nysf(j) = j * nny |
---|
| 409 | nynf(j) = ( j + 1 ) * nny - 1 |
---|
| 410 | ENDDO |
---|
| 411 | |
---|
| 412 | ! |
---|
| 413 | !-- Local array bounds of the respective PEs |
---|
[1003] | 414 | nxl = nxlf(pcoord(1)) |
---|
| 415 | nxr = nxrf(pcoord(1)) |
---|
| 416 | nys = nysf(pcoord(2)) |
---|
| 417 | nyn = nynf(pcoord(2)) |
---|
| 418 | nzb = 0 |
---|
| 419 | nzt = nz |
---|
| 420 | nnz = nz |
---|
[1] | 421 | |
---|
| 422 | ! |
---|
[707] | 423 | !-- Set switches to define if the PE is situated at the border of the virtual |
---|
| 424 | !-- processor grid |
---|
| 425 | IF ( nxl == 0 ) left_border_pe = .TRUE. |
---|
| 426 | IF ( nxr == nx ) right_border_pe = .TRUE. |
---|
| 427 | IF ( nys == 0 ) south_border_pe = .TRUE. |
---|
| 428 | IF ( nyn == ny ) north_border_pe = .TRUE. |
---|
| 429 | |
---|
| 430 | ! |
---|
[1] | 431 | !-- Calculate array bounds and gridpoint numbers for the transposed arrays |
---|
| 432 | !-- (needed in the pressure solver) |
---|
| 433 | !-- For the transposed arrays, cyclic boundaries as well as top and bottom |
---|
| 434 | !-- boundaries are omitted, because they are obstructive to the transposition |
---|
| 435 | |
---|
| 436 | ! |
---|
| 437 | !-- 1. transposition z --> x |
---|
[1001] | 438 | !-- This transposition is not neccessary in case of a 1d-decomposition along x |
---|
[2938] | 439 | IF ( psolver == 'poisfft' .OR. calculate_spectra ) THEN |
---|
[1304] | 440 | |
---|
[1922] | 441 | IF ( pdims(2) /= 1 ) THEN |
---|
| 442 | IF ( MOD( nz , pdims(1) ) /= 0 ) THEN |
---|
[4264] | 443 | WRITE( message_string, * ) 'transposition z --> x:& ', & |
---|
| 444 | 'nz=',nz,' is not an integral multiple ',& |
---|
| 445 | 'of pdims(1)=',pdims(1) |
---|
[1922] | 446 | CALL message( 'init_pegrid', 'PA0230', 1, 2, 0, 6, 0 ) |
---|
| 447 | ENDIF |
---|
[1] | 448 | ENDIF |
---|
[1922] | 449 | |
---|
| 450 | nys_x = nys |
---|
| 451 | nyn_x = nyn |
---|
| 452 | nny_x = nny |
---|
| 453 | nnz_x = nz / pdims(1) |
---|
| 454 | nzb_x = 1 + myidx * nnz_x |
---|
| 455 | nzt_x = ( myidx + 1 ) * nnz_x |
---|
| 456 | sendrecvcount_zx = nnx * nny * nnz_x |
---|
| 457 | |
---|
[1] | 458 | ENDIF |
---|
| 459 | |
---|
[1922] | 460 | |
---|
[4264] | 461 | IF ( psolver == 'poisfft' ) THEN |
---|
[1] | 462 | ! |
---|
[1922] | 463 | !-- 2. transposition x --> y |
---|
| 464 | IF ( MOD( nx+1 , pdims(2) ) /= 0 ) THEN |
---|
[4264] | 465 | WRITE( message_string, * ) 'transposition x --> y:& ', & |
---|
| 466 | 'nx+1=',nx+1,' is not an integral ', & |
---|
| 467 | 'multiple of pdims(2)=',pdims(2) |
---|
[1922] | 468 | CALL message( 'init_pegrid', 'PA0231', 1, 2, 0, 6, 0 ) |
---|
| 469 | ENDIF |
---|
[1] | 470 | |
---|
[1922] | 471 | nnz_y = nnz_x |
---|
| 472 | nzb_y = nzb_x |
---|
| 473 | nzt_y = nzt_x |
---|
| 474 | nnx_y = (nx+1) / pdims(2) |
---|
| 475 | nxl_y = myidy * nnx_y |
---|
| 476 | nxr_y = ( myidy + 1 ) * nnx_y - 1 |
---|
| 477 | sendrecvcount_xy = nnx_y * nny_x * nnz_y |
---|
[1] | 478 | ! |
---|
[4264] | 479 | !-- 3. transposition y --> z |
---|
[1922] | 480 | !-- (ELSE: x --> y in case of 1D-decomposition along x) |
---|
| 481 | nxl_z = nxl_y |
---|
| 482 | nxr_z = nxr_y |
---|
| 483 | nny_z = (ny+1) / pdims(1) |
---|
| 484 | nys_z = myidx * nny_z |
---|
| 485 | nyn_z = ( myidx + 1 ) * nny_z - 1 |
---|
| 486 | sendrecvcount_yz = nnx_y * nny_z * nnz_y |
---|
[1304] | 487 | |
---|
[1922] | 488 | IF ( pdims(2) /= 1 ) THEN |
---|
[1] | 489 | ! |
---|
[1922] | 490 | !-- y --> z |
---|
| 491 | !-- This transposition is not neccessary in case of a 1d-decomposition |
---|
| 492 | !-- along x, except that the uptream-spline method is switched on |
---|
| 493 | IF ( MOD( ny+1 , pdims(1) ) /= 0 ) THEN |
---|
[4264] | 494 | WRITE( message_string, * ) 'transposition y --> z:& ', & |
---|
| 495 | 'ny+1=',ny+1,' is not an integral ', & |
---|
| 496 | 'multiple of pdims(1)=',pdims(1) |
---|
[1922] | 497 | CALL message( 'init_pegrid', 'PA0232', 1, 2, 0, 6, 0 ) |
---|
| 498 | ENDIF |
---|
[1] | 499 | |
---|
[1922] | 500 | ELSE |
---|
[1] | 501 | ! |
---|
[1922] | 502 | !-- x --> y |
---|
| 503 | !-- This condition must be fulfilled for a 1D-decomposition along x |
---|
| 504 | IF ( MOD( ny+1 , pdims(1) ) /= 0 ) THEN |
---|
[4264] | 505 | WRITE( message_string, * ) 'transposition x --> y:& ', & |
---|
| 506 | 'ny+1=',ny+1,' is not an integral ', & |
---|
| 507 | 'multiple of pdims(1)=',pdims(1) |
---|
[1922] | 508 | CALL message( 'init_pegrid', 'PA0233', 1, 2, 0, 6, 0 ) |
---|
| 509 | ENDIF |
---|
| 510 | |
---|
[1] | 511 | ENDIF |
---|
| 512 | |
---|
| 513 | ENDIF |
---|
| 514 | |
---|
| 515 | ! |
---|
| 516 | !-- Indices for direct transpositions z --> y (used for calculating spectra) |
---|
[1922] | 517 | IF ( calculate_spectra ) THEN |
---|
[1003] | 518 | IF ( MOD( nz, pdims(2) ) /= 0 ) THEN |
---|
[1922] | 519 | WRITE( message_string, * ) 'direct transposition z --> y (needed ', & |
---|
[4264] | 520 | 'for spectra):& nz=',nz,' is not an ', & |
---|
| 521 | 'integral multiple of pdims(2)=',pdims(2) |
---|
[254] | 522 | CALL message( 'init_pegrid', 'PA0234', 1, 2, 0, 6, 0 ) |
---|
[1] | 523 | ELSE |
---|
[1003] | 524 | nxl_yd = nxl |
---|
| 525 | nxr_yd = nxr |
---|
| 526 | nzb_yd = 1 + myidy * ( nz / pdims(2) ) |
---|
| 527 | nzt_yd = ( myidy + 1 ) * ( nz / pdims(2) ) |
---|
| 528 | sendrecvcount_zyd = nnx * nny * ( nz / pdims(2) ) |
---|
[1] | 529 | ENDIF |
---|
| 530 | ENDIF |
---|
| 531 | |
---|
[1922] | 532 | IF ( psolver == 'poisfft' .OR. calculate_spectra ) THEN |
---|
[1] | 533 | ! |
---|
[1922] | 534 | !-- Indices for direct transpositions y --> x |
---|
| 535 | !-- (they are only possible in case of a 1d-decomposition along x) |
---|
| 536 | IF ( pdims(2) == 1 ) THEN |
---|
| 537 | nny_x = nny / pdims(1) |
---|
| 538 | nys_x = myid * nny_x |
---|
| 539 | nyn_x = ( myid + 1 ) * nny_x - 1 |
---|
| 540 | nzb_x = 1 |
---|
| 541 | nzt_x = nz |
---|
| 542 | sendrecvcount_xy = nnx * nny_x * nz |
---|
| 543 | ENDIF |
---|
| 544 | |
---|
[1] | 545 | ENDIF |
---|
| 546 | |
---|
[1922] | 547 | IF ( psolver == 'poisfft' ) THEN |
---|
[1] | 548 | ! |
---|
[1922] | 549 | !-- Indices for direct transpositions x --> y |
---|
| 550 | !-- (they are only possible in case of a 1d-decomposition along y) |
---|
| 551 | IF ( pdims(1) == 1 ) THEN |
---|
| 552 | nnx_y = nnx / pdims(2) |
---|
| 553 | nxl_y = myid * nnx_y |
---|
| 554 | nxr_y = ( myid + 1 ) * nnx_y - 1 |
---|
| 555 | nzb_y = 1 |
---|
| 556 | nzt_y = nz |
---|
| 557 | sendrecvcount_xy = nnx_y * nny * nz |
---|
| 558 | ENDIF |
---|
| 559 | |
---|
[1] | 560 | ENDIF |
---|
| 561 | |
---|
| 562 | ! |
---|
| 563 | !-- Arrays for storing the array bounds are needed any more |
---|
| 564 | DEALLOCATE( nxlf , nxrf , nynf , nysf ) |
---|
| 565 | |
---|
[807] | 566 | |
---|
[145] | 567 | ! |
---|
| 568 | !-- Collect index bounds from other PEs (to be written to restart file later) |
---|
| 569 | ALLOCATE( hor_index_bounds(4,0:numprocs-1) ) |
---|
| 570 | |
---|
| 571 | IF ( myid == 0 ) THEN |
---|
| 572 | |
---|
| 573 | hor_index_bounds(1,0) = nxl |
---|
| 574 | hor_index_bounds(2,0) = nxr |
---|
| 575 | hor_index_bounds(3,0) = nys |
---|
| 576 | hor_index_bounds(4,0) = nyn |
---|
| 577 | |
---|
| 578 | ! |
---|
| 579 | !-- Receive data from all other PEs |
---|
| 580 | DO i = 1, numprocs-1 |
---|
| 581 | CALL MPI_RECV( ibuf, 4, MPI_INTEGER, i, MPI_ANY_TAG, comm2d, status, & |
---|
| 582 | ierr ) |
---|
| 583 | hor_index_bounds(:,i) = ibuf(1:4) |
---|
| 584 | ENDDO |
---|
| 585 | |
---|
| 586 | ELSE |
---|
| 587 | ! |
---|
| 588 | !-- Send index bounds to PE0 |
---|
| 589 | ibuf(1) = nxl |
---|
| 590 | ibuf(2) = nxr |
---|
| 591 | ibuf(3) = nys |
---|
| 592 | ibuf(4) = nyn |
---|
| 593 | CALL MPI_SEND( ibuf, 4, MPI_INTEGER, 0, myid, comm2d, ierr ) |
---|
| 594 | |
---|
| 595 | ENDIF |
---|
| 596 | |
---|
[807] | 597 | |
---|
[1] | 598 | #if defined( __print ) |
---|
| 599 | ! |
---|
| 600 | !-- Control output |
---|
| 601 | IF ( myid == 0 ) THEN |
---|
| 602 | PRINT*, '*** processor topology ***' |
---|
| 603 | PRINT*, ' ' |
---|
| 604 | PRINT*, 'myid pcoord left right south north idx idy nxl: nxr',& |
---|
| 605 | &' nys: nyn' |
---|
| 606 | PRINT*, '------------------------------------------------------------',& |
---|
| 607 | &'-----------' |
---|
| 608 | WRITE (*,1000) 0, pcoord(1), pcoord(2), pleft, pright, psouth, pnorth, & |
---|
| 609 | myidx, myidy, nxl, nxr, nys, nyn |
---|
| 610 | 1000 FORMAT (I4,2X,'(',I3,',',I3,')',3X,I4,2X,I4,3X,I4,2X,I4,2X,I3,1X,I3, & |
---|
| 611 | 2(2X,I4,':',I4)) |
---|
| 612 | |
---|
| 613 | ! |
---|
[108] | 614 | !-- Receive data from the other PEs |
---|
[1] | 615 | DO i = 1,numprocs-1 |
---|
| 616 | CALL MPI_RECV( ibuf, 12, MPI_INTEGER, i, MPI_ANY_TAG, comm2d, status, & |
---|
| 617 | ierr ) |
---|
| 618 | WRITE (*,1000) i, ( ibuf(j) , j = 1,12 ) |
---|
| 619 | ENDDO |
---|
| 620 | ELSE |
---|
| 621 | |
---|
| 622 | ! |
---|
| 623 | !-- Send data to PE0 |
---|
| 624 | ibuf(1) = pcoord(1); ibuf(2) = pcoord(2); ibuf(3) = pleft |
---|
| 625 | ibuf(4) = pright; ibuf(5) = psouth; ibuf(6) = pnorth; ibuf(7) = myidx |
---|
| 626 | ibuf(8) = myidy; ibuf(9) = nxl; ibuf(10) = nxr; ibuf(11) = nys |
---|
| 627 | ibuf(12) = nyn |
---|
| 628 | CALL MPI_SEND( ibuf, 12, MPI_INTEGER, 0, myid, comm2d, ierr ) |
---|
| 629 | ENDIF |
---|
| 630 | #endif |
---|
| 631 | |
---|
[667] | 632 | ! |
---|
[709] | 633 | !-- Determine the number of ghost point layers |
---|
[3347] | 634 | IF ( scalar_advec == 'ws-scheme' .OR. & |
---|
[3999] | 635 | momentum_advec == 'ws-scheme' .OR. nested_run ) THEN |
---|
[667] | 636 | nbgp = 3 |
---|
| 637 | ELSE |
---|
| 638 | nbgp = 1 |
---|
[709] | 639 | ENDIF |
---|
[667] | 640 | |
---|
[102] | 641 | ! |
---|
[4241] | 642 | !-- Check that the number of computational grid points is not smaller than the number of |
---|
| 643 | !-- ghost points. |
---|
| 644 | IF ( nnx < nbgp ) THEN |
---|
| 645 | WRITE( message_string, * ) 'number of subdomain grid points along x (', nnx, ') is smaller',& |
---|
| 646 | 'than the number of ghost points (', nbgp, ')' |
---|
| 647 | CALL message( 'init_pegrid', 'PA0682', 1, 2, 0, 6, 0 ) |
---|
| 648 | ENDIF |
---|
| 649 | IF ( nny < nbgp ) THEN |
---|
| 650 | WRITE( message_string, * ) 'number of subdomain grid points along y (', nny, ') is smaller',& |
---|
| 651 | 'than the number of ghost points (', nbgp, ')' |
---|
| 652 | CALL message( 'init_pegrid', 'PA0683', 1, 2, 0, 6, 0 ) |
---|
| 653 | ENDIF |
---|
| 654 | |
---|
| 655 | ! |
---|
[709] | 656 | !-- Create a new MPI derived datatype for the exchange of surface (xy) data, |
---|
| 657 | !-- which is needed for coupled atmosphere-ocean runs. |
---|
| 658 | !-- First, calculate number of grid points of an xy-plane. |
---|
[667] | 659 | ngp_xy = ( nxr - nxl + 1 + 2 * nbgp ) * ( nyn - nys + 1 + 2 * nbgp ) |
---|
[102] | 660 | CALL MPI_TYPE_VECTOR( ngp_xy, 1, nzt-nzb+2, MPI_REAL, type_xy, ierr ) |
---|
| 661 | CALL MPI_TYPE_COMMIT( type_xy, ierr ) |
---|
[667] | 662 | |
---|
[4564] | 663 | IF ( TRIM( coupling_mode ) /= 'uncoupled' ) THEN |
---|
[667] | 664 | |
---|
| 665 | ! |
---|
| 666 | !-- Pass the number of grid points of the atmosphere model to |
---|
| 667 | !-- the ocean model and vice versa |
---|
| 668 | IF ( coupling_mode == 'atmosphere_to_ocean' ) THEN |
---|
| 669 | |
---|
| 670 | nx_a = nx |
---|
| 671 | ny_a = ny |
---|
| 672 | |
---|
[709] | 673 | IF ( myid == 0 ) THEN |
---|
| 674 | |
---|
| 675 | CALL MPI_SEND( nx_a, 1, MPI_INTEGER, numprocs, 1, comm_inter, & |
---|
| 676 | ierr ) |
---|
| 677 | CALL MPI_SEND( ny_a, 1, MPI_INTEGER, numprocs, 2, comm_inter, & |
---|
| 678 | ierr ) |
---|
| 679 | CALL MPI_SEND( pdims, 2, MPI_INTEGER, numprocs, 3, comm_inter, & |
---|
| 680 | ierr ) |
---|
| 681 | CALL MPI_RECV( nx_o, 1, MPI_INTEGER, numprocs, 4, comm_inter, & |
---|
| 682 | status, ierr ) |
---|
| 683 | CALL MPI_RECV( ny_o, 1, MPI_INTEGER, numprocs, 5, comm_inter, & |
---|
| 684 | status, ierr ) |
---|
| 685 | CALL MPI_RECV( pdims_remote, 2, MPI_INTEGER, numprocs, 6, & |
---|
[667] | 686 | comm_inter, status, ierr ) |
---|
| 687 | ENDIF |
---|
| 688 | |
---|
[709] | 689 | CALL MPI_BCAST( nx_o, 1, MPI_INTEGER, 0, comm2d, ierr ) |
---|
| 690 | CALL MPI_BCAST( ny_o, 1, MPI_INTEGER, 0, comm2d, ierr ) |
---|
| 691 | CALL MPI_BCAST( pdims_remote, 2, MPI_INTEGER, 0, comm2d, ierr ) |
---|
[667] | 692 | |
---|
| 693 | ELSEIF ( coupling_mode == 'ocean_to_atmosphere' ) THEN |
---|
| 694 | |
---|
| 695 | nx_o = nx |
---|
| 696 | ny_o = ny |
---|
| 697 | |
---|
| 698 | IF ( myid == 0 ) THEN |
---|
[709] | 699 | |
---|
| 700 | CALL MPI_RECV( nx_a, 1, MPI_INTEGER, 0, 1, comm_inter, status, & |
---|
| 701 | ierr ) |
---|
| 702 | CALL MPI_RECV( ny_a, 1, MPI_INTEGER, 0, 2, comm_inter, status, & |
---|
| 703 | ierr ) |
---|
| 704 | CALL MPI_RECV( pdims_remote, 2, MPI_INTEGER, 0, 3, comm_inter, & |
---|
| 705 | status, ierr ) |
---|
| 706 | CALL MPI_SEND( nx_o, 1, MPI_INTEGER, 0, 4, comm_inter, ierr ) |
---|
| 707 | CALL MPI_SEND( ny_o, 1, MPI_INTEGER, 0, 5, comm_inter, ierr ) |
---|
| 708 | CALL MPI_SEND( pdims, 2, MPI_INTEGER, 0, 6, comm_inter, ierr ) |
---|
[667] | 709 | ENDIF |
---|
| 710 | |
---|
| 711 | CALL MPI_BCAST( nx_a, 1, MPI_INTEGER, 0, comm2d, ierr) |
---|
| 712 | CALL MPI_BCAST( ny_a, 1, MPI_INTEGER, 0, comm2d, ierr) |
---|
| 713 | CALL MPI_BCAST( pdims_remote, 2, MPI_INTEGER, 0, comm2d, ierr) |
---|
| 714 | |
---|
| 715 | ENDIF |
---|
| 716 | |
---|
[709] | 717 | ngp_a = ( nx_a+1 + 2 * nbgp ) * ( ny_a+1 + 2 * nbgp ) |
---|
| 718 | ngp_o = ( nx_o+1 + 2 * nbgp ) * ( ny_o+1 + 2 * nbgp ) |
---|
[667] | 719 | |
---|
| 720 | ! |
---|
[709] | 721 | !-- Determine if the horizontal grid and the number of PEs in ocean and |
---|
| 722 | !-- atmosphere is same or not |
---|
| 723 | IF ( nx_o == nx_a .AND. ny_o == ny_a .AND. & |
---|
[667] | 724 | pdims(1) == pdims_remote(1) .AND. pdims(2) == pdims_remote(2) ) & |
---|
| 725 | THEN |
---|
| 726 | coupling_topology = 0 |
---|
| 727 | ELSE |
---|
| 728 | coupling_topology = 1 |
---|
| 729 | ENDIF |
---|
| 730 | |
---|
| 731 | ! |
---|
| 732 | !-- Determine the target PEs for the exchange between ocean and |
---|
| 733 | !-- atmosphere (comm2d) |
---|
[709] | 734 | IF ( coupling_topology == 0 ) THEN |
---|
| 735 | ! |
---|
| 736 | !-- In case of identical topologies, every atmosphere PE has exactly one |
---|
| 737 | !-- ocean PE counterpart and vice versa |
---|
| 738 | IF ( TRIM( coupling_mode ) == 'atmosphere_to_ocean' ) THEN |
---|
[667] | 739 | target_id = myid + numprocs |
---|
| 740 | ELSE |
---|
| 741 | target_id = myid |
---|
| 742 | ENDIF |
---|
| 743 | |
---|
| 744 | ELSE |
---|
| 745 | ! |
---|
| 746 | !-- In case of nonequivalent topology in ocean and atmosphere only for |
---|
| 747 | !-- PE0 in ocean and PE0 in atmosphere a target_id is needed, since |
---|
[709] | 748 | !-- data echxchange between ocean and atmosphere will be done only |
---|
| 749 | !-- between these PEs. |
---|
| 750 | IF ( myid == 0 ) THEN |
---|
| 751 | |
---|
| 752 | IF ( TRIM( coupling_mode ) == 'atmosphere_to_ocean' ) THEN |
---|
[667] | 753 | target_id = numprocs |
---|
| 754 | ELSE |
---|
| 755 | target_id = 0 |
---|
| 756 | ENDIF |
---|
[709] | 757 | |
---|
[667] | 758 | ENDIF |
---|
[709] | 759 | |
---|
[667] | 760 | ENDIF |
---|
| 761 | |
---|
| 762 | ENDIF |
---|
| 763 | |
---|
[1] | 764 | #else |
---|
| 765 | |
---|
| 766 | ! |
---|
| 767 | !-- Array bounds when running on a single PE (respectively a non-parallel |
---|
| 768 | !-- machine) |
---|
[1003] | 769 | nxl = 0 |
---|
| 770 | nxr = nx |
---|
| 771 | nnx = nxr - nxl + 1 |
---|
| 772 | nys = 0 |
---|
| 773 | nyn = ny |
---|
| 774 | nny = nyn - nys + 1 |
---|
| 775 | nzb = 0 |
---|
| 776 | nzt = nz |
---|
| 777 | nnz = nz |
---|
[1] | 778 | |
---|
[145] | 779 | ALLOCATE( hor_index_bounds(4,0:0) ) |
---|
| 780 | hor_index_bounds(1,0) = nxl |
---|
| 781 | hor_index_bounds(2,0) = nxr |
---|
| 782 | hor_index_bounds(3,0) = nys |
---|
| 783 | hor_index_bounds(4,0) = nyn |
---|
| 784 | |
---|
[1] | 785 | ! |
---|
| 786 | !-- Array bounds for the pressure solver (in the parallel code, these bounds |
---|
| 787 | !-- are the ones for the transposed arrays) |
---|
[1003] | 788 | nys_x = nys |
---|
| 789 | nyn_x = nyn |
---|
| 790 | nzb_x = nzb + 1 |
---|
| 791 | nzt_x = nzt |
---|
[1] | 792 | |
---|
[1003] | 793 | nxl_y = nxl |
---|
| 794 | nxr_y = nxr |
---|
| 795 | nzb_y = nzb + 1 |
---|
| 796 | nzt_y = nzt |
---|
[1] | 797 | |
---|
[1003] | 798 | nxl_z = nxl |
---|
| 799 | nxr_z = nxr |
---|
| 800 | nys_z = nys |
---|
| 801 | nyn_z = nyn |
---|
[1] | 802 | |
---|
| 803 | #endif |
---|
| 804 | |
---|
| 805 | ! |
---|
| 806 | !-- Calculate number of grid levels necessary for the multigrid poisson solver |
---|
| 807 | !-- as well as the gridpoint indices on each level |
---|
[1575] | 808 | IF ( psolver(1:9) == 'multigrid' ) THEN |
---|
[1] | 809 | |
---|
| 810 | ! |
---|
| 811 | !-- First calculate number of possible grid levels for the subdomains |
---|
| 812 | mg_levels_x = 1 |
---|
| 813 | mg_levels_y = 1 |
---|
| 814 | mg_levels_z = 1 |
---|
| 815 | |
---|
| 816 | i = nnx |
---|
| 817 | DO WHILE ( MOD( i, 2 ) == 0 .AND. i /= 2 ) |
---|
| 818 | i = i / 2 |
---|
| 819 | mg_levels_x = mg_levels_x + 1 |
---|
| 820 | ENDDO |
---|
| 821 | |
---|
| 822 | j = nny |
---|
| 823 | DO WHILE ( MOD( j, 2 ) == 0 .AND. j /= 2 ) |
---|
| 824 | j = j / 2 |
---|
| 825 | mg_levels_y = mg_levels_y + 1 |
---|
| 826 | ENDDO |
---|
| 827 | |
---|
[181] | 828 | k = nz ! do not use nnz because it might be > nz due to transposition |
---|
| 829 | ! requirements |
---|
[1] | 830 | DO WHILE ( MOD( k, 2 ) == 0 .AND. k /= 2 ) |
---|
| 831 | k = k / 2 |
---|
| 832 | mg_levels_z = mg_levels_z + 1 |
---|
| 833 | ENDDO |
---|
[2197] | 834 | ! |
---|
| 835 | !-- The optimized MG-solver does not allow odd values for nz at the coarsest |
---|
| 836 | !-- grid level |
---|
| 837 | IF ( TRIM( psolver ) /= 'multigrid_noopt' ) THEN |
---|
| 838 | IF ( MOD( k, 2 ) /= 0 ) mg_levels_z = mg_levels_z - 1 |
---|
[3057] | 839 | ! |
---|
| 840 | !-- An odd value of nz does not work. The finest level must have an even |
---|
| 841 | !-- value. |
---|
| 842 | IF ( mg_levels_z == 0 ) THEN |
---|
| 843 | message_string = 'optimized multigrid method requires nz to be even' |
---|
[3058] | 844 | CALL message( 'init_pegrid', 'PA0495', 1, 2, 0, 6, 0 ) |
---|
[3057] | 845 | ENDIF |
---|
[2197] | 846 | ENDIF |
---|
[1] | 847 | |
---|
| 848 | maximum_grid_level = MIN( mg_levels_x, mg_levels_y, mg_levels_z ) |
---|
[3890] | 849 | ! |
---|
| 850 | !-- Check if subdomain sizes prevents any coarsening. |
---|
| 851 | !-- This case, the maximum number of grid levels is 1, i.e. effectively |
---|
| 852 | !-- a Gauss-Seidel scheme is applied rather than a multigrid approach. |
---|
[3897] | 853 | !-- Give a warning in this case. |
---|
| 854 | IF ( maximum_grid_level == 1 .AND. mg_switch_to_pe0_level == -1 ) THEN |
---|
| 855 | message_string = 'No grid coarsening possible, multigrid ' // & |
---|
| 856 | 'approach effectively reduces to a Gauss-Seidel ' //& |
---|
| 857 | 'scheme.' |
---|
| 858 | |
---|
| 859 | CALL message( 'poismg', 'PA0648', 0, 1, 0, 6, 0 ) |
---|
[3890] | 860 | ENDIF |
---|
[1] | 861 | |
---|
| 862 | ! |
---|
| 863 | !-- Find out, if the total domain allows more levels. These additional |
---|
[709] | 864 | !-- levels are identically processed on all PEs. |
---|
[197] | 865 | IF ( numprocs > 1 .AND. mg_switch_to_pe0_level /= -1 ) THEN |
---|
[709] | 866 | |
---|
[1] | 867 | IF ( mg_levels_z > MIN( mg_levels_x, mg_levels_y ) ) THEN |
---|
[709] | 868 | |
---|
[1] | 869 | mg_switch_to_pe0_level_l = maximum_grid_level |
---|
| 870 | |
---|
| 871 | mg_levels_x = 1 |
---|
| 872 | mg_levels_y = 1 |
---|
| 873 | |
---|
| 874 | i = nx+1 |
---|
| 875 | DO WHILE ( MOD( i, 2 ) == 0 .AND. i /= 2 ) |
---|
| 876 | i = i / 2 |
---|
| 877 | mg_levels_x = mg_levels_x + 1 |
---|
| 878 | ENDDO |
---|
| 879 | |
---|
| 880 | j = ny+1 |
---|
| 881 | DO WHILE ( MOD( j, 2 ) == 0 .AND. j /= 2 ) |
---|
| 882 | j = j / 2 |
---|
| 883 | mg_levels_y = mg_levels_y + 1 |
---|
| 884 | ENDDO |
---|
| 885 | |
---|
| 886 | maximum_grid_level_l = MIN( mg_levels_x, mg_levels_y, mg_levels_z ) |
---|
| 887 | |
---|
| 888 | IF ( maximum_grid_level_l > mg_switch_to_pe0_level_l ) THEN |
---|
| 889 | mg_switch_to_pe0_level_l = maximum_grid_level_l - & |
---|
| 890 | mg_switch_to_pe0_level_l + 1 |
---|
| 891 | ELSE |
---|
| 892 | mg_switch_to_pe0_level_l = 0 |
---|
| 893 | ENDIF |
---|
[709] | 894 | |
---|
[1] | 895 | ELSE |
---|
[3057] | 896 | |
---|
[1] | 897 | mg_switch_to_pe0_level_l = 0 |
---|
| 898 | maximum_grid_level_l = maximum_grid_level |
---|
[709] | 899 | |
---|
[1] | 900 | ENDIF |
---|
| 901 | |
---|
| 902 | ! |
---|
| 903 | !-- Use switch level calculated above only if it is not pre-defined |
---|
| 904 | !-- by user |
---|
| 905 | IF ( mg_switch_to_pe0_level == 0 ) THEN |
---|
| 906 | IF ( mg_switch_to_pe0_level_l /= 0 ) THEN |
---|
| 907 | mg_switch_to_pe0_level = mg_switch_to_pe0_level_l |
---|
| 908 | maximum_grid_level = maximum_grid_level_l |
---|
| 909 | ENDIF |
---|
| 910 | |
---|
| 911 | ELSE |
---|
| 912 | ! |
---|
| 913 | !-- Check pre-defined value and reset to default, if neccessary |
---|
[3045] | 914 | IF ( mg_switch_to_pe0_level < mg_switch_to_pe0_level_l .OR. & |
---|
[1] | 915 | mg_switch_to_pe0_level >= maximum_grid_level_l ) THEN |
---|
[3045] | 916 | message_string = 'mg_switch_to_pe0_level ' // & |
---|
[2271] | 917 | 'out of range and reset to 0' |
---|
[254] | 918 | CALL message( 'init_pegrid', 'PA0235', 0, 1, 0, 6, 0 ) |
---|
[1] | 919 | mg_switch_to_pe0_level = 0 |
---|
| 920 | ELSE |
---|
| 921 | ! |
---|
| 922 | !-- Use the largest number of possible levels anyway and recalculate |
---|
| 923 | !-- the switch level to this largest number of possible values |
---|
| 924 | maximum_grid_level = maximum_grid_level_l |
---|
| 925 | |
---|
| 926 | ENDIF |
---|
[709] | 927 | |
---|
[1] | 928 | ENDIF |
---|
| 929 | |
---|
| 930 | ENDIF |
---|
| 931 | |
---|
[1056] | 932 | ALLOCATE( grid_level_count(maximum_grid_level), & |
---|
| 933 | nxl_mg(0:maximum_grid_level), nxr_mg(0:maximum_grid_level), & |
---|
| 934 | nyn_mg(0:maximum_grid_level), nys_mg(0:maximum_grid_level), & |
---|
| 935 | nzt_mg(0:maximum_grid_level) ) |
---|
[1] | 936 | |
---|
| 937 | grid_level_count = 0 |
---|
[1056] | 938 | ! |
---|
| 939 | !-- Index zero required as dummy due to definition of arrays f2 and p2 in |
---|
| 940 | !-- recursive subroutine next_mg_level |
---|
| 941 | nxl_mg(0) = 0; nxr_mg(0) = 0; nyn_mg(0) = 0; nys_mg(0) = 0; nzt_mg(0) = 0 |
---|
[778] | 942 | |
---|
[1] | 943 | nxl_l = nxl; nxr_l = nxr; nys_l = nys; nyn_l = nyn; nzt_l = nzt |
---|
| 944 | |
---|
| 945 | DO i = maximum_grid_level, 1 , -1 |
---|
| 946 | |
---|
| 947 | IF ( i == mg_switch_to_pe0_level ) THEN |
---|
[1804] | 948 | #if defined( __parallel ) |
---|
[1] | 949 | ! |
---|
| 950 | !-- Save the grid size of the subdomain at the switch level, because |
---|
| 951 | !-- it is needed in poismg. |
---|
| 952 | ind(1) = nxl_l; ind(2) = nxr_l |
---|
| 953 | ind(3) = nys_l; ind(4) = nyn_l |
---|
| 954 | ind(5) = nzt_l |
---|
| 955 | ALLOCATE( ind_all(5*numprocs), mg_loc_ind(5,0:numprocs-1) ) |
---|
| 956 | CALL MPI_ALLGATHER( ind, 5, MPI_INTEGER, ind_all, 5, & |
---|
| 957 | MPI_INTEGER, comm2d, ierr ) |
---|
| 958 | DO j = 0, numprocs-1 |
---|
| 959 | DO k = 1, 5 |
---|
| 960 | mg_loc_ind(k,j) = ind_all(k+j*5) |
---|
| 961 | ENDDO |
---|
| 962 | ENDDO |
---|
| 963 | DEALLOCATE( ind_all ) |
---|
| 964 | ! |
---|
[709] | 965 | !-- Calculate the grid size of the total domain |
---|
[1] | 966 | nxr_l = ( nxr_l-nxl_l+1 ) * pdims(1) - 1 |
---|
| 967 | nxl_l = 0 |
---|
| 968 | nyn_l = ( nyn_l-nys_l+1 ) * pdims(2) - 1 |
---|
| 969 | nys_l = 0 |
---|
| 970 | ! |
---|
| 971 | !-- The size of this gathered array must not be larger than the |
---|
| 972 | !-- array tend, which is used in the multigrid scheme as a temporary |
---|
[778] | 973 | !-- array. Therefore the subdomain size of an PE is calculated and |
---|
| 974 | !-- the size of the gathered grid. These values are used in |
---|
| 975 | !-- routines pres and poismg |
---|
| 976 | subdomain_size = ( nxr - nxl + 2 * nbgp + 1 ) * & |
---|
| 977 | ( nyn - nys + 2 * nbgp + 1 ) * ( nzt - nzb + 2 ) |
---|
[3045] | 978 | gathered_size = ( nxr_l - nxl_l + 3 ) * ( nyn_l - nys_l + 3 ) * & |
---|
[1] | 979 | ( nzt_l - nzb + 2 ) |
---|
| 980 | |
---|
[1804] | 981 | #else |
---|
[3045] | 982 | message_string = 'multigrid gather/scatter impossible ' // & |
---|
[1] | 983 | 'in non parallel mode' |
---|
[254] | 984 | CALL message( 'init_pegrid', 'PA0237', 1, 2, 0, 6, 0 ) |
---|
[1] | 985 | #endif |
---|
| 986 | ENDIF |
---|
| 987 | |
---|
| 988 | nxl_mg(i) = nxl_l |
---|
| 989 | nxr_mg(i) = nxr_l |
---|
| 990 | nys_mg(i) = nys_l |
---|
| 991 | nyn_mg(i) = nyn_l |
---|
| 992 | nzt_mg(i) = nzt_l |
---|
| 993 | |
---|
| 994 | nxl_l = nxl_l / 2 |
---|
| 995 | nxr_l = nxr_l / 2 |
---|
| 996 | nys_l = nys_l / 2 |
---|
| 997 | nyn_l = nyn_l / 2 |
---|
| 998 | nzt_l = nzt_l / 2 |
---|
[3890] | 999 | |
---|
[1] | 1000 | ENDDO |
---|
| 1001 | |
---|
[780] | 1002 | ! |
---|
[3045] | 1003 | !-- Temporary problem: Currently calculation of maxerror in routine poismg crashes |
---|
[780] | 1004 | !-- if grid data are collected on PE0 already on the finest grid level. |
---|
| 1005 | !-- To be solved later. |
---|
| 1006 | IF ( maximum_grid_level == mg_switch_to_pe0_level ) THEN |
---|
| 1007 | message_string = 'grid coarsening on subdomain level cannot be performed' |
---|
| 1008 | CALL message( 'poismg', 'PA0236', 1, 2, 0, 6, 0 ) |
---|
| 1009 | ENDIF |
---|
| 1010 | |
---|
[1] | 1011 | ELSE |
---|
| 1012 | |
---|
[667] | 1013 | maximum_grid_level = 0 |
---|
[1] | 1014 | |
---|
| 1015 | ENDIF |
---|
| 1016 | |
---|
[722] | 1017 | ! |
---|
| 1018 | !-- Default level 0 tells exchange_horiz that all ghost planes have to be |
---|
| 1019 | !-- exchanged. grid_level is adjusted in poismg, where only one ghost plane |
---|
| 1020 | !-- is required. |
---|
| 1021 | grid_level = 0 |
---|
[1] | 1022 | |
---|
[1804] | 1023 | #if defined( __parallel ) |
---|
[1] | 1024 | ! |
---|
| 1025 | !-- Gridpoint number for the exchange of ghost points (y-line for 2D-arrays) |
---|
[667] | 1026 | ngp_y = nyn - nys + 1 + 2 * nbgp |
---|
[1] | 1027 | |
---|
| 1028 | ! |
---|
[709] | 1029 | !-- Define new MPI derived datatypes for the exchange of ghost points in |
---|
| 1030 | !-- x- and y-direction for 2D-arrays (line) |
---|
[1968] | 1031 | CALL MPI_TYPE_VECTOR( nxr-nxl+1+2*nbgp, nbgp, ngp_y, MPI_REAL, type_x, & |
---|
[709] | 1032 | ierr ) |
---|
[1] | 1033 | CALL MPI_TYPE_COMMIT( type_x, ierr ) |
---|
| 1034 | |
---|
[667] | 1035 | CALL MPI_TYPE_VECTOR( nbgp, ngp_y, ngp_y, MPI_REAL, type_y, ierr ) |
---|
| 1036 | CALL MPI_TYPE_COMMIT( type_y, ierr ) |
---|
[1968] | 1037 | ! |
---|
| 1038 | !-- Define new MPI derived datatypes for the exchange of ghost points in |
---|
[3542] | 1039 | !-- x- and y-direction for 2D-INTEGER arrays (line) - on normal grid. |
---|
| 1040 | !-- Define types for 32-bit and 8-bit Integer. The 8-bit Integer are only |
---|
| 1041 | !-- required on normal grid, while 32-bit Integer may be also required on |
---|
| 1042 | !-- coarser grid level in case of multigrid solver. |
---|
| 1043 | ! |
---|
| 1044 | !-- 8-bit Integer |
---|
| 1045 | CALL MPI_TYPE_VECTOR( nxr-nxl+1+2*nbgp, nbgp, ngp_y, MPI_BYTE, & |
---|
| 1046 | type_x_byte, ierr ) |
---|
| 1047 | CALL MPI_TYPE_COMMIT( type_x_byte, ierr ) |
---|
| 1048 | |
---|
| 1049 | CALL MPI_TYPE_VECTOR( nbgp, ngp_y, ngp_y, MPI_BYTE, & |
---|
| 1050 | type_y_byte, ierr ) |
---|
| 1051 | CALL MPI_TYPE_COMMIT( type_y_byte, ierr ) |
---|
| 1052 | ! |
---|
| 1053 | !-- 32-bit Integer |
---|
[1968] | 1054 | ALLOCATE( type_x_int(0:maximum_grid_level), & |
---|
| 1055 | type_y_int(0:maximum_grid_level) ) |
---|
[3542] | 1056 | |
---|
[1968] | 1057 | CALL MPI_TYPE_VECTOR( nxr-nxl+1+2*nbgp, nbgp, ngp_y, MPI_INTEGER, & |
---|
| 1058 | type_x_int(0), ierr ) |
---|
| 1059 | CALL MPI_TYPE_COMMIT( type_x_int(0), ierr ) |
---|
[667] | 1060 | |
---|
[1968] | 1061 | CALL MPI_TYPE_VECTOR( nbgp, ngp_y, ngp_y, MPI_INTEGER, type_y_int(0), ierr ) |
---|
| 1062 | CALL MPI_TYPE_COMMIT( type_y_int(0), ierr ) |
---|
[1] | 1063 | ! |
---|
| 1064 | !-- Calculate gridpoint numbers for the exchange of ghost points along x |
---|
| 1065 | !-- (yz-plane for 3D-arrays) and define MPI derived data type(s) for the |
---|
| 1066 | !-- exchange of ghost points in y-direction (xz-plane). |
---|
| 1067 | !-- Do these calculations for the model grid and (if necessary) also |
---|
| 1068 | !-- for the coarser grid levels used in the multigrid method |
---|
[2696] | 1069 | ALLOCATE ( ngp_xz(0:maximum_grid_level), & |
---|
| 1070 | ngp_xz_int(0:maximum_grid_level), & |
---|
| 1071 | ngp_yz(0:maximum_grid_level), & |
---|
| 1072 | ngp_yz_int(0:maximum_grid_level), & |
---|
| 1073 | type_xz(0:maximum_grid_level), & |
---|
| 1074 | type_xz_int(0:maximum_grid_level), & |
---|
| 1075 | type_yz(0:maximum_grid_level), & |
---|
| 1076 | type_yz_int(0:maximum_grid_level) ) |
---|
[1] | 1077 | |
---|
| 1078 | nxl_l = nxl; nxr_l = nxr; nys_l = nys; nyn_l = nyn; nzb_l = nzb; nzt_l = nzt |
---|
[709] | 1079 | |
---|
[667] | 1080 | ! |
---|
| 1081 | !-- Discern between the model grid, which needs nbgp ghost points and |
---|
| 1082 | !-- grid levels for the multigrid scheme. In the latter case only one |
---|
| 1083 | !-- ghost point is necessary. |
---|
[709] | 1084 | !-- First definition of MPI-datatypes for exchange of ghost layers on normal |
---|
[667] | 1085 | !-- grid. The following loop is needed for data exchange in poismg.f90. |
---|
| 1086 | ! |
---|
| 1087 | !-- Determine number of grid points of yz-layer for exchange |
---|
| 1088 | ngp_yz(0) = (nzt - nzb + 2) * (nyn - nys + 1 + 2 * nbgp) |
---|
[709] | 1089 | |
---|
[667] | 1090 | ! |
---|
[709] | 1091 | !-- Define an MPI-datatype for the exchange of left/right boundaries. |
---|
| 1092 | !-- Although data are contiguous in physical memory (which does not |
---|
| 1093 | !-- necessarily require an MPI-derived datatype), the data exchange between |
---|
| 1094 | !-- left and right PE's using the MPI-derived type is 10% faster than without. |
---|
[667] | 1095 | CALL MPI_TYPE_VECTOR( nxr-nxl+1+2*nbgp, nbgp*(nzt-nzb+2), ngp_yz(0), & |
---|
[709] | 1096 | MPI_REAL, type_xz(0), ierr ) |
---|
[667] | 1097 | CALL MPI_TYPE_COMMIT( type_xz(0), ierr ) |
---|
[1] | 1098 | |
---|
[709] | 1099 | CALL MPI_TYPE_VECTOR( nbgp, ngp_yz(0), ngp_yz(0), MPI_REAL, type_yz(0), & |
---|
| 1100 | ierr ) |
---|
[667] | 1101 | CALL MPI_TYPE_COMMIT( type_yz(0), ierr ) |
---|
[709] | 1102 | |
---|
[667] | 1103 | ! |
---|
[2696] | 1104 | !-- Define data types for exchange of 3D Integer arrays. |
---|
| 1105 | ngp_yz_int(0) = (nzt - nzb + 2) * (nyn - nys + 1 + 2 * nbgp) |
---|
| 1106 | |
---|
| 1107 | CALL MPI_TYPE_VECTOR( nxr-nxl+1+2*nbgp, nbgp*(nzt-nzb+2), ngp_yz_int(0), & |
---|
| 1108 | MPI_INTEGER, type_xz_int(0), ierr ) |
---|
| 1109 | CALL MPI_TYPE_COMMIT( type_xz_int(0), ierr ) |
---|
| 1110 | |
---|
| 1111 | CALL MPI_TYPE_VECTOR( nbgp, ngp_yz_int(0), ngp_yz_int(0), MPI_INTEGER, & |
---|
| 1112 | type_yz_int(0), ierr ) |
---|
| 1113 | CALL MPI_TYPE_COMMIT( type_yz_int(0), ierr ) |
---|
| 1114 | |
---|
| 1115 | ! |
---|
[709] | 1116 | !-- Definition of MPI-datatypes for multigrid method (coarser level grids) |
---|
[1575] | 1117 | IF ( psolver(1:9) == 'multigrid' ) THEN |
---|
[667] | 1118 | ! |
---|
[709] | 1119 | !-- Definition of MPI-datatyoe as above, but only 1 ghost level is used |
---|
| 1120 | DO i = maximum_grid_level, 1 , -1 |
---|
[1968] | 1121 | ! |
---|
[2696] | 1122 | !-- For 3D-exchange on different multigrid level, one ghost point for |
---|
| 1123 | !-- REAL arrays, two ghost points for INTEGER arrays |
---|
[1575] | 1124 | ngp_xz(i) = (nzt_l - nzb_l + 2) * (nxr_l - nxl_l + 3) |
---|
[667] | 1125 | ngp_yz(i) = (nzt_l - nzb_l + 2) * (nyn_l - nys_l + 3) |
---|
| 1126 | |
---|
[2696] | 1127 | ngp_xz_int(i) = (nzt_l - nzb_l + 2) * (nxr_l - nxl_l + 3) |
---|
| 1128 | ngp_yz_int(i) = (nzt_l - nzb_l + 2) * (nyn_l - nys_l + 3) |
---|
| 1129 | ! |
---|
| 1130 | !-- MPI data type for REAL arrays, for xz-layers |
---|
| 1131 | CALL MPI_TYPE_VECTOR( nxr_l-nxl_l+3, nzt_l-nzb_l+2, ngp_yz(i), & |
---|
[709] | 1132 | MPI_REAL, type_xz(i), ierr ) |
---|
[667] | 1133 | CALL MPI_TYPE_COMMIT( type_xz(i), ierr ) |
---|
[1] | 1134 | |
---|
[2696] | 1135 | ! |
---|
| 1136 | !-- MPI data type for INTEGER arrays, for xz-layers |
---|
| 1137 | CALL MPI_TYPE_VECTOR( nxr_l-nxl_l+3, nzt_l-nzb_l+2, ngp_yz_int(i), & |
---|
| 1138 | MPI_INTEGER, type_xz_int(i), ierr ) |
---|
| 1139 | CALL MPI_TYPE_COMMIT( type_xz_int(i), ierr ) |
---|
| 1140 | |
---|
| 1141 | ! |
---|
| 1142 | !-- MPI data type for REAL arrays, for yz-layers |
---|
[709] | 1143 | CALL MPI_TYPE_VECTOR( 1, ngp_yz(i), ngp_yz(i), MPI_REAL, type_yz(i), & |
---|
| 1144 | ierr ) |
---|
[667] | 1145 | CALL MPI_TYPE_COMMIT( type_yz(i), ierr ) |
---|
[2696] | 1146 | ! |
---|
| 1147 | !-- MPI data type for INTEGER arrays, for yz-layers |
---|
| 1148 | CALL MPI_TYPE_VECTOR( 1, ngp_yz_int(i), ngp_yz_int(i), MPI_INTEGER, & |
---|
| 1149 | type_yz_int(i), ierr ) |
---|
| 1150 | CALL MPI_TYPE_COMMIT( type_yz_int(i), ierr ) |
---|
[667] | 1151 | |
---|
[1968] | 1152 | |
---|
| 1153 | !-- For 2D-exchange of INTEGER arrays on coarser grid level, where 2 ghost |
---|
[3542] | 1154 | !-- points need to be exchanged. Only required for 32-bit Integer arrays. |
---|
| 1155 | CALL MPI_TYPE_VECTOR( nxr_l-nxl_l+5, 2, nyn_l-nys_l+5, MPI_INTEGER, & |
---|
[1968] | 1156 | type_x_int(i), ierr ) |
---|
| 1157 | CALL MPI_TYPE_COMMIT( type_x_int(i), ierr ) |
---|
| 1158 | |
---|
| 1159 | |
---|
[3542] | 1160 | CALL MPI_TYPE_VECTOR( 2, nyn_l-nys_l+5, nyn_l-nys_l+5, MPI_INTEGER, & |
---|
[1968] | 1161 | type_y_int(i), ierr ) |
---|
| 1162 | CALL MPI_TYPE_COMMIT( type_y_int(i), ierr ) |
---|
| 1163 | |
---|
[667] | 1164 | nxl_l = nxl_l / 2 |
---|
| 1165 | nxr_l = nxr_l / 2 |
---|
| 1166 | nys_l = nys_l / 2 |
---|
| 1167 | nyn_l = nyn_l / 2 |
---|
| 1168 | nzt_l = nzt_l / 2 |
---|
[709] | 1169 | |
---|
[667] | 1170 | ENDDO |
---|
[709] | 1171 | |
---|
| 1172 | ENDIF |
---|
[1677] | 1173 | |
---|
[1] | 1174 | #endif |
---|
| 1175 | |
---|
[1804] | 1176 | #if defined( __parallel ) |
---|
[1] | 1177 | ! |
---|
[1933] | 1178 | !-- Setting of flags for inflow/outflow/nesting conditions. |
---|
[1964] | 1179 | IF ( pleft == MPI_PROC_NULL ) THEN |
---|
[3182] | 1180 | IF ( bc_lr == 'dirichlet/radiation' .OR. bc_lr == 'nested' .OR. & |
---|
| 1181 | bc_lr == 'nesting_offline' ) THEN |
---|
| 1182 | bc_dirichlet_l = .TRUE. |
---|
[1964] | 1183 | ELSEIF ( bc_lr == 'radiation/dirichlet' ) THEN |
---|
[3182] | 1184 | bc_radiation_l = .TRUE. |
---|
[1] | 1185 | ENDIF |
---|
| 1186 | ENDIF |
---|
[1933] | 1187 | |
---|
[1964] | 1188 | IF ( pright == MPI_PROC_NULL ) THEN |
---|
| 1189 | IF ( bc_lr == 'dirichlet/radiation' ) THEN |
---|
[3182] | 1190 | bc_radiation_r = .TRUE. |
---|
| 1191 | ELSEIF ( bc_lr == 'radiation/dirichlet' .OR. bc_lr == 'nested' .OR. & |
---|
| 1192 | bc_lr == 'nesting_offline' ) THEN |
---|
| 1193 | bc_dirichlet_r = .TRUE. |
---|
[1] | 1194 | ENDIF |
---|
| 1195 | ENDIF |
---|
| 1196 | |
---|
[1964] | 1197 | IF ( psouth == MPI_PROC_NULL ) THEN |
---|
| 1198 | IF ( bc_ns == 'dirichlet/radiation' ) THEN |
---|
[3182] | 1199 | bc_radiation_s = .TRUE. |
---|
| 1200 | ELSEIF ( bc_ns == 'radiation/dirichlet' .OR. bc_ns == 'nested' .OR. & |
---|
| 1201 | bc_ns == 'nesting_offline' ) THEN |
---|
| 1202 | bc_dirichlet_s = .TRUE. |
---|
[1] | 1203 | ENDIF |
---|
| 1204 | ENDIF |
---|
| 1205 | |
---|
[1964] | 1206 | IF ( pnorth == MPI_PROC_NULL ) THEN |
---|
[3182] | 1207 | IF ( bc_ns == 'dirichlet/radiation' .OR. bc_ns == 'nested' .OR. & |
---|
| 1208 | bc_ns == 'nesting_offline' ) THEN |
---|
| 1209 | bc_dirichlet_n = .TRUE. |
---|
[1964] | 1210 | ELSEIF ( bc_ns == 'radiation/dirichlet' ) THEN |
---|
[3182] | 1211 | bc_radiation_n = .TRUE. |
---|
[1] | 1212 | ENDIF |
---|
| 1213 | ENDIF |
---|
[2938] | 1214 | ! |
---|
| 1215 | !-- In case of synthetic turbulence geneartor determine ids. |
---|
| 1216 | !-- Please note, if no forcing or nesting is applied, the generator is applied |
---|
| 1217 | !-- only at the left lateral boundary. |
---|
| 1218 | IF ( use_syn_turb_gen ) THEN |
---|
[3182] | 1219 | IF ( bc_dirichlet_l ) THEN |
---|
[2938] | 1220 | id_stg_left_l = myidx |
---|
| 1221 | ELSE |
---|
| 1222 | id_stg_left_l = 0 |
---|
| 1223 | ENDIF |
---|
[3182] | 1224 | IF ( bc_dirichlet_r ) THEN |
---|
[2938] | 1225 | id_stg_right_l = myidx |
---|
| 1226 | ELSE |
---|
| 1227 | id_stg_right_l = 0 |
---|
| 1228 | ENDIF |
---|
[3182] | 1229 | IF ( bc_dirichlet_s ) THEN |
---|
[2938] | 1230 | id_stg_south_l = myidy |
---|
| 1231 | ELSE |
---|
| 1232 | id_stg_south_l = 0 |
---|
| 1233 | ENDIF |
---|
[3182] | 1234 | IF ( bc_dirichlet_n ) THEN |
---|
[2938] | 1235 | id_stg_north_l = myidy |
---|
| 1236 | ELSE |
---|
| 1237 | id_stg_north_l = 0 |
---|
| 1238 | ENDIF |
---|
[1968] | 1239 | |
---|
[2938] | 1240 | IF ( collective_wait ) CALL MPI_BARRIER( comm2d, ierr ) |
---|
| 1241 | CALL MPI_ALLREDUCE( id_stg_left_l, id_stg_left, 1, MPI_INTEGER, & |
---|
| 1242 | MPI_SUM, comm1dx, ierr ) |
---|
| 1243 | |
---|
| 1244 | IF ( collective_wait ) CALL MPI_BARRIER( comm2d, ierr ) |
---|
| 1245 | CALL MPI_ALLREDUCE( id_stg_right_l, id_stg_right, 1, MPI_INTEGER, & |
---|
| 1246 | MPI_SUM, comm1dx, ierr ) |
---|
| 1247 | |
---|
| 1248 | IF ( collective_wait ) CALL MPI_BARRIER( comm2d, ierr ) |
---|
| 1249 | CALL MPI_ALLREDUCE( id_stg_south_l, id_stg_south, 1, MPI_INTEGER, & |
---|
| 1250 | MPI_SUM, comm1dy, ierr ) |
---|
| 1251 | |
---|
| 1252 | IF ( collective_wait ) CALL MPI_BARRIER( comm2d, ierr ) |
---|
| 1253 | CALL MPI_ALLREDUCE( id_stg_north_l, id_stg_north, 1, MPI_INTEGER, & |
---|
| 1254 | MPI_SUM, comm1dy, ierr ) |
---|
| 1255 | |
---|
| 1256 | ENDIF |
---|
| 1257 | |
---|
[151] | 1258 | ! |
---|
| 1259 | !-- Broadcast the id of the inflow PE |
---|
[3182] | 1260 | IF ( bc_dirichlet_l ) THEN |
---|
[163] | 1261 | id_inflow_l = myidx |
---|
[151] | 1262 | ELSE |
---|
| 1263 | id_inflow_l = 0 |
---|
| 1264 | ENDIF |
---|
[622] | 1265 | IF ( collective_wait ) CALL MPI_BARRIER( comm2d, ierr ) |
---|
[151] | 1266 | CALL MPI_ALLREDUCE( id_inflow_l, id_inflow, 1, MPI_INTEGER, MPI_SUM, & |
---|
| 1267 | comm1dx, ierr ) |
---|
| 1268 | |
---|
[163] | 1269 | ! |
---|
| 1270 | !-- Broadcast the id of the recycling plane |
---|
| 1271 | !-- WARNING: needs to be adjusted in case of inflows other than from left side! |
---|
[3884] | 1272 | IF ( turbulent_inflow ) THEN |
---|
| 1273 | |
---|
[4045] | 1274 | IF ( NINT( recycling_width / dx, KIND=idp ) >= nxl .AND. & |
---|
| 1275 | NINT( recycling_width / dx, KIND=idp ) <= nxr ) THEN |
---|
[3884] | 1276 | id_recycling_l = myidx |
---|
| 1277 | ELSE |
---|
| 1278 | id_recycling_l = 0 |
---|
| 1279 | ENDIF |
---|
| 1280 | IF ( collective_wait ) CALL MPI_BARRIER( comm2d, ierr ) |
---|
| 1281 | CALL MPI_ALLREDUCE( id_recycling_l, id_recycling, 1, MPI_INTEGER, MPI_SUM, & |
---|
| 1282 | comm1dx, ierr ) |
---|
| 1283 | |
---|
[163] | 1284 | ENDIF |
---|
| 1285 | |
---|
[2050] | 1286 | ! |
---|
| 1287 | !-- Broadcast the id of the outflow PE and outflow-source plane |
---|
| 1288 | IF ( turbulent_outflow ) THEN |
---|
| 1289 | |
---|
[3182] | 1290 | IF ( bc_radiation_r ) THEN |
---|
[2050] | 1291 | id_outflow_l = myidx |
---|
| 1292 | ELSE |
---|
| 1293 | id_outflow_l = 0 |
---|
| 1294 | ENDIF |
---|
| 1295 | IF ( collective_wait ) CALL MPI_BARRIER( comm2d, ierr ) |
---|
| 1296 | CALL MPI_ALLREDUCE( id_outflow_l, id_outflow, 1, MPI_INTEGER, MPI_SUM, & |
---|
| 1297 | comm1dx, ierr ) |
---|
| 1298 | |
---|
| 1299 | IF ( NINT( outflow_source_plane / dx ) >= nxl .AND. & |
---|
| 1300 | NINT( outflow_source_plane / dx ) <= nxr ) THEN |
---|
| 1301 | id_outflow_source_l = myidx |
---|
| 1302 | ELSE |
---|
| 1303 | id_outflow_source_l = 0 |
---|
| 1304 | ENDIF |
---|
| 1305 | IF ( collective_wait ) CALL MPI_BARRIER( comm2d, ierr ) |
---|
| 1306 | CALL MPI_ALLREDUCE( id_outflow_source_l, id_outflow_source, 1, & |
---|
| 1307 | MPI_INTEGER, MPI_SUM, comm1dx, ierr ) |
---|
| 1308 | |
---|
| 1309 | ENDIF |
---|
| 1310 | |
---|
[3885] | 1311 | CALL location_message( 'creating virtual PE grids + MPI derived data types', 'finished' ) |
---|
[1384] | 1312 | |
---|
[1804] | 1313 | #else |
---|
[1159] | 1314 | IF ( bc_lr == 'dirichlet/radiation' ) THEN |
---|
[3182] | 1315 | bc_dirichlet_l = .TRUE. |
---|
| 1316 | bc_radiation_r = .TRUE. |
---|
[1159] | 1317 | ELSEIF ( bc_lr == 'radiation/dirichlet' ) THEN |
---|
[3182] | 1318 | bc_radiation_l = .TRUE. |
---|
| 1319 | bc_dirichlet_r = .TRUE. |
---|
[1] | 1320 | ENDIF |
---|
| 1321 | |
---|
[1159] | 1322 | IF ( bc_ns == 'dirichlet/radiation' ) THEN |
---|
[3182] | 1323 | bc_dirichlet_n = .TRUE. |
---|
| 1324 | bc_radiation_s = .TRUE. |
---|
[1159] | 1325 | ELSEIF ( bc_ns == 'radiation/dirichlet' ) THEN |
---|
[3182] | 1326 | bc_radiation_n = .TRUE. |
---|
| 1327 | bc_dirichlet_s = .TRUE. |
---|
[1] | 1328 | ENDIF |
---|
| 1329 | #endif |
---|
[807] | 1330 | |
---|
[106] | 1331 | ! |
---|
[978] | 1332 | !-- At the inflow or outflow, u or v, respectively, have to be calculated for |
---|
| 1333 | !-- one more grid point. |
---|
[3182] | 1334 | IF ( bc_dirichlet_l .OR. bc_radiation_l ) THEN |
---|
[106] | 1335 | nxlu = nxl + 1 |
---|
| 1336 | ELSE |
---|
| 1337 | nxlu = nxl |
---|
| 1338 | ENDIF |
---|
[3182] | 1339 | IF ( bc_dirichlet_s .OR. bc_radiation_s ) THEN |
---|
[106] | 1340 | nysv = nys + 1 |
---|
| 1341 | ELSE |
---|
| 1342 | nysv = nys |
---|
| 1343 | ENDIF |
---|
[1] | 1344 | |
---|
[114] | 1345 | ! |
---|
| 1346 | !-- Allocate wall flag arrays used in the multigrid solver |
---|
[1575] | 1347 | IF ( psolver(1:9) == 'multigrid' ) THEN |
---|
[114] | 1348 | |
---|
| 1349 | DO i = maximum_grid_level, 1, -1 |
---|
| 1350 | |
---|
| 1351 | SELECT CASE ( i ) |
---|
| 1352 | |
---|
| 1353 | CASE ( 1 ) |
---|
| 1354 | ALLOCATE( wall_flags_1(nzb:nzt_mg(i)+1, & |
---|
| 1355 | nys_mg(i)-1:nyn_mg(i)+1, & |
---|
| 1356 | nxl_mg(i)-1:nxr_mg(i)+1) ) |
---|
| 1357 | |
---|
| 1358 | CASE ( 2 ) |
---|
| 1359 | ALLOCATE( wall_flags_2(nzb:nzt_mg(i)+1, & |
---|
| 1360 | nys_mg(i)-1:nyn_mg(i)+1, & |
---|
| 1361 | nxl_mg(i)-1:nxr_mg(i)+1) ) |
---|
| 1362 | |
---|
| 1363 | CASE ( 3 ) |
---|
| 1364 | ALLOCATE( wall_flags_3(nzb:nzt_mg(i)+1, & |
---|
| 1365 | nys_mg(i)-1:nyn_mg(i)+1, & |
---|
| 1366 | nxl_mg(i)-1:nxr_mg(i)+1) ) |
---|
| 1367 | |
---|
| 1368 | CASE ( 4 ) |
---|
| 1369 | ALLOCATE( wall_flags_4(nzb:nzt_mg(i)+1, & |
---|
| 1370 | nys_mg(i)-1:nyn_mg(i)+1, & |
---|
| 1371 | nxl_mg(i)-1:nxr_mg(i)+1) ) |
---|
| 1372 | |
---|
| 1373 | CASE ( 5 ) |
---|
| 1374 | ALLOCATE( wall_flags_5(nzb:nzt_mg(i)+1, & |
---|
| 1375 | nys_mg(i)-1:nyn_mg(i)+1, & |
---|
| 1376 | nxl_mg(i)-1:nxr_mg(i)+1) ) |
---|
| 1377 | |
---|
| 1378 | CASE ( 6 ) |
---|
| 1379 | ALLOCATE( wall_flags_6(nzb:nzt_mg(i)+1, & |
---|
| 1380 | nys_mg(i)-1:nyn_mg(i)+1, & |
---|
| 1381 | nxl_mg(i)-1:nxr_mg(i)+1) ) |
---|
| 1382 | |
---|
| 1383 | CASE ( 7 ) |
---|
| 1384 | ALLOCATE( wall_flags_7(nzb:nzt_mg(i)+1, & |
---|
| 1385 | nys_mg(i)-1:nyn_mg(i)+1, & |
---|
| 1386 | nxl_mg(i)-1:nxr_mg(i)+1) ) |
---|
| 1387 | |
---|
| 1388 | CASE ( 8 ) |
---|
| 1389 | ALLOCATE( wall_flags_8(nzb:nzt_mg(i)+1, & |
---|
| 1390 | nys_mg(i)-1:nyn_mg(i)+1, & |
---|
| 1391 | nxl_mg(i)-1:nxr_mg(i)+1) ) |
---|
| 1392 | |
---|
| 1393 | CASE ( 9 ) |
---|
| 1394 | ALLOCATE( wall_flags_9(nzb:nzt_mg(i)+1, & |
---|
| 1395 | nys_mg(i)-1:nyn_mg(i)+1, & |
---|
| 1396 | nxl_mg(i)-1:nxr_mg(i)+1) ) |
---|
| 1397 | |
---|
| 1398 | CASE ( 10 ) |
---|
| 1399 | ALLOCATE( wall_flags_10(nzb:nzt_mg(i)+1, & |
---|
| 1400 | nys_mg(i)-1:nyn_mg(i)+1, & |
---|
| 1401 | nxl_mg(i)-1:nxr_mg(i)+1) ) |
---|
| 1402 | |
---|
| 1403 | CASE DEFAULT |
---|
[254] | 1404 | message_string = 'more than 10 multigrid levels' |
---|
| 1405 | CALL message( 'init_pegrid', 'PA0238', 1, 2, 0, 6, 0 ) |
---|
[114] | 1406 | |
---|
| 1407 | END SELECT |
---|
| 1408 | |
---|
| 1409 | ENDDO |
---|
| 1410 | |
---|
| 1411 | ENDIF |
---|
| 1412 | |
---|
[1] | 1413 | END SUBROUTINE init_pegrid |
---|