]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGLF/FORWARD/analysis2/AliFMDSharingFilter.cxx
Improvements
[u/mrichter/AliRoot.git] / PWGLF / FORWARD / analysis2 / AliFMDSharingFilter.cxx
1 //
2 // Class to do the sharing correction.  That is, a filter that merges 
3 // adjacent strip signals presumably originating from a single particle 
4 // that impinges on the detector in such a way that it deposite energy 
5 // into two or more strips. 
6 //
7 // Input: 
8 //    - AliESDFMD object  - from reconstruction
9 //
10 // Output: 
11 //    - AliESDFMD object  - copy of input, but with signals merged 
12 //
13 // Corrections used: 
14 //    - AliFMDCorrELossFit
15 //
16 // Histograms: 
17 //    - For each ring (FMD1i, FMD2i, FMD2o, FMD3i, FMD3o) the distribution of 
18 //      signals before and after the filter.  
19 //    - For each ring (see above), an array of distributions of number of 
20 //      hit strips for each vertex bin (if enabled - see Init method)
21 // 
22 //
23 //
24 #include "AliFMDSharingFilter.h"
25 #include "AliFMDStripIndex.h"
26 #include <AliESDFMD.h>
27 #include <TAxis.h>
28 #include <TList.h>
29 #include <TH1.h>
30 #include <TMath.h>
31 #include "AliForwardCorrectionManager.h"
32 #include "AliFMDCorrELossFit.h"
33 #include <AliLog.h>
34 #include <TROOT.h>
35 #include <THStack.h>
36 #include <TParameter.h>
37 #include <iostream>
38 #include <iomanip>
39
40 ClassImp(AliFMDSharingFilter)
41 #if 0
42 ; // This is for Emacs
43 #endif 
44
45 #define DBG(L,M)                                                        \
46   do { if (L>fDebug)break; std::cout << (M) << std::flush;} while(false) 
47 #define DBGL(L,M)                                                       \
48   do { if (L>fDebug)break; std::cout << (M) << std::endl;} while(false) 
49                               
50
51
52 //____________________________________________________________________
53 AliFMDSharingFilter::AliFMDSharingFilter()
54   : TNamed(), 
55     fRingHistos(),
56     fCorrectAngles(kFALSE), 
57     // fSummed(0),
58     fHighCuts(0),
59     fLowCuts(0),
60     // fOper(0),
61     fDebug(0),
62     fZeroSharedHitsBelowThreshold(false),
63     fLCuts(),
64     fHCuts(),
65     fUseSimpleMerging(false),
66     fThreeStripSharing(true),
67     fRecalculateEta(false),
68     // fExtraDead(0),
69     fXtraDead(0),
70     fInvalidIsEmpty(false),
71     fMergingDisabled(false),
72     fIgnoreESDForAngleCorrection(false)
73 {
74   // 
75   // Default Constructor - do not use 
76   //
77   DGUARD(fDebug,1, "Default CTOR for AliFMDSharingFilter");
78 }
79
80 //____________________________________________________________________
81 AliFMDSharingFilter::AliFMDSharingFilter(const char* title)
82   : TNamed("fmdSharingFilter", title), 
83     fRingHistos(), 
84     fCorrectAngles(kFALSE), 
85     // fSummed(0),
86     fHighCuts(0),
87     fLowCuts(0),
88     // fOper(0),
89     fDebug(0),
90     fZeroSharedHitsBelowThreshold(false),
91     fLCuts(),
92     fHCuts(),
93     fUseSimpleMerging(false),
94     fThreeStripSharing(true),
95     fRecalculateEta(false), 
96     // fExtraDead(51200),
97     fXtraDead(AliFMDStripIndex::Pack(3,'O',19,511)+1),
98     fInvalidIsEmpty(false),
99     fMergingDisabled(false),
100     fIgnoreESDForAngleCorrection(false)
101 {
102   // 
103   // Constructor 
104   // 
105   // Parameters:
106   //    title Title of object  - not significant 
107   //
108   DGUARD(fDebug,1, "Named CTOR for AliFMDSharingFilter: %s", title);
109   fRingHistos.SetName(GetName());
110   fRingHistos.SetOwner();
111   
112   fRingHistos.Add(new RingHistos(1, 'I'));
113   fRingHistos.Add(new RingHistos(2, 'I'));
114   fRingHistos.Add(new RingHistos(2, 'O'));
115   fRingHistos.Add(new RingHistos(3, 'I'));
116   fRingHistos.Add(new RingHistos(3, 'O'));
117
118   fHCuts.Set(AliFMDMultCuts::kLandauSigmaWidth, 1);
119   fLCuts.Set(AliFMDMultCuts::kFixed, .15);
120
121   // fExtraDead.Reset(-1);
122 }
123
124 //____________________________________________________________________
125 AliFMDSharingFilter::~AliFMDSharingFilter()
126 {
127   // 
128   // Destructor
129   //
130   DGUARD(fDebug,3, "DTOR for AliFMDSharingFilter");
131   // fRingHistos.Delete();
132 }
133
134 //____________________________________________________________________
135 AliFMDSharingFilter::RingHistos*
136 AliFMDSharingFilter::GetRingHistos(UShort_t d, Char_t r) const
137 {
138   // 
139   // Get the ring histogram container 
140   // 
141   // Parameters:
142   //    d Detector
143   //    r Ring 
144   // 
145   // Return:
146   //    Ring histogram container 
147   //
148   Int_t idx = -1;
149   switch (d) { 
150   case 1: idx = 0; break;
151   case 2: idx = 1 + (r == 'I' || r == 'i' ? 0 : 1); break;
152   case 3: idx = 3 + (r == 'I' || r == 'i' ? 0 : 1); break;
153   }
154   if (idx < 0 || idx >= fRingHistos.GetEntries()) return 0;
155   
156   return static_cast<RingHistos*>(fRingHistos.At(idx));
157 }
158
159 //____________________________________________________________________
160 void
161 AliFMDSharingFilter::AddDead(UShort_t d, Char_t r, UShort_t s, UShort_t t)
162 {
163   if (d < 1 || d > 3) {
164     Warning("AddDead", "Invalid detector FMD%d", d);
165     return;
166   }
167   Bool_t inner = (r == 'I' || r == 'i');
168   if (d == 1 && !inner) { 
169     Warning("AddDead", "Invalid ring FMD%d%c", d, r);
170     return;
171   }
172   if ((inner && s >= 20) || (!inner && s >= 40)) { 
173     Warning("AddDead", "Invalid sector FMD%d%c[%02d]", d, r, s);
174     return;
175   }
176   if ((inner && t >= 512) || (!inner && t >= 256)) { 
177     Warning("AddDead", "Invalid strip FMD%d%c[%02d,%03d]", d, r, s, t);
178     return;
179   }
180     
181   Int_t id = AliFMDStripIndex::Pack(d, r, s, t);
182   // Int_t i  = 0;
183   fXtraDead.SetBitNumber(id, true);
184 #if 0
185   for (i = 0; i < fExtraDead.GetSize(); i++) {
186     Int_t j = fExtraDead.At(i);
187     if (j == id) return; // Already there 
188     if (j <  0) break; // Free slot 
189   }
190   if (i >= fExtraDead.GetSize()) { 
191     Warning("AddDead", "No free slot to add FMD%d%c[%02d,%03d] at", 
192             d, r, s, t);
193     return;
194   }
195   fExtraDead[i] = id;
196 #endif
197 }
198 //____________________________________________________________________
199 void
200 AliFMDSharingFilter::AddDeadRegion(UShort_t d,  Char_t r, 
201                                    UShort_t s1, UShort_t s2, 
202                                    UShort_t t1, UShort_t t2)
203 {
204   // Add a dead region spanning from FMD<d><r>[<s1>,<t1>] to 
205   // FMD<d><r>[<s2>,<t2>] (both inclusive)
206   for (Int_t s = s1; s <= s2; s++) 
207     for (Int_t t = t1; t <= t2; t++) 
208       AddDead(d, r, s, t);
209 }
210 //____________________________________________________________________
211 void
212 AliFMDSharingFilter::AddDead(const Char_t* script)
213 {
214   if (!script || script[0] == '\0') return;
215   
216   gROOT->Macro(Form("%s((AliFMDSharingFilter*)%p);", script, this));
217 }
218
219 //____________________________________________________________________
220 Bool_t
221 AliFMDSharingFilter::IsDead(UShort_t d, Char_t r, UShort_t s, UShort_t t) const
222 {
223   Int_t id = AliFMDStripIndex::Pack(d, r, s, t);
224   return fXtraDead.TestBitNumber(id); 
225 #if 0
226   for (Int_t i = 0; i < fExtraDead.GetSize(); i++) {
227     Int_t j = fExtraDead.At(i);
228     if (j == id) {
229       //Info("IsDead", "FMD%d%c[%02d,%03d] marked as dead here", d, r, s, t);
230       return true;
231     }
232     if (j < 0) break; // High water mark 
233   }
234   return false;
235 #endif
236 }
237 //____________________________________________________________________
238 void
239 AliFMDSharingFilter::SetupForData(const TAxis& axis)
240 {
241   // Initialise - called on first event
242   DGUARD(fDebug,1, "Initialize for AliFMDSharingFilter");
243   AliForwardCorrectionManager& fcm  = AliForwardCorrectionManager::Instance();
244   const AliFMDCorrELossFit*    fits = fcm.GetELossFit();
245
246   // Compactify the xtra dead bits 
247   fXtraDead.Compact();
248
249   // Get the high cut.  The high cut is defined as the 
250   // most-probably-value peak found from the energy distributions, minus 
251   // 2 times the width of the corresponding Landau.
252   
253   TAxis eAxis(axis.GetNbins(),
254               axis.GetXmin(),
255               axis.GetXmax());
256   if(fits) 
257     eAxis.Set(fits->GetEtaAxis().GetNbins(), 
258               fits->GetEtaAxis().GetXmin(),
259               fits->GetEtaAxis().GetXmax());
260
261   UShort_t nEta = eAxis.GetNbins();
262         
263   fHighCuts->SetBins(nEta, eAxis.GetXmin(), eAxis.GetXmax(), 5, .5, 5.5);
264   fHighCuts->GetYaxis()->SetBinLabel(1, "FMD1i");
265   fHighCuts->GetYaxis()->SetBinLabel(2, "FMD2i");
266   fHighCuts->GetYaxis()->SetBinLabel(3, "FMD2o");
267   fHighCuts->GetYaxis()->SetBinLabel(4, "FMD3i");
268   fHighCuts->GetYaxis()->SetBinLabel(5, "FMD3o");
269
270   fLowCuts->SetBins(nEta, eAxis.GetXmin(), eAxis.GetXmax(), 5, .5, 5.5);
271   fLowCuts->GetYaxis()->SetBinLabel(1, "FMD1i");
272   fLowCuts->GetYaxis()->SetBinLabel(2, "FMD2i");
273   fLowCuts->GetYaxis()->SetBinLabel(3, "FMD2o");
274   fLowCuts->GetYaxis()->SetBinLabel(4, "FMD3i");
275   fLowCuts->GetYaxis()->SetBinLabel(5, "FMD3o");
276
277   // Cache our cuts in histograms 
278   fLCuts.FillHistogram(fLowCuts);
279   fHCuts.FillHistogram(fHighCuts);
280 }
281
282 //____________________________________________________________________
283 #define ETA2COS(ETA)                                            \
284   TMath::Cos(2*TMath::ATan(TMath::Exp(-TMath::Abs(ETA))))
285
286 Bool_t
287 AliFMDSharingFilter::Filter(const AliESDFMD& input, 
288                             Bool_t           /*lowFlux*/,
289                             AliESDFMD&       output, 
290                             Double_t         zvtx)
291 {
292   // 
293   // Filter the input AliESDFMD object
294   // 
295   // Parameters:
296   //    input     Input 
297   //    lowFlux   If this is a low-flux event 
298   //    output    Output AliESDFMD object 
299   // 
300   // Return:
301   //    True on success, false otherwise 
302   //
303   DGUARD(fDebug,1, "Filter event in AliFMDSharingFilter");
304   output.Clear();
305   TIter    next(&fRingHistos);
306   RingHistos* o      = 0;
307   while ((o = static_cast<RingHistos*>(next()))) o->Clear();
308
309   Int_t nSingle    = 0;
310   Int_t nDouble    = 0;
311   Int_t nTriple    = 0;
312
313   for(UShort_t d = 1; d <= 3; d++) {
314     Int_t nRings = (d == 1 ? 1 : 2);
315     for (UShort_t q = 0; q < nRings; q++) {
316       Char_t      r      = (q == 0 ? 'I' : 'O');
317       UShort_t    nsec   = (q == 0 ?  20 :  40);
318       UShort_t    nstr   = (q == 0 ? 512 : 256);
319       RingHistos* histos = GetRingHistos(d, r);
320       
321       for(UShort_t s = 0; s < nsec;  s++) {     
322         // `used' flags if the _current_ strip was used by _previous_ 
323         // iteration. 
324         Bool_t   used            = kFALSE;
325         // `eTotal' contains the current sum of merged signals so far 
326         Double_t eTotal          = -1;
327         // Int_t    nDistanceBefore = -1;
328         // Int_t    nDistanceAfter  = -1;
329         // `twoLow' flags if we saw two consequtive strips with a 
330         // signal between the two cuts. 
331         Bool_t   twoLow          = kFALSE;
332         Int_t    nStripsAboveCut = 0;
333         
334         for(UShort_t t = 0; t < nstr; t++) {
335           // nDistanceBefore++;
336           // nDistanceAfter++;
337
338           output.SetMultiplicity(d,r,s,t,0.);
339           Float_t mult         = SignalInStrip(input,d,r,s,t);
340           Float_t multNext     = (t<nstr-1) ? SignalInStrip(input,d,r,s,t+1) :0;
341           Float_t multNextNext = (t<nstr-2) ? SignalInStrip(input,d,r,s,t+2) :0;
342           if (multNext     ==  AliESDFMD::kInvalidMult) multNext     = 0;
343           if (multNextNext ==  AliESDFMD::kInvalidMult) multNextNext = 0;
344           if(!fThreeStripSharing) multNextNext = 0;
345
346           // Get the pseudo-rapidity 
347           Double_t eta = input.Eta(d,r,s,t);
348           Double_t phi = input.Phi(d,r,s,t) * TMath::Pi() / 180.;
349           if (s == 0) output.SetEta(d,r,s,t,eta);
350           
351           if(fRecalculateEta) { 
352             Double_t etaOld  = eta;
353             Double_t etaCalc = AliForwardUtil::GetEtaFromStrip(d,r,s,t,zvtx);
354             eta              = etaCalc;
355             
356             if (mult > 0 && mult != AliESDFMD::kInvalidMult ) {
357               Double_t cosOld =  ETA2COS(etaOld);
358               Double_t cosNew =  ETA2COS(etaCalc);
359               Double_t corr   =  cosNew / cosOld;
360               mult            *= corr;
361               multNext        *= corr;
362               multNextNext    *= corr;
363             }
364           } // Recalculate eta 
365
366           // Special case for pre revision 43611 AliFMDReconstructor.
367           // If fInvalidIsEmpty and we get an invalid signal from the
368           // ESD, then we need to set this signal to zero.  Note, dead
369           // strips added in the ForwardAODConfig.C file are not
370           // effected by this, and we can use that to set the proper
371           // dead strips.
372           if (mult == AliESDFMD::kInvalidMult && fInvalidIsEmpty) 
373             mult = 0;
374
375           // Keep dead-channel information - either from the ESD (but
376           // see above for older data) or from the settings in the
377           // ForwardAODConfig.C file.
378           if (mult == AliESDFMD::kInvalidMult || IsDead(d,r,s,t)) {
379             output.SetMultiplicity(d,r,s,t,AliESDFMD::kInvalidMult);
380             histos->fBefore->Fill(-1);
381             mult = AliESDFMD::kInvalidMult;
382           }
383           
384           if (mult != AliESDFMD::kInvalidMult) 
385             // Always fill the ESD sum histogram 
386             histos->fSumESD->Fill(eta, phi, mult);
387
388           // If no signal or dead strip, go on. 
389           if (mult == AliESDFMD::kInvalidMult || mult == 0) {
390             if (mult == 0) histos->fSum->Fill(eta,phi,mult);
391             // Flush a possible signal 
392             if (eTotal > 0 && t > 0) 
393               output.SetMultiplicity(d,r,s,t-1,eTotal);
394             // Reset states so we do not try to merge over a dead strip. 
395             eTotal = -1;
396             used   = false;
397             twoLow = false;
398             if (t > 0)  
399               histos->fNConsecutive->Fill(nStripsAboveCut);
400             if (mult == AliESDFMD::kInvalidMult)
401               // Why not fill immidiately here? 
402               nStripsAboveCut = -1;
403             else
404               // Why not fill immidiately here? 
405               nStripsAboveCut = 0;      
406             continue;
407           }
408
409           // Fill the diagnostics histogram 
410           histos->fBefore->Fill(mult);
411
412           Double_t mergedEnergy = mult;
413           // it seems to me that this logic could be condensed a bit
414           if(mult > GetLowCut(d, r, eta)) {               
415             if(nStripsAboveCut < 1) {
416               if(t > 0)
417                 histos->fNConsecutive->Fill(nStripsAboveCut);
418               nStripsAboveCut=0;
419             }
420             nStripsAboveCut++;
421           }     
422           else {
423             if (t > 0)
424               histos->fNConsecutive->Fill(nStripsAboveCut);
425             nStripsAboveCut=0;
426           }             
427
428           if (!fMergingDisabled) {
429             mergedEnergy = 0;
430
431             // The current sum
432             Float_t etot = 0;
433           
434             // Fill in neighbor information
435             if (t < nstr-1) histos->fNeighborsBefore->Fill(mult,multNext);
436
437             Bool_t thisValid = mult     > GetLowCut(d, r, eta);
438             Bool_t nextValid = multNext > GetLowCut(d, r, eta);
439             Bool_t thisSmall = mult     < GetHighCut(d, r, eta ,false);
440             Bool_t nextSmall = multNext < GetHighCut(d, r, eta ,false);
441           
442             // If this strips signal is above the high cut, reset distance
443             // if (!thisSmall) {
444             //    histos->fDistanceBefore->Fill(nDistanceBefore);
445             //    nDistanceBefore = -1;
446             // }
447           
448             // If the total signal in the past 1 or 2 strips are non-zero
449             // we need to check 
450             if (eTotal > 0) {
451               // Here, we have already flagged one strip as a candidate 
452             
453               // If 3-strip merging is enabled, then check the next 
454               // strip to see that it falls within cut, or if we have 
455               // two low signals 
456               if (fThreeStripSharing && nextValid && (nextSmall || twoLow)) {
457                 eTotal = eTotal + multNext;
458                 used = kTRUE;
459                 histos->fTriple->Fill(eTotal);
460                 nTriple++;
461                 twoLow = kFALSE;
462               }
463               // Otherwise, we got a double hit before, and that 
464               // should be stored. 
465               else {
466                 used = kFALSE;
467                 histos->fDouble->Fill(eTotal);
468                 nDouble++;
469               }
470               // Store energy loss and reset sum 
471               etot   = eTotal;
472               eTotal = -1;
473             } // if (eTotal>0)
474             else {
475               // If we have no current sum 
476             
477               // Check if this is marked as used, and if so, continue
478               if (used) {used = kFALSE; continue; }
479             
480               // If the signal is abvoe the cut, set current
481               if (thisValid) etot = mult;
482             
483               // If the signal is abiove the cut, and so is the next 
484               // signal and either of them are below the high cut, 
485               if (thisValid  && nextValid  && (thisSmall || nextSmall)) {
486               
487                 // If this is below the high cut, and the next is too, then 
488                 // we have two low signals 
489                 if (thisSmall && nextSmall) twoLow = kTRUE;
490               
491                 // If this signal is bigger than the next, and the 
492                 // one after that is below the low-cut, then update 
493                 // the sum
494                 if (mult>multNext && multNextNext < GetLowCut(d, r, eta)) {
495                   etot = mult + multNext;
496                   used = kTRUE;
497                   histos->fDouble->Fill(etot);
498                   nDouble++;
499                 }
500                 // Otherwise, we may need to merge with a third strip
501                 else {
502                   etot   = 0;
503                   eTotal = mult + multNext;
504                 }
505               }
506               // This is a signle hit 
507               else if(etot > 0) {
508                 histos->fSingle->Fill(etot);
509                 histos->fSinglePerStrip->Fill(etot,t);
510                 nSingle++;
511               }
512             } // else if (etotal >= 0)
513           
514             mergedEnergy = etot;
515             // if (mergedEnergy > GetHighCut(d, r, eta ,false)) {
516             //   histos->fDistanceAfter->Fill(nDistanceAfter);
517             //   nDistanceAfter    = -1;
518             // }
519             //if(mult>0 && multNext >0)
520             //  std::cout<<mult<<"  "<<multNext<<"  "<<mergedEnergy<<std::endl;
521           } // if (!fMergingDisabled)
522
523           if (!fCorrectAngles)
524             mergedEnergy = AngleCorrect(mergedEnergy, eta);
525           // if (mergedEnergy > 0) histos->Incr();
526           
527           if (t != 0) 
528             histos->fNeighborsAfter->Fill(output.Multiplicity(d,r,s,t-1), 
529                                           mergedEnergy);
530           histos->fBeforeAfter->Fill(mult, mergedEnergy);
531           if(mergedEnergy > 0)
532             histos->fAfter->Fill(mergedEnergy);
533           histos->fSum->Fill(eta,phi,mergedEnergy);
534           
535           output.SetMultiplicity(d,r,s,t,mergedEnergy);
536         } // for strip
537         histos->fNConsecutive->Fill(nStripsAboveCut); // fill the last sector 
538       } // for sector
539     } // for ring 
540   } // for detector
541   DMSG(fDebug, 3,"single=%9d, double=%9d, triple=%9d", 
542        nSingle, nDouble, nTriple);
543   next.Reset();
544   // while ((o = static_cast<RingHistos*>(next()))) o->Finish();
545
546   return kTRUE;
547 }
548
549 //_____________________________________________________________________
550 Double_t 
551 AliFMDSharingFilter::SignalInStrip(const AliESDFMD& input, 
552                                    UShort_t         d,
553                                    Char_t           r,
554                                    UShort_t         s,
555                                    UShort_t         t) const
556 {
557   // 
558   // Get the signal in a strip 
559   // 
560   // Parameters:
561   //    fmd   ESD object
562   //    d     Detector
563   //    r     Ring 
564   //    s     Sector 
565   //    t     Strip
566   // 
567   // Return:
568   //    The energy signal 
569   //
570   Double_t mult = input.Multiplicity(d,r,s,t);
571   // In case of 
572   //  - bad value (invalid or 0) 
573   //  - we want angle corrected and data is 
574   //  - we don't want angle corrected and data isn't 
575   // just return read value  
576   if (mult == AliESDFMD::kInvalidMult               || 
577       mult == 0                                     ||
578       (fCorrectAngles && (fIgnoreESDForAngleCorrection || input.IsAngleCorrected())) || 
579       (!fCorrectAngles && !fIgnoreESDForAngleCorrection && !input.IsAngleCorrected()))
580     return mult;
581
582   // If we want angle corrected data, correct it, 
583   // otherwise de-correct it 
584   if (fCorrectAngles) mult = AngleCorrect(mult, input.Eta(d,r,s,t));
585   else                mult = DeAngleCorrect(mult, input.Eta(d,r,s,t));
586   return mult;
587 }
588
589 namespace {
590   Double_t Rng2Cut(UShort_t d, Char_t r, Double_t eta, TH2* h) {
591     Double_t ret = 1024;
592     Int_t ybin = 0;                                                     
593     switch(d) {                                                         
594     case 1: ybin = 1; break;                                            
595     case 2: ybin = (r=='i' || r=='I') ? 2 : 3; break;                   
596     case 3: ybin = (r=='i' || r=='I') ? 4 : 5; break;                   
597     default: return ret;
598     }                                                                   
599     Int_t xbin = h->GetXaxis()->FindBin(eta);                           
600     if (xbin < 1 && xbin > h->GetXaxis()->GetNbins()) return ret;
601     ret = h->GetBinContent(xbin,ybin);                                  
602     return ret;
603   }
604 }
605     
606   //_____________________________________________________________________
607 Double_t 
608 AliFMDSharingFilter::GetLowCut(UShort_t d, Char_t r, Double_t eta) const
609 {
610   //
611   // Get the low cut.  Normally, the low cut is taken to be the lower
612   // value of the fit range used when generating the energy loss fits.
613   // However, if fLowCut is set (using SetLowCit) to a value greater
614   // than 0, then that value is used.
615   //
616   return Rng2Cut(d, r, eta, fLowCuts);
617   // return fLCuts.GetMultCut(d,r,eta,false);
618 }
619                         
620 //_____________________________________________________________________
621 Double_t 
622 AliFMDSharingFilter::GetHighCut(UShort_t d, Char_t r, 
623                                 Double_t eta, Bool_t /*errors*/) const
624 {
625   //
626   // Get the high cut.  The high cut is defined as the 
627   // most-probably-value peak found from the energy distributions, minus 
628   // 2 times the width of the corresponding Landau.
629   //
630   return Rng2Cut(d, r, eta, fHighCuts);
631   // return fHCuts.GetMultCut(d,r,eta,errors); 
632 }
633
634 //____________________________________________________________________
635 Double_t
636 AliFMDSharingFilter::AngleCorrect(Double_t mult, Double_t eta) const
637 {
638   // 
639   // Angle correct the signal 
640   // 
641   // Parameters:
642   //    mult Angle Un-corrected Signal 
643   //    eta  Pseudo-rapidity 
644   // 
645   // Return:
646   //    Angle corrected signal 
647   //
648   Double_t theta =  2 * TMath::ATan(TMath::Exp(-eta));
649   if (eta < 0) theta -= TMath::Pi();
650   return mult * TMath::Cos(theta);
651 }
652 //____________________________________________________________________
653 Double_t
654 AliFMDSharingFilter::DeAngleCorrect(Double_t mult, Double_t eta) const
655 {
656   // 
657   // Angle de-correct the signal 
658   // 
659   // Parameters:
660   //    mult Angle corrected Signal 
661   //    eta  Pseudo-rapidity 
662   // 
663   // Return:
664   //    Angle un-corrected signal 
665   //
666   Double_t theta =  2 * TMath::ATan(TMath::Exp(-eta));
667   if (eta < 0) theta -= TMath::Pi();
668   return mult / TMath::Cos(theta);
669 }
670
671 //____________________________________________________________________
672 void
673 AliFMDSharingFilter::Terminate(const TList* dir, TList* output, Int_t nEvents)
674 {
675   // 
676   // Scale the histograms to the total number of events 
677   // 
678   // Parameters:
679   //    dir     Where the output is 
680   //    nEvents Number of events 
681   //
682   DGUARD(fDebug,1, "Scale histograms in AliFMDSharingFilter");
683   if (nEvents <= 0) return;
684   TList* d = static_cast<TList*>(dir->FindObject(GetName()));
685   if (!d) return;
686
687   TList* out = new TList;
688   out->SetName(d->GetName());
689   out->SetOwner();
690
691   TParameter<int>* nFiles = 
692     static_cast<TParameter<int>*>(d->FindObject("nFiles"));
693
694   TH2* lowCuts  = static_cast<TH2*>(d->FindObject("lowCuts"));
695   TH2* highCuts = static_cast<TH2*>(d->FindObject("highCuts"));
696   if (lowCuts && nFiles) {
697     lowCuts->Scale(1. / nFiles->GetVal());
698     out->Add(lowCuts->Clone());
699   }
700   else 
701     AliWarning("low cuts histogram not found in input list");
702   if (highCuts && nFiles) {
703     highCuts->Scale(1. / nFiles->GetVal());
704     out->Add(highCuts->Clone());
705   }
706   else 
707     AliWarning("high cuts histogram not found in input list");
708   
709   TIter    next(&fRingHistos);
710   RingHistos* o = 0;
711   THStack* sums = new THStack("sums", "Sum of ring signals");
712   THStack* sumsESD = new THStack("sumsESD", "Sum of ring ESD signals");
713   while ((o = static_cast<RingHistos*>(next()))) {
714     o->Terminate(d, nEvents);
715     if (!o->fSum) { 
716       Warning("Terminate", "No sum histogram found for ring %s", o->GetName());
717       continue;
718     }
719     TH1D* sum = o->fSum->ProjectionX(o->GetName(), 1, o->fSum->GetNbinsY(),"e");
720     sum->Scale(1., "width");
721     sum->SetTitle(o->GetName());
722     sum->SetDirectory(0);
723     sum->SetYTitle("#sum_{c} #Delta/#Delta_{mip}");
724     sums->Add(sum);
725
726     sum = o->fSumESD->ProjectionX(o->GetName(), 1, o->fSumESD->GetNbinsY(),"e");
727     sum->Scale(1., "width");
728     sum->SetTitle(o->GetName());
729     sum->SetDirectory(0);
730     sum->SetYTitle("#sum_{s} #Delta/#Delta_{mip}");
731     sumsESD->Add(sum);
732   }
733   out->Add(sums);
734   out->Add(sumsESD);
735   output->Add(out);
736 }
737
738 //____________________________________________________________________
739 void
740 AliFMDSharingFilter::CreateOutputObjects(TList* dir)
741 {
742   // 
743   // Define the output histograms.  These are put in a sub list of the
744   // passed list.   The histograms are merged before the parent task calls 
745   // AliAnalysisTaskSE::Terminate 
746   // 
747   // Parameters:
748   //    dir Directory to add to 
749   //
750   DGUARD(fDebug,1, "Define output in AliFMDSharingFilter");
751   TList* d = new TList;
752   d->SetName(GetName());
753   dir->Add(d);
754
755 #if 0
756   fSummed = new TH2I("operations", "Strip operations", 
757                      kMergedInto, kNone-.5, kMergedInto+.5,
758                      51201, -.5, 51200.5);
759   fSummed->SetXTitle("Operation");
760   fSummed->SetYTitle("# of strips");
761   fSummed->SetZTitle("Events");
762   fSummed->GetXaxis()->SetBinLabel(kNone,            "None");
763   fSummed->GetXaxis()->SetBinLabel(kCandidate,       "Candidate");
764   fSummed->GetXaxis()->SetBinLabel(kMergedWithOther, "Merged w/other");
765   fSummed->GetXaxis()->SetBinLabel(kMergedInto,      "Merged into");
766   fSummed->SetDirectory(0);
767   d->Add(fSummed);
768 #endif
769
770   fHighCuts = new TH2D("highCuts", "High cuts used", 1,0,1, 1,0,1);
771   fHighCuts->SetXTitle("#eta");
772   fHighCuts->SetDirectory(0);
773   d->Add(fHighCuts);
774
775   fLowCuts = new TH2D("lowCuts", "Low cuts used", 1,0,1, 1,0,1);
776   fLowCuts->SetXTitle("#eta");
777   fLowCuts->SetDirectory(0);
778   d->Add(fLowCuts);
779
780   // d->Add(lowCut);
781   // d->Add(nXi);
782   // d->Add(sigma);
783   d->Add(AliForwardUtil::MakeParameter("angle", fCorrectAngles));
784   d->Add(AliForwardUtil::MakeParameter("lowSignal", 
785                                        fZeroSharedHitsBelowThreshold));
786   d->Add(AliForwardUtil::MakeParameter("simple", fUseSimpleMerging));
787   d->Add(AliForwardUtil::MakeParameter("sumThree", fThreeStripSharing));
788   d->Add(AliForwardUtil::MakeParameter("disabled", fMergingDisabled));
789   TParameter<int>* nFiles = new TParameter<int>("nFiles", 1);
790   nFiles->SetMergeMode('+');
791   d->Add(nFiles);
792   
793   TObjArray* extraDead = new TObjArray;
794   extraDead->SetOwner();
795   extraDead->SetName("extraDead");
796 #if 0
797   for (Int_t i = 0; i < fExtraDead.GetSize(); i++) { 
798     if (fExtraDead.At(i) < 0) break;
799     UShort_t dd, s, t;
800     Char_t  r;
801     Int_t   id = fExtraDead.At(i);
802     AliFMDStripIndex::Unpack(id, dd, r, s, t);
803     extraDead->Add(AliForwardUtil::MakeParameter(Form("FMD%d%c[%02d,%03d]",
804                                                       dd, r, s, t), id));
805   }
806 #endif
807   fXtraDead.Compact();
808   UInt_t firstBit = fXtraDead.FirstSetBit();
809   UInt_t nBits    = fXtraDead.GetNbits();
810   for (UInt_t i = firstBit; i < nBits; i++) {
811     if (!fXtraDead.TestBitNumber(i)) continue;
812     UShort_t dd, s, t;
813     Char_t  r;
814     AliFMDStripIndex::Unpack(i, dd, r, s, t);
815     extraDead->Add(AliForwardUtil::MakeParameter(Form("FMD%d%c[%02d,%03d]",
816                                                       dd, r, s, t), 
817                                                  UShort_t(i)));
818   }
819   // d->Add(&fXtraDead);
820   d->Add(extraDead);
821   fLCuts.Output(d,"lCuts");
822   fHCuts.Output(d,"hCuts");
823
824   TIter    next(&fRingHistos);
825   RingHistos* o = 0;
826   while ((o = static_cast<RingHistos*>(next()))) {
827     o->CreateOutputObjects(d);
828   }
829 }
830 #define PF(N,V,...)                                     \
831   AliForwardUtil::PrintField(N,V, ## __VA_ARGS__)
832 #define PFB(N,FLAG)                             \
833   do {                                                                  \
834     AliForwardUtil::PrintName(N);                                       \
835     std::cout << std::boolalpha << (FLAG) << std::noboolalpha << std::endl; \
836   } while(false)
837 #define PFV(N,VALUE)                                    \
838   do {                                                  \
839     AliForwardUtil::PrintName(N);                       \
840     std::cout << (VALUE) << std::endl; } while(false)
841
842 //____________________________________________________________________
843 void
844 AliFMDSharingFilter::Print(Option_t* /*option*/) const
845 {
846   // 
847   // Print information
848   // 
849   // Parameters:
850   //    option Not used 
851   //
852   AliForwardUtil::PrintTask(*this);
853   gROOT->IncreaseDirLevel();
854
855   PFB("Use corrected angles",  fCorrectAngles);
856   PFB("Zero below threshold",  fZeroSharedHitsBelowThreshold);
857   PFB("Use simple sharing",    fUseSimpleMerging);
858   PFB("Consider invalid null", fInvalidIsEmpty);
859   PFB("Allow 3 strip merging", fThreeStripSharing);
860   PF("Low cuts",        "");
861   fLCuts.Print();
862   PF("High cuts",       "");
863   fHCuts.Print();
864   gROOT->DecreaseDirLevel();
865 }
866   
867 //====================================================================
868 AliFMDSharingFilter::RingHistos::RingHistos()
869   : AliForwardUtil::RingHistos(), 
870     fBefore(0), 
871     fAfter(0), 
872     fSingle(0),
873     fDouble(0),
874     fTriple(0),
875     fSinglePerStrip(0),
876     // fDistanceBefore(0),
877     // fDistanceAfter(0),
878     fBeforeAfter(0),
879     fNeighborsBefore(0),
880     fNeighborsAfter(0),
881     fSumESD(0),
882     fSum(0),
883     fNConsecutive(0)    
884      // ,
885     // fHits(0),
886     // fNHits(0)
887 {
888   // 
889   // Default CTOR
890   //
891
892 }
893
894 //____________________________________________________________________
895 AliFMDSharingFilter::RingHistos::RingHistos(UShort_t d, Char_t r)
896   : AliForwardUtil::RingHistos(d,r), 
897     fBefore(0), 
898     fAfter(0),
899     fSingle(0),
900     fDouble(0),
901     fTriple(0),    
902     fSinglePerStrip(0),
903     // fDistanceBefore(0),
904     // fDistanceAfter(0),
905     fBeforeAfter(0),
906     fNeighborsBefore(0),
907     fNeighborsAfter(0),
908     fSumESD(0),
909     fSum(0),
910     fNConsecutive(0)    
911      //,
912     // fHits(0),
913     // fNHits(0)
914 {
915   // 
916   // Constructor
917   // 
918   // Parameters:
919   //    d detector
920   //    r ring 
921   //
922   fBefore = new TH1D("esdEloss", Form("Energy loss in %s (reconstruction)", 
923                                       GetName()), 640, -1, 15);
924   fBefore->SetXTitle("#Delta E/#Delta E_{mip}");
925   fBefore->SetYTitle("P(#Delta E/#Delta E_{mip})");
926   fBefore->SetFillColor(Color());
927   fBefore->SetFillStyle(3001);
928   fBefore->SetLineColor(kBlack);
929   fBefore->SetLineStyle(2);
930   fBefore->SetDirectory(0);
931
932   fAfter  = static_cast<TH1D*>(fBefore->Clone("anaEloss"));
933   fAfter->SetTitle(Form("Energy loss in %s (sharing corrected)", GetName()));
934   fAfter->SetFillColor(Color()+2);
935   fAfter->SetLineStyle(1);
936   fAfter->SetDirectory(0);
937   
938   fSingle = new TH1D("singleEloss", "Energy loss (single strips)", 
939                      600, 0, 15);
940   fSingle->SetXTitle("#Delta/#Delta_{mip}");
941   fSingle->SetYTitle("P(#Delta/#Delta_{mip})");
942   fSingle->SetFillColor(Color());
943   fSingle->SetFillStyle(3001);
944   fSingle->SetLineColor(kBlack);
945   fSingle->SetLineStyle(2);
946   fSingle->SetDirectory(0);
947
948   fDouble = static_cast<TH1D*>(fSingle->Clone("doubleEloss"));
949   fDouble->SetTitle("Energy loss (two strips)");
950   fDouble->SetFillColor(Color()+1);
951   fDouble->SetDirectory(0);
952   
953   fTriple = static_cast<TH1D*>(fSingle->Clone("tripleEloss"));
954   fTriple->SetTitle("Energy loss (three strips)"); 
955   fTriple->SetFillColor(Color()+2);
956   fTriple->SetDirectory(0);
957   
958   //Int_t nBinsForInner = (r == 'I' ? 32 : 16);
959   Int_t nBinsForInner = (r == 'I' ? 512 : 256);
960   Int_t nStrips       = (r == 'I' ? 512 : 256);
961   
962   fSinglePerStrip = new TH2D("singlePerStrip", "SinglePerStrip", 
963                              600,0,15, nBinsForInner,0,nStrips);
964   fSinglePerStrip->SetXTitle("#Delta/#Delta_{mip}");
965   fSinglePerStrip->SetYTitle("Strip #");
966   fSinglePerStrip->SetZTitle("Counts");
967   fSinglePerStrip->SetDirectory(0);
968
969 #if 0
970   fDistanceBefore = new TH1D("distanceBefore", "Distance before sharing", 
971                              nStrips , 0,nStrips );
972   fDistanceBefore->SetXTitle("Distance");
973   fDistanceBefore->SetYTitle("Counts");
974   fDistanceBefore->SetFillColor(kGreen+2);
975   fDistanceBefore->SetFillStyle(3001);
976   fDistanceBefore->SetLineColor(kBlack);
977   fDistanceBefore->SetLineStyle(2);
978   fDistanceBefore->SetDirectory(0);
979
980   fDistanceAfter = static_cast<TH1D*>(fDistanceBefore->Clone("distanceAfter"));
981   fDistanceAfter->SetTitle("Distance after sharing"); 
982   fDistanceAfter->SetFillColor(kGreen+1);
983   fDistanceAfter->SetDirectory(0);
984 #endif
985   
986   Double_t max = 15;
987   Double_t min = -1;
988   Int_t    n   = int((max-min) / (max / 300));
989   fBeforeAfter = new TH2D("beforeAfter", "Before and after correlation", 
990                           n, min, max, n, min, max);
991   fBeforeAfter->SetXTitle("#Delta E/#Delta E_{mip} before");
992   fBeforeAfter->SetYTitle("#Delta E/#Delta E_{mip} after");
993   fBeforeAfter->SetZTitle("Correlation");
994   fBeforeAfter->SetDirectory(0);
995
996   fNeighborsBefore = static_cast<TH2D*>(fBeforeAfter->Clone("neighborsBefore"));
997   fNeighborsBefore->SetTitle("Correlation of neighbors before");
998   fNeighborsBefore->SetXTitle("#Delta E_{i}/#Delta E_{mip}");
999   fNeighborsBefore->SetYTitle("#Delta E_{i+1}/#Delta E_{mip}");
1000   fNeighborsBefore->SetDirectory(0);
1001   
1002   fNeighborsAfter = 
1003     static_cast<TH2D*>(fNeighborsBefore->Clone("neighborsAfter"));
1004   fNeighborsAfter->SetTitle("Correlation of neighbors after");
1005   fNeighborsAfter->SetDirectory(0);
1006
1007   fSumESD = new TH2D("summedESD", "Summed ESD signal", 200, -4, 6, 
1008                   NSector(), 0, 2*TMath::Pi());
1009   fSumESD->SetDirectory(0);
1010   fSumESD->Sumw2();
1011   fSumESD->SetMarkerColor(Color());
1012   // fSum->SetFillColor(Color());
1013   fSumESD->SetXTitle("#eta");
1014   fSumESD->SetYTitle("#varphi [radians]");
1015   fSumESD->SetZTitle("#sum_{strip} #Delta/#Delta_{mip}(#eta,#varphi) ");
1016
1017   fSum = static_cast<TH2D*>(fSumESD->Clone("summed"));
1018   fSum->SetTitle("Summed cluster signal");
1019   fSum->SetZTitle("#sum_{cluster} #Delta/#Delta_{mip}(#eta,#varphi) ");
1020   fSum->SetDirectory(0);
1021  
1022   // Perhaps we need to ensure that this histogram has enough range to
1023   // accommondate all possible ranges - that is, from -1 to the number
1024   // of strips in this ring(-type) - i.e., NStrips().  Perhaps the
1025   // axis should be defined with increasin bin size - e.g.,
1026   //
1027   //   -1.5,-.5,.5,1.5,...,100.5,128.5,192.5,...,NStrips()
1028   // 
1029   fNConsecutive = new TH1D("nConsecutive","# consecutive strips above low cut",
1030                            201,-1.5,199.5);
1031   fNConsecutive->SetXTitle("N_{strips}");
1032   fNConsecutive->SetYTitle("N_{entries}");
1033   fNConsecutive->SetFillColor(kYellow+2);
1034   fNConsecutive->SetFillStyle(3001); 
1035   fNConsecutive->SetDirectory(0);
1036   
1037  
1038 #if 0  
1039   fHits = new TH1D("hits", "Number of hits", 200, 0, 200000);
1040   fHits->SetDirectory(0);
1041   fHits->SetXTitle("# of hits");
1042   fHits->SetFillColor(kGreen+1);
1043 #endif
1044 }
1045 //____________________________________________________________________
1046 AliFMDSharingFilter::RingHistos::RingHistos(const RingHistos& o)
1047   : AliForwardUtil::RingHistos(o), 
1048     fBefore(o.fBefore), 
1049     fAfter(o.fAfter),
1050     fSingle(o.fSingle),
1051     fDouble(o.fDouble),
1052     fTriple(o.fTriple),
1053     fSinglePerStrip(o.fSinglePerStrip),
1054     // fDistanceBefore(o.fDistanceBefore),
1055     // fDistanceAfter(o.fDistanceAfter),    
1056     fBeforeAfter(o.fBeforeAfter),
1057     fNeighborsBefore(o.fNeighborsBefore),
1058     fNeighborsAfter(o.fNeighborsAfter),
1059     fSumESD(o.fSumESD), //,
1060     fSum(o.fSum),
1061     fNConsecutive(o.fNConsecutive)
1062      //,
1063     // fHits(o.fHits),
1064     // fNHits(o.fNHits)
1065 {
1066   // 
1067   // Copy constructor 
1068   // 
1069   // Parameters:
1070   //    o Object to copy from 
1071   //
1072 }
1073
1074 //____________________________________________________________________
1075 AliFMDSharingFilter::RingHistos&
1076 AliFMDSharingFilter::RingHistos::operator=(const RingHistos& o)
1077 {
1078   // 
1079   // Assignment operator 
1080   // 
1081   // Parameters:
1082   //    o Object to assign from 
1083   // 
1084   // Return:
1085   //    Reference to this 
1086   //
1087   if (&o == this) return *this;
1088   AliForwardUtil::RingHistos::operator=(o);
1089   fDet = o.fDet;
1090   fRing = o.fRing;
1091   
1092   if (fBefore)         delete  fBefore;
1093   if (fAfter)          delete  fAfter;
1094   if (fSingle)         delete  fSingle;
1095   if (fDouble)         delete  fDouble;
1096   if (fTriple)         delete  fTriple;
1097   if (fSinglePerStrip) delete  fSinglePerStrip;
1098   if (fNConsecutive)   delete  fNConsecutive;
1099   // if (fDistanceBefore) delete fDistanceBefore;
1100   // if (fDistanceAfter)  delete fDistanceAfter;
1101   // if (fHits)                delete fHits;
1102   
1103   
1104   fBefore          = static_cast<TH1D*>(o.fBefore->Clone());
1105   fAfter           = static_cast<TH1D*>(o.fAfter->Clone());
1106   fSingle          = static_cast<TH1D*>(o.fSingle->Clone());
1107   fDouble          = static_cast<TH1D*>(o.fDouble->Clone());
1108   fTriple          = static_cast<TH1D*>(o.fTriple->Clone());
1109   fSinglePerStrip  = static_cast<TH2D*>(o.fSinglePerStrip->Clone());
1110   // fDistanceBefore  = static_cast<TH1D*>(o.fDistanceBefore->Clone());
1111   // fDistanceAfter   = static_cast<TH1D*>(o.fDistanceAfter->Clone());
1112   fBeforeAfter     = static_cast<TH2D*>(o.fBeforeAfter->Clone());
1113   fNeighborsBefore = static_cast<TH2D*>(o.fNeighborsBefore->Clone());
1114   fNeighborsAfter  = static_cast<TH2D*>(o.fNeighborsAfter->Clone());
1115   // fHits            = static_cast<TH1D*>(o.fHits->Clone());
1116   fSumESD          = static_cast<TH2D*>(o.fSumESD->Clone());
1117   fSum             = static_cast<TH2D*>(o.fSum->Clone());
1118   fNConsecutive    = static_cast<TH1D*>(o.fNConsecutive->Clone());
1119
1120   return *this;
1121 }
1122 //____________________________________________________________________
1123 AliFMDSharingFilter::RingHistos::~RingHistos()
1124 {
1125   // 
1126   // Destructor 
1127   //
1128 }
1129 #if 0
1130 //____________________________________________________________________
1131 void
1132 AliFMDSharingFilter::RingHistos::Finish()
1133 {
1134   // 
1135   // Finish off 
1136   // 
1137   //
1138   // fHits->Fill(fNHits);
1139 }
1140 #endif
1141 //____________________________________________________________________
1142 void
1143 AliFMDSharingFilter::RingHistos::Terminate(const TList* dir, Int_t nEvents)
1144 {
1145   // 
1146   // Scale the histograms to the total number of events 
1147   // 
1148   // Parameters:
1149   //    nEvents Number of events 
1150   //    dir     Where the output is 
1151   //
1152   TList* l = GetOutputList(dir);
1153   if (!l) return; 
1154
1155 #if 0
1156   TH1D* before = static_cast<TH1D*>(l->FindObject("esdEloss"));
1157   TH1D* after  = static_cast<TH1D*>(l->FindObject("anaEloss"));
1158   if (before) before->Scale(1./nEvents);
1159   if (after)  after->Scale(1./nEvents);
1160 #endif
1161
1162   TH2D* summed = static_cast<TH2D*>(l->FindObject("summed"));
1163   if (summed) summed->Scale(1./nEvents);
1164   fSum = summed;
1165
1166   TH2D* summedESD = static_cast<TH2D*>(l->FindObject("summedESD"));
1167   if (summedESD) summedESD->Scale(1./nEvents);
1168   fSumESD = summedESD;
1169
1170   TH1D* consecutive = static_cast<TH1D*>(l->FindObject("nConsecutive"));
1171   if (consecutive) consecutive->Scale(1./nEvents);
1172   fNConsecutive= consecutive;
1173 }
1174
1175 //____________________________________________________________________
1176 void
1177 AliFMDSharingFilter::RingHistos::CreateOutputObjects(TList* dir)
1178 {
1179   // 
1180   // Make output 
1181   // 
1182   // Parameters:
1183   //    dir where to store 
1184   //
1185   TList* d = DefineOutputList(dir);
1186
1187   d->Add(fBefore);
1188   d->Add(fAfter);
1189   d->Add(fSingle);
1190   d->Add(fDouble);
1191   d->Add(fTriple);
1192   d->Add(fSinglePerStrip);
1193   // d->Add(fDistanceBefore);
1194   // d->Add(fDistanceAfter);
1195   d->Add(fBeforeAfter);
1196   d->Add(fNeighborsBefore);
1197   d->Add(fNeighborsAfter);
1198   // d->Add(fHits);
1199   d->Add(fSumESD);
1200   d->Add(fSum);
1201   d->Add(fNConsecutive);
1202
1203   // Removed to avoid doubly adding the list which destroys 
1204   // the merging
1205   //dir->Add(d);
1206 }
1207
1208 //____________________________________________________________________
1209 //
1210 // EOF
1211 //