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