psmile_mg_cells_gauss2.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_mg_cells_gauss2
00009 !
00010 ! !INTERFACE:
00011       subroutine psmile_mg_cells_gauss2 ( grid_id,    &
00012                       search_grid_type, found, loc,   &
00013                       loc_fnd_shape, control, ierror )
00014 !
00015 ! !USES:
00016 !
00017       use PRISM_constants
00018 !
00019       use PSMILe, dummy_interface => PSMILe_mg_cells_gauss2
00020 #ifdef DEBUG_TRACE
00021       use psmile_debug_trace
00022 #endif
00023 
00024       Implicit none
00025 !
00026 ! !INPUT PARAMETERS:
00027 !
00028       Integer, Intent (In)          :: grid_id
00029 
00030       Integer, Intent (In)          :: search_grid_type
00031 
00032       Integer, Intent (In)          :: loc_fnd_shape (2, ndim_3d)
00033 
00034 !     Dimension of loc and found
00035 
00036       Integer, Intent (In)          :: control (2, ndim_3d)
00037 
00038 !     Input  value: index loc_fnd_shape in "coords" to be searched
00039 !     Output value: index loc_fnd_shape in "coords" for which a cell was found
00040 !
00041 ! !INPUT/OUTPUT PARAMETERS:
00042 !
00043       Integer, Intent (InOut)         :: found ( loc_fnd_shape(1,1):loc_fnd_shape(2,1), 
00044                                                  loc_fnd_shape(1,2):loc_fnd_shape(2,2), 
00045                                                  loc_fnd_shape(1,3):loc_fnd_shape(2,3) )
00046 
00047 !     Finest level number on which a grid cell was found for point i,j,k.
00048 !     Input :
00049 !     Level number < -1: Point was not found and
00050 !                        and last level number was (-found(i,j,k))
00051 !     Level number = -(nlev+1): Never found
00052 !     Output :
00053 !     found(i,j,k) = +1: Point (i,j,k) is located on a point
00054 !                        of the method grid.
00055 !     found(i,j,k) = -1: Point (i,j,k) is located in a cell
00056 !                        of the method grid.
00057 !     found(,j,k) = -(nlev+2): Not located in the method grid.
00058 
00059       Integer, Intent (InOut)         :: loc ( loc_fnd_shape(1,1):loc_fnd_shape(2,1), 
00060                                                loc_fnd_shape(1,2):loc_fnd_shape(2,2), 
00061                                                loc_fnd_shape(1,3):loc_fnd_shape(2,3) )
00062 !
00063 !     Indices of the grid cell in which the point was found.
00064 !     The indices are relative to "shape".
00065 !     Note: The indices are indices in the Gauss Reduced Grid.
00066 !
00067 ! !OUTPUT PARAMETERS:
00068 !
00069       Integer, Intent (Out)           :: ierror
00070 
00071 !     Returns the error code of PSMILE_mg_cells_gauss2;
00072 !             ierror = 0 : No error
00073 !             ierror > 0 : Severe error
00074 !
00075 ! !DEFINED PARAMETERS:
00076 !
00077 !  val_coupler = Code for locations which should to be transferred
00078 !                to the coupler process.
00079 !
00080       Integer, Parameter              :: val_coupler = -1
00081 !
00082 !
00083 !     ... for levels
00084 !
00085       Integer, Parameter              :: lev = 1
00086 !
00087 ! !LOCAL VARIABLES
00088 !
00089       Integer                         :: i, j, k
00090 !
00091       Integer                         :: not_found 
00092 !
00093 !     ... for error parameters
00094 !
00095       Integer, Parameter              :: nerrp = 2
00096       Integer                         :: ierrp (nerrp)
00097 !
00098 #ifdef DEBUG_TRACE
00099 !
00100 !     ... for debugging
00101 !
00102       Integer                         :: m1, m2, m3, m4
00103 #endif
00104 !
00105 ! !DESCRIPTION:
00106 !
00107 ! At this point psmile_mg_final_gauss2_dble has already determined the exact
00108 ! source cell for each target point.
00109 ! This routine adjusts the found array according to the interpolation
00110 ! required for each target point. It is designed for grids of type
00111 ! "psmile_gaussreduced_regvrt".
00112 !
00113 !
00114 ! !REVISION HISTORY:
00115 !
00116 !   Date      Programmer   Description
00117 ! ----------  ----------   -----------
00118 ! 03.07.05    H. Ritzdorf  created
00119 ! 15.03.2011  M. Hanke     revisited
00120 !
00121 !EOP
00122 !----------------------------------------------------------------------
00123 !
00124 !  $Id: psmile_mg_cells_gauss2.F90 3023 2011-03-17 10:21:20Z hanke $
00125 !  $Author: hanke $
00126 !
00127    Character(len=len_cvs_string), save :: mycvs = 
00128   '$Id: psmile_mg_cells_gauss2.F90 3023 2011-03-17 10:21:20Z hanke $'
00129 !
00130 !----------------------------------------------------------------------
00131 !  Initialization
00132 !----------------------------------------------------------------------
00133 !
00134       ierror = 0
00135 !
00136       not_found = - (Grids(grid_id)%nlev + 1)
00137 !
00138 #ifdef VERBOSE
00139       print 9990, trim(ch_id), lev, control
00140       call psmile_flushstd
00141 #endif /* VERBOSE */
00142 
00143 !
00144 !----------------------------------------------------------------------
00145 !  Assertions
00146 !----------------------------------------------------------------------
00147 !
00148 #ifdef PRISM_ASSERTION
00149       do k = control(1,3), control(2,3)
00150          do j = control(1,2), control (2,2)
00151 !cdir vector
00152             do i = control(1,1), control (2,1)
00153                if ( found(i,j,k) == lev .and. &
00154                     ( loc(i,j,k) < grids(grid_id)%grid_shape(1,1) .or. &
00155                       loc(i,j,k) > grids(grid_id)%grid_shape(2,1)) ) exit
00156             end do
00157 
00158             if (i <= control (2,1)) then
00159                print *, 'i,j,k', i,j,k, ', loc', loc(i,j,k)
00160                print *, 'grids(grid_id)%grid_shape', grids(grid_id)%grid_shape(:,1)
00161                call psmile_assert (__FILE__, __LINE__, &
00162                     "Incorrect gauss index found")
00163             endif
00164          end do
00165       end do
00166 #endif
00167 !
00168 !----------------------------------------------------------------------
00169 !  Actual work
00170 !----------------------------------------------------------------------
00171 !
00172 !
00173 !===> ... Ensure that cell is marked as found if one corner was found 
00174 !         The lower left index will serve as cell index from now on.
00175 !
00176       if ( search_grid_type /= PRISM_Gaussreduced_regvrt ) then
00177 
00178          do k = loc_fnd_shape(1,3), loc_fnd_shape(2,3) ! only 2D search
00179             do j = loc_fnd_shape(1,2), loc_fnd_shape(2,2)-1
00180                do i = loc_fnd_shape(1,1), loc_fnd_shape(2,1)-1
00181 
00182                   if ( found(i  ,j  ,k) /= lev .and. &
00183                        (found(i+1,j  ,k) == lev .or. &
00184                         found(i  ,j+1,k) == lev .or. &
00185                         found(i+1,j+1,k) == lev ) ) then
00186 
00187                        found(i,j,k) = 1
00188 
00189                        if      ( abs(found(i+1,j  ,k)) == lev ) then
00190 
00191                           loc(i,j,k) = loc(i+1,j  ,k)
00192 
00193                        else if ( abs(found(i  ,j+1,k)) == lev ) then
00194 
00195                           loc(i,j,k) = loc(i  ,j+1,k)
00196 
00197                        else if ( abs(found(i+1,j+1,k)) == lev ) then
00198 
00199                           loc(i,j,k) = loc(i+1,j+1,k)
00200 
00201                        endif
00202 
00203                   endif
00204                enddo ! i
00205             enddo ! j
00206          enddo ! k
00207 
00208          ! set the points of the top row to not found, because we
00209          ! are searching for cells that are represented by the lower
00210          ! left corner
00211          do k = loc_fnd_shape(1,3), loc_fnd_shape(2,3)
00212             j = loc_fnd_shape(2,2)
00213             do i = loc_fnd_shape(1,1), loc_fnd_shape(2,1)
00214                found(i,j,k) = not_found
00215             enddo
00216          enddo
00217          ! set the points of the right most column 
00218          do k = loc_fnd_shape(1,3), loc_fnd_shape(2,3)
00219             i = loc_fnd_shape(2,1)
00220             do j = loc_fnd_shape(1,2), loc_fnd_shape(2,2)
00221                found(i,j,k) = not_found
00222             enddo
00223          enddo
00224 
00225       endif ! ( search_grid_type /= PRISM_Gaussreduced_regvrt )
00226 !
00227 !===> ... Mark cells that need to be send to the coupler for interpolation
00228 !         by "val_coupler"
00229 !
00230       where (found(:,:,:) == lev)
00231          found(:,:,:) = val_coupler
00232       endwhere
00233 !
00234 !-------------------------------------------------------------------------------
00235 !===> All done
00236 !-------------------------------------------------------------------------------
00237 !
00238 
00239 #ifdef DEBUG_TRACE
00240       if ( loc_fnd_shape(1,1) <= ictl_ind(1) .and. loc_fnd_shape(2,1) >= ictl_ind(1) .and. &
00241            loc_fnd_shape(1,2) <= ictl_ind(2) .and. loc_fnd_shape(2,2) >= ictl_ind(2)) then 
00242          print 8990, ictl_ind (1:ndim_2d), &
00243               found (ictl_ind(1),ictl_ind(2),control(1,3)), &
00244               loc   (ictl_ind(1),ictl_ind(2),control(1,3))
00245          m1 = 1
00246          m3 = 0
00247          do m2 = grids(grid_id)%grid_shape(1,1), loc (ictl_ind(1),ictl_ind(2),control(1,3))
00248             m4 = Grids(grid_id)%partition(m1,1)
00249             m3 = m3 + 1
00250             if ( m3 == Grids(grid_id)%extent(m1,1) ) then
00251                m1 = m1 + 1
00252                m3 = 0
00253             endif
00254          enddo
00255          print *, '       corresponding to loc global index ', m4+m3
00256 
00257 8990     format (1x, '### psmile_mg_cells_gauss2: ictl_ind', 2i5, &
00258               '; found, loc ', i3, 1x, i8)
00259       endif
00260 #endif /* DEBUG_TRACE */
00261 
00262 #ifdef VERBOSE
00263       print 9980, trim(ch_id), lev
00264 
00265       call psmile_flushstd
00266 #endif /* VERBOSE */
00267 !
00268 !  Formats:
00269 !
00270 9990  format (1x, a, ': psmile_mg_cells_gauss2: level', i3, &
00271            ', control', 6i6)
00272 9980  format (1x, a, ': psmile_mg_cells_gauss2: eof, level', i3)
00273 
00274       end subroutine psmile_mg_cells_gauss2

Generated on 18 Mar 2011 for Oasis4 by  doxygen 1.6.1