]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPHIC/ConfigTPHIC.C
Adding macros to create Calibration objects
[u/mrichter/AliRoot.git] / TPHIC / ConfigTPHIC.C
1 static Int_t    eventsPerRun = 100;
2 enum PprGeo_t 
3 {
4     kHoles, kNoHoles
5 };
6 static PprGeo_t geo = kHoles;
7
8 void Config()
9 {
10     // 7-DEC-2000 09:00
11     // Switch on Transition Radiation simulation. 6/12/00 18:00
12     // iZDC=1  7/12/00 09:00
13     // ThetaRange is (0., 180.). It was (0.28,179.72) 7/12/00 09:00
14     // Theta range given through pseudorapidity limits 22/6/2001
15
16     // Set Random Number seed
17     // gRandom->SetSeed(12345);
18
19
20    // libraries required by geant321
21     gSystem->Load("libgeant321");
22
23     new     TGeant3("C++ Interface to Geant3");
24
25     if (!gSystem->Getenv("CONFIG_FILE"))
26     {
27         TFile  *rootfile = new TFile("galice.root", "recreate");
28
29         rootfile->SetCompressionLevel(2);
30     }
31
32     TGeant3 *geant3 = (TGeant3 *) gMC;
33
34     //
35     // Set External decayer
36     TVirtualMCDecayer *decayer = new AliDecayerPythia();
37
38     decayer->SetForceDecay(kAll);
39     decayer->Init();
40     gMC->SetExternalDecayer(decayer);
41     //
42     //
43     //=======================================================================
44     // ******* GEANT STEERING parameters FOR ALICE SIMULATION *******
45     geant3->SetTRIG(1);         //Number of events to be processed 
46     geant3->SetSWIT(4, 1);
47     geant3->SetDEBU(0, 0, 1);
48     //geant3->SetSWIT(2,2);
49     geant3->SetDCAY(1);
50     geant3->SetPAIR(1);
51     geant3->SetCOMP(1);
52     geant3->SetPHOT(1);
53     geant3->SetPFIS(0);
54     geant3->SetDRAY(0);
55     geant3->SetANNI(1);
56     geant3->SetBREM(1);
57     geant3->SetMUNU(1);
58     geant3->SetCKOV(1);
59     geant3->SetHADR(1);         //Select pure GEANH (HADR 1) or GEANH/NUCRIN (HADR 3)
60     geant3->SetLOSS(2);
61     geant3->SetMULS(1);
62     geant3->SetRAYL(1);
63     geant3->SetAUTO(1);         //Select automatic STMIN etc... calc. (AUTO 1) or manual (AUTO 0)
64     geant3->SetABAN(0);         //Restore 3.16 behaviour for abandoned tracks
65     geant3->SetOPTI(2);         //Select optimisation level for GEANT geometry searches (0,1,2)
66     geant3->SetERAN(5.e-7);
67
68     Float_t cut = 1.e-3;        // 1MeV cut by default
69     Float_t tofmax = 1.e10;
70
71     //             GAM ELEC NHAD CHAD MUON EBREM MUHAB EDEL MUDEL MUPA TOFMAX
72     geant3->SetCUTS(cut, cut, cut, cut, cut, cut, cut, cut, cut, cut,
73                     tofmax);
74     //
75     //=======================================================================
76     // ************* STEERING parameters FOR ALICE SIMULATION **************
77     // --- Specify event type to be tracked through the ALICE setup
78     // --- All positions are in cm, angles in degrees, and P and E in GeV
79
80     Int_t process=1;
81     GenTPHIC(process);
82
83     // 
84     // Activate this line if you want the vertex smearing to happen
85     // track by track
86     //
87     //gener->SetVertexSmear(perTrack); 
88     // Field (L3 0.4 T)
89     TGeoGlobalMagField::Instance()->SetField(new AliMagF("Maps","Maps", 1., 1., AliMagF::k5kG));
90     rootfile->cd();
91
92
93     Int_t   iABSO  =  0;
94     Int_t   iDIPO  =  0;
95     Int_t   iFMD   =  0;
96     Int_t   iFRAME =  0;
97     Int_t   iHALL  =  0;
98     Int_t   iITS   =  0;
99     Int_t   iMAG   =  0;
100     Int_t   iMUON  =  0;
101     Int_t   iPHOS  =  0;
102     Int_t   iPIPE  =  0;
103     Int_t   iPMD   =  0;
104     Int_t   iHMPID  =  0;
105     Int_t   iSHIL  =  0;
106     Int_t   iT0 =  0;
107     Int_t   iTOF   =  0;
108     Int_t   iTPC   =  0;
109     Int_t   iTRD   =  0;
110     Int_t   iZDC   =  0;
111     Int_t   iEMCAL =  0;
112     Int_t   iACORDE   =  0;
113     Int_t   iVZERO =  0;
114
115     //=================== Alice BODY parameters =============================
116     AliBODY *BODY = new AliBODY("BODY", "Alice envelop");
117
118     if (iMAG)
119     {
120         //=================== MAG parameters ============================
121         // --- Start with Magnet since detector layouts may be depending ---
122         // --- on the selected Magnet dimensions ---
123         AliMAG *MAG = new AliMAG("MAG", "Magnet");
124     }
125
126
127     if (iABSO)
128     {
129         //=================== ABSO parameters ============================
130         AliABSO *ABSO = new AliABSOv0("ABSO", "Muon Absorber");
131     }
132
133     if (iDIPO)
134     {
135         //=================== DIPO parameters ============================
136
137         AliDIPO *DIPO = new AliDIPOv2("DIPO", "Dipole version 2");
138     }
139
140     if (iHALL)
141     {
142         //=================== HALL parameters ============================
143
144         AliHALL *HALL = new AliHALL("HALL", "Alice Hall");
145     }
146
147
148     if (iFRAME)
149     {
150         //=================== FRAME parameters ============================
151
152         AliFRAMEv2 *FRAME = new AliFRAMEv2("FRAME", "Space Frame");
153         if (geo == kHoles) {
154             FRAME->SetHoles(1);
155         } else {
156             FRAME->SetHoles(0);
157         }
158     }
159
160     if (iSHIL)
161     {
162         //=================== SHIL parameters ============================
163
164         AliSHIL *SHIL = new AliSHILv2("SHIL", "Shielding Version 2");
165     }
166
167
168     if (iPIPE)
169     {
170         //=================== PIPE parameters ============================
171
172         AliPIPE *PIPE = new AliPIPEv0("PIPE", "Beam Pipe");
173     }
174  
175     if(iITS) {
176
177     //=================== ITS parameters ============================
178     //
179     // As the innermost detector in ALICE, the Inner Tracking System "impacts" on
180     // almost all other detectors. This involves the fact that the ITS geometry
181     // still has several options to be followed in parallel in order to determine
182     // the best set-up which minimizes the induced background. All the geometries
183     // available to date are described in the following. Read carefully the comments
184     // and use the default version (the only one uncommented) unless you are making
185     // comparisons and you know what you are doing. In this case just uncomment the
186     // ITS geometry you want to use and run Aliroot.
187     //
188     // Detailed geometries:         
189     //
190     //
191     //AliITS *ITS  = new AliITSv5symm("ITS","Updated ITS TDR detailed version with symmetric services");
192     //
193     //AliITS *ITS  = new AliITSv5asymm("ITS","Updates ITS TDR detailed version with asymmetric services");
194     //
195         AliITSvPPRasymm *ITS  = new AliITSvPPRasymm("ITS","New ITS PPR detailed version with asymmetric services");
196         ITS->SetMinorVersion(2);                                         // don't touch this parameter if you're not an ITS developer
197         ITS->SetReadDet(kFALSE);                                         // don't touch this parameter if you're not an ITS developer
198     //    ITS->SetWriteDet("$ALICE_ROOT/ITS/ITSgeometry_vPPRasymm2.det");  // don't touch this parameter if you're not an ITS developer
199         ITS->SetThicknessDet1(200.);   // detector thickness on layer 1 must be in the range [100,300]
200         ITS->SetThicknessDet2(200.);   // detector thickness on layer 2 must be in the range [100,300]
201         ITS->SetThicknessChip1(200.);  // chip thickness on layer 1 must be in the range [150,300]
202         ITS->SetThicknessChip2(200.);  // chip thickness on layer 2 must be in the range [150,300]
203         ITS->SetRails(0);            // 1 --> rails in ; 0 --> rails out
204         ITS->SetCoolingFluid(1);   // 1 --> water ; 0 --> freon
205         //
206     //AliITSvPPRsymm *ITS  = new AliITSvPPRsymm("ITS","New ITS PPR detailed version with symmetric services");
207     //ITS->SetMinorVersion(2);                                       // don't touch this parameter if you're not an ITS developer
208     //ITS->SetReadDet(kFALSE);                                       // don't touch this parameter if you're not an ITS developer
209     //ITS->SetWriteDet("$ALICE_ROOT/ITS/ITSgeometry_vPPRsymm2.det"); // don't touch this parameter if you're not an ITS developer
210     //ITS->SetThicknessDet1(200.);   // detector thickness on layer 1 must be in the range [100,300]
211     //ITS->SetThicknessDet2(200.);   // detector thickness on layer 2 must be in the range [100,300]
212     //ITS->SetThicknessChip1(200.);  // chip thickness on layer 1 must be in the range [150,300]
213     //ITS->SetThicknessChip2(200.);  // chip thickness on layer 2 must be in the range [150,300]
214     //ITS->SetRails(0);              // 1 --> rails in ; 0 --> rails out
215     //ITS->SetCoolingFluid(1);       // 1 --> water ; 0 --> freon
216     //
217     //
218     // Coarse geometries (warning: no hits are produced with these coarse geometries and they unuseful 
219     // for reconstruction !):
220     //                                                     
221     //
222     //AliITSvPPRcoarseasymm *ITS  = new AliITSvPPRcoarseasymm("ITS","New ITS PPR coarse version with asymmetric services");
223     //ITS->SetRails(0);                // 1 --> rails in ; 0 --> rails out
224     //ITS->SetSupportMaterial(0);      // 0 --> Copper ; 1 --> Aluminum ; 2 --> Carbon
225     //
226     //AliITS *ITS  = new AliITSvPPRcoarsesymm("ITS","New ITS PPR coarse version with symmetric services");
227     //ITS->SetRails(0);                // 1 --> rails in ; 0 --> rails out
228     //ITS->SetSupportMaterial(0);      // 0 --> Copper ; 1 --> Aluminum ; 2 --> Carbon
229     //                      
230     //
231     //
232     // Geant3 <-> EUCLID conversion
233     // ============================
234     //
235     // SetEUCLID is a flag to output (=1) or not to output (=0) both geometry and
236     // media to two ASCII files (called by default ITSgeometry.euc and
237     // ITSgeometry.tme) in a format understandable to the CAD system EUCLID.
238     // The default (=0) means that you dont want to use this facility.
239     //
240         ITS->SetEUCLID(0);  
241     }
242
243     if (iTPC)
244     {
245         //============================ TPC parameters ================================
246         // --- This allows the user to specify sectors for the SLOW (TPC geometry 2)
247         // --- Simulator. SecAL (SecAU) <0 means that ALL lower (upper)
248         // --- sectors are specified, any value other than that requires at least one 
249         // --- sector (lower or upper)to be specified!
250         // --- Reminder: sectors 1-24 are lower sectors (1-12 -> z>0, 13-24 -> z<0)
251         // ---           sectors 25-72 are the upper ones (25-48 -> z>0, 49-72 -> z<0)
252         // --- SecLows - number of lower sectors specified (up to 6)
253         // --- SecUps - number of upper sectors specified (up to 12)
254         // --- Sens - sensitive strips for the Slow Simulator !!!
255         // --- This does NOT work if all S or L-sectors are specified, i.e.
256         // --- if SecAL or SecAU < 0
257         //
258         //
259         //-----------------------------------------------------------------------------
260
261         //  gROOT->LoadMacro("SetTPCParam.C");
262         //  AliTPCParam *param = SetTPCParam();
263         AliTPC *TPC = new AliTPCv2("TPC", "Default");
264
265         // All sectors included 
266         TPC->SetSecAL(-1);
267         TPC->SetSecAU(-1);
268
269     }
270
271
272     if (iTOF) {
273         if (geo == kHoles) {
274         //=================== TOF parameters ============================
275             AliTOF *TOF = new AliTOFv2FHoles("TOF", "TOF with Holes");
276         } else {
277             AliTOF *TOF = new AliTOFv4T0("TOF", "normal TOF");
278         }
279     }
280
281
282     if (iHMPID)
283     {
284         //=================== HMPID parameters ===========================
285         AliHMPID *HMPID = new AliHMPIDv3("HMPID", "normal HMPID");
286
287     }
288
289
290     if (iZDC)
291     {
292         //=================== ZDC parameters ============================
293
294         AliZDC *ZDC = new AliZDCv2("ZDC", "normal ZDC");
295     }
296
297     if (iTRD)
298     {
299         //=================== TRD parameters ============================
300
301         AliTRD *TRD = new AliTRDv1("TRD", "TRD slow simulator");
302
303         // Select the gas mixture (0: 97% Xe + 3% isobutane, 1: 90% Xe + 10% CO2)
304         TRD->SetGasMix(1);
305         if (geo == kHoles) {
306             // With hole in front of PHOS
307             TRD->SetPHOShole();
308             // With hole in front of HMPID
309             TRD->SetHMPIDhole();
310         }
311             // Switch on TR
312             AliTRDsim *TRDsim = TRD->CreateTR();
313     }
314
315     if (iFMD)
316     {
317         //=================== FMD parameters ============================
318         AliFMD *FMD = new AliFMDv1("FMD", "normal FMD");
319         FMD->SetRingsSi1(256);
320         FMD->SetRingsSi2(128);
321         FMD->SetSectorsSi1(20);
322         FMD->SetSectorsSi2(40);      
323    }
324
325     if (iMUON)
326     {
327         //=================== MUON parameters ===========================
328
329         AliMUON *MUON = new AliMUONv1("MUON", "default");
330     }
331     //=================== PHOS parameters ===========================
332
333     if (iPHOS)
334     {
335         AliPHOS *PHOS = new AliPHOSv1("PHOS", "IHEP");
336     }
337
338
339     if (iPMD)
340     {
341         //=================== PMD parameters ============================
342         AliPMD *PMD = new AliPMDv1("PMD", "normal PMD");
343     }
344
345     if (iT0)
346     {
347         //=================== T0 parameters ============================
348         AliT0 *T0 = new AliT0v1("T0", "T0 Detector");
349     }
350
351     if (iEMCAL)
352     {
353         //=================== EMCAL parameters ============================
354         AliEMCAL *EMCAL = new AliEMCALv2("EMCAL", "EMCAL_COMPLETE");
355     }
356
357      if (iACORDE)
358     {
359         //=================== ACORDE parameters ============================
360         AliACORDE *ACORDE = new AliACORDEv1("ACORDE", "normal ACORDE");
361     }
362
363      if (iVZERO)
364     {
365         //=================== ACORDE parameters ============================
366         AliVZERO *VZERO = new AliVZEROv2("VZERO", "normal VZERO");
367     }
368
369 }
370
371 Float_t EtaToTheta(Float_t arg){
372   return (180./TMath::Pi())*2.*atan(exp(-arg));
373 }
374
375 //-----------------------------------------------------------------------------
376 GenTPHIC(Int_t process)
377 {
378   AliGenTPHIC *gener = new AliGenTPHIC();
379
380   gener->SetMomentumRange(0, 999);
381   gener->SetPhiRange(0, 360);
382   Float_t thmin = EtaToTheta(8);   // theta min. <---> eta max
383   Float_t thmax = EtaToTheta(-8);  // theta max. <---> eta min 
384   gener->SetThetaRange(thmin,thmax);
385   gener->SetOrigin(0, 0, 0);  //vertex position
386   gener->SetSigma(0, 0, 0);   //Sigma in (X,Y,Z) (cm) on IP position
387
388   // Set TPHIC-specific parameters
389   gener->SetProcess(process);
390   if      (process == 1) { // gamma gamma -> X
391   }
392   else if (process == 2) { // gamma gamma -> eta_c
393     gener->SetKfOnium(441);
394     gener->SetGGwidthOnium(7.4e-06);
395   }
396   else if (process == 3) { // gamma gamma -> mu+ mu-
397     gener->SetKfFermion = 13;
398   }
399   else if (process == 4) { // gamma gamma -> W+ W-
400     gener->SetMggRange(70.,200.);
401     gener->SetYggRange(-5.,5.);
402     gener->SetLumFunName("lum_ca_70_200.dat");
403     gener->SetLumFunFlag(-1);
404   }
405   else if (process == 5) { // gamma gamma -> chi_1+ chi_1- (not implemented)
406   }
407   else if (process == 5) { // gamma gamma -> rho0 rho0
408     gener->SetKfVmesons(113,113);
409   }
410
411   gener->SetDebug(1);
412   gener->SetEventListRange(1,10);
413   gener->Init();
414 }