]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/trigger/AliHLTTriggerBarrelGeomMultiplicity.cxx
- addign classes to the build system
[u/mrichter/AliRoot.git] / HLT / trigger / AliHLTTriggerBarrelGeomMultiplicity.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: Matthias Richter <Matthias.Richter@ift.uib.no>        *
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 /// @file   AliHLTTriggerBarrelGeomMultiplicity.cxx
18 /// @author Oystein Djuvsland
19 /// @date   2009-10-08
20 /// @brief  HLT trigger component for charged particle multiplicity in
21 ///         the central barrel.
22
23 // see header file for class documentation
24 // or
25 // refer to README to build package
26 // or
27 // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt
28
29 #include "AliHLTTriggerBarrelGeomMultiplicity.h"
30 #include "AliESDEvent.h"
31 #include "AliHLTTriggerDecision.h"
32 #include "AliHLTDomainEntry.h"
33 #include "AliHLTGlobalBarrelTrack.h"
34 #include "TObjArray.h"
35 #include "TObjString.h"
36
37 /** ROOT macro for the implementation of ROOT specific class methods */
38 ClassImp(AliHLTTriggerBarrelGeomMultiplicity)
39
40 AliHLTTriggerBarrelGeomMultiplicity::AliHLTTriggerBarrelGeomMultiplicity()
41   : AliHLTTrigger()
42   , fMinTracks(1)
43   , fSolenoidBz(0.0)
44 {
45   // see header file for class documentation
46   // or
47   // refer to README to build package
48   // or
49   // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt
50
51 }
52
53 const char* AliHLTTriggerBarrelGeomMultiplicity::fgkOCDBEntry="HLT/ConfigHLT/BarrelGeomMultiplicityTrigger";
54
55 AliHLTTriggerBarrelGeomMultiplicity::~AliHLTTriggerBarrelGeomMultiplicity()
56 {
57   // see header file for class documentation
58 }
59
60 const char* AliHLTTriggerBarrelGeomMultiplicity::GetTriggerName() const
61 {
62   // see header file for class documentation
63   return "BarrelGeomMultiplicityTrigger";
64 }
65
66 AliHLTComponent* AliHLTTriggerBarrelGeomMultiplicity::Spawn()
67 {
68   // see header file for class documentation
69   return new AliHLTTriggerBarrelGeomMultiplicity;
70 }
71
72 int AliHLTTriggerBarrelGeomMultiplicity::DoTrigger()
73 {
74   // see header file for class documentation
75   int iResult=0;
76   int numberOfTracks=-1;
77
78   // try the ESD as input
79   const TObject* obj = GetFirstInputObject(kAliHLTAllDataTypes, "AliESDEvent");
80   AliESDEvent* esd = dynamic_cast<AliESDEvent*>(const_cast<TObject*>(obj));
81   TString description;
82   TString ptcut,tdca,ldca,dcaref,op1st,op2nd;
83   if (esd != NULL) {
84     numberOfTracks=0;
85     esd->GetStdContent();
86     
87     for (Int_t i = 0; i < esd->GetNumberOfTracks(); i++) {
88       if (CheckCondition(esd->GetTrack(i), esd->GetMagneticField())) numberOfTracks++;
89     }
90   }
91
92   // try the AliHLTExternal track data as input
93   if (iResult>=0 && numberOfTracks<0) {
94     for (const AliHLTComponentBlockData* pBlock=GetFirstInputBlock(kAliHLTDataTypeTrack);
95          pBlock!=NULL; pBlock=GetNextInputBlock()) {
96       if (numberOfTracks<0) numberOfTracks=0;
97       vector<AliHLTGlobalBarrelTrack> tracks;
98       if ((iResult=AliHLTGlobalBarrelTrack::ConvertTrackDataArray(reinterpret_cast<const AliHLTTracksData*>(pBlock->fPtr), pBlock->fSize, tracks))>0) {
99         for (vector<AliHLTGlobalBarrelTrack>::iterator element=tracks.begin();
100              element!=tracks.end(); element++) {
101           if (CheckCondition(&(*element), fSolenoidBz)) numberOfTracks++;
102         }
103       } else if (iResult<0) {
104         HLTError("can not extract tracks from data block of type %s (specification %08x) of size %d: error %d", 
105                  DataType2Text(pBlock->fDataType).c_str(), pBlock->fSpecification, pBlock->fSize, iResult);
106       }
107     }
108   }
109   if (numberOfTracks>=fMinTracks) {
110
111     ApplyTrigger();
112
113     description.Form("Event contains %d track(s) with ", numberOfTracks);
114     SetDescription(description.Data());
115     // Enable the central detectors for readout.
116     GetReadoutList().Enable(
117                             AliHLTReadoutList::kITSSPD |
118                             AliHLTReadoutList::kITSSDD |
119                             AliHLTReadoutList::kITSSSD |
120                             AliHLTReadoutList::kTPC |
121                             AliHLTReadoutList::kTRD |
122                             AliHLTReadoutList::kTOF |
123                             AliHLTReadoutList::kHMPID |
124                             AliHLTReadoutList::kPHOS
125                             );
126       // Add the available HLT information for readout too.
127       GetTriggerDomain().Add("CLUSTERS", "TPC ");
128       TriggerEvent(true);
129       return 0;
130     }
131     description.Form("No tracks matching the tresholds found in the central barrel (min tracks %d, ", fMinTracks);
132     description+=ptcut;
133     description+=op1st;
134     description+=ldca;
135     description+=op2nd;
136     description+=tdca;
137     description+=dcaref;
138     description+=")";
139   } else {
140     description.Form("No input blocks found");
141   }
142   SetDescription(description.Data());
143   TriggerEvent(false);
144   return iResult;
145 }
146
147 template<class T>
148 bool AliHLTTriggerBarrelGeomMultiplicity::CheckCondition(T* track, float b)
149 {
150   // see header file for class documentation
151   if (!track) return false;
152
153   ret = IsInDetectors(track, b);
154
155   return ret;
156
157 }
158
159 template<class T>
160 bool AliHLTTriggerBarrelGeomMultiplicity::IsInDetector(T* track, b)
161 {
162
163   for(Int_t i = 0; i < fDetectorList->GetEntries(); i++)
164     {
165       AliHLTTriggerDetectorGeom *det = static_cast<AliHLTTriggerDetectorGeom*>(fDetectorList.At(i))
166       Double_t trackPoint[3];
167
168       det->GetInitialPoint(trackPoint);
169
170       bool ret = track->Intersect(trackPoint, det->NormVector(), b);
171
172       if(track->Eta() >= det->EtaMin() && 
173          track->Eta() <= det->EtaMax() &&
174          track->Phi() >= det->PhiMin() &&
175          track->Phi() <= det->PhiMax())
176         {
177           return true;
178         }
179     }
180   return false;
181 }
182
183 int AliHLTTriggerBarrelGeomMultiplicity::DoInit(int argc, const char** argv)
184 {
185   // see header file for class documentation
186
187   // first configure the default
188   int iResult=0;
189   iResult=ConfigureFromCDBTObjString(kAliHLTCDBSolenoidBz);
190   if (iResult>=0) iResult=ConfigureFromCDBTObjString(fgkOCDBEntry);
191
192   // configure from the command line parameters if specified
193   if (iResult>=0 && argc>0)
194     iResult=ConfigureFromArgumentString(argc, argv);
195   return iResult;
196 }
197
198 int AliHLTTriggerBarrelGeomMultiplicity::DoDeinit()
199 {
200   // see header file for class documentation
201   return 0;
202 }
203
204 int AliHLTTriggerBarrelGeomMultiplicity::Reconfigure(const char* cdbEntry, const char* /*chainId*/)
205 {
206   // see header file for class documentation
207
208   // configure from the specified antry or the default one
209   const char* entry=cdbEntry;
210   if (!entry || entry[0]==0) {
211     ConfigureFromCDBTObjString(kAliHLTCDBSolenoidBz);
212     entry=fgkOCDBEntry;
213   }
214
215   return ConfigureFromCDBTObjString(entry);
216 }
217
218 int AliHLTTriggerBarrelGeomMultiplicity::ReadPreprocessorValues(const char* /*modules*/)
219 {
220   // see header file for class documentation
221
222   // TODO 2009-09-10: implementation
223   // for the moment very quick, just reload the magnetic field
224   return ConfigureFromCDBTObjString(kAliHLTCDBSolenoidBz);
225 }
226
227 int AliHLTTriggerBarrelGeomMultiplicity::ScanConfigurationArgument(int argc, const char** argv)
228 {
229   // see header file for class documentation
230   if (argc<=0) return 0;
231   int i=0;
232   TString argument=argv[i];
233
234   // -maxpt
235   if (argument.CompareTo("-maxpt")==0) {
236     if (++i>=argc) return -EPROTO;
237     argument=argv[i];
238     fPtMax=argument.Atof();
239     return 2;
240   }    
241
242   // -minpt
243   if (argument.CompareTo("-minpt")==0) {
244     if (++i>=argc) return -EPROTO;
245     argument=argv[i];
246     fPtMin=argument.Atof();
247     return 2;
248   }    
249
250   // -mintracks
251   if (argument.CompareTo("-mintracks")==0) {
252     if (++i>=argc) return -EPROTO;
253     argument=argv[i];
254     fMinTracks=argument.Atoi();
255     return 2;
256   }    
257
258   // -dca-reference
259   // reference point for the transverse and longitudinal dca cut
260   if (argument.CompareTo("-dca-reference")==0) {
261     if (++i>=argc) return -EPROTO;
262     argument=argv[i];
263     // scan x,y,z
264     TObjArray* pTokens=argument.Tokenize("'");
265     if (pTokens) {
266       for (int c=0; c<pTokens->GetEntriesFast() && c<fgkDCAReferenceSize; c++) {
267         argument=((TObjString*)pTokens->At(c))->GetString();
268         fDCAReference[i]=argument.Atof();
269       }
270       delete pTokens;
271     }
272     return 2;
273   }
274
275   // -min-ldca
276   // minimum longitudinal dca to reference point
277   if (argument.CompareTo("-min-ldca")==0) {
278     if (++i>=argc) return -EPROTO;
279     argument=argv[i];
280     fMinLDca=argument.Atof();
281     return 2;
282   }
283   
284   // -max-ldca
285   // maximum longitudinal dca to reference point
286   if (argument.CompareTo("-max-ldca")==0) {
287     if (++i>=argc) return -EPROTO;
288     argument=argv[i];
289     fMaxLDca=argument.Atof();
290     return 2;
291   }
292
293   // -min-tdca
294   // minimum transverse dca to reference point
295   if (argument.CompareTo("-min-tdca")==0) {
296     if (++i>=argc) return -EPROTO;
297     argument=argv[i];
298     fMinTDca=argument.Atof();
299     return 2;
300   }
301   
302   // -max-tdca
303   // maximum transverse dca to reference point
304   if (argument.CompareTo("-max-tdca")==0) {
305     if (++i>=argc) return -EPROTO;
306     argument=argv[i];
307     fMaxTDca=argument.Atof();
308     return 2;
309   }
310
311   // -solenoidBz
312   if (argument.CompareTo("-solenoidBz")==0) {
313     if (++i>=argc) return -EPROTO;
314     argument=argv[i];
315     fSolenoidBz=argument.Atof();
316     return 2;
317   }
318
319   // unknown argument
320   return -EINVAL;
321 }