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 <jochen@thaeder.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 AliHLTJETFastJetComponent.cxx
21 @author Jochen Thaeder <jochen@thaeder.de>
22 @brief Component to run the FastJet jetfinder
29 #include "AliHLTJETFastJetComponent.h"
32 #include "TObjString.h"
36 /** ROOT macro for the implementation of ROOT specific class methods */
37 ClassImp(AliHLTJETFastJetComponent)
40 * ---------------------------------------------------------------------------------
41 * Constructor / Destructor
42 * ---------------------------------------------------------------------------------
45 // #################################################################################
46 AliHLTJETFastJetComponent::AliHLTJETFastJetComponent()
51 fJetReaderHeader(NULL),
55 // see header file for class documentation
57 // refer to README to build package
59 // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt
62 // #################################################################################
63 AliHLTJETFastJetComponent::~AliHLTJETFastJetComponent() {
64 // see header file for class documentation
69 * ---------------------------------------------------------------------------------
70 * Public functions to implement AliHLTComponent's interface.
71 * These functions are required for the registration process
72 * ---------------------------------------------------------------------------------
75 // #################################################################################
76 const Char_t* AliHLTJETFastJetComponent::GetComponentID() {
77 // see header file for class documentation
78 return "JETFastJetFinder";
81 // #################################################################################
82 void AliHLTJETFastJetComponent::GetInputDataTypes( vector<AliHLTComponentDataType>& list) {
83 // see header file for class documentation
85 list.push_back( kAliHLTDataTypeMCObject|kAliHLTDataOriginHLT );
86 list.push_back( kAliHLTDataTypeESDObject|kAliHLTDataOriginOffline );
87 list.push_back( kAliHLTDataTypeESDObject|kAliHLTDataOriginHLT );
90 // #################################################################################
91 AliHLTComponentDataType AliHLTJETFastJetComponent::GetOutputDataType() {
92 // see header file for class documentation
93 return (kAliHLTDataTypeJet|kAliHLTDataOriginHLT);
96 // #################################################################################
97 void AliHLTJETFastJetComponent::GetOutputDataSize( ULong_t& constBase, Double_t& inputMultiplier ) {
98 // see header file for class documentation
101 inputMultiplier = 0.3;
104 // #################################################################################
105 AliHLTComponent* AliHLTJETFastJetComponent::Spawn() {
106 // see header file for class documentation
107 return new AliHLTJETFastJetComponent();
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 * ---------------------------------------------------------------------------------
118 // #################################################################################
119 Int_t AliHLTJETFastJetComponent::DoInit( Int_t argc, const Char_t** argv ) {
120 // see header file for class documentation
122 if ( fJetFinder || fJetHeader || fJetReader || fJetReaderHeader ||
123 fTrackCuts || fJetCuts || fJets )
126 // ---------------------------------------------------------------------
128 // ---------------------------------------------------------------------
130 TString comment = "HLT FastJet interface";
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;
137 // ---------------------------------------------------------------------
139 // ---------------------------------------------------------------------
142 Int_t bMissingParam=0;
146 // -- Loop over all arguments
147 for ( Int_t iter = 0; iter<argc && iResult>=0; iter++) {
150 if (argument.IsNull())
154 if ( !argument.CompareTo("-algorithm") ) {
155 if ((bMissingParam=(++iter>=argc))) break;
157 TString parameter(argv[iter]);
158 parameter.Remove(TString::kLeading, ' ');
160 if ( !parameter.CompareTo("Kt") ) {
161 algorithm = AliHLTJETBase::kKt;
164 comment += parameter;
167 else if ( !parameter.CompareTo("AntiKt") ) {
168 algorithm = AliHLTJETBase::kAntiKt;
171 comment += parameter;
175 HLTError("Wrong parameter %s for argument %s.", parameter.Data(), argument.Data());
181 else if ( !argument.CompareTo("-coneRadius") ) {
182 if ((bMissingParam=(++iter>=argc))) break;
184 TString parameter(argv[iter]);
185 parameter.Remove(TString::kLeading, ' ');
187 if ( parameter.IsFloat() ) {
188 coneRadius = parameter.Atof();
191 comment += parameter;
195 HLTError("Wrong parameter %s for argument %s.", parameter.Data(), argument.Data());
201 else if ( !argument.CompareTo("-trackCutMinPt") ) {
202 if ((bMissingParam=(++iter>=argc))) break;
204 TString parameter(argv[iter]);
205 parameter.Remove(TString::kLeading, ' ');
207 if ( parameter.IsFloat() ) {
208 trackCutMinPt = parameter.Atof();
211 comment += parameter;
215 HLTError("Wrong parameter %s for argument %s.", parameter.Data(), argument.Data());
221 else if ( !argument.CompareTo("-jetCutMinEt") ) {
222 if ((bMissingParam=(++iter>=argc))) break;
224 TString parameter(argv[iter]);
225 parameter.Remove(TString::kLeading, ' ');
227 if ( parameter.IsFloat() ) {
228 jetCutMinEt = parameter.Atof();
231 comment += parameter;
235 HLTError("Wrong parameter %s for argument %s.", parameter.Data(), argument.Data());
240 // -- Argument not known
242 HLTError("Unknown argument %s.", argument.Data());
245 } // for ( Int iter = 0; iter<argc && iResult>=0; iter++) {
247 // -- Check if parameter is missing
248 if ( bMissingParam ) {
249 HLTError("Missing parameter for argument %s.", argument.Data());
256 // ---------------------------------------------------------------------
258 // ---------------------------------------------------------------------
259 if ( ! (fTrackCuts = new AliHLTJETTrackCuts()) ) {
260 HLTError("Error instantiating track cuts");
264 fTrackCuts->SetChargedOnly( kTRUE );
265 fTrackCuts->SetMinPt( trackCutMinPt );
267 // ---------------------------------------------------------------------
269 // ---------------------------------------------------------------------
270 if ( ! (fJetCuts = new AliHLTJETJetCuts()) ) {
271 HLTError("Error instantiating jet cuts");
275 fJetCuts->SetMinEt( jetCutMinEt );
277 // ---------------------------------------------------------------------
278 // -- Jet Reader Header
279 // ---------------------------------------------------------------------
280 if ( ! (fJetReaderHeader = new AliHLTJETReaderHeader()) ) {
281 HLTError("Error instantiating jet reader header");
286 fJetReaderHeader->SetJetAlgorithm( algorithm );
288 // Set prt to track cuts
289 fJetReaderHeader->SetTrackCuts( fTrackCuts );
291 // Set Eta min/max and Phi min/max
292 fJetReaderHeader->SetFiducialEta( -0.9, 0.9) ;
293 fJetReaderHeader->SetFiducialPhi( 0.0, TMath::TwoPi() ) ;
296 fJetReaderHeader->SetConeRadius(coneRadius);
298 // ---------------------------------------------------------------------
300 // ---------------------------------------------------------------------
301 if ( ! (fJetReader = new AliHLTJETReader()) ) {
302 HLTError("Error instantiating jet reader");
306 fJetReader->SetReaderHeader(fJetReaderHeader);
308 // ---------------------------------------------------------------------
310 // ---------------------------------------------------------------------
311 if ( ! (fJets = new AliHLTJets()) ) {
312 HLTError("Error instantiating jet container");
316 fJets->SetComment(comment);
318 // ---------------------------------------------------------------------
320 // ---------------------------------------------------------------------
321 if ( ! (fJetHeader = new AliHLTJETFastJetHeader()) ) {
322 HLTError("Error instantiating fastjet header");
326 fJetHeader->SetReaderHeader(fJetReaderHeader);
327 fJetHeader->SetJetCuts(fJetCuts);
329 // ---------------------------------------------------------------------
331 // ---------------------------------------------------------------------
332 if ( ! (fJetFinder = new AliHLTJETFastJetFinder()) ) {
333 HLTError("Error instantiating fastjet finder");
337 fJetFinder->SetJetHeader(fJetHeader);
338 fJetFinder->SetJetReader(fJetReader);
339 fJetFinder->SetOutputJets(fJets);
341 // ---------------------------------------------------------------------
342 // -- Initialize Jet Finder
343 // ---------------------------------------------------------------------
344 if ( (fJetFinder->Initialize()) ) {
345 HLTError("Error initializing fastjet finder");
352 // #################################################################################
353 Int_t AliHLTJETFastJetComponent::DoDeinit() {
354 // see header file for class documentation
368 if ( fJetReaderHeader )
369 delete fJetReaderHeader;
370 fJetReaderHeader = NULL;
383 // #################################################################################
384 Int_t AliHLTJETFastJetComponent::DoEvent( const AliHLTComponentEventData& /*evtData*/,
385 AliHLTComponentTriggerData& /*trigData*/ ) {
386 // see header file for class documentation
390 const TObject* iter = NULL;
394 if ( GetFirstInputObject(kAliHLTDataTypeSOR) && !iResult ) {
395 HLTInfo("On-line SOR Event");
398 // -- ADD MC Object -- On-line
399 // ------------------------------
400 for ( iter=GetFirstInputObject(kAliHLTDataTypeMCObject|kAliHLTDataOriginHLT);
401 iter != NULL && !iResult; iter=GetNextInputObject() ) {
403 // -- Set automatic MC usage, --> needed in off-line
404 fJetReaderHeader->SetUseMC(kTRUE);
406 // -- Set input event
407 fJetReader->SetInputEvent( NULL, NULL, const_cast<TObject*>(iter) );
409 // -- Fill vector with MC
410 if ( ! fJetReader->FillVectorHLTMC() ) {
411 HLTError("Error filling vector.");
412 iResult = -EINPROGRESS;
417 if ( ! fJetFinder->ProcessHLTEvent() ) {
418 HLTError("Error processing fastjet event.");
419 iResult = -EINPROGRESS;
425 PushBack(fJets, kAliHLTDataTypeJet|kAliHLTDataOriginHLT, GetSpecification());
429 // -- ADD ESD Object -- Off-line
430 // -------------------------------
431 for ( iter=GetFirstInputObject(kAliHLTDataTypeESDObject|kAliHLTDataOriginOffline);
432 iter != NULL && !iResult; iter=GetNextInputObject() ) {
434 // -- Set automatic MC usage, --> needed in off-line
435 fJetReaderHeader->SetUseMC(kFALSE);
437 // -- Set input event
438 fJetReader->SetInputEvent( const_cast<TObject*>(iter), NULL, NULL );
440 // -- Fill vector with ESD
441 if ( ! fJetReader->FillVectorESD() ) {
442 HLTError("Error filling vector.");
448 if ( ! fJetFinder->ProcessHLTEvent() ) {
449 HLTError("Error processing fastjet event.");
450 iResult = -EINPROGRESS;
456 PushBack(fJets, kAliHLTDataTypeJet|kAliHLTDataOriginHLT, GetSpecification());
460 // -- ADD ESD Object -- On-line
461 // ------------------------------
462 for ( iter=GetFirstInputObject(kAliHLTDataTypeESDObject|kAliHLTDataOriginHLT);
463 iter != NULL && !iResult; iter=GetNextInputObject() ) {
465 // -- Set automatic MC usage, --> needed in off-line
466 fJetReaderHeader->SetUseMC(kFALSE);
468 // -- Set input event
469 fJetReader->SetInputEvent( const_cast<TObject*>(iter), NULL, NULL );
471 // -- Fill vector with ESD
472 if ( ! fJetReader->FillVectorESD() ) {
473 HLTError("Error filling vector.");
479 if ( ! fJetFinder->ProcessHLTEvent() ) {
480 HLTError("Error processing fastjet event.");
481 iResult = -EINPROGRESS;
487 PushBack(fJets, kAliHLTDataTypeJet|kAliHLTDataOriginHLT, GetSpecification());