KaliVeda
Toolkit for HIC analysis
KVFAZIA.cpp
1 //Created by KVClassFactory on Tue Jan 27 11:37:39 2015
2 //Author: ,,,
3 
4 #include "KVFAZIA.h"
5 #include "KVGeoImport.h"
6 #include "KVSignal.h"
7 #include "KVFAZIADetector.h"
8 #include "KVGroup.h"
9 #include "KVFAZIABlock.h"
10 #include "KVDetectorEvent.h"
11 #include "KVTarget.h"
12 #include "TSystem.h"
13 #include "KVDataSet.h"
14 #include "KVConfig.h"
15 #include "KVFAZIAIDSiPSA.h"
16 #include "KVFAZIAIDCsI.h"
17 #include "KVFAZIAIDSiCsI.h"
18 #include "KVFAZIAIDSiSi.h"
19 
20 #include "TGeoCompositeShape.h"
21 #include "TGeoEltu.h"
22 
23 #include <KVFAZIAIDSiSiCsI.h>
24 #include <KVReconstructedNucleus.h>
25 
26 #ifdef WITH_MFM
27 #include "MFMFaziaFrame.h"
28 #endif
29 
30 #ifdef WITH_PROTOBUF
31 #include "FzEventSet.pb.h"
32 #include "KVFzDataReader.h"
33 #endif
34 
36 
37 // BEGIN_HTML <!--
39 /* -->
40 <h2>KVFAZIA</h2>
41 <h4>Base class for description of the FAZIA set up</h4>
42 <!-- */
43 // --> END_HTML
45 KVFAZIA* gFazia;
46 
47 static Char_t const* const FzDataType_str[] = { "QH1", "I1", "QL1", "Q2", "I2", "Q3", "ADC", "UNK" };
48 static Char_t const* const FzDetector_str[] = { "SI1", "SI1", "SI1", "SI2", "SI2", "CSI" };
49 
50 
53 
54 KVFAZIA::KVFAZIA(const Char_t* title)
55  : KVMultiDetArray("FAZIA", title)
56 {
57  // Default constructor
59  gFazia = this;
60  fDetectorLabels = "";
61  fSignalTypes = "QL1,I1,QH1,Q2,I2,Q3";
64 
65  // values of trapezoidal filter rise time set in the fpgas to be linked with a database...
66  fQH1risetime = GetSetupParameter("QH1.FPGARiseTime");
67  fQ2risetime = GetSetupParameter("Q2.FPGARiseTime");
68  fQ3slowrisetime = GetSetupParameter("Q3.slow.FPGARiseTime");
69  fQ3fastrisetime = GetSetupParameter("Q3.fast.FPGARiseTime");
70 
71  Info("KVFAZIA", "fpga shapers: %lf %lf %lf %lf", fQH1risetime, fQ2risetime, fQ3slowrisetime, fQ3fastrisetime);
72 }
73 
74 
75 
77 
79 {
80 
81  Double_t lval = -1;
82  if (gDataSet) lval = gDataSet->GetDataSetEnv(parname, 0.0);
83  else lval = gEnv->GetValue(parname, 0.0);
84  return lval;
85 }
86 
87 
88 
91 
93 {
94  // Overrides base method in order to set the value of the trigger bit pattern for the event
96  if (l.HasIntParameter("FAZIA.TRIGPAT")) SetTriggerPattern(l.GetIntValue("FAZIA.TRIGPAT"));
97  else SetTriggerPattern(0);// FAZIA absent from event
98 }
99 
100 
101 
104 
106 {
107  // Override base method in order to read FAZIA trigger for each run
108 
111 }
112 
113 
114 
128 
130 {
131  // Returns the symbolic name for the principal DAQ trigger used for the current run
132  // e.g. 'Mult2', 'Mult1/100', etc. (see SetTriggerPatternsForDataSet()).
133  //
134  // This can be used to test if the actual DAQ trigger for an event was consistent
135  // with the principal trigger by doing:
136  //
137  //~~~~{.cpp}
138  // if( gFazia->GetTrigger().IsTrigger( gFazia->GetTriggerForCurrentRun() ) )
139  // {
140  // ===ok in this case trigger is consistent==
141  // }
142  //~~~~
143 
144  if (gExpDB) {
145  if (gExpDB->GetTable("FAZIA.Triggers")) {
146  auto rundb = gExpDB->GetDBRun(GetCurrentRunNumber());
147  if (rundb) {
148  auto links = rundb->GetLinks("FAZIA.Triggers");
149  if (links && links->GetEntries())
150  return links->First()->GetName();
151  }
152  }
153  }
154  return "";
155 }
156 
157 
158 
166 
168 {
169  // Set identification telescope objects for FAZIA geometry by hand.
170  //
171  // - CSI
172  // - SI1-SI2-CSI (name: ID_SI_CSI_xxxx)
173  // - SI1-SI2
174  // - SI1
175 
177  TIter next_traj(GetTrajectories());
178  KVGeoDNTrajectory* traj;
179  while ((traj = (KVGeoDNTrajectory*)next_traj())) { // loop over all trajectories
180 
181  traj->IterateBackFrom(); // from closest-in to furthest-out detector
182 
184  KVDetector* csi{nullptr}, *si2{nullptr}, *si1{nullptr};
185  while ((N = traj->GetNextNode())) {
186  auto d = N->GetDetector();
187  if (d->IsLabelled("SI1") && d->IsOK()) si1 = d;
188  else if (d->IsLabelled("SI2") && d->IsOK()) si2 = d;
189  else if (d->IsLabelled("CSI") && d->IsOK()) csi = d;
190  }
191 
192  auto add_telescope = [ = ](KVIDTelescope * idt, KVGroup * g) {
193  idt->SetGroup(g);
194  traj->AccessIDTelescopeList()->Add(idt);
195  fIDTelescopes->Add(idt);
196  };
197 
198  if (csi) {
199  // add CSI, SI2-CSI and SI1-SI2-CSI telescopes
200  KVIDTelescope* idt = new KVFAZIAIDCsI;
201  idt->AddDetector(csi);
202  add_telescope(idt, csi->GetGroup());
203  if (si2) {
204  idt = new KVFAZIAIDSiCsI;
205  idt->AddDetector(si2);
206  idt->AddDetector(csi);
207  add_telescope(idt, csi->GetGroup());
208  if (si1) {
209  idt = new KVFAZIAIDSiSiCsI;
210  idt->AddDetector(si1);
211  idt->AddDetector(si2);
212  idt->AddDetector(csi);
213  idt->SetName(Form("ID_SI_CSI_%d", csi->GetIndex()));
214  add_telescope(idt, csi->GetGroup());
215  }
216  }
217  }
218  if (si2) {
219  auto idt = new KVFAZIAIDSiPSA;
220  idt->AddDetector(si2);
221  add_telescope(idt, si2->GetGroup());
222  }
223  if (si1 && si2) {
224  auto idt = new KVFAZIAIDSiSi;
225  idt->AddDetector(si1);
226  idt->AddDetector(si2);
227  add_telescope(idt, si2->GetGroup());
228  }
229  if (si1) {
230  auto idt = new KVFAZIAIDSiPSA;
231  idt->AddDetector(si1);
232  add_telescope(idt, si1->GetGroup());
233  }
234  }
235 }
236 
237 
238 
239 
242 
244 {
245  // Destructor
246 
247  if (gFazia == this) gFazia = nullptr;
248 }
249 
250 
251 
253 
255 {
256  if (fDetectorLabels == "") fDetectorLabels += label;
257  else if (!fDetectorLabels.Contains(label)) fDetectorLabels += Form(",%s", label);
258 }
259 
260 
261 
265 
267 {
268  // Look for the geometry object <-> detector name correspondance file in the dataset directory
269  // If not found, we create it
270 
271  fCorrespondanceFile = "";
272  if (gDataSet) fCorrespondanceFile = gDataSet->GetFullPathToDataSetFile(Form("%s.names", ClassName()));
273  if (fCorrespondanceFile != "") return;
274 
275 #ifdef WITH_GNU_INSTALL
277 #else
278  fCorrespondanceFile.Form("%s/%s.names", gDataSet->GetDataSetDir(), ClassName());
279 #endif
280  Info("GenerateCorrespondanceFile", "Creation de %s", fCorrespondanceFile.Data());
281  KVEnv env;
282 
283  fDetectorLabels = "SI1,SI2,CSI";
284 
285  SetNameOfDetectors(env);
286  if (env.GetTable() && env.GetTable()->GetEntries() > 0) {
287  env.AddCommentLine(Form("Automatic generated file by %s::GenerateCorrespondanceFile", ClassName()));
288  env.AddCommentLine("Make link between geometric ROOT objects and detector names");
290  }
291  fDetectorLabels = "";
292 }
293 
294 
295 
303 
305 {
306  //define the format of detectors name
307  // label-index
308  // where index = block*100+quartet*10+telescope
309  // example :
310  // SI1-123 is the Silicon 1 of the block 1, the quartet 2 and the telescope 3
311  //
312 
313  for (Int_t bb = fStartingBlockNumber; bb < fNblocks; bb += 1) {
314  for (Int_t qq = 1; qq <= 4; qq += 1) {
315  for (Int_t tt = 1; tt <= 4; tt += 1) {
316  fDetectorLabels.Begin(",");
317  while (!fDetectorLabels.End()) {
318  KVString sdet = fDetectorLabels.Next();
319  env.SetValue(
320  Form("BLOCK_%d_QUARTET_%d_%s-T%d", bb, qq, sdet.Data(), tt),
321  Form("%s-%d", sdet.Data(), bb * 100 + qq * 10 + tt)
322  );
323  }
324  }
325  }
326  }
327 }
328 
329 
330 
334 
336 {
337  // Finalise description of array performing all operations which require ROOT
338  // geometry to be closed
339 
342  imp.SetDetectorPlugin(GetDataSetEnv(GetDataSet(), "FAZIADetectorPlugin", "FAZIADetector"));
344  // any additional structure name formatting definitions
346  imp.AddAcceptedDetectorName("SI1-");
347  imp.AddAcceptedDetectorName("SI2-");
348  imp.AddAcceptedDetectorName("CSI-");
349 
350  // the following parameters are optimized for a 12-block compact
351  // geometry placed at 80cm with rings 1-5 of INDRA removed.
352  // make sure that the expected number of detectors get imported!
354 
356 
358  SetBit(kIsBuilt);
359 }
360 
361 
362 
365 
367 {
368  //Called by the Build method
369  AbstractMethod("GetGeometryParameters");
370 }
371 
372 
373 
393 
395 {
396  // Read and set up definitions of trigger patterns for this dataset.
397  // These should be given by variables such as:
398  //
399  //~~~~
400  // +[dataset].FAZIA.TriggerPatterns: [name1]
401  // [dataset].FAZIA.TriggerPattern.[name1]: [value1]
402  // +[dataset].FAZIA.TriggerPatterns: [name2]
403  // [dataset].FAZIA.TriggerPattern.[name2]: [value2]
404  //~~~~
405  //
406  // where [name*]='Mult1','Mult1/100','Mult2', etc. (see KVFAZIATrigger for known trigger patterns).
407  //
408  // and [value*] is the value of the corresponding bit pattern, e.g. if bit '3' (0b100) corresponds to
409  // 'Mult2' (i.e. multiplicity >= 2) then this would give
410  //
411  //~~~~
412  // [dataset].FAZIA.TriggerPattern.Mult2: 4
413  //~~~~
414 
415  KVString patterns = GetDataSetEnv(dataset, "FAZIA.TriggerPatterns", "");
416  if (patterns.Length()) {
417  patterns.Begin(" ");
418  while (!patterns.End()) {
419  auto pattern = patterns.Next(kTRUE);
420  uint16_t val = (uint16_t)GetDataSetEnv(dataset, Form("FAZIA.TriggerPattern.%s", pattern.Data()), 0.);
421  fTrigger.AddTriggerPattern(pattern, val);
422  }
423  }
424 }
425 
426 
427 
430 
432 {
433  //Called by the Build method
434  Info("BuildFAZIA", "to be defined in child class ...");
435 
436 }
437 
438 
439 
441 
443 {
444 
445  KVMaterial target_holder_mat("Al");
446  new TGeoBBox("TARGET_FRAME", 3., 3., 0.1 / 2.);
447  new TGeoEltu("TARGET_HOLE", 2., 2., 0.1 / 2.);
448  TGeoCompositeShape* cs = new TGeoCompositeShape("TARGET_FRAME", "TARGET_FRAME - TARGET_HOLE");
449  TGeoVolume* target_frame = new TGeoVolume("TARGET_FRAME", cs, target_holder_mat.GetGeoMedium());
450  gGeoManager->GetTopVolume()->AddNode(target_frame, 1);
451 
452  KVTarget* T = GetTarget();
453  if (T) {
454  KVMaterial* targMat = (KVMaterial*)T->GetLayers()->First();
455  TGeoVolume* target = gGeoManager->MakeEltu("TARGET", targMat->GetGeoMedium(), 2., 2., targMat->GetThickness() / 2.);
457  }
458 }
459 
460 
461 
464 
466 {
467  // Build the FAZIA array
470 
472 
473  BuildFAZIA();
474 
475  if (fBuildTarget)
476  BuildTarget();
477 
478  if (fCloseGeometryNow) {
481  }
482 
484 }
485 
486 
487 
499 
501 {
502  // First step in event reconstruction based on current status of detectors in array.
503  //
504  // Fills the given KVDetectorEvent with the list of all groups which have fired.
505  // i.e. loop over all groups of the array and test whether KVGroup::Fired() returns true or false.
506  //
507  // This can be made more efficient if the detectors which were hit in the event are already known:
508  // then their list should be given to argument dets
509  //
510  // If the list of fired detectors dets is not given, we use the internal fFiredACQParams list
511  // which is filled with all hit detectors when raw data is treated in treat_event()
512 
513  if (!fHandledRawData) {
514  //Info("GetDetectorEvent","i didnt handle any data...");
515  return;
516  }
517  if (!dets || !dets->GetEntries()) {
518  if (fFiredDetectors.GetEntries()) {
519  dets = &fFiredDetectors;
520  //Info("GetDetectorEvent", "using internal list");
521  }
522  }
523  if (dets && dets->GetEntries()) {
524  TIter next_det(dets);
525 
526  KVDetector* det = 0;
527  while ((det = (KVDetector*)next_det())) {
528 
529  if (det->GetGroup()->Fired()) detev->AddGroup(det->GetGroup());
530 
531  }
532  }
533  else {
534  //Info("GetDetectorEvent", "Calling base method");
536  }
537 }
538 
539 
540 
545 
547 {
548  // Protected method, called when required to fill fDetList with pointers to
549  // the detectors whose names are stored in fDetNames.
550  // Also set all raw data values in the detectors.
551 
552  KVFAZIADetector* det = 0;
553 
554  DetList->Clear();
555  DetNames.Begin("/");
556  while (!DetNames.End()) {
557  KVString sdet = DetNames.Next(kTRUE);
558  det = (KVFAZIADetector*)GetDetector(sdet.Data());
559  if (!det) {
561  }
562 
563  if (det) {
564  DetList->Add(det);
565  // read and set from the particle's parameter list any relevant detector signal values
566  // each signal is stored with a name "[detname].[signal name]"
567  // except GTTag and DetTag which have the same value for all detectors of the same telescope
568  // and so are only stored once with name "DetTag" or "GTTag".
569 
570  TIter it(&det->GetListOfDetectorSignals());
571  KVDetectorSignal* ds;
572  while ((ds = (KVDetectorSignal*)it())) {
573  if (ds->IsRaw() && !ds->IsExpression())
574  // only look for raw data, excluding any expressions based only on raw data
575  {
576  TString pname = Form("%s.%s", det->GetName(), ds->GetName());
577  if (rnuc->GetParameters()->HasParameter(pname))
578  ds->SetValue(rnuc->GetParameters()->GetDoubleValue(pname));
579  }
580  }
581  if (rnuc->GetParameters()->HasParameter("GTTag"))
582  det->SetGTTag(rnuc->GetParameters()->GetIntValue("GTTag"));
583  if (rnuc->GetParameters()->HasParameter("DetTag"))
584  det->SetDetTag(rnuc->GetParameters()->GetIntValue("DetTag"));
585  }
586  }
587 }
588 
589 
590 
593 
595 {
596  // Specialized group reconstructor for FAZIA
597 
598  KVGroupReconstructor* gr(nullptr);
599  if (GetGroup(g->GetName())) { // make sure group belongs to us
601  }
602  return gr;
603 }
604 
605 
606 
608 
610 {
611 
612  TString sname;
613  if (bb == 4) {
614  if (fDataSet == "FAZIASYM") {
615  sname.Form("%s-RUTH", FzDataType_str[idsig]);
616  }
617  else {
618  sname.Form("%s-%d", FzDataType_str[idsig], 100 * bb + 10 * qq + tt);
619  }
620  }
621  else if (bb == 6) {
622  if (fDataSet == "FAZIAPRE") {
623  sname.Form("%s-RUTH", FzDataType_str[idsig]);
624  }
625  else {
626  sname.Form("%s-%d", FzDataType_str[idsig], 100 * bb + 10 * qq + tt);
627  }
628  }
629  else {
630  sname.Form("%s-%d", FzDataType_str[idsig], 100 * bb + 10 * qq + tt);
631  }
632  return sname;
633 
634 }
635 
636 
637 #ifdef WITH_PROTOBUF
638 
640 
642 {
643  return treat_event(((KVFzDataReader&)R).get_fazia_event());
644 }
645 
646 
647 
649 
651 {
652  Int_t value = (val << 2);
653  value >>= 2;
654  Double_t dval = -1.;
655  switch (sigid) {
656  case DAQ::FzData_FzDataType_QH1:
657  if (eid == 0) dval = value / (fQH1risetime * 1e3 / 10.);
658  break;
659  case DAQ::FzData_FzDataType_I1:
660  break;
661  case DAQ::FzData_FzDataType_QL1:
662  break;
663  case DAQ::FzData_FzDataType_Q2:
664  if (eid == 0) dval = value / (fQ2risetime * 1e3 / 10.);
665  break;
666  case DAQ::FzData_FzDataType_I2:
667  break;
668  case DAQ::FzData_FzDataType_Q3:
669  if (eid == 0) dval = value / (fQ3slowrisetime * 1e3 / 10.);
670  if (eid == 1) dval = value / (fQ3fastrisetime * 1e3 / 10.);
671  break;
672  }
673  return dval;
674 }
675 
676 
677 
680 
681 Bool_t KVFAZIA::treat_event(const DAQ::FzEvent& e)
682 {
683  // Read raw data for an event
684 
685  Bool_t good = kTRUE;
686 
687  //get info from trigger
688  int ts = e.trinfo_size();
689  double dt = 0;
690  uint64_t tot = 0;
691  for (Int_t tr = ts - 1; tr >= 0; tr--) {
692  const DAQ::FzTrigInfo& rdtrinfo = e.trinfo(tr);
693  uint64_t triggervalue = rdtrinfo.value();
694  if (tr == ts - 5) {
695  fReconParameters.SetValue("FAZIA.TRIGPAT", (int)triggervalue);
696  SetTriggerPattern((uint16_t)triggervalue);
697  }
698  else if (tr == ts - 6) fReconParameters.SetValue64bit("FAZIA.EC", ((triggervalue << 12) + e.ec()));
699  else if (tr == ts - 8) {
700  dt = triggervalue / (1e6);
701  fReconParameters.SetValue("FAZIA.TRIGRATE.DT", dt);
702  double faziats = fReconParameters.GetValue64bit("FAZIA.TS") * 1.e-8;
703  if (oldfaziats > 0.) fReconParameters.SetValue("FAZIA.CENTRUM.DT", faziats - oldfaziats);
704  oldfaziats = faziats;
705  }
706  else if (tr == ts - 9) fReconParameters.SetValue("FAZIA.TRIGRATE.EXT", 1.*triggervalue / dt);
707  else if (tr == ts - 10) fReconParameters.SetValue("FAZIA.TRIGRATE.MAN", 1.*triggervalue / dt);
708  else if (tr == ts - 11) fReconParameters.SetValue(Form("FAZIA.TRIGRATE.PAT%d", tr - 2), 1.*triggervalue / dt);
709  else if (tr == ts - 12) fReconParameters.SetValue(Form("FAZIA.TRIGRATE.PAT%d", tr - 2), 1.*triggervalue / dt);
710  else if (tr == ts - 13) fReconParameters.SetValue(Form("FAZIA.TRIGRATE.PAT%d", tr - 2), 1.*triggervalue / dt);
711  else if (tr == ts - 14) fReconParameters.SetValue(Form("FAZIA.TRIGRATE.PAT%d", tr - 2), 1.*triggervalue / dt);
712  else if (tr == ts - 15) fReconParameters.SetValue(Form("FAZIA.TRIGRATE.PAT%d", tr - 2), 1.*triggervalue / dt);
713  else if (tr == ts - 16) fReconParameters.SetValue(Form("FAZIA.TRIGRATE.PAT%d", tr - 2), 1.*triggervalue / dt);
714  else if (tr == ts - 17) fReconParameters.SetValue(Form("FAZIA.TRIGRATE.PAT%d", tr - 2), 1.*triggervalue / dt);
715  else if (tr == ts - 18) fReconParameters.SetValue(Form("FAZIA.TRIGRATE.PAT%d", tr - 2), 1.*triggervalue / dt);
716  else if (tr == ts - 19) {
717  fReconParameters.SetValue("FAZIA.TRIGRATE.TOT", 1.*triggervalue / dt);
718  tot = triggervalue;
719  }
720  else if (tr == ts - 20) {
721  fReconParameters.SetValue("FAZIA.TRIGRATE.VAL", 1.*triggervalue / dt);
722  fReconParameters.SetValue("FAZIA.DEADTIME", 100.*(1. - 1.*triggervalue / tot));
723  }
724  else {}
725  }
726 
727 // if (tot > 0) fReconParameters.ls();
728 
729  for (int b = 0; b < e.block_size(); ++b) {
730 
731  // check block errors
732  if (e.block(b).len_error() || e.block(b).crc_error() || (!good)) {
733  //Warning("treat_event", "BLOCK LEN OR CRC ERROR B%03d", e.block(b).blkid());
734  good = kFALSE;
735  break; //stop iteration on blocks
736  }
737  int fIdBlk = e.block(b).blkid();
738 
739  for (int f = 0; f < e.block(b).fee_size(); ++f) {
740 
741  const DAQ::FzFee& rdfee = e.block(b).fee(f);
742 
743  for (int h = 0; h < rdfee.hit_size(); ++h) {
744 
745  const DAQ::FzHit& rdhit = rdfee.hit(h);
746  // check fee errors
747  if (rdfee.len_error() || rdfee.crc_error() || (!good)) {
748  Warning("treat_event", "FEE LEN OR CRC ERROR B%03d-FE%d", e.block(b).blkid(), rdfee.feeid());
749  good = kFALSE;
750  break; //stop iteration on hits
751  }
752  int fIdFee = rdhit.feeid();
753  int fIdTel = rdhit.telid();
754 
755  for (Int_t mm = 0; mm < rdhit.data_size(); mm++) {
756  const DAQ::FzData& rdata = rdhit.data(mm);
757  int fIdSignal = rdata.type();
758 
759  int DetTag = rdhit.dettag();
760  int GTTag = rdhit.gttag();
761  if (DetTag >= 16384 && GTTag < 16384) GTTag += 32768;
762 
763  //on decompile le HIT
764  int fIdQuartet = fQuartet[fIdFee][fIdTel];
765  int fIdTelescope = fTelescope[fIdFee][fIdTel];
766 
767  KVFAZIADetector* det = (KVFAZIADetector*)GetDetector(Form("%s-%d", FzDetector_str[fIdSignal], 100 * fIdBlk + 10 * fIdQuartet + fIdTelescope));
768  if (!det) {
769  Error("treat_event", "No detector %s-%d found in FAZIA geometry...", FzDetector_str[fIdSignal], 100 * fIdBlk + 10 * fIdQuartet + fIdTelescope);
770  continue;
771  }
772  det->SetDetTag(DetTag);
773  det->SetGTTag(GTTag);
774 
775  if (!rdata.has_energy() && !rdata.has_waveform()) {
776  if (FzDataType_str[fIdSignal] != "I2") Warning("treat_event", "[NO DATA] [%s %s]", det->GetName(), FzDataType_str[fIdSignal]);
777  continue;
778  }
779 
780  if (rdata.has_energy()) {
781  const DAQ::Energy& ren = rdata.energy();
782  for (Int_t ee = 0; ee < ren.value_size(); ee++) {
783  Double_t energy = TreatEnergy(fIdSignal, ee, ren.value(ee));
784  fFiredSignals.Add(det->SetFPGAEnergy(fIdSignal, ee, energy));
785  }
786  fFiredDetectors.Add(det);
787  }
788  if (rdata.has_baseline()) {
789  Float_t bl = rdata.baseline();
790  fFiredSignals.Add(det->SetBaseLine(fIdSignal, bl));
791  }
792  if (rdata.has_waveform()) {
793  const DAQ::Waveform& rwf = rdata.waveform();
794  Int_t supp;
795 
796  if (fIdSignal <= 5) {
797  TString sname = GetSignalName(fIdBlk, fIdQuartet, fIdTelescope, fIdSignal);//QH1-123 etc.
798  if (sname == "")
799  Warning("treat_event", "signal name is empty !!! blk=%d qua=%d tel=%d\n", fIdBlk, fIdQuartet, fIdTelescope);
800 
801  TGraph sig(rwf.sample_size());
802 
803  for (Int_t nn = 0; nn < rwf.sample_size(); nn++) {
804  if (fIdSignal != DAQ::FzData::ADC) {
805  if (rwf.sample(nn) > 8191) {
806  supp = rwf.sample(nn) | 0xFFFFC000;
807  }
808  else {
809  supp = rwf.sample(nn);
810  }
811  }
812  else {
813  supp = rwf.sample(nn);
814  }
815  sig.SetPoint(nn, nn, supp);
816  }
817  det->SetSignal(&sig, sname);
818  }
819  else {
820  if (fIdSignal > 5)
821  Warning("treat_event", "datatype %d>5 - taille = %d\n", fIdSignal, rwf.sample_size());
822  }
823  }
824  }
825  }
826  }
827  }
828 
829 // cout << "good=" << good << endl;
830 // fFPGAParameters.ls();
831 // fSignals.ls();
832 
833  return good;
834 }
835 
836 #endif
837 
838 #ifdef WITH_MFM
839 
844 
846 {
847  // Treatment of raw data in MFM frames with type MFM_FAZIA_FRAME_TYPE
848  // The timestamp is extracted from the frame header and added to fReconParameters
849  // in a 64 bit parameter with name "FAZIA.TS"
850 
851  if (f.GetFrameType() != MFM_FAZIA_FRAME_TYPE) return kFALSE;
852  fReconParameters.SetValue64bit("FAZIA.TS", f.GetTimeStamp());
853 
854 #ifdef WITH_PROTOBUF
855  DAQ::FzEventSet fazia_set;
856  DAQ::FzEvent fazia_event;
857  // Parse protobuf data in MFM frame
858  if (fazia_set.ParseFromArray(f.GetPointUserData(), ((MFMFaziaFrame&)f).GetEventSize())) {
859  // Parsed an event set
860  if (fazia_set.ev_size() > 1) {
861  Warning("handle_raw_data_event_mfmframe",
862  "Got a FzEventSet from data: cannot handle multiple events at once!");
863  return kFALSE;
864  }
865  return treat_event(fazia_set.ev(0));
866  }
867  else if (fazia_event.ParseFromArray(f.GetPointUserData(), ((MFMFaziaFrame&)f).GetEventSize())) {
868  // Parsed an event
869  return treat_event(fazia_event);
870  }
871 #endif
872  return kTRUE;
873 }
874 
875 #endif
876 
877 
881 
883 {
884  // set up correspondence between FPGA number/FEE number (from acquisition)
885  // and Quartet/Telescope numbers
886 
887  TString DataFilePath;
888  if (!KVBase::SearchKVFile("ElecDetLink.env", DataFilePath, "data")) {
889  Error("CreateCorrespondence", "ElecDetLink.env not found");
890  return;
891  }
892  KVEnv DetLink;
893  DetLink.ReadFile(DataFilePath, kEnvUser);
894  for (int t = 1; t <= 4; t++) {
895  for (int q = 1; q <= 4; q++) {
896  TString elec = DetLink.GetValue(Form("T%1d-Q%1d", t, q), " ");
897  if (!elec.IsWhitespace()) {
898  int fee, fpga;
899  sscanf(elec.Data(), "FPGA%d-FE%d", &fpga, &fee);
900  fQuartet[fee][fpga] = q;
901  fTelescope[fee][fpga] = t;
902  }
903  else {
904  Error("CreateCorrespondence", "Problem reading FAZIA ElecDetLink.env file : T%1d-Q%1d = %s", t, q, elec.Data());
905  }
906  }
907  }
908 }
909 
910 
911 
926 
928 {
929  // Read a file containing runlists for each principal trigger used during an experiment
930  //
931  // The file should be in TEnv format like so:
932  //
933  //~~~~
934  // Mult1: 100-122,541-1938
935  // Mult2: 91-765
936  //~~~~
937  //
938  // where each trigger pattern name must be known and declared to occur during the dataset
939  // (see SetTriggerPatternsForDataSet()) and the list of runs are given using KVNumberList syntax.
940  //
941  // The data is added to the database in a table 'FAZIA.Triggers'.
942 
943  TString fullpath;
944  if (!db->FindCalibFile("Triggers", fullpath)) return;
945 
946  Info("ReadTriggerPatterns()", "Reading FAZIA triggers used during runs...");
947  auto trigs = db->AddTable("FAZIA.Triggers", "Principal triggers used by FAZIA");
948 
949  KVDBRecord* dbrec = 0;
950  TEnv env;
951  TEnvRec* rec = 0;
952  env.ReadFile(fullpath.Data(), kEnvAll);
953  TIter it(env.GetTable());
954 
955  while ((rec = (TEnvRec*)it.Next())) {
956  KVString srec(rec->GetName());
957  KVNumberList nl(rec->GetValue());
958  dbrec = new KVDBRecord(rec->GetName(), "FAZIA Trigger");
959  dbrec->AddKey("Runs", "List of Runs");
960  trigs->AddRecord(dbrec);
961  db->LinkRecordToRunRange(dbrec, nl);
962  }
963 }
964 
965 
966 
969 
971 {
972  // Set the FAZIA-specific general identification code for the given telescope
973 
974  if (idt->InheritsFrom(KVFAZIAIDSiPSA::Class())) {
975  if (idt->GetDetector(1)->IsLabelled("SI1")) idt->SetIDCode(IDCodes::ID_SI1_PSA);
976  else idt->SetIDCode(IDCodes::ID_SI2_PSA);
977  }
978  else if (idt->InheritsFrom(KVFAZIAIDSiSi::Class())) idt->SetIDCode(IDCodes::ID_SI1_SI2);
979  else if (idt->InheritsFrom(KVFAZIAIDSiCsI::Class())) idt->SetIDCode(IDCodes::ID_SI2_CSI);
980  else if (idt->InheritsFrom(KVFAZIAIDSiSiCsI::Class())) {
981  // ID code for these telescopes depends on what detectors the grid's VARY uses
982  auto labsY = idt->GetDetectorLabelsForGridCoord("y");
983  if (labsY.Contains("SI1") && labsY.Contains("SI2"))
984  idt->SetIDCode(IDCodes::ID_SI12_CSI);
985  else if (labsY.Contains("SI1"))
986  idt->SetIDCode(IDCodes::ID_SI1_CSI);
987  else if (labsY.Contains("SI2"))
988  idt->SetIDCode(IDCodes::ID_SI2_CSI);
989  }
990  else if (idt->InheritsFrom(KVFAZIAIDCsI::Class())) idt->SetIDCode(IDCodes::ID_CSI_PSA);
991  else {
992  Error("SetIDCodeForIDTelescope", "Request for telescope name=%s of unknown class=%s",
993  idt->GetName(), idt->IsA()->GetName());
994  }
995 }
996 
997 
int Int_t
unsigned int UInt_t
#define d(i)
#define f(i)
#define e(i)
bool Bool_t
char Char_t
float Float_t
constexpr Bool_t kFALSE
double Double_t
constexpr Bool_t kTRUE
R__EXTERN TEnv * gEnv
kEnvUser
kEnvAll
#define N
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t b
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t Int_t Int_t Window_t TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t winding char text const char depth char const char Int_t count const char ColorStruct_t color const char Pixmap_t Pixmap_t PictureAttributes_t attr const char char ret_data h unsigned char height h Atom_t Int_t ULong_t ULong_t unsigned char prop_list Atom_t Atom_t target
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void value
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t g
R__EXTERN TGeoManager * gGeoManager
float * q
char * Form(const char *fmt,...)
static const Char_t * WorkingDirectory()
Definition: KVBase.h:181
Bool_t IsLabelled(const Char_t *l) const
Definition: KVBase.h:207
static const Char_t * GetDataSetEnv(const Char_t *dataset, const Char_t *type, const Char_t *defval)
Definition: KVBase.cpp:1619
static Bool_t SearchKVFile(const Char_t *name, TString &fullpath, const Char_t *kvsubdir="")
Definition: KVBase.cpp:538
Record folder for the database.
Definition: KVDBRecord.h:43
virtual Bool_t AddKey(KVDBKey *key, Bool_t check=kTRUE)
Definition: KVDBRecord.cpp:65
virtual KVRList * GetLinks(const Char_t *key) const
Returns the list of records linked to this record in table "key".
Definition: KVDBRecord.cpp:206
virtual KVDBTable * GetTable(const Char_t *table) const
Definition: KVDataBase.h:159
virtual Bool_t AddTable(KVDBTable *table)
Definition: KVDataBase.cpp:84
const Char_t * GetDataSetDir() const
Definition: KVDataSet.cpp:729
const Char_t * GetDataSetEnv(const Char_t *type, const Char_t *defval="") const
Definition: KVDataSet.cpp:767
TString GetFullPathToDataSetFile(const Char_t *filename)
Definition: KVDataSet.cpp:1898
List of hit groups in a multidetector array.
void AddGroup(KVGroup *grp)
Base class for output signal data produced by a detector.
virtual Bool_t IsExpression() const
virtual void SetValue(Double_t x)
virtual Bool_t IsRaw() const
Base class for detector geometry description.
Definition: KVDetector.h:160
KVGroup * GetGroup() const
const KVSeqCollection & GetListOfDetectorSignals() const
Definition: KVDetector.h:786
Extension of TEnv to allow the writing of comments in the file.
Definition: KVEnv.h:17
void AddCommentLine(const Char_t *line)
Definition: KVEnv.cpp:104
virtual Int_t WriteFile(const char *fname, EEnvLevel level=kEnvAll)
Definition: KVEnv.cpp:174
Base class to describe database of an experiment ,,.
Definition: KVExpDB.h:20
KVDBRun * GetDBRun(Int_t number) const
Definition: KVExpDB.h:76
virtual void LinkRecordToRunRange(KVDBRecord *rec, UInt_t first_run, UInt_t last_run)
Definition: KVExpDB.cpp:87
Bool_t FindCalibFile(const Char_t *type, TString &fullpath, const TString &array_name="") const
Definition: KVExpDB.cpp:487
Base class for FAZIA detectors.
void SetSignal(TGraph *signal, const Char_t *signal_name)
KVDetectorSignal * SetBaseLine(int sigid, Float_t baseline)
void SetGTTag(Int_t t)
void SetDetTag(Int_t t)
static const Char_t * GetNewName(KVString oldname)
KVDetectorSignal * SetFPGAEnergy(int sigid, Int_t idx, Double_t energy)
Telescope for FAZIA identification using SI1 and/or SI2 with CSI.
void AddTriggerPattern(const TString &name, uint16_t value)
Description of a FAZIA detector geometry.
Definition: KVFAZIA.h:33
Double_t fImport_Xorg
for geometry import
Definition: KVFAZIA.h:50
Double_t fImport_ThetaMin
for geometry import
Definition: KVFAZIA.h:46
Double_t fQ2risetime
Definition: KVFAZIA.h:58
KVGroupReconstructor * GetReconstructorForGroup(const KVGroup *) const
Specialized group reconstructor for FAZIA.
Definition: KVFAZIA.cpp:594
Double_t fQH1risetime
values of trapezoidal filter rise time set in the fpgas defined in .kvrootrc
Definition: KVFAZIA.h:57
Double_t TreatEnergy(Int_t sigid, Int_t eid, UInt_t val)
Definition: KVFAZIA.cpp:650
Double_t fQ3slowrisetime
Definition: KVFAZIA.h:59
void SetTriggerPatternsForDataSet(const TString &dataset)
Definition: KVFAZIA.cpp:394
void FillDetectorList(KVReconstructedNucleus *rnuc, KVHashList *DetList, const KVString &DetNames)
Definition: KVFAZIA.cpp:546
void SetGeometryImportParameters(Double_t dt=0.25, Double_t dp=1.0, Double_t tmin=2., Double_t pmin=0, Double_t tmax=20., Double_t pmax=360., Double_t xorg=0, Double_t yorg=0, Double_t zorg=0)
Definition: KVFAZIA.h:238
int fQuartet[8][2]
quartet number from #FEE and #FPGA
Definition: KVFAZIA.h:53
Double_t GetSetupParameter(const Char_t *parname)
Definition: KVFAZIA.cpp:78
Double_t fImport_dTheta
for geometry import
Definition: KVFAZIA.h:44
void SetIDCodeForIDTelescope(KVIDTelescope *) const
Set the FAZIA-specific general identification code for the given telescope.
Definition: KVFAZIA.cpp:970
std::string GetTriggerForCurrentRun() const
Definition: KVFAZIA.cpp:129
void AddDetectorLabel(const Char_t *label)
Definition: KVFAZIA.cpp:254
Bool_t handle_raw_data_event_protobuf(KVProtobufDataReader &)
Definition: KVFAZIA.cpp:641
virtual void BuildTarget()
Definition: KVFAZIA.cpp:442
Int_t fStartingBlockNumber
Definition: KVFAZIA.h:39
void DeduceIdentificationTelescopesFromGeometry()
Definition: KVFAZIA.cpp:167
TString fCorrespondanceFile
Bool_t fBuildTarget; //kTRUE to include target frame in the geometry.
Definition: KVFAZIA.h:41
Double_t fImport_PhiMax
for geometry import
Definition: KVFAZIA.h:49
virtual void SetNameOfDetectors(KVEnv &env)
Definition: KVFAZIA.cpp:304
Double_t fImport_Zorg
for geometry import
Definition: KVFAZIA.h:52
virtual void BuildFAZIA()
methods to be implemented in child classes
Definition: KVFAZIA.cpp:431
KVFAZIA(const Char_t *title="")
Default constructor.
Definition: KVFAZIA.cpp:54
void ReadTriggerPatterns(KVExpDB *db)
Definition: KVFAZIA.cpp:927
virtual void DefineStructureFormats(KVGeoImport &)
Definition: KVFAZIA.h:83
void GenerateCorrespondanceFile()
Definition: KVFAZIA.cpp:266
Double_t fQ3fastrisetime
Definition: KVFAZIA.h:60
KVString fSignalTypes
Definition: KVFAZIA.h:43
void PerformClosedROOTGeometryOperations()
Definition: KVFAZIA.cpp:335
Bool_t treat_event(const DAQ::FzEvent &)
Read raw data for an event.
Definition: KVFAZIA.cpp:681
virtual ~KVFAZIA()
Destructor.
Definition: KVFAZIA.cpp:243
Double_t fImport_dPhi
for geometry import
Definition: KVFAZIA.h:45
Double_t fImport_ThetaMax
for geometry import
Definition: KVFAZIA.h:47
int fTelescope[8][2]
telescope number from #FEE and #FPGA
Definition: KVFAZIA.h:54
KVFAZIATrigger fTrigger
trigger pattern read from data for each event
Definition: KVFAZIA.h:63
virtual void MakeCalibrationTables(KVExpDB *)
Override base method in order to read FAZIA trigger for each run.
Definition: KVFAZIA.cpp:105
Bool_t handle_raw_data_event_mfmframe(const MFMCommonFrame &)
Definition: KVFAZIA.cpp:845
virtual void GetGeometryParameters()
Called by the Build method.
Definition: KVFAZIA.cpp:366
void SetTriggerPattern(uint16_t fp)
Definition: KVFAZIA.h:69
double oldfaziats
dummy ts to control trigger info transmission rate
Definition: KVFAZIA.h:66
TString GetSignalName(Int_t bb, Int_t qq, Int_t tt, Int_t idsig)
Definition: KVFAZIA.cpp:609
Double_t fImport_PhiMin
for geometry import
Definition: KVFAZIA.h:48
Int_t fNblocks
number of blocks
Definition: KVFAZIA.h:38
void GetDetectorEvent(KVDetectorEvent *detev, const TSeqCollection *dets)
Definition: KVFAZIA.cpp:500
virtual void Build(Int_t=-1)
Build the FAZIA array.
Definition: KVFAZIA.cpp:465
Double_t fImport_Yorg
for geometry import
Definition: KVFAZIA.h:51
KVString fDetectorLabels
Definition: KVFAZIA.h:42
void CreateCorrespondence()
Definition: KVFAZIA.cpp:882
virtual void SetRawDataFromReconEvent(KVNameValueList &)
Overrides base method in order to set the value of the trigger bit pattern for the event.
Definition: KVFAZIA.cpp:92
Handle FAZIA protobuf-format raw data files.
Path taken by particles through multidetector geometry.
KVGeoDetectorNode * GetNextNode() const
KVSeqCollection * AccessIDTelescopeList()
void IterateBackFrom(const KVGeoDetectorNode *node0=nullptr) const
Information on relative positions of detectors & particle trajectories.
Import detector array described by ROOT geometry and set up corresponding KVMultiDetArray object.
Definition: KVGeoImport.h:68
void AddAcceptedDetectorName(const char *name)
void SetOrigin(double x, double y, double z)
Definition: KVGeoImport.h:101
void SetDetectorPlugin(const TString &name)
Definition: KVGeoImport.h:95
void ImportGeometry(Double_t dTheta=0.1, Double_t dPhi=1.0, Double_t ThetaMin=0.0, Double_t PhiMin=0.0, Double_t ThetaMax=180.0, Double_t PhiMax=360.0)
Definition: KVGeoImport.cpp:99
void SetNameCorrespondanceList(const Char_t *)
virtual Bool_t Fired(Option_t *opt="any") const
KVDetector * GetDetector(const Char_t *name) const
Return detector in this structure with given name.
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:20
Extended version of ROOT THashList.
Definition: KVHashList.h:29
Base class for all detectors or associations of detectors in array which can identify charged particl...
Definition: KVIDTelescope.h:84
void SetGroup(KVGroup *kvg)
KVDetector * GetDetector(UInt_t n) const
virtual void SetIDCode(UShort_t c)
KVGroup * GetGroup() const
virtual void AddDetector(KVDetector *d)
KVString GetDetectorLabelsForGridCoord(const KVString &axis) const
Description of physical materials used to construct detectors & targets; interface to range tables.
Definition: KVMaterial.h:94
virtual Double_t GetThickness() const
Definition: KVMaterial.cpp:487
virtual TGeoMedium * GetGeoMedium(const Char_t *="")
static KVIonRangeTable * GetRangeTable()
Definition: KVMaterial.cpp:166
Base class for describing the geometry of a detector array.
virtual void GetDetectorEvent(KVDetectorEvent *detev, const TSeqCollection *fired_params=0)
static Bool_t fCloseGeometryNow
virtual KVGroup * GetGroup(const Char_t *) const
Return pointer to group with name.
void CreateGeoManager(Double_t dx=500, Double_t dy=500, Double_t dz=500)
Bool_t fHandledRawData
set to true if multidetector handles data in last call to HandleRawData
KVTarget * GetTarget()
KVSeqCollection * fIDTelescopes
deltaE-E telescopes in groups
TString GetDataSet() const
virtual void MakeCalibrationTables(KVExpDB *)
UInt_t GetCurrentRunNumber() const
virtual void SetIdentifications()
const TSeqCollection * GetTrajectories() const
static Bool_t fBuildTarget
virtual void SetRawDataFromReconEvent(KVNameValueList &)
KVUniqueNameList fFiredDetectors
list of fired detectors after reading raw data event
KVNameValueList fReconParameters
general purpose list of parameters for storing information on data reconstruction
TString fDataSet
name of associated dataset, used with MakeMultiDetector()
KVUnownedList fFiredSignals
list of fired signals after reading raw data event
Handles lists of named parameters with different types, a list of KVNamedParameter objects.
Int_t GetIntValue(const Char_t *name) const
Double_t GetDoubleValue(const Char_t *name) const
void SetValue(const Char_t *name, value_type value)
void SetValue64bit(const Char_t *name, ULong64_t)
ULong64_t GetValue64bit(const Char_t *name) const
Bool_t HasParameter(const Char_t *name) const
Strings used to represent a set of ranges of values.
Definition: KVNumberList.h:85
KVNameValueList * GetParameters() const
Definition: KVParticle.h:815
Read Google Protobuf DAQ files.
Nuclei reconstructed from data measured by a detector array .
virtual void Clear(Option_t *option="")
virtual void Add(TObject *obj)
virtual void Delete(Option_t *option="")
Extension of ROOT TString class which allows backwards compatibility with ROOT v3....
Definition: KVString.h:73
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
Calculation/correction of energy losses of particles through an experimental target.
Definition: KVTarget.h:127
virtual void Add(TObject *obj)
virtual Int_t GetEntries() const
THashList * GetTable() const
virtual const char * GetValue(const char *name, const char *dflt) const
virtual Int_t ReadFile(const char *fname, EEnvLevel level)
virtual void SetValue(const char *name, const char *value, EEnvLevel level=kEnvChange, const char *type=nullptr)
void CloseGeometry(Option_t *option="d")
TGeoVolume * GetTopVolume() const
TGeoVolume * MakeEltu(const char *name, TGeoMedium *medium, Double_t a, Double_t b, Double_t dz)
virtual TGeoNode * AddNode(TGeoVolume *vol, Int_t copy_no, TGeoMatrix *mat=nullptr, Option_t *option="")
virtual void SetPoint(Int_t i, Double_t x, Double_t y)
TObject * Next()
const char * GetName() const override
static TClass * Class()
virtual void SetName(const char *name)
TClass * IsA() const override
void AbstractMethod(const char *method) const
virtual const char * GetName() const
void SetBit(UInt_t f)
virtual const char * ClassName() const
virtual void Warning(const char *method, const char *msgfmt,...) const
virtual Bool_t InheritsFrom(const char *classname) const
virtual void Error(const char *method, const char *msgfmt,...) const
virtual void Info(const char *method, const char *msgfmt,...) const
TObject * First() const override
Ssiz_t Length() const
const char * Data() const
Bool_t IsWhitespace() const
void Form(const char *fmt,...)
Bool_t Contains(const char *pat, ECaseCompare cmp=kExact) const
TGraphErrors * gr
TH1 * h
double T(double x)
rec
constexpr Double_t R()
TLine l
auto * tt
ClassImp(TPyArg)