2 #if !defined(__CINT__) || defined(__MAKECINT__)
8 #include "TClonesArray.h"
10 #include "AliModule.h"
11 #include "AliFRAMEv1.h"
14 #include "AliTRDhit.h"
15 #include "AliTRDsim.h"
16 #include "AliTRDsimple.h"
17 #include "AliTRDsimpleGen.h"
18 #include "AliTRDsimpleMC.h"
19 #include "AliTRDgeometry.h"
20 #include "AliTRDdigitizer.h"
21 #include "AliTRDdigitsManager.h"
22 #include "AliTRDdataArrayI.h"
26 void AliTRDrunSimple()
29 ///////////////////////////////////////////////////////////////////////////////
31 // Macro to run the simplified simulator
33 ///////////////////////////////////////////////////////////////////////////////
35 //_____________________________________________________________________________
37 // Initialization part
38 //_____________________________________________________________________________
41 // The simple simulator
42 AliTRDsimple *simple = new AliTRDsimple();
45 // Initialize a dummy frame so that the TRD works
46 new AliFRAMEv1("FRAME","Space Frame");
48 // Initialize the TRD detector
49 AliTRDv1 *trd = new AliTRDv1("TRD","TRD slow simulator");
50 trd->SetHitTypeStandard();
52 // Needed for some material properties
53 trd->CreateMaterials();
55 // Select the gas mixture (0: 97% Xe + 3% isobutane, 1: 85% Xe + 15% CO2)
58 // Define the parameter object
59 // If no external parameter object is defined,
60 // default parameter will be used
61 AliTRDparameter *parameter = new AliTRDparameter("TRDparameter"
62 ,"TRD parameter class");
64 parameter->SetADCthreshold(0);
65 parameter->SetNTimeBin(30); // The number of timebins
66 parameter->SetExpandTimeBin(5,15); // The additional timebins
68 // Switch on TR production
69 AliTRDsim *trdsim = trd->CreateTR();
70 trdsim->SetNFoils(100);
71 trdsim->SetFoilThick(0.0013);
72 trdsim->SetGapThick(0.0500);
74 // Initialize the TRD object
78 AliTRDdigitizer *digitizer = new AliTRDdigitizer("TRDdigitizer","Digitizer class");
80 digitizer->SetParameter(parameter);
81 digitizer->SetSimple();
82 digitizer->InitDetector();
84 // The event generator
85 AliTRDsimpleGen *generator = simple->GetGenerator();
86 generator->SetMomentum(3.0,3.0);
87 generator->SetPdg(11); // Electrons
89 //_____________________________________________________________________________
92 //_____________________________________________________________________________
95 Int_t timeMax = parameter->GetTimeTotal();
96 Int_t adcRange = ((Int_t) parameter->GetADCoutRange());
98 TH1F *hQ = new TH1F("hQ" ,"Charge per hit (all)" ,100,0.0,1000.0);
99 TH1F *hQdedx = new TH1F("hQdedx","Charge per hit (dedx)",100,0.0,1000.0);
100 TH1F *hQtr = new TH1F("hQtr ","Charge per hit (tr)" ,100,0.0,1000.0);
102 TProfile *hQX = new TProfile("hQX" ,"Charge per hit vs X (all)" ,35,0.0,3.5);
103 TProfile *hQXdedx = new TProfile("hQXdedx","Charge per hit vs X (dedx)",35,0.0,3.5);
104 TProfile *hQXtr = new TProfile("hQXtr ","Charge per hit vs X (tr)" ,35,0.0,3.5);
106 TH1F *hNstep = new TH1F("hNstep","Number of steps / cm" , 151,-0.5, 150.5);
107 TH1F *hNel = new TH1F("hNel" ,"Number of electrons / cm",1001,-0.5,1000.5);
109 TH1F *hAmp = new TH1F("hAmp","Amplitude of the digits"
110 ,adcRange+1,-0.5,((Float_t) adcRange)+0.5);
111 TProfile *hAmpTime = new TProfile("hAmpTime","Amplitude vs timebin"
112 ,timeMax,-0.5,((Float_t) timeMax)-0.5);
114 //_____________________________________________________________________________
117 //_____________________________________________________________________________
120 // Number of events (tracks)
121 Int_t nEvent = 10000;
123 Float_t x0 = parameter->GetTime0(0) - AliTRDgeometry::DrThick();
125 TClonesArray *hitsArray = trd->Hits();
127 for (Int_t iEvent = 0; iEvent < nEvent; iEvent++) {
129 if (!(iEvent % 100) && (iEvent)) {
130 printf("Event no. %d\n",iEvent);
133 // Generate the hits for one track
134 simple->ProcessEvent(iEvent);
137 Float_t nElectrons = 0.0;
138 for (Int_t iHit = 0; iHit < hitsArray->GetEntries(); iHit++) {
140 AliTRDhit *hit = (AliTRDhit *) hitsArray->UncheckedAt(iHit);
141 Int_t charge = TMath::Abs(hit->GetCharge());
142 Float_t x = hit->X() - x0;
145 if (hit->FromDrift() ||
146 hit->FromAmplification()) {
147 hQdedx->Fill(charge);
148 hQXdedx->Fill(x,charge);
149 nElectrons += charge;
151 if (hit->FromTRphoton()) {
153 hQXtr->Fill(x,charge);
158 hNstep->Fill(((AliTRDsimpleMC *) gMC)->GetNStep() / 3.5);
159 hNel->Fill(nElectrons / 3.5);
162 digitizer->MakeDigits();
164 // Analyze the digits
167 AliTRDdigitsManager *digitsManager = digitizer->Digits();
168 AliTRDdataArrayI *digitsArray = digitsManager->GetDigits(12);
169 for (Int_t time = 0; time < timeMax; time++) {
171 Int_t amp = digitsArray->GetDataUnchecked(row,col,time);
173 hAmp->Fill((Double_t) amp);
174 hAmpTime->Fill((Double_t) time, (Double_t) amp);
181 //_____________________________________________________________________________
184 //_____________________________________________________________________________
187 TCanvas *cHit = new TCanvas("cHit","Hits",50,50,800,600);
193 hQ->SetXTitle("Q (a.u.)");
194 hQ->SetYTitle("Entries");
196 hQdedx->SetLineColor(3);
197 hQdedx->Draw("SAME");
198 hQtr->SetLineColor(2);
203 hQX->SetLineColor(4);
204 hQX->SetXTitle("x (cm)");
205 hQX->SetYTitle("<Q> (a.u.)");
207 hQXdedx->SetLineColor(3);
208 hQXdedx->Draw("SAMEHIST");
212 hQXtr->SetLineColor(2);
213 hQXtr->SetXTitle("x (cm)");
214 hQXtr->SetYTitle("<Q> (a.u.)");
218 hNstep->SetLineColor(4);
219 hNstep->SetXTitle("N_{step} / cm");
220 hNstep->SetYTitle("Entries");
224 hNel->SetLineColor(4);
225 hNel->SetXTitle("N_{electron} / cm");
226 hNel->SetYTitle("Entries");
229 TCanvas *cDigit = new TCanvas("cDigit","Digits",50,50,600,400);
234 hAmp->SetLineColor(2);
235 hAmp->SetXTitle("ADC channels");
236 hAmp->SetYTitle("Entries");
240 hAmpTime->SetLineColor(2);
241 hAmpTime->SetXTitle("Timebin");
242 hAmpTime->SetYTitle("<ADC channels>");
243 hAmpTime->Draw("HIST");
245 //_____________________________________________________________________________
247 // Save the histograms
248 //_____________________________________________________________________________
251 TFile *fileOut = new TFile("simple.root","RECREATE");