source: CPL/oasis3-mct_5.0/examples/spoc/spoc_communication/read_grid.F90 @ 6328

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

First import of oasis3-mct_5.0 (from oasis git server, branch OASIS3-MCT_5.0)

File size: 3.7 KB
Line 
1  !****************************************************************************************
2  SUBROUTINE read_grid (nlon, nlat, nc, id_begi, id_begj, id_lon, id_lat, &
3                              data_filename, w_unit,                       &
4                              dda_lon, dda_lat, dda_clo, dda_cla, dda_srf, ida_mask)
5  !**************************************************************************************
6  !
7  USE netcdf
8  IMPLICIT NONE
9  !
10  INTEGER, INTENT(in)             :: nlon, nlat, nc, id_begi, id_begj, id_lon, id_lat
11  CHARACTER(len=*), INTENT(in)   :: data_filename
12  INTEGER, INTENT(in)             :: w_unit
13  DOUBLE PRECISION, DIMENSION(id_lon, id_lat), INTENT(out)       :: dda_lon, dda_lat, dda_srf
14  DOUBLE PRECISION, DIMENSION(id_lon, id_lat, nc), INTENT(out)   :: dda_clo, dda_cla
15  INTEGER, DIMENSION(id_lon, id_lat), INTENT(out)                :: ida_mask
16  !
17  INTEGER :: il_file_id, il_lon_id, il_lat_id, il_clo_id, il_cla_id, il_srf_id, il_msk_id
18  !               
19  INTEGER,  DIMENSION(3)          :: ila_dim, ila_st
20  !
21#ifdef _DEBUG
22  WRITE(w_unit,*) 'Starting read_grid'
23  CALL flush(w_unit)
24#endif
25  CALL hdlerr (NF90_OPEN(data_filename, NF90_NOWRITE, il_file_id), __LINE__ )
26  !
27  !**************************************************************************************
28  !
29  CALL hdlerr( NF90_INQ_VARID(il_file_id, 'lon' , il_lon_id), __LINE__ )
30  CALL hdlerr( NF90_INQ_VARID(il_file_id, 'lat' , il_lat_id), __LINE__ )
31  CALL hdlerr( NF90_INQ_VARID(il_file_id, 'clo' , il_clo_id), __LINE__ )
32  CALL hdlerr( NF90_INQ_VARID(il_file_id, 'cla' , il_cla_id), __LINE__ )
33  CALL hdlerr( NF90_INQ_VARID(il_file_id, 'srf' , il_srf_id), __LINE__ )
34  CALL hdlerr( NF90_INQ_VARID(il_file_id, 'imask' , il_msk_id), __LINE__ )
35  !
36  CALL flush(w_unit)
37  ila_st(1) = id_begi
38  ila_st(2) = id_begj
39  ila_st(3) = 1
40  !
41  ila_dim(1) = id_lon
42  ila_dim(2) = id_lat
43  ila_dim(3) = nc
44  !
45  CALL hdlerr( NF90_GET_VAR (il_file_id, il_lon_id, dda_lon, ila_st(1:2), ila_dim(1:2)), __LINE__ )
46  CALL hdlerr( NF90_GET_VAR (il_file_id, il_lat_id, dda_lat, ila_st(1:2), ila_dim(1:2)), __LINE__ )
47#ifdef _DEBUG
48  WRITE(w_unit,*) 'Local grid longitudes and latitudes read from file'
49  CALL flush(w_unit)
50#endif
51  !
52  CALL hdlerr( NF90_GET_VAR (il_file_id, il_clo_id, dda_clo, ila_st(1:3), ila_dim(1:3)), __LINE__ )
53  CALL hdlerr( NF90_GET_VAR (il_file_id, il_cla_id, dda_cla, ila_st(1:3), ila_dim(1:3)), __LINE__ )
54#ifdef _DEBUG
55  WRITE(w_unit,*) 'Local grid corner longitudes and latitudes read from file'
56  CALL flush(w_unit)
57#endif
58  !
59  CALL hdlerr( NF90_GET_VAR (il_file_id, il_srf_id, dda_srf, ila_st(1:2), ila_dim(1:2)), __LINE__ )
60  CALL hdlerr( NF90_GET_VAR (il_file_id, il_msk_id, ida_mask, ila_st(1:2), ila_dim(1:2)), __LINE__ )
61  !
62#ifdef _DEBUG
63  WRITE(w_unit,*) 'Local grid areas and mask read from file'
64  CALL flush(w_unit)
65#endif
66  !
67  CALL hdlerr( NF90_CLOSE(il_file_id), __LINE__ )
68  !
69  ! OASIS3 mask convention (1=masked, 0=not masked) is opposite to usual one)
70  !
71  WHERE (ida_mask == 0) ; ida_mask = 1; ELSEWHERE ; ida_mask = 0; END WHERE
72  !
73#ifdef _DEBUG
74  WRITE(w_unit,*) 'End of routine read_grid'
75  CALL flush(w_unit)
76#endif
77  !
78END SUBROUTINE read_grid
79!
80!*********************************************************************************
81SUBROUTINE hdlerr(istatus, line)
82  !*********************************************************************************
83  use netcdf
84  implicit none
85  !
86  INCLUDE 'mpif.h'
87  !
88  ! Check for error message from NetCDF call
89  !
90  integer, intent(in) :: istatus, line
91  integer             :: ierror
92  !
93  IF (istatus .NE. NF90_NOERR) THEN
94      write ( * , * ) 'NetCDF problem at line',line
95      write ( * , * ) 'Stopped '
96      call MPI_Abort ( MPI_COMM_WORLD, 1, ierror )
97  ENDIF
98  !
99  RETURN
100END SUBROUTINE hdlerr
Note: See TracBrowser for help on using the repository browser.