]>
Commit | Line | Data |
---|---|---|
333b9f18 | 1 | // One can use the configuration macro in compiled mode by |
2 | // root [0] gSystem->Load("libgeant321"); | |
3 | // root [0] gSystem->SetIncludePath("-I$ROOTSYS/include -I$ALICE_ROOT/include\ | |
4 | // -I$ALICE_ROOT -I$ALICE/geant3/TGeant3"); | |
5 | // root [0] .x grun.C(1,"Config_PDC07.C++") | |
6 | ||
7 | #if !defined(__CINT__) || defined(__MAKECINT__) | |
8 | #include <Riostream.h> | |
9 | #include <TRandom.h> | |
10 | #include <TDatime.h> | |
11 | #include <TSystem.h> | |
12 | #include <TVirtualMC.h> | |
13 | #include <TGeant3TGeo.h> | |
14 | #include "EVGEN/AliGenCocktail.h" | |
15 | #include "EVGEN/AliGenParam.h" | |
16 | #include "EVGEN/AliGenMUONlib.h" | |
17 | #include "STEER/AliRunLoader.h" | |
18 | #include "STEER/AliRun.h" | |
19 | #include "STEER/AliConfig.h" | |
20 | #include "PYTHIA6/AliDecayerPythia.h" | |
21 | #include "PYTHIA6/AliGenPythia.h" | |
99c7d495 | 22 | #include "STEER/AliMagF.h" |
333b9f18 | 23 | #include "STRUCT/AliBODY.h" |
24 | #include "STRUCT/AliMAG.h" | |
25 | #include "STRUCT/AliABSOv3.h" | |
26 | #include "STRUCT/AliDIPOv3.h" | |
27 | #include "STRUCT/AliHALLv3.h" | |
28 | #include "STRUCT/AliFRAMEv2.h" | |
29 | #include "STRUCT/AliSHILv3.h" | |
30 | #include "STRUCT/AliPIPEv3.h" | |
31 | #include "ITS/AliITSv11Hybrid.h" | |
32 | #include "TPC/AliTPCv2.h" | |
33 | #include "TOF/AliTOFv6T0.h" | |
34 | #include "HMPID/AliHMPIDv3.h" | |
35 | #include "ZDC/AliZDCv3.h" | |
36 | #include "TRD/AliTRDv1.h" | |
37 | #include "TRD/AliTRDgeometry.h" | |
38 | #include "FMD/AliFMDv1.h" | |
39 | #include "MUON/AliMUONv1.h" | |
40 | #include "PHOS/AliPHOSv1.h" | |
41 | #include "PHOS/AliPHOSSimParam.h" | |
42 | #include "PMD/AliPMDv1.h" | |
43 | #include "T0/AliT0v1.h" | |
44 | #include "EMCAL/AliEMCALv2.h" | |
45 | #include "ACORDE/AliACORDEv1.h" | |
46 | #include "VZERO/AliVZEROv7.h" | |
47 | #endif | |
48 | ||
49 | ||
50 | enum PDC07Proc_t | |
51 | { | |
52 | //--- Heavy Flavour Production --- | |
53 | kCharmPbPb5500, kCharmpPb8800, kCharmpp14000, kCharmpp14000wmi, | |
54 | kD0PbPb5500, kD0pPb8800, kD0pp14000, | |
55 | kDPlusPbPb5500, kDPluspPb8800, kDPluspp14000, | |
56 | kBeautyPbPb5500, kBeautypPb8800, kBeautypp14000, kBeautypp14000wmi, | |
57 | // -- Pythia Mb | |
58 | kPyMbNoHvq, kPyOmegaPlus, kPyOmegaMinus, kPyJetJet, kBeautyppwmiJet,kPyBJetEMCAL, | |
59 | kPyGammaJetPHOS, kPyJetJetPHOS, kPyJetJetPHOSv2, kPyGammaBremsPHOS, | |
60 | kPyGammaJetEMCAL, kPyJetJetEMCAL, kPyGammaBremsEMCAL, kRunMax | |
61 | }; | |
62 | ||
63 | const char * pprRunName[] = { | |
64 | "kCharmPbPb5500", "kCharmpPb8800", "kCharmpp14000", "kCharmpp14000wmi", | |
65 | "kD0PbPb5500", "kD0pPb8800", "kD0pp14000", | |
66 | "kDPlusPbPb5500", "kDPluspPb8800", "kDPluspp14000", | |
67 | "kBeautyPbPb5500", "kBeautypPb8800", "kBeautypp14000", "kBeautypp14000wmi", | |
68 | "kPyMbNoHvq", "kPyOmegaPlus", "kPyOmegaMinus", "kPyJetJet", "kBeautyppwmiJet", "kPyBJetEMCAL", | |
69 | "kPyGammaJetPHOS", "kPyJetJetPHOS", "kPyJetJetPHOSv2", "kPyGammaBremsPHOS", | |
70 | "kPyGammaJetEMCAL", "kPyJetJetEMCAL", "kPyGammaBremsEMCAL" | |
71 | }; | |
72 | ||
73 | ||
74 | //--- Decay Mode --- | |
75 | enum DecayHvFl_t | |
76 | { | |
77 | kNature, kHadr, kSemiEl, kSemiMu | |
78 | }; | |
79 | //--- Magnetic Field --- | |
80 | enum Mag_t | |
81 | { | |
82 | kNoField, k5kG, kFieldMax | |
83 | }; | |
84 | //--- Trigger config --- | |
85 | enum TrigConf_t | |
86 | { | |
87 | kDefaultPPTrig, kDefaultPbPbTrig | |
88 | }; | |
89 | ||
90 | const char * TrigConfName[] = { | |
91 | "p-p","Pb-Pb" | |
92 | }; | |
93 | ||
94 | //--- Functions --- | |
95 | class AliGenPythia ; | |
96 | AliGenPythia *PythiaHVQ(PDC07Proc_t proc); | |
97 | AliGenPythia *PythiaHard(PDC07Proc_t proc); | |
98 | AliGenerator *MbCocktail(); | |
99 | AliGenerator *PyMbTriggered(Int_t pdg); | |
100 | void ProcessEnvironmentVars(); | |
101 | ||
102 | // This part for configuration | |
103 | static DecayHvFl_t decHvFl = kNature; | |
104 | static Mag_t mag = k5kG; | |
105 | static TrigConf_t trig = kDefaultPPTrig; // default pp trigger configuration | |
106 | static Int_t runNumber = 0; | |
107 | static Int_t eventNumber = 0; | |
108 | //========================// | |
109 | // Set Random Number seed // | |
110 | //========================// | |
111 | TDatime dt; | |
112 | static UInt_t seed = dt.Get(); | |
113 | ||
114 | // nEvts = -1 : you get 1 QQbar pair and all the fragmentation and | |
115 | // decay chain | |
116 | // nEvts = N>0 : you get N charm / beauty Hadrons | |
117 | Int_t nEvts = -1; | |
118 | // stars = kTRUE : all heavy resonances and their decay stored | |
119 | // = kFALSE: only final heavy hadrons and their decays stored | |
120 | Bool_t stars = kTRUE; | |
121 | ||
122 | // To be used only with kCharmpp14000wmi and kBeautypp14000wmi | |
123 | // To get a "reasonable" agreement with MNR results, events have to be | |
124 | // generated with the minimum ptHard set to 2.76 GeV. | |
125 | // To get a "perfect" agreement with MNR results, events have to be | |
126 | // generated in four ptHard bins with the following relative | |
127 | // normalizations: | |
128 | // CHARM | |
129 | // 2.76-3 GeV: 25% | |
130 | // 3-4 GeV: 40% | |
131 | // 4-8 GeV: 29% | |
132 | // >8 GeV: 6% | |
133 | // BEAUTY | |
134 | // 2.76-4 GeV: 5% | |
135 | // 4-6 GeV: 31% | |
136 | // 6-8 GeV: 28% | |
137 | // >8 GeV: 36% | |
138 | ||
139 | ||
140 | Float_t eCMS=14000; | |
141 | //PDC07Proc_t proc = kPyJetJet; | |
142 | //PDC07Proc_t proc = kPyGammaJetEMCAL; | |
143 | PDC07Proc_t proc = kPyBJetEMCAL; | |
144 | Int_t nPtHardGJBins = 11; | |
145 | Int_t ptHardGJLo[] = {5 ,10,20,30,40,50,60,70,80,90 ,100}; | |
146 | Int_t ptHardGJHi[] = {10,20,30,40,50,60,70,80,90,100,-1}; | |
147 | Int_t nPtHardJJBins = 17; | |
148 | Int_t ptHardJJLo[] = {2 ,12,16,20,24,29,35,41,50,60,72,86 ,104,124,149,179,215}; | |
149 | Int_t ptHardJJHi[] = {12,16,20,24,29,35,41,50,60,72,86,104,124,149,179,215,258}; | |
150 | Float_t ptHardMin = 10.; | |
151 | Float_t ptHardMax = 20.; | |
152 | Float_t ptGammaPi0Min = 0.; | |
153 | Int_t iquenching = 0; | |
154 | Float_t qhat = 20.; //for iquenching equal 1 and 4 | |
155 | //Float_t medlength = 6.; // for iquenching equal to 4 | |
156 | Int_t year = 0; //ALICE geometry configuration | |
157 | Bool_t pi0FragPhotonInDet = kFALSE ; | |
158 | ||
159 | // Comment line | |
160 | static TString comment; | |
161 | ||
162 | void Config() | |
163 | { | |
164 | ||
165 | // Get settings from environment variables | |
166 | ProcessEnvironmentVars(); | |
167 | ||
168 | gSystem->Load("liblhapdf.so"); // Parton density functions | |
169 | gSystem->Load("libEGPythia6.so"); // TGenerator interface | |
170 | if(iquenching !=4) gSystem->Load("libpythia6.so"); // Pythia | |
171 | else gSystem->Load("libqpythia.so"); // qpythia | |
172 | gSystem->Load("libAliPythia6.so"); // ALICE specific implementations | |
173 | ||
174 | // libraries required by geant321 | |
175 | #if defined(__CINT__) | |
176 | gSystem->Load("liblhapdf"); | |
177 | gSystem->Load("libEGPythia6"); | |
178 | if(iquenching !=4) gSystem->Load("libpythia6.so"); // Pythia | |
179 | else gSystem->Load("libqpythia.so"); // qpythia | |
180 | gSystem->Load("libAliPythia6"); | |
181 | gSystem->Load("libgeant321"); | |
182 | #endif | |
183 | ||
184 | new TGeant3TGeo("C++ Interface to Geant3"); | |
185 | ||
186 | // Output every 100 tracks | |
187 | ((TGeant3*)gMC)->SetSWIT(4,100); | |
188 | ||
189 | //======================================================================= | |
190 | ||
191 | // Run loader | |
192 | AliRunLoader* rl=0x0; | |
193 | rl = AliRunLoader::Open("galice.root", | |
194 | AliConfig::GetDefaultEventFolderName(), | |
195 | "recreate"); | |
196 | if (rl == 0x0) | |
197 | { | |
198 | gAlice->Fatal("Config.C","Can not instatiate the Run Loader"); | |
199 | return; | |
200 | } | |
201 | rl->SetCompressionLevel(2); | |
202 | rl->SetNumberOfEventsPerFile(1000); | |
203 | gAlice->SetRunLoader(rl); | |
204 | ||
205 | // Set the trigger configuration | |
28da60d3 | 206 | AliSimulation::Instance()->SetTriggerConfig(TrigConfName[trig]); |
333b9f18 | 207 | cout<<"Trigger configuration is set to "<<TrigConfName[trig]<<endl; |
208 | ||
209 | // | |
210 | //======================================================================= | |
211 | // ************* STEERING parameters FOR ALICE SIMULATION ************** | |
212 | // --- Specify event type to be tracked through the ALICE setup | |
213 | // --- All positions are in cm, angles in degrees, and P and E in GeV | |
214 | ||
215 | ||
216 | gMC->SetProcess("DCAY",1); | |
217 | gMC->SetProcess("PAIR",1); | |
218 | gMC->SetProcess("COMP",1); | |
219 | gMC->SetProcess("PHOT",1); | |
220 | gMC->SetProcess("PFIS",0); | |
221 | gMC->SetProcess("DRAY",0); | |
222 | gMC->SetProcess("ANNI",1); | |
223 | gMC->SetProcess("BREM",1); | |
224 | gMC->SetProcess("MUNU",1); | |
225 | gMC->SetProcess("CKOV",1); | |
226 | gMC->SetProcess("HADR",1); | |
227 | gMC->SetProcess("LOSS",2); | |
228 | gMC->SetProcess("MULS",1); | |
229 | gMC->SetProcess("RAYL",1); | |
230 | ||
231 | Float_t cut = 1.e-3; // 1MeV cut by default | |
232 | Float_t tofmax = 1.e10; | |
233 | ||
234 | gMC->SetCut("CUTGAM", cut); | |
235 | gMC->SetCut("CUTELE", cut); | |
236 | gMC->SetCut("CUTNEU", cut); | |
237 | gMC->SetCut("CUTHAD", cut); | |
238 | gMC->SetCut("CUTMUO", cut); | |
239 | gMC->SetCut("BCUTE", cut); | |
240 | gMC->SetCut("BCUTM", cut); | |
241 | gMC->SetCut("DCUTE", cut); | |
242 | gMC->SetCut("DCUTM", cut); | |
243 | gMC->SetCut("PPCUTM", cut); | |
244 | gMC->SetCut("TOFMAX", tofmax); | |
245 | ||
246 | //======================// | |
247 | // Set External decayer // | |
248 | //======================// | |
249 | TVirtualMCDecayer* decayer = new AliDecayerPythia(); | |
250 | // DECAYS | |
251 | // | |
252 | switch(decHvFl) { | |
253 | case kNature: | |
254 | decayer->SetForceDecay(kNeutralPion); | |
255 | break; | |
256 | case kHadr: | |
257 | decayer->SetForceDecay(kHadronicD); | |
258 | break; | |
259 | case kSemiEl: | |
260 | decayer->SetForceDecay(kSemiElectronic); | |
261 | break; | |
262 | case kSemiMu: | |
263 | decayer->SetForceDecay(kSemiMuonic); | |
264 | break; | |
265 | } | |
266 | decayer->Init(); | |
267 | gMC->SetExternalDecayer(decayer); | |
268 | //if(proc == kPyJetJetPHOSv2) // in this case we decay the pi0 | |
269 | // decayer->SetForceDecay(kNeutralPion); | |
270 | ||
271 | //=========================// | |
272 | // Generator Configuration // | |
273 | //=========================// | |
274 | AliGenPythia* gener = 0x0; | |
275 | ||
276 | if (proc <= kBeautypp14000wmi) { | |
277 | AliGenPythia *pythia = PythiaHVQ(proc); | |
278 | // FeedDown option | |
279 | pythia->SetFeedDownHigherFamily(kFALSE); | |
280 | // Stack filling option | |
281 | if(!stars) pythia->SetStackFillOpt(AliGenPythia::kParentSelection); | |
282 | // Set Count mode | |
283 | if(nEvts>0) pythia->SetCountMode(AliGenPythia::kCountParents); | |
284 | // | |
285 | // DECAYS | |
286 | // | |
287 | switch(decHvFl) { | |
288 | case kNature: | |
289 | pythia->SetForceDecay(kNeutralPion); | |
290 | break; | |
291 | case kHadr: | |
292 | pythia->SetForceDecay(kHadronicD); | |
293 | break; | |
294 | case kSemiEl: | |
295 | pythia->SetForceDecay(kSemiElectronic); | |
296 | break; | |
297 | case kSemiMu: | |
298 | pythia->SetForceDecay(kSemiMuonic); | |
299 | break; | |
300 | } | |
301 | // | |
302 | // GEOM & KINE CUTS | |
303 | // | |
304 | pythia->SetMomentumRange(0,99999999); | |
305 | pythia->SetPhiRange(0., 360.); | |
306 | pythia->SetThetaRange(0,180); | |
307 | switch(ycut) { | |
308 | case kFull: | |
309 | pythia->SetYRange(-12,12); | |
310 | break; | |
311 | case kBarrel: | |
312 | pythia->SetYRange(-2,2); | |
313 | break; | |
314 | case kMuonArm: | |
315 | pythia->SetYRange(1,6); | |
316 | break; | |
317 | } | |
318 | gener = pythia; | |
319 | } else if (proc == kPyMbNoHvq) { | |
320 | gener = MbCocktail(); | |
321 | } else if (proc == kPyOmegaMinus) { | |
322 | gener = PyMbTriggered(3334); | |
323 | } else if (proc == kPyOmegaPlus) { | |
324 | gener = PyMbTriggered(-3334); | |
325 | } else if (proc <= kPyGammaBremsEMCAL) { | |
326 | AliGenPythia *pythia = PythiaHard(proc); | |
327 | // FeedDown option | |
328 | pythia->SetFeedDownHigherFamily(kFALSE); | |
329 | // Set Count mode | |
330 | if(nEvts>0) pythia->SetCountMode(AliGenPythia::kCountParents); | |
331 | ||
332 | // | |
333 | // GEOM & KINE CUTS | |
334 | // | |
335 | pythia->SetMomentumRange(0,99999999); | |
336 | // if(proc == kPyJetJetPHOSv2) | |
337 | // pythia->SetForceDecay(kNeutralPion); | |
338 | // else | |
339 | // pythia->SetForceDecay(kAll); | |
340 | ||
341 | pythia->SetForceDecay(kNeutralPion); | |
342 | pythia->SetPycellParameters(2., 274, 432, 0., 4., 5., 1.0); | |
343 | pythia->SetGluonRadiation(1,1); //ISR and FRS effects on. | |
344 | pythia->SetPtKick(3.); // set the intrinsic kt to about 5.3 GeV/c | |
345 | gener = pythia; | |
346 | } | |
347 | ||
348 | //gener->SetTrackingFlag(1); | |
349 | ||
350 | //Quenching | |
351 | gener->SetQuench(iquenching); | |
352 | Float_t k = 6e5*(qhat/1.7) ; //qhat=1.7, k = 6e5, default value | |
353 | if(iquenching == 1)//quenching weights | |
354 | AliPythia::Instance()->InitQuenching(0.,0.1,k,0,0.95,6); | |
355 | else if(iquenching == 4)//q-pythia | |
356 | gener->SetQhat(k); // !transport coefficient (GeV2/fm) | |
357 | ||
358 | ||
359 | // PRIMARY VERTEX | |
360 | ||
361 | gener->SetOrigin(0., 0., 0.); // vertex position | |
362 | ||
363 | // Size of the interaction diamond | |
364 | // Longitudinal | |
365 | Float_t sigmaz = 7.55 / TMath::Sqrt(2.); // [cm] 14 TeV (p-p), 5.5 TeV (Pb-Pb) | |
366 | if(eCMS == 10000) sigmaz = 5.4 / TMath::Sqrt(2.); // [cm] 10 TeV | |
367 | if(eCMS == 900) sigmaz = 10.5 / TMath::Sqrt(2.); // [cm] 0.9 TeV | |
368 | cout<<">>> SigmaZ x srqt(2.) "<< sigmaz * TMath::Sqrt(2.) << endl; | |
369 | ||
370 | // Transverse | |
371 | Float_t betast = 10; // beta* [m] | |
372 | Float_t eps = 3.75e-6; // emittance [m] | |
373 | Float_t gamma = eCMS / 2.0 / 0.938272; // relativistic gamma [1] | |
374 | Float_t sigmaxy = TMath::Sqrt(eps * betast / gamma) / TMath::Sqrt(2.) * 100.; // [cm] | |
375 | printf("\n \n Diamond size x-y: %10.3e z: %10.3e\n \n", sigmaxy, sigmaz); | |
376 | ||
377 | gener->SetSigma(sigmaxy, sigmaxy, sigmaz); // Sigma in (X,Y,Z) (cm) on IP position | |
378 | gener->SetCutVertexZ(3.); // Truncate at 3 sigma | |
379 | gener->SetVertexSmear(kPerEvent); | |
380 | ||
381 | gener->Init(); | |
382 | ||
383 | ||
384 | // FIELD | |
385 | // | |
99c7d495 | 386 | |
387 | AliMagF* field = 0x0; | |
333b9f18 | 388 | if (mag == kNoField) { |
389 | comment = comment.Append(" | L3 field 0.0 T"); | |
4642ac4b | 390 | field = new AliMagF("Maps","Maps", 0., 0., AliMagF::k5kGUniform); |
333b9f18 | 391 | } else if (mag == k5kG) { |
392 | comment = comment.Append(" | L3 field 0.5 T"); | |
4642ac4b | 393 | field = new AliMagF("Maps","Maps", -1., -1., AliMagF::k5kG); |
333b9f18 | 394 | } |
395 | printf("\n \n Comment: %s \n \n", comment.Data()); | |
99c7d495 | 396 | TGeoGlobalMagField::Instance()->SetField(field); |
397 | ||
333b9f18 | 398 | rl->CdGAFile(); |
399 | gAlice->SetField(field); | |
400 | ||
401 | ||
402 | Int_t iABSO = 1; | |
403 | Int_t iACORDE = 0; | |
404 | Int_t iDIPO = 1; | |
405 | Int_t iEMCAL = 1; | |
406 | Int_t iFMD = 1; | |
407 | Int_t iFRAME = 1; | |
408 | Int_t iHALL = 1; | |
409 | Int_t iITS = 1; | |
410 | Int_t iMAG = 1; | |
411 | Int_t iMUON = 1; | |
412 | Int_t iPHOS = 1; | |
413 | Int_t iPIPE = 1; | |
414 | Int_t iPMD = 1; | |
415 | Int_t iHMPID = 1; | |
416 | Int_t iSHIL = 1; | |
417 | Int_t iT0 = 1; | |
418 | Int_t iTOF = 1; | |
419 | Int_t iTPC = 1; | |
420 | Int_t iTRD = 1; | |
421 | Int_t iVZERO = 1; | |
422 | Int_t iZDC = 1; | |
423 | ||
424 | ||
425 | //=================== Alice BODY parameters ============================= | |
426 | AliBODY *BODY = new AliBODY("BODY", "Alice envelop"); | |
427 | ||
428 | ||
429 | if (iMAG) | |
430 | { | |
431 | //=================== MAG parameters ============================ | |
432 | // --- Start with Magnet since detector layouts may be depending --- | |
433 | // --- on the selected Magnet dimensions --- | |
434 | AliMAG *MAG = new AliMAG("MAG", "Magnet"); | |
435 | } | |
436 | ||
437 | ||
438 | if (iABSO) | |
439 | { | |
440 | //=================== ABSO parameters ============================ | |
441 | AliABSO *ABSO = new AliABSOv3("ABSO", "Muon Absorber"); | |
442 | } | |
443 | ||
444 | if (iDIPO) | |
445 | { | |
446 | //=================== DIPO parameters ============================ | |
447 | ||
448 | AliDIPO *DIPO = new AliDIPOv3("DIPO", "Dipole version 3"); | |
449 | } | |
450 | ||
451 | if (iHALL) | |
452 | { | |
453 | //=================== HALL parameters ============================ | |
454 | ||
455 | AliHALL *HALL = new AliHALLv3("HALL", "Alice Hall"); | |
456 | } | |
457 | ||
458 | ||
459 | if (iFRAME) | |
460 | { | |
461 | //=================== FRAME parameters ============================ | |
462 | ||
463 | AliFRAMEv2 *FRAME = new AliFRAMEv2("FRAME", "Space Frame"); | |
464 | FRAME->SetHoles(1); | |
465 | } | |
466 | ||
467 | if (iSHIL) | |
468 | { | |
469 | //=================== SHIL parameters ============================ | |
470 | ||
471 | AliSHIL *SHIL = new AliSHILv3("SHIL", "Shielding Version 3"); | |
472 | } | |
473 | ||
474 | ||
475 | if (iPIPE) | |
476 | { | |
477 | //=================== PIPE parameters ============================ | |
478 | ||
479 | AliPIPE *PIPE = new AliPIPEv3("PIPE", "Beam Pipe"); | |
480 | } | |
481 | ||
482 | if (iITS) | |
483 | { | |
484 | //=================== ITS parameters ============================ | |
485 | ||
486 | AliITS *ITS = new AliITSv11Hybrid("ITS","ITS v11Hybrid"); | |
487 | } | |
488 | ||
489 | if (iTPC) | |
490 | { | |
491 | //============================ TPC parameters ===================== | |
492 | ||
493 | AliTPC *TPC = new AliTPCv2("TPC", "Default"); | |
494 | } | |
495 | ||
496 | ||
497 | if (iTOF) { | |
498 | //=================== TOF parameters ============================ | |
499 | ||
500 | AliTOF *TOF = new AliTOFv6T0("TOF", "normal TOF"); | |
501 | } | |
502 | ||
503 | ||
504 | if (iHMPID) | |
505 | { | |
506 | //=================== HMPID parameters =========================== | |
507 | ||
508 | AliHMPID *HMPID = new AliHMPIDv3("HMPID", "normal HMPID"); | |
509 | } | |
510 | ||
511 | ||
512 | if (iZDC) | |
513 | { | |
514 | //=================== ZDC parameters ============================ | |
515 | ||
516 | AliZDC *ZDC = new AliZDCv3("ZDC", "normal ZDC"); | |
517 | } | |
518 | ||
519 | if (iTRD) | |
520 | { | |
521 | //=================== TRD parameters ============================ | |
522 | ||
523 | AliTRD *TRD = new AliTRDv1("TRD", "TRD slow simulator"); | |
524 | AliTRDgeometry *geoTRD = TRD->GetGeometry(); | |
525 | if(year == 2009){ | |
526 | // Partial geometry: modules at 0,8,9,17 | |
527 | // Partial geometry: modules at 1,7,10,16 expected for 2009 | |
528 | // starting at 3h in positive direction | |
529 | ||
530 | //geoTRD->SetSMstatus(1,0); | |
531 | //geoTRD->SetSMstatus(2,0); | |
532 | //geoTRD->SetSMstatus(3,0); | |
533 | //geoTRD->SetSMstatus(4,0); | |
534 | //geoTRD->SetSMstatus(5,0); | |
535 | //geoTRD->SetSMstatus(6,0); | |
536 | //geoTRD->SetSMstatus(7,0); | |
537 | ||
538 | ||
539 | //geoTRD->SetSMstatus(10,0); | |
540 | //geoTRD->SetSMstatus(11,0); | |
541 | //geoTRD->SetSMstatus(12,0); | |
542 | //geoTRD->SetSMstatus(13,0); | |
543 | //geoTRD->SetSMstatus(14,0); | |
544 | //geoTRD->SetSMstatus(15,0); | |
545 | //geoTRD->SetSMstatus(16,0); | |
546 | ||
547 | geoTRD->SetSMstatus(2,0); | |
548 | geoTRD->SetSMstatus(3,0); | |
549 | geoTRD->SetSMstatus(4,0); | |
550 | geoTRD->SetSMstatus(5,0); | |
551 | geoTRD->SetSMstatus(6,0); | |
552 | ||
553 | geoTRD->SetSMstatus(11,0); | |
554 | geoTRD->SetSMstatus(12,0); | |
555 | geoTRD->SetSMstatus(13,0); | |
556 | geoTRD->SetSMstatus(14,0); | |
557 | geoTRD->SetSMstatus(15,0); | |
558 | } | |
559 | } | |
560 | ||
561 | if (iFMD) | |
562 | { | |
563 | //=================== FMD parameters ============================ | |
564 | AliFMD *FMD = new AliFMDv1("FMD", "normal FMD"); | |
565 | } | |
566 | ||
567 | if (iMUON) | |
568 | { | |
569 | //=================== MUON parameters =========================== | |
570 | // New MUONv1 version (geometry defined via builders) | |
571 | AliMUON *MUON = new AliMUONv1("MUON", "default"); | |
572 | } | |
573 | //=================== PHOS parameters =========================== | |
574 | ||
575 | if (iPHOS) | |
576 | { | |
577 | AliPHOS *PHOS = new AliPHOSv1("PHOS", "IHEP"); | |
578 | } | |
579 | ||
580 | ||
581 | if (iPMD) | |
582 | { | |
583 | //=================== PMD parameters ============================ | |
584 | AliPMD *PMD = new AliPMDv1("PMD", "normal PMD"); | |
585 | } | |
586 | ||
587 | if (iT0) | |
588 | { | |
589 | //=================== T0 parameters ============================ | |
590 | AliT0 *T0 = new AliT0v1("T0", "T0 Detector"); | |
591 | } | |
592 | ||
593 | if (iEMCAL) | |
594 | { | |
595 | //=================== EMCAL parameters ============================ | |
596 | AliEMCAL *EMCAL = new AliEMCALv2("EMCAL", "EMCAL_COMPLETE"); | |
597 | } | |
598 | ||
599 | if (iACORDE) | |
600 | { | |
601 | //=================== ACORDE parameters ============================ | |
602 | AliACORDE *ACORDE = new AliACORDEv1("ACORDE", "normal ACORDE"); | |
603 | } | |
604 | ||
605 | if (iVZERO) | |
606 | { | |
607 | //=================== VZERO parameters ============================ | |
608 | AliVZERO *VZERO = new AliVZEROv7("VZERO", "normal VZERO"); | |
609 | } | |
610 | ||
611 | AliLog::Message(AliLog::kInfo, "End of Config", "Config.C", "Config.C", "Config()"," Config.C", __LINE__); | |
612 | } | |
613 | ||
614 | // PYTHIA | |
615 | ||
616 | AliGenPythia *PythiaHVQ(PDC07Proc_t proc) { | |
617 | //*******************************************************************// | |
618 | // Configuration file for charm / beauty generation with PYTHIA // | |
619 | // // | |
620 | // The parameters have been tuned in order to reproduce the inclusive// | |
621 | // heavy quark pt distribution given by the NLO pQCD calculation by // | |
622 | // Mangano, Nason and Ridolfi. // | |
623 | // // | |
624 | // For details and for the NORMALIZATION of the yields see: // | |
625 | // N.Carrer and A.Dainese, // | |
626 | // "Charm and beauty production at the LHC", // | |
627 | // ALICE-INT-2003-019, [arXiv:hep-ph/0311225]; // | |
628 | // PPR Chapter 6.6, CERN/LHCC 2005-030 (2005). // | |
629 | //*******************************************************************// | |
630 | AliGenPythia * gener = 0x0; | |
631 | ||
632 | switch(proc) { | |
633 | case kCharmPbPb5500: | |
634 | comment = comment.Append(" Charm in Pb-Pb at 5.5 TeV"); | |
635 | gener = new AliGenPythia(nEvts); | |
636 | gener->SetProcess(kPyCharmPbPbMNR); | |
637 | gener->SetStrucFunc(kCTEQ4L); | |
638 | gener->SetPtHard(2.1,-1.0); | |
639 | gener->SetEnergyCMS(5500.); | |
640 | gener->SetNuclei(208,208); | |
641 | break; | |
642 | case kCharmpPb8800: | |
643 | comment = comment.Append(" Charm in p-Pb at 8.8 TeV"); | |
644 | gener = new AliGenPythia(nEvts); | |
645 | gener->SetProcess(kPyCharmpPbMNR); | |
646 | gener->SetStrucFunc(kCTEQ4L); | |
647 | gener->SetPtHard(2.1,-1.0); | |
648 | gener->SetEnergyCMS(8800.); | |
649 | gener->SetProjectile("P",1,1); | |
650 | gener->SetTarget("Pb",208,82); | |
651 | break; | |
652 | case kCharmpp14000: | |
653 | comment = comment.Append(" Charm in pp at 14 TeV"); | |
654 | gener = new AliGenPythia(nEvts); | |
655 | gener->SetProcess(kPyCharmppMNR); | |
656 | gener->SetStrucFunc(kCTEQ4L); | |
657 | gener->SetPtHard(2.1,-1.0); | |
658 | gener->SetEnergyCMS(14000.); | |
659 | break; | |
660 | case kCharmpp14000wmi: | |
661 | comment = comment.Append(" Charm in pp at 14 TeV with mult. interactions"); | |
662 | gener = new AliGenPythia(-1); | |
663 | gener->SetProcess(kPyCharmppMNRwmi); | |
664 | gener->SetStrucFunc(kCTEQ5L); | |
665 | gener->SetPtHard(ptHardMin,ptHardMax); | |
666 | gener->SetEnergyCMS(14000.); | |
667 | break; | |
668 | case kD0PbPb5500: | |
669 | comment = comment.Append(" D0 in Pb-Pb at 5.5 TeV"); | |
670 | gener = new AliGenPythia(nEvts); | |
671 | gener->SetProcess(kPyD0PbPbMNR); | |
672 | gener->SetStrucFunc(kCTEQ4L); | |
673 | gener->SetPtHard(2.1,-1.0); | |
674 | gener->SetEnergyCMS(5500.); | |
675 | gener->SetNuclei(208,208); | |
676 | break; | |
677 | case kD0pPb8800: | |
678 | comment = comment.Append(" D0 in p-Pb at 8.8 TeV"); | |
679 | gener = new AliGenPythia(nEvts); | |
680 | gener->SetProcess(kPyD0pPbMNR); | |
681 | gener->SetStrucFunc(kCTEQ4L); | |
682 | gener->SetPtHard(2.1,-1.0); | |
683 | gener->SetEnergyCMS(8800.); | |
684 | gener->SetProjectile("P",1,1); | |
685 | gener->SetTarget("Pb",208,82); | |
686 | break; | |
687 | case kD0pp14000: | |
688 | comment = comment.Append(" D0 in pp at 14 TeV"); | |
689 | gener = new AliGenPythia(nEvts); | |
690 | gener->SetProcess(kPyD0ppMNR); | |
691 | gener->SetStrucFunc(kCTEQ4L); | |
692 | gener->SetPtHard(2.1,-1.0); | |
693 | gener->SetEnergyCMS(14000.); | |
694 | break; | |
695 | case kDPlusPbPb5500: | |
696 | comment = comment.Append(" DPlus in Pb-Pb at 5.5 TeV"); | |
697 | gener = new AliGenPythia(nEvts); | |
698 | gener->SetProcess(kPyDPlusPbPbMNR); | |
699 | gener->SetStrucFunc(kCTEQ4L); | |
700 | gener->SetPtHard(2.1,-1.0); | |
701 | gener->SetEnergyCMS(5500.); | |
702 | gener->SetNuclei(208,208); | |
703 | break; | |
704 | case kDPluspPb8800: | |
705 | comment = comment.Append(" DPlus in p-Pb at 8.8 TeV"); | |
706 | gener = new AliGenPythia(nEvts); | |
707 | gener->SetProcess(kPyDPluspPbMNR); | |
708 | gener->SetStrucFunc(kCTEQ4L); | |
709 | gener->SetPtHard(2.1,-1.0); | |
710 | gener->SetEnergyCMS(8800.); | |
711 | gener->SetProjectile("P",1,1); | |
712 | gener->SetTarget("Pb",208,82); | |
713 | break; | |
714 | case kDPluspp14000: | |
715 | comment = comment.Append(" DPlus in pp at 14 TeV"); | |
716 | gener = new AliGenPythia(nEvts); | |
717 | gener->SetProcess(kPyDPlusppMNR); | |
718 | gener->SetStrucFunc(kCTEQ4L); | |
719 | gener->SetPtHard(2.1,-1.0); | |
720 | gener->SetEnergyCMS(14000.); | |
721 | break; | |
722 | case kBeautyPbPb5500: | |
723 | comment = comment.Append(" Beauty in Pb-Pb at 5.5 TeV"); | |
724 | gener = new AliGenPythia(nEvts); | |
725 | gener->SetProcess(kPyBeautyPbPbMNR); | |
726 | gener->SetStrucFunc(kCTEQ4L); | |
727 | gener->SetPtHard(2.75,-1.0); | |
728 | gener->SetEnergyCMS(5500.); | |
729 | gener->SetNuclei(208,208); | |
730 | break; | |
731 | case kBeautypPb8800: | |
732 | comment = comment.Append(" Beauty in p-Pb at 8.8 TeV"); | |
733 | gener = new AliGenPythia(nEvts); | |
734 | gener->SetProcess(kPyBeautypPbMNR); | |
735 | gener->SetStrucFunc(kCTEQ4L); | |
736 | gener->SetPtHard(2.75,-1.0); | |
737 | gener->SetEnergyCMS(8800.); | |
738 | gener->SetProjectile("P",1,1); | |
739 | gener->SetTarget("Pb",208,82); | |
740 | break; | |
741 | case kBeautypp14000: | |
742 | comment = comment.Append(" Beauty in pp at 14 TeV"); | |
743 | gener = new AliGenPythia(nEvts); | |
744 | gener->SetProcess(kPyBeautyppMNR); | |
745 | gener->SetStrucFunc(kCTEQ4L); | |
746 | gener->SetPtHard(2.75,-1.0); | |
747 | gener->SetEnergyCMS(14000.); | |
748 | break; | |
749 | case kBeautypp14000wmi: | |
750 | comment = comment.Append(" Beauty in pp at 14 TeV with mult. interactions"); | |
751 | gener = new AliGenPythia(-1); | |
752 | gener->SetProcess(kPyBeautyppMNRwmi); | |
753 | gener->SetStrucFunc(kCTEQ5L); | |
754 | gener->SetPtHard(ptHardMin,ptHardMax); | |
755 | gener->SetEnergyCMS(14000.); | |
756 | break; | |
757 | } | |
758 | ||
759 | return gener; | |
760 | } | |
761 | ||
762 | AliGenPythia *PythiaHard(PDC07Proc_t proc) { | |
763 | //*******************************************************************// | |
764 | // Configuration file for hard QCD processes generation with PYTHIA // | |
765 | // // | |
766 | //*******************************************************************// | |
767 | AliGenPythia * gener = 0x0; | |
768 | ||
769 | switch(proc) { | |
770 | ||
771 | case kPyJetJet: | |
772 | comment = comment.Append(" pp->jet + jet over at 14 TeV, no restriction"); | |
773 | AliGenPythia * gener = new AliGenPythia(nEvts); | |
774 | gener->SetEnergyCMS(eCMS);// Centre of mass energy | |
775 | gener->SetProcess(kPyJets);// Process type | |
776 | gener->SetJetEtaRange(-1.5, 1.5);// Final state kinematic cuts | |
777 | gener->SetJetPhiRange(0., 360.); | |
778 | gener->SetJetEtRange(10., 1000.); | |
779 | gener->SetPtHard(ptHardMin, ptHardMax);// Pt transfer of the hard scattering | |
780 | gener->SetStrucFunc(kCTEQ4L); | |
781 | break; | |
782 | case kBeautyppwmiJet: | |
783 | comment = comment.Append(" Beauty Jets, no acceptance cuts"); | |
784 | AliGenPythia *gener = new AliGenPythia(-1); | |
785 | gener->SetProcess(kPyBeautyppMNRwmi); | |
786 | gener->SetStrucFunc(kCTEQ4L); | |
787 | //Jet specific stuff below | |
788 | gener->SetPtHard(ptHardMin,ptHardMax); | |
789 | gener->SetJetEtaRange(-1,1); //used like that in LHC08d10, does nothing? | |
790 | gener->SetJetPhiRange(60.,210.);//used like that in LHC08d10, does nothing? | |
791 | gener->SetJetEtRange(10.,1000.);//used like that in LHC08d10,does nothing | |
792 | //jet specific stuff above | |
793 | gener->SetEnergyCMS(eCMS); | |
794 | gener->SetForceDecay(kSemiElectronic); | |
795 | gener->SetYRange(-2,2); | |
796 | gener->SetFeedDownHigherFamily(kFALSE); | |
797 | gener->SetCountMode(AliGenPythia::kCountParents); | |
798 | break; | |
799 | case kPyBJetEMCAL: | |
800 | comment = comment.Append(" Beauty Jets over EMCAL"); | |
801 | AliGenPythia *gener = new AliGenPythia(-1); | |
802 | gener->SetProcess(kPyBeautyJets); | |
803 | gener->SetStrucFunc(kCTEQ4L); | |
804 | //Jet specific stuff below | |
805 | gener->SetPtHard(ptHardMin,ptHardMax); | |
806 | gener->SetJetEtaRange(-1,1); | |
807 | gener->SetJetPhiRange(60.,210.); | |
808 | gener->SetJetEtRange(10.,1000.); | |
809 | //jet specific stuff above | |
810 | gener->SetEnergyCMS(eCMS); | |
811 | gener->SetForceDecay(kSemiElectronic); | |
812 | gener->SetYRange(-2,2); | |
813 | gener->SetFeedDownHigherFamily(kFALSE); | |
814 | gener->SetCountMode(AliGenPythia::kCountParents); | |
815 | gener->SetElectronMinPt(1.); | |
816 | gener->SetElectronInEMCAL(kTRUE); | |
817 | break; | |
818 | case kPyGammaJetPHOS: | |
819 | comment = comment.Append(" pp->jet + gamma over PHOS"); | |
820 | gener = new AliGenPythia(nEvts); | |
821 | gener->SetEnergyCMS(eCMS); | |
822 | gener->SetProcess(kPyDirectGamma); | |
823 | gener->SetStrucFunc(kCTEQ4L); | |
824 | gener->SetPtHard(ptHardMin,ptHardMax); | |
825 | //gener->SetYHard(-1.0,1.0); | |
826 | gener->SetGammaEtaRange(-0.13,0.13); | |
827 | gener->SetGammaPhiRange(219.,321.);//Over 5 modules +-2 deg | |
828 | break; | |
829 | case kPyJetJetPHOS: | |
830 | comment = comment.Append(" pp->jet + jet over PHOS"); | |
831 | gener = new AliGenPythia(nEvts); | |
832 | gener->SetEnergyCMS(eCMS); | |
833 | gener->SetProcess(kPyJets); | |
834 | gener->SetStrucFunc(kCTEQ4L); | |
835 | gener->SetPtHard(ptHardMin,ptHardMax); | |
836 | //gener->SetYHard(-1.0,1.0); | |
837 | gener->SetJetEtaRange(-1.,1.); | |
838 | gener->SetJetPhiRange(200.,340.); | |
839 | gener->SetJetEtRange(10., 1000.); | |
840 | gener->SetPi0InPHOS(pi0FragPhotonInDet); | |
841 | gener->SetFragPhotonOrPi0MinPt(ptGammaPi0Min); | |
842 | ||
843 | printf("\n \n Event generator: Minimum pT of particle in calorimeter %f \n \n", ptGammaPi0Min); | |
844 | break; | |
845 | case kPyGammaBremsPHOS: | |
846 | comment = comment.Append(" pp->jet + jet+bremsphoton over PHOS at 14 TeV"); | |
847 | gener = new AliGenPythia(nEvts); | |
848 | gener->SetEnergyCMS(eCMS); | |
849 | gener->SetProcess(kPyJets); | |
850 | gener->SetStrucFunc(kCTEQ4L); | |
851 | gener->SetPtHard(ptHardMin,ptHardMax); | |
852 | //gener->SetYHard(-1.0,1.0); | |
853 | gener->SetJetEtaRange(-1.,1.); | |
854 | gener->SetJetPhiRange(200.,340.); | |
855 | gener->SetFragPhotonInPHOS(pi0FragPhotonInDet); | |
856 | gener->SetFragPhotonOrPi0MinPt(ptGammaPi0Min); | |
857 | printf("\n \n Event generator: Minimum pT of particle in calorimeter %f \n \n", ptGammaPi0Min); | |
858 | break; | |
859 | case kPyJetJetPHOSv2: | |
860 | comment = comment.Append(" pp->jet + jet over PHOS version2 "); | |
861 | gener = new AliGenPythia(nEvts); | |
862 | gener->SetEnergyCMS(eCMS); | |
863 | gener->SetProcess(kPyJets); | |
864 | gener->SetStrucFunc(kCTEQ4L); | |
865 | gener->SetPtHard(ptHardMin,ptHardMax); | |
866 | //gener->SetYHard(-1.0,1.0); | |
867 | gener->SetJetEtaRange(-1.,1.); | |
868 | gener->SetJetPhiRange(200.,340.); | |
869 | //gener->SetPi0InPHOS(kTRUE); | |
870 | gener->SetPhotonInPHOSeta(pi0FragPhotonInDet); | |
871 | gener->SetPhotonMinPt(ptGammaPi0Min); | |
872 | gener->SetForceDecay(kAll); | |
873 | break; | |
874 | case kPyGammaJetEMCAL: | |
875 | comment = comment.Append(" pp->jet + gamma over EMCAL at 14 TeV"); | |
876 | gener = new AliGenPythia(nEvts); | |
877 | gener->SetEnergyCMS(eCMS); | |
878 | gener->SetProcess(kPyDirectGamma); | |
879 | gener->SetStrucFunc(kCTEQ4L); | |
880 | gener->SetPtHard(ptHardMin,ptHardMax); | |
881 | //gener->SetYHard(-1.0,1.0); | |
882 | gener->SetGammaEtaRange(-0.701,0.701); | |
883 | gener->SetGammaPhiRange(79.,191.);//Over 6 supermodules +-2 deg | |
884 | break; | |
885 | case kPyJetJetEMCAL: | |
886 | comment = comment.Append(" pp->jet + jet over EMCAL at 14 TeV"); | |
887 | gener = new AliGenPythia(nEvts); | |
888 | gener->SetEnergyCMS(eCMS); | |
889 | gener->SetProcess(kPyJets); | |
890 | gener->SetStrucFunc(kCTEQ4L); | |
891 | gener->SetPtHard(ptHardMin,ptHardMax); | |
892 | //gener->SetYHard(-1.0,1.0); | |
893 | gener->SetJetEtaRange(-1,1); | |
894 | gener->SetJetPhiRange(60.,210.); | |
895 | gener->SetJetEtRange(10., 1000.); | |
896 | gener->SetPi0InEMCAL(pi0FragPhotonInDet); | |
897 | gener->SetFragPhotonOrPi0MinPt(ptGammaPi0Min); | |
898 | printf("\n \n Event generator: Minimum pT of particle in calorimeter %f \n \n", ptGammaPi0Min); | |
899 | break; | |
900 | case kPyGammaBremsEMCAL: | |
901 | comment = comment.Append(" pp->jet + jet+bremsphoton over EMCAL at 14 TeV"); | |
902 | gener = new AliGenPythia(nEvts); | |
903 | gener->SetEnergyCMS(eCMS); | |
904 | gener->SetProcess(kPyJets); | |
905 | gener->SetStrucFunc(kCTEQ4L); | |
906 | gener->SetPtHard(ptHardMin,ptHardMax); | |
907 | //gener->SetYHard(-1.0,1.0); | |
908 | gener->SetJetEtaRange(-1,1); | |
909 | gener->SetJetPhiRange(60.,210.); | |
910 | gener->SetJetEtRange(10., 1000.); | |
911 | gener->SetFragPhotonInEMCAL(pi0FragPhotonInDet); | |
912 | gener->SetFragPhotonOrPi0MinPt(ptGammaPi0Min); | |
913 | ||
914 | printf("\n \n Event generator: Minimum pT of particle in calorimeter %f \n \n", ptGammaPi0Min); | |
915 | break; | |
916 | ||
917 | } | |
918 | ||
919 | return gener; | |
920 | } | |
921 | ||
922 | AliGenerator* MbCocktail() | |
923 | { | |
924 | comment = comment.Append(" pp at 14 TeV: Pythia low-pt, no heavy quarks + J/Psi from parameterisation"); | |
925 | AliGenCocktail * gener = new AliGenCocktail(); | |
926 | gener->UsePerEventRates(); | |
927 | ||
928 | // Pythia | |
929 | ||
930 | AliGenPythia* pythia = new AliGenPythia(-1); | |
931 | pythia->SetMomentumRange(0, 999999.); | |
932 | pythia->SetThetaRange(0., 180.); | |
933 | pythia->SetYRange(-12.,12.); | |
934 | pythia->SetPtRange(0,1000.); | |
935 | pythia->SetProcess(kPyMb); | |
936 | pythia->SetEnergyCMS(14000.); | |
937 | pythia->SwitchHFOff(); | |
938 | ||
939 | // J/Psi parameterisation | |
940 | ||
941 | AliGenParam* jpsi = new AliGenParam(1, AliGenMUONlib::kJpsi, "CDF scaled", "Jpsi"); | |
942 | jpsi->SetPtRange(0.,100.); | |
943 | jpsi->SetYRange(-8., 8.); | |
944 | jpsi->SetPhiRange(0., 360.); | |
945 | jpsi->SetForceDecay(kAll); | |
946 | ||
947 | gener->AddGenerator(pythia, "Pythia", 1.); | |
948 | gener->AddGenerator(jpsi, "J/Psi", 8.e-4); | |
949 | ||
950 | return gener; | |
951 | } | |
952 | ||
953 | AliGenerator* PyMbTriggered(Int_t pdg) | |
954 | { | |
955 | AliGenPythia* pythia = new AliGenPythia(-1); | |
956 | pythia->SetMomentumRange(0, 999999.); | |
957 | pythia->SetThetaRange(0., 180.); | |
958 | pythia->SetYRange(-12.,12.); | |
959 | pythia->SetPtRange(0,1000.); | |
960 | pythia->SetProcess(kPyMb); | |
961 | pythia->SetEnergyCMS(14000.); | |
962 | pythia->SetTriggerParticle(pdg, 0.9); | |
963 | return pythia; | |
964 | } | |
965 | ||
966 | void ProcessEnvironmentVars() | |
967 | { | |
968 | cout << "######################################" << endl; | |
969 | cout << "## Processing environment variables ##" << endl; | |
970 | cout << "######################################" << endl; | |
971 | ||
972 | // Run Number | |
973 | if (gSystem->Getenv("RUN")) { | |
974 | runNumber = atoi(gSystem->Getenv("RUN")); | |
975 | } | |
976 | //cout<<"Run number "<<runNumber<<endl; | |
977 | ||
978 | // Event Number | |
979 | if (gSystem->Getenv("EVENT")) { | |
980 | eventNumber = atoi(gSystem->Getenv("EVENT")); | |
981 | } | |
982 | //cout<<"Event number "<<eventNumber<<endl; | |
983 | ||
984 | // Random Number seed | |
985 | if (gSystem->Getenv("CONFIG_SEED")) { | |
986 | seed = atoi(gSystem->Getenv("CONFIG_SEED")); | |
987 | } | |
988 | else if(gSystem->Getenv("EVENT") && gSystem->Getenv("RUN")){ | |
989 | seed = runNumber * 100000 + eventNumber; | |
990 | } | |
991 | ||
992 | gRandom->SetSeed(seed); | |
993 | ||
994 | cout<<"////////////////////////////////////////////////////////////////////////////////////"<<endl; | |
995 | cout<<"Seed for random number generation= "<< seed<<" "<< gRandom->GetSeed()<<endl; | |
996 | cout<<"////////////////////////////////////////////////////////////////////////////////////"<<endl; | |
997 | ||
998 | // Run type | |
999 | if (gSystem->Getenv("DC_RUN_TYPE")) { | |
1000 | cout<<"run type "<<gSystem->Getenv("DC_RUN_TYPE")<<endl; | |
1001 | for (Int_t iRun = 0; iRun < kRunMax; iRun++) { | |
1002 | if (strcmp(gSystem->Getenv("DC_RUN_TYPE"), pprRunName[iRun])==0) { | |
1003 | proc = (PDC07Proc_t)iRun; | |
1004 | } | |
1005 | } | |
1006 | } | |
1007 | else | |
1008 | cout << "Environment variable DC_RUN_TYPE is not defined" << endl; | |
1009 | if (gSystem->Getenv("YEAR")) | |
1010 | year = atoi(gSystem->Getenv("YEAR")); | |
1011 | if (gSystem->Getenv("ECMS")) | |
1012 | eCMS = atof(gSystem->Getenv("ECMS")); | |
1013 | if (gSystem->Getenv("TRIGGER")) | |
1014 | trig = atoi(gSystem->Getenv("TRIGGER")); | |
1015 | if ( gSystem->Getenv("PI0GAMMAINDET") ) | |
1016 | pi0FragPhotonInDet = atoi(gSystem->Getenv("PI0GAMMAINDET")); | |
1017 | if (gSystem->Getenv("PTGAMMAPI0MIN")) | |
1018 | ptGammaPi0Min = atof(gSystem->Getenv("PTGAMMAPI0MIN")); | |
1019 | if (gSystem->Getenv("QUENCHING")) | |
1020 | iquenching = atof(gSystem->Getenv("QUENCHING")); | |
1021 | if (gSystem->Getenv("QHAT")) | |
1022 | qhat = atof(gSystem->Getenv("QHAT")); | |
1023 | //if (gSystem->Getenv("MEDLENGHT")) | |
1024 | // medlength = atof(gSystem->Getenv("MEDLENGHT")); | |
1025 | if (gSystem->Getenv("PTHARDMIN")) | |
1026 | ptHardMin = atof(gSystem->Getenv("PTHARDMIN")); | |
1027 | if (gSystem->Getenv("PTHARDMAX")) | |
1028 | ptHardMax = atof(gSystem->Getenv("PTHARDMAX")); | |
1029 | ||
1030 | if(gSystem->Getenv("PTHARDMIN") && gSystem->Getenv("PTHARDMAX")){ | |
1031 | if( strncmp(pprRunName[proc],"kPyJetJet",9)==0 || strncmp(pprRunName[proc],"kPyBJetEMCAL",11)==0 ){ | |
1032 | ptHardMin = ptHardJJLo[runNumber%nPtHardJJBins]; | |
1033 | ptHardMax = ptHardJJHi[runNumber%nPtHardJJBins]; | |
1034 | } | |
1035 | else if( strncmp(pprRunName[proc],"kPyGammaJet",11)==0){ | |
1036 | ptHardMin = ptHardGJLo[runNumber%nPtHardGJBins]; | |
1037 | ptHardMax = ptHardGJHi[runNumber%nPtHardGJBins]; | |
1038 | } | |
1039 | } | |
1040 | ||
1041 | cout<<">> Run type set to "<<pprRunName[proc]<<endl; | |
1042 | cout<<">> Collision energy set to "<<eCMS <<endl; | |
1043 | cout<<">> Collisions trigger type "<<trig<<" "<<TrigConfName[trig]<<endl; | |
1044 | cout<<">> ptHard limits: "<<ptHardMin<<" to " <<ptHardMax<<" GeV"<<endl; | |
1045 | if(pi0FragPhotonInDet) cout<<">> pt gamma/pi0 threshold "<< ptGammaPi0Min<<" GeV "<<endl; | |
1046 | if(iquenching) cout<<">> quenching on? "<< iquenching<<" qhat "<<qhat//<<" medlength "<<medlength | |
1047 | <<endl; | |
1048 | cout<<">> Year geometry : "<<year<<endl; | |
1049 | cout << "######################################" << endl; | |
1050 | cout << "######################################" << endl; | |
1051 | } |