KaliVeda  1.13/01
Heavy-Ion Analysis Toolkit
KVINDRAGroupReconstructor.cpp
Go to the documentation of this file.
1 //Created by KVClassFactory on Wed Feb 21 13:42:47 2018
2 //Author: John Frankland,,,
3 
5 #include <KVIDGCsI.h>
6 #include <KVIDINDRACsI.h>
7 #include <KVINDRADetector.h>
8 
10 
12 
13 
27 
29 {
30  // \param traj trajectory currently being scanned
31  // \param node current detector on trajectory to test
32  // \returns pointer to a new reconstructed particle added to this group's event; nullptr if nothing is to be done
33  //
34  // Specialised particle reconstruction for INDRA data.
35  //
36  // If the fired detector in question is a CsI we check, if identification is available, whether
37  // this corresponds to a 'gamma'. If so we count it (event parameter "INDRA_GAMMA_MULT")
38  // and add the name of the detector to the parameter "INDRA_GAMMA_DETS"
39  // but do not begin the reconstruction of a particle. This allows to continue along the trajectory and
40  // directly reconstruct any charged particle which may stop in the Si detector in coincidence with
41  // a 'gamma' in the CsI.
42 
43  if (node->GetDetector()->IsType("CSI")) {
44 
45  if (node->GetDetector()->Fired(GetPartSeedCond())) {
46 
47  ++nfireddets;
49  if (idt) {
50 
52  if (idt->IsReadyForID()) {
53 
54  idt->Identify(&idr);
55  if (idr.IDOK && (idr.IDcode == KVINDRA::IDCodes::ID_GAMMA || idr.IDquality == KVIDGCsI::kICODE10)) { // gamma in CsI
56  GetEventFragment()->GetParameters()->IncrementValue("INDRA_GAMMA_MULT", 1);
57  GetEventFragment()->GetParameters()->IncrementValue("INDRA_GAMMA_DETS", node->GetName());
58  node->GetDetector()->SetAnalysed();
59  return nullptr;
60  }
61 
62  // if we arrive here, CSI identification for the particle has been performed and
63  // the result is not a gamma (which are rejected; no particle is reconstructed).
64  // as the coordinates in the identification map are randomized (KVACQParamSignal),
65  // we do not want to perform a second identification attempt in KVINDRAGroupReconstructor::IdentifyParticle:
66  // the results may not be the same, and for particles close to the threshold a charged particle
67  // identified here may be subsequently identified as a gamma in a second identification attempt,
68  // leading to the incoherency that there are gamma particles (IDcode=KVINDRA::IDCodes::ID_GAMMA)
69  // in the final reconstructed event.
70  // therefore in this case we add a new particle to the event, and copy the results of this
71  // first identification into the particle, with IDR->IDattempted=true, so that in
72  // KVINDRAGroupReconstructor::IdentifyParticle the CSI identification will not be redone.
73 
74  auto new_part = GetEventFragment()->AddParticle();
75  *(new_part->GetIdentificationResult(1)) = idr;
76  return new_part;
77  }
78  }
79  return GetEventFragment()->AddParticle();
80  }
81  return nullptr;
82  }
84 }
85 
86 
87 
89 
91 {
93  for (KVReconstructedEvent::Iterator it = GetEventFragment()->begin(); it != GetEventFragment()->end(); ++it) {
94 
95  KVReconstructedNucleus* d = it.get_pointer();
96 
97  if (!d->IsIdentified()) d->SetIDCode(KVINDRA::IDCodes::NO_IDENTIFICATION); // unidentifiable particle
98  else {
99  if (d->GetIDCode() == KVINDRA::IDCodes::ID_CSI_PSA && d->GetZ() == 0) {
100  std::cout << "\n\n";
101  std::cout << " GAMMA\n\n";
102  d->Print();
103  std::cout << "\n\n";
104  auto traj = (KVGeoDNTrajectory*)d->GetStoppingDetector()->GetNode()->GetTrajectories()->First();
105  KVIDINDRACsI* idt = (KVIDINDRACsI*)traj->GetIDTelescopes()->FindObjectByType(CSI_ID_TYPE);
106  std::cout << idt << std::endl;
107  std::cout << "Type: " << d->GetStoppingDetector()->GetType() << std::endl;
108  std::cout << "Fired: " << d->GetStoppingDetector()->Fired(GetPartSeedCond()) << std::endl;
109  std::cout << "ReadyforID: " << d->GetIdentifyingTelescope()->IsReadyForID() << std::endl;
110  std::cout << "\n\n";
111  exit(1);
112  }
113  }
114  }
115 }
116 
117 
118 
137 
139 {
140  // INDRA-specific particle identification.
141  // Here we attribute the Veda6-style general identification codes depending on the
142  // result of KVReconstructedNucleus::Identify and the subcodes from the different
143  // identification algorithms:
144  // If the particle's mass A was NOT measured, we make sure that it is calculated
145  // from the measured Z using the mass formula defined by default
146  //
147  //IDENTIFIED PARTICLES
148  //Identified particles with ID code = 2 with subcodes 4 & 5
149  //(masse hors limite superieure/inferieure) are relabelled
150  //with kIDCode10 (identification entre les lignes CsI)
151  //
152  //UNIDENTIFIED PARTICLES
153  //Unidentified particles receive the general ID code for non-identified particles (kIDCode14)
154  //EXCEPT if their identification in CsI gave subcodes 6 or 7
155  //(Zmin) then they are relabelled "Identified" with IDcode = 9 (ident. incomplete dans CsI ou Phoswich (Z.min))
156  //Their "identifying" telescope is set to the CsI ID telescope
157 
159  // if a successful identification has occured, identifying_telescope & partID are now set
160  // and the particle has IsIdentified()=true
161 
162  // INDRA coherency treatment
163  PART.SetParameter("Coherent", kTRUE);
164  PART.SetParameter("Pileup", kFALSE);
165  PART.SetParameter("UseFullChIoEnergyForCalib", kTRUE);
166  Bool_t ok = DoCoherencyAnalysis(PART);
167 
168  if (ok) { // identification may change here due to coherency analysis
169  PART.SetIsIdentified();
171  } // if not ok, do we need to unset any previously identified particle?
172 
173  if (PART.IsIdentified()) {
174 
175  /******* IDENTIFIED PARTICLES *******/
176  if (partID.IsType(CSI_ID_TYPE)) { /**** CSI IDENTIFICATION ****/
177 
178  //Identified particles with ID code = 2 with subcodes 4 & 5
179  //(masse hors limite superieure/inferieure) are relabelled
180  //with kIDCode10 (identification entre les lignes CsI)
181 
182  Int_t grid_code = partID.IDquality;
183  if (grid_code == KVIDGCsI::kICODE4 || grid_code == KVIDGCsI::kICODE5) {
184  partID.IDcode = KVINDRA::IDCodes::ID_CSI_MASS_OUT_OF_RANGE;
186  }
187 
188  }
189 
190  }
191  else {
192 
193  /******* UNIDENTIFIED PARTICLES *******/
194 
195  /*** general ID code for non-identified particles ***/
196  PART.SetIDCode(KVINDRA::IDCodes::NO_IDENTIFICATION);
197  auto csirl = id_by_type.find(CSI_ID_TYPE.Data());
198  if (csirl != id_by_type.end()) {
199  //Particles remaining unidentified are checked: if their identification in CsI R-L gave subcodes 6 or 7
200  //(Zmin) then they are relabelled "Identified" with IDcode = 9 (ident. incomplete dans CsI ou Phoswich (Z.min))
201  //
202  //Particles with ID code = 2 with subcodes 4 & 5 (masse hors limite superieure/inferieure) are relabelled
203  //with kIDCode10 (identification entre les lignes CsI)
204  //
205  //Their "identifying" telescope is set to the CsI ID telescope
206  if (csirl->second->IDattempted) {
207  if (csirl->second->IDquality == KVIDGCsI::kICODE6 || csirl->second->IDquality == KVIDGCsI::kICODE7) {
208  PART.SetIsIdentified();
209  csirl->second->IDcode = KVINDRA::IDCodes::ID_CSI_FRAGMENT;
210  partID = *(csirl->second);
213  }
214  else if (csirl->second->IDquality == KVIDGCsI::kICODE4 || csirl->second->IDquality == KVIDGCsI::kICODE5) {
215  PART.SetIsIdentified();
216  csirl->second->IDcode = KVINDRA::IDCodes::ID_CSI_MASS_OUT_OF_RANGE;
217  partID = *(csirl->second);
220  }
221  }
222  }
223 
224  }
225 }
226 
227 
228 
235 
237 {
238  // calculate fEChIo from residual energy
239  //
240  // returns kFALSE if it doesn't work, and sets particle bad calibration status
241  //
242  // returns kTRUE if it works, and sets calib status to SOME_ENERGY_LOSSES_CALCULATED
243 
244  Double_t e0 = theChio->GetDeltaEFromERes(n->GetZ(), n->GetA(), ERES);
247  if (fEChIo <= 0) {
249  fEChIo = 0;
250  return kFALSE;
251  }
253  n->SetECode(KVINDRA::ECodes::SOME_ENERGY_LOSSES_CALCULATED);
254  return kTRUE;
255 }
256 
257 
258 
264 
266 {
267  // Calculate and set the energy of a (previously identified) reconstructed particle
268  //
269  // This is only possible for correctly identified particles.
270  // We exclude IDCODE9 particles (Zmin in CsI-RL)
271 
272  fEChIo = fESi = fECsI = 0;
273 
274  print_part = false;
275 
277 
278  if (PART->GetIDCode() != KVINDRA::IDCodes::ID_CSI_FRAGMENT && PART->GetIDCode() != KVINDRA::IDCodes::ID_CSI_MASS_OUT_OF_RANGE) {
279  // this status may be modified depending on what happens in DoCalibration
280  SetCalibrationStatus(*PART, KVINDRA::ECodes::NORMAL_CALIBRATION);
281  DoCalibration(PART);
282  }
283 
284  PART->SetParameter("INDRA.ECHIO", fEChIo);
285  PART->SetParameter("INDRA.ESI", fESi);
286  PART->SetParameter("INDRA.ECSI", fECsI);
287  //add correction for target energy loss - moving charged particles only!
288  Double_t E_targ = 0.;
289  if (PART->GetZ() && PART->GetEnergy() > 0) {
290  E_targ = GetTargetEnergyLossCorrection(PART);
291  PART->SetTargetEnergyLoss(E_targ);
292  }
293  Double_t E_tot = PART->GetEnergy() + E_targ;
294  PART->SetEnergy(E_tot);
295  // set particle momentum from telescope dimensions (random)
297  CheckCsIEnergy(PART);
298 
299  //if(print_part) PART->Print();
300 }
301 
302 
303 
309 
311 {
312  // Beryllium-8 = 2 alpha particles of same energy
313  // We halve the total light output of the CsI to calculate the energy of 1 alpha
314  // Then multiply resulting energy by 2
315  // Note: fECsI is -ve, because energy is calculated not measured
316 
317  auto csi = GetCsI(n);
318  Double_t half_light = csi->GetDetectorSignalValue("TotLight") * 0.5;
319  KVNucleus tmp(2, 4);
320  double ecsi = 2.*csi->GetCorrectedEnergy(&tmp, half_light, kFALSE);
321  if (ecsi > 0) {
322  SetCalibrationStatus(*n, KVINDRA::ECodes::SOME_ENERGY_LOSSES_CALCULATED);
323  // calculated energy returned as negative value
324  return -ecsi;
325  }
327  return 0;
328 }
329 
330 
331 
338 
340 {
341  // Check calculated CsI energy loss of particle.
342  // If it is greater than the maximum theoretical energy loss
343  // (depending on the length of CsI, the Z & A of the particle)
344  // we set the energy calibration code to 3 (historical VEDA code
345  // for particles with E_csi > E_max_csi)
346 
347  KVDetector* csi = GetCsI(n);
348  if (csi && (n->GetZ() > 0) && (n->GetZ() < 3) && (csi->GetDetectorSignalValue("Energy", Form("Z=%d,A=%d", n->GetZ(), n->GetA())) > csi->GetMaxDeltaE(n->GetZ(), n->GetA()))) {
349  n->SetECode(KVINDRA::ECodes::WARNING_CSI_MAX_ENERGY);
350  }
351 }
352 
353 
354 /*
355  The following is a copy of the old
356 void KVINDRAReconEvent::SecondaryAnalyseGroup(KVGroup* grp)
357  method from branch 1.10
358 
359 void KVINDRAGroupReconstructor:PerformSecondaryAnalysis()
360 {
361  // Perform identifications and calibrations of particles not included
362  // in first round (methods IdentifyEvent() and CalibrateEvent()).
363  //
364  // Here we treat particles with GetStatus()==KVReconstructedNucleus::kStatusOKafterSub
365  // after subtracting the energy losses of all previously calibrated particles in group from the
366  // measured energy losses in the detectors they crossed.
367 
368  // loop over al identified & calibrated particles in group and subtract calculated
369  // energy losses from all detectors
370  KVINDRAReconNuc* nuc;
371  TList sixparts;
372  TIter parts(grp->GetParticles());
373  while ((nuc = (KVINDRAReconNuc*)parts())) {
374  if (nuc->IsIdentified() && nuc->IsCalibrated()) {
375  nuc->SubtractEnergyFromAllDetectors();
376  // reconstruct particles from pile-up in silicon detectors revealed by coherency CsIR/L - SiCsI
377  if (nuc->IsSiPileup() && nuc->GetSi()->GetEnergy() > 0.1) {
378  KVINDRAReconNuc* SIX = AddParticle();
379  SIX->Reconstruct(nuc->GetSi());
380  sixparts.Add(SIX);
381  }
382  // reconstruct particles from pile-up in si75 detectors revealed by coherency
383  if (nuc->IsSi75Pileup()) {
384  KVINDRAReconNuc* SIX = AddParticle();
385  SIX->Reconstruct(nuc->GetSi75());
386  sixparts.Add(SIX);
387  }
388  // reconstruct particles from pile-up in sili detectors revealed by coherency
389  if (nuc->IsSiLiPileup()) {
390  KVINDRAReconNuc* SIX = AddParticle();
391  SIX->Reconstruct(nuc->GetSiLi());
392  sixparts.Add(SIX);
393  }
394 
395  // reconstruct particles from pile-up in ChIo detectors revealed by coherency CsIR/L - ChIoCsI
396  if (nuc->IsChIoPileup() && nuc->GetChIo()->GetEnergy() > 1.0) {
397  KVINDRAReconNuc* SIX = AddParticle();
398  SIX->Reconstruct(nuc->GetChIo());
399  sixparts.Add(SIX);
400  }
401  }
402  }
403  // reanalyse group
404  KVReconstructedNucleus::AnalyseParticlesInGroup(grp);
405 
406  Int_t nident = 0; //number of particles identified in each step
407  if (sixparts.GetEntries()) { // identify any particles added by coherency CsIR/L - SiCsI
408  KVINDRAReconNuc* SIX;
409  TIter nextsix(&sixparts);
410  while ((SIX = (KVINDRAReconNuc*)nextsix())) {
411  if (SIX->GetStatus() == KVReconstructedNucleus::kStatusOK) {
412  SIX->Identify();
413  if (SIX->IsIdentified()) {
414  nident++;
415  if (SIX->GetCodes().TestIDCode(kIDCode5)) SIX->SetIDCode(kIDCode7);
416  else SIX->SetIDCode(kIDCode6);
417  SIX->Calibrate();
418  if (SIX->IsCalibrated()) SIX->SubtractEnergyFromAllDetectors();
419  } else {
420  // failure of ChIo-Si identification: particle stopped in ChIo ?
421  // estimation of Z (minimum) from energy loss (if detector is calibrated)
422  UInt_t zmin = ((KVDetector*)SIX->GetDetectorList()->Last())->FindZmin(-1., SIX->GetMassFormula());
423  if (zmin) {
424  SIX->SetZ(zmin);
425  SIX->SetIsIdentified();
426  SIX->SetIDCode(kIDCode7);
427  // "Identifying" telescope is taken from list of ID telescopes
428  // to which stopping detector belongs
429  SIX->SetIdentifyingTelescope((KVIDTelescope*)SIX->GetStoppingDetector()->GetIDTelescopes()->Last());
430  SIX->Calibrate();
431  }
432  }
433  }
434  }
435  }
436  if (nident) { // newly-identified particles may change status of others in group
437  // reanalyse group
438  KVReconstructedNucleus::AnalyseParticlesInGroup(grp);
439  nident = 0;
440  }
441 
442  TIter parts2(grp->GetParticles()); // list may have changed if we have added particles
443  // identify & calibrate any remaining particles with status=KVReconstructedNucleus::kStatusOK
444  while ((nuc = (KVINDRAReconNuc*)parts2())) {
445  if (!nuc->IsIdentified() && nuc->GetStatus() == KVReconstructedNucleus::kStatusOK) {
446  nuc->ResetNSegDet();
447  nuc->Identify();
448  if (nuc->IsIdentified()) {
449  nident++;
450  nuc->Calibrate();
451  if (nuc->IsCalibrated()) nuc->SubtractEnergyFromAllDetectors();
452  }
453  }
454  }
455  if (nident) { // newly-identified particles may change status of others in group
456  // reanalyse group
457  KVReconstructedNucleus::AnalyseParticlesInGroup(grp);
458  nident = 0;
459  }
460 
461  // any kStatusOKafterShare particles ?
462  TList shareChIo;
463  parts2.Reset();
464  while ((nuc = (KVINDRAReconNuc*)parts2())) {
465  if (!nuc->IsIdentified() && nuc->GetStatus() == KVReconstructedNucleus::kStatusOKafterShare) {
466  shareChIo.Add(nuc);
467  }
468  }
469  Int_t nshares = shareChIo.GetEntries();
470  if (nshares) {
471  KVChIo* theChIo = ((KVINDRAReconNuc*)shareChIo.At(0))->GetChIo();
472  if (theChIo && nshares > 1) {
473  // divide chio energy equally
474  Double_t Eshare = theChIo->GetEnergyLoss() / nshares;
475  theChIo->SetEnergyLoss(Eshare);
476  // modify PG and GG of ChIo according to new energy loss
477  Double_t volts = theChIo->GetVoltsFromEnergy(Eshare);
478  Double_t GG = theChIo->GetCanalGGFromVolts(volts);
479  Double_t PG = theChIo->GetCanalPGFromVolts(volts);
480  theChIo->GetACQParam("PG")->SetData(TMath::Min(4095, (Int_t)PG));
481  theChIo->GetACQParam("GG")->SetData(TMath::Min(4095, (Int_t)GG));
482  }
483  // now try to identify
484  TIter nextSh(&shareChIo);
485  while ((nuc = (KVINDRAReconNuc*)nextSh())) {
486  nuc->SetNSegDet(10);
487  nuc->Identify();
488  if (nuc->IsIdentified()) {
489  nuc->SetIDCode(kIDCode8);
490  nuc->Calibrate();
491  }
492  }
493  }
494 
495  // any remaining stopped in first stage particles ?
496  parts2.Reset();
497  while ((nuc = (KVINDRAReconNuc*)parts2())) {
498  if (!nuc->IsIdentified() && nuc->GetStatus() == KVReconstructedNucleus::kStatusStopFirstStage) {
499  // estimation of Z (minimum) from energy loss (if detector is calibrated)
500  UInt_t zmin = ((KVDetector*)nuc->GetDetectorList()->Last())->FindZmin(-1., nuc->GetMassFormula());
501  if (zmin) {
502  nuc->SetZ(zmin);
503  nuc->SetIsIdentified();
504  nuc->SetIDCode(kIDCode5);
505  // "Identifying" telescope is taken from list of ID telescopes
506  // to which stopping detector belongs
507  nuc->SetIdentifyingTelescope((KVIDTelescope*)nuc->GetStoppingDetector()->GetIDTelescopes()->Last());
508  nuc->Calibrate();
509  }
510  }
511  }
512 }
513 */
514 
int Int_t
ClassImp(KVPartitionList) void KVPartitionList
Initialisation.
#define d(i)
const Bool_t kFALSE
bool Bool_t
double Double_t
const Bool_t kTRUE
char * Form(const char *fmt,...)
virtual Bool_t IsType(const Char_t *typ) const
Definition: KVBase.h:184
Base class for detector geometry description.
Definition: KVDetector.h:159
virtual Double_t GetMaxDeltaE(Int_t Z, Int_t A)
virtual Bool_t Fired(Option_t *opt="any") const
Definition: KVDetector.h:450
Double_t GetDetectorSignalValue(const KVString &type, const KVNameValueList &params="") const
Definition: KVDetector.h:492
virtual void SetEResAfterDetector(Double_t e)
Definition: KVDetector.h:629
virtual Double_t GetDeltaEFromERes(Int_t Z, Int_t A, Double_t Eres)
void SetAnalysed(Bool_t b=kTRUE)
Definition: KVDetector.h:446
virtual Double_t GetCorrectedEnergy(KVNucleus *, Double_t e=-1., Bool_t transmission=kTRUE)
Definition: KVDetector.cpp:810
KVNameValueList * GetParameters() const
Definition: KVEvent.h:203
Path taken by particles through multidetector geometry.
const KVSeqCollection * GetIDTelescopes() const
Information on relative positions of detectors & particle trajectories.
const Char_t * GetName() const
Name of node is same as name of associated detector.
KVDetector * GetDetector() const
virtual KVReconstructedNucleus * ReconstructTrajectory(const KVGeoDNTrajectory *traj, const KVGeoDetectorNode *node)
KVIDTelescope * identifying_telescope
telescope which identified current particle
std::unordered_map< std::string, KVIdentificationResult * > id_by_type
identification results by type for current particle
KVReconstructedEvent * GetEventFragment() const
int nfireddets
number of fired detectors in group for current event
virtual void IdentifyParticle(KVReconstructedNucleus &PART)
Double_t GetTargetEnergyLossCorrection(KVReconstructedNucleus *ion)
void SetCalibrationStatus(KVReconstructedNucleus &PART, UShort_t code)
TString GetPartSeedCond() const
KVIdentificationResult partID
identification to be applied to current particle
virtual Bool_t Identify(KVIdentificationResult *, Double_t x=-1., Double_t y=-1.)
Base class for all detectors or associations of detectors in array which can identify charged particl...
Definition: KVIDTelescope.h:83
virtual Bool_t IsReadyForID()
Reconstruct particles in INDRA groups.
void SetNoCalibrationStatus(KVReconstructedNucleus *n)
void IdentifyParticle(KVReconstructedNucleus &PART)
Bool_t CalculateChIoDEFromResidualEnergy(KVReconstructedNucleus *n, Double_t ERES)
double DoBeryllium8Calibration(KVReconstructedNucleus *n)
virtual bool DoCoherencyAnalysis(KVReconstructedNucleus &)=0
void CalibrateParticle(KVReconstructedNucleus *PART)
void CheckCsIEnergy(KVReconstructedNucleus *n)
KVDetector * GetCsI(KVReconstructedNucleus *n)
virtual void DoCalibration(KVReconstructedNucleus *)
void SetBadCalibrationStatus(KVReconstructedNucleus *n)
KVReconstructedNucleus * ReconstructTrajectory(const KVGeoDNTrajectory *traj, const KVGeoDetectorNode *node)
KVDetector * theChio
the ChIo of the group
Full result of one attempted particle identification.
Bool_t IDOK
general quality of identification, =kTRUE if acceptable identification made
Int_t IDquality
specific quality code returned by identification procedure
Int_t IDcode
a general identification code for this type of identification
void IncrementValue(const Char_t *name, value_type value)
Description of properties and kinematics of atomic nuclei.
Definition: KVNucleus.h:125
Int_t GetZ() const
Return the number of proton / atomic number.
Definition: KVNucleus.cpp:770
Double_t GetEnergy() const
Definition: KVParticle.h:623
void SetParameter(const Char_t *name, ValType value) const
Definition: KVParticle.h:820
void SetEnergy(Double_t e)
Definition: KVParticle.h:601
Nuclei reconstructed from data measured by a detector array .
const KVReconNucTrajectory * GetReconstructionTrajectory() const
void SetIdentification(KVIdentificationResult *, KVIDTelescope *)
virtual Int_t GetIDCode() const
virtual void SetTargetEnergyLoss(Double_t e)
virtual void GetAnglesFromReconstructionTrajectory(Option_t *opt="random")
virtual TObject * FindObjectByType(const Char_t *) const
Iterator end() const
Particle * AddParticle()
const char * Data() const
const Int_t n
Double_t Abs(Double_t d)