]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG/muon/AliMuonPairCuts.cxx
Updates to run with deltas (L. Cunqueiro)
[u/mrichter/AliRoot.git] / PWG / muon / AliMuonPairCuts.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
7  * Permission to use, copy, modify and distribute this software and its   *
8  * documentation strictly for non-commercial purposes is hereby granted   *
9  * without fee, provided that the above copyright notice appears in all   *
10  * copies and that both the copyright notice and this permission notice   *
11  * appear in the supporting documentation. The authors make no claims     *
12  * about the suitability of this software for any purpose. It is          *
13  * provided "as is" without express or implied warranty.                  *
14  **************************************************************************/
15
16 #include "AliMuonPairCuts.h"
17
18 #include "TMath.h"
19 #include "TList.h"
20 #include "TLorentzVector.h"
21
22 #include "AliLog.h"
23 #include "AliVParticle.h"
24
25 /// \cond CLASSIMP
26 ClassImp(AliMuonPairCuts) // Class implementation in ROOT context
27 /// \endcond
28
29 //________________________________________________________________________
30 AliMuonPairCuts::AliMuonPairCuts() :
31   AliAnalysisCuts(),
32   fMuonTrackCuts()
33 {
34   /// Default ctor.
35 }
36
37 //________________________________________________________________________
38 AliMuonPairCuts::AliMuonPairCuts(const char* name, const char* title, Bool_t isESD) :
39   AliAnalysisCuts(name, title),
40   fMuonTrackCuts("muonTrackCutInPair","muonTrackCutInPair", isESD)
41 {
42   /// Constructor
43   SetDefaultFilterMask();
44 }
45
46
47 //________________________________________________________________________
48 AliMuonPairCuts::AliMuonPairCuts(const char* name, const char* title, const AliMuonTrackCuts& trackCuts) :
49 AliAnalysisCuts(name, title),
50 fMuonTrackCuts(trackCuts)
51 {
52   /// Test Constructor
53   SetDefaultFilterMask();
54 }
55
56
57
58 //________________________________________________________________________
59 AliMuonPairCuts::AliMuonPairCuts(const AliMuonPairCuts& obj) :
60   AliAnalysisCuts(obj),
61   fMuonTrackCuts(obj.fMuonTrackCuts)
62 {
63   /// Copy constructor
64 }
65
66
67 //________________________________________________________________________
68 AliMuonPairCuts::~AliMuonPairCuts()
69 {
70   /// Destructor
71 }
72
73 //________________________________________________________________________
74 Bool_t AliMuonPairCuts::SetRun( Int_t runNumber )
75 {
76   /// Get parameters from OADB for runNumber
77   return fMuonTrackCuts.SetRun(runNumber);
78 }
79
80
81 //________________________________________________________________________
82 Bool_t AliMuonPairCuts::IsSelected( TObject* /*obj*/ )
83 {
84   /// Not implemented
85   AliError("Requires a list of two AliVParticle");
86   return kFALSE;
87 }
88
89 //________________________________________________________________________
90 Bool_t AliMuonPairCuts::IsSelected( TList* list )
91 {
92   /// Pair is selected
93   UInt_t filterMask = GetFilterMask();
94   UInt_t selectionMask = GetSelectionMask(list);
95   
96   return ( ( selectionMask & filterMask ) == filterMask );
97 }
98
99
100 //________________________________________________________________________
101 UInt_t AliMuonPairCuts::GetSelectionMask( const TObject* obj )
102 {
103   /// Get selection mask (overloaded function)
104   const TList* list = static_cast<const TList*> ( obj );
105   if ( list->GetEntries() < 2 ) {
106     AliError("Requires a list of two AliVParticle");
107     return 0;
108   }
109
110   return GetSelectionMask(list->At(0), list->At(1));
111
112 }
113
114
115 //________________________________________________________________________
116 UInt_t AliMuonPairCuts::GetSelectionMask( const TObject* track1, const TObject* track2 )
117 {
118   /// Get selection mask from AliVParticles
119   
120   UInt_t selectionMask = 0;
121     
122   UInt_t maskTrack1 = fMuonTrackCuts.GetSelectionMask(track1);
123   UInt_t maskTrack2 = fMuonTrackCuts.GetSelectionMask(track2);
124   
125   UInt_t maskAND = maskTrack1 & maskTrack2;
126   UInt_t maskOR = maskTrack1 | maskTrack2;
127   if ( maskAND & AliMuonTrackCuts::kMuEta ) selectionMask |= kBothMuEta;
128   if ( maskAND & AliMuonTrackCuts::kMuThetaAbs ) selectionMask |= kBothMuThetaAbs;
129   if ( maskAND & AliMuonTrackCuts::kMuPdca ) selectionMask |= kBothMuPdca;
130   if ( maskAND & AliMuonTrackCuts::kMuTrackChiSquare ) selectionMask |= kBothMuTrackChiSquare;
131   if ( maskAND & AliMuonTrackCuts::kMuMatchApt ) selectionMask |= kBothMuMatchApt;
132   if ( maskAND & AliMuonTrackCuts::kMuMatchLpt ) selectionMask |= kBothMuMatchLpt;
133   if ( maskAND & AliMuonTrackCuts::kMuMatchHpt ) selectionMask |= kBothMuMatchHpt;
134   if ( maskOR & AliMuonTrackCuts::kMuMatchApt ) selectionMask |= kOneMuMatchApt;
135   if ( maskOR & AliMuonTrackCuts::kMuMatchLpt ) selectionMask |= kOneMuMatchLpt;
136   if ( maskOR & AliMuonTrackCuts::kMuMatchHpt ) selectionMask |= kOneMuMatchHpt;
137   
138   TLorentzVector vec[2];
139   Double_t chargeProduct = 1.;
140   for ( Int_t itrack=0; itrack<2; ++itrack ) {
141     const AliVParticle* track = static_cast<const AliVParticle*> ( ( itrack == 0 ) ? track1 : track2 );
142     chargeProduct *= track->Charge();
143     Double_t trackP = track->P();
144     Double_t energy = TMath::Sqrt(trackP*trackP + MuonMass2());
145     vec[itrack].SetPxPyPzE(track->Px(), track->Py(), track->Pz(), energy);
146   }
147   
148   if ( chargeProduct < 0. ) selectionMask |= kDimuUnlikeSign;
149   
150   TLorentzVector vecPair = vec[0] + vec[1];
151   Double_t rapidity = vecPair.Rapidity();
152   if ( rapidity > -4. && rapidity < -2.5 ) selectionMask |= kDimuRapidity; 
153   
154   return selectionMask;
155 }
156
157
158 //_____________________________________________________________________________
159 Double_t AliMuonPairCuts::MuonMass2() const
160 {
161   /// A usefull constant
162   static Double_t m2 = 1.11636129640000012e-02; // using a constant here as the line below is a problem for CINT...
163   return m2;
164 }
165
166
167 //________________________________________________________________________
168 void AliMuonPairCuts::SetDefaultFilterMask ()
169 {
170   /// Standard cuts for muon pair
171   SetFilterMask ( kBothMuEta | kBothMuThetaAbs | kBothMuMatchLpt | kDimuUnlikeSign | kDimuRapidity );
172 }  
173
174 //________________________________________________________________________
175 void AliMuonPairCuts::SetIsMC ( Bool_t isMC )
176 {
177   /// Set Is MC
178   fMuonTrackCuts.SetIsMC(isMC);
179 }
180
181 //________________________________________________________________________
182 void AliMuonPairCuts::Print(Option_t* option) const
183 {
184   //
185   /// Print info
186   //
187   
188   TString sopt(option);
189   sopt.ToLower();
190   if ( sopt.IsNull() || sopt.Contains("*") || sopt.Contains("all") ) sopt += " pair track";
191   if ( sopt.Contains("pair") ) {
192     UInt_t filterMask = GetFilterMask();
193     printf(" *** Muon pair filter mask: *** \n");
194     printf("  0x%x\n", filterMask);
195     if ( filterMask & kBothMuEta ) printf("  mu1 && mu2 pass eta cut\n");
196     if ( filterMask & kBothMuThetaAbs ) printf("  mu1 && mu2 pass theta_abs cut\n");
197     if ( filterMask & kBothMuPdca ) printf("  mu1 && mu2 pass pxDCA cut\n");
198     if ( filterMask & kBothMuTrackChiSquare ) printf("  Chi2 cut on track\n");
199     if ( filterMask & kBothMuMatchApt ) printf("  mu1 && mu2 match Apt\n");
200     if ( filterMask & kBothMuMatchLpt ) printf("  mu1 && mu2 match Lpt\n");
201     if ( filterMask & kBothMuMatchHpt ) printf("  mu1 && mu2 match Hpt\n");
202     if ( filterMask & kBothMuMatchApt ) printf("  mu1 || mu2 match Apt\n");
203     if ( filterMask & kBothMuMatchLpt ) printf("  mu1 || mu2 match Lpt\n");
204     if ( filterMask & kBothMuMatchHpt ) printf("  mu1 || mu2 match Hpt\n");
205     if ( filterMask & kDimuUnlikeSign ) printf("  Unlike sign\n");
206     if ( filterMask & kDimuRapidity ) printf("  -4 < y_{mumu} < -2.5\n");
207     printf(" ******************** \n");
208   }
209   
210   if ( sopt.Contains("track") )
211     fMuonTrackCuts.Print(sopt.Data());
212 }