BmnRoot
Loading...
Searching...
No Matches
MpdMCStack.cxx
Go to the documentation of this file.
1// -------------------------------------------------------------------------
2// ----- MpdMCStack source file -----
3// ----- Created 09/10/08 by M. Al-Turany -----
4// -------------------------------------------------------------------------
5
6#include "MpdMCStack.h"
7#ifdef BMNROOT
8#include "CbmMCTrack.h"
9#else
10#include "FairMCTrack.h"
11#endif
12#include "FairRootManager.h"
13#include "FairLogger.h"
14
15#include "TGeoTrack.h"
16#include "TClonesArray.h"
17#include "TEveManager.h"
18#include "TParticlePDG.h"
19#include "TDatabasePDG.h"
20
21#include <iostream>
22using namespace std;
23
24
25// ----- Default constructor -------------------------------------------
28
29// ----- Standard constructor ------------------------------------------
30MpdMCStack::MpdMCStack(const char* name, Int_t iVerbose)
31 : FairTask(name, iVerbose),
32 fEveTrList(new TObjArray(16))
33{}
34
35// ----- Destructor ----------------------------------------------------
38
39// -------------------------------------------------------------------------
40InitStatus MpdMCStack::Init()
41{
42 if (fVerbose > 1) cout<<"MpdMCStack::Init()"<<endl;
43
44 FairRootManager* fManager = FairRootManager::Instance();
45
46 fTrackList = (TClonesArray*) fManager->GetObject("MCTrack");
47 if (fTrackList == 0)
48 {
49 LOG(error)<<"MpdMCStack::Init() branch "<<GetName()<<" not found! Task will be deactivated";
50 SetActive(kFALSE);
51 }
52 if (fVerbose > 2)
53 cout<<"MpdMCStack::Init() get track list"<<fTrackList<<endl;
54
56 if (fVerbose > 2) cout<<"MpdMCStack::Init() get instance of MpdEventManager"<<endl;
57
58 fPro = new FairGeanePro();
59 gMC3 = (TGeant3*) gMC;
60
61 x1[0] = 0; x1[1] = 0; x1[2] = 0;
62 p1[0] = 0; p1[1] = 0; p1[2] = 0;
63 x2[0] = 0; x2[1] = 0; x2[2] = 0;
64 p2[0] = 0; p2[1] = 0; p2[2] = 0;
65 for (Int_t i = 0; i < 15; i++)
66 ein[i] = 0;
67
68 //Float_t length[1];
69 //length[0] = 100.0;
70 //gMC3->Eufill(1, ein, length);
71
72 fTrajFilter = FairTrajFilter::Instance();
73
74 if (IsActive())
75 return kSUCCESS;
76 else
77 return kERROR;
78}
79
80// -------------------------------------------------------------------------
81void MpdMCStack::Exec(Option_t* /*option*/)
82{
83 if (!IsActive())
84 return;
85
86 if (fVerbose > 1) cout<<"MpdMCStack::Exec"<<endl;
87 Reset();
88
89 for (Int_t i = 0; i < fTrackList->GetEntriesFast(); i++)
90 {
91 if (fVerbose > 2) cout<<"MpdMCStack::Exec "<<i<<endl;
92
93#ifdef BMNROOT
94 CbmMCTrack* tr = (CbmMCTrack*) fTrackList->At(i);
95#else
96 FairMCTrack* tr = (FairMCTrack*) fTrackList->At(i);
97#endif
98
99 TVector3 Ptot;
100 tr->GetMomentum(Ptot);
101 Int_t MotherId =tr->GetMotherId();
102 TVector3 Vertex;
103 tr->GetStartVertex(Vertex);
104
105 Double_t time = tr->GetStartT()*1e-09;
106 x1[0] = Vertex.x(); x1[1] = Vertex.y(); x1[2] = Vertex.z();
107 p1[0] = Ptot.Px(); p1[1] = Ptot.Py(); p1[2] = Ptot.Pz();
108
109 //TParticle* P = (TParticle*) tr->GetParticle();
110 TParticlePDG* fParticlePDG = TDatabasePDG::Instance()->GetParticle(tr->GetPdgCode());
111
112 Double_t mass = 0, ene = 0;
113 if (fParticlePDG)
114 mass = fParticlePDG->Mass();
115 if (mass >= 0)
116 ene = TMath::Sqrt(mass*mass + Ptot.Mag2());
117
118 //TParticle(Int_t pdg, Int_t status, Int_t mother1, Int_t mother2, Int_t daughter1, Int_t daughter2, Double_t px, Double_t py, Double_t pz, Double_t etot, Double_t vx, Double_t vy, Double_t vz, Double_t time)
119 TParticle* P = new TParticle(tr->GetPdgCode(), i, MotherId, -1, -1, -1, Ptot.Px(), Ptot.Py(),Ptot.Pz(),ene, Vertex.x(), Vertex.z(), Vertex.z(), time);
120 //Int_t particle_mother_id = P->GetMother(0);
121 Float_t particle_energy = P->Energy();
122 fEventManager->ExpandEnergyLimits(particle_energy);
123
124 if (fVerbose > 2)
125 cout<<"MinEnergyLimit "<<fEventManager->GetMinEnergyLimit()<<" MaxEnergyLimit "<<fEventManager->GetMaxEnergyLimit()<<endl;
126
127 if ((fEventManager->IsPriOnly() && (P->GetMother(0) > -1)))
128 continue;
129 if ((fEventManager->fCurrentPDG.size() != 0) && (fEventManager->fCurrentPDG.find(tr->GetPdgCode()) == fEventManager->fCurrentPDG.end()))
130 continue;
131 if (fVerbose > 2)
132 cout<<"PEnergy "<<particle_energy<<" Min "<<fEventManager->GetMinEnergyCut()<<" Max "<<fEventManager->GetMaxEnergyCut()<<endl;
133 if ((particle_energy < fEventManager->GetMinEnergyCut()) || (particle_energy > fEventManager->GetMaxEnergyCut()))
134 continue;
135
136 //Int_t Np = tr->GetNpoints();
137 fTrList = GetTrGroup(P);
138 TEveTrack* track= new TEveTrack(P, tr->GetPdgCode(), fTrPr);
139 track->SetLineColor(fEventManager->Color(tr->GetPdgCode()));
140
141 // PROPAGATION
142 // Int_t GeantCode = TDatabasePDG::Instance()->ConvertPdgToGeant3(tr->GetPdgCode());
143 // gMC3->Ertrak(x1, p1, x2, p2,G eantCode, "L");
144 fPro->SetDestinationLength(100.0);
145 fPro->Propagate(x1, p1, x2, p2, tr->GetPdgCode());
146 TGeoTrack* tr1 = fTrajFilter->GetCurrentTrk();
147
148 const Double_t* point;
149 Int_t Np = tr1->GetNpoints();
150 for (Int_t n = 0; n < Np; n++)
151 {
152 point = tr1->GetPoint(n);
153 track->SetPoint(n, point[0], point[1], point[2]);
154 TEveVector pos = TEveVector(point[0], point[1], point[2]);
155
156 TEvePathMark* path = new TEvePathMark();
157 path->fV = pos ;
158 path->fTime = point[3];
159 if (n == 0)
160 {
161 TEveVector Mom = TEveVector(P->Px(), P->Py(),P->Pz());
162 path->fP = Mom;
163 }
164
165 track->AddPathMark(*path);
166 if (fVerbose > 3) cout<<"Path marker added "<<path<<endl;
167 }
168
169 fTrList->AddElement(track);
170 if (fVerbose > 3) cout<<"Track added "<<track->GetName()<<endl;
171 }
172
173 //for (Int_t i = 0; i < fEveTrList->GetEntriesFast(); i++)
174 //{
175 // TEveTrackList* TrListIn = (TEveTrackList*) fEveTrList->At(i);
176 // TrListIn->FindMomentumLimits(TrListIn, kFALSE);
177 //}
178
179 //gEve->Redraw3D(kFALSE);
180}
181
182// -------------------------------------------------------------------------
185
186// -------------------------------------------------------------------------
189
190// -------------------------------------------------------------------------
192{
193 for (Int_t i = 0; i < fEveTrList->GetEntriesFast(); i++)
194 {
195 TEveTrackList* ele = (TEveTrackList*) fEveTrList->At(i);
196 gEve->RemoveElement(ele, fEventManager);
197 }
198 fEveTrList->Clear();
199}
200
201TEveTrackList* MpdMCStack::GetTrGroup(TParticle* P)
202{
203 fTrList = 0;
204 for (Int_t i = 0; i < fEveTrList->GetEntriesFast(); i++)
205 {
206 TEveTrackList* TrListIn = (TEveTrackList*) fEveTrList->At(i);
207
208 if (strcmp(TrListIn->GetName(), P->GetName()) == 0)
209 {
210 fTrList= TrListIn;
211 break;
212 }
213 }
214
215 if (fTrList == 0)
216 {
217 fTrPr = new TEveTrackPropagator();
218 fTrList = new TEveTrackList(P->GetName(), fTrPr);
219 fTrList->SetMainColor(fEventManager->Color(P->GetPdgCode()));
220 fTrList->SetRnrLine(kTRUE);
221
222 fEveTrList->Add(fTrList);
223 gEve->AddElement(fTrList, fEventManager);
224 }
225
226 return fTrList;
227}
int i
Definition P4_F32vec4.h:22
void GetMomentum(TVector3 &momentum) const
Definition CbmMCTrack.h:139
Int_t GetMotherId() const
Definition CbmMCTrack.h:57
Double_t GetStartT() const
Definition CbmMCTrack.h:64
void GetStartVertex(TVector3 &vertex) const
Definition CbmMCTrack.h:149
Int_t GetPdgCode() const
Definition CbmMCTrack.h:56
virtual void ExpandEnergyLimits(Float_t check_energy, float add_multiplier=1)
virtual Bool_t IsPriOnly()
static MpdEventManager * Instance()
virtual Float_t GetMinEnergyLimit()
unordered_set< Int_t > fCurrentPDG
virtual Int_t Color(Int_t pdg)
virtual Float_t GetMinEnergyCut()
virtual Float_t GetMaxEnergyCut()
virtual Float_t GetMaxEnergyLimit()
Float_t ein[15]
Definition MpdMCStack.h:53
TEveTrackList * GetTrGroup(TParticle *P)
TEveTrackList * fTrList
Definition MpdMCStack.h:51
virtual ~MpdMCStack()
TGeant3 * gMC3
Definition MpdMCStack.h:54
FairGeanePro * fPro
Definition MpdMCStack.h:55
TObjArray * fEveTrList
Definition MpdMCStack.h:49
virtual InitStatus Init()
virtual void Exec(Option_t *option)
FairTrajFilter * fTrajFilter
Definition MpdMCStack.h:56
TClonesArray * fTrackList
Definition MpdMCStack.h:48
virtual void Finish()
Float_t x2[3]
Definition MpdMCStack.h:53
TEveTrackPropagator * fTrPr
Definition MpdMCStack.h:50
Float_t p1[3]
Definition MpdMCStack.h:53
virtual void SetParContainers()
Float_t p2[3]
Definition MpdMCStack.h:53
MpdEventManager * fEventManager
Definition MpdMCStack.h:46
Float_t x1[3]
Definition MpdMCStack.h:53
STL namespace.