]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/global/FlatESDConverter.C
give fileName of friengs to FlatESDConverter
[u/mrichter/AliRoot.git] / HLT / global / FlatESDConverter.C
1 /**
2  * >> Testing Macro to fill FlatESDEvent from ALiESDEvent <<
3  **
4  * Primary Authors : Sergey Gorbunov, Jochen Thaeder, Chiara Zampolli
5  *
6  * Usage:
7  *  aliroot -b -l -q LoadLibs.C FlatESDConverter.C++
8  *
9  **************************************************************************/
10
11 #if !defined(__CINT__) || defined(__MAKECINT__)
12 #include "AliESDEvent.h"
13 #include "AliESD.h"
14 #include "AliESDfriend.h"
15 #include <TFile.h>
16 #include <TTree.h>
17 #include <TSystem.h>
18 #include "./AliFlatESDEvent.h"
19 #include "./AliFlatESDTrack.h"
20 #include "./AliFlatTPCCluster.h"
21 #include "./AliFlatExternalTrackParam.h"
22 #include "Riostream.h"
23 #endif   
24
25 void FlatESDConverter(const char* filename="AliESDs.root", const char* filenameFriends="AliESDfriends.root",const char* filenameOut="out.dat", Bool_t useESDFriends = kTRUE, Bool_t useHLTtree = kFALSE) {
26   // -- Convert AliESDEvent to AliFlatESDEvent
27
28   ofstream outFile(Form("%s",filenameOut), std::ifstream::binary | std::ifstream::out);
29   //ofstream outFile("outFlatESD.dat");
30
31   TFile *file    = new TFile(Form("%s", filename));
32         
33         TTree *esdTree = useHLTtree? dynamic_cast<TTree*>(file->Get("HLTesdTree")) : dynamic_cast<TTree*>(file->Get("esdTree"));
34
35   
36   // -- Connect ESD
37   AliESDEvent *esd = new AliESDEvent;
38   esd->ReadFromTree(esdTree);
39   
40   // -- Connect ESD friend
41   AliESDfriend *esdFriend = NULL; 
42   if (useESDFriends && !esdTree->FindBranch("ESDfriend.")) {
43     esdTree->AddFriend("esdFriendTree", Form("%s", filenameFriends));
44     esdTree->SetBranchStatus("ESDfriend.", 1);
45
46     esdFriend = dynamic_cast<AliESDfriend*>((const_cast<AliESDEvent*>(esd))->FindListObject("AliESDfriend"));
47     if (esdFriend)
48       esdTree->SetBranchAddress("ESDfriend.", &esdFriend);
49   } // if (!esdTree->FindBranch("ESDfriend.")) {
50   
51   AliFlatESDEvent *flatEsd = NULL;
52
53   // -- Event Loop
54   for (Int_t idxEvent = 0; idxEvent < esdTree->GetEntries(); idxEvent++) {
55     esdTree->GetEntry(idxEvent);
56     // -- Book memory for AliFlatESDEvent
57    // -- TEST >>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
58
59     Byte_t *mem = new Byte_t[AliFlatESDEvent::EstimateSize(esd, useESDFriends)];
60
61     flatEsd = reinterpret_cast<AliFlatESDEvent*>(mem);    
62
63     // -- Fill AliFlatESDEvent
64     flatEsd->Fill(esd, useESDFriends);  
65  #if 0
66      Printf("TEST: Event %d || Tracks %d | FRIEND Tracks %d || estimated size %llu || sizeof(AliFlatESDEvent) %llu", 
67            idxEvent, esd->GetNumberOfTracks(), esdFriend->GetNumberOfTracks(), 
68            AliFlatESDEvent::EstimateSize(esd, useESDFriends), flatEsd->GetSize());
69
70
71     AliFlatESDTrack *track = flatEsd->GetTracks();
72     for (Int_t idxTrack = 0; idxTrack < flatEsd->GetNumberOfTracks(); ++idxTrack) {      
73        AliESDtrack *esdTrack = esd->GetTrack(idxTrack);
74        
75        if (track && !esdTrack) {
76          Printf("ERROR THIS SHOULD NOT HAPPEN AT ALL !!! TRACK %d HAS NO ESD TRACK!!!", idxTrack);
77          return;
78        }
79
80       if (track) {
81         AliFlatExternalTrackParam* exp1 = track->GetTrackParamCp();
82         AliFlatExternalTrackParam* exp2 = track->GetTrackParamIp();
83         AliFlatExternalTrackParam* exp3 = track->GetTrackParamTPCInner();
84         AliFlatExternalTrackParam* exp4 = track->GetTrackParamOp();
85
86         Float_t alphaFLAT[4] = {0., 0., 0., 0.};
87         if (exp1) alphaFLAT[0] = exp1->GetAlpha();
88         if (exp2) alphaFLAT[1] = exp2->GetAlpha();
89         if (exp3) alphaFLAT[2] = exp3->GetAlpha();
90         if (exp4) alphaFLAT[3] = exp4->GetAlpha();
91         
92         Float_t alphaOLD[4]  = {0., 0., 0., 0.};
93         if (esdTrack->GetConstrainedParam())  alphaOLD[0] = esdTrack->GetConstrainedParam()->GetAlpha();
94         if (esdTrack->GetInnerParam())        alphaOLD[1] = esdTrack->GetInnerParam()->GetAlpha();
95         if (esdTrack->GetTPCInnerParam())     alphaOLD[2] = esdTrack->GetTPCInnerParam()->GetAlpha();
96         if (esdTrack->GetOuterParam())        alphaOLD[3] = esdTrack->GetOuterParam()->GetAlpha();
97
98         Printf("  TEST: FlatTrack %d > FlatExternalTrackParam > %p %p %p %p", idxTrack, exp1, exp2, exp3, exp4);
99         Printf("  TEST: FlatTrack %d > Alpha %f %f %f %f", idxTrack, alphaFLAT[0], alphaFLAT[1], alphaFLAT[2], alphaFLAT[3]);
100         Printf("  TEST: Old Track %d > Alpha %f %f %f %f", idxTrack, alphaOLD[0], alphaOLD[1], alphaOLD[2], alphaOLD[3]);
101         Printf("  TEST: Diff      %d > Alpha %f %f %f %f", idxTrack, 
102                alphaFLAT[0]-alphaOLD[0], alphaFLAT[1]-alphaOLD[1], alphaFLAT[2]-alphaOLD[2], alphaFLAT[3]-alphaOLD[3]);
103
104         Printf("  TEST: FlatTrack %d has %d FlatClusters", idxTrack, track->GetNumberOfTPCClusters());
105
106         for (Int_t idxCluster = 0; idxCluster < track->GetNumberOfTPCClusters(); ++idxCluster)
107           Printf("    TEST: FlatTrack %d > FlatCluster %d has row %d", idxTrack, idxCluster, track->GetTPCCluster(idxCluster).GetPadRow());
108       }
109       track = track->GetNextTrack();
110     }
111     // -- <<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< TEST
112 #endif
113
114     outFile.write(reinterpret_cast<char*>(mem), flatEsd->GetSize());
115
116     delete[] mem; 
117
118   } // for (Int_t idxEvent = 1; idxEvent < 2; idxEvent++) {
119
120   outFile.close();
121
122   return;
123 }