FairRoot
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
NewDetector.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 "NewDetector.h"
9 
10 #include "FairGeoInterface.h"
11 #include "FairGeoLoader.h"
12 #include "FairGeoNode.h"
13 #include "FairGeoVolume.h"
14 #include "FairRootManager.h"
15 #include "FairRun.h"
16 #include "FairRuntimeDb.h"
17 #include "FairVolume.h"
18 #include "MyProjDetectorList.h"
19 #include "MyProjStack.h"
20 #include "NewDetectorGeo.h"
21 #include "NewDetectorGeoPar.h"
22 #include "NewDetectorPoint.h"
23 
24 #include <TClonesArray.h>
25 #include <TGeoBBox.h>
26 #include <TGeoCompositeShape.h>
27 #include <TGeoManager.h>
28 #include <TGeoMaterial.h>
29 #include <TGeoMedium.h>
30 #include <TGeoTube.h>
31 #include <TVirtualMC.h>
32 #include <iostream>
33 using std::cout;
34 using std::endl;
35 
37  : FairDetector("NewDetector", kTRUE, kNewDetector)
38  , fTrackID(-1)
39  , fVolumeID(-1)
40  , fPos()
41  , fMom()
42  , fTime(-1.)
43  , fLength(-1.)
44  , fELoss(-1)
45  , fNewDetectorPointCollection(new TClonesArray("NewDetectorPoint"))
46 {}
47 
48 NewDetector::NewDetector(const char* name, Bool_t active)
49  : FairDetector(name, active, kNewDetector)
50  , fTrackID(-1)
51  , fVolumeID(-1)
52  , fPos()
53  , fMom()
54  , fTime(-1.)
55  , fLength(-1.)
56  , fELoss(-1)
57  , fNewDetectorPointCollection(new TClonesArray("NewDetectorPoint"))
58 {}
59 
61  : FairDetector(right)
62  , fTrackID(-1)
63  , fVolumeID(-1)
64  , fPos()
65  , fMom()
66  , fTime(-1.)
67  , fLength(-1.)
68  , fELoss(-1)
69  , fNewDetectorPointCollection(new TClonesArray("NewDetectorPoint"))
70 {}
71 
73 {
74  if (fNewDetectorPointCollection) {
75  fNewDetectorPointCollection->Delete();
76  delete fNewDetectorPointCollection;
77  }
78 }
79 
81 {
87  DefineSensitiveVolumes();
88 
91  NewDetectorGeoPar* par = (NewDetectorGeoPar*)(rtdb->getContainer("NewDetectorGeoPar"));
92 }
93 
95 {
98  // Set parameters at entrance of volume. Reset ELoss.
99  if (TVirtualMC::GetMC()->IsTrackEntering()) {
100  fELoss = 0.;
101  fTime = TVirtualMC::GetMC()->TrackTime() * 1.0e09;
102  fLength = TVirtualMC::GetMC()->TrackLength();
103  TVirtualMC::GetMC()->TrackPosition(fPos);
104  TVirtualMC::GetMC()->TrackMomentum(fMom);
105  }
106 
107  // Sum energy loss for all steps in the active volume
108  fELoss += TVirtualMC::GetMC()->Edep();
109 
110  // Create NewDetectorPoint at exit of active volume
111  if (TVirtualMC::GetMC()->IsTrackExiting() || TVirtualMC::GetMC()->IsTrackStop()
112  || TVirtualMC::GetMC()->IsTrackDisappeared()) {
113  fTrackID = TVirtualMC::GetMC()->GetStack()->GetCurrentTrackNumber();
114  fVolumeID = vol->getMCid();
115  if (fELoss == 0.) {
116  return kFALSE;
117  }
118  AddHit(fTrackID,
119  fVolumeID,
120  TVector3(fPos.X(), fPos.Y(), fPos.Z()),
121  TVector3(fMom.Px(), fMom.Py(), fMom.Pz()),
122  fTime,
123  fLength,
124  fELoss);
125 
126  // Increment number of NewDetector det points in TParticle
127  MyProjStack* stack = (MyProjStack*)TVirtualMC::GetMC()->GetStack();
128  stack->AddPoint(kNewDetector);
129  }
130 
131  return kTRUE;
132 }
133 
135 {
136 
137  LOG(info) << "NewDetector: " << fNewDetectorPointCollection->GetEntriesFast() << " points registered in this event";
138 
139  fNewDetectorPointCollection->Clear();
140 }
141 
143 {
144 
151  if (!gMC->IsMT()) {
152  FairRootManager::Instance()->Register("NewDetectorPoint", "NewDetector", fNewDetectorPointCollection, kTRUE);
153  } else {
154  FairRootManager::Instance()->RegisterAny("NewDetectorPoint", fNewDetectorPointCollection, kTRUE);
155  }
156 }
157 
158 TClonesArray* NewDetector::GetCollection(Int_t iColl) const
159 {
160  if (iColl == 0) {
161  return fNewDetectorPointCollection;
162  } else {
163  return NULL;
164  }
165 }
166 
167 void NewDetector::Reset() { fNewDetectorPointCollection->Clear(); }
168 
170 {
171  TGeoVolume* top = gGeoManager->GetTopVolume();
172  TGeoMedium* Si = gGeoManager->GetMedium("Si");
173  TGeoMedium* Carbon = gGeoManager->GetMedium("C");
174 
175  if (Si == 0) {
176  TGeoMaterial* matSi = new TGeoMaterial("Si", 28.0855, 14, 2.33);
177  Si = new TGeoMedium("Si", 2, matSi);
178  }
179  if (Carbon == 0) {
180  TGeoMaterial* matCarbon = new TGeoMaterial("C", 12.011, 6.0, 2.265);
181  Carbon = new TGeoMedium("C", 3, matCarbon);
182  }
183 
184  TGeoVolume* det1 = gGeoManager->MakeTubs("Det1", Si, 5, 80, 0.1, 0, 360);
185  TGeoRotation r1;
186  r1.SetAngles(0, 0, 0);
187  TGeoTranslation t1(0, 0, 0);
188  TGeoCombiTrans c1(t1, r1);
189  TGeoHMatrix* h1 = new TGeoHMatrix(c1);
190  top->AddNode(det1, 1, h1);
191  det1->SetLineColor(kGreen);
192  AddSensitiveVolume(det1);
193 
194  TGeoVolume* passive1 = gGeoManager->MakeTubs("Pass1", Si, 5, 120, 10, 0, 360);
195  TGeoRotation rp1;
196  rp1.SetAngles(0, 0, 0);
197  TGeoTranslation tp1(0, 0, 20);
198  TGeoCombiTrans cp1(tp1, rp1);
199  TGeoHMatrix* hp1 = new TGeoHMatrix(cp1);
200  top->AddNode(passive1, 1, hp1);
201  passive1->SetLineColor(kRed);
202 
203  TGeoVolume* det2 = gGeoManager->MakeTubs("Det2", Si, 5, 150, 0.1, 0, 360);
204  TGeoRotation r2;
205  r2.SetAngles(0, 0, 0);
206  TGeoTranslation t2(0, 0, 70);
207  TGeoCombiTrans c2(t2, r2);
208  TGeoHMatrix* h2 = new TGeoHMatrix(c2);
209  top->AddNode(det2, 1, h2);
210  det2->SetLineColor(kGreen);
211  AddSensitiveVolume(det2);
212 
213  TGeoVolume* det3 = gGeoManager->MakeTubs("Det3", Si, 5, 150, 0.1, 0, 360);
214  TGeoRotation r3;
215  r3.SetAngles(0, 0, 0);
216  TGeoTranslation t3(0, 0, 150);
217  TGeoCombiTrans c3(t3, r3);
218  TGeoHMatrix* h3 = new TGeoHMatrix(c3);
219  top->AddNode(det3, 1, h3);
220  det3->SetLineColor(kGreen);
221  AddSensitiveVolume(det3);
222 }
223 
225  Int_t detID,
226  TVector3 pos,
227  TVector3 mom,
228  Double_t time,
229  Double_t length,
230  Double_t eLoss)
231 {
232  TClonesArray& clref = *fNewDetectorPointCollection;
233  Int_t size = clref.GetEntriesFast();
234  return new (clref[size]) NewDetectorPoint(trackID, detID, pos, mom, time, length, eLoss);
235 }
236 
237 FairModule* NewDetector::CloneModule() const { return new NewDetector(*this); }
238 
239 void NewDetector::DefineSensitiveVolumes()
240 {
241  TObjArray* volumes = gGeoManager->GetListOfVolumes();
242  TIter next(volumes);
243  TGeoVolume* volume;
244  while ((volume = static_cast<TGeoVolume*>(next()))) {
245  if (IsSensitive(volume->GetName())) {
246  LOG(debug2) << "Sensitive Volume " << volume->GetName();
247  AddSensitiveVolume(volume);
248  }
249  }
250 }
251 
252 Bool_t NewDetector::IsSensitive(const std::string& name)
253 {
254  if (name.find("Det") != std::string::npos) {
255  return kTRUE;
256  }
257  return kFALSE;
258 }
259 
void AddSensitiveVolume(TGeoVolume *v)
Definition: FairModule.cxx:294
virtual void see Tutorial4 for examples LOG(warn)<< "This function is deprecated. Use FairAlignmentHandler instead
list of container factories
Definition: FairRuntimeDb.h:24
static FairRun * Instance()
Definition: FairRun.cxx:31
virtual FairModule * CloneModule() const
static FairRootManager * Instance()
virtual void Register()
ClassImp(FairEventBuilder)
virtual ~NewDetector()
Definition: NewDetector.cxx:72
virtual Bool_t ProcessHits(FairVolume *v=0)
Definition: NewDetector.cxx:94
virtual void EndOfEvent()
FairParSet * getContainer(const Text_t *)
virtual void Initialize()
virtual void Initialize()
Definition: NewDetector.cxx:80
FairMQExParamsParOne * par
virtual void Reset()
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 Bool_t IsSensitive(const std::string &name)
virtual TClonesArray * GetCollection(Int_t iColl) const
void AddPoint(DetectorId iDet)
NewDetectorPoint * AddHit(Int_t trackID, Int_t detID, TVector3 pos, TVector3 mom, Double_t time, Double_t length, Double_t eLoss)
void RegisterAny(const char *name, T *&obj, Bool_t toFile)
void ConstructGeometry()