00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011 subroutine psmile_enddef_comp (comp_id, global_comp_id, Number_of_Grids, &
00012 comp_info, ierror)
00013
00014
00015
00016
00017 use PRISM_constants
00018
00019 use psmile_grid, only : psmile_transform_extent_cyclic, &
00020 max_num_trans_parts, code_no_trans
00021
00022 use PSMILe, dummy_interface => PSMILe_enddef_comp
00023
00024 implicit none
00025
00026
00027
00028 Integer, Intent (In) :: comp_id
00029
00030
00031
00032 Integer, Intent (In) :: global_comp_id
00033
00034
00035
00036 Integer, Intent (In) :: Number_of_Grids
00037
00038
00039
00040
00041
00042
00043 Type (Enddef_comp), Intent (Out) :: comp_info
00044
00045
00046
00047
00048
00049
00050
00051
00052
00053
00054
00055
00056 integer, Intent (Out) :: ierror
00057
00058
00059
00060
00061
00062
00063
00064 integer :: i, j, n
00065
00066
00067
00068
00069 Real (PSMILe_float_kind) :: extent (1:2, 1:ndim_3d)
00070 Real (PSMILe_float_kind) :: parts (1:2, 1:ndim_3d,
00071 1:Number_of_Grids*max_num_trans_parts)
00072 Integer :: extent_id (Number_of_Grids)
00073 Integer :: grid_type (Number_of_Grids)
00074 Integer :: part_info (nd_extent_infos,
00075 1:Number_of_Grids*max_num_trans_parts)
00076 Integer :: tr_codes (max_num_trans_parts)
00077
00078 Integer :: nparts, n_trans
00079
00080
00081
00082 Integer :: grid_id
00083
00084
00085
00086 integer :: comm, comm_size
00087 integer :: n_total
00088 integer :: count (Comps(comp_id)%size)
00089 integer :: disp (Comps(comp_id)%size)
00090
00091
00092
00093 integer, parameter :: nerrp = 2
00094 integer :: ierrp (nerrp)
00095
00096
00097
00098
00099
00100
00101
00102
00103
00104
00105
00106
00107
00108
00109
00110
00111
00112
00113 Character(len=len_cvs_string), save :: mycvs =
00114 '$Id: psmile_enddef_comp.F90 2687 2010-10-28 15:15:52Z coquart $'
00115
00116
00117
00118 #ifdef VERBOSE
00119 print 9990, trim(ch_id), comp_id
00120
00121 call psmile_flushstd
00122 #endif /* VERBOSE */
00123
00124 #ifdef DEBUG
00125 print*, 'In psmile_enddef_comp global_comp_id :',global_comp_id
00126 print*, 'In psmile_enddef_comp Comps(comp_id)%n_grids :',Comps(comp_id)%n_grids
00127 print*, 'In psmile_enddef_comp Number_of_Grids :',Number_of_Grids
00128 call psmile_flushstd
00129 #endif
00130
00131 #ifdef PRISM_ASSERTION
00132 if (Number_of_Grids < 0 .or. &
00133 Number_of_Grids > Comps(comp_id)%n_grids) then
00134 print *, "Number of Grids", Number_of_Grids, Comps(comp_id)%n_grids
00135 call psmile_assert ( __FILE__, __LINE__, 'Incorrect number of grids')
00136 endif
00137 #endif /* PRISM_ASSERTION */
00138
00139
00140
00141 ierror = 0
00142
00143 comm = Comps(comp_id)%comm
00144 comm_size = Comps(comp_id)%size
00145
00146 comp_info%comp_id = comp_id
00147 comp_info%global_comp_id = global_comp_id
00148 comp_info%size = comm_size
00149
00150
00151
00152
00153
00154
00155
00156
00157
00158
00159
00160
00161
00162
00163
00164
00165
00166
00167
00168
00169
00170
00171 n = 0
00172 do grid_id = 1, Number_of_Grids_allocated
00173
00174 if (Grids(grid_id)%status /= PSMILe_status_free .and. &
00175 Grids(grid_id)%comp_id == comp_id .and. &
00176 Grids(grid_id)%used_for_coupling) then
00177
00178 n = n + 1
00179 extent_id (n) = grid_id
00180 grid_type (n) = Grids(grid_id)%grid_type
00181
00182 endif
00183
00184 if (Appl%stand_alone .and. &
00185 Grids(grid_id)%status /= PSMILe_status_free .and. &
00186 Grids(grid_id)%comp_id == comp_id .and. &
00187 Grids(grid_id)%used_for_io) then
00188
00189 n = n + 1
00190 extent_id (n) = grid_id
00191 grid_type (n) = Grids(grid_id)%grid_type
00192
00193 endif
00194
00195 enddo
00196
00197 #ifdef PRISM_ASSERTION
00198 if (n /= Number_of_Grids) then
00199 write (*, 9970) n, Number_of_Grids
00200 call psmile_assert ( __FILE__, __LINE__, 'n /= Number_of_Grids')
00201 endif
00202 #endif /* PRISM_ASSERTION */
00203
00204 nparts = 0
00205
00206 do n = 1, Number_of_Grids
00207 if (grid_type (n) == PRISM_Gridless) then
00208 grid_id = extent_id (n)
00209
00210 if (associated (Grids(grid_id)%partition)) then
00211 n_trans = size(Grids(grid_id)%partition(:,1))
00212 do i = 1, ndim_3d
00213 do j = 1, n_trans
00214 parts (1, i, nparts+j) = Grids(grid_id)%partition (j, i) + 1
00215 parts (2, i, nparts+j) = Grids(grid_id)%partition (j, i) + &
00216 Grids(grid_id)%extent(j,i)
00217 end do
00218 end do
00219 else
00220 parts (1:2, 1:ndim_3d, nparts+1) = &
00221 Grids(grid_id)%grid_shape (1:2, 1:ndim_3d)
00222 n_trans = 1
00223 endif
00224
00225 part_info (3, nparts+1:nparts+n_trans) = code_no_trans
00226
00227 else
00228 call psmile_get_grid_extent (extent_id (n), extent, ierror)
00229 if ( ierror > 0 ) return
00230
00231 call psmile_transform_extent_cyclic (grid_type (n), extent, &
00232 parts (1,1,nparts+1), tr_codes, n_trans, ierror)
00233 if ( ierror > 0 ) return
00234
00235 part_info (3, nparts+1:nparts+n_trans) = tr_codes (1:n_trans)
00236 endif
00237
00238 part_info (1, nparts+1:nparts+n_trans) = extent_id (n)
00239 part_info (2, nparts+1:nparts+n_trans) = grid_type (n)
00240 part_info (4, nparts+1:nparts+n_trans) = &
00241 Grids(extent_id(n))%global_grid_id
00242
00243 nparts = nparts + n_trans
00244 end do
00245
00246 #ifdef DEBUG
00247 print*, 'In psmile_enddef_comp nparts:',nparts
00248 call psmile_flushstd
00249 #endif
00250
00251
00252
00253 Allocate (comp_info%Number_of_grids_vector(1:comm_size), STAT = ierror)
00254 if ( ierror > 0 ) then
00255 ierrp (1) = ierror
00256 ierrp (2) = comm_size
00257 call psmile_error ( PRISM_Error_Alloc, 'Number_of_grids_vector', &
00258 ierrp, 2, __FILE__, __LINE__ )
00259 return
00260 endif
00261
00262
00263
00264 Allocate (comp_info%psmile_ranks(1:comm_size), STAT = ierror)
00265 if ( ierror > 0 ) then
00266 ierrp (1) = ierror
00267 ierrp (2) = comm_size
00268 call psmile_error ( PRISM_Error_Alloc, 'psmile_ranks', &
00269 ierrp, 2, __FILE__, __LINE__ )
00270 return
00271 endif
00272
00273
00274
00275 call MPI_Allgather (nparts, 1, MPI_INTEGER, &
00276 comp_info%Number_of_Grids_Vector, 1, MPI_INTEGER, &
00277 comm, ierror)
00278 if (ierror /= MPI_SUCCESS) then
00279 ierrp (1) = ierror
00280 ierror = PRISM_Error_MPI
00281
00282 call psmile_error ( ierror, 'MPI_AllGather', &
00283 ierrp, 1, __FILE__, __LINE__ )
00284 return
00285 endif
00286
00287 #ifdef DEBUG
00288 print*, 'In psmile_enddef_comp after all proc of the comp get the infos:'
00289 print*, 'comp_info%Number_of_Grids_Vector :',comp_info%Number_of_Grids_Vector
00290 call psmile_flushstd
00291 #endif
00292
00293 n_total = SUM (comp_info%Number_of_Grids_Vector)
00294
00295
00296
00297
00298
00299 if ( n_total < 1 ) then
00300 print *, trim(ch_id), &
00301 ' : psmile_enddef_comp: No grids of fields defined for component'
00302 call psmile_abort
00303 endif
00304
00305
00306
00307 call MPI_Allgather (psmile_rank, 1, MPI_INTEGER, &
00308 comp_info%psmile_ranks, 1, MPI_INTEGER, &
00309 comm, ierror)
00310 if (ierror /= MPI_SUCCESS) then
00311 ierrp (1) = ierror
00312 ierror = PRISM_Error_MPI
00313
00314 call psmile_error ( ierror, 'MPI_AllGather', &
00315 ierrp, 1, __FILE__, __LINE__ )
00316 return
00317 endif
00318
00319 #ifdef DEBUG
00320 print*, 'comp_info%psmile_ranks :',comp_info%psmile_ranks
00321 call psmile_flushstd
00322 #endif
00323
00324
00325
00326
00327
00328 if (n_total > 0) then
00329
00330 Allocate (comp_info%all_extents(1:2, 1:ndim_3d, 1:n_total), &
00331 STAT = ierror)
00332 if ( ierror > 0 ) then
00333 ierrp (1) = ierror
00334 ierrp (2) = 2 * ndim_3d * n_total
00335 call psmile_error ( PRISM_Error_Alloc, 'all_extents', &
00336 ierrp, 2, __FILE__, __LINE__ )
00337 return
00338 endif
00339
00340
00341
00342 count (:) = comp_info%Number_of_Grids_Vector (:) * (2*ndim_3d)
00343
00344 disp (1) = 0
00345
00346 do i = 2, comm_size
00347 disp (i) = disp (i-1) + count (i-1)
00348 enddo
00349
00350 call MPI_Allgatherv (parts, nparts*(2*ndim_3d), PSMILe_float_datatype, &
00351 comp_info%all_extents, count, disp, PSMILe_float_datatype, &
00352 comm, ierror)
00353 if (ierror /= MPI_SUCCESS) then
00354 ierrp (1) = ierror
00355 ierror = PRISM_Error_MPI
00356
00357 call psmile_error ( ierror, 'MPI_AllGatherv', &
00358 ierrp, 1, __FILE__, __LINE__ )
00359 return
00360 endif
00361
00362 #ifdef DEBUG
00363 print*, 'comp_info%all_extents :',comp_info%all_extents
00364 call psmile_flushstd
00365 #endif
00366
00367
00368
00369 Allocate (comp_info%all_extent_infos(nd_extent_infos, 1:n_total), &
00370 STAT = ierror)
00371 if ( ierror > 0 ) then
00372 ierrp (1) = ierror
00373 ierrp (2) = n_total
00374 call psmile_error ( PRISM_Error_Alloc, 'all_extent_infos', &
00375 ierrp, 2, __FILE__, __LINE__ )
00376 return
00377 endif
00378
00379
00380
00381
00382 count (:) = comp_info%Number_of_Grids_Vector (:) * nd_extent_infos
00383
00384 disp (1) = 0
00385
00386 do i = 2, comm_size
00387 disp (i) = disp (i-1) + count (i-1)
00388 enddo
00389
00390 call MPI_Allgatherv (part_info, nparts*nd_extent_infos, MPI_INTEGER, &
00391 comp_info%all_extent_infos, count, disp, MPI_INTEGER, &
00392 comm, ierror)
00393 if (ierror /= MPI_SUCCESS) then
00394 ierrp (1) = ierror
00395 ierror = PRISM_Error_MPI
00396
00397 call psmile_error ( ierror, 'MPI_AllGatherv', &
00398 ierrp, 1, __FILE__, __LINE__ )
00399 return
00400 endif
00401
00402
00403
00404 call psmile_enddef_comp_periodic (comp_id, extent_id, Number_of_grids, &
00405 ierror)
00406 if ( ierror > 0 ) return
00407
00408 endif
00409
00410
00411
00412 #ifdef VERBOSE
00413 print 9980, trim(ch_id), ierror
00414
00415 call psmile_flushstd
00416 #endif /* VERBOSE */
00417
00418
00419
00420 9990 format (1x, a, ': psmile_enddef_comp: comp_id =', i3)
00421 9980 format (1x, a, ': psmile_enddef_comp: eof, ierror =', i3)
00422
00423 #ifdef PRISM_ASSERTION
00424 9970 format (/1x, 'psmile_enddef_comp: inconsistent number of grids:', &
00425 'n = ', i7, '; Number_of_Grids =', i7)
00426 #endif /* PRISM_ASSERTION */
00427
00428 end subroutine PSMILe_enddef_comp