]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/TRD/AliHLTTRDEsdWriterComponent.cxx
updating HLT TRD ESDWriter follow reorganization of the friends in AliESDs.root ...
[u/mrichter/AliRoot.git] / HLT / TRD / AliHLTTRDEsdWriterComponent.cxx
1 /**************************************************************************
2  * This file is property of and copyright by the ALICE HLT Project        * 
3  * ALICE Experiment at CERN, All rights reserved.                         *
4  *                                                                        *
5  * Primary Authors: Matthias Richter <Matthias.Richter@ift.uib.no>        *
6  *                  for The ALICE HLT Project.                            *
7  *                                                                        *
8  * Permission to use, copy, modify and distribute this software and its   *
9  * documentation strictly for non-commercial purposes is hereby granted   *
10  * without fee, provided that the above copyright notice appears in all   *
11  * copies and that both the copyright notice and this permission notice   *
12  * appear in the supporting documentation. The authors make no claims     *
13  * about the suitability of this software for any purpose. It is          *
14  * provided "as is" without express or implied warranty.                  *
15  **************************************************************************/
16
17 /** @file   AliHLTTRDEsdWriterComponent.cxx
18     @author Mateusz Ploskon
19     @date   
20     @brief  Writer component to store tracks of the HLT TRD
21
22                                                                           */
23 #include "AliHLTTRDEsdWriterComponent.h"
24 #include "AliHLTTRDDefinitions.h"
25 #include "AliHLTTRDUtils.h"
26 #include "AliHLTCTPData.h"
27 #include "AliESDEvent.h"
28 #include "AliESDfriend.h"
29 #include "AliESDtrack.h"
30 #include "AliTRDtrackV1.h"
31 #include "TTree.h"
32 #include "TFile.h"
33
34 /** ROOT macro for the implementation of ROOT specific class methods */
35 ClassImp(AliHLTTRDEsdWriterComponent)
36
37 AliHLTTRDEsdWriterComponent::AliHLTTRDEsdWriterComponent()
38 :AliHLTRootFileWriterComponent()
39   ,fTree(NULL)
40   ,fFrTree(NULL)
41   ,fESD(NULL)
42   ,fESDfriend(NULL)
43   ,fFile(NULL)
44   ,fFrFile(NULL)
45   ,fTracksArray(NULL)
46 {
47   // see header file for class documentation
48   // or
49   // refer to README to build package
50   // or
51   // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt
52 }
53
54 AliHLTTRDEsdWriterComponent::AliHLTTRDEsdWriterComponent(const AliHLTTRDEsdWriterComponent&)
55   :AliHLTRootFileWriterComponent()
56   ,fTree(NULL)
57   ,fFrTree(NULL)
58   ,fESD(NULL)
59   ,fESDfriend(NULL)
60   ,fFile(NULL)
61   ,fTracksArray(NULL)
62 {
63 }
64
65 void AliHLTTRDEsdWriterComponent::GetInputDataTypes( vector<AliHLTComponent_DataType>& list)
66 {
67   // Get the list of input data  
68   list.clear(); // We do not have any requirements for our input data type(s).
69   list.push_back( AliHLTTRDDefinitions::fgkTracksDataType );
70   list.push_back( AliHLTTRDDefinitions::fgkHiLvlTracksDataType );
71 }
72
73 AliHLTTRDEsdWriterComponent::~AliHLTTRDEsdWriterComponent()
74 {
75   // see header file for class documentation
76 }
77
78 int AliHLTTRDEsdWriterComponent::InitWriter()
79 {
80   // see header file for class documentation
81   
82   fFile = new TFile("AliHLTTRDESDs.root", "recreate");
83   fESD = new AliESDEvent;
84   fESD->CreateStdContent();
85   fTree = new TTree("esdTree", "Tree with HLT::TRD ESD objects");
86   fESD->WriteToTree(fTree);
87   fFrFile = new TFile("AliHLTTRDESDfriends.root", "recreate");
88   fESDfriend = new AliESDfriend();
89   fFrTree = new TTree("esdFriendTree", "Tree with HLT::TRD ESD Friend objects");
90   fFrTree->Branch("ESDfriend.","AliESDfriend", &fESDfriend);
91   fESD->AddObject(fESDfriend);
92   fFile->cd();
93   fTree->GetUserInfo()->Add(fESD);
94   fTracksArray = new TClonesArray("AliTRDtrackV1");
95
96   SetupCTPData();
97
98   return 0;
99 }
100
101 int AliHLTTRDEsdWriterComponent::CloseWriter()
102 {
103   // see header file for class documentation
104
105   //fTree->Print();
106   fFile->cd();
107   fTree->Write(fTree->GetName(),TObject::kOverwrite);
108   fFile->Write();
109   fFile->Close();
110   fFrFile->cd();
111   fFrTree->Write(fFrTree->GetName(),TObject::kOverwrite);
112   fFrFile->Write();
113   fFrFile->Close();
114   delete fFile; fFile=0;
115   delete fFrFile; fFrFile=0;
116   //delete fTree;
117   delete fTracksArray; fTracksArray=0;
118
119   return AliHLTRootFileWriterComponent::CloseWriter();
120 }
121
122 int AliHLTTRDEsdWriterComponent::DumpEvent( const AliHLTComponentEventData& /*evtData*/,
123                                             const AliHLTComponentBlockData* /*blocks*/, 
124                                             AliHLTComponentTriggerData& trigData )
125 {
126   TClonesArray* TCAarray[18] = {0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0};
127   Int_t usedEntries = 0;
128   Int_t blockOrObject = 0;
129
130   for (const AliHLTComponentBlockData* pBlock=GetFirstInputBlock(AliHLTTRDDefinitions::fgkTracksDataType); pBlock; pBlock=GetNextInputBlock()) 
131     {
132       TCAarray[0] = fTracksArray;
133       AliHLTTRDUtils::ReadTracks(TCAarray[0], pBlock->fPtr, pBlock->fSize);
134       usedEntries = 1;
135       blockOrObject = -1;
136     }
137
138   for(const TObject *iter = GetFirstInputObject(AliHLTTRDDefinitions::fgkHiLvlTracksDataType); iter; iter = GetNextInputObject()) 
139     {
140       if(blockOrObject<0){
141         HLTError("You may not mix high level and low level!");
142         return -1;
143       }
144
145       TCAarray[usedEntries] = dynamic_cast<TClonesArray*>(const_cast<TObject*>(iter));
146       if(!TCAarray[usedEntries])continue;
147       usedEntries++;
148       blockOrObject = 1;
149     }
150
151   if(!blockOrObject)
152     return 0;
153
154   fESD->Reset(); 
155   fESD->SetMagneticField(GetBz());
156   fESD->SetRunNumber(GetRunNo());
157   fESD->SetPeriodNumber(GetPeriodNumber());
158   fESD->SetOrbitNumber(GetOrbitNumber());
159   fESD->SetBunchCrossNumber(GetBunchCrossNumber());
160   fESD->SetTimeStamp(GetTimeStamp());
161   fESD->SetEventType(7);
162
163   const AliHLTCTPData* pCTPData=CTPData();
164   if (pCTPData) {
165     AliHLTUInt64_t mask=pCTPData->ActiveTriggers(trigData);
166     for (int index=0; index<gkNCTPTriggerClasses; index++) {
167       if ((mask&((AliHLTUInt64_t)0x1<<index)) == 0) continue;
168       fESD->SetTriggerClass(pCTPData->Name(index), index);
169     }
170     fESD->SetTriggerMask(mask);
171   }
172   
173   for(int i=0; i<usedEntries; i++){
174     const TClonesArray* inArr = TCAarray[i];
175     for(int ii=0; ii<inArr->GetEntriesFast(); ii++){
176       AliTRDtrackV1* inV1 = (AliTRDtrackV1*)inArr->UncheckedAt(ii);
177       AliESDtrack *esdTrack = new AliESDtrack();
178       esdTrack->UpdateTrackParams(inV1, AliESDtrack::kTRDout);
179       esdTrack->SetLabel(inV1->GetLabel());
180       inV1->UpdateESDtrack(esdTrack);
181       AliTRDtrackV1 *calibTrack = new AliTRDtrackV1(*inV1);
182       calibTrack->SetOwner();
183       esdTrack->AddCalibObject(calibTrack);
184       fESD->AddTrack(esdTrack);
185       delete esdTrack;
186     }
187   }
188   
189   fESD->GetESDfriend(fESDfriend);
190   Int_t nb = fTree->Fill();
191   HLTInfo("Tree filled with %i bytes\n", nb);  
192   nb = fFrTree->Fill();
193   HLTInfo("FrTree filled with %i bytes\n", nb);
194   fESD->Reset();
195   fESDfriend->~AliESDfriend();
196   new (fESDfriend) AliESDfriend();
197
198   if(blockOrObject<0){
199     TCAarray[0]->Delete();
200   }
201
202   return 0;
203 }