Versions Compared

Key

  • This line was added.
  • This line was removed.
  • Formatting was changed.

Include Page
PageMenuBegin
PageMenuBegin
Table of Contents
Include Page
PageMenuEnd
PageMenuEnd

Introduction

This document offers a brief overview of a high-level interface to the psana analysis framework for user applications written in the Python programming language. The interface is informally known as "interactive psana*. The tool's use is not solely limited to the interactive analysis scenarios. It also allows its users to benefit from a rich set of services of the core framework while retaining a full control over an iteration in data sets (runs, files, etc.). This combination makes it possible for the interactive exploration and (if needed) visualization of the experimental data. Note that by the later we always mean data files in the XTC or HDF5 formats produced at the LCLS DAQ or Data Management system. We also suggest visiting the Glossary of Terms psana - Interactive API which is found in the end of the document.

...

The final comment, before we'll proceed to the practical steps, is that a reader of the document isn't required to be fully familiar with the batch framework. Those areas where such knowledge would be needed are expected to be covered by the document. Though, we still encourage our users to spend some time to get an overview of the Data Analysis Tools we provide at PCDS. That's because many problems in doing the data analysis can be solved by the batch version of psana in a more efficient and natural way. These two flavors of the framework are not meant to compete with each other, they are designed to complement each other to cover a broader spectrum of analysis scenarios.

Test data

Include Page
Tutorial test data
Tutorial test data

Setting up the analysis environment

...

  1. obtain and properly configure your UNIX account at PCDS. Specific instructions can be found in the Account Setup section of the Analysis Workbook.
  2. select the latest analysis release by running the following command before launching the Python interpreter:

    Code Block
    bgColor#F7F7ED
    
    % sit_setup ana-current
    

The later command should be done just once per session. It will initialize all relevant environment variables, including PATH, LD_LIBRARY_PATH, PYTHONPATH and some others. This will also give you an access to an appropriate version of the Python interpreter and the corresponding Python modules. We recommend using ipython. The following example illustrates how to launch the interpreter and test if the interactive psana module is available in the session environment. When everything is set up correctly one should see:

Code Block
bgColor#F7F7ED

% ipython

Python 2.7.2 (default, Jan 14 2013, 21:09:22)
Type "copyright", "credits" or "license" for more information.

IPython 0.13.1 -- An enhanced Interactive Python.
?         -> Introduction and overview of IPython's features.
%quickref -> Quick reference.
help      -> Python's own help system.
object?   -> Details about 'object', use 'object??' for extra details.

In [1]: import psana

In [2]: psana.
Display all 108 possibilities? (y or n)
psana.Acqiris            psana.Gsc16ai            psana.ndarray_float32_1  psana.ndarray_int32_5    psana.ndarray_uint32_3
psana.Andor              psana.Imp                psana.ndarray_float32_2  psana.ndarray_int32_6    psana.ndarray_uint32_4
psana.Bld                psana.Ipimb              psana.ndarray_float32_3  psana.ndarray_int64_1    psana.ndarray_uint32_5
psana.BldInfo            psana.Lusi               psana.ndarray_float32_4  psana.ndarray_int64_2    psana.ndarray_uint32_6
..

...

Here is the psana version of the traditional "Hello World!" program:

Code Block
bgColor#F7F7ED

from psana import *

dataset_name = "exp=CXI/cxitut13:run=22"
ds = DataSource(dataset_name)

for num,evt in enumerate(ds.events()):
    id = evt.get(EventId)
    print "Event #",num," has id:",id

The code of this example does nothing but scanning through all events of the data set and reporting identifiers (class EventId) of each event (class Event). The identifiers () encapsulate a number of attributes, including: timestamp, fiducials, etc. And here is how the output should look like:

Code Block
bgColor#F7F7ED

Event # 0  has id: XtcEventId(run=22, time=2013-01-18 17:58:53.047288760-08, fiducials=19404, ticks=331022, vector=1)
Event # 1  has id: XtcEventId(run=22, time=2013-01-18 17:58:53.055622526-08, fiducials=19407, ticks=330630, vector=2)
Event # 2  has id: XtcEventId(run=22, time=2013-01-18 17:58:53.063956294-08, fiducials=19410, ticks=329468, vector=3)
..

Now let's go through the example's code line-by-line and see what it does at each step:

  1. importing all definitions from the psanamodule into the global namespace. This includes functions, classes and types:

    Code Block
    bgColor#F7F7ED
    
    from psana import *
    
  2. opening a data set and checking if it exists/available to your process:

    Code Block
    bgColor#F7F7ED
    
    dataset_name = "exp=CXI/cxitut13:run=22"
    ds = DataSource(dataset_name)
    
  3. iterating over all events in the data set and obtaining an event identification object for each event:

    Code Block
    bgColor#F7F7ED
    
    for num,evt in enumerate(ds.events()):
        id = evt.get(EventId)
    

...

The data set string encodes various parameters, some of which are needed to locate data files, while others would affect the behavior of the file reader. The general syntax of the string is:

Code Block
bgColor#F7F7ED

par[=val][:par[=val][...]

These are some of the parameters which are supported in psana:

  • experiment name (which may optionally contain the name of an instrument)

    Code Block
    bgColor#F7F7ED
    
    exp=cxi12313
    exp=CXI/cxi12313
    
  • run number specification (can be a single run, a range of runs, a series of runs, or a combination of all above):

    Code Block
    bgColor#F7F7ED
    
    run=1
    run=10-20
    run=1,2,3,4
    run=1,20-20,31,41
    
  • file type (presently the default type is 'xtc')

    Code Block
    bgColor#F7F7ED
    
    xtc
    h5
    
  • a random stream (if a value is omitted) or a specific stream. Note this option only makes a sense for XTC files:

    Code Block
    bgColor#F7F7ED
    
    one-stream
    one-stream=2
    
  • allow reading from live files while they're still being recorded (by the DAQ or by the Data Migration service). Note that this feature is only available when running psanaat PCDS, in all other cases the option will be ignored:

    Code Block
    bgColor#F7F7ED
    
    live
    

Putting all together one would see a data set specification which would tell the framework to read data of stream #2 from XTC files of run 41 while these files were sill being recorded (by the DAQ, or the data migration service, or by another process of the same user, etc.):

Code Block
bgColor#F7F7ED

exp=CXI/cxi12313:run=41:xtc:one-stream=2:live

...

In this section we're going to focus on an event object to see how to get various information from it. Let's begin with an example where we're fetching and plotting an image captured at the Princeton camera (which is one of the detectors available at the XCS instrument). In this example we won't be iterating over all events. Only the first event will be considered:

Code Block
bgColor#F7F7ED

from psana import *

ds = DataSource('exp=XCS/xcstut13:run=15')
src = Source('DetInfo(XcsBeamline.0:Princeton.0)')

itr = ds.events()
evt = itr.next()
frame = evt.get(Princeton.FrameV1, src)

import matplotlib.pyplot as plt

plt.figure('Princeton Camera')
plt.imshow(frame.data())
plt.show()

...

The type key is a real Python type known to the framework. All types which would be recognized by a particular version of the framework can be obtained by the 'dot' operator of the ipython interpreter as shown below:

Code Block
bgColor#F7F7ED

% ipython

In [1]: import psana

In [2]: psana.
Display all 108 possibilities? (y or n)
psana.Acqiris            psana.Gsc16ai            psana.ndarray_float32_1  psana.ndarray_int32_5    psana.ndarray_uint32_3
psana.Andor              psana.Imp                psana.ndarray_float32_2  psana.ndarray_int32_6    psana.ndarray_uint32_4
psana.Bld                psana.Ipimb              psana.ndarray_float32_3  psana.ndarray_int64_1    psana.ndarray_uint32_5
psana.BldInfo            psana.Lusi               psana.ndarray_float32_4  psana.ndarray_int64_2    psana.ndarray_uint32_6
psana.Camera             psana.OceanOptics        psana.ndarray_float32_5  psana.ndarray_int64_3    psana.ndarray_uint64_1
psana.ControlData        psana.Opal1k             psana.ndarray_float32_6  psana.ndarray_int64_4    psana.ndarray_uint64_2
psana.CsPad              psana.Orca               psana.ndarray_float64_1  psana.ndarray_int64_5    psana.ndarray_uint64_3
psana.CsPad2x2           psana.PNCCD              psana.ndarray_float64_2  psana.ndarray_int64_6    psana.ndarray_uint64_4
psana.DataSource         psana.PSAna              psana.ndarray_float64_3  psana.ndarray_int8_1     psana.ndarray_uint64_5
..

Note that some of those entries aren't really types. Besides, they may be nested Python modules providing more types like this one:

Code Block
bgColor#F7F7ED

In [3]: psana.Princeton.
psana.Princeton.Config    psana.Princeton.ConfigV3  psana.Princeton.Frame     psana.Princeton.Info
psana.Princeton.ConfigV1  psana.Princeton.ConfigV4  psana.Princeton.FrameV1   psana.Princeton.InfoV1
psana.Princeton.ConfigV2  psana.Princeton.ConfigV5  psana.Princeton.FrameV2

...

The second (source) key should be constructed using a special class called Source which also exported by the psana module:

Code Block
bgColor#F7F7ED

src = Source('<object address string>')

...

The following example demonstrates how to dump a catalog of event components:

Code Block
bgColor#F7F7ED

In [8]: evt.keys()
Out[8]:
[EventKey(type=psana.EvrData.DataV3, src='DetInfo(NoDetector.0:Evr.0)'),
 EventKey(type=psana.Camera.FrameV1, src='DetInfo(CxiDg1.0:Tm6740.0)'),
 EventKey(type=psana.Camera.FrameV1, src='DetInfo(CxiDg2.0:Tm6740.0)'),
 EventKey(type=psana.CsPad.DataV1, src='DetInfo(CxiDs1.0:Cspad.0)'),
 EventKey(type=psana.Bld.BldDataEBeamV3, src='BldInfo(EBeam)'),
 EventKey(type=psana.Bld.BldDataFEEGasDetEnergy, src='BldInfo(FEEGasDetEnergy)'),
 EventKey(type=psana.EventId),
 EventKey(type=None)]

In [9]:

The information reported by the method would give us an idea what's found in the event, and how to obtain those components using the get() method. When psana has been imported into the global namespace, the mapping is:

Code Block
bgColor#F7F7ED

from psana import *
...
obj1 = evt.get( EvrData.DataV3,             Source('DetInfo(NoDetector.0:Evr.0)'))
obj2 = evt.get( Camera.FrameV1,             Source('DetInfo(CxiDg1.0:Tm6740.0)'))
obj3 = evt.get( Camera.FrameV1,             Source('DetInfo(CxiDg2.0:Tm6740.0)'))
obj4 = evt.get( CsPad.DataV1,               Source('DetInfo(CxiDs1.0:Cspad.0)'))
obj5 = evt.get( Bld.BldDataEBeamV3,         Source('BldInfo(EBeam)'))
obj6 = evt.get( Bld.BldDataFEEGasDetEnergy, Source('BldInfo(FEEGasDetEnergy)'))
obj7 = evt.get( EventId)

...

For those would like to build some automation in discovering which components and of what kind exist in the event there is another option. A user can iterate over the list of key elements to examine their attributes. Each such element would encapsulate a type, a source and a key (string) of the corresponding event component. Consider the following example:

Code Block
bgColor#F7F7ED

for k in evt.keys():
    print "type:   ", k.type()
    print "source: ", k.src()
    print "key:    ", k.key()
    ..

The type() method will return one of those Python type objects which were already mentioned in section "Four forms of the get() method". If a user is looking for a key which has a specific type then the following type comparison can be used:

Code Block
bgColor#F7F7ED

for k in evt.keys():
    if k.type() == Camera.FrameV1:
        print "got Camera.FrameV1"

...

Since APIs of those specific source classes differ one from another then a user would need to obtain the final type using the following technique:

Code Block
bgColor#F7F7ED

for k in evt.keys():
    src = k.src()
    src_type = type(src)
    if   src_type == DetInfo  : print "detector name: ", src.detName(), " device name:   ", src.devName(), ...
    elif src_type == BldInfo  : print "detector type: ", src.type(),    " detector name: ", src.detName(), ...
    elif src_type == ProcInfo : print "IP address:    ", src.ipAddr(),  " process id:    ", src.processId()

...

And the last method of the event API will return a run number. This information may be useful for data sets spanning across many runs:

Code Block
bgColor#F7F7ED

ds = DataSource("exp=CXI/cxitut13:run=22,23,24,25)
for evt in ds.events():
    print "run: ", evt.run()
    ..

...

The event environment encapsulate a broad spectrum of data and services which have various origins. This environment is needed to evaluate or process event data in a proper context. Some of these data may have different life cycles than events. Other parts of this information (such as calibrations) may not even come directly from the input data stream (the DAQ system). The information is available through a special object which is obtained by calling the data set object's method env():

Code Block
bgColor#F7F7ED

dsname = 'exp=MEC/mec70813:run=35'
ds = DataSource(dsname)

env = ds.env()
env.
env.calibDir     env.configStore  env.expNum       env.fwkName      env.hmgr         env.jobName
env.calibStore   env.epicsStore   env.experiment   env.getConfig    env.instrument   env.subprocess

...

This is a sample analysis session illustrating a result of calling the methods:

Code Block
bgColor#F7F7ED

ds  = DataSource('exp=CXI/cxitut13:run=22')
env = ds.env()

print '   framework name:',env.fwkName()
print '         job name:',env.jobName()
print '       instrument:',env.instrument()
print '    experiment id:',env.expNum()
print '  experiment name:',env.experiment()
print 'subprocess number:',env.subprocess()

This will produce an an output which will look like this:

Code Block
bgColor#F7F7ED

   framework name: psana
         job name: cxitut13:run=22
       instrument: CXI
    experiment id: 304
  experiment name: cxitut13
subprocess number: 0

...

Consider the following example in which the first try to list configuration keys will result in an empty dictionary:

Code Block
bgColor#F7F7ED

ds = DataSource('exp=CXI/cxitut13:run=22')
configStore = ds.env().configStore()
configStore.keys()

[]

The second attempt made after getting to the first event will produce some meaningful output:

Code Block
bgColor#F7F7ED

itr = ds.events()
evt = itr.next()
configStore.keys()

[EventKey(type=psana.ControlData.ConfigV2, src='ProcInfo(0.0.0.0, pid=17877)'),
 EventKey(type=None, src='ProcInfo(0.0.0.0, pid=17877)'),
 EventKey(type=psana.EvrData.ConfigV7, src='DetInfo(NoDetector.0:Evr.0)'),
 EventKey(type=None, src='DetInfo(NoDetector.0:Evr.0)'),
 EventKey(type=psana.EvrData.ConfigV7, src='DetInfo(NoDetector.0:Evr.1)'),
 EventKey(type=None, src='DetInfo(NoDetector.0:Evr.1)'),
 EventKey(type=psana.EvrData.ConfigV7, src='DetInfo(NoDetector.0:Evr.2)'),
 EventKey(type=None, src='DetInfo(NoDetector.0:Evr.2)'),
 EventKey(type=psana.Epics.ConfigV1, src='DetInfo(EpicsArch.0:NoDevice.0)'),
 EventKey(type=None, src='DetInfo(EpicsArch.0:NoDevice.0)'),
 EventKey(type=psana.Epics.ConfigV1, src='DetInfo(EpicsArch.0:NoDevice.1)'),
 EventKey(type=None, src='DetInfo(EpicsArch.0:NoDevice.1)'),
 EventKey(type=psana.Acqiris.ConfigV1, src='DetInfo(CxiEndstation.0:Acqiris.0)'),
 EventKey(type=None, src='DetInfo(CxiEndstation.0:Acqiris.0)'),
 EventKey(type=psana.Ipimb.ConfigV2, src='DetInfo(CxiEndstation.0:Ipimb.0)'),
 EventKey(type=psana.Lusi.IpmFexConfigV2, src='DetInfo(CxiEndstation.0:Ipimb.0)'),
 EventKey(type=None, src='DetInfo(CxiEndstation.0:Ipimb.0)'),
 EventKey(type=None, src='DetInfo(CxiEndstation.0:Ipimb.0)'),
 EventKey(type=psana.Camera.FrameFexConfigV1, src='DetInfo(CxiEndstation.0:Opal4000.1)'),
 EventKey(type=psana.Opal1k.ConfigV1, src='DetInfo(CxiEndstation.0:Opal4000.1)'),
 EventKey(type=None, src='DetInfo(CxiEndstation.0:Opal4000.1)'),
 EventKey(type=None, src='DetInfo(CxiEndstation.0:Opal4000.1)'),
 EventKey(type=psana.Ipimb.ConfigV2, src='DetInfo(CxiDg1.0:Ipimb.0)'),
 EventKey(type=psana.Lusi.IpmFexConfigV2, src='DetInfo(CxiDg1.0:Ipimb.0)'),
 ..

The configuration keys have a structure which is reminiscent to the one of the event keys. This is illustrated below:

Code Block
bgColor#F7F7ED

In [140]: configStore.
configStore.get   configStore.keys

In [142]: keys = configStore.keys()

In [143]: keys[4]
Out[143]: EventKey(type=psana.EvrData.ConfigV7, src='DetInfo(NoDetector.0:Evr.1)')

In [144]: keys[4].type()
Out[144]: psana.EvrData.ConfigV7


In [145]: obj = configStore.get(EvrData.ConfigV7,Source('DetInfo(NoDetector.0:Evr.1)'))

In [146]: obj.
obj.eventcodes   obj.neventcodes  obj.noutputs     obj.npulses      obj.output_maps  obj.pulses       obj.seq_config

In [147]: type(obj)
Out[147]: psana.EvrData.ConfigV7

...

ControlPV is a configuration data which is updated on every step (steps are explained later in the document when discussing various ways of iterating over events in a data set). Like any other configuration data it is accessible through the environment object. Here is an example of getting controlPV data:

Code Block
bgColor#F7F7ED

from psana import *

ds = DataSource('exp=xpp66613:run=300:h5')

for step in ds.steps():
    control = ds.env().configStore().get(ControlData.Config)
    print [(c.name(), c.value()) for c in control.pvControls()]

This will result in the following output:

Code Block
bgColor#F7F7ED

[('lxt_ttc', -1.9981747581849466e-12)]
[('lxt_ttc', -1.8004943676675365e-12)]
[('lxt_ttc', -1.6001426205245978e-12)]
[('lxt_ttc', -1.3997908733800049e-12)]
[('lxt_ttc', -1.199439126235412e-12)]
[('lxt_ttc', -9.990873790924733e-13)]
[('lxt_ttc', -7.987356319478803e-13)]
[('lxt_ttc', -5.983838848049417e-13)]
[('lxt_ttc', -3.9803213766034873e-13)]
[('lxt_ttc', -2.0035174714459297e-13)]
[('lxt_ttc', 0.0)]
[('lxt_ttc', 2.0035174714459297e-13)]
[('lxt_ttc', 4.007034942875316e-13)]
[('lxt_ttc', 6.010552414321246e-13)]
[('lxt_ttc', 8.014069885767175e-13)]
..

...

All EPICS variables can be accessed through the EpicsStore object of the environment:

Code Block
bgColor#F7F7ED

ds = DataSource(dsname)
epics = ds.env().epicsStore()
epics.
epics.alias   epics.aliases  epics.getPV   epics.names
epics.pvName  epics.pvNames  epics.status  epics.value

The store interface allows:

  • obtaining official PV names of the EPICS variables (method namespvNames()).  These names are often obscure.
  • obtaining alias names which are known to the store (method aliases()).  Aliases are experiment-dependent simpler names specified by the users at the time the data is taken.
  • method names() shows both the aliases and the pvNames
  • return the alias for checking if there is an alias name for specified PV name (method alias()), or vs vice-versa (method pvName())
  • obtaining values of PVs (method value()).  This method can accept either a name or an alias.
  • obtaining descriptor objects for PVs (method getPV()).  This is useful for understanding the type/shape of the epics data, for example.
Note
titleHow to properly use the EPICS store

It's important to understand that the contents of the store is relevant to the most recent event obtained from a data set. This means that:

  • the store will be empty after a data set has been just open and no single event has been fetched from it
  • the contents of the store will change from one event to the other one

Therefore it's up to a user code to implement a correct logic for fetching events and EPICS variables to ensure that they're properly synchronized.

Here is a code which would be tracking values of some PV for all events and reporting events when the value of the PV changes. Notice that values of this (which is also true for most EPICS variables) won't change at each event:

Code Block
bgColor#F7F7ED

ds = DataSource(dsname)
epics = ds.env().epicsStore()
prev_val = None
for i, evt in enumerate(ds.events()):
    val = epics.value('LAS:FS0:ACOU:amp_rf1_17_2:rd')
    if val != prev_val:
        print "%6d:" % i, val
        prev_val = val


     0: 2016
   725: 2024
   845: 2019
   966: 2014
  1932: 2021
  2053: 2020
  2174: 2016
  3381: 2020
  3502: 2024
  4830: 2018
  6279: 2019
  6400: 2018
  7728: 2023
  7849: 2019
  9177: 2016
 10626: 2021
 10747: 2022
 12075: 2016
 13524: 2020
 ..

...

The Histogram Manager is the only public (user-level) service which is implemented in the current version of the framework. A reference to the manager can be obtain using:

Code Block
bgColor#F7F7ED

ds = DataSource('exp=CXI/cxitut13:run=22')
hist_manager = ds.env().hmgr()

...

The previous examples have already demonstrated the very basic technique for finding all events in a data set. The interactive psana has actually more elaborate ways of browsing through the data:

  • iterating over all events of a dataset:

    Code Block
    bgColor#F7F7ED
    
    ds = DataSource(...)
    for evt in ds.events():
        ...
    
  • iterating over all runs of a dataset, then iterating over all events of each run:

    Code Block
    bgColor#F7F7ED
    
    ds = DataSource('exp=...:run=12,13,14')
    for run in ds.runs():
        for evt in run.events():
            ...
    
  • iterating over all runs of a dataset, then iterating over all steps of each runs, then iterating over all events of each step:

    Code Block
    bgColor#F7F7ED
    
    ds = DataSource('exp=...:run=12,13,14')
    for run in ds.runs():
        for step in run.steps():
            for evt in step.events():
                ...
    
  • iterating directly over all steps of a dataset, then iterating over all events of each step:

    Code Block
    bgColor#F7F7ED
    
    ds = DataSource('exp=...:run=12,13,14')
    for step in ds.steps():
        for evt in step.events():
            ...
    

...

External modules are activated in the framework by mean of a specially prepared configuration file which has to be given to the framework before opening a data set. Otherwise the file won't make any effect. Here is this example:

Code Block
bgColor#F7F7ED

cfg = "/reg/g/psdm/tutorials/cxi/cspad_imaging/frame_reco.cfg"
setConfigFile(cfg)
ds = DataSource('exp=CXI/cxitut13:run=22')

...

In this section we're going to explore in a little bit more details the effect of the configuration file which was used in the previous. First, let's have a look at the contents of that file:

Code Block
bgColor#F7F7ED

[psana]
verbose         = 0
modules         = CSPadPixCoords.CSPadImageProducer

[CSPadPixCoords.CSPadImageProducer]

source          = CxiDs1.0:Cspad.0
typeGroupName   = CsPad::CalibV1
key             =
imgkey          = reconstructed
tiltIsApplied   = true
print_bits      = 0

This configuration refers to some real module (CSPadImageProducer from the OFFLINE Analysis package CSPadPixCoords) which will be doing the geometric reconstruction of the full CSPad image. This module is a part of any latest analysis releases. The main effect of the module is that it will extend each event by an additional component which otherwise wouldn't be present in the event (look for the first one in the output below):

Code Block
bgColor#F7F7ED

evt.keys()
..
EventKey(type=psana.ndarray_int16_2, src='DetInfo(CxiDs1.0:Cspad.0)', key='reconstructed'),
..
EventKey(type=psana.CsPad.DataV2, src='DetInfo(CxiDsd.0:Cspad.0)'),
..

That component has an additional key='reconstructed' (please, refer back to the section on the Four forms of the get() method for more information on that parameter). An object returned with that key has a type of a 2D array of 16-bit elements representing CSPad pixels. With a little bit of help from Matplotlib one can easily turn this into an image. The next example will illustrate how to show both raw (unprocessed) and reconstructed version of the CSPad image:

Code Block
bgColor#F7F7ED

import matplotlib.pyplot as plt

plt.figure()
plt.ion()
plt.show()

cspad = evt.get(CsPad.DataV2,Source('DetInfo(CxiDs1.0:Cspad.0)'))
a = []
for i in range(0,4):
    quad = cspad.quads(i)
    d = quad.data()
    a.append(np.vstack([d[i] for i in range(0,8)]))

frame_raw = np.hstack(a)

frame_reconstructed = evt.get(ndarray_int16_2,Source('DetInfo(CxiDs1.0:Cspad.0)'),'reconstructed')

plt.imshow(frame_raw)
plt.clim(850,1200)

plt.imshow(frame_reconstructed)
plt.clim(850,1200)

...

The only caveat with making multiple calls to the DataSource() method is that the very last call to the setConfigFile() operation will affect any instances of the framework open with DataSource(). In other words, an order in which functions setConfigFile() and DataSource() does matter. Consider the following example:

Code Block
bgColor#F7F7ED

ds1 = DataSource(dsname1)   ## no configuration is assumed

setConfigFile('c1.cfg')

ds2 = DataSource(dsname2)   ## the dataset will be processed with c1.cfg
ds3 = DataSource(dsname3)   ## the dataset will be also processed with c1.cfg

setConfigFile('c2.cfg')

ds4 = DataSource(dsname4)   ## the dataset will be processed with c2.cfg
...

Re-opening a dataset is no different from opening multiple different data sets:

Code Block
bgColor#F7F7ED

ds = DataSource(dsname)

for evt in ds.events():
    ...

ds = DataSource(dsname)     ## this is a fresh dataset object in which
                            ## all iterators are poised to the very first event (run, step)

...

Unlike its batch version, the interactive psana allows multiple events to be present at a time within the process memory. Consider the following extreme scenario:

Code Block
bgColor#F7F7ED

ds = DataSource(dsname)

events = [evt for evt in ds.events()]

## Now all events are in memory, hence they can be addressed directly
## from the list

evt = events[123]
print evt.run()

However, if you aren't so lucky, and a machine doesn't have enough memory then you would see the run-time exception:

Code Block
bgColor#F7F7ED

RuntimeError: St9bad_alloc

The second thing to worry about is to make sure all relevant environment data are properly handled. In particular, the event environment objects should be obtained and stored along with each cached event at the same time the event object is retrieved from a dataset and before moving to the next event. The problem was already mentioned in this document when discussing the EPICS Store. As an illustration, let's suppose we need to compare images stored in each pair of consecutive events for all events in a run, and to do the comparison we also need the corresponding values of some PV. In that case the correct algorithm may look like this:

Code Block
bgColor#F7F7ED

ds = DataSource('exp=XCS/xcstut13:run=15')
epics = ds.env().epicsStore()

prev_evt = None
prev_pv  = None
for evt in ds.events()
    pv = epics.value('LAS:FS0:ACOU:amp_rf1_17_2:rd')
    if prev_evt  is not None:
        ...  ## compare (prev_evt,pv) vs (evt,pv)
    prev_evt = evt
    prev_pv  = pv

...

Consider the following example:

Code Block
bgColor#F7F7ED

ds = DataSource('exp=XCS/xcstut13:run=15')

for evt in ds.events():
    frame = evt.get(Princeton.FrameV1, Source('DetInfo(XcsBeamline.0:Princeton.0)'))
    ..

One problem with that code is that it would construct the component address object (using Source()) at every single step of the iteration. This may cost you some slowdown in your application's performance due to an overhead of parsing the address string and constructing a new source object at each step of the iteration. That won't be necessary, and the example can be easily rewritten like:

Code Block
bgColor#F7F7ED

ds = DataSource('exp=XCS/xcstut13:run=15')
src = Source('DetInfo(XcsBeamline.0:Princeton.0)')

for evt in ds.events():
    frame = evt.get(Princeton.FrameV1, src)
    ..

...

The internal implementation of psana won't fully construct components of an event unless they are requested by a user's code. This will make the following code less efficient when fetching only those components which are actually needed by an application:

Code Block
bgColor#F7F7ED

ds = DataSource('exp=XCS/xcstut13:run=15')
for evt in ds.events():
    dets = [evt.get(k.type(),k.src()) for k in evt.keys() if (k.type() != EventId) and (k.type() != None)]
    # now all component objects are stored in the list
    ...

...

Reading images for all events of a dataset is better to be done with the XTC format:

Code Block
bgColor#F7F7ED

ds  = DataSource('exp=XCS/xcstut13:run=15')
src = Source('DetInfo(XcsBeamline.0:Princeton.0)')

for evt in ds.events():
    frame = evt.get(Princeton.FrameV1, src)
    ..

Reading EPICS variables for all events of a dataset is better to be done with the HDF5 format:

Code Block
bgColor#F7F7ED

ds = DataSource('exp=XCS/xcstut13:run=15:h5')
epics = ds.env().epicsStore()

for evt in ds.events()
    pv = epics.value('LAS:FS0:ACOU:amp_rf1_17_2:rd')
    ...

Iterating over the first event of each step is better to be done with the HDF5 format. This is actually a very good example where the direct access to events helps to avoid unnecessary I/O operations when jumping between steps. And obviously, this benefit will only be seen in those datasets which have many steps:

Code Block
bgColor#F7F7ED

ds = DataSource('exp=XCS/xcstut13:run=15')
for step in ds.steps():
    itr = step.events()
    evt = itr.next()   ## do something about this event and then proceed
    ..                 ## to the next step

...