psmile_range_subgrid_3d_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_3d_dble
00008 !
00009 ! !INTERFACE:
00010 
00011       subroutine psmile_range_subgrid_3d_dble (           &
00012                     array1, array2, array3,               &
00013                             idlow, idhigh, jdlow, jdhigh, &
00014                             kdlow, kdhigh,                &
00015                             nbr_corners,                  &
00016                     shape, rinter, inter, ierror)
00017 !
00018 ! !USES:
00019 !
00020       use PRISM_constants
00021 !
00022       use PSMILe, dummy_interface => PSMILe_Range_Subgrid_3d_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)             :: kdlow, kdhigh
00031       integer, Intent (In)             :: nbr_corners
00032 
00033 !     Dimensions of "array1" and "array2"
00034 
00035       double precision, Intent (In)    :: array1 (idlow:idhigh, 
00036                                                   jdlow:jdhigh, 
00037                                                   kdlow:kdhigh, 
00038                                                   nbr_corners)
00039       double precision, Intent (In)    :: array2 (idlow:idhigh, 
00040                                                   jdlow:jdhigh, 
00041                                                   kdlow:kdhigh, 
00042                                                   nbr_corners)
00043       double precision, Intent (In)    :: array3 (idlow:idhigh, 
00044                                                   jdlow:jdhigh, 
00045                                                   kdlow:kdhigh, 
00046                                                   nbr_corners)
00047 
00048 !     Fully dimensioned arrays for which the range of 3d-intersection
00049 !     should be computed.
00050 
00051       integer, Intent (In)             :: shape (2, ndim_3d)
00052 
00053 !     Definintion of the subgrid
00054 
00055       double precision, Intent (In)    :: rinter (2, ndim_3d)
00056 
00057 !     Specifies the intersection for which the range should be found.
00058 !
00059 ! !OUTPUT PARAMETERS:
00060 !
00061       Integer, Intent (Out)            :: inter (2, ndim_3d)
00062 
00063 !     Returns the range of the intersection
00064 
00065       integer, Intent (Out)            :: ierror
00066 
00067 !     Returns the error code of PSMILe_Range_Subgrid_3d_dble;
00068 !             ierror = 0 : No error
00069 !             ierror > 0 : Severe error
00070 !
00071 ! !DEFINED PARAMETERS:
00072 !
00073       logical, parameter :: all = .true. ! Brut force
00074 !
00075 ! !LOCAL VARIABLES
00076 !
00077       integer            :: i, j, k
00078 !
00079       logical            :: same_i, same_ij
00080 !
00081 ! !DESCRIPTION:
00082 !
00083 ! Subroutine "PSMILe_Range_Subgrid_3d_dble" computes the range of
00084 ! of grid contained in array~(ibeg:iend,~jbeg:jend,~nbr_corners) which contains
00085 ! the intersection "rinter"
00086 !
00087 !
00088 ! !REVISION HISTORY:
00089 !   Date      Programmer   Description
00090 ! ----------  ----------   -----------
00091 ! 01.12.03    H. Ritzdorf  created
00092 !
00093 !EOP
00094 !----------------------------------------------------------------------
00095 !
00096 ! $Id: psmile_range_subgrid_3d_dble.F90 2325 2010-04-21 15:00:07Z valcke $
00097 ! $Author: valcke $
00098 !
00099    Character(len=len_cvs_string), save :: mycvs = 
00100        '$Id: psmile_range_subgrid_3d_dble.F90 2325 2010-04-21 15:00:07Z valcke $'
00101 !
00102 !----------------------------------------------------------------------
00103 
00104 #ifdef VERBOSE
00105       print *, trim(ch_id), ': PSMILe_Range_Subgrid_3d_dble'
00106 
00107       call psmile_flushstd
00108 #endif /* VERBOSE */
00109 !
00110 !  Initialization
00111 !
00112       ierror = 0
00113 !
00114       inter (1, :) = shape (2, :) + 1
00115       inter (2, :) = shape (1, :) - 1
00116 !
00117       same_i  = shape(1,1) == idlow .and. shape(2,1) == idhigh
00118       same_ij = shape(1,2) == jdlow .and. shape(2,2) == jdhigh .and. same_i
00119 !
00120       if (same_ij .or. all) then
00121 !
00122          do k = shape (1,3), shape (2,3)
00123             do j = jdlow, jdhigh
00124 !cdir vector
00125                do i = idlow, idhigh
00126                if (maxval(array1(i,j,k,:)) >= rinter (1,1) .and. &
00127                    minval(array1(i,j,k,:)) <= rinter (2,1) .and. &
00128                    maxval(array2(i,j,k,:)) >= rinter (1,2) .and. &
00129                    minval(array2(i,j,k,:)) <= rinter (2,2) .and. &
00130                    maxval(array3(i,j,k,:)) >= rinter (1,3) .and. &
00131                    minval(array3(i,j,k,:)) <= rinter (2,3)) then
00132                   inter (1, 1) = min (inter (1,1), i)
00133                   inter (2, 1) = max (inter (2,1), i)
00134                   inter (1, 2) = min (inter (1,2), j)
00135                   inter (2, 2) = max (inter (2,2), j)
00136                   inter (1, 3) = min (inter (1,3), k)
00137                   inter (2, 3) = max (inter (2,3), k)
00138                endif
00139                end do
00140             end do
00141          end do
00142 !
00143          if (all .and. .not. same_ij) then
00144             inter (1, :) = max (inter (1,:), shape (1, :))
00145             inter (2, :) = min (inter (2,:), shape (2, :))
00146          endif
00147 !
00148       else if (same_i) then
00149 !
00150          do k = shape (1,3), shape (2,3)
00151             do j = shape (1,2), shape (2,2)
00152 !cdir vector
00153                do i = idlow, idhigh
00154                if (maxval(array1(i,j,k,:)) >= rinter (1,1) .and. &
00155                    minval(array1(i,j,k,:)) <= rinter (2,1) .and. &
00156                    maxval(array2(i,j,k,:)) >= rinter (1,2) .and. &
00157                    minval(array2(i,j,k,:)) <= rinter (2,2) .and. &
00158                    maxval(array3(i,j,k,:)) >= rinter (1,3) .and. &
00159                    minval(array3(i,j,k,:)) <= rinter (2,3)) then
00160                   inter (1, 1) = min (inter (1,1), i)
00161                   inter (2, 1) = max (inter (2,1), i)
00162                   inter (1, 2) = min (inter (1,2), j)
00163                   inter (2, 2) = max (inter (2,2), j)
00164                   inter (1, 3) = min (inter (1,3), k)
00165                   inter (2, 3) = max (inter (2,3), k)
00166                endif
00167                end do
00168             end do
00169          end do
00170 !
00171       else
00172 !
00173          do k = shape (1,3), shape (2,3)
00174             do j = shape (1,2), shape (2,2)
00175 !cdir vector
00176                do i = shape (1,1), shape (2,1)
00177                if (maxval(array1(i,j,k,:)) >= rinter (1,1) .and. &
00178                    minval(array1(i,j,k,:)) <= rinter (2,1) .and. &
00179                    maxval(array2(i,j,k,:)) >= rinter (1,2) .and. &
00180                    minval(array2(i,j,k,:)) <= rinter (2,2) .and. &
00181                    maxval(array3(i,j,k,:)) >= rinter (1,3) .and. &
00182                    minval(array3(i,j,k,:)) <= rinter (2,3)) then
00183                   inter (1, 1) = min (inter (1,1), i)
00184                   inter (2, 1) = max (inter (2,1), i)
00185                   inter (1, 2) = min (inter (1,2), j)
00186                   inter (2, 2) = max (inter (2,2), j)
00187                   inter (1, 3) = min (inter (1,3), k)
00188                   inter (2, 3) = max (inter (2,3), k)
00189                endif
00190                end do
00191             end do
00192          end do
00193       endif
00194 !
00195 !===> All done
00196 !
00197 #ifdef VERBOSE
00198       print *, trim(ch_id), ': PSMILe_Range_Subgrid_3d_dble eof', &
00199                             ': ierror =', ierror
00200 !                         , '; inter  =', inter
00201 
00202       call psmile_flushstd
00203 #endif /* VERBOSE */
00204 !
00205       end subroutine PSMILe_Range_Subgrid_3d_dble

Generated on 18 Mar 2011 for Oasis4 by  doxygen 1.6.1