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

Last change on this file since 8780 was 8780, checked in by steve, 11 years ago

Some changes to allow netcdf4 use

File size: 7.7 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 Scientific.IO.NetCDF import NetCDFFile
45    #return NetCDFFile(file_name, netcdf_mode)
46
47    from netCDF4 import Dataset
48    return Dataset(file_name, netcdf_mode, format='NETCDF3_64BIT')
49
50    #return Dataset(file_name, netcdf_mode, format='NETCDF3_CLASSIC')
51
52
53    # COMMENT SR; Can't use scipy.io.netcdf as we can't append to
54    # a file as of 2013/03/26
55    #from scipy.io.netcdf import netcdf_file
56    #return netcdf_file(file_name, netcdf_mode, version=2)
57
58
59
60
61class Write_nc:
62    """Write an nc file.
63
64    Note, this should be checked to meet cdc netcdf conventions for gridded
65    data. http://www.cdc.noaa.gov/cdc/conventions/cdc_netcdf_standard.shtml
66    """
67
68    def __init__(self,
69                 quantity_name,
70                 file_name,
71                 time_step_count,
72                 time_step,
73                 lon,
74                 lat):
75        """Instantiate a Write_nc instance (NetCDF file writer).
76
77        time_step_count is the number of time steps.
78        time_step is the time step size
79
80        pre-condition: quantity_name must be 'HA', 'UA'or 'VA'.
81        """
82
83        self.quantity_name = quantity_name
84        quantity_units = {'HA':'CENTIMETERS',
85                          'UA':'CENTIMETERS/SECOND',
86                          'VA':'CENTIMETERS/SECOND'}
87
88        multiplier_dic = {'HA':100.0,   # To convert from m to cm
89                          'UA':100.0,   #             and m/s to cm/sec
90                          'VA':-100.0}  # MUX files have positive x in the
91                                        # Southern direction.  This corrects
92                                        # for it, when writing nc files.
93
94        self.quantity_multiplier =  multiplier_dic[self.quantity_name]
95
96        #self.file_name = file_name
97        self.time_step_count = time_step_count
98        self.time_step = time_step
99
100        # NetCDF file definition
101        self.outfile = NetCDFFile(file_name, netcdf_mode_w)
102        outfile = self.outfile
103
104        #Create new file
105        nc_lon_lat_header(outfile, lon, lat)
106
107        # TIME
108        outfile.createDimension(time_name, None)
109        outfile.createVariable(time_name, precision, (time_name,))
110
111        #QUANTITY
112        outfile.createVariable(self.quantity_name, precision,
113                               (time_name, lat_name, lon_name))
114        outfile.variables[self.quantity_name].missing_value = -1.e+034
115        outfile.variables[self.quantity_name].units = \
116                                 quantity_units[self.quantity_name]
117        outfile.variables[lon_name][:]= ensure_numeric(lon)
118        outfile.variables[lat_name][:]= ensure_numeric(lat)
119
120        #Assume no one will be wanting to read this, while we are writing
121        #outfile.close()
122
123    def store_timestep(self, quantity_slice):
124        """Write a time slice of quantity info
125
126        quantity_slice is the data to be stored at this time step
127        """
128
129        # Get the variables
130        time = self.outfile.variables[time_name]
131        quantity = self.outfile.variables[self.quantity_name]
132
133        # get index oflice to write
134        i = len(time)
135
136        #Store time
137        time[i] = i * self.time_step    #self.domain.time
138        quantity[i,:] = quantity_slice * self.quantity_multiplier
139
140    def close(self):
141        """ Close file underlying the class instance. """
142        self.outfile.close()
143
144
145
146def write_elevation_nc(file_out, lon, lat, depth_vector):
147    """Write an nc elevation file."""
148
149    # NetCDF file definition
150    outfile = NetCDFFile(file_out, netcdf_mode_w)
151
152    #Create new file
153    nc_lon_lat_header(outfile, lon, lat)
154
155    # ELEVATION
156    zname = 'ELEVATION'
157    outfile.createVariable(zname, precision, (lat_name, lon_name))
158    outfile.variables[zname].units = 'CENTIMETERS'
159    outfile.variables[zname].missing_value = -1.e+034
160
161    outfile.variables[lon_name][:] = ensure_numeric(lon)
162    outfile.variables[lat_name][:] = ensure_numeric(lat)
163
164    depth = num.reshape(depth_vector, (len(lat), len(lon)))
165    outfile.variables[zname][:] = depth
166
167    outfile.close()
168
169
170def nc_lon_lat_header(outfile, lon, lat):
171    """Write lat/lon headers to a NetCDF file.
172
173    outfile is the netcdf file handle.
174    lon - a list/array of the longitudes
175    lat - a list/array of the latitudes
176    """
177
178    outfile.institution = 'Geoscience Australia'
179    outfile.description = 'Converted from URS binary C'
180
181    # Longitude
182    outfile.createDimension(lon_name, len(lon))
183    outfile.createVariable(lon_name, precision, (lon_name,))
184    outfile.variables[lon_name].point_spacing = 'uneven'
185    outfile.variables[lon_name].units = 'degrees_east'
186    outfile.variables[lon_name].assignValue(lon)
187
188    # Latitude
189    outfile.createDimension(lat_name, len(lat))
190    outfile.createVariable(lat_name, precision, (lat_name,))
191    outfile.variables[lat_name].point_spacing = 'uneven'
192    outfile.variables[lat_name].units = 'degrees_north'
193    outfile.variables[lat_name].assignValue(lat)
194
195
196
197
198def filter_netcdf(filename1, filename2, first=0, last=None, step=1):
199    """Filter data file, selecting timesteps first:step:last.
200   
201    Read netcdf filename1, pick timesteps first:step:last and save to
202    nettcdf file filename2
203    """
204
205    from Scientific.IO.NetCDF import NetCDFFile
206
207    # Get NetCDF
208    infile = NetCDFFile(filename1, netcdf_mode_r)  #Open existing file for read
209    outfile = NetCDFFile(filename2, netcdf_mode_w)  #Open new file
210
211    # Copy dimensions
212    for d in infile.dimensions:
213        outfile.createDimension(d, infile.dimensions[d])
214
215    # Copy variable definitions
216    for name in infile.variables:
217        var = infile.variables[name]
218        outfile.createVariable(name, var.dtype.char, var.dimensions)
219
220    # Copy the static variables
221    for name in infile.variables:
222        if name == 'time' or name == 'stage':
223            pass
224        else:
225            outfile.variables[name][:] = infile.variables[name][:]
226
227    # Copy selected timesteps
228    time = infile.variables['time']
229    stage = infile.variables['stage']
230
231    newtime = outfile.variables['time']
232    newstage = outfile.variables['stage']
233
234    if last is None:
235        last = len(time)
236
237    selection = range(first, last, step)
238    for i, j in enumerate(selection):
239        log.critical('Copying timestep %d of %d (%f)'
240                     % (j, last-first, time[j]))
241        newtime[i] = time[j]
242        newstage[i,:] = stage[j,:]
243
244    # Close
245    infile.close()
246    outfile.close()
247
Note: See TracBrowser for help on using the repository browser.