source: trunk/anuga_core/source/anuga/file/netcdf.py @ 8813

Last change on this file since 8813 was 8812, checked in by steve, 12 years ago

Fixed unit test in test_urs2sww.py to work with netcdf4 on windows

File size: 8.2 KB
Line 
1""" This module is responsible for loading and saving NetCDF NC files
2"""
3
4import numpy as num
5
6from anuga.coordinate_transforms.redfearn import \
7     convert_from_latlon_to_utm
8from anuga.config import minimum_storable_height as \
9     default_minimum_storable_height
10from anuga.config import netcdf_mode_r, netcdf_mode_w
11from anuga.config import netcdf_float, netcdf_int
12from anuga.utilities.numerical_tools import ensure_numeric,  mean
13
14import anuga.utilities.log as log
15
16
17# Definitions of various NetCDF dimension names, etc.
18lon_name = 'LON'
19lat_name = 'LAT'
20time_name = 'TIME'
21precision = netcdf_float # So if we want to change the precision its done here
22
23
24
25def NetCDFFile(file_name, netcdf_mode=netcdf_mode_r):
26    """Wrapper to isolate changes of the netcdf libray.
27
28    In theory we should be able to change over to NetCDF4 via this
29    wrapper, by ensuring the interface to the NetCDF library isthe same as the
30    the old Scientific.IO.NetCDF library.
31
32    There is a difference between extracting dimensions. We have used the following
33    to cover netcdf4 and scientific python
34
35    try: # works with netcdf4
36        number_of_timesteps = len(fid.dimensions['number_of_timesteps'])
37        number_of_points = len(fid.dimensions['number_of_points'])
38    except: # works with Scientific.IO.NetCDF
39        number_of_timesteps = fid.dimensions['number_of_timesteps']
40        number_of_points = fid.dimensions['number_of_points']
41   
42    """
43
44    from netCDF4 import Dataset
45    if netcdf_mode == 'wl' :
46        return Dataset(file_name, 'w', format='NETCDF3_64BIT')
47    else:
48        return Dataset(file_name, netcdf_mode, format='NETCDF3_64BIT')
49
50
51#    try:
52#        from Scientific.IO.NetCDF import NetCDFFile
53#        return NetCDFFile(file_name, netcdf_mode)
54#    except:
55#        from netCDF4 import Dataset
56#        if netcdf_mode == 'wl' :
57#            return Dataset(file_name, 'w', format='NETCDF3_64BIT')
58#        else:
59#            return Dataset(file_name, netcdf_mode, format='NETCDF3_64BIT')
60
61#    from netCDF4 import Dataset
62#    return Dataset(file_name, netcdf_mode, format='NETCDF3_64BIT')
63
64    #return Dataset(file_name, netcdf_mode, format='NETCDF3_CLASSIC')
65
66
67    # COMMENT SR; Can't use scipy.io.netcdf as we can't append to
68    # a file as of 2013/03/26
69    #from scipy.io.netcdf import netcdf_file
70    #return netcdf_file(file_name, netcdf_mode, version=2)
71
72
73
74
75class Write_nc:
76    """Write an nc file.
77
78    Note, this should be checked to meet cdc netcdf conventions for gridded
79    data. http://www.cdc.noaa.gov/cdc/conventions/cdc_netcdf_standard.shtml
80    """
81
82    def __init__(self,
83                 quantity_name,
84                 file_name,
85                 time_step_count,
86                 time_step,
87                 lon,
88                 lat):
89        """Instantiate a Write_nc instance (NetCDF file writer).
90
91        time_step_count is the number of time steps.
92        time_step is the time step size
93
94        pre-condition: quantity_name must be 'HA', 'UA'or 'VA'.
95        """
96
97        self.quantity_name = quantity_name
98        quantity_units = {'HA':'CENTIMETERS',
99                          'UA':'CENTIMETERS/SECOND',
100                          'VA':'CENTIMETERS/SECOND'}
101
102        multiplier_dic = {'HA':100.0,   # To convert from m to cm
103                          'UA':100.0,   #             and m/s to cm/sec
104                          'VA':-100.0}  # MUX files have positive x in the
105                                        # Southern direction.  This corrects
106                                        # for it, when writing nc files.
107
108        self.quantity_multiplier =  multiplier_dic[self.quantity_name]
109
110        #self.file_name = file_name
111        self.time_step_count = time_step_count
112        self.time_step = time_step
113
114        # NetCDF file definition
115        self.outfile = NetCDFFile(file_name, netcdf_mode_w)
116        outfile = self.outfile
117
118        #Create new file
119        nc_lon_lat_header(outfile, lon, lat)
120
121        # TIME
122        outfile.createDimension(time_name, None)
123        outfile.createVariable(time_name, precision, (time_name,))
124
125        #QUANTITY
126        outfile.createVariable(self.quantity_name, precision,
127                               (time_name, lat_name, lon_name))
128        outfile.variables[self.quantity_name].missing_value = -1.e+034
129        outfile.variables[self.quantity_name].units = \
130                                 quantity_units[self.quantity_name]
131        outfile.variables[lon_name][:]= ensure_numeric(lon)
132        outfile.variables[lat_name][:]= ensure_numeric(lat)
133
134        #Assume no one will be wanting to read this, while we are writing
135        #outfile.close()
136
137    def store_timestep(self, quantity_slice):
138        """Write a time slice of quantity info
139
140        quantity_slice is the data to be stored at this time step
141        """
142
143        # Get the variables
144        time = self.outfile.variables[time_name]
145        quantity = self.outfile.variables[self.quantity_name]
146
147        # get index oflice to write
148        i = len(time)
149
150        #Store time
151        time[i] = i * self.time_step    #self.domain.time
152        quantity[i,:] = quantity_slice * self.quantity_multiplier
153
154    def close(self):
155        """ Close file underlying the class instance. """
156        self.outfile.close()
157
158
159
160def write_elevation_nc(file_out, lon, lat, depth_vector):
161    """Write an nc elevation file."""
162
163    # NetCDF file definition
164    outfile = NetCDFFile(file_out, netcdf_mode_w)
165
166    #Create new file
167    nc_lon_lat_header(outfile, lon, lat)
168
169    # ELEVATION
170    zname = 'ELEVATION'
171    outfile.createVariable(zname, precision, (lat_name, lon_name))
172    outfile.variables[zname].units = 'CENTIMETERS'
173    outfile.variables[zname].missing_value = -1.e+034
174
175    outfile.variables[lon_name][:] = ensure_numeric(lon)
176    outfile.variables[lat_name][:] = ensure_numeric(lat)
177
178    depth = num.reshape(depth_vector, (len(lat), len(lon)))
179    outfile.variables[zname][:] = depth
180
181    outfile.close()
182
183
184def nc_lon_lat_header(outfile, lon, lat):
185    """Write lat/lon headers to a NetCDF file.
186
187    outfile is the netcdf file handle.
188    lon - a list/array of the longitudes
189    lat - a list/array of the latitudes
190    """
191
192    outfile.institution = 'Geoscience Australia'
193    outfile.description = 'Converted from URS binary C'
194
195    # Longitude
196    outfile.createDimension(lon_name, len(lon))
197    outfile.createVariable(lon_name, precision, (lon_name,))
198    outfile.variables[lon_name].point_spacing = 'uneven'
199    outfile.variables[lon_name].units = 'degrees_east'
200    outfile.variables[lon_name].assignValue(lon)
201
202    # Latitude
203    outfile.createDimension(lat_name, len(lat))
204    outfile.createVariable(lat_name, precision, (lat_name,))
205    outfile.variables[lat_name].point_spacing = 'uneven'
206    outfile.variables[lat_name].units = 'degrees_north'
207    outfile.variables[lat_name].assignValue(lat)
208
209
210
211
212def filter_netcdf(filename1, filename2, first=0, last=None, step=1):
213    """Filter data file, selecting timesteps first:step:last.
214   
215    Read netcdf filename1, pick timesteps first:step:last and save to
216    nettcdf file filename2
217    """
218
219    from Scientific.IO.NetCDF import NetCDFFile
220
221    # Get NetCDF
222    infile = NetCDFFile(filename1, netcdf_mode_r)  #Open existing file for read
223    outfile = NetCDFFile(filename2, netcdf_mode_w)  #Open new file
224
225    # Copy dimensions
226    for d in infile.dimensions:
227        outfile.createDimension(d, infile.dimensions[d])
228
229    # Copy variable definitions
230    for name in infile.variables:
231        var = infile.variables[name]
232        outfile.createVariable(name, var.dtype.char, var.dimensions)
233
234    # Copy the static variables
235    for name in infile.variables:
236        if name == 'time' or name == 'stage':
237            pass
238        else:
239            outfile.variables[name][:] = infile.variables[name][:]
240
241    # Copy selected timesteps
242    time = infile.variables['time']
243    stage = infile.variables['stage']
244
245    newtime = outfile.variables['time']
246    newstage = outfile.variables['stage']
247
248    if last is None:
249        last = len(time)
250
251    selection = range(first, last, step)
252    for i, j in enumerate(selection):
253        log.critical('Copying timestep %d of %d (%f)'
254                     % (j, last-first, time[j]))
255        newtime[i] = time[j]
256        newstage[i,:] = stage[j,:]
257
258    # Close
259    infile.close()
260    outfile.close()
261
Note: See TracBrowser for help on using the repository browser.