BmnRoot
Loading...
Searching...
No Matches
BmnCSC.cxx
Go to the documentation of this file.
1#include "BmnCSC.h"
2
3#include "BmnCSCPoint.h"
4#include "CbmStack.h"
5#include "FairGeoMedia.h"
6#include "FairGeoRootBuilder.h"
7#include "FairRootManager.h"
8#include "FairVolume.h"
9#include "TFile.h"
10#include "TGDMLParse.h"
11#include "TGeoManager.h"
12#include "TParticle.h"
13#include "TParticlePDG.h"
14#include "TRegexp.h"
15#include "TVirtualMC.h"
16
17#include <iostream>
18
19using namespace std;
20
22 : FairDetector("CSC", kTRUE)
23{
24 fPointCollection = new TClonesArray("BmnCSCPoint");
25 fPosIndex = 0;
26 fVerboseLevel = 1;
27 ResetParameters();
28}
29
30BmnCSC::BmnCSC(const char* name, Bool_t active)
31 : FairDetector(name, active)
32{
33 fPointCollection = new TClonesArray("BmnCSCPoint");
34 fPosIndex = 0;
35 fVerboseLevel = 1;
36 ResetParameters();
37}
38
40{
41 if (fPointCollection) {
42 fPointCollection->Delete();
43 delete fPointCollection;
44 }
45}
46
47//------------------------------------------------------------------------------
48
49Bool_t BmnCSC::ProcessHits(FairVolume* vol)
50{
51
52 // Set parameters at entrance of volume. Reset ELoss.
53 if (gMC->IsTrackEntering()) {
54
55 ResetParameters();
56 fELoss = 0.;
57 fTime = gMC->TrackTime() * 1.0e09;
58 fLength = gMC->TrackLength();
59 fIsPrimary = 0;
60 fCharge = -1;
61 fPdgId = 0;
62
63 TLorentzVector PosIn;
64 gMC->TrackPosition(PosIn);
65 fPosIn.SetXYZ(PosIn.X(), PosIn.Y(), PosIn.Z());
66
67 TLorentzVector MomIn;
68 gMC->TrackMomentum(MomIn);
69 fMomIn.SetXYZ(MomIn.X(), MomIn.Y(), MomIn.Z());
70
71 TParticle* part = 0;
72 part = gMC->GetStack()->GetCurrentTrack();
73 if (part) {
74 fIsPrimary = (Int_t)part->IsPrimary();
75 fCharge = (Int_t)part->GetPDG()->Charge();
76 fPdgId = (Int_t)part->GetPdgCode();
77 }
78
79 fVolumeID = vol->getMCid();
80
81 fTrackID = gMC->GetStack()->GetCurrentTrackNumber();
82 }
83
84 // Sum energy loss for all steps in the active volume
85 fELoss += gMC->Edep();
86
87 // Create BmnCSCPoint at EXIT of active volume;
88 if ((gMC->IsTrackExiting() || gMC->IsTrackStop() || gMC->IsTrackDisappeared()) && fELoss > 0) {
89
90 TLorentzVector PosOut;
91 gMC->TrackPosition(PosOut);
92 fPosOut.SetXYZ(PosOut.X(), PosOut.Y(), PosOut.Z());
93
94 TLorentzVector MomOut;
95 gMC->TrackMomentum(MomOut);
96 fMomOut.SetXYZ(MomOut.X(), MomOut.Y(), MomOut.Z());
97
98 // correction step to avoid the seg. violation error due to invalid memory access
99 TVector3 diff_pos = fPosIn - fPosOut;
100
101 if (diff_pos.Mag() < 0.001)
102 return kFALSE; // ignore points produced by tracks with zero length inside the current sens. volume
103 if (fMomOut.Mag() == 0)
104 return kFALSE; // ignore points produced by tracks with zero momentum inside the current sens. volume
105
106 TVector3 corr_step = fMomOut;
107 corr_step.SetMag(0.001); // 10 um
108 TVector3 pos = fPosOut;
109 fPosOut = pos - corr_step;
110 gGeoManager->FindNode(fPosOut[0], fPosOut[1], fPosOut[2]);
111
112 // Determine station and module numbers for the current hit ------------
113 Int_t stationNum = -1; // current station number (default)
114 Int_t moduleNum = -1; // current module number (default)
115
116 TGeoVolume* currentVolume = gGeoManager->GetCurrentVolume();
117 TString currentVolumeName = currentVolume->GetName();
118 TRegexp expr = "^CSCSensor_module[0-9]+_station[0-9]+.*";
119 if (currentVolumeName.Contains(expr)) {
120 TRegexp mod_expt = "module[0-9]+";
121 TRegexp stat_expt = "station[0-9]+";
122
123 moduleNum = TString(TString(currentVolumeName(mod_expt))(TRegexp("[0-9]+"))).Atoi();
124 stationNum = TString(TString(currentVolumeName(stat_expt))(TRegexp("[0-9]+"))).Atoi();
125 }
126
127 if (stationNum == -1 || moduleNum == -1)
128 return kFALSE; // check if the current point has incorrect indices
129
130 /*BmnCSCPoint *p = */ AddHit(fTrackID, fVolumeID, fPosIn, fPosOut, fMomIn, fMomOut, fTime, fLength, fELoss,
131 fIsPrimary, fCharge, fPdgId, stationNum, moduleNum);
132
133 ((CbmStack*)gMC->GetStack())->AddPoint(kCSC);
134 }
135
136 return kTRUE;
137}
138
140{
141 if (fVerboseLevel)
142 Print("");
143 fPointCollection->Clear();
144 fPosIndex = 0;
145}
146
148{
149 FairRootManager::Instance()->Register("CSCPoint", "CSC", fPointCollection, kTRUE);
150}
151
152TClonesArray* BmnCSC::GetCollection(Int_t iColl) const
153{
154 if (iColl == 0)
155 return fPointCollection;
156 return NULL;
157}
158
159void BmnCSC::Print(Option_t*) const
160{
161 Int_t nHits = fPointCollection->GetEntriesFast();
162 cout << "-I- BmnCSC: " << nHits << " points registered in this event." << endl;
163
164 if (fVerboseLevel > 1) {
165 for (Int_t i = 0; i < nHits; i++) {
166 (*fPointCollection)[i]->Print();
167 }
168 }
169}
170
172{
173 fPointCollection->Clear();
174 ResetParameters();
175}
176
177void BmnCSC::CopyClones(TClonesArray* cl1, TClonesArray* cl2, Int_t offset)
178{
179 Int_t nEntries = cl1->GetEntriesFast();
180 cout << "-I- BmnCSC: " << nEntries << " entries to add." << endl;
181 TClonesArray& clref = *cl2;
182 BmnCSCPoint* oldpoint = NULL;
183
184 for (Int_t i = 0; i < nEntries; i++) {
185 oldpoint = (BmnCSCPoint*)cl1->At(i);
186 Int_t index = oldpoint->GetTrackID() + offset;
187 oldpoint->SetTrackID(index);
188 new (clref[fPosIndex]) BmnCSCPoint(*oldpoint);
189 fPosIndex++;
190 }
191
192 cout << "-I- BmnCSC: " << cl2->GetEntriesFast() << " merged entries." << endl;
193}
194
196{
197 TString fileName = GetGeometryFileName();
198
199 if (fileName.EndsWith(".root")) {
200 LOG(info) << "Constructing CSC geometry from ROOT file " << fileName.Data();
201 ConstructRootGeometry();
202 }
203
204 else if (fileName.EndsWith(".gdml"))
205 {
206 LOG(info) << "Constructing CSC geometry from GDML file " << fileName.Data();
207 ConstructGDMLGeometry(nullptr);
208 }
209
210 else
211 {
212 LOG(fatal) << "Geometry format of CSC file " << fileName.Data() << " not supported.";
213 }
214}
215
216// ----- ConstructGDMLGeometry ---------------------------------------------
218{
219 TFile* old = gFile;
220 TGDMLParse parser;
221 TGeoVolume* gdmlTop;
222
223 // Before importing GDML
224 Int_t maxInd = gGeoManager->GetListOfMedia()->GetEntries() - 1;
225
226 gdmlTop = parser.GDMLReadFile(GetGeometryFileName());
227
228 // Cheating - reassigning media indices after GDML import (need to fix this in TGDMLParse class!!!)
229 // for (Int_t i=0; i<gGeoManager->GetListOfMedia()->GetEntries(); i++)
230 // gGeoManager->GetListOfMedia()->At(i)->Dump();
231 // After importing GDML
232 Int_t j = gGeoManager->GetListOfMedia()->GetEntries() - 1;
233 Int_t curId;
234 TGeoMedium* m;
235 do {
236 m = (TGeoMedium*)gGeoManager->GetListOfMedia()->At(j);
237 curId = m->GetId();
238 m->SetId(curId + maxInd);
239 j--;
240 } while (curId > 1);
241
242 // LOG(debug) << "====================================================================";
243 // for (Int_t i=0; i<gGeoManager->GetListOfMedia()->GetEntries(); i++)
244 // gGeoManager->GetListOfMedia()->At(i)->Dump();
245
246 Int_t newMaxInd = gGeoManager->GetListOfMedia()->GetEntries() - 1;
247
248 gGeoManager->GetTopVolume()->AddNode(gdmlTop, 1, 0);
249 ExpandNodeForGdml(gGeoManager->GetTopVolume()->GetNode(gGeoManager->GetTopVolume()->GetNdaughters() - 1));
250
251 for (Int_t k = maxInd + 1; k < newMaxInd + 1; k++) {
252 TGeoMedium* medToDel = (TGeoMedium*)(gGeoManager->GetListOfMedia()->At(maxInd + 1));
253 LOG(debug) << " removing media " << medToDel->GetName() << " with id " << medToDel->GetId() << " (k=" << k
254 << ")";
255 gGeoManager->GetListOfMedia()->Remove(medToDel);
256 }
257 gGeoManager->SetAllIndex();
258
259 gFile = old;
260}
261
262void BmnCSC::ExpandNodeForGdml(TGeoNode* node)
263{
264 LOG(debug) << "----------------------------------------- ExpandNodeForGdml for node " << node->GetName();
265
266 TGeoVolume* curVol = node->GetVolume();
267
268 LOG(debug) << " volume: " << curVol->GetName();
269
270 if (curVol->IsAssembly()) {
271 LOG(debug) << " skipping volume-assembly";
272 } else {
273 TGeoMedium* curMed = curVol->GetMedium();
274 TGeoMaterial* curMat = curVol->GetMaterial();
275 TGeoMedium* curMedInGeoManager = gGeoManager->GetMedium(curMed->GetName());
276 TGeoMaterial* curMatOfMedInGeoManager = curMedInGeoManager->GetMaterial();
277 TGeoMaterial* curMatInGeoManager = gGeoManager->GetMaterial(curMat->GetName());
278
279 // Current medium and material assigned to the volume from GDML
280 LOG(debug2) << " curMed\t\t\t\t" << curMed << "\t" << curMed->GetName() << "\t" << curMed->GetId();
281 LOG(debug2) << " curMat\t\t\t\t" << curMat << "\t" << curMat->GetName() << "\t" << curMat->GetIndex();
282
283 // Medium and material found in the gGeoManager - either the pre-loaded one or one from GDML
284 LOG(debug2) << " curMedInGeoManager\t\t" << curMedInGeoManager << "\t" << curMedInGeoManager->GetName()
285 << "\t" << curMedInGeoManager->GetId();
286 LOG(debug2) << " curMatOfMedInGeoManager\t\t" << curMatOfMedInGeoManager << "\t"
287 << curMatOfMedInGeoManager->GetName() << "\t" << curMatOfMedInGeoManager->GetIndex();
288 LOG(debug2) << " curMatInGeoManager\t\t" << curMatInGeoManager << "\t" << curMatInGeoManager->GetName()
289 << "\t" << curMatInGeoManager->GetIndex();
290
291 TString matName = curMat->GetName();
292 TString medName = curMed->GetName();
293
294 if (curMed->GetId() != curMedInGeoManager->GetId()) {
295 if (fFixedMedia.find(medName) == fFixedMedia.end()) {
296 LOG(debug) << " Medium needs to be fixed";
297 fFixedMedia[medName] = curMedInGeoManager;
298 Int_t ind = curMat->GetIndex();
299 gGeoManager->RemoveMaterial(ind);
300 LOG(debug) << " removing material " << curMat->GetName() << " with index " << ind;
301 for (Int_t i = ind; i < gGeoManager->GetListOfMaterials()->GetEntries(); i++) {
302 TGeoMaterial* m = (TGeoMaterial*)gGeoManager->GetListOfMaterials()->At(i);
303 m->SetIndex(m->GetIndex() - 1);
304 }
305
306 LOG(debug) << " Medium fixed";
307 } else {
308 LOG(debug) << " Already fixed medium found in the list ";
309 }
310 } else {
311 if (fFixedMedia.find(medName) == fFixedMedia.end()) {
312 LOG(debug) << " There is no correct medium in the memory yet";
313
314 FairGeoLoader* geoLoad = FairGeoLoader::Instance();
315 FairGeoInterface* geoFace = geoLoad->getGeoInterface();
316 FairGeoMedia* geoMediaBase = geoFace->getMedia();
317 FairGeoBuilder* geobuild = geoLoad->getGeoBuilder();
318
319 FairGeoMedium* curMedInGeo = geoMediaBase->getMedium(medName);
320 if (curMedInGeo == 0) {
321 LOG(fatal) << " Media not found in Geo file: " << medName;
325 } else {
326 LOG(debug) << " Found media in Geo file" << medName;
327 /*Int_t nmed = */ geobuild->createMedium(curMedInGeo);
328 fFixedMedia[medName] = (TGeoMedium*)gGeoManager->GetListOfMedia()->Last();
329 gGeoManager->RemoveMaterial(curMatOfMedInGeoManager->GetIndex());
330 LOG(debug) << " removing material " << curMatOfMedInGeoManager->GetName() << " with index "
331 << curMatOfMedInGeoManager->GetIndex();
332 for (Int_t i = curMatOfMedInGeoManager->GetIndex();
333 i < gGeoManager->GetListOfMaterials()->GetEntries(); i++)
334 {
335 TGeoMaterial* m = (TGeoMaterial*)gGeoManager->GetListOfMaterials()->At(i);
336 m->SetIndex(m->GetIndex() - 1);
337 }
338 }
339
340 if (curMedInGeo->getSensitivityFlag()) {
341 LOG(debug) << " Adding sensitive " << curVol->GetName();
342 AddSensitiveVolume(curVol);
343 }
344 } else {
345 LOG(debug) << " Already fixed medium found in the list";
346 LOG(debug) << "!!! Sensitivity: " << fFixedMedia[medName]->GetParam(0);
347 if (fFixedMedia[medName]->GetParam(0) == 1) {
348 LOG(debug) << " Adding sensitive " << curVol->GetName();
349 AddSensitiveVolume(curVol);
350 }
351 }
352 }
353
354 curVol->SetMedium(fFixedMedia[medName]);
355 gGeoManager->SetAllIndex();
356
357 // gGeoManager->GetListOfMaterials()->Print();
358 // gGeoManager->GetListOfMedia()->Print();
359 }
360
362 if (curVol->GetNdaughters() != 0) {
363 TObjArray* NodeChildList = curVol->GetNodes();
364 TGeoNode* curNodeChild;
365 for (Int_t j = 0; j < NodeChildList->GetEntriesFast(); j++) {
366 curNodeChild = (TGeoNode*)NodeChildList->At(j);
367 ExpandNodeForGdml(curNodeChild);
368 }
369 }
370}
371
372//--- Check if Sensitive -------------------------------------------------------
373Bool_t BmnCSC::CheckIfSensitive(std::string name)
374{
375 TString tsname = name;
376 if (tsname.Contains("CSCSensor")) {
377 return kTRUE;
378 }
379 return kFALSE;
380}
381//------------------------------------------------------------------------------
382
383//------------------------------------------------------------------------------
384BmnCSCPoint* BmnCSC::AddHit(Int_t trackID,
385 Int_t detID,
386 TVector3 posIn,
387 TVector3 posOut,
388 TVector3 momIn,
389 TVector3 momOut,
390 Double_t time,
391 Double_t length,
392 Double_t eLoss,
393 Int_t isPrimary,
394 Double_t charge,
395 Int_t pdgId,
396 Int_t stationNum,
397 Int_t moduleNum)
398{
399
400 TClonesArray& clref = *fPointCollection;
401 Int_t size = clref.GetEntriesFast();
402 // std::cout << "ELoss: " << eLoss << "\n";
403 new (clref[size])
404 BmnCSCPoint(trackID, detID, posIn, posOut, momIn, momOut, time, length, eLoss, isPrimary, charge, pdgId);
405
406 BmnCSCPoint* hit = (BmnCSCPoint*)clref.At(clref.GetEntriesFast() - 1);
407 hit->SetStation(stationNum);
408 hit->SetModule(moduleNum);
409 return hit;
410}
411//------------------------------------------------------------------------------
int i
Definition P4_F32vec4.h:22
__m128 m
Definition P4_F32vec4.h:27
@ kCSC
void SetStation(Int_t station)
Definition BmnCSCPoint.h:58
void SetModule(Int_t module)
Definition BmnCSCPoint.h:59
virtual Bool_t ProcessHits(FairVolume *vol=0)
Definition BmnCSC.cxx:49
virtual void Register()
Definition BmnCSC.cxx:147
virtual void Reset()
Definition BmnCSC.cxx:171
virtual void ConstructGeometry()
Definition BmnCSC.cxx:195
virtual TClonesArray * GetCollection(Int_t iColl) const
Definition BmnCSC.cxx:152
void ExpandNodeForGdml(TGeoNode *node)
Definition BmnCSC.cxx:262
virtual void CopyClones(TClonesArray *cl1, TClonesArray *cl2, Int_t offset)
Definition BmnCSC.cxx:177
virtual ~BmnCSC()
Definition BmnCSC.cxx:39
virtual void EndOfEvent()
Definition BmnCSC.cxx:139
virtual Bool_t CheckIfSensitive(std::string name)
Definition BmnCSC.cxx:373
BmnCSC()
Definition BmnCSC.cxx:21
virtual void Print(Option_t *) const
Definition BmnCSC.cxx:159
virtual void ConstructGDMLGeometry(TGeoMatrix *)
Definition BmnCSC.cxx:217
STL namespace.