psmile_reducedgrid_map.F90

Go to the documentation of this file.
00001 !-----------------------------------------------------------------------
00002 ! Copyright 2006-2010, CERFACS, Toulouse, France.
00003 ! Copyright 2006-2010, SGI Germany, Munich, Germany.
00004 ! Copyright 2006-2010, NEC Europe Ltd., London, UK.
00005 ! All rights reserved. Use is subject to OASIS4 license terms.
00006 !-----------------------------------------------------------------------
00007 !BOP
00008 !
00009 ! !ROUTINE: psmile_reducedgrid_map
00010 !
00011 ! !INTERFACE:
00012 
00013       subroutine psmile_reducedgrid_map ( grid_id, nbr_latitudes, &
00014                  nbr_points_per_lat, ierror )
00015 !
00016 ! !USES:
00017 !
00018       use PRISM
00019 !
00020       use PSMILe, dummy => psmile_reducedgrid_map
00021 !
00022       use psmile_grid_reduced_gauss, only : psmile_init_gauss_data
00023 
00024       implicit none
00025 !
00026 ! !INPUT PARAMETERS:
00027 !
00028       integer, Intent (In)                :: grid_id
00029 
00030 !     Specifies the handle to the grid information for which this call is.
00031 
00032       integer, Intent (In)                :: nbr_latitudes
00033 
00034 !     Number of latitudes of the Gaussian reduced grid.
00035 
00036       integer, Intent (In)                :: nbr_points_per_lat(nbr_latitudes)
00037 
00038 !     Number of longitude points per latitude circle
00039 !
00040       integer, Intent (Out)               :: ierror
00041 !
00042 !     Returns the error code of psmile_reducedgrid_map;
00043 !             ierror = 0 : No error
00044 !             ierror > 0 : Severe error
00045 !
00046 !
00047 ! !LOCAL VARIABLES
00048 !
00049       integer, parameter :: nerrp = 2
00050       integer            :: ierrp (nerrp)
00051 !
00052 ! !DESCRIPTION:
00053 !
00054 !  Subroutine "psmile_reducedgrid_map" receive mapping information for the
00055 !  1 d coordinate arrays.
00056 !
00057 !
00058 !
00059 ! !REVISION HISTORY:
00060 !   Date      Programmer   Description
00061 ! ----------  ----------   -----------
00062 ! 05.06.27    R. Redler    created
00063 !
00064 !EOP
00065 !----------------------------------------------------------------------
00066 !
00067 ! $Id: psmile_reducedgrid_map.F90 2917 2011-01-27 09:24:20Z hanke $
00068 ! $Author: hanke $
00069 !
00070   Character(len=len_cvs_string), save :: mycvs = 
00071       '$Id: psmile_reducedgrid_map.F90 2917 2011-01-27 09:24:20Z hanke $'
00072 !
00073 !----------------------------------------------------------------------
00074 
00075 #ifdef VERBOSE
00076       print *, trim(ch_id), ': psmile_reducedgrid_map: grid_id =', grid_id
00077 
00078       call psmile_flushstd
00079 #endif /* VERBOSE */
00080 
00081 !------------------------------------------------------------------------
00082 !  1st Initialization
00083 !------------------------------------------------------------------------
00084 
00085       ierror = 0
00086 
00087       call psmile_init_gauss_data (grid_id, nbr_points_per_lat)
00088 
00089 !------------------------------------------------------------------------
00090 !  2nd Control Grid id
00091 !------------------------------------------------------------------------
00092 
00093       if (grid_id < 1 .or. &
00094           grid_id > Number_of_Grids_allocated) then
00095 
00096          ierrp (1) = grid_id
00097          ierrp (2) = Number_of_Grids_allocated
00098 
00099          ierror = PRISM_Error_Arg
00100 
00101          call psmile_error ( ierror, 'grid_id', &
00102                              ierrp, 2, __FILE__, __LINE__ )
00103          return
00104       endif
00105 
00106       if (Grids(grid_id)%status == PSMILe_status_free) then
00107 
00108          ierrp (1) = grid_id
00109 
00110          ierror = PRISM_Error_Arg
00111 
00112          call psmile_error ( PRISM_Error_Arg, 'grid_id (not active)', &
00113                              ierrp, 1, __FILE__, __LINE__ )
00114          return
00115       endif
00116 
00117 !------------------------------------------------------------------------
00118 !  3rd Allocate memory
00119 !------------------------------------------------------------------------
00120 
00121       Allocate(Grids(grid_id)%nbr_points_per_lat(nbr_latitudes), STAT=ierror)
00122 
00123       if (ierror > 0) then
00124          ierrp (1) = ierror
00125          ierrp (2) = nbr_latitudes
00126          ierror = PRISM_Error_Alloc
00127 
00128          call psmile_error ( ierror, 'nbr_points_per_lat', &
00129                              ierrp, 2, __FILE__, __LINE__ )
00130          return
00131       endif
00132 
00133 !------------------------------------------------------------------------
00134 !  4th Store Data
00135 !------------------------------------------------------------------------
00136 
00137       Grids(grid_id)%nbr_latitudes = nbr_latitudes
00138       Grids(grid_id)%nbr_points_per_lat(:) = nbr_points_per_lat(:)
00139 
00140 !------------------------------------------------------------------------
00141 !  5th Epilogue
00142 !------------------------------------------------------------------------
00143 
00144 #ifdef VERBOSE
00145       print *, trim(ch_id), ': psmile_reducedgrid_map: eof ierror =', &
00146                ierror, '; grid_id', grid_id
00147 
00148       call psmile_flushstd
00149 #endif /* VERBOSE */
00150 !
00151       end subroutine psmile_reducedgrid_map

Generated on 18 Mar 2011 for Oasis4 by  doxygen 1.6.1