KaliVeda  1.12/06
Heavy-Ion Analysis Toolkit
KVDetector.cpp
Go to the documentation of this file.
1 #include "KVDetector.h"
2 #include "Riostream.h"
3 #include "TROOT.h"
4 #include "TRandom3.h"
5 #include "TBuffer.h"
6 #include "KVUnits.h"
7 #include "KVTelescope.h"
8 #include "KVGroup.h"
9 #include "KVCalibrator.h"
10 #include "KVACQParam.h"
11 #include "TPluginManager.h"
12 #include "TObjString.h"
13 #include "TClass.h"
14 /*** geometry ***/
15 #include "TGeoVolume.h"
16 #include "TGeoManager.h"
17 #include "TGeoMatrix.h"
18 #include "TGeoBBox.h"
19 #include "TGeoArb8.h"
20 #include "TGeoTube.h"
21 #include "KVACQParamSignal.h"
22 #include "KVCalibratedSignal.h"
25 
26 #include <TGeoPhysicalNode.h>
27 #include <TGraph.h>
28 //#include <KVGeoDNTrajectory.h>
29 
30 using namespace std;
31 
33 
34 
36 
37 
40 
42 {
43  //default initialisations
44  fCalibrators = 0;
45  fACQParams = 0;
46  fParticles = 0;
47  fSegment = 0;
48  fGain = 1.;
49  fCalWarning = 0;
50  fAbsorbers = new KVList;
51  fActiveLayer = nullptr;
52  fIDTelescopes = new KVList(kFALSE);
53  fIDTelescopes->SetCleanup(kTRUE);
54  fIDTelAlign = new KVList(kFALSE);
55  fIDTelAlign->SetCleanup(kTRUE);
56  fIDTele4Ident = 0;
57  fIdentP = fUnidentP = 0;
58  fTotThickness = 0.;
59  fDepthInTelescope = 0.;
60  fFiredMask.Set("0");
61  fELossF = fEResF = fRangeF = 0;
62  fEResforEinc = -1.;
63  fAlignedDetectors[0] = 0;
64  fAlignedDetectors[1] = 0;
65  fSimMode = kFALSE;
66  fPresent = kTRUE;
67  fDetecting = kTRUE;
68  fParentStrucList.SetCleanup();
69  fSingleLayer = kTRUE;
70  fNode.SetDetector(this);
71  SetKVDetectorFiredACQParameterListFormatString();
72  // detector owns any signals which are defined for it
73  fDetSignals.SetOwner();
74  // adding a new signal with the same name as an existing one
75  // will delete the existing signal and replace it
76  fDetSignals.ReplaceObjects();
77 }
78 
79 
80 
86 
88 {
89  // Set the value of fKVDetectorFiredACQParameterListFormatString to be
90  // [classname].Fired.ACQParameterList.%s
91  // where [classname] is the name of the class of whatever object
92  // is calling this method
93 
94  fKVDetectorFiredACQParameterListFormatString.Form("%s.Fired.ACQParameterList.",
95  ClassName());
96  fKVDetectorFiredACQParameterListFormatString += "%s";
97 }
98 
99 
100 
103 
105 {
106 //default ctor
107  init();
108  fDetCounter++;
109  SetName(Form("Det_%d", fDetCounter));
110 }
111 
112 
113 
116 
118  const Float_t thick): KVMaterial()
119 {
120  // Create a new detector of a given material and thickness in centimetres (default value = 0.0)
121 
122  init();
123  SetType("DET");
124  fDetCounter++;
125  SetName(Form("Det_%d", fDetCounter));
126  AddAbsorber(new KVMaterial(type, thick));
127 }
128 
129 
130 
133 
135 {
136 //copy ctor
137  init();
138  fDetCounter++;
139  SetName(Form("Det_%d", fDetCounter));
140 #if ROOT_VERSION_CODE >= ROOT_VERSION(3,4,0)
141  obj.Copy(*this);
142 #else
143  ((KVDetector&) obj).Copy(*this);
144 #endif
145 }
146 
147 
148 #if ROOT_VERSION_CODE >= ROOT_VERSION(3,4,0)
149 
154 
155 void KVDetector::Copy(TObject& obj) const
156 #else
157 void KVDetector::Copy(TObject& obj)
158 #endif
159 {
160  //copy 'this' to 'obj'
161  //The structure of the detector is copied, with new cloned objects for
162  //each absorber layer. The active layer is set in the new detector.
163 
164  TIter next(fAbsorbers);
165  KVMaterial* mat;
166  while ((mat = (KVMaterial*) next())) {
167  ((KVDetector&) obj).AddAbsorber((KVMaterial*) mat->Clone());
168  }
169  //set active layer
170  Int_t in_actif = fAbsorbers->IndexOf(fActiveLayer);
171  ((KVDetector&) obj).SetActiveLayer(((KVDetector&)obj).GetAbsorber(in_actif));
172 }
173 
174 
175 
176 
178 
179 KVDetector::~KVDetector()
180 {
181  fIDTelescopes->Clear();
185  delete fAbsorbers;
190  fIDTelAlign->Clear();
193  if (fAlignedDetectors[0]) fAlignedDetectors[0]->Clear("nodelete");
195  if (fAlignedDetectors[1]) fAlignedDetectors[1]->Clear("nodelete");
197 }
198 
199 
200 
205 
207 {
208  // Set material of active layer.
209  // If no absorbers have been added to the detector, create and add
210  // one (active layer by default)
211 
212  if (!GetActiveLayer())
214  else
216 }
217 
218 
219 
231 
233 {
234  //Calculate the energy loss of a charged particle traversing the detector,
235  //the particle is slowed down, it is added to the list of all particles hitting the
236  //detector. The apparent energy loss of the particle in the active layer of the
237  //detector is set.
238  //Do nothing if particle has zero (or -ve) energy.
239  //
240  //If the optional argument 'norm' is given, it is supposed to be a vector
241  //normal to the detector, oriented from the origin towards the detector.
242  //In this case the effective thicknesses of the detector's absorbers 'seen' by the particle
243  //depending on its direction of motion is used for the calculation.
244 
245  if (kvp->GetKE() <= 0.)
246  return;
247 
248  AddHit(kvp); //add nucleus to list of particles hitting detector in the event
249  //set flag to say that particle has been slowed down
250  kvp->SetIsDetected();
251  //If this is the first absorber that the particle crosses, we set a "reminder" of its
252  //initial energy
253  if (!kvp->GetPInitial())
254  kvp->SetE0();
255 
256  std::vector<Double_t> thickness;
257  if (norm) {
258  // modify thicknesses of all layers according to orientation,
259  // and store original thicknesses in array
260  TVector3 p = kvp->GetMomentum();
261  thickness.reserve(fAbsorbers->GetEntries());
262  KVMaterial* abs;
263  int i = 0;
264  TIter next(fAbsorbers);
265  while ((abs = (KVMaterial*) next())) {
266  thickness[i++] = abs->GetThickness();
267  Double_t T = abs->GetEffectiveThickness((*norm), p);
268  abs->SetThickness(T);
269  }
270  }
271  Double_t eloss = GetTotalDeltaE(kvp->GetZ(), kvp->GetA(), kvp->GetEnergy());
272  Double_t dE = GetDeltaE(kvp->GetZ(), kvp->GetA(), kvp->GetEnergy());
273  if (norm) {
274  // reset thicknesses of absorbers
275  KVMaterial* abs;
276  int i = 0;
277  TIter next(fAbsorbers);
278  while ((abs = (KVMaterial*) next())) {
279  abs->SetThickness(thickness[i++]);
280  }
281  }
282  Double_t epart = kvp->GetEnergy() - eloss;
283  if (epart < 1e-3) {
284  //printf("%s, pb d arrondi on met l energie de la particule a 0\n",GetName());
285  epart = 0.0;
286  }
287  kvp->SetEnergy(epart);
288  Double_t eloss_old = GetEnergyLoss();
289  SetEnergyLoss(eloss_old + dE);
290 }
291 
292 
293 
294 
304 
306 {
307  //Calculate the total energy loss of a charged particle traversing the detector.
308  //This does not affect the "stored" energy loss value of the detector,
309  //nor its ACQData, nor the energy of the particle.
310  //
311  //If the optional argument 'norm' is given, it is supposed to be a vector
312  //normal to the detector, oriented from the origin towards the detector.
313  //In this case the effective thicknesses of the detector's absorbers 'seen' by the particle
314  //depending on its direction of motion is used for the calculation.
315 
316  std::vector<Double_t> thickness;
317  if (norm) {
318  // modify thicknesses of all layers according to orientation,
319  // and store original thicknesses in array
320  TVector3 p = kvp->GetMomentum();
321  thickness.reserve(fAbsorbers->GetEntries());
322  KVMaterial* abs;
323  int i = 0;
324  TIter next(fAbsorbers);
325  while ((abs = (KVMaterial*) next())) {
326  thickness[i++] = abs->GetThickness();
327  Double_t T = abs->GetEffectiveThickness((*norm), p);
328  abs->SetThickness(T);
329  }
330  }
331  Double_t eloss = GetTotalDeltaE(kvp->GetZ(), kvp->GetA(), kvp->GetEnergy());
332  if (norm) {
333  // reset thicknesses of absorbers
334  KVMaterial* abs;
335  int i = 0;
336  TIter next(fAbsorbers);
337  while ((abs = (KVMaterial*) next())) {
338  abs->SetThickness(thickness[i++]);
339  }
340  }
341  return eloss;
342 }
343 
344 
345 
346 
356 
358 {
359  //Calculate the energy of particle 'kvn' before its passage through the detector,
360  //based on the current kinetic energy, Z & A of nucleus 'kvn', supposed to be
361  //after passing through the detector.
362  //
363  //If the optional argument 'norm' is given, it is supposed to be a vector
364  //normal to the detector, oriented from the origin towards the detector.
365  //In this case the effective thicknesses of the detector's absorbers 'seen' by the particle
366  //depending on its direction of motion is used for the calculation.
367 
368  Double_t Einc = 0.0;
369  //make 'clone' of particle
370  KVNucleus clone_part(kvp->GetZ(), kvp->GetA());
371  clone_part.SetMomentum(kvp->GetMomentum());
372  //composite detector - calculate losses in all layers
373  KVMaterial* abs;
374  TIter next(fAbsorbers, kIterBackward); //work through list backwards
375  while ((abs = (KVMaterial*) next())) {
376 
377  //calculate energy of particle before current absorber
378  Einc = abs->GetParticleEIncFromERes(&clone_part, norm);
379 
380  //set new energy of particle
381  clone_part.SetKE(Einc);
382  }
383  return Einc;
384 }
385 
386 
387 
388 
392 
393 void KVDetector::Print(Option_t* opt) const
394 {
395  //Print info on this detector
396  //if option="data" the energy loss and coder channel data are displayed
397 
398  if (!strcmp(opt, "data")) {
399  cout << ((KVDetector*) this)->
400  GetName() << " -- E=" << ((KVDetector*) this)->
401  GetEnergy();
402  cout << " ";
403  TIter next(fACQParams);
404  KVACQParam* acq;
405  while ((acq = (KVACQParam*) next())) {
406  cout << acq->GetName() << "=" << (Short_t) acq->
407  GetCoderData() << "/" << TMath::Nint(acq->GetPedestal());
408  cout << " ";
409  }
411  cout << "(Belongs to an unidentified particle)";
412  cout << endl;
413  }
414  else if (!strcmp(opt, "all")) {
415  //Give full details of detector
416  //
417  TString option(" ");
418  cout << option << ClassName() << " : " << ((KVDetector*) this)->
419  GetName() << endl;
420  //composite detector - print all layers
421  KVMaterial* abs;
422  TIter next(fAbsorbers);
423  while ((abs = (KVMaterial*) next())) {
424  if (GetActiveLayer() == abs)
425  cout << " ### ACTIVE LAYER ### " << endl;
426  cout << option << option;
427  abs->Print();
428  if (GetActiveLayer() == abs)
429  cout << " #################### " << endl;
430  }
431  cout << option << "Gain: " << GetGain() << endl;
432  if (fParticles) {
433  cout << option << " --- Detected particles:" << endl;
434  fParticles->Print();
435  }
436  if (fIDTelescopes) {
437  cout << option << " --- Detector belongs to the following Identification Telescopes:" << endl;
438  fIDTelescopes->ls();
439  }
440  if (fIDTelAlign) {
441  cout << option << " --- Identification Telescopes in front of detector:" << endl;
442  fIDTelAlign->ls();
443  }
444  if (fIDTele4Ident) {
445  cout << option << " --- Identification Telescopes used to identify particles stopping in this detector:" << endl;
446  fIDTele4Ident->ls();
447  }
448  }
449  else {
450  //just print name
451  cout << ClassName() << ": " << ((KVDetector*) this)->
452  GetName() << endl;
453  }
454 }
455 
456 
457 
458 
459 
470 
472 {
473  // This method is called by KVASMultiDetArray::MakeListOfDetectors
474  // after the array geometry has been defined (i.e. all detectors have
475  // been placed in the array). The string returned by this method
476  // is used to set the name of the detector.
477  //
478  // Override this method in child classes in order to define a naming
479  // convention for specific detectors of the array.
480  //
481  // By default we return the same name as KVDetector::GetName
482 
483  fFName = GetName();
484  return fFName.Data();
485 }
486 
487 
488 
509 
511 {
512  // Associate a calibration with this detector (the object will be deleted by the detector)
513  //
514  // This will add a new signal to the list of the detector's signals
515  //
516  // Also sets calibrator's name to [detname]_[caltype]
517  //
518  // If the input signal required by the calibrator is not defined for the detector,
519  // this method returns kFALSE and the calibrator will be deleted.
520  //
521  // The (optional) KVNameValueList argument can be used to pass any extra parameters/options.
522  // For example, if it contains a parameter `ZRange`:
523  //
524  //~~~~~~~~~~~~~~~
525  // ZRange=1-10
526  //~~~~~~~~~~~~~~~
527  //
528  // then the calibrator will be handled by a KVZDependentCalibratedSignal (handles several
529  // calibrators which provide the same output signal, each one is used for a specific range
530  // of atomic numbers)
531 
532  if (!cal) return kFALSE;
533  // check for input signal
535  if (!in) {
536  Warning("AddCalibrator", "%s : input signal %s not found for calibrator %s. No output signal created.",
537  GetName(), cal->GetInputSignalType().Data(), cal->GetType());
538  delete cal;
539  return kFALSE;
540  }
541  if (!fCalibrators)
542  fCalibrators = new KVList();
543 
544  cal->SetDetector(this);
545  fCalibrators->Add(cal);
546  cal->SetName(Form("%s_%s", GetName(), cal->GetType()));
547 
548  if (cal->GetOutputSignalType() == "") {
549  Warning("AddCalibrator", "%s : output signal not defined for calibrator %s. No output signal created.",
550  GetName(), cal->GetType());
551  }
552  else {
553  KVDetectorSignal* new_cal_sig(nullptr);
554  if (opts.HasParameter("ZRange")) {
555  // If 'ZRange' parameter is given we need to find the KVZDependentCalibratedSignal
556  // and add a new signal to it.
558  if (!sig) sig = new_cal_sig = new KVZDependentCalibratedSignal(in, cal->GetOutputSignalType());
559  dynamic_cast<KVZDependentCalibratedSignal*>(sig)->AddSignal(
560  new KVCalibratedSignal(in, cal), opts.GetStringValue("ZRange")
561  );
562  }
563  else {
564  new_cal_sig = new KVCalibratedSignal(in, cal);
565  }
566  if (new_cal_sig) fDetSignals.Add(new_cal_sig);
567  }
568  return kTRUE;
569 }
570 
571 
572 
583 
585 {
586  // Replace calibrator of given type with the given calibrator object
587  // The calibrator object should not be shared with any other detectors: it now belongs
588  // to this detector, which will delete it when necessary.
589  // If an exising calibrator with the same type is already defined, it will be
590  // deleted and removed from the detector's calibrator list
591  //
592  // Returns kFALSE in case of problems.
593  //
594  // The (optional) KVNameValueList argument can be used to pass any extra parameters/options.
595 
596  KVCalibrator* old_cal = GetCalibrator(type);
597  if (old_cal) {
598  fCalibrators->Remove(old_cal);
600  delete old_cal;
601  }
602  return AddCalibrator(cal, opts);
603 }
604 
605 
606 
611 
613 {
614  // A detector is considered to be calibrated if it has
615  // a signal "Energy" available and if depending on the supplied parameters
616  // this signal can be calculated
617 
618  KVCalibratedSignal* e_sig = dynamic_cast<KVCalibratedSignal*>(GetDetectorSignal("Energy"));
619  return (e_sig && e_sig->IsAvailableFor(params) && IsOK());
620 }
621 
622 
623 
626 
628 {
629  // Add given acquisition parameter to this detector.
630 
631  if (!fACQParams) {
632  fACQParams = new KVList;
633  fACQParams->SetName(Form("List of acquisition parameters for detector %s", GetName()));
634  }
635  par->SetDetector(this);
636  fACQParams->Add(par);
637  fDetSignals.Add(new KVACQParamSignal(par));
638 }
639 
640 
641 
645 
647 {
648  // Look for acquisition parameter with given name in list
649  // of parameters associated with this detector.
650 
651  if (!fACQParams) {
652  return 0;
653  }
654  return ((KVACQParam*) fACQParams->FindObject(name));
655 }
656 
657 
658 
665 
667 {
668  // Access acquisition data value associated to parameter with given name.
669  // Returns value as a floating-point number which is the raw channel number
670  // read from the coder plus a random number in the range [-0.5,+0.5].
671  // If the detector has no DAQ parameter of the given type,
672  // or if the raw channel number = 0, the value returned is -1.
673 
674  KVACQParam* par = GetACQParam(name);
675  return (par ? par->GetData() : -1.);
676 }
677 
678 
679 
682 
684 {
685  // Access pedestal value associated to parameter with given name.
686 
687  KVACQParam* par = GetACQParam(name);
688  return (par ? par->GetPedestal() : 0);
689 }
690 
691 
692 
695 
696 void KVDetector::SetPedestal(const Char_t* name, Float_t ped)
697 {
698  // Set value of pedestal associated to parameter with given name.
699 
700  KVACQParam* par = GetACQParam(name);
701  if (par) {
702  par->SetPedestal(ped);
703  }
704 }
705 
706 
707 
708 
712 
714 {
715  //Set energy loss(es) etc. to zero
716  //If opt="N" we do not reset acquisition parameters/raw detector signals
717 
719  fIdentP = fUnidentP = 0;
722  if (strncmp(opt, "N", 1)) {
723  if (fACQParams) {
724  TIter next(fACQParams);
725  KVACQParam* par;
726  while ((par = (KVACQParam*) next())) {
727  par->Clear();
728  }
729  }
731  KVDetectorSignal* ds;
732  while ((ds = (KVDetectorSignal*)it())) {
733  if (ds->IsRaw() && !ds->IsExpression()) ds->Reset();
734  }
735  }
736  ClearHits();
737  //reset all layers in detector
738  KVMaterial* mat;
739  TIter next(fAbsorbers);
740  while ((mat = (KVMaterial*) next())) {
741  mat->Clear();
742  }
743  fEResforEinc = -1.;
744 }
745 
746 
747 
752 
754 {
755  // Add a layer of absorber material to the detector
756  // By default, the first layer added is set as the "Active" layer.
757  // Call SetActiveLayer to change this.
758  fAbsorbers->Add(mat);
759  if (!TestBit(kActiveSet))
760  SetActiveLayer(mat);
761  if (fAbsorbers->GetSize() > 1) fSingleLayer = kFALSE;
762 }
763 
764 
765 
768 
770 {
771  //Returns pointer to the i-th absorber in the detector (i=0 first absorber, i=1 second, etc.)
772 
773  if (i < 0) {
774  //special case of reading old detectors
775  //no warning
776  return 0;
777  }
778  if (!fAbsorbers) {
779  Error("GetAbsorber", "No absorber list");
780  return 0;
781  }
782  if (fAbsorbers->GetSize() < 1) {
783  Error("GetAbsorber", "Nothing in absorber list");
784  return 0;
785  }
786  if (i >= fAbsorbers->GetSize()) {
787  Error("GetAbsorber", "i=%d but only %d absorbers in list", i,
788  fAbsorbers->GetSize());
789  return 0;
790  }
791 
792  return (KVMaterial*) fAbsorbers->At(i);
793 }
794 
795 
796 
801 
803 {
804  //Attribute acquisition parameters to this detector.
805  //This method does nothing; it should be overridden in child classes to attribute
806  //parameters specific to each detector.
807  ;
808 }
809 
810 
811 
815 
817 {
818  // Used when a calibrator object is removed or replaced
819  // We remove and delete the corresponding output signal from the list of detector signals
820 
821  if (K->GetOutputSignalType() != "") {
822  KVDetectorSignal* ds = GetDetectorSignal(K->GetOutputSignalType());
823  if (ds) {
824  fDetSignals.Remove(ds);
825  delete ds;
826  }
827  }
828 }
829 
830 
831 
837 
839 {
840  // Removes all calibrations associated to this detector: in other words, we delete all
841  // the KVCalibrator objects in list fCalibrators.
842  //
843  // We also destroy all signals provided by these calibrators
844 
845  if (fCalibrators) {
846  KVCalibrator* K;
847  TIter it(fCalibrators);
848  while ((K = (KVCalibrator*)it())) {
850  }
851  fCalibrators->Delete();
852  }
853 }
854 
855 
856 
857 
860 
862 {
863  //Add ID telescope to list of telescopes to which detector belongs
864  fIDTelescopes->Add(idt);
865 }
866 
867 
868 
874 
876 {
877  // Return list of all ID telescopes containing detectors placed in front of
878  // this one.
879 
880  // temporary kludge during transition to trajectory-based reconstruction
881  // ROOT-geometry-based detectors will not have fIDTelAlign filled
882  if (ROOTGeo() && !fIDTelAlign->GetEntries()) {
884  (KVGeoDNTrajectory*)GetNode()->GetTrajectories()->First(),
885  GetNode()
886  );
887  if (Rtr) fIDTelAlign->AddAll(Rtr->GetIDTelescopes());
888  }
889  return fIDTelAlign;
890 }
891 
892 
893 
894 
899 
901 {
902  //Returns list of identification telescopes to be used in order to try to identify
903  //particles stopping in this detector. This is the same as GetAlignedIDTelescopes
904  //but only including the telescopes of which this detector is a member.
905  if (fIDTele4Ident) return fIDTele4Ident;
906  if (!fIDTelescopes || !fIDTelAlign) return 0;
907  fIDTele4Ident = new TList;
909  TObject* idt;
910  while ((idt = next())) {
911  if (fIDTelescopes->FindObject(idt)) fIDTele4Ident->Add(idt);
912  }
913  return fIDTele4Ident;
914 }
915 
916 
917 
918 
919 
948 
950 {
951  // Returns the total energy loss in the detector for a given nucleus
952  // including inactive absorber layers.
953  // e = energy loss in active layer (if not given, we use current value)
954  // transmission = kTRUE (default): the particle is assumed to emerge with
955  // a non-zero residual energy Eres after the detector.
956  // = kFALSE: the particle is assumed to stop in the detector.
957  //
958  // WARNING: if transmission=kTRUE, and if the residual energy after the detector
959  // is known (i.e. measured in a detector placed after this one), you should
960  // first call
961  // SetEResAfterDetector(Eres);
962  // before calling this method. Otherwise, especially for heavy ions, the
963  // correction may be false for particles which are just above the punch-through energy.
964  //
965  // WARNING 2: if measured energy loss in detector active layer is greater than
966  // maximum possible theoretical value for given nucleus' Z & A, this may be because
967  // the A was not measured but calculated from Z and hence could be false, or perhaps
968  // there was an (undetected) pile-up of two or more particles in the detector.
969  // In this case we return the uncorrected energy measured in the active layer
970  // and we add the following parameters to the particle (in nuc->GetParameters()):
971  //
972  // GetCorrectedEnergy.Warning = 1
973  // GetCorrectedEnergy.Detector = [name]
974  // GetCorrectedEnergy.MeasuredDE = [value]
975  // GetCorrectedEnergy.MaxDE = [value]
976  // GetCorrectedEnergy.Transmission = 0 or 1
977  // GetCorrectedEnergy.ERES = [value]
978 
979  Int_t z = nuc->GetZ();
980  Int_t a = nuc->GetA();
981 
982  if (e < 0.) e = GetEnergy();
983  if (e <= 0) {
985  return 0;
986  }
987 
988  enum SolType solution = kEmax;
989  if (!transmission) solution = kEmin;
990 
991  // check that apparent energy loss in detector is compatible with a & z
992  Double_t maxDE = GetMaxDeltaE(z, a);
993  Double_t EINC, ERES = GetEResAfterDetector();
994  if (e > maxDE) {
995  nuc->GetParameters()->SetValue("GetCorrectedEnergy.Warning", 1);
996  nuc->GetParameters()->SetValue("GetCorrectedEnergy.Detector", GetName());
997  nuc->GetParameters()->SetValue("GetCorrectedEnergy.MeasuredDE", e);
998  nuc->GetParameters()->SetValue("GetCorrectedEnergy.MaxDE", maxDE);
999  nuc->GetParameters()->SetValue("GetCorrectedEnergy.Transmission", (Int_t)transmission);
1000  nuc->GetParameters()->SetValue("GetCorrectedEnergy.ERES", ERES);
1001  return e;
1002 
1003  }
1004  if (transmission && ERES > 0.) {
1005  // if residual energy is known we use it to calculate EINC.
1006  // if EINC < max of dE curve, we change solution
1007  EINC = GetIncidentEnergyFromERes(z, a, ERES);
1008  if (EINC < GetEIncOfMaxDeltaE(z, a)) solution = kEmin;
1009  // we could keep the EINC value calculated using ERES, but then
1010  // the corrected dE of this detector would not depend on the
1011  // measured dE !
1012  }
1013  EINC = GetIncidentEnergy(z, a, e, solution);
1014  if (EINC < 0) {
1015  SetEResAfterDetector(-1.);
1016  return EINC;
1017  }
1018  ERES = GetERes(z, a, EINC);
1019 
1020  SetEResAfterDetector(-1.);
1021  return (EINC - ERES);
1022 }
1023 
1024 
1025 
1026 
1044 
1046 {
1047  //For particles which stop in the first stage of an identification telescope,
1048  //we can at least estimate a minimum Z value based on the energy lost in this
1049  //detector.
1050  //
1051  //This is based on the KVMaterial::GetMaxDeltaE method, giving the maximum
1052  //energy loss in the active layer of the detector for a given nucleus (A,Z).
1053  //
1054  //The "Zmin" is the Z of the particle which gives a maximum energy loss just greater
1055  //than that measured in the detector. Particles with Z<Zmin could not lose as much
1056  //energy and so are excluded.
1057  //
1058  //If ELOSS is not given, we use the current value of GetEnergy()
1059  //Use 'mass_formula' to change the formula used to calculate the A of the nucleus
1060  //from its Z. Default is valley of stability value. (see KVNucleus::GetAFromZ).
1061  //
1062  //If the value of ELOSS or GetEnergy() is <=0 we return Zmin=0
1063 
1064  ELOSS = (ELOSS < 0. ? GetEnergy() : ELOSS);
1065  if (ELOSS <= 0) return 0;
1066 
1067  UInt_t z = 40;
1068  UInt_t zmin, zmax;
1069  zmin = 1;
1070  zmax = 92;
1071  Double_t difference;
1072  UInt_t last_positive_difference_z = 1;
1073  KVNucleus particle;
1074  if (mass_formula > -1)
1075  particle.SetMassFormula((UChar_t)mass_formula);
1076 
1077  difference = 0.;
1078 
1079  while (zmax > zmin + 1) {
1080 
1081  particle.SetZ(z);
1082 
1083  difference = GetMaxDeltaE(z, particle.GetA()) - ELOSS;
1084  //if difference < 0 the z is too small
1085  if (difference < 0.0) {
1086 
1087  zmin = z;
1088  z += (UInt_t)((zmax - z) / 2 + 0.5);
1089 
1090  }
1091  else {
1092 
1093  zmax = z;
1094  last_positive_difference_z = z;
1095  z -= (UInt_t)((z - zmin) / 2 + 0.5);
1096 
1097  }
1098  }
1099  if (difference < 0) {
1100  //if the last z we tried actually made the difference become negative again
1101  //(i.e. z too small) we return the last z which gave a +ve difference
1102  z = last_positive_difference_z;
1103  }
1104  return z;
1105 }
1106 
1107 
1108 
1109 
1118 
1120 {
1121  // Calculates energy loss (in MeV) in active layer of detector, taking into account preceding layers
1122  //
1123  // Arguments are:
1124  // x[0] is incident energy in MeV
1125  // Parameters are:
1126  // par[0] Z of ion
1127  // par[1] A of ion
1128 
1129  Double_t e = x[0];
1130  TIter next(fAbsorbers);
1131  KVMaterial* mat;
1132  mat = (KVMaterial*)next();
1133  while (fActiveLayer != mat) {
1134  // calculate energy losses in absorbers before active layer
1135  e = mat->GetERes(par[0], par[1], e); //residual energy after layer
1136  if (e <= 0.)
1137  return 0.; // return 0 if particle stops in layers before active layer
1138  mat = (KVMaterial*)next();
1139  }
1140  //calculate energy loss in active layer
1141  return mat->GetDeltaE(par[0], par[1], e);
1142 }
1143 
1144 
1145 
1146 
1156 
1158 {
1159  // Calculates range (in centimetres) of ions in detector as a function of incident energy (in MeV),
1160  // taking into account all layers of the detector.
1161  //
1162  // Arguments are:
1163  // x[0] = incident energy in MeV
1164  // Parameters are:
1165  // par[0] = Z of ion
1166  // par[1] = A of ion
1167 
1168  Double_t Einc = x[0];
1169  Int_t Z = (Int_t)par[0];
1170  Int_t A = (Int_t)par[1];
1171  Double_t range = 0.;
1172  TIter next(fAbsorbers);
1173  KVMaterial* mat = (KVMaterial*)next();
1174  if (!mat) return 0.;
1175  do {
1176  // range in this layer
1177  Double_t this_range = mat->GetLinearRange(Z, A, Einc);
1178  KVMaterial* next_mat = (KVMaterial*)next();
1179  if (this_range > mat->GetThickness()) {
1180  // particle traverses layer.
1181  if (next_mat)
1182  range += mat->GetThickness();
1183  else // if this is last layer, the range continues to increase beyond size of detector
1184  range += this_range;
1185  // calculate incident energy for next layer (i.e. residual energy after this layer)
1186  Einc = mat->GetERes(Z, A, Einc);
1187  }
1188  else {
1189  // particle stops in layer
1190  range += this_range;
1191  return range;
1192  }
1193  mat = next_mat;
1194  }
1195  while (mat);
1196  // particle traverses detector
1197  return range;
1198 }
1199 
1200 
1201 
1202 
1212 
1214 {
1215  // Calculates residual energy (in MeV) of particle after traversing all layers of detector.
1216  // Returned value is -1000 if particle stops in one of the layers of the detector.
1217  //
1218  // Arguments are:
1219  // x[0] is incident energy in MeV
1220  // Parameters are:
1221  // par[0] Z of ion
1222  // par[1] A of ion
1223 
1224  Double_t e = x[0];
1225  TIter next(fAbsorbers);
1226  KVMaterial* mat;
1227  while ((mat = (KVMaterial*)next())) {
1228  Double_t eres = mat->GetERes(par[0], par[1], e); //residual energy after layer
1229  if (eres <= 0.)
1230  return -1000.; // return -1000 if particle stops in layers before active layer
1231  e = eres;
1232  }
1233  return e;
1234 }
1235 
1236 
1237 
1238 
1253 
1255 {
1256  //Static function which will create an instance of the KVDetector-derived class
1257  //corresponding to 'name'
1258  //These are defined as 'Plugin' objects in the file $KVROOT/KVFiles/.kvrootrc :
1259  // [name_of_dataset].detector_type
1260  // detector_type
1261  // To use the dataset-dependent plugin, call this method with
1262  // name = "[name_of_dataset].detector_type"
1263  // If not, the default plugin will be used
1264  //first we check if there is a special plugin for the DataSet
1265  //if not we take the default one
1266  //
1267  //'thickness' is passed as argument to the constructor for the detector plugin
1268 
1269  //check and load plugin library
1270  TPluginHandler* ph = nullptr;
1271  KVString nom(name);
1272  if (!nom.Contains(".") && !(ph = LoadPlugin("KVDetector", name))) return nullptr;
1273  if (nom.Contains(".")) {
1274  // name format like [dataset].[det_type]
1275  // in case dataset name contains "." we parse string to find detector type: assumed after the last "."
1276  nom.RBegin(".");
1277  KVString det_type = nom.RNext();
1278  if (!(ph = LoadPlugin("KVDetector", name))) {
1279  if (!(ph = LoadPlugin("KVDetector", det_type))) {
1280  return nullptr;
1281  }
1282  }
1283  }
1284 
1285  //execute constructor/macro for detector
1286  return (KVDetector*) ph->ExecPlugin(1, thickness);
1287 }
1288 
1289 
1290 
1291 
1303 
1304 void KVDetector::GetVerticesInOwnFrame(TVector3* corners, Double_t depth, Double_t layer_thickness)
1305 {
1306  // This will fill the array corners[8] with the coordinates of the vertices of the
1307  // front (corners[0],...,corners[3]) & back (corners[4],...,corners[7]) sides of the volume
1308  // representing either a single absorber or the whole detector.
1309  //
1310  // depth = depth of detector/absorber inside the KVTelescope it belongs to (in centimetres)
1311  // layer_thickness = thickness of absorber/detector (in centimetres)
1312  //
1313  // Positioning information is taken from the KVTelescope to which this detector
1314  // belongs; if this is not the case, nothing will be done.
1315 
1316  // relative distance to back of detector
1317  Double_t dist_to_back = depth + layer_thickness;
1318 
1319  // get coordinates
1320  KVTelescope* fTelescope = (KVTelescope*)GetParentStructure("TELESCOPE");
1321  if (fTelescope) {
1322  fTelescope->GetCornerCoordinatesInOwnFrame(corners, depth);
1323  fTelescope->GetCornerCoordinatesInOwnFrame(&corners[4], dist_to_back);
1324  }
1325  else {
1326  GetCornerCoordinatesInOwnFrame(corners, depth);
1327  GetCornerCoordinatesInOwnFrame(&corners[4], dist_to_back);
1328  }
1329 }
1330 
1331 
1332 
1333 
1348 
1350 {
1351  // Construct a TGeoVolume shape to represent this detector in the current
1352  // geometry managed by gGeoManager.
1353  //
1354  // Making the volume requires:
1355  // - to construct a 'mother' volume (TGeoArb8) defined by the (theta-min/max, phi-min/max)
1356  // and the total thickness of the detector (all layers)
1357  // - to construct and position volumes (TGeoArb8) for each layer of the detector inside the
1358  // 'mother' volume. Each layer is represented by a TGeoArb8 whose two parallel faces correspond
1359  // to the front and back surfaces of the layer.
1360  //
1361  // If the detector is composed of a single absorber, we do not create a superfluous "mother" volume.
1362  //
1363  // gGeoManager must point to current instance of geometry manager.
1364 
1365  if (!gGeoManager) return NULL;
1366 
1367  if (fTotThickness == 0) GetTotalThicknessInCM(); // calculate sum of absorber thicknesses in centimetres
1368  // get from telescope info on relative depth of detector i.e. depth inside telescope
1369 
1370  KVTelescope* fTelescope = (KVTelescope*)GetParentStructure("TELESCOPE");
1371  if (fTelescope && fDepthInTelescope == 0)
1372  fDepthInTelescope = fTelescope->GetDepthInCM(fTelescope->GetDetectorRank(this));
1373 
1374  TVector3 coords[8];
1375  Double_t vertices[16];
1376  Double_t tot_thick_det = fTotThickness;
1377  TGeoMedium* med;
1378  TGeoVolume* mother_vol = 0;
1379 
1380  Bool_t multi_layer = fAbsorbers->GetSize() > 1;
1381 
1382  if (multi_layer) {
1383  mother_vol = gGeoManager->MakeVolumeAssembly(Form("%s_DET", GetName()));
1384  }
1385 
1386  /**** BUILD & ADD ABSORBER VOLUMES ****/
1387  TIter next(fAbsorbers);
1388  KVMaterial* abs;
1389  Double_t depth_in_det = 0.;
1390  Int_t no_abs = 1;
1391  while ((abs = (KVMaterial*)next())) {
1392  // get medium for absorber
1393  med = abs->GetGeoMedium();
1394  Double_t thick = abs->GetThickness();
1395  // do not include layers with zero thickness
1396  if (thick == 0.0) continue;
1397  GetVerticesInOwnFrame(coords, fDepthInTelescope + depth_in_det, thick);
1398  for (int i = 0; i < 8; i++) {
1399  vertices[2 * i] = coords[i].X();
1400  vertices[2 * i + 1] = coords[i].Y();
1401  }
1402  Double_t dz = thick / 2.;
1403  TString vol_name;
1404  if (abs == GetActiveLayer()) vol_name = GetName();
1405  else vol_name = Form("%s_%d_%s", GetName(), no_abs, abs->GetType());
1406  TGeoVolume* vol =
1407  gGeoManager->MakeArb8(vol_name.Data(), med, dz, vertices);
1408  vol->SetLineColor(med->GetMaterial()->GetDefaultColor());
1409  if (multi_layer) {
1410  /*** position absorber in mother ***/
1411  Double_t trans_z = -tot_thick_det / 2. + depth_in_det + dz; // (note: reference is CENTRE of absorber)
1412  TGeoTranslation* tr = new TGeoTranslation(0., 0., trans_z);
1413  mother_vol->AddNode(vol, 1, tr);
1414  }
1415  else {
1416  // single absorber: mother is absorber is detector is mother is ...
1417  mother_vol = vol;
1418  }
1419  depth_in_det += thick;
1420  no_abs++;
1421  }
1422 
1423  return mother_vol;
1424 }
1425 
1426 
1427 
1428 
1441 
1443 {
1444  // Construct and position a TGeoVolume shape to represent this detector in the current
1445  // geometry managed by gGeoManager.
1446  //
1447  // Adding the detector to the geometry requires:
1448  // - to construct a 'mother' volume (TGeoArb8) defined by the (theta-min/max, phi-min/max)
1449  // and the total thickness of the detector (all layers)
1450  // - to construct and position volumes (TGeoArb8) for each layer of the detector inside the
1451  // 'mother' volume
1452  // - to position 'mother' inside the top-level geometry
1453  //
1454  // gGeoManager must point to current instance of geometry manager.
1455 
1456  if (!gGeoManager) return;
1457 
1458  // get volume for detector
1459  TGeoVolume* vol = GetGeoVolume();
1460 
1461  // rotate detector to orientation corresponding to (theta,phi)
1462  Double_t theta = GetTheta();
1463  Double_t phi = GetPhi();
1464  TGeoRotation rot1, rot2;
1465  rot2.SetAngles(phi + 90., theta, 0.);
1466  rot1.SetAngles(-90., 0., 0.);
1467  // distance to detector centre = distance to telescope + depth of detector in telescope + half total thickness of detector
1468  KVTelescope* fTelescope = (KVTelescope*)GetParentStructure("TELESCOPE");
1469  Double_t dist_telescope = (fTelescope ? fTelescope->GetDistance() : 0.);
1470  Double_t dist = dist_telescope + fDepthInTelescope + fTotThickness / 2.;
1471  // translate detector to correct distance from target (note: reference is CENTRE of detector)
1472  Double_t trans_z = dist;
1473 
1474  TGeoTranslation tran(0, 0, trans_z);
1475  TGeoHMatrix h = rot2 * tran * rot1;
1476  TGeoHMatrix* ph = new TGeoHMatrix(h);
1477 
1478  // add detector volume to geometry
1479  gGeoManager->GetTopVolume()->AddNode(vol, 1, ph);
1480 }
1481 
1482 
1483 
1508 
1510 {
1511  // Set bitmask used to determine which acquisition parameters are
1512  // taken into account by KVDetector::Fired based on the environment variables
1513  // [dataset].KVACQParam.[par name].Working: NO
1514  // [dataset].[classname].Fired.ACQParameterList.[type]: PG,GG,T
1515  // where [classname]=KVDetector by default, or the name of some class
1516  // derived from KVDetector which calls the method KVDetector::SetKVDetectorFiredACQParameterListFormatString()
1517  // in its constructor.
1518  // The first allows to define certain acquisition parameters as not functioning;
1519  // they will not be taken into account.
1520  // The second allows to "fine-tune" what is meant by "all" or "any" acquisition parameters
1521  // (i.e. when using Fired("all"), Fired("any"), Fired("Pall", etc.).
1522  // For each detector type, give a comma-separated list of the acquisition
1523  // parameter types to be taken into account in the KVDetector::Fired method.
1524  // Only those parameters which appear in the list will be considered:
1525  // then "all" means => all parameters in the list
1526  // and "any" means => any of the parameters in the list
1527  // These lists are read during construction of multidetector arrays (KVMultiDetArray::Build),
1528  // the method KVMultiDetArray::SetACQParams uses them to define a mask for each detector
1529  // of the array.
1530  // Bits are set/reset in the order of the acquisition parameter list of the detector.
1531  // If no variable [dataset].[classname].Fired.ACQParameterList.[type] exists,
1532  // we set a bitmask authorizing all acquisition parameters of the detector, e.g.
1533  // if the detector has 3 acquisition parameters the bitmask will be "111"
1534 
1535  TObjArray* toks = lpar.Tokenize(",");
1536  TIter next(fACQParams);
1537  Bool_t no_variable_defined = (toks->GetEntries() == 0);
1538  KVACQParam* par;
1539  Int_t id = 0;
1540  while ((par = (KVACQParam*)next())) {
1541  if (!par->IsWorking()) fFiredMask.ResetBit(id); // ignore non-working parameters
1542  else {
1543  if (no_variable_defined || toks->FindObject(par->GetType())) fFiredMask.SetBit(id);
1544  else fFiredMask.ResetBit(id);
1545  }
1546  id++;
1547  }
1548  delete toks;
1549 }
1550 
1551 
1552 
1554 
1556 {
1557  cout << "(" << v.X() << "," << v.Y() << "," << v.Z() << ")";
1558 };
1559 
1560 
1561 
1564 
1566 {
1567  // Return surface area of first layer of detector in cm2.
1568 
1569  if (ROOTGeo()) {
1570  if (!GetEntranceWindow().GetShape()->InheritsFrom("TGeoArb8")) {
1571  // simple shape area
1572  return GetEntranceWindow().GetShape()->GetFacetArea(1);
1573  }
1574  // Monte Carlo calculation for TGeoArb8 shapes
1575  return GetEntranceWindow().GetSurfaceArea();
1576  }
1577 
1578  KVTelescope* fTelescope = (KVTelescope*)GetParentStructure("TELESCOPE");
1579  if (fTelescope && fDepthInTelescope == 0)
1580  fDepthInTelescope = fTelescope->GetDepthInCM(fTelescope->GetDetectorRank(this));
1581 
1582  TVector3 coords[4];
1583 
1584  if (fTelescope) fTelescope->GetCornerCoordinates(coords, fDepthInTelescope);
1585  else GetCornerCoordinates(coords, 0);
1586  cout << "DETECTOR COORDINATES (in cm):" << endl;
1587  cout << "=================================" << endl;
1588  cout << " A : ";
1589  printvec(coords[0]);
1590  cout << endl;
1591  cout << " B : ";
1592  printvec(coords[1]);
1593  cout << endl;
1594  cout << " C : ";
1595  printvec(coords[2]);
1596  cout << endl;
1597  cout << " D : ";
1598  printvec(coords[3]);
1599  cout << endl;
1600 
1601  cout << "DETECTOR DIMENSIONS (in cm):" << endl;
1602  cout << "================================" << endl;
1603  Double_t c = (coords[0] - coords[1]).Mag();
1604  Double_t b = (coords[1] - coords[2]).Mag();
1605  Double_t d = (coords[2] - coords[3]).Mag();
1606  Double_t a = (coords[0] - coords[3]).Mag();
1607  cout << " AB = " << c << endl;
1608  cout << " BC = " << b << endl;
1609  cout << " CD = " << d << endl;
1610  cout << " AD = " << a << endl;
1611 
1612  cout << "DETECTOR SURFACE AREA = ";
1613  Double_t surf = pow((a + b), 2.0) * (a - b + 2.0 * c) * (b - a + 2.0 * c);
1614  surf = sqrt(surf) / 400.0;
1615  cout << surf << " cm2" << endl;
1616 
1617  return surf;
1618 }
1619 
1620 
1621 
1626 
1628 {
1629  // Return pointer toTF1 giving residual energy after detector as function of incident energy,
1630  // for a given nucleus (Z,A).
1631  // The TF1::fNpx parameter is taken from environment variable KVDetector.ResidualEnergy.Npx
1632 
1633  if (!fEResF) {
1634  fEResF = new TF1(Form("KVDetector:%s:ERes", GetName()), this, &KVDetector::EResDet,
1635  0., 1.e+04, 2, "KVDetector", "EResDet");
1636  fEResF->SetNpx(gEnv->GetValue("KVDetector.ResidualEnergy.Npx", 20));
1637  }
1639  fEResF->SetRange(0., GetSmallestEmaxValid(Z, A));
1640  fEResF->SetTitle(Form("Residual energy [MeV] after detector %s for Z=%d A=%d", GetName(), Z, A));
1641 
1642  return fEResF;
1643 }
1644 
1645 
1646 
1651 
1653 {
1654  // Return pointer toTF1 giving range (in centimetres) in detector as function of incident energy,
1655  // for a given nucleus (Z,A).
1656  // The TF1::fNpx parameter is taken from environment variable KVDetector.Range.Npx
1657 
1658  if (!fRangeF) {
1659  fRangeF = new TF1(Form("KVDetector:%s:Range", GetName()), this, &KVDetector::RangeDet,
1660  0., 1.e+04, 2, "KVDetector", "RangeDet");
1661  fRangeF->SetNpx(gEnv->GetValue("KVDetector.Range.Npx", 20));
1662  }
1665  fRangeF->SetTitle(Form("Range [cm] in detector %s for Z=%d A=%d", GetName(), Z, A));
1666 
1667  return fRangeF;
1668 }
1669 
1670 
1671 
1676 
1678 {
1679  // Return pointer to TF1 giving energy loss in active layer of detector as function of incident energy,
1680  // for a given nucleus (Z,A).
1681  // The TF1::fNpx parameter is taken from environment variable KVDetector.EnergyLoss.Npx
1682 
1683  if (!fELossF) {
1684  fELossF = new TF1(Form("KVDetector:%s:ELossActive", GetName()), this, &KVDetector::ELossActive,
1685  0., 1.e+04, 2, "KVDetector", "ELossActive");
1686  fELossF->SetNpx(gEnv->GetValue("KVDetector.EnergyLoss.Npx", 20));
1687  }
1690  fELossF->SetTitle(Form("Energy loss [MeV] in detector %s for Z=%d A=%d", GetName(), Z, A));
1691  return fELossF;
1692 }
1693 
1694 
1695 
1700 
1702 {
1703  // Overrides KVMaterial::GetEIncOfMaxDeltaE
1704  // Returns incident energy corresponding to maximum energy loss in the
1705  // active layer of the detector, for a given nucleus.
1706 
1707  return GetELossFunction(Z, A)->GetMaximumX();
1708 }
1709 
1710 
1711 
1716 
1718 {
1719  // Overrides KVMaterial::GetMaxDeltaE
1720  // Returns maximum energy loss in the
1721  // active layer of the detector, for a given nucleus.
1722 
1723  return GetELossFunction(Z, A)->GetMaximum();
1724 }
1725 
1726 
1727 
1732 
1734 {
1735  // Overrides KVMaterial::GetDeltaE
1736  // Returns energy loss of given nucleus in the active layer of the detector.
1737 
1738  // optimization for single-layer detectors
1739  if (fSingleLayer) {
1740  return fActiveLayer->GetDeltaE(Z, A, Einc);
1741  }
1742  return GetELossFunction(Z, A)->Eval(Einc);
1743 }
1744 
1745 
1746 
1750 
1752 {
1753  // Returns calculated total energy loss of ion in ALL layers of the detector.
1754  // This is just (Einc - GetERes(Z,A,Einc))
1755 
1756  return Einc - GetERes(Z, A, Einc);
1757 }
1758 
1759 
1760 
1765 
1767 {
1768  // Overrides KVMaterial::GetERes
1769  // Returns residual energy of given nucleus after the detector.
1770  // Returns 0 if Einc<=0
1771 
1772  if (Einc <= 0.) return 0.;
1773  Double_t eres = GetEResFunction(Z, A)->Eval(Einc);
1774  // Eres function returns -1000 when particle stops in detector,
1775  // in order for function inversion (GetEIncFromEres) to work
1776  if (eres < 0.) eres = 0.;
1777  return eres;
1778 }
1779 
1780 
1781 
1804 
1806 {
1807  // Overrides KVMaterial::GetIncidentEnergy
1808  // Returns incident energy corresponding to energy loss delta_e in active layer of detector for a given nucleus.
1809  // If delta_e is not given, the current energy loss in the active layer is used.
1810  //
1811  // By default the solution corresponding to the highest incident energy is returned
1812  // This is the solution found for Einc greater than the maximum of the dE(Einc) curve.
1813  // If you want the low energy solution set SolType = KVIonRangeTable::kEmin.
1814  //
1815  // WARNING: calculating the incident energy of a particle using only the dE in a detector
1816  // is ambiguous, as in general (and especially for very heavy ions) the maximum of the dE
1817  // curve occurs for Einc greater than the punch-through energy, therefore it is not always
1818  // true to assume that if the particle does not stop in the detector the required solution
1819  // is that for type=KVIonRangeTable::kEmax. For a range of energies between punch-through
1820  // and dE_max, the required solution is still that for type=KVIonRangeTable::kEmin.
1821  // If the residual energy of the particle is unknown, there is no way to know which is the
1822  // correct solution.
1823  //
1824  // WARNING 2
1825  // If the given energy loss in the active layer is greater than the maximum theoretical dE
1826  // for given Z & A, (dE > GetMaxDeltaE(Z,A)) then we return a NEGATIVE incident energy
1827  // corresponding to the maximum, GetEIncOfMaxDeltaE(Z,A)
1828 
1829  if (Z < 1) return 0.;
1830 
1831  Double_t DE = (delta_e > 0 ? delta_e : GetEnergyLoss());
1832 
1833  // If the given energy loss in the active layer is greater than the maximum theoretical dE
1834  // for given Z & A, (dE > GetMaxDeltaE(Z,A)) then we return a NEGATIVE incident energy
1835  // corresponding to the maximum, GetEIncOfMaxDeltaE(Z,A)
1836  if (DE > GetMaxDeltaE(Z, A)) return -GetEIncOfMaxDeltaE(Z, A);
1837 
1838  TF1* dE = GetELossFunction(Z, A);
1839  Double_t e1, e2;
1840  dE->GetRange(e1, e2);
1841  switch (type) {
1842  case kEmin:
1843  e2 = GetEIncOfMaxDeltaE(Z, A);
1844  break;
1845  case kEmax:
1846  e1 = GetEIncOfMaxDeltaE(Z, A);
1847  break;
1848  }
1849  int status;
1850  Double_t EINC = ProtectedGetX(dE, DE, status, e1, e2);
1851  if (status != 0) {
1852 // Warning("GetIncidentEnergy",
1853 // "In %s : Called for Z=%d A=%d with dE=%f sol_type %d",
1854 // GetName(), Z, A, DE, type);
1855 // Warning("GetIncidentEnergy",
1856 // "Max DeltaE in this case is %f MeV",
1857 // GetMaxDeltaE(Z, A));
1858 // Warning("GetIncidentEnergy",
1859 // "Between Einc limits [%f,%f] no solution found",
1860 // e1,e2);
1861 // Warning("GetIncidentEnergy",
1862 // "Returned value is -Einc for %s dE, %f",
1863 // (status>0 ? "maximum" : "minimum"), EINC);
1864  return -EINC;
1865  }
1866  return EINC;
1867 }
1868 
1869 
1870 
1876 
1878 {
1879  // Overrides KVMaterial::GetDeltaEFromERes
1880  //
1881  // Calculate energy loss in active layer of detGetAlignedDetector for nucleus (Z,A)
1882  // having a residual kinetic energy Eres (MeV)
1883 
1884  if (Z < 1 || Eres <= 0.) return 0.;
1885  Double_t Einc = GetIncidentEnergyFromERes(Z, A, Eres);
1886  if (Einc <= 0.) return 0.;
1887  return GetELossFunction(Z, A)->Eval(Einc);
1888 }
1889 
1890 
1891 
1898 
1900 {
1901  // Overrides KVMaterial::GetIncidentEnergyFromERes
1902  //
1903  // Calculate incident energy of nucleus from residual energy.
1904  //
1905  // Returns -1 if Eres is out of defined range of values
1906 
1907  if (Z < 1 || Eres <= 0.) return 0.;
1908  //return GetEResFunction(Z, A)->GetX(Eres);
1909  Int_t status;
1910  Double_t einc = KVBase::ProtectedGetX(GetEResFunction(Z, A), Eres, status);
1911  if (status != 0) {
1912  // problem with inversion - value out of defined range of function
1913  return -1;
1914  }
1915  return einc;
1916 }
1917 
1918 
1919 
1923 
1925 {
1926  // Returns the smallest maximum energy for which range tables are valid
1927  // for all absorbers in the detector, and given ion (Z,A)
1928 
1929  Double_t maxmin = -1.;
1930  TIter next(fAbsorbers);
1931  KVMaterial* abs;
1932  while ((abs = (KVMaterial*)next())) {
1933  if (maxmin < 0.) maxmin = abs->GetEmaxValid(Z, A);
1934  else {
1935  if (abs->GetEmaxValid(Z, A) < maxmin) maxmin = abs->GetEmaxValid(Z, A);
1936  }
1937  }
1938  return maxmin;
1939 }
1940 
1941 
1942 
1960 
1962 {
1963  // Create detector from text file in 'TEnv' format.
1964  //
1965  // Example:
1966  // ========
1967  //
1968  // Layer: Gold
1969  // Gold.Material: Au
1970  // Gold.AreaDensity: 200.*KVUnits::ug
1971  // +Layer: Gas1
1972  // Gas1.Material: C3F8
1973  // Gas1.Thickness: 5.*KVUnits::cm
1974  // Gas1.Pressure: 50.*KVUnits::mbar
1975  // Gas1.Active: yes
1976  // +Layer: Si1
1977  // Si1.Material: Si
1978  // Si1.Thickness: 300.*KVUnits::um
1979 
1980  TEnv fEnvFile(envrc);
1981 
1982  KVString layers(fEnvFile.GetValue("Layer", ""));
1983  layers.Begin(" ");
1984  while (!layers.End()) {
1985  KVString layer = layers.Next();
1986  KVString mat = fEnvFile.GetValue(Form("%s.Material", layer.Data()), "");
1987  KVString tS = fEnvFile.GetValue(Form("%s.Thickness", layer.Data()), "");
1988  KVString pS = fEnvFile.GetValue(Form("%s.Pressure", layer.Data()), "");
1989  KVString dS = fEnvFile.GetValue(Form("%s.AreaDensity", layer.Data()), "");
1990  Double_t thick, dens, press;
1991  thick = dens = press = 0.;
1992  KVMaterial* M = 0;
1993  if (pS != "" && tS != "") {
1994  press = (Double_t)gROOT->ProcessLineFast(Form("%s*1.e+12", pS.Data()));
1995  press /= 1.e+12;
1996  thick = (Double_t)gROOT->ProcessLineFast(Form("%s*1.e+12", tS.Data()));
1997  thick /= 1.e+12;
1998  M = new KVMaterial(mat.Data(), thick, press);
1999  }
2000  else if (tS != "") {
2001  thick = (Double_t)gROOT->ProcessLineFast(Form("%s*1.e+12", tS.Data()));
2002  thick /= 1.e+12;
2003  M = new KVMaterial(mat.Data(), thick);
2004  }
2005  else if (dS != "") {
2006  dens = (Double_t)gROOT->ProcessLineFast(Form("%s*1.e+12", dS.Data()));
2007  dens /= 1.e+12;
2008  M = new KVMaterial(dens, mat.Data());
2009  }
2010  if (M) {
2011  AddAbsorber(M);
2012  if (fEnvFile.GetValue(Form("%s.Active", layer.Data()), kFALSE)) SetActiveLayer(M);
2013  }
2014  }
2015 }
2016 
2017 
2018 
2019 
2037 
2039 {
2040  // Returns list of detectors (including this one) which are in geometrical aligment
2041  // with respect to the target position (assuming this detector is part of a multidetector
2042  // array described by KVMultiDetArray).
2043  //
2044  // By default the list is in the order starting from this detector and going towards the target
2045  // (direction=KVGroup::kBackwards).
2046  // Call with argument direction=KVGroup::kForwards to have the list of detectors in the order
2047  // "seen" by a particle flying out from the target and arriving in this detector.
2048  //
2049  // If this detector is not part of a KVMultiDetArray (i.e. we have no information on
2050  // its geometrical relation to other detectors), we return 0x0.
2051  //
2052  // The list pointers are stored in member variable fAlignedDetectors[] for rapid retrieval,
2053  // the lists will be deleted with this detector.
2054  //
2055  // See KVGroup::GetAlignedDetectors for more details.
2056 
2057  if (!GetGroup() || direction > 1) return 0x0;
2058  if (fAlignedDetectors[direction]) return fAlignedDetectors[direction];
2059  return (fAlignedDetectors[direction] = GetGroup()->GetAlignedDetectors(this, direction));
2060 }
2061 
2062 
2063 
2064 
2066 
2068 {
2069  if (!GetGroup() || direction > 1) return;
2070  if (fAlignedDetectors[direction]) fAlignedDetectors[direction] = 0;
2071 }
2072 
2073 
2074 
2078 
2080 {
2081  // WARNING: SAME AS KVDetector::GetLinearRange
2082  // Only linear range in centimetres is calculated for detectors!
2083  return GetLinearRange(Z, A, Einc);
2084 }
2085 
2086 
2087 
2093 
2095 {
2096  // Returns range of ion in centimetres in this detector,
2097  // taking into account all layers.
2098  // Note that for Einc > punch through energy, this range is no longer correct
2099  // (but still > total thickness of detector).
2100  return GetRangeFunction(Z, A)->Eval(Einc);
2101 }
2102 
2103 
2104 
2108 
2110 {
2111  // Returns energy (in MeV) necessary for ion (Z,A) to punch through all
2112  // layers of this detector
2113 
2114  if (fSingleLayer) {
2115  // Optimize calculation time for single-layer detector
2116  return fActiveLayer->GetPunchThroughEnergy(Z, A);
2117  }
2118  return GetRangeFunction(Z, A)->GetX(GetTotalThicknessInCM());
2119 }
2120 
2121 
2122 
2123 
2128 
2130 {
2131  // Creates and fills a TGraph with the punch through energy in MeV vs. Z for the given detector,
2132  // for Z=1-92. The mass of each nucleus is calculated according to the given mass formula
2133  // (see KVNucleus).
2134 
2135  TGraph* punch = new TGraph(92);
2136  punch->SetName(Form("KVDetpunchthrough_%s_mass%d", GetName(), massform));
2137  punch->SetTitle(Form("Simple Punch-through %s (MeV) (mass formula %d)", GetName(), massform));
2138  KVNucleus nuc;
2139  nuc.SetMassFormula(massform);
2140  for (int Z = 1; Z <= 92; Z++) {
2141  nuc.SetZ(Z);
2142  punch->SetPoint(Z - 1, Z, GetPunchThroughEnergy(nuc.GetZ(), nuc.GetA()));
2143  }
2144  return punch;
2145 }
2146 
2147 
2148 
2153 
2155 {
2156  // Creates and fills a TGraph with the punch through energy in MeV/nucleon vs. Z for the given detector,
2157  // for Z=1-92. The mass of each nucleus is calculated according to the given mass formula
2158  // (see KVNucleus).
2159 
2160  TGraph* punch = new TGraph(92);
2161  punch->SetName(Form("KVDetpunchthroughEsurA_%s_mass%d", GetName(), massform));
2162  punch->SetTitle(Form("Simple Punch-through %s (AMeV) (mass formula %d)", GetName(), massform));
2163  KVNucleus nuc;
2164  nuc.SetMassFormula(massform);
2165  for (int Z = 1; Z <= 92; Z++) {
2166  nuc.SetZ(Z);
2167  punch->SetPoint(Z - 1, Z, GetPunchThroughEnergy(nuc.GetZ(), nuc.GetA()) / nuc.GetA());
2168  }
2169  return punch;
2170 }
2171 
2172 
2173 
2174 
2176 
2178 {
2179  return (KVGroup*)GetParentStructure("GROUP");
2180 }
2181 
2182 
2183 
2185 
2187 {
2188  return (GetGroup() ? GetGroup()->GetNumber() : 0);
2189 }
2190 
2191 
2192 
2194 
2196 {
2197  fParentStrucList.Add(elem);
2198 }
2199 
2200 
2201 
2203 
2205 {
2206  fParentStrucList.Remove(elem);
2207 }
2208 
2209 
2210 
2214 
2216 {
2217  // Get parent geometry structure element of given type.
2218  // Give unique name of structure if more than one element of same type is possible.
2219  KVGeoStrucElement* el = 0;
2220  if (strcmp(name, "")) {
2222  el = (KVGeoStrucElement*)strucs->FindObject(name);
2223  delete strucs;
2224  }
2225  else
2227  return el;
2228 }
2229 
2230 
2231 /* KVGeoDNTrajectory* KVDetector::GetTrajectoryForReconstruction()
2232 {
2233  // Return pointer to trajectory to be used for reconstruction of a
2234  // particle stopping in this detector.
2235  // If only one trajectory going forwards from this detector exists,
2236  // it is returned by default.
2237  // If there are more than one possible trajectories, a choice is made:
2238  // - we choose the trajectory with the least 'unfired' detectors
2239 
2240 
2243 
2244  if(!GetNode()->GetNTraj()) return 0x0;
2245  Int_t ntrajfor = GetNode()->GetNTrajForwards();
2246  if(ntrajfor>1){
2247  KVGeoDNTrajectory* tr = GetNode()->GetForwardTrajectoryWithLeastUnfiredDetectors();
2248  if(!tr) tr = GetNode()->GetForwardTrajectoryWithMostFiredDetectors();
2249  return tr;
2250  }
2251  else if(ntrajfor==1){
2252  return (KVGeoDNTrajectory*)GetNode()->GetForwardTrajectories()->First();
2253  }
2254  return (KVGeoDNTrajectory*)GetNode()->GetTrajectories()->First();
2255 }
2256  */
2258 {
2259  // Set ROOT geometry global matrix transformation to coordinate frame of active layer volume
2260  SetMatrix(m);
2261 }
2262 
2263 
2264 
2267 
2269 {
2270  // Set ROOT geometry shape of active layer volume
2271  SetShape(s);
2272 }
2273 
2274 
2275 
2278 
2280 {
2281  // Set ROOT geometry global matrix transformation to coordinate frame of entrance window
2283 }
2284 
2285 
2286 
2289 
2291 {
2292  // Set ROOT geometry shape of entrance window
2294 }
2295 
2296 
2297 
2305 
2307 {
2308  // Overrides KVMaterial::SetThickness
2309  //
2310  // If ROOT geometry is defined, we modify the DZ thickness of the volume representing
2311  // this detector in accordance with the new thickness.
2312  //
2313  // This is only implemented for single-layer detectors with a simple shape.
2314 
2315  if (ROOTGeo() && fSingleLayer) {
2318  TGeoBBox* shape = (TGeoBBox*)pn->GetShape();
2319  TGeoShape* newshape = nullptr;
2320  // bad kludge - is there no better way to clone a shape and change its dZ?
2321  if (shape->IsA() == TGeoBBox::Class()) {
2322  newshape = new TGeoBBox(shape->GetDX(), shape->GetDY(), 0.5 * thick);
2323  }
2324  else if (shape->IsA() == TGeoTube::Class()) {
2325  TGeoTube* oldtube = static_cast<TGeoTube*>(shape);
2326  newshape = new TGeoTube(oldtube->GetRmin(), oldtube->GetRmax(), 0.5 * thick);
2327  }
2328  else {
2329  Error("SetThickness", "No implementation for %s (%s)", shape->IsA()->GetName(), GetName());
2330  }
2331  if (newshape) {
2332  pn->Align(nullptr, newshape);
2333  Info("SetThickness", "Modified ROOT geometry for %s: new thickness=%g cm", GetName(), thick);
2335  }
2336  }
2337  KVMaterial::SetThickness(thick);
2338 }
2339 
2340 
2341 
2347 
2349 {
2350  // Return kTRUE if the two detectors have the same internal structure, i.e.
2351  // - the same number of absorber layers
2352  // - in the same order
2353  // - with the same material & thickness
2354 
2355  int nabs = GetNumberOfAbsorberLayers();
2356  if (other->GetNumberOfAbsorberLayers() != nabs) return kFALSE;
2357  bool same = true;
2358  for (int iabs = 0; iabs < nabs; ++iabs) {
2359  KVMaterial* this_abs = GetAbsorber(iabs);
2360  KVMaterial* other_abs = other->GetAbsorber(iabs);
2361  if (!this_abs->IsType(other_abs->GetType())
2362  || this_abs->GetMass() != other_abs->GetMass()
2363  || this_abs->GetThickness() != other_abs->GetThickness())
2364  same = false;
2365  }
2366  return same;
2367 }
2368 
2369 
2370 
2379 
2381 {
2382  // Add a new KVDetectorSignalExpression to this detector:
2383  //
2384  // 'type' will be the name of the new signal.
2385  //
2386  // '_expr' is a mathematical expression using any of the known signals of the detector.
2387  //
2388  // If the expression is not valid, no signal will be created and method returns kFALSE.
2389 
2391  if (!ds->IsValid()) {
2392  delete ds;
2393  ds = nullptr;
2394  }
2395  else
2396  AddDetectorSignal(ds);
2397  return (ds != nullptr);
2398 }
2399 
2400 
int Int_t
unsigned int UInt_t
void printvec(TVector3 &v)
KVIonRangeTableMaterial * M
ClassImp(KVPartitionList) void KVPartitionList
Initialisation.
#define SafeDelete(p)
#define d(i)
#define b(i)
#define c(i)
#define e(i)
unsigned char UChar_t
char Char_t
const Bool_t kFALSE
bool Bool_t
short Short_t
double Double_t
float Float_t
const Bool_t kTRUE
const char Option_t
const Bool_t kIterBackward
R__EXTERN TEnv * gEnv
int type
R__EXTERN TGeoManager * gGeoManager
double pow(double, double)
#define gROOT
char * Form(const char *fmt,...)
void ResetBit(UChar_t)
Definition: Binary_t.h:516
void SetBit(UChar_t)
Definition: Binary_t.h:498
Wrapper signal for KVACQParam objects.
GANIL VXI/VME acquisition parameter.
Definition: KVACQParam.h:15
void SetPedestal(Float_t ped)
Definition: KVACQParam.h:98
Bool_t IsWorking() const
Definition: KVACQParam.h:126
Float_t GetPedestal() const
Definition: KVACQParam.h:102
Double_t GetData() const
Definition: KVACQParam.h:72
void SetDetector(KVDetector *kd)
Definition: KVACQParam.h:34
void Clear(Option_t *="")
Clear object properties : name, type/title, number, label.
Definition: KVACQParam.h:44
virtual Bool_t IsType(const Char_t *typ) const
Definition: KVBase.h:178
const Char_t * GetType() const
Definition: KVBase.h:170
virtual void SetType(const Char_t *str)
Definition: KVBase.h:166
Double_t ProtectedGetX(const TF1 *func, Double_t val, int &status, Double_t xmin=0.0, Double_t xmax=0.0) const
Definition: KVBase.cpp:1654
static TPluginHandler * LoadPlugin(const Char_t *base, const Char_t *uri="0")
Definition: KVBase.cpp:756
UInt_t GetNumber() const
Definition: KVBase.h:213
Calibrated output from detector.
virtual Bool_t IsAvailableFor(const KVNameValueList &params) const
Base class for all detector calibrations.
Definition: KVCalibrator.h:75
TString GetInputSignalType() const
Definition: KVCalibrator.h:200
void SetDetector(KVDetector *d)
Definition: KVCalibrator.h:208
TString GetOutputSignalType() const
Definition: KVCalibrator.h:204
Detector output from a mathematical combination of other signals.
Output signal data produced by a detector.
virtual Bool_t IsExpression() const
virtual Bool_t IsRaw() const
virtual void Reset()
Base class for detector geometry description, interface to energy-loss calculations.
Definition: KVDetector.h:121
virtual KVDetectorSignal * GetDetectorSignal(const TString &type) const
Definition: KVDetector.h:463
static KVDetector * MakeDetector(const Char_t *name, Float_t thick)
virtual const Char_t * GetArrayName()
Definition: KVDetector.cpp:471
virtual Bool_t IsOK() const
Definition: KVDetector.h:628
void SetMatrix(const TGeoHMatrix *m)
Definition: KVDetector.h:147
KVPosition fEWPosition
position of entrance window i.e. first volume in detector geometry
Definition: KVDetector.h:124
Bool_t AddDetectorSignalExpression(const TString &type, const KVString &_expr)
virtual Double_t GetMaxDeltaE(Int_t Z, Int_t A)
@ kIdentifiedParticle
Definition: KVDetector.h:138
@ kUnidentifiedParticle
Definition: KVDetector.h:137
KVUniqueNameList fParentStrucList
list of geometry structures which directly contain this detector
Definition: KVDetector.h:125
Int_t GetNumberOfAbsorberLayers() const
Definition: KVDetector.h:264
virtual void GetVerticesInOwnFrame(TVector3 *, Double_t, Double_t)
virtual Double_t GetERes(Int_t Z, Int_t A, Double_t Einc)
virtual Double_t GetELostByParticle(KVNucleus *, TVector3 *norm=0)
Definition: KVDetector.cpp:305
Double_t ELossActive(Double_t *x, Double_t *par)
TList * GetTelescopesForIdentification()
Definition: KVDetector.cpp:900
virtual void AddIDTelescope(TObject *idt)
Add ID telescope to list of telescopes to which detector belongs.
Definition: KVDetector.cpp:861
KVList * GetAlignedIDTelescopes()
Definition: KVDetector.cpp:875
KVList * fCalibrators
list of associated calibrator objects
Definition: KVDetector.h:197
KVGeoStrucElement * GetParentStructure(const Char_t *type, const Char_t *name="") const
KVGroup * GetGroup() const
void AddDetectorSignal(KVDetectorSignal *ds)
Definition: KVDetector.h:736
virtual Double_t GetRange(Int_t Z, Int_t A, Double_t Einc)
virtual void SetEnergyLoss(Double_t e) const
Definition: KVDetector.h:332
Bool_t ReplaceCalibrator(const Char_t *type, KVCalibrator *cal, const KVNameValueList &opts="")
Definition: KVDetector.cpp:584
virtual void SetPedestal(const Char_t *, Float_t)
Set value of pedestal associated to parameter with given name.
Definition: KVDetector.cpp:696
virtual Double_t GetTotalDeltaE(Int_t Z, Int_t A, Double_t Einc)
virtual TGeoVolume * GetGeoVolume()
virtual Double_t GetEnergy() const
Definition: KVDetector.h:308
virtual TF1 * GetEResFunction(Int_t Z, Int_t A)
Binary8_t fFiredMask
bitmask used by Fired to determine which parameters to take into account
Definition: KVDetector.h:208
virtual Double_t GetEIncOfMaxDeltaE(Int_t Z, Int_t A)
virtual TGraph * DrawPunchThroughEsurAVsZ(Int_t massform=KVNucleus::kBetaMass)
static Int_t fDetCounter
Definition: KVDetector.h:127
void AddAbsorber(KVMaterial *)
Definition: KVDetector.cpp:753
virtual Double_t GetEnergyLoss() const
Definition: KVDetector.h:328
virtual Double_t GetEResAfterDetector() const
Definition: KVDetector.h:580
virtual void SetACQParams()
Definition: KVDetector.cpp:802
KVList * fIDTelescopes
list of ID telescopes to which detector belongs
Definition: KVDetector.h:129
void SetShape(TGeoBBox *s)
Definition: KVDetector.h:151
void AddParentStructure(KVGeoStrucElement *elem)
Double_t GetTheta() const
Definition: KVDetector.h:690
Double_t fEResforEinc
used by GetIncidentEnergy & GetCorrectedEnergy
Definition: KVDetector.h:218
KVList * fIDTelAlign
list of ID telescopes made of this detector and all aligned detectors placed in front of it
Definition: KVDetector.h:130
virtual Int_t FindZmin(Double_t ELOSS=-1., Char_t mass_formula=-1)
TF1 * fELossF
parametric function dE in active layer vs. incident energy
Definition: KVDetector.h:214
TList * fAlignedDetectors[2]
stores lists of aligned detectors in both directions
Definition: KVDetector.h:219
virtual void Clear(Option_t *opt="")
Definition: KVDetector.cpp:713
const KVPosition & GetEntranceWindow() const
Definition: KVDetector.h:647
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:769
void ClearHits()
Definition: KVDetector.h:399
void remove_signal_for_calibrator(KVCalibrator *K)
Definition: KVDetector.cpp:816
KVList * fAbsorbers
list of absorbers making up the detector
Definition: KVDetector.h:200
TList * fIDTele4Ident
list of ID telescopes used for particle ID
Definition: KVDetector.h:131
Double_t GetGain() const
Definition: KVDetector.h:789
TF1 * fEResF
parametric function Eres residual energy after all layers of detector
Definition: KVDetector.h:215
virtual Double_t GetSmallestEmaxValid(Int_t Z, Int_t A)
virtual Float_t GetACQData(const Char_t *) const
Definition: KVDetector.cpp:666
Bool_t HasSameStructureAs(const KVDetector *) const
void SetThickness(Double_t thick)
virtual void AddToGeometry()
KVList * fACQParams
list of raw data parameters read from coders
Definition: KVDetector.h:198
UInt_t GetGroupNumber()
virtual Double_t GetIncidentEnergy(Int_t Z, Int_t A, Double_t delta_e=-1.0, enum SolType type=kEmax)
virtual void ReadDefinitionFromFile(const Char_t *)
Double_t GetTotalThicknessInCM()
Definition: KVDetector.h:274
Bool_t IsCalibrated() const
Definition: KVDetector.h:359
void AddACQParam(KVACQParam *)
Add given acquisition parameter to this detector.
Definition: KVDetector.cpp:627
virtual TF1 * GetELossFunction(Int_t Z, Int_t A)
void SetActiveLayer(KVMaterial *actif)
Definition: KVDetector.h:244
KVDetector()
default ctor
Definition: KVDetector.cpp:104
Bool_t BelongsToUnidentifiedParticle() const
Definition: KVDetector.h:498
Double_t RangeDet(Double_t *x, Double_t *par)
KVMaterial * GetActiveLayer() const
Definition: KVDetector.h:249
virtual void RemoveCalibrators()
Definition: KVDetector.cpp:838
virtual Double_t GetEntranceWindowSurfaceArea()
Return surface area of first layer of detector in cm2.
virtual TGraph * DrawPunchThroughEnergyVsZ(Int_t massform=KVNucleus::kBetaMass)
KVUniqueNameList fDetSignals
list of signals associated with detector
Definition: KVDetector.h:190
Double_t GetPhi() const
Definition: KVDetector.h:708
virtual Double_t GetPunchThroughEnergy(Int_t Z, Int_t A)
Int_t fUnidentP
temporary counters, determine state of identified/unidentified particle flags
Definition: KVDetector.h:142
TF1 * fRangeF
parametric function range of particles in detector
Definition: KVDetector.h:216
virtual void SetFiredBitmask(KVString &)
virtual Double_t GetParticleEIncFromERes(KVNucleus *, TVector3 *norm=0)
Definition: KVDetector.cpp:357
virtual void SetEResAfterDetector(Double_t e)
Definition: KVDetector.h:576
void SetActiveLayerMatrix(const TGeoHMatrix *)
virtual void Copy(TObject &obj) const
Definition: KVDetector.cpp:155
virtual KVACQParam * GetACQParam(const Char_t *) const
Definition: KVDetector.cpp:646
void RemoveParentStructure(KVGeoStrucElement *elem)
const KVSeqCollection & GetListOfDetectorSignals() const
Definition: KVDetector.h:732
virtual Float_t GetPedestal(const Char_t *) const
Access pedestal value associated to parameter with given name.
Definition: KVDetector.cpp:683
KVList * fParticles
list of particles hitting detector in an event
Definition: KVDetector.h:199
virtual TList * GetAlignedDetectors(UInt_t direction=1)
KVGeoDetectorNode * GetNode()
Definition: KVDetector.h:285
void SetEntranceWindowMatrix(const TGeoHMatrix *)
Set ROOT geometry global matrix transformation to coordinate frame of entrance window.
virtual Double_t GetDeltaE(Int_t Z, Int_t A, Double_t Einc)
Double_t EResDet(Double_t *x, Double_t *par)
Double_t fDepthInTelescope
used to store depth of detector in parent telescope
Definition: KVDetector.h:206
KVMaterial * fActiveLayer
The active absorber in the detector.
Definition: KVDetector.h:128
void AddHit(KVNucleus *part)
Definition: KVDetector.h:374
KVCalibrator * GetCalibrator(const Char_t *name, const Char_t *type) const
Definition: KVDetector.h:757
virtual Double_t GetIncidentEnergyFromERes(Int_t Z, Int_t A, Double_t Eres)
Int_t fIdentP
temporary counters, determine state of identified/unidentified particle flags
Definition: KVDetector.h:141
void SetKVDetectorFiredACQParameterListFormatString()
Definition: KVDetector.cpp:87
virtual Double_t GetDeltaEFromERes(Int_t Z, Int_t A, Double_t Eres)
void SetEntranceWindowShape(TGeoBBox *)
Set ROOT geometry shape of entrance window.
void ResetAlignedDetectors(UInt_t direction=1)
virtual void DetectParticle(KVNucleus *, TVector3 *norm=0)
Definition: KVDetector.cpp:232
void SetActiveLayerShape(TGeoBBox *)
Set ROOT geometry shape of active layer volume.
void SetAnalysed(Bool_t b=kTRUE)
Definition: KVDetector.h:416
TString fFName
dynamically generated full name of detector
Definition: KVDetector.h:196
virtual void Print(Option_t *option="") const
Definition: KVDetector.cpp:393
Bool_t fSingleLayer
=kTRUE if detector has a single absorber layer
Definition: KVDetector.h:225
virtual Double_t GetLinearRange(Int_t Z, Int_t A, Double_t Einc)
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:949
void init()
default initialisations
Definition: KVDetector.cpp:41
Bool_t AddCalibrator(KVCalibrator *cal, const KVNameValueList &opts="")
Definition: KVDetector.cpp:510
Double_t fTotThickness
used to store value calculated by GetTotalThicknessInCM
Definition: KVDetector.h:205
TGeoBBox * GetShape() const
Definition: KVDetector.h:159
virtual void SetMaterial(const Char_t *type)
Definition: KVDetector.cpp:206
Path taken by particles through multidetector geometry.
const KVSeqCollection * GetIDTelescopes() const
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:19
const KVGeoDNTrajectory * GetTrajectoryForReconstruction(const KVGeoDNTrajectory *t, const KVGeoDetectorNode *n) const
Definition: KVGroup.h:103
Extended TList class which owns its objects by default.
Definition: KVList.h:27
Description of physical materials used to construct detectors; interface to range tables.
Definition: KVMaterial.h:41
virtual void SetThickness(Double_t thick)
Definition: KVMaterial.cpp:354
virtual Double_t GetThickness() const
Definition: KVMaterial.cpp:380
virtual TGeoMedium * GetGeoMedium(const Char_t *="")
Definition: KVMaterial.cpp:971
virtual void Print(Option_t *option="") const
Show information on this material.
Definition: KVMaterial.cpp:571
virtual Double_t GetEmaxValid(Int_t Z, Int_t A)
virtual Double_t GetPunchThroughEnergy(Int_t Z, Int_t A)
virtual void SetMaterial(const Char_t *type)
Definition: KVMaterial.cpp:191
virtual Double_t GetERes(Int_t Z, Int_t A, Double_t Einc)
Definition: KVMaterial.cpp:816
virtual Double_t GetDeltaE(Int_t Z, Int_t A, Double_t Einc)
Definition: KVMaterial.cpp:662
virtual void Clear(Option_t *opt="")
Reset absorber - set energy lost by particles to zero.
Definition: KVMaterial.cpp:879
KVMaterial()
default ctor
Definition: KVMaterial.cpp:70
Double_t GetEffectiveThickness(TVector3 &norm, TVector3 &direction)
Definition: KVMaterial.cpp:533
virtual Double_t GetParticleEIncFromERes(KVNucleus *, TVector3 *norm=0)
Definition: KVMaterial.cpp:632
virtual Double_t GetLinearRange(Int_t Z, Int_t A, Double_t Einc)
Definition: KVMaterial.cpp:699
Double_t GetMass() const
Returns atomic mass of material. Will be isotopic mass if set.
Definition: KVMaterial.cpp:252
Handles lists of named parameters with different types, a list of KVNamedParameter objects.
void SetValue(const Char_t *name, value_type value)
const Char_t * GetStringValue(const Char_t *name) const
Bool_t HasParameter(const Char_t *name) const
Description of properties and kinematics of atomic nuclei.
Definition: KVNucleus.h:125
Int_t GetA() const
Definition: KVNucleus.cpp:799
void SetZ(Int_t z, Char_t mt=-1)
Definition: KVNucleus.cpp:704
void SetMassFormula(UChar_t mt)
Definition: KVNucleus.h:347
Int_t GetZ() const
Return the number of proton / atomic number.
Definition: KVNucleus.cpp:770
TVector3 * GetPInitial() const
Definition: KVParticle.h:687
TVector3 GetMomentum() const
Definition: KVParticle.h:565
KVNameValueList * GetParameters() const
Definition: KVParticle.h:735
void SetMomentum(const TVector3 &v)
Definition: KVParticle.h:535
Double_t GetEnergy() const
Definition: KVParticle.h:582
void SetKE(Double_t ecin)
Definition: KVParticle.cpp:230
void SetE0(TVector3 *e=0)
Definition: KVParticle.h:676
void SetIsDetected()
Definition: KVParticle.h:691
Double_t GetKE() const
Definition: KVParticle.h:575
void SetEnergy(Double_t e)
Definition: KVParticle.h:560
Base class used for handling geometry in a multidetector array.
Definition: KVPosition.h:90
virtual void SetShape(TGeoBBox *)
Definition: KVPosition.cpp:743
Bool_t ROOTGeo() const
Returns kTRUE if ROOT geometry is used, kFALSE if not.
Definition: KVPosition.cpp:598
void GetCornerCoordinates(TVector3 *, Double_t=0)
Definition: KVPosition.cpp:496
virtual Double_t GetDistance(void) const
Definition: KVPosition.h:189
virtual Double_t GetSurfaceArea(int npoints=100000) const
Definition: KVPosition.cpp:974
virtual void SetMatrix(const TGeoHMatrix *)
Definition: KVPosition.cpp:702
virtual TGeoBBox * GetShape() const
Definition: KVPosition.cpp:790
void GetCornerCoordinatesInOwnFrame(TVector3 *, Double_t=0)
Definition: KVPosition.cpp:544
KaliVeda extensions to ROOT collection classes.
virtual void Clear(Option_t *option="")
virtual Int_t GetSize() const
virtual TObject * At(Int_t idx) const
KVSeqCollection * GetSubListWithType(const Char_t *retvalue) const
virtual TObject * FindObjectByType(const Char_t *) const
virtual void Add(TObject *obj)
virtual TObject * Remove(TObject *obj)
Remove object from list.
virtual void Delete(Option_t *option="")
virtual TObject * FindObject(const char *name) const
Extension of ROOT TString class which allows backwards compatibility with ROOT v3....
Definition: KVString.h:72
void Begin(TString delim) const
Definition: KVString.cpp:562
void RBegin(TString delim) const
Definition: KVString.cpp:741
Bool_t End() const
Definition: KVString.cpp:625
KVString Next(Bool_t strip_whitespace=kFALSE) const
Definition: KVString.cpp:675
KVString RNext(Bool_t strip_whitespace=kFALSE) const
Definition: KVString.cpp:814
Associates two detectors placed one behind the other.
Definition: KVTelescope.h:35
Int_t GetDetectorRank(KVDetector *kvd)
returns position (1=front, 2=next, etc.) detector in the telescope structure
Double_t GetDepthInCM(Int_t ndet) const
Definition: KVTelescope.h:73
virtual void Add(TObject *obj)
Handle several calibrations valid for different Z ranges.
virtual void ls(Option_t *option="") const
virtual void Print(Option_t *option, const char *wildcard, Int_t recurse=1) const
void SetName(const char *name)
virtual void AddAll(const TCollection *col)
virtual Int_t GetEntries() const
virtual const char * GetValue(const char *name, const char *dflt) const
virtual void SetTitle(const char *title="")
virtual void SetRange(Double_t xmin, Double_t xmax)
virtual void SetNpx(Int_t npx=100)
virtual Double_t GetX(Double_t y, Double_t xmin=0, Double_t xmax=0, Double_t epsilon=1.E-10, Int_t maxiter=100, Bool_t logx=false) const
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
virtual Double_t GetDX() const
virtual Double_t GetFacetArea(Int_t index=0) const
virtual Double_t GetDY() const
TGeoVolume * MakeArb8(const char *name, TGeoMedium *medium, Double_t dz, Double_t *vertices=0)
void RefreshPhysicalNodes(Bool_t lock=kTRUE)
TObjArray * GetListOfPhysicalNodes()
TGeoVolumeAssembly * MakeVolumeAssembly(const char *name)
TGeoVolume * GetTopVolume() const
TGeoPhysicalNode * MakePhysicalNode(const char *path=0)
virtual Int_t GetDefaultColor() const
TGeoMaterial * GetMaterial() const
Bool_t Align(TGeoMatrix *newmat=0, TGeoShape *newshape=0, Bool_t check=kFALSE, Double_t ovlp=0.001)
TGeoShape * GetShape(Int_t level=-1) const
void SetAngles(Double_t phi, Double_t theta, Double_t psi)
virtual const char * GetName() const
virtual Double_t GetRmin() const
virtual Double_t GetRmax() const
virtual void SetLineColor(Color_t lcolor)
virtual void AddNode(TGeoVolume *vol, Int_t copy_no, TGeoMatrix *mat=0, Option_t *option="")
virtual void SetPoint(Int_t i, Double_t x, Double_t y)
virtual void SetName(const char *name="")
virtual void SetTitle(const char *title="")
virtual void Add(TObject *obj)
virtual void Clear(Option_t *option="")
virtual const char * GetName() const
virtual void SetName(const char *name)
virtual TObject * Clone(const char *newname="") const
Int_t GetEntries() const
virtual TObject * FindObject(const char *name) const
R__ALWAYS_INLINE Bool_t TestBit(UInt_t f) const
virtual const char * ClassName() const
virtual void Warning(const char *method, const char *msgfmt,...) const
virtual Bool_t InheritsFrom(const char *classname) const
virtual void Error(const char *method, const char *msgfmt,...) const
void ResetBit(UInt_t f)
virtual void Info(const char *method, const char *msgfmt,...) const
Longptr_t ExecPlugin(int nargs, const T &... params)
TObjArray * Tokenize(const TString &delim) const
const char * Data() const
Bool_t Contains(const char *pat, ECaseCompare cmp=kExact) const
Double_t Y() const
Double_t X() const
RooCmdArg ClassName(const char *name)
VecExpr< UnaryOp< Sqrt< T >, SVector< T, D >, T >, T, D > sqrt(const SVector< T, D > &rhs)
T Mag(const SVector< T, D > &rhs)
Double_t x[n]
gr SetName("gr")
TH1 * h
const long double s
Definition: KVUnits.h:94
const long double m
Definition: KVUnits.h:70
double T(double x)
double dist(AxisAngle const &r1, AxisAngle const &r2)
Int_t Nint(T x)
constexpr Double_t K()
v
auto * a