KaliVeda  1.13/01
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 (ParVal) delete [] ParVal;
48  if (ParNum) delete [] ParNum;
49  if (fGanilData) delete fGanilData;
50 }
51 
52 
53 
56 
58 {
59  //default initialisations
60 
62  fUserTree = nullptr;
63  ParVal = nullptr;
64  ParNum = nullptr;
66  fGanilData = nullptr;
67 }
68 
69 
70 
137 
139 {
140  // To fill a TTree with the data in the current file, create a TTree:
141  // TFile* file = new TFile("run1.root","recreate");
142  // TTree* T = new TTree("Run1", "Raw data for Run1");
143  // and then call this method: SetUserTree(T)
144  // If you read all events of the file, the TTree will be automatically filled
145  // with data :
146  // while( runfile->GetNextEvent() ) ;
147  //
148  // Two different TTree structures are available, depending on the option string:
149  //
150  // opt = "arrays": [default]
151  //
152  // The TTree will have the following structure:
153  //
154  // *Br 0 :NbParFired : NbParFired/I = number of fired parameters in event
155  // *............................................................................*
156  // *Br 1 :ParNum : ParNum[NbParFired]/i = array of indices of fired parameters
157  // *............................................................................*
158  // *Br 2 :ParVal : ParVal[NbParFired]/s = array of values of fired parameters
159  //
160  // This structure is the fastest to fill and produces the smallest file sizes.
161  // In order to be able to directly access the parameters as if option "leaves" were used
162  // (i.e. one branch/leaf for each parameter), we add two aliases for each parameter to
163  // the tree:
164  // PARNAME = value of parameter if present in event
165  // PARNAME_M = number of times parameter appears in event
166  // Assuming that each parameter only appears at most once in each event, i.e. PARNAME_M=0 or 1,
167  // then
168  // root[0] T->Draw("PARNAME", "PARNAME_M")
169  // will histogram the value of PARNAME for each event in which it is present.
170  // (if the selection condition "PARNAME_M" is not used, the histogram will also be filled with a 0
171  // for each event in which PARNAME does not appear).
172  // N.B. the PARNAME alias is in fact the sum of the values of PARNAME in each event.
173  // If PARNAME_M>1 in some events, it is not the individual values but their sum which will
174  // be histogrammed in this case.
175  //
176  // Thus, if the data file has parameters called "PAR_1" and "PAR_2",
177  // the following command will work
178  //
179  // root[0] T->Draw("PAR_1:PAR_2", "PAR_1_M&&PAR_2_M", "col")
180  //
181  // even though no branches "PAR_1" or "PAR_2" exist.
182  //
183  //
184  //
185  // opt = "leaves":
186  //
187  // The TTree will have a branch/leaf for each parameter. This option is slower and produces
188  // larger files.
189  //
190  // If the option string contains both "arrays" and "leaves", then both structures will be used
191  // (in this case there is a high redundancy, as each parameter is stored twice).
192  //
193  // The full list of parameters is stored in a TObjArray in the list returned by TTree::GetUserInfo().
194  // Each parameter is represented by a TNamed object.
195  // In order to retrieve the name of the parameter with index 674 (e.g. taken from branch ParNum),
196  // do:
197  // TObjArray* parlist = (TObjArray*) T->GetUserInfo()->FindObject("ParameterList");
198  // cout << "Par 674 name = " << (*parlist)[674]->GetName() << endl;
199  //
200  //
201  // Automatic creation & filling of Scalers TTree
202  //
203  // give an option string containing "scalers", i.e. "leaves,scalers", or "ARRAYS+SCALERS", etc.
204  // a TTree with name 'Scalers' will be created, all scaler buffers will be written in it.
205 
206 
207  TString option = opt;
208  option.ToUpper();
209  make_arrays = option.Contains("ARRAYS");
210  make_leaves = option.Contains("LEAVES");
211  Bool_t make_scalers = option.Contains("SCALERS");
212  if (make_scalers) {
214  }
215 
216  fUserTree = T;
217  if (make_arrays) {
218  Int_t maxParFired = GetRawDataParameters().GetEntries();
219  ParVal = new UShort_t[maxParFired];
220  ParNum = new UInt_t[maxParFired];
221  fUserTree->Branch("NbParFired", &NbParFired, "NbParFired/I");
222  fUserTree->Branch("ParNum", ParNum, "ParNum[NbParFired]/i");
223  fUserTree->Branch("ParVal", ParVal, "ParVal[NbParFired]/s");
224  }
225  if (make_leaves) {
226  TIter next_rawpar(&GetRawDataParameters());
227  KVEBYEDAT_ACQParam* acqpar;
228  while ((acqpar = (KVEBYEDAT_ACQParam*)next_rawpar())) {
229  TString leaf;
230  leaf.Form("%s/S", acqpar->GetName());
231  // for parameters with <=8 bits only use 1 byte for storage
232  if (acqpar->GetNbBits() <= 8) leaf += "1";
233  fUserTree->Branch(acqpar->GetName(), *(acqpar->ConnectData()), leaf.Data());
234  }
235  }
236 
237 #if ROOT_VERSION_CODE > ROOT_VERSION(5,25,4)
238 #if ROOT_VERSION_CODE < ROOT_VERSION(5,26,1)
239  // The TTree::OptimizeBaskets mechanism is disabled, as for ROOT versions < 5.26/00b
240  // this lead to a memory leak
242 #endif
243 #endif
244 
245  // add list of parameter names in fUserTree->GetUserInfos()
246  // and if option="arrays" add aliases for each parameter & its multiplicity
247 
248  // TObjArray has to be as big as the largest parameter number in the list
249  // of raw data parameters. So first loop over parameters to find max param number.
250  UInt_t maxpar = 0;
251  TIter next(&GetRawDataParameters());
252  KVEBYEDAT_ACQParam* par;
253  while ((par = (KVEBYEDAT_ACQParam*)next())) if (par->GetNumber() > maxpar) maxpar = par->GetNumber();
254 
255  TObjArray* parlist = new TObjArray(maxpar, 1);
256  parlist->SetName("ParameterList");
257 
258  next.Reset();
259  while ((par = (KVEBYEDAT_ACQParam*)next())) {
260  parlist->AddAt(new TNamed(par->GetName(), Form("index=%d", par->GetNumber())), par->GetNumber());
261  if (make_arrays) {
262  fUserTree->SetAlias(par->GetName(), Form("Sum$((ParNum==%d)*ParVal)", par->GetNumber()));
263  fUserTree->SetAlias(Form("%s_M", par->GetName()), Form("Sum$(ParNum==%d)", par->GetNumber()));
264  }
265  }
266  fUserTree->GetUserInfo()->Add(parlist);
267 }
268 
269 
270 
282 
284 {
285  //Open and initialise a GANIL data file for reading.
286  //By default, scaler buffers are ignored.
287  //This can be changed by changing the value of
288  //
289  // KVGANILDataReader.ScalerBuffersManagement: kSkipScaler
290  //
291  //If file cannot be opened, this object will be made a zombie. Do not use.
292  //To test if file is open, use IsOpen().
293  //The basename of the file (excluding any path) can be obtained using GetName()
294  //The full pathname of the file can be obtained using GetTitle()
295 
296  if (fGanilData) delete fGanilData;
297  fGanilData = NewGanTapeInterface(dataset);
298 
301  SetTitle(file);
302 
303  // handling scaler buffers
304  TString what_scale = gEnv->GetValue("KVGANILDataReader.ScalerBuffersManagement", "kSkipScaler");
305  what_scale.Prepend("GTGanilData::");
306  Long_t ws = gInterpreter->ProcessLine(what_scale);
308 
309  fGanilData->Open();
310 
311  //test whether file has been opened
312  if (!fGanilData->IsOpen()) {
313  //if initial attempt fails, we try to open the file as a 'raw' TFile
314  //This may work when we are attempting to open a remote file i.e.
315  //via rfio, and a subsequent attempt to open the file using the GanTape
316  //routines may then succeed.
317  TUrl rawtfile(file, kTRUE);
318  rawtfile.SetOptions("filetype=raw");
319  TFile* rawfile = TFile::Open(rawtfile.GetUrl());
320  if (!rawfile) {
321  Error("OpenFile", "File cannot be opened: %s",
322  file);
323  MakeZombie();
324  return;
325  }
326  //TFile::Open managed to open file! Try again...
327  delete rawfile;
328  fGanilData->Open();
329  if (!fGanilData->IsOpen()) {
330  //failed again ??!
331  Error("OpenFile", "File cannot be opened: %s",
332  file);
333  MakeZombie();
334  return;
335  }
336  }
337 }
338 
339 
340 
341 
349 
351 {
352  //Generate all required KVEBYEDAT_ACQParam objects corresponding to parameters in file
353  //
354  //fParameters is filled with a KVEBYEDAT_ACQParam for every acquisition parameter in the file.
355  //
356  //To access the full list of data parameters in the file after this method has been
357  //called (i.e. after the file is opened), use GetRawDataParameters().
358 
360  KVEBYEDAT_ACQParam* par;
361  GTDataPar* daq_par;
362  while ((daq_par = (GTDataPar*) next())) {//loop over all parameters
363  par = new KVEBYEDAT_ACQParam(daq_par->GetName());
364  par->SetNumber(daq_par->Index());
365  par->SetNbBits(daq_par->Bits());
366  fParameters.Add(par);
367  fGanilData->Connect(daq_par->Index(), par->ConnectData());
368  fGanilData->ConnectFired(daq_par->Index(), par->ConnectFired());
369  }
370 }
371 
372 
373 
374 
382 
384 {
385  // Read next event in raw data file.
386  // Returns false if no event found (end of file).
387  // The list of all fired acquisition parameters is filled, and can be retrieved with
388  // GetFiredDataParameters().
389  // If SetUserTree(TTree*) has been called, the TTree is filled with the values of all
390  // parameters in this event.
391 
392  Bool_t ok = fGanilData->Next();
394  if (fUserTree) {
395  if (make_arrays) {
397  TIter next(&fFired);
398  KVEBYEDAT_ACQParam* par;
399  int i = 0;
400  while ((par = (KVEBYEDAT_ACQParam*)next())) {
401  ParVal[i] = par->GetData();
402  ParNum[i] = par->GetNumber();
403  i++;
404  }
405  }
406  fUserTree->Fill();
407  }
408  return ok;
409 }
410 
411 
412 
413 
418 
420 {
421  // Creates and returns new instance of class derived from GTGanilData used to read GANIL acquisition data
422  // Actual class is determined by Plugin.GTGanilData in environment files and name of dataset.
423 
424  //check and load plugin library
425  TPluginHandler* ph = LoadPlugin("GTGanilData", dataset);
426  if (!ph) {
427  // default: GTGanilData
428  return (new GTGanilData);
429  }
430 
431  //execute constructor/macro for plugin
432  return ((GTGanilData*) ph->ExecPlugin(0));
433 }
434 
435 
436 
437 
439 
441 {
442  return fGanilData;
443 }
444 
445 
446 
447 
451 
453 {
454  //Static method, used by KVDataSet::Open
455  //The name of the dataset is passed in the option string
456  return new KVGANILDataReader(filename, dataset);
457 }
458 
459 
460 
461 
464 
466 {
467  // clears and then fills list fFired with all fired acquisition parameters in event
468  fFired.Clear();
469  TIter next(&fParameters);
470  KVEBYEDAT_ACQParam* par;
471  while ((par = (KVEBYEDAT_ACQParam*)next())) if (par->Fired()) fFired.Add(par);
472 }
473 
474 
475 
477 
479 {
480  return fGanilData->IsScalerBuffer();
481 }
482 
483 
484 
486 
488 {
489  return fGanilData->GetScalers()->GetNbChannel();
490 }
491 
492 
493 
494 
496 
498 {
499  return fGanilData->GetScalers()->GetScalerPtr(index)->GetCount();
500 }
501 
502 
503 
504 
506 
508 {
509  return fGanilData->GetScalers()->GetScalerPtr(index)->GetStatus();
510 }
511 
512 
513 
514 
516 
518 {
519  return fGanilData->GetEventCount();
520 }
521 
522 
523 
527 
529 {
530  // Return run number of file being read.
531  // Only call when file has been successfully opened.
532  return fGanilData->GetRunNumber();
533 }
534 
535 
536 
537 
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 ConnectFired(const Int_t index, Bool_t **p) const
Connect a pointer to a data to the defined index in the Data Array.
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
virtual void SetNumber(UInt_t num)
Definition: KVBase.h:215
static TPluginHandler * LoadPlugin(const Char_t *base, const Char_t *uri="0")
Definition: KVBase.cpp:793
UInt_t GetNumber() const
Definition: KVBase.h:219
GANIL VXI/VME 16 bit (maximum) EBYEDAT acquisition parameter.
Bool_t Fired() const
UShort_t GetData() const
UShort_t ** ConnectData()
UChar_t GetNbBits() const
void SetNbBits(UChar_t n)
Reads GANIL acquisition files (EBYEDAT)
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")
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
const KVSeqCollection & GetRawDataParameters() const
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()
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
KVUniqueNameList fFired
list of fired parameters in one event
KVUniqueNameList fParameters
list of all data parameters contained in file
virtual void SetOwner(Bool_t enable=kTRUE)
virtual void Clear(Option_t *option="")
virtual void Add(TObject *obj)
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()