KaliVeda
Toolkit for HIC analysis
KVMultiDetArray.cpp
1 //Created by KVClassFactory on Tue Apr 16 09:45:50 2013
2 //Author: John Frankland,,,
3 
4 #include "KVMultiDetArray.h"
5 #include "KVDetector.h"
6 #include "KVDetectorEvent.h"
7 #include "KVReconstructedEvent.h"
8 #include "KVReconstructedNucleus.h"
9 #include "KVRList.h"
10 #include "KVNucleus.h"
11 #include "KVGroup.h"
12 #include "KVMaterial.h"
13 #include "KVTarget.h"
14 #include "KVIDTelescope.h"
15 #include <KVString.h>
16 #include <TObjString.h>
17 #include <TObjArray.h>
18 #include <KVIDGridManager.h>
19 #include <KVDataSetManager.h>
20 #include <KVUpDater.h>
21 #include "TPluginManager.h"
22 #include "KVDataSet.h"
23 #include "TGeoManager.h"
24 #include "KVHashList.h"
25 #include "KVNameValueList.h"
26 #include "KVUniqueNameList.h"
27 #include "KVIonRangeTable.h"
28 #include "KVRangeTableGeoNavigator.h"
29 #include <KVDataAnalyser.h>
30 #include <KVNamedParameter.h>
31 #include <KVCalibrator.h>
32 #include <KVDBParameterSet.h>
33 #ifdef WITH_RSQLITE
34 #include <KVSQLROOTIDGridManager.h>
35 #endif
36 #ifdef WITH_OPENGL
37 #include <TGLViewer.h>
38 #include <TVirtualPad.h>
39 #endif
40 #ifdef WITH_BUILTIN_GRU
41 #include "KVGANILDataReader.h"
42 #else
43 #include "KVRawDataReader.h"
44 #endif
45 #ifdef WITH_MFM
46 #include "KVMFMDataFileReader.h"
47 #include "MFMEbyedatFrame.h"
48 #endif
49 #ifdef WITH_PROTOBUF
50 #include "KVProtobufDataReader.h"
51 #endif
52 using namespace std;
53 
54 KVMultiDetArray* gMultiDetArray = nullptr;
55 
61 
63 
64 
65 
69  : KVGeoStrucElement(), fTrajectories(kTRUE)
70 {
71  // Default constructor
72  init();
73  gMultiDetArray = this;
74 }
75 
76 
77 
80 
82  : KVGeoStrucElement(name, type), fTrajectories(kTRUE)
83 {
84  // Constructor with name and optional type
85  init();
86  gMultiDetArray = this;
87 }
88 
89 
90 
102 
104 {
105  //Basic initialisation called by constructor.
106  //Creates detectors list fDetectors,
107  //groups list fGroups, identification telescopes list
108  //fIDTelescopes
109  //
110  //The fGroups & fIDTelescopes lists contain objects owned by the multidetector array,
111  //but which may be deleted by other objects (or as a result of the deletion of other
112  //objects: i.e. if all the detectors in a group are deleted, the group itself is destroyed).
113  //We use the ROOT automatic garbage collection to make sure that any object deleted
114  //elsewhere is removed automatically from these lists.
115 
116  fIDTelescopes = new KVHashList();
117  fIDTelescopes->SetOwner(kTRUE); // owns its objects
119 
120  fHitGroups = 0;
121 
122  fTarget = 0;
123  fCurrentRun = 0;
124 
126  fCalibStatusDets = 0;
127  fSimMode = kFALSE;
128 
130 
131  fUpDater = 0;
132 
133  if (!gIDGridManager) new KVIDGridManager;
134 
135  // all trajectories belong to us
137 
138  //all detectors belong to us
140 
141  fRawDataReader = nullptr;
142  fHandledRawData = false;
143 
144  // any extra raw data signals created when reading data belong to us
146 }
147 
148 
149 
152 
154 {
155  //destroy (delete) the MDA and all the associated structure, detectors etc.
156 
158  //destroy all identification telescopes
161  delete fIDTelescopes;
162  }
163  fIDTelescopes = 0;
164 
165  if (gMultiDetArray == this)
166  gMultiDetArray = nullptr;
167 
168  if (fStatusIDTelescopes) {
170  delete fStatusIDTelescopes;
172  }
173  if (fCalibStatusDets) {
175  delete fCalibStatusDets;
176  fCalibStatusDets = 0;
177  }
178 
180 }
181 
182 
183 
184 
185 
187 
189 {
190 
191 
192 
193 }
194 
195 
196 
197 
221 
223 {
224  // Create one or more KVIDTelescope particle-identification objects from the two detectors
225  //
226  // The different ID telescopes are defined as 'Plugin' objects in the file $KVROOT/KVFiles/.kvrootrc :
227  // # The KVMultiDetArray::GetIDTelescopes(KVDetector*de, KVDetector*e) method uses these plugins to
228  // # create KVIDTelescope instances adapted to the specific array geometry and detector types.
229  // # For each pair of detectors we look for a plugin with one of the following names:
230  // # [name_of_dataset].array_name.de_detector_type[de detector thickness]-e_detector_type[de detector thickness]
231  // # Each characteristic in [] brackets may or may not be present in the name; first we test for names
232  // # with these characteristics, then all combinations where one or other of the characteristics is not present.
233  // # In addition, we first test all combinations which begin with [name_of_dataset].
234  // # The first plugin found in this order will be used.
235  // # In addition, if for one of the two detectors there is a plugin called
236  // # [name_of_dataset].array_name.de_detector_type[de detector thickness]
237  // # [name_of_dataset].array_name.e_detector_type[e detector thickness]
238  // # then we add also an instance of this 1-detector identification telescope.
239  //
240  // This method is called by DeduceIdentificationTelescopesFromGeometry
241  // in order to set up all ID telescopes of the array.
242  //
243  // \returns number of ID telescopes created
244  // \param[in] de node corresponding to \f$\Delta E\f$ detector; may be nullptr
245  // \param[in] e residual energy detector; always points to valid object
246 
247  assert(e);
248 
249  Int_t ntels = 0;
250 
251  if (fDataSet == "" && gDataSet) fDataSet = gDataSet->GetName();
252 
253  //look for ID telescopes starting from furthest from target
254  if (e->IsOK()) ntels += try_all_singleID_telescopes(e, list);
255  if (de && (de->GetDetector() != e)) {
256  if (e->IsOK() && de->GetDetector()->IsOK()) ntels += try_all_doubleID_telescopes(de->GetDetector(), e, list);
257  if (de->GetDetector()->IsOK()) ntels += try_all_singleID_telescopes(de->GetDetector(), list);
258  }
259 
260  return ntels;
261 }
262 
263 
264 
282 
284 {
285  // Attempt to find a plugin KVIDTelescope class for making a single-detector
286  // ID telescope from detector *d
287  // We look for plugins with the following signatures (uri):
288  //
289  // [array name].[type]
290  // [array_name].[type][thickness]
291  //
292  // where 'type' is the type of the detector in UPPER or lowercase letters
293  // 'thickness' is the nearest-integer thickness of the detector as returned by d->GetThickness()
294  // In addition, if a dataset is set (gDataSet!=nullptr) we try also for dataset-specific
295  // plugins:
296  //
297  // [dataset].[array name].[type]
298  // [dataset].[array name].[type][thickness]
299  //
300  // Returns number of generated telescopes
301 
302  TString uri = Form("%s.%s", GetName(), d->GetType());
303  Int_t ntels = 0;
304  if (!(ntels += try_upper_and_lower_singleIDtelescope(uri, d, l))) {
305  Int_t d_thick = TMath::Nint(d->GetThickness());
306  uri += d_thick;
307  ntels += try_upper_and_lower_singleIDtelescope(uri, d, l);
308  }
309  return ntels;
310 }
311 
312 
313 
314 
337 
339 {
340  // Attempt to find a plugin KVIDTelescope class for making an ID telescope from detectors de & e.
341  // We look for plugins with the following signatures (uri):
342  //
343  // [array name].[de-type]-[e-type]
344  // [array name].[de-type][thickness]-[e-type]
345  // [array name].[de-type]-[e-type][thickness]
346  // [array name].[de-type][thickness]-[e-type][thickness]
347  //
348  // where 'type' is the type of the detector in UPPER or lowercase letters
349  // 'thickness' is the nearest-integer thickness of the detector.
350  // In addition, if a dataset is set (gDataSet!=nullptr) we try also for dataset-specific
351  // plugins:
352  //
353  // [dataset].[array name].[de-type][thickness]-[e-type][thickness]
354  // [dataset].[array name].[de-type][thickness]-[e-type]
355  // [dataset].[array name].[de-type]-[e-type][thickness]
356  // [dataset].[array name].[de-type]-[e-type]
357  //
358  // if no plugin is found, we return a KVIDTelescope base class object
359  //
360  // Returns 1 (we always generate exactly one telescope)
361 
362  TString de_type = de->GetType();
363  TString e_type = e->GetType();
364  TString de_thick = Form("%d", TMath::Nint(de->GetThickness()));
365  TString e_thick = Form("%d", TMath::Nint(e->GetThickness()));
366 
367  TString uri = de_type + de_thick + "-" + e_type + e_thick;
368  uri.Prepend(Form("%s.", GetName()));
369  if (try_upper_and_lower_doubleIDtelescope(uri, de, e, l)) return 1;
370 
371  uri = de_type + de_thick + "-" + e_type;
372  uri.Prepend(Form("%s.", GetName()));
373  if (try_upper_and_lower_doubleIDtelescope(uri, de, e, l)) return 1;
374 
375  uri = de_type + "-" + e_type + e_thick;
376  uri.Prepend(Form("%s.", GetName()));
377  if (try_upper_and_lower_doubleIDtelescope(uri, de, e, l)) return 1;
378 
379  uri = de_type + "-" + e_type;
380  uri.Prepend(Form("%s.", GetName()));
381  if (try_upper_and_lower_doubleIDtelescope(uri, de, e, l)) return 1;
382 
383  // default id telescope object
384  KVIDTelescope* idt = new KVIDTelescope;
385  uri = de_type + "-" + e_type;
386  idt->SetLabel(uri);
387  // set type as "[DEtype]-[Etype]" where "[*type]" is type of DE/E detectors
388  idt->SetType(Form("%s-%s", de_type.Data(), e_type.Data()));
389  set_up_telescope(de, e, idt, l);
390 
391  return 1;
392 }
393 
394 
395 
403 
405 {
406  // Attempt to find a plugin KVIDTelescope class for making a single-detector
407  // ID telescope from detector *d with the given signature/uri
408  // Both original & all-upper-case versions of uri are tried.
409  // uri is tried both with & without prepended dataset name (if set)
410  // Returns true if successful (the new ID telescope will be added to internal
411  // list fIDTelescopes and also to TCollection* l)
412 
413  if (try_a_singleIDtelescope(uri, d, l)) return true;
414  uri.ToUpper();
415  return try_a_singleIDtelescope(uri, d, l);
416 }
417 
418 
419 
426 
428 {
429  // Attempt to find a plugin KVIDTelescope class for making an ID telescope with the given signature/uri
430  // Both original & all-upper-case versions of uri are tried.
431  // uri is tried both with & without prepended dataset name (if set)
432  // Returns true if successful (the new ID telescope will be added to internal
433  // list fIDTelescopes and also to TCollection* l)
434 
435  if (try_a_doubleIDtelescope(uri, de, e, l)) return true;
436  uri.ToUpper();
437  return try_a_doubleIDtelescope(uri, de, e, l);
438 }
439 
440 
441 
450 
452 {
453  // Attempt to find a plugin KVIDTelescope class for making a single-detector
454  // ID telescope from detector *d with the given signature/uri
455  // Both original & all-upper-case versions of uri are tried.
456  // uri is tried both with & without prepended dataset name (if set)
457  // Returns true if successful (the new ID telescope will be added to internal
458  // list fIDTelescopes and also to TCollection* l)
459 
460  // dataset-specific version takes precedence over default
461  TString duri = uri;
462  if (gDataSet) {
463  // try with dataset name
464  duri.Prepend(Form("%s.", fDataSet.Data()));
465  KVIDTelescope* idt;
466  if ((idt = KVIDTelescope::MakeIDTelescope(duri))) {
468  return true;
469  }
470  }
471 
472  // look for default version
473  KVIDTelescope* idt;
474  if ((idt = KVIDTelescope::MakeIDTelescope(uri))) {
476  return true;
477  }
478 
479  return false;
480 }
481 
482 
483 
490 
492 {
493  // Attempt to find a plugin KVIDTelescope class for making an ID telescope with the given signature/uri
494  // uri is tried both with & without prepended dataset name (if set)
495  // Returns true if successful (the new ID telescope will be added to internal
496  // list fIDTelescopes and also to TCollection* l)
497 
498  // dataset-specific version takes precedence over default
499  TString duri = uri;
500  if (gDataSet) {
501  // try with dataset name
502  duri.Prepend(Form("%s.", fDataSet.Data()));
503  KVIDTelescope* idt;
504  if ((idt = KVIDTelescope::MakeIDTelescope(duri))) {
505  set_up_telescope(de, e, idt, l);
506  return true;
507  }
508  }
509  // look for default version
510  KVIDTelescope* idt;
511  if ((idt = KVIDTelescope::MakeIDTelescope(uri))) {
512  set_up_telescope(de, e, idt, l);
513  return true;
514  }
515 
516  return false;
517 }
518 
519 
520 
523 
525 {
526  // Set up detectors in de-e identification telescope and add to fIDTelescopes and to l
527 
528  idt->AddDetector(de);
529  idt->AddDetector(e);
530  if (de->GetGroup()) {
531  idt->SetGroup(de->GetGroup());
532  }
533  else {
534  idt->SetGroup(e->GetGroup());
535  }
536  // if telescope already exists, we delete this new version and add a reference to
537  // the original into list l
539  if (p) {
540  l->Add(p);
541  delete idt;
542  }
543  else {
544  fIDTelescopes->Add(idt);
545  l->Add(idt);
546  }
547 }
548 
549 
550 
553 
555 {
556  // Set up detector in single-stage identification telescope and add to fIDTelescopes and to l
557 
558  idt->AddDetector(det);
559  idt->SetGroup(det->GetGroup());
560  // if telescope already exists, we delete this new version and add a reference to
561  // the original into list l
563  if (p) {
564  l->Add(p);
565  delete idt;
566  }
567  else {
568  fIDTelescopes->Add(idt);
569  l->Add(idt);
570  }
571 }
572 
573 
574 
575 
579 
581 {
582  //Number groups according to position in list fGroups and set fGr counter to the number
583  //of groups in the list
584  Int_t fGr = 0;
585  KVGroup* g = 0;
586  KVSeqCollection* fGroups = GetStructures()->GetSubListWithType("GROUP");
587  TIter next(fGroups);
588  while ((g = (KVGroup*) next())) {
589  g->SetNumber(++fGr);
590  }
591  delete fGroups;
592 }
593 
594 
595 
598 
600 {
601  //Return pointer to DeltaE-E ID Telescope with "name"
602 
604 }
605 
606 
607 
612 
614 {
615  //Reset all groups (lists of detected particles etc.)
616  //and detectors in groups (energy losses, ACQparams etc. etc.)
617  //and the target if there is one
618 
619  unique_ptr<KVSeqCollection> fGroups(GetStructures()->GetSubListWithType("GROUP"));
620 
621  TIter next(fGroups.get());
622  KVGroup* grp;
623  while ((grp = (KVGroup*) next())) {
624  grp->Reset();
625  }
626  if (GetTarget())
627  GetTarget()->Clear();
628 }
629 
630 
631 
632 
638 
640 {
641  // Sets up calibrators for all detectors with a defined calibration for run
642  // Set parameters for all detectors with links to table "Calibrations" for run
643  // If 'myname' is given, we look in "myname.Calibrations"
644 
645  //Reset all calibrators of all detectors first
646  TIter next(GetDetectors());
647  KVDetector* kvd;
648  while ((kvd = (KVDetector*) next())) kvd->RemoveCalibrators();
649 
650  TString tabname = (myname != "" ? Form("%s.Calibrations", myname.Data()) : "Calibrations");
651  //Info("SetCalibratorParameters", "For array %s in table %s", GetName(), tabname.Data());
652  KVRList* run_links = r->GetLinks(tabname);
653  if (run_links) Info("SetCalibratorParameters", "Found %d calibrations for this run", run_links->GetEntries());
654  else {
655 // Warning("SetCalibratorParameters", "Got no links for %s", tabname.Data());
656 // r->GetKeys()->ls();
657  return;
658  }
659  TIter nxt_link(run_links);
660  KVDBParameterSet* dbps;
661  while ((dbps = (KVDBParameterSet*)nxt_link())) {
662 
663  KVDetector* det = GetDetector(dbps->GetName());
664  if (!det) {
665  Warning("SetCalibratorParameters", "Got parameters for unknown detector: %s", dbps->GetName());
666  continue;
667  }
668 
669  KVNameValueList class_options;
670  KVString clop;
671  if (dbps->HasParameter("CalibOptions")) clop = dbps->GetStringParameter("CalibOptions");
672  if (clop != "") {
673  clop.Begin(",");
674  while (!clop.End()) {
675  KVString clopp = clop.Next(true);
676  clopp.Begin("=");
677  KVString par(clopp.Next(true)), val(clopp.Next(true));
678  class_options.SetValue(par, val);
679  }
680  }
682  cal->SetType(dbps->GetTitle());
683  if (clop != "") {
684  try {
685  cal->SetOptions(class_options);
686  }
687  catch (std::exception& e) {
688  Error("SetCalibratorParameters",
689  "Problem for %s [%s] : %s", det->GetName(), cal->GetType(), e.what());
690  delete cal;
691  continue;
692  }
693  }
694  cal->SetInputSignalType(dbps->GetStringParameter("SignalIn"));
695  cal->SetOutputSignalType(dbps->GetStringParameter("SignalOut"));
696  if (!det->AddCalibrator(cal, dbps->GetParameters())) {
697  // Calibrator invalid - probably input signal is not defined for detector
698  // N.B. 'cal' deleted by KVDetector::AddCalibrator
699  continue;
700  }
701 
702  if (dbps->GetParamNumber() > cal->GetNumberParams()) {
703  Warning("SetCalibratorParameters", "Wrong number of parameters (%d) for calibrator %s for detector %s : should be %d",
704  dbps->GetParamNumber(), dbps->GetTitle(), dbps->GetName(), cal->GetNumberParams());
705  dbps->Print();
706  continue;
707  }
708  for (int i = 0; i < dbps->GetParamNumber(); ++i) {
709  if (i >= cal->GetNumberParams())
710  cal->SetParameter(i, 0);
711  else
712  cal->SetParameter(i, dbps->GetParameter(i));
713  }
714  cal->SetStatus(true);
715  }
716 }
717 
718 
719 
737 
739 {
740  // Detector signals corresponding to raw data acquisition parameters are typically only created and
741  // added to detectors when some raw data has been read including those parameters.
742  //
743  // However, when reconstructing data, we may define identification matrices or calibration formulae
744  // which use these signals before starting to read data. Therefore we need to know beforehand what
745  // detector signals are expected to be available once data has been read.
746  //
747  // These are defined, according to detector types, by variables of the form
748  //
749  //~~~~
750  // [dataset].[array].[detector-type].ExpectedDetectorSignals: [comma-separated list of signal names]
751  //~~~~
752  //
753  // where [dataset] is optionally used to provide dataset-specific definitions.
754  //
755  // Here we add a detector signal of each expected type to each detector of the array
756 
757  TIter it(GetDetectors());
758  KVDetector* det;
759  while ((det = (KVDetector*)it())) {
760  auto s = KVBase::GetDataSetEnv<KVString>(fDataSet.Data(), Form("%s.%s.ExpectedDetectorSignals", GetName(), det->GetType()));
761  if (s.IsNull()) continue;
762  s.Begin(",");
763  while (!s.End()) {
764  det->AddDetectorSignal(s.Next(kTRUE));
765  }
766  }
767 }
768 
769 
770 
781 
783 {
784  // First step in event reconstruction based on current status of detectors in array.
785  // Fills the given KVDetectorEvent with the list of all groups which have fired.
786  // i.e. loop over all groups of the array and test whether KVGroup::Fired() returns true or false.
787  //
788  // If the list of fired detectors 'fired_dets' is given, then we use this list
789  // to find the associated groups. If not given, or if it is empty, we may use the internal fFiredDetectors list.
790  //
791  // Call method detev->Clear() before reading another event in order to reset all of the hit groups
792  // (including all detectors etc.) and emptying the list.
793 
794  if (!fired_dets || !fired_dets->GetEntries()) {
795  if (fFiredDetectors.GetEntries()) fired_dets = &fFiredDetectors;
796  }
797  if (fired_dets && fired_dets->GetEntries()) {
798  // list of fired detectorsgiven
799  TIter next_det(fired_dets);
800  KVDetector* det = 0;
801  KVGroup* grp = 0;
802  while ((det = (KVDetector*)next_det())) {
803  if ((grp = det->GetGroup()) && grp->GetParents()->Contains(this)) detev->AddGroup(grp);
804  }
805  }
806  else {
807  //loop over groups
808  unique_ptr<KVSeqCollection> fGroups(GetStructures()->GetSubListWithType("GROUP"));
809 
810  TIter next_grp(fGroups.get());
811  KVGroup* grp;
812  while ((grp = (KVGroup*) next_grp())) {
813  if (grp->Fired()) {
814  //if (!fHitGroups->FindObject(grp))
815  // grp->Print();
816  //add new group to list of hit groups
817  detev->AddGroup(grp);
818  }
819  }
820  }
821 }
822 
823 
824 
841 
843 {
844  //Static function which will create and 'Build' the multidetector object corresponding to
845  //a given run of dataset 'dataset_name'. Any previously existing multidetector will be
846  //deleted.
847  //We first activate the given dataset if not already done
848  //
849  //Multidetector arrays are defined as 'Plugin' objects in the file $KVROOT/KVFiles/.kvrootrc :
850  //
851  //Plugin.KVMultiDet: [dataset_name] [classname] [library] "[constructor]()"
852  //
853  //The constructors/macros are always without arguments
854  //
855  //Dataset name is stored in fDataSet
856 
857  // store the run number (if given) so that if the dataset needs to update its database
858  // (which requires building the multidetector by calling MakeMultiDetector()) then it
859  // will build the correct geometry for the run
860  if (run > 0) fMakeMultiDetectorRunNumber = run;
862 
863  if (gDataSetManager && (!gDataSet || (gDataSet != gDataSetManager->GetDataSet(dataset_name)))) {
864  printf("Info in <KVMultiDetArray::MakeMultiDetector>: Changing dataset\n");
865  gDataSetManager->GetDataSet(dataset_name)->cd();
866  }
867 
868  if (gMultiDetArray && gMultiDetArray->GetDataSet() != dataset_name) {
869  printf("Info in <KVMultiDetArray::MakeMultiDetector>: Deleting existing array %s\n", gMultiDetArray->GetName());
870  if (gIDGridManager) {
871  delete gIDGridManager;
872  gIDGridManager = nullptr;
873  }
874  delete gMultiDetArray;
875  gMultiDetArray = nullptr;
876  }
877 
878 
879  // Creation of database when dataset is selected for first time may
880  // include creation of multidetector array (by calling this method)
881  KVMultiDetArray* mda = nullptr;
882  if (!gMultiDetArray) {
883  TPluginHandler* ph;
884  if (!(ph = LoadPlugin(classname.Data(), dataset_name)))
885  return nullptr;
886 
887  //execute constructor/macro for multidetector - assumed without arguments
888  mda = (KVMultiDetArray*) ph->ExecPlugin(0);
889  mda->fDataSet = dataset_name;
890  mda->Build(run);
891  // set dataset-dependent lists of acceptable ID/E codes for reconstructed nuclei
892  auto codes = GetDataSetEnv<KVString>(dataset_name, Form("%s.ReconstructedNuclei.AcceptIDCodes", mda->GetName()));
893  if (codes != "") mda->fAcceptIDCodes.Set(codes);
894  codes = GetDataSetEnv<KVString>(dataset_name, Form("%s.ReconstructedNuclei.AcceptECodes", mda->GetName()));
895  if (codes != "") mda->fAcceptECodes.Set(codes);
896  // set dataset-dependent condition for seeding reconstructed nuclei
897  mda->SetPartSeedCond(GetDataSetEnv<KVString>(dataset_name, Form("%s.ReconstructedNuclei.ParticleSeedCond", mda->GetName())));
898 #ifdef WITH_RSQLITE
899  // save contents of grid manager in an SQL-ROOT database if not already done
900  if (fMakeMultiDetectorSetParameters && !gIDGridManager->IsSQLROOT()) {
901  // 'if(fMakeMultiDetectorSetParameters...' ensures that we are not currently building a subarray
902  // of a KVExpSetUp. Grids can only be saved once in a single file after all have been read.
903  TString filepath;
904  if (is_gnuinstall()) {
905  // GNU-style install: use working directory $HOME/.kaliveda
906  filepath = GetWORKDIRFilePath(gDataSet->GetName());
907  }
908  else
909  filepath = gDataSet->GetDataSetDir();
910 
911  filepath += "/idgrids_DB";
912 
913  int n_grids_to_write = gIDGridManager->GetGrids()->GetEntries();
914 
915  // we retrieve the index multiplier for the dataset in case it is >1,
916  // we use it to "correct" the runlists for the grids
917  auto index_multiplier = gDataSet->GetDataSetEnv("DataSet.RunFileIndexMultiplier.raw", 1);
918 
919  // we also set the runlist for any grids which don't have one, to the runlist of the entire dataset.
920  // in other words, grids with no runlist are considered valid for all runs
921  if (n_grids_to_write) {
922  KVSQLROOTFile f(filepath, "recreate");
923 
924  printf("Info in <KVMultiDetArray::MakeMultiDetector>: Saving %d grids in SQL-ROOT database file %s\n",
925  n_grids_to_write, filepath.Data());
926 
927  TIter it(gIDGridManager->GetGrids());
928  KVIDGraph* gr;
929  if (!db) db = gExpDB;
930  while ((gr = (KVIDGraph*)it())) {
931  // self-consistently modify runlist if necessary
932  auto rl = gr->GetRuns().GetListDividedBy(index_multiplier);
933  if (rl.IsEmpty() && db) rl = db->GetRunList();
934  gr->SetRuns(rl);
935  f.WriteObject(gr, {
936  {"IDLabel", gr->GetIDTelescopeLabel()},
937  {"IDTelescopes", gr->GetParameters()->GetStringValue("IDTelescopes")},
938  {"VarX", gr->GetVarX()},
939  {"VarY", gr->GetVarY()},
940  {"Runlist", gr->GetRunList()},
941  }
942  );
943  --n_grids_to_write;
944  if (!(n_grids_to_write % 1000)) printf("Info in <KVMultiDetArray::MakeMultiDetector>: ...%d grids left...\n",
945  n_grids_to_write);
946 
947  }
948  }
949  }
950 #endif
951  }
952  else {
953  mda = gMultiDetArray;
954  }
955  // set parameters if required & allowed & not done yet
957  return mda;
958 }
959 
960 
961 
973 
975 {
976  // Return pointer to KVUpDater defined by dataset for this multidetector, the class used
977  // is defined as a plugin like this:
978  //
979  // # Plugin.KVUpDater: name_of_dataset name_of_class name_of_plugin_library constructor_to_call
980  //
981  // However, if a dataset defines a variable like this:
982  //
983  // [dataset].ExpSetUp.Updater.[multidetector name]: [name_of_dataset for plugin]
984  //
985  // then we use the updater plugin defined for the given dataset
986 
987  if (!fUpDater) {
988  auto alt_updater = KVBase::GetDataSetEnv<KVString>(fDataSet, Form("ExpSetUp.Updater.%s", GetName()));
989  if (alt_updater != "") fUpDater = KVUpDater::MakeUpDater(alt_updater, this);
991  }
992  Info("GetUpDater", "updater class for dataset %s: %s", fDataSet.Data(), fUpDater->IsA()->GetName());
993  return fUpDater;
994 }
995 
996 
997 
998 
1007 
1008 void KVMultiDetArray::SetParameters(UInt_t run, Bool_t physics_parameters_only)
1009 {
1010  // Set run-dependent parameters of the array.
1011  //
1012  // if physics_parameters_only==false, identification and calibration parameters are set.
1013  // if physics_parameters_only==true, just the minimum necessary for physics analysis of reduced data are set.
1014  //
1015  // This can only be done if gDataSet has been set i.e. a dataset has been chosen,
1016  // otherwise this just has the effect of setting the current run number
1017 
1018  fCurrentRun = run;
1019  KVDataSet* ds = gDataSet;
1020  if (!ds) {
1021  if (gDataSetManager)
1022  ds = gDataSetManager->GetDataSet(fDataSet.Data());
1023  }
1024  if (ds) {
1025  GetUpDater()->SetParameters(run, physics_parameters_only);
1026  SetBit(kParamsSet);
1027  }
1028 }
1029 
1030 
1031 
1032 
1037 
1039 {
1040  //Set identification parameters for run.
1041  //This can only be done if gDataSet has been set i.e. a dataset has been chosen
1042  //Otherwise this just has the effect of setting the current run number
1043 
1044  fCurrentRun = run;
1045  KVDataSet* ds = gDataSet;
1046  if (!ds) {
1047  if (gDataSetManager)
1048  ds = gDataSetManager->GetDataSet(fDataSet.Data());
1049  }
1050  if (ds) {
1053  }
1054 }
1055 
1056 
1057 
1058 
1063 
1065 {
1066  //Set calibration parameters for run.
1067  //This can only be done if gDataSet has been set i.e. a dataset has been chosen
1068  //Otherwise this just has the effect of setting the current run number
1069 
1070  fCurrentRun = run;
1071  KVDataSet* ds = gDataSet;
1072  if (!ds) {
1073  if (gDataSetManager)
1074  ds = gDataSetManager->GetDataSet(fDataSet.Data());
1075  }
1076  if (ds) {
1079  }
1080 }
1081 
1082 
1083 
1084 
1101 
1103 {
1104  //Initialisation of all ACTIVE identification telescopes in the array, i.e. those appearing in a line
1105  //in the .kvrootrc file such as this:
1106  //
1107  //# [dataset name].ActiveIdentifications: [type1] [type2] ...
1108  //
1109  //The 'types' given correspond to the value given by KVIDTelescope::GetLabel(), these are the
1110  //identifiers used to retrieve the different plugin classes in GetIDTelescopes(KVDetector*,KVDetector*,KVList*).
1111  //
1112  //For each type of identification in the list, we retrieve the first identification telescope with this
1113  //label from the list of all ID telescopes, in order to call its KVIDTelescope::SetIdentificationParameters() method.
1114  //This method (when rederived in child classes of KVIDTelescope) initialises the identification objects
1115  //for ALL of the ID telescopes of the same type (class) in the array.
1116  //
1117  //Note that, in general, the parameters of the identifications for a given run are not
1118  //set until SetParameters or SetRunIdentificationParameters is called.
1119 
1120 #ifdef WITH_RSQLITE
1121  // if a KVSQLROOTFile has been filled with the dataset identification parameters we use it.
1122  if (gIDGridManager->IsSQLROOT()) {
1123  // nothing more to do
1124  return;
1125  }
1126  else if (!gDataSet->DataBaseUpdateInProgress())
1127  // the update of the database may have been caused by 1 or more identification files being modified
1128  // (even if they are not part of the database). in this case we read in all the grids again even
1129  // if it has already been done before.
1130  {
1131  TString filepath;
1132  if (is_gnuinstall()) {
1133  // GNU-style install: use working directory $HOME/.kaliveda
1134  filepath = GetWORKDIRFilePath(gDataSet->GetName());
1135  }
1136  else
1137  filepath = gDataSet->GetDataSetDir();
1138 
1139  filepath += "/idgrids_DB";
1140 
1141  if (SearchKVFile(filepath.Data(), filepath)) {
1142  delete gIDGridManager;
1143  gIDGridManager = new KVSQLROOTIDGridManager(filepath);
1144  return;
1145  }
1146  }
1147 #endif
1148 
1149  auto id_labels = GetDataSetEnv<KVString>(fDataSet, "ActiveIdentifications");
1150  if (id_labels == "" || (gDataSet && !gDataSet->HasCalibIdentInfos())) {
1151  Info("SetIdentifications", "No active identifications");
1152  return;
1153  }
1154  //split list of labels
1155  id_labels.Begin(" ");
1156  int ok(0);
1157  //loop over labels/identification 'types'
1158  while (!id_labels.End()) {
1159 
1160  //get first telescope in list with right label
1162  //set ID parameters for all telescopes of this 'type'
1163  if (idt) {
1164  Info("SetIdentifications", "Initialising %s identifications...", idt->GetLabel());
1165  if (idt->SetIdentificationParameters(this))
1166  Info("SetIdentifications", "OK");
1167  ++ok;
1168  }
1169 
1170  }
1171  if (!ok) {
1172  // None of the labels in the list correspond to telescopes in the array
1173  Warning("SetIdentfications", "No telescopes found with labels given in %s.ActiveIdentifications list: %s",
1174  gDataSet->GetName(), id_labels.Data());
1175  }
1176 }
1177 
1178 
1179 
1180 
1187 
1189 {
1190  // Calls Initialize() method of each identification telescope (see KVIDTelescope
1191  // and derived classes) and sets the general identification code defined for each
1192  // telescope.
1193  //
1194  // Calling this method is essential before identification of particles is attempted.
1195 
1196  TIter next(fIDTelescopes);
1197  KVIDTelescope* idt;
1198  while ((idt = (KVIDTelescope*)next())) {
1199  idt->Initialize();
1201  }
1202 }
1203 
1204 
1205 
1216 
1218 {
1219  // Read all identification grids from the file and add them to the IDGridManager object
1220  // used by this array. This method sets up the links between each grid and the
1221  // IDtelescope(s) it is to be used for, unlike calling
1222  //
1223  // gIDGridManager->ReadAsciiFile(grids)
1224  //
1225  // which does not.
1226  //
1227  // Returns kFALSE if there is a problem reading the file
1228 
1229  if (gIDGridManager->ReadAsciiFile(grids)) {
1230  TIter next(gIDGridManager->GetLastReadGrids());
1231  KVIDGraph* gr;
1232  while ((gr = (KVIDGraph*)next())) FillListOfIDTelescopes(gr);
1233  return kTRUE;
1234  }
1235  return kFALSE;
1236 }
1237 
1238 
1239 
1240 
1244 
1246 {
1247  // Print full status report on ID telescopes in array, using informations stored in
1248  // fStatusIDTelescopes (see GetStatusOfIDTelescopes).
1249 
1250  cout << endl << "-----STATUS OF IDENTIFICATION TELESCOPES";
1251  if (GetCurrentRunNumber()) cout << " FOR RUN "
1252  << GetCurrentRunNumber();
1253  cout << "------" << endl << endl;
1254  //get list of active telescopes
1255  KVString id_labels;
1256  if (gDataSet) id_labels = gDataSet->GetDataSetEnv<KVString>("ActiveIdentifications");
1257  else {
1258  auto typelist = GetIDTelescopeTypes();
1259  TIter it(&typelist);
1260  TObjString* type;
1261  while ((type = (TObjString*)it())) {
1262  if (id_labels == "") id_labels += type->GetString().Data();
1263  else {
1264  id_labels += Form(" %s", type->GetString().Data());
1265  }
1266  }
1267  }
1268  if (id_labels == "") {
1269  cout << " *** No active identifications *** " << endl;
1270  return;
1271  }
1272  // iterate over labels
1273  unique_ptr<TObjArray> toks(id_labels.Tokenize(' '));
1274 
1275  //update status infos
1277 
1278  TIter next_type(fStatusIDTelescopes);
1279  TList* id_type_list = 0;
1280  while ((id_type_list = (TList*)next_type())) {
1281 
1282  cout << " *** " << id_type_list->GetName() << " Identifications -------------------" << endl;
1283  if (!toks->FindObject(id_type_list->GetName())) {
1284  cout << " [NOT ACTIVE]" << endl;
1285  }
1286  TList* ok_list = (TList*)id_type_list->FindObject("OK");
1287  TList* notok_list = (TList*)id_type_list->FindObject("NOT OK");
1288  TList* print_list = 0;
1289  Int_t Nok = ok_list->GetEntries();
1290  Int_t Notok = notok_list->GetEntries();
1291  if (Nok && Notok) {
1292  if (Nok < Notok) print_list = ok_list;
1293  else print_list = notok_list;
1294  }
1295  if (Nok && (!Notok)) cout << " ALL telescopes are OK" << endl;
1296  else if (Notok && (!Nok)) cout << " NO telescopes are OK" << endl;
1297  else {
1298  cout << " " << ok_list->GetEntries() << " telescopes are OK, "
1299  << notok_list->GetEntries() << " telescopes are NOT OK" << endl;
1300  cout << " " << print_list->GetName() << " :" << endl;
1301  TIter it(print_list);
1302  TObject* ob = it();
1303  cout << ob->GetName();
1304  while ((ob = it())) cout << "," << ob->GetName();
1305  cout << endl;
1306  }
1307  cout << endl;
1308 
1309  }
1310 }
1311 
1312 
1313 
1314 
1322 
1324 {
1325  // Fill and return pointer to list fStatusIDTelescopes which contains
1326  // a list for each type of ID telescope in the array, each list contains a list
1327  // "OK" with the ID telescopes which have IsReadyForID()=kTRUE, and
1328  // a list "NOT OK" with the others.
1329  //
1330  // The returned TList object must not be deleted (it belongs to the KVMultiDetArray).
1331 
1332  if (!fStatusIDTelescopes) {
1333  fStatusIDTelescopes = new TList;
1335  }
1336  else {
1338  }
1340  TIter next(fIDTelescopes);
1341  KVIDTelescope* idt = 0;
1342  while ((idt = (KVIDTelescope*)next())) {
1343 
1344  TString id_type = idt->GetLabel();
1345  TList* id_type_list = (TList*)fStatusIDTelescopes->FindObject(id_type.Data());
1346  if (!id_type_list) {
1347  id_type_list = new TList;
1348  id_type_list->SetOwner(kTRUE);
1349  id_type_list->SetName(id_type.Data());
1350  fStatusIDTelescopes->Add(id_type_list);
1351  id_type_list->Add(new TList);
1352  ((TList*)id_type_list->At(0))->SetName("OK");
1353  id_type_list->Add(new TList);
1354  ((TList*)id_type_list->At(1))->SetName("NOT OK");
1355  }
1356  if (idt->IsReadyForID())
1357  ((TList*)id_type_list->FindObject("OK"))->Add(idt);
1358  else
1359  ((TList*)id_type_list->FindObject("NOT OK"))->Add(idt);
1360  }
1361  return fStatusIDTelescopes;
1362 }
1363 
1364 
1365 
1366 
1370 
1372 {
1373  // Create, fill and return pointer to a list of TObjString containing the name of each type
1374  // of ID telescope (actually the label) in the array.
1375 
1376  KVUniqueNameList type_list(kTRUE);
1377  type_list.SetOwner();
1378  if (!fIDTelescopes || !fIDTelescopes->GetEntries()) return type_list;
1379  TIter next(fIDTelescopes);
1380  KVIDTelescope* idt = 0;
1381  while ((idt = (KVIDTelescope*)next())) {
1382  type_list.Add(new TObjString(idt->GetLabel()));
1383  }
1384  return type_list;
1385 }
1386 
1387 
1388 
1389 
1397 
1399 {
1400  // Create, fill and return pointer to a list of KVIDTelescopes with
1401  // the given type (label) in the array.
1402  // WARNING! - check pointer is not zero (we return NULL if ID telescopes
1403  // list is not defined or empty)
1404  //
1405  // Delete the KVList after use (it does not own the KVIDTelescopes).
1406 
1407  if (!fIDTelescopes || !fIDTelescopes->GetEntries()) return NULL;
1409 }
1410 
1411 
1412 
1413 
1421 
1423 {
1424  // Fill and return pointer to list fCalibStatusDets which contains
1425  // a list for each type of detector in the array, each list contains a list
1426  // "OK" with the detectors which are calibrated, and
1427  // a list "NOT OK" with the others.
1428  //
1429  // The returned TList object must not be deleted (it belongs to the KVMultiDetArray).
1430 
1431  if (!fCalibStatusDets) {
1432  fCalibStatusDets = new TList;
1434  }
1435  else {
1437  }
1438  if (!GetDetectors()->GetEntries()) return fCalibStatusDets;
1439  TIter next(GetDetectors());
1440  KVDetector* det = 0;
1441  while ((det = (KVDetector*)next())) {
1442 
1443  TString type = det->GetType();
1444  TList* type_list = (TList*)fCalibStatusDets->FindObject(type.Data());
1445  if (!type_list) {
1446  type_list = new TList;
1447  type_list->SetOwner(kTRUE);
1448  type_list->SetName(type.Data());
1449  fCalibStatusDets->Add(type_list);
1450  type_list->Add(new TList);
1451  ((TList*)type_list->At(0))->SetName("OK");
1452  type_list->Add(new TList);
1453  ((TList*)type_list->At(1))->SetName("NOT OK");
1454  }
1455  if (det->IsCalibrated())
1456  ((TList*)type_list->FindObject("OK"))->Add(det);
1457  else
1458  ((TList*)type_list->FindObject("NOT OK"))->Add(det);
1459  }
1460  return fCalibStatusDets;
1461 }
1462 
1463 
1464 
1465 
1469 
1471 {
1472  // Print full status report on calibration of detectors in array, using informations stored in
1473  // fCalibStatusDets (see GetCalibrationStatusOfDetectors).
1474 
1475  if (!GetCurrentRunNumber()) {
1476  Info("PrintCalibStatusOfDetectors", "Cannot know status without knowing RUN NUMBER");
1477  return;
1478  }
1479 
1480  cout << endl << "-----------STATUS OF CALIBRATIONS FOR RUN "
1481  << GetCurrentRunNumber() << "------------" << endl << endl;
1482 
1483  //update status infos
1485 
1486  TIter next_type(fCalibStatusDets);
1487  TList* id_type_list = 0;
1488  while ((id_type_list = (TList*)next_type())) {
1489 
1490  cout << " *** " << id_type_list->GetName() << " Detectors -------------------" << endl;
1491  TList* ok_list = (TList*)id_type_list->FindObject("OK");
1492  TList* notok_list = (TList*)id_type_list->FindObject("NOT OK");
1493  TList* print_list = 0;
1494  Int_t Nok = ok_list->GetEntries();
1495  Int_t Notok = notok_list->GetEntries();
1496  if (Nok && Notok) {
1497  if (Nok < Notok) print_list = ok_list;
1498  else print_list = notok_list;
1499  }
1500  if (Nok && (!Notok)) cout << " ALL calibrations are OK" << endl;
1501  else if (Notok && (!Nok)) cout << " NO calibrations are OK" << endl;
1502  else {
1503  cout << " " << ok_list->GetEntries() << " calibrations are OK, "
1504  << notok_list->GetEntries() << " calibrations are NOT OK" << endl;
1505  cout << " " << print_list->GetName() << " :" << endl;
1506  TIter it(print_list);
1507  TObject* ob = it();
1508  cout << ob->GetName();
1509  while ((ob = it())) cout << "," << ob->GetName();
1510  cout << endl;
1511  }
1512  cout << endl;
1513 
1514  }
1515 }
1516 
1517 
1518 
1519 
1535 
1537 {
1538  // Calculate the energy loss in the current target of the multidetector
1539  // for the reconstructed charged particle 'ion', assuming that the current
1540  // energy and momentum of this particle correspond to its state on
1541  // leaving the target.
1542  //
1543  // WARNING: for this correction to work, the target must be in the right 'state':
1544  //
1545  // gMultiDetArray->GetTarget()->SetIncoming(kFALSE);
1546  // gMultiDetArray->GetTarget()->SetOutgoing(kTRUE);
1547  //
1548  // (see KVTarget::GetParticleEIncFromERes).
1549  //
1550  // The returned value is the energy lost in the target in MeV.
1551  // The energy/momentum of 'ion' are not affected.
1552 
1553  if (fTarget && ion) return (fTarget->GetParticleEIncFromERes(ion) - ion->GetEnergy());
1554  return 0;
1555 }
1556 
1557 
1558 
1559 
1562 
1564 {
1565  // Return pointer to the (ROOT) geometry of the array.
1566  return gGeoManager;
1567 }
1568 
1569 
1570 
1572 
1574 {
1575  return fNavigator.get();
1576 }
1577 
1578 
1579 
1598 
1600 {
1601  // Actual thicknesses of detectors can be given in one or more files associated with a dataset/multidetector.
1602  //
1603  // We look for the first of the files
1604  //~~~~
1605  //[array_name].DetectorThicknessFiles.dat
1606  //DetectorThicknessFiles.dat
1607  //~~~~
1608  // which exists for the current dataset, and if found we read each file listed in it.
1609  //
1610  // Otherwise we look for the first of the files
1611  //~~~~
1612  //[array_name].DetectorThicknesses.dat
1613  //DetectorThicknesses.dat
1614  //~~~~
1615  // which exists for the current dataset, and if found we read the thicknesses from it.
1616  //
1617  // See set_detector_thicknesses() for details of the format of the individual thickness files.
1618 
1619  auto find_a_file = [this](const TString & _base_filename) {
1620  TString filename = TString(GetName()) + "." + _base_filename;
1621  TString fullpath;
1624 
1625  filename = _base_filename;
1628 
1629  return TString();
1630  };
1631 
1632  // look for list of files
1633  auto fullpath = find_a_file("DetectorThicknessFiles.dat");
1634  if (!fullpath.IsNull()) {
1635  KVFileReader fr;
1636  fr.OpenFileToRead(fullpath);
1637  while (fr.IsOK()) {
1638  fr.ReadLine(0);
1639  if (fr.GetCurrentLine().BeginsWith("#") || fr.GetCurrentLine() == "") continue;
1641  }
1642  return;
1643  }
1644 
1645  // look for single file
1646  fullpath = find_a_file("DetectorThicknesses.dat");
1647  if (!fullpath.IsNull()) set_detector_thicknesses(fullpath);
1648 }
1649 
1650 
1651 
1670 
1672 {
1673  // Use file given by fullpath to set the real thicknesses of the detectors.
1674  // Any detector which is not in the file will be left with its nominal thickness.
1675  //
1676  // EXAMPLE FILE:
1677  //
1678  //# thickness of detector DET01 in default units
1679  //DET01: 56.4627
1680  //
1681  //# DET03 has several layers
1682  //DET03.Abs0: 61.34
1683  //DET03.Abs1: 205.62
1684  //
1685  // !!! WARNING !!!
1686  // Single-layer detectors: The units are those defined by default for the detector's
1687  // Get/SetThickness methods.
1688  // Multi-layer: Each layer is a KVMaterial object. The thickness MUST be given in centimetres
1689  // (default thickness unit for KVMaterial).
1690 
1691  TEnv thickdat;
1692  if (thickdat.ReadFile(fullpath, kEnvUser) != 0) {
1693  Error("SetDetectorThicknesses", "Problem opening file %s", fullpath.Data());
1694  return;
1695  }
1696  Info("SetDetectorThicknesses", "Setting thicknesses of detectors from file %s", fullpath.Data());
1697  TIter next(GetDetectors());
1698  KVDetector* det;
1699  while ((det = (KVDetector*)next())) {
1700  if (thickdat.Defined(det->GetName())) {
1701  // simple single layer detector
1702  Double_t thick = thickdat.GetValue(det->GetName(), 0.0);
1703  det->SetThickness(thick);
1704  //Info("SetDetectorThicknesses", "Set thickness of %s to %f", det->GetName(), thick);
1705  }
1706  else {
1707  Char_t i = 0;
1708  TString absname;
1709  absname.Form("%s.Abs%d", det->GetName(), (Int_t)i);
1710  if (thickdat.Defined(absname.Data())) {
1711  // detector with several layers
1712  KVMaterial* abs = 0;
1713  while ((abs = det->GetAbsorber(i))) {
1714  Double_t thick = thickdat.GetValue(absname.Data(), 0.0);
1715  abs->SetThickness(thick);
1716  //Info("SetDetectorThicknesses", "Set thickness of %s.Abs%d to %f", det->GetName(), (Int_t)i, thick);
1717  i++;
1718  absname.Form("%s.Abs%d", det->GetName(), (Int_t)i);
1719  if (!thickdat.Defined(absname.Data())) break;
1720  }
1721  }
1722  }
1723  }
1724 }
1725 
1726 
1727 
1736 
1738 {
1739  // Called by KVGeoImport::ImportGeometry
1740  //
1741  // Creates KVRangeTableGeoNavigator for calculating energy losses of
1742  // particles propagated through array.
1743  //
1744  // If no name and/or title are defined for the array, the name and title of the TGeoManager
1745  // object will be used for the array.
1746 
1747  if (!strcmp(GetName(), "")) SetName(g->GetName());
1748  if (!strcmp(GetTitle(), "")) SetTitle(g->GetTitle());
1750 }
1751 
1752 
1753 
1782 
1784 {
1785  // Method for positioning volumes in detector geometries
1786  //
1787  // Given:
1788  //
1789  // distance [cm] = distance from target (origin) to the CENTRE of the volume in position
1790  // theta [deg] = polar angle of vector from target to centre of volume in position
1791  // phi [deg] = azimuthal angle of vector
1792  //
1793  // this method generates the matrix which is required to position the volume as required
1794  // while also turning the volume so that the side nearest the target (i.e. the entrance
1795  // window of the detector) remains perpendicular to the vector joining the origin and
1796  // the centre of the volume.
1797  //
1798  // If required, a further translation can be given which will be applied to the volume after
1799  // it has been placed with the required orientation at the nominal distance. This can be used
1800  // e.g. for detector misalignment, when detectors are in a structure which guarantees their line
1801  // of sight to be orthogonal to their surface at a nominal distance, but the nominal distance
1802  // is not respected.
1803  //
1804  // Example of use:
1805  //
1806  //~~~~~~~~~~~~
1807  // TGeoVolume* vol;// volume to be positioned
1808  // double depth = vol->GetShape()->GetDZ(); // half-width of volume in direction of target
1809  // // place front of volume at 100cm, with theta=45 deg. and phi=60 deg.
1810  // gGeoManager->GetTopVolume()->AddNode(vol, 1, KVMultiDetArray::GetVolumePositioningMatrix(100+depth,45,60));
1811  //~~~~~~~~~~~~
1812 
1813  TGeoRotation rot;
1814  TGeoTranslation trans;
1815  rot.SetAngles(phi + 90, theta, -90) ;
1816  trans.SetDz(distance) ;
1817  TGeoHMatrix h;
1818  if (postTrans) h = (*postTrans) * rot * trans ;
1819  else h = rot * trans;
1820  TGeoHMatrix* ph = new TGeoHMatrix(h);
1821  return ph;
1822 }
1823 
1824 
1825 
1826 
1830 
1832 {
1833  // For each grid which is valid for this run, we call the KVIDTelescope::SetIDGrid method
1834  // of each associated ID telescope.
1835 
1836  if (gIDGridManager->IsSQLROOT()) {
1837  gIDGridManager->LoadGridsForRun(run);
1838  }
1839 
1840  TIter next(gIDGridManager->GetGrids());
1841  KVIDGraph* gr = 0;
1842  while ((gr = (KVIDGraph*) next())) {
1843  if (!gIDGridManager->IsSQLROOT()) {
1844  if (!gr->GetRuns().Contains((Int_t) run))
1845  continue;
1846  }
1847  else
1849 
1850  TIter nxtid(gr->GetIDTelescopes());
1851  KVIDTelescope* idt;
1852  while ((idt = (KVIDTelescope*) nxtid())) {
1853  idt->SetIDGrid(gr);
1854  }
1855  }
1856 }
1857 
1858 
1862 
1864 {
1865  // Fill list of ID telescopes with which this grid is associated
1866  // from list of names read from ascii file.
1867 
1868  gr->ClearListOfTelescopes();
1869  if (gr->GetParameters()->HasParameter("IDTelescopes")) {
1870  KVString tel_list = gr->GetParameters()->GetStringValue("IDTelescopes");
1871  tel_list.Begin(",");
1872  while (!tel_list.End()) {
1873  TString tel_name = tel_list.Next();
1874  KVIDTelescope* idt = GetIDTelescope(tel_name.Data()) ;
1875  if (idt) gr->AddIDTelescope(idt);
1876  }
1877  }
1878 }
1879 
1880 
1881 
1891 
1893 {
1894  // Use OpenGL viewer to view multidetector geometry
1895  //
1896  // If option="tracks" we draw any tracks corresponding to the last simulated
1897  // event whose detection was simulated with DetectEvent
1898  // If option="tracks:[numberlist]" with a list of numbers,
1899  // it will be interpreted as a KVNumberList containing the Z of tracks to be drawn
1900  // e.g. option="tracks:1-92" draw only tracks with 1<=Z<=92 (no neutrons)
1901  // option="tracks:2" draw only helium isotopes
1902 
1903  GetGeometry()->GetTopVolume()->Draw("ogl");
1904  KVString opt(option);
1905  opt.Begin(":");
1906  if (opt.Next() == "tracks") {
1907  if (!opt.End()) {
1908  KVNumberList zlist(opt.Next());
1909  GetNavigator()->DrawTracks(&zlist);
1910  }
1911  else
1912  GetNavigator()->DrawTracks();
1913  }
1914 #ifdef WITH_OPENGL
1915  TGLViewer* view = (TGLViewer*)gPad->GetViewer3D();
1918  view->SetSmoothLines(kTRUE);
1919  view->SetSmoothPoints(kTRUE);
1920 #endif
1921 }
1922 
1923 
1924 
1926 
1928 {
1929  fNavigator.reset(geo);
1930 }
1931 
1932 
1933 
1938 
1940 {
1941  // Create TH2F histograms for all IDTelescopes of the array
1942  // They will be added to the list
1943  // histograms will have resolution of dimension*dimension
1944 
1946  KVIDTelescope* idt;
1947  while ((idt = (KVIDTelescope*)it())) {
1948  TString name(idt->GetName());
1949  name.ReplaceAll("-", "_");
1950  list->Add(new TH2F(name, Form("Hits in %s", idt->GetName()), dimension, 0., 0., dimension, 0., 0.));
1951  }
1952 }
1953 
1954 
1955 
1958 
1960 {
1961  // Fill TH2F histograms for all IDTelescopes of the array
1962 
1964  KVIDTelescope* idt;
1965  while ((idt = (KVIDTelescope*)it())) {
1966  TString name(idt->GetName());
1967  name.ReplaceAll("-", "_");
1968  TH2F* h = (TH2F*)list->FindObject(name);
1969  if (h) h->Fill(idt->GetIDMapX(), idt->GetIDMapY());
1970  }
1971 }
1972 
1973 
1974 
1977 
1979 {
1980  // Modify the transparency of detector volumes in OpenGL view
1981 
1982  TIter itV(GetGeometry()->GetListOfVolumes());
1983  TGeoVolume* vol;
1984  while ((vol = (TGeoVolume*)itV())) vol->SetTransparency(t);
1985 }
1986 
1987 
1988 
1992 
1994 {
1995  // Calculate all possible (sub-)trajectories
1996  // for particle reconstruction (GetReconTrajectories())
1997 
1998  unique_ptr<KVSeqCollection> groups(GetStructureTypeList("GROUP"));
1999  TIter it(groups.get());
2000  KVGroup* group;
2001  Int_t ntr = 0;
2002  Info("CalculateReconstructionTrajectories", "Calculating trajectories for particle reconstruction:");
2004  std::cout << "\xd" << " -- calculated " << ntr << " reconstruction trajectories" << std::flush;
2005  while ((group = (KVGroup*)it())) {
2006  ntr += group->CalculateReconstructionTrajectories();
2008  std::cout << "\xd" << " -- calculated " << ntr << " reconstruction trajectories" << std::flush;
2009  }
2011  std::cout << " -- calculated " << ntr << " reconstruction trajectories" << std::endl;
2012  else
2013  std::cout << std::endl;
2014 }
2015 
2016 
2017 
2022 
2024 {
2025  // Track over all possible particle trajectories calling GetIDTelescopes(KVDetector*,KVDetector*)
2026  // for each pair of (present & functioning) detectors.
2027  // This will create all possible KVIDTelescope identification objects and put them in list fIDTelescopes
2028 
2029  fIDTelescopes->Delete();
2030  TIter next_traj(GetTrajectories());
2031  KVGeoDNTrajectory* traj;
2032  Int_t count = 0;
2033  Info("DeduceIdentificationTelescopesFromGeometry", "Calculating...");
2035  std::cout << "\xd" << " -- created " << count << " telescopes" << std::flush;
2036  while ((traj = (KVGeoDNTrajectory*)next_traj())) { // loop over all trajectories
2037 
2038  traj->IterateFrom(); // from furthest-out to closest-in detector
2039 
2041  while ((N = traj->GetNextNode())) {
2042  KVGeoDetectorNode* Nplus1 = traj->GetNodeInFront(N);
2043  count += GetIDTelescopes(Nplus1, N->GetDetector(), traj->AccessIDTelescopeList());
2045  std::cout << "\xd" << " -- created " << count << " telescopes" << std::flush;
2046  }
2047  }
2049  std::cout << " -- created " << count << " telescopes" << std::endl;
2050  else
2051  std::cout << std::endl;
2052 }
2053 
2054 
2055 
2059 
2061 {
2062  // Eliminate any trajectories which are just sub-trajectories of others
2063  // For each trajectory in list fTrajectories, we add a reference to the trajectory to each node on the trajectory
2064 
2065  TIter it(&fTrajectories);
2066  KVGeoDNTrajectory* tr;
2067  KVList duplicates;
2068  // look for duplicate sub-trajectories
2069  while ((tr = (KVGeoDNTrajectory*)it())) {
2070  int len_tr = tr->GetN();
2071  TIter it2(&fTrajectories);
2072  KVGeoDNTrajectory* tr2;
2073  while ((tr2 = (KVGeoDNTrajectory*)it2())) {
2074  if ((tr2 != tr) && (len_tr < tr2->GetN()) && (tr2->ContainsPath(tr))) {
2075  duplicates.Add(tr);
2076  break;
2077  }
2078  }
2079  }
2080  // remove duplicates
2081  if (duplicates.GetEntries()) {
2082  TIter it_dup(&duplicates);
2083  while ((tr = (KVGeoDNTrajectory*)it_dup())) {
2084  fTrajectories.Remove(tr);
2085  }
2086  Info("AssociateTrajectoriesAndNodes", "Removed %d duplicated sub-trajectories", duplicates.GetEntries());
2087  }
2088  Info("AssociateTrajectoriesAndNodes", "Calculated %d particle trajectories", fTrajectories.GetEntries());
2089  it.Reset();
2090  while ((tr = (KVGeoDNTrajectory*)it())) {
2091  tr->AddToNodes();
2092  }
2093 }
2094 
2095 
2096 
2098 
2100 {
2101  if (N->GetNTraj() > 1) {
2102  if (!multitraj_nodes.FindObject(N)) { // look for any detectors which are on multiple trajectories
2103  //cout << "multitraj node found: " << N->GetName() << " (" << N->GetNTraj() << ")" << endl;
2104  multitraj_nodes.Add(N);
2105  TIter tr(N->GetTrajectories());
2106  KVGeoDNTrajectory* traj;
2107  while ((traj = (KVGeoDNTrajectory*)tr())) { // for each trajectory associated with detector
2108  if (tried_trajectories.FindObject(traj)) continue; // trajectory already used
2109  tried_trajectories.Add(traj);
2110  traj->IterateFrom();
2111  KVGeoDetectorNode* node;
2112  while ((node = traj->GetNextNode())) { // store names of all detectors on trajectory
2113  detectors_of_group.Add(node);
2114  RecursiveTrajectoryClustering(node, tried_trajectories, multitraj_nodes, detectors_of_group);
2115  }
2116  }
2117  }
2118  }
2119  else if (N->GetNTraj() == 1) {
2120  // single-trajectory node.
2121  // work along trajectory adding nodes to group
2122  KVGeoDNTrajectory* traj = (KVGeoDNTrajectory*)N->GetTrajectories()->First();
2123  if (tried_trajectories.FindObject(traj)) return; // trajectory already used
2124  tried_trajectories.Add(traj);
2125  traj->IterateFrom();
2126  KVGeoDetectorNode* node;
2127  while ((node = traj->GetNextNode())) { // store names of all detectors on trajectory
2128  detectors_of_group.Add(node);
2129  RecursiveTrajectoryClustering(node, tried_trajectories, multitraj_nodes, detectors_of_group);
2130  }
2131  }
2132  else {
2133  // orphan node? single-detector array?
2134  detectors_of_group.Add(N);
2135  }
2136 }
2137 
2138 
2139 
2148 
2150 {
2151  // Create and return pointer to new KVGroupReconstructor for reconstructing particles
2152  // in the given group. Returns nullptr if group is not part of this array.
2153  //
2154  // Plugins for specific arrays can be defined as plugins using the name of the array:
2155  // +Plugin.KVGroupReconstructor: my_array my_group_reconstructor my_lib "my_group_reconstructor()"
2156  //
2157  // If we are in 'SimMode', the default reconstructor is a KVFilterGroupReconstructor
2158 
2159  KVGroupReconstructor* gr(nullptr);
2160  if (GetGroup(g->GetName())) {
2161  // look for plugin
2163  if (!gr) {
2164  if (IsSimMode()) return KVGroupReconstructor::Factory("Filter", g);
2165  else
2166  gr = new KVGroupReconstructor(g);
2167  }
2168  }
2169  return gr;
2170 }
2171 
2172 
2173 
2177 
2179 {
2180  // Call this method just after opening a raw data file in order to perform any
2181  // necessary initialisations, depending on the type of data
2182 
2183 #ifdef WITH_BUILTIN_GRU
2184  if (r->GetDataFormat() == "EBYEDAT")
2185  dynamic_cast<KVGANILDataReader*>(r)->ConnectRawDataParameters();
2186 #endif
2187 }
2188 
2189 
2190 
2195 
2197 {
2198  // Deduce the "groups" in the array from the trajectories
2199  // Any trajectories with 1 or more common detectors define a group.
2200  // The group is constituted of all detectors belonging to the trajectories of the group.
2201 
2202  Info("DeduceGroupsFromTrajectories", "Deducing groups of detectors from trajectories");
2203  Int_t number_of_groups = 0;
2204  TIter next_det(GetDetectors());
2205  unique_ptr<KVSeqCollection> stl(GetStructureTypeList("GROUP"));
2206  if (stl.get() && stl->GetEntries()) {
2207  Info("DeduceGroupsFromTrajectories", "Deleting existing %d groups in array", stl->GetEntries());
2208  ClearStructures("GROUP");
2209  Info("DeduceGroupsFromTrajectories", "Done");
2210  }
2211  KVDetector* det;
2212  KVUniqueNameList tried_trajectories;//avoid double-counting/infinite loops
2213  KVUniqueNameList multitraj_nodes;//avoid double-counting/infinite loops
2214  while ((det = (KVDetector*) next_det())) {
2215  if (det->GetGroup()) continue; // group assignment already done
2216  KVUniqueNameList detectors_of_group;
2217  RecursiveTrajectoryClustering(det->GetNode(), tried_trajectories, multitraj_nodes, detectors_of_group);
2218  if (!detectors_of_group.GetEntries()) continue;
2219  KVGroup* Group = new KVGroup;
2220  Group->SetNumber(++number_of_groups);
2221  Add(Group);
2222  TIter next_node(&detectors_of_group);
2224  while ((d = (KVGeoDetectorNode*)next_node())) Group->Add(d->GetDetector());
2225  }
2226  TIter tr(&fTrajectories);
2227  KVGeoDNTrajectory* t;
2228  Info("DeduceGroupsFromTrajectories", "Filling group trajectory lists");
2229  while ((t = (KVGeoDNTrajectory*)tr())) {
2230  if (t->GetNodeAt(0)->GetDetector()->GetGroup())
2231  t->GetNodeAt(0)->GetDetector()->GetGroup()->AddTrajectory(t);
2232  else {
2233  t->Print();
2234  t->GetNodeAt(0)->GetDetector()->Print();
2235  }
2236  }
2237 }
2238 
2239 
2240 
2244 
2246 {
2247  // Called when required to fill KVReconstructedNucleus::fDetList with pointers to
2248  // the detectors whose names are stored in KVReconstructedNucleus::fDetNames.
2249 
2250  DetList->Clear();
2251  DetNames.Begin("/");
2252  while (!DetNames.End()) {
2253  KVDetector* det = GetDetector(DetNames.Next(kTRUE));
2254  if (det) DetList->Add(det);
2255  }
2256 }
2257 
2258 
2259 
2279 
2281 {
2282  // Set status of particle by comparing its identification/calibration codes
2283  // with those set as acceptable in fAcceptIDCodes and fAcceptECodes.
2284  // The status can be tested with method KVReconstructedNucles::IsOK().
2285  //
2286  // The default lists are defined in variables of the form
2287  //```
2288  // [DataSet].[name].ReconstructedNuclei.AcceptIDCodes: [list]
2289  // [DataSet].[name].ReconstructedNuclei.AcceptECodes: [list]
2290  //```
2291  // where:
2292  // + `DataSet` is an optional dataset name for dataset-specific lists
2293  // + `name` is the name of the multidetector array
2294  // + `list` is a numeric list (KVNumberList format)
2295  //
2296  // The default lists can be overridden using methods AcceptIDCodes(), AcceptECodes(),
2297  // AcceptAllIDCodes() and AcceptAllECodes().
2298  //
2299  // If either list is empty, no selection is made for the corresponding code
2300 
2301  Bool_t ok = kTRUE;
2303  if (!fAcceptECodes.IsEmpty()) ok = ok && fAcceptECodes.Contains(NUC->GetECode());
2304  NUC->SetIsOK(ok);
2305 }
2306 
2307 
2308 #ifdef WITH_BUILTIN_GRU
2309 
2312 
2314 {
2315  // General method for reading raw data in old GANIL ebyedat format
2316  AbstractMethod("handle_raw_data_event_ebyedat");
2317  return kFALSE;
2318 }
2319 
2320 #endif
2321 
2322 
2325 
2327 {
2328  // reset acquisition parameters etc. before reading new raw data event
2329 
2332  fHandledRawData = false;
2333  // reset fired signals
2334  TIter nxt(&fFiredSignals);
2335  KVDetectorSignal* ds;
2336  while ((ds = (KVDetectorSignal*)nxt())) {
2337  ds->SetFired(false);
2338  ds->SetValue(0);
2339  }
2340  fFiredSignals.Clear();
2341 }
2342 
2343 
2344 
2352 
2354 {
2355  // Perform any operations to finalise the description of the multidetector
2356  // which can only be done once the geometry is closed, e.g. use KVGeoImport
2357  // to set up nodes, trajectories, detectors, idtelescopes, etc.
2358  // This has to be kept separate for use with KVExpSetUp which first fills
2359  // a single ROOT geometry with all component KVMultiDetArray geometries,
2360  // then closes the geometry only when all have been built.
2361 }
2362 
2363 
2364 
2367 
2369 {
2370  // Copy any parameters in fReconParameters in to the reconstructed event parameter list
2371  e->GetParameters()->Concatenate(fReconParameters);
2372 }
2373 
2374 
2375 
2388 
2390 {
2391  // values of fired raw data signals (acquisition parameters) from last read raw event
2392  // are copied to the fReconParameters list of parameters to be stored with the
2393  // reconstructed event.
2394  //
2395  // the format for each signal is:
2396  //
2397  // ACQPAR.[array].[detector].[signal]
2398  // ACQPAR.[array].[signal]
2399  //
2400  // in the first case for signals associated with detectors, in the latter case signals
2401  // which are not associated with a detector
2402 
2403  TIter it(GetFiredSignals());
2404  KVDetectorSignal* o;
2405  while ((o = (KVDetectorSignal*)it())) {
2406  fReconParameters.SetValue(Form("ACQPAR.%s.%s", GetName(), o->GetFullName().Data()), o->GetValue());
2407  }
2408 }
2409 
2410 
2411 
2420 
2422 {
2423  // Update array according to last event read using the KVRawDataReader object
2424  // (it is assumed that KVRawDataReader::GetNextEvent() was called before calling this method)
2425  //
2426  // Return kTRUE if raw data was treated
2427  //
2428  // All fired acquisition parameters are written in the fReconParameters list,
2429  // ready to be copied to the reconstructed event
2430 
2431  fRawDataReader = rawdata;
2433  if (rawdata->GetDataFormat() == "MFM") {
2434 #ifdef WITH_MFM
2436 #endif
2437  }
2438  else if (rawdata->GetDataFormat() == "PROTOBUF") {
2439 #ifdef WITH_PROTOBUF
2441 #endif
2442  }
2443  else if (rawdata->GetDataFormat() == "EBYEDAT") {
2444 #ifdef WITH_BUILTIN_GRU
2446 #endif
2447  }
2448  if (fHandledRawData) {
2450  }
2451  return fHandledRawData;
2452 }
2453 
2454 
2455 #ifdef WITH_MFM
2456 
2465 
2467 {
2468  // Update array according to last event read from MFM buffer
2469  // (it is assumed that MFMBufferReader::ReadNextFrame() was called before calling this method)
2470  //
2471  // Return kTRUE if raw data was treated
2472  //
2473  // All fired acquisition parameters are written in the fReconParameters list,
2474  // ready to be copied to the reconstructed event
2475 
2477  bool ok = false;
2478  ok = handle_raw_data_event_mfmfile(bufrdr);
2479  if (ok) {
2481  }
2482  return ok;
2483 }
2484 
2485 #endif
2486 
2487 
2498 
2500 {
2501  // Given a pointer to a detector (may be nullptr) try to set the data sig_data in the associated signal
2502  // of the given type, sig_typ.
2503  //
2504  // If the signal does not exist for the detector, it will be created.
2505  //
2506  // If detector pointer is null, the signal will be looked for in fExtraRawDataSignals
2507  // and created if necessary.
2508  //
2509  // All fired detectors and signals are added to the lists fFiredDetectors and fFiredSignals.
2510 
2511  KVDetectorSignal* det_signal = nullptr;
2512  if (detector) {
2513  det_signal = detector->GetDetectorSignal(sig_type);
2514  if (!det_signal) {
2515  det_signal = detector->AddDetectorSignal(sig_type);
2516  }
2517  fFiredDetectors.Add(detector);
2518  }
2519  else {
2520  // raw data not associated with a detector
2521  TString sig_name;
2522  if (detname != "") sig_name = Form("%s.%s", detname.Data(), sig_type.Data());
2523  else sig_name = sig_type;
2524  det_signal = fExtraRawDataSignals.get_object<KVDetectorSignal>(sig_name);
2525  if (!det_signal) {
2526  det_signal = new KVDetectorSignal(sig_name);
2527  fExtraRawDataSignals.Add(det_signal);
2528  }
2529  }
2530  if (det_signal) {
2531  det_signal->SetValue(sig_data);
2532  det_signal->SetFired();
2533  fFiredSignals.Add(det_signal);
2534  }
2535 }
2536 
2537 
2538 
2544 
2546 {
2547  // Take values 'ACQPAR.[array_name].[detname].[signal]' or 'ACQPAR.[array_name].[signal]'
2548  // in the parameter list and use them to set values of raw acquisition parameters.
2549  //
2550  // Any detector signals which don't already exist will be created
2551 
2552  prepare_to_handle_new_raw_data(); // clear previous fired parameters/detectors
2553  int N = l.GetNpar();
2554  for (int i = 0; i < N; ++i) {
2555  KVNamedParameter* np = l.GetParameter(i);
2556 
2557  KVString name(np->GetName());
2558  if (name.BeginsWith("ACQPAR")) {
2559  // 3 '.' => 4 values means associated detector
2560  // 2 '.' => 3 values means no detector
2561  int dots = name.GetNValues(".");
2562  bool with_det = (dots == 4);
2563  assert(with_det || (dots == 3)); // sanity check
2564  name.Begin(".");
2565  name.Next(); // "ACQPAR"
2566  if (name.Next() != GetName()) continue; // check name of array - somebody else's parameter ?
2567  KVString det_name;
2568  KVString sig_type;
2569  KVDetector* det = nullptr;
2570  if (with_det) {
2571  det_name = name.Next();
2572  sig_type = name.Next();
2573  det = GetDetector(det_name);
2574  }
2575  else {
2576  sig_type = name.Next();
2577  }
2578  add_and_set_detector_signal(det, det_name, np->GetDouble(), sig_type);
2579  }
2580  }
2581 }
2582 
2583 
2584 
2603 
2605 {
2606  // We first look for following files with the name given by
2607  //
2608  // [dataset].[name].OoODetectors: [name.OoODetectors.dat]
2609  // which should contain the runlists for each malfunctioning detector.
2610  // If found we add to the experiment database a table '[name].OoO Detectors' where [name] is the name of this array.
2611  //
2612  // [dataset].[name].AbsentDetectors: [name.AbsentDetectors.dat]
2613  // which should contain the runlists for each absent detector.
2614  // If found we add to the experiment database a table '[name].Absent Detectors' where [name] is the name of this array.
2615  //
2616  // Then we look for a file with the name given by
2617  //
2618  // [dataset].[name].CalibrationFiles: [CalibrationFiles.dat]
2619  //
2620  // which should contain the names of files to read with each type of calibration
2621  // If found we add to the experiment database a table '[name].Calibrations' where [name] is the name of this array,
2622  // containing all calibrations as KVDBParameterSet objects with the name of the detector concerned.
2623  db->SetDBType(Form("%sDB", GetName()));
2624  ReadOoODetectors(db);
2625  ReadAbsentDetectors(db);
2627 }
2628 
2629 
2630 
2632 
2634 {
2635  TString basic_name = db->GetCalibFileName(keyw);
2636  if (basic_name == "") {
2637  Info(meth, "No name found for \"%s\" file", keyw);
2638  return "";
2639  }
2640  Info(meth, "Search for %s for dataset %s ...", basic_name.Data(), fDataSet.Data());
2641  TString fp;
2642  SearchKVFile(basic_name.Data(), fp, fDataSet);
2643  if (fp == "") {
2644  Info(meth, "\tNo file found ...");
2645  }
2646  return fp;
2647 }
2648 
2649 
2650 
2652 
2653 unique_ptr<KVFileReader> KVMultiDetArray::GetKVFileReader(KVExpDB* db, const Char_t* meth, const Char_t* keyw)
2654 {
2655 
2656  TString fp = GetFileName(db, meth, keyw);
2657  if (fp == "")
2658  return unique_ptr<KVFileReader>();
2659 
2660  unique_ptr<KVFileReader> fr(new KVFileReader());
2661  if (!fr->OpenFileToRead(fp.Data())) {
2662  Error(meth, "Error in opening file %s", fp.Data());
2663  fr.reset(nullptr);
2664  }
2665  else
2666  Info(meth, "Reading %s file", fp.Data());
2667  return fr;
2668 }
2669 
2670 
2671 
2673 
2675 {
2676 
2677  unique_ptr<KVFileReader> fr = GetKVFileReader(db, "ReadCalibrationFiles()", "CalibrationFiles");
2678  if (!fr.get())
2679  return;
2680 
2681  KVDBTable* calib_table = db->AddTable(Form("%s.Calibrations", GetName()), Form("Calibrations for %s", GetName()));
2682  while (fr->IsOK()) {
2683  fr->ReadLine(0);
2684  if (fr->GetCurrentLine().BeginsWith("#") || fr->GetCurrentLine() == "") {}
2685  else {
2686  ReadCalibFile(fr->GetCurrentLine().Data(), db, calib_table);
2687  }
2688  }
2689  fr->CloseFile();
2690 }
2691 
2692 
2693 
2736 
2737 void KVMultiDetArray::ReadCalibFile(const Char_t* filename, KVExpDB* db, KVDBTable* calib_table)
2738 {
2739  // Read a calibration file with the format
2740  //
2741  //~~~~~~~~~~~~~
2742  // RunList: 1546-7485
2743  // SignalIn: PG
2744  // SignalOut: Volts
2745  // CalibType: ChannelVolt
2746  // CalibOptions: func=pol3,min=0,max=1
2747  // ZRange: 2-92
2748  // [detector1]: 0.0,0.261829,0.0
2749  // [detector2]: 0.1 0.539535 1.2
2750  //~~~~~~~~~~~~~
2751  //
2752  //The `[RunList]` is optional: if not given, the calibration will be applied to all runs in the database.
2753  //
2754  //If different parameters are required for different sets of runs, they should be written in different
2755  //files (all of which are listed in `CalibrationFiles.dat` or `[array].CalibrationFiles.dat`).
2756  //
2757  //Numerical parameter values for each detector can be given separated by commas, whitespace or tab characters.
2758  //
2759  //The `[CalibClass]`, if given, must correspond to a KVCalibrator plugin name. The list of plugin names and the corresponding
2760  //classes can be retrieved with
2761  //
2762  //~~~~~~~~~~~
2763  //KVBase::GetListOfPlugins("KVCalibrator")
2764  //KVBase::GetListOfPluginURIs("KVCalibrator")
2765  //~~~~~~~~~~~
2766  //
2767  //KVCalibrator objects are added to detectors as required by the contents of calibration files.
2768  //If any detector has an existing calibrator of type `[CalibType]` which is not of the given class
2769  //it will be replaced with a new calibrator corresponding to the plugin.
2770  //
2771  //The `[CalibOptions]` is optional: list in `[CalibOptions]` will be used
2772  //to complete set-up of any new calibrator objects by calling the KVCalibrator::SetOptions()
2773  //method.
2774  //
2775  //`[CalibOptions]` should hold a comma-separated list of `parameter=value` pairs which will be used
2776  //to fill a KVNameValueList for the method call. See the KVCalibrator::SetOptions() method.
2777  //
2778  //`[ZRange]` is an option if several calibrations need to be used to provide the same signal
2779  //for certain detectors depending on the atomic number Z of the particle detected.
2780 
2781 
2782  TString fullpath = "";
2783  if (!SearchKVFile(filename, fullpath, fDataSet)) {
2784  Info("ReadCalibFile", "%s does not exist or not found", filename);
2785  return;
2786  }
2787 
2788  Info("ReadCalibFile", "file : %s found", fullpath.Data());
2789  TEnv env;
2790  env.ReadFile(fullpath, kEnvAll);
2791 
2792  // read options from file
2793  KVNameValueList options;
2794  KVString opt_list = "RunList SignalIn SignalOut CalibType CalibClass CalibOptions ZRange";
2795  opt_list.Begin(" ");
2796  while (!opt_list.End()) {
2797  KVString opt = opt_list.Next();
2798  KVString opt_val = env.GetValue(opt, "");
2799  opt_val.Remove(TString::kBoth, ' ');
2800  options.SetValue(opt, opt_val.Data());
2801  }
2802  // check for stupid spellnig mitskaes
2803  if (TString(env.GetValue("Runlist", "")) != "") {
2804  Warning("ReadCalibFile", "Calibration has 'Runlist' parameter (ignored): %s, did you mean 'RunList'?", env.GetValue("Runlist", ""));
2805  }
2806 
2807  if (options.GetTStringValue("SignalIn") == "") {
2808  Error("ReadCalibFile", "No input signal defined : SignalIn");
2809  return;
2810  }
2811  if (options.GetTStringValue("SignalOut") == "") {
2812  Error("ReadCalibFile", "No output signal defined : SignalOut");
2813  return;
2814  }
2815  if (options.GetTStringValue("CalibType") == "") {
2816  Error("ReadCalibFile", "No calibration type defined : CalibType");
2817  return;
2818  }
2819  Bool_t check_class(options.GetTStringValue("CalibClass") != "");
2820  TString calibrator_class;
2821  if (check_class) {
2822  TPluginHandler* ph = LoadPlugin("KVCalibrator", options.GetStringValue("CalibClass"));
2823  if (ph) calibrator_class = ph->GetClass();
2824  else {
2825  Error("ReadCalibFile", "No calibrator plugin of type %s", options.GetStringValue("CalibClass"));
2826  return;
2827  }
2828  }
2829 
2830  KVString clop;
2831  if (options.HasParameter("CalibOptions")) clop = options.GetStringValue("CalibOptions");
2832 
2833  KVString zrange;
2834  if (options.HasParameter("ZRange")) zrange = options.GetStringValue("ZRange");
2835 
2836  KVNumberList run_list = db->GetRunList();
2837  if (options.GetTStringValue("RunList") != "") {
2838  run_list.Set(options.GetTStringValue("RunList"));
2839  Info("ReadCalibFile", "Calibration used for runs %s", run_list.AsString());
2840  }
2841  else {
2842  Info("ReadCalibFile", "Calibration used for all runs in database");
2843  }
2844 
2845  TIter next(env.GetTable());
2846  TEnvRec* rec = 0;
2847  KVDBParameterSet* par = 0;
2848 
2849  while ((rec = (TEnvRec*)next())) {
2850 
2851  TString sname(rec->GetName());
2852  KVDetector* det = GetDetector(sname);
2853  if (!det) continue;
2854 
2855  KVString lval(rec->GetValue());
2856  // valid delimiters for numerical parameters are: space - comma - tab
2857  // here we add ':' because if there is space in the TEnv file between the key name and the ':'
2858  // i.e. if we have
2859  // SI2-11 : 0.1,0.3
2860  // then the ':' becomes part of the value: TEnv::GetValue("SI2-11","") will give ": 0.1,0.3"
2861  // instead of
2862  // SI2-11: 0.1,0.3
2863  // => TEnv::GetValue("SI2-11","") gives "0.1,0.3"
2864  TString param_delimiters = ":, \t";
2865  par = new KVDBParameterSet(sname.Data(), options.GetStringValue("CalibType"), lval.GetNValues(param_delimiters));
2866  par->SetParameter("SignalIn", options.GetStringValue("SignalIn"));
2867  par->SetParameter("SignalOut", options.GetStringValue("SignalOut"));
2868  // put infos on required calibrator class into database so that it can be replaced
2869  // as needed in SetCalibratorParameters
2870  par->SetParameter("CalibClass", options.GetStringValue("CalibClass"));
2871  if (clop != "") par->SetParameter("CalibOptions", clop);
2872  if (zrange != "") par->SetParameter("ZRange", zrange);
2873  Int_t np = 0;
2874  lval.Begin(param_delimiters);
2875  while (!lval.End()) {
2876  par->SetParameter(np++, lval.Next().Atof());
2877  }
2878  calib_table->AddRecord(par);
2879  db->LinkRecordToRunRange(par, run_list);
2880 
2881  }
2882 }
2883 
2884 
2885 #ifdef WITH_MFM
2886 
2892 
2894 {
2895  // Update array according to last event read using the KVMFMDataFileReader object
2896  // (it is assumed that KVRawDataReader::GetNextEvent() was called before calling this method)
2897  //
2898  // Return kTRUE if raw data was treated
2899 
2900  if (mfmreader.IsFrameReadMerge()) {
2901  return handle_raw_data_event_mfmmergeframe(mfmreader.GetMergeManager());
2902  }
2903  else {
2904  return handle_raw_data_event_mfmframe(mfmreader.GetFrameRead());
2905  }
2906  return kFALSE;
2907 }
2908 
2909 
2910 
2914 
2915 Bool_t KVMultiDetArray::handle_raw_data_event_mfmmergeframe(const MFMMergeFrameManager& mergeframe)
2916 {
2917  // Method used to handle merged MFM frames
2918  // We call handle_raw_data_event_mfmframe() for each frame contained in the merge
2919 
2920  Bool_t ok = false;
2921  while (mergeframe.ReadNextFrame()) {
2922  Bool_t me = handle_raw_data_event_mfmframe(mergeframe.GetFrameRead());
2923  ok = (ok || me);
2924  }
2925  return ok;
2926 }
2927 
2928 
2929 
2940 
2942 {
2943  // Method used to treat raw data in MFM format read by KVMFMDataFileReader
2944  //
2945  // Here we dispatch two types of frame - MFMEbyedatFrame & MFMMesytecMDPPFrame -
2946  // to specific methods - handle_raw_data_event_mfmframe_ebyedat() and
2947  // handle_raw_data_event_mfmframe_mesytec_mdpp()
2948  // which need to be implemented in child classes for specific arrays which
2949  // use these data formats.
2950  //
2951  // Return kTRUE if raw data was treated
2952 #ifdef WITH_MESYTEC
2953  if (mfmframe.GetFrameType() == MFM_MESYTEC_FRAME_TYPE)
2954  return handle_raw_data_event_mfmframe_mesytec_mdpp((const MFMMesytecMDPPFrame&)mfmframe);
2955 #endif
2956  if (mfmframe.GetFrameType() == MFM_EBY_EN_FRAME_TYPE
2957  || mfmframe.GetFrameType() == MFM_EBY_TS_FRAME_TYPE
2958  || mfmframe.GetFrameType() == MFM_EBY_EN_TS_FRAME_TYPE)
2959  return handle_raw_data_event_mfmframe_ebyedat((const MFMEbyedatFrame&)mfmframe);
2960 
2961  return kFALSE;
2962 }
2963 
2964 
2965 
2968 
2970 {
2971  // Read a raw data event from a EBYEDAT MFM Frame.
2972 
2973  AbstractMethod("handle_raw_data_event_mfmframe_ebyedat");
2974  return kFALSE;
2975 }
2976 
2977 
2978 #ifdef WITH_MESYTEC
2979 
2982 
2984 {
2985  // Read a raw data event from a Mesytec MFM Frame.
2986 
2987  AbstractMethod("handle_raw_data_event_mfmframe_mesytec_mdpp");
2988  return kFALSE;
2989 }
2990 
2991 #endif
2992 #endif
2993 
2994 #ifdef WITH_PROTOBUF
2995 
2997 
2999 {
3000  AbstractMethod("handle_raw_data_event_protobuf");
3001  return kFALSE;
3002 }
3003 
3004 #endif
3005 
3006 
3009 
3011 {
3012  // For each IDtelescope in array, calculate an identification grid
3013 
3014  TIter nxtid(GetListOfIDTelescopes());
3015  KVIDTelescope* idt;
3016  while ((idt = (KVIDTelescope*) nxtid())) {
3017  idt->CalculateDeltaE_EGrid("1-92", 0, 20);
3018  }
3019 }
3020 
3021 
3022 
3027 
3029 {
3030  // Sets status of detectors (KVDetector::IsPresent() and KVDetector::IsWorking()) for a given run of a dataset.
3031  //
3032  // If 'myname' is given, we look in database table "myname.OoODets"
3033 
3034  KVRList* absdet = (myname != "" ? kvrun->GetLinks(Form("%s.Absent Detectors", myname.Data())) : kvrun->GetLinks("Absent Detectors"));
3035  KVRList* ooodet = (myname != "" ? kvrun->GetLinks(Form("%s.OoO Detectors", myname.Data())) : kvrun->GetLinks("OoO Detectors"));
3036 
3037  TIter next(GetDetectors());
3038  KVDetector* det;
3039 
3040  Int_t ndet_absent = 0;
3041  Int_t ndet_ooo = 0;
3042  TString absent_dets, ooo_dets;
3043 
3044  while ((det = (KVDetector*)next())) {
3045  //Test de la presence ou non du detecteur
3046  if (!absdet) {
3047  det->SetPresent();
3048  }
3049  else {
3050  if (absdet->FindObject(det->GetName())) {
3051  det->SetPresent(kFALSE);
3052  if (ndet_absent) absent_dets += ",";
3053  absent_dets += det->GetName();
3054  ndet_absent += 1;
3055  }
3056  else {
3057  det->SetPresent();
3058  }
3059  }
3060  if (det->IsPresent()) {
3061  //Test du bon fonctionnement ou non du detecteur
3062  if (!ooodet) {
3063  det->SetDetecting();
3064  }
3065  else {
3066  if (ooodet->FindObject(det->GetName())) {
3067  det->SetDetecting(kFALSE);
3068  if (ndet_ooo) ooo_dets += ",";
3069  ooo_dets += det->GetName();
3070  ndet_ooo += 1;
3071  }
3072  else {
3073  det->SetDetecting();
3074  }
3075  }
3076  }
3077  }
3078 
3079  if (ndet_absent) Info("CheckStatusOfDetectors", "%d detectors absent during run : %s", ndet_absent, absent_dets.Data());
3080  else Info("CheckStatusOfDetectors", "All detectors present during run");
3081  if (ndet_ooo) Info("CheckStatusOfDetectors", "%d detectors malfunctioned during run : %s", ndet_ooo, ooo_dets.Data());
3082  else Info("CheckStatusOfDetectors", "All detectors functioning during run");
3083 }
3084 
3085 
3086 
3088 
3090 {
3091  CheckStatusOfDetectors(dbr, myname);
3092 }
3093 
3094 
3095 
3110 
3112 {
3113  // Read a file containing runlists for each temporarily non-functioning detector.
3114  //
3115  // The file should be in TEnv format like so:
3116  //
3117  //~~~~
3118  // DET_1: 100-122,541-1938
3119  // DET_2,DET_3: 91-765
3120  //~~~~
3121  //
3122  // i.e. more than one detector can be associated with the same runs (comma-separated list of
3123  // detector names) and the list of runs are given using KVNumberList syntax.
3124  //
3125  // The data is added to the database in a table '[name].OoO Detectors' with the name of this array.
3126 
3127  TString fullpath;
3128  if (!db->FindCalibFile("OoODet", fullpath, GetName())) return;
3129 
3130  Info("ReadOoODetectors()", "Reading lists of out-of-order detectors...");
3131  auto fOoODet = db->AddTable(Form("%s.OoO Detectors", GetName()), "Name of out of order detectors");
3132 
3133  KVDBRecord* dbrec = 0;
3134  TEnv env;
3135  TEnvRec* rec = 0;
3136  env.ReadFile(fullpath.Data(), kEnvAll);
3137  TIter it(env.GetTable());
3138 
3139  while ((rec = (TEnvRec*)it.Next())) {
3140  KVString srec(rec->GetName());
3141  KVNumberList nl(rec->GetValue());
3142  if (srec.Contains(",")) {
3143  srec.Begin(",");
3144  while (!srec.End()) {
3145  dbrec = new KVDBRecord(srec.Next(kTRUE), "OoO Detector");
3146  dbrec->AddKey("Runs", "List of Runs");
3147  fOoODet->AddRecord(dbrec);
3148  db->LinkRecordToRunRange(dbrec, nl);
3149  }
3150  }
3151  else {
3152  dbrec = new KVDBRecord(rec->GetName(), "OoO Detector");
3153  dbrec->AddKey("Runs", "List of Runs");
3154  fOoODet->AddRecord(dbrec);
3155  db->LinkRecordToRunRange(dbrec, nl);
3156  }
3157  }
3158 }
3159 
3160 
3175 
3177 {
3178  // Read a file containing runlists for each temporarily absent/dismounted detector.
3179  //
3180  // The file should be in TEnv format like so:
3181  //
3182  //~~~~
3183  // DET_1: 100-122,541-1938
3184  // DET_2,DET_3: 91-765
3185  //~~~~
3186  //
3187  // i.e. more than one detector can be associated with the same runs (comma-separated list of
3188  // detector names) and the list of runs are given using KVNumberList syntax.
3189  //
3190  // The data is added to the database in a table '[name].Absent Detectors' with the name of this array.
3191 
3192  TString fullpath;
3193  if (!db->FindCalibFile("AbsentDet", fullpath, GetName())) return;
3194 
3195  Info("ReadAbsentDetectors()", "Reading lists of absent/dismounted detectors... file=[%s]", fullpath.Data());
3196  auto fAbsDet = db->AddTable(Form("%s.Absent Detectors", GetName()), "Name of out of order detectors");
3197 
3198  KVDBRecord* dbrec = 0;
3199  TEnv env;
3200  TEnvRec* rec = 0;
3201  env.ReadFile(fullpath.Data(), kEnvAll);
3202  TIter it(env.GetTable());
3203 
3204  while ((rec = (TEnvRec*)it.Next())) {
3205  KVString srec(rec->GetName());
3206  KVNumberList nl(rec->GetValue());
3207  if (srec.Contains(",")) {
3208  srec.Begin(",");
3209  while (!srec.End()) {
3210  dbrec = new KVDBRecord(srec.Next(kTRUE), "Absent Detector");
3211  dbrec->AddKey("Runs", "List of Runs");
3212  fAbsDet->AddRecord(dbrec);
3213  db->LinkRecordToRunRange(dbrec, nl);
3214  }
3215  }
3216  else {
3217  dbrec = new KVDBRecord(rec->GetName(), "Absent Detector");
3218  dbrec->AddKey("Runs", "List of Runs");
3219  fAbsDet->AddRecord(dbrec);
3220  db->LinkRecordToRunRange(dbrec, nl);
3221  }
3222  }
3223 }
3224 
3225 
int Int_t
unsigned int UInt_t
ROOT::R::TRInterface & r
#define SafeDelete(p)
#define d(i)
#define f(i)
#define e(i)
bool Bool_t
unsigned short UShort_t
char Char_t
constexpr Bool_t kFALSE
double Double_t
constexpr Bool_t kTRUE
const char Option_t
kEnvUser
kEnvAll
#define N
winID h TVirtualViewer3D TVirtualGLPainter p
Option_t Option_t option
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 filename
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 np
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 Atom_t Time_t UChar_t Atom_t typelist
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
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 Atom_t Time_t type
char name[80]
R__EXTERN TGeoManager * gGeoManager
char * Form(const char *fmt,...)
#define gPad
void SetLabel(const Char_t *lab)
Definition: KVBase.h:194
virtual const Char_t * GetType() const
Definition: KVBase.h:176
void Error(const char *method, const char *msgfmt,...) const override
Definition: KVBase.cpp:1683
static const Char_t * GetWORKDIRFilePath(const Char_t *namefile="")
Definition: KVBase.cpp:121
const Char_t * GetLabel() const
Definition: KVBase.h:198
static bool is_gnuinstall()
Definition: KVBase.h:275
void Print(Option_t *option="") const override
Definition: KVBase.cpp:413
virtual void SetType(const Char_t *str)
Definition: KVBase.h:172
static Bool_t SearchKVFile(const Char_t *name, TString &fullpath, const Char_t *kvsubdir="")
Definition: KVBase.cpp:541
static TPluginHandler * LoadPlugin(const Char_t *base, const Char_t *uri="0")
Definition: KVBase.cpp:796
void Warning(const char *method, const char *msgfmt,...) const override
Definition: KVBase.cpp:1670
Base class for all detector calibrations.
Definition: KVCalibrator.h:99
Int_t GetNumberParams() const
Definition: KVCalibrator.h:153
virtual void SetOptions(const KVNameValueList &)
void SetParameter(int i, Double_t par_val) const
Definition: KVCalibrator.h:158
void SetStatus(Bool_t ready)
Definition: KVCalibrator.h:166
void SetOutputSignalType(const TString &type)
Definition: KVCalibrator.h:220
void SetInputSignalType(const TString &type)
Definition: KVCalibrator.h:216
static KVCalibrator * MakeCalibrator(const Char_t *type)
To store calibration parameters in a database ,.
TString GetStringParameter(const TString &name) const
Double_t GetParameter(UShort_t i=0) const
Bool_t HasParameter(const TString &name) const
void Print(Option_t *option="") const override
Int_t GetParamNumber() const
void SetParameter(UShort_t i, Double_t val)
const KVNameValueList & GetParameters() const
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
Description of an experimental run in database ,,.
Definition: KVDBRun.h:41
Table in a database.
Definition: KVDBTable.h:34
virtual Bool_t AddRecord(KVDBRecord *add)
Definition: KVDBTable.cpp:74
static Bool_t IsRunningBatchAnalysis()
virtual Bool_t AddTable(KVDBTable *table)
Definition: KVDataBase.cpp:84
KVDataSet * GetDataSet(Int_t) const
Return pointer to DataSet using index in list of all datasets, index>=0.
Manage an experimental dataset corresponding to a given experiment or campaign.
Definition: KVDataSet.h:146
TString GetFullPathToDataSetFile(const Char_t *filename) const
Definition: KVDataSet.cpp:2067
const Char_t * GetDataSetDir() const
Definition: KVDataSet.cpp:762
Bool_t HasCalibIdentInfos() const
Definition: KVDataSet.h:380
void cd() const
Definition: KVDataSet.cpp:784
Bool_t DataBaseUpdateInProgress() const
Definition: KVDataSet.h:186
ValType GetDataSetEnv(const Char_t *type, const ValType &defval={}) const
Definition: KVDataSet.h:269
static Bool_t FindDataSetFile(const TString &dataset, const Char_t *filename)
Definition: KVDataSet.cpp:2101
List of hit groups in a multidetector array.
void AddGroup(KVGroup *grp)
Base class for output signal data produced by a detector.
virtual void SetValue(Double_t x)
virtual Double_t GetValue(const KVNameValueList &params="") const
TString GetFullName() const
void SetFired(Bool_t yes=true)
Base class for detector geometry description, interface to energy-loss calculations.
Definition: KVDetector.h:160
void SetThickness(Double_t thick) override
KVGroup * GetGroup() const
void AddDetectorSignal(KVDetectorSignal *ds)
Definition: KVDetector.h:278
void Print(Option_t *option="") const override
Definition: KVDetector.cpp:357
KVMaterial * GetAbsorber(Int_t i) const
Returns pointer to the i-th absorber in the detector (i=0 first absorber, i=1 second,...
Definition: KVDetector.cpp:624
virtual KVDetectorSignal * GetDetectorSignal(const KVString &type) const
Definition: KVDetector.h:546
Bool_t IsCalibrated() const
Definition: KVDetector.h:407
void SetDetecting(Bool_t yes=kTRUE)
Definition: KVDetector.h:674
virtual void RemoveCalibrators()
Definition: KVDetector.cpp:702
KVGeoDetectorNode * GetNode()
Definition: KVDetector.h:346
virtual Bool_t IsPresent() const
Definition: KVDetector.h:660
Bool_t AddCalibrator(KVCalibrator *cal, const KVNameValueList &opts="")
Definition: KVDetector.cpp:440
void SetPresent(Bool_t yes=kTRUE)
Definition: KVDetector.h:665
Base class to describe database of an experiment ,,.
Definition: KVExpDB.h:61
void SetDBType(const TString &s)
Definition: KVExpDB.h:118
TString GetCalibFileName(const Char_t *type) const
Definition: KVExpDB.h:188
const KVNumberList & GetRunList() const
Definition: KVExpDB.h:156
virtual void LinkRecordToRunRange(KVDBRecord *rec, UInt_t first_run, UInt_t last_run)
Definition: KVExpDB.cpp:285
Bool_t FindCalibFile(const Char_t *type, TString &fullpath, const TString &array_name="") const
Definition: KVExpDB.cpp:667
Handle reading columns of numeric data in text files.
Definition: KVFileReader.h:121
KVString GetCurrentLine()
Definition: KVFileReader.h:320
ReadStatus ReadLine(const KVString &pattern="")
Definition: KVFileReader.h:243
Bool_t IsOK()
Definition: KVFileReader.h:231
Bool_t OpenFileToRead(const KVString &filename)
Definition: KVFileReader.h:210
Reads GANIL acquisition files (EBYEDAT)
Path taken by particles through multidetector geometry.
KVGeoDetectorNode * GetNextNode() const
KVSeqCollection * AccessIDTelescopeList()
void AddToNodes()
Add reference to this trajectory to all nodes on it.
void IterateFrom(const KVGeoDetectorNode *node0=nullptr) const
KVGeoDetectorNode * GetNodeInFront(const KVGeoDetectorNode *n) const
Bool_t ContainsPath(const KVGeoDNTrajectory *other) const
KVGeoDetectorNode * GetNodeAt(Int_t i) const
Information on relative positions of detectors & particle trajectories.
KVDetector * GetDetector() const
Base class for propagation of particles through array geometry.
void DrawTracks(KVNumberList *=nullptr)
Base class describing elements of array geometry.
virtual Bool_t Fired(Option_t *opt="any") const
void SetOwnsDetectors(Bool_t yes=kTRUE)
virtual KVDetector * GetDetector(const Char_t *name) const
Return detector in this structure with given name.
const KVSeqCollection * GetDetectors() const
const KVSeqCollection * GetStructures() const
virtual void Add(KVBase *)
KVSeqCollection * GetStructureTypeList(const Char_t *type) const
void ClearStructures(const Char_t *type="")
const KVSeqCollection * GetParents() const
Base class for particle reconstruction in one group of a detector array.
static KVGroupReconstructor * Factory(const TString &plugin="", const KVGroup *g=nullptr)
Group of detectors which can be treated independently of all others in array.
Definition: KVGroup.h:20
void SetNumber(UInt_t num) override
Definition: KVGroup.h:38
void AddTrajectory(KVGeoDNTrajectory *t)
Definition: KVGroup.h:91
void Reset(Option_t *opt="")
Definition: KVGroup.h:45
Extended version of ROOT THashList.
Definition: KVHashList.h:29
Base class for particle identification in a 2D map.
Definition: KVIDGraph.h:32
Handles a stock of identification grids to be used by one or more identification telescopes.
virtual void LoadGridsForRun(UInt_t)
KVSeqCollection * GetGrids()
Bool_t ReadAsciiFile(const Char_t *filename)
const TList * GetLastReadGrids() const
virtual bool IsSQLROOT() const
Base class for all detectors or associations of detectors in array which can identify charged particl...
Definition: KVIDTelescope.h:85
virtual Double_t GetIDMapY(Option_t *opt="")
virtual Bool_t IsReadyForID()
KVIDGrid * CalculateDeltaE_EGrid(const KVNameValueList &AperZ, Int_t npoints=30, Double_t xfactor=1.)
void SetGroup(KVGroup *kvg)
virtual Double_t GetIDMapX(Option_t *opt="")
static KVIDTelescope * MakeIDTelescope(const Char_t *name)
virtual Bool_t SetIdentificationParameters(const KVMultiDetArray *)
void SetIDGrid(KVIDGraph *)
virtual void AddDetector(KVDetector *d)
virtual void Initialize(void)
Extended TList class which owns its objects by default.
Definition: KVList.h:22
Read MFM format acquisition data.
Description of physical materials used to construct detectors & targets; interface to range tables.
Definition: KVMaterial.h:90
static KVIonRangeTable * GetRangeTable()
Definition: KVMaterial.cpp:158
Base class for describing the geometry of a detector array.
const KVSeqCollection * GetFiredSignals() const
KVNumberList fAcceptECodes
list of acceptable calibration codes for reconstructed nuclei
bool try_a_singleIDtelescope(TString uri, KVDetector *d, TCollection *l)
KVUniqueNameList fExtraRawDataSignals
any signals read from raw data not associated with a detector
KVSeqCollection * GetListOfIDTelescopes() const
void FillListOfIDTelescopes(KVIDGraph *gr) const
virtual Bool_t handle_raw_data_event_mfmframe_ebyedat(const MFMEbyedatFrame &)
Read a raw data event from a EBYEDAT MFM Frame.
KVNumberList fAcceptIDCodes
list of acceptable identification codes for reconstructed nuclei
virtual void GetDetectorEvent(KVDetectorEvent *detev, const TSeqCollection *fired_params=0)
virtual Bool_t handle_raw_data_event_protobuf(KVProtobufDataReader &)
void MakeHistogramsForAllIDTelescopes(KVSeqCollection *list, Int_t dimension=100)
TList * GetStatusOfIDTelescopes()
int try_all_singleID_telescopes(KVDetector *d, TCollection *l)
static Bool_t fCloseGeometryNow
void Draw(Option_t *option="") override
virtual void DeduceIdentificationTelescopesFromGeometry()
static KVMultiDetArray * MakeMultiDetector(const Char_t *dataset_name, Int_t run=-1, TString classname="KVMultiDetArray", KVExpDB *db=nullptr)
std::unique_ptr< KVFileReader > GetKVFileReader(KVExpDB *db, const Char_t *meth, const Char_t *keyw)
virtual Double_t GetTargetEnergyLossCorrection(KVReconstructedNucleus *)
virtual Bool_t handle_raw_data_event_mfmmergeframe(const MFMMergeFrameManager &)
virtual Bool_t handle_raw_data_event_mfmfile(MFMBufferReader &)
void ReadOoODetectors(KVExpDB *db)
Bool_t fHandledRawData
set to true if multidetector handles data in last call to HandleRawData
KVSeqCollection * GetIDTelescopesWithType(const Char_t *type)
TList * fCalibStatusDets
used by GetStatusIDTelescopes
KVDetectorEvent * fHitGroups
list of hit groups in simulation
void ReadAbsentDetectors(KVExpDB *db)
KVTarget * GetTarget()
void RecursiveTrajectoryClustering(KVGeoDetectorNode *N, KVUniqueNameList &tried_trajectories, KVUniqueNameList &multitraj_nodes, KVUniqueNameList &detectors_of_group)
KVSeqCollection * fIDTelescopes
deltaE-E telescopes in groups
UInt_t fCurrentRun
Number of the current run used to call SetParameters.
static Bool_t fMakeMultiDetectorSetParameters
void Clear(Option_t *opt="") override
void CalculateIdentificationGrids()
For each IDtelescope in array, calculate an identification grid.
void prepare_to_handle_new_raw_data()
reset acquisition parameters etc. before reading new raw data event
int try_all_doubleID_telescopes(KVDetector *de, KVDetector *e, TCollection *l)
Bool_t IsSimMode() const
void FillHistogramsForAllIDTelescopes(KVSeqCollection *list)
Fill TH2F histograms for all IDTelescopes of the array.
virtual void set_detector_thicknesses(const TString &)
Bool_t ReadGridsFromAsciiFile(const Char_t *) const
virtual void InitialiseRawDataReading(KVRawDataReader *)
virtual void SetExpectedDetectorSignalNames()
virtual void SetParameters(UInt_t n, Bool_t physics_parameters_only=kFALSE)
virtual void SetDetectorParametersForRun(KVDBRun *, const TString &="")
virtual void copy_fired_parameters_to_recon_param_list()
TString GetDataSet() const
KVGroup * GetGroup(const Char_t *name) const
static TGeoHMatrix * GetVolumePositioningMatrix(Double_t distance, Double_t theta, Double_t phi, TGeoTranslation *postTrans=nullptr)
KVGeoNavigator * GetNavigator() const
KVUpDater * fUpDater
used to set parameters for multidetector
void SetGeometry(TGeoManager *)
virtual void MakeCalibrationTables(KVExpDB *)
virtual ~KVMultiDetArray()
destroy (delete) the MDA and all the associated structure, detectors etc.
virtual void SetRunIdentificationParameters(UShort_t n)
Bool_t HandleRawDataBuffer(MFMBufferReader &)
TString GetFileName(KVExpDB *, const Char_t *meth, const Char_t *keyw)
void CheckStatusOfDetectors(KVDBRun *, const TString &="")
virtual void PerformClosedROOTGeometryOperations()
UInt_t GetCurrentRunNumber() const
void SetDetectorTransparency(Char_t)
Modify the transparency of detector volumes in OpenGL view.
virtual KVGroupReconstructor * GetReconstructorForGroup(const KVGroup *) const
virtual void SetIdentifications()
static Bool_t fMakeMultiDetectorPhysicsParametersOnly
TList * fStatusIDTelescopes
used by GetStatusIDTelescopes
Bool_t fSimMode
=kTRUE in "simulation mode" (use for calculating response to simulated events)
virtual void Build(Int_t run=-1)
virtual void FillDetectorList(KVReconstructedNucleus *rnuc, KVHashList *DetList, const KVString &DetNames)
virtual Bool_t handle_raw_data_event_ebyedat(KVGANILDataReader &)
General method for reading raw data in old GANIL ebyedat format.
const TSeqCollection * GetTrajectories() const
static Bool_t fBuildTarget
KVTarget * fTarget
target used in experiment
virtual void set_up_telescope(KVDetector *de, KVDetector *e, KVIDTelescope *idt, TCollection *l)
Set up detectors in de-e identification telescope and add to fIDTelescopes and to l.
TList * GetCalibrationStatusOfDetectors()
bool try_upper_and_lower_singleIDtelescope(TString uri, KVDetector *d, TCollection *l)
virtual Bool_t handle_raw_data_event_mfmframe(const MFMCommonFrame &)
TGeoManager * GetGeometry() const
Return pointer to the (ROOT) geometry of the array.
KVRawDataReader * fRawDataReader
last raw data reader object used in call to HandleRawData
virtual void SetRawDataFromReconEvent(KVNameValueList &)
void add_and_set_detector_signal(KVDetector *det, KVString detname, Double_t sig_data, KVString sig_type)
static Int_t fMakeMultiDetectorRunNumber
virtual Bool_t HandleRawDataEvent(KVRawDataReader *)
virtual void InitializeIDTelescopes()
KVUniqueNameList fFiredDetectors
list of fired detectors after reading raw data event
virtual void RenumberGroups()
virtual void SetRunCalibrationParameters(UShort_t n)
void PrintStatusOfIDTelescopes()
bool try_upper_and_lower_doubleIDtelescope(TString uri, KVDetector *de, KVDetector *e, TCollection *l)
void PrintCalibStatusOfDetectors()
void ReadCalibrationFiles(KVExpDB *db)
virtual void SetCalibratorParameters(KVDBRun *, const TString &="")
Int_t fFilterType
type of filtering (used by DetectEvent)
virtual void SetPartSeedCond(const Char_t *cond)
std::unique_ptr< KVGeoNavigator > fNavigator
for propagating particles through array geometry
void SetGridsInTelescopes(UInt_t run)
KVIDTelescope * GetIDTelescope(const Char_t *name) const
Return pointer to DeltaE-E ID Telescope with "name".
void CalculateReconstructionTrajectories()
KVUniqueNameList fTrajectories
list of all possible trajectories through detectors of array
virtual void SetIDCodeForIDTelescope(KVIDTelescope *) const
virtual void SetReconParametersInEvent(KVReconstructedEvent *) const
Copy any parameters in fReconParameters in to the reconstructed event parameter list.
KVNameValueList fReconParameters
general purpose list of parameters for storing information on data reconstruction
TString fDataSet
name of associated dataset, used with MakeMultiDetector()
void AssociateTrajectoriesAndNodes()
void ReadCalibFile(const Char_t *filename, KVExpDB *db, KVDBTable *calib_table)
virtual void AcceptParticleForAnalysis(KVReconstructedNucleus *) const
Int_t GetIDTelescopes(KVGeoDetectorNode *, KVDetector *, TCollection *list)
void SetNavigator(KVGeoNavigator *geo)
KVMultiDetArray()
Default constructor.
KVUpDater * GetUpDater()
KVUniqueNameList GetIDTelescopeTypes()
void DeduceGroupsFromTrajectories()
virtual void set_up_single_stage_telescope(KVDetector *det, KVIDTelescope *idt, TCollection *l)
Set up detector in single-stage identification telescope and add to fIDTelescopes and to l.
bool try_a_doubleIDtelescope(TString uri, KVDetector *de, KVDetector *e, TCollection *l)
virtual Bool_t handle_raw_data_event_mfmframe_mesytec_mdpp(const MFMMesytecMDPPFrame &)
Read a raw data event from a Mesytec MFM Frame.
KVUnownedList fFiredSignals
list of fired signals after reading raw data event
Handles lists of named parameters with different types, a list of KVNamedParameter objects.
void SetValue(const Char_t *name, value_type value)
void Clear(Option_t *opt="") override
const Char_t * GetStringValue(const Char_t *name) const
Bool_t HasParameter(const Char_t *name) const
TString GetTStringValue(const Char_t *name) const
A generic named parameter storing values of different types.
Strings used to represent a set of ranges of values.
Definition: KVNumberList.h:85
Bool_t Contains(Int_t val) const
returns kTRUE if the value 'val' is contained in the ranges defined by the number list
const Char_t * AsString(Int_t maxchars=0) const
void Set(const TString &l)
Definition: KVNumberList.h:135
Bool_t IsEmpty() const
Definition: KVNumberList.h:175
void SetIsOK(Bool_t flag=kTRUE)
Definition: KVParticle.cpp:371
Double_t GetEnergy() const
Definition: KVParticle.h:624
Read Google Protobuf DAQ files.
Wrapper for TRefArray adding some functionality.
Definition: KVRList.h:37
KVBase * FindObject(const Char_t *name, const Char_t *type) const
Definition: KVRList.cpp:106
Propagate particles through array geometry calculating energy losses.
Abstract base class for reading raw (DAQ) data.
virtual TString GetDataFormat() const =0
Event containing KVReconstructedNucleus nuclei reconstructed from hits in detectors.
Nuclei reconstructed from data measured by a detector array .
virtual Int_t GetECode() const
virtual Int_t GetIDCode() const
Combine ROOT file containing objects with SQLite database with info on the objects.
Definition: KVSQLROOTFile.h:76
ID grid manager using KVSQLROOTFile backend.
KaliVeda extensions to ROOT collection classes.
virtual TObject * FindObjectByLabel(const Char_t *) const
T * get_object(const TString &name) const
TObject * Remove(TObject *obj) override
Remove object from list.
void Add(TObject *obj) override
TObject * FindObject(const char *name) const override
KVSeqCollection * GetSubListWithType(const Char_t *retvalue) const
void Clear(Option_t *option="") override
virtual void SetCleanup(Bool_t enable=kTRUE)
void SetOwner(Bool_t enable=kTRUE) override
void Delete(Option_t *option="") override
KVSeqCollection * GetSubListWithLabel(const Char_t *retvalue) const
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
Int_t GetNValues(TString delim) const
Definition: KVString.cpp:886
Double_t GetParticleEIncFromERes(KVNucleus *, TVector3 *norm=0) override
Definition: KVTarget.cpp:958
void Clear(Option_t *opt="") override
Definition: KVTarget.cpp:748
Optimised list in which named objects can only be placed once.
void Add(TObject *obj) override
Abstract class implementing necessary methods for setting multidetector parameters for each run of th...
Definition: KVUpDater.h:25
virtual void SetParameters(UInt_t, Bool_t physics_parameters_only=kFALSE)
Definition: KVUpDater.cpp:91
static KVUpDater * MakeUpDater(const Char_t *uri, KVMultiDetArray *)
Definition: KVUpDater.cpp:59
virtual void SetCalibrationParameters(UInt_t)
Set calibration parameters for this run.
Definition: KVUpDater.cpp:198
virtual void SetIdentificationParameters(UInt_t)
Definition: KVUpDater.cpp:157
void SetName(const char *name)
const char * GetName() const override
virtual Int_t GetEntries() const
virtual void SetOwner(Bool_t enable=kTRUE)
Bool_t Contains(const char *name) const
THashList * GetTable() const
virtual const char * GetValue(const char *name, const char *dflt) const
virtual Int_t ReadFile(const char *fname, EEnvLevel level)
Bool_t Defined(const char *name) const
void SetStyle(Short_t st)
void SetCurrentCamera(ECameraType camera)
void SetSmoothPoints(Bool_t s)
void SetSmoothLines(Bool_t s)
TGeoVolume * GetTopVolume() const
void SetAngles(Double_t phi, Double_t theta, Double_t psi)
void SetDz(Double_t dz) override
void Draw(Option_t *option="") override
void SetTransparency(Char_t transparency=0)
virtual Int_t Fill(const char *name, Double_t w)
TObject * Next()
void Reset()
TObject * FindObject(const char *name) const override
void Add(TObject *obj) override
void Delete(Option_t *option="") override
TObject * At(Int_t idx) const override
virtual void SetTitle(const char *title="")
const char * GetName() const override
const char * GetTitle() const override
virtual void SetName(const char *name)
void AbstractMethod(const char *method) const
virtual const char * GetName() const
void SetBit(UInt_t f)
R__ALWAYS_INLINE Bool_t TestBit(UInt_t f) const
virtual void Info(const char *method, const char *msgfmt,...) const
const char * GetClass() const
Longptr_t ExecPlugin(int nargs)
Int_t GetEntries() const override
Double_t Atof() const
const char * Data() const
void ToUpper()
TObjArray * Tokenize(const TString &delim) const
Bool_t BeginsWith(const char *s, ECaseCompare cmp=kExact) const
TString & Prepend(char c, Ssiz_t rep=1)
void Form(const char *fmt,...)
TString & Remove(EStripType s, char c)
Bool_t Contains(const char *pat, ECaseCompare cmp=kExact) const
RVec< PromoteType< T > > abs(const RVec< T > &v)
TGraphErrors * gr
TH1 * h
void init()
rec
Int_t Nint(T x)
TLine l
ClassImp(TPyArg)