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