]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/trigger/AliHLTTriggerEmcalElectron.cxx
flat friends update
[u/mrichter/AliRoot.git] / HLT / trigger / AliHLTTriggerEmcalElectron.cxx
1 // $Id: AliHLTTriggerEmcalElectron.cxx 50471 2011-07-07 09:50:47Z fronchet $
2 //**************************************************************************
3 //* This file is property of and copyright by the ALICE HLT Project        * 
4 //* ALICE Experiment at CERN, All rights reserved.                         *
5 //*                                                                        *
6 //* Primary Authors: marcelfigueredo@gmail.com                             *
7 //*                  for The ALICE HLT Project.                            *
8 //*                                                                        *
9 //* Permission to use, copy, modify and distribute this software and its   *
10 //* documentation strictly for non-commercial purposes is hereby granted   *
11 //* without fee, provided that the above copyright notice appears in all   *
12 //* copies and that both the copyright notice and this permission notice   *
13 //* appear in the supporting documentation. The authors make no claims     *
14 //* about the suitability of this software for any purpose. It is          *
15 //* provided "as is" without express or implied warranty.                  *
16 //**************************************************************************
17
18 #include "AliHLTTriggerEmcalElectron.h"
19 #include "AliESDEvent.h"
20 #include "AliVCluster.h" 
21 #include "AliESDCaloCluster.h"
22 #include "AliHLTTriggerDecision.h"
23 #include "AliHLTDomainEntry.h"
24 #include "AliHLTCaloClusterReader.h"
25 #include "AliHLTCaloClusterDataStruct.h"
26 #include "TRefArray.h"
27 #include "TString.h"
28 #include "TMap.h"
29 #include "AliESDtrack.h"
30 #include "AliHLTScalars.h"
31 /** ROOT macro for the implementation of ROOT specific class methods */
32 ClassImp(AliHLTTriggerEmcalElectron)
33
34 AliHLTTriggerEmcalElectron::AliHLTTriggerEmcalElectron() : 
35   AliHLTTrigger(),
36   fEThreshold(0.0),
37   fEoverPThreshold(0.),
38   fEoverPLimit(0.),
39   fdEta(1.),
40   fdPhi(1.),
41   
42 //   fClustersRefs(NULL),
43 //   fClusterReader(NULL),
44   fOCDBEntry("HLT/ConfigHLT/EmcalElectronTrigger"),
45 //   fOCDBEntry(""), 
46   fDetector("EMCAL"),
47   fInputDataType(),
48   fMakeStats(kFALSE)
49 {
50   // see header file for class documentation
51   // or
52   // refer to README to build package
53   // or
54   // visit http://web.ift.uib.no/~kjeks/doc/alice-hlts
55
56   //if ( fMakeStats ) AliHLTScalars scalars;
57   //AliHLTScalars scalars;
58 //   fClusterReader = new AliHLTCaloClusterReader();
59 //   fClustersRefs = new TRefArray();
60
61 }
62
63
64 AliHLTTriggerEmcalElectron::~AliHLTTriggerEmcalElectron() {
65   // see header file for class documentation
66 //   if (fClusterReader)
67 //     delete fClusterReader;
68 //   fClusterReader = NULL;
69 // 
70 //   if(fClustersRefs)
71 //     delete fClustersRefs;
72 //   fClustersRefs = NULL;
73 }
74
75 //Trigger name
76 const char* AliHLTTriggerEmcalElectron::GetTriggerName() const {
77   // see header file for class documentation
78   return "EmcalElectronTrigger";
79 }
80
81 AliHLTComponent* AliHLTTriggerEmcalElectron::Spawn() {
82   // see header file for class documentation
83   return new AliHLTTriggerEmcalElectron;
84 }
85
86 Int_t AliHLTTriggerEmcalElectron::DoTrigger() {
87   // see header file for class documentation
88   
89   Int_t iResult = 0;
90   AliHLTScalars scalars;
91
92   if ( GetFirstInputBlock( kAliHLTDataTypeSOR ) || GetFirstInputBlock( kAliHLTDataTypeEOR ) )
93     return 0;
94
95   //Try the caloclusterstruct input
96
97 /*  
98   for (const AliHLTComponentBlockData* pBlock=GetFirstInputBlock(fInputDataType); pBlock!=NULL; pBlock=GetNextInputBlock()) {
99     AliHLTCaloClusterHeaderStruct *caloClusterHeader = reinterpret_cast<AliHLTCaloClusterHeaderStruct*>(pBlock->fPtr);
100     fClusterReader->SetMemory(caloClusterHeader);
101     
102     AliHLTCaloClusterDataStruct * caloClusterStruct;
103     while( (caloClusterStruct = fClusterReader->NextCluster()) != 0) {
104       if (TriggerOnEoverP(caloClusterStruct)) {
105         return iResult;
106       }
107     }
108   }*/
109
110   //Try the ESD input
111   const TObject* obj = GetFirstInputObject(kAliHLTAllDataTypes, "AliESDEvent");
112   AliESDEvent* esd = dynamic_cast<AliESDEvent*>(const_cast<TObject*>(obj));
113   
114   if (esd != NULL) {
115     esd->GetStdContent();
116
117 //     Int_t ncc = GetClustersFromEsd(esd, fClustersRefs); //marcel test
118     Int_t ncc = esd->GetNumberOfCaloClusters();  
119 //     Int_t ncc = GetClustersFromEsd(esd, fClustersRefs); 
120     
121     for (Int_t i = 0; i < ncc ; i++) {
122       
123 //       AliESDCaloCluster * cluster = static_cast<AliESDCaloCluster*>(fClustersRefs->At(i));
124       AliVCluster *cluster = (AliVCluster*) esd->GetCaloCluster(i); //MARCEL test
125       if(TriggerOnEoverP(cluster,esd)) { 
126         return iResult;
127       }
128     }
129   }
130
131   // If we got to this point then we did not find any clusters with E > fEThreshold
132   // generate negative trigger decision
133   TString description;
134   description.Form("No %s clusters corresponding to Energy >  %.02f GeV and %.02f < E/P < %.02f", fDetector.Data(),fEThreshold,fEoverPThreshold,fEoverPLimit);
135   SetDescription(description.Data());
136   TriggerEvent(false);
137   
138   if ( fMakeStats ) iResult = PushBack(&scalars, kAliHLTDataTypeEventStatistics|kAliHLTDataOriginHLT);
139   return iResult;
140
141 }
142
143
144 template <class T>
145 Bool_t AliHLTTriggerEmcalElectron::TriggerOnEoverP(T* cluster,AliESDEvent *esd) {
146   
147   AliHLTScalars scalars;
148   
149   if (cluster->E() > fEThreshold) {    
150
151       Int_t trackindex=cluster->GetTrackMatchedIndex(); 
152       if(trackindex<0)return kFALSE;
153
154       Double_t dEta=cluster->GetTrackDz();
155       Double_t dPhi=cluster->GetTrackDx();
156
157       AliESDtrack* track = esd->GetTrack(trackindex);
158       if(!track)return kFALSE;     
159       Double_t EoverP=cluster->E()/track->P(); 
160         
161       if ( fMakeStats ) {
162         Double_t deltaR=TMath::Sqrt(dEta*dEta+dPhi*dPhi);
163         scalars.Add("dR","Residuals dR of track matching", deltaR);
164         scalars.Add("dEta","Residuals dEta of track matching", dEta);
165         scalars.Add("dPhi","Residuals dPhi of track matching", dPhi);
166         }
167       
168       if(TMath::Abs(dEta)>fdEta)return kFALSE;
169       if(TMath::Abs(dPhi)>fdPhi)return kFALSE;
170
171       if ( fMakeStats ) {  
172         scalars.Add("TracksPt","TPC tracks pT", track->Pt());
173         scalars.Add("ClusterEn","Cluster Energy",cluster->E());
174         scalars.Add("EoverP","EoverP for matched tracks",EoverP);
175         }
176             
177       if(EoverP<fEoverPThreshold)return kFALSE;
178       if(EoverP>fEoverPLimit)return kFALSE; 
179  
180       //We have a cluster satisfying trigger criteria
181     TString description;
182     description.Form("Event contains at least one %s cluster with energy greater than %.02f corresponding to  %.02f < E/P < %.02f, residuals: dPhi < %.02f dEta < %.02f", fDetector.Data(),fEThreshold,fEoverPThreshold,fEoverPLimit,fdPhi,fdEta);
183     SetDescription(description.Data());
184     
185     // Enable the detectors for readout.
186
187     //GetReadoutList().Enable(AliHLTReadoutList::kPHOS);
188 //     SetCaloReadoutList("EMCAL");  //FR
189 //     SetCaloReadoutList();  //FR
190      
191     // Add the available HLT information for readout too.
192     GetTriggerDomain().Add(kAliHLTAnyDataTypeID, fDetector.Data());
193     
194     //Set trigger decision
195     TriggerEvent(kTRUE);
196     
197     return kTRUE;
198   } 
199
200
201   return kFALSE;
202
203 }
204
205 int AliHLTTriggerEmcalElectron::DoInit(int argc, const char** argv) {
206   // see header file for class documentation
207 //   return 0;// marcel test
208   // first configure the default
209   int iResult=ConfigureFromCDBTObjString(fOCDBEntry);
210
211   // configure from the command line parameters if specified
212   if (iResult>=0 && argc>0) {
213     iResult=ConfigureFromArgumentString(argc, argv);
214     HLTImportant("Trigger threshold set from argument string:  %.02f GeV: E/P:%.02f and %.02f, residuals: dPhi:%.02f dEta:%.02f", fDetector.Data(),fEThreshold,fEoverPThreshold,fEoverPLimit,fdPhi,fdEta);   
215   } else if ( iResult >=0 ) {
216     HLTImportant("Trigger threshold set from OCDB database entry:  %.02f GeV: E/P:%.02f and %.02f, residuals: dPhi:%.02f dEta:%.02f", fDetector.Data(),fEThreshold,fEoverPThreshold,fEoverPLimit,fdPhi,fdEta);
217   }
218   return iResult;
219 }
220
221 int AliHLTTriggerEmcalElectron::DoDeinit() {
222
223   // see header file for class documentation
224  
225   return 0;
226 }
227
228 int AliHLTTriggerEmcalElectron::Reconfigure(const char* cdbEntry, const char* /*chainId*/) {
229   // see header file for class documentation
230
231   // configure from the specified antry or the default one
232   const char* entry=cdbEntry;
233   if (!entry || entry[0]==0) entry=fOCDBEntry;
234
235   return ConfigureFromCDBTObjString(entry);
236 }
237
238 int AliHLTTriggerEmcalElectron::ScanConfigurationArgument(int argc, const char** argv) {
239   // see header file for class documentation
240   if (argc<=0) return 0;
241   int i=0;
242   TString argument=argv[i];
243
244   // -maxpt
245   if (argument.CompareTo("-energy")==0) {
246     if (++i>=argc) return -EPROTO;
247     argument=argv[i];
248     fEThreshold=argument.Atof(); // 
249     return 2;
250   }    
251
252   if (argument.CompareTo("-minEoverP")==0) {
253     if (++i>=argc) return -EPROTO;
254     argument=argv[i];
255     fEoverPThreshold=argument.Atof(); // 
256     return 2;
257   } 
258   
259     if (argument.CompareTo("-maxEoverP")==0) {
260     if (++i>=argc) return -EPROTO;
261     argument=argv[i];
262     fEoverPLimit=argument.Atof(); // 
263     return 2;
264   }
265   
266     if (argument.CompareTo("-dEta")==0) {
267     if (++i>=argc) return -EPROTO;
268     argument=argv[i];
269     fdEta=argument.Atof(); // 
270     return 2;
271   } 
272   
273     if (argument.CompareTo("-dPhi")==0) {
274     if (++i>=argc) return -EPROTO;
275     argument=argv[i];
276     fdPhi=argument.Atof(); // 
277     return 2;
278   } 
279   
280     if (argument.CompareTo("-makestats")==0) {
281     fMakeStats = kTRUE;
282     return 2;
283   }
284   
285 // unknown argument
286   return -EINVAL;
287 }
288
289 //______________________________________________________________
290
291 void AliHLTTriggerEmcalElectron::GetOutputDataTypes(AliHLTComponentDataTypeList &list) const {
292   // return the output data types generated
293   
294   list.push_back(kAliHLTDataTypeTriggerDecision);
295   list.push_back(kAliHLTDataTypeEventStatistics|kAliHLTDataOriginHLT);
296 }
297 //_______________________________________________________________
298 // 
299 void AliHLTTriggerEmcalElectron::GetOutputDataSize(unsigned long& constBase, double& inputMultiplier) {
300   // see header file for documentation
301   constBase = sizeof(AliHLTTriggerDecision) + sizeof(AliHLTDomainEntry)*14;
302   inputMultiplier = 1;
303 }
304
305
306 void AliHLTTriggerEmcalElectron::GetOCDBObjectDescription( TMap* const targetMap) {
307   
308   // Get a list of OCDB object description.
309   if (!targetMap) return;
310   targetMap->Add(new TObjString(fOCDBEntry),
311                  new TObjString(Form("%s threshold trigger OCDB object", fDetector.Data()) ) 
312                  );
313 }
314
315