source: sasview/src/sas/sascalc/file_converter/nxcansas_writer.py @ 0836f77

ESS_GUIESS_GUI_DocsESS_GUI_batch_fittingESS_GUI_bumps_abstractionESS_GUI_iss1116ESS_GUI_iss879ESS_GUI_iss959ESS_GUI_openclESS_GUI_orderingESS_GUI_sync_sascalccostrafo411magnetic_scattrelease-4.1.1release-4.1.2release-4.2.2ticket-1009ticket-1094-headlessticket-1242-2d-resolutionticket-1243ticket-1249ticket885unittest-saveload
Last change on this file since 0836f77 was 0836f77, checked in by lewis, 8 years ago

Write source metadata to NXcanSAS

  • Property mode set to 100644
File size: 12.6 KB
Line 
1"""
2    NXcanSAS 1/2D data reader for writing HDF5 formatted NXcanSAS files.
3"""
4
5import h5py
6import numpy as np
7import re
8import os
9
10from sas.sascalc.dataloader.readers.cansas_reader_HDF5 import Reader as Cansas2Reader
11from sas.sascalc.dataloader.data_info import Data1D, Data2D
12
13class NXcanSASWriter(Cansas2Reader):
14    """
15    A class for writing in NXcanSAS data files. Any number of data sets may be
16    written to the file. Currently 1D and 2D SAS data sets are supported
17
18    NXcanSAS spec: http://download.nexusformat.org/sphinx/classes/contributed_definitions/NXcanSAS.html
19
20    :Dependencies:
21        The NXcanSAS writer requires h5py => v2.5.0 or later.
22    """
23
24    def write(self, dataset, filename):
25        """
26        Write an array of Data1d or Data2D objects to an NXcanSAS file, as
27        one SASEntry with multiple SASData elements. The metadata of the first
28        elememt in the array will be written as the SASentry metadata
29        (detector, instrument, sample, etc).
30
31        :param dataset: A list of Data1D or Data2D objects to write
32        :param filename: Where to write the NXcanSAS file
33        """
34
35        def _h5_string(string):
36            """
37            Convert a string to a numpy string in a numpy array. This way it is
38            written to the HDF5 file as a fixed length ASCII string and is
39            compatible with the Reader read() method.
40            """
41            if not isinstance(string, str):
42                string = str(string)
43
44            return np.array([np.string_(string)])
45
46        def _write_h5_string(entry, value, key):
47            entry[key] = _h5_string(value)
48
49        def _h5_float(x):
50            if not (isinstance(x, list)):
51                x = [x]
52            return np.array(x, dtype=np.float32)
53
54        def _write_h5_float(entry, value, key):
55            entry.create_dataset(key, data=_h5_float(value))
56
57        def _write_h5_vector(entry, vector, names=['x_position', 'y_position'],
58            units=None, write_fn=_write_h5_string):
59            """
60            Write a vector to an h5 entry
61
62            :param entry: The H5Py entry to write to
63            :param vector: The Vector to write
64            :param names: What to call the x,y and z components of the vector
65                when writing to the H5Py entry
66            :param units: The units of the vector (optional)
67            :param write_fn: A function to convert the value to the required
68                format and write it to the H5Py entry, of the form
69                f(entry, value, name) (optional)
70            """
71            if len(names) < 2:
72                raise ValueError("Length of names must be >= 2.")
73
74            if vector.x is not None:
75                write_fn(entry, vector.x, names[0])
76                if units is not None:
77                    entry[names[0]].attrs['units'] = units
78            if vector.y is not None:
79                write_fn(entry, vector.y, names[1])
80                if units is not None:
81                    entry[names[1]].attrs['units'] = units
82            if len(names) == 3 and vector.z is not None:
83                write_fn(entry, vector.z, names[2])
84                if units is not None:
85                    entry[names[2]].attrs['units'] = units
86
87        valid_data = all([issubclass(d.__class__, (Data1D, Data2D)) for d in dataset])
88        if not valid_data:
89            raise ValueError("All entries of dataset must be Data1D or Data2D objects")
90
91        # Get run name and number from first Data object
92        data_info = dataset[0]
93        run_number = ''
94        run_name = ''
95        if len(data_info.run) > 0:
96            run_number = data_info.run[0]
97            if len(data_info.run_name) > 0:
98                run_name = data_info.run_name[run_number]
99
100        f = h5py.File(filename, 'w')
101        sasentry = f.create_group('sasentry01')
102        sasentry['definition'] = _h5_string('NXcanSAS')
103        sasentry['run'] = _h5_string(run_number)
104        sasentry['run'].attrs['name'] = run_name
105        sasentry['title'] = _h5_string(data_info.title)
106        sasentry.attrs['canSAS_class'] = 'SASentry'
107        sasentry.attrs['version'] = '1.0'
108
109        i = 1
110
111        for data_obj in dataset:
112            data_entry = sasentry.create_group("sasdata{0:0=2d}".format(i))
113            data_entry.attrs['canSAS_class'] = 'SASdata'
114            if isinstance(data_obj, Data1D):
115                self._write_1d_data(data_obj, data_entry)
116            elif isinstance(data_obj, Data2D):
117                self._write_2d_data(data_obj, data_entry)
118            i += 1
119
120        data_info = dataset[0]
121        # Sample metadata
122        sample_entry = sasentry.create_group('sassample')
123        sample_entry.attrs['canSAS_class'] = 'SASsample'
124        sample_entry['ID'] = _h5_string(data_info.sample.name)
125        sample_attrs = ['thickness', 'temperature', 'transmission']
126        for key in sample_attrs:
127            if getattr(data_info.sample, key) is not None:
128                sample_entry.create_dataset(key,
129                    data=_h5_float(getattr(data_info.sample, key)))
130        _write_h5_vector(sample_entry, data_info.sample.position)
131        # NXcanSAS doesn't save information about pitch, only roll
132        # and yaw. The _write_h5_vector method writes vector.y, but we
133        # need to write vector.z for yaw
134        data_info.sample.orientation.y = data_info.sample.orientation.z
135        _write_h5_vector(sample_entry, data_info.sample.orientation,
136            names=['polar_angle', 'azimuthal_angle'])
137        if data_info.sample.details is not None\
138            and data_info.sample.details != []:
139            details = reduce(lambda x,y: x + "\n" + y, data_info.sample.details)
140            sample_entry['details'] = _h5_string(details)
141
142        # Instrumment metadata
143        instrument_entry = sasentry.create_group('sasinstrument')
144        instrument_entry.attrs['canSAS_class'] = 'SASinstrument'
145        instrument_entry['name'] = _h5_string(data_info.instrument)
146
147        # Source metadata
148        source_entry = instrument_entry.create_group('sassource')
149        source_entry.attrs['canSAS_class'] = 'SASsource'
150        if data_info.source.radiation is None:
151            source_entry['radiation'] = _h5_string('neutron')
152        else:
153            source_entry['radiation'] = _h5_string(data_info.source.radiation)
154        if data_info.source.beam_shape is not None:
155            source_entry['beam_shape'] = _h5_string(data_info.source.beam_shape)
156        wavelength_keys = { 'wavelength': 'incident_wavelength',
157            'wavelength_min':'wavelength_min',
158            'wavelength_max': 'wavelength_max',
159            'wavelength_spread': 'incident_wavelength_spread' }
160        for sasname, nxname in wavelength_keys.iteritems():
161            value = getattr(data_info.source, sasname)
162            units = getattr(data_info.source, sasname + '_unit')
163            if value is not None:
164                source_entry[nxname] = _h5_float(value)
165                source_entry[nxname].attrs['units'] = units
166        _write_h5_vector(source_entry, data_info.source.beam_size,
167            names=['beam_size_x', 'beam_size_y'],
168            units=data_info.source.beam_size_unit, write_fn=_write_h5_float)
169
170
171        # Collimation metadata
172        if len(data_info.collimation) > 0:
173            i = 1
174            for coll_info in data_info.collimation:
175                collimation_entry = instrument_entry.create_group(
176                    'sascollimation{0:0=2d}'.format(i))
177                collimation_entry.attrs['canSAS_class'] = 'SAScollimation'
178                if coll_info.length is not None:
179                    _write_h5_float(collimation_entry, coll_info.length, 'SDD')
180                    collimation_entry['SDD'].attrs['units'] = coll_info.length_unit
181                if coll_info.name is not None:
182                    collimation_entry['name'] = _h5_string(coll_info.name)
183        else:
184            # Create a blank one - at least 1 set of collimation metadata
185            # required by format
186            collimation_entry = instrument_entry.create_group('sascollimation01')
187
188        # Detector metadata
189        if len(data_info.detector) > 0:
190            i = 1
191            for det_info in data_info.detector:
192                detector_entry = instrument_entry.create_group(
193                    'sasdetector{0:0=2d}'.format(i))
194                detector_entry.attrs['canSAS_class'] = 'SASdetector'
195                if det_info.distance is not None:
196                    _write_h5_float(detector_entry, det_info.distance, 'SDD')
197                    detector_entry['SDD'].attrs['units'] = det_info.distance_unit
198                if det_info.name is not None:
199                    detector_entry['name'] = _h5_string(det_info.name)
200                else:
201                    detector_entry['name'] = _h5_string('')
202                if det_info.slit_length is not None:
203                    _write_h5_float(detector_entry, det_info.slit_length, 'slit_length')
204                    detector_entry['slit_length'].attrs['units'] = det_info.slit_length_unit
205                _write_h5_vector(detector_entry, det_info.offset)
206                # NXcanSAS doesn't save information about pitch, only roll
207                # and yaw. The _write_h5_vector method writes vector.y, but we
208                # need to write vector.z for yaw
209                det_info.orientation.y = det_info.orientation.z
210                _write_h5_vector(detector_entry, det_info.orientation,
211                    names=['polar_angle', 'azimuthal_angle'])
212                _write_h5_vector(detector_entry, det_info.beam_center,
213                    names=['beam_center_x', 'beam_center_y'],
214                    write_fn=_write_h5_float, units=det_info.beam_center_unit)
215                _write_h5_vector(detector_entry, det_info.pixel_size,
216                    names=['x_pixel_size', 'y_pixel_size'],
217                    write_fn=_write_h5_float, units=det_info.pixel_size_unit)
218
219                i += 1
220        else:
221            # Create a blank one - at least 1 detector required by format
222            detector_entry = instrument_entry.create_group('sasdetector01')
223            detector_entry.attrs['canSAS_class'] = 'SASdetector'
224            detector_entry.attrs['name'] = ''
225
226        # TODO: implement writing SASnote
227        i = 1
228        note_entry = sasentry.create_group('sasnote{0:0=2d}'.format(i))
229        note_entry.attrs['canSAS_class'] = 'SASnote'
230
231        f.close()
232
233    def _write_1d_data(self, data_obj, data_entry):
234        """
235        Writes the contents of a Data1D object to a SASdata h5py Group
236
237        :param data_obj: A Data1D object to write to the file
238        :param data_entry: A h5py Group object representing the SASdata
239        """
240        data_entry.attrs['signal'] = 'I'
241        data_entry.attrs['I_axes'] = 'Q'
242        data_entry.attrs['I_uncertainties'] = 'Idev'
243        data_entry.attrs['Q_indicies'] = 0
244
245        dI = data_obj.dy
246        if dI is None:
247            dI = np.zeros((data_obj.y.shape))
248
249        data_entry.create_dataset('Q', data=data_obj.x)
250        data_entry.create_dataset('I', data=data_obj.y)
251        data_entry.create_dataset('Idev', data=dI)
252
253    def _write_2d_data(self, data, data_entry):
254        """
255        Writes the contents of a Data2D object to a SASdata h5py Group
256
257        :param data: A Data2D object to write to the file
258        :param data_entry: A h5py Group object representing the SASdata
259        """
260        data_entry.attrs['signal'] = 'I'
261        data_entry.attrs['I_axes'] = 'Q,Q'
262        data_entry.attrs['I_uncertainties'] = 'Idev'
263        data_entry.attrs['Q_indicies'] = [0,1]
264
265        (n_rows, n_cols) = (len(data.y_bins), len(data.x_bins))
266
267        if n_rows == 0 and n_cols == 0:
268            # Calculate rows and columns, assuming detector is square
269            # Same logic as used in PlotPanel.py _get_bins
270            n_cols = int(np.floor(np.sqrt(len(data.qy_data))))
271            n_rows = int(np.floor(len(data.qy_data) / n_cols))
272
273            if n_rows * n_cols != len(data.qy_data):
274                raise ValueError("Unable to calculate dimensions of 2D data")
275
276        I = np.reshape(data.data, (n_rows, n_cols))
277        dI = np.zeros((n_rows, n_cols))
278        if not all(data.err_data == [None]):
279            dI = np.reshape(data.err_data, (n_rows, n_cols))
280        qx =  np.reshape(data.qx_data, (n_rows, n_cols))
281        qy = np.reshape(data.qy_data, (n_rows, n_cols))
282
283        I_entry = data_entry.create_dataset('I', data=I)
284        I_entry.attrs['units'] = data.I_unit
285        Qx_entry = data_entry.create_dataset('Qx', data=qx)
286        Qx_entry.attrs['units'] = data.Q_unit
287        Qy_entry = data_entry.create_dataset('Qy', data=qy)
288        Qy_entry.attrs['units'] = data.Q_unit
289        Idev_entry = data_entry.create_dataset('Idev', data=dI)
290        Idev_entry.attrs['units'] = data.I_unit
Note: See TracBrowser for help on using the repository browser.