4 //**************************************************************************
5 //* This file is property of and copyright by the ALICE HLT Project *
6 //* ALICE Experiment at CERN, All rights reserved. *
8 //* Primary Authors: Jochen Thaeder <thaeder@kip.uni-heidelberg.de> *
9 //* for The ALICE HLT Project. *
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 //**************************************************************************
20 /** @file AliHLTJETConeJetComponent.cxx
21 @author Jochen Thaeder <thaeder@kip.uni-heidelberg.de>
23 @brief Component to run the ConeJet jetfinder
34 #include "AliHLTJETConeJetComponent.h"
37 #include "TObjString.h"
39 /** ROOT macro for the implementation of ROOT specific class methods */
40 ClassImp(AliHLTJETConeJetComponent)
43 * ---------------------------------------------------------------------------------
44 * Constructor / Destructor
45 * ---------------------------------------------------------------------------------
48 // #################################################################################
49 AliHLTJETConeJetComponent::AliHLTJETConeJetComponent() :
54 fJetReaderHeader(NULL),
58 // see header file for class documentation
60 // refer to README to build package
62 // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt
66 // #################################################################################
67 AliHLTJETConeJetComponent::~AliHLTJETConeJetComponent() {
68 // see header file for class documentation
72 * ---------------------------------------------------------------------------------
73 * Public functions to implement AliHLTComponent's interface.
74 * These functions are required for the registration process
75 * ---------------------------------------------------------------------------------
78 // #################################################################################
79 const Char_t* AliHLTJETConeJetComponent::GetComponentID() {
80 // see header file for class documentation
81 return "JETConeJetFinder";
84 // #################################################################################
85 void AliHLTJETConeJetComponent::GetInputDataTypes( vector<AliHLTComponentDataType>& list) {
86 // see header file for class documentation
88 list.push_back( kAliHLTDataTypeMCObject|kAliHLTDataOriginHLT );
89 list.push_back( kAliHLTDataTypeESDObject|kAliHLTDataOriginOffline );
90 list.push_back( kAliHLTDataTypeESDObject|kAliHLTDataOriginHLT );
93 // #################################################################################
94 AliHLTComponentDataType AliHLTJETConeJetComponent::GetOutputDataType() {
95 // see header file for class documentation
96 return (kAliHLTDataTypeJet|kAliHLTDataOriginHLT);
99 // #################################################################################
100 void AliHLTJETConeJetComponent::GetOutputDataSize( ULong_t& constBase, Double_t& inputMultiplier ) {
101 // see header file for class documentation
104 inputMultiplier = 0.3;
107 // #################################################################################
108 AliHLTComponent* AliHLTJETConeJetComponent::Spawn() {
109 // see header file for class documentation
110 return new AliHLTJETConeJetComponent();
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 * ---------------------------------------------------------------------------------
121 // #################################################################################
122 Int_t AliHLTJETConeJetComponent::DoInit( Int_t argc, const Char_t** argv ) {
123 // see header file for class documentation
125 if ( fJetFinder || fJetHeader || fJetReader || fJetReaderHeader ||
126 fTrackCuts || fSeedCuts || fJetCuts || fJets )
129 // ---------------------------------------------------------------------
131 // ---------------------------------------------------------------------
133 TString comment = "HLT Fast Fixed Seeded Cone finder ";
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;
142 // ---------------------------------------------------------------------
144 // ---------------------------------------------------------------------
147 Int_t bMissingParam=0;
151 // -- Loop over all arguments
152 for ( Int_t iter = 0; iter<argc && iResult>=0; iter++) {
155 if (argument.IsNull())
159 if ( !argument.CompareTo("-algorithm") ) {
160 if ((bMissingParam=(++iter>=argc))) break;
162 TString parameter(argv[iter]);
163 parameter.Remove(TString::kLeading, ' ');
165 if ( !parameter.CompareTo("FSCSquareCell") ) {
166 algorithm = AliHLTJETBase::kFFSCSquareCell;
169 comment += parameter;
172 else if ( !parameter.CompareTo("FSCRadiusCell") ) {
173 algorithm = AliHLTJETBase::kFFSCRadiusCell;
176 comment += parameter;
180 HLTError("Wrong parameter %s for argument %s.", parameter.Data(), argument.Data());
186 else if ( !argument.CompareTo("-leading") ) {
187 if ((bMissingParam=(++iter>=argc))) break;
189 TString parameter(argv[iter]);
190 parameter.Remove(TString::kLeading, ' ');
192 if ( !parameter.CompareTo("0") ) {
196 comment += parameter;
199 else if ( !parameter.CompareTo("1") ) {
203 comment += parameter;
207 HLTError("Wrong parameter %s for argument %s.", parameter.Data(), argument.Data());
213 else if ( !argument.CompareTo("-coneRadius") ) {
214 if ((bMissingParam=(++iter>=argc))) break;
216 TString parameter(argv[iter]);
217 parameter.Remove(TString::kLeading, ' ');
219 if ( parameter.IsFloat() ) {
220 coneRadius = parameter.Atof();
223 comment += parameter;
227 HLTError("Wrong parameter %s for argument %s.", parameter.Data(), argument.Data());
233 else if ( !argument.CompareTo("-trackCutMinPt") ) {
234 if ((bMissingParam=(++iter>=argc))) break;
236 TString parameter(argv[iter]);
237 parameter.Remove(TString::kLeading, ' ');
239 if ( parameter.IsFloat() ) {
240 trackCutMinPt = parameter.Atof();
243 comment += parameter;
247 HLTError("Wrong parameter %s for argument %s.", parameter.Data(), argument.Data());
253 else if ( !argument.CompareTo("-seedCutMinPt") ) {
254 if ((bMissingParam=(++iter>=argc))) break;
256 TString parameter(argv[iter]);
257 parameter.Remove(TString::kLeading, ' ');
259 if ( parameter.IsFloat() ) {
260 seedCutMinPt = parameter.Atof();
263 comment += parameter;
267 HLTError("Wrong parameter %s for argument %s.", parameter.Data(), argument.Data());
273 else if ( !argument.CompareTo("-jetCutMinEt") ) {
274 if ((bMissingParam=(++iter>=argc))) break;
276 TString parameter(argv[iter]);
277 parameter.Remove(TString::kLeading, ' ');
279 if ( parameter.IsFloat() ) {
280 jetCutMinEt = parameter.Atof();
283 comment += parameter;
287 HLTError("Wrong parameter %s for argument %s.", parameter.Data(), argument.Data());
292 // -- Argument not known
294 HLTError("Unknown argument %s.", argument.Data());
297 } // for ( Int iter = 0; iter<argc && iResult>=0; iter++) {
299 // -- Check if parameter is missing
300 if ( bMissingParam ) {
301 HLTError("Missing parameter for argument %s.", argument.Data());
308 // ---------------------------------------------------------------------
310 // ---------------------------------------------------------------------
311 if ( ! (fTrackCuts = new AliHLTJETTrackCuts()) ) {
312 HLTError("Error instantiating track cuts");
316 fTrackCuts->SetChargedOnly( kTRUE );
317 fTrackCuts->SetMinPt( trackCutMinPt );
319 // ---------------------------------------------------------------------
321 // ---------------------------------------------------------------------
322 if ( ! (fSeedCuts = new AliHLTJETConeSeedCuts()) ) {
323 HLTError("Error instantiating seed cuts");
327 fSeedCuts->SetMinPt( seedCutMinPt );
329 // ---------------------------------------------------------------------
331 // ---------------------------------------------------------------------
332 if ( ! (fJetCuts = new AliHLTJETJetCuts()) ) {
333 HLTError("Error instantiating jet cuts");
337 fJetCuts->SetMinEt( jetCutMinEt );
339 // ---------------------------------------------------------------------
340 // -- Jet Reader Header
341 // ---------------------------------------------------------------------
342 if ( ! (fJetReaderHeader = new AliHLTJETReaderHeader()) ) {
343 HLTError("Error instantiating jet reader header");
348 fJetReaderHeader->SetJetAlgorithm( algorithm );
350 // Set prt to track cuts
351 fJetReaderHeader->SetTrackCuts( fTrackCuts );
352 fJetReaderHeader->SetSeedCuts( fSeedCuts );
354 // Set Eta min/max and Phi min/max
355 fJetReaderHeader->SetFiducialEta( -0.9, 0.9) ;
356 fJetReaderHeader->SetFiducialPhi( 0.0, TMath::TwoPi() ) ;
359 fJetReaderHeader->SetGridEtaBinning( 0.05 );
360 fJetReaderHeader->SetGridPhiBinning( 0.05 );
363 fJetReaderHeader->SetConeRadius(coneRadius);
365 // ---------------------------------------------------------------------
367 // ---------------------------------------------------------------------
368 if ( ! (fJetReader = new AliHLTJETReader()) ) {
369 HLTError("Error instantiating jet reader");
373 fJetReader->SetReaderHeader(fJetReaderHeader);
375 // ---------------------------------------------------------------------
377 // ---------------------------------------------------------------------
378 if ( ! (fJets = new AliHLTJets()) ) {
379 HLTError("Error instantiating jet container");
383 fJets->SetComment(comment);
385 // ---------------------------------------------------------------------
387 // ---------------------------------------------------------------------
388 if ( ! (fJetHeader = new AliHLTJETConeHeader()) ) {
389 HLTError("Error instantiating cone jet header");
393 fJetHeader->SetJetCuts(fJetCuts);
394 fJetHeader->SetUseLeading(leading);
396 // ---------------------------------------------------------------------
398 // ---------------------------------------------------------------------
399 if ( ! (fJetFinder = new AliHLTJETConeFinder()) ) {
400 HLTError("Error instantiating jet finder");
404 fJetFinder->SetJetHeader(fJetHeader);
405 fJetFinder->SetJetReader(fJetReader);
406 fJetFinder->SetOutputJets(fJets);
408 // ---------------------------------------------------------------------
409 // -- Initialize Jet Finder
410 // ---------------------------------------------------------------------
411 if ( (fJetFinder->Initialize()) ) {
412 HLTError("Error initializing cone jet finder");
419 // #################################################################################
420 Int_t AliHLTJETConeJetComponent::DoDeinit() {
421 // see header file for class documentation
435 if ( fJetReaderHeader )
436 delete fJetReaderHeader;
437 fJetReaderHeader = NULL;
458 // #################################################################################
459 Int_t AliHLTJETConeJetComponent::DoEvent( const AliHLTComponentEventData& /*evtData*/,
460 AliHLTComponentTriggerData& /*trigData*/ ) {
461 // see header file for class documentation
465 const TObject* iter = NULL;
469 if ( GetFirstInputObject(kAliHLTDataTypeSOR) && !iResult ) {
470 HLTInfo("On-line SOR Event");
473 // -- ADD MC Object -- On-line
474 // ------------------------------
475 for ( iter=GetFirstInputObject(kAliHLTDataTypeMCObject|kAliHLTDataOriginHLT);
476 iter != NULL && !iResult; iter=GetNextInputObject() ) {
478 // -- Set automatic MC usage, --> needed in off-line
479 fJetReaderHeader->SetUseMC(kTRUE);
481 // -- Set input event
482 fJetReader->SetInputEvent( NULL, NULL, const_cast<TObject*>(iter) );
484 // -- Fill grid with MC
485 if ( ! fJetReader->FillGridHLTMC() ) {
486 HLTError("Error filling grid.");
487 iResult = -EINPROGRESS;
492 if ( ! fJetFinder->ProcessHLTEvent() ) {
493 HLTError("Error processing cone event.");
494 iResult = -EINPROGRESS;
500 PushBack(fJets, kAliHLTDataTypeJet|kAliHLTDataOriginHLT, GetSpecification());
504 // -- ADD ESD Object -- Off-line
505 // -------------------------------
506 for ( iter=GetFirstInputObject(kAliHLTDataTypeESDObject|kAliHLTDataOriginOffline);
507 iter != NULL && !iResult; iter=GetNextInputObject() ) {
509 // -- Set automatic MC usage, --> needed in off-line
510 fJetReaderHeader->SetUseMC(kFALSE);
512 // -- Set input event
513 fJetReader->SetInputEvent( const_cast<TObject*>(iter), NULL, NULL );
515 // -- Fill grid with ESD
516 if ( ! fJetReader->FillGridESD() ) {
517 HLTError("Error filling grid.");
523 if ( ! fJetFinder->ProcessHLTEvent() ) {
524 HLTError("Error processing cone event.");
525 iResult = -EINPROGRESS;
531 PushBack(fJets, kAliHLTDataTypeJet|kAliHLTDataOriginHLT, GetSpecification());
535 // -- ADD ESD Object -- On-line
536 // ------------------------------
537 for ( iter=GetFirstInputObject(kAliHLTDataTypeESDObject|kAliHLTDataOriginHLT);
538 iter != NULL && !iResult; iter=GetNextInputObject() ) {
540 // -- Set automatic MC usage, --> needed in off-line
541 fJetReaderHeader->SetUseMC(kFALSE);
543 // -- Set input event
544 fJetReader->SetInputEvent( const_cast<TObject*>(iter), NULL, NULL );
546 // -- Fill grid with ESD
547 if ( ! fJetReader->FillGridESD() ) {
548 HLTError("Error filling grid.");
554 if ( ! fJetFinder->ProcessHLTEvent() ) {
555 HLTError("Error processing cone event.");
556 iResult = -EINPROGRESS;
562 PushBack(fJets, kAliHLTDataTypeJet|kAliHLTDataOriginHLT, GetSpecification());