]> git.uio.no Git - u/mrichter/AliRoot.git/blame - HLT/MUON/OfflineInterface/AliHLTMUONESDMaker.cxx
First version of the analysis QA
[u/mrichter/AliRoot.git] / HLT / MUON / OfflineInterface / AliHLTMUONESDMaker.cxx
CommitLineData
649ab027 1/**************************************************************************
2 * This file is property of and copyright by the ALICE HLT Project *
3 * All rights reserved. *
4 * *
5 * Primary Authors: *
6 * Artur Szostak <artursz@iafrica.com> *
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/* $Id: $ */
18
19///
20/// @file AliHLTMUONESDMaker.cxx
21/// @author Artur Szostak <artursz@iafrica.com>
22/// @date 30 June 2008
23/// @brief Implementation of the AliHLTMUONESDMaker component.
24///
25/// The ESD maker component converts dHLT raw internal reconstructed information
26/// into AliESDEvent objects.
27///
28
29#include "AliHLTMUONESDMaker.h"
30#include "AliHLTMUONEvent.h"
31#include "AliHLTMUONConstants.h"
32#include "AliHLTMUONUtils.h"
33#include "AliHLTMUONRecHit.h"
34#include "AliHLTMUONTriggerRecord.h"
35#include "AliHLTMUONMansoTrack.h"
36#include "AliHLTMUONDecision.h"
37#include "AliMUONConstants.h"
38#include "AliMUONVCluster.h"
39#include "AliESDEvent.h"
40#include "AliESDRun.h"
41#include "AliESDMuonTrack.h"
42#include "AliESDMuonCluster.h"
43#include "TClonesArray.h"
44#include <cmath>
45#include <cassert>
46
47
48ClassImp(AliHLTMUONESDMaker);
49
50
51AliHLTMUONESDMaker::AliHLTMUONESDMaker() :
52 AliHLTMUONProcessor(),
53 fWarnForUnexpecedBlock(false),
54 fMakeMinimalESD(false)
55{
56 /// Default constructor.
57}
58
59
60AliHLTMUONESDMaker::~AliHLTMUONESDMaker()
61{
62 /// Default destructor.
63}
64
65
66int AliHLTMUONESDMaker::DoInit(int argc, const char** argv)
67{
68 /// Inherited from AliHLTComponent.
69 /// Parses the command line parameters and initialises the component.
70
71 HLTInfo("Initialising dHLT ESD maker component.");
72
73 fWarnForUnexpecedBlock = false;
74 fMakeMinimalESD = false;
75
76 for (int i = 0; i < argc; i++)
77 {
78 if (strcmp(argv[i], "-make_minimal_esd") == 0)
79 {
80 fMakeMinimalESD = true;
81 continue;
82 }
83
84 if (strcmp(argv[i], "-warn_on_unexpected_block") == 0)
85 {
86 fWarnForUnexpecedBlock = true;
87 continue;
88 }
89
90 HLTError("Unknown option '%s'.", argv[i]);
91 return -EINVAL;
92 }
93
94 return 0;
95}
96
97
98int AliHLTMUONESDMaker::DoDeinit()
99{
100 /// Inherited from AliHLTComponent. Performs a cleanup of the component.
101
102 HLTInfo("Deinitialising dHLT ESD maker component.");
103 return 0;
104}
105
106
107const char* AliHLTMUONESDMaker::GetComponentID()
108{
109 /// Inherited from AliHLTComponent. Returns the component ID.
110
111 return AliHLTMUONConstants::ESDMakerId();
112}
113
114
115AliHLTComponentDataType AliHLTMUONESDMaker::GetOutputDataType()
116{
117 /// Inherited from AliHLTComponent.
118 /// Returns the ESD object data type with MUON origin.
119
120 return AliHLTMUONConstants::ESDDataType();
121}
122
123
124void AliHLTMUONESDMaker::GetInputDataTypes(
125 vector<AliHLTComponentDataType>& list
126 )
127{
128 /// Inherited from AliHLTProcessor.
129 /// Returns the list of expected input data types.
130
131 list.push_back(AliHLTMUONConstants::TriggerRecordsBlockDataType());
132 list.push_back(AliHLTMUONConstants::MansoTracksBlockDataType());
133}
134
135
136void AliHLTMUONESDMaker::GetOutputDataSize(
137 unsigned long& constBase, double& inputMultiplier
138 )
139{
140 /// Inherited from AliHLTComponent.
141 /// Returns an estimate of the expected output data size.
142
94135365 143 constBase = sizeof(AliESDEvent) + 1024*1024; // The extra 1 MByte is for auxilary objects created in AliESDEvent.
649ab027 144 inputMultiplier = 10;
145}
146
147
148AliHLTComponent* AliHLTMUONESDMaker::Spawn()
149{
150 /// Inherited from AliHLTComponent. Creates a new object instance.
151
152 return new AliHLTMUONESDMaker();
153}
154
155
156int AliHLTMUONESDMaker::DoEvent(
157 const AliHLTComponentEventData& /*evtData*/,
158 AliHLTComponentTriggerData& /*trigData*/
159 )
160{
161 /// Inherited from AliHLTProcessor. Processes the new event data.
162
163 AliESDEvent event;
164 AliHLTUInt32_t clusterIndex = 0; // for the cluster unique ID.
165
166 // Create and fill in the standard ESD objects or just create the muon
167 // tracks array if so requested.
168 if (fMakeMinimalESD)
169 {
170 TClonesArray* muonTracks = new TClonesArray("AliESDMuonTrack",0);
171 muonTracks->SetName("MuonTracks");
172 event.AddObject(muonTracks);
173 event.GetStdContent();
174 }
175 else
176 {
177 event.CreateStdContent();
178 event.SetRunNumber(GetRunNo());
179 }
180
181 const AliHLTComponentBlockData* block = NULL;
182 AliHLTUInt32_t specification = 0; // Contains the output data block spec bits.
183 std::vector<const AliHLTMUONTriggerRecordStruct*> triggerRecords;
184
185 // First process the blocks of trigger records. We simply mark each trigger
186 // record in the triggerRecords array.
187 for (int i = 0; i < GetNumberOfInputBlocks(); i++)
188 {
189 block = GetInputBlock(i);
190 assert( block != NULL );
191
192 HLTDebug("Handling block: %u, with fDataType = '%s', fPtr = %p and fSize = %u bytes.",
193 i, DataType2Text(block->fDataType).c_str(), block->fPtr, block->fSize
194 );
195
196 if (block->fDataType == AliHLTMUONConstants::TriggerRecordsBlockDataType())
197 {
198 specification |= block->fSpecification;
3240b3ce 199 AliHLTMUONTriggerRecordsBlockReader inblock(block->fPtr, block->fSize);
649ab027 200 if (not BlockStructureOk(inblock)) continue;
201
202 for (AliHLTUInt32_t n = 0; n < inblock.Nentries(); n++)
203 {
204 triggerRecords.push_back(&inblock[n]);
205 }
206 }
207 else
208 {
209 if (block->fDataType != AliHLTMUONConstants::MansoTracksBlockDataType())
210 {
211 // Log a message indicating that we got a data block that we
212 // do not know how to handle.
213 if (fWarnForUnexpecedBlock)
214 HLTWarning("Received a data block of a type we cannot handle: '%s', spec: 0x%X",
215 DataType2Text(block->fDataType).c_str(), block->fSpecification
216 );
217 else
218 HLTDebug("Received a data block of a type we cannot handle: '%s', spec: 0x%X",
219 DataType2Text(block->fDataType).c_str(), block->fSpecification
220 );
221 }
222 }
223 }
224
225 // Now we can look for tracks to add. We needed the ROOT trigger records
226 // and reco hits created before we can create track objects.
227 for (block = GetFirstInputBlock(AliHLTMUONConstants::MansoTracksBlockDataType());
228 block != NULL;
229 block = GetNextInputBlock()
230 )
231 {
232 specification |= block->fSpecification;
3240b3ce 233 AliHLTMUONMansoTracksBlockReader inblock(block->fPtr, block->fSize);
649ab027 234 if (not BlockStructureOk(inblock)) continue;
235
236 for (AliHLTUInt32_t n = 0; n < inblock.Nentries(); n++)
237 {
238 const AliHLTMUONMansoTrackStruct& t = inblock[n];
239 AliESDMuonTrack muTrack;
240
241 AliHLTMUONParticleSign sign;
242 bool hitset[4];
243 AliHLTMUONUtils::UnpackMansoTrackFlags(
244 t.fFlags, sign, hitset
245 );
246
247 double signVal = 0;
248 switch (sign)
249 {
250 case kSignMinus: signVal = +1.; break;
251 case kSignUnknown: signVal = 0.; break;
252 case kSignPlus: signVal = -1.; break;
253 default:
254 HLTWarning("Got a track with an invalid sign value: %d", int(sign));
255 }
256
257 TVector3 mom(t.fPx, t.fPy, t.fPz);
258 if (mom.Mag() != 0)
259 muTrack.SetInverseBendingMomentum(signVal/mom.Mag());
260 else
261 muTrack.SetInverseBendingMomentum(0.);
262 muTrack.SetThetaX(atan2(t.fPx, t.fPz));
263 muTrack.SetThetaY(atan2(t.fPy, t.fPz));
264 muTrack.SetZ(0.);
265 muTrack.SetBendingCoor(0.);
266 muTrack.SetNonBendingCoor(0.);
267
268 // The Manso algorithm assumes the information at the
269 // Distance of Closest Approach and chamber 1 is the same
270 // as the vertex.
271 if (mom.Mag() != 0)
272 muTrack.SetInverseBendingMomentumAtDCA(1./mom.Mag());
273 else
274 muTrack.SetInverseBendingMomentumAtDCA(0.);
275 muTrack.SetThetaXAtDCA(atan2(t.fPx, t.fPz));
276 muTrack.SetThetaYAtDCA(atan2(t.fPy, t.fPz));
277 muTrack.SetBendingCoorAtDCA(0.);
278 muTrack.SetNonBendingCoorAtDCA(0.);
279
280 if (mom.Mag() != 0)
281 muTrack.SetInverseBendingMomentumUncorrected(1./mom.Mag());
282 else
283 muTrack.SetInverseBendingMomentumUncorrected(0.);
284 muTrack.SetThetaXUncorrected(atan2(t.fPx, t.fPz));
285 muTrack.SetThetaYUncorrected(atan2(t.fPy, t.fPz));
286 muTrack.SetZUncorrected(0.);
287 muTrack.SetBendingCoorUncorrected(0.);
288 muTrack.SetNonBendingCoorUncorrected(0.);
289
290 muTrack.SetChi2(t.fChi2);
291
292 // Fill in the track hit points.
293 Int_t nHits = 0;
294 for (int i = 0; i < 4; i++)
295 {
296 if (not hitset[i]) continue;
297
298 AliHLTUInt8_t chamber;
299 AliHLTUInt16_t detElemId;
300 AliHLTMUONUtils::UnpackRecHitFlags(t.fHit[i].fFlags, chamber, detElemId);
301
302 AliESDMuonCluster cluster;
303 cluster.SetUniqueID(AliMUONVCluster::BuildUniqueID(chamber, detElemId, clusterIndex++));
304 cluster.SetXYZ(t.fHit[i].fX, t.fHit[i].fY, t.fHit[i].fZ);
305 cluster.SetErrXY( // Use nominal values.
306 AliMUONConstants::DefaultNonBendingReso(),
307 AliMUONConstants::DefaultBendingReso()
308 );
309 cluster.SetCharge(-1.); // Indicate no total charge calculated.
310 cluster.SetChi2(-1.); // Indicate no fit made.
311 muTrack.AddCluster(cluster);
312 nHits++;
313 muTrack.AddInMuonClusterMap(i+6);
314 }
315
316 // Find the corresponding trigger record.
317 const AliHLTMUONTriggerRecordStruct* trigrec = NULL;
318 for (size_t k = 0; k < triggerRecords.size(); k++)
319 {
320 if (triggerRecords[k]->fId == t.fTrigRec)
321 {
322 trigrec = triggerRecords[k];
323 break;
324 }
325 }
326 // If the trigger record was found then fill its hit information also.
327 if (trigrec != NULL)
328 {
329 AliHLTMUONParticleSign trsign;
330 bool trhitset[4];
331 AliHLTMUONUtils::UnpackTriggerRecordFlags(
332 trigrec->fFlags, trsign, trhitset
333 );
334
335 for (int i = 0; i < 4; i++)
336 {
337 if (not trhitset[i]) continue;
338
339 AliHLTUInt8_t chamber;
340 AliHLTUInt16_t detElemId;
341 AliHLTMUONUtils::UnpackRecHitFlags(trigrec->fHit[i].fFlags, chamber, detElemId);
342
343 AliESDMuonCluster cluster;
344 cluster.SetUniqueID(AliMUONVCluster::BuildUniqueID(chamber, detElemId, clusterIndex++));
345 cluster.SetXYZ(
346 trigrec->fHit[i].fX,
347 trigrec->fHit[i].fY,
348 trigrec->fHit[i].fZ
349 );
350 cluster.SetErrXY( // Use nominal values.
351 AliMUONConstants::TriggerNonBendingReso(),
352 AliMUONConstants::TriggerBendingReso()
353 );
354 cluster.SetCharge(-1.); // Indicate no total charge calculated.
355 cluster.SetChi2(-1.); // Indicate no fit made.
356 muTrack.AddCluster(cluster);
357 nHits++;
358 muTrack.AddInMuonClusterMap(i+10);
359 }
360 }
361 else
362 {
363 HLTWarning("Trigger record (ID = %d) not found for track ID = %d.",
364 t.fTrigRec, t.fId
365 );
366 }
367
368 muTrack.SetNHit(nHits);
369 event.AddMuonTrack(&muTrack);
370 }
371 }
372
373 PushBack(&event, AliHLTMUONConstants::ESDDataType(), specification);
374
375 return 0;
376}
377