New URL for NEMO forge!   http://forge.nemo-ocean.eu

Since March 2022 along with NEMO 4.2 release, the code development moved to a self-hosted GitLab.
This present forge is now archived and remained online for history.
readcoordmesh.f90 in branches/2014/dev_r4650_UKMO14.5_SST_BIAS_CORRECTION/NEMOGCM/TOOLS/SECTIONS_DIADCT/src – NEMO

source: branches/2014/dev_r4650_UKMO14.5_SST_BIAS_CORRECTION/NEMOGCM/TOOLS/SECTIONS_DIADCT/src/readcoordmesh.f90 @ 5967

Last change on this file since 5967 was 5967, checked in by timgraham, 8 years ago

Reset keywords before merging with head of trunk

  • Property svn:executable set to *
  • Property svn:keywords set to Id
File size: 5.9 KB
Line 
1MODULE readcoordmesh
2  !!=====================================================================
3  !!                       ***  MODULE  readcoordmesh  ***
4  !!
5  !! History: 2011: Clement Bricaud, Mercator-Ocean
6  !!
7  !!=====================================================================
8  !! * Modules used
9  USE netcdf
10  USE declarations
11
12  IMPLICIT NONE
13  PRIVATE
14 
15  PUBLIC  read_coord_mesh
16
17CONTAINS
18
19  SUBROUTINE read_coord_mesh
20  !!---------------------------------------------------------------------
21  !!              ***  ROUTINE coord_mesh_read  ***
22  !!
23  !! ** Purpose :   Read a coordinate file and a meshmask file in NetCDF format
24  !!
25  !!---------------------------------------------------------------------     
26  PRINT*,'              '
27  PRINT*,'READ COORDINATES AND MESHMASK'
28  PRINT*,'-----------------------------'
29
30  ! Get coordinates dimensions
31  CALL getdim(cdfile="coordinates.nc")
32
33  !Allocate coordinates array with domain size
34  ALLOCATE(glamt(jpi,jpj)) ; ALLOCATE(gphit(jpi,jpj))
35  ALLOCATE(glamf(jpi,jpj)) ; ALLOCATE(gphif(jpi,jpj))
36  ALLOCATE(e1t(jpi,jpj)  )   
37
38  !Read glamt
39  CALL read_ncdf(cdfile="coordinates.nc",cdvar="glamt",ksize=(/jpi,jpj,1,1/),ptab=glamt)
40
41  !Read gphit
42  CALL read_ncdf(cdfile="coordinates.nc",cdvar="gphit",ksize=(/jpi,jpj,1,1/),ptab=gphit)
43
44  !Read glamf
45  CALL read_ncdf(cdfile="coordinates.nc",cdvar="glamf",ksize=(/jpi,jpj,1,1/),ptab=glamf)
46
47  !Read gphif
48  CALL read_ncdf(cdfile="coordinates.nc",cdvar="gphif",ksize=(/jpi,jpj,1,1/),ptab=gphif)
49
50  !Read e1t
51  CALL read_ncdf(cdfile="coordinates.nc",cdvar="e1t",ksize=(/jpi,jpj,1,1/),ptab=e1t)
52
53  END SUBROUTINE read_coord_mesh
54
55  SUBROUTINE getdim(cdfile)
56  !!----------------------------------------------------------------------
57  !!              ***  ROUTINE getdim  ***
58  !!
59  !! ** Purpose :   get dimensions of a netcdf file
60  !!
61  !!----------------------------------------------------------------------
62  !! * Arguments
63  CHARACTER(*),INTENT(IN):: cdfile
64
65  !! * Local declarations
66  INTEGER           :: ncid                 ! file unit
67  INTEGER           :: idims                ! number of dimensions
68  INTEGER           :: istatus, id_var      ! dummy variable
69  CHARACTER(len=30) :: clname               ! dimension name   
70  INTEGER, ALLOCATABLE,DIMENSION(:) :: ndim ! dimension value
71  !!----------------------------------------------------------------------
72
73  !Open file
74  istatus=NF90_OPEN(TRIM(cdfile),nf90_nowrite,ncid)
75
76  IF( istatus/= NF90_NOERR )THEN
77     PRINT*,TRIM(cdfile),' not found.stop ' ; STOP
78  ELSE
79 
80     ! read number of dimensions   
81     istatus=NF90_INQUIRE(ncid,ndimensions=idims)
82
83     ALLOCATE( ndim(idims) )
84
85     ! read each dimension
86     PRINT*,'     File dimensions: '
87     DO id_var=1,idims
88        istatus=NF90_Inquire_Dimension(ncid,id_var,clname,ndim(id_var))
89        PRINT*,'       ',id_var,clname,ndim(id_var)
90     ENDDO
91
92     !close
93     istatus=NF90_CLOSE( ncid )             
94     IF( istatus/=NF90_NOERR )THEN
95        PRINT*,'Problem for closing ',TRIM(cdfile);STOP
96     ELSE
97        PRINT*,'     close ',TRIM(cdfile),' OK'
98     ENDIF
99
100  ENDIF
101
102  !domain dimensions
103  jpi = ndim(1)
104  jpj = ndim(2)
105
106  DEALLOCATE( ndim )
107  END SUBROUTINE getdim
108
109  SUBROUTINE read_ncdf(cdfile,cdvar,ksize,ptab,kstart,kcount)
110  !!----------------------------------------------------------------------
111  !!              ***  ROUTINE coord_mesh_read  ***
112  !!
113  !! ** Purpose :   Read a coordinate and a meshmask file in NetCDF format
114  !!
115  !!----------------------------------------------------------------------
116  !! * Arguments
117  CHARACTER(*),        INTENT(IN)                                    :: cdfile
118  CHARACTER(*),        INTENT(IN)                                    :: cdvar
119  INTEGER,DIMENSION(4),INTENT(IN)                                    :: ksize
120  INTEGER,DIMENSION(4),INTENT(IN),OPTIONAL                           :: kstart,kcount
121  REAL(wp),DIMENSION(ksize(1),ksize(2),ksize(3),ksize(4)),INTENT(OUT):: ptab
122
123  !! * Local declarations
124  INTEGER                                 ::istatus,ncid,id_var,len
125  CHARACTER(len=30) :: clname , cdvar2
126
127  INTEGER :: idims
128  INTEGER,DIMENSION(3)::idimids 
129  !!----------------------------------------------------------------------
130  ptab=0.
131  PRINT*,'read ',TRIM(cdvar),' in ',TRIM(cdfile)
132 
133  !OPEN
134  !----
135  istatus=NF90_OPEN(TRIM(cdfile),nf90_nowrite,ncid)
136  IF( istatus/= NF90_NOERR )THEN
137     PRINT*,TRIM(cdfile),' not found.stop ' ; STOP
138  ENDIF
139
140  !READ
141  !----
142  !search variable
143  istatus=NF90_INQ_VARID (ncid,TRIM(cdvar),id_var)
144  IF( istatus/=NF90_NOERR )THEN
145      PRINT*,TRIM(cdvar),' not found in ',TRIM(cdfile),'.stop';STOP
146  ENDIF
147
148  !get variable
149  !------------
150  istatus=nf90_inquire_variable(ncid,id_var, cdvar2, ndims=idims, dimids=idimids)
151  IF ( PRESENT(kstart) .AND. PRESENT(kcount) )THEN 
152      istatus=NF90_GET_VAR(ncid,id_var,ptab,start=kstart,count=kcount)
153  ELSE
154      istatus=NF90_GET_VAR(ncid,id_var,ptab)
155  ENDIF
156
157  CALL ERR_HDL(istatus)
158
159  IF( istatus/=NF90_NOERR )THEN
160           PRINT*,'Problem for reading ',TRIM(cdvar),' in ',TRIM(cdfile); STOP
161  ELSE
162      PRINT*,'     read ',TRIM(cdvar),' OK'
163  ENDIF
164
165  !CLOSE
166  !-----
167  istatus=NF90_CLOSE( ncid )
168  IF( istatus/=NF90_NOERR )THEN
169      PRINT*,'Problem for closing ',TRIM(cdfile);stop
170  ELSE
171      PRINT*,'     close ',TRIM(cdfile),' OK'
172  ENDIF
173
174
175  END SUBROUTINE read_ncdf
176
177  SUBROUTINE ERR_HDL(kstatus)
178  !! ----------------------------------------------------------
179  !!   ***  SUBROUTINE err_hdl
180  !!
181  !!   ** Purpose :  Error handle for NetCDF routine.
182  !!          Stop if kstatus indicates error conditions.
183  !!
184  !! History :
185  !!     Original: J.M. Molines (01/99)
186  !!
187  !! -----------------------------------------------------------
188  INTEGER, INTENT(in) ::  kstatus
189
190  !! -----------------------------------------------------------
191  IF( kstatus /=  NF90_NOERR ) THEN
192     PRINT *, 'ERROR in NETCDF routine, status=',kstatus
193     PRINT *,NF90_STRERROR(kstatus)
194     STOP
195  END IF
196
197  END SUBROUTINE ERR_HDL
198
199END MODULE readcoordmesh
Note: See TracBrowser for help on using the repository browser.