psmile_mg_first_subgrid_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_First_subgrid_1d_dble
00008 !
00009 ! !INTERFACE:
00010 
00011       subroutine psmile_mg_first_subgrid_1d_dble (       &
00012                     array, idlow, idhigh, nbr_corners,   &
00013                            range,                        &
00014                     chmin, chmax, midp,                  &
00015                     levdim, ierror)
00016 !
00017 ! !USES:
00018 !
00019       use PRISM_constants
00020 !
00021       use PSMILe, dummy_interface => PSMILe_MG_First_subgrid_1d_dble
00022 
00023       implicit none
00024 !
00025 ! !INPUT PARAMETERS:
00026 !
00027       Integer, Intent (In)             :: idlow, idhigh, nbr_corners
00028 
00029 !     Dimensions of "array"
00030 !
00031       Double Precision, Intent (In)    :: array (idlow:idhigh, nbr_corners)
00032 
00033 !     Fully dimensioned Array for which the first level
00034 !     should be computed.
00035 
00036       Integer, Intent (In)             :: range (2)
00037 
00038 !     Definintion of the subgrid
00039 
00040       Integer, Intent (In)             :: levdim
00041 
00042 !     Dimension of the first MG level
00043 !
00044 ! !OUTPUT PARAMETERS:
00045 !
00046       Double Precision, Intent (Out)   :: chmin (0:levdim)
00047 !
00048 !     Minimum extension of bounding box of each coordinate cell.
00049 !
00050       Double Precision, Intent (Out)   :: chmax (0:levdim)
00051 !
00052 !     Maximum extension of bounding box of each coordinate cell.
00053 !
00054       Double Precision, Intent (Out)   :: midp  (0:levdim)
00055 !
00056 !     Mid point of bounding box of each coordinate cell.
00057 !
00058       Integer, Intent (Out)            :: ierror
00059 
00060 !     Returns the error code of PSMILe_MG_First_subgrid_1d_dble
00061 !             ierror = 0 : No error
00062 !             ierror > 0 : Severe error
00063 !
00064 ! !DEFINED PARAMETERS:
00065 !
00066       Double Precision, Parameter    :: epsp1 = 1.0d0 + 5d-6
00067 !
00068 ! !LOCAL VARIABLES
00069 !
00070       integer                        :: i, ibeg, iend
00071 !
00072 ! !DESCRIPTION:
00073 !
00074 ! Subroutine "PSMILe_MG_First_subgrid_1d_dble" computes the bounding
00075 ! boxes and the mid point for the corner coordinates ("control volumes")
00076 ! for a 1-dimensional subgrid. The bounding boxes and mid points are
00077 ! used in the Multigrid search.
00078 !
00079 ! !REVISION HISTORY:
00080 !
00081 !   Date      Programmer   Description
00082 ! ----------  ----------   -----------
00083 ! 03.06.25    H. Ritzdorf  created
00084 !
00085 !EOP
00086 !----------------------------------------------------------------------
00087 !
00088 !  $Id: psmile_mg_first_subgrid_1d_dble.F90 2959 2011-02-16 13:28:28Z hanke $
00089 !  $Author: hanke $
00090 !
00091    Character(len=len_cvs_string), save :: mycvs = 
00092        '$Id: psmile_mg_first_subgrid_1d_dble.F90 2959 2011-02-16 13:28:28Z hanke $'
00093 !
00094 !----------------------------------------------------------------------
00095 !
00096 #ifdef VERBOSE
00097       print 9990, trim(ch_id)
00098 
00099       call psmile_flushstd
00100 #endif /* VERBOSE */
00101 !
00102 !  Initialization
00103 !
00104       ierror = 0
00105 
00106 !
00107       ibeg = range (1)
00108       iend = min (range(1)+levdim, range (2))
00109 !
00110 !     ... Compute minimum extension of grid cell
00111 !
00112 !cdir vector
00113       do i = ibeg, iend
00114          chmin (i-range(1)) = min (array (i, 1), array (i, 2))
00115       enddo
00116 !
00117 !     ... Compute maximum extension of grid cell
00118 !
00119 !cdir vector
00120       do i = ibeg, iend
00121          chmax (i-range(1)) = max (array (i, 1), array (i, 2))
00122       enddo
00123 !
00124 !     TODO: Attaching grid cells which interfere at the cyclic 
00125 !           boundary should be stored and searched in order to
00126 !           improve the MG search.
00127 !           Type (Cyclic_cell_information) cyclic_info
00128 !           wieder hochholen !
00129 !
00130 !     ... Compute mid point of grid cell
00131 !
00132 !cdir vector
00133       do i = ibeg, iend
00134          midp (i-range(1))  = (array (i, 1) + array (i, 2)) * 0.5
00135       enddo
00136 !
00137 !     ... Create dummy cells if necessary
00138 !
00139       if (iend < range(1)+levdim) then
00140          ibeg = ibeg - range(1)
00141          iend = iend - range(1)
00142 
00143          chmin (iend+1:levdim) = chmax (iend) * epsp1
00144          chmax (iend+1:levdim) = chmax (iend)
00145          midp  (iend+1:levdim) = midp  (iend)
00146       endif
00147 !
00148 !===> All done
00149 !
00150 #ifdef DEBUGX
00151       do i = range(1), range(2)
00152          print *, 'i, array', i, array(i,1), array(i,2)
00153       end do
00154 !
00155       do i = 0, levdim
00156          print *, 'i, chmin, chmax', i, chmin(i), chmax(i)
00157       end do
00158 #endif
00159 
00160 #ifdef VERBOSE
00161 !     print *, 'chmin', chmin(0), chmin (levdim)
00162 !     print *, 'chmax', chmax(0), chmax (levdim)
00163 !     print *, 'midp ', midp (0), midp  (levdim)
00164 
00165       print 9980, trim(ch_id), ierror
00166 
00167       call psmile_flushstd
00168 #endif /* VERBOSE */
00169 !
00170 !  Formats:
00171 !
00172 9990  format (1x, a, ': psmile_mg_first_subgrid_1d_dble')
00173 9980  format (1x, a, ': psmile_mg_first_subgrid_1d_dble eof: ierror =', i4)
00174 !
00175       end subroutine PSMILe_MG_First_subgrid_1d_dble

Generated on 18 Mar 2011 for Oasis4 by  doxygen 1.6.1