]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGPP/MUON/dep/AliAnalysisTaskMuonResolution.h
Merge branch 'TPCdev' of https://git.cern.ch/reps/AliRoot into TPCdev
[u/mrichter/AliRoot.git] / PWGPP / MUON / dep / AliAnalysisTaskMuonResolution.h
1 #ifndef ALIANALYSISTASKMUONRESOLUTION_H
2 #define ALIANALYSISTASKMUONRESOLUTION_H
3 /* Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4  * See cxx source for full Copyright notice                               */
5
6 /* $Id$ */ 
7
8 /// \ingroup muondep
9 /// \class AliAnalysisTaskMuonResolution
10 /// \brief Muon spectrometer resolution
11 //Author: Philippe Pillot - SUBATECH Nantes
12
13 #include <TString.h>
14 #include <TMatrixD.h>
15 #include <TF1.h>
16 #include "AliAnalysisTaskSE.h"
17 #include "AliMuonTrackCuts.h"
18
19 class TH1;
20 class TH2;
21 class TGraphErrors;
22 class TObjArray;
23 class TList;
24 class AliMUONTrack;
25 class AliMUONTrackParam;
26 class AliMUONGeometryTransformer;
27 class AliMUONVCluster;
28
29 class AliAnalysisTaskMuonResolution : public AliAnalysisTaskSE {
30 public:
31   
32   AliAnalysisTaskMuonResolution();
33   AliAnalysisTaskMuonResolution(const char *name);
34   virtual ~AliAnalysisTaskMuonResolution();
35   
36   /// Set location of the default OCDB storage (if not set use "raw://")
37   void SetDefaultStorage(const char* ocdbPath) { fDefaultStorage = ocdbPath; }
38   
39   void SetStartingResolution(Int_t chId, Double_t valNB, Double_t valB);
40   void SetStartingResolution(Double_t valNB[10], Double_t valB[10]);
41   void GetStartingResolution(Double_t valNB[10], Double_t valB[10]) const;
42   
43   void SetHalfChShift(Int_t hchId, Double_t valNB, Double_t valB);
44   void SetHalfChShift(Double_t valNB[20], Double_t valB[20]);
45   void GetHalfChShift(Double_t valNB[20], Double_t valB[20]) const;
46   void ShiftHalfCh(Bool_t flag = kTRUE) { fShiftHalfCh = flag; }
47   void PrintHalfChShift(Bool_t flag = kTRUE) { fPrintHalfChShift = flag; }
48   
49   void SetDEShift(Int_t iDE, Double_t valNB, Double_t valB);
50   void SetDEShift(Double_t valNB[200], Double_t valB[200]);
51   void GetDEShift(Double_t valNB[200], Double_t valB[200]) const;
52   void ShiftDE(Bool_t flag = kTRUE) { fShiftDE = flag; }
53   void PrintDEShift(Bool_t flag = kTRUE) { fPrintDEShift = flag; }
54   
55   /// set the minimum momentum value of the tracks used to compute the resolution
56   void SetMinMomentum(Double_t val) { fMinMomentum = val; }
57   
58   /// set the flag to use only tracks passing the physics selection
59   void SelectPhysics(Bool_t flag = kTRUE) {fSelectPhysics = flag;}
60   
61   /// set the flag to use only tracks matched with trigger or not
62   void MatchTrigger(Bool_t flag = kTRUE) { fMatchTrig = flag; }
63   
64   /// set the flag to use only tracks passing the acceptance cuts (Rabs, eta)
65   void ApplyAccCut(Bool_t flag = kTRUE) { fApplyAccCut = flag; }
66   
67   // set standard cuts to select tracks to be considered
68   void SetMuonTrackCuts(AliMuonTrackCuts &trackCuts);
69   
70   /// Select events belonging to at least one of the trigger classes selected by the mask to fill histograms:
71   /// - if the physics selection is used, apply the mask to the trigger word returned by the physics selection
72   /// - if not, apply the mask to the trigger word built by looking for triggers listed in "fSelectTriggerClass"
73   void SelectTrigger(Bool_t flag = kTRUE, UInt_t mask = AliVEvent::kMUON) {fSelectTrigger = flag; fTriggerMask = mask;}
74   
75   /// set the extrapolation mode to get the track parameters and covariances at a given cluster:
76   /// 0 = extrapolate from the closest cluster; 1 = extrapolate from the previous cluster except between stations 2-3-4
77   void SetExtrapMode(Int_t val) { fExtrapMode = val; }
78   
79   /// set the flag to add or not the systematic shifts of the residuals to the resolution
80   void CorrectForSystematics(Bool_t flag = kTRUE) { fCorrectForSystematics = flag; }
81   
82   /// Set the OCDB path to the alignment file used in the reco (if not set use default storage)
83   void SetAlignStorage(const char* ocdbPath) { fNewAlignStorage = ocdbPath; }
84   
85   void ReAlign(const char* oldAlignStorage = 0x0, const char* newAlignStorage = "");
86   
87   /// return the list of summary canvases
88   TObjArray* GetCanvases() {return fCanvases;}
89   
90   /// set the flag to show the progression bar
91   void ShowProgressBar(Bool_t flag = kTRUE) {fShowProgressBar = flag;}
92   
93   /// set the flag to print the cluster resolution per chamber/DE
94   void PrintClusterRes(Bool_t perCh = kTRUE, Bool_t perDE = kFALSE) {fPrintClResPerCh = perCh; fPrintClResPerDE = perDE;}
95   
96   void FitResiduals(Bool_t flag = kTRUE);
97   
98   /// set the flag to remove mono-cathod clusters (either considering all the pads or only the ones directly below)
99   void RemoveMonoCathodClusters(Bool_t flag = kTRUE, Bool_t checkAllPads = kTRUE) {fRemoveMonoCathCl = flag; fCheckAllPads = checkAllPads;}
100   
101   /// set the flag to improve the track before measuring the resolution
102   void ImproveTracks(Bool_t flag = kTRUE) {fImproveTracks = flag;}
103   
104   virtual void   UserCreateOutputObjects();
105   virtual void   UserExec(Option_t *);
106   virtual void   NotifyRun();
107   virtual void   Terminate(Option_t *);
108   
109 private:
110   
111   /// Not implemented
112   AliAnalysisTaskMuonResolution(const AliAnalysisTaskMuonResolution& rhs);
113   /// Not implemented
114   AliAnalysisTaskMuonResolution& operator = (const AliAnalysisTaskMuonResolution& rhs);
115   
116   void ModifyClusters(AliMUONTrack& track);
117   void Zoom(TH1* h, Double_t fractionCut = 0.01);
118   void ZoomLeft(TH1* h, Double_t fractionCut = 0.02);
119   void ZoomRight(TH1* h, Double_t fractionCut = 0.02);
120   void GetMeanRMS(TH1* h, Double_t& mean, Double_t& meanErr,Double_t& rms, Double_t& rmsErr, TGraphErrors* gMean = 0x0,
121                   TGraphErrors* gRMS = 0x0, Int_t i = 0, Double_t x = 0, Bool_t zoom = kTRUE, Bool_t enableFit = kTRUE);
122   void FillMeanSigmaClusterVsX(const TH2* hIn, const TH2* hOut, TGraphErrors* gMean, TGraphErrors* gSigma);
123   void Cov2CovP(const AliMUONTrackParam &param, TMatrixD &covP);
124   UInt_t BuildTriggerWord(const TString& FiredTriggerClasses);
125   void CheckPads(AliMUONVCluster *cl, Bool_t &hasBending, Bool_t &hasNonBending) const;
126   void CheckPadsBelow(AliMUONVCluster *cl, Bool_t &hasBending, Bool_t &hasNonBending) const;
127   
128 private:
129   
130   enum eResiduals {
131     kResidualPerChClusterIn             = 0,  ///< cluster-track residual-X/Y distribution per chamber (cluster attached to the track)
132     kResidualPerChClusterOut            = 2,  ///< cluster-track residual-X/Y distribution per chamber (cluster not attached to the track)
133     kTrackResPerCh                      = 4,  ///< track resolution-X/Y per chamber
134     kMCSPerCh                           = 6,  ///< MCS X/Y-dispersion of extrapolated track per chamber
135     kClusterRes2PerCh                   = 8,  ///< cluster X/Y-resolution per chamber
136     kResidualPerDEClusterIn             = 10, ///< cluster-track residual-X/Y distribution per DE (cluster attached to the track)
137     kResidualPerDEClusterOut            = 12, ///< cluster-track residual-X/Y distribution per DE (cluster not attached to the track)
138     kTrackResPerDE                      = 14, ///< track resolution-X/Y per DE
139     kMCSPerDE                           = 16, ///< MCS X/Y-dispersion of extrapolated track per DE
140     kResidualPerHalfChClusterIn         = 18, ///< cluster-track residual-X/Y distribution per half chamber (cluster attached to the track)
141     kResidualPerHalfChClusterOut        = 20, ///< cluster-track residual-X/Y distribution per half chamber (cluster not attached to the track)
142     kTrackResPerHalfCh                  = 22, ///< track resolution-X/Y per half chamber
143     kMCSPerHalfCh                       = 24, ///< MCS X/Y-dispersion of extrapolated track per half chamber
144     kLocalChi2PerCh                     = 26, ///< local chi2-X/Y/total distribution per chamber
145     kLocalChi2PerDE                     = 29  ///< local chi2-X/Y/total distribution per DE
146   };
147   
148   enum eResidualsVsP {
149     kResidualInChVsPClusterIn           = 0,  ///< cluster-track residual-X/Y distribution in chamber i versus momentum (cluster attached to the track)
150     kResidualInChVsPClusterOut          = 20, ///< cluster-track residual-X/Y distribution in chamber i versus momentum (cluster not attached to the track)
151     kResidualVsPClusterIn               = 40, ///< cluster-track residual-X/Y distribution integrated over chambers versus momentum (cluster attached to the track)
152     kResidualVsPClusterOut              = 42  ///< cluster-track residual-X/Y distribution integrated over chambers versus momentum (cluster not attached to the track)
153   };
154   
155   enum eResidualsVsCent {
156     kResidualInChVsCentClusterIn        = 0,  ///< cluster-track residual-X/Y distribution in chamber i versus centrality (cluster attached to the track)
157     kResidualInChVsCentClusterOut       = 20, ///< cluster-track residual-X/Y distribution in chamber i versus centrality (cluster not attached to the track)
158     kResidualVsCentClusterIn            = 40, ///< cluster-track residual-X/Y distribution integrated over chambers versus centrality (cluster attached to the track)
159     kResidualVsCentClusterOut           = 42  ///< cluster-track residual-X/Y distribution integrated over chambers versus centrality (cluster not attached to the track)
160   };
161   
162   enum eResidualsVsAngle {
163     kResidualInChVsAngleClusterIn       = 0,  ///< cluster-track residual-X/Y distribution in chamber i versus track angle in X/Y direction (cluster attached to the track)
164     kResidualInChVsAngleClusterOut      = 20, ///< cluster-track residual-X/Y distribution in chamber i versus track angle in X/Y direction (cluster not attached to the track)
165     kResidualVsAngleClusterIn           = 40, ///< cluster-track residual-X/Y distribution integrated over chambers versus track angle in X/Y direction (cluster attached to the track)
166     kResidualVsAngleClusterOut          = 42, ///< cluster-track residual-X/Y distribution integrated over chambers versus track angle in X/Y direction (cluster not attached to the track)
167     kResidualInHalfChVsAngleClusterIn   = 44  ///< cluster-track residual-X/Y distribution in half-chamber i versus track angle in X/Y direction (cluster attached to the track)
168   };
169   
170   enum eLocalChi2 {
171     kLocalChi2PerChMean                 = 0,  ///< local chi2-X/Y/total per chamber: mean
172     kLocalChi2PerDEMean                 = 3   ///< local chi2-X/Y/total per DE: mean
173   };
174   
175   enum eChamberRes {
176     kResidualPerChMeanClusterIn         = 0,  ///< cluster-track residual-X/Y per chamber: mean (cluster in)
177     kResidualPerChMeanClusterOut        = 2,  ///< cluster-track residual-X/Y per chamber: mean (cluster out)
178     kResidualPerChSigmaClusterIn        = 4,  ///< cluster-track residual-X/Y per chamber: sigma (cluster in)
179     kResidualPerChSigmaClusterOut       = 6,  ///< cluster-track residual-X/Y per chamber: sigma (cluster out)
180     kResidualPerChDispersionClusterOut  = 8,  ///< cluster-track residual-X/Y per chamber: dispersion (cluster out)
181     kCombinedResidualPerChSigma         = 10, ///< combined cluster-track residual-X/Y per chamber
182     kTrackResPerChMean                  = 12, ///< track X/Y-resolution per chamber
183     kMCSPerChMean                       = 14, ///< MCS X/Y-dispersion of extrapolated track per chamber
184     kClusterResPerCh                    = 16, ///< cluster X/Y-resolution per chamber
185     kCalcClusterResPerCh                = 18, ///< calculated cluster X/Y-resolution per chamber
186     kResidualPerDEMeanClusterIn         = 20, ///< cluster-track residual-X/Y per DE: mean (cluster in)
187     kResidualPerDEMeanClusterOut        = 22, ///< cluster-track residual-X/Y per DE: mean (cluster out)
188     kCombinedResidualPerDESigma         = 24, ///< combined cluster-track residual-X/Y per DE
189     kClusterResPerDE                    = 26, ///< cluster X/Y-resolution per DE
190     kResidualPerHalfChMeanClusterIn     = 28, ///< cluster-track residual-X/Y per half chamber: mean (cluster in)
191     kResidualPerHalfChMeanClusterOut    = 30, ///< cluster-track residual-X/Y per half chamber: mean (cluster out)
192     kCombinedResidualPerHalfChSigma     = 32, ///< combined cluster-track residual-X/Y per half chamber
193     kClusterResPerHalfCh                = 34, ///< cluster X/Y-resolution per half chamber
194     kResidualMeanClusterInVsP           = 36, ///< cluster-track residual-X/Y per chamber versus momentum: mean (cluster in)
195     kCombinedResidualSigmaVsP           = 38, ///< cluster X/Y-resolution per chamber versus momentum
196     kCombinedResidualAllChSigmaVsP      = 40, ///< cluster X/Y-resolution integrated over chambers versus momentum
197     kResidualMeanClusterInVsCent        = 42, ///< cluster-track residual-X/Y per chamber versus centrality: mean (cluster in)
198     kCombinedResidualSigmaVsCent        = 44, ///< cluster X/Y-resolution per chamber versus centrality
199     kCombinedResidualAllChSigmaVsCent   = 46, ///< cluster X/Y-resolution integrated over chambers versus centrality
200     kResidualMeanClusterInVsAngle       = 48, ///< cluster-track residual-X/Y per chamber versus track angle in X/Y direction: mean (cluster in)
201     kCombinedResidualSigmaVsAngle       = 50, ///< cluster X/Y-resolution per chamber versus track angle in X/Y direction
202     kCombinedResidualAllChSigmaVsAngle  = 52, ///< cluster X/Y-resolution integrated over chambers versus track angle in X/Y direction
203     kHChResidualMeanClusterInVsAngle    = 54  ///< cluster-track residual-X/Y per half-chamber versus track angle in X/Y direction: mean (cluster in)
204   };
205   
206   enum eTrackRes {
207     kUncorrPRes                         = 0,  ///< muon momentum reconstructed resolution at first cluster vs p
208     kPRes                               = 1,  ///< muon momentum reconstructed resolution at vertex vs p
209     kUncorrPtRes                        = 2,  ///< muon transverse momentum reconstructed resolution at first cluster vs p
210     kPtRes                              = 3,  ///< muon transverse momentum reconstructed resolution at vertex vs p
211     kUncorrSlopeRes                     = 4,  ///< muon slope-X/Y reconstructed resolution at first cluster vs p
212     kSlopeRes                           = 6   ///< muon slope-X/Y reconstructed resolution at vertex vs p
213   };
214   
215   enum eCanvases {
216     kResPerCh                           = 0,  ///< summary canvas
217     kResPerDE                           = 1,  ///< summary canvas
218     kResPerHalfCh                       = 2,  ///< summary canvas
219     kResPerChVsP                        = 3,  ///< summary canvas
220     kResPerChVsCent                     = 4,  ///< summary canvas
221     kResPerChVsAngle                    = 5,  ///< summary canvas
222     kShiftPerChVsP                      = 6,  ///< summary canvas
223     kShiftPerChVsCent                   = 7,  ///< summary canvas
224     kShiftPerChVsAngle                  = 8,  ///< summary canvas
225     kDetailResPerCh                     = 9,  ///< summary canvas
226     kDetailResPerHalfCh                 = 13, ///< summary canvas
227     kDetailResPerDE                     = 17  ///< summary canvas
228   };
229   
230   static const Int_t fgkMinEntries; ///< minimum number of entries needed to compute resolution
231   
232   TObjArray*  fResiduals;       //!< List of residual histos
233   TObjArray*  fResidualsVsP;    //!< List of residual vs. p histos
234   TObjArray*  fResidualsVsCent; //!< List of residual vs. centrality histos
235   TObjArray*  fResidualsVsAngle;//!< List of residual vs. track angle histos
236   TObjArray*  fLocalChi2;       //!< List of plots related to local chi2 per chamber/DE
237   TObjArray*  fChamberRes;      //!< List of plots related to chamber/DE resolution
238   TObjArray*  fTrackRes;        //!< List of plots related to track resolution (p, pT, ...)
239   TObjArray*  fCanvases;        //!< List of canvases summarizing the results
240   TObjArray*  fTmpHists;        //!< List of temporary histograms
241   
242   Double_t fClusterResNB[10];   ///< cluster resolution in non-bending direction
243   Double_t fClusterResB[10];    ///< cluster resolution in bending direction
244   
245   Double_t fHalfChShiftNB[20];  ///< half-chamber deplacements in non-bending direction
246   Double_t fHalfChShiftB[20];   ///< half-chamber deplacements in bending direction
247   
248   Double_t fDEShiftNB[200];     ///< DE deplacements in non-bending direction
249   Double_t fDEShiftB[200];      ///< DE deplacements in bending direction
250   
251   TString  fDefaultStorage;        ///< location of the default OCDB storage
252   Int_t    fNEvents;               //!< number of processed events
253   Bool_t   fShowProgressBar;       ///< show the progression bar
254   Bool_t   fPrintClResPerCh;       ///< print the cluster resolution per chamber
255   Bool_t   fPrintClResPerDE;       ///< print the cluster resolution per DE
256   TF1*     fGaus;                  ///< gaussian function to fit the residuals
257   Double_t fMinMomentum;           ///< use only tracks with momentum higher than this value
258   Bool_t   fSelectPhysics;         ///< use only tracks passing the physics selection
259   Bool_t   fMatchTrig;             ///< use only tracks matched with trigger
260   Bool_t   fApplyAccCut;           ///< use only tracks passing the acceptance cuts (Rabs, eta)
261   Bool_t   fSelectTrigger;         ///< use only tracks passing the trigger selection
262   UInt_t   fTriggerMask;           ///< trigger mask to be used when selecting tracks
263   Int_t    fExtrapMode;            ///< extrapolation mode to get the track parameters and covariances at a given cluster
264   Bool_t   fCorrectForSystematics; ///< add or not the systematic shifts of the residuals to the resolution
265   Bool_t   fRemoveMonoCathCl;      ///< remove or not the mono-cathod clusters
266   Bool_t   fCheckAllPads;          ///< use all pads or only the ones directly below the cluster to look for mono-cathods
267   Bool_t   fImproveTracks;         ///< flag telling whether to improve or not the track before measuring the resolution
268   Bool_t   fShiftHalfCh;           ///< flag telling wether to displace half-chambers by fHalfChShift(N)B[i] or not
269   Bool_t   fPrintHalfChShift;      ///< print the half-chamber displacements
270   Bool_t   fShiftDE;               ///< flag telling wether to displace DEs by fDEShift(N)B[i] or not
271   Bool_t   fPrintDEShift;          ///< print the DE displacements
272   Bool_t   fOCDBLoaded;            //!< flag telling if the OCDB has been properly loaded or not
273   Int_t    fNDE;                   //!< total number of DE
274   Int_t    fDEIndices[1100];       //!< index of DE in histograms refered by ID
275   Int_t    fDEIds[200];            //!< ID of DE refered by index in histograms
276   Bool_t   fReAlign;               ///< flag telling whether to re-align the spectrometer or not before computing resolution
277   TString  fOldAlignStorage;       ///< location of the OCDB storage where to find old MUON/Align/Data (use the default one if empty)
278   TString  fNewAlignStorage;       ///< location of the OCDB storage where to find new MUON/Align/Data (use the default one if empty)
279   AliMUONGeometryTransformer* fOldGeoTransformer; //!< geometry transformer used to recontruct the present data
280   AliMUONGeometryTransformer* fNewGeoTransformer; //!< new geometry transformer containing the new alignment to be applied
281   
282   TList* fSelectTriggerClass; //!< list of trigger class that can be selected to fill histograms
283   
284   AliMuonTrackCuts* fMuonTrackCuts; ///< cuts to select tracks to be considered
285   
286   ClassDef(AliAnalysisTaskMuonResolution, 4); // chamber resolution analysis
287 };
288
289 //________________________________________________________________________
290 inline void AliAnalysisTaskMuonResolution::SetStartingResolution(Int_t chId, Double_t valNB, Double_t valB)
291 {
292   /// set chamber non-bending and bending resolutions
293   if (chId < 0 || chId >= 10) return;
294   fClusterResNB[chId] = valNB;
295   fClusterResB[chId] = valB;
296 }
297
298 //________________________________________________________________________
299 inline void AliAnalysisTaskMuonResolution::SetStartingResolution(Double_t valNB[10], Double_t valB[10])
300 {
301   /// set chambers non-bending and bending resolutions
302   for (Int_t i = 0; i < 10; i++) {
303     fClusterResNB[i] = valNB[i];
304     fClusterResB[i] = valB[i];
305   }
306 }
307
308 //________________________________________________________________________
309 inline void AliAnalysisTaskMuonResolution::GetStartingResolution(Double_t valNB[10], Double_t valB[10]) const
310 {
311   /// set chambers non-bending and bending resolutions
312   for (Int_t i = 0; i < 10; i++) {
313     valNB[i] = fClusterResNB[i];
314     valB[i] = fClusterResB[i];
315   }
316 }
317
318 //________________________________________________________________________
319 inline void AliAnalysisTaskMuonResolution::SetHalfChShift(Int_t hchId, Double_t valNB, Double_t valB)
320 {
321   /// set chamber non-bending and bending resolutions
322   if (hchId < 0 || hchId >= 20) return;
323   fHalfChShiftNB[hchId] = valNB;
324   fHalfChShiftB[hchId] = valB;
325 }
326
327 //________________________________________________________________________
328 inline void AliAnalysisTaskMuonResolution::SetHalfChShift(Double_t valNB[20], Double_t valB[20])
329 {
330   /// set chambers non-bending and bending resolutions
331   for (Int_t i = 0; i < 20; i++) {
332     fHalfChShiftNB[i] = valNB[i];
333     fHalfChShiftB[i] = valB[i];
334   }
335 }
336
337 //________________________________________________________________________
338 inline void AliAnalysisTaskMuonResolution::GetHalfChShift(Double_t valNB[20], Double_t valB[20]) const
339 {
340   /// set chambers non-bending and bending resolutions
341   for (Int_t i = 0; i < 20; i++) {
342     valNB[i] = fHalfChShiftNB[i];
343     valB[i] = fHalfChShiftB[i];
344   }
345 }
346
347 //________________________________________________________________________
348 inline void AliAnalysisTaskMuonResolution::SetDEShift(Int_t iDE, Double_t valNB, Double_t valB)
349 {
350   /// set chamber non-bending and bending resolutions
351   if (iDE < 0 || iDE >= 200) return;
352   fDEShiftNB[iDE] = valNB;
353   fDEShiftB[iDE] = valB;
354 }
355
356 //________________________________________________________________________
357 inline void AliAnalysisTaskMuonResolution::SetDEShift(Double_t valNB[200], Double_t valB[200])
358 {
359   /// set chambers non-bending and bending resolutions
360   for (Int_t i = 0; i < 200; i++) {
361     fDEShiftNB[i] = valNB[i];
362     fDEShiftB[i] = valB[i];
363   }
364 }
365
366 //________________________________________________________________________
367 inline void AliAnalysisTaskMuonResolution::GetDEShift(Double_t valNB[200], Double_t valB[200]) const
368 {
369   /// set chambers non-bending and bending resolutions
370   for (Int_t i = 0; i < 200; i++) {
371     valNB[i] = fDEShiftNB[i];
372     valB[i] = fDEShiftB[i];
373   }
374 }
375
376 //________________________________________________________________________
377 inline void AliAnalysisTaskMuonResolution::ReAlign(const char* oldAlignStorage, const char* newAlignStorage)
378 {
379   /// Set the flag to activate the re-alignment and the specific storage where to find the old/new alignment data.
380   /// If oldAlignStorage = 0x0 we assume the spectrometer was not aligned before (default geometry)
381   /// If old(new)AlignStorage = "" we assume the old(new) alignment data are in the default storage
382   if (oldAlignStorage) fOldAlignStorage = oldAlignStorage;
383   else fOldAlignStorage = "none";
384   fNewAlignStorage = newAlignStorage;
385   fReAlign = kTRUE;
386 }
387
388 //________________________________________________________________________
389 inline void AliAnalysisTaskMuonResolution::FitResiduals(Bool_t flag)
390 {
391   /// set gaussian function to fit the residual distribution to extract the mean and the dispersion.
392   /// if not set: take the mean and the RMS of the distribution
393   delete fGaus;
394   if (flag) fGaus = new TF1("fGaus","gaus");
395   else fGaus = NULL;
396 }
397
398 //________________________________________________________________________
399 inline void AliAnalysisTaskMuonResolution::SetMuonTrackCuts(AliMuonTrackCuts &trackCuts)
400 {
401   /// set standard cuts to select tracks to be considered
402   delete fMuonTrackCuts;
403   fMuonTrackCuts = new AliMuonTrackCuts(trackCuts);
404 }
405
406 #endif
407