source: CONFIG_DEVT/IPSLCM6.5_work_ENSEMBLES/oasis3-mct/examples/toy_configuration_components_C/read_all_data.F90 @ 5725

Last change on this file since 5725 was 5725, checked in by aclsce, 3 years ago

Added new oasis3-MCT version to be used to handle ensembles simulations with XIOS.

File size: 16.5 KB
Line 
1MODULE read_all_data
2  !
3  USE netcdf
4  IMPLICIT NONE
5  !
6  !
7  CONTAINS
8!****************************************************************************************
9SUBROUTINE read_dimgrid(nlon, nlat, data_filename, cl_grd, w_unit, FILE_Debug)
10  !**************************************************************************************
11  !
12  INTEGER, INTENT(in)               :: w_unit, FILE_Debug
13  CHARACTER(len=*), INTENT(in)      :: data_filename
14  CHARACTER(len=4)                  :: cl_grd ! name of the grid
15  !               
16  INTEGER, INTENT(out)       :: nlon, nlat
17  !
18  ! Local variables to the routine
19  INTEGER                  :: il_file_id, il_lon_id, &
20                              il_lat_id, &
21                              lon_dims, lat_dims
22
23  INTEGER, DIMENSION(NF90_MAX_VAR_DIMS) :: lon_dims_ids, lat_dims_ids, &
24                                           lon_dims_len, lat_dims_len
25
26  logical                    :: exists
27  INTEGER                    :: i 
28  CHARACTER(len=8)           :: cl_nam ! cl_grd+.lon,+.lat ...
29  character(len=*),parameter :: subname = '(read_dimgrid)'
30  !
31  IF (FILE_Debug >= 2) THEN
32      WRITE(w_unit,*) 'Data_filename :',TRIM(data_filename)
33      CALL FLUSH(w_unit)
34  ENDIF
35  !
36  ! Dimensions
37  !
38  ! Check if file exists before open it
39  inquire(file=TRIM(data_filename),exist=exists)
40  if (exists .eqv. .FALSE. ) then
41     write(w_unit,*) 'File ',trim(data_filename),' does not exists'
42     call routine_model_abort(w_unit,__FILE__,__LINE__,subname)
43  endif
44
45  CALL hdlerr(NF90_OPEN(TRIM(data_filename), NF90_NOWRITE, il_file_id), w_unit, subname, __FILE__,__LINE__ )
46  !
47  cl_nam=TRIM(cl_grd)//".lon" 
48  IF (FILE_Debug >= 2) THEN
49      WRITE(w_unit,*) 'Longitudes :',TRIM(cl_nam)
50      CALL FLUSH(w_unit)
51  ENDIF
52  CALL hdlerr( NF90_INQ_VARID(il_file_id, TRIM(cl_nam),  il_lon_id), w_unit, subname, __FILE__, __LINE__ )
53
54  CALL hdlerr( NF90_INQUIRE_VARIABLE(il_file_id, varid=il_lon_id, ndims=lon_dims, dimids=lon_dims_ids), w_unit, subname, __FILE__, __LINE__ )
55  !
56  ! The value lon_dims_len(i) is obtained thanks to the lon_dims_ids ID already obtained from the file
57  DO i=1,lon_dims
58    CALL hdlerr( NF90_INQUIRE_DIMENSION(ncid=il_file_id,dimid=lon_dims_ids(i),len=lon_dims_len(i)), w_unit, subname, __FILE__, __LINE__ )
59  ENDDO
60  !
61  nlon=lon_dims_len(1)
62  nlat=lon_dims_len(2)
63  !++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
64  !
65  cl_nam=TRIM(cl_grd)//".lat" 
66  IF (FILE_Debug >= 2) THEN
67      WRITE(w_unit,*) 'Latitudes :',TRIM(cl_nam)
68      CALL FLUSH(w_unit)
69  ENDIF
70  CALL hdlerr( NF90_INQ_VARID(il_file_id, TRIM(cl_nam),  il_lat_id), w_unit, subname, __FILE__,   __LINE__ )
71  !
72  CALL hdlerr( NF90_INQUIRE_VARIABLE(ncid=il_file_id, varid=il_lat_id, ndims=lat_dims, dimids=lat_dims_ids), w_unit, subname, __FILE__, __LINE__ )
73  !
74  ! The value lat_dims_len(i) is obtained thanks to the lat_dims_ids ID already obtained from the file
75  DO i=1,lat_dims
76    CALL hdlerr( NF90_INQUIRE_DIMENSION(ncid=il_file_id,dimid=lat_dims_ids(i),len=lat_dims_len(i)), w_unit, subname, __FILE__, __LINE__ )
77  ENDDO
78  !
79  IF ( (lat_dims_len(1) .NE. lon_dims_len(1)).OR.(lat_dims_len(2) .NE. lon_dims_len(2)) ) THEN
80      WRITE(w_unit,*) 'Problem model1 in read_dimgrid'
81      WRITE(w_unit,*) 'Dimensions of the latitude are not the same as the ones of the longitude'
82      CALL flush(w_unit)
83      STOP
84  ENDIF
85  !
86  CALL hdlerr(NF90_CLOSE(il_file_id),  w_unit, subname, __FILE__, __LINE__ )
87  !
88  IF (FILE_Debug >= 2) THEN
89      WRITE(w_unit,*) 'Reading input file ',TRIM(data_filename)
90      WRITE(w_unit,*) 'Global dimensions nlon= ',nlon,' nlat= ',nlat
91      CALL FLUSH(w_unit)
92  ENDIF
93  !
94  !
95END SUBROUTINE read_dimgrid
96
97  !****************************************************************************************
98  SUBROUTINE read_grid(nlon, nlat, id_begi, id_begj, id_lon, id_lat, &
99                        cl_grd, data_filename, w_unit, FILE_Debug,  &
100                        dda_lon, dda_lat)
101  !**************************************************************************************
102  !
103  INTEGER, PARAMETER :: wp = SELECTED_REAL_KIND(12,307) ! double
104  !
105  INTEGER, INTENT(in)            :: w_unit, FILE_Debug                                   
106  INTEGER, INTENT(in)            :: nlon, nlat, id_begi, id_begj, id_lon, id_lat
107  CHARACTER(len=*), INTENT(in)   :: data_filename
108  CHARACTER(len=*), INTENT(in)   :: cl_grd ! name of the grid
109  !
110  REAL (kind=wp), DIMENSION(id_lon, id_lat), INTENT(out)  :: dda_lon, dda_lat
111  !
112  ! Local variables to the routine
113  CHARACTER(len=8)           :: cl_nam ! cl_grd+.lon,+.lat ...
114  character(len=*),parameter :: subname = '(read_grid)'
115  logical                    :: exists
116  INTEGER                    :: il_file_id, il_lon_id, il_lat_id 
117  INTEGER,  DIMENSION(2)     :: ila_dim
118  INTEGER,  DIMENSION(2)     :: ila_st
119  !
120  !
121  IF (FILE_Debug >= 2) THEN
122      WRITE(w_unit,*) 'Data_filename :',data_filename
123      CALL FLUSH(w_unit)
124  ENDIF
125  !
126  ! Check if file exists before open it
127  inquire(file=TRIM(data_filename),exist=exists)
128  if (exists .eqv. .FALSE. ) then
129     write(w_unit,*) 'File ',TRIM(data_filename),' does not exists'
130     call routine_model_abort(w_unit,__FILE__,__LINE__,subname)
131  endif
132
133  CALL hdlerr(NF90_OPEN(TRIM(data_filename), NF90_NOWRITE, il_file_id), w_unit, subname, __FILE__,__LINE__ )
134  !
135  cl_nam=TRIM(cl_grd)//".lon" 
136  IF (FILE_Debug >= 2) THEN
137      WRITE(w_unit,*) 'Longitudes :',TRIM(cl_nam)
138      CALL FLUSH(w_unit)
139  ENDIF
140  CALL hdlerr( NF90_INQ_VARID(il_file_id, TRIM(cl_nam),  il_lon_id),  w_unit, subname, __FILE__, __LINE__ )
141
142  cl_nam=TRIM(cl_grd)//".lat" 
143  IF (FILE_Debug >= 2) THEN
144      WRITE(w_unit,*) 'Latitudes :',TRIM(cl_nam)
145      CALL FLUSH(w_unit)
146  ENDIF
147  CALL hdlerr( NF90_INQ_VARID(il_file_id, TRIM(cl_nam),  il_lat_id),  w_unit, subname, __FILE__, __LINE__ )
148  !
149  ila_st(1) = id_begi
150  ila_st(2) = id_begj
151  !
152  ila_dim(1) = id_lon
153  ila_dim(2) = id_lat
154  !
155  ! Data
156  !
157  CALL hdlerr( NF90_GET_VAR (il_file_id, il_lon_id, dda_lon, &
158     ila_st, ila_dim),  w_unit, subname, __FILE__, __LINE__ )
159  IF (FILE_Debug >= 2) THEN
160      WRITE(w_unit,*) '... local grid longitudes reading done'
161      CALL FLUSH(w_unit)
162  ENDIF
163  !
164  CALL hdlerr( NF90_GET_VAR (il_file_id, il_lat_id, dda_lat, &
165     ila_st, ila_dim),  w_unit, subname, __FILE__, __LINE__ )
166  IF (FILE_Debug >= 2) THEN
167      WRITE(w_unit,*) '... local grid latitudes reading done'
168      CALL FLUSH(w_unit)
169  ENDIF
170  !
171  !
172  CALL hdlerr( NF90_CLOSE(il_file_id),  w_unit, subname, __FILE__, __LINE__ )
173  !
174  IF (FILE_Debug >= 2) THEN
175      WRITE(w_unit,*) 'End of routine read_grid'
176      CALL FLUSH(w_unit)
177  ENDIF
178  !
179END SUBROUTINE read_grid
180
181!****************************************************************************************
182LOGICAL FUNCTION inquire_mask(cl_grd, data_filename, w_unit, FILE_Debug)
183!**************************************************************************************
184   !
185   INTEGER, INTENT(in)           :: w_unit, FILE_Debug         
186   CHARACTER(len=*), INTENT(in)  :: data_filename
187   CHARACTER(len=*), INTENT(in)  :: cl_grd ! name of the source grid
188
189   ! local variables to the routine
190   CHARACTER(len=8)          :: cl_nam ! cl_grd+.msk
191   character(len=*),parameter :: subname = '(inquire_mask)'
192   INTEGER                    :: il_file_id, il_msk_id
193   logical                    :: exists 
194   !
195   ! Check if file exists before open it
196   inquire(file=TRIM(data_filename),exist=exists)
197   if (exists .eqv. .FALSE. ) then
198      write(w_unit,*) 'File ',TRIM(data_filename),' does not exists'
199      call routine_model_abort(w_unit,__FILE__,__LINE__,subname)
200   endif
201
202   CALL hdlerr(NF90_OPEN(data_filename, NF90_NOWRITE, il_file_id), w_unit, subname, __FILE__,__LINE__ )
203   !
204   !
205   cl_nam=TRIM(cl_grd)//".msk" 
206   inquire_mask =  NF90_INQ_VARID(il_file_id, TRIM(cl_nam),  il_msk_id) == NF90_NOERR
207   !
208   CALL hdlerr( NF90_CLOSE(il_file_id),  w_unit, subname, __FILE__, __LINE__ )
209   !
210   IF (FILE_Debug >= 2) THEN
211      WRITE(w_unit,*) 'End of function inquire_mask for ',TRIM(cl_grd)
212      CALL FLUSH(w_unit)
213   ENDIF
214   !
215END FUNCTION inquire_mask
216
217!****************************************************************************************
218LOGICAL FUNCTION inquire_frac(cl_grd, data_filename, w_unit, FILE_Debug)
219!**************************************************************************************
220   !
221   INTEGER, INTENT(in)           :: w_unit, FILE_Debug         
222   CHARACTER(len=*), INTENT(in)  :: data_filename
223   CHARACTER(len=*), INTENT(in)  :: cl_grd ! name of the source grid
224
225   ! local variables to the routine
226   CHARACTER(len=8)          :: cl_nam ! cl_grd+.frc
227   character(len=*),parameter :: subname = '(inquire_mask)'
228   INTEGER                    :: il_file_id, il_frc_id
229   logical                    :: exists 
230   !
231   ! Check if file exists before open it
232   inquire(file=TRIM(data_filename),exist=exists)
233   if (exists .eqv. .FALSE. ) then
234      write(w_unit,*) 'File ',TRIM(data_filename),' does not exists'
235      call routine_model_abort(w_unit,__FILE__,__LINE__,subname)
236   endif
237
238   CALL hdlerr(NF90_OPEN(data_filename, NF90_NOWRITE, il_file_id), w_unit, subname, __FILE__,__LINE__ )
239   !
240   !
241   cl_nam=TRIM(cl_grd)//".frc" 
242   inquire_frac =  NF90_INQ_VARID(il_file_id, TRIM(cl_nam),  il_frc_id) == NF90_NOERR
243   !
244   CALL hdlerr( NF90_CLOSE(il_file_id),  w_unit, subname, __FILE__, __LINE__ )
245   !
246   IF (FILE_Debug >= 2) THEN
247      WRITE(w_unit,*) 'End of function inquire_frac for ',TRIM(cl_grd)
248      CALL FLUSH(w_unit)
249   ENDIF
250   !
251END FUNCTION inquire_frac
252
253  !****************************************************************************************
254  SUBROUTINE read_mask(nlon, nlat, id_begi, id_begj, id_lon, id_lat, &
255                        cl_grd, data_filename, w_unit, FILE_Debug,  &
256                        ida_msk)
257  !**************************************************************************************
258  !
259  INTEGER, PARAMETER :: wp = SELECTED_REAL_KIND(12,307) ! double
260  !
261  INTEGER, INTENT(in)            :: w_unit, FILE_Debug                                   
262  INTEGER, INTENT(in)            :: nlon, nlat, id_begi, id_begj, id_lon, id_lat
263  CHARACTER(len=*), INTENT(in)   :: data_filename
264  CHARACTER(len=*), INTENT(in)   :: cl_grd ! name of the grid
265  !
266  INTEGER, DIMENSION(id_lon, id_lat), INTENT(out)      :: ida_msk
267  !REAL (kind=wp), DIMENSION(id_lon, id_lat), optional  :: ida_frc
268  !
269  ! Local variables to the routine
270  CHARACTER(len=8)           :: cl_nam ! cl_grd+.lon,+.lat ...
271  character(len=*),parameter :: subname = '(read_mask)'
272  logical                    :: exists, frc_exists
273  INTEGER                    :: il_file_id, il_msk_id
274  INTEGER,  DIMENSION(2)     :: ila_dim
275  INTEGER,  DIMENSION(2)     :: ila_st
276  !
277  !
278  IF (FILE_Debug >= 2) THEN
279      WRITE(w_unit,*) 'Data_filename :',TRIM(data_filename)
280      CALL FLUSH(w_unit)
281  ENDIF
282  !
283  ! Check if file exists before open it
284  inquire(file=TRIM(data_filename),exist=exists)
285  if (exists .eqv. .FALSE. ) then
286     write(w_unit,*) 'File ',TRIM(data_filename),' does not exists'
287     call routine_model_abort(w_unit,__FILE__,__LINE__,subname)
288  endif
289
290  CALL hdlerr(NF90_OPEN(data_filename, NF90_NOWRITE, il_file_id), w_unit, subname, __FILE__,__LINE__ )
291  !
292  !
293  cl_nam=TRIM(cl_grd)//".msk" 
294  CALL hdlerr( NF90_INQ_VARID(il_file_id, TRIM(cl_nam),  il_msk_id),  w_unit, subname, __FILE__, __LINE__ )
295  !
296  ila_st(1) = id_begi
297  ila_st(2) = id_begj
298  !
299  ila_dim(1) = id_lon
300  ila_dim(2) = id_lat
301  !
302  ! Data
303  !
304  CALL hdlerr( NF90_GET_VAR (il_file_id, il_msk_id, ida_msk, &
305     ila_st, ila_dim),  w_unit, subname, __FILE__, __LINE__ )
306  IF (FILE_Debug >= 2) THEN
307      WRITE(w_unit,*) '... local mask reading done'
308      CALL FLUSH(w_unit)
309  ENDIF
310  !
311  !
312  CALL hdlerr( NF90_CLOSE(il_file_id),  w_unit, subname, __FILE__, __LINE__ )
313  !
314  IF (FILE_Debug >= 2) THEN
315      WRITE(w_unit,*) 'End of routine read_mask and eventually fraction of mask over ocean'
316      CALL FLUSH(w_unit)
317  ENDIF
318  !
319END SUBROUTINE read_mask
320
321  !****************************************************************************************
322  SUBROUTINE read_frac(nlon, nlat, id_begi, id_begj, id_lon, id_lat, &
323                        cl_grd, data_filename, w_unit, FILE_Debug,  &
324                        dda_frc)
325  !****************************************************************************************
326  !
327  INTEGER, PARAMETER :: wp = SELECTED_REAL_KIND(12,307) ! double
328  !
329  INTEGER, INTENT(in)            :: w_unit, FILE_Debug                                   
330  INTEGER, INTENT(in)            :: nlon, nlat, id_begi, id_begj, id_lon, id_lat
331  CHARACTER(len=*), INTENT(in)   :: data_filename  ! masks.nc
332  CHARACTER(len=*), INTENT(in)   :: cl_grd ! name of the grid
333  !
334  REAL (kind=wp), DIMENSION(id_lon, id_lat), INTENT(out)  :: dda_frc
335  !
336  ! Local variables to the routine
337  CHARACTER(len=8)           :: cl_nam ! cl_grd+.frc
338  character(len=*),parameter :: subname = '(read_frac)'
339  logical                    :: exists
340  INTEGER                    :: il_file_id, il_frc_id 
341  INTEGER,  DIMENSION(2)     :: ila_dim
342  INTEGER,  DIMENSION(2)     :: ila_st
343  !
344  !
345  IF (FILE_Debug >= 2) THEN
346      WRITE(w_unit,*) 'Data_filename :',TRIM(data_filename)
347      CALL FLUSH(w_unit)
348  ENDIF
349  !
350  !
351  ! Check if file exists before open it
352  inquire(file=TRIM(data_filename),exist=exists)
353  if (exists .eqv. .FALSE. ) then
354     write(w_unit,*) 'File ',TRIM(data_filename),' does not exists'
355     call routine_model_abort(w_unit,__FILE__,__LINE__,subname)
356  endif
357
358  CALL hdlerr(NF90_OPEN(TRIM(data_filename), NF90_NOWRITE, il_file_id), w_unit, subname, __FILE__,__LINE__ )
359  !
360  cl_nam=TRIM(cl_grd)//".frc" 
361  CALL hdlerr( NF90_INQ_VARID(il_file_id, TRIM(cl_nam), il_frc_id),  w_unit, subname, __FILE__, __LINE__ )
362  !
363  ila_st(1) = id_begi
364  ila_st(2) = id_begj
365  !
366  ila_dim(1) = id_lon
367  ila_dim(2) = id_lat
368  !
369  ! Data
370  !
371  CALL hdlerr( NF90_GET_VAR (il_file_id, il_frc_id, dda_frc, &
372     ila_st, ila_dim),  w_unit, subname, __FILE__, __LINE__ )
373  IF (FILE_Debug >= 2) THEN
374      WRITE(w_unit,*) '... local area reading done'
375      CALL FLUSH(w_unit)
376  ENDIF
377  !
378  CALL hdlerr( NF90_CLOSE(il_file_id),  w_unit, subname, __FILE__, __LINE__ )
379  !
380  IF (FILE_Debug >= 2) THEN
381      WRITE(w_unit,*) 'End of routine read_frac'
382      CALL FLUSH(w_unit)
383  ENDIF
384  !
385END SUBROUTINE read_frac
386
387  !****************************************************************************************
388  SUBROUTINE read_area (nlon, nlat, id_begi, id_begj, id_lon, id_lat, &
389                        cl_grd, data_filename, w_unit, FILE_Debug,  &
390                        dda_srf)
391  !****************************************************************************************
392  !
393  INTEGER, PARAMETER :: wp = SELECTED_REAL_KIND(12,307) ! double
394  !
395  INTEGER, INTENT(in)            :: w_unit, FILE_Debug                                   
396  INTEGER, INTENT(in)            :: nlon, nlat, id_begi, id_begj, id_lon, id_lat
397  CHARACTER(len=*), INTENT(in)   :: data_filename
398  CHARACTER(len=*), INTENT(in)   :: cl_grd ! name of the grid
399  !
400  REAL (kind=wp), DIMENSION(id_lon, id_lat), INTENT(out)  :: dda_srf
401  !
402  ! Local variables to the routine
403  CHARACTER(len=8)           :: cl_nam ! cl_grd+.srf
404  character(len=*),parameter :: subname = '(read_area)'
405  logical                    :: exists
406  INTEGER                    :: il_file_id, il_srf_id 
407  INTEGER,  DIMENSION(2)     :: ila_dim
408  INTEGER,  DIMENSION(2)     :: ila_st
409  !
410  !
411  IF (FILE_Debug >= 2) THEN
412      WRITE(w_unit,*) 'Data_filename :',TRIM(data_filename)
413      CALL FLUSH(w_unit)
414  ENDIF
415  !
416  !
417  ! Check if file exists before open it
418  inquire(file=TRIM(data_filename),exist=exists)
419  if (exists .eqv. .FALSE. ) then
420     write(w_unit,*) 'File ',TRIM(data_filename),' does not exists'
421     call routine_model_abort(w_unit,__FILE__,__LINE__,subname)
422  endif
423
424  CALL hdlerr(NF90_OPEN(TRIM(data_filename), NF90_NOWRITE, il_file_id), w_unit, subname, __FILE__,__LINE__ )
425  !
426  cl_nam=TRIM(cl_grd)//".srf" 
427  CALL hdlerr( NF90_INQ_VARID(il_file_id, TRIM(cl_nam), il_srf_id),  w_unit, subname, __FILE__, __LINE__ )
428  !
429  ila_st(1) = id_begi
430  ila_st(2) = id_begj
431  !
432  ila_dim(1) = id_lon
433  ila_dim(2) = id_lat
434  !
435  ! Data
436  !
437  CALL hdlerr( NF90_GET_VAR (il_file_id, il_srf_id, dda_srf, &
438     ila_st, ila_dim),  w_unit, subname, __FILE__, __LINE__ )
439  IF (FILE_Debug >= 2) THEN
440      WRITE(w_unit,*) '... local area reading done'
441      CALL FLUSH(w_unit)
442  ENDIF
443  !
444  CALL hdlerr( NF90_CLOSE(il_file_id),  w_unit, subname, __FILE__, __LINE__ )
445  !
446  IF (FILE_Debug >= 2) THEN
447      WRITE(w_unit,*) 'End of routine read_area'
448      CALL FLUSH(w_unit)
449  ENDIF
450  !
451END SUBROUTINE read_area
452
453!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
454END MODULE read_all_data
Note: See TracBrowser for help on using the repository browser.