CbmRoot
CbmLitTrackingGeometryConstructor.cxx
Go to the documentation of this file.
1 
8 #include "CbmUtils.h"
10 
11 #include "CbmHistManager.h"
12 #include "FairRunAna.h"
13 #include "FairRuntimeDb.h"
14 
15 #include "TDirectory.h"
16 #include "TFile.h"
17 #include "TGeoArb8.h"
18 #include "TGeoBBox.h"
19 #include "TGeoCone.h"
20 #include "TGeoManager.h"
21 #include "TGeoMatrix.h"
22 #include "TGeoMedium.h"
23 #include "TGeoPgon.h"
24 #include "TGeoShape.h"
25 #include "TGeoSphere.h"
26 #include "TGeoVolume.h"
27 #include "TProfile2D.h"
28 #include "TRegexp.h"
29 
30 #include <algorithm>
31 #include <iostream>
32 #include <limits>
33 #include <map>
34 #include <set>
35 #include <sstream>
36 
37 using std::cout;
38 using std::map;
39 using std::set;
40 
42  : fGeo(NULL)
43  , fNofTrdStations(-1)
44  , fNofMuchStations(-1)
45  , fNofMvdStations(-1)
46  , fNofStsStations(-1)
47  , fDet() {
48  fGeo = gGeoManager;
50 }
51 
53 
56  static CbmLitTrackingGeometryConstructor instance;
57  return &instance;
58 }
59 
62  GetMuchLayout(layout);
63 }
64 
67  GetMuchLayout(layout);
68 }
69 
70 template<class T>
73  cout << "Getting MUCH layout for parallel version of tracking...\n";
74 
75  // Read file with TProfile2D containing silicon equivalent of the material
76  TString parDir =
77  TString(gSystem->Getenv("VMCWORKDIR")) + TString("/parameters");
78  TString matBudgetFile = parDir + "/littrack/much_v12b.silicon.root";
79  TFile* oldFile = gFile;
80  TDirectory* oldDirectory = gDirectory;
81  TFile* file = new TFile(matBudgetFile, "READ");
82  CbmHistManager hm;
83  hm.ReadFromFile(file);
84 
85  CbmLitFieldGridCreator gridCreator;
86 
87  TGeoNode* much = nullptr;
88  TObjArray* nodes = gGeoManager->GetTopVolume()->GetNodes();
89  for (Int_t iNode = 0; iNode < nodes->GetEntriesFast(); iNode++) {
90  TGeoNode* node = (TGeoNode*) nodes->At(iNode);
91  if (TString(node->GetName())
92  .Contains("much", TString::kIgnoreCase)) { // Top MUCH node
93  much = node;
94  break;
95  }
96  }
97  assert(much);
98  TObjArray* muchNodes = much->GetNodes();
99  Int_t currentStation = 0;
100  Int_t currentLayer = 0;
101  for (Int_t iMuchNode = 0; iMuchNode < muchNodes->GetEntriesFast();
102  iMuchNode++) {
103  const TGeoNode* muchNode =
104  static_cast<const TGeoNode*>(muchNodes->At(iMuchNode));
105 
106  if (TString(muchNode->GetName()).Contains("station")) {
107  TObjArray* layerNodes = muchNode->GetNodes();
108 
109  for (Int_t iLayerNode = 0; iLayerNode < layerNodes->GetEntriesFast();
110  iLayerNode++) {
111  const TGeoNode* layer =
112  static_cast<const TGeoNode*>(layerNodes->At(iLayerNode));
113 
114  Double_t z = much->GetMatrix()->GetTranslation()[2]
115  + muchNode->GetMatrix()->GetTranslation()[2]
116  + layer->GetMatrix()->GetTranslation()[2];
117 
118  // Convert material for this station
119  TProfile2D* profile =
120  (iLayerNode == 0)
121  ? hm.P2("hrl_ThicknessSilicon_MuchAbsorber_"
122  + Cbm::ToString<Int_t>(currentStation) + "_P2")
123  : hm.P2("hrl_ThicknessSilicon_Much_"
124  + Cbm::ToString<Int_t>(currentLayer) + "_P2");
125  //profile->Rebin2D(200, 200);
127  ConvertTProfile2DToLitMaterialGrid(profile, &material);
128 
130  vs.SetMaterial(material);
131 
133  station.AddVirtualStation(vs);
134 
135  layout.AddStation(station);
136 
137 
138  currentLayer++;
139  }
140  currentStation++;
141  }
142  }
143 
144  gFile = oldFile;
145  gDirectory = oldDirectory;
146  file->Close();
147  delete file;
148 
149  cout << layout;
150  cout << "Finish getting MUCH layout for parallel version of tracking\n";
151 }
152 
155  GetTrdLayout(layout);
156 }
157 
160  GetTrdLayout(layout);
161 }
162 
163 template<class T>
166  cout << "Getting TRD layout for parallel version of tracking...\n";
167 
168  // Read file with TProfile2D containing silicon equivalent of the material
169  TString parDir =
170  TString(gSystem->Getenv("VMCWORKDIR")) + TString("/parameters");
171  TString matBudgetFile = parDir + "/littrack/trd_v13p_3e.silicon.root";
172  TFile* oldFile = gFile;
173  TDirectory* oldDirectory = gDirectory;
174  TFile* file = new TFile(matBudgetFile, "READ");
175  CbmHistManager hm;
176  hm.ReadFromFile(file);
177 
178  CbmLitFieldGridCreator gridCreator;
179 
180  Double_t startZPosition = 100.; // Last STS station
181  // Double_t detectorZPosition = 400.; // We want to extrapolate up to here using virtual stations
182  Double_t dZ = 10.; // Distance between neighboring virtual stations
183 
184  lit::parallel::LitMaterialGrid richMaterial;
185  GetRichMaterial(&richMaterial);
186 
187  // Virtual stations
188  Int_t nofVirtualStations = 31;
189  for (Int_t iStation = 0; iStation < nofVirtualStations; iStation++) {
190  Double_t z = startZPosition + iStation * dZ;
191  lit::parallel::LitFieldGrid fieldGrid;
192  gridCreator.CreateGrid(z, fieldGrid);
193 
195  station.SetZ(z);
196  station.SetField(fieldGrid);
197  if (iStation == 10) station.SetMaterial(richMaterial);
198 
199  layout.AddVirtualStation(station);
200  }
201 
202  // Stations
203  Int_t nofStations = GetNofTrdStations();
204  for (Int_t iStation = 0; iStation < nofStations; iStation++) {
205  // Convert material for this station
206  TProfile2D* profile = hm.P2("hrl_ThicknessSilicon_Trd_"
207  + Cbm::ToString<Int_t>(iStation) + "_P2");
208  //profile->Rebin2D(200, 200);
210  ConvertTProfile2DToLitMaterialGrid(profile, &material);
211 
213  vs.SetMaterial(material);
214 
216  station.AddVirtualStation(vs);
217 
218  layout.AddStation(station);
219  }
220 
221  gFile = oldFile;
222  gDirectory = oldDirectory;
223  file->Close();
224  delete file;
225 
226  cout << layout;
227  cout << "Finish getting TRD layout for parallel version of tracking\n";
228 }
229 
231  lit::parallel::LitMaterialGrid* material) {
232  if (!fDet.GetDet(ECbmModuleId::kRich)) return;
233  // Read file with TProfile2D containing silicon equivalent of the material
234  TString parDir =
235  TString(gSystem->Getenv("VMCWORKDIR")) + TString("/parameters");
236  TString matBudgetFile = parDir + "/littrack/rich_v08a.silicon.root";
237  TFile* oldFile = gFile;
238  TDirectory* oldDirectory = gDirectory;
239  TFile* file = new TFile(matBudgetFile, "READ");
240  CbmHistManager hm;
241  hm.ReadFromFile(file);
242 
243  TProfile2D* profile = hm.P2("hrl_ThicknessSilicon_Rich_P2");
244  // profile->Rebin2D(200, 200);
245  ConvertTProfile2DToLitMaterialGrid(profile, material, 3.);
246 
247  gFile = oldFile;
248  gDirectory = oldDirectory;
249  file->Close();
250  delete file;
251 }
252 
254  const TProfile2D* profile,
256  Double_t maximumValue) {
257  // Int_t nofBinsX = profile->GetNbinsX();
258  // Int_t nofBinsY = profile->GetNbinsY();
259  // vector<vector<fscal> >material(nofBinsX);
260  // for (Int_t i = 0; i < nofBinsX; i++) material[i].resize(nofBinsY);
261  // for (Int_t iX = 1; iX <= nofBinsX; iX++) {
262  // for (Int_t iY = 1; iY <= nofBinsY; iY++) {
263  // Double_t content = profile->GetBinContent(iX, iY);
264  // if (maximumValue > 0 && content > maximumValue) content = maximumValue;
265  // material[iX - 1][iY - 1] = content;
266  // }
267  // }
268  // Double_t xmin = profile->GetXaxis()->GetXmin();
269  // Double_t xmax = profile->GetXaxis()->GetXmax();
270  // Double_t ymin = profile->GetYaxis()->GetXmin();
271  // Double_t ymax = profile->GetYaxis()->GetXmax();
272  // grid->SetMaterial(material, xmin, xmax, ymin, ymax, nofBinsX, nofBinsY);
273 
274  Int_t nofBinsX = profile->GetNbinsX();
275  Int_t nofBinsY = profile->GetNbinsY();
276  Int_t minShrinkBinX = std::numeric_limits<Double_t>::max();
277  Int_t maxShrinkBinX = std::numeric_limits<Double_t>::min();
278  Int_t minShrinkBinY = std::numeric_limits<Double_t>::max();
279  Int_t maxShrinkBinY = std::numeric_limits<Double_t>::min();
280  Bool_t isSet = false;
281  for (Int_t iBinX = 1; iBinX <= nofBinsX; iBinX++) {
282  for (Int_t iBinY = 1; iBinY <= nofBinsY; iBinY++) {
283  Double_t content = profile->GetBinContent(iBinX, iBinY);
284  if (content != 0.) {
285  minShrinkBinX = std::min(iBinX, minShrinkBinX);
286  maxShrinkBinX = std::max(iBinX, maxShrinkBinX);
287  minShrinkBinY = std::min(iBinY, minShrinkBinY);
288  maxShrinkBinY = std::max(iBinY, maxShrinkBinY);
289  isSet = true;
290  }
291  }
292  }
293 
294  Int_t nofShrinkBinsX = maxShrinkBinX - minShrinkBinX + 1;
295  Int_t nofShrinkBinsY = maxShrinkBinY - minShrinkBinY + 1;
296  vector<vector<fscal>> material(nofShrinkBinsX);
297  for (Int_t i = 0; i < nofShrinkBinsX; i++)
298  material[i].resize(nofShrinkBinsY);
299 
300  for (Int_t iX = minShrinkBinX; iX <= maxShrinkBinX; iX++) {
301  for (Int_t iY = minShrinkBinY; iY <= maxShrinkBinY; iY++) {
302  Double_t content = profile->GetBinContent(iX, iY);
303  if (maximumValue > 0 && content > maximumValue) content = maximumValue;
304  material[iX - minShrinkBinX][iY - minShrinkBinY] = content;
305  }
306  }
307  Double_t xmin = profile->GetXaxis()->GetBinLowEdge(minShrinkBinX);
308  Double_t xmax = profile->GetXaxis()->GetBinUpEdge(maxShrinkBinX);
309  Double_t ymin = profile->GetYaxis()->GetBinLowEdge(minShrinkBinY);
310  Double_t ymax = profile->GetYaxis()->GetBinUpEdge(maxShrinkBinY);
311  grid->SetMaterial(
312  material, xmin, xmax, ymin, ymax, nofShrinkBinsX, nofShrinkBinsY);
313 }
314 
316  static Bool_t firstTime = true;
317  if (firstTime) {
318  Int_t layerCounter = 0;
319  TObjArray* topNodes = fGeo->GetTopNode()->GetNodes();
320  for (Int_t iTopNode = 0; iTopNode < topNodes->GetEntriesFast();
321  iTopNode++) {
322  TGeoNode* topNode = static_cast<TGeoNode*>(topNodes->At(iTopNode));
323  if (TString(topNode->GetName()).Contains("trd")) {
324  TObjArray* layers = topNode->GetNodes();
325  for (Int_t iLayer = 0; iLayer < layers->GetEntriesFast(); iLayer++) {
326  TGeoNode* layer = static_cast<TGeoNode*>(layers->At(iLayer));
327  if (TString(layer->GetName()).Contains("layer")) { layerCounter++; }
328  }
329  }
330  }
331  fNofTrdStations = layerCounter;
332  firstTime = false;
333  }
334  return fNofTrdStations;
335 }
336 
338  static Bool_t firstTime = true;
339  if (firstTime) {
340  fNofMuchStations = 0;
341  TGeoNode* much = nullptr;
342  TObjArray* nodes = gGeoManager->GetTopVolume()->GetNodes();
343  for (Int_t iNode = 0; iNode < nodes->GetEntriesFast(); iNode++) {
344  TGeoNode* node = (TGeoNode*) nodes->At(iNode);
345  if (TString(node->GetName())
346  .Contains("much", TString::kIgnoreCase)) { // Top MUCH node
347  much = node;
348  break;
349  }
350  }
351  if (NULL == much) { // No MUCH detector return 0 stations
352  firstTime = false;
353  return fNofMuchStations;
354  }
355  TObjArray* muchNodes = much->GetNodes();
356  for (Int_t iMuchNode = 0; iMuchNode < muchNodes->GetEntriesFast();
357  iMuchNode++) {
358  const TGeoNode* muchNode =
359  static_cast<const TGeoNode*>(muchNodes->At(iMuchNode));
360  if (TString(muchNode->GetName()).Contains("station")) {
361  // old structure defined in ASCII geometry
362  if (TString(muchNode->GetName()).Contains("muchstation")) {
363  TObjArray* layerNodes = muchNode->GetNodes();
364  fNofMuchStations += layerNodes->GetEntriesFast();
365  } else {
366  // New structure defined in ROOT geometry files
367  TObjArray* muchSubNodes = muchNode->GetNodes();
368  for (Int_t iMuchSubNode = 0;
369  iMuchSubNode < muchSubNodes->GetEntriesFast();
370  iMuchSubNode++) {
371  const TGeoNode* muchSubNode =
372  static_cast<const TGeoNode*>(muchSubNodes->At(iMuchSubNode));
373  TObjArray* layerNodes = muchSubNode->GetNodes();
374  fNofMuchStations += layerNodes->GetEntriesFast();
375  }
376  }
377  }
378  }
379  firstTime = false;
380  }
381  return fNofMuchStations;
382 }
383 
385  static Bool_t firstTime = true;
386  if (firstTime) {
387  fNofMuchAbsorbers = 0;
388  const TGeoNode* much = static_cast<const TGeoNode*>(
389  fGeo->GetTopNode()->GetNodes()->FindObject("much_0"));
390  if (NULL == much) { // No MUCH detector return 0 stations
391  firstTime = false;
392  return fNofMuchAbsorbers;
393  }
394  TObjArray* muchNodes = much->GetNodes();
395  for (Int_t iMuchNode = 0; iMuchNode < muchNodes->GetEntriesFast();
396  iMuchNode++) {
397  const TGeoNode* muchNode =
398  static_cast<const TGeoNode*>(muchNodes->At(iMuchNode));
399  if (TString(muchNode->GetName()).Contains("absorber")) {
400  // old structure defined in ASCII geometry
401  if (TString(muchNode->GetName()).Contains("muchabsorber")) {
403  } else {
404  // New structure defined in ROOT geometry files
405  TObjArray* muchSubNodes = muchNode->GetNodes();
406  fNofMuchAbsorbers += muchSubNodes->GetEntriesFast();
407  }
408  }
409  }
410  firstTime = false;
411  }
412  return fNofMuchAbsorbers;
413 }
414 
417 }
418 
420  static Bool_t firstTime = true;
421  if (firstTime) {
422  fNofMvdStations = 0;
423  TGeoNode* topNode = gGeoManager->GetTopVolume()->FindNode("pipevac1_0");
424  if (topNode) {
425  TObjArray* nodes = topNode->GetNodes();
426  for (Int_t iNode = 0; iNode < nodes->GetEntriesFast(); iNode++) {
427  TGeoNode* node = (TGeoNode*) nodes->At(iNode);
428  if (TString(node->GetName()).Contains("mvdstation")) {
429  fNofMvdStations++;
430  }
431  }
432  } else {
433  TObjArray* nodes = gGeoManager->GetTopNode()->GetNodes();
434  for (Int_t iNode = 0; iNode < nodes->GetEntriesFast(); iNode++) {
435  TGeoNode* node = (TGeoNode*) nodes->At(iNode);
436  TString nodeName = node->GetName();
437  nodeName.ToLower();
438  if (nodeName.Contains("pipe")) {
439  TObjArray* nodes2 = node->GetNodes();
440  for (Int_t iiNode = 0; iiNode < nodes2->GetEntriesFast(); iiNode++) {
441  TGeoNode* node2 = (TGeoNode*) nodes2->At(iiNode);
442  TString nodeName2 = node2->GetName();
443  nodeName2.ToLower();
444  if (nodeName2.Contains("pipevac1")) {
445  // check if there is a mvd in the pipevac
446  TObjArray* nodes3 = node2->GetNodes();
447  for (Int_t iiiNode = 0; iiiNode < nodes3->GetEntriesFast();
448  iiiNode++) {
449  TGeoNode* node3 = (TGeoNode*) nodes3->At(iiiNode);
450  TString nodeName3 = node3->GetName();
451  nodeName3.ToLower();
452  if (nodeName3.Contains("mvd")) {
453  //fNofMvdStations = node3->GetNodes()->GetEntriesFast();
454  //break;
455  // Fix for mvd_v14a.
456  // Count number of daughter nodes which contain MVDDo
457  TObjArray* nodes4 = node3->GetNodes();
458  for (Int_t iiiiNode = 0; iiiiNode < nodes4->GetEntriesFast();
459  iiiiNode++) {
460  TGeoNode* node4 = (TGeoNode*) nodes4->At(iiiiNode);
461  TString nodeName4 = node4->GetName();
462  nodeName4.ToLower();
463  if (nodeName4.Contains("mvd")) { fNofMvdStations++; }
464  }
465  }
466  }
467  break;
468  }
469  }
470  break;
471  }
472  }
473  }
474  firstTime = false;
475  }
476  cout << "Number of MVD station: " << fNofMvdStations << std::endl;
477  return fNofMvdStations;
478 }
479 
481  static Bool_t firstTime = true;
482  if (firstTime) {
483  fNofStsStations = 0;
484  TObjArray* nodes = gGeoManager->GetTopVolume()->GetNodes();
485  for (Int_t iNode = 0; iNode < nodes->GetEntriesFast(); iNode++) {
486  TGeoNode* node = (TGeoNode*) nodes->At(iNode);
487  if (TString(node->GetName()).Contains("STS_")) { // Top STS node
488  fNofStsStations = node->GetNodes()->GetEntriesFast();
489  break;
490  }
491  }
492  firstTime = false;
493  }
494  return fNofStsStations;
495 }
496 
498  Int_t station,
499  Int_t layer) {
500  static Bool_t firstTime = true;
501  static vector<Int_t> sumOfLayers;
502  if (firstTime) {
503  fNofMuchStations = 0;
504  const TGeoNode* much = static_cast<const TGeoNode*>(
505  fGeo->GetTopNode()->GetNodes()->FindObject("much_0"));
506  if (NULL == much) { // No MUCH detector return 0
507  firstTime = false;
508  return 0;
509  }
510  map<Int_t, Int_t> nofLayersPerStation;
511  TObjArray* muchNodes = much->GetNodes();
512  for (Int_t iMuchNode = 0; iMuchNode < muchNodes->GetEntriesFast();
513  iMuchNode++) {
514  const TGeoNode* muchNode =
515  static_cast<const TGeoNode*>(muchNodes->At(iMuchNode));
516  if (TString(muchNode->GetName()).Contains("station")) {
517  // old structure defined in ASCII geometry
518  if (TString(muchNode->GetName()).Contains("muchstation")) {
519  Int_t station =
520  std::atoi(string(muchNode->GetName() + 11, 2).c_str()) - 1;
521  TObjArray* layerNodes = muchNode->GetNodes();
522  nofLayersPerStation[station] = layerNodes->GetEntriesFast();
523  } else {
524  // New structure defined in ROOT geometry files
525  TObjArray* muchSubNodes = muchNode->GetNodes();
526  for (Int_t iMuchSubNode = 0;
527  iMuchSubNode < muchSubNodes->GetEntriesFast();
528  iMuchSubNode++) {
529  const TGeoNode* muchSubNode =
530  static_cast<const TGeoNode*>(muchSubNodes->At(iMuchSubNode));
531  Int_t station =
532  std::atoi(string(muchSubNode->GetName() + 11, 2).c_str()) - 1;
533  TObjArray* layerNodes = muchSubNode->GetNodes();
534  nofLayersPerStation[station] += layerNodes->GetEntriesFast();
535  }
536  }
537  }
538  }
539  map<Int_t, Int_t>::const_iterator it;
540  Int_t sum = 0;
541  sumOfLayers.push_back(0);
542  for (it = nofLayersPerStation.begin(); it != nofLayersPerStation.end();
543  it++) {
544  sum += (*it).second;
545  sumOfLayers.push_back(sum);
546  }
547  firstTime = false;
548  }
549  return sumOfLayers[station] + layer;
550 }
lit::parallel::LitStation::AddVirtualStation
void AddVirtualStation(const LitVirtualStation< T > &virtualStation)
Add virtual station to detector layout.
Definition: LitStation.h:47
lit::parallel::LitDetectorLayout::AddStation
void AddStation(const LitStation< T > &station)
Add station to detector layout.
Definition: LitDetectorLayout.h:49
CbmLitTrackingGeometryConstructor::ConvertTProfile2DToLitMaterialGrid
void ConvertTProfile2DToLitMaterialGrid(const TProfile2D *profile, lit::parallel::LitMaterialGrid *grid, Double_t maximumValue=0)
Definition: CbmLitTrackingGeometryConstructor.cxx:253
lit::parallel::LitVirtualStation::SetZ
void SetZ(T z)
Definition: LitVirtualStation.h:48
CbmLitTrackingGeometryConstructor::GetMuchLayoutVec
void GetMuchLayoutVec(lit::parallel::LitDetectorLayoutVec &layout)
Return MUCH detector layout for parallel MUCH tracking in SIMD format.
Definition: CbmLitTrackingGeometryConstructor.cxx:60
lit::parallel::LitMaterialGrid
Class stores a grid of material thickness in silicon equivalent.
Definition: LitMaterialGrid.h:37
CbmHistManager::ReadFromFile
void ReadFromFile(TFile *file)
Read histograms from file.
Definition: core/base/CbmHistManager.cxx:110
CbmLitTrackingGeometryConstructor.h
Tracking geometry constructor.
i
int i
Definition: L1/vectors/P4_F32vec4.h:25
CbmLitTrackingGeometryConstructor::GetNofMvdStations
Int_t GetNofMvdStations()
Return number of stations in MVD.
Definition: CbmLitTrackingGeometryConstructor.cxx:419
CbmLitTrackingGeometryConstructor::GetNofMuchStations
Int_t GetNofMuchStations()
Return number of stations in MUCH.
Definition: CbmLitTrackingGeometryConstructor.cxx:337
lit::parallel::LitVirtualStation::SetMaterial
void SetMaterial(const LitMaterialGrid &material)
Definition: LitVirtualStation.h:44
lit::parallel::LitStation
Detector station.
Definition: LitStation.h:31
CbmLitTrackingGeometryConstructor
Definition: CbmLitTrackingGeometryConstructor.h:22
lit::parallel::LitMaterialGrid::SetMaterial
void SetMaterial(const vector< vector< fscal >> &material, fscal xmin, fscal xmax, fscal ymin, fscal ymax, int nofBinsX, int nofBinsY)
Returns Z position of the grid.
Definition: LitMaterialGrid.h:67
lit::parallel::LitVirtualStation
Virtual detector station which stores information needed for track propagation.
Definition: LitVirtualStation.h:31
CbmLitFieldGridCreator
Definition: CbmLitFieldGridCreator.h:17
CbmLitTrackingGeometryConstructor::fNofStsStations
Int_t fNofStsStations
Definition: CbmLitTrackingGeometryConstructor.h:150
CbmHistManager.h
Histogram manager.
CbmLitTrackingGeometryConstructor::GetTrdLayoutScal
void GetTrdLayoutScal(lit::parallel::LitDetectorLayoutScal &layout)
Return TRD detector layout for TRD parallel tracking in scalar format.
Definition: CbmLitTrackingGeometryConstructor.cxx:158
lit::parallel::LitDetectorLayout
Represents detector layout.
Definition: LitDetectorLayout.h:33
min
friend F32vec4 min(const F32vec4 &a, const F32vec4 &b)
Definition: L1/vectors/P4_F32vec4.h:33
CbmLitTrackingGeometryConstructor::GetNofStsStations
Int_t GetNofStsStations()
Return number of stations in STS.
Definition: CbmLitTrackingGeometryConstructor.cxx:480
CbmLitTrackingGeometryConstructor::GetTrdLayout
void GetTrdLayout(lit::parallel::LitDetectorLayout< T > &layout)
Return TRD detector layout for TRD parallel tracking.
Definition: CbmLitTrackingGeometryConstructor.cxx:164
lit::parallel::LitFieldGrid
Class stores a grid of magnetic field values in XY slice at Z position.
Definition: LitFieldGrid.h:44
lit::parallel::LitVirtualStation::SetField
void SetField(const LitFieldGrid &field)
Definition: LitVirtualStation.h:47
CbmHistManager
Histogram manager.
Definition: CbmHistManager.h:41
CbmLitTrackingGeometryConstructor::fNofMuchAbsorbers
Int_t fNofMuchAbsorbers
Definition: CbmLitTrackingGeometryConstructor.h:151
CbmLitTrackingGeometryConstructor::GetNofTrdStations
Int_t GetNofTrdStations()
Return number of stations in TRD.
Definition: CbmLitTrackingGeometryConstructor.cxx:315
CbmLitTrackingGeometryConstructor::GetNofMuchTrdStations
Int_t GetNofMuchTrdStations()
Return number of stations in MUCH + TRD.
Definition: CbmLitTrackingGeometryConstructor.cxx:415
CbmLitTrackingGeometryConstructor::CbmLitTrackingGeometryConstructor
CbmLitTrackingGeometryConstructor()
Constructor. Constructor is protected since singleton pattern is used. Pointer to object is returned ...
Definition: CbmLitTrackingGeometryConstructor.cxx:41
CbmLitDetectorSetup::GetDet
bool GetDet(ECbmModuleId detId) const
Return detector presence in setup.
Definition: CbmLitDetectorSetup.cxx:27
CbmLitTrackingGeometryConstructor::ConvertMuchToAbsoluteStationNr
Int_t ConvertMuchToAbsoluteStationNr(Int_t station, Int_t layer)
Definition: CbmLitTrackingGeometryConstructor.cxx:497
ECbmModuleId::kRich
@ kRich
Ring-Imaging Cherenkov Detector.
CbmLitTrackingGeometryConstructor::fGeo
TGeoManager * fGeo
Definition: CbmLitTrackingGeometryConstructor.h:146
CbmLitTrackingGeometryConstructor::GetMuchLayoutScal
void GetMuchLayoutScal(lit::parallel::LitDetectorLayoutScal &layout)
Return MUCH detector layout for parallel MUCH tracking in scalar format.
Definition: CbmLitTrackingGeometryConstructor.cxx:65
CbmLitFieldGridCreator.h
Class creates grid with magnetic field values at a certain Z position.
CbmLitTrackingGeometryConstructor::GetNofMuchAbsorbers
Int_t GetNofMuchAbsorbers()
Return number of MUCH absorbers.
Definition: CbmLitTrackingGeometryConstructor.cxx:384
CbmUtils.h
CbmLitTrackingGeometryConstructor::fNofMvdStations
Int_t fNofMvdStations
Definition: CbmLitTrackingGeometryConstructor.h:149
CbmLitTrackingGeometryConstructor::GetMuchLayout
void GetMuchLayout(lit::parallel::LitDetectorLayout< T > &layout)
Return MUCH detector layout for parallel MUCH tracking.
Definition: CbmLitTrackingGeometryConstructor.cxx:71
CbmLitTrackingGeometryConstructor::~CbmLitTrackingGeometryConstructor
virtual ~CbmLitTrackingGeometryConstructor()
Destructor.
Definition: CbmLitTrackingGeometryConstructor.cxx:52
CbmLitTrackingGeometryConstructor::fNofTrdStations
Int_t fNofTrdStations
Definition: CbmLitTrackingGeometryConstructor.h:147
CbmHistManager::P2
TProfile2D * P2(const std::string &name) const
Return pointer to TH2 histogram.
Definition: CbmHistManager.h:283
CbmLitTrackingGeometryConstructor::GetTrdLayoutVec
void GetTrdLayoutVec(lit::parallel::LitDetectorLayoutVec &layout)
Return TRD detector layout for TRD parallel tracking in SIMD format.
Definition: CbmLitTrackingGeometryConstructor.cxx:153
lit::parallel::LitDetectorLayout::AddVirtualStation
void AddVirtualStation(const LitVirtualStation< T > &virtualStation)
Add virtual station to detector layout.
Definition: LitDetectorLayout.h:57
max
friend F32vec4 max(const F32vec4 &a, const F32vec4 &b)
Definition: L1/vectors/P4_F32vec4.h:36
CbmLitTrackingGeometryConstructor::GetRichMaterial
void GetRichMaterial(lit::parallel::LitMaterialGrid *material)
Definition: CbmLitTrackingGeometryConstructor.cxx:230
CbmLitTrackingGeometryConstructor::fDet
CbmLitDetectorSetup fDet
Definition: CbmLitTrackingGeometryConstructor.h:152
CbmLitTrackingGeometryConstructor::Instance
static CbmLitTrackingGeometryConstructor * Instance()
Return pointer to singleton object.
Definition: CbmLitTrackingGeometryConstructor.cxx:55
CbmLitTrackingGeometryConstructor::fNofMuchStations
Int_t fNofMuchStations
Definition: CbmLitTrackingGeometryConstructor.h:148
CbmLitFieldGridCreator::CreateGrid
void CreateGrid(fscal Z, lit::parallel::LitFieldGrid &grid)
Main function which creates grid with magnetic field values in (X, Y) slice.
Definition: CbmLitFieldGridCreator.cxx:23
CbmLitDetectorSetup::DetermineSetup
void DetermineSetup()
Determines detector presence using TGeoManager.
Definition: CbmLitDetectorSetup.cxx:79