| SURFEX v7.3
   
    General documentation of Surfex | 
00001 !################################################################ 00002 SUBROUTINE READ_NAM_GRID_GAUSS(HPROGRAM,KGRID_PAR,KL,PGRID_PAR) 00003 !################################################################ 00004 ! 00005 !!**** *READ_NAM_GRID_GAUSS* - routine to read in namelist the horizontal grid 00006 !! 00007 !! PURPOSE 00008 !! ------- 00009 !! 00010 !!** METHOD 00011 !! ------ 00012 !! 00013 !! EXTERNAL 00014 !! -------- 00015 !! 00016 !! 00017 !! IMPLICIT ARGUMENTS 00018 !! ------------------ 00019 !! 00020 !! REFERENCE 00021 !! --------- 00022 !! 00023 !! 00024 !! AUTHOR 00025 !! ------ 00026 !! V. Masson *Meteo France* 00027 !! 00028 !! MODIFICATIONS 00029 !! ------------- 00030 !! Original 01/2004 00031 !! B. Decharme 2008 Comput and save the Mesh size 00032 !------------------------------------------------------------------------------- 00033 ! 00034 !* 0. DECLARATIONS 00035 ! ------------ 00036 ! 00037 USE MODD_CSTS, ONLY : XPI 00038 ! 00039 USE MODE_POS_SURF 00040 ! 00041 USE MODI_OPEN_NAMELIST 00042 USE MODI_CLOSE_NAMELIST 00043 USE MODI_GET_LUOUT 00044 ! 00045 USE MODE_GRIDTYPE_GAUSS 00046 ! 00047 USE EGGANGLES , ONLY : P_ASIN 00048 ! 00049 ! 00050 USE YOMHOOK ,ONLY : LHOOK, DR_HOOK 00051 USE PARKIND1 ,ONLY : JPRB 00052 ! 00053 USE MODI_ABOR1_SFX 00054 ! 00055 IMPLICIT NONE 00056 ! 00057 !* 0.1 Declarations of arguments 00058 ! ------------------------- 00059 ! 00060 CHARACTER(LEN=6), INTENT(IN) :: HPROGRAM ! calling program 00061 INTEGER, INTENT(INOUT) :: KGRID_PAR ! size of PGRID_PAR 00062 INTEGER, INTENT(OUT) :: KL ! number of points 00063 REAL, DIMENSION(KGRID_PAR), INTENT(OUT) :: PGRID_PAR ! parameters defining this grid 00064 ! 00065 !* 0.2 Declarations of local variables 00066 ! ------------------------------- 00067 ! 00068 INTEGER :: ILUOUT ! output listing logical unit 00069 INTEGER :: ILUNAM ! namelist file logical unit 00070 REAL, DIMENSION(:), ALLOCATABLE :: ZLAT_XY ! pseudo-latitudes 00071 REAL, DIMENSION(:), ALLOCATABLE :: ZLON_XY ! pseudo-longitudes 00072 REAL, DIMENSION(:), ALLOCATABLE :: ZLAT ! latitudes 00073 REAL, DIMENSION(:), ALLOCATABLE :: ZLON ! longitudes 00074 REAL, DIMENSION(:), ALLOCATABLE :: ZMESH_SIZE ! Mesh size 00075 ! 00076 !* 0.3 Declarations of namelist 00077 ! ------------------------ 00078 ! 00079 INTEGER :: NDGLG ! number of pseudo-latitudes 00080 REAL :: RMUCEN ! sine of the latitude of the rotated pole 00081 REAL :: RLOCEN ! longitude of the rotated pole (radian) 00082 REAL :: RSTRET ! stretching factor (must be greater than or equal to 1) 00083 INTEGER, DIMENSION(1000) :: NRGRI ! number of pseudo-longitudes on each 00084 ! pseudo-latitude circle on pseau 00085 ! northern hemisphere (starting from 00086 ! the rotated pole) 00087 ! 00088 REAL :: ZLAPO ! latitude of the rotated pole (deg) 00089 REAL :: ZLOPO ! longitude of the rotated pole (deg) 00090 REAL :: ZCODIL ! stretching factor (must be greater than or equal to 1) 00091 INTEGER :: ITYP ! type of transform (0 --> no rotation, 1 otherwise) 00092 INTEGER :: INLATI ! number of latitudes 00093 INTEGER, DIMENSION(:), ALLOCATABLE :: INLOPA ! number of pseudo-longitudes on each 00094 ! pseudo-latitude circle 00095 00096 INTEGER :: JSTGLO 00097 00098 ! 00099 REAL, DIMENSION(:), POINTER :: ZGRID_PAR 00100 ! 00101 LOGICAL :: GFOUND 00102 REAL(KIND=JPRB) :: ZHOOK_HANDLE 00103 ! 00104 NAMELIST/NAMDIM/NDGLG 00105 NAMELIST/NAMGEM/RMUCEN, RLOCEN, RSTRET 00106 NAMELIST/NAMRGRI/NRGRI 00107 ! 00108 !------------------------------------------------------------------------------ 00109 ! 00110 !* 1. Default values 00111 ! 00112 IF (LHOOK) CALL DR_HOOK('READ_NAM_GRID_GAUSS',0,ZHOOK_HANDLE) 00113 NDGLG = 0 00114 RMUCEN = 1. 00115 RLOCEN = XPI 00116 RSTRET = 1. 00117 ! 00118 NRGRI(:) = 0 00119 !------------------------------------------------------------------------------ 00120 ! 00121 !* 2. opening of namelist 00122 ! 00123 CALL GET_LUOUT(HPROGRAM,ILUOUT) 00124 ! 00125 CALL OPEN_NAMELIST(HPROGRAM,ILUNAM) 00126 ! 00127 !--------------------------------------------------------------------------- 00128 ! 00129 !* 3. Reading of projection parameters 00130 ! -------------------------------- 00131 ! 00132 CALL POSNAM(ILUNAM,'NAMGEM',GFOUND,ILUOUT) 00133 IF (GFOUND) READ(UNIT=ILUNAM,NML=NAMGEM) 00134 ! 00135 IF (RSTRET<1.) THEN 00136 WRITE(ILUOUT,*) '****************************************************' 00137 WRITE(ILUOUT,*) 'stretching factor RSTRET for the Gaussian grid' 00138 WRITE(ILUOUT,*) 'definition must be greater than or equal to 1' 00139 WRITE(ILUOUT,*) 'You have set RSTRET=', RSTRET 00140 WRITE(ILUOUT,*) 'Please modify the value of RSTRET in namelist NAMGEM' 00141 WRITE(ILUOUT,*) '****************************************************' 00142 CALL ABOR1_SFX('READ_NAM_GRID_GAUSS: STRETCHING FACTOR MUST BE >= 1.') 00143 END IF 00144 ! 00145 ZLAPO = 180. / XPI * P_ASIN(RMUCEN) 00146 ZLOPO = 180. / XPI * RLOCEN 00147 ! 00148 ZCODIL = RSTRET 00149 ! 00150 !--------------------------------------------------------------------------- 00151 ! 00152 !* 4. Reading parameters of the grid 00153 ! ------------------------------ 00154 ! 00155 CALL POSNAM(ILUNAM,'NAMDIM',GFOUND,ILUOUT) 00156 IF (GFOUND) READ(UNIT=ILUNAM,NML=NAMDIM) 00157 CALL POSNAM(ILUNAM,'NAMRGRI',GFOUND,ILUOUT) 00158 IF (GFOUND) READ(UNIT=ILUNAM,NML=NAMRGRI) 00159 ! 00160 INLATI = NDGLG 00161 ALLOCATE(INLOPA(INLATI)) 00162 INLOPA(1:INLATI/2) = NRGRI(1:INLATI/2) 00163 INLOPA(INLATI/2+1:INLATI) = NRGRI(INLATI/2:1:-1) 00164 ! 00165 !--------------------------------------------------------------------------- 00166 ! 00167 !* 5. Computes pseudo-latitudes and pseudo-longitudes of all points 00168 ! ------------------------------------------------------------- 00169 ! 00170 !* number of points 00171 KL = SUM(INLOPA) 00172 00173 ! 00174 !* type of transform 00175 IF (ZLAPO>89.99 .AND. ABS(ZLOPO)<0.00001) THEN 00176 ITYP=0 00177 ELSE 00178 ITYP=1 00179 ENDIF 00180 ! 00181 ALLOCATE(ZLAT_XY(KL)) 00182 ALLOCATE(ZLON_XY(KL)) 00183 00184 CALL COMP_GRIDTYPE_GAUSS(INLATI,INLOPA,KL,ITYP,ZLAT_XY,ZLON_XY) 00185 00186 ! 00187 !--------------------------------------------------------------------------- 00188 ! 00189 !* 6. Computes latitudes and longitudes 00190 ! --------------------------------- 00191 ! 00192 !* all points are used 00193 ALLOCATE(ZLAT(KL)) 00194 ALLOCATE(ZLON(KL)) 00195 CALL LATLON_GAUSS(ZLON_XY,ZLAT_XY,KL,ZLOPO,ZLAPO,ZCODIL,ZLON,ZLAT) 00196 ! 00197 !--------------------------------------------------------------------------- 00198 ! 00199 !* 7. Computes mesh size 00200 ! --------------------------------- 00201 ! 00202 ALLOCATE(ZMESH_SIZE(KL)) 00203 ! 00204 CALL MESH_SIZE_GAUSS(KL,INLATI,INLOPA,ZLAPO,ZLOPO,ZCODIL,& 00205 ZLAT_XY,ZLON,ZLAT,ZMESH_SIZE) 00206 ! 00207 !--------------------------------------------------------------------------- 00208 CALL CLOSE_NAMELIST(HPROGRAM,ILUNAM) 00209 !--------------------------------------------------------------------------- 00210 ! 00211 !* 8. All this information stored into pointer PGRID_PAR 00212 ! -------------------------------------------------- 00213 ! 00214 CALL PUT_GRIDTYPE_GAUSS(ZGRID_PAR,INLATI,ZLAPO,ZLOPO,ZCODIL,INLOPA, & 00215 KL,ZLAT,ZLON,ZLAT_XY,ZLON_XY,ZMESH_SIZE ) 00216 ! 00217 DEALLOCATE(ZLAT) 00218 DEALLOCATE(ZLON) 00219 DEALLOCATE(ZLAT_XY) 00220 DEALLOCATE(ZLON_XY) 00221 DEALLOCATE(INLOPA) 00222 DEALLOCATE(ZMESH_SIZE) 00223 !--------------------------------------------------------------------------- 00224 ! 00225 !* 1st call : initializes dimension 00226 ! 00227 IF (KGRID_PAR==0) THEN 00228 KGRID_PAR = SIZE(ZGRID_PAR) 00229 ! 00230 ELSE 00231 ! 00232 !* 2nd call : initializes grid array 00233 ! 00234 PGRID_PAR(:) = 0. 00235 PGRID_PAR(:) = ZGRID_PAR 00236 END IF 00237 ! 00238 DEALLOCATE(ZGRID_PAR) 00239 IF (LHOOK) CALL DR_HOOK('READ_NAM_GRID_GAUSS',1,ZHOOK_HANDLE) 00240 00241 ! 00242 !--------------------------------------------------------------------------- 00243 ! 00244 END SUBROUTINE READ_NAM_GRID_GAUSS
 1.8.0
 1.8.0