psmile_get_mask_handle.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_Get_Mask_Handle
00008 !
00009 ! !INTERFACE:
00010 
00011       subroutine psmile_get_mask_handle (maskId, ierror)
00012 
00013 !
00014 ! !USES:
00015 !
00016       use PRISM_constants
00017       use PSMILe, dummy_interface => PSMILe_Get_mask_handle
00018       implicit none
00019 !
00020 ! !OUTPUT PARAMETERS:
00021 !
00022       integer,Intent(Out) :: maskId
00023 
00024 !     Handle to the mask information created.
00025 
00026       integer, Intent (Out)               :: ierror
00027 
00028 !     Returns the error code of PSMILe_Get_Mask_Handle;
00029 !             ierror = 0 : No error
00030 !             ierror > 0 : Severe error
00031 !
00032 !
00033 ! !LOCAL VARIABLES
00034 !
00035       integer             :: i, new_dim
00036 
00037       integer, parameter  :: nerrp = 2
00038       integer             :: ierrp (nerrp)
00039 
00040       type (Mask), Dimension(:), Pointer :: New_Masks
00041 !
00042 !
00043 ! !DESCRIPTION:
00044 !
00045 !  Subroutine "PSMILe_Get_mask_handle" returns a mask handle.
00046 !
00047 !
00048 ! !REVISION HISTORY:
00049 !
00050 !   Date      Programmer   Description
00051 ! ----------  ----------   -----------
00052 ! 01.12.03    R. Redler    created
00053 !
00054 !EOP
00055 !----------------------------------------------------------------------
00056 !
00057 !  $Id: psmile_get_mask_handle.F90 2325 2010-04-21 15:00:07Z valcke $
00058 !  $Autor$
00059 !
00060    Character(len=len_cvs_string), save :: mycvs = 
00061        '$Id: psmile_get_mask_handle.F90 2325 2010-04-21 15:00:07Z valcke $'
00062 !
00063 !----------------------------------------------------------------------
00064 !
00065 !   Initialization
00066 !
00067       ierror = 0
00068 !
00069 !   Search for a free mask index
00070 !
00071       do i = 1, Number_of_Masks_allocated
00072          if (Masks(i)%status == PSMILe_status_free) exit
00073       end do
00074 !
00075       if (i <= Number_of_Masks_allocated) then
00076          maskId = i
00077          Masks(maskId)%status = PSMILe_status_defined
00078          return
00079       endif
00080 !
00081 !   Allocate new Masks vector and copy old vector
00082 !
00083       new_dim = Number_of_Masks_allocated + 8
00084 !
00085       Allocate (New_Masks (new_dim), STAT = ierror)
00086       if (ierror > 0) then
00087          ierrp (1) = ierror
00088          ierrp (2) = new_dim
00089          call psmile_error ( PRISM_Error_Alloc, 'New_Masks', &
00090                              ierrp, 2, __FILE__, __LINE__ )
00091          return
00092       endif
00093 
00094       New_Masks(Number_of_Masks_allocated+1:new_dim)%status = &
00095          PSMILe_status_free
00096 !
00097       do i = Number_of_Masks_allocated+1, new_dim
00098         Nullify ( New_Masks(i)%mask_array )
00099       enddo
00100 !
00101       if (Number_of_Masks_allocated > 0) then
00102          New_Masks (1:Number_of_Masks_allocated) = &
00103              Masks (1:Number_of_Masks_allocated)
00104 !
00105 !   De-allocate Masks vector
00106 !
00107          Deallocate (Masks, STAT = ierror)
00108          if (ierror > 0) then
00109             ierrp (1) = ierror
00110             call psmile_error ( PRISM_Error_Dealloc, 'Masks', &
00111                                 ierrp, 1, __FILE__, __LINE__ )
00112             return
00113          endif
00114       endif
00115 !
00116 !   Update Number of Masks
00117 !
00118       Masks => New_Masks
00119 !
00120       maskId = Number_of_Masks_allocated + 1
00121       Number_of_Masks_allocated = new_dim
00122 !
00123       Masks(maskId)%status = PSMILe_status_defined
00124 !
00125       end subroutine PSMILe_Get_Mask_handle

Generated on 18 Mar 2011 for Oasis4 by  doxygen 1.6.1