]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/trigger/AliHLTTriggerFastJet.cxx
jet trigger component and jet wrapper for EMCAL and relative package files
[u/mrichter/AliRoot.git] / HLT / trigger / AliHLTTriggerFastJet.cxx
1 //**************************************************************************                        
2 //* This file is property of and copyright by the ALICE HLT Project        *                        
3 //* ALICE Experiment at CERN, All rights reserved.                         *                        
4 //*                                                                        *                        
5 //* Primary Authors: leonidas.xaplanteris@gmail.com                        *                        
6 //*                  for The ALICE HLT Project.                            *                        
7 //*                                                                        *                        
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 //**************************************************************************   
16
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"
26 #include "TString.h"
27 #include "TMap.h"
28
29 #include <fastjet/PseudoJet.hh>
30 #include <fastjet/JetDefinition.hh>
31 #include <fastjet/AreaDefinition.hh>
32
33 #include "AliFJWrapper.h"
34
35
36 ClassImp(AliHLTTriggerFastJet)
37
38 AliHLTTriggerFastJet::AliHLTTriggerFastJet(TString detector) :
39 AliHLTTrigger(),
40   fEThreshold(0.0),
41   fDetector(detector),
42   fFastJetWrapper(NULL),
43   fClustersRefs(NULL),
44   fOCDBEntry(""),
45   fOffset(10000),
46   fInputDataType()
47 {
48   fFastJetWrapper = new AliFJWrapper("FastJet","FastJet");
49
50   fClustersRefs = new TRefArray();
51
52 }
53 //_____________________________________________________________
54
55 AliHLTTriggerFastJet::~AliHLTTriggerFastJet() {
56
57   if (fFastJetWrapper)
58     delete fFastJetWrapper;
59   fFastJetWrapper = NULL;
60
61   if (fClustersRefs)
62     delete fClustersRefs;
63   fClustersRefs = NULL;
64
65 }
66 //_____________________________________________________________
67
68 int AliHLTTriggerFastJet::DoInit(int argc, const char** argv) {
69   
70   int iResult = ConfigureFromCDBTObjString(fOCDBEntry);
71   
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 );
77   }
78   return iResult;
79 }
80 //______________________________________________________________
81
82 int AliHLTTriggerFastJet::DoDeInit() {
83   
84   return 0;
85
86 }
87 //______________________________________________________________
88
89 Int_t AliHLTTriggerFastJet::DoTrigger() {
90
91   Int_t iResult = 0;
92
93   if ( GetFirstInputBlock( kAliHLTDataTypeSOR ) || GetFirstInputBlock( kAliHLTDataTypeEOR ) )
94     return 0;
95
96   const TObject *obj = GetFirstInputObject( kAliHLTAllDataTypes, "AliESDEvent" );
97   AliESDEvent   *esd = dynamic_cast<AliESDEvent*>(const_cast<TObject*>(obj));
98   
99   // check for MatchedTrack to avoid double counting
100   
101
102   if ( esd != NULL ) {
103     esd->GetStdContent();
104     
105     Double_t vertex[3] = {0,0,0};
106     esd->GetVertex()->GetXYZ(vertex);
107     TLorentzVector gamma;
108     
109     // -- add VZERO
110     AliESDVZERO *v0  = esd->GetVZEROData();
111     Float_t MtotVOA = v0->GetMTotV0A();
112     Float_t MtotVOC = v0->GetMTotV0A();
113
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());
119       if ( tpctrack )
120         fFastJetWrapper->AddInputVector(tpctrack->Px(),tpctrack->Py(),tpctrack->Pz(),tpctrack->P());
121       delete tpctrack;
122     }
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);
129         fOffset++;
130         gamma.Clear();
131       }
132     }
133     
134     
135   }
136   
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);
145   
146   double median = 0.0;
147   double sigma  = 0.0;
148
149   fFastJetWrapper->Run();
150   fFastJetWrapper->GetMedianAndSigma(median,sigma);
151
152   if ( TriggerOnJet(fFastJetWrapper->GetSubtractedJetsPts(median)) )
153        return iResult;
154   
155   TString description;
156   description.Form(" No %s jets with energy >  %.02f GeV found! ", fDetector.Data(), fEThreshold);
157   SetDescription(description.Data());
158   TriggerEvent(kFALSE);
159   
160   return iResult;
161
162 }
163 //______________________________________________________________
164
165 template <class T>
166 Bool_t AliHLTTriggerFastJet::TriggerOnJet(T Jet) {
167
168   for ( unsigned int ij = 0; ij<Jet.size(); ij++ ) {
169     if ( Jet[ij] > fEThreshold ) {
170       TString description;
171       description.Form(" Event contains at least one %s jet with energy > %.02f GeV! ", fDetector.Data(), fEThreshold);
172       SetDescription(description.Data());
173       
174       // Enable the detectors for readout
175       GetReadoutList().Enable(AliHLTReadoutList::kEMCAL |
176                               AliHLTReadoutList::kTPC);
177       
178       // Add the available HLT info for readout
179       GetTriggerDomain().Add(kAliHLTAnyDataTypeID, fDetector.Data());
180       
181       // Set trigger desicion
182       TriggerEvent(kTRUE);
183       
184       return(kTRUE);
185     }
186   }
187   
188   return kFALSE;
189
190 }
191 //______________________________________________________________
192
193 int AliHLTTriggerFastJet::Reconfigure(const char* cdbEntry, const char* /*chainId*/) {
194   
195   const char* entry = cdbEntry;
196   if ( !entry || entry[0]==0 ) entry=fOCDBEntry;
197
198   return ConfigureFromCDBTObjString(entry);
199
200 }
201
202 //______________________________________________________________
203
204 int AliHLTTriggerFastJet::ScanConfigurationArgument(int argc, const char** argv) {
205
206   if ( argc<=0 ) return 0;
207   int i = 0;
208   TString argument = argv[i];
209
210   if ( argument.CompareTo("-energy") == 0 ) {
211     if (++i>=argc) return -EPROTO;
212     argument = argv[i];
213     fEThreshold=argument.Atof();
214     return 2;
215   }
216   
217   return -EINVAL;
218
219 }
220 //______________________________________________________________
221
222 void AliHLTTriggerFastJet::GetOutputDataSize(unsigned long &constBase, double &inputMultiplier) {
223   
224   constBase = sizeof(AliHLTTriggerDecision) + sizeof(AliHLTDomainEntry)*14;
225   inputMultiplier = 1;
226
227 }
228 //______________________________________________________________
229
230 void AliHLTTriggerFastJet::GetOCDBObjectDescription(TMap* const targetMap) {
231   
232   if ( !targetMap ) return;
233   targetMap->Add(new TObjString(fOCDBEntry),
234                  new TObjString(Form("%s threshold trigger OCDB object", fDetector.Data()) ) 
235                  );
236
237 }
238 //_______________________________________________________________