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