00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011 subroutine psmile_send_req_coords_dble (msg_intersections, &
00012 dest, tag, ierror)
00013
00014
00015
00016 use PRISM_constants
00017
00018 use PSMILe, dummy_interface => PSMILe_Send_req_coords_dble
00019
00020 implicit none
00021
00022
00023
00024 Type (enddef_msg_intersections), Intent (In) :: msg_intersections
00025
00026
00027
00028 Integer, Intent (In) :: dest
00029
00030
00031
00032
00033 Integer, Intent (In) :: tag
00034
00035
00036
00037
00038
00039 Integer, Intent (Out) :: ierror
00040
00041
00042
00043
00044
00045
00046
00047
00048
00049 Integer :: grid_id, method_id
00050 Integer :: field_id, mask_id
00051 Type (Coords_Block), Pointer :: coords_pointer
00052
00053
00054
00055 Integer :: subarray_type
00056 Integer :: sizes (ndim_3d)
00057 Integer :: subsizes (ndim_3d)
00058 Integer :: starts (ndim_3d)
00059
00060
00061
00062 Integer :: inter (2, ndim_3d)
00063
00064 Integer :: i, ipart, npart
00065
00066
00067
00068 Integer, parameter :: nerrp = 3
00069 Integer :: ierrp (nerrp)
00070
00071
00072
00073
00074
00075
00076
00077
00078
00079
00080
00081
00082
00083
00084
00085
00086
00087
00088
00089
00090 Character(len=len_cvs_string), save :: mycvs =
00091 '$Id: psmile_send_req_coords_dble.F90 2787 2010-11-29 16:51:32Z hanke $'
00092
00093
00094
00095
00096
00097 ierror = 0
00098
00099 npart = msg_intersections%num_parts
00100 grid_id = msg_intersections%tgt_grid_id
00101 method_id = msg_intersections%first_tgt_method_id
00102 field_id = msg_intersections%first_tgt_var_id
00103
00104
00105 mask_id = msg_intersections%tgt_mask_id
00106
00107 #ifdef VERBOSE
00108 print 9970, trim(ch_id)
00109 print 9990, trim(ch_id), grid_id, dest
00110
00111 call psmile_flushstd
00112 #endif /* VERBOSE */
00113
00114 coords_pointer => Methods(method_id)%coords_pointer
00115
00116 #ifdef PRISM_ASSERTION
00117
00118
00119
00120 do i = 1, ndim_3d
00121 if (.not. Associated(coords_pointer%coords_dble(1)%vector) ) then
00122 print *, trim(ch_id), 'method id, i', method_id, i
00123 call psmile_assert ( __FILE__, __LINE__, &
00124 'Pointer coords_dble(i)%vector is not set')
00125 endif
00126 end do
00127
00128 if (field_id < 1 .or. field_id > Number_of_Fields_allocated) then
00129 print *, trim(ch_id), 'method id, field_id', method_id, field_id
00130 call psmile_assert ( __FILE__, __LINE__, &
00131 'Field id is out of range')
00132 endif
00133
00134 if (mask_id /= PRISM_UNDEFINED .and. &
00135 (mask_id < 1 .or. mask_id > Number_of_Masks_allocated)) then
00136 print *, trim(ch_id), 'method id, field_id, mask_id', &
00137 method_id, field_id, mask_id
00138 call psmile_assert ( __FILE__, __LINE__, &
00139 'Mask id is out of range')
00140 endif
00141
00142 #endif /* PRISM_ASSERTION */
00143
00144 do ipart = 1, npart
00145
00146 inter = msg_intersections%intersections(ipart)%intersection
00147
00148 select case ( Grids(grid_id)%grid_type )
00149
00150
00151
00152
00153
00154
00155 case (PRISM_Reglonlatvrt)
00156
00157 do i = 1, ndim_3d
00158 call MPI_Send (coords_pointer%coords_dble(i)%vector( &
00159 inter(1,i)-coords_pointer%coords_shape(1,i)+1), &
00160 inter(2,i)-inter(1,i)+1, MPI_DOUBLE_PRECISION, &
00161 dest, tag, comm_psmile, ierror)
00162 if (ierror /= MPI_SUCCESS) then
00163 ierrp (1) = ierror
00164 ierrp (2) = dest
00165 ierrp (3) = tag
00166 ierror = PRISM_Error_Send
00167
00168 call psmile_error (ierror, 'MPI_Send', &
00169 ierrp, 3, __FILE__, __LINE__ )
00170 return
00171 endif
00172 end do
00173
00174
00175
00176
00177
00178
00179 case (PRISM_Irrlonlat_regvrt)
00180
00181 sizes (1) = coords_pointer%coords_shape(2,1) - &
00182 coords_pointer%coords_shape(1,1) + 1
00183 sizes (2) = coords_pointer%coords_shape(2,2) - &
00184 coords_pointer%coords_shape(1,2) + 1
00185
00186 subsizes (1) = inter(2,1) - inter (1,1) + 1
00187 subsizes (2) = inter(2,2) - inter (1,2) + 1
00188
00189 starts (1) = inter (1,1) - coords_pointer%coords_shape(1,1)
00190 starts (2) = inter (1,2) - coords_pointer%coords_shape(1,2)
00191
00192 #if defined ( PRISM_with_MPI1 )
00193 call psmile_type_create_subarray (ndim_2d, sizes, subsizes, starts, &
00194 MPI_DOUBLE_PRECISION, subarray_type, ierror)
00195 if (ierror /= MPI_SUCCESS) then
00196 ierrp (1) = ierror
00197 ierror = PRISM_Error_MPI
00198
00199 call psmile_error (ierror, 'PSMILe_Type_create_subarry', &
00200 ierrp, 1, __FILE__, __LINE__ )
00201 return
00202 endif
00203 #else
00204 call MPI_Type_create_subarray (ndim_2d, sizes, subsizes, starts, &
00205 MPI_ORDER_FORTRAN, MPI_DOUBLE_PRECISION, &
00206 subarray_type, ierror)
00207 if (ierror /= MPI_SUCCESS) then
00208 ierrp (1) = ierror
00209 ierror = PRISM_Error_MPI
00210
00211 call psmile_error (ierror, 'MPI_Type_create_subarry', &
00212 ierrp, 1, __FILE__, __LINE__ )
00213 return
00214 endif
00215 #endif
00216
00217 call MPI_Type_commit (subarray_type, ierror)
00218 if (ierror /= MPI_SUCCESS) then
00219 ierrp (1) = ierror
00220 ierror = PRISM_Error_MPI
00221
00222 call psmile_error ( ierror, 'MPI_Type_commit', &
00223 ierrp, 1, __FILE__, __LINE__ )
00224 return
00225 endif
00226
00227 do i = 1, ndim_2d
00228 call MPI_Send (coords_pointer%coords_dble(i)%vector, 1, &
00229 subarray_type, dest, tag, comm_psmile, ierror)
00230 if (ierror /= MPI_SUCCESS) then
00231 ierrp (1) = ierror
00232 ierrp (2) = dest
00233 ierrp (3) = tag
00234 ierror = PRISM_Error_Send
00235
00236 call psmile_error (ierror, 'MPI_Send', &
00237 ierrp, 3, __FILE__, __LINE__ )
00238 return
00239 endif
00240 end do
00241
00242 call MPI_Send (coords_pointer%coords_dble(3)%vector(inter(1,3)- &
00243 coords_pointer%coords_shape(1,3)+1), &
00244 inter(2,3)-inter(1,3)+1, MPI_DOUBLE_PRECISION, &
00245 dest, tag, comm_psmile, ierror)
00246 if (ierror /= MPI_SUCCESS) then
00247 ierrp (1) = ierror
00248 ierrp (2) = dest
00249 ierrp (3) = tag
00250 ierror = PRISM_Error_Send
00251
00252 call psmile_error (ierror, 'MPI_Send', &
00253 ierrp, 3, __FILE__, __LINE__ )
00254 return
00255 endif
00256 #ifdef DEBUGX
00257 print *, 'psmile_send_req_coords_dble', &
00258 coords_pointer%coords_dble(3)%vector(inter(1,3)- &
00259 coords_pointer%coords_shape(1,3)+1), &
00260 coords_pointer%coords_dble(3)%vector(inter(1,3)- &
00261 coords_pointer%coords_shape(1,3)+2)
00262 print *, 'psmile_send_req_coords_dble', inter(1,3), inter(2,3), &
00263 coords_pointer%coords_shape(1,3), &
00264 coords_pointer%coords_shape(2,3)
00265 #endif
00266
00267 call MPI_type_free (subarray_type, ierror)
00268 if (ierror /= MPI_SUCCESS) then
00269 ierrp (1) = ierror
00270 ierror = PRISM_Error_MPI
00271
00272 call psmile_error (ierror, 'MPI_Type_free', &
00273 ierrp, 1, __FILE__, __LINE__ )
00274 return
00275 endif
00276
00277
00278
00279
00280
00281 case (PRISM_Irrlonlatvrt)
00282
00283 sizes (:) = coords_pointer%coords_shape(2,:) - &
00284 coords_pointer%coords_shape(1,:) + 1
00285
00286 subsizes (:) = inter(2,:) - inter (1,:) + 1
00287
00288 starts (:) = inter (1,:) - coords_pointer%coords_shape(1,:)
00289
00290 #if defined ( PRISM_with_MPI1 )
00291 call psmile_type_create_subarray (ndim_3d, sizes, subsizes, starts, &
00292 MPI_DOUBLE_PRECISION, subarray_type, ierror)
00293 if (ierror /= MPI_SUCCESS) then
00294 ierrp (1) = ierror
00295 ierror = PRISM_Error_MPI
00296
00297 call psmile_error (ierror, 'PSMILe_Type_create_subarry', &
00298 ierrp, 1, __FILE__, __LINE__ )
00299 return
00300 endif
00301 #else
00302 call MPI_Type_create_subarray (ndim_3d, sizes, subsizes, starts, &
00303 MPI_ORDER_FORTRAN, MPI_DOUBLE_PRECISION, &
00304 subarray_type, ierror)
00305 if (ierror /= MPI_SUCCESS) then
00306 ierrp (1) = ierror
00307 ierror = PRISM_Error_MPI
00308
00309 call psmile_error (ierror, 'MPI_Type_create_subarry', &
00310 ierrp, 1, __FILE__, __LINE__ )
00311 return
00312 endif
00313 #endif
00314
00315 call MPI_Type_commit (subarray_type, ierror)
00316 if (ierror /= MPI_SUCCESS) then
00317 ierrp (1) = ierror
00318 ierror = PRISM_Error_MPI
00319
00320 call psmile_error ( ierror, 'MPI_Type_commit', &
00321 ierrp, 1, __FILE__, __LINE__ )
00322 return
00323 endif
00324
00325 do i = 1, ndim_3d
00326 call MPI_Send (coords_pointer%coords_dble(i)%vector, 1, &
00327 subarray_type, dest, tag, comm_psmile, ierror)
00328 if (ierror /= MPI_SUCCESS) then
00329 ierrp (1) = ierror
00330 ierrp (2) = dest
00331 ierrp (3) = tag
00332 ierror = PRISM_Error_Send
00333
00334 call psmile_error (ierror, 'MPI_Send', &
00335 ierrp, 3, __FILE__, __LINE__ )
00336 return
00337 endif
00338 end do
00339
00340 call MPI_type_free (subarray_type, ierror)
00341 if (ierror /= MPI_SUCCESS) then
00342 ierrp (1) = ierror
00343 ierror = PRISM_Error_MPI
00344
00345 call psmile_error (ierror, 'MPI_Type_free', &
00346 ierrp, 1, __FILE__, __LINE__ )
00347 return
00348 endif
00349
00350
00351
00352
00353
00354 case (PRISM_Gaussreduced_regvrt)
00355
00356 do i = 1, 2
00357
00358 call MPI_Send (coords_pointer%coords_dble(i)%vector( &
00359 inter(1,1)-coords_pointer%coords_shape(1,1)+1), &
00360 inter(2,1)-inter(1,1)+1, MPI_DOUBLE_PRECISION, &
00361 dest, tag, comm_psmile, ierror)
00362 if (ierror /= MPI_SUCCESS) then
00363 ierrp (1) = ierror
00364 ierrp (2) = dest
00365 ierrp (3) = tag
00366 ierror = PRISM_Error_Send
00367
00368 call psmile_error (ierror, 'MPI_Send', &
00369 ierrp, 3, __FILE__, __LINE__ )
00370 return
00371 endif
00372 end do
00373
00374 call MPI_Send (coords_pointer%coords_dble(3)%vector( &
00375 inter(1,3)-coords_pointer%coords_shape(1,3)+1), &
00376 inter(2,3)-inter(1,3)+1, MPI_DOUBLE_PRECISION, &
00377 dest, tag, comm_psmile, ierror)
00378 if (ierror /= MPI_SUCCESS) then
00379 ierrp (1) = ierror
00380 ierrp (2) = dest
00381 ierrp (3) = tag
00382 ierror = PRISM_Error_Send
00383
00384 call psmile_error (ierror, 'MPI_Send', &
00385 ierrp, 3, __FILE__, __LINE__ )
00386 return
00387 endif
00388
00389
00390
00391
00392
00393 case DEFAULT
00394
00395 ierrp (1) = Grids(grid_id)%grid_type
00396
00397 ierror = PRISM_Error_Internal
00398
00399 call psmile_error ( ierror, 'unsupported grid generation type', &
00400 ierrp, 1, __FILE__, __LINE__ )
00401
00402 end select
00403
00404
00405
00406 if (mask_id /= PRISM_UNDEFINED) then
00407 sizes (:) = Masks(mask_id)%mask_shape(2,:) - &
00408 Masks(mask_id)%mask_shape(1,:) + 1
00409
00410 subsizes (:) = inter(2,:) - inter (1,:) + 1
00411
00412 starts (:) = inter (1,:) - Masks(mask_id)%mask_shape(1,:)
00413
00414 #if defined ( PRISM_with_MPI1 )
00415 call psmile_type_create_subarray (ndim_3d, sizes, subsizes, starts, &
00416 MPI_LOGICAL, subarray_type, ierror)
00417 if (ierror /= MPI_SUCCESS) then
00418 ierrp (1) = ierror
00419 ierror = PRISM_Error_MPI
00420
00421 call psmile_error (ierror, 'PSMILe_Type_create_subarry', &
00422 ierrp, 1, __FILE__, __LINE__ )
00423 return
00424 endif
00425 #else
00426 call MPI_Type_create_subarray (ndim_3d, sizes, subsizes, starts, &
00427 MPI_ORDER_FORTRAN, MPI_LOGICAL, &
00428 subarray_type, ierror)
00429 if (ierror /= MPI_SUCCESS) then
00430 ierrp (1) = ierror
00431 ierror = PRISM_Error_MPI
00432
00433 call psmile_error (ierror, 'MPI_Type_create_subarry', &
00434 ierrp, 1, __FILE__, __LINE__ )
00435 return
00436 endif
00437 #endif
00438
00439 call MPI_Type_commit (subarray_type, ierror)
00440 if (ierror /= MPI_SUCCESS) then
00441 ierrp (1) = ierror
00442 ierror = PRISM_Error_MPI
00443
00444 call psmile_error ( ierror, 'MPI_Type_commit', &
00445 ierrp, 1, __FILE__, __LINE__ )
00446 return
00447 endif
00448
00449 call MPI_Send (Masks(mask_id)%mask_array, 1, subarray_type, &
00450 dest, tag, comm_psmile, ierror)
00451 if (ierror /= MPI_SUCCESS) then
00452 ierrp (1) = ierror
00453 ierrp (2) = dest
00454 ierrp (3) = tag
00455 ierror = PRISM_Error_Send
00456
00457 call psmile_error (ierror, 'MPI_Send', &
00458 ierrp, 3, __FILE__, __LINE__ )
00459 return
00460 endif
00461
00462 call MPI_type_free (subarray_type, ierror)
00463 if (ierror /= MPI_SUCCESS) then
00464 ierrp (1) = ierror
00465 ierror = PRISM_Error_MPI
00466
00467 call psmile_error (ierror, 'MPI_Type_free', &
00468 ierrp, 1, __FILE__, __LINE__ )
00469 return
00470 endif
00471
00472 endif
00473
00474 end do
00475
00476
00477
00478 #ifdef VERBOSE
00479 print 9980, trim(ch_id), grid_id, ierror
00480
00481 call psmile_flushstd
00482 #endif /* VERBOSE */
00483
00484
00485
00486 9990 format (1x, a, ': PSMILe_Send_req_coords_dble: grid_id =', i3, &
00487 '; dest =', i4)
00488 9980 format (1x, a, ': PSMILe_Send_req_coords_dble: end grid_id =', i3, &
00489 '; ierror =', i4)
00490 9970 format (1x, a, ': PSMILe_Send_req_coords_dble: start')
00491
00492 end subroutine PSMILe_Send_req_coords_dble