ATLAS Offline Software
Loading...
Searching...
No Matches
LayerMaterialInspector.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2024 CERN for the benefit of the ATLAS collaboration
3*/
4
6// LayerMaterialInspector.cxx, (c) ATLAS Detector software
8
14#include "TrkGeometry/Layer.h"
15#include "TrkSurfaces/Surface.h"
19// Root
20#include "TH2F.h"
21#include "TTree.h"
22#include "TString.h"
23// Gaudi
24#include "GaudiKernel/ITHistSvc.h"
25
26// constructor
27Trk::LayerMaterialInspector::LayerMaterialInspector(const std::string& t, const std::string& n, const IInterface* p) :
29 m_treeFolder("/val/")
30{}
31
32// destructor
34= default;
35
37{
38 return StatusCode::SUCCESS;
39}
40
41StatusCode Trk::LayerMaterialInspector::processNode(const Trk::Layer& lay, size_t) const
42{
43
45 TString folderId = tvol ? tvol->volumeName() : std::string("Unknown");
46 TString folderName = m_treeFolder+folderId;
47 folderName.ReplaceAll("::","_");
48
49 // prepare the information
50 int lIndex = lay.layerIndex().value();
51 const Amg::Vector3D& lCenter = lay.surfaceRepresentation().center();
53
54 // skip navigation layers
55 if (!lMaterial) return StatusCode::SUCCESS;
56
57 int lType = 0; // 1 - cylinder, 2 - disc
58 float dim0 = 0.;
59 float dim1 = 0.;
60
61 int bins0 = 1;
62 int bins1 = 1;
63
64 float centerX = lCenter.x();
65 float centerY = lCenter.y();
66 float centerZ = lCenter.z();
67
68 // cylinder bounds
70 lType = 1;
71 // cylinder bounds
72 const Trk::CylinderBounds* cb = dynamic_cast<const Trk::CylinderBounds*>(&(lay.surfaceRepresentation().bounds()));
73 if (cb){
74 dim0 = cb->r();
75 dim1 = cb->halflengthZ();
76 }
77 } else if ( lay.surfaceRepresentation().type() == Trk::SurfaceType::Disc ) {
78 lType = 2;
79 // disc bounds
80 const Trk::DiscBounds* db = dynamic_cast<const Trk::DiscBounds*>(&(lay.surfaceRepresentation().bounds()));
81 if (db){
82 dim0 = db->rMin();
83 dim1 = db->rMax();
84 }
85 }
86
87 if (lMaterial){
88 const Trk::BinUtility* lBinUtility = lMaterial->binUtility();
89 bins0 = lBinUtility ? lBinUtility->max(0)+1 : 1;
90 bins1 = lBinUtility ? lBinUtility->max(1)+1 : 1;
91 }
92
93
94 // -------------------------------------------------------------------------------------
95 TString hName = ( lType == 1 ) ? "CylinderLayer_" : "DiscLayer_";
96 hName += int(lIndex);
97
98 TString pXo = "_pX0";
99 TString info = "_Information";
100
101 TH2F* lMaterialHist = lType == 1 ?
102 new TH2F(hName+pXo, hName, bins0, -M_PI*dim0, M_PI*dim0, bins1, -dim1, dim1) :
103 new TH2F(hName+pXo, hName, bins0, dim0, dim1, bins1, -M_PI, M_PI);
104
105 TTree* lTreeInformation = new TTree(hName+info,hName);
106 lTreeInformation->Branch("LayerCenterX", &centerX, "cX/F");
107 lTreeInformation->Branch("LayerCenterY", &centerY, "cY/F");
108 lTreeInformation->Branch("LayerCenterZ", &centerZ, "cZ/F");
109 lTreeInformation->Branch("LayerIndex", &lIndex, "idx/I");
110 if (lType == 1){
111 lTreeInformation->Branch("Radius", &dim0, "r/F");
112 lTreeInformation->Branch("HalflengthZ", &dim1, "hz/F");
113 lTreeInformation->Branch("BinsRPhi", &bins0, "bin0/I");
114 lTreeInformation->Branch("BinsHalfZ", &bins1, "bin1/I");
115 } else if (lType == 2){
116 lTreeInformation->Branch("InnerRadius", &dim0, "rmin/F");
117 lTreeInformation->Branch("OuterRadius", &dim1, "rmax/F");
118 lTreeInformation->Branch("BinsR", &bins0, "bin0/I");
119 lTreeInformation->Branch("BinsPhi", &bins1, "bin1/I");
120 }
121
122 // --------------- registration for histogram and tree ---------------------------------
123 TString regHistName = folderName+"/"+hName+"/"+hName+pXo;
124 TString regTreeName = folderName+"/"+hName+"/"+hName+info;
125
126 SmartIF<ITHistSvc> tHistSvc{service("THistSvc")};
127 ATH_CHECK( tHistSvc.isValid() );
128
129 if ( (tHistSvc->regTree(std::string(regTreeName.Data()),lTreeInformation)).isFailure() ){
130 ATH_MSG_WARNING("Could not register TTree with name " << hName);
131 }
132 lTreeInformation->Fill();
133
134 if ( (tHistSvc->regHist(std::string(regHistName.Data()),lMaterialHist)).isFailure() ){
135 ATH_MSG_WARNING("Could not register THist with name " << hName);
136 }
137
138 if (lMaterial){
139 // fill the material into the bins
140 for ( int ib0 = 0; ib0 < bins0; ++ib0 )
141 for ( int ib1 = 0; ib1 < bins1; ++ib1 ){
142 // get the material
143 const Trk::MaterialProperties* mps = lMaterial->material(size_t(ib0),size_t(ib1));
144 if (mps)
145 lMaterialHist->SetBinContent(ib0+1,ib1+1,mps->thicknessInX0());
146 }
147 }
148
149 return StatusCode::SUCCESS;
150}
151
153{
154 return StatusCode::SUCCESS;
155}
156
157
158
159
#define M_PI
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_WARNING(x)
A generic symmetric BinUtility, for fully symmetric binning in terms of binning grid and binning type...
Definition BinUtility.h:39
size_t max(size_t ba=0) const
First bin maximal value.
Definition BinUtility.h:212
Bounds for a cylindrical Surface.
virtual double r() const override final
This method returns the radius.
double halflengthZ() const
This method returns the halflengthZ.
Class to describe the bounds for a planar DiscSurface.
Definition DiscBounds.h:44
int value() const
layerIndex expressed in an integer
Definition LayerIndex.h:71
virtual ~LayerMaterialInspector()
Destructor.
virtual StatusCode processNode(const TrackingVolume &tvol, size_t level=0) const
Processor Action to work on TrackingVolumes - the level is for the hierachy tree.
LayerMaterialInspector(const std::string &, const std::string &, const IInterface *)
Constructor.
This virtual base class encapsulates the logics to build pre/post/full update material for Layer stru...
virtual const BinUtility * binUtility() const =0
Return the BinUtility.
virtual const MaterialProperties * material(size_t ib0, size_t ib1) const =0
Direct access via bins to the MaterialProperties.
Base Class for a Detector Layer in the Tracking realm.
Definition Layer.h:72
virtual const Surface & surfaceRepresentation() const =0
Transforms the layer into a Surface representation for extrapolation.
const LayerMaterialProperties * layerMaterialProperties() const
getting the LayerMaterialProperties including full/pre/post update
const LayerIndex & layerIndex() const
get the layerIndex
const TrackingVolume * enclosingTrackingVolume() const
get the confining TrackingVolume
Material with information about thickness of material.
float thicknessInX0() const
Return the radiationlength fraction.
RecursiveGeometryProcessor(const std::string &, const std::string &, const IInterface *)
Constructor.
Abstract Base Class for tracking surfaces.
virtual constexpr SurfaceType type() const =0
Returns the Surface type to avoid dynamic casts.
const Amg::Vector3D & center() const
Returns the center position of the Surface.
Full Volume description used in Tracking, it inherits from Volume to get the geometrical structure,...
const std::string & volumeName() const
Returns the VolumeName - for debug reason, might be depreciated later.
Eigen::Matrix< double, 3, 1 > Vector3D
Ensure that the ATLAS eigen extensions are properly loaded.