00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011 subroutine psmile_search_donor_2d_dble ( &
00012 found, locations, len, search, ipart, &
00013 grid_id, first_method_id, var_id, tol, ierror)
00014
00015
00016
00017 use PRISM_constants
00018
00019 use PSMILe, dummy_interface => PSMILe_Search_donor_2d_dble
00020
00021 implicit none
00022
00023
00024
00025 Integer, Intent (In) :: len
00026
00027
00028
00029 Type (Enddef_search) :: search
00030
00031
00032
00033 Integer, Intent (In) :: ipart
00034
00035
00036
00037 Integer, Intent (In) :: grid_id
00038
00039
00040
00041 Integer, Intent (In) :: first_method_id
00042
00043
00044
00045 Integer, Intent (In) :: var_id
00046
00047
00048
00049
00050 Double Precision, Intent (In) :: tol
00051
00052
00053
00054
00055
00056
00057
00058 Integer, Intent (InOut) :: found (len)
00059
00060
00061
00062
00063 Integer, Intent (InOut) :: locations (ndim_2d, len)
00064
00065
00066
00067
00068
00069
00070 integer, Intent (Out) :: ierror
00071
00072
00073
00074
00075
00076
00077
00078 Type (Corner_Block), Pointer :: corner_pointer
00079 Integer :: i
00080
00081
00082
00083 Integer :: ibeg, iend
00084 Integer :: control (2, ndim_3d)
00085
00086 Double Precision, Pointer :: array1 (:)
00087 Double Precision, Pointer :: array2 (:)
00088 Integer :: shape (2, ndim_3d)
00089 Integer :: range (2, ndim_3d)
00090 Integer :: j, len1
00091
00092
00093
00094 Integer :: lev, levc, nlev
00095 Integer :: ijkinc (ndim_3d), ijkcoa (ndim_2d)
00096 Type(Grid), Pointer :: grid_info
00097
00098
00099
00100 Integer, parameter :: nerrp = 2
00101 Integer :: ierrp (nerrp)
00102
00103
00104
00105
00106
00107
00108
00109
00110
00111
00112
00113
00114
00115
00116
00117
00118
00119
00120 Character(len=len_cvs_string), save :: mycvs =
00121 '$Id: psmile_search_donor_2d_dble.F90 2687 2010-10-28 15:15:52Z coquart $'
00122
00123
00124
00125
00126
00127 #ifdef VERBOSE
00128 print 9990, trim(ch_id), Grids(grid_id)%comp_id, search%sender
00129
00130 call psmile_flushstd
00131 #endif /* VERBOSE */
00132
00133 ierror = 0
00134 grid_info => Grids(grid_id)
00135
00136 #ifdef PRISM_ASSERTION
00137 #endif
00138
00139
00140
00141
00142
00143 if (search%grid_type == PRISM_Reglonlatvrt) then
00144
00145
00146
00147 Allocate (array1(len), STAT = ierror)
00148 if ( ierror > 0 ) then
00149 ierrp (1) = ierror
00150 ierrp (2) = len
00151 ierror = PRISM_Error_Alloc
00152 call psmile_error ( ierror, 'array1', &
00153 ierrp, 2, __FILE__, __LINE__ )
00154 return
00155 endif
00156
00157 Allocate (array2(len), STAT = ierror)
00158 if ( ierror > 0 ) then
00159 ierrp (1) = ierror
00160 ierrp (2) = len
00161 ierror = PRISM_Error_Alloc
00162 call psmile_error ( ierror, 'array2', &
00163 ierrp, 2, __FILE__, __LINE__ )
00164 return
00165 endif
00166
00167 shape (:, 1:ndim_2d) = search%range(:, 1:ndim_2d, ipart)
00168 shape (:, ndim_3d) = 1
00169
00170 range = shape
00171
00172 len1 = shape(2,1)-shape(1,1) + 1
00173 #ifdef PRISM_ASSERTION
00174 if (search%dims(1, ipart) /= len1) then
00175 call psmile_assert (__FILE__, __LINE__, &
00176 "dim(1) != len1")
00177 endif
00178 #endif
00179
00180 do j = 1, shape(2,2)-shape(1,2) + 1
00181 array1 ((j-1)*len1+1:j*len1) = search%search_dble(1,ipart)%vector(1:len1)
00182 array2 ((j-1)*len1+1:j*len1) = search%search_dble(2,ipart)%vector(j)
00183 end do
00184
00185 else
00186
00187 array1 => search%search_dble(1,ipart)%vector
00188 array2 => search%search_dble(2,ipart)%vector
00189
00190 if (search%grid_type == PRISM_Irrlonlat_regvrt) then
00191 shape(:, 1:ndim_2d) = search%shape(:, 1:ndim_2d, ipart)
00192 shape(:, ndim_3d) = 1
00193
00194 range(:, 1:ndim_2d) = search%range(:, 1:ndim_2d, ipart)
00195 range(:, ndim_3d) = 1
00196
00197 else if (search%grid_type == PRISM_Gaussreduced_regvrt) then
00198 shape(:, 1:ndim_2d) = search%shape(:, 1:ndim_2d, ipart)
00199 shape(:, ndim_3d) = 1
00200
00201 range(:, 1:ndim_2d) = search%range(:, 1:ndim_2d, ipart)
00202 range(:, ndim_3d) = 1
00203
00204 else
00205 shape(:, :) = search%shape(:, :, ipart)
00206 range(:, :) = search%range(:, :, ipart)
00207 endif
00208 endif
00209
00210
00211
00212
00213
00214 nlev = grid_info%nlev
00215 lev = nlev
00216
00217 ibeg = 1
00218 iend = len
00219
00220 #ifdef PRISM_ASSERTION
00221 if (grid_info%mg_infos(lev)%levdim(1) /= 0 .or. &
00222 grid_info%mg_infos(lev)%levdim(2) /= 0) then
00223
00224 call psmile_assert (__FILE__, __LINE__, &
00225 "coarsest level dim != 0")
00226 endif
00227 #endif
00228
00229 call psmile_mg_coarse_2d_dble (lev, &
00230 grid_info%mg_infos(lev)%double_arrays%chmin, &
00231 grid_info%mg_infos(lev)%double_arrays%chmax, &
00232 found, locations, array1, array2, &
00233 ibeg, iend)
00234
00235 if (ibeg > iend) then
00236
00237 #ifdef VERBOSE
00238 print 9980, trim(ch_id), grid_info%comp_id, search%sender, ierror
00239
00240 call psmile_flushstd
00241 #endif /* VERBOSE */
00242 return
00243 endif
00244
00245 control = range
00246
00247
00248
00249
00250
00251
00252
00253
00254
00255
00256
00257
00258
00259
00260
00261
00262
00263 ijkinc (:) = 1
00264
00265 do lev = nlev-1, 1, -1
00266
00267 levc = lev + 1
00268
00269 ijkcoa (:) = 2
00270
00271 do i = 1, ndim_2d
00272 if (grid_info%mg_infos(lev )%levdim(i) == &
00273 grid_info%mg_infos(levc)%levdim(i)) then
00274 ijkcoa (i) = 1
00275 endif
00276 enddo
00277
00278
00279
00280 call psmile_mg_next_level_2d_dble ( grid_id, lev, nlev, &
00281 grid_info%mg_infos(lev)%double_arrays%chmin(1)%vector, &
00282 grid_info%mg_infos(lev)%double_arrays%chmin(2)%vector, &
00283 grid_info%mg_infos(lev)%double_arrays%chmax(1)%vector, &
00284 grid_info%mg_infos(lev)%double_arrays%chmax(2)%vector, &
00285 grid_info%mg_infos(lev)%double_arrays%midp(1)%vector, &
00286 grid_info%mg_infos(lev)%double_arrays%midp(2)%vector, &
00287 grid_info%mg_infos(lev)%levdim, &
00288 found, locations, range, &
00289 array1, array2, shape, &
00290 control, ijkinc, ijkcoa, ierror)
00291
00292 if (ierror > 0) return
00293
00294
00295 enddo
00296
00297
00298
00299
00300 do i = 1, len
00301 locations (1:ndim_2d, i) = locations (1:ndim_2d, i) &
00302 + grid_info%ijk0 (1:ndim_2d)
00303 end do
00304
00305 #define SEARCH_ON_FINAL_GRID
00306 #ifdef SEARCH_ON_FINAL_GRID
00307
00308
00309
00310
00311 corner_pointer => Grids(grid_id)%corner_pointer
00312
00313 if (corner_pointer%nbr_corners == 8) then
00314
00315
00316 call psmile_mg_final_2d_dble (grid_id, nlev, &
00317 grid_info%mg_infos(1)%double_arrays%chmin(1)%vector, &
00318 grid_info%mg_infos(1)%double_arrays%chmin(2)%vector, &
00319 grid_info%mg_infos(1)%double_arrays%chmax(1)%vector, &
00320 grid_info%mg_infos(1)%double_arrays%chmax(2)%vector, &
00321 grid_info%mg_infos(1)%double_arrays%midp(1)%vector, &
00322 grid_info%mg_infos(1)%double_arrays%midp(2)%vector, &
00323 grid_info%mg_infos(1)%levdim, &
00324 found, locations, range, &
00325 array1, array2, shape, control, &
00326 corner_pointer%corners_dble(1)%vector, &
00327 corner_pointer%corners_dble(2)%vector, &
00328 corner_pointer%corner_shape, &
00329 corner_pointer%nbr_corners/2, &
00330 tol, ierror)
00331 if (ierror > 0) return
00332 endif
00333 #endif
00334
00335
00336
00337
00338
00339 if (search%grid_type == PRISM_Reglonlatvrt) then
00340 Deallocate (array2)
00341 Deallocate (array1)
00342 endif
00343
00344 #ifdef VERBOSE
00345 print 9980, trim(ch_id), grid_info%comp_id, search%sender, ierror
00346
00347 call psmile_flushstd
00348 #endif /* VERBOSE */
00349
00350
00351
00352 9990 format (1x, a, ': PSMILe_Search_donor_2d_dble: comp_id =', i3, &
00353 '; sender =', i4)
00354 9980 format (1x, a, ': PSMILe_Search_donor_2d_dble: comp_id =', i3, &
00355 '; eof sender =', i3, ', ierror =', i4)
00356
00357 end subroutine PSMILe_Search_donor_2d_dble