FairRoot
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
FairEveGeoTracks.cxx
Go to the documentation of this file.
1 /********************************************************************************
2  * Copyright (C) 2020 GSI Helmholtzzentrum fuer Schwerionenforschung GmbH *
3  * *
4  * This software is distributed under the terms of the *
5  * GNU Lesser General Public Licence (LGPL) version 3, *
6  * copied verbatim in the file "LICENSE" *
7  ********************************************************************************/
8 /*
9  * FairGeoTracks.cxx
10  *
11  * Created on: 23 kwi 2020
12  * Author: Daniel Wielanek
13  * E-mail: daniel.wielanek@gmail.com
14  * Warsaw University of Technology, Faculty of Physics
15  */
16 
17 #include "FairEveGeoTracks.h"
18 
19 #include "FairEveTrack.h" // for FairEveTrack
20 #include "FairEventManager.h" // for FairEventManager
21 #include "FairGetEventTime.h"
22 #include "FairRootManager.h" // for FairRootManager
23 
24 #include <TBranch.h>
25 #include <TClonesArray.h> // for TClonesArray
26 #include <TEveManager.h> // for TEveManager, gEve
27 #include <TEveTrack.h> // for TEveTrackList
28 #include <TGeoTrack.h> // for TGeoTrack
29 #include <TParticle.h> // for TParticle
30 #include <TString.h> // for Form
31 #include <TVector3.h> // for TVector3
32 #include <iostream>
33 #include <limits> // for numeric_limits
34 
36  : FairEveTracks()
37  , fContainer(nullptr)
38  , fShowPrimary(kTRUE)
39  , fShowSecondary(kTRUE)
40  , fUsePdg(kFALSE)
41  , fPdgCut(0)
42  , fTMin(0)
43  , fTMax(std::numeric_limits<Double_t>::max())
44 {
45  SetElementNameTitle("FairGeoTracks", "FairGeoTracks");
46 }
47 
49 {
51  fContainer = (TClonesArray *)mngr->GetObject("GeoTracks");
52  if (fContainer == nullptr)
53  return kFATAL;
54  fBranch = mngr->GetInTree()->GetBranch("GeoTracks");
56  return FairEveTracks::Init();
57 }
58 
60 {
61  TGeoTrack *tr = (TGeoTrack *)fContainer->UncheckedAt(id);
62  if (!CheckCuts(tr))
63  return;
64  TParticle *p = (TParticle *)tr->GetParticle();
65  Color_t color = GetEventManager()->Color(p->GetPdgCode());
66  TEveTrackList *trList = FindTrackGroup(p->GetName(), color);
67 
68  FairEveTrack *track = new FairEveTrack(p, p->GetPdgCode(), trList->GetPropagator());
69  track->SetElementTitle(Form("p={%4.3f,%4.3f,%4.3f}", p->Px(), p->Py(), p->Pz()));
70  track->SetMainColor(color);
71  Double_t x, y, z, t;
72  tr->GetPoint(0, x, y, z, t);
73  TVector3 pos(x, y, z);
74  TVector3 mom(p->Px(), p->Py(), p->Pz());
75  track->SetFirstPoint(mom, pos);
76  for (int i = 1; i < tr->GetNpoints(); i++) {
77  tr->GetPoint(i, x, y, z, t);
78  pos.SetXYZ(x, y, z);
79  track->SetNextPoint(pos);
80  }
81  track->SetRnrLine(kTRUE);
82  trList->AddElement(track);
83 }
84 
85 void FairEveGeoTracks::DrawAnimatedTrack(TGeoTrack *tr, double t0)
86 {
87  const Double_t timeScale = 1E+9;
88  if (!CheckCuts(tr))
89  return;
90  if (tr->GetNpoints() < 2)
91  return; // not enought points
92  if (tr->GetPoint(tr->GetNpoints() - 1)[3] * timeScale + t0 < fTMin)
93  return; // last point before tmin
94  if (tr->GetPoint(0)[3] * timeScale + t0 > fTMax)
95  return; // first point after tmax
96  TParticle *p = (TParticle *)tr->GetParticle();
97  Color_t color = GetEventManager()->Color(p->GetPdgCode());
98  TEveTrackList *trList = FindTrackGroup(p->GetName(), color);
99  FairEveTrack *track = new FairEveTrack(p, p->GetPdgCode(), trList->GetPropagator());
100  track->SetElementTitle(Form("p={%4.3f,%4.3f,%4.3f}", p->Px(), p->Py(), p->Pz()));
101  track->SetMainColor(color);
102  Double_t x, y, z, t; // currentPoint
103  Double_t xp, yp, zp, tp; // previousPoint
104  bool firstPoint = true;
105  bool previousPoint = false;
106  for (int i = 0; i < tr->GetNpoints(); i++) {
107  tr->GetPoint(i, x, y, z, t);
108  t = t * timeScale + t0;
109  TVector3 pos(x, y, z);
110  if (t < fTMin) {
111  xp = x;
112  yp = y;
113  zp = z;
114  tp = t;
115  previousPoint = true;
116  } else if (t > fTMin && t < fTMax) {
117  if (firstPoint && previousPoint) {
118  Double_t dT = (fTMin - tp) / (t - tp);
119  Double_t dx = x - xp;
120  Double_t dy = y - yp;
121  Double_t dz = z - zp;
122  pos.SetXYZ(xp + dx * dT, yp + dy * dT, zp + dz * dT);
123  }
124  if (firstPoint) {
125  TVector3 mom(p->Px(), p->Py(), p->Pz());
126  track->SetFirstPoint(mom, pos);
127  firstPoint = false;
128  } else {
129  track->SetNextPoint(pos);
130  }
131  xp = x;
132  yp = y;
133  zp = z;
134  tp = t;
135  }
136  if (t > fTMax) { // outside of time limits
137  Double_t dT = (fTMax - tp) / (t - tp);
138  Double_t dx = x - xp;
139  Double_t dy = y - yp;
140  Double_t dz = z - zp;
141  pos.SetXYZ(xp + dx * dT, yp + dy * dT, zp + dz * dT);
142  track->SetNextPoint(pos);
143  break;
144  }
145  }
146  track->SetRnrLine(kTRUE);
147  trList->AddElement(track);
148 }
149 
151 {
152  bool useGeoTrackHandler = false;
153  if (FairRunAna::Instance()->IsTimeStamp() && fBranch != nullptr) {
154  Double_t simTime = FairEventManager::Instance()->GetEvtTime();
155  std::pair<int, double> evt = FairGetEventTime::Instance().GetEvent(simTime);
156 
157  if (evt.first < 0) {
158  fContainer->Clear();
159  } else {
160  fBranch->GetEvent(evt.first);
161  }
162  if (FairEventManager::Instance()->GetClearHandler() == kTRUE) {
163  fGeoTrackHandler.Reset();
164  }
165  if (evt.first > -1)
166  fGeoTrackHandler.FillTClonesArray(fContainer, evt.first, evt.second, simTime);
167  useGeoTrackHandler = true;
168  }
169  Int_t nTracks = 0;
170  if (useGeoTrackHandler) {
171  nTracks = fGeoTrackHandler.GetData().size();
172  } else {
173  nTracks = fContainer->GetEntriesFast();
174  }
175  RemoveElements();
176  FairEventManager::Instance()->GetTimeLimits(fTMin, fTMax);
177  if (fTMin > fTMax) { // wrong time limits draw entire tracks
178  for (int iTrack = 0; iTrack < nTracks; iTrack++) {
179  DrawTrack(iTrack);
180  }
181  } else {
182  for (int iTrack = 0; iTrack < nTracks; iTrack++) {
183  TGeoTrack *track;
184  double t0 = .0;
185  if (useGeoTrackHandler) {
186  track = fGeoTrackHandler.GetData()[iTrack].first;
187  t0 = fGeoTrackHandler.GetData()[iTrack].second;
188  } else {
189  track = (TGeoTrack *)fContainer->At(iTrack);
190  }
191  DrawAnimatedTrack(track, t0);
192  }
193  }
194  gEve->Redraw3D(kFALSE);
195 }
196 
197 Bool_t FairEveGeoTracks::CheckCuts(TGeoTrack *tr)
198 {
199  TParticle *p = (TParticle *)tr->GetParticle();
200  if (UsePtCut()) {
201  if (p->Pt() < GetPtMin())
202  return kFALSE;
203  if (p->Pt() > GetPtMax())
204  return kFALSE;
205  }
206  if (UseEtaCut()) {
207  if (p->Eta() < GetEtaMin())
208  return kFALSE;
209  if (p->Eta() > GetEtaMax())
210  return kFALSE;
211  }
212  if (UseEnergyCut()) {
213  if (p->Energy() < GetEnergyMin())
214  return kFALSE;
215  if (p->Energy() > GetEnergyMax())
216  return kFALSE;
217  }
218  if (fUsePdg) {
219  if (p->GetPdgCode() != fPdgCut)
220  return kFALSE;
221  }
222  if (fShowPrimary && p->IsPrimary())
223  return kTRUE;
224  if (fShowSecondary && !p->IsPrimary())
225  return kTRUE;
226  return kFALSE;
227 }
228 
std::vector< std::pair< T *, double > > & GetData()
virtual InitStatus Init()
void DrawTrack(Int_t id)
static FairGetEventTime & Instance()
InitStatus
Definition: FairTask.h:33
FairEventManager * GetEventManager() const
Definition: FairEveTracks.h:58
Bool_t CheckCuts(TGeoTrack *tr)
Float_t GetEvtTime()
current time in ns to display in the event display. Either set value or event time taken from FairRoo...
Double_t GetEtaMax() const
Definition: FairEveTracks.h:52
Bool_t UseEnergyCut() const
Definition: FairEveTracks.h:57
static FairRootManager * Instance()
void FillTClonesArray(TClonesArray *inputArray, int evtIndex, double t0Event, double t0Current)
Bool_t UsePtCut() const
Definition: FairEveTracks.h:55
Double_t GetEnergyMin() const
Definition: FairEveTracks.h:53
TObject * GetObject(const char *BrName)
std::pair< int, double > GetEvent(double simTime) const
Double_t GetEnergyMax() const
Definition: FairEveTracks.h:54
static FairRunAna * Instance()
Definition: FairRunAna.cxx:61
virtual Int_t Color(Int_t pdg)
void GetTimeLimits(Double_t &min, Double_t &max)
void DrawAnimatedTrack(TGeoTrack *tr, double t0=0)
Bool_t UseEtaCut() const
Definition: FairEveTracks.h:56
static FairEventManager * Instance()
Double_t GetPtMax() const
Definition: FairEveTracks.h:50
virtual InitStatus Init()
Double_t GetEtaMin() const
Definition: FairEveTracks.h:51
TEveTrackList * FindTrackGroup(TString groupName, Color_t color)
Double_t GetPtMin() const
Definition: FairEveTracks.h:49