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