@@ -29,19 +29,24 @@ class NCFile(FileBase):
2929 mode : str
3030 ``r``, ``w`` or ``a`` for read, write or append. Default is ``a``.
3131 clobber : bool, optional
32+ If True (default), opening a file with mode='w' will clobber an
33+ existing file with the same name. If False, an exception will be
34+ raised if a file with the same name already exists.
35+ kw : dict, optional
36+ Optional additional keyword arguments used when creating the backend
37+ file.
3238
3339 Note
3440 ----
3541 Each class instance creates one unique NetCDF4-file, with one step-counter.
3642 It is possible to store multiple fields in each file, but all snapshots of
3743 the fields must be taken at the same time. If you want one field stored
3844 every 10th timestep and another every 20th timestep, then use two different
39- class instances and as such two NetCDF4-files .
45+ class instances with two different filenames ``ncname`` .
4046 """
4147 def __init__ (self , ncname , domain = None , mode = 'a' , clobber = True , ** kw ):
42- FileBase .__init__ (self , domain = domain , ** kw )
48+ FileBase .__init__ (self , ncname , domain = domain )
4349 from netCDF4 import Dataset
44- self .filename = ncname
4550 # netCDF4 does not seem to handle 'a' if the file does not already exist
4651 if mode == 'a' and not os .path .exists (ncname ):
4752 mode = 'w'
@@ -51,11 +56,9 @@ def __init__(self, ncname, domain=None, mode='a', clobber=True, **kw):
5156 if not 'time' in self .f .variables :
5257 self .f .createDimension ('time' , None )
5358 self .f .createVariable ('time' , np .float , ('time' ))
54-
5559 self .close ()
5660
57- def _check_domain (self , write_domain , field ):
58- """Check dimensions of domain and write to file if missing"""
61+ def _check_domain (self , group , field ):
5962 N = field .global_shape [field .rank :]
6063 if self .domain is None :
6164 self .domain = []
@@ -92,22 +95,59 @@ def _check_domain(self, write_domain, field):
9295 def backend ():
9396 return 'netcdf4'
9497
95- def open (self ):
98+ def open (self , mode = 'r+' ):
9699 from netCDF4 import Dataset
97- self .f = Dataset (self .filename , mode = 'r+' , parallel = True , comm = comm )
100+ self .f = Dataset (self .filename , mode = mode , parallel = True , comm = comm )
98101
99102 def write (self , step , fields , ** kw ):
100- """Write snapshot step of ``fields`` to NetCDF4 file
103+ """Write snapshot `` step`` of ``fields`` to NetCDF4 file
101104
102105 Parameters
103106 ----------
104107 step : int
105- Index of snapshot
108+ Index of snapshot.
106109 fields : dict
107110 The fields to be dumped to file. (key, value) pairs are group name
108111 and either arrays or 2-tuples, respectively. The arrays are complete
109112 arrays to be stored, whereas 2-tuples are arrays with associated
110113 *global* slices.
114+ as_scalar : boolean, optional
115+ Whether to store rank > 0 arrays as scalars. Default is False.
116+
117+ Example
118+ -------
119+ >>> from mpi4py import MPI
120+ >>> from mpi4py_fft import PFFT, NCFile, newDistArray
121+ >>> comm = MPI.COMM_WORLD
122+ >>> T = PFFT(comm, (15, 16, 17))
123+ >>> u = newDistArray(T, forward_output=False, val=1)
124+ >>> v = newDistArray(T, forward_output=False, val=2)
125+ >>> f = NCFile('ncfilename.nc', mode='w')
126+ >>> f.write(0, {'u': [u, (u, [slice(None), 4, slice(None)])],
127+ ... 'v': [v, (v, [slice(None), 5, 5])]})
128+ >>> f.write(1, {'u': [u, (u, [slice(None), 4, slice(None)])],
129+ ... 'v': [v, (v, [slice(None), 5, 5])]})
130+ >>> f.close()
131+
132+ This stores the following datasets to the file ``ncfilename.nc``.
133+ Using in a terminal 'ncdump -h ncfilename.nc', one gets::
134+
135+ netcdf ncfilename {
136+ dimensions:
137+ time = UNLIMITED ; // (2 currently)
138+ x = 15 ;
139+ y = 16 ;
140+ z = 17 ;
141+ variables:
142+ double time(time) ;
143+ double x(x) ;
144+ double y(y) ;
145+ double z(z) ;
146+ double u(time, x, y, z) ;
147+ double u_slice_4_slice(time, x, z) ;
148+ double v(time, x, y, z) ;
149+ double v_slice_5_5(time, x) ;
150+ }
111151
112152 """
113153 self .open ()
@@ -122,17 +162,6 @@ def write(self, step, fields, **kw):
122162 self .close ()
123163
124164 def read (self , u , name , ** kw ):
125- """Read into array ``u``
126-
127- Parameters
128- ----------
129- u : array
130- The array to read into.
131- name : str
132- Name of array to be read.
133- step : int, optional
134- Index of field to be read. Default is 0.
135- """
136165 step = kw .get ('step' , 0 )
137166 self .open ()
138167 s = u .local_slice ()
@@ -141,10 +170,9 @@ def read(self, u, name, **kw):
141170 self .close ()
142171
143172 def _write_slice_step (self , name , step , slices , field , ** kw ):
144- assert name not in self .dims
173+ assert name not in self .dims # Crashes if user tries to name fields x, y, z, .
145174 rank = field .rank
146- slices = (slice (None ),)* rank + tuple (slices )
147- slices = list (slices )
175+ slices = list ((slice (None ),)* rank + tuple (slices ))
148176 slname = self ._get_slice_name (slices [rank :])
149177 s = field .local_slice ()
150178 slices , inside = self ._get_local_slices (slices , s )
@@ -168,7 +196,7 @@ def _write_slice_step(self, name, step, slices, field, **kw):
168196 self .f .sync ()
169197
170198 def _write_group (self , name , u , step , ** kw ):
171- assert name not in self .dims
199+ assert name not in self .dims # Crashes if user tries to name fields x, y, z, .
172200 s = u .local_slice ()
173201 if name not in self .f .variables :
174202 h = self .f .createVariable (name , u .dtype , self .dims )
0 commit comments