psmile_mg_method_2d_dble.F90

Go to the documentation of this file.
00001 !-----------------------------------------------------------------------
00002 ! Copyright 2006-2010, NEC Europe Ltd., London, UK.
00003 ! Copyright 2010, DKRZ, Hamburg, Germany.
00004 ! All rights reserved. Use is subject to OASIS4 license terms.
00005 !-----------------------------------------------------------------------
00006 !BOP
00007 !
00008 ! !ROUTINE: PSMILe_MG_method_2d_dble
00009 !
00010 ! !INTERFACE:
00011 
00012 subroutine psmile_mg_method_2d_dble (                    &
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 ! !USES:
00022 !
00023   use prism_constants
00024   use psmile, dummy_interface => PSMILe_mg_method_2d_dble
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 ! !INPUT PARAMETERS:
00035 !
00036 
00037   Type (Enddef_comp), Intent (In) :: comp_info
00038 
00039 ! Info on the component in which the donor should be searched.
00040 
00041   Integer, Intent (In)            :: nlev
00042 
00043 ! Total number of levels
00044 
00045   Integer, Intent (In)            :: range (2, ndim_3d)
00046 
00047 ! Dimension of loc and found
00048 
00049   Integer, Intent (In)            :: search_shape (2, ndim_3d)
00050 
00051 ! Dimension of coordinate arrays "coord1" and "coord2"
00052 ! which contain the coordinates to be searched.
00053 
00054   Double Precision, 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   Double Precision, 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 ! Coordinates to be searched
00062 
00063   Integer, Intent (In)            :: control (2, ndim_3d)
00064 
00065 ! Index range in "coords" to be searched
00066 
00067   Integer, Intent (In)            :: coords_shape (2, ndim_2d)
00068 
00069 ! Dimension of coordinates (method) x_coords, ...
00070 
00071   Double Precision, Intent (In)   :: x_coords(coords_shape(1,1):coords_shape(2,1), 
00072                                               coords_shape(1,2):coords_shape(2,2))
00073 
00074   Double Precision, Intent (In)   :: y_coords(coords_shape(1,1):coords_shape(2,1), 
00075                                               coords_shape(1,2):coords_shape(2,2))
00076 
00077 ! Coordinates of the method
00078 
00079   Integer, Intent (In)            :: grid_valid_shape (2, ndim_2d)
00080 
00081 ! Specifies the valid block shape for the "ndim_2d"-dimensional block
00082 ! Dimension of chmin, chmax and midp
00083 
00084   Logical,          Intent (In)   :: cyclic (ndim_2d)
00085 
00086 ! Is this a (locally) cyclic coordinate in corresponding direction
00087 
00088   Double Precision, Intent (In)   :: period (ndim_2d)
00089 !
00090 ! Period of a cyclic coordinate.
00091 !
00092   Double Precision, 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   Double Precision, 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 ! Minimum of grid coordinates per grid cell
00099 ! chmin (i, j) contains the minimal value of bounding box of
00100 !  method cell ((i,j), (i+1,j), (i,j+1), (i+1,j+1))
00101 !
00102 ! The first box
00103 !     chmin (grid_valid_shape(1,1)-1,j) contains the minimal value of
00104 !                                   the bounding box of virtual cell
00105 !       (array(grid_valid_shape(1,1), j), array(grid_valid_shape(1,1), j+1),
00106 !        bounding box of corner control volume of (grid_valid_shape(1,1),j)
00107 !
00108 ! The last box
00109 !     chmin (grid_valid_shape(2,1),  *) contains the minimal value of
00110 !                                   the bounding box of virtual cell
00111 !       (array(grid_valid_shape(2,1), j), array(grid_valid_shape(2,1), j+1),
00112 !        bounding box of corner control volume of (grid_valid_shape(2,1),j)
00113 !
00114 ! Corresponding values are contained for other virtual cells at block
00115 ! border.
00116 !
00117   Double Precision, 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   Double Precision, 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 ! Maximum of grid coordinates per grid cell
00124 
00125   Double Precision, Intent (In)   :: tol
00126 !
00127 ! Absolute tolerance for search of "identical" points TOL >= 0.0
00128 !
00129 ! !INPUT/OUTPUT PARAMETERS:
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 ! Finest level number on which a grid cell was found for point i,j,k.
00136 ! Input :
00137 !     Level number < -1: Point was not found and
00138 !                        and last level number was (-found(i,j,k))
00139 !     Level number = -(nlev+1): Never found
00140 ! Output :
00141 !     found(i,j,k) = +1: Point (i,j,k) is located on a point
00142 !                        of the method grid.
00143 !     found(i,j,k) = -1: Point (i,j,k) is located in a cell
00144 !                        of the method grid.
00145 !     Otherwise, not contained in the method grid.
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 ! Indices of the grid cell in which the point was found.
00153 ! The indices are relative to "shape".
00154 !
00155 ! !OUTPUT PARAMETERS:
00156 !
00157   Integer, Intent (Out)           :: ierror
00158 
00159 !     Returns the error code of PSMILE_mg_method_2d_dble;
00160 !             ierror = 0 : No error
00161 !             ierror > 0 : Severe error
00162 !
00163 ! !DEFINED PARAMETERS:
00164 !
00165 !  val_direct   = Code for locations which should to be directly transferred
00166 !                 to the destination process.
00167 !  val_coupler  = Code for locations which should to be transferred
00168 !                 to the coupler process.
00169 !  val_notfound = Code for locations which was not found
00170 !
00171   Integer, Parameter              :: nrefd = 5 * 5
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 !     ... for levels
00181 !
00182   Integer, Parameter              :: lev = 1
00183 !
00184 ! !LOCAL VARIABLES
00185 !
00186 !     ... for locations searched
00187 !
00188   Integer                         :: i, j, k, l
00189   Integer                         :: ibegl, jpart, n
00190   Integer                         :: curr_loc_list (temp_len, ndim_2d)
00191 !
00192 !     ... For locations controlled
00193 !
00194 
00195   Integer, parameter              :: ijkctl (ndim_2d, nrefd) = reshape (              
00196                                           (/-1,-1,   0,-1,   1,-1,                    
00197                                             -1, 0,   0, 0,   1, 0,                    
00198                                             -1, 1,   0, 1,   1, 1,                    
00199                                             -2,-2,  -2,-1,  -2, 0,  -2, 1,  -2, 2,    
00200                                             -1,-2,  -1, 2,                            
00201                                              0,-2,   0, 2,                            
00202                                              1,-2,   1, 2,                            
00203                                              2,-2,   2,-1,   2, 0,   2, 1,   2, 2 /), 
00204                                              (/ndim_2d, nrefd/))
00205   Integer                         :: ijkdst (temp_len, ndim_2d, nrefd)
00206   Integer                         :: irange (temp_len, 2, ndim_2d)
00207 
00208   Integer                         :: ii, jj, nref (temp_len)
00209 !
00210 !     ... For cyclic cells
00211 !
00212   Logical                         :: cyclic_req
00213 
00214   Double Precision                :: period2 (ndim_2d)
00215 !
00216 !     ... Local reference point and edges
00217 !
00218   Type (point_dble)               :: xy_source(temp_len, 4)
00219   Type (point_dble)               :: xy_source_point(temp_len)
00220   Type (point_dble), target       :: xy_target(temp_len)
00221   Type (point_dble), target       :: xy_target_rotated(temp_len)
00222   Type p_xy_target
00223      Type (point_dble), pointer   :: p
00224   end Type p_xy_target
00225   Type(p_xy_target)               :: xy_target_pointer(temp_len)
00226 !
00227 !     ...
00228 !
00229   Integer                         :: num_curr_locs, nc_full, np, nc_part
00230   Integer                         :: list (temp_len*2), full_list (temp_len)
00231   Integer                         :: part_list (temp_len)
00232   Logical                         :: fnd (temp_len)
00233 
00234 #ifdef DEBUG
00235 !     ... for locations searched
00236 !
00237   Integer                         :: nfnd (0:3)
00238 #endif /* DEBUG */
00239 !
00240 ! !DESCRIPTION:
00241 !
00242 ! Subroutine "PSMILe_mg_method_2d_dble" searches the corresponding points
00243 ! of the method (grid) for the subgrid coords by the sending process.
00244 !
00245 !
00246 ! !REVISION HISTORY:
00247 !
00248 !   Date      Programmer   Description
00249 ! ----------  ----------   -----------
00250 ! 03.07.21    H. Ritzdorf  created
00251 ! 10.09.10    M. Hanke     major revision to address problem at poles, date line, and 0-Meridian
00252 !             R. Redler
00253 !
00254 !EOP
00255 !----------------------------------------------------------------------
00256 !
00257 !  $Id: psmile_mg_method_2d_dble.F90 3248 2011-06-23 13:03:19Z coquart $
00258 !  $Author: coquart $
00259 !
00260   Character(len=len_cvs_string), save :: mycvs = 
00261        '$Id: psmile_mg_method_2d_dble.F90 3248 2011-06-23 13:03:19Z coquart $'
00262 !
00263 !----------------------------------------------------------------------
00264 !
00265 !  Initialization
00266 !
00267   ierror = 0
00268   val_notfound = -(nlev+1)
00269 
00270 #ifdef VERBOSE
00271   print 9990, trim(ch_id), lev, control
00272 
00273   call psmile_flushstd
00274 #endif /* VERBOSE */
00275 
00276 #ifdef DEBUG_TRACE
00277   if ( range(1,1) <= ictl_ind(1) .and. range(2,1) >= ictl_ind(1) .and. &
00278        range(1,2) <= ictl_ind(2) .and. range(2,2) >= ictl_ind(2)) then
00279      print 8990, trim(ch_id), ictl_ind (1:ndim_2d), &
00280           found ( ictl_ind(1),ictl_ind(2),control(1,3)), &
00281           loc (:, ictl_ind(1),ictl_ind(2),control(1,3))
00282   endif
00283 #endif /* DEBUG_TRACE */
00284   !
00285 #ifdef PRISM_ASSERTION
00286   if (tol < 0.0) then
00287      call psmile_assert (__FILE__, __LINE__, &
00288           "tol must be >= 0.0")
00289   endif
00290 
00291   if ( grid_valid_shape(1,1)-1 < coords_shape (1,1) .or. &
00292        grid_valid_shape(2,1)+1 > coords_shape (2,1) .or. &
00293        grid_valid_shape(1,2)-1 < coords_shape (1,2) .or. &
00294        grid_valid_shape(2,2)+1 > coords_shape (2,2)) then
00295 
00296      print *, 'grid_valid_shape:', grid_valid_shape
00297      print *, 'coords_shape    :', coords_shape
00298      call psmile_assert (__FILE__, __LINE__, &
00299           "insufficient coords_shape")
00300   endif
00301 #endif
00302   !
00303   ! create temporary memory (control (2,1)-control(1,1)+1)*ndtry
00304   ! and do n =1,nref before i-loop ?!?
00305   !
00306   cyclic_req = cyclic (1) .or. cyclic (2)
00307   period2 (:) = period (:) * 0.5
00308   !
00309 #ifdef DEBUG
00310   nfnd (:) = 0
00311 #endif /* DEBUG */
00312   !
00313   !===>
00314   !
00315 
00316   do k = control(1,3), control(2,3)
00317      do j = control(1,2), control (2,2)
00318 
00319 !
00320 !-----------------------------------------------------------------------
00321 !             Look in the fine parts of the coarse grid cell
00322 !             loc(:, i,j,k)
00323 !
00324 !            TODO: List over 1-dimensional Vector (i,j,k)
00325 !-----------------------------------------------------------------------
00326 !
00327         ibegl = control (1, 1)
00328 !
00329 !===> ... Look for the next cells to be controlled
00330 !         Store indices of cells in list "list"
00331 !
00332 loop_i: do while (ibegl <= control(2,1))
00333 
00334            num_curr_locs = 0
00335            i = ibegl
00336 
00337            do while (num_curr_locs < temp_len .and. i <= control(2,1))
00338               ibegl = i
00339 !cdir vector
00340               do i = ibegl, min (ibegl+(temp_len*2-1)-num_curr_locs, control(2,1))
00341                  if (found(i, j, k) == lev) then
00342                     num_curr_locs = num_curr_locs + 1
00343                     list (num_curr_locs) = i
00344                  endif
00345               end do
00346            end do
00347 
00348            if (num_curr_locs == 0) exit loop_i ! exit while loop (ibegl < ...)
00349 
00350            if (num_curr_locs < temp_len) then
00351               ibegl = control(2,1) + 1
00352            else  if (num_curr_locs > temp_len) then
00353               num_curr_locs = temp_len
00354               ibegl = list (temp_len+1)
00355            else
00356               ibegl = list (temp_len) + 1
00357            endif
00358 
00359            if (cyclic_req) then
00360               nc_full = num_curr_locs
00361               full_list (1:num_curr_locs) = list (1:num_curr_locs)
00362            else
00363               nc_full = 0
00364            endif
00365 
00366 !
00367 !===> ... Look in the real cell (curr_loc_list(i, 1),curr_loc_list(i, 2))
00368 !
00369 ! get the current locations
00370            curr_loc_list (1:num_curr_locs, 1) = loc (1, list(1:num_curr_locs), j, k)
00371            curr_loc_list (1:num_curr_locs, 2) = loc (2, list(1:num_curr_locs), j, k)
00372 
00373 #ifdef PRISM_ASSERTION
00374 !cdir vector
00375            do i = 1, num_curr_locs
00376               if ( curr_loc_list(i,1) < coords_shape(1,1) .or. &
00377                    curr_loc_list(i,1) > coords_shape(2,1) .or. &
00378                    curr_loc_list(i,2) < coords_shape(1,2) .or. &
00379                    curr_loc_list(i,2) > coords_shape(2,2) ) exit
00380            end do
00381 
00382            if (i <= num_curr_locs) then
00383               print *, 'i, curr_loc_list', i, curr_loc_list(i, :)
00384               print *, 'coords_shape', coords_shape
00385               call psmile_assert (__FILE__, __LINE__, &
00386                    'wrong index')
00387            endif
00388 #endif
00389 !cdir vector
00390            ! for each location in the current list
00391            do i = 1, num_curr_locs
00392               ! get the range which contains the neighbours of the current location.
00393               ! the range is trunctated at the grid boundaries
00394               irange (i, 1, :) = max (curr_loc_list(i,:)-2, grid_valid_shape (1, :))
00395               irange (i, 2, :) = min (curr_loc_list(i,:)+2, grid_valid_shape (2, :))
00396            end do
00397 
00398            ! compute the number of available neighbours + the actual location
00399            nref (1:num_curr_locs) = &
00400                (irange (1:num_curr_locs, 2,1) - irange (1:num_curr_locs, 1,1) + 1) * &
00401                (irange (1:num_curr_locs, 2,2) - irange (1:num_curr_locs, 1,2) + 1)
00402 
00403 #ifdef PRISM_ASSERTION
00404            if (minval (nref(1:num_curr_locs)) <= 0) then
00405               call psmile_assert (__FILE__, __LINE__, &
00406                    'nref <= 0')
00407            endif
00408 #endif
00409 !cdir vector
00410            do i = 1, num_curr_locs
00411               ! if the all (9) neighbours were available (no trunctation at grid boundary)
00412               if (nref(i) == nrefd) then
00413                  do n = 1, nrefd
00414                     ! generate list of all neighbours + current location
00415                     ijkdst (i, :, n) = curr_loc_list (i, :) + ijkctl (:, n)
00416                  end do
00417               else
00418                  ! generate list of all neighbours + current location
00419                  n = 0
00420                  do jj = irange (i, 1, 2), irange (i, 2, 2)
00421                     do ii = irange (i, 1, 1), irange (i, 2, 1)
00422                        n = n + 1
00423                        ijkdst (i, 1, n) = ii
00424                        ijkdst (i, 2, n) = jj
00425                     end do ! ii
00426                  end do ! jj
00427               endif
00428 
00429 #ifdef DEBUG_TRACE
00430               if (j == ictl_ind(2)) then
00431                  if (list(i) == ictl_ind(1)) then
00432                     print *, 'ictl: indices going to be checked:'
00433 
00434                     do n = 1, nref(n)
00435                         print *, n, "ijkdst (i, :, n)", ijkdst (i, :, n)
00436                     end do
00437                  endif
00438               end if
00439 #endif /* DEBUG_TRACE */
00440            end do
00441            ! get the target coords for all locations currently processed
00442 !cdir vector
00443            do i = 1, num_curr_locs
00444               xy_target(i)%x = coords1(list(i),j,k)
00445               xy_target(i)%y = coords2(list(i),j,k)
00446            enddo
00447 
00448            nc_part = 0
00449 
00450 !
00451 !===> ... Check for matching points
00452 !
00453 !           do n = 1, nrefd
00454            do n = 1, 9 ! this is a little optimisation, because only the first
00455                        ! 9 can be matching points (see definition of ijkctl)
00456 !cdir vector
00457               do i = 1, num_curr_locs
00458                  if (n <= nref(i)) then
00459                     xy_source_point%x = x_coords (ijkdst(i,1,n), ijkdst (i,2,n))
00460                     xy_source_point%y = y_coords (ijkdst(i,1,n), ijkdst (i,2,n))
00461                     ! if the source cell is close to the pole
00462                     ! transform source and target point to the equator
00463                     if ( y_coords (ijkdst(i,1,n), ijkdst (i,2,n)) >  pole_threshold .or. &
00464                          y_coords (ijkdst(i,1,n), ijkdst (i,2,n)) < -pole_threshold ) then
00465 
00466                        call psmile_transrot_dble (                     &
00467                              x_coords (ijkdst(i,1,n), ijkdst (i,2,n)), &
00468                              y_coords (ijkdst(i,1,n), ijkdst (i,2,n)), &
00469                              xy_source_point(i)%x, xy_source_point(i)%y )
00470 
00471                        xy_target_rotated(i) = xy_target(i)
00472 
00473                        call psmile_transrot_dble (                        &
00474                              coords1(list(i),j,k), coords2(list(i),j,k),   &
00475                              xy_target_rotated(i)%x, xy_target_rotated(i)%y )
00476 
00477                        xy_target_pointer(i)%p => xy_target_rotated(i)
00478                     else
00479                        xy_source_point(i)%x = x_coords (ijkdst(i,1,n), ijkdst (i,2,n))
00480                        xy_source_point(i)%y = y_coords (ijkdst(i,1,n), ijkdst (i,2,n))
00481                        xy_target_pointer(i)%p => xy_target(i)
00482                     endif
00483 
00484                     ! if source and target point match
00485                     if (abs(xy_target_pointer(i)%p%x - xy_source_point(i)%x) < tol .and. &
00486                         abs(xy_target_pointer(i)%p%y - xy_source_point(i)%y) < tol) then
00487 !not needed            found (list(i),j,k)  = val_direct ! == lev == 1
00488                        loc (:, list(i),j,k) = ijkdst (i, :, n)
00489                        nc_part = nc_part + 1
00490                        part_list (nc_part) = i
00491                     endif
00492                  end if
00493               end do
00494            end do
00495 
00496            ! if matching points were found
00497            if (nc_part > 0) then
00498 
00499 #ifdef DEBUG
00500               nfnd (1) = nfnd (1) + nc_part
00501 #endif /* DEBUG */
00502 
00503               ! if all target cells in the current list have a matching source cell
00504               ! get the next list of target points
00505               if (num_curr_locs == nc_part) cycle loop_i
00506 !
00507 !===> ...... Remove indices of "part_list" from "list"
00508 !            ??? Vectorisierung mittels hilfsvektors ?
00509 !
00510               do np = 0, nc_part-1
00511                  i = part_list (nc_part - np)
00512 
00513                  if (i /= num_curr_locs-np) then
00514                     list (i) = list (num_curr_locs-np)
00515                     nref (i) = nref (num_curr_locs-np)
00516                     xy_target(i) = xy_target(num_curr_locs-np)
00517 
00518                     curr_loc_list (i, :) = curr_loc_list (num_curr_locs-np, :)
00519 
00520                     ijkdst (i, :, :) = ijkdst (num_curr_locs-np, :, :)
00521                     irange (i, :, :) = irange (num_curr_locs-np, :, :)
00522                  endif
00523               end do
00524 
00525               num_curr_locs = num_curr_locs - nc_part
00526            endif
00527 !
00528 !-------------------------------------------------------------------------------
00529 !     ... Points are not located on the grid point of the method.
00530 !         Look for the correct method cell.
00531 !-------------------------------------------------------------------------------
00532 !
00533 !cdir vector
00534            do i = 1, num_curr_locs
00535               found (list(i), j, k) = val_coupler
00536            end do
00537 !
00538 !===> ....... Get cell range "irange" for search of points
00539 !             which are located at the lower boundary
00540 !             (irange (i,1,:) must be extended) (see virtual cells)
00541 !
00542            nc_part = 0
00543 !cdir vector
00544            do i = 1, num_curr_locs
00545               if ( curr_loc_list(i,1) == grid_valid_shape(1,1) .or. &
00546                    curr_loc_list(i,2) == grid_valid_shape(1,2)) then
00547                  nc_part = nc_part + 1
00548                  part_list (nc_part) = i
00549               endif
00550            end do
00551            !
00552 !cdir vector
00553            do np = 1, nc_part
00554               i = part_list (np)
00555               !
00556               ! This will place some locations out of the grid_valid_shape. This must therefore be
00557               !  considered during the ASSERTION check in psmile_neigh_trili_irreg2
00558               !                                           psmile_neigh_near_irreg2
00559               !                                           psmile_neigh_tricu_irreg2, etc.
00560               !
00561               irange (i,1, :) = max (curr_loc_list(i,:)-2, grid_valid_shape(1,:)-1)
00562               irange (i,2, :) = min (curr_loc_list(i,:)+2, grid_valid_shape(2,:))
00563 
00564               nref(i) = (irange (i,2,1) - irange (i,1,1) + 1) * &
00565                         (irange (i,2,2) - irange (i,1,2) + 1)
00566 
00567               if (nref(i) == nrefd) then
00568                  do n = 1, nrefd
00569                     ijkdst (i, :, n) = curr_loc_list (i, :) + ijkctl (:, n)
00570                  end do
00571               else
00572                  n = 0
00573                  do jj = irange (i,1, 2), irange (i,2, 2)
00574                     do ii = irange (i,1, 1), irange (i,2, 1)
00575                        n = n + 1
00576                        ijkdst (i, 1, n) = ii
00577                        ijkdst (i, 2, n) = jj
00578                     end do ! ii
00579                  end do ! jj
00580               endif
00581 #ifdef DEBUG_TRACE
00582               if (j == ictl_ind(2)) then
00583                  if (list(i) == ictl_ind(1)) then
00584                     print *, 'ictl: new indices going to be checked:'
00585 
00586                     do n = 1, nref(n)
00587                         print *, n, "ijkdst (i, :, n)", ijkdst (i, :, n)
00588                     end do
00589                  endif
00590               end if
00591 #endif /* DEBUG_TRACE */
00592            end do
00593 
00594 #ifdef PRISM_ASSERTION
00595 !cdir vector
00596            do np = 1, nc_part
00597               if (nref(part_list (np)) <= 0) exit
00598            end do
00599 
00600            if (np <= nc_part) then
00601               call psmile_assert (__FILE__, __LINE__, &
00602                    'nref <= 0')
00603            endif
00604 #endif
00605 !
00606 !===> ... Find correct source method cell
00607 !
00608            fnd (1:num_curr_locs) = .false.
00609 
00610            do n = 1, nrefd
00611 !cdir vector
00612               do i = 1, num_curr_locs
00613 
00614                  ! This is the dead of vectorisation!
00615                  ! split loop over i into several smaller loops.
00616                  ! MoHa: once this runs we will have a look at
00617                  ! some performance measurments and determine
00618                  ! whether this is required or not
00619 
00620                  if ( n <= nref (i) .and. .not. fnd(i) ) then
00621 
00622                     ! move target longitudes into the common local longitude range
00623                     do while ( xy_target(i)%x < chmax1(ijkdst(i, 1,n),ijkdst(i, 2,n)) - period(1) )
00624                        xy_target(i)%x =  xy_target(i)%x + period(1)
00625                     enddo
00626 
00627                     do while ( xy_target(i)%x > chmin1(ijkdst(i, 1,n),ijkdst(i, 2,n)) + period(1) )
00628                        xy_target(i)%x =  xy_target(i)%x - period(1)
00629                     enddo
00630 
00631                     ! if the target coord is within the bounding box of the current cell
00632                     ! -> get the respective cell
00633                     if ( xy_target(i)%x >= chmin1(ijkdst(i,1,n), ijkdst(i,2,n)) - tol .and. &
00634                          xy_target(i)%y >= chmin2(ijkdst(i,1,n), ijkdst(i,2,n)) - tol .and. &
00635                          xy_target(i)%x <= chmax1(ijkdst(i,1,n), ijkdst(i,2,n)) + tol .and. &
00636                          xy_target(i)%y <= chmax2(ijkdst(i,1,n), ijkdst(i,2,n)) + tol) then
00637 
00638 #ifdef DEBUG_TRACE
00639                        if (j == ictl_ind(2)) then
00640                           if (list(i) == ictl_ind(1)) then
00641                              print *, 'ictl: is within bounding box:'
00642                              PRINT*, 'n :',n
00643                              print *, 'ch1', chmin1(ijkdst(i,1,n), ijkdst(i,2,n)), &
00644                                                chmax1(ijkdst(i,1,n), ijkdst(i,2,n))
00645                              print *, 'ch2', chmin2(ijkdst(i,1,n), ijkdst(i,2,n)), &
00646                                                chmax2(ijkdst(i,1,n), ijkdst(i,2,n))
00647                           endif
00648                        end if
00649 #endif /* DEBUG_TRACE */
00650 
00651                        ! if the bounding box is close to the pole
00652                        if ( chmin2(ijkdst(i,1,n),ijkdst(i,2,n)) >  pole_threshold .or. &
00653                             chmax2(ijkdst(i,1,n),ijkdst(i,2,n)) < -pole_threshold ) then
00654 
00655                           ! translate source and target grid coordintates to the equator
00656                           ! to do: make pole_threshold dependent on hor. resolution
00657                           call psmile_transrot_dble (                       &
00658                                x_coords (ijkdst(i,1,n)  , ijkdst(i,2,n)  ), &
00659                                y_coords (ijkdst(i,1,n)  , ijkdst(i,2,n)  ), &
00660                                xy_source(i, 1)%x, xy_source(i, 1)%y )
00661 
00662                           call psmile_transrot_dble (                       &
00663                                x_coords (ijkdst(i,1,n)+1, ijkdst(i,2,n)  ), &
00664                                y_coords (ijkdst(i,1,n)+1, ijkdst(i,2,n)  ), &
00665                                xy_source(i, 2)%x, xy_source(i, 2)%y )
00666 
00667                           call psmile_transrot_dble (                       &
00668                                x_coords (ijkdst(i,1,n)+1, ijkdst(i,2,n)+1), &
00669                                y_coords (ijkdst(i,1,n)+1, ijkdst(i,2,n)+1), &
00670                                xy_source(i, 3)%x, xy_source(i, 3)%y )
00671 
00672                           call psmile_transrot_dble (                       &
00673                                x_coords (ijkdst(i,1,n)  , ijkdst(i,2,n)+1), &
00674                                y_coords (ijkdst(i,1,n)  , ijkdst(i,2,n)+1), &
00675                                xy_source(i, 4)%x, xy_source(i, 4)%y )
00676 
00677                           call psmile_transrot_dble (                        &
00678                                coords1(list(i),j,k), coords2(list(i),j,k),   &
00679                                xy_target_rotated(i)%x, xy_target_rotated(i)%y )
00680 
00681                           ! the transformation might cause some bad coordinates
00682                           call psmile_transform_cell_cyclic ( xy_source(i, :)%x, period(1), ierror )
00683 
00684                           ! move target longitudes into the common local longitude range
00685                           do while ( xy_target_rotated(i)%x < minval(xy_source(i, :)%x) - period2(1) )
00686                              xy_target_rotated(i)%x =  xy_target_rotated(i)%x + period(1)
00687                           enddo
00688 
00689                           do while ( xy_target_rotated(i)%x > maxval(xy_source(i, :)%x) + period2(1) )
00690                              xy_target_rotated(i)%x =  xy_target_rotated(i)%x - period(1)
00691                           enddo
00692 
00693                           xy_target_pointer(i)%p => xy_target_rotated(i)
00694 
00695 #ifdef DEBUG_TRACE
00696                           if (j == ictl_ind(2)) then
00697                              if (list(i) == ictl_ind(1)) then
00698 
00699                                 print *, 'ictl: point and cell are rotated to equator'
00700                                 PRINT*, 'n :',n
00701                                 print "('  point coordinates before: ', 2f8.2)", xy_target(i)
00702                                 print "('  cell  coordinates before: ', 8f8.2)", &
00703                                    x_coords (ijkdst(i,1,n)  , ijkdst(i,2,n)  ),  &
00704                                    y_coords (ijkdst(i,1,n)  , ijkdst(i,2,n)  ),  &
00705                                    x_coords (ijkdst(i,1,n)+1, ijkdst(i,2,n)  ),  &
00706                                    y_coords (ijkdst(i,1,n)+1, ijkdst(i,2,n)  ),  &
00707                                    x_coords (ijkdst(i,1,n)+1, ijkdst(i,2,n)+1),  &
00708                                    y_coords (ijkdst(i,1,n)+1, ijkdst(i,2,n)+1),  &
00709                                    x_coords (ijkdst(i,1,n)  , ijkdst(i,2,n)+1),  &
00710                                    y_coords (ijkdst(i,1,n)  , ijkdst(i,2,n)+1)
00711                              endif
00712                           end if
00713 #endif /* DEBUG_TRACE */
00714                        else
00715 
00716                           xy_source(i, 1)%x = x_coords (ijkdst(i,1,n)  , ijkdst(i,2,n)  )
00717                           xy_source(i, 1)%y = y_coords (ijkdst(i,1,n)  , ijkdst(i,2,n)  )
00718 
00719                           xy_source(i, 2)%x = x_coords (ijkdst(i,1,n)+1, ijkdst(i,2,n)  )
00720                           xy_source(i, 2)%y = y_coords (ijkdst(i,1,n)+1, ijkdst(i,2,n)  )
00721 
00722                           xy_source(i, 3)%x = x_coords (ijkdst(i,1,n)+1, ijkdst(i,2,n)+1)
00723                           xy_source(i, 3)%y = y_coords (ijkdst(i,1,n)+1, ijkdst(i,2,n)+1)
00724 
00725                           xy_source(i, 4)%x = x_coords (ijkdst(i,1,n)  , ijkdst(i,2,n)+1)
00726                           xy_source(i, 4)%y = y_coords (ijkdst(i,1,n)  , ijkdst(i,2,n)+1)
00727 
00728                           xy_target_pointer(i)%p => xy_target(i)
00729 
00730                           ! if any corner of the current cell is outside its own bounding box, then
00731                           ! the cell requires transforming
00732                           if (any(xy_source(i, :)%x < chmin1(ijkdst(i,1,n),ijkdst(i,2,n))) .or. &
00733                               any(xy_source(i, :)%x > chmax1(ijkdst(i,1,n),ijkdst(i,2,n)))) then
00734                              call psmile_transform_cell_cyclic ( xy_source(i, :)%x, period(1), ierror )
00735                           endif
00736 
00737                        endif
00738 
00739                        ! check whether the point is within the cell
00740                        if (point_is_in_quadrangle(xy_target_pointer(i)%p, xy_source(i,:))) then
00741                           loc(1,list(i),j,k) = ijkdst(i,1,n)
00742                           loc(2,list(i),j,k) = ijkdst(i,2,n)
00743                           fnd (i) = .true.
00744                        endif
00745 #ifdef DEBUG_TRACE
00746                        if (j == ictl_ind(2)) then
00747                           if (list(i) == ictl_ind(1)) then
00748 
00749                              print "('  point coordinates: ', 2f8.2)", xy_target_pointer(i)%p
00750                              print "('  cell  coordinates: ', 8f8.2)", xy_source(i,:)
00751 
00752                              if (fnd (i)) then
00753                                 print *, 'point is within this cell'
00754                              else
00755                                 print *, 'point is not within this cell'
00756                              endif
00757                           endif
00758                        end if
00759                     else
00760                        if (j == ictl_ind(2)) then
00761                           if (list(i) == ictl_ind(1)) then
00762                              print *, 'ictl: is not within bounding box:'
00763                              PRINT*, 'n :',n
00764                              print *, 'ch1', chmin1(ijkdst(i,1,n), ijkdst(i,2,n)), &
00765                                                chmax1(ijkdst(i,1,n), ijkdst(i,2,n))
00766                              print *, 'ch2', chmin2(ijkdst(i,1,n), ijkdst(i,2,n)), &
00767                                                chmax2(ijkdst(i,1,n), ijkdst(i,2,n))
00768                           endif
00769                        end if
00770 #endif /* DEBUG_TRACE */
00771                     endif
00772                  endif
00773               end do ! i
00774            end do ! n
00775            
00776 ! MoHa: This loop marks target cells as not found that where found on the source corner
00777 ! grid but not on the method grid. In order for the nearest neighbour search to be
00778 ! activated for these cells the found array needs to be set to val_coupler for these
00779 ! points. That is why this loop is deactivated for now. (see ticket #109)
00780 #if 0
00781            ! check if there are point which did not find a matching cell
00782            ! mark them as not found
00783            do n = 1, nrefd
00784 !cdir vector
00785               do i = 1, num_curr_locs
00786                  if ( n <= nref (i) .and. .not. fnd(i) ) then
00787 #ifdef DEBUG_TRACE
00788                      if (list(i) == ictl_ind(1) .and. j == ictl_ind(2) .and. &
00789                          found(list(i), j, k) /= val_notfound) then
00790 
00791                         print *, "ictl: The point was in a corner cell but not in any method cell."
00792                      endif
00793 #endif
00794                     found(list(i), j, k) = val_notfound
00795                  endif
00796               enddo
00797            enddo
00798 #endif
00799 
00800 #ifdef DEBUG_TRACE
00801            if (j == ictl_ind(2)) then
00802               do i = 1, num_curr_locs
00803                  if (list(i) == ictl_ind(1)) exit
00804               end do
00805               !
00806               if (i <= num_curr_locs) then
00807                  print *, 'mg_method: ictl_ind', ictl_ind(1:ndim_2d), &
00808                       ': ', fnd(i), coords1(list(i),j,k), coords2(list(i),j,k)
00809                  print *, ' .... found, i, j: ', fnd(i)
00810                  if ( fnd(i) ) then
00811                    print *, ' .... found, i, j: ', loc(1:2,list(i),j,k)
00812                  else
00813                    print *, ' did not find new i, j!'
00814                  endif
00815               endif
00816            end if
00817 #endif /* DEBUG_TRACE */
00818 !
00819 !-----------------------------------------------------------------------------
00820 !     ... All points supposed to be sent to the coupler are examined
00821 !     ... Is there an periodic index
00822 !-----------------------------------------------------------------------------
00823 !
00824 
00825            do l = 1, ndim_2d
00826               if (cyclic(l)) then
00827                  !
00828                  ! This is an cyclic/periodic index;
00829                  ! Move points, which are located in virtual cell 0
00830                  ! into cell "grid_valid_shape(2,l)"
00831                  !
00832                  do jpart = 1, num_curr_locs, 5000
00833 !cdir vector loopcnt=5000
00834                     do i = jpart, min (num_curr_locs, jpart+5000-1)
00835                        if (loc (l, list(i),j,k) < grid_valid_shape (1,l) &
00836                            .and. fnd(i)) then
00837                            loc (l, list(i),j,k) = grid_valid_shape (2,l)
00838                        endif
00839                     enddo
00840 
00841                  enddo ! jpart
00842               endif
00843            end do ! l
00844 !
00845 !===> ... Get next list of points for first logical index
00846 !
00847         end do loop_i ! while
00848      end do ! j
00849   end do ! k
00850 
00851 !
00852 !===> All done
00853 !
00854 #ifdef DEBUG_TRACE
00855   if ( range(1,1) <= ictl_ind(1) .and. range(2,1) >= ictl_ind(1) .and. &
00856        range(1,2) <= ictl_ind(2) .and. range(2,2) >= ictl_ind(2)) then
00857 
00858      print *, 'coords', coords1(ictl_ind(1),ictl_ind(2),control(1,3)), &
00859           coords2(ictl_ind(1),ictl_ind(2),control(1,3))
00860 
00861      print 8990,  trim(ch_id), ictl_ind(1:ndim_2d), &
00862           found ( ictl_ind(1),ictl_ind(2),control(1,3)), &
00863           loc (:, ictl_ind(1),ictl_ind(2),control(1,3))
00864   endif
00865 
00866 8990 format (1x, a, ': psmile_mg_method_2d_dble: ictl_ind', 2i5, &
00867           '; found', i3, ' loc ', 2i6)
00868 #endif /* DEBUG_TRACE */
00869 
00870 #ifdef DEBUG
00871   print 9970, trim(ch_id), lev, nfnd
00872 9970 format (1x, a, ': psmile_mg_method_2d_dble: lev =', i3, &
00873        ': not :', i5, ', dir: ', i6, ', fnd :', i6, i5)
00874 
00875 #ifdef DEBUG_PRINT
00876   call psmile_print_found_dble (comp_info, lev,           &
00877        found, loc, range,                          &
00878        coords1, coords2, coords1,                  &
00879        search_shape, search_shape, search_shape,   &
00880        control, ndim_2d, "psmile_mg_method_2d_dble: eof")
00881 
00882 #endif /* DEBUG_PRINT */
00883 #endif /* DEBUG */
00884 
00885 #ifdef VERBOSE
00886   print 9980, trim(ch_id), lev
00887 
00888   call psmile_flushstd
00889 #endif /* VERBOSE */
00890 !
00891 !  Formats:
00892 !
00893 
00894 9990 format (1x, a, ': psmile_mg_method_2d_dble: level', i3, &
00895              ', control', 6i6)
00896 9980 format (1x, a, ': psmile_mg_method_2d_dble: eof, level', i3)
00897 
00898 end subroutine psmile_mg_method_2d_dble

Generated on 1 Dec 2011 for Oasis4 by  doxygen 1.6.1