psmile_mg_search_2d_real.F90

Go to the documentation of this file.
00001 !-----------------------------------------------------------------------
00002 ! Copyright 2007-2010, NEC Europe Ltd., London, UK.
00003 ! All rights reserved. Use is subject to OASIS4 license terms.
00004 !-----------------------------------------------------------------------
00005 !BOP
00006 !
00007 ! !ROUTINE: PSMILe_mg_search_2d_real
00008 !
00009 ! !INTERFACE:
00010 
00011       subroutine psmile_mg_search_2d_real (grid_id,         &
00012                         found, locations, len, search_data, &
00013                         ipart, tol, ierror)
00014 !
00015 ! !USES:
00016 !
00017       use PRISM_constants
00018 !
00019       use PSMILe, dummy_interface => PSMILe_mg_search_2d_real
00020 
00021       implicit none
00022 !
00023 ! !INPUT PARAMETERS:
00024 !
00025       Integer, Intent (In)            :: grid_id
00026 !
00027 !     Grid id
00028 !
00029       Integer, Intent (In)            :: len
00030 
00031 !     Number of coords to be searched
00032 !
00033       Type (Enddef_search_data)       :: search_data
00034 
00035 !     Info's on coordinates to be searched
00036 
00037       Integer, Intent (In)            :: ipart
00038 
00039 !     Number of the partition to be searched
00040 
00041       Real, Intent (In)               :: tol
00042 !
00043 !     Absolute tolerance for search of "identical" points
00044 !     TOL >= 0.0
00045 !
00046 !
00047 ! !INPUT/OUTPUT PARAMETERS:
00048 !
00049       Integer, Intent (InOut)         :: found (len)
00050 
00051 !     Finest level number on which a grid cell was found for point I.
00052 !     Level number = -(nlev+1): Never found (input value)
00053 
00054       Integer, Intent (InOut)         :: locations (ndim_2d, len)
00055 
00056 !     Indices of the grid cell in which the point was found.
00057 !     Assumed input value locations (:, len) = 0
00058 !
00059 ! !OUTPUT PARAMETERS:
00060 !
00061       integer, Intent (Out)               :: ierror
00062 
00063 !     Returns the error code of PSMILe_mg_search_2d_real;
00064 !             ierror = 0 : No error
00065 !             ierror > 0 : Severe error
00066 !
00067 ! !LOCAL VARIABLES
00068 !
00069       Type (Corner_Block), Pointer :: corner_pointer
00070       Integer                      :: i
00071 !
00072 !     ... for locations searched
00073 !
00074       Integer                      :: ibeg, iend
00075       Integer                      :: control (2, ndim_3d)
00076 !
00077       Real, Pointer                :: tgt_coords_x (:)
00078       Real, Pointer                :: tgt_coords_y (:)
00079       Integer                      :: shape (2, ndim_3d)
00080       Integer                      :: range (2, ndim_3d)
00081       Integer                      :: j, stride
00082 !
00083 !     ... for levels
00084 !
00085       Integer                      :: lev, levc, nlev
00086       Integer                      :: ijkinc (ndim_3d), ijkcoa (ndim_2d)
00087       Type(Grid), Pointer          :: grid_info
00088 !
00089 !     ... for error parameters
00090 !
00091       Integer, parameter           :: nerrp = 2
00092       Integer                      :: ierrp (nerrp)
00093 !
00094 ! !DESCRIPTION:
00095 !
00096 ! Subroutine "PSMILe_mg_search_2d_real" searches the donor
00097 ! for the subgrid sent by the sending process.
00098 !
00099 !
00100 ! !REVISION HISTORY:
00101 !   Date      Programmer   Description
00102 ! ----------  ----------   -----------
00103 ! 03.07.21    H. Ritzdorf  created
00104 !
00105 !EOP
00106 !----------------------------------------------------------------------
00107 !
00108 ! $Id: psmile_mg_search_2d_real.F90 3175 2011-05-04 10:48:06Z hanke $
00109 ! $Author: hanke $
00110 !
00111    Character(len=len_cvs_string), save :: mycvs = 
00112        '$Id: psmile_mg_search_2d_real.F90 3175 2011-05-04 10:48:06Z hanke $'
00113 !
00114 !----------------------------------------------------------------------
00115 !
00116 !  Initialization
00117 !
00118 #ifdef VERBOSE
00119       print 9990, trim(ch_id), Grids(grid_id)%comp_id
00120 
00121       call psmile_flushstd
00122 #endif /* VERBOSE */
00123 
00124       ierror = 0
00125       grid_info => Grids(grid_id)
00126 !
00127 !-----------------------------------------------------------------------
00128 !     Generate 2d array if necessary
00129 !-----------------------------------------------------------------------
00130 !
00131       if (search_data%grid_type == PRISM_Reglonlatvrt) then
00132 !        ... Dimension in 1st and 2nd direction is only 1-dimensional
00133 !        ... Generate 2d arrays
00134 !
00135          Allocate (tgt_coords_x(len), tgt_coords_y(len), STAT = ierror)
00136          if ( ierror > 0 ) then
00137             ierrp (1) = ierror
00138             ierrp (2) = 2* len
00139             ierror = PRISM_Error_Alloc
00140             call psmile_error ( ierror, 'tgt_coords_x, tgt_coords_y', &
00141                                 ierrp, 2, __FILE__, __LINE__ )
00142             return
00143          endif
00144 !
00145          shape (:, 1:ndim_2d) = search_data%range(:, 1:ndim_2d, ipart)
00146          shape (:,   ndim_3d) = 1
00147 !
00148          range = shape
00149 
00150          stride = shape(2,1)-shape(1,1) + 1
00151 #ifdef PRISM_ASSERTION
00152          if (search_data%dim_size(1, ipart) /= stride) then
00153             call psmile_assert (__FILE__, __LINE__, &
00154                                 "dim(1) != stride")
00155          endif
00156 #endif
00157 !
00158          do j = 1, shape(2,2)-shape(1,2) + 1
00159             tgt_coords_x ((j-1)*stride+1:j*stride) = search_data%search_dble(1,ipart)%vector(1:stride)
00160             tgt_coords_y ((j-1)*stride+1:j*stride) = search_data%search_dble(2,ipart)%vector(j)
00161          end do
00162 
00163       else
00164 
00165          tgt_coords_x => search_data%search_real(1,ipart)%vector
00166          tgt_coords_y => search_data%search_real(2,ipart)%vector
00167 !
00168          if (search_data%grid_type == PRISM_Irrlonlat_regvrt .or. &
00169              search_data%grid_type == PRISM_Gaussreduced_regvrt) then
00170 
00171             shape(:, 1:ndim_2d) = search_data%shape(:, 1:ndim_2d, ipart)
00172             shape(:,   ndim_3d) = 1
00173 
00174             range(:, 1:ndim_2d) = search_data%range(:, 1:ndim_2d, ipart)
00175             range(:,   ndim_3d) = 1
00176          else
00177             shape(:, :) = search_data%shape(:, :, ipart)
00178             range(:, :) = search_data%range(:, :, ipart)
00179          endif
00180       endif
00181 !
00182 !-----------------------------------------------------------------------
00183 !     Coarsest level nlev
00184 !-----------------------------------------------------------------------
00185 !
00186       nlev = grid_info%nlev
00187       lev = nlev
00188 !
00189       ibeg = 1
00190       iend = len
00191 
00192 #ifdef PRISM_ASSERTION
00193       if (grid_info%mg_infos(lev)%levdim(1) /= 0 .or. &
00194           grid_info%mg_infos(lev)%levdim(2) /= 0) then
00195 
00196          call psmile_assert (__FILE__, __LINE__, &
00197                              "coarsest level dim != 0")
00198       endif
00199 #endif
00200 !
00201       call psmile_mg_coarse_2d_real (lev, &
00202                      grid_info%mg_infos(lev)%real_arrays%chmin,    &
00203                      grid_info%mg_infos(lev)%real_arrays%chmax,    &
00204                      found, locations, tgt_coords_x, tgt_coords_y, &
00205                      ibeg, iend)
00206 !
00207       if (ibeg > iend) then
00208 
00209 #ifdef VERBOSE
00210          print 9980, trim(ch_id), grid_info%comp_id, ierror
00211 
00212          call psmile_flushstd
00213 #endif /* VERBOSE */
00214          return
00215       endif
00216 
00217       control = range
00218 !
00219 !===> Multiple grids
00220 !
00221 !     ijkinc (1) = max (1, (iend-ibeg)/2)
00222 !     ijkinc (2) = max (1, (jend-jbeg)/2)
00223 !
00224       ijkinc (:) = 1
00225 !
00226       do lev = nlev-1, 1, -1
00227 !
00228          levc = lev + 1
00229 !
00230          ijkcoa (:) = 2
00231 !
00232          do i = 1, ndim_2d
00233             if (grid_info%mg_infos(lev )%levdim(i) == &
00234                 grid_info%mg_infos(levc)%levdim(i)) then
00235                ijkcoa (i) = 1
00236             endif
00237          enddo
00238 !
00239 ! chmin, chmax, midp als 4d array alloziieren !
00240 !
00241          call psmile_mg_next_level_2d_real ( grid_id, lev, nlev,     &
00242                 grid_info%mg_infos(lev)%real_arrays%chmin(1)%vector, &
00243                 grid_info%mg_infos(lev)%real_arrays%chmin(2)%vector, &
00244                 grid_info%mg_infos(lev)%real_arrays%chmax(1)%vector, &
00245                 grid_info%mg_infos(lev)%real_arrays%chmax(2)%vector, &
00246                 grid_info%mg_infos(lev)%real_arrays%midp(1)%vector,  &
00247                 grid_info%mg_infos(lev)%real_arrays%midp(2)%vector,  &
00248                 grid_info%mg_infos(lev)%levdim,                      &
00249                 found, locations, range,                             &
00250                 tgt_coords_x, tgt_coords_y, shape,                   &
00251                 control, ijkinc, ijkcoa, ierror)
00252 
00253          if (ierror > 0) return
00254 !
00255 !        ijkinc (:) = max (1, ijkinc(:)/ijkcoa(:))
00256       enddo ! lev
00257 !
00258 !   Transfer  locations from mg locations (range (0:levdim(:)))
00259 !   into user locations (i.e. locations in corner/method dimensions)
00260 !
00261       do i = 1, len
00262          locations (1:ndim_2d, i) = locations (1:ndim_2d, i)  &
00263                                   + grid_info%ijk0 (1:ndim_2d)
00264       end do
00265 !
00266 #define SEARCH_ON_FINAL_GRID
00267 #ifdef SEARCH_ON_FINAL_GRID
00268 !
00269 !   Get final locations controlling the real cells
00270 !   Probably not efficient
00271 !
00272       corner_pointer => Grids(grid_id)%corner_pointer
00273 !
00274       if (corner_pointer%nbr_corners == 8) then
00275 !        4 corners in 2 direction, 2 corners in vertical direction
00276 !
00277          call psmile_mg_final_2d_real (grid_id, nlev,                 &
00278                    grid_info%mg_infos(1)%real_arrays%chmin(1)%vector, &
00279                    grid_info%mg_infos(1)%real_arrays%chmin(2)%vector, &
00280                    grid_info%mg_infos(1)%real_arrays%chmax(1)%vector, &
00281                    grid_info%mg_infos(1)%real_arrays%chmax(2)%vector, &
00282                    grid_info%mg_infos(1)%real_arrays%midp(1)%vector,  &
00283                    grid_info%mg_infos(1)%real_arrays%midp(2)%vector,  &
00284                    grid_info%mg_infos(1)%levdim,                      &
00285                    found, locations, range,                           &
00286                    tgt_coords_x, tgt_coords_y, shape, control,        &
00287                    corner_pointer%corners_real(1)%vector,             &
00288                    corner_pointer%corners_real(2)%vector,             &
00289                    corner_pointer%corner_shape,                       &
00290                    corner_pointer%nbr_corners/2,                      &
00291                    tol, ierror)
00292          if (ierror > 0) return
00293       endif
00294 #endif
00295 !
00296 !-----------------------------------------------------------------------
00297 !     All done; Free memory allocated
00298 !-----------------------------------------------------------------------
00299 !
00300       if (search_data%grid_type == PRISM_Reglonlatvrt) then
00301          Deallocate (tgt_coords_y)
00302          Deallocate (tgt_coords_x)
00303       endif
00304 !
00305 #ifdef VERBOSE
00306       print 9980, trim(ch_id), grid_info%comp_id, ierror
00307 
00308       call psmile_flushstd
00309 #endif /* VERBOSE */
00310 !
00311 !  Formats:
00312 !
00313 9990 format (1x, a, ': PSMILe_mg_search_2d_real: comp_id =', i3)
00314 9980 format (1x, a, ': PSMILe_mg_search_2d_real: comp_id =', i3, '; eof ierror =', i4)
00315 
00316       end subroutine PSMILe_mg_search_2d_real

Generated on 1 Dec 2011 for Oasis4 by  doxygen 1.6.1