KaliVeda  1.13/01
Heavy-Ion Analysis Toolkit
KVINDRA.cpp
Go to the documentation of this file.
1 /***************************************************************************
2 $Id: KVINDRA.cpp,v 1.68 2009/01/21 10:05:51 franklan Exp $
3  kvindra.cpp - description
4  -------------------
5  begin : Mon May 20 2002
6  copyright : (C) 2002 by J.D. Frankland
7  email : frankland@ganil.fr
8  ***************************************************************************/
9 
10 /***************************************************************************
11  * *
12  * This program is free software; you can redistribute it and/or modify *
13  * it under the terms of the GNU General Public License as published by *
14  * the Free Software Foundation; either version 2 of the License, or *
15  * (at your option) any later version. *
16  * *
17  ***************************************************************************/
18 
19 #include "Riostream.h"
20 #include "KVINDRA.h"
21 #include "KVMaterial.h"
22 #include "KVDetector.h"
23 #include "KVTelescope.h"
24 #include "KVIDTelescope.h"
25 #include "KVIDSiCsI.h"
26 #include "KVIDINDRACsI.h"
27 #include "KVIDChIoSi.h"
28 #include "KVIDChIoCsI.h"
29 #include "KVIDChIoSi75.h"
30 #include "KVIDSi75SiLi.h"
31 #include "KVIDSiLiCsI.h"
32 #include "KVIDPhoswich.h"
33 #include "KVLayer.h"
34 #include "KVRing.h"
35 #include "KVGroup.h"
36 #include "KVSilicon.h"
37 #include "KVChIo.h"
38 #include "KVCsI.h"
39 #include "KVPhoswich.h"
40 #include "KVNucleus.h"
41 #include "KVDetectorEvent.h"
42 #include "KVNucleusEvent.h"
43 #include "KVINDRAReconEvent.h"
44 #include "TEnv.h"
45 #include "KVFileReader.h"
46 #include "INDRAGeometryBuilder.h"
47 #include <KVASGroup.h>
48 #include "TROOT.h"
49 #include "KVGeoNavigator.h"
50 #ifdef WITH_BUILTIN_GRU
51 #include <KVGANILDataReader.h>
52 #endif
53 #include <KVGeoImport.h>
55 #include <KVRawDataReader.h>
56 #ifdef WITH_MFM
57 #include "KVMFMDataFileReader.h"
58 #include "MFMEbyedatFrame.h"
59 #ifdef WITH_MESYTEC
60 #include "MFMMesytecMDPPFrame.h"
61 #include "mesytec_buffer_reader.h"
62 #endif
63 #endif
65 #include <cassert>
66 
67 using namespace std;
68 
70 
71 //Use this static array to translate EBaseIndra_type
72 //signal type to a string giving the signal type
73 
74 
77  "", //dummy for index=0
78  "GG", "PG", "T",
79  "GG", "PG", "T",
80  "R", "L", "T",
81  "GG", "PG", "T",
82  "GG", "PG", "T"
83 };
84 
85 
86 #define KV_DEBUG 1
87 
89 
93 
95 {
96  //Default constructor
97  //Set up lists of ChIo, Si, CsI, Phoswich
98  fChIo = new KVHashList();
99  fChIo->SetCleanup(kTRUE);
100  fSi = new KVHashList();
101  fSi->SetCleanup(kTRUE);
102  fCsI = new KVHashList();
103  fCsI->SetCleanup(kTRUE);
104  fPhoswich = new KVHashList();
105  fPhoswich->SetCleanup(kTRUE);
106  fTrigger = 0;
107  gIndra = this;
108  fPHDSet = kFALSE;
109  fSelecteur = nullptr;
110  // unlike a normal multidetector, INDRA does not own its detectors
111  // they are owned by the TELESCOPE structures
112  SetOwnsDetectors(kFALSE);
113  fMesytecData = kFALSE;
114  fEbyedatData = kFALSE;
115 }
116 
117 
118 
122 
123 KVINDRA::~KVINDRA()
124 {
125  //Delets lists of ChIo, Si, etc.
126  //Resets global gIndra pointer.
127 
128  if (fChIo && fChIo->TestBit(kNotDeleted)) {
129  fChIo->Clear();
130  delete fChIo;
131  }
132  fChIo = 0;
133  if (fSi && fSi->TestBit(kNotDeleted)) {
134  fSi->Clear();
135  delete fSi;
136  }
137  fSi = 0;
138  if (fCsI && fCsI->TestBit(kNotDeleted)) {
139  fCsI->Clear();
140  delete fCsI;
141  }
142  fCsI = 0;
143  if (fPhoswich && fPhoswich->TestBit(kNotDeleted)) {
144  fPhoswich->Clear();
145  delete fPhoswich;
146  }
147  fPhoswich = 0;
148  if (fSelecteur) delete fSelecteur;
149  gIndra = 0;
150 }
151 
152 
153 
154 
168 
170 {
171  // Construction of INDRA detector array.
172  //
173  // Uses infos in file
174  // $KVROOT/KVFiles/data/indra_struct.[dataset].env
175  // or
176  // $KVROOT/KVFiles/data/indra_struct.env
177  //
178  // if no dataset-specific file found
179  //
180  // Alternatively, by defining the variable
181  //
182  // [dataset].INDRA.StructureFile: [path to file]
183 
184 
185  TString path = Form("indra-struct.%s.env", fDataSet.Data());
186  TString path2 = GetDataSetEnv(fDataSet, "INDRA.StructureFile", "");
187  if (path2 == "") {
188  if (!SearchKVFile(path.Data(), path2, "data")) {
189  path = "indra-struct.env";
190  SearchKVFile(path.Data(), path2, "data");
191  }
192  }
193  fStrucInfos.ReadFile(path2, kEnvAll);
194 
195  KVString lruns = fStrucInfos.GetValue("AddOnForRuns", "");
196  //test if additional geometrical specification exists
197  if (lruns != "") {
198  lruns.Begin(",");
199  while (!lruns.End()) {
200  KVString sruns = lruns.Next();
201  KVNumberList nlr(sruns.Data());
202  //the current run needs specific geometry
203  if (nlr.Contains(fCurrentRun)) {
204  path = fStrucInfos.GetValue(sruns.Data(), "");
205  Info("BuildGeometry", "Additional geometry for run=%d in file #%s#", fCurrentRun, path.Data());
206  SearchKVFile(path.Data(), path2, "data");
207  if (path2 == "") {
208  Warning("BuildGeometry", "fichier %s inconnu", path.Data());
209  } else {
210  fStrucInfos.ReadFile(path2, kEnvChange);
211  }
212  }
213  }
214  }
215 
216  SetName(fStrucInfos.GetValue("INDRA.Name", ""));
217  SetTitle(fStrucInfos.GetValue("INDRA.Title", ""));
218 
219  KVString layers = fStrucInfos.GetValue("INDRA.Layers", "");
220  layers.Begin(" ");
221  while (!layers.End()) BuildLayer(layers.Next());
222 }
223 
224 
225 
228 
229 void KVINDRA::BuildLayer(const Char_t* name)
230 {
231  // Build layer 'name' with infos in file "$KVROOT/KVFiles/data/indra-struct.[dataset].env"
232 
233  KVLayer* layer = new KVLayer;
234  Int_t number = fStrucInfos.GetValue(Form("INDRA.Layer.%s.Number", name), 0);
235  layer->SetName(name);
236  layer->SetNumber(number);
237  Add(layer);
238  Info("BuildLayer", "Building layer %d:%s", number, name);
239  KVNumberList rings = fStrucInfos.GetValue(Form("INDRA.Layer.%s.Rings", name), "");
240  rings.Begin();
241  while (!rings.End()) {
242  Int_t ring_num = rings.Next();
243  TString prefix = Form("INDRA.Layer.%s.Ring.%d", name, ring_num);
244  KVRing* ring = BuildRing(ring_num, prefix);
245  layer->Add(ring);
246  }
247 }
248 
249 
250 
253 
254 KVRing* KVINDRA::BuildRing(Int_t number, const Char_t* prefix)
255 {
256  // Build ring with infos in file "$KVROOT/KVFiles/data/indra-struct.[dataset].env"
257 
258  KVRing* ring = new KVRing;
259  Info("BuildRing", "Building ring %d (%s)", number, prefix);
260 
261  Double_t thmin = fStrucInfos.GetValue(Form("%s.ThMin", prefix), 0.0);
262  Double_t thmax = fStrucInfos.GetValue(Form("%s.ThMax", prefix), 0.0);
263  Double_t phi_0 = fStrucInfos.GetValue(Form("%s.Phi0", prefix), 0.0);
264  Int_t ntel = fStrucInfos.GetValue(Form("%s.Ntel", prefix), 0);
265  Int_t step = fStrucInfos.GetValue(Form("%s.Step", prefix), 1);
266  Int_t num_first = fStrucInfos.GetValue(Form("%s.Nfirst", prefix), 1);
267  Double_t dphi = fStrucInfos.GetValue(Form("%s.Dphi", prefix), -1.0);
268  KVNumberList absent = fStrucInfos.GetValue(Form("%s.Absent", prefix), "");
269 
270  ring->SetPolarMinMax(thmin, thmax);
271  KVINDRATelescope* kvt = BuildTelescope(prefix, num_first);
272  Double_t phi_min = phi_0 - 0.5 * (kvt->GetAzimuthalWidth());
273  phi_min = (phi_min < 0. ? phi_min + 360. : phi_min);
274  dphi = (dphi < 0. ? 360. / ntel : dphi);
275  Double_t phi_max =
276  phi_0 + ((ntel - 1.0) * dphi) + 0.5 * (kvt->GetAzimuthalWidth());
277  phi_max = (phi_max > 360. ? phi_max - 360. : phi_max);
278  ring->SetAzimuthalMinMax(phi_min, phi_max);
279  ring->SetNumber(number);
280 
281  Double_t phi = phi_0;
282  UInt_t counter = num_first;
283  for (int i = 0; i < ntel; i++) {
284  kvt->SetPolarMinMax(thmin, thmax);
285  kvt->SetPhi(phi);
286  phi = (phi + dphi > 360. ? (phi + dphi - 360.) : (phi + dphi));
287  kvt->SetNumber(counter);
288  if (!absent.Contains(counter)) ring->Add(kvt);
289  else delete kvt;
290  counter += step;
291  if (i + 1 < ntel) kvt = BuildTelescope(prefix, counter);
292  }
293  ring->SetName(Form("RING_%d", number));
294  return ring;
295 }
296 
297 
298 
302 
304 {
305  // Build telescope from infos in file "$KVROOT/KVFiles/data/indra-struct.[dataset].env"
306 
307  //Info("BuildTelescope", "Building telescope %s",name);
309  KVString telescopes = fStrucInfos.GetValue(Form("%s.Telescope", prefix), "");
310  telescopes.Begin(" ");
311  KVString name;
312  while (!telescopes.End()) {
313  name = telescopes.Next();
314  KVNumberList mods = fStrucInfos.GetValue(Form("%s.Telescope.%s.Modules", prefix, name.Data()), "1-100");
315  if (mods.Contains(module)) break;
316  }
317 
318  KVString detectors = fStrucInfos.GetValue(Form("INDRA.Telescope.%s.Detectors", name.Data()), "");
319  Double_t aziwidth = fStrucInfos.GetValue(Form("INDRA.Telescope.%s.AziWidth", name.Data()), 1.0);
320  Double_t dist = fStrucInfos.GetValue(Form("%s.Dist", prefix), 0.0);
321  kvt->SetDistance(dist);
322  detectors.Begin(" ");
323  Int_t idet = 1;
324  while (!detectors.End()) {
325  KVString detector = detectors.Next();
326  detector.Begin("()");
327  KVString dettype = detector.Next();
328  Double_t detthick = detector.Next().Atof();
329  KVDetector* det = KVDetector::MakeDetector(Form("%s.%s", fDataSet.Data(), dettype.Data()), detthick);
330  kvt->Add(det);
331  Double_t depth = fStrucInfos.GetValue(Form("INDRA.Telescope.%s.Depth.%s", name.Data(), dettype.Data()), 0.);
332  kvt->SetDepth(idet, depth);
333  idet++;
334  }
335  kvt->SetAzimuthalWidth(aziwidth);
336  return kvt;
337 }
338 
339 
340 
349 
351 {
352  // Kludge to make INDRA ROOT geometry work like any other
353  // Normally ID telescopes are deduced from successive detectors on the different trajectories
354  // Each trajectory then possesses its list of ID telescopes
355  // These lists are then used when reconstruction trajectories are calculated
356  //
357  // When INDRA ROOT geometry is used, trajectory lists need to be filled
358  // by hand before reconstruction trajectories are calculated
359 
360  TIter it_traj(GetTrajectories());
361  KVGeoDNTrajectory* tr;
362  while ((tr = (KVGeoDNTrajectory*)it_traj())) {
363  // get first detector on trajectory
364  KVDetector* det = tr->GetNodeAt(0)->GetDetector();
365  // list of 'aligned' id telescopes
366  TIter it_idt(det->GetAlignedIDTelescopes());
367  KVIDTelescope* idt;
368  while ((idt = (KVIDTelescope*)it_idt())) {
369  if (tr->ContainsAll(idt->GetDetectors())) { // all detectors on trajectory
370  if (idt->GetDetectors()->GetEntries() > 1) {
371  // make sure dE and E detector are consecutive on trajectory
372  // i.e. dE has to be immediately in front of E on the trajectory
373  if (tr->GetNodeInFront(idt->GetDetector(2)->GetNode()) == idt->GetDetector(1)->GetNode())
374  tr->AccessIDTelescopeList()->Add(idt);
375  } else
376  tr->AccessIDTelescopeList()->Add(idt);//single-detector telescope
377  }
378  }
379  }
380 }
381 
382 
383 
384 
388 
390 {
391  // Overrides KVASMultiDetArray::Build
392  // Correspondance between CsI detectors and pin lasers is set up if known.
393 
394  if (run != -1) {
395  fCurrentRun = run;
396  }
397  BuildGeometry();
398 
399  MakeListOfDetectors();
400 
401  SetGroupsAndIDTelescopes();
402 
403  SetNamesOfIDTelescopes();
404 
405  SetDetectorThicknesses();
406 
407  SetIdentifications();
408 
409  //set flag to say Build() was called
410  SetBit(kIsBuilt);
411 
412  SetPinLasersForCsI();
413 
414  SetExpectedDetectorSignalNames();
415 
416  // link EBYEDAT parameter names to detectors
417  TIter next_det(GetDetectors());
418  KVDetector* det;
419  while ((det = (KVDetector*)next_det())) {
420  TString det_name(det->GetName());
421  if (det_name.BeginsWith("CI") || det_name.BeginsWith("SI")) {
422  fEbyedatParamDetMap.SetValue(det_name + "_PG", det_name);
423  fEbyedatParamDetMap.SetValue(det_name + "_GG", det_name);
424  } else if (det_name.BeginsWith("CSI")) {
425  fEbyedatParamDetMap.SetValue(det_name + "_R", det_name);
426  fEbyedatParamDetMap.SetValue(det_name + "_L", det_name);
427  }
428  // handle change of "time marker" parameter type for certain experiments
429  TString time_marker_type = gDataSet->GetDataSetEnv(Form("KVACQParam.%s.T", det_name.Data()), "T");
430  time_marker_type.Prepend("_");
431  fEbyedatParamDetMap.SetValue(det_name + time_marker_type, det_name);
432  }
433 }
434 
435 
436 //void KVINDRA::SetArrayACQParams()
437 //{
438 // // Overrides KVASMultiDetArray::SetArrayACQParams() in order to
439 // // add the following acquisition parameters which are not associated to a detector:
440 // //
441 // // STAT_EVE
442 // // R_DEC
443 // // CONFIG
444 // // PILA_01_PG
445 // // PILA_01_GG
446 // // PILA_02_PG
447 // // PILA_02_GG
448 // // PILA_03_PG
449 // // PILA_03_GG
450 // // PILA_04_PG
451 // // PILA_04_GG
452 // // PILA_05_PG
453 // // PILA_05_GG
454 // // PILA_06_PG
455 // // PILA_06_GG
456 // // PILA_07_PG
457 // // PILA_07_GG
458 // // PILA_08_PG
459 // // PILA_08_GG
460 // // SI_PIN1_PG
461 // // SI_PIN1_GG
462 // // SI_PIN2_PG
463 // // SI_PIN2_GG
464 // //
465 // // We also create and initialize the KVINDRATriggerInfo object (fSelecteur) used to read
466 // // the status of the DAQ trigger event by event (access through GetTriggerInfo()).
467 
468 // KVACQParam* ste = new KVACQParam("STAT_EVE");
469 // AddArrayACQParam(ste);
470 // KVACQParam* dec = new KVACQParam("R_DEC");
471 // AddArrayACQParam(dec);
472 // KVACQParam* conf = new KVACQParam("CONFIG");
473 // AddArrayACQParam(conf);
474 // fSelecteur = new KVINDRATriggerInfo;
475 // fSelecteur->SetSTAT_EVE_PAR(ste);
476 // fSelecteur->SetR_DEC_PAR(dec);
477 // fSelecteur->SetVXCONFIG_PAR(conf);
478 
479 // AddArrayACQParam(new KVACQParam("PILA_01_PG"));
480 // AddArrayACQParam(new KVACQParam("PILA_01_GG"));
481 // AddArrayACQParam(new KVACQParam("PILA_02_PG"));
482 // AddArrayACQParam(new KVACQParam("PILA_02_GG"));
483 // AddArrayACQParam(new KVACQParam("PILA_03_PG"));
484 // AddArrayACQParam(new KVACQParam("PILA_03_GG"));
485 // AddArrayACQParam(new KVACQParam("PILA_04_PG"));
486 // AddArrayACQParam(new KVACQParam("PILA_04_GG"));
487 // AddArrayACQParam(new KVACQParam("PILA_05_PG"));
488 // AddArrayACQParam(new KVACQParam("PILA_05_GG"));
489 // AddArrayACQParam(new KVACQParam("PILA_06_PG"));
490 // AddArrayACQParam(new KVACQParam("PILA_06_GG"));
491 // AddArrayACQParam(new KVACQParam("PILA_07_PG"));
492 // AddArrayACQParam(new KVACQParam("PILA_07_GG"));
493 // AddArrayACQParam(new KVACQParam("PILA_08_PG"));
494 // AddArrayACQParam(new KVACQParam("PILA_08_GG"));
495 // AddArrayACQParam(new KVACQParam("SI_PIN1_PG"));
496 // AddArrayACQParam(new KVACQParam("SI_PIN1_GG"));
497 // AddArrayACQParam(new KVACQParam("SI_PIN2_PG"));
498 // AddArrayACQParam(new KVACQParam("SI_PIN2_GG"));
499 //}
500 
501 
502 
505 
507 {
508  //Overrides KVASMultiDetArray method to add FillListsOfDetectorsByType()
509 
511  FillListsOfDetectorsByType();
512 }
513 
514 
515 
518 
520 {
521  //Fill lists of ChIo, Si, CsI and phoswich
522 
523  fChIo->Clear();
524  fSi->Clear();
525  fCsI->Clear();
526  fPhoswich->Clear();
527  TIter next_det(GetDetectors());
528  KVDetector* kvd;
529  while ((kvd = (KVDetector*) next_det())) {
530  kvd->SetNameOfArray("INDRA");
531  if (kvd->InheritsFrom("KVChIo")) {
532  fChIo->Add(kvd);
533  }
534  if (kvd->InheritsFrom("KVSilicon")) {
535  fSi->Add(kvd);
536  }
537  if (kvd->InheritsFrom("KVCsI")) {
538  fCsI->Add(kvd);
539  }
540  if (kvd->InheritsFrom("KVPhoswich")) {
541  fPhoswich->Add(kvd);
542  }
543  }
544 }
545 
546 
547 
548 
554 
556 {
557  //Find groups of telescopes in angular alignment placed on different layers.
558  //List is in fGroups.
559  //Also creates all ID telescopes in array and stores them in fIDTelescopes.
560  //Any previous groups/idtelescopes are deleted beforehand.
561  CalculateGroupsFromGeometry();
562 
563  //now read list of groups and create list of ID telescopes
564  CreateIDTelescopesInGroups();
565 }
566 
567 
568 
572 
574 {
575 
576  //Returns a pointer to the Ionisation Chamber placed directly in front of the
577  //detector "detname". If no ChIo is present, a null pointer is returned.
578 
579  KVINDRADetector* kvd;
580  kvd = dynamic_cast<KVINDRADetector*>(GetDetector(detname));
581  return (kvd ? (KVChIo*)kvd->GetChIo() : NULL);
582 }
583 
584 
585 
586 
590 
592 {
593  //Return pointer to layer in INDRA structure corresponding to ionisation
594  //chambers.
595 
596  return (KVLayer*)GetStructure("LAYER", "CHIO");
597 }
598 
599 
600 
601 
605 
607 {
608  // Define multiplicity trigger used for acquisition and filter.
609  // Events with multipicity >= trig are OK.
610 
611  fTrigger = trig;
612 }
613 
614 
615 
641 
643 {
644  //Find a detector based on the old BaseIndra type definitions:
645  //
646 //enum EBaseIndra_type {
647 // ChIo_GG=1,
648 // ChIo_PG,//=2
649 // ChIo_T,//=3
650 // Si_GG,//=4
651 // Si_PG,//=5
652 // Si_T,//=6
653 // CsI_R,//=7
654 // CsI_L,//=8
655 // CsI_T,//=9
656 // Si75_GG,//=10
657 // Si75_PG,//=11
658 // Si75_T,//=12
659 // SiLi_GG,//=13
660 // SiLi_PG,//=14
661 // SiLi_T//=15
662 //};
663 //enum EBaseIndra_typePhos {
664 // Phos_R=1,
665 // Phos_L,//=2
666 // Phos_T,//=3
667 //};
668 
669  KVINDRADetector* det = 0;
670 
671  char nom_det[20];
672  if (cou == 1) {
673  if (type >= Phos_R && type <= Phos_T) {
674  sprintf(nom_det, "PHOS_%02d", mod);
675  det =
676  (KVINDRADetector*) GetListOfPhoswich()->FindObjectByName(nom_det);
677  }
678  }
679  if (det)
680  return det;
681  else if (type >= ChIo_GG && type <= ChIo_T) {
682  sprintf(nom_det, "CI_%02d%02d", cou, mod);
683  det = (KVINDRADetector*) GetListOfChIo()->FindObject(nom_det);
684  return det;
685  } else if (type >= Si_GG && type <= Si_T) {
686  sprintf(nom_det, "SI_%02d%02d", cou, mod);
687  det = (KVINDRADetector*) GetListOfSi()->FindObject(nom_det);
688  return det;
689  } else if (type >= SiLi_GG && type <= SiLi_T) {
690  sprintf(nom_det, "SILI_%02d", cou);
691  det = (KVINDRADetector*) GetListOfSi()->FindObject(nom_det);
692  return det;
693  } else if (type >= Si75_GG && type <= Si75_T) {
694  sprintf(nom_det, "SI75_%02d", cou);
695  det = (KVINDRADetector*) GetListOfSi()->FindObject(nom_det);
696  return det;
697  } else if (type >= CsI_R && type <= CsI_T) {
698  sprintf(nom_det, "CSI_%02d%02d", cou, mod);
699  det = (KVINDRADetector*) GetListOfCsI()->FindObject(nom_det);
700  return det;
701  }
702  return 0;
703 }
704 
705 
706 
707 
711 
713 {
714  //Override KVASMultiDetArray method for special case of "etalon" modules:
715  //we need to add ChIo-CsI identification telescope by hand
716 
718 
719  if (de->InheritsFrom("KVSiLi") && e->InheritsFrom("KVCsI")) {
720  KVChIo* chio = (KVChIo*)((KVINDRADetector*)e)->GetChIo();
721  if (chio) {
722  n += KVASMultiDetArray::GetIDTelescopes(chio, e, idtels);
723  }
724  }
725  return n;
726 }
727 
728 
729 
734 
736 {
737  // Change default names of ID telescopes to INDRA standard
738  //
739  // This method also sets the types of the ID telescopes
740 
741  TIter it(GetListOfIDTelescopes());
742  KVIDTelescope* idt;
743  const TString etalon_numbers[] = {"1002", "1102", "1202", "1304", "1403", "1503", "1602", "1702"};
744  while ((idt = (KVIDTelescope*)it())) {
745  KVString N = idt->GetDetector(1)->GetName();
746  N.Begin("_");
747  KVString de_type = N.Next();
748  KVString de_number = N.Next();
749  if (idt->GetSize() == 1) {
750  if (de_type == "PHOS") {
751  idt->SetName(Form("%s_R_L_%s", de_type.Data(), de_number.Data()));
752  idt->SetType(Form("%s_R_L", de_type.Data()));
753  } else {
754  // CsI identification telescopes are called either CSI_R_L_RRMM (before 2021)
755  // or CSI_RRMM (for datasets from 2021 onwards)
756  TString csi_id_name_fmt = KVBase::GetDataSetEnv(fDataSet, "INDRA.CSI.IDTelescopeNameFormat", "CSI_%s");
757  idt->SetName(Form(csi_id_name_fmt, de_number.Data()));
758  // similarly, the type of CSI identifications was CSI_R_L before 2021, CSI from 2021 onwards
759  if (csi_id_name_fmt.Contains("R_L")) {
760  idt->SetType("CSI_R_L");
762  } else idt->SetType("CSI");
763  }
764  } else {
765  N = idt->GetDetector(2)->GetName();
766  N.Begin("_");
767  KVString e_type = N.Next();
768  KVString e_number = N.Next();
769  if (e_type == "SILI" || e_type == "SI75") {
770  e_number = etalon_numbers[dynamic_cast<KVINDRADetector*>(idt->GetDetector(2))->GetRingNumber() - 10];
771  }
772  idt->SetName(Form("%s_%s_%s", de_type.Data(), e_type.Data(), e_number.Data()));
773  idt->SetType(Form("%s_%s", de_type.Data(), e_type.Data()));
774  }
775  }
776  // changed names of all id telescope objects in hashlist: must rehash
777  dynamic_cast<KVHashList*>(GetListOfIDTelescopes())->Rehash();
778 }
779 
780 
781 
785 
787 {
788  // Finalise the ROOT geometry description by performing operations which can
789  // only be done once the geometry is closed
790  CreateROOTGeometry();
791 }
792 
793 
794 
796 
797 void KVINDRA::handle_ebyedat_raw_data_parameter(const char* param_name, uint16_t val)
798 {
799  KVString detname;
800  KVDetector* det = nullptr;
801  KVString sig_type;
802  // use look-up table parname => detectorname
803  if (fEbyedatParamDetMap.HasParameter(param_name)) {
804  det = GetDetector(fEbyedatParamDetMap.GetStringValue(param_name));
805  detname = det->GetName();
806  // parameter type given by whatever comes after final '_'
807  std::string lab(param_name);
808  sig_type = lab.substr(lab.rfind('_') + 1);
809  } else { // no detector associated to parameter
810  sig_type = param_name;
811  }
812  add_and_set_detector_signal(det, detname, val, sig_type);
813 }
814 
815 
816 
824 
826 {
827  // Put values of all fired detector signals into the parameter list which will be copied
828  // into the reconstructed event.
829  //
830  // If Mesytec data is being read, we set the parameter "INDRA.MESYTEC", or if Ebyedat data, "INDRA.EBYEDAT".
831  // Note that in old reconstructed data these parameters did not exist, therefore in their
832  // absence EBYEDAT is assumed.
833 
834  if (fMesytecData) fReconParameters.SetValue("INDRA.MESYTEC", kTRUE);
835  else if (fEbyedatData) fReconParameters.SetValue("INDRA.EBYEDAT", kTRUE);
837 }
838 
839 
840 #ifdef WITH_BUILTIN_GRU
841 
845 
847 {
848  // Set raw data in detectors/array coming from a GANIL EBYEDAT format
849  // acquisition file.
850 
851  fEbyedatData = kTRUE;
852  // loop over fired data parameters
853  TIter nxt_frd(&rdr.GetFiredDataParameters());
854  KVEBYEDAT_ACQParam* eby_par;
855  while ((eby_par = (KVEBYEDAT_ACQParam*)nxt_frd())) {
856  handle_ebyedat_raw_data_parameter(eby_par->GetName(), eby_par->GetData());
857  }
858 
859  return kTRUE;
860 }
861 
862 #endif
863 
864 #ifdef WITH_MFM
865 
881 
883 {
884  // General method for reading raw data in MFM-encapsulated ebyedat format
885  // Fills list of hit acquisition parameters.
886  // Returns kTRUE if at least one parameter belonging to the array is present.
887  //
888  // Any unknown parameters in the event (i.e. ones for which no KVACQParam object
889  // has been defined) are written in the fReconParameters list with names
890  // "ACQPAR.[array name].[parameter name]"
891  //
892  // Retrieve CENTRUM timestamp from data if present.
893  // It will be added to fReconParameters as a 64-bit value "INDRA.TS" (if != 0)
894  // Event number is retrieved and stored as "INDRA.EN" (if != 0)
895  // Any parameter which appears as [name] and [name]_UP is an unsigned 32-bit value
896  // split into two 16-bit words. We replace the two parameters with a 64-bit
897  // value (to hold correctly all unsigned 32-bit values) with [name].
898 
899  fEbyedatData = kTRUE;
900 
901  uint16_t val;
902  std::string lab;
903 
904  for (int i = 0; i < f.GetNbItems(); ++i) {
905  f.GetDataItem(i, lab, val);
906  handle_ebyedat_raw_data_parameter(lab.c_str(), val);
907  }
908 
909  ULong64_t ts = f.GetCENTRUMTimestamp();
910  if (ts != 0) fReconParameters.SetValue64bit("INDRA.TS", ts);
911  ULong64_t en = f.GetEventNumber();
912  if (en != 0) fReconParameters.SetValue64bit("INDRA.EN", en);
913  int npars = fReconParameters.GetNpar();
914  std::vector<TString> names;
915  for (int i = 0; i < npars; ++i) {
916  TString name = fReconParameters.GetParameter(i)->GetName();
917  if (name.EndsWith("_UP")) names.push_back(name);
918  }
919  if (names.size()) {
920  for (std::vector<TString>::iterator it = names.begin(); it != names.end(); ++it) {
921  TString name = (*it);
922  name.Remove(name.Index("_UP"), 3);
923  TString name_up = (*it);
924  ULong64_t par_up = fReconParameters.GetIntValue(name_up);
925  ULong64_t par = fReconParameters.GetIntValue(name);
926  UInt_t par32 = (par_up << 16) + par;
927  fReconParameters.RemoveParameter(name_up);
928  fReconParameters.RemoveParameter(name);
929  fReconParameters.SetValue64bit(name, par32);
930  }
931  }
932  return kTRUE;
933 }
934 
935 #endif
936 
937 
940 
942 {
943  // Set the INDRA-specific general identification code for the given telescope
944 
945  if (idt->InheritsFrom(KVIDPhoswich::Class())) idt->SetIDCode(IDCodes::ID_PHOSWICH);
946  else if (idt->InheritsFrom(KVIDINDRACsI::Class())) idt->SetIDCode(IDCodes::ID_CSI_PSA);
947  else if (idt->InheritsFrom(KVIDSiCsI::Class())) idt->SetIDCode(IDCodes::ID_SI_CSI);
948  else if (idt->InheritsFrom(KVIDChIoSi::Class())) idt->SetIDCode(IDCodes::ID_CI_SI);
949  else if (idt->InheritsFrom(KVIDChIoCsI::Class())) idt->SetIDCode(IDCodes::ID_CI_CSI);
950  else if (idt->InheritsFrom(KVIDChIoSi75::Class())) idt->SetIDCode(IDCodes::ID_CI_SI75);
951  else if (idt->InheritsFrom(KVIDSi75SiLi::Class())) idt->SetIDCode(IDCodes::ID_SI75_SILI);
952  else if (idt->InheritsFrom(KVIDChIoSi75::Class())) idt->SetIDCode(IDCodes::ID_CI_SI75);
953  else {
954  Error("SetIDCodeForIDTelescope", "Request for telescope name=%s of unknown class=%s",
955  idt->GetName(), idt->IsA()->GetName());
956  }
957 }
958 
959 
960 
986 
988 {
989  // Sets the KVCsI::fPinLaser member of each CsI detector with the number of the
990  // pin laser associated for the stability control of these detectors.
991  //
992  // We look for a file with the following format:
993  //
994  // CSI_0101 1
995  // CSI_0102 1
996  // CSI_0103 1
997  // CSI_0104 1
998  // etc.
999  //
1000  // i.e. 'name of CsI detector' 'number of pin laser (1-8)'
1001  // Comment lines must begin with '#'
1002  //
1003  // The default name of this file is defined in .kvrootrc by
1004  //
1005  // INDRADB.CsIPinCorr: CsI_PILA.dat
1006  //
1007  // Dataset-specific version can be specified:
1008  //
1009  // INDRA_e999.INDRADB.CsIPinCorr: CorrCsIPin_2054.dat
1010  //
1011  // This file should be in the directory corresponding to the current dataset,
1012  // i.e. in $KVROOT/KVFiles/name_of_dataset
1013 
1014  ifstream pila_file;
1015  if (gDataSet->OpenDataSetFile(gDataSet->GetDataSetEnv("INDRADB.CsIPinCorr", ""), pila_file)) {
1016 
1017  Info("SetPinLasersForCsI", "Setting correspondance CsI-PinLaser using file %s.",
1018  gDataSet->GetDataSetEnv("INDRADB.CsIPinCorr", ""));
1019  // read file, set correspondance
1020  KVString line;
1021  line.ReadLine(pila_file);
1022  while (pila_file.good()) {
1023  if (!line.BeginsWith("#")) {
1024  line.Begin(" ");
1025  KVString detname = line.Next(kTRUE);
1026  KVCsI* det = (KVCsI*)GetDetector(detname.Data());
1027  Int_t pila = line.Next(kTRUE).Atoi();
1028  if (det) {
1029  det->SetPinLaser(pila);
1030  }
1031  }
1032  line.ReadLine(pila_file);
1033  }
1034  pila_file.close();
1035  } else {
1036  Info("SetPinLasersForCsI", "File %s not found. Correspondance Csi-PinLaser is unknown.",
1037  gDataSet->GetDataSetEnv("CsIPinCorr", ""));
1038  }
1039 }
1040 
1041 
1042 
1043 //void KVINDRA::LinkToCodeurs()
1044 //{
1045 
1046 // // Link detectors with electronic modules
1047 // // for the moment only QDC for Si and ChIo are implemented
1048 // // This information is accessible via KVINDRADetector::GetNumeroCodeur()
1049 // // To be active one has to put in the dataset directory
1050 // // a file name Codeurs.dat containing the name of the file for the concerned type
1051 // // of electronic module
1052 // // for example see INDRA_e613 dataset
1053 // // [dataset name].INDRADB.Codeurs: ...
1054 
1055 
1056 // KVFileReader flist;
1057 // TString fp;
1058 // if (!KVBase::SearchKVFile(gDataSet->GetDataSetEnv("INDRADB.Codeurs", ""), fp, gDataSet->GetName())) {
1059 // Warning("LinkToCodeurs", "Fichier %s, inconnu au bataillon", gDataSet->GetDataSetEnv("INDRADB.Codeurs", ""));
1060 // return;
1061 // }
1062 
1063 // if (!flist.OpenFileToRead(fp.Data())) {
1064 // //Error("ReadGainList","Error opening file named %s",fp.Data());
1065 // return;
1066 // }
1067 // Info("LinkToCodeurs()", "Reading correspondance Codeur-Detecteur ...");
1068 
1069 // TEnv* env = 0;
1070 // KVINDRADetector* idet = 0;
1071 // while (flist.IsOK()) {
1072 // flist.ReadLine(NULL);
1073 // KVString file = flist.GetCurrentLine();
1074 // if (file != "") {
1075 // if (KVBase::SearchKVFile(file.Data(), fp, gDataSet->GetName())) {
1076 // env = new TEnv();
1077 // env->ReadFile(fp.Data(), kEnvAll);
1078 // TEnvRec* rec = 0;
1079 // TObjArray* toks = 0;
1080 // TIter it(env->GetTable());
1081 // while ((rec = (TEnvRec*)it.Next())) {
1082 // if (!strcmp(rec->GetName(), "type")) {
1083 // Info("LinkToCodeurs", "Module type %s", rec->GetValue());
1084 // }
1085 // else {
1086 // toks = TString(rec->GetValue()).Tokenize(",");
1087 // for (Int_t ii = 0; ii < toks->GetEntries(); ii += 1) {
1088 // idet = (KVINDRADetector*)gIndra->GetDetector(((TObjString*)toks->At(ii))->GetString().Data());
1089 // if (idet)
1090 // idet->SetNumeroCodeur(TString(rec->GetName()).Atoi());
1091 // }
1092 // delete toks;
1093 // }
1094 // }
1095 // delete env;
1096 // }
1097 // }
1098 // }
1099 
1100 // flist.CloseFile();
1101 
1102 
1103 //}
1104 
1105 
1106 
1107 
1113 
1114 void KVINDRA::GetDetectorEvent(KVDetectorEvent* detev, const TSeqCollection* fired_dets)
1115 {
1116  // Overrides KVASMultiDetArray::GetDetectorEvent.
1117  // If the list of fired detectors is given (meaning we are reading raw data)
1118  // then we check that what we have read is in fact an INDRA event
1119  // (see KVINDRATriggerInfo::IsINDRAEvent()) : if not, we do not try to find the hit groups.
1120 
1121  if (((fired_dets && fired_dets->GetEntries()) || fFiredDetectors.GetEntries())
1122  && (GetTriggerInfo() && !GetTriggerInfo()->IsINDRAEvent())) return;
1123  KVASMultiDetArray::GetDetectorEvent(detev, fired_dets);
1124 }
1125 
1126 
1127 
1128 //void KVINDRA::SetGGtoPGConversionFactors()
1129 //{
1130 // // Sets the parameters for linear conversion of silicon & ChIo coder values
1131 // // between GG and PG, using the following formula:
1132 // //
1133 // // PG = alpha + beta*(GG - GG_0) + PG_0
1134 // //
1135 // // where GG_0 and PG_0 are respectively GG and PG pedestals
1136 // //
1137 // // We look for the file whose name is given by the .kvrootrc variable
1138 // // [dataset].INDRADB.GGtoPGFactors:
1139 // // or by default
1140 // // INDRADB.GGtoPGFactors:
1141 // // and expect to find in it a line for each detector of the form:
1142 // // Det_Name alpha beta
1143 // // Comments in the file can be written on lines beginning with the character '#'
1144 
1145 // ifstream datfile;
1146 // if (!gDataSet->OpenDataSetFile(gDataSet->GetDataSetEnv("INDRADB.GGtoPGFactors", ""), datfile)) {
1147 
1148 // Info("SetGGtoPGConversionFactors", "Cannot open file with parameters for conversion (%s).",
1149 // gDataSet->GetDataSetEnv("INDRADB.GGtoPGFactors", ""));
1150 // return;
1151 // }
1152 // else {
1153 // Info("SetGGtoPGConversionFactors", "Reading parameters from file %s",
1154 // gDataSet->GetDataSetEnv("INDRADB.GGtoPGFactors", ""));
1155 
1156 // Char_t detname[30];
1157 // Double_t a, b;
1158 // TString aline;
1159 // aline.ReadLine(datfile);
1160 // while (datfile.good()) {
1161 
1162 // if (aline[0] != '#') { //skip comments
1163 
1164 // sscanf(aline.Data(), "%s %lf %lf", detname, &a, &b);
1165 // KVINDRADetector* det = (KVINDRADetector*)GetDetector(detname);
1166 // if (!det) {
1167 // //no detector found with cou, mod and type
1168 // Error("SetGGtoPGConversionFactors", "Unknown detector : %s", detname);
1169 // }
1170 // else {
1171 // det->SetGGtoPGConversionFactors(a, b);
1172 // //Info("SetGGtoPGConversionFactors", "%s : PG = %f + %f * GG", detname, a, b);
1173 // }
1174 // }
1175 // aline.ReadLine(datfile);
1176 // } //while( datfile.good()
1177 // datfile.close();
1178 // }
1179 //}
1180 
1181 
1182 
1193 
1195 {
1196  // Overrides KVASMultiDetArray::CreateGeoManager in order to use INDRAGeometryBuilder
1197  // which builds the TGeo representation of INDRA using the Y. Huguet CAO data.
1198  //
1199  // The optional arguments (dx,dy,dz) are the half-lengths in centimetres of the "world"/"top" volume
1200  // into which all the detectors of the array are placed. This should be big enough so that all detectors
1201  // fit in. The default values of 500 give a "world" which is a cube 1000cmx1000cmx1000cm (with sides
1202  // going from -500cm to +500cm on each axis).
1203  //
1204  // If closegeo=kFALSE we leave the geometry open for other structures to be added.
1205 
1206  if (!IsBuilt()) {
1207  Error("CreateROOTGeometry", "gIndra has to be build first");
1208  return;
1209  }
1210  if (!GetNavigator()) {
1211  //Error("CreateROOTGeometry","No existing navigator"); return;
1213  GetNavigator()->SetNameCorrespondanceList("INDRA.names");
1214  }
1215 
1216  // set up shape & matrix pointers in detectors
1217  Info("CreateROOTGeometry", "Scanning geometry shapes and matrices...");
1219  gimp.SetNameCorrespondanceList("INDRA.names");
1220  KVNucleusEvent evt;
1221  KVNucleus* nuc = evt.AddNucleus();
1222  nuc->SetZAandE(1, 1, 1);
1223  KVINDRADetector* det;
1224  TIter next(GetDetectors());
1225  Int_t nrootgeo = 0;
1226  while ((det = (KVINDRADetector*)next())) {
1227  nuc->SetTheta(det->GetTheta());
1228  nuc->SetPhi(det->GetPhi());
1229  gimp.SetLastDetector(0);
1230  gimp.PropagateEvent(&evt);
1231  if (!(det->ROOTGeo())) {
1232  Info("CreateROOTGeometry", "Volume checking for %s", det->GetName());
1233  Double_t theta0 = det->GetTheta();
1234  Double_t phi0 = det->GetPhi();
1235  for (Double_t TH = theta0 - 0.5; TH <= theta0 + 0.5; TH += 0.1) {
1236  for (Double_t PH = phi0 - 10; PH <= phi0 + 10; PH += 1) {
1237  nuc->SetTheta(TH);
1238  nuc->SetPhi(PH);
1239  gimp.SetLastDetector(0);
1240  gimp.PropagateEvent(&evt);
1241  if (det->ROOTGeo()) break;
1242  }
1243  if (det->ROOTGeo()) break;
1244  }
1245  }
1246  if (!(det->ROOTGeo())) {
1247  Info("CreateROOTGeometry", "Volume checking failed for : %s", det->GetName());
1248  }
1249  // check etalon trajectories (if etalons are present)
1250  if (det->ROOTGeo() && det->GetRingNumber() > 9) {
1251  if (GetDetector(Form("SI75_%d", det->GetRingNumber())) || GetDetector(Form("SILI_%d", det->GetRingNumber()))) {
1252  if ((det->IsCalled("CSI_1002") || det->IsCalled("CSI_1102")
1253  || det->IsCalled("CSI_1202") || det->IsCalled("CSI_1304")
1254  || det->IsCalled("CSI_1403") || det->IsCalled("CSI_1503")
1255  || det->IsCalled("CSI_1602") || det->IsCalled("CSI_1702"))
1256  && det->GetNode()->GetNDetsInFront() < 2) {
1257  Info("CreateROOTGeometry", "Trajectory checking for %s", det->GetName());
1258  Double_t theta0 = det->GetTheta();
1259  Double_t phi0 = det->GetPhi();
1260  for (Double_t TH = theta0 - 0.5; TH <= theta0 + 0.5; TH += 0.1) {
1261  for (Double_t PH = phi0 - 10; PH <= phi0 + 10; PH += 1) {
1262  nuc->SetTheta(TH);
1263  nuc->SetPhi(PH);
1264  gimp.SetLastDetector(0);
1265  gimp.PropagateEvent(&evt);
1266  if (det->GetNode()->GetNDetsInFront() == 2) break;
1267  }
1268  if (det->GetNode()->GetNDetsInFront() == 2) break;
1269  }
1270  }
1271  }
1272  }
1273  nrootgeo += (det->ROOTGeo());
1274  }
1275 
1276  Info("CreateROOTGeometry", "ROOT geometry initialised for %d/%d detectors", nrootgeo, GetDetectors()->GetEntries());
1277 
1278  // Set up trajectories
1279  TIter it(GetDetectors());
1280  KVDetector* d;
1281  while ((d = (KVDetector*)it())) d->GetNode()->RehashLists();// make sure detector nodes are correct
1282  AssociateTrajectoriesAndNodes();
1283  DeduceGroupsFromTrajectories();
1284  FillTrajectoryIDTelescopeLists();
1285  CalculateReconstructionTrajectories();
1286  GetNavigator()->AbsorbDetectorPaths(&gimp);
1287 }
1288 
1289 
1290 
1294 
1296 {
1297  // Override base class method
1298  // If ROOT geometry is requested but has not been built, we create it
1299 
1300  if (on) {
1301  CreateGeoManager();
1303  igb.Build(kFALSE, fCloseGeometryNow);
1304  if (fCloseGeometryNow) PerformClosedROOTGeometryOperations();
1305  } else {
1307  }
1308 }
1309 
1310 
1311 
1315 
1317 {
1318  // Set minimum OK multiplicity for (reconstructed) event
1319  // This is the multiplicity trigger used for the current run (if known)
1320 
1321  if (GetTrigger() > 0) e->SetMinimumOKMultiplicity(GetTrigger());
1322 }
1323 
1324 
1325 
1331 
1333 {
1334  // Special INDRA group reconstructors:
1335  // KVINDRAForwardGroupReconstructor rings 1-9
1336  // KVINDRABackwardGroupReconstructor rings 10-17
1337  // KVINDRAEtalonGroupReconstructor for groups with etalon telescopes
1338 
1339  KVGroupReconstructor* gr(nullptr);
1340  if (GetGroup(g->GetName())) { // make sure group belongs to us
1341  // etalons ?
1342  if (g->GetDetectorByType("SILI") || g->GetDetectorByType("SI75"))
1343  gr = KVGroupReconstructor::Factory("INDRA.etalon");
1344  else {
1345  KVINDRADetector* id = (KVINDRADetector*)g->GetDetectors()->First();
1346  if (id->GetRingNumber() < 10) gr = KVGroupReconstructor::Factory("INDRA.forward");
1347  else gr = KVGroupReconstructor::Factory("INDRA.backward");
1348  }
1349  }
1350  return gr;
1351 }
1352 
1353 
1354 
1357 
1359 {
1360  // If "INDRA.EN" parameter has been set, we use it to set the event number
1361 
1363  if (GetReconParameters().HasValue64bit("INDRA.EN")) e->SetNumber(GetReconParameters().GetValue64bit("INDRA.EN"));
1364 }
1365 
1366 
1367 
1373 
1375 {
1376  // Overrides base method in KVMultiDetArray.
1377  //
1378  // If we are reading old reconstructed data with EBYEDAT parameters, we need special
1379  // treatment to decode the detector name & signal type.
1380 
1381  if (!l.GetBoolValue("INDRA.MESYTEC") && !l.GetBoolValue("INDRA.EBYEDAT")) {
1382  // assume we are reading old recon data with parameters in list such as "ACQPAR.INDRA.SI_0201_PG"
1383  fMesytecData = kFALSE;
1384  fEbyedatData = kTRUE;
1385  prepare_to_handle_new_raw_data();
1386 
1387  int N = l.GetNpar();
1388  for (int i = 0; i < N; ++i) {
1389  KVNamedParameter* np = l.GetParameter(i);
1390 
1391  KVString name(np->GetName());
1392  if (name.BeginsWith("ACQPAR")) {
1393  // 2 '.' => 3 values
1394  int dots = name.GetNValues(".");
1395  assert(dots == 3); // sanity check
1396  name.Begin(".");
1397  name.Next(); // "ACQPAR"
1398  if (name.Next() != "INDRA") continue; // check name of array - somebody else's parameter ?
1399  std::string parname(name.Next());
1400  KVString detname;
1401  KVString sig_type;
1402  KVDetector* det = nullptr;
1403  if (parname.find('_') != kNPOS) {
1404  // try to split into '[detname]_[partype]'
1405  detname = parname.substr(0, parname.rfind('_'));
1406  sig_type = parname.substr(parname.rfind('_') + 1);
1407  det = GetDetector(detname);
1408  } else {
1409  sig_type = parname;
1410  }
1411  add_and_set_detector_signal(det, detname, np->GetDouble(), sig_type);
1412  }
1413  }
1414  } else {
1415  fMesytecData = l.GetBoolValue("INDRA.MESYTEC");
1416  fEbyedatData = l.GetBoolValue("INDRA.EBYEDAT");
1418  }
1419 }
1420 
1421 
1422 
1426 
1428 {
1429  // Call this method just after opening a raw data file in order to perform any
1430  // necessary initialisations, depending on the type of data
1431 
1433 
1434  fMesytecData = kFALSE;
1435  fEbyedatData = kFALSE;
1436 #ifdef WITH_MESYTEC
1437 #ifdef WITH_MFM
1438  if (r->GetDataFormat() == "MFM") {
1439  TString crate = gDataSet->GetDataSetEnv("MesytecCrateConf", "");
1440  if (crate != "") {
1441  Info("InitialiseRawDataReading", "Setting Mesytec crate config for run %d, raw run number=%d",
1442  GetCurrentRunNumber(), r->GetRunNumberReadFromFile());
1443 
1444  TString channel_conf_list = gDataSet->GetDataSetEnv("MesytecChannelsConf", "");
1445  if (channel_conf_list != "") {
1446  TEnv conf_list;
1447  if (conf_list.ReadFile(gDataSet->GetFullPathToDataSetFile(channel_conf_list), kEnvUser) == -1) {
1448  Fatal("InitialiseRawDataReading",
1449  "Could not open file containing names of Mesytec channel config for each run: %s",
1450  gDataSet->GetFullPathToDataSetFile(channel_conf_list).Data());
1451  }
1452  KVString files = conf_list.GetValue("File", "");
1453  files.Begin(" ");
1454  while (!files.End()) {
1455  KVString next_file = files.Next(kTRUE);
1456  KVNumberList runlist = conf_list.GetValue(Form("%s.RunList", next_file.Data()), "");
1457  if (runlist.Contains(r->GetRunNumberReadFromFile())) {
1458  // found file for this run
1459  Info("InitialiseRawDataReading", "Using file %s for raw run number %d",
1460  next_file.Data(), r->GetRunNumberReadFromFile());
1461  dynamic_cast<KVMFMDataFileReader*>(r)->InitialiseMesytecConfig(
1463  gDataSet->GetFullPathToDataSetFile(next_file).Data()
1464  );
1465  }
1466  }
1467  } else {
1468  Fatal("InitialiseRawDataReading",
1469  "Name of file containing names of Mesytec channel config for each run should be defined in dataset variable %s.MesytecChannelsConf",
1470  GetDataSet().Data());
1471  }
1472  }
1473  }
1474 #endif
1475 #endif
1476 }
1477 
1478 
1479 #ifdef WITH_MFM
1480 
1487 
1488 #ifdef WITH_MESYTEC
1489 Bool_t KVINDRA::handle_raw_data_event_mfmframe_mesytec_mdpp(const MFMMesytecMDPPFrame& f)
1490 {
1491  // Read a raw data event from a Mesytec MFM Frame.
1492  //
1493  // All data is transferred to KVDetectorSignal objects, either associated to detectors of the
1494  // array (if they can be identified), either associated more globally with the array/event
1495  // itself. The latter are created as needed and go into the fExtraRawDataSignals list.
1496 
1497  if (!fMesytecData) {
1498  // first time we read data, we tweak the parameter names
1499  mesytec::module::set_data_type_alias("qdc_long", "TotLight"); // => CsI
1500  mesytec::module::set_data_type_alias("qdc_short", "R"); // => like old CsI 'R' (rapide) parameter
1501  mesytec::module::set_data_type_alias("tdc", "T"); // => like old 'marqueur de temps' parameters
1502  mesytec::module::set_data_type_alias("adc", "ADC"); // => Si or ChIo signals
1503 
1504  fMesytecData = kTRUE;
1505  }
1506  fReconParameters.SetValue64bit("INDRA.TS", f.GetTimeStamp());
1507  auto mfmfilereader = dynamic_cast<KVMFMDataFileReader*>(fRawDataReader);
1508  mfmfilereader->GetMesytecBufferReader().read_event_in_buffer(
1509  (const uint8_t*)f.GetPointUserData(), f.GetBlobSize(),
1510  [ = ](const mesytec::mdpp::event & evt) {
1511  auto& setup = mfmfilereader->GetMesytecBufferReader().get_setup();
1512  // loop over module data in event, set data in detectors when possible
1513  for (auto& mdat : evt.modules) {
1514  auto mod_id = mdat.module_id;
1515  auto& current_module = setup.get_module(mod_id);
1516  if (current_module.is_mdpp_module()) {
1517  // data for detectors from MDPP module
1518  for (auto& voie : mdat.data) {
1519  KVString detname(setup.get_detector(mod_id, voie.channel));
1520  KVString sig_type(voie.data_type);
1521  Double_t sig_data = voie.data;
1522  add_and_set_detector_signal(GetDetector(detname), detname, sig_data, sig_type);
1523  }
1524  } else if (current_module.is_mvlc_scaler()) {
1525  // data from 64 bit scalers in MVLC
1526  //
1527  // 4 data words of 16 bits (least significant is first word) => 64 bit scaler
1528  if (mdat.data.size() == 4) {
1529  ULong64_t x = 0;
1530  int i = 0;
1531  for (auto& d : mdat.data) x += ((uint64_t)d.data_word) << (16 * (i++));
1532  fReconParameters.SetValue64bit(current_module.name.c_str(), x);
1533  }
1534 // else {
1535 // Warning("handle_raw_data_event_mfmframe_mesytec_mdpp",
1536 // "Scaler data %s is corrupt : %d words instead of 4",
1537 // current_module.name.c_str(), mdat.data.size());
1538 // }
1539  }
1540  }
1541  },
1542  f.GetRevision()
1543  );
1544  return kTRUE;
1545 }
1546 
1547 #endif
1548 #endif
1549 
int Int_t
unsigned int UInt_t
KVDataSet * gDataSet
Definition: KVDataSet.cpp:29
KVINDRA * gIndra
Definition: KVINDRA.cpp:88
@ SiLi_GG
Definition: KVINDRA.h:54
@ Si75_GG
Definition: KVINDRA.h:51
@ ChIo_T
Definition: KVINDRA.h:44
@ SiLi_T
Definition: KVINDRA.h:56
@ CsI_T
Definition: KVINDRA.h:50
@ ChIo_GG
Definition: KVINDRA.h:42
@ CsI_R
Definition: KVINDRA.h:48
@ Si75_T
Definition: KVINDRA.h:53
@ Si_GG
Definition: KVINDRA.h:45
@ Si_T
Definition: KVINDRA.h:47
@ Phos_T
Definition: KVINDRA.h:61
@ Phos_R
Definition: KVINDRA.h:59
ClassImp(KVPartitionList) void KVPartitionList
Initialisation.
ROOT::R::TRInterface & r
#define d(i)
#define f(i)
#define e(i)
const Ssiz_t kNPOS
unsigned char UChar_t
char Char_t
const Bool_t kFALSE
bool Bool_t
double Double_t
const Bool_t kTRUE
kEnvUser
kEnvChange
kEnvAll
#define N
int type
R__EXTERN TGeoManager * gGeoManager
char * Form(const char *fmt,...)
Build INDRA geometry from Huguet CAO infos.
void Build(Bool_t withTarget=kTRUE, Bool_t closeGeometry=kTRUE)
virtual void SetNumber(UInt_t num)
Definition: KVBase.h:215
virtual void SetType(const Char_t *str)
Definition: KVBase.h:172
static const Char_t * GetDataSetEnv(const Char_t *dataset, const Char_t *type, const Char_t *defval)
Definition: KVBase.cpp:1619
virtual Bool_t IsCalled(const Char_t *name) const
Definition: KVBase.h:189
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
void SetPinLaser(Int_t n)
Definition: KVCsI.h:33
const Char_t * GetDataSetEnv(const Char_t *type, const Char_t *defval="") const
Definition: KVDataSet.cpp:758
TString GetFullPathToDataSetFile(const Char_t *filename)
Definition: KVDataSet.cpp:1887
Bool_t OpenDataSetFile(const Char_t *filename, std::ifstream &file)
Definition: KVDataSet.cpp:1850
Base class for detector geometry description.
Definition: KVDetector.h:159
void SetNameOfArray(const TString &n)
Definition: KVDetector.h:775
static KVDetector * MakeDetector(const Char_t *name, Float_t thick)
KVList * GetAlignedIDTelescopes()
Definition: KVDetector.cpp:736
Double_t GetTheta() const
Definition: KVDetector.h:743
Double_t GetPhi() const
Definition: KVDetector.h:761
KVGeoDetectorNode * GetNode()
Definition: KVDetector.h:325
GANIL VXI/VME 16 bit (maximum) EBYEDAT acquisition parameter.
UShort_t GetData() const
Abstract base class container for multi-particle events.
Definition: KVEvent.h:66
Reads GANIL acquisition files (EBYEDAT)
const KVSeqCollection & GetFiredDataParameters() const
Path taken by particles through multidetector geometry.
Bool_t ContainsAll(const TCollection *l) const
KVSeqCollection * AccessIDTelescopeList()
KVGeoDetectorNode * GetNodeInFront(const KVGeoDetectorNode *n) const
KVGeoDetectorNode * GetNodeAt(Int_t i) const
KVDetector * GetDetector() const
Int_t GetNDetsInFront() const
Returns number of detectors directly in front of this one.
Import detector array described by ROOT geometry and set up corresponding KVMultiDetArray object.
Definition: KVGeoImport.h:67
void SetLastDetector(KVDetector *)
void PropagateEvent(KVEvent *, TVector3 *TheOrigin=0)
void SetNameCorrespondanceList(const Char_t *)
virtual void Add(KVBase *)
Base class for particle reconstruction in one group of a detector array.
static KVGroupReconstructor * Factory(const TString &plugin="")
Group of detectors which can be treated independently of all others in array.
Definition: KVGroup.h:19
Extended version of ROOT THashList.
Definition: KVHashList.h:28
Base class for all detectors or associations of detectors in array which can identify charged particl...
Definition: KVIDTelescope.h:83
KVDetector * GetDetector(UInt_t n) const
virtual void SetIDCode(UShort_t c)
const KVList * GetDetectors() const
UInt_t GetSize() const
Base class for detectors of INDRA array.
UInt_t GetRingNumber() const
KVINDRADetector * GetChIo() const
Base class for telescopes in INDRA array.
INDRA multidetector array geometry.
Definition: KVINDRA.h:70
virtual void Build(Int_t run=-1)
Definition: KVINDRA.cpp:389
void InitialiseRawDataReading(KVRawDataReader *)
Definition: KVINDRA.cpp:1427
KVINDRA()
Definition: KVINDRA.cpp:94
virtual KVChIo * GetChIoOf(const Char_t *detname)
Definition: KVINDRA.cpp:573
void handle_ebyedat_raw_data_parameter(const char *param_name, uint16_t val)
Definition: KVINDRA.cpp:797
KVINDRATelescope * BuildTelescope(const Char_t *prefix, Int_t mod)
Definition: KVINDRA.cpp:303
virtual Bool_t handle_raw_data_event_ebyedat(KVGANILDataReader &)
Definition: KVINDRA.cpp:846
void CreateROOTGeometry()
Definition: KVINDRA.cpp:1194
void SetRawDataFromReconEvent(KVNameValueList &)
Definition: KVINDRA.cpp:1374
void SetPinLasersForCsI()
Definition: KVINDRA.cpp:987
void copy_fired_parameters_to_recon_param_list()
Definition: KVINDRA.cpp:825
KVGroupReconstructor * GetReconstructorForGroup(const KVGroup *) const
Definition: KVINDRA.cpp:1332
virtual void SetGroupsAndIDTelescopes()
Definition: KVINDRA.cpp:555
Int_t GetIDTelescopes(KVDetector *, KVDetector *, TCollection *)
Definition: KVINDRA.cpp:712
KVLayer * GetChIoLayer()
Definition: KVINDRA.cpp:591
virtual void SetROOTGeometry(Bool_t on=kTRUE)
Definition: KVINDRA.cpp:1295
static Char_t SignalTypes[16][3]
Use this static array to translate EBaseIndra_type signal type to a string giving the signal type.
Definition: KVINDRA.h:73
void FillListsOfDetectorsByType()
Fill lists of ChIo, Si, CsI and phoswich.
Definition: KVINDRA.cpp:519
virtual void MakeListOfDetectors()
Overrides KVASMultiDetArray method to add FillListsOfDetectorsByType()
Definition: KVINDRA.cpp:506
void FillTrajectoryIDTelescopeLists()
Definition: KVINDRA.cpp:350
void SetIDCodeForIDTelescope(KVIDTelescope *) const
Set the INDRA-specific general identification code for the given telescope.
Definition: KVINDRA.cpp:941
void PerformClosedROOTGeometryOperations()
Definition: KVINDRA.cpp:786
void SetReconParametersInEvent(KVReconstructedEvent *) const
If "INDRA.EN" parameter has been set, we use it to set the event number.
Definition: KVINDRA.cpp:1358
void SetTrigger(UChar_t trig)
Definition: KVINDRA.cpp:606
Bool_t handle_raw_data_event_mfmframe_ebyedat(const MFMEbyedatFrame &)
Definition: KVINDRA.cpp:882
KVRing * BuildRing(Int_t number, const Char_t *prefix)
Build ring with infos in file "$KVROOT/KVFiles/data/indra-struct.[dataset].env".
Definition: KVINDRA.cpp:254
virtual KVINDRADetector * GetDetectorByType(UInt_t cou, UInt_t mod, UInt_t type) const
Definition: KVINDRA.cpp:642
void SetNamesOfIDTelescopes() const
Definition: KVINDRA.cpp:735
virtual void GetDetectorEvent(KVDetectorEvent *detev, const TSeqCollection *fired_dets=0)
Definition: KVINDRA.cpp:1114
void BuildLayer(const Char_t *name)
Build layer 'name' with infos in file "$KVROOT/KVFiles/data/indra-struct.[dataset]....
Definition: KVINDRA.cpp:229
virtual void BuildGeometry()
Definition: KVINDRA.cpp:169
void SetMinimumOKMultiplicity(KVEvent *) const
Definition: KVINDRA.cpp:1316
Set of detectors at a similar distance from target (obsolete)
Definition: KVLayer.h:32
Read MFM format acquisition data.
static KVIonRangeTable * GetRangeTable()
Definition: KVMaterial.cpp:166
virtual void GetDetectorEvent(KVDetectorEvent *detev, const TSeqCollection *fired_params=0)
virtual void InitialiseRawDataReading(KVRawDataReader *)
virtual void copy_fired_parameters_to_recon_param_list()
virtual Int_t GetIDTelescopes(KVDetector *, KVDetector *, TCollection *list)
virtual void SetRawDataFromReconEvent(KVNameValueList &)
virtual void SetROOTGeometry(Bool_t on=kTRUE)
virtual void SetReconParametersInEvent(KVReconstructedEvent *) const
Copy any parameters in fReconParameters in to the reconstructed event parameter list.
Handles lists of named parameters with different types, a list of KVNamedParameter objects.
A generic named parameter storing values of different types.
Double_t GetDouble() const
An event container for KVNucleus objects.
Description of properties and kinematics of atomic nuclei.
Definition: KVNucleus.h:125
void SetZAandE(Int_t z, Int_t a, Double_t ekin)
Set atomic number, mass number, and kinetic energy in MeV.
Definition: KVNucleus.cpp:733
Strings used to represent a set of ranges of values.
Definition: KVNumberList.h:83
Bool_t Contains(Int_t val) const
returns kTRUE if the value 'val' is contained in the ranges defined by the number list
Bool_t End(void) const
Definition: KVNumberList.h:197
void Begin(void) const
Int_t Next(void) const
void SetTheta(Double_t theta)
Definition: KVParticle.h:695
void SetPhi(Double_t phi)
Definition: KVParticle.h:699
virtual void SetAzimuthalMinMax(Double_t min, Double_t max)
Set min and max azimuthal angles and calculate (mean) phi.
Definition: KVPosition.cpp:216
Bool_t ROOTGeo() const
Returns kTRUE if ROOT geometry is used, kFALSE if not.
Definition: KVPosition.cpp:598
void SetDistance(Double_t d)
Definition: KVPosition.h:185
Double_t GetAzimuthalWidth(Double_t phmin=-1., Double_t phimax=-1.) const
Definition: KVPosition.cpp:612
void SetPhi(Double_t p)
Definition: KVPosition.h:181
virtual void SetPolarMinMax(Double_t min, Double_t max)
Set min and max polar angles and calculate (mean) theta.
Definition: KVPosition.cpp:171
Propagate particles through array geometry calculating energy losses.
Abstract base class for reading raw (DAQ) data.
Event containing KVReconstructedNucleus nuclei reconstructed from hits in detectors.
Ring in INDRA array (obsolete)
Definition: KVRing.h:19
void Add(KVBase *)
Only KVTelescope-derived structures can be placed in a KVRing.
Definition: KVRing.cpp:59
virtual void Add(TObject *obj)
Extension of ROOT TString class which allows backwards compatibility with ROOT v3....
Definition: KVString.h:72
void Begin(TString delim) const
Definition: KVString.cpp:565
Bool_t End() const
Definition: KVString.cpp:634
KVString Next(Bool_t strip_whitespace=kFALSE) const
Definition: KVString.cpp:695
Int_t GetNValues(TString delim) const
Definition: KVString.cpp:886
void SetDepth(UInt_t ndet, Float_t depth)
set the depth of detector number ndet(=1,2,...) in mm.
void Add(KVBase *element)
Definition: KVTelescope.cpp:65
virtual void SetPolarMinMax(Double_t min, Double_t max)
Set min and max polar angles and calculate (mean) theta.
Definition: KVTelescope.h:99
virtual void SetAzimuthalWidth(Double_t aw)
Definition: KVTelescope.h:104
virtual Int_t GetEntries() const
virtual const char * GetValue(const char *name, const char *dflt) const
virtual Int_t ReadFile(const char *fname, EEnvLevel level)
virtual const char * GetName() const
virtual void SetName(const char *name)
virtual TObject * FindObject(const char *name) const
virtual Bool_t InheritsFrom(const char *classname) const
Bool_t EndsWith(const char *pat, ECaseCompare cmp=kExact) const
Double_t Atof() const
Bool_t BeginsWith(const char *s, ECaseCompare cmp=kExact) const
const char * Data() const
TString & Prepend(char c, Ssiz_t rep=1)
TString & Remove(EStripType s, char c)
Bool_t Contains(const char *pat, ECaseCompare cmp=kExact) const
Ssiz_t Index(const char *pat, Ssiz_t i=0, ECaseCompare cmp=kExact) const
unsigned long long ULong64_t
TLine * line
Double_t x[n]
const Int_t n
gr SetName("gr")
TGraphErrors * gr
const long double g
masses
Definition: KVUnits.h:72
void Add(RHist< DIMENSIONS, PRECISION, STAT_TO... > &to, const RHist< DIMENSIONS, PRECISION, STAT_FROM... > &from)
double dist(AxisAngle const &r1, AxisAngle const &r2)
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),...)
auto * l