1 /**************************************************************************
2 * Copyright(c) 1998-2009, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
7 * Permission to use, copy, modify and distribute this software and its *
8 * documentation strictly for non-commercial purposes is hereby granted *
9 * without fee, provided that the above copyright notice appears in all *
10 * copies and that both the copyright notice and this permission notice *
11 * appear in the supporting documentation. The authors make no claims *
12 * about the suitability of this software for any purpose. It is *
13 * provided "as is" without express or implied warranty. *
14 **************************************************************************/
16 //*************************************************************************
17 // Class AliAnalysisTaskVertexESD
18 // AliAnalysisTask to extract from ESD the information for the analysis
19 // of the primary vertex reconstruction efficiency and resolution
20 // (for MC events) and distributions (for real data). Three vertices:
25 // Author: A.Dainese, andrea.dainese@pd.infn.it
26 //*************************************************************************
32 #include <TClonesArray.h>
33 #include <TObjArray.h>
38 #include "AliAnalysisTask.h"
39 #include "AliAnalysisManager.h"
41 #include "AliMultiplicity.h"
42 #include "AliVertexerTracks.h"
43 #include "AliESDtrack.h"
44 #include "AliExternalTrackParam.h"
45 #include "AliESDVertex.h"
46 #include "AliESDEvent.h"
47 #include "AliESDInputHandler.h"
48 #include "AliTrackReference.h"
49 //#include "AliTriggerAnalysis.h"
51 #include "AliMCEventHandler.h"
52 #include "AliMCEvent.h"
56 #include "AliGenEventHeader.h"
57 #include "AliAnalysisTaskVertexESD.h"
60 ClassImp(AliAnalysisTaskVertexESD)
62 //________________________________________________________________________
63 AliAnalysisTaskVertexESD::AliAnalysisTaskVertexESD(const char *name) :
64 AliAnalysisTaskSE(name),
65 fCheckEventType(kTRUE),
68 fRecoVtxITSTPC(kFALSE),
69 fRecoVtxITSTPCHalfEvent(kFALSE),
70 fOnlyITSTPCTracks(kFALSE),
71 fOnlyITSSATracks(kFALSE),
73 fFillTreeBeamSpot(kFALSE),
91 // Define input and output slots here
92 // Output slot #0 writes into a TList container
93 DefineOutput(1, TList::Class()); //My private output
95 //________________________________________________________________________
96 AliAnalysisTaskVertexESD::~AliAnalysisTaskVertexESD()
100 // histograms are in the output list and deleted when the output
101 // list is deleted by the TSelector dtor
110 //________________________________________________________________________
111 void AliAnalysisTaskVertexESD::UserCreateOutputObjects()
116 // Several histograms are more conveniently managed in a TList
120 fNtupleVertexESD = new TNtuple("fNtupleVertexESD","vertices","run:tstamp:bunchcross:triggered:dndygen:xtrue:ytrue:ztrue:xSPD:xerrSPD:ySPD:yerrSPD:zSPD:zerrSPD:ntrksSPD:SPD3D:dphiSPD:xTPC:xerrTPC:yTPC:yerrTPC:zTPC:zerrTPC:ntrksTPC:constrTPC:xTRK:xerrTRK:yTRK:yerrTRK:zTRK:zerrTRK:ntrksTRK:constrTRK:ntrklets:nESDtracks:nITSrefit5or6:nTPCin:nTPCinEta09:SPD0cls:xSPDp:xerrSPDp:ySPDp:yerrSPDp:zSPDp:zerrSPDp:ntrksSPDp:xTPCnc:xerrTPCnc:yTPCnc:yerrTPCnc:zTPCnc:zerrTPCnc:ntrksTPCnc:xTPCc:xerrTPCc:yTPCc:yerrTPCc:zTPCc:zerrTPCc:ntrksTPCc:xTRKnc:xerrTRKnc:yTRKnc:yerrTRKnc:zTRKnc:zerrTRKnc:ntrksTRKnc:xTRKc:xerrTRKc:yTRKc:yerrTRKc:zTRKc:zerrTRKc:ntrksTRKc:deltaxTRKnc:deltayTRKnc:deltazTRKnc:deltaxerrTRKnc:deltayerrTRKnc:deltazerrTRKnc:ntrksEvenTRKnc:ntrksOddTRKnc");
122 fOutput->Add(fNtupleVertexESD);
124 fhSPDVertexX = new TH1F("fhSPDVertexX","SPDVertex x; x vertex [cm]; events",200,-1,1);
125 fOutput->Add(fhSPDVertexX);
126 fhSPDVertexY = new TH1F("fhSPDVertexY","SPDVertex y; y vertex [cm]; events",200,-1,1);
127 fOutput->Add(fhSPDVertexY);
128 fhSPDVertexZ = new TH1F("fhSPDVertexZ","SPDVertex z; z vertex [cm]; events",200,-20,20);
129 fOutput->Add(fhSPDVertexZ);
130 fhTRKVertexX = new TH1F("fhTRKVertexX","TRKVertex x; x vertex [cm]; events",200,-1,1);
131 fOutput->Add(fhTRKVertexX);
132 fhTRKVertexY = new TH1F("fhTRKVertexY","TRKVertex y; y vertex [cm]; events",200,-1,1);
133 fOutput->Add(fhTRKVertexY);
134 fhTRKVertexZ = new TH1F("fhTRKVertexZ","TRKVertex z; z vertex [cm]; events",200,-20,20);
135 fOutput->Add(fhTRKVertexZ);
136 fhTPCVertexX = new TH1F("fhTPCVertexX","TPCVertex x; x vertex [cm]; events",200,-3,3);
137 fOutput->Add(fhTPCVertexX);
138 fhTPCVertexY = new TH1F("fhTPCVertexY","TPCVertex y; y vertex [cm]; events",200,-3,3);
139 fOutput->Add(fhTPCVertexY);
140 fhTPCVertexZ = new TH1F("fhTPCVertexZ","TPCVertex z; z vertex [cm]; events",200,-20,20);
141 fOutput->Add(fhTPCVertexZ);
143 fhTrackRefs = new TH2F("fhTrackRefs","Track references; x; y",1000,-4,4,1000,-4,4);
144 fOutput->Add(fhTrackRefs);
146 fTreeBeamSpot = new TTree("fTreeBeamSpot", "BeamSpotTree");
147 UShort_t triggered, ntrkletsS, ntrksTRKnc;
149 Float_t cetTimeLHC,xTRKnc, yTRKnc, zTRKnc;
150 fTreeBeamSpot->Branch("run", &run, "run/i");
151 fTreeBeamSpot->Branch("cetTimeLHC", &cetTimeLHC, "cetTimeLHC/F");
152 fTreeBeamSpot->Branch("bx", &bx, "bx/i");
153 fTreeBeamSpot->Branch("triggered", &triggered, "triggered/s");
154 fTreeBeamSpot->Branch("ntrkletsS", &ntrkletsS, "ntrkletsS/s");
155 fTreeBeamSpot->Branch("xTRKnc", &xTRKnc, "xTRKnc/F");
156 fTreeBeamSpot->Branch("yTRKnc", &yTRKnc, "yTRKnc/F");
157 fTreeBeamSpot->Branch("zTRKnc", &zTRKnc, "zTRKnc/F");
158 fTreeBeamSpot->Branch("ntrksTRKnc", &ntrksTRKnc, "ntrksTRKnc/s");
159 fOutput->Add(fTreeBeamSpot);
161 PostData(1, fOutput);
166 //________________________________________________________________________
167 void AliAnalysisTaskVertexESD::UserExec(Option_t *)
170 // Called for each event
173 Printf("ERROR: fESD not available");
176 AliESDEvent* esdE = (AliESDEvent*) InputEvent();
178 if(fCheckEventType && (esdE->GetEventType())!=7) return;
182 mcVertex[0]=9999.; mcVertex[1]=9999.; mcVertex[2]=9999.;
183 Float_t dNchdy=-999.;
184 // *********** MC info ***************
186 AliMCEventHandler* eventHandler = dynamic_cast<AliMCEventHandler*> (AliAnalysisManager::GetAnalysisManager()->GetMCtruthEventHandler());
188 Printf("ERROR: Could not retrieve MC event handler");
192 AliMCEvent* mcEvent = eventHandler->MCEvent();
194 Printf("ERROR: Could not retrieve MC event");
198 AliStack* stack = mcEvent->Stack();
200 AliDebug(AliLog::kError, "Stack not available");
204 AliHeader* header = mcEvent->Header();
206 AliDebug(AliLog::kError, "Header not available");
209 AliGenEventHeader* genHeader = header->GenEventHeader();
210 genHeader->PrimaryVertex(mcVertex);
212 Int_t ngenpart = (Int_t)stack->GetNtrack();
213 //printf("# generated particles = %d\n",ngenpart);
215 for(Int_t ip=0; ip<ngenpart; ip++) {
216 TParticle* part = (TParticle*)stack->Particle(ip);
217 // keep only electorns, muons, pions, kaons and protons
218 Int_t apdg = TMath::Abs(part->GetPdgCode());
219 if(apdg!=11 && apdg!=13 && apdg!=211 && apdg!=321 && apdg!=2212) continue;
220 // reject secondaries
221 if(TMath::Sqrt((part->Vx()-mcVertex[0])*(part->Vx()-mcVertex[0])+(part->Vy()-mcVertex[1])*(part->Vy()-mcVertex[1]))>0.0010) continue;
222 // reject incoming protons
223 Double_t energy = part->Energy();
224 if(energy>900.) continue;
225 Double_t pz = part->Pz();
226 Double_t y = 0.5*TMath::Log((energy+pz+1.e-13)/(energy-pz+1.e-13));
227 if(TMath::Abs(y)<1.0) dNchdy += 0.5; // count 1/2 of particles in |y|<1
229 TClonesArray *trefs=0;
230 Int_t ntrefs = mcEvent->GetParticleAndTR(ip,part,trefs);
231 if(ntrefs<0) continue;
232 for(Int_t iref=0; iref<ntrefs; iref++) {
233 AliTrackReference *tref = (AliTrackReference*)trefs->At(iref);
234 if(tref->R()>10.) continue;
235 fhTrackRefs->Fill(tref->X(),tref->Y());
238 //printf("# primary particles = %7.1f\n",dNchdy);
240 // *********** MC info ***************
244 //ULong64_t triggerMask;
245 //ULong64_t spdFO = (1 << 14);
246 //ULong64_t v0left = (1 << 10);
247 //ULong64_t v0right = (1 << 11);
249 //triggerMask=esdE->GetTriggerMask();
250 // MB1: SPDFO || V0L || V0R
251 //Bool_t eventTriggered = (triggerMask & spdFO || ((triggerMask & v0left) || (triggerMask & v0right)));
253 //triggerMask & spdFO && ((triggerMask&v0left) || (triggerMask&v0right))
254 //Bool_t eventTriggered = (triggerMask & spdFO);
256 //static AliTriggerAnalysis* triggerAnalysis = new AliTriggerAnalysis();
257 Bool_t eventTriggered = 0;//triggerAnalysis->IsTriggerFired(esdE, AliTriggerAnalysis::kSPDGFO /*| AliTriggerAnalysis::kOfflineFlag*/);
259 // use response of AliPhysicsSelection
260 eventTriggered = ((AliInputEventHandler*)(AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler()))->IsEventSelected();
263 Int_t nInFile = esdE->GetEventNumberInFile();
265 const AliMultiplicity *alimult = esdE->GetMultiplicity();
266 Int_t ntrklets=0,spd0cls=0;
268 ntrklets = alimult->GetNumberOfTracklets();
270 for(Int_t l=0;l<alimult->GetNumberOfTracklets();l++){
271 if(alimult->GetDeltaPhi(l)<-9998.) ntrklets--;
273 spd0cls = alimult->GetNumberOfSingleClusters()+ntrklets;
277 UShort_t triggered, ntrkletsS, ntrksTRKnc;
279 Float_t cetTimeLHC,xTRKnc, yTRKnc, zTRKnc;
280 fTreeBeamSpot->SetBranchAddress("run", &run);
281 fTreeBeamSpot->SetBranchAddress("cetTimeLHC", &cetTimeLHC);
282 fTreeBeamSpot->SetBranchAddress("bx", &bx);
283 fTreeBeamSpot->SetBranchAddress("triggered", &triggered);
284 fTreeBeamSpot->SetBranchAddress("ntrkletsS", &ntrkletsS);
285 fTreeBeamSpot->SetBranchAddress("xTRKnc", &xTRKnc);
286 fTreeBeamSpot->SetBranchAddress("yTRKnc", &yTRKnc);
287 fTreeBeamSpot->SetBranchAddress("zTRKnc", &zTRKnc);
288 fTreeBeamSpot->SetBranchAddress("ntrksTRKnc", &ntrksTRKnc);
291 Double_t tstamp = esdE->GetTimeStamp();
292 Float_t cetTime =(tstamp-1262304000.)+7200.;
294 Float_t cetTime1h =(tstamp-1262304000.)+3600.;
296 cetTimeLHC = (Float_t)cetTime1h;
298 Int_t ntracks = esdE->GetNumberOfTracks();
299 Int_t nITS5or6=0,nTPCin=0,nTPCinEta09=0;
300 //printf("Tracks # = %d\n",esdE->GetNumberOfTracks());
301 for(Int_t itr=0; itr<ntracks; itr++) {
302 AliESDtrack *t = esdE->GetTrack(itr);
303 if(t->GetNcls(0)>=5) nITS5or6++;
304 Double_t z0; t->GetZAt(0,esdE->GetMagneticField(),z0);
305 if(t->GetNcls(1)>0 && TMath::Abs(t->GetD(0,0,esdE->GetMagneticField()))<2.8 && TMath::Abs(z0)<20) {
307 if(TMath::Abs(t->GetTgl())<1.5) nTPCinEta09++;
312 const AliESDVertex *spdv=esdE->GetPrimaryVertexSPD();
313 const AliESDVertex *spdvp=esdE->GetPileupVertexSPD(0);
314 const AliESDVertex *tpcv=esdE->GetPrimaryVertexTPC();
315 const AliESDVertex *trkv=esdE->GetPrimaryVertexTracks();
320 if(spdv->GetNContributors()>0) {
321 TString title=spdv->GetTitle();
322 if(title.Contains("3D")) {
323 fhSPDVertexX->Fill(spdv->GetXv());
324 fhSPDVertexY->Fill(spdv->GetYv());
326 fhSPDVertexZ->Fill(spdv->GetZv());
331 if(trkv->GetNContributors()>0) {
332 fhTRKVertexX->Fill(trkv->GetXv());
333 fhTRKVertexY->Fill(trkv->GetYv());
334 fhTRKVertexZ->Fill(trkv->GetZv());
339 if(tpcv->GetNContributors()>0) {
340 fhTPCVertexX->Fill(tpcv->GetXv());
341 fhTPCVertexY->Fill(tpcv->GetYv());
342 fhTPCVertexZ->Fill(tpcv->GetZv());
349 Float_t expile=-999.;
350 Float_t eypile=-999.;
351 Float_t ezpile=-999.;
354 if(esdE->GetNumberOfPileupVerticesSPD()>0 && spdvp){
355 xpile=spdvp->GetXv();
356 expile=spdvp->GetXRes();
357 ypile=spdvp->GetYv();
358 eypile=spdvp->GetYRes();
359 zpile=spdvp->GetZv();
360 ezpile=spdvp->GetZRes();
361 ntrkspile=spdvp->GetNContributors();
366 Float_t xnt[82]; for(Int_t iii=0;iii<isize;iii++) xnt[iii]=0.;
369 Float_t secnt[9]; for(Int_t iii=0;iii<isizeSecNt;iii++) secnt[iii]=0.;
374 xnt[index++]=(Float_t)esdE->GetRunNumber();
375 secnt[indexSecNt++]=(Float_t)esdE->GetRunNumber();
376 run = (Int_t)esdE->GetRunNumber();
378 xnt[index++]=cetTime; //(Float_t)esdE->GetTimeStamp();
379 //secnt[indexSecNt++]=cetTime;
380 secnt[indexSecNt++]=cetTime1h;
382 xnt[index++]=(Float_t)esdE->GetBunchCrossNumber();
383 secnt[indexSecNt++]=(Float_t)esdE->GetBunchCrossNumber();
384 bx = (Int_t)esdE->GetBunchCrossNumber();
386 xnt[index++]=(eventTriggered ? 1. : 0.);
387 secnt[indexSecNt++]=(eventTriggered ? 1. : 0.);
388 triggered = (UShort_t)(eventTriggered ? 1 : 0);
390 xnt[index++]=(Float_t)dNchdy;
392 xnt[index++]=mcVertex[0];
393 xnt[index++]=mcVertex[1];
394 xnt[index++]=mcVertex[2];
396 xnt[index++]=spdv->GetXv();
397 xnt[index++]=spdv->GetXRes();
398 xnt[index++]=spdv->GetYv();
399 xnt[index++]=spdv->GetYRes();
400 xnt[index++]=spdv->GetZv();
401 xnt[index++]=spdv->GetZRes();
402 xnt[index++]=spdv->GetNContributors();
403 TString spdtitle = spdv->GetTitle();
404 xnt[index++]=(spdtitle.Contains("vertexer: 3D") ? 1. : 0.);
405 xnt[index++]=spdv->GetDispersion();
407 xnt[index++]=tpcv->GetXv();
408 xnt[index++]=tpcv->GetXRes();
409 xnt[index++]=tpcv->GetYv();
410 xnt[index++]=tpcv->GetYRes();
411 xnt[index++]=tpcv->GetZv();
412 xnt[index++]=tpcv->GetZRes();
413 xnt[index++]=tpcv->GetNContributors();
414 TString tpctitle = tpcv->GetTitle();
415 xnt[index++]=(tpctitle.Contains("WithConstraint") ? 1. : 0.);
417 xnt[index++]=trkv->GetXv();
418 xnt[index++]=trkv->GetXRes();
419 xnt[index++]=trkv->GetYv();
420 xnt[index++]=trkv->GetYRes();
421 xnt[index++]=trkv->GetZv();
422 xnt[index++]=trkv->GetZRes();
423 xnt[index++]=trkv->GetNContributors();// tpccontrorig;
424 TString trktitle = trkv->GetTitle();
425 xnt[index++]=(trktitle.Contains("WithConstraint") ? 1. : 0.);
427 xnt[index++]=float(ntrklets);
428 secnt[indexSecNt++]=float(ntrklets);
429 ntrkletsS = (UShort_t)ntrklets;
431 xnt[index++]=float(ntracks);
432 xnt[index++]=float(nITS5or6);
434 xnt[index++]=float(nTPCin);
435 xnt[index++]=float(nTPCinEta09);
437 xnt[index++]=spd0cls;
445 xnt[index++]=(Float_t)ntrkspile;
449 // add recalculated vertices TRK and TPC
452 AliESDVertex *tpcvnc = ReconstructPrimaryVertexTPC(kFALSE);
453 xnt[index++]=tpcvnc->GetXv();
454 xnt[index++]=tpcvnc->GetXRes();
455 xnt[index++]=tpcvnc->GetYv();
456 xnt[index++]=tpcvnc->GetYRes();
457 xnt[index++]=tpcvnc->GetZv();
458 xnt[index++]=tpcvnc->GetZRes();
459 xnt[index++]=tpcvnc->GetNContributors();
460 delete tpcvnc; tpcvnc=0;
462 AliESDVertex *tpcvc = ReconstructPrimaryVertexTPC(kTRUE);
463 xnt[index++]=tpcvc->GetXv();
464 xnt[index++]=tpcvc->GetXRes();
465 xnt[index++]=tpcvc->GetYv();
466 xnt[index++]=tpcvc->GetYRes();
467 xnt[index++]=tpcvc->GetZv();
468 xnt[index++]=tpcvc->GetZRes();
469 xnt[index++]=tpcvc->GetNContributors();
470 delete tpcvc; tpcvc=0;
474 AliESDVertex *trkvnc = ReconstructPrimaryVertexITSTPC(kFALSE);
475 xnt[index++]=trkvnc->GetXv();
476 xnt[index++]=trkvnc->GetXRes();
477 xnt[index++]=trkvnc->GetYv();
478 xnt[index++]=trkvnc->GetYRes();
479 xnt[index++]=trkvnc->GetZv();
480 xnt[index++]=trkvnc->GetZRes();
481 xnt[index++]=trkvnc->GetNContributors();
483 secnt[indexSecNt++]=trkvnc->GetXv();
484 secnt[indexSecNt++]=trkvnc->GetYv();
485 secnt[indexSecNt++]=trkvnc->GetZv();
486 secnt[indexSecNt++]=trkvnc->GetNContributors();
488 xTRKnc = (Float_t)trkvnc->GetXv();
489 yTRKnc = (Float_t)trkvnc->GetYv();
490 zTRKnc = (Float_t)trkvnc->GetZv();
491 ntrksTRKnc = (trkvnc->GetNContributors()<0 ? 0 : (UShort_t)trkvnc->GetNContributors());
493 delete trkvnc; trkvnc=0;
496 AliESDVertex *trkvc = ReconstructPrimaryVertexITSTPC(kTRUE);
497 xnt[index++]=trkvc->GetXv();
498 xnt[index++]=trkvc->GetXRes();
499 xnt[index++]=trkvc->GetYv();
500 xnt[index++]=trkvc->GetYRes();
501 xnt[index++]=trkvc->GetZv();
502 xnt[index++]=trkvc->GetZRes();
503 xnt[index++]=trkvc->GetNContributors();
504 delete trkvc; trkvc=0;
507 if(fRecoVtxITSTPCHalfEvent) {
508 AliESDVertex *trkvncodd = ReconstructPrimaryVertexITSTPC(kFALSE,1);
509 AliESDVertex *trkvnceven = ReconstructPrimaryVertexITSTPC(kFALSE,2);
510 if(trkvncodd->GetNContributors()>0 &&
511 trkvnceven->GetNContributors()>0) {
512 xnt[index++]=trkvncodd->GetXv()-trkvnceven->GetXv();
513 xnt[index++]=trkvncodd->GetYv()-trkvnceven->GetYv();
514 xnt[index++]=trkvncodd->GetZv()-trkvnceven->GetZv();
515 xnt[index++]=TMath::Sqrt(trkvncodd->GetXRes()*trkvncodd->GetXRes()+trkvnceven->GetXRes()*trkvnceven->GetXRes());
516 xnt[index++]=TMath::Sqrt(trkvncodd->GetYRes()*trkvncodd->GetYRes()+trkvnceven->GetYRes()*trkvnceven->GetYRes());
517 xnt[index++]=TMath::Sqrt(trkvncodd->GetZRes()*trkvncodd->GetZRes()+trkvnceven->GetZRes()*trkvnceven->GetZRes());
518 xnt[index++]=trkvnceven->GetNContributors();
519 xnt[index++]=trkvncodd->GetNContributors();
530 delete trkvncodd; trkvncodd=0;
531 delete trkvnceven; trkvnceven=0;
535 if(index>isize) printf("AliAnalysisTaskVertexESD: ERROR, index!=isize\n");
536 if(fFillNtuple) fNtupleVertexESD->Fill(xnt);
538 if(indexSecNt>isizeSecNt) printf("AliAnalysisTaskVertexESD: ERROR, indexSecNt!=isizeSecNt\n");
540 // if(indexTree>isizeTree) printf("AliAnalysisTaskVertexESD: ERROR, indexTree!=isizeTree\n");
541 // only every second event (to reduce output size)
542 if(fFillTreeBeamSpot && (nInFile%2 == 0)) fTreeBeamSpot->Fill();
545 // Post the data already here
546 PostData(1, fOutput);
551 //________________________________________________________________________
552 void AliAnalysisTaskVertexESD::Terminate(Option_t *)
554 // Draw result to the screen
555 // Called once at the end of the query
556 fOutput = dynamic_cast<TList*> (GetOutputData(1));
558 Printf("ERROR: fOutput not available");
562 if (!fNtupleVertexESD){
563 Printf("ERROR: fNtuple not available");
567 fNtupleVertexESD = dynamic_cast<TNtuple*>(fOutput->FindObject("fNtupleVertexESD"));
573 //_________________________________________________________________________
574 AliESDVertex* AliAnalysisTaskVertexESD::ReconstructPrimaryVertexTPC(Bool_t constr) const {
575 // On the fly reco of TPC vertex from ESD
576 AliESDEvent* evt = (AliESDEvent*) fInputEvent;
577 AliVertexerTracks vertexer(evt->GetMagneticField());
578 vertexer.SetTPCMode(); // defaults
579 Float_t diamondcovxy[3]; evt->GetDiamondCovXY(diamondcovxy);
580 Double_t pos[3]={evt->GetDiamondX(),evt->GetDiamondY(),0};
581 Double_t cov[6]={diamondcovxy[0],diamondcovxy[1],diamondcovxy[2],0.,0.,10.*10.};
582 AliESDVertex *initVertex = new AliESDVertex(pos,cov,1.,1);
583 vertexer.SetVtxStart(initVertex);
585 if(!constr) vertexer.SetConstraintOff();
587 return vertexer.FindPrimaryVertex(evt);
590 //_________________________________________________________________________
591 AliESDVertex* AliAnalysisTaskVertexESD::ReconstructPrimaryVertexITSTPC(Bool_t constr,Int_t mode) const {
592 // On the fly reco of ITS+TPC vertex from ESD
593 // mode=0 use all tracks
594 // mode=1 use odd-number tracks
595 // mode=2 use even-number tracks
597 AliESDEvent* evt = (AliESDEvent*) fInputEvent;
598 AliVertexerTracks vertexer(evt->GetMagneticField());
599 vertexer.SetITSMode(); // defaults
600 vertexer.SetMinClusters(4); // default is 5
601 //vertexer.SetITSpureSA();
602 Float_t diamondcovxy[3]; evt->GetDiamondCovXY(diamondcovxy);
603 Double_t pos[3]={evt->GetDiamondX(),evt->GetDiamondY(),0};
604 Double_t cov[6]={diamondcovxy[0],diamondcovxy[1],diamondcovxy[2],0.,0.,10.*10.};
605 AliESDVertex *initVertex = new AliESDVertex(pos,cov,1.,1);
606 vertexer.SetVtxStart(initVertex);
608 if(!constr) vertexer.SetConstraintOff();
610 // use only ITS-TPC or only ITS-SA tracks
611 if(fOnlyITSTPCTracks || fOnlyITSSATracks || mode>0) {
613 Int_t *skip = new Int_t[evt->GetNumberOfTracks()];
614 for(Int_t itr=0;itr<evt->GetNumberOfTracks(); itr++) {
615 AliESDtrack* track = evt->GetTrack(itr);
616 if(fOnlyITSTPCTracks && track->GetNcls(1)==0) { // skip ITSSA
620 if(fOnlyITSSATracks && track->GetNcls(1)>0) { // skip ITSTPC
624 if(mode==1 && itr%2==0) skip[iskip++]=itr;
625 if(mode==2 && itr%2==1) skip[iskip++]=itr;
627 vertexer.SetSkipTracks(iskip,skip);
628 delete [] skip; skip=NULL;
631 return vertexer.FindPrimaryVertex(evt);