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

Generated on 18 Mar 2011 for Oasis4 by  doxygen 1.6.1