psmile_tricu_3d_extra_off.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_Tricu_3d_extra_off
00008 !
00009 ! !INTERFACE:
00010 
00011       subroutine psmile_tricu_3d_extra_off (comp_info, search, &
00012                         mask_array, mask_shape, mask_available,    &
00013                         ibuf, len_item, n_send, num_neigh,         &
00014                         grid_valid_shape, cyclic, ierror)
00015 !
00016 ! !USES:
00017 !
00018       use PRISM
00019 !
00020       use PSMILe, dummy_interface => PSMILe_Tricu_3d_extra_off
00021 
00022       implicit none
00023 !
00024 ! !INPUT PARAMETERS:
00025 !
00026       Type (Enddef_comp),          Intent (In)    :: comp_info
00027 !
00028 !     Info on the component in which the donor cells
00029 !     should be searched.
00030 
00031       Integer, Intent (In)            :: mask_shape (2, ndim_3d)
00032 
00033 !     Dimension of (method) mask array
00034 
00035       Logical, Intent (In)            :: mask_array (mask_shape (1,1): 
00036                                                      mask_shape (2,1), 
00037                                                      mask_shape (1,2): 
00038                                                      mask_shape (2,2), 
00039                                                      mask_shape (1,3): 
00040                                                      mask_shape (2,3))
00041 !     Mask of the method
00042 
00043       Logical, Intent (In)            :: mask_available
00044 
00045 !     Is mask specified in array "mask_array" ?
00046 !     mask_available = .false. : Mask is not available
00047 !                      .true.  : Mask is     available
00048 
00049       Integer,                     Intent (In)    :: len_item
00050 !
00051 !     Number of pieces of information in integer buffer "ibuf".
00052 !        1:3: (i,j,k) index of lower left corner
00053 !        4  : Index in
00054 !        5  : Code for points required
00055 !
00056       Integer,                     Intent (In)    :: n_send
00057 !
00058 !     Number of extra points sent
00059 
00060       Integer,                     Intent (In)    :: num_neigh
00061 !
00062 !     (Maximal) Number of interpolations bases required.
00063 !     Need to be the first "num_neigh" bases in "ijkstd" (see below).
00064 !
00065       Integer,                     Intent (In)    :: grid_valid_shape (2, 
00066                                                              ndim_3d)
00067 !
00068 !     Specifies the valid block shape for the "ndim_3d"-dimensional block
00069 !     (without halo/overlap regions) with
00070 !
00071 !     grid_valid_shape (1,i) = Lowest  "i"-th index of block
00072 !     grid_valid_shape (2,i) = Highest "i"-th index of block
00073 !
00074       Logical,            Intent (In) :: cyclic (ndim_3d)
00075 !    
00076 !     Does the block have cyclic coordinates in the corresponding direction ?
00077 !
00078 ! !INPUT/OUTPUT PARAMETERS:
00079 !
00080       Type (Enddef_global_search), Intent (InOut) :: search
00081 
00082 !     Data on the points to be searched
00083 
00084       Integer,                     Intent (InOut) :: ibuf (len_item, n_send)
00085 !
00086 !     Integer buffer containing the packed integer data.
00087 !     On return, ibuf (5, *) contains the code of the points found.
00088 !
00089 ! !OUTPUT PARAMETERS:
00090 !
00091       Integer,                     Intent (Out)   :: ierror
00092 
00093 !     Returns the error code of PSMILe_Search_donor_extra_off;
00094 !             ierror = 0 : No error
00095 !             ierror > 0 : Severe error
00096 !
00097 ! !LOCAL PARAMETERS:
00098 !
00099 ! n_corners_3d = Number of interpolation bases in 3d case
00100 !
00101 ! masked_out   = Hash value for a point which is masked out
00102 !
00103       Integer, Parameter           :: n_corners_3d = 4 * 4 * 4
00104 !
00105       Integer, Parameter           :: masked_out = 0
00106 !
00107 ! !LOCAL VARIABLES
00108 !
00109 !     ... loop variables
00110 !
00111       Integer                      :: i, j, n
00112 !
00113 !     ... 
00114 !
00115       Integer                      :: index0
00116       Integer                      :: length (ndim_3d)
00117 !
00118 !     ... For searched points of control volume
00119 !
00120       Integer                      :: code, n_corners, nlocs
00121 !
00122       Integer                      :: ijkstd (ndim_3d, n_corners_3d)
00123       Integer                      :: ijkctl (ndim_3d, n_corners_3d)
00124       Integer, Allocatable         :: locs (:, :)
00125       Integer, Allocatable         :: hash (:)
00126       Logical, Allocatable         :: ijkmsk (:, :)
00127 !
00128 !     ... for error parameters
00129 !
00130       Integer, parameter           :: nerrp = 1
00131       Integer                      :: ierrp (nerrp)
00132 !
00133 ! !DESCRIPTION:
00134 !
00135 ! Subroutine "PSMILe_Tricu_3d_extra_off" performs the additional (global)
00136 ! search for
00137 !
00138 ! TODO: PSMILe_Trili_3d_extra_off und PSMILe_Tricu_3d_extra_off vereinigen
00139 !       locs alloziieren und ausrechnen und danach gleiche Routine
00140 !       verwenden
00141 !       
00142 !
00143 ! !REVISION HISTORY:
00144 !
00145 !   Date      Programmer   Description
00146 ! ----------  ----------   -----------
00147 ! 03.07.05    H. Ritzdorf  created
00148 !
00149 !EOP
00150 !----------------------------------------------------------------------
00151 !
00152 ! $Id: psmile_tricu_3d_extra_off.F90 2082 2009-10-21 13:31:19Z hanke $
00153 ! $Author: hanke $
00154 !
00155    Character(len=len_cvs_string), save :: mycvs = 
00156        '$Id: psmile_tricu_3d_extra_off.F90 2082 2009-10-21 13:31:19Z hanke $'
00157 !
00158 !----------------------------------------------------------------------
00159 !
00160 !  Relative indices of corner
00161 !
00162       data ((ijkstd (i, n), i=1,ndim_3d), n = 1, n_corners_3d)            &
00163 !
00164                             /  -1,-1, 0,   0,-1, 0,   1,-1, 0,   2,-1, 0, &
00165                                -1, 0, 0,   0, 0, 0,   1, 0, 0,   2, 0, 0, &
00166                                -1, 1, 0,   0, 1, 0,   1, 1, 0,   2, 1, 0, &
00167                                -1, 2, 0,   0, 2, 0,   1, 2, 0,   2, 2, 0, &
00168 !
00169                                -1,-1,-1,   0,-1,-1,   1,-1,-1,   2,-1,-1, &
00170                                -1, 0,-1,   0, 0,-1,   1, 0,-1,   2, 0,-1, &
00171                                -1, 1,-1,   0, 1,-1,   1, 1,-1,   2, 1,-1, &
00172                                -1, 2,-1,   0, 2,-1,   1, 2,-1,   2, 2,-1, &
00173 !
00174                                -1,-1, 1,   0,-1, 1,   1,-1, 1,   2,-1, 1, &
00175                                -1, 0, 1,   0, 0, 1,   1, 0, 1,   2, 0, 1, &
00176                                -1, 1, 1,   0, 1, 1,   1, 1, 1,   2, 1, 1, &
00177                                -1, 2, 1,   0, 2, 1,   1, 2, 1,   2, 2, 1, &
00178 !
00179                                -1,-1, 2,   0,-1, 2,   1,-1, 2,   2,-1, 2, &
00180                                -1, 0, 2,   0, 0, 2,   1, 0, 2,   2, 0, 2, &
00181                                -1, 1, 2,   0, 1, 2,   1, 1, 2,   2, 1, 2, &
00182                                -1, 2, 2,   0, 2, 2,   1, 2, 2,   2, 2, 2 /
00183 !
00184 !----------------------------------------------------------------------
00185 !
00186 !  Initialization
00187 !
00188 #ifdef VERBOSE
00189       print 9990, trim(ch_id), comp_info%comp_id
00190 
00191       call psmile_flushstd
00192 #endif /* VERBOSE */
00193 !
00194       ierror = 0
00195       n_corners = num_neigh
00196 !
00197       length (1:ndim_3d) = grid_valid_shape(2,1:ndim_3d) - &
00198                            grid_valid_shape(1,1:ndim_3d) + 1
00199 !
00200 #ifdef PRISM_ASSERTION
00201       if (num_neigh < 1 .or. num_neigh > n_corners_3d) then
00202          print *, 'num_neigh', num_neigh, n_corners_3d
00203          call psmile_assert ( __FILE__, __LINE__, &
00204                 'Invalid value for num_neigh (number of interpolation bases')
00205       endif
00206 #endif
00207 !
00208       if (n_corners > 31) then
00209          ierror = PRISM_Error_internal
00210          ierrp (1) = n_corners
00211 
00212          call psmile_error ( ierror, & 
00213                             'Number of corners too large', &
00214                              ierrp, 1, __FILE__, __LINE__ )
00215          return
00216       endif
00217 !
00218 !===> Allocate mask "ijkmsk" for points to be searched.
00219 !     This corresponds to the decoded code "code".
00220 !
00221       Allocate (ijkmsk (n_send, n_corners), stat = ierror)
00222       if (ierror /= 0) then
00223          ierrp (1) = n_send * n_corners
00224 
00225          ierror = PRISM_Error_Alloc
00226 
00227          call psmile_error ( ierror, 'ijkmsk', ierrp, 1, &
00228                  __FILE__, __LINE__ )
00229          return
00230       endif
00231 !
00232 !===> Decode code for searched points and
00233 !     create list of locations "locs" of points required
00234 !
00235       code = 1
00236 !
00237          do j = 1, n_corners
00238 !cdir vector
00239             do i = 1, n_send
00240             ijkmsk (i,j) = IAND (code, ibuf (5, i)) /= 0
00241             end do
00242 
00243          code = code * 2
00244          end do
00245 !
00246       nlocs = Count (ijkmsk(:,:))
00247 !
00248 #ifdef PRISM_ASSERTION
00249       if (nlocs <= 0) then
00250          call psmile_assert ( __FILE__, __LINE__, &
00251                              'Number of points to be searched == 0 ?!?')
00252       endif
00253 #endif
00254 !
00255 #ifdef DEBUGX
00256       print *, 'code for locations searched', n_send
00257          do i = 1, n_send
00258          print *, 'code', ibuf (5, i), ijkmsk(i, :)
00259          end do
00260 #endif
00261 !
00262       Allocate (locs (ndim_3d, nlocs), hash (nlocs), stat = ierror)
00263       if (ierror /= 0) then
00264          ierrp (1) = (ndim_3d+1) * nlocs
00265 
00266          ierror = PRISM_Error_Alloc
00267 
00268          call psmile_error ( ierror, 'locs, hash', ierrp, 1, &
00269                  __FILE__, __LINE__ )
00270          return
00271       endif
00272 !
00273 !===> For special cases (2d Hyperplanes)
00274 !     ??? Braucht man hier nicht ein spezielles Flag,
00275 !         dass das globale Gebiet eine Hyperplane ist ?
00276 !
00277       ijkctl = ijkstd
00278 !
00279          do j = 1, ndim_3d
00280          if (length(j) == 1) ijkctl (j, :) = 0
00281          end do
00282 !
00283       n = 0
00284          do j = 1, n_corners
00285 !cdir vector
00286             do i = 1, n_send
00287             if (ijkmsk (i,j)) then
00288                n = n + 1
00289                locs (1, n) = ibuf (1, i) + ijkctl (1, j)
00290                locs (2, n) = ibuf (2, i) + ijkctl (2, j)
00291                locs (3, n) = ibuf (3, i) + ijkctl (3, j)
00292             end if
00293             end do ! i
00294          end do ! j
00295 !
00296 #ifdef PRISM_ASSERTION
00297       if (n /= nlocs) then
00298          print *, 'n, nlocs', n, nlocs
00299          call psmile_assert ( __FILE__, __LINE__, &
00300                 'Inconsistent values of n and nlocs; must be same value')
00301       endif
00302 #endif
00303 !
00304 !===> Set indices for cyclic coordinate directions
00305 !
00306 !rr This is not save. The locs(i,n) may become zero, and in this
00307 !rr case the fuction below will keep it at zero and will not apply
00308 !rr the cyclic condition correction, and the locs will be moved out
00309 !rr of the valid shape. Therefore an extra test is necessary to cover
00310 !rr this special case. 
00311 !rr
00312          do i = 1, ndim_3d
00313          if (cyclic(i) .and. length(i) > 1) then
00314             index0 = grid_valid_shape(1,i)
00315 !
00316 !cdir vector
00317                do n = 1, nlocs
00318                locs (i, n) = index0 + mod (locs(i, n) - index0, length(i))
00319                end do
00320 
00321                do n = 1, nlocs
00322                   if ( locs(i, n) < grid_valid_shape(1,i) ) &
00323                        locs(i, n) = locs(i, n) + length(i)
00324                enddo
00325          endif
00326          end do
00327 !
00328 !===> Create compact list of points to be sent and 
00329 !     return hash value "hash" (index in fortran order) for each location
00330 !
00331       call psmile_hash_extra (search, locs, hash, nlocs,                &
00332               mask_array, mask_shape, mask_available, grid_valid_shape, &
00333               ierror)
00334       if (ierror > 0) return
00335 !
00336 #ifdef DEBUGX
00337       print *, 'hash values and locs ', grid_valid_shape
00338       n = 0
00339          do j = 1, n_corners
00340 !cdir vector
00341             do i = 1, n_send
00342             if (ijkmsk (i,j)) then
00343                n = n + 1
00344                print *, 'i, j, n, hash(n), locs(:,n) ', &
00345                          i, j, n, hash(n), locs(:,n)
00346             end if
00347             end do ! i
00348          end do ! j
00349 #endif
00350 !
00351 !===> Save coded mask of points found in "ibuf"
00352 !     ??? Ob das so clever ist ?
00353 !
00354       if (search%n_liste > 0) then
00355          ibuf (5, 1:n_send) = 0
00356 !
00357          n = 0
00358          code = 1
00359             do j = 1, n_corners
00360 !cdir vector
00361                do i = 1, n_send
00362                if (ijkmsk (i,j)) then
00363 
00364                   n = n + 1
00365 !
00366                   if (hash (n) >= 0) then
00367 !
00368 !===> ...... Required "j-th" interpolation base was found
00369 !            (may be masked out (hash(n) == masked_out).
00370 !            Corresponding value is stored in position "index_found(n)".
00371 !
00372                      ibuf (5, i) = ibuf (5, i) + code
00373                   endif
00374 !
00375                endif
00376                end do ! i
00377 
00378             code = code * 2
00379             end do ! j
00380 
00381 #ifdef PRISM_ASSERTION
00382          if (n /= nlocs) then
00383             print *, 'nlocs, n', nlocs, n
00384             call psmile_assert ( __FILE__, __LINE__, &
00385                                 'n /= nlocs')
00386          endif
00387 #endif
00388 !
00389 #ifdef DEBUGX
00390         print *, 'code for locations found', n_send
00391             do i = 1, n_send
00392             print *, 'code', i, ibuf (5, i), ijkmsk(i, :)
00393             end do
00394 #endif 
00395 
00396       endif
00397 !
00398 !===> Deallocate memory allocated
00399 !
00400       Deallocate (locs, hash, ijkmsk)
00401 !
00402 !===> All done
00403 !
00404 #ifdef VERBOSE
00405       print 9980, trim(ch_id), ierror, search%n_found, n_send
00406 
00407       call psmile_flushstd
00408 #endif /* VERBOSE */
00409 !
00410       return
00411 !
00412 !  Formats:
00413 !
00414 #ifdef VERBOSE
00415 9990 format (1x, a, ': psmile_tricu_3d_extra_off: comp_id =', i3, &
00416                     '; sender =', i4)
00417 9980 format (1x, a, ': psmile_tricu_3d_extra_off: eof', &
00418                     ', ierror =', i4, ', n_found =', i5, i5)
00419 #endif /* VERBOSE */
00420 
00421       end subroutine PSMILe_Tricu_3d_extra_off

Generated on 18 Mar 2011 for Oasis4 by  doxygen 1.6.1