]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STARLIGHT/test/Config.C
updated STARLIGHT to r191 of http://starlight.hepforge.org/svn/trunk
[u/mrichter/AliRoot.git] / STARLIGHT / test / Config.C
1 //
2 // Configuration for p-p Diffraction studies
3 // STARLIGHT pp 7000
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   kStarlight, kRunMax,
55 };
56
57 const char * pprRunName[] = {
58   "kStarlight"
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 void ProcessEnvironmentVars();
72
73 // Generator, field, beam energy
74 static PDC06Proc_t   proc         = kStarlight;
75 static Mag_t         mag          = k5kG;
76 static Float_t       energy       = 7000; // energy in CMS
77 static Float_t       pBeamEnergy  = 3500.; // energy p-Beam
78 static Int_t         runNumber    = 0; 
79 static Int_t         eventOffset  = 0;
80
81 //========================//
82 // Set Random Number seed //
83 //========================//
84 TDatime dt;
85 static UInt_t seed    = dt.Get();
86
87 // Comment line
88 static TString comment;
89
90 void Config()
91 {
92   // Get settings from environment variables
93   ProcessEnvironmentVars();
94
95   gRandom->SetSeed(seed);
96   cerr<<"Seed for random number generation= "<<seed<<endl;
97
98   // Libraries required by geant321
99 #if defined(__CINT__)
100   gSystem->Load("liblhapdf");      // Parton density functions
101   gSystem->Load("libEGPythia6");   // TGenerator interface
102 //   if (proc == kPythia6 || proc == kPhojet) {
103 //     gSystem->Load("libpythia6");        // Pythia 6.2
104 //   } else {
105      gSystem->Load("libpythia6.4.21");   // Pythia 6.4
106 //   }
107   gSystem->Load("libAliPythia6");  // ALICE specific implementations
108   gSystem->Load("libgeant321");
109 #endif
110
111   new TGeant3TGeo("C++ Interface to Geant3");
112
113   //=======================================================================
114   //  Create the output file
115
116
117   AliRunLoader* rl=0x0;
118
119   cout<<"Config.C: Creating Run Loader ..."<<endl;
120   rl = AliRunLoader::Open("galice.root",
121                           AliConfig::GetDefaultEventFolderName(),
122                           "recreate");
123   if (rl == 0x0)
124     {
125       gAlice->Fatal("Config.C","Can not instatiate the Run Loader");
126       return;
127     }
128   rl->SetCompressionLevel(2);
129   rl->SetNumberOfEventsPerFile(1000);
130 //   rl->MakeHeader();
131 //   rl->MakeStack();
132   gAlice->SetRunLoader(rl);
133   // gAlice->SetGeometryFromFile("geometry.root");
134   // gAlice->SetGeometryFromCDB();
135
136   // Set the trigger configuration: proton-proton
137   AliSimulation::Instance()->SetTriggerConfig("p-p"); 
138
139   //
140   //=======================================================================
141   // ************* STEERING parameters FOR ALICE SIMULATION **************
142   // --- Specify event type to be tracked through the ALICE setup
143   // --- All positions are in cm, angles in degrees, and P and E in GeV
144
145   gMC->SetProcess("DCAY",1);
146   gMC->SetProcess("PAIR",1);
147   gMC->SetProcess("COMP",1);
148   gMC->SetProcess("PHOT",1);
149   gMC->SetProcess("PFIS",0);
150   gMC->SetProcess("DRAY",0);
151   gMC->SetProcess("ANNI",1);
152   gMC->SetProcess("BREM",1);
153   gMC->SetProcess("MUNU",1);
154   gMC->SetProcess("CKOV",1);
155   gMC->SetProcess("HADR",1);
156   gMC->SetProcess("LOSS",2);
157   gMC->SetProcess("MULS",1);
158   gMC->SetProcess("RAYL",1);
159
160   Float_t cut = 1.e-3;        // 1MeV cut by default
161   Float_t tofmax = 1.e10;
162
163   gMC->SetCut("CUTGAM", cut);
164   gMC->SetCut("CUTELE", cut);
165   gMC->SetCut("CUTNEU", cut);
166   gMC->SetCut("CUTHAD", cut);
167   gMC->SetCut("CUTMUO", cut);
168   gMC->SetCut("BCUTE",  cut);
169   gMC->SetCut("BCUTM",  cut);
170   gMC->SetCut("DCUTE",  cut);
171   gMC->SetCut("DCUTM",  cut);
172   gMC->SetCut("PPCUTM", cut);
173   gMC->SetCut("TOFMAX", tofmax);
174
175   //======================//
176   // Set External decayer //
177   //======================//
178   TVirtualMCDecayer* decayer = new AliDecayerPythia();
179   decayer->SetForceDecay(kAll);
180
181   decayer->Init();
182   gMC->SetExternalDecayer(decayer);
183
184   //=========================//
185   // Generator Configuration //
186   //=========================//
187   AliGenerator* gener = 0x0;
188
189   if (proc == kStarlight) {
190     TString ap="../../build/lib/tgt_linux/";
191     gSystem->Load(ap+"libStarLight.so");
192     gSystem->Load(ap+"libAliStarLight.so");
193
194     AliGenStarLight* sl = new AliGenStarLight(1000*1000);
195
196     sl->SetParameter("BEAM_1_Z     =    1    #Z of projectile");
197     sl->SetParameter("BEAM_1_A     =    1    #A of projectile");
198     sl->SetParameter("BEAM_2_Z     =   82    #Z of target");
199     sl->SetParameter("BEAM_2_A     =  208    #A of target");
200     sl->SetParameter("BEAM_1_GAMMA = 4264.4  #Gamma of the colliding ions");
201     sl->SetParameter("BEAM_2_GAMMA = 1682.4  #Gamma of the colliding ions");
202     sl->SetParameter("W_MAX        =   15    #Max value of w");
203     sl->SetParameter("W_MIN        =    -1   #Min value of w");
204     sl->SetParameter("W_N_BINS     =   40    #Bins i w");
205     sl->SetParameter("RAP_MAX      =    8.   #max y");
206     sl->SetParameter("RAP_N_BINS   =   80    #Bins i y");
207     sl->SetParameter("CUT_PT       =    0    #Cut in pT? 0 = (no, 1 = yes)");
208     sl->SetParameter("PT_MIN       =    1.0  #Minimum pT in GeV");
209     sl->SetParameter("PT_MAX       =    3.0  #Maximum pT in GeV");
210     sl->SetParameter("CUT_ETA      =    0    #Cut in pseudorapidity? (0 = no, 1 = yes)");
211     sl->SetParameter("ETA_MIN      =  -10    #Minimum pseudorapidity");
212     sl->SetParameter("ETA_MAX      =   10    #Maximum pseudorapidity");
213     sl->SetParameter("PROD_MODE    =    4    #gg or gP switch (1 = 2-photon, 2 = coherent vector meson (narrow), 3 = coherent vector meson (wide), # 4 = incoherent vector meson, 5 = A+A DPMJet single, 6 = A+A DPMJet double, 7 = p+A DPMJet single, 8 = p+A Pythia single )");
214     sl->SetParameter("PROD_PID     =  913    #Channel of interest (not relevant for photonuclear processes)");
215     sl->SetParameter("RND_SEED     = 34533   #Random number seed");
216     sl->SetParameter("BREAKUP_MODE  =   5    #Controls the nuclear breakup");
217     sl->SetParameter("INTERFERENCE  =   0    #Interference (0 = off, 1 = on)");
218     sl->SetParameter("IF_STRENGTH   =   1.   #% of intefernce (0.0 - 0.1)");
219     sl->SetParameter("COHERENT      =   1    #Coherent=1,Incoherent=0");
220     sl->SetParameter("INCO_FACTOR   =   1.   #percentage of incoherence");
221     sl->SetParameter("INT_PT_MAX    =   0.24 #Maximum pt considered, when interference is turned on");
222     sl->SetParameter("INT_PT_N_BINS = 120    #Number of pt bins when interference is turned on");
223
224     sl->SetRapidityMotherRange(-1., 1.);
225
226     sl->Init();
227     sl->GetTStarLight()->PrintInputs(std::cout);
228     
229     AliGenerator* gener = sl;
230     gener->SetOrigin(0, 0, 0);    // vertex position
231     gener->SetSigma(0, 0, 5.3);   // Sigma in (X,Y,Z) (cm) on IP position
232     gener->SetCutVertexZ(1.);     // Truncate at 1 sigma
233     gener->SetVertexSmear(kPerEvent); 
234     gener->SetTrackingFlag(1);
235 //     gener->SetYRange(-0.9, 0.9);
236     gener->Init();
237   }
238
239   //
240   //
241   // Size of the interaction diamond
242   // Longitudinal
243   Float_t sigmaz  = 5.4 / TMath::Sqrt(2.); // [cm]
244   //   if (energy == 900)
245   //sigmaz  = 10.5 / TMath::Sqrt(2.); // [cm]
246   //sigmaz = 3.7;
247   if (energy == 7000)
248     sigmaz  = 6.3 / TMath::Sqrt(2.); // [cm]
249   
250   //
251   // Transverse
252   Float_t betast  = 10;                 // beta* [m]
253   if (runNumber >= 117048) betast = 2.;
254   if (runNumber >  122375) betast = 3.5; // starting with fill 1179
255
256   Float_t eps     = 5.e-6;            // emittance [m]
257   Float_t gamma   = energy / 2.0 / 0.938272;  // relativistic gamma [1]
258   Float_t sigmaxy = TMath::Sqrt(eps * betast / gamma) / TMath::Sqrt(2.) * 100.;  // [cm]
259   printf("\n \n Diamond size x-y: %10.3e z: %10.3e\n \n", sigmaxy, sigmaz);
260
261   gener->SetSigma(sigmaxy, sigmaxy, sigmaz);      // Sigma in (X,Y,Z) (cm) on IP position
262   gener->SetVertexSmear(kPerEvent);
263   gener->Init();
264
265   printf("\n \n Comment: %s \n \n", comment.Data());
266
267   //!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
268   TGeoGlobalMagField::Instance()->SetField(new AliMagF("Maps","Maps", -1., -1., AliMagF::k5kG,
269                                                        AliMagF::kBeamTypepp, pBeamEnergy));
270
271   rl->CdGAFile();
272
273   Int_t iABSO  = 1;
274   Int_t iACORDE= 0;
275   Int_t iDIPO  = 1;
276   Int_t iEMCAL = 1;
277   Int_t iFMD   = 1;
278   Int_t iFRAME = 1;
279   Int_t iHALL  = 1;
280   Int_t iITS   = 1;
281   Int_t iMAG   = 1;
282   Int_t iMUON  = 1;
283   Int_t iPHOS  = 1;
284   Int_t iPIPE  = 1;
285   Int_t iPMD   = 1;
286   Int_t iHMPID = 1;
287   Int_t iSHIL  = 1;
288   Int_t iT0    = 1;
289   Int_t iTOF   = 1;
290   Int_t iTPC   = 1;
291   Int_t iTRD   = 1;
292   Int_t iVZERO = 1;
293   Int_t iZDC   = 1;
294
295
296   //=================== Alice BODY parameters =============================
297   AliBODY *BODY = new AliBODY("BODY", "Alice envelop");
298
299
300   if (iMAG)
301     {
302       //=================== MAG parameters ============================
303       // --- Start with Magnet since detector layouts may be depending ---
304       // --- on the selected Magnet dimensions ---
305       AliMAG *MAG = new AliMAG("MAG", "Magnet");
306     }
307
308
309   if (iABSO)
310     {
311       //=================== ABSO parameters ============================
312       AliABSO *ABSO = new AliABSOv3("ABSO", "Muon Absorber");
313     }
314
315   if (iDIPO)
316     {
317       //=================== DIPO parameters ============================
318
319       AliDIPO *DIPO = new AliDIPOv3("DIPO", "Dipole version 3");
320     }
321
322   if (iHALL)
323     {
324       //=================== HALL parameters ============================
325
326       AliHALL *HALL = new AliHALLv3("HALL", "Alice Hall");
327     }
328
329
330   if (iFRAME)
331     {
332       //=================== FRAME parameters ============================
333
334       AliFRAMEv2 *FRAME = new AliFRAMEv2("FRAME", "Space Frame");
335       FRAME->SetHoles(1);
336     }
337
338   if (iSHIL)
339     {
340       //=================== SHIL parameters ============================
341
342       AliSHIL *SHIL = new AliSHILv3("SHIL", "Shielding Version 3");
343     }
344
345
346   if (iPIPE)
347     {
348       //=================== PIPE parameters ============================
349
350       AliPIPE *PIPE = new AliPIPEv3("PIPE", "Beam Pipe");
351     }
352
353   if (iITS)
354     {
355       //=================== ITS parameters ============================
356
357       AliITS *ITS  = new AliITSv11("ITS","ITS v11");
358     }
359
360   if (iTPC)
361     {
362       //============================ TPC parameters =====================
363
364       AliTPC *TPC = new AliTPCv2("TPC", "Default");
365     }
366
367
368   if (iTOF) {
369     //=================== TOF parameters ============================
370
371     AliTOF *TOF = new AliTOFv6T0("TOF", "normal TOF");
372   }
373
374
375   if (iHMPID)
376     {
377       //=================== HMPID parameters ===========================
378
379       AliHMPID *HMPID = new AliHMPIDv3("HMPID", "normal HMPID");
380
381     }
382
383
384   if (iZDC)
385     {
386       //=================== ZDC parameters ============================
387
388       AliZDC *ZDC = new AliZDCv3("ZDC", "normal ZDC");
389     }
390
391   if (iTRD)
392     {
393       //=================== TRD parameters ============================
394
395       AliTRD *TRD = new AliTRDv1("TRD", "TRD slow simulator");
396       AliTRDgeometry *geoTRD = TRD->GetGeometry();
397       // Partial geometry: modules at 0,1,7,8,9,16,17
398       // starting at 3h in positive direction
399       geoTRD->SetSMstatus(2,0);
400       geoTRD->SetSMstatus(3,0);
401       geoTRD->SetSMstatus(4,0);
402       geoTRD->SetSMstatus(5,0);
403       geoTRD->SetSMstatus(6,0);
404       geoTRD->SetSMstatus(11,0);
405       geoTRD->SetSMstatus(12,0);
406       geoTRD->SetSMstatus(13,0);
407       geoTRD->SetSMstatus(14,0);
408       geoTRD->SetSMstatus(15,0);
409       geoTRD->SetSMstatus(16,0);
410     }
411
412   if (iFMD)
413     {
414       //=================== FMD parameters ============================
415
416       AliFMD *FMD = new AliFMDv1("FMD", "normal FMD");
417     }
418
419   if (iMUON)
420     {
421       //=================== MUON parameters ===========================
422       // New MUONv1 version (geometry defined via builders)
423
424       AliMUON *MUON = new AliMUONv1("MUON", "default");
425     }
426
427   if (iPHOS)
428     {
429       //=================== PHOS parameters ===========================
430
431       AliPHOS *PHOS = new AliPHOSv1("PHOS", "noCPV_Modules123");
432
433     }
434
435
436   if (iPMD)
437     {
438       //=================== PMD parameters ============================
439
440       AliPMD *PMD = new AliPMDv1("PMD", "normal PMD");
441     }
442
443   if (iT0)
444     {
445       //=================== T0 parameters ============================
446       AliT0 *T0 = new AliT0v1("T0", "T0 Detector");
447     }
448
449   if (iEMCAL)
450     {
451       //=================== EMCAL parameters ============================
452
453       AliEMCAL *EMCAL = new AliEMCALv2("EMCAL", "EMCAL_FIRSTYEAR");
454     }
455
456   if (iACORDE)
457     {
458       //=================== ACORDE parameters ============================
459
460       AliACORDE *ACORDE = new AliACORDEv1("ACORDE", "normal ACORDE");
461     }
462
463   if (iVZERO)
464     {
465       //=================== ACORDE parameters ============================
466
467       AliVZERO *VZERO = new AliVZEROv7("VZERO", "normal VZERO");
468     }
469 }
470
471 void ProcessEnvironmentVars()
472 {
473   // Run type
474   if (gSystem->Getenv("CONFIG_RUN_TYPE")) {
475     for (Int_t iRun = 0; iRun < kRunMax; iRun++) {
476       if (strcmp(gSystem->Getenv("CONFIG_RUN_TYPE"), pprRunName[iRun])==0) {
477         proc = (PDC06Proc_t)iRun;
478         cout<<"Run type set to "<<pprRunName[iRun]<<endl;
479       }
480     }
481   }
482   // Field
483   if (gSystem->Getenv("CONFIG_FIELD")) {
484     for (Int_t iField = 0; iField < kFieldMax; iField++) {
485       if (strcmp(gSystem->Getenv("CONFIG_FIELD"), pprField[iField])==0) {
486         mag = (Mag_t)iField;
487         cout<<"Field set to "<<pprField[iField]<<endl;
488       }
489     }
490   }
491   // Energy
492   if (gSystem->Getenv("CONFIG_ENERGY")) {
493     energy = atoi(gSystem->Getenv("CONFIG_ENERGY"));
494     cout<<"Energy set to "<<energy<<" GeV"<<endl;
495   }
496   // Random Number seed
497   if (gSystem->Getenv("CONFIG_SEED")) {
498     seed = atoi(gSystem->Getenv("CONFIG_SEED"));
499   }
500   // Run number
501   if (gSystem->Getenv("DC_RUN")) {
502     runNumber = atoi(gSystem->Getenv("DC_RUN"));
503   }
504   // Event offset
505   if (gSystem->Getenv("DC_EVENT")) {
506     eventOffset = 400*(atoi(gSystem->Getenv("DC_EVENT")) - 1); // 400 events per job (->sim.C)
507     std::cout << "eventOffset= " << eventOffset << std::endl;
508   }
509 }
510