CFEL - ASG Software Suite  2.5.0
CASS
roottree_converter.cpp
Go to the documentation of this file.
1 // Copyright (C) 2011 Lutz Foucar
2 
3 /**
4  * @file roottree_converter.cpp file contains definition of processor 2001
5  *
6  * @author Lutz Foucar
7  */
8 
9 #include <sstream>
10 #include <string>
11 #include <vector>
12 #include <stdexcept>
13 #include <iostream>
14 
15 #include <TObject.h>
16 #include <TFile.h>
17 #include <TTree.h>
18 #include <TDirectory.h>
19 #include <TH1.h>
20 #include <TH2.h>
21 
22 #include <QtCore/QDateTime>
23 #include <QtCore/QString>
24 #include <QtCore/QStringList>
25 
26 #include "roottree_converter.h"
27 #include "result.hpp"
28 #include "cass_settings.h"
29 #include "cass_event.h"
30 #include "convenience_functions.h"
31 #include "delayline_detector.h"
32 #include "tree_structure.h"
33 #include "machine_device.hpp"
34 #include "rootfile_helper.h"
35 #include "log.h"
36 #include "processor_manager.h"
37 
38 
39 using namespace cass;
40 using namespace std;
41 using namespace ACQIRIS;
42 using namespace MachineData;
43 
44 namespace cass
45 {
46 /** typedef for easier code */
47 typedef HelperAcqirisDetectors::helperinstancesmap_t::key_type detectorkey_t;
48 
49 /** typedef for easier code */
50 typedef list<detectorkey_t> dlddetectors_t;
51 
52 /** typedef for easier code */
53 typedef pair<DelaylineDetector::particles_t::key_type,
54  HelperAcqirisDetectors::helperinstancesmap_t::key_type> particleskey_t;
55 
56 /** typedef for easier code */
57 typedef list<particleskey_t> particleslist_t;
58 
59 /** load the settings of all acqiris detectors defined in .ini file
60  *
61  * @author Lutz Foucar
62  */
64 {
65  CASSSettings s;
66  s.beginGroup("AcqirisDetectors");
67  QStringList detectorNamesList(s.childGroups());
68  QStringList::const_iterator detName(detectorNamesList.begin());
69  for (; detName != detectorNamesList.end(); ++detName)
70  HelperAcqirisDetectors::instance(detName->toStdString())->loadSettings();
71 }
72 
73 /** check whether the key points to a delayline detector
74  *
75  * @return true when key points to a delayline detector
76  * @param detkey key for the possible delayline detector
77  *
78  * @author Lutz Foucar
79  */
80 bool isDLD(const detectorkey_t &detkey)
81 {
82  HelperAcqirisDetectors::shared_pointer dethelp (HelperAcqirisDetectors::instance(detkey));
83  return (dethelp->detectortype() == Delayline);
84 }
85 
86 /** convert a qstring to the key in the list of detectors
87  *
88  * will convert the requested delayline detector described in the qstring to the
89  * delayline detector key. Before returning the key, check whether the requested
90  * detector is realy a delayline detector and whether it is defined in the ini
91  * file.
92  *
93  * @return the key to the requested detector
94  * @param qstr the QString to convert to the key
95  *
96  * @author Lutz Foucar
97  */
98 detectorkey_t qstring2detector(const QString & qstr)
99 {
100  detectorkey_t dld(qstr.toStdString());
101  if (!isDLD(dld))
102  throw invalid_argument("pp2001::loadSettings(): Error detector '" + dld +
103  "' is not a Delaylinedetector.");
104  CASSSettings s;
105  s.beginGroup("AcqirisDetectors");
106  QStringList detectorNamesList(s.childGroups());
107  if (!detectorNamesList.contains(qstr))
108  throw invalid_argument("pp2001::loadSettings(): Error detector '" + dld +
109  "' is not defined.");
110  return dld;
111 }
112 
113 /** convert qstring to the particle key pair
114  *
115  * first try to find the detector that the particle belongs to. If found which
116  * detector the particle belongs to then add this as second part of the returned
117  * key.
118  *
119  * @return the key pair containing the particle name and the key to the detector
120  * that the particle belongs to
121  * @param qstr the QString of the particle that is requested
122  *
123  * @author Lutz Foucar
124  */
125 particleskey_t qstring2particle(const QString & qstr)
126 {
127  particleskey_t particlekey(make_pair(qstr.toStdString(),""));
128  const HelperAcqirisDetectors::helperinstancesmap_t & knownDetectors(HelperAcqirisDetectors::instances());
129  HelperAcqirisDetectors::helperinstancesmap_t::const_iterator det(knownDetectors.begin());
130  for (;det != knownDetectors.end(); ++det)
131  {
132  if (isDLD(det->first))
133  {
134  HelperAcqirisDetectors &dethelp(*HelperAcqirisDetectors::instance(det->first));
135  const DetectorBackend &detback((dethelp.detector()));
136  const DelaylineDetector &dld(dynamic_cast<const DelaylineDetector&>(detback));
137  const DelaylineDetector::particles_t &particles(dld.particles());
138  DelaylineDetector::particles_t::const_iterator particle(particles.find(particlekey.first));
139  if (particle != particles.end())
140  {
141  particlekey.second = det->first;
142  return particlekey;
143  }
144  }
145  }
146  throw invalid_argument("pp2001::loadSettings(): Error particle '" + particlekey.first +
147  "' is not part of any Delaylinedetector.");
148  return particlekey;
149 }
150 
151 /** copy map values to map
152  *
153  * will copy each key pair to the map of the tree structure.
154  *
155  * @param first iterator to the first element to copy
156  * @param last const iterator to one past the last element to copy
157  * @param dest the destination where the elements will be copied to
158  *
159  * @author Lutz Foucar
160  */
161 void copyMapValues(map<string,double>::const_iterator first,
162  map<string,double>::const_iterator last,
163  treehit_t& dest)
164 {
165  while(first != last)
166  {
167  dest[(*first).first] = (*first).second;
168  ++first;
169  }
170 }
171 }//end namespace cass
172 
173 pp2001::pp2001(const name_t &name, std::string filename)
174  : Processor(name),
175  _rootfile(ROOTFileHelper::create(filename)),
176  _tree(new TTree("CASSData","Selected preprocessed data from the CASSEvent")),
177  _treestructure_ptr(&_treestructure),
178  _machinestructure_ptr(&_machinestructure),
179  _eventstatusstructure_ptr(&_eventstatusstructure),
180  _ppstructure_ptr(&_ppstructure)
181 {
182  if (!_rootfile)
183  throw invalid_argument("pp2001 (" + name + "): '" + filename +
184  "' could not be opened! Maybe deleting the file helps.");
185  loadSettings(0);
186 }
187 
189 {
190  throw logic_error("pp2001::result: '"+name()+"' should never be called");
191 }
192 
194 {
195  CASSSettings settings;
196  settings.beginGroup("Processor");
198  setupGeneral();
199  QStringList pps(settings.value("Processors").toStringList());
200  QStringList::const_iterator ppname(pps.begin());
201  bool allDepsAreThere(true);
202  for (ppname = pps.begin(); ppname != pps.constEnd(); ++ppname)
203  {
204  shared_pointer pp(setupDependency("",ppname->toStdString()));
205  allDepsAreThere = pp && allDepsAreThere;
206  if (pp && pp->result().dim() != 0)
207  throw invalid_argument("pp2001 (" + name() + "): Processor '" + pp->name() +
208  "' is not a 0d histogram.");
209  _pps.push_back(pp);
210  }
211  bool ret (setupCondition());
212  if (!(ret && allDepsAreThere))
213  {
214  _pps.clear();
215  return;
216  }
217  loadAllDets();
218  QStringList detectors(settings.value("Detectors").toStringList());
219  _detectors.resize(detectors.size());
220  transform(detectors.begin(),detectors.end(),_detectors.begin(),qstring2detector);
221  QStringList particles(settings.value("Particles").toStringList());
222  _particles.resize(particles.size());
223  transform(particles.begin(),particles.end(),_particles.begin(),qstring2particle);
224  if (_tree->FindBranch("DLDetectorData") == 0)
225  if (!_detectors.empty() || !_particles.empty())
226  _tree->Branch("DLDetectorData","map<string,vector<map<string,double> > >",&_treestructure_ptr);
227  if (_tree->FindBranch("EvendId") == 0)
228  _tree->Branch("EvendId",&_eventid,"id/l");
229  if (settings.value("MachineData",false).toBool())
230  if (_tree->FindBranch("MachineData") == 0)
231  _tree->Branch("MachineData","map<string,double>",&_machinestructure_ptr);
232  if (settings.value("EventStatus",false).toBool())
233  if (_tree->FindBranch("EventStatus") == 0)
234  _tree->Branch("EventStatus","vector<bool>",&_eventstatusstructure_ptr);
235  if (!_pps.empty())
236  if (_tree->FindBranch("Processors") == 0)
237  _tree->Branch("Processors","map<string,double>",&_ppstructure_ptr);
238 
239  _hide = true;
240  string output("Processor '" + name() + "' will write the hits of detectors: ");
241  dlddetectors_t::const_iterator detectorsIt(_detectors.begin());
242  for (;detectorsIt!=_detectors.end();++detectorsIt)
243  output += ("'" + (*detectorsIt) + "', ");
244  output += " and particles: ";
245  particleslist_t::const_iterator particle(_particles.begin());
246  for (;particle != _particles.end();++particle)
247  output += ("'" + particle->first + "(" + particle->second + ")', ");
248  output += (" to rootfile '" + string(_rootfile->GetName()) + "'. Condition is '" +
249  _condition->name() + "'");
250  Log::add(Log::INFO,output);
251 }
252 
254 {
255  _tree->Write();
257 }
258 
260 {
261  if (!_condition->result(evt.id()).isTrue())
262  return;
263  QMutexLocker locker(&_lock);
264 
265  _eventid = evt.id();
266  dlddetectors_t::const_iterator detector(_detectors.begin());
267  dlddetectors_t::const_iterator detectorEnd(_detectors.end());
268  for (;detector != detectorEnd; ++detector)
269  {
270  DetectorBackend &rawdet(
271  HelperAcqirisDetectors::instance(*detector)->detector(evt));
272  DelaylineDetector &det(dynamic_cast<DelaylineDetector&>(rawdet));
273 
275  treedet.clear();
276  detectorHits_t::iterator hit(det.hits().begin());
277  for (; hit != det.hits().end(); ++hit)
278  {
279  treehit_t treehit;
280  detectorHit_t &hitvalues(*hit);
281 // copyMapValues(hitvalues.begin(),hitvalues.end(),treehit);
282  treehit["x"] = hitvalues[x];
283  treehit["y"] = hitvalues[y];
284  treehit["t"] = hitvalues[t];
285  treehit["method"] = hitvalues[method];
286  treedet.push_back(treehit);
287  }
288  }
289  particleslist_t::const_iterator particle(_particles.begin());
290  particleslist_t::const_iterator particleEnd(_particles.end());
291  for (;particle !=particleEnd;++particle)
292  {
293  DetectorBackend &rawdet(
294  HelperAcqirisDetectors::instance(particle->second)->detector(evt));
295  DelaylineDetector &det(dynamic_cast<DelaylineDetector&>(rawdet));
296 
297  treedetector_t &treeparticle(_treestructure[particle->first]);
298  treeparticle.clear();
299  particleHits_t & hits(det.particles()[particle->first].hits());
300  particleHits_t::iterator hit(hits.begin());
301  particleHits_t::iterator hitEnd(hits.end());
302  for (; hit != hitEnd; ++hit)
303  {
304  treehit_t treehit;
305  particleHit_t &hitvalues(*hit);
306 // copyMapValues(hitvalues.begin(),hitvalues.end(),treehit);
307  treehit["px"] = hitvalues[px];
308  treehit["py"] = hitvalues[py];
309  treehit["pz"] = hitvalues[pz];
310  treehit["x_mm"] = hitvalues[x_mm];
311  treehit["y_mm"] = hitvalues[y_mm];
312  treehit["tof_ns"] = hitvalues[tof_ns];
313  treehit["xCor_mm"] = hitvalues[xCor_mm];
314  treehit["yCor_mm"] = hitvalues[yCor_mm];
315  treehit["xCorScal_mm"] = hitvalues[xCorScal_mm];
316  treehit["yCorScal_mm"] = hitvalues[yCorScal_mm];
317  treehit["xCorScalRot_mm"] = hitvalues[xCorScalRot_mm];
318  treehit["yCorScalRot_mm"] = hitvalues[yCorScalRot_mm];
319  treehit["tofCor_ns"] = hitvalues[tofCor_ns];
320  treehit["roh"] = hitvalues[roh];
321  treehit["theta"] = hitvalues[theta];
322  treehit["phi"] = hitvalues[phi];
323  treehit["e_au"] = hitvalues[e_au];
324  treehit["e_eV"] = hitvalues[e_eV];
325  treeparticle.push_back(treehit);
326  }
327  }
328  const Device &machinedata
329  (dynamic_cast<const Device&>(*(evt.devices().find(cass::CASSEvent::MachineData)->second)));
330  copyMapValues(machinedata.BeamlineData().begin(), machinedata.BeamlineData().end(), _machinestructure);
331  copyMapValues(machinedata.EpicsData().begin(), machinedata.EpicsData().end(), _machinestructure);
332  _eventstatusstructure.resize(machinedata.EvrData().size());
333  copy(machinedata.EvrData().begin(),machinedata.EvrData().end(),_eventstatusstructure.begin());
334 
335  /** copy the values of each 0d Processor into the processor structure */
336  std::list<shared_pointer>::const_iterator ProcessorsIt(_pps.begin());
337  std::list<shared_pointer>::const_iterator ProcessorsEnd(_pps.end());
338  for (;ProcessorsIt != ProcessorsEnd;++ProcessorsIt)
339  {
340  Processor &pp(*(*ProcessorsIt));
341  const result_t &val(pp.result(_eventid));
342  QReadLocker lock(&val.lock);
343  float value(val.getValue());
344  _ppstructure[pp.name()] = value;
345  }
346  _tree->Fill();
347 }
348 
std::vector< double > detectorHit_t
define a detector hit
eventStatus_t * _eventstatusstructure_ptr
pointer to the event status structure
Event to store all LCLS Data.
Definition: cass_event.h:32
uint64_t _eventid
copy of the event id
HelperAcqirisDetectors::helperinstancesmap_t::key_type detectorkey_t
typedef for easier code
pp2001(const name_t &, std::string)
Construct processor for converting histograms to root histograms.
virtual const result_t & result(const CASSEvent::id_t eventid=0)
overwrite the retrieval of an histogram
const name_t name() const
retrieve the name of this processor
Definition: processor.h:167
file contains declaration of the CASSEvent
Settings for CASS.
Definition: cass_settings.h:30
uint64_t id_t
define the id type
Definition: cass_event.h:52
define which of the hitfinders defined above will be used as hit
bool _hide
flag to tell whether this pp should be hidden in the dropdown list
Definition: processor.h:262
machinestructure_t * _machinestructure_ptr
pointer to the machine structure defined above
STL namespace.
Average out the iShit status to get the avererage hits
void loadAllDets()
load the settings of all acqiris detectors defined in .ini file
std::map< std::string, double > treehit_t
a hit is just a map of hit values identified by a string
bool isDLD(const detectorkey_t &detkey)
check whether the key points to a delayline detector
result classes
simple pixel detectors
list< detectorkey_t > dlddetectors_t
typedef for easier code
file contains declaration of processor 2001
TTree * _tree
the root tree to fill
static void add(Level level, const std::string &line)
add a string to the log
Definition: log.cpp:31
file contains the classes that describe a delayline detector.
TFile * _rootfile
the root file
std::vector< double > particleHit_t
define a particle hit
fromStdString(const std::string &str)
QMutex _lock
a lock to make the process reentrant
particleskey_t qstring2particle(const QString &qstr)
convert qstring to the particle key pair
treestructure_t _treestructure
structure that should be written to tree
base class for processors.
Definition: processor.h:39
shared_pointer setupDependency(const std::string &depVarName, const name_t &name="")
setup the dependecy.
Definition: processor.cpp:114
devices_t & devices()
setters
Definition: cass_event.h:66
file contains declaration of classes and functions that help other processors to do their job...
root file creation
contains the manager for the processors
id_t & id()
setters
Definition: cass_event.h:64
list< particleskey_t > particleslist_t
typedef for easier code
std::list< ACQIRIS::HelperAcqirisDetectors::helperinstancesmap_t::key_type > _detectors
list of detectors who's hits should be filled into the tree
contains singleton definition for creating root files
value(const QString &key, const QVariant &defaultValue=QVariant()
std::list< shared_pointer > _pps
container for all 0d Processors that should be added to the tree
virtual void processEvent(const CASSEvent &)
only a stub does nothing, but needs to be there because its pure virtual in base class ...
static void close(TFile *rootfile)
close root file
void setupGeneral()
general setup of the processor
Definition: processor.cpp:85
std::vector< treehit_t > treedetector_t
lots of hits are a delayline detector
pair< DelaylineDetector::particles_t::key_type, HelperAcqirisDetectors::helperinstancesmap_t::key_type > particleskey_t
typedef for easier code
definitions of a machine device
void copyMapValues(map< string, double >::const_iterator first, map< string, double >::const_iterator last, treehit_t &dest)
copy map values to map
void loadSettings(CASSSettings &s, CFDParameters &p, uint32_t &instrument, size_t &channelNbr)
implementation of loading settings for both CFD classes
Definition: cfd.cpp:210
virtual void loadSettings(size_t)
load the settings of this pp
file contains specialized class that do the settings for cass
ppstructure_t _ppstructure
0d processor structure
shared_pointer _condition
pointer to the processor that will contain the condition
Definition: processor.h:277
eventStatus_t _eventstatusstructure
event status structure
Electron detector
Definition: hdf5-input.ini:62
bool setupCondition(bool defaultConditionType=true)
setup the condition.
Definition: processor.cpp:94
std::string name_t
define the name type
Definition: processor.h:46
contains a logger for cass
std::tr1::shared_ptr< Processor > shared_pointer
a shared pointer of this
Definition: processor.h:43
virtual void aboutToQuit()
dump all histogram to a root file just before quitting
std::vector< particleHit_t > particleHits_t
define container for all particle hits
check if there is some light in the chamber based upon the GMD value
detectorkey_t qstring2detector(const QString &qstr)
convert a qstring to the key in the list of detectors
defining structures for the root tree
beginGroup(const QString &prefix)
machinestructure_t _machinestructure
machine data structure
treestructure_t * _treestructure_ptr
pointer to the above structure (needed by the tree)
ppstructure_t * _ppstructure_ptr
pointer to the 0d processor structure
virtual const result_t & result(const CASSEvent::id_t eventid=0)
retrieve a result for a given id.
Definition: processor.cpp:54
std::list< std::pair< ACQIRIS::DelaylineDetector::particles_t::key_type, ACQIRIS::HelperAcqirisDetectors::helperinstancesmap_t::key_type > > _particles
list of particles whos hits should be added to the tree