er  dev
ERDecay10Heto8He.cxx
1 /********************************************************************************
2  * Copyright (C) Joint Institute for Nuclear Research *
3  * *
4  * This software is distributed under the terms of the *
5  * GNU Lesser General Public Licence version 3 (LGPL) version 3, *
6  * copied verbatim in the file "LICENSE" *
7  ********************************************************************************/
8 
9 #include "ERDecay10Heto8He.h"
10 
11 #include <iostream>
12 #include <string>
13 #include <sstream>
14 using namespace std;
15 
16 #include "TVirtualMC.h"
17 #include "TLorentzVector.h"
18 #include "TMCProcess.h"
19 #include "TRandom.h"
20 #include "TVector.h"
21 
22 #include "FairRunSim.h"
23 #include "FairLogger.h"
24 
25 #include "ERDecayMCEventHeader.h"
26 #include "ER10Heto8HeEventHeader.h"
27 #include "ERMCEventHeader.h"
28 
29 ERDecay10Heto8He::ERDecay10Heto8He():
30  ERDecay("10Heto8He"),
31  fDecayFinish(kFALSE),
32  fTargetReactZ(0.),
33  fMinStep(0.01),
34  f8He(NULL),
35  f10He(NULL),
36  f3H(NULL),
37  fn(NULL),
38  f1H(NULL),
39  f10HeMass(0.),
40  fIs10HeUserMassSet(false),
41  fIs10HeExcitationSet(false),
42  fADInput(NULL),
43  fADFunction(NULL),
44  fDecayFilePath(""),
45  fDecayFileFinished(kFALSE),
46  fDecayFileCurrentEvent(0)
47 {
48  fRnd = new TRandom3();
49  // fRnd->SetSeed();
50  fRnd2 = new TRandom3();
51  fRnd2->SetSeed();
52  fReactionPhaseSpace = new TGenPhaseSpace();
53  fDecayPhaseSpace = new TGenPhaseSpace();
54  FairRunSim* run = FairRunSim::Instance();
55  fUnstableIon10He = new FairIon("10He", 2, 10, 2);
56  // fIon1H = new FairIon("1H", 1, 1, 1);
57  run->AddNewIon(fUnstableIon10He);
58  // run->AddNewIon(fIon1H);
59 
60  fLv1H = new TLorentzVector();
61  fLv10He = new TLorentzVector();
62 }
63 
64 // // //-------------------------------------------------------------------------------------------------
65 ERDecay10Heto8He::~ERDecay10Heto8He() {
66  if (fDecayFile.is_open())
67  fDecayFile.close();
68  if (fDecayFilePath == ""){ // LV from TGenPhaseSpace will be deleted in TGenPhaseSpace
69  if (fLv8Hed){
70  delete fLv8Hed;
71  delete fLvn1;
72  delete fLvn2;
73  // delete fLvn3;
74  // delete fLvn4;
75  }
76  }
77 }
78 
79 //-------------------------------------------------------------------------------------------------
80 void ERDecay10Heto8He::SetHe10Exitation(Double_t excMean, Double_t fwhm, Double_t distibWeight) {
81  f10HeExcitationMean.push_back(excMean);
82  f10HeExcitationSigma.push_back(fwhm / 2.355);
83  if (!fIs10HeExcitationSet) {
84  f10HeExcitationWeight.push_back(distibWeight);
85  fIs10HeExcitationSet = true;
86  return ;
87  }
88  f10HeExcitationWeight.push_back(f10HeExcitationWeight.back() + distibWeight);
89 }
90 
91 //-------------------------------------------------------------------------------------------------
92 Bool_t ERDecay10Heto8He::Init() {
93 
94  cout << "Decayer Init." << endl;
95 
96  f8He = TDatabasePDG::Instance()->GetParticle("8He");
97  if ( ! f8He ) {
98  std::cerr << "-W- ERDecay10Heto8He: Ion 8He not found in database!" << endl;
99  return kFALSE;
100  }
101 
102  f3H = TDatabasePDG::Instance()->GetParticle("Triton");
103  if ( ! f3H ) {
104  std::cerr << "-W- ERDecay10Heto8He: Ion Triton not found in database!" << endl;
105  return kFALSE;
106  }
107 
108 
109  fn = TDatabasePDG::Instance()->GetParticle("neutron");
110  if ( ! fn ) {
111  std::cerr << "-W- ERDecay10Heto8He: Particle neutron not found in database!" << endl;
112  return kFALSE;
113  }
114 
115  f1H = TDatabasePDG::Instance()->GetParticle("proton");
116  if ( ! f1H ) {
117  std::cerr << "-W- ERDecay10Heto8He: Particle proton not found in database!" << endl;
118  return kFALSE;
119  }
120 
121  if (fIs10HeUserMassSet) {
122  fUnstableIon10He->SetMass(f10HeMass / .931494028);
123  } else {
124  f10HeMass = f10He->Mass(); // if user mass is not defined in ERDecay10Heto8He::SetHe10Mass() than get a GEANT mass
125  }
126  CalculateTargetParameters();
127 
128 
129  printf("\n++++++++++++\nList of particles\n++++++++++++\n");
130  f8He->Print();
131  printf("mass 8He %f\n",f8He->Mass());
132  printf("mass 10He %f\n",f10HeMass);
133  printf("mass p %f\n",f1H->Mass());
134  printf("mass n %f\n",fn->Mass());
135  f3H->Print();
136  f1H->Print();
137  fn->Print();
138 
139 
140  if (fDecayFilePath != ""){
141  LOG(INFO) << "Use decay kinematics from external text file" << FairLogger::endl;
142  fDecayFile.open(fDecayFilePath.Data());
143  if (!fDecayFile.is_open())
144  LOG(FATAL) << "Can`t open decay file " << fDecayFilePath << FairLogger::endl;
145  //Пропускаем шапку файла
146  std::string header;
147  std::getline(fDecayFile,header);
148 
149  fLv8Hed = new TLorentzVector();
150  fLvn1 = new TLorentzVector();
151  fLvn2 = new TLorentzVector();
152  // fLvn3 = new TLorentzVector();
153  // fLvn4 = new TLorentzVector();
154  }
155 
156 
157  // f6He = TDatabasePDG::Instance()->GetParticle("6He");
158  // if ( ! f6He ) {
159  // std::cerr << "-W- ERDecay10Heto8He: Ion 6He not found in database!" << endl;
160  // return kFALSE;
161  // }
162 
163  // printf("8He PDG = %i\n",f8He->PdgCode());
164  // printf("6He PDG = %i\n",f6He->PdgCode());
165  // printf("triton PDG = %i\n",f3H->PdgCode());
166  // printf("neutron PDG = %i\n",fn->PdgCode());
167  // printf("proton PDG = %i\n",f1H->PdgCode());
168  // printf("10He PDG = %i\n",fUnstableIon10He->PdgCode());
169 
170  return kTRUE;
171 }
172 
173 //-------------------------------------------------------------------------------------------------
174 Bool_t ERDecay10Heto8He::Stepping() {
175  if(!fDecayFinish && gMC->TrackPid() == 1000020080
176  && TString(gMC->CurrentVolName()).Contains(GetInteractionVolumeName()))
177  {
178  if (!fIsInterationPointFound) {
179  if (!FindInteractionPoint()) {
180  fDecayFinish = kTRUE;
181  return kTRUE;
182  } else {
183  fDistanceFromEntrance = 0;
184  }
185  }
186  gMC->SetMaxStep(fMinStep);
187  TLorentzVector curPos;
188  gMC->TrackPosition(curPos);
189  Double_t trackStep = gMC->TrackStep();
190  fDistanceFromEntrance += trackStep;
191  // std::cout << "Track step: " << fDistanceFromEntrance << "; Diff " << (curPos.Vect() - fInputPoint).Mag() << endl;
192  // std::cout << "Track step: " << fDistanceFromEntrance << endl;
193  if (fDistanceFromEntrance > fDistanceToInteractPoint) {
194  // std::cout << "Start reation in target. Defined pos: " << fDistanceToInteractPoint << ", current pos: " << curPos.Z() << endl;
195 
196  // 8He + 2H → 3He + 7H
197  // 8He + 3H -> 1H + 10He
198  TLorentzVector lv8Heb;
199  gMC->TrackMomentum(lv8Heb);
200 
201  if (lv8Heb.P() == 0) { // temporary fix of bug with zero kinetic energy
202  return kTRUE;
203  }
204 
205  TLorentzVector lv3H(0., 0., 0., f3H->Mass());
206  TLorentzVector lvReaction;
207  lvReaction = lv8Heb + lv3H;
208 
209  const TVector3 boost = lvReaction.BoostVector(); //Get Pcm 3 vector
210  Double_t ECM = 0;
211  TLorentzVector lv8HebCM, lv3HCM;
212  lv8HebCM = lv8Heb;
213  lv3HCM = lv3H;
214  lv8HebCM.Boost(-boost);
215  lv3HCM.Boost(-boost);
216  ECM = lv8HebCM(3) + lv3HCM(3);
217 
218  Int_t reactionHappen = kFALSE;
219 
220  Double_t decay10HeMass;
221  Int_t reactionAttempsCounter = 0;
222  Double_t excitation = 0; // excitation energy
223  while (reactionHappen==kFALSE) { // while reaction condition is not fullfilled
224  decay10HeMass = f10HeMass;
225  if (fIs10HeExcitationSet) {
226  Double_t randWeight = gRandom->Uniform(0., f10HeExcitationWeight.back());
227  Int_t distribNum = 0;
228  // choose distribution by weight
229  for (; distribNum < f10HeExcitationWeight.size(); distribNum++) {
230  if (randWeight < f10HeExcitationWeight[distribNum]) {
231  break;
232  }
233  }
234  excitation = gRandom->Gaus(f10HeExcitationMean[distribNum], f10HeExcitationSigma[distribNum]);
235  fUnstableIon10He->SetExcEnergy(excitation);
236  }
237  decay10HeMass += excitation;
238  if((ECM - f1H->Mass() - decay10HeMass) > 0) { // выход из цикла while для PhaseGenerator
239  reactionHappen = kTRUE;
240  LOG(DEBUG) << "[ERDecayEXP1811] Reaction is happen" << endl;
241  }
242  reactionAttempsCounter++;
243  if (reactionAttempsCounter > 1000){
244  LOG(DEBUG) << "[ERDecayEXP1811] Reaction is forbidden for this CM energy" << endl;
245  fDecayFinish = kTRUE;
246  return kTRUE;
247  }
248  }
249 
250  ReactionPhaseGenerator(ECM, decay10HeMass);
251  fLv10He->Boost(boost);
252  fLv1H->Boost(boost);
253 
254  //7H → f3H + n +n +n +n
255  //10He -> 8He + n + n
256  if (!DecayPhaseGenerator(excitation)){
257  fDecayFinish = kTRUE;
258  return kTRUE;
259  }
260 
261  Int_t He8bTrackNb, He10TrackNb,H1TrackNb, He8dTrackNb, n1TrackNb, n2TrackNb;
262 
263  He8bTrackNb = gMC->GetStack()->GetCurrentTrackNumber();
264  // std::cout << "He8TrackNb " << He8TrackNb << std::endl;
265  /* //???????????
266  gMC->GetStack()->PushTrack(1, He8TrackNb, f7H->PdgCode(),
267  fLv7H->Px(), fLv7H->Py(), fLv7H->Pz(),
268  fLv7H->E(), curPos.X(), curPos.Y(), curPos.Z(),
269  gMC->TrackTime(), 0., 0., 0.,
270  kPDecay, H7TrackNb, decay7HMass, 0);*/
271 
272  gMC->GetStack()->PushTrack(1, He8bTrackNb, f1H->PdgCode(),
273  fLv1H->Px(), fLv1H->Py(), fLv1H->Pz(),
274  fLv1H->E(), curPos.X(), curPos.Y(), curPos.Z(),
275  gMC->TrackTime(), 0., 0., 0.,
276  kPDecay, H1TrackNb, f1H->Mass(), 0);
277  gMC->GetStack()->PushTrack(1, He8bTrackNb, f8He->PdgCode(),
278  fLv8Hed->Px(), fLv8Hed->Py(), fLv8Hed->Pz(),
279  fLv8Hed->E(), curPos.X(), curPos.Y(), curPos.Z(),
280  gMC->TrackTime(), 0., 0., 0.,
281  kPDecay, He8dTrackNb, f8He->Mass(), 0);
282  gMC->GetStack()->PushTrack(1, He8bTrackNb, fn->PdgCode(),
283  fLvn1->Px(),fLvn1->Py(),fLvn1->Pz(),
284  fLvn1->E(), curPos.X(), curPos.Y(), curPos.Z(),
285  gMC->TrackTime(), 0., 0., 0.,
286  kPDecay, n1TrackNb, fn->Mass(), 0);
287  gMC->GetStack()->PushTrack(1, He8bTrackNb, fn->PdgCode(),
288  fLvn2->Px(),fLvn2->Py(),fLvn2->Pz(),
289  fLvn2->E(), curPos.X(), curPos.Y(), curPos.Z(),
290  gMC->TrackTime(), 0., 0., 0.,
291  kPDecay, n2TrackNb, fn->Mass(), 0);
292  gMC->StopTrack();
293  fDecayFinish = kTRUE;
294  gMC->SetMaxStep(100.);
295 
296  FairRunSim* run = FairRunSim::Instance();
297  if (TString(run->GetMCEventHeader()->ClassName()).Contains("ERDecayMCEventHeader")){
298  ERDecayMCEventHeader* header = (ERDecayMCEventHeader*)run->GetMCEventHeader();
299  header->SetDecayPos(curPos.Vect());
300  header->SetInputIon(He8bTrackNb);
301  header->AddOutputParticle(He10TrackNb);
302  header->AddOutputParticle(H1TrackNb);
303  header->AddOutputParticle(He8dTrackNb);
304  header->AddOutputParticle(n1TrackNb);
305  header->AddOutputParticle(n2TrackNb);
306  }
307  if (TString(run->GetMCEventHeader()->ClassName()).Contains("ER10Heto8HeEventHeader")){
308  ER10Heto8HeEventHeader* header = (ER10Heto8HeEventHeader*)run->GetMCEventHeader();
309  header->SetData(curPos.Vect(), lv8Heb, lv3H, *fLv1H, *fLv8Hed, *fLv10He, *fLvn1, *fLvn2);
310  header->SetTrigger(1);
311  }
312  }
313  }
314  return kTRUE;
315 }
316 
317 //-------------------------------------------------------------------------------------------------
318 void ERDecay10Heto8He::BeginEvent() {
319  fDecayFinish = kFALSE;
320  fIsInterationPointFound = kFALSE;
321  fTargetReactZ = fRnd->Uniform(-fTargetThickness / 2, fTargetThickness / 2);
322  FairRunSim* run = FairRunSim::Instance();
323 }
324 
325 //-------------------------------------------------------------------------------------------------
326 void ERDecay10Heto8He::FinishEvent() {
327  FairRunSim* run = FairRunSim::Instance();
328  if (TString(run->GetMCEventHeader()->ClassName()).Contains("ERDecayMCEventHeader")){
329  ERDecayMCEventHeader* header = (ERDecayMCEventHeader*)run->GetMCEventHeader();
330  header->Clear();
331  }
332  if (TString(run->GetMCEventHeader()->ClassName()).Contains("ER10Heto8HeEventHeader")){
333  ER10Heto8HeEventHeader* header = (ER10Heto8HeEventHeader*)run->GetMCEventHeader();
334  header->Clear();
335  }
336 }
337 
338 //-------------------------------------------------------------------------------------------------
339 void ERDecay10Heto8He::ReactionPhaseGenerator(Double_t Ecm, Double_t he10Mass) {
340  Double_t m1 = he10Mass;
341  Double_t m2 = f1H->Mass();
342 
343  // Energy of 1-st particle in cm.
344  // total energy of the first particle is calculated as
345  Double_t E1 = 0.5 * (Ecm * Ecm + m1 * m1 - m2 * m2) / Ecm;
346 
347  //Impulse in CM
348  Double_t Pcm = TMath::Sqrt(E1 * E1 - m1 * m1);
349  //Generate angles of particles in CM
350  Double_t thetaCM;
351  if(!fADInput) { // if file with angular distribution isn't setted than isotropic distribution is generated
352 // thetaCM = TMath::ACos(gRandom->Uniform(-1, 1));
353 // thetaCM = TMath::ACos(gRandom->Uniform(-1, -0.965925826));
354 
355  thetaCM = TMath::ACos(gRandom->Uniform(0.984807753,1)); //вырезаю задние углы для попадания в протонный детектор
356  } else {
357  thetaCM = fADFunction->GetRandom(fThetaMin, fThetaMax)*TMath::DegToRad();
358  }
359  Double_t phi = gRandom->Uniform(0., 2. * TMath::Pi());
360  TVector3 Pcmv;
361  Pcmv.SetMagThetaPhi(Pcm, thetaCM, phi);
362 
363  fLv10He->SetXYZM(0., 0., 0., 0.);
364  fLv1H->SetXYZM(0., 0., 0., 0.);
365  fLv10He->SetXYZM(Pcmv(0), Pcmv(1), Pcmv(2), m1);
366  fLv1H->SetXYZM(-Pcmv(0), -Pcmv(1), -Pcmv(2), m2);
367 }
368 
369 //-------------------------------------------------------------------------------------------------
370 Bool_t ERDecay10Heto8He::DecayPhaseGenerator(const Double_t excitation) {
371  if (fDecayFilePath == ""){ // if decay file not defined, per morm decay using phase space
372  // Double_t decayMasses[5];
373  Double_t decayMasses[3];
374  decayMasses[0] = f8He->Mass();
375  decayMasses[1] = fn->Mass();
376  decayMasses[2] = fn->Mass();
377 
378  // printf("8He = %.10f\t n = %.10f\n",decayMasses[0],decayMasses[1]);
379  // printf("+++ fLv10He %f\t%f\n",fLv10He->E()-fLv10He->Mag(),fLv10He->Mag()-decayMasses[0]-decayMasses[1]-decayMasses[2]);
380  fDecayPhaseSpace->SetDecay(*fLv10He, 3, decayMasses);
381  fDecayPhaseSpace->Generate();
382  fLv8Hed = fDecayPhaseSpace->GetDecay(0);
383  fLvn1 = fDecayPhaseSpace->GetDecay(1);
384  fLvn2 = fDecayPhaseSpace->GetDecay(2);
385  return kTRUE;
386  }
387  if (fDecayFile.eof()){
388  LOG(ERROR) << "Decay file finished! There are no more events in file " << fDecayFilePath
389  << " to be processed." << FairLogger::endl;
390  return kFALSE;
391  }
392  std::string event_line;
393  std::getline(fDecayFile,event_line);
394  std::istringstream iss(event_line);
395  std::vector<std::string> outputs_components((std::istream_iterator<std::string>(iss)),
396  std::istream_iterator<std::string>());
397  // if (outputs_components.size() < 5*3){
398  if (outputs_components.size() < 3*3){
399  LOG(ERROR) << "Wrong components number in raw in decay file!" << FairLogger::endl;
400  return kFALSE;
401  }
402  // Fill momentum vectors in CM.
403  TVector3 p8Hed(std::stod(outputs_components[0]),std::stod(outputs_components[1]),
404  std::stod(outputs_components[2]));
405  TVector3 pn1(std::stod(outputs_components[3]),std::stod(outputs_components[4]),
406  std::stod(outputs_components[5]));
407  TVector3 pn2(std::stod(outputs_components[6]),std::stod(outputs_components[7]),
408  std::stod(outputs_components[8]));
409 
410  // TVector3 pn1(std::stod(outputs_components[0]),std::stod(outputs_components[1]),
411  // std::stod(outputs_components[2]));
412  // TVector3 pn2(std::stod(outputs_components[3]),std::stod(outputs_components[4]),
413  // std::stod(outputs_components[5]));
414  // TVector3 p8Hed(std::stod(outputs_components[12]),std::stod(outputs_components[13]),
415  // std::stod(outputs_components[14]));
416 
417  // Apply scale factor
418  const auto excitationScale = excitation > 0. ? sqrt(excitation / fDecayFileExcitation) : 1.;
419  // printf("\n\n++++++ %f ++++++\n",excitationScale);
420  const auto MeV2GeV = 1./1000.;
421  const auto scale = excitationScale * MeV2GeV;
422  pn1 *= scale;
423  pn2 *= scale;
424  p8Hed *= scale;
425  const auto fill_output_lorentz_vectors_in_lab =
426  [this](TLorentzVector* lv, const TVector3& p, const Double_t mass) {
427  lv->SetXYZM(p.X(), p.Y(), p.Z(), mass);
428  lv->Boost(fLv10He->BoostVector());
429  };
430  fill_output_lorentz_vectors_in_lab(fLvn1, pn1, fn->Mass());
431  fill_output_lorentz_vectors_in_lab(fLvn2, pn2, fn->Mass());
432  fill_output_lorentz_vectors_in_lab(fLv8Hed, p8Hed, f8He->Mass());
433  return kTRUE;
434 }
435 
436 //-------------------------------------------------------------------------------------------------
437 Double_t ERDecay10Heto8He::ADEvaluate(Double_t *x, Double_t *p) {
438  if (fADInput->IsZombie()) {
439  Error("ERDecay10Heto8He::ADEvaluate", "AD input was not loaded");
440  return -1;
441  }
442  // on each step of creating distribution function returns interpolated value of input data
443  return fADInput->Eval(x[0]);
444 }
445 
446 //-------------------------------------------------------------------------------------------------
448  TString ADFilePath = gSystem->Getenv("VMCWORKDIR");
449  ADFilePath += "/input/" + ADFile;
450  std::ifstream f;
451  f.open(ADFilePath.Data());
452  if (!f.is_open()) {
453  LOG(FATAL) << "Can't open file " << ADFilePath << FairLogger::endl;
454  }
455  Int_t nPoints = std::count(std::istreambuf_iterator<char>(f),
456  std::istreambuf_iterator<char>(), '\n');
457  f.seekg(0, std::ios::beg);
458  TVectorD tet(nPoints);
459  TVectorD sigma(nPoints);
460  LOG(DEBUG2) << "nPoints = " << nPoints << FairLogger::endl;
461  Int_t i = 0;
462  while (!f.eof()) {
463  // Костыль
464  if (i == nPoints) break;
465  f >> tet(i) >> sigma(i);
466  LOG(DEBUG2) << i << ": " << tet(i) << "\t" << sigma(i) << FairLogger::endl;
467  i++;
468  }
469  fADInput = new TGraph(tet, sigma);
470  if (fADInput->GetN() <= 0) { //if there are no points in input file
471  LOG(INFO) << "ERDecay10Heto8He::SetAngularDistribution: "
472  << "Too few inputs for creation of AD function!" << FairLogger::endl;
473  return;
474  }
475  Double_t* angle = fADInput->GetX(); // get first column variables that contains number of point
476 
477  // Creation of angular distribution function using class member function.
478  // Constructor divides interval (0; fADInput->GetN()-1) into grid.
479  // On each step of grid it calls ADEvaluate() to get interpolated values of input data.
480  fThetaMin = angle[0];
481  fThetaMax = angle[fADInput->GetN()-1];
482  fADFunction = new TF1("angDistr", this, &ERDecay10Heto8He::ADEvaluate,
483  fThetaMin, fThetaMax, 0, "ERDecay10Heto8He", "ADEvaluate");
484  fADFunction->Eval(1.);
485 }
486 //-------------------------------------------------------------------------------------------------
487 // void ERDecay10Heto8He::PrintMass(){
488  // printf("\n++++++++++++\nList of particles\n++++++++++++\n");
489  // f8He->Print();
490  // f10He->Print();
491  // f3H->Print();
492  // f1H->Print();
493  // fn->Print();
494  // // return;
495 // }
496 //-------------------------------------------------------------------------------------------------
497 ClassImp(ERDecay10Heto8He)
TString GetInteractionVolumeName()
Method returns interaction volume name.
Definition: ERDecay.h:45
void ReactionPhaseGenerator(Double_t Ecm, Double_t he10Mass)
Body reaction in phase space approach.
void SetAngularDistribution(TString ADfile)
Sets distribution is contained in file.
Bool_t FindInteractionPoint()
Definition: ERDecay.cxx:122