source: CONFIG_DEVT/IPSLCM6.5_work_ENSEMBLES/oasis3-mct/examples/test_interpolation/read_all_data.F90 @ 5863

Last change on this file since 5863 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: 11.7 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  USE netcdf
12  IMPLICIT NONE
13  !
14  INTEGER                  :: i,j,k,w_unit,FILE_Debug
15  !
16  INTEGER                  :: il_file_id,il_grid_id,il_lon_id, &
17     il_lat_id,il_indice_id, &
18     lon_dims,lat_dims,imask_dims
19  !
20  INTEGER, DIMENSION(NF90_MAX_VAR_DIMS) :: lon_dims_ids,lat_dims_ids,&
21     imask_dims_ids,lon_dims_len,&
22     lat_dims_len,imask_dims_len 
23  !               
24  INTEGER, INTENT(out)       :: nlon,nlat
25  !
26  logical                    :: exists
27  CHARACTER(len=30)          :: data_filename
28  CHARACTER(len=4)           :: cl_grd ! name of the source grid
29  CHARACTER(len=8)           :: cl_nam ! cl_grd+.lon,+.lat ...
30  character(len=*),parameter :: subname = '(read_dimgrid)'
31  !
32  IF (FILE_Debug >= 2) THEN
33      WRITE(w_unit,*) 'Data_filename :',data_filename
34      CALL FLUSH(w_unit)
35  ENDIF
36  !
37  ! Dimensions
38  !
39  ! Check if file exists before open it
40  inquire(file=trim(data_filename),exist=exists)
41  if (exists .eqv. .FALSE. ) then
42     write(w_unit,*) 'File ',trim(data_filename),' does not exists'
43     call routine_model_abort(w_unit,__FILE__,__LINE__,subname)
44  endif
45
46  CALL hdlerr(NF90_OPEN(data_filename, NF90_NOWRITE, il_file_id), w_unit, subname, __FILE__,__LINE__ )
47  !
48  cl_nam=cl_grd//".lon" 
49  IF (FILE_Debug >= 2) THEN
50      WRITE(w_unit,*) 'Longitudes :',cl_nam
51      CALL FLUSH(w_unit)
52  ENDIF
53  CALL hdlerr( NF90_INQ_VARID(il_file_id, cl_nam,  il_lon_id), w_unit, subname, __FILE__, __LINE__ )
54  cl_nam=cl_grd//".lat" 
55  IF (FILE_Debug >= 2) THEN
56      WRITE(w_unit,*) 'Latitudes :',cl_nam
57      CALL FLUSH(w_unit)
58  ENDIF
59  CALL hdlerr( NF90_INQ_VARID(il_file_id, cl_nam,  il_lat_id), w_unit, subname, __FILE__,   __LINE__ )
60
61  CALL hdlerr( NF90_INQUIRE_VARIABLE(il_file_id, varid=il_lon_id, ndims=lon_dims, dimids=lon_dims_ids), w_unit, subname, __FILE__, __LINE__ )
62  !
63  ! The value lon_dims_len(i) is obtained thanks to the lon_dims_ids ID already obtained from the file
64  DO i=1,lon_dims
65    CALL hdlerr( NF90_INQUIRE_DIMENSION(ncid=il_file_id,dimid=lon_dims_ids(i),len=lon_dims_len(i)), w_unit, subname, __FILE__, __LINE__ )
66  ENDDO
67  !
68  nlon=lon_dims_len(1)
69  nlat=lon_dims_len(2)
70  !++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
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 ',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,&
99                        data_filename, cl_grd, w_unit, FILE_Debug, &
100                        gridlon,gridlat)
101  !**************************************************************************************
102  !
103  INTEGER, PARAMETER :: wp = SELECTED_REAL_KIND(12,307) ! double
104  !
105  INTEGER                    :: i,j,k,w_unit,FILE_Debug
106  !
107  INTEGER                    :: il_file_id,il_lon_id, il_lat_id                                     
108  !               
109  INTEGER, INTENT(in)        :: nlon,nlat
110  !
111  logical                    :: exists
112  CHARACTER(len=30)          :: data_filename
113  CHARACTER(len=4)           :: cl_grd ! name of the source grid
114  CHARACTER(len=8)           :: cl_nam ! cl_grd+.lon,+.lat ...
115  character(len=*),parameter :: subname = '(read_grid)'
116  !
117  INTEGER,  DIMENSION(2)   :: ila_dim
118  INTEGER,  DIMENSION(2)   :: ila_what
119  !
120  REAL (kind=wp), DIMENSION(nlon,nlat)  :: gridlon,gridlat
121  !
122  !
123  IF (FILE_Debug >= 2) THEN
124      WRITE(w_unit,*) 'Data_filename :',data_filename
125      CALL FLUSH(w_unit)
126  ENDIF
127  !
128  ! Check if file exists before open it
129  inquire(file=trim(data_filename),exist=exists)
130  if (exists .eqv. .FALSE. ) then
131     write(w_unit,*) 'File ',trim(data_filename),' does not exists'
132     call routine_model_abort(w_unit,__FILE__,__LINE__,subname)
133  endif
134
135  CALL hdlerr(NF90_OPEN(data_filename, NF90_NOWRITE, il_file_id), w_unit, subname, __FILE__,__LINE__ )
136  !
137  !
138  cl_nam=cl_grd//".lon" 
139  IF (FILE_Debug >= 2) THEN
140      WRITE(w_unit,*) 'Longitudes :',cl_nam
141      CALL FLUSH(w_unit)
142  ENDIF
143  CALL hdlerr( NF90_INQ_VARID(il_file_id, cl_nam,  il_lon_id),  w_unit, subname, __FILE__, __LINE__ )
144  cl_nam=cl_grd//".lat" 
145  IF (FILE_Debug >= 2) THEN
146      WRITE(w_unit,*) 'Latitudes :',cl_nam
147      CALL FLUSH(w_unit)
148  ENDIF
149  CALL hdlerr( NF90_INQ_VARID(il_file_id, cl_nam,  il_lat_id),  w_unit, subname, __FILE__, __LINE__ )
150  !
151  ila_what(:)=1
152  !
153  ila_dim(1)=nlon
154  ila_dim(2)=nlat
155  !
156  ! Data
157  !
158  CALL hdlerr( NF90_GET_VAR (il_file_id, il_lon_id, gridlon, &
159     ila_what, ila_dim),  w_unit, subname, __FILE__, __LINE__ )
160  IF (FILE_Debug >= 2) THEN
161      WRITE(w_unit,*) '... global grid longitudes reading done'
162      CALL FLUSH(w_unit)
163  ENDIF
164  !
165  CALL hdlerr( NF90_GET_VAR (il_file_id, il_lat_id, gridlat, &
166     ila_what, ila_dim),  w_unit, subname, __FILE__, __LINE__ )
167  IF (FILE_Debug >= 2) THEN
168      WRITE(w_unit,*) '... global grid latitudes reading done'
169      CALL FLUSH(w_unit)
170  ENDIF
171  !
172  !
173  CALL hdlerr( NF90_CLOSE(il_file_id),  w_unit, subname, __FILE__, __LINE__ )
174  !
175  IF (FILE_Debug >= 2) THEN
176      WRITE(w_unit,*) 'End of routine read_grid'
177      CALL FLUSH(w_unit)
178  ENDIF
179  !
180END SUBROUTINE read_grid
181
182!****************************************************************************************
183LOGICAL FUNCTION inquire_mask (data_filename, cl_grd, w_unit, FILE_Debug)
184!**************************************************************************************
185   !
186   INTEGER                    :: w_unit,FILE_Debug
187   !
188   INTEGER                    :: il_file_id, il_indice_id
189   !
190   logical                    :: exists               
191   CHARACTER(len=30)          :: data_filename
192   CHARACTER(len=*)           :: cl_grd ! name of the source grid
193   CHARACTER(len=18)          :: cl_nam ! cl_grd+.msk
194   character(len=*),parameter :: subname = '(inquire_mask)'
195   !
196   ! Check if file exists before open it
197   inquire(file=trim(data_filename),exist=exists)
198   if (exists .eqv. .FALSE. ) then
199      write(w_unit,*) 'File ',trim(data_filename),' does not exists'
200      call routine_model_abort(w_unit,__FILE__,__LINE__,subname)
201   endif
202
203   CALL hdlerr(NF90_OPEN(data_filename, NF90_NOWRITE, il_file_id), w_unit, subname, __FILE__,__LINE__ )
204   !
205   !
206   cl_nam=TRIM(cl_grd)//".msk" 
207   inquire_mask =  NF90_INQ_VARID(il_file_id, TRIM(cl_nam),  il_indice_id) == NF90_NOERR
208   !
209   CALL hdlerr( NF90_CLOSE(il_file_id),  w_unit, subname, __FILE__, __LINE__ )
210   !
211   IF (FILE_Debug >= 2) THEN
212      WRITE(w_unit,*) 'End of function inquire_mask for ',TRIM(cl_grd)
213      CALL FLUSH(w_unit)
214   ENDIF
215   !
216END FUNCTION inquire_mask
217
218  !****************************************************************************************
219  SUBROUTINE read_mask (nlon, nlat, &
220                        data_filename, cl_grd, w_unit, FILE_Debug, &
221                        indice_mask)
222  !**************************************************************************************
223  !
224  INTEGER                    :: i,j,k,w_unit,FILE_Debug
225  !
226  INTEGER                    :: il_file_id, il_indice_id
227  !               
228  INTEGER, INTENT(in)        :: nlon,nlat
229  !
230  logical                    :: exists
231  CHARACTER(len=30)          :: data_filename
232  CHARACTER(len=*)           :: cl_grd ! name of the source grid
233  CHARACTER(len=18)          :: cl_nam ! cl_grd+.lon,+.lat ...
234  character(len=*),parameter :: subname = '(read_mask)'
235  !
236  INTEGER,  DIMENSION(2)   :: ila_dim
237  INTEGER,  DIMENSION(2)   :: ila_what
238  !
239  INTEGER, DIMENSION(nlon,nlat)  :: indice_mask
240  !
241  ! Check if file exists before open it
242  inquire(file=trim(data_filename),exist=exists)
243  if (exists .eqv. .FALSE. ) then
244     write(w_unit,*) 'File ',trim(data_filename),' does not exists'
245     call routine_model_abort(w_unit,__FILE__,__LINE__,subname)
246  endif
247
248  CALL hdlerr(NF90_OPEN(data_filename, NF90_NOWRITE, il_file_id), w_unit, subname, __FILE__,__LINE__ )
249  !
250  !
251  cl_nam=TRIM(cl_grd)//".msk" 
252  CALL hdlerr( NF90_INQ_VARID(il_file_id, TRIM(cl_nam),  il_indice_id),  w_unit, subname, __FILE__, __LINE__ )
253  !
254  ila_what(:)=1
255  !
256  ila_dim(1)=nlon
257  ila_dim(2)=nlat
258  !
259  ! Data
260  !
261  CALL hdlerr( NF90_GET_VAR (il_file_id, il_indice_id, indice_mask, &
262     ila_what, ila_dim),  w_unit, subname, __FILE__, __LINE__ )
263  IF (FILE_Debug >= 2) THEN
264      WRITE(w_unit,*) '... global grid mask reading done'
265      CALL FLUSH(w_unit)
266  ENDIF
267  !
268  CALL hdlerr( NF90_CLOSE(il_file_id),  w_unit, subname, __FILE__, __LINE__ )
269  !
270  IF (FILE_Debug >= 2) THEN
271      WRITE(w_unit,*) 'End of routine read_mask'
272      CALL FLUSH(w_unit)
273  ENDIF
274  !
275END SUBROUTINE read_mask
276
277  !****************************************************************************************
278  SUBROUTINE read_area (nlon,nlat, &
279                        data_filename, cl_grd, w_unit, FILE_Debug, &
280                        gridsrf)
281  !****************************************************************************************
282  !
283  INTEGER, PARAMETER :: wp = SELECTED_REAL_KIND(12,307) ! double
284  !
285  INTEGER                    :: i,j,k,w_unit,FILE_Debug
286  !
287  INTEGER                    :: il_file_id,  il_srf_id
288  !               
289  INTEGER, INTENT(in)        :: nlon,nlat
290  !
291  logical                    :: exists
292  CHARACTER(len=30)          :: data_filename
293  CHARACTER(len=4)           :: cl_grd ! name of the source grid
294  CHARACTER(len=8)           :: cl_nam ! cl_grd+.lon,+.lat ...
295  character(len=*),parameter :: subname = '(read_area)'
296  !
297  INTEGER,  DIMENSION(2)   :: ila_dim
298  INTEGER,  DIMENSION(2)   :: ila_what
299  !
300  REAL (kind=wp), DIMENSION(nlon,nlat)  :: gridsrf
301  !
302  !
303  ! Check if file exists before open it
304  inquire(file=trim(data_filename),exist=exists)
305  if (exists .eqv. .FALSE. ) then
306     write(w_unit,*) 'File ',trim(data_filename),' does not exists'
307     call routine_model_abort(w_unit,__FILE__,__LINE__,subname)
308  endif
309
310  CALL hdlerr(NF90_OPEN(data_filename, NF90_NOWRITE, il_file_id), w_unit, subname, __FILE__,__LINE__ )
311  !
312  cl_nam=cl_grd//".srf" 
313  CALL hdlerr( NF90_INQ_VARID(il_file_id, cl_nam, il_srf_id),  w_unit, subname, __FILE__, __LINE__ )
314  !
315  ila_what(:)=1
316  !
317  ila_dim(1)=nlon
318  ila_dim(2)=nlat
319  !
320  ! Data
321  !
322  CALL hdlerr( NF90_GET_VAR (il_file_id, il_srf_id, gridsrf, &
323     ila_what, ila_dim),  w_unit, subname, __FILE__, __LINE__ )
324  IF (FILE_Debug >= 2) THEN
325      WRITE(w_unit,*) '... global grid mask reading done'
326      CALL FLUSH(w_unit)
327  ENDIF
328  !
329  CALL hdlerr( NF90_CLOSE(il_file_id),  w_unit, subname, __FILE__, __LINE__ )
330  !
331  IF (FILE_Debug >= 2) THEN
332      WRITE(w_unit,*) 'End of routine read_mask'
333      CALL FLUSH(w_unit)
334  ENDIF
335  !
336END SUBROUTINE read_area
337
338!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
339END MODULE read_all_data
Note: See TracBrowser for help on using the repository browser.