psmile_search_donor_1d_dble.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_Search_donor_1d_dble
00008 !
00009 ! !INTERFACE:
00010 
00011       subroutine psmile_search_donor_1d_dble (grid_id, idim, &
00012                         found, locations, coords, len, tol, ierror)
00013 !
00014 ! !USES:
00015 !
00016       use PRISM_constants
00017 !
00018       use PSMILe, dummy_interface => PSMILe_Search_donor_1d_dble
00019 
00020       implicit none
00021 !
00022 ! !INPUT PARAMETERS:
00023 !
00024       Integer, Intent (In)            :: grid_id
00025 
00026 !     Info on the component in which the donor
00027 !     should be searched.
00028 
00029       Integer, Intent (In)            :: idim
00030 
00031 !     Dimension in which should be searched
00032 
00033       Integer, Intent (In)            :: len
00034 
00035 !     Number of coords to be searched
00036 
00037       Double Precision, Intent (In)   :: coords (len)
00038 
00039 !     Coordinates to be searched
00040 
00041       Double Precision, Intent (In)   :: tol
00042 !
00043 !     Absolute tolerance for search of "identical" points
00044 !     TOL >= 0.0
00045 !
00046 ! !INPUT/OUTPUT PARAMETERS:
00047 !
00048       Integer, Intent (InOut)         :: found (len)
00049 
00050 !     Finest level number on which a grid cell was found for point I.
00051 !     Level number < -1: Point was not found and
00052 !                        and last level number was (-found(i,j,k))
00053 !     Level number = (-nlev+1): Never found (input value)
00054 
00055       Integer, Intent (InOut)         :: locations (len)
00056 
00057 !     Indices of the grid cell in which the point was found.
00058 !     Assumed input value locations (:) = 0
00059 !
00060 ! !OUTPUT PARAMETERS:
00061 !
00062       integer, Intent (Out)               :: ierror
00063 
00064 !     Returns the error code of PSMILe_Search_donor_1d_dble;
00065 !             ierror = 0 : No error
00066 !             ierror > 0 : Severe error
00067 !
00068 ! !LOCAL VARIABLES
00069 !
00070 !     ... for locations searched
00071 !
00072       Integer                      :: ibeg, iend
00073 !     Integer                      :: shape   (2, ndim_3d)
00074 !     Integer                      :: range   (2, ndim_3d)
00075 !     Integer                      :: control (2, ndim_3d)
00076 !
00077 !     ... for levels
00078 !
00079       Integer                      :: lev, levc, nlev
00080       Integer                      :: ijkinc (1), ijkcoa (1)
00081       Type(Grid), Pointer          :: grid_info
00082 !
00083 !     ... for error parameters
00084 !
00085 !     Integer, parameter           :: nerrp = 1
00086 !     Integer                      :: ierrp (nerrp)
00087 !
00088 ! !DESCRIPTION:
00089 !
00090 ! Subroutine "PSMILe_Search_donor_1d_dble" searches the donor
00091 ! for the subgrid sent by the sending process.
00092 !
00093 !
00094 ! !REVISION HISTORY:
00095 !
00096 !   Date      Programmer   Description
00097 ! ----------  ----------   -----------
00098 ! 03.07.21    H. Ritzdorf  created
00099 !
00100 !EOP
00101 !----------------------------------------------------------------------
00102 !
00103 ! $Id: psmile_search_donor_1d_dble.F90 2820 2010-12-10 09:20:01Z hanke $
00104 ! $Author: hanke $
00105 !
00106    Character(len=len_cvs_string), save :: mycvs = 
00107        '$Id: psmile_search_donor_1d_dble.F90 2820 2010-12-10 09:20:01Z hanke $'
00108 !
00109 !----------------------------------------------------------------------
00110 !
00111 !  Initialization
00112 !
00113 #ifdef VERBOSE
00114       print 9990, trim(ch_id), Grids(grid_id)%comp_id
00115 
00116       call psmile_flushstd
00117 #endif /* VERBOSE */
00118 !
00119       ierror = 0
00120       grid_info => Grids(grid_id)
00121 !
00122 !-----------------------------------------------------------------------
00123 !     Define shapes and ranges of coordinates to be searched
00124 !-----------------------------------------------------------------------
00125 !
00126 #ifdef NEW
00127 das wird gebraucht, wenn das source und das destination gitter auf dem
00128 gleichen prozess liegen und die search daten nicht gesendet werden,
00129 weil man dann nicht einfach ueber die coords rueberloopen kann.
00130 Trotzdem sollte man hier den Spezialfall haben,
00131 dass die search%daten zusammenliegen !
00132 Also zusaetzliche neue Routinen, die das machen und den NEW code 
00133 verwenden.
00134 !
00135       if (search%grid_type == PRISM_Reglonlatvrt .or. &
00136          (search%grid_type == PRISM_Irrlonlat_Regvrt .and idim == ndim_3d)) then
00137 !        ... Dimension in 1st, 2nd and 3rd direction is only 1-dimensional
00138 
00139          shape = 1
00140          range = 1
00141 
00142          shape (:, idim) = search%shape(:, idim, ipart)
00143          range (:, idim) = search%range(:, idim, ipart)
00144 
00145       else if (search%grid_type == PRISM_Gaussreduced_Regvrt) then
00146 !        ... Dimension of aux grid in 1st, 2nd and 3rd direction is only 1-dimensional
00147 
00148          shape = 1
00149          range = 1
00150 
00151          shape (:, idim) = search%shape(:, idim, ipart)
00152          range (:, idim) = search%range(:, idim, ipart)
00153 
00154       else if (search%grid_type == PRISM_Irrlonlat_Regvrt) then
00155 !        ... Dimension in 1st and 2nd direction is 2-dimensional
00156 !        ... Dimension in 3rd         direction is 1-dimensional
00157 
00158          shape (:, 1:ndim_2d) = search%shape(:, 1:ndim_2d, ipart)
00159          range (:, 1:ndim_2d) = search%range(:, 1:ndim_2d, ipart)
00160 
00161          shape (:, ndim_3d) = 1
00162          range (:, ndim_3d) = 1
00163 
00164       else
00165 !        ... Dimension in all direction is 3-dimensional
00166 
00167          shape (:, :) = search%shape(:, :, ipart)
00168          range (:, :) = search%range(:, :, ipart)
00169 
00170       endif
00171 !
00172       control = range
00173 #endif
00174 !
00175 !-----------------------------------------------------------------------
00176 !     Coarsest level nlev
00177 !-----------------------------------------------------------------------
00178 !
00179 !   ??? Hier fragt sich, ob nicht das suchen in einer einfachen Liste Schneller
00180 !   ??? ist. Das muss man mal alles genau austesten.
00181 !
00182       nlev = grid_info%nlev
00183       lev = nlev
00184 !
00185       ibeg = 1
00186       iend = len
00187 !
00188 #ifdef PRISM_ASSERTION
00189       if ( nlev <= 0 ) then
00190          call psmile_assert (__FILE__, __LINE__, "nlev <= 0")
00191       endif
00192 #endif
00193 !
00194       call psmile_mg_coarse_1d_dble (lev, &
00195                      grid_info%mg_infos(lev)%double_arrays%chmin(idim)%vector, &
00196                      grid_info%mg_infos(lev)%double_arrays%chmax(idim)%vector, &
00197                      found, locations, coords, &
00198                      ibeg, iend)
00199 !
00200       if (ibeg > iend) then
00201 #ifdef VERBOSE
00202          print 9980, trim(ch_id), grid_info%comp_id, ierror
00203 
00204          call psmile_flushstd
00205 #endif /* VERBOSE */
00206          return
00207       endif
00208 !
00209 !===> Multiple grids
00210 !
00211 !     ijkinc(1) = max (1, (iend-ibeg)/2)
00212       ijkinc(1) = 1
00213 !
00214          do lev = nlev-1, 1, -1
00215 !
00216          levc = lev + 1
00217 !
00218          ijkcoa (1) = 2
00219 !
00220          if (grid_info%mg_infos(lev)%levdim(idim) == &
00221              grid_info%mg_infos(levc)%levdim(idim)) then
00222 !           skippen ?!?!
00223             ijkcoa (1) = 1
00224          endif
00225 !
00226          call psmile_mg_next_level_1d_dble ( grid_id, idim, lev, nlev,    &
00227                 grid_info%mg_infos(lev)%double_arrays%chmin(idim)%vector, &
00228                 grid_info%mg_infos(lev)%double_arrays%chmax(idim)%vector, &
00229                 grid_info%mg_infos(lev)%double_arrays%midp(idim)%vector,  &
00230                 grid_info%mg_infos(lev)%levdim(idim),                   &
00231                 found,  locations, coords, ibeg, iend,                  &
00232                         ijkinc(1), ijkcoa(1), ierror)
00233          if (ierror > 0) return
00234 !
00235 !        ijkinc (1) = max (1, ijkinc(1)/ijkcoa(1))
00236          enddo ! lev
00237 !
00238 !   Compute locations
00239 !
00240       locations (ibeg:iend) = locations (ibeg:iend) + grid_info%ijk0 (idim)
00241 !
00242 !===> All done
00243 !
00244 #ifdef VERBOSE
00245       print 9980, trim(ch_id), grid_info%comp_id, ierror
00246 
00247       call psmile_flushstd
00248 #endif /* VERBOSE */
00249 !
00250 !  Formats:
00251 !
00252 9990 format (1x, a, ': psmile_search_donor_1d_dble: comp_id =', i3)
00253 9980 format (1x, a, ': psmile_search_donor_1d_dble: comp_id =', i3, &
00254                     '; eof ierror =', i4)
00255 
00256       end subroutine PSMILe_Search_donor_1d_dble

Generated on 18 Mar 2011 for Oasis4 by  doxygen 1.6.1