psmile_search_donor_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_Search_donor_extra_off
00008 !
00009 ! !INTERFACE:
00010 
00011       subroutine psmile_search_donor_extra_off (comp_info, &
00012                         search, var_id, tol, ierror)
00013 !
00014 ! !USES:
00015 !
00016       use PRISM
00017 !
00018       use PSMILe, dummy_interface => PSMILe_Search_donor_extra_off
00019 
00020       implicit none
00021 !
00022 ! !INPUT PARAMETERS:
00023 !
00024       Type (Enddef_comp),          Intent (In)    :: comp_info
00025 !
00026 !     Info on the component in which the donor cells
00027 !     should be searched.
00028 
00029       Integer,                     Intent (In)    :: var_id
00030 
00031 !     Handle to the grid function
00032 
00033       Double Precision,            Intent (In)    :: tol
00034 
00035 !     Absolute tolerance for search of "identical" points
00036 !     TOL >= 0.0
00037 
00038 !
00039 ! !INPUT/OUTPUT PARAMETERS:
00040 !
00041       Type (Enddef_global_search), Intent (InOut) :: search
00042 
00043 !     Data on the points to be searched
00044 !
00045 ! !OUTPUT PARAMETERS:
00046 !
00047       Integer,                     Intent (Out)   :: ierror
00048 
00049 !     Returns the error code of PSMILe_Search_donor_extra_off;
00050 !             ierror = 0 : No error
00051 !             ierror > 0 : Severe error
00052 !
00053 ! !LOCAL PARAMETERS
00054 !
00055 ! nd_dist  = Dummy parameter for routine "PSMILe_Return_extra_off_...."
00056 ! nb_extra = Dummy parameter for routine "PSMILe_Return_extra_off_...."
00057 !
00058       Integer, Parameter            :: nd_dist  = 1
00059       Integer, Parameter            :: nb_extra = 1
00060 !
00061 ! !LOCAL VARIABLES
00062 !
00063 !     ... loop variables
00064 !
00065       Integer                       :: i, j
00066 !
00067 !     ... field pointer
00068 !
00069       Type (Gridfunction), Pointer  :: field
00070       Type (Coords_Block), Pointer  :: coords_pointer
00071 !
00072 !     ... Grid pointer
00073 !
00074       Type (Grid), Pointer          :: grid_info
00075       Integer                       :: datatype, grid_id
00076       Integer                       :: method_id
00077 !
00078       Type (Corner_Block), Pointer  :: corner_pointer
00079 !
00080 !     ... for indices of points to be searched
00081 !
00082 ! len_item = Length of a data item per point in search%ibuf
00083 !
00084       Integer                       :: n_send, len_item, len_rtem
00085       Integer                       :: num_neigh
00086       Integer                       :: shift (ndim_3d)
00087 !
00088 !     ... for masks
00089 !         dummy_mask_array is a dummy array for to be
00090 !         transferred to interpolation routines
00091 !         Note: The target attributes (see !rr) were removed
00092 !               since there problems with a compiler and array bound checking
00093 !
00094       Integer                       :: mask_id
00095       Logical                       :: src_mask_available
00096 !
00097       Logical, Target               :: dummy_mask_array (1)
00098 !rr   Integer, target               :: dummy_mask_shape (2, ndim_3d)
00099       Integer                       :: dummy_mask_shape (2, ndim_3d)
00100 !
00101       Logical, Pointer              :: mask_array (:) 
00102 !rr   Integer, Pointer              :: mask_shape (:, :)
00103       Integer                       :: mask_shape (2, ndim_3d)
00104 !
00105 !     ... dummy distance
00106 !
00107       Integer                       :: df (nd_dist)
00108       Double Precision              :: ddist (nd_dist, nb_extra)
00109       Real                          :: rdist (nd_dist, nb_extra)
00110 !
00111 !     ... for error parameters
00112 !
00113       Integer, Parameter            :: nerrp = 1
00114       Integer                       :: ierrp (nerrp)
00115 !
00116 ! !DESCRIPTION:
00117 !
00118 ! Subroutine "PSMILe_Search_donor_extra_off" performs the additional (global)
00119 ! search for coordinates sent by the requesting process if offset are 
00120 ! available for the current grid and for the requesting coordinates.
00121 !
00122 ! !REVISION HISTORY:
00123 !
00124 !   Date      Programmer   Description
00125 ! ----------  ----------   -----------
00126 ! 03.07.05    H. Ritzdorf  created
00127 !
00128 !EOP
00129 !----------------------------------------------------------------------
00130 !
00131 ! $Id: psmile_search_donor_extra_off.F90 2967 2011-02-18 09:59:46Z hanke $
00132 ! $Author: hanke $
00133 !
00134    Character(len=len_cvs_string), save :: mycvs = 
00135        '$Id: psmile_search_donor_extra_off.F90 2967 2011-02-18 09:59:46Z hanke $'
00136 !----------------------------------------------------------------------
00137 !
00138 !  Initialization
00139 !
00140 #ifdef VERBOSE
00141       print 9990, trim(ch_id), comp_info%comp_id, search%sender
00142 
00143       call psmile_flushstd
00144 #endif /* VERBOSE */
00145 !
00146       ierror = 0
00147 !
00148       field     => Fields (var_id)
00149       method_id = field%method_id
00150       mask_id   = field%mask_id
00151 
00152       grid_id        =  Methods(method_id)%grid_id
00153       coords_pointer => Methods(method_id)%coords_pointer
00154 !
00155       grid_info      => Grids (grid_id)
00156       corner_pointer => grid_info%corner_pointer
00157 !
00158       datatype = corner_pointer%corner_datatype
00159 !
00160       src_mask_available = mask_id /= PRISM_UNDEFINED
00161 !
00162       if (src_mask_available) then
00163          mask_array => Masks(mask_id)%mask_array
00164 !rr      mask_shape => Masks(mask_id)%mask_shape
00165          mask_shape =  Masks(mask_id)%mask_shape
00166       else
00167          mask_array => dummy_mask_array
00168 !rr      mask_shape => dummy_mask_shape
00169          mask_shape =  dummy_mask_shape
00170       endif
00171 !
00172 #ifdef PRISM_ASSERTION
00173 !
00174 !===> Internal control
00175 !
00176       if (.not. Associated (grid_info%partition)) then
00177           print *, trim(ch_id), ": comp id, grid_id =", &
00178                    comp_info%comp_id, grid_id
00179           call psmile_assert ( __FILE__, __LINE__, &
00180                               'No partition info associated')
00181       endif
00182 !
00183 !     if (.not. Associated (grid_info%extent)) then
00184 !         print *, trim(ch_id), ": comp id, grid_id =", &
00185 !                  comp_info%comp_id, grid_id
00186 !         call PSMILe_Assert ( __FILE__, __LINE__, &
00187 !                             'No extent info associated')
00188 !     endif
00189 !
00190 #endif
00191 
00192 #ifdef VERBOSE
00193       if (.not. Associated (grid_info%extent)) then
00194          print 9970, trim(ch_id),                            &
00195                      grid_info%partition(1,1:ndim_3d),  &
00196                      grid_info%grid_shape(2,1:ndim_3d)- &
00197                      grid_info%grid_shape(1,1:ndim_3d)+1, &
00198                      search%msg_extra(12:14)
00199       else
00200          print 9970, trim(ch_id),                           &
00201                      grid_info%partition(1,1:ndim_3d), &
00202                      grid_info%extent(1,1:ndim_3d),    &
00203                      search%msg_extra(12:14)
00204       endif
00205 
00206       call psmile_flushstd
00207 #endif /* VERBOSE */
00208 !
00209 !  n_send = Number of points to be searched
00210 !
00211       n_send    = search%msg_extra (8)
00212       len_item  = search%msg_extra (9)
00213       len_rtem  = search%msg_extra (10)
00214       num_neigh = search%msg_extra (16)
00215 !
00216 !===> Set output parameter
00217 !
00218       search%n_liste = 0
00219       search%n_found = 0
00220 !
00221 !  Transfer global grid indices from remote grid
00222 !      into local  grid indices of  current grid
00223 !
00224 !  Note: The grid indices are not transformed for 
00225 !        PRISM_Gaussreduced_regvrt grids, since the
00226 !        local grid indices has to be computed using 
00227 !        coordinates.
00228 !
00229 #ifdef DEBUG
00230       print *, 'shift in extra', shift
00231       print *, 'partition', grid_info%partition (1, 1:ndim_3d)
00232       print *, 'extent', grid_info%extent (1, 1:ndim_3d)
00233       print *, 'grid_shape     ', grid_info%grid_shape(1, 1:ndim_3d)
00234       print *, 'grid_shape ende', grid_info%grid_shape(2, 1:ndim_3d)
00235       print *, 'cyclic', grid_info%cyclic
00236       print *, 'periodic', grid_info%periodic
00237       print *, 'len_periodic', grid_info%len_periodic
00238 #endif
00239 !
00240       if (grid_info%grid_type /= PRISM_Gaussreduced_regvrt) then
00241 
00242          shift (1:ndim_3d) = grid_info%grid_shape(1, 1:ndim_3d) - &
00243                              grid_info%partition (1, 1:ndim_3d) - 1
00244 !
00245 !cdir vector
00246          do i = 1, n_send
00247             search%ibuf ((i-1)*len_item+1) = &
00248             search%ibuf ((i-1)*len_item+1) + shift (1)
00249             search%ibuf ((i-1)*len_item+2) = &
00250             search%ibuf ((i-1)*len_item+2) + shift (2)
00251             search%ibuf ((i-1)*len_item+3) = &
00252             search%ibuf ((i-1)*len_item+3) + shift (3)
00253          end do
00254 !
00255          do j = 1, ndim_3d
00256             if (grid_info%periodic(j) == PSMILE_true) then
00257 !cdir vector
00258                do i = 1, n_send
00259                   if (search%ibuf ((i-1)*len_item+j) < grid_info%grid_shape(1,j)-1) then
00260                      search%ibuf ((i-1)*len_item+j) = &
00261                      search%ibuf ((i-1)*len_item+j) + grid_info%len_periodic(j)
00262                   else if (search%ibuf ((i-1)*len_item+j) > grid_info%grid_shape(2,j)+1) then
00263                      search%ibuf ((i-1)*len_item+j) = &
00264                      search%ibuf ((i-1)*len_item+j) - grid_info%len_periodic(j)
00265                   endif
00266                end do
00267             endif
00268          end do ! j = 1, ndim_3d
00269       endif
00270 !
00271 #ifdef DEBUGX
00272       print *, trim(ch_id), ': search%msg_extra(1:4)', search%msg_extra(1:4)
00273 !
00274       print *, 'psmile_search_donor_extra_off: grid shape', grid_info%grid_shape
00275       print *, 'srcloc (1:3), indices (j), len_item', len_item
00276       do i = 1, n_send
00277          print 8900, search%ibuf ((i-1)*len_item+1:i*len_item)
00278          print *, search%dbuf ((i-1)*len_rtem+1:i*len_rtem)
00279       end do
00280 #endif
00281 !
00282 !----------------------------------------------------------------------
00283 !    Depending on the interpolation mode
00284 !----------------------------------------------------------------------
00285 !
00286       select case (search%msg_extra(1))
00287 
00288       case (PSMILe_trilinear, PSMILe_bilinear, PSMILe_linear)
00289 
00290           if (grid_info%grid_type == PRISM_Gaussreduced_regvrt) then
00291           
00292              call psmile_trili_gauss2_extra (search, grid_id,     &
00293                       mask_array, mask_shape, src_mask_available, &
00294                       search%ibuf, len_item, n_send, num_neigh,   &
00295                       ierror)
00296 
00297           else
00298              call psmile_trili_3d_extra_off (comp_info, search,       &
00299                        grid_id,                                       &
00300                        mask_array, mask_shape, src_mask_available,    &
00301                        search%ibuf, len_item, n_send, num_neigh,      &
00302                        ierror)
00303           endif
00304 
00305       case (PSMILe_bicubic)
00306 
00307           if (grid_info%grid_type == PRISM_Gaussreduced_regvrt) then
00308 
00309              call psmile_tricu_gauss2_extra (search, grid_id,     &
00310                       mask_array, mask_shape, src_mask_available, &
00311                       search%ibuf, len_item, n_send, num_neigh,   &
00312                       ierror)
00313           else
00314 
00315              call psmile_tricu_3d_extra_off (comp_info, search,         &
00316                          mask_array, mask_shape, src_mask_available,    &
00317                          search%ibuf, len_item, n_send, num_neigh,      &
00318                          grid_info%grid_shape, grid_info%cyclic, ierror)
00319 
00320           endif
00321        
00322 !     Error: unsupported interpolation method
00323 
00324       case default
00325          ierrp (1) = search%msg_extra(1)
00326          ierror = PRISM_Error_Internal
00327 
00328          call psmile_error ( ierror, 'unsupported 3d interpolation method', &
00329                              ierrp, 1, __FILE__, __LINE__ )
00330          return
00331 
00332       end select
00333 !
00334 !----------------------------------------------------------------------
00335 !    Return info to the destination process
00336 !----------------------------------------------------------------------
00337 !
00338       if (datatype == MPI_REAL) then
00339          call psmile_return_extra_off_real (comp_info, search, var_id, &
00340                      df, rdist, nd_dist, nb_extra, ierror)
00341 
00342       else if (datatype == MPI_DOUBLE_PRECISION) then
00343          call psmile_return_extra_off_dble (comp_info, search, var_id, &
00344                      df, ddist, nd_dist, nb_extra, ierror)
00345 
00346 #if defined ( PRISM_QUAD_TYPE )
00347       else if (datatype == MPI_REAL16) the
00348          call psmile_return_extra_off_dble (comp_info, search, var_id, &
00349                      df, qdist, nd_dist, nb_extra, ierror)
00350 #endif
00351       endif
00352       if (ierror > 0) return
00353 !
00354 !===> All done
00355 !
00356 #ifdef VERBOSE
00357       print 9980, trim(ch_id), grid_id, search%sender, ierror
00358 
00359       call psmile_flushstd
00360 #endif /* VERBOSE */
00361 !
00362       return
00363 !
00364 !  Formats:
00365 !
00366 #ifdef VERBOSE
00367 9990 format (1x, a, ': psmile_search_donor_extra_off: comp_id =', i3, &
00368                     '; sender =', i4)
00369 9980 format (1x, a, ': psmile_search_donor_extra_off: eof grid id =', i3, &
00370                     '; sender =', i3, ', ierror =', i4)
00371 
00372 9970 format (1x, a, ': psmile_search_donor_extra_off:', &
00373              1x, 'offset local  ', 3i7, &
00374             /1x, 'extent        ', 3i7, &
00375             /1x, 'offset remote ', 3i7)
00376 #endif /* VERBOSE */
00377 8900 format (1x, 3i6, 5(';', i6) )
00378 
00379       end subroutine PSMILe_Search_donor_extra_off

Generated on 18 Mar 2011 for Oasis4 by  doxygen 1.6.1