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 
491 bool KVMultiDetArray::try_a_doubleIDtelescope(TString uri, KVDetector* de, KVDetector* e, TCollection* l)
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 
524 void KVMultiDetArray::set_up_telescope(KVDetector* de, KVDetector* e, KVIDTelescope* idt, TCollection* l)
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  KVString s = KVBase::GetDataSetEnv(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  KVString codes = GetDataSetEnv(dataset_name, Form("%s.ReconstructedNuclei.AcceptIDCodes", mda->GetName()), "");
893  if (codes != "") mda->fAcceptIDCodes.Set(codes);
894  codes = GetDataSetEnv(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(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  if (n_grids_to_write) {
920  KVSQLROOTFile f(filepath, "recreate");
921 
922  printf("Info in <KVMultiDetArray::MakeMultiDetector>: Saving %d grids in SQL-ROOT database file %s\n",
923  n_grids_to_write, filepath.Data());
924 
925  TIter it(gIDGridManager->GetGrids());
926  KVIDGraph* gr;
927  while ((gr = (KVIDGraph*)it())) {
928  // self-consistently modify runlist if necessary
929  auto rl = gr->GetRuns().GetListDividedBy(index_multiplier);
930  gr->SetRuns(rl);
931  f.WriteObject(gr, {
932  {"IDLabel", gr->GetIDTelescopeLabel()},
933  {"IDTelescopes", gr->GetParameters()->GetStringValue("IDTelescopes")},
934  {"VarX", gr->GetVarX()},
935  {"VarY", gr->GetVarY()},
936  {"Runlist", gr->GetRunList()},
937  }
938  );
939  --n_grids_to_write;
940  if (!(n_grids_to_write % 1000)) printf("Info in <KVMultiDetArray::MakeMultiDetector>: ...%d grids left...\n",
941  n_grids_to_write);
942 
943  }
944  }
945  }
946 #endif
947  }
948  else {
949  mda = gMultiDetArray;
950  }
951  // set parameters if required & allowed & not done yet
953  return mda;
954 }
955 
956 
957 
969 
971 {
972  // Return pointer to KVUpDater defined by dataset for this multidetector, the class used
973  // is defined as a plugin like this:
974  //
975  // # Plugin.KVUpDater: name_of_dataset name_of_class name_of_plugin_library constructor_to_call
976  //
977  // However, if a dataset defines a variable like this:
978  //
979  // [dataset].ExpSetUp.Updater.[multidetector name]: [name_of_dataset for plugin]
980  //
981  // then we use the updater plugin defined for the given dataset
982 
983  if (!fUpDater) {
984  KVString alt_updater = KVBase::GetDataSetEnv(fDataSet, Form("ExpSetUp.Updater.%s", GetName()), "");
985  if (alt_updater != "") fUpDater = KVUpDater::MakeUpDater(alt_updater, this);
987  }
988  Info("GetUpDater", "updater class for dataset %s: %s", fDataSet.Data(), fUpDater->IsA()->GetName());
989  return fUpDater;
990 }
991 
992 
993 
994 
1003 
1004 void KVMultiDetArray::SetParameters(UInt_t run, Bool_t physics_parameters_only)
1005 {
1006  // Set run-dependent parameters of the array.
1007  //
1008  // if physics_parameters_only==false, identification and calibration parameters are set.
1009  // if physics_parameters_only==true, just the minimum necessary for physics analysis of reduced data are set.
1010  //
1011  // This can only be done if gDataSet has been set i.e. a dataset has been chosen,
1012  // otherwise this just has the effect of setting the current run number
1013 
1014  fCurrentRun = run;
1015  KVDataSet* ds = gDataSet;
1016  if (!ds) {
1017  if (gDataSetManager)
1018  ds = gDataSetManager->GetDataSet(fDataSet.Data());
1019  }
1020  if (ds) {
1021  GetUpDater()->SetParameters(run, physics_parameters_only);
1022  SetBit(kParamsSet);
1023  }
1024 }
1025 
1026 
1027 
1028 
1033 
1035 {
1036  //Set identification parameters for run.
1037  //This can only be done if gDataSet has been set i.e. a dataset has been chosen
1038  //Otherwise this just has the effect of setting the current run number
1039 
1040  fCurrentRun = run;
1041  KVDataSet* ds = gDataSet;
1042  if (!ds) {
1043  if (gDataSetManager)
1044  ds = gDataSetManager->GetDataSet(fDataSet.Data());
1045  }
1046  if (ds) {
1049  }
1050 }
1051 
1052 
1053 
1054 
1059 
1061 {
1062  //Set calibration parameters for run.
1063  //This can only be done if gDataSet has been set i.e. a dataset has been chosen
1064  //Otherwise this just has the effect of setting the current run number
1065 
1066  fCurrentRun = run;
1067  KVDataSet* ds = gDataSet;
1068  if (!ds) {
1069  if (gDataSetManager)
1070  ds = gDataSetManager->GetDataSet(fDataSet.Data());
1071  }
1072  if (ds) {
1075  }
1076 }
1077 
1078 
1079 
1080 
1097 
1099 {
1100  //Initialisation of all ACTIVE identification telescopes in the array, i.e. those appearing in a line
1101  //in the .kvrootrc file such as this:
1102  //
1103  //# [dataset name].ActiveIdentifications: [type1] [type2] ...
1104  //
1105  //The 'types' given correspond to the value given by KVIDTelescope::GetLabel(), these are the
1106  //identifiers used to retrieve the different plugin classes in GetIDTelescopes(KVDetector*,KVDetector*,KVList*).
1107  //
1108  //For each type of identification in the list, we retrieve the first identification telescope with this
1109  //label from the list of all ID telescopes, in order to call its KVIDTelescope::SetIdentificationParameters() method.
1110  //This method (when rederived in child classes of KVIDTelescope) initialises the identification objects
1111  //for ALL of the ID telescopes of the same type (class) in the array.
1112  //
1113  //Note that, in general, the parameters of the identifications for a given run are not
1114  //set until SetParameters or SetRunIdentificationParameters is called.
1115 
1116 #ifdef WITH_RSQLITE
1117  // if a KVSQLROOTFile has been filled with the dataset identification parameters we use it.
1118  if (gIDGridManager->IsSQLROOT()) {
1119  // nothing more to do
1120  return;
1121  }
1122  else if (!gDataSet->DataBaseUpdateInProgress())
1123  // the update of the database may have been caused by 1 or more identification files being modified
1124  // (even if they are not part of the database). in this case we read in all the grids again even
1125  // if it has already been done before.
1126  {
1127  TString filepath;
1128  if (is_gnuinstall()) {
1129  // GNU-style install: use working directory $HOME/.kaliveda
1130  filepath = GetWORKDIRFilePath(gDataSet->GetName());
1131  }
1132  else
1133  filepath = gDataSet->GetDataSetDir();
1134 
1135  filepath += "/idgrids_DB";
1136 
1137  if (SearchKVFile(filepath.Data(), filepath)) {
1138  delete gIDGridManager;
1139  gIDGridManager = new KVSQLROOTIDGridManager(filepath);
1140  return;
1141  }
1142  }
1143 #endif
1144 
1145  KVString id_labels = GetDataSetEnv(fDataSet, "ActiveIdentifications", "");
1146  if (id_labels == "" || (gDataSet && !gDataSet->HasCalibIdentInfos())) {
1147  Info("SetIdentifications", "No active identifications");
1148  return;
1149  }
1150  //split list of labels
1151  id_labels.Begin(" ");
1152  int ok(0);
1153  //loop over labels/identification 'types'
1154  while (!id_labels.End()) {
1155 
1156  //get first telescope in list with right label
1158  //set ID parameters for all telescopes of this 'type'
1159  if (idt) {
1160  Info("SetIdentifications", "Initialising %s identifications...", idt->GetLabel());
1161  if (idt->SetIdentificationParameters(this))
1162  Info("SetIdentifications", "OK");
1163  ++ok;
1164  }
1165 
1166  }
1167  if (!ok) {
1168  // None of the labels in the list correspond to telescopes in the array
1169  Warning("SetIdentfications", "No telescopes found with labels given in %s.ActiveIdentifications list: %s",
1170  gDataSet->GetName(), id_labels.Data());
1171  }
1172 }
1173 
1174 
1175 
1176 
1183 
1185 {
1186  // Calls Initialize() method of each identification telescope (see KVIDTelescope
1187  // and derived classes) and sets the general identification code defined for each
1188  // telescope.
1189  //
1190  // Calling this method is essential before identification of particles is attempted.
1191 
1192  TIter next(fIDTelescopes);
1193  KVIDTelescope* idt;
1194  while ((idt = (KVIDTelescope*)next())) {
1195  idt->Initialize();
1197  }
1198 }
1199 
1200 
1201 
1212 
1214 {
1215  // Read all identification grids from the file and add them to the IDGridManager object
1216  // used by this array. This method sets up the links between each grid and the
1217  // IDtelescope(s) it is to be used for, unlike calling
1218  //
1219  // gIDGridManager->ReadAsciiFile(grids)
1220  //
1221  // which does not.
1222  //
1223  // Returns kFALSE if there is a problem reading the file
1224 
1225  if (gIDGridManager->ReadAsciiFile(grids)) {
1226  TIter next(gIDGridManager->GetLastReadGrids());
1227  KVIDGraph* gr;
1228  while ((gr = (KVIDGraph*)next())) FillListOfIDTelescopes(gr);
1229  return kTRUE;
1230  }
1231  return kFALSE;
1232 }
1233 
1234 
1235 
1236 
1240 
1242 {
1243  // Print full status report on ID telescopes in array, using informations stored in
1244  // fStatusIDTelescopes (see GetStatusOfIDTelescopes).
1245 
1246  cout << endl << "-----STATUS OF IDENTIFICATION TELESCOPES";
1247  if (GetCurrentRunNumber()) cout << " FOR RUN "
1248  << GetCurrentRunNumber();
1249  cout << "------" << endl << endl;
1250  //get list of active telescopes
1251  KVString id_labels;
1252  if (gDataSet) id_labels = gDataSet->GetDataSetEnv("ActiveIdentifications");
1253  else {
1254  unique_ptr<KVUniqueNameList> typelist(GetIDTelescopeTypes());
1255  TIter it(typelist.get());
1256  TObjString* type;
1257  while ((type = (TObjString*)it())) {
1258  if (id_labels == "") id_labels += type->GetString().Data();
1259  else {
1260  id_labels += Form(" %s", type->GetString().Data());
1261  }
1262  }
1263  }
1264  if (id_labels == "") {
1265  cout << " *** No active identifications *** " << endl;
1266  return;
1267  }
1268  // iterate over labels
1269  unique_ptr<TObjArray> toks(id_labels.Tokenize(' '));
1270 
1271  //update status infos
1273 
1274  TIter next_type(fStatusIDTelescopes);
1275  TList* id_type_list = 0;
1276  while ((id_type_list = (TList*)next_type())) {
1277 
1278  cout << " *** " << id_type_list->GetName() << " Identifications -------------------" << endl;
1279  if (!toks->FindObject(id_type_list->GetName())) {
1280  cout << " [NOT ACTIVE]" << endl;
1281  }
1282  TList* ok_list = (TList*)id_type_list->FindObject("OK");
1283  TList* notok_list = (TList*)id_type_list->FindObject("NOT OK");
1284  TList* print_list = 0;
1285  Int_t Nok = ok_list->GetEntries();
1286  Int_t Notok = notok_list->GetEntries();
1287  if (Nok && Notok) {
1288  if (Nok < Notok) print_list = ok_list;
1289  else print_list = notok_list;
1290  }
1291  if (Nok && (!Notok)) cout << " ALL telescopes are OK" << endl;
1292  else if (Notok && (!Nok)) cout << " NO telescopes are OK" << endl;
1293  else {
1294  cout << " " << ok_list->GetEntries() << " telescopes are OK, "
1295  << notok_list->GetEntries() << " telescopes are NOT OK" << endl;
1296  cout << " " << print_list->GetName() << " :" << endl;
1297  TIter it(print_list);
1298  TObject* ob = it();
1299  cout << ob->GetName();
1300  while ((ob = it())) cout << "," << ob->GetName();
1301  cout << endl;
1302  }
1303  cout << endl;
1304 
1305  }
1306 }
1307 
1308 
1309 
1310 
1318 
1320 {
1321  // Fill and return pointer to list fStatusIDTelescopes which contains
1322  // a list for each type of ID telescope in the array, each list contains a list
1323  // "OK" with the ID telescopes which have IsReadyForID()=kTRUE, and
1324  // a list "NOT OK" with the others.
1325  //
1326  // The returned TList object must not be deleted (it belongs to the KVMultiDetArray).
1327 
1328  if (!fStatusIDTelescopes) {
1329  fStatusIDTelescopes = new TList;
1331  }
1332  else {
1334  }
1336  TIter next(fIDTelescopes);
1337  KVIDTelescope* idt = 0;
1338  while ((idt = (KVIDTelescope*)next())) {
1339 
1340  TString id_type = idt->GetLabel();
1341  TList* id_type_list = (TList*)fStatusIDTelescopes->FindObject(id_type.Data());
1342  if (!id_type_list) {
1343  id_type_list = new TList;
1344  id_type_list->SetOwner(kTRUE);
1345  id_type_list->SetName(id_type.Data());
1346  fStatusIDTelescopes->Add(id_type_list);
1347  id_type_list->Add(new TList);
1348  ((TList*)id_type_list->At(0))->SetName("OK");
1349  id_type_list->Add(new TList);
1350  ((TList*)id_type_list->At(1))->SetName("NOT OK");
1351  }
1352  if (idt->IsReadyForID())
1353  ((TList*)id_type_list->FindObject("OK"))->Add(idt);
1354  else
1355  ((TList*)id_type_list->FindObject("NOT OK"))->Add(idt);
1356  }
1357  return fStatusIDTelescopes;
1358 }
1359 
1360 
1361 
1362 
1368 
1370 {
1371  // Create, fill and return pointer to a list of TObjString containing the name of each type
1372  // of ID telescope (actually the label) in the array.
1373  //
1374  // Delete the list after use (it owns the TObjString objects)
1375 
1376  KVUniqueNameList* type_list = new KVUniqueNameList(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  {
1621  TString filename = TString(GetName()) + "." + _base_filename;
1622  TString fullpath;
1625 
1626  filename = _base_filename;
1629 
1630  return TString();
1631  };
1632 
1633  // look for list of files
1634  auto fullpath = find_a_file("DetectorThicknessFiles.dat");
1635  if(!fullpath.IsNull())
1636  {
1637  KVFileReader fr;
1638  fr.OpenFileToRead(fullpath);
1639  while (fr.IsOK()) {
1640  fr.ReadLine(0);
1641  if (fr.GetCurrentLine().BeginsWith("#") || fr.GetCurrentLine() == "") continue;
1643  }
1644  return;
1645  }
1646 
1647  // look for single file
1648  fullpath = find_a_file("DetectorThicknesses.dat");
1649  if(!fullpath.IsNull()) set_detector_thicknesses(fullpath);
1650 }
1651 
1652 
1653 
1672 
1674 {
1675  // Use file given by fullpath to set the real thicknesses of the detectors.
1676  // Any detector which is not in the file will be left with its nominal thickness.
1677  //
1678  // EXAMPLE FILE:
1679  //
1680  //# thickness of detector DET01 in default units
1681  //DET01: 56.4627
1682  //
1683  //# DET03 has several layers
1684  //DET03.Abs0: 61.34
1685  //DET03.Abs1: 205.62
1686  //
1687  // !!! WARNING !!!
1688  // Single-layer detectors: The units are those defined by default for the detector's
1689  // Get/SetThickness methods.
1690  // Multi-layer: Each layer is a KVMaterial object. The thickness MUST be given in centimetres
1691  // (default thickness unit for KVMaterial).
1692 
1693  TEnv thickdat;
1694  if (thickdat.ReadFile(fullpath, kEnvUser) != 0) {
1695  Error("SetDetectorThicknesses", "Problem opening file %s", fullpath.Data());
1696  return;
1697  }
1698  Info("SetDetectorThicknesses", "Setting thicknesses of detectors from file %s", fullpath.Data());
1699  TIter next(GetDetectors());
1700  KVDetector* det;
1701  while ((det = (KVDetector*)next())) {
1702  if (thickdat.Defined(det->GetName())) {
1703  // simple single layer detector
1704  Double_t thick = thickdat.GetValue(det->GetName(), 0.0);
1705  det->SetThickness(thick);
1706  //Info("SetDetectorThicknesses", "Set thickness of %s to %f", det->GetName(), thick);
1707  }
1708  else {
1709  Char_t i = 0;
1710  TString absname;
1711  absname.Form("%s.Abs%d", det->GetName(), (Int_t)i);
1712  if (thickdat.Defined(absname.Data())) {
1713  // detector with several layers
1714  KVMaterial* abs = 0;
1715  while ((abs = det->GetAbsorber(i))) {
1716  Double_t thick = thickdat.GetValue(absname.Data(), 0.0);
1717  abs->SetThickness(thick);
1718  //Info("SetDetectorThicknesses", "Set thickness of %s.Abs%d to %f", det->GetName(), (Int_t)i, thick);
1719  i++;
1720  absname.Form("%s.Abs%d", det->GetName(), (Int_t)i);
1721  if (!thickdat.Defined(absname.Data())) break;
1722  }
1723  }
1724  }
1725  }
1726 }
1727 
1728 
1729 
1738 
1740 {
1741  // Called by KVGeoImport::ImportGeometry
1742  //
1743  // Creates KVRangeTableGeoNavigator for calculating energy losses of
1744  // particles propagated through array.
1745  //
1746  // If no name and/or title are defined for the array, the name and title of the TGeoManager
1747  // object will be used for the array.
1748 
1749  if (!strcmp(GetName(), "")) SetName(g->GetName());
1750  if (!strcmp(GetTitle(), "")) SetTitle(g->GetTitle());
1752 }
1753 
1754 
1755 
1784 
1786 {
1787  // Method for positioning volumes in detector geometries
1788  //
1789  // Given:
1790  //
1791  // distance [cm] = distance from target (origin) to the CENTRE of the volume in position
1792  // theta [deg] = polar angle of vector from target to centre of volume in position
1793  // phi [deg] = azimuthal angle of vector
1794  //
1795  // this method generates the matrix which is required to position the volume as required
1796  // while also turning the volume so that the side nearest the target (i.e. the entrance
1797  // window of the detector) remains perpendicular to the vector joining the origin and
1798  // the centre of the volume.
1799  //
1800  // If required, a further translation can be given which will be applied to the volume after
1801  // it has been placed with the required orientation at the nominal distance. This can be used
1802  // e.g. for detector misalignment, when detectors are in a structure which guarantees their line
1803  // of sight to be orthogonal to their surface at a nominal distance, but the nominal distance
1804  // is not respected.
1805  //
1806  // Example of use:
1807  //
1808  //~~~~~~~~~~~~
1809  // TGeoVolume* vol;// volume to be positioned
1810  // double depth = vol->GetShape()->GetDZ(); // half-width of volume in direction of target
1811  // // place front of volume at 100cm, with theta=45 deg. and phi=60 deg.
1812  // gGeoManager->GetTopVolume()->AddNode(vol, 1, KVMultiDetArray::GetVolumePositioningMatrix(100+depth,45,60));
1813  //~~~~~~~~~~~~
1814 
1815  TGeoRotation rot;
1816  TGeoTranslation trans;
1817  rot.SetAngles(phi + 90, theta, -90) ;
1818  trans.SetDz(distance) ;
1819  TGeoHMatrix h;
1820  if (postTrans) h = (*postTrans) * rot * trans ;
1821  else h = rot * trans;
1822  TGeoHMatrix* ph = new TGeoHMatrix(h);
1823  return ph;
1824 }
1825 
1826 
1827 
1828 
1832 
1834 {
1835  // For each grid which is valid for this run, we call the KVIDTelescope::SetIDGrid method
1836  // of each associated ID telescope.
1837 
1838  if (gIDGridManager->IsSQLROOT()) {
1839  gIDGridManager->LoadGridsForRun(run);
1840  }
1841 
1842  TIter next(gIDGridManager->GetGrids());
1843  KVIDGraph* gr = 0;
1844  while ((gr = (KVIDGraph*) next())) {
1845  if (!gIDGridManager->IsSQLROOT()) {
1846  if (!gr->GetRuns().Contains((Int_t) run))
1847  continue;
1848  }
1849  else
1851 
1852  TIter nxtid(gr->GetIDTelescopes());
1853  KVIDTelescope* idt;
1854  while ((idt = (KVIDTelescope*) nxtid())) {
1855  idt->SetIDGrid(gr);
1856  }
1857  }
1858 }
1859 
1860 
1864 
1866 {
1867  // Fill list of ID telescopes with which this grid is associated
1868  // from list of names read from ascii file.
1869 
1870  gr->ClearListOfTelescopes();
1871  if (gr->GetParameters()->HasParameter("IDTelescopes")) {
1872  KVString tel_list = gr->GetParameters()->GetStringValue("IDTelescopes");
1873  tel_list.Begin(",");
1874  while (!tel_list.End()) {
1875  TString tel_name = tel_list.Next();
1876  KVIDTelescope* idt = GetIDTelescope(tel_name.Data()) ;
1877  if (idt) gr->AddIDTelescope(idt);
1878  }
1879  }
1880 }
1881 
1882 
1883 
1893 
1895 {
1896  // Use OpenGL viewer to view multidetector geometry
1897  //
1898  // If option="tracks" we draw any tracks corresponding to the last simulated
1899  // event whose detection was simulated with DetectEvent
1900  // If option="tracks:[numberlist]" with a list of numbers,
1901  // it will be interpreted as a KVNumberList containing the Z of tracks to be drawn
1902  // e.g. option="tracks:1-92" draw only tracks with 1<=Z<=92 (no neutrons)
1903  // option="tracks:2" draw only helium isotopes
1904 
1905  GetGeometry()->GetTopVolume()->Draw("ogl");
1906  KVString opt(option);
1907  opt.Begin(":");
1908  if (opt.Next() == "tracks") {
1909  if (!opt.End()) {
1910  KVNumberList zlist(opt.Next());
1911  GetNavigator()->DrawTracks(&zlist);
1912  }
1913  else
1914  GetNavigator()->DrawTracks();
1915  }
1916 #ifdef WITH_OPENGL
1917  TGLViewer* view = (TGLViewer*)gPad->GetViewer3D();
1920  view->SetSmoothLines(kTRUE);
1921  view->SetSmoothPoints(kTRUE);
1922 #endif
1923 }
1924 
1925 
1926 
1928 
1930 {
1931  fNavigator.reset(geo);
1932 }
1933 
1934 
1935 
1940 
1942 {
1943  // Create TH2F histograms for all IDTelescopes of the array
1944  // They will be added to the list
1945  // histograms will have resolution of dimension*dimension
1946 
1948  KVIDTelescope* idt;
1949  while ((idt = (KVIDTelescope*)it())) {
1950  TString name(idt->GetName());
1951  name.ReplaceAll("-", "_");
1952  list->Add(new TH2F(name, Form("Hits in %s", idt->GetName()), dimension, 0., 0., dimension, 0., 0.));
1953  }
1954 }
1955 
1956 
1957 
1960 
1962 {
1963  // Fill TH2F histograms for all IDTelescopes of the array
1964 
1966  KVIDTelescope* idt;
1967  while ((idt = (KVIDTelescope*)it())) {
1968  TString name(idt->GetName());
1969  name.ReplaceAll("-", "_");
1970  TH2F* h = (TH2F*)list->FindObject(name);
1971  if (h) h->Fill(idt->GetIDMapX(), idt->GetIDMapY());
1972  }
1973 }
1974 
1975 
1976 
1979 
1981 {
1982  // Modify the transparency of detector volumes in OpenGL view
1983 
1984  TIter itV(GetGeometry()->GetListOfVolumes());
1985  TGeoVolume* vol;
1986  while ((vol = (TGeoVolume*)itV())) vol->SetTransparency(t);
1987 }
1988 
1989 
1990 
1994 
1996 {
1997  // Calculate all possible (sub-)trajectories
1998  // for particle reconstruction (GetReconTrajectories())
1999 
2000  unique_ptr<KVSeqCollection> groups(GetStructureTypeList("GROUP"));
2001  TIter it(groups.get());
2002  KVGroup* group;
2003  Int_t ntr = 0;
2004  Info("CalculateReconstructionTrajectories", "Calculating trajectories for particle reconstruction:");
2006  std::cout << "\xd" << " -- calculated " << ntr << " reconstruction trajectories" << std::flush;
2007  while ((group = (KVGroup*)it())) {
2008  ntr += group->CalculateReconstructionTrajectories();
2010  std::cout << "\xd" << " -- calculated " << ntr << " reconstruction trajectories" << std::flush;
2011  }
2013  std::cout << " -- calculated " << ntr << " reconstruction trajectories" << std::endl;
2014  else
2015  std::cout << std::endl;
2016 }
2017 
2018 
2019 
2024 
2026 {
2027  // Track over all possible particle trajectories calling GetIDTelescopes(KVDetector*,KVDetector*)
2028  // for each pair of (present & functioning) detectors.
2029  // This will create all possible KVIDTelescope identification objects and put them in list fIDTelescopes
2030 
2031  fIDTelescopes->Delete();
2032  TIter next_traj(GetTrajectories());
2033  KVGeoDNTrajectory* traj;
2034  Int_t count = 0;
2035  Info("DeduceIdentificationTelescopesFromGeometry", "Calculating...");
2037  std::cout << "\xd" << " -- created " << count << " telescopes" << std::flush;
2038  while ((traj = (KVGeoDNTrajectory*)next_traj())) { // loop over all trajectories
2039 
2040  traj->IterateFrom(); // from furthest-out to closest-in detector
2041 
2043  while ((N = traj->GetNextNode())) {
2044  KVGeoDetectorNode* Nplus1 = traj->GetNodeInFront(N);
2045  count += GetIDTelescopes(Nplus1, N->GetDetector(), traj->AccessIDTelescopeList());
2047  std::cout << "\xd" << " -- created " << count << " telescopes" << std::flush;
2048  }
2049  }
2051  std::cout << " -- created " << count << " telescopes" << std::endl;
2052  else
2053  std::cout << std::endl;
2054 }
2055 
2056 
2057 
2061 
2063 {
2064  // Eliminate any trajectories which are just sub-trajectories of others
2065  // For each trajectory in list fTrajectories, we add a reference to the trajectory to each node on the trajectory
2066 
2067  TIter it(&fTrajectories);
2068  KVGeoDNTrajectory* tr;
2069  KVList duplicates;
2070  // look for duplicate sub-trajectories
2071  while ((tr = (KVGeoDNTrajectory*)it())) {
2072  int len_tr = tr->GetN();
2073  TIter it2(&fTrajectories);
2074  KVGeoDNTrajectory* tr2;
2075  while ((tr2 = (KVGeoDNTrajectory*)it2())) {
2076  if ((tr2 != tr) && (len_tr < tr2->GetN()) && (tr2->ContainsPath(tr))) {
2077  duplicates.Add(tr);
2078  break;
2079  }
2080  }
2081  }
2082  // remove duplicates
2083  if (duplicates.GetEntries()) {
2084  TIter it_dup(&duplicates);
2085  while ((tr = (KVGeoDNTrajectory*)it_dup())) {
2086  fTrajectories.Remove(tr);
2087  }
2088  Info("AssociateTrajectoriesAndNodes", "Removed %d duplicated sub-trajectories", duplicates.GetEntries());
2089  }
2090  Info("AssociateTrajectoriesAndNodes", "Calculated %d particle trajectories", fTrajectories.GetEntries());
2091  it.Reset();
2092  while ((tr = (KVGeoDNTrajectory*)it())) {
2093  tr->AddToNodes();
2094  }
2095 }
2096 
2097 
2098 
2100 
2102 {
2103  if (N->GetNTraj() > 1) {
2104  if (!multitraj_nodes.FindObject(N)) { // look for any detectors which are on multiple trajectories
2105  //cout << "multitraj node found: " << N->GetName() << " (" << N->GetNTraj() << ")" << endl;
2106  multitraj_nodes.Add(N);
2107  TIter tr(N->GetTrajectories());
2108  KVGeoDNTrajectory* traj;
2109  while ((traj = (KVGeoDNTrajectory*)tr())) { // for each trajectory associated with detector
2110  if (tried_trajectories.FindObject(traj)) continue; // trajectory already used
2111  tried_trajectories.Add(traj);
2112  traj->IterateFrom();
2113  KVGeoDetectorNode* node;
2114  while ((node = traj->GetNextNode())) { // store names of all detectors on trajectory
2115  detectors_of_group.Add(node);
2116  RecursiveTrajectoryClustering(node, tried_trajectories, multitraj_nodes, detectors_of_group);
2117  }
2118  }
2119  }
2120  }
2121  else if (N->GetNTraj() == 1) {
2122  // single-trajectory node.
2123  // work along trajectory adding nodes to group
2124  KVGeoDNTrajectory* traj = (KVGeoDNTrajectory*)N->GetTrajectories()->First();
2125  if (tried_trajectories.FindObject(traj)) return; // trajectory already used
2126  tried_trajectories.Add(traj);
2127  traj->IterateFrom();
2128  KVGeoDetectorNode* node;
2129  while ((node = traj->GetNextNode())) { // store names of all detectors on trajectory
2130  detectors_of_group.Add(node);
2131  RecursiveTrajectoryClustering(node, tried_trajectories, multitraj_nodes, detectors_of_group);
2132  }
2133  }
2134  else {
2135  // orphan node? single-detector array?
2136  detectors_of_group.Add(N);
2137  }
2138 }
2139 
2140 
2141 
2150 
2152 {
2153  // Create and return pointer to new KVGroupReconstructor for reconstructing particles
2154  // in the given group. Returns nullptr if group is not part of this array.
2155  //
2156  // Plugins for specific arrays can be defined as plugins using the name of the array:
2157  // +Plugin.KVGroupReconstructor: my_array my_group_reconstructor my_lib "my_group_reconstructor()"
2158  //
2159  // If we are in 'SimMode', the default reconstructor is a KVFilterGroupReconstructor
2160 
2161  KVGroupReconstructor* gr(nullptr);
2162  if (GetGroup(g->GetName())) {
2163  // look for plugin
2165  if (!gr) {
2166  if (IsSimMode()) return KVGroupReconstructor::Factory("Filter", g);
2167  else
2168  gr = new KVGroupReconstructor(g);
2169  }
2170  }
2171  return gr;
2172 }
2173 
2174 
2175 
2179 
2181 {
2182  // Call this method just after opening a raw data file in order to perform any
2183  // necessary initialisations, depending on the type of data
2184 
2185 #ifdef WITH_BUILTIN_GRU
2186  if (r->GetDataFormat() == "EBYEDAT")
2187  dynamic_cast<KVGANILDataReader*>(r)->ConnectRawDataParameters();
2188 #endif
2189 }
2190 
2191 
2192 
2197 
2199 {
2200  // Deduce the "groups" in the array from the trajectories
2201  // Any trajectories with 1 or more common detectors define a group.
2202  // The group is constituted of all detectors belonging to the trajectories of the group.
2203 
2204  Info("DeduceGroupsFromTrajectories", "Deducing groups of detectors from trajectories");
2205  Int_t number_of_groups = 0;
2206  TIter next_det(GetDetectors());
2207  unique_ptr<KVSeqCollection> stl(GetStructureTypeList("GROUP"));
2208  if (stl.get() && stl->GetEntries()) {
2209  Info("DeduceGroupsFromTrajectories", "Deleting existing %d groups in array", stl->GetEntries());
2210  ClearStructures("GROUP");
2211  Info("DeduceGroupsFromTrajectories", "Done");
2212  }
2213  KVDetector* det;
2214  KVUniqueNameList tried_trajectories;//avoid double-counting/infinite loops
2215  KVUniqueNameList multitraj_nodes;//avoid double-counting/infinite loops
2216  while ((det = (KVDetector*) next_det())) {
2217  if (det->GetGroup()) continue; // group assignment already done
2218  KVUniqueNameList detectors_of_group;
2219  RecursiveTrajectoryClustering(det->GetNode(), tried_trajectories, multitraj_nodes, detectors_of_group);
2220  if (!detectors_of_group.GetEntries()) continue;
2221  KVGroup* Group = new KVGroup;
2222  Group->SetNumber(++number_of_groups);
2223  Add(Group);
2224  TIter next_node(&detectors_of_group);
2226  while ((d = (KVGeoDetectorNode*)next_node())) Group->Add(d->GetDetector());
2227  }
2228  TIter tr(&fTrajectories);
2229  KVGeoDNTrajectory* t;
2230  Info("DeduceGroupsFromTrajectories", "Filling group trajectory lists");
2231  while ((t = (KVGeoDNTrajectory*)tr())) {
2232  if (t->GetNodeAt(0)->GetDetector()->GetGroup())
2233  t->GetNodeAt(0)->GetDetector()->GetGroup()->AddTrajectory(t);
2234  else {
2235  t->Print();
2236  t->GetNodeAt(0)->GetDetector()->Print();
2237  }
2238  }
2239 }
2240 
2241 
2242 
2246 
2248 {
2249  // Called when required to fill KVReconstructedNucleus::fDetList with pointers to
2250  // the detectors whose names are stored in KVReconstructedNucleus::fDetNames.
2251 
2252  DetList->Clear();
2253  DetNames.Begin("/");
2254  while (!DetNames.End()) {
2255  KVDetector* det = GetDetector(DetNames.Next(kTRUE));
2256  if (det) DetList->Add(det);
2257  }
2258 }
2259 
2260 
2261 
2281 
2283 {
2284  // Set status of particle by comparing its identification/calibration codes
2285  // with those set as acceptable in fAcceptIDCodes and fAcceptECodes.
2286  // The status can be tested with method KVReconstructedNucles::IsOK().
2287  //
2288  // The default lists are defined in variables of the form
2289  //```
2290  // [DataSet].[name].ReconstructedNuclei.AcceptIDCodes: [list]
2291  // [DataSet].[name].ReconstructedNuclei.AcceptECodes: [list]
2292  //```
2293  // where:
2294  // + `DataSet` is an optional dataset name for dataset-specific lists
2295  // + `name` is the name of the multidetector array
2296  // + `list` is a numeric list (KVNumberList format)
2297  //
2298  // The default lists can be overridden using methods AcceptIDCodes(), AcceptECodes(),
2299  // AcceptAllIDCodes() and AcceptAllECodes().
2300  //
2301  // If either list is empty, no selection is made for the corresponding code
2302 
2303  Bool_t ok = kTRUE;
2305  if (!fAcceptECodes.IsEmpty()) ok = ok && fAcceptECodes.Contains(NUC->GetECode());
2306  NUC->SetIsOK(ok);
2307 }
2308 
2309 
2310 #ifdef WITH_BUILTIN_GRU
2311 
2314 
2316 {
2317  // General method for reading raw data in old GANIL ebyedat format
2318  AbstractMethod("handle_raw_data_event_ebyedat");
2319  return kFALSE;
2320 }
2321 
2322 #endif
2323 
2324 
2327 
2329 {
2330  // reset acquisition parameters etc. before reading new raw data event
2331 
2334  fHandledRawData = false;
2335  // reset fired signals
2336  TIter nxt(&fFiredSignals);
2337  KVDetectorSignal* ds;
2338  while ((ds = (KVDetectorSignal*)nxt())) {
2339  ds->SetFired(false);
2340  ds->SetValue(0);
2341  }
2342  fFiredSignals.Clear();
2343 }
2344 
2345 
2346 
2354 
2356 {
2357  // Perform any operations to finalise the description of the multidetector
2358  // which can only be done once the geometry is closed, e.g. use KVGeoImport
2359  // to set up nodes, trajectories, detectors, idtelescopes, etc.
2360  // This has to be kept separate for use with KVExpSetUp which first fills
2361  // a single ROOT geometry with all component KVMultiDetArray geometries,
2362  // then closes the geometry only when all have been built.
2363 }
2364 
2365 
2366 
2369 
2371 {
2372  // Copy any parameters in fReconParameters in to the reconstructed event parameter list
2373  e->GetParameters()->Concatenate(fReconParameters);
2374 }
2375 
2376 
2377 
2390 
2392 {
2393  // values of fired raw data signals (acquisition parameters) from last read raw event
2394  // are copied to the fReconParameters list of parameters to be stored with the
2395  // reconstructed event.
2396  //
2397  // the format for each signal is:
2398  //
2399  // ACQPAR.[array].[detector].[signal]
2400  // ACQPAR.[array].[signal]
2401  //
2402  // in the first case for signals associated with detectors, in the latter case signals
2403  // which are not associated with a detector
2404 
2405  TIter it(GetFiredSignals());
2406  KVDetectorSignal* o;
2407  while ((o = (KVDetectorSignal*)it())) {
2408  fReconParameters.SetValue(Form("ACQPAR.%s.%s", GetName(), o->GetFullName().Data()), o->GetValue());
2409  }
2410 }
2411 
2412 
2413 
2422 
2424 {
2425  // Update array according to last event read using the KVRawDataReader object
2426  // (it is assumed that KVRawDataReader::GetNextEvent() was called before calling this method)
2427  //
2428  // Return kTRUE if raw data was treated
2429  //
2430  // All fired acquisition parameters are written in the fReconParameters list,
2431  // ready to be copied to the reconstructed event
2432 
2433  fRawDataReader = rawdata;
2435  if (rawdata->GetDataFormat() == "MFM") {
2436 #ifdef WITH_MFM
2438 #endif
2439  }
2440  else if (rawdata->GetDataFormat() == "PROTOBUF") {
2441 #ifdef WITH_PROTOBUF
2443 #endif
2444  }
2445  else if (rawdata->GetDataFormat() == "EBYEDAT") {
2446 #ifdef WITH_BUILTIN_GRU
2448 #endif
2449  }
2450  if (fHandledRawData) {
2452  }
2453  return fHandledRawData;
2454 }
2455 
2456 
2457 #ifdef WITH_MFM
2458 
2467 
2469 {
2470  // Update array according to last event read from MFM buffer
2471  // (it is assumed that MFMBufferReader::ReadNextFrame() was called before calling this method)
2472  //
2473  // Return kTRUE if raw data was treated
2474  //
2475  // All fired acquisition parameters are written in the fReconParameters list,
2476  // ready to be copied to the reconstructed event
2477 
2479  bool ok = false;
2480  ok = handle_raw_data_event_mfmfile(bufrdr);
2481  if (ok) {
2483  }
2484  return ok;
2485 }
2486 
2487 #endif
2488 
2489 
2500 
2501 void KVMultiDetArray::add_and_set_detector_signal(KVDetector* detector, KVString detname, Double_t sig_data, KVString sig_type)
2502 {
2503  // Given a pointer to a detector (may be nullptr) try to set the data sig_data in the associated signal
2504  // of the given type, sig_typ.
2505  //
2506  // If the signal does not exist for the detector, it will be created.
2507  //
2508  // If detector pointer is null, the signal will be looked for in fExtraRawDataSignals
2509  // and created if necessary.
2510  //
2511  // All fired detectors and signals are added to the lists fFiredDetectors and fFiredSignals.
2512 
2513  KVDetectorSignal* det_signal = nullptr;
2514  if (detector) {
2515  det_signal = detector->GetDetectorSignal(sig_type);
2516  if (!det_signal) {
2517  det_signal = detector->AddDetectorSignal(sig_type);
2518  }
2519  fFiredDetectors.Add(detector);
2520  }
2521  else {
2522  // raw data not associated with a detector
2523  TString sig_name;
2524  if (detname != "") sig_name = Form("%s.%s", detname.Data(), sig_type.Data());
2525  else sig_name = sig_type;
2526  det_signal = fExtraRawDataSignals.get_object<KVDetectorSignal>(sig_name);
2527  if (!det_signal) {
2528  det_signal = new KVDetectorSignal(sig_name);
2529  fExtraRawDataSignals.Add(det_signal);
2530  }
2531  }
2532  if (det_signal) {
2533  det_signal->SetValue(sig_data);
2534  det_signal->SetFired();
2535  fFiredSignals.Add(det_signal);
2536  }
2537 }
2538 
2539 
2540 
2546 
2548 {
2549  // Take values 'ACQPAR.[array_name].[detname].[signal]' or 'ACQPAR.[array_name].[signal]'
2550  // in the parameter list and use them to set values of raw acquisition parameters.
2551  //
2552  // Any detector signals which don't already exist will be created
2553 
2554  prepare_to_handle_new_raw_data(); // clear previous fired parameters/detectors
2555  int N = l.GetNpar();
2556  for (int i = 0; i < N; ++i) {
2557  KVNamedParameter* np = l.GetParameter(i);
2558 
2559  KVString name(np->GetName());
2560  if (name.BeginsWith("ACQPAR")) {
2561  // 3 '.' => 4 values means associated detector
2562  // 2 '.' => 3 values means no detector
2563  int dots = name.GetNValues(".");
2564  bool with_det = (dots == 4);
2565  assert(with_det || (dots == 3)); // sanity check
2566  name.Begin(".");
2567  name.Next(); // "ACQPAR"
2568  if (name.Next() != GetName()) continue; // check name of array - somebody else's parameter ?
2569  KVString det_name;
2570  KVString sig_type;
2571  KVDetector* det = nullptr;
2572  if (with_det) {
2573  det_name = name.Next();
2574  sig_type = name.Next();
2575  det = GetDetector(det_name);
2576  }
2577  else {
2578  sig_type = name.Next();
2579  }
2580  add_and_set_detector_signal(det, det_name, np->GetDouble(), sig_type);
2581  }
2582  }
2583 }
2584 
2585 
2586 
2605 
2607 {
2608  // We first look for following files with the name given by
2609  //
2610  // [dataset].[name].OoODetectors: [name.OoODetectors.dat]
2611  // which should contain the runlists for each malfunctioning detector.
2612  // If found we add to the experiment database a table '[name].OoO Detectors' where [name] is the name of this array.
2613  //
2614  // [dataset].[name].AbsentDetectors: [name.AbsentDetectors.dat]
2615  // which should contain the runlists for each absent detector.
2616  // If found we add to the experiment database a table '[name].Absent Detectors' where [name] is the name of this array.
2617  //
2618  // Then we look for a file with the name given by
2619  //
2620  // [dataset].[name].CalibrationFiles: [CalibrationFiles.dat]
2621  //
2622  // which should contain the names of files to read with each type of calibration
2623  // If found we add to the experiment database a table '[name].Calibrations' where [name] is the name of this array,
2624  // containing all calibrations as KVDBParameterSet objects with the name of the detector concerned.
2625  db->SetDBType(Form("%sDB", GetName()));
2626  ReadOoODetectors(db);
2627  ReadAbsentDetectors(db);
2629 }
2630 
2631 
2632 
2634 
2636 {
2637  TString basic_name = db->GetCalibFileName(keyw);
2638  if (basic_name == "") {
2639  Info(meth, "No name found for \"%s\" file", keyw);
2640  return "";
2641  }
2642  Info(meth, "Search for %s for dataset %s ...", basic_name.Data(), fDataSet.Data());
2643  TString fp;
2644  SearchKVFile(basic_name.Data(), fp, fDataSet);
2645  if (fp == "") {
2646  Info(meth, "\tNo file found ...");
2647  }
2648  return fp;
2649 }
2650 
2651 
2652 
2654 
2655 unique_ptr<KVFileReader> KVMultiDetArray::GetKVFileReader(KVExpDB* db, const Char_t* meth, const Char_t* keyw)
2656 {
2657 
2658  TString fp = GetFileName(db, meth, keyw);
2659  if (fp == "")
2660  return unique_ptr<KVFileReader>();
2661 
2662  unique_ptr<KVFileReader> fr(new KVFileReader());
2663  if (!fr->OpenFileToRead(fp.Data())) {
2664  Error(meth, "Error in opening file %s", fp.Data());
2665  fr.reset(nullptr);
2666  }
2667  else
2668  Info(meth, "Reading %s file", fp.Data());
2669  return fr;
2670 }
2671 
2672 
2673 
2675 
2677 {
2678 
2679  unique_ptr<KVFileReader> fr = GetKVFileReader(db, "ReadCalibrationFiles()", "CalibrationFiles");
2680  if (!fr.get())
2681  return;
2682 
2683  KVDBTable* calib_table = db->AddTable(Form("%s.Calibrations", GetName()), Form("Calibrations for %s", GetName()));
2684  while (fr->IsOK()) {
2685  fr->ReadLine(0);
2686  if (fr->GetCurrentLine().BeginsWith("#") || fr->GetCurrentLine() == "") {}
2687  else {
2688  ReadCalibFile(fr->GetCurrentLine().Data(), db, calib_table);
2689  }
2690  }
2691  fr->CloseFile();
2692 }
2693 
2694 
2695 
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  //The `[CalibClass]`, if given, must correspond to a KVCalibrator plugin name. The list of plugin names and the corresponding
2758  //classes can be retrieved with
2759  //
2760  //~~~~~~~~~~~
2761  //KVBase::GetListOfPlugins("KVCalibrator")
2762  //KVBase::GetListOfPluginURIs("KVCalibrator")
2763  //~~~~~~~~~~~
2764  //
2765  //KVCalibrator objects are added to detectors as required by the contents of calibration files.
2766  //If any detector has an existing calibrator of type `[CalibType]` which is not of the given class
2767  //it will be replaced with a new calibrator corresponding to the plugin.
2768  //
2769  //The `[CalibOptions]` is optional: list in `[CalibOptions]` will be used
2770  //to complete set-up of any new calibrator objects by calling the KVCalibrator::SetOptions()
2771  //method.
2772  //
2773  //`[CalibOptions]` should hold a comma-separated list of `parameter=value` pairs which will be used
2774  //to fill a KVNameValueList for the method call. See the KVCalibrator::SetOptions() method.
2775  //
2776  //`[ZRange]` is an option if several calibrations need to be used to provide the same signal
2777  //for certain detectors depending on the atomic number Z of the particle detected.
2778 
2779 
2780  TString fullpath = "";
2781  if (!SearchKVFile(filename, fullpath, fDataSet)) {
2782  Info("ReadCalibFile", "%s does not exist or not found", filename);
2783  return;
2784  }
2785 
2786  Info("ReadCalibFile", "file : %s found", fullpath.Data());
2787  TEnv env;
2788  env.ReadFile(fullpath, kEnvAll);
2789 
2790  // read options from file
2791  KVNameValueList options;
2792  KVString opt_list = "RunList SignalIn SignalOut CalibType CalibClass CalibOptions ZRange";
2793  opt_list.Begin(" ");
2794  while (!opt_list.End()) {
2795  KVString opt = opt_list.Next();
2796  KVString opt_val = env.GetValue(opt, "");
2797  opt_val.Remove(TString::kBoth, ' ');
2798  options.SetValue(opt, opt_val.Data());
2799  }
2800  // check for stupid spellnig mitskaes
2801  if (TString(env.GetValue("Runlist", "")) != "") {
2802  Warning("ReadCalibFile", "Calibration has 'Runlist' parameter (ignored): %s, did you mean 'RunList'?", env.GetValue("Runlist", ""));
2803  }
2804 
2805  if (options.GetTStringValue("SignalIn") == "") {
2806  Error("ReadCalibFile", "No input signal defined : SignalIn");
2807  return;
2808  }
2809  if (options.GetTStringValue("SignalOut") == "") {
2810  Error("ReadCalibFile", "No output signal defined : SignalOut");
2811  return;
2812  }
2813  if (options.GetTStringValue("CalibType") == "") {
2814  Error("ReadCalibFile", "No calibration type defined : CalibType");
2815  return;
2816  }
2817  Bool_t check_class(options.GetTStringValue("CalibClass") != "");
2818  TString calibrator_class;
2819  if (check_class) {
2820  TPluginHandler* ph = LoadPlugin("KVCalibrator", options.GetStringValue("CalibClass"));
2821  if (ph) calibrator_class = ph->GetClass();
2822  else {
2823  Error("ReadCalibFile", "No calibrator plugin of type %s", options.GetStringValue("CalibClass"));
2824  return;
2825  }
2826  }
2827 
2828  KVString clop;
2829  if (options.HasParameter("CalibOptions")) clop = options.GetStringValue("CalibOptions");
2830 
2831  KVString zrange;
2832  if (options.HasParameter("ZRange")) zrange = options.GetStringValue("ZRange");
2833 
2834  KVNumberList run_list = db->GetRunList();
2835  if (options.GetTStringValue("RunList") != "") {
2836  run_list.Set(options.GetTStringValue("RunList"));
2837  Info("ReadCalibFile", "Calibration used for runs %s", run_list.AsString());
2838  }
2839  else {
2840  Info("ReadCalibFile", "Calibration used for all runs in database");
2841  }
2842 
2843  TIter next(env.GetTable());
2844  TEnvRec* rec = 0;
2845  KVDBParameterSet* par = 0;
2846 
2847  while ((rec = (TEnvRec*)next())) {
2848 
2849  TString sname(rec->GetName());
2850  KVDetector* det = GetDetector(sname);
2851  if (!det) continue;
2852 
2853  KVString lval(rec->GetValue());
2854  par = new KVDBParameterSet(sname.Data(), options.GetStringValue("CalibType"), lval.GetNValues(","));
2855  par->SetParameter("SignalIn", options.GetStringValue("SignalIn"));
2856  par->SetParameter("SignalOut", options.GetStringValue("SignalOut"));
2857  // put infos on required calibrator class into database so that it can be replaced
2858  // as needed in SetCalibratorParameters
2859  par->SetParameter("CalibClass", options.GetStringValue("CalibClass"));
2860  if (clop != "") par->SetParameter("CalibOptions", clop);
2861  if (zrange != "") par->SetParameter("ZRange", zrange);
2862  Int_t np = 0;
2863  lval.Begin(",");
2864  while (!lval.End()) {
2865  par->SetParameter(np++, lval.Next().Atof());
2866  }
2867  calib_table->AddRecord(par);
2868  db->LinkRecordToRunRange(par, run_list);
2869 
2870  }
2871 }
2872 
2873 
2874 #ifdef WITH_MFM
2875 
2881 
2883 {
2884  // Update array according to last event read using the KVMFMDataFileReader object
2885  // (it is assumed that KVRawDataReader::GetNextEvent() was called before calling this method)
2886  //
2887  // Return kTRUE if raw data was treated
2888 
2889  if (mfmreader.IsFrameReadMerge()) {
2890  return handle_raw_data_event_mfmmergeframe(mfmreader.GetMergeManager());
2891  }
2892  else {
2893  return handle_raw_data_event_mfmframe(mfmreader.GetFrameRead());
2894  }
2895  return kFALSE;
2896 }
2897 
2898 
2899 
2903 
2904 Bool_t KVMultiDetArray::handle_raw_data_event_mfmmergeframe(const MFMMergeFrameManager& mergeframe)
2905 {
2906  // Method used to handle merged MFM frames
2907  // We call handle_raw_data_event_mfmframe() for each frame contained in the merge
2908 
2909  Bool_t ok = false;
2910  while (mergeframe.ReadNextFrame()) {
2911  Bool_t me = handle_raw_data_event_mfmframe(mergeframe.GetFrameRead());
2912  ok = (ok || me);
2913  }
2914  return ok;
2915 }
2916 
2917 
2918 
2929 
2931 {
2932  // Method used to treat raw data in MFM format read by KVMFMDataFileReader
2933  //
2934  // Here we dispatch two types of frame - MFMEbyedatFrame & MFMMesytecMDPPFrame -
2935  // to specific methods - handle_raw_data_event_mfmframe_ebyedat() and
2936  // handle_raw_data_event_mfmframe_mesytec_mdpp()
2937  // which need to be implemented in child classes for specific arrays which
2938  // use these data formats.
2939  //
2940  // Return kTRUE if raw data was treated
2941 #ifdef WITH_MESYTEC
2942  if (mfmframe.GetFrameType() == MFM_MESYTEC_MDPP_FRAME_TYPE)
2943  return handle_raw_data_event_mfmframe_mesytec_mdpp((const MFMMesytecMDPPFrame&)mfmframe);
2944 #endif
2945  if (mfmframe.GetFrameType() == MFM_EBY_EN_FRAME_TYPE
2946  || mfmframe.GetFrameType() == MFM_EBY_TS_FRAME_TYPE
2947  || mfmframe.GetFrameType() == MFM_EBY_EN_TS_FRAME_TYPE)
2948  return handle_raw_data_event_mfmframe_ebyedat((const MFMEbyedatFrame&)mfmframe);
2949 
2950  return kFALSE;
2951 }
2952 
2953 
2954 
2957 
2959 {
2960  // Read a raw data event from a EBYEDAT MFM Frame.
2961 
2962  AbstractMethod("handle_raw_data_event_mfmframe_ebyedat");
2963  return kFALSE;
2964 }
2965 
2966 
2967 #ifdef WITH_MESYTEC
2968 
2971 
2973 {
2974  // Read a raw data event from a Mesytec MFM Frame.
2975 
2976  AbstractMethod("handle_raw_data_event_mfmframe_mesytec_mdpp");
2977  return kFALSE;
2978 }
2979 
2980 #endif
2981 #endif
2982 
2983 #ifdef WITH_PROTOBUF
2984 
2986 
2988 {
2989  AbstractMethod("handle_raw_data_event_protobuf");
2990  return kFALSE;
2991 }
2992 
2993 #endif
2994 
2995 
2998 
3000 {
3001  // For each IDtelescope in array, calculate an identification grid
3002 
3003  TIter nxtid(GetListOfIDTelescopes());
3004  KVIDTelescope* idt;
3005  while ((idt = (KVIDTelescope*) nxtid())) {
3006  idt->CalculateDeltaE_EGrid("1-92", 0, 20);
3007  }
3008 }
3009 
3010 
3011 
3016 
3018 {
3019  // Sets status of detectors (KVDetector::IsPresent() and KVDetector::IsWorking()) for a given run of a dataset.
3020  //
3021  // If 'myname' is given, we look in database table "myname.OoODets"
3022 
3023  KVRList* absdet = (myname != "" ? kvrun->GetLinks(Form("%s.Absent Detectors", myname.Data())) : kvrun->GetLinks("Absent Detectors"));
3024  KVRList* ooodet = (myname != "" ? kvrun->GetLinks(Form("%s.OoO Detectors", myname.Data())) : kvrun->GetLinks("OoO Detectors"));
3025 
3026  TIter next(GetDetectors());
3027  KVDetector* det;
3028 
3029  Int_t ndet_absent = 0;
3030  Int_t ndet_ooo = 0;
3031  TString absent_dets, ooo_dets;
3032 
3033  while ((det = (KVDetector*)next())) {
3034  //Test de la presence ou non du detecteur
3035  if (!absdet) {
3036  det->SetPresent();
3037  }
3038  else {
3039  if (absdet->FindObject(det->GetName())) {
3040  det->SetPresent(kFALSE);
3041  if (ndet_absent) absent_dets += ",";
3042  absent_dets += det->GetName();
3043  ndet_absent += 1;
3044  }
3045  else {
3046  det->SetPresent();
3047  }
3048  }
3049  if (det->IsPresent()) {
3050  //Test du bon fonctionnement ou non du detecteur
3051  if (!ooodet) {
3052  det->SetDetecting();
3053  }
3054  else {
3055  if (ooodet->FindObject(det->GetName())) {
3056  det->SetDetecting(kFALSE);
3057  if (ndet_ooo) ooo_dets += ",";
3058  ooo_dets += det->GetName();
3059  ndet_ooo += 1;
3060  }
3061  else {
3062  det->SetDetecting();
3063  }
3064  }
3065  }
3066  }
3067 
3068  if (ndet_absent) Info("CheckStatusOfDetectors", "%d detectors absent during run : %s", ndet_absent, absent_dets.Data());
3069  else Info("CheckStatusOfDetectors", "All detectors present during run");
3070  if (ndet_ooo) Info("CheckStatusOfDetectors", "%d detectors malfunctioned during run : %s", ndet_ooo, ooo_dets.Data());
3071  else Info("CheckStatusOfDetectors", "All detectors functioning during run");
3072 }
3073 
3074 
3075 
3077 
3079 {
3080  CheckStatusOfDetectors(dbr, myname);
3081 }
3082 
3083 
3084 
3099 
3101 {
3102  // Read a file containing runlists for each temporarily non-functioning detector.
3103  //
3104  // The file should be in TEnv format like so:
3105  //
3106  //~~~~
3107  // DET_1: 100-122,541-1938
3108  // DET_2,DET_3: 91-765
3109  //~~~~
3110  //
3111  // i.e. more than one detector can be associated with the same runs (comma-separated list of
3112  // detector names) and the list of runs are given using KVNumberList syntax.
3113  //
3114  // The data is added to the database in a table '[name].OoO Detectors' with the name of this array.
3115 
3116  TString fullpath;
3117  if (!db->FindCalibFile("OoODet", fullpath, GetName())) return;
3118 
3119  Info("ReadOoODetectors()", "Reading lists of out-of-order detectors...");
3120  auto fOoODet = db->AddTable(Form("%s.OoO Detectors", GetName()), "Name of out of order detectors");
3121 
3122  KVDBRecord* dbrec = 0;
3123  TEnv env;
3124  TEnvRec* rec = 0;
3125  env.ReadFile(fullpath.Data(), kEnvAll);
3126  TIter it(env.GetTable());
3127 
3128  while ((rec = (TEnvRec*)it.Next())) {
3129  KVString srec(rec->GetName());
3130  KVNumberList nl(rec->GetValue());
3131  if (srec.Contains(",")) {
3132  srec.Begin(",");
3133  while (!srec.End()) {
3134  dbrec = new KVDBRecord(srec.Next(kTRUE), "OoO Detector");
3135  dbrec->AddKey("Runs", "List of Runs");
3136  fOoODet->AddRecord(dbrec);
3137  db->LinkRecordToRunRange(dbrec, nl);
3138  }
3139  }
3140  else {
3141  dbrec = new KVDBRecord(rec->GetName(), "OoO Detector");
3142  dbrec->AddKey("Runs", "List of Runs");
3143  fOoODet->AddRecord(dbrec);
3144  db->LinkRecordToRunRange(dbrec, nl);
3145  }
3146  }
3147 }
3148 
3149 
3164 
3166 {
3167  // Read a file containing runlists for each temporarily absent/dismounted detector.
3168  //
3169  // The file should be in TEnv format like so:
3170  //
3171  //~~~~
3172  // DET_1: 100-122,541-1938
3173  // DET_2,DET_3: 91-765
3174  //~~~~
3175  //
3176  // i.e. more than one detector can be associated with the same runs (comma-separated list of
3177  // detector names) and the list of runs are given using KVNumberList syntax.
3178  //
3179  // The data is added to the database in a table '[name].Absent Detectors' with the name of this array.
3180 
3181  TString fullpath;
3182  if (!db->FindCalibFile("AbsentDet", fullpath, GetName())) return;
3183 
3184  Info("ReadAbsentDetectors()", "Reading lists of absent/dismounted detectors... file=[%s]", fullpath.Data());
3185  auto fAbsDet = db->AddTable(Form("%s.Absent Detectors", GetName()), "Name of out of order detectors");
3186 
3187  KVDBRecord* dbrec = 0;
3188  TEnv env;
3189  TEnvRec* rec = 0;
3190  env.ReadFile(fullpath.Data(), kEnvAll);
3191  TIter it(env.GetTable());
3192 
3193  while ((rec = (TEnvRec*)it.Next())) {
3194  KVString srec(rec->GetName());
3195  KVNumberList nl(rec->GetValue());
3196  if (srec.Contains(",")) {
3197  srec.Begin(",");
3198  while (!srec.End()) {
3199  dbrec = new KVDBRecord(srec.Next(kTRUE), "Absent Detector");
3200  dbrec->AddKey("Runs", "List of Runs");
3201  fAbsDet->AddRecord(dbrec);
3202  db->LinkRecordToRunRange(dbrec, nl);
3203  }
3204  }
3205  else {
3206  dbrec = new KVDBRecord(rec->GetName(), "Absent Detector");
3207  dbrec->AddKey("Runs", "List of Runs");
3208  fAbsDet->AddRecord(dbrec);
3209  db->LinkRecordToRunRange(dbrec, nl);
3210  }
3211  }
3212 }
3213 
3214 
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
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:274
void Print(Option_t *option="") const override
Definition: KVBase.cpp:389
virtual void SetType(const Char_t *str)
Definition: KVBase.h:172
static ValType GetDataSetEnv(const KVString &dataset, const KVString &type, const ValType &defval)
Definition: KVBase.h:304
static Bool_t SearchKVFile(const Char_t *name, TString &fullpath, const Char_t *kvsubdir="")
Definition: KVBase.cpp:517
static TPluginHandler * LoadPlugin(const Char_t *base, const Char_t *uri="0")
Definition: KVBase.cpp:772
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 ,.
Bool_t HasParameter(const Char_t *name) const
TString GetStringParameter(const TString &name) const
Double_t GetParameter(UShort_t i=0) 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:40
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:36
const Char_t * GetDataSetDir() const
Definition: KVDataSet.cpp:746
Bool_t HasCalibIdentInfos() const
Definition: KVDataSet.h:229
KVString GetDataSetEnv(const Char_t *type, const Char_t *defval="") const
Definition: KVDataSet.cpp:784
TString GetFullPathToDataSetFile(const Char_t *filename)
Definition: KVDataSet.cpp:1926
void cd() const
Definition: KVDataSet.cpp:762
Bool_t DataBaseUpdateInProgress() const
Definition: KVDataSet.h:74
static Bool_t FindDataSetFile(const TString &dataset, const Char_t *filename)
Definition: KVDataSet.cpp:1960
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 to describe database of an experiment ,,.
Definition: KVExpDB.h:20
void SetDBType(const TString &s)
Definition: KVExpDB.h:63
TString GetCalibFileName(const Char_t *type) const
Definition: KVExpDB.h:133
const KVNumberList & GetRunList() const
Definition: KVExpDB.h:101
virtual void LinkRecordToRunRange(KVDBRecord *rec, UInt_t first_run, UInt_t last_run)
Definition: KVExpDB.cpp:90
Bool_t FindCalibFile(const Char_t *type, TString &fullpath, const TString &array_name="") const
Definition: KVExpDB.cpp:486
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:19
void SetNumber(UInt_t num) override
Definition: KVGroup.h:38
void Reset(Option_t *opt="")
Definition: KVGroup.cpp:63
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.
KVList * GetGrids()
virtual void LoadGridsForRun(UInt_t)
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:84
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:28
Read MFM format acquisition data.
Description of physical materials used to construct detectors & targets; interface to range tables.
Definition: KVMaterial.h:89
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
KVUniqueNameList * GetIDTelescopeTypes()
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()
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.
static KVMultiDetArray * MakeMultiDetector(const Char_t *dataset_name, Int_t run=-1, TString classname="KVMultiDetArray")
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()
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 Warning(const char *method, const char *msgfmt,...) const
virtual void Error(const char *method, const char *msgfmt,...) 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)