psmile_bbcells_2d_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_bbcells_2d_real
00008 !
00009 ! !INTERFACE:
00010 
00011       subroutine psmile_bbcells_2d_real (           &
00012                       coords, coords_shape, sub_range, &
00013                       corner_shape,                    &
00014                       chmin, chmax, midp, levdim,      &
00015                       period, ierror)
00016 !
00017 ! !USES:
00018 !
00019       use PRISM_constants
00020 !
00021       use psmile_grid, only : psmile_transform_cell_cyclic
00022       use psmile, dummy_interface => psmile_bbcells_2d_real
00023 
00024       Implicit none
00025 !
00026 ! !INPUT PARAMETERS:
00027 !
00028       Integer,           Intent (In)  :: coords_shape (2, ndim_2d)
00029 
00030 !     Dimension of (method/grid) coordinate "coords"
00031 
00032       Real,              Intent (In)  :: 
00033          coords(coords_shape(1,1):coords_shape(2,1), 
00034                 coords_shape(1,2):coords_shape(2,2))
00035 
00036 !     Coordinates of the method
00037 
00038       Integer,           Intent (In)  :: sub_range (2, ndim_2d)
00039 
00040 !     Definintion of the subgrid for which the (bounding) cells should
00041 !     be computed.
00042 
00043       Integer,           Intent (In)  :: corner_shape (2, ndim_2d)
00044 
00045 !     Dimension of the array "corners" containing the corner points
00046 !     of the control volume (for a single coordinate).
00047 
00048       Integer,           Intent (In)  :: levdim (ndim_2d)
00049 
00050 !     Dimension of output arrays "chmin", "chmax" and "midp".
00051 !     For testing.
00052 !
00053       Real,              Intent (In)  :: period
00054 !
00055 !     Type (Split_Information), Intent (InOut) :: split_info
00056 !
00057 ! !OUTPUT PARAMETERS:
00058 !
00059       Real,              Intent (Out) :: chmin (sub_range(1,1)-1:sub_range(2,1), 
00060                                                 sub_range(1,2)-1:sub_range(2,2))
00061 !
00062 !     Minimum value of bounding box of method cell
00063 !
00064       Real,              Intent (Out) :: chmax (sub_range(1,1)-1:sub_range(2,1), 
00065                                                 sub_range(1,2)-1:sub_range(2,2))
00066 !
00067 !     Maximum value of bounding box of method cell
00068 !
00069       Real,              Intent (Out) :: midp  (sub_range(1,1)-1:sub_range(2,1), 
00070                                                 sub_range(1,2)-1:sub_range(2,2))
00071 !
00072 !     Mid point of bounding box of method cell
00073 !
00074       Integer,           Intent (Out) :: ierror
00075 
00076 !     Returns the error code of psmile_bbcells_2d_real
00077 !             ierror = 0 : No error
00078 !             ierror > 0 : Severe error
00079 !
00080 ! !DEFINED PARAMETERS:
00081 !
00082 ! r_nbr = 1 / Number of corners per grid cell.
00083 !
00084       Integer,          Parameter  :: n_corners = 4
00085       Real,             Parameter  :: r_nbr = 0.25d0
00086 !
00087 ! !LOCAL VARIABLES
00088 !
00089       Integer                      :: i, ibeg, iend
00090       Integer                      :: j, jbeg, jend
00091 !
00092 !     Auxiliary fields for cyclic method cells
00093 !
00094       Real                         :: cell (n_corners)
00095 !
00096       Real                         :: period2
00097 !
00098 ! Error parameters
00099 !
00100       Integer, Parameter           :: nerrp = 2
00101       Integer                      :: ierrp (nerrp)
00102 !
00103 !
00104 ! !DESCRIPTION:
00105 !
00106 ! Subroutine "psmile_bbcells_2d_real" computes the bounding box for
00107 ! each cell in (sub_range(1,1):sub_range(2,1)-1, sub_range(1,2):sub_range(2,2)-1)
00108 ! (i.e. witin the interior of the 2-dimensional block)
00109 ! of the 2-dimensional coordinate defined by "array".
00110 !
00111 ! chmin (i, j) returns the minimal value of bounding box of
00112 !              method cell ((i,j), (i+1,j), (i,j+1), (i+1,j+1))
00113 ! chmax (i, j) returns the maximal value of bounding box of
00114 !              method cell ((i,j), (i+1,j), (i,j+1), (i+1,j+1))
00115 ! midp  (i, j) returns the mid point of
00116 !              method cell ((i,j), (i+1,j), (i,j+1), (i+1,j+1))
00117 !
00118 ! !REVISION HISTORY:
00119 !
00120 !   Date      Programmer   Description
00121 ! ----------  ----------   -----------
00122 ! 05.05.09    H. Ritzdorf  created
00123 ! 08.09.10    M. Hanke     simplified handling of cyclic coordinate ranges
00124 !
00125 !EOP
00126 !----------------------------------------------------------------------
00127 !
00128 !  $Id: psmile_bbcells_2d_real.F90 2687 2010-10-28 15:15:52Z coquart $
00129 !  $Author: coquart $
00130 !
00131    Character(len=len_cvs_string), save :: mycvs = 
00132        '$Id: psmile_bbcells_2d_real.F90 2687 2010-10-28 15:15:52Z coquart $'
00133 !
00134 !----------------------------------------------------------------------
00135 !
00136 #ifdef VERBOSE
00137       print 9990, trim(ch_id)
00138 
00139       call psmile_flushstd
00140 #endif /* VERBOSE */
00141 
00142 #ifdef PRISM_ASSERTION
00143       if (sub_range(1,1) < corner_shape(1,1) .or. &
00144           sub_range(2,1) > corner_shape(2,1) .or. &
00145           sub_range(1,2) < corner_shape(1,2) .or. &
00146           sub_range(2,2) > corner_shape(2,2)) then
00147          print *, trim(ch_id), "psmile_bbcells_2d_real: sub_range", sub_range
00148          print *, trim(ch_id), "psmile_bbcells_2d_real: corner_shape", &
00149                   corner_shape
00150 
00151          call psmile_assert (__FILE__, __LINE__, &
00152               "range must be in corner_shape")
00153       endif
00154 !
00155       if (minval(levdim(:)-(sub_range(2,:)-sub_range(1,:)+2)) < 0) then
00156          print *, trim(ch_id), "psmile_bbcells_2d_real: levdim", levdim
00157          print *, trim(ch_id), "psmile_bbcells_2d_real: sub_range ", sub_range
00158 
00159          call psmile_assert (__FILE__, __LINE__, &
00160               "levdim < sub_range(2,:)-sub_range(1,:)+2")
00161       endif
00162 #endif
00163 !
00164 !  Initialization
00165 !
00166       ierror = 0
00167 !
00168       ibeg = sub_range (1,1) - 1
00169       iend = sub_range (2,1) + 1
00170 !
00171       jbeg = sub_range (1,2) - 1
00172       jend = sub_range (2,2) + 1
00173 !
00174       period2 = period * 0.5d0
00175 !
00176 !=======================================================================
00177 !     Special case: 1-dimensional
00178 !=======================================================================
00179 !
00180 !     MoHa: does the 1d case require special handling of cyclic
00181 !           coordinate ranges?
00182 !
00183       if (sub_range(1,2) == sub_range (2,2)) then
00184 !
00185          j = jbeg
00186 
00187 !  ...... Set bounding box of method cells
00188 !
00189 !  ...... Compute values using the method coordinates
00190 !
00191 !cdir vector
00192          do i = ibeg, iend-1
00193             chmin (i, j) = min (coords(i,j), coords (i+1,j))
00194          enddo
00195 !
00196 !cdir vector
00197          do i = ibeg, iend-1
00198             chmax (i, j) = max (coords(i,j), coords (i+1,j))
00199          enddo
00200 !
00201 !cdir vector
00202          do i = ibeg, iend-1
00203             midp (i, j) = (coords(i,j) + coords (i+1,j)) * 0.5
00204          enddo
00205 !
00206       else if (sub_range(1,1) == sub_range(2,1)) then
00207 !
00208          i = ibeg
00209 
00210 !  ...... Set bounding box of method cells
00211 !
00212 !cdir vector
00213          do j = jbeg, jend-1
00214             chmin (i, j) = min (coords(i,j), coords (i,j+1))
00215          enddo
00216 !
00217 !cdir vector
00218          do j = jbeg, jend-1
00219             chmax (i, j) = max (coords(i,j), coords (i,j+1))
00220          enddo
00221 !
00222 !cdir vector
00223          do j = jbeg, jend-1
00224             midp (i, j) = (coords(i,j) + coords (i,j+1)) * 0.5
00225          enddo
00226       else
00227 !
00228 !=======================================================================
00229 !     Regular 2d array
00230 !=======================================================================
00231 !
00232 !     Compute minimum extension of method cell
00233 !
00234 
00235          do j = jbeg, jend-1
00236 !cdir vector
00237             do i = ibeg, iend-1
00238                chmin (i, j) = min (coords(i,j  ), coords (i+1,j  ), &
00239                                  coords(i,j+1), coords (i+1,j+1))
00240             enddo
00241 !
00242 !     Compute maximum extension of method cell
00243 !
00244 !cdir vector
00245             do i = ibeg, iend-1
00246                chmax (i, j) = max (coords(i,j  ), coords (i+1,j  ), &
00247                                  coords(i,j+1), coords (i+1,j+1))
00248             enddo
00249 !
00250 !     Compute mid point of method cell
00251 !
00252 !cdir vector
00253             do i = ibeg, iend-1
00254                midp (i, j)  = (coords(i,j  ) + coords (i+1,j  ) + &
00255                               coords(i,j+1) + coords (i+1,j+1)) * r_nbr
00256             enddo
00257          enddo
00258 !
00259 !     Apply handling for cyclic coordinate ranges
00260 !
00261          do j = jbeg, jend-1
00262 !cdir vector
00263             do i = ibeg, iend-1
00264                ! if the current cell might require special handling
00265                if (chmax(i,j)-chmin(i,j) > period2) then
00266 
00267                   ! get the cell
00268                   cell(1) = coords(i  ,j  ); cell(2) = coords(i+1,j)
00269                   cell(3) = coords(i+1,j+1); cell(4) = coords(i,j+1)
00270 
00271                   ! apply handling for cyclic coordinate ranges
00272                   ! (in case this handling was not required this
00273                   ! routine does nothing)
00274                   call psmile_transform_cell_cyclic (cell, period, ierror)
00275 
00276 
00277                   ! compute new min max and mid point values
00278                   chmin (i, j) = minval(cell)
00279                   chmax (i, j) = maxval(cell)
00280                   midp  (i, j) = sum(cell)  * r_nbr
00281                end if
00282             end do
00283          end do
00284       endif
00285 !
00286 !===> All done
00287 !
00288 #ifdef VERBOSE
00289       print 9980, trim(ch_id), ierror
00290 
00291       call psmile_flushstd
00292 #endif /* VERBOSE */
00293 !
00294 !  Formats:
00295 !
00296 9990 format (1x, a, ': psmile_bbcells_2d_real')
00297 9980 format (1x, a, ': psmile_bbcells_2d_real: eof ierror =', i3)
00298 !
00299       end subroutine psmile_bbcells_2d_real

Generated on 18 Mar 2011 for Oasis4 by  doxygen 1.6.1