]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
fix flow calculation in scalar product
authorsnelling <snelling@f7af4fe6-9843-0410-8265-dc069ae4e863>
Tue, 3 Mar 2009 17:37:17 +0000 (17:37 +0000)
committersnelling <snelling@f7af4fe6-9843-0410-8265-dc069ae4e863>
Tue, 3 Mar 2009 17:37:17 +0000 (17:37 +0000)
PWG2/FLOW/AliFlowCommon/AliFlowAnalysisWithScalarProduct.cxx
PWG2/FLOW/AliFlowCommon/AliFlowAnalysisWithScalarProduct.h
PWG2/FLOW/AliFlowTasks/AliAnalysisTaskScalarProduct.cxx

index 37ddb29b07d45a987f379e55cef1b680d9d881af..37f24ca9f0f862936a9bd05999f72aee4586e4d8 100644 (file)
@@ -51,10 +51,7 @@ ClassImp(AliFlowAnalysisWithScalarProduct)
    fHistProUQetaPOI(NULL),
    fHistProUQPtRP(NULL),
    fHistProUQPtPOI(NULL),
-   fHistProVetaRP(NULL),
-   fHistProVetaPOI(NULL),
-   fHistProVPtRP(NULL),
-   fHistProVPtPOI(NULL),
+   fHistProQaQb(NULL),
    fCommonHists(NULL),
    fCommonHistsRes(NULL)
 {
@@ -126,25 +123,9 @@ void AliFlowAnalysisWithScalarProduct::Init() {
   fHistProUQPtPOI->SetYTitle("<uQ>");
   fHistList->Add(fHistProUQPtPOI);
 
-  fHistProVetaRP = new TProfile("Flow_VetaRP_SP","Flow_VetaRP_SP",iNbinsEta,dEtaMin,dEtaMax);
-  fHistProVetaRP->SetXTitle("{eta}");
-  fHistProVetaRP->SetYTitle("v_{2}");
-  fHistList->Add(fHistProVetaRP);
-
-  fHistProVetaPOI = new TProfile("Flow_VetaPOI_SP","Flow_VetaPOI_SP",iNbinsEta,dEtaMin,dEtaMax);
-  fHistProVetaPOI->SetXTitle("{eta}");
-  fHistProVetaPOI->SetYTitle("v_{2}");
-  fHistList->Add(fHistProVetaPOI);
-
-  fHistProVPtRP = new TProfile("Flow_VPtRP_SP","Flow_VPtRP_SP",iNbinsPt,dPtMin,dPtMax);
-  fHistProVPtRP->SetXTitle("p_t (GeV)");
-  fHistProVPtRP->SetYTitle("v_{2}");
-  fHistList->Add(fHistProVPtRP);
-
-  fHistProVPtPOI = new TProfile("Flow_VPtPOI_SP","Flow_VPtPOI_SP",iNbinsPt,dPtMin,dPtMax);
-  fHistProVPtPOI->SetXTitle("p_t (GeV)");
-  fHistProVPtPOI->SetYTitle("v_{2}");
-  fHistList->Add(fHistProVPtPOI);
+  fHistProQaQb = new TProfile("Flow_QaQb_SP","Flow_QaQb_SP", 1, -0.5, 0.5);
+  fHistProQaQb->SetYTitle("<QaQb>");
+  fHistList->Add(fHistProQaQb);
 
   fCommonHists = new AliFlowCommonHist("AliFlowCommonHistSP");
   fHistList->Add(fCommonHists);
@@ -169,9 +150,8 @@ void AliFlowAnalysisWithScalarProduct::Make(AliFlowEventSimple* anEvent) {
     //get Q vectors for the subevents
     AliFlowVector vQa = anEvent->GetQsub(-0.9,0.);  //Random subevents?
     AliFlowVector vQb = anEvent->GetQsub(0.,0.9);
-    //CHECK THIS!
-    Double_t dQaQb = vQa*vQb; //CHECK this!
-    dQaQb = 2*TMath::Sqrt(dQaQb);
+    Double_t dQaQb = vQa*vQb; 
+    fHistProQaQb -> Fill(0.,dQaQb);    
                 
     //loop over the tracks of the event
     AliFlowTrackSimple*   pTrack = NULL; 
@@ -200,22 +180,16 @@ void AliFlowAnalysisWithScalarProduct::Make(AliFlowEventSimple* anEvent) {
 
          //dUQ = scalar product of vU and vQm
          Double_t dUQ = vU * vQm;
-         Double_t dV = 0.;
-         if (dQaQb != 0.) {dV = dUQ/dQaQb;}
          Double_t dPt = pTrack->Pt();
          Double_t dEta = pTrack->Eta();
          //fill the profile histograms
          if (pTrack->UseForIntegratedFlow()) {
            fHistProUQetaRP -> Fill(dEta,dUQ);
            fHistProUQPtRP -> Fill(dPt,dUQ);
-           fHistProVetaRP -> Fill(dEta,dV);
-           fHistProVPtRP -> Fill(dPt,dV);
          }
          if (pTrack->UseForDifferentialFlow()) {
            fHistProUQetaPOI -> Fill(dEta,dUQ);
            fHistProUQPtPOI -> Fill(dPt,dUQ);
-           fHistProVetaPOI -> Fill(dEta,dV);
-           fHistProVPtPOI -> Fill(dPt,dV);
          }  
        }//track selected
       }//loop over tracks
@@ -228,7 +202,7 @@ void AliFlowAnalysisWithScalarProduct::Make(AliFlowEventSimple* anEvent) {
   //--------------------------------------------------------------------    
 void AliFlowAnalysisWithScalarProduct::Finish() {
    
-  //fill the commonhist results and calculate integrated flow
+  //calculate flow and fill the AliFlowCommonHistResults
   if (fDebug) cout<<"AliFlowAnalysisWithScalarProduct::Finish()"<<endl;
 
   cout<<"*************************************"<<endl;
@@ -236,85 +210,132 @@ void AliFlowAnalysisWithScalarProduct::Finish() {
   cout<<"      Integrated flow from           "<<endl;
   cout<<"         Scalar product              "<<endl;
   cout<<endl;
-
-  //copy content of profile into TH1D and fill the AliFlowCommonHistResults
+  
   Int_t iNbinsPt  = AliFlowCommonConstants::GetNbinsPt();
   Int_t iNbinsEta = AliFlowCommonConstants::GetNbinsEta();
-
-  //v as a function of eta for RP selection
-  for(Int_t b=0;b<iNbinsEta;b++) {
-    Double_t dv2pro  = fHistProVetaRP->GetBinContent(b);
-    Double_t dv2err  = fHistProVetaRP->GetBinError(b); //copy error for now
-    //fill TH1D
-    fCommonHistsRes->FillDifferentialFlowEtaRP(b, dv2pro, dv2err); 
-  } //loop over bins b
+  Double_t dQaQbAv  = fHistProQaQb->GetBinContent(1); //average over events
+  Double_t dQaQbErr = fHistProQaQb->GetBinError(1);
+  if (dQaQbAv <= 0.){
+    //set v to -0
+    fCommonHistsRes->FillIntegratedFlowRP(-0.,0.);
+    fCommonHistsRes->FillIntegratedFlow(-0.,0.);
+    cout<<"dV(RP) = -0. +- 0."<<endl;
+    fCommonHistsRes->FillIntegratedFlowPOI(-0.,0.);
+    cout<<"dV(POI) = -0. +- 0."<<endl;
+  } else {
+    Double_t dQaQbSqrt = 2*TMath::Sqrt(dQaQbAv);
+    Double_t dQaQbSqrtErr = (dQaQbSqrt/2)*(dQaQbErr/dQaQbAv);
+    Double_t dQaQbSqrtErrRel = dQaQbSqrtErr/dQaQbSqrt;
+    Double_t dQaQbSqrtErrRel2 = dQaQbSqrtErrRel*dQaQbSqrtErrRel;
+
+    //v as a function of eta for RP selection
+    for(Int_t b=0;b<iNbinsEta;b++) {
+      Double_t duQpro = fHistProUQetaRP->GetBinContent(b);
+      Double_t duQerr = fHistProUQetaRP->GetBinError(b); //copy error for now
+      Double_t duQerrRel = 0.;
+      if (duQpro != 0.) {duQerrRel = duQerr/duQpro;}
+      Double_t duQerrRel2 = duQerrRel*duQerrRel;
+
+      Double_t dv2pro     = duQpro/dQaQbSqrt;
+      Double_t dv2errRel2 = duQerrRel2 + dQaQbSqrtErrRel2;
+      Double_t dv2errRel  = TMath::Sqrt(dv2errRel2);
+      Double_t dv2err     = dv2pro*dv2errRel; 
+      //fill TH1D
+      fCommonHistsRes->FillDifferentialFlowEtaRP(b, dv2pro, dv2err); 
+    } //loop over bins b
     
-  //v as a function of eta for POI selection
-  for(Int_t b=0;b<iNbinsEta;b++) {
-    Double_t dv2pro  = fHistProVetaPOI->GetBinContent(b);
-    Double_t dv2err  = fHistProVetaPOI->GetBinError(b); //copy error for now
-    //fill TH1D
-    fCommonHistsRes->FillDifferentialFlowEtaPOI(b, dv2pro, dv2err); 
-  } //loop over bins b
+    //v as a function of eta for POI selection
+    for(Int_t b=0;b<iNbinsEta;b++) {
+      Double_t duQpro = fHistProUQetaPOI->GetBinContent(b);
+      Double_t duQerr = fHistProUQetaPOI->GetBinError(b); //copy error for now
+      Double_t duQerrRel = 0.;
+      if (duQpro != 0.) {duQerrRel = duQerr/duQpro;}
+      Double_t duQerrRel2 = duQerrRel*duQerrRel;
+
+      Double_t dv2pro     = duQpro/dQaQbSqrt;
+      Double_t dv2errRel2 = duQerrRel2 + dQaQbSqrtErrRel2;
+      Double_t dv2errRel  = TMath::Sqrt(dv2errRel2);
+      Double_t dv2err     = dv2pro*dv2errRel; 
+      //fill TH1D
+      fCommonHistsRes->FillDifferentialFlowEtaPOI(b, dv2pro, dv2err); 
+    } //loop over bins b
     
-  //v as a function of Pt for RP selection
-  TH1F* fHistPtInt = fCommonHists->GetHistPtInt(); //for calculating integrated flow
-  Double_t dVRP = 0.;
-  Double_t dSum = 0.;
-  Double_t dErrV =0.;
-
-  for(Int_t b=0;b<iNbinsPt;b++) {
-    Double_t dv2pro  = fHistProVPtRP->GetBinContent(b);
-    Double_t dv2err  = fHistProVPtRP->GetBinError(b); //copy error for now
-    //fill TH1D
-    fCommonHistsRes->FillDifferentialFlowPtRP(b, dv2pro, dv2err);
-    //calculate integrated flow for RP selection
-    if (fHistPtInt){
-      Double_t dYieldPt = fHistPtInt->GetBinContent(b);
-      dVRP += dv2pro*dYieldPt;
-      dSum +=dYieldPt;
-      dErrV += dYieldPt*dYieldPt*dv2err*dv2err;
-    } else { cout<<"fHistPtDiff is NULL"<<endl; }
-  } //loop over bins b
-
-  if (dSum != 0.) {
-    dVRP /= dSum; //the pt distribution should be normalised
-    dErrV /= (dSum*dSum);
-    dErrV = TMath::Sqrt(dErrV);
-  }
-  fCommonHistsRes->FillIntegratedFlowRP(dVRP,dErrV);
-
-  cout<<"dV(RP) = "<<dVRP<<" +- "<<dErrV<<endl;
+    //v as a function of Pt for RP selection
+    TH1F* fHistPtInt = fCommonHists->GetHistPtInt(); //for calculating integrated flow
+    Double_t dVRP = 0.;
+    Double_t dSum = 0.;
+    Double_t dErrV =0.;
+
+    for(Int_t b=0;b<iNbinsPt;b++) {
+      Double_t duQpro = fHistProUQPtRP->GetBinContent(b);
+      Double_t duQerr = fHistProUQPtRP->GetBinError(b); //copy error for now
+      Double_t duQerrRel = 0.;
+      if (duQpro != 0.) {duQerrRel = duQerr/duQpro;}
+      Double_t duQerrRel2 = duQerrRel*duQerrRel;
+
+      Double_t dv2pro     = duQpro/dQaQbSqrt;
+      Double_t dv2errRel2 = duQerrRel2 + dQaQbSqrtErrRel2;
+      Double_t dv2errRel  = TMath::Sqrt(dv2errRel2);
+      Double_t dv2err     = dv2pro*dv2errRel; 
+      //fill TH1D
+      fCommonHistsRes->FillDifferentialFlowPtRP(b, dv2pro, dv2err);
+      //calculate integrated flow for RP selection
+      if (fHistPtInt){
+       Double_t dYieldPt = fHistPtInt->GetBinContent(b);
+       dVRP += dv2pro*dYieldPt;
+       dSum +=dYieldPt;
+       dErrV += dYieldPt*dYieldPt*dv2err*dv2err;
+      } else { cout<<"fHistPtDiff is NULL"<<endl; }
+    } //loop over bins b
+
+    if (dSum != 0.) {
+      dVRP /= dSum; //the pt distribution should be normalised
+      dErrV /= (dSum*dSum);
+      dErrV = TMath::Sqrt(dErrV);
+    }
+    fCommonHistsRes->FillIntegratedFlowRP(dVRP,dErrV);
+    fCommonHistsRes->FillIntegratedFlow(dVRP,dErrV);
+
+    cout<<"dV(RP) = "<<dVRP<<" +- "<<dErrV<<endl;
        
-  //v as a function of Pt for POI selection 
-  TH1F* fHistPtDiff = fCommonHists->GetHistPtDiff(); //for calculating integrated flow
-  Double_t dVPOI = 0.;
-  dSum = 0.;
-  dErrV =0.;
+    //v as a function of Pt for POI selection 
+    TH1F* fHistPtDiff = fCommonHists->GetHistPtDiff(); //for calculating integrated flow
+    Double_t dVPOI = 0.;
+    dSum = 0.;
+    dErrV =0.;
   
-  for(Int_t b=0;b<iNbinsPt;b++) {
-    Double_t dv2pro = fHistProVPtPOI->GetBinContent(b);
-    Double_t dv2err = fHistProVPtPOI->GetBinError(b); //copy error for now
-    //fill TH1D
-    fCommonHistsRes->FillDifferentialFlowPtPOI(b, dv2pro, dv2err); 
-
-    //calculate integrated flow for POI selection
-    if (fHistPtDiff){
-      Double_t dYieldPt = fHistPtDiff->GetBinContent(b);
-      dVPOI += dv2pro*dYieldPt;
-      dSum +=dYieldPt;
-      dErrV += dYieldPt*dYieldPt*dv2err*dv2err;
-    } else { cout<<"fHistPtDiff is NULL"<<endl; }
-  } //loop over bins b
-
-  if (dSum != 0.) {
-    dVPOI /= dSum; //the pt distribution should be normalised
-    dErrV /= (dSum*dSum);
-    dErrV = TMath::Sqrt(dErrV);
+    for(Int_t b=0;b<iNbinsPt;b++) {
+      Double_t duQpro = fHistProUQPtPOI->GetBinContent(b);
+      Double_t duQerr = fHistProUQPtPOI->GetBinError(b); //copy error for now
+      Double_t duQerrRel = 0.;
+      if (duQpro != 0.) {duQerrRel = duQerr/duQpro;}
+      Double_t duQerrRel2 = duQerrRel*duQerrRel;
+
+      Double_t dv2pro     = duQpro/dQaQbSqrt;
+      Double_t dv2errRel2 = duQerrRel2 + dQaQbSqrtErrRel2;
+      Double_t dv2errRel  = TMath::Sqrt(dv2errRel2);
+      Double_t dv2err     = dv2pro*dv2errRel; 
+      //fill TH1D
+      fCommonHistsRes->FillDifferentialFlowPtPOI(b, dv2pro, dv2err); 
+
+      //calculate integrated flow for POI selection
+      if (fHistPtDiff){
+       Double_t dYieldPt = fHistPtDiff->GetBinContent(b);
+       dVPOI += dv2pro*dYieldPt;
+       dSum +=dYieldPt;
+       dErrV += dYieldPt*dYieldPt*dv2err*dv2err;
+      } else { cout<<"fHistPtDiff is NULL"<<endl; }
+    } //loop over bins b
+
+    if (dSum != 0.) {
+      dVPOI /= dSum; //the pt distribution should be normalised
+      dErrV /= (dSum*dSum);
+      dErrV = TMath::Sqrt(dErrV);
+    }
+    fCommonHistsRes->FillIntegratedFlowPOI(dVPOI,dErrV);
+
+    cout<<"dV(POI) = "<<dVPOI<<" +- "<<dErrV<<endl;
   }
-  fCommonHistsRes->FillIntegratedFlowPOI(dVPOI,dErrV);
-
-  cout<<"dV(POI) = "<<dVPOI<<" +- "<<dErrV<<endl;
   cout<<endl;
   cout<<"*************************************"<<endl;
   cout<<"*************************************"<<endl;           
index b3c3887e1a6c1c435ddfafc7c442f7d1634a55c7..379e99fd10c07d89c9ba9ae962d8fefc32cb2fc2 100644 (file)
@@ -37,13 +37,13 @@ class AliFlowAnalysisWithScalarProduct {
    void    Make(AliFlowEventSimple* anEvent);        //calculates variables and fills histograms
    void    Finish();                                 //saves histograms
    void    WriteHistograms(TString* outputFileName); //writes histograms locally
-   void    WriteHistograms(TString outputFileName); //writes histograms locally
+   void    WriteHistograms(TString outputFileName);  //writes histograms locally
 
-   void      SetDebug(Bool_t kt)            { this->fDebug = kt ; }
-   Bool_t    GetDebug() const               { return this->fDebug ; }
+   void      SetDebug(Bool_t kt)   { this->fDebug = kt ; }
+   Bool_t    GetDebug() const      { return this->fDebug ; }
 
    // Output 
-   TList*   GetHistList() const             { return this->fHistList ; }     // Gets output histogram list
+   TList*   GetHistList() const    { return this->fHistList ; }     // Gets output histogram list
   
    TProfile* GetHistProUQetaRP() const {return this->fHistProUQetaRP;}   
    void      SetHistProUQetaRP(TProfile* const aHistProUQetaRP) {this->fHistProUQetaRP = aHistProUQetaRP;}
@@ -53,14 +53,8 @@ class AliFlowAnalysisWithScalarProduct {
    void      SetHistProUQPtRP(TProfile* const aHistProUQPtRP) {this->fHistProUQPtRP = aHistProUQPtRP;}
    TProfile* GetHistProUQPtPOI() const {return this->fHistProUQPtPOI;}
    void      SetHistProUQPtPOI(TProfile* const aHistProUQPtPOI) {this->fHistProUQPtPOI = aHistProUQPtPOI;}
-   TProfile* GetHistProVetaRP() const {return this->fHistProVetaRP;}
-   void      SetHistProVetaRP(TProfile* const aHistProVetaRP) {this->fHistProVetaRP = aHistProVetaRP;}
-   TProfile* GetHistProVetaPOI() const {return this->fHistProVetaPOI;} 
-   void      SetHistProVetaPOI(TProfile* const aHistProVetaPOI) {this->fHistProVetaPOI = aHistProVetaPOI;}
-   TProfile* GetHistProVPtRP() const {return this->fHistProVPtRP;} 
-   void      SetHistProVPtRP(TProfile* const aHistProVPtRP) {this->fHistProVPtRP = aHistProVPtRP;}
-   TProfile* GetHistProVPtPOI() const {return this->fHistProVPtPOI;} 
-   void      SetHistProVPtPOI(TProfile* const aHistProVPtPOI) {this->fHistProVPtPOI = aHistProVPtPOI;}
+   TProfile* GetHistProQaQb() const {return this->fHistProQaQb;}
+   void      SetHistProQaQb(TProfile* const aHistProQaQb) {this->fHistProQaQb = aHistProQaQb;}
 
    AliFlowCommonHist* GetCommonHists() const {return this->fCommonHists; }
    void SetCommonHists(AliFlowCommonHist* const someCommonHists) {this->fCommonHists = someCommonHists; }
@@ -73,18 +67,15 @@ class AliFlowAnalysisWithScalarProduct {
    AliFlowAnalysisWithScalarProduct(const AliFlowAnalysisWithScalarProduct& anAnalysis);            //copy constructor
    AliFlowAnalysisWithScalarProduct& operator=(const AliFlowAnalysisWithScalarProduct& anAnalysis); //assignment operator 
       
-   Int_t              fEventNumber;      // event counter
-   Bool_t             fDebug ;           // flag for analysis: more print statements
-
-   TList*             fHistList;         //list to hold all output histograms  
-   TProfile*          fHistProUQetaRP;   //uQ(eta) for RP
-   TProfile*          fHistProUQetaPOI;  //uQ(eta) for POI
-   TProfile*          fHistProUQPtRP;    //uQ(pt) for RP
-   TProfile*          fHistProUQPtPOI;   //uQ(pt) for POI
-   TProfile*          fHistProVetaRP;    //v2(eta) for RP
-   TProfile*          fHistProVetaPOI;   //v2(eta) for POI
-   TProfile*          fHistProVPtRP;     //v2(pt) for RP
-   TProfile*          fHistProVPtPOI;    //v2(pt) for POI
+   Int_t      fEventNumber;      // event counter
+   Bool_t     fDebug ;           // flag for analysis: more print statements
+
+   TList*     fHistList;         //list to hold all output histograms  
+   TProfile*  fHistProUQetaRP;   //uQ(eta) for RP
+   TProfile*  fHistProUQetaPOI;  //uQ(eta) for POI
+   TProfile*  fHistProUQPtRP;    //uQ(pt) for RP
+   TProfile*  fHistProUQPtPOI;   //uQ(pt) for POI
+   TProfile*  fHistProQaQb;      //average of QaQb (for event plane resolution)
    AliFlowCommonHist*        fCommonHists;    //control histograms
    AliFlowCommonHistResults* fCommonHistsRes; //results histograms
 
index ed808f264910f2637278a68e2403023176868973..402a7146ad62aa76a8e36fe63bd7389d2775557b 100644 (file)
@@ -140,18 +140,20 @@ void AliAnalysisTaskScalarProduct::Terminate(Option_t *)
       (fListHistos->FindObject("AliFlowCommonHistSP"));
     AliFlowCommonHistResults *pCommonHistResults = dynamic_cast<AliFlowCommonHistResults*> 
       (fListHistos->FindObject("AliFlowCommonHistResultsSP"));
-    TProfile* pHistProVetaRP  = dynamic_cast<TProfile*>(fListHistos->FindObject("Flow_VetaRP_SP"));
-    TProfile* pHistProVetaPOI = dynamic_cast<TProfile*>(fListHistos->FindObject("Flow_VetaPOI_SP"));
-    TProfile* pHistProVPtRP   = dynamic_cast<TProfile*>(fListHistos->FindObject("Flow_VPtRP_SP"));
-    TProfile* pHistProVPtPOI  = dynamic_cast<TProfile*>(fListHistos->FindObject("Flow_VPtPOI_SP"));
-
-    if (pCommonHist && pCommonHistResults) {
+    TProfile* pHistProQaQb     = dynamic_cast<TProfile*>(fListHistos->FindObject("Flow_QaQb_SP"));
+    TProfile* pHistProUQetaRP  = dynamic_cast<TProfile*>(fListHistos->FindObject("Flow_UQetaRP_SP"));
+    TProfile* pHistProUQetaPOI = dynamic_cast<TProfile*>(fListHistos->FindObject("Flow_UQetaPOI_SP"));
+    TProfile* pHistProUQPtRP   = dynamic_cast<TProfile*>(fListHistos->FindObject("Flow_UQPtRP_SP"));
+    TProfile* pHistProUQPtPOI  = dynamic_cast<TProfile*>(fListHistos->FindObject("Flow_UQPtPOI_SP"));
+    if (pCommonHist && pCommonHistResults && pHistProQaQb && pHistProUQetaRP
+       && pHistProUQetaPOI && pHistProUQPtRP && pHistProUQPtPOI) {
       fSPTerm -> SetCommonHists(pCommonHist);
       fSPTerm -> SetCommonHistsRes(pCommonHistResults);
-      fSPTerm -> SetHistProVetaRP(pHistProVetaRP);
-      fSPTerm -> SetHistProVetaPOI(pHistProVetaPOI);
-      fSPTerm -> SetHistProVPtRP(pHistProVPtRP);
-      fSPTerm -> SetHistProVPtPOI(pHistProVPtPOI);
+      fSPTerm -> SetHistProQaQb(pHistProQaQb);
+      fSPTerm -> SetHistProUQetaRP(pHistProUQetaRP);
+      fSPTerm -> SetHistProUQetaPOI(pHistProUQetaPOI);
+      fSPTerm -> SetHistProUQPtRP(pHistProUQPtRP);
+      fSPTerm -> SetHistProUQPtPOI(pHistProUQPtPOI);
       fSPTerm -> Finish();
       PostData(0,fListHistos);
     }