psmile_set_mask.F90

Go to the documentation of this file.
00001 !-----------------------------------------------------------------------
00002 ! Copyright 2007-2010, CERFACS, Toulouse, France.
00003 ! Copyright 2007-2010, SGI Germany, Munich, Germany.
00004 ! Copyright 2007-2010, NEC Europe Ltd., London, UK.
00005 ! All rights reserved. Use is subject to OASIS4 license terms.
00006 !-----------------------------------------------------------------------
00007 !BOP
00008 !
00009 ! !ROUTINE: psmile_set_mask
00010 !
00011 ! !INTERFACE:
00012 
00013       subroutine psmile_set_mask ( mask_id, grid_id, mask_actual_shape, &
00014                                            mask_array, new_mask, ierror )
00015 !
00016 ! !USES:
00017 !
00018       use PRISM
00019       use PSMILe, dummy => psmile_set_mask
00020 
00021       implicit none
00022 !
00023 ! !INPUT PARAMETERS:
00024 !
00025       integer, Intent (In)                :: grid_id
00026 
00027 !     Returned handle from prism_def_grid to the grid information
00028 
00029       logical, Intent (In)                :: mask_array (*)
00030 
00031 !     Array containg the mask data
00032 
00033       integer, Intent (In)                :: mask_actual_shape (1:2, *)
00034 
00035 !     Specifies the shape for the "n_dim"-dimensional mask
00036 !     (including halo/overlap regions) with
00037 
00038 !     mask_actual_shape (1,     1) = Lowest  first   dimension of mask
00039 !     mask_actual_shape (2,     1) = Highest first   dimension of mask
00040 !     mask_actual_shape (1,     2) = Lowest  second  dimension of mask
00041 !     mask_actual_shape (2,     2) = Highest second  dimension of mask
00042 !        ...
00043 !     mask_actual_shape (1,     i) = Lowest  i-th    dimension of mask
00044 !     mask_actual_shape (2,     i) = Highest i-th    dimension of mask
00045 !     mask_actual_shape (1, n_dim) = Lowest  n_dim-th dimension of mask
00046 !     mask_actual_shape (2, n_dim) = Highest n_dim-th dimension of mask
00047 
00048 !     where n_dim corresponds to value of input parameter of n_dim
00049 !     of corresponding subroutine call prism_def_grid.
00050 
00051       logical, Intent (In)                :: new_mask
00052 
00053 !     indicator whether a new mask_id was requested (.true.) or whether an
00054 !     old mask is going to be redefined (.false.)
00055 !
00056 ! !INPUT/OUTPUT PARAMETERS:
00057 !
00058       integer, Intent (InOut)             :: mask_id
00059 
00060 !     Returns the handle to the mask information if "new_mask = .true."
00061 !     is passed, otherwise, "mask_id" specifies the handle to the mask data
00062 !     to be changed.
00063 !
00064 ! !OUTPUT PARAMETERS:
00065 !
00066       integer, Intent (Out)               :: ierror
00067 
00068 !     Returns the error code of psmile_set_mask;
00069 !             ierror = 0 : No error
00070 !             ierror > 0 : Severe error
00071 !
00072 ! !LOCAL VARIABLES
00073 !
00074       integer                      :: mask_shape (1:2, ndim_3d)
00075 
00076       integer(kind=int64)          :: len_alloc
00077       integer                      :: i, n_dim
00078 !
00079       integer, parameter           :: nerrp = 2
00080       integer                      :: ierrp (nerrp)
00081 !
00082 ! !DESCRIPTION:
00083 !
00084 !  This routine defines a mask. The indivudual masks will be attached
00085 !  to a variable and a particular set of points in the call to
00086 !  prism_def_var via by using the mask_id.
00087 !
00088 !
00089 ! !REVISION HISTORY:
00090 !   Date      Programmer   Description
00091 ! ----------  ----------   -----------
00092 ! 01.12.03    H. Ritzdorf  created
00093 !
00094 !EOP
00095 !----------------------------------------------------------------------
00096 !
00097 ! $Id: psmile_set_mask.F90 2981 2011-02-24 09:44:15Z hanke $
00098 ! $Author: hanke $
00099 !
00100   Character(len=len_cvs_string), save :: mycvs = 
00101       '$Id: psmile_set_mask.F90 2981 2011-02-24 09:44:15Z hanke $'
00102 !
00103 !----------------------------------------------------------------------
00104 
00105 #ifdef VERBOSE
00106       print *, trim(ch_id), ': psmile_set_mask: mask_id', ' new_mask ', new_mask
00107 
00108       call psmile_flushstd
00109 #endif /* VERBOSE */
00110 
00111 !
00112 !  Initialization
00113 !
00114       ierror = 0
00115 !
00116 !  Control input arguments
00117 !
00118       if ( .not. new_mask ) then
00119 !
00120          if (mask_id < 1 .or. &
00121              mask_id > Number_of_Masks_allocated ) then
00122             ierrp (1) = mask_id
00123             ierrp (2) = Number_of_Masks_allocated
00124 
00125             ierror = PRISM_Error_Arg
00126 
00127             call psmile_error ( ierror, 'mask_id', &
00128                                 ierrp, 2, __FILE__, __LINE__ )
00129             return
00130          endif
00131 !
00132          if (Masks(mask_id)%mask_type == PSMILe_VectorMask) then
00133             ierrp (1) = mask_id
00134 
00135             ierror = PRISM_Error_Arg
00136 
00137             call psmile_error ( PRISM_Error_Arg, 'Use prism_set_vectormask for this mask_id', &
00138                                 ierrp, 1, __FILE__, __LINE__ )
00139             return
00140          endif
00141 !
00142          if (Masks(mask_id)%status == PSMILe_status_free) then
00143             ierrp (1) = mask_id
00144 
00145             ierror = PRISM_Error_Arg
00146 
00147             call psmile_error ( PRISM_Error_Arg, 'mask_id (not active)', &
00148                                 ierrp, 1, __FILE__, __LINE__ )
00149             return
00150          endif
00151 
00152 #ifdef VERBOSE
00153          print *, trim(ch_id), ': psmile_set_mask: overwriting old mask for mask_id ', &
00154                               mask_id
00155          call psmile_flushstd
00156 #endif /* VERBOSE */
00157 
00158       endif
00159 !
00160       if (Grids(grid_id)%grid_structure == PSMILe_Grid_Block) then
00161 
00162          n_dim = Grids(grid_id)%n_dim
00163 
00164          if ( Grids(grid_id)%grid_type == PRISM_Gaussreduced_regvrt ) then
00165            mask_shape(1:2,1) = mask_actual_shape(1:2,1)
00166            mask_shape(1:2,2) = 1
00167            mask_shape(1:2,3) = mask_actual_shape(1:2,2)
00168          else
00169            do i = 1, n_dim
00170               mask_shape(1:2,i) = mask_actual_shape(1:2,i)
00171            enddo
00172            do i = n_dim+1, ndim_3d
00173               mask_shape(1:2,i) = 1
00174            enddo
00175          endif
00176 !
00177          do i = 1, n_dim
00178 
00179             if (mask_shape(1,i) .gt. Grids(grid_id)%grid_shape(1,i) .or. &
00180                 mask_shape(2,i) .lt. Grids(grid_id)%grid_shape(2,i)) exit
00181          enddo
00182 
00183          if (i <= n_dim) then
00184 
00185              ierror = PRISM_Error_Arglist
00186 
00187              call psmile_error ( PRISM_Error_Arg, 'mask_actual_shape', &
00188                                  mask_actual_shape, n_dim*2, &
00189                                  __FILE__, __LINE__ )
00190              return
00191          endif
00192 
00193       else
00194 !
00195 !  3-dimensional mask for an invalid grid type
00196 !
00197          ierrp (1) = grid_id
00198 
00199          ierror = PRISM_Error_Arg
00200 
00201          call psmile_error ( PRISM_Error_Arg, 'grid_id (not for a 3-dimensional mask)', &
00202                              ierrp, 1, __FILE__, __LINE__ )
00203          return
00204       endif
00205 !
00206 !  Store Mask data
00207 !
00208       if (new_mask) then
00209 
00210 #ifdef VERBOSE
00211          print *, trim(ch_id), ': psmile_set_mask: getting new mask_id'
00212          call psmile_flushstd
00213 #endif /* VERBOSE */
00214 
00215          call psmile_get_mask_handle (mask_id, ierror)
00216          if (ierror > 0) return
00217 
00218       else if (Associated(Masks(mask_id)%mask_array)) then
00219 
00220          Deallocate (Masks(mask_id)%mask_array, STAT = ierror)
00221 
00222          if (ierror > 0) then
00223             ierrp (1) = ierror
00224             ierror = PRISM_Error_Dealloc
00225             call psmile_error ( ierror, 'mask_array', &
00226                                 ierrp, 1, __FILE__, __LINE__ )
00227             return
00228          endif
00229 
00230       endif
00231 !
00232 !  Get size of mask_array
00233 !
00234       len_alloc = 1
00235       do i = 1, n_dim
00236          len_alloc = len_alloc * (mask_shape(2,i) - mask_shape(1,i)+1)
00237       enddo
00238 !
00239 !  Allocate
00240 !
00241       Allocate (Masks(mask_id)%mask_array(1:len_alloc), STAT = ierror)
00242 
00243       if ( ierror > 0 ) then
00244          ierrp (1) = ierror
00245          ierrp (2) = len_alloc
00246 
00247          ierror = PRISM_Error_Alloc
00248          call psmile_error ( ierror, 'mask_array', &
00249                              ierrp, 2, __FILE__, __LINE__ )
00250          return
00251       endif
00252 !
00253 !     Copy data to internal mask_array
00254 !
00255       Masks(mask_id)%mask_array(1:len_alloc) = mask_array (1:len_alloc)
00256 !
00257 !     Save or initialize Type elements
00258 !
00259       Masks(mask_id)%mask_shape              = 1
00260       Masks(mask_id)%mask_shape(1:2,1:n_dim) = mask_shape(1:2,1:n_dim)
00261       Masks(mask_id)%grid_id                 = grid_id
00262       Masks(mask_id)%mask_type               = PSMILe_PointMask
00263 !
00264 !     Mark status of mask as changed
00265 !
00266       Masks(mask_id)%status = PSMILe_status_defined
00267 
00268 #ifdef VERBOSE
00269       print *, trim(ch_id), ': psmile_set_mask: eof ierror, mask_id =', ierror, mask_id
00270 
00271       call psmile_flushstd
00272 #endif /* VERBOSE */
00273 !
00274       end subroutine psmile_set_mask

Generated on 18 Mar 2011 for Oasis4 by  doxygen 1.6.1