00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011 subroutine psmile_mg_method_1d_real ( &
00012 comp_info, nlev, &
00013 found, loc, range, &
00014 coords1, &
00015 shape, control, &
00016 method_id, &
00017 x_coords, &
00018 coords_shape, &
00019 grid_valid_shape, cyclic, &
00020 chmin, chmax, &
00021 tol, ierror)
00022
00023
00024
00025 use PRISM_constants
00026
00027 use PSMILe, dummy_interface => PSMILe_mg_method_1d_real
00028
00029 implicit none
00030
00031
00032
00033 Type (Enddef_comp), Intent (In) :: comp_info
00034
00035
00036
00037
00038 Integer, Intent (In) :: nlev
00039
00040
00041
00042 Integer, Intent (In) :: range (2, ndim_3d)
00043
00044
00045
00046 Integer, Intent (In) :: shape (2, ndim_3d)
00047
00048
00049
00050 Integer, Intent (InOut) :: found (range(1,1):range(2,1),
00051 range(1,2):range(2,2),
00052 range(1,3):range(2,3))
00053
00054
00055
00056
00057
00058
00059
00060
00061
00062
00063
00064
00065
00066 Integer, Intent (InOut) :: loc (range(1,1):range(2,1),
00067 range(1,2):range(2,2),
00068 range(1,3):range(2,3))
00069
00070
00071
00072
00073 Real, Intent (In) :: coords1 (shape(1,1):shape(2,1),
00074 shape(1,2):shape(2,2),
00075 shape(1,3):shape(2,3))
00076
00077
00078
00079 Integer, Intent (In) :: control (2, ndim_3d)
00080
00081
00082
00083
00084 Integer, Intent (In) :: method_id
00085
00086
00087
00088 Integer, Intent (In) :: coords_shape (2)
00089
00090
00091
00092 Real, Intent (In) :: x_coords(coords_shape(1):coords_shape(2))
00093
00094
00095
00096 Integer, Intent (In) :: grid_valid_shape (2)
00097
00098
00099
00100
00101 Logical, Intent (In) :: cyclic
00102
00103
00104
00105 Real, Intent (In) :: chmin (grid_valid_shape(1)-1:
00106 grid_valid_shape(2))
00107
00108
00109
00110
00111
00112
00113
00114
00115
00116
00117
00118
00119
00120
00121
00122
00123
00124
00125
00126 Real, Intent (In) :: chmax (grid_valid_shape(1)-1:
00127 grid_valid_shape(2))
00128
00129
00130
00131 Real, Intent (In) :: tol
00132
00133
00134
00135
00136
00137
00138 integer, Intent (Out) :: ierror
00139
00140
00141
00142
00143
00144
00145
00146
00147
00148
00149
00150
00151
00152 Integer, Parameter :: nrefd = 3
00153
00154 Integer, Parameter :: val_direct = 1
00155 Integer, Parameter :: val_coupler = -1
00156
00157
00158
00159 Integer, Parameter :: lev = 1
00160
00161
00162
00163
00164
00165 Integer :: i, j, k
00166 Integer :: ibegl, iendl, jpart, n
00167 Integer :: ic
00168
00169
00170
00171 Real :: dist (nrefd)
00172 Integer :: irange (2)
00173
00174 Integer :: nref
00175 Integer :: nmin (1)
00176
00177
00178
00179 Integer :: ii
00180
00181 #ifdef DEBUG
00182
00183
00184 Integer :: nfnd (0:3)
00185 #endif /* DEBUG */
00186
00187
00188
00189
00190
00191
00192
00193
00194
00195
00196
00197
00198
00199
00200
00201
00202
00203
00204
00205
00206 Character(len=len_cvs_string), save :: mycvs =
00207 '$Id: psmile_mg_method_1d_real.F90 2325 2010-04-21 15:00:07Z valcke $'
00208
00209
00210
00211
00212 ierror = 0
00213
00214 #ifdef VERBOSE
00215 print 9990, trim(ch_id), lev, control
00216
00217 call psmile_flushstd
00218 #endif /* VERBOSE */
00219
00220 #ifdef PRISM_ASSERTION
00221 if (tol < 0.0) then
00222 call psmile_assert (__FILE__, __LINE__, &
00223 "tol must be >= 0.0")
00224 endif
00225
00226 #ifdef FOO
00227
00228
00229 if (range(1) < shape (1) .or. shape (2) < range (2)) then
00230 print *, "range", range
00231 print *, "shape", shape
00232 call psmile_assert (__FILE__, __LINE__, &
00233 "range must be within shape")
00234 endif
00235
00236
00237
00238 if (control(1) < range (1) .or. range (2) < control (2)) then
00239 print *, "control", control
00240 print *, "range ", range
00241 call psmile_assert (__FILE__, __LINE__, &
00242 "control must be within range")
00243 endif
00244 #endif /* FOO */
00245 #endif
00246
00247
00248
00249
00250
00251
00252
00253
00254
00255 #ifdef DEBUG
00256 nfnd (:) = 0
00257 #endif /* DEBUG */
00258
00259 #ifdef DEBUGX
00260 print *, 'begin ---'
00261 do k = control (1, 3), control (2,3)
00262 do j = control(1,2), control (2,2)
00263 do i = control(1,1), control (2,1)
00264 print *, 'i,j,k,coord, loc, found', i,j,k, coords1 (i,j,k), &
00265 loc(i,j,k), found (i,j,k)
00266 end do
00267 end do
00268 end do
00269 #endif
00270
00271
00272
00273 do k = control (1, 3), control (2,3)
00274 do j = control(1,2), control (2,2)
00275
00276
00277
00278
00279
00280
00281
00282
00283 ibegl = control (1, 1)
00284
00285
00286
00287 do while (ibegl <= control(2,1))
00288
00289 do i = ibegl, control(2,1)
00290 if (found(i, j, k) == lev) exit
00291 end do
00292
00293 if (i > control(2,1)) exit
00294 ibegl = i
00295
00296
00297 do i = ibegl, control(2,1)
00298 if (found(i, j, k) /= lev) exit
00299 end do
00300
00301 iendl = i - 1
00302
00303
00304
00305 #ifdef PRISM_ASSERTION
00306 do i = ibegl, iendl
00307 if (loc (i,j,k) < coords_shape(1) .or. &
00308 loc (i,j,k) > coords_shape(2)) exit
00309 end do
00310
00311 if (i <= iendl) then
00312 print *, "i,j,k, ic, coords_shape", &
00313 i,j,k, loc (i,j,k), coords_shape
00314 print *, 'control', control
00315 print *, 'range ', range
00316 call psmile_assert (__FILE__, __LINE__, &
00317 'wrong index')
00318 endif
00319 #endif
00320
00321 do jpart = 0, (iendl-ibegl)/5000
00322
00323 do i = ibegl+jpart*5000, min (iendl, ibegl+(jpart+1)*5000 - 1)
00324
00325 ic = loc (i,j,k)
00326
00327 irange (1) = max (ic-1, grid_valid_shape (1))
00328 irange (2) = min (ic+1, grid_valid_shape (2))
00329
00330
00331
00332 do n = irange(1), irange (2)
00333 dist (n-irange(1)+1) = abs (x_coords (n)-coords1(i,j,k))
00334 end do
00335
00336
00337
00338 nref = irange (2) - irange (1) + 1
00339 #ifdef PRISM_ASSERTION
00340 if (nref <= 0) then
00341 call psmile_assert (__FILE__, __LINE__, &
00342 'nref <= 0')
00343 endif
00344 #endif
00345 nmin = MINLOC (dist(1:nref))
00346 #ifdef MINLOCFIX
00347 if (nmin(1).eq.0) nmin=1
00348 #endif /* MINLOCFIX */
00349
00350 loc (i,j,k) = irange (1) + nmin(1) - 1
00351
00352 #ifdef DEBUG2
00353 print *, "psmile_mg_method_1d_real: i, j, k --- ", i,j,k
00354 print *, "psmile_mg_method_1d_real loc ", loc(i,j,k), nmin(1), dist(nmin(1))
00355 print *, "psmile_mg_method_1d_real ic ", ic
00356 print *, "psmile_mg_method_1d_real irange", irange
00357 print *, "psmile_mg_method_1d_real search", coords1(i,j,k)
00358 print *, "psmile_mg_method_1d_real source", x_coords (irange(1):irange(2))
00359 #endif
00360
00361 if (dist(nmin(1)) .lt. tol) then
00362
00363
00364
00365
00366
00367
00368 #ifdef DEBUG
00369 nfnd (1) = nfnd (1) + 1
00370 #endif /* DEBUG */
00371 else
00372 found (i,j,k) = val_coupler
00373 #ifdef DEBUG
00374 nfnd (2) = nfnd (2) + 1
00375 #endif /* DEBUG */
00376
00377 if (coords1(i,j,k) < chmin (loc(i,j,k)) .or. &
00378 coords1(i,j,k) > chmax (loc(i,j,k))) then
00379
00380 irange (1) = max (ic-1, grid_valid_shape (1)-1)
00381
00382 do ii = irange (1), irange (2)
00383 if (chmin(ii) <= coords1(i,j,k) .and. &
00384 chmax(ii) >= coords1(i,j,k)) exit
00385 end do
00386
00387 if (ii <= irange(2)) then
00388
00389 loc (i,j,k) = ii
00390 else
00391
00392
00393 #ifdef DEBUG
00394 nfnd (3) = nfnd (3) + 1
00395 nfnd (2) = nfnd (2) - 1
00396 #endif /* DEBUG */
00397
00398
00399 endif
00400 endif
00401 endif
00402
00403 end do
00404 end do
00405
00406
00407
00408 if (cyclic) then
00409
00410
00411
00412
00413
00414
00415 do i = ibegl, iendl
00416 if (loc (i,j,k) < grid_valid_shape (1)) then
00417 loc (i,j,k) = grid_valid_shape (2)
00418 endif
00419 enddo
00420 #ifdef REMOVED
00421 else
00422
00423
00424
00425
00426
00427
00428
00429 do i = ibegl, iendl
00430 loc (i,j,k) = max (loc(i,j,k), grid_valid_shape (1))
00431 enddo
00432 #endif
00433
00434 endif
00435
00436
00437
00438 ibegl = iendl + 2
00439 end do
00440
00441 end do
00442 end do
00443
00444
00445
00446 #ifdef DEBUGX
00447 print *, 'end ---'
00448 do k = control (1, 3), control (2,3)
00449 do j = control(1,2), control (2,2)
00450 do i = control(1,1), control (2,1)
00451 print *, 'i,j,k, coords, loc, found', i,j,k, coords1 (i,j,k), &
00452 loc(i,j,k), found (i,j,k)
00453 end do
00454 end do
00455 end do
00456 #endif
00457
00458 #ifdef DEBUG
00459 print 9970, trim(ch_id), lev, nfnd
00460 9970 format (1x, a, ': psmile_mg_method_1d_real: lev =', i3, &
00461 ': not :', i5, ', dir: ', i6, ', fnd :', i6, i5)
00462 #endif /* DEBUG */
00463
00464 #ifdef VERBOSE
00465 print 9980, trim(ch_id), lev
00466
00467 call psmile_flushstd
00468 #endif /* VERBOSE */
00469
00470
00471
00472
00473 9990 format (1x, a, ': psmile_mg_method_1d_real: level', i3, &
00474 ', control', 6i6)
00475 9980 format (1x, a, ': psmile_mg_method_1d_real: eof, level', i3)
00476
00477 end subroutine psmile_mg_method_1d_real