SURFEX  V8_0
Surfex V8_0 release
 All Classes Files Functions Variables
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 (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 !
43 !
45 !
46 USE modd_pgd_grid, ONLY : xmeshlength
47 !
48 USE modi_get_luout
49 USE modi_open_namelist
50 USE modi_close_namelist
51 USE modi_open_file
52 USE modi_close_file
53 USE modi_readhead
54 USE modi_ini_ssowork
55 USE modi_pt_by_pt_treatment
56 !
57 !
58 USE yomhook ,ONLY : lhook, dr_hook
59 USE parkind1 ,ONLY : jprb
60 !
61 USE modi_abor1_sfx
62 !
63 IMPLICIT NONE
64 !
65 !* 0.1 Declaration of arguments
66 ! ------------------------
67 !
68 !
69 TYPE(surf_atm_sso_t), INTENT(INOUT) :: uss
70 !
71  CHARACTER(LEN=6), INTENT(IN) :: hprogram ! Type of program
72  CHARACTER(LEN=6), INTENT(IN) :: hscheme ! Scheme treated
73  CHARACTER(LEN=6), INTENT(IN) :: hsubroutine ! Name of the subroutine to call
74  CHARACTER(LEN=28), INTENT(IN) :: hfilename ! Name of the field file.
75 !
76 !* 0.2 Declaration of local variables read in the data file head
77 ! ---------------------------------------------------------
78 !
79 REAL :: zglblatmin ! minimum latitude of data box in the file
80 REAL :: zglblonmin ! minimum longitude of data box in the file
81 REAL :: zglblatmax ! maximum latitude of data box in the file
82 REAL :: zglblonmax ! maximum longitude of data box in the file
83 INTEGER :: inbline ! number of latitude rows (number of lines
84 INTEGER :: inbcol ! number of longitude rows (number of columns)
85 REAL :: znodata ! value below which data are not considered
86 !
87 !* 0.3 Declaration of local variables
88 ! ------------------------------
89 !
90 INTEGER :: ifile ! logical units
91 INTEGER :: iluout ! output listing logical unit
92 INTEGER :: ierr ! return codes
93 !
94 INTEGER :: jloop ! loop index
95  CHARACTER(LEN=100):: ystring ! string
96 !
97 REAL :: zdlat ! latitude mesh in the data file
98 REAL :: zdlon ! longitude mesh in the data file
99 INTEGER :: jline ! index of line
100 INTEGER :: jcol ! index of column
101 !
102 REAL, DIMENSION(:), ALLOCATABLE :: zvalue ! value of a record of data points
103 REAL, DIMENSION(:), POINTER :: zlat ! latitude of data points
104 REAL, DIMENSION(:), POINTER :: zlon ! longitude of data points
105 REAL(KIND=JPRB) :: zhook_handle
106 !
107 !
108 !----------------------------------------------------------------------------
109 !
110 IF (lhook) CALL dr_hook('READ_LATLON',0,zhook_handle)
111  CALL get_luout(hprogram,iluout)
112 !
113 !* 1. Openning of header
114 ! ------------------
115 !
116  CALL open_namelist(hprogram,ifile,hfilename)
117 !
118 !----------------------------------------------------------------------------
119 !
120 !* 2. Reading of the global field
121 ! ---------------------------
122 !
123 !* 2.1 Head of data file
124 ! -----------------
125 !
126  CALL readhead(ifile,zglblatmin,zglblatmax,zglblonmin,zglblonmax, &
127  inbline,inbcol,znodata,zdlat,zdlon,zlat,zlon,ierr)
128 IF (ierr/=0) THEN
129  CALL abor1_sfx('READ_LATLON: PROBLEM IN FILE HEADER')
130 END IF
131 !
132 !* 2.2 Closing of header
133 ! -----------------
134 !
135  CALL close_namelist(hprogram,ifile)
136 !
137 !----------------------------------------------------------------------------
138 !
139 !* 3. Adapt subgrid mesh to input file resolution
140 ! -------------------------------------------
141 !
142 IF (hsubroutine=='A_OROG') CALL ini_ssowork(xmeshlength,zdlat,zdlon)
143 !
144 !----------------------------------------------------------------------------
145 !
146 !* 4. Openning of file
147 ! ----------------
148 !
149  CALL open_file(hprogram,ifile,hfilename,'FORMATTED',haction='READ')
150 DO jloop=1,8
151  READ(ifile,fmt='(A100)') ystring
152 END DO
153 !
154 !----------------------------------------------------------------------------
155 !
156 !* 5. Allocation of array containing the data
157 ! ---------------------------------------
158 !
159 ALLOCATE(zvalue(inbcol))
160 !
161 !----------------------------------------------------------------------------
162 !
163 !* 6. Loop on lines
164 ! -------------
165 !
166 DO jline=1,inbline
167 !
168 !----------------------------------------------------------------------------
169 !
170 !* 7. Reading in the file
171 ! -------------------
172 !
173  READ(ifile,fmt=*) zvalue(:)
174 !
175 !
176 !----------------------------------------------------------------------------
177 !
178 !* 8. Loop on columns
179 ! ---------------
180 !
181  DO jcol=1,inbcol
182 
183 !-------------------------------------------------------------------------------
184 !
185 !* 9. value not valid
186 ! ---------------
187 !
188  IF (abs(zvalue(jcol)-znodata)<=1.e-10) cycle
189 !
190 !-------------------------------------------------------------------------------
191 !
192 !* 10. Call to the adequate subroutine (point by point treatment)
193 ! ----------------------------------------------------------
194 !
195  CALL pt_by_pt_treatment(uss, &
196  iluout,zlat(jline:jline),zlon(jcol:jcol),zvalue(jcol:jcol),&
197  hsubroutine )
198 !
199 !-------------------------------------------------------------------------------
200  END DO
201 !-------------------------------------------------------------------------------
202 END DO
203 !-------------------------------------------------------------------------------
204 !
205 !* 11. deallocations
206 ! -------------
207 !
208 DEALLOCATE(zlat)
209 DEALLOCATE(zlon)
210 DEALLOCATE(zvalue)
211 !
212 !-------------------------------------------------------------------------------
213 !
214 !* 12. closes the file
215 ! ---------------
216 !
217  CALL close_file(hprogram,ifile)
218 IF (lhook) CALL dr_hook('READ_LATLON',1,zhook_handle)
219 !
220 !-------------------------------------------------------------------------------
221 !
222 END SUBROUTINE read_latlon
subroutine read_latlon(USS, HPROGRAM, HSCHEME, HSUBROUTINE, HFILENAME)
Definition: read_latlon.F90:6
subroutine pt_by_pt_treatment(USS, KLUOUT, PLAT, PLON, PVALUE, HSUBROUTINE, KNBLINES, PNODATA)
subroutine abor1_sfx(YTEXT)
Definition: abor1_sfx.F90:6
subroutine readhead(KGLB, PGLBLATMIN, PGLBLATMAX, PGLBLONMIN, PGLBLONMAX, KNBLAT, KNBLON, PCUTVAL, PDLAT, PDLON, PLAT, PLON, KERR)
Definition: readhead.F90:6
subroutine ini_ssowork(PMESHLENGTH, PDLAT, PDLON)
Definition: ini_ssowork.F90:6
subroutine close_namelist(HPROGRAM, KLUNAM)
subroutine close_file(HPROGRAM, KUNIT)
Definition: close_file.F90:6
subroutine get_luout(HPROGRAM, KLUOUT)
Definition: get_luout.F90:6
subroutine open_namelist(HPROGRAM, KLUNAM, HFILE)
subroutine open_file(HPROGRAM, KUNIT, HFILE, HFORM, HACTION, HACCESS, KRECL)
Definition: open_file.F90:6