SURFEX  V8_0
Surfex V8_0 release
 All Classes Files Functions Variables
mode_read_netcdf_mercator.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 !! Modified 09/2013 : S. Senesi : adapt READ_NETCDF_SST to read 2D fields other than SST
7 !!
9 !!!=============================================================================
10 !-------------------------------------------------------------------------------
11 !
12 !
13 USE modi_abor1_sfx
14 !
15 USE yomhook ,ONLY : lhook, dr_hook
16 USE parkind1 ,ONLY : jprb
17 !
18  CONTAINS
19 !-------------------------------------------------------------------
20 !-------------------------------------------------------------------
21 ! ####################
22  SUBROUTINE handle_err_mer(status,line)
23 ! ####################
24 IMPLICIT NONE
25 INTEGER, INTENT(IN) :: status
26  CHARACTER(LEN=80), INTENT(IN) :: line
27 REAL(KIND=JPRB) :: zhook_handle
28 !
29 include 'netcdf.inc'
30 !
31 IF (lhook) CALL dr_hook('MODE_READ_NETCDF_MERCATOR:HANDLE_ERR_MER',0,zhook_handle)
32 IF (status /= nf_noerr) THEN
33  CALL abor1_sfx('MODE_READ_NETCDF_MERCATOR: HANDLE_ERR_MER')
34 END IF
35 IF (lhook) CALL dr_hook('MODE_READ_NETCDF_MERCATOR:HANDLE_ERR_MER',1,zhook_handle)
36 END SUBROUTINE handle_err_mer
37 !-------------------------------------------------------------------
38 !-------------------------------------------------------------------
39 ! ####################
40  SUBROUTINE get1dcdf(KCDF_ID,IDVAR,PMISSVALUE,PVALU1D)
41 ! ####################
42 !
43 IMPLICIT NONE
44 !
45 INTEGER,INTENT(IN) :: kcdf_id !netcdf file identifiant
46 INTEGER,INTENT(IN) :: idvar !variable to read identifiant
47 REAL, INTENT(OUT) :: pmissvalue !undefined value
48 REAL,DIMENSION(:),INTENT(OUT) :: pvalu1d !value array
49 !
50 integer :: status
51 character(len=80) :: haction
52 integer,save :: ndims=1
53 integer :: kvartype
54 integer,DIMENSION(:),ALLOCATABLE :: nvardimid,nvardimlen
55 character(len=80),DIMENSION(:),ALLOCATABLE :: nvardimnam
56 integer :: jloop
57 integer :: ngatts
58 character(len=80),DIMENSION(:),ALLOCATABLE :: hname
59 REAL,DIMENSION(:),ALLOCATABLE :: zvalu1d !value array
60 REAL(KIND=JPRB) :: zhook_handle
61 !
62 include 'netcdf.inc'
63 !
64 IF (lhook) CALL dr_hook('MODE_READ_NETCDF_MERCATOR:GET1DCDF',0,zhook_handle)
65 pmissvalue=-9999.9
66 ALLOCATE(nvardimid(ndims))
67 ALLOCATE(nvardimlen(ndims))
68 ALLOCATE(nvardimnam(ndims))
69 nvardimid(:)=0
70 nvardimlen(:)=0
71 nvardimnam(:)=' '
72 !
73 haction='get variable type'
74 status=nf_inq_vartype(kcdf_id,idvar,kvartype)
75 if (status/=nf_noerr) CALL handle_err_mer(status,haction)
76 !write(0,*) 'variable type = ',KVARTYPE
77 !
78 haction='get variable dimensions name'
79 status=nf_inq_dimname(kcdf_id,idvar,nvardimnam(ndims))
80 if (status/=nf_noerr) CALL handle_err_mer(status,haction)
81 !
82 haction='get variable dimensions length'
83 status=nf_inq_dimlen(kcdf_id,idvar,nvardimlen(ndims))
84 if (status/=nf_noerr) CALL handle_err_mer(status,haction)
85 !write(0,*) 'variable dimension ',NDIMS,' named ',NVARDIMNAM(NDIMS),&
86 ! &'has a length of',NVARDIMLEN(NDIMS)
87 !!
88 haction='get attributs'
89 status=nf_inq_varnatts(kcdf_id,idvar,ngatts)
90 if (status/=nf_noerr) CALL handle_err_mer(status,haction)
91 !write(0,*) 'number of attributes = ',NGATTS
92 allocate(hname(1:ngatts))
93 !
94 ALLOCATE(zvalu1d(1:nvardimlen(ndims)))
95 zvalu1d=0.
96 !
97 IF (kvartype>=5) then
98  haction='get variable values (1D)'
99  status=nf_get_var_double(kcdf_id,idvar,zvalu1d(:))
100  if (status/=nf_noerr) CALL handle_err_mer(status,haction)
101 ENDIF
102 !
103 pvalu1d(:)=zvalu1d(:)
104 !
105 IF (ALLOCATED(zvalu1d )) DEALLOCATE(zvalu1d)
106 IF (lhook) CALL dr_hook('MODE_READ_NETCDF_MERCATOR:GET1DCDF',1,zhook_handle)
107 !
108 END SUBROUTINE get1dcdf
109 !-------------------------------------------------------------------
110 !-------------------------------------------------------------------
111 ! ####################
112  SUBROUTINE get2dcdf(KCDF_ID,IDVAR,PDIM1,HDIM1NAME,PDIM2,HDIM2NAME,&
113  pmissvalue,pvalu2d)
114 ! ####################
115 USE modd_surf_par, ONLY : xundef
116 !
117 IMPLICIT NONE
118 !
119 INTEGER,INTENT(IN) :: kcdf_id !netcdf file identifiant
120 INTEGER,INTENT(IN) :: idvar !variable to read identifiant
121 REAL,DIMENSION(:),INTENT(OUT) :: pdim1,pdim2 !dimensions for PVALU2D array
122  CHARACTER(len=80),INTENT(OUT) :: hdim1name,hdim2name !dimensions names
123 REAL, INTENT(OUT) :: pmissvalue
124 REAL,DIMENSION(:,:),INTENT(OUT) :: pvalu2d !value array
125 !
126 integer :: status
127 character(len=80) :: haction
128 integer,save :: ndims=2
129 integer :: kvartype
130 integer,DIMENSION(:),ALLOCATABLE :: nvardimid,nvardimlen
131 character(len=80),DIMENSION(:),ALLOCATABLE :: nvardimnam
132 integer :: jloop2, jloop, j1, j2
133 integer :: ngatts
134 character(len=80),DIMENSION(:),ALLOCATABLE :: hname
135 real :: zmiss1,zmiss2
136 real :: zscfa, zoffs
137 REAL,DIMENSION(:,:),ALLOCATABLE :: zvalu2d !value array
138 INTEGER,DIMENSION(:,:),ALLOCATABLE :: ivalu2d
139 REAL(KIND=JPRB) :: zhook_handle
140 !
141 include 'netcdf.inc'
142 !
143 IF (lhook) CALL dr_hook('MODE_READ_NETCDF_MERCATOR:GET2DCDF',0,zhook_handle)
144 pmissvalue=-9999.9
145 ALLOCATE(nvardimid(ndims))
146 ALLOCATE(nvardimlen(ndims))
147 ALLOCATE(nvardimnam(ndims))
148 nvardimid(:)=0
149 nvardimlen(:)=0
150 nvardimnam(:)=' '
151 !
152 haction='get variable type'
153 status=nf_inq_vartype(kcdf_id,idvar,kvartype)
154 if (status/=nf_noerr) CALL handle_err_mer(status,haction)
155 !write(0,*) 'variable type = ',KVARTYPE
156 !
157 haction='get variable dimensions identifiant'
158 status=nf_inq_vardimid(kcdf_id,idvar,nvardimid)
159 if (status/=nf_noerr) CALL handle_err_mer(status,haction)
160 !
161 haction='get attributs'
162 status=nf_inq_varnatts(kcdf_id,idvar,ngatts)
163 if (status/=nf_noerr) CALL handle_err_mer(status,haction)
164 !write(0,*) 'number of attributes = ',NGATTS
165 allocate(hname(1:ngatts))
166 !
167 zscfa=1.
168 zoffs=0.
169 DO jloop=1,ngatts
170  status=nf_inq_attname(kcdf_id,idvar,jloop,hname(jloop))
171  if (status/=nf_noerr) CALL handle_err_mer(status,haction)
172  !write(0,*) 'attributes names = ', hname(JLOOP)
173  if (trim(hname(jloop))=='missing_value') then
174  !write(0,*) 'missing value search '
175  haction='get missing value'
176  status=nf_get_att_double(kcdf_id,idvar,"missing_value",pmissvalue)
177  if (status/=nf_noerr) CALL handle_err_mer(status,haction)
178  !write(0,*) 'missing value = ',PMISSVALUE
179  else
180  if (trim(hname(jloop))=='_FillValue') then
181  !write(0,*) 'missing value found '
182  haction='get _FillValue'
183  status=nf_get_att_double(kcdf_id,idvar,"_FillValue",pmissvalue)
184  if (status/=nf_noerr) CALL handle_err_mer(status,haction)
185  !write(0,*) 'missing value = ',PMISSVALUE
186  endif
187  endif
188  if (trim(hname(jloop))=='scale_factor') then
189  !write(0,*) 'missing value found '
190  haction='get scale factor'
191  status=nf_get_att_double(kcdf_id,idvar,"scale_factor",zscfa)
192  if (status/=nf_noerr) CALL handle_err_mer(status,haction)
193  !write(0,*) 'missing value = ',PMISSVALUE
194  endif
195  if (trim(hname(jloop))=='add_offset') then
196  !write(0,*) 'missing value found '
197  haction='get offset'
198  status=nf_get_att_double(kcdf_id,idvar,"add_offset",zoffs)
199  if (status/=nf_noerr) CALL handle_err_mer(status,haction)
200  !write(0,*) 'missing value = ',PMISSVALUE
201  endif
202 ENDDO
203 !
204 !
205 DO jloop2=1,ndims
206  haction='get variable dimensions name'
207  status=nf_inq_dimname(kcdf_id,nvardimid(jloop2),nvardimnam(jloop2))
208  if (status/=nf_noerr) CALL handle_err_mer(status,haction)
209  haction='get variable dimensions length'
210  status=nf_inq_dimlen(kcdf_id,nvardimid(jloop2),nvardimlen(jloop2))
211  if (status/=nf_noerr) CALL handle_err_mer(status,haction)
212  !write(0,*) 'variable dimension ',JLOOP2,' named ',NVARDIMNAM(JLOOP2),&
213  ! &'has a length of',NVARDIMLEN(JLOOP2)
214 ENDDO
215 !
216 IF (kvartype>=5) then
217  ALLOCATE(zvalu2d(1:nvardimlen(1),1:nvardimlen(2)))
218  zvalu2d=0.
219  haction='get variable values (2D)'
220  status=nf_get_var_double(kcdf_id,idvar,zvalu2d(:,:))
221  if (status/=nf_noerr) CALL handle_err_mer(status,haction)
222 ELSE
223  ALLOCATE(ivalu2d(1:nvardimlen(1),1:nvardimlen(2)))
224  ivalu2d=0.
225  haction='get variable values (2D)'
226  status=nf_get_var_int(kcdf_id,idvar,ivalu2d(:,:))
227  if (status/=nf_noerr) CALL handle_err_mer(status,haction)
228 ENDIF
229 pvalu2d(:,:)=xundef
230 DO j1=1,nvardimlen(1)
231  DO j2=1,nvardimlen(2)
232  IF (kvartype>=5) THEN
233  IF (zvalu2d(j1,j2)/=pmissvalue) pvalu2d(j1,j2)=zvalu2d(j1,j2)*zscfa+zoffs
234  ELSE
235  IF (zvalu2d(j1,j2)/=pmissvalue) pvalu2d(j1,j2)=ivalu2d(j1,j2)*zscfa+zoffs
236  ENDIF
237  ENDDO
238 ENDDO
239 !
240  CALL get1dcdf(kcdf_id,nvardimid(1),zmiss1,pdim1)
241  CALL get1dcdf(kcdf_id,nvardimid(2),zmiss2,pdim2)
242 hdim1name=nvardimnam(1)
243 hdim2name=nvardimnam(2)
244 IF (ALLOCATED(zvalu2d )) DEALLOCATE(zvalu2d)
245 IF (lhook) CALL dr_hook('MODE_READ_NETCDF_MERCATOR:GET2DCDF',1,zhook_handle)
246 !
247 END SUBROUTINE get2dcdf
248 !-------------------------------------------------------------------
249 !-------------------------------------------------------------------
250 ! ####################
251  SUBROUTINE get3dcdf(KCDF_ID,IDVAR,PDIM1,HDIM1NAME,PDIM2,HDIM2NAME,&
252  pdim3,hdim3name,pmissvalue,pvalu3d)
253 ! ####################
254 USE modd_surf_par, ONLY : xundef
255 !
256 IMPLICIT NONE
257 !
258 INTEGER,INTENT(IN) :: kcdf_id !netcdf file identifiant
259 INTEGER,INTENT(IN) :: idvar !variable to read identifiant
260 REAL,DIMENSION(:),INTENT(OUT) :: pdim1,pdim2,pdim3 !dimensions for PVALU2D array
261  CHARACTER(len=80),INTENT(OUT) :: hdim1name,hdim2name,hdim3name !dimensions names
262 REAL, INTENT(OUT) :: pmissvalue
263 REAL,DIMENSION(:,:,:),INTENT(OUT) :: pvalu3d !value array
264 !
265 integer :: status
266 character(len=80) :: haction
267 integer,save :: ndims=3
268 integer :: kvartype
269 integer,DIMENSION(:),ALLOCATABLE :: nvardimid,nvardimlen
270 character(len=80),DIMENSION(:),ALLOCATABLE :: nvardimnam
271 integer :: jloop2, jloop
272 integer :: j1,j2,j3
273 integer :: ngatts
274 character(len=80),DIMENSION(:),ALLOCATABLE :: hname
275 real :: zmiss1,zmiss2,zmiss3
276 real :: zscfa, zoffs
277 REAL,DIMENSION(:,:,:),ALLOCATABLE :: zvalu3d !value array
278 INTEGER,DIMENSION(:,:,:),ALLOCATABLE :: ivalu3d
279 REAL(KIND=JPRB) :: zhook_handle
280 !
281 include 'netcdf.inc'
282 !
283 IF (lhook) CALL dr_hook('MODE_READ_NETCDF_MERCATOR:GET3DCDF',0,zhook_handle)
284 pmissvalue=-9999.9
285 ALLOCATE(nvardimid(ndims))
286 ALLOCATE(nvardimlen(ndims))
287 ALLOCATE(nvardimnam(ndims))
288 nvardimid(:)=0
289 nvardimlen(:)=0
290 nvardimnam(:)=' '
291 !
292 haction='get variable type'
293 status=nf_inq_vartype(kcdf_id,idvar,kvartype)
294 if (status/=nf_noerr) CALL handle_err_mer(status,haction)
295 !write(0,*) 'variable type = ',KVARTYPE
296 !
297 haction='get variable dimensions identifiant'
298 status=nf_inq_vardimid(kcdf_id,idvar,nvardimid)
299 if (status/=nf_noerr) CALL handle_err_mer(status,haction)
300 !write(0,*) 'variable dimension identifiant ',NVARDIMID
301 !
302 haction='get attributs'
303 status=nf_inq_varnatts(kcdf_id,idvar,ngatts)
304 if (status/=nf_noerr) CALL handle_err_mer(status,haction)
305 !write(0,*) 'number of attributes = ',NGATTS
306 allocate(hname(1:ngatts))
307 !
308 zscfa=1.
309 zoffs=0.
310 DO jloop=1,ngatts
311  status=nf_inq_attname(kcdf_id,idvar,jloop,hname(jloop))
312  if (status/=nf_noerr) CALL handle_err_mer(status,haction)
313  !write(0,*) 'attributes names = ', hname(JLOOP)
314  if (trim(hname(jloop))=='missing_value') then
315  !write(0,*) 'missing value found '
316  haction='get missing value'
317  status=nf_get_att_double(kcdf_id,idvar,"missing_value",pmissvalue)
318  if (status/=nf_noerr) CALL handle_err_mer(status,haction)
319  !write(0,*) 'missing value = ',PMISSVALUE
320  else
321  if (trim(hname(jloop))=='_FillValue') then
322  !write(0,*) 'missing value found '
323  haction='get _FillValue'
324  status=nf_get_att_double(kcdf_id,idvar,"_FillValue",pmissvalue)
325  if (status/=nf_noerr) CALL handle_err_mer(status,haction)
326  !write(0,*) 'missing value = ',PMISSVALUE
327  endif
328  endif
329  if (trim(hname(jloop))=='scale_factor') then
330  !write(0,*) 'missing value found '
331  haction='get scale factor'
332  status=nf_get_att_double(kcdf_id,idvar,"scale_factor",zscfa)
333  if (status/=nf_noerr) CALL handle_err_mer(status,haction)
334  !write(0,*) 'missing value = ',PMISSVALUE
335  endif
336  if (trim(hname(jloop))=='add_offset') then
337  !write(0,*) 'missing value found '
338  haction='get offset'
339  status=nf_get_att_double(kcdf_id,idvar,"add_offset",zoffs)
340  if (status/=nf_noerr) CALL handle_err_mer(status,haction)
341  !write(0,*) 'missing value = ',PMISSVALUE
342  endif
343 ENDDO
344 !
345 !
346 DO jloop2=1,ndims
347  haction='get variable dimensions name'
348  status=nf_inq_dimname(kcdf_id,nvardimid(jloop2),nvardimnam(jloop2))
349  if (status/=nf_noerr) CALL handle_err_mer(status,haction)
350  haction='get variable dimensions length'
351  status=nf_inq_dimlen(kcdf_id,nvardimid(jloop2),nvardimlen(jloop2))
352  if (status/=nf_noerr) CALL handle_err_mer(status,haction)
353  !write(0,*) 'variable dimension ',JLOOP2,' named ',NVARDIMNAM(JLOOP2),&
354  ! &'has a length of',NVARDIMLEN(JLOOP2)
355 ENDDO
356 !
357 IF (kvartype>=5) then
358  ALLOCATE(zvalu3d(1:nvardimlen(1),1:nvardimlen(2),1:nvardimlen(3)))
359  zvalu3d=0.
360  haction='get variable values (3D)'
361  status=nf_get_var_double(kcdf_id,idvar,zvalu3d(:,:,:))
362  if (status/=nf_noerr) CALL handle_err_mer(status,haction)
363 ELSE
364  ALLOCATE(ivalu3d(1:nvardimlen(1),1:nvardimlen(2),1:nvardimlen(3)))
365  ivalu3d=0.
366  haction='get variable values (3D)'
367  status=nf_get_var_int(kcdf_id,idvar,ivalu3d(:,:,:))
368  if (status/=nf_noerr) CALL handle_err_mer(status,haction)
369 ENDIF
370 !
371 pvalu3d(:,:,:)=xundef
372 DO j1=1,nvardimlen(1)
373  DO j2=1,nvardimlen(2)
374  DO j3=1,nvardimlen(3)
375  IF (kvartype>=5) THEN
376  IF (zvalu3d(j1,j2,j3)/=pmissvalue) pvalu3d(j1,j2,j3)=zvalu3d(j1,j2,j3)*zscfa+zoffs
377  ELSE
378  IF (ivalu3d(j1,j2,j3)/=pmissvalue) pvalu3d(j1,j2,j3)=ivalu3d(j1,j2,j3)*zscfa+zoffs
379  ENDIF
380  ENDDO
381  ENDDO
382 ENDDO
383 !
384  CALL get1dcdf(kcdf_id,nvardimid(1),zmiss1,pdim1)
385  CALL get1dcdf(kcdf_id,nvardimid(2),zmiss2,pdim2)
386  CALL get1dcdf(kcdf_id,nvardimid(3),zmiss3,pdim3)
387 hdim1name=nvardimnam(1)
388 hdim2name=nvardimnam(2)
389 hdim3name=nvardimnam(3)
390 IF (ALLOCATED(zvalu3d )) DEALLOCATE(zvalu3d)
391 IF (lhook) CALL dr_hook('MODE_READ_NETCDF_MERCATOR:GET3DCDF',1,zhook_handle)
392 !
393 END SUBROUTINE get3dcdf
394 !--------------------------------------------------------------------
395 !-------------------------------------------------------------------
396 !------------------------------------------------------------------------------
397 !==============================================================================
398 ! ####################
399  SUBROUTINE read_dim_cdf(HFILENAME,HNCVARNAME,KDIM)
400 ! ####################
401 !
402 IMPLICIT NONE
403 !
404  CHARACTER(LEN=28), INTENT(IN) :: hfilename ! Name of the field file.
405  CHARACTER(LEN=28), INTENT(IN) :: hncvarname ! Name of variable to read in netcdf file
406 INTEGER, INTENT(OUT):: kdim ! value of dimension to get
407 !
408 integer :: status
409 integer :: kcdf_id
410 integer :: nbvars
411 character(len=80) :: haction
412 character(len=80),DIMENSION(:),ALLOCATABLE :: varname
413 integer ::jloop1,jloop
414 integer ::id_vartoget,id_vartoget1,id_vartoget2
415 integer ::nvardims
416 integer,DIMENSION(2) ::nlen2d
417 REAL(KIND=JPRB) :: zhook_handle
418 !
419 include 'netcdf.inc'
420 !
421 !* 1. Open the netcdf file
422 ! --------------------
423 IF (lhook) CALL dr_hook('MODE_READ_NETCDF_MERCATOR:READ_DIM_CDF',0,zhook_handle)
424 haction='open netcdf'
425 status=nf_open(hfilename,nf_nowrite,kcdf_id)
426 !write(0,*) 'identifiant de ',HFILENAME,'=',kcdf_id
427 if (status/=nf_noerr) then
428  CALL handle_err_mer(status,haction)
429 !else
430 ! write(0,*) 'netcdf file opened: ',HFILENAME
431 endif
432 !
433 !-----------
434 !
435 !* 2. get the number of variables in netcdf file
436 ! ------------------------------------------
437 haction='get number of variables'
438 status=nf_inq_nvars(kcdf_id,nbvars)
439 if (status/=nf_noerr) CALL handle_err_mer(status,haction)
440 !write(0,*) 'nb vars', NBVARS
441 ALLOCATE(varname(nbvars))
442 !
443 !-----------
444 !
445 !* 3. get the variables names in netcdf file
446 ! --------------------------------------
447 id_vartoget1=0
448 id_vartoget2=0
449 DO jloop1=1,nbvars
450  haction='get variables names'
451  status=nf_inq_varname(kcdf_id,jloop1,varname(jloop1))
452  if (status/=nf_noerr) CALL handle_err_mer(status,haction)
453  !write(0,*) 'var',JLOOP1,' name: ',VARNAME(JLOOP1)
454  if (varname(jloop1)==hncvarname) then
455  !write(0,*) 'var',JLOOP1,' corresponding to variable required'
456  id_vartoget1=jloop1
457  endif
458  if (varname(jloop1)/=hncvarname) then
459  if((lgt(trim(varname(jloop1)),trim(hncvarname))).AND.&
460  (scan(trim(varname(jloop1)),trim(hncvarname))==1)) then
461  !write(0,*) 'var',JLOOP1,VARNAME(JLOOP1),' could correspond to variable required ?'
462  !write(0,*) HNCVARNAME,' is variable required; only ',VARNAME(JLOOP1),' found'
463  id_vartoget2=jloop1
464  endif
465  endif
466 ENDDO
467 if (id_vartoget1/=0) then
468  id_vartoget=id_vartoget1
469 else
470  id_vartoget=id_vartoget2
471 endif
472 if (id_vartoget==0) then
473  haction='close netcdf'
474  status=nf_close(kcdf_id)
475  if (status/=nf_noerr) CALL handle_err_mer(status,haction)
476  CALL abor1_sfx('MODE_READ_NETCDF_MERCATOR: READ_DIM_CDF')
477 endif
478 !-----------
479 !
480 !* 4. get the total dimension of HNCVARNAME
481 ! -------------------------------------
482 !
483 ! 4.1 get the variable dimensions number
484 ! -----------------------------------
485 !
486 haction='get variable dimensions number'
487 status=nf_inq_varndims(kcdf_id,id_vartoget,nvardims)
488 if (status/=nf_noerr) CALL handle_err_mer(status,haction)
489 !write(0,*) 'variable dimensions number = ',NVARDIMS
490 !
491 ! 4.2 get the variable dimensions length
492 ! ----------------------------------
493 SELECT CASE (nvardims)
494 !CAS 1D
495  CASE (1)
496  haction='get variable dimensions length'
497  status=nf_inq_dimlen(kcdf_id,id_vartoget,kdim)
498  if (status/=nf_noerr) CALL handle_err_mer(status,haction)
499 !
500 !CAS 2D
501  CASE (2)
502  kdim=1
503  DO jloop=1,nvardims
504  haction='get variable dimensions length'
505  status=nf_inq_dimlen(kcdf_id,id_vartoget,nlen2d(jloop))
506  if (status/=nf_noerr) CALL handle_err_mer(status,haction)
507  kdim=kdim*nlen2d(jloop)
508  ENDDO
509 END SELECT
510 !-----------
511 !* 10. Close the netcdf file
512 ! ---------------------
513 haction='close netcdf'
514 status=nf_close(kcdf_id)
515 if (status/=nf_noerr) CALL handle_err_mer(status,haction)
516 !write(0,*) 'OK: netcdf file closed: ',HFILENAME
517 !
518 !-----------
519 !* 11. Deallocate
520 ! ----------
521 IF (ALLOCATED(varname )) DEALLOCATE(varname)
522 IF (lhook) CALL dr_hook('MODE_READ_NETCDF_MERCATOR:READ_DIM_CDF',1,zhook_handle)
523 !
524 END SUBROUTINE read_dim_cdf
525 !-------------------------------------------------------------------
526 !-------------------------------------------------------------------
527 ! ####################
528  SUBROUTINE prep_netcdf_grid(HFILENAME,HNCVARNAME)
529 ! ####################
530 !
532 USE modd_surf_par
533 !
534 IMPLICIT NONE
535 !
536  CHARACTER(LEN=28), INTENT(IN) :: hfilename ! Name of the field file.
537  CHARACTER(LEN=28), INTENT(IN) :: hncvarname ! Name of variable to read in netcdf file
538 !
539 integer :: status
540 integer :: kcdf_id
541 integer :: nbvars
542 character(len=80) :: haction
543 character(len=80),DIMENSION(:),ALLOCATABLE :: varname
544 integer,DIMENSION(:),ALLOCATABLE :: nvardimid
545 integer ::jloop1,jloop
546 integer ::id_vartoget,id_vartoget1,id_vartoget2
547 integer ::nvardims
548 integer,DIMENSION(3) ::ndimlen
549 character(LEN=80),DIMENSION(3) :: ndimnam
550 integer :: idim
551 integer :: inlon
552 real :: zzlamiss,zzlomiss
553 REAL(KIND=JPRB) :: zhook_handle
554 !
555 include 'netcdf.inc'
556 !
557 IF (lhook) CALL dr_hook('MODE_READ_NETCDF_MERCATOR:PREP_NETCDF_GRID',0,zhook_handle)
558 ninlat =-nundef
559 nindepth=-nundef
560 nilength=-nundef
561 !
562 xilat1=xundef
563 xilon1=xundef
564 xilat2=xundef
565 xilon2=xundef
566 !* 1. Open the netcdf file
567 ! --------------------
568 haction='open netcdf'
569 status=nf_open(hfilename,nf_nowrite,kcdf_id)
570 !write(0,*) 'identifiant de ',HFILENAME,'=',kcdf_id
571 if (status/=nf_noerr) then
572  CALL handle_err_mer(status,haction)
573 !else
574 ! write(0,*) 'netcdf file opened: ',HFILENAME
575 endif
576 !
577 !-----------
578 !
579 !* 2. get the number of variables in netcdf file
580 ! ------------------------------------------
581 haction='get number of variables'
582 status=nf_inq_nvars(kcdf_id,nbvars)
583 if (status/=nf_noerr) CALL handle_err_mer(status,haction)
584 !write(0,*) 'nb vars', NBVARS
585 ALLOCATE(varname(nbvars))
586 !
587 !-----------
588 !
589 !* 3. get the variables names in netcdf file
590 ! --------------------------------------
591 id_vartoget1=0
592 id_vartoget2=0
593 DO jloop1=1,nbvars
594  haction='get variables names'
595  status=nf_inq_varname(kcdf_id,jloop1,varname(jloop1))
596  if (status/=nf_noerr) CALL handle_err_mer(status,haction)
597  !write(0,*) 'var',JLOOP1,' name: ',VARNAME(JLOOP1)
598  if (varname(jloop1)==hncvarname) then
599  !write(0,*) 'var',JLOOP1,' corresponding to variable required'
600  id_vartoget1=jloop1
601  endif
602  if (varname(jloop1)/=hncvarname) then
603  if((lgt(trim(varname(jloop1)),trim(hncvarname))).AND.&
604  (scan(trim(varname(jloop1)),trim(hncvarname))==1)) then
605  !write(0,*) 'var',JLOOP1,VARNAME(JLOOP1),' could correspond to variable required ?'
606  !write(0,*) HNCVARNAME,' is variable required; only ',VARNAME(JLOOP1),' found'
607  id_vartoget2=jloop1
608  endif
609  endif
610 ENDDO
611 if (id_vartoget1/=0) then
612  id_vartoget=id_vartoget1
613 else
614  id_vartoget=id_vartoget2
615 endif
616 if (id_vartoget==0) then
617  haction='close netcdf'
618  status=nf_close(kcdf_id)
619  if (status/=nf_noerr) CALL handle_err_mer(status,haction)
620  IF (lhook) CALL dr_hook('MODE_READ_NETCDF_MERCATOR:PREP_NETCDF_GRID',1,zhook_handle)
621  RETURN
622 endif
623 nilength=0
624 !-----------
625 !
626 !* 4. get the total dimension of HNCVARNAME
627 ! -------------------------------------
628 !
629 ! 4.1 get the variable dimensions number
630 ! -----------------------------------
631 !
632 haction='get variable dimensions number'
633 status=nf_inq_varndims(kcdf_id,id_vartoget,nvardims)
634 if (status/=nf_noerr) CALL handle_err_mer(status,haction)
635 !write(0,*) 'variable dimensions number = ',NVARDIMS
636 ALLOCATE(nvardimid(nvardims))
637 haction='get variable dimensions identifiant'
638 status=nf_inq_vardimid(kcdf_id,id_vartoget,nvardimid)
639 if (status/=nf_noerr) CALL handle_err_mer(status,haction)
640 !
641 ! 4.2 get the variable dimensions length
642 ! ----------------------------------
643 SELECT CASE (nvardims)
644 !CAS 1D
645  CASE (1)
646  haction='get variable dimensions length'
647  status=nf_inq_dimlen(kcdf_id,id_vartoget,idim)
648  if (status/=nf_noerr) CALL handle_err_mer(status,haction)
649 !
650 !CAS 2D,3D
651  CASE (2,3)
652  DO jloop=1,nvardims
653  haction='get variable dimensions length'
654  status=nf_inq_dimlen(kcdf_id,nvardimid(jloop),ndimlen(jloop))
655  if (status/=nf_noerr) CALL handle_err_mer(status,haction)
656  haction='get variable dimensions names'
657  status=nf_inq_dimname(kcdf_id,nvardimid(jloop),ndimnam(jloop))
658  if (status/=nf_noerr) CALL handle_err_mer(status,haction)
659  if ((ndimnam(jloop)=='lat').OR.(ndimnam(jloop)=='latitude')) then
660  ninlat=ndimlen(jloop)
661  if (.not.allocated(xilatarray)) allocate(xilatarray(ndimlen(jloop)))
662  if (.not.allocated(ninlon)) allocate(ninlon(ninlat))
663  CALL get1dcdf(kcdf_id,nvardimid(jloop),zzlamiss,xilatarray(:))
664  endif
665  if ((ndimnam(jloop)=='lon').OR.(ndimnam(jloop)=='longitude')) then
666  inlon=ndimlen(jloop)
667  if (.not.allocated(xilonarray)) allocate(xilonarray(ndimlen(jloop)))
668  CALL get1dcdf(kcdf_id,nvardimid(jloop),zzlomiss,xilonarray(:))
669  endif
670  if (ndimnam(jloop)=='depth') nindepth=ndimlen(jloop)
671  ENDDO
672  ninlon(:)=inlon
673 END SELECT
674 !-----------
675 !* 10. Close the netcdf file
676 ! ---------------------
677 haction='close netcdf'
678 status=nf_close(kcdf_id)
679 if (status/=nf_noerr) CALL handle_err_mer(status,haction)
680 !write(0,*) 'OK: netcdf file closed: ',HFILENAME
681 !
682 !-----------
683 !GRID PARAM FOR HORIBL_SURF
684 DO jloop1=1,ninlat
685  nilength = nilength + ninlon(jloop1)
686 ENDDO
687 xilat1=xilatarray(1)
688 xilon1=xilonarray(1)
689 xilat2=xilatarray(SIZE(xilatarray))
690 xilon2=xilonarray(SIZE(xilonarray))
691 !
692 !* 11. Deallocate
693 ! ----------
694 IF (ALLOCATED(varname )) DEALLOCATE(varname)
695 IF (lhook) CALL dr_hook('MODE_READ_NETCDF_MERCATOR:PREP_NETCDF_GRID',1,zhook_handle)
696 !
697 END SUBROUTINE prep_netcdf_grid
698 !------------------------------------------------------------------------------
699 !==============================================================================
700 ! ####################
701  SUBROUTINE read_z1d_cdf(HFILENAME,HNCVARNAME,PVAL)
702 ! ####################
703 !
704 IMPLICIT NONE
705 !
706  CHARACTER(LEN=28), INTENT(IN) :: hfilename ! Name of the field file.
707  CHARACTER(LEN=28), INTENT(IN) :: hncvarname ! Name of variable to read in netcdf file
708 REAL, DIMENSION(:), INTENT(OUT) :: pval ! value to get
709 !
710 integer :: status
711 integer :: kcdf_id
712 integer :: nbvars
713 character(len=80) :: haction
714 character(len=80),DIMENSION(:),ALLOCATABLE :: varname
715 integer ::jloop1
716 integer ::id_vartoget,id_vartoget1,id_vartoget2
717 integer ::nvardims
718 integer ::nlen
719 real,DIMENSION(:),ALLOCATABLE :: zvalu
720 real :: zmiss
721 REAL(KIND=JPRB) :: zhook_handle
722 !
723 include 'netcdf.inc'
724 !
725 !* 1. Open the netcdf file
726 ! --------------------
727 IF (lhook) CALL dr_hook('MODE_READ_NETCDF_MERCATOR:READ_Z1D_CDF',0,zhook_handle)
728 status=-9999
729 kcdf_id=-9999
730 haction='open netcdf'
731 status=nf_open(hfilename,nf_nowrite,kcdf_id)
732 if (status/=nf_noerr) then
733  CALL handle_err_mer(status,haction)
734 endif
735 !-----------
736 !* 2. get the number of variables in netcdf file
737 ! ------------------------------------------
738 haction='get number of variables'
739 status=nf_inq_nvars(kcdf_id,nbvars)
740 if (status/=nf_noerr) CALL handle_err_mer(status,haction)
741 !write(0,*) 'nb vars', NBVARS
742 ALLOCATE(varname(nbvars))
743 !-----------
744 !* 3. get the variables names in netcdf file
745 ! --------------------------------------
746 id_vartoget1=0
747 id_vartoget2=0
748 DO jloop1=1,nbvars
749  haction='get variables names'
750  status=nf_inq_varname(kcdf_id,jloop1,varname(jloop1))
751  if (status/=nf_noerr) CALL handle_err_mer(status,haction)
752  if (varname(jloop1)==hncvarname) then
753  id_vartoget1=jloop1
754  endif
755  if (varname(jloop1)/=hncvarname) then
756  if((lgt(trim(varname(jloop1)),trim(hncvarname))).AND.&
757  (scan(trim(varname(jloop1)),trim(hncvarname))==1)) then
758  id_vartoget2=jloop1
759  endif
760  endif
761 ENDDO
762 if (id_vartoget1/=0) then
763  id_vartoget=id_vartoget1
764 else
765  id_vartoget=id_vartoget2
766 endif
767 if (id_vartoget==0) then
768  haction='close netcdf'
769  status=nf_close(kcdf_id)
770  if (status/=nf_noerr) CALL handle_err_mer(status,haction)
771  CALL abor1_sfx('MODE_READ_NETCDF_MERCATOR: READ_Z1D_CDF')
772 endif
773 !-----------
774 !* 4. get the variable in netcdf file
775 ! -------------------------------
776 ! 4.1 get the variable dimensions number
777 ! -----------------------------------
778 haction='get variable dimensions number'
779 status=nf_inq_varndims(kcdf_id,id_vartoget,nvardims)
780 if (status/=nf_noerr) CALL handle_err_mer(status,haction)
781 !
782 ! 4.2 get the variable dimensions length and values
783 ! ----------------------------------------------
784 SELECT CASE (nvardims)
785 !CAS 1D
786  CASE (1)
787  haction='get variable dimensions length'
788  status=nf_inq_dimlen(kcdf_id,id_vartoget,nlen)
789  if (status/=nf_noerr) CALL handle_err_mer(status,haction)
790  ALLOCATE(zvalu(nlen))
791  !write(0,*) 'call GET1DCDF'
792  CALL get1dcdf(kcdf_id,id_vartoget,zmiss,zvalu)
793  pval(:)=zvalu(:)
794 !CAS 2D
795  CASE (2)
796  write(0,*) 'YOU ARE TRYING TO READ A 2D FIELD FOR :', trim(hncvarname)
797  CALL abor1_sfx('MODE_READ_NETCDF_MERCATOR: READ_Z1D_CDF')
798 END SELECT
799 !-----------
800 !* 5. Close the netcdf file
801 ! ---------------------
802 haction='close netcdf'
803 status=nf_close(kcdf_id)
804 if (status/=nf_noerr) CALL handle_err_mer(status,haction)
805 !-----------
806 !* 6. Deallocate
807 ! ----------
808 IF (ALLOCATED(varname )) DEALLOCATE(varname)
809 IF (ALLOCATED(zvalu )) DEALLOCATE(zvalu )
810 !!
811 IF (lhook) CALL dr_hook('MODE_READ_NETCDF_MERCATOR:READ_Z1D_CDF',1,zhook_handle)
812 END SUBROUTINE read_z1d_cdf
813 !------------------------------------------------------------------------------
814 !==============================================================================
815 ! ####################
816  SUBROUTINE read_latlonval_cdf(HFILENAME,HNCVARNAME,PLON,PLAT,PVAL)
817 ! ####################
818 !
819 IMPLICIT NONE
820 !
821  CHARACTER(LEN=28), INTENT(IN) :: hfilename ! Name of the field file.
822  CHARACTER(LEN=28), INTENT(IN) :: hncvarname ! Name of variable to read in netcdf file
823 REAL, DIMENSION(:), INTENT(OUT) :: plon,plat ! Longitudes/latitudes in netcdf file
824 REAL, DIMENSION(:), INTENT(OUT) :: pval ! value to get
825 !
826 integer :: status
827 integer :: kcdf_id
828 integer :: nbvars
829 character(len=80) :: haction
830 character(len=80),DIMENSION(:),ALLOCATABLE :: varname
831 integer ::jloop1,jdim1,jdim2,jloop
832 integer ::id_vartoget,id_vartoget1,id_vartoget2
833 integer ::nvardims
834 integer ::nlen
835 integer,DIMENSION(2) ::nlen2d
836 integer,DIMENSION(:),ALLOCATABLE :: nvardimid,nvardimlen
837 character(len=80),DIMENSION(:),ALLOCATABLE :: nvardimnam
838 real,DIMENSION(:),ALLOCATABLE :: zvalu
839 real,DIMENSION(:,:),ALLOCATABLE :: zvalu2d
840 real :: zmiss
841 real,DIMENSION(:),ALLOCATABLE :: zdim1
842 real,DIMENSION(:),ALLOCATABLE :: zdim2
843 character(len=80) :: ydim1name,ydim2name
844 integer :: ilonfound,ilatfound, iarg
845 REAL(KIND=JPRB) :: zhook_handle
846 !
847 include 'netcdf.inc'
848 !
849 !
850 !
851 !* 1. Open the netcdf file
852 ! --------------------
853 IF (lhook) CALL dr_hook('MODE_READ_NETCDF_MERCATOR:READ_LATLONVAL_CDF',0,zhook_handle)
854 status=-9999
855 kcdf_id=-9999
856 haction='open netcdf'
857 status=nf_open(hfilename,nf_nowrite,kcdf_id)
858 !write(0,*) 'status=',status
859 !write(0,*) 'identifiant de ',HFILENAME,'=',kcdf_id
860 if (status/=nf_noerr) then
861  CALL handle_err_mer(status,haction)
862 !else
863 ! write(0,*) 'netcdf file opened: ',HFILENAME
864 endif
865 !
866 !-----------
867 !
868 !* 2. get the number of variables in netcdf file
869 ! ------------------------------------------
870 haction='get number of variables'
871 status=nf_inq_nvars(kcdf_id,nbvars)
872 if (status/=nf_noerr) CALL handle_err_mer(status,haction)
873 !write(0,*) 'nb vars', NBVARS
874 ALLOCATE(varname(nbvars))
875 !
876 !-----------
877 !
878 !* 3. get the variables names in netcdf file
879 ! --------------------------------------
880 id_vartoget1=0
881 id_vartoget2=0
882 DO jloop1=1,nbvars
883  haction='get variables names'
884  status=nf_inq_varname(kcdf_id,jloop1,varname(jloop1))
885  if (status/=nf_noerr) CALL handle_err_mer(status,haction)
886  !write(0,*) 'var',JLOOP1,' name: ',VARNAME(JLOOP1)
887  if (varname(jloop1)==hncvarname) then
888  !write(0,*) 'var',JLOOP1,' corresponding to variable required'
889  id_vartoget1=jloop1
890  endif
891  if (varname(jloop1)/=hncvarname) then
892  if((lgt(trim(varname(jloop1)),trim(hncvarname))).AND.&
893  (scan(trim(varname(jloop1)),trim(hncvarname))==1)) then
894  !write(0,*) 'var',JLOOP1,VARNAME(JLOOP1),' could correspond to variable required ?'
895  !write(0,*) HNCVARNAME,' is variable required; only ',VARNAME(JLOOP1),' found'
896  id_vartoget2=jloop1
897  endif
898  endif
899 ENDDO
900 if (id_vartoget1/=0) then
901  id_vartoget=id_vartoget1
902 else
903  id_vartoget=id_vartoget2
904 endif
905 if (id_vartoget==0) then
906  haction='close netcdf'
907  status=nf_close(kcdf_id)
908  if (status/=nf_noerr) CALL handle_err_mer(status,haction)
909  CALL abor1_sfx('MODE_READ_NETCDF_MERCATOR: READ_LATLONVAL_CDF')
910 endif
911 !-----------
912 !
913 !* 4. get the variable in netcdf file
914 ! -------------------------------
915 !
916 ! 4.1 get the variable dimensions number
917 ! -----------------------------------
918 !
919 haction='get variable dimensions number'
920 status=nf_inq_varndims(kcdf_id,id_vartoget,nvardims)
921 if (status/=nf_noerr) CALL handle_err_mer(status,haction)
922 !write(0,*) 'variable dimensions number = ',NVARDIMS
923 !
924 ! 4.2 get the variable dimensions length and values
925 ! ----------------------------------------------
926 SELECT CASE (nvardims)
927 !CAS 1D
928  CASE (1)
929  haction='get variable dimensions length'
930  status=nf_inq_dimlen(kcdf_id,id_vartoget,nlen)
931  if (status/=nf_noerr) CALL handle_err_mer(status,haction)
932  ALLOCATE(zvalu(nlen))
933  !write(0,*) 'call GET1DCDF'
934  CALL get1dcdf(kcdf_id,id_vartoget,zmiss,zvalu)
935  pval(:)=zvalu(:)
936 !CAS 2D
937  CASE (2)
938  DO jloop=1,nvardims
939  haction='get variable dimensions length'
940  status=nf_inq_dimlen(kcdf_id,id_vartoget,nlen2d(jloop))
941  if (status/=nf_noerr) CALL handle_err_mer(status,haction)
942  ENDDO
943  ALLOCATE(zvalu2d(nlen2d(1),nlen2d(2)))
944  ALLOCATE(zdim1(nlen2d(1)))
945  ALLOCATE(zdim2(nlen2d(2)))
946  !write(0,*) 'call GET2DCDF'
947  CALL get2dcdf(kcdf_id,id_vartoget,zdim1,ydim1name,zdim2,ydim2name,&
948  zmiss,zvalu2d)
949  !write(0,*) 'YDIM1NAME: ',YDIM1NAME
950  !write(0,*) 'YDIM2NAME: ',YDIM2NAME
951  if ((ydim1name=='lon').OR.(ydim1name=='longitude')) ilonfound=1
952  if ((ydim2name=='lon').OR.(ydim2name=='longitude')) ilonfound=2
953  if ((ydim1name=='lat').OR.(ydim1name=='latitude')) ilatfound=1
954  if ((ydim2name=='lat').OR.(ydim2name=='latitude')) ilatfound=2
955  iarg=0
956 !
957 ! 4.3 complete arrays
958 ! ---------------
959  IF ((ilonfound==1).AND.(ilatfound==2)) then
960  !write(0,*) 'ILONFOUND',ILONFOUND,'ILATFOUND',ILATFOUND
961  DO jdim1=1,SIZE(zdim1)
962  DO jdim2=1,SIZE(zdim2)
963  iarg=iarg+1
964  pval(iarg)=zvalu2d(jdim1,jdim2)
965  plon(iarg)=zdim1(jdim1)
966  plat(iarg)=zdim2(jdim2)
967  ENDDO
968  ENDDO
969  ELSEIF ((ilonfound==2).AND.(ilatfound==1)) then
970  !write(0,*) 'ILONFOUND',ILONFOUND,'ILATFOUND',ILATFOUND
971  DO jdim1=1,SIZE(zdim1)
972  DO jdim2=1,SIZE(zdim2)
973  iarg=iarg+1
974  pval(iarg)=zvalu2d(jdim1,jdim2)
975  plat(iarg)=zdim1(jdim1)
976  plon(iarg)=zdim2(jdim2)
977  ENDDO
978  ENDDO
979  ELSE
980  write(0,*) '*****WARNING*****: incompatible dimensions to lat/lon/value arrays'
981  ENDIF
982 !
983 END SELECT
984 !
985 !
986 !-----------
987 !* 10. Close the netcdf file
988 ! ---------------------
989 haction='close netcdf'
990 status=nf_close(kcdf_id)
991 if (status/=nf_noerr) CALL handle_err_mer(status,haction)
992 !write(0,*) 'OK: netcdf file closed: ',HFILENAME
993 !
994 !-----------
995 !* 11. Deallocate
996 ! ----------
997 IF (ALLOCATED(varname )) DEALLOCATE(varname)
998 IF (ALLOCATED(zvalu )) DEALLOCATE(zvalu )
999 IF (ALLOCATED(zvalu2d )) DEALLOCATE(zvalu2d)
1000 IF (ALLOCATED(zdim1 )) DEALLOCATE(zdim1 )
1001 IF (ALLOCATED(zdim2 )) DEALLOCATE(zdim2 )
1002 !
1003 !
1004 IF (ALLOCATED(nvardimid )) DEALLOCATE(nvardimid )
1005 IF (ALLOCATED(nvardimnam )) DEALLOCATE(nvardimnam)
1006 IF (ALLOCATED(nvardimlen )) DEALLOCATE(nvardimlen)
1007 IF (lhook) CALL dr_hook('MODE_READ_NETCDF_MERCATOR:READ_LATLONVAL_CDF',1,zhook_handle)
1008 END SUBROUTINE read_latlonval_cdf
1009 !------------------------------------------------------------------------------
1010 !==============================================================================
1011 ! ####################
1012  SUBROUTINE read_latlondepval_cdf(HFILENAME,HNCVARNAME,PLON,PLAT,PDEP,PVAL)
1013 ! ####################
1014 !
1015 IMPLICIT NONE
1016 !
1017  CHARACTER(LEN=28), INTENT(IN) :: hfilename ! Name of the field file.
1018  CHARACTER(LEN=28), INTENT(IN) :: hncvarname ! Name of variable to read in netcdf file
1019 REAL, DIMENSION(:), INTENT(OUT) :: plon,plat ! Longitudes/latitudes in netcdf file
1020 REAL, DIMENSION(:), INTENT(OUT) :: pdep ! depth in netcdf file
1021 REAL, DIMENSION(:,:), INTENT(OUT) :: pval ! value to get
1022 !
1023 integer :: status
1024 integer :: kcdf_id
1025 integer :: nbvars
1026 character(len=80) :: haction
1027 character(len=80),DIMENSION(:),ALLOCATABLE :: varname
1028 integer ::jloop1,jdim1,jdim2,jdim3,jloop
1029 !integer ::JLOOP2,JLOOP
1030 integer ::id_vartoget,id_vartoget1,id_vartoget2
1031 integer ::nvardims
1032 integer,DIMENSION(3) ::nlen3d
1033 integer,DIMENSION(:),ALLOCATABLE :: nvardimid,nvardimlen
1034 character(len=80),DIMENSION(:),ALLOCATABLE :: nvardimnam
1035 real,DIMENSION(:,:,:),ALLOCATABLE :: zvalu3d
1036 real :: zmiss
1037 real,DIMENSION(:),ALLOCATABLE :: zdim1
1038 real,DIMENSION(:),ALLOCATABLE :: zdim2
1039 real,DIMENSION(:),ALLOCATABLE :: zdim3
1040 character(len=80) :: ydim1name,ydim2name,ydim3name
1041 integer :: ilonfound,ilatfound,idepfound
1042 integer :: iarg
1043 REAL(KIND=JPRB) :: zhook_handle
1044 !
1045 include 'netcdf.inc'
1046 !
1047 !
1048 !
1049 !* 1. Open the netcdf file
1050 ! --------------------
1051 IF (lhook) CALL dr_hook('MODE_READ_NETCDF_MERCATOR:READ_LATLONDEPVAL_CDF',0,zhook_handle)
1052 haction='open netcdf'
1053 status=nf_open(hfilename,nf_nowrite,kcdf_id)
1054 if (status/=nf_noerr) CALL handle_err_mer(status,haction)
1055 !write(0,*) 'netcdf file opened: ',HFILENAME
1056 !
1057 !-----------
1058 !
1059 !* 2. get the number of variables in netcdf file
1060 ! ------------------------------------------
1061 haction='get number of variables'
1062 status=nf_inq_nvars(kcdf_id,nbvars)
1063 if (status/=nf_noerr) CALL handle_err_mer(status,haction)
1064 !write(0,*) 'nb vars', NBVARS
1065 ALLOCATE(varname(nbvars))
1066 !
1067 !-----------
1068 !
1069 !* 3. get the variables names in netcdf file
1070 ! --------------------------------------
1071 id_vartoget1=0
1072 id_vartoget2=0
1073 DO jloop1=1,nbvars
1074  haction='get variables names'
1075  status=nf_inq_varname(kcdf_id,jloop1,varname(jloop1))
1076  if (status/=nf_noerr) CALL handle_err_mer(status,haction)
1077  !write(0,*) 'var',JLOOP1,' name: ',VARNAME(JLOOP1)
1078  if (varname(jloop1)==hncvarname) then
1079  !write(0,*) 'var',JLOOP1,' corresponding to variable required'
1080  id_vartoget1=jloop1
1081  endif
1082  if (varname(jloop1)/=hncvarname) then
1083  if((lgt(trim(varname(jloop1)),trim(hncvarname))).AND.&
1084  (scan(trim(varname(jloop1)),trim(hncvarname))==1)) then
1085  !write(0,*) 'var',JLOOP1,VARNAME(JLOOP1),' could correspond to variable required ?'
1086  !write(0,*) HNCVARNAME,' is variable required; only ',VARNAME(JLOOP1),' found'
1087  id_vartoget2=jloop1
1088  endif
1089  endif
1090 ENDDO
1091 if (id_vartoget1/=0) then
1092  id_vartoget=id_vartoget1
1093 else
1094  id_vartoget=id_vartoget2
1095 endif
1096 if (id_vartoget==0) then
1097  haction='close netcdf'
1098  status=nf_close(kcdf_id)
1099  if (status/=nf_noerr) CALL handle_err_mer(status,haction)
1100  CALL abor1_sfx('MODE_READ_NETCDF_MERCATOR: READ_LATLONDEPVAL_CDF')
1101 endif
1102 !-----------
1103 !
1104 !* 4. get the variable in netcdf file
1105 ! -------------------------------
1106 !
1107 ! 4.1 get the variable dimensions number
1108 ! -----------------------------------
1109 !
1110 haction='get variable dimensions number'
1111 status=nf_inq_varndims(kcdf_id,id_vartoget,nvardims)
1112 if (status/=nf_noerr) CALL handle_err_mer(status,haction)
1113 !write(0,*) 'variable dimensions number = ',NVARDIMS
1114 ALLOCATE(nvardimid(nvardims))
1115 haction='get variable dimensions identifiant'
1116 status=nf_inq_vardimid(kcdf_id,id_vartoget,nvardimid)
1117 if (status/=nf_noerr) CALL handle_err_mer(status,haction)
1118 !
1119 !
1120 ! 4.2 get the variable dimensions length and values
1121 ! ----------------------------------------------
1122 SELECT CASE (nvardims)
1123 !CAS 1D, 2D
1124  CASE (1,2)
1125  write(0,*) '********************************************'
1126  write(0,*) '* number of dimension to low: ',nvardims,' *'
1127  write(0,*) '* you need a 3-dimension variable *'
1128  write(0,*) '********************************************'
1129 !CAS 3D
1130  CASE (3)
1131  DO jloop=1,nvardims
1132  haction='get variable dimensions length'
1133  status=nf_inq_dimlen(kcdf_id,nvardimid(jloop),nlen3d(jloop))
1134  if (status/=nf_noerr) CALL handle_err_mer(status,haction)
1135  ENDDO
1136  ALLOCATE(zvalu3d(nlen3d(1),nlen3d(2),nlen3d(3)))
1137  ALLOCATE(zdim1(nlen3d(1)))
1138  ALLOCATE(zdim2(nlen3d(2)))
1139  ALLOCATE(zdim3(nlen3d(3)))
1140  !write(0,*) 'call GET3DCDF'
1141  CALL get3dcdf(kcdf_id,id_vartoget,zdim1,ydim1name,zdim2,ydim2name,&
1142  zdim3,ydim3name,zmiss,zvalu3d)
1143  !write(0,*) 'YDIM1NAME: ',YDIM1NAME
1144  !write(0,*) 'YDIM2NAME: ',YDIM2NAME
1145  !write(0,*) 'YDIM3NAME: ',YDIM3NAME
1146  if ((ydim1name=='lon').OR.(ydim1name=='longitude')) ilonfound=1
1147  if ((ydim2name=='lon').OR.(ydim2name=='longitude')) ilonfound=2
1148  if ((ydim3name=='lon').OR.(ydim3name=='longitude')) ilonfound=3
1149  if ((ydim1name=='lat').OR.(ydim1name=='latitude')) ilatfound=1
1150  if ((ydim2name=='lat').OR.(ydim2name=='latitude')) ilatfound=2
1151  if ((ydim3name=='lat').OR.(ydim3name=='latitude')) ilatfound=3
1152  if (ydim1name=='depth') idepfound=1
1153  if (ydim2name=='depth') idepfound=2
1154  if (ydim3name=='depth') idepfound=3
1155  iarg=0
1156  !write(0,*) 'ILONFOUND',ILONFOUND,'ILATFOUND',ILATFOUND,'IDEPFOUND',IDEPFOUND
1157 !!
1158 !! 4.3 complete arrays
1159 !! ---------------
1160  IF ((ilonfound==1).AND.(ilatfound==2).AND.(idepfound==3)) then
1161  !write(0,*) 'SIZE LON=',SIZE(ZDIM1),'SIZE LAT',SIZE(ZDIM2),'SIZE DEP',SIZE(ZDIM3)
1162  !write(0,*) 'SIZE PLON=',SIZE(PLON),'SIZE PLAT',SIZE(PLAT),'SIZE PDEP',SIZE(PDEP)
1163  pdep(:)=zdim3(:)
1164  DO jdim2=1,SIZE(zdim2)
1165  DO jdim1=1,SIZE(zdim1)
1166  iarg=iarg+1
1167  plon(iarg)=zdim1(jdim1)
1168  plat(iarg)=zdim2(jdim2)
1169  DO jdim3=1,SIZE(zdim3)
1170  pval(iarg,jdim3)=zvalu3d(jdim1,jdim2,jdim3)
1171  ENDDO
1172  ENDDO
1173  ENDDO
1174  !write(0,*) 'END complete arrays'
1175 !
1176  ELSEIF ((ilonfound==2).AND.(ilatfound==1).AND.(idepfound==3)) then
1177  pdep(:)=zdim3(:)
1178  DO jdim1=1,SIZE(zdim1)
1179  DO jdim2=1,SIZE(zdim2)
1180  iarg=iarg+1
1181  plon(iarg)=zdim2(jdim2)
1182  plat(iarg)=zdim1(jdim1)
1183  DO jdim3=1,SIZE(zdim3)
1184  pval(iarg,jdim3)=zvalu3d(jdim1,jdim2,jdim3)
1185  ENDDO
1186  ENDDO
1187  ENDDO
1188 !
1189  ELSEIF ((ilonfound==1).AND.(ilatfound==3).AND.(idepfound==2)) then
1190  pdep(:)=zdim2(:)
1191  DO jdim3=1,SIZE(zdim3)
1192  DO jdim1=1,SIZE(zdim1)
1193  iarg=iarg+1
1194  plon(iarg)=zdim1(jdim1)
1195  plat(iarg)=zdim3(jdim3)
1196  DO jdim2=1,SIZE(zdim2)
1197  pval(iarg,jdim2)=zvalu3d(jdim1,jdim2,jdim3)
1198  ENDDO
1199  ENDDO
1200  ENDDO
1201 !
1202  ELSEIF ((ilatfound==1).AND.(ilonfound==3).AND.(idepfound==2)) then
1203  pdep(:)=zdim2(:)
1204  DO jdim1=1,SIZE(zdim1)
1205  DO jdim3=1,SIZE(zdim3)
1206  iarg=iarg+1
1207  plon(iarg)=zdim3(jdim3)
1208  plat(iarg)=zdim1(jdim1)
1209  DO jdim2=1,SIZE(zdim2)
1210  pval(iarg,jdim2)=zvalu3d(jdim1,jdim2,jdim3)
1211  ENDDO
1212  ENDDO
1213  ENDDO
1214 !
1215  ELSEIF ((ilonfound==2).AND.(ilatfound==3).AND.(idepfound==1)) then
1216  pdep(:)=zdim1(:)
1217  DO jdim3=1,SIZE(zdim3)
1218  DO jdim2=1,SIZE(zdim2)
1219  iarg=iarg+1
1220  plon(iarg)=zdim2(jdim2)
1221  plat(iarg)=zdim3(jdim3)
1222  DO jdim1=1,SIZE(zdim1)
1223  pval(iarg,jdim1)=zvalu3d(jdim1,jdim2,jdim3)
1224  ENDDO
1225  ENDDO
1226  ENDDO
1227 !
1228  ELSEIF ((ilatfound==2).AND.(ilonfound==3).AND.(idepfound==1)) then
1229  pdep(:)=zdim1(:)
1230  DO jdim2=1,SIZE(zdim2)
1231  DO jdim3=1,SIZE(zdim3)
1232  iarg=iarg+1
1233  plon(iarg)=zdim3(jdim3)
1234  plat(iarg)=zdim2(jdim2)
1235  DO jdim1=1,SIZE(zdim1)
1236  pval(iarg,jdim1)=zvalu3d(jdim1,jdim2,jdim3)
1237  ENDDO
1238  ENDDO
1239  ENDDO
1240 !
1241  ELSE
1242  write(0,*) '*****WARNING*****: incompatible dimensions to lat/lon/value arrays'
1243  ENDIF
1244 !
1245 END SELECT
1246 !
1247 !-----------
1248 !* 10. Close the netcdf file
1249 ! ---------------------
1250 haction='close netcdf'
1251 !write(0,*) HACTION
1252 status=nf_close(kcdf_id)
1253 if (status/=nf_noerr) CALL handle_err_mer(status,haction)
1254 !write(0,*) 'OK: netcdf file closed: ',HFILENAME
1255 !
1256 !-----------
1257 !* 11. Deallocate
1258 ! ----------
1259 IF (ALLOCATED(varname )) DEALLOCATE(varname)
1260 IF (ALLOCATED(zvalu3d )) DEALLOCATE(zvalu3d)
1261 IF (ALLOCATED(zdim1 )) DEALLOCATE(zdim1 )
1262 IF (ALLOCATED(zdim2 )) DEALLOCATE(zdim2 )
1263 IF (ALLOCATED(zdim3 )) DEALLOCATE(zdim3 )
1264 !
1265 !
1266 IF (ALLOCATED(nvardimid )) DEALLOCATE(nvardimid )
1267 IF (ALLOCATED(nvardimnam )) DEALLOCATE(nvardimnam)
1268 IF (ALLOCATED(nvardimlen )) DEALLOCATE(nvardimlen)
1269 201 FORMAT(4(3x,f10.4))
1270 IF (lhook) CALL dr_hook('MODE_READ_NETCDF_MERCATOR:READ_LATLONDEPVAL_CDF',1,zhook_handle)
1271 END SUBROUTINE read_latlondepval_cdf
1272 !------------------------------------------------------------------------------
1273 !==============================================================================
1274 ! ####################
1275  SUBROUTINE read_netcdf_sst(HFILENAME,HNCVARNAME,PFIELD)
1276 ! ####################
1277 !
1278 USE modd_grid_latlonregul, ONLY : nindepth,nilength
1279 USE modd_surf_par, ONLY : xundef
1280 USE modd_csts, ONLY : xtt
1281 USE modd_prep, ONLY : cinterp_type
1282 !
1283 IMPLICIT NONE
1284 !
1285  CHARACTER(LEN=28), INTENT(IN) :: hfilename ! Name of the field file.
1286  CHARACTER(LEN=28), INTENT(IN) :: hncvarname ! Name of variable to read in netcdf file
1287 REAL, POINTER,DIMENSION(:) :: pfield ! value to get
1288 !
1289 REAL,DIMENSION(:), ALLOCATABLE :: zlati
1290 REAL,DIMENSION(:), ALLOCATABLE :: zlong
1291 REAL,TARGET,DIMENSION(:,:), ALLOCATABLE :: zvalue
1292 REAL,DIMENSION(:), ALLOCATABLE :: zdepth
1293 REAL,TARGET,DIMENSION(:), ALLOCATABLE :: zval
1294 integer :: jloop
1295 !PLM
1296 REAL :: zundef=999.
1297 REAL(KIND=JPRB) :: zhook_handle
1298 !
1299 !
1300 include 'netcdf.inc'
1301 !
1302 IF (lhook) CALL dr_hook('MODE_READ_NETCDF_MERCATOR:READ_NETCDF_SST',0,zhook_handle)
1303 IF (nilength<0) then
1304  ALLOCATE(pfield(1))
1305  pfield(:)=xundef
1306  cinterp_type='UNIF ' !!prescribed uniform field
1307 ELSEIF (nindepth<0) THEN
1308  ALLOCATE(zlati(nilength) )
1309  ALLOCATE(zlong(nilength) )
1310  ALLOCATE(zval(nilength) )
1311  CALL read_latlonval_cdf(hfilename,hncvarname,zlong,zlati,zval)
1312  ALLOCATE(pfield(nilength))
1313  pfield(:)=xundef
1314  pfield(:) = zval(:)
1315  IF (trim(hncvarname) == 'temperature') THEN
1316  WHERE (zval(:)/=zundef .AND. zval(:)<100.) pfield(:)=pfield(:)+xtt
1317  ENDIF
1318  cinterp_type='HORIBL' !!interpolation from gaussian, legendre or regular grid
1319 ! !!CINGRID_TYPE='GAUSS ' ou ='AROME '
1320 ! !!CINGRID_TYPE='LATLON '
1321 ELSE
1322  ALLOCATE(zvalue(nilength,nindepth))
1323  ALLOCATE(zlati(nilength) )
1324  ALLOCATE(zlong(nilength) )
1325  ALLOCATE(zdepth(nindepth))
1326 !
1327  CALL read_latlondepval_cdf(hfilename,hncvarname,zlong,zlati,zdepth,zvalue)
1328 !
1329  ALLOCATE(pfield(nilength))
1330  pfield(:)=xundef
1331  pfield(:)=zvalue(:,1)
1332  IF (trim(hncvarname) == 'temperature') THEN
1333  WHERE (zvalue(:,1)/=zundef .AND. zvalue(:,1)<100.) pfield(:)=pfield(:)+xtt
1334  ENDIF
1335  cinterp_type='HORIBL' !!interpolation from gaussian, legendre or regular grid
1336 ! !!CINGRID_TYPE='GAUSS ' ou ='AROME '
1337 ! !!CINGRID_TYPE='LATLON '
1338 ENDIF
1339 !
1340 IF (ALLOCATED(zvalue )) DEALLOCATE(zvalue )
1341 IF (ALLOCATED(zlong )) DEALLOCATE(zlong )
1342 IF (ALLOCATED(zlati )) DEALLOCATE(zlati )
1343 IF (ALLOCATED(zdepth )) DEALLOCATE(zdepth )
1344 IF (ALLOCATED(zval )) DEALLOCATE(zval )
1345 !
1346 202 FORMAT(3(3x,f10.4))
1347 IF (lhook) CALL dr_hook('MODE_READ_NETCDF_MERCATOR:READ_NETCDF_SST',1,zhook_handle)
1348 !
1349 END SUBROUTINE read_netcdf_sst
1350 !------------------------------------------------------------------------------
1351 !==============================================================================
1352 ! ####################
1353  SUBROUTINE read_netcdf_zs_sea(HFILENAME,HNCVARNAME,PFIELD)
1354 ! ####################
1355 !
1356 USE modd_grid_latlonregul, ONLY : ninlat,ninlon,nindepth,nilength
1357 USE modd_prep, ONLY : cinterp_type
1358 !
1359 IMPLICIT NONE
1360 !
1361  CHARACTER(LEN=28), INTENT(IN) :: hfilename ! Name of the field file.
1362  CHARACTER(LEN=28), INTENT(IN) :: hncvarname ! Name of variable to read in netcdf file
1363 REAL, POINTER, DIMENSION(:) :: pfield ! value to get
1364 !
1365 REAL,DIMENSION(:), ALLOCATABLE :: zlati
1366 REAL,DIMENSION(:), ALLOCATABLE :: zlong
1367 REAL,TARGET, DIMENSION(:), ALLOCATABLE:: zvalue
1368 integer :: jloop
1369 REAL(KIND=JPRB) :: zhook_handle
1370 !
1371 !
1372 include 'netcdf.inc'
1373 !
1374 IF (lhook) CALL dr_hook('MODE_READ_NETCDF_MERCATOR:READ_NETCDF_ZS_SEA',0,zhook_handle)
1375 if(nindepth>0) then
1376  !write(0,*) '*****warning*****',HNCVARNAME,' is a 3D field'
1377  ALLOCATE(pfield(1))
1378  pfield(:)=0.
1379  cinterp_type='UNIF ' !!prescribed uniform field
1380 elseif(nilength>0) then
1381  ALLOCATE(zvalue(nilength))
1382  ALLOCATE(zlati(nilength) )
1383  ALLOCATE(zlong(nilength) )
1384 !
1385  CALL read_latlonval_cdf(hfilename,hncvarname,zlong,zlati,zvalue)
1386  ALLOCATE(pfield(nilength))
1387  pfield(:)=zvalue(:)
1388  cinterp_type='HORIBL' !!interpolation from gaussian, legendre or regular grid
1389 ! !!CINGRID_TYPE='GAUSS ' ou ='AROME '
1390 ! !!CINGRID_TYPE='LATLON '
1391 else
1392  ALLOCATE(pfield(1))
1393  pfield(:)=0.
1394  cinterp_type='UNIF ' !!prescribed uniform field
1395 endif
1396 !
1397 IF (ALLOCATED(zvalue )) DEALLOCATE(zvalue )
1398 IF (ALLOCATED(zlong )) DEALLOCATE(zlong )
1399 IF (ALLOCATED(zlati )) DEALLOCATE(zlati )
1400 !
1401 202 FORMAT(3(3x,f10.4))
1402 IF (lhook) CALL dr_hook('MODE_READ_NETCDF_MERCATOR:READ_NETCDF_ZS_SEA',1,zhook_handle)
1403 !
1404 END SUBROUTINE read_netcdf_zs_sea
1405 !------------------------------------------------------------------------------
1406 !==============================================================================
1407 END MODULE mode_read_netcdf_mercator
subroutine get1dcdf(KCDF_ID, IDVAR, PMISSVALUE, PVALU1D)
subroutine get3dcdf(KCDF_ID, IDVAR, PDIM1, HDIM1NAME, PDIM2, HDIM2NAME, PDIM3, HDIM3NAME, PMISSVALUE, PVALU3D)
subroutine read_netcdf_sst(HFILENAME, HNCVARNAME, PFIELD)
subroutine abor1_sfx(YTEXT)
Definition: abor1_sfx.F90:6
subroutine read_latlondepval_cdf(HFILENAME, HNCVARNAME, PLON, PLAT, PDEP, PVAL)
subroutine read_dim_cdf(HFILENAME, HNCVARNAME, KDIM)
subroutine prep_netcdf_grid(HFILENAME, HNCVARNAME)
subroutine read_latlonval_cdf(HFILENAME, HNCVARNAME, PLON, PLAT, PVAL)
subroutine read_z1d_cdf(HFILENAME, HNCVARNAME, PVAL)
subroutine handle_err_mer(status, line)
subroutine get2dcdf(KCDF_ID, IDVAR, PDIM1, HDIM1NAME, PDIM2, HDIM2NAME, PMISSVALUE, PVALU2D)
subroutine read_netcdf_zs_sea(HFILENAME, HNCVARNAME, PFIELD)