]> git.uio.no Git - u/mrichter/AliRoot.git/blob - macros/Config.C
Updated version of the Geant4 configuration (Andrei)
[u/mrichter/AliRoot.git] / macros / Config.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/AliGenHIJINGpara.h"
19 #include "STEER/AliMagFMaps.h"
20 #include "STRUCT/AliBODY.h"
21 #include "STRUCT/AliMAG.h"
22 #include "STRUCT/AliABSOv3.h"
23 #include "STRUCT/AliDIPOv3.h"
24 #include "STRUCT/AliHALLv3.h"
25 #include "STRUCT/AliFRAMEv2.h"
26 #include "STRUCT/AliSHILv3.h"
27 #include "STRUCT/AliPIPEv3.h"
28 #include "ITS/AliITSvPPRasymmFMD.h"
29 #include "TPC/AliTPCv2.h"
30 #include "TOF/AliTOFv6T0.h"
31 #include "HMPID/AliHMPIDv1.h"
32 #include "ZDC/AliZDCv2.h"
33 #include "TRD/AliTRDv1.h"
34 #include "FMD/AliFMDv1.h"
35 #include "MUON/AliMUONv1.h"
36 #include "PHOS/AliPHOSv1.h"
37 #include "PMD/AliPMDv1.h"
38 #include "T0/AliT0v1.h"
39 #include "EMCAL/AliEMCALv2.h"
40 #include "ACORDE/AliACORDEv0.h"
41 #include "VZERO/AliVZEROv7.h"
42 #endif
43
44 Float_t EtaToTheta(Float_t arg);
45
46 void Config()
47 {
48     // ThetaRange is (0., 180.). It was (0.28,179.72) 7/12/00 09:00
49     // Theta range given through pseudorapidity limits 22/6/2001
50
51     // Set Random Number seed
52   gRandom->SetSeed(123456); // Set 0 to use the currecnt time
53   AliLog::Message(AliLog::kInfo, Form("Seed for random number generation = %d",gRandom->GetSeed()), "Config.C", "Config.C", "Config()","Config.C", __LINE__);
54
55
56    // libraries required by geant321
57 #if defined(__CINT__)
58     gSystem->Load("libgeant321");
59 #endif
60
61     new     TGeant3TGeo("C++ Interface to Geant3");
62
63     AliRunLoader* rl=0x0;
64
65     AliLog::Message(AliLog::kInfo, "Creating Run Loader", "Config.C", "Config.C", "Config()"," Config.C", __LINE__);
66
67     rl = AliRunLoader::Open("galice.root",
68                             AliConfig::GetDefaultEventFolderName(),
69                             "recreate");
70     if (rl == 0x0)
71       {
72         gAlice->Fatal("Config.C","Can not instatiate the Run Loader");
73         return;
74       }
75     rl->SetCompressionLevel(2);
76     rl->SetNumberOfEventsPerFile(3);
77     gAlice->SetRunLoader(rl);
78
79     // Set the trigger configuration
80     gAlice->SetTriggerDescriptor("Pb-Pb");
81     cout<<"Trigger configuration is set to  Pb-Pb"<<endl;
82
83     //
84     // Set External decayer
85     TVirtualMCDecayer *decayer = new AliDecayerPythia();
86
87     decayer->SetForceDecay(kAll);
88     decayer->Init();
89     gMC->SetExternalDecayer(decayer);
90     //=======================================================================
91     // ************* STEERING parameters FOR ALICE SIMULATION **************
92     // --- Specify event type to be tracked through the ALICE setup
93     // --- All positions are in cm, angles in degrees, and P and E in GeV
94
95
96     gMC->SetProcess("DCAY",1);
97     gMC->SetProcess("PAIR",1);
98     gMC->SetProcess("COMP",1);
99     gMC->SetProcess("PHOT",1);
100     gMC->SetProcess("PFIS",0);
101     gMC->SetProcess("DRAY",0);
102     gMC->SetProcess("ANNI",1);
103     gMC->SetProcess("BREM",1);
104     gMC->SetProcess("MUNU",1);
105     gMC->SetProcess("CKOV",1);
106     gMC->SetProcess("HADR",1);
107     gMC->SetProcess("LOSS",2);
108     gMC->SetProcess("MULS",1);
109     gMC->SetProcess("RAYL",1);
110
111     Float_t cut = 1.e-3;        // 1MeV cut by default
112     Float_t tofmax = 1.e10;
113
114     gMC->SetCut("CUTGAM", cut);
115     gMC->SetCut("CUTELE", cut);
116     gMC->SetCut("CUTNEU", cut);
117     gMC->SetCut("CUTHAD", cut);
118     gMC->SetCut("CUTMUO", cut);
119     gMC->SetCut("BCUTE",  cut); 
120     gMC->SetCut("BCUTM",  cut); 
121     gMC->SetCut("DCUTE",  cut); 
122     gMC->SetCut("DCUTM",  cut); 
123     gMC->SetCut("PPCUTM", cut);
124     gMC->SetCut("TOFMAX", tofmax); 
125
126
127     int     nParticles = 100;
128     if (gSystem->Getenv("CONFIG_NPARTICLES"))
129     {
130         nParticles = atoi(gSystem->Getenv("CONFIG_NPARTICLES"));
131     }
132
133     AliGenCocktail *gener = new AliGenCocktail();
134     gener->SetPhiRange(0, 360);
135     // Set pseudorapidity range from -8 to 8.
136     Float_t thmin = EtaToTheta(8);   // theta min. <---> eta max
137     Float_t thmax = EtaToTheta(-8);  // theta max. <---> eta min 
138     gener->SetThetaRange(thmin,thmax);
139     gener->SetOrigin(0, 0, 0);  //vertex position
140     gener->SetSigma(0, 0, 0);   //Sigma in (X,Y,Z) (cm) on IP position
141
142     AliGenHIJINGpara *hijingparam = new AliGenHIJINGpara(nParticles);
143     hijingparam->SetMomentumRange(0.2, 999);
144     gener->AddGenerator(hijingparam,"HIJING PARAM",1);
145
146 //    AliGenBox *genbox = new AliGenBox(nParticles);
147 //    genbox->SetPart(22);
148 //    genbox->SetPtRange(0.3, 10.00);
149 //    gener->AddGenerator(genbox,"GENBOX GAMMA for PHOS",1);
150     gener->Init();
151
152
153     // 
154     // Activate this line if you want the vertex smearing to happen
155     // track by track
156     //
157     //gener->SetVertexSmear(perTrack); 
158     // Field (L3 0.4 T)
159     AliMagFMaps* field = new AliMagFMaps("Maps","Maps", 2, 1., 10., 1);
160     gAlice->SetField(field);    
161
162
163     Int_t   iABSO  =  1;
164     Int_t   iDIPO  =  1;
165     Int_t   iFMD   =  1;
166     Int_t   iFRAME =  1;
167     Int_t   iHALL  =  1;
168     Int_t   iITS   =  1;
169     Int_t   iMAG   =  1;
170     Int_t   iMUON  =  1;
171     Int_t   iPHOS  =  1;
172     Int_t   iPIPE  =  1;
173     Int_t   iPMD   =  1;
174     Int_t   iHMPID  =  1;
175     Int_t   iSHIL  =  1;
176     Int_t   iT0 =  1;
177     Int_t   iTOF   =  1;
178     Int_t   iTPC   =  1;
179     Int_t   iTRD   =  1;
180     Int_t   iZDC   =  1;
181     Int_t   iEMCAL =  1;
182     Int_t   iACORDE   =  0;
183     Int_t   iVZERO =  1;
184     rl->CdGAFile();
185     //=================== Alice BODY parameters =============================
186     AliBODY *BODY = new AliBODY("BODY", "Alice envelop");
187
188     if (iMAG)
189     {
190         //=================== MAG parameters ============================
191         // --- Start with Magnet since detector layouts may be depending ---
192         // --- on the selected Magnet dimensions ---
193         AliMAG *MAG = new AliMAG("MAG", "Magnet");
194     }
195
196
197     if (iABSO)
198     {
199         //=================== ABSO parameters ============================
200         AliABSO *ABSO = new AliABSOv3("ABSO", "Muon Absorber");
201     }
202
203     if (iDIPO)
204     {
205         //=================== DIPO parameters ============================
206
207         AliDIPO *DIPO = new AliDIPOv3("DIPO", "Dipole version 3");
208     }
209
210     if (iHALL)
211     {
212         //=================== HALL parameters ============================
213
214         AliHALL *HALL = new AliHALLv3("HALL", "Alice Hall");
215     }
216
217
218     if (iFRAME)
219     {
220         //=================== FRAME parameters ============================
221
222         AliFRAMEv2 *FRAME = new AliFRAMEv2("FRAME", "Space Frame");
223     }
224
225     if (iSHIL)
226     {
227         //=================== SHIL parameters ============================
228
229         AliSHIL *SHIL = new AliSHILv3("SHIL", "Shielding Version 3");
230     }
231
232
233     if (iPIPE)
234     {
235         //=================== PIPE parameters ============================
236
237         AliPIPE *PIPE = new AliPIPEv3("PIPE", "Beam Pipe");
238     }
239  
240     if(iITS) {
241
242     //=================== ITS parameters ============================
243     //
244     // As the innermost detector in ALICE, the Inner Tracking System "impacts" on
245     // almost all other detectors. This involves the fact that the ITS geometry
246     // still has several options to be followed in parallel in order to determine
247     // the best set-up which minimizes the induced background. All the geometries
248     // available to date are described in the following. Read carefully the comments
249     // and use the default version (the only one uncommented) unless you are making
250     // comparisons and you know what you are doing. In this case just uncomment the
251     // ITS geometry you want to use and run Aliroot.
252     //
253     // Detailed geometries:         
254     //
255     //
256     //
257         AliITSvPPRasymmFMD *ITS  = new AliITSvPPRasymmFMD("ITS","ITS PPR detailed version with asymmetric services");
258         ITS->SetMinorVersion(2);  // don't touch this parameter if you're not an ITS developer
259         ITS->SetReadDet(kFALSE);          // don't touch this parameter if you're not an ITS developer
260         //    ITS->SetWriteDet("$ALICE_ROOT/ITS/ITSgeometry_vPPRasymm2.det");  // don't touch this parameter if you're not an ITS developer
261         ITS->SetThicknessDet1(200.);   // detector thickness on layer 1 must be in the range [100,300]
262         ITS->SetThicknessDet2(200.);   // detector thickness on layer 2 must be in the range [100,300]
263         ITS->SetThicknessChip1(150.);  // chip thickness on layer 1 must be in the range [150,300]
264         ITS->SetThicknessChip2(150.);  // chip thickness on layer 2 must be in the range [150,300]
265         ITS->SetRails(0);              // 1 --> rails in ; 0 --> rails out
266         ITS->SetCoolingFluid(1);       // 1 --> water ; 0 --> freon
267
268  
269     //
270     // Coarse geometries (warning: no hits are produced with these coarse geometries and they unuseful 
271     // for reconstruction !):
272     //                                                     
273     //
274     //AliITSvPPRcoarseasymm *ITS  = new AliITSvPPRcoarseasymm("ITS","New ITS PPR coarse version with asymmetric services");
275     //ITS->SetRails(0);                // 1 --> rails in ; 0 --> rails out
276     //ITS->SetSupportMaterial(0);      // 0 --> Copper ; 1 --> Aluminum ; 2 --> Carbon
277     //
278     //AliITS *ITS  = new AliITSvPPRcoarsesymm("ITS","New ITS PPR coarse version with symmetric services");
279     //ITS->SetRails(0);                // 1 --> rails in ; 0 --> rails out
280     //ITS->SetSupportMaterial(0);      // 0 --> Copper ; 1 --> Aluminum ; 2 --> Carbon
281     //                      
282     //
283     //
284     // Geant3 <-> EUCLID conversion
285     // ============================
286     //
287     // SetEUCLID is a flag to output (=1) or not to output (=0) both geometry and
288     // media to two ASCII files (called by default ITSgeometry.euc and
289     // ITSgeometry.tme) in a format understandable to the CAD system EUCLID.
290     // The default (=0) means that you dont want to use this facility.
291     //
292      ITS->SetEUCLID(0);  
293     }
294
295     if (iTPC)
296     {
297         //============================ TPC parameters ===================
298         AliTPC *TPC = new AliTPCv2("TPC", "Default");
299     }
300
301
302     if (iTOF) {
303         //=================== TOF parameters ============================
304         AliTOF *TOF = new AliTOFv6T0("TOF", "normal TOF");
305     }
306
307
308     if (iHMPID)
309     {
310         //=================== HMPID parameters ===========================
311         AliHMPID *HMPID = new AliHMPIDv1("HMPID", "normal HMPID");
312
313     }
314
315
316     if (iZDC)
317     {
318         //=================== ZDC parameters ============================
319
320         AliZDC *ZDC = new AliZDCv2("ZDC", "normal ZDC");
321     }
322
323     if (iTRD)
324     {
325         //=================== TRD parameters ============================
326
327         AliTRD *TRD = new AliTRDv1("TRD", "TRD slow simulator");
328     }
329
330     if (iFMD)
331     {
332         //=================== FMD parameters ============================
333         AliFMD *FMD = new AliFMDv1("FMD", "normal FMD");
334    }
335
336     if (iMUON)
337     {
338         //=================== MUON parameters ===========================
339         // New MUONv1 version (geometry defined via builders)
340         AliMUON *MUON = new AliMUONv1("MUON", "default");
341     }
342     //=================== PHOS parameters ===========================
343
344     if (iPHOS)
345     {
346         AliPHOS *PHOS = new AliPHOSv1("PHOS", "IHEP");
347     }
348
349
350     if (iPMD)
351     {
352         //=================== PMD parameters ============================
353         AliPMD *PMD = new AliPMDv1("PMD", "normal PMD");
354     }
355
356     if (iT0)
357     {
358         //=================== T0 parameters ============================
359         AliT0 *T0 = new AliT0v1("T0", "T0 Detector");
360     }
361
362     if (iEMCAL)
363     {
364         //=================== EMCAL parameters ============================
365         AliEMCAL *EMCAL = new AliEMCALv2("EMCAL", "SHISH_77_TRD1_2X2_FINAL_110DEG");
366     }
367
368      if (iACORDE)
369     {
370         //=================== ACORDE parameters ============================
371         AliACORDE *ACORDE = new AliACORDEv0("ACORDE", "normal ACORDE");
372     }
373
374      if (iVZERO)
375     {
376         //=================== ACORDE parameters ============================
377         AliVZERO *VZERO = new AliVZEROv7("VZERO", "normal VZERO");
378     }
379
380      AliLog::Message(AliLog::kInfo, "End of Config", "Config.C", "Config.C", "Config()"," Config.C", __LINE__);
381
382 }
383
384 Float_t EtaToTheta(Float_t arg){
385   return (180./TMath::Pi())*2.*atan(exp(-arg));
386 }