KaliVeda  1.12/06
Heavy-Ion Analysis Toolkit
KVGANILDataReader.cpp
Go to the documentation of this file.
1 //Created by KVClassFactory on Wed Sep 23 16:07:38 2009
2 //Author: Chbihi
3 
4 #include "KVGANILDataReader.h"
5 #include "GTGanilData.h"
6 #include "GTOneScaler.h"
7 #include "TSystem.h"
8 #include "TUrl.h"
9 #include "TPluginManager.h"
10 #include "RVersion.h"
11 #include <TInterpreter.h>
12 #include "TEnv.h"
13 
15 
16 
17 
18 
27 {
28  //Open and initialise a GANIL data file for reading.
29  //By default, Scaler buffers are ignored.
30  //If file cannot be opened, this object will be made a zombie. Do not use.
31  //To test if file is open, use IsOpen().
32  //The basename of the file (excluding any path) can be obtained using GetName()
33  //The full pathname of the file can be obtained using GetTitle()
34 
35  init();
36  OpenFile(file, dataset);
37 }
38 
39 
40 
43 
45 {
46  // Destructor
47  if (fGanilData) {
48  delete fGanilData;
49  fGanilData = 0;
50  }
51  fParameters->Clear();
52  delete fParameters;
53  if (fExtParams) {
54  fExtParams->Delete("slow");
55  delete fExtParams;
56  }
57  delete fFired;
58  if (ParVal) delete [] ParVal;
59  if (ParNum) delete [] ParNum;
60 }
61 
62 
63 
66 
68 {
69  //default initialisations
70 
71  fExtParams = 0;
72  fParameters = new KVHashList;
74  fGanilData = 0;
75  fUserTree = 0;
76  fFired = new KVHashList;
77  ParVal = 0;
78  ParNum = 0;
80 }
81 
82 
83 
150 
152 {
153  // To fill a TTree with the data in the current file, create a TTree:
154  // TFile* file = new TFile("run1.root","recreate");
155  // TTree* T = new TTree("Run1", "Raw data for Run1");
156  // and then call this method: SetUserTree(T)
157  // If you read all events of the file, the TTree will be automatically filled
158  // with data :
159  // while( runfile->GetNextEvent() ) ;
160  //
161  // Two different TTree structures are available, depending on the option string:
162  //
163  // opt = "arrays": [default]
164  //
165  // The TTree will have the following structure:
166  //
167  // *Br 0 :NbParFired : NbParFired/I = number of fired parameters in event
168  // *............................................................................*
169  // *Br 1 :ParNum : ParNum[NbParFired]/i = array of indices of fired parameters
170  // *............................................................................*
171  // *Br 2 :ParVal : ParVal[NbParFired]/s = array of values of fired parameters
172  //
173  // This structure is the fastest to fill and produces the smallest file sizes.
174  // In order to be able to directly access the parameters as if option "leaves" were used
175  // (i.e. one branch/leaf for each parameter), we add two aliases for each parameter to
176  // the tree:
177  // PARNAME = value of parameter if present in event
178  // PARNAME_M = number of times parameter appears in event
179  // Assuming that each parameter only appears at most once in each event, i.e. PARNAME_M=0 or 1,
180  // then
181  // root[0] T->Draw("PARNAME", "PARNAME_M")
182  // will histogram the value of PARNAME for each event in which it is present.
183  // (if the selection condition "PARNAME_M" is not used, the histogram will also be filled with a 0
184  // for each event in which PARNAME does not appear).
185  // N.B. the PARNAME alias is in fact the sum of the values of PARNAME in each event.
186  // If PARNAME_M>1 in some events, it is not the individual values but their sum which will
187  // be histogrammed in this case.
188  //
189  // Thus, if the data file has parameters called "PAR_1" and "PAR_2",
190  // the following command will work
191  //
192  // root[0] T->Draw("PAR_1:PAR_2", "PAR_1_M&&PAR_2_M", "col")
193  //
194  // even though no branches "PAR_1" or "PAR_2" exist.
195  //
196  //
197  //
198  // opt = "leaves":
199  //
200  // The TTree will have a branch/leaf for each parameter. This option is slower and produces
201  // larger files.
202  //
203  // If the option string contains both "arrays" and "leaves", then both structures will be used
204  // (in this case there is a high redundancy, as each parameter is stored twice).
205  //
206  // The full list of parameters is stored in a TObjArray in the list returned by TTree::GetUserInfo().
207  // Each parameter is represented by a TNamed object.
208  // In order to retrieve the name of the parameter with index 674 (e.g. taken from branch ParNum),
209  // do:
210  // TObjArray* parlist = (TObjArray*) T->GetUserInfo()->FindObject("ParameterList");
211  // cout << "Par 674 name = " << (*parlist)[674]->GetName() << endl;
212  //
213  //
214  // Automatic creation & filling of Scalers TTree
215  //
216  // give an option string containing "scalers", i.e. "leaves,scalers", or "ARRAYS+SCALERS", etc.
217  // a TTree with name 'Scalers' will be created, all scaler buffers will be written in it.
218 
219 
220  TString option = opt;
221  option.ToUpper();
222  make_arrays = option.Contains("ARRAYS");
223  make_leaves = option.Contains("LEAVES");
224  Bool_t make_scalers = option.Contains("SCALERS");
225  if (make_scalers) {
227  }
228 
229  fUserTree = T;
230  if (make_arrays) {
231  Int_t maxParFired = GetRawDataParameters()->GetEntries();
232  ParVal = new UShort_t[maxParFired];
233  ParNum = new UInt_t[maxParFired];
234  fUserTree->Branch("NbParFired", &NbParFired, "NbParFired/I");
235  fUserTree->Branch("ParNum", ParNum, "ParNum[NbParFired]/i");
236  fUserTree->Branch("ParVal", ParVal, "ParVal[NbParFired]/s");
237  }
238  if (make_leaves) {
239  TIter next_rawpar(GetRawDataParameters());
240  KVACQParam* acqpar;
241  while ((acqpar = (KVACQParam*)next_rawpar())) {
242  TString leaf;
243  leaf.Form("%s/S", acqpar->GetName());
244  // for parameters with <=8 bits only use 1 byte for storage
245  if (acqpar->GetNbBits() <= 8) leaf += "1";
246  fUserTree->Branch(acqpar->GetName(), *(acqpar->ConnectData()), leaf.Data());
247  }
248  }
249 
250 #if ROOT_VERSION_CODE > ROOT_VERSION(5,25,4)
251 #if ROOT_VERSION_CODE < ROOT_VERSION(5,26,1)
252  // The TTree::OptimizeBaskets mechanism is disabled, as for ROOT versions < 5.26/00b
253  // this lead to a memory leak
255 #endif
256 #endif
257 
258  // add list of parameter names in fUserTree->GetUserInfos()
259  // and if option="arrays" add aliases for each parameter & its multiplicity
260 
261  // TObjArray has to be as big as the largest parameter number in the list
262  // of raw data parameters. So first loop over parameters to find max param number.
263  UInt_t maxpar = 0;
264  TIter next(GetRawDataParameters());
265  KVACQParam* par;
266  while ((par = (KVACQParam*)next())) if (par->GetNumber() > maxpar) maxpar = par->GetNumber();
267 
268  TObjArray* parlist = new TObjArray(maxpar, 1);
269  parlist->SetName("ParameterList");
270 
271  next.Reset();
272  while ((par = (KVACQParam*)next())) {
273  parlist->AddAt(new TNamed(par->GetName(), Form("index=%d", par->GetNumber())), par->GetNumber());
274  if (make_arrays) {
275  fUserTree->SetAlias(par->GetName(), Form("Sum$((ParNum==%d)*ParVal)", par->GetNumber()));
276  fUserTree->SetAlias(Form("%s_M", par->GetName()), Form("Sum$(ParNum==%d)", par->GetNumber()));
277  }
278  }
279  fUserTree->GetUserInfo()->Add(parlist);
280 }
281 
282 
283 
295 
297 {
298  //Open and initialise a GANIL data file for reading.
299  //By default, scaler buffers are ignored.
300  //This can be changed by changing the value of
301  //
302  // KVGANILDataReader.ScalerBuffersManagement: kSkipScaler
303  //
304  //If file cannot be opened, this object will be made a zombie. Do not use.
305  //To test if file is open, use IsOpen().
306  //The basename of the file (excluding any path) can be obtained using GetName()
307  //The full pathname of the file can be obtained using GetTitle()
308 
309  if (fGanilData) {
310  delete fGanilData;
311  fGanilData = 0;
312  }
313 
314  fGanilData = NewGanTapeInterface(dataset);
315 
318  SetTitle(file);
319 
320  // handling scaler buffers
321  TString what_scale = gEnv->GetValue("KVGANILDataReader.ScalerBuffersManagement", "kSkipScaler");
322  what_scale.Prepend("GTGanilData::");
323  Long_t ws = gInterpreter->ProcessLine(what_scale);
325 
326  fGanilData->Open();
327 
328  //test whether file has been opened
329  if (!fGanilData->IsOpen()) {
330  //if initial attempt fails, we try to open the file as a 'raw' TFile
331  //This may work when we are attempting to open a remote file i.e.
332  //via rfio, and a subsequent attempt to open the file using the GanTape
333  //routines may then succeed.
334  TUrl rawtfile(file, kTRUE);
335  rawtfile.SetOptions("filetype=raw");
336  TFile* rawfile = TFile::Open(rawtfile.GetUrl());
337  if (!rawfile) {
338  Error("OpenFile", "File cannot be opened: %s",
339  file);
340  MakeZombie();
341  return;
342  }
343  //TFile::Open managed to open file! Try again...
344  delete rawfile;
345  fGanilData->Open();
346  if (!fGanilData->IsOpen()) {
347  //failed again ??!
348  Error("OpenFile", "File cannot be opened: %s",
349  file);
350  MakeZombie();
351  return;
352  }
353  }
354 }
355 
356 
357 
358 
373 
375 {
376  //Call this method in order to initialise links between data in file and KVACQParam
377  //objects associated with detectors if defined. Call with no argument in order to
378  //generate all required KVACQParam objects corresponding to parameters in file
379  //(case where no KVMultiDetArray is defined).
380  //
381  //fParameters is filled with a KVACQParam for every acquisition parameter in the file.
382  //The list contains all KVACQParams defined for the detectors of a multidetector array.
383  //For any parameters for which no KVACQParam is found in the list we create new KVACQParam
384  //objects which will be deleted with this KVGANILDataReader (these objects can be accessed from the list
385  //returned by GetUnknownParameters()).
386  //
387  //To access the full list of data parameters in the file after this method has been
388  //called (i.e. after the file is opened), use GetRawDataParameters().
389 
391  KVACQParam* par;
392  GTDataPar* daq_par;
393  while ((daq_par = (GTDataPar*) next())) {//loop over all parameters
394  par = CheckACQParam(list_acq_params, daq_par->GetName());
395  fGanilData->Connect(par->GetName(), par->ConnectData());
396  par->SetNumber(daq_par->Index());
397  par->SetNbBits(daq_par->Bits());
398  fParameters->Add(par);
399  }
400 }
401 
402 
403 
404 
410 
411 KVACQParam* KVGANILDataReader::CheckACQParam(const TSeqCollection* list_acq_params, const Char_t* par_name)
412 {
413  //Check the named acquisition parameter.
414  //We look for a corresponding parameter in the list of acq params belonging to
415  //the current KVMultiDetArray (if one exists).
416  //If none is found, we create a new acq param which is added to the list of "unknown parameters"
417 
418  KVACQParam* par;
419  if (!list_acq_params || !(par = (KVACQParam*)list_acq_params->FindObject(par_name))) {
420  //create new unknown parameter
421  par = new KVACQParam;
422  par->SetName(par_name);
423  if (!fExtParams) {
424  fExtParams = new KVHashList;
426  }
427  fExtParams->Add(par);
428  }
429  return par;
430 }
431 
432 
433 
434 
442 
444 {
445  // Read next event in raw data file.
446  // Returns false if no event found (end of file).
447  // The list of all fired acquisition parameters is filled, and can be retrieved with
448  // GetFiredDataParameters().
449  // If SetUserTree(TTree*) has been called, the TTree is filled with the values of all
450  // parameters in this event.
451 
452  Bool_t ok = fGanilData->Next();
454  if (fUserTree) {
455  if (make_arrays) {
457  TIter next(fFired);
458  KVACQParam* par;
459  int i = 0;
460  while ((par = (KVACQParam*)next())) {
461  ParVal[i] = par->GetCoderData();
462  ParNum[i] = par->GetNumber();
463  i++;
464  }
465  }
466  fUserTree->Fill();
467  }
468  return ok;
469 }
470 
471 
472 
473 
478 
480 {
481  // Creates and returns new instance of class derived from GTGanilData used to read GANIL acquisition data
482  // Actual class is determined by Plugin.GTGanilData in environment files and name of dataset.
483 
484  //check and load plugin library
485  TPluginHandler* ph = LoadPlugin("GTGanilData", dataset);
486  if (!ph) {
487  // default: GTGanilData
488  return (new GTGanilData);
489  }
490 
491  //execute constructor/macro for plugin
492  return ((GTGanilData*) ph->ExecPlugin(0));
493 }
494 
495 
496 
497 
499 
501 {
502  return fGanilData;
503 }
504 
505 
506 
507 
511 
513 {
514  //Static method, used by KVDataSet::Open
515  //The name of the dataset is passed in the option string
516  return new KVGANILDataReader(filename, dataset);
517 }
518 
519 
520 
521 
524 
526 {
527  // clears and then fills list fFired with all fired acquisition parameters in event
528  fFired->Clear();
529  TIter next(fParameters);
530  KVACQParam* par;
531  while ((par = (KVACQParam*)next())) if (par->Fired()) fFired->Add(par);
532 }
533 
534 
535 
537 
539 {
540  return fGanilData->IsScalerBuffer();
541 }
542 
543 
544 
546 
548 {
549  return fGanilData->GetScalers()->GetNbChannel();
550 }
551 
552 
553 
554 
556 
558 {
559  return fGanilData->GetScalers()->GetScalerPtr(index)->GetCount();
560 }
561 
562 
563 
564 
566 
568 {
569  return fGanilData->GetScalers()->GetScalerPtr(index)->GetStatus();
570 }
571 
572 
573 
574 
576 
578 {
579  return fGanilData->GetEventCount();
580 }
581 
582 
583 
587 
589 {
590  // Return run number of file being read.
591  // Only call when file has been successfully opened.
592  return fGanilData->GetRunNumber();
593 }
594 
595 
596 
597 
int Int_t
unsigned int UInt_t
long Long_t
ClassImp(KVPartitionList) void KVPartitionList
Initialisation.
unsigned short UShort_t
char Char_t
const Bool_t kFALSE
bool Bool_t
const Bool_t kTRUE
const char Option_t
R__EXTERN TEnv * gEnv
#define gInterpreter
char * Form(const char *fmt,...)
R__EXTERN TSystem * gSystem
Classes handling VME/VXI parameters in GANIL DAQ data.
int Bits(void) const
int Index(void) const
Read GANIL formatted tapes or files.
Definition: GTGanilData.h:63
void Connect(const Int_t index, UShort_t **p) const
Connect a pointer to a data to the defined index in the Data Array.
Int_t GetRunNumber(void) const
Returns current run number.
const TList * GetListOfDataParameters() const
Definition: GTGanilData.h:112
void SetFileName(const TString filename)
Definition: GTGanilData.cpp:90
Bool_t IsOpen(void) const
bool IsScalerBuffer(void) const
Definition: GTGanilData.h:90
Int_t GetEventCount() const
Definition: GTGanilData.h:116
GTScalers * GetScalers(void)
Definition: GTGanilData.h:94
void SetScalerBuffersManagement(const ScalerWhat_t sc)
bool Next(void)
void Open(void)
ScalerWhat_t
What to do with scaler buffer.
Definition: GTGanilData.h:66
@ kAutoWriteScaler
have to take care of it.
Definition: GTGanilData.h:71
UInt_t GetStatus(void) const
Definition: GTOneScaler.h:45
UInt_t GetCount(void) const
Definition: GTOneScaler.h:41
const GTOneScaler * GetScalerPtr(Int_t index) const
Return a constant pointer to the Scaler Index.
Definition: GTScalers.cpp:98
Int_t GetNbChannel(void) const
Definition: GTScalers.h:41
GANIL VXI/VME acquisition parameter.
Definition: KVACQParam.h:15
void SetNbBits(UChar_t n)
Definition: KVACQParam.h:136
UShort_t ** ConnectData()
Definition: KVACQParam.h:51
Short_t GetCoderData() const
Definition: KVACQParam.h:64
UChar_t GetNbBits() const
Definition: KVACQParam.h:140
Bool_t Fired(Option_t *what="") const
Definition: KVACQParam.h:85
virtual void SetNumber(UInt_t num)
Definition: KVBase.h:209
static TPluginHandler * LoadPlugin(const Char_t *base, const Char_t *uri="0")
Definition: KVBase.cpp:756
UInt_t GetNumber() const
Definition: KVBase.h:213
Reads GANIL acquisition files (EBYEDAT)
virtual KVACQParam * CheckACQParam(const TSeqCollection *, const Char_t *)
Int_t GetEventCount() const
static KVGANILDataReader * Open(const Char_t *filename, Option_t *opt="")
void init()
default initialisations
virtual void SetUserTree(TTree *, Option_t *="arrays")
void ConnectRawDataParameters(const TSeqCollection *list_acq_params=nullptr)
Bool_t HasScalerBuffer() const
void FillFiredParameterList()
clears and then fills list fFired with all fired acquisition parameters in event
Int_t GetNumberOfScalers() const
Int_t GetScalerStatus(Int_t index) const
virtual GTGanilData * NewGanTapeInterface(Option_t *dataset)
virtual Bool_t GetNextEvent()
Int_t GetRunNumberReadFromFile() const
GTGanilData * fGanilData
object used to read GANIL acquisition file
virtual GTGanilData * GetGanTapeInterface()
const KVSeqCollection * GetRawDataParameters() const
void OpenFile(const Char_t *, Option_t *dataset)
TTree * fUserTree
user TTree to fill with data
virtual ~KVGANILDataReader()
Destructor.
UInt_t GetScalerCount(Int_t index) const
KVHashList * fFired
list of fired parameters in one event
KVHashList * fExtParams
list of data parameters in file not defined by gMultiDetArray
KVHashList * fParameters
list of all data parameters contained in file
Extended version of ROOT THashList.
Definition: KVHashList.h:28
virtual void SetOwner(Bool_t enable=kTRUE)
virtual void Clear(Option_t *option="")
virtual void SetCleanup(Bool_t enable=kTRUE)
virtual void Add(TObject *obj)
virtual void Delete(Option_t *option="")
virtual TObject * FindObject(const char *name) const
void SetName(const char *name)
virtual Int_t GetEntries() const
virtual const char * GetValue(const char *name, const char *dflt) const
static TFile * Open(const char *name, Option_t *option="", const char *ftitle="", Int_t compress=ROOT::RCompressionSetting::EDefaults::kUseCompiledDefault, Int_t netopt=0)
void Reset()
virtual void Add(TObject *obj)
virtual const char * GetName() const
virtual void SetTitle(const char *title="")
virtual void SetName(const char *name)
virtual void AddAt(TObject *obj, Int_t idx)
virtual void Error(const char *method, const char *msgfmt,...) const
void MakeZombie()
Longptr_t ExecPlugin(int nargs, const T &... params)
void ToUpper()
const char * Data() const
TString & Prepend(char c, Ssiz_t rep=1)
void Form(const char *fmt,...)
Bool_t Contains(const char *pat, ECaseCompare cmp=kExact) const
virtual const char * BaseName(const char *pathname)
virtual Int_t Fill()
virtual Bool_t SetAlias(const char *aliasName, const char *aliasFormula)
virtual TList * GetUserInfo()
virtual void SetAutoFlush(Long64_t autof=-30000000)
virtual Int_t Branch(const char *folder, Int_t bufsize=32000, Int_t splitlevel=99)
const char * GetUrl(Bool_t withDeflt=kFALSE) const
void SetOptions(const char *opt)
double T(double x)
TFile * OpenFile(const TString &fin)
void ws()