]> git.uio.no Git - u/mrichter/AliRoot.git/blob - doc/aliroot-primer/scripts/Config.C
Fixes for bug #52499: Field polarities inconsistiency
[u/mrichter/AliRoot.git] / doc / aliroot-primer / scripts / Config.C
1 // Function converting pseudorapidity
2 // interval to polar angle interval. It is used to set 
3 // the limits in the generator
4 Float_t EtaToTheta(Float_t arg){
5   return (180./TMath::Pi())*2.*atan(exp(-arg));
6 }
7
8 // Set Random Number seed using the current time
9 TDatime dat;
10 static UInt_t sseed = dat.Get();
11
12 void Config()
13 {
14   gRandom->SetSeed(sseed);
15   cout<<"Seed for random number generation= "<<gRandom->GetSeed()<<endl; 
16     
17   // Load GEANT 3 library. It has to be in LD_LIBRARY_PATH
18   gSystem->Load("libgeant321");
19
20   // Instantiation of the particle transport package. gMC is set internaly
21   new TGeant3TGeo("C++ Interface to Geant3");
22
23   // Create run loader and set some properties
24   AliRunLoader* rl =  AliRunLoader::Open("galice.root",
25                                          AliConfig::GetDefaultEventFolderName(),
26                                          "recreate");
27   if (!rl) Fatal("Config.C","Can not instatiate the Run Loader");
28   rl->SetCompressionLevel(2);
29   rl->SetNumberOfEventsPerFile(3);
30
31   // Register the run loader in gAlice
32   gAlice->SetRunLoader(rl);
33
34   // Set external decayer
35   TVirtualMCDecayer *decayer = new AliDecayerPythia();
36   decayer->SetForceDecay(kAll); // kAll means no specific decay is forced
37   decayer->Init();
38
39   // Register the external decayer in the transport package
40   gMC->SetExternalDecayer(decayer);
41
42   // STEERING parameters FOR ALICE SIMULATION
43   // Specify event type to be transported through the ALICE setup
44   // All positions are in cm, angles in degrees, and P and E in GeV
45   // For the details see the GEANT 3 manual
46
47   // Switch on/off the physics processes (global)
48   // Please consult the file data/galice.cuts for detector
49   // specific settings, i.e. DRAY
50   gMC->SetProcess("DCAY",1); // Particle decay
51   gMC->SetProcess("PAIR",1); // Pair production
52   gMC->SetProcess("COMP",1); // Compton scattering
53   gMC->SetProcess("PHOT",1); // Photo effect
54   gMC->SetProcess("PFIS",0); // Photo fission
55   gMC->SetProcess("DRAY",0); // Delta rays
56   gMC->SetProcess("ANNI",1); // Positron annihilation
57   gMC->SetProcess("BREM",1); // Bremstrahlung
58   gMC->SetProcess("MUNU",1); // Muon nuclear interactions
59   gMC->SetProcess("CKOV",1); // Cerenkov production
60   gMC->SetProcess("HADR",1); // Hadronic interactions
61   gMC->SetProcess("LOSS",2); // Energy loss (2=complete fluct.)
62   gMC->SetProcess("MULS",1); // Multiple scattering
63   gMC->SetProcess("RAYL",1); // Rayleigh scattering
64   
65   // Set the transport package cuts
66   Float_t cut = 1.e-3;        // 1MeV cut by default
67   Float_t tofmax = 1.e10;
68
69   gMC->SetCut("CUTGAM", cut); // Cut for gammas
70   gMC->SetCut("CUTELE", cut); // Cut for electrons
71   gMC->SetCut("CUTNEU", cut); // Cut for neutral hadrons
72   gMC->SetCut("CUTHAD", cut); // Cut for charged hadrons
73   gMC->SetCut("CUTMUO", cut); // Cut for muons
74   gMC->SetCut("BCUTE",  cut); // Cut for electron brems.
75   gMC->SetCut("BCUTM",  cut); // Cut for muon brems.
76   gMC->SetCut("DCUTE",  cut); // Cut for electron delta-rays
77   gMC->SetCut("DCUTM",  cut); // Cut for muon delta-rays
78   gMC->SetCut("PPCUTM", cut); // Cut for e+e- pairs by muons
79   gMC->SetCut("TOFMAX", tofmax); // Time of flight cut
80   
81   // Set up the particle generation
82
83   // AliGenCocktail permits to combine several different generators
84   AliGenCocktail *gener = new AliGenCocktail();
85
86   // The phi range is always inside 0-360
87   gener->SetPhiRange(0, 360);
88
89   // Set pseudorapidity range from -8 to 8.
90   Float_t thmin = EtaToTheta(8);   // theta min. <---> eta max
91   Float_t thmax = EtaToTheta(-8);  // theta max. <---> eta min 
92   gener->SetThetaRange(thmin,thmax);
93
94   gener->SetOrigin(0, 0, 0);  // vertex position
95   gener->SetSigma(0, 0, 5.3);   // Sigma in (X,Y,Z) (cm) on IP position
96   gener->SetCutVertexZ(1.);     // Truncate at 1 sigma
97   gener->SetVertexSmear(kPerEvent); 
98
99   // First cocktail component: 100 ``background'' particles 
100   AliGenHIJINGpara *hijingparam = new AliGenHIJINGpara(100);
101   hijingparam->SetMomentumRange(0.2, 999);
102   gener->AddGenerator(hijingparam,"HIJING PARAM",1);
103
104   // Second cocktail component: one gamma in PHOS direction
105   AliGenBox *genbox = new AliGenBox(1);
106   genbox->SetMomentumRange(10,11.);
107   genbox->SetPhiRange(270.5,270.7);
108   genbox->SetThetaRange(90.5,90.7);
109   genbox->SetPart(kGamma);
110   gener->AddGenerator(genbox,"GENBOX GAMMA for PHOS",1);
111
112   gener->Init(); // Initialization of the coctail generator
113
114   // Field (the last parameter is 1 => L3 0.4 T)
115   TGeoGlobalMagField::Instance()->SetField(new AliMagF("Maps","Maps", 2, -1., -1., 10., AliMagF::k5kG));
116
117   // Make sure the current ROOT directory is in galice.root 
118   rl->CdGAFile();
119
120   // Build the setup and set some detector parameters
121
122   // ALICE BODY parameters. BODY is always present
123   AliBODY *BODY = new AliBODY("BODY", "ALICE envelop");
124
125   // Start with Magnet since detector layouts may be depending
126   // on the selected Magnet dimensions
127   AliMAG *MAG = new AliMAG("MAG", "Magnet");
128
129   AliABSO *ABSO = new AliABSOv0("ABSO", "Muon Absorber");       // Absorber
130
131   AliDIPO *DIPO = new AliDIPOv2("DIPO", "Dipole version 2");    // Dipole magnet
132
133   AliHALL *HALL = new AliHALL("HALL", "ALICE Hall");            // Hall
134     
135   AliFRAMEv2 *FRAME = new AliFRAMEv2("FRAME", "Space Frame");   // Space frame
136
137   AliSHIL *SHIL = new AliSHILv2("SHIL", "Shielding Version 2"); // Shielding
138
139   AliPIPE *PIPE = new AliPIPEv0("PIPE", "Beam Pipe");           // Beam pipe
140     
141   // ITS parameters
142   AliITSvPPRasymmFMD *ITS  = new AliITSvPPRasymmFMD("ITS",
143                "ITS PPR detailed version with asymmetric services");
144   ITS->SetMinorVersion(2);      // don't change it if you're not an ITS developer
145   ITS->SetReadDet(kFALSE);      // don't change it if you're not an ITS developer
146   ITS->SetThicknessDet1(200.);  // detector thickness on layer 1:[100,300] mkm
147   ITS->SetThicknessDet2(200.);  // detector thickness on layer 2:[100,300] mkm
148   ITS->SetThicknessChip1(150.); // chip thickness on layer 1: [150,300] mkm
149   ITS->SetThicknessChip2(150.); // chip thickness on layer 2: [150,300]
150   ITS->SetRails(0);             // 1 --> rails in ; 0 --> rails out
151   ITS->SetCoolingFluid(1);      // 1 --> water ; 0 --> freon
152   ITS->SetEUCLID(0);            // no output for the EUCLID CAD system 
153
154   
155   AliTPC *TPC = new AliTPCv2("TPC", "Default");                 // TPC
156
157   AliTOF *TOF = new AliTOFv5T0("TOF", "normal TOF");            // TOF
158
159   AliHMPID *HMPID = new AliHMPIDv1("HMPID", "normal HMPID");         // HMPID
160
161   AliZDC *ZDC = new AliZDCv2("ZDC", "normal ZDC");              // ZDC
162
163   AliTRD *TRD = new AliTRDv1("TRD", "TRD slow simulator");      // TRD
164
165   AliFMD *FMD = new AliFMDv1("FMD", "normal FMD");              // FMD
166
167   AliMUON *MUON = new AliMUONv1("MUON", "default");             // MUON
168
169   AliPHOS *PHOS = new AliPHOSv1("PHOS", "IHEP");                // PHOS
170
171   AliPMD *PMD = new AliPMDv1("PMD", "normal PMD");              // PMD
172
173   AliT0 *T0 = new AliT0v1("T0", "T0 Detector");  // T0
174
175   // EMCAL
176   AliEMCAL *EMCAL = new AliEMCALv2("EMCAL", "EMCAL_COMPLETE");
177
178   AliVZERO *VZERO = new AliVZEROv7("VZERO", "normal VZERO");    // VZERO
179 }