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

Last change on this file since 4775 was 4775, checked in by aclsce, 5 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: 3.6 KB
Line 
1  !****************************************************************************************
2  SUBROUTINE read_grid (nlon,nlat,corners_ij_lus, &
3                                       data_filename, w_unit, &
4                                       gridlon,gridlat, &
5                                       gridclo,gridcla, &
6                                       gridsrf, &
7                                       indice_mask)
8  !**************************************************************************************
9  !
10  USE netcdf
11  IMPLICIT NONE
12  !
13  INTEGER                  :: i,j,k,w_unit
14  !
15  INTEGER                  :: il_file_id,il_grid_id,il_lon_id, &
16                              il_lat_id,il_clo_id,il_cla_id,il_srf_id,il_indice_id, &
17                              lon_dims,lat_dims,clo_dims,cla_dims,&
18                              imask_dims
19  !
20  INTEGER, DIMENSION(NF90_MAX_VAR_DIMS) :: lon_dims_ids,lat_dims_ids,clo_dims_ids,&
21                                           cla_dims_ids,imask_dims_ids,lon_dims_len,&
22                                           lat_dims_len,clo_dims_len,cla_dims_len,&
23                                           imask_dims_len 
24  !               
25  INTEGER, INTENT(in)     :: nlon,nlat,corners_ij_lus
26  !
27  CHARACTER(len=30)        :: data_filename
28  !
29  INTEGER,  DIMENSION(3)   :: ila_dim
30  INTEGER,  DIMENSION(3)   :: ila_corners,ila_what
31  !
32  DOUBLE PRECISION, DIMENSION(nlon,nlat)                :: gridlon,gridlat,gridsrf
33  DOUBLE PRECISION, DIMENSION(nlon,nlat,corners_ij_lus) :: gridclo,gridcla
34  INTEGER, DIMENSION(nlon,nlat)                      :: indice_mask
35  !
36  !
37  ! Dimensions
38  !
39  CALL hdlerr(NF90_OPEN(data_filename, NF90_NOWRITE, il_file_id), __LINE__ )
40  !
41  !
42  CALL hdlerr( NF90_INQ_VARID(il_file_id, 'lon' , il_lon_id), __LINE__ )
43  CALL hdlerr( NF90_INQ_VARID(il_file_id, 'lat' , il_lat_id), __LINE__ )
44  CALL hdlerr( NF90_INQ_VARID(il_file_id, 'clo' , il_clo_id), __LINE__ )
45  CALL hdlerr( NF90_INQ_VARID(il_file_id, 'cla' , il_cla_id), __LINE__ )
46  CALL hdlerr( NF90_INQ_VARID(il_file_id, 'srf' , il_srf_id), __LINE__ )
47  CALL hdlerr( NF90_INQ_VARID(il_file_id, 'imask' , il_indice_id), __LINE__ )
48  !
49  ila_what(:)=1
50  !
51  ila_dim(1)=nlon
52  ila_dim(2)=nlat
53  ila_dim(3)=1
54  !
55  ila_corners(1)=nlon
56  ila_corners(2)=nlat
57  ila_corners(3)=corners_ij_lus
58  !
59  ! Data
60  !
61  CALL hdlerr( NF90_OPEN(data_filename, NF90_NOWRITE, il_file_id), __LINE__ )
62  !
63  CALL hdlerr( NF90_GET_VAR (il_file_id, il_lon_id, gridlon, &
64     ila_what(1:2), ila_dim(1:2)), __LINE__ )
65  WRITE(w_unit,*) 'Global grid longitudes reading done'
66  CALL flush(w_unit)
67  !
68  CALL hdlerr( NF90_GET_VAR (il_file_id, il_lat_id, gridlat, &
69     ila_what(1:2), ila_dim(1:2)), __LINE__ )
70  WRITE(w_unit,*) 'Global grid latitudes reading done'
71  CALL flush(w_unit)
72  !
73  CALL hdlerr( NF90_GET_VAR(il_file_id, il_clo_id, gridclo, &
74     ila_what, ila_corners), __LINE__ )
75  WRITE(w_unit,*) 'Global grid longitude corners reading done'
76  CALL flush(w_unit)
77  !
78  CALL hdlerr( NF90_GET_VAR (il_file_id, il_cla_id, gridcla, &
79     ila_what, ila_corners), __LINE__ )
80  WRITE(w_unit,*) 'Global grid latitude corners reading done'
81  CALL flush(w_unit)
82  !
83  CALL hdlerr( NF90_GET_VAR (il_file_id, il_srf_id, gridsrf, &
84     ila_what(1:2), ila_dim(1:2)), __LINE__ )
85  WRITE(w_unit,*) 'Global grid surfaces reading done'
86  CALL flush(w_unit)
87  !
88  CALL hdlerr( NF90_GET_VAR (il_file_id, il_indice_id, indice_mask, &
89     ila_what, ila_dim), __LINE__ )
90  WRITE(w_unit,*) 'Global grid mask reading done'
91  CALL flush(w_unit)
92  !
93  CALL hdlerr( NF90_CLOSE(il_file_id), __LINE__ )
94  !
95  WRITE(w_unit,*) 'End of routine read_grid'
96  CALL flush(w_unit)
97  !
98END SUBROUTINE read_grid
Note: See TracBrowser for help on using the repository browser.