psmile_tricu_gauss2_extra.F90

Go to the documentation of this file.
00001 !-----------------------------------------------------------------------
00002 ! Copyright 2006-2010, NEC Europe Ltd., London, UK.
00003 ! Copyright 2011, DKRZ, Hamburg, Germany.
00004 ! All rights reserved. Use is subject to OASIS4 license terms.
00005 !-----------------------------------------------------------------------
00006 !BOP
00007 !
00008 ! !ROUTINE: psmile_tricu_gauss2_extra
00009 !
00010 ! !INTERFACE:
00011 
00012       subroutine psmile_tricu_gauss2_extra (search, grid_id, mask_array, &
00013                         mask_shape, mask_available, ibuf, len_item,      &
00014                         n_send, num_neigh, ierror)
00015 !
00016 ! !USES:
00017 !
00018       use PRISM
00019       use PSMILe, dummy_interface => psmile_tricu_gauss2_extra
00020       use psmile_grid_reduced_gauss
00021 
00022       Implicit none
00023 !
00024 ! !INPUT PARAMETERS:
00025 !
00026 
00027       Integer,                      Intent (In)    :: grid_id
00028 
00029 
00030       Integer,                      Intent (In)    :: mask_shape (2, ndim_3d)
00031 
00032 !     Dimension of (method) mask array
00033 
00034       Logical,                      Intent (In)    ::     
00035           mask_array ( mask_shape (1,1):mask_shape (2,1), 
00036                        mask_shape (1,2):mask_shape (2,2), 
00037                        mask_shape (1,3):mask_shape (2,3) )
00038 !     Mask of the method
00039 
00040       Logical,                     Intent (In)     :: mask_available
00041 
00042 !     Is mask specified in array "mask_array" ?
00043 !     mask_available = .false. : Mask is not available
00044 !                      .true.  : Mask is     available
00045 !
00046 
00047       Integer,                     Intent (In)     :: len_item
00048 !
00049 !     Number of pieces of information in integer buffer "ibuf".
00050 !        1:3: (i,j,k) index of lower left corner
00051 !        4  : Index in
00052 !        5  : Code for points required
00053 !
00054       Integer,                     Intent (In)    :: n_send
00055 !
00056 !     Number of extra points sent
00057 
00058       Integer,                     Intent (In)    :: num_neigh
00059 !
00060 !     (Maximal) Number of interpolations bases required..
00061 !
00062 !
00063 ! !INPUT/OUTPUT PARAMETERS:
00064 !
00065       Type (Enddef_global_search), Intent (InOut) :: search
00066 
00067 !     Data on the points to be searched
00068 
00069       Integer,                     Intent (InOut) :: ibuf (len_item, n_send)
00070 !
00071 !     Integer buffer containing the packed integer data.
00072 !     On return, ibuf (5, *) contains the code of the points found.
00073 !
00074 ! !OUTPUT PARAMETERS:
00075 !
00076       Integer,                     Intent (Out)   :: ierror
00077 
00078 !     Returns the error code of PSMILe_Search_donor_extra;
00079 !             ierror = 0 : No error
00080 !             ierror > 0 : Severe error
00081 !
00082 ! !LOCAL PARAMETERS:
00083 !
00084 ! n_corners_3d = Number of interpolation bases in 3d case
00085 !
00086 ! masked_out   = Hash value for a point which is masked out
00087 !
00088       Integer, Parameter           :: n_corners_3d = 4 * 4 * 4
00089 !
00090       Integer, Parameter           :: masked_out = 0
00091 !
00092 ! !LOCAL VARIABLES
00093 !
00094 !     ... loop variables
00095 !
00096       Integer                      :: j, n, i, k
00097 !
00098 !     ... For Grids
00099 !
00100       Type (Grid), Pointer         :: grid_info
00101       Integer,     Pointer         :: grid_valid_shape (:, :)
00102 !
00103 !     ... For Gaussreduced Grids
00104 !
00105       Logical                      :: virtual_cell_available
00106 !
00107 !     ... For searched points of control volume
00108 !
00109       Integer                      :: code, n_corners, outside
00110 !
00111       Integer, Allocatable         :: locs (:, :)
00112       Integer, Allocatable         :: neighbours_3d_extra (:, :, :)
00113       Integer, Allocatable         :: hash (:)
00114 !
00115 !
00116 !
00117       Logical                      :: missing_points(num_neigh, n_send)
00118       Integer                      :: decode_mask(num_neigh)
00119       Integer                      :: num_missing_points
00120 
00121       integer, parameter           :: shift_3d(ndim_3d, 9) = 
00122                                           reshape ((/-1,-1,0, 0,-1,0, +1,-1,0,
00123                                                      -1, 0,0, 0, 0,0, +1, 0,0,
00124                                                      -1,+1,0, 0,+1,0, +1,+1,0/), 
00125                                                    (/ndim_3d,9/))
00126 !
00127 !     ... for error parameters
00128 !
00129       Integer, parameter           :: nerrp = 1
00130       Integer                      :: ierrp (nerrp)
00131 !
00132 #ifdef DEBUG_TRACE
00133       Integer                      :: ictl_req, ind_loc (n_corners_3d)
00134 #endif
00135 !
00136 ! !DESCRIPTION:
00137 !
00138 ! Subroutine "psmile_tricu_gauss2_extra" performs the
00139 ! additional (global) search for grids of type "PRISM_Gaussreduced_regvrt".
00140 !
00141 ! !REVISION HISTORY:
00142 !
00143 !   Date      Programmer   Description
00144 ! ----------  ----------   -----------
00145 ! 03.07.05    H. Ritzdorf  created as psmile_tricu_gauss2_extra_dble.F90
00146 ! 27.01.11    M. Hanke     simplified code and corrected shifting, removed difference
00147 !                          between dble and real implementation
00148 !
00149 !EOP
00150 !----------------------------------------------------------------------
00151 !
00152 ! $Id: psmile_tricu_gauss2_extra.F90 3023 2011-03-17 10:21:20Z hanke $
00153 ! $Author: hanke $
00154 !
00155    Character(len=len_cvs_string), save :: mycvs = 
00156        '$Id: psmile_tricu_gauss2_extra.F90 3023 2011-03-17 10:21:20Z hanke $'
00157 !
00158 !----------------------------------------------------------------------
00159 !
00160 !  Initialization
00161 !
00162 #ifdef VERBOSE
00163       print 9990, trim(ch_id), grids(grid_id)%comp_id
00164 
00165       call psmile_flushstd
00166 #endif /* VERBOSE */
00167 !
00168       ierror = 0
00169       n_corners = num_neigh
00170 
00171 !
00172 !===> Get grid info
00173 !
00174 !  grid_valid_shape (2, ndim_3d) : Specifies the valid block shape
00175 !                                  for the "ndim_3d"-dimensional block
00176 !                                  (without halo/overlap regions)
00177 !
00178       grid_info => Grids(grid_id)
00179 !
00180       grid_valid_shape => grid_info%grid_shape
00181 !
00182       virtual_cell_available = len_item > ndim_3d + 2 ! unschoen
00183 !
00184       if (grid_info%grid_type /= PRISM_Gaussreduced_Regvrt) then
00185          print *, "grid_type", grid_info%grid_type
00186          call psmile_assert (__FILE__, __LINE__, &
00187                    "This routine is not designed for this grid type")
00188       endif
00189 !
00190       outside = grid_valid_shape (1, 1) - 11
00191 !
00192 #ifdef PRISM_ASSERTION
00193 
00194       if (num_neigh < 1 .or. num_neigh > n_corners_3d) then
00195          print *, 'num_neigh', num_neigh, n_corners_3d
00196          call psmile_assert ( __FILE__, __LINE__, &
00197                 'Invalid value for num_neigh (number of interpolation bases')
00198       endif
00199 !
00200 !===> Control input values contained in "ibuf"
00201 !
00202 !cdir vector
00203       do n = 1, n_send
00204          if (ibuf (2, n) < 1 .or. ibuf (2, n) > grid_info%nbr_latitudes) exit
00205       end do
00206 !
00207       if (n <= n_send) then
00208          print *, "i,j,k, gnbr_lats", ibuf(1:3,n), grid_info%nbr_latitudes
00209          call psmile_assert (__FILE__, __LINE__, &
00210                    "second (virtual) index out of range")
00211       endif
00212 !
00213 !cdir vector
00214       do n = 1, n_send
00215          if (ibuf (1, n) < 1 .or. &
00216              ibuf (1, n) > grid_info%nbr_points_per_lat(ibuf(2,n))) exit
00217       end do
00218 !
00219       if (n <= n_send) then
00220          print *, "i,j,k, gpoints_per_lat", &
00221                   ibuf(1:3,n), grid_info%nbr_points_per_lat(ibuf(2,n))
00222          call psmile_assert (__FILE__, __LINE__, &
00223                    "first (virtual) index out of range")
00224       endif
00225 #endif
00226 
00227 #ifdef DEBUG_TRACE
00228 !cdir vector
00229       do ictl_req = 1, n_send
00230          if (ibuf (4, ictl_req) == 1) exit
00231       end do
00232 #endif
00233 
00234       ! 1st we rebuilt the neighbours_3d (global 3d representation) array
00235       ! for the source locations received
00236 
00237       allocate (neighbours_3d_extra (ndim_3d, n_send, 16), stat=ierror)
00238       if (ierror /= 0) then
00239          call psmile_error ( PRISM_Error_Alloc, 'locs', (/ndim_3d * n_send * 16/), 1, &
00240                  __FILE__, __LINE__ )
00241          return
00242       endif
00243 
00244       if (virtual_cell_available) then
00245          ! for all locations
00246          do n = 1, n_send
00247 
00248             if (ibuf(6, n) == 0) then
00249 
00250                neighbours_3d_extra(:,n,1:16) = &
00251                   psmile_gauss_get_bicu_stencil (grid_id, ibuf(1:3,n))
00252 
00253             else
00254 
00255                neighbours_3d_extra(:,n,1:16) = &
00256                   psmile_gauss_shift_bicu_stencil (grid_id, ibuf(1:3,n), &
00257                                                                                 shift_3d(:,ibuf(6, n)))
00258 
00259             endif
00260          enddo ! n
00261       else
00262          neighbours_3d_extra(:,n,1:16) = psmile_gauss_get_bicu_stencil (grid_id, ibuf(1:3,n))
00263       endif
00264 !
00265 !  decode which points are missing
00266 !
00267       do i = 1, num_neigh
00268          decode_mask(i) = ishft (1, i - 1)
00269       enddo
00270 
00271       do n = 1, n_send
00272          missing_points(:, n) = iand (decode_mask(:), ibuf(5, n)) /= 0
00273       enddo
00274 !
00275 !  how many points are missing?
00276 !
00277       num_missing_points = count (missing_points)
00278 #ifdef PRISM_ASSERTION
00279       if (num_missing_points <= 0) then
00280          call psmile_assert ( __FILE__, __LINE__, &
00281                              'Number of points to be searched == 0 ?!?')
00282       endif
00283 #endif
00284 
00285 #ifdef DEBUGX
00286       print *, 'code for locations to be searched (n_send =', num_missing_points, ')'
00287       do n = 1, n_send
00288          print *, n, 'srcloc =', ibuf(1:3,n), 'code =', ibuf (5, n)
00289          print *, '   decoded     =', missing_points(:,n)
00290       end do
00291 #endif
00292 !
00293 #ifdef DEBUG_TRACE
00294       if (ictl_req <= n_send) then
00295          print *, 'code for indices_req', ibuf (4, ictl_req)
00296          print *, 'ictl_req', ictl_req, 'srcloc', ibuf(1:3,ictl_req), &
00297                   'code', ibuf (5, ictl_req), missing_points(:, ictl_req)
00298       endif
00299 #endif
00300 !
00301 !  Allocate memory for the indices of the missing points
00302 !
00303       Allocate (locs (ndim_3d, num_missing_points), stat = ierror)
00304       if (ierror /= 0) then
00305          call psmile_error ( PRISM_Error_Alloc, 'locs', (/ndim_3d * num_missing_points/), 1, &
00306                  __FILE__, __LINE__ )
00307          return
00308       endif
00309 
00310       k = 0
00311 
00312       do i = 1, n_corners
00313 !cdir vector
00314          do n = 1, n_send
00315             if (missing_points (i,n)) then
00316 
00317                k = k + 1
00318 
00319                locs (:, k) = neighbours_3d_extra(:, n, iand (i-1, 15)+1)
00320 
00321                if (grid_valid_shape(2,3) - grid_valid_shape(1,3) >= ishft(n_corners, -4) - 1 .and. &
00322                    n_corners > 16) then
00323                   locs (3, k) = locs (3, k) + ishft (i-1, -4)
00324                endif
00325             end if
00326          end do ! n
00327       end do ! i
00328 !
00329 #ifdef PRISM_ASSERTION
00330       if (k /= num_missing_points) then
00331          print *, 'k, num_missing_points', k, num_missing_points
00332          call psmile_assert ( __FILE__, __LINE__, &
00333                 'Inconsistent values of k and num_missing_points; must be same value')
00334       endif
00335 #endif
00336 
00337 #ifdef DEBUG_TRACE
00338       if (ictl_req <= n_send) then
00339          ind_loc (:) = 0
00340          n = 0
00341          do j = 1, n_corners
00342 !cdir vector
00343             do i = 1, n_send
00344                if (missing_points (j,i)) then
00345                   n = n + 1
00346                   if (i == ictl_req) ind_loc (j) = n
00347                end if
00348             end do ! i
00349          end do ! j
00350 !
00351          print *, 'locs searched for indices_req', ibuf (4, ictl_req)
00352          do j = 1, n_corners
00353             n = ind_loc (j)
00354             if (n > 0) &
00355             print *, 'ictl_req', ictl_req, ', n ', n, &
00356                            ': locs', locs (1:ndim_3d, n)
00357          end do
00358       endif
00359 #endif
00360 
00361 !
00362 !===> Transform indices from global (virtual) 3d indices
00363 !                       into local  (virtual) 3d indices
00364 !
00365       locs(1,:) = psmile_gauss_3d_to_local_1d (grid_id, locs(:,:), &
00366                                                num_missing_points, outside)
00367       locs(2,:) = 1
00368       !locs(3,:) = locs(3,:)
00369 !
00370 !===> Create compact list of points to be sent and
00371 !     return hash value "hash" (index in fortran order) for each location
00372 !
00373       Allocate (hash (num_missing_points), stat = ierror)
00374       if (ierror /= 0) then
00375          ierrp (1) = num_missing_points
00376 
00377          ierror = PRISM_Error_Alloc
00378 
00379          call psmile_error ( ierror, 'hash', ierrp, 1, &
00380                  __FILE__, __LINE__ )
00381          return
00382       endif
00383 !
00384       call psmile_hash_extra (search, locs, hash, num_missing_points,                &
00385               mask_array, mask_shape, mask_available, grid_valid_shape, &
00386               ierror)
00387       if (ierror > 0) return
00388 !
00389 #ifdef DEBUGX
00390       print *, 'hash values and locs ', grid_valid_shape
00391       n = 0
00392       do j = 1, n_corners
00393 !cdir vector
00394          do i = 1, n_send
00395             if (missing_points (j,i)) then
00396                n = n + 1
00397                print *, 'i, j, n, hash(n), locs(:,n) ', &
00398                          i, j, n, hash(n), locs(:,n)
00399             end if
00400          end do ! i
00401       end do ! j
00402 #endif
00403 !
00404 !===> Save coded mask of points found in "ibuf"
00405 !     ??? Ob das so clever ist ?
00406 !
00407       if (search%n_liste > 0) then
00408          ibuf (5, 1:n_send) = 0
00409 !
00410          n = 0
00411          code = 1
00412          do j = 1, n_corners
00413 !cdir vector
00414             do i = 1, n_send
00415                if (missing_points (j,i)) then
00416 
00417                   n = n + 1
00418 !
00419                   if (hash (n) >= 0) then
00420 !
00421 !===> ...... Required "j-th" interpolation base was found
00422 !            (may be masked out (hash(n) == masked_out).
00423 !            Corresponding value is stored in position "index_found(n)".
00424 !
00425                      ibuf (5, i) = ibuf (5, i) + code
00426                   endif
00427 !
00428                endif
00429             end do ! i
00430 
00431             code = code * 2
00432          end do ! j
00433 
00434 #ifdef PRISM_ASSERTION
00435          if (n /= num_missing_points) then
00436             print *, 'num_missing_points, n', num_missing_points, n
00437             call psmile_assert ( __FILE__, __LINE__, &
00438                                 'n /= num_missing_points')
00439          endif
00440 #endif
00441 !
00442 #ifdef DEBUGX
00443          print *, 'code for locations found', n_send
00444          do i = 1, n_send
00445             print *, 'code', i, ibuf (5, i), missing_points(:, i)
00446          end do
00447 #endif
00448 !
00449 #ifdef DEBUG_TRACE
00450          if (ictl_req <= n_send) then
00451             print *, 'code for locations found for indices_req', ibuf (4, ictl_req)
00452             print *, 'ictl_req', ictl_req, ': code', ibuf (5, ictl_req), missing_points(:, ictl_req)
00453          endif
00454 #endif
00455 
00456       endif
00457 !
00458 !===> Deallocate memory allocated
00459 !
00460       Deallocate (locs, hash, neighbours_3d_extra)
00461 !
00462 !===> All done
00463 !
00464 #ifdef VERBOSE
00465       print 9980, trim(ch_id), ierror, search%n_found, n_send
00466 
00467       call psmile_flushstd
00468 #endif /* VERBOSE */
00469 !
00470       return
00471 !
00472 !  Formats:
00473 !
00474 #ifdef VERBOSE
00475 9990 format (1x, a, ': psmile_tricu_gauss2_extra: comp_id =', i3)
00476 9980 format (1x, a, ': psmile_tricu_gauss2_extra: eof', &
00477                     ', ierror =', i4, ', n_found =', i8, i8)
00478 #endif /* VERBOSE */
00479 
00480       end subroutine psmile_tricu_gauss2_extra

Generated on 18 Mar 2011 for Oasis4 by  doxygen 1.6.1