]>
Commit | Line | Data |
---|---|---|
1 | // | |
2 | // Configuration for the Physics Data Challenge 2006 | |
3 | // | |
4 | ||
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 | ||
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> | |
17 | #include <TGeant3TGeo.h> | |
18 | #include "EVGEN/AliGenCocktail.h" | |
19 | #include "EVGEN/AliGenParam.h" | |
20 | #include "EVGEN/AliGenMUONlib.h" | |
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" | |
35 | #include "ITS/AliITSgeom.h" | |
36 | #include "ITS/AliITSvPPRasymmFMD.h" | |
37 | #include "TPC/AliTPCv2.h" | |
38 | #include "TOF/AliTOFv6T0.h" | |
39 | #include "HMPID/AliHMPIDv1.h" | |
40 | #include "ZDC/AliZDCv2.h" | |
41 | #include "TRD/AliTRDv1.h" | |
42 | #include "FMD/AliFMDv1.h" | |
43 | #include "MUON/AliMUONv1.h" | |
44 | #include "PHOS/AliPHOSv1.h" | |
45 | #include "PMD/AliPMDv1.h" | |
46 | #include "T0/AliT0v1.h" | |
47 | #include "EMCAL/AliEMCALv2.h" | |
48 | #include "ACORDE/AliACORDEv0.h" | |
49 | #include "VZERO/AliVZEROv7.h" | |
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 | |
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" | |
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 | }; | |
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 | ||
99 | //--- Functions --- | |
100 | AliGenPythia *PythiaHVQ(PDC06Proc_t proc); | |
101 | AliGenerator *MbCocktail(); | |
102 | AliGenerator *PyMbTriggered(Int_t pdg); | |
103 | void ProcessEnvironmentVars(); | |
104 | ||
105 | // This part for configuration | |
106 | static PDC06Proc_t proc = kPyMbNoHvq; | |
107 | static DecayHvFl_t decHvFl = kNature; | |
108 | static YCut_t ycut = kFull; | |
109 | static Mag_t mag = k5kG; | |
110 | static TrigConf_t trig = kDefaultPPTrig; // default pp trigger configuration | |
111 | //========================// | |
112 | // Set Random Number seed // | |
113 | //========================// | |
114 | TDatime dt; | |
115 | static UInt_t seed = dt.Get(); | |
116 | ||
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 | ||
151 | ||
152 | // Get settings from environment variables | |
153 | ProcessEnvironmentVars(); | |
154 | ||
155 | gRandom->SetSeed(seed); | |
156 | cerr<<"Seed for random number generation= "<<seed<<endl; | |
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); | |
183 | ||
184 | // Set the trigger configuration | |
185 | gAlice->SetTriggerDescriptor(TrigConfName[trig]); | |
186 | cout<<"Trigger configuration is set to "<<TrigConfName[trig]<<endl; | |
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 | //=========================// | |
253 | AliGenerator* gener = 0x0; | |
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 | |
325 | gener->SetCutVertexZ(3.); // Truncate at 3 sigma | |
326 | gener->SetVertexSmear(kPerEvent); | |
327 | ||
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; | |
349 | Int_t iACORDE = 0; | |
350 | Int_t iDIPO = 1; | |
351 | Int_t iEMCAL = 1; | |
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; | |
358 | Int_t iPHOS = 1; | |
359 | Int_t iPIPE = 1; | |
360 | Int_t iPMD = 1; | |
361 | Int_t iHMPID = 1; | |
362 | Int_t iSHIL = 1; | |
363 | Int_t iT0 = 1; | |
364 | Int_t iTOF = 1; | |
365 | Int_t iTPC = 1; | |
366 | Int_t iTRD = 1; | |
367 | Int_t iVZERO = 1; | |
368 | Int_t iZDC = 1; | |
369 | ||
370 | ||
371 | //=================== Alice BODY parameters ============================= | |
372 | AliBODY *BODY = new AliBODY("BODY", "Alice envelop"); | |
373 | ||
374 | ||
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 | } | |
382 | ||
383 | ||
384 | if (iABSO) | |
385 | { | |
386 | //=================== ABSO parameters ============================ | |
387 | AliABSO *ABSO = new AliABSOv0("ABSO", "Muon Absorber"); | |
388 | } | |
389 | ||
390 | if (iDIPO) | |
391 | { | |
392 | //=================== DIPO parameters ============================ | |
393 | ||
394 | AliDIPO *DIPO = new AliDIPOv2("DIPO", "Dipole version 2"); | |
395 | } | |
396 | ||
397 | if (iHALL) | |
398 | { | |
399 | //=================== HALL parameters ============================ | |
400 | ||
401 | AliHALL *HALL = new AliHALL("HALL", "Alice Hall"); | |
402 | } | |
403 | ||
404 | ||
405 | if (iFRAME) | |
406 | { | |
407 | //=================== FRAME parameters ============================ | |
408 | ||
409 | AliFRAMEv2 *FRAME = new AliFRAMEv2("FRAME", "Space Frame"); | |
410 | } | |
411 | ||
412 | if (iSHIL) | |
413 | { | |
414 | //=================== SHIL parameters ============================ | |
415 | ||
416 | AliSHIL *SHIL = new AliSHILv2("SHIL", "Shielding Version 2"); | |
417 | } | |
418 | ||
419 | ||
420 | if (iPIPE) | |
421 | { | |
422 | //=================== PIPE parameters ============================ | |
423 | ||
424 | AliPIPE *PIPE = new AliPIPEv0("PIPE", "Beam Pipe"); | |
425 | } | |
426 | ||
427 | if(iITS) { | |
428 | ||
429 | //=================== ITS parameters ============================ | |
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 | // | |
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(kFALSE); // 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(150.); // chip thickness on layer 1 must be in the range [150,300] | |
454 | ITS->SetThicknessChip2(150.); // 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 | ||
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"); | |
463 | //ITS->SetRails(0); // 1 --> rails in ; 0 --> rails out | |
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"); | |
467 | //ITS->SetRails(0); // 1 --> rails in ; 0 --> rails out | |
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 | // | |
480 | ITS->SetEUCLID(0); | |
481 | } | |
482 | ||
483 | if (iTPC) | |
484 | { | |
485 | //============================ TPC parameters ===================== | |
486 | AliTPC *TPC = new AliTPCv2("TPC", "Default"); | |
487 | } | |
488 | ||
489 | ||
490 | if (iTOF) { | |
491 | //=================== TOF parameters ============================ | |
492 | AliTOF *TOF = new AliTOFv6T0("TOF", "normal TOF"); | |
493 | // Partial geometry: modules at 2,3,4,6,7,11,12,14,15,16 | |
494 | // starting at 6h in positive direction | |
495 | // Int_t TOFSectors[18]={-1,-1,0,0,0,-1,0,0,-1,-1,-1,0,0,-1,0,0,0,0}; | |
496 | // Partial geometry: modules at 1,2,6,7,9,10,11,12,15,16,17 | |
497 | // (ALICE numbering convention) | |
498 | Int_t TOFSectors[18]={-1,0,0,-1,-1,-1,0,0,-1,0,0,0,0,-1,-1,0,0,0}; | |
499 | TOF->SetTOFSectors(TOFSectors); | |
500 | } | |
501 | ||
502 | ||
503 | if (iHMPID) | |
504 | { | |
505 | //=================== HMPID parameters =========================== | |
506 | AliHMPID *HMPID = new AliHMPIDv1("HMPID", "normal HMPID"); | |
507 | ||
508 | } | |
509 | ||
510 | ||
511 | if (iZDC) | |
512 | { | |
513 | //=================== ZDC parameters ============================ | |
514 | ||
515 | AliZDC *ZDC = new AliZDCv2("ZDC", "normal ZDC"); | |
516 | } | |
517 | ||
518 | if (iTRD) | |
519 | { | |
520 | //=================== TRD parameters ============================ | |
521 | ||
522 | AliTRD *TRD = new AliTRDv1("TRD", "TRD slow simulator"); | |
523 | AliTRDgeometry *geoTRD = TRD->GetGeometry(); | |
524 | // Partial geometry: modules at 2,3,4,6,11,12,14,15 | |
525 | // starting at 6h in positive direction | |
526 | geoTRD->SetSMstatus(0,0); | |
527 | geoTRD->SetSMstatus(1,0); | |
528 | geoTRD->SetSMstatus(5,0); | |
529 | geoTRD->SetSMstatus(7,0); | |
530 | geoTRD->SetSMstatus(8,0); | |
531 | geoTRD->SetSMstatus(9,0); | |
532 | geoTRD->SetSMstatus(10,0); | |
533 | geoTRD->SetSMstatus(13,0); | |
534 | geoTRD->SetSMstatus(16,0); | |
535 | geoTRD->SetSMstatus(17,0); | |
536 | } | |
537 | ||
538 | if (iFMD) | |
539 | { | |
540 | //=================== FMD parameters ============================ | |
541 | AliFMD *FMD = new AliFMDv1("FMD", "normal FMD"); | |
542 | } | |
543 | ||
544 | if (iMUON) | |
545 | { | |
546 | //=================== MUON parameters =========================== | |
547 | // New MUONv1 version (geometry defined via builders) | |
548 | AliMUON *MUON = new AliMUONv1("MUON", "default"); | |
549 | } | |
550 | //=================== PHOS parameters =========================== | |
551 | ||
552 | if (iPHOS) | |
553 | { | |
554 | AliPHOS *PHOS = new AliPHOSv1("PHOS", "IHEP"); | |
555 | } | |
556 | ||
557 | ||
558 | if (iPMD) | |
559 | { | |
560 | //=================== PMD parameters ============================ | |
561 | AliPMD *PMD = new AliPMDv1("PMD", "normal PMD"); | |
562 | } | |
563 | ||
564 | if (iT0) | |
565 | { | |
566 | //=================== T0 parameters ============================ | |
567 | AliT0 *T0 = new AliT0v1("T0", "T0 Detector"); | |
568 | } | |
569 | ||
570 | if (iEMCAL) | |
571 | { | |
572 | //=================== EMCAL parameters ============================ | |
573 | AliEMCAL *EMCAL = new AliEMCALv2("EMCAL", "SHISH_77_TRD1_2X2_FINAL_110DEG"); | |
574 | } | |
575 | ||
576 | if (iACORDE) | |
577 | { | |
578 | //=================== ACORDE parameters ============================ | |
579 | AliACORDE *ACORDE = new AliACORDEv0("ACORDE", "normal ACORDE"); | |
580 | } | |
581 | ||
582 | if (iVZERO) | |
583 | { | |
584 | //=================== VZERO parameters ============================ | |
585 | AliVZERO *VZERO = new AliVZEROv7("VZERO", "normal VZERO"); | |
586 | } | |
587 | } | |
588 | // | |
589 | // PYTHIA | |
590 | // | |
591 | AliGenPythia *PythiaHVQ(PDC06Proc_t proc) { | |
592 | //*******************************************************************// | |
593 | // Configuration file for charm / beauty generation with PYTHIA // | |
594 | // // | |
595 | // The parameters have been tuned in order to reproduce the inclusive// | |
596 | // heavy quark pt distribution given by the NLO pQCD calculation by // | |
597 | // Mangano, Nason and Ridolfi. // | |
598 | // // | |
599 | // For details and for the NORMALIZATION of the yields see: // | |
600 | // N.Carrer and A.Dainese, // | |
601 | // "Charm and beauty production at the LHC", // | |
602 | // ALICE-INT-2003-019, [arXiv:hep-ph/0311225]; // | |
603 | // PPR Chapter 6.6, CERN/LHCC 2005-030 (2005). // | |
604 | //*******************************************************************// | |
605 | AliGenPythia * gener = 0x0; | |
606 | ||
607 | switch(proc) { | |
608 | case kCharmPbPb5500: | |
609 | comment = comment.Append(" Charm in Pb-Pb at 5.5 TeV"); | |
610 | gener = new AliGenPythia(nEvts); | |
611 | gener->SetProcess(kPyCharmPbPbMNR); | |
612 | gener->SetStrucFunc(kCTEQ4L); | |
613 | gener->SetPtHard(2.1,-1.0); | |
614 | gener->SetEnergyCMS(5500.); | |
615 | gener->SetNuclei(208,208); | |
616 | break; | |
617 | case kCharmpPb8800: | |
618 | comment = comment.Append(" Charm in p-Pb at 8.8 TeV"); | |
619 | gener = new AliGenPythia(nEvts); | |
620 | gener->SetProcess(kPyCharmpPbMNR); | |
621 | gener->SetStrucFunc(kCTEQ4L); | |
622 | gener->SetPtHard(2.1,-1.0); | |
623 | gener->SetEnergyCMS(8800.); | |
624 | gener->SetProjectile("P",1,1); | |
625 | gener->SetTarget("Pb",208,82); | |
626 | break; | |
627 | case kCharmpp14000: | |
628 | comment = comment.Append(" Charm in pp at 14 TeV"); | |
629 | gener = new AliGenPythia(nEvts); | |
630 | gener->SetProcess(kPyCharmppMNR); | |
631 | gener->SetStrucFunc(kCTEQ4L); | |
632 | gener->SetPtHard(2.1,-1.0); | |
633 | gener->SetEnergyCMS(14000.); | |
634 | break; | |
635 | case kCharmpp14000wmi: | |
636 | comment = comment.Append(" Charm in pp at 14 TeV with mult. interactions"); | |
637 | gener = new AliGenPythia(-1); | |
638 | gener->SetProcess(kPyCharmppMNRwmi); | |
639 | gener->SetStrucFunc(kCTEQ5L); | |
640 | gener->SetPtHard(ptHardMin,ptHardMax); | |
641 | gener->SetEnergyCMS(14000.); | |
642 | break; | |
643 | case kD0PbPb5500: | |
644 | comment = comment.Append(" D0 in Pb-Pb at 5.5 TeV"); | |
645 | gener = new AliGenPythia(nEvts); | |
646 | gener->SetProcess(kPyD0PbPbMNR); | |
647 | gener->SetStrucFunc(kCTEQ4L); | |
648 | gener->SetPtHard(2.1,-1.0); | |
649 | gener->SetEnergyCMS(5500.); | |
650 | gener->SetNuclei(208,208); | |
651 | break; | |
652 | case kD0pPb8800: | |
653 | comment = comment.Append(" D0 in p-Pb at 8.8 TeV"); | |
654 | gener = new AliGenPythia(nEvts); | |
655 | gener->SetProcess(kPyD0pPbMNR); | |
656 | gener->SetStrucFunc(kCTEQ4L); | |
657 | gener->SetPtHard(2.1,-1.0); | |
658 | gener->SetEnergyCMS(8800.); | |
659 | gener->SetProjectile("P",1,1); | |
660 | gener->SetTarget("Pb",208,82); | |
661 | break; | |
662 | case kD0pp14000: | |
663 | comment = comment.Append(" D0 in pp at 14 TeV"); | |
664 | gener = new AliGenPythia(nEvts); | |
665 | gener->SetProcess(kPyD0ppMNR); | |
666 | gener->SetStrucFunc(kCTEQ4L); | |
667 | gener->SetPtHard(2.1,-1.0); | |
668 | gener->SetEnergyCMS(14000.); | |
669 | break; | |
670 | case kDPlusPbPb5500: | |
671 | comment = comment.Append(" DPlus in Pb-Pb at 5.5 TeV"); | |
672 | gener = new AliGenPythia(nEvts); | |
673 | gener->SetProcess(kPyDPlusPbPbMNR); | |
674 | gener->SetStrucFunc(kCTEQ4L); | |
675 | gener->SetPtHard(2.1,-1.0); | |
676 | gener->SetEnergyCMS(5500.); | |
677 | gener->SetNuclei(208,208); | |
678 | break; | |
679 | case kDPluspPb8800: | |
680 | comment = comment.Append(" DPlus in p-Pb at 8.8 TeV"); | |
681 | gener = new AliGenPythia(nEvts); | |
682 | gener->SetProcess(kPyDPluspPbMNR); | |
683 | gener->SetStrucFunc(kCTEQ4L); | |
684 | gener->SetPtHard(2.1,-1.0); | |
685 | gener->SetEnergyCMS(8800.); | |
686 | gener->SetProjectile("P",1,1); | |
687 | gener->SetTarget("Pb",208,82); | |
688 | break; | |
689 | case kDPluspp14000: | |
690 | comment = comment.Append(" DPlus in pp at 14 TeV"); | |
691 | gener = new AliGenPythia(nEvts); | |
692 | gener->SetProcess(kPyDPlusppMNR); | |
693 | gener->SetStrucFunc(kCTEQ4L); | |
694 | gener->SetPtHard(2.1,-1.0); | |
695 | gener->SetEnergyCMS(14000.); | |
696 | break; | |
697 | case kBeautyPbPb5500: | |
698 | comment = comment.Append(" Beauty in Pb-Pb at 5.5 TeV"); | |
699 | gener = new AliGenPythia(nEvts); | |
700 | gener->SetProcess(kPyBeautyPbPbMNR); | |
701 | gener->SetStrucFunc(kCTEQ4L); | |
702 | gener->SetPtHard(2.75,-1.0); | |
703 | gener->SetEnergyCMS(5500.); | |
704 | gener->SetNuclei(208,208); | |
705 | break; | |
706 | case kBeautypPb8800: | |
707 | comment = comment.Append(" Beauty in p-Pb at 8.8 TeV"); | |
708 | gener = new AliGenPythia(nEvts); | |
709 | gener->SetProcess(kPyBeautypPbMNR); | |
710 | gener->SetStrucFunc(kCTEQ4L); | |
711 | gener->SetPtHard(2.75,-1.0); | |
712 | gener->SetEnergyCMS(8800.); | |
713 | gener->SetProjectile("P",1,1); | |
714 | gener->SetTarget("Pb",208,82); | |
715 | break; | |
716 | case kBeautypp14000: | |
717 | comment = comment.Append(" Beauty in pp at 14 TeV"); | |
718 | gener = new AliGenPythia(nEvts); | |
719 | gener->SetProcess(kPyBeautyppMNR); | |
720 | gener->SetStrucFunc(kCTEQ4L); | |
721 | gener->SetPtHard(2.75,-1.0); | |
722 | gener->SetEnergyCMS(14000.); | |
723 | break; | |
724 | case kBeautypp14000wmi: | |
725 | comment = comment.Append(" Beauty in pp at 14 TeV with mult. interactions"); | |
726 | gener = new AliGenPythia(-1); | |
727 | gener->SetProcess(kPyBeautyppMNRwmi); | |
728 | gener->SetStrucFunc(kCTEQ5L); | |
729 | gener->SetPtHard(ptHardMin,ptHardMax); | |
730 | gener->SetEnergyCMS(14000.); | |
731 | break; | |
732 | } | |
733 | ||
734 | return gener; | |
735 | } | |
736 | ||
737 | AliGenerator* MbCocktail() | |
738 | { | |
739 | comment = comment.Append(" pp at 14 TeV: Pythia low-pt, no heavy quarks + J/Psi from parameterisation"); | |
740 | AliGenCocktail * gener = new AliGenCocktail(); | |
741 | gener->UsePerEventRates(); | |
742 | ||
743 | // | |
744 | // Pythia | |
745 | AliGenPythia* pythia = new AliGenPythia(-1); | |
746 | pythia->SetMomentumRange(0, 999999.); | |
747 | pythia->SetThetaRange(0., 180.); | |
748 | pythia->SetYRange(-12.,12.); | |
749 | pythia->SetPtRange(0,1000.); | |
750 | pythia->SetProcess(kPyMb); | |
751 | pythia->SetEnergyCMS(14000.); | |
752 | pythia->SwitchHFOff(); | |
753 | ||
754 | // | |
755 | // J/Psi parameterisation | |
756 | ||
757 | AliGenParam* jpsi = new AliGenParam(1, AliGenMUONlib::kJpsi, "CDF scaled", "Jpsi"); | |
758 | jpsi->SetPtRange(0.,100.); | |
759 | jpsi->SetYRange(-8., 8.); | |
760 | jpsi->SetPhiRange(0., 360.); | |
761 | jpsi->SetForceDecay(kAll); | |
762 | // | |
763 | // | |
764 | gener->AddGenerator(jpsi, "J/Psi", 8.e-4); | |
765 | gener->AddGenerator(pythia, "Pythia", 1.); | |
766 | ||
767 | ||
768 | return gener; | |
769 | } | |
770 | ||
771 | AliGenerator* PyMbTriggered(Int_t pdg) | |
772 | { | |
773 | AliGenPythia* pythia = new AliGenPythia(-1); | |
774 | pythia->SetMomentumRange(0, 999999.); | |
775 | pythia->SetThetaRange(0., 180.); | |
776 | pythia->SetYRange(-12.,12.); | |
777 | pythia->SetPtRange(0,1000.); | |
778 | pythia->SetProcess(kPyMb); | |
779 | pythia->SetEnergyCMS(14000.); | |
780 | pythia->SetTriggerParticle(pdg, 0.9); | |
781 | return pythia; | |
782 | } | |
783 | ||
784 | void ProcessEnvironmentVars() | |
785 | { | |
786 | // Run type | |
787 | if (gSystem->Getenv("CONFIG_RUN_TYPE")) { | |
788 | for (Int_t iRun = 0; iRun < kRunMax; iRun++) { | |
789 | if (strcmp(gSystem->Getenv("CONFIG_RUN_TYPE"), pprRunName[iRun])==0) { | |
790 | proc = (PDC06Proc_t)iRun; | |
791 | cout<<"Run type set to "<<pprRunName[iRun]<<endl; | |
792 | } | |
793 | } | |
794 | } | |
795 | ||
796 | // Random Number seed | |
797 | if (gSystem->Getenv("CONFIG_SEED")) { | |
798 | seed = atoi(gSystem->Getenv("CONFIG_SEED")); | |
799 | } | |
800 | } | |
801 | ||
802 | ||
803 |