]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
From Francesco Noferini: important change in the STEER and ANALYSIS dir that allows...
authoragheata <agheata@f7af4fe6-9843-0410-8265-dc069ae4e863>
Fri, 6 Jul 2012 13:59:16 +0000 (13:59 +0000)
committeragheata <agheata@f7af4fe6-9843-0410-8265-dc069ae4e863>
Fri, 6 Jul 2012 13:59:16 +0000 (13:59 +0000)
16 files changed:
ANALYSIS/AliAnalysisTaskPIDCombined.cxx
ANALYSIS/AliAnalysisTaskPIDCombined.h
ANALYSIS/AliAnalysisTaskPIDResponse.cxx
ANALYSIS/AliAnalysisTaskPIDResponse.h
STEER/AOD/AliAODTrack.cxx
STEER/AOD/AliAODTrack.h
STEER/AOD/AliAODpidUtil.cxx
STEER/AOD/AliAODpidUtil.h
STEER/ESD/AliESDpid.cxx
STEER/ESD/AliESDpid.h
STEER/ESD/AliESDtrack.cxx
STEER/ESD/AliESDtrack.h
STEER/STEERBase/AliPIDResponse.cxx
STEER/STEERBase/AliPIDResponse.h
STEER/STEERBase/AliTPCPIDResponse.cxx
STEER/STEERBase/AliVTrack.h

index fe85f36c3e893bc22c8324e6cfbd18018382fc07..465b1cec426fe9d785e73ed2fd0be2d70ce0408f 100644 (file)
@@ -64,7 +64,9 @@ AliAnalysisTaskPIDCombined::AliAnalysisTaskPIDCombined() :
   fPIDResponse(0x0),\r
   fPIDCombined(0x0),\r
   fTrackCuts(0x0),\r
-  fTrackFilter(0x0)\r
+  fTrackFilter(0x0),\r
+  fDeDx(NULL),
+  fDeDxTuned(NULL)
 {\r
   //\r
   // Constructor\r
@@ -89,7 +91,9 @@ AliAnalysisTaskPIDCombined::AliAnalysisTaskPIDCombined(const char *name) :
   fPIDResponse(0x0),\r
   fPIDCombined(0x0),\r
   fTrackCuts(0x0),\r
-  fTrackFilter(0x0)\r
+  fTrackFilter(0x0),\r
+  fDeDx(NULL),
+  fDeDxTuned(NULL)
 {\r
   //\r
   // Constructor\r
@@ -221,7 +225,11 @@ void AliAnalysisTaskPIDCombined::UserCreateOutputObjects()
   }\r
 \r
 \r
-\r
+  fDeDx = new TH2D("hDeDx",";p_{TPC};dE/dx (a.u.)",500,0,5,500,0,500);
+  fHistList.Add(fDeDx);
+  fDeDxTuned = new TH2D("hDeDxTuned",";p_{TPC};dE/dx (a.u.)",500,0,5,500,0,500);
+  fHistList.Add(fDeDxTuned);
+
   fHistList.SetOwner();\r
   PostData(1,&fHistList);\r
 \r
@@ -301,6 +309,10 @@ void AliAnalysisTaskPIDCombined::UserExec(Option_t *)
 \r
       }\r
 \r
+      fPIDResponse->GetTPCsignalTunedOnData(track);
+
+      fDeDx->Fill(mom,track->GetTPCsignal());
+      fDeDxTuned->Fill(mom,track->GetTPCsignalTunedOnData());
 \r
     }\r
   }\r
index e4fd958731d1422c169cd63d893037671c2bea3d..8b930359090835f0600906d06cf66715f309ee77 100644 (file)
@@ -61,6 +61,10 @@ private:
   AliPIDCombined       *fPIDCombined;     //! combined PID object\r
   AliESDtrackCuts      *fTrackCuts;            //! track selection\r
   AliAnalysisFilter    *fTrackFilter;         //! track filter\r
+
+  TH2D *fDeDx;                              //! histo with the dedx
+  TH2D *fDeDxTuned;                         //! histo to check the dedx tuning in MC
+
 \r
   AliAnalysisTaskPIDCombined(const AliAnalysisTaskPIDCombined &c);\r
   AliAnalysisTaskPIDCombined& operator= (const AliAnalysisTaskPIDCombined &c);\r
index 4e90a9e37164e7e84a4b20d6d5c8eb4bc3474c1f..853c67caa7804889beb68cbd2544bf06be272656 100644 (file)
@@ -38,7 +38,9 @@ fOADBPath(),
 fPIDResponse(0x0),
 fRun(0),
 fOldRun(0),
-fRecoPass(0)
+fRecoPass(0),
+fIsTunedOnData(kFALSE),
+fRecoPassTuned(0)  
 {
   //
   // Dummy constructor
@@ -53,7 +55,9 @@ fOADBPath(),
 fPIDResponse(0x0),
 fRun(0),
 fOldRun(0),
-fRecoPass(0)
+fRecoPass(0),
+fIsTunedOnData(kFALSE),
+fRecoPassTuned(0)
 {
   //
   // Default constructor
@@ -90,6 +94,8 @@ void AliAnalysisTaskPIDResponse::UserCreateOutputObjects()
 
   fPIDResponse->SetOADBPath(AliAnalysisManager::GetOADBPath());
   if (!fOADBPath.IsNull()) fPIDResponse->SetOADBPath(fOADBPath.Data());
+
+  if(fIsTunedOnData) fPIDResponse->SetTunedOnData(kTRUE,fRecoPassTuned);
 }
 
 //______________________________________________________________________________
index 65c6b27ff996529ac5ca704a859f966ba52b71c6..2894ae7a551e5d1370f7df5453c6dc7827226109 100644 (file)
@@ -40,6 +40,7 @@ public:
 
   void SetOADBPath(const char* path) {fOADBPath=path;}
   const char* GetOADBPath() const { return fOADBPath.Data(); }
+  void SetTuneOnData(Bool_t flag,Int_t recopass){fIsTunedOnData=flag;fRecoPassTuned=recopass;};
 
 private:
   Bool_t fIsMC;                        //  If we run on MC data
@@ -49,6 +50,9 @@ private:
   Int_t   fRun;                        //! current run number
   Int_t   fOldRun;                     //! current run number
   Int_t   fRecoPass;                   //! reconstruction pass
+
+  Bool_t  fIsTunedOnData;              // flag to tune MC on data (dE/dx)
+  Int_t   fRecoPassTuned;              // Reco pass tuned on data for MC
   
   //
   void SetRecoInfo();
index 0362a984b025722d26c5183f147b0498142831b1..72bccf988ee4f9a5607ec33213dd60f9e606de38 100644 (file)
@@ -55,6 +55,7 @@ AliAODTrack::AliAODTrack() :
   fProdVertex(NULL),
   fTrackPhiOnEMCal(-999),
   fTrackEtaOnEMCal(-999),
+  fTPCsignalTuned(0),
   fAODEvent(NULL)
 {
   // default constructor
@@ -106,6 +107,7 @@ AliAODTrack::AliAODTrack(Short_t id,
   fProdVertex(prodVertex),
   fTrackPhiOnEMCal(-999),
   fTrackEtaOnEMCal(-999),
+  fTPCsignalTuned(0),
   fAODEvent(NULL)
 {
   // constructor
@@ -161,6 +163,7 @@ AliAODTrack::AliAODTrack(Short_t id,
   fProdVertex(prodVertex),
   fTrackPhiOnEMCal(-999),
   fTrackEtaOnEMCal(-999),
+  fTPCsignalTuned(0),
   fAODEvent(NULL)
 {
   // constructor
@@ -210,6 +213,7 @@ AliAODTrack::AliAODTrack(const AliAODTrack& trk) :
   fProdVertex(trk.fProdVertex),
   fTrackPhiOnEMCal(trk.fTrackPhiOnEMCal),
   fTrackEtaOnEMCal(trk.fTrackEtaOnEMCal),
+  fTPCsignalTuned(trk.fTPCsignalTuned),
   fAODEvent(trk.fAODEvent)
 {
   // Copy constructor
@@ -257,6 +261,7 @@ AliAODTrack& AliAODTrack::operator=(const AliAODTrack& trk)
     fCaloIndex         = trk.fCaloIndex;
     fTrackPhiOnEMCal   = trk.fTrackPhiOnEMCal;
     fTrackEtaOnEMCal   = trk.fTrackEtaOnEMCal;
+    fTPCsignalTuned    = trk.fTPCsignalTuned;
 
     delete fCovMatrix;
     if(trk.fCovMatrix) fCovMatrix=new AliAODRedCov<6>(*trk.fCovMatrix);
index 06f6d74bfe086278ef94677eace144d050254769..62dc4c089bb8c0417d9cdc758307fbd076c80c25 100644 (file)
@@ -271,6 +271,8 @@ class AliAODTrack : public AliVTrack {
   //pid signal interface
   Double_t  GetITSsignal()       const { return fDetPid?fDetPid->GetITSsignal():0.;    }
   Double_t  GetTPCsignal()       const { return fDetPid?fDetPid->GetTPCsignal():0.;    }
+  Double_t  GetTPCsignalTunedOnData() const { return fTPCsignalTuned;}
+  void      SetTPCsignalTunedOnData(Double_t signal) {fTPCsignalTuned = signal;}
   UShort_t  GetTPCsignalN()      const { return fDetPid?fDetPid->GetTPCsignalN():0;    }
   virtual AliTPCdEdxInfo* GetTPCdEdxInfo() const {return fDetPid?fDetPid->GetTPCdEdxInfo():0;}
   Double_t  GetTPCmomentum()     const { return fDetPid?fDetPid->GetTPCmomentum():0.;  }
@@ -284,7 +286,7 @@ class AliAODTrack : public AliVTrack {
   UChar_t   GetTRDntrackletsPID() const;
   void      GetHMPIDpid(Double_t */*p*/) const { return; } // TODO: To be implemented properly with the new HMPID object
 
-  const AliAODEvent* GetAODEvent(){return fAODEvent;}
+  const AliAODEvent* GetAODEvent() const {return fAODEvent;}
   void SetAODEvent(const AliAODEvent* ptr){fAODEvent = ptr;}
 
   AliAODPid    *GetDetPid() const { return fDetPid; }
@@ -411,9 +413,11 @@ class AliAODTrack : public AliVTrack {
   Double_t      fTrackPhiOnEMCal;   // phi of track after being propagated to 430cm
   Double_t      fTrackEtaOnEMCal;   // eta of track after being propagated to 430cm
 
+  Double_t      fTPCsignalTuned;    //! TPC signal tuned on data when using MC
+
   const AliAODEvent* fAODEvent;     //! 
 
-  ClassDef(AliAODTrack, 17);
+  ClassDef(AliAODTrack, 18);
 };
 
 inline Bool_t  AliAODTrack::IsPrimaryCandidate() const
index a005b67289dd10f1b26a9e96392bfe867c472426..34c022fbeb3bc23d20f0c5abc9d04845c59d1405 100644 (file)
@@ -22,6 +22,7 @@
 //      Origin: Rosa Romita, GSI, r.romita@gsi.de
 //-----------------------------------------------------------------
 
+#include "TRandom.h"
 #include "AliLog.h"
 #include "AliPID.h"
 #include "AliAODpidUtil.h"
@@ -30,6 +31,8 @@
 #include "AliAODPid.h"
 #include "AliTRDPIDResponse.h"
 #include "AliESDtrack.h"
+#include "AliAODMCHeader.h"
+#include "AliAODMCParticle.h"
 
 ClassImp(AliAODpidUtil)
 
@@ -65,6 +68,72 @@ ClassImp(AliAODpidUtil)
   return 0;
 }
 //_________________________________________________________________________
+Float_t AliAODpidUtil::GetTPCsignalTunedOnData(const AliVTrack *t) const {
+    AliAODTrack *track = (AliAODTrack *) t;
+    Float_t dedx = track->GetTPCsignalTunedOnData();
+    if(dedx > 0) return dedx;
+
+    Double_t mom = t->GetTPCmomentum();
+
+    dedx = t->GetTPCsignal();
+    track->SetTPCsignalTunedOnData(dedx);
+
+    if(dedx < 20) return dedx;
+
+    
+    AliPID::EParticleType type = AliPID::kPion;
+    
+    AliAODMCHeader *mcHeader = dynamic_cast<AliAODMCHeader*>(track->GetAODEvent()->GetList()->FindObject(AliAODMCHeader::StdBranchName()));
+    if (mcHeader) {
+       TClonesArray *mcArray = (TClonesArray*)track->GetAODEvent()->GetList()->FindObject(AliAODMCParticle::StdBranchName());
+       
+       Bool_t kGood = kTRUE;
+       
+       Int_t iS = TMath::Abs(((AliAODMCParticle*)mcArray->At(TMath::Abs(t->GetLabel())))->GetPdgCode());
+       if(iS==AliPID::ParticleCode(AliPID::kElectron)){
+           type = AliPID::kElectron;
+       }
+       else if(iS==AliPID::ParticleCode(AliPID::kMuon)){
+           type = AliPID::kMuon;
+       }
+       else if(iS==AliPID::ParticleCode(AliPID::kPion)){
+           type = AliPID::kPion;
+       }
+       else if(iS==AliPID::ParticleCode(AliPID::kKaon)){
+           type = AliPID::kKaon;
+       }
+       else if(iS==AliPID::ParticleCode(AliPID::kProton)){
+           type = AliPID::kProton;
+       }
+       else if(iS==AliPID::ParticleCode(AliPID::kDeuteron)){ // d
+           type = AliPID::kDeuteron;
+       }
+       else if(iS==AliPID::ParticleCode(AliPID::kTriton)){ // t
+           type = AliPID::kTriton;
+       }
+       else if(iS==AliPID::ParticleCode(AliPID::kHe3)){ // 3He
+           type = AliPID::kHe3;
+       }
+       else if(iS==AliPID::ParticleCode(AliPID::kAlpha)){ // 4He
+           type = AliPID::kAlpha;
+       }
+       else
+           kGood = kFALSE;
+
+       if(kGood){
+           Double_t bethe=fTPCResponse.GetExpectedSignal(mom,type);
+           Double_t sigma=fTPCResponse.GetExpectedSigma(mom,t->GetTPCsignalN(),type);
+           dedx = gRandom->Gaus(bethe,sigma);
+
+           if(iS == AliPID::ParticleCode(AliPID::kHe3) || iS == AliPID::ParticleCode(AliPID::kAlpha)) dedx *= 5;
+       }
+
+    }
+
+    track->SetTPCsignalTunedOnData(dedx);
+    return dedx;
+}
+//_________________________________________________________________________
 void AliAODpidUtil::MakeTPCPID(const AliAODTrack *track,Double_t *p) const
 {
   //
index fb8f200dc4061e65f7a87bf6d77d11847ab474af..7f3718a508021ed779eda7acef82b67981185216 100644 (file)
@@ -41,6 +41,8 @@ public:
 
   virtual Float_t NumberOfSigmasTOF(const AliVParticle *vtrack, AliPID::EParticleType type) const;
   
+  Float_t GetTPCsignalTunedOnData(const AliVTrack *t) const;
+
 private:
   
   ClassDef(AliAODpidUtil,3)  // PID calculation class
index 41b0aa4e5c9f41c2a296fe05968dcdacc73f25c3..8bf4d56f4e4c680ca778896e378f7a3be4f53e23 100644 (file)
 #include "TArrayI.h"
 #include "TArrayF.h"
 
+#include "TRandom.h"
 #include "AliLog.h"
 #include "AliPID.h"
 #include "AliTOFHeader.h"
 #include "AliESDpid.h"
 #include "AliESDEvent.h"
 #include "AliESDtrack.h"
+#include "AliMCEvent.h"
+#include "AliMCEventHandler.h"
+#include "AliAnalysisManager.h"
+
 
 ClassImp(AliESDpid)
 
@@ -65,6 +70,73 @@ Int_t AliESDpid::MakePID(AliESDEvent *event, Bool_t TPConly, Float_t timeZeroTOF
   return 0;
 }
 //_________________________________________________________________________
+Float_t AliESDpid::GetTPCsignalTunedOnData(const AliVTrack *t) const {
+    AliESDtrack *track = (AliESDtrack *) t;
+    Float_t dedx = track->GetTPCsignalTunedOnData();
+
+    if(dedx > 0) return dedx;
+
+    Double_t mom = t->GetTPCmomentum();
+    dedx = t->GetTPCsignal();
+    track->SetTPCsignalTunedOnData(dedx);
+    if(dedx < 20) return dedx;
+
+    AliPID::EParticleType type = AliPID::kPion;
+
+    AliMCEventHandler* eventHandler=NULL;
+    eventHandler = dynamic_cast<AliMCEventHandler*> (AliAnalysisManager::GetAnalysisManager()->GetMCtruthEventHandler());
+    if (eventHandler) {
+       AliMCEvent* mcEvent = eventHandler->MCEvent();
+       if(mcEvent){
+           Bool_t kGood = kTRUE;
+           AliMCParticle *MCpart = (AliMCParticle *) mcEvent->GetTrack(TMath::Abs(t->GetLabel()));
+           TParticle *part = MCpart->Particle();
+           
+           Int_t iS = TMath::Abs(part->GetPdgCode());
+
+           if(iS==AliPID::ParticleCode(AliPID::kElectron)){
+               type = AliPID::kElectron;
+           }
+           else if(iS==AliPID::ParticleCode(AliPID::kMuon)){
+               type = AliPID::kMuon;
+           }
+           else if(iS==AliPID::ParticleCode(AliPID::kPion)){
+               type = AliPID::kPion;
+           }
+           else if(iS==AliPID::ParticleCode(AliPID::kKaon)){
+               type = AliPID::kKaon;
+           }
+           else if(iS==AliPID::ParticleCode(AliPID::kProton)){
+               type = AliPID::kProton;
+           }
+           else if(iS==AliPID::ParticleCode(AliPID::kDeuteron)){ // d
+               type = AliPID::kDeuteron;
+           }
+           else if(iS==AliPID::ParticleCode(AliPID::kTriton)){ // t
+               type = AliPID::kTriton;
+           }
+           else if(iS==AliPID::ParticleCode(AliPID::kHe3)){ // 3He
+               type = AliPID::kHe3;
+           }
+           else if(iS==AliPID::ParticleCode(AliPID::kAlpha)){ // 4He
+               type = AliPID::kAlpha;
+           }
+           else
+               kGood = kFALSE;
+
+           if(kGood){
+               Double_t bethe=fTPCResponse.GetExpectedSignal(mom,type);
+               Double_t sigma=fTPCResponse.GetExpectedSigma(mom,t->GetTPCsignalN(),type);
+               dedx = gRandom->Gaus(bethe,sigma);
+               if(iS == AliPID::ParticleCode(AliPID::kHe3) || iS == AliPID::ParticleCode(AliPID::kAlpha)) dedx *= 5;
+           }
+       }
+    }
+
+    track->SetTPCsignalTunedOnData(dedx);
+    return dedx;
+}
+//_________________________________________________________________________
 void AliESDpid::MakeTPCPID(AliESDtrack *track) const
 {
   //
index 0dcaa9e2499653b6b82ccd2601177715267d456b..50c8043bd0056f408948f197d5741f4721fbe8f6 100644 (file)
@@ -43,6 +43,8 @@ public:
   void SetNMaxSigmaTOFTPCMismatch(Float_t range) {fRangeTOFMismatch=range;}
   Float_t GetNMaxSigmaTOFTPCMismatch() const {return fRangeTOFMismatch;}
 
+  Float_t GetTPCsignalTunedOnData(const AliVTrack *t) const;
+
 private:
   Float_t           fRangeTOFMismatch; // nSigma max for TOF matching with TPC
 
index 3aa9b6ac7b65b82e710b8840c2676a3c6423c3fd..0a148c2b643d834a643ba579f19b5015c1f166f1 100644 (file)
@@ -200,6 +200,7 @@ AliESDtrack::AliESDtrack() :
   fGlobalChi2(0),
   fITSsignal(0),
   fTPCsignal(0),
+  fTPCsignalTuned(0),
   fTPCsignalS(0),
   fTPCdEdxInfo(0),
   fTRDsignal(0),
@@ -314,6 +315,7 @@ AliESDtrack::AliESDtrack(const AliESDtrack& track):
   fGlobalChi2(track.fGlobalChi2),
   fITSsignal(track.fITSsignal),
   fTPCsignal(track.fTPCsignal),
+  fTPCsignalTuned(track.fTPCsignalTuned),
   fTPCsignalS(track.fTPCsignalS),
   fTPCdEdxInfo(0),
   fTRDsignal(track.fTRDsignal),
@@ -440,6 +442,7 @@ AliESDtrack::AliESDtrack(const AliVTrack *track) :
   fGlobalChi2(0),
   fITSsignal(0),
   fTPCsignal(0),
+  fTPCsignalTuned(0),
   fTPCsignalS(0),
   fTPCdEdxInfo(0),
   fTRDsignal(0),
@@ -586,6 +589,7 @@ AliESDtrack::AliESDtrack(TParticle * part) :
   fGlobalChi2(0),
   fITSsignal(0),
   fTPCsignal(0),
+  fTPCsignalTuned(0),
   fTPCsignalS(0),
   fTPCdEdxInfo(0),
   fTRDsignal(0),
@@ -916,6 +920,7 @@ AliESDtrack &AliESDtrack::operator=(const AliESDtrack &source){
   fITSsignal  = source.fITSsignal;     
   for (Int_t i=0;i<4;i++) {fITSdEdxSamples[i]=source.fITSdEdxSamples[i];}
   fTPCsignal  = source.fTPCsignal;     
+  fTPCsignalTuned  = source.fTPCsignalTuned;
   fTPCsignalS = source.fTPCsignalS;    
   for(int i = 0; i< 4;++i){
     fTPCPoints[i] = source.fTPCPoints[i];  
@@ -1065,6 +1070,7 @@ Bool_t AliESDtrack::FillTPCOnlyTrack(AliESDtrack &track){
   track.fTPCchi2 = fTPCchi2; 
   track.fTPCchi2Iter1 = fTPCchi2Iter1; 
   track.fTPCsignal = fTPCsignal;
+  track.fTPCsignalTuned = fTPCsignalTuned;
   track.fTPCsignalS = fTPCsignalS;
   for(int i = 0;i<4;++i)track.fTPCPoints[i] = fTPCPoints[i];
 
@@ -1144,6 +1150,7 @@ void AliESDtrack::MakeMiniESDtrack(){
   fTPCClusterMap = 0;  
   fTPCSharedMap = 0;  
   fTPCsignal= 0;      
+  fTPCsignalTuned= 0;
   fTPCsignalS= 0;      
   fTPCsignalN= 0;      
   for (Int_t i=0;i<AliPID::kSPECIES;i++) fTPCr[i]=0; 
index efd80cf25f85a341395fb9547c23cac79e99b5da..83218a6c0ef52185a229d60bd4ecefedac0bb8e6 100644 (file)
@@ -218,10 +218,14 @@ public:
   void    SetTPCsignal(Float_t signal, Float_t sigma, UChar_t npoints){ 
      fTPCsignal = signal; fTPCsignalS = sigma; fTPCsignalN = npoints;
   }
+  void    SetTPCsignalTunedOnData(Float_t signal){
+      fTPCsignalTuned = signal;
+  }
   void  SetTPCdEdxInfo(AliTPCdEdxInfo * dEdxInfo); 
 
   AliTPCdEdxInfo * GetTPCdEdxInfo() const {return fTPCdEdxInfo;}
   Double_t GetTPCsignal() const {return fTPCsignal;}
+  Double_t GetTPCsignalTunedOnData() const {return fTPCsignalTuned;}
   Double_t GetTPCsignalSigma() const {return fTPCsignalS;}
   UShort_t GetTPCsignalN() const {return fTPCsignalN;}
   Double_t GetTPCmomentum() const {return fIp?fIp->GetP():GetP();}
@@ -477,6 +481,7 @@ protected:
   Double32_t  fITSdEdxSamples[4]; // [0.,0.,10] ITS dE/dx samples
 
   Double32_t  fTPCsignal;        // [0.,0.,10] detector's PID signal
+  Double32_t  fTPCsignalTuned;   //! [0.,0.,10] detector's PID signal tuned on data when using MC
   Double32_t  fTPCsignalS;       // [0.,0.,10] RMS of dEdx measurement
   AliTPCdEdxInfo * fTPCdEdxInfo; // object containing dE/dx information for different pad regions
   Double32_t  fTPCPoints[4];     // [0.,0.,10] TPC points -first, max. dens, last and max density
@@ -536,7 +541,7 @@ protected:
   static bool fgkOnlineMode; //! indicate the online mode to skip some of the functionality
 
   AliESDtrack & operator=(const AliESDtrack & );
-  ClassDef(AliESDtrack,65)  //ESDtrack 
+  ClassDef(AliESDtrack,66)  //ESDtrack 
 };
 
 
index dd339cda9c96108e4dfc4189736d2b6fb4853382..2148be25ae90c93a2c326cb4f987645caeb54b57 100644 (file)
@@ -72,7 +72,8 @@ fTOFtail(1.1),
 fTOFPIDParams(0x0),
 fEMCALPIDParams(0x0),
 fCurrentEvent(0x0),
-fCurrCentrality(0.0)
+fCurrCentrality(0.0),
+fTuneMConData(kFALSE)
 {
   //
   // default ctor
@@ -123,7 +124,8 @@ fTOFtail(1.1),
 fTOFPIDParams(0x0),
 fEMCALPIDParams(0x0),
 fCurrentEvent(0x0),
-fCurrCentrality(0.0)
+fCurrCentrality(0.0),
+fTuneMConData(kFALSE)
 {
   //
   // copy ctor
@@ -371,6 +373,8 @@ AliPIDResponse::EDetPidStatus AliPIDResponse::ComputeTPCProbability  (const AliV
   Double_t dedx=track->GetTPCsignal();
   Bool_t mismatch=kTRUE/*, heavy=kTRUE*/;
 
+  if(fTuneMConData) dedx = this->GetTPCsignalTunedOnData(track);
+
   for (Int_t j=0; j<AliPID::kSPECIES; j++) {
     AliPID::EParticleType type=AliPID::EParticleType(j);
     Double_t bethe=fTPCResponse.GetExpectedSignal(mom,type);
@@ -750,8 +754,9 @@ void AliPIDResponse::SetTPCParametrisation()
   TString datatype="DATA";
   //in case of mc fRecoPass is per default 1
   if (fIsMC) {
-    datatype="MC";
-    fRecoPass=1;
+      datatype="MC";
+      if(!fTuneMConData) datatype="MC";
+      fRecoPass=1;
   }
   
   //
@@ -763,7 +768,7 @@ void AliPIDResponse::SetTPCParametrisation()
 
   // period
   TString period=fLHCperiod;
-  if (fIsMC) period=fMCperiodTPC;
+  if (fIsMC && !fTuneMConData) period=fMCperiodTPC;
 
   AliInfo(Form("Searching splines for: %s %s PASS%d %s",datatype.Data(),period.Data(),fRecoPass,fBeamType.Data()));
   Bool_t found=kFALSE;
@@ -771,13 +776,15 @@ void AliPIDResponse::SetTPCParametrisation()
   //set the new PID splines
   //
   if (fArrPidResponseMaster){
+    Int_t recopass = fRecoPass;
+    if(fTuneMConData) recopass = fRecoPassUser;
     TObject *grAll=0x0;
     //for MC don't use period information
 //     if (fIsMC) period="[A-Z0-9]*";
     //for MC use MC period information
 //pattern for the default entry (valid for all particles)
-    TPRegexp reg(Form("TSPLINE3_%s_([A-Z]*)_%s_PASS%d_%s_MEAN",datatype.Data(),period.Data(),fRecoPass,fBeamType.Data()));
-    
+    TPRegexp reg(Form("TSPLINE3_%s_([A-Z]*)_%s_PASS%d_%s_MEAN",datatype.Data(),period.Data(),recopass,fBeamType.Data()));
+   
     //loop over entries and filter them
     for (Int_t iresp=0; iresp<fArrPidResponseMaster->GetEntriesFast();++iresp){
       TObject *responseFunction=fArrPidResponseMaster->At(iresp);
index 3c7530e1e772a0807114c5db287e6eec6b1dec4f..df581906260d810e5bb425e90c4bf859083e00dc 100644 (file)
@@ -100,6 +100,10 @@ public:
   void SetTOFtail(Float_t tail=1.1){if(tail > 0) fTOFtail=tail; else printf("TOF tail should be greater than 0 (nothing done)\n");};
   void SetTOFResponse(AliVEvent *vevent,EStartTimeType_t option);
 
+  virtual Float_t GetTPCsignalTunedOnData(const AliVTrack *t) const {return t->GetTPCsignal();};
+  Bool_t IsTunedOnData() const {return fTuneMConData;};
+  void SetTunedOnData(Bool_t flag=kTRUE,Int_t recoPass=0){fTuneMConData = flag; if(recoPass>0) fRecoPassUser = recoPass;};
+
   AliPIDResponse(const AliPIDResponse &other);
   AliPIDResponse& operator=(const AliPIDResponse &other);
 
@@ -144,6 +148,8 @@ private:
 
   Float_t fCurrCentrality;             //! current centrality
   
+  Bool_t fTuneMConData;                // switch to force the MC to be similar to data (dE/dx)
+
   void ExecNewRun();
   
   //
@@ -173,7 +179,7 @@ private:
   //
   void SetRecoInfo();
   
-  ClassDef(AliPIDResponse, 7);  //PID response handling
+  ClassDef(AliPIDResponse, 8);  //PID response handling
 };
 
 inline Float_t AliPIDResponse::NumberOfSigmasTPC(const AliVParticle *vtrack, AliPID::EParticleType type) const {
@@ -182,6 +188,8 @@ inline Float_t AliPIDResponse::NumberOfSigmasTPC(const AliVParticle *vtrack, Ali
   Double_t sig  = track->GetTPCsignal();
   UInt_t   sigN = track->GetTPCsignalN();
 
+  if(fTuneMConData) sig = GetTPCsignalTunedOnData(track);
+
   Double_t nSigma = -999.;
   if (sigN>0) nSigma=fTPCResponse.GetNumberOfSigmas(mom,sig,sigN,type);
 
index 8dc4ee73683a908c1e20581f75638eb50e8cdfa7..93811c2bc0e23d0d8ab62ecc13f07ff464cb8683 100644 (file)
@@ -125,7 +125,7 @@ Double_t AliTPCPIDResponse::GetExpectedSignal(const Float_t mom,
   // assigned clusters and/or the track dip angle, for example.  
   //
   
-  Double_t mass=AliPID::ParticleMass(n);
+  Double_t mass=AliPID::ParticleMassZ(n);
   if (!fUseDatabase||fResponseFunctions.GetEntriesFast()>AliPID::kUnknown) return Bethe(mom/mass);
   //
   TSpline3 * responseFunction = (TSpline3 *) fResponseFunctions.UncheckedAt(n);
index 80ab1a1a29badec982395d26c8b4688d4a00bf9c..65f2bcae545ca501d618b8d736f8cfd8ce3ad10f 100644 (file)
@@ -78,6 +78,7 @@ public:
 \r
   virtual Double_t  GetITSsignal()       const {return 0.;}\r
   virtual Double_t  GetTPCsignal()       const {return 0.;}\r
+  virtual Double_t  GetTPCsignalTunedOnData() const {return 0.;}
   virtual UShort_t  GetTPCsignalN()      const {return 0 ;}\r
   virtual Double_t  GetTPCmomentum()     const {return 0.;}\r
   virtual Double_t  GetTOFsignal()       const {return 0.;}\r