3 #include "AliMultiplicityTask.h"
14 #include <TParticle.h>
17 #include <TObjString.h>
21 #include <AliESDVertex.h>
22 #include <AliESDEvent.h>
24 #include <AliHeader.h>
25 #include <AliGenEventHeader.h>
26 #include <AliMultiplicity.h>
27 #include <AliAnalysisManager.h>
28 #include <AliMCEventHandler.h>
29 #include <AliMCEvent.h>
30 #include <AliESDInputHandler.h>
32 #include "AliESDtrackCuts.h"
33 #include "AliPWG0Helper.h"
34 #include "multiplicity/AliMultiplicityCorrection.h"
35 #include "AliCorrection.h"
36 #include "AliCorrectionMatrix3D.h"
37 #include "AliPhysicsSelection.h"
38 #include "AliTriggerAnalysis.h"
40 ClassImp(AliMultiplicityTask)
42 AliMultiplicityTask::AliMultiplicityTask(const char* opt) :
43 AliAnalysisTask("AliMultiplicityTask", ""),
46 fAnalysisMode((AliPWG0Helper::AnalysisMode) (AliPWG0Helper::kSPD | AliPWG0Helper::kFieldOn)),
47 fTrigger(AliTriggerAnalysis::kMB1),
49 fDiffTreatment(AliPWG0Helper::kMCFlags),
54 fSystSkipParticles(kFALSE),
55 fSelectProcessType(0),
66 // Constructor. Initialization of pointers
69 for (Int_t i = 0; i<8; i++)
70 fParticleCorrection[i] = 0;
72 for (Int_t i=0; i<3; i++)
75 // Define input and output slots here
76 DefineInput(0, TChain::Class());
77 DefineOutput(0, TList::Class());
79 if (fOption.Contains("only-process-type-nd"))
80 fSelectProcessType = 1;
82 if (fOption.Contains("only-process-type-sd"))
83 fSelectProcessType = 2;
85 if (fOption.Contains("only-process-type-dd"))
86 fSelectProcessType = 3;
88 if (fSelectProcessType != 0)
89 AliInfo(Form("WARNING: Systematic study enabled. Only considering process type %d", fSelectProcessType));
92 AliMultiplicityTask::~AliMultiplicityTask()
98 // histograms are in the output list and deleted when the output
99 // list is deleted by the TSelector dtor
107 //________________________________________________________________________
108 void AliMultiplicityTask::ConnectInputData(Option_t *)
113 Printf("AliMultiplicityTask::ConnectInputData called");
115 TTree* tree = dynamic_cast<TTree*> (GetInputData(0));
117 Printf("ERROR: Could not read tree from input slot 0");
119 // Disable all branches and enable only the needed ones
121 tree->SetBranchStatus("*", 0);
123 tree->SetBranchStatus("AliESDHeader*", 1);
124 tree->SetBranchStatus("*Vertex*", 1);
126 if (fAnalysisMode & AliPWG0Helper::kSPD) {
127 tree->SetBranchStatus("AliMultiplicity*", 1);
130 if (fAnalysisMode & AliPWG0Helper::kTPC || fAnalysisMode & AliPWG0Helper::kTPCITS) {
131 //AliESDtrackCuts::EnableNeededBranches(tree);
132 tree->SetBranchStatus("Tracks*", 1);
136 AliESDInputHandler *esdH = dynamic_cast<AliESDInputHandler*> (AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler());
139 Printf("ERROR: Could not get ESDInputHandler");
141 fESD = esdH->GetEvent();
144 // disable info messages of AliMCEvent (per event)
145 AliLog::SetClassDebugLevel("AliMCEvent", AliLog::kWarning - AliLog::kDebug + 1);
148 void AliMultiplicityTask::CreateOutputObjects()
150 // create result objects and add to output list
155 fMultiplicity = new AliMultiplicityCorrection("Multiplicity", "Multiplicity");
156 fOutput->Add(fMultiplicity);
158 fdNdpT = new TH1F("fdNdpT", "fdNdpT", 1000, 0, 10);
160 fOutput->Add(fdNdpT);
162 if (fOption.Contains("particle-efficiency"))
163 for (Int_t i = 0; i<8; i++)
165 fParticleCorrection[i] = new AliCorrection(Form("correction_%d", i), Form("correction_%d", i));
166 fOutput->Add(fParticleCorrection[i]);
169 if (fOption.Contains("pt-spectrum-hist"))
171 TFile* file = TFile::Open("ptspectrum_fit.root");
174 TString subStr(fOption(fOption.Index("pt-spectrum")+17, 3));
175 TString histName(Form("ptspectrum_%s", subStr.Data()));
176 AliInfo(Form("Pt-Spectrum modification. Using %s.", histName.Data()));
177 fPtSpectrum = (TH1D*) file->Get(histName);
179 AliError("Histogram not found");
182 AliError("Could not open ptspectrum_fit.root. Pt Spectrum study could not be enabled.");
185 if (fOption.Contains("pt-spectrum-func"))
189 Printf("Using function for systematic p_t study");
193 Printf("ERROR: Could not find function for systematic p_t study");
194 fPtSpectrum = new TH1D("fPtSpectrum", "fPtSpectrum", 1, 0, 100);
195 fPtSpectrum->SetBinContent(1, 1);
200 Printf("WARNING: Systematic study enabled. Pt spectrum will be modified");
202 if (fOption.Contains("particle-species")) {
203 fParticleSpecies = new TNtuple("fParticleSpecies", "fParticleSpecies", "vtx:Pi_True:K_True:p_True:o_True:Pi_Rec:K_Rec:p_Rec:o_Rec:nolabel:doublePrim:doubleCount");
204 fOutput->Add(fParticleSpecies);
207 fTemp1 = new TH2F("fTemp1", "fTemp1", 100, -0.5, 99.5, 100, -0.5, 99.5);
208 fOutput->Add(fTemp1);
210 for (Int_t i=0; i<3; i++)
212 fEta[i] = new TH1F(Form("fEta_%d", i), ";#eta", 100, -2, 2);
213 fOutput->Add(fEta[i]);
216 fVertex = new TH3F("vertex_check", "vertex_check;x;y;z", 100, -1, 1, 100, -1, 1, 100, -30, 30);
217 fOutput->Add(fVertex);
219 fEtaPhi = new TH2F("fEtaPhi", "fEtaPhi;#eta;#phi in rad.;count", 80, -4, 4, 18*20, 0, 2 * TMath::Pi());
220 fOutput->Add(fEtaPhi);
222 // TODO set seed for random generator
225 void AliMultiplicityTask::Exec(Option_t*)
229 // Check prerequisites
232 AliDebug(AliLog::kError, "ESD not available");
236 AliInputEventHandler* inputHandler = (AliInputEventHandler*) AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler();
239 Printf("ERROR: Could not receive input handler");
243 Bool_t eventTriggered = inputHandler->IsEventSelected();
245 static AliTriggerAnalysis* triggerAnalysis = 0;
246 if (!triggerAnalysis)
248 AliPhysicsSelection* physicsSelection = dynamic_cast<AliPhysicsSelection*> (inputHandler->GetEventSelection());
249 triggerAnalysis = physicsSelection->GetTriggerAnalysis();
252 eventTriggered = triggerAnalysis->IsTriggerFired(fESD, fTrigger);
254 const AliESDVertex* vtxESD = AliPWG0Helper::GetVertex(fESD, fAnalysisMode);
255 if (vtxESD && !AliPWG0Helper::TestVertex(vtxESD, fAnalysisMode))
258 // remove vertices outside +- 15 cm
259 if (vtxESD && TMath::Abs(vtxESD->GetZv()) > 15)
262 Bool_t eventVertex = (vtxESD != 0);
268 fVertex->Fill(vtxESD->GetXv(), vtxESD->GetYv(), vtxESD->GetZv());
271 // post the data already here
272 PostData(0, fOutput);
274 //const Float_t kPtCut = 0.3;
276 // create list of (label, eta) tuples
277 Int_t inputCount = 0;
280 if (fAnalysisMode & AliPWG0Helper::kSPD)
285 const AliMultiplicity* mult = fESD->GetMultiplicity();
288 AliDebug(AliLog::kError, "AliMultiplicity not available");
292 labelArr = new Int_t[mult->GetNumberOfTracklets()];
293 etaArr = new Float_t[mult->GetNumberOfTracklets()];
295 Bool_t foundInEta10 = kFALSE;
297 // get multiplicity from ITS tracklets
298 for (Int_t i=0; i<mult->GetNumberOfTracklets(); ++i)
300 //printf("%d %f %f %f\n", i, mult->GetTheta(i), mult->GetPhi(i), mult->GetDeltaPhi(i));
302 Float_t deltaPhi = mult->GetDeltaPhi(i);
304 if (fDeltaPhiCut > 0 && TMath::Abs(deltaPhi) > fDeltaPhiCut)
307 if (fSystSkipParticles && gRandom->Uniform() < (0.0153))
309 Printf("Skipped tracklet!");
313 etaArr[inputCount] = mult->GetEta(i);
314 if (mult->GetLabel(i, 0) == mult->GetLabel(i, 1))
316 labelArr[inputCount] = mult->GetLabel(i, 0);
319 labelArr[inputCount] = -1;
321 for (Int_t j=0; j<3; j++)
323 if (vtx[2] > fMultiplicity->GetVertexBegin(j) && vtx[2] < fMultiplicity->GetVertexEnd(j))
324 fEta[j]->Fill(etaArr[inputCount]);
327 if (vtxESD && TMath::Abs(vtx[2]) < 10)
328 fEtaPhi->Fill(etaArr[inputCount], mult->GetPhi(i));
330 // we have to repeat the trigger here, because the tracklet might have been kicked out fSystSkipParticles
331 if (TMath::Abs(etaArr[inputCount]) < 1)
332 foundInEta10 = kTRUE;
337 if (fSystSkipParticles && (fTrigger & AliTriggerAnalysis::kOneParticle) && !foundInEta10)
338 eventTriggered = kFALSE;
341 else if (fAnalysisMode & AliPWG0Helper::kTPC || fAnalysisMode & AliPWG0Helper::kTPCITS)
345 AliDebug(AliLog::kError, "fESDTrackCuts not available");
351 // get multiplicity from ESD tracks
352 TObjArray* list = fEsdTrackCuts->GetAcceptedTracks(fESD, (fAnalysisMode & AliPWG0Helper::kTPC));
353 Int_t nGoodTracks = list->GetEntries();
355 labelArr = new Int_t[nGoodTracks];
356 etaArr = new Float_t[nGoodTracks];
358 // loop over esd tracks
359 for (Int_t i=0; i<nGoodTracks; i++)
361 AliESDtrack* esdTrack = dynamic_cast<AliESDtrack*> (list->At(i));
364 AliDebug(AliLog::kError, Form("ERROR: Could not retrieve track %d.", i));
368 if (esdTrack->Pt() < 0.15)
371 Float_t d0z0[2],covd0z0[3];
372 esdTrack->GetImpactParameters(d0z0,covd0z0);
373 Float_t sigma= 0.0050+0.0060/TMath::Power(esdTrack->Pt(),0.9);
374 Float_t d0max = 7.*sigma;
375 if (TMath::Abs(d0z0[0]) > d0max)
378 if (vtxESD && TMath::Abs(vtx[2]) < 10)
379 fEtaPhi->Fill(esdTrack->Eta(), esdTrack->Phi());
381 etaArr[inputCount] = esdTrack->Eta();
382 labelArr[inputCount] = TMath::Abs(esdTrack->GetLabel());
390 // eta range for nMCTracksSpecies, nESDTracksSpecies
391 Float_t etaRange = 1.0;
392 // switch (fAnalysisMode) {
393 // case AliPWG0Helper::kInvalid: break;
394 // case AliPWG0Helper::kTPC:
395 // case AliPWG0Helper::kTPCITS:
396 // etaRange = 1.0; break;
397 // case AliPWG0Helper::kSPD: etaRange = 1.0; break;
401 if (!fReadMC) // Processing of ESD information
403 Int_t nESDTracks05 = 0;
404 Int_t nESDTracks10 = 0;
405 Int_t nESDTracks14 = 0;
407 for (Int_t i=0; i<inputCount; ++i)
409 Float_t eta = etaArr[i];
411 if (TMath::Abs(eta) < 0.5)
414 if (TMath::Abs(eta) < 1.0)
417 if (TMath::Abs(eta) < 1.3)
421 //if (nESDTracks05 >= 20 || nESDTracks10 >= 30 || nESDTracks14 >= 32)
422 // Printf("File: %s, IEV: %d, TRG: ---, Orbit: 0x%x, Period: %d, BC: %d; Tracks: %d %d %d", ((TTree*) GetInputData(0))->GetCurrentFile()->GetName(), fESD->GetEventNumberInFile(), fESD->GetOrbitNumber(),fESD->GetPeriodNumber(),fESD->GetBunchCrossNumber(), nESDTracks05, nESDTracks10, nESDTracks14);
425 fMultiplicity->FillTriggeredEvent(nESDTracks05, nESDTracks10, nESDTracks14);
427 if (eventTriggered && eventVertex)
428 fMultiplicity->FillMeasured(vtx[2], nESDTracks05, nESDTracks10, nESDTracks14);
430 else if (fReadMC) // Processing of MC information
432 AliMCEventHandler* eventHandler = dynamic_cast<AliMCEventHandler*> (AliAnalysisManager::GetAnalysisManager()->GetMCtruthEventHandler());
434 Printf("ERROR: Could not retrieve MC event handler");
438 AliMCEvent* mcEvent = eventHandler->MCEvent();
440 Printf("ERROR: Could not retrieve MC event");
444 AliStack* stack = mcEvent->Stack();
447 AliDebug(AliLog::kError, "Stack not available");
451 AliHeader* header = mcEvent->Header();
454 AliDebug(AliLog::kError, "Header not available");
460 Printf("WARNING: Replacing vertex by MC vertex. This is for systematical checks only.");
462 AliGenEventHeader* genHeader = header->GenEventHeader();
464 genHeader->PrimaryVertex(vtxMC);
469 // get process information
470 AliPWG0Helper::MCProcessType processType = AliPWG0Helper::GetEventProcessType(fESD, header, stack, fDiffTreatment);
472 Bool_t processEvent = kTRUE;
473 if (fSelectProcessType > 0)
475 processEvent = kFALSE;
478 if (fSelectProcessType == 1 && processType == AliPWG0Helper::kND)
479 processEvent = kTRUE;
481 // single diffractive
482 if (fSelectProcessType == 2 && processType == AliPWG0Helper::kSD)
483 processEvent = kTRUE;
485 // double diffractive
486 if (fSelectProcessType == 3 && processType == AliPWG0Helper::kDD)
487 processEvent = kTRUE;
490 Printf("Skipping this event, because it is not of the requested process type (%d)", (Int_t) processType);
496 AliGenEventHeader* genHeader = header->GenEventHeader();
498 genHeader->PrimaryVertex(vtxMC);
500 // get number of tracks from MC
501 // loop over mc particles
502 Int_t nPrim = stack->GetNprimary();
503 Int_t nMCPart = stack->GetNtrack();
505 // tracks in different eta ranges
506 Int_t nMCTracks05 = 0;
507 Int_t nMCTracks10 = 0;
508 Int_t nMCTracks14 = 0;
509 Int_t nMCTracksAll = 0;
511 // tracks per particle species, in |eta| < 2 (systematic study)
512 Int_t nMCTracksSpecies[4]; // (pi, K, p, other)
513 for (Int_t i = 0; i<4; ++i)
514 nMCTracksSpecies[i] = 0;
516 for (Int_t iMc = 0; iMc < nPrim; ++iMc)
518 AliDebug(AliLog::kDebug+1, Form("MC Loop: Processing particle %d.", iMc));
520 TParticle* particle = stack->Particle(iMc);
524 AliDebug(AliLog::kError, Form("UNEXPECTED: particle with label %d not found in stack (mc loop).", iMc));
528 Bool_t debug = kFALSE;
529 if (AliPWG0Helper::IsPrimaryCharged(particle, nPrim, debug) == kFALSE)
531 //printf("%d) DROPPED ", iMC);
536 //printf("%d) OK ", iMC);
539 //if (particle->Pt() < kPtCut)
542 Int_t particleWeight = 1;
544 Float_t pt = particle->Pt();
546 // in case of systematic study, weight according to the change of the pt spectrum
547 // it cannot be just multiplied because we cannot count "half of a particle"
548 // instead a random generator decides if the particle is counted twice (if value > 1)
549 // or not (if value < 0)
552 Int_t bin = fPtSpectrum->FindBin(pt);
553 if (bin > 0 && bin <= fPtSpectrum->GetNbinsX())
555 Float_t factor = fPtSpectrum->GetBinContent(bin);
558 Float_t random = gRandom->Uniform();
559 if (factor > 1 && random < factor - 1)
563 else if (factor < 1 && random < 1 - factor)
569 //Printf("MC weight is: %d", particleWeight);
571 if (TMath::Abs(particle->Eta()) < 0.5)
572 nMCTracks05 += particleWeight;
574 if (TMath::Abs(particle->Eta()) < 1.0)
575 nMCTracks10 += particleWeight;
577 if (TMath::Abs(particle->Eta()) < 1.3)
578 nMCTracks14 += particleWeight;
580 nMCTracksAll += particleWeight;
582 if (particle->Pt() > 0 && TMath::Abs(particle->Eta()) < 1.0)
583 fdNdpT->Fill(particle->Pt(), 1.0 / particle->Pt());
585 if (fParticleCorrection[0] || fParticleSpecies)
588 switch (TMath::Abs(particle->GetPdgCode()))
590 case 211: id = 0; break;
591 case 321: id = 1; break;
592 case 2212: id = 2; break;
593 default: id = 3; break;
596 if (TMath::Abs(particle->Eta()) < etaRange)
597 nMCTracksSpecies[id]++;
599 if (fParticleCorrection[id])
600 fParticleCorrection[id]->GetTrackCorrection()->FillGene(vtxMC[2], particle->Eta(), particle->Pt());
602 } // end of mc particle
604 fMultiplicity->FillGenerated(vtxMC[2], eventTriggered, eventVertex, processType, (Int_t) nMCTracks05, (Int_t) nMCTracks10, (Int_t) nMCTracks14, (Int_t) nMCTracksAll);
607 Int_t nESDTracks05 = 0;
608 Int_t nESDTracks10 = 0;
609 Int_t nESDTracks14 = 0;
611 // tracks per particle species, in |eta| < 2 (systematic study)
612 Int_t nESDTracksSpecies[7]; // (pi, K, p, other, nolabel, doublecount_prim, doublecount_all)
613 for (Int_t i = 0; i<7; ++i)
614 nESDTracksSpecies[i] = 0;
616 Bool_t* foundPrimaries = new Bool_t[nPrim]; // to prevent double counting
617 for (Int_t i=0; i<nPrim; i++)
618 foundPrimaries[i] = kFALSE;
620 Bool_t* foundPrimaries2 = new Bool_t[nPrim]; // to prevent double counting
621 for (Int_t i=0; i<nPrim; i++)
622 foundPrimaries2[i] = kFALSE;
624 Bool_t* foundTracks = new Bool_t[nMCPart]; // to prevent double counting
625 for (Int_t i=0; i<nMCPart; i++)
626 foundTracks[i] = kFALSE;
628 for (Int_t i=0; i<inputCount; ++i)
630 Float_t eta = etaArr[i];
631 Int_t label = labelArr[i];
633 Int_t particleWeight = 1;
635 // in case of systematic study, weight according to the change of the pt spectrum
638 TParticle* mother = 0;
640 // preserve label for later
641 Int_t labelCopy = label;
643 labelCopy = AliPWG0Helper::FindPrimaryMotherLabel(stack, labelCopy);
645 mother = stack->Particle(labelCopy);
647 // in case of pt study we do not count particles w/o label, because they cannot be scaled
651 // it cannot be just multiplied because we cannot count "half of a particle"
652 // instead a random generator decides if the particle is counted twice (if value > 1)
653 // or not (if value < 0)
654 Int_t bin = fPtSpectrum->FindBin(mother->Pt());
655 if (bin > 0 && bin <= fPtSpectrum->GetNbinsX())
657 Float_t factor = fPtSpectrum->GetBinContent(bin);
660 Float_t random = gRandom->Uniform();
661 if (factor > 1 && random < factor - 1)
665 else if (factor < 1 && random < 1 - factor)
671 //Printf("ESD weight is: %d", particleWeight);
673 if (TMath::Abs(eta) < 0.5)
674 nESDTracks05 += particleWeight;
676 if (TMath::Abs(eta) < 1.0)
677 nESDTracks10 += particleWeight;
679 if (TMath::Abs(eta) < 1.3)
680 nESDTracks14 += particleWeight;
682 if (fParticleSpecies)
684 Int_t motherLabel = -1;
685 TParticle* mother = 0;
689 motherLabel = AliPWG0Helper::FindPrimaryMotherLabel(stack, label);
690 if (motherLabel >= 0)
691 mother = stack->Particle(motherLabel);
695 // count tracks that did not have a label
696 if (TMath::Abs(eta) < etaRange)
697 nESDTracksSpecies[4]++;
701 // get particle type (pion, proton, kaon, other) of mother
703 switch (TMath::Abs(mother->GetPdgCode()))
705 case 211: idMother = 0; break;
706 case 321: idMother = 1; break;
707 case 2212: idMother = 2; break;
708 default: idMother = 3; break;
711 // double counting is ok for particle ratio study
712 if (TMath::Abs(eta) < etaRange)
713 nESDTracksSpecies[idMother]++;
715 // double counting is not ok for efficiency study
717 // check if we already counted this particle, this way distinguishes double counted particles (bug/artefact in tracking) or double counted primaries due to secondaries (physics)
718 if (foundTracks[label])
720 if (TMath::Abs(eta) < etaRange)
721 nESDTracksSpecies[6]++;
725 foundTracks[label] = kTRUE;
727 // particle (primary) already counted?
728 if (foundPrimaries[motherLabel])
730 if (TMath::Abs(eta) < etaRange)
731 nESDTracksSpecies[5]++;
734 foundPrimaries[motherLabel] = kTRUE;
739 if (fParticleCorrection[0])
741 if (label >= 0 && stack->IsPhysicalPrimary(label))
743 TParticle* particle = stack->Particle(label);
745 // get particle type (pion, proton, kaon, other)
747 switch (TMath::Abs(particle->GetPdgCode()))
749 case 211: id = 0; break;
750 case 321: id = 1; break;
751 case 2212: id = 2; break;
752 default: id = 3; break;
755 // todo check if values are not completely off??
757 // particle (primary) already counted?
758 if (!foundPrimaries2[label])
760 foundPrimaries2[label] = kTRUE;
761 fParticleCorrection[id]->GetTrackCorrection()->FillMeas(vtxMC[2], particle->Eta(), particle->Pt());
767 if (fParticleCorrection[0])
769 // if the particle decays/stops before this radius we do not see it
770 // 8cm larger than SPD layer 2
771 // 123cm TPC radius where a track has about 50 clusters (cut limit)
772 const Float_t endRadius = (fAnalysisMode & AliPWG0Helper::kSPD) ? 8. : 123;
774 // loop over all primaries that have not been found
775 for (Int_t i=0; i<nPrim; i++)
778 if (foundPrimaries2[i])
781 TParticle* particle = 0;
782 TClonesArray* trackrefs = 0;
783 mcEvent->GetParticleAndTR(i, particle, trackrefs);
785 // true primary and charged
786 if (!AliPWG0Helper::IsPrimaryCharged(particle, nPrim))
789 //skip particles with larger |eta| than 3, to keep the log clean, is anyway not included in correction map
790 if (TMath::Abs(particle->Eta()) > 3)
793 // skipping checking of process type of daughter: Neither kPBrem, kPDeltaRay nor kPCerenkov should appear in the event generation
795 // get particle type (pion, proton, kaon, other)
797 switch (TMath::Abs(particle->GetPdgCode()))
799 case 211: id = 4; break;
800 case 321: id = 5; break;
801 case 2212: id = 6; break;
802 default: id = 7; break;
805 if (!fParticleCorrection[id])
808 // get last track reference
809 AliTrackReference* trackref = dynamic_cast<AliTrackReference*> (trackrefs->Last());
813 Printf("ERROR: Could not get trackref of %d (count %d)", i, trackrefs->GetEntries());
818 // particle in tracking volume long enough...
819 if (trackref->R() > endRadius)
822 if (particle->GetLastDaughter() >= 0)
824 Int_t uID = stack->Particle(particle->GetLastDaughter())->GetUniqueID();
825 //if (uID != kPBrem && uID != kPDeltaRay && uID < kPCerenkov)
830 Printf("Particle %d (%s) decayed at %f, daugher uniqueID: %d:", i, particle->GetName(), trackref->R(), uID);
832 Printf("Daugthers:");
833 for (Int_t d = particle->GetFirstDaughter(); d <= particle->GetLastDaughter(); d++)
834 stack->Particle(d)->Print();
837 fParticleCorrection[id]->GetTrackCorrection()->FillGene(vtxMC[2], particle->Eta(), particle->Pt());
842 if (trackref->DetectorId() == -1)
845 Printf("Particle %d stopped at %f:", i, trackref->R());
849 fParticleCorrection[id]->GetTrackCorrection()->FillMeas(vtxMC[2], particle->Eta(), particle->Pt());
853 Printf("Particle %d simply not tracked", i);
859 delete[] foundTracks;
860 delete[] foundPrimaries;
861 delete[] foundPrimaries2;
863 // if ((Int_t) nMCTracks14 > 10 && nESDTracks14 <= 3)
865 // TTree* tree = dynamic_cast<TTree*> (GetInputData(0));
866 // printf("WARNING: Event %lld %s (vtx-z = %f, recon: %f, contrib: %d, res: %f) has %d generated and %d reconstructed...\n", tree->GetReadEntry(), tree->GetCurrentFile()->GetName(), vtxMC[2], vtx[2], vtxESD->GetNContributors(), vtxESD->GetZRes(), nMCTracks14, nESDTracks14);
871 fMultiplicity->FillTriggeredEvent(nESDTracks05, nESDTracks10, nESDTracks14);
872 fMultiplicity->FillNoVertexEvent(vtxMC[2], eventVertex, nMCTracks05, nMCTracks10, nMCTracks14, nESDTracks05, nESDTracks10, nESDTracks14);
875 // if (nESDTracks05 == 0)
876 // fTemp1->Fill(nMCTracks05, nESDTracks05);
880 if (eventTriggered && eventVertex)
882 // fill response matrix using vtxMC (best guess)
883 fMultiplicity->FillCorrection(vtxMC[2], nMCTracks05, nMCTracks10, nMCTracks14, nMCTracksAll, nESDTracks05, nESDTracks10, nESDTracks14);
885 fMultiplicity->FillMeasured(vtx[2], nESDTracks05, nESDTracks10, nESDTracks14);
887 if (fParticleSpecies)
888 fParticleSpecies->Fill(vtxMC[2], nMCTracksSpecies[0], nMCTracksSpecies[1], nMCTracksSpecies[2], nMCTracksSpecies[3], nESDTracksSpecies[0], nESDTracksSpecies[1], nESDTracksSpecies[2], nESDTracksSpecies[3], nESDTracksSpecies[4], nESDTracksSpecies[5], nESDTracksSpecies[6]);
899 void AliMultiplicityTask::Terminate(Option_t *)
901 // The Terminate() function is the last function to be called during
902 // a query. It always runs on the client, it can be used to present
903 // the results graphically or save the results to file.
905 fOutput = dynamic_cast<TList*> (GetOutputData(0));
907 Printf("ERROR: fOutput not available");
911 fMultiplicity = dynamic_cast<AliMultiplicityCorrection*> (fOutput->FindObject("Multiplicity"));
912 for (Int_t i = 0; i < 8; ++i)
913 fParticleCorrection[i] = dynamic_cast<AliCorrection*> (fOutput->FindObject(Form("correction_%d", i)));
914 fParticleSpecies = dynamic_cast<TNtuple*> (fOutput->FindObject("fParticleSpecies"));
916 fdNdpT = dynamic_cast<TH1*> (fOutput->FindObject("fdNdpT"));
920 AliDebug(AliLog::kError, Form("ERROR: Histograms not available %p", (void*) fMultiplicity));
924 TString fileName("multiplicity");
925 if (fSelectProcessType == 1)
927 if (fSelectProcessType == 2)
929 if (fSelectProcessType == 3)
932 TFile* file = TFile::Open(fileName, "RECREATE");
934 fMultiplicity->SaveHistograms();
935 for (Int_t i = 0; i < 8; ++i)
936 if (fParticleCorrection[i])
937 fParticleCorrection[i]->SaveHistograms();
938 if (fParticleSpecies)
939 fParticleSpecies->Write();
943 fTemp1 = dynamic_cast<TH1*> (fOutput->FindObject("fTemp1"));
947 fEtaPhi = dynamic_cast<TH2F*> (fOutput->FindObject("fEtaPhi"));
951 fVertex = dynamic_cast<TH3F*> (fOutput->FindObject("vertex_check"));
955 for (Int_t i=0; i<3; i++)
957 fEta[i] = dynamic_cast<TH1*> (fOutput->FindObject(Form("fEta_%d", i)));
961 Float_t events = fMultiplicity->GetMultiplicityESD(i)->Integral(1, fMultiplicity->GetMultiplicityESD(i)->GetNbinsX());
963 fEta[i]->Scale(1.0 / events);
964 fEta[i]->Scale(1.0 / fEta[i]->GetXaxis()->GetBinWidth(1));
969 TObjString option(fOption);
974 Printf("Written result to %s", fileName.Data());