KaliVeda  1.13/01
Heavy-Ion Analysis Toolkit
KVINDRAPulserDataTree.cpp
Go to the documentation of this file.
1 /*
2 $Id: KVINDRAPulserDataTree.cpp,v 1.7 2009/03/27 16:42:58 franklan Exp $
3 $Revision: 1.7 $
4 $Date: 2009/03/27 16:42:58 $
5 */
6 
7 //Created by KVClassFactory on Wed Jan 21 11:56:26 2009
8 //Author: franklan
9 
10 #include "KVINDRAPulserDataTree.h"
11 #include "KVDataSet.h"
12 #include "KVINDRA.h"
13 #include "KVINDRADBRun.h"
14 
15 using namespace std;
16 
18 
19 
20 
21 
22 
26 {
27  // Default constructor
28  fArb = 0;
29  fRun = 0;
30  fVal = 0;
31  fIndex = 0;
32  fRunlist = 0;
33 }
34 
35 
36 
39 
41 {
42  // Destructor
43  if (fVal) delete [] fVal;
44  if (fIndex) {
45  fIndex->Delete();
46  delete fIndex;
47  }
48 }
49 
50 
51 
112 
114 {
115  // Create and fill tree with pulser data.
116  // We look for the following two directories:
117  //
118  // $KVROOT/KVFiles/name_of_dataset/gene_detecteurs
119  // $KVROOT/KVFiles/name_of_dataset/gene_pins
120  //
121  // If not found, we look for the following compressed archives:
122  //
123  // $KVROOT/KVFiles/name_of_dataset/gene_detecteurs.tgz
124  // $KVROOT/KVFiles/name_of_dataset/gene_pins.tgz
125  //
126  // and if found, uncompress them ('tar -zxf').
127  // [[ N.B. in this case the extracted directories will be deleted after reading,
128  // [[ ensuring that if new archives are supplied,
129  // [[ we always use the latest versions of files.
130  //
131  // The default names of these directories are defined in .kvrootrc by:
132  //
133  // KVINDRAPulserDataTree.GeneDetDir: gene_detecteurs
134  // KVINDRAPulserDataTree.GenePinDir: gene_pins
135  //
136  // Dataset-dependent alternatives can be defined using:
137  //
138  // dataset_name.KVINDRAPulserDataTree.GeneDetDir: dataset_specific_value
139  //
140  // The first directory (gene_detecteurs) must contain 1 file per run with names like:
141  //
142  // run8820.gene
143  //
144  // These files are generated using example analysis class GetGeneMean
145  // (see Examples). They contain the mean value of every acquisition parameter
146  // associated with a detector in the run:
147  //
148  // CI_0201_GG 3095.28
149  // CI_0201_PG 275.626
150  // CI_0203_GG 2863.66
151  // CI_0203_PG 263.308
152  // CI_0205_GG 3042.83
153  // etc.
154  //
155  // i.e. 'name of acquisition parameter' 'mean value for run'
156  //
157  // The second directory (gene_pins) must contain 1 file per run with names like:
158  //
159  // run8820.genepin
160  // OR run8820.laserpin
161  // OR run8820.genelaserpin
162  //
163  // These files are generated using example analysis class GetGeneMeanPin
164  // (see Examples). They contain the mean values of the acquisition parameters
165  // associated with pin diodes in the run:
166  //
167  // PILA_01_PG_gene 792.616
168  // PILA_01_PG_laser 1747.18
169  // PILA_01_GG_gene 0
170  // etc.
171  //
172  // We create a TTree with 1 branch for each detector acquisition parameter.
173  // For PILA and SI_PIN parameters, we create a 'gene' and a 'laser' branch for each.
174 
175  fGeneDir = new KVTarArchive(GetDirectoryName("GeneDetDir"), gDataSet->GetDataSetDir());
176  fPinDir = new KVTarArchive(GetDirectoryName("GenePinDir"), gDataSet->GetDataSetDir());
177  if (!fGeneDir->IsOK() && !fPinDir->IsOK()) {
178  Info("Build", "No data available to build pulser data tree");
179  return;
180  }
181  CreateTree();
182  ReadData();
183  delete fGeneDir;
184  delete fPinDir;
185 }
186 
187 
188 
195 
197 {
198  // Returns the name of the directory defined by the .kvrootrc environment variable
199  //
200  // KVINDRAPulserDataTree.[dirvar]
201  // OR
202  // dataset_name.KVINDRAPulserDataTree.[dirvar]:
203 
204  TString search, datasetenv;
205  datasetenv.Form("KVINDRAPulserDataTree.%s", dirvar);
206  search = gDataSet->GetDataSetEnv(datasetenv.Data(), "");
207  if (search == "") {
208  Error("GetDirectoryName", "%s is not defined for dataset %s. Check .kvrootrc files.",
209  datasetenv.Data(), gDataSet->GetName());
210  }
211  return search;
212 }
213 
214 
215 
223 
225 {
226  // Create new TTree with
227  // 1 branch 'Run' with run number
228  // 1 branch for each acquisition parameter of every detector (except time markers)
229  // 2 branches for each 'PILA_...' or 'SI_PIN...' parameter, suffixed with '_laser' and '_gene'
230  //
231  // NB if multidetector has not been built, it will be built by this method
232 
233  fArb = new TTree("PulserData", "Created by KVINDRAPulserDataTree");
234  fArb->SetDirectory(0);
235 
236  fArb->Branch("Run", &fRun, "Run/I");
237 
239 
240  const KVSeqCollection* acq_pars = nullptr;//gIndra->GetACQParams();
241 
242  fTab_siz = acq_pars->GetEntries() + 20;
243  fVal = new Float_t[fTab_siz];
244  fIndex = new THashTable(20, 5);
245  fIndex->SetOwner(kTRUE);
246 
247  TIter nxtACQ(acq_pars);
248  KVEBYEDAT_ACQParam* ap = 0;
249  Int_t ap_num = 0;
250  KVBase* iob;
251  while ((ap = (KVEBYEDAT_ACQParam*)nxtACQ())) {
252 
253  TString ap_name(ap->GetName());
254  TString ap_type(ap->GetType());
255  if (ap_name.BeginsWith("PILA") || ap_name.BeginsWith("SI_PIN")) {
256  ap_name += "_laser";
257  iob = new KVBase(ap_name.Data());
258  iob->SetNumber(ap_num);
259  fIndex->Add(iob);
260  fArb->Branch(ap_name.Data(), &fVal[ap_num++], Form("%s/F", ap_name.Data()));
261  ap_name.Form("%s%s", ap->GetName(), "_gene");
262  iob = new KVBase(ap_name.Data());
263  iob->SetNumber(ap_num);
264  fIndex->Add(iob);
265  fArb->Branch(ap_name.Data(), &fVal[ap_num++], Form("%s/F", ap_name.Data()));
266  }
267  else if (ap_type != "T") {
268  iob = new KVBase(ap_name.Data());
269  iob->SetNumber(ap_num);
270  fIndex->Add(iob);
271  fArb->Branch(ap_name.Data(), &fVal[ap_num++], Form("%s/F", ap_name.Data()));
272  }
273  if (ap_num > fTab_siz - 2) {
274  Error("CreateTree",
275  "Number of branches to create is greater than estimated (%d). Not all parameters can be treated.",
276  fTab_siz);
277  return;
278  }
279 
280  }
281  //keep number of used 'slots' in array
282  fTab_siz = ap_num;
283 }
284 
285 
286 
289 
291 {
292  // Read data in one file
293 
294  KVString line;
295  line.ReadLine(fin);
296  while (fin.good()) {
297  if (!line.BeginsWith("#")) {
298  line.Begin(" ");
299  KVString br_name = line.Next(kTRUE);
300  Int_t index = GetIndex(br_name.Data());
301  fVal[index] = line.Next(kTRUE).Atof();
302  }
303  line.ReadLine(fin);
304  }
305  fin.close();
306 }
307 
308 
309 
312 
314 {
315  // Read data for one run, fill tree
316 
317  UChar_t msg = 0;
318  fRun = run;
319  if (fGeneDir->IsOK()) {
320  ifstream f;
321  if (OpenGeneData(run, f)) ReadFile(f);
322  else msg = msg | 1;
323  }
324  if (fPinDir->IsOK()) {
325  ifstream f;
326  if (OpenPinData(run, f)) ReadFile(f);
327  else msg = msg | 2;
328  }
329  fArb->Fill();
330  return msg;
331 }
332 
333 
334 
339 
341 {
342  // Open gene data for one run
343  // We look for file 'runXXXX.gene' in the directory given by
344  // environment variable KVINDRAPulserDataTree.GeneDetDir.
345 
346  TString fname;
347  fname.Form("/run%d.gene", run);
348  fname.Prepend(gDataSet->GetDataSetEnv("KVINDRAPulserDataTree.GeneDetDir", ""));
349  return gDataSet->OpenDataSetFile(fname.Data(), f);
350 }
351 
352 
353 
362 
364 {
365  // Open pin data for one run
366  // We look for one of the following files in the directory given by
367  // environment variable KVINDRAPulserDataTree.GenePinDir:
368  //
369  // runXXXX.genepin
370  // OR runXXXX.laserpin
371  // OR runXXXX.genelaserpin
372 
373  TString fname;
374  fname.Form("/run%d.genepin", run);
375  TString pindir(gDataSet->GetDataSetEnv("KVINDRAPulserDataTree.GenePinDir", ""));
376  fname.Prepend(pindir);
377  if (gDataSet->OpenDataSetFile(fname.Data(), f)) return kTRUE;
378  fname.Form("/run%d.laserpin", run);
379  fname.Prepend(pindir);
380  if (gDataSet->OpenDataSetFile(fname.Data(), f)) return kTRUE;
381  fname.Form("/run%d.genelaserpin", run);
382  fname.Prepend(pindir);
383  return gDataSet->OpenDataSetFile(fname.Data(), f);
384 }
385 
386 
387 
390 
392 {
393  // Read data for every run in dataset
394 
395  if (!fRunlist) {
396  Error("ReadData", "Must set list of runs first using SetRunList(TList*)");
397  return;
398  }
399  Info("ReadData", "Reading pulser and laser data for all runs");
400  TIter Nxt_r(fRunlist);
401  KVINDRADBRun* run = 0;
402  KVNumberList missing1, missing2;
403  while ((run = (KVINDRADBRun*)Nxt_r())) {
404  Int_t run_num = run->GetNumber();
405  //reset all array members to -1
406  for (int i = 0; i < fTab_siz; i++) fVal[i] = -1.0;
407  std::cout << "\rInfo in <KVINDRAPulserDataTree::ReadData>: Reading data for run " << run_num << std::flush;
408  UChar_t msg = ReadData(run_num);
409  if (msg & 1) missing1.Add(run_num);
410  if (msg & 2) missing2.Add(run_num);
411  }
412  std::cout << std::endl;
413  if (missing1.GetEntries())
414  Warning("ReadData", "Missing file 'run[run_num].gene' for runs: %s", missing1.GetList());
415  if (missing2.GetEntries())
416  Warning("ReadData", "Missing file 'run[run_num].[gene][laser]pin' for runs: %s", missing2.GetList());
417 }
418 
419 
420 
423 
425 {
426  // Read pulser data tree from file
427 
428  fArb = (TTree*)file->Get("PulserData");
429  if (fArb) {
430  //disable all branches except Run number
431  fArb->SetBranchStatus("*", 0);
432  fArb->SetBranchStatus("Run", 1);
433  }
434 }
435 
436 
437 
441 
443 {
444  // Write pulser data tree in file
445  // We build and index based on the Run number and store it in the tree.
446 
447  if (fArb) {
448  fArb->SetDirectory(file);
449  fArb->BuildIndex("Run");
450  fArb->Write();
451  }
452 }
453 
454 
455 
465 
467 {
468  // Return mean value of pulser/laser for given parameter and run.
469  // For detectors, param should be name of an acquisition parameter
470  // e.g. CI_0201_PG, CSI_1301_L, etc.
471  // For pin laser diodes, param should be name of associated acquisition parameter
472  // with either '_laser' or '_gene' appended
473  // e.g. PILA_05_PG_laser, SI_PIN1_PG_gene
474  //
475  // Returns -1.0 if no data available for this parameter/run.
476 
477  if (!fArb) return -1.0;
478 
479  //find corresponding branch
480  TBranch* br = fArb->GetBranch(param);
481  if (!br) {
482  //no branch found - wrong name given ?
483  Error("GetMean", "No branch found with name %s", param);
484  return -1.0;
485  }
486  //enable branch
487  fArb->SetBranchStatus(param, 1);
488  //connect variable to branch
489  Float_t value = -1.0;
490  br->SetAddress(&value);
491  //read entry corresponding to run
492  Int_t bytes = fArb->GetEntryWithIndex(run);
493  if (bytes < 0) {
494  //unknown run number
495  Error("GetMean", "Unknown run %d", run);
496  return -1.0;
497  }
498  //disable branch
499  fArb->SetBranchStatus(param, 0);
500 
501  return value;
502 }
503 
504 
int Int_t
KVDataSet * gDataSet
Definition: KVDataSet.cpp:29
KVINDRA * gIndra
Definition: KVINDRA.cpp:88
ClassImp(KVPartitionList) void KVPartitionList
Initialisation.
#define f(i)
unsigned char UChar_t
char Char_t
bool Bool_t
float Float_t
const Bool_t kTRUE
char * Form(const char *fmt,...)
Base class for KaliVeda framework.
Definition: KVBase.h:141
virtual const Char_t * GetType() const
Definition: KVBase.h:176
virtual void SetNumber(UInt_t num)
Definition: KVBase.h:215
virtual Int_t GetNumber() const
Definition: KVDBRecord.h:72
const Char_t * GetDataSetDir() const
Definition: KVDataSet.cpp:720
const Char_t * GetDataSetEnv(const Char_t *type, const Char_t *defval="") const
Definition: KVDataSet.cpp:758
Bool_t OpenDataSetFile(const Char_t *filename, std::ifstream &file)
Definition: KVDataSet.cpp:1850
GANIL VXI/VME 16 bit (maximum) EBYEDAT acquisition parameter.
Database entry for each run of an INDRA experiment.
Definition: KVINDRADBRun.h:29
Handles TTree with mean pulser data for every run.
void ReadFile(std::ifstream &)
Read data in one file.
Bool_t OpenGeneData(Int_t, std::ifstream &)
virtual ~KVINDRAPulserDataTree()
Destructor.
Float_t GetMean(const Char_t *, Int_t)
TString GetDirectoryName(const Char_t *)
void ReadData()
Read data for every run in dataset.
void ReadTree(TFile *)
Read pulser data tree from file.
Bool_t OpenPinData(Int_t, std::ifstream &)
static KVMultiDetArray * MakeMultiDetector(const Char_t *dataset_name, Int_t run=-1, TString classname="KVMultiDetArray")
Strings used to represent a set of ranges of values.
Definition: KVNumberList.h:83
const Char_t * GetList() const
Int_t GetEntries() const
Definition: KVNumberList.h:169
void Add(Int_t)
Add value 'n' to the list.
KaliVeda extensions to ROOT collection classes.
Extension of ROOT TString class which allows backwards compatibility with ROOT v3....
Definition: KVString.h:72
Handles directories stored in .tgz archive files.
Definition: KVTarArchive.h:21
virtual void SetAddress(void *add)
virtual Int_t GetEntries() const
virtual const char * GetName() const
Bool_t BeginsWith(const char *s, ECaseCompare cmp=kExact) const
const char * Data() const
TString & Prepend(char c, Ssiz_t rep=1)
void Form(const char *fmt,...)
TLine * line
void Info(const char *location, const char *va_(fmt),...)
void Error(const char *location, const char *va_(fmt),...)
void Warning(const char *location, const char *va_(fmt),...)
size_t fIndex