/************************************************************************* * Copyright(c) 1998-2008, ALICE Experiment at CERN, All rights reserved. * * * * Author: The ALICE Off-line Project. * * Contributors are mentioned in the code where appropriate. * * * * Permission to use, copy, modify and distribute this software and its * * documentation strictly for non-commercial purposes is hereby granted * * without fee, provided that the above copyright notice appears in all * * copies and that both the copyright notice and this permission notice * * appear in the supporting documentation. The authors make no claims * * about the suitability of this software for any purpose. It is * * provided "as is" without express or implied warranty. * **************************************************************************/ /******************************** * estimating reference flow by * * fitting q-distribution * * * * author: Ante Bilandzic * * (abilandzic@gmail.com) * * * * based on the macro written * * by Sergei Voloshin * *******************************/ #define AliFlowAnalysisWithFittingQDistribution_cxx #include "Riostream.h" #include "AliFlowCommonConstants.h" #include "AliFlowCommonHist.h" #include "AliFlowCommonHistResults.h" #include "TChain.h" #include "TFile.h" #include "TList.h" #include "TF1.h" #include "TH2D.h" #include "TParticle.h" #include "TProfile.h" #include "AliFlowEventSimple.h" #include "AliFlowTrackSimple.h" #include "AliFlowAnalysisWithFittingQDistribution.h" class TH1; class TGraph; class TPave; class TLatex; class TMarker; class TObjArray; class TList; class TSystem; class TROOT; class AliFlowVector; class TVector; //================================================================================================================ using std::endl; using std::cout; ClassImp(AliFlowAnalysisWithFittingQDistribution) AliFlowAnalysisWithFittingQDistribution::AliFlowAnalysisWithFittingQDistribution(): fHistList(NULL), fBookOnlyBasicCCH(kTRUE), fCommonHists(NULL), fCommonHistsResults(NULL), fnBinsPhi(0), fPhiMin(0), fPhiMax(0), fPhiBinWidth(0), fnBinsPt(0), fPtMin(0), fPtMax(0), fPtBinWidth(0), fnBinsEta(0), fEtaMin(0), fEtaMax(0), fEtaBinWidth(0), fHarmonic(2), fAnalysisLabel(NULL), fMultiplicityIs(AliFlowCommonConstants::kRP), fWeightsList(NULL), fUsePhiWeights(kFALSE), fUsePtWeights(kFALSE), fUseEtaWeights(kFALSE), fUseParticleWeights(NULL), fPhiWeights(NULL), fPtWeights(NULL), fEtaWeights(NULL), fSumOfParticleWeights(NULL), fqDistribution(NULL), fqMin(0.), fqMax(100.), fqNbins(10000), fStoreqDistributionVsMult(kFALSE), fqDistributionVsMult(NULL), fMinMult(0.), fMaxMult(10000.), fnBinsMult(1000), fFittingParameters(NULL), fTreshold(5.), fvStart(0.05), fvMin(0.0), fvMax(0.25), fSigma2Start(0.75), fSigma2Min(0.5), fSigma2Max(2.5), fFinalResultIsFromSigma2Fitted(kTRUE), fPrintOnTheScreen(kTRUE), fDoFit(kTRUE), fExactNoRPs(0) { // constructor // base list to hold all output objects: fHistList = new TList(); fHistList->SetName("cobjFQD"); fHistList->SetOwner(kTRUE); // analysis label; fAnalysisLabel = new TString(); // list to hold histograms with phi, pt and eta weights: fWeightsList = new TList(); // initialize all arrays: this->InitializeArrays(); } // end of constructor //================================================================================================================ AliFlowAnalysisWithFittingQDistribution::~AliFlowAnalysisWithFittingQDistribution() { // destructor delete fHistList; } // end of destructor //================================================================================================================ void AliFlowAnalysisWithFittingQDistribution::Init() { // Access constants and book everything. //save old value and prevent histograms from being added to directory //to avoid name clashes in case multiple analaysis objects are used //in an analysis Bool_t oldHistAddStatus = TH1::AddDirectoryStatus(); TH1::AddDirectory(kFALSE); // access constants: this->AccessConstants(); // booking: this->BookCommonHistograms(); this->BookAndFillWeightsHistograms(); this->BookEverythingForDistributions(); // store fitting parameters: this->StoreFittingParameters(); // nest lists: fWeightsList->SetName("Weights"); fWeightsList->SetOwner(kTRUE); fHistList->Add(fWeightsList); // set harmonic in common control histograms (to be improved (should I do this somewhere else?)): (fCommonHists->GetHarmonic())->Fill(0.5,fHarmonic); //restore old status TH1::AddDirectory(oldHistAddStatus); } // end of void AliFlowAnalysisWithFittingQDistribution::Init() //================================================================================================================ void AliFlowAnalysisWithFittingQDistribution::Make(AliFlowEventSimple* anEvent) { // Loop over data only in this method. // a) Check all pointers used in this method; // b) Fill the common control histograms; // c) Loop over data and calculate Q-vector and sum of particle weights; // d) Fill the histogram for q-distribution and sum of particle weights. // a) Check all pointers used in this method: this->CheckPointersUsedInMake(); // b) fill the common control histograms: fCommonHists->FillControlHistograms(anEvent); // c) Loop over data and fill histogram for q-distribution: Double_t dPhi = 0.; // azimuthal angle in the laboratory frame Double_t dPt = 0.; // transverse momentum Double_t dEta = 0.; // pseudorapidity Double_t wPhi = 1.; // phi weight Double_t wPt = 1.; // pt weight Double_t wEta = 1.; // eta weight Double_t dReQ = 0.; // real part of Q-vector Double_t dImQ = 0.; // imaginary part of Q-vector Int_t nCounterNoRPs = 0; // needed only for shuffling Int_t n = fHarmonic; // shortcut for the harmonic Double_t dSumOfParticleWeights = 0.; // when particle weights are not used dSumOfParticleWeights is equal to multiplicity AliFlowTrackSimple *aftsTrack = NULL; Int_t nPrim = anEvent->NumberOfTracks(); // nPrim = total number of primary tracks if(fExactNoRPs > 0 && anEvent->GetNumberOfRPs() 0 && nCounterNoRPs>fExactNoRPs){continue;} // needed only for shuffling aftsTrack=anEvent->GetTrack(i); if(aftsTrack) { if(!(aftsTrack->InRPSelection())){continue;} // consider only tracks which are RPs nCounterNoRPs++; dPhi = aftsTrack->Phi(); dPt = aftsTrack->Pt(); dEta = aftsTrack->Eta(); if(fUsePhiWeights && fPhiWeights && fnBinsPhi) // determine phi weight for this particle: { wPhi = fPhiWeights->GetBinContent(1+(Int_t)(TMath::Floor(dPhi*fnBinsPhi/TMath::TwoPi()))); } if(fUsePtWeights && fPtWeights && fPtBinWidth) // determine pt weight for this particle: { wPt = fPtWeights->GetBinContent(1+(Int_t)(TMath::Floor((dPt-fPtMin)/fPtBinWidth))); } if(fUseEtaWeights && fEtaWeights && fEtaBinWidth) // determine eta weight for this particle: { wEta = fEtaWeights->GetBinContent(1+(Int_t)(TMath::Floor((dEta-fEtaMin)/fEtaBinWidth))); } // Calculate real and imaginary part of Q-vector and sum of particle weights for this event: // Q-vector: dReQ += wPhi*wPt*wEta*TMath::Cos(n*dPhi); dImQ += wPhi*wPt*wEta*TMath::Sin(n*dPhi); // sum of particle weights: dSumOfParticleWeights += wPhi*wPt*wEta; // when particle weights are not used this sum gives # of RPs, i.e. multiplicity } // end of if(aftsTrack) } // end of for(Int_t i=0;i 1.) { q = pow(dReQ*dReQ+dImQ*dImQ,0.5)/pow(dSumOfParticleWeights,0.5); fqDistribution->Fill(q,1.); fSumOfParticleWeights->Fill(dSumOfParticleWeights,1.); if(fStoreqDistributionVsMult) { // Multiplicity bin of an event (relevant for all histos vs M): Double_t dMultiplicityBin = 0.; if(fMultiplicityIs==AliFlowCommonConstants::kRP) { Int_t nNumberOfRPsEBE = anEvent->GetNumberOfRPs(); // number of RPs (i.e. number of reference particles) dMultiplicityBin = nNumberOfRPsEBE+0.5; } else if(fMultiplicityIs==AliFlowCommonConstants::kExternal) { Int_t nReferenceMultiplicityEBE = anEvent->GetReferenceMultiplicity(); // reference multiplicity for current event dMultiplicityBin = nReferenceMultiplicityEBE+0.5; } else if(fMultiplicityIs==AliFlowCommonConstants::kPOI) { Int_t nNumberOfPOIsEBE = anEvent->GetNumberOfPOIs(); // number of POIs (i.e. number of particles of interest) dMultiplicityBin = nNumberOfPOIsEBE+0.5; } // Fill qDist vs M: fqDistributionVsMult->Fill(dMultiplicityBin,q); } // end of if(fStoreqDistributionVsMult) } // end of if(dSumOfParticleWeights > 1.) } // end of Make() //================================================================================================================ void AliFlowAnalysisWithFittingQDistribution::Finish(Bool_t doFit) { // Calculate the final results. // a) Check all pointers used in this method; // b) Acces common constants; // c) Access fitting paremeters; // d) Do the final fit of q-distribution; // e) Fill common hist results; // f) Print on the screen the final results. // a) Check all pointers used in this method: this->CheckPointersUsedInFinish(); // b) Access common constants: this->AccessConstants(); // c) Access fitting paremeters: this->AccessFittingParameters(); // d) Do the final fit of q-distribution: if(doFit && fDoFit) { this->DoFit(kFALSE); // sigma^2 not fitted (fixed to 0.5) this->DoFit(kTRUE); // sigma^2 fitted } // e) Fill common hist results (by default fill with results obtained with sigma^2 fitted, // alternatively use a setter SetFinalResultIsFromSigma2Fitted(kFALSE)): this->FillCommonHistResults(fFinalResultIsFromSigma2Fitted); // f) Print on the screen the final results: if(fPrintOnTheScreen) this->PrintOnTheScreen(); } // end of void AliFlowAnalysisWithFittingQDistribution::Finish(Bool_t doFit) //================================================================================================================ void AliFlowAnalysisWithFittingQDistribution::CheckPointersUsedInMake() { // Check all pointers used in method Make(). if(!fCommonHists) { cout<Data(); AliFlowCommonHist *commonHist = dynamic_cast(outputListHistos->FindObject(commonHistName.Data())); if(commonHist) this->SetCommonHists(commonHist); TString commonHistResName = "AliFlowCommonHistResultsFQD"; commonHistResName += fAnalysisLabel->Data(); AliFlowCommonHistResults *commonHistRes = dynamic_cast (outputListHistos->FindObject(commonHistResName.Data())); if(commonHistRes) this->SetCommonHistsResults(commonHistRes); // 2.) weights: TList *weightsList = dynamic_cast(outputListHistos->FindObject("Weights")); if(weightsList){this->SetWeightsList(weightsList);} Bool_t bUsePhiWeights = kFALSE; Bool_t bUsePtWeights = kFALSE; Bool_t bUseEtaWeights = kFALSE; TString fUseParticleWeightsName = "fUseParticleWeightsFQD"; fUseParticleWeightsName += fAnalysisLabel->Data(); TProfile *useParticleWeights = NULL; if(weightsList) { useParticleWeights = dynamic_cast(weightsList->FindObject(fUseParticleWeightsName.Data())); } if(useParticleWeights) { this->SetUseParticleWeights(useParticleWeights); bUsePhiWeights = (Int_t)useParticleWeights->GetBinContent(1); this->SetUsePhiWeights(bUsePhiWeights); bUsePtWeights = (Int_t)useParticleWeights->GetBinContent(2); this->SetUsePtWeights(bUsePtWeights); bUseEtaWeights = (Int_t)useParticleWeights->GetBinContent(3); this->SetUseEtaWeights(bUseEtaWeights); } // 3.) distributions and 4.) final results of fitting: TString sigmaFlag[2] = {"#sigma^{2} not fitted","#sigma^{2} fitted"}; // q-distribution: TString qDistributionName = "fqDistribution"; qDistributionName += fAnalysisLabel->Data(); // sum of particle weights: TString sumOfParticleWeightsName = "fSumOfParticleWeights"; sumOfParticleWeightsName += fAnalysisLabel->Data(); // final results for reference flow: TString intFlowName = "fIntFlowFQD"; intFlowName += fAnalysisLabel->Data(); // sigma^2: TString sigma2Name = "fSigma2"; sigma2Name += fAnalysisLabel->Data(); // chi^2: TString chi2Name = "fChi2"; chi2Name += fAnalysisLabel->Data(); // fitting function: TString fittingFunctionName = "fFittingFunction"; fittingFunctionName += fAnalysisLabel->Data(); TH1D *qDistribution = NULL; TH1D *sumOfParticleWeights = NULL; TH1D *intFlow[2] = {NULL}; TH1D *sigma2[2] = {NULL}; TH1D *chi2[2] = {NULL}; TF1 *fittingFunction[2] = {NULL}; // q-distribution: qDistribution = dynamic_cast(outputListHistos->FindObject(Form("%s",qDistributionName.Data()))); if(qDistribution) { this->SetqDistribution(qDistribution); } else { cout<<"WARNING: qDistribution is NULL in AFAWFQD::GOH() !!!!"<(outputListHistos->FindObject(Form("%s",sumOfParticleWeightsName.Data()))); if(sumOfParticleWeights) { this->SetSumOfParticleWeights(sumOfParticleWeights); } else { cout<<"WARNING: sumOfParticleWeights is NULL in AFAWFQD::GOH() !!!!"<(outputListHistos->FindObject(Form("%s, %s",intFlowName.Data(),sigmaFlag[f].Data()))); if(intFlow[f]) { this->SetIntFlow(intFlow[f],f); } else { cout<<"WARNING: intFlow[f] is NULL in AFAWFQD::GOH() !!!!"<SetSigma2(sigma2[f],f); } else { cout<<"WARNING: sigma2[f] is NULL in AFAWFQD::GOH() !!!!"<SetChi2(chi2[f],f); } else { cout<<"WARNING: chi2[f] is NULL in AFAWFQD::GOH() !!!!"<SetFittingFunction(fittingFunction[f],f); } else { cout<<"WARNING: fittingFunction[f] is NULL in AFAWFQD::GOH() !!!!"<Data(); TProfile *fittingParameters = NULL; fittingParameters = dynamic_cast(outputListHistos->FindObject(fittingParametersName.Data())); if(fittingParameters) { this->SetFittingParameters(fittingParameters); } else { cout<<"WARNING:fittingParameters is NULL in AFAWFQD::GOH() !!!!"<Data(),"RECREATE"); //output->WriteObject(fHistList, "cobjFQD","SingleKey"); fHistList->SetName("cobjFQD"); fHistList->SetOwner(kTRUE); fHistList->Write(fHistList->GetName(), TObject::kSingleKey); delete output; } //================================================================================================================ void AliFlowAnalysisWithFittingQDistribution::WriteHistograms(TString outputFileName) { //store the final results in output .root file TFile *output = new TFile(outputFileName.Data(),"RECREATE"); //output->WriteObject(fHistList, "cobjFQD","SingleKey"); fHistList->SetName("cobjFQD"); fHistList->SetOwner(kTRUE); fHistList->Write(fHistList->GetName(), TObject::kSingleKey); delete output; } //================================================================================================================ void AliFlowAnalysisWithFittingQDistribution::WriteHistograms(TDirectoryFile *outputFileName) { //store the final results in output .root file fHistList->SetName("cobjFQD"); fHistList->SetOwner(kTRUE); outputFileName->Add(fHistList); outputFileName->Write(outputFileName->GetName(), TObject::kSingleKey); } //================================================================================================================ void AliFlowAnalysisWithFittingQDistribution::InitializeArrays() { // Initialize all arrays. for(Int_t s2F=0;s2F<2;s2F++) // sigma^2 not fitted (0) or fitted (1) { fIntFlow[s2F] = NULL; fSigma2[s2F] = NULL; fChi2[s2F] = NULL; fFittingFunction[s2F] = NULL; } } // end of void AliFlowAnalysisWithFittingQDistribution::InitializeArrays() //================================================================================================================ void AliFlowAnalysisWithFittingQDistribution::BookCommonHistograms() { // Book common histograms. // common control histogram: TString commonHistName = "AliFlowCommonHistFQD"; commonHistName += fAnalysisLabel->Data(); fCommonHists = new AliFlowCommonHist(commonHistName.Data(),commonHistName.Data(),fBookOnlyBasicCCH); fHistList->Add(fCommonHists); // common histograms for final results: TString commonHistResName = "AliFlowCommonHistResultsFQD"; commonHistResName += fAnalysisLabel->Data(); fCommonHistsResults = new AliFlowCommonHistResults(commonHistResName.Data(),"",fHarmonic); fHistList->Add(fCommonHistsResults); } // end of void AliFlowAnalysisWithFittingQDistribution::BookCommonHistograms() //================================================================================================================ void AliFlowAnalysisWithFittingQDistribution::BookAndFillWeightsHistograms() { // Book and fill histograms which hold phi, pt and eta weights. if(!fWeightsList) { cout<<"WARNING: fWeightsList is NULL in AFAWFQD::BAFWH() !!!!"<Data(); fUseParticleWeights = new TProfile(fUseParticleWeightsName.Data(),"0 = particle weight not used, 1 = particle weight used ",3,0,3); fUseParticleWeights->SetLabelSize(0.08); (fUseParticleWeights->GetXaxis())->SetBinLabel(1,"w_{#phi}"); (fUseParticleWeights->GetXaxis())->SetBinLabel(2,"w_{p_{T}}"); (fUseParticleWeights->GetXaxis())->SetBinLabel(3,"w_{#eta}"); fUseParticleWeights->Fill(0.5,(Int_t)fUsePhiWeights); fUseParticleWeights->Fill(1.5,(Int_t)fUsePtWeights); fUseParticleWeights->Fill(2.5,(Int_t)fUseEtaWeights); fWeightsList->Add(fUseParticleWeights); // phi weights: if(fUsePhiWeights) { if(fWeightsList->FindObject("phi_weights")) { fPhiWeights = dynamic_cast(fWeightsList->FindObject("phi_weights")); if(!fPhiWeights){printf("\n WARNING (FQD): !fPhiWeights !!!!\n");exit(0);} if(TMath::Abs(fPhiWeights->GetBinWidth(1)-fPhiBinWidth)>pow(10.,-6.)) { cout<FindObject(\"phi_weights\") is NULL in AFAWFQD::BAFWH() !!!!"<FindObject("pt_weights")) { fPtWeights = dynamic_cast(fWeightsList->FindObject("pt_weights")); if(!fPtWeights){printf("\n WARNING (FQD): !fPtWeights !!!!\n");exit(0);} if(TMath::Abs(fPtWeights->GetBinWidth(1)-fPtBinWidth)>pow(10.,-6.)) { cout<FindObject(\"pt_weights\") is NULL in AFAWFQD::BAFWH() !!!!"<FindObject("eta_weights")) { fEtaWeights = dynamic_cast(fWeightsList->FindObject("eta_weights")); if(!fEtaWeights){printf("\n WARNING (FQD): !fEtaWeights !!!!\n");exit(0);} if(TMath::Abs(fEtaWeights->GetBinWidth(1)-fEtaBinWidth)>pow(10.,-6.)) { cout<FindObject(\"eta_weights\") is NULL in AFAWFQD::BAFWH() !!!!"<GetNbinsPhi(); fPhiMin = AliFlowCommonConstants::GetMaster()->GetPhiMin(); fPhiMax = AliFlowCommonConstants::GetMaster()->GetPhiMax(); if(fnBinsPhi) fPhiBinWidth = (fPhiMax-fPhiMin)/fnBinsPhi; fnBinsPt = AliFlowCommonConstants::GetMaster()->GetNbinsPt(); fPtMin = AliFlowCommonConstants::GetMaster()->GetPtMin(); fPtMax = AliFlowCommonConstants::GetMaster()->GetPtMax(); if(fnBinsPt) fPtBinWidth = (fPtMax-fPtMin)/fnBinsPt; fnBinsEta = AliFlowCommonConstants::GetMaster()->GetNbinsEta(); fEtaMin = AliFlowCommonConstants::GetMaster()->GetEtaMin(); fEtaMax = AliFlowCommonConstants::GetMaster()->GetEtaMax(); if(fnBinsEta) fEtaBinWidth = (fEtaMax-fEtaMin)/fnBinsEta; } // end of void AliFlowAnalysisWithFittingQDistribution::AccessConstants() //================================================================================================================================ void AliFlowAnalysisWithFittingQDistribution::BookEverythingForDistributions() { // Book all objects relevant for fitting of q-distributions. // Flags: TString sigmaFlag[2] = {"#sigma^{2} not fitted","#sigma^{2} fitted"}; // q-distribution: TString fqDistributionName = "fqDistribution"; fqDistributionName += fAnalysisLabel->Data(); fqDistribution = new TH1D(Form("%s",fqDistributionName.Data()),"q-distribution",fqNbins,fqMin,fqMax); fqDistribution->SetXTitle(Form("q_{%d}=|Q_{%d}|/#sqrt{M}",fHarmonic,fHarmonic)); fqDistribution->SetYTitle("Counts"); fHistList->Add(fqDistribution); // q-distribution vs multiplicity: if(fStoreqDistributionVsMult) { TString fqDistributionVsMultName = "fqDistributionVsMult"; fqDistributionVsMultName += fAnalysisLabel->Data(); fqDistributionVsMult = new TH2D(Form("%s",fqDistributionVsMultName.Data()),"q-distribution vs M",fnBinsMult,fMinMult,fMaxMult,fqNbins,fqMin,fqMax); fqDistributionVsMult->GetYaxis()->SetTitle(Form("q_{%d}=|Q_{%d}|/#sqrt{M}",fHarmonic,fHarmonic)); if(fMultiplicityIs==AliFlowCommonConstants::kRP) { fqDistributionVsMult->GetXaxis()->SetTitle("# RPs"); } else if(fMultiplicityIs==AliFlowCommonConstants::kExternal) { fqDistributionVsMult->GetXaxis()->SetTitle("Reference multiplicity (from ESD)"); } else if(fMultiplicityIs==AliFlowCommonConstants::kPOI) { fqDistributionVsMult->GetXaxis()->SetTitle("# POIs"); } fHistList->Add(fqDistributionVsMult); } // end of if(fStoreqDistributionVsMult) // Sum of particle weights: TString fSumOfParticleWeightsName = "fSumOfParticleWeights"; fSumOfParticleWeightsName += fAnalysisLabel->Data(); fSumOfParticleWeights = new TH1D(Form("%s",fSumOfParticleWeightsName.Data()),"Sum of particle weights",10000,0,100000); // (to be improved - harwired limits and number of bins) fSumOfParticleWeights->SetXTitle("#sum_{i=1}^{N} w_{i}"); fSumOfParticleWeights->GetXaxis()->SetTitleSize(0.03); fSumOfParticleWeights->SetYTitle("Counts"); fHistList->Add(fSumOfParticleWeights); // Final results for integrated flow: TString fIntFlowName = "fIntFlowFQD"; fIntFlowName += fAnalysisLabel->Data(); // sigma^2: TString fSigma2Name = "fSigma2"; fSigma2Name += fAnalysisLabel->Data(); // chi^2: TString fChi2Name = "fChi2"; fChi2Name += fAnalysisLabel->Data(); // fitting function: TString fittingFunctionName = "fFittingFunction"; fittingFunctionName += fAnalysisLabel->Data(); for(Int_t f=0;f<2;f++) // sigma^2 not fitted (0) or fitted (1) { // final results for integrated flow: fIntFlow[f] = new TH1D(Form("%s, %s",fIntFlowName.Data(),sigmaFlag[f].Data()),Form("Reference Flow (%s)",sigmaFlag[f].Data()),1,0,1); fIntFlow[f]->SetLabelSize(0.08); fHistList->Add(fIntFlow[f]); // sigma^2: fSigma2[f] = new TH1D(Form("%s, %s",fSigma2Name.Data(),sigmaFlag[f].Data()),Form("#sigma^{2} (%s)",sigmaFlag[f].Data()),1,0,1); fSigma2[f]->SetLabelSize(0.08); (fSigma2[f]->GetXaxis())->SetBinLabel(1,"#sigma^{2}"); fHistList->Add(fSigma2[f]); // chi^2: fChi2[f] = new TH1D(Form("%s, %s",fChi2Name.Data(),sigmaFlag[f].Data()),Form("#chi^{2} (%s)",sigmaFlag[f].Data()),1,0,1); fChi2[f]->SetLabelSize(0.08); (fChi2[f]->GetXaxis())->SetLabelOffset(0.01); (fChi2[f]->GetXaxis())->SetBinLabel(1,"#chi^{2}"); fHistList->Add(fChi2[f]); // fitting functions: fFittingFunction[f] = new TF1(Form("%s, %s",fittingFunctionName.Data(),sigmaFlag[f].Data()),"[2]*(x/[1])*exp(-(x*x+[0]*[0])/(2.*[1]))*TMath::BesselI0(x*[0]/[1])"); fHistList->Add(fFittingFunction[f]); } // end of for(Int_t f=0;f<2;f++) // sigma^2 not fitted or fitted // Book profile fFittingParameters which will hold all fitting parameters: TString fFittingParametersName = "fFittingParameters"; fFittingParametersName += fAnalysisLabel->Data(); fFittingParameters = new TProfile(fFittingParametersName.Data(),"Parameters for fitting q-distribution",13,0,13); fFittingParameters->SetLabelSize(0.05); fFittingParameters->GetXaxis()->SetBinLabel(1,"treshold"); fFittingParameters->GetXaxis()->SetBinLabel(2,Form("starting v_{%d}",fHarmonic)); fFittingParameters->GetXaxis()->SetBinLabel(3,Form("min. v_{%d}",fHarmonic)); fFittingParameters->GetXaxis()->SetBinLabel(4,Form("max. v_{%d}",fHarmonic)); fFittingParameters->GetXaxis()->SetBinLabel(5,"starting #sigma^{2}"); fFittingParameters->GetXaxis()->SetBinLabel(6,"min. #sigma^{2}"); fFittingParameters->GetXaxis()->SetBinLabel(7,"max. #sigma^{2}"); fFittingParameters->GetXaxis()->SetBinLabel(8,"Final result is from #sigma^{2} fitted?"); fFittingParameters->GetXaxis()->SetBinLabel(9,"Print results on the screen?"); fFittingParameters->GetXaxis()->SetBinLabel(10,"fStoreqDistributionVsMult"); fFittingParameters->GetXaxis()->SetBinLabel(11,"fMultiplicityIs"); fFittingParameters->GetXaxis()->SetBinLabel(12,"fDoFit"); fFittingParameters->GetXaxis()->SetBinLabel(13,"fExactNoRPs"); fHistList->Add(fFittingParameters); } // end of void AliFlowAnalysisWithFittingQDistribution::BookEverythingForDistributions() //================================================================================================================================ void AliFlowAnalysisWithFittingQDistribution::DoFit(Bool_t sigma2Fitted) { // Do the final fit of q-distribution. Int_t s2F = (Int_t)(sigma2Fitted); // shortcut Double_t AvM = fSumOfParticleWeights->GetMean(1); // average multiplicity //Int_t nEvts = (Int_t)fSumOfParticleWeights->GetEntries(); // number of events: // Start fitting from the bin with at least fTreshold entries, // finish fitting at the bin with at least fTreshold entries: Int_t binMin = fqDistribution->FindFirstBinAbove(fTreshold); Int_t binMax = fqDistribution->FindLastBinAbove(fTreshold); Double_t binWidth = fqDistribution->GetBinWidth(4); // assuming that all bins have the same width if(TMath::Abs(binWidth) < 1.e-44) { cout<GetBinContent(b); } Double_t norm = binWidth*ent; // norm (assuming that all bins have the same width) // Fitting function: fFittingFunction[s2F]->SetRange(qmin,qmax); fFittingFunction[s2F]->SetParNames("v*sqrt{sum of particle weights}","sigma^2","norm"); fFittingFunction[s2F]->SetParameters(fvStart*pow(AvM,0.5),fSigma2Start,norm); fFittingFunction[s2F]->SetParLimits(0,fvMin*pow(AvM,0.5),fvMax*pow(AvM,0.5)); if(s2F == 0) { fFittingFunction[s2F]->FixParameter(1,0.5); } else { fFittingFunction[s2F]->SetParLimits(1,fSigma2Min,fSigma2Max); } fFittingFunction[s2F]->FixParameter(2,norm); // Fitting is done here: // Remark: do fit only if # of entries > 50 - this is only a pragmatics fix to avoid TMinuit crash (to be improved) if(ent > 50) { fqDistribution->Fit(fFittingFunction[s2F]->GetName(),"NQ","",qmin,qmax); } // Final results: Double_t v = 0.; // reference flow Double_t vError = 0.; // error of reference flow Double_t sigma2 = 0.; // sigma^2 Double_t sigma2Error = 0.; // error of sigma^2 Double_t chi2 = 0; // chi^2 from Minuit // Reference flow: if(AvM) { v = fFittingFunction[s2F]->GetParameter(0)/pow(AvM,0.5); vError = fFittingFunction[s2F]->GetParError(0)/pow(AvM,0.5); fIntFlow[s2F]->SetBinContent(1,v); // s2F is shortcut for "sigma^2 fitted" fIntFlow[s2F]->SetBinError(1,vError); // s2F is shortcut for "sigma^2 fitted" } else { cout<SetBinContent(1,sigma2); fSigma2[0]->SetBinError(1,0.); } else // sigma^2 fitted { sigma2 = fFittingFunction[s2F]->GetParameter(1); sigma2Error = fFittingFunction[s2F]->GetParError(1); fSigma2[1]->SetBinContent(1,sigma2); fSigma2[1]->SetBinError(1,sigma2Error); } // chi^2: chi2 = fFittingFunction[s2F]->GetChisquare(); fChi2[s2F]->SetBinContent(1,chi2); } // end of void AliFlowAnalysisWithFittingQDistribution::DoFit() //================================================================================================================================ void AliFlowAnalysisWithFittingQDistribution::FillCommonHistResults(Bool_t sigma2Fitted) { // Fill common result histogram for reference flow and resolution. // Remark: by default the result obtained with sigma^2 fitted is being stored. // In order to store the result obtained with sigma^2 fixed use a setter SetFinalResultIsFromSigma2Fitted(kFALSE). Int_t s2F = (Int_t)sigma2Fitted; // shortcut // Reference flow: Double_t v = fIntFlow[s2F]->GetBinContent(1); Double_t vError = fIntFlow[s2F]->GetBinError(1); fCommonHistsResults->FillIntegratedFlow(v,vError); // Resolution: Double_t AvM = fSumOfParticleWeights->GetMean(1); Double_t chi2 = AvM*pow(v,2.); // chi^2 if(chi2>=0.) { fCommonHistsResults->FillChi(pow(chi2,0.5)); } } // end of void AliFlowAnalysisWithFittingQDistribution::FillCommonHistResultsIntFlow(Bool_t useParticleWeights, Bool_t sigma2NotFixed) //================================================================================================================================ void AliFlowAnalysisWithFittingQDistribution::PrintOnTheScreen() { // Print the final results on the screen. // shortcut for the harmonic: Int_t n = 0; if(fCommonHists->GetHarmonic()) { n = (Int_t)(fCommonHists->GetHarmonic())->GetBinContent(1); } // printing: cout<GetBinError(1)<GetBinContent(1)-fvMin)<1.e-10 || TMath::Abs(fIntFlow[0]->GetBinContent(1)-fvMax)<1.e-10) { cout<<" WARNING: value of v_"<GetBinError(1)<GetBinError(1)<GetBinContent(1)-fvMin)<1.e-10 || TMath::Abs(fIntFlow[0]->GetBinContent(1)-fvMax)<1.e-10) { cout<<" WARNING: value of v_"<GetBinContent(1)-fSigma2Min)<1.e-10 || TMath::Abs(fSigma2[1]->GetBinContent(1)-fSigma2Max)<1.e-10) { cout<<" WARNING: value of sigma^2 is on the boundary"<Reset(); fFittingParameters->Fill(0.5,fTreshold); fFittingParameters->Fill(1.5,fvStart); fFittingParameters->Fill(2.5,fvMin); fFittingParameters->Fill(3.5,fvMax); fFittingParameters->Fill(4.5,fSigma2Start); fFittingParameters->Fill(5.5,fSigma2Min); fFittingParameters->Fill(6.5,fSigma2Max); fFittingParameters->Fill(7.5,fFinalResultIsFromSigma2Fitted); fFittingParameters->Fill(8.5,fPrintOnTheScreen); fFittingParameters->Fill(9.5,fStoreqDistributionVsMult); // which multiplicity was used: if(fMultiplicityIs==AliFlowCommonConstants::kRP) // # of Reference Particles { fFittingParameters->Fill(10.5,0); // 0 = # of Reference Particles } else if(fMultiplicityIs==AliFlowCommonConstants::kExternal) { fFittingParameters->Fill(10.5,1); // 1 = ref. mult. from ESD } else if(fMultiplicityIs==AliFlowCommonConstants::kPOI) { fFittingParameters->Fill(10.5,2); // 2 = # of Particles of Interest } fFittingParameters->Fill(11.5,fDoFit); fFittingParameters->Fill(12.5,fExactNoRPs); } // end of void AliFlowAnalysisWithFittingQDistribution::StoreFittingParameters() //================================================================================================================================ void AliFlowAnalysisWithFittingQDistribution::AccessFittingParameters() { // Access fitting parameters for the fit of q-distribution. fTreshold = fFittingParameters->GetBinContent(1); fvStart = fFittingParameters->GetBinContent(2); fvMin = fFittingParameters->GetBinContent(3); fvMax = fFittingParameters->GetBinContent(4); fSigma2Start = fFittingParameters->GetBinContent(5); fSigma2Min = fFittingParameters->GetBinContent(6); fSigma2Max = fFittingParameters->GetBinContent(7); fFinalResultIsFromSigma2Fitted = (Bool_t)fFittingParameters->GetBinContent(8); fPrintOnTheScreen = (Bool_t)fFittingParameters->GetBinContent(9); fStoreqDistributionVsMult = (Bool_t)fFittingParameters->GetBinContent(10); fDoFit = (Bool_t)fFittingParameters->GetBinContent(12); fExactNoRPs = (Int_t)fFittingParameters->GetBinContent(13); } // end of void AliFlowAnalysisWithFittingQDistribution::AccessFittingParameters()