psmile_control_cell_gauss2_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_Control_cell_gauss2_real
00008 !
00009 ! !INTERFACE:
00010 
00011       subroutine psmile_control_cell_gauss2_real (grid_id,           &
00012                       ic, nc, list, j, k,                            &
00013                       coords1,  coords2,  shape,                     &
00014                       corners1, corners2, corner_shape, nbr_corners, &
00015                       tol, fnd, ierror)
00016 !
00017 ! !USES:
00018 !
00019       use PRISM_constants
00020 !
00021       use PSMILe, dummy_interface => PSMILe_Control_cell_gauss2_real
00022 #ifdef DEBUG_TRACE
00023       use psmile_debug_trace
00024 #endif
00025 
00026       Implicit none
00027 !
00028 ! !INPUT PARAMETERS:
00029 !
00030       Integer, Intent (In)            :: grid_id
00031 !
00032 !     Grid Id (for info degraded cells)
00033 !
00034       Integer, Intent (In)            :: nc
00035 !
00036 !     Number of cells (and points) to be controlled
00037 !
00038       Integer, Intent (In)            :: ic (nc)
00039 !
00040 !     Indices of the (corner) cells to be controlled
00041 !     Note: These indices are GaussReduced indices
00042 !
00043       Integer, Intent (In)            :: list (nc)
00044 !
00045 !     I-Indices of the points in to be controlled
00046 !
00047       Integer, Intent (In)            :: j
00048 !
00049 !     J-Index of the points in to be controlled
00050 !
00051       Integer, Intent (In)            :: k
00052 !
00053 !     K-Index of the points in to be controlled
00054 !     Note: list(:), j, k are indices of points to be searched !
00055 !
00056       Integer, Intent (In)            :: shape (2, ndim_3d)
00057 
00058 !     Dimension of coordinate arrays "coord1", "coord2"
00059 !     which contain the coordinates to be searched.
00060 
00061       Real, Intent (In)               ::   
00062          coords1 (shape(1,1):shape(2,1), shape(1,2):shape(2,2), 
00063                   shape(1,3):shape(2,3))
00064 
00065       Real, Intent (In)               ::   
00066          coords2 (shape(1,1):shape(2,1), shape(1,2):shape(2,2), 
00067                   shape(1,3):shape(2,3))
00068 
00069 !     Coordinates to be searched
00070 
00071       Integer, Intent (In)            :: corner_shape (2, ndim_2d)
00072 
00073 !     Dimension of corner coordinate arrays
00074 
00075       Integer, Intent (In)            :: nbr_corners
00076 
00077 !     Number of corners
00078 !     Note: This routine assumes nbr_corners == 2.
00079 
00080       Real,   Intent (In)             ::                 
00081          corners1 ( corner_shape(1,1):corner_shape(2,1), nbr_corners)
00082 
00083 !     Coordinates of corner points of control volume (Longitude)
00084 
00085       Real,   Intent (In)             ::                 
00086          corners2 ( corner_shape(1,1):corner_shape(2,1), nbr_corners)
00087 
00088 !     Coordinates of corner points of control volume (Latitude)
00089 
00090       Real, Intent (In)               :: tol
00091 !
00092 !     Absolute tolerance for search of "identical" points
00093 !     TOL >= 0.0
00094 !
00095 ! !OUTPUT PARAMETERS:
00096 !
00097       Logical, Intent (Out)           :: fnd (nc)
00098 !
00099 !     Was the i-th point of "list" found within cell "ic"
00100 !
00101       Integer, Intent (Out)           :: ierror
00102 
00103 !     Returns the error code of PSMILe_Control_cell_gauss2_real;
00104 !             ierror = 0 : No error
00105 !             ierror > 0 : Severe error
00106 !
00107 ! !DEFINED PARAMETERS:
00108 !
00109 !  NTET = Number of triangles used to determine whether
00110 !         point is within the grid cell
00111 !  EPS  = Tolerance for determinate of deformed grid cells
00112 !
00113       Integer, Parameter              :: ntet = 3
00114 !
00115       Real, Parameter                 :: eps = 1.0d-9
00116       Real, Parameter                 :: rmxeps = 1.0d30
00117       Real, Parameter                 :: dzero = 0.0
00118 !
00119 ! !LOCAL VARIABLES
00120 !
00121 !     ... for tolerances
00122 !
00123       Real                            :: tol1
00124 !
00125 !     ... for locations controlled
00126 !
00127       Integer                         :: i, jj, m, n
00128 !
00129 !     ... Logical reference point and edges (constant)
00130 ! 
00131       Integer, Save                   :: ijkedg (2, ntet)
00132 !
00133 !     ... Reference point and edges and distances
00134 ! 
00135       Real                            :: xyzpnt (nc, ntet, ndim_2d), 
00136                                          xyzref (nc, ntet+1, ndim_2d)
00137       Real                            :: xyzedg (nc, ntet, ndim_2d, 2)
00138       Real                            :: d1, d2
00139       Real                            :: det (nc, ntet)
00140       Real                            :: dnorm
00141 
00142 !
00143 ! !DESCRIPTION:
00144 !
00145 ! Subroutine "PSMILe_Control_cell_gauss2_real" controls whether the point
00146 ! "(list(i),j,k)" is contained in the corner cell "ic(i)".
00147 ! The routine assumes that the number of corners for the 2-dimensional
00148 ! cell is 4 (2 in longitude and 2 in latitude direction).
00149 !
00150 ! Subroutine "PSMILe_Control_cell_gauss2_real" is designed for grids of
00151 ! type "PRISM_Gaussreduced_regvrt".
00152 !
00153 !
00154 ! !REVISION HISTORY:
00155 !
00156 !   Date      Programmer   Description
00157 ! ----------  ----------   -----------
00158 ! 27.05.08    H. Ritzdorf  created
00159 !
00160 !EOP
00161 !----------------------------------------------------------------------
00162 !
00163 !  $Id: psmile_control_cell_gauss2_real.F90 2927 2011-01-28 14:04:12Z hanke $
00164 !  $Author: hanke $
00165 !
00166    Character(len=len_cvs_string), save :: mycvs = 
00167        '$Id: psmile_control_cell_gauss2_real.F90 2927 2011-01-28 14:04:12Z hanke $'
00168 !
00169 !----------------------------------------------------------------------
00170 !
00171 !  Relative indices of triangles
00172 !
00173 !     The cell is split into 3 triangles
00174 !
00175 !                   4-------3
00176 !                   !       !
00177 !                   !       !
00178 !            Corner 1 ----- 2
00179 !
00180       data ((ijkedg (jj, n), jj = 1, 2), n = 1, ntet) &
00181            / 2, 4,  3, 1,  4, 2/
00182 !   Corner    1      2      3
00183 !   Edges of corner 1 are connections to Corner 2 and Corner 4
00184 !
00185 !----------------------------------------------------------------------
00186 !
00187 !  Initialization
00188 !
00189 #ifdef VERBOSE
00190       print 9990, trim(ch_id), nc
00191 
00192       call psmile_flushstd
00193 #endif /* VERBOSE */
00194 !
00195 #ifdef PRISM_ASSERTION
00196       if (tol < 0.0) then
00197          call psmile_assert (__FILE__, __LINE__, &
00198                              "tol must be >= 0.0")
00199       endif
00200 !
00201 !cdir vector
00202          do i = 1, nc
00203          if (ic(i) < corner_shape(1,1) .or. ic(i) > corner_shape(2,1) ) then
00204             exit
00205          endif
00206          end do
00207 !
00208       if (i <= nc) then
00209          print *, 'i, ic', i, ic(i)
00210          print *, 'corner_shape', corner_shape
00211          call psmile_assert (__FILE__, __LINE__, &
00212                              'wrong index')
00213       endif
00214 #endif
00215 !
00216       ierror = 0
00217 !
00218 !===> Compute reference points, corners and edges of transformed
00219 !     triangles.
00220 !
00221       xyzref (:, 1, 1) = corners1 (ic(:), 1)
00222       xyzref (:, 1, 2) = corners2 (ic(:), 1)
00223 !
00224       xyzref (:, 2, 1) = corners1 (ic(:), 2)
00225       xyzref (:, 2, 2) = xyzref (:, 1, 2)
00226 !
00227       xyzref (:, 3, 1) = xyzref (:, 2, 1)
00228       xyzref (:, 3, 2) = corners2 (ic(:), 2)
00229 !
00230       xyzref (:, 4, 1) = xyzref (:, 1, 1)
00231       xyzref (:, 4, 2) = xyzref (:, 3, 2)
00232 !
00233 #ifdef DEBUG_TRACE
00234       do i = 1, nc
00235         if ( list(i) == ictl_ind(1) .and. j == ictl_ind(2) ) then
00236            print *, 'Search coords ', &
00237                coords1(list(i),j,shape(1,3)), coords2(list(i),j,shape(1,3))
00238            print *, 'Corners: '
00239            print *, '    C1 : ', xyzref (i, 1, 1), xyzref (i, 1, 2)
00240            print *, '    C2 : ', xyzref (i, 2, 1), xyzref (i, 2, 2)
00241            print *, '    C3 : ', xyzref (i, 3, 1), xyzref (i, 3, 2)
00242            print *, '    C4 : ', xyzref (i, 4, 1), xyzref (i, 4, 2)
00243         endif
00244       enddo
00245 #endif
00246         do n = 1, ntet
00247         xyzpnt (:, n, 1) = coords1 (list(:), j, k) - xyzref (:, n, 1)
00248         xyzpnt (:, n, 2) = coords2 (list(:), j, k) - xyzref (:, n, 2)
00249 !
00250            do m = 1, 2
00251 !cdir vector
00252               do i = 1, nc
00253               xyzedg (i, n, 1, m) = xyzref (i, ijkedg(m, n), 1) &
00254                                   - xyzref (i, n, 1)
00255 
00256               xyzedg (i, n, 2, m) = xyzref (i, ijkedg(m, n), 2) &
00257                                   - xyzref (i, n, 2)
00258               end do
00259            end do
00260         end do ! n
00261 !
00262 !===> Solve bilinear system; compute spat product
00263 !
00264 !   (xyzedg(*,1), xyzedg(*,2)) (x) = xyzpnt
00265 !
00266       det(:, :) = xyzedg(:, :, 1,1)*xyzedg(:, :, 2,2) &
00267                 - xyzedg(:, :, 1,2)*xyzedg(:, :, 2,1)
00268 !
00269 #define ALLOW_DEGRADED_CELL
00270 #ifdef ALLOW_DEGRADED_CELL
00271 !
00272 ! ??? ist das koscher ?? wird d1(n), d2(n), d3(n) == 0
00273 !
00274          do n = 1, ntet
00275 !cdir vector
00276             do i = 1, nc
00277             dnorm = abs(xyzedg(i, n,1,1)) + abs(xyzedg(i, n,2,1)) + &
00278                     abs(xyzedg(i, n,1,2)) + abs(xyzedg(i, n,2,2))
00279 
00280             if (abs(det(i, n)) < eps*dnorm*dnorm .or. dnorm == dzero) then
00281                det (i, n) = sign (rmxeps, det(i, n))
00282             else
00283                det (i, n) = 1.0 / det (i, n)
00284             endif
00285             end do ! i
00286          end do ! n
00287 
00288 #else 
00289 
00290 loop_n:  do n = 1, ntet
00291 !cdir vector
00292             do i = 1, nc
00293             dnorm = abs(xyzedg(i, n,1,1)) + abs(xyzedg(i, n,2,1)) + &
00294                     abs(xyzedg(i, n,1,2)) + abs(xyzedg(i, n,2,2))
00295 
00296             if (abs(det(i, n)) < eps*dnorm*dnorm .or. &
00297                 dnorm(i,n) == dzero) then
00298                det (i, n) = dnorm * 1e-25
00299 ! eigentlich exit loop_n und print_degraded cell
00300 ! ist bei vektorisierung verschwunden
00301 ! sollte man in part_list suchen
00302             endif
00303             end do ! i
00304          end do ! n
00305 !
00306 #ifdef PRINT_DEGRADED_CELL
00307       if (n < ntet) then
00308 !
00309 !===> cell is degraded
00310 !
00311          print 9960, trim(ch_id), grid_id, ic (i), &
00312                (xyzref (i, n, 1), xyzref (i, n, 2), n = 1, 4)
00313 !
00314             do n = 1, ntet
00315             if (abs(det(i, n)) < eps) then
00316                print 9950, n, det (i, n), dnorm, &
00317                  xyzedg (i, n,1, 1), xyzedg (i, n,1, 2), &
00318                  xyzedg (i, n,2, 1), xyzedg (i, n,2, 2)
00319             endif
00320             end do
00321 !
00322          nicht gefunden fnd (i) = .false. versucht, durch
00323          ganz kleine Determinante zu loesen
00324       endif
00325 #endif /* PRINT_DEGRADED_CELL */
00326 !
00327       det (:, :) = 1.0 / det (:, :)
00328 !
00329 #endif /* ALLOW_DEGRADED_CELL */
00330 !
00331 !===> Return mask for points "list(i)"
00332 !     which were found in the specific cell (ic(i,*), j,k)
00333 !
00334       tol1 = tol + 1.0
00335       fnd (:) = .false.
00336 !
00337          do n = 1, ntet
00338 !cdir vector
00339             do i = 1, nc
00340 
00341             d1 = (xyzpnt(i, n,1)*xyzedg(i, n,2,2)- &
00342                   xyzpnt(i, n,2)*xyzedg(i, n,1,2)) * det (i, n)
00343             d2 = (xyzpnt(i, n,2)*xyzedg(i, n,1,1)- &
00344                   xyzpnt(i, n,1)*xyzedg(i, n,2,1)) * det (i, n)
00345 !
00346             fnd (i) = fnd (i) .or. &
00347                (min(d1, d2) >= -tol .and. &
00348                 max(d1,dzero) + max(d2,dzero) <= tol1)
00349             end do ! nc
00350          end do ! ntet
00351 !
00352 !===> All done
00353 !
00354 #ifdef VERBOSE
00355       print 9980, trim(ch_id), count(fnd(:))
00356 
00357       call psmile_flushstd
00358 #endif /* VERBOSE */
00359 !
00360 !  Formats:
00361 !
00362 
00363 9990 format (1x, a, ': psmile_control_cell_gauss2_real: nc', i5)
00364 9980 format (1x, a, ': psmile_control_cell_gauss2_real: eof, suc', i5)
00365 
00366 #ifdef PRINT_DEGRADED_CELL
00367 9960 format (1x, a, ': psmile_control_cell_gauss2_real: grid id =', i5, &
00368             /1x,  ' Grid cell (', i4,
00369                  ') is degraded! Cell coordinates are:', 1p,      &
00370             4 (/1x, '(', 1 (e17.9, ','), e17.9, ')'))
00371 9950 format (1x, i1, '-th spawn; det =', 1p, e17.9, ' (', e17.9, &
00372              ') of', 2 (/1x, 2e17.9, ','))
00373 #endif
00374 
00375       end subroutine psmile_control_cell_gauss2_real

Generated on 18 Mar 2011 for Oasis4 by  doxygen 1.6.1