BmnRoot
Loading...
Searching...
No Matches
BmnScWall.cxx
Go to the documentation of this file.
1#include "BmnScWall.h"
2
3#include "BmnScWallGeo.h"
4#include "BmnScWallPoint.h"
5#include "CbmStack.h"
6#include "FairGeoInterface.h"
7#include "FairGeoLoader.h"
8#include "FairGeoNode.h"
9#include "FairGeoRootBuilder.h"
10#include "FairMCPoint.h"
11#include "FairRootManager.h"
12#include "FairVolume.h"
13#include "TClonesArray.h"
14#include "TGeoArb8.h"
15#include "TGeoMCGeometry.h"
16#include "TGeoManager.h"
17#include "TLorentzVector.h"
18#include "TParticle.h"
19#include "TVector3.h"
20#include "TVirtualMC.h"
21
22#include <iostream>
23// add on for debug
24// #include "FairGeoG3Builder.h"
25#include "FairRun.h"
26#include "FairRuntimeDb.h"
27#include "TObjArray.h"
28#include "TParticlePDG.h"
29
30// ----- Default constructor -------------------------------------------
32{
33 fScWallCollection = new TClonesArray("BmnScWallPoint");
34 volDetector = 0;
35 fEventID = -1;
36 fCellSmallCutVolId = 0;
37 fCellSmallTrapVolId = 0;
38
39 fCellSmallCutVolId_10mm = 0;
40 fCellSmallTrapVolId_10mm = 0;
41 fCellSmallCutVolId_20mm = 0;
42 fCellSmallTrapVolId_20mm = 0;
43
44 fCellLargeCutVolId = 0;
45 fCellLargeTrapVolId = 0;
46 fVerboseLevel = 1;
47}
48
49// ----- Standard constructor ------------------------------------------
50BmnScWall::BmnScWall(const char* name, Bool_t active)
51 : FairDetector(name, active)
52{
53 fScWallCollection = new TClonesArray("BmnScWallPoint");
54 volDetector = 0;
55 fEventID = -1;
56 fCellSmallCutVolId = 0;
57 fCellSmallTrapVolId = 0;
58 fCellSmallCutVolId_10mm = 0;
59 fCellSmallTrapVolId_10mm = 0;
60 fCellSmallCutVolId_20mm = 0;
61 fCellSmallTrapVolId_20mm = 0;
62 fCellLargeCutVolId = 0;
63 fCellLargeTrapVolId = 0;
64 fVerboseLevel = 1;
65}
66
67// ----- Destructor ----------------------------------------------------
69{
70 if (fScWallCollection) {
71 fScWallCollection->Delete();
72 delete fScWallCollection;
73 }
74}
75
76// ----- Public method Intialize ---------------------------------------
78{
79 // Init function
80
81 FairDetector::Initialize();
82 // FairRun* sim = FairRun::Instance();
83 // FairRuntimeDb* rtdb = sim->GetRuntimeDb();
84
85 fCellSmallCutVolId = gMC->VolId(BmnScWallGeoPar::scwallSensitiveCell_name_small_cutted);
86 fCellSmallCutVolId_10mm = gMC->VolId(BmnScWallGeoPar::scwallSensitiveCell_name_small_cutted_10mm);
87 fCellSmallCutVolId_20mm = gMC->VolId(BmnScWallGeoPar::scwallSensitiveCell_name_small_cutted_20mm);
88
89 fCellSmallTrapVolId = gMC->VolId(BmnScWallGeoPar::scwallSensitiveCell_name_small_trap);
90 fCellSmallTrapVolId_10mm = gMC->VolId(BmnScWallGeoPar::scwallSensitiveCell_name_small_trap_10mm);
91 fCellSmallTrapVolId_20mm = gMC->VolId(BmnScWallGeoPar::scwallSensitiveCell_name_small_trap_20mm);
92
93 fCellLargeCutVolId = gMC->VolId(BmnScWallGeoPar::scwallSensitiveCell_name_large_cutted);
94 fCellLargeTrapVolId = gMC->VolId(BmnScWallGeoPar::scwallSensitiveCell_name_large_trap);
95}
96
97//_____________________________________________________________________________
99{
100 // Returns the hit for the specified slice.
101 // ---
102
103 return (BmnScWallPoint*)fScWallCollection->At(i);
104}
105
106//_____________________________________________________________________________
107BmnScWallPoint* BmnScWall::GetHit(Int_t slice, Int_t cell) const
108{
109 // Returns the hit for the specified vsc and module.
110 // ---
111
112 BmnScWallPoint* hit;
113 Int_t nofHits = fScWallCollection->GetEntriesFast();
114 for (Int_t i = 0; i < nofHits; i++) {
115 hit = GetHit(i);
116 if (hit->GetCopy() == slice && hit->GetCopyMother() == cell)
117 return hit;
118 }
119
120 return 0;
121}
122
123// ----- Public method ProcessHits --------------------------------------
124Bool_t BmnScWall::ProcessHits(FairVolume* vol)
125{
128 Int_t copyNoCellSmallCut = 0;
129 Int_t copyNoCellSmallCut_10mm = 0;
130 Int_t copyNoCellSmallCut_20mm = 0;
131 Int_t copyNoCellSmallTrap = 0;
132 Int_t copyNoCellSmallTrap_10mm = 0;
133 Int_t copyNoCellSmallTrap_20mm = 0;
134 Int_t copyNoCellLargeCut = 0;
135 Int_t copyNoCellLargeTrap = 0;
136 Int_t copyNoSLICECom = 0, copyNoCELLCom = 0;
137
138 TLorentzVector tPos1, tMom1;
139 TLorentzVector tPos, tMom;
140
141 Int_t slice = 0, cell = 0;
142
143 Double_t time = 0;
144 Double_t length = 0;
145
146 // TParticle* part;
147
148 Double_t QCF = 1; // quenching for Birk
149 Double_t BirkConst = 12.6; // 0.126 mm/MeV for polystyrene
150 // 0.126 *(0.1/0.001) = 12.6 cm/GeV
151 //(0.126 mm/MeV - from Wikipedia, 0.07943mm/MeV in Geant4)
152
153 if (gMC->CurrentVolID(copyNoCellSmallCut) != fCellSmallCutVolId
154 && gMC->CurrentVolID(copyNoCellSmallTrap) != fCellSmallTrapVolId
155 && gMC->CurrentVolID(copyNoCellLargeCut) != fCellLargeCutVolId
156 && gMC->CurrentVolID(copyNoCellLargeTrap) != fCellLargeTrapVolId
157 && gMC->CurrentVolID(copyNoCellSmallCut_10mm) != fCellSmallCutVolId_10mm
158 && gMC->CurrentVolID(copyNoCellSmallCut_20mm) != fCellSmallCutVolId_20mm
159 && gMC->CurrentVolID(copyNoCellSmallTrap_10mm) != fCellSmallTrapVolId_10mm
160 && gMC->CurrentVolID(copyNoCellSmallTrap_20mm) != fCellSmallTrapVolId_20mm)
161 {
162 return kFALSE;
163 }
164
165 Int_t ivol = vol->getMCid();
166
167 if (gMC->CurrentVolID(copyNoCellSmallCut) == fCellSmallCutVolId
168 || gMC->CurrentVolID(copyNoCellSmallTrap) == fCellSmallTrapVolId
169 || gMC->CurrentVolID(copyNoCellLargeCut) == fCellLargeCutVolId
170 || gMC->CurrentVolID(copyNoCellLargeTrap) != fCellLargeTrapVolId
171 || gMC->CurrentVolID(copyNoCellSmallCut_10mm) == fCellSmallCutVolId_10mm
172 || gMC->CurrentVolID(copyNoCellSmallCut_20mm) == fCellSmallCutVolId_20mm
173 || gMC->CurrentVolID(copyNoCellSmallTrap_10mm) == fCellSmallTrapVolId_10mm
174 || gMC->CurrentVolID(copyNoCellSmallTrap_20mm) == fCellSmallTrapVolId_20mm)
175 {
176
177 gMC->CurrentVolOffID(0, slice); // scwall
178 gMC->CurrentVolOffID(1, cell);
179 copyNoSLICECom = slice;
180 copyNoCELLCom = cell;
181 }
182
183 // cout <<"2_1_0 " << ivol <<" " <<gMC->CurrentVolOffName(2) << " " <<gMC->CurrentVolOffName(1) << " " <<
184 // gMC->CurrentVolOffName(0) <<std::endl;
185
186 if (gMC->IsTrackEntering()) {
187
188 ResetParameters();
189
190 fELoss = 0.;
191 time = 0.;
192 length = 0.;
193
194 TLorentzVector PosIn;
195 gMC->TrackPosition(PosIn);
196 fPos.SetXYZ(PosIn.X(), PosIn.Y(), PosIn.Z());
197
198 TLorentzVector MomIn;
199 gMC->TrackMomentum(MomIn);
200 fMom.SetXYZ(MomIn.Px(), MomIn.Py(), MomIn.Pz());
201
202 fTime = gMC->TrackTime() * 1.0e09;
203 fLength = gMC->TrackLength();
204 // part = gMC->GetStack()->GetCurrentTrack();
205
206 // std::cout << "ENTER MpdZdc::ProcessHits: TrackID:" <<part->GetPdgCode() <<" " <<fELoss <<" " <<
207 // gMC->CurrentVolPath() << " " << PosIn.Z() <<" " <<MomIn.Pz() <<std::endl;
208 } // if (gMC->IsTrackEntering())
209
210 if (gMC->IsTrackInside()) {
211
212 gMC->TrackPosition(tPos);
213 gMC->TrackMomentum(tMom);
214 length += gMC->TrackStep();
215
216 // fELoss +=gMC->Edep();
217 // Birk corrections
218 if (gMC->TrackStep() > 0)
219 QCF = 1. + (BirkConst / gMC->TrackStep()) * gMC->Edep();
220 else
221 QCF = 1;
222 fELoss += (gMC->Edep()) / QCF;
223
224 time += gMC->TrackTime() * 1.0e09; // nsec
225
226 if (gMC->IsTrackStop() || gMC->IsTrackDisappeared()) {
227
228 // part = gMC->GetStack()->GetCurrentTrack();
229 // Double_t charge = part->GetPDG()->Charge() / 3.;
230
231 // Create BmnZdcPoint
232 fTrackID = gMC->GetStack()->GetCurrentTrackNumber();
233 // time = gMC->TrackTime() * 1.0e09;
234 // length = gMC->TrackLength();
235 // gMC->TrackPosition(tPos);
236 // gMC->TrackMomentum(tMom);
237 // gMC->CurrentVolOffID(2, module);
238 // gMC->CurrentVolOffID(1, slice);
239
240 if (fELoss > 0) {
241 // std::cout << "INSIDE MpdZdc::ProcessHits: TrackID:" <<part->GetPdgCode() <<" " << fTrackID << " "
242 // <<fELoss <<" " << gMC->CurrentVolPath() << " " << tPos.Z() <<" " <<tMom.Pz() <<" " << ivol <<" "
243 // <<gMC->CurrentVolOffName(2) << " " <<gMC->CurrentVolOffName(1) << " " << gMC->CurrentVolOffName(0)
244 // <<std::endl; std::cout << "STOP_DIS MpdZdc::ProcessHits: TrackID:" <<part->GetPdgCode() <<" " <<
245 // fTrackID << " " <<fELoss <<" " << gMC->CurrentVolPath() << " " << tPos.Z() <<" " <<tMom.Pz()
246 // <<std::endl;
247 if (copyNoSLICECom == slice && copyNoCELLCom == cell) {
248 if (!GetHit(slice, cell)) {
249 AddHit(fTrackID, ivol, slice, cell, TVector3(tPos.X(), tPos.Y(), tPos.Z()),
250 TVector3(tMom.Px(), tMom.Py(), tMom.Pz()), time, length, fELoss);
251 } else {
252 GetHit(slice, cell)
253 ->AddCELL(fTrackID, ivol, slice, cell, TVector3(tPos.X(), tPos.Y(), tPos.Z()),
254 TVector3(tMom.Px(), tMom.Py(), tMom.Pz()), time, length, fELoss);
255 }
256 } // if(copyNoSLICECom==slice && copyNoCELLCom==cell)
257
258 } // if(fELoss>0)
259 } // if ( gMC->IsTrackStop() || gMC->IsTrackDisappeared() )
260 } // if ( gMC->IsTrackInside())
261
262 if (gMC->IsTrackExiting()) {
263 // Create BmnZdcPoint
264 fTrackID = gMC->GetStack()->GetCurrentTrackNumber();
265 time += gMC->TrackTime() * 1.0e09;
266 length += gMC->TrackLength();
267
268 // fELoss +=gMC->Edep();
269 // Birk corrections
270 if (gMC->TrackStep() > 0)
271 QCF = 1. + (BirkConst / gMC->TrackStep()) * gMC->Edep();
272 else
273 QCF = 1;
274 fELoss += (gMC->Edep()) / QCF;
275
276 gMC->TrackPosition(tPos);
277 gMC->TrackMomentum(tMom);
278
279 if (fELoss > 0) {
280 // part = gMC->GetStack()->GetCurrentTrack();
281 // std::cout << "EXITING MpdZdc::ProcessHits: TrackID:" <<part->GetPdgCode() <<" " << fTrackID << " "
282 // <<fELoss <<" " << gMC->CurrentVolPath() << " " << tPos.Z() <<" " <<tMom.Pz() <<std::endl;
283 if (copyNoSLICECom == slice && copyNoCELLCom == cell) {
284 if (!GetHit(slice, cell)) {
285 AddHit(fTrackID, ivol, slice, cell, TVector3(tPos.X(), tPos.Y(), tPos.Z()),
286 TVector3(tMom.Px(), tMom.Py(), tMom.Pz()), time, length, fELoss);
287 } else {
288 GetHit(slice, cell)
289 ->AddCELL(fTrackID, ivol, slice, cell, TVector3(tPos.X(), tPos.Y(), tPos.Z()),
290 TVector3(tMom.Px(), tMom.Py(), tMom.Pz()), time, length, fELoss);
291 }
292 } // if(copyNoVTYVECCom==copyNoVTYVEC
293 } // if(fELoss>0)
294 } // if ( gMC->IsTrackExiting()) {
295
296 Int_t points = gMC->GetStack()->GetCurrentTrack()->GetMother(1);
297 points = (points & (~(1 << 30))) | (1 << 30);
298
299 gMC->GetStack()->GetCurrentTrack()->SetMother(1, points);
300
301 ((CbmStack*)gMC->GetStack())->AddPoint(kSCWALL);
302
303 /*//S.Merts
304 fELoss += gMC->Edep();
305 if (gMC->IsTrackExiting() || gMC->IsTrackStop() || gMC->IsTrackDisappeared()) {
306 fTrackID = gMC->GetStack()->GetCurrentTrackNumber();
307 fVolumeID = vol->getMCid();
308 if (fELoss == 0.) {
309 return kFALSE;
310 }
311 AddHit(fTrackID, fVolumeID, fPos, fMom, fTime, fLength, fELoss);
312
313 ((CbmStack*)gMC->GetStack())->AddPoint(kSCWALL);
314 }
315 */
316
317 return kTRUE;
318}
319
320// ----- Public method EndOfEvent -----------------------------------------
322{
323 if (fVerboseLevel)
324 Print("");
325 Reset();
326}
327
328// ----- Public method Register -------------------------------------------
330{
331 FairRootManager::Instance()->Register("ScWallPoint", "ScWall", fScWallCollection, kTRUE);
332}
333
334// ----- Public method GetCollection --------------------------------------
335TClonesArray* BmnScWall::GetCollection(Int_t iColl) const
336{
337 if (iColl == 0)
338 return fScWallCollection;
339 return NULL;
340}
341
342// ----- Public method Print ----------------------------------------------
343void BmnScWall::Print(Option_t*) const
344{
345 Int_t nHits = fScWallCollection->GetEntriesFast();
346 cout << "-I- BmnScWall: " << nHits << " points registered in this event." << endl;
347
348 if (fVerboseLevel > 1)
349 for (Int_t i = 0; i < nHits; i++)
350 (*fScWallCollection)[i]->Print();
351}
352
353// ----- Public method Reset ----------------------------------------------
355{
356 fScWallCollection->Delete();
357}
358
359//------ ConstructGeometry ----------------------------------------------
361{
362 TString fileName = GetGeometryFileName();
363 if (fileName.EndsWith(".root")) {
364 LOG(info) << "Constructing ScWall geometry from ROOT file " << fileName.Data();
365 ConstructRootGeometry();
366 } else
367 LOG(fatal) << "Geometry format of ScWall file " << fileName.Data() << " not supported.";
368}
369
370//------ CheckIfSensitive -----------------------------------------------
371Bool_t BmnScWall::CheckIfSensitive(std::string name)
372{
373 TString tsname = name;
374 if (tsname.Contains("sens"))
375 return kTRUE;
376 return kFALSE;
377}
378
379// ----- Private method AddHit --------------------------------------------
381 Int_t detID,
382 Int_t copyNo,
383 Int_t copyNoMother,
384 TVector3 pos,
385 TVector3 mom,
386 Double_t time,
387 Double_t length,
388 Double_t eLoss)
389{
390 TClonesArray& clref = *fScWallCollection;
391 Int_t size = clref.GetEntriesFast();
392 return new (clref[size]) BmnScWallPoint(trackID, detID, copyNo, copyNoMother, pos, mom, time, length, eLoss);
393}
int i
Definition P4_F32vec4.h:22
@ kSCWALL
static const TString scwallSensitiveCell_name_small_trap_20mm
static const TString scwallSensitiveCell_name_small_trap
static const TString scwallSensitiveCell_name_small_cutted_10mm
static const TString scwallSensitiveCell_name_small_cutted_20mm
static const TString scwallSensitiveCell_name_small_trap_10mm
static const TString scwallSensitiveCell_name_large_cutted
static const TString scwallSensitiveCell_name_small_cutted
static const TString scwallSensitiveCell_name_large_trap
Short_t GetCopy() const
void AddCELL(Int_t trackID, Int_t detID, Int_t idslice, Int_t idcell, TVector3 pos, TVector3 mom, Double_t dt, Double_t dl, Double_t de)
Short_t GetCopyMother() const
BmnScWallPoint * AddHit(Int_t trackID, Int_t detID, Int_t copyNo, Int_t copyNoMother, TVector3 pos, TVector3 mom, Double_t tof, Double_t length, Double_t eLoss)
virtual Bool_t ProcessHits(FairVolume *vol=0)
virtual Bool_t CheckIfSensitive(std::string name)
virtual void Reset()
BmnScWallPoint * GetHit(Int_t i) const
Definition BmnScWall.cxx:98
virtual void ConstructGeometry()
virtual TClonesArray * GetCollection(Int_t iColl) const
virtual void EndOfEvent()
virtual void Print(Option_t *) const
virtual void Register()
virtual ~BmnScWall()
Definition BmnScWall.cxx:68
virtual void Initialize()
Definition BmnScWall.cxx:77