psmile_mg_control_cell_2d_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_MG_Control_cell_2d_real
00008 !
00009 ! !INTERFACE:
00010 
00011       subroutine psmile_mg_control_cell_2d_real (          &
00012                       chmin1, chmin2,                      &
00013                       chmax1, chmax2,                      &
00014                       midp1,  midp2,           levdim,     &
00015                       ijk, xyz, nold, ignore, wide, found, newijk)
00016 !
00017 ! !USES:
00018 !
00019       use PRISM_constants
00020 !
00021       use PSMILe, dummy_interface => PSMILe_mg_control_cell_2d_real
00022 
00023       implicit none
00024 !
00025 ! !INPUT PARAMETERS:
00026 !
00027       Integer, Intent (In)                        :: levdim (ndim_2d)
00028 
00029 !     Dimensions of chmin, chmax and midp
00030 
00031       Real, Intent (In)                           :: chmin1 (0:levdim(1), 
00032                                                              0:levdim(2))
00033       Real, Intent (In)                           :: chmin2 (0:levdim(1), 
00034                                                              0:levdim(2))
00035 
00036 !     Minimum of grid coordinates per grid cell
00037 
00038       Real, Intent (In)                           :: chmax1 (0:levdim(1), 
00039                                                              0:levdim(2))
00040       Real, Intent (In)                           :: chmax2 (0:levdim(1), 
00041                                                              0:levdim(2))
00042 
00043 !     Maximum of grid coordinates per grid cell
00044 
00045       Real, Intent (In)                           ::  midp1 (0:levdim(1), 
00046                                                              0:levdim(2))
00047       Real, Intent (In)                           ::  midp2 (0:levdim(1), 
00048                                                              0:levdim(2))
00049 
00050 !     Mid point of the cell
00051 
00052       Integer, Intent (In)                        :: ijk (ndim_2d)
00053 
00054 !     Index of the cell for which the cell parts should be
00055 !     controlled.
00056 
00057       Real, Intent (In)                           :: xyz (ndim_2d)
00058 
00059 !     Corrdinates of the point for which is searched
00060 
00061       Integer, Intent (In)                        :: nold
00062 
00063 !     Old number of the cell part where the point XYZ was found last.
00064 !     And all other parts "i" with 
00065 !                DIST(i) < DIST(NOLD) were already controlled.
00066 !     NOLD = 0: Start first search in this cell.
00067 !
00068       Integer, Intent (In)                        :: ignore (ndim_2d)
00069 !
00070 !     Cell indices to be ignored.
00071 !
00072       logical, Intent (In)                        :: wide
00073 
00074 !     Wide control ?
00075 !
00076 ! !OUTPUT PARAMETERS:
00077 !
00078       integer, Intent (Out)                       :: found
00079 
00080 !     Was the point xyz found ?
00081 !     found = 0: Not found
00082 !     found > 0: Part of the cell where point XYZ was found.
00083 
00084       integer, Intent (Out)                       :: newijk (ndim_2d)
00085 
00086 !     Indices on level LEV if found > 0 is returned.
00087 !
00088 ! !DEFINED PARAMETERS:
00089 !
00090       Real, Parameter                             :: remax = 1.0e20
00091 !
00092       Integer, Parameter                          :: ndtrys = 4
00093       Integer, Parameter                          :: ndtryw = 16
00094 !
00095 ! !LOCAL VARIABLES
00096 !
00097       Real                                        :: dist (ndtryw)
00098       Integer                                     :: i, n, ndtry
00099       Integer                                     :: nmin (1)
00100 !
00101       Integer, Pointer                            :: ijkctl (:, :)
00102       Integer, Target                             :: ijkctls (ndim_2d, ndtrys)
00103       Integer, Target                             :: ijkctlw (ndim_2d, ndtryw)
00104       Integer                                     :: ijkdst (ndim_2d, ndtryw)
00105       Logical                                     :: within (ndtryw)
00106 !
00107 #ifdef DEBUGX
00108       logical control
00109       common /huhu1/ control
00110 #endif
00111 !
00112 ! !DESCRIPTION:
00113 !
00114 ! Subroutine "PSMILe_mg_control_cell_2d_real" searches the point XYZ
00115 ! in the n the neihbourhood of a cell. The cell which is controlled
00116 ! starts in IJK. All cells with (i) in
00117 !      (ijk(:)-1:ijk(:)+2) if wide is true.
00118 !      (ijk(:)-1:ijk(:)+1) otherwise
00119 ! are controlled.
00120 !
00121 ! ?!? Vielleicht koennte man alle Punkte zurueckgeben, die man
00122 !     gefunden hat.
00123 !
00124 !
00125 ! !REVISION HISTORY:
00126 !
00127 !   Date      Programmer   Description
00128 ! ----------  ----------   -----------
00129 ! 03.07.21    H. Ritzdorf  created
00130 !
00131 !EOP
00132 !----------------------------------------------------------------------
00133 !
00134 !  $Id: psmile_mg_control_cell_2d_real.F90 2325 2010-04-21 15:00:07Z valcke $
00135 !  $Author: valcke $
00136 !
00137    Character(len=len_cvs_string), save :: mycvs = 
00138        '$Id: psmile_mg_control_cell_2d_real.F90 2325 2010-04-21 15:00:07Z valcke $'
00139 !
00140 !----------------------------------------------------------------------
00141 !
00142 !  Relative control indices for sub-cells of coarse grid cells
00143 !
00144       data ((ijkctls (i, n), i=1,ndim_2d), n = 1,ndtrys) &
00145                              / 0, 0,   1, 0, &
00146                                0, 1,   1, 1/
00147 !
00148 !  Relative control indices for sub-cells of coarse grid cells
00149 !
00150 !  27 (-1:1, -1:1, -1:1) Zellen sollte auch reichen
00151 !  Die reihenfolge der ersten 8 muss so bleiben !
00152 !
00153       data ((ijkctlw (i, n), i=1,ndim_2d), n = 1,ndtryw) &
00154                              / 0, 0,   1, 0,   0, 1,   1, 1, &
00155                               -1,-1,   0,-1,   1,-1,   2,-1, &
00156                               -1, 0,                   2, 0, &
00157                               -1, 1,                   2, 1, &
00158                               -1, 2,   0, 2,   1, 2,   2, 2/
00159 !
00160 !----------------------------------------------------------------------
00161 !
00162 !  Initialization
00163 !
00164 #ifdef VERBOSE
00165       print 9990, trim(ch_id)
00166 
00167       call psmile_flushstd
00168 #endif /* VERBOSE */
00169 !
00170 !===> Compute indices to be controlled
00171 !
00172       if (wide) then
00173          ndtry = ndtryw
00174          ijkctl => ijkctlw
00175       else
00176          ndtry = ndtrys
00177          ijkctl => ijkctls
00178       endif
00179 !
00180 #ifdef PRISM_ASSERTION
00181       if (nold < 0 .or. nold > ndtry) then
00182          call psmile_assert (__FILE__, __LINE__, &
00183                              'wrong nold')
00184       endif
00185 !
00186       if (minval(ijk(:)) < 0 .or.  minval(levdim(:)-ijk(:)) < 0) then
00187          call psmile_assert (__FILE__, __LINE__, 'wrong ijk')
00188       endif
00189 #endif
00190 !
00191 !===> Compute indices to be controlled
00192 !
00193         do n = 1, ndtry
00194         ijkdst (:, n) = ijk(:) + ijkctl(:,n)
00195         end do
00196 !
00197 !===> Control index range
00198 !
00199         do n = 1, ndtry
00200         within (n) = ijkdst (1,n) >= 0 .and. ijkdst(1,n) <= levdim(1) .and. &
00201                      ijkdst (2,n) >= 0 .and. ijkdst(2,n) <= levdim(2)
00202         end do
00203 #ifdef DEBUGX
00204 if (control) print *, 'within a:', within (1:ndtry)
00205 #endif
00206 !
00207 !===> Control location of point XYZ in the boxes
00208 !
00209 !  ??? ist (sent-chmin)*(sent-chmax) <= 0.0 schneller ???
00210 !
00211 !cdir vector
00212         do n = 1, ndtry
00213         if (within(n)) then
00214            within (n) = xyz(1) >= chmin1(ijkdst(1,n), ijkdst(2,n)) .and. &
00215                         xyz(2) >= chmin2(ijkdst(1,n), ijkdst(2,n)) .and. &
00216                         xyz(1) <= chmax1(ijkdst(1,n), ijkdst(2,n)) .and. &
00217                         xyz(2) <= chmax2(ijkdst(1,n), ijkdst(2,n))
00218         endif
00219         end do
00220 #ifdef DEBUGX
00221 if (control) print *, 'within b:', within (1:ndtry)
00222 #endif
00223 !
00224 !===> Compute distances to the cell midpoints
00225 !
00226 !cdir vector
00227         do n = 1, ndtry
00228         if (within(n)) then
00229             dist (n) = (xyz(1) - midp1(ijkdst(1,n), ijkdst(2,n)))**2 &
00230                      + (xyz(2) - midp2(ijkdst(1,n), ijkdst(2,n)))**2
00231 
00232         else
00233 
00234             dist (n) = remax
00235 
00236         endif
00237         end do
00238 !
00239 !===> Set old values to remax
00240 !
00241       if (nold > 0) then
00242 
00243 #ifdef PRISM_ASSERTIONX
00244          if (dist(nold) .eq. remax) then
00245             write (*, 9970) xyz
00246 9970  format (1x, 1p, 3e18.9)
00247 !
00248             call psmile_assert (__FILE__, __LINE__, &
00249                                 'incorrect nold')
00250          endif
00251 #endif
00252 !
00253 !===> Exclude all old distances
00254 !
00255             do n = 1, nold-1
00256             if (dist(n) <= dist(nold)) dist (n) = remax
00257             end do
00258 !
00259             do n = nold+1, ndtry
00260             if (dist(n) < dist(nold)) dist (n) = remax
00261             end do
00262 !
00263          dist (nold) = remax
00264       endif
00265 !
00266 !===> Look for the minimum distance
00267 !
00268 100   continue
00269 !
00270       nmin = MINLOC (dist(1:ndtry))
00271 #ifdef MINLOCFIX
00272       if (nmin(1).eq.0) nmin=1
00273 #endif /* MINLOCFIX */
00274 !
00275 !===> Was a valid cell found ?
00276 !
00277       if (dist(nmin(1)) == remax) then
00278          found = 0
00279       else
00280 !
00281          if (ijkdst(1, nmin(1)) == ignore (1) .and. &
00282              ijkdst(2, nmin(1)) == ignore (2)) then
00283 !           Cell has to be  ignored
00284 
00285             dist (nmin(1)) = remax
00286             go to 100
00287          endif
00288 !
00289          found = nmin (1)
00290 !
00291          newijk (:) = ijkdst (:, nmin(1))
00292 
00293 #ifdef PRISM_ASSERTION
00294          if ( xyz(1) < chmin1(newijk(1), newijk(2)) .or. &
00295               xyz(1) > chmax1(newijk(1), newijk(2)) .or. &
00296               xyz(2) < chmin2(newijk(1), newijk(2)) .or. &
00297               xyz(2) > chmax2(newijk(1), newijk(2))) then
00298             print *, 'xyz', xyz
00299             print *, 'chmin', chmin1(newijk(1), newijk(2)), &
00300                               chmin2(newijk(1), newijk(2))
00301             print *, 'chmax', chmax1(newijk(1), newijk(2)), &
00302                               chmax2(newijk(1), newijk(2))
00303 
00304             call psmile_assert (__FILE__, __LINE__, &
00305                       'not contained in computed cell')
00306         endif
00307 #endif
00308       endif
00309 !
00310 !===> All done
00311 !
00312 #ifdef VERBOSE
00313       print 9980, trim(ch_id)
00314 
00315       call psmile_flushstd
00316 #endif /* VERBOSE */
00317 !
00318       return
00319 !
00320 !  Formats:
00321 !
00322 9990 format (1x, a, ': PSMILe_mg_control_cell_2d_real:')
00323 9980 format (1x, a, ': PSMILe_mg_control_cell_2d_real: eof')
00324 
00325       end subroutine PSMILe_mg_control_cell_2d_real

Generated on 18 Mar 2011 for Oasis4 by  doxygen 1.6.1