]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG/muon/AliMuonTrackCuts.cxx
1) Adding class AliAnalysisMuonUtility which contains static methods allowing to...
[u/mrichter/AliRoot.git] / PWG / muon / AliMuonTrackCuts.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 "AliMuonTrackCuts.h"
17
18 #include "TMath.h"
19 #include "TList.h"
20 #include "TArrayD.h"
21 #include "TArrayI.h"
22 #include "TObjArray.h"
23 #include "TObjString.h"
24 #include "TFile.h"
25 #include "TParameter.h"
26 #include "TKey.h"
27 #include "TVector3.h"
28 #include "TSystem.h"
29
30 #include "AliLog.h"
31 #include "AliVParticle.h"
32 #include "AliESDMuonTrack.h"
33 #include "AliAODTrack.h"
34 #include "AliAnalysisManager.h"
35 #include "AliAnalysisMuonUtility.h"
36
37 /// \cond CLASSIMP
38 ClassImp(AliMuonTrackCuts) // Class implementation in ROOT context
39 /// \endcond
40
41
42 //________________________________________________________________________
43 AliMuonTrackCuts::AliMuonTrackCuts() :
44   AliAnalysisCuts(),
45   fIsMC(kFALSE),
46   fUseCustomParam(kFALSE),
47   fSharpPtCut(kFALSE),
48   fParameters(TArrayD(kNParameters))
49 {
50   /// Default ctor.
51   fParameters.Reset();
52 }
53
54 //________________________________________________________________________
55 AliMuonTrackCuts::AliMuonTrackCuts(const char* name, const char* title ) :
56 AliAnalysisCuts(name, title),
57   fIsMC(kFALSE),
58   fUseCustomParam(kFALSE),
59   fSharpPtCut(kFALSE),
60   fParameters(TArrayD(kNParameters))
61 {
62   /// Constructor
63   fParameters.Reset();
64   SetDefaultFilterMask();
65 }
66
67
68 //________________________________________________________________________
69 AliMuonTrackCuts::AliMuonTrackCuts(const AliMuonTrackCuts& obj) :
70   AliAnalysisCuts(obj),
71   fIsMC(obj.fIsMC),
72   fUseCustomParam(obj.fUseCustomParam),
73   fSharpPtCut(obj.fSharpPtCut),
74   fParameters(obj.fParameters)
75 {
76   /// Copy constructor
77 }
78
79
80 //________________________________________________________________________
81 AliMuonTrackCuts& AliMuonTrackCuts::operator=(const AliMuonTrackCuts& obj)
82 {
83   /// Assignment operator
84   if ( this != &obj ) { 
85     AliAnalysisCuts::operator=(obj);
86     fIsMC = obj.fIsMC;
87     fUseCustomParam = obj.fUseCustomParam;
88     fSharpPtCut = obj.fSharpPtCut;
89     fParameters = obj.fParameters;
90   }
91   return *this;
92 }
93
94
95 //________________________________________________________________________
96 AliMuonTrackCuts::~AliMuonTrackCuts()
97 {
98   /// Destructor
99 }
100
101 //________________________________________________________________________
102 Bool_t AliMuonTrackCuts::RunMatchesRange( Int_t runNumber, const Char_t* objName ) const
103 {
104   /// Check if the object contains the run
105   TString sname(objName);
106   TObjArray* array = sname.Tokenize("_");
107   array->SetOwner();
108   Int_t runRange[2] = { -1, -1 };
109   if ( array->GetEntries() >= 3 ) {
110     for ( Int_t irun=0; irun<2; ++irun ) {
111       TString currRun = array->At(irun+1)->GetName();
112       if ( currRun.IsDigit() ) runRange[irun] = currRun.Atoi();
113     }
114   }
115   delete array;
116   return ( runNumber >= runRange[0] && runNumber <= runRange[1]);
117 }
118
119 //________________________________________________________________________
120 void AliMuonTrackCuts::SetUseCustomParam( Bool_t useCustomParam, Int_t runNumber  )
121 {
122   /// Flag to select custom parameters
123   /// It first searches the default parameters in OADB
124   /// then disables the access to the OADB
125   /// and allows to manually modify parameters
126
127   if ( ! fUseCustomParam && useCustomParam ) SetRun(runNumber);
128   fUseCustomParam = useCustomParam;
129 }
130
131 //________________________________________________________________________
132 Bool_t AliMuonTrackCuts::SetRun( Int_t runNumber )
133 {
134   /// Get parameters from OADB for runNumber
135   
136   if ( fUseCustomParam ) return kFALSE;
137   return StreamParameters(runNumber, -1);
138 }
139
140
141 //________________________________________________________________________
142 Bool_t AliMuonTrackCuts::StreamParameters( Int_t runNumber,  Int_t runMax )
143 {
144   if ( runMax > 0 ) { // Stream to OADB
145     if ( ! fUseCustomParam ) {
146       AliError("Users are not allowed to update OADB. Use SetUseCustomParam() instead");
147       return kFALSE;
148     }
149   }
150
151   TString filename = Form("%s/PWG/MUON/MuonTrackCuts.root",AliAnalysisManager::GetOADBPath());
152   if ( fIsMC ) filename.ReplaceAll(".root", "_MC.root");
153
154   TString parNames[kNParameters];
155   parNames[kMeanDcaX]       = "MeanDcaX";
156   parNames[kMeanDcaY]       = "MeanDcaY";
157   parNames[kMeanDcaZ]       = "MeanDcaZ";
158   parNames[kMeanPCorr23]    = "MeanPCorr23";
159   parNames[kMeanPCorr310]   = "MeanPCorr310";
160   parNames[kSigmaPdca23]    = "SigmaPdca23";
161   parNames[kSigmaPdca310]   = "SigmaPdca310";
162   parNames[kNSigmaPdcaCut]  = "NSigmaPdcaCut";
163   parNames[kChi2NormCut]    = "Chi2NormCut";
164   parNames[kRelPResolution] = "RelPResolution";
165   parNames[kSharpPtApt]     = "SharpPtApt";
166   parNames[kSharpPtLpt]     = "SharpPtLpt";
167   parNames[kSharpPtHpt]     = "SharpPtHpt";
168
169   TObjArray* paramList = 0x0;
170
171   if ( runMax < 0 ) { // Get from OADB
172     TFile* file = TFile::Open(filename.Data(), "READ");
173     if ( ! file ) {
174       AliError(Form("OADB file %s not found!", filename.Data()));
175       return kFALSE;
176     }
177
178     TList* listOfKeys = file->GetListOfKeys();
179     TIter next(listOfKeys);
180     TObject* key = 0x0;
181     Bool_t foundMatch = kFALSE;
182     TObject* defaultObj = 0x0;
183     while ( ( key = next() ) ) {
184       TString objName = key->GetName();
185       objName.ToUpper();
186       if ( RunMatchesRange(runNumber, objName.Data()) ) {
187         paramList = static_cast<TObjArray*>(file->Get(key->GetName()));
188         foundMatch = kTRUE;
189         break;
190       }
191       if ( objName.Contains("DEFAULT") ) defaultObj = file->Get(key->GetName());
192     }
193
194     if ( ! foundMatch ) {
195       AliWarning("Run number not found in OADB: using default");
196       if ( defaultObj ) paramList = static_cast<TObjArray*>(defaultObj);
197       else {
198         file->Close();
199         AliError("Default parameters not found in OADB!");
200         return kFALSE;
201       }
202     }
203
204     AliInfo(Form("Required run %i. Param. set: %s", runNumber, paramList->GetName()));
205
206     for ( Int_t ipar=0; ipar<kNParameters; ++ipar ) {
207       TParameter<Double_t>* param = static_cast<TParameter<Double_t>*>(paramList->FindObject(parNames[ipar].Data()));
208       if ( ! param ) {
209         AliWarning(Form("Parameter %s not set", parNames[ipar].Data()));
210         continue;
211       }
212       fParameters[ipar] = param->GetVal();
213     }
214
215     file->Close();
216   }
217   else {
218     if ( ! paramList ) {
219       paramList = new TObjArray(kNParameters);
220       paramList->SetOwner();
221     }
222     for ( Int_t ipar=0; ipar<kNParameters; ++ipar ) {
223       TParameter<Double_t>* param= new TParameter<Double_t>(parNames[ipar].Data(), fParameters[ipar]);
224       paramList->AddAt(param, ipar);
225     }
226
227     TString paramListName = "MuonCuts_";
228     paramListName += ( runNumber < 0 ) ? "Default" : Form("%i_%i",runNumber, runMax);
229     AliInfo(Form("Adding %s to file %s", paramListName.Data(), filename.Data()));
230     paramList->SetName(paramListName.Data());
231     TFile* file = TFile::Open(filename.Data(), "UPDATE");
232     paramList->Write(paramListName.Data(), TObject::kSingleKey);
233     file->Close();
234     delete paramList;
235   }
236   return kTRUE;
237 }
238
239 //________________________________________________________________________
240 Bool_t AliMuonTrackCuts::SetParameter(Int_t iparam, Double_t value)
241 {
242   /// Set parameter
243   if ( fUseCustomParam ) {
244     fParameters.SetAt(value, iparam);
245     return kTRUE;
246   }
247
248   AliWarning("Parameters automatically taken from OADB. If you want to use with custom parameters, use SetUseCustomParam()");
249   return kFALSE;
250 }
251
252
253 //________________________________________________________________________
254 Bool_t AliMuonTrackCuts::IsSelected( TObject* obj )
255 {
256   /// Track is selected
257   UInt_t filterMask = GetFilterMask();
258   UInt_t selectionMask = GetSelectionMask(obj);
259   
260   return ( ( selectionMask & filterMask ) == filterMask );
261 }
262
263
264 //________________________________________________________________________
265 UInt_t AliMuonTrackCuts::GetSelectionMask( const TObject* obj )
266 {
267   /// Get selection mask
268   
269   const AliVParticle* track = static_cast<const AliVParticle*> ( obj );
270   
271   UInt_t selectionMask = 0;
272
273   if ( ! AliAnalysisMuonUtility::IsMuonTrack(track) ) return selectionMask;
274
275   Double_t eta = track->Eta();
276   if ( eta > -4. && eta < -2.5 ) selectionMask |= kMuEta;
277
278   Double_t thetaAbsEndDeg = AliAnalysisMuonUtility::GetThetaAbsDeg(track);
279
280   if ( thetaAbsEndDeg > 2. && thetaAbsEndDeg < 10. ) selectionMask |= kMuThetaAbs;
281
282   Int_t matchTrig = AliAnalysisMuonUtility::GetMatchTrigger(track);
283   Int_t cutLevel[3] = {kMuMatchApt, kMuMatchLpt, kMuMatchHpt};
284   Double_t pt = track->Pt();
285   for ( Int_t ilevel=0; ilevel<3; ilevel++ ) {
286     if ( matchTrig < ilevel+1 ) break;
287     if ( fSharpPtCut && pt < GetSharpPtCut(ilevel) ) break;
288     selectionMask |= cutLevel[ilevel];
289   }
290
291   if ( AliAnalysisMuonUtility::GetChi2perNDFtracker(track) < GetChi2NormCut() ) selectionMask |= kMuTrackChiSquare;
292
293   TVector3 dcaAtVz = GetCorrectedDCA(track);
294   Double_t pTotMean = GetAverageMomentum(track);
295
296   Double_t pDca = pTotMean * dcaAtVz.Mag();
297     
298   Double_t pTot = track->P();
299   Double_t rAbsEnd = AliAnalysisMuonUtility::GetRabs(track);
300   
301   // Momentum resolution only
302   // The cut depends on the momentum resolution. In particular:
303   // Sigma_pDCA = Sqrt( Sigma_pDCA_measured^2 + Sigma_p/p * pDCA)
304   // The relative momentum distribution Sigma_p/p is estimated from the track Delta_p,
305   // and from the error on sagitta Delta_s:
306   // Delta_p/p = p*Delta_s/(1+p*Delta_s)
307   // A cut at N sigmas requres a cut in N Delta_s, so
308   // Sigma_p/p(N) = p*N*Delta_s/(1+p*N*Delta_s)
309   //Double_t pResolutionEffect = pDca * pTot * GetRelPResolution() / ( 1. + GetNSigmaPdca() * GetRelPResolution()*pTot );
310   
311   //Double_t pResolutionEffect = 0.4 * pTot;  // Values used in 2010 data
312   //Double_t pResolutionEffect = 0.32 * pTot; // Values in 2011
313   //Double_t sigmaPdca = GetSigmaPdca(rAbsEnd);
314   //Double_t sigmaPdcaWithRes = TMath::Sqrt( sigmaPdca*sigmaPdca + pResolutionEffect*pResolutionEffect );
315   
316   // Momentum resolution and slope resolution 
317   // Due to the momentum resolution, the measured momentum is biased
318   // Since we want to keep as much signal as possible, we want to avoid
319   // that a measured pxDCA is rejected since the momentum is overestimated
320   // p_true = p_meas - Delta_p
321   // p_true = p_meas - N*Delta_s*p_meas / (1+n*Delta_s*p_meas)
322   // Hence:
323   // p_true x DCA < N * Sigma_pDCA_meas =>
324   // p_meas x DCA < N * Sigma_pDCA_meas / ( 1 - N*Delta_s*p_meas / (1+n*Delta_s*p_meas))
325   // Finally the cut value has to be summed in quadrature with the error on DCA,
326   // which is given by the slope resolution
327   // p_meas x DCA < N * Sqrt( ( Sigma_pDCA_meas / ( 1 - N*Delta_s*p_meas / (1+n*Delta_s*p_meas)) )^2 + (distance * sigma_slope * p_meas )^2)
328   Double_t nrp = GetNSigmaPdca() * GetRelPResolution() * pTot;
329   Double_t pResolutionEffect = GetSigmaPdca(rAbsEnd) / ( 1. - nrp / ( 1. + nrp ) );
330   Double_t slopeResolutionEffect = 535. * GetSlopeResolution() * pTot;
331   
332   Double_t sigmaPdcaWithRes = TMath::Sqrt( pResolutionEffect*pResolutionEffect + slopeResolutionEffect*slopeResolutionEffect );
333   
334   if ( pDca < GetNSigmaPdca() * sigmaPdcaWithRes ) selectionMask |= kMuPdca;
335   
336   AliDebug(1, Form("Selection mask 0x%x\n", selectionMask));
337
338   return selectionMask;
339 }
340
341
342 //________________________________________________________________________
343 Bool_t AliMuonTrackCuts::IsSelected( TList* /* list */)
344 {
345   /// Not implemented
346   AliError("Function not implemented: Use IsSelected(TObject*)");
347   return kFALSE;
348 }
349
350
351 //________________________________________________________________________
352 Int_t AliMuonTrackCuts::GetThetaAbsBin ( Double_t rAtAbsEnd ) const
353 {
354   /// Get theta abs bin
355   Double_t thetaAbsEndDeg = TMath::ATan( rAtAbsEnd / 505. ) * TMath::RadToDeg();
356   Int_t thetaAbsBin = ( thetaAbsEndDeg < 3. ) ? kThetaAbs23 : kThetaAbs310;
357   return thetaAbsBin;
358 }
359
360
361 //________________________________________________________________________
362 TVector3 AliMuonTrackCuts::GetCorrectedDCA ( const AliVParticle* track ) const
363 {
364   /// Get corrected DCA
365   
366   TVector3 vertex(AliAnalysisMuonUtility::GetXatVertex(track), AliAnalysisMuonUtility::GetYatVertex(track), AliAnalysisMuonUtility::GetZatVertex(track));
367   
368   TVector3 dcaTrack(AliAnalysisMuonUtility::GetXatDCA(track), AliAnalysisMuonUtility::GetYatDCA(track), AliAnalysisMuonUtility::GetZatDCA(track));
369   
370   TVector3 dcaAtVz = dcaTrack - vertex - GetMeanDCA();
371
372   return dcaAtVz;
373 }
374
375 //________________________________________________________________________
376 Double_t AliMuonTrackCuts::GetAverageMomentum ( const AliVParticle* track ) const
377 {
378   /// Get average momentum before and after the absorber
379   
380   Double_t pTotMean = 0.;
381   Double_t pTot = track->P();
382   //if ( isESD ) pTotMean = 0.5 * ( pTot + ((AliESDMuonTrack*)track)->PUncorrected() );
383   if ( ! AliAnalysisMuonUtility::IsAODTrack(track) ) pTotMean = ((AliESDMuonTrack*)track)->PUncorrected(); // Increased stability if using uncorrected value
384   else {
385     pTotMean = pTot - GetMeanPCorr(((AliAODTrack*)track)->GetRAtAbsorberEnd());
386   }
387
388   return pTotMean;
389 }
390
391
392 //________________________________________________________________________
393 void AliMuonTrackCuts::SetMeanDCA ( Double_t xAtDca, Double_t yAtDca, Double_t zAtDca )
394 {
395     /// Set mean DCA from track
396   SetParameter(kMeanDcaX, xAtDca);
397   SetParameter(kMeanDcaY, yAtDca);
398   SetParameter(kMeanDcaZ, zAtDca);
399 }
400
401 //________________________________________________________________________
402 TVector3 AliMuonTrackCuts::GetMeanDCA () const
403
404     /// Get mean DCA from track
405   return TVector3(fParameters[kMeanDcaX], fParameters[kMeanDcaY], fParameters[kMeanDcaZ]);
406 }
407
408 //________________________________________________________________________
409 void AliMuonTrackCuts::SetMeanPCorr ( Double_t pCorrThetaAbs23, Double_t pCorrThetaAbs310 )
410 {
411   /// Set mean p correction
412   SetParameter(kMeanPCorr23, pCorrThetaAbs23);
413   SetParameter(kMeanPCorr310, pCorrThetaAbs310);
414 }
415
416 //________________________________________________________________________
417 Double_t AliMuonTrackCuts::GetMeanPCorr ( Double_t rAtAbsEnd ) const
418 {
419   /// Get mean p correction
420   return fParameters[kMeanPCorr23+GetThetaAbsBin(rAtAbsEnd)];
421 }
422
423 //________________________________________________________________________
424 void AliMuonTrackCuts::SetSigmaPdca ( Double_t sigmaThetaAbs23, Double_t sigmaThetaAbs310 )
425
426   /// Set sigma pdca
427   SetParameter(kSigmaPdca23, sigmaThetaAbs23);
428   SetParameter(kSigmaPdca310, sigmaThetaAbs310);
429 }
430
431 //________________________________________________________________________
432 Double_t AliMuonTrackCuts::GetSigmaPdca ( Double_t rAtAbsEnd ) const
433
434   /// Get mean pdca
435   return fParameters[kSigmaPdca23+GetThetaAbsBin(rAtAbsEnd)];
436 }
437
438 //________________________________________________________________________
439 void AliMuonTrackCuts::SetNSigmaPdca ( Double_t nSigmas )
440
441   /// Set N sigma pdca cut
442   SetParameter(kNSigmaPdcaCut, nSigmas);
443 }
444
445 //________________________________________________________________________
446 Double_t AliMuonTrackCuts::GetNSigmaPdca () const
447 {
448   /// Get N sigma pdca cut
449   return fParameters[kNSigmaPdcaCut];
450 }
451
452 //________________________________________________________________________
453 void AliMuonTrackCuts::SetChi2NormCut ( Double_t chi2normCut )
454 {
455   /// Set cut on normalized chi2 of tracks
456   SetParameter(kChi2NormCut,chi2normCut);
457 }
458
459 //________________________________________________________________________
460 Double_t AliMuonTrackCuts::GetChi2NormCut () const
461 {
462   /// Get cut on normalized chi2 of tracks
463   return fParameters[kChi2NormCut];
464 }
465
466 //________________________________________________________________________
467 void AliMuonTrackCuts::SetRelPResolution ( Double_t relPResolution )
468 {
469   /// Set relative momentum resolution
470   SetParameter(kRelPResolution,relPResolution);
471 }
472
473 //________________________________________________________________________
474 Double_t AliMuonTrackCuts::GetRelPResolution () const
475 {
476   /// Get relative momentum resolution
477   return fParameters[kRelPResolution];
478 }
479
480
481 //________________________________________________________________________
482 void AliMuonTrackCuts::SetSlopeResolution ( Double_t slopeResolution )
483 {
484   /// Set slope resolution
485   SetParameter(kSlopeResolution,slopeResolution);
486 }
487
488 //________________________________________________________________________
489 Double_t AliMuonTrackCuts::GetSlopeResolution () const
490 {
491   /// Get slope resolution
492   return fParameters[kSlopeResolution];
493 }
494
495 //________________________________________________________________________
496 void AliMuonTrackCuts::SetSharpPtCut ( Double_t valueApt, Double_t valueLpt, Double_t valueHpt  )
497 {
498   /// Set sharp tracker cut matching the trigger level
499
500   SetParameter(kSharpPtApt, valueApt);
501   SetParameter(kSharpPtLpt, valueLpt);
502   SetParameter(kSharpPtHpt, valueHpt);
503 }
504
505 //________________________________________________________________________
506 Double_t AliMuonTrackCuts::GetSharpPtCut ( Int_t trigPtCut, Bool_t warn ) const
507 {
508   /// Get sharp tracker cut matching the trigger level
509   /// trigPtCut can be 0 (Apt), 1 (Lpt) or 2 (Hpt)
510   if ( trigPtCut < 0 || trigPtCut > 2 ) {
511     if ( warn ) AliError("Allowed values for trigPtCut are 0 (Apt), 1 (Lpt), 2 (Hpt)");
512     return 0.;
513   }
514   Int_t ipar = kSharpPtApt + trigPtCut;
515   return fParameters[ipar];
516 }
517
518 //________________________________________________________________________
519 void AliMuonTrackCuts::SetDefaultFilterMask ()
520 {
521   /// Standard cuts for single muon
522   SetFilterMask ( kMuEta | kMuThetaAbs | kMuPdca | kMuMatchApt );
523 }
524
525 //________________________________________________________________________
526 Bool_t AliMuonTrackCuts::TrackPtCutMatchTrigClass ( const AliVParticle* track, const TArrayI ptCutFromClass) const
527 {
528   /// Check if track passes the trigger pt cut level used in the trigger class
529   Int_t matchTrig = AliAnalysisMuonUtility::GetMatchTrigger(track);
530   Bool_t matchTrackerPt = kTRUE;
531   if ( IsApplySharpPtCutInMatching() ) {
532     matchTrackerPt = ( track->Pt() >= GetSharpPtCut(ptCutFromClass[0]-1,kFALSE) );
533   }
534   Bool_t passCut = ( ( matchTrig >= ptCutFromClass[0] ) && matchTrackerPt );
535   AliDebug(1,Form("Class matchTrig %i %i  trackMatchTrig %i  trackPt %g (required %i)  passCut %i", ptCutFromClass[0], ptCutFromClass[1], matchTrig, track->Pt(), IsApplySharpPtCutInMatching(),passCut));
536   return passCut;
537 }
538
539
540 //________________________________________________________________________
541 void AliMuonTrackCuts::Print(Option_t* option) const
542 {
543   //
544   /// Print info
545   //
546   TString sopt(option);
547   sopt.ToLower();
548   if ( sopt.IsNull() || sopt.Contains("*") || sopt.Contains("all") ) sopt = "mask param";
549   UInt_t filterMask = GetFilterMask();
550   Int_t cutLevel[3] = {kMuMatchApt, kMuMatchLpt, kMuMatchHpt};
551   TString cutLevelName[3] = {"Apt", "Lpt", "Hpt"};
552   if ( sopt.Contains("mask") ) {
553     printf(" *** Muon track filter mask: *** \n");
554     printf("  0x%x\n", filterMask);
555     if ( filterMask & kMuEta ) printf("  -4 < eta < -2.5\n");
556     if ( filterMask & kMuThetaAbs ) printf("  2 < theta_abs < 10 deg\n");
557     if ( filterMask & kMuPdca ) printf("  pxDCA cut\n");
558     for ( Int_t ilevel=0; ilevel<3; ilevel++ ) {
559       if ( filterMask & cutLevel[ilevel] ) {
560         printf("  match %s", cutLevelName[ilevel].Data());
561         if ( fSharpPtCut ) printf(" && sharp pt from tracker");
562         printf("\n");
563       }
564     }
565     if ( filterMask & kMuTrackChiSquare ) printf("  Chi2 cut on track\n");
566     printf(" ******************** \n");
567   }
568   if ( sopt.Contains("param") ) {
569     printf(" *** Muon track parameter summary: ***\n");
570     printf("  Mean vertex DCA: (%g, %g, %g)\n", fParameters[kMeanDcaX], fParameters[kMeanDcaY], fParameters[kMeanDcaZ]);
571     printf("  Mean p correction (GeV/c): theta2-3 = %g  theta3-10 = %g\n", fParameters[kMeanPCorr23], fParameters[kMeanPCorr310]);
572     printf("  Sigma p x DCA (cm x GeV/c): theta2-3 = %g  theta3-10 = %g\n", fParameters[kSigmaPdca23], fParameters[kSigmaPdca310]);
573     printf("  Cut p x DCA in units of sigma: %g\n", fParameters[kNSigmaPdcaCut]);
574     printf("  Cut on track chi2/NDF: %g\n", fParameters[kChi2NormCut]);
575     printf("  Momentum resolution: %g\n", fParameters[kRelPResolution]);
576     printf("  Slope resolution: %g\n", fParameters[kSlopeResolution]);
577     printf("  Sharp pt cut: %g (Apt)  %g (Lpt)  %g (Hpt)\n", fParameters[kSharpPtApt], fParameters[kSharpPtLpt], fParameters[kSharpPtHpt]);
578     printf(" ********************************\n");
579   }
580 }