]> git.uio.no Git - u/mrichter/AliRoot.git/blame - TRD/ConfigPID.C
new schema for config
[u/mrichter/AliRoot.git] / TRD / ConfigPID.C
CommitLineData
defc9040 1// Macro to generate e, pi, mu, K, p 200 each with box generator
2// One can use the configuration macro in compiled mode by
3// root [0] gSystem->Load("libgeant321");
4// root [0] gSystem->SetIncludePath("-I$ROOTSYS/include -I$ALICE_ROOT/include\
5// -I$ALICE_ROOT -I$ALICE/geant3/TGeant3");
6// root [0] .x grun.C(1,"ConfigPPR.C++")
7// Prashant Shukla (shukla@physi.uni-heidelberg.de)
8
9#if !defined(__CINT__) || defined(__MAKECINT__)
10#include <Riostream.h>
11#include <TRandom.h>
12#include <TSystem.h>
13#include <TVirtualMC.h>
14#include <TGeant3TGeo.h>
15#include <TPDGCode.h>
16#include <TF1.h>
17#include "STEER/AliRunLoader.h"
18#include "STEER/AliRun.h"
19#include "STEER/AliConfig.h"
20#include "STEER/AliGenerator.h"
21#include "STEER/AliLog.h"
22#include "PYTHIA6/AliDecayerPythia.h"
23#include "EVGEN/AliGenHIJINGpara.h"
24#include "THijing/AliGenHijing.h"
25#include "EVGEN/AliGenCocktail.h"
26#include "EVGEN/AliGenSlowNucleons.h"
27#include "EVGEN/AliSlowNucleonModelExp.h"
28#include "EVGEN/AliGenParam.h"
29#include "EVGEN/AliGenMUONlib.h"
30#include "EVGEN/AliGenSTRANGElib.h"
31#include "EVGEN/AliGenMUONCocktail.h"
32#include "EVGEN/AliGenCocktail.h"
33#include "EVGEN/AliGenGeVSim.h"
34#include "EVGEN/AliGeVSimParticle.h"
35#include "PYTHIA6/AliGenPythia.h"
36#include "STEER/AliMagFMaps.h"
37#include "STRUCT/AliBODY.h"
38#include "STRUCT/AliMAG.h"
39#include "STRUCT/AliABSOv0.h"
40#include "STRUCT/AliDIPOv2.h"
41#include "STRUCT/AliHALL.h"
42#include "STRUCT/AliFRAMEv2.h"
43#include "STRUCT/AliSHILv2.h"
44#include "STRUCT/AliPIPEv0.h"
45#include "ITS/AliITSvPPRasymmFMD.h"
46#include "TPC/AliTPCv2.h"
47#include "TOF/AliTOFv4T0.h"
48#include "RICH/AliRICHv1.h"
49#include "ZDC/AliZDCv2.h"
50#include "TRD/AliTRDv1.h"
51#include "FMD/AliFMDv1.h"
52#include "MUON/AliMUONv1.h"
53#include "MUON/AliMUONSt1GeometryBuilderV2.h"
54#include "MUON/AliMUONSt2GeometryBuilderV2.h"
55#include "MUON/AliMUONSlatGeometryBuilder.h"
56#include "MUON/AliMUONTriggerGeometryBuilder.h"
57#include "PHOS/AliPHOSv1.h"
58#include "PMD/AliPMDv1.h"
59#include "START/AliSTARTv1.h"
60#include "EMCAL/AliEMCALv1.h"
61#include "CRT/AliCRTv0.h"
62#include "VZERO/AliVZEROv5.h"
63#endif
64
65enum PprRun_t
66{
67 test50,
68 kParam_8000, kParam_4000, kParam_2000,
69 kHijing_cent1, kHijing_cent2,
70 kHijing_per1, kHijing_per2, kHijing_per3, kHijing_per4, kHijing_per5,
71 kCocktailTRD, kpieTRD
72};
73
74const char* pprRunName[] = {
75 "test50",
76 "kParam_8000", "kParam_4000", "kParam_2000",
77 "kHijing_cent1", "kHijing_cent2",
78 "kHijing_per1", "kHijing_per2", "kHijing_per3", "kHijing_per4",
79 "kHijing_per5",
80 "kCocktailTRD", "kpieTRD"
81};
82
83enum PprGeo_t
84{
85 kHoles, kNoHoles
86};
87
88enum PprRad_t
89{
90 kGluonRadiation, kNoGluonRadiation
91};
92
93enum PprMag_t
94{
95 k2kG, k4kG, k5kG
96};
97
98
99// This part for configuration
100//static PprRun_t srun = test50;
101static PprRun_t srun = kpieTRD;
102static PprGeo_t sgeo = kNoHoles;
103static PprRad_t srad = kGluonRadiation;
104static PprMag_t smag = k5kG;
105static Int_t sseed = 0; //Set 0 to use the current time
106
107// Comment line
108static TString comment;
109
110// Functions
111Float_t EtaToTheta(Float_t arg);
112AliGenerator* GeneratorFactory(PprRun_t srun);
113AliGenHijing* HijingStandard();
114void ProcessEnvironmentVars();
115
116void Config()
117{
118 // ThetaRange is (0., 180.). It was (0.28,179.72) 7/12/00 09:00
119 // Theta range given through pseudorapidity limits 22/6/2001
120
121 // Get settings from environment variables
122 ProcessEnvironmentVars();
123
124 // Set Random Number seed
125 gRandom->SetSeed(sseed);
126 cout<<"Seed for random number generation= "<<gRandom->GetSeed()<<endl;
127
128
129 // libraries required by geant321
130#if defined(__CINT__)
131 gSystem->Load("libgeant321");
132#endif
133
134 new TGeant3TGeo("C++ Interface to Geant3");
135
136 AliRunLoader* rl=0x0;
137
138 cout<<"Config.C: Creating Run Loader ..."<<endl;
139 rl = AliRunLoader::Open("galice.root",
140 AliConfig::GetDefaultEventFolderName(),
141 "recreate");
142 if (rl == 0x0)
143 {
144 gAlice->Fatal("Config.C","Can not instatiate the Run Loader");
145 return;
146 }
147 rl->SetCompressionLevel(2);
148 rl->SetNumberOfEventsPerFile(3);
149 gAlice->SetRunLoader(rl);
150
151 //
152 // Set External decayer
153 AliDecayer *decayer = new AliDecayerPythia();
154
155 decayer->SetForceDecay(kAll);
156
157 decayer->Init();
158 gMC->SetExternalDecayer(decayer);
159 //
160 //
161 //=======================================================================
162 //
163 //=======================================================================
164 // ************* STEERING parameters FOR ALICE SIMULATION **************
165 // --- Specify event type to be tracked through the ALICE setup
166 // --- All positions are in cm, angles in degrees, and P and E in GeV
167
168 gMC->SetProcess("DCAY",1);
169 gMC->SetProcess("PAIR",1);
170 gMC->SetProcess("COMP",1);
171 gMC->SetProcess("PHOT",1);
172 gMC->SetProcess("PFIS",0);
173 gMC->SetProcess("DRAY",0);
174 gMC->SetProcess("ANNI",1);
175 gMC->SetProcess("BREM",1);
176 gMC->SetProcess("MUNU",1);
177 gMC->SetProcess("CKOV",1);
178 gMC->SetProcess("HADR",1);
179 gMC->SetProcess("LOSS",2);
180 gMC->SetProcess("MULS",1);
181 gMC->SetProcess("RAYL",1);
182
183 Float_t cut = 1.e-3; // 1MeV cut by default
184 Float_t tofmax = 1.e10;
185
186 gMC->SetCut("CUTGAM", cut);
187 gMC->SetCut("CUTELE", cut);
188 gMC->SetCut("CUTNEU", cut);
189 gMC->SetCut("CUTHAD", cut);
190 gMC->SetCut("CUTMUO", cut);
191 gMC->SetCut("BCUTE", cut);
192 gMC->SetCut("BCUTM", cut);
193 gMC->SetCut("DCUTE", cut);
194 gMC->SetCut("DCUTM", cut);
195 gMC->SetCut("PPCUTM", cut);
196 gMC->SetCut("TOFMAX", tofmax);
197
198 // Debug and log level
199 // AliLog::SetGlobalDebugLevel(0);
200 // AliLog::SetGlobalLogLevel(AliLog::kError);
201
202 // Generator Configuration
203 AliGenerator* gener = GeneratorFactory(srun);
204 gener->SetOrigin(0, 0, 0); // vertex position
205 gener->SetSigma(0, 0, 5.3); // Sigma in (X,Y,Z) (cm) on IP position
206 gener->SetCutVertexZ(1.); // Truncate at 1 sigma
207 gener->SetVertexSmear(kPerEvent);
208 gener->SetTrackingFlag(1);
209 gener->Init();
210
211 if (smag == k2kG) {
212 comment = comment.Append(" | L3 field 0.2 T");
213 } else if (smag == k4kG) {
214 comment = comment.Append(" | L3 field 0.4 T");
215 } else if (smag == k5kG) {
216 comment = comment.Append(" | L3 field 0.5 T");
217 }
218
219
220 if (srad == kGluonRadiation)
221 {
222 comment = comment.Append(" | Gluon Radiation On");
223
224 } else {
225 comment = comment.Append(" | Gluon Radiation Off");
226 }
227
228 if (sgeo == kHoles)
229 {
230 comment = comment.Append(" | Holes for PHOS/RICH");
231
232 } else {
233 comment = comment.Append(" | No holes for PHOS/RICH");
234 }
235
236 printf("\n \n Comment: %s \n \n", comment.Data());
237
238
239// Field (L3 0.4 T)
240 AliMagFMaps* field = new AliMagFMaps("Maps","Maps", 2, 1., 10., smag);
241 // AliMagFMaps* field = new AliMagFMaps("Maps","Maps", 2, 0, 10., smag); // B = 0
242 field->SetL3ConstField(0); //Using const. field in the barrel
243 rl->CdGAFile();
244 gAlice->SetField(field);
245//
246 Int_t iABSO = 1;
247 Int_t iDIPO = 1;
248 Int_t iFMD = 1;
249 Int_t iFRAME = 1;
250 Int_t iHALL = 1;
251 Int_t iITS = 1;
252 Int_t iMAG = 1;
253 Int_t iMUON = 0;
254 Int_t iPHOS = 1;
255 Int_t iPIPE = 1;
256 Int_t iPMD = 0;
257 Int_t iRICH = 0;
258 Int_t iSHIL = 1;
259 Int_t iSTART = 1;
260 Int_t iTOF = 0;
261 Int_t iTPC = 1;
262 Int_t iTRD = 1;
263 Int_t iZDC = 1;
264 Int_t iEMCAL = 0;
265 Int_t iVZERO = 1;
266 Int_t iCRT = 0;
267
268 //=================== Alice BODY parameters =============================
269 AliBODY *BODY = new AliBODY("BODY", "Alice envelop");
270
271
272 if (iMAG)
273 {
274 //=================== MAG parameters ============================
275 // --- Start with Magnet since detector layouts may be depending ---
276 // --- on the selected Magnet dimensions ---
277 AliMAG *MAG = new AliMAG("MAG", "Magnet");
278 }
279
280
281 if (iABSO)
282 {
283 //=================== ABSO parameters ============================
284 AliABSO *ABSO = new AliABSOv0("ABSO", "Muon Absorber");
285 }
286
287 if (iDIPO)
288 {
289 //=================== DIPO parameters ============================
290
291 AliDIPO *DIPO = new AliDIPOv2("DIPO", "Dipole version 2");
292 }
293
294 if (iHALL)
295 {
296 //=================== HALL parameters ============================
297
298 AliHALL *HALL = new AliHALL("HALL", "Alice Hall");
299 }
300
301
302 if (iFRAME)
303 {
304 //=================== FRAME parameters ============================
305
306 AliFRAMEv2 *FRAME = new AliFRAMEv2("FRAME", "Space Frame");
307 if (sgeo == kHoles) {
308 FRAME->SetHoles(1);
309 } else {
310 FRAME->SetHoles(0);
311 }
312 }
313
314 if (iSHIL)
315 {
316 //=================== SHIL parameters ============================
317
318 AliSHIL *SHIL = new AliSHILv2("SHIL", "Shielding Version 2");
319 }
320
321
322 if (iPIPE)
323 {
324 //=================== PIPE parameters ============================
325
326 AliPIPE *PIPE = new AliPIPEv0("PIPE", "Beam Pipe");
327 }
328
329 if(iITS) {
330
331 //=================== ITS parameters ============================
332 //
333 // As the innermost detector in ALICE, the Inner Tracking System "impacts" on
334 // almost all other detectors. This involves the fact that the ITS geometry
335 // still has several options to be followed in parallel in order to determine
336 // the best set-up which minimizes the induced background. All the geometries
337 // available to date are described in the following. Read carefully the comments
338 // and use the default version (the only one uncommented) unless you are making
339 // comparisons and you know what you are doing. In this case just uncomment the
340 // ITS geometry you want to use and run Aliroot.
341 //
342 // Detailed geometries:
343 //
344 //
345 //AliITS *ITS = new AliITSv5symm("ITS","Updated ITS TDR detailed version with symmetric services");
346 //
347 //AliITS *ITS = new AliITSv5asymm("ITS","Updates ITS TDR detailed version with asymmetric services");
348 //
349 AliITSvPPRasymmFMD *ITS = new AliITSvPPRasymmFMD("ITS","New ITS PPR detailed version with asymmetric services");
350 ITS->SetMinorVersion(2); // don't touch this parameter if you're not an ITS developer
351 ITS->SetReadDet(kTRUE); // don't touch this parameter if you're not an ITS developer
352 // ITS->SetWriteDet("$ALICE_ROOT/ITS/ITSgeometry_vPPRasymm2.det"); // don't touch this parameter if you're not an ITS developer
353 ITS->SetThicknessDet1(200.); // detector thickness on layer 1 must be in the range [100,300]
354 ITS->SetThicknessDet2(200.); // detector thickness on layer 2 must be in the range [100,300]
355 ITS->SetThicknessChip1(200.); // chip thickness on layer 1 must be in the range [150,300]
356 ITS->SetThicknessChip2(200.); // chip thickness on layer 2 must be in the range [150,300]
357 ITS->SetRails(0); // 1 --> rails in ; 0 --> rails out
358 ITS->SetCoolingFluid(1); // 1 --> water ; 0 --> freon
359
360 // Coarse geometries (warning: no hits are produced with these coarse geometries and they unuseful
361 // for reconstruction !):
362 //
363 //
364 //AliITSvPPRcoarseasymm *ITS = new AliITSvPPRcoarseasymm("ITS","New ITS PPR coarse version with asymmetric services");
365 //ITS->SetRails(0); // 1 --> rails in ; 0 --> rails out
366 //ITS->SetSupportMaterial(0); // 0 --> Copper ; 1 --> Aluminum ; 2 --> Carbon
367 //
368 //AliITS *ITS = new AliITSvPPRcoarsesymm("ITS","New ITS PPR coarse version with symmetric services");
369 //ITS->SetRails(0); // 1 --> rails in ; 0 --> rails out
370 //ITS->SetSupportMaterial(0); // 0 --> Copper ; 1 --> Aluminum ; 2 --> Carbon
371 //
372 //
373 //
374 // Geant3 <-> EUCLID conversion
375 // ============================
376 //
377 // SetEUCLID is a flag to output (=1) or not to output (=0) both geometry and
378 // media to two ASCII files (called by default ITSgeometry.euc and
379 // ITSgeometry.tme) in a format understandable to the CAD system EUCLID.
380 // The default (=0) means that you dont want to use this facility.
381 //
382 ITS->SetEUCLID(0);
383 }
384
385 if (iTPC)
386 {
387 //============================ TPC parameters ================================
388 // --- This allows the user to specify sectors for the SLOW (TPC geometry 2)
389 // --- Simulator. SecAL (SecAU) <0 means that ALL lower (upper)
390 // --- sectors are specified, any value other than that requires at least one
391 // --- sector (lower or upper)to be specified!
392 // --- Reminder: sectors 1-24 are lower sectors (1-12 -> z>0, 13-24 -> z<0)
393 // --- sectors 25-72 are the upper ones (25-48 -> z>0, 49-72 -> z<0)
394 // --- SecLows - number of lower sectors specified (up to 6)
395 // --- SecUps - number of upper sectors specified (up to 12)
396 // --- Sens - sensitive strips for the Slow Simulator !!!
397 // --- This does NOT work if all S or L-sectors are specified, i.e.
398 // --- if SecAL or SecAU < 0
399 //
400 //
401 //-----------------------------------------------------------------------------
402
403 // gROOT->LoadMacro("SetTPCParam.C");
404 // AliTPCParam *param = SetTPCParam();
405 AliTPC *TPC = new AliTPCv2("TPC", "Default");
406
407 // All sectors included
408 TPC->SetSecAL(-1);
409 TPC->SetSecAU(-1);
410
411 }
412
413
414 if (iTOF) {
415 //=================== TOF parameters ============================
416 AliTOF *TOF = new AliTOFv4T0("TOF", "normal TOF");
417 }
418
419
420 if (iRICH)
421 {
422 //=================== RICH parameters ===========================
423 AliRICH *RICH = new AliRICHv1("RICH", "normal RICH");
424
425 }
426
427
428 if (iZDC)
429 {
430 //=================== ZDC parameters ============================
431
432 AliZDC *ZDC = new AliZDCv2("ZDC", "normal ZDC");
433 }
434
435 if (iTRD)
436 {
437 //=================== TRD parameters ============================
438
439 AliTRD *TRD = new AliTRDv1("TRD", "TRD slow simulator");
440
441 // Select the gas mixture (0: 97% Xe + 3% isobutane, 1: 90% Xe + 10% CO2)
442 TRD->SetGasMix(1);
443 if (sgeo == kHoles) {
444 // With hole in front of PHOS
445 TRD->SetPHOShole();
446 // With hole in front of RICH
447 TRD->SetRICHhole();
448 }
449 // Switch on TR
450 AliTRDsim *TRDsim = TRD->CreateTR();
451 }
452
453 if (iFMD)
454 {
455 //=================== FMD parameters ============================
456 AliFMD *FMD = new AliFMDv1("FMD", "normal FMD");
457 }
458
459 if (iMUON)
460 {
461 //=================== MUON parameters ===========================
462 // New MUONv1 version (geometry defined via builders)
463 AliMUON *MUON = new AliMUONv1("MUON", "default");
464 ((AliMUONv1*)MUON)->SetStepManagerVersionDE(true);
465 MUON->AddGeometryBuilder(new AliMUONSt1GeometryBuilderV2(MUON));
466 MUON->AddGeometryBuilder(new AliMUONSt2GeometryBuilderV2(MUON));
467 MUON->AddGeometryBuilder(new AliMUONSlatGeometryBuilder(MUON));
468 MUON->AddGeometryBuilder(new AliMUONTriggerGeometryBuilder(MUON));
469 }
470 //=================== PHOS parameters ===========================
471
472 if (iPHOS)
473 {
474 AliPHOS *PHOS = new AliPHOSv1("PHOS", "IHEP");
475 }
476
477
478 if (iPMD)
479 {
480 //=================== PMD parameters ============================
481 AliPMD *PMD = new AliPMDv1("PMD", "normal PMD");
482 }
483
484 if (iSTART)
485 {
486 //=================== START parameters ============================
487 AliSTART *START = new AliSTARTv1("START", "START Detector");
488 }
489
490 if (iEMCAL)
491 {
492 //=================== EMCAL parameters ============================
493 AliEMCAL *EMCAL = new AliEMCALv1("EMCAL", "EMCAL_55_25");
494 }
495
496 if (iCRT)
497 {
498 //=================== CRT parameters ============================
499 AliCRT *CRT = new AliCRTv0("CRT", "normal ACORDE");
500 }
501
502 if (iVZERO)
503 {
504 //=================== CRT parameters ============================
505 AliVZERO *VZERO = new AliVZEROv5("VZERO", "normal VZERO");
506 }
507
508
509}
510
511Float_t EtaToTheta(Float_t arg){
512 return (180./TMath::Pi())*2.*atan(exp(-arg));
513}
514
515
516
517AliGenerator* GeneratorFactory(PprRun_t srun) {
518 Int_t isw = 3;
519 if (srad == kNoGluonRadiation) isw = 0;
520
521
522 AliGenerator * gGener = 0x0;
523 switch (srun) {
524 case test50:
525 {
526 comment = comment.Append(":HIJINGparam test 50 particles");
527 AliGenHIJINGpara *gener = new AliGenHIJINGpara(50);
528 gener->SetMomentumRange(0, 999999.);
529 gener->SetPhiRange(0., 360.);
530 // Set pseudorapidity range from -8 to 8.
531 Float_t thmin = EtaToTheta(8); // theta min. <---> eta max
532 Float_t thmax = EtaToTheta(-8); // theta max. <---> eta min
533 gener->SetThetaRange(thmin,thmax);
534 gGener=gener;
535 }
536 break;
537 case kParam_8000:
538 {
539 comment = comment.Append(":HIJINGparam N=8000");
540 AliGenHIJINGpara *gener = new AliGenHIJINGpara(86030);
541 gener->SetMomentumRange(0, 999999.);
542 gener->SetPhiRange(0., 360.);
543 // Set pseudorapidity range from -8 to 8.
544 Float_t thmin = EtaToTheta(8); // theta min. <---> eta max
545 Float_t thmax = EtaToTheta(-8); // theta max. <---> eta min
546 gener->SetThetaRange(thmin,thmax);
547 gGener=gener;
548 }
549 break;
550 case kParam_4000:
551 {
552 comment = comment.Append("HIJINGparam N=4000");
553 AliGenHIJINGpara *gener = new AliGenHIJINGpara(43015);
554 gener->SetMomentumRange(0, 999999.);
555 gener->SetPhiRange(0., 360.);
556 // Set pseudorapidity range from -8 to 8.
557 Float_t thmin = EtaToTheta(8); // theta min. <---> eta max
558 Float_t thmax = EtaToTheta(-8); // theta max. <---> eta min
559 gener->SetThetaRange(thmin,thmax);
560 gGener=gener;
561 }
562 break;
563 case kParam_2000:
564 {
565 comment = comment.Append("HIJINGparam N=2000");
566 AliGenHIJINGpara *gener = new AliGenHIJINGpara(21507);
567 gener->SetMomentumRange(0, 999999.);
568 gener->SetPhiRange(0., 360.);
569 // Set pseudorapidity range from -8 to 8.
570 Float_t thmin = EtaToTheta(8); // theta min. <---> eta max
571 Float_t thmax = EtaToTheta(-8); // theta max. <---> eta min
572 gener->SetThetaRange(thmin,thmax);
573 gGener=gener;
574 }
575 break;
576//
577// Hijing Central
578//
579 case kHijing_cent1:
580 {
581 comment = comment.Append("HIJING cent1");
582 AliGenHijing *gener = HijingStandard();
583// impact parameter range
584 gener->SetImpactParameterRange(0., 5.);
585 gGener=gener;
586 }
587 break;
588 case kHijing_cent2:
589 {
590 comment = comment.Append("HIJING cent2");
591 AliGenHijing *gener = HijingStandard();
592// impact parameter range
593 gener->SetImpactParameterRange(0., 2.);
594 gGener=gener;
595 }
596 break;
597//
598// Hijing Peripheral
599//
600 case kHijing_per1:
601 {
602 comment = comment.Append("HIJING per1");
603 AliGenHijing *gener = HijingStandard();
604// impact parameter range
605 gener->SetImpactParameterRange(5., 8.6);
606 gGener=gener;
607 }
608 break;
609 case kHijing_per2:
610 {
611 comment = comment.Append("HIJING per2");
612 AliGenHijing *gener = HijingStandard();
613// impact parameter range
614 gener->SetImpactParameterRange(8.6, 11.2);
615 gGener=gener;
616 }
617 break;
618 case kHijing_per3:
619 {
620 comment = comment.Append("HIJING per3");
621 AliGenHijing *gener = HijingStandard();
622// impact parameter range
623 gener->SetImpactParameterRange(11.2, 13.2);
624 gGener=gener;
625 }
626 break;
627 case kHijing_per4:
628 {
629 comment = comment.Append("HIJING per4");
630 AliGenHijing *gener = HijingStandard();
631// impact parameter range
632 gener->SetImpactParameterRange(13.2, 15.);
633 gGener=gener;
634 }
635 break;
636 case kHijing_per5:
637 {
638 comment = comment.Append("HIJING per5");
639 AliGenHijing *gener = HijingStandard();
640// impact parameter range
641 gener->SetImpactParameterRange(15., 100.);
642 gGener=gener;
643 }
644 break;
645 case kCocktailTRD:
646 {
647 comment = comment.Append(" Cocktail for TRD at 5.5 TeV");
648 AliGenCocktail *gener = new AliGenCocktail();
649
650 AliGenParam *phi = new AliGenParam(10,
651 new AliGenMUONlib(),
652 AliGenMUONlib::kPhi,
653 "Vogt PbPb");
654
655 phi->SetPtRange(0, 100);
656 phi->SetYRange(-1., +1.);
657 phi->SetForceDecay(kDiElectron);
658
659 AliGenParam *omega = new AliGenParam(10,
660 new AliGenMUONlib(),
661 AliGenMUONlib::kOmega,
662 "Vogt PbPb");
663
664 omega->SetPtRange(0, 100);
665 omega->SetYRange(-1., +1.);
666 omega->SetForceDecay(kDiElectron);
667
668 AliGenParam *jpsi = new AliGenParam(10,
669 new AliGenMUONlib(),
670 AliGenMUONlib::kJpsiFamily,
671 "Vogt PbPb");
672
673 jpsi->SetPtRange(0, 100);
674 jpsi->SetYRange(-1., +1.);
675 jpsi->SetForceDecay(kDiElectron);
676
677 AliGenParam *ups = new AliGenParam(10,
678 new AliGenMUONlib(),
679 AliGenMUONlib::kUpsilonFamily,
680 "Vogt PbPb");
681 ups->SetPtRange(0, 100);
682 ups->SetYRange(-1., +1.);
683 ups->SetForceDecay(kDiElectron);
684
685 AliGenParam *charm = new AliGenParam(10,
686 new AliGenMUONlib(),
687 AliGenMUONlib::kCharm,
688 "central");
689 charm->SetPtRange(0, 100);
690 charm->SetYRange(-1.5, +1.5);
691 charm->SetForceDecay(kSemiElectronic);
692
693
694 AliGenParam *beauty = new AliGenParam(10,
695 new AliGenMUONlib(),
696 AliGenMUONlib::kBeauty,
697 "central");
698 beauty->SetPtRange(0, 100);
699 beauty->SetYRange(-1.5, +1.5);
700 beauty->SetForceDecay(kSemiElectronic);
701
702 AliGenParam *beautyJ = new AliGenParam(10,
703 new AliGenMUONlib(),
704 AliGenMUONlib::kBeauty,
705 "central");
706 beautyJ->SetPtRange(0, 100);
707 beautyJ->SetYRange(-1.5, +1.5);
708 beautyJ->SetForceDecay(kBJpsiDiElectron);
709
710 gener->AddGenerator(phi,"Phi",1);
711 gener->AddGenerator(omega,"Omega",1);
712 gener->AddGenerator(jpsi,"J/psi",1);
713 gener->AddGenerator(ups,"Upsilon",1);
714 gener->AddGenerator(charm,"Charm",1);
715 gener->AddGenerator(beauty,"Beauty",1);
716 gener->AddGenerator(beautyJ,"J/Psi from Beauty",1);
717 gGener=gener;
718 }
719 break;
720 case kpieTRD:
721 {
722 comment = comment.Append("e, pi for TRD");
723 AliGenCocktail *gener = new AliGenCocktail();
724 Double_t momen=2.0;
725 Int_t Npart=200;
726 AliGenBox *electron = new AliGenBox(Npart);
727 electron->SetMomentumRange(momen,momen);
728 electron->SetPhiRange(0,360);
729 Float_t thmin = EtaToTheta(.9); // theta min. <---> eta max
730 Float_t thmax = EtaToTheta(-.9); // theta max. <---> eta min
731 electron->SetThetaRange(thmin,thmax);
732 electron->SetOrigin(0,0,0); //vertex position
733 electron->SetSigma(0,0,0); //Sigma in (X,Y,Z) (cm) on IP position
734 electron->SetPart(11); //GEANT particle type
735
736 AliGenBox *pion = new AliGenBox(Npart);
737 pion->SetMomentumRange(momen,momen);
738 pion->SetPhiRange(0,360);
739 pion->SetThetaRange(thmin,thmax);
740 pion->SetOrigin(0,0,0); //vertex position
741 pion->SetSigma(0,0,0); //Sigma in (X,Y,Z) (cm) on IP position
742 pion->SetPart(211); //GEANT particle type
743
744 AliGenBox *muon = new AliGenBox(Npart);
745 muon->SetMomentumRange(momen,momen);
746 muon->SetPhiRange(0,360);
747 muon->SetThetaRange(thmin,thmax);
748 muon->SetOrigin(0,0,0); //vertex position
749 muon->SetSigma(0,0,0); //Sigma in (X,Y,Z) (cm) on IP position
750 muon->SetPart(13); //GEANT particle type
751
752 AliGenBox *kaon = new AliGenBox(Npart);
753 kaon->SetMomentumRange(momen,momen);
754 kaon->SetPhiRange(0,360);
755 kaon->SetThetaRange(thmin,thmax);
756 kaon->SetOrigin(0,0,0); //vertex position
757 kaon->SetSigma(0,0,0); //Sigma in (X,Y,Z) (cm) on IP position
758 kaon->SetPart(321); //GEANT particle type
759
760 AliGenBox *proton = new AliGenBox(Npart);
761 proton->SetMomentumRange(momen,momen);
762 proton->SetPhiRange(0,360);
763 proton->SetThetaRange(thmin,thmax);
764 proton->SetOrigin(0,0,0); //vertex position
765 proton->SetSigma(0,0,0); //Sigma in (X,Y,Z) (cm) on IP position
766 proton->SetPart(2212); //GEANT particle type
767
768 gener->AddGenerator(electron,"electron",1);
769 gener->AddGenerator(pion,"pion",1);
770 gener->AddGenerator(muon,"muon",1);
771 gener->AddGenerator(kaon,"kaon",1);
772 gener->AddGenerator(proton,"proton",1);
773 gGener=gener;
774 }
775 break;
776 default: break;
777 }
778 return gGener;
779}
780
781void ProcessEnvironmentVars()
782{
783 // Run type
784 if (gSystem->Getenv("CONFIG_RUN_TYPE")) {
785 for (Int_t iRun = 0; iRun < kRunMax; iRun++) {
786 if (strcmp(gSystem->Getenv("CONFIG_RUN_TYPE"), pprRunName[iRun])==0) {
787 srun = (PprRun_t)iRun;
788 cout<<"Run type set to "<<pprRunName[iRun]<<endl;
789 }
790 }
791 }
792
793 // Random Number seed
794 if (gSystem->Getenv("CONFIG_SEED")) {
795 sseed = atoi(gSystem->Getenv("CONFIG_SEED"));
796 }
797}