source: CPL/oasis3-mct/branches/OASIS3-MCT_5.0_branch/examples/regrid_environment/src/read_all_data.F90 @ 6331

Last change on this file since 6331 was 6331, checked in by aclsce, 17 months ago

Moved oasis-mct_5.0 in oasis3-mct/branches directory.

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