source: sasview/src/sas/sascalc/file_converter/nxcansas_writer.py @ 9ddffe0

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 9ddffe0 was 9ddffe0, checked in by lewis, 8 years ago

Correctly write orientation vectors

  • Property mode set to 100644
File size: 10.9 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
50        def _h5_float(x):
51            if not (isinstance(x, list)):
52                x = [x]
53            return np.array(x, dtype=np.float32)
54
55        def _write_h5_float(entry, value, key):
56            entry.create_dataset(key, data=_h5_float(value))
57
58        def _write_h5_vector(entry, vector, names=['x_position', 'y_position'],
59            units=None, write_fn=_write_h5_string):
60            if len(names) < 2:
61                raise ValueError("Length of names must be >= 2.")
62
63            if vector.x is not None:
64                write_fn(entry, vector.x, names[0])
65                if units is not None:
66                    entry[names[0]].attrs['units'] = units
67            if vector.y is not None:
68                write_fn(entry, vector.y, names[1])
69                if units is not None:
70                    entry[names[1]].attrs['units'] = units
71            if len(names) == 3 and vector.z is not None:
72                write_fn(entry, vector.z, names[2])
73                if units is not None:
74                    entry[names[2]].attrs['units'] = units
75
76        valid_data = all([issubclass(d.__class__, (Data1D, Data2D)) for d in dataset])
77        if not valid_data:
78            raise ValueError("All entries of dataset must be Data1D or Data2D objects")
79
80        # Get run name and number from first Data object
81        data_info = dataset[0]
82        run_number = ''
83        run_name = ''
84        if len(data_info.run) > 0:
85            run_number = data_info.run[0]
86            if len(data_info.run_name) > 0:
87                run_name = data_info.run_name[run_number]
88
89        f = h5py.File(filename, 'w')
90        sasentry = f.create_group('sasentry01')
91        sasentry['definition'] = _h5_string('NXcanSAS')
92        sasentry['run'] = _h5_string(run_number)
93        sasentry['run'].attrs['name'] = run_name
94        sasentry['title'] = _h5_string(data_info.title)
95        sasentry.attrs['canSAS_class'] = 'SASentry'
96        sasentry.attrs['version'] = '1.0'
97
98        i = 1
99
100        for data_obj in dataset:
101            data_entry = sasentry.create_group("sasdata{0:0=2d}".format(i))
102            data_entry.attrs['canSAS_class'] = 'SASdata'
103            if isinstance(data_obj, Data1D):
104                self._write_1d_data(data_obj, data_entry)
105            elif isinstance(data_obj, Data2D):
106                self._write_2d_data(data_obj, data_entry)
107            i += 1
108
109        data_info = dataset[0]
110        # Sample metadata
111        sample_entry = sasentry.create_group('sassample')
112        sample_entry.attrs['canSAS_class'] = 'SASsample'
113        sample_entry['ID'] = _h5_string(data_info.sample.name)
114        sample_attrs = ['thickness', 'temperature', 'transmission']
115        for key in sample_attrs:
116            if getattr(data_info.sample, key) is not None:
117                sample_entry.create_dataset(key,
118                    data=_h5_float(getattr(data_info.sample, key)))
119        _write_h5_vector(sample_entry, data_info.sample.position)
120        # NXcanSAS doesn't save information about pitch, only roll
121        # and yaw. The _write_h5_vector method writes vector.y, but we
122        # need to write vector.z for yaw
123        data_info.sample.orientation.y = data_info.orientation.z
124        _write_h5_vector(sample_entry, data_info.sample.orientation,
125            names=['polar_angle', 'azimuthal_angle'])
126
127        # Instrumment metadata
128        instrument_entry = sasentry.create_group('sasinstrument')
129        instrument_entry.attrs['canSAS_class'] = 'SASinstrument'
130        instrument_entry['name'] = _h5_string(data_info.instrument)
131
132        # Source metadata
133        source_entry = instrument_entry.create_group('sassource')
134        source_entry.attrs['canSAS_class'] = 'SASsource'
135        if data_info.source.radiation is None:
136            source_entry['radiation'] = _h5_string('neutron')
137        else:
138            source_entry['radiation'] = _h5_string(data_info.source.radiation)
139
140        # Collimation metadata
141        if len(data_info.collimation) > 0:
142            i = 1
143            for coll_info in data_info.collimation:
144                collimation_entry = instrument_entry.create_group(
145                    'sascollimation{0:0=2d}'.format(i))
146                collimation_entry.attrs['canSAS_class'] = 'SAScollimation'
147                if coll_info.length is not None:
148                    _write_h5_float(collimation_entry, coll_info.length, 'SDD')
149                    collimation_entry['SDD'].attrs['units'] = coll_info.length_unit
150                if coll_info.name is not None:
151                    collimation_entry['name'] = _h5_string(coll_info.name)
152        else:
153            # Create a blank one - at least 1 set of collimation metadata
154            # required by format
155            collimation_entry = instrument_entry.create_group('sascollimation01')
156
157        # Detector metadata
158        if len(data_info.detector) > 0:
159            i = 1
160            for det_info in data_info.detector:
161                detector_entry = instrument_entry.create_group(
162                    'sasdetector{0:0=2d}'.format(i))
163                detector_entry.attrs['canSAS_class'] = 'SASdetector'
164                if det_info.distance is not None:
165                    _write_h5_float(detector_entry, det_info.distance, 'SDD')
166                    detector_entry['SDD'].attrs['units'] = det_info.distance_unit
167                if det_info.name is not None:
168                    detector_entry['name'] = _h5_string(det_info.name)
169                else:
170                    detector_entry['name'] = _h5_string('')
171                if det_info.slit_length is not None:
172                    _write_h5_float(detector_entry, det_info.slit_length, 'slit_length')
173                    detector_entry['slit_length'].attrs['units'] = det_info.slit_length_unit
174                _write_h5_vector(detector_entry, det_info.offset)
175                # NXcanSAS doesn't save information about pitch, only roll
176                # and yaw. The _write_h5_vector method writes vector.y, but we
177                # need to write vector.z for yaw
178                det_info.orientation.y = det_info.orientation.z
179                _write_h5_vector(detector_entry, det_info.orientation,
180                    names=['polar_angle', 'azimuthal_angle'])
181                _write_h5_vector(detector_entry, det_info.beam_center,
182                    names=['beam_center_x', 'beam_center_y'],
183                    write_fn=_write_h5_float, units=det_info.beam_center_unit)
184                _write_h5_vector(detector_entry, det_info.pixel_size,
185                    names=['x_pixel_size', 'y_pixel_size'],
186                    write_fn=_write_h5_float, units=det_info.pixel_size_unit)
187
188                i += 1
189        else:
190            # Create a blank one - at least 1 detector required by format
191            detector_entry = instrument_entry.create_group('sasdetector01')
192            detector_entry.attrs['canSAS_class'] = 'SASdetector'
193            detector_entry.attrs['name'] = ''
194
195        # TODO: implement writing SASnote
196        i = 1
197        note_entry = sasentry.create_group('sasnote{0:0=2d}'.format(i))
198        note_entry.attrs['canSAS_class'] = 'SASnote'
199
200        f.close()
201
202    def _write_1d_data(self, data_obj, data_entry):
203        """
204        Writes the contents of a Data1D object to a SASdata h5py Group
205
206        :param data_obj: A Data1D object to write to the file
207        :param data_entry: A h5py Group object representing the SASdata
208        """
209        data_entry.attrs['signal'] = 'I'
210        data_entry.attrs['I_axes'] = 'Q'
211        data_entry.attrs['I_uncertainties'] = 'Idev'
212        data_entry.attrs['Q_indicies'] = 0
213
214        dI = data_obj.dy
215        if dI is None:
216            dI = np.zeros((data_obj.y.shape))
217
218        data_entry.create_dataset('Q', data=data_obj.x)
219        data_entry.create_dataset('I', data=data_obj.y)
220        data_entry.create_dataset('Idev', data=dI)
221
222    def _write_2d_data(self, data, data_entry):
223        """
224        Writes the contents of a Data2D object to a SASdata h5py Group
225
226        :param data: A Data2D object to write to the file
227        :param data_entry: A h5py Group object representing the SASdata
228        """
229        data_entry.attrs['signal'] = 'I'
230        data_entry.attrs['I_axes'] = 'Q,Q'
231        data_entry.attrs['I_uncertainties'] = 'Idev'
232        data_entry.attrs['Q_indicies'] = [0,1]
233
234        (n_rows, n_cols) = (len(data.y_bins), len(data.x_bins))
235
236        if n_rows == 0 and n_cols == 0:
237            # Calculate rows and columns, assuming detector is square
238            # Same logic as used in PlotPanel.py _get_bins
239            n_cols = int(np.floor(np.sqrt(len(data.qy_data))))
240            n_rows = int(np.floor(len(data.qy_data) / n_cols))
241
242            if n_rows * n_cols != len(data.qy_data):
243                raise ValueError("Unable to calculate dimensions of 2D data")
244
245        I = np.reshape(data.data, (n_rows, n_cols))
246        dI = np.zeros((n_rows, n_cols))
247        if not all(data.err_data == [None]):
248            dI = np.reshape(data.err_data, (n_rows, n_cols))
249        qx =  np.reshape(data.qx_data, (n_rows, n_cols))
250        qy = np.reshape(data.qy_data, (n_rows, n_cols))
251
252        I_entry = data_entry.create_dataset('I', data=I)
253        I_entry.attrs['units'] = data.I_unit
254        Qx_entry = data_entry.create_dataset('Qx', data=qx)
255        Qx_entry.attrs['units'] = data.Q_unit
256        Qy_entry = data_entry.create_dataset('Qy', data=qy)
257        Qy_entry.attrs['units'] = data.Q_unit
258        Idev_entry = data_entry.create_dataset('Idev', data=dI)
259        Idev_entry.attrs['units'] = data.I_unit
Note: See TracBrowser for help on using the repository browser.