]> git.uio.no Git - u/mrichter/AliRoot.git/blob - macros/Config_AliGenCosmicsParam.C
### files: AliTPCTempMap.h (.cxx)
[u/mrichter/AliRoot.git] / macros / Config_AliGenCosmicsParam.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 <TRandom.h>
10 #include <TSystem.h>
11 #include <TVirtualMC.h>
12 #include <TGeant3TGeo.h>
13 #include "STEER/AliRunLoader.h"
14 #include "STEER/AliRun.h"
15 #include "STEER/AliConfig.h"
16 #include "PYTHIA6/AliDecayerPythia.h"
17 #include "EVGEN/AliGenCocktail.h"
18 #include "EVGEN/AliGenCosmicsParam.h"
19 #include "EVGEN/AliGenHIJINGpara.h"
20 #include "STEER/AliMagFMaps.h"
21 #include "STRUCT/AliBODY.h"
22 #include "STRUCT/AliMAG.h"
23 #include "STRUCT/AliABSOv0.h"
24 #include "STRUCT/AliDIPOv2.h"
25 #include "STRUCT/AliHALL.h"
26 #include "STRUCT/AliFRAMEv2.h"
27 #include "STRUCT/AliSHILv2.h"
28 #include "STRUCT/AliPIPEv0.h"
29 #include "ITS/AliITSv11Hybrid.h"
30 #include "TPC/AliTPCv2.h"
31 #include "TOF/AliTOFv5T0.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 "START/AliSTARTv1.h"
40 #include "EMCAL/AliEMCALv2.h"
41 #include "ACORDE/AliACORDEv1.h"
42 #include "VZERO/AliVZEROv7.h"
43 #endif
44
45 //--- Magnetic Field ---
46 enum Mag_t
47 {
48     k2kG, k4kG, k5kG
49 };
50
51 Float_t EtaToTheta(Float_t arg);
52
53 static Mag_t         mag      = k5kG; 
54
55 void Config()
56 {
57     // ThetaRange is (0., 180.). It was (0.28,179.72) 7/12/00 09:00
58     // Theta range given through pseudorapidity limits 22/6/2001
59
60     // Set Random Number seed
61   TDatime dt;
62   UInt_t curtime=dt.Get();
63   UInt_t procid=gSystem->GetPid();
64   UInt_t seed=curtime-procid;
65   if (gSystem->Getenv("envevno")) {
66     seed=atoi(gSystem->Getenv("envevno"));
67     seed += 1000; // 0 e' il seed dall'orologio....
68     printf("...taking seed as event number + 1000...\n");
69   }
70   printf("...setting seed as %d...\n",seed);
71   gRandom->SetSeed(seed);
72   printf("Seed for random number generation = %d \n",gRandom->GetSeed());
73
74
75   //  AliLog::Message(AliLog::kInfo, Form("Seed for random number generation = %d",gRandom->GetSeed()), "Config.C", "Config.C", "Config()","Config.C", __LINE__);
76
77
78    // libraries required by geant321
79 #if defined(__CINT__)
80     gSystem->Load("libgeant321");
81 #endif
82
83     new     TGeant3TGeo("C++ Interface to Geant3");
84
85     if(!AliCDBManager::Instance()->IsDefaultStorageSet()){
86       AliCDBManager::Instance()->SetDefaultStorage("local://$ALICE_ROOT");
87       AliCDBManager::Instance()->SetRun(0);
88     }
89
90     AliRunLoader* rl=0x0;
91
92     AliLog::Message(AliLog::kInfo, "Creating Run Loader", "Config.C", "Config.C", "Config()"," Config.C", __LINE__);
93
94     rl = AliRunLoader::Open("galice.root",
95                             AliConfig::GetDefaultEventFolderName(),
96                             "recreate");
97     if (rl == 0x0)
98       {
99         gAlice->Fatal("Config.C","Can not instatiate the Run Loader");
100         return;
101       }
102     rl->SetCompressionLevel(2);
103     rl->SetNumberOfEventsPerFile(1000);
104     gAlice->SetRunLoader(rl);
105     // gAlice->SetGeometryFromFile("geometry.root");
106     // gAlice->SetGeometryFromCDB();
107
108     // Set the trigger configuration
109     //    gAlice->SetTriggerDescriptor("Pb-Pb");
110     gAlice->SetTriggerDescriptor("ACORDE");
111     cout<<"Trigger configuration is set to  ACORDE"<<endl;
112
113     //
114     // Set External decayer
115     TVirtualMCDecayer *decayer = new AliDecayerPythia();
116
117     decayer->SetForceDecay(kAll);
118     decayer->Init();
119     gMC->SetExternalDecayer(decayer);
120     //=======================================================================
121     // ************* STEERING parameters FOR ALICE SIMULATION **************
122     // --- Specify event type to be tracked through the ALICE setup
123     // --- All positions are in cm, angles in degrees, and P and E in GeV
124
125
126     gMC->SetProcess("DCAY",1);
127     gMC->SetProcess("PAIR",1);
128     gMC->SetProcess("COMP",1);
129     gMC->SetProcess("PHOT",1);
130     gMC->SetProcess("PFIS",0);
131     gMC->SetProcess("DRAY",0);
132     gMC->SetProcess("ANNI",1);
133     gMC->SetProcess("BREM",1);
134     gMC->SetProcess("MUNU",1);
135     gMC->SetProcess("CKOV",1);
136     gMC->SetProcess("HADR",1);
137     gMC->SetProcess("LOSS",2);
138     gMC->SetProcess("MULS",1);
139     gMC->SetProcess("RAYL",1);
140
141     Float_t cut = 1.e-3;        // 1MeV cut by default
142     Float_t tofmax = 1.e10;
143
144     gMC->SetCut("CUTGAM", cut);
145     gMC->SetCut("CUTELE", cut);
146     gMC->SetCut("CUTNEU", cut);
147     gMC->SetCut("CUTHAD", cut);
148     gMC->SetCut("CUTMUO", cut);
149     gMC->SetCut("BCUTE",  cut); 
150     gMC->SetCut("BCUTM",  cut); 
151     gMC->SetCut("DCUTE",  cut); 
152     gMC->SetCut("DCUTM",  cut); 
153     gMC->SetCut("PPCUTM", cut);
154     gMC->SetCut("TOFMAX", tofmax); 
155
156
157     //  AliGenCosmicsParam
158     AliGenCosmicsParam *gener = new AliGenCosmicsParam();
159     gener->SetParamACORDE();
160     gener->SetYOrigin(260.); // warning: just above TPC, no TOF, no ACORDE
161     gener->SetMomentumRange(8.,1000.);
162     gener->SetMaxAngleWRTVertical(3.1415/4.);
163     gener->SetInSPDinner(); // "acceptance trigger"
164     gener->SetBkG(0.); // needed for "acceptance trigger"
165     gener->Init();
166     gGener = gener;
167
168     // MAGNETIC FIELD IN THE BARREL
169     // AliMagFMaps* field = new AliMagFMaps("Maps","Maps", 2, 1., 10., mag); // FIELD
170         AliMagFMaps* field = new AliMagFMaps("Maps","Maps", 2, 0., 10., mag); // NO FIELD
171     field->SetL3ConstField(0); //Using const. field in the barrel
172     rl->CdGAFile();
173     gAlice->SetField(field);    
174
175
176     Int_t   iABSO  =  1;
177     Int_t   iDIPO  =  0;
178     Int_t   iFMD   =  0;
179     Int_t   iFRAME =  1;
180     Int_t   iHALL  =  1;
181     Int_t   iITS   =  1;
182     Int_t   iMAG   =  1;
183     Int_t   iMUON  =  0;
184     Int_t   iPHOS  =  0;
185     Int_t   iPIPE  =  1;
186     Int_t   iPMD   =  0;
187     Int_t   iHMPID =  0;
188     Int_t   iSHIL  =  1;
189     Int_t   iSTART =  0;
190     Int_t   iTOF   =  0;
191     Int_t   iTPC   =  1;
192     Int_t   iTRD   =  0;
193     Int_t   iZDC   =  0;
194     Int_t   iEMCAL =  0;
195     Int_t   iACORDE   =  1;
196     Int_t   iVZERO =  0;
197     rl->CdGAFile();
198     //=================== Alice BODY parameters =============================
199     AliBODY *BODY = new AliBODY("BODY", "Alice envelop");
200
201     if (iMAG)
202     {
203         //=================== MAG parameters ============================
204         // --- Start with Magnet since detector layouts may be depending ---
205         // --- on the selected Magnet dimensions ---
206         AliMAG *MAG = new AliMAG("MAG", "Magnet");
207     }
208
209
210     if (iABSO)
211     {
212         //=================== ABSO parameters ============================
213         AliABSO *ABSO = new AliABSOv0("ABSO", "Muon Absorber");
214     }
215
216     if (iDIPO)
217     {
218         //=================== DIPO parameters ============================
219
220         AliDIPO *DIPO = new AliDIPOv2("DIPO", "Dipole version 2");
221     }
222
223     if (iHALL)
224     {
225         //=================== HALL parameters ============================
226
227         AliHALL *HALL = new AliHALL("HALL", "Alice Hall");
228     }
229
230
231     if (iFRAME)
232     {
233         //=================== FRAME parameters ============================
234
235         AliFRAMEv2 *FRAME = new AliFRAMEv2("FRAME", "Space Frame");
236         FRAME->SetHoles(1);
237     }
238
239     if (iSHIL)
240     {
241         //=================== SHIL parameters ============================
242
243         AliSHIL *SHIL = new AliSHILv2("SHIL", "Shielding Version 2");
244     }
245
246
247     if (iPIPE)
248     {
249         //=================== PIPE parameters ============================
250
251         AliPIPE *PIPE = new AliPIPEv0("PIPE", "Beam Pipe");
252     }
253  
254     if (iITS)
255     {
256         //=================== ITS parameters ============================
257
258         AliITS *ITS  = new AliITSv11Hybrid("ITS","ITS v11Hybrid");
259     }
260
261     if (iTPC)
262     {
263         //============================ TPC parameters ===================
264         AliTPC *TPC = new AliTPCv2("TPC", "Default");
265     }
266
267
268     if (iTOF) {
269         //=================== TOF parameters ============================
270         AliTOF *TOF = new AliTOFv5T0("TOF", "normal TOF");
271     }
272
273
274     if (iHMPID)
275     {
276         //=================== HMPID parameters ===========================
277         AliHMPID *HMPID = new AliHMPIDv3("HMPID", "normal HMPID");
278
279     }
280
281
282     if (iZDC)
283     {
284         //=================== ZDC parameters ============================
285
286         AliZDC *ZDC = new AliZDCv3("ZDC", "normal ZDC");
287     }
288
289     if (iTRD)
290     {
291         //=================== TRD parameters ============================
292
293         AliTRD *TRD = new AliTRDv1("TRD", "TRD slow simulator");
294     }
295
296     if (iFMD)
297     {
298         //=================== FMD parameters ============================
299         AliFMD *FMD = new AliFMDv1("FMD", "normal FMD");
300    }
301
302     if (iMUON)
303     {
304         //=================== MUON parameters ===========================
305         // New MUONv1 version (geometry defined via builders)
306         AliMUON *MUON = new AliMUONv1("MUON", "default");
307     }
308     //=================== PHOS parameters ===========================
309
310     if (iPHOS)
311     {
312         AliPHOS *PHOS = new AliPHOSv1("PHOS", "IHEP");
313     }
314
315
316     if (iPMD)
317     {
318         //=================== PMD parameters ============================
319         AliPMD *PMD = new AliPMDv1("PMD", "normal PMD");
320     }
321
322     if (iSTART)
323     {
324         //=================== START parameters ============================
325         AliSTART *START = new AliSTARTv1("START", "START Detector");
326     }
327
328     if (iEMCAL)
329     {
330         //=================== EMCAL parameters ============================
331         AliEMCAL *EMCAL = new AliEMCALv2("EMCAL", "EMCAL_COMPLETE");
332     }
333
334      if (iACORDE)
335     {
336         //=================== ACORDE parameters ============================
337         AliACORDE *ACORDE = new AliACORDEv1("ACORDE", "normal ACORDE");
338     }
339
340      if (iVZERO)
341     {
342         //=================== ACORDE parameters ============================
343         AliVZERO *VZERO = new AliVZEROv7("VZERO", "normal VZERO");
344     }
345
346      AliLog::Message(AliLog::kInfo, "End of Config", "Config.C", "Config.C", "Config()"," Config.C", __LINE__);
347
348 }
349
350 Float_t EtaToTheta(Float_t arg){
351   return (180./TMath::Pi())*2.*atan(exp(-arg));
352 }