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