]> git.uio.no Git - u/mrichter/AliRoot.git/blame - TPC/testMC/ConfigHPT1.C
Use PDG code incude
[u/mrichter/AliRoot.git] / TPC / testMC / ConfigHPT1.C
CommitLineData
75c56124 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(0); // 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
141 AliGenCocktail *gener = new AliGenCocktail();
142 gener->SetPhiRange(0, 360);
143 // Set pseudorapidity range from -8 to 8.
144 Float_t thmin = EtaToTheta(1.2); // theta min. <---> eta max
145 Float_t thmax = EtaToTheta(-1.2); // theta max. <---> eta min
146 gener->SetThetaRange(thmin,thmax);
147 gener->SetOrigin(0, 0, 0); //vertex position
148 gener->SetSigma(0, 0, 0); //Sigma in (X,Y,Z) (cm) on IP position
149
150 // AliGenHIJINGpara *hijingparam = new AliGenHIJINGpara(nParticles);
151 //hijingparam->SetMomentumRange(0.2, 999);
152 //gener->AddGenerator(hijingparam,"HIJING PARAM",1);
153 //
154 //
155 Int_t nParticles2 =1;
156 // mu
157 AliGenBox *genboxMP = new AliGenBox(nParticles2);
158 genboxMP->SetPart(13);
159 genboxMP->SetPtRange(10, 100.00);
160
161
162 gener->AddGenerator(genboxMP,"GENBOX",1);
163
164 gener->Init();
165
166
167
168 //
169 // Activate this line if you want the vertex smearing to happen
170 // track by track
171 //
172 //gener->SetVertexSmear(perTrack);
173 // Field (L3 0.4 T)
174 AliMagFMaps* field = new AliMagFMaps("Maps","Maps", 2, 1., 10., 1);
175 gAlice->SetField(field);
176
177
178 Int_t iABSO = 1;
179 Int_t iDIPO = 1;
180 Int_t iFMD = 0;
181 Int_t iFRAME = 1;
182 Int_t iHALL = 1;
183 Int_t iITS = 1;
184 Int_t iMAG = 1;
185 Int_t iMUON = 0;
186 Int_t iPHOS = 0;
187 Int_t iPIPE = 1;
188 Int_t iPMD = 0;
189 Int_t iHMPID = 0;
190 Int_t iSHIL = 1;
191 Int_t iT0 = 1;
192 Int_t iTOF = 1;
193 Int_t iTPC = 1;
194 Int_t iTRD = 1;
195 Int_t iZDC = 1;
196 Int_t iEMCAL = 0;
197 Int_t iACORDE = 0;
198 Int_t iVZERO = 0;
199 rl->CdGAFile();
200 //=================== Alice BODY parameters =============================
201 AliBODY *BODY = new AliBODY("BODY", "Alice envelop");
202
203 if (iMAG)
204 {
205 //=================== MAG parameters ============================
206 // --- Start with Magnet since detector layouts may be depending ---
207 // --- on the selected Magnet dimensions ---
208 AliMAG *MAG = new AliMAG("MAG", "Magnet");
209 }
210
211
212 if (iABSO)
213 {
214 //=================== ABSO parameters ============================
215 AliABSO *ABSO = new AliABSOv3("ABSO", "Muon Absorber");
216 }
217
218 if (iDIPO)
219 {
220 //=================== DIPO parameters ============================
221
222 AliDIPO *DIPO = new AliDIPOv3("DIPO", "Dipole version 3");
223 }
224
225 if (iHALL)
226 {
227 //=================== HALL parameters ============================
228
229 AliHALL *HALL = new AliHALLv3("HALL", "Alice Hall");
230 }
231
232
233 if (iFRAME)
234 {
235 //=================== FRAME parameters ============================
236
237 AliFRAMEv2 *FRAME = new AliFRAMEv2("FRAME", "Space Frame");
238 }
239
240 if (iSHIL)
241 {
242 //=================== SHIL parameters ============================
243
244 AliSHIL *SHIL = new AliSHILv3("SHIL", "Shielding Version 3");
245 }
246
247
248 if (iPIPE)
249 {
250 //=================== PIPE parameters ============================
251
252 AliPIPE *PIPE = new AliPIPEv3("PIPE", "Beam Pipe");
253 }
254
255 if (iITS)
256 {
257 //=================== ITS parameters ============================
258 AliITSvPPRasymmFMD *ITS = new AliITSvPPRasymmFMD("ITS",
259 "ITS PPR detailed version with asymmetric services");
260 }
261
262 if (iTPC)
263 {
264 //============================ TPC parameters ===================
265 AliTPC *TPC = new AliTPCv2("TPC", "Default");
266 }
267
268
269 if (iTOF) {
270 //=================== TOF parameters ============================
271 AliTOF *TOF = new AliTOFv6T0("TOF", "normal TOF");
272 }
273
274
275 if (iHMPID)
276 {
277 //=================== HMPID parameters ===========================
278 AliHMPID *HMPID = new AliHMPIDv2("HMPID", "normal HMPID");
279
280 }
281
282
283 if (iZDC)
284 {
285 //=================== ZDC parameters ============================
286
287 AliZDC *ZDC = new AliZDCv2("ZDC", "normal ZDC");
288 }
289
290 if (iTRD)
291 {
292 //=================== TRD parameters ============================
293
294 AliTRD *TRD = new AliTRDv1("TRD", "TRD slow simulator");
295 }
296
297 if (iFMD)
298 {
299 //=================== FMD parameters ============================
300 AliFMD *FMD = new AliFMDv1("FMD", "normal FMD");
301 }
302
303 if (iMUON)
304 {
305 //=================== MUON parameters ===========================
306 // New MUONv1 version (geometry defined via builders)
307 AliMUON *MUON = new AliMUONv1("MUON", "default");
308 }
309 //=================== PHOS parameters ===========================
310
311 if (iPHOS)
312 {
313 AliPHOS *PHOS = new AliPHOSv1("PHOS", "IHEP");
314 }
315
316
317 if (iPMD)
318 {
319 //=================== PMD parameters ============================
320 AliPMD *PMD = new AliPMDv1("PMD", "normal PMD");
321 }
322
323 if (iT0)
324 {
325 //=================== T0 parameters ============================
326 AliT0 *T0 = new AliT0v1("T0", "T0 Detector");
327 }
328
329 if (iEMCAL)
330 {
331 //=================== EMCAL parameters ============================
332 AliEMCAL *EMCAL = new AliEMCALv2("EMCAL", "SHISH_77_TRD1_2X2_FINAL_110DEG");
333 }
334
335 if (iACORDE)
336 {
337 //=================== ACORDE parameters ============================
338 AliACORDE *ACORDE = new AliACORDEv0("ACORDE", "normal ACORDE");
339 }
340
341 if (iVZERO)
342 {
343 //=================== ACORDE parameters ============================
344 AliVZERO *VZERO = new AliVZEROv7("VZERO", "normal VZERO");
345 }
346
347 AliLog::Message(AliLog::kInfo, "End of Config", "Config.C", "Config.C", "Config()"," Config.C", __LINE__);
348
349}
350
351Float_t EtaToTheta(Float_t arg){
352 return (180./TMath::Pi())*2.*atan(exp(-arg));
353}