00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012 subroutine psmile_mg_method_2d_real ( &
00013 comp_info, nlev, found, loc, range, &
00014 coords1, coords2, search_shape, control, &
00015 x_coords, y_coords, coords_shape, &
00016 grid_valid_shape, &
00017 cyclic, period, &
00018 chmin1, chmin2, chmax1, chmax2, &
00019 tol, ierror)
00020
00021
00022
00023 use prism_constants
00024 use psmile, dummy_interface => PSMILe_mg_method_2d_real
00025 use psmile_geometry, only : point_is_in_quadrangle
00026 use psmile_grid, only : pole_threshold, psmile_transform_cell_cyclic
00027 #ifdef DEBUG_TRACE
00028 use psmile_debug_trace
00029 #endif
00030
00031 implicit none
00032
00033
00034
00035
00036
00037 Type (Enddef_comp), Intent (In) :: comp_info
00038
00039
00040
00041 Integer, Intent (In) :: nlev
00042
00043
00044
00045 Integer, Intent (In) :: range (2, ndim_3d)
00046
00047
00048
00049 Integer, Intent (In) :: search_shape (2, ndim_3d)
00050
00051
00052
00053
00054 Real, Intent (In) :: coords1(search_shape(1,1):search_shape(2,1),
00055 search_shape(1,2):search_shape(2,2),
00056 search_shape(1,3):search_shape(2,3))
00057
00058 Real, Intent (In) :: coords2(search_shape(1,1):search_shape(2,1),
00059 search_shape(1,2):search_shape(2,2),
00060 search_shape(1,3):search_shape(2,3))
00061
00062
00063 Integer, Intent (In) :: control (2, ndim_3d)
00064
00065
00066
00067 Integer, Intent (In) :: coords_shape (2, ndim_2d)
00068
00069
00070
00071 Real, Intent (In) :: x_coords(coords_shape(1,1):coords_shape(2,1),
00072 coords_shape(1,2):coords_shape(2,2))
00073
00074 Real, Intent (In) :: y_coords(coords_shape(1,1):coords_shape(2,1),
00075 coords_shape(1,2):coords_shape(2,2))
00076
00077
00078
00079 Integer, Intent (In) :: grid_valid_shape (2, ndim_2d)
00080
00081
00082
00083
00084 Logical, Intent (In) :: cyclic (ndim_2d)
00085
00086
00087
00088 Real, Intent (In) :: period (ndim_2d)
00089
00090
00091
00092 Real, Intent (In) :: chmin1 (grid_valid_shape(1,1)-1:grid_valid_shape(2,1),
00093 grid_valid_shape(1,2)-1:grid_valid_shape(2,2))
00094
00095 Real, Intent (In) :: chmin2 (grid_valid_shape(1,1)-1:grid_valid_shape(2,1),
00096 grid_valid_shape(1,2)-1:grid_valid_shape(2,2))
00097
00098
00099
00100
00101
00102
00103
00104
00105
00106
00107
00108
00109
00110
00111
00112
00113
00114
00115
00116
00117 Real, Intent (In) :: chmax1 (grid_valid_shape(1,1)-1:grid_valid_shape(2,1),
00118 grid_valid_shape(1,2)-1:grid_valid_shape(2,2))
00119 Real, Intent (In) :: chmax2 (grid_valid_shape(1,1)-1:grid_valid_shape(2,1),
00120 grid_valid_shape(1,2)-1:grid_valid_shape(2,2))
00121
00122
00123
00124
00125 Real, Intent (In) :: tol
00126
00127
00128
00129
00130
00131 Integer, Intent (InOut) :: found (range(1,1):range(2,1),
00132 range(1,2):range(2,2),
00133 range(1,3):range(2,3))
00134
00135
00136
00137
00138
00139
00140
00141
00142
00143
00144
00145
00146
00147 Integer, Intent (InOut) :: loc (ndim_2d,
00148 range(1,1):range(2,1),
00149 range(1,2):range(2,2),
00150 range(1,3):range(2,3))
00151
00152
00153
00154
00155
00156
00157 Integer, Intent (Out) :: ierror
00158
00159
00160
00161
00162
00163
00164
00165
00166
00167
00168
00169
00170
00171 Integer, Parameter :: nrefd = 3 * 3
00172 Integer, Parameter :: n_corners = 4
00173
00174 Integer, Parameter :: val_direct = 1
00175 Integer, Parameter :: val_coupler = -1
00176 Integer :: val_notfound
00177
00178 Integer, Parameter :: temp_len = 256
00179
00180
00181
00182 Integer, Parameter :: lev = 1
00183
00184
00185
00186
00187
00188 Integer :: i, j, k, l
00189 Integer :: ibegl, jpart, n
00190 Integer :: curr_loc_list (temp_len, ndim_2d)
00191
00192
00193
00194 Integer :: ijkctl (ndim_2d, nrefd)
00195 Integer :: ijkdst (temp_len, ndim_2d, nrefd)
00196 Integer :: irange (temp_len, 2, ndim_2d)
00197
00198 Integer :: ii, jj, nref (temp_len)
00199
00200
00201
00202 Logical :: cyclic_req
00203
00204 Real :: period2 (ndim_2d)
00205
00206
00207
00208 Type (point_real) :: xy_source(temp_len, 4)
00209 Type (point_real) :: xy_source_point(temp_len)
00210 Type (point_real), target :: xy_target(temp_len)
00211 Type (point_real), target :: xy_target_rotated(temp_len)
00212 Type p_xy_target
00213 Type (point_real), pointer :: p
00214 end Type p_xy_target
00215 Type(p_xy_target) :: xy_target_pointer(temp_len)
00216
00217
00218
00219 Integer :: num_curr_locs, nc_full, np, nc_part
00220 Integer :: list (temp_len*2), full_list (temp_len)
00221 Integer :: part_list (temp_len)
00222 Logical :: fnd (temp_len)
00223
00224 #ifdef DEBUG
00225
00226
00227 Integer :: nfnd (0:3)
00228 #endif /* DEBUG */
00229
00230
00231
00232
00233
00234
00235
00236
00237
00238
00239
00240
00241
00242
00243
00244
00245
00246
00247
00248
00249
00250 Character(len=len_cvs_string), save :: mycvs =
00251 '$Id: psmile_mg_method_2d_real.F90 2927 2011-01-28 14:04:12Z hanke $'
00252
00253
00254
00255
00256
00257 data ((ijkctl (i, n), i=1,ndim_2d), n = 1,nrefd) &
00258 / -1,-1, 0,-1, 1,-1, &
00259 -1, 0, 0, 0, 1, 0, &
00260 -1, 1, 0, 1, 1, 1/
00261
00262
00263
00264
00265 ierror = 0
00266 val_notfound = -(nlev+1)
00267
00268 #ifdef VERBOSE
00269 print 9990, trim(ch_id), lev, control
00270
00271 call psmile_flushstd
00272 #endif /* VERBOSE */
00273
00274 #ifdef DEBUG_TRACE
00275 if ( range(1,1) <= ictl_ind(1) .and. range(2,1) >= ictl_ind(1) .and. &
00276 range(1,2) <= ictl_ind(2) .and. range(2,2) >= ictl_ind(2)) then
00277 print 8990, trim(ch_id), ictl_ind (1:ndim_2d), &
00278 found ( ictl_ind(1),ictl_ind(2),control(1,3)), &
00279 loc (:, ictl_ind(1),ictl_ind(2),control(1,3))
00280 endif
00281 #endif /* DEBUG_TRACE */
00282
00283 #ifdef PRISM_ASSERTION
00284 if (tol < 0.0) then
00285 call psmile_assert (__FILE__, __LINE__, &
00286 "tol must be >= 0.0")
00287 endif
00288
00289 if ( grid_valid_shape(1,1)-1 < coords_shape (1,1) .or. &
00290 grid_valid_shape(2,1)+1 > coords_shape (2,1) .or. &
00291 grid_valid_shape(1,2)-1 < coords_shape (1,2) .or. &
00292 grid_valid_shape(2,2)+1 > coords_shape (2,2)) then
00293
00294 print *, 'grid_valid_shape:', grid_valid_shape
00295 print *, 'coords_shape :', coords_shape
00296 call psmile_assert (__FILE__, __LINE__, &
00297 "insufficient coords_shape")
00298 endif
00299 #endif
00300
00301
00302
00303
00304 cyclic_req = cyclic (1) .or. cyclic (2)
00305 period2 (:) = period (:) * 0.5
00306
00307 #ifdef DEBUG
00308 nfnd (:) = 0
00309 #endif /* DEBUG */
00310
00311
00312
00313
00314 do k = control(1,3), control(2,3)
00315 do j = control(1,2), control (2,2)
00316
00317
00318
00319
00320
00321
00322
00323
00324
00325 ibegl = control (1, 1)
00326
00327
00328
00329
00330 loop_i: do while (ibegl <= control(2,1))
00331
00332 num_curr_locs = 0
00333 i = ibegl
00334
00335 do while (num_curr_locs < temp_len .and. i <= control(2,1))
00336 ibegl = i
00337
00338 do i = ibegl, min (ibegl+(temp_len*2-1)-num_curr_locs, control(2,1))
00339 if (found(i, j, k) == lev) then
00340 num_curr_locs = num_curr_locs + 1
00341 list (num_curr_locs) = i
00342 endif
00343 end do
00344 end do
00345
00346 if (num_curr_locs == 0) exit loop_i
00347
00348 if (num_curr_locs < temp_len) then
00349 ibegl = control(2,1) + 1
00350 else if (num_curr_locs > temp_len) then
00351 num_curr_locs = temp_len
00352 ibegl = list (temp_len+1)
00353 else
00354 ibegl = list (temp_len) + 1
00355 endif
00356
00357 if (cyclic_req) then
00358 nc_full = num_curr_locs
00359 full_list (1:num_curr_locs) = list (1:num_curr_locs)
00360 else
00361 nc_full = 0
00362 endif
00363
00364
00365
00366
00367
00368 curr_loc_list (1:num_curr_locs, 1) = loc (1, list(1:num_curr_locs), j, k)
00369 curr_loc_list (1:num_curr_locs, 2) = loc (2, list(1:num_curr_locs), j, k)
00370
00371 #ifdef PRISM_ASSERTION
00372
00373 do i = 1, num_curr_locs
00374 if ( curr_loc_list(i,1) < coords_shape(1,1) .or. &
00375 curr_loc_list(i,1) > coords_shape(2,1) .or. &
00376 curr_loc_list(i,2) < coords_shape(1,2) .or. &
00377 curr_loc_list(i,2) > coords_shape(2,2) ) exit
00378 end do
00379
00380 if (i <= num_curr_locs) then
00381 print *, 'i, curr_loc_list', i, curr_loc_list(i, :)
00382 print *, 'coords_shape', coords_shape
00383 call psmile_assert (__FILE__, __LINE__, &
00384 'wrong index')
00385 endif
00386 #endif
00387
00388
00389 do i = 1, num_curr_locs
00390
00391
00392 irange (i, 1, :) = max (curr_loc_list(i,:)-1, grid_valid_shape (1, :))
00393 irange (i, 2, :) = min (curr_loc_list(i,:)+1, grid_valid_shape (2, :))
00394 end do
00395
00396
00397 nref (1:num_curr_locs) = &
00398 (irange (1:num_curr_locs, 2,1) - irange (1:num_curr_locs, 1,1) + 1) * &
00399 (irange (1:num_curr_locs, 2,2) - irange (1:num_curr_locs, 1,2) + 1)
00400
00401 #ifdef PRISM_ASSERTION
00402 if (minval (nref(1:num_curr_locs)) <= 0) then
00403 call psmile_assert (__FILE__, __LINE__, &
00404 'nref <= 0')
00405 endif
00406 #endif
00407
00408 do i = 1, num_curr_locs
00409
00410 if (nref(i) == nrefd) then
00411 do n = 1, nrefd
00412
00413 ijkdst (i, :, n) = curr_loc_list (i, :) + ijkctl (:, n)
00414 end do
00415 else
00416
00417 n = 0
00418 do jj = irange (i, 1, 2), irange (i, 2, 2)
00419 do ii = irange (i, 1, 1), irange (i, 2, 1)
00420 n = n + 1
00421 ijkdst (i, 1, n) = ii
00422 ijkdst (i, 2, n) = jj
00423 end do
00424 end do
00425 endif
00426 end do
00427
00428
00429
00430 do i = 1, num_curr_locs
00431 xy_target(i)%x = coords1(list(i),j,k)
00432 xy_target(i)%y = coords2(list(i),j,k)
00433 enddo
00434
00435 nc_part = 0
00436
00437
00438
00439
00440 do n = 1, nrefd
00441
00442 do i = 1, num_curr_locs
00443 if (n <= nref(i)) then
00444 xy_source_point%x = x_coords (ijkdst(i,1,n), ijkdst (i,2,n))
00445 xy_source_point%y = y_coords (ijkdst(i,1,n), ijkdst (i,2,n))
00446
00447
00448 if ( y_coords (ijkdst(i,1,n), ijkdst (i,2,n)) > pole_threshold .or. &
00449 y_coords (ijkdst(i,1,n), ijkdst (i,2,n)) < -pole_threshold ) then
00450
00451 call psmile_transrot_real ( &
00452 x_coords (ijkdst(i,1,n), ijkdst (i,2,n)), &
00453 y_coords (ijkdst(i,1,n), ijkdst (i,2,n)), &
00454 xy_source_point(i)%x, xy_source_point(i)%y )
00455
00456 xy_target_rotated(i) = xy_target(i)
00457
00458 call psmile_transrot_real ( &
00459 coords1(list(i),j,k), coords2(list(i),j,k), &
00460 xy_target_rotated(i)%x, xy_target_rotated(i)%y )
00461
00462 xy_target_pointer(i)%p => xy_target_rotated(i)
00463 else
00464 xy_source_point(i)%x = x_coords (ijkdst(i,1,n), ijkdst (i,2,n))
00465 xy_source_point(i)%y = y_coords (ijkdst(i,1,n), ijkdst (i,2,n))
00466 xy_target_pointer(i)%p => xy_target(i)
00467 endif
00468
00469
00470 if (abs(xy_target_pointer(i)%p%x - xy_source_point(i)%x) < tol .and. &
00471 abs(xy_target_pointer(i)%p%y - xy_source_point(i)%y) < tol) then
00472
00473 loc (:, list(i),j,k) = ijkdst (i, :, n)
00474 nc_part = nc_part + 1
00475 part_list (nc_part) = i
00476 endif
00477 end if
00478 end do
00479 end do
00480
00481
00482 if (nc_part > 0) then
00483
00484 #ifdef DEBUG
00485 nfnd (1) = nfnd (1) + nc_part
00486 #endif /* DEBUG */
00487
00488
00489
00490 if (num_curr_locs == nc_part) cycle loop_i
00491
00492
00493
00494
00495 do np = 0, nc_part-1
00496 i = part_list (nc_part - np)
00497
00498 if (i /= num_curr_locs-np) then
00499 list (i) = list (num_curr_locs-np)
00500 nref (i) = nref (num_curr_locs-np)
00501 xy_target(i) = xy_target(num_curr_locs-np)
00502
00503 curr_loc_list (i, :) = curr_loc_list (num_curr_locs-np, :)
00504
00505 ijkdst (i, :, :) = ijkdst (num_curr_locs-np, :, :)
00506 irange (i, :, :) = irange (num_curr_locs-np, :, :)
00507 endif
00508 end do
00509
00510 num_curr_locs = num_curr_locs - nc_part
00511 endif
00512
00513
00514
00515
00516
00517
00518
00519 do i = 1, num_curr_locs
00520 found (list(i), j, k) = val_coupler
00521 end do
00522
00523
00524
00525
00526
00527 nc_part = 0
00528
00529 do i = 1, num_curr_locs
00530 if ( curr_loc_list(i,1) == grid_valid_shape(1,1) .or. &
00531 curr_loc_list(i,2) == grid_valid_shape(1,2)) then
00532 nc_part = nc_part + 1
00533 part_list (nc_part) = i
00534 endif
00535 end do
00536
00537
00538 do np = 1, nc_part
00539 i = part_list (np)
00540
00541
00542
00543
00544
00545
00546 irange (i,1, :) = max (curr_loc_list(i,:)-1, grid_valid_shape(1,:)-1)
00547 irange (i,2, :) = min (curr_loc_list(i,:)+1, grid_valid_shape(2,:))
00548
00549 nref(i) = (irange (i,2,1) - irange (i,1,1) + 1) * &
00550 (irange (i,2,2) - irange (i,1,2) + 1)
00551
00552 if (nref(i) == nrefd) then
00553 do n = 1, nrefd
00554 ijkdst (i, :, n) = curr_loc_list (i, :) + ijkctl (:, n)
00555 end do
00556 else
00557 n = 0
00558 do jj = irange (i,1, 2), irange (i,2, 2)
00559 do ii = irange (i,1, 1), irange (i,2, 1)
00560 n = n + 1
00561 ijkdst (i, 1, n) = ii
00562 ijkdst (i, 2, n) = jj
00563 end do
00564 end do
00565 endif
00566 end do
00567
00568 #ifdef PRISM_ASSERTION
00569
00570 do np = 1, nc_part
00571 if (nref(part_list (np)) <= 0) exit
00572 end do
00573
00574 if (np <= nc_part) then
00575 call psmile_assert (__FILE__, __LINE__, &
00576 'nref <= 0')
00577 endif
00578 #endif
00579
00580
00581
00582 fnd (1:num_curr_locs) = .false.
00583
00584 do n = 1, nrefd
00585
00586 do i = 1, num_curr_locs
00587
00588
00589
00590
00591
00592
00593
00594 if ( n <= nref (i) .and. .not. fnd(i) ) then
00595
00596
00597 do while ( xy_target(i)%x < chmax1(ijkdst(i, 1,n),ijkdst(i, 2,n)) - period(1) )
00598 xy_target(i)%x = xy_target(i)%x + period(1)
00599 enddo
00600
00601 do while ( xy_target(i)%x > chmin1(ijkdst(i, 1,n),ijkdst(i, 2,n)) + period(1) )
00602 xy_target(i)%x = xy_target(i)%x - period(1)
00603 enddo
00604
00605
00606
00607 if ( xy_target(i)%x >= chmin1(ijkdst(i,1,n), ijkdst(i,2,n)) .and. &
00608 xy_target(i)%y >= chmin2(ijkdst(i,1,n), ijkdst(i,2,n)) .and. &
00609 xy_target(i)%x <= chmax1(ijkdst(i,1,n), ijkdst(i,2,n)) .and. &
00610 xy_target(i)%y <= chmax2(ijkdst(i,1,n), ijkdst(i,2,n))) then
00611
00612 #ifdef DEBUG_TRACE
00613 if (j == ictl_ind(2)) then
00614 if (list(i) == ictl_ind(1)) then
00615 print *, 'is within bounding box:'
00616 print *, 'ch1', chmin1(ijkdst(i,1,n), ijkdst(i,2,n)), &
00617 chmax1(ijkdst(i,1,n), ijkdst(i,2,n))
00618 print *, 'ch2', chmin2(ijkdst(i,1,n), ijkdst(i,2,n)), &
00619 chmax2(ijkdst(i,1,n), ijkdst(i,2,n))
00620 endif
00621 end if
00622 #endif /* DEBUG_TRACE */
00623
00624
00625 if ( chmin2(ijkdst(i,1,n),ijkdst(i,2,n)) > pole_threshold .or. &
00626 chmax2(ijkdst(i,1,n),ijkdst(i,2,n)) < -pole_threshold ) then
00627
00628
00629
00630 call psmile_transrot_real ( &
00631 x_coords (ijkdst(i,1,n) , ijkdst(i,2,n) ), &
00632 y_coords (ijkdst(i,1,n) , ijkdst(i,2,n) ), &
00633 xy_source(i, 1)%x, xy_source(i, 1)%y )
00634
00635 call psmile_transrot_real ( &
00636 x_coords (ijkdst(i,1,n)+1, ijkdst(i,2,n) ), &
00637 y_coords (ijkdst(i,1,n)+1, ijkdst(i,2,n) ), &
00638 xy_source(i, 2)%x, xy_source(i, 2)%y )
00639
00640 call psmile_transrot_real ( &
00641 x_coords (ijkdst(i,1,n)+1, ijkdst(i,2,n)+1), &
00642 y_coords (ijkdst(i,1,n)+1, ijkdst(i,2,n)+1), &
00643 xy_source(i, 3)%x, xy_source(i, 3)%y )
00644
00645 call psmile_transrot_real ( &
00646 x_coords (ijkdst(i,1,n) , ijkdst(i,2,n)+1), &
00647 y_coords (ijkdst(i,1,n) , ijkdst(i,2,n)+1), &
00648 xy_source(i, 4)%x, xy_source(i, 4)%y )
00649
00650 call psmile_transrot_real ( &
00651 coords1(list(i),j,k), coords2(list(i),j,k), &
00652 xy_target_rotated(i)%x, xy_target_rotated(i)%y )
00653
00654
00655 call psmile_transform_cell_cyclic ( xy_source(i, :)%x, period(1), ierror )
00656
00657
00658 do while ( xy_target_rotated(i)%x < minval(xy_source(i, :)%x) - period2(1) )
00659 xy_target_rotated(i)%x = xy_target_rotated(i)%x + period(1)
00660 enddo
00661
00662 do while ( xy_target_rotated(i)%x > maxval(xy_source(i, :)%x) + period2(1) )
00663 xy_target_rotated(i)%x = xy_target_rotated(i)%x - period(1)
00664 enddo
00665
00666 xy_target_pointer(i)%p => xy_target_rotated(i)
00667
00668 #ifdef DEBUG_TRACE
00669 if (j == ictl_ind(2)) then
00670 if (list(i) == ictl_ind(1)) then
00671
00672 print *, 'point and cell are rotated to equator'
00673 print "(' point coordinates before: ', 2f8.2)", xy_target(i)
00674 print "(' cell coordinates before: ', 8f8.2)", &
00675 x_coords (ijkdst(i,1,n) , ijkdst(i,2,n) ), &
00676 y_coords (ijkdst(i,1,n) , ijkdst(i,2,n) ), &
00677 x_coords (ijkdst(i,1,n)+1, ijkdst(i,2,n) ), &
00678 y_coords (ijkdst(i,1,n)+1, ijkdst(i,2,n) ), &
00679 x_coords (ijkdst(i,1,n)+1, ijkdst(i,2,n)+1), &
00680 y_coords (ijkdst(i,1,n)+1, ijkdst(i,2,n)+1), &
00681 x_coords (ijkdst(i,1,n) , ijkdst(i,2,n)+1), &
00682 y_coords (ijkdst(i,1,n) , ijkdst(i,2,n)+1)
00683 endif
00684 end if
00685 #endif /* DEBUG_TRACE */
00686 else
00687
00688 xy_source(i, 1)%x = x_coords (ijkdst(i,1,n) , ijkdst(i,2,n) )
00689 xy_source(i, 1)%y = y_coords (ijkdst(i,1,n) , ijkdst(i,2,n) )
00690
00691 xy_source(i, 2)%x = x_coords (ijkdst(i,1,n)+1, ijkdst(i,2,n) )
00692 xy_source(i, 2)%y = y_coords (ijkdst(i,1,n)+1, ijkdst(i,2,n) )
00693
00694 xy_source(i, 3)%x = x_coords (ijkdst(i,1,n)+1, ijkdst(i,2,n)+1)
00695 xy_source(i, 3)%y = y_coords (ijkdst(i,1,n)+1, ijkdst(i,2,n)+1)
00696
00697 xy_source(i, 4)%x = x_coords (ijkdst(i,1,n) , ijkdst(i,2,n)+1)
00698 xy_source(i, 4)%y = y_coords (ijkdst(i,1,n) , ijkdst(i,2,n)+1)
00699
00700 xy_target_pointer(i)%p => xy_target(i)
00701
00702
00703
00704 if (any(xy_source(i, :)%x < chmin1(ijkdst(i,1,n),ijkdst(i,2,n))) .or. &
00705 any(xy_source(i, :)%x > chmax1(ijkdst(i,1,n),ijkdst(i,2,n)))) then
00706 call psmile_transform_cell_cyclic ( xy_source(i, :)%x, period(1), ierror )
00707 endif
00708
00709 endif
00710
00711
00712 if (point_is_in_quadrangle(xy_target_pointer(i)%p, xy_source(i,:))) then
00713 loc(1,list(i),j,k) = ijkdst(i,1,n)
00714 loc(2,list(i),j,k) = ijkdst(i,2,n)
00715 fnd (i) = .true.
00716 endif
00717 #ifdef DEBUG_TRACE
00718 if (j == ictl_ind(2)) then
00719 if (list(i) == ictl_ind(1)) then
00720
00721 print "(' point coordinates: ', 2f8.2)", xy_target_pointer(i)%p
00722 print "(' cell coordinates: ', 8f8.2)", xy_source(i,:)
00723
00724 if (fnd (i)) then
00725 print *, 'point is within this cell'
00726 else
00727 print *, 'point is not within this cell'
00728 endif
00729 endif
00730 end if
00731 else
00732 if (j == ictl_ind(2)) then
00733 if (list(i) == ictl_ind(1)) then
00734 print *, 'is not within bounding box:'
00735 print *, 'ch1', chmin1(ijkdst(i,1,n), ijkdst(i,2,n)), &
00736 chmax1(ijkdst(i,1,n), ijkdst(i,2,n))
00737 print *, 'ch2', chmin2(ijkdst(i,1,n), ijkdst(i,2,n)), &
00738 chmax2(ijkdst(i,1,n), ijkdst(i,2,n))
00739 endif
00740 end if
00741 #endif /* DEBUG_TRACE */
00742 endif
00743 endif
00744 end do
00745 end do
00746
00747
00748
00749 do n = 1, nrefd
00750
00751 do i = 1, num_curr_locs
00752 if ( n <= nref (i) .and. .not. fnd(i) ) then
00753 found(list(i), j, k) = val_notfound
00754 endif
00755 enddo
00756 enddo
00757
00758 #ifdef DEBUG_TRACE
00759 if (j == ictl_ind(2)) then
00760 do i = 1, num_curr_locs
00761 if (list(i) == ictl_ind(1)) exit
00762 end do
00763
00764 if (i <= num_curr_locs) then
00765 print *, 'mg_method: ictl_ind', ictl_ind(1:ndim_2d), &
00766 ': ', fnd(i), coords1(list(i),j,k), coords2(list(i),j,k)
00767 print *, ' .... found, i, j: ', fnd(i)
00768 if ( fnd(i) ) then
00769 print *, ' .... found, i, j: ', loc(1:2,list(i),j,k)
00770 else
00771 print *, ' did not find new i, j!'
00772 endif
00773 endif
00774 end if
00775 #endif /* DEBUG_TRACE */
00776
00777
00778
00779
00780
00781
00782
00783 do l = 1, ndim_2d
00784 if (cyclic(l)) then
00785
00786
00787
00788
00789
00790 do jpart = 1, num_curr_locs, 5000
00791
00792 do i = jpart, min (num_curr_locs, jpart+5000-1)
00793 if (loc (l, list(i),j,k) < grid_valid_shape (1,l) &
00794 .and. fnd(i)) then
00795 loc (l, list(i),j,k) = grid_valid_shape (2,l)
00796 endif
00797 enddo
00798
00799 enddo
00800 endif
00801 end do
00802
00803
00804
00805 end do loop_i
00806 end do
00807 end do
00808
00809
00810
00811
00812 #ifdef DEBUG_TRACE
00813 if ( range(1,1) <= ictl_ind(1) .and. range(2,1) >= ictl_ind(1) .and. &
00814 range(1,2) <= ictl_ind(2) .and. range(2,2) >= ictl_ind(2)) then
00815
00816 print *, 'coords', coords1(ictl_ind(1),ictl_ind(2),control(1,3)), &
00817 coords2(ictl_ind(1),ictl_ind(2),control(1,3))
00818
00819 print 8990, trim(ch_id), ictl_ind(1:ndim_2d), &
00820 found ( ictl_ind(1),ictl_ind(2),control(1,3)), &
00821 loc (:, ictl_ind(1),ictl_ind(2),control(1,3))
00822 endif
00823
00824 8990 format (1x, a, ': psmile_mg_method_2d_real: ictl_ind', 2i5, &
00825 '; found', i3, ' loc ', 2i6)
00826 #endif /* DEBUG_TRACE */
00827
00828 #ifdef DEBUG
00829 print 9970, trim(ch_id), lev, nfnd
00830 9970 format (1x, a, ': psmile_mg_method_2d_real: lev =', i3, &
00831 ': not :', i5, ', dir: ', i6, ', fnd :', i6, i5)
00832
00833 #ifdef DEBUG_PRINT
00834 call psmile_print_found_real (comp_info, lev, &
00835 found, loc, range, &
00836 coords1, coords2, coords1, &
00837 search_shape, search_shape, search_shape, &
00838 control, ndim_2d, "psmile_mg_method_2d_real: eof")
00839
00840 #endif /* DEBUG_PRINT */
00841 #endif /* DEBUG */
00842
00843 #ifdef VERBOSE
00844 print 9980, trim(ch_id), lev
00845
00846 call psmile_flushstd
00847 #endif /* VERBOSE */
00848
00849
00850
00851
00852 9990 format (1x, a, ': psmile_mg_method_2d_real: level', i3, &
00853 ', control', 6i6)
00854 9980 format (1x, a, ': psmile_mg_method_2d_real: eof, level', i3)
00855
00856 end subroutine psmile_mg_method_2d_real