]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TFluka/macro/FlukaConfig.C
.so cleanup: more gSystem->Load()
[u/mrichter/AliRoot.git] / TFluka / macro / FlukaConfig.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/AliITSv11.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 Float_t EtaToTheta(Float_t arg);
46 void    LoadPythia();
47
48
49 void Config()
50 {
51     // ThetaRange is (0., 180.). It was (0.28,179.72) 7/12/00 09:00
52     // Theta range given through pseudorapidity limits 22/6/2001
53
54     // Set Random Number seed
55   //gRandom->SetSeed(123456); // Set 0 to use the current time
56   
57   AliLog::Message(AliLog::kInfo, Form("Seed for random number generation = %d",gRandom->GetSeed()), "Config.C", "Config.C", "Config()","Config.C", __LINE__);
58
59   // Load Pythia libraries
60   LoadPythia();
61   // Libraries required by geant321
62   Bool_t isFluka = kTRUE;
63     if (isFluka) {
64       gSystem->Load("libGeom");
65       cout << "\t* Loading TFluka..." << endl;  
66       gSystem->Load("libfluka");    
67       
68       cout << "\t* Instantiating TFluka..." << endl;
69       new  TFluka("C++ Interface to Fluka", 0/*verbositylevel*/);
70     }
71     else {
72       cout << "\t* Loading Geant3..." << endl;  
73       gSystem->Load("libgeant321");
74       
75       cout << "\t* Instantiating Geant3TGeo..." << endl;
76       new     TGeant3TGeo("C++ Interface to Geant3");
77     }
78     AliRunLoader* rl=0x0;
79
80     AliLog::Message(AliLog::kInfo, "Creating Run Loader", "Config.C", "Config.C", "Config()"," Config.C", __LINE__);
81
82     rl = AliRunLoader::Open("galice.root",
83                             AliConfig::GetDefaultEventFolderName(),
84                             "recreate");
85     if (rl == 0x0)
86       {
87         gAlice->Fatal("Config.C","Can not instatiate the Run Loader");
88         return;
89       }
90     rl->SetCompressionLevel(2);
91     rl->SetNumberOfEventsPerFile(3);
92     gAlice->SetRunLoader(rl);
93     
94     // gAlice->SetGeometryFromFile("geometry.root");
95
96     // Uncomment if you want to load geometry from OCDB!   >>>>
97 /*    
98     if(!AliCDBManager::Instance()->IsDefaultStorageSet()){
99          cout << "#####################################################" << endl;
100          cout << "#                                                   #" << endl;
101          cout << "#     WARNING: CDB DEFAULT STORAGE NOT SET !!!      #" << endl;
102          cout << "#     SETTING IT TO local://$ALICE_ROOT/OCDB !!!         #" << endl;
103          cout << "#                                                   #" << endl;
104          cout << "#####################################################" << endl;
105           
106          AliCDBManager::Instance()->SetDefaultStorage("local://$ALICE_ROOT/OCDB");
107     }
108     
109     if(AliCDBManager::Instance()->GetRun() < 0){
110          cout << "#####################################################" << endl;
111          cout << "#                                                   #" << endl;
112          cout << "#     WARNING: RUN NUMBER NOT SET !!!               #" << endl;
113          cout << "#     SETTING IT TO 0 !!!                           #" << endl;
114          cout << "#                                                   #" << endl;
115          cout << "#####################################################" << endl;
116           
117          AliCDBManager::Instance()->SetRun(0);
118     }
119     gAlice->SetGeometryFromCDB();
120 */
121     // Uncomment if you want to load geometry from OCDB!   <<<<
122
123     // Set the trigger configuration
124     AliSimulation::Instance()->SetTriggerConfig("Pb-Pb");
125     cout<<"Trigger configuration is set to  Pb-Pb"<<endl;
126
127     //
128     // Set External decayer
129     TVirtualMCDecayer *decayer = new AliDecayerPythia();
130
131     decayer->SetForceDecay(kAll);
132     decayer->Init();
133     gMC->SetExternalDecayer(decayer);
134     //=======================================================================
135     // ************* STEERING parameters FOR ALICE SIMULATION **************
136     // --- Specify event type to be tracked through the ALICE setup
137     // --- All positions are in cm, angles in degrees, and P and E in GeV
138
139
140     gMC->SetProcess("DCAY",1);
141     gMC->SetProcess("PAIR",1);
142     gMC->SetProcess("COMP",1);
143     gMC->SetProcess("PHOT",1);
144     gMC->SetProcess("PFIS",0);
145     gMC->SetProcess("DRAY",0);
146     gMC->SetProcess("ANNI",1);
147     gMC->SetProcess("BREM",1);
148     gMC->SetProcess("MUNU",1);
149     gMC->SetProcess("CKOV",1);
150     gMC->SetProcess("HADR",1);
151     gMC->SetProcess("LOSS",2);
152     gMC->SetProcess("MULS",1);
153     gMC->SetProcess("RAYL",1);
154
155     Float_t cut = 1.e-3;        // 1MeV cut by default
156     Float_t tofmax = 1.e10;
157
158     gMC->SetCut("CUTGAM", cut);
159     gMC->SetCut("CUTELE", cut);
160     gMC->SetCut("CUTNEU", cut);
161     gMC->SetCut("CUTHAD", cut);
162     gMC->SetCut("CUTMUO", cut);
163     gMC->SetCut("BCUTE",  cut); 
164     gMC->SetCut("BCUTM",  cut); 
165     gMC->SetCut("DCUTE",  cut); 
166     gMC->SetCut("DCUTM",  cut); 
167     gMC->SetCut("PPCUTM", cut);
168     gMC->SetCut("TOFMAX", tofmax); 
169
170
171     int     nParticles = 100;
172     if (gSystem->Getenv("CONFIG_NPARTICLES"))
173     {
174         nParticles = atoi(gSystem->Getenv("CONFIG_NPARTICLES"));
175     }
176
177     AliGenCocktail *gener = new AliGenCocktail();
178     gener->SetPhiRange(0, 360);
179     // Set pseudorapidity range from -8 to 8.
180     Float_t thmin = EtaToTheta(8);   // theta min. <---> eta max
181     Float_t thmax = EtaToTheta(-8);  // theta max. <---> eta min 
182     gener->SetThetaRange(thmin,thmax);
183     gener->SetOrigin(0, 0, 0);  //vertex position
184     gener->SetSigma(0, 0, 0);   //Sigma in (X,Y,Z) (cm) on IP position
185
186     AliGenHIJINGpara *hijingparam = new AliGenHIJINGpara(nParticles);
187     hijingparam->SetMomentumRange(0.2, 999);
188     gener->AddGenerator(hijingparam,"HIJING PARAM",1);
189
190 //    AliGenBox *genbox = new AliGenBox(nParticles);
191 //    genbox->SetPart(kGamma);
192 //    genbox->SetPtRange(0.3, 10.00);
193 //    gener->AddGenerator(genbox,"GENBOX GAMMA for PHOS",1);
194     gener->Init();
195
196
197     // 
198     // Activate this line if you want the vertex smearing to happen
199     // track by track
200     //
201     //gener->SetVertexSmear(perTrack); 
202     // Field (L3 0.4 T)
203     TGeoGlobalMagField::Instance()->SetField(new AliMagF("Maps","Maps", -1., -1., AliMagF::k5kG));
204
205     Int_t   iABSO   = 1;
206     Int_t   iDIPO   = 1;
207     Int_t   iFMD    = 1;
208     Int_t   iFRAME  = 1;
209     Int_t   iHALL   = 1;
210     Int_t   iITS    = 1;
211     Int_t   iMAG    = 1;
212     Int_t   iMUON   = 1;
213     Int_t   iPHOS   = 1;
214     Int_t   iPIPE   = 1;
215     Int_t   iPMD    = 1;
216     Int_t   iHMPID  = 1;
217     Int_t   iSHIL   = 1;
218     Int_t   iT0     = 1;
219     Int_t   iTOF    = 1;
220     Int_t   iTPC    = 1;
221     Int_t   iTRD    = 1;
222     Int_t   iZDC    = 1;
223     Int_t   iEMCAL  = 1;
224     Int_t   iACORDE = 1;
225     Int_t   iVZERO  = 1;
226     rl->CdGAFile();
227     //=================== Alice BODY parameters =============================
228     AliBODY *BODY = new AliBODY("BODY", "Alice envelop");
229
230     if (iMAG)
231     {
232         //=================== MAG parameters ============================
233         // --- Start with Magnet since detector layouts may be depending ---
234         // --- on the selected Magnet dimensions ---
235         AliMAG *MAG = new AliMAG("MAG", "Magnet");
236     }
237
238
239     if (iABSO)
240     {
241         //=================== ABSO parameters ============================
242         AliABSO *ABSO = new AliABSOv3("ABSO", "Muon Absorber");
243     }
244
245     if (iDIPO)
246     {
247         //=================== DIPO parameters ============================
248
249         AliDIPO *DIPO = new AliDIPOv3("DIPO", "Dipole version 3");
250     }
251
252     if (iHALL)
253     {
254         //=================== HALL parameters ============================
255
256         AliHALL *HALL = new AliHALLv3("HALL", "Alice Hall");
257     }
258
259
260     if (iFRAME)
261     {
262         //=================== FRAME parameters ============================
263
264         AliFRAMEv2 *FRAME = new AliFRAMEv2("FRAME", "Space Frame");
265         FRAME->SetHoles(1);
266     }
267
268     if (iSHIL)
269     {
270         //=================== SHIL parameters ============================
271
272         AliSHIL *SHIL = new AliSHILv3("SHIL", "Shielding Version 3");
273     }
274
275
276     if (iPIPE)
277     {
278         //=================== PIPE parameters ============================
279
280         AliPIPE *PIPE = new AliPIPEv3("PIPE", "Beam Pipe");
281     }
282  
283     if (iITS)
284     {
285         //=================== ITS parameters ============================
286
287         AliITS *ITS  = new AliITSv11("ITS","ITS v11");
288     }
289
290     if (iTPC)
291     {
292         //============================ TPC parameters ===================
293         AliTPC *TPC = new AliTPCv2("TPC", "Default");
294         TPC->SetPrimaryIonisation();
295     }
296
297
298     if (iTOF) {
299         //=================== TOF parameters ============================
300         AliTOF *TOF = new AliTOFv6T0("TOF", "normal TOF");
301     }
302
303
304     if (iHMPID)
305     {
306         //=================== HMPID parameters ===========================
307         AliHMPID *HMPID = new AliHMPIDv3("HMPID", "normal HMPID");
308
309     }
310
311
312     if (iZDC)
313     {
314         //=================== ZDC parameters ============================
315
316         AliZDC *ZDC = new AliZDCv3("ZDC", "normal ZDC");
317     }
318
319     if (iTRD)
320     {
321         //=================== TRD parameters ============================
322
323         AliTRD *TRD = new AliTRDv1("TRD", "TRD slow simulator");
324     }
325
326     if (iFMD)
327     {
328         //=================== FMD parameters ============================
329         AliFMD *FMD = new AliFMDv1("FMD", "normal FMD");
330    }
331
332     if (iMUON)
333     {
334         //=================== MUON parameters ===========================
335         // New MUONv1 version (geometry defined via builders)
336         AliMUON *MUON = new AliMUONv1("MUON", "default");
337     }
338     //=================== PHOS parameters ===========================
339
340     if (iPHOS)
341     {
342         AliPHOS *PHOS = new AliPHOSv1("PHOS", "IHEP");
343     }
344
345
346     if (iPMD)
347     {
348         //=================== PMD parameters ============================
349         AliPMD *PMD = new AliPMDv1("PMD", "normal PMD");
350     }
351
352     if (iT0)
353     {
354         //=================== T0 parameters ============================
355         AliT0 *T0 = new AliT0v1("T0", "T0 Detector");
356     }
357
358     if (iEMCAL)
359     {
360         //=================== EMCAL parameters ============================
361         AliEMCAL *EMCAL = new AliEMCALv2("EMCAL", "EMCAL_COMPLETE");
362     }
363
364      if (iACORDE)
365     {
366         //=================== ACORDE parameters ============================
367         AliACORDE *ACORDE = new AliACORDEv1("ACORDE", "normal ACORDE");
368     }
369
370      if (iVZERO)
371     {
372         //=================== VZERO parameters ============================
373         AliVZERO *VZERO = new AliVZEROv7("VZERO", "normal VZERO");
374     }
375
376      AliLog::Message(AliLog::kInfo, "End of Config", "Config.C", "Config.C", "Config()"," Config.C", __LINE__);
377
378 }
379
380 Float_t EtaToTheta(Float_t arg){
381   return (180./TMath::Pi())*2.*atan(exp(-arg));
382 }
383
384
385 void LoadPythia()
386 {
387     // Load Pythia related libraries
388     gSystem->Load("liblhapdf");      // Parton density functions
389     gSystem->Load("libEGPythia6");   // TGenerator interface
390     gSystem->Load("libpythia6");     // Pythia
391     gSystem->Load("libAliPythia6");  // ALICE specific implementations
392 }