BmnRoot
Loading...
Searching...
No Matches
BmnBd.cxx
Go to the documentation of this file.
1#include "BmnBd.h"
2#include "BmnBdGeo.h"
3#include "BmnBdGeoPar.h"
4#include "CbmStack.h"
5
6#include "FairGeoInterface.h"
7#include "FairGeoLoader.h"
8#include "FairGeoNode.h"
9#include "FairRootManager.h"
10#include "FairVolume.h"
11#include "FairRuntimeDb.h"
12#include "FairRun.h"
13
14#include "TClonesArray.h"
15#include "TLorentzVector.h"
16#include "TParticlePDG.h"
17#include "TVector3.h"
18#include "TObjArray.h"
19
20#include <iostream>
21
22// ----- Default constructor -------------------------------------------
24 fBdCollection = new TClonesArray("BmnBdPoint");
25 volDetector = 0;
26 fPosIndex = 0;
27 //fpreflag = 0;
28 //fpostflag = 0;
29 //fEventID = -1;
30 fVerboseLevel = 1;
31}
32
33// ----- Standard constructor ------------------------------------------
34BmnBd::BmnBd(const char* name, Bool_t active)
35 : FairDetector(name, active) {
36 fBdCollection = new TClonesArray("BmnBdPoint");
37 fPosIndex = 0;
38 volDetector = 0;
39 //fpreflag = 0;
40 //fpostflag = 0;
41 //fEventID=-1;
42 fVerboseLevel = 1;
43}
44
45// ----- Destructor ----------------------------------------------------
47 if (fBdCollection) {
48 fBdCollection->Delete();
49 delete fBdCollection;
50 }
51}
52
53// ----- Public method Initialize ---------------------------------------
55 FairDetector::Initialize();
56 //FairRun* sim = FairRun::Instance();
57 //FairRuntimeDb* rtdb = sim->GetRuntimeDb();
58}
59
60// ------ Begin of the event ----------------------------------------
62{
63}
64
65// ----- Public method ProcessHits --------------------------------------
66Bool_t BmnBd::ProcessHits(FairVolume* vol) {
67 // if (TMath::Abs(gMC->TrackCharge()) <= 0) return kFALSE;
68
69 static const Double_t dP = 1.032 ;
70 static const Double_t BirkC1 = 0.013/dP;
71 static const Double_t BirkC2 = 9.6e-6/(dP * dP);
72
73 static Double_t lightYield;
74 static Double_t timeIn;
75 static Double_t timeOut;
76 static Double_t lengthtrack;
77 //static Double_t de_dx;
78
79 Int_t ivol = vol->getMCid();
80 //TLorentzVector tPos1, tMom1;
81 //TLorentzVector tPos, tMom;
82
83 //#define EDEBUG
84#ifdef EDEBUG
85 static Int_t lEDEBUGcounter = 0;
86 if (lEDEBUGcounter < 1)
87 std::cout << "EDEBUG-- BmnBd::ProcessHits: entered" << gMC->CurrentVolPath() << endl;
88#endif
89
90 if(gMC->IsTrackEntering())
91 {
92 ResetParameters();
93
94 fELoss = 0.;
95 fIsPrimary = 0;
96 fCharge = -1;
97 fPdgId = 0;
98 lightYield = 0.;
99 lengthtrack = 0.;
100
101 TLorentzVector PosIn;
102 gMC->TrackPosition(PosIn);
103 fPosIn.SetXYZ(PosIn.X(), PosIn.Y(), PosIn.Z());
104
105 TLorentzVector MomIn;
106 gMC->TrackMomentum(MomIn);
107 fMomIn.SetXYZ(MomIn.Px(), MomIn.Py(), MomIn.Pz());
108
109 timeIn = gMC->TrackTime() * 1.0e09;
110
111#ifdef EDEBUG
112 //gMC->TrackPosition(tPos1);
113 //gMC->TrackMomentum(tMom1);
114#endif
115
116 TParticle* part = 0;
117 part = gMC->GetStack()->GetCurrentTrack();
118
119 if (part) {
120 fIsPrimary = (Int_t)part->IsPrimary();
121 fCharge = (Int_t)part->GetPDG()->Charge() / 3.;
122 fPdgId = (Int_t)part->GetPdgCode();
123 }
124 }
125
126 if (gMC->TrackCharge() != 0)
127 {
128 Double_t BirkC1Mod = 0;
129
130 if (TMath::Abs(gMC->TrackCharge()) >= 2)
131 BirkC1Mod = BirkC1 * 7.2/12.6;
132 else
133 BirkC1Mod = BirkC1;
134
135 if (gMC->TrackStep() > 0)
136 {
137 Double_t dedxcm = gMC -> Edep() * 1000./gMC -> TrackStep(); //[MeV/cm]
138 Double_t curLightYield = gMC -> Edep() * 1000./(1. + BirkC1Mod * dedxcm + BirkC2 * dedxcm * dedxcm);
139 curLightYield /= 1000.; //[GeV]
140 lightYield += curLightYield;
141 }
142 }
143
144 Double_t eLoss = gMC->Edep();
145 if (eLoss > 0)
146 fELoss += eLoss;
147
148 Double_t slengthtrack = gMC-> TrackStep();
149 lengthtrack += slengthtrack;
150
151 if ((gMC->IsTrackExiting() || gMC->IsTrackStop() || gMC->IsTrackDisappeared()) && (fELoss > 0))
152 {
153 TLorentzVector PosOut;
154 gMC->TrackPosition(PosOut);
155 fPosOut.SetXYZ(PosOut.X(), PosOut.Y(), PosOut.Z());
156
157 TLorentzVector MomOut;
158 gMC->TrackMomentum(MomOut);
159 fMomOut.SetXYZ(MomOut.Px(), MomOut.Py(), MomOut.Pz());
160
161 timeOut = gMC->TrackTime() * 1.0e09;
162/*
163 TArrayI processesID;
164
165 ofstream fout("process.txt", ios::app);
166
167 //if (gMC->TrackPid() == 11 && gMC->Etot() <= 0.00006){
168
169 gMC->StepProcesses(processesID);
170 fout << gMC->TrackPid() << " " << gMC->Edep() << " " ;
171 for (Int_t i = 0; i<processesID.GetSize(); i++){
172
173 fout << TMCProcessName[processesID[i]] << " ";
174 }
175 fout << endl;
176
177 //fout<<" "<<gMC->TrackPid() <<" " <<gMC->Edep()<<" "<< TMCProcessName[processesID[i]]<<endl;
178 fout.close();
179*/
180
181#ifndef EDEBUG
182 if (fELoss == 0. ) return kFALSE;
183#endif
184
185 fTrackID = gMC->GetStack()->GetCurrentTrackNumber();
186
187 Double_t time = gMC->TrackTime() * 1.0e09;
188 Double_t length = gMC->TrackLength();
189
190 Int_t copyNo;
191 /*Int_t ivol1 = */gMC->CurrentVolID(copyNo);
192
193 if (copyNo < 41)
194 {
195 Int_t iCell;
196 gMC->CurrentVolOffID(1, iCell);
197
198#ifdef EDEBUG
199 if (lEDEBUGcounter < 100)
200 {
201 std::cout << "EDEBUG-- BmnBd::ProcessHits: TrackID:" << fTrackID <<
202 // " ELoss: " << fELoss <<
203 // " particle: " << (part->GetName()) <<
204 " " << gMC->CurrentVolPath() << " " << tPos.Z() <<
205 // " " << (gMC->GetStack()->GetCurrentTrack()->GetMother(1)) <<
206 " " << ivol << "=="<< gMC->CurrentVolID(copyNo) << ","<< copyNo <<
207 " " << gMC->CurrentVolOffID(1,iCell) << " " << iCell <<
208 " " << gMC->CurrentVolOffName(1) << " " << gMC->CurrentVolOffName(0) <<
209 // " " << vol->getRealName() << " " << gMC->CurrentVolPath() <<
210 // " ivol,iCell,copyNo= " << ivol << ","<< iCell << ","<< copyNo <<
211 // " " << vol->getRealName() << " "<< gMC->CurrentVolName() << " "<< gMC->CurrentVolPath() <<
212 // " " << ivol << ","<< vol->getVolumeId() << " : "<< gMC->CurrentVolID(copyNo) << ","<< copyNo <<
213 // " "<< gMC->CurrentVolOffName(2) << " "<< gMC->CurrentVolOffName(3) <<
214 std::endl;
215 lEDEBUGcounter++;
216 }
217
218 AddHit(fTrackID, ivol, copyNo, fPosIn, fPosOut, fMomIn, fMomOut,
219 time, length, fELoss, fIsPrimary, fCharge, fPdgId, lightYield, timeIn, timeOut, lengthtrack);
220#else
221 AddHit(fTrackID, ivol, copyNo, fPosIn, fPosOut, fMomIn, fMomOut,
222 time, length, fELoss, fIsPrimary, fCharge, fPdgId, lightYield, timeIn, timeOut, lengthtrack);
223#endif
224
226// Int_t nBdPoints = (points & (1<<30)) >> 30;
227// nBdPoints ++;
228// if (nBdPoints > 1) nBdPoints = 1;
229// points = ( points & ( ~ (1<<30) ) ) | (nBdPoints << 30);
232
233 ((CbmStack*)gMC->GetStack())->AddPoint(kBD);
234 }
235 }
236
237 return kTRUE;
238
239 #undef EDEBUG
240}
241
242// ----- Public method EndOfEvent -----------------------------------------
244 if (fVerboseLevel) Print("");
245 Reset();
246}
247
248// ----- Public method Register -------------------------------------------
250 FairRootManager::Instance()->Register("BdPoint","BD", fBdCollection, kTRUE);
251}
252
253// ----- Public method GetCollection --------------------------------------
254TClonesArray* BmnBd::GetCollection(Int_t iColl) const {
255 if (iColl == 0) return fBdCollection;
256
257 return NULL;
258}
259
260// ----- Public method Print ----------------------------------------------
261void BmnBd::Print(Option_t*) const {
262 Int_t nHits = fBdCollection->GetEntriesFast();
263 cout << "-I- BmnBd: " << nHits << " points registered in this event." << endl;
264
265 if (fVerboseLevel > 1)
266 for (Int_t i=0; i < nHits; i++) (*fBdCollection)[i]->Print();
267}
268
269// ----- Public method Reset ----------------------------------------------
271 fBdCollection->Delete();
272 fPosIndex = 0;
273}
274
275// guarda in FairRootManager::CopyClones
276// ----- Public method CopyClones -----------------------------------------
277void BmnBd::CopyClones(TClonesArray* cl1, TClonesArray* cl2, Int_t offset ) {
278 Int_t nEntries = cl1->GetEntriesFast();
279 TClonesArray& clref = *cl2;
280 BmnBdPoint* oldpoint = NULL;
281 for (Int_t i = 0; i < nEntries; i++)
282 {
283 oldpoint = (BmnBdPoint*) cl1->At(i);
284 Int_t index = oldpoint->GetTrackID() + offset;
285 oldpoint->SetTrackID(index);
286 new (clref[fPosIndex]) BmnBdPoint(*oldpoint);
287 fPosIndex++;
288 }
289 cout << " -I- BmnBd: " << cl2->GetEntriesFast() << " merged entries." << endl;
290}
291
292//--------------------------------------------------------------------------------------------------------------------------------------
294{
295 TString fileName = GetGeometryFileName();
296 if (fileName.EndsWith(".root"))
297 {
298 LOG(info) << "Constructing BD geometry from ROOT file" << fileName.Data();
299 ConstructRootGeometry();
300 }
301 else if ( fileName.EndsWith(".geo") )
302 {
303 LOG(info) << "Constructing BD geometry from ASCII file" << fileName.Data();
305 }
306 else LOG(fatal) << "Geometry format of BD file " << fileName.Data() << " not supported.";
307}
308
310{
311 FairGeoLoader* geoLoad = FairGeoLoader::Instance();
312 FairGeoInterface* geoFace = geoLoad->getGeoInterface();
313 BmnBdGeo* BdGeo = new BmnBdGeo();
314 BdGeo->setGeomFile(GetGeometryFileName());
315 geoFace->addGeoModule(BdGeo);
316
317 Bool_t rc = geoFace->readSet(BdGeo);
318 if (rc) BdGeo->create(geoLoad->getGeoBuilder());
319 TList* volList = BdGeo->getListOfVolumes();
320
321 // store geo parameter
322 FairRun *fRun = FairRun::Instance();
323 FairRuntimeDb *rtdb = FairRun::Instance()->GetRuntimeDb();
324 BmnBdGeoPar* par = (BmnBdGeoPar*)(rtdb->getContainer("BmnBdGeoPar"));
325 TObjArray *fSensNodes = par->GetGeoSensitiveNodes();
326 TObjArray *fPassNodes = par->GetGeoPassiveNodes();
327
328 TListIter iter(volList);
329 FairGeoNode *node = nullptr;
330 FairGeoVolume *aVol = nullptr;
331
332 while ( (node = (FairGeoNode*)iter.Next()) )
333 {
334 aVol = dynamic_cast<FairGeoVolume*>(node);
335 cout<<aVol->GetName()<<endl;
336
337 if (node->isSensitive()) fSensNodes->AddLast(aVol);
338 else fPassNodes->AddLast(aVol);
339 }
340
341 par->setChanged();
342 par->setInputVersion(fRun->GetRunId(),1);
343 ProcessNodes( volList );
344}
345
346//--------------------------------------------------------------------------------------------------------------------------------------
347Bool_t BmnBd::CheckIfSensitive(std::string name)
348{
349 TString tsname = name;
350 if (tsname.Contains("Active")) return kTRUE;
351
352 return kFALSE;
353}
354
355// ----- Private method AddHit --------------------------------------------
356BmnBdPoint* BmnBd::AddHit(Int_t trackID, Int_t detID, Int_t copyNo, TVector3 posIn, TVector3 posOut,
357 TVector3 momIn, TVector3 momOut, Double_t time, Double_t length, Double_t eLoss, Bool_t isPrimary,
358 Double_t charge, Int_t pdgId, Double_t lightYield, Double_t timeIn, Double_t timeOut, Double_t lengthtrack)
359{
360 TClonesArray& clref = *fBdCollection;
361 Int_t size = clref.GetEntriesFast();
362 return new(clref[size]) BmnBdPoint(trackID, detID, copyNo, posIn, posOut, momIn, momOut, time, length, eLoss, isPrimary,
363 charge, pdgId, lightYield, timeIn, timeOut, lengthtrack);
364}
int i
Definition P4_F32vec4.h:22
@ kBD
TObjArray * GetGeoPassiveNodes()
Definition BmnBdGeoPar.h:18
TObjArray * GetGeoSensitiveNodes()
Definition BmnBdGeoPar.h:17
virtual void BeginEvent()
Definition BmnBd.cxx:61
virtual void CopyClones(TClonesArray *cl1, TClonesArray *cl2, Int_t offset)
Definition BmnBd.cxx:277
BmnBdPoint * AddHit(Int_t trackID, Int_t detID, Int_t copyNo, TVector3 posIn, TVector3 posOut, TVector3 momIn, TVector3 momOut, Double_t tof, Double_t length, Double_t eLoss, Bool_t isPrimary, Double_t charge, Int_t pdgId, Double_t lightYield, Double_t timeIn, Double_t timeOut, Double_t lengthtrack)
Definition BmnBd.cxx:356
virtual TClonesArray * GetCollection(Int_t iColl) const
Definition BmnBd.cxx:254
virtual void Print(Option_t *) const
Definition BmnBd.cxx:261
virtual void Initialize()
Definition BmnBd.cxx:54
BmnBd()
Definition BmnBd.cxx:23
virtual void Register()
Definition BmnBd.cxx:249
virtual void ConstructAsciiGeometry()
Definition BmnBd.cxx:309
virtual void Reset()
Definition BmnBd.cxx:270
virtual Bool_t CheckIfSensitive(std::string name)
Definition BmnBd.cxx:347
virtual void EndOfEvent()
Definition BmnBd.cxx:243
virtual Bool_t ProcessHits(FairVolume *vol=0)
Definition BmnBd.cxx:66
virtual ~BmnBd()
Definition BmnBd.cxx:46
virtual void ConstructGeometry()
Definition BmnBd.cxx:293