CbmRoot
CbmRich.cxx
Go to the documentation of this file.
1 #include "CbmRich.h"
2 
3 #include "CbmRichPoint.h"
4 
5 #include "CbmGeometryUtils.h"
6 
7 #include "CbmStack.h"
8 #include "FairGeoInterface.h"
9 #include "FairGeoLoader.h"
10 #include "FairGeoNode.h"
11 #include "FairGeoRootBuilder.h"
12 #include "FairRootManager.h"
13 #include "FairRun.h"
14 #include "FairRuntimeDb.h"
15 #include "FairVolume.h"
16 
17 #include "Math/Interpolator.h"
18 #include "TClonesArray.h"
19 #include "TGeoMCGeometry.h"
20 #include "TGeoManager.h"
21 #include "TGeoMatrix.h"
22 #include "TGeoMedium.h"
23 #include "TGeoNode.h"
24 #include "TLorentzVector.h"
25 #include "TObjArray.h"
26 #include "TParticle.h"
27 #include "TParticlePDG.h"
28 #include "TVirtualMC.h"
29 #include <TFile.h>
30 
31 #include "FairGeoBuilder.h"
32 #include "FairGeoMedia.h"
33 
34 #include "TGDMLParse.h"
35 
36 
37 #include <iostream>
38 
39 using std::cout;
40 using std::endl;
41 
42 std::map<TString, TGeoMedium*> CbmRich::fFixedMedia;
43 
45  : FairDetector("RICH", kTRUE, ToIntegralType(ECbmModuleId::kRich))
46  , fPosIndex(0)
47  , fRegisterPhotonsOnSensitivePlane(true)
48  , fRichPoints(new TClonesArray("CbmRichPoint"))
49  , fRichRefPlanePoints(new TClonesArray("CbmRichPoint"))
50  , fRichMirrorPoints(new TClonesArray("CbmRichPoint"))
51  , fRotation(NULL)
52  , fPositionRotation(NULL) {
53  fVerboseLevel = 1;
54 }
55 
56 CbmRich::CbmRich(const char* name,
57  Bool_t active,
58  Double_t px,
59  Double_t py,
60  Double_t pz,
61  Double_t rx,
62  Double_t ry,
63  Double_t rz)
64  : FairDetector(name, active, ToIntegralType(ECbmModuleId::kRich))
65  , fPosIndex(0)
66  , fRegisterPhotonsOnSensitivePlane(true)
67  , fRichPoints(new TClonesArray("CbmRichPoint"))
68  , fRichRefPlanePoints(new TClonesArray("CbmRichPoint"))
69  , fRichMirrorPoints(new TClonesArray("CbmRichPoint"))
70  , fRotation(new TGeoRotation("", rx, ry, rz))
71  , fPositionRotation(new TGeoCombiTrans(px, py, pz, fRotation)) {
72  fVerboseLevel = 1;
73 }
74 
76  if (NULL != fRichPoints) {
77  fRichPoints->Delete();
78  delete fRichPoints;
79  }
80  if (NULL != fRichRefPlanePoints) {
81  fRichRefPlanePoints->Delete();
82  delete fRichRefPlanePoints;
83  }
84 
85  if (NULL != fRichMirrorPoints) {
86  fRichMirrorPoints->Delete();
87  delete fRichMirrorPoints;
88  }
89 }
90 
91 void CbmRich::Initialize() { FairDetector::Initialize(); }
92 
93 Bool_t CbmRich::CheckIfSensitive(std::string name) {
94  //return true;
95  TString volName = name;
96  if (volName.Contains("pmt_pixel") || volName.Contains("sens_plane"))
97  return kTRUE;
98  // mirrors
99  if (volName.Contains("mirror_tile_type")) return kTRUE;
100  return kFALSE;
101 }
102 
103 
104 Bool_t CbmRich::ProcessHits(FairVolume* vol) {
105  Int_t pdgCode = gMC->TrackPid();
106  Int_t iVol = vol->getMCid();
107  TString volName = TString(vol->GetName());
108 
109  if (volName.Contains("pmt_pixel")) {
110  if (gMC->IsTrackEntering()) {
111  // cout << gGeoManager->GetPath() << endl;
112  TParticle* part = gMC->GetStack()->GetCurrentTrack();
113  Double_t charge = part->GetPDG()->Charge() / 3.;
114  Int_t trackID = gMC->GetStack()->GetCurrentTrackNumber();
115  Double_t time = gMC->TrackTime() * 1.0e09;
116  Double_t length = gMC->TrackLength();
117  Double_t eLoss = gMC->Edep();
118 
119  TLorentzVector tPos, tMom;
120  gMC->TrackPosition(tPos);
121  gMC->TrackMomentum(tMom);
122 
123  if (pdgCode == 50000050) { // Cherenkovs only
124 
125  AddHit(trackID,
126  iVol,
127  TVector3(tPos.X(), tPos.Y(), tPos.Z()),
128  TVector3(tMom.Px(), tMom.Py(), tMom.Pz()),
129  time,
130  length,
131  eLoss);
132 
133  // Increment number of RichPoints for this track
134  CbmStack* stack = (CbmStack*) gMC->GetStack();
136  return kTRUE;
137  } else {
138  if (charge == 0.) {
139  return kFALSE; // no neutrals
140  } else { // charged particles
141  AddHit(trackID,
142  iVol,
143  TVector3(tPos.X(), tPos.Y(), tPos.Z()),
144  TVector3(tMom.Px(), tMom.Py(), tMom.Pz()),
145  time,
146  length,
147  eLoss);
148 
149  // Increment number of RichPoints for this track
150  CbmStack* stack = (CbmStack*) gMC->GetStack();
152  return kTRUE;
153  }
154  }
155  }
156  }
157 
158  // Treat imaginary plane in front of the mirrors: Only charged particles at entrance
159  if (volName.Contains("sens_plane")) {
160  // cout << volName << endl;
161  // Collecting points of tracks and imaginary plane intersection
162  if (gMC->IsTrackEntering()) {
163  TParticle* part = gMC->GetStack()->GetCurrentTrack();
164  Double_t charge = part->GetPDG()->Charge() / 3.;
165 
166  if ((fRegisterPhotonsOnSensitivePlane && pdgCode == 50000050)
167  || (!fRegisterPhotonsOnSensitivePlane && charge != 0.)) {
168  Int_t trackID = gMC->GetStack()->GetCurrentTrackNumber();
169  Double_t time = gMC->TrackTime() * 1.0e09;
170  Double_t length = gMC->TrackLength();
171  Double_t eLoss = gMC->Edep();
172  TLorentzVector tPos, tMom;
173 
174  gMC->TrackPosition(tPos);
175  gMC->TrackMomentum(tMom);
176 
177  AddRefPlaneHit(trackID,
178  iVol,
179  TVector3(tPos.X(), tPos.Y(), tPos.Z()),
180  TVector3(tMom.Px(), tMom.Py(), tMom.Pz()),
181  time,
182  length,
183  eLoss);
184 
185  //Increment number of RefPlanePoints for this track
186  CbmStack* stack = (CbmStack*) gMC->GetStack();
188  return kTRUE;
189  } else {
190  return kFALSE;
191  }
192  }
193  }
194 
195  // Treat mirror points
196  Bool_t isMirror = (volName.Contains("mirror_tile_type"));
197  if (isMirror) {
198 
199  // Collecting points
200  if (gMC->IsTrackEntering()) {
201  TParticle* part = gMC->GetStack()->GetCurrentTrack();
202  Double_t charge = part->GetPDG()->Charge() / 3.;
203  if (charge != 0.) {
204 
205  Int_t trackID = gMC->GetStack()->GetCurrentTrackNumber();
206 
207  Double_t time = gMC->TrackTime() * 1.0e09;
208  Double_t length = gMC->TrackLength();
209  Double_t eLoss = gMC->Edep();
210  TLorentzVector tPos, tMom;
211 
212  gMC->TrackPosition(tPos);
213  gMC->TrackMomentum(tMom);
214 
215  AddMirrorHit(trackID,
216  iVol,
217  TVector3(tPos.X(), tPos.Y(), tPos.Z()),
218  TVector3(tMom.Px(), tMom.Py(), tMom.Pz()),
219  time,
220  length,
221  eLoss);
222  return kTRUE;
223  }
224  }
225  }
226 
227  return kFALSE;
228 }
229 
231  if (fVerboseLevel) Print("");
232  Reset();
233 }
234 
236  FairRootManager::Instance()->Register(
237  "RichPoint", "Rich", fRichPoints, kTRUE);
238  FairRootManager::Instance()->Register(
239  "RefPlanePoint", "RichRefPlane", fRichRefPlanePoints, kTRUE);
240  FairRootManager::Instance()->Register(
241  "RichMirrorPoint", "RichMirror", fRichMirrorPoints, kFALSE);
242 }
243 
244 TClonesArray* CbmRich::GetCollection(Int_t iColl) const {
245  if (iColl == 0) return fRichPoints;
246  if (iColl == 1) return fRichRefPlanePoints;
247  if (iColl == 2) return fRichMirrorPoints;
248  return NULL;
249 }
250 
251 void CbmRich::Print(Option_t*) const {
252  Int_t nHits = fRichPoints->GetEntriesFast();
253  LOG(info) << fName << ": " << nHits << " points registered in this event.";
254 
255  if (fVerboseLevel > 1)
256  for (Int_t i = 0; i < nHits; i++)
257  (*fRichPoints)[i]->Print();
258 }
259 
261  fRichPoints->Delete();
262  fRichRefPlanePoints->Delete();
263  fRichMirrorPoints->Delete();
264  fPosIndex = 0;
265 }
266 
267 void CbmRich::CopyClones(TClonesArray* cl1, TClonesArray* cl2, Int_t offset) {
268  Int_t nEntries = cl1->GetEntriesFast();
269  LOG(info) << "CbmRich: " << nEntries << " entries to add.";
270  TClonesArray& clref = *cl2;
271  CbmRichPoint* oldpoint = NULL;
272  for (Int_t i = 0; i < nEntries; i++) {
273  oldpoint = (CbmRichPoint*) cl1->At(i);
274  Int_t index = oldpoint->GetTrackID() + offset;
275  oldpoint->SetTrackID(index);
276  new (clref[fPosIndex]) CbmRichPoint(*oldpoint);
277  fPosIndex++;
278  }
279  LOG(info) << "CbmRich: " << cl2->GetEntriesFast() << " merged entries.";
280 }
281 
283  LOG(info) << "CbmRich::ConstructOpGeometry()";
285 }
286 
288  TString fileName = GetGeometryFileName();
289  if (fileName.EndsWith(".root")) {
291  LOG(info) << "Importing RICH geometry from ROOT file " << fgeoName.Data();
292  TGeoCombiTrans*
293  fCombiTrans {};
294  Cbm::GeometryUtils::ImportRootGeometry(fgeoName, this, fCombiTrans);
295  } else {
296  LOG(info) << "Constructing RICH geometry from ROOT file "
297  << fgeoName.Data();
298  FairModule::ConstructRootGeometry();
299  }
300  } else if (fileName.EndsWith(".gdml")) {
301  LOG(info) << "Constructing RICH geometry from GDML file "
302  << fileName.Data();
304  } else {
305  LOG(fatal) << "Geometry format of RICH file " << fileName.Data()
306  << " not supported. Only ROOT and GDML formats are supported.";
307  }
308 }
309 
311  FairRun* run = FairRun::Instance();
312  Bool_t isGeant4 = std::string(run->GetName()) == "TGeant4";
313  LOG(info) << "CbmRich::SetRichGlassPropertiesForGeant4() Transport engine:"
314  << std::string(run->GetName());
315 
316  // for GEANT3 simulation no need to implement reflectivity
317  if (!isGeant4) {
318  LOG(info) << "CbmRich::SetRichGlassPropertiesForGeant4() fIsGeant4 is "
319  "false. No need to set RICH glass properties. Return.";
320  return;
321  } else {
322  LOG(info) << "CbmRich::SetRichGlassPropertiesForGeant4() fIsGeant4 is "
323  "true. RICH glass properties will be set.";
324  }
325 
326  std::vector<Double_t> energyAr = {
327  6.1997, 5.9045, 5.6361, 5.391, 5.1664, 4.9598, 4.769, 4.5924, 4.4284,
328  4.2757, 4.1331, 3.9998, 3.8748, 3.7574, 3.6469, 3.5427, 3.4443, 3.3512,
329  3.263, 3.1793, 3.0998, 3.0242, 2.9522, 2.8836, 2.818, 2.7554, 2.6955,
330  2.6382, 2.5832, 2.5305, 2.4799, 2.4313, 2.3845, 2.3395, 2.2962, 2.2544,
331  2.2142, 2.1753, 2.1378, 2.1016, 2.0666, 2.0327, 1.9999, 1.9682, 1.9374,
332  1.9076, 1.8787, 1.8507, 1.8234, 1.797, 1.7713, 1.7464, 1.7221, 1.6985,
333  1.6756, 1.6533, 1.6315, 1.6103, 1.5897, 1.5695, 1.5499};
334  // transform to GeV
335  for (size_t i = 0; i < energyAr.size(); i++) {
336  energyAr[i] = 1.e-9 * energyAr[i];
337  }
338 
339  std::vector<Double_t> reflectivityAr = {
340  0.22529, 0.1862, 0.15901, 0.14713, 0.13963, 0.13898, 0.13762, 0.13622,
341  0.13868, 0.13951, 0.14031, 0.14073, 0.1409, 0.13977, 0.14205, 0.14072,
342  0.1396, 0.13826, 0.13665, 0.13513, 0.13463, 0.13287, 0.13182, 0.13084,
343  0.12824, 0.12601, 0.12622, 0.12681, 0.12193, 0.12011, 0.12109, 0.11908,
344  0.11526, 0.11364, 0.11385, 0.12015, 0.11935, 0.11712, 0.1208, 0.12021,
345  0.11909, 0.11783, 0.12257, 0.1215, 0.12199, 0.12494, 0.13101, 0.12089,
346  0.12284, 0.12569, 0.13136, 0.13307, 0.13705, 0.13844, 0.13753, 0.14416,
347  0.14761, 0.14953, 0.15218, 0.15315, 0.15719};
348 
349  for (size_t i = 0; i < reflectivityAr.size(); i++) {
350  reflectivityAr[i] = 1. - reflectivityAr[i];
351  }
352 
353  gMC->DefineOpSurface(
354  "RICHglassSurf", kGlisur, kDielectric_metal, kPolished, 0.0);
355  gMC->SetMaterialProperty("RICHglassSurf",
356  "REFLECTIVITY",
357  energyAr.size(),
358  energyAr.data(),
359  reflectivityAr.data());
360 
361  for (int i = 0; i < 10; i++) {
362  if (gGeoManager->FindVolumeFast(
363  ("mirror_tile_type" + std::to_string(i)).c_str())
364  == nullptr)
365  continue;
366  gMC->SetSkinSurface(("RICHglassSkin" + std::to_string(i)).c_str(),
367  ("mirror_tile_type" + std::to_string(i)).c_str(),
368  "RICHglassSurf");
369  }
370 }
371 
372 void CbmRich::ConstructGdmlGeometry(TGeoMatrix* geoMatrix) {
373  TFile* old = gFile;
374  TGDMLParse parser;
375  TGeoVolume* gdmlTop;
376 
377  // Before importing GDML
378  Int_t maxInd = gGeoManager->GetListOfMedia()->GetEntries() - 1;
379 
380  gdmlTop = parser.GDMLReadFile(GetGeometryFileName());
381 
382  // Cheating - reassigning media indices after GDML import (need to fix this in TGDMLParse class!!!)
383  // for (Int_t i=0; i<gGeoManager->GetListOfMedia()->GetEntries(); i++)
384  // gGeoManager->GetListOfMedia()->At(i)->Dump();
385  // After importing GDML
386  Int_t j = gGeoManager->GetListOfMedia()->GetEntries() - 1;
387  Int_t curId;
388  TGeoMedium* m;
389  do {
390  m = (TGeoMedium*) gGeoManager->GetListOfMedia()->At(j);
391  curId = m->GetId();
392  m->SetId(curId + maxInd);
393  j--;
394  } while (curId > 1);
395  // LOG(debug) << "====================================================================";
396  // for (Int_t i=0; i<gGeoManager->GetListOfMedia()->GetEntries(); i++)
397  // gGeoManager->GetListOfMedia()->At(i)->Dump();
398 
399  Int_t newMaxInd = gGeoManager->GetListOfMedia()->GetEntries() - 1;
400 
401  gGeoManager->GetTopVolume()->AddNode(gdmlTop, 1, geoMatrix);
402  ExpandNodeForGdml(gGeoManager->GetTopVolume()->GetNode(
403  gGeoManager->GetTopVolume()->GetNdaughters() - 1));
404 
405  for (Int_t k = maxInd + 1; k < newMaxInd + 1; k++) {
406  TGeoMedium* medToDel =
407  (TGeoMedium*) (gGeoManager->GetListOfMedia()->At(maxInd + 1));
408  LOG(debug) << " removing media " << medToDel->GetName() << " with id "
409  << medToDel->GetId() << " (k=" << k << ")";
410  gGeoManager->GetListOfMedia()->Remove(medToDel);
411  }
412  gGeoManager->SetAllIndex();
413 
414  gFile = old;
415 }
416 
417 void CbmRich::ExpandNodeForGdml(TGeoNode* node) {
418  LOG(debug)
419  << "----------------------------------------- ExpandNodeForGdml for node "
420  << node->GetName();
421 
422  TGeoVolume* curVol = node->GetVolume();
423 
424  LOG(debug) << " volume: " << curVol->GetName();
425 
426  if (curVol->IsAssembly()) {
427  LOG(debug) << " skipping volume-assembly";
428  } else {
429 
430  TGeoMedium* curMed = curVol->GetMedium();
431  TGeoMaterial* curMat = curVol->GetMaterial();
432  TGeoMedium* curMedInGeoManager = gGeoManager->GetMedium(curMed->GetName());
433  TGeoMaterial* curMatOfMedInGeoManager = curMedInGeoManager->GetMaterial();
434  TGeoMaterial* curMatInGeoManager =
435  gGeoManager->GetMaterial(curMat->GetName());
436 
437  // Current medium and material assigned to the volume from GDML
438  LOG(debug2) << " curMed\t\t\t\t" << curMed << "\t" << curMed->GetName()
439  << "\t" << curMed->GetId();
440  LOG(debug2) << " curMat\t\t\t\t" << curMat << "\t" << curMat->GetName()
441  << "\t" << curMat->GetIndex();
442 
443  // Medium and material found in the gGeoManager - either the pre-loaded one or one from GDML
444  LOG(debug2) << " curMedInGeoManager\t\t" << curMedInGeoManager << "\t"
445  << curMedInGeoManager->GetName() << "\t"
446  << curMedInGeoManager->GetId();
447  LOG(debug2) << " curMatOfMedInGeoManager\t\t" << curMatOfMedInGeoManager
448  << "\t" << curMatOfMedInGeoManager->GetName() << "\t"
449  << curMatOfMedInGeoManager->GetIndex();
450  LOG(debug2) << " curMatInGeoManager\t\t" << curMatInGeoManager << "\t"
451  << curMatInGeoManager->GetName() << "\t"
452  << curMatInGeoManager->GetIndex();
453 
454  TString matName = curMat->GetName();
455  TString medName = curMed->GetName();
456 
457  if (curMed->GetId() != curMedInGeoManager->GetId()) {
458  if (fFixedMedia.find(medName) == fFixedMedia.end()) {
459  LOG(debug) << " Medium needs to be fixed";
460  fFixedMedia[medName] = curMedInGeoManager;
461  Int_t ind = curMat->GetIndex();
462  gGeoManager->RemoveMaterial(ind);
463  LOG(debug) << " removing material " << curMat->GetName()
464  << " with index " << ind;
465  for (Int_t i = ind; i < gGeoManager->GetListOfMaterials()->GetEntries();
466  i++) {
467  TGeoMaterial* m =
468  (TGeoMaterial*) gGeoManager->GetListOfMaterials()->At(i);
469  m->SetIndex(m->GetIndex() - 1);
470  }
471 
472  LOG(debug) << " Medium fixed";
473  } else {
474  LOG(debug) << " Already fixed medium found in the list ";
475  }
476 
477  } else {
478  if (fFixedMedia.find(medName) == fFixedMedia.end()) {
479  LOG(debug) << " There is no correct medium in the memory yet";
480 
481  FairGeoLoader* geoLoad = FairGeoLoader::Instance();
482  FairGeoInterface* geoFace = geoLoad->getGeoInterface();
483  FairGeoMedia* geoMediaBase = geoFace->getMedia();
484  // FairGeoBuilder* geobuild = geoLoad->getGeoBuilder();
485 
486  FairGeoMedium* curMedInGeo = geoMediaBase->getMedium(medName);
487  if (curMedInGeo == 0) {
488  LOG(fatal) << medName << " Media not found in Geo file.";
492  } else {
493  LOG(debug) << " Found media in Geo file" << medName;
494  // Int_t nmed = geobuild->createMedium(curMedInGeo);
495  fFixedMedia[medName] =
496  (TGeoMedium*) gGeoManager->GetListOfMedia()->Last();
497  gGeoManager->RemoveMaterial(curMatOfMedInGeoManager->GetIndex());
498  LOG(debug) << " removing material "
499  << curMatOfMedInGeoManager->GetName() << " with index "
500  << curMatOfMedInGeoManager->GetIndex();
501  for (Int_t i = curMatOfMedInGeoManager->GetIndex();
502  i < gGeoManager->GetListOfMaterials()->GetEntries();
503  i++) {
504  TGeoMaterial* m =
505  (TGeoMaterial*) gGeoManager->GetListOfMaterials()->At(i);
506  m->SetIndex(m->GetIndex() - 1);
507  }
508  }
509 
510  if (curMedInGeo->getSensitivityFlag()) {
511  LOG(debug) << " Adding sensitive " << curVol->GetName();
512  AddSensitiveVolume(curVol);
513  }
514 
515  } else {
516  LOG(debug) << " Already fixed medium found in the list";
517  LOG(debug) << "!!! Sensitivity: " << fFixedMedia[medName]->GetParam(0);
518  if (fFixedMedia[medName]->GetParam(0) == 1) {
519  LOG(debug) << " Adding sensitive " << curVol->GetName();
520  AddSensitiveVolume(curVol);
521  }
522  }
523  }
524 
525  curVol->SetMedium(fFixedMedia[medName]);
526  gGeoManager->SetAllIndex();
527 
528  // gGeoManager->GetListOfMaterials()->Print();
529  // gGeoManager->GetListOfMedia()->Print();
530  }
531 
533  if (curVol->GetNdaughters() != 0) {
534  TObjArray* NodeChildList = curVol->GetNodes();
535  TGeoNode* curNodeChild;
536  for (Int_t j = 0; j < NodeChildList->GetEntriesFast(); j++) {
537  curNodeChild = (TGeoNode*) NodeChildList->At(j);
538  ExpandNodeForGdml(curNodeChild);
539  }
540  }
541 }
542 
544  Int_t detID,
545  TVector3 pos,
546  TVector3 mom,
547  Double_t time,
548  Double_t length,
549  Double_t eLoss) {
550  TClonesArray& clref = *fRichPoints;
551  Int_t size = clref.GetEntriesFast();
552  return new (clref[size])
553  CbmRichPoint(trackID, detID, pos, mom, time, length, eLoss);
554 }
555 
557  Int_t detID,
558  TVector3 pos,
559  TVector3 mom,
560  Double_t time,
561  Double_t length,
562  Double_t eLoss) {
563  TClonesArray& clref = *fRichRefPlanePoints;
564  Int_t tsize = clref.GetEntriesFast();
565  return new (clref[tsize])
566  CbmRichPoint(trackID, detID, pos, mom, time, length, eLoss);
567 }
568 
570  Int_t detID,
571  TVector3 pos,
572  TVector3 mom,
573  Double_t time,
574  Double_t length,
575  Double_t eLoss) {
576  TClonesArray& clref = *fRichMirrorPoints;
577  Int_t tsize = clref.GetEntriesFast();
578  return new (clref[tsize])
579  CbmRichPoint(trackID, detID, pos, mom, time, length, eLoss);
580 }
581 
CbmRich::Initialize
virtual void Initialize()
Initialize detector. Stores volume IDs for RICH detector and mirror.
Definition: CbmRich.cxx:91
CbmRichPoint.h
CbmRich::GetCollection
virtual TClonesArray * GetCollection(Int_t iColl) const
Return hit collection.
Definition: CbmRich.cxx:244
CbmRich::ProcessHits
virtual Bool_t ProcessHits(FairVolume *vol=0)
Defines the action to be taken when a step is inside the active volume. Creates CbmRichPoints and Cbm...
Definition: CbmRich.cxx:104
CbmRich::fRichRefPlanePoints
TClonesArray * fRichRefPlanePoints
Definition: CbmRich.h:183
i
int i
Definition: L1/vectors/P4_F32vec4.h:25
CbmRich
Defines the active detector RICH. Constructs the geometry and creates MCPoints.
Definition: CbmRich.h:38
ECbmModuleId
ECbmModuleId
Definition: CbmDefs.h:33
CbmRich::fRichPoints
TClonesArray * fRichPoints
Definition: CbmRich.h:182
CbmRich::CopyClones
virtual void CopyClones(TClonesArray *cl1, TClonesArray *cl2, Int_t offset)
Copies the hit collection with a given track index offset.
Definition: CbmRich.cxx:267
CbmRich::AddRefPlaneHit
CbmRichPoint * AddRefPlaneHit(Int_t trackID, Int_t detID, TVector3 pos, TVector3 mom, Double_t time, Double_t length, Double_t eLoss)
Adds a RichRefPlanePoint to the TClonesArray.
Definition: CbmRich.cxx:556
CbmRich::AddHit
CbmRichPoint * AddHit(Int_t trackID, Int_t detID, TVector3 pos, TVector3 mom, Double_t time, Double_t length, Double_t eLoss)
Adds a RichPoint to the TClonesArray.
Definition: CbmRich.cxx:543
CbmStack.h
CbmRich::Reset
virtual void Reset()
Clears the hit collection.
Definition: CbmRich.cxx:260
CbmRich::fFixedMedia
static std::map< TString, TGeoMedium * > fFixedMedia
Definition: CbmRich.h:189
CbmStack::AddPoint
void AddPoint(ECbmModuleId iDet)
Definition: CbmStack.cxx:383
CbmRich::ExpandNodeForGdml
void ExpandNodeForGdml(TGeoNode *node)
Assign materials by taking description from medoa.geo and not from GDML for a certain node.
Definition: CbmRich.cxx:417
Cbm::GeometryUtils::IsNewGeometryFile
Bool_t IsNewGeometryFile(TString &filename)
Definition: CbmGeometryUtils.cxx:133
CbmRich::CheckIfSensitive
virtual Bool_t CheckIfSensitive(std::string name)
Definition: CbmRich.cxx:93
CbmRich::fRichMirrorPoints
TClonesArray * fRichMirrorPoints
Definition: CbmRich.h:184
CbmRich.h
Defines the active detector RICH. Constructs the geometry and creates MCPoints.
ECbmModuleId::kRich
@ kRich
Ring-Imaging Cherenkov Detector.
ECbmModuleId::kRef
@ kRef
Reference plane.
ClassImp
ClassImp(CbmConverterManager) InitStatus CbmConverterManager
Definition: CbmConverterManager.cxx:12
CbmRich::Print
virtual void Print(Option_t *) const
Screen output of hit collection.
Definition: CbmRich.cxx:251
CbmRich::EndOfEvent
virtual void EndOfEvent()
If verbosity level is set, print hit collection at the end of the event and resets it afterwards.
Definition: CbmRich.cxx:230
CbmRich::SetRichGlassPropertiesForGeant4
void SetRichGlassPropertiesForGeant4()
Set Cherenkov propeties for RICH mirror.
Definition: CbmRich.cxx:310
CbmRich::fRegisterPhotonsOnSensitivePlane
Bool_t fRegisterPhotonsOnSensitivePlane
Definition: CbmRich.h:180
CbmRich::ConstructGdmlGeometry
void ConstructGdmlGeometry(TGeoMatrix *geoMatrix)
Construct geometry from GDML file.
Definition: CbmRich.cxx:372
ToIntegralType
constexpr auto ToIntegralType(T enumerator) -> typename std::underlying_type< T >::type
Definition: CbmDefs.h:24
m
__m128 m
Definition: L1/vectors/P4_F32vec4.h:26
CbmRich::ConstructGeometry
virtual void ConstructGeometry()
Construct geometry. Currently ROOT and ASCII formats are supported. The concrete method for geometry ...
Definition: CbmRich.cxx:287
pos
TVector3 pos
Definition: CbmMvdSensorDigiToHitTask.cxx:60
Cbm::GeometryUtils::ImportRootGeometry
void ImportRootGeometry(TString &filename, FairModule *mod, TGeoMatrix *mat)
Definition: CbmGeometryUtils.cxx:140
CbmRich::fPositionRotation
TGeoCombiTrans * fPositionRotation
Definition: CbmRich.h:192
CbmRich::fPosIndex
Int_t fPosIndex
Definition: CbmRich.h:177
CbmRich::~CbmRich
virtual ~CbmRich()
Destructor.
Definition: CbmRich.cxx:75
CbmRich::Register
virtual void Register()
Registers the hit collection in the ROOT manager.
Definition: CbmRich.cxx:235
CbmStack
Definition: CbmStack.h:45
CbmRich::CbmRich
CbmRich()
Default constructor.
Definition: CbmRich.cxx:44
CbmRichPoint
Definition: CbmRichPoint.h:24
CbmRich::AddMirrorHit
CbmRichPoint * AddMirrorHit(Int_t trackID, Int_t detID, TVector3 pos, TVector3 mom, Double_t time, Double_t length, Double_t eLoss)
Adds a RichMirrorPoint to the TClonesArray.
Definition: CbmRich.cxx:569
CbmRich::ConstructOpGeometry
void ConstructOpGeometry()
Put some optical properties.
Definition: CbmRich.cxx:282
CbmGeometryUtils.h