]> git.uio.no Git - u/mrichter/AliRoot.git/blame_incremental - macros/Config.C
TPC default configuration
[u/mrichter/AliRoot.git] / macros / Config.C
... / ...
CommitLineData
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/AliHMPIDv2.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
44Float_t EtaToTheta(Float_t arg);
45
46void 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 if(!AliCDBManager::Instance()->IsDefaultStorageSet()){
64 AliCDBManager::Instance()->SetDefaultStorage("local://$ALICE_ROOT");
65 AliCDBManager::Instance()->SetRun(0);
66 }
67
68 AliRunLoader* rl=0x0;
69
70 AliLog::Message(AliLog::kInfo, "Creating Run Loader", "Config.C", "Config.C", "Config()"," Config.C", __LINE__);
71
72 rl = AliRunLoader::Open("galice.root",
73 AliConfig::GetDefaultEventFolderName(),
74 "recreate");
75 if (rl == 0x0)
76 {
77 gAlice->Fatal("Config.C","Can not instatiate the Run Loader");
78 return;
79 }
80 rl->SetCompressionLevel(2);
81 rl->SetNumberOfEventsPerFile(3);
82 gAlice->SetRunLoader(rl);
83 // gAlice->SetGeometryFromFile("geometry.root");
84 // gAlice->SetGeometryFromCDB();
85
86 // Set the trigger configuration
87 gAlice->SetTriggerDescriptor("Pb-Pb");
88 cout<<"Trigger configuration is set to Pb-Pb"<<endl;
89
90 //
91 // Set External decayer
92 TVirtualMCDecayer *decayer = new AliDecayerPythia();
93
94 decayer->SetForceDecay(kAll);
95 decayer->Init();
96 gMC->SetExternalDecayer(decayer);
97 //=======================================================================
98 // ************* STEERING parameters FOR ALICE SIMULATION **************
99 // --- Specify event type to be tracked through the ALICE setup
100 // --- All positions are in cm, angles in degrees, and P and E in GeV
101
102
103 gMC->SetProcess("DCAY",1);
104 gMC->SetProcess("PAIR",1);
105 gMC->SetProcess("COMP",1);
106 gMC->SetProcess("PHOT",1);
107 gMC->SetProcess("PFIS",0);
108 gMC->SetProcess("DRAY",0);
109 gMC->SetProcess("ANNI",1);
110 gMC->SetProcess("BREM",1);
111 gMC->SetProcess("MUNU",1);
112 gMC->SetProcess("CKOV",1);
113 gMC->SetProcess("HADR",1);
114 gMC->SetProcess("LOSS",2);
115 gMC->SetProcess("MULS",1);
116 gMC->SetProcess("RAYL",1);
117
118 Float_t cut = 1.e-3; // 1MeV cut by default
119 Float_t tofmax = 1.e10;
120
121 gMC->SetCut("CUTGAM", cut);
122 gMC->SetCut("CUTELE", cut);
123 gMC->SetCut("CUTNEU", cut);
124 gMC->SetCut("CUTHAD", cut);
125 gMC->SetCut("CUTMUO", cut);
126 gMC->SetCut("BCUTE", cut);
127 gMC->SetCut("BCUTM", cut);
128 gMC->SetCut("DCUTE", cut);
129 gMC->SetCut("DCUTM", cut);
130 gMC->SetCut("PPCUTM", cut);
131 gMC->SetCut("TOFMAX", tofmax);
132
133
134 int nParticles = 100;
135 if (gSystem->Getenv("CONFIG_NPARTICLES"))
136 {
137 nParticles = atoi(gSystem->Getenv("CONFIG_NPARTICLES"));
138 }
139
140 AliGenCocktail *gener = new AliGenCocktail();
141 gener->SetPhiRange(0, 360);
142 // Set pseudorapidity range from -8 to 8.
143 Float_t thmin = EtaToTheta(8); // theta min. <---> eta max
144 Float_t thmax = EtaToTheta(-8); // theta max. <---> eta min
145 gener->SetThetaRange(thmin,thmax);
146 gener->SetOrigin(0, 0, 0); //vertex position
147 gener->SetSigma(0, 0, 0); //Sigma in (X,Y,Z) (cm) on IP position
148
149 AliGenHIJINGpara *hijingparam = new AliGenHIJINGpara(nParticles);
150 hijingparam->SetMomentumRange(0.2, 999);
151 gener->AddGenerator(hijingparam,"HIJING PARAM",1);
152
153// AliGenBox *genbox = new AliGenBox(nParticles);
154// genbox->SetPart(22);
155// genbox->SetPtRange(0.3, 10.00);
156// gener->AddGenerator(genbox,"GENBOX GAMMA for PHOS",1);
157 gener->Init();
158
159
160 //
161 // Activate this line if you want the vertex smearing to happen
162 // track by track
163 //
164 //gener->SetVertexSmear(perTrack);
165 // Field (L3 0.4 T)
166 AliMagFMaps* field = new AliMagFMaps("Maps","Maps", 2, 1., 10., 1);
167 gAlice->SetField(field);
168
169
170 Int_t iABSO = 1;
171 Int_t iDIPO = 1;
172 Int_t iFMD = 1;
173 Int_t iFRAME = 1;
174 Int_t iHALL = 1;
175 Int_t iITS = 1;
176 Int_t iMAG = 1;
177 Int_t iMUON = 1;
178 Int_t iPHOS = 1;
179 Int_t iPIPE = 1;
180 Int_t iPMD = 1;
181 Int_t iHMPID = 1;
182 Int_t iSHIL = 1;
183 Int_t iT0 = 1;
184 Int_t iTOF = 1;
185 Int_t iTPC = 1;
186 Int_t iTRD = 1;
187 Int_t iZDC = 1;
188 Int_t iEMCAL = 1;
189 Int_t iACORDE = 0;
190 Int_t iVZERO = 1;
191 rl->CdGAFile();
192 //=================== Alice BODY parameters =============================
193 AliBODY *BODY = new AliBODY("BODY", "Alice envelop");
194
195 if (iMAG)
196 {
197 //=================== MAG parameters ============================
198 // --- Start with Magnet since detector layouts may be depending ---
199 // --- on the selected Magnet dimensions ---
200 AliMAG *MAG = new AliMAG("MAG", "Magnet");
201 }
202
203
204 if (iABSO)
205 {
206 //=================== ABSO parameters ============================
207 AliABSO *ABSO = new AliABSOv3("ABSO", "Muon Absorber");
208 }
209
210 if (iDIPO)
211 {
212 //=================== DIPO parameters ============================
213
214 AliDIPO *DIPO = new AliDIPOv3("DIPO", "Dipole version 3");
215 }
216
217 if (iHALL)
218 {
219 //=================== HALL parameters ============================
220
221 AliHALL *HALL = new AliHALLv3("HALL", "Alice Hall");
222 }
223
224
225 if (iFRAME)
226 {
227 //=================== FRAME parameters ============================
228
229 AliFRAMEv2 *FRAME = new AliFRAMEv2("FRAME", "Space Frame");
230 }
231
232 if (iSHIL)
233 {
234 //=================== SHIL parameters ============================
235
236 AliSHIL *SHIL = new AliSHILv3("SHIL", "Shielding Version 3");
237 }
238
239
240 if (iPIPE)
241 {
242 //=================== PIPE parameters ============================
243
244 AliPIPE *PIPE = new AliPIPEv3("PIPE", "Beam Pipe");
245 }
246
247 if (iITS)
248 {
249 //=================== ITS parameters ============================
250
251 AliITSvPPRasymmFMD *ITS = new AliITSvPPRasymmFMD("ITS","ITS PPR detailed version with asymmetric services");
252 }
253
254 if (iTPC)
255 {
256 //============================ TPC parameters ===================
257 AliTPC *TPC = new AliTPCv2("TPC", "Default");
258 }
259
260
261 if (iTOF) {
262 //=================== TOF parameters ============================
263 AliTOF *TOF = new AliTOFv6T0("TOF", "normal TOF");
264 }
265
266
267 if (iHMPID)
268 {
269 //=================== HMPID parameters ===========================
270 AliHMPID *HMPID = new AliHMPIDv2("HMPID", "normal HMPID");
271
272 }
273
274
275 if (iZDC)
276 {
277 //=================== ZDC parameters ============================
278
279 AliZDC *ZDC = new AliZDCv2("ZDC", "normal ZDC");
280 }
281
282 if (iTRD)
283 {
284 //=================== TRD parameters ============================
285
286 AliTRD *TRD = new AliTRDv1("TRD", "TRD slow simulator");
287 }
288
289 if (iFMD)
290 {
291 //=================== FMD parameters ============================
292 AliFMD *FMD = new AliFMDv1("FMD", "normal FMD");
293 }
294
295 if (iMUON)
296 {
297 //=================== MUON parameters ===========================
298 // New MUONv1 version (geometry defined via builders)
299 AliMUON *MUON = new AliMUONv1("MUON", "default");
300 }
301 //=================== PHOS parameters ===========================
302
303 if (iPHOS)
304 {
305 AliPHOS *PHOS = new AliPHOSv1("PHOS", "IHEP");
306 }
307
308
309 if (iPMD)
310 {
311 //=================== PMD parameters ============================
312 AliPMD *PMD = new AliPMDv1("PMD", "normal PMD");
313 }
314
315 if (iT0)
316 {
317 //=================== T0 parameters ============================
318 AliT0 *T0 = new AliT0v1("T0", "T0 Detector");
319 }
320
321 if (iEMCAL)
322 {
323 //=================== EMCAL parameters ============================
324 AliEMCAL *EMCAL = new AliEMCALv2("EMCAL", "SHISH_77_TRD1_2X2_FINAL_110DEG");
325 }
326
327 if (iACORDE)
328 {
329 //=================== ACORDE parameters ============================
330 AliACORDE *ACORDE = new AliACORDEv0("ACORDE", "normal ACORDE");
331 }
332
333 if (iVZERO)
334 {
335 //=================== ACORDE parameters ============================
336 AliVZERO *VZERO = new AliVZEROv7("VZERO", "normal VZERO");
337 }
338
339 AliLog::Message(AliLog::kInfo, "End of Config", "Config.C", "Config.C", "Config()"," Config.C", __LINE__);
340
341}
342
343Float_t EtaToTheta(Float_t arg){
344 return (180./TMath::Pi())*2.*atan(exp(-arg));
345}