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

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

Changing to netcdf.py to run scientific python netcdf in preference to netcdf4

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