er  dev
ERDecayEXP1803.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 "ERDecayEXP1803.h"
10 
11 #include <iostream>
12 using namespace std;
13 
14 #include "TVirtualMC.h"
15 #include "TLorentzVector.h"
16 #include "TMCProcess.h"
17 #include "TRandom.h"
18 
19 #include "FairRunSim.h"
20 #include "FairLogger.h"
21 
22 #include "ERDecayMCEventHeader.h"
23 #include "ERMCEventHeader.h"
24 
25 ERDecayEXP1803::ERDecayEXP1803():
26  ERDecay("EXP1803"),
27  fDecayFinish(kFALSE),
28  fTargetReactZ(0.),
29  fMinStep(0.01),
30  f6He(NULL),
31  f2H (NULL),
32  f3He(NULL),
33  f5H (NULL),
34  f3H (NULL),
35  fn (NULL),
36  fIon3He(NULL),
37  f5HMass(0.),
38  fIs5HUserMassSet(false),
39  fIs5HExcitationSet(false),
40  fADInput(NULL),
41  fADFunction(NULL)
42 {
43  fRnd = new TRandom3();
44  // fRnd->SetSeed();
45  fRnd2 = new TRandom3();
46  fRnd2->SetSeed();
47  fReactionPhaseSpace = new TGenPhaseSpace();
48  fDecayPhaseSpace = new TGenPhaseSpace();
49  FairRunSim* run = FairRunSim::Instance();
50  fUnstableIon5H = new FairIon("5H", 1, 5, 1);
51  fIon3He = new FairIon("3He", 2, 3, 2);
52  run->AddNewIon(fUnstableIon5H);
53  run->AddNewIon(fIon3He);
54 
55  fLv5H = new TLorentzVector();
56  fLv3He = new TLorentzVector();
57 
58  cout << "ERDecayEXP1803 constructed." << endl;
59 }
60 
61 //-------------------------------------------------------------------------------------------------
62 ERDecayEXP1803::~ERDecayEXP1803() {
63 }
64 
65 //-------------------------------------------------------------------------------------------------
66 void ERDecayEXP1803::SetH5Exitation(Double_t excMean, Double_t fwhm, Double_t distibWeight) {
67  f5HExcitationMean.push_back(excMean);
68  f5HExcitationSigma.push_back(fwhm / 2.355);
69  if (!fIs5HExcitationSet) {
70  f5HExcitationWeight.push_back(distibWeight);
71  fIs5HExcitationSet = true;
72  return ;
73  }
74  f5HExcitationWeight.push_back(f5HExcitationWeight.back() + distibWeight);
75 }
76 
77 //-------------------------------------------------------------------------------------------------
78 Bool_t ERDecayEXP1803::Init() {
79 
80  cout << "Decayer Init." << endl;
81 
82  f6He = TDatabasePDG::Instance()->GetParticle("6He");
83  if ( ! f6He ) {
84  std::cerr << "-W- ERDecayEXP1803: Ion 6He not found in database!" << endl;
85  return kFALSE;
86  }
87 
88  f2H = TDatabasePDG::Instance()->GetParticle("Deuteron");
89  if ( ! f2H ) {
90  std::cerr << "-W- ERDecayEXP1803: Ion Deuteron not found in database!" << endl;
91  return kFALSE;
92  }
93 
94  f5H = TDatabasePDG::Instance()->GetParticle("5H");
95  if ( ! f5H ) {
96  std::cerr << "-W- ERDecayEXP1803: Ion 5H not found in database!" << endl;
97  return kFALSE;
98  }
99 
100  f3He = TDatabasePDG::Instance()->GetParticle(fIon3He->GetName());
101  if ( ! f3He ) {
102  std::cerr << "-W- ERDecayEXP1803: Ion 3He not found in database!" << endl;
103  return kFALSE;
104  }
105 
106  f3H = TDatabasePDG::Instance()->GetParticle("Triton");
107  if ( ! f3H ) {
108  std::cerr << "-W- ERDecayEXP1803: Ion Triton not found in database!" << endl;
109  return kFALSE;
110  }
111 
112  fn = TDatabasePDG::Instance()->GetParticle("neutron");
113  if ( ! fn ) {
114  std::cerr << "-W- ERDecayEXP1803: Particle neutron not found in database!" << endl;
115  return kFALSE;
116  }
117  if (fIs5HUserMassSet) {
118  fUnstableIon5H->SetMass(f5HMass / .931494028);
119  } else {
120  f5HMass = f5H->Mass(); // if user mass is not defined in ERDecayEXP1803::SetH5Mass() than get a GEANT mass
121  }
122  // CalculateTargetParameters();
123  return kTRUE;
124 }
125 
126 //-------------------------------------------------------------------------------------------------
127 Bool_t ERDecayEXP1803::Stepping() {
128  if(!fDecayFinish && gMC->TrackPid() == 1000020060
129  && TString(gMC->CurrentVolName()).Contains(GetInteractionVolumeName()))
130  {
131  if (!fIsInterationPointFound) {
132  if (!FindInteractionPoint()) {
133  fDecayFinish = kTRUE;
134  return kFALSE;
135  } else {
136  fDistanceFromEntrance = 0;
137  }
138  }
139  gMC->SetMaxStep(fMinStep);
140  TLorentzVector curPos;
141  gMC->TrackPosition(curPos);
142  Double_t trackStep = gMC->TrackStep();
143  fDistanceFromEntrance += trackStep;
144  // std::cout << "Track step: " << fDistanceFromEntrance << "; Diff " << (curPos.Vect() - fInputPoint).Mag() << endl;
145  // std::cout << "Track step: " << fDistanceFromEntrance << endl;
146  if (fDistanceFromEntrance > fDistanceToInteractPoint) {
147  // std::cout << "Start reation in target. Defined pos: " << fDistanceToInteractPoint << ", current pos: " << curPos.Z() << endl;
148  // 6He + 2H → 3He + 5H
149  TLorentzVector lv6He;
150  gMC->TrackMomentum(lv6He);
151 
152  if (lv6He.P() == 0) { // temporary fix of bug with zero kinetic energy
153  return kTRUE;
154  }
155 
156  TLorentzVector lv2H(0., 0., 0., f2H->Mass());
157  TLorentzVector lvReaction;
158  lvReaction = lv6He + lv2H;
159 
160  const TVector3 boost = lvReaction.BoostVector(); //Get Pcm 3 vector
161  Double_t ECM = 0;
162  TLorentzVector lv6HeCM, lv2HCM;
163  lv6HeCM = lv6He;
164  lv2HCM = lv2H;
165  lv6HeCM.Boost(-boost);
166  lv2HCM.Boost(-boost);
167  ECM = lv6HeCM(3) + lv2HCM(3);
168 
169  Int_t decayHappen = kFALSE;
170 
171  Double_t decay5HMass;
172  while (decayHappen==kFALSE) { // while decay condition is not fullfilled
173  decay5HMass = f5HMass;
174  Double_t excitation = 0; // excitation energy
175  if (fIs5HExcitationSet) {
176  Double_t randWeight = gRandom->Uniform(0., f5HExcitationWeight.back());
177  Int_t distribNum = 0;
178  // choose distribution by weight
179  for (; distribNum < f5HExcitationWeight.size(); distribNum++) {
180  if (randWeight < f5HExcitationWeight[distribNum]) {
181  break;
182  }
183  }
184  excitation = gRandom->Gaus(f5HExcitationMean[distribNum], f5HExcitationSigma[distribNum]);
185  fUnstableIon5H->SetExcEnergy(excitation);
186  }
187  decay5HMass += excitation;
188  if((ECM - f3He->Mass() - decay5HMass) > 0) { // выход из цикла while для PhaseGenerator
189  decayHappen = kTRUE;
190  }
191  }
192 
193  PhaseGenerator(ECM, decay5HMass);
194  fLv5H->Boost(boost);
195  fLv3He->Boost(boost);
196 
197  //5H → f3H + n +n.
198  Double_t decayMasses[3];
199  decayMasses[0] = f3H->Mass();
200  decayMasses[1] = fn->Mass();
201  decayMasses[2] = fn->Mass();
202 
203  fDecayPhaseSpace->SetDecay(*fLv5H, 3, decayMasses);
204  fDecayPhaseSpace->Generate();
205 
206  TLorentzVector *lv3H = fDecayPhaseSpace->GetDecay(0);
207  TLorentzVector *lvn1 = fDecayPhaseSpace->GetDecay(1);
208  TLorentzVector *lvn2 = fDecayPhaseSpace->GetDecay(2);
209 
210 
211  Int_t He6TrackNb, H5TrackNb, He3TrackNb, H3TrackNb, n1TrackNb, n2TrackNb;
212 
213  He6TrackNb = gMC->GetStack()->GetCurrentTrackNumber();
214  // std::cout << "He6TrackNb " << He6TrackNb << std::endl;
215  /*
216  gMC->GetStack()->PushTrack(1, He6TrackNb, f5H->PdgCode(),
217  fLv5H->Px(), fLv5H->Py(), fLv5H->Pz(),
218  fLv5H->E(), curPos.X(), curPos.Y(), curPos.Z(),
219  gMC->TrackTime(), 0., 0., 0.,
220  kPDecay, H5TrackNb, decay5HMass, 0);*/
221 
222  gMC->GetStack()->PushTrack(1, He6TrackNb, f3He->PdgCode(),
223  fLv3He->Px(), fLv3He->Py(), fLv3He->Pz(),
224  fLv3He->E(), curPos.X(), curPos.Y(), curPos.Z(),
225  gMC->TrackTime(), 0., 0., 0.,
226  kPDecay, He3TrackNb, f3He->Mass(), 0);
227  gMC->GetStack()->PushTrack(1, He6TrackNb, f3H->PdgCode(),
228  lv3H->Px(), lv3H->Py(), lv3H->Pz(),
229  lv3H->E(), curPos.X(), curPos.Y(), curPos.Z(),
230  gMC->TrackTime(), 0., 0., 0.,
231  kPDecay, H3TrackNb, f3H->Mass(), 0);
232  gMC->GetStack()->PushTrack(1, He6TrackNb, fn->PdgCode(),
233  lvn1->Px(),lvn1->Py(),lvn1->Pz(),
234  lvn1->E(), curPos.X(), curPos.Y(), curPos.Z(),
235  gMC->TrackTime(), 0., 0., 0.,
236  kPDecay, n1TrackNb, fn->Mass(), 0);
237  gMC->GetStack()->PushTrack(1, He6TrackNb, fn->PdgCode(),
238  lvn2->Px(),lvn2->Py(),lvn2->Pz(),
239  lvn2->E(), curPos.X(), curPos.Y(), curPos.Z(),
240  gMC->TrackTime(), 0., 0., 0.,
241  kPDecay, n2TrackNb, fn->Mass(), 0);
242  gMC->StopTrack();
243  fDecayFinish = kTRUE;
244  gMC->SetMaxStep(100.);
245 
246  FairRunSim* run = FairRunSim::Instance();
247  if (TString(run->GetMCEventHeader()->ClassName()).Contains("ERDecayMCEventHeader")){
248  ERDecayMCEventHeader* header = (ERDecayMCEventHeader*)run->GetMCEventHeader();
249  header->SetDecayPos(curPos.Vect());
250  header->SetInputIon(He6TrackNb);
251  header->AddOutputParticle(H5TrackNb);
252  header->AddOutputParticle(He3TrackNb);
253  header->AddOutputParticle(H3TrackNb);
254  header->AddOutputParticle(n1TrackNb);
255  header->AddOutputParticle(n2TrackNb);
256  }
257  }
258  }
259  return kTRUE;
260 }
261 
262 //-------------------------------------------------------------------------------------------------
263 void ERDecayEXP1803::BeginEvent() {
264  fDecayFinish = kFALSE;
265  fIsInterationPointFound = kFALSE;
266  fTargetReactZ = fRnd->Uniform(-fTargetThickness / 2, fTargetThickness / 2);
267  FairRunSim* run = FairRunSim::Instance();
268 }
269 
270 //-------------------------------------------------------------------------------------------------
271 void ERDecayEXP1803::FinishEvent() {
272  FairRunSim* run = FairRunSim::Instance();
273  if (TString(run->GetMCEventHeader()->ClassName()).Contains("ERDecayMCEventHeader")){
274  ERDecayMCEventHeader* header = (ERDecayMCEventHeader*)run->GetMCEventHeader();
275  header->Clear();
276  }
277 }
278 
279 //-------------------------------------------------------------------------------------------------
280 void ERDecayEXP1803::PhaseGenerator(Double_t Ecm, Double_t h5Mass) {
281  Double_t m1 = h5Mass;
282  Double_t m2 = f3He->Mass();
283 
284  // Energy of 1-st particle in cm.
285  //total energy of the first particle is calculated as
286  Double_t E1 = 0.5 * (Ecm * Ecm + m1 * m1 - m2 * m2) / Ecm;
287 
288  //Impulse in CM
289  Double_t Pcm = TMath::Sqrt(E1 * E1 - m1 * m1);
290 
291  //Generate angles of particles in CM
292  Double_t thetaCM;
293  if(fADInput == NULL) { // if file with angular distribution isn't setted than isotropic distribution is generated
294  thetaCM = TMath::ACos(gRandom->Uniform(-1, 1));
295  }
296  else {
297  thetaCM = fADFunction->GetRandom(1., fADInput->GetN()-1)*TMath::DegToRad();
298  }
299  Double_t phi = gRandom->Uniform(0., 2. * TMath::Pi());
300 
301  TVector3 Pcmv;
302  Pcmv.SetMagThetaPhi(Pcm, thetaCM, phi);
303 
304  fLv5H->SetXYZM(0., 0., 0., 0.);
305  fLv3He->SetXYZM(0., 0., 0., 0.);
306  fLv5H->SetXYZM(Pcmv(0), Pcmv(1), Pcmv(2), m1);
307  fLv3He->SetXYZM(-Pcmv(0), -Pcmv(1), -Pcmv(2), m2);
308 }
309 
310 //-------------------------------------------------------------------------------------------------
311 Double_t ERDecayEXP1803::ADEvaluate(Double_t *x, Double_t *p) {
312  if (fADInput->IsZombie()) {
313  Error("ERDecayEXP1803::ADEvaluate", "AD input was not loaded");
314  return -1;
315  }
316  // on each step of creating distribution function returns interpolated value of input data
317  return fADInput->Eval(x[0]);
318 }
319 
320 //-------------------------------------------------------------------------------------------------
322  TString ADFilePath = gSystem->Getenv("VMCWORKDIR");
323  ADFilePath += "/input/generators/" + ADFile;
324 
325  fADInput = new TGraph(ADFilePath, "%lg %*lg %lg"); // TGraph object is used for reading data from file with distribution
326 
327  if (fADInput->GetN() <= 0) { //if there are no points in input file
328  LOG(INFO) << "ERDecayEXP1803::SetAngularDistribution: "
329  << "Too few inputs for creation of AD function!" << FairLogger::endl;
330  return;
331  }
332  Double_t* angle = fADInput->GetX(); // get first column variables that contains number of point
333 
334  // Creation of angular distribution function using class member function.
335  // Constructor divides interval (0; fADInput->GetN()-1) into grid.
336  // On each step of grid it calls ADEvaluate() to get interpolated values of input data.
337  fADFunction = new TF1("angDistr", this, &ERDecayEXP1803::ADEvaluate,
338  angle[0], angle[fADInput->GetN()-1], 0, "ERDecayEXP1803", "ADEvaluate");
339 }
340 //-------------------------------------------------------------------------------------------------
341 ClassImp(ERDecayEXP1803)
void PhaseGenerator(Double_t Ecm, Double_t h5Mass)
Body decay in phase space approach.
void SetAngularDistribution(TString ADfile)
Sets distribution is contained in file.
Double_t ADEvaluate(Double_t *x, Double_t *p)
function describing AD (angular distribution) of binary reaction