]> git.uio.no Git - u/mrichter/AliRoot.git/blame - macros/ConfigRaw2015.C
CMake: removing qpythia from the depedencies
[u/mrichter/AliRoot.git] / macros / ConfigRaw2015.C
CommitLineData
d080385b 1//
2// Configuration file for generating 2015 geometry
3//
4
5// To use this macro in compiled mode one first would need to do:
6// root [0] gSystem->Load("libgeant321");
7// root [1] gSystem->SetIncludePath("-I$ROOTSYS/include -I$ALICE_ROOT/include -Igeant3_src/TGeant3");
8//
9// However it is meant to be called in the following sequence:
10// root [0] AliCDBManager* cdb = AliCDBManager::Instance();
11// root [0] cdb->SetDefaultStorage("local://$ALICE_ROOT/OCDB");
12// root [0] cdb->SetRun(200000);
13// root [0] gROOT->LoadMacro("ConfigRaw2015.C");
14// root [0] gInterpreter->ProcessLine(gAlice->GetConfigFunction());
15// root [0] gAlice->GetMCApp()->Init();
16
17#if !defined(__CINT__) || defined(__MAKECINT__)
18#include <Riostream.h>
19#include <TRandom.h>
20#include <TDatime.h>
21#include <TSystem.h>
22#include <TVirtualMC.h>
23#include <TGeoGlobalMagField.h>
24#include <TGeant3TGeo.h>
25#include "STEER/STEER/AliRunLoader.h"
26#include "STEER/STEER/AliRun.h"
27#include "STEER/STEER/AliConfig.h"
28#include "STEER/STEER/AliSimulation.h"
29#include "PYTHIA6/AliPythia6/AliDecayerPythia.h"
30#include "PYTHIA6/AliPythia6/AliGenPythia.h"
31#include "TDPMjet/AliGenDPMjet.h"
32#include "STEER/STEERBase/AliMagF.h"
33#include "STRUCT/AliBODY.h"
34#include "STRUCT/AliMAG.h"
35#include "STRUCT/AliABSOv3.h"
36#include "STRUCT/AliDIPOv3.h"
37#include "STRUCT/AliHALLv3.h"
38#include "STRUCT/AliFRAMEv2.h"
39#include "STRUCT/AliSHILv3.h"
40#include "STRUCT/AliPIPEv3.h"
41#include "ITS/ITSsim/AliITSv11.h"
42#include "TPC/TPCsim/AliTPCv2.h"
43#include "TOF/TOFsim/AliTOFv6T0.h"
44#include "HMPID/HMPIDsim/AliHMPIDv3.h"
45#include "ZDC/ZDCsim/AliZDCv4.h"
46#include "TRD/TRDsim/AliTRDv1.h"
47#include "TRD/TRDbase/AliTRDgeometry.h"
48#include "FMD/FMDsim/AliFMDv1.h"
49#include "MUON/MUONsim/AliMUONv1.h"
50#include "PHOS/PHOSsim/AliPHOSv1.h"
51#include "PHOS/PHOSbase/AliPHOSSimParam.h"
52#include "PMD/PMDsim/AliPMDv1.h"
53#include "T0/T0sim/AliT0v1.h"
54#include "EMCAL/EMCALsim/AliEMCALv2.h"
55#include "ACORDE/ACORDEsim/AliACORDEv1.h"
56#include "VZERO/VZEROsim/AliVZEROv7.h"
57#include "AD/ADsim/AliADv1.h"
58#endif
59
60enum PDC06Proc_t
61{
62 kPythia6, kPythia6D6T, kPhojet, kRunMax
63};
64
65const char * pprRunName[] = {
66 "kPythia6", "kPythia6D6T", "kPhojet"
67};
68
69enum Mag_t
70{
71 kNoField, k5kG, kFieldMax
72};
73
74const char * pprField[] = {
75 "kNoField", "k5kG"
76};
77
78//--- Functions ---
79class AliGenPythia;
80AliGenerator *MbPythia();
81AliGenerator *MbPythiaTuneD6T();
82AliGenerator *MbPhojet();
83void ProcessEnvironmentVars();
84
85// Geterator, field, beam energy
86static PDC06Proc_t proc = kPhojet;
87static Mag_t mag = k5kG;
88static Float_t energy = 10000; // energy in CMS
89//========================//
90// Set Random Number seed //
91//========================//
92TDatime dt;
93static UInt_t seed = dt.Get();
94
95// Comment line
96static TString comment;
97
98void Config()
99{
100
101
102 // Get settings from environment variables
103 ProcessEnvironmentVars();
104
105 gRandom->SetSeed(seed);
106 cerr<<"Seed for random number generation= "<<seed<<endl;
107
108 // Libraries required by geant321
109#if defined(__CINT__)
110 gSystem->Load("liblhapdf"); // Parton density functions
111 gSystem->Load("libEGPythia6"); // TGenerator interface
112 if (proc != kPythia6D6T) {
113 gSystem->Load("libpythia6"); // Pythia 6.2
114 } else {
115 gSystem->Load("libqpythia"); // Pythia 6.4
116 }
117 gSystem->Load("libAliPythia6"); // ALICE specific implementations
118 gSystem->Load("libgeant321");
119#endif
120
121 new TGeant3TGeo("C++ Interface to Geant3");
122
123 //=======================================================================
124 // Create the output file
125
126
127 AliRunLoader* rl=0x0;
128
129 cout<<"Config.C: Creating Run Loader ..."<<endl;
130 rl = AliRunLoader::Open("galice.root",
131 AliConfig::GetDefaultEventFolderName(),
132 "recreate");
133 if (rl == 0x0)
134 {
135 gAlice->Fatal("Config.C","Can not instatiate the Run Loader");
136 return;
137 }
138 rl->SetCompressionLevel(2);
139 rl->SetNumberOfEventsPerFile(1000);
140 gAlice->SetRunLoader(rl);
141 // gAlice->SetGeometryFromFile("geometry.root");
142 // gAlice->SetGeometryFromCDB();
143
144 // Set the trigger configuration: proton-proton
145 AliSimulation::Instance()->SetTriggerConfig("p-p");
146
147 //
148 //=======================================================================
149 // ************* STEERING parameters FOR ALICE SIMULATION **************
150 // --- Specify event type to be tracked through the ALICE setup
151 // --- All positions are in cm, angles in degrees, and P and E in GeV
152
153
154 gMC->SetProcess("DCAY",1);
155 gMC->SetProcess("PAIR",1);
156 gMC->SetProcess("COMP",1);
157 gMC->SetProcess("PHOT",1);
158 gMC->SetProcess("PFIS",0);
159 gMC->SetProcess("DRAY",0);
160 gMC->SetProcess("ANNI",1);
161 gMC->SetProcess("BREM",1);
162 gMC->SetProcess("MUNU",1);
163 gMC->SetProcess("CKOV",1);
164 gMC->SetProcess("HADR",1);
165 gMC->SetProcess("LOSS",2);
166 gMC->SetProcess("MULS",1);
167 gMC->SetProcess("RAYL",1);
168
169 Float_t cut = 1.e-3; // 1MeV cut by default
170 Float_t tofmax = 1.e10;
171
172 gMC->SetCut("CUTGAM", cut);
173 gMC->SetCut("CUTELE", cut);
174 gMC->SetCut("CUTNEU", cut);
175 gMC->SetCut("CUTHAD", cut);
176 gMC->SetCut("CUTMUO", cut);
177 gMC->SetCut("BCUTE", cut);
178 gMC->SetCut("BCUTM", cut);
179 gMC->SetCut("DCUTE", cut);
180 gMC->SetCut("DCUTM", cut);
181 gMC->SetCut("PPCUTM", cut);
182 gMC->SetCut("TOFMAX", tofmax);
183
184
185
186
187 //======================//
188 // Set External decayer //
189 //======================//
190 TVirtualMCDecayer* decayer = new AliDecayerPythia();
191 decayer->SetForceDecay(kAll);
192 decayer->Init();
193 gMC->SetExternalDecayer(decayer);
194
195 //=========================//
196 // Generator Configuration //
197 //=========================//
198 AliGenerator* gener = 0x0;
199
200 if (proc == kPythia6) {
201 gener = MbPythia();
202 } else if (proc == kPythia6D6T) {
203 gener = MbPythiaTuneD6T();
204 } else if (proc == kPhojet) {
205 gener = MbPhojet();
206 }
207
208
209
210 // PRIMARY VERTEX
211 //
212 gener->SetOrigin(0., 0., 0.); // vertex position
213 //
214 //
215 // Size of the interaction diamond
216 // Longitudinal
217 Float_t sigmaz = 5.4 / TMath::Sqrt(2.); // [cm]
218 if (energy == 900)
219 sigmaz = 10.5 / TMath::Sqrt(2.); // [cm]
220 //
221 // Transverse
222 Float_t betast = 10; // beta* [m]
223 Float_t eps = 3.75e-6; // emittance [m]
224 Float_t gamma = energy / 2.0 / 0.938272; // relativistic gamma [1]
225 Float_t sigmaxy = TMath::Sqrt(eps * betast / gamma) / TMath::Sqrt(2.) * 100.; // [cm]
226 printf("\n \n Diamond size x-y: %10.3e z: %10.3e\n \n", sigmaxy, sigmaz);
227
228 gener->SetSigma(sigmaxy, sigmaxy, sigmaz); // Sigma in (X,Y,Z) (cm) on IP position
229 gener->SetCutVertexZ(3.); // Truncate at 3 sigma
230 gener->SetVertexSmear(kPerEvent);
231
232 gener->Init();
233
234 // FIELD
235 //
236 AliMagF* field = 0x0;
237 if (mag == kNoField) {
238 comment = comment.Append(" | L3 field 0.0 T");
239 field = new AliMagF("Maps","Maps", 0., 0., AliMagF::k5kGUniform,AliMagF::kBeamTypepp, energy/2.0);
240 } else if (mag == k5kG) {
241 comment = comment.Append(" | L3 field 0.5 T");
242 field = new AliMagF("Maps","Maps", -1., -1., AliMagF::k5kG, AliMagF::kBeamTypepp, energy/2.0);
243 }
244 printf("\n \n Comment: %s \n \n", comment.Data());
245
246 TGeoGlobalMagField::Instance()->SetField(field);
247
248 rl->CdGAFile();
249
250 Int_t iABSO = 1;
2990093d 251 Int_t iACORDE= 1;
d080385b 252 Int_t iAD = 0;
253 Int_t iDIPO = 1;
254 Int_t iEMCAL = 1;
255 Int_t iFMD = 1;
256 Int_t iFRAME = 1;
257 Int_t iHALL = 1;
258 Int_t iITS = 1;
259 Int_t iMAG = 1;
260 Int_t iMUON = 1;
261 Int_t iPHOS = 1;
262 Int_t iPIPE = 1;
263 Int_t iPMD = 1;
264 Int_t iHMPID = 1;
265 Int_t iSHIL = 1;
266 Int_t iT0 = 1;
267 Int_t iTOF = 1;
268 Int_t iTPC = 1;
269 Int_t iTRD = 1;
270 Int_t iVZERO = 1;
271 Int_t iZDC = 1;
272
273
274 //=================== Alice BODY parameters =============================
275 AliBODY *BODY = new AliBODY("BODY", "Alice envelop");
276
277
278 if (iMAG)
279 {
280 //=================== MAG parameters ============================
281 // --- Start with Magnet since detector layouts may be depending ---
282 // --- on the selected Magnet dimensions ---
283 AliMAG *MAG = new AliMAG("MAG", "Magnet");
284 }
285
286
287 if (iABSO)
288 {
289 //=================== ABSO parameters ============================
290 AliABSO *ABSO = new AliABSOv3("ABSO", "Muon Absorber");
291 }
292
293 if (iDIPO)
294 {
295 //=================== DIPO parameters ============================
296
297 AliDIPO *DIPO = new AliDIPOv3("DIPO", "Dipole version 3");
298 }
299
300 if (iHALL)
301 {
302 //=================== HALL parameters ============================
303
304 AliHALL *HALL = new AliHALLv3("HALL", "Alice Hall");
305 }
306
307
308 if (iFRAME)
309 {
310 //=================== FRAME parameters ============================
311
312 AliFRAMEv2 *FRAME = new AliFRAMEv2("FRAME", "Space Frame");
313 FRAME->SetHoles(1);
314 }
315
316 if (iSHIL)
317 {
318 //=================== SHIL parameters ============================
319
320 AliSHIL *SHIL = new AliSHILv3("SHIL", "Shielding Version 3");
321 }
322
323
324 if (iPIPE)
325 {
326 //=================== PIPE parameters ============================
327
328 AliPIPE *PIPE = new AliPIPEv3("PIPE", "Beam Pipe");
329 }
330
331 if (iITS)
332 {
333 //=================== ITS parameters ============================
334
335 AliITS *ITS = new AliITSv11("ITS","ITS v11");
336 }
337
338 if (iTPC)
339 {
340 //============================ TPC parameters =====================
341
342 AliTPC *TPC = new AliTPCv2("TPC", "Ne-CO2");
343 }
344
345
346 if (iTOF) {
347 //=================== TOF parameters ============================
348
349 AliTOF *TOF = new AliTOFv6T0("TOF", "normal TOF");
350 }
351
352
353 if (iHMPID)
354 {
355 //=================== HMPID parameters ===========================
356
357 AliHMPID *HMPID = new AliHMPIDv3("HMPID", "normal HMPID");
358
359 }
360
361
362 if (iZDC)
363 {
364 //=================== ZDC parameters ============================
365
366 AliZDC *ZDC = new AliZDCv4("ZDC", "normal ZDC");
367 }
368
369 if (iTRD)
370 {
371 //=================== TRD parameters ============================
372
373 AliTRD *TRD = new AliTRDv1("TRD", "TRD slow simulator");
374 // SM indexing starts at 3h in positive direction. We assume all SM
375 // are installed (i.e. all SM missing last year: 4, 5, 12, 13, 14)
376 // AliTRDgeometry *geoTRD = TRD->GetGeometry();
377 // geoTRD->SetSMstatus(4,0);
378 }
379
380 if (iFMD)
381 {
382 //=================== FMD parameters ============================
383
384 AliFMD *FMD = new AliFMDv1("FMD", "normal FMD");
385 }
386
387 if (iMUON)
388 {
389 //=================== MUON parameters ===========================
390 // New MUONv1 version (geometry defined via builders)
391
392 AliMUON *MUON = new AliMUONv1("MUON", "default");
393 }
394
395 if (iPHOS)
396 {
397 //=================== PHOS parameters ===========================
398
399 AliPHOS *PHOS = new AliPHOSv1("PHOS", "Modules1234");
400 }
401
402
403 if (iPMD)
404 {
405 //=================== PMD parameters ============================
406
407 AliPMD *PMD = new AliPMDv1("PMD", "normal PMD");
408 }
409
410 if (iT0)
411 {
412 //=================== T0 parameters ============================
413 AliT0 *T0 = new AliT0v1("T0", "T0 Detector");
414 }
415
416 if (iEMCAL)
417 {
418 //=================== EMCAL parameters ============================
419
420 AliEMCAL* EMCAL = new AliEMCALv2("EMCAL", "EMCAL_COMPLETE12SMV1_DCAL_8SM", kFALSE);
421 }
422
423 if (iACORDE)
424 {
425 //=================== ACORDE parameters ============================
426
427 AliACORDE *ACORDE = new AliACORDEv1("ACORDE", "normal ACORDE");
428 }
429
430 if (iVZERO)
431 {
432 //=================== VZERO parameters ============================
433
434 AliVZERO *VZERO = new AliVZEROv7("VZERO", "normal VZERO");
435 }
436
437 if (iAD)
438 {
439 //=================== AD parameters ============================
440
441 AliAD *AD = new AliADv1("AD", "The AD geometry");
442 }
443}
444//
445// PYTHIA
446//
447
448AliGenerator* MbPythia()
449{
450 comment = comment.Append(" pp at 14 TeV: Pythia low-pt");
451//
452// Pythia
453 AliGenPythia* pythia = new AliGenPythia(-1);
454 pythia->SetMomentumRange(0, 999999.);
455 pythia->SetThetaRange(0., 180.);
456 pythia->SetYRange(-12.,12.);
457 pythia->SetPtRange(0,1000.);
458 pythia->SetProcess(kPyMb);
459 pythia->SetEnergyCMS(energy);
460
461 return pythia;
462}
463
464AliGenerator* MbPythiaTuneD6T()
465{
466 comment = comment.Append(" pp at 14 TeV: Pythia low-pt");
467//
468// Pythia
469 AliGenPythia* pythia = new AliGenPythia(-1);
470 pythia->SetMomentumRange(0, 999999.);
471 pythia->SetThetaRange(0., 180.);
472 pythia->SetYRange(-12.,12.);
473 pythia->SetPtRange(0,1000.);
474 pythia->SetProcess(kPyMb);
475 pythia->SetEnergyCMS(energy);
476// Tune
477// 109 D6T : Rick Field's CDF Tune D6T (NB: needs CTEQ6L pdfs externally)
478 pythia->SetTune(109);
479 pythia->SetStrucFunc(kCTEQ6l);
480//
481 return pythia;
482}
483
484AliGenerator* MbPhojet()
485{
486 comment = comment.Append(" pp at 14 TeV: Pythia low-pt");
487//
488// DPMJET
489#if defined(__CINT__)
490 gSystem->Load("libdpmjet"); // Parton density functions
491 gSystem->Load("libTDPMjet"); // Parton density functions
492#endif
493 AliGenDPMjet* dpmjet = new AliGenDPMjet(-1);
494 dpmjet->SetMomentumRange(0, 999999.);
495 dpmjet->SetThetaRange(0., 180.);
496 dpmjet->SetYRange(-12.,12.);
497 dpmjet->SetPtRange(0,1000.);
498 dpmjet->SetProcess(kDpmMb);
499 dpmjet->SetEnergyCMS(energy);
500
501 return dpmjet;
502}
503
504void ProcessEnvironmentVars()
505{
506 // Run type
507 if (gSystem->Getenv("CONFIG_RUN_TYPE")) {
508 for (Int_t iRun = 0; iRun < kRunMax; iRun++) {
509 if (strcmp(gSystem->Getenv("CONFIG_RUN_TYPE"), pprRunName[iRun])==0) {
510 proc = (PDC06Proc_t)iRun;
511 cout<<"Run type set to "<<pprRunName[iRun]<<endl;
512 }
513 }
514 }
515
516 // Field
517 if (gSystem->Getenv("CONFIG_FIELD")) {
518 for (Int_t iField = 0; iField < kFieldMax; iField++) {
519 if (strcmp(gSystem->Getenv("CONFIG_FIELD"), pprField[iField])==0) {
520 mag = (Mag_t)iField;
521 cout<<"Field set to "<<pprField[iField]<<endl;
522 }
523 }
524 }
525
526 // Energy
527 if (gSystem->Getenv("CONFIG_ENERGY")) {
528 energy = atoi(gSystem->Getenv("CONFIG_ENERGY"));
529 cout<<"Energy set to "<<energy<<" GeV"<<endl;
530 }
531
532 // Random Number seed
533 if (gSystem->Getenv("CONFIG_SEED")) {
534 seed = atoi(gSystem->Getenv("CONFIG_SEED"));
535 }
536}