psmile_neigh_nearx_sub_reg_dble.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_nearx_sub_reg_dble
00008 !
00009 ! !INTERFACE:
00010 
00011       subroutine psmile_neigh_nearx_sub_reg_dble (           &
00012                            grid_id,                          &
00013                            x_coords, y_coords, z_coords,     &
00014                            coords_shape,                     &
00015                            mask_array, mask_shape,  mask_available,    &
00016                            sin_values_lon, cos_values_lon,   &
00017                            sin_values_lat, cos_values_lat,   &
00018                            grid_valid_shape, search_mode,    &
00019                            neighbors_3d, num_neigh, nloc,    &
00020                            extra_search, ijk,                &
00021                            sin_search, cos_search, z_search, &
00022                            jbeg, jend, width, ierror)
00023 !
00024 ! !USES:
00025 !
00026       use PRISM_constants
00027 !
00028       use PSMILe, dummy_interface => PSMILe_Neigh_nearx_sub_reg_dble
00029 
00030       implicit none
00031 
00032 ! !DEFINED PARAMETERS:
00033 !
00034 !  NREF_3D = Number of locations in 3d to be controlled
00035 !          = 4 ** ndim_3d
00036 !  NREF_2D = Number of locations in 2d to be controlled
00037 !          = 4 ** ndim_2d
00038 !  NREF_1D = Number of locations in 1d to be controlled
00039 !          = 4 ** ndim_1d
00040 !
00041 !  lon   = Index of Longitudes in arrays "sin_values" and "cos_values"
00042 !  lat   = Index of Latitudes  in arrays "sin_values" and "cos_values"
00043 !  dble_earth = Earth radius in meter
00044 !
00045       Integer, Parameter              :: nref_3d = 4 * 4 * 4
00046       Integer, Parameter              :: nref_2d = 4 * 4
00047       Integer, Parameter              :: nref_1d = 4
00048       Integer, Parameter              :: lon = 1
00049       Integer, Parameter              :: lat = 2
00050       Integer, Parameter              :: maxlen = 4
00051 !
00052       Double Precision, Parameter     :: dble_earth = 6400000.0
00053 !
00054       Double Precision, Parameter     :: tol = 1d-6 ! muss argument werden
00055       Double Precision, Parameter     :: eps = 1d-6
00056       Double Precision, Parameter     :: acosp1 =  1.0d0
00057       Double Precision, Parameter     :: acosm1 = -1.0d0
00058 !
00059 ! !INPUT PARAMETERS:
00060 !
00061       Integer, Intent (In)            :: grid_id
00062 
00063 !     Info on the component in which the donor cells
00064 !     should be searched.
00065 
00066       Integer, Intent (In)            :: nloc
00067 !
00068 !     Total number of locations to be transferred
00069 !
00070       Integer, Intent (In)            :: coords_shape (2, ndim_3d)
00071 
00072 !     Dimension of coordinates (method) x_coords, ...
00073 
00074       Double Precision, Intent (In)   :: x_coords(coords_shape(1,1): 
00075                                                   coords_shape(2,1))
00076       Double Precision, Intent (In)   :: y_coords(coords_shape(1,2): 
00077                                                   coords_shape(2,2))
00078       Double Precision, Intent (In)   :: z_coords(coords_shape(1,3): 
00079                                                   coords_shape(2,3))
00080 
00081 !     Coordinates of the method
00082 
00083       Integer, Intent (In)            :: mask_shape (2, ndim_3d)
00084 
00085 !     Dimension of (method) mask array
00086 
00087       Logical, Intent (In)            :: mask_array (mask_shape (1,1): 
00088                                                      mask_shape (2,1), 
00089                                                      mask_shape (1,2): 
00090                                                      mask_shape (2,2), 
00091                                                      mask_shape (1,3): 
00092                                                      mask_shape (2,3))
00093 !     Mask of the method
00094 
00095       Logical, Intent (In)            :: mask_available
00096 
00097 !     Is mask specified in array "mask_array" ?
00098 !     mask_available = .false. : Mask is not available
00099 !                      .true.  : Mask is     available
00100 
00101       Integer, Intent (In)            :: grid_valid_shape (2, ndim_3d)
00102 
00103 !     Specifies the valid block shape for the "ndim_3d"-dimensional block
00104 
00105       Double Precision, Intent (In)   :: sin_values_lon (grid_valid_shape(1,1): 
00106                                                          grid_valid_shape(2,1))
00107       Double Precision, Intent (In)   :: sin_values_lat (grid_valid_shape(1,2): 
00108                                                          grid_valid_shape(2,2))
00109 !
00110 !     Sin values of Longitudes and Latitudes (x_coords, y_coords)
00111 !
00112       Double Precision, Intent (In)   :: cos_values_lon (grid_valid_shape(1,1): 
00113                                                          grid_valid_shape(2,1))
00114       Double Precision, Intent (In)   :: cos_values_lat (grid_valid_shape(1,2): 
00115                                                          grid_valid_shape(2,2))
00116 !    
00117 !     Cos values of Longitudes and Latitudes (x_coords, y_coords)
00118 !
00119       Integer, Intent (In)            :: search_mode
00120 !
00121 !     Specifies the search mode for nearest neigbours with
00122 !        search_mode = 3 : Full 3d-search
00123 !        search_mode = 2 : Search in 2d-hyperplane (1st 2 direction lon, lat)
00124 !        search_mode = 1 : Search in 1d-plane
00125 
00126       Integer, Intent (In)            :: num_neigh
00127 !
00128 !     Last dimension of neighbors array "neighbors_3d" and
00129 !     number of neighbors to be searched.
00130 !
00131       Integer,          Intent (In)   :: jbeg, jend
00132 !
00133 !     Dimensioning of arrays
00134 !
00135       Integer,          Intent (In)   :: ijk (ndim_3d, jbeg:jend)
00136 !
00137 !     Coarse grid index on level nlev-2 (for dimensions specified 
00138 !     by search mode). Coarsening factor had to be 4.
00139 
00140       Double Precision, Intent (In)   :: sin_search (jbeg:jend, lat)
00141 !
00142 !     Sin values of the points for which the nearest neighbours
00143 !     should be searched.
00144 !
00145       Double Precision, Intent (In)   :: cos_search (jbeg:jend, lat)
00146 !
00147 !     Cos values of the points for which the nearest neighbours
00148 !     should be searched.
00149 !
00150       Double Precision, Intent (In)   :: z_search (jbeg:jend)
00151 !
00152 !     Z values of the points for which the nearest neighbours
00153 !     should be searched.
00154 !
00155       Integer, Intent (In)            :: width (ndim_3d)
00156 !
00157       Type (Extra_search_info), Intent(In) :: extra_search
00158 !
00159 !     Number of locations where
00160 !       (*) required mask values were not true
00161 !
00162 ! !OUTPUT PARAMETERS:
00163 !
00164       Integer, Intent (Out)           :: neighbors_3d (ndim_3d, nloc, num_neigh)
00165 
00166 !     Indices of neighbor locations (interpolation bases)
00167 !
00168       Integer, Intent (Out)           :: ierror
00169 
00170 !     Returns the error code of PSMILE_Neigh_nearx_sub_reg_dble;
00171 !             ierror = 0 : No error
00172 !             ierror > 0 : Severe error
00173 !
00174 ! !LOCAL VARIABLES
00175 !
00176 !     ... For locations controlled
00177 !
00178       Integer                         :: nrefv (1:ndim_3d)
00179       Integer                         :: dim1 (2)
00180       Logical                         :: global_search
00181       Double Precision, Pointer       :: dist_dble (:, :)
00182 !
00183 !     ... (1d) Indices in (entire) data to be extra searched
00184 !
00185       Integer                         :: i, k
00186       Integer                         :: j, jpart
00187       Integer                         :: leni, lenj, lenk, lenij
00188 !
00189       Integer, Pointer                :: indices (:)
00190 !
00191 !     ... 3d Indices in local data to be extra searched
00192 !
00193       Integer                         :: ii, jj, kk
00194       Integer                         :: n, nref, nrefd
00195 !
00196       Integer                         :: ijkdst (ndim_3d, nref_3d)
00197       Integer                         :: ijkctl (ndim_3d, nref_3d)
00198       Integer                         :: irange (2, ndim_3d)
00199 !
00200 !     ... auxiliary arrays
00201 !
00202       Integer                         :: itemp (ndim_3d)
00203 !
00204 !     ... for mask values
00205 !
00206       Logical                         :: mask_dist (nref_3d)
00207 !
00208 !     ... for nearest neighbour control
00209 !
00210       Integer                         :: nmin (1)
00211       Integer                         :: nfound
00212 !
00213       Double Precision                :: dble_huge
00214       Double Precision                :: distance (nref_3d), dist, fac, val
00215 
00216       Integer                         :: nlev
00217       Type (Extra_search_nn)          :: nn_srch
00218       Type (Grid), Pointer            :: grid_info
00219       Type (Extra_search_dble)        :: arrays
00220       Type (Enddef_mg_double), Pointer :: my_arrays
00221 !
00222 !     ... for error handling
00223 !
00224       Integer, Parameter              :: nerrp = 2
00225       Integer                         :: ierrp (nerrp)
00226 !
00227 !
00228 ! !DESCRIPTION:
00229 !
00230 ! Subroutine "PSMILe_Neigh_nearx_sub_reg_dble" searches the next "num_neigh"
00231 ! nearest neighbours on the method-grid (x_coords, y_coords, z_coords)
00232 ! for the extra indices specified in "extra_search%indices".
00233 !
00234 ! There are 3 different of search modes supported by the routine:
00235 ! search_mode = 3: The nearest neighbours are searched are searched
00236 !                  within the 3d region.
00237 ! search_mode = 2: The nearest neighbours are searched only in the
00238 !                  2d hyperplane (lon- and lat-direction).
00239 ! search_mode = 1: The nearest neighbours are searched only in the
00240 !                  1d plane.
00241 !
00242 ! Subroutine "PSMILe_Neigh_nearx_sub_reg_dble" is designed for grids of
00243 ! of type "PRISM_Reglonlatvrt".
00244 !
00245 ! !REVISION HISTORY:
00246 !
00247 !   Date      Programmer   Description
00248 ! ----------  ----------   -----------
00249 ! 11.11.04    H. Ritzdorf  created
00250 !
00251 !EOP
00252 !----------------------------------------------------------------------
00253 !
00254 !  $Id: psmile_neigh_nearx_sub_reg_dble.F90 2325 2010-04-21 15:00:07Z valcke $
00255 !  $Author: valcke $
00256 !
00257    Character(len=len_cvs_string), save :: mycvs = 
00258        '$Id: psmile_neigh_nearx_sub_reg_dble.F90 2325 2010-04-21 15:00:07Z valcke $'
00259 !
00260 !----------------------------------------------------------------------
00261 !
00262 !  Relative control indices 
00263 !      
00264 #ifdef OLD
00265       Integer, Parameter              :: nrefd = 3 * 3 * 3
00266 !
00267       data ((ijkctl (i, n), i=1,ndim_3d), n = 1,nref_3d) &
00268                             /  0, 0, 0,                       &
00269                                1, 0, 0,   0, 1, 0,   0, 0, 1, &
00270                               -1, 0, 0,   0,-1, 0,   0, 0,-1, &
00271                               -1,-1, 0,   1,-1, 0,            &
00272                               -1, 1, 0,   1, 1, 0,            &
00273                               -1, 0,-1,   1, 0,-1,            &
00274                               -1, 0, 1,   1, 0, 1,            &
00275                                0,-1,-1,   0, 1,-1,            &
00276                                0,-1, 1,   0, 1, 1,            &
00277                               -1,-1,-1,   1,-1,-1,            &
00278                               -1, 1,-1,   1, 1,-1, &
00279                               -1,-1, 1,   1,-1, 1, &
00280                               -1, 1, 1,   1, 1, 1/
00281 #endif
00282 ! schoener waere eine liste, die weniger sortierung (s.o.) nachher verlangt
00283 ! veielleicht sollte man immer annehmen, dass ijk in der "Mitte" liegt
00284 !
00285       data ((ijkctl (i, n), i=1,ndim_3d), n = 1,nref_3d) &
00286             /  0, 0, 0,   1, 0, 0,   2, 0, 0,   3, 0, 0, &
00287                0, 1, 0,   1, 1, 0,   2, 1, 0,   3, 1, 0, &
00288                0, 2, 0,   1, 2, 0,   2, 2, 0,   3, 2, 0, &
00289                0, 3, 0,   1, 3, 0,   2, 3, 0,   3, 3, 0, &
00290                0, 0, 1,   1, 0, 1,   2, 0, 1,   3, 0, 1, &
00291                0, 1, 1,   1, 1, 1,   2, 1, 1,   3, 1, 1, &
00292                0, 2, 1,   1, 2, 1,   2, 2, 1,   3, 2, 1, &
00293                0, 3, 1,   1, 3, 1,   2, 3, 1,   3, 3, 1, &
00294                0, 0, 2,   1, 0, 2,   2, 0, 2,   3, 0, 2, &
00295                0, 1, 2,   1, 1, 2,   2, 1, 2,   3, 1, 2, &
00296                0, 2, 2,   1, 2, 2,   2, 2, 2,   3, 2, 2, &
00297                0, 3, 2,   1, 3, 2,   2, 3, 2,   3, 3, 2, &
00298                0, 0, 3,   1, 0, 3,   2, 0, 3,   3, 0, 3, &
00299                0, 1, 3,   1, 1, 3,   2, 1, 3,   3, 1, 3, &
00300                0, 2, 3,   1, 2, 3,   2, 2, 3,   3, 2, 3, &
00301                0, 3, 3,   1, 3, 3,   2, 3, 3,   3, 3, 3/
00302 !
00303       data nrefv /nref_1d, nref_2d, nref_3d/
00304 !                            
00305 !----------------------------------------------------------------------
00306 !
00307 !  Initialization
00308 !
00309 #ifdef VERBOSE
00310       print 9990, trim(ch_id)
00311 
00312       call psmile_flushstd
00313 #endif /* VERBOSE */
00314 !
00315       ierror = 0
00316 
00317       grid_info => Grids(grid_id)
00318       indices   => extra_search%indices
00319       global_search = extra_search%global_search
00320 !
00321       dble_huge = huge (distance(1))
00322 !
00323 !===> Allocate vector to store the distances found if required
00324 !
00325       if (global_search) then
00326          dist_dble => extra_search%dist_dble
00327          dim1 (1) = 1
00328          dim1 (2) = extra_search%n_extra
00329       else
00330          dim1 (1) = jbeg
00331          dim1 (2) = jend
00332 !
00333          Allocate (dist_dble (jbeg:jend, num_neigh), stat = ierror)
00334 
00335          if (ierror /= 0) then
00336             ierrp (1) = ierror
00337             ierrp (2) = num_neigh * (jend-jbeg+1)
00338 
00339             ierror = PRISM_Error_Alloc
00340             call psmile_error ( ierror, 'dist_dble', ierrp, 2, &
00341                 __FILE__, __LINE__ ) 
00342             return
00343          endif
00344       endif
00345 !
00346 #ifdef PRISM_ASERTION
00347       if (global_search) then
00348          if (.not. Associated (extra_search%dist_dble)) then
00349             call psmile_assert (__FILE__, &
00350                 __LINE__, 'extra_search%dist_dble was not allocated')
00351          endif
00352       endif
00353 #endif
00354 !
00355 !===> Get number of reference points
00356 !
00357       nrefd = nrefv (search_mode)
00358 !
00359 !===> Compute distances for each point in (jbeg:jend)
00360 !
00361       do jpart = jbeg, jend, 5000
00362 
00363 !cdir vector loopcnt=5000
00364          do j = jpart, min(jend, jpart+5000-1)
00365 
00366          i = indices(j)
00367 !
00368 !===> ... Get standard range to be controlled
00369 !
00370          irange (1, :) = max (ijk(:, j),          grid_valid_shape (1, :))
00371          irange (2, :) = min (ijk(:, j)+width(:), grid_valid_shape (2, :))
00372 !
00373          nref = (irange (2,1) - irange (1,1) + 1) * &
00374                 (irange (2,2) - irange (1,2) + 1) * &
00375                 (irange (2,3) - irange (1,3) + 1)
00376 !
00377 #ifdef PRISM_ASSERTION
00378 !
00379          if (nref > nrefd) then
00380             print *, trim(ch_id), " nref, nredfd", nref, nrefd
00381             call psmile_assert (__FILE__, &
00382         __LINE__, 'nref > nrefd ')
00383          endif
00384 #endif
00385 !
00386 !===> ... Set indices to be controlled
00387 !
00388          if (nref == nrefd) then
00389                do n = 1, nrefd
00390                ijkdst (:, n) = ijk (:, j) + ijkctl (:, n)
00391                end do
00392                if ( Grids(grid_id)%grid_type == PRISM_Gaussreduced_regvrt ) ijkdst (2, :) = 1
00393          else
00394 ! 
00395             leni  =  irange (2, 1) - irange (1, 1) + 1
00396             lenij = (irange (2, 2) - irange (1, 2) + 1) * leni 
00397 !
00398                do kk = 1, irange (2,3)-irange(1,3) + 1
00399                n = (kk-1)*lenij
00400                   do jj = 1, irange (2,2)-irange(1,2) + 1
00401 
00402 !cdir vector
00403                      do ii = 1, leni
00404                      ijkdst (1, n+ii) = irange(1,1) + ii - 1
00405                      end do ! ii
00406 
00407                   ijkdst (2, n+1:n+leni) = irange (1,2) + jj - 1
00408                   n = n + leni
00409                   end do ! jj
00410 
00411                ijkdst (3, (kk-1)*lenij+1:kk*lenij) = irange (1,3) + kk - 1
00412                end do ! kk
00413                if ( Grids(grid_id)%grid_type == PRISM_Gaussreduced_regvrt ) ijkdst (2, :) = 1
00414          endif
00415 !
00416 !====> ... Get Mask values if necessary
00417 !
00418          if (mask_available) then
00419 !cdir vector
00420                do n = 1, nref
00421                mask_dist (n) = mask_array (ijkdst (1,n), ijkdst (2,n), &
00422                                            ijkdst (3,n))
00423                end do
00424 !
00425             n = count (mask_dist(1:nref))
00426 !
00427             if (n == 0) then
00428 !              no value available
00429 !
00430                nref = 0
00431                go to 10
00432             endif
00433 !
00434             if (n < nref) then
00435                n = 0
00436                   do ii = 1, nref
00437                   if (mask_dist (ii)) then
00438                       n = n + 1
00439                       if (n /= ii) ijkdst (:, n) = ijkdst (:, ii)
00440                   endif
00441                   end do
00442                nref = n
00443                if ( Grids(grid_id)%grid_type == PRISM_Gaussreduced_regvrt ) ijkdst (2, :) = 1
00444             endif
00445          endif
00446 !
00447 !===> ... Compute distances to points
00448 !
00449          fac = dble_earth + z_search (j)
00450 !
00451 !cdir vector
00452             do n = 1, nref
00453             if ( Grids(grid_id)%grid_type == PRISM_Gaussreduced_regvrt ) then
00454             val   = (  sin_values_lat(ijkdst (1,n)) * sin_search(j, lat)    &
00455                     +  cos_values_lat(ijkdst (1,n)) * cos_search(j, lat)    &
00456                     * (cos_values_lon(ijkdst (1,n)) * cos_search(j, lon)    &
00457                       +sin_values_lon(ijkdst (1,n)) * sin_search(j, lon)) )
00458             else
00459             val   = (  sin_values_lat(ijkdst (2,n)) * sin_search(j, lat)    &
00460                     +  cos_values_lat(ijkdst (2,n)) * cos_search(j, lat)    &
00461                     * (cos_values_lon(ijkdst (1,n)) * cos_search(j, lon)    &
00462                       +sin_values_lon(ijkdst (1,n)) * sin_search(j, lon)) )
00463             endif
00464 #ifdef PRISM_ASSERTION
00465             if (val < acosm1 - eps .or. val > acosp1 + eps) then
00466                print *, i, j, val
00467                call psmile_assert (__FILE__, __LINE__, &
00468                                    'invalid value for acos')
00469             endif
00470 #endif
00471             val = min (acosp1, val)
00472             val = max (acosm1, val)
00473 
00474             dist = fac * acos (val)
00475 
00476             distance (n) = dist*dist +                                 &
00477                            ( z_search (j)-z_coords(ijkdst (3,n)))**2
00478 #ifdef DEBUG2
00479             if (i == 2) then
00480                print *, 'n, distance, dist, diff', n, distance (n), dist, &
00481                         z_search (j)-z_coords(ijkdst (3,n))
00482             endif
00483 #endif
00484             end do ! n
00485 !
00486          nmin = minloc (distance(1:nref))
00487 #ifdef MINLOCFIX
00488          if (nmin(1).eq.0) nmin=1
00489 #endif /* MINLOCFIX */
00490 !
00491 #ifdef DEBUG2
00492          if (i == 2) then
00493             print *, 'psmile_neigh_nearx_sub_reg_dble: ictl =', i, j, k
00494             print *, 'ijkdst and distance:'
00495                do n = 1, nref
00496                print *, ijkdst (:, n), distance (n)
00497                end do
00498          endif
00499 #endif
00500 !
00501 !===> ... Special case: Minimal distance == 0
00502 !         Move other neighbours to "outside" of region.
00503 !
00504          if (distance(nmin(1)) <= tol) then
00505             neighbors_3d (:, i, 1) = ijkdst (:, nmin(1))
00506 !
00507                do n = 2, num_neigh
00508                neighbors_3d (:, i, n) = grid_valid_shape (1, :) - 1
00509                end do
00510 !
00511             dist_dble (j, 1:num_neigh) = 0.0
00512 !
00513             cycle
00514          endif
00515 !
00516 !===> ... Sort distances (silly sort)
00517 !         Da muss noch etwas besseres her !
00518 !
00519 10       continue
00520 !
00521          nfound = min (num_neigh, nref)
00522          if (nfound < num_neigh) then
00523             distance(nfound+1:num_neigh) = dble_huge
00524 
00525                do n = nfound+1, num_neigh
00526                ijkdst (:, n) = grid_valid_shape (1, :) - 1
00527                end do
00528          endif
00529 !
00530             do n = 1, nfound
00531 !
00532             nmin = minloc (distance(n:nref))
00533 #ifdef MINLOCFIX
00534             if (nmin(1).eq.0) nmin=1
00535 #endif /* MINLOCFIX */
00536 !
00537             if (nmin(1) /= 1) then
00538                jj = n + nmin (1) - 1
00539 !
00540                itemp  (:)    = ijkdst (:, n)
00541                ijkdst (:, n) = ijkdst (:, jj)
00542                ijkdst (:, jj) = itemp  (:)
00543 !
00544                dist = distance (n)
00545                distance (n) = distance (jj)
00546                distance (jj) = dist
00547             endif
00548             end do ! n
00549 !
00550 #ifdef PRISM_ASSERTION
00551             do n = 1, nfound-1
00552             if (distance(n) > distance (n+1)) exit
00553             end do
00554 !
00555          if (n <= nfound-1) then
00556             print *, 'n =' , n
00557             print *, 'distance (n)  ', distance (n)
00558             print *, 'distance (n+1)', distance (n+1)
00559             print *, 'distance', distance (1:nfound)
00560             call psmile_assert (__FILE__, __LINE__, &
00561                                 'incorrect sort')
00562          endif
00563 #endif
00564 !
00565 !===> ... End of local search;
00566 !         fill neighbors_3d and
00567 !         store distances for global search if necessary
00568 !
00569 ! Brauche ich die sqrt ?
00570 !
00571          dist_dble (j, 1:nfound) = sqrt (distance (1:nfound))
00572          if (nfound < num_neigh) then
00573             dist_dble (j, nfound+1:num_neigh) = dble_huge
00574 
00575 ! das sollte man nicht machen (fuer nlev > 2)
00576 ! um mg suche einen startpunkt zu geben !
00577 !
00578                do n = nfound+1, num_neigh
00579                ijkdst (:, n) = grid_valid_shape (1, :) - 1
00580                end do
00581          endif
00582 !
00583             do n = 1, num_neigh
00584             neighbors_3d (:, i, n) = ijkdst (:, n)
00585             end do
00586 !
00587          end do ! j
00588       end do ! jpart
00589 !
00590 !=======================================================================
00591 !     Use simple search
00592 !=======================================================================
00593 !
00594       if ( Grids(grid_id)%grid_type == PRISM_Gaussreduced_regvrt ) then
00595 
00596          call psmile_srch_nneigh_gauss2_dble (grid_id, nn_srch, arrays, &
00597                              search_mode, nref_3d,                      &
00598                              sin_values_lon, cos_values_lon,            &
00599                              sin_values_lat, cos_values_lat,            &
00600                              grid_valid_shape,                          &
00601                              z_coords, coords_shape,                    &
00602                              neighbors_3d, nloc, num_neigh,             &
00603                              sin_search, cos_search, z_search,          &
00604                              dist_dble, dim1, extra_search, jbeg, jend, &
00605                              mask_array, mask_shape, mask_available,    &
00606                              tol, ierror)
00607          if (ierror > 0) return
00608 
00609       else
00610 !
00611 !=======================================================================
00612 !     Use multigrid info for search in other regions
00613 !=======================================================================
00614 !
00615          if ( grid_info%nlev > 2 ) then
00616 !
00617 !===> Allocate memory and calculate distances on the coarsest level
00618 !     to begin the search for all points.
00619 !
00620             nlev = max(grid_info%nlev-3,2)
00621             nn_srch%leni =  grid_info%mg_infos(nlev)%levdim(1) + 1
00622             nn_srch%lenj =  grid_info%mg_infos(nlev)%levdim(2) + 1
00623             nn_srch%lenk =  grid_info%mg_infos(nlev)%levdim(3) + 1
00624 
00625             leni = nn_srch%leni
00626             lenj = nn_srch%lenj
00627             lenk = nn_srch%lenk
00628 
00629             ii = 0
00630 
00631             allocate(arrays%radius(leni*lenj*lenk), STAT=jj); ii = ii + jj
00632 
00633             allocate(arrays%sin_ccorner_lon(1)%vector(leni), &
00634                  arrays%cos_ccorner_lon(1)%vector(leni), &
00635                  arrays%sin_ccorner_lat(1)%vector(lenj), &
00636                  arrays%cos_ccorner_lat(1)%vector(lenj), STAT=jj); ii = ii + jj
00637 
00638             allocate(arrays%sin_cmidp_lon%vector(leni), &
00639                  arrays%cos_cmidp_lon%vector(leni), &
00640                  arrays%sin_cmidp_lat%vector(lenj), &
00641                  arrays%cos_cmidp_lat%vector(lenj), STAT=jj); ii = ii + jj
00642 
00643             Allocate(arrays%sin_fmidp_lon%vector(maxlen), &
00644                  arrays%cos_fmidp_lon%vector(maxlen), &
00645                  arrays%sin_fmidp_lat%vector(maxlen), &
00646                  arrays%cos_fmidp_lat%vector(maxlen), STAT=jj); ii = ii + jj
00647 
00648             if ( ii /= 0 ) then
00649                ierrp (1) = ierror
00650                ierrp (2) = ii
00651 
00652                ierror = PRISM_Error_Alloc
00653                call psmile_error ( ierror, 'arrays%...%vector', ierrp, 2, &
00654                     __FILE__, __LINE__ )
00655                return
00656 
00657             endif
00658 
00659             my_arrays => grid_info%mg_infos(nlev)%double_arrays
00660 
00661             arrays%sin_ccorner_lon(1)%vector(1:leni) = sin(my_arrays%chmin(1)%vector(1:leni)*dble_deg2rad)
00662             arrays%cos_ccorner_lon(1)%vector(1:leni) = cos(my_arrays%chmin(1)%vector(1:leni)*dble_deg2rad)
00663             arrays%sin_ccorner_lat(1)%vector(1:lenj) = sin(my_arrays%chmin(2)%vector(1:lenj)*dble_deg2rad)
00664             arrays%cos_ccorner_lat(1)%vector(1:lenj) = cos(my_arrays%chmin(2)%vector(1:lenj)*dble_deg2rad)
00665 
00666             arrays%sin_cmidp_lon%vector(1:leni) = sin(my_arrays%midp(1)%vector(1:leni)*dble_deg2rad)
00667             arrays%cos_cmidp_lon%vector(1:leni) = cos(my_arrays%midp(1)%vector(1:leni)*dble_deg2rad)
00668             arrays%sin_cmidp_lat%vector(1:lenj) = sin(my_arrays%midp(2)%vector(1:lenj)*dble_deg2rad)
00669             arrays%cos_cmidp_lat%vector(1:lenj) = cos(my_arrays%midp(2)%vector(1:lenj)*dble_deg2rad)
00670 
00671             n = 0
00672 
00673             do k = 1, lenk
00674                fac = dble_earth + my_arrays%midp(3)%vector(k)
00675                do j = 1, lenj
00676                   do i = 1, leni
00677                      n = n+1
00678 
00679                      val = ( arrays%sin_ccorner_lat(1)%vector(j) * arrays%sin_cmidp_lat%vector(j)    &
00680                           +  arrays%cos_ccorner_lat(1)%vector(j) * arrays%cos_cmidp_lat%vector(j)    &
00681                           * (arrays%cos_ccorner_lon(1)%vector(i) * arrays%cos_cmidp_lon%vector(i)    &
00682                           +  arrays%sin_ccorner_lon(1)%vector(i) * arrays%sin_cmidp_lon%vector(i)) )
00683 
00684                      val = min (acosp1, val)
00685                      val = max (acosm1, val)
00686 
00687                      dist = fac * acos (val)
00688 
00689                      if ( search_mode == 3 ) then
00690                         arrays%radius(n) = sqrt(dist*dist+(my_arrays%midp(3)%vector(k)- &
00691                              my_arrays%chmin(3)%vector(k))**2)
00692                      else if ( search_mode == 2 ) then
00693                         arrays%radius(n) = abs(dist)
00694                      endif
00695 
00696                   enddo
00697                enddo
00698             enddo
00699 !
00700 !===> ... Use multigrid info for search in orther regions
00701 !
00702             call psmile_mg_srch_nneigh_reg_dble (grid_id, nn_srch,      &
00703                       arrays, search_mode, nref_3d, grid_valid_shape,   &
00704                       neighbors_3d, nloc, num_neigh,                    &
00705                       sin_search, cos_search, z_search,                 &
00706                       dist_dble, dim1, indices, jbeg, jend,             &
00707                       mask_array, mask_shape, mask_available,           &
00708                  tol, ierror)
00709             if (ierror > 0) return
00710 !
00711 !===> ... Deallocate memory
00712 !
00713             deallocate(arrays%sin_ccorner_lon(1)%vector, &
00714                  arrays%cos_ccorner_lon(1)%vector, &
00715                  arrays%sin_ccorner_lat(1)%vector, &
00716                  arrays%cos_ccorner_lat(1)%vector)
00717 
00718             deallocate(arrays%sin_cmidp_lon%vector, &
00719                  arrays%cos_cmidp_lon%vector, &
00720                  arrays%sin_cmidp_lat%vector, &
00721                  arrays%cos_cmidp_lat%vector)
00722 
00723             deallocate(arrays%sin_fmidp_lon%vector, &
00724                  arrays%cos_fmidp_lon%vector, &
00725                  arrays%sin_fmidp_lat%vector, &
00726                  arrays%cos_fmidp_lat%vector)
00727 
00728             Deallocate(arrays%radius)
00729          endif
00730 !
00731          if (.not. global_search) then
00732             Deallocate (dist_dble)
00733          endif
00734 
00735       endif
00736 !
00737 !===> All done
00738 !
00739 #ifdef VERBOSE
00740       print 9980, trim(ch_id), ierror
00741 
00742       call psmile_flushstd
00743 #endif /* VERBOSE */
00744 !
00745 !  Formats:
00746 !
00747 9990 format (1x, a, ': psmile_neigh_nearx_sub_reg_dble')
00748 9980 format (1x, a, ': psmile_neigh_nearx_sub_reg_dble: eof, ierror =', i3)
00749 
00750       end subroutine PSMILe_Neigh_nearx_sub_reg_dble

Generated on 18 Mar 2011 for Oasis4 by  doxygen 1.6.1