KaliVeda
Toolkit for HIC analysis
KVDetector.cpp
1 #include "KVDetector.h"
2 #include "TROOT.h"
3 #include "KVGroup.h"
4 #include "KVCalibrator.h"
5 #include "TPluginManager.h"
6 #include "TClass.h"
7 /*** geometry ***/
8 #include "TGeoVolume.h"
9 #include "TGeoManager.h"
10 #include "TGeoMatrix.h"
11 #include "TGeoBBox.h"
12 #include "TGeoArb8.h"
13 #include "TGeoTube.h"
14 #include "KVCalibratedSignal.h"
15 #include "KVDetectorSignalTrace.h"
16 #include "KVDetectorSignalExpression.h"
17 #include "KVZDependentCalibratedSignal.h"
18 #include "KVMaterialStack.h"
19 #include <TGeoPhysicalNode.h>
20 #include <TGraph.h>
21 
22 using namespace std;
23 
25 
27 
28 
31 
33 {
34  //default initialisations
35  fCalibrators = nullptr;
36  fParticles.SetCleanup();
37  fActiveLayer = nullptr;
38  fIdentP = fUnidentP = 0;
39  fELossF = fEResF = fRangeF = nullptr;
40  fEResforEinc = -1.;
41  fSimMode = kFALSE;
42  fPresent = kTRUE;
43  fDetecting = kTRUE;
44  fParentStrucList.SetCleanup();
45  fSingleLayer = kTRUE;
46  fNode.SetDetector(this);
47  // detector owns any signals which are defined for it
48  fDetSignals.SetOwner();
49  // adding a new signal with the same name as an existing one
50  // will delete the existing signal and replace it
51  fDetSignals.ReplaceObjects();
52 }
53 
54 
55 
58 
60 {
61 //default ctor
62  init();
63  fDetCounter++;
64  SetName(Form("Det_%d", fDetCounter));
65 }
66 
67 
68 
71 
73  const Float_t thick): KVMaterial()
74 {
75  // Create a new detector of a given material and thickness in centimetres (default value = 0.0)
76 
77  init();
78  SetType("DET");
79  fDetCounter++;
80  SetName(Form("Det_%d", fDetCounter));
81  AddAbsorber(new KVMaterial(type, thick));
82 }
83 
84 
85 
90 
91 KVDetector::KVDetector(const Char_t* gas, const Double_t thick, const Double_t pressure, const Double_t temperature)
92 {
93  // Create gaseous dteector with given type, linear thickness in cm, pressure in Torr, and temperature in degrees C (default value 19°C).
94  //
95  // \note this just defines some gas, the 'detector' has no windows!
96 
97  init();
98  SetType("DET");
99  fDetCounter++;
100  SetName(Form("Det_%d", fDetCounter));
101  AddAbsorber(new KVMaterial(gas, thick, pressure, temperature));
102 }
103 
104 
105 
108 
110 {
111 //copy ctor
112  init();
113  fDetCounter++;
114  SetName(Form("Det_%d", fDetCounter));
115 #if ROOT_VERSION_CODE >= ROOT_VERSION(3,4,0)
116  obj.Copy(*this);
117 #else
118  ((KVDetector&) obj).Copy(*this);
119 #endif
120 }
121 
122 
123 #if ROOT_VERSION_CODE >= ROOT_VERSION(3,4,0)
124 
129 
130 void KVDetector::Copy(TObject& obj) const
131 #else
132 void KVDetector::Copy(TObject& obj)
133 #endif
134 {
135  //copy 'this' to 'obj'
136  //The structure of the detector is copied, with new cloned objects for
137  //each absorber layer. The active layer is set in the new detector.
138 
139  TIter next(&fAbsorbers);
140  KVMaterial* mat;
141  while ((mat = (KVMaterial*) next())) {
142  ((KVDetector&) obj).AddAbsorber((KVMaterial*) mat->Clone());
143  }
144  //set active layer
145  Int_t in_actif = fAbsorbers.IndexOf(fActiveLayer);
146  ((KVDetector&) obj).SetActiveLayer(((KVDetector&)obj).GetAbsorber(in_actif));
147 }
148 
149 
150 
151 
153 
154 KVDetector::~KVDetector()
155 {
160 }
161 
162 
163 
168 
170 {
171  // Set material of active layer.
172  // If no absorbers have been added to the detector, create and add
173  // one (active layer by default)
174 
175  if (!GetActiveLayer())
177  else
179 }
180 
181 
182 
194 
196 {
197  //Calculate the energy loss of a charged particle traversing the detector,
198  //the particle is slowed down, it is added to the list of all particles hitting the
199  //detector. The apparent energy loss of the particle in the active layer of the
200  //detector is set.
201  //Do nothing if particle has zero (or -ve) energy.
202  //
203  //If the optional argument 'norm' is given, it is supposed to be a vector
204  //normal to the detector, oriented from the origin towards the detector.
205  //In this case the effective thicknesses of the detector's absorbers 'seen' by the particle
206  //depending on its direction of motion is used for the calculation.
207 
208  if (kvp->GetKE() <= 0.)
209  return;
210 
211  AddHit(kvp); //add nucleus to list of particles hitting detector in the event
212  //set flag to say that particle has been slowed down
213  kvp->SetIsDetected();
214  //If this is the first absorber that the particle crosses, we set a "reminder" of its
215  //initial energy
216  if (!kvp->GetPInitial())
217  kvp->SetE0();
218 
219  std::vector<Double_t> thickness;
220  if (norm) {
221  // modify thicknesses of all layers according to orientation,
222  // and store original thicknesses in array
223  TVector3 p = kvp->GetMomentum();
224  thickness.reserve(fAbsorbers.GetEntries());
225  KVMaterial* abs;
226  int i = 0;
227  TIter next(&fAbsorbers);
228  while ((abs = (KVMaterial*) next())) {
229  thickness[i++] = abs->GetThickness();
230  Double_t T = abs->GetEffectiveThickness((*norm), p);
231  abs->SetThickness(T);
232  }
233  }
234  Double_t eloss = GetTotalDeltaE(kvp->GetZ(), kvp->GetA(), kvp->GetEnergy());
235  Double_t dE = GetDeltaE(kvp->GetZ(), kvp->GetA(), kvp->GetEnergy());
236  if (norm) {
237  // reset thicknesses of absorbers
238  KVMaterial* abs;
239  int i = 0;
240  TIter next(&fAbsorbers);
241  while ((abs = (KVMaterial*) next())) {
242  abs->SetThickness(thickness[i++]);
243  }
244  }
245  Double_t epart = kvp->GetEnergy() - eloss;
246  if (epart < 1e-3) {
247  //printf("%s, pb d arrondi on met l energie de la particule a 0\n",GetName());
248  epart = 0.0;
249  }
250  kvp->SetEnergy(epart);
251  Double_t eloss_old = GetEnergyLoss();
252  SetEnergyLoss(eloss_old + dE);
253 }
254 
255 
256 
257 
267 
269 {
270  //Calculate the total energy loss of a charged particle traversing the detector.
271  //This does not affect the "stored" energy loss value of the detector,
272  //nor its ACQData, nor the energy of the particle.
273  //
274  //If the optional argument 'norm' is given, it is supposed to be a vector
275  //normal to the detector, oriented from the origin towards the detector.
276  //In this case the effective thicknesses of the detector's absorbers 'seen' by the particle
277  //depending on its direction of motion is used for the calculation.
278 
279  std::vector<Double_t> thickness;
280  if (norm) {
281  // modify thicknesses of all layers according to orientation,
282  // and store original thicknesses in array
283  TVector3 p = kvp->GetMomentum();
284  thickness.reserve(fAbsorbers.GetEntries());
285  KVMaterial* abs;
286  int i = 0;
287  TIter next(&fAbsorbers);
288  while ((abs = (KVMaterial*) next())) {
289  thickness[i++] = abs->GetThickness();
290  Double_t T = abs->GetEffectiveThickness((*norm), p);
291  abs->SetThickness(T);
292  }
293  }
294  Double_t eloss = GetTotalDeltaE(kvp->GetZ(), kvp->GetA(), kvp->GetEnergy());
295  if (norm) {
296  // reset thicknesses of absorbers
297  KVMaterial* abs;
298  int i = 0;
299  TIter next(&fAbsorbers);
300  while ((abs = (KVMaterial*) next())) {
301  abs->SetThickness(thickness[i++]);
302  }
303  }
304  return eloss;
305 }
306 
307 
308 
309 
319 
321 {
322  //Calculate the energy of particle 'kvn' before its passage through the detector,
323  //based on the current kinetic energy, Z & A of nucleus 'kvn', supposed to be
324  //after passing through the detector.
325  //
326  //If the optional argument 'norm' is given, it is supposed to be a vector
327  //normal to the detector, oriented from the origin towards the detector.
328  //In this case the effective thicknesses of the detector's absorbers 'seen' by the particle
329  //depending on its direction of motion is used for the calculation.
330 
331  Double_t Einc = 0.0;
332  //make 'clone' of particle
333  KVNucleus clone_part(kvp->GetZ(), kvp->GetA());
334  clone_part.SetMomentum(kvp->GetMomentum());
335  //composite detector - calculate losses in all layers
336  KVMaterial* abs;
337  TIter next(&fAbsorbers, kIterBackward); //work through list backwards
338  while ((abs = (KVMaterial*) next())) {
339 
340  //calculate energy of particle before current absorber
341  Einc = abs->GetParticleEIncFromERes(&clone_part, norm);
342 
343  //set new energy of particle
344  clone_part.SetKE(Einc);
345  }
346  return Einc;
347 }
348 
349 
350 
351 
355 
356 void KVDetector::Print(Option_t* opt) const
357 {
358  //Print info on this detector
359  //if option="data" the energy loss and raw data are displayed
360 
361  if (!strcmp(opt, "data")) {
362  cout << Form("%10s -- ", GetName());
363  // print raw data signals
364  KVDetectorSignal* sig;
366  while ((sig = (KVDetectorSignal*)it())) {
367  if (sig->IsRaw()) {
368  std::cout << sig->GetName() << "=" << sig->GetValue() << " | ";
369  }
370  }
371  // print calibrated data
372  it.Reset();
373  while ((sig = (KVDetectorSignal*)it())) {
374  if (!sig->IsRaw() && !sig->GetValueNeedsExtraParameters()) {
375  std::cout << sig->GetName() << "=" << sig->GetValue() << " | ";
376  }
377  }
378  cout << " ";
380  cout << "(Belongs to an unidentified particle)";
381  cout << endl;
382  }
383  else if (!strcmp(opt, "all")) {
384  //Give full details of detector
385  //
386  TString option(" ");
387  cout << option << ClassName() << " : " << ((KVDetector*) this)->
388  GetName() << endl;
389  //composite detector - print all layers
390  KVMaterial* abs;
391  TIter next(&fAbsorbers);
392  while ((abs = (KVMaterial*) next())) {
393  if (GetActiveLayer() == abs)
394  cout << " ### ACTIVE LAYER ### " << endl;
395  cout << option << option;
396  abs->Print();
397  if (GetActiveLayer() == abs)
398  cout << " #################### " << endl;
399  }
400  if (fParticles.GetEntries()) {
401  cout << option << " --- Detected particles:" << endl;
402  fParticles.Print();
403  }
404  }
405  else {
406  //just print name
407  cout << ClassName() << ": " << ((KVDetector*) this)->
408  GetName() << endl;
409  }
410 }
411 
412 
413 
438 
440 {
441  // Associate a calibration with this detector.
442  //
443  // This will add a new signal to the list of the detector's signals.
444  //
445  // Also sets calibrator's name to `[detname]_[caltype]` where `caltype` is the type of the KVCalibration object.
446  //
447  // \param[in] cal pointer to KVCalibrator object (must be on heap, i.e. created with new: detector handles deletion)
448  // \param[in] opts
449  // \parblock
450  // can be used to pass any extra parameters/options needed by the calibrator.
451  // For example, if it contains a parameter `ZRange`:
452  //
453  //~~~~~~~~~~~~~~~
454  // ZRange=1-10
455  //~~~~~~~~~~~~~~~
456  //
457  // then the calibrator will be handled by a KVZDependentCalibratedSignal (handles several
458  // calibrators which provide the same output signal, each one is used for a specific range
459  // of atomic numbers)
460  // \endparblock
461  //
462  // \returns kFALSE in case of problems (non-existent input signal for calibrator, output signal not defined for calibrator),
463  // otherwise kTRUE
464 
465  if (!cal) return kFALSE;
466  // look for input signal
468  // input signal must exist
469  if (!in) {
470  Error("AddCalibrator", "Detector:%s has no signal %s, required as input by calibrator %s to provide output %s",
471  GetName(), cal->GetInputSignalType().Data(), cal->GetType(), cal->GetOutputSignalType().Data());
472  return kFALSE;
473  }
474  if (!fCalibrators)
475  fCalibrators = new KVList();
476 
477  cal->SetDetector(this);
478  fCalibrators->Add(cal);
479  cal->SetName(Form("%s_%s", GetName(), cal->GetType()));
480 
481  if (cal->GetOutputSignalType() == "") {
482  Warning("AddCalibrator", "%s : output signal not defined for calibrator %s. No output signal created.",
483  GetName(), cal->GetType());
484  }
485  else {
486  KVDetectorSignal* new_cal_sig(nullptr);
487  if (opts.HasParameter("ZRange")) {
488  // If 'ZRange' parameter is given we need to find the KVZDependentCalibratedSignal
489  // and add a new signal to it.
491  if (!sig) {
492  new_cal_sig = sig = new KVZDependentCalibratedSignal(in, cal->GetOutputSignalType());
493  }
494  dynamic_cast<KVZDependentCalibratedSignal*>(sig)->AddSignal(new KVCalibratedSignal(in, cal), opts.GetStringValue("ZRange"));
495  }
496  else {
497  new_cal_sig = new KVCalibratedSignal(in, cal);
498  }
499  if (new_cal_sig) fDetSignals.Add(new_cal_sig);
500  }
501  return kTRUE;
502 }
503 
504 
505 
516 
518 {
519  // Replace calibrator of given type with the given calibrator object
520  // The calibrator object should not be shared with any other detectors: it now belongs
521  // to this detector, which will delete it when necessary.
522  // If an exising calibrator with the same type is already defined, it will be
523  // deleted and removed from the detector's calibrator list
524  //
525  // Returns kFALSE in case of problems.
526  //
527  // The (optional) KVNameValueList argument can be used to pass any extra parameters/options.
528 
529  KVCalibrator* old_cal = GetCalibrator(type);
530  if (old_cal) {
531  fCalibrators->Remove(old_cal);
533  delete old_cal;
534  }
535  return AddCalibrator(cal, opts);
536 }
537 
538 
539 
544 
546 {
547  // A detector is considered to be calibrated if it has
548  // a signal "Energy" available and if depending on the supplied parameters
549  // this signal can be calculated
550 
551  KVCalibratedSignal* e_sig = dynamic_cast<KVCalibratedSignal*>(GetDetectorSignal("Energy"));
552  return (e_sig && e_sig->IsAvailableFor(params) && IsOK());
553 }
554 
555 
556 
557 
564 
566 {
567  //Set energy loss(es) etc. to zero
568  //
569  //If opt="N":
570  // + we do not reset acquisition parameters/raw detector signals
571  // + in SimMode we do not reset energy losses (this is called before reconstruction happens)
572 
573  bool option_N = (strncmp(opt, "N", 1) == 0);
574  bool sim_mode_option_N = option_N && IsSimMode();
575 
576  if (!sim_mode_option_N) KVMaterial::Clear(opt);
578  fIdentP = fUnidentP = 0;
581  if (!option_N) {
583  KVDetectorSignal* ds;
584  while ((ds = (KVDetectorSignal*)it())) {
585  if (ds->IsRaw() && !ds->IsExpression()) ds->Reset();
586  }
587  }
588  ClearHits();
589  if (!sim_mode_option_N) {
590  //reset all layers in detector
591  KVMaterial* mat;
592  TIter next(&fAbsorbers);
593  while ((mat = (KVMaterial*) next())) {
594  mat->Clear();
595  }
596  }
597  fEResforEinc = -1.;
598 }
599 
600 
601 
606 
608 {
609  // Add a layer of absorber material to the detector
610  // By default, the first layer added is set as the "Active" layer.
611  // Call SetActiveLayer to change this.
612  fAbsorbers.Add(mat);
613  if (!TestBit(kActiveSet))
614  SetActiveLayer(mat);
615  if (fAbsorbers.GetSize() > 1) fSingleLayer = kFALSE;
616 }
617 
618 
619 
622 
624 {
625  //Returns pointer to the i-th absorber in the detector (i=0 first absorber, i=1 second, etc.)
626 
627  if (!fAbsorbers.GetEntries()) {
628  Error("GetAbsorber", "No absorbers defined for detector");
629  return nullptr;
630  }
631  if (i < 0 || i >= fAbsorbers.GetEntries()) {
632  Error("GetAbsorber", "i=%d but valid values are 0-%d [number of absorbers in list]", i,
634  return nullptr;
635  }
636  return (KVMaterial*) fAbsorbers.At(i);
637 }
638 
639 
640 
645 
647 {
648  // Completely reset the KVDetector as if it had just been created by a call to the default constructor:
649  // + removes all absorber layers
650  // + removes all calibrators/detector signals
651 
652  Clear();
655  fAbsorbers.Clear();
659  fCalibrators = nullptr;
660  fActiveLayer = nullptr;
662  fIdentP = fUnidentP = 0;
663  fELossF = fEResF = fRangeF = nullptr;
664  fEResforEinc = -1.;
665  fSimMode = kFALSE;
666  fPresent = kTRUE;
667  fDetecting = kTRUE;
669 }
670 
671 
672 
676 
678 {
679  // Used when a calibrator object is removed or replaced
680  // We remove and delete the corresponding output signal from the list of detector signals
681 
682  if (K->GetOutputSignalType() != "") {
683  KVDetectorSignal* ds = GetDetectorSignal(K->GetOutputSignalType());
684  if (ds) {
685  fDetSignals.Remove(ds);
686  delete ds;
687  }
688  }
689 }
690 
691 
692 
698 
700 {
701  // Removes all calibrations associated to this detector: in other words, we delete all
702  // the KVCalibrator objects in list fCalibrators.
703  //
704  // We also destroy all signals provided by these calibrators
705 
706  if (fCalibrators) {
707  KVCalibrator* K;
708  TIter it(fCalibrators);
709  while ((K = (KVCalibrator*)it())) {
711  }
712  fCalibrators->Delete();
713  }
714 }
715 
716 
717 
718 
748 
750 {
751  // Returns the total energy loss in the detector for a given nucleus
752  // including inactive absorber layers.
753  // e = energy loss in active layer (if not given, we use current value)
754  // transmission = kTRUE (default): the particle is assumed to emerge with
755  // a non-zero residual energy Eres after the detector.
756  // = kFALSE: the particle is assumed to stop in the detector.
757  //
758  // WARNING: if transmission=kTRUE, and if the residual energy after the detector
759  // is known (i.e. measured in a detector placed after this one), you should
760  // first call
761  // SetEResAfterDetector(Eres);
762  // before calling this method. Otherwise, especially for heavy ions, the
763  // correction may be false for particles which are just above the punch-through energy.
764  //
765  // WARNING 2: if measured energy loss in detector active layer is greater than
766  // maximum possible theoretical value for given nucleus' Z & A, this may be because
767  // the A was not measured but calculated from Z and hence could be false, or perhaps
768  // there was an (undetected) pile-up of two or more particles in the detector.
769  // In this case we return the corrected energy
770  // corresponding to the maximum theoretical energy loss in the active layer
771  // and we add the following parameters to the particle (in nuc->GetParameters()):
772  //
773  // GetCorrectedEnergy.Warning = 1
774  // GetCorrectedEnergy.Detector = [name]
775  // GetCorrectedEnergy.MeasuredDE = [value]
776  // GetCorrectedEnergy.MaxDE = [value]
777  // GetCorrectedEnergy.Transmission = 0 or 1
778  // GetCorrectedEnergy.ERES = [value]
779 
780  Int_t z = nuc->GetZ();
781  Int_t a = nuc->GetA();
782 
783  if (e < 0.) e = GetEnergy();
784  if (e <= 0) {
786  return 0;
787  }
788 
789  // check if calling this function for a previous detector in the particle's trajectory
790  // led to using an effective incident angle for the particle: if so, use it here
791  if (nuc->GetParameters()->HasDoubleParameter("GetCorrectedEnergy.IncidentAngle")) {
792  KVMaterialStack stack(this);
793  stack.SetIncidenceAngle(nuc->GetParameters()->GetDoubleValue("GetCorrectedEnergy.IncidentAngle"));
794  // check that apparent energy loss in detector is compatible with a & z
795  Double_t maxDE = stack.GetMaxDeltaE(z, a);
796  if (e > maxDE) {
797  auto inc_angle = stack.GetMinimumIncidentAngleForDEMax(z, a, e);
798  stack.SetIncidenceAngle(inc_angle);
799  nuc->GetParameters()->SetValue(Form("%s.GetCorrectedEnergy.Warning", GetName()), 1);
800  nuc->GetParameters()->SetValue(Form("%s.GetCorrectedEnergy.MeasuredDE", GetName()), e);
801  nuc->GetParameters()->SetValue(Form("%s.GetCorrectedEnergy.MaxDE", GetName()), maxDE);
802  nuc->GetParameters()->SetValue(Form("%s.GetCorrectedEnergy.Transmission", GetName()), (Int_t)transmission);
803  nuc->GetParameters()->SetValue(Form("%s.GetCorrectedEnergy.ERES", GetName()), GetEResAfterDetector());
804  nuc->GetParameters()->SetValue(Form("%s.GetCorrectedEnergy.IncidentAngle", GetName()), inc_angle);
805  nuc->GetParameters()->SetValue("GetCorrectedEnergy.IncidentAngle", inc_angle);
806  }
807  return get_corrected_energy(&stack, nuc, e, transmission);
808  }
809  // check that apparent energy loss in detector is compatible with a & z
810  Double_t maxDE = GetMaxDeltaE(z, a);
811  if (e > maxDE) {
812  KVMaterialStack stack(this);
813  auto inc_angle = stack.GetMinimumIncidentAngleForDEMax(z, a, e);
814  stack.SetIncidenceAngle(inc_angle);
815  nuc->GetParameters()->SetValue(Form("%s.GetCorrectedEnergy.Warning", GetName()), 1);
816  nuc->GetParameters()->SetValue(Form("%s.GetCorrectedEnergy.MeasuredDE", GetName()), e);
817  nuc->GetParameters()->SetValue(Form("%s.GetCorrectedEnergy.MaxDE", GetName()), maxDE);
818  nuc->GetParameters()->SetValue(Form("%s.GetCorrectedEnergy.Transmission", GetName()), (Int_t)transmission);
819  nuc->GetParameters()->SetValue(Form("%s.GetCorrectedEnergy.ERES", GetName()), GetEResAfterDetector());
820  nuc->GetParameters()->SetValue(Form("%s.GetCorrectedEnergy.IncidentAngle", GetName()), inc_angle);
821  nuc->GetParameters()->SetValue("GetCorrectedEnergy.IncidentAngle", inc_angle);
822  return get_corrected_energy(&stack, nuc, e, transmission);
823  }
824  return get_corrected_energy(this, nuc, e, transmission);
825 }
826 
827 
828 
829 
847 
849 {
850  //For particles which stop in the first stage of an identification telescope,
851  //we can at least estimate a minimum Z value based on the energy lost in this
852  //detector.
853  //
854  //This is based on the KVMaterial::GetMaxDeltaE method, giving the maximum
855  //energy loss in the active layer of the detector for a given nucleus (A,Z).
856  //
857  //The "Zmin" is the Z of the particle which gives a maximum energy loss just greater
858  //than that measured in the detector. Particles with Z<Zmin could not lose as much
859  //energy and so are excluded.
860  //
861  //If ELOSS is not given, we use the current value of GetEnergy()
862  //Use 'mass_formula' to change the formula used to calculate the A of the nucleus
863  //from its Z. Default is valley of stability value. (see KVNucleus::GetAFromZ).
864  //
865  //If the value of ELOSS or GetEnergy() is <=0 we return Zmin=0
866 
867  ELOSS = (ELOSS < 0. ? GetEnergy() : ELOSS);
868  if (ELOSS <= 0) return 0;
869 
870  UInt_t z = 40;
871  UInt_t zmin, zmax;
872  zmin = 1;
873  zmax = 92;
874  KVNucleus particle(zmin);
875 
876  if (mass_formula > -1)
877  particle.SetMassFormula((UChar_t)mass_formula);
878 
879  if (GetMaxDeltaE(particle.GetZ(), particle.GetA()) > ELOSS)
880  return zmin; // Z=1 solution
881 
882  do {
883 
884  particle.SetZ(z);
885 
886  auto difference = GetMaxDeltaE(z, particle.GetA()) - ELOSS;
887  //if difference < 0 the z is too small
888  if (difference < 0.0) {
889 
890  zmin = z;
891  z += (UInt_t)((zmax - z) / 2 + 0.5);
892 
893  }
894  else {
895 
896  zmax = z;
897  z -= (UInt_t)((z - zmin) / 2 + 0.5);
898 
899  }
900  }
901  while (zmax > (zmin + 1));
902 
903  return zmax;
904 }
905 
906 
907 
908 
917 
919 {
920  // Calculates energy loss (in MeV) in active layer of detector, taking into account preceding layers
921  //
922  // Arguments are:
923  // x[0] is incident energy in MeV
924  // Parameters are:
925  // par[0] Z of ion
926  // par[1] A of ion
927 
928  Double_t e = x[0];
929  TIter next(&fAbsorbers);
930  KVMaterial* mat;
931  mat = (KVMaterial*)next();
932  while (fActiveLayer != mat) {
933  // calculate energy losses in absorbers before active layer
934  e = mat->GetERes(par[0], par[1], e); //residual energy after layer
935  if (e <= 0.)
936  return 0.; // return 0 if particle stops in layers before active layer
937  mat = (KVMaterial*)next();
938  }
939  //calculate energy loss in active layer
940  return mat->GetDeltaE(par[0], par[1], e);
941 }
942 
943 
944 
945 
955 
957 {
958  // Calculates range (in centimetres) of ions in detector as a function of incident energy (in MeV),
959  // taking into account all layers of the detector.
960  //
961  // Arguments are:
962  // x[0] = incident energy in MeV
963  // Parameters are:
964  // par[0] = Z of ion
965  // par[1] = A of ion
966 
967  Double_t Einc = x[0];
968  Int_t Z = (Int_t)par[0];
969  Int_t A = (Int_t)par[1];
970  Double_t range = 0.;
971  TIter next(&fAbsorbers);
972  KVMaterial* mat = (KVMaterial*)next();
973  if (!mat) return 0.;
974  do {
975  // range in this layer
976  Double_t this_range = mat->GetLinearRange(Z, A, Einc);
977  KVMaterial* next_mat = (KVMaterial*)next();
978  if (this_range > mat->GetThickness()) {
979  // particle traverses layer.
980  if (next_mat)
981  range += mat->GetThickness();
982  else // if this is last layer, the range continues to increase beyond size of detector
983  range += this_range;
984  // calculate incident energy for next layer (i.e. residual energy after this layer)
985  Einc = mat->GetERes(Z, A, Einc);
986  }
987  else {
988  // particle stops in layer
989  range += this_range;
990  return range;
991  }
992  mat = next_mat;
993  }
994  while (mat);
995  // particle traverses detector
996  return range;
997 }
998 
999 
1000 
1001 
1011 
1013 {
1014  // Calculates residual energy (in MeV) of particle after traversing all layers of detector.
1015  // Returned value is -1000 if particle stops in one of the layers of the detector.
1016  //
1017  // Arguments are:
1018  // x[0] is incident energy in MeV
1019  // Parameters are:
1020  // par[0] Z of ion
1021  // par[1] A of ion
1022 
1023  Double_t e = x[0];
1024  TIter next(&fAbsorbers);
1025  KVMaterial* mat;
1026  while ((mat = (KVMaterial*)next())) {
1027  Double_t eres = mat->GetERes(par[0], par[1], e); //residual energy after layer
1028  if (eres <= 0.)
1029  return -1000.; // return -1000 if particle stops in layers before active layer
1030  e = eres;
1031  }
1032  return e;
1033 }
1034 
1035 
1036 
1037 
1052 
1054 {
1055  //Static function which will create an instance of the KVDetector-derived class
1056  //corresponding to 'name'
1057  //These are defined as 'Plugin' objects in the file $KVROOT/KVFiles/.kvrootrc :
1058  // [name_of_dataset].detector_type
1059  // detector_type
1060  // To use the dataset-dependent plugin, call this method with
1061  // name = "[name_of_dataset].detector_type"
1062  // If not, the default plugin will be used
1063  //first we check if there is a special plugin for the DataSet
1064  //if not we take the default one
1065  //
1066  //'thickness' is passed as argument to the constructor for the detector plugin
1067 
1068  //check and load plugin library
1069  TPluginHandler* ph = nullptr;
1070  KVString nom(name);
1071  if (!nom.Contains(".") && !(ph = LoadPlugin("KVDetector", name))) return nullptr;
1072  if (nom.Contains(".")) {
1073  // name format like [dataset].[det_type]
1074  // in case dataset name contains "." we parse string to find detector type: assumed after the last "."
1075  nom.RBegin(".");
1076  KVString det_type = nom.RNext();
1077  if (!(ph = LoadPlugin("KVDetector", name))) {
1078  if (!(ph = LoadPlugin("KVDetector", det_type))) {
1079  return nullptr;
1080  }
1081  }
1082  }
1083 
1084  //execute constructor/macro for detector
1085  return (KVDetector*) ph->ExecPlugin(1, thickness);
1086 }
1087 
1088 
1089 
1092 
1094 {
1095  // Return surface area of first layer of detector in cm2.
1096 
1097  if (!GetEntranceWindow().GetShape()->InheritsFrom("TGeoArb8")) {
1098  // simple shape area
1099  return GetEntranceWindow().GetShape()->GetFacetArea(1);
1100  }
1101  // Monte Carlo calculation for TGeoArb8 shapes
1102  return GetEntranceWindow().GetSurfaceArea();
1103 }
1104 
1105 
1106 
1111 
1113 {
1114  // Return pointer toTF1 giving residual energy after detector as function of incident energy,
1115  // for a given nucleus (Z,A).
1116  // The TF1::fNpx parameter is taken from environment variable KVDetector.ResidualEnergy.Npx
1117 
1118  if (!fEResF) {
1119  fEResF = new TF1(Form("KVDetector:%s:ERes", GetName()), this, &KVDetector::EResDet,
1120  0., 1.e+04, 2, "KVDetector", "EResDet");
1121  fEResF->SetNpx(gEnv->GetValue("KVDetector.ResidualEnergy.Npx", 20));
1122  }
1124  fEResF->SetRange(0., GetSmallestEmaxValid(Z, A));
1125  fEResF->SetTitle(Form("Residual energy [MeV] after detector %s for Z=%d A=%d", GetName(), Z, A));
1126 
1127  return fEResF;
1128 }
1129 
1130 
1131 
1136 
1138 {
1139  // Return pointer toTF1 giving range (in centimetres) in detector as function of incident energy,
1140  // for a given nucleus (Z,A).
1141  // The TF1::fNpx parameter is taken from environment variable KVDetector.Range.Npx
1142 
1143  if (!fRangeF) {
1144  fRangeF = new TF1(Form("KVDetector:%s:Range", GetName()), this, &KVDetector::RangeDet,
1145  0., 1.e+04, 2, "KVDetector", "RangeDet");
1146  fRangeF->SetNpx(gEnv->GetValue("KVDetector.Range.Npx", 20));
1147  }
1150  fRangeF->SetTitle(Form("Range [cm] in detector %s for Z=%d A=%d", GetName(), Z, A));
1151 
1152  return fRangeF;
1153 }
1154 
1155 
1156 
1161 
1163 {
1164  // Return pointer to TF1 giving energy loss in active layer of detector as function of incident energy,
1165  // for a given nucleus (Z,A).
1166  // The TF1::fNpx parameter is taken from environment variable KVDetector.EnergyLoss.Npx
1167 
1168  if (!fELossF) {
1169  fELossF = new TF1(Form("KVDetector:%s:ELossActive", GetName()), this, &KVDetector::ELossActive,
1170  0., 1.e+04, 2, "KVDetector", "ELossActive");
1171  fELossF->SetNpx(gEnv->GetValue("KVDetector.EnergyLoss.Npx", 20));
1172  }
1175  fELossF->SetTitle(Form("Energy loss [MeV] in detector %s for Z=%d A=%d", GetName(), Z, A));
1176  return fELossF;
1177 }
1178 
1179 
1180 
1185 
1187 {
1188  // Overrides KVMaterial::GetEIncOfMaxDeltaE
1189  // Returns incident energy corresponding to maximum energy loss in the
1190  // active layer of the detector, for a given nucleus.
1191 
1192  return GetELossFunction(Z, A)->GetMaximumX();
1193 }
1194 
1195 
1196 
1201 
1203 {
1204  // Overrides KVMaterial::GetMaxDeltaE
1205  // Returns maximum energy loss in the
1206  // active layer of the detector, for a given nucleus.
1207 
1208  return GetELossFunction(Z, A)->GetMaximum();
1209 }
1210 
1211 
1212 
1217 
1219 {
1220  // Overrides KVMaterial::GetDeltaE
1221  // Returns energy loss of given nucleus in the active layer of the detector.
1222 
1223  // optimization for single-layer detectors
1224  if (fSingleLayer) {
1225  return fActiveLayer->GetDeltaE(Z, A, Einc);
1226  }
1227  return GetELossFunction(Z, A)->Eval(Einc);
1228 }
1229 
1230 
1231 
1235 
1237 {
1238  // Returns calculated total energy loss of ion in ALL layers of the detector.
1239  // This is just (Einc - GetERes(Z,A,Einc))
1240 
1241  return Einc - GetERes(Z, A, Einc);
1242 }
1243 
1244 
1245 
1250 
1252 {
1253  // Overrides KVMaterial::GetERes
1254  // Returns residual energy of given nucleus after the detector.
1255  // Returns 0 if Einc<=0
1256 
1257  if (Einc <= 0.) return 0.;
1258  Double_t eres = GetEResFunction(Z, A)->Eval(Einc);
1259  // Eres function returns -1000 when particle stops in detector,
1260  // in order for function inversion (GetEIncFromEres) to work
1261  if (eres < 0.) eres = 0.;
1262  return eres;
1263 }
1264 
1265 
1266 
1289 
1291 {
1292  // Overrides KVMaterial::GetIncidentEnergy
1293  // Returns incident energy corresponding to energy loss delta_e in active layer of detector for a given nucleus.
1294  // If delta_e is not given, the current energy loss in the active layer is used.
1295  //
1296  // By default the solution corresponding to the highest incident energy is returned
1297  // This is the solution found for Einc greater than the maximum of the dE(Einc) curve.
1298  // If you want the low energy solution set SolType = KVIonRangeTable::kEmin.
1299  //
1300  // WARNING: calculating the incident energy of a particle using only the dE in a detector
1301  // is ambiguous, as in general (and especially for very heavy ions) the maximum of the dE
1302  // curve occurs for Einc greater than the punch-through energy, therefore it is not always
1303  // true to assume that if the particle does not stop in the detector the required solution
1304  // is that for type=KVIonRangeTable::kEmax. For a range of energies between punch-through
1305  // and dE_max, the required solution is still that for type=KVIonRangeTable::kEmin.
1306  // If the residual energy of the particle is unknown, there is no way to know which is the
1307  // correct solution.
1308  //
1309  // WARNING 2
1310  // If the given energy loss in the active layer is greater than the maximum theoretical dE
1311  // for given Z & A, (dE > GetMaxDeltaE(Z,A)) then we return a NEGATIVE incident energy
1312  // corresponding to the maximum, GetEIncOfMaxDeltaE(Z,A)
1313 
1314  if (Z < 1) return 0.;
1315 
1316  Double_t DE = (delta_e > 0 ? delta_e : GetEnergyLoss());
1317 
1318  // If the given energy loss in the active layer is greater than the maximum theoretical dE
1319  // for given Z & A, (dE > GetMaxDeltaE(Z,A)) then we return a NEGATIVE incident energy
1320  // corresponding to the maximum, GetEIncOfMaxDeltaE(Z,A)
1321  if (DE > GetMaxDeltaE(Z, A)) return -GetEIncOfMaxDeltaE(Z, A);
1322 
1323  TF1* dE = GetELossFunction(Z, A);
1324  Double_t e1, e2;
1325  dE->GetRange(e1, e2);
1326  switch (type) {
1327  case kEmin:
1328  e2 = GetEIncOfMaxDeltaE(Z, A);
1329  break;
1330  case kEmax:
1331  e1 = GetEIncOfMaxDeltaE(Z, A);
1332  break;
1333  }
1334  KVBase::GetX_status status;
1335  Double_t EINC = ProtectedGetX(dE, DE, status, e1, e2);
1337 // Warning("GetIncidentEnergy",
1338 // "In %s : Called for Z=%d A=%d with dE=%f sol_type %d",
1339 // GetName(), Z, A, DE, type);
1340 // Warning("GetIncidentEnergy",
1341 // "Max DeltaE in this case is %f MeV",
1342 // GetMaxDeltaE(Z, A));
1343 // Warning("GetIncidentEnergy",
1344 // "Between Einc limits [%f,%f] no solution found",
1345 // e1,e2);
1346 // Warning("GetIncidentEnergy",
1347 // "Returned value is -Einc for %s dE, %f",
1348 // (status>0 ? "maximum" : "minimum"), EINC);
1349  return -EINC;
1350  }
1351  return EINC;
1352 }
1353 
1354 
1355 
1361 
1363 {
1364  // Overrides KVMaterial::GetDeltaEFromERes
1365  //
1366  // Calculate energy loss in active layer of detGetAlignedDetector for nucleus (Z,A)
1367  // having a residual kinetic energy Eres (MeV)
1368 
1369  if (Z < 1 || Eres <= 0.) return 0.;
1370  Double_t Einc = GetIncidentEnergyFromERes(Z, A, Eres);
1371  if (Einc <= 0.) return 0.;
1372  return GetELossFunction(Z, A)->Eval(Einc);
1373 }
1374 
1375 
1376 
1383 
1385 {
1386  // Overrides KVMaterial::GetIncidentEnergyFromERes
1387  //
1388  // Calculate incident energy of nucleus from residual energy.
1389  //
1390  // Returns -1 if Eres is out of defined range of values
1391 
1392  if (Z < 1 || Eres <= 0.) return 0.;
1393  //return GetEResFunction(Z, A)->GetX(Eres);
1394  KVBase::GetX_status status;
1395  Double_t einc = KVBase::ProtectedGetX(GetEResFunction(Z, A), Eres, status);
1396  if (status == KVBase::GetX_status::above_maximum || status == KVBase::GetX_status::below_minimum){// problem with inversion - value out of defined range of function
1397  return -1;
1398  }
1399  return einc;
1400 }
1401 
1402 
1403 
1407 
1409 {
1410  // Returns the smallest maximum energy for which range tables are valid
1411  // for all absorbers in the detector, and given ion (Z,A)
1412 
1413  Double_t maxmin = -1.;
1414  TIter next(&fAbsorbers);
1415  KVMaterial* abs;
1416  while ((abs = (KVMaterial*)next())) {
1417  if (maxmin < 0.) maxmin = abs->GetEmaxValid(Z, A);
1418  else {
1419  if (abs->GetEmaxValid(Z, A) < maxmin) maxmin = abs->GetEmaxValid(Z, A);
1420  }
1421  }
1422  return maxmin;
1423 }
1424 
1425 
1426 
1444 
1446 {
1447  // Create detector from text file in 'TEnv' format.
1448  //
1449  // Example:
1450  // ========
1451  //
1452  // Layer: Gold
1453  // Gold.Material: Au
1454  // Gold.AreaDensity: 200.*KVUnits::ug
1455  // +Layer: Gas1
1456  // Gas1.Material: C3F8
1457  // Gas1.Thickness: 5.*KVUnits::cm
1458  // Gas1.Pressure: 50.*KVUnits::mbar
1459  // Gas1.Active: yes
1460  // +Layer: Si1
1461  // Si1.Material: Si
1462  // Si1.Thickness: 300.*KVUnits::um
1463 
1464  TEnv fEnvFile(envrc);
1465 
1466  KVString layers(fEnvFile.GetValue("Layer", ""));
1467  layers.Begin(" ");
1468  while (!layers.End()) {
1469  KVString layer = layers.Next();
1470  KVString mat = fEnvFile.GetValue(Form("%s.Material", layer.Data()), "");
1471  KVString tS = fEnvFile.GetValue(Form("%s.Thickness", layer.Data()), "");
1472  KVString pS = fEnvFile.GetValue(Form("%s.Pressure", layer.Data()), "");
1473  KVString dS = fEnvFile.GetValue(Form("%s.AreaDensity", layer.Data()), "");
1474  Double_t thick, dens, press;
1475  thick = dens = press = 0.;
1476  KVMaterial* M = 0;
1477  if (pS != "" && tS != "") {
1478  press = (Double_t)gROOT->ProcessLineFast(Form("%s*1.e+12", pS.Data()));
1479  press /= 1.e+12;
1480  thick = (Double_t)gROOT->ProcessLineFast(Form("%s*1.e+12", tS.Data()));
1481  thick /= 1.e+12;
1482  M = new KVMaterial(mat.Data(), thick, press);
1483  }
1484  else if (tS != "") {
1485  thick = (Double_t)gROOT->ProcessLineFast(Form("%s*1.e+12", tS.Data()));
1486  thick /= 1.e+12;
1487  M = new KVMaterial(mat.Data(), thick);
1488  }
1489  else if (dS != "") {
1490  dens = (Double_t)gROOT->ProcessLineFast(Form("%s*1.e+12", dS.Data()));
1491  dens /= 1.e+12;
1492  M = new KVMaterial(dens, mat.Data());
1493  }
1494  if (M) {
1495  AddAbsorber(M);
1496  if (fEnvFile.GetValue(Form("%s.Active", layer.Data()), kFALSE)) SetActiveLayer(M);
1497  }
1498  }
1499 }
1500 
1501 
1502 
1506 
1508 {
1509  // WARNING: SAME AS KVDetector::GetLinearRange
1510  // Only linear range in centimetres is calculated for detectors!
1511  return GetLinearRange(Z, A, Einc);
1512 }
1513 
1514 
1515 
1521 
1523 {
1524  // Returns range of ion in centimetres in this detector,
1525  // taking into account all layers.
1526  // Note that for Einc > punch through energy, this range is no longer correct
1527  // (but still > total thickness of detector).
1528  return GetRangeFunction(Z, A)->Eval(Einc);
1529 }
1530 
1531 
1532 
1536 
1538 {
1539  // Returns energy (in MeV) necessary for ion (Z,A) to punch through all
1540  // layers of this detector
1541 
1542  if (fSingleLayer) {
1543  // Optimize calculation time for single-layer detector
1544  return fActiveLayer->GetPunchThroughEnergy(Z, A);
1545  }
1546  //return GetRangeFunction(Z, A)->GetX(GetTotalThicknessInCM());
1547  KVBase::GetX_status status;
1548  return ProtectedGetX(GetRangeFunction(Z, A),GetTotalThicknessInCM(),status);
1549 }
1550 
1551 
1552 
1553 
1558 
1560 {
1561  // Creates and fills a KVDrawable<TGraph> with the punch through energy in MeV vs. Z for the given detector,
1562  // for Z=1-92. The mass of each nucleus is calculated according to the given mass formula
1563  // (see KVNucleus).
1564 
1565  TGraph* punch = new TGraph(92);
1566  punch->SetName(Form("KVDetpunchthrough_%s_mass%d", GetName(), massform));
1567  punch->SetTitle(Form("Simple Punch-through %s (MeV) (mass formula %d)", GetName(), massform));
1568  KVNucleus nuc;
1569  nuc.SetMassFormula(massform);
1570  for (int Z = 1; Z <= 92; Z++) {
1571  nuc.SetZ(Z);
1572  punch->SetPoint(Z - 1, Z, GetPunchThroughEnergy(nuc.GetZ(), nuc.GetA()));
1573  }
1574  return punch;
1575 }
1576 
1577 
1578 
1583 
1585 {
1586  // Creates and fills a KVDrawable<TGraph> with the punch through energy in MeV/nucleon vs. Z for the given detector,
1587  // for Z=1-92. The mass of each nucleus is calculated according to the given mass formula
1588  // (see KVNucleus).
1589 
1590  TGraph* punch = new TGraph(92);
1591  punch->SetName(Form("KVDetpunchthroughEsurA_%s_mass%d", GetName(), massform));
1592  punch->SetTitle(Form("Simple Punch-through %s (AMeV) (mass formula %d)", GetName(), massform));
1593  KVNucleus nuc;
1594  nuc.SetMassFormula(massform);
1595  for (int Z = 1; Z <= 92; Z++) {
1596  nuc.SetZ(Z);
1597  punch->SetPoint(Z - 1, Z, GetPunchThroughEnergy(nuc.GetZ(), nuc.GetA()) / nuc.GetA());
1598  }
1599  return punch;
1600 }
1601 
1602 
1603 
1604 
1606 
1608 {
1609  return (KVGroup*)GetParentStructure("GROUP");
1610 }
1611 
1612 
1613 
1615 
1617 {
1618  return (GetGroup() ? GetGroup()->GetNumber() : 0);
1619 }
1620 
1621 
1622 
1624 
1626 {
1627  fParentStrucList.Add(elem);
1628 }
1629 
1630 
1631 
1633 
1635 {
1636  fParentStrucList.Remove(elem);
1637 }
1638 
1639 
1640 
1644 
1646 {
1647  // Get parent geometry structure element of given type.
1648  // Give unique name of structure if more than one element of same type is possible.
1649  KVGeoStrucElement* el = 0;
1650  if (strcmp(name, "")) {
1652  el = (KVGeoStrucElement*)strucs->FindObject(name);
1653  delete strucs;
1654  }
1655  else
1657  return el;
1658 }
1659 
1660 
1661 
1664 
1666 {
1667  // Set ROOT geometry global matrix transformation to coordinate frame of active layer volume
1668  SetMatrix(m);
1669 }
1670 
1671 
1672 
1675 
1677 {
1678  // Set ROOT geometry shape of active layer volume
1679  SetShape(s);
1680 }
1681 
1682 
1683 
1686 
1688 {
1689  // Set ROOT geometry global matrix transformation to coordinate frame of entrance window
1691 }
1692 
1693 
1694 
1697 
1699 {
1700  // Set ROOT geometry shape of entrance window
1701  fEWPosition.SetShape(s);
1702 }
1703 
1704 
1705 
1716 
1718 {
1719  // Overrides KVMaterial::SetThickness
1720  //
1721  // If ROOT geometry is defined, we modify the DZ thickness of the volume representing
1722  // this detector in accordance with the new thickness.
1723  //
1724  // This is only implemented for single-layer detectors with the following shapes:
1725  // - TGeoBBox (rectangular box)
1726  // - TGeoTube (tube)
1727  // - TGeoArb8 (arbitrary trapezoid with less than 8 vertices standing on two parallel planes perpendicular to Z axis)
1728 
1729  if (ROOTGeo() && fSingleLayer) {
1732  TGeoBBox* shape = (TGeoBBox*)pn->GetShape();
1733  TGeoShape* newshape = nullptr;
1734  // bad kludge - is there no better way to clone a shape and change its dZ?
1735  if (shape->IsA() == TGeoBBox::Class()) {
1736  newshape = new TGeoBBox(shape->GetDX(), shape->GetDY(), 0.5 * thick);
1737  }
1738  else if (shape->IsA() == TGeoTube::Class()) {
1739  TGeoTube* oldtube = static_cast<TGeoTube*>(shape);
1740  newshape = new TGeoTube(oldtube->GetRmin(), oldtube->GetRmax(), 0.5 * thick);
1741  }
1742  else if (shape->IsA() == TGeoArb8::Class()) {
1743  TGeoArb8* oldtube = static_cast<TGeoArb8*>(shape);
1744  auto vert = oldtube->GetVertices();
1745  newshape = new TGeoArb8(0.5 * thick, vert);
1746  }
1747  else {
1748  Error("SetThickness", "No implementation for %s (%s)", shape->IsA()->GetName(), GetName());
1749  }
1750  if (newshape) {
1751  pn->Align(nullptr, newshape);
1753  }
1754  }
1755  KVMaterial::SetThickness(thick);
1756 }
1757 
1758 
1759 
1765 
1767 {
1768  // Return kTRUE if the two detectors have the same internal structure, i.e.
1769  // - the same number of absorber layers
1770  // - in the same order
1771  // - with the same material & thickness
1772 
1773  int nabs = GetNumberOfAbsorberLayers();
1774  if (other->GetNumberOfAbsorberLayers() != nabs) return kFALSE;
1775  bool same = true;
1776  for (int iabs = 0; iabs < nabs; ++iabs) {
1777  KVMaterial* this_abs = GetAbsorber(iabs);
1778  KVMaterial* other_abs = other->GetAbsorber(iabs);
1779  if (!this_abs->IsType(other_abs->GetType())
1780  || this_abs->GetMass() != other_abs->GetMass()
1781  || this_abs->GetThickness() != other_abs->GetThickness())
1782  same = false;
1783  }
1784  return same;
1785 }
1786 
1787 
1788 
1796 
1798 {
1799  // Add a new signal to the list of detector's signals.
1800  // \param[in] type define the name of the signal to add
1801  // \returns pointer to the new signal object
1802  // \note do not `delete`{.cpp} the signal object: the detector handles deletion
1803  //
1804  // If the type is "SignalTrace", we add a special KVDetectorSignalTrace object
1805 
1806  if(type=="SignalTrace")
1807  return AddDetectorSignalTrace();
1808  auto signal = new KVDetectorSignal(type, this);
1809  fDetSignals.Add(signal);
1810  return signal;
1811 }
1812 
1813 
1814 
1820 
1822 {
1823  // Add a new signal trace (with type "SignalTrace") to the list of detector's signals.
1824  //
1825  // \returns pointer to the new signal object
1826  // \note do not `delete`{.cpp} the signal object: the detector handles deletion
1827 
1828  auto signal = new KVDetectorSignalTrace(this);
1829  fDetSignals.Add(signal);
1830  return signal;
1831 }
1832 
1833 
1834 
1842 
1844 {
1845  // Add a new KVDetectorSignalExpression to this detector
1846  //
1847  // \param[in] type the name/type of the new signal
1848  // \param[in] _expr mathematical expression using any of the known signals of the detector
1849  //
1850  // \note If the expression is not valid, no signal will be created and method returns kFALSE.
1851 
1853  if (!ds->IsValid()) {
1854  delete ds;
1855  ds = nullptr;
1856  }
1857  else
1858  AddDetectorSignal(ds);
1859  return (ds != nullptr);
1860 }
1861 
1862 
1863 
1872 
1874 {
1875  // \param[in] P pressure in [Torr]
1876  //
1877  // For a gaseous detector, set/change the pressure of the active gas layer.
1878  //
1879  // For ROOT geometries, we change the medium/material of the corresponding node in the geometry
1880  // in order to reflect the change in pressure (gases are represented by different media/materials
1881  // for each pressure/temperature) so that it will be taken into account for example when filtering simulated data.
1882 
1884  if (ROOTGeo()) {
1885  // find node in geometry
1888  // use new medium reflecting change in pressure
1890  }
1891 }
1892 
1893 
1894 
1903 
1905 {
1906  // \param[in] T pressure in [deg. C]
1907  //
1908  // For a gaseous detector, set/change the temperature of the active gas layer.
1909  //
1910  // For ROOT geometries, we change the medium/material of the corresponding node in the geometry
1911  // in order to reflect the change in temperature (gases are represented by different media/materials
1912  // for each pressure/temperature) so that it will be taken into account for example when filtering simulated data.
1913 
1915  if (ROOTGeo()) {
1916  // find node in geometry
1919  // use new medium reflecting change in pressure
1921  }
1922 }
1923 
1924 
int Int_t
unsigned int UInt_t
#define SafeDelete(p)
#define e(i)
bool Bool_t
unsigned char UChar_t
char Char_t
float Float_t
constexpr Bool_t kFALSE
double Double_t
constexpr Bool_t kTRUE
const char Option_t
const Bool_t kIterBackward
R__EXTERN TEnv * gEnv
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 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
#define gROOT
char * Form(const char *fmt,...)
virtual const Char_t * GetType() const
Definition: KVBase.h:176
GetX_status
Definition: KVBase.h:320
void Error(const char *method, const char *msgfmt,...) const override
Definition: KVBase.cpp:1709
static Double_t ProtectedGetX(const TF1 *func, Double_t val, GetX_status &status, std::optional< Double_t > xmin={}, std::optional< Double_t > xmax={})
Definition: KVBase.cpp:1623
virtual Bool_t IsType(const Char_t *typ) const
Definition: KVBase.h:184
virtual void SetType(const Char_t *str)
Definition: KVBase.h:172
static TPluginHandler * LoadPlugin(const Char_t *base, const Char_t *uri="0")
Definition: KVBase.cpp:796
void Warning(const char *method, const char *msgfmt,...) const override
Definition: KVBase.cpp:1696
UInt_t GetNumber() const
Definition: KVBase.h:219
Output signal from detector obtained by calibration.
Bool_t IsAvailableFor(const KVNameValueList &params) const override
Base class for all detector calibrations.
Definition: KVCalibrator.h:99
virtual void SetDetector(KVDetector *d)
Definition: KVCalibrator.h:236
TString GetInputSignalType() const
Definition: KVCalibrator.h:228
TString GetOutputSignalType() const
Definition: KVCalibrator.h:232
Signal output from a mathematical combination of other signals.
Detector signal waveform.
Base class for output signal data produced by a detector.
virtual Bool_t IsExpression() const
virtual Bool_t IsRaw() const
virtual Bool_t GetValueNeedsExtraParameters() const
virtual Double_t GetValue(const KVNameValueList &params="") const
virtual void Reset()
Base class for detector geometry description, interface to energy-loss calculations.
Definition: KVDetector.h:160
static KVDetector * MakeDetector(const Char_t *name, Float_t thick)
void SetThickness(Double_t thick) override
virtual Bool_t IsSimMode() const
Definition: KVDetector.h:643
virtual Bool_t IsOK() const
Definition: KVDetector.h:672
KVPosition fEWPosition
position of entrance window i.e. first volume in detector geometry
Definition: KVDetector.h:163
KVMaterial * GetActiveLayer() const override
Definition: KVDetector.h:313
KVUniqueNameList fParentStrucList
list of geometry structures which directly contain this detector
Definition: KVDetector.h:164
Int_t GetNumberOfAbsorberLayers() const
Definition: KVDetector.h:328
void SetMaterial(const Char_t *type) override
Definition: KVDetector.cpp:169
Double_t ELossActive(Double_t *x, Double_t *par)
Definition: KVDetector.cpp:918
virtual KVDrawable< TGraph > DrawPunchThroughEnergyVsZ(Int_t massform=KVNucleus::kBetaMass)
void SetTemperature(Double_t T) override
KVList * fCalibrators
list of associated calibrator objects
Definition: KVDetector.h:258
KVGeoStrucElement * GetParentStructure(const Char_t *type, const Char_t *name="") const
KVGroup * GetGroup() const
void AddDetectorSignal(KVDetectorSignal *ds)
Definition: KVDetector.h:278
Double_t GetRange(Int_t Z, Int_t A, Double_t Einc) override
Double_t GetIncidentEnergyFromERes(Int_t Z, Int_t A, Double_t Eres) override
Bool_t fDetecting
=kTRUE if detector is "detecting", =kFALSE if not
Definition: KVDetector.h:274
Bool_t ReplaceCalibrator(const Char_t *type, KVCalibrator *cal, const KVNameValueList &opts="")
Definition: KVDetector.cpp:517
virtual Double_t GetTotalDeltaE(Int_t Z, Int_t A, Double_t Einc)
virtual Double_t GetEnergy() const
Definition: KVDetector.h:365
virtual TF1 * GetEResFunction(Int_t Z, Int_t A)
static Int_t fDetCounter
Definition: KVDetector.h:166
void AddAbsorber(KVMaterial *)
Definition: KVDetector.cpp:607
KVDetectorSignalTrace * AddDetectorSignalTrace()
virtual Double_t GetEResAfterDetector() const
Definition: KVDetector.h:627
Double_t GetIncidentEnergy(Int_t Z, Int_t A, Double_t delta_e=-1.0, enum SolType type=kEmax) override
Double_t GetPunchThroughEnergy(Int_t Z, Int_t A) override
void AddParentStructure(KVGeoStrucElement *elem)
Double_t GetDeltaE(Int_t Z, Int_t A, Double_t Einc, Double_t=0.) override
Double_t fEResforEinc
used by GetIncidentEnergy & GetCorrectedEnergy
Definition: KVDetector.h:270
void Print(Option_t *option="") const override
Definition: KVDetector.cpp:356
virtual KVDrawable< TGraph > DrawPunchThroughEsurAVsZ(Int_t massform=KVNucleus::kBetaMass)
virtual Int_t FindZmin(Double_t ELOSS=-1., Char_t mass_formula=-1)
Definition: KVDetector.cpp:848
TF1 * fELossF
parametric function dE in active layer vs. incident energy
Definition: KVDetector.h:266
Double_t GetLinearRange(Int_t Z, Int_t A, Double_t Einc) override
const KVPosition & GetEntranceWindow() const
Definition: KVDetector.h:689
KVMaterial * GetAbsorber(Int_t i) const
Returns pointer to the i-th absorber in the detector (i=0 first absorber, i=1 second,...
Definition: KVDetector.cpp:623
@ kIdentifiedParticle
Definition: KVDetector.h:174
@ kUnidentifiedParticle
Definition: KVDetector.h:173
void ClearHits()
Definition: KVDetector.h:438
void SetEnergyLoss(Double_t e) const override
Definition: KVDetector.h:389
void remove_signal_for_calibrator(KVCalibrator *K)
Definition: KVDetector.cpp:677
TF1 * fEResF
parametric function Eres residual energy after all layers of detector
Definition: KVDetector.h:267
Bool_t HasSameStructureAs(const KVDetector *) const
UInt_t GetGroupNumber()
virtual KVDetectorSignal * GetDetectorSignal(const KVString &type) const
Definition: KVDetector.h:539
virtual void ReadDefinitionFromFile(const Char_t *)
Bool_t IsCalibrated() const
Definition: KVDetector.h:406
Double_t GetEnergyLoss() const override
Definition: KVDetector.h:385
void SetShape(TGeoBBox *s) override
Definition: KVDetector.h:187
void Copy(TObject &obj) const override
Definition: KVDetector.cpp:130
virtual TF1 * GetELossFunction(Int_t Z, Int_t A)
KVUnownedList fParticles
list of particles hitting detector in an event
Definition: KVDetector.h:259
void SetActiveLayer(KVMaterial *actif)
Definition: KVDetector.h:303
KVDetector()
default ctor
Definition: KVDetector.cpp:59
Bool_t BelongsToUnidentifiedParticle() const
Definition: KVDetector.h:567
void DetectParticle(KVNucleus *, TVector3 *norm=0) override
Definition: KVDetector.cpp:195
Double_t RangeDet(Double_t *x, Double_t *par)
Definition: KVDetector.cpp:956
Double_t GetParticleEIncFromERes(KVNucleus *, TVector3 *norm=0) override
Definition: KVDetector.cpp:320
virtual void RemoveCalibrators()
Definition: KVDetector.cpp:699
virtual Double_t GetEntranceWindowSurfaceArea()
Return surface area of first layer of detector in cm2.
TGeoBBox * GetShape() const override
Definition: KVDetector.h:195
KVUniqueNameList fDetSignals
list of signals associated with detector
Definition: KVDetector.h:224
Double_t get_corrected_energy(AbsorberStack *stack, KVNucleus *nuc, Double_t e, Bool_t transmission)
Definition: KVDetector.h:230
Int_t fUnidentP
temporary counters, determine state of identified/unidentified particle flags
Definition: KVDetector.h:178
TF1 * fRangeF
parametric function range of particles in detector
Definition: KVDetector.h:268
Bool_t AddDetectorSignalExpression(const KVString &type, const KVString &_expr)
virtual void SetEResAfterDetector(Double_t e)
Definition: KVDetector.h:623
void SetActiveLayerMatrix(const TGeoHMatrix *)
Set ROOT geometry global matrix transformation to coordinate frame of active layer volume.
Double_t GetERes(Int_t Z, Int_t A, Double_t Einc, Double_t=0.) override
Double_t GetMaxDeltaE(Int_t Z, Int_t A) override
void RemoveParentStructure(KVGeoStrucElement *elem)
virtual Double_t GetSmallestEmaxValid(Int_t Z, Int_t A) const
const KVSeqCollection & GetListOfDetectorSignals() const
Definition: KVDetector.h:784
KVGeoDetectorNode * GetNode()
Definition: KVDetector.h:345
Bool_t fSimMode
=kTRUE when using to simulate detector response, =kFALSE when analysing data
Definition: KVDetector.h:272
Double_t GetDeltaEFromERes(Int_t Z, Int_t A, Double_t Eres) override
void SetMatrix(const TGeoHMatrix *m) override
Definition: KVDetector.h:183
void SetEntranceWindowMatrix(const TGeoHMatrix *)
Set ROOT geometry global matrix transformation to coordinate frame of entrance window.
Double_t EResDet(Double_t *x, Double_t *par)
KVMaterial * fActiveLayer
The active absorber in the detector.
Definition: KVDetector.h:167
void AddHit(KVNucleus *part)
Definition: KVDetector.h:420
void RemoveAllAbsorbers()
Definition: KVDetector.cpp:646
Double_t GetEIncOfMaxDeltaE(Int_t Z, Int_t A) override
Double_t GetTotalThicknessInCM() const
Definition: KVDetector.h:334
KVCalibrator * GetCalibrator(const Char_t *name, const Char_t *type) const
Definition: KVDetector.h:830
Int_t fIdentP
temporary counters, determine state of identified/unidentified particle flags
Definition: KVDetector.h:177
void SetPressure(Double_t P) override
Double_t GetELostByParticle(KVNucleus *, TVector3 *norm=0) override
Definition: KVDetector.cpp:268
void SetEntranceWindowShape(TGeoBBox *)
Set ROOT geometry shape of entrance window.
Bool_t fPresent
=kTRUE if detector is present, =kFALSE if it has been removed
Definition: KVDetector.h:273
void SetActiveLayerShape(TGeoBBox *)
Set ROOT geometry shape of active layer volume.
void SetAnalysed(Bool_t b=kTRUE)
Definition: KVDetector.h:453
Bool_t fSingleLayer
=kTRUE if detector has a single absorber layer
Definition: KVDetector.h:276
virtual TF1 * GetRangeFunction(Int_t Z, Int_t A)
virtual Double_t GetCorrectedEnergy(KVNucleus *, Double_t e=-1., Bool_t transmission=kTRUE)
Definition: KVDetector.cpp:749
void init()
default initialisations
Definition: KVDetector.cpp:32
void Clear(Option_t *opt="") override
Definition: KVDetector.cpp:565
KVList fAbsorbers
list of absorbers making up the detector
Definition: KVDetector.h:260
Bool_t AddCalibrator(KVCalibrator *cal, const KVNameValueList &opts="")
Definition: KVDetector.cpp:439
Simple wrapper for objects which can be drawn (graphs, histograms)
Definition: KVDrawable.h:29
const Char_t * GetFullPathToNode() const
Base class describing elements of array geometry.
Group of detectors which can be treated independently of all others in array.
Definition: KVGroup.h:20
Extended TList class which owns its objects by default.
Definition: KVList.h:22
A stack of materials in which successive energy losses of charged particles can be calculated ,...
Double_t GetMinimumIncidentAngleForDEMax(Int_t Z, Int_t A, Double_t Emax)
void SetIncidenceAngle(double psi)
Double_t GetMaxDeltaE(Int_t Z, Int_t A)
Description of physical materials used to construct detectors & targets; interface to range tables.
Definition: KVMaterial.h:90
virtual void SetPressure(Double_t)
Definition: KVMaterial.cpp:580
virtual void SetTemperature(Double_t)
Definition: KVMaterial.cpp:646
virtual void SetThickness(Double_t thick)
Definition: KVMaterial.cpp:451
virtual Double_t GetThickness() const
Definition: KVMaterial.cpp:484
virtual TGeoMedium * GetGeoMedium(const Char_t *="")
void Clear(Option_t *opt="") override
Reset absorber - set stored energy lost by particles in absorber to zero.
virtual Double_t GetPunchThroughEnergy(Int_t Z, Int_t A)
virtual void SetMaterial(const Char_t *type)
Definition: KVMaterial.cpp:217
virtual Double_t GetERes(Int_t Z, Int_t A, Double_t Einc, Double_t dx=0.)
KVMaterial()
default ctor
Definition: KVMaterial.cpp:73
virtual Double_t GetDeltaE(Int_t Z, Int_t A, Double_t Einc, Double_t dx=0.)
Definition: KVMaterial.cpp:833
virtual Double_t GetLinearRange(Int_t Z, Int_t A, Double_t Einc)
Definition: KVMaterial.cpp:901
Double_t GetMass() const
Definition: KVMaterial.cpp:302
Handles lists of named parameters with different types, a list of KVNamedParameter objects.
Double_t GetDoubleValue(const Char_t *name) const
void SetValue(const Char_t *name, value_type value)
const Char_t * GetStringValue(const Char_t *name) const
Bool_t HasDoubleParameter(const Char_t *name) const
Bool_t HasParameter(const Char_t *name) const
Description of properties and kinematics of atomic nuclei.
Definition: KVNucleus.h:123
Int_t GetA() const
Definition: KVNucleus.cpp:796
void SetZ(Int_t z, Char_t mt=-1)
Definition: KVNucleus.cpp:700
void SetMassFormula(UChar_t mt)
Definition: KVNucleus.h:341
Int_t GetZ() const
Return the number of proton / atomic number.
Definition: KVNucleus.cpp:767
TVector3 * GetPInitial() const
Definition: KVParticle.h:729
TVector3 GetMomentum() const
Definition: KVParticle.h:607
KVNameValueList * GetParameters() const
Definition: KVParticle.h:818
Double_t GetEnergy() const
Definition: KVParticle.h:624
void SetKE(Double_t ecin)
Definition: KVParticle.cpp:246
void SetMomentum(const TVector3 *v)
Definition: KVParticle.h:545
void SetE0(TVector3 *e=0)
Definition: KVParticle.h:718
void SetIsDetected()
Definition: KVParticle.h:733
Double_t GetKE() const
Definition: KVParticle.h:617
void SetEnergy(Double_t e)
Definition: KVParticle.h:602
Base class used for handling geometry in a multidetector array.
Definition: KVPosition.h:91
virtual void SetShape(TGeoBBox *)
Definition: KVPosition.cpp:758
virtual Double_t GetSurfaceArea(int npoints=100000) const
Definition: KVPosition.cpp:989
virtual void SetMatrix(const TGeoHMatrix *)
Definition: KVPosition.cpp:717
virtual TGeoBBox * GetShape() const
Definition: KVPosition.cpp:805
virtual Bool_t ROOTGeo() const
Definition: KVPosition.h:210
KaliVeda extensions to ROOT collection classes.
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
Int_t GetSize() const override
void Clear(Option_t *option="") override
virtual TObject * FindObjectByType(const Char_t *) const
TObject * At(Int_t idx) const override
void Delete(Option_t *option="") override
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
void RBegin(TString delim) const
Definition: KVString.cpp:768
Bool_t End() const
Definition: KVString.cpp:634
KVString Next(Bool_t strip_whitespace=kFALSE) const
Definition: KVString.cpp:695
KVString RNext(Bool_t strip_whitespace=kFALSE) const
Definition: KVString.cpp:841
void Add(TObject *obj) override
Handle several calibrations valid for different Z ranges.
virtual void Print(Option_t *option, const char *wildcard, Int_t recurse=1) const
virtual Int_t GetEntries() const
virtual const char * GetValue(const char *name, const char *dflt) const
virtual void SetRange(Double_t xmin, Double_t xmax)
void SetTitle(const char *title="") override
virtual void SetNpx(Int_t npx=100)
virtual void GetRange(Double_t &xmin, Double_t &xmax) const
virtual void SetParameters(const Double_t *params)
virtual Double_t Eval(Double_t x, Double_t y=0, Double_t z=0, Double_t t=0) const
virtual Double_t GetMaximum(Double_t xmin=0, Double_t xmax=0, Double_t epsilon=1.E-10, Int_t maxiter=100, Bool_t logx=false) const
virtual Double_t GetMaximumX(Double_t xmin=0, Double_t xmax=0, Double_t epsilon=1.E-10, Int_t maxiter=100, Bool_t logx=false) const
Double_t * GetVertices()
static TClass * Class()
virtual Double_t GetDX() const
virtual Double_t GetFacetArea(Int_t index=0) const
virtual Double_t GetDY() const
static TClass * Class()
TClass * IsA() const override
void RefreshPhysicalNodes(Bool_t lock=kTRUE)
TGeoPhysicalNode * MakePhysicalNode(const char *path=nullptr)
TObjArray * GetListOfPhysicalNodes()
Bool_t Align(TGeoMatrix *newmat=nullptr, TGeoShape *newshape=nullptr, Bool_t check=kFALSE, Double_t ovlp=0.001)
TGeoVolume * GetVolume(Int_t level=-1) const
TGeoShape * GetShape(Int_t level=-1) const
virtual Double_t GetRmin() const
static TClass * Class()
virtual Double_t GetRmax() const
virtual void SetMedium(TGeoMedium *medium)
virtual void SetPoint(Int_t i, Double_t x, Double_t y)
void SetName(const char *name="") override
void SetTitle(const char *title="") override
void Reset()
TObject * Clone(const char *newname="") const override
const char * GetName() const override
virtual void SetName(const char *name)
TObject * FindObject(const char *name) const override
R__ALWAYS_INLINE Bool_t TestBit(UInt_t f) const
virtual const char * ClassName() const
virtual Bool_t InheritsFrom(const char *classname) const
void ResetBit(UInt_t f)
Longptr_t ExecPlugin(int nargs)
const char * Data() const
Bool_t Contains(const char *pat, ECaseCompare cmp=kExact) const
RVec< PromoteType< T > > abs(const RVec< T > &v)
Double_t x[n]
gr SetName("gr")
double T(double x)
void init()
constexpr Double_t K()
TMarker m
TArc a
ClassImp(TPyArg)