psmile_compact_locations.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_compact_locations
00008 !
00009 ! !INTERFACE:
00010 
00011       subroutine psmile_compact_locations ( grid_id, search, ndim, found, ierror )
00012 !
00013 ! !USES:
00014 !
00015       use PRISM_constants
00016 !
00017       use PSMILe, dummy_interface => PSMILe_compact_locations
00018 
00019       implicit none
00020 !
00021 ! !INPUT PARAMETERS:
00022 !
00023       Integer, Intent (In)                  :: grid_id
00024 
00025 !     Pointer to the grid information
00026 
00027       Type (Enddef_search), Intent (In)     :: search
00028 
00029 !     Info's on locations to be processed
00030 
00031       Integer, Intent (In)                  :: ndim
00032 
00033 !     Dimension of found array
00034 
00035       Type (integer_vector), Intent (InOut) :: found (search%npart, ndim)
00036 !
00037 !     Finest level number on which a grid cell was found in the
00038 !     method grid for point I.
00039 !     found()%vector (i) =  1 = val_direct :
00040 !                               Point was found on a point of the method grid.
00041 !     found()%vector (i) = -1 = val_coupler:
00042 !                               Point was found in a cell  of the method grid
00043 !     Level number < -1       : Last level number on which point was found.
00044 !     Level number = -(nlev+1): Never found (input value)
00045 !
00046       Integer, Intent (Out)              :: ierror
00047 
00048 !     Returns the error code of PSMILe_MG_found_loc_to_3d;
00049 !             ierror = 0 : No error
00050 !             ierror > 0 : Severe error
00051 !
00052 ! !LOCAL VARIABLES
00053 !
00054       Integer, Parameter           :: indl = 1
00055       Integer, Parameter           :: indz = 2
00056       Integer                      :: compact_type
00057 
00058       Integer                      :: i, j, k, ipart, npart
00059 
00060       Integer                      :: nlev
00061       Integer                      :: ncpl
00062       Integer                      :: idx, idxr, idxu, idxur
00063       Integer                      :: sizei, sizej, sizek
00064 
00065       Integer                      :: range_3d(2,ndim_3d)
00066       Integer                      :: range   (2,ndim_3d,search%npart)
00067 
00068       Integer, Allocatable         :: index_3d(:,:,:)
00069       Integer, Allocatable         :: index_2d(:,:)
00070 !
00071 !     ... for error handling
00072 !
00073       Integer, Parameter           :: nerrp = 2
00074       Integer                      :: ierrp (nerrp)
00075 !
00076 ! !DESCRIPTION:
00077 !
00078 ! Subroutine "PSMILe_compact_locations" eliminates target locations which
00079 !   occure more than once by marking redundant points as not found.
00080 !
00081 ! !REVISION HISTORY:
00082 !   Date      Programmer   Description
00083 ! ----------  ----------   -----------
00084 ! 07.03.07    R. Redler  created
00085 !
00086 !EOP
00087 !----------------------------------------------------------------------
00088 !
00089 ! $Id: psmile_compact_locations.F90 2788 2010-11-30 14:34:07Z hanke $
00090 ! $Author: hanke $
00091 !
00092    Character(len=len_cvs_string), save :: mycvs = 
00093      '$Id: psmile_compact_locations.F90 2788 2010-11-30 14:34:07Z hanke $'
00094 !
00095 !----------------------------------------------------------------------
00096 !
00097 !  Initialization
00098 !
00099 #ifdef VERBOSE
00100       print 9990, trim(ch_id)
00101 
00102       call psmile_flushstd
00103 #endif /* VERBOSE */
00104 
00105       ierror = 0
00106       npart  = search%npart
00107       compact_type = search%msg_intersections%requires_conserv_remap
00108 
00109       nlev = Grids(grid_id)%nlev
00110 
00111       range = search%range
00112 
00113       ! determines maximum extent of part ranges
00114       do i = 1, ndim_3d
00115         range_3d(1,i) = minval(range(1,i,:) )
00116         range_3d(2,i) = maxval(range(2,i,:) )
00117       enddo
00118 
00119       ! does some preprocessing of the part ranges
00120       do ipart = 1, npart
00121          if ( compact_type == PSMILe_conserv2D ) then
00122 
00123              if ( Grids(grid_id)%grid_type == PRISM_Gaussreduced_regvrt .and. &
00124                   search%grid_type /= PRISM_Gaussreduced_regvrt .or. &
00125                   Grids(grid_id)%grid_type == PRISM_Reglonlatvrt) then
00126 
00127                   ! the lat column/row in 2d is omitted
00128                   ! for reglonlatvrt grids the only "the lower left" corners are
00129                   ! of interest, therefore we ignore the upper right ones
00130                   ! (does this also apply to irrlonlat*?)
00131                   range(2,1,ipart) = range(2,1,ipart)-1
00132                   range(2,2,ipart) = range(2,2,ipart)-1
00133              endif
00134 
00135          else if ( compact_type == PSMILe_conserv3D ) then
00136 
00137             range(2,1,ipart) = range(2,1,ipart)-1
00138             ! does this really only apply to PRISM_Gaussreduced_regvrt?
00139             if ( search%grid_type /= PRISM_Gaussreduced_regvrt ) &
00140             range(2,2,ipart) = range(2,2,ipart)-1
00141             range(2,3,ipart) = range(2,3,ipart)-1
00142 
00143          endif
00144       enddo
00145 
00146       select case ( Grids(grid_id)%grid_type )
00147 
00148       case ( PRISM_Reglonlatvrt, PRISM_Irrlonlatvrt )
00149 
00150          Allocate ( index_3d(range_3d(1,1):range_3d(2,1) , &
00151                              range_3d(1,2):range_3d(2,2) , &
00152                              range_3d(1,3):range_3d(2,3)), STAT = ierror)
00153 
00154          if ( ierror > 0 ) then
00155             ierrp (1) = ierror
00156             ierrp (2) = (range_3d(2,1) - range_3d(1,1) + 1) * &
00157                         (range_3d(2,2) - range_3d(1,2) + 1) * &
00158                         (range_3d(2,3) - range_3d(1,3) + 1)
00159             ierror = PRISM_Error_Alloc
00160             call psmile_error ( ierror, 'index_3d', &
00161                  ierrp, 2, __FILE__, __LINE__ )
00162             return
00163          endif
00164 
00165          index_3d = 0
00166 
00167          do ipart = 1, npart
00168 
00169             sizei = search%range(2,1,ipart) - search%range(1,1,ipart) + 1
00170             sizej = search%range(2,2,ipart) - search%range(1,2,ipart) + 1
00171             sizek = search%range(2,3,ipart) - search%range(1,3,ipart) + 1
00172 
00173             do k = range(1,3,ipart), range(2,3,ipart)
00174                do j = range(1,2,ipart), range(2,2,ipart)
00175                   do i = range(1,1,ipart), range(2,1,ipart)
00176                      idx = ( k - search%range(1,3,ipart) ) * sizei * sizej +  &
00177                            ( j - search%range(1,2,ipart) ) * sizei         +  &
00178                              i - search%range(1,1,ipart) + 1
00179 
00180                      if ( index_3d (i,j,k) == 0 ) then
00181                         if ( abs(found(ipart,1)%vector(idx)) == 1 ) index_3d (i,j,k) = 1
00182                      else
00183                         found(ipart,1)%vector(idx) = -(nlev+1)
00184                      endif
00185 
00186                   enddo
00187                enddo
00188             enddo
00189 
00190          enddo
00191 
00192          Deallocate ( index_3d, STAT = ierror)
00193 
00194          if ( ierror > 0 ) then
00195             ierrp (1) = ierror
00196             ierrp (2) = (range_3d(2,1) - range_3d(1,1) + 1) * &
00197                         (range_3d(2,2) - range_3d(1,2) + 1) * &
00198                         (range_3d(2,3) - range_3d(1,3) + 1)
00199             ierror = PRISM_Error_Dealloc
00200             call psmile_error ( ierror, 'index_3d', &
00201                  ierrp, 2, __FILE__, __LINE__ )
00202             return
00203          endif
00204 
00205       case ( PRISM_Irrlonlat_regvrt, PRISM_Gaussreduced_regvrt )
00206 
00207          Allocate ( index_2d(range_3d(1,1):range_3d(2,1) , &
00208                              range_3d(1,2):range_3d(2,2)), STAT = ierror)
00209 
00210          if ( ierror > 0 ) then
00211             ierrp (1) = ierror
00212             ierrp (2) = (range_3d(2,1) - range_3d(1,1) + 1) * &
00213                         (range_3d(2,2) - range_3d(1,2) + 1)
00214             ierror = PRISM_Error_Alloc
00215             call psmile_error ( ierror, 'index_2d', &
00216                  ierrp, 2, __FILE__, __LINE__ )
00217             return
00218          endif
00219 
00220          index_2d = 0
00221 
00222          do ipart = 1, npart
00223 
00224             sizei = search%range(2,1,ipart) - search%range(1,1,ipart) + 1
00225 
00226             ncpl = 0
00227 
00228             do j = range(1,2,ipart), range(2,2,ipart)
00229                do i = range(1,1,ipart), range(2,1,ipart)
00230                   idx = ( j - search%range(1,2,ipart) ) * sizei + &
00231                           i - search%range(1,1,ipart) + 1
00232 
00233                   if ( index_2d (i,j) == 0 ) then
00234                      if ( abs(found(ipart,indl)%vector(idx)) == 1 ) then
00235                         index_2d (i,j) = 1
00236                         ncpl = ncpl + 1
00237                      endif
00238                   else
00239                      found(ipart,indl)%vector(idx) = -(nlev+1)
00240                   endif
00241                enddo
00242             enddo
00243 
00244             if ( ncpl == 0 ) then
00245                do k = search%range(1,3,ipart), search%range(2,3,ipart)
00246                   found(ipart,indz)%vector(k-search%range(1,3,ipart)+1) = -(nlev+1)
00247                enddo
00248             endif
00249 
00250          enddo
00251 
00252          Deallocate ( index_2d, STAT = ierror)
00253 
00254          if ( ierror > 0 ) then
00255             ierrp (1) = ierror
00256             ierrp (2) = (range_3d(2,1) - range_3d(1,1) + 1) * &
00257                         (range_3d(2,2) - range_3d(1,2) + 1)
00258             ierror = PRISM_Error_Dealloc
00259             call psmile_error ( ierror, 'index_2d', &
00260                  ierrp, 2, __FILE__, __LINE__ )
00261             return
00262          endif
00263 
00264       case DEFAULT
00265 
00266          ierrp (1) = Grids(grid_id)%grid_type
00267          ierror = PRISM_Error_Internal
00268          call psmile_error ( ierror, 'unsupported grid type', &
00269               ierrp, 1, __FILE__, __LINE__ )
00270 
00271       end select
00272 !
00273 !-----------------------------------------------------------------------
00274 !     All done
00275 !-----------------------------------------------------------------------
00276 !
00277 #ifdef VERBOSE
00278       print 9980, trim(ch_id), ierror
00279 
00280       call psmile_flushstd
00281 #endif /* VERBOSE */
00282 !
00283 !  Formats:
00284 !
00285 9990  format (1x, a, ': psmile_compact_locations: ')
00286 9980  format (1x, a, ': psmile_compact_locations: eof ierror =', i4)
00287 !
00288       end subroutine PSMILe_compact_locations

Generated on 18 Mar 2011 for Oasis4 by  doxygen 1.6.1