9 #include "ERDecay10Heto8He.h" 16 #include "TVirtualMC.h" 17 #include "TLorentzVector.h" 18 #include "TMCProcess.h" 22 #include "FairRunSim.h" 23 #include "FairLogger.h" 25 #include "ERDecayMCEventHeader.h" 26 #include "ER10Heto8HeEventHeader.h" 27 #include "ERMCEventHeader.h" 29 ERDecay10Heto8He::ERDecay10Heto8He():
40 fIs10HeUserMassSet(false),
41 fIs10HeExcitationSet(false),
45 fDecayFileFinished(kFALSE),
46 fDecayFileCurrentEvent(0)
48 fRnd =
new TRandom3();
50 fRnd2 =
new TRandom3();
52 fReactionPhaseSpace =
new TGenPhaseSpace();
53 fDecayPhaseSpace =
new TGenPhaseSpace();
54 FairRunSim* run = FairRunSim::Instance();
55 fUnstableIon10He =
new FairIon(
"10He", 2, 10, 2);
57 run->AddNewIon(fUnstableIon10He);
60 fLv1H =
new TLorentzVector();
61 fLv10He =
new TLorentzVector();
65 ERDecay10Heto8He::~ERDecay10Heto8He() {
66 if (fDecayFile.is_open())
68 if (fDecayFilePath ==
""){
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;
88 f10HeExcitationWeight.push_back(f10HeExcitationWeight.back() + distibWeight);
92 Bool_t ERDecay10Heto8He::Init() {
94 cout <<
"Decayer Init." << endl;
96 f8He = TDatabasePDG::Instance()->GetParticle(
"8He");
98 std::cerr <<
"-W- ERDecay10Heto8He: Ion 8He not found in database!" << endl;
102 f3H = TDatabasePDG::Instance()->GetParticle(
"Triton");
104 std::cerr <<
"-W- ERDecay10Heto8He: Ion Triton not found in database!" << endl;
109 fn = TDatabasePDG::Instance()->GetParticle(
"neutron");
111 std::cerr <<
"-W- ERDecay10Heto8He: Particle neutron not found in database!" << endl;
115 f1H = TDatabasePDG::Instance()->GetParticle(
"proton");
117 std::cerr <<
"-W- ERDecay10Heto8He: Particle proton not found in database!" << endl;
121 if (fIs10HeUserMassSet) {
122 fUnstableIon10He->SetMass(f10HeMass / .931494028);
124 f10HeMass = f10He->Mass();
126 CalculateTargetParameters();
129 printf(
"\n++++++++++++\nList of particles\n++++++++++++\n");
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());
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;
147 std::getline(fDecayFile,header);
149 fLv8Hed =
new TLorentzVector();
150 fLvn1 =
new TLorentzVector();
151 fLvn2 =
new TLorentzVector();
174 Bool_t ERDecay10Heto8He::Stepping() {
175 if(!fDecayFinish && gMC->TrackPid() == 1000020080
178 if (!fIsInterationPointFound) {
180 fDecayFinish = kTRUE;
183 fDistanceFromEntrance = 0;
186 gMC->SetMaxStep(fMinStep);
187 TLorentzVector curPos;
188 gMC->TrackPosition(curPos);
189 Double_t trackStep = gMC->TrackStep();
190 fDistanceFromEntrance += trackStep;
193 if (fDistanceFromEntrance > fDistanceToInteractPoint) {
198 TLorentzVector lv8Heb;
199 gMC->TrackMomentum(lv8Heb);
201 if (lv8Heb.P() == 0) {
205 TLorentzVector lv3H(0., 0., 0., f3H->Mass());
206 TLorentzVector lvReaction;
207 lvReaction = lv8Heb + lv3H;
209 const TVector3 boost = lvReaction.BoostVector();
211 TLorentzVector lv8HebCM, lv3HCM;
214 lv8HebCM.Boost(-boost);
215 lv3HCM.Boost(-boost);
216 ECM = lv8HebCM(3) + lv3HCM(3);
218 Int_t reactionHappen = kFALSE;
220 Double_t decay10HeMass;
221 Int_t reactionAttempsCounter = 0;
222 Double_t excitation = 0;
223 while (reactionHappen==kFALSE) {
224 decay10HeMass = f10HeMass;
225 if (fIs10HeExcitationSet) {
226 Double_t randWeight = gRandom->Uniform(0., f10HeExcitationWeight.back());
227 Int_t distribNum = 0;
229 for (; distribNum < f10HeExcitationWeight.size(); distribNum++) {
230 if (randWeight < f10HeExcitationWeight[distribNum]) {
234 excitation = gRandom->Gaus(f10HeExcitationMean[distribNum], f10HeExcitationSigma[distribNum]);
235 fUnstableIon10He->SetExcEnergy(excitation);
237 decay10HeMass += excitation;
238 if((ECM - f1H->Mass() - decay10HeMass) > 0) {
239 reactionHappen = kTRUE;
240 LOG(DEBUG) <<
"[ERDecayEXP1811] Reaction is happen" << endl;
242 reactionAttempsCounter++;
243 if (reactionAttempsCounter > 1000){
244 LOG(DEBUG) <<
"[ERDecayEXP1811] Reaction is forbidden for this CM energy" << endl;
245 fDecayFinish = kTRUE;
250 ReactionPhaseGenerator(ECM, decay10HeMass);
251 fLv10He->Boost(boost);
256 if (!DecayPhaseGenerator(excitation)){
257 fDecayFinish = kTRUE;
261 Int_t He8bTrackNb, He10TrackNb,H1TrackNb, He8dTrackNb, n1TrackNb, n2TrackNb;
263 He8bTrackNb = gMC->GetStack()->GetCurrentTrackNumber();
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);
293 fDecayFinish = kTRUE;
294 gMC->SetMaxStep(100.);
296 FairRunSim* run = FairRunSim::Instance();
297 if (TString(run->GetMCEventHeader()->ClassName()).Contains(
"ERDecayMCEventHeader")){
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);
307 if (TString(run->GetMCEventHeader()->ClassName()).Contains(
"ER10Heto8HeEventHeader")){
309 header->SetData(curPos.Vect(), lv8Heb, lv3H, *fLv1H, *fLv8Hed, *fLv10He, *fLvn1, *fLvn2);
310 header->SetTrigger(1);
318 void ERDecay10Heto8He::BeginEvent() {
319 fDecayFinish = kFALSE;
320 fIsInterationPointFound = kFALSE;
321 fTargetReactZ = fRnd->Uniform(-fTargetThickness / 2, fTargetThickness / 2);
322 FairRunSim* run = FairRunSim::Instance();
326 void ERDecay10Heto8He::FinishEvent() {
327 FairRunSim* run = FairRunSim::Instance();
328 if (TString(run->GetMCEventHeader()->ClassName()).Contains(
"ERDecayMCEventHeader")){
332 if (TString(run->GetMCEventHeader()->ClassName()).Contains(
"ER10Heto8HeEventHeader")){
340 Double_t m1 = he10Mass;
341 Double_t m2 = f1H->Mass();
345 Double_t E1 = 0.5 * (Ecm * Ecm + m1 * m1 - m2 * m2) / Ecm;
348 Double_t Pcm = TMath::Sqrt(E1 * E1 - m1 * m1);
355 thetaCM = TMath::ACos(gRandom->Uniform(0.984807753,1));
357 thetaCM = fADFunction->GetRandom(fThetaMin, fThetaMax)*TMath::DegToRad();
359 Double_t phi = gRandom->Uniform(0., 2. * TMath::Pi());
361 Pcmv.SetMagThetaPhi(Pcm, thetaCM, phi);
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);
370 Bool_t ERDecay10Heto8He::DecayPhaseGenerator(
const Double_t excitation) {
371 if (fDecayFilePath ==
""){
373 Double_t decayMasses[3];
374 decayMasses[0] = f8He->Mass();
375 decayMasses[1] = fn->Mass();
376 decayMasses[2] = fn->Mass();
380 fDecayPhaseSpace->SetDecay(*fLv10He, 3, decayMasses);
381 fDecayPhaseSpace->Generate();
382 fLv8Hed = fDecayPhaseSpace->GetDecay(0);
383 fLvn1 = fDecayPhaseSpace->GetDecay(1);
384 fLvn2 = fDecayPhaseSpace->GetDecay(2);
387 if (fDecayFile.eof()){
388 LOG(ERROR) <<
"Decay file finished! There are no more events in file " << fDecayFilePath
389 <<
" to be processed." << FairLogger::endl;
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>());
398 if (outputs_components.size() < 3*3){
399 LOG(ERROR) <<
"Wrong components number in raw in decay file!" << FairLogger::endl;
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]));
418 const auto excitationScale = excitation > 0. ? sqrt(excitation / fDecayFileExcitation) : 1.;
420 const auto MeV2GeV = 1./1000.;
421 const auto scale = excitationScale * MeV2GeV;
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());
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());
437 Double_t ERDecay10Heto8He::ADEvaluate(Double_t *x, Double_t *p) {
438 if (fADInput->IsZombie()) {
439 Error(
"ERDecay10Heto8He::ADEvaluate",
"AD input was not loaded");
443 return fADInput->Eval(x[0]);
448 TString ADFilePath = gSystem->Getenv(
"VMCWORKDIR");
449 ADFilePath +=
"/input/" + ADFile;
451 f.open(ADFilePath.Data());
453 LOG(FATAL) <<
"Can't open file " << ADFilePath << FairLogger::endl;
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;
464 if (i == nPoints)
break;
465 f >> tet(i) >> sigma(i);
466 LOG(DEBUG2) << i <<
": " << tet(i) <<
"\t" << sigma(i) << FairLogger::endl;
469 fADInput =
new TGraph(tet, sigma);
470 if (fADInput->GetN() <= 0) {
471 LOG(INFO) <<
"ERDecay10Heto8He::SetAngularDistribution: " 472 <<
"Too few inputs for creation of AD function!" << FairLogger::endl;
475 Double_t* angle = fADInput->GetX();
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.);
TString GetInteractionVolumeName()
Method returns interaction volume name.
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()