psmile_neigh_nearx_3d_reg_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_nearx_3d_reg_real
00008 !
00009 ! !INTERFACE:
00010 
00011       subroutine psmile_neigh_nearx_3d_reg_real (          &
00012                            grid_id,                        &
00013                            coords1, coords2, coords3,      &
00014                            x_coords, y_coords, z_coords,   &
00015                            coords_shape,                   &
00016                            mask_array, mask_shape,  mask_available,    &
00017                            sin_values_lon, cos_values_lon, &
00018                            sin_values_lat, cos_values_lat, &
00019                            grid_valid_shape, search_mode,  &
00020                            srclocs, nlocs, nloc, nprev,    &
00021                            neighbors_3d, num_neigh,        &
00022                            extra_search, ierror)
00023 !
00024 ! !USES:
00025 !
00026       use PRISM_constants
00027 !
00028       use PSMILe, dummy_interface => PSMILe_Neigh_nearx_3d_reg_real
00029 
00030       implicit none
00031 !
00032 ! !INPUT PARAMETERS:
00033 !
00034       Integer, Intent (In)            :: grid_id
00035 
00036 !     Link to the info on the component in which the donor cells
00037 !     should be searched.
00038 
00039       Integer, Intent (In)            :: nloc
00040 !
00041 !     Total number of locations to be transferred
00042 !
00043       Integer, Intent (In)            :: nlocs (ndim_3d)
00044 !
00045 !     Total  Number of locations to be transferred
00046 
00047       Integer, Intent (In)            :: nprev
00048 !
00049 !     Number of locations for which neighbours were already searched.
00050 !
00051       Type (integer_vector), Intent (In) :: srclocs (ndim_3d)
00052 
00053 !     Indices of the grid cell in which the point was found.
00054 !
00055       Real, Intent (In)               :: coords1 (nlocs(1),nlocs(2),nlocs(3))
00056       Real, Intent (In)               :: coords2 (nlocs(1),nlocs(2),nlocs(3))
00057       Real, Intent (In)               :: coords3 (nlocs(1),nlocs(2),nlocs(3))
00058 
00059 !     Coordinates of the target points (which were searched and found)
00060 !     Note: Target grids are full arrays which were generated 
00061 !           out of partial arrays.
00062 
00063       Integer, Intent (In)            :: coords_shape (2, ndim_3d)
00064 
00065 !     Dimension of coordinates (method) x_coords, ...
00066 
00067       Real, Intent (In)               :: x_coords(coords_shape(1,1): 
00068                                                   coords_shape(2,1))
00069       Real, Intent (In)               :: y_coords(coords_shape(1,2): 
00070                                                   coords_shape(2,2))
00071       Real, Intent (In)               :: z_coords(coords_shape(1,3): 
00072                                                   coords_shape(2,3))
00073 
00074 !     Coordinates of the method
00075 
00076       Integer, Intent (In)            :: mask_shape (2, ndim_3d)
00077 
00078 !     Dimension of (method) mask array
00079 
00080       Logical, Intent (In)            :: mask_array (mask_shape (1,1): 
00081                                                      mask_shape (2,1), 
00082                                                      mask_shape (1,2): 
00083                                                      mask_shape (2,2), 
00084                                                      mask_shape (1,3): 
00085                                                      mask_shape (2,3))
00086 !     Mask of the method
00087 
00088       Logical, Intent (In)            :: mask_available
00089 
00090 !     Is mask specified in array "mask_array" ?
00091 !     mask_available = .false. : Mask is not available
00092 !                      .true.  : Mask is     available
00093 
00094       Integer, Intent (In)            :: grid_valid_shape (2, ndim_3d)
00095 
00096 !     Specifies the valid block shape for the "ndim_3d"-dimensional block
00097 
00098       Real, Intent (In)               :: sin_values_lon (grid_valid_shape(1,1): 
00099                                                          grid_valid_shape(2,1))
00100       Real, Intent (In)               :: sin_values_lat (grid_valid_shape(1,2): 
00101                                                          grid_valid_shape(2,2))
00102 !
00103 !     Sin values of Longitudes and Latitudes (x_coords, y_coords)
00104 !
00105       Real, Intent (In)               :: cos_values_lon (grid_valid_shape(1,1): 
00106                                                          grid_valid_shape(2,1))
00107       Real, Intent (In)               :: cos_values_lat (grid_valid_shape(1,2): 
00108                                                          grid_valid_shape(2,2))
00109 !    
00110 !     Cos values of Longitudes and Latitudes (x_coords, y_coords)
00111 !
00112       Integer, Intent (In)            :: search_mode
00113 !
00114 !     Specifies the search mode for nearest neigbours with
00115 !        search_mode = 3 : Full 3d-search
00116 !        search_mode = 2 : Search in 2d-hyperplane (1st 2 direction lon, lat)
00117 !        search_mode = 1 : Search in 1d-plane
00118 !
00119       Integer, Intent (In)            :: num_neigh
00120 !
00121 !     Last dimension of neighbors array "neighbors_3d" and
00122 !     number of neighbors to be searched.
00123 !
00124 ! !INPUT/OUTPUT PARAMETERS:
00125 !
00126       Type (Extra_search_info)        :: extra_search
00127 !
00128 !     Number of locations where
00129 !       (*) required mask values were not true
00130 !
00131 ! !OUTPUT PARAMETERS:
00132 !
00133       Integer, Intent (Out)           :: neighbors_3d (ndim_3d, nloc, num_neigh)
00134 
00135 !     Indices of neighbor locations (interpolation bases)
00136 !
00137       Integer, Intent (Out)           :: ierror
00138 
00139 !     Returns the error code of PSMILE_Neigh_nearx_3d_reg_real;
00140 !             ierror = 0 : No error
00141 !             ierror > 0 : Severe error
00142 !
00143 ! !DEFINED PARAMETERS:
00144 !
00145 !  lon   = Index of Longitudes in arrays "sin_values" and "cos_values"
00146 !  lat   = Index of Latitudes  in arrays "sin_values" and "cos_values"
00147 !
00148       Integer, Parameter              :: lon = 1
00149       Integer, Parameter              :: lat = 2
00150 !
00151 ! !LOCAL VARIABLES
00152 !
00153 !     ... For locations controlled
00154 !
00155       Real, Pointer                   :: sin_search (:, :)
00156       Real, Pointer                   :: cos_search (:, :)
00157       Real, Pointer                   :: z_search (:)
00158 !
00159 !     ... (1d) Indices in (entire) data to be extra searched
00160 !
00161 #ifdef PRISM_ASSERTION
00162       Integer                         :: i
00163 #endif
00164       Integer                         :: j
00165       Integer                         :: ibeg, iend
00166       Integer                         :: jbeg, jend
00167       Integer                         :: lenij
00168 !
00169       Integer, Pointer                :: indices (:)
00170       Integer                         :: nprev1
00171 !
00172 !     ... 3d Indices in local data to be extra searched
00173 !
00174       Integer                         :: ijk0 (ndim_3d)
00175       Integer, Pointer                :: ijk (:, :)
00176       Integer                         :: width (ndim_3d)
00177 !
00178 !     ... for error handling
00179 !
00180       Integer, Parameter              :: nerrp = 2
00181       Integer                         :: ierrp (nerrp)
00182 !
00183 !
00184 ! !DESCRIPTION:
00185 !
00186 ! Subroutine "PSMILe_Neigh_nearx_3d_reg_real" searches the next "num_neigh"
00187 ! nearest neighbours on the method-grid (x_coords, y_coords, z_coords)
00188 ! for the extra indices specified in "extra_search%indices".
00189 !
00190 ! There are 3 different of search modes supported by the routine:
00191 ! search_mode = 3: The nearest neighbours are searched are searched
00192 !                  within the 3d region.
00193 ! search_mode = 2: The nearest neighbours are searched only in the
00194 !                  2d hyperplane (lon- and lat-direction).
00195 ! search_mode = 1: The nearest neighbours are searched only in the
00196 !                  1d plane.
00197 !
00198 ! Subroutine "PSMILe_Neigh_nearx_3d_reg_real" is designed
00199 ! (*) for source grids of type "PRISM_Reglonlatvrt" and
00200 ! (*) source locations "srcloc" which are stored as a 3 single 1d vectors
00201 !     containing the 3d indices of the cells found for the target points
00202 !     "coords" to be searched.
00203 !
00204 ! !REVISION HISTORY:
00205 !
00206 !   Date      Programmer   Description
00207 ! ----------  ----------   -----------
00208 ! 03.07.21    H. Ritzdorf  created
00209 !
00210 !EOP
00211 !----------------------------------------------------------------------
00212 !
00213 !  $Id: psmile_neigh_nearx_3d_reg_real.F90 2325 2010-04-21 15:00:07Z valcke $
00214 !  $Author: valcke $
00215 !
00216    Character(len=len_cvs_string), save :: mycvs = 
00217        '$Id: psmile_neigh_nearx_3d_reg_real.F90 2325 2010-04-21 15:00:07Z valcke $'
00218 !                            
00219 !----------------------------------------------------------------------
00220 !
00221 !  Initialization
00222 !
00223 #ifdef VERBOSE
00224       print 9990, trim(ch_id)
00225 
00226       call psmile_flushstd
00227 #endif /* VERBOSE */
00228 !
00229       ierror = 0
00230 !
00231       nprev1 = nprev + 1
00232 !
00233       ibeg = nprev + 1
00234       iend = nprev + nlocs(1)*nlocs(2)*nlocs(3)
00235 !
00236       indices => extra_search%indices
00237 !
00238 !===> Get range in indices to be controlled
00239 !     Note: indices are locations in the full 3d cube and
00240 !           contains the previous numbers
00241 !           TODO: DAS SOLLTE in lokale Indices geaendert werden
00242 !
00243          do jbeg = 1, extra_search%n_extra
00244          if (indices(jbeg) >= ibeg) exit
00245          end do
00246 !
00247       if (jbeg > extra_search%n_extra) return
00248 !
00249          do jend = jbeg, extra_search%n_extra
00250          if (indices(jend) > iend) exit
00251          end do
00252       jend = jend - 1
00253 !
00254       if (jbeg > jend) return
00255 !
00256       ibeg = indices(jbeg)
00257       iend = indices(jend)
00258 !
00259 #ifdef PRISM_ASSERTION
00260 !
00261 !===> Internal control
00262 !
00263       if (search_mode < 1 .or. search_mode > ndim_3d) then
00264          print *, trim(ch_id), " search_mode ", search_mode
00265          call psmile_assert (__FILE__, __LINE__, &
00266                              'search_mode is out of range ')
00267       endif
00268 !
00269       if (nloc < nprev + nlocs (1) * nlocs (2) * nlocs (3)) then
00270          print *, trim(ch_id), " nloc, nprev, nlocs ", nloc, nprev, nlocs
00271          call psmile_assert (__FILE__, __LINE__, &
00272                              'nloc < nprev + PRODUCT (nlocs) ')
00273       endif
00274 !
00275 !===> Are the locations within the correct shape ?
00276 !     Note: the scrloc's are source locations in the
00277 !           method grid which has an virtual cell 0.
00278 !
00279          do j = 1, ndim_3d
00280 !cdir vector
00281             do i = 1, nlocs (j)
00282 !
00283             if (srclocs(j)%vector(i) < grid_valid_shape(1,j)-1 .or. &
00284                 srclocs(j)%vector(i) > grid_valid_shape(2,j)) exit
00285             end do
00286 !
00287          if (i <= nlocs(j)) then
00288             print *, "Incorrect index in j, i", j, i, srclocs(j)%vector(i)
00289             call psmile_assert (__FILE__, __LINE__, &
00290                                 'wrong index in srcloc')
00291          endif
00292          end do
00293 #endif
00294 !
00295 !-----------------------------------------------------------------------
00296 !     Allocate and set temporary array which is used to transform
00297 !     the coordinates of points 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(jbeg:jend, lat), &
00307                 cos_search(jbeg:jend, lat), &
00308                   z_search(jbeg:jend),      STAT = ierror)
00309 
00310       if ( ierror > 0 ) then
00311          ierrp (1) = ierror
00312          ierrp (2) = (jend-jbeg+1) * (lat * 2 + 1)
00313 
00314          ierror = PRISM_Error_Alloc
00315          call psmile_error ( ierror, 'sin_search/cos_search/z_search', &
00316                              ierrp, 2, __FILE__, __LINE__ )
00317          return
00318       endif
00319 !
00320       Allocate (ijk(ndim_3d, jbeg:jend), STAT = ierror)
00321 
00322       if ( ierror > 0 ) then
00323          ierrp (1) = ierror
00324          ierrp (2) = (jend-jbeg+1) * ndim_3d
00325 
00326          ierror = PRISM_Error_Alloc
00327          call psmile_error ( ierror, 'ijk', &
00328                              ierrp, 2, __FILE__, __LINE__ )
00329          return
00330       endif
00331 !
00332 !===> Compute 3d Indices
00333 !
00334       lenij = nlocs(1) * nlocs (2)
00335 !
00336       ijk (3, jbeg:jend) =    (indices(jbeg:jend)-nprev1)/lenij + 1
00337       ijk (2, jbeg:jend) =   ((indices(jbeg:jend)-nprev1) - &
00338                                 (ijk(3,jbeg:jend)-1)*lenij)/nlocs(1) + 1
00339       ijk (1, jbeg:jend) = mod(indices(jbeg:jend)-nprev1, nlocs(1)) + 1
00340 
00341 #ifdef PRISM_ASSERTION
00342 !cdir vector
00343          do j = jbeg, jend
00344          if (ijk(1,j) < 1 .or. ijk(1,j) > nlocs(1) .or. &
00345              ijk(2,j) < 1 .or. ijk(2,j) > nlocs(2) .or. &
00346              ijk(3,j) < 1 .or. ijk(3,j) > nlocs(3)) exit
00347          end do
00348 
00349       if (j <= jend) then
00350          print *, "Incorrect index in ijk", j, jbeg, jend, ijk(:, j)
00351          call psmile_assert (__FILE__, __LINE__, &
00352                                 'wrong index in ijk')
00353       endif
00354 #endif
00355 !
00356 !===> Compute values required to compute distances
00357 !     Note: Target grids are full arrays which were generated 
00358 !           out of partial arrays.
00359 !
00360 !cdir vector
00361          do j = jbeg, jend
00362          sin_search (j, lon) = coords1 (ijk(1,j), ijk(2,j), ijk(3,j)) &
00363                              * real_deg2rad
00364 !        sin_search (j, lon) = coords1 (ijk(1,j), 1, 1) * real_deg2rad
00365          end do
00366 !
00367 !cdir vector
00368          do j = jbeg, jend
00369          sin_search (j, lat) = coords2 (ijk(1,j), ijk(2,j), ijk(3,j)) &
00370                              * real_deg2rad
00371 !        sin_search (j, lon) = coords1 (1, ijk(2,j), 1) * real_deg2rad
00372          end do
00373 !
00374       cos_search = cos (sin_search)
00375       sin_search = sin (sin_search)
00376 !
00377 !cdir vector
00378          do j = jbeg, jend
00379          z_search (j) = coords3 (ijk(1,j), ijk(2,j), ijk(3,j))
00380 !        z_search (j) = coords3 (1, 1, ijk(3,j))
00381          end do
00382 !
00383 !-----------------------------------------------------------------------
00384 !     Compute coarse grid indices
00385 !-----------------------------------------------------------------------
00386 !
00387       ijk0 = Grids(grid_id)%ijk0
00388 !
00389 !===> Get coarse grid index on level nlev-2.
00390 !     By coarsening with a faktor of 4
00391 !
00392       width (1:search_mode) = 3
00393       width (search_mode+1:ndim_3d) = 0
00394 !
00395       ijk (1, jbeg:jend) = ijk0(1) + ((srclocs(1)%vector(ijk(1, jbeg:jend)) - &
00396                                        ijk0(1)) / 4) * 4
00397 !
00398       if (search_mode >= 2) then
00399          ijk (2, jbeg:jend) = ijk0(2) + ((srclocs(2)%vector(ijk(2, jbeg:jend)) - &
00400                                          ijk0(2)) / 4) * 4
00401       else
00402          ijk (2, jbeg:jend) = srclocs(2)%vector(ijk(2, jbeg:jend))
00403       endif
00404 !
00405       if (search_mode >= 3) then
00406          ijk (3, jbeg:jend) = ijk0(3) + ((srclocs(3)%vector(ijk(3, jbeg:jend)) - &
00407                                          ijk0(3)) / 4) * 4
00408       else
00409          ijk (3, jbeg:jend) = srclocs(3)%vector(ijk(3, jbeg:jend))
00410       endif
00411 !
00412 !===> Compute nearest neigbours for all points
00413 !     specified in "indices(jbeg:jend)"
00414 !
00415       call psmile_neigh_nearx_sub_reg_real (               &
00416                            grid_id,                        &
00417                            x_coords, y_coords, z_coords,   &
00418                            coords_shape,                   &
00419                            mask_array, mask_shape,  mask_available,    &
00420                            sin_values_lon, cos_values_lon, &
00421                            sin_values_lat, cos_values_lat, &
00422                            grid_valid_shape, search_mode,  &
00423                            neighbors_3d, num_neigh, nloc,  &
00424                            extra_search,                   &
00425                            ijk, sin_search, cos_search,    &
00426                            z_search, jbeg, jend, width, ierror)
00427 !
00428 !===> Deallocate memory
00429 !
00430       Deallocate (ijk)
00431       Deallocate (cos_search, sin_search, z_search)
00432 !
00433 !===> All done
00434 !
00435 #ifdef VERBOSE
00436       print 9980, trim(ch_id), ierror
00437 
00438       call psmile_flushstd
00439 #endif /* VERBOSE */
00440 !
00441 !  Formats:
00442 !
00443 9990 format (1x, a, ': psmile_neigh_nearx_3d_reg_real')
00444 9980 format (1x, a, ': psmile_neigh_nearx_3d_reg_real: eof, ierror =', i3)
00445 
00446       end subroutine PSMILe_Neigh_nearx_3d_reg_real

Generated on 18 Mar 2011 for Oasis4 by  doxygen 1.6.1