FairRoot
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
FairTestDetector.cxx
Go to the documentation of this file.
1 /********************************************************************************
2  * Copyright (C) 2014 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 #include "FairTestDetector.h"
9 
10 #include "FairDetectorList.h" // for DetectorId::kTutDet
11 #include "FairLink.h" // for FairLink
12 #include "FairRootManager.h" // for FairRootManager
13 #include "FairRun.h" // for FairRun
14 #include "FairRuntimeDb.h" // for FairRuntimeDb
15 #include "FairStack.h" // for FairStack
16 #include "FairTestDetectorGeo.h" // for FairTestDetectorGeo
17 #include "FairTestDetectorGeoPar.h" // for FairTestDetectorGeoPar
18 #include "FairTestDetectorPoint.h" // for FairTestDetectorPoint
19 #include "FairVolume.h" // for FairVolume
20 
21 #include <TClonesArray.h> // for TClonesArray
22 #include <TVirtualMC.h> // for TVirtualMC
23 #include <TVirtualMCStack.h> // for TVirtualMCStack
24 
26  : FairDetector("FairTestDetector", kTRUE, kTutDet)
27  , fTrackID(-1)
28  , fVolumeID(-1)
29  , fPos()
30  , fMom()
31  , fPosOut()
32  , fMomOut()
33  , fTime(-1.)
34  , fLength(-1.)
35  , fELoss(-1)
36  , fEventNr(0)
37  , fFairTestDetectorPointCollection(new TClonesArray("FairTestDetectorPoint"))
38 {}
39 
40 FairTestDetector::FairTestDetector(const char* name, Bool_t active)
41  : FairDetector(name, active, kTutDet)
42  , fTrackID(-1)
43  , fVolumeID(-1)
44  , fPos()
45  , fMom()
46  , fPosOut()
47  , fMomOut()
48  , fTime(-1.)
49  , fLength(-1.)
50  , fELoss(-1)
51  , fEventNr(0)
52  , fFairTestDetectorPointCollection(new TClonesArray("FairTestDetectorPoint"))
53 {}
54 
56 {
57  if (fFairTestDetectorPointCollection) {
58  fFairTestDetectorPointCollection->Delete();
59  delete fFairTestDetectorPointCollection;
60  }
61 }
62 
64 {
67  rtdb->getContainer("FairTestDetectorGeoPar");
68 }
69 
71 {
74  // Set parameters at entrance of volume. Reset ELoss.
75  if (TVirtualMC::GetMC()->IsTrackEntering()) {
76  fELoss = 0.;
77  fTime = TVirtualMC::GetMC()->TrackTime() * 1.0e09;
78  fLength = TVirtualMC::GetMC()->TrackLength();
79  TVirtualMC::GetMC()->TrackPosition(fPos);
80  TVirtualMC::GetMC()->TrackMomentum(fMom);
81  }
82 
83  // Sum energy loss for all steps in the active volume
84  fELoss += TVirtualMC::GetMC()->Edep();
85 
86  // Create FairTestDetectorPoint at exit of active volume
87  if (TVirtualMC::GetMC()->IsTrackExiting() || TVirtualMC::GetMC()->IsTrackStop()
88  || TVirtualMC::GetMC()->IsTrackDisappeared()) {
89  fTrackID = TVirtualMC::GetMC()->GetStack()->GetCurrentTrackNumber();
90  fVolumeID = vol->getMCid();
91  TVirtualMC::GetMC()->TrackPosition(fPosOut);
92  TVirtualMC::GetMC()->TrackMomentum(fMomOut);
93  if (fELoss == 0.) {
94  return kFALSE;
95  }
96  AddHit(fTrackID,
97  fVolumeID,
98  TVector3(fPos.X(), fPos.Y(), fPos.Z()),
99  TVector3(fMom.Px(), fMom.Py(), fMom.Pz()),
100  TVector3(fPosOut.X(), fPosOut.Y(), fPosOut.Z()),
101  TVector3(fMomOut.Px(), fMomOut.Py(), fMomOut.Pz()),
102  fTime,
103  fLength,
104  fELoss);
105 
106  // Increment number of FairTestDetector det points in TParticle
107  FairStack* stack = static_cast<FairStack*>(TVirtualMC::GetMC()->GetStack());
108  stack->AddPoint(kTutDet);
109  }
110 
111  return kTRUE;
112 }
113 
115 {
116 
117  fFairTestDetectorPointCollection->Clear();
118  fEventNr++;
119 }
120 
122 {
123 
131  "FairTestDetectorPoint", "FairTestDetector", fFairTestDetectorPointCollection, kTRUE);
132 }
133 
134 TClonesArray* FairTestDetector::GetCollection(Int_t iColl) const
135 {
136  if (iColl == 0) {
137  return fFairTestDetectorPointCollection;
138  } else {
139  return nullptr;
140  }
141 }
142 
143 void FairTestDetector::Reset() { fFairTestDetectorPointCollection->Clear(); }
144 
146 {
152  ConstructASCIIGeometry<FairTestDetectorGeo, FairTestDetectorGeoPar>(Geo, "FairTestDetectorGeoPar");
153 }
154 
156  Int_t detID,
157  TVector3 pos,
158  TVector3 mom,
159  TVector3 posOut,
160  TVector3 momOut,
161  Double_t time,
162  Double_t length,
163  Double_t eLoss)
164 {
165  TClonesArray& clref = *fFairTestDetectorPointCollection;
166  Int_t size = clref.GetEntriesFast();
167  FairTestDetectorPoint* myPoint =
168  new (clref[size]) FairTestDetectorPoint(trackID, detID, pos, mom, posOut, momOut, time, length, eLoss);
169  myPoint->SetLink(FairLink(-1, fEventNr, FairRootManager::Instance()->GetBranchId("MCTrack"), trackID));
170  return myPoint;
171 }
172 
virtual Bool_t ProcessHits(FairVolume *v=0)
list of container factories
Definition: FairRuntimeDb.h:24
virtual void Initialize()
virtual void Register()
void SetLink(FairLink link)
Sets the Links with a single FairLink.
virtual void EndOfEvent()
static FairRun * Instance()
Definition: FairRun.cxx:31
static FairRootManager * Instance()
ClassImp(FairEventBuilder)
FairParSet * getContainer(const Text_t *)
virtual void Initialize()
virtual void Reset()
virtual ~FairTestDetector()
Int_t getMCid()
Definition: FairVolume.h:57
FairRuntimeDb * GetRuntimeDb(void)
Definition: FairRun.h:80
void Register(const char *name, const char *Foldername, TNamed *obj, Bool_t toFile)
virtual TClonesArray * GetCollection(Int_t iColl) const
FairTestDetectorPoint * AddHit(Int_t trackID, Int_t detID, TVector3 pos, TVector3 mom, TVector3 posOut, TVector3 momOut, Double_t time, Double_t length, Double_t eLoss)
void AddPoint(DetectorId iDet)
Definition: FairStack.cxx:344