1 | """ This module is responsible for loading and saving NetCDF NC files |
---|
2 | """ |
---|
3 | |
---|
4 | import numpy as num |
---|
5 | |
---|
6 | from Scientific.IO.NetCDF import NetCDFFile |
---|
7 | |
---|
8 | from anuga.coordinate_transforms.redfearn import \ |
---|
9 | convert_from_latlon_to_utm |
---|
10 | from anuga.config import minimum_storable_height as \ |
---|
11 | default_minimum_storable_height |
---|
12 | from anuga.config import netcdf_mode_r, netcdf_mode_w |
---|
13 | from anuga.config import netcdf_float, netcdf_int |
---|
14 | from anuga.utilities.numerical_tools import ensure_numeric, mean |
---|
15 | |
---|
16 | import anuga.utilities.log as log |
---|
17 | |
---|
18 | |
---|
19 | # Definitions of various NetCDF dimension names, etc. |
---|
20 | lon_name = 'LON' |
---|
21 | lat_name = 'LAT' |
---|
22 | time_name = 'TIME' |
---|
23 | precision = netcdf_float # So if we want to change the precision its done here |
---|
24 | |
---|
25 | class Write_nc: |
---|
26 | """Write an nc file. |
---|
27 | |
---|
28 | Note, this should be checked to meet cdc netcdf conventions for gridded |
---|
29 | data. http://www.cdc.noaa.gov/cdc/conventions/cdc_netcdf_standard.shtml |
---|
30 | """ |
---|
31 | |
---|
32 | def __init__(self, |
---|
33 | quantity_name, |
---|
34 | file_name, |
---|
35 | time_step_count, |
---|
36 | time_step, |
---|
37 | lon, |
---|
38 | lat): |
---|
39 | """Instantiate a Write_nc instance (NetCDF file writer). |
---|
40 | |
---|
41 | time_step_count is the number of time steps. |
---|
42 | time_step is the time step size |
---|
43 | |
---|
44 | pre-condition: quantity_name must be 'HA', 'UA'or 'VA'. |
---|
45 | """ |
---|
46 | |
---|
47 | self.quantity_name = quantity_name |
---|
48 | quantity_units = {'HA':'CENTIMETERS', |
---|
49 | 'UA':'CENTIMETERS/SECOND', |
---|
50 | 'VA':'CENTIMETERS/SECOND'} |
---|
51 | |
---|
52 | multiplier_dic = {'HA':100.0, # To convert from m to cm |
---|
53 | 'UA':100.0, # and m/s to cm/sec |
---|
54 | 'VA':-100.0} # MUX files have positive x in the |
---|
55 | # Southern direction. This corrects |
---|
56 | # for it, when writing nc files. |
---|
57 | |
---|
58 | self.quantity_multiplier = multiplier_dic[self.quantity_name] |
---|
59 | |
---|
60 | #self.file_name = file_name |
---|
61 | self.time_step_count = time_step_count |
---|
62 | self.time_step = time_step |
---|
63 | |
---|
64 | # NetCDF file definition |
---|
65 | self.outfile = NetCDFFile(file_name, netcdf_mode_w) |
---|
66 | outfile = self.outfile |
---|
67 | |
---|
68 | #Create new file |
---|
69 | nc_lon_lat_header(outfile, lon, lat) |
---|
70 | |
---|
71 | # TIME |
---|
72 | outfile.createDimension(time_name, None) |
---|
73 | outfile.createVariable(time_name, precision, (time_name,)) |
---|
74 | |
---|
75 | #QUANTITY |
---|
76 | outfile.createVariable(self.quantity_name, precision, |
---|
77 | (time_name, lat_name, lon_name)) |
---|
78 | outfile.variables[self.quantity_name].missing_value = -1.e+034 |
---|
79 | outfile.variables[self.quantity_name].units = \ |
---|
80 | quantity_units[self.quantity_name] |
---|
81 | outfile.variables[lon_name][:]= ensure_numeric(lon) |
---|
82 | outfile.variables[lat_name][:]= ensure_numeric(lat) |
---|
83 | |
---|
84 | #Assume no one will be wanting to read this, while we are writing |
---|
85 | #outfile.close() |
---|
86 | |
---|
87 | def store_timestep(self, quantity_slice): |
---|
88 | """Write a time slice of quantity info |
---|
89 | |
---|
90 | quantity_slice is the data to be stored at this time step |
---|
91 | """ |
---|
92 | |
---|
93 | # Get the variables |
---|
94 | time = self.outfile.variables[time_name] |
---|
95 | quantity = self.outfile.variables[self.quantity_name] |
---|
96 | |
---|
97 | # get index oflice to write |
---|
98 | i = len(time) |
---|
99 | |
---|
100 | #Store time |
---|
101 | time[i] = i * self.time_step #self.domain.time |
---|
102 | quantity[i,:] = quantity_slice * self.quantity_multiplier |
---|
103 | |
---|
104 | def close(self): |
---|
105 | """ Close file underlying the class instance. """ |
---|
106 | self.outfile.close() |
---|
107 | |
---|
108 | |
---|
109 | |
---|
110 | def write_elevation_nc(file_out, lon, lat, depth_vector): |
---|
111 | """Write an nc elevation file.""" |
---|
112 | |
---|
113 | # NetCDF file definition |
---|
114 | outfile = NetCDFFile(file_out, netcdf_mode_w) |
---|
115 | |
---|
116 | #Create new file |
---|
117 | nc_lon_lat_header(outfile, lon, lat) |
---|
118 | |
---|
119 | # ELEVATION |
---|
120 | zname = 'ELEVATION' |
---|
121 | outfile.createVariable(zname, precision, (lat_name, lon_name)) |
---|
122 | outfile.variables[zname].units = 'CENTIMETERS' |
---|
123 | outfile.variables[zname].missing_value = -1.e+034 |
---|
124 | |
---|
125 | outfile.variables[lon_name][:] = ensure_numeric(lon) |
---|
126 | outfile.variables[lat_name][:] = ensure_numeric(lat) |
---|
127 | |
---|
128 | depth = num.reshape(depth_vector, (len(lat), len(lon))) |
---|
129 | outfile.variables[zname][:] = depth |
---|
130 | |
---|
131 | outfile.close() |
---|
132 | |
---|
133 | |
---|
134 | def nc_lon_lat_header(outfile, lon, lat): |
---|
135 | """Write lat/lon headers to a NetCDF file. |
---|
136 | |
---|
137 | outfile is the netcdf file handle. |
---|
138 | lon - a list/array of the longitudes |
---|
139 | lat - a list/array of the latitudes |
---|
140 | """ |
---|
141 | |
---|
142 | outfile.institution = 'Geoscience Australia' |
---|
143 | outfile.description = 'Converted from URS binary C' |
---|
144 | |
---|
145 | # Longitude |
---|
146 | outfile.createDimension(lon_name, len(lon)) |
---|
147 | outfile.createVariable(lon_name, precision, (lon_name,)) |
---|
148 | outfile.variables[lon_name].point_spacing = 'uneven' |
---|
149 | outfile.variables[lon_name].units = 'degrees_east' |
---|
150 | outfile.variables[lon_name].assignValue(lon) |
---|
151 | |
---|
152 | # Latitude |
---|
153 | outfile.createDimension(lat_name, len(lat)) |
---|
154 | outfile.createVariable(lat_name, precision, (lat_name,)) |
---|
155 | outfile.variables[lat_name].point_spacing = 'uneven' |
---|
156 | outfile.variables[lat_name].units = 'degrees_north' |
---|
157 | outfile.variables[lat_name].assignValue(lat) |
---|
158 | |
---|
159 | |
---|
160 | |
---|
161 | |
---|
162 | def filter_netcdf(filename1, filename2, first=0, last=None, step=1): |
---|
163 | """Filter data file, selecting timesteps first:step:last. |
---|
164 | |
---|
165 | Read netcdf filename1, pick timesteps first:step:last and save to |
---|
166 | nettcdf file filename2 |
---|
167 | """ |
---|
168 | |
---|
169 | from Scientific.IO.NetCDF import NetCDFFile |
---|
170 | |
---|
171 | # Get NetCDF |
---|
172 | infile = NetCDFFile(filename1, netcdf_mode_r) #Open existing file for read |
---|
173 | outfile = NetCDFFile(filename2, netcdf_mode_w) #Open new file |
---|
174 | |
---|
175 | # Copy dimensions |
---|
176 | for d in infile.dimensions: |
---|
177 | outfile.createDimension(d, infile.dimensions[d]) |
---|
178 | |
---|
179 | # Copy variable definitions |
---|
180 | for name in infile.variables: |
---|
181 | var = infile.variables[name] |
---|
182 | outfile.createVariable(name, var.dtype.char, var.dimensions) |
---|
183 | |
---|
184 | # Copy the static variables |
---|
185 | for name in infile.variables: |
---|
186 | if name == 'time' or name == 'stage': |
---|
187 | pass |
---|
188 | else: |
---|
189 | outfile.variables[name][:] = infile.variables[name][:] |
---|
190 | |
---|
191 | # Copy selected timesteps |
---|
192 | time = infile.variables['time'] |
---|
193 | stage = infile.variables['stage'] |
---|
194 | |
---|
195 | newtime = outfile.variables['time'] |
---|
196 | newstage = outfile.variables['stage'] |
---|
197 | |
---|
198 | if last is None: |
---|
199 | last = len(time) |
---|
200 | |
---|
201 | selection = range(first, last, step) |
---|
202 | for i, j in enumerate(selection): |
---|
203 | log.critical('Copying timestep %d of %d (%f)' |
---|
204 | % (j, last-first, time[j])) |
---|
205 | newtime[i] = time[j] |
---|
206 | newstage[i,:] = stage[j,:] |
---|
207 | |
---|
208 | # Close |
---|
209 | infile.close() |
---|
210 | outfile.close() |
---|
211 | |
---|