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