1 //**************************************************************************
2 //* This file is property of and copyright by the ALICE HLT Project *
3 //* ALICE Experiment at CERN, All rights reserved. *
5 //* Primary Authors: leonidas.xaplanteris@gmail.com *
6 //* for The ALICE HLT Project. *
8 //* Permission to use, copy, modify and distribute this software and its *
9 //* documentation strictly for non-commercial purposes is hereby granted *
10 //* without fee, provided that the above copyright notice appears in all *
11 //* copies and that both the copyright notice and this permission notice *
12 //* appear in the supporting documentation. The authors make no claims *
13 //* about the suitability of this software for any purpose. It is *
14 //* provided "as is" without express or implied warranty. *
15 //**************************************************************************
17 #include "AliHLTTriggerFastJet.h"
18 #include "AliESDEvent.h"
19 #include "AliESDtrackCuts.h"
20 #include "AliESDCaloCluster.h"
21 #include "AliESDVZERO.h"
22 #include "AliHLTTriggerDecision.h"
23 #include "AliHLTDomainEntry.h"
24 #include "TLorentzVector.h"
25 #include "TRefArray.h"
29 #include <fastjet/PseudoJet.hh>
30 #include <fastjet/JetDefinition.hh>
31 #include <fastjet/AreaDefinition.hh>
33 #include "AliFJWrapper.h"
36 ClassImp(AliHLTTriggerFastJet)
38 AliHLTTriggerFastJet::AliHLTTriggerFastJet(TString detector) :
42 fFastJetWrapper(NULL),
48 fFastJetWrapper = new AliFJWrapper("FastJet","FastJet");
50 fClustersRefs = new TRefArray();
53 //_____________________________________________________________
55 AliHLTTriggerFastJet::~AliHLTTriggerFastJet() {
58 delete fFastJetWrapper;
59 fFastJetWrapper = NULL;
66 //_____________________________________________________________
68 int AliHLTTriggerFastJet::DoInit(int argc, const char** argv) {
70 int iResult = ConfigureFromCDBTObjString(fOCDBEntry);
72 if ( iResult>=0 && argc>0 ) {
73 iResult = ConfigureFromArgumentString(argc, argv);
74 HLTImportant("Trigger threshold set from argument string : %.02f GeV:", fEThreshold );
75 } else if ( iResult>=0 ) {
76 HLTImportant("Trigger threshold set from OCDB database entry : %.02f Gev:", fEThreshold );
80 //______________________________________________________________
82 int AliHLTTriggerFastJet::DoDeInit() {
87 //______________________________________________________________
89 Int_t AliHLTTriggerFastJet::DoTrigger() {
93 if ( GetFirstInputBlock( kAliHLTDataTypeSOR ) || GetFirstInputBlock( kAliHLTDataTypeEOR ) )
96 const TObject *obj = GetFirstInputObject( kAliHLTAllDataTypes, "AliESDEvent" );
97 AliESDEvent *esd = dynamic_cast<AliESDEvent*>(const_cast<TObject*>(obj));
99 // check for MatchedTrack to avoid double counting
103 esd->GetStdContent();
105 Double_t vertex[3] = {0,0,0};
106 esd->GetVertex()->GetXYZ(vertex);
107 TLorentzVector gamma;
110 AliESDVZERO *v0 = esd->GetVZEROData();
111 Float_t MtotVOA = v0->GetMTotV0A();
112 Float_t MtotVOC = v0->GetMTotV0A();
114 // add the tracks to FastJet Wrapper
115 for ( Int_t i = 0; i<esd->GetNumberOfTracks(); i++ ) {
116 AliESDtrack *esdtrack = esd->GetTrack(i);
117 fFastJetWrapper->AddInputVector(esdtrack->Px(),esdtrack->Py(),esdtrack->Pz(),esdtrack->P());
118 AliESDtrack *tpctrack = AliESDtrackCuts::GetTPCOnlyTrack(esd,esdtrack->GetID());
120 fFastJetWrapper->AddInputVector(tpctrack->Px(),tpctrack->Py(),tpctrack->Pz(),tpctrack->P());
123 // add the Calo Clusters to FastJet Wrapper
124 for ( Int_t j = 0; j<GetClustersFromEsd(esd,fClustersRefs); j++) {
125 AliESDCaloCluster *cluster = static_cast<AliESDCaloCluster*>(fClustersRefs->At(j));
126 if ( cluster->IsEMCAL() ) {
127 cluster->GetMomentum(gamma,vertex);
128 fFastJetWrapper->AddInputVector(gamma.Px(),gamma.Py(),gamma.Pz(),gamma.P(),fOffset);
137 fFastJetWrapper->SetStrategy(fastjet::Best);
138 fFastJetWrapper->SetAlgorithm(fastjet::antikt_algorithm);
139 fFastJetWrapper->SetRecombScheme(fastjet::BIpt_scheme);
140 fFastJetWrapper->SetAreaType(fastjet::active_area);
141 fFastJetWrapper->SetNRepeats(1);
142 fFastJetWrapper->SetGhostArea(0.01);
143 fFastJetWrapper->SetMaxRap(0.9);
144 fFastJetWrapper->SetR(0.4);
149 fFastJetWrapper->Run();
150 fFastJetWrapper->GetMedianAndSigma(median,sigma);
152 if ( TriggerOnJet(fFastJetWrapper->GetSubtractedJetsPts(median)) )
156 description.Form(" No %s jets with energy > %.02f GeV found! ", fDetector.Data(), fEThreshold);
157 SetDescription(description.Data());
158 TriggerEvent(kFALSE);
163 //______________________________________________________________
166 Bool_t AliHLTTriggerFastJet::TriggerOnJet(T Jet) {
168 for ( unsigned int ij = 0; ij<Jet.size(); ij++ ) {
169 if ( Jet[ij] > fEThreshold ) {
171 description.Form(" Event contains at least one %s jet with energy > %.02f GeV! ", fDetector.Data(), fEThreshold);
172 SetDescription(description.Data());
174 // Enable the detectors for readout
175 GetReadoutList().Enable(AliHLTReadoutList::kEMCAL |
176 AliHLTReadoutList::kTPC);
178 // Add the available HLT info for readout
179 GetTriggerDomain().Add(kAliHLTAnyDataTypeID, fDetector.Data());
181 // Set trigger desicion
191 //______________________________________________________________
193 int AliHLTTriggerFastJet::Reconfigure(const char* cdbEntry, const char* /*chainId*/) {
195 const char* entry = cdbEntry;
196 if ( !entry || entry[0]==0 ) entry=fOCDBEntry;
198 return ConfigureFromCDBTObjString(entry);
202 //______________________________________________________________
204 int AliHLTTriggerFastJet::ScanConfigurationArgument(int argc, const char** argv) {
206 if ( argc<=0 ) return 0;
208 TString argument = argv[i];
210 if ( argument.CompareTo("-energy") == 0 ) {
211 if (++i>=argc) return -EPROTO;
213 fEThreshold=argument.Atof();
220 //______________________________________________________________
222 void AliHLTTriggerFastJet::GetOutputDataSize(unsigned long &constBase, double &inputMultiplier) {
224 constBase = sizeof(AliHLTTriggerDecision) + sizeof(AliHLTDomainEntry)*14;
228 //______________________________________________________________
230 void AliHLTTriggerFastJet::GetOCDBObjectDescription(TMap* const targetMap) {
232 if ( !targetMap ) return;
233 targetMap->Add(new TObjString(fOCDBEntry),
234 new TObjString(Form("%s threshold trigger OCDB object", fDetector.Data()) )
238 //_______________________________________________________________