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