psmile_get_grid_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_grid_handle
00008 !
00009 ! !INTERFACE:
00010 
00011       subroutine psmile_get_grid_handle (grid_id, ierror)
00012 !
00013 ! !USES:
00014 !
00015       use PRISM_constants
00016       use PSMILe, dummy_interface => PSMILe_Get_grid_handle
00017 
00018       implicit none
00019 !
00020 ! !OUTPUT PARAMETERS:
00021 !
00022       integer, Intent (Out)               :: grid_id
00023 
00024 !     Returns the handle to the grid information created.
00025 
00026       integer, Intent (Out)               :: ierror
00027 
00028 !     Returns the error code of PSMILe_Get_grid_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 (Grid), Dimension(:), Pointer :: New_Grids
00041 !
00042 !
00043 ! !DESCRIPTION:
00044 !
00045 ! Subroutine "PSMILe_Get_grid_handle" returns a grid handle.
00046 !
00047 !
00048 ! !REVISION HISTORY:
00049 !
00050 !   Date      Programmer    Description
00051 ! ----------  -----------   -----------
00052 ! 01.12.03    H. Ritzdorf   created
00053 !
00054 !EOP
00055 !----------------------------------------------------------------------
00056 !
00057 ! $Id: psmile_get_grid_handle.F90 2325 2010-04-21 15:00:07Z valcke $
00058 ! $Author: valcke $
00059 !
00060    Character(len=len_cvs_string), save :: mycvs = 
00061        '$Id: psmile_get_grid_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 grid index
00070 !
00071          do i = 1, Number_of_Grids_allocated
00072          if (Grids(i)%status == PSMILe_status_free) exit
00073          end do
00074 !
00075       if (i <= Number_of_Grids_allocated) then
00076          grid_id = i
00077          Grids(grid_id)%status = PSMILe_status_defined
00078          return
00079       endif
00080 !
00081 !   Allocate new Grids vector, initialize and copy old vector
00082 !
00083       new_dim = Number_of_Grids_allocated + 8
00084 !
00085       Allocate (New_Grids (new_dim), STAT = ierror)
00086       if (ierror > 0) then
00087          ierrp (1) = ierror
00088          ierrp (2) = new_dim
00089          ierror = PRISM_Error_Alloc
00090 
00091          call psmile_error ( ierror, 'New_Grids', &
00092                              ierrp, 2, __FILE__, __LINE__ )
00093          return
00094       endif
00095 !
00096       New_Grids (1:Number_of_Grids_allocated) = &
00097           Grids (1:Number_of_Grids_allocated)
00098 
00099       New_Grids(Number_of_Grids_allocated+1:new_dim)%status = &
00100           PSMILe_status_free
00101 
00102       New_Grids(Number_of_Grids_allocated+1:new_dim)%nlev = 0
00103 
00104       do i = Number_of_Grids_allocated+1, new_dim
00105         New_Grids(i)%len_periodic (:) = 0
00106         New_Grids(i)%nbr_halo_segments = 0
00107       end do
00108 
00109       do i = Number_of_Grids_allocated+1, new_dim
00110         Nullify ( New_Grids(i)%corner_pointer )
00111         Nullify ( New_Grids(i)%partition )
00112         Nullify ( New_Grids(i)%extent )
00113         Nullify ( New_Grids(i)%nbr_points_per_lat )
00114         Nullify ( New_Grids(i)%mg_infos )
00115         Nullify ( New_Grids(i)%send_list )
00116         Nullify ( New_Grids(i)%recv_list )
00117         Nullify ( New_Grids(i)%get_list )
00118         Nullify ( New_Grids(i)%put_list )
00119         Nullify ( New_Grids(i)%remote_index )
00120         Nullify ( New_Grids(i)%star )
00121         Nullify ( New_Grids(i)%face )
00122         Nullify ( New_Grids(i)%global_beg )
00123         Nullify ( New_Grids(i)%global_end )
00124         Nullify ( New_Grids(i)%l2g )
00125         Nullify ( New_Grids(i)%g2l )
00126         Nullify ( New_Grids(i)%halo )
00127       enddo
00128 !
00129 !   De-allocate Grids vector
00130 !
00131       Deallocate (Grids, STAT = ierror)
00132       if (ierror > 0) then
00133          ierrp (1) = ierror
00134          ierror = PRISM_Error_Dealloc
00135 
00136          call psmile_error ( ierror, 'Grids', &
00137                              ierrp, 1, __FILE__, __LINE__ )
00138          return
00139       endif
00140 !
00141 !   Update Number of Grids
00142 !
00143       Grids => New_Grids
00144 !
00145       grid_id = Number_of_Grids_allocated + 1
00146       Number_of_Grids_allocated = new_dim
00147 !
00148       Grids(grid_id)%status = PSMILe_status_defined
00149 !
00150       end subroutine PSMILe_Get_grid_handle

Generated on 18 Mar 2011 for Oasis4 by  doxygen 1.6.1