]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/ZDC/AliHLTZDCESDRecoComponent.cxx
delete the temporary tree of clusters instead of Reset().
[u/mrichter/AliRoot.git] / HLT / ZDC / AliHLTZDCESDRecoComponent.cxx
1 //-*- Mode: C++ -*-
2 // $Id: AliHLTZDCESDRecoComponent.cxx $
3 /**************************************************************************
4  * This file is property of and copyright by the ALICE HLT Project        * 
5  * ALICE Experiment at CERN, All rights reserved.                         *
6  *                                                                        *
7  * Primary Authors: Chiara Oppedisano <Chiara.Oppedisano@to.infn.it>      *
8  *                  for The ALICE HLT Project.                            *
9  *                                                                        *
10  * Permission to use, copy, modify and distribute this software and its   *
11  * documentation strictly for non-commercial purposes is hereby granted   *
12  * without fee, provided that the above copyright notice appears in all   *
13  * copies and that both the copyright notice and this permission notice   *
14  * appear in the supporting documentation. The authors make no claims     *
15  * about the suitability of this software for any purpose. It is          *
16  * provided "as is" without express or implied warranty.                  *
17  **************************************************************************/
18
19 /** @file    AliHLTZDCESDRecoComponent.cxx
20     @author  Chiara Oppedisano <Chiara.Oppedisano@to.infn.it>
21     @brief   ZDC reconstruction component
22 */
23
24 #if __GNUC__>= 3
25 using namespace std;
26 #endif
27
28 #include "AliHLTErrorGuard.h"
29 #include "AliHLTSystem.h"
30 #include "AliHLTZDCESDRecoComponent.h"
31 #include "AliHLTDefinitions.h"
32 #include "AliRawReaderMemory.h"
33 #include "AliGeomManager.h"
34 #include "AliGRPObject.h"
35 #include "AliZDCReconstructor.h"
36 #include "AliZDCRecoParam.h"
37 #include "AliZDCRecoParampp.h"
38 #include "AliZDCRecoParamPbPb.h"
39 #include "AliESDEvent.h"
40 #include "AliESDZDC.h"
41 #include "AliCDBManager.h"
42 #include "AliHLTDataTypes.h"
43 #include <cstdlib>
44 #include <cerrno>
45
46 ClassImp(AliHLTZDCESDRecoComponent)
47     
48 //_____________________________________________________________________
49 AliHLTZDCESDRecoComponent::AliHLTZDCESDRecoComponent()
50   : AliHLTProcessor(),    
51     fRawReader(NULL),
52     fReconstructor(NULL)
53 {
54 }
55
56 //_____________________________________________________________________
57 AliHLTZDCESDRecoComponent::~AliHLTZDCESDRecoComponent()
58 {
59     // see header file for class documentation
60 }
61
62 //_____________________________________________________________________
63 const char* AliHLTZDCESDRecoComponent::GetComponentID()
64 {
65     // see header file for class documentation
66     return "ZDCESDReco"; // The ID of this component
67 }
68
69 //_____________________________________________________________________
70 void AliHLTZDCESDRecoComponent::GetInputDataTypes( vector<AliHLTComponentDataType>& list)
71 {
72     // see header file for class documentation
73       list.push_back(kAliHLTDataTypeDDLRaw|kAliHLTDataOriginZDC);
74 }
75
76 //_____________________________________________________________________
77 AliHLTComponentDataType AliHLTZDCESDRecoComponent::GetOutputDataType()
78 {
79       return kAliHLTDataTypeESDContent|kAliHLTDataOriginZDC;
80 }
81
82 //_____________________________________________________________________
83 void AliHLTZDCESDRecoComponent::GetOutputDataSize( unsigned long& constBase, double& inputMultiplier )
84 {
85     // see header file for class documentation
86     constBase = 1000;
87     inputMultiplier = 0;
88 }
89
90 //_____________________________________________________________________
91 AliHLTComponent* AliHLTZDCESDRecoComponent::Spawn()
92 {
93     // see header file for class documentation
94     return new AliHLTZDCESDRecoComponent;
95 }
96
97 //_____________________________________________________________________
98 int AliHLTZDCESDRecoComponent::DoInit( int argc, const char** argv )
99 {
100   // see header file for class documentation
101   int iResult=0;
102
103   // -- Load GeomManager
104   if(AliGeomManager::GetGeometry()==NULL){
105     AliGeomManager::LoadGeometry();
106   }
107   
108   // read configuration object : HLT/ConfigZDC/
109   TString cdbPath="HLT/ConfigZDC/";
110   cdbPath+=GetComponentID();
111   iResult=ConfigureFromCDBTObjString(cdbPath);
112
113   // init stage 3: read the component arguments
114   if (iResult>=0) {
115     iResult=ConfigureFromArgumentString(argc, argv);
116   }
117
118   // init stage 1: default values for all data members
119   
120   TObject* pOCDBEntry = LoadAndExtractOCDBObject("GRP/GRP/Data");
121   AliGRPObject* pGRP = pOCDBEntry?dynamic_cast<AliGRPObject*>(pOCDBEntry):NULL;
122   Float_t  beamEnergy=0.;
123   if(pGRP) beamEnergy = pGRP->GetBeamEnergy(); 
124
125   TString beamType="";
126   if(pGRP) beamType = pGRP->GetBeamType();
127
128   //!!! set the beam type and the energy explicitly
129
130   beamType="A-A"; 
131   beamEnergy = 2760.;
132
133   // implement the component initialization
134   do {
135     if (iResult<0) break;
136
137     fReconstructor = new AliZDCReconstructor;
138     if (!fReconstructor) {
139       iResult=-ENOMEM;
140       break;
141     }
142
143     fRawReader = new AliRawReaderMemory;
144     if (!fRawReader) {
145       iResult=-ENOMEM;
146       break;
147     }
148
149     // implement further initialization
150   } while (0);
151
152   if (iResult<0) {
153     // implement cleanup
154
155     if (fRawReader) delete fRawReader;
156     fRawReader = NULL;
157     
158     if (fReconstructor) delete fReconstructor;
159     fReconstructor = NULL;
160   }
161
162   if (iResult>=0) {
163
164    if((beamType.CompareTo("p-p"))==0) fReconstructor->SetRecoParam(AliZDCRecoParampp::GetLowFluxParam());
165     else if((beamType.CompareTo("A-A"))==0) fReconstructor->SetRecoParam(AliZDCRecoParamPbPb::GetHighFluxParam(beamEnergy));
166     else HLTWarning(" Beam type not known by ZDC!");
167   
168
169     fReconstructor->Init(beamType, beamEnergy/2);
170   }
171
172   return iResult;
173 }
174
175 //_____________________________________________________________________
176 int AliHLTZDCESDRecoComponent::ScanConfigurationArgument(int /*argc*/, const char** argv)
177 {
178   // Scan configuration arguments
179   // Return the number of processed arguments
180   //        -EPROTO if argument format error (e.g. number expected but not found)
181   //
182   // The AliHLTComponent base class implements a parsing loop for argument strings and
183   // arrays of strings which is invoked by ConfigureFromArgumentString/ConfigureFromCDBTObjString
184   // The component needs to implement ScanConfigurationArgument in order to decode the arguments.
185
186   int i=0;
187   TString argument=argv[i];
188
189   if (argument.IsNull()) return 0;
190
191   return 0;
192 }
193
194 //_____________________________________________________________________
195 int AliHLTZDCESDRecoComponent::DoDeinit()
196 {
197     if(fRawReader) delete fRawReader;
198     if(fReconstructor) delete fReconstructor;
199     return 0;
200 }
201
202 //_____________________________________________________________________
203 int AliHLTZDCESDRecoComponent::DoEvent(const AliHLTComponentEventData& /*evtData*/,
204                                        AliHLTComponentTriggerData& /*trigData*/)
205 {
206   // event processing function
207   int iResult=0;
208   
209   // check if this is a data event, there are a couple of special events
210   // which should be ignored for normal processing
211   if (!IsDataEvent()) return 0;
212   
213   // get ZDC raw input data block and set up the rawreader
214   const AliHLTComponentBlockData* pBlock = GetFirstInputBlock(kAliHLTDataTypeDDLRaw|kAliHLTDataOriginZDC);
215   if (!pBlock) {
216     HLTInfo("No ZDC input block !!!");
217     return 0;
218   }
219
220   // add input block to raw reader
221   if (!fRawReader->SetMemory((UChar_t*) pBlock->fPtr, pBlock->fSize )){
222     HLTError("Could not add buffer of data block  %s, 0x%08x to rawreader",
223              DataType2Text(pBlock->fDataType).c_str(), pBlock->fSpecification);
224     iResult = -1;
225   }
226   
227   TTree *clusterTree = new TTree();
228   if (!clusterTree) {
229     iResult=-ENOMEM;
230   }
231
232   if (iResult >= 0) {
233
234     // set ZDC EquipmentID
235     fRawReader->SetEquipmentID(3840);
236
237     fReconstructor->Reconstruct(fRawReader, clusterTree);
238
239     fReconstructor->FillZDCintoESD(clusterTree, (AliESDEvent *) NULL);
240
241     // send AliESDZDC
242     PushBack(static_cast<TObject*>(fReconstructor->GetZDCESDData()), 
243              kAliHLTDataTypeESDContent|kAliHLTDataOriginZDC,0);
244    
245   }
246
247   delete clusterTree;
248
249   // clear the rawreader
250   fRawReader->ClearBuffers();    
251   
252   return iResult;
253 }
254
255
256 //_____________________________________________________________________
257 int AliHLTZDCESDRecoComponent::Reconfigure(const char* cdbEntry, const char* chainId)
258 {
259   // reconfigure the component from the specified CDB entry, or default CDB entry
260   // function is invoked by the framework if a reconfigure command was received.
261   // 
262   int iResult=0;
263   TString cdbPath;
264   if (cdbEntry) {
265     cdbPath=cdbEntry;
266   } else {
267     cdbPath="HLT/ConfigZDC/";
268     cdbPath+=GetComponentID();
269   }
270   AliInfoClass(Form("reconfigure '%s' from entry %s%s", chainId, cdbPath.Data(), cdbEntry?"":" (default)"));
271   iResult=ConfigureFromCDBTObjString(cdbPath);
272
273   return iResult;
274 }
275
276
277 //_____________________________________________________________________
278 int AliHLTZDCESDRecoComponent::ReadPreprocessorValues(const char* /*modules*/)
279 {
280   // see header file for class documentation
281   ALIHLTERRORGUARD(5, "ReadPreProcessorValues not implemented for this component");
282   return 0;
283 }