]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/trigger/AliHLTTriggerBarrelGeomMultiplicity.cxx
Modifications to the geometry trigger classes. The trigger has now been tested with...
[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 "TObjArray.h"
40 #include "AliCDBEntry.h"
41 #include "AliCDBManager.h"
42 #include "TFile.h"
43 #include "AliHLTTrigger.h"
44
45 /** ROOT macro for the implementation of ROOT specific class methods */
46 ClassImp(AliHLTTriggerBarrelGeomMultiplicity)
47
48 AliHLTTriggerBarrelGeomMultiplicity::AliHLTTriggerBarrelGeomMultiplicity()
49   : AliHLTTrigger()
50   , fSolenoidBz(0)
51   , fMinTracks(1)
52   , fDetectorArray(0)
53   , fTriggerDecisionPars(0)
54   , fTriggerName(0)
55   , fOCDBEntry(0)
56 {
57   // see header file for class documentation
58   // or
59   // refer to README to build package
60   // or
61   // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt
62   fDetectorArray = new TObjArray(1);
63
64   fDetectorArray = new TObjArray;
65
66 }
67
68 AliHLTTriggerBarrelGeomMultiplicity::~AliHLTTriggerBarrelGeomMultiplicity()
69 {
70   // see header file for class documentation
71 }
72
73 const char* AliHLTTriggerBarrelGeomMultiplicity::GetTriggerName() const 
74 {
75   // see header file for class documentation
76   //  const char* name = fTriggerName;
77   //  if(name) return name;
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 n  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 template<class T>
177 bool AliHLTTriggerBarrelGeomMultiplicity::CheckCondition(T* track, float b)
178 {
179   bool ret = false;
180
181   // see header file for class documentation
182   if (!track) return false;
183
184   ret = IsInDetectors(track, b);
185
186   return ret;
187
188 }
189
190 template<class T>
191 bool AliHLTTriggerBarrelGeomMultiplicity::IsInDetectors(T* track, float b)
192 {
193   // See header file for class documentation  
194   for(Int_t i = 0; i < fDetectorArray->GetEntries(); i++)
195     {
196       AliHLTTriggerDetectorGeom *det = static_cast<AliHLTTriggerDetectorGeom*>(fDetectorArray->At(i));
197
198       Double_t trackPoint[3];
199       Double_t normVector[3];
200
201       det->GetInitialPoint(trackPoint);
202       det->GetNormVector(normVector);
203
204       bool ret = track->Intersect(trackPoint, normVector, b);
205
206       if(ret)
207         {
208           if(det->IsInDetector(trackPoint)) return true;
209         }
210     }
211   return false;
212 }
213
214 int AliHLTTriggerBarrelGeomMultiplicity::DoInit(int argc, const char** argv)
215 {
216   // see header file for class documentation
217
218   // first configure the default
219   int iResult=0;
220
221   if (iResult>=0 && argc>0)
222     iResult=ConfigureFromArgumentString(argc, argv);
223
224   if (!fTriggerDecisionPars) {
225     HLTError("decision parameter not initialized");
226     iResult=-ENODEV;
227   }
228
229   return iResult;
230 }
231
232 int AliHLTTriggerBarrelGeomMultiplicity::DoDeinit()
233  {
234   // see header file for class documentation
235   return 0;
236 }
237
238 int AliHLTTriggerBarrelGeomMultiplicity::ReadPreprocessorValues(const char* /*modules*/)
239 {
240     // see header file for function documentation
241   int nDetectorGeoms=0;
242
243   // TODO 2009-10-10: implementation
244   // for the moment very quick, just reload the magnetic field
245   return ConfigureFromCDBTObjString(kAliHLTCDBSolenoidBz);
246 }
247
248 int AliHLTTriggerBarrelGeomMultiplicity::GetDetectorGeomsFromCDBObject(const char *cdbEntry, const char* chainId)
249 {
250     // see header file for function documentation
251   int nDetectorGeoms=0;
252
253   if(fDetectorArray)
254     {
255       fDetectorArray->Clear();
256     }
257   else
258     {
259       fDetectorArray = new TObjArray();
260     }
261   
262   const char *path = cdbEntry;
263
264   if(!path) path = fOCDBEntry;
265   
266   if(path)
267     {
268       //     const char* chainId=GetChainId();
269       HLTInfo("configure from entry %s, chain id %s", path, (chainId!=NULL && chainId[0]!=0)?chainId:"<none>");
270       AliCDBEntry *pEntry = AliCDBManager::Instance()->Get(path/*,GetRunNo()*/);
271       if (pEntry) 
272         {
273           TObjArray* pArr=dynamic_cast<TObjArray*>(pEntry->GetObject());
274           if (pArr) 
275             {
276
277               for(int i = 0; i < pArr->GetEntries(); i++)
278                 {
279                   if(!strcmp(pArr->At(i)->ClassName(), "AliHLTTriggerDecisionParameters"))
280                     {
281                       fTriggerDecisionPars = dynamic_cast<AliHLTTriggerDecisionParameters*>(pArr->At(i));
282                     }
283                   else if(pArr->At(i)->InheritsFrom("AliHLTTriggerDetectorGeom"))
284                     {
285                       fDetectorArray->AddLast(dynamic_cast<AliHLTTriggerDetectorGeom*>(pArr->At(i)));
286                       nDetectorGeoms++;
287                       HLTDebug("received detector geometry of type %s", pArr->At(i)->ClassName());
288                     }
289                   else
290                     {
291                       HLTWarning("Unknown object of type %s in configuration object", pArr->At(i)->ClassName());
292                     }
293                 }
294             } 
295           else 
296             {
297               HLTError("configuration object \"%s\" has wrong type, required TObjArray", path);
298               nDetectorGeoms=-EINVAL;
299             }
300         }
301       else 
302         {
303           HLTError("can not fetch object \"%s\" from OCDB", path);
304           nDetectorGeoms=-ENOENT;
305         }
306     }
307
308   HLTInfo("received %d detector geometries", nDetectorGeoms);
309
310   return nDetectorGeoms;
311 }
312
313 int AliHLTTriggerBarrelGeomMultiplicity::GetDetectorGeomsFromFile(const char *filename)
314 {
315     // see header file for function documentation
316   int nDetectorGeoms=0;
317
318   if(fDetectorArray)
319     {
320       fDetectorArray->Clear();
321     }
322   else
323     {
324       fDetectorArray = new TObjArray();
325     }
326   
327
328   if (filename) 
329     {
330       TFile *geomfile = TFile::Open(filename, "READ");
331       
332       if(geomfile)
333         {
334           HLTInfo("configure from file \"%s\"", filename);
335           TObjArray* pArr=dynamic_cast<TObjArray*>(geomfile->Get("GeomConf"));
336           if (pArr) 
337             {
338
339               for(int i = 0; i < pArr->GetEntries(); i++)
340                 {
341                   if(!strcmp(pArr->At(i)->ClassName(), "AliHLTTriggerDecisionParameters"))
342                     {
343                       fTriggerDecisionPars = dynamic_cast<AliHLTTriggerDecisionParameters*>(pArr->At(i));
344                     }
345                   else if(pArr->At(i)->InheritsFrom("AliHLTTriggerDetectorGeom"))
346                     {
347                       fDetectorArray->AddLast(dynamic_cast<AliHLTTriggerDetectorGeom*>(pArr->At(i)));
348                       nDetectorGeoms++;
349                       HLTDebug("received detector geometry of type %s", pArr->At(i)->ClassName());
350                     }
351                   else
352                     {
353                       HLTWarning("Unknown object of type %s in configuration object", pArr->At(i)->ClassName());
354                     }
355                 }
356             } 
357           else 
358             {
359               HLTError("configuration object has wrong type, required TObjArray");
360               nDetectorGeoms=-EINVAL;
361             }
362           } 
363       else 
364         {
365           HLTError("could not open file \"%s\"", filename);
366           nDetectorGeoms=-ENOENT;
367         }
368     }
369   else
370     {
371       HLTError("ROOT file name not specified");
372     }
373   HLTInfo("received %d detector geometries", nDetectorGeoms);
374
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