]> git.uio.no Git - u/mrichter/AliRoot.git/blob - test/fpprod/Config.C
The present commit corresponds to an important change in the way the
[u/mrichter/AliRoot.git] / test / fpprod / 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/AliMagFCheb.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   gAlice->SetTriggerDescriptor("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   AliGenerator* gener = 0x0;
186   
187   if (proc == kPythia6) {
188       gener = MbPythia();
189   } else if (proc == kPhojet) {
190       gener = MbPhojet();
191   }
192   
193   
194
195   // PRIMARY VERTEX
196   //
197   gener->SetOrigin(0., 0., 0.);    // vertex position
198   //
199   //
200   // Size of the interaction diamond
201   // Longitudinal
202   Float_t sigmaz  = 5.4 / TMath::Sqrt(2.); // [cm]
203   if (energy == 900)
204     sigmaz  = 10.5 / TMath::Sqrt(2.); // [cm]
205   //
206   // Transverse
207   Float_t betast  = 10;                 // beta* [m]
208   Float_t eps     = 3.75e-6;            // emittance [m]
209   Float_t gamma   = energy / 2.0 / 0.938272;  // relativistic gamma [1]
210   Float_t sigmaxy = TMath::Sqrt(eps * betast / gamma) / TMath::Sqrt(2.) * 100.;  // [cm]
211   printf("\n \n Diamond size x-y: %10.3e z: %10.3e\n \n", sigmaxy, sigmaz);
212     
213   gener->SetSigma(sigmaxy, sigmaxy, sigmaz);      // Sigma in (X,Y,Z) (cm) on IP position
214   gener->SetCutVertexZ(3.);        // Truncate at 3 sigma
215   gener->SetVertexSmear(kPerEvent);
216
217   gener->Init();
218
219   // FIELD
220   //
221   AliMagF* field = 0x0;
222   if (mag == kNoField) {
223     comment = comment.Append(" | L3 field 0.0 T");
224     field = new AliMagF("Maps","Maps", 2, 0., 0., 10., AliMagF::k2kG);
225   } else if (mag == k5kG) {
226     comment = comment.Append(" | L3 field 0.5 T");
227     field = new AliMagF("Maps","Maps", 2, 1., 1., 10., AliMagF::k5kG);
228   }
229   printf("\n \n Comment: %s \n \n", comment.Data());
230
231   TGeoGlobalMagField::Instance()->SetField(field);
232     
233   rl->CdGAFile();
234   
235   Int_t iABSO  = 1;
236   Int_t iACORDE= 0;
237   Int_t iDIPO  = 1;
238   Int_t iEMCAL = 0;
239   Int_t iFMD   = 1;
240   Int_t iFRAME = 1;
241   Int_t iHALL  = 1;
242   Int_t iITS   = 1;
243   Int_t iMAG   = 1;
244   Int_t iMUON  = 1;
245   Int_t iPHOS  = 1;
246   Int_t iPIPE  = 1;
247   Int_t iPMD   = 0;
248   Int_t iHMPID = 1;
249   Int_t iSHIL  = 1;
250   Int_t iT0    = 1;
251   Int_t iTOF   = 1;
252   Int_t iTPC   = 1;
253   Int_t iTRD   = 1;
254   Int_t iVZERO = 1;
255   Int_t iZDC   = 1;
256   
257
258     //=================== Alice BODY parameters =============================
259     AliBODY *BODY = new AliBODY("BODY", "Alice envelop");
260
261
262     if (iMAG)
263     {
264         //=================== MAG parameters ============================
265         // --- Start with Magnet since detector layouts may be depending ---
266         // --- on the selected Magnet dimensions ---
267         AliMAG *MAG = new AliMAG("MAG", "Magnet");
268     }
269
270
271     if (iABSO)
272     {
273         //=================== ABSO parameters ============================
274         AliABSO *ABSO = new AliABSOv3("ABSO", "Muon Absorber");
275     }
276
277     if (iDIPO)
278     {
279         //=================== DIPO parameters ============================
280
281         AliDIPO *DIPO = new AliDIPOv3("DIPO", "Dipole version 3");
282     }
283
284     if (iHALL)
285     {
286         //=================== HALL parameters ============================
287
288         AliHALL *HALL = new AliHALLv3("HALL", "Alice Hall");
289     }
290
291
292     if (iFRAME)
293     {
294         //=================== FRAME parameters ============================
295
296         AliFRAMEv2 *FRAME = new AliFRAMEv2("FRAME", "Space Frame");
297         FRAME->SetHoles(1);
298     }
299
300     if (iSHIL)
301     {
302         //=================== SHIL parameters ============================
303
304         AliSHIL *SHIL = new AliSHILv3("SHIL", "Shielding Version 3");
305     }
306
307
308     if (iPIPE)
309     {
310         //=================== PIPE parameters ============================
311
312         AliPIPE *PIPE = new AliPIPEv3("PIPE", "Beam Pipe");
313     }
314  
315     if (iITS)
316     {
317         //=================== ITS parameters ============================
318
319         AliITS *ITS  = new AliITSv11Hybrid("ITS","ITS v11Hybrid");
320     }
321
322     if (iTPC)
323     {
324       //============================ TPC parameters =====================
325
326         AliTPC *TPC = new AliTPCv2("TPC", "Default");
327     }
328
329
330     if (iTOF) {
331         //=================== TOF parameters ============================
332
333         AliTOF *TOF = new AliTOFv6T0("TOF", "normal TOF");
334     }
335
336
337     if (iHMPID)
338     {
339         //=================== HMPID parameters ===========================
340
341         AliHMPID *HMPID = new AliHMPIDv3("HMPID", "normal HMPID");
342
343     }
344
345
346     if (iZDC)
347     {
348         //=================== ZDC parameters ============================
349
350         AliZDC *ZDC = new AliZDCv3("ZDC", "normal ZDC");
351     }
352
353     if (iTRD)
354     {
355         //=================== TRD parameters ============================
356
357         AliTRD *TRD = new AliTRDv1("TRD", "TRD slow simulator");
358         AliTRDgeometry *geoTRD = TRD->GetGeometry();
359         // Partial geometry: modules at 0,8,9,17
360         // starting at 3h in positive direction
361         geoTRD->SetSMstatus(1,0);
362         geoTRD->SetSMstatus(2,0);
363         geoTRD->SetSMstatus(3,0);
364         geoTRD->SetSMstatus(4,0);
365         geoTRD->SetSMstatus(5,0);
366         geoTRD->SetSMstatus(6,0);
367         geoTRD->SetSMstatus(7,0);
368         geoTRD->SetSMstatus(10,0);
369         geoTRD->SetSMstatus(11,0);
370         geoTRD->SetSMstatus(12,0);
371         geoTRD->SetSMstatus(13,0);
372         geoTRD->SetSMstatus(14,0);
373         geoTRD->SetSMstatus(15,0);
374         geoTRD->SetSMstatus(16,0);
375     }
376
377     if (iFMD)
378     {
379         //=================== FMD parameters ============================
380
381         AliFMD *FMD = new AliFMDv1("FMD", "normal FMD");
382    }
383
384     if (iMUON)
385     {
386         //=================== MUON parameters ===========================
387         // New MUONv1 version (geometry defined via builders)
388
389         AliMUON *MUON = new AliMUONv1("MUON", "default");
390     }
391
392     if (iPHOS)
393     {
394         //=================== PHOS parameters ===========================
395
396         AliPHOS *PHOS = new AliPHOSv1("PHOS", "IHEP");
397         //Set simulation parameters different from the default ones.
398         AliPHOSSimParam* simEmc = AliPHOSSimParam::GetInstance() ;
399   
400         // APD noise of warm (+20C) PHOS:
401         // a2 = a1*(Y1/Y2)*(M1/M2), where a1 = 0.012 is APD noise at -25C,
402         // Y1 = 4.3 photo-electrons/MeV, Y2 = 1.7 p.e/MeV - light yields at -25C and +20C,
403         // M1 = 50, M2 = 50 - APD gain factors chosen for t1 = -25C and t2 = +20C,
404         // Y = MeanLightYield*APDEfficiency.
405
406         Float_t apdNoise = 0.012*2.5; 
407         simEmc->SetAPDNoise(apdNoise);
408
409         //Raw Light Yield at +20C
410         simEmc->SetMeanLightYield(18800);
411
412         //ADC channel width at +18C.
413         simEmc->SetADCchannelW(0.0125);
414     }
415
416
417     if (iPMD)
418     {
419         //=================== PMD parameters ============================
420
421         AliPMD *PMD = new AliPMDv1("PMD", "normal PMD");
422     }
423
424     if (iT0)
425     {
426         //=================== T0 parameters ============================
427         AliT0 *T0 = new AliT0v1("T0", "T0 Detector");
428     }
429
430     if (iEMCAL)
431     {
432         //=================== EMCAL parameters ============================
433
434         AliEMCAL *EMCAL = new AliEMCALv2("EMCAL", "EMCAL_COMPLETE");
435     }
436
437      if (iACORDE)
438     {
439         //=================== ACORDE parameters ============================
440
441         AliACORDE *ACORDE = new AliACORDEv1("ACORDE", "normal ACORDE");
442     }
443
444      if (iVZERO)
445     {
446         //=================== ACORDE parameters ============================
447
448         AliVZERO *VZERO = new AliVZEROv7("VZERO", "normal VZERO");
449     }
450 }
451 //
452 //           PYTHIA
453 //
454
455 AliGenerator* MbPythia()
456 {
457       comment = comment.Append(" pp at 14 TeV: Pythia low-pt");
458 //
459 //    Pythia
460       AliGenPythia* pythia = new AliGenPythia(-1); 
461       pythia->SetMomentumRange(0, 999999.);
462       pythia->SetThetaRange(0., 180.);
463       pythia->SetYRange(-12.,12.);
464       pythia->SetPtRange(0,1000.);
465       pythia->SetProcess(kPyMb);
466       pythia->SetEnergyCMS(energy);
467       
468       return pythia;
469 }
470
471 AliGenerator* MbPhojet()
472 {
473       comment = comment.Append(" pp at 14 TeV: Pythia low-pt");
474 //
475 //    DPMJET
476 #if defined(__CINT__)
477   gSystem->Load("libdpmjet");      // Parton density functions
478   gSystem->Load("libTDPMjet");      // Parton density functions
479 #endif
480       AliGenDPMjet* dpmjet = new AliGenDPMjet(-1); 
481       dpmjet->SetMomentumRange(0, 999999.);
482       dpmjet->SetThetaRange(0., 180.);
483       dpmjet->SetYRange(-12.,12.);
484       dpmjet->SetPtRange(0,1000.);
485       dpmjet->SetProcess(kDpmMb);
486       dpmjet->SetEnergyCMS(energy);
487
488       return dpmjet;
489 }
490
491 void ProcessEnvironmentVars()
492 {
493     // Run type
494     if (gSystem->Getenv("CONFIG_RUN_TYPE")) {
495       for (Int_t iRun = 0; iRun < kRunMax; iRun++) {
496         if (strcmp(gSystem->Getenv("CONFIG_RUN_TYPE"), pprRunName[iRun])==0) {
497           proc = (PDC06Proc_t)iRun;
498           cout<<"Run type set to "<<pprRunName[iRun]<<endl;
499         }
500       }
501     }
502
503     // Field
504     if (gSystem->Getenv("CONFIG_FIELD")) {
505       for (Int_t iField = 0; iField < kFieldMax; iField++) {
506         if (strcmp(gSystem->Getenv("CONFIG_FIELD"), pprField[iField])==0) {
507           mag = (Mag_t)iField;
508           cout<<"Field set to "<<pprField[iField]<<endl;
509         }
510       }
511     }
512
513     // Energy
514     if (gSystem->Getenv("CONFIG_ENERGY")) {
515       energy = atoi(gSystem->Getenv("CONFIG_ENERGY"));
516       cout<<"Energy set to "<<energy<<" GeV"<<endl;
517     }
518
519     // Random Number seed
520     if (gSystem->Getenv("CONFIG_SEED")) {
521       seed = atoi(gSystem->Getenv("CONFIG_SEED"));
522     }
523 }