]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/TPCLib/AliHLTTPCdEdxMonitoringComponent.cxx
334197367ccfa40c3e313e0d3343a53186ae198e
[u/mrichter/AliRoot.git] / HLT / TPCLib / AliHLTTPCdEdxMonitoringComponent.cxx
1 // $Id$
2
3 //**************************************************************************
4 //* This file is property of and copyright by the ALICE HLT Project        * 
5 //* ALICE Experiment at CERN, All rights reserved.                         *
6 //*                                                                        *
7 //* Primary Authors: Per-Ivar Lønne <perivarlonne@gmail.com>               *
8 //*                  for The ALICE HLT Project.                            *
9 //*                                                                        *
10 //* Permission to use, copy, modify and distribute this software and its   *
11 //* documentation strictly for non-commercial purposes is hereby granted   *
12 //* without fee, provided that the above copyright notice appears in all   *
13 //* copies and that both the copyright notice and this permission notice   *
14 //* appear in the supporting documentation. The authors make no claims     *
15 //* about the suitability of this software for any purpose. It is          *
16 //* provided "as is" without express or implied warranty.                  *
17 //**************************************************************************/
18
19 /// @file   AliHLTTPCdEdxMonitoringComponent.cxx
20 /// @author Per-Ivar Lønne, Jochen Thaeder, Matthias Richter, Alexander Kalweit
21 /// @date   21.08.2011
22 /// @brief  Component for reading ESD from chain and produce a dEdx monitoring plot
23 ///
24
25 #if __GNUC__>= 3
26 using namespace std;
27 #endif
28
29
30 #include "TSystem.h"
31 #include "AliESDtrackCuts.h"
32 #include <AliHLTDAQ.h>
33 #include "AliESDEvent.h"
34 #include "AliESDtrack.h"
35 #include "AliLog.h"
36 #include "TMap.h"
37 #include "TObjString.h"
38 #include "AliHLTDataTypes.h"
39 #include "AliHLTComponentBenchmark.h"
40 #include "TROOT.h"
41 #include "TStyle.h"
42 #include "TColor.h"
43 #include "TH2F.h"
44
45 #include "AliHLTTPCdEdxMonitoringComponent.h"
46
47 /** ROOT macro for the implementation of ROOT specific class methods */
48 ClassImp( AliHLTTPCdEdxMonitoringComponent )
49   
50 AliHLTTPCdEdxMonitoringComponent::AliHLTTPCdEdxMonitoringComponent() 
51   :
52   AliHLTProcessor(),
53   fESDTrackCuts(NULL),
54   fHist(),
55   fxbins(),
56   fxmin(),
57   fxmax(),
58   fybins(),
59   fymin(),
60   fymax()
61 {
62   //Constructor
63 }
64
65 AliHLTTPCdEdxMonitoringComponent::~AliHLTTPCdEdxMonitoringComponent()
66 {
67   //Destructor
68 }
69
70 // ######################################################################### //
71 const char* AliHLTTPCdEdxMonitoringComponent::GetComponentID()
72 {
73   // get component id
74   return "TPCdEdxMonitoring";
75 }
76
77 // ######################################################################### //
78 void AliHLTTPCdEdxMonitoringComponent::GetInputDataTypes( vector<AliHLTComponentDataType>& list)
79 {
80   // get list of input data types
81   list.push_back(kAliHLTDataTypeESDObject|kAliHLTDataOriginAny);
82 }
83
84 // ######################################################################### //
85 AliHLTComponentDataType AliHLTTPCdEdxMonitoringComponent::GetOutputDataType()
86 {
87   // get the output data size
88   return kAliHLTDataTypeHistogram|kAliHLTDataOriginHLT;
89 }
90
91 // ######################################################################### //
92 void AliHLTTPCdEdxMonitoringComponent::GetOutputDataSize( unsigned long& constBase, Double_t& inputMultiplier )
93 {
94   // get output size estimator
95   constBase = 4096;
96   inputMultiplier = 0;
97 }
98
99 // ######################################################################### //
100 void AliHLTTPCdEdxMonitoringComponent::GetOCDBObjectDescription( TMap* const targetMap)
101 {
102   // Get a list of OCDB object description.
103   // The list of objects is provided in a TMap
104   // - key: complete OCDB path, e.g. GRP/GRP/Data
105   // - value: short description why the object is needed
106   // Key and value objects created inside this class go into ownership of
107   // target TMap.
108   
109     if (!targetMap) return;
110     targetMap->Add(new TObjString("HLT/ConfigTPC/TPCdEdxMonitoring"),
111     new TObjString("configuration object"));
112     targetMap->Add(new TObjString("GRP/GRP/Data"),
113     new TObjString("GRP object"));
114   
115 }
116
117 // ######################################################################### //
118 AliHLTComponent* AliHLTTPCdEdxMonitoringComponent::Spawn()
119 {
120   // Spawn function, return new class instance
121   return new AliHLTTPCdEdxMonitoringComponent;
122 }
123
124 // ######################################################################### //
125 Int_t AliHLTTPCdEdxMonitoringComponent::DoInit( Int_t argc, const char** argv )
126 {
127   // init the component
128   Int_t iResult=0;
129   fESDTrackCuts = new AliESDtrackCuts("AliESDtrackCuts","HLT");
130
131   fxbins=Int_t(300);
132   fxmin=Double_t(-2);
133   fxmax=Double_t(2);
134   fybins=Int_t(500);
135   fymin=Double_t(0);
136   fymax=Double_t(500);
137   
138 if (!fESDTrackCuts)
139     {
140       iResult=-ENOMEM;
141     }
142
143   if (iResult<0)
144     {
145       if (fESDTrackCuts)
146         delete fESDTrackCuts;
147       fESDTrackCuts = NULL;
148     }      
149        
150   if (iResult>=0) 
151     {
152       SetDefaultConfiguration();
153       TString cdbPath="HLT/ConfigTPC/"; 
154       cdbPath+=GetComponentID();
155       iResult=ConfigureFromCDBTObjString(cdbPath);
156     
157   
158       if (iResult>=0) 
159         {
160           iResult=ConfigureFromArgumentString(argc, argv);      
161         }
162     }
163
164   if (iResult>=0) {
165     HLTInfo("ESD track cuts : %s",fESDTrackCuts->GetTitle() );
166   }
167   fHist = new TH2F("hHLT", "HLT", fxbins, fxmin, fxmax, fybins, fymin, fymax);
168   fHist->GetXaxis()->SetTitle("momentum/charge #frac{p}{z}  (GeV/c)");
169   fHist->GetYaxis()->SetTitle("dE/dx in TPC (a.u.)");
170   Plotstyle();
171   return iResult;
172 }
173   
174
175 // ######################################################################### //
176 Int_t AliHLTTPCdEdxMonitoringComponent::DoDeinit()
177 {
178   // component cleanup, delete all instances of helper classes here
179   if (fESDTrackCuts)
180     delete fESDTrackCuts;
181   fESDTrackCuts = NULL;
182   if (fHist)
183     delete fHist;
184   fHist=NULL;
185   return 0;
186 }
187
188
189 // ######################################################################### //
190 Int_t AliHLTTPCdEdxMonitoringComponent::DoEvent(const AliHLTComponentEventData& /*evtData*/,
191                                               AliHLTComponentTriggerData& /*trigData*/)
192 {
193   // event processing function
194   float sig;
195   Double_t ptot;
196
197   // check if this is a data event, there are a couple of special events
198   // which should be ignored for normal processing
199   if (!IsDataEvent()) return 0;
200
201   const TObject* obj = GetFirstInputObject(kAliHLTAllDataTypes, "AliESDEvent");
202
203   // input objects are not supposed to be changed by the component, so they
204   // are defined const. However, the implementation of AliESDEvent does not
205   // support this and we need the const_cast
206   AliESDEvent* esd = dynamic_cast<AliESDEvent*>(const_cast<TObject*>(obj));
207   if (esd != NULL) {
208     AliInfoClass(Form("==================== event %3d ================================", GetEventCount()));
209     esd->GetStdContent();
210     for (Int_t i = 0; i < esd->GetNumberOfTracks(); i++) {
211       AliESDtrack* track = esd->GetTrack(i);
212       AliInfoClass(Form("-------------------- track %3d --------------------------------", i));
213       track->Print("");
214       sig=track->GetTPCsignal();
215       if(!fESDTrackCuts->AcceptTrack(track)) continue;
216       if (!track->GetInnerParam()) continue;
217       ptot = track->GetInnerParam()->GetP()*track->GetSign();
218       fHist->Fill(ptot, sig);
219     }
220   }
221   // publish the histogram
222   PushBack(dynamic_cast<TObject*>(fHist), kAliHLTDataTypeHistogram);
223
224   return 0;
225 }
226
227 // ######################################################################### //
228 Int_t AliHLTTPCdEdxMonitoringComponent::ScanConfigurationArgument(Int_t argc, const char** argv)
229 {
230   // Scan configuration arguments
231   // Return the number of processed arguments
232   //        -EPROTO if argument format error (e.g. number expected but not found)
233   //
234   // The AliHLTComponent base class implements a parsing loop for argument strings and
235   // arrays of strings which is invoked by ConfigureFromArgumentString/ConfigureFromCDBTObjString
236   // The component needs to implement ScanConfigurationArgument in order to decode the arguments.
237   
238   if (argc<=0) return 0;
239   Int_t ii =0;
240   TString argument=argv[ii];
241   
242   if (argument.IsNull()) return 0;
243   
244   if( !fESDTrackCuts){
245     HLTError("No ESD track cuts availible");
246     return -ENOMEM;
247   }
248   
249
250   //**********************************//
251   //        Histogram Binning         //
252   //**********************************//
253
254
255   // -xbins
256   if (argument.CompareTo("-xbins")==0) 
257     {
258       if (++ii>=argc) return -EPROTO;
259       argument=argv[ii];
260
261       fxbins = argument.Atoi();
262
263       return 2;
264     }
265
266   // -xmin
267   if (argument.CompareTo("-xmin")==0) 
268     {
269       if (++ii>=argc) return -EPROTO;
270       argument=argv[ii];
271
272       fxmin = argument.Atoi();
273
274       return 2;
275     }
276
277   if (argument.CompareTo("-xmax")==0) 
278     {
279       if (++ii>=argc) return -EPROTO;
280       argument=argv[ii];
281
282       fxmax = argument.Atoi();
283
284       return 2;
285     }
286
287   // -xbins
288   if (argument.CompareTo("-ybins")==0) 
289     {
290       if (++ii>=argc) return -EPROTO;
291       argument=argv[ii];
292
293       fybins = argument.Atoi();
294
295       return 2;
296     }
297
298   // -xmin
299   if (argument.CompareTo("-ymin")==0) 
300     {
301       if (++ii>=argc) return -EPROTO;
302       argument=argv[ii];
303
304       fymin = argument.Atoi();
305
306       return 2;
307     }
308
309   if (argument.CompareTo("-ymax")==0) 
310     {
311       if (++ii>=argc) return -EPROTO;
312       argument=argv[ii];
313
314       fymax = argument.Atoi();
315
316       return 2;
317     }
318
319
320   //**********************************//
321   //           Track Cuts             //
322   //**********************************//
323
324  // -maxpt
325   if (argument.CompareTo("-maxpt")==0) {
326     if (++ii>=argc) return -EPROTO;
327     argument=argv[ii];
328
329     Float_t minPt, maxPt;
330     fESDTrackCuts->GetPtRange(minPt,maxPt);
331     maxPt = argument.Atof(); 
332     fESDTrackCuts->SetPtRange(minPt,maxPt);
333
334     TString title = fESDTrackCuts->GetTitle();
335     if (!title.CompareTo("No track cuts")) title = "";
336     else title += " && ";
337     title += Form("p_t < %f", maxPt);
338     fESDTrackCuts->SetTitle(title);
339     return 2;
340   }    
341
342   // -minpt
343   if (argument.CompareTo("-minpt")==0) {
344     if (++ii>=argc) return -EPROTO;
345     argument=argv[ii];
346
347     Float_t minPt, maxPt;
348     fESDTrackCuts->GetPtRange(minPt,maxPt);
349     minPt = argument.Atof(); 
350     fESDTrackCuts->SetPtRange(minPt,maxPt);
351
352     TString title = fESDTrackCuts->GetTitle();
353     if (!title.CompareTo("No track cuts")) title = "";
354     else title += " && ";
355     title += Form("p_t > %f", minPt);
356     fESDTrackCuts->SetTitle(title);
357     return 2;
358   }    
359
360   // -min-ldca
361   // minimum longitudinal dca to vertex
362   if (argument.CompareTo("-min-ldca")==0) {
363     if (++ii>=argc) return -EPROTO;
364     argument=argv[ii];
365
366     fESDTrackCuts->SetMinDCAToVertexZ(argument.Atof());
367     TString title = fESDTrackCuts->GetTitle();
368     if (!title.CompareTo("No track cuts")) title = "";
369     else title += " && ";
370     title += Form("DCAz > %f", argument.Atof());
371     fESDTrackCuts->SetTitle(title);
372     return 2;
373   }
374   
375   // -max-ldca
376   // maximum longitudinal dca to vertex
377   if (argument.CompareTo("-max-ldca")==0) {
378     if (++ii>=argc) return -EPROTO;
379     argument=argv[ii];
380
381     fESDTrackCuts->SetMaxDCAToVertexZ(argument.Atof());
382     TString title = fESDTrackCuts->GetTitle();
383     if (!title.CompareTo("No track cuts")) title = "";
384     else title += " && ";
385     title += Form("DCAz < %f", argument.Atof());
386     fESDTrackCuts->SetTitle(title);
387     return 2;
388   }
389
390   // -min-tdca
391   // minimum transverse dca to vertex
392   if (argument.CompareTo("-min-tdca")==0) {
393     if (++ii>=argc) return -EPROTO;
394     argument=argv[ii];
395
396     fESDTrackCuts->SetMinDCAToVertexXY(argument.Atof());
397     TString title = fESDTrackCuts->GetTitle();
398     if (!title.CompareTo("No track cuts")) title = "";
399     else title += " && ";
400     title += Form("DCAr > %f", argument.Atof());
401     fESDTrackCuts->SetTitle(title);
402     return 2;
403   }
404   
405   // -max-tdca
406   // maximum transverse dca to vertex
407   if (argument.CompareTo("-max-tdca")==0) {
408     if (++ii>=argc) return -EPROTO;
409     argument=argv[ii];
410
411     fESDTrackCuts->SetMaxDCAToVertexXY(argument.Atof());
412     TString title = fESDTrackCuts->GetTitle();
413     if (!title.CompareTo("No track cuts")) title = "";
414     else title += " && ";
415     title += Form("DCAr < %f", argument.Atof());
416     fESDTrackCuts->SetTitle(title);
417     return 2;
418   }
419
420   // -etarange
421   // +/- eta 
422   if (argument.CompareTo("-etarange")==0) {
423     if (++ii>=argc) return -EPROTO;
424     argument=argv[ii];
425     Float_t eta = argument.Atof();
426
427     fESDTrackCuts->SetEtaRange(-eta,eta);     
428     TString title = fESDTrackCuts->GetTitle();
429     if (!title.CompareTo("No track cuts")) title = "";
430     else title += " && ";
431     title += Form("Eta[%f,%f]", argument.Atof(),argument.Atof());
432     fESDTrackCuts->SetTitle(title);
433     return 2;
434   }
435
436   // -minNClsTPC
437   // minimum clusters in TPC
438   if (argument.CompareTo("-minNClsTPC")==0) {
439     if (++ii>=argc) return -EPROTO;
440     argument=argv[ii];
441
442     Int_t ncls;
443     ncls = Int_t(argument.Atof()); 
444     fESDTrackCuts->SetMinNClustersTPC(ncls);  
445
446     TString title = fESDTrackCuts->GetTitle();
447     if (!title.CompareTo("No track cuts")) title = "";
448     else title += " && ";
449     title += Form("minNClsTPC < %i", ncls);
450     fESDTrackCuts->SetTitle(title);
451     return 2;
452   }    
453
454
455
456   // unknown argument
457   return -EINVAL;
458   
459   return 0;
460   
461 }
462
463 // ######################################################################### //
464 Int_t AliHLTTPCdEdxMonitoringComponent::Reconfigure(const char* cdbEntry, const char* chainId)
465 {
466   // reconfigure the component from the specified CDB entry, or default CDB entry
467   HLTInfo("reconfigure '%s' from entry %s", chainId, cdbEntry);
468
469   Int_t iResult=0;
470   TString cdbPath;
471   if (cdbEntry)
472     {
473       cdbPath=cdbEntry;
474     } else {
475       cdbPath="HLT/ConfigTPC/";
476       cdbPath+=GetComponentID();
477     }
478
479   iResult=ConfigureFromCDBTObjString(cdbPath); //// Or use return 0, and skip this line?
480
481   return iResult;
482 }
483
484 // ######################################################################### //
485 Int_t AliHLTTPCdEdxMonitoringComponent::ReadPreprocessorValues(const char* modules)
486 {
487   // read the preprocessor values for the detectors in the modules list
488   Int_t iResult=0;
489   TString detectors(modules!=NULL?modules:"");
490   HLTInfo("read preprocessor values for detector(s): %s", detectors.IsNull()?"none":detectors.Data());
491   return iResult; //Done differently in AliHLTMultiplicityCorrelationsComponent...
492 }
493
494
495
496 // ######################################################################### //
497 void AliHLTTPCdEdxMonitoringComponent::SetDefaultConfiguration() 
498 {
499   if (fESDTrackCuts)
500     {
501       //fESDTrackCuts = AliESDtrackCuts::GetStandardITSTPCTrackCuts2010(kTRUE);
502       fESDTrackCuts->SetEtaRange(-0.8,+0.8);
503       fESDTrackCuts->SetPtRange(0.15,1e10);
504       fESDTrackCuts->SetMaxCovDiagonalElements(2, 2, 0.5, 0.5, 2);  // BEWARE STANDARD VALUES ARE: 2, 2, 0.5, 0.5, 2
505       fESDTrackCuts->SetMaxNsigmaToVertex(3);
506       fESDTrackCuts->SetRequireSigmaToVertex(kTRUE);
507       fESDTrackCuts->SetAcceptKinkDaughters(kFALSE);
508       fESDTrackCuts->SetMinNClustersTPC(70);
509       fESDTrackCuts->SetMaxChi2PerClusterTPC(4);
510       fESDTrackCuts->SetMaxDCAToVertexXY(3);
511       fESDTrackCuts->SetMaxDCAToVertexZ(3);
512       fESDTrackCuts->SetRequireTPCRefit(kTRUE);
513       //fESDTrackCuts->SetRequireITSRefit(kTRUE); //Kills HLT simulated reconstructions?
514       fESDTrackCuts->SetClusterRequirementITS(AliESDtrackCuts::kSPD, AliESDtrackCuts::kAny); //TEMPORARY <-> REMOVE
515       fESDTrackCuts->SetMinNClustersITS(3);
516       
517     }
518  fxbins=300;
519  fxmin=-2;
520  fxmax=2;
521  fybins=500;  
522  fymin=0;
523  fymax=500;
524   return;
525 }
526
527
528
529 // ######################################################################### //
530 void AliHLTTPCdEdxMonitoringComponent::Plotstyle()
531 {
532   gROOT->SetStyle("Plain");
533   
534   gStyle->SetCanvasBorderMode(0);
535   gStyle->SetCanvasColor(10);
536   gStyle->SetCanvasDefH(550);
537   gStyle->SetCanvasDefW(575);
538   gStyle->SetPadBorderMode(0);
539   gStyle->SetPadColor(10);
540   gStyle->SetPadTickX(1);
541   gStyle->SetPadTickY(1);
542   gStyle->SetStatColor(10);
543   gStyle->SetFrameFillColor(10);
544   gStyle->SetPalette(1,0);
545   
546   //
547   
548   //gStyle->SetStatX(0.7);
549   //gStyle->SetStatW(0.2);
550   //gStyle->SetLabelOffset(1.2);
551   //gStyle->SetLabelFont(72);
552   //gStyle->SetLabelSize(0.6);
553   //gStyle->SetTitleOffset(1.2);
554   gStyle->SetTitleFontSize(0.04);
555   
556   
557   gStyle->SetOptStat(10);
558   gStyle->SetLineWidth(2);
559   gStyle->SetMarkerSize(1.0);
560   gStyle->SetTextSize(0.04);
561   gStyle->SetTitleSize(0.04,"xyz");
562   gStyle->SetLabelSize(0.04,"xyz");
563   gStyle->SetLabelOffset(0.02,"xyz");
564   gStyle->SetLabelFont(42,"xyz");
565   gStyle->SetTitleOffset(1.3,"x");
566   gStyle->SetTitleOffset(1.6,"y");
567   gStyle->SetTitleOffset(1.6,"z");
568   gStyle->SetTitleFont(42,"xyz");
569   gStyle->SetTitleColor(1,"xyz");
570   //gStyle->SetPadTopMargin(0.1);
571   //gStyle->SetPadRightMargin(0.1);
572   gStyle->SetPadBottomMargin(0.2);
573   gStyle->SetPadLeftMargin(0.2);
574   
575   
576   const Int_t NCont=255;
577   const Int_t NRGBs = 5;
578   Double_t stops[NRGBs] = { 0.00, 0.34, 0.61, 0.84, 1.00 };
579   Double_t red[NRGBs]   = { 0.00, 0.00, 0.87, 1.00, 0.51 };
580   Double_t green[NRGBs] = { 0.00, 0.81, 1.00, 0.20, 0.00 };
581   Double_t blue[NRGBs]  = { 0.51, 1.00, 0.12, 0.00, 0.00 };
582   TColor::CreateGradientColorTable(NRGBs, stops, red, green, blue, NCont);
583   gStyle->SetNumberContours(NCont); 
584 }