]> git.uio.no Git - u/mrichter/AliRoot.git/blob - prod/LHC08d6/Config.C
Working on the electron cut
[u/mrichter/AliRoot.git] / prod / LHC08d6 / Config.C
1 // One can use the configuration macro in compiled mode by
2 // root [0] gSystem->Load("libgeant321");
3 // root [0] gSystem->SetIncludePath("-I$ROOTSYS/include -I$ALICE_ROOT/include\
4 //                   -I$ALICE_ROOT -I$ALICE/geant3/TGeant3");
5 // root [0] .x grun.C(1,"Config.C++")
6
7 #if !defined(__CINT__) || defined(__MAKECINT__)
8 #include <Riostream.h>
9 #include <TPDGCode.h>
10 #include <TRandom.h>
11 #include <TSystem.h>
12 #include <TVirtualMC.h>
13 #include <TGeant3TGeo.h>
14 #include "STEER/AliRunLoader.h"
15 #include "STEER/AliRun.h"
16 #include "STEER/AliConfig.h"
17 #include "PYTHIA6/AliDecayerPythia.h"
18 #include "EVGEN/AliGenCocktail.h"
19 #include "EVGEN/AliGenHIJINGpara.h"
20 #include "STEER/AliMagF.h"
21 #include "STRUCT/AliBODY.h"
22 #include "STRUCT/AliMAG.h"
23 #include "STRUCT/AliABSOv3.h"
24 #include "STRUCT/AliDIPOv3.h"
25 #include "STRUCT/AliHALLv3.h"
26 #include "STRUCT/AliFRAMEv2.h"
27 #include "STRUCT/AliSHILv3.h"
28 #include "STRUCT/AliPIPEv3.h"
29 #include "ITS/AliITSv11Hybrid.h"
30 #include "TPC/AliTPCv2.h"
31 #include "TOF/AliTOFv6T0.h"
32 #include "HMPID/AliHMPIDv3.h"
33 #include "ZDC/AliZDCv3.h"
34 #include "TRD/AliTRDv1.h"
35 #include "FMD/AliFMDv1.h"
36 #include "MUON/AliMUONv1.h"
37 #include "PHOS/AliPHOSv1.h"
38 #include "PMD/AliPMDv1.h"
39 #include "T0/AliT0v1.h"
40 #include "EMCAL/AliEMCALv2.h"
41 #include "ACORDE/AliACORDEv1.h"
42 #include "VZERO/AliVZEROv7.h"
43 #endif
44
45 //--- Functions ---
46 void ProcessEnvironmentVars();
47 Float_t EtaToTheta(Float_t arg);
48 void    LoadPythia();
49
50 //--- Trigger config ---
51 enum TrigConf_t
52 {
53     kDefaultPPTrig, kDefaultPbPbTrig
54 };
55 const char * TrigConfName[] = {
56     "p-p","Pb-Pb"
57 };
58
59 // This part for configuration
60 static AliMagF::BMap_t         mag         = AliMagF::k5kG;
61 static TrigConf_t    trig        = kDefaultPbPbTrig; // default pp trigger configuration
62 static Int_t         runNumber   = 0;
63 static Int_t         eventNumber = 0;
64
65 //========================//
66 // Set Random Number seed //
67 //========================//
68 TDatime dt;
69 static UInt_t seed    = dt.Get();
70 static Int_t   runNumber= 0;
71 static Float_t bMin = 0.;
72 static Float_t bMax = 5.;
73 static UInt_t  quench = 1;
74 static UInt_t  shad = 1;
75 static Float_t etaMin = -8.0;
76 static Float_t etaMax = 8.0;
77 static Float_t phiMin = 0.;
78 static Float_t phiMax = 360.;
79 static TString comment;
80
81 void Config()
82 {
83   // ThetaRange is (0., 180.). It was (0.28,179.72) 7/12/00 09:00
84   // Theta range given through pseudorapidity limits 22/6/2001
85
86   // Set Random Number seed
87   //gRandom->SetSeed(123456); // Set 0 to use the current time
88
89   AliLog::Message(AliLog::kInfo, Form("Seed for random number generation = %d",gRandom->GetSeed()), "Config.C", "Config.C", "Config()","Config.C", __LINE__);
90
91   // Get settings from environment variables
92   ProcessEnvironmentVars();
93
94   // Load Pythia libraries
95   LoadPythia();
96   // Libraries required by geant321
97 #if defined(__CINT__)
98   gSystem->Load("libgeant321");
99 #endif
100
101   new     TGeant3TGeo("C++ Interface to Geant3");
102
103   // Output every 100 tracks
104   ((TGeant3*)gMC)->SetSWIT(4,100);
105
106   AliRunLoader* rl=0x0;
107
108     AliLog::Message(AliLog::kInfo, "Creating Run Loader", "Config.C", "Config.C", "Config()"," Config.C", __LINE__);
109
110     rl = AliRunLoader::Open("galice.root",
111                             AliConfig::GetDefaultEventFolderName(),
112                             "recreate");
113     if (rl == 0x0)
114       {
115         gAlice->Fatal("Config.C","Can not instatiate the Run Loader");
116         return;
117       }
118     rl->SetCompressionLevel(2);
119     rl->SetNumberOfEventsPerFile(2);
120     gAlice->SetRunLoader(rl);
121
122     // gAlice->SetGeometryFromFile("geometry.root");
123
124     // Uncomment if you want to load geometry from OCDB!   >>>>
125 /*
126     if(!AliCDBManager::Instance()->IsDefaultStorageSet()){
127          cout << "#####################################################" << endl;
128          cout << "#                                                   #" << endl;
129          cout << "#     WARNING: CDB DEFAULT STORAGE NOT SET !!!      #" << endl;
130          cout << "#     SETTING IT TO local://$ALICE_ROOT/OCDB !!!         #" << endl;
131          cout << "#                                                   #" << endl;
132          cout << "#####################################################" << endl;
133
134          AliCDBManager::Instance()->SetDefaultStorage("local://$ALICE_ROOT/OCDB");
135     }
136
137     if(AliCDBManager::Instance()->GetRun() < 0){
138          cout << "#####################################################" << endl;
139          cout << "#                                                   #" << endl;
140          cout << "#     WARNING: RUN NUMBER NOT SET !!!               #" << endl;
141          cout << "#     SETTING IT TO 0 !!!                           #" << endl;
142          cout << "#                                                   #" << endl;
143          cout << "#####################################################" << endl;
144
145          AliCDBManager::Instance()->SetRun(0);
146     }
147     gAlice->SetGeometryFromCDB();
148 */
149     // Uncomment if you want to load geometry from OCDB!   <<<<
150
151     // Set the trigger configuration
152     AliSimulation::Instance()->SetTriggerConfig(TrigConfName[trig]);
153     cout<<"Trigger configuration is set to  "<<TrigConfName[trig]<<endl;
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     AliGenHijing *gener = new AliGenHijing(-1);
201     // centre of mass energy
202     gener->SetEnergyCMS(5500);
203     // reference frame
204     gener->SetReferenceFrame("CMS     ");
205     // projectile
206     gener->SetProjectile("A       ", 208, 82);
207     gener->SetTarget    ("A       ", 208, 82);
208     // impact parameter range
209     gener->SetImpactParameterRange(bMin, bMax); // bMin = 0 - bMax = 3
210     // evaluate cross section before run
211     gener->SetEvaluate(0);
212     // tell hijing to keep the full parent child chain
213     gener->KeepFullEvent();
214     // enable jet quenching
215     gener->SetJetQuenching(quench); // 1
216     // enable shadowing
217     gener->SetShadowing(shad); // 1
218     // neutral pion and heavy particle decays switched off
219     gener->SetDecaysOff(1);
220     // Don't track spectators
221     gener->SetSpectators(0);
222     // trigger
223     //  gener->SetTrigger(0);
224     // kinematic selection
225     gener->SetSelectAll(0);
226     // momentum range
227     gener->SetMomentumRange(0,999);
228     // No resytriction on phi, theta
229     Float_t thmin = EtaToTheta(etaMax); // Theta min <---> eta max 2.
230     Float_t thmax = EtaToTheta(etaMin); // Theta max <---> eta min -2.
231     gener->SetPhiRange(phiMin, phiMax); // 0 - 360
232     gener->SetThetaRange(thmin,thmax);
233     // PRIMARY VERTEX
234     gener->SetOrigin(0, 0, 0);  //vertex position
235 //    gener->SetSigma(0, 0, 5.3);   //Sigma in (X,Y,Z) (cm) on IP position
236 //    gener->SetCutVertexZ(3.);        // Truncate at 3 sigma
237 //    gener->SetVertexSmear(kPerEvent);
238
239     // Size of the interaction diamond
240     // Longitudinal
241     Float_t sigmaz  = 7.55 / TMath::Sqrt(2.); // [cm]
242     // Transverse
243     Float_t betast  = 10;                 // beta* [m]
244     Float_t eps     = 3.75e-6;            // emittance [m]
245 //    Float_t gamma   = 7000. / 0.938272;   // relativistic gamma [1]
246     Float_t gamma   = 2750. / 0.938272;   // relativistic gamma [1]
247     Float_t sigmaxy = TMath::Sqrt(eps * betast / gamma) / TMath::Sqrt(2.) * 100.;  // [cm]
248     printf("\n \n Diamond size x-y: %10.3e z: %10.3e\n \n", sigmaxy, sigmaz);
249
250     gener->SetSigma(sigmaxy, sigmaxy, sigmaz);      // Sigma in (X,Y,Z) (cm) on IP position
251     gener->SetCutVertexZ(3.);        // Truncate at 3 sigma
252     gener->SetVertexSmear(kPerEvent);
253
254     //
255     // Activate this line if you want the vertex smearing to happen
256     // track by track
257     //
258     //gener->SetVertexSmear(perTrack);
259
260     gener->Init();
261
262     // Field (L3 0.5 T)
263     TGeoGlobalMagField::Instance()->SetField(new AliMagF("Maps","Maps", -1., -1., mag));
264
265     Int_t   iABSO   = 1;
266     Int_t   iDIPO   = 1;
267     Int_t   iFMD    = 1;
268     Int_t   iFRAME  = 1;
269     Int_t   iHALL   = 1;
270     Int_t   iITS    = 1;
271     Int_t   iMAG    = 1;
272     Int_t   iMUON   = 1;
273     Int_t   iPHOS   = 1;
274     Int_t   iPIPE   = 1;
275     Int_t   iPMD    = 1;
276     Int_t   iHMPID  = 1;
277     Int_t   iSHIL   = 1;
278     Int_t   iT0     = 1;
279     Int_t   iTOF    = 1;
280     Int_t   iTPC    = 1;
281     Int_t   iTRD    = 1;
282     Int_t   iZDC    = 1;
283     Int_t   iEMCAL  = 1;
284     Int_t   iACORDE = 1;
285     Int_t   iVZERO  = 1;
286     rl->CdGAFile();
287     //=================== Alice BODY parameters =============================
288     AliBODY *BODY = new AliBODY("BODY", "Alice envelop");
289
290     if (iMAG)
291     {
292         //=================== MAG parameters ============================
293         // --- Start with Magnet since detector layouts may be depending ---
294         // --- on the selected Magnet dimensions ---
295         AliMAG *MAG = new AliMAG("MAG", "Magnet");
296     }
297
298
299     if (iABSO)
300     {
301         //=================== ABSO parameters ============================
302         AliABSO *ABSO = new AliABSOv3("ABSO", "Muon Absorber");
303     }
304
305     if (iDIPO)
306     {
307         //=================== DIPO parameters ============================
308
309         AliDIPO *DIPO = new AliDIPOv3("DIPO", "Dipole version 3");
310     }
311
312     if (iHALL)
313     {
314         //=================== HALL parameters ============================
315
316         AliHALL *HALL = new AliHALLv3("HALL", "Alice Hall");
317     }
318
319
320     if (iFRAME)
321     {
322         //=================== FRAME parameters ============================
323
324         AliFRAMEv2 *FRAME = new AliFRAMEv2("FRAME", "Space Frame");
325         FRAME->SetHoles(1);
326     }
327
328     if (iSHIL)
329     {
330         //=================== SHIL parameters ============================
331
332         AliSHIL *SHIL = new AliSHILv3("SHIL", "Shielding Version 3");
333     }
334
335
336     if (iPIPE)
337     {
338         //=================== PIPE parameters ============================
339
340         AliPIPE *PIPE = new AliPIPEv3("PIPE", "Beam Pipe");
341     }
342
343     if (iITS)
344     {
345         //=================== ITS parameters ============================
346
347         AliITS *ITS  = new AliITSv11Hybrid("ITS","ITS v11Hybrid");
348     }
349
350     if (iTPC)
351     {
352         //============================ TPC parameters ===================
353         AliTPC *TPC = new AliTPCv2("TPC", "Default");
354     }
355
356
357     if (iTOF) {
358         //=================== TOF parameters ============================
359         AliTOF *TOF = new AliTOFv6T0("TOF", "normal TOF");
360     }
361
362
363     if (iHMPID)
364     {
365         //=================== HMPID parameters ===========================
366         AliHMPID *HMPID = new AliHMPIDv3("HMPID", "normal HMPID");
367
368     }
369
370
371     if (iZDC)
372     {
373         //=================== ZDC parameters ============================
374
375         AliZDC *ZDC = new AliZDCv3("ZDC", "normal ZDC");
376     }
377
378     if (iTRD)
379     {
380         //=================== TRD parameters ============================
381
382         AliTRD *TRD = new AliTRDv1("TRD", "TRD slow simulator");
383     }
384
385     if (iFMD)
386     {
387         //=================== FMD parameters ============================
388         AliFMD *FMD = new AliFMDv1("FMD", "normal FMD");
389    }
390
391     if (iMUON)
392     {
393         //=================== MUON parameters ===========================
394         // New MUONv1 version (geometry defined via builders)
395         AliMUON *MUON = new AliMUONv1("MUON", "default");
396     }
397     //=================== PHOS parameters ===========================
398
399     if (iPHOS)
400     {
401         AliPHOS *PHOS = new AliPHOSv1("PHOS", "IHEP");
402     }
403
404
405     if (iPMD)
406     {
407         //=================== PMD parameters ============================
408         AliPMD *PMD = new AliPMDv1("PMD", "normal PMD");
409     }
410
411     if (iT0)
412     {
413         //=================== T0 parameters ============================
414         AliT0 *T0 = new AliT0v1("T0", "T0 Detector");
415     }
416
417     if (iEMCAL)
418     {
419         //=================== EMCAL parameters ============================
420         AliEMCAL *EMCAL = new AliEMCALv2("EMCAL", "EMCAL_COMPLETE");
421     }
422
423      if (iACORDE)
424     {
425         //=================== ACORDE parameters ============================
426         AliACORDE *ACORDE = new AliACORDEv1("ACORDE", "normal ACORDE");
427     }
428
429      if (iVZERO)
430     {
431         //=================== VZERO parameters ============================
432         AliVZERO *VZERO = new AliVZEROv7("VZERO", "normal VZERO");
433     }
434
435      AliLog::Message(AliLog::kInfo, "End of Config", "Config.C", "Config.C", "Config()"," Config.C", __LINE__);
436
437 }
438
439 Float_t EtaToTheta(Float_t arg){
440   return (180./TMath::Pi())*2.*atan(exp(-arg));
441 }
442
443 void ProcessEnvironmentVars()
444 {
445     cout << "Processing environment variables" << endl;
446     // Random Number seed
447     if (gSystem->Getenv("CONFIG_SEED")) {
448       seed = atoi(gSystem->Getenv("CONFIG_SEED"));
449     }
450
451     gRandom->SetSeed(seed);
452     cout<<"Seed for random number generation= "<<seed<<endl;
453
454     // Run Number
455     if (gSystem->Getenv("DC_RUN")) {
456       runNumber = atoi(gSystem->Getenv("DC_RUN"));
457     }
458     cout<<"Run number "<<runNumber<<endl;
459
460     // Impact param
461     if (gSystem->Getenv("CONFIG_BMIN")) {
462       bMin = atof(gSystem->Getenv("CONFIG_BMIN"));
463     }
464
465     if (gSystem->Getenv("CONFIG_BMAX")) {
466       bMax = atof(gSystem->Getenv("CONFIG_BMAX"));
467     }
468     cout<<"Impact parameter in ["<<bMin<<","<<bMax<<"]"<<endl;
469
470     // Quenching scenario
471     if (gSystem->Getenv("QUENCH")) {
472       quench = atoi(gSystem->Getenv("QUENCH"));
473     }
474     if(quench==1) cout<<"With quenching "<<endl;
475     else cout<<"Without quenching "<<endl;
476
477     // Quenching scenario
478     if (gSystem->Getenv("QHAT")) {
479       shad = atoi(gSystem->Getenv("QHAT"));
480     }
481     if(shad==1) cout<<"With shadowing "<<endl;
482     else cout<<"Without shadowing "<<endl;
483
484     // Acceptance param
485     if (gSystem->Getenv("CONFIG_ETAMIN")) {
486       etaMin = atof(gSystem->Getenv("CONFIG_ETAMIN"));
487     }
488
489     if (gSystem->Getenv("CONFIG_ETAMAX")) {
490       etaMax = atof(gSystem->Getenv("CONFIG_ETAMAX"));
491     }
492     cout<<"Eta acceptance ["<<etaMin<<","<<etaMax<<"]"<<endl;
493
494     if (gSystem->Getenv("CONFIG_PHIMIN")) {
495       phiMin = atof(gSystem->Getenv("CONFIG_PHIMIN"));
496     }
497
498     if (gSystem->Getenv("CONFIG_PHIMAX")) {
499       phiMax = atof(gSystem->Getenv("CONFIG_PHIMAX"));
500     }
501     cout<<"Phi acceptance ["<<phiMin<<","<<phiMax<<"]"<<endl;
502
503
504 }
505
506 void LoadPythia()
507 {
508     // Load Pythia related libraries
509     gSystem->Load("liblhapdf.so");      // Parton density functions
510     gSystem->Load("libEGPythia6.so");   // TGenerator interface
511     gSystem->Load("libpythia6.so");     // Pythia
512     gSystem->Load("libAliPythia6.so");  // ALICE specific implementations
513 }