source: CPL/oasis3-mct/branches/OASIS3-MCT_2.0_branch/examples/test_interpolation/read_all_data.F90 @ 4775

Last change on this file since 4775 was 4775, checked in by aclsce, 4 years ago
  • Imported oasis3-mct from Cerfacs svn server (not suppotred anymore).

The version has been extracted from https://oasis3mct.cerfacs.fr/svn/branches/OASIS3-MCT_2.0_branch/oasis3-mct@1818

File size: 8.6 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  CHARACTER(len=30)        :: data_filename
27  CHARACTER(len=4)         :: cl_grd ! name of the source grid
28  CHARACTER(len=8)         :: cl_nam ! cl_grd+.lon,+.lat ...
29  !
30  IF (FILE_Debug >= 2) THEN
31      WRITE(w_unit,*) 'Data_filename :',data_filename
32      CALL FLUSH(w_unit)
33  ENDIF
34  !
35  ! Dimensions
36  !
37  CALL hdlerr(NF90_OPEN(data_filename, NF90_NOWRITE, il_file_id), __LINE__ )
38  !
39  cl_nam=cl_grd//".lon" 
40  IF (FILE_Debug >= 2) THEN
41      WRITE(w_unit,*) 'Longitudes :',cl_nam
42      CALL FLUSH(w_unit)
43  ENDIF
44  CALL hdlerr( NF90_INQ_VARID(il_file_id, cl_nam,  il_lon_id),    __LINE__ )
45  cl_nam=cl_grd//".lat" 
46  IF (FILE_Debug >= 2) THEN
47      WRITE(w_unit,*) 'Latitudes :',cl_nam
48      CALL FLUSH(w_unit)
49  ENDIF
50  CALL hdlerr( NF90_INQ_VARID(il_file_id, cl_nam,  il_lat_id),    __LINE__ )
51
52  CALL hdlerr( NF90_INQUIRE_VARIABLE(il_file_id, varid=il_lon_id, ndims=lon_dims, dimids=lon_dims_ids), __LINE__ )
53  !
54  ! The value lon_dims_len(i) is obtained thanks to the lon_dims_ids ID already obtained from the file
55  DO i=1,lon_dims
56    CALL hdlerr( NF90_INQUIRE_DIMENSION(ncid=il_file_id,dimid=lon_dims_ids(i),len=lon_dims_len(i)), __LINE__ )
57  ENDDO
58  !
59  nlon=lon_dims_len(1)
60  nlat=lon_dims_len(2)
61  !++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
62  !
63  CALL hdlerr( NF90_INQUIRE_VARIABLE(ncid=il_file_id, varid=il_lat_id, ndims=lat_dims, dimids=lat_dims_ids), __LINE__ )
64  !
65  ! The value lat_dims_len(i) is obtained thanks to the lat_dims_ids ID already obtained from the file
66  DO i=1,lat_dims
67    CALL hdlerr( NF90_INQUIRE_DIMENSION(ncid=il_file_id,dimid=lat_dims_ids(i),len=lat_dims_len(i)), __LINE__ )
68  ENDDO
69  !
70  IF ( (lat_dims_len(1) .NE. lon_dims_len(1)).OR.(lat_dims_len(2) .NE. lon_dims_len(2)) ) THEN
71      WRITE(w_unit,*) 'Problem model1 in read_dimgrid'
72      WRITE(w_unit,*) 'Dimensions of the latitude are not the same as the ones of the longitude'
73      CALL flush(w_unit)
74      STOP
75  ENDIF
76  !
77  CALL hdlerr(NF90_CLOSE(il_file_id), __LINE__ )
78  !
79  IF (FILE_Debug >= 2) THEN
80      WRITE(w_unit,*) 'Reading input file ',data_filename
81      WRITE(w_unit,*) 'Global dimensions nlon=',nlon,' nlat=',nlat
82      CALL FLUSH(w_unit)
83  ENDIF
84  !
85  !
86END SUBROUTINE read_dimgrid
87
88  !****************************************************************************************
89  SUBROUTINE read_grid (nlon,nlat,&
90                        data_filename, cl_grd, w_unit, FILE_Debug, &
91                        gridlon,gridlat)
92  !**************************************************************************************
93  !
94#ifdef NO_USE_DOUBLE_PRECISION
95  INTEGER, PARAMETER :: wp = SELECTED_REAL_KIND(6,37)   ! real
96#elif defined USE_DOUBLE_PRECISION
97  INTEGER, PARAMETER :: wp = SELECTED_REAL_KIND(12,307) ! double
98#endif
99  !
100  INTEGER                  :: i,j,k,w_unit,FILE_Debug
101  !
102  INTEGER                  :: il_file_id,il_lon_id, il_lat_id                                     
103  !               
104  INTEGER, INTENT(in)     :: nlon,nlat
105  !
106  CHARACTER(len=30)        :: data_filename
107  CHARACTER(len=4)         :: cl_grd ! name of the source grid
108  CHARACTER(len=8)         :: cl_nam ! cl_grd+.lon,+.lat ...
109  !
110  INTEGER,  DIMENSION(2)   :: ila_dim
111  INTEGER,  DIMENSION(2)   :: ila_what
112  !
113  REAL (kind=wp), DIMENSION(nlon,nlat)  :: gridlon,gridlat
114  !
115  !
116  IF (FILE_Debug >= 2) THEN
117      WRITE(w_unit,*) 'Data_filename :',data_filename
118      CALL FLUSH(w_unit)
119  ENDIF
120  !
121  ! Dimensions
122  !
123  CALL hdlerr(NF90_OPEN(data_filename, NF90_NOWRITE, il_file_id), __LINE__ )
124  !
125  cl_nam=cl_grd//".lon" 
126  IF (FILE_Debug >= 2) THEN
127      WRITE(w_unit,*) 'Longitudes :',cl_nam
128      CALL FLUSH(w_unit)
129  ENDIF
130  CALL hdlerr( NF90_INQ_VARID(il_file_id, cl_nam,  il_lon_id),    __LINE__ )
131  cl_nam=cl_grd//".lat" 
132  IF (FILE_Debug >= 2) THEN
133      WRITE(w_unit,*) 'Latitudes :',cl_nam
134      CALL FLUSH(w_unit)
135  ENDIF
136  CALL hdlerr( NF90_INQ_VARID(il_file_id, cl_nam,  il_lat_id),    __LINE__ )
137  !
138  ila_what(:)=1
139  !
140  ila_dim(1)=nlon
141  ila_dim(2)=nlat
142  !
143  ! Data
144  !
145  !
146  CALL hdlerr( NF90_GET_VAR (il_file_id, il_lon_id, gridlon, &
147     ila_what, ila_dim), __LINE__ )
148  IF (FILE_Debug >= 2) THEN
149      WRITE(w_unit,*) 'Global grid longitudes reading done'
150      CALL FLUSH(w_unit)
151  ENDIF
152  !
153  CALL hdlerr( NF90_GET_VAR (il_file_id, il_lat_id, gridlat, &
154     ila_what, ila_dim), __LINE__ )
155  IF (FILE_Debug >= 2) THEN
156      WRITE(w_unit,*) 'Global grid latitudes reading done'
157      CALL FLUSH(w_unit)
158  ENDIF
159  !
160  !
161  CALL hdlerr( NF90_CLOSE(il_file_id), __LINE__ )
162  !
163  IF (FILE_Debug >= 2) THEN
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_mask (nlon, nlat, &
172                        data_filename, cl_grd, w_unit, FILE_Debug, &
173                        indice_mask)
174  !**************************************************************************************
175  !
176  INTEGER                  :: i,j,k,w_unit,FILE_Debug
177  !
178  INTEGER                  :: il_file_id, il_indice_id
179  !               
180  INTEGER, INTENT(in)     :: nlon,nlat
181  !
182  CHARACTER(len=30)        :: data_filename
183  CHARACTER(len=4)         :: cl_grd ! name of the source grid
184  CHARACTER(len=8)         :: cl_nam ! cl_grd+.lon,+.lat ...
185  !
186  INTEGER,  DIMENSION(2)   :: ila_dim
187  INTEGER,  DIMENSION(2)   :: ila_what
188  !
189  INTEGER, DIMENSION(nlon,nlat)  :: indice_mask
190  !
191  !
192  ! Dimensions
193  !
194  CALL hdlerr(NF90_OPEN(data_filename, NF90_NOWRITE, il_file_id), __LINE__ )
195  !
196  cl_nam=cl_grd//".msk" 
197  CALL hdlerr( NF90_INQ_VARID(il_file_id, cl_nam,  il_indice_id),    __LINE__ )
198  !
199  ila_what(:)=1
200  !
201  ila_dim(1)=nlon
202  ila_dim(2)=nlat
203  !
204  ! Data
205  !
206  CALL hdlerr( NF90_GET_VAR (il_file_id, il_indice_id, indice_mask, &
207     ila_what, ila_dim), __LINE__ )
208  IF (FILE_Debug >= 2) THEN
209      WRITE(w_unit,*) 'Global grid mask reading done'
210      CALL FLUSH(w_unit)
211  ENDIF
212  !
213  CALL hdlerr( NF90_CLOSE(il_file_id), __LINE__ )
214  !
215  IF (FILE_Debug >= 2) THEN
216      WRITE(w_unit,*) 'End of routine read_mask'
217      CALL FLUSH(w_unit)
218  ENDIF
219  !
220END SUBROUTINE read_mask
221
222  !****************************************************************************************
223  SUBROUTINE read_area (nlon,nlat, &
224                        data_filename, cl_grd, w_unit, FILE_Debug, &
225                        gridsrf)
226  !****************************************************************************************
227  !
228#ifdef NO_USE_DOUBLE_PRECISION
229  INTEGER, PARAMETER :: wp = SELECTED_REAL_KIND(6,37)   ! real
230#elif defined USE_DOUBLE_PRECISION
231  INTEGER, PARAMETER :: wp = SELECTED_REAL_KIND(12,307) ! double
232#endif
233  !
234  INTEGER                  :: i,j,k,w_unit,FILE_Debug
235  !
236  INTEGER                  :: il_file_id,  il_srf_id
237  !               
238  INTEGER, INTENT(in)     :: nlon,nlat
239  !
240  CHARACTER(len=30)        :: data_filename
241  CHARACTER(len=4)         :: cl_grd ! name of the source grid
242  CHARACTER(len=8)         :: cl_nam ! cl_grd+.lon,+.lat ...
243  !
244  INTEGER,  DIMENSION(2)   :: ila_dim
245  INTEGER,  DIMENSION(2)   :: ila_what
246  !
247  REAL (kind=wp), DIMENSION(nlon,nlat)  :: gridsrf
248  !
249  !
250  ! Dimensions
251  !
252  CALL hdlerr(NF90_OPEN(data_filename, NF90_NOWRITE, il_file_id), __LINE__ )
253  !
254  cl_nam=cl_grd//".srf" 
255  CALL hdlerr( NF90_INQ_VARID(il_file_id, cl_nam, il_srf_id),    __LINE__ )
256  !
257  ila_what(:)=1
258  !
259  ila_dim(1)=nlon
260  ila_dim(2)=nlat
261  !
262  ! Data
263  !
264  CALL hdlerr( NF90_GET_VAR (il_file_id, il_srf_id, gridsrf, &
265     ila_what, ila_dim), __LINE__ )
266  IF (FILE_Debug >= 2) THEN
267      WRITE(w_unit,*) 'Global grid mask reading done'
268      CALL FLUSH(w_unit)
269  ENDIF
270  !
271  CALL hdlerr( NF90_CLOSE(il_file_id), __LINE__ )
272  !
273  IF (FILE_Debug >= 2) THEN
274      WRITE(w_unit,*) 'End of routine read_mask'
275      CALL FLUSH(w_unit)
276  ENDIF
277  !
278END SUBROUTINE read_area
279
280!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
281END MODULE read_all_data
Note: See TracBrowser for help on using the repository browser.