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