psmile_range_subgrid_2d_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_Range_Subgrid_2d_dble
00008 !
00009 ! !INTERFACE:
00010 
00011       subroutine psmile_range_subgrid_2d_dble (           &
00012                     array1, array2,                       &
00013                             idlow, idhigh, jdlow, jdhigh, &
00014                             nbr_corners,                  &
00015                             ibeg, iend, jbeg, jend,       &
00016                     rinter, inter, ierror)
00017 !
00018 ! !USES:
00019 !
00020       use PRISM_constants
00021 !
00022       use PSMILe, dummy_interface => PSMILe_Range_Subgrid_2d_dble
00023 
00024       implicit none
00025 !
00026 ! !INPUT PARAMETERS:
00027 !
00028       integer, Intent (In)             :: idlow, idhigh
00029       integer, Intent (In)             :: jdlow, jdhigh
00030       integer, Intent (In)             :: nbr_corners
00031 
00032 !     Dimensions of "array1" and "array2"
00033 
00034       double precision, Intent (In)    :: array1 (idlow:idhigh, 
00035                                                   jdlow:jdhigh, 
00036                                                   nbr_corners)
00037       double precision, Intent (In)    :: array2 (idlow:idhigh, 
00038                                                   jdlow:jdhigh, 
00039                                                   nbr_corners)
00040 
00041 !     Fully dimensioned arrays for which the range of 2d-intersection
00042 !     should be computed.
00043 
00044       integer, Intent (In)             :: ibeg, iend, jbeg, jend
00045 
00046 !     Definintion of the subgrid
00047 
00048       double precision, Intent (In)    :: rinter (2, 2)
00049 
00050 !     Specifies the intersection for which the range should be found.
00051 !
00052 ! !OUTPUT PARAMETERS:
00053 !
00054       Integer, Intent (Out)            :: inter (2, 2)
00055 
00056 !     Returns the range of the intersection
00057 
00058       integer, Intent (Out)            :: ierror
00059 
00060 !     Returns the error code of PSMILe_Range_Subgrid_2d_dble;
00061 !             ierror = 0 : No error
00062 !             ierror > 0 : Severe error
00063 !
00064 ! !DEFINED PARAMETERS:
00065 !
00066       logical, parameter :: all = .true. ! Brut force
00067 !
00068 ! !LOCAL VARIABLES
00069 !
00070       integer            :: i, j
00071 !
00072       logical            :: same
00073 !
00074 ! !DESCRIPTION:
00075 !
00076 ! Subroutine "PSMILe_Range_Subgrid_2d_dble" computes the range of
00077 ! of grid contained in array~(ibeg:iend,~jbeg:jend,~nbr_corners) which contains
00078 ! the intersection "rinter"
00079 !
00080 !
00081 ! !REVISION HISTORY:
00082 !
00083 !   Date      Programmer   Description
00084 ! ----------  ----------   -----------
00085 ! 01.12.03    H. Ritzdorf  created
00086 !
00087 !EOP
00088 !----------------------------------------------------------------------
00089 !
00090 ! $Id: psmile_range_subgrid_2d_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_range_subgrid_2d_dble.F90 2325 2010-04-21 15:00:07Z valcke $'
00095 !
00096 !----------------------------------------------------------------------
00097 
00098 #ifdef VERBOSE
00099       print *, trim(ch_id), ': PSMILe_Range_Subgrid_2d_dble'
00100 
00101       call psmile_flushstd
00102 #endif /* VERBOSE */
00103 !
00104 !  Initialization
00105 !
00106       ierror = 0
00107 !
00108       inter (1, 1) = iend + 1
00109       inter (2, 1) = ibeg - 1
00110       inter (1, 2) = jend + 1
00111       inter (2, 2) = jbeg - 1
00112 !
00113       same = ibeg == idlow .and. iend == idhigh
00114 !
00115       if (same .or. all) then
00116 !
00117          do j = jbeg, jend
00118 !cdir vector
00119             do i = idlow, idhigh
00120             if (maxval(array1(i,j,:)) >= rinter (1,1) .and. &
00121                 minval(array1(i,j,:)) <= rinter (2,1) .and. &
00122                 maxval(array2(i,j,:)) >= rinter (1,2) .and. &
00123                 minval(array2(i,j,:)) <= rinter (2,2)) then
00124                inter (1, 1) = min (inter (1,1), i)
00125                inter (2, 1) = max (inter (2,1), i)
00126                inter (1, 2) = min (inter (1,2), j)
00127                inter (2, 2) = max (inter (2,2), j)
00128             endif
00129             end do
00130          end do
00131 !
00132          if (all .and. .not. same) then
00133             inter (1, 1) = max (inter (1,1), ibeg)
00134             inter (2, 1) = min (inter (2,1), iend)
00135          endif
00136 !
00137       else
00138 !
00139          do j = jbeg, jend
00140 !cdir vector
00141             do i = ibeg, iend
00142             if (maxval(array1(i,j,:)) >= rinter (1,1) .and. &
00143                 minval(array1(i,j,:)) <= rinter (2,1) .and. &
00144                 maxval(array2(i,j,:)) >= rinter (1,2) .and. &
00145                 minval(array2(i,j,:)) <= rinter (2,2)) then
00146                inter (1, 1) = min (inter (1,1), i)
00147                inter (2, 1) = max (inter (2,1), i)
00148                inter (1, 2) = min (inter (1,2), j)
00149                inter (2, 2) = max (inter (2,2), j)
00150             endif
00151             end do
00152          end do
00153       endif
00154 !
00155 !===> All done
00156 !
00157 #ifdef VERBOSE
00158       print *, trim(ch_id), ': PSMILe_Range_Subgrid_2d_dble eof', &
00159                             ': ierror =', ierror &
00160                           , '; inter  =', inter
00161 
00162       call psmile_flushstd
00163 #endif /* VERBOSE */
00164 !
00165       end subroutine PSMILe_Range_Subgrid_2d_dble

Generated on 18 Mar 2011 for Oasis4 by  doxygen 1.6.1