psmile_neigh_tricu_gauss2_irreg.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_neigh_tricu_gauss2_irreg
00009 !
00010 ! !INTERFACE:
00011 
00012       subroutine psmile_neigh_tricu_gauss2_irreg (                &
00013                grid_id, grid_valid_shape, interpolation_mode,     &
00014                srcloc, srclocz, virtual_cell, nlocs, nloc, nprev, &
00015                virtual_cell_available, neighbors_3d, num_neigh,   &
00016                neigh_bascule, ierror)
00017 !
00018 ! !USES:
00019 !
00020       use PRISM_constants
00021       use PSMILe, dummy_interface => PSMILe_Neigh_tricu_gauss2_irreg
00022       use psmile_grid_reduced_gauss
00023 #ifdef DEBUG_TRACE
00024       use psmile_debug_trace
00025 #endif
00026 
00027       implicit none
00028 !
00029 ! !INPUT PARAMETERS:
00030 !
00031       Integer, Intent (In)            :: nloc
00032 !
00033 !     Total number of locations to be transferred
00034 !
00035       Integer, Intent (In)            :: nlocs (ndim_2d)
00036 !
00037 !     Number of locations in each direction
00038 !
00039       Integer, Intent (In)            :: srcloc  (ndim_2d, nlocs (1))
00040       Integer, Intent (In)            :: srclocz (nlocs (2))
00041 !
00042 !     Indices of the (method) grid cell in which the point was found.
00043 
00044       Integer, Intent (In)            :: nprev
00045 !
00046 !     Number of locations for which neighbours were already searched.
00047 !
00048       Integer, Intent (In)            :: virtual_cell (nlocs(1))
00049 
00050 !     Code whether destination point is located in virtual cell with
00051 !
00052 !     virtual_cell (n) = 0 : Not located in a virtual cell
00053 !     virtual_cell (n) > 0 :     Located in a virtual cell
00054 !                      = 1 : (-1,-1) direction
00055 !                      = 2 : ( 0,-1) direction
00056 !                      = 3 : ( 1,-1) direction
00057 !                      = 4 : (-1, 0) direction
00058 !                      = 6 : ( 1, 0) direction
00059 !                      = 7 : (-1, 1) direction
00060 !                      = 8 : ( 0, 1) direction
00061 !                      = 9 : ( 1, 1) direction
00062 
00063       Integer, Intent (In)            :: grid_id
00064 
00065 !     Grid handle
00066 
00067       Integer, Intent (In)            :: grid_valid_shape (2, ndim_3d)
00068 
00069 !     Specifies the valid block shape for the "ndim_3d"-dimensional block
00070 !
00071       Integer,            Intent (In) :: interpolation_mode
00072 !
00073 !     Specifies the interpolation mode with
00074 !        interpolation_mode = 3 : Tri-cubic interpolation
00075 !                                 (determines 64 corners) not supported
00076 !        interpolation_mode = 2 : Bilcubic interpolation
00077 !                                 (determines 16 corners)
00078 !
00079       Integer,            Intent (In) :: num_neigh
00080 !
00081 !     Last dimension of neighbors array "neighbors_3d"
00082 !
00083       Logical,            Intent (In) :: virtual_cell_available
00084 !
00085 !     Is the virtual cell info available ?
00086 !
00087 ! !OUTPUT PARAMETERS:
00088 !
00089 !     Location of indices for neighbors_3d
00090 !
00091 !           13   14   15   16
00092 !
00093 !     ^      9   10   11   12
00094 !     |
00095 !     |      5    6    7    8
00096 !
00097 !     j      1    2    3    4
00098 !
00099 !       i ---->
00100 !
00101       Integer, Intent (Out)           :: neighbors_3d (ndim_3d, nloc, num_neigh)
00102 
00103 !     Indices of neighbor locations (interpolation bases)
00104 
00105       Integer, Intent (Out)           :: neigh_bascule(nloc)
00106 
00107 !     Indicator for special index ordering
00108 
00109       Integer, Intent (Out)           :: ierror
00110 
00111 !     Returns the error code of PSMILE_Neigh_tricu_gauss2_irreg;
00112 !             ierror = 0 : No error
00113 !             ierror > 0 : Severe error
00114 !
00115 ! !DEFINED PARAMETERS:
00116 !
00117       integer, parameter              :: shift_3d(ndim_3d, 9) = 
00118                                              reshape ((/-1,-1,0, 0,-1,0, +1,-1,0,
00119                                                         -1, 0,0, 0, 0,0, +1, 0,0,
00120                                                         -1,+1,0, 0,+1,0, +1,+1,0/), 
00121                                                       (/ndim_3d,9/))
00122 !
00123 ! !LOCAL VARIABLES
00124 !
00125 !     ... for locations searched
00126 !
00127       Integer                         :: i, n
00128 !
00129 !     ... for error handling
00130 !
00131       Integer, Parameter              :: nerrp = 2
00132       Integer                         :: ierrp (nerrp)
00133 
00134 #ifdef DEBUG_TRACE
00135 !
00136 !     ... for debugging
00137 !
00138       integer                         :: ictl_neigh_global_1d(16), ictl_neigh_local_1d(16)
00139 #endif
00140 !
00141 ! !DESCRIPTION:
00142 !
00143 ! Subroutine "psmile_neigh_tricu_gauss2_irreg" generates the stencil for the cubic
00144 ! interpolation and stores the local 1d coordinates the stencil in neighbors_3d.
00145 ! In addition it stores the information whether a stencil covers the pole or
00146 ! not in neigh_bascule.
00147 !
00148 ! The subroutine is designed for (source) grids of type "PRISM_Gaussreduced_regvert"
00149 !  and target grids of type "PRISM_Irrlonlat_regvrt".
00150 !
00151 ! !REVISION HISTORY:
00152 !
00153 !   Date      Programmer   Description
00154 ! ----------  ----------   -----------
00155 ! 03.07.21    H. Ritzdorf  created for linear interpolation
00156 ! 05.05.09    R. Redler    modified for cubic interpolation
00157 ! 31.01.11    M. Hanke     rewritten based on psmile_neigh_tricu_gauss2.F90
00158 !
00159 !EOP
00160 !----------------------------------------------------------------------
00161 !
00162 !  $Id: psmile_neigh_tricu_gauss2_irreg.F90 2968 2011-02-18 15:26:34Z hanke $
00163 !  $Author: hanke $
00164 !
00165    Character(len=len_cvs_string), save :: mycvs = 
00166        '$Id: psmile_neigh_tricu_gauss2_irreg.F90 2968 2011-02-18 15:26:34Z hanke $'
00167 !
00168 !----------------------------------------------------------------------
00169 
00170 #ifdef VERBOSE
00171       print 9990, trim(ch_id)
00172       call psmile_flushstd
00173 #endif /* VERBOSE */
00174 
00175 !
00176 !===> inital checks
00177 !
00178 #ifdef PRISM_ASSERTION
00179 !
00180 !===> Control interpolation mode
00181 !
00182       if (interpolation_mode /= 2) then
00183          ierrp (1) = interpolation_mode
00184          ierror = PRISM_Error_Internal
00185 
00186          call psmile_error ( ierror, 'unsupported interpolation mode in psmile_neigh_tricu_gauss2_irreg', &
00187                              ierrp, 1, __FILE__, __LINE__ )
00188          return
00189       endif
00190 !
00191 !===> Is the dimension of array "neighbors_3d" sufficient
00192 !
00193       if (num_neigh < 16) then
00194          print 9970, trim(ch_id), num_neigh, 16
00195          call psmile_assert (__FILE__, __LINE__, &
00196                              'Number of neighbors is insufficient')
00197       endif
00198 !
00199 !===>
00200 !
00201       if (nloc < nprev + nlocs (1) * nlocs (2) ) then
00202          print *, trim(ch_id), " nloc, nprev, nlocs ", nloc, nprev, nlocs
00203          call psmile_assert (__FILE__, __LINE__, &
00204                              'nloc < nprev + PRODUCT (nlocs) ')
00205       endif
00206 !
00207 !===> Are the locations within the correct shape ?
00208 !
00209 !cdir vector
00210          do i = 1, nlocs (1)
00211             if ( abs(srcloc(1,i)) < grid_valid_shape(1,1) .or. &
00212                  abs(srcloc(1,i)) > grid_valid_shape(2,1) ) exit
00213          end do
00214 
00215          if (i <= nlocs(1)) then
00216             print *, "grid valid shape is ", grid_valid_shape
00217             print *, "Incorrect index in srcloc, i", i, srcloc(1,i)
00218             call psmile_assert (__FILE__, __LINE__, 'wrong index')
00219          endif
00220 
00221 !cdir vector
00222          do i = 1, nlocs (2)
00223             if ( srclocz(i) < grid_valid_shape(1,3)-1 .or. &
00224                  srclocz(i) > grid_valid_shape(2,3)) exit
00225          end do
00226 !
00227          if (i <= nlocs(2)) then
00228             print *, "Incorrect index in srclocs(2), i", i, srclocz(i), &
00229                      grid_valid_shape(1,3)-1, grid_valid_shape(2,3)
00230             call psmile_assert (__FILE__, __LINE__, 'wrong index')
00231          endif
00232 #endif
00233 !
00234 !===> Initialization
00235 !
00236       ierror = 0
00237 !
00238 !==> Actual work
00239 !
00240       ! 1st generate the 3d stencil
00241 
00242       if (virtual_cell_available) then
00243          ! for all locations
00244          do i = 1, nlocs(1)
00245 
00246             n = i + nprev
00247 
00248             if (virtual_cell(i) == 0) then
00249 
00250                neighbors_3d(:,n,1:16) = psmile_gauss_get_bicu_stencil (grid_id, &
00251                                           psmile_gauss_local_1d_to_3d (grid_id, abs(srcloc(1,i))))
00252 
00253             else
00254 
00255                neighbors_3d(:,n,1:16) = psmile_gauss_shift_bicu_stencil (grid_id,                  &
00256                                           psmile_gauss_local_1d_to_3d (grid_id, abs(srcloc(1,i))), &
00257                                           shift_3d(:,virtual_cell(i)))
00258 
00259             endif
00260          enddo ! í
00261       else
00262          do i = 1, nlocs(1)
00263 
00264             n = i + nprev
00265 
00266             neighbors_3d(:,n,1:16) = psmile_gauss_get_bicu_stencil (grid_id, &
00267                                        psmile_gauss_local_1d_to_3d (grid_id, abs(srcloc(1,i))))
00268          enddo ! i
00269       endif
00270 
00271 #ifdef DEBUG_TRACE
00272       if (ictl > nprev .and. ictl <= nprev + nlocs(1)) then
00273          print *, "neighbors_3d (3d):"
00274          do i = 1, 16
00275             print *, ' neighbors_3d (:,ictl,i) ', i, ':', neighbors_3d (:, ictl, i)
00276          enddo
00277       endif
00278 #endif
00279 
00280       ! 2nd check for stencils that go over the pole
00281 
00282       where (neighbors_3d(2,nprev+1:nprev+nlocs(1),6) >= grids(grid_id)%nbr_latitudes-1 .or. &
00283              neighbors_3d(2,nprev+1:nprev+nlocs(1),6) == 1)
00284 
00285          neigh_bascule(nprev+1:nprev+nlocs(1)) = 1
00286 
00287       else where
00288 
00289          neigh_bascule(nprev+1:nprev+nlocs(1)) = 0
00290 
00291       end where
00292 
00293 #ifdef DEBUG_TRACE
00294       if (ictl > nprev .and. ictl <= nprev + nlocs(1)) then
00295          print *, "neigh_bascule(ictl)", neigh_bascule(ictl)
00296       endif
00297 #endif
00298 
00299       ! 3rd transform the 3d coordinates into local 1d ones
00300 
00301       do i = 1, 16
00302 
00303          neighbors_3d(1,nprev+1:nprev+nlocs(1),i) = &
00304             psmile_gauss_3d_to_global_1d (grid_id, neighbors_3d(:,nprev+1:nprev+nlocs(1),i), nlocs(1))
00305          neighbors_3d(2,nprev+1:nprev+nlocs(1),i) = 1
00306 #ifdef DEBUG_TRACE
00307          if (ictl > nprev .and. ictl <= nprev + nlocs(1)) then
00308             ictl_neigh_global_1d(i) = neighbors_3d (1, ictl, i)
00309          endif
00310 #endif
00311 
00312          neighbors_3d(1,nprev+1:nprev+nlocs(1),i) = &
00313             psmile_gauss_1d_global_to_local (grid_id, neighbors_3d(1,nprev+1:nprev+nlocs(1),i), nlocs(1), 0)
00314 
00315 #ifdef DEBUG_TRACE
00316          if (ictl > nprev .and. ictl <= nprev + nlocs(1)) then
00317             ictl_neigh_local_1d(i) = neighbors_3d (1, ictl, i)
00318          endif
00319 #endif
00320       enddo
00321 
00322 #ifdef DEBUG_TRACE
00323       if (ictl > nprev .and. ictl <= nprev + nlocs(1)) then
00324          print *, "neighbors_3d (1d global):"
00325          print *, ' neighbors_3d (1,ictl,1:16) ', ictl_neigh_global_1d
00326          print *, "neighbors_3d (1d local):"
00327          print *, ' neighbors_3d (1,ictl,1:16) ', ictl_neigh_local_1d
00328       endif
00329 #endif
00330 !
00331       if ( num_neigh == 32 ) then
00332 !cdir vector
00333          do n = 1 + nprev, nprev + nlocs(1)
00334             neighbors_3d (1:2, n, 17:32) = neighbors_3d (1:2, n, 1:16)
00335          enddo
00336       endif
00337 !
00338 !===> Duplicate first two index for third direction
00339 !
00340       do i = 2, nlocs (2)
00341          neighbors_3d (1:2, nprev+(i-1)*nlocs(1)+1:nprev+i*nlocs(1),:) = &
00342          neighbors_3d (1:2, nprev+               1:nprev+  nlocs(1),:)
00343       end do
00344 !
00345 !===> Set third index
00346 !
00347       do i = 1, nlocs (2)
00348          neighbors_3d (3, nprev+(i-1)*nlocs(1)+1:nprev+i*nlocs(1), :) = srclocz (i)
00349       end do
00350 !
00351 !===> ... K indices of neighbour indices
00352 !
00353       if (grid_valid_shape(2,3) - grid_valid_shape(1,3) > 0 &
00354           .and. num_neigh == 32) then
00355          do n = 17, 32
00356                neighbors_3d (3, nprev+1:nprev+nlocs(1)*nlocs(2),n) = &
00357                   neighbors_3d (3, nprev+1:nprev+nlocs(1)*nlocs(2), 1) + 1
00358          enddo
00359       endif
00360 !
00361 !===> All done
00362 !
00363 #ifdef VERBOSE
00364       print 9980, trim(ch_id)
00365 
00366       call psmile_flushstd
00367 #endif /* VERBOSE */
00368 !
00369 !  Formats:
00370 !
00371 9990 format (1x, a, ': psmile_neigh_tricu_gauss2_irreg')
00372 9980 format (1x, a, ': psmile_neigh_tricu_gauss2_irreg: eof')
00373 
00374 #ifdef PRISM_ASSERTION
00375 9970 format (1x, a, ': number of neighbors', i2, '; expected at least', i2)
00376 #endif
00377 
00378       end subroutine psmile_neigh_tricu_gauss2_irreg

Generated on 18 Mar 2011 for Oasis4 by  doxygen 1.6.1