2 // Class to do the sharing correction for MC data.
5 // - AliESDFMD object - from reconstruction
10 // - AliESDFMD object - copy of input, but with signals merged
16 // - For each ring (FMD1i, FMD2i, FMD2o, FMD3i, FMD3o) the distribution of
17 // signals before and after the filter.
18 // - For each ring (see above), an array of distributions of number of
19 // hit strips for each vertex bin (if enabled - see Init method)
21 #include "AliFMDMCSharingFilter.h"
22 #include <AliESDFMD.h>
23 #include <AliMCEvent.h>
24 #include <AliTrackReference.h>
30 #include "AliFMDStripIndex.h"
36 ClassImp(AliFMDMCSharingFilter)
38 ; // This is for Emacs
42 //____________________________________________________________________
43 AliFMDMCSharingFilter::AliFMDMCSharingFilter(const char* title)
44 : AliFMDSharingFilter(title),
56 // title Title of object - not significant
58 fFMD1i = new TH2D("FMD1i_corr", "Merged vs MC", 21, -.5, 20.5, 300, 0, 15);
59 fFMD2i = new TH2D("FMD2i_corr", "Merged vs MC", 21, -.5, 20.5, 300, 0, 15);
60 fFMD2o = new TH2D("FMD2o_corr", "Merged vs MC", 21, -.5, 20.5, 300, 0, 15);
61 fFMD3i = new TH2D("FMD3i_corr", "Merged vs MC", 21, -.5, 20.5, 300, 0, 15);
62 fFMD3o = new TH2D("FMD3o_corr", "Merged vs MC", 21, -.5, 20.5, 300, 0, 15);
63 fFMD1i->SetYTitle("#Delta E/#Delta_{mip} (ESD)");
64 fFMD1i->SetXTitle("Hits (MC)");
65 fFMD2i->SetYTitle("#Delta E/#Delta_{mip} (ESD)");
66 fFMD2i->SetXTitle("Hits (MC)");
67 fFMD2o->SetYTitle("#Delta E/#Delta_{mip} (ESD)");
68 fFMD2o->SetXTitle("Hits (MC)");
69 fFMD3i->SetYTitle("#Delta E/#Delta_{mip} (ESD)");
70 fFMD3i->SetXTitle("Hits (MC)");
71 fFMD3o->SetYTitle("#Delta E/#Delta_{mip} (ESD)");
72 fFMD3o->SetXTitle("Hits (MC)");
73 fFMD1i->SetDirectory(0);
74 fFMD2i->SetDirectory(0);
75 fFMD2o->SetDirectory(0);
76 fFMD3i->SetDirectory(0);
77 fFMD3o->SetDirectory(0);
78 fSumEta = new TH1D("mcSumEta", "MC INEL Truth", 200, -4, 6);
79 fSumEta->SetXTitle("#eta");
80 fSumEta->SetYTitle("dN_{ch}/d#eta");
81 fSumEta->SetDirectory(0);
83 fSumEta->SetMarkerColor(kOrange+2);
84 fSumEta->SetMarkerStyle(22);
85 fSumEta->SetFillColor(0);
86 fSumEta->SetFillStyle(0);
90 //____________________________________________________________________
91 AliFMDMCSharingFilter::AliFMDMCSharingFilter(const AliFMDMCSharingFilter& o)
92 : AliFMDSharingFilter(o),
104 // o Object to copy from
108 //____________________________________________________________________
109 AliFMDMCSharingFilter::~AliFMDMCSharingFilter()
114 if (fFMD1i) delete fFMD1i;
115 if (fFMD2i) delete fFMD2i;
116 if (fFMD2o) delete fFMD2o;
117 if (fFMD3i) delete fFMD3i;
118 if (fFMD3o) delete fFMD3o;
119 if (fSumEta) delete fSumEta;
122 //____________________________________________________________________
123 AliFMDMCSharingFilter&
124 AliFMDMCSharingFilter::operator=(const AliFMDMCSharingFilter& o)
127 // Assignment operator
130 // o Object to assign from
135 AliFMDSharingFilter::operator=(o);
139 //____________________________________________________________________
141 AliFMDMCSharingFilter::StoreParticle(UShort_t d,
145 AliESDFMD& output) const
148 // Store a particle hit in FMD<i>dr</i>[<i>s,t</i>] in @a output
155 // output Output ESD object
157 Double_t old = output.Multiplicity(d,r,s,t);
158 if (old == AliESDFMD::kInvalidMult) old = 0;
159 output.SetMultiplicity(d,r,s,t,old+1);
162 //____________________________________________________________________
164 AliFMDMCSharingFilter::FilterMC(const AliESDFMD& input,
165 const AliMCEvent& event,
171 // Filter the input kinematics and track references, using
172 // some of the ESD information
175 // input Input ESD event
176 // event Input MC event
177 // vz Vertex position
178 // output Output ESD-like object
179 // primary Per-event histogram of primaries
182 // True on succes, false otherwise
186 // Copy eta values to output
187 for (UShort_t ed = 1; ed <= 3; ed++) {
188 UShort_t nq = (ed == 1 ? 1 : 2);
189 for (UShort_t eq = 0; eq < nq; eq++) {
190 Char_t er = (eq == 0 ? 'I' : 'O');
191 UShort_t ns = (eq == 0 ? 20 : 40);
192 UShort_t nt = (eq == 0 ? 512 : 256);
193 for (UShort_t es = 0; es < ns; es++)
194 for (UShort_t et = 0; et < nt; et++)
195 output.SetEta(ed, er, es, et, input.Eta(ed, er, es, et));
198 AliStack* stack = const_cast<AliMCEvent&>(event).Stack();
199 Int_t nTracks = event.GetNumberOfTracks();
200 for (Int_t iTr = 0; iTr < nTracks; iTr++) {
201 AliMCParticle* particle =
202 static_cast<AliMCParticle*>(event.GetTrack(iTr));
204 // Check the returned particle
205 if (!particle) continue;
207 // Check if this charged and a primary
208 Bool_t isCharged = particle->Charge() != 0;
209 if (!isCharged) continue;
210 Bool_t isPrimary = stack->IsPhysicalPrimary(iTr);
212 // Fill 'dn/deta' histogram
214 fSumEta->Fill(particle->Eta());
215 primary->Fill(particle->Eta(), particle->Phi());
218 Int_t nTrRef = particle->GetNumberOfTrackReferences();
221 UShort_t oD = 0, oS = 1024, oT = 1024;
223 for (Int_t iTrRef = 0; iTrRef < nTrRef; iTrRef++) {
224 AliTrackReference* ref = particle->GetTrackReference(iTrRef);
229 // Check that we hit an FMD element
230 if (ref->DetectorId() != AliTrackReference::kFMD)
233 // Get the detector coordinates
236 AliFMDStripIndex::Unpack(ref->UserId(), d, r, s, t);
237 if (oD > 0 && oR != '\0' && (d != oD || r != oR)) {
239 StoreParticle(oD, oR, oS, oT, output);
247 // The longest passage is determined through the angle
248 Double_t x = ref->X();
249 Double_t y = ref->Y();
250 Double_t z = ref->Z()-vz;
251 Double_t rr = TMath::Sqrt(x*x+y*y);
252 Double_t theta= TMath::ATan2(rr,z);
253 Double_t ang = TMath::Abs(TMath::Pi()-theta);
258 } // Loop over track references
259 if (longest < 0) continue;
261 // Get the reference corresponding to the longest path through the detector
262 AliTrackReference* ref = particle->GetTrackReference(longest);
264 // Get the detector coordinates
267 AliFMDStripIndex::Unpack(ref->UserId(), d, r, s, t);
269 StoreParticle(d,r,s,t,output);
270 } // Loop over tracks
274 //____________________________________________________________________
276 AliFMDMCSharingFilter::CompareResults(const AliESDFMD& esd,
280 // Compare the result of merging to the monte-carlo truth. This
281 // fills the correlation histograms
284 // esd ESD after sharing correction
288 // Copy eta values to output
289 for (UShort_t d = 1; d <= 3; d++) {
290 UShort_t nq = (d == 1 ? 1 : 2);
291 for (UShort_t q = 0; q < nq; q++) {
292 Char_t r = (q == 0 ? 'I' : 'O');
293 UShort_t ns = (q == 0 ? 20 : 40);
294 UShort_t nt = (q == 0 ? 512 : 256);
297 case 1: co = fFMD1i; break;
298 case 2: co = (q == 0 ? fFMD2i : fFMD2o); break;
299 case 3: co = (q == 0 ? fFMD3i : fFMD3o); break;
302 for (UShort_t s = 0; s < ns; s++) {
303 for (UShort_t t = 0; t < nt; t++) {
304 Float_t mEsd = esd.Multiplicity(d, r, s, t);
305 Float_t mMc = mc.Multiplicity(d, r, s, t);
314 //____________________________________________________________________
316 AliFMDMCSharingFilter::DefineOutput(TList* dir)
319 // Define the output histograms. These are put in a sub list of the
320 // passed list. The histograms are merged before the parent task calls
321 // AliAnalysisTaskSE::Terminate
324 // dir Directory to add to
326 AliFMDSharingFilter::DefineOutput(dir);
327 TList* d = static_cast<TList*>(dir->FindObject(GetName()));
328 TList* cd = new TList;
329 cd->SetName("esd_mc_comparion");
339 //____________________________________________________________________
341 AliFMDMCSharingFilter::ScaleHistograms(TList* dir, Int_t nEvents)
344 // Scale the histograms to the total number of events
347 // dir Where the output is
348 // nEvents Number of events
350 AliFMDSharingFilter::ScaleHistograms(dir, nEvents);
351 TH1D* sumEta = static_cast<TH1D*>(dir->FindObject("mcSumEta"));
353 AliWarning(Form("No mcSumEta histogram found in output list"));
356 sumEta->Scale(1. / nEvents, "width");
359 //____________________________________________________________________
361 AliFMDMCSharingFilter::Print(Option_t* option) const
369 char ind[gROOT->GetDirLevel()+1];
370 for (Int_t i = 0; i < gROOT->GetDirLevel(); i++) ind[i] = ' ';
371 ind[gROOT->GetDirLevel()] = '\0';
372 AliFMDSharingFilter::Print(option);
375 //____________________________________________________________________