psmile_def_partition.F90

Go to the documentation of this file.
00001 !
00002 !-----------------------------------------------------------------------
00003 ! Copyright 2007-2010, CERFACS, Toulouse, France.
00004 ! Copyright 2007-2010, SGI Germany, Munich, Germany.
00005 ! Copyright 2007-2010, NEC Europe Ltd., London, UK.
00006 ! All rights reserved. Use is subject to OASIS4 license terms.
00007 !-----------------------------------------------------------------------
00008 !BOP
00009 !
00010 ! !ROUTINE: psmile_def_partition
00011 !
00012 ! !INTERFACE:
00013 !
00014       subroutine psmile_def_partition ( grid_id, nbr_blocks, &
00015                                  partition_array, extent_array, ierror )
00016 !
00017 ! !USES:
00018 !
00019       use PRISM
00020       use psmile_grid, only : get_size_of_shape
00021       use PSMILe, dummy => psmile_def_partition
00022 
00023       implicit none
00024 !
00025 ! !INPUT PARAMETERS:
00026 
00027       Integer, Intent (In)                :: grid_id
00028 
00029 !     Specifies handle to the grid information created by routine prism_def_grid:
00030 
00031       Integer, Intent (In)                :: nbr_blocks
00032 
00033 !     number of blocks in the global index space, covered by the grid_valid_shape domain
00034 
00035       Integer, Intent (In)                :: partition_array(1:nbr_blocks,*)
00036 
00037 !     Vector containing for each block the offset for all dimensions in the
00038 !     global index space. The rank of partition_array depends on grid_type.
00039 !
00040 !RVSV  Proposal:
00041 
00042       Integer, Intent (In)                :: extent_array(1:nbr_blocks,*)
00043 !     For later use:
00044 !     Vector containing for each block the extent for all dimensions in the
00045 !     global index space. The rank of extent_array depends on grid_type.
00046 !     Needed if local blocks do not cover the full local partition, means
00047 !     nbr_blocks > 1.
00048 !     The extent_array is referenced only if nbr_blocks > 1. If the values
00049 !     are equal PRISM_UNDEFINED they are ignored.     
00050       
00051 !
00052 ! !OUTPUT PARAMETERS:
00053 !
00054       Integer, Intent (Out)               :: ierror
00055 
00056 !     Returns the error code of psmile_def_partition;
00057 !             ierror = 0 : No error
00058 !             ierror > 0 : Severe error
00059 !
00060 ! !LOCAL VARIABLES
00061 !
00062       integer            :: size_of_valid_shape(2)
00063 
00064       integer, parameter :: nerrp = 2
00065       integer            :: ierrp (nerrp)
00066 
00067       integer            :: i, len, lenz
00068 !
00069 ! !DESCRIPTION:
00070 !
00071 ! Subroutine "psmile_def_partition" defines the partition of the block
00072 !       indices relative to the global index space.
00073 !
00074 !
00075 ! !REVISION HISTORY:
00076 !   Date      Programmer   Description
00077 ! ----------  ----------   -----------
00078 ! 02.09.02    R. Redler    created
00079 ! 04.07.07    WP3a         interface revised.
00080 !
00081 !
00082 !EOP
00083 !----------------------------------------------------------------------
00084 !
00085 ! $Id: psmile_def_partition.F90 2820 2010-12-10 09:20:01Z hanke $
00086 ! $Author: hanke $
00087 !
00088   Character(len=len_cvs_string), save :: mycvs = 
00089       '$Id: psmile_def_partition.F90 2820 2010-12-10 09:20:01Z hanke $'
00090 !
00091 !----------------------------------------------------------------------
00092 !
00093 #ifdef VERBOSE
00094       print *, trim(ch_id), ': psmile_def_partition: grid_id', grid_id
00095 
00096       call psmile_flushstd
00097 #endif /* VERBOSE */
00098 !----------------------------------------------------------------------
00099 !  1st Initialization
00100 !----------------------------------------------------------------------
00101 
00102       ierror = 0
00103 
00104 !----------------------------------------------------------------------
00105 !  2nd Control Grid id
00106 !----------------------------------------------------------------------
00107 
00108       if (grid_id < 1 .or. &
00109           grid_id > Number_of_Grids_allocated) then
00110 
00111          ierrp (1) = grid_id
00112          ierrp (2) = Number_of_Grids_allocated
00113 
00114          ierror = PRISM_Error_Arg
00115 
00116          call psmile_error ( ierror, 'grid_id is not defined', &
00117                              ierrp, 2, __FILE__, __LINE__ )
00118          return
00119       endif
00120 !
00121       if (Grids(grid_id)%status == PSMILe_status_free) then
00122 
00123          ierrp (1) = grid_id
00124 
00125          ierror = PRISM_Error_Arg
00126 
00127          call psmile_error ( PRISM_Error_Arg, 'grid_id (not active)', &
00128                              ierrp, 1, __FILE__, __LINE__ )
00129          return
00130       endif
00131 
00132 !----------------------------------------------------------------------
00133 !  3rd Determine expected length of array partition_array and allocate
00134 !----------------------------------------------------------------------
00135 
00136       size_of_valid_shape = get_size_of_shape (Grids(grid_id)%grid_type)
00137       ! reduced gauss grids get special treatment
00138       if (Grids(grid_id)%grid_type == PRISM_Gaussreduced_regvrt .or. &
00139           Grids(grid_id)%grid_type == PRISM_Gaussreduced_sigmavrt) &
00140                size_of_valid_shape(2) = size_of_valid_shape(2) + 1
00141       len = size_of_valid_shape(2)
00142 
00143       Allocate ( Grids(grid_id)%partition(1:nbr_blocks,1:len), STAT = ierror )
00144 
00145       if (ierror > 0) then
00146          ierrp (1) = ierror
00147          ierrp (2) = 1
00148 
00149          ierror = PRISM_Error_Alloc
00150 
00151          call psmile_error ( ierror, 'Grids(grid_id)%partition', &
00152                              ierrp, 2, __FILE__, __LINE__ )
00153          return
00154       endif
00155 
00156 !----------------------------------------------------------------------
00157 !  4th Store data
00158 !----------------------------------------------------------------------
00159       
00160       Allocate ( Grids(grid_id)%extent(1:nbr_blocks,1:len), STAT = ierror )
00161 
00162       if (ierror > 0) then
00163          ierrp (1) = ierror
00164          ierrp (2) = 1
00165 
00166          ierror = PRISM_Error_Alloc
00167 
00168          call psmile_error ( ierror, 'Grids(grid_id)%extent', &
00169               ierrp, 2, __FILE__, __LINE__ )
00170          return
00171       endif
00172 
00173       if ( Grids(grid_id)%grid_type == PRISM_Gaussreduced_regvrt ) then
00174 
00175 !        Control partition and extension
00176 !
00177 !        Rene: Check implications for IO of Gaussgrids as the partitions
00178 !              and extents are used there as well
00179 !
00180 !cdir vector
00181          do i = 1, nbr_blocks-1
00182             if (max(partition_array (i,  1), partition_array (i+1,1)) <  &
00183                 min(partition_array (i  ,1)+extent_array (i  ,1),        &
00184                     partition_array (i+1,1)+extent_array (i+1,1)) .and.  &
00185                 max(partition_array (i,  2), partition_array (i+1,2)) <  &
00186                 min(partition_array (i  ,2)+extent_array (i  ,2),        &
00187                     partition_array (i+1,2)+extent_array (i+1,2))) exit
00188          end do
00189 !
00190          if (i < nbr_blocks) then
00191             ierror = PRISM_Error_Arg
00192             ierrp (1) = grid_id
00193             ierrp (2) = i
00194 !
00195             write (*, 9960) trim(ch_id), i,                               &
00196                             partition_array(i   ,1:2), extent_array(i,   1:2), &
00197                             partition_array(i+1, 1:2), extent_array(i+1, 1:2)
00198             call psmile_error ( ierror, 'Error in definition of partition: partiton(i)+extent(i) has to be <= partition(i+1)', &
00199                                 ierrp, 2, __FILE__, __LINE__ )
00200             
00201          endif
00202 
00203 !        Store partition and extents as 3-dimensional arrays
00204 
00205          lenz = Grids(grid_id)%grid_shape(2,3) - Grids(grid_id)%grid_shape(1,3)+1
00206 
00207          do i = 1, nbr_blocks
00208             if ( extent_array   (i,2) /= lenz ) print 9970, trim(ch_id), lenz
00209             if ( partition_array(i,2) /= 0    ) print 9980, trim(ch_id)
00210          enddo
00211 
00212          Grids(grid_id)%extent   (1:nbr_blocks,1) = extent_array   (1:nbr_blocks,1)
00213          Grids(grid_id)%partition(1:nbr_blocks,1) = partition_array(1:nbr_blocks,1)
00214 
00215          Grids(grid_id)%extent   (1:nbr_blocks,2) = 1
00216          Grids(grid_id)%partition(1:nbr_blocks,2) = 0
00217 
00218 #ifdef TODO
00219          Grids(grid_id)%extent   (1:nbr_blocks,3) = extent_array   (1:nbr_blocks,2)
00220          Grids(grid_id)%partition(1:nbr_blocks,3) = partition_array(1:nbr_blocks,2)
00221 #else
00222          Grids(grid_id)%extent   (1:nbr_blocks,3) = lenz
00223          Grids(grid_id)%partition(1:nbr_blocks,3) = 0
00224 #endif
00225 
00226       else
00227 
00228          Grids(grid_id)%partition(1:nbr_blocks,1:len) = partition_array(1:nbr_blocks,1:len)
00229          Grids(grid_id)%extent(1:nbr_blocks,1:len) = extent_array(1:nbr_blocks,1:len)
00230 
00231       endif
00232 
00233 !----------------------------------------------------------------------
00234 !  Epilogue
00235 !----------------------------------------------------------------------
00236 
00237 #ifdef VERBOSE
00238       print *, trim(ch_id), ': psmile_def_partition: eof ierror =', ierror
00239 
00240       call psmile_flushstd
00241 #endif /* VERBOSE */
00242 
00243 !  Formats
00244 
00245 9960  format (/1x, a, ': Error in partition number', i5, ':', &
00246               /1x, 'partition (i,   :) =', i7, ',', i7,   &
00247               /1x, 'extent    (i,   :) =', i7, ',', i7,   &
00248               /1x, 'partition (i+1, :) =', i7, ',', i7,   &
00249               /1x, 'extent    (i+1, :) =', i7, ',', i7)
00250 9970  format (/1x, a, ': #### WARNING in psmile_def_partition:', &
00251               /10x, 'Partitioning in the vertical not supported ', &
00252                     'for Gaussreduced grids!' &
00253               /10x, 'We reset the vertical extent to ', i4)
00254 
00255 9980  format (/1x, a, ': #### WARNING in psmile_def_partition:', &
00256               /10x, 'Partitioning in the vertical not supported ', &
00257                     'for Gaussreduced grids!' &
00258               /10x, 'We reset the vertical offset to zero')
00259       end subroutine psmile_def_partition

Generated on 18 Mar 2011 for Oasis4 by  doxygen 1.6.1