psmile_neigh_nearest_3d_real.F90

Go to the documentation of this file.
00001 !-----------------------------------------------------------------------
00002 ! Copyright 2006-2010, NEC Europe Ltd., London, UK.
00003 ! All rights reserved. Use is subject to OASIS4 license terms.
00004 !-----------------------------------------------------------------------
00005 !BOP
00006 !
00007 ! !ROUTINE: PSMILe_Neigh_nearest_3d_real
00008 !
00009 ! !INTERFACE:
00010 
00011       subroutine psmile_neigh_nearest_3d_real (grid_id,    &
00012                       coords1, coords2, coords3,           &
00013                       x_coords, y_coords, z_coords, coords_shape, &
00014                       mask_array, mask_shape,  mask_available,    &
00015                       sin_values, cos_values,              &
00016                       grid_valid_shape,                    &
00017                       srcloc, nloc, nprev, nsearch,        &
00018                       neighbors_3d, num_neigh,             &
00019                       extra_search, ierror)
00020 !
00021 ! !USES:
00022 !
00023       use PRISM_constants
00024 !
00025       use PSMILe, dummy_interface => PSMILe_Neigh_nearest_3d_real
00026 
00027       Implicit none
00028 !
00029 ! !INPUT PARAMETERS:
00030 !
00031       Integer, Intent (In)            :: grid_id
00032 
00033 !     Link to the info on the component in which the donor
00034 !     should be searched.
00035 
00036       Integer, Intent (In)            :: nloc
00037 !
00038 !     Total  Number of locations to be transferred
00039 
00040       Integer, Intent (In)            :: nprev
00041 !
00042 !     Number of locations for which neighbours were already searched.
00043 !
00044       Integer, Intent (In)            :: nsearch
00045 !
00046 !     Number of locations for which neighbours should be searched in this
00047 !     run.
00048 !
00049       Integer, Intent (In)            :: srcloc (ndim_3d, nloc)
00050 !
00051 !     Locations to be transferred
00052 
00053 !     Indices of the grid cell in which the point was found.
00054 !     The indices are relative to shape.
00055 
00056       Real, Intent (In)               :: coords1 (nloc)
00057       Real, Intent (In)               :: coords2 (nloc)
00058       Real, Intent (In)               :: coords3 (nloc)
00059 
00060 !     Coordinates to be searched (and found)
00061 
00062       Integer, Intent (In)            :: coords_shape (2, ndim_3d)
00063 
00064 !     Dimension of coordinates (method) x_coords, ...
00065 
00066       Real, Intent (In)               :: 
00067          x_coords( coords_shape(1,1):coords_shape(2,1), 
00068                    coords_shape(1,2):coords_shape(2,2), 
00069                    coords_shape(1,3):coords_shape(2,3))
00070       Real, Intent (In)               :: 
00071          y_coords( coords_shape(1,1):coords_shape(2,1), 
00072                    coords_shape(1,2):coords_shape(2,2), 
00073                    coords_shape(1,3):coords_shape(2,3))
00074       Real, Intent (In)               :: 
00075          z_coords( coords_shape(1,1):coords_shape(2,1), 
00076                    coords_shape(1,2):coords_shape(2,2), 
00077                    coords_shape(1,3):coords_shape(2,3))
00078 
00079 !     Coordinates of the method
00080 
00081       Integer, Intent (In)            :: mask_shape (2, ndim_3d)
00082 
00083 !     Dimension of (method) mask array
00084 
00085       Logical, Intent (In)            :: 
00086          mask_array ( mask_shape (1,1):mask_shape (2,1), 
00087                       mask_shape (1,2):mask_shape (2,2), 
00088                       mask_shape (1,3):mask_shape (2,3))
00089 !     Mask of the method
00090 
00091       Logical, Intent (In)            :: mask_available
00092 
00093 !     Is mask specified in array "mask_array" ?
00094 !     mask_available = .false. : Mask is not available
00095 !                      .true.  : Mask is     available
00096 
00097       Integer, Intent (In)            :: grid_valid_shape (2, ndim_3d)
00098 
00099 !     Specifies the valid block shape for the "ndim_3d"-dimensional block
00100 
00101       Real, Intent (In)               :: 
00102          sin_values ( grid_valid_shape(1,1):grid_valid_shape(2,1), 
00103                       grid_valid_shape(1,2):grid_valid_shape(2,2), 
00104                       grid_valid_shape(1,3):grid_valid_shape(2,3), 2)
00105 !
00106 !     Sin values of Longitudes and Latitudes (x_coords, y_coords)
00107 !
00108       Real, Intent (In)               :: 
00109          cos_values ( grid_valid_shape(1,1):grid_valid_shape(2,1), 
00110                       grid_valid_shape(1,2):grid_valid_shape(2,2), 
00111                       grid_valid_shape(1,3):grid_valid_shape(2,3), 2)
00112 !
00113 !     Cos values of Longitudes and Latitudes (x_coords, y_coords)
00114 !
00115       Integer, Intent (In)            :: num_neigh
00116 !
00117 !     Last dimension of neighbors array "neighbors_3d" and
00118 !     number of neighbors to be searched.
00119 !
00120 ! !INPUT/OUTPUT PARAMETERS:
00121 !
00122       Type (Extra_search_info)        :: extra_search
00123 !
00124 !     Number of locations where
00125 !       (*) required mask values were not true
00126 !
00127 ! !OUTPUT PARAMETERS:
00128 !
00129       Integer, Intent (Out)           :: neighbors_3d (ndim_3d, nloc, num_neigh)
00130 
00131 !     Indices of neighbor locations (interpolation bases)
00132 !
00133       integer, Intent (Out)           :: ierror
00134 
00135 !     Returns the error code of PSMILE_Neigh_nearest_3d_real;
00136 !             ierror = 0 : No error
00137 !             ierror > 0 : Severe error
00138 !
00139 ! !DEFINED PARAMETERS:
00140 !
00141 !  NREFD = Number of locations to be controlled
00142 !        = 4 ** ndim_3d
00143 !
00144 !  lon   = Index of Longitudes in arrays "sin_values" and "cos_values"
00145 !  lat   = Index of Latitudes  in arrays "sin_values" and "cos_values"
00146 !  REAL_EARTH = Earth radius
00147 !
00148       Integer, Parameter              :: nrefd = 4 * 4 * 4
00149       Integer, Parameter              :: lon = 1
00150       Integer, Parameter              :: lat = 2
00151 !
00152       Real, Parameter                 :: real_earth = 6400.0
00153 !
00154       Real, Parameter                 :: tol = 1d-6 ! muss argument werden
00155       Real, Parameter                 :: eps = 1d-6
00156       Real, Parameter                 :: acosp1 =  1.0
00157       Real, Parameter                 :: acosm1 = -1.0
00158 !
00159 ! !LOCAL VARIABLES
00160 !
00161 !     ... for locations searched
00162 !
00163       Integer                         :: i, ibeg, iend, n
00164 !
00165 !     ... For locations controlled
00166 !
00167       Real, Pointer                   :: sin_search (:, :)
00168       Real, Pointer                   :: cos_search (:, :)
00169 !
00170       Integer                         :: ijkctl (ndim_3d, nrefd)
00171       Integer                         :: ijkdst (ndim_3d, nrefd)
00172       Integer                         :: irange (2, ndim_3d)
00173       Integer                         :: itemp (ndim_3d)
00174 !
00175       Logical                         :: mask_dist (nrefd)
00176       Real                            :: distance (nrefd), dist, fac, val
00177       Real                            :: dist_control
00178 !
00179       Integer                         :: ijk (ndim_3d), ijk0 (ndim_3d)
00180       Integer                         :: leni, lenij
00181       Integer                         :: ii, jj, kk, nref
00182 !
00183       Integer                         :: nfound
00184       Integer                         :: nmin (1)
00185 !
00186 !     ... for error handling
00187 !
00188       Real                            :: real_huge
00189 !
00190 !     ... for error handling
00191 !
00192       Integer, Parameter              :: nerrp = 2
00193       Integer                         :: ierrp (nerrp)
00194 !
00195 !
00196 ! !DESCRIPTION:
00197 !
00198 ! Subroutine "PSMILe_Neigh_nearest_3d_real" searches the next "num_neigh"
00199 ! neighbours on the method-grid (x_coords, y_coords, z_coords)
00200 ! for the subgrid coords sent by the destination process.
00201 !
00202 !
00203 ! !REVISION HISTORY:
00204 !
00205 !   Date      Programmer   Description
00206 ! ----------  ----------   -----------
00207 ! 03.07.21    H. Ritzdorf  created
00208 !
00209 !EOP
00210 !----------------------------------------------------------------------
00211 !
00212 !  $Id: psmile_neigh_nearest_3d_real.F90 2325 2010-04-21 15:00:07Z valcke $
00213 !  $Author: valcke $
00214 !
00215    Character(len=len_cvs_string), save :: mycvs = 
00216        '$Id: psmile_neigh_nearest_3d_real.F90 2325 2010-04-21 15:00:07Z valcke $'
00217 !
00218 !----------------------------------------------------------------------
00219 ! 
00220 !  Relative control indices 
00221 !      
00222 #ifdef OLD
00223       Integer, Parameter              :: nrefd = 3 * 3 * 3
00224 !
00225       data ((ijkctl (i, n), i=1,ndim_3d), n = 1,nrefd) &
00226                             /  0, 0, 0,                       &
00227                                1, 0, 0,   0, 1, 0,   0, 0, 1, &
00228                               -1, 0, 0,   0,-1, 0,   0, 0,-1, &
00229                               -1,-1, 0,   1,-1, 0,            &
00230                               -1, 1, 0,   1, 1, 0,            &
00231                               -1, 0,-1,   1, 0,-1,            &
00232                               -1, 0, 1,   1, 0, 1,            &
00233                                0,-1,-1,   0, 1,-1,            &
00234                                0,-1, 1,   0, 1, 1,            &
00235                               -1,-1,-1,   1,-1,-1,            &
00236                               -1, 1,-1,   1, 1,-1,            &
00237                               -1,-1, 1,   1,-1, 1,            &
00238                               -1, 1, 1,   1, 1, 1/
00239 #endif
00240 ! schoener waere eine liste, die weniger sortierung (s.o.) nachher verlangt
00241 ! veielleicht sollte man immer annehmen, dass ijk in der "Mitte" liegt
00242 !
00243       data ((ijkctl (i, n), i=1,ndim_3d), n = 1,nrefd) &
00244             /  0, 0, 0,   1, 0, 0,   2, 0, 0,   3, 0, 0, &
00245                0, 1, 0,   1, 1, 0,   2, 1, 0,   3, 1, 0, &
00246                0, 2, 0,   1, 2, 0,   2, 2, 0,   3, 2, 0, &
00247                0, 3, 0,   1, 3, 0,   2, 3, 0,   3, 3, 0, &
00248                0, 0, 1,   1, 0, 1,   2, 0, 1,   3, 0, 1, &
00249                0, 1, 1,   1, 1, 1,   2, 1, 1,   3, 1, 1, &
00250                0, 2, 1,   1, 2, 1,   2, 2, 1,   3, 2, 1, &
00251                0, 3, 1,   1, 3, 1,   2, 3, 1,   3, 3, 1, &
00252                0, 0, 2,   1, 0, 2,   2, 0, 2,   3, 0, 2, &
00253                0, 1, 2,   1, 1, 2,   2, 1, 2,   3, 1, 2, &
00254                0, 2, 2,   1, 2, 2,   2, 2, 2,   3, 2, 2, &
00255                0, 3, 2,   1, 3, 2,   2, 3, 2,   3, 3, 2, &
00256                0, 0, 3,   1, 0, 3,   2, 0, 3,   3, 0, 3, &
00257                0, 1, 3,   1, 1, 3,   2, 1, 3,   3, 1, 3, &
00258                0, 2, 3,   1, 2, 3,   2, 2, 3,   3, 2, 3, &
00259                0, 3, 3,   1, 3, 3,   2, 3, 3,   3, 3, 3/
00260 
00261 !----------------------------------------------------------------------
00262 !
00263 !  Initialization
00264 !
00265 #ifdef VERBOSE
00266       print 9990, trim(ch_id)
00267 
00268       call psmile_flushstd
00269 #endif /* VERBOSE */
00270 !
00271       ibeg = nprev + 1
00272       iend = nprev + nsearch
00273 !
00274 #ifdef PRISM_ASSERTION
00275 !cdir vector
00276          do i = ibeg, iend
00277 !
00278          if (srcloc(1,i) < coords_shape(1,1) .or. &
00279              srcloc(1,i) > coords_shape(2,1) .or. &
00280              srcloc(2,i) < coords_shape(1,2) .or. &
00281              srcloc(2,i) > coords_shape(2,2) .or. &
00282              srcloc(3,i) < coords_shape(1,3) .or. &
00283              srcloc(3,i) > coords_shape(2,3)) exit
00284          end do
00285 !
00286       if (i <= iend) then
00287          print *, i, srcloc(:, i), coords_shape
00288          call psmile_assert (__FILE__, __LINE__, &
00289                              'wrong index')
00290       endif
00291 #endif
00292 !
00293       real_huge = huge (distance(1))
00294 !
00295 !-----------------------------------------------------------------------
00296 !     Allocate and set temporary array which is used to transform
00297 !     the coordinates of poinst to be searched.
00298 !     (*) Longitudes and Latitudes are transformed
00299 !         from degrees into radients.
00300 !     (*) The z-values are currently not transformed
00301 !         ??? Whats about PRISM_Irrlonlat_sigmavrt, ...
00302 !
00303 !     Is it possible to save the values, in order to be reused.
00304 !-----------------------------------------------------------------------
00305 !
00306       Allocate (sin_search(ibeg:iend, lat), STAT = ierror)
00307 
00308       if ( ierror > 0 ) then
00309          ierrp (1) = ierror
00310          ierrp (2) = nsearch * lat
00311 
00312          ierror = PRISM_Error_Alloc
00313          call psmile_error ( ierror, 'sin_search', &
00314                              ierrp, 2, __FILE__, __LINE__ )
00315          return
00316       endif
00317 
00318       Allocate (cos_search(ibeg:iend, lat), STAT = ierror)
00319 
00320       if ( ierror > 0 ) then
00321          ierrp (1) = ierror
00322          ierrp (2) = nsearch * lat
00323 
00324          ierror = PRISM_Error_Alloc
00325          call psmile_error ( ierror, 'cos_search', &
00326                              ierrp, 2, __FILE__, __LINE__ )
00327          return
00328       endif
00329 !
00330       sin_search (ibeg:iend, lon) = coords1 (ibeg:iend) * real_deg2rad
00331       sin_search (ibeg:iend, lat) = coords2 (ibeg:iend) * real_deg2rad
00332 !
00333       cos_search = cos (sin_search)
00334       sin_search = sin (sin_search)
00335 !
00336 !     z_search(ibeg:iend) = coords3 (ibeg:iend)
00337 !
00338 !-----------------------------------------------------------------------
00339 !     Compute distances
00340 !
00341 !     dist = (real_earth + z_search)*ACOS ( SIN(lat_coords)*SIN(lat_search) +
00342 !                                           COS(lat_coords)*COS(lat_search) *
00343 !                                          (COS(lon_coords)*COS(lon_search) +
00344 !                                           SIN(lon_coords)*SIN(lon_search)) )
00345 !     dist = SQRT(dist**2 + (z_search-z_coords) ** 2)
00346 !
00347 !     with xxx_search = Coordinates to be searched
00348 !     and  xxx_coords = Grid (Methods) coordinates (transformed)
00349 !
00350 ! Suchen mit 
00351 !                    nlev = Grids(grid_id)%nlev
00352 !                    Grids(grid_id)%mg_infos(lev)%real_arrays%chmin, &
00353 !                    Grids(grid_id)%mg_infos(lev)%real_arrays%chmax, &
00354 !
00355 !-----------------------------------------------------------------------
00356 !
00357       ijk0 = Grids(grid_id)%ijk0
00358 !
00359          do i = ibeg, iend
00360 !
00361 !===> Get coarse grid index on level nlev-2
00362 !
00363          ijk = ijk0 + ((srcloc(:,i) - ijk0) / 4) * 4
00364 !
00365 !===> Get standard range to be controlled
00366 !
00367          irange (1, :) = max (ijk,   grid_valid_shape (1, :))
00368          irange (2, :) = min (ijk+3, grid_valid_shape (2, :))
00369 !
00370          nref = (irange (2,1) - irange (1,1) + 1) * &
00371                 (irange (2,2) - irange (1,2) + 1) * &
00372                 (irange (2,3) - irange (1,3) + 1)
00373 !
00374          if (nref == nrefd) then
00375                do n = 1, nrefd
00376                ijkdst (:, n) = ijk + ijkctl (:, n)
00377                end do
00378          else
00379 ! 
00380             leni  =  irange (2, 1) - irange (1, 1) + 1
00381             lenij = (irange (2, 2) - irange (1, 2) + 1) * leni 
00382 !
00383             do kk = 1, irange (2,3)-irange(1,3) + 1
00384             n = (kk-1)*lenij
00385                do jj = 1, irange (2,2)-irange(1,2) + 1
00386 
00387                   do ii = 1, leni
00388                   ijkdst (1, n+ii) = irange(1,1) + ii - 1
00389                   end do ! ii
00390 
00391                ijkdst (2, n+1:n+leni) = irange (1,2) + jj - 1
00392                n = n + leni
00393                end do ! jj
00394 
00395             ijkdst (3, (kk-1)*lenij+1:kk*lenij) = irange (1,3) + kk - 1
00396             end do ! kk
00397          endif
00398 !
00399 !====> ... Get Mask values if necessary
00400 !
00401          if (mask_available) then
00402 !cdir vector
00403                do n = 1, nref
00404                mask_dist (n) = mask_array (ijkdst (1,n), ijkdst (2,n), &
00405                                            ijkdst (3,n))
00406                end do
00407 !
00408             n = count (mask_dist(1:nref))
00409 !
00410             if (n == 0) then
00411 !              no value available
00412 !
00413                nref = 0
00414                go to 10
00415             endif
00416 !
00417             if (n < nref) then
00418                n = 0
00419                   do ii = 1, nref
00420                   if (mask_dist (ii)) then
00421                       n = n + 1
00422                       if (n /= ii) ijkdst (:, n) = ijkdst (:, ii)
00423                   endif
00424                   end do
00425                nref = n
00426             endif
00427          endif
00428 !
00429 !===> ... Compute distances to points
00430 !
00431          fac = real_earth + coords3 (i)
00432 !cdir vector
00433             do n = 1, nref
00434             val   = (  sin_values(ijkdst (1,n),                               &
00435                                   ijkdst (2,n),                               &
00436                                   ijkdst (3,n), lat) * sin_search(i, lat)     &
00437                     +  cos_values(ijkdst (1,n),                               &
00438                                   ijkdst (2,n),                               &
00439                                   ijkdst (3,n), lat) * cos_search(i, lat)     &
00440                     * (cos_values(ijkdst (1,n),                               &
00441                                   ijkdst (2,n),                               &
00442                                   ijkdst (3,n), lon) * cos_search(i, lon)     &
00443                       +sin_values(ijkdst (1,n),                               &
00444                                   ijkdst (2,n),                               &
00445                                   ijkdst (3,n), lon) * sin_search(i, lon)) )
00446 #ifdef PRISM_ASSERTION
00447             if (val < acosm1 - eps .or. val > acosp1 + eps) then
00448                print *, i, val
00449                call psmile_assert (__FILE__, __LINE__, &
00450                                    'invalid value for acos')
00451             endif
00452 #endif
00453             val = min (acosp1, val)
00454             val = max (acosm1, val)
00455 
00456             dist = fac * acos (val)
00457 
00458             distance (n) = dist*dist +                                 &
00459                            ( coords3 (i)-z_coords(ijkdst (1,n),       &
00460                                                    ijkdst (2,n),       &
00461                                                    ijkdst (3,n)))**2
00462 #ifdef DEBUG2
00463          if (i == 2) then
00464             print *, 'dist, diff', dist, &
00465                coords3 (i)-z_coords(ijkdst (1,n), ijkdst (2,n), ijkdst (3,n))
00466          endif
00467 #endif
00468             end do ! n
00469 !
00470          nmin = minloc (distance(1:nref))
00471 #ifdef MINLOCFIX
00472          if (nmin(1).eq.0) nmin=1
00473 #endif /* MINLOCFIX */
00474 !
00475 #ifdef DEBUG2
00476          if (i == 2) then
00477             print *, 'psmile_neigh_nearest_3d_real: ictl =', i
00478             print *, 'ijkdst and distance:'
00479                do n = 1, nref
00480                print *, ijkdst (:, n), distance (n)
00481                end do
00482          endif
00483 #endif
00484 !
00485 !===> ... Special case: Minimal distance == 0
00486 !
00487          if (distance(nmin(1)) <= tol) then
00488             neighbors_3d (:, i, 1) = ijkdst (:, nmin(1))
00489 !
00490                do n = 2, num_neigh
00491                neighbors_3d (:, i, n) = grid_valid_shape (1, :) - 1
00492                end do
00493 !
00494             cycle
00495          endif
00496 !
00497 !===> ... Sort distances (silly sort)
00498 !         Da muss noch etwas besseres her !
00499 !
00500 10       continue
00501 !
00502          nfound = min (num_neigh, nref)
00503          if (nfound < num_neigh) then
00504             distance(nfound+1:num_neigh) = real_huge
00505 
00506                do n = nfound+1, num_neigh
00507                ijkdst (:, n) = grid_valid_shape (1, :) - 1
00508                end do
00509          endif
00510 !
00511          do n = 1, nfound
00512 !
00513             nmin = minloc (distance(n:nref))
00514 #ifdef MINLOCFIX
00515             if (nmin(1).eq.0) nmin=1
00516 #endif /* MINLOCFIX */
00517 !
00518             if (nmin(1) /= 1) then
00519                jj = n + nmin (1) - 1
00520 !
00521                itemp  (:)    = ijkdst (:, n)
00522                ijkdst (:, n) = ijkdst (:, jj)
00523                ijkdst (:, jj) = itemp  (:)
00524 !
00525                dist = distance (n)
00526                distance (n) = distance (jj)
00527                distance (jj) = dist
00528             endif
00529          end do
00530 !
00531 #ifdef PRISM_ASSERTION
00532             do n = 1, nfound-1
00533             if (distance(n) > distance (n+1)) exit
00534             end do
00535 !
00536          if (n <= nfound-1) then
00537             print *, 'n =' , n
00538             print *, 'distance (n)  ', distance (n)
00539             print *, 'distance (n+1)', distance (n+1)
00540             print *, 'distance', distance (1:nfound)
00541             call psmile_assert (__FILE__, __LINE__, &
00542                                 'incorrect sort')
00543          endif
00544 #endif
00545 !
00546 !===> ... Use multigrid info for search in order regions
00547 !
00548          dist_control = distance (num_neigh)
00549 !
00550 #ifdef TODO
00551          call psmile_mg_search_nneigh_3d_real (distance, ijkdst, num_neigh, &
00552                      grid_id, ierror)
00553 #endif /* TODO */
00554 !
00555 !===> ... End of search; fill neighbors_3d
00556 !
00557             do n = 1, num_neigh
00558             neighbors_3d (:, i, n) = ijkdst (:, n)
00559             end do
00560 !
00561          end do ! i
00562 !
00563 !===> Deallocate memory
00564 !
00565       Deallocate (cos_search)
00566       Deallocate (sin_search)
00567 !
00568 !===> All done
00569 !
00570 #ifdef VERBOSE
00571       print 9980, trim(ch_id), ierror
00572 
00573       call psmile_flushstd
00574 #endif /* VERBOSE */
00575 !
00576 !  Formats:
00577 !
00578 9990 format (1x, a, ': psmile_neigh_nearest_3d_real')
00579 9980 format (1x, a, ': psmile_neigh_nearest_3d_real: eof, ierror =', i3)
00580 
00581       end subroutine PSMILe_Neigh_nearest_3d_real

Generated on 18 Mar 2011 for Oasis4 by  doxygen 1.6.1