1 /**************************************************************************
2 * Copyright(c) 1998-2007, 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 **************************************************************************/
18 //-------------------------------------------------------------------------
20 // Author: Markus Oldenburg, CERN
21 //-------------------------------------------------------------------------
26 #include <TFriendElement.h>
27 #include <TProcessID.h>
28 #include <TCollection.h>
29 #include "Riostream.h"
30 #include "AliAODEvent.h"
31 #include "AliAODHeader.h"
32 #include "AliAODTrack.h"
33 #include "AliAODDimuon.h"
37 // definition of std AOD member names
38 const char* AliAODEvent::fAODListName[kAODListN] = {"header",
63 //______________________________________________________________________________
64 AliAODEvent::AliAODEvent() :
93 // default constructor
94 if (TClass::IsCallingNew() != TClass::kDummyNew) fAODObjects = new TList();
97 //______________________________________________________________________________
98 AliAODEvent::AliAODEvent(const AliAODEvent& aod):
100 fAODObjects(new TList()),
103 fHeader(new AliAODHeader(*aod.fHeader)),
104 fTracks(new TClonesArray(*aod.fTracks)),
105 fVertices(new TClonesArray(*aod.fVertices)),
106 fV0s(new TClonesArray(*aod.fV0s)),
107 fCascades(new TClonesArray(*aod.fCascades)),
108 fTracklets(new AliAODTracklets(*aod.fTracklets)),
109 fJets(new TClonesArray(*aod.fJets)),
110 fEmcalCells(new AliAODCaloCells(*aod.fEmcalCells)),
111 fPhosCells(new AliAODCaloCells(*aod.fPhosCells)),
112 fCaloClusters(new TClonesArray(*aod.fCaloClusters)),
113 fEMCALTrigger(new AliAODCaloTrigger(*aod.fEMCALTrigger)),
114 fPHOSTrigger(new AliAODCaloTrigger(*aod.fPHOSTrigger)),
115 fFmdClusters(new TClonesArray(*aod.fFmdClusters)),
116 fPmdClusters(new TClonesArray(*aod.fPmdClusters)),
117 fHMPIDrings(new TClonesArray(*aod.fHMPIDrings)),
118 fDimuons(new TClonesArray(*aod.fDimuons)),
119 fAODTZERO(new AliAODTZERO(*aod.fAODTZERO)),
120 fAODVZERO(new AliAODVZERO(*aod.fAODVZERO)),
121 fAODZDC(new AliAODZDC(*aod.fAODZDC)),
122 fTOFHeader(new AliTOFHeader(*aod.fTOFHeader))
124 ,fAODMFT(new AliAODMFT(*aod.fAODMFT))
130 AddObject(fVertices);
132 AddObject(fCascades);
133 AddObject(fTracklets);
135 AddObject(fEmcalCells);
136 AddObject(fPhosCells);
137 AddObject(fCaloClusters);
138 AddObject(fEMCALTrigger);
139 AddObject(fPHOSTrigger);
140 AddObject(fFmdClusters);
141 AddObject(fPmdClusters);
142 AddObject(fHMPIDrings);
144 AddObject(fAODTZERO);
145 AddObject(fAODVZERO);
147 AddObject(fTOFHeader);
149 AddObject(fAODVZERO);
151 fConnected = aod.fConnected;
156 //______________________________________________________________________________
157 AliAODEvent & AliAODEvent::operator=(const AliAODEvent& aod) {
159 // Assignment operator
161 if(&aod == this) return *this;
162 AliVEvent::operator=(aod);
164 // This assumes that the list is already created
165 // and that the virtual void Copy(Tobject&) function
166 // is correctly implemented in the derived class
167 // otherwise only TObject::Copy() will be used
169 if((fAODObjects->GetSize()==0)&&(aod.fAODObjects->GetSize()>=kAODListN)){
170 // We cover the case that we do not yet have the
171 // standard content but the source has it
175 // Here we have the standard content without user additions, but the content is
176 // not matching the aod source.
178 // Iterate the list of source objects
179 TIter next(aod.GetList());
182 while ((its = next())) {
183 name = its->GetName();
184 // Check if we have this object type in out list
185 TObject *mine = fAODObjects->FindObject(name);
187 // We have to create the same type of object.
188 TClass* pClass=TClass::GetClass(its->ClassName());
190 AliWarning(Form("Can not find class description for entry %s (%s)\n",
191 its->ClassName(), name.Data()));
194 mine=(TObject*)pClass->New();
196 // not in this: can be added to list
197 AliWarning(Form("%s:%d Could not find %s for copying \n",
198 (char*)__FILE__,__LINE__,name.Data()));
201 if(mine->InheritsFrom("TNamed")) {
202 ((TNamed*)mine)->SetName(name);
203 } else if(mine->InheritsFrom("TCollection")){
204 if(mine->InheritsFrom("TClonesArray")) {
205 TClonesArray *itscl = dynamic_cast<TClonesArray*>(its);
207 AliWarning(Form("Class description for entry %s (%s) not TClonesArray\n",
208 its->ClassName(), name.Data()));
212 dynamic_cast<TClonesArray*>(mine)->SetClass(itscl->GetClass(), itscl->GetSize());
214 dynamic_cast<TCollection*>(mine)->SetName(name);
216 AliDebug(1, Form("adding object %s of type %s", mine->GetName(), mine->ClassName()));
219 // Now we have an object of the same type and name, but different content.
220 if(!its->InheritsFrom("TCollection")){
221 // simple objects (do they have a Copy method that calls operator= ?)
223 } else if (its->InheritsFrom("TClonesArray")) {
224 // Create or expand the tclonesarray pointers
225 // so we can directly copy to the object
226 TClonesArray *its_tca = (TClonesArray*)its;
227 TClonesArray *mine_tca = (TClonesArray*)mine;
228 // this leaves the capacity of the TClonesArray the same
229 // except for a factor of 2 increase when size > capacity
230 // does not release any memory occupied by the tca
231 Int_t its_entries = its_tca->GetEntriesFast();
232 mine_tca->ExpandCreate(its_entries);
233 for(int i=0; i<its_entries; i++){
235 TObject *mine_tca_obj = mine_tca->At(i);
236 TObject *its_tca_obj = its_tca->At(i);
237 // no need to delete first
238 // pointers within the class should be handled by Copy()...
239 // Can there be Empty slots?
240 its_tca_obj->Copy(*mine_tca_obj);
243 AliWarning(Form("%s:%d cannot copy TCollection \n",
244 (char*)__FILE__,__LINE__));
247 fConnected = aod.fConnected;
252 //______________________________________________________________________________
253 AliAODEvent::~AliAODEvent()
259 // fAODObjects->Delete("slow");
264 //______________________________________________________________________________
265 void AliAODEvent::AddObject(TObject* obj)
267 // Add an object to the list of objects.
268 // Please be aware that in order to increase performance you should
269 // refrain from using TObjArrays (if possible). Use TClonesArrays, instead.
271 // if ( !fAODObjects ) {
272 // fAODObjects = new TList();
273 // fAODObjects->SetOwner();
275 if ( !fAODObjects->FindObject(obj) )
277 fAODObjects->AddLast(obj);
281 //______________________________________________________________________________
282 void AliAODEvent::RemoveObject(TObject* obj)
284 // Removes an object from the list of objects.
286 fAODObjects->Remove(obj);
289 //______________________________________________________________________________
290 TObject *AliAODEvent::FindListObject(const char *objName) const
292 // Return the pointer to the object with the given name.
294 return fAODObjects->FindObject(objName);
297 //______________________________________________________________________________
298 void AliAODEvent::CreateStdContent()
300 // create the standard AOD content and set pointers
302 // create standard objects and add them to the TList of objects
303 AddObject(new AliAODHeader());
304 AddObject(new TClonesArray("AliAODTrack", 0));
305 AddObject(new TClonesArray("AliAODVertex", 0));
306 AddObject(new TClonesArray("AliAODv0", 0));
307 AddObject(new TClonesArray("AliAODcascade", 0));
308 AddObject(new AliAODTracklets());
309 AddObject(new TClonesArray("AliAODJet", 0));
310 AddObject(new AliAODCaloCells());
311 AddObject(new AliAODCaloCells());
312 AddObject(new TClonesArray("AliAODCaloCluster", 0));
313 AddObject(new AliAODCaloTrigger()); // EMCAL
314 AddObject(new AliAODCaloTrigger()); // PHOS
315 AddObject(new TClonesArray("AliAODFmdCluster", 0));
316 AddObject(new TClonesArray("AliAODPmdCluster", 0));
317 AddObject(new TClonesArray("AliAODHMPIDrings", 0));
318 AddObject(new TClonesArray("AliAODDimuon", 0));
319 AddObject(new AliAODTZERO());
320 AddObject(new AliAODVZERO());
321 AddObject(new AliAODZDC());
322 AddObject(new AliTOFHeader());
324 AddObject(new AliAODMFT());
329 // read back pointers
335 void AliAODEvent::MakeEntriesReferencable()
337 // Make all entries referencable in a subsequent process
339 TIter next(fAODObjects);
341 while ((obj = next()))
343 if(obj->InheritsFrom("TCollection"))
345 AssignIDtoCollection((TCollection*)obj);
350 //______________________________________________________________________________
351 void AliAODEvent::SetStdNames()
353 // introduce the standard naming
355 if(fAODObjects->GetEntries()==kAODListN){
356 for(int i = 0;i < fAODObjects->GetEntries();i++){
357 TObject *fObj = fAODObjects->At(i);
358 if(fObj->InheritsFrom("TNamed")){
359 ((TNamed*)fObj)->SetName(fAODListName[i]);
361 else if(fObj->InheritsFrom("TClonesArray")){
362 ((TClonesArray*)fObj)->SetName(fAODListName[i]);
367 printf("%s:%d SetStdNames() Wrong number of Std Entries \n",(char*)__FILE__,__LINE__);
371 void AliAODEvent::CreateStdFolders()
373 // Create the standard folder structure
374 if(fAODFolder)delete fAODFolder;
375 fAODFolder = gROOT->GetRootFolder()->AddFolder("AOD", "AOD");
376 if(fAODObjects->GetEntries()==kAODListN){
377 for(int i = 0;i < fAODObjects->GetEntries();i++){
378 TObject *fObj = fAODObjects->At(i);
379 if(fObj->InheritsFrom("TClonesArray")){
380 fAODFolder->AddFolder(fAODListName[i], fAODListName[i], (TCollection*) fObj);
382 fAODFolder->AddFolder(fAODListName[i], fAODListName[i], 0);
387 printf("%s:%d CreateStdFolders() Wrong number of Std Entries \n",(char*)__FILE__,__LINE__);
391 //______________________________________________________________________________
392 void AliAODEvent::GetStdContent()
394 // set pointers for standard content
396 fHeader = (AliAODHeader*)fAODObjects->FindObject("header");
397 fTracks = (TClonesArray*)fAODObjects->FindObject("tracks");
398 fVertices = (TClonesArray*)fAODObjects->FindObject("vertices");
399 fV0s = (TClonesArray*)fAODObjects->FindObject("v0s");
400 fCascades = (TClonesArray*)fAODObjects->FindObject("cascades");
401 fTracklets = (AliAODTracklets*)fAODObjects->FindObject("tracklets");
402 fJets = (TClonesArray*)fAODObjects->FindObject("jets");
403 fEmcalCells = (AliAODCaloCells*)fAODObjects->FindObject("emcalCells");
404 fPhosCells = (AliAODCaloCells*)fAODObjects->FindObject("phosCells");
405 fCaloClusters = (TClonesArray*)fAODObjects->FindObject("caloClusters");
406 fEMCALTrigger = (AliAODCaloTrigger*)fAODObjects->FindObject("emcalTrigger");
407 fPHOSTrigger = (AliAODCaloTrigger*)fAODObjects->FindObject("phosTrigger");
408 fFmdClusters = (TClonesArray*)fAODObjects->FindObject("fmdClusters");
409 fPmdClusters = (TClonesArray*)fAODObjects->FindObject("pmdClusters");
410 fHMPIDrings = (TClonesArray*)fAODObjects->FindObject("hmpidRings");
411 fDimuons = (TClonesArray*)fAODObjects->FindObject("dimuons");
412 fAODTZERO = (AliAODTZERO*)fAODObjects->FindObject("AliAODTZERO");
413 fAODVZERO = (AliAODVZERO*)fAODObjects->FindObject("AliAODVZERO");
414 fAODZDC = (AliAODZDC*)fAODObjects->FindObject("AliAODZDC");
415 fTOFHeader = (AliTOFHeader*)fAODObjects->FindObject("AliTOFHeader");
417 fAODMFT = (AliAODMFT*)fAODObjects->FindObject("AliAODMFT");
421 //______________________________________________________________________________
422 void AliAODEvent::ResetStd(Int_t trkArrSize,
425 Int_t cascadeArrSize,
430 Int_t hmpidRingsSize,
434 // deletes content of standard arrays and resets size
437 if (trkArrSize > fTracks->GetSize())
438 fTracks->Expand(trkArrSize);
442 if (vtxArrSize > fVertices->GetSize())
443 fVertices->Expand(vtxArrSize);
447 if (v0ArrSize > fV0s->GetSize())
448 fV0s->Expand(v0ArrSize);
452 if (cascadeArrSize > fCascades->GetSize())
453 fCascades->Expand(cascadeArrSize);
457 if (jetSize > fJets->GetSize())
458 fJets->Expand(jetSize);
461 fCaloClusters->Delete();
462 if (caloClusSize > fCaloClusters->GetSize())
463 fCaloClusters->Expand(caloClusSize);
466 fFmdClusters->Delete();
467 if (fmdClusSize > fFmdClusters->GetSize())
468 fFmdClusters->Expand(fmdClusSize);
471 fPmdClusters->Delete();
472 if (pmdClusSize > fPmdClusters->GetSize())
473 fPmdClusters->Expand(pmdClusSize);
476 fHMPIDrings->Delete();
477 if (hmpidRingsSize > fHMPIDrings->GetSize())
478 fHMPIDrings->Expand(hmpidRingsSize);
482 if (dimuonArrSize > fDimuons->GetSize())
483 fDimuons->Expand(dimuonArrSize);
486 fTracklets->DeleteContainer();
488 fPhosCells->DeleteContainer();
490 fEmcalCells->DeleteContainer();
493 fEMCALTrigger->DeAllocate();
495 fPHOSTrigger->DeAllocate();
499 void AliAODEvent::ClearStd()
501 // clears the standard arrays
507 fVertices ->Delete();
511 fCascades ->Delete();
513 fTracklets ->DeleteContainer();
517 fEmcalCells ->DeleteContainer();
519 fPhosCells ->DeleteContainer();
521 fCaloClusters ->Delete();
523 fFmdClusters ->Clear();
525 fPmdClusters ->Clear();
527 fHMPIDrings ->Clear();
532 fEMCALTrigger->DeAllocate();
534 fPHOSTrigger->DeAllocate();
537 //_________________________________________________________________
538 Int_t AliAODEvent::GetPHOSClusters(TRefArray *clusters) const
540 // fills the provided TRefArray with all found phos clusters
544 AliAODCaloCluster *cl = 0;
545 Bool_t first = kTRUE;
546 for (Int_t i = 0; i < GetNumberOfCaloClusters() ; i++) {
547 if ( (cl = GetCaloCluster(i)) ) {
550 new (clusters) TRefArray(TProcessID::GetProcessWithUID(cl));
554 //printf("IsPHOS cluster %d, E %2.3f Size: %d \n",i,cl->E(),clusters->GetEntriesFast());
558 return clusters->GetEntriesFast();
561 //_________________________________________________________________
562 Int_t AliAODEvent::GetEMCALClusters(TRefArray *clusters) const
564 // fills the provided TRefArray with all found emcal clusters
567 AliAODCaloCluster *cl = 0;
568 Bool_t first = kTRUE;
569 for (Int_t i = 0; i < GetNumberOfCaloClusters(); i++) {
570 if ( (cl = GetCaloCluster(i)) ) {
573 new (clusters) TRefArray(TProcessID::GetProcessWithUID(cl));
577 //printf("IsEMCal cluster %d, E %2.3f Size: %d \n",i,cl->E(),clusters->GetEntriesFast());
581 return clusters->GetEntriesFast();
585 //______________________________________________________________________________
586 Int_t AliAODEvent::GetMuonTracks(TRefArray *muonTracks) const
588 // fills the provided TRefArray with all found muon tracks
592 AliAODTrack *track = 0;
593 for (Int_t iTrack = 0; iTrack < GetNTracks(); iTrack++) {
594 track = GetTrack(iTrack);
595 if (track->IsMuonTrack()) {
596 muonTracks->Add(track);
600 return muonTracks->GetEntriesFast();
604 //______________________________________________________________________________
605 Int_t AliAODEvent::GetNumberOfMuonTracks() const
607 // get number of muon tracks
609 for (Int_t iTrack = 0; iTrack < GetNTracks(); iTrack++) {
610 if ((GetTrack(iTrack))->IsMuonTrack()) {
618 //______________________________________________________________________________
619 void AliAODEvent::ReadFromTree(TTree *tree, Option_t* opt /*= ""*/)
621 // Connects aod event to tree
624 AliWarning("Zero Pointer to Tree \n");
628 if(!tree->GetTree())tree->LoadTree(0);
630 // Try to find AliAODEvent
631 AliAODEvent *aodEvent = 0;
632 aodEvent = (AliAODEvent*)tree->GetTree()->GetUserInfo()->FindObject("AliAODEvent");
634 // This event is connected to the tree by definition, just say so
635 aodEvent->SetConnected();
636 // Check if already connected to tree
637 TList* connectedList = (TList*) (tree->GetUserInfo()->FindObject("AODObjectsConnectedToTree"));
638 if (connectedList && (!strcmp(opt, "reconnect"))) {
639 // If connected use the connected list of objects
640 if (fAODObjects != connectedList) {
642 fAODObjects = connectedList;
649 // prevent a memory leak when reading back the TList
650 // if (!(strcmp(opt, "reconnect"))) fAODObjects->Delete();
652 // create a new TList from the UserInfo TList...
653 // copy constructor does not work...
654 // fAODObjects = (TList*)(aodEvent->GetList()->Clone());
655 fAODObjects = (TList*)aodEvent->GetList();
656 fAODObjects->SetOwner(kTRUE);
657 if(fAODObjects->GetEntries()<kAODListN)
659 AliWarning(Form("AliAODEvent::ReadFromTree() TList contains less than the standard contents %d < %d"
660 " That might be fine though (at least for filtered AODs)",fAODObjects->GetEntries(),kAODListN));
663 // Let's find out whether we have friends
664 TList* friendL = tree->GetTree()->GetListOfFriends();
669 while ((fe = (TFriendElement*)next())){
670 aodEvent = (AliAODEvent*)(fe->GetTree()->GetUserInfo()->FindObject("AliAODEvent"));
672 printf("No UserInfo on tree \n");
675 // TList* objL = (TList*)(aodEvent->GetList()->Clone());
676 TList* objL = (TList*)aodEvent->GetList();
677 printf("Get list of object from tree %d !!\n", objL->GetEntries());
678 TIter nextobject(objL);
680 while((obj = nextobject()))
682 printf("Adding object from friend %s !\n", obj->GetName());
683 fAODObjects->Add(obj);
684 } // object "branch" loop
688 // set the branch addresses
689 TIter next(fAODObjects);
691 while((el=(TNamed*)next())){
692 TString bname(el->GetName());
693 // check if branch exists under this Name
694 TBranch *br = tree->GetTree()->GetBranch(bname.Data());
696 tree->SetBranchAddress(bname.Data(),fAODObjects->GetObjectRef(el));
698 br = tree->GetBranch(Form("%s.",bname.Data()));
700 tree->SetBranchAddress(Form("%s.",bname.Data()),fAODObjects->GetObjectRef(el));
703 printf("%s %d AliAODEvent::ReadFromTree() No Branch found with Name %s. \n",
704 (char*)__FILE__,__LINE__,bname.Data());
709 // when reading back we are not owner of the list
710 // must not delete it
711 fAODObjects->SetOwner(kTRUE);
712 fAODObjects->SetName("AODObjectsConnectedToTree");
713 // we are not owner of the list objects
714 // must not delete it
715 tree->GetUserInfo()->Add(fAODObjects);
719 // we can't get the list from the user data, create standard content
720 // and set it by hand
722 TIter next(fAODObjects);
724 while((el=(TNamed*)next())){
725 TString bname(el->GetName());
726 tree->SetBranchAddress(bname.Data(),fAODObjects->GetObjectRef(el));
729 // when reading back we are not owner of the list
730 // must not delete it
731 fAODObjects->SetOwner(kTRUE);
734 //______________________________________________________________________________
735 Int_t AliAODEvent::GetNumberOfPileupVerticesSPD() const{
736 // count number of SPD pileup vertices
737 Int_t nVertices=GetNumberOfVertices();
738 Int_t nPileupVertices=0;
739 for(Int_t iVert=0; iVert<nVertices; iVert++){
740 AliAODVertex *v=GetVertex(iVert);
741 if(v->GetType()==AliAODVertex::kPileupSPD) nPileupVertices++;
743 return nPileupVertices;
745 //______________________________________________________________________________
746 Int_t AliAODEvent::GetNumberOfPileupVerticesTracks() const{
747 // count number of track pileup vertices
748 Int_t nVertices=GetNumberOfVertices();
749 Int_t nPileupVertices=0;
750 for(Int_t iVert=0; iVert<nVertices; iVert++){
751 AliAODVertex *v=GetVertex(iVert);
752 if(v->GetType()==AliAODVertex::kPileupTracks) nPileupVertices++;
754 return nPileupVertices;
756 //______________________________________________________________________________
757 AliAODVertex* AliAODEvent::GetPrimaryVertexSPD() const{
758 // Get SPD primary vertex
759 Int_t nVertices=GetNumberOfVertices();
760 for(Int_t iVert=0; iVert<nVertices; iVert++){
761 AliAODVertex *v=GetVertex(iVert);
762 if(v->GetType()==AliAODVertex::kMainSPD) return v;
766 //______________________________________________________________________________
767 AliAODVertex* AliAODEvent::GetPileupVertexSPD(Int_t iV) const{
768 // Get pile-up vertex iV
769 Int_t nVertices=GetNumberOfVertices();
771 for(Int_t iVert=0; iVert<nVertices; iVert++){
772 AliAODVertex *v=GetVertex(iVert);
773 if(v->GetType()==AliAODVertex::kPileupSPD){
774 if(counter==iV) return v;
780 //______________________________________________________________________________
781 AliAODVertex* AliAODEvent::GetPileupVertexTracks(Int_t iV) const{
782 // Get pile-up vertex iV
783 Int_t nVertices=GetNumberOfVertices();
785 for(Int_t iVert=0; iVert<nVertices; iVert++){
786 AliAODVertex *v=GetVertex(iVert);
787 if(v->GetType()==AliAODVertex::kPileupTracks){
788 if(counter==iV) return v;
794 //______________________________________________________________________________
795 Bool_t AliAODEvent::IsPileupFromSPD(Int_t minContributors,
797 Double_t nSigmaZdist,
798 Double_t nSigmaDiamXY,
799 Double_t nSigmaDiamZ) const{
801 // This function checks if there was a pile up
802 // reconstructed with SPD
804 AliAODVertex *mainV=GetPrimaryVertexSPD();
805 if(!mainV) return kFALSE;
806 Int_t nc1=mainV->GetNContributors();
807 if(nc1<1) return kFALSE;
808 Int_t nPileVert=GetNumberOfPileupVerticesSPD();
809 if(nPileVert==0) return kFALSE;
810 Int_t nVertices=GetNumberOfVertices();
812 for(Int_t iVert=0; iVert<nVertices; iVert++){
813 AliAODVertex *pv=GetVertex(iVert);
814 if(pv->GetType()!=AliAODVertex::kPileupSPD) continue;
815 Int_t nc2=pv->GetNContributors();
816 if(nc2>=minContributors){
817 Double_t z1=mainV->GetZ();
818 Double_t z2=pv->GetZ();
819 Double_t distZ=TMath::Abs(z2-z1);
820 Double_t distZdiam=TMath::Abs(z2-GetDiamondZ());
821 Double_t cutZdiam=nSigmaDiamZ*TMath::Sqrt(GetSigma2DiamondZ());
822 if(GetSigma2DiamondZ()<0.0001)cutZdiam=99999.; //protection for missing z diamond information
823 if(distZ>minZdist && distZdiam<cutZdiam){
824 Double_t x2=pv->GetX();
825 Double_t y2=pv->GetY();
826 Double_t distXdiam=TMath::Abs(x2-GetDiamondX());
827 Double_t distYdiam=TMath::Abs(y2-GetDiamondY());
828 Double_t cov1[6],cov2[6];
829 mainV->GetCovarianceMatrix(cov1);
830 pv->GetCovarianceMatrix(cov2);
831 Double_t errxDist=TMath::Sqrt(cov2[0]+GetSigma2DiamondX());
832 Double_t erryDist=TMath::Sqrt(cov2[2]+GetSigma2DiamondY());
833 Double_t errzDist=TMath::Sqrt(cov1[5]+cov2[5]);
834 Double_t cutXdiam=nSigmaDiamXY*errxDist;
835 if(GetSigma2DiamondX()<0.0001)cutXdiam=99999.; //protection for missing diamond information
836 Double_t cutYdiam=nSigmaDiamXY*erryDist;
837 if(GetSigma2DiamondY()<0.0001)cutYdiam=99999.; //protection for missing diamond information
838 if( (distXdiam<cutXdiam) && (distYdiam<cutYdiam) && (distZ>nSigmaZdist*errzDist) ){
847 //______________________________________________________________________________
848 void AliAODEvent::Print(Option_t *) const
850 // Print the names of the all branches
851 TIter next(fAODObjects);
853 Printf(">>>>> AOD Content <<<<<");
854 while((el=(TNamed*)next())){
855 Printf(">> %s ",el->GetName());
857 Printf(">>>>> <<<<<");
862 void AliAODEvent::AssignIDtoCollection(const TCollection* col)
864 // Static method which assigns a ID to each object in a collection
865 // In this way the objects are marked as referenced and written with
866 // an ID. This has the advantage that TRefs to this objects can be
867 // written by a subsequent process.
870 while ((obj = next()))
871 TProcessID::AssignID(obj);
874 Bool_t AliAODEvent::IsPileupFromSPDInMultBins() const {
875 Int_t nTracklets=GetTracklets()->GetNumberOfTracklets();
876 if(nTracklets<20) return IsPileupFromSPD(3,0.8);
877 else if(nTracklets<50) return IsPileupFromSPD(4,0.8);
878 else return IsPileupFromSPD(5,0.8);
881 void AliAODEvent::Reset()
884 // Std content + Non std content
888 if(fAODObjects->GetSize()>kAODListN){
889 // we have non std content
890 // this also covers aodfriends
891 for(int i = kAODListN;i < fAODObjects->GetSize();++i){
892 TObject *pObject = fAODObjects->At(i);
894 if(pObject->InheritsFrom(TClonesArray::Class())){
895 ((TClonesArray*)pObject)->Delete();
897 else if(!pObject->InheritsFrom(TCollection::Class())){
898 TClass *pClass = TClass::GetClass(pObject->ClassName());
899 if (pClass && pClass->GetListOfMethods()->FindObject("Clear")) {
900 AliDebug(1, Form("Clear for object %s class %s", pObject->GetName(), pObject->ClassName()));
904 AliDebug(1, Form("ResetWithPlacementNew for object %s class %s", pObject->GetName(), pObject->ClassName()));
905 Long_t dtoronly = TObject::GetDtorOnly();
906 TObject::SetDtorOnly(pObject);
908 pClass->New(pObject);
909 TObject::SetDtorOnly((void*)dtoronly);
913 AliWarning(Form("No reset for %s \n",
914 pObject->ClassName()));
920 Float_t AliAODEvent::GetVZEROEqMultiplicity(Int_t i) const
922 // Get VZERO Multiplicity for channel i
923 // Themethod uses the equalization factors
924 // stored in the ESD-run object in order to
925 // get equal multiplicities within a VZERO rins (1/8 of VZERO)
926 if (!fAODVZERO || !fHeader) return -1;
929 Float_t factorSum = 0;
930 for(Int_t j = 8*ring; j < (8*ring+8); ++j) {
931 factorSum += fHeader->GetVZEROEqFactors(j);
933 Float_t factor = fHeader->GetVZEROEqFactors(i)*8./factorSum;
935 return (fAODVZERO->GetMultiplicity(i)/factor);
938 //------------------------------------------------------------
939 void AliAODEvent::SetTOFHeader(const AliTOFHeader *header)
942 // Set the TOF event_time
947 //fTOFHeader->SetName(fgkESDListName[kTOFHeader]);
950 // for analysis of reconstructed events
951 // when this information is not avaliable
952 fTOFHeader = new AliTOFHeader(*header);
953 //AddObject(fTOFHeader);
957 //------------------------------------------------------------
958 AliAODHMPIDrings *AliAODEvent::GetHMPIDringForTrackID(Int_t trackID) const
961 // Returns the HMPID object if any for a given track ID
965 for(Int_t ien = 0 ; ien < GetNHMPIDrings(); ien++)
967 if( GetHMPIDring(ien)->GetHmpTrkID() == trackID ) return GetHMPIDring(ien);
972 //------------------------------------------------------------
973 Int_t AliAODEvent::GetNHMPIDrings() const
976 // If there is a list of HMPID rings in the given AOD event, return their number
978 if ( fHMPIDrings) return fHMPIDrings->GetEntriesFast();
981 //------------------------------------------------------------
982 AliAODHMPIDrings *AliAODEvent::GetHMPIDring(Int_t nRings) const
985 // If there is a list of HMPID rings in the given AOD event, return corresponding ring
988 if( (AliAODHMPIDrings*)fHMPIDrings->UncheckedAt(nRings) ) {
989 return (AliAODHMPIDrings*)fHMPIDrings->UncheckedAt(nRings);
995 //------------------------------------------------------------