X-Git-Url: http://git.uio.no/git/?a=blobdiff_plain;f=PWG0%2FAliPWG0Helper.cxx;h=66ff63255ed5fa4e639a7392432bfc93e2a09e49;hb=afabc52f1bc55966b235138435c07b5040b9ccff;hp=13e01024122930435e726d6d90e8a0387dca93e0;hpb=29771dc84287fe10feadbd49f5235ed82ed0b968;p=u%2Fmrichter%2FAliRoot.git diff --git a/PWG0/AliPWG0Helper.cxx b/PWG0/AliPWG0Helper.cxx index 13e01024122..66ff63255ed 100644 --- a/PWG0/AliPWG0Helper.cxx +++ b/PWG0/AliPWG0Helper.cxx @@ -5,58 +5,128 @@ #include #include #include +#include #include +#include +#include +#include +#include + +#include +#include #include #include +#include #include +#include +#include + +#include +#include +#include +#include +#include //____________________________________________________________________ ClassImp(AliPWG0Helper) +Int_t AliPWG0Helper::fgLastProcessType = -1; + //____________________________________________________________________ -Bool_t AliPWG0Helper::IsEventTriggered(AliESD* aEsd) +Bool_t AliPWG0Helper::TestVertex(const AliESDVertex* vertex, AnalysisMode analysisMode, Bool_t debug) { - // check if the event was triggered - // - // this function needs the branch fTriggerMask - // - // MB should be - // ITS_SPD_GFO_L0 : 32 - // VZERO_OR_LEFT : 1 - // VZERO_OR_RIGHT : 2 + // Checks if a vertex meets the needed quality criteria - ULong64_t triggerMask = aEsd->GetTriggerMask(); + Float_t requiredZResolution = -1; + if (analysisMode & kSPD || analysisMode & kTPCITS) + { + requiredZResolution = 0.1; + } + else if (analysisMode & kTPC) + requiredZResolution = 10.; - if (triggerMask&32 && ((triggerMask&1) || (triggerMask&2))) - return kTRUE; + // check resolution + Double_t zRes = vertex->GetZRes(); + + if (zRes > requiredZResolution) { + if (debug) + Printf("AliPWG0Helper::TestVertex: Resolution too poor %f (required: %f", zRes, requiredZResolution); + return kFALSE; + } - return kFALSE; + return kTRUE; } //____________________________________________________________________ -Bool_t AliPWG0Helper::IsVertexReconstructed(AliESD* aEsd) +const AliESDVertex* AliPWG0Helper::GetVertex(AliESDEvent* aEsd, AnalysisMode analysisMode, Bool_t debug, Bool_t bRedoTPC) { - // checks if the vertex is reasonable + // Get the vertex from the ESD and returns it if the vertex is valid // - // this function needs the branches fSPDVertex* + // Second argument decides which vertex is used (this selects + // also the quality criteria that are applied) + const AliESDVertex* vertex = 0; + if (analysisMode & kSPD || analysisMode & kTPCITS) + { + vertex = aEsd->GetPrimaryVertexSPD(); + if (debug) + Printf("AliPWG0Helper::GetVertex: Returning SPD vertex"); + } + else if (analysisMode & kTPC) + { + if(bRedoTPC){ + if (debug) + Printf("AliPWG0Helper::GetVertex: Redoing vertex"); + Double_t kBz = aEsd->GetMagneticField(); + AliVertexerTracks vertexer(kBz); + vertexer.SetTPCMode(); + AliESDVertex *vTPC = vertexer.FindPrimaryVertex(aEsd); + aEsd->SetPrimaryVertexTPC(vTPC); + for (Int_t i=0; iGetNumberOfTracks(); i++) { + AliESDtrack *t = aEsd->GetTrack(i); + t->RelateToVertexTPC(vTPC, kBz, kVeryBig); + } + delete vTPC; + } - const AliESDVertex* vtxESD = aEsd->GetVertex(); + vertex = aEsd->GetPrimaryVertexTPC(); + if (debug) + Printf("AliPWG0Helper::GetVertex: Returning vertex from tracks"); + } + else + Printf("AliPWG0Helper::GetVertex: ERROR: Invalid second argument %d", analysisMode); - // the vertex should be reconstructed - if (strcmp(vtxESD->GetName(), "default")==0) - return kFALSE; + if (!vertex) { + if (debug) + Printf("AliPWG0Helper::GetVertex: No vertex found in ESD"); + return 0; + } - Double_t vtx_res[3]; - vtx_res[0] = vtxESD->GetXRes(); - vtx_res[1] = vtxESD->GetYRes(); - vtx_res[2] = vtxESD->GetZRes(); + // check Ncontributors + if (vertex->GetNContributors() <= 0) { + if (debug){ + Printf("AliPWG0Helper::GetVertex: NContributors() <= 0: %d",vertex->GetNContributors()); + Printf("AliPWG0Helper::GetVertex: NIndices(): %d",vertex->GetNIndices()); + vertex->Print(); + } + return 0; + } - if (vtx_res[2]==0 || vtx_res[2]>0.1) - return kFALSE; + // check resolution + Double_t zRes = vertex->GetZRes(); + if (zRes == 0) { + Printf("AliPWG0Helper::GetVertex: UNEXPECTED: resolution is 0."); + return 0; + } - return kTRUE; + if (debug) + { + Printf("AliPWG0Helper::GetVertex: Returning valid vertex: %s", vertex->GetTitle()); + vertex->Print(); + } + + return vertex; } //____________________________________________________________________ @@ -68,6 +138,8 @@ Bool_t AliPWG0Helper::IsPrimaryCharged(TParticle* aParticle, Int_t aTotalPrimari // // This function or a equivalent should be available in some common place of AliRoot // + // WARNING: Call this function only for particles that are among the particles from the event generator! + // --> stack->Particle(id) with id < stack->GetNprimary() // if the particle has a daughter primary, we do not want to count it if (aParticle->GetFirstDaughter() != -1 && aParticle->GetFirstDaughter() < aTotalPrimaries) @@ -78,6 +150,7 @@ Bool_t AliPWG0Helper::IsPrimaryCharged(TParticle* aParticle, Int_t aTotalPrimari } Int_t pdgCode = TMath::Abs(aParticle->GetPdgCode()); + // skip quarks and gluon if (pdgCode <= 10 || pdgCode == 21) @@ -87,10 +160,17 @@ Bool_t AliPWG0Helper::IsPrimaryCharged(TParticle* aParticle, Int_t aTotalPrimari return kFALSE; } + Int_t status = aParticle->GetStatusCode(); + // skip non final state particles.. + if(status!=1){ + if (adebug) + printf("Dropping particle because it is not a final state particle.\n"); + return kFALSE; + } + if (strcmp(aParticle->GetName(),"XXX") == 0) { - if (adebug) - printf("WARNING: There is a particle named XXX.\n"); + Printf("WARNING: There is a particle named XXX (pdg code %d).", pdgCode); return kFALSE; } @@ -98,8 +178,7 @@ Bool_t AliPWG0Helper::IsPrimaryCharged(TParticle* aParticle, Int_t aTotalPrimari if (strcmp(pdgPart->ParticleClass(),"Unknown") == 0) { - if (adebug) - printf("WARNING: There is a particle with an unknown particle class (pdg code %d).\n", pdgCode); + Printf("WARNING: There is a particle with an unknown particle class (pdg code %d).", pdgCode); return kFALSE; } @@ -176,7 +255,7 @@ void AliPWG0Helper::CreateDividedProjections(TH3* hist, TH3* hist2, const char* TH1* division = dynamic_cast (proj->Clone(Form("%s_div_%s", proj->GetName(), proj2->GetName()))); //printf("doing axis: %s, x axis has %d %d bins, min %f %f max %f %f\n", axis, division->GetNbinsX(), proj2->GetNbinsX(), division->GetXaxis()->GetBinLowEdge(1), proj2->GetXaxis()->GetBinLowEdge(1), division->GetXaxis()->GetBinUpEdge(division->GetNbinsX()), proj2->GetXaxis()->GetBinUpEdge(proj2->GetNbinsX())); //printf("doing axis: %s, y axis has %d %d bins, min %f %f max %f %f\n", axis, division->GetNbinsY(), proj2->GetNbinsY(), division->GetYaxis()->GetBinLowEdge(1), proj2->GetYaxis()->GetBinLowEdge(1), division->GetYaxis()->GetBinUpEdge(division->GetNbinsY()), proj2->GetYaxis()->GetBinUpEdge(proj2->GetNbinsY())); - division->Divide(proj2); + division->Divide(proj, proj2, 1, 1, "B"); division->SetTitle(Form("%s divided %s", proj->GetTitle(), proj2->GetTitle())); if (putErrors) @@ -222,3 +301,237 @@ const char* AliPWG0Helper::GetAxisTitle(TH3* hist, const char axis) return 0; } + + +AliPWG0Helper::MCProcessType AliPWG0Helper::GetPythiaEventProcessType(AliGenEventHeader* aHeader, Bool_t adebug) { + + AliGenPythiaEventHeader* pythiaGenHeader = dynamic_cast(aHeader); + + if (!pythiaGenHeader) { + printf("AliPWG0Helper::GetProcessType : Unknown gen Header type). \n"); + return kInvalidProcess; + } + + + Int_t pythiaType = pythiaGenHeader->ProcessType(); + fgLastProcessType = pythiaType; + MCProcessType globalType = kInvalidProcess; + + + if (adebug) { + printf("AliPWG0Helper::GetProcessType : Pythia process type found: %d \n",pythiaType); + } + + + if(pythiaType==92||pythiaType==93){ + globalType = kSD; + } + else if(pythiaType==94){ + globalType = kDD; + } + //else if(pythiaType != 91){ // also exclude elastic to be sure... CKB?? + else { + globalType = kND; + } + return globalType; +} + + +AliPWG0Helper::MCProcessType AliPWG0Helper::GetDPMjetEventProcessType(AliGenEventHeader* aHeader, Bool_t adebug) { + // + // get the process type of the event. + // + + // can only read pythia headers, either directly or from cocktalil header + AliGenDPMjetEventHeader* dpmJetGenHeader = dynamic_cast(aHeader); + + if (!dpmJetGenHeader) { + printf("AliPWG0Helper::GetDPMjetProcessType : Unknown header type (not DPMjet or). \n"); + return kInvalidProcess; + } + + Int_t dpmJetType = dpmJetGenHeader->ProcessType(); + fgLastProcessType = dpmJetType; + MCProcessType globalType = kInvalidProcess; + + + if (adebug) { + printf("AliPWG0Helper::GetDPMJetProcessType : DPMJet process type found: %d \n",dpmJetType); + } + + + if(dpmJetType == 1){ // this is explicitly inelastic + globalType = kND; + } + else if(dpmJetType==5||dpmJetType==6){ + globalType = kSD; + } + else if(dpmJetType==7||dpmJetType==4){// DD and double pomeron + globalType = kDD; + } + return globalType; +} + + +AliPWG0Helper::MCProcessType AliPWG0Helper::GetEventProcessType(AliHeader* aHeader, Bool_t adebug) { + // + // get the process type of the event. + // + + + // Check for simple headers first + + AliGenPythiaEventHeader* pythiaGenHeader = dynamic_cast(aHeader->GenEventHeader()); + if (pythiaGenHeader) { + return GetPythiaEventProcessType(pythiaGenHeader,adebug); + } + + AliGenDPMjetEventHeader* dpmJetGenHeader = dynamic_cast(aHeader->GenEventHeader()); + if (dpmJetGenHeader) { + return GetDPMjetEventProcessType(dpmJetGenHeader,adebug); + } + + + // check for cocktail + + AliGenCocktailEventHeader* genCocktailHeader = dynamic_cast(aHeader->GenEventHeader()); + if (!genCocktailHeader) { + printf("AliPWG0Helper::GetProcessType : Unknown header type (not Pythia or Cocktail). \n"); + return kInvalidProcess; + } + + TList* headerList = genCocktailHeader->GetHeaders(); + if (!headerList) { + return kInvalidProcess; + } + + for (Int_t i=0; iGetEntries(); i++) { + + pythiaGenHeader = dynamic_cast(headerList->At(i)); + if (pythiaGenHeader) { + return GetPythiaEventProcessType(pythiaGenHeader,adebug); + } + + dpmJetGenHeader = dynamic_cast(headerList->At(i)); + if (dpmJetGenHeader) { + return GetDPMjetEventProcessType(dpmJetGenHeader,adebug); + } + } + return kInvalidProcess; +} + + + +//____________________________________________________________________ +TParticle* AliPWG0Helper::FindPrimaryMother(AliStack* stack, Int_t label) +{ + // + // Finds the first mother among the primary particles of the particle identified by