"README.md" did not exist on "31debfd2c951ef1fec4d7cecb6512855a625c1f0"
Newer
Older
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
//==========================================================================
// AIDA Detector description implementation
//--------------------------------------------------------------------------
// Copyright (C) Organisation europeenne pour la Recherche nucleaire (CERN)
// All rights reserved.
//
// For the licensing terms see $DD4hepINSTALL/LICENSE.
// For the list of contributors see $DD4hepINSTALL/doc/CREDITS.
//
//==========================================================================
//
// Compute the material budget in terms of integrated radiation and ineraction
// lengths as seen from the IP for various subdetector regions given in
// a simple steering file. Creates a root file w/ histograms and dumps numbers
// the screen.
//
// Author : F.Gaede, DESY
//
//==========================================================================
#include "TError.h"
// Framework include files
#include "DD4hep/Detector.h"
#include "DD4hep/DetType.h"
#include "DD4hep/Printout.h"
#include "DDRec/MaterialManager.h"
// #include "TGeoVolume.h"
// #include "TGeoManager.h"
// #include "TGeoNode.h"
#include "TFile.h"
#include "TH1F.h"
#include <fstream>
#include "main.h"
using namespace dd4hep;
using namespace dd4hep::rec;
void dumpExampleSteering() ;
namespace {
struct SDetHelper{
std::string name ;
double r0 ;
double r1 ;
double z0 ;
double z1 ;
TH1F* hx ;
TH1F* hl ;
};
/// comput a point on a cylinder (r,z) for a given direction (theta,phi) from the IP
Vector3D pointOnCylinder( double theta, double r, double z, double phi){
double theta0 = atan2( r, z ) ;
Vector3D v = ( theta > theta0 ?
Vector3D( r , phi , r / tan( theta ) , Vector3D::cylindrical ) :
Vector3D( z * tan( theta ) , phi , z , Vector3D::cylindrical ) ) ;
return v ;
}
}
int main_wrapper(int argc, char** argv) {
struct Handler {
Handler() { SetErrorHandler(Handler::print); }
static void print(int level, Bool_t abort, const char *location, const char *msg) {
if ( level > kInfo || abort ) ::printf("%s: %s\n", location, msg);
}
static void usage() {
std::cout << " usage: materialBudget compact.xml steering.txt" << std::endl
<< " -> create histograms with the material budget as seen from the IP within fixed ranges of (rmin, rmax, zmin. zmax)" << std::endl
<< " see example steering file for details ..." << std::endl
<< " -x : dump example steering file " << std::endl
<< std::endl;
exit(1);
}
} _handler;
if( argc == 2 ){
std::string dummy = argv[1] ;
if( dummy == "-x" )
dumpExampleSteering() ;
else
Handler::usage();
}
if( argc != 3 )
Handler::usage();
// ================== parse steering file ========================================================
std::string compactFile = argv[1];
std::string steeringFile = argv[2];
int nbins = 90 ;
double phi0 = M_PI / 2. ;
std::string outFileName("material_budget.root") ;
std::vector<SDetHelper> subdets ;
std::ifstream infile(steeringFile );
std::string line;
while (std::getline(infile, line))
{
// ignore empty lines and comments
if( line.empty() || line.find("#") == 0 )
continue ;
std::istringstream iss(line);
std::string token ;
iss >> token ;
if( token == "nbins" ){
iss >> nbins ;
}
else if( token == "rootfile" ){
iss >> outFileName ;
}
else if( token == "phi" ){
iss >> phi0 ;
phi0 = phi0 / 180. * M_PI ;
}
else if( token == "subdet" ){
SDetHelper det ;
iss >> det.name >> det.r0 >> det.z0 >> det.r1 >> det.z1 ;
subdets.push_back( det );
}
if ( !iss.eof() || iss.fail() ){
std::cout << " ERROR parsing line : " << line << std::endl ;
exit(1) ;
}
}
// =================================================================================================
setPrintLevel(WARNING);
Detector& description = Detector::getInstance();
description.fromXML(compactFile ) ;
//----- open root file and book histograms
TFile* rootFile = new TFile(outFileName.c_str(),"recreate");
for(auto& det : subdets){
std::string hxn(det.name), hxnn(det.name) ;
std::string hln(det.name), hlnn(det.name) ;
hxn += "x0" ;
hxnn += " integrated X0 vs -theta" ;
det.hx = new TH1F( hxn.c_str(), hxnn.c_str(), nbins, -90. , 0. ) ;
hln += "lambda" ;
hlnn += " integrated int. lengths vs -theta" ;
det.hl = new TH1F( hln.c_str(), hlnn.c_str(), nbins, -90. , 0. ) ;
}
//-------------------------
Volume world = description.world().volume() ;
MaterialManager matMgr( world ) ;
double dTheta = 0.5*M_PI/nbins; // bin size
std::cout << "====================================================================================================" << std::endl ;
std::cout << "theta:f/" ;
for(auto& det : subdets){ std::cout << det.name << "_x0:f/" << det.name << "_lam:f/" ; }
std::cout << std::endl ;
for(int i=0 ; i< nbins ;++i){
double theta = (0.5+i)*dTheta ;
std::cout << std::scientific << theta << " " ;
for(auto& det : subdets){
Vector3D p0 = pointOnCylinder( theta, det.r0 , det.z0 , phi0 ) ;// double theta, double r, double z, double phi)
Vector3D p1 = pointOnCylinder( theta, det.r1 , det.z1 , phi0 ) ;// double theta, double r, double z, double phi)
const MaterialVec& materials = matMgr.materialsBetween(p0, p1);
double sum_x0(0.), sum_lambda(0.),path_length(0.);
for( auto amat : materials ){
TGeoMaterial* mat = amat.first->GetMaterial();
double length = amat.second;
double nx0 = length / mat->GetRadLen();
sum_x0 += nx0;
double nLambda = length / mat->GetIntLen();
sum_lambda += nLambda;
path_length += length;
}
det.hx->Fill( -theta/M_PI*180. , sum_x0 ) ;
det.hl->Fill( -theta/M_PI*180. , sum_lambda ) ;
std::cout << std::scientific << sum_x0 << " " << sum_lambda << " " ; // << path_length ;
}
std::cout << std::endl ;
}
std::cout << "====================================================================================================" << std::endl ;
rootFile->Write();
rootFile->Close();
return 0;
}
void dumpExampleSteering(){
std::cout << "# Example steering file for materialBudget (taken from ILD_l5_v02)" << std::endl ;
std::cout << std::endl ;
std::cout << "# root output file" << std::endl ;
std::cout << "rootfile material_budget.root" << std::endl ;
std::cout << "# number of bins for polar angle (default 90)" << std::endl ;
std::cout << "nbins 90" << std::endl ;
std::cout << std::endl ;
std::cout << "# phi direction in deg (default: 90./y-axis)" << std::endl ;
std::cout << "phi 90." << std::endl ;
std::cout << std::endl ;
std::cout << "# names and subdetector ranges given in [rmin,zmin,rmax,zmax] - e.g. for ILD_l5_vo2 (run dumpdetector -d to get numbers... ) " << std::endl ;
std::cout << std::endl ;
std::cout << "subdet vxd 0. 0. 6.549392e+00 1.450000e+01" << std::endl ;
std::cout << "subdet sitftd 0. 0. 3.024000e+01 2.211800e+02" << std::endl ;
std::cout << "subdet tpc 0. 0. 1.692100e+02 2.225000e+02" << std::endl ;
std::cout << "subdet outtpc 0. 0. 1.769800e+02 2.350000e+02" << std::endl ;
std::cout << "subdet set 0. 0. 1.775200e+02 2.300000e+02" << std::endl ;
exit(0);
}