00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011 subroutine psmile_get_next_field (comp_info, search, &
00012 field_list, n_vars, n_vars_ret, var_id, ierror)
00013
00014
00015
00016 use PRISM_constants
00017
00018 use PSMILe, dummy_interface => PSMILe_Get_next_field
00019
00020 implicit none
00021
00022
00023
00024
00025
00026 Type (Enddef_comp), Intent (In) :: comp_info
00027
00028
00029
00030
00031 Type (Enddef_search), Intent (InOut) :: search
00032
00033
00034
00035 Integer, Intent (In) :: n_vars
00036
00037
00038
00039 Integer, Intent (InOut) :: n_vars_ret
00040
00041
00042
00043 Integer, Intent (In) :: field_list (nd_field_list, n_vars)
00044
00045
00046
00047
00048
00049 Integer, Intent (Out) :: var_id
00050
00051
00052
00053
00054 Integer, Intent (Out) :: ierror
00055
00056
00057
00058
00059
00060
00061
00062 logical, parameter :: new_search = .false.
00063
00064
00065
00066
00067
00068 Integer :: npart, n
00069
00070
00071
00072 Integer, Allocatable :: recv_req (:)
00073 Integer :: save_lreq (2:3)
00074 Integer :: index
00075 Integer :: status (MPI_STATUS_SIZE)
00076 #ifndef PRISM_with_MPI2
00077 Integer, Allocatable :: statuses (:, :)
00078 #endif
00079
00080
00081
00082 Integer, parameter :: nerrp = 3
00083 Integer :: ierrp (nerrp)
00084
00085
00086
00087
00088
00089
00090
00091
00092
00093
00094
00095
00096
00097
00098
00099
00100
00101
00102
00103
00104 Character(len=len_cvs_string), save :: mycvs =
00105 '$Id: psmile_get_next_field.F90 2788 2010-11-30 14:34:07Z hanke $'
00106
00107
00108
00109
00110
00111 #ifdef VERBOSE
00112 print 9990, trim(ch_id)
00113
00114 call psmile_flushstd
00115 #endif /* VERBOSE */
00116
00117 ierror = 0
00118
00119 npart = search%npart
00120
00121 n_vars_ret = n_vars_ret + 1
00122
00123
00124
00125
00126
00127
00128
00129
00130
00131 #ifdef PRISM_ASSERTION
00132 if ( search%npart /= search%msg_intersections%num_parts ) then
00133 call psmile_assert ( __FILE__, __LINE__, &
00134 ' search%npart /= search%msg_intersections%num_parts ')
00135 endif
00136 #endif
00137 if ( search%msg_intersections%requires_conserv_remap == PSMILe_conserv2D .and. &
00138 field_list (6, n_vars_ret) /= PSMILe_conserv2D ) then
00139 #ifdef PRISM_ASSERTION
00140 if ( search%grid_type == PRISM_Gaussreduced_regvrt) then
00141 call psmile_assert ( __FILE__, __LINE__, &
00142 'Case is not yet supported for PRISM_Gaussreduced_regvrt.')
00143 endif
00144 #endif
00145 if ( search%grid_type /= PRISM_Gaussreduced_regvrt ) then
00146 do n = 1, npart
00147 search%msg_intersections%intersections(n)%intersection(2,1:2) = &
00148 search%msg_intersections%intersections(n)%intersection(2,1:2) - 1
00149 enddo
00150 endif
00151
00152 else if ( search%msg_intersections%requires_conserv_remap == PSMILe_conserv3D .and. &
00153 field_list (6, n_vars_ret) /= PSMILe_conserv3D ) then
00154 #ifdef PRISM_ASSERTION
00155 if ( search%grid_type == PRISM_Gaussreduced_regvrt) then
00156 call psmile_assert ( __FILE__, __LINE__, &
00157 'Case is not yet supported for PRISM_Gaussreduced_regvrt.')
00158 endif
00159 #endif
00160 if ( search%grid_type /= PRISM_Gaussreduced_regvrt ) then
00161 do n = 1, npart
00162 search%msg_intersections%intersections(n)%intersection(2,:) = &
00163 search%msg_intersections%intersections(n)%intersection(2,:) - 1
00164 enddo
00165 endif
00166
00167 endif
00168
00169 search%msg_intersections%first_tgt_method_id = field_list (1, n_vars_ret)
00170 search%msg_intersections%first_tgt_var_id = field_list (2, n_vars_ret)
00171 search%msg_intersections%tgt_mask_id = field_list (3, n_vars_ret)
00172 search%msg_intersections%transient_in_id = field_list (4, n_vars_ret)
00173 search%msg_intersections%transient_out_id = field_list (5, n_vars_ret)
00174 search%msg_intersections%requires_conserv_remap = field_list (6, n_vars_ret)
00175
00176
00177
00178
00179
00180
00181 call psmile_pack_msg_intersections(search%msg_intersections, search%msgint)
00182
00183 call psmile_bsend (search%msgint, nd_msgint, MPI_INTEGER, &
00184 search%sender, reqtag, comm_psmile, ierror)
00185
00186 if (ierror /= MPI_SUCCESS) then
00187 ierrp (1) = ierror
00188 ierrp (2) = search%sender
00189 ierrp (3) = reqtag
00190 ierror = PRISM_Error_Send
00191
00192 call psmile_error (ierror, 'MPI_Send', &
00193 ierrp, 3, __FILE__, __LINE__ )
00194 return
00195 endif
00196
00197
00198
00199 Allocate (recv_req ((ndim_3d+2)*npart), stat = ierror)
00200 if ( ierror /= 0 ) then
00201 ierrp (1) = (ndim_3d+1) * npart
00202 ierror = PRISM_Error_Alloc
00203 call psmile_error ( ierror, 'recv_req', &
00204 ierrp, 1, __FILE__, __LINE__ )
00205 return
00206 endif
00207
00208 recv_req = MPI_REQUEST_NULL
00209
00210 call psmile_recv_req_subgrid (search%msg_intersections, &
00211 search%sender, grdtag, search, &
00212 recv_req, new_search, ierror)
00213 if (ierror > 0) return
00214
00215
00216
00217
00218 call psmile_find_corr_field (comp_info, search, &
00219 var_id, ierror)
00220 if (ierror > 0) return
00221
00222 #ifdef PRISM_ASSERTION
00223 if (paction%lrequest (3) /= MPI_REQUEST_NULL) then
00224 call psmile_assert ( __FILE__, __LINE__, &
00225 'Why is request for grid data active ?')
00226
00227 endif
00228 #endif
00229
00230
00231
00232
00233
00234
00235
00236
00237 save_lreq (2:3) = paction%lrequest (2:3)
00238
00239 paction%lrequest (2) = MPI_REQUEST_NULL
00240 paction%lrequest (3) = recv_req(1)
00241
00242 index = 0
00243 #ifdef DEBUGX
00244 do n = 1, paction%nreq
00245 print *, ' Request list waitany ', n, paction%lrequest(n)
00246 enddo
00247 #endif
00248 do while (index /= 3)
00249
00250 call MPI_Waitany (paction%nreq, paction%lrequest, index, status, ierror)
00251
00252 if ( ierror /= MPI_SUCCESS ) then
00253 ierrp (1) = ierror
00254 ierrp (2) = status (MPI_SOURCE)
00255 ierrp (3) = status (MPI_TAG)
00256
00257 ierror = PRISM_Error_MPI
00258
00259 call psmile_error ( ierror, 'MPI_Waitany', &
00260 ierrp, 3, __FILE__, __LINE__ )
00261 return
00262 endif
00263
00264 #ifdef PRISM_ASSERTION
00265 if (index == MPI_UNDEFINED) then
00266 call psmile_assert ( __FILE__, __LINE__, &
00267 'request list is empty')
00268 endif
00269 #endif
00270
00271 if (index /= 3) then
00272 call psmile_enddef_action (search, index, status, ierror)
00273 if (ierror > 0) return
00274 endif
00275 end do
00276
00277
00278
00279 #ifdef PRISM_ASSERTION
00280 if (paction%lrequest (2) /= MPI_REQUEST_NULL .or. &
00281 paction%lrequest (3) /= MPI_REQUEST_NULL) then
00282 print *, 'request: ', paction%lrequest (2:3)
00283 call psmile_assert ( __FILE__, __LINE__, &
00284 'Illegal request stored')
00285
00286 endif
00287 #endif
00288
00289 paction%lrequest (2:3) = save_lreq (2:3)
00290
00291
00292
00293 #ifdef PRISM_with_MPI2
00294 call MPI_Waitall ((ndim_3d+2)*npart-1, recv_req(2:), &
00295 MPI_STATUSES_IGNORE, ierror)
00296 #else
00297 Allocate (statuses (MPI_STATUS_SIZE,npart*(ndim_3d+2)-1), stat = ierror)
00298 if ( ierror /= 0 ) then
00299 ierrp (1) = MPI_STATUS_SIZE * (npart*(ndim_3d+2)-1)
00300 ierror = PRISM_Error_Alloc
00301 call psmile_error ( ierror, 'statuses', &
00302 ierrp, 1, __FILE__, __LINE__ )
00303 return
00304 endif
00305
00306 #ifdef DEBUGX
00307 do n = 2, (ndim_3d+2)*npart
00308 print *, ' Request list waitall ', n, recv_req(n)
00309 enddo
00310 #endif
00311 call MPI_Waitall ((ndim_3d+2)*npart-1, recv_req(2:), &
00312 statuses, ierror)
00313 #endif
00314
00315 if ( ierror /= MPI_SUCCESS ) then
00316 ierrp (1) = ierror
00317
00318 ierror = PRISM_Error_MPI
00319
00320 call psmile_error ( ierror, 'MPI_Waitall', &
00321 ierrp, 1, __FILE__, __LINE__ )
00322 return
00323 endif
00324
00325 #ifndef PRISM_with_MPI2
00326 Deallocate (statuses)
00327 #endif
00328 Deallocate (recv_req)
00329
00330
00331
00332
00333 call psmile_transform_coords (comp_info, search, ierror)
00334 if (ierror > 0) return
00335
00336
00337
00338 #ifdef VERBOSE
00339 print 9980, trim(ch_id), ierror
00340
00341 call psmile_flushstd
00342 #endif /* VERBOSE */
00343
00344
00345
00346 #ifdef VERBOSE
00347
00348 9990 format (1x, a, ': psmile_get_next_field:')
00349 9980 format (1x, a, ': psmile_get_next_field: eof ierror =', i3)
00350
00351 #endif /* VERBOSE */
00352
00353 end subroutine PSMILe_Get_next_field