er  dev
ERDecayEXP1811.h
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 #ifndef ERDecayEXP1811_H
10 #define ERDecayEXP1811_H
11 
12 #include <vector>
13 #include <fstream>
14 
15 #include "TGraph.h"
16 #include "TF1.h"
17 #include "TRandom3.h"
18 #include "TRandom2.h"
19 #include "TGenPhaseSpace.h"
20 #include "TLorentzVector.h"
21 
22 #include "FairIon.h"
23 
24 #include "ERDecay.h"
25 
26 class ERDecayEXP1811 : public ERDecay {
27 
28 public:
30  ~ERDecayEXP1811();
31 
32  /*Modifiers*/
33  void SetMinStep(Double_t minStep) {fMinStep = minStep;}
34  void SetTargetThickness(Double_t targetThickness) {fTargetThickness = targetThickness;}
35  void SetH7Mass(Double_t mass) {f7HMass = mass; fIs7HUserMassSet = true;}
36  void SetH7Exitation(Double_t excMean, Double_t fwhm, Double_t distibWeight);
37  void SetDecayFile(const TString& filePath, Double_t excitationEnergyInFile /*[GeV]*/){ fDecayFilePath = filePath; }
38 
42  void SetAngularDistribution(TString ADfile);
43 
44 public:
45  Bool_t Init();
46  Bool_t Stepping();
47 
48  void BeginEvent();
49  void FinishEvent();
50 
51 private:
56  void ReactionPhaseGenerator(Double_t Ecm, Double_t h7Mass);
57 
58  Bool_t DecayPhaseGenerator(Double_t excitation);
59 
60  std::vector<TLorentzVector> ReadDecayEvent();
61 
62 private:
63  TRandom3 *fRnd;
64  TRandom3 *fRnd2;
65 
66  TParticlePDG *f8He;
67  TParticlePDG *f2H;
68  TParticlePDG *f3He;
69  TParticlePDG *f7H;
70  TParticlePDG *f3H;
71  TParticlePDG *fn;
72 
73  TLorentzVector *fLv3He;
74  TLorentzVector *fLv7H;
75  TLorentzVector *fLv3H;
76  TLorentzVector *fLvn1;
77  TLorentzVector *fLvn2;
78  TLorentzVector *fLvn3;
79  TLorentzVector *fLvn4;
80 
81  FairIon *fIon3He;
82  FairIon *fUnstableIon7H;
83 
84  TGenPhaseSpace *fReactionPhaseSpace;
85  TGenPhaseSpace *fDecayPhaseSpace;
86  Double_t fTargetReactZ;
87  Double_t fMinStep;
88  Double_t fTargetThickness;
89  Bool_t fDecayFinish;
90 
91  std::vector<Double_t> f7HExcitationMean;
92  std::vector<Double_t> f7HExcitationSigma;
93  std::vector<Double_t> f7HExcitationWeight;
94 
95  Double_t f7HMass;
96  Bool_t fIs7HUserMassSet;
97  Bool_t fIs7HExcitationSet;
98 
99  TString fDecayFilePath;
100  Double_t fDecayFileExcitation = 1. /*[GeV]*/;
101  Bool_t fDecayFileFinished;
102  Int_t fDecayFileCurrentEvent;
103  std::ifstream fDecayFile;
104 
105  TGraph *fADInput = nullptr;
106  TF1 *fADFunction = nullptr;
107  Double_t fThetaMin = 0.;
108  Double_t fThetaMax = 0.;
109 
110  //ADEvaluate function is necessary for TF1 constructor
111  Double_t ADEvaluate(Double_t *x, Double_t *p);
112 
113  ClassDef(ERDecayEXP1811,1)
114 };
115 
116 #endif
void SetAngularDistribution(TString ADfile)
Sets distribution is contained in file.
Double_t fThetaMin
function describing AD (angular distribution) of binary reaction
TF1 * fADFunction
distribution (angular distribution) graph containing AD input
void ReactionPhaseGenerator(Double_t Ecm, Double_t h7Mass)
Body reaction in phase space approach.