]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/JET/cone/AliHLTJETConeJetComponent.cxx
* Added fast interface to fastjet
[u/mrichter/AliRoot.git] / HLT / JET / cone / AliHLTJETConeJetComponent.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 <thaeder@kip.uni-heidelberg.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   AliHLTJETConeJetComponent.cxx
21     @author Jochen Thaeder <thaeder@kip.uni-heidelberg.de>
22     @date   
23     @brief  Component to run the ConeJet jetfinder
24 */
25
26 #if __GNUC__>= 3
27 using namespace std;
28 #endif
29
30 #include <cstdlib>
31 #include <cerrno>
32 #include <sys/time.h>
33
34 #include "AliHLTJETConeJetComponent.h" 
35
36 #include "TString.h"
37 #include "TObjString.h"
38
39 /** ROOT macro for the implementation of ROOT specific class methods */
40 ClassImp(AliHLTJETConeJetComponent)
41
42 /*
43  * ---------------------------------------------------------------------------------
44  *                            Constructor / Destructor
45  * ---------------------------------------------------------------------------------
46  */
47
48 // #################################################################################
49 AliHLTJETConeJetComponent::AliHLTJETConeJetComponent() :
50   fJetFinder(NULL),
51   fJetHeader(NULL),
52   fSeedCuts(NULL),
53   fJetReader(NULL),
54   fJetReaderHeader(NULL),
55   fTrackCuts(NULL),
56   fJetCuts(NULL),
57   fJets(NULL) {
58   // see header file for class documentation
59   // or
60   // refer to README to build package
61   // or
62   // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt
63
64 }
65
66 // #################################################################################
67 AliHLTJETConeJetComponent::~AliHLTJETConeJetComponent() {
68   // see header file for class documentation
69 }
70
71 /*
72  * ---------------------------------------------------------------------------------
73  * Public functions to implement AliHLTComponent's interface.
74  * These functions are required for the registration process
75  * ---------------------------------------------------------------------------------
76  */
77
78 // #################################################################################
79 const Char_t* AliHLTJETConeJetComponent::GetComponentID() {
80   // see header file for class documentation
81   return "JETConeJetFinder";
82 }
83
84 // #################################################################################
85 void AliHLTJETConeJetComponent::GetInputDataTypes( vector<AliHLTComponentDataType>& list) {
86   // see header file for class documentation
87   list.clear(); 
88   list.push_back( kAliHLTDataTypeMCObject|kAliHLTDataOriginHLT );
89   list.push_back( kAliHLTDataTypeESDObject|kAliHLTDataOriginOffline );
90   list.push_back( kAliHLTDataTypeESDObject|kAliHLTDataOriginHLT );
91 }
92
93 // #################################################################################
94 AliHLTComponentDataType AliHLTJETConeJetComponent::GetOutputDataType() {
95   // see header file for class documentation
96   return (kAliHLTDataTypeJet|kAliHLTDataOriginHLT);
97 }
98
99 // #################################################################################
100 void AliHLTJETConeJetComponent::GetOutputDataSize( ULong_t& constBase, Double_t& inputMultiplier ) {
101   // see header file for class documentation
102
103   constBase = 1000;
104   inputMultiplier = 0.3;
105 }
106
107 // #################################################################################
108 AliHLTComponent* AliHLTJETConeJetComponent::Spawn() {
109   // see header file for class documentation
110   return new AliHLTJETConeJetComponent();
111 }
112
113 /*
114  * ---------------------------------------------------------------------------------
115  * Protected functions to implement AliHLTComponent's interface.
116  * These functions provide initialization as well as the actual processing
117  * capabilities of the component. 
118  * ---------------------------------------------------------------------------------
119  */
120
121 // #################################################################################
122 Int_t AliHLTJETConeJetComponent::DoInit( Int_t argc, const Char_t** argv ) {
123   // see header file for class documentation
124
125   if ( fJetFinder || fJetHeader || fJetReader || fJetReaderHeader || 
126        fTrackCuts || fSeedCuts || fJetCuts || fJets )
127     return -EINPROGRESS;
128
129   // ---------------------------------------------------------------------
130   // -- Defaults
131   // ---------------------------------------------------------------------
132
133   TString comment       = "HLT Fast Fixed Seeded Cone finder ";
134
135   AliHLTJETBase::JetAlgorithmType_t algorithm = AliHLTJETBase::kFFSCSquareCell; 
136   Bool_t  leading       = kFALSE;
137   Float_t coneRadius    =  0.4;
138   Float_t trackCutMinPt =  1.0;
139   Float_t seedCutMinPt  =  5.0;
140   Float_t jetCutMinEt   = 15.0;
141
142   // ---------------------------------------------------------------------
143   // -- Get Arguments
144   // ---------------------------------------------------------------------
145
146   Int_t iResult = 0;
147   Int_t bMissingParam=0;
148   
149   TString argument="";
150
151   // -- Loop over all arguments
152   for ( Int_t iter = 0; iter<argc && iResult>=0; iter++) {
153     argument=argv[iter];
154
155     if (argument.IsNull()) 
156       continue;
157
158     // -- algorithm
159     if ( !argument.CompareTo("-algorithm") ) {
160       if ((bMissingParam=(++iter>=argc))) break;
161       
162       TString parameter(argv[iter]);
163       parameter.Remove(TString::kLeading, ' ');
164       
165       if ( !parameter.CompareTo("FSCSquareCell") ) {
166         algorithm = AliHLTJETBase::kFFSCSquareCell;
167         comment += argument;
168         comment += " ";
169         comment += parameter;
170         comment += ' ';
171       }
172       else if ( !parameter.CompareTo("FSCRadiusCell") ) {
173         algorithm = AliHLTJETBase::kFFSCRadiusCell;
174         comment += argument;
175         comment += " ";
176         comment += parameter;
177         comment += ' ';
178       }
179       else {
180         HLTError("Wrong parameter %s for argument %s.", parameter.Data(), argument.Data());
181         iResult=-EINVAL;
182       }
183     } 
184
185     // -- leading
186     else if ( !argument.CompareTo("-leading") ) {
187       if ((bMissingParam=(++iter>=argc))) break;
188       
189       TString parameter(argv[iter]);
190       parameter.Remove(TString::kLeading, ' ');
191       
192       if ( !parameter.CompareTo("0") ) {
193         leading = kFALSE;
194         comment += argument;
195         comment += " ";
196         comment += parameter;
197         comment += ' ';
198       }
199       else if ( !parameter.CompareTo("1") ) {
200         leading = kTRUE;
201         comment += argument;
202         comment += " ";
203         comment += parameter;
204         comment += ' ';
205       }
206       else {
207         HLTError("Wrong parameter %s for argument %s.", parameter.Data(), argument.Data());
208         iResult=-EINVAL;
209       }
210     } 
211
212     // -- coneRadius
213     else if ( !argument.CompareTo("-coneRadius") ) {
214       if ((bMissingParam=(++iter>=argc))) break;
215       
216       TString parameter(argv[iter]);
217       parameter.Remove(TString::kLeading, ' ');
218       
219       if ( parameter.IsFloat() ) {
220         coneRadius = parameter.Atof();
221         comment += argument;
222         comment += " ";
223         comment += parameter;
224         comment += ' ';
225       }
226       else {
227         HLTError("Wrong parameter %s for argument %s.", parameter.Data(), argument.Data());
228         iResult=-EINVAL;
229       }
230     } 
231
232     // -- trackCutMinPt
233     else if ( !argument.CompareTo("-trackCutMinPt") ) {
234       if ((bMissingParam=(++iter>=argc))) break;
235       
236       TString parameter(argv[iter]);
237       parameter.Remove(TString::kLeading, ' ');
238       
239       if ( parameter.IsFloat() ) {
240         trackCutMinPt = parameter.Atof();
241         comment += argument;
242         comment += " ";
243         comment += parameter;
244         comment += ' ';
245       }
246       else {
247         HLTError("Wrong parameter %s for argument %s.", parameter.Data(), argument.Data());
248         iResult=-EINVAL;
249       }
250     } 
251
252     // -- seedCutMinPt
253     else if ( !argument.CompareTo("-seedCutMinPt") ) {
254       if ((bMissingParam=(++iter>=argc))) break;
255       
256       TString parameter(argv[iter]);
257       parameter.Remove(TString::kLeading, ' ');
258       
259       if ( parameter.IsFloat() ) {
260         seedCutMinPt = parameter.Atof();
261         comment += argument;
262         comment += " ";
263         comment += parameter;
264         comment += ' ';
265       }
266       else {
267         HLTError("Wrong parameter %s for argument %s.", parameter.Data(), argument.Data());
268         iResult=-EINVAL;
269       }
270     } 
271
272     // -- jetCutMinEt
273     else if ( !argument.CompareTo("-jetCutMinEt") ) {
274       if ((bMissingParam=(++iter>=argc))) break;
275       
276       TString parameter(argv[iter]);
277       parameter.Remove(TString::kLeading, ' ');
278       
279       if ( parameter.IsFloat() ) {
280         jetCutMinEt = parameter.Atof();
281         comment += argument;
282         comment += " ";
283         comment += parameter;
284         comment += ' ';
285       }
286       else {
287         HLTError("Wrong parameter %s for argument %s.", parameter.Data(), argument.Data());
288         iResult=-EINVAL;
289       }
290     } 
291     
292     // -- Argument not known
293     else {
294       HLTError("Unknown argument %s.", argument.Data());
295       iResult = -EINVAL;
296     }
297   } // for ( Int iter = 0; iter<argc && iResult>=0; iter++) {
298   
299   // -- Check if parameter is missing
300   if ( bMissingParam ) {
301     HLTError("Missing parameter for argument %s.", argument.Data());
302     iResult=-EINVAL;
303   }
304
305   if (iResult)
306     return iResult;
307
308   // ---------------------------------------------------------------------
309   // -- Jet Track Cuts
310   // ---------------------------------------------------------------------
311   if ( ! (fTrackCuts = new AliHLTJETTrackCuts()) ) {
312     HLTError("Error instantiating track cuts");
313     return -EINPROGRESS;
314   }
315
316   fTrackCuts->SetChargedOnly( kTRUE );
317   fTrackCuts->SetMinPt( trackCutMinPt );
318
319   // ---------------------------------------------------------------------
320   // -- Jet Seed Cuts
321   // ---------------------------------------------------------------------
322   if ( ! (fSeedCuts = new AliHLTJETConeSeedCuts()) ) {
323     HLTError("Error instantiating seed cuts");
324     return -EINPROGRESS;
325   }
326
327   fSeedCuts->SetMinPt( seedCutMinPt );
328
329   // ---------------------------------------------------------------------
330   // -- Jet Jet Cuts
331   // ---------------------------------------------------------------------
332   if ( ! (fJetCuts = new AliHLTJETJetCuts()) ) {
333     HLTError("Error instantiating jet cuts");
334     return -EINPROGRESS;
335   }
336
337   fJetCuts->SetMinEt( jetCutMinEt );
338
339   // ---------------------------------------------------------------------
340   // -- Jet Reader Header
341   // ---------------------------------------------------------------------
342   if ( ! (fJetReaderHeader = new AliHLTJETReaderHeader()) ) {
343     HLTError("Error instantiating jet reader header");
344     return -EINPROGRESS;
345   }
346   
347   // Set Algorithm
348   fJetReaderHeader->SetJetAlgorithm( algorithm );
349
350   // Set prt to track cuts
351   fJetReaderHeader->SetTrackCuts( fTrackCuts );
352   fJetReaderHeader->SetSeedCuts( fSeedCuts );
353
354   // Set Eta min/max and Phi min/max
355   fJetReaderHeader->SetFiducialEta( -0.9, 0.9) ;
356   fJetReaderHeader->SetFiducialPhi(  0.0, TMath::TwoPi() ) ;
357
358   // Set grid binning
359   fJetReaderHeader->SetGridEtaBinning( 0.05 );
360   fJetReaderHeader->SetGridPhiBinning( 0.05 );
361  
362   // Set cone radius
363   fJetReaderHeader->SetConeRadius(coneRadius);
364
365   // ---------------------------------------------------------------------
366   // -- Jet Reader
367   // ---------------------------------------------------------------------
368   if ( ! (fJetReader = new AliHLTJETReader()) ) {
369     HLTError("Error instantiating jet reader");
370     return -EINPROGRESS;
371   }
372
373   fJetReader->SetReaderHeader(fJetReaderHeader);
374
375   // ---------------------------------------------------------------------
376   // -- Jet Container
377   // ---------------------------------------------------------------------
378   if ( ! (fJets = new AliHLTJets()) ) {
379     HLTError("Error instantiating jet container");
380     return -EINPROGRESS;
381   }
382
383   fJets->SetComment(comment);
384
385   // ---------------------------------------------------------------------
386   // -- Jet Header
387   // ---------------------------------------------------------------------
388   if ( ! (fJetHeader = new AliHLTJETConeHeader()) ) {
389     HLTError("Error instantiating cone jet header");
390     return -EINPROGRESS;
391   }
392
393   fJetHeader->SetJetCuts(fJetCuts);
394   fJetHeader->SetUseLeading(leading);
395
396   // ---------------------------------------------------------------------
397   // -- Jet Finder
398   // ---------------------------------------------------------------------
399   if ( ! (fJetFinder = new AliHLTJETConeFinder()) ) {
400     HLTError("Error instantiating jet finder");
401     return -EINPROGRESS;
402   }
403
404   fJetFinder->SetJetHeader(fJetHeader);
405   fJetFinder->SetJetReader(fJetReader);
406   fJetFinder->SetOutputJets(fJets);
407
408   // ---------------------------------------------------------------------
409   // -- Initialize Jet Finder
410   // ---------------------------------------------------------------------
411   if ( (fJetFinder->Initialize()) ) {
412     HLTError("Error initializing cone jet finder");
413     return -EINPROGRESS;
414   }
415
416   return 0;
417 }
418
419 // #################################################################################
420 Int_t AliHLTJETConeJetComponent::DoDeinit() {
421   // see header file for class documentation
422
423   if ( fJetFinder )
424     delete fJetFinder;
425   fJetFinder = NULL;
426
427   if ( fJetHeader )
428     delete fJetHeader;
429   fJetHeader = NULL;
430
431   if ( fJetReader )
432     delete fJetReader;
433   fJetReader = NULL;
434
435   if ( fJetReaderHeader )
436     delete fJetReaderHeader;
437   fJetReaderHeader = NULL;
438
439   if ( fTrackCuts )
440     delete fTrackCuts;
441   fTrackCuts = NULL;
442
443   if ( fSeedCuts )
444     delete fSeedCuts;
445   fSeedCuts = NULL;
446
447   if ( fJetCuts )
448     delete fJetCuts;
449   fJetCuts = NULL;
450
451   if ( fJets )
452     delete fJets;
453   fJets = NULL;
454   
455   return 0;
456 }
457
458 // #################################################################################
459 Int_t AliHLTJETConeJetComponent::DoEvent( const AliHLTComponentEventData& /*evtData*/,
460                                           AliHLTComponentTriggerData& /*trigData*/ ) {
461   // see header file for class documentation
462
463   Int_t iResult = 0;
464
465   const TObject* iter = NULL;
466
467   // -- Start-Of-Run
468   // -----------------
469   if ( GetFirstInputObject(kAliHLTDataTypeSOR) && !iResult ) {
470     HLTInfo("On-line SOR Event");
471   }
472
473   // -- ADD MC Object -- On-line
474   // ------------------------------
475   for ( iter=GetFirstInputObject(kAliHLTDataTypeMCObject|kAliHLTDataOriginHLT); 
476         iter != NULL && !iResult; iter=GetNextInputObject() ) {
477
478     // -- Set automatic MC usage, --> needed in off-line
479     fJetReaderHeader->SetUseMC(kTRUE);
480
481     // -- Set input event
482     fJetReader->SetInputEvent( NULL, NULL, const_cast<TObject*>(iter) );    
483
484     // -- Fill grid with MC
485     if ( ! fJetReader->FillGridHLTMC() ) {
486       HLTError("Error filling grid.");
487       iResult = -EINPROGRESS;
488     }
489
490     // -- Find jets
491     if ( !iResult) {
492       if ( ! fJetFinder->ProcessHLTEvent() ) {
493         HLTError("Error processing cone event.");
494         iResult = -EINPROGRESS;
495       }
496     }
497     
498     // -- PushBack
499     if ( !iResult) {
500       PushBack(fJets, kAliHLTDataTypeJet|kAliHLTDataOriginHLT, GetSpecification());
501     }
502   }
503
504   // -- ADD ESD Object -- Off-line
505   // -------------------------------
506   for ( iter=GetFirstInputObject(kAliHLTDataTypeESDObject|kAliHLTDataOriginOffline); 
507         iter != NULL && !iResult; iter=GetNextInputObject() ) {
508
509     // -- Set automatic MC usage, --> needed in off-line
510     fJetReaderHeader->SetUseMC(kFALSE);
511
512     // -- Set input event
513     fJetReader->SetInputEvent( const_cast<TObject*>(iter), NULL, NULL );    
514   
515     // -- Fill grid with ESD
516     if ( ! fJetReader->FillGridESD() ) {
517       HLTError("Error filling grid.");
518       iResult = -1;  
519     }
520
521     // -- Find jets
522     if ( !iResult) {
523       if ( ! fJetFinder->ProcessHLTEvent() ) {
524         HLTError("Error processing cone event.");
525         iResult = -EINPROGRESS;
526       }
527     }
528
529     // -- PushBack
530     if ( !iResult) {
531       PushBack(fJets, kAliHLTDataTypeJet|kAliHLTDataOriginHLT, GetSpecification());
532     }
533   }
534
535   // -- ADD ESD Object -- On-line
536   // ------------------------------
537   for ( iter=GetFirstInputObject(kAliHLTDataTypeESDObject|kAliHLTDataOriginHLT); 
538         iter != NULL && !iResult; iter=GetNextInputObject() ) {
539
540     // -- Set automatic MC usage, --> needed in off-line
541     fJetReaderHeader->SetUseMC(kFALSE);
542
543     // -- Set input event
544     fJetReader->SetInputEvent( const_cast<TObject*>(iter), NULL, NULL );    
545
546     // -- Fill grid with ESD
547     if ( ! fJetReader->FillGridESD() ) {
548       HLTError("Error filling grid.");
549       iResult = -1;  
550     }
551
552     // -- Find jets
553     if ( !iResult) {
554       if ( ! fJetFinder->ProcessHLTEvent() ) {
555         HLTError("Error processing cone event.");
556         iResult = -EINPROGRESS;
557       }
558     }
559
560     // -- PushBack
561     if ( !iResult) {
562       PushBack(fJets, kAliHLTDataTypeJet|kAliHLTDataOriginHLT, GetSpecification());
563     }
564   }
565
566   return iResult;
567 }