psmile_mg_ctrl_subgrid_3d_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_ctrl_subgrid_3d_real
00008 !
00009 ! !INTERFACE:
00010 
00011       subroutine psmile_mg_ctrl_subgrid_3d_real (array,  &
00012                     corner_shape, nbr_corners, range,    &
00013                     chmin, chmax, levdim,                &
00014                     period, grid_id, ind, ierror)
00015 !
00016 ! !USES:
00017 !
00018       use PRISM_constants
00019 !
00020       use PSMILe, dummy_interface => PSMILe_MG_ctrl_subgrid_3d_real
00021 
00022       Implicit none
00023 !
00024 ! !INPUT PARAMETERS:
00025 !
00026       Integer, Intent (In)             :: corner_shape (2, ndim_3d)
00027       Integer, Intent (In)             :: nbr_corners
00028 
00029 !     Dimensions of "array"
00030 !
00031       Real, Intent (In)                :: array ( 
00032                                            corner_shape(1,1):corner_shape(2,1), 
00033                                            corner_shape(1,2):corner_shape(2,2), 
00034                                            corner_shape(1,3):corner_shape(2,3), 
00035                                            nbr_corners)
00036 
00037 !     Fully dimensioned array "array" for which the first MG level
00038 !     was be computed.
00039 
00040       Integer, Intent (In)             :: range (2, ndim_3d)
00041 
00042 !     Definintion of the subgrid
00043 
00044       Integer, Intent (In)             :: levdim (ndim_3d)
00045 
00046 !     Dimensions of the first MG level
00047 
00048       Real, Intent (In)                :: chmin (0:levdim(1), 0:levdim(2), 
00049                                                  0:levdim(3))
00050 !
00051 !     Minimum value of bounding box of coordinate array "array".
00052 !
00053       Real, Intent (In)                :: chmax (0:levdim(1), 0:levdim(2), 
00054                                                  0:levdim(3))
00055 !
00056 !     Maximum value of bounding box of coordinate array "array".
00057 !
00058       Real, Intent (In)                :: period
00059 !
00060 !     Period of cyclic coordinate
00061 
00062       Integer, Intent (In)             :: grid_id
00063 
00064 !     Grid id
00065 
00066       Integer, Intent (In)             :: ind
00067 
00068 !     Coordinate index of the grid
00069 !
00070 ! !OUTPUT PARAMETERS:
00071 !
00072       Integer, Intent (Out)            :: ierror
00073 
00074 !     Returns the error code of PSMILe_MG_ctrl_subgrid_3d_real
00075 !             ierror = 0 : No error
00076 !             ierror > 0 : Severe error
00077 !
00078 ! !LOCAL VARIABLES
00079 !
00080       Integer                          :: i, j, k, kbeg
00081       Integer                          :: ijkend (ndim_3d)
00082       Real                             :: dist_cyclic
00083 !
00084 ! !DESCRIPTION:
00085 !
00086 ! Subroutine "PSMILe_MG_ctrl_subgrid_3d_real" controls whether
00087 ! the "control volumes" defined by the corner arrays "array" fulfill the
00088 ! internal requirements of the MG search.
00089 !
00090 ! !REVISION HISTORY:
00091 !
00092 !   Date      Programmer   Description
00093 ! ----------  ----------   -----------
00094 ! 03.06.25    H. Ritzdorf  created
00095 !
00096 !EOP
00097 !----------------------------------------------------------------------
00098 !
00099 !  $Id: psmile_mg_ctrl_subgrid_3d_real.F90 2325 2010-04-21 15:00:07Z valcke $
00100 !  $Author: valcke $
00101 !
00102    Character(len=len_cvs_string), save :: mycvs = 
00103        '$Id: psmile_mg_ctrl_subgrid_3d_real.F90 2325 2010-04-21 15:00:07Z valcke $'
00104 !
00105 !----------------------------------------------------------------------
00106 !
00107 #ifdef VERBOSE
00108       print 9990, trim(ch_id)
00109 
00110       call psmile_flushstd
00111 #endif /* VERBOSE */
00112 !
00113 !  Initialization
00114 !
00115       ierror = 0
00116 !
00117       ijkend(:) = min (range(1,:)+levdim(:), range (2,:)) - range (1,:)
00118 !
00119 !     Control cells and warn user
00120 !
00121       dist_cyclic = period * 0.5
00122 
00123 loopk:   do k = 0, ijkend (3)
00124             do j = 0, ijkend (2)
00125 !cdir vector
00126                do i = 0, ijkend (1)
00127                if (chmax(i,j,k) - chmin (i,j,k) > dist_cyclic) exit loopk
00128                enddo
00129             enddo
00130          enddo loopk
00131 !
00132       if (k <= ijkend (3)) then
00133 !
00134 !     ... A cell violating the internal requirements was probably found
00135 !         Warn user !
00136 !
00137          print 9970, trim(ch_id), grid_id, trim(Grids(grid_id)%grid_name), &
00138                      name_coord(ind)
00139 !
00140          kbeg = k
00141             do k = kbeg, ijkend (3)
00142                do j = 0, ijkend (2)
00143                   do i = 0, ijkend (1)
00144                   if (chmax(i,j,k) - chmin (i,j,k) > dist_cyclic) then
00145                      ierror = ierror - 1
00146                      print 9960, i+range(1,1), j+range(1,2), k+range(1,3), &
00147                           array (i+range(1,1), j+range(1,2), k+range(1,3), :)
00148                   end if
00149                   end do
00150                end do
00151             end do
00152 !
00153          print *
00154       endif
00155 !
00156 !===> All done
00157 !
00158 #ifdef VERBOSE
00159       if ( ierror == 0 ) then
00160          print 9980, trim(ch_id), ierror
00161       else
00162          print 9981, trim(ch_id), abs(ierror)     
00163       endif
00164 
00165       call psmile_flushstd
00166 #endif /* VERBOSE */
00167 !
00168 !  Formats:
00169 !
00170 9990  format (1x, a, ': psmile_mg_ctrl_subgrid_3d_real')
00171 9980  format (1x, a, ': psmile_mg_ctrl_subgrid_3d_real eof: ierror =', i4)
00172 9981  format (1x, a, ': psmile_mg_ctrl_subgrid_3d_real eof: issue warning for', i8, ' cells.')
00173 
00174 9970  format (/1x, a, ': #### WARNING in psmile_mg_ctrl_subgrid_3d_real:', &
00175               /10x, 'Cells with incorrect (periodic) corner coordinates ', &
00176                     'probably found !' &
00177               /10x, 'Grid id ', i4, '; name:', 1x, a, &
00178               /10x, 'Cell indices of ', a, ' coordinate direction:')
00179 9960  format (1x, 'Indices(', i7, ',', i7, ',', i7, '); coordinates', &
00180               (1x, 4f16.9))
00181 !
00182       end subroutine PSMILe_MG_ctrl_subgrid_3d_real

Generated on 18 Mar 2011 for Oasis4 by  doxygen 1.6.1