9 #include "ERDecayEXP1803.h" 14 #include "TVirtualMC.h" 15 #include "TLorentzVector.h" 16 #include "TMCProcess.h" 19 #include "FairRunSim.h" 20 #include "FairLogger.h" 22 #include "ERDecayMCEventHeader.h" 23 #include "ERMCEventHeader.h" 25 ERDecayEXP1803::ERDecayEXP1803():
38 fIs5HUserMassSet(false),
39 fIs5HExcitationSet(false),
43 fRnd =
new TRandom3();
45 fRnd2 =
new TRandom3();
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);
55 fLv5H =
new TLorentzVector();
56 fLv3He =
new TLorentzVector();
58 cout <<
"ERDecayEXP1803 constructed." << endl;
62 ERDecayEXP1803::~ERDecayEXP1803() {
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;
74 f5HExcitationWeight.push_back(f5HExcitationWeight.back() + distibWeight);
78 Bool_t ERDecayEXP1803::Init() {
80 cout <<
"Decayer Init." << endl;
82 f6He = TDatabasePDG::Instance()->GetParticle(
"6He");
84 std::cerr <<
"-W- ERDecayEXP1803: Ion 6He not found in database!" << endl;
88 f2H = TDatabasePDG::Instance()->GetParticle(
"Deuteron");
90 std::cerr <<
"-W- ERDecayEXP1803: Ion Deuteron not found in database!" << endl;
94 f5H = TDatabasePDG::Instance()->GetParticle(
"5H");
96 std::cerr <<
"-W- ERDecayEXP1803: Ion 5H not found in database!" << endl;
100 f3He = TDatabasePDG::Instance()->GetParticle(fIon3He->GetName());
102 std::cerr <<
"-W- ERDecayEXP1803: Ion 3He not found in database!" << endl;
106 f3H = TDatabasePDG::Instance()->GetParticle(
"Triton");
108 std::cerr <<
"-W- ERDecayEXP1803: Ion Triton not found in database!" << endl;
112 fn = TDatabasePDG::Instance()->GetParticle(
"neutron");
114 std::cerr <<
"-W- ERDecayEXP1803: Particle neutron not found in database!" << endl;
117 if (fIs5HUserMassSet) {
118 fUnstableIon5H->SetMass(f5HMass / .931494028);
120 f5HMass = f5H->Mass();
127 Bool_t ERDecayEXP1803::Stepping() {
128 if(!fDecayFinish && gMC->TrackPid() == 1000020060
129 && TString(gMC->CurrentVolName()).Contains(GetInteractionVolumeName()))
131 if (!fIsInterationPointFound) {
132 if (!FindInteractionPoint()) {
133 fDecayFinish = kTRUE;
136 fDistanceFromEntrance = 0;
139 gMC->SetMaxStep(fMinStep);
140 TLorentzVector curPos;
141 gMC->TrackPosition(curPos);
142 Double_t trackStep = gMC->TrackStep();
143 fDistanceFromEntrance += trackStep;
146 if (fDistanceFromEntrance > fDistanceToInteractPoint) {
149 TLorentzVector lv6He;
150 gMC->TrackMomentum(lv6He);
152 if (lv6He.P() == 0) {
156 TLorentzVector lv2H(0., 0., 0., f2H->Mass());
157 TLorentzVector lvReaction;
158 lvReaction = lv6He + lv2H;
160 const TVector3 boost = lvReaction.BoostVector();
162 TLorentzVector lv6HeCM, lv2HCM;
165 lv6HeCM.Boost(-boost);
166 lv2HCM.Boost(-boost);
167 ECM = lv6HeCM(3) + lv2HCM(3);
169 Int_t decayHappen = kFALSE;
171 Double_t decay5HMass;
172 while (decayHappen==kFALSE) {
173 decay5HMass = f5HMass;
174 Double_t excitation = 0;
175 if (fIs5HExcitationSet) {
176 Double_t randWeight = gRandom->Uniform(0., f5HExcitationWeight.back());
177 Int_t distribNum = 0;
179 for (; distribNum < f5HExcitationWeight.size(); distribNum++) {
180 if (randWeight < f5HExcitationWeight[distribNum]) {
184 excitation = gRandom->Gaus(f5HExcitationMean[distribNum], f5HExcitationSigma[distribNum]);
185 fUnstableIon5H->SetExcEnergy(excitation);
187 decay5HMass += excitation;
188 if((ECM - f3He->Mass() - decay5HMass) > 0) {
193 PhaseGenerator(ECM, decay5HMass);
195 fLv3He->Boost(boost);
198 Double_t decayMasses[3];
199 decayMasses[0] = f3H->Mass();
200 decayMasses[1] = fn->Mass();
201 decayMasses[2] = fn->Mass();
203 fDecayPhaseSpace->SetDecay(*fLv5H, 3, decayMasses);
204 fDecayPhaseSpace->Generate();
206 TLorentzVector *lv3H = fDecayPhaseSpace->GetDecay(0);
207 TLorentzVector *lvn1 = fDecayPhaseSpace->GetDecay(1);
208 TLorentzVector *lvn2 = fDecayPhaseSpace->GetDecay(2);
211 Int_t He6TrackNb, H5TrackNb, He3TrackNb, H3TrackNb, n1TrackNb, n2TrackNb;
213 He6TrackNb = gMC->GetStack()->GetCurrentTrackNumber();
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);
243 fDecayFinish = kTRUE;
244 gMC->SetMaxStep(100.);
246 FairRunSim* run = FairRunSim::Instance();
247 if (TString(run->GetMCEventHeader()->ClassName()).Contains(
"ERDecayMCEventHeader")){
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);
263 void ERDecayEXP1803::BeginEvent() {
264 fDecayFinish = kFALSE;
265 fIsInterationPointFound = kFALSE;
266 fTargetReactZ = fRnd->Uniform(-fTargetThickness / 2, fTargetThickness / 2);
267 FairRunSim* run = FairRunSim::Instance();
271 void ERDecayEXP1803::FinishEvent() {
272 FairRunSim* run = FairRunSim::Instance();
273 if (TString(run->GetMCEventHeader()->ClassName()).Contains(
"ERDecayMCEventHeader")){
281 Double_t m1 = h5Mass;
282 Double_t m2 = f3He->Mass();
286 Double_t E1 = 0.5 * (Ecm * Ecm + m1 * m1 - m2 * m2) / Ecm;
289 Double_t Pcm = TMath::Sqrt(E1 * E1 - m1 * m1);
293 if(fADInput == NULL) {
294 thetaCM = TMath::ACos(gRandom->Uniform(-1, 1));
297 thetaCM = fADFunction->GetRandom(1., fADInput->GetN()-1)*TMath::DegToRad();
299 Double_t phi = gRandom->Uniform(0., 2. * TMath::Pi());
302 Pcmv.SetMagThetaPhi(Pcm, thetaCM, phi);
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);
312 if (fADInput->IsZombie()) {
313 Error(
"ERDecayEXP1803::ADEvaluate",
"AD input was not loaded");
317 return fADInput->Eval(x[0]);
322 TString ADFilePath = gSystem->Getenv(
"VMCWORKDIR");
323 ADFilePath +=
"/input/generators/" + ADFile;
325 fADInput =
new TGraph(ADFilePath,
"%lg %*lg %lg");
327 if (fADInput->GetN() <= 0) {
328 LOG(INFO) <<
"ERDecayEXP1803::SetAngularDistribution: " 329 <<
"Too few inputs for creation of AD function!" << FairLogger::endl;
332 Double_t* angle = fADInput->GetX();
338 angle[0], angle[fADInput->GetN()-1], 0,
"ERDecayEXP1803",
"ADEvaluate");
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