source: trunk/bindings/python/nxs/tree.py @ 1822

Revision 1819, 118.2 KB checked in by Ray Osborn, 3 months ago (diff)

Refs #322: Mostly tidying things up.

  • Updated exceptions to use preferred Python syntax and removed some redundant TypeError traps.
  • Replaced spaces in NXobject names by underscores.
  • Removed check for the dimension products being negative due to overflows.
  • In log plots, the minimum value of 0.1 is only used for integer data, where it is a reasonable replacement for values of 0. In floating point data, the default matplotlib behaviour is restored (but can be overridden by using the zmin keyword argument).
Line 
1#!/usr/bin/env python
2# This program is public domain
3# Author: Paul Kienzle, Ray Osborn
4
5"""
6NeXus data as Python trees
7==========================
8The `nexus.tree` modules are designed to accomplish two goals:
9
10    1. To provide convenient access to existing data contained in NeXus files.
11    2. To enable new NeXus data to be created and manipulated interactively.
12
13These goals are achieved by mapping hierarchical NeXus data structures directly
14into python objects, which either represent NeXus groups or NeXus fields.
15Entries in a group are referenced much like fields in a class are referenced in
16python. The entire data hierarchy can be referenced at any time, whether the
17NeXus data has been loaded in from an existing NeXus file or created dynamically
18within the python session. This provides a much more natural scripting interface
19to NeXus data than the directory model of the `nexus.napi` interface.
20
21Example 1: Loading a NeXus file
22-------------------------------
23The following commands loads NeXus data from a file, displays (some of) the
24contents as a tree, and then accesses individual data items.
25
26    >>> from nexpy.api import nexus as nx
27    >>> a=nx.load('sns/data/ARCS_7326.nxs')
28    >>> print a.tree
29    root:NXroot
30      @HDF5_Version = 1.8.2
31      @NeXus_version = 4.2.1
32      @file_name = ARCS_7326.nxs
33      @file_time = 2010-05-05T01:59:25-05:00
34      entry:NXentry
35        data:NXdata
36          data = float32(631x461x4x825)
37            @axes = rotation_angle:tilt_angle:sample_angle:time_of_flight
38            @signal = 1
39          rotation_angle = float32(632)
40            @units = degree
41          sample_angle = [ 210.  215.  220.  225.  230.]
42            @units = degree
43          tilt_angle = float32(462)
44            @units = degree
45          time_of_flight = float32(826)
46            @units = microsecond
47        run_number = 7326
48        sample:NXsample
49          pulse_time = 2854.94747365
50            @units = microsecond
51    .
52    .
53    .
54    >>> a.entry.run_number
55    NXfield(7326)
56
57So the tree returned from load() has an entry for each group, field and
58attribute.  You can traverse the hierarchy using the names of the groups.  For
59example, tree.entry.instrument.detector.distance is an example of a field
60containing the distance to each pixel in the detector. Entries can also be
61referenced by NXclass name, such as tree.NXentry[0].instrument. Since there may
62be multiple entries of the same NeXus class, the NXclass attribute returns a
63(possibly empty) list.
64
65The load() and save() functions are implemented using the class
66`nexus.tree.NeXusTree`, a subclass of `nexus.napi.NeXus` which allows all the
67usual API functions.
68
69Example 2: Creating a NeXus file dynamically
70--------------------------------------------
71The second example shows how to create NeXus data dynamically and saves it to a
72file. The data are first created as Numpy arrays
73
74    >>> import numpy as np
75    >>> x=y=np.linspace(0,2*np.pi,101)
76    >>> X,Y=np.meshgrid(y,x)
77    >>> z=np.sin(X)*np.sin(Y)
78
79Then a NeXus data groups are created and the data inserted to produce a
80NeXus-compliant structure that can be saved to a file.
81
82    >>> root=nx.NXroot(NXentry())
83    >>> print root.tree
84    root:NXroot
85      entry:NXentry
86    >>> root.entry.data=nx.NXdata(z,[x,y])
87
88Additional metadata can be inserted before saving the data to a file.
89
90    >>> root.entry.sample=nx.NXsample()
91    >>> root.entry.sample.temperature = 40.0
92    >>> root.entry.sample.temperature.units = 'K'
93    >>> root.save('example.nxs')
94
95NXfield objects have much of the functionality of Numpy arrays. They may be used
96in simple arithmetic expressions with other NXfields, Numpy arrays or scalar
97values and will be cast as ndarray objects if used as arguments in Numpy
98modules.
99
100    >>> x=nx.NXfield(np.linspace(0,10.0,11))
101    >>> x
102    NXfield([  0.   1.   2. ...,   8.   9.  10.])
103    >>> x + 10
104    NXfield([ 10.  11.  12. ...,  18.  19.  20.])
105    >>> sin(x)
106    array([ 0.        ,  0.84147098,  0.90929743, ...,  0.98935825,
107        0.41211849, -0.54402111])
108
109If the arithmetic operation is assigned to a NeXus group attribute, it will be
110automatically cast as a valid NXfield object with the type and shape determined
111by the Numpy array type and shape.
112
113    >>> entry.data.result = sin(x)
114    >>> entry.data.result
115    NXfield([ 0.          0.84147098  0.90929743 ...,  0.98935825  0.41211849
116     -0.54402111])
117    >>> entry.data.result.dtype, entry.data.result.shape
118    (dtype('float64'), (11,))
119
120NeXus Objects
121-------------
122Properties of the entry in the tree are referenced by attributes that depend
123on the object type, different nx attributes may be available.
124
125Objects (class NXobject) have attributes shared by both groups and fields::
126    * nxname   object name
127    * nxclass  object class for groups, 'NXfield' for fields
128    * nxgroup  group containing the entry, or None for the root
129    * attrs    dictionary of NeXus attributes for the object
130
131Groups (class NXgroup) have attributes for accessing children::
132    * entries  dictionary of entries within the group
133    * component('nxclass')  return group entries of a particular class
134    * dir()    print the list of entries in the group
135    * tree     return the list of entries and subentries in the group
136    * plot()   plot signal and axes for the group, if available
137
138Fields (class NXfield) have attributes for accessing data:
139    * shape    dimensions of data in the field
140    * dtype    data type
141    * nxdata   data in the field
142
143Linked fields or groups (class NXlink) have attributes for accessing the link::
144    * nxlink   reference to the linked field or group
145
146NeXus attributes (class NXattr) have a type and a value only::
147    * dtype    attribute type
148    * nxdata   attribute data
149
150There is a subclass of NXgroup for each group class defined by the NeXus standard,
151so it is possible to create an NXgroup of NeXus class NXsample directly using:
152
153    >>> sample = NXsample()
154
155The default group name will be the class name following the 'NX', so the above
156group will have an nxname of 'sample'. However, this is overridden by the
157attribute name when it is assigned as a group attribute, e.g.,
158
159    >>> entry.sample1 = NXsample()
160    >>> entry.sample1.nxname
161    sample1
162
163You can traverse the tree by component class instead of component name. Since
164there may be multiple components of the same class in one group you will need to
165specify which one to use.  For example,
166
167    tree.NXentry[0].NXinstrument[0].NXdetector[0].distance
168
169references the first detector of the first instrument of the first entry.
170Unfortunately, there is no guarantee regarding the order of the entries, and it
171may vary from call to call, so this is mainly useful in iterative searches.
172
173
174Unit Conversion
175---------------
176Data can be stored in the NeXus file in a variety of units, depending on which
177facility is storing the file.  This makes life difficult for reduction and
178analysis programs which must know the units they are working with.  Our solution
179to this problem is to allow the reader to retrieve data from the file in
180particular units.  For example, if detector distance is stored in the file using
181millimeters you can retrieve them in meters using::
182
183    entry.instrument.detector.distance.convert('m')
184
185See `nexus.unit` for more details on the unit formats supported.
186
187Reading and Writing Slabs
188-------------------------
189The slab interface to field data works by opening the file handle and keeping it
190open as long as the slab interface is needed.  This is done in python 2.5 using
191the with statement.  Once the context is entered, get() and put() methods on the
192object allow you to read and write data a slab at a time.  For example::
193
194    # Read a Ni x Nj x Nk array one vector at a time
195    with root.NXentry[0].data.data as slab:
196        Ni,Nj,Nk = slab.shape
197        size = [1,1,Nk]
198        for i in range(Ni):
199            for j in range(Nj):
200                value = slab.get([i,j,0],size)
201
202The equivalent can be done in Python 2.4 and lower using the context
203functions __enter__ and __exit__::
204
205    slab = data.slab.__enter__()
206    ... do the slab functions ...
207    data.slab.__exit__()
208
209Plotting NeXus data
210-------------------
211There is a plot() method for groups that automatically looks for 'signal' and
212'axes' attributes within the group in order to determine what to plot. These are
213defined by the 'nxsignal' and 'nxaxes' properties of the group. This means that
214the method will determine whether the plot should be one- or two- dimensional.
215For higher than two dimensions, only the top slice is plotted by default.
216
217The plot method accepts as arguments the standard matplotlib.pyplot.plot format
218strings to customize one-dimensional plots, axis and scale limits, and will
219transmit keyword arguments to the matplotlib plotting methods.
220
221    >>> a=nx.load('chopper.nxs')
222    >>> a.entry.monitor1.plot()
223    >>> a.entry.monitor2.plot('r+', xmax=2600)
224   
225It is possible to plot over the existing figure with the oplot() method and to
226plot with logarithmic intensity scales with the logplot() method. The x- and
227y-axes can also be rendered logarithmically using the logx and logy keywards.
228
229Although the plot() method uses matplotlib by default to plot the data, you can replace
230this with your own plotter by setting nexus.NXgroup._plotter to your own plotter
231class.  The plotter class has one method::
232
233    plot(signal, axes, entry, title, format, **opts)
234
235where signal is the field containing the data, axes are the fields listing the
236signal sample points, entry is file/path within the file to the data group and
237title is the title of the group or the parent NXentry, if available.
238"""
239from __future__ import with_statement
240from copy import copy, deepcopy
241
242import numpy as np
243import napi
244from napi import NeXusError
245
246#Memory in MB
247NX_MEMORY = 500
248
249__all__ = ['NeXusTree', 'NXobject', 'NXfield', 'NXgroup', 'NXattr',
250           'NX_MEMORY', 'setmemory', 'load', 'save', 'tree', 'centers',
251           'NXlink', 'NXlinkfield', 'NXlinkgroup', 'SDS', 'NXlinkdata']
252
253#List of defined base classes (later added to __all__)
254_nxclasses = ['NXroot', 'NXentry', 'NXsubentry', 'NXdata', 'NXmonitor',
255              'NXlog', 'NXsample', 'NXinstrument', 'NXaperture', 'NXattenuator',
256              'NXbeam', 'NXbeam_stop', 'NXbending_magnet', 'NXcharacterization',
257              'NXcollection', 'NXcollimator', 'NXcrystal', 'NXdetector',
258              'NXdisk_chopper', 'NXenvironment', 'NXevent_data',
259              'NXfermi_chopper', 'NXfilter', 'NXflipper', 'NXgeometry',
260              'NXguide', 'NXinsertion_device', 'NXmirror', 'NXmoderator',
261              'NXmonochromator', 'NXnote', 'NXorientation', 'NXparameter',
262              'NXpolarizer', 'NXpositioner', 'NXprocess', 'NXsensor', 'NXshape',
263              'NXsource', 'NXtranslation', 'NXuser', 'NXvelocity_selector', 'NXtree']
264
265np.set_printoptions(threshold=5)
266
267class NeXusTree(napi.NeXus):
268
269    """
270    Structure-based interface to the NeXus file API.
271
272    Usage::
273
274      file = NeXusTree(filename, ['r','rw','w'])
275        - open the NeXus file
276      root = file.readfile()
277        - read the structure of the NeXus file.  This returns a NeXus tree.
278      file.writefile(root)
279        - write a NeXus tree to the file.
280      data = file.readpath(path)
281        - read data from a particular path
282
283    Example::
284
285      nx = NeXusTree('REF_L_1346.nxs','r')
286      tree = nx.readfile()
287      for entry in tree.NXentry:
288          process(entry)
289      copy = NeXusTree('modified.nxs','w')
290      copy.writefile(tree)
291
292    Note that the large datasets are not loaded immediately.  Instead, the
293    when the data set is requested, the file is reopened, the data read, and
294    the file closed again.  open/close are available for when we want to
295    read/write slabs without the overhead of moving the file cursor each time.
296    The NXdata objects in the returned tree hold the object values.
297    """
298
299    def readfile(self):
300        """
301        Read the NeXus file structure from the file and return a tree of NXobjects.
302
303        Large datasets are not read until they are needed.
304        """
305        self.open()
306        self.openpath("/")
307        root = self._readgroup()
308        self.close()
309        root._group = None
310        # Resolve links (not necessary now that link is set as a property)
311        #self._readlinks(root, root)
312        root._file = self
313        return root
314
315    def writefile(self, tree):
316        """
317        Write the NeXus file structure to a file.
318
319        The file is assumed to start empty. Updating individual objects can be
320        done using the napi interface, with nx.handle as the nexus file handle.
321        """
322        self.open()
323        links = []
324        for entry in tree.entries.values():
325            links += self._writegroup(entry, path="")
326        self._writelinks(links)
327        self.close()
328
329    def readpath(self, path):
330        """
331        Return the data on a particular file path.
332
333        Returns a numpy array containing the data, a python scalar, or a
334        string depending on the shape and storage class.
335        """
336        self.open()
337        self.openpath(path)
338        try:
339            return self.getdata()
340        except ValueError:
341            return None
342
343    def _readdata(self, name):
344        """
345        Read a data object and return it as an NXfield or NXlink.
346        """
347        # Finally some data, but don't read it if it is big
348        # Instead record the location, type and size
349        self.opendata(name)
350        attrs={}
351        attrs = self.getattrs()
352        if 'target' in attrs and attrs['target'] != self.path:
353            # This is a linked dataset; don't try to load it.
354            data = NXlinkfield(target=attrs['target'], name=name)
355        else:
356            dims,type = self.getinfo()
357            #Read in the data if it's not too large
358            if np.prod(dims) < 1000:# i.e., less than 1k dims
359                try:
360                    value = self.getdata()
361                except ValueError:
362                    value = None
363            else:
364                value = None
365            data = NXfield(value=value,name=name,dtype=type,shape=dims,attrs=attrs)
366        self.closedata()
367        data._infile = data._saved = data._changed = True
368        return data
369
370    # These are groups that HDFView explicitly skips
371    _skipgroups = ['CDF0.0','_HDF_CHK_TBL_','Attr0.0','RIG0.0','RI0.0',
372                   'RIATTR0.0N','RIATTR0.0C']
373
374    def _readchildren(self,n):
375        children = {}
376        for _item in range(n):
377            name,nxclass = self.getnextentry()
378            if nxclass in self._skipgroups:
379                pass # Skip known bogus classes
380            elif nxclass == 'SDS': # NXgetnextentry returns 'SDS' as the class for NXfields
381                children[name] = self._readdata(name)
382            else:
383                self.opengroup(name,nxclass)
384                children[name] = self._readgroup()
385                self.closegroup()
386        return children
387
388    def _readgroup(self):
389        """
390        Read the currently open group and return it as an NXgroup.
391        """
392        n,name,nxclass = self.getgroupinfo()
393        attrs = {}
394        attrs = self.getattrs()
395        if 'target' in attrs and attrs['target'] != self.path:
396            # This is a linked group; don't try to load it.
397            group = NXlinkgroup(target=attrs['target'], name=name)
398        else:
399            children = self._readchildren(n)
400            # If we are subclassed with a handler for the particular
401            # NXentry class name use that constructor for the group
402            # rather than the generic NXgroup class.
403            group = NXgroup(nxclass=nxclass,name=name,attrs=attrs,entries=children)
404            # Build chain back structure
405            for obj in children.values():
406                obj._group = group
407        group._infile = group._saved = group._changed = True
408        return group
409
410    def _readlinks(self, root, group):
411        """
412        Convert linked objects into direct references.
413        """
414        for entry in group.entries.values():
415            if isinstance(entry, NXlink):
416                link = root
417                try:
418                    for level in entry._target[1:].split('/'):
419                        link = getattr(link,level)
420                    entry.nxlink = link
421                except AttributeError:
422                    pass
423            elif isinstance(entry, NXgroup):
424                self._readlinks(root, entry)
425
426    def _writeattrs(self, attrs):
427        """
428        Return the attributes for the currently open group/data.
429
430        If no group or data object is open, the file attributes are returned.
431        """
432        for name,pair in attrs.iteritems():
433            self.putattr(name,pair.nxdata,pair.dtype)
434
435    def _writedata(self, data, path):
436        """
437        Write the given data to a file.
438
439        NXlinks cannot be written until the linked group is created, so
440        this routine returns the set of links that need to be written.
441        Call writelinks on the list.
442        """
443
444        path = path + "/" + data.nxname
445
446        # If the data is linked then
447        if hasattr(data,'_target'):
448            return [(path, data._target)]
449
450        shape = data.shape
451        if shape == (): shape = (1,)
452
453        #If the array size is too large, their product needs a long integer
454        if np.prod(shape) > 10000:
455            # Compress the fastest moving dimension of large datasets
456            slab_dims = np.ones(len(shape),'i')
457            if shape[-1] < 100000:
458                slab_dims[-1] = shape[-1]
459            else:
460                slab_dims[-1] = 100000
461            self.compmakedata(data.nxname, data.dtype, shape, 'lzw', slab_dims)
462        else:
463            # Don't use compression for small datasets
464            try:
465                self.makedata(data.nxname, data.dtype, shape)
466            except StandardError,errortype:
467                print "Error in tree, makedata: ",errortype
468
469        self.opendata(data.nxname)
470        self._writeattrs(data.attrs)
471        value = data.nxdata
472        if value is not None:
473            self.putdata(data.nxdata)
474        self.closedata()
475        return []
476
477    def _writegroup(self, group, path):
478        """
479        Write the given group structure, including the data.
480
481        NXlinks cannot be written until the linked group is created, so
482        this routine returns the set of links that need to be written.
483        Call writelinks on the list.
484        """
485        path = path + "/" + group.nxname
486
487        links = []
488        self.makegroup(group.nxname, group.nxclass)
489        self.opengroup(group.nxname, group.nxclass)
490        self._writeattrs(group.attrs)
491        if hasattr(group, '_target'):
492            links += [(path, group._target)]
493        for child in group.entries.values():
494            if child.nxclass == 'NXfield':
495                links += self._writedata(child,path)
496            elif hasattr(child,'_target'):
497                links += [(path+"/"+child.nxname,child._target)]
498            else:
499                links += self._writegroup(child,path)
500        self.closegroup()
501        return links
502
503    def _writelinks(self, links):
504        """
505        Create links within the NeXus file.
506
507        THese are defined by the set of pairs returned by _writegroup.
508        """
509        gid = {}
510
511        # identify targets
512        for path,target in links:
513            gid[target] = None
514
515        # find gids for targets
516        for target in gid.iterkeys():
517            self.openpath(target)
518            # Can't tell from the name if we are linking to a group or
519            # to a dataset, so cheat and rely on getdataID to signal
520            # an error if we are not within a group.
521            try:
522                gid[target] = self.getdataID()
523            except NeXusError:
524                gid[target] = self.getgroupID()
525
526        # link sources to targets
527        for path,target in links:
528            if path != target:
529                # ignore self-links
530                parent = "/".join(path.split("/")[:-1])
531                self.openpath(parent)
532                self.makelink(gid[target])
533
534
535def _readaxes(axes):
536    """
537    Return a list of axis names stored in the 'axes' attribute.
538
539    The delimiter separating each axis can be white space, a comma, or a colon.
540    """
541    import re
542    sep=re.compile('[\[]*(\s*,*:*)+[\]]*')
543    return filter(lambda x: len(x)>0, sep.split(axes))
544
545
546class AttrDict(dict):
547
548    """
549    A dictionary class to assign all attributes to the NXattr class.
550    """
551
552    def __setitem__(self, key, value):
553        if isinstance(value, NXattr):
554            dict.__setitem__(self, key, value)
555        else:
556            dict.__setitem__(self, key, NXattr(value))
557
558
559class NXattr(object):
560
561    """
562    Class for NeXus attributes of a NXfield or NXgroup object.
563
564    This class is only used for NeXus attributes that are stored in a
565    NeXus file and helps to distinguish them from Python attributes.
566    There are two Python attributes for each NeXus attribute.
567
568    Python Attributes
569    -----------------
570    nxdata : string, Numpy scalar, or Numpy ndarray
571        The value of the NeXus attribute.
572    dtype : string
573        The data type of the NeXus attribute. This is set to 'char' for
574        a string attribute or the string of the corresponding Numpy data type
575        for a numeric attribute.
576
577    NeXus Attributes
578    ----------------
579    NeXus attributes are stored in the 'attrs' dictionary of the parent object,
580    NXfield or NXgroup, but can often be referenced or assigned using the
581    attribute name as if it were an object attribute.
582
583    For example, after assigning the NXfield, the following three attribute
584    assignments are all equivalent::
585
586        >>> entry.sample.temperature = NXfield(40.0)
587        >>> entry.sample.temperature.attrs['units'] = 'K'
588        >>> entry.sample.temperature.units = NXattr('K')
589        >>> entry.sample.temperature.units = 'K'
590
591    The third version above is only allowed for NXfield attributes and is
592    not allowed if the attribute has the same name as one of the following
593    internally defined attributes, i.e.,
594
595    ['entries', 'attrs', 'dtype','shape']
596
597    or if the attribute name begins with 'nx' or '_'. It is only possible to
598    reference attributes with one of the proscribed names using the 'attrs'
599    dictionary.
600
601    """
602
603    def __init__(self,value=None,dtype=''):
604        if isinstance(value, NXattr):
605            self._data,self._dtype = value.nxdata,value.dtype
606        elif dtype:
607            if dtype in np.typeDict:
608                self._data,self._dtype = np.__dict__[dtype](value),dtype
609            elif dtype == 'char':
610                self._data,self._dtype = str(value),dtype
611            else:
612                raise NeXusError("Invalid data type")
613        else:
614            if isinstance(value, str):
615                self._data,self._dtype = str(value), 'char'
616            elif value is not None:
617                if isinstance(value, NXobject):
618                    raise NeXusError("A data attribute cannot be a NXfield or NXgroup")
619                else:
620                    self._data = np.array(value)
621                self._dtype = self._data.dtype.name
622                if self._data.size == 1:
623                    self._data = np.__dict__[self._dtype](self._data)
624            else:
625                self._data,self._dtype = None, 'char'
626
627    def __str__(self):
628        return str(self.nxdata)
629
630    def __repr__(self):
631        if str(self.dtype) == 'char':
632            return "NXattr('%s')"%self.nxdata
633        else:
634            return "NXattr(%s)"%self.nxdata
635
636    def __eq__(self, other):
637        """
638        Return true if the value of the attribute is the same as the other.
639        """
640        if isinstance(other, NXattr):
641            return self.nxdata == other.nxdata
642        else:
643            return self.nxdata == other
644
645    def _getdata(self):
646        """
647        Return the attribute value.
648        """
649        return self._data
650
651    def _getdtype(self):
652        return self._dtype
653
654    nxdata = property(_getdata,doc="The attribute values")
655    dtype = property(_getdtype, "Data type of NeXus attribute")
656
657_npattrs = filter(lambda x: not x.startswith('_'), np.ndarray.__dict__.keys())
658
659class NXobject(object):
660
661    """
662    Abstract base class for elements in NeXus files.
663
664    The object has a subclass of NXfield, NXgroup, or one of the NXgroup
665    subclasses. Child nodes should be accessible directly as object attributes.
666    Constructors for NXobject objects are defined by either the NXfield or
667    NXgroup classes.
668
669    Python Attributes
670    -----------------
671    nxclass : string
672        The class of the NXobject. NXobjects can have class NXfield, NXgroup, or
673        be one of the NXgroup subclasses.
674    nxname : string
675        The name of the NXobject. Since it is possible to reference the same
676        Python object multiple times, this is not necessarily the same as the
677        object name. However, if the object is part of a NeXus tree, this will
678        be the attribute name within the tree.
679    nxgroup : NXgroup
680        The parent group containing this object within a NeXus tree. If the
681        object is not part of any NeXus tree, it will be set to None.
682    nxpath : string
683        The path to this object with respect to the root of the NeXus tree. For
684        NeXus data read from a file, this will be a group of class NXroot, but
685        if the NeXus tree was defined interactively, it can be any valid
686        NXgroup.
687    nxroot : NXgroup
688        The root object of the NeXus tree containing this object. For
689        NeXus data read from a file, this will be a group of class NXroot, but
690        if the NeXus tree was defined interactively, it can be any valid
691        NXgroup.
692    nxfile : NeXusTree
693        The file handle of the root object of the NeXus tree containing this
694        object.
695    filename : string
696        The file name of NeXus object's tree file handle.
697    attrs : dict
698        A dictionary of the NeXus object's attributes.
699
700    Methods
701    -------
702    dir(self, attrs=False, recursive=False):
703        Print the group directory.
704
705        The directory is a list of NeXus objects within this group, either NeXus
706        groups or NXfield data. If 'attrs' is True, NXfield attributes are
707        displayed. If 'recursive' is True, the contents of child groups are also
708        displayed.
709
710    tree:
711        Return the object's tree as a string.
712
713        It invokes the 'dir' method with both 'attrs' and 'recursive'
714        set to True. Note that this is defined as a property attribute and
715        does not require parentheses.
716
717    save(self, filename, format='w5')
718        Save the NeXus group into a file
719
720        The object is wrapped in an NXroot group (with name 'root') and an
721        NXentry group (with name 'entry'), if necessary, in order to produce
722        a valid NeXus file.
723
724    """
725
726    _class = "unknown"
727    _name = "unknown"
728    _group = None
729    _file = None
730    _infile = False
731    _saved = False
732    _changed = True
733
734    def __str__(self):
735        return "%s:%s"%(self.nxclass,self.nxname)
736
737    def __repr__(self):
738        return "NXobject('%s','%s')"%(self.nxclass,self.nxname)
739
740    def _setattrs(self, attrs):
741        for k,v in attrs.items():
742            self._attrs[k] = v
743
744    def _str_name(self,indent=0):
745        if self.nxclass == 'NXfield':
746            return " "*indent+self.nxname
747        else:
748            return " "*indent+self.nxname+':'+self.nxclass
749
750    def _str_value(self,indent=0):
751        return ""
752
753    def _str_attrs(self,indent=0):
754        names = self.attrs.keys()
755        names.sort()
756        result = []
757        for k in names:
758            result.append(" "*indent+"@%s = %s"%(k,self.attrs[k].nxdata))
759        return "\n".join(result)
760
761    def _str_tree(self,indent=0,attrs=False,recursive=False):
762        """
763        Print current object and possibly children.
764        """
765        result = [self._str_name(indent=indent)]
766        if attrs and self.attrs:
767            result.append(self._str_attrs(indent=indent+2))
768        # Print children
769        entries = self.entries
770        if entries:
771            names = entries.keys()
772            names.sort()
773            if recursive:
774                for k in names:
775                    result.append(entries[k]._str_tree(indent=indent+2,
776                                                       attrs=attrs, recursive=True))
777            else:
778                for k in names:
779                    result.append(entries[k]._str_name(indent=indent+2))
780        result
781        return "\n".join(result)
782
783    def walk(self):
784        if False: yield
785
786    def dir(self,attrs=False,recursive=False):
787        """
788        Print the object directory.
789
790        The directory is a list of NeXus objects within this object, either
791        NeXus groups or NXfields. If 'attrs' is True, NXfield attributes are
792        displayed. If 'recursive' is True, the contents of child groups are
793        also displayed.
794        """
795        print self._str_tree(attrs=attrs,recursive=recursive)
796
797    @property
798    def tree(self):
799        """
800        Return the directory tree as a string.
801
802        The tree contains all child objects of this object and their children.
803        It invokes the 'dir' method with both 'attrs' and 'recursive' set
804        to True.
805        """
806        return self._str_tree(attrs=True,recursive=True)
807
808    def __enter__(self):
809        """
810        Open the datapath for reading or writing.
811
812        Note: the results are undefined if you try accessing
813        more than one slab at a time.  Don't nest your
814        "with data" statements!
815        """
816        self._close_on_exit = not self.nxfile.isopen
817        self.nxfile.open() # Force file open even if closed
818        self.nxfile.openpath(self.nxpath)
819        self._incontext = True
820        return self.nxfile
821
822    def __exit__(self, type, value, traceback):
823        """
824        Close the file associated with the data.
825        """
826        self._incontext = False
827        if self._close_on_exit:
828            self.nxfile.close()
829
830    def save(self, filename=None, format='w5'):
831        """
832        Save the NeXus object to a data file.
833
834        An error is raised if the object is an NXroot group from an external file
835        that has been opened as readonly and no file name is specified.
836
837        The object is wrapped in an NXroot group (with name 'root') and an
838        NXentry group (with name 'entry'), if necessary, in order to produce
839        a valid NeXus file.
840       
841        Example
842        -------
843        >>> data = NXdata(sin(x), x)
844        >>> data.save('file.nxs')
845        >>> print data.nxroot.tree
846        root:NXroot
847          @HDF5_Version = 1.8.2
848          @NeXus_version = 4.2.1
849          @file_name = file.nxs
850          @file_time = 2012-01-20T13:14:49-06:00
851          entry:NXentry
852            data:NXdata
853              axis1 = float64(101)
854              signal = float64(101)
855                @axes = axis1
856                @signal = 1             
857        >>> root.entry.data.axis1.units = 'meV'
858        >>> root.save()
859        """
860        if filename:
861            if self.nxclass == "NXroot":
862                root = self
863            elif self.nxclass == "NXentry":
864                root = NXroot(self)
865            else:
866                root = NXroot(NXentry(self))
867            if root.nxfile: root.nxfile.close()
868            file = NeXusTree(filename, format)
869            file.writefile(root)
870            file.close()
871            root._file = NeXusTree(filename, 'rw')
872            root._setattrs(root._file.getattrs())
873            for node in root.walk():
874                node._infile = node._saved = True
875           
876        elif self.nxfile:
877            for entry in self.nxroot.values():
878                entry.write()
879
880        else:
881            raise NeXusError("No output file specified")
882
883    @property
884    def infile(self):
885        """
886        Returns True if the object has been created in a NeXus file.
887        """
888        return self._infile
889   
890    @property
891    def saved(self):
892        """
893        Returns True if the object has been saved to a file.
894        """
895        return self._saved
896
897    @property
898    def changed(self):
899        """
900        Returns True if the object has been changed.
901       
902        This property is for use by external scripts that need to track
903        which NeXus objects have been changed.
904        """
905        return self._changed
906   
907    def set_unchanged(self, recursive=False):
908        """
909        Set an object's change status to unchanged.
910        """
911        if recursive:
912            for node in self.walk():
913                node._changed = False
914        else:
915            self._changed = False
916   
917    def _getclass(self):
918        return self._class
919
920    def _getname(self):
921        return self._name
922
923    def _setname(self, value):
924        self._name = str(value)
925
926    def _getgroup(self):
927        return self._group
928
929    def _getpath(self):
930        if self.nxgroup is None:
931            return ""
932        elif isinstance(self.nxgroup, NXroot):
933            return "/" + self.nxname
934        else:
935            return self.nxgroup._getpath()+"/"+self.nxname
936
937    def _getroot(self):
938        if self.nxgroup is None:
939            return self
940        elif isinstance(self.nxgroup, NXroot):
941            return self.nxgroup
942        else:
943            return self.nxgroup._getroot()
944
945    def _getfile(self):
946        return self.nxroot._file
947
948    def _getfilename(self):
949        return self.nxroot._file.filename
950
951    def _getattrs(self):
952        return self._attrs
953
954    nxclass = property(_getclass, doc="Class of NeXus object")
955    nxname = property(_getname, _setname, doc="Name of NeXus object")
956    nxgroup = property(_getgroup, doc="Parent group of NeXus object")
957    nxpath = property(_getpath, doc="Path to NeXus object")
958    nxroot = property(_getroot, doc="Root group of NeXus object's tree")
959    nxfile = property(_getfile, doc="File handle of NeXus object's tree")
960    attrs = property(_getattrs, doc="NeXus attributes for an object")
961
962
963class NXfield(NXobject):
964
965    """
966    A NeXus data field.
967
968    This is a subclass of NXobject that contains scalar, array, or string data
969    and associated NeXus attributes.
970
971    NXfield(value=None, name='unknown', dtype='', shape=[], attrs={}, file=None,
972            path=None, group=None, **attr)
973
974    Input Parameters
975    ----------------
976    value : scalar value, Numpy array, or string
977        The numerical or string value of the NXfield, which is directly
978        accessible as the NXfield attribute 'nxdata'.
979    name : string
980        The name of the NXfield, which is directly accessible as the NXfield
981        attribute 'name'. If the NXfield is initialized as the attribute of a
982        parent object, the name is automatically set to the name of this
983        attribute.
984    dtype : string
985        The data type of the NXfield value, which is directly accessible as the
986        NXfield attribute 'dtype'. Valid input types correspond to standard
987        Numpy data types, using names defined by the NeXus API,
988        i.e., 'float32' 'float64'
989              'int8' 'int16' 'int32' 'int64'
990              'uint8' 'uint16' 'uint32' 'uint64'
991              'char'
992        If the data type is not specified, then it is determined automatically
993        by the data type of the 'value' parameter.
994    shape : list of ints
995        The dimensions of the NXfield data, which is accessible as the NXfield
996        attribute 'shape'. This corresponds to the shape of the Numpy array.
997        Scalars (numeric or string) are stored as Numpy zero-rank arrays,
998        for which shape=[].
999    attrs : dict
1000        A dictionary containing NXfield attributes. The dictionary values should
1001        all have class NXattr.
1002    file : filename
1003        The file from which the NXfield has been read.
1004    path : string
1005        The path to this object with respect to the root of the NeXus tree,
1006        using the convention for unix file paths.
1007    group : NXgroup or subclass of NXgroup
1008        The parent NeXus object. If the NXfield is initialized as the attribute
1009        of a parent group, this attribute is automatically set to the parent group.
1010
1011    Python Attributes
1012    -----------------
1013    nxclass : 'NXfield'
1014        The class of the NXobject.
1015    nxname : string
1016        The name of the NXfield. Since it is possible to reference the same
1017        Python object multiple times, this is not necessarily the same as the
1018        object name. However, if the field is part of a NeXus tree, this will
1019        be the attribute name within the tree.
1020    nxgroup : NXgroup
1021        The parent group containing this field within a NeXus tree. If the
1022        field is not part of any NeXus tree, it will be set to None.
1023    dtype : string or Numpy dtype
1024        The data type of the NXfield value. If the NXfield has been initialized
1025        but the data values have not been read in or defined, this is a string.
1026        Otherwise, it is set to the equivalent Numpy dtype.
1027    shape : list or tuple of ints
1028        The dimensions of the NXfield data. If the NXfield has been initialized
1029        but the data values have not been read in or defined, this is a list of
1030        ints. Otherwise, it is set to the equivalent Numpy shape, which is a
1031        tuple. Scalars (numeric or string) are stored as Numpy zero-rank arrays,
1032        for which shape=().
1033    attrs : dict
1034        A dictionary of all the NeXus attributes associated with the field.
1035        These are objects with class NXattr.
1036    nxdata : scalar, Numpy array or string
1037        The data value of the NXfield. This is normally initialized using the
1038        'value' parameter (see above). If the NeXus data is contained
1039        in a file and the size of the NXfield array is too large to be stored
1040        in memory, the value is not read in until this attribute is directly
1041        accessed. Even then, if there is insufficient memory, a value of None
1042        will be returned. In this case, the NXfield array should be read as a
1043        series of smaller slabs using 'get'.
1044    nxdata_as('units') : scalar value or Numpy array
1045        If the NXfield 'units' attribute has been set, the data values, stored
1046        in 'nxdata', are returned after conversion to the specified units.
1047    nxpath : string
1048        The path to this object with respect to the root of the NeXus tree. For
1049        NeXus data read from a file, this will be a group of class NXroot, but
1050        if the NeXus tree was defined interactively, it can be any valid
1051        NXgroup.
1052    nxroot : NXgroup
1053        The root object of the NeXus tree containing this object. For
1054        NeXus data read from a file, this will be a group of class NXroot, but
1055        if the NeXus tree was defined interactively, it can be any valid
1056        NXgroup.
1057
1058    NeXus Attributes
1059    ----------------
1060    NeXus attributes are stored in the 'attrs' dictionary of the NXfield, but
1061    can usually be assigned or referenced as if they are Python attributes, as
1062    long as the attribute name is not the same as one of those listed above.
1063    This is to simplify typing in an interactive session and should not cause
1064    any problems because there is no name clash with attributes so far defined
1065    within the NeXus standard. When writing modules, it is recommended that the
1066    attributes always be referenced using the 'attrs' dictionary if there is
1067    any doubt.
1068
1069    a) Assigning a NeXus attribute
1070
1071    In the example below, after assigning the NXfield, the following three
1072    NeXus attribute assignments are all equivalent:
1073
1074        >>> entry.sample.temperature = NXfield(40.0)
1075        >>> entry.sample.temperature.attrs['units'] = 'K'
1076        >>> entry.sample.temperature.units = NXattr('K')
1077        >>> entry.sample.temperature.units = 'K'
1078
1079    b) Referencing a NeXus attribute
1080
1081    If the name of the NeXus attribute is not the same as any of the Python
1082    attributes listed above, or one of the methods listed below, or any of the
1083    attributes defined for Numpy arrays, they can be referenced as if they were
1084    a Python attribute of the NXfield. However, it is only possible to reference
1085    attributes with one of the proscribed names using the 'attrs' dictionary.
1086
1087        >>> entry.sample.temperature.tree = 10.0
1088        >>> entry.sample.temperature.tree
1089        temperature = 40.0
1090          @tree = 10.0
1091          @units = K
1092        >>> entry.sample.temperature.attrs['tree']
1093        NXattr(10.0)
1094
1095    Numerical Operations on NXfields
1096    --------------------------------
1097    NXfields usually consist of arrays of numeric data with associated
1098    meta-data, the NeXus attributes. The exception is when they contain
1099    character strings. This makes them similar to Numpy arrays, and this module
1100    allows the use of NXfields in numerical operations in the same way as Numpy
1101    ndarrays. NXfields are technically not a sub-class of the ndarray class, but
1102    most Numpy operations work on NXfields, returning either another NXfield or,
1103    in some cases, an ndarray that can easily be converted to an NXfield.
1104
1105        >>> x = NXfield((1.0,2.0,3.0,4.0))
1106        >>> print x+1
1107        [ 2.  3.  4.  5.]
1108        >>> print 2*x
1109        [ 2.  4.  6.  8.]
1110        >>> print x/2
1111        [ 0.5  1.   1.5  2. ]
1112        >>> print x**2
1113        [  1.   4.   9.  16.]
1114        >>> print x.reshape((2,2))
1115        [[ 1.  2.]
1116         [ 3.  4.]]
1117        >>> y = NXfield((0.5,1.5,2.5,3.5))
1118        >>> x+y
1119        NXfield(name=x,value=[ 1.5  3.5  5.5  7.5])
1120        >>> x*y
1121        NXfield(name=x,value=[  0.5   3.    7.5  14. ])
1122        >>> (x+y).shape
1123        (4,)
1124        >>> (x+y).dtype
1125        dtype('float64')
1126
1127    All these operations return valid NXfield objects containing the same
1128    attributes as the first NXobject in the expression. The 'reshape' and
1129    'transpose' methods also return NXfield objects.
1130
1131    It is possible to use the standard slice syntax.
1132
1133        >>> x=NXfield(np.linspace(0,10,11))
1134        >>> x
1135        NXfield([  0.   1.   2. ...,   8.   9.  10.])
1136        >>> x[2:5]
1137        NXfield([ 2.  3.  4.])
1138
1139    In addition, it is possible to use floating point numbers as the slice
1140    indices. If one of the indices is not integer, both indices are used to
1141    extract elements in the array with values between the two index values.
1142
1143        >>> x=NXfield(np.linspace(0,100.,11))
1144        >>> x
1145        NXfield([   0.   10.   20. ...,   80.   90.  100.])
1146        >>> x[20.:50.]
1147        NXfield([ 20.  30.  40.  50.])
1148
1149    The standard Numpy ndarray attributes and methods will also work with
1150    NXfields, but will return scalars or Numpy arrays.
1151
1152        >>> x.size
1153        4
1154        >>> x.sum()
1155        10.0
1156        >>> x.max()
1157        4.0
1158        >>> x.mean()
1159        2.5
1160        >>> x.var()
1161        1.25
1162        >>> x.reshape((2,2)).sum(1)
1163        array([ 3.,  7.])
1164
1165    Finally, NXfields are cast as ndarrays for operations that require them.
1166    The returned value will be the same as for the equivalent ndarray
1167    operation, e.g.,
1168
1169    >>> np.sin(x)
1170    array([ 0.84147098,  0.90929743,  0.14112001, -0.7568025 ])
1171    >>> np.sqrt(x)
1172    array([ 1.        ,  1.41421356,  1.73205081,  2.        ])
1173
1174    Methods
1175    -------
1176    dir(self, attrs=False):
1177        Print the NXfield specification.
1178
1179        This outputs the name, dimensions and data type of the NXfield.
1180        If 'attrs' is True, NXfield attributes are displayed.
1181
1182    tree:
1183        Return the NXfield's tree.
1184
1185        It invokes the 'dir' method with both 'attrs' and 'recursive'
1186        set to True. Note that this is defined as a property attribute and
1187        does not require parentheses.
1188
1189
1190    save(self, filename, format='w5')
1191        Save the NXfield into a file wrapped in a NXroot group and NXentry group
1192        with default names. This is equivalent to
1193
1194        >>> NXroot(NXentry(NXfield(...))).save(filename)
1195
1196    Examples
1197    --------
1198    >>> x = NXfield(np.linspace(0,2*np.pi,101), units='degree')
1199    >>> phi = x.nxdata_as(units='radian')
1200    >>> y = NXfield(np.sin(phi))
1201
1202    # Read a Ni x Nj x Nk array one vector at a time
1203    >>> with root.NXentry[0].data.data as slab:
1204            Ni,Nj,Nk = slab.shape
1205            size = [1,1,Nk]
1206            for i in range(Ni):
1207                for j in range(Nj):
1208                    value = slab.get([i,j,0],size)
1209
1210    """
1211
1212    def __init__(self, value=None, name='field', dtype=None, shape=(), group=None,
1213                 attrs={}, **attr):
1214        if isinstance(value, list) or isinstance(value, tuple):
1215            value = np.array(value)
1216        self._value = value
1217        self._class = 'NXfield'
1218        self._name = name.replace(' ','_')
1219        self._group = group
1220        self._dtype = dtype
1221        if dtype == 'char':
1222            self._dtype = 'char'
1223        elif dtype in np.typeDict:
1224            self._dtype = np.dtype(dtype)
1225        elif dtype:
1226            raise NeXusError("Invalid data type: %s" % dtype)
1227        self._shape = tuple(shape)
1228        # Append extra keywords to the attribute list
1229        self._attrs = AttrDict()
1230        for key in attr.keys():
1231            attrs[key] = attr[key]
1232        # Convert NeXus attributes to python attributes
1233        self._setattrs(attrs)
1234        if 'units' in attrs:
1235            units = attrs['units']
1236        else:
1237            units = None
1238        self._incontext = False
1239        del attrs
1240        if value is not None and dtype == 'char': value = str(value)
1241        self._setdata(value)
1242        self._saved = False
1243        self._changed = True
1244
1245    def __repr__(self):
1246        if self._value is not None:
1247            if str(self.dtype) == 'char':
1248                return "NXfield('%s')" % str(self)
1249            else:
1250                return "NXfield(%s)" % self._str_value()
1251        else:
1252            return "NXfield(dtype=%s,shape=%s)" % (self.dtype,self.shape)
1253
1254    def __getattr__(self, name):
1255        """
1256        Enable standard numpy ndarray attributes if not otherwise defined.
1257        """
1258        if name in _npattrs:
1259            return self.nxdata.__getattribute__(name)
1260        elif name in self.attrs:
1261            return self.attrs[name].nxdata
1262        raise KeyError(name+" not in "+self.nxname)
1263
1264    def __setattr__(self, name, value):
1265        """
1266        Add an attribute to the NXfield 'attrs' dictionary unless the attribute
1267        name starts with 'nx' or '_', or unless it is one of the standard Python
1268        attributes for the NXfield class.
1269        """
1270        if name.startswith('_') or name.startswith('nx'):
1271            object.__setattr__(self, name, value)
1272        elif isinstance(value, NXattr):
1273            self._attrs[name] = value
1274            self._saved = False
1275            self._changed = True
1276        else:
1277            self._attrs[name] = NXattr(value)
1278            self._saved = False
1279            self._changed = True
1280
1281    def __getitem__(self, index):
1282        """
1283        Returns a slice from the NXfield.
1284
1285        In most cases, the slice values are applied to the NXfield nxdata array
1286        and returned within an NXfield object with the same metadata. However,
1287        if the array is one-dimensional and the index start and stop values
1288        are real, the nxdata array is returned with values between those limits.
1289        This is to allow axis arrays to be limited by their actual value. This
1290        real-space slicing should only be used on monotonically increasing (or
1291        decreasing) one-dimensional arrays.
1292        """
1293        if isinstance(index, slice) and \
1294           (isinstance(index.start, float) or isinstance(index.stop, float)):
1295            index = slice(self.index(index.start), self.index(index.stop,max=True)+1)
1296        if self._value is not None:
1297            result = self.nxdata.__getitem__(index)
1298        else:
1299            offset = np.zeros(len(self.shape),dtype=int)
1300            size = np.array(self.shape)
1301            if isinstance(index, int):
1302                offset[0] = index
1303                size[0] = 1
1304            else:
1305                if isinstance(index, slice): index = [index]
1306                i = 0
1307                for ind in index:
1308                    if isinstance(ind, int):
1309                        offset[i] = ind
1310                        size[i] = 1
1311                    else:
1312                        if ind.start: offset[i] = ind.start
1313                        if ind.stop: size[i] = ind.stop - offset[i]
1314                    i = i + 1
1315            try:
1316                result = self.get(offset, size)
1317            except ValueError:
1318                result = self.nxdata.__getitem__(index)
1319        return NXfield(result, name=self.nxname, attrs=self.attrs)
1320
1321    def __setitem__(self, index, value):
1322        """
1323        Assign a slice to the NXfield.
1324        """
1325        if self._value is not None:
1326            self.nxdata[index] = value
1327            self._saved = False
1328            self._changed = True
1329        else:
1330            raise NeXusError("NXfield dataspace not yet allocated")
1331
1332    def __deepcopy__(self, memo):
1333        dpcpy = self.__class__()
1334        memo[id(self)] = dpcpy
1335        dpcpy._value = copy(self._value)
1336        dpcpy._name = copy(self.nxname)
1337        dpcpy._dtype = copy(self.dtype)
1338        dpcpy._shape = copy(self.shape)
1339        for k, v in self.attrs.items():
1340            dpcpy.attrs[k] = copy(v)
1341        return dpcpy
1342
1343    def __len__(self):
1344        """
1345        Return the length of the NXfield data.
1346        """
1347        return np.prod(self.shape)
1348
1349    def index(self, value, max=False):
1350        """
1351        Return the index of the NXfield nxdata array that is greater than or equal to the value.
1352
1353        If max, then return the index that is less than or equal to the value.
1354        This should only be used on one-dimensional monotonically increasing arrays.
1355        """
1356        if max:
1357            return len(self.nxdata)-len(self.nxdata[self.nxdata>=value])
1358        else:
1359            return len(self.nxdata[self.nxdata<value])
1360
1361    def __array__(self):
1362        """
1363        Cast the NXfield as an array when it is expected by numpy
1364        """
1365        return self.nxdata
1366
1367    def __eq__(self, other):
1368        """
1369        Return true if the values of the NXfield are the same.
1370        """
1371        if isinstance(other, NXfield):
1372            if isinstance(self.nxdata, np.ndarray) and isinstance(other.nxdata, np.ndarray):
1373                return all(self.nxdata == other.nxdata)
1374            else:
1375                return self.nxdata == other.nxdata
1376        else:
1377            return False
1378
1379    def __ne__(self, other):
1380        """
1381        Return true if the values of the NXfield are not the same.
1382        """
1383        if isinstance(other, NXfield):
1384            if isinstance(self.nxdata, np.ndarray) and isinstance(other.nxdata, np.ndarray):
1385                return any(self.nxdata != other.nxdata)
1386            else:
1387                return self.nxdata != other.nxdata
1388        else:
1389            return True
1390
1391    def __add__(self, other):
1392        """
1393        Return the sum of the NXfield and another NXfield or number.
1394        """
1395        if isinstance(other, NXfield):
1396            return NXfield(value=self.nxdata+other.nxdata, name=self.nxname,
1397                           attrs=self.attrs)
1398        else:
1399            return NXfield(value=self.nxdata+other, name=self.nxname,
1400                           attrs=self.attrs)
1401
1402    def __radd__(self, other):
1403        """
1404        Return the sum of the NXfield and another NXfield or number.
1405
1406        This variant makes __add__ commutative.
1407        """
1408        return self.__add__(other)
1409
1410    def __sub__(self, other):
1411        """
1412        Return the NXfield with the subtraction of another NXfield or number.
1413        """
1414        if isinstance(other, NXfield):
1415            return NXfield(value=self.nxdata-other.nxdata, name=self.nxname,
1416                           attrs=self.attrs)
1417        else:
1418            return NXfield(value=self.nxdata-other, name=self.nxname,
1419                           attrs=self.attrs)
1420
1421    def __mul__(self, other):
1422        """
1423        Return the product of the NXfield and another NXfield or number.
1424        """
1425        if isinstance(other, NXfield):
1426            return NXfield(value=self.nxdata*other.nxdata, name=self.nxname,
1427                           attrs=self.attrs)
1428        else:
1429            return NXfield(value=self.nxdata*other, name=self.nxname,
1430                          attrs=self.attrs)
1431
1432    def __rmul__(self, other):
1433        """
1434        Return the product of the NXfield and another NXfield or number.
1435
1436        This variant makes __mul__ commutative.
1437        """
1438        return self.__mul__(other)
1439
1440    def __div__(self, other):
1441        """
1442        Return the NXfield divided by another NXfield or number.
1443        """
1444        if isinstance(other, NXfield):
1445            return NXfield(value=self.nxdata/other.nxdata, name=self.nxname,
1446                           attrs=self.attrs)
1447        else:
1448            return NXfield(value=self.nxdata/other, name=self.nxname,
1449                           attrs=self.attrs)
1450
1451    def __rdiv__(self, other):
1452        """
1453        Return the inverse of the NXfield divided by another NXfield or number.
1454        """
1455        if isinstance(other, NXfield):
1456            return NXfield(value=other.nxdata/self.nxdata, name=self.nxname,
1457                           attrs=self.attrs)
1458        else:
1459            return NXfield(value=other/self.nxdata, name=self.nxname,
1460                           attrs=self.attrs)
1461
1462    def __pow__(self, power):
1463        """
1464        Return the NXfield raised to the specified power.
1465        """
1466        return NXfield(value=pow(self.nxdata,power), name=self.nxname,
1467                       attrs=self.attrs)
1468
1469    def reshape(self, shape):
1470        """
1471        Returns an NXfield with the specified shape.
1472        """
1473        return NXfield(value=self.nxdata.reshape(shape), name=self.nxname,
1474                       attrs=self.attrs)
1475
1476    def transpose(self):
1477        """
1478        Returns an NXfield containing the transpose of the data array.
1479        """
1480        return NXfield(value=self.nxdata.transpose(), name=self.nxname,
1481                       attrs=self.attrs)
1482
1483    @property
1484    def T(self):
1485        return self.transpose()
1486
1487    def centers(self):
1488        """
1489        Returns an NXfield with the centers of a single axis
1490        assuming it contains bin boundaries.
1491        """
1492        return NXfield((self.nxdata[:-1]+self.nxdata[1:])/2,
1493                        name=self.nxname,attrs=self.attrs)
1494
1495    def read(self):
1496        """
1497        Read the NXfield, including attributes, from the NeXus file.
1498
1499        The data values are read provided they do not exceed NX_MEMORY. In that
1500        case, the data have to be read in as slabs using the get method.
1501        """
1502        if self.nxfile:
1503            with self as path:
1504                self._setattrs(path.getattrs())
1505                shape, dtype = path.getinfo()
1506                if dtype == 'char':
1507                    self._value = path.getdata()
1508                elif np.prod(shape) * np.dtype(dtype).itemsize <= NX_MEMORY*1024*1024:
1509                    self._value = path.getdata()
1510                else:
1511                    raise MemoryError('Data size larger than NX_MEMORY=%s MB' % NX_MEMORY)
1512                self._shape = tuple(shape)
1513                self._dtype = dtype
1514                if dtype == 'char':
1515                    self._dtype = 'char'
1516                elif dtype in np.typeDict:
1517                    self._dtype = np.dtype(dtype)
1518                self._infile = self._saved = self._changed = True
1519        else:
1520            raise IOError("Data is not attached to a file")
1521
1522    def write(self):
1523        """
1524        Write the NXfield, including attributes, to the NeXus file.
1525        """
1526        if self.nxfile:
1527            if self.nxfile.mode == napi.ACC_READ:
1528                raise NeXusError("NeXus file is readonly")
1529            if not self.infile:
1530                shape = self.shape
1531                if shape == (): shape = (1,)
1532                with self.nxgroup as path:
1533                    if np.prod(shape) > 10000:
1534                    # Compress the fastest moving dimension of large datasets
1535                        slab_dims = np.ones(len(shape),'i')
1536                        if shape[-1] < 100000:
1537                            slab_dims[-1] = shape[-1]
1538                        else:
1539                            slab_dims[-1] = 100000
1540                        path.compmakedata(self.nxname, self.dtype, shape, 'lzw', 
1541                                          slab_dims)
1542                    else:
1543                    # Don't use compression for small datasets
1544                        path.makedata(self.nxname, self.dtype, shape)
1545                self._infile = True
1546            if not self.saved:           
1547                with self as path:
1548                    path._writeattrs(self.attrs)
1549                    value = self.nxdata
1550                    if value is not None:
1551                        path.putdata(value)
1552                self._saved = True
1553        else:
1554            raise IOError("Data is not attached to a file")
1555
1556    def get(self, offset, size):
1557        """
1558        Return a slab from the data array.
1559
1560        Offsets are 0-origin. Shape can be inferred from the data.
1561        Offset and shape must each have one entry per dimension.
1562
1563        Corresponds to NXgetslab(handle,data,offset,shape)
1564        """
1565        if self.nxfile:
1566            with self as path:
1567                value = path.getslab(offset,size)
1568                return value
1569        else:
1570            raise IOError("Data is not attached to a file")
1571
1572    def put(self, data, offset, refresh=True):
1573        """
1574        Put a slab into the data array.
1575
1576        Offsets are 0-origin.  Shape can be inferred from the data.
1577        Offset and shape must each have one entry per dimension.
1578
1579        Corresponds to NXputslab(handle,data,offset,shape)
1580        """
1581        if self.nxfile:
1582            if self.nxfile.mode == napi.ACC_READ:
1583                raise NeXusError("NeXus file is readonly")
1584            with self as path:
1585                if isinstance(data, NXfield):
1586                    path.putslab(data.nxdata.astype(self.dtype), offset, data.shape)
1587                else:
1588                    data = np.array(data)
1589                    path.putslab(data.astype(self.dtype), offset, data.shape)
1590            if refresh: self.read()
1591        else:
1592            raise IOError("Data is not attached to a file")
1593
1594    def add(self, data, offset, refresh=True):
1595        """
1596        Add a slab into the data array.
1597
1598        Calls get to read in existing data before adding the value
1599        and calling put. It assumes that the two sets of data have
1600        compatible data types.
1601        """
1602        if isinstance(data, NXfield):
1603            value = self.get(offset, data.shape)
1604            self.put(data.nxdata.astype(self.dtype)+value, offset)
1605        else:
1606            value = self.get(offset, data.shape)
1607            self.put(data.astype(self.dtype)+value, offset)
1608        if refresh: self.refresh()
1609
1610    def refresh(self):
1611        """
1612        Rereads the data from the file.
1613
1614        If put has been called, then nxdata is no longer synchronized with the
1615        file making a refresh necessary. This will only be performed if nxdata
1616        already stores the data.
1617        """
1618        if self._value is not None:
1619            if self.nxfile:
1620                self._value = self.nxfile.readpath(self.nxpath)
1621                self._infile = self._saved = True
1622            else:
1623                raise IOError("Data is not attached to a file")
1624
1625    def convert(self, units=""):
1626        """
1627        Return the data in particular units.
1628        """
1629        try:
1630            import units
1631        except ImportError:
1632            raise NeXusError("No conversion utility available")
1633        if self._value is not None:
1634            return self._converter(self._value,units)
1635        else:
1636            return None
1637
1638    def __str__(self):
1639        """
1640        If value is loaded, return the value as a string.  If value is
1641        not loaded, return the empty string.  Only the first view values
1642        for large arrays will be printed.
1643        """
1644        if self._value is not None:
1645            return str(self._value)
1646        return ""
1647
1648    def _str_value(self,indent=0):
1649        v = str(self)
1650        if '\n' in v:
1651            v = '\n'.join([(" "*indent)+s for s in v.split('\n')])
1652        return v
1653
1654    def _str_tree(self,indent=0,attrs=False,recursive=False):
1655        dims = 'x'.join([str(n) for n in self.shape])
1656        s = str(self)
1657        if '\n' in s or s == "":
1658            s = "%s(%s)"%(self.dtype, dims)
1659        v=[" "*indent + "%s = %s"%(self.nxname, s)]
1660        if attrs and self.attrs: v.append(self._str_attrs(indent=indent+2))
1661        return "\n".join(v)
1662
1663    def walk(self):
1664        yield self
1665
1666    def _getaxes(self):
1667        """
1668        Return a list of NXfields containing axes.
1669
1670        Only works if the NXfield has the 'axes' attribute
1671        """
1672        try:
1673            return [getattr(self.nxgroup,name) for name in _readaxes(self.axes)]
1674        except KeyError:
1675            return None
1676
1677    def _getdata(self):
1678        """
1679        Return the data if it is not larger than NX_MEMORY.
1680        """
1681        if self._value is None:
1682            if self.nxfile:
1683                if str(self.dtype) == 'char':
1684                    self._value = self.nxfile.readpath(self.nxpath)
1685                elif np.prod(self.shape) * np.dtype(self.dtype).itemsize <= NX_MEMORY*1024*1024:
1686                    self._value = self.nxfile.readpath(self.nxpath)
1687                else:
1688                    raise MemoryError('Data size larger than NX_MEMORY=%s MB' % NX_MEMORY)
1689                self._saved = True
1690            else:
1691                return None
1692
1693        return self._value
1694
1695    def _setdata(self, value):
1696        if value is not None:
1697            if str(self._dtype) == 'char' or isinstance(value,str):
1698                self._value = str(value)
1699                self._shape = (len(self._value),)
1700                self._dtype = 'char'
1701            else:
1702                if self.dtype in np.typeDict:
1703                    self._value = np.array(value,self.dtype)
1704                else:
1705                    self._value = np.array(value)
1706                self._shape = self._value.shape
1707                self._dtype = self._value.dtype
1708            self._saved = False
1709            self._changed = True
1710       
1711    def _getdtype(self):
1712        return self._dtype
1713
1714    def _getshape(self):
1715        return self._shape
1716
1717    def _getsize(self):
1718        return len(self)
1719
1720    nxdata = property(_getdata,_setdata,doc="The data values")
1721    nxaxes = property(_getaxes,doc="The plotting axes")
1722    dtype = property(_getdtype,doc="Data type of NeXus field")
1723    shape = property(_getshape,doc="Shape of NeXus field")
1724    size = property(_getsize,doc="Size of NeXus field")
1725
1726SDS = NXfield # For backward compatibility
1727
1728def _fixaxes(signal, axes):
1729    """
1730    Remove length-one dimensions from plottable data
1731    """
1732    shape = list(signal.shape)
1733    while 1 in shape: shape.remove(1)
1734    newaxes = []
1735    for axis in axes:
1736        if axis.size > 1: newaxes.append(axis)
1737    return signal.nxdata.view().reshape(shape), newaxes
1738
1739class PylabPlotter(object):
1740
1741    """
1742    Matplotlib plotter class for NeXus data.
1743    """
1744
1745    def plot(self, signal, axes, title, errors, fmt, 
1746             xmin, xmax, ymin, ymax, zmin, zmax, **opts):
1747        """
1748        Plot the data entry.
1749
1750        Raises NeXusError if the data cannot be plotted.
1751        """
1752        try:
1753            import matplotlib.pyplot as plt
1754        except ImportError:
1755            raise NeXusError("Default plotting package (matplotlib) not available.")
1756
1757        over = False
1758        if "over" in opts.keys():
1759            if opts["over"]: over = True
1760            del opts["over"]
1761
1762        log = logx = logy = False
1763        if "log" in opts.keys():
1764            if opts["log"]: log = True
1765            del opts["log"]
1766        if "logy" in opts.keys():
1767            if opts["logy"]: logy = True
1768            del opts["logy"]
1769        if "logx" in opts.keys():
1770            if opts["logx"]: logx = True
1771            del opts["logx"]
1772
1773        if over:
1774            plt.autoscale(enable=False)
1775        else:
1776            plt.autoscale(enable=True)
1777            plt.clf()
1778
1779        # Provide a new view of the data if there is a dimension of length 1
1780        if 1 in signal.shape:
1781            data, axes = _fixaxes(signal, axes)
1782        else:
1783            data = signal.nxdata
1784
1785        # Find the centers of the bins for histogrammed data
1786        axis_data = centers(data, axes)
1787
1788        #One-dimensional Plot
1789        if len(data.shape) == 1:
1790            plt.ioff()
1791            if hasattr(signal, 'units'):
1792                if not errors and signal.units == 'counts':
1793                    errors = NXfield(np.sqrt(data))
1794            if errors:
1795                ebars = errors.nxdata
1796                plt.errorbar(axis_data[0], data, ebars, fmt=fmt, **opts)
1797            else:
1798                plt.plot(axis_data[0], data, fmt, **opts)
1799            if not over:
1800                ax = plt.gca()
1801                xlo, xhi = ax.set_xlim(auto=True)       
1802                ylo, yhi = ax.set_ylim(auto=True)               
1803                if xmin: xlo = xmin
1804                if xmax: xhi = xmax
1805                ax.set_xlim(xlo, xhi)
1806                if ymin: ylo = ymin
1807                if ymax: yhi = ymax
1808                ax.set_ylim(ylo, yhi)
1809                if logx: ax.set_xscale('symlog')
1810                if log or logy: ax.set_yscale('symlog')
1811                plt.xlabel(label(axes[0]))
1812                plt.ylabel(label(signal))
1813                plt.title(title)
1814            plt.ion()
1815
1816        #Two dimensional plot
1817        else:
1818            from matplotlib.image import NonUniformImage
1819            from matplotlib.colors import LogNorm
1820
1821            if len(data.shape) > 2:
1822                slab = [slice(None), slice(None)]
1823                for _dim in data.shape[2:]:
1824                    slab.append(0)
1825                data = data[slab].view().reshape(data.shape[:2])
1826                print "Warning: Only the top 2D slice of the data is plotted"
1827
1828            x = axis_data[0]
1829            y = axis_data[1]
1830            if not zmin: zmin = np.min(data)
1831            if not zmax: zmax = np.max(data)
1832            z = np.clip(data,zmin,zmax).T
1833           
1834            if log:
1835                opts["norm"] = LogNorm()
1836                if z.min() <= 0 and np.issubdtype(z[0,0],int):
1837                    z = np.clip(z,0.1,zmax)
1838
1839            ax = plt.gca()
1840            extent = (x[0],x[-1],y[0],y[-1])
1841            im = NonUniformImage(ax, extent=extent, origin=None, **opts)
1842            im.set_data(x,y,z)
1843            ax.images.append(im)
1844            xlo, xhi = ax.set_xlim(x[0],x[-1])
1845            ylo, yhi = ax.set_ylim(y[0],y[-1])
1846            if xmin: 
1847                xlo = xmin
1848            else:
1849                xlo = x[0]
1850            if xmax: 
1851                xhi = xmax
1852            else:
1853                xhi = x[-1]
1854            if ymin: 
1855                yhi = ymin
1856            else:
1857                yhi = y[0]
1858            if ymax: 
1859                yhi = ymax
1860            else:
1861                yhi = y[-1]
1862            ax.set_xlim(xlo, xhi)
1863            ax.set_ylim(ylo, yhi)
1864            plt.xlabel(label(axes[0]))
1865            plt.ylabel(label(axes[1]))
1866            plt.title(title)
1867            plt.colorbar(im)
1868
1869        plt.gcf().canvas.draw_idle()
1870
1871    @staticmethod
1872    def show():
1873        import matplotlib.pyplot as plt
1874        plt.show()   
1875
1876
1877class NXgroup(NXobject):
1878
1879    """
1880    A NeXus group object.
1881
1882    This is a subclass of NXobject and is the base class for the specific
1883    NeXus group classes, e.g., NXentry, NXsample, NXdata.
1884
1885    NXgroup(*items, **opts)
1886
1887    Parameters
1888    ----------
1889    The NXgroup parameters consist of a list of positional and/or keyword
1890    arguments.
1891
1892    Positional Arguments : These must be valid NeXus objects, either an NXfield
1893    or a NeXus group. These are added without modification as children of this
1894    group.
1895
1896    Keyword Arguments : Apart from a list of special keywords shown below,
1897    keyword arguments are used to add children to the group using the keywords
1898    as attribute names. The values can either be valid NXfields or NXgroups,
1899    in which case the 'name' attribute is changed to the keyword, or they
1900    can be numerical or string data, which are converted to NXfield objects.
1901
1902    Special Keyword Arguments:
1903
1904    name : string
1905        The name of the NXgroup, which is directly accessible as the NXgroup
1906        attribute 'name'. If the NXgroup is initialized as the attribute of
1907        a parent group, the name is automatically set to the name of this
1908        attribute. If 'nxclass' is specified and has the usual prefix 'NX',
1909        the default name is the class name without this prefix.
1910    nxclass : string
1911        The class of the NXgroup.
1912    entries : dict
1913        A dictionary containing a list of group entries. This is an
1914        alternative way of adding group entries to the use of keyword
1915        arguments.
1916    file : filename
1917        The file from which the NXfield has been read.
1918    path : string
1919        The path to this object with respect to the root of the NeXus tree,
1920        using the convention for unix file paths.
1921    group : NXobject (NXgroup or subclass of NXgroup)
1922        The parent NeXus group, which is accessible as the group attribute
1923        'group'. If the group is initialized as the attribute of
1924        a parent group, this is set to the parent group.
1925
1926    Python Attributes
1927    -----------------
1928    nxclass : string
1929        The class of the NXobject.
1930    nxname : string
1931        The name of the NXfield.
1932    entries : dictionary
1933        A dictionary of all the NeXus objects contained within an NXgroup.
1934    attrs : dictionary
1935        A dictionary of all the NeXus attributes, i.e., attribute with class NXattr.
1936    entries : dictionary
1937        A dictionary of all the NeXus objects contained within the group.
1938    attrs : dictionary
1939        A dictionary of all the group's NeXus attributes, which all have the
1940        class NXattr.
1941    nxpath : string
1942        The path to this object with respect to the root of the NeXus tree. For
1943        NeXus data read from a file, this will be a group of class NXroot, but
1944        if the NeXus tree was defined interactively, it can be any valid
1945        NXgroup.
1946    nxroot : NXgroup
1947        The root object of the NeXus tree containing this object. For
1948        NeXus data read from a file, this will be a group of class NXroot, but
1949        if the NeXus tree was defined interactively, it can be any valid
1950        NXgroup.
1951
1952    NeXus Group Entries
1953    -------------------
1954    Just as in a NeXus file, NeXus groups can contain either data or other
1955    groups, represented by NXfield and NXgroup objects respectively. To
1956    distinguish them from regular Python attributes, all NeXus objects are
1957    stored in the 'entries' dictionary of the NXgroup. However, they can usually
1958    be assigned or referenced as if they are Python attributes, i.e., using the
1959    dictionary name directly as the group attribute name, as long as this name
1960    is not the same as one of the Python attributes defined above or as one of
1961    the NXfield Python attributes.
1962
1963    a) Assigning a NeXus object to a NeXus group
1964
1965    In the example below, after assigning the NXgroup, the following three
1966    NeXus object assignments to entry.sample are all equivalent:
1967
1968        >>> entry.sample = NXsample()
1969        >>> entry.sample['temperature'] = NXfield(40.0)
1970        >>> entry.sample.temperature = NXfield(40.0)
1971        >>> entry.sample.temperature = 40.0
1972        >>> entry.sample.temperature
1973        NXfield(40.0)
1974
1975    If the assigned value is not a valid NXobject, then it is cast as an NXfield
1976    with a type determined from the Python data type.
1977
1978        >>> entry.sample.temperature = 40.0
1979        >>> entry.sample.temperature
1980        NXfield(40.0)
1981        >>> entry.data.data.x=np.linspace(0,10,11).astype('float32')
1982        >>> entry.data.data.x
1983        NXfield([  0.   1.   2. ...,   8.   9.  10.])
1984
1985    b) Referencing a NeXus object in a NeXus group
1986
1987    If the name of the NeXus object is not the same as any of the Python
1988    attributes listed above, or the methods listed below, they can be referenced
1989    as if they were a Python attribute of the NXgroup. However, it is only possible
1990    to reference attributes with one of the proscribed names using the group
1991    dictionary, i.e.,
1992
1993        >>> entry.sample.tree = 100.0
1994        >>> print entry.sample.tree
1995        sample:NXsample
1996          tree = 100.0
1997        >>> entry.sample['tree']
1998        NXfield(100.0)
1999
2000    For this reason, it is recommended to use the group dictionary to reference
2001    all group objects within Python scripts.
2002
2003    NeXus Attributes
2004    ----------------
2005    NeXus attributes are not currently used much with NXgroups, except for the
2006    root group, which has a number of global attributes to store the file name,
2007    file creation time, and NeXus and HDF version numbers. However, the
2008    mechanism described for NXfields works here as well. All NeXus attributes
2009    are stored in the 'attrs' dictionary of the NXgroup, but can be referenced
2010    as if they are Python attributes as long as there is no name clash.
2011
2012        >>> entry.sample.temperature = 40.0
2013        >>> entry.sample.attrs['tree'] = 10.0
2014        >>> print entry.sample.tree
2015        sample:NXsample
2016          @tree = 10.0
2017          temperature = 40.0
2018        >>> entry.sample.attrs['tree']
2019        NXattr(10.0)
2020
2021    Methods
2022    -------
2023    insert(self, NXobject, name='unknown'):
2024        Insert a valid NXobject (NXfield or NXgroup) into the group.
2025
2026        If NXobject has a 'name' attribute and the 'name' keyword is not given,
2027        then the object is inserted with the NXobject name.
2028
2029    makelink(self, NXobject):
2030        Add the NXobject to the group entries as a link (NXlink).
2031
2032    dir(self, attrs=False, recursive=False):
2033        Print the group directory.
2034
2035        The directory is a list of NeXus objects within this group, either NeXus
2036        groups or NXfield data. If 'attrs' is True, NXfield attributes are
2037        displayed. If 'recursive' is True, the contents of child groups are also
2038        displayed.
2039
2040    tree:
2041        Return the group tree.
2042
2043        It invokes the 'dir' method with both 'attrs' and 'recursive'
2044        set to True.
2045
2046    save(self, filename, format='w5')
2047        Save the NeXus group into a file
2048
2049        The object is wrapped in an NXroot group (with name 'root') and an
2050        NXentry group (with name 'entry'), if necessary, in order to produce
2051        a valid NeXus file.
2052
2053    Examples
2054    --------
2055    >>> x = NXfield(np.linspace(0,2*np.pi,101), units='degree')
2056    >>> entry = NXgroup(x, name='entry', nxclass='NXentry')
2057    >>> entry.sample = NXgroup(temperature=NXfield(40.0,units='K'),
2058                               nxclass='NXsample')
2059    >>> print entry.sample.tree
2060    sample:NXsample
2061      temperature = 40.0
2062        @units = K
2063
2064    Note: All the currently defined NeXus classes are defined as subclasses
2065          of the NXgroup class. It is recommended that these are used
2066          directly, so that the above examples become:
2067
2068    >>> entry = NXentry(x)
2069    >>> entry.sample = NXsample(temperature=NXfield(40.0,units='K'))
2070
2071    or
2072
2073    >>> entry.sample.temperature = 40.0
2074    >>> entry.sample.temperature.units='K'
2075
2076    """
2077
2078    # Plotter to use for plot calls
2079    _plotter = PylabPlotter()
2080
2081    def __init__(self, *items, **opts):
2082        if "name" in opts.keys():
2083            self._name = opts["name"].replace(' ','_')
2084            del opts["name"]
2085        self._entries = {}
2086        if "entries" in opts.keys():
2087            for k,v in opts["entries"].items():
2088                setattr(self, k, v)
2089            del opts["entries"]
2090        self._attrs = AttrDict()
2091        if "attrs" in opts.keys():
2092            self._setattrs(opts["attrs"])
2093            del opts["attrs"]
2094        if "nxclass" in opts.keys():
2095            self._class = opts["nxclass"]
2096            del opts["nxclass"]
2097        if "group" in opts.keys():
2098            self._group = opts["group"]
2099            del opts["group"]
2100        for k,v in opts.items():
2101            setattr(self, k, v)
2102        if self.nxclass.startswith("NX"):
2103            if self.nxname == "unknown": self._name = self.nxclass[2:]
2104            try: # If one exists, set the class to a valid NXgroup subclass
2105                self.__class__ = globals()[self.nxclass]
2106            except KeyError:
2107                pass
2108        for item in items:
2109            try:
2110                setattr(self, item.nxname, item)
2111            except AttributeError:
2112                raise NeXusError("Non-keyword arguments must be valid NXobjects")
2113        self._saved = False
2114        self._changed = True
2115
2116#    def __cmp__(self, other):
2117#        """Sort groups by their distances or names."""
2118#        try:
2119#            return cmp(self.distance, other.distance)
2120#        except KeyError:
2121#            return cmp(self.nxname, other.nxname)
2122
2123    def __repr__(self):
2124        return "%s('%s')" % (self.__class__.__name__,self.nxname)
2125
2126    def _str_value(self,indent=0):
2127        return ""
2128
2129    def walk(self):
2130        yield self
2131        for node in self.entries.values():
2132            for child in node.walk():
2133                yield child
2134
2135    def __getattr__(self, key):
2136        """
2137        Provide direct access to groups via nxclass name.
2138        """
2139        if key.startswith('NX'):
2140            return self.component(key)
2141        elif key in self.entries:
2142            return self.entries[key]
2143        elif key in self.attrs:
2144            return self.attrs[key].nxdata
2145        raise KeyError(key+" not in "+self.nxclass+":"+self.nxname)
2146
2147    def __setattr__(self, name, value):
2148        """
2149        Set an attribute as an object or regular Python attribute.
2150
2151        It is assumed that attributes starting with 'nx' or '_' are regular
2152        Python attributes. All other attributes are converted to valid NXobjects,
2153        with class NXfield, NXgroup, or a sub-class of NXgroup, depending on the
2154        assigned value.
2155
2156        The internal value of the attribute name, i.e., 'name', is set to the
2157        attribute name used in the assignment.  The parent group of the
2158        attribute, i.e., 'group', is set to the parent group of the attribute.
2159
2160        If the assigned value is a numerical (scalar or array) or string object,
2161        it is converted to an object of class NXfield, whose attribute, 'nxdata',
2162        is set to the assigned value.
2163        """
2164        if name.startswith('_') or name.startswith('nx'):
2165            object.__setattr__(self, name, value)
2166        elif isinstance(value, NXattr):
2167            self._attrs[name] = value
2168            self._saved = False
2169            self._changed = True
2170        else:
2171            self[name] = value
2172
2173    def __getitem__(self, index):
2174        """
2175        Returns a slice from the NXgroup nxsignal attribute (if it exists) as
2176        a new NXdata group, if the index is a slice object.
2177
2178        In most cases, the slice values are applied to the NXfield nxdata array
2179        and returned within an NXfield object with the same metadata. However,
2180        if the array is one-dimensional and the index start and stop values
2181        are real, the nxdata array is returned with values between the limits
2182        set by those axis values.
2183        This is to allow axis arrays to be limited by their actual value. This
2184        real-space slicing should only be used on monotonically increasing (or
2185        decreasing) one-dimensional arrays.
2186        """
2187        if isinstance(index, str): #i.e., requesting a dictionary value
2188            return self._entries[index]
2189
2190        if not self.nxsignal:
2191            raise NeXusError("No plottable signal")
2192        if not hasattr(self,"nxclass"):
2193            raise NeXusError("Indexing not allowed for groups of unknown class")
2194        if isinstance(index, int):
2195            axes = self.nxaxes
2196            axes[0] = axes[0][index]
2197            result = NXdata(self.nxsignal[index], axes)
2198            if self.nxerrors: result.errors = self.errors[index]
2199        elif isinstance(index, slice):
2200            axes = self.nxaxes
2201            axes[0] = axes[0][index]
2202            if isinstance(index.start, float) or isinstance(index.stop, float):
2203                index = slice(self.nxaxes[0].index(index.start),
2204                              self.nxaxes[0].index(index.stop,max=True)+1)
2205                result = NXdata(self.nxsignal[index], axes)
2206                if self.nxerrors: result.errors = self.errors[index]
2207            else:
2208                result = NXdata(self.nxsignal[index], axes)
2209                if self.nxerrors: result.errors = self.errors[index]
2210        else:
2211            i = 0
2212            slices = []
2213            axes = self.nxaxes
2214            for ind in index:
2215                axes[i] = axes[i][ind]
2216                if isinstance(ind, slice) and \
2217                   (isinstance(ind.start, float) or isinstance(ind.stop, float)):
2218                    slices.append(slice(self.nxaxes[i].index(ind.start),
2219                                        self.nxaxes[i].index(ind.stop)))
2220                else:
2221                    slices.append(ind)
2222                i = i + 1
2223            result = NXdata(self.nxsignal.__getitem__(tuple(slices)), axes)
2224            if self.nxerrors: result.errors = self.errors.__getitem__(tuple(slices))
2225        axes = []
2226        for axis in result.nxaxes:
2227            if len(axis) > 1: axes.append(axis)
2228        result.nxsignal.axes = ":".join([axis.nxname for axis in axes])
2229        if self.nxtitle:
2230            result.title = self.nxtitle
2231        return result
2232
2233    def __setitem__(self, key, value):
2234        """
2235        Adds or modifies an item in the NeXus group.
2236        """
2237        if key in self.entries: 
2238            infile = self._entries[key]._infile
2239            if isinstance(self._entries[key], NXlink):
2240                if self._entries[key].nxlink:
2241                    setattr(self._entries[key].nxlink.nxgroup, key, value)
2242                return
2243            attrs = self._entries[key].attrs
2244        else:
2245            infile = None
2246            attrs = {}
2247        if isinstance(value, NXlink):
2248            self._entries[key] = value
2249        elif isinstance(value, NXobject):
2250            if value.nxgroup is not None:
2251                memo = {}
2252                value = deepcopy(value, memo)
2253                value._attrs = copy(value._attrs)
2254            value._group = self
2255            value._name = key
2256            self._entries[key] = value
2257        else:
2258            self._entries[key] = NXfield(value=value, name=key, group=self, attrs=attrs)
2259        if infile is not None: self[key]._infile = infile
2260        self._changed = True
2261
2262    def __deepcopy__(self, memo):
2263        dpcpy = self.__class__()
2264        memo[id(self)] = dpcpy
2265        for k,v in self.items():
2266            if isinstance(v, NXgroup):
2267                dpcpy[k] = deepcopy(v, memo)
2268            else:
2269                dpcpy[k] = copy(v)
2270        for k, v in self.attrs.items():
2271            dpcpy.attrs[k] = copy(v)
2272        return dpcpy
2273
2274    def keys(self):
2275        """
2276        Returns the names of NeXus objects in the group.
2277        """
2278        return self._entries.keys()
2279
2280    def values(self):
2281        """
2282        Returns the values of NeXus objects in the group.
2283        """
2284        return self._entries.values()
2285
2286    def items(self):
2287        """
2288        Returns a list of the NeXus objects in the group as (key,value) pairs.
2289        """
2290        return self._entries.items()
2291
2292    def has_key(self, name):
2293        """
2294        Returns true if the NeXus object with the specified name is in the group.
2295        """
2296        return self._entries.has_key(name)
2297
2298    def insert(self, value, name='unknown'):
2299        """
2300        Adds an attribute to the group.
2301
2302        If it is not a valid NeXus object (NXfield or NXgroup), the attribute
2303        is converted to an NXfield.
2304        """
2305        if isinstance(value, NXobject):
2306            if name == 'unknown': name = value.nxname
2307            if name in self._entries:
2308                raise NeXusError("'%s' already exists in group" % name)
2309            value._group = self
2310            self._entries[name] = value
2311        else:
2312            self._entries[name] = NXfield(value=value, name=name, group=self)
2313
2314    def makelink(self, target):
2315        """
2316        Creates a linked NXobject within the group.
2317
2318        All attributes are inherited from the parent object including the name
2319        """
2320        if isinstance(target, NXobject):
2321            self[target.nxname] = NXlink(target=target, group=self)
2322        else:
2323            raise NeXusError("Link target must be an NXobject")
2324
2325    def read(self):
2326        """
2327        Read the NXgroup and all its children from the NeXus file.
2328        """
2329        if self.nxfile:
2330            with self as path:
2331                n, nxname, nxclass = path.getgroupinfo()
2332                if nxclass != self.nxclass:
2333                    raise NeXusError("The NeXus group class does not match the file")
2334                self._setattrs(path.getattrs())
2335                entries = path.entries()
2336            for name,nxclass in entries:
2337                path = self.nxpath + '/' + name
2338                if nxclass == 'SDS':
2339                    attrs = self.nxfile.getattrs()
2340                    if 'target' in attrs and attrs['target'] != path:
2341                        self._entries[name] = NXlinkfield(target=attrs['target'])           
2342                    else:
2343                        self._entries[name] = NXfield(name=name)
2344                else:
2345                    attrs = self.nxfile.getattrs()
2346                    if 'target' in attrs and attrs['target'] != path:
2347                        self._entries[name] = NXlinkgroup(name=name,
2348                                                          target=attrs['target'])
2349                    else:
2350                        self._entries[name] = NXgroup(nxclass=nxclass)
2351                self._entries[name]._group = self
2352            #Make sure non-linked variables are processed first.
2353            for entry in self._entries.values():
2354                for node in entry.walk():
2355                    if not isinstance(node, NXlink): node.read()
2356            for entry in self._entries.values():
2357                for node in entry.walk():
2358                    if isinstance(node, NXlink): node.read()
2359            self._infile = self._saved = self._changed = True
2360        else:
2361            raise IOError("Data is not attached to a file")
2362
2363    def write(self):
2364        """
2365        Write the NXgroup, including its children, to the NeXus file.
2366        """
2367        if self.nxfile:
2368            if self.nxfile.mode == napi.ACC_READ:
2369                raise NeXusError("NeXus file is readonly")
2370            if not self.infile:
2371                with self.nxgroup as path:
2372                    path.makegroup(self.nxname, self.nxclass)
2373                self._infile = True
2374            with self as path:
2375                path._writeattrs(self.attrs)
2376                for entry in self.walk():
2377                    if entry is not self: entry.write()
2378                self._infile = self._saved = True
2379        else:
2380            raise IOError("Group is not attached to a file")
2381
2382    def sum(self, axis=None):
2383        """
2384        Return the sum of the NXdata group using the Numpy sum method
2385        on the NXdata signal.
2386
2387        The result contains a copy of all the metadata contained in
2388        the NXdata group.
2389        """
2390        if not self.nxsignal:
2391            raise NeXusError("No signal to sum")
2392        if not hasattr(self,"nxclass"):
2393            raise NeXusError("Summing not allowed for groups of unknown class")
2394        if axis is None:
2395            return self.nxsignal.sum()
2396        else:
2397            signal = NXfield(self.nxsignal.sum(axis), name=self.nxsignal.nxname)
2398            axes = self.nxaxes
2399            summedaxis = axes.pop(axis)
2400            units = ""
2401            if hasattr(summedaxis, "units"): units = summedaxis.units
2402            signal.long_name = "Integral from %s to %s %s" % \
2403                               (summedaxis[0], summedaxis[-1], units)
2404            average = NXfield(0.5*(summedaxis.nxdata[0]+summedaxis.nxdata[-1]), 
2405                                   name=summedaxis.nxname)
2406            if units: average.units = units
2407            result = NXdata(signal, axes, average)
2408            if self.nxerrors:
2409                errors = np.sqrt((self.nxerrors.nxdata**2).sum(axis))
2410                result.errors = NXfield(errors, name="errors")
2411            if self.nxtitle:
2412                result.title = self.nxtitle
2413            return result
2414
2415    def moment(self, order=1):
2416        """
2417        Return an NXfield containing the moments of the NXdata group
2418        assuming the signal is one-dimensional.
2419
2420        Currently, only the first moment has been defined. Eventually, the
2421        order of the moment will be defined by the 'order' parameter.
2422        """
2423        if not self.nxsignal:
2424            raise NeXusError("No signal to calculate")
2425        elif len(self.nxsignal.shape) > 1:
2426            raise NeXusError("Operation only possible on one-dimensional signals")
2427        elif order > 1:
2428            raise NeXusError("Higher moments not yet implemented")
2429        if not hasattr(self,"nxclass"):
2430            raise NeXusError("Operation not allowed for groups of unknown class")
2431        return (centers(self.nxsignal,self.nxaxes)*self.nxsignal).sum() \
2432                /self.nxsignal.sum()
2433
2434    def component(self, nxclass):
2435        """
2436        Find all child objects that have a particular class.
2437        """
2438        return [E for _name,E in self.entries.items() if E.nxclass==nxclass]
2439
2440    def signals(self):
2441        """
2442        Return a dictionary of NXfield's containing signal data.
2443
2444        The key is the value of the signal attribute.
2445        """
2446        signals = {}
2447        for obj in self.entries.values():
2448            if 'signal' in obj.attrs:
2449                signals[obj.nxsignal.nxdata] = obj
2450        return signals
2451
2452    def _signal(self):
2453        """
2454        Return the NXfield containing the signal data.
2455        """
2456        for obj in self.entries.values():
2457            if 'signal' in obj.attrs and str(obj.signal) == '1':
2458#                if isinstance(self[obj.nxname],NXlink):
2459#                    return self[obj.nxname].nxlink
2460#                else:
2461                return self[obj.nxname]
2462        return None
2463   
2464    def _set_signal(self, signal):
2465        """
2466        Setter for the signal attribute.
2467       
2468        The argument should be a valid NXfield within the group.
2469        """
2470        self[signal.nxname].signal = NXattr(1)
2471
2472    def _axes(self):
2473        """
2474        Return a list of NXfields containing the axes.
2475        """
2476        try:
2477            return [getattr(self,name) for name in _readaxes(self.nxsignal.axes)]
2478        except KeyError:
2479            axes = {}
2480            for obj in self.entries:
2481                if 'axis' in getattr(self,obj).attrs:
2482                    axes[getattr(self,obj).axis] = getattr(self,obj)
2483            return [axes[key] for key in sorted(axes.keys())]
2484
2485    def _set_axes(self, axes):
2486        """
2487        Setter for the signal attribute.
2488       
2489        The argument should be a list of valid NXfields within the group.
2490        """
2491        if not isinstance(axes, list):
2492            axes = [axes]
2493        self.nxsignal.axes = NXattr(":".join([axis.nxname for axis in axes]))
2494
2495    def _errors(self):
2496        """
2497        Return the NXfield containing the signal errors.
2498        """
2499        try:
2500            return self.entries['errors']
2501        except KeyError:
2502            return None
2503
2504    def _title(self):
2505        """
2506        Return the title as a string.
2507
2508        If there is no title attribute in the string, the parent
2509        NXentry group in the group's path is searched.
2510        """
2511        title = self.nxpath
2512        if 'title' in self.entries:
2513            return str(self.title)
2514        elif self.nxgroup:
2515            if 'title' in self.nxgroup.entries:
2516                return str(self.nxgroup.title)
2517        return self.nxpath
2518
2519    def _getentries(self):
2520        return self._entries
2521
2522    nxsignal = property(_signal, _set_signal, "Signal NXfield within group")
2523    nxaxes = property(_axes, _set_axes, "List of axes within group")
2524    nxerrors = property(_errors, "Errors NXfield within group")
2525    nxtitle = property(_title, "Title for group plot")
2526    entries = property(_getentries,doc="NeXus objects within group")
2527
2528    def plot(self, fmt='bo', xmin=None, xmax=None, ymin=None, ymax=None,
2529             zmin=None, zmax=None, **opts):
2530        """
2531        Plot data contained within the group.
2532
2533        The format argument is used to set the color and type of the
2534        markers or lines for one-dimensional plots, using the standard
2535        matplotlib syntax. The default is set to blue circles. All
2536        keyword arguments accepted by matplotlib.pyplot.plot can be
2537        used to customize the plot.
2538       
2539        In addition to the matplotlib keyword arguments, the following
2540        are defined:
2541       
2542            log = True     - plot the intensity on a log scale
2543            logy = True    - plot the y-axis on a log scale
2544            logx = True    - plot the x-axis on a log scale
2545            over = True    - plot on the current figure
2546
2547        Raises NeXusError if the data could not be plotted.
2548        """
2549
2550        group = self
2551        if self.nxclass == "NXroot":
2552            group = group.NXentry[0]
2553        if group.nxclass == "NXentry":
2554            try:
2555                group = group.NXdata[0]
2556            except IndexError:
2557                raise NeXusError('No NXdata group found')
2558
2559        # Find a plottable signal
2560        signal = group.nxsignal
2561        if not signal:
2562            raise NeXusError('No plottable signal defined')
2563
2564        # Find errors
2565        errors= group.nxerrors
2566
2567        # Find the associated axes
2568        axes = group.nxaxes
2569
2570        # Construct title
2571        title = group.nxtitle
2572
2573        # Plot with the available plotter
2574        group._plotter.plot(signal, axes, title, errors, fmt, 
2575                            xmin, xmax, ymin, ymax, zmin, zmax, **opts)
2576   
2577    def oplot(self, fmt='bo', **opts):
2578        """
2579        Plot the data contained within the group over the current figure.
2580        """
2581        self.plot(fmt=fmt, over=True, **opts)
2582
2583    def logplot(self, fmt='bo', xmin=None, xmax=None, ymin=None, ymax=None,
2584                zmin=None, zmax=None, **opts):
2585        """
2586        Plot the data intensity contained within the group on a log scale.
2587        """
2588        self.plot(fmt=fmt, log=True,
2589                  xmin=xmin, xmax=xmax, ymin=ymin, ymax=ymax,
2590                  zmin=zmin, zmax=zmax, **opts)
2591
2592class NXlink(NXobject):
2593
2594    """
2595    Class for NeXus linked objects.
2596
2597    The real object will be accessible by following the link attribute.
2598    """
2599
2600    _class = "NXlink"
2601
2602    def __init__(self, target=None, name='link', group=None):
2603        self._group = group
2604        self._class = "NXlink"
2605        if isinstance(target, NXobject):
2606            self._name = target.nxname
2607            self._target = target.nxpath
2608            self.nxlink.attrs["target"] = target.nxpath
2609            if target.nxclass == "NXlink":
2610                raise NeXusError("Cannot link to another NXlink object")
2611            elif target.nxclass == "NXfield":
2612                self.__class__ = NXlinkfield
2613            else:
2614                self.__class__ = NXlinkgroup
2615        else:
2616            self._name = name
2617            self._target = target
2618
2619    def __getattr__(self, key):
2620        try:
2621            try:
2622                return self.nxlink.__dict__[key]
2623            except KeyError:
2624                return self.nxlink.__getattr__(key)
2625        except KeyError:
2626            raise KeyError((key+" not in %s" % self._target))
2627
2628    def __setattr__(self, name, value):
2629        if name.startswith('_')  or name.startswith('nx'):
2630            object.__setattr__(self, name, value)
2631        elif self.nxlink:
2632            self.nxlink.__setattr__(name, value)
2633
2634    def __repr__(self):
2635        return "NXlink('%s')"%(self._target)
2636
2637    def __str__(self):
2638        return str(self.nxlink)
2639
2640    def _str_tree(self, indent=0, attrs=False, recursive=False):
2641        if self.nxlink:
2642            return self.nxlink._str_tree(indent, attrs, recursive)
2643        else:
2644            return " "*indent+self.nxname+' -> '+self._target
2645
2646    def _getlink(self):
2647        link = self.nxroot
2648        if link:
2649            try:
2650                for level in self._target[1:].split('/'):
2651                    link = link.entries[level]
2652                return link
2653            except AttributeError:
2654                return None
2655        else:
2656            return None
2657
2658    def _getattrs(self):
2659        return self.nxlink.attrs
2660
2661    nxlink = property(_getlink, "Linked object")
2662    attrs = property(_getattrs,doc="NeXus attributes for object")
2663
2664    def read(self):
2665        """
2666        Read the linked NXobject.
2667        """
2668        self.nxlink.read()
2669        self._infile = self._saved = self._changed = True
2670
2671
2672class NXlinkfield(NXlink, NXfield):
2673
2674    """
2675    Class for a NeXus linked field.
2676
2677    The real field will be accessible by following the link attribute.
2678    """
2679
2680    def write(self):
2681        """
2682        Write the linked NXfield.
2683        """
2684        self.nxlink.write()
2685        if not self.infile:
2686            with self.nxlink as path:
2687                target = path.getdataID()
2688            with self.nxgroup as path:
2689                path.makelink(target)
2690            self._infile = self._saved = True
2691
2692    def get(self, offset, size):
2693        """
2694        Get a slab from the data array.
2695
2696        Offsets are 0-origin.  Shape can be inferred from the data.
2697        Offset and shape must each have one entry per dimension.
2698
2699        This operation should be performed in a "with group.data"
2700        conext.
2701
2702        Raises ValueError cannot convert units.
2703
2704        Corresponds to NXgetslab(handle,data,offset,shape)
2705        """
2706        if self.nxfile:
2707            with self.nxlink as path:
2708                value = path.getslab(offset,size)
2709        else:
2710            raise IOError("Data is not attached to a file")
2711
2712NXlinkdata = NXlinkfield # For backward compatibility
2713
2714class NXlinkgroup(NXlink, NXgroup):
2715
2716    """
2717    Class for a NeXus linked group.
2718
2719    The real group will be accessible by following the link attribute.
2720    """
2721
2722    def write(self):
2723        """
2724        Write the linked NXgroup.
2725        """
2726        self.nxlink.write()
2727        if not self.infile:
2728            with self.nxlink as path:
2729                target = path.getgroupID()
2730            with self.nxgroup as path:
2731                path.makelink(target)
2732            self._infile = self._saved = True
2733
2734    def _getentries(self):
2735        return self.nxlink.entries
2736
2737    entries = property(_getentries,doc="NeXus objects within group")
2738
2739
2740class NXentry(NXgroup):
2741
2742    """
2743    NXentry group. This is a subclass of the NXgroup class.
2744
2745    Each NXdata and NXmonitor object of the same name will be added
2746    together, raising an NeXusError if any of the groups do not exist
2747    in both NXentry groups or if any of the NXdata additions fail.
2748    The resulting NXentry group contains a copy of all the other metadata
2749    contained in the first group. Note that other extensible data, such
2750    as the run duration, are not currently added together.
2751
2752    See the NXgroup documentation for more details.
2753    """
2754
2755    def __init__(self, *items, **opts):
2756        self._class = "NXentry"
2757        NXgroup.__init__(self, *items, **opts)
2758
2759    def __add__(self, other):
2760        """
2761        Add two NXentry objects
2762        """
2763        result = NXentry(entries=self.entries, attrs=self.attrs)
2764        try:
2765            names = [group.nxname for group in self.component("NXdata")]
2766            for name in names:
2767                if isinstance(other.entries[name], NXdata):
2768                    result.entries[name] = self.entries[name] + other.entries[name]
2769                else:
2770                    raise KeyError
2771            names = [group.nxname for group in self.component("NXmonitor")]
2772            for name in names:
2773                if isinstance(other.entries[name], NXmonitor):
2774                    result.entries[name] = self.entries[name] + other.entries[name]
2775                else:
2776                    raise KeyError
2777            return result
2778        except KeyError:
2779            raise NeXusError("Inconsistency between two NXentry groups")
2780
2781    def __sub__(self, other):
2782        """
2783        Subtract two NXentry objects
2784        """
2785        result = NXentry(entries=self.entries, attrs=self.attrs)
2786        try:
2787            names = [group.nxname for group in self.component("NXdata")]
2788            for name in names:
2789                if isinstance(other.entries[name], NXdata):
2790                    result.entries[name] = self.entries[name] - other.entries[name]
2791                else:
2792                    raise KeyError
2793            names = [group.nxname for group in self.component("NXmonitor")]
2794            for name in names:
2795                if isinstance(other.entries[name], NXmonitor):
2796                    result.entries[name] = self.entries[name] - other.entries[name]
2797                else:
2798                    raise KeyError
2799            return result
2800        except KeyError:
2801            raise NeXusError("Inconsistency between two NXentry groups")
2802
2803
2804class NXsubentry(NXentry):
2805
2806    """
2807    NXsubentry group. This is a subclass of the NXsubentry class.
2808
2809    See the NXgroup documentation for more details.
2810    """
2811
2812    def __init__(self, *items, **opts):
2813        self._class = "NXsubentry"
2814        NXgroup.__init__(self, *items, **opts)
2815
2816
2817class NXdata(NXgroup):
2818
2819    """
2820    NXdata group. This is a subclass of the NXgroup class.
2821
2822    The constructor assumes that the first argument contains the signal and
2823    the second contains either the axis, for one-dimensional data, or a list
2824    of axes, for multidimensional data. These arguments can either be NXfield
2825    objects or Numpy arrays, which are converted to NXfield objects with default
2826    names. Alternatively, the signal and axes NXfields can be defined using the
2827    'nxsignal' and 'nxaxes' properties. See the examples below.
2828   
2829    Various arithmetic operations (addition, subtraction, multiplication,
2830    and division) have been defined for combining NXdata groups with other
2831    NXdata groups, Numpy arrays, or constants, raising a NeXusError if the
2832    shapes don't match. Data errors are propagated in quadrature if
2833    they are defined, i.e., if the 'nexerrors' attribute is not None,
2834
2835    Attributes
2836    ----------
2837    nxsignal : The NXfield containing the attribute 'signal' with value 1
2838    nxaxes   : A list of NXfields containing the signal axes
2839    nxerrors : The NXfield containing the errors
2840
2841    Methods
2842    -------
2843    plot(self, fmt, over=False, log=False, logy=False, logx=False, **opts)
2844        Plot the NXdata group using the defined signal and axes. Valid
2845        Matplotlib parameters, specifying markers, colors, etc, can be
2846        specified using format argument or through keyword arguments.
2847
2848    logplot(self, fmt, over=False, logy=False, logx=False, **opts)
2849        Plot the NXdata group using the defined signal and axes with
2850        the intensity plotted on a logarithmic scale. In one-dimensional
2851        plots, this is the y-axis. In two-dimensional plots, it is the
2852        color scale.
2853
2854    oplot(self, fmt, **opts)
2855        Plot the NXdata group using the defined signal and axes over
2856        the current plot.
2857
2858    moment(self, order=1)
2859        Calculate moments of the NXdata group. This assumes that the
2860        signal is one-dimenional. Currently, only the first moment is
2861        implemented.
2862
2863    Examples
2864    --------
2865    There are three methods of creating valid NXdata groups with the
2866    signal and axes NXfields defined according to the NeXus standard.
2867   
2868    1) Create the NXdata group with Numpy arrays that will be assigned
2869       default names.
2870       
2871    >>> x = np.linspace(0, 2*np.pi, 101)
2872    >>> line = NXdata(sin(x), x)
2873    data:NXdata
2874      signal = float64(101)
2875        @axes = x
2876        @signal = 1
2877      axis1 = float64(101)
2878     
2879    2) Create the NXdata group with NXfields that have their internal
2880       names already assigned.
2881
2882    >>> x = NXfield(linspace(0,2*pi,101), name='x')
2883    >>> y = NXfield(linspace(0,2*pi,101), name='y')   
2884    >>> X, Y = np.meshgrid(x, y)
2885    >>> z = NXfield(sin(X) * sin(Y), name='z')
2886    >>> entry = NXentry()
2887    >>> entry.grid = NXdata(z, (x, y))
2888    >>> grid.tree()
2889    entry:NXentry
2890      grid:NXdata
2891        x = float64(101)
2892        y = float64(101)
2893        z = float64(101x101)
2894          @axes = x:y
2895          @signal = 1
2896
2897    3) Create the NXdata group with keyword arguments defining the names
2898       and set the signal and axes using the nxsignal and nxaxes properties.
2899
2900    >>> x = linspace(0,2*pi,101)
2901    >>> y = linspace(0,2*pi,101) 
2902    >>> X, Y = np.meshgrid(x, y)
2903    >>> z = sin(X) * sin(Y)
2904    >>> entry = NXentry()
2905    >>> entry.grid = NXdata(z=sin(X)*sin(Y), x=x, y=y)
2906    >>> entry.grid.nxsignal = entry.grid.z
2907    >>> entry.grid.nxaxes = [entry.grid.x.entry.grid.y]
2908    >>> grid.tree()
2909    entry:NXentry
2910      grid:NXdata
2911        x = float64(101)
2912        y = float64(101)
2913        z = float64(101x101)
2914          @axes = x:y
2915          @signal = 1
2916    """
2917
2918    def __init__(self, signal=None, axes=None, *items, **opts):
2919        self._class = "NXdata"
2920        NXgroup.__init__(self, *items, **opts)
2921        if signal is not None:
2922            if isinstance(signal,NXfield):
2923                if signal.nxname == "unknown": signal.nxname = "signal"
2924                if "signal" not in signal.attrs: signal.signal = 1
2925                self[signal.nxname] = signal
2926                signalname = signal.nxname
2927            else:
2928                self["signal"] = signal
2929                self["signal"].signal = 1
2930                signalname = "signal"
2931            if axes is not None:
2932                if not isinstance(axes,tuple) and not isinstance(axes,list):
2933                    axes = [axes]
2934                axisnames = {}
2935                i = 0
2936                for axis in axes:
2937                    i = i + 1
2938                    if isinstance(axis,NXfield):
2939                        if axis._name == "unknown": axis._name = "axis%s" % i
2940                        self[axis.nxname] = axis
2941                        axisnames[i] = axis.nxname
2942                    else:
2943                        axisname = "axis%s" % i
2944                        self[axisname] = axis
2945                        axisnames[i] = axisname
2946                self[signalname].axes = ":".join(axisnames.values())
2947
2948    def __add__(self, other):
2949        """
2950        Define a method for adding a NXdata group to another NXdata group
2951        or to a number. Only the signal data is affected.
2952
2953        The result contains a copy of all the metadata contained in
2954        the first NXdata group. The module checks that the dimensions are
2955        compatible, but does not check that the NXfield names or values are
2956        identical. This is so that spelling variations or rounding errors
2957        do not make the operation fail. However, it is up to the user to
2958        ensure that the results make sense.
2959        """
2960        result = NXdata(entries=self.entries, attrs=self.attrs)
2961        if isinstance(other, NXdata):
2962            if self.nxsignal and self.nxsignal.shape == other.nxsignal.shape:
2963                result.entries[self.nxsignal.nxname] = self.nxsignal + other.nxsignal
2964                if self.nxerrors:
2965                    if other.nxerrors:
2966                        result.errors = np.sqrt(self.errors.nxdata**2+other.errors.nxdata**2)
2967                    else:
2968                        result.errors = self.errors
2969                return result
2970        elif isinstance(other, NXgroup):
2971            raise NeXusError("Cannot add two arbitrary groups")
2972        else:
2973            result.entries[self.nxsignal.nxname] = self.nxsignal + other
2974            result.entries[self.nxsignal.nxname].nxname = self.nxsignal.nxname
2975            return result
2976
2977    def __sub__(self, other):
2978        """
2979        Define a method for subtracting a NXdata group or a number from
2980        the NXdata group. Only the signal data is affected.
2981
2982        The result contains a copy of all the metadata contained in
2983        the first NXdata group. The module checks that the dimensions are
2984        compatible, but does not check that the NXfield names or values are
2985        identical. This is so that spelling variations or rounding errors
2986        do not make the operation fail. However, it is up to the user to
2987        ensure that the results make sense.
2988        """
2989        result = NXdata(entries=self.entries, attrs=self.attrs)
2990        if isinstance(other, NXdata):
2991            if self.nxsignal and self.nxsignal.shape == other.nxsignal.shape:
2992                result.entries[self.nxsignal.nxname] = self.nxsignal - other.nxsignal
2993                if self.nxerrors:
2994                    if other.nxerrors:
2995                        result.errors = np.sqrt(self.errors.nxdata**2+other.errors.nxdata**2)
2996                    else:
2997                        result.errors = self.errors
2998                return result
2999        elif isinstance(other, NXgroup):
3000            raise NeXusError("Cannot subtract two arbitrary groups")
3001        else:
3002            result.entries[self.nxsignal.nxname] = self.nxsignal - other
3003            result.entries[self.nxsignal.nxname].nxname = self.nxsignal.nxname
3004            return result
3005
3006    def __mul__(self, other):
3007        """
3008        Define a method for multiplying the NXdata group with a NXdata group
3009        or a number. Only the signal data is affected.
3010
3011        The result contains a copy of all the metadata contained in
3012        the first NXdata group. The module checks that the dimensions are
3013        compatible, but does not check that the NXfield names or values are
3014        identical. This is so that spelling variations or rounding errors
3015        do not make the operation fail. However, it is up to the user to
3016        ensure that the results make sense.
3017        """
3018        result = NXdata(entries=self.entries, attrs=self.attrs)
3019        if isinstance(other, NXdata):
3020
3021            # error here signal not defined in this scope
3022            #if self.nxsignal and signal.shape == other.nxsignal.shape:
3023            if self.nxsignal and self.nxsignal.shape == other.nxsignal.shape:
3024                result.entries[self.nxsignal.nxname] = self.nxsignal * other.nxsignal
3025                if self.nxerrors:
3026                    if other.nxerrors:
3027                        result.errors = np.sqrt((self.errors.nxdata*other.nxsignal.nxdata)**2+
3028                                                (other.errors.nxdata*self.nxsignal.nxdata)**2)
3029                    else:
3030                        result.errors = self.errors
3031                return result
3032        elif isinstance(other, NXgroup):
3033            raise NeXusError("Cannot multiply two arbitrary groups")
3034        else:
3035            result.entries[self.nxsignal.nxname] = self.nxsignal * other
3036            result.entries[self.nxsignal.nxname].nxname = self.nxsignal.nxname
3037            if self.nxerrors:
3038                result.errors = self.errors * other
3039            return result
3040
3041    def __rmul__(self, other):
3042        """
3043        Define a method for multiplying NXdata groups.
3044
3045        This variant makes __mul__ commutative.
3046        """
3047        return self.__mul__(other)
3048
3049    def __div__(self, other):
3050        """
3051        Define a method for dividing the NXdata group by a NXdata group
3052        or a number. Only the signal data is affected.
3053
3054        The result contains a copy of all the metadata contained in
3055        the first NXdata group. The module checks that the dimensions are
3056        compatible, but does not check that the NXfield names or values are
3057        identical. This is so that spelling variations or rounding errors
3058        do not make the operation fail. However, it is up to the user to
3059        ensure that the results make sense.
3060        """
3061        result = NXdata(entries=self.entries, attrs=self.attrs)
3062        if isinstance(other, NXdata):
3063            if self.nxsignal and self.nxsignal.shape == other.nxsignal.shape:
3064                # error here, signal and othersignal not defined here
3065                #result.entries[self.nxsignal.nxname] = signal / othersignal
3066                result.entries[self.nxsignal.nxname] = self.nxsignal / other.nxsignal
3067                resultvalues = result.entries[self.nxsignal.nxname].nxdata
3068                if self.nxerrors:
3069                    if other.nxerrors:
3070                        result.errors = (np.sqrt(self.errors.nxdata**2 +
3071                                         (resultvalues*other.errors.nxdata)**2)
3072                                         / other.nxsignal)
3073                    else:
3074                        result.errors = self.errors
3075                return result
3076        elif isinstance(other, NXgroup):
3077            raise NeXusError("Cannot divide two arbitrary groups")
3078        else:
3079            result.entries[self.nxsignal.nxname] = self.nxsignal / other
3080            result.entries[self.nxsignal.nxname].nxname = self.nxsignal.nxname
3081            if self.nxerrors: result.errors = self.errors / other
3082            return result
3083
3084
3085class NXmonitor(NXdata):
3086
3087    """
3088    NXmonitor group. This is a subclass of the NXdata class.
3089
3090    See the NXdata and NXgroup documentation for more details.
3091    """
3092
3093    def __init__(self, signal=None, axes=(), *items, **opts):
3094        NXdata.__init__(self, signal=signal, axes=axes, *items, **opts)
3095        self._class = "NXmonitor"
3096        if "name" not in opts.keys():
3097            self._name = "monitor"
3098
3099
3100class NXlog(NXgroup):
3101
3102    """
3103    NXlog group. This is a subclass of the NXgroup class.
3104
3105    Methods
3106    -------
3107    plot(self, **opts)
3108        Plot the logged values against the elapsed time. Valid
3109        Matplotlib parameters, specifying markers, colors, etc, can be
3110        specified using the 'opts' dictionary.
3111
3112    See the NXgroup documentation for more details.
3113    """
3114
3115    def __init__(self, *items, **opts):
3116        self._class = "NXlog"
3117        NXgroup.__init__(self, *items, **opts)
3118
3119    def plot(self, **opts):
3120        axis = [self.time]
3121        title = "%s Log" % self.value.nxname.upper()
3122        self._plotter.plot(self.value, axis, title, **opts)
3123
3124
3125#-------------------------------------------------------------------------
3126#Add remaining base classes as subclasses of NXgroup and append to __all__
3127
3128for _class in _nxclasses:
3129    if _class not in globals():
3130        docstring = """
3131                    %s group. This is a subclass of the NXgroup class.
3132
3133                    See the NXgroup documentation for more details.
3134                    """ % _class
3135        globals()[_class]=type(_class, (NXgroup,),
3136                               {'_class':_class,'__doc__':docstring})
3137    __all__.append(_class)
3138
3139#-------------------------------------------------------------------------
3140
3141
3142def centers(signal, axes):
3143    """
3144    Return the centers of the axes.
3145
3146    This works regardless if the axes contain bin boundaries or centers.
3147    """
3148    def findc(axis, dimlen):
3149        if axis.shape[0] == dimlen+1:
3150            return (axis.nxdata[:-1] + axis.nxdata[1:])/2
3151        else:
3152            assert axis.shape[0] == dimlen
3153            return axis.nxdata
3154    return [findc(a,signal.shape[i]) for i,a in enumerate(axes)]
3155
3156def setmemory(value):
3157    """
3158    Set the memory limit for data arrays (in MB).
3159    """
3160    global NX_MEMORY
3161    NX_MEMORY = value
3162
3163def label(field):
3164    """
3165    Return a label for a data field suitable for use on a graph axis.
3166    """
3167    if hasattr(field,'long_name'):
3168        return field.long_name
3169    elif hasattr(field,'units'):
3170        return "%s (%s)"%(field.nxname,field.units)
3171    else:
3172        return field.nxname
3173
3174# File level operations
3175def load(filename, mode='r'):
3176    """
3177    Read a NeXus file returning a tree of objects.
3178
3179    This is aliased to 'read' because of potential name clashes with Numpy
3180    """
3181    file = NeXusTree(filename,mode)
3182    tree = file.readfile()
3183    file.close()
3184    return tree
3185
3186#Definition for when there are name clashes with Numpy
3187nxload = load
3188__all__.append('nxload')
3189
3190def save(filename, group, format='w5'):
3191    """
3192    Write a NeXus file from a tree of objects.
3193    """
3194    if group.nxclass == "NXroot":
3195        tree = group
3196    elif group.nxclass == "NXentry":
3197        tree = NXroot(group)
3198    else:
3199        tree = NXroot(NXentry(group))
3200    file = NeXusTree(filename, format)
3201    file.writefile(tree)
3202    file.close()
3203
3204def tree(file):
3205    """
3206    Read and summarize the named NeXus file.
3207    """
3208    nxfile = load(file)
3209    nxfile.tree
3210
3211def demo(argv):
3212    """
3213    Process a list of command line commands.
3214
3215    'argv' should contain program name, command, arguments, where command is one
3216    of the following:
3217        copy fromfile.nxs tofile.nxs
3218        ls f1.nxs f2.nxs ...
3219    """
3220    if len(argv) > 1:
3221        op = argv[1]
3222    else:
3223        op = 'help'
3224    if op == 'ls':
3225        for f in argv[2:]: dir(f)
3226    elif op == 'copy' and len(argv)==4:
3227        tree = load(argv[2])
3228        save(argv[3], tree)
3229    elif op == 'plot' and len(argv)==4:
3230        tree = load(argv[2])
3231        for entry in argv[3].split('.'):
3232            tree = getattr(tree,entry)
3233        tree.plot()
3234        tree._plotter.show()
3235
3236    else:
3237        usage = """
3238usage: %s cmd [args]
3239    copy fromfile.nxs tofile.nxs
3240    ls *.nxs
3241    plot file.nxs entry.data
3242        """%(argv[0],)
3243        print usage
3244
3245
3246if __name__ == "__main__":
3247    import sys
3248    demo(sys.argv)
3249
Note: See TracBrowser for help on using the repository browser.