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

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

Some changes to netcdf config

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