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