]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/trigger/AliHLTTriggerBarrelGeomMultiplicity.cxx
- fixing compilation warning
[u/mrichter/AliRoot.git] / HLT / trigger / AliHLTTriggerBarrelGeomMultiplicity.cxx
1 // $Id$
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: Oystein Djuvsland                                     *
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 /// @file   AliHLTTriggerBarrelGeomMultiplicity.cxx
19 /// @author Oystein Djuvsland
20 /// @date   2009-10-08
21 /// @brief  HLT trigger component for charged particle multiplicity 
22 ///         within a geometrical acceptance in the central barrel.
23
24 // see header file for class documentation
25 // or
26 // refer to README to build package
27 // or
28 // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt
29
30 #include "AliHLTTriggerBarrelGeomMultiplicity.h"
31 #include "AliHLTTriggerDetectorGeom.h"
32 #include "AliHLTTriggerDecisionParameters.h"
33 #include "AliESDEvent.h"
34 #include "AliHLTTriggerDecision.h"
35 #include "AliHLTDomainEntry.h"
36 #include "AliHLTGlobalBarrelTrack.h"
37 #include "TObjArray.h"
38 #include "TObjString.h"
39 #include "AliCDBEntry.h"
40 #include "AliCDBManager.h"
41 #include "TFile.h"
42 #include "AliHLTTrigger.h"
43
44 /** ROOT macro for the implementation of ROOT specific class methods */
45 ClassImp(AliHLTTriggerBarrelGeomMultiplicity)
46
47 AliHLTTriggerBarrelGeomMultiplicity::AliHLTTriggerBarrelGeomMultiplicity()
48   : AliHLTTrigger()
49   , fSolenoidBz(0)
50   , fMinTracks(1)
51   , fDetectorArray(0)
52   , fTriggerDecisionPars(0)
53   , fTriggerName(0)
54   , fOCDBEntry(0)
55 {
56   // see header file for class documentation
57   // or
58   // refer to README to build package
59   // or
60   // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt
61   fDetectorArray = new TObjArray(1);
62
63   fDetectorArray = new TObjArray(1);
64
65   fDetectorArray = new TObjArray;
66
67 }
68
69 AliHLTTriggerBarrelGeomMultiplicity::~AliHLTTriggerBarrelGeomMultiplicity()
70 {
71   // see header file for class documentation
72 }
73
74 const char* AliHLTTriggerBarrelGeomMultiplicity::GetTriggerName() const 
75 {
76   // see header file for class documentation
77
78   return "BarrelGeomMultiplicityTrigger";
79 }
80
81 AliHLTComponent* AliHLTTriggerBarrelGeomMultiplicity::Spawn()
82 {
83   // see header file for class documentation
84   return new AliHLTTriggerBarrelGeomMultiplicity;
85 }
86
87 int AliHLTTriggerBarrelGeomMultiplicity::Reconfigure(const char *cdbEntry, const char *chainId)
88 {
89   // see header file for class documentation
90
91   // configure from the specified entry or the default
92   const char* entry=cdbEntry;
93
94   if (!entry)
95     {
96       HLTDebug("No CDB path specified");
97       entry = fOCDBEntry; 
98     }
99
100   return GetDetectorGeomsFromCDBObject(entry, chainId);
101
102
103 int AliHLTTriggerBarrelGeomMultiplicity::DoTrigger()
104 {
105   // see header file for class documentation
106   int iResult=0;
107   int numberOfTracks=-1;
108
109   if (!fTriggerDecisionPars) {
110     iResult=-ENODEV;
111   }
112
113   // try the ESD as input
114   const TObject* obj = GetFirstInputObject(kAliHLTAllDataTypes, "AliESDEvent");
115   AliESDEvent* esd = dynamic_cast<AliESDEvent*>(const_cast<TObject*>(obj));
116   TString description;
117
118   if (esd != NULL) 
119     {
120       numberOfTracks=0;
121       esd->GetStdContent();
122       for (Int_t i = 0; i < esd->GetNumberOfTracks(); i++) 
123         {
124           if (CheckCondition(esd->GetTrack(i), esd->GetMagneticField())) numberOfTracks++;
125         }
126     }
127
128   // try the AliHLTExternal track data as input
129   if (iResult>=0 && numberOfTracks<0) 
130     {
131       for (const AliHLTComponentBlockData* pBlock=GetFirstInputBlock(kAliHLTDataTypeTrack);
132            pBlock!=NULL; pBlock=GetNextInputBlock()) 
133         {
134           if (numberOfTracks<0) numberOfTracks=0;
135           vector<AliHLTGlobalBarrelTrack> tracks;
136           if ((iResult=AliHLTGlobalBarrelTrack::ConvertTrackDataArray(reinterpret_cast<const AliHLTTracksData*>(pBlock->fPtr), pBlock->fSize, tracks))>0) 
137             {
138               for (vector<AliHLTGlobalBarrelTrack>::iterator element=tracks.begin();
139                    element!=tracks.end(); element++) 
140                 {
141                   if (CheckCondition(&(*element), fSolenoidBz)) numberOfTracks++;
142                 }
143             } 
144           else if (iResult<0) 
145             {
146               HLTError("can not extract tracks from data block of type %s (specification %08x) of size %d: error %d", 
147                        DataType2Text(pBlock->fDataType).c_str(), pBlock->fSpecification, pBlock->fSize, iResult);
148             }
149         }
150     }
151
152   bool condition=false;
153   description="Geometrical conditions not matched";
154   AliHLTReadoutList readout;
155
156   if (numberOfTracks>=fMinTracks) 
157     {
158       condition=true;
159       description=fTriggerDecisionPars->GetDescription();
160       readout=fTriggerDecisionPars->GetReadoutListParameter();
161       HLTDebug("Geometrical acceptance trigger %s triggered", fTriggerDecisionPars->GetTriggerName().Data());
162     }
163
164   AliHLTTriggerDecision decision(
165                                  condition,
166                                  fTriggerDecisionPars->GetTriggerName().Data(),
167                                  AliHLTTriggerDomain(readout),
168                                  description.Data()
169                                  );
170   TriggerEvent(&decision, kAliHLTDataTypeTObject|kAliHLTDataOriginOut);
171
172   return iResult;
173
174 }
175
176
177 template<class T>
178 bool AliHLTTriggerBarrelGeomMultiplicity::CheckCondition(T* track, float b)
179 {
180
181   bool ret = false;
182
183   // see header file for class documentation
184   if (!track) return false;
185
186   ret = IsInDetectors(track, b);
187
188   return ret;
189
190 }
191
192 template<class T>
193 bool AliHLTTriggerBarrelGeomMultiplicity::IsInDetectors(T* track, float b)
194 {
195   // See header file for class documentation  
196   for(Int_t i = 0; i < fDetectorArray->GetEntries(); i++)
197     {
198       AliHLTTriggerDetectorGeom *det = static_cast<AliHLTTriggerDetectorGeom*>(fDetectorArray->At(i));
199
200       Double_t trackPoint[3];
201       Double_t normVector[3];
202
203       det->GetInitialPoint(trackPoint);
204       det->GetNormVector(normVector);
205
206       bool ret = track->Intersect(trackPoint, normVector, b);
207
208       if(ret)
209         {
210           if(det->IsInDetector(trackPoint)) return true;
211         }
212     }
213   return false;
214 }
215
216 int AliHLTTriggerBarrelGeomMultiplicity::DoInit(int argc, const char** argv)
217 {
218   // see header file for class documentation
219
220   // first configure the default
221   int iResult=0;
222
223   if (iResult>=0 && argc>0)
224     iResult=ConfigureFromArgumentString(argc, argv);
225
226   if (!fTriggerDecisionPars) {
227     HLTError("decision parameter not initialized");
228     iResult=-ENODEV;
229   }
230
231   return iResult;
232 }
233
234 int AliHLTTriggerBarrelGeomMultiplicity::DoDeinit()
235  {
236   // see header file for class documentation
237   return 0;
238 }
239
240 int AliHLTTriggerBarrelGeomMultiplicity::ReadPreprocessorValues(const char* /*modules*/)
241 {
242     // see header file for function documentation
243
244   // TODO 2009-10-10: implementation
245   // for the moment very quick, just reload the magnetic field
246   return ConfigureFromCDBTObjString(kAliHLTCDBSolenoidBz);
247 }
248
249 int AliHLTTriggerBarrelGeomMultiplicity::GetDetectorGeomsFromCDBObject(const char *cdbEntry, const char* chainId)
250 {
251     // see header file for function documentation
252   int nDetectorGeoms=0;
253
254   if(fDetectorArray)
255     {
256       fDetectorArray->Clear();
257     }
258   else
259     {
260       fDetectorArray = new TObjArray();
261     }
262   
263   const char *path = cdbEntry;
264
265   if(!path) path = fOCDBEntry;
266   
267   if(path)
268     {
269       //     const char* chainId=GetChainId();
270       HLTInfo("configure from entry %s, chain id %s", path, (chainId!=NULL && chainId[0]!=0)?chainId:"<none>");
271       AliCDBEntry *pEntry = AliCDBManager::Instance()->Get(path/*,GetRunNo()*/);
272       if (pEntry) 
273         {
274           TObjArray* pArr=dynamic_cast<TObjArray*>(pEntry->GetObject());
275           if (pArr) 
276             {
277
278               for(int i = 0; i < pArr->GetEntries(); i++)
279                 {
280                   if(!strcmp(pArr->At(i)->ClassName(), "AliHLTTriggerDecisionParameters"))
281                     {
282                       fTriggerDecisionPars = dynamic_cast<AliHLTTriggerDecisionParameters*>(pArr->At(i));
283                     }
284                   else if(pArr->At(i)->InheritsFrom("AliHLTTriggerDetectorGeom"))
285                     {
286                       fDetectorArray->AddLast(dynamic_cast<AliHLTTriggerDetectorGeom*>(pArr->At(i)));
287                       nDetectorGeoms++;
288                       HLTDebug("received detector geometry of type %s", pArr->At(i)->ClassName());
289                     }
290                   else
291                     {
292                       HLTWarning("Unknown object of type %s in configuration object", pArr->At(i)->ClassName());
293                     }
294                 }
295             } 
296           else 
297             {
298               HLTError("configuration object \"%s\" has wrong type, required TObjArray", path);
299               nDetectorGeoms=-EINVAL;
300             }
301         }
302       else 
303         {
304           HLTError("can not fetch object \"%s\" from OCDB", path);
305           nDetectorGeoms=-ENOENT;
306         }
307     }
308
309   HLTInfo("received %d detector geometries", nDetectorGeoms);
310
311   return nDetectorGeoms;
312 }
313
314 int AliHLTTriggerBarrelGeomMultiplicity::GetDetectorGeomsFromFile(const char *filename)
315 {
316     // see header file for function documentation
317   int nDetectorGeoms=0;
318
319   if(fDetectorArray)
320     {
321       fDetectorArray->Clear();
322     }
323   else
324     {
325       fDetectorArray = new TObjArray();
326     }
327   
328
329   if (filename) 
330     {
331       TFile *geomfile = TFile::Open(filename, "READ");
332       
333       if(geomfile)
334         {
335           HLTInfo("configure from file \"%s\"", filename);
336           TObjArray* pArr=dynamic_cast<TObjArray*>(geomfile->Get("GeomConf"));
337           if (pArr) 
338             {
339
340               for(int i = 0; i < pArr->GetEntries(); i++)
341                 {
342                   if(!strcmp(pArr->At(i)->ClassName(), "AliHLTTriggerDecisionParameters"))
343                     {
344                       fTriggerDecisionPars = dynamic_cast<AliHLTTriggerDecisionParameters*>(pArr->At(i));
345                     }
346                   else if(pArr->At(i)->InheritsFrom("AliHLTTriggerDetectorGeom"))
347                     {
348                       fDetectorArray->AddLast(dynamic_cast<AliHLTTriggerDetectorGeom*>(pArr->At(i)));
349                       nDetectorGeoms++;
350                       HLTDebug("received detector geometry of type %s", pArr->At(i)->ClassName());
351                     }
352                   else
353                     {
354                       HLTWarning("Unknown object of type %s in configuration object", pArr->At(i)->ClassName());
355                     }
356                 }
357             } 
358           else 
359             {
360               HLTError("configuration object has wrong type, required TObjArray");
361               nDetectorGeoms=-EINVAL;
362             }
363           } 
364       else 
365         {
366           HLTError("could not open file \"%s\"", filename);
367           nDetectorGeoms=-ENOENT;
368         }
369     }
370   else
371     {
372       HLTError("ROOT file name not specified");
373     }
374   HLTInfo("received %d detector geometries", nDetectorGeoms);
375
376   return nDetectorGeoms;
377 }
378
379 int AliHLTTriggerBarrelGeomMultiplicity::ScanConfigurationArgument(int argc, const char** argv)
380 {
381   // See header file for class documentation
382   if (argc<=0) return 0;
383   int i=0;
384   TString argument=argv[i];
385
386   if (argument.CompareTo("-geomfile")==0) 
387     {
388       if (++i>=argc) return -EPROTO;
389     
390       GetDetectorGeomsFromFile(argv[i]);
391     
392       return 2;
393     }    
394
395   if (argument.CompareTo("-triggername")==0) 
396     {
397       if (++i>=argc) return -EPROTO;
398       
399       fTriggerName = new char[128];
400       sprintf(fTriggerName, argv[i]);
401       
402       fOCDBEntry = fTriggerName;
403
404       return 2;
405   }    
406   return 0;
407 }
408