SURFEX v8.1
General documentation of Surfex
read_latlon.F90
Go to the documentation of this file.
1 !SFX_LIC Copyright 1994-2014 CNRS, Meteo-France and Universite Paul Sabatier
2 !SFX_LIC This is part of the SURFEX software governed by the CeCILL-C licence
3 !SFX_LIC version 1. See LICENSE, CeCILL-C_V1-en.txt and CeCILL-C_V1-fr.txt
4 !SFX_LIC for details. version 1.
5 ! #########
6  SUBROUTINE read_latlon (UG, U, USS, &
7  HPROGRAM,HSCHEME,HSUBROUTINE,HFILENAME)
8 ! #########################################################
9 !
10 !!**** *READ_LATLON* reads a latlon file and call treatment subroutine
11 !!
12 !! PURPOSE
13 !! -------
14 !!
15 !! METHOD
16 !! ------
17 !!
18 !! EXTERNAL
19 !! --------
20 !!
21 !! IMPLICIT ARGUMENTS
22 !! ------------------
23 !!
24 !! REFERENCE
25 !! ---------
26 !!
27 !! AUTHOR
28 !! ------
29 !!
30 !! V. Masson Meteo-France
31 !!
32 !! MODIFICATION
33 !! ------------
34 !!
35 !! Original 11/09/95
36 !!
37 !----------------------------------------------------------------------------
38 !
39 !* 0. DECLARATION
40 ! -----------
41 !
42 !
44 USE modd_surf_atm_n, ONLY : surf_atm_t
45 USE modd_sso_n, ONLY : sso_t
46 !
47 USE modd_pgd_grid, ONLY : xmeshlength
48 !
49 USE modi_get_luout
50 USE modi_open_namelist
51 USE modi_close_namelist
52 USE modi_open_file
53 USE modi_close_file
54 USE modi_readhead
55 USE modi_ini_ssowork
56 USE modi_pt_by_pt_treatment
57 !
58 !
59 USE yomhook ,ONLY : lhook, dr_hook
60 USE parkind1 ,ONLY : jprb
61 !
62 USE modi_abor1_sfx
63 !
64 IMPLICIT NONE
65 !
66 !* 0.1 Declaration of arguments
67 ! ------------------------
68 !
69 TYPE(surf_atm_grid_t), INTENT(INOUT) :: UG
70 TYPE(surf_atm_t), INTENT(INOUT) :: U
71 TYPE(sso_t), INTENT(INOUT) :: USS
72 !
73  CHARACTER(LEN=6), INTENT(IN) :: HPROGRAM ! Type of program
74  CHARACTER(LEN=6), INTENT(IN) :: HSCHEME ! Scheme treated
75  CHARACTER(LEN=6), INTENT(IN) :: HSUBROUTINE ! Name of the subroutine to call
76  CHARACTER(LEN=28), INTENT(IN) :: HFILENAME ! Name of the field file.
77 !
78 !* 0.2 Declaration of local variables read in the data file head
79 ! ---------------------------------------------------------
80 !
81 REAL :: ZGLBLATMIN ! minimum latitude of data box in the file
82 REAL :: ZGLBLONMIN ! minimum longitude of data box in the file
83 REAL :: ZGLBLATMAX ! maximum latitude of data box in the file
84 REAL :: ZGLBLONMAX ! maximum longitude of data box in the file
85 INTEGER :: INBLINE ! number of latitude rows (number of lines
86 INTEGER :: INBCOL ! number of longitude rows (number of columns)
87 REAL :: ZNODATA ! value below which data are not considered
88 !
89 !* 0.3 Declaration of local variables
90 ! ------------------------------
91 !
92 INTEGER :: IFILE ! logical units
93 INTEGER :: ILUOUT ! output listing logical unit
94 INTEGER :: IERR ! return codes
95 !
96 INTEGER :: JLOOP, IFACT ! loop index
97  CHARACTER(LEN=100):: YSTRING ! string
98 !
99 REAL :: ZDLAT ! latitude mesh in the data file
100 REAL :: ZDLON ! longitude mesh in the data file
101 INTEGER :: JLINE ! index of line
102 INTEGER :: JCOL ! index of column
103 !
104 REAL, DIMENSION(:), ALLOCATABLE :: ZVALUE ! value of a record of data points
105 REAL, DIMENSION(:), POINTER :: ZLAT ! latitude of data points
106 REAL, DIMENSION(:), POINTER :: ZLON ! longitude of data points
107 LOGICAL :: GCOMPRESS
108 !
109 REAL(KIND=JPRB) :: ZHOOK_HANDLE
110 !
111 !----------------------------------------------------------------------------
112 !
113 IF (lhook) CALL dr_hook('READ_LATLON',0,zhook_handle)
114  CALL get_luout(hprogram,iluout)
115 !
116 !* 1. Openning of header
117 ! ------------------
118 !
119  CALL open_namelist(hprogram,ifile,hfilename)
120 !
121 !----------------------------------------------------------------------------
122 !
123 !* 2. Reading of the global field
124 ! ---------------------------
125 !
126 !* 2.1 Head of data file
127 ! -----------------
128 !
129  CALL readhead(ifile,zglblatmin,zglblatmax,zglblonmin,zglblonmax, &
130  inbline,inbcol,znodata,zdlat,zdlon,zlat,zlon,ierr,ifact,&
131  gcompress)
132 IF (ierr/=0) THEN
133  CALL abor1_sfx('READ_LATLON: PROBLEM IN FILE HEADER')
134 END IF
135 !
136 !* 2.2 Closing of header
137 ! -----------------
138 !
139  CALL close_namelist(hprogram,ifile)
140 !
141 !----------------------------------------------------------------------------
142 !
143 !* 3. Adapt subgrid mesh to input file resolution
144 ! -------------------------------------------
145 !
146 IF (hsubroutine=='A_OROG') CALL ini_ssowork(xmeshlength,zdlat,zdlon)
147 !
148 !----------------------------------------------------------------------------
149 !
150 !* 4. Openning of file
151 ! ----------------
152 !
153  CALL open_file(hprogram,ifile,hfilename,'FORMATTED',haction='READ')
154 DO jloop=1,8
155  READ(ifile,fmt='(A100)') ystring
156 END DO
157 !
158 !----------------------------------------------------------------------------
159 !
160 !* 5. Allocation of array containing the data
161 ! ---------------------------------------
162 !
163 ALLOCATE(zvalue(inbcol))
164 !
165 !----------------------------------------------------------------------------
166 !
167 !* 6. Loop on lines
168 ! -------------
169 !
170 DO jline=1,inbline
171 !
172 !----------------------------------------------------------------------------
173 !
174 !* 7. Reading in the file
175 ! -------------------
176 !
177  READ(ifile,fmt=*) zvalue(:)
178 !
179 !
180 !----------------------------------------------------------------------------
181 !
182 !* 8. Loop on columns
183 ! ---------------
184 !
185  DO jcol=1,inbcol
186 
187 !-------------------------------------------------------------------------------
188 !
189 !* 9. value not valid
190 ! ---------------
191 !
192  IF (abs(zvalue(jcol)-znodata)<=1.e-10) cycle
193 !
194 !-------------------------------------------------------------------------------
195 !
196 !* 10. Call to the adequate subroutine (point by point treatment)
197 ! ----------------------------------------------------------
198 !
199  zvalue(:) = zvalue(:) / float(ifact)
200 
201  CALL pt_by_pt_treatment(ug, u, uss, &
202  iluout,zlat(jline:jline),zlon(jcol:jcol),zvalue(jcol:jcol),&
203  hsubroutine )
204 !
205 !-------------------------------------------------------------------------------
206  END DO
207 !-------------------------------------------------------------------------------
208 END DO
209 !-------------------------------------------------------------------------------
210 !
211 !* 11. deallocations
212 ! -------------
213 !
214 DEALLOCATE(zlat)
215 DEALLOCATE(zlon)
216 DEALLOCATE(zvalue)
217 !
218 !-------------------------------------------------------------------------------
219 !
220 !* 12. closes the file
221 ! ---------------
222 !
223  CALL close_file(hprogram,ifile)
224 IF (lhook) CALL dr_hook('READ_LATLON',1,zhook_handle)
225 !
226 !-------------------------------------------------------------------------------
227 !
228 END SUBROUTINE read_latlon
subroutine open_file(HPROGRAM, KUNIT, HFILE, HFORM, HACTION, HACCESS, KR
Definition: open_file.F90:7
subroutine abor1_sfx(YTEXT)
Definition: abor1_sfx.F90:7
subroutine read_latlon(UG, U, USS, HPROGRAM, HSCHEME, HSUBROUTINE, HFILENAME)
Definition: read_latlon.F90:8
integer, parameter jprb
Definition: parkind1.F90:32
subroutine ini_ssowork(PMESHLENGTH, PDLAT, PDLON)
Definition: ini_ssowork.F90:7
subroutine close_namelist(HPROGRAM, KLUNAM)
subroutine close_file(HPROGRAM, KUNIT)
Definition: close_file.F90:7
subroutine pt_by_pt_treatment(UG, U, USS, KLUOUT, PLAT, PLON, PVALUE, HSUBROUTINE
subroutine get_luout(HPROGRAM, KLUOUT)
Definition: get_luout.F90:7
logical lhook
Definition: yomhook.F90:15
subroutine open_namelist(HPROGRAM, KLUNAM, HFILE)
subroutine readhead(KGLB, PGLBLATMIN, PGLBLATMAX, PGLBLONMIN, PGLBLONM
Definition: readhead.F90:7