KaliVeda  1.12/06
Heavy-Ion Analysis Toolkit
KVINDRADstToRootTransfert.cpp
Go to the documentation of this file.
1 //Author: Eric Bonnet
2 
4 #include "KVClassFactory.h"
5 #include "KVBatchSystem.h"
6 
7 #include "TROOT.h"
8 #include "TSystem.h"
9 #include "TObjArray.h"
10 #include "Riostream.h"
11 #include "TFile.h"
12 #include "TTree.h"
13 #include "TSQLResult.h"
14 #include "TSQLRow.h"
15 #include "KVINDRA.h"
16 #include "KVINDRACodes.h"
17 #include "KVINDRAReconNuc.h"
18 #include "KVINDRAReconEvent.h"
19 #include "KVChIo.h"
20 #include "KVSilicon.h"
21 #include "KVCsI.h"
22 #include "KVPhoswich.h"
23 #include "KVDetector.h"
24 #include "KVIDTelescope.h"
25 #include "KVDBRun.h"
27 #include "KVDataSetManager.h"
28 #include "KVDataRepository.h"
29 #include "KVBatchSystem.h"
30 #include "KVGANILDataReader.h"
31 #include "KVINDRADB.h"
32 
34 
35 using namespace std;
36 
38 
39 
40 
41 
45 {
46  //Default constructor
47 }
48 
49 
50 
53 
55 {
56  //Destructor
57 }
58 
59 
60 
62 
64 {
65  Info("InitRun", "ds InitRun");
66 
67  KVString dst_file = gDataSet->GetFullPathToRunfile("dst", fRunNumber);
68  Info("InitRun", "dst file %s", dst_file.Data());
69 
70  Info("InitRun", "%s %s", gDataSet->GetDataPathSubdir(), gDataSet->GetRunfileName("dst", fRunNumber));
71  gDataRepository->CopyFileFromRepository(gDataSet, "dst", gDataSet->GetRunfileName("dst", fRunNumber), ".");
72  gROOT->ProcessLine(".! ls -lhrt *");
73 
74  Info("InitRun", "Starting analysis of run %d", fRunNumber);
75 
76  TDatime now1;
77  Info("InitRun", "Debut lecture DST %s", now1.AsString());
78  DefineSHELLVariables();
79 
80 
81  ReadDST();
82 
83  Info("InitRun", "Bilan ressource apres ReadDST");
84  ProcInfo_t pid;
85  if (gSystem->GetProcInfo(&pid) == 0) {
86  TString du = gSystem->GetFromPipe("du -hs");
87  TObjArray* toks = du.Tokenize(" .");
88  TString disk = ((TObjString*)toks->At(0))->String();
89  delete toks;
90  cout << " ------------- Process infos -------------" << endl;
91  printf(" CpuUser = %f s. VirtMem = %f MB DiskUsed = %s\n",
92  pid.fCpuUser, pid.fMemVirtual / 1024., disk.Data());
93  }
94  else {
95  Warning("InitRun", "pas d info disponible sur le bilan ressource");
96  }
97 
98 
99  TDatime now2;
100  Info("InitRun", "Fin lecture DST %s", now2.AsString());
101 
102  KVString inst;
103  inst.Form(".! rm %s", gDataSet->GetRunfileName("dst", fRunNumber));
104  gROOT->ProcessLine(inst.Data());
105 
106  inst.Form(".! rm xrun%d", fRunNumber);
107  gROOT->ProcessLine(inst.Data());
108 
109 }
110 
111 
112 
114 
116 {
117 
118  KVString inst;
119  inst.Form(".! %s", KVBase::GetBINDIRFilePath("veda2root_kali"));
120  gROOT->ProcessLine(inst.Data());
121 
122 }
123 
124 
125 
137 
139 {
140 
141  // Convert DST file for 1 run to KaliVeda ROOT format.
142  // By default the new ROOT file will be written in the same data repository as the DST file we are reading.
143  // This can be changed by setting the environment variable(s):
144  //
145  // DSTtoROOT.DataAnalysisTask.OutputRepository: [name of repository]
146  // [name of dataset].DSTtoROOT.DataAnalysisTask.OutputRepository: [name of repository]
147  //
148  // If no value is set for the current dataset (second variable), the value of the
149  // first variable will be used. If neither is defined, the new file will be written in the same repository as
150  // the DST file (if possible, i.e. if repository is not remote).
151 
152  ifstream f_in("list_of_files");
153  KVString line = "";
154 
155  Int_t nfiles, necrit, rn_verif;
156  nfiles = necrit = rn_verif = 0;
157  TObjArray* toks = 0;
158  while (f_in.good()) {
159  line.ReadLine(f_in);
160  if (line != "") {
161  toks = line.Tokenize("=");
162  Int_t val = ((TObjString*)toks->At(1))->GetString().Atoi();
163 
164  if (line.BeginsWith("file_number=")) nfiles = val;
165  else if (line.BeginsWith("write_evts=")) necrit = val;
166  else if (line.BeginsWith("run_number=")) rn_verif = val;
167  else {}
168  delete toks;
169  }
170  else {}
171  }
172  f_in.close();
173 
174  if (rn_verif != fRunNumber) {
175  Error("ProcessRun", "Le numero de run inscrit dans le fichier list_of_files (%d) est different du run courant (%d)", rn_verif, fRunNumber);
176  Error("ProcessRun", "\t->analysis is going to stop");
177 
178  TString inst;
179  inst.Form(".! rm list_of_files");
180  gROOT->ProcessLine(inst.Data());
181 
182  inst.Form(".! rm arbre_root_*.txt");
183  gROOT->ProcessLine(inst.Data());
184 
185  return;
186  }
187  Info("ProcessRun", "Lecture de list_of_files :\n - nbre de fichiers a lire %d\n - nbre d evts a lire %d\n", nfiles, necrit);
188 
189  TDatime now1;
190  Info("ProcessRun", "Debut Conversion au format du run %d ROOT %s", fRunNumber, now1.AsString());
191 
193 
197 
198  // Create new ROOT file with TTree for converted events.
199 
200  // get dataset to which we must associate new run
201  KVDataSet* OutputDataset =
203 
204  TFile* fi = OutputDataset->NewRunfile("root", fRunNumber);
205 
206  KVString stit;
207  stit.Form("%s : %s : root events converted from DST",
208  gIndraDB->GetRun(fRunNumber)->GetName(), gIndraDB->GetRun(fRunNumber)->GetTitle());
209  //tree for reconstructed events
210  data_tree = new TTree("ReconstructedEvents", stit.Data());
211 #if ROOT_VERSION_CODE > ROOT_VERSION(5,25,4)
212 #if ROOT_VERSION_CODE < ROOT_VERSION(5,26,1)
213  // The TTree::OptimizeBaskets mechanism is disabled, as for ROOT versions < 5.26/00b
214  // this lead to a memory leak
215  data_tree->SetAutoFlush(0);
216 #endif
217 #endif
218  //leaves for reconstructed events
219  KVEvent::MakeEventBranch(data_tree, "INDRAReconEvent", "KVINDRAReconEvent", evt);
220 
221  //tree for raw data
222  rawtree = new TTree("RawData", Form("%s : %s : raw data",
223  gIndraDB->GetRun(fRunNumber)->GetName(), gIndraDB->GetRun(fRunNumber)->GetTitle()));
224  rawtree->Branch("RunNumber", &fRunNumber, "RunNumber/I");
225  rawtree->Branch("EventNumber", &fEventNumber, "EventNumber/I");
226 
227  // the format of the raw data tree must be the same as that generated by KVGanilDataReader::SetUserTree
228  // with the option "arrays".
229  // we depend on it in KVINDRAReconDataAnalyser in order to read the raw data and set the detector acquisition parameters.
230  // the event numbers of the two trees must be synchronised
231  TString raw_opt = "arrays";
232 
233  KVGANILDataReader* raw_data = gDataSet->OpenRunfile<KVGANILDataReader>("raw", fRunNumber);
234  if (!raw_data) {
235  Error("ProcessRun", "No RAW data file available for this run. ABORT!!");
236  exit(1);
237  }
238  raw_data->SetUserTree(rawtree, raw_opt.Data());
239  Info("InitRun", "Created raw data tree (%s : %s). Format: %s",
240  rawtree->GetName(), rawtree->GetTitle(), raw_opt.Data());
241  // we use the KVGANILDataReader just in order to create the required TTree branches
242  // and initialize the acquisition parameter objects with the index given by the acquisition
243  // we can now delete it
244  delete raw_data;
245  // get list of all acquisition parameters
246  params = gIndra->GetACQParams();
247  // internal data member pointers are pointing to internal array of GTGanilData object
248  // used by now deleted KVGANILDataReader. need to set them to internal fData member.
249  TIter nextparam(params);
250  KVACQParam* acqpar;
251  while ((acqpar = (KVACQParam*)nextparam())) acqpar->UseInternalDataMember();
252  // reset addresses of branches which were pointing to internal members of deleted KVGANILDataReader.
253  rawtree->SetBranchAddress("NbParFired", &NbParFired);
254  rawtree->SetBranchAddress("ParVal", ParVal);
255  rawtree->SetBranchAddress("ParNum", ParNum);
256 
257  events_good = events_read = 0;
258 
259  if (camp2) {
260  Info("ProcessRun", "Special treatment INDRA 2nd campaign data:");
261  Info("ProcessRun", "Particles detected in phoswich ring 1 have general ID code=4");
262  Info("ProcessRun", "After translation, they will have Veda ID code=2 (like 1st campaign)");
263  }
264 
265  Info("ProcessRun", "Bilan ressource avant lecture des %d fichiers ascii", nfiles);
266  ProcInfo_t pid;
267  if (gSystem->GetProcInfo(&pid) == 0) {
268  TString du = gSystem->GetFromPipe("du -hs");
269  TObjArray* toks = du.Tokenize(" .");
270  TString disk = ((TObjString*)toks->At(0))->String();
271  delete toks;
272  cout << " ------------- Process infos -------------" << endl;
273  printf(" CpuUser = %f s. VirtMem = %f MB DiskUsed = %s\n",
274  pid.fCpuUser, pid.fMemVirtual / 1024., disk.Data());
275  }
276 
277 
278  KVString inst;
279  fEventNumber = 1;
280  for (Int_t nf = 1; nf <= nfiles; nf += 1) {
281 
282  stit.Form("arbre_root_%d.txt", nf);
283  ifstream f_data(stit.Data());
284  events_in_file = 0;
285  while (f_data.good()) {
286  lire_evt(f_data, evt);
287  evt->Clear();
288 
289  if (events_read % 10000 == 0 && events_read > 0) {
290  cout << " +++ " << events_read << " events processed +++ " << endl;
291  ProcInfo_t pid;
292  if (gSystem->GetProcInfo(&pid) == 0) {
293  TString du = gSystem->GetFromPipe("du -hs");
294  TObjArray* toks = du.Tokenize(" .");
295  TString disk = ((TObjString*)toks->At(0))->String();
296  delete toks;
297  cout << " ------------- Process infos -------------" << endl;
298  printf(" CpuUser = %f s. VirtMem = %f MB DiskUsed = %s\n",
299  pid.fCpuUser, pid.fMemVirtual / 1024., disk.Data());
300  }
301  }
302  }
303  f_data.close();
304 
305  inst.Form(".! rm arbre_root_%d.txt", nf);
306 
307  Info("ProcessRun", "Bilan ressource apres lecture du fichier numero %d/%d", nf, nfiles);
308  ProcInfo_t pid;
309  if (gSystem->GetProcInfo(&pid) == 0) {
310  TString du = gSystem->GetFromPipe("du -hs");
311  TObjArray* toks = du.Tokenize(" .");
312  TString disk = ((TObjString*)toks->At(0))->String();
313  delete toks;
314  cout << " ------------- Process infos -------------" << endl;
315  printf(" CpuUser = %f s. VirtMem = %f MB DiskUsed = %s\n",
316  pid.fCpuUser, pid.fMemVirtual / 1024., disk.Data());
317  }
318 
319  gROOT->ProcessLine(inst.Data());
320  }
321 
323  //check number of events against scaler info for this run
324  Double_t check_events = (1.*events_read) / ((Double_t)necrit);
325  cout << "Nb. evts read from files: " << events_read << endl;
326  cout << "Nb. evts written to TTree: " << events_good << endl;
327  cout << "Nb. evts in FORTRAN file: " << necrit << endl;
328  if (check_events < 0.99) cout << "ERROR: number of events read <99% number in FORTRAN file" << endl;
329  }
330 
331  //At the end of each run we:
332  // write the tree into the new file
333  // add the file to the repository
334 
335  fi->cd();
336 
337  gDataAnalyser->WriteBatchInfo(data_tree);
338  data_tree->Write(); //write tree to file
339  rawtree->Write();
340 
341  //add new file to repository
342  OutputDataset->CommitRunfile("root", fRunNumber, fi);
343 
344  TDatime now2;
345  Info("ProcessRun", "Fin Conversion format ROOT %s", now2.AsString());
346 
347 
348  inst.Form(".! rm list_of_files");
349  gROOT->ProcessLine(inst.Data());
350 
351 }
352 
353 
354 
356 
358 {
359 
360  Info("EndRun", "ds EndRun");
361 
362 
363 }
364 
365 
366 
368 
370 {
371 
372  if (gDataSet != GetDataSet()) GetDataSet()->cd();
373  SetCampagneNumber();
374 
375  if (fCampNumber == -1) return;
376 
377  //InitAnalysis();
378  Info("SubmitTask", "Liste de runs : %s", GetRunList().AsString());
379  //loop over runs
380  Info("SubmitTask", "RunningInLaunchDirectory : %d", Int_t(RunningInLaunchDirectory()));
381 
382  GetRunList().Begin();
383  while (!GetRunList().End()) {
384  Info("SubmitTask", "%s", GetBatchName());
385  fRunNumber = GetRunList().Next();
386  Info("SubmitTask", "Traitement du Run %d", fRunNumber);
387  InitRun();
388  ProcessRun();
389  EndRun();
390  }
391 
392  //EndAnalysis();
393 
394 }
395 
396 
397 
399 
401 {
402 
403 
404  Info("DefineSHELLVariables", "campagne %d", fCampNumber);
405 
406  KVString shell_var;
407 
408  shell_var.Form("%d", fCampNumber);
409  gSystem->Setenv("CAMPAGNE", shell_var.Data());
410 
411 // shell_var.Form("%s",gSystem->ExpandPathName("$PWD"));
412 // gSystem->Setenv("DIR_PERSO",shell_var.Data());
413 
414  shell_var.Form("%s/lib", gSystem->ExpandPathName("$KVROOT"));
415  gSystem->Setenv("REP_LIBS", shell_var.Data());
416 
417  shell_var.Form("%s/src/faire_arbre_c%d.f", gSystem->ExpandPathName("$KVROOT"), fCampNumber);
418  gSystem->Setenv("FOR_PERSO", shell_var.Data());
419 
420  shell_var.Form("run%d", fRunNumber);
421  gSystem->Setenv("RUN_PREFIX", shell_var.Data());
422 
423  //gSystem->Setenv("DST_FILE",gDataSet->GetRunfileName("dst",fRunNumber));
424 
425  camp1 = camp2 = camp4 = kFALSE;
426 
427  if (fCampNumber == 1) {
428 
429  shell_var.Form("%s/veda%d/for", gSystem->ExpandPathName("$THRONG_DIR"), fCampNumber);
430  gSystem->Setenv("VEDA_FOR", shell_var.Data());
431 
432  shell_var.Form("%s/veda%d/files/", gSystem->ExpandPathName("$THRONG_DIR"), fCampNumber);
433  gSystem->Setenv("VEDA_DATA", shell_var.Data());
434 
435  camp1 = kTRUE;
436 
437  }
438  else if (fCampNumber == 2) {
439 
440  shell_var.Form("%s/veda%d/for%d", gSystem->ExpandPathName("$THRONG_DIR"), fCampNumber, fCampNumber);
441  gSystem->Setenv("VEDA_FOR", shell_var.Data());
442 
443  shell_var.Form("%s/veda%d/data/", gSystem->ExpandPathName("$THRONG_DIR"), fCampNumber);
444  gSystem->Setenv("VEDA_DATA", shell_var.Data());
445 
446  camp2 = kTRUE;
447 
448  }
449  else if (fCampNumber == 4) {
450 
451  shell_var.Form("%s/veda%d/for/include", gSystem->ExpandPathName("$THRONG_DIR"), fCampNumber);
452  gSystem->Setenv("VEDA_FOR", shell_var.Data());
453 
454  shell_var.Form("%s/veda%d/data/", gSystem->ExpandPathName("$THRONG_DIR"), fCampNumber);
455  gSystem->Setenv("VEDA_DATA", shell_var.Data());
456 
457  camp4 = kTRUE;
458 
459  }
460  else {}
461 
462 
463 
464 }
465 
466 
478 
480 {
481 //code = 2 : part. ident. dans CsI
482  //for 1st campaign: code=2 for phoswich (ring 1) also
483 //code = 9 : ident. incomplete dans CsI ou Phoswich (Z.min)
484 //code =10 : ident. "entre les lignes" dans CsI
485 //set coder values of detectors concerned by a particle with ID code 2
486 //Les codes 9 et 10 sont des codes qui ont ete ajoutes pour donner
487 //au code 2 la notion d'identification sans ambiguites.
488 
489  //returns pointer of the detector for which
490  //KVINDRAReconNuc::Reconstruct must be called.
491  //cout << "Code2and9and10" << endl;
492  if (ring == 1) {
493  if (!camp4) {
494  //phoswich
495  KVPhoswich* phos =
496  (KVPhoswich*)gIndra->GetDetectorByType(ring, mod, Phos_R);
497  if (!phos) return 0;
498  phos->SetEnergy(de1);
499  phos->GetACQParam("R")->SetData((UShort_t)canal[Phos_R]);
500  phos->GetACQParam("L")->SetData((UShort_t)canal[Phos_L]);
501  phos->GetACQParam("T")->SetData((UShort_t)mt);
502  identifying_telescope = gIndra->GetIDTelescope(Form("PHOS_R_L_%02d", mod));
503  return phos;
504  }
505  else { // campagne 4: GSI Campaign
506  KVCsI* csi = (KVCsI*) gIndra->GetDetectorByType(ring, mod, CsI_R);
507  if (!csi) return 0;
508  KVSilicon* si = (KVSilicon*) gIndra->GetDetectorByType(ring, mod, Si_GG);
509  if (si)si->SetEnergy(de2);
510  csi->SetEnergy(de3);
511  if (si) {
512  si->GetACQParam("GG")->SetData((UShort_t)canal[Si_GG]);
513  si->GetACQParam("PG")->SetData((UShort_t)canal[Si_PG]);
514  si->GetACQParam("T")->SetData((UShort_t)canal[Si_T]);
515  }
516  csi->GetACQParam("R")->SetData((UShort_t)canal[CsI_R]);
517  csi->GetACQParam("L")->SetData((UShort_t)canal[CsI_L]);
518  csi->GetACQParam("T")->SetData((UShort_t)mt);
519  identifying_telescope = gIndra->GetIDTelescope(Form("CSI_R_L_%02d%02d", ring, mod));
520  return csi;
521  }
522  }
523  else if (ring >= 2 && ring <= 9) {
524  //chio-si-csi
525  KVCsI* csi = (KVCsI*) gIndra->GetDetectorByType(ring, mod, CsI_R);
526  if (!csi) return 0;
527  KVSilicon* si = (KVSilicon*) gIndra->GetDetectorByType(ring, mod, Si_GG);
528  KVChIo* chio = (KVChIo*)csi->GetChIo();
529  if (!chio) return 0;
530  chio->SetEnergy(de1);
531  if (si)si->SetEnergy(de2);
532  csi->SetEnergy(de3);
533  chio->GetACQParam("GG")->SetData((UShort_t)canal[ChIo_GG]);
534  chio->GetACQParam("PG")->SetData((UShort_t)canal[ChIo_PG]);
535  chio->GetACQParam("T")->SetData((UShort_t)canal[ChIo_T]);
536  if (si) {
537  si->GetACQParam("GG")->SetData((UShort_t)canal[Si_GG]);
538  si->GetACQParam("PG")->SetData((UShort_t)canal[Si_PG]);
539  si->GetACQParam("T")->SetData((UShort_t)canal[Si_T]);
540  }
541  csi->GetACQParam("R")->SetData((UShort_t)canal[CsI_R]);
542  csi->GetACQParam("L")->SetData((UShort_t)canal[CsI_L]);
543  csi->GetACQParam("T")->SetData((UShort_t)mt);
544  identifying_telescope = gIndra->GetIDTelescope(Form("CSI_R_L_%02d%02d", ring, mod));
545  return csi;
546  }
547  else if (ring > 9) {
548  //chio-csi or chio-si75-sili-csi
549  KVCsI* csi = (KVCsI*) gIndra->GetDetectorByType(ring, mod, CsI_R);
550  if (!csi) return 0;
551  KVChIo* chio = (KVChIo*)csi->GetChIo();
552  if (!chio) return 0;
553  KVSi75* si75 = 0;
554  KVSiLi* sili = 0;
555  if (canal[Si75_GG] > 0) {
556  si75 = (KVSi75*)gIndra->GetDetectorByType(ring, mod, Si75_GG);
557  }
558  if (canal[SiLi_GG] > 0) {
559  sili = (KVSiLi*)gIndra->GetDetectorByType(ring, mod, SiLi_GG);
560  }
561  chio->SetEnergy(de1);
562  if (si75)si75->SetEnergy(de4);
563  if (sili)sili->SetEnergy(de5);
564  csi->SetEnergy(de3);
565  chio->GetACQParam("GG")->SetData((UShort_t)canal[ChIo_GG]);
566  chio->GetACQParam("PG")->SetData((UShort_t)canal[ChIo_PG]);
567  chio->GetACQParam("T")->SetData((UShort_t)canal[ChIo_T]);
568  if (si75) {
569  si75->GetACQParam("GG")->SetData((UShort_t)canal[Si75_GG]);
570  si75->GetACQParam("PG")->SetData((UShort_t)canal[Si75_PG]);
571  si75->GetACQParam("T")->SetData((UShort_t)canal[Si75_T]);
572  }
573  if (sili) {
574  sili->GetACQParam("GG")->SetData((UShort_t)canal[SiLi_GG]);
575  sili->GetACQParam("PG")->SetData((UShort_t)canal[SiLi_PG]);
576  sili->GetACQParam("T")->SetData((UShort_t)canal[SiLi_T]);
577  }
578  csi->GetACQParam("R")->SetData((UShort_t)canal[CsI_R]);
579  csi->GetACQParam("L")->SetData((UShort_t)canal[CsI_L]);
580  csi->GetACQParam("T")->SetData((UShort_t)mt);
581  identifying_telescope = gIndra->GetIDTelescope(Form("CSI_R_L_%02d%02d", ring, mod));
582  return csi;
583  }
584  return 0;
585 }
586 
587 
596 
598 {
599 // cout << "Code3" << endl;
600 //code = 3 : fragment identifie dans Si-CsI
601 // ou fragment ident. dans Si75-SiLi ou SiLi-CsI
602  //set coder values of detectors concerned by a particle with ID code 3
603  //returns pointer of the detector for which
604  //KVINDRAReconNuc::Reconstruct must be called.
605 
606  // First Ring if GSI Campaign
607  if (ring == 1 && camp4) {
608  KVCsI* csi = (KVCsI*) gIndra->GetDetectorByType(ring, mod, CsI_R);
609  if (!csi) return 0;
610  KVSilicon* si = (KVSilicon*) gIndra->GetDetectorByType(ring, mod, Si_GG);
611  if (!si) return 0;
612  si->SetEnergy(de2);
613  csi->SetEnergy(de3);
614  si->GetACQParam("GG")->SetData((UShort_t)canal[Si_GG]);
615  si->GetACQParam("PG")->SetData((UShort_t)canal[Si_PG]);
616  si->GetACQParam("T")->SetData((UShort_t)canal[Si_T]);
617  csi->GetACQParam("R")->SetData((UShort_t)canal[CsI_R]);
618  csi->GetACQParam("L")->SetData((UShort_t)canal[CsI_L]);
619  csi->GetACQParam("T")->SetData((UShort_t)mt);
620  identifying_telescope = gIndra->GetIDTelescope(Form("SI_CSI_%02d%02d", ring, mod));
621  return csi;
622  }
623  else if (ring >= 2 && ring <= 9) {
624  //chio-si-csi
625  KVCsI* csi = (KVCsI*) gIndra->GetDetectorByType(ring, mod, CsI_R);
626  if (!csi) return 0;
627  KVSilicon* si = (KVSilicon*) gIndra->GetDetectorByType(ring, mod, Si_GG);
628  if (!si) return 0;
629  KVChIo* chio = (KVChIo*)csi->GetChIo();
630  if (!chio) return 0;
631  chio->SetEnergy(de1);
632  chio->GetACQParam("GG")->SetData((UShort_t)canal[ChIo_GG]);
633  chio->GetACQParam("PG")->SetData((UShort_t)canal[ChIo_PG]);
634  chio->GetACQParam("T")->SetData((UShort_t)canal[ChIo_T]);
635  si->SetEnergy(de2);
636  csi->SetEnergy(de3);
637  si->GetACQParam("GG")->SetData((UShort_t)canal[Si_GG]);
638  si->GetACQParam("PG")->SetData((UShort_t)canal[Si_PG]);
639  si->GetACQParam("T")->SetData((UShort_t)canal[Si_T]);
640  csi->GetACQParam("R")->SetData((UShort_t)canal[CsI_R]);
641  csi->GetACQParam("L")->SetData((UShort_t)canal[CsI_L]);
642  csi->GetACQParam("T")->SetData((UShort_t)mt);
643  identifying_telescope = gIndra->GetIDTelescope(Form("SI_CSI_%02d%02d", ring, mod));
644  return csi;
645  }
646  else if (ring >= 10) {
647  KVCsI* csi = (KVCsI*) gIndra->GetDetectorByType(ring, mod, CsI_R);
648  if (!csi) return 0;
649  KVChIo* chio = (KVChIo*)csi->GetChIo();
650  if (!chio) return 0;
651  KVSi75* si75 = (KVSi75*)gIndra->GetDetectorByType(ring, mod, Si75_GG);
652  if (!si75) return 0;
653  KVSiLi* sili = (KVSiLi*)gIndra->GetDetectorByType(ring, mod, SiLi_GG);
654  if (!sili) return 0;
655  chio->SetEnergy(de1);
656  si75->SetEnergy(de4);
657  sili->SetEnergy(de5);
658  chio->GetACQParam("GG")->SetData((UShort_t)canal[ChIo_GG]);
659  chio->GetACQParam("PG")->SetData((UShort_t)canal[ChIo_PG]);
660  chio->GetACQParam("T")->SetData((UShort_t)canal[ChIo_T]);
661  si75->GetACQParam("GG")->SetData((UShort_t)canal[Si75_GG]);
662  si75->GetACQParam("PG")->SetData((UShort_t)canal[Si75_PG]);
663  si75->GetACQParam("T")->SetData((UShort_t)canal[Si75_T]);
664  sili->GetACQParam("GG")->SetData((UShort_t)canal[SiLi_GG]);
665  sili->GetACQParam("PG")->SetData((UShort_t)canal[SiLi_PG]);
666  sili->GetACQParam("T")->SetData((UShort_t)canal[SiLi_T]);
667  if (canal[CsI_R] > 0.) {
668  csi->SetEnergy(de3);
669  csi->GetACQParam("R")->SetData((UShort_t)canal[CsI_R]);
670  csi->GetACQParam("L")->SetData((UShort_t)canal[CsI_L]);
671  csi->GetACQParam("T")->SetData((UShort_t)mt);
672  identifying_telescope = gIndra->GetIDTelescope(Form("SILI_CSI_%02d%02d", ring, mod));
673  return csi;
674  }
675  else {
676  identifying_telescope = gIndra->GetIDTelescope(Form("SI75_SILI_%02d%02d", ring, mod));
677  return sili;
678  }
679  }
680  return 0;
681 }
682 
683 
693 
695 {
696 // cout << "Code0" << endl;
697 //code = 0 : gamma
698 // ** Mention speciale pour les Gammas :
699 //
700 // Pour ceux-ci, on ne stocke que certaines infos disponibles :
701 //
702 // - Signal codeur CsI/Phoswich rapide --> de3(i) canaux en negatif
703 // - Marqueur de temps --> Mt(i)
704 
705  if (ring == 1 && !camp4) {
706  //phoswich
707  KVPhoswich* phos = (KVPhoswich*)gIndra->GetDetectorByType(ring, mod, Phos_R);
708  if (!phos) return 0;
709  phos->GetACQParam("R")->SetData((UShort_t)TMath::Abs(de3));
710  phos->GetACQParam("T")->SetData((UShort_t)mt);
711  identifying_telescope = gIndra->GetIDTelescope(Form("PHOS_R_L_%02d", mod));
712  return phos;
713  }
714  KVCsI* csi = (KVCsI*) gIndra->GetDetectorByType(ring, mod, CsI_R);
715  if (csi) {
716  csi->GetACQParam("T")->SetData((UShort_t)mt);
717  csi->GetACQParam("R")->SetData((UShort_t)TMath::Abs(de3));
718  }
719  identifying_telescope = gIndra->GetIDTelescope(Form("CSI_R_L_%02d%02d", ring, mod));
720  return csi;
721 }
722 
723 
727 
729 {
730 // cout << "Code1" << endl;
731 //code = 1 : neutron (seulement couronnes 2 a 9)
732 
733  if (ring < 2 || ring > 9) return 0;
734 
735  KVCsI* csi = (KVCsI*) gIndra->GetDetectorByType(ring, mod, CsI_R);
736  if (!csi) return 0;
737  csi->SetEnergy(de3);
738  csi->GetACQParam("R")->SetData((UShort_t)canal[CsI_R]);
739  csi->GetACQParam("L")->SetData((UShort_t)canal[CsI_L]);
740  csi->GetACQParam("T")->SetData((UShort_t)mt);
741  identifying_telescope = gIndra->GetIDTelescope(Form("CSI_R_L_%02d%02d", ring, mod));
742  return csi;
743 }
744 
745 
768 
770 {
771 // code = 4 : fragment identifie dans ChIo-Si
772 // ou fragment ou part. identifie dans ChIo-Si75
773  //for 2nd campaign: code=4 => phoswich (ring 1) also
774 // code = 5 : fragment non-ident. (montee de Bragg)(Z min.)
775 // ou stoppe dans Chio (Z min)
776 //code = 6 : frag. cree par la coherence CsI (stoppe ds Si.)
777 // En general, lorsqu'il n'y a pas "multiple hit" dans un module
778 // INDRA, canal contient les valeurs bruts de codeurs. Dans le
779 // cas d'un "multiple hit" avec creation de particule , canal
780 // contient les valeurs estimees de la contribution CHIO ou SI
781 // obtenues par la coherence.
782 // code = 8 : multiple comptage dans ChIo avec arret
783 // - dans les siliciums (cour. 2-9)
784 // - dans les CsI (cour. 10-17)
785 // Nota: Pour le code 8 la contribution ChIo est partagee egalement entre
786 // toutes les particules qui s'arretent dans les modules suivants.
787 // Le Z individuel est donc surement faux mais cela permet
788 // d'avoir une assez bonne estimation du Z.total.
789 //
790 // 09/2012 pour la 2e campagne presence de code14 qui semble etre
791 // des codes 4 donc on les traite comme tels
792 
793  if (ring == 1 && camp2) {
794  //phoswich
795  KVPhoswich* phos =
796  (KVPhoswich*)gIndra->GetDetectorByType(ring, mod, Phos_R);
797  if (!phos) return 0;
798  phos->SetEnergy(de1);
799  phos->GetACQParam("R")->SetData((UShort_t)canal[Phos_R]);
800  phos->GetACQParam("L")->SetData((UShort_t)canal[Phos_L]);
801  phos->GetACQParam("T")->SetData((UShort_t)mt);
802  identifying_telescope = gIndra->GetIDTelescope(Form("PHOS_R_L_%02d", mod));
803  return phos;
804  }
805  else if (ring >= 2 && ring <= 9) {
806  //chio-si
807  KVSilicon* si = (KVSilicon*) gIndra->GetDetectorByType(ring, mod, Si_GG);
808  if (!si) return 0;
809  KVChIo* chio = (KVChIo*)si->GetChIo();
810  if (!chio) return 0;
811  chio->SetEnergy(de1);
812  chio->GetACQParam("GG")->SetData((UShort_t)canal[ChIo_GG]);
813  chio->GetACQParam("PG")->SetData((UShort_t)canal[ChIo_PG]);
814  chio->GetACQParam("T")->SetData((UShort_t)canal[ChIo_T]);
815  if (de2 > 0.0) {
816  si->SetEnergy(de2);
817  si->GetACQParam("GG")->SetData((UShort_t)canal[Si_GG]);
818  si->GetACQParam("PG")->SetData((UShort_t)canal[Si_PG]);
819  si->GetACQParam("T")->SetData((UShort_t)canal[Si_T]);
820  identifying_telescope = gIndra->GetIDTelescope(Form("CI_SI_%02d%02d", ring, mod));
821  return si;
822  }
823  else
824  return chio;
825  }
826  else if (ring >= 10) {
827  //chio-csi or chio-si75
828  KVCsI* csi = (KVCsI*) gIndra->GetDetectorByType(ring, mod, CsI_R);
829  if (!csi) return 0;
830  KVChIo* chio = (KVChIo*)csi->GetChIo();
831  if (!chio) return 0;
832  KVSi75* si75 = (KVSi75*)gIndra->GetDetectorByType(ring, mod, Si75_GG);
833  if (!si75) return 0;
834  chio->SetEnergy(de1);
835  chio->GetACQParam("GG")->SetData((UShort_t)canal[ChIo_GG]);
836  chio->GetACQParam("PG")->SetData((UShort_t)canal[ChIo_PG]);
837  chio->GetACQParam("T")->SetData((UShort_t)canal[ChIo_T]);
838  if (canal[CsI_R] > 0.) {
839  csi->SetEnergy(de3);
840  csi->GetACQParam("R")->SetData((UShort_t)canal[CsI_R]);
841  csi->GetACQParam("L")->SetData((UShort_t)canal[CsI_L]);
842  csi->GetACQParam("T")->SetData((UShort_t)mt);
843  identifying_telescope = gIndra->GetIDTelescope(Form("CI_CSI_%02d%02d", ring, mod));
844  return csi;
845  }
846  else if (de4 > 0.0) {
847  si75->SetEnergy(de4);
848  si75->GetACQParam("GG")->SetData((UShort_t)canal[Si75_GG]);
849  si75->GetACQParam("PG")->SetData((UShort_t)canal[Si75_PG]);
850  si75->GetACQParam("T")->SetData((UShort_t)canal[Si75_T]);
851  identifying_telescope = gIndra->GetIDTelescope(Form("CI_SI75_%02d%02d", ring, mod));
852  return si75;
853  }
854  else
855  return chio;
856  }
857  return 0;
858 }
859 
860 
864 
866 {
867 // cout << "Code7" << endl;
868 
869 //code = 7 : frag. cree par la coherence ChIo (stoppe ds ChIo)
870 
871  KVChIo* chio = (KVChIo*)gIndra->GetDetectorByType(ring, mod, ChIo_GG);
872  if (!chio) {
873  //if chio not found, may be because particule has (ring,mod) corresponding
874  //to a Si/CsI, and not those of a ChIo module
875  //try CsI with given ring number and module
876  KVCsI* csi = (KVCsI*)gIndra->GetDetectorByType(ring, mod, CsI_R);
877  if (csi) chio = (KVChIo*)csi->GetChIo();
878  if (!chio) return 0;
879  }
880  chio->SetEnergy(de1);
881  chio->GetACQParam("GG")->SetData((UShort_t)canal[ChIo_GG]);
882  chio->GetACQParam("PG")->SetData((UShort_t)canal[ChIo_PG]);
883  chio->GetACQParam("T")->SetData((UShort_t)canal[ChIo_T]);
884  return chio;
885 }
886 
887 
892 
894 {
895  //
896  //Lire un evenement
897  //
898 
899  FNMETHOD CodeFunc[] = {
911  };
912 
913  TString buff;
914  Int_t mul, num_ev;
915 
916  buff.ReadLine(f_in);
917  num_ev = buff.Atoi();
918 
919  if (f_in.good()) {
920  events_in_file++; //total number events in each file
921  events_read++; //total events read
922  events_good++;
923 
924  // Instead of using the DST event number (which may not have a continuous sequence, e.g. 1,3,4,7,12,...)
925  // we use the sequential event counter fEventNumber (=1,2,3,...) which is also used as the event number
926  // of the raw data tree which we are filling with one entry for each DST event using the informations
927  // on raw data (canal[..]) written in the DST files. This is to ensure correct synchronisation
928  //evt->SetNumber(num_ev);
929  evt->SetNumber(fEventNumber);
930 
931  f_in >> mul;
932  for (Int_t i = 0; i < mul; i++) {
933 
934  Int_t z, a, code, icou, imod, ecode, ncans;
935  Double_t z_indra, a_indra, ener;
936  for (Int_t gj = 0; gj < 16; gj++) canal[gj] = 0.0;
937  f_in >> a;
938  f_in >> z;
939  f_in >> z_indra;
940  f_in >> a_indra;
941 //Check energy is not NAN
942 // f_in >> ener;
943  buff.ReadLine(f_in);
944  if (buff.IsFloat()) {
945  ener = buff.Atof();
946  }
947  else {
948  ener = 0.;
949  cout << "BAD: event#" << num_ev << " z=" << z << " a=" << a << " / Non-numeric value read for ener. Set E=0" << endl;
950  cout << "(read: " << buff.Data() << ")" << endl;
951  }
952  f_in >> icou;
953  f_in >> imod;
954  f_in >> code;
955  f_in >> ecode;
956 
957 // f_in >> de_mylar;
958  buff.ReadLine(f_in);
959  if (buff.IsFloat()) de_mylar = buff.Atof();
960  else {
961  de_mylar = 0.;
962  cout << "BAD: event#" << num_ev << " z=" << z << " a=" << a << " / Non-numeric value read for de_mylar. Set de_mylar=0" << endl;
963  cout << "(read: " << buff.Data() << ")" << endl;
964  }
965  de_mylar *= 1.0; //just to avoid "variable set but not used" warning
966 
967 // f_in >> de1;
968  buff.ReadLine(f_in);
969  if (buff.IsFloat()) de1 = buff.Atof();
970  else {
971  de1 = 0.;
972  cout << "BAD: event#" << num_ev << " z=" << z << " a=" << a << " / Non-numeric value read for de1. Set de1=0" << endl;
973  cout << "(read: " << buff.Data() << ")" << endl;
974  }
975 
976 // f_in >> de2;
977  buff.ReadLine(f_in);
978  if (buff.IsFloat()) de2 = buff.Atof();
979  else {
980  de2 = 0.;
981  cout << "BAD: event#" << num_ev << " z=" << z << " a=" << a << " / Non-numeric value read for de2. Set de2=0" << endl;
982  cout << "(read: " << buff.Data() << ")" << endl;
983  }
984 
985 // f_in >> de3;
986  buff.ReadLine(f_in);
987  if (buff.IsFloat()) de3 = buff.Atof();
988  else {
989  de3 = 0.;
990  cout << "BAD: event#" << num_ev << " z=" << z << " a=" << a << " / Non-numeric value read for de3. Set de3=0" << endl;
991  cout << "(read: " << buff.Data() << ")" << endl;
992  }
993 
994 // f_in >> de4;
995  buff.ReadLine(f_in);
996  if (buff.IsFloat()) de4 = buff.Atof();
997  else {
998  de4 = 0.;
999  cout << "BAD: event#" << num_ev << " z=" << z << " a=" << a << " / Non-numeric value read for de4. Set de4=0" << endl;
1000  cout << "(read: " << buff.Data() << ")" << endl;
1001  }
1002 
1003 // f_in >> de5;
1004  buff.ReadLine(f_in);
1005  if (buff.IsFloat()) de5 = buff.Atof();
1006  else {
1007  de5 = 0.;
1008  cout << "BAD: event#" << num_ev << " z=" << z << " a=" << a << " / Non-numeric value read for de5. Set de5=0" << endl;
1009  cout << "(read: " << buff.Data() << ")" << endl;
1010  }
1011 
1012  f_in >> mt;
1013 
1014  if (!icou) {
1015  cout << "BAD: event#" << num_ev << " particle#" << i << " Z=" << z << " code=" << code << endl;
1016  }
1017 
1018  Int_t idf;
1019  for (Int_t mk = 1; mk <= 4; mk++) {
1020  f_in >> idf;
1021  code_idf[mk] = idf;
1022  }
1023 
1024  f_in >> ncans;
1025  if (ncans) {
1026  Int_t kk;
1027  Double_t chan;
1028  for (Int_t mk = 1; mk <= ncans; mk++) {
1029  f_in >> kk >> chan;
1030  canal[kk] = chan;
1031  }
1032  }
1033 
1034  if ((code < 11 || (camp2 && code == 14)) && icou > 0 && imod > 0 && icou < 18 && imod < 25) {
1035  //work out where particle stopped
1036  //set up corresponding detectors with energy losses and raw data (canal)
1037  //set stopping detector's marqueur de temps to the MT of the particle
1038  //(the other Mt are not written on the DST's).
1039  //then use KVINDRAReconNuc::Reconstruct to set up list of hit detectors,
1040  //hit group and telescope for this particle
1041  //
1042  //pour la 2e campagne presence de code14 qui semble etre
1043  //des codes 4 donc on les traite comme tels
1044 
1045  KVDetector* det_stop = 0;
1046  identifying_telescope = 0;
1047  if (code == 14) det_stop = (this->*CodeFunc[4])(icou, imod);
1048  else det_stop = (this->*CodeFunc[code])(icou, imod);
1049 
1050  if (!det_stop) {
1051  cout << "det_stop=0: event#" << num_ev << " particle#" << i << " icou=" << icou << " imod=" << imod <<
1052  " z=" << z << " a=" << a << " code=" << code << " ecode=" << ecode << endl;
1053  }
1054  else {
1055  KVINDRAReconNuc* tmp = evt->AddParticle();
1056  if (!tmp) {
1057  Fatal("lire_evt", "KVINDRAReconEvent::AddParticle() returns NULL pointer.");
1058  gSystem->Exit(1);
1059  }
1060  det_stop->GetACQParam("T")->SetData((UShort_t)mt);
1061  tmp->Reconstruct(det_stop);
1062  if (identifying_telescope) {
1063  tmp->SetIdentifyingTelescope(identifying_telescope);
1064  tmp->SetIsIdentified();
1065  tmp->SetZMeasured(kTRUE);
1066  }
1067  if (ecode > 0)
1068  tmp->SetIsCalibrated();
1069  //modify ID code for phoswich particles 2nd campaign
1070  //we put VEDA ID code=2 just like for 1st campaign
1071  //if(camp2 && icou==1 && code==4) code=2;
1074  tmp->SetZ(z);
1075  tmp->SetRealZ(z_indra);
1076  tmp->SetA(TMath::Abs(a));
1077  if (a > 0 && (code == 2 || code == 3)) {
1078  tmp->SetAMeasured(kTRUE);
1079  tmp->SetRealA(a_indra);
1080  }
1081  else {
1082  tmp->SetAMeasured(kFALSE);
1083  // reset Z in order to force calculation of mass according to default mass formula of KVINDRAReconNuc
1084  tmp->SetZ(z);
1085  }
1086  if (code == 1) {
1087  //set Z and realZ = 0 in case of 'neutron'
1088  tmp->SetZ(0);
1089  tmp->SetRealZ(0);
1090  tmp->SetA(1);
1091  tmp->SetRealA(0);
1092  }
1093  if (code == 5) {
1094  tmp->SetZMeasured(kFALSE);
1095  }
1096  tmp->SetEnergy(ener);
1097 
1098  tmp->SetEnergyChIo(de1 + de_mylar);
1099  tmp->SetEnergySi(de2);
1100  tmp->SetEnergyCsI(de3);
1101  tmp->SetEnergySi75(de4);
1102  tmp->SetEnergySiLi(de5);
1103 
1104 
1105  }//if(det_stop)
1106  }//if code<11
1107  }//for(int i=0;
1108 
1109  //write event to tree
1110  FillRawTree();
1111  data_tree->Fill();
1112  fEventNumber++;
1113  }//if(f_in.good
1114 
1115 }
1116 
1117 
1118 
1123 
1125 {
1126  //Prints list of available runs, sorted according to multiplicity
1127  //trigger, for selected dataset, data type/analysis task, and system
1128  //Returns list containing all run numbers
1129 
1130  KVNumberList all_runs =
1131  GetDataSet()->GetRunList(datatype.Data(), GetSystem());
1132  KVINDRADBRun* dbrun;
1133 
1134  //first read list and find what triggers are available
1135  vector<int> triggers;
1136  all_runs.Begin();
1137  while (!all_runs.End()) {
1138  dbrun = (KVINDRADBRun*)GetDataSet()->GetDataBase()->GetDBRun(all_runs.Next());
1139  if (triggers.size() == 0
1140  || std::find(triggers.begin(), triggers.end(), dbrun->GetTrigger()) != triggers.end()) {
1141  triggers.push_back(dbrun->GetTrigger());
1142  }
1143  }
1144  //sort triggers in ascending order
1145  std::sort(triggers.begin(), triggers.end());
1146 
1147  for (std::vector<int>::iterator it = triggers.begin(); it != triggers.end(); ++it) {
1148  cout << " ---> Trigger M>" << *it << endl;
1149  all_runs.Begin();
1150  while (!all_runs.End()) {
1151  dbrun = (KVINDRADBRun*)GetDataSet()->GetDataBase()->GetDBRun(all_runs.Next());
1152  if (dbrun->GetTrigger() == *it) {
1153  cout << " " << Form("%4d", dbrun->GetNumber());
1154  cout << Form("\t(%9llu events)", dbrun->GetEvents());
1155  cout << "\t[File written: " << dbrun->GetDatime().
1156  AsString() << "]";
1157  if (dbrun->GetComments())
1158  cout << "\t" << dbrun->GetComments();
1159  cout << endl;
1160  }
1161  }
1162  cout << endl;
1163  }
1164  return all_runs;
1165 }
1166 
1167 
1168 
1172 
1174 {
1175  // Put values of all "fired" acquisition parameters into the arrays ParNum/ParVal
1176  // and store them in the raw data tree
1177 
1178  NbParFired = 0;
1179  TIter next(params);
1180  KVACQParam* par = 0;
1181  while ((par = (KVACQParam*)next())) {
1182  if (par->Fired()) {
1183  ParVal[NbParFired] = par->GetCoderData();
1184  ParNum[NbParFired++] = par->GetNumber();
1185  }
1186  }
1187  rawtree->Fill();
1188 }
1189 
1190 
int Int_t
KVDataAnalyser * gDataAnalyser
KVDataRepositoryManager * gDataRepositoryManager
KVDataRepository * gDataRepository
KVDataSet * gDataSet
Definition: KVDataSet.cpp:30
KVINDRADB * gIndraDB
Definition: KVINDRADB.cpp:38
KVDetector *(KVINDRADstToRootTransfert::* FNMETHOD)(int, int)
KVINDRA * gIndra
Definition: KVINDRA.cpp:78
@ Si_PG
Definition: KVINDRA.h:47
@ ChIo_PG
Definition: KVINDRA.h:44
@ SiLi_GG
Definition: KVINDRA.h:55
@ CsI_L
Definition: KVINDRA.h:50
@ Si75_GG
Definition: KVINDRA.h:52
@ ChIo_T
Definition: KVINDRA.h:45
@ SiLi_T
Definition: KVINDRA.h:57
@ SiLi_PG
Definition: KVINDRA.h:56
@ ChIo_GG
Definition: KVINDRA.h:43
@ CsI_R
Definition: KVINDRA.h:49
@ Si75_T
Definition: KVINDRA.h:54
@ Si_GG
Definition: KVINDRA.h:46
@ Si_T
Definition: KVINDRA.h:48
@ Si75_PG
Definition: KVINDRA.h:53
@ Phos_R
Definition: KVINDRA.h:60
@ Phos_L
Definition: KVINDRA.h:61
ClassImp(KVPartitionList) void KVPartitionList
Initialisation.
unsigned short UShort_t
const Bool_t kFALSE
double Double_t
const Bool_t kTRUE
#define gROOT
char * Form(const char *fmt,...)
R__EXTERN TSystem * gSystem
GANIL VXI/VME acquisition parameter.
Definition: KVACQParam.h:15
Short_t GetCoderData() const
Definition: KVACQParam.h:64
void UseInternalDataMember()
Definition: KVACQParam.h:144
Bool_t Fired(Option_t *what="") const
Definition: KVACQParam.h:85
void SetData(UShort_t val)
Definition: KVACQParam.h:58
static const Char_t * GetBINDIRFilePath(const Char_t *namefile="")
Definition: KVBase.cpp:118
virtual void SetNumber(UInt_t num)
Definition: KVBase.h:209
UInt_t GetNumber() const
Definition: KVBase.h:213
Ionisation chamber detectors of the INDRA multidetector array.
Definition: KVChIo.h:29
CsI(Tl) scintillation detectors of the INDRA multidetector array.
Definition: KVCsI.h:15
virtual Int_t GetNumber() const
Definition: KVDBRecord.h:72
ULong64_t GetEvents() const
Definition: KVDBRun.h:133
KVDBSystem * GetSystem() const
Definition: KVDBRun.cpp:242
const TDatime & GetDatime() const
Definition: KVDBRun.h:113
Int_t GetTrigger() const
Definition: KVDBRun.h:103
const Char_t * GetComments() const
Definition: KVDBRun.h:146
void WriteBatchInfo(TTree *)
KVDataSet * GetDataSet(const Char_t *repository, const Char_t *dataset) const
Return pointer to named dataset in the given repository.
virtual void CopyFileFromRepository(const KVDataSet *dataset, const Char_t *datatype, const Char_t *filename, const Char_t *destination)
Manage an experimental dataset corresponding to a given experiment or campaign.
Definition: KVDataSet.h:207
TString GetOutputRepository(const Char_t *taskname) const
Definition: KVDataSet.cpp:2057
FileType * OpenRunfile(const Char_t *type, Int_t run)
Definition: KVDataSet.h:333
const Char_t * GetRunfileName(const Char_t *type, Int_t run) const
Definition: KVDataSet.cpp:913
virtual const Char_t * GetDataPathSubdir() const
Returns name of top-level directory in data repository used to store data files for this dataset.
Definition: KVDataSet.h:261
TString GetFullPathToRunfile(const Char_t *type, Int_t run) const
Definition: KVDataSet.cpp:887
void CommitRunfile(const Char_t *type, Int_t run, TFile *file)
Definition: KVDataSet.cpp:1311
TFile * NewRunfile(const Char_t *type, Int_t run)
Definition: KVDataSet.cpp:1057
Base class for detector geometry description, interface to energy-loss calculations.
Definition: KVDetector.h:121
virtual void SetEnergy(Double_t e) const
Definition: KVDetector.h:320
virtual KVACQParam * GetACQParam(const Char_t *) const
Definition: KVDetector.cpp:646
virtual void Clear(Option_t *opt="")
Definition: KVEvent.cpp:200
static void MakeEventBranch(TTree *tree, const TString &branchname, const TString &classname, T &event, Int_t bufsize=10000000)
Definition: KVEvent.h:468
Reads GANIL acquisition files (EBYEDAT)
virtual void SetUserTree(TTree *, Option_t *="arrays")
static UChar_t VedaECodeToBitmask(UChar_t veda_e_code)
static UShort_t VedaIDCodeToBitmask(UChar_t veda_id_code)
Database entry for each run of an INDRA experiment.
Definition: KVINDRADBRun.h:29
KVINDRADBRun * GetRun(Int_t run) const
Definition: KVINDRADB.h:132
virtual KVACQParam * GetACQParam(const Char_t *) const
KVINDRADetector * GetChIo() const
Conversion of INDRA DST file to KaliVeda ROOT format.
KVDetector * Code1(int ring, int mod)
KVDetector * Code2and9and10(int ring, int mod)
virtual ~KVINDRADstToRootTransfert()
Destructor.
KVDetector * Code0(int ring, int mod)
void lire_evt(std::ifstream &f_in, KVINDRAReconEvent *evt)
virtual KVNumberList PrintAvailableRuns(KVString &datatype)
KVDetector * Code4and5and6and8(int ring, int mod)
KVDetector * Code3(int ring, int mod)
KVDetector * Code7(int ring, int mod)
Event reconstructed from energy losses in INDRA multidetector.
KVINDRAReconNuc * AddParticle()
Wrapper for KVEvent::GetNextParticle casting result to KVINDRAReconNuc*.
Nuclei reconstructed from data measured in the INDRA array.
void SetEnergyCsI(Float_t val)
void SetEnergySi(Float_t val)
void SetEnergyChIo(Float_t val)
void SetEnergySi75(Float_t val)
void SetEnergySiLi(Float_t val)
virtual KVINDRADetector * GetDetectorByType(UInt_t cou, UInt_t mod, UInt_t type) const
Definition: KVINDRA.cpp:651
const KVSeqCollection * GetACQParams() const
UInt_t GetCurrentRunNumber() const
static KVMultiDetArray * MakeMultiDetector(const Char_t *dataset_name, Int_t run=-1, TString classname="KVMultiDetArray")
KVIDTelescope * GetIDTelescope(const Char_t *name) const
Return pointer to DeltaE-E ID Telescope with "name".
void SetA(Int_t a)
Definition: KVNucleus.cpp:655
void SetZ(Int_t z, Char_t mt=-1)
Definition: KVNucleus.cpp:704
Strings used to represent a set of ranges of values.
Definition: KVNumberList.h:83
Bool_t End(void) const
Definition: KVNumberList.h:196
void Begin(void) const
Int_t Next(void) const
void SetEnergy(Double_t e)
Definition: KVParticle.h:560
Phoswich detector in the INDRA array (first 3 campaigns)
Definition: KVPhoswich.h:35
void SetEnergy(Double_t e) const
Set energy lost in both layers.
Definition: KVPhoswich.cpp:92
virtual void Reconstruct(KVDetector *kvd)
virtual void SetAMeasured(Bool_t yes=kTRUE)
virtual void SetZMeasured(Bool_t yes=kTRUE)
void SetIdentifyingTelescope(KVIDTelescope *i)
80um silicon detector for INDRA etalon telescopes
Definition: KVSilicon.h:60
2mm + 40um dead zone Si(Li) detector for INDRA etalon telescopes
Definition: KVSilicon.h:77
Silicon detectors of the INDRA array.
Definition: KVSilicon.h:18
Extension of ROOT TString class which allows backwards compatibility with ROOT v3....
Definition: KVString.h:72
const char * AsString() const
Bool_t cd(const char *path=nullptr) override
virtual const char * GetName() const
virtual void SetTitle(const char *title="")
virtual const char * GetTitle() const
virtual void SetName(const char *name)
TObject * At(Int_t idx) const
Int_t Atoi() const
Double_t Atof() const
Bool_t IsFloat() const
TObjArray * Tokenize(const TString &delim) const
const char * Data() const
void Form(const char *fmt,...)
std::istream & ReadLine(std::istream &str, Bool_t skipWhite=kTRUE)
virtual int GetProcInfo(ProcInfo_t *info) const
virtual void Exit(int code, Bool_t mode=kTRUE)
virtual TString GetFromPipe(const char *command)
virtual char * ExpandPathName(const char *path)
virtual void Setenv(const char *name, const char *value)
TLine * line
void Info(const char *location, const char *va_(fmt),...)
void Error(const char *location, const char *va_(fmt),...)
void Fatal(const char *location, const char *va_(fmt),...)
void Warning(const char *location, const char *va_(fmt),...)
void End()
Double_t Abs(Double_t d)
const char * String
Long_t fMemVirtual
Float_t fCpuUser
auto * a