KaliVeda
Toolkit for HIC analysis
KVReconstructedNucleus.cpp
1 #include "KVReconstructedNucleus.h"
2 #include "KVIDTelescope.h"
3 #include "KVGroup.h"
4 #include "KVMultiDetArray.h"
5 
6 using namespace std;
7 
9 
10 
13 
15 {
16  //default initialisation
17  fReconTraj = nullptr;
18  fRealZ = fRealA = 0.;
19  fDetNames = "/";
20  fIDTelName = "";
21  fIDTelescope = nullptr;
22  fNSegDet = 0;
23  fAnalStatus = 99;
24  fTargetEnergyLoss = 0;
25  ResetBit(kIsIdentified);
26  ResetBit(kIsCalibrated);
27  ResetBit(kCoherency);
28  ResetBit(kZMeasured);
29  ResetBit(kAMeasured);
30 }
31 
32 
33 
34 
37 
38 KVReconstructedNucleus::KVReconstructedNucleus() : fIDResults("KVIdentificationResult", 5)
39 {
40  //default ctor.
41  init();
42 }
43 
44 
45 
48 
50  : KVNucleus(), fIDResults("KVIdentificationResult", 5)
51 {
52  //copy ctor
53  init();
54 #if ROOT_VERSION_CODE >= ROOT_VERSION(3,4,0)
55  obj.Copy(*this);
56 #else
57  ((KVReconstructedNucleus&) obj).Copy(*this);
58 #endif
59 }
60 
61 
62 
63 
67 
69 {
70  // Stream an object of class KVReconstructedNucleus.
71  //Customized streamer.
72 
73  UInt_t R__s, R__c;
74  if (R__b.IsReading()) {
75  Version_t R__v = R__b.ReadVersion(&R__s, &R__c);
76  if (R__v < 17) {
77  // Before v17, fIDResults was a static array: KVIdentificationResult fIDresults[5]
78  // We convert this to the new TClonesArray format
79  KVNucleus::Streamer(R__b);
80  fDetNames.Streamer(R__b);
81  fIDTelName.Streamer(R__b);
82  if (R__v >= 16) {
83  R__b >> fNSegDet;
84  R__b >> fAnalStatus;
85  }
86  R__b >> fRealZ;
87  R__b >> fRealA;
88  R__b >> fTargetEnergyLoss;
89  int R__i;
90  KVIdentificationResult id_array[5];
91  for (R__i = 0; R__i < 5; R__i++) {
92  id_array[R__i].Streamer(R__b);
93  id_array[R__i].Copy(*GetIdentificationResult(R__i + 1));
94  }
95  R__b.CheckByteCount(R__s, R__c, KVReconstructedNucleus::IsA());
96  }
97  else
98  R__b.ReadClassBuffer(KVReconstructedNucleus::Class(), this, R__v, R__s, R__c);
99  // if the multidetector object exists, update some informations
100  // concerning the detectors etc. hit by this particle
101  if (gMultiDetArray) {
102  fReconTraj = nullptr;
104  fIDTelescope = nullptr;
105  if (fIDTelName != "") fIDTelescope = gMultiDetArray->GetIDTelescope(fIDTelName.Data());
106  if (fReconTraj) {
107  if (R__v < 16) fNSegDet = fReconTraj->GetNumberOfIndependentIdentifications(); // fNSegDet/fAnalStatus non-persistent before v.16
108  }
109  else {
110  TIter next_det(&fDetList);
111  KVDetector* det;
112  while ((det = (KVDetector*)next_det())) {
113  det->AddHit(this);
114  if (det->IsDetecting()) { //to be coherent with AddDetector() method
115  //modify detector's counters depending on particle's identification state
116  if (IsIdentified())
117  det->IncrementIdentifiedParticles();
118  else
119  det->IncrementUnidentifiedParticles();
120  }
121  }
122  }
123  }
124  }
125  else {
127  }
128 }
129 
130 
131 
132 
134 
136 {
137  switch (GetStatus()) {
138  case kStatusOK:
139  cout <<
140  "Particle alone in group, or identification independently of other particles in group is directly possible." << endl;
141  break;
142 
143  case kStatusOKafterSub:
144  cout <<
145  "Particle reconstructed after identification of others in group and subtraction of their calculated energy losses in common detectors."
146  << endl;
147  break;
148 
149  case kStatusOKafterShare:
150  cout <<
151  "Particle identification estimated after arbitrary sharing of energy lost in common detectors between several reconstructed particles."
152  << endl;
153  break;
154 
156  cout <<
157  "Particle stopped in first stage of telescope. Estimation of minimum Z."
158  << endl;
159  break;
160 
161  case kStatusPileupDE:
162  cout <<
163  "Undetectable pile-up in first member of identifying telesscope (apparent status=OK). Would lead to incorrect identification by DE-E method (Z and/or A overestimated)."
164  << endl;
165  break;
166 
167  case kStatusPileupGhost:
168  cout <<
169  "Undetectable ghost particle in filtered simulation. Another particle passed through all of the same detectors (pile-up)."
170  << endl;
171  break;
172 
173 
174  default:
175  cout << GetStatus() << endl;
176  break;
177  }
178 }
179 
180 
181 
184 
186 {
187  // Returns kTRUE if particle was detected in array with given name
188  return GetArrayName() == name;
189 }
190 
191 
192 
195 
197 {
198  // Returns name of array particle was detected in (if known)
199  if (GetParameters()->HasStringParameter("ARRAY")) return GetParameters()->GetStringValue("ARRAY");
200  if (GetStoppingDetector()) return GetStoppingDetector()->GetGroup()->GetParentStructure<KVMultiDetArray>()->GetName();
201  return "";
202 }
203 
204 
205 
207 
209 {
210 
211  int ndets = GetNumDet();
212  if (ndets) {
213 
214  for (int i = ndets - 1; i >= 0; i--) {
215  KVDetector* det = GetDetector(i);
216  if (det) det->Print("data");
217  }
218  for (int i = 1; i <= GetNumberOfIdentificationResults(); i++) {
220  if (idr->IDattempted) idr->Print();
221  }
222  }
223  if (GetStoppingDetector()) cout << "STOPPED IN : " <<
224  GetStoppingDetector()->GetName() << endl;
225  if (IsIdentified()) {
226  if (GetIdentifyingTelescope()) cout << "IDENTIFIED IN : " <<
227  GetIdentifyingTelescope()->GetName() << endl;
228  cout << " =======> ";
229  cout << " Z=" << GetZ() << " A=" << GetA();
230  if (IsAMeasured()) cout << " Areal=" << GetRealA();
231  else cout << " Zreal=" << GetRealZ();
232  }
233  else {
234  cout << "(unidentified)" << endl;
235  }
236  if (IsCalibrated()) {
237  cout << " Total Energy = " << GetEnergy() << " MeV, Theta=" << GetTheta() << " Phi=" << GetPhi() << endl;
238  cout << " Target energy loss correction : " << GetTargetEnergyLoss() << " MeV" << endl;
239  }
240  else {
241  cout << "(uncalibrated)" << endl;
242  }
243  cout << "RECONSTRUCTION STATUS : " << endl;
245  if (fReconTraj) fReconTraj->ls();
246  if (GetParameters()->GetNpar()) GetParameters()->Print();
247 }
248 
249 
250 
251 #if ROOT_VERSION_CODE >= ROOT_VERSION(3,4,0)
252 
257 
259 #else
261 #endif
262 {
263  //
264  //Copy this to obj
265  //
266  KVNucleus::Copy(obj);
268  robj.SetIdentifyingTelescope(GetIdentifyingTelescope());
269  robj.fDetNames = fDetNames;
270  robj.SetRealZ(GetRealZ());
271  robj.SetRealA(GetRealA());
272  robj.SetTargetEnergyLoss(GetTargetEnergyLoss());
273  robj.fReconTraj = fReconTraj;
274  robj.fIDTelescope = fIDTelescope;
275  robj.fNSegDet = fNSegDet;
276  robj.fAnalStatus = fAnalStatus;
277  // copy id results
278  Int_t nidres = GetNumberOfIdentificationResults();
279  for (int i = 1; i <= nidres; ++i) {
280  GetIdentificationResult(i)->Copy(*robj.GetIdentificationResult(i));
281  }
282  robj.SetBit(kIsIdentified, TestBit(kIsIdentified));
283  robj.SetBit(kIsCalibrated, TestBit(kIsCalibrated));
284  robj.SetBit(kCoherency, TestBit(kCoherency));
285  robj.SetBit(kZMeasured, TestBit(kZMeasured));
286  robj.SetBit(kAMeasured, TestBit(kAMeasured));
287 }
288 
289 
290 
295 
297 {
298  // Copy all characteristics of 'other' and also change all references to
299  // 'other' to references to 'this' (i.e. in detectors hit by particle).
300  // 'other' will not be fully valid after this operation (shouldn't be used further)
301 
302  other->Copy(*this);
303  KVGeoDetectorNode* node;
304  const KVReconNucTrajectory* traj = other->GetReconstructionTrajectory();
305  traj->IterateFrom();
306  while ((node = traj->GetNextNode())) {
307  KVDetector* d = node->GetDetector();
308  d->GetHits()->Remove(const_cast<KVReconstructedNucleus*>(other));
309  d->GetHits()->Add(this);
310  }
311 }
312 
313 
314 
318 
320 {
321  // Reset nucleus. Calls KVNucleus::Clear.
322  // if opt!="N": Calls KVGroup::Reset for the group where it was reconstructed.
323 
324  KVNucleus::Clear(opt);
325  if (GetGroup() && strncmp(opt, "N", 1))
326  GetGroup()->Reset();
327  fDetList.Clear();
328  fIDResults.Clear("C");
329  init();
330 }
331 
332 
333 
337 
339 {
340  // this method modifies the reconstructed trajectory
341  // probably called during the identification coherency check
342  fReconTraj = t;
343  fDetNames = t->GetPathString();//to store/retrieve reconstruction trajectory on file
344 }
345 
346 
347 
350 
352 {
353  // Method called in initial reconstruction of particle.
354 
355  fReconTraj = t;
357  t->AddUnidentifiedParticle(0);//add 1 unidentified particle to each detector on trajectory
358  fDetNames = t->GetPathString();//to store/retrieve reconstruction trajectory on file
359 }
360 
361 
362 
375 
377 {
378  // Calculate angles theta and phi for nucleus based on the detectors on its reconstruction trajectory.
379  // The momentum is set using these angles, its mass and its kinetic energy.
380  //
381  // The detector with the smallest solid angle along the trajectory is the one which defines the
382  // angles for the reconstructed particle.
383  //
384  // The (optional) option string can be "random" or "mean":
385  //
386  // If "random" (default) the angles are drawn at random between the over the surface of the detector.
387  //
388  // If "mean" the (theta,phi) position of the centre of the detector is used to fix the nucleus' direction.
389 
390 
391  if (GetEnergy() <= 0.0)//don't try if particle has no correctly defined energy
392  return;
394  return;
395 
396  KVDetector* angle_det = nullptr;
397  double small_solid = 1.e+09;
400  while ((n = GetReconstructionTrajectory()->GetNextNode())) {
401  if (n->GetDetector()->GetSolidAngle() < small_solid) {
402  angle_det = n->GetDetector();
403  small_solid = n->GetDetector()->GetSolidAngle();
404  }
405  }
406  if (!strcmp(opt, "random")) {
407  //random angles
408  TVector3 dir = angle_det->GetRandomDirection("random");
409  SetMomentum(GetEnergy(), dir);
410  }
411  else {
412  //middle of telescope
413  TVector3 dir = angle_det->GetDirection();
414  SetMomentum(GetEnergy(), dir);
415  }
416 }
417 
418 
419 
423 
425 {
426  // Called by Streamer when reading in data
427  // The fDetNames string is used to associate the particle with its reconstruction trajectory
428 
429  if (gMultiDetArray) {
430  // for old data, detnames was written as "/DET_1/DET_2/..."
431  // trajectory paths are written as "DET_1/DET_2/..."
432  if (fDetNames[0] == '/') fDetNames.Remove(0, 1);
433  fDetNames.Begin("/");
434  KVString n_stop_det = fDetNames.Next();
435  KVDetector* stop_det = gMultiDetArray->GetDetector(n_stop_det);
436  if (stop_det) {
437  KVGroup* gr = stop_det->GetGroup();
438  gr->AddHit(this);//for backwards compatibility, group coherency analysis performed in KVReconstructedEvent::Streamer (old data)
439  fReconTraj = (KVReconNucTrajectory*)gr->FindReconTraj(fDetNames);
440  if (fReconTraj) {
443  while ((n = fReconTraj->GetNextNode())) {
444  n->GetDetector()->AddHit(this);
445  if (IsIdentified()) n->GetDetector()->IncrementIdentifiedParticles();
446  else n->GetDetector()->IncrementUnidentifiedParticles();
447  }
448  }
449  }
450  }
451 }
452 
453 
454 
460 
462 {
463  // Change the particle's reconstruction trajectory to a different one starting from the
464  // same stopping detector (and therefore in the same group).
465  //
466  // trajectory paths are written as "DET_1/DET_2/..."
467 
469  if (!new_traj) {
470  Error("ReplaceReconTraj", "Trajectory %s not found - meant to replace %s", traj_name.Data(), fReconTraj->GetTitle());
471  return;
472  }
473  // first iterate over existing trajectory and reset detectors
476  while ((n = fReconTraj->GetNextNode())) {
477  n->GetDetector()->RemoveHit(this);
478  if (IsIdentified()) n->GetDetector()->IncrementIdentifiedParticles(-1);
479  else n->GetDetector()->IncrementUnidentifiedParticles(-1);
480  }
481  // set new trajectory
482  fReconTraj = new_traj;
484  while ((n = fReconTraj->GetNextNode())) {
485  n->GetDetector()->AddHit(this);
486  if (IsIdentified()) n->GetDetector()->IncrementIdentifiedParticles();
487  else n->GetDetector()->IncrementUnidentifiedParticles();
488  }
489 }
490 
491 
492 
500 
502 {
503  // Set identification of nucleus from informations in identification result object
504  // The mass (A) information in KVIdentificationResult is only used if the mass
505  // was measured as part of the identification. Otherwise the nucleus' mass formula
506  // will be used to calculate A from the measured Z.
507  //
508  // The identifying telescope is set to idt.
509 
511  SetIDCode(idr->IDcode);
512  SetZMeasured(idr->Zident);
513  SetAMeasured(idr->Aident);
514  SetZ(idr->Z);
515  if (idr->Aident) {
516  SetA(idr->A);
517  SetRealA(idr->PID);
518  }
519  else {
520  SetRealZ(idr->PID);
521  }
522 }
523 
524 
525 
527 
529 {
530  printf(" A:%6s", GetParameters()->GetStringValue("ARRAY"));
531  if (GetStoppingDetector()) printf(" D:%10s", GetStoppingDetector()->GetName());
532  printf(" IDCODE=%2d", GetIDCode());
533  if (IsIdentified()) {
534  if (GetIdentifyingTelescope()) printf(" ID:%15s", GetIdentifyingTelescope()->GetName());
535  if (IsZMeasured()) printf(" Z=%2d", GetZ());
536  else printf(" ");
537  if (IsAMeasured()) printf(" A=%3d : ", GetA());
538  else printf(" : ");
539  if (IsCalibrated()) printf(" E=%g MeV", GetEnergy());
540  if (GetParameters()->IsValue("Coherent", false)) printf("/not coherent/");
541  if (GetParameters()->IsValue("Pileup", true)) printf("/pileup/");
542  }
543  printf("\n");
544 }
545 
546 
int Int_t
unsigned int UInt_t
#define d(i)
bool Bool_t
short Version_t
const char Option_t
char name[80]
KVGeoDetectorNode * GetNextNode() const
KVString GetPathString() const
void AddUnidentifiedParticle(int modify_identified=-1) const
void IterateFrom(const KVGeoDetectorNode *node0=nullptr) const
Information on relative positions of detectors & particle trajectories.
KVDetector * GetDetector() const
virtual KVDetector * GetDetector(const Char_t *name) const
Return detector in this structure with given name.
Group of detectors which can be treated independently of all others in array.
Definition: KVGroup.h:19
void Reset(Option_t *opt="")
Definition: KVGroup.cpp:63
const KVGeoDNTrajectory * FindReconTraj(const KVString &path)
Definition: KVGroup.h:92
Base class for all detectors or associations of detectors in array which can identify charged particl...
Definition: KVIDTelescope.h:84
Full result of one attempted particle identification.
Bool_t IDattempted
=kTRUE if identification was attempted
void Print(Option_t *opt="") const override
Bool_t Aident
= kTRUE if A of particle established
Double_t PID
= "real" Z if Zident==kTRUE and Aident==kFALSE, "real" A if Zident==Aident==kTRUE
Int_t A
A of particle found (if Aident==kTRUE)
Int_t Z
Z of particle found (if Zident==kTRUE)
Int_t IDcode
a general identification code for this type of identification
void Copy(TObject &) const override
Copy this to obj.
Bool_t Zident
=kTRUE if Z of particle established
Base class for describing the geometry of a detector array.
KVIDTelescope * GetIDTelescope(const Char_t *name) const
Return pointer to DeltaE-E ID Telescope with "name".
const Char_t * GetStringValue(const Char_t *name) const
void Print(Option_t *opt="") const override
Description of properties and kinematics of atomic nuclei.
Definition: KVNucleus.h:123
void Copy(TObject &) const override
Copy this KVNucleus into the KVNucleus object referenced by "obj".
Definition: KVNucleus.cpp:827
void Clear(Option_t *opt="") override
Definition: KVNucleus.cpp:291
Int_t GetA() const
Definition: KVNucleus.cpp:792
void SetA(Int_t a)
Definition: KVNucleus.cpp:647
void SetZ(Int_t z, Char_t mt=-1)
Definition: KVNucleus.cpp:696
Int_t GetZ() const
Return the number of proton / atomic number.
Definition: KVNucleus.cpp:763
KVNameValueList * GetParameters() const
Definition: KVParticle.h:818
Double_t GetTheta() const
Definition: KVParticle.h:683
Double_t GetEnergy() const
Definition: KVParticle.h:624
const Char_t * GetName() const override
return the field fName
Definition: KVParticle.cpp:455
Double_t GetPhi() const
Definition: KVParticle.h:691
void SetMomentum(const TVector3 *v)
Definition: KVParticle.h:545
Path through detector array used to reconstruct detected particle.
Int_t GetNumberOfIndependentIdentifications() const
void ls(Option_t *="") const override
Nuclei reconstructed from data measured by a detector array .
virtual void SetIDCode(UShort_t s)
void ls(Option_t *="") const override
KVDetector * GetDetector(const TString &label) const
virtual Bool_t IsZMeasured() const
virtual void SetAMeasured(Bool_t yes=kTRUE)
Float_t fRealZ
Z returned by identification routine.
void SetReconstructionTrajectory(const KVReconNucTrajectory *t)
Method called in initial reconstruction of particle.
virtual Double_t GetTargetEnergyLoss() const
void CopyAndMoveReferences(const KVReconstructedNucleus *)
@ kStatusOKafterShare
of energy losses of other particles in the same group which have Status=0
@ kStatusStopFirstStage
(arbitrarily) between this and the other particle(s) with Status=2
@ kStatusPileupDE
telescope; a minimum Z could be estimated from the measured energy loss.
Bool_t InArray(const TString &) const
Returns kTRUE if particle was detected in array with given name.
const KVReconNucTrajectory * fReconTraj
trajectory used to reconstruct particle
Int_t GetNumberOfIdentificationResults() const
KVIdentificationResult * GetIdentificationResult(Int_t i)
KVHashList fDetList
non-persistent list of pointers to detectors
const KVReconNucTrajectory * GetReconstructionTrajectory() const
Float_t fRealA
A returned by identification routine.
void SetIdentification(KVIdentificationResult *, KVIDTelescope *)
void Copy(TObject &) const override
void Clear(Option_t *option="") override
virtual Int_t GetIDCode() const
TClonesArray fIDResults
results of every identification attempt made for this nucleus, in order of the ID telescopes used
KVDetector * GetStoppingDetector() const
void SetDetector(int i, KVDetector *);
virtual void SetTargetEnergyLoss(Double_t e)
KVString fDetNames
list of names of detectors through which particle passed
void ModifyReconstructionTrajectory(const KVReconNucTrajectory *t)
virtual void GetAnglesFromReconstructionTrajectory(Option_t *opt="random")
void init()
default initialisation
KVString fIDTelName
name of identification telescope which identified this particle (if any)
TString GetArrayName() const
Returns name of array particle was detected in (if known)
void Print(Option_t *option="") const override
virtual void SetZMeasured(Bool_t yes=kTRUE)
Int_t fNSegDet
number of segmented/independent detectors hit by particle
void ReplaceReconTraj(const TString &traj_name)
Double_t fTargetEnergyLoss
calculated energy lost in target
KVIDTelescope * GetIdentifyingTelescope() const
void SetIdentifyingTelescope(KVIDTelescope *i)
KVIDTelescope * fIDTelescope
non-persistent pointer to identification telescope
virtual Bool_t IsAMeasured() const
Int_t fAnalStatus
status of particle after analysis of reconstructed event
void Clear(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
KVString Next(Bool_t strip_whitespace=kFALSE) const
Definition: KVString.cpp:695
virtual Version_t ReadVersion(UInt_t *start=nullptr, UInt_t *bcnt=nullptr, const TClass *cl=nullptr)=0
virtual Int_t CheckByteCount(UInt_t startpos, UInt_t bcnt, const char *classname)=0
virtual Int_t ReadClassBuffer(const TClass *cl, void *pointer, const TClass *onfile_class=nullptr)=0
Bool_t IsReading() const
virtual Int_t WriteClassBuffer(const TClass *cl, void *pointer)=0
void Clear(Option_t *option="") override
TClass * IsA() const override
void Streamer(TBuffer &) override
static TClass * Class()
const char * GetName() const override
void Streamer(TBuffer &) override
const char * GetTitle() const override
void SetBit(UInt_t f)
virtual void Error(const char *method, const char *msgfmt,...) const
const char * Data() const
virtual void Streamer(TBuffer &)
TString & Remove(EStripType s, char c)
const Int_t n
TGraphErrors * gr
ClassImp(TPyArg)