psmile_mg_prev_levels_1d_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_Prev_levels_1d_dble
00008 !
00009 ! !INTERFACE:
00010 
00011       subroutine psmile_mg_prev_levels_1d_dble (grid_id, idim, lev, nlev, &
00012                            lstijk, xyz, found, newijk)
00013 !
00014 ! !USES:
00015 !
00016       use PRISM_constants
00017 !
00018       use PSMILe, dummy_interface => PSMILe_mg_prev_levels_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)            :: lev
00034 
00035 !     Level number of the next grid
00036 
00037       Integer, Intent (In)            :: nlev
00038 
00039 !     Total number of levels
00040 
00041       Integer, Intent (In)            :: lstijk
00042 
00043 !     Last indices of cell on coarse level LEV+1
00044 
00045       Double Precision, Intent (In)   :: xyz
00046 
00047 !     Coordinates of the point which is searched
00048 !
00049 ! !OUTPUT PARAMETERS:
00050 !
00051       integer, Intent (Out)           :: found
00052 
00053 !     Was the point xyz found ?
00054 !     found = 0: Not found
00055 !     found > 0: Part of the cell where point XYZ was found.
00056 !
00057       integer, Intent (Out)           :: newijk
00058 !
00059 !     Indices on level LEV if found > 0 is returned.
00060 !
00061 ! !DEFINED PARAMETERS:
00062 !
00063       Logical, parameter              :: wide = .true.
00064       Logical, parameter              :: small = .false.
00065 !
00066 ! !LOCAL VARIABLES
00067 !
00068 !     ... for levels
00069 
00070       Integer                         :: levc, levd, levdbg, levu, levubg
00071       Integer                         :: nold (nlev), ijk (nlev)
00072       Logical                         :: all (nlev)
00073       Type (Grid), Pointer            :: grid_info
00074 !
00075 ! !DESCRIPTION:
00076 !
00077 ! Subroutine "PSMILe_mg_prev_levels_1d_dble" searches for the donor cell
00078 ! of point "sent".  The point was last found on level LEVC in donor cell (IC).
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_prev_levels_1d_dble.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_prev_levels_1d_dble.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
00102 
00103       call psmile_flushstd
00104 #endif /* VERBOSE */
00105 !
00106       grid_info => Grids(grid_id)
00107 
00108       levc = lev + 1
00109 
00110 #ifdef PRISM_ASSERTION
00111       if (levc >= nlev) then
00112          call psmile_assert (__FILE__, __LINE__, 'incorrect level')
00113       endif
00114 !
00115       if (lstijk < 0 .or. &
00116           lstijk > grid_info%mg_infos(levc)%levdim(idim)) then
00117          print *, lstijk, grid_info%mg_infos(levc)%levdim(idim)
00118          call psmile_assert (__FILE__, __LINE__, 'wrong lstijk')
00119       endif
00120 #endif
00121 !
00122       ijk (levc) = lstijk
00123 !
00124          do levd = levc, nlev-1
00125          all  (levd) = .false.
00126          nold (levd) = 1
00127          end do
00128 !
00129 !===> ... Compute indices on the coarser levels
00130 !
00131          do levd = levc+1, nlev-1
00132          if (grid_info%mg_infos(levd)%levdim(idim) /= &
00133              grid_info%mg_infos(levd-1)%levdim(idim) ) then
00134             ijk (levd) = ijk (levd-1) / 2
00135          else
00136             ijk (levd) = ijk (levd-1)
00137          endif
00138          end do
00139 !
00140 !---------------------------------------------------------------------
00141 !     Go down in the level hierarchie
00142 !---------------------------------------------------------------------
00143 !
00144       levdbg = levc
00145 !
00146 100   continue
00147 !
00148          do levd = levdbg, nlev-1
00149 !
00150 #ifdef PRISM_ASSERTION
00151          if (ijk(levd) < 0 .or. &
00152              ijk(levd) > grid_info%mg_infos(levd)%levdim(idim)) then
00153             print *, levd, ijk(levd), grid_info%mg_infos(levd)%levdim(idim)
00154             print *, lstijk, levc
00155             call psmile_assert (__FILE__, __LINE__, 'wrong ijk(levd)')
00156          endif
00157 #endif /* PRISM_ASSERTION */
00158 !
00159          call psmile_mg_control_cell_1d_dble ( &
00160              grid_info%mg_infos(levd)%double_arrays%chmin(idim)%vector, &
00161              grid_info%mg_infos(levd)%double_arrays%chmax(idim)%vector, &
00162              grid_info%mg_infos(levd)%double_arrays%midp (idim)%vector, &
00163              grid_info%mg_infos(levd)%levdim(idim), &
00164              ijk(levd), xyz, nold(levd), all(levd), wide, found, newijk)
00165 !
00166          if (found > 0) go to 1990
00167 !
00168          end do
00169 !
00170 !===> Not found
00171 !
00172 
00173 #ifdef VERBOSE
00174       print 9980, trim(ch_id), lev
00175 
00176       call psmile_flushstd
00177 #endif /* VERBOSE */
00178 
00179       return
00180 !
00181 !---------------------------------------------------------------------
00182 !     Up again; startpoint is newijk on level levd
00183 !---------------------------------------------------------------------
00184 !
00185 1990  nold (levd) = found
00186       all (levd) = .true.
00187 !
00188       levdbg = levd
00189       levubg = levd - 1
00190 !
00191       nold (levubg) = 0
00192       all  (levubg) = .true.
00193 !
00194       if (grid_info%mg_infos(levubg)%levdim(idim) > &
00195           grid_info%mg_infos(levubg+1)%levdim(idim)) then
00196          ijk (levubg) = min (newijk*2, grid_info%mg_infos(levubg)%levdim(idim))
00197       else
00198          ijk (levubg) = newijk
00199       endif
00200 !
00201 1995  continue
00202 !
00203          do levu = levubg, lev, -1
00204          if (levu .ne. levubg) then
00205             if (grid_info%mg_infos(levu)%levdim(idim) > &
00206                 grid_info%mg_infos(levu+1)%levdim(idim)) then
00207                ijk (levu) = min (newijk*2, grid_info%mg_infos(levu)%levdim(idim))
00208             else
00209                ijk (levu) = newijk
00210             endif
00211          endif
00212 !
00213 #ifdef PRISM_ASSERTION
00214          if (ijk(levu) < 0 .or. &
00215              ijk(levu) > grid_info%mg_infos(levu)%levdim(idim)) then
00216             print *, levu, ijk(levu), grid_info%mg_infos(levu)%levdim(idim)
00217             call psmile_assert (__FILE__, __LINE__, 'wrong ijk(levu)')
00218          endif
00219 #endif /* PRISM_ASSERTION */
00220 !
00221          call psmile_mg_control_cell_1d_dble ( &
00222              grid_info%mg_infos(levu)%double_arrays%chmin(idim)%vector, &
00223              grid_info%mg_infos(levu)%double_arrays%chmax(idim)%vector, &
00224              grid_info%mg_infos(levu)%double_arrays%midp (idim)%vector, &
00225              grid_info%mg_infos(levu)%levdim(idim),                   &
00226              ijk(levu), xyz, nold(levu), all(levu), small, found, newijk)
00227 !
00228          if (found .eq. 0) then
00229             if (levu .lt. levd-1) then
00230 !
00231 !===> ... go back to the previous level; and try again
00232 !
00233                levubg = levu + 1
00234                go to 1995
00235             else
00236 !
00237 !===> ... go down
00238 !
00239                go to 100
00240             endif
00241 !
00242          else if (levu > lev) then
00243             all (levu) = .true.
00244             nold (levu) = found
00245 !
00246             nold (levu-1) = 0
00247             all (levu-1) = .true.
00248          endif
00249          end do ! levu
00250 !
00251 !===> All done
00252 !
00253 #ifdef VERBOSE
00254       print 9980, trim(ch_id), lev
00255 
00256       call psmile_flushstd
00257 #endif /* VERBOSE */
00258 !
00259 !  Formats:
00260 !
00261 9990 format (1x, a, ': psmile_mg_prev_levels_1d_dble: level', i3)
00262 9980 format (1x, a, ': psmile_mg_prev_levels_1d_dble: eof, level', i3)
00263 
00264       end subroutine PSMILe_mg_prev_levels_1d_dble

Generated on 18 Mar 2011 for Oasis4 by  doxygen 1.6.1