e3d91248 |
1 | // |
2 | // Configuration for the Physics Data Challenge 2006 |
3 | // |
4 | |
13217d97 |
5 | // One can use the configuration macro in compiled mode by |
6 | // root [0] gSystem->Load("libgeant321"); |
7 | // root [0] gSystem->SetIncludePath("-I$ROOTSYS/include -I$ALICE_ROOT/include\ |
8 | // -I$ALICE_ROOT -I$ALICE/geant3/TGeant3"); |
9 | // root [0] .x grun.C(1,"Config_PDC06.C++") |
10 | |
e3d91248 |
11 | #if !defined(__CINT__) || defined(__MAKECINT__) |
12 | #include <Riostream.h> |
13 | #include <TRandom.h> |
14 | #include <TDatime.h> |
15 | #include <TSystem.h> |
16 | #include <TVirtualMC.h> |
13217d97 |
17 | #include <TGeant3TGeo.h> |
18 | #include "EVGEN/AliGenCocktail.h" |
19 | #include "EVGEN/AliGenParam.h" |
20 | #include "EVGEN/AliGenMUONlib.h" |
e3d91248 |
21 | #include "STEER/AliRunLoader.h" |
22 | #include "STEER/AliRun.h" |
23 | #include "STEER/AliConfig.h" |
24 | #include "PYTHIA6/AliDecayerPythia.h" |
25 | #include "PYTHIA6/AliGenPythia.h" |
26 | #include "STEER/AliMagFMaps.h" |
27 | #include "STRUCT/AliBODY.h" |
28 | #include "STRUCT/AliMAG.h" |
29 | #include "STRUCT/AliABSOv0.h" |
30 | #include "STRUCT/AliDIPOv2.h" |
31 | #include "STRUCT/AliHALL.h" |
32 | #include "STRUCT/AliFRAMEv2.h" |
33 | #include "STRUCT/AliSHILv2.h" |
34 | #include "STRUCT/AliPIPEv0.h" |
13217d97 |
35 | #include "ITS/AliITSgeom.h" |
e3d91248 |
36 | #include "ITS/AliITSvPPRasymmFMD.h" |
37 | #include "TPC/AliTPCv2.h" |
13217d97 |
38 | #include "TOF/AliTOFv5T0.h" |
e3d91248 |
39 | #include "RICH/AliRICHv1.h" |
13217d97 |
40 | #include "ZDC/AliZDCv2.h" |
e3d91248 |
41 | #include "TRD/AliTRDv1.h" |
13217d97 |
42 | #include "FMD/AliFMDv1.h" |
e3d91248 |
43 | #include "MUON/AliMUONv1.h" |
44 | #include "PHOS/AliPHOSv1.h" |
45 | #include "PMD/AliPMDv1.h" |
46 | #include "START/AliSTARTv1.h" |
13217d97 |
47 | #include "EMCAL/AliEMCALv2.h" |
48 | #include "CRT/AliCRTv0.h" |
4a2f6442 |
49 | #include "VZERO/AliVZEROv7.h" |
e3d91248 |
50 | #endif |
51 | |
52 | |
53 | enum PDC06Proc_t |
54 | { |
55 | //--- Heavy Flavour Production --- |
56 | kCharmPbPb5500, kCharmpPb8800, kCharmpp14000, kCharmpp14000wmi, |
57 | kD0PbPb5500, kD0pPb8800, kD0pp14000, |
58 | kDPlusPbPb5500, kDPluspPb8800, kDPluspp14000, |
59 | kBeautyPbPb5500, kBeautypPb8800, kBeautypp14000, kBeautypp14000wmi, |
60 | // -- Pythia Mb |
13217d97 |
61 | kPyMbNoHvq, kPyOmegaPlus, kPyOmegaMinus, kRunMax |
62 | }; |
63 | |
64 | const char * pprRunName[] = { |
65 | "kCharmPbPb5500", "kCharmpPb8800", "kCharmpp14000", "kCharmpp14000wmi", |
66 | "kD0PbPb5500", "kD0pPb8800", "kD0pp14000", |
67 | "kDPlusPbPb5500", "kDPluspPb8800", "kDPluspp14000", |
68 | "kBeautyPbPb5500", "kBeautypPb8800", "kBeautypp14000", "kBeautypp14000wmi", |
69 | "kPyMbNoHvq", "kPyOmegaPlus", "kPyOmegaMinus" |
e3d91248 |
70 | }; |
71 | |
72 | |
73 | //--- Decay Mode --- |
74 | enum DecayHvFl_t |
75 | { |
76 | kNature, kHadr, kSemiEl, kSemiMu |
77 | }; |
78 | //--- Rapidity Cut --- |
79 | enum YCut_t |
80 | { |
81 | kFull, kBarrel, kMuonArm |
82 | }; |
83 | //--- Magnetic Field --- |
84 | enum Mag_t |
85 | { |
86 | k2kG, k4kG, k5kG |
87 | }; |
6521ae06 |
88 | |
89 | //--- Trigger config --- |
90 | enum TrigConf_t |
91 | { |
92 | kDefaultPPTrig, kDefaultPbPbTrig |
93 | }; |
94 | |
95 | const char * TrigConfName[] = { |
96 | "p-p","Pb-Pb" |
97 | }; |
98 | |
e3d91248 |
99 | //--- Functions --- |
100 | AliGenPythia *PythiaHVQ(PDC06Proc_t proc); |
101 | AliGenerator *MbCocktail(); |
102 | AliGenerator *PyMbTriggered(Int_t pdg); |
13217d97 |
103 | void ProcessEnvironmentVars(); |
e3d91248 |
104 | |
105 | // This part for configuration |
a7e9beb2 |
106 | static PDC06Proc_t proc = kPyOmegaPlus; |
e3d91248 |
107 | static DecayHvFl_t decHvFl = kNature; |
108 | static YCut_t ycut = kFull; |
109 | static Mag_t mag = k5kG; |
6521ae06 |
110 | static TrigConf_t trig = kDefaultPPTrig; // default pp trigger configuration |
13217d97 |
111 | //========================// |
112 | // Set Random Number seed // |
113 | //========================// |
114 | TDatime dt; |
115 | static UInt_t seed = dt.Get(); |
116 | |
e3d91248 |
117 | // nEvts = -1 : you get 1 QQbar pair and all the fragmentation and |
118 | // decay chain |
119 | // nEvts = N>0 : you get N charm / beauty Hadrons |
120 | Int_t nEvts = -1; |
121 | // stars = kTRUE : all heavy resonances and their decay stored |
122 | // = kFALSE: only final heavy hadrons and their decays stored |
123 | Bool_t stars = kTRUE; |
124 | |
125 | // To be used only with kCharmppMNRwmi and kBeautyppMNRwmi |
126 | // To get a "reasonable" agreement with MNR results, events have to be |
127 | // generated with the minimum ptHard set to 2.76 GeV. |
128 | // To get a "perfect" agreement with MNR results, events have to be |
129 | // generated in four ptHard bins with the following relative |
130 | // normalizations: |
131 | // CHARM |
132 | // 2.76-3 GeV: 25% |
133 | // 3-4 GeV: 40% |
134 | // 4-8 GeV: 29% |
135 | // >8 GeV: 6% |
136 | // BEAUTY |
137 | // 2.76-4 GeV: 5% |
138 | // 4-6 GeV: 31% |
139 | // 6-8 GeV: 28% |
140 | // >8 GeV: 36% |
141 | Float_t ptHardMin = 2.76; |
142 | Float_t ptHardMax = -1.; |
143 | |
144 | |
145 | // Comment line |
146 | static TString comment; |
147 | |
148 | void Config() |
149 | { |
150 | |
13217d97 |
151 | |
152 | // Get settings from environment variables |
153 | ProcessEnvironmentVars(); |
154 | |
155 | gRandom->SetSeed(seed); |
156 | cerr<<"Seed for random number generation= "<<seed<<endl; |
e3d91248 |
157 | |
158 | // libraries required by geant321 |
159 | #if defined(__CINT__) |
160 | gSystem->Load("libgeant321"); |
161 | #endif |
162 | |
163 | new TGeant3TGeo("C++ Interface to Geant3"); |
164 | |
165 | //======================================================================= |
166 | // Create the output file |
167 | |
168 | |
169 | AliRunLoader* rl=0x0; |
170 | |
171 | cout<<"Config.C: Creating Run Loader ..."<<endl; |
172 | rl = AliRunLoader::Open("galice.root", |
173 | AliConfig::GetDefaultEventFolderName(), |
174 | "recreate"); |
175 | if (rl == 0x0) |
176 | { |
177 | gAlice->Fatal("Config.C","Can not instatiate the Run Loader"); |
178 | return; |
179 | } |
180 | rl->SetCompressionLevel(2); |
181 | rl->SetNumberOfEventsPerFile(1000); |
182 | gAlice->SetRunLoader(rl); |
6521ae06 |
183 | |
184 | // Set the trigger configuration |
185 | gAlice->SetTriggerDescriptor(TrigConfName[trig]); |
186 | cout<<"Trigger configuration is set to "<<TrigConfName[trig]<<endl; |
e3d91248 |
187 | |
188 | // |
189 | //======================================================================= |
190 | // ************* STEERING parameters FOR ALICE SIMULATION ************** |
191 | // --- Specify event type to be tracked through the ALICE setup |
192 | // --- All positions are in cm, angles in degrees, and P and E in GeV |
193 | |
194 | |
195 | gMC->SetProcess("DCAY",1); |
196 | gMC->SetProcess("PAIR",1); |
197 | gMC->SetProcess("COMP",1); |
198 | gMC->SetProcess("PHOT",1); |
199 | gMC->SetProcess("PFIS",0); |
200 | gMC->SetProcess("DRAY",0); |
201 | gMC->SetProcess("ANNI",1); |
202 | gMC->SetProcess("BREM",1); |
203 | gMC->SetProcess("MUNU",1); |
204 | gMC->SetProcess("CKOV",1); |
205 | gMC->SetProcess("HADR",1); |
206 | gMC->SetProcess("LOSS",2); |
207 | gMC->SetProcess("MULS",1); |
208 | gMC->SetProcess("RAYL",1); |
209 | |
210 | Float_t cut = 1.e-3; // 1MeV cut by default |
211 | Float_t tofmax = 1.e10; |
212 | |
213 | gMC->SetCut("CUTGAM", cut); |
214 | gMC->SetCut("CUTELE", cut); |
215 | gMC->SetCut("CUTNEU", cut); |
216 | gMC->SetCut("CUTHAD", cut); |
217 | gMC->SetCut("CUTMUO", cut); |
218 | gMC->SetCut("BCUTE", cut); |
219 | gMC->SetCut("BCUTM", cut); |
220 | gMC->SetCut("DCUTE", cut); |
221 | gMC->SetCut("DCUTM", cut); |
222 | gMC->SetCut("PPCUTM", cut); |
223 | gMC->SetCut("TOFMAX", tofmax); |
224 | |
225 | |
226 | |
227 | |
228 | // Set External decayer // |
229 | //======================// |
230 | TVirtualMCDecayer* decayer = new AliDecayerPythia(); |
231 | // DECAYS |
232 | // |
233 | switch(decHvFl) { |
234 | case kNature: |
235 | decayer->SetForceDecay(kAll); |
236 | break; |
237 | case kHadr: |
238 | decayer->SetForceDecay(kHadronicD); |
239 | break; |
240 | case kSemiEl: |
241 | decayer->SetForceDecay(kSemiElectronic); |
242 | break; |
243 | case kSemiMu: |
244 | decayer->SetForceDecay(kSemiMuonic); |
245 | break; |
246 | } |
247 | decayer->Init(); |
248 | gMC->SetExternalDecayer(decayer); |
249 | |
250 | //=========================// |
251 | // Generator Configuration // |
252 | //=========================// |
13217d97 |
253 | AliGenerator* gener = 0x0; |
e3d91248 |
254 | |
255 | if (proc <= kBeautypp14000wmi) { |
256 | AliGenPythia *pythia = PythiaHVQ(proc); |
257 | // FeedDown option |
258 | pythia->SetFeedDownHigherFamily(kFALSE); |
259 | // Stack filling option |
260 | if(!stars) pythia->SetStackFillOpt(AliGenPythia::kParentSelection); |
261 | // Set Count mode |
262 | if(nEvts>0) pythia->SetCountMode(AliGenPythia::kCountParents); |
263 | // |
264 | // DECAYS |
265 | // |
266 | switch(decHvFl) { |
267 | case kNature: |
268 | pythia->SetForceDecay(kAll); |
269 | break; |
270 | case kHadr: |
271 | pythia->SetForceDecay(kHadronicD); |
272 | break; |
273 | case kSemiEl: |
274 | pythia->SetForceDecay(kSemiElectronic); |
275 | break; |
276 | case kSemiMu: |
277 | pythia->SetForceDecay(kSemiMuonic); |
278 | break; |
279 | } |
280 | // |
281 | // GEOM & KINE CUTS |
282 | // |
283 | pythia->SetMomentumRange(0,99999999); |
284 | pythia->SetPhiRange(0., 360.); |
285 | pythia->SetThetaRange(0,180); |
286 | switch(ycut) { |
287 | case kFull: |
288 | pythia->SetYRange(-999,999); |
289 | break; |
290 | case kBarrel: |
291 | pythia->SetYRange(-2,2); |
292 | break; |
293 | case kMuonArm: |
294 | pythia->SetYRange(1,6); |
295 | break; |
296 | } |
297 | gener = pythia; |
298 | } else if (proc == kPyMbNoHvq) { |
299 | gener = MbCocktail(); |
300 | } else if (proc == kPyOmegaMinus) { |
301 | gener = PyMbTriggered(3334); |
302 | } else if (proc == kPyOmegaPlus) { |
303 | gener = PyMbTriggered(-3334); |
304 | } |
305 | |
306 | |
307 | |
308 | // PRIMARY VERTEX |
309 | // |
310 | gener->SetOrigin(0., 0., 0.); // vertex position |
311 | // |
312 | // |
313 | // Size of the interaction diamond |
314 | // Longitudinal |
315 | Float_t sigmaz = 7.55 / TMath::Sqrt(2.); // [cm] |
316 | // |
317 | // Transverse |
318 | Float_t betast = 10; // beta* [m] |
319 | Float_t eps = 3.75e-6; // emittance [m] |
320 | Float_t gamma = 7000. / 0.938272; // relativistic gamma [1] |
321 | Float_t sigmaxy = TMath::Sqrt(eps * betast / gamma) / TMath::Sqrt(2.) * 100.; // [cm] |
322 | printf("\n \n Diamond size x-y: %10.3e z: %10.3e\n \n", sigmaxy, sigmaz); |
323 | |
324 | gener->SetSigma(sigmaxy, sigmaxy, sigmaz); // Sigma in (X,Y,Z) (cm) on IP position |
550c2171 |
325 | gener->SetCutVertexZ(3.); // Truncate at 3 sigma |
e3d91248 |
326 | gener->SetVertexSmear(kPerEvent); |
327 | |
e3d91248 |
328 | gener->Init(); |
329 | |
330 | // FIELD |
331 | // |
332 | if (mag == k2kG) { |
333 | comment = comment.Append(" | L3 field 0.2 T"); |
334 | } else if (mag == k4kG) { |
335 | comment = comment.Append(" | L3 field 0.4 T"); |
336 | } else if (mag == k5kG) { |
337 | comment = comment.Append(" | L3 field 0.5 T"); |
338 | } |
339 | printf("\n \n Comment: %s \n \n", comment.Data()); |
340 | |
341 | AliMagFMaps* field = new AliMagFMaps("Maps","Maps", 2, 1., 10., mag); |
342 | field->SetL3ConstField(0); //Using const. field in the barrel |
343 | rl->CdGAFile(); |
344 | gAlice->SetField(field); |
345 | |
346 | |
347 | |
348 | Int_t iABSO = 1; |
4a2f6442 |
349 | Int_t iCRT = 0; |
e3d91248 |
350 | Int_t iDIPO = 1; |
b9f5ae9d |
351 | Int_t iEMCAL = 1; |
e3d91248 |
352 | Int_t iFMD = 1; |
353 | Int_t iFRAME = 1; |
354 | Int_t iHALL = 1; |
355 | Int_t iITS = 1; |
356 | Int_t iMAG = 1; |
357 | Int_t iMUON = 1; |
b9f5ae9d |
358 | Int_t iPHOS = 1; |
e3d91248 |
359 | Int_t iPIPE = 1; |
360 | Int_t iPMD = 1; |
361 | Int_t iRICH = 1; |
362 | Int_t iSHIL = 1; |
363 | Int_t iSTART = 1; |
b9f5ae9d |
364 | Int_t iTOF = 1; |
e3d91248 |
365 | Int_t iTPC = 1; |
366 | Int_t iTRD = 1; |
b9f5ae9d |
367 | Int_t iVZERO = 1; |
e3d91248 |
368 | Int_t iZDC = 1; |
b9f5ae9d |
369 | |
e3d91248 |
370 | |
b9f5ae9d |
371 | //=================== Alice BODY parameters ============================= |
372 | AliBODY *BODY = new AliBODY("BODY", "Alice envelop"); |
e3d91248 |
373 | |
e3d91248 |
374 | |
b9f5ae9d |
375 | if (iMAG) |
376 | { |
377 | //=================== MAG parameters ============================ |
378 | // --- Start with Magnet since detector layouts may be depending --- |
379 | // --- on the selected Magnet dimensions --- |
380 | AliMAG *MAG = new AliMAG("MAG", "Magnet"); |
381 | } |
e3d91248 |
382 | |
e3d91248 |
383 | |
b9f5ae9d |
384 | if (iABSO) |
385 | { |
386 | //=================== ABSO parameters ============================ |
387 | AliABSO *ABSO = new AliABSOv0("ABSO", "Muon Absorber"); |
388 | } |
e3d91248 |
389 | |
b9f5ae9d |
390 | if (iDIPO) |
391 | { |
392 | //=================== DIPO parameters ============================ |
e3d91248 |
393 | |
b9f5ae9d |
394 | AliDIPO *DIPO = new AliDIPOv2("DIPO", "Dipole version 2"); |
395 | } |
e3d91248 |
396 | |
b9f5ae9d |
397 | if (iHALL) |
398 | { |
399 | //=================== HALL parameters ============================ |
e3d91248 |
400 | |
b9f5ae9d |
401 | AliHALL *HALL = new AliHALL("HALL", "Alice Hall"); |
402 | } |
e3d91248 |
403 | |
e3d91248 |
404 | |
b9f5ae9d |
405 | if (iFRAME) |
406 | { |
407 | //=================== FRAME parameters ============================ |
e3d91248 |
408 | |
b9f5ae9d |
409 | AliFRAMEv2 *FRAME = new AliFRAMEv2("FRAME", "Space Frame"); |
b9f5ae9d |
410 | } |
e3d91248 |
411 | |
b9f5ae9d |
412 | if (iSHIL) |
413 | { |
414 | //=================== SHIL parameters ============================ |
e3d91248 |
415 | |
b9f5ae9d |
416 | AliSHIL *SHIL = new AliSHILv2("SHIL", "Shielding Version 2"); |
417 | } |
e3d91248 |
418 | |
e3d91248 |
419 | |
b9f5ae9d |
420 | if (iPIPE) |
421 | { |
422 | //=================== PIPE parameters ============================ |
e3d91248 |
423 | |
b9f5ae9d |
424 | AliPIPE *PIPE = new AliPIPEv0("PIPE", "Beam Pipe"); |
425 | } |
426 | |
427 | if(iITS) { |
e3d91248 |
428 | |
b9f5ae9d |
429 | //=================== ITS parameters ============================ |
e3d91248 |
430 | // |
431 | // As the innermost detector in ALICE, the Inner Tracking System "impacts" on |
432 | // almost all other detectors. This involves the fact that the ITS geometry |
433 | // still has several options to be followed in parallel in order to determine |
434 | // the best set-up which minimizes the induced background. All the geometries |
435 | // available to date are described in the following. Read carefully the comments |
436 | // and use the default version (the only one uncommented) unless you are making |
437 | // comparisons and you know what you are doing. In this case just uncomment the |
438 | // ITS geometry you want to use and run Aliroot. |
439 | // |
440 | // Detailed geometries: |
441 | // |
442 | // |
443 | //AliITS *ITS = new AliITSv5symm("ITS","Updated ITS TDR detailed version with symmetric services"); |
444 | // |
445 | //AliITS *ITS = new AliITSv5asymm("ITS","Updates ITS TDR detailed version with asymmetric services"); |
446 | // |
b9f5ae9d |
447 | AliITSvPPRasymmFMD *ITS = new AliITSvPPRasymmFMD("ITS","New ITS PPR detailed version with asymmetric services"); |
448 | ITS->SetMinorVersion(2); // don't touch this parameter if you're not an ITS developer |
449 | ITS->SetReadDet(kTRUE); // don't touch this parameter if you're not an ITS developer |
450 | // ITS->SetWriteDet("$ALICE_ROOT/ITS/ITSgeometry_vPPRasymm2.det"); // don't touch this parameter if you're not an ITS developer |
451 | ITS->SetThicknessDet1(200.); // detector thickness on layer 1 must be in the range [100,300] |
452 | ITS->SetThicknessDet2(200.); // detector thickness on layer 2 must be in the range [100,300] |
453 | ITS->SetThicknessChip1(200.); // chip thickness on layer 1 must be in the range [150,300] |
454 | ITS->SetThicknessChip2(200.); // chip thickness on layer 2 must be in the range [150,300] |
455 | ITS->SetRails(0); // 1 --> rails in ; 0 --> rails out |
456 | ITS->SetCoolingFluid(1); // 1 --> water ; 0 --> freon |
457 | |
e3d91248 |
458 | // Coarse geometries (warning: no hits are produced with these coarse geometries and they unuseful |
459 | // for reconstruction !): |
460 | // |
461 | // |
462 | //AliITSvPPRcoarseasymm *ITS = new AliITSvPPRcoarseasymm("ITS","New ITS PPR coarse version with asymmetric services"); |
b9f5ae9d |
463 | //ITS->SetRails(0); // 1 --> rails in ; 0 --> rails out |
e3d91248 |
464 | //ITS->SetSupportMaterial(0); // 0 --> Copper ; 1 --> Aluminum ; 2 --> Carbon |
465 | // |
466 | //AliITS *ITS = new AliITSvPPRcoarsesymm("ITS","New ITS PPR coarse version with symmetric services"); |
b9f5ae9d |
467 | //ITS->SetRails(0); // 1 --> rails in ; 0 --> rails out |
e3d91248 |
468 | //ITS->SetSupportMaterial(0); // 0 --> Copper ; 1 --> Aluminum ; 2 --> Carbon |
469 | // |
470 | // |
471 | // |
472 | // Geant3 <-> EUCLID conversion |
473 | // ============================ |
474 | // |
475 | // SetEUCLID is a flag to output (=1) or not to output (=0) both geometry and |
476 | // media to two ASCII files (called by default ITSgeometry.euc and |
477 | // ITSgeometry.tme) in a format understandable to the CAD system EUCLID. |
478 | // The default (=0) means that you dont want to use this facility. |
479 | // |
b9f5ae9d |
480 | ITS->SetEUCLID(0); |
481 | } |
e3d91248 |
482 | |
b9f5ae9d |
483 | if (iTPC) |
484 | { |
485 | //============================ TPC parameters ===================== |
486 | AliTPC *TPC = new AliTPCv2("TPC", "Default"); |
487 | } |
e3d91248 |
488 | |
e3d91248 |
489 | |
b9f5ae9d |
490 | if (iTOF) { |
491 | //=================== TOF parameters ============================ |
492 | AliTOF *TOF = new AliTOFv5T0("TOF", "normal TOF"); |
493 | } |
e3d91248 |
494 | |
e3d91248 |
495 | |
b9f5ae9d |
496 | if (iRICH) |
497 | { |
498 | //=================== RICH parameters =========================== |
499 | AliRICH *RICH = new AliRICHv1("RICH", "normal RICH"); |
e3d91248 |
500 | |
b9f5ae9d |
501 | } |
e3d91248 |
502 | |
e3d91248 |
503 | |
b9f5ae9d |
504 | if (iZDC) |
505 | { |
506 | //=================== ZDC parameters ============================ |
e3d91248 |
507 | |
b9f5ae9d |
508 | AliZDC *ZDC = new AliZDCv2("ZDC", "normal ZDC"); |
509 | } |
e3d91248 |
510 | |
b9f5ae9d |
511 | if (iTRD) |
512 | { |
513 | //=================== TRD parameters ============================ |
e3d91248 |
514 | |
b9f5ae9d |
515 | AliTRD *TRD = new AliTRDv1("TRD", "TRD slow simulator"); |
b9f5ae9d |
516 | } |
e3d91248 |
517 | |
b9f5ae9d |
518 | if (iFMD) |
519 | { |
520 | //=================== FMD parameters ============================ |
521 | AliFMD *FMD = new AliFMDv1("FMD", "normal FMD"); |
522 | } |
e3d91248 |
523 | |
b9f5ae9d |
524 | if (iMUON) |
525 | { |
526 | //=================== MUON parameters =========================== |
527 | // New MUONv1 version (geometry defined via builders) |
528 | AliMUON *MUON = new AliMUONv1("MUON", "default"); |
529 | } |
530 | //=================== PHOS parameters =========================== |
e3d91248 |
531 | |
b9f5ae9d |
532 | if (iPHOS) |
533 | { |
534 | AliPHOS *PHOS = new AliPHOSv1("PHOS", "IHEP"); |
535 | } |
e3d91248 |
536 | |
537 | |
b9f5ae9d |
538 | if (iPMD) |
539 | { |
540 | //=================== PMD parameters ============================ |
541 | AliPMD *PMD = new AliPMDv1("PMD", "normal PMD"); |
542 | } |
e3d91248 |
543 | |
b9f5ae9d |
544 | if (iSTART) |
545 | { |
546 | //=================== START parameters ============================ |
547 | AliSTART *START = new AliSTARTv1("START", "START Detector"); |
548 | } |
e3d91248 |
549 | |
b9f5ae9d |
550 | if (iEMCAL) |
551 | { |
552 | //=================== EMCAL parameters ============================ |
0588179e |
553 | AliEMCAL *EMCAL = new AliEMCALv2("EMCAL", "SHISH_77_TRD1_2X2_FINAL_110DEG"); |
b9f5ae9d |
554 | } |
e3d91248 |
555 | |
b9f5ae9d |
556 | if (iCRT) |
557 | { |
558 | //=================== CRT parameters ============================ |
559 | AliCRT *CRT = new AliCRTv0("CRT", "normal ACORDE"); |
560 | } |
e3d91248 |
561 | |
b9f5ae9d |
562 | if (iVZERO) |
563 | { |
564 | //=================== CRT parameters ============================ |
4a2f6442 |
565 | AliVZERO *VZERO = new AliVZEROv7("VZERO", "normal VZERO"); |
b9f5ae9d |
566 | } |
e3d91248 |
567 | } |
568 | // |
569 | // PYTHIA |
570 | // |
571 | AliGenPythia *PythiaHVQ(PDC06Proc_t proc) { |
572 | //*******************************************************************// |
573 | // Configuration file for charm / beauty generation with PYTHIA // |
574 | // // |
575 | // The parameters have been tuned in order to reproduce the inclusive// |
576 | // heavy quark pt distribution given by the NLO pQCD calculation by // |
577 | // Mangano, Nason and Ridolfi. // |
578 | // // |
579 | // For details and for the NORMALIZATION of the yields see: // |
580 | // N.Carrer and A.Dainese, // |
581 | // "Charm and beauty production at the LHC", // |
582 | // ALICE-INT-2003-019, [arXiv:hep-ph/0311225]; // |
583 | // PPR Chapter 6.6, CERN/LHCC 2005-030 (2005). // |
584 | //*******************************************************************// |
585 | AliGenPythia * gener = 0x0; |
586 | |
587 | switch(proc) { |
588 | case kCharmPbPb5500: |
589 | comment = comment.Append(" Charm in Pb-Pb at 5.5 TeV"); |
590 | gener = new AliGenPythia(nEvts); |
591 | gener->SetProcess(kPyCharmPbPbMNR); |
592 | gener->SetStrucFunc(kCTEQ4L); |
593 | gener->SetPtHard(2.1,-1.0); |
594 | gener->SetEnergyCMS(5500.); |
595 | gener->SetNuclei(208,208); |
596 | break; |
597 | case kCharmpPb8800: |
598 | comment = comment.Append(" Charm in p-Pb at 8.8 TeV"); |
599 | gener = new AliGenPythia(nEvts); |
600 | gener->SetProcess(kPyCharmpPbMNR); |
601 | gener->SetStrucFunc(kCTEQ4L); |
602 | gener->SetPtHard(2.1,-1.0); |
603 | gener->SetEnergyCMS(8800.); |
604 | gener->SetProjectile("P",1,1); |
605 | gener->SetTarget("Pb",208,82); |
606 | break; |
607 | case kCharmpp14000: |
608 | comment = comment.Append(" Charm in pp at 14 TeV"); |
609 | gener = new AliGenPythia(nEvts); |
610 | gener->SetProcess(kPyCharmppMNR); |
611 | gener->SetStrucFunc(kCTEQ4L); |
612 | gener->SetPtHard(2.1,-1.0); |
613 | gener->SetEnergyCMS(14000.); |
614 | break; |
615 | case kCharmpp14000wmi: |
616 | comment = comment.Append(" Charm in pp at 14 TeV with mult. interactions"); |
617 | gener = new AliGenPythia(-1); |
618 | gener->SetProcess(kPyCharmppMNRwmi); |
619 | gener->SetStrucFunc(kCTEQ5L); |
620 | gener->SetPtHard(ptHardMin,ptHardMax); |
621 | gener->SetEnergyCMS(14000.); |
622 | break; |
623 | case kD0PbPb5500: |
624 | comment = comment.Append(" D0 in Pb-Pb at 5.5 TeV"); |
625 | gener = new AliGenPythia(nEvts); |
626 | gener->SetProcess(kPyD0PbPbMNR); |
627 | gener->SetStrucFunc(kCTEQ4L); |
628 | gener->SetPtHard(2.1,-1.0); |
629 | gener->SetEnergyCMS(5500.); |
630 | gener->SetNuclei(208,208); |
631 | break; |
632 | case kD0pPb8800: |
633 | comment = comment.Append(" D0 in p-Pb at 8.8 TeV"); |
634 | gener = new AliGenPythia(nEvts); |
635 | gener->SetProcess(kPyD0pPbMNR); |
636 | gener->SetStrucFunc(kCTEQ4L); |
637 | gener->SetPtHard(2.1,-1.0); |
638 | gener->SetEnergyCMS(8800.); |
639 | gener->SetProjectile("P",1,1); |
640 | gener->SetTarget("Pb",208,82); |
641 | break; |
642 | case kD0pp14000: |
643 | comment = comment.Append(" D0 in pp at 14 TeV"); |
644 | gener = new AliGenPythia(nEvts); |
645 | gener->SetProcess(kPyD0ppMNR); |
646 | gener->SetStrucFunc(kCTEQ4L); |
647 | gener->SetPtHard(2.1,-1.0); |
648 | gener->SetEnergyCMS(14000.); |
649 | break; |
650 | case kDPlusPbPb5500: |
651 | comment = comment.Append(" DPlus in Pb-Pb at 5.5 TeV"); |
652 | gener = new AliGenPythia(nEvts); |
653 | gener->SetProcess(kPyDPlusPbPbMNR); |
654 | gener->SetStrucFunc(kCTEQ4L); |
655 | gener->SetPtHard(2.1,-1.0); |
656 | gener->SetEnergyCMS(5500.); |
657 | gener->SetNuclei(208,208); |
658 | break; |
659 | case kDPluspPb8800: |
660 | comment = comment.Append(" DPlus in p-Pb at 8.8 TeV"); |
661 | gener = new AliGenPythia(nEvts); |
662 | gener->SetProcess(kPyDPluspPbMNR); |
663 | gener->SetStrucFunc(kCTEQ4L); |
664 | gener->SetPtHard(2.1,-1.0); |
665 | gener->SetEnergyCMS(8800.); |
666 | gener->SetProjectile("P",1,1); |
667 | gener->SetTarget("Pb",208,82); |
668 | break; |
669 | case kDPluspp14000: |
670 | comment = comment.Append(" DPlus in pp at 14 TeV"); |
671 | gener = new AliGenPythia(nEvts); |
672 | gener->SetProcess(kPyDPlusppMNR); |
673 | gener->SetStrucFunc(kCTEQ4L); |
674 | gener->SetPtHard(2.1,-1.0); |
675 | gener->SetEnergyCMS(14000.); |
676 | break; |
677 | case kBeautyPbPb5500: |
678 | comment = comment.Append(" Beauty in Pb-Pb at 5.5 TeV"); |
679 | gener = new AliGenPythia(nEvts); |
680 | gener->SetProcess(kPyBeautyPbPbMNR); |
681 | gener->SetStrucFunc(kCTEQ4L); |
682 | gener->SetPtHard(2.75,-1.0); |
683 | gener->SetEnergyCMS(5500.); |
684 | gener->SetNuclei(208,208); |
685 | break; |
686 | case kBeautypPb8800: |
687 | comment = comment.Append(" Beauty in p-Pb at 8.8 TeV"); |
688 | gener = new AliGenPythia(nEvts); |
689 | gener->SetProcess(kPyBeautypPbMNR); |
690 | gener->SetStrucFunc(kCTEQ4L); |
691 | gener->SetPtHard(2.75,-1.0); |
692 | gener->SetEnergyCMS(8800.); |
693 | gener->SetProjectile("P",1,1); |
694 | gener->SetTarget("Pb",208,82); |
695 | break; |
696 | case kBeautypp14000: |
697 | comment = comment.Append(" Beauty in pp at 14 TeV"); |
698 | gener = new AliGenPythia(nEvts); |
699 | gener->SetProcess(kPyBeautyppMNR); |
700 | gener->SetStrucFunc(kCTEQ4L); |
701 | gener->SetPtHard(2.75,-1.0); |
702 | gener->SetEnergyCMS(14000.); |
703 | break; |
704 | case kBeautypp14000wmi: |
705 | comment = comment.Append(" Beauty in pp at 14 TeV with mult. interactions"); |
706 | gener = new AliGenPythia(-1); |
707 | gener->SetProcess(kPyBeautyppMNRwmi); |
708 | gener->SetStrucFunc(kCTEQ5L); |
709 | gener->SetPtHard(ptHardMin,ptHardMax); |
710 | gener->SetEnergyCMS(14000.); |
711 | break; |
712 | } |
713 | |
714 | return gener; |
715 | } |
716 | |
717 | AliGenerator* MbCocktail() |
718 | { |
719 | comment = comment.Append(" pp at 14 TeV: Pythia low-pt, no heavy quarks + J/Psi from parameterisation"); |
13217d97 |
720 | AliGenCocktail * gener = new AliGenCocktail(); |
e3d91248 |
721 | gener->UsePerEventRates(); |
e3d91248 |
722 | |
723 | // |
724 | // Pythia |
725 | AliGenPythia* pythia = new AliGenPythia(-1); |
726 | pythia->SetMomentumRange(0, 999999.); |
727 | pythia->SetThetaRange(0., 180.); |
728 | pythia->SetYRange(-12.,12.); |
729 | pythia->SetPtRange(0,1000.); |
730 | pythia->SetProcess(kPyMb); |
731 | pythia->SetEnergyCMS(14000.); |
732 | pythia->SwitchHFOff(); |
733 | |
734 | // |
735 | // J/Psi parameterisation |
736 | |
737 | AliGenParam* jpsi = new AliGenParam(1, AliGenMUONlib::kJpsi, "CDF scaled", "Jpsi"); |
738 | jpsi->SetPtRange(0.,100.); |
739 | jpsi->SetYRange(-8., 8.); |
740 | jpsi->SetPhiRange(0., 360.); |
741 | jpsi->SetForceDecay(kAll); |
742 | // |
743 | // |
744 | gener->AddGenerator(pythia, "Pythia", 1.); |
745 | gener->AddGenerator(jpsi, "J/Psi", 8.e-4); |
746 | |
747 | return gener; |
748 | } |
749 | |
750 | AliGenerator* PyMbTriggered(Int_t pdg) |
751 | { |
752 | AliGenPythia* pythia = new AliGenPythia(-1); |
753 | pythia->SetMomentumRange(0, 999999.); |
754 | pythia->SetThetaRange(0., 180.); |
755 | pythia->SetYRange(-12.,12.); |
756 | pythia->SetPtRange(0,1000.); |
757 | pythia->SetProcess(kPyMb); |
758 | pythia->SetEnergyCMS(14000.); |
759 | pythia->SetTriggerParticle(pdg, 0.9); |
760 | return pythia; |
761 | } |
762 | |
13217d97 |
763 | void ProcessEnvironmentVars() |
764 | { |
765 | // Run type |
766 | if (gSystem->Getenv("CONFIG_RUN_TYPE")) { |
767 | for (Int_t iRun = 0; iRun < kRunMax; iRun++) { |
768 | if (strcmp(gSystem->Getenv("CONFIG_RUN_TYPE"), pprRunName[iRun])==0) { |
769 | proc = (PDC06Proc_t)iRun; |
770 | cout<<"Run type set to "<<pprRunName[iRun]<<endl; |
771 | } |
772 | } |
773 | } |
774 | |
775 | // Random Number seed |
776 | if (gSystem->Getenv("CONFIG_SEED")) { |
777 | seed = atoi(gSystem->Getenv("CONFIG_SEED")); |
778 | } |
779 | } |
e3d91248 |
780 | |
781 | |
782 | |