psmile_mg_coarse_1d_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_coarse_1d_real
00008 !
00009 ! !INTERFACE:
00010 
00011       subroutine psmile_mg_coarse_1d_real (lev, chmin, chmax,        &
00012                                            found, locations, coords, &
00013                                            ibeg, iend)
00014 !
00015 ! !USES:
00016 !
00017       use PRISM_constants
00018 !
00019       use PSMILe, dummy_interface => PSMILe_mg_coarse_1d_real
00020 
00021       implicit none
00022 !
00023 ! !INPUT/OUTPUT PARAMETERS:
00024 !
00025       Integer, Intent (InOut)         :: ibeg
00026 
00027 !     Input  value: First index in "coords" to be searched
00028 !     Output value: First index in "coords" for which a cell was found
00029 
00030       Integer, Intent (InOut)         :: iend
00031 
00032 !     Input  value: Last  index in "coords" to be searched
00033 !                   Dimension of "found" and "location".
00034 !     Output value: Last  index in "coords" for which a cell was found
00035 
00036       Integer, Intent (InOut)         :: found (iend)
00037 
00038 !     Finest level number on which a grid cell was found for point I.
00039 !     Level number = -(nlev+1): Never found (input value)
00040 
00041       Integer, Intent (InOut)         :: locations (iend)
00042 
00043 !     Indices of the grid cell in which the point was found.
00044 !     Assumed input value locations (:, :) = 0
00045 !
00046 !
00047 ! !INPUT PARAMETERS:
00048 !
00049       Integer, Intent (In)            :: lev
00050 
00051 !     Level number of the coarsest grid
00052 
00053       Real, Intent (In)               :: chmin (0:0)
00054 
00055 !     Minimum of grid coordinates per grid cell
00056 
00057       Real, Intent (In)               :: chmax (0:0)
00058 
00059 !     Maximum of grid coordinates per grid cell
00060 
00061       Real, Intent (In)               :: coords (iend)
00062 !
00063 !     Coordinates to be searched
00064 !
00065 ! !LOCAL VARIABLES
00066 !
00067 !     ... for locations searched
00068 !
00069       Integer                      :: i, i2
00070 #ifdef VERBOSE
00071       Integer                      :: n
00072 #endif
00073 !
00074 ! !DESCRIPTION:
00075 !
00076 ! Subroutine "PSMILe_mg_coarse_1d_real" searches the donor cells
00077 ! on the coarsest MG grid (only 1 cell) for the subgrid sent by the
00078 ! sending process.
00079 !
00080 !
00081 ! !REVISION HISTORY:
00082 !
00083 !   Date      Programmer   Description
00084 ! ----------  ----------   -----------
00085 ! 03.07.21    H. Ritzdorf  created
00086 !
00087 !EOP
00088 !----------------------------------------------------------------------
00089 !
00090 !  $Id: psmile_mg_coarse_1d_real.F90 2325 2010-04-21 15:00:07Z valcke $
00091 !  $Author: valcke $
00092 !
00093    Character(len=len_cvs_string), save :: mycvs = 
00094        '$Id: psmile_mg_coarse_1d_real.F90 2325 2010-04-21 15:00:07Z valcke $'
00095 !
00096 !----------------------------------------------------------------------
00097 !
00098 !  Initialization
00099 !
00100 #ifdef VERBOSE
00101       print 9990, trim(ch_id), lev, ibeg, iend
00102 
00103       call psmile_flushstd
00104 #endif /* VERBOSE */
00105 !
00106 #ifdef PRISM_ASSERTION
00107       if (chmin(0) > chmax (0)) then
00108          call psmile_assert ( __FILE__, __LINE__, &
00109                               'incorrect coarsest box')
00110       endif
00111 #endif
00112 !
00113 !===> Is point I in the coarsest box ?
00114 !
00115 !cdir vector
00116          do i = ibeg, iend
00117          if (chmin(0) <= coords(i) .and. coords(i) <= chmax(0)) then
00118             found (i) = lev
00119          endif
00120          end do
00121 !
00122 !===> Update ibeg, iend
00123 !
00124 !cdir vector
00125          do i = ibeg, iend
00126          if (found (i) == lev) exit
00127          enddo
00128 !
00129       if (i > iend) then
00130          iend = ibeg - 1
00131 !
00132       else
00133          ibeg = i
00134 !
00135          if (found(iend) /= lev) then
00136 !cdir vector
00137                do i = ibeg, iend
00138                if (found(i) == lev) then
00139                   i2 = i
00140                endif
00141                enddo
00142 !
00143             iend = i2
00144 
00145          endif
00146       endif
00147 !
00148 !===> All done
00149 !
00150 #ifdef VERBOSE
00151 !
00152 !===> Get number of points found
00153 !
00154 #ifdef DONT_USE_COUNT
00155       n = 0
00156 !cdir vector
00157          do i = ibeg, iend
00158          if (found(i) == lev) n = n + 1
00159          done
00160 #else
00161       n = count (found(ibeg:iend) == lev)
00162 #endif
00163 !
00164       print 9980, trim(ch_id), lev, n
00165 
00166       call psmile_flushstd
00167 #endif /* VERBOSE */
00168 !
00169 !  Formats:
00170 !
00171 9990 format (1x, a, ': PSMILe_mg_coarse_1d_real: level', i3, &
00172                     ', ibeg, iend', 2i6)
00173 9980 format (1x, a, ': PSMILe_mg_coarse_1d_real: eof, level', i3, &
00174                     ', found =', i9)
00175 
00176       end subroutine PSMILe_mg_coarse_1d_real

Generated on 18 Mar 2011 for Oasis4 by  doxygen 1.6.1