]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/JET/fastjet/AliHLTJETFastJetComponent.cxx
Adding the correction for wrongly simulated pid TOF signal
[u/mrichter/AliRoot.git] / HLT / JET / fastjet / AliHLTJETFastJetComponent.cxx
1 //-*- Mode: C++ -*-
2 // $Id: $
3
4 //**************************************************************************
5 //* This file is property of and copyright by the ALICE HLT Project        * 
6 //* ALICE Experiment at CERN, All rights reserved.                         *
7 //*                                                                        *
8 //* Primary Authors: Jochen Thaeder <jochen@thaeder.de>                    *
9 //*                  for The ALICE HLT Project.                            *
10 //*                                                                        *
11 //* Permission to use, copy, modify and distribute this software and its   *
12 //* documentation strictly for non-commercial purposes is hereby granted   *
13 //* without fee, provided that the above copyright notice appears in all   *
14 //* copies and that both the copyright notice and this permission notice   *
15 //* appear in the supporting documentation. The authors make no claims     *
16 //* about the suitability of this software for any purpose. It is          *
17 //* provided "as is" without express or implied warranty.                  *
18 //**************************************************************************
19
20 /** @file   AliHLTJETFastJetComponent.cxx
21     @author Jochen Thaeder <jochen@thaeder.de>
22     @brief  Component to run the FastJet jetfinder
23 */
24
25 #include <cstdlib>
26 #include <cerrno>
27 #include <sys/time.h>
28
29 #include "AliHLTJETFastJetComponent.h" 
30
31 #include "TString.h"
32 #include "TObjString.h"
33
34 using namespace std;
35
36 /** ROOT macro for the implementation of ROOT specific class methods */
37 ClassImp(AliHLTJETFastJetComponent)
38
39 /*
40  * ---------------------------------------------------------------------------------
41  *                            Constructor / Destructor
42  * ---------------------------------------------------------------------------------
43  */
44
45 // #################################################################################
46 AliHLTJETFastJetComponent::AliHLTJETFastJetComponent() 
47   :
48   fJetFinder(NULL),
49   fJetHeader(NULL),
50   fJetReader(NULL),
51   fJetReaderHeader(NULL),
52   fTrackCuts(NULL),
53   fJetCuts(NULL),
54   fJets(NULL) {
55   // see header file for class documentation
56   // or
57   // refer to README to build package
58   // or
59   // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt
60 }
61
62 // #################################################################################
63 AliHLTJETFastJetComponent::~AliHLTJETFastJetComponent() {
64   // see header file for class documentation
65
66 }
67
68 /*
69  * ---------------------------------------------------------------------------------
70  * Public functions to implement AliHLTComponent's interface.
71  * These functions are required for the registration process
72  * ---------------------------------------------------------------------------------
73  */
74
75 // #################################################################################
76 const Char_t* AliHLTJETFastJetComponent::GetComponentID() {
77   // see header file for class documentation
78   return "JETFastJetFinder";
79 }
80
81 // #################################################################################
82 void AliHLTJETFastJetComponent::GetInputDataTypes( vector<AliHLTComponentDataType>& list) {
83   // see header file for class documentation
84   list.clear(); 
85   list.push_back( kAliHLTDataTypeMCObject|kAliHLTDataOriginHLT );
86   list.push_back( kAliHLTDataTypeESDObject|kAliHLTDataOriginOffline );
87   list.push_back( kAliHLTDataTypeESDObject|kAliHLTDataOriginHLT );
88 }
89
90 // #################################################################################
91 AliHLTComponentDataType AliHLTJETFastJetComponent::GetOutputDataType() {
92   // see header file for class documentation
93   return (kAliHLTDataTypeJet|kAliHLTDataOriginHLT);
94 }
95
96 // #################################################################################
97 void AliHLTJETFastJetComponent::GetOutputDataSize( ULong_t& constBase, Double_t& inputMultiplier ) {
98   // see header file for class documentation
99
100   constBase = 1000;
101   inputMultiplier = 0.3;
102 }
103
104 // #################################################################################
105 AliHLTComponent* AliHLTJETFastJetComponent::Spawn() {
106   // see header file for class documentation
107   return new AliHLTJETFastJetComponent();
108 }
109
110 /*
111  * ---------------------------------------------------------------------------------
112  * Protected functions to implement AliHLTComponent's interface.
113  * These functions provide initialization as well as the actual processing
114  * capabilities of the component. 
115  * ---------------------------------------------------------------------------------
116  */
117
118 // #################################################################################
119 Int_t AliHLTJETFastJetComponent::DoInit( Int_t argc, const Char_t** argv ) {
120   // see header file for class documentation
121
122   if ( fJetFinder || fJetHeader || fJetReader || fJetReaderHeader || 
123        fTrackCuts || fJetCuts || fJets )
124     return -EINPROGRESS;
125
126   // ---------------------------------------------------------------------
127   // -- Defaults
128   // ---------------------------------------------------------------------
129
130   TString comment       = "HLT FastJet interface";
131
132   AliHLTJETBase::JetAlgorithmType_t algorithm = AliHLTJETBase::kKt; 
133   Float_t coneRadius    =  0.4;
134   Float_t trackCutMinPt =  1.0;
135   Float_t jetCutMinEt   = 15.0;
136
137   // ---------------------------------------------------------------------
138   // -- Get Arguments
139   // ---------------------------------------------------------------------
140
141   Int_t iResult = 0;
142   Int_t bMissingParam=0;
143   
144   TString argument="";
145
146   // -- Loop over all arguments
147   for ( Int_t iter = 0; iter<argc && iResult>=0; iter++) {
148     argument=argv[iter];
149
150     if (argument.IsNull()) 
151       continue;
152
153     // -- algorithm
154     if ( !argument.CompareTo("-algorithm") ) {
155       if ((bMissingParam=(++iter>=argc))) break;
156       
157       TString parameter(argv[iter]);
158       parameter.Remove(TString::kLeading, ' ');
159       
160       if ( !parameter.CompareTo("Kt") ) {
161         algorithm = AliHLTJETBase::kKt;
162         comment += argument;
163         comment += " ";
164         comment += parameter;
165         comment += ' ';
166       }
167       else if ( !parameter.CompareTo("AntiKt") ) {
168         algorithm = AliHLTJETBase::kAntiKt;
169         comment += argument;
170         comment += " ";
171         comment += parameter;
172         comment += ' ';
173       }
174       else {
175         HLTError("Wrong parameter %s for argument %s.", parameter.Data(), argument.Data());
176         iResult=-EINVAL;
177       }
178     } 
179
180     // -- coneRadius
181     else if ( !argument.CompareTo("-coneRadius") ) {
182       if ((bMissingParam=(++iter>=argc))) break;
183       
184       TString parameter(argv[iter]);
185       parameter.Remove(TString::kLeading, ' ');
186       
187       if ( parameter.IsFloat() ) {
188         coneRadius = parameter.Atof();
189         comment += argument;
190         comment += " ";
191         comment += parameter;
192         comment += ' ';
193       }
194       else {
195         HLTError("Wrong parameter %s for argument %s.", parameter.Data(), argument.Data());
196         iResult=-EINVAL;
197       }
198     } 
199
200     // -- trackCutMinPt
201     else if ( !argument.CompareTo("-trackCutMinPt") ) {
202       if ((bMissingParam=(++iter>=argc))) break;
203       
204       TString parameter(argv[iter]);
205       parameter.Remove(TString::kLeading, ' ');
206       
207       if ( parameter.IsFloat() ) {
208         trackCutMinPt = parameter.Atof();
209         comment += argument;
210         comment += " ";
211         comment += parameter;
212         comment += ' ';
213       }
214       else {
215         HLTError("Wrong parameter %s for argument %s.", parameter.Data(), argument.Data());
216         iResult=-EINVAL;
217       }
218     } 
219
220     // -- jetCutMinEt
221     else if ( !argument.CompareTo("-jetCutMinEt") ) {
222       if ((bMissingParam=(++iter>=argc))) break;
223       
224       TString parameter(argv[iter]);
225       parameter.Remove(TString::kLeading, ' ');
226       
227       if ( parameter.IsFloat() ) {
228         jetCutMinEt = parameter.Atof();
229         comment += argument;
230         comment += " ";
231         comment += parameter;
232         comment += ' ';
233       }
234       else {
235         HLTError("Wrong parameter %s for argument %s.", parameter.Data(), argument.Data());
236         iResult=-EINVAL;
237       }
238     } 
239     
240     // -- Argument not known
241     else {
242       HLTError("Unknown argument %s.", argument.Data());
243       iResult = -EINVAL;
244     }
245   } // for ( Int iter = 0; iter<argc && iResult>=0; iter++) {
246   
247   // -- Check if parameter is missing
248   if ( bMissingParam ) {
249     HLTError("Missing parameter for argument %s.", argument.Data());
250     iResult=-EINVAL;
251   }
252
253   if (iResult)
254     return iResult;
255
256   // ---------------------------------------------------------------------
257   // -- Jet Track Cuts
258   // ---------------------------------------------------------------------
259   if ( ! (fTrackCuts = new AliHLTJETTrackCuts()) ) {
260     HLTError("Error instantiating track cuts");
261     return -EINPROGRESS;
262   }
263
264   fTrackCuts->SetChargedOnly( kTRUE );
265   fTrackCuts->SetMinPt( trackCutMinPt );
266
267   // ---------------------------------------------------------------------
268   // -- Jet Jet Cuts
269   // ---------------------------------------------------------------------
270   if ( ! (fJetCuts = new AliHLTJETJetCuts()) ) {
271     HLTError("Error instantiating jet cuts");
272     return -EINPROGRESS;
273   }
274
275   fJetCuts->SetMinEt( jetCutMinEt );
276
277   // ---------------------------------------------------------------------
278   // -- Jet Reader Header
279   // ---------------------------------------------------------------------
280   if ( ! (fJetReaderHeader = new AliHLTJETReaderHeader()) ) {
281     HLTError("Error instantiating jet reader header");
282     return -EINPROGRESS;
283   }
284   
285   // Set Algorithm
286   fJetReaderHeader->SetJetAlgorithm( algorithm );
287
288   // Set prt to track cuts
289   fJetReaderHeader->SetTrackCuts( fTrackCuts );
290
291   // Set Eta min/max and Phi min/max
292   fJetReaderHeader->SetFiducialEta( -0.9, 0.9) ;
293   fJetReaderHeader->SetFiducialPhi(  0.0, TMath::TwoPi() ) ;
294  
295   // Set cone radius
296   fJetReaderHeader->SetConeRadius(coneRadius);
297
298   // ---------------------------------------------------------------------
299   // -- Jet Reader
300   // ---------------------------------------------------------------------
301   if ( ! (fJetReader = new AliHLTJETReader()) ) {
302     HLTError("Error instantiating jet reader");
303     return -EINPROGRESS;
304   }
305
306   fJetReader->SetReaderHeader(fJetReaderHeader);
307
308   // ---------------------------------------------------------------------
309   // -- Jet Container
310   // ---------------------------------------------------------------------
311   if ( ! (fJets = new AliHLTJets()) ) {
312     HLTError("Error instantiating jet container");
313     return -EINPROGRESS;
314   }
315
316   fJets->SetComment(comment);
317
318   // ---------------------------------------------------------------------
319   // -- Jet Header
320   // ---------------------------------------------------------------------
321   if ( ! (fJetHeader = new AliHLTJETFastJetHeader()) ) {
322     HLTError("Error instantiating fastjet header");
323     return -EINPROGRESS;
324   }
325
326   fJetHeader->SetReaderHeader(fJetReaderHeader);
327   fJetHeader->SetJetCuts(fJetCuts);
328
329   // ---------------------------------------------------------------------
330   // -- Jet Finder
331   // ---------------------------------------------------------------------
332   if ( ! (fJetFinder = new AliHLTJETFastJetFinder()) ) {
333     HLTError("Error instantiating fastjet finder");
334     return -EINPROGRESS;
335   }
336
337   fJetFinder->SetJetHeader(fJetHeader);
338   fJetFinder->SetJetReader(fJetReader);
339   fJetFinder->SetOutputJets(fJets);
340
341   // ---------------------------------------------------------------------
342   // -- Initialize Jet Finder
343   // ---------------------------------------------------------------------
344   if ( (fJetFinder->Initialize()) ) {
345     HLTError("Error initializing fastjet finder");
346     return -EINPROGRESS;
347   }
348
349   return 0;
350 }
351
352 // #################################################################################
353 Int_t AliHLTJETFastJetComponent::DoDeinit() {
354   // see header file for class documentation
355
356  if ( fJetFinder )
357     delete fJetFinder;
358   fJetFinder = NULL;
359
360   if ( fJetHeader )
361     delete fJetHeader;
362   fJetHeader = NULL;
363
364   if ( fJetReader )
365     delete fJetReader;
366   fJetReader = NULL;
367
368   if ( fJetReaderHeader )
369     delete fJetReaderHeader;
370   fJetReaderHeader = NULL;
371
372   if ( fJetCuts )
373     delete fJetCuts;
374   fJetCuts = NULL;
375
376   if ( fJets )
377     delete fJets;
378   fJets = NULL;
379
380   return 0;
381 }
382
383 // #################################################################################
384 Int_t AliHLTJETFastJetComponent::DoEvent( const AliHLTComponentEventData& /*evtData*/,
385                                           AliHLTComponentTriggerData& /*trigData*/ ) {
386   // see header file for class documentation
387
388   Int_t iResult = 0;
389
390   const TObject* iter = NULL;
391
392   // -- Start-Of-Run
393   // -----------------
394   if ( GetFirstInputObject(kAliHLTDataTypeSOR) && !iResult ) {
395     HLTInfo("On-line SOR Event");
396   }
397
398   // -- ADD MC Object -- On-line
399   // ------------------------------
400   for ( iter=GetFirstInputObject(kAliHLTDataTypeMCObject|kAliHLTDataOriginHLT); 
401         iter != NULL && !iResult; iter=GetNextInputObject() ) {
402
403     // -- Set automatic MC usage, --> needed in off-line
404     fJetReaderHeader->SetUseMC(kTRUE);
405
406     // -- Set input event
407     fJetReader->SetInputEvent( NULL, NULL, const_cast<TObject*>(iter) );    
408
409     // -- Fill vector with MC
410     if ( ! fJetReader->FillVectorHLTMC() ) {
411       HLTError("Error filling vector.");
412       iResult = -EINPROGRESS;
413     }
414
415     // -- Find jets
416     if ( !iResult) {
417       if ( ! fJetFinder->ProcessHLTEvent() ) {
418         HLTError("Error processing fastjet event.");
419         iResult = -EINPROGRESS;
420       }
421     }
422     
423     // -- PushBack
424     if ( !iResult) {
425       PushBack(fJets, kAliHLTDataTypeJet|kAliHLTDataOriginHLT, GetSpecification());
426     }
427   }
428
429   // -- ADD ESD Object -- Off-line
430   // -------------------------------
431   for ( iter=GetFirstInputObject(kAliHLTDataTypeESDObject|kAliHLTDataOriginOffline); 
432         iter != NULL && !iResult; iter=GetNextInputObject() ) {
433
434     // -- Set automatic MC usage, --> needed in off-line
435     fJetReaderHeader->SetUseMC(kFALSE);
436
437     // -- Set input event
438     fJetReader->SetInputEvent( const_cast<TObject*>(iter), NULL, NULL );    
439   
440     // -- Fill vector with ESD
441     if ( ! fJetReader->FillVectorESD() ) {
442       HLTError("Error filling vector.");
443       iResult = -1;  
444     }
445
446     // -- Find jets
447     if ( !iResult) {
448       if ( ! fJetFinder->ProcessHLTEvent() ) {
449         HLTError("Error processing fastjet event.");
450         iResult = -EINPROGRESS;
451       }
452     }
453
454     // -- PushBack
455     if ( !iResult) {
456       PushBack(fJets, kAliHLTDataTypeJet|kAliHLTDataOriginHLT, GetSpecification());
457     }
458   }
459
460   // -- ADD ESD Object -- On-line
461   // ------------------------------
462   for ( iter=GetFirstInputObject(kAliHLTDataTypeESDObject|kAliHLTDataOriginHLT); 
463         iter != NULL && !iResult; iter=GetNextInputObject() ) {
464
465     // -- Set automatic MC usage, --> needed in off-line
466     fJetReaderHeader->SetUseMC(kFALSE);
467
468     // -- Set input event
469     fJetReader->SetInputEvent( const_cast<TObject*>(iter), NULL, NULL );    
470
471     // -- Fill vector with ESD
472     if ( ! fJetReader->FillVectorESD() ) {
473       HLTError("Error filling vector.");
474       iResult = -1;  
475     }
476
477     // -- Find jets
478     if ( !iResult) {
479       if ( ! fJetFinder->ProcessHLTEvent() ) {
480         HLTError("Error processing fastjet event.");
481         iResult = -EINPROGRESS;
482       }
483     }
484
485     // -- PushBack
486     if ( !iResult) {
487       PushBack(fJets, kAliHLTDataTypeJet|kAliHLTDataOriginHLT, GetSpecification());
488     }
489   }
490
491   return iResult;
492 }