1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
16 #include "AliMuonPairCuts.h"
20 #include "TLorentzVector.h"
23 #include "AliVParticle.h"
26 ClassImp(AliMuonPairCuts) // Class implementation in ROOT context
29 //________________________________________________________________________
30 AliMuonPairCuts::AliMuonPairCuts() :
37 //________________________________________________________________________
38 AliMuonPairCuts::AliMuonPairCuts(const char* name, const char* title, Bool_t isESD) :
39 AliAnalysisCuts(name, title),
40 fMuonTrackCuts("muonTrackCutInPair","muonTrackCutInPair", isESD)
43 SetDefaultFilterMask();
47 //________________________________________________________________________
48 AliMuonPairCuts::AliMuonPairCuts(const char* name, const char* title, const AliMuonTrackCuts& trackCuts) :
49 AliAnalysisCuts(name, title),
50 fMuonTrackCuts(trackCuts)
53 SetDefaultFilterMask();
58 //________________________________________________________________________
59 AliMuonPairCuts::AliMuonPairCuts(const AliMuonPairCuts& obj) :
61 fMuonTrackCuts(obj.fMuonTrackCuts)
67 //________________________________________________________________________
68 AliMuonPairCuts& AliMuonPairCuts::operator=(const AliMuonPairCuts& obj)
70 /// Assignment operator
72 AliAnalysisCuts::operator=(obj);
73 fMuonTrackCuts = obj.fMuonTrackCuts;
79 //________________________________________________________________________
80 AliMuonPairCuts::~AliMuonPairCuts()
85 //________________________________________________________________________
86 Bool_t AliMuonPairCuts::SetRun( Int_t runNumber )
88 /// Get parameters from OADB for runNumber
89 return fMuonTrackCuts.SetRun(runNumber);
93 //________________________________________________________________________
94 Bool_t AliMuonPairCuts::IsSelected( TObject* /*obj*/ )
97 AliError("Requires a list of two AliVParticle");
101 //________________________________________________________________________
102 Bool_t AliMuonPairCuts::IsSelected( TList* list )
105 UInt_t filterMask = GetFilterMask();
106 UInt_t selectionMask = GetSelectionMask(list);
108 return ( ( selectionMask & filterMask ) == filterMask );
112 //________________________________________________________________________
113 UInt_t AliMuonPairCuts::GetSelectionMask( const TObject* obj )
115 /// Get selection mask (overloaded function)
116 const TList* list = static_cast<const TList*> ( obj );
117 if ( list->GetEntries() < 2 ) {
118 AliError("Requires a list of two AliVParticle");
122 return GetSelectionMask(list->At(0), list->At(1));
127 //________________________________________________________________________
128 UInt_t AliMuonPairCuts::GetSelectionMask( const TObject* track1, const TObject* track2 )
130 /// Get selection mask from AliVParticles
132 UInt_t selectionMask = 0;
134 UInt_t maskTrack1 = fMuonTrackCuts.GetSelectionMask(track1);
135 UInt_t maskTrack2 = fMuonTrackCuts.GetSelectionMask(track2);
137 UInt_t maskAND = maskTrack1 & maskTrack2;
138 UInt_t maskOR = maskTrack1 | maskTrack2;
139 if ( maskAND & AliMuonTrackCuts::kMuEta ) selectionMask |= kBothMuEta;
140 if ( maskAND & AliMuonTrackCuts::kMuThetaAbs ) selectionMask |= kBothMuThetaAbs;
141 if ( maskAND & AliMuonTrackCuts::kMuPdca ) selectionMask |= kBothMuPdca;
142 if ( maskAND & AliMuonTrackCuts::kMuTrackChiSquare ) selectionMask |= kBothMuTrackChiSquare;
143 if ( maskAND & AliMuonTrackCuts::kMuMatchApt ) selectionMask |= kBothMuMatchApt;
144 if ( maskAND & AliMuonTrackCuts::kMuMatchLpt ) selectionMask |= kBothMuMatchLpt;
145 if ( maskAND & AliMuonTrackCuts::kMuMatchHpt ) selectionMask |= kBothMuMatchHpt;
146 if ( maskOR & AliMuonTrackCuts::kMuMatchApt ) selectionMask |= kOneMuMatchApt;
147 if ( maskOR & AliMuonTrackCuts::kMuMatchLpt ) selectionMask |= kOneMuMatchLpt;
148 if ( maskOR & AliMuonTrackCuts::kMuMatchHpt ) selectionMask |= kOneMuMatchHpt;
150 TLorentzVector vec[2];
151 Double_t chargeProduct = 1.;
152 for ( Int_t itrack=0; itrack<2; ++itrack ) {
153 const AliVParticle* track = static_cast<const AliVParticle*> ( ( itrack == 0 ) ? track1 : track2 );
154 chargeProduct *= track->Charge();
155 Double_t trackP = track->P();
156 Double_t energy = TMath::Sqrt(trackP*trackP + MuonMass2());
157 vec[itrack].SetPxPyPzE(track->Px(), track->Py(), track->Pz(), energy);
160 if ( chargeProduct < 0. ) selectionMask |= kDimuUnlikeSign;
162 TLorentzVector vecPair = vec[0] + vec[1];
163 Double_t rapidity = vecPair.Rapidity();
164 if ( rapidity > -4. && rapidity < -2.5 ) selectionMask |= kDimuRapidity;
166 return selectionMask;
170 //_____________________________________________________________________________
171 Double_t AliMuonPairCuts::MuonMass2() const
173 /// A usefull constant
174 static Double_t m2 = 1.11636129640000012e-02; // using a constant here as the line below is a problem for CINT...
179 //________________________________________________________________________
180 void AliMuonPairCuts::SetDefaultFilterMask ()
182 /// Standard cuts for muon pair
183 SetFilterMask ( kBothMuEta | kBothMuThetaAbs | kBothMuMatchLpt | kDimuUnlikeSign | kDimuRapidity );
186 //________________________________________________________________________
187 void AliMuonPairCuts::SetIsMC ( Bool_t isMC )
190 fMuonTrackCuts.SetIsMC(isMC);
193 //________________________________________________________________________
194 void AliMuonPairCuts::Print(Option_t* option) const
200 TString sopt(option);
202 if ( sopt.IsNull() || sopt.Contains("*") || sopt.Contains("all") ) sopt += " pair track";
203 if ( sopt.Contains("pair") ) {
204 UInt_t filterMask = GetFilterMask();
205 printf(" *** Muon pair filter mask: *** \n");
206 printf(" 0x%x\n", filterMask);
207 if ( filterMask & kBothMuEta ) printf(" mu1 && mu2 pass eta cut\n");
208 if ( filterMask & kBothMuThetaAbs ) printf(" mu1 && mu2 pass theta_abs cut\n");
209 if ( filterMask & kBothMuPdca ) printf(" mu1 && mu2 pass pxDCA cut\n");
210 if ( filterMask & kBothMuTrackChiSquare ) printf(" Chi2 cut on track\n");
211 if ( filterMask & kBothMuMatchApt ) printf(" mu1 && mu2 match Apt\n");
212 if ( filterMask & kBothMuMatchLpt ) printf(" mu1 && mu2 match Lpt\n");
213 if ( filterMask & kBothMuMatchHpt ) printf(" mu1 && mu2 match Hpt\n");
214 if ( filterMask & kBothMuMatchApt ) printf(" mu1 || mu2 match Apt\n");
215 if ( filterMask & kBothMuMatchLpt ) printf(" mu1 || mu2 match Lpt\n");
216 if ( filterMask & kBothMuMatchHpt ) printf(" mu1 || mu2 match Hpt\n");
217 if ( filterMask & kDimuUnlikeSign ) printf(" Unlike sign\n");
218 if ( filterMask & kDimuRapidity ) printf(" -4 < y_{mumu} < -2.5\n");
219 printf(" ******************** \n");
222 if ( sopt.Contains("track") )
223 fMuonTrackCuts.Print(sopt.Data());