forked from alisw/AliRoot
-
Notifications
You must be signed in to change notification settings - Fork 0
/
AliGenHalo.cxx
433 lines (378 loc) · 9.47 KB
/
AliGenHalo.cxx
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
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
/**************************************************************************
* Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
* *
* Author: The ALICE Off-line Project. *
* Contributors are mentioned in the code where appropriate. *
* *
* Permission to use, copy, modify and distribute this software and its *
* documentation strictly for non-commercial purposes is hereby granted *
* without fee, provided that the above copyright notice appears in all *
* copies and that both the copyright notice and this permission notice *
* appear in the supporting documentation. The authors make no claims *
* about the suitability of this software for any purpose. It is *
* provided "as is" without express or implied warranty. *
**************************************************************************/
/* $Id$ */
// Read beam halo background particles from a boundary source
// Boundary source is in the LHCb format http://www.hep.manchester.ac.uk/u/robert/LHC_backgrounds/Note-MIBStudies.pdf
// and has been provided by Robert Appleby
// Author: [email protected]
#include <stdlib.h>
#include <TDatabasePDG.h>
#include <TCanvas.h>
#include <TGraph.h>
#include <TPDGCode.h>
#include <TSystem.h>
#include "AliGenHalo.h"
#include "AliGenEventHeader.h"
#include "AliRun.h"
#include "AliLog.h"
ClassImp(AliGenHalo)
AliGenHalo::AliGenHalo()
:AliGenerator(-1),
fFile(0),
fFileName(0),
fSide(1),
fRunPeriod(kY3D90),
fTimePerEvent(1.e-4),
fNskip(0),
fZ1(0),
fZ2(0),
fG1(0),
fG2(0),
fGPASize(0),
fLossID(0),
fLossA(0),
fPdg(0),
fLossT0(0),
fLossZ(0),
fLossW(0),
fXS(0),
fYS(0),
fZS(0),
fDX(0),
fDY(0),
fEkin(0),
fTS(0),
fWS(0)
{
// Constructor
fName = "Halo";
fTitle = "Halo from LHC Beam";
//
// Read all particles
fNpart = -1;
SetAnalog(0);
}
AliGenHalo::AliGenHalo(Int_t npart)
:AliGenerator(npart),
fFile(0),
fFileName(0),
fSide(1),
fRunPeriod(kY3D90),
fTimePerEvent(1.e-4),
fNskip(0),
fZ1(0),
fZ2(0),
fG1(0),
fG2(0),
fGPASize(0),
fLossID(0),
fLossA(0),
fPdg(0),
fLossT0(0),
fLossZ(0),
fLossW(0),
fXS(0),
fYS(0),
fZS(0),
fDX(0),
fDY(0),
fEkin(0),
fTS(0),
fWS(0)
{
// Constructor
fName = "Halo";
fTitle= "Halo from LHC Beam";
//
fNpart = npart;
//
SetAnalog(0);
}
//____________________________________________________________
AliGenHalo::~AliGenHalo()
{
// Destructor
}
//____________________________________________________________
void AliGenHalo::Init()
{
// Initialisation
fFile = fopen(fFileName,"r");
Int_t ir = 0;
if (fFile) {
printf("\n File %s opened for reading, %p ! \n ", fFileName.Data(), (void*)fFile);
} else {
printf("\n Opening of file %s failed, %p ! \n ", fFileName.Data(), (void*)fFile);
return;
}
if (fNskip > 0) {
// Skip the first fNskip events
SkipEvents();
}
//
//
//
// Read file with gas pressure values
char *name = 0;
if (fRunPeriod < 5) {
name = gSystem->ExpandPathName("$(ALICE_ROOT)/LHC/gasPressure.dat" );
fGPASize = 21;
fG1 = new Float_t[fGPASize];
fG2 = new Float_t[fGPASize];
fZ1 = new Float_t[fGPASize];
fZ2 = new Float_t[fGPASize];
} else if (fRunPeriod == 5) {
name = gSystem->ExpandPathName("$(ALICE_ROOT)/LHC/pressure_2003_startup.dat");
fGPASize = 18853;
fG1 = new Float_t[fGPASize];
fZ1 = new Float_t[fGPASize];
} else if (fRunPeriod ==6) {
name = gSystem->ExpandPathName("$(ALICE_ROOT)/LHC/pressure_2003_conditioned.dat");
fGPASize = 12719;
fG1 = new Float_t[fGPASize];
fZ1 = new Float_t[fGPASize];
} else {
Fatal("Init()", "No gas pressure file for given run period !");
}
FILE* file = 0;
if (name) file = fopen(name, "r");
if (!file) {
AliError("No gas pressure file");
return;
}
Float_t z;
Int_t i;
Float_t p[5];
const Float_t kCrossSection = 0.094e-28; // m^2
const Float_t kFlux = 1.e11 / 25.e-9; // protons/s
Float_t pFlux[5] = {0.2, 0.2, 0.3, 0.3, 1.0};
if (fRunPeriod < 5) {
//
// Ring 1
//
for (i = 0; i < fGPASize; i++)
{
ir = fscanf(file, "%f %f %f %f %f %f", &z, &p[0], &p[1], &p[2] , &p[3], &p[4]);
if (ir == 0) break;
fG1[i] = p[fRunPeriod];
if (i > 0) {
fZ1[i] = fZ1[i-1] + z;
} else {
fZ1[i] = 20.;
}
}
//
// Ring 2
//
for (i = 0; i < fGPASize; i++)
{
ir = fscanf(file, "%f %f %f %f %f %f", &z, &p[0], &p[1], &p[2] , &p[3], &p[4]);
if (ir == 0) break;
fG2[i] = p[fRunPeriod];
if (i > 0) {
fZ2[i] = fZ2[i-1] + z;
} else {
fZ2[i] = 20.;
}
}
//
// Interaction rates
//
for (i = 0; i < fGPASize; i++)
{
fG1[i] = fG1[i] * kCrossSection * pFlux[fRunPeriod] * kFlux; // 1/m/s
fG2[i] = fG2[i] * kCrossSection * pFlux[fRunPeriod] * kFlux; // 1/m/s
}
} else {
for (i = 0; i < fGPASize; i++)
{
ir = fscanf(file, "%f %e %e %e %e %e", &z, &p[0], &p[1], &p[2], &p[3], &p[4]);
if (ir == 0) break;
z /= 1000.;
fG1[i] = p[4] * kCrossSection * kFlux; // 1/m/s
// 1/3 of nominal intensity at startup
if (fRunPeriod == kLHCPR674Startup) fG1[i] /= 3.;
fZ1[i] = z;
}
}
//
// Transform into interaction rates
//
Float_t sum1 = 0.;
Float_t sum2 = 0.;
for (Int_t iz = 0; iz < 300; iz++) {
Float_t zpos = 20. + iz * 1.;
zpos *= 100;
Float_t wgt1 = GasPressureWeight( zpos);
Float_t wgt2 = GasPressureWeight(-zpos);
sum1 += wgt1;
sum2 += wgt2;
}
sum1/=250.;
sum2/=250.;
printf("\n %f %f \n \n", sum1, sum2);
delete file;
}
//____________________________________________________________
void AliGenHalo::Generate()
{
// Generate by reading particles from input file
Float_t polar[3]= {0., 0., 0.};
Float_t origin[3];
Float_t p[3], p0;
Float_t tz, txy;
Float_t mass;
//
Int_t nt;
static Bool_t first = kTRUE;
static Int_t oldID = -1;
//
if (first && (fNskip == 0)) ReadNextParticle();
first = kFALSE;
oldID = fLossID;
Int_t np = 0;
while(1) {
// Push particle to stack
mass = TDatabasePDG::Instance()->GetParticle(fPdg)->Mass();
p0 = TMath::Sqrt(fEkin * fEkin + 2. * mass * fEkin);
txy = TMath::Sqrt(fDX * fDX + fDY * fDY);
if (txy > 1.) {
tz = 0.;
} else {
tz = - TMath::Sqrt(1. - txy);
}
p[0] = p0 * fDX;
p[1] = p0 * fDY;
p[2] = p0 * tz;
origin[0] = fXS;
origin[1] = fYS;
origin[2] = 1950.;
PushTrack(fTrackIt , -1, fPdg , p, origin, polar, fTS - 1950./2.9979e10, kPNoProcess, nt, fWS);
np++;
Int_t nc = ReadNextParticle();
if (fLossID != oldID || nc == 0) {
oldID = fLossID;
break;
}
}
SetHighWaterMark(nt);
AliGenEventHeader* header = new AliGenEventHeader("HALO");
header->SetNProduced(np);
// Passes header either to the container or to gAlice
if (fContainer) {
header->SetName(fName);
fContainer->AddHeader(header);
} else {
gAlice->SetGenEventHeader(header);
}
}
Float_t AliGenHalo::GasPressureWeight(Float_t zPrimary)
{
//
// Return z-dependent gasspressure weight = interaction rate [1/m/s].
//
Float_t weight = 0.;
zPrimary /= 100.; // m
if (fRunPeriod < 5) {
Float_t zAbs = TMath::Abs(zPrimary);
if (zPrimary > 0.)
{
if (zAbs > fZ1[20]) {
weight = 2.e4;
} else {
for (Int_t i = 1; i < 21; i++) {
if (zAbs < fZ1[i]) {
weight = fG1[i];
break;
}
}
}
} else {
if (zAbs > fZ2[20]) {
weight = 2.e4;
} else {
for (Int_t i = 1; i < 21; i++) {
if (zAbs < fZ2[i]) {
weight = fG2[i];
break;
}
}
}
}
} else {
Int_t index = TMath::BinarySearch(fGPASize, fZ1, zPrimary);
weight = fG1[index];
}
return weight;
}
void AliGenHalo::Draw(Option_t *)
{
// Draws the gas pressure distribution
Float_t z[400];
Float_t p[400];
for (Int_t i = 0; i < 400; i++)
{
z[i] = -20000. + Float_t(i) * 100;
p[i] = GasPressureWeight(z[i]);
}
TGraph* gr = new TGraph(400, z, p);
TCanvas* c1 = new TCanvas("c1","Canvas 1",400,10,600,700);
c1->cd();
gr->Draw("AL");
}
Int_t AliGenHalo::ReadNextParticle()
{
// Read the next particle from the file
Int_t ncols = fscanf(fFile,"%d %f %f %d %f %d %f %f %f %f %f %f %f %f",
&fLossID, &fLossT0, &fLossZ, &fLossA, &fLossW, &fPdg, &fXS, &fYS, &fZS, &fDX, &fDY, &fEkin, &fTS, &fWS);
fLossZ /= 10.;
fXS /= 10.;
fYS /= 10.;
fZS /= 10.;
fTS *= 1.e-9;
return (ncols);
}
void AliGenHalo::SkipEvents()
{
//
// Skip the first fNskip events
Int_t skip = fNskip;
Int_t oldID = -1;
while (skip >= 0)
{
ReadNextParticle();
if (oldID != fLossID) {
oldID = fLossID;
skip--;
}
}
}
void AliGenHalo::CountEvents()
{
// Count total number of events
Int_t nev = 0;
Int_t oldID = -1;
Int_t nc = 1;
while (nc != -1)
{
nc = ReadNextParticle();
if (oldID != fLossID) {
oldID = fLossID;
nev++;
printf("Number of events %10d %10d \n", nev, oldID);
}
}
rewind(fFile);
}