3 //**************************************************************************
4 //* This file is property of and copyright by the ALICE HLT Project *
5 //* ALICE Experiment at CERN, All rights reserved. *
7 //* Primary Authors: Per-Ivar Lønne <perivarlonne@gmail.com> *
8 //* for The ALICE HLT Project. *
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 //**************************************************************************/
19 /// @file AliHLTTPCdEdxMonitoringComponent.cxx
20 /// @author Per-Ivar Lønne, Jochen Thaeder, Matthias Richter, Alexander Kalweit
22 /// @brief Component for reading ESD from chain and produce a dEdx monitoring plot
31 #include "AliESDtrackCuts.h"
32 #include <AliHLTDAQ.h>
33 #include "AliESDEvent.h"
34 #include "AliESDtrack.h"
37 #include "TObjString.h"
38 #include "AliHLTDataTypes.h"
39 #include "AliHLTComponentBenchmark.h"
45 #include "AliHLTTPCdEdxMonitoringComponent.h"
47 /** ROOT macro for the implementation of ROOT specific class methods */
48 ClassImp( AliHLTTPCdEdxMonitoringComponent )
50 AliHLTTPCdEdxMonitoringComponent::AliHLTTPCdEdxMonitoringComponent()
65 AliHLTTPCdEdxMonitoringComponent::~AliHLTTPCdEdxMonitoringComponent()
70 // ######################################################################### //
71 const char* AliHLTTPCdEdxMonitoringComponent::GetComponentID()
74 return "TPCdEdxMonitoring";
77 // ######################################################################### //
78 void AliHLTTPCdEdxMonitoringComponent::GetInputDataTypes( vector<AliHLTComponentDataType>& list)
80 // get list of input data types
81 list.push_back(kAliHLTDataTypeESDObject|kAliHLTDataOriginAny);
84 // ######################################################################### //
85 AliHLTComponentDataType AliHLTTPCdEdxMonitoringComponent::GetOutputDataType()
87 // get the output data size
88 return kAliHLTDataTypeHistogram|kAliHLTDataOriginHLT;
91 // ######################################################################### //
92 void AliHLTTPCdEdxMonitoringComponent::GetOutputDataSize( unsigned long& constBase, Double_t& inputMultiplier )
94 // get output size estimator
99 // ######################################################################### //
100 void AliHLTTPCdEdxMonitoringComponent::GetOCDBObjectDescription( TMap* const targetMap)
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
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"));
117 // ######################################################################### //
118 AliHLTComponent* AliHLTTPCdEdxMonitoringComponent::Spawn()
120 // Spawn function, return new class instance
121 return new AliHLTTPCdEdxMonitoringComponent;
124 // ######################################################################### //
125 Int_t AliHLTTPCdEdxMonitoringComponent::DoInit( Int_t argc, const char** argv )
127 // init the component
129 fESDTrackCuts = new AliESDtrackCuts("AliESDtrackCuts","HLT");
146 delete fESDTrackCuts;
147 fESDTrackCuts = NULL;
152 SetDefaultConfiguration();
153 TString cdbPath="HLT/ConfigTPC/";
154 cdbPath+=GetComponentID();
155 iResult=ConfigureFromCDBTObjString(cdbPath);
160 iResult=ConfigureFromArgumentString(argc, argv);
165 HLTInfo("ESD track cuts : %s",fESDTrackCuts->GetTitle() );
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.)");
175 // ######################################################################### //
176 Int_t AliHLTTPCdEdxMonitoringComponent::DoDeinit()
178 // component cleanup, delete all instances of helper classes here
180 delete fESDTrackCuts;
181 fESDTrackCuts = NULL;
189 // ######################################################################### //
190 Int_t AliHLTTPCdEdxMonitoringComponent::DoEvent(const AliHLTComponentEventData& /*evtData*/,
191 AliHLTComponentTriggerData& /*trigData*/)
193 // event processing function
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;
201 const TObject* obj = GetFirstInputObject(kAliHLTAllDataTypes, "AliESDEvent");
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));
208 esd->GetStdContent();
209 for (Int_t i = 0; i < esd->GetNumberOfTracks(); i++) {
210 AliESDtrack* track = esd->GetTrack(i);
211 sig=track->GetTPCsignal();
212 if(!fESDTrackCuts->AcceptTrack(track)) continue;
213 if (!track->GetInnerParam()) continue;
214 ptot = track->GetInnerParam()->GetP()*track->GetSign();
215 fHist->Fill(ptot, sig);
218 // publish the histogram
219 PushBack(fHist, kAliHLTDataTypeHistogram | kAliHLTDataOriginHLT);
224 // ######################################################################### //
225 Int_t AliHLTTPCdEdxMonitoringComponent::ScanConfigurationArgument(Int_t argc, const char** argv)
227 // Scan configuration arguments
228 // Return the number of processed arguments
229 // -EPROTO if argument format error (e.g. number expected but not found)
231 // The AliHLTComponent base class implements a parsing loop for argument strings and
232 // arrays of strings which is invoked by ConfigureFromArgumentString/ConfigureFromCDBTObjString
233 // The component needs to implement ScanConfigurationArgument in order to decode the arguments.
235 if (argc<=0) return 0;
237 TString argument=argv[ii];
239 if (argument.IsNull()) return 0;
242 HLTError("No ESD track cuts availible");
247 //**********************************//
248 // Histogram Binning //
249 //**********************************//
253 if (argument.CompareTo("-xbins")==0)
255 if (++ii>=argc) return -EPROTO;
258 fxbins = argument.Atoi();
264 if (argument.CompareTo("-xmin")==0)
266 if (++ii>=argc) return -EPROTO;
269 fxmin = argument.Atoi();
274 if (argument.CompareTo("-xmax")==0)
276 if (++ii>=argc) return -EPROTO;
279 fxmax = argument.Atoi();
285 if (argument.CompareTo("-ybins")==0)
287 if (++ii>=argc) return -EPROTO;
290 fybins = argument.Atoi();
296 if (argument.CompareTo("-ymin")==0)
298 if (++ii>=argc) return -EPROTO;
301 fymin = argument.Atoi();
306 if (argument.CompareTo("-ymax")==0)
308 if (++ii>=argc) return -EPROTO;
311 fymax = argument.Atoi();
317 //**********************************//
319 //**********************************//
322 if (argument.CompareTo("-maxpt")==0) {
323 if (++ii>=argc) return -EPROTO;
326 Float_t minPt, maxPt;
327 fESDTrackCuts->GetPtRange(minPt,maxPt);
328 maxPt = argument.Atof();
329 fESDTrackCuts->SetPtRange(minPt,maxPt);
331 TString title = fESDTrackCuts->GetTitle();
332 if (!title.CompareTo("No track cuts")) title = "";
333 else title += " && ";
334 title += Form("p_t < %f", maxPt);
335 fESDTrackCuts->SetTitle(title);
340 if (argument.CompareTo("-minpt")==0) {
341 if (++ii>=argc) return -EPROTO;
344 Float_t minPt, maxPt;
345 fESDTrackCuts->GetPtRange(minPt,maxPt);
346 minPt = argument.Atof();
347 fESDTrackCuts->SetPtRange(minPt,maxPt);
349 TString title = fESDTrackCuts->GetTitle();
350 if (!title.CompareTo("No track cuts")) title = "";
351 else title += " && ";
352 title += Form("p_t > %f", minPt);
353 fESDTrackCuts->SetTitle(title);
358 // minimum longitudinal dca to vertex
359 if (argument.CompareTo("-min-ldca")==0) {
360 if (++ii>=argc) return -EPROTO;
363 fESDTrackCuts->SetMinDCAToVertexZ(argument.Atof());
364 TString title = fESDTrackCuts->GetTitle();
365 if (!title.CompareTo("No track cuts")) title = "";
366 else title += " && ";
367 title += Form("DCAz > %f", argument.Atof());
368 fESDTrackCuts->SetTitle(title);
373 // maximum longitudinal dca to vertex
374 if (argument.CompareTo("-max-ldca")==0) {
375 if (++ii>=argc) return -EPROTO;
378 fESDTrackCuts->SetMaxDCAToVertexZ(argument.Atof());
379 TString title = fESDTrackCuts->GetTitle();
380 if (!title.CompareTo("No track cuts")) title = "";
381 else title += " && ";
382 title += Form("DCAz < %f", argument.Atof());
383 fESDTrackCuts->SetTitle(title);
388 // minimum transverse dca to vertex
389 if (argument.CompareTo("-min-tdca")==0) {
390 if (++ii>=argc) return -EPROTO;
393 fESDTrackCuts->SetMinDCAToVertexXY(argument.Atof());
394 TString title = fESDTrackCuts->GetTitle();
395 if (!title.CompareTo("No track cuts")) title = "";
396 else title += " && ";
397 title += Form("DCAr > %f", argument.Atof());
398 fESDTrackCuts->SetTitle(title);
403 // maximum transverse dca to vertex
404 if (argument.CompareTo("-max-tdca")==0) {
405 if (++ii>=argc) return -EPROTO;
408 fESDTrackCuts->SetMaxDCAToVertexXY(argument.Atof());
409 TString title = fESDTrackCuts->GetTitle();
410 if (!title.CompareTo("No track cuts")) title = "";
411 else title += " && ";
412 title += Form("DCAr < %f", argument.Atof());
413 fESDTrackCuts->SetTitle(title);
419 if (argument.CompareTo("-etarange")==0) {
420 if (++ii>=argc) return -EPROTO;
422 Float_t eta = argument.Atof();
424 fESDTrackCuts->SetEtaRange(-eta,eta);
425 TString title = fESDTrackCuts->GetTitle();
426 if (!title.CompareTo("No track cuts")) title = "";
427 else title += " && ";
428 title += Form("Eta[%f,%f]", argument.Atof(),argument.Atof());
429 fESDTrackCuts->SetTitle(title);
434 // minimum clusters in TPC
435 if (argument.CompareTo("-minNClsTPC")==0) {
436 if (++ii>=argc) return -EPROTO;
440 ncls = Int_t(argument.Atof());
441 fESDTrackCuts->SetMinNClustersTPC(ncls);
443 TString title = fESDTrackCuts->GetTitle();
444 if (!title.CompareTo("No track cuts")) title = "";
445 else title += " && ";
446 title += Form("minNClsTPC < %i", ncls);
447 fESDTrackCuts->SetTitle(title);
460 // ######################################################################### //
461 Int_t AliHLTTPCdEdxMonitoringComponent::Reconfigure(const char* cdbEntry, const char* chainId)
463 // reconfigure the component from the specified CDB entry, or default CDB entry
464 HLTInfo("reconfigure '%s' from entry %s", chainId, cdbEntry);
472 cdbPath="HLT/ConfigTPC/";
473 cdbPath+=GetComponentID();
476 iResult=ConfigureFromCDBTObjString(cdbPath); //// Or use return 0, and skip this line?
481 // ######################################################################### //
482 Int_t AliHLTTPCdEdxMonitoringComponent::ReadPreprocessorValues(const char* modules)
484 // read the preprocessor values for the detectors in the modules list
486 TString detectors(modules!=NULL?modules:"");
487 HLTInfo("read preprocessor values for detector(s): %s", detectors.IsNull()?"none":detectors.Data());
488 return iResult; //Done differently in AliHLTMultiplicityCorrelationsComponent...
493 // ######################################################################### //
494 void AliHLTTPCdEdxMonitoringComponent::SetDefaultConfiguration()
498 //fESDTrackCuts = AliESDtrackCuts::GetStandardITSTPCTrackCuts2010(kTRUE);
499 fESDTrackCuts->SetEtaRange(-0.8,+0.8);
500 fESDTrackCuts->SetPtRange(0.15,1e10);
501 fESDTrackCuts->SetMaxCovDiagonalElements(2, 2, 0.5, 0.5, 2); // BEWARE STANDARD VALUES ARE: 2, 2, 0.5, 0.5, 2
502 fESDTrackCuts->SetMaxNsigmaToVertex(3);
503 fESDTrackCuts->SetRequireSigmaToVertex(kTRUE);
504 fESDTrackCuts->SetAcceptKinkDaughters(kFALSE);
505 fESDTrackCuts->SetMinNClustersTPC(70);
506 fESDTrackCuts->SetMaxChi2PerClusterTPC(4);
507 fESDTrackCuts->SetMaxDCAToVertexXY(3);
508 fESDTrackCuts->SetMaxDCAToVertexZ(3);
509 fESDTrackCuts->SetRequireTPCRefit(kTRUE);
510 //fESDTrackCuts->SetRequireITSRefit(kTRUE); //Kills HLT simulated reconstructions?
511 fESDTrackCuts->SetClusterRequirementITS(AliESDtrackCuts::kSPD, AliESDtrackCuts::kAny); //TEMPORARY <-> REMOVE
512 fESDTrackCuts->SetMinNClustersITS(3);
526 // ######################################################################### //
527 void AliHLTTPCdEdxMonitoringComponent::Plotstyle()
529 gROOT->SetStyle("Plain");
531 gStyle->SetCanvasBorderMode(0);
532 gStyle->SetCanvasColor(10);
533 gStyle->SetCanvasDefH(550);
534 gStyle->SetCanvasDefW(575);
535 gStyle->SetPadBorderMode(0);
536 gStyle->SetPadColor(10);
537 gStyle->SetPadTickX(1);
538 gStyle->SetPadTickY(1);
539 gStyle->SetStatColor(10);
540 gStyle->SetFrameFillColor(10);
541 gStyle->SetPalette(1,0);
545 //gStyle->SetStatX(0.7);
546 //gStyle->SetStatW(0.2);
547 //gStyle->SetLabelOffset(1.2);
548 //gStyle->SetLabelFont(72);
549 //gStyle->SetLabelSize(0.6);
550 //gStyle->SetTitleOffset(1.2);
551 gStyle->SetTitleFontSize(0.04);
554 gStyle->SetOptStat(10);
555 gStyle->SetLineWidth(2);
556 gStyle->SetMarkerSize(1.0);
557 gStyle->SetTextSize(0.04);
558 gStyle->SetTitleSize(0.04,"xyz");
559 gStyle->SetLabelSize(0.04,"xyz");
560 gStyle->SetLabelOffset(0.02,"xyz");
561 gStyle->SetLabelFont(42,"xyz");
562 gStyle->SetTitleOffset(1.3,"x");
563 gStyle->SetTitleOffset(1.6,"y");
564 gStyle->SetTitleOffset(1.6,"z");
565 gStyle->SetTitleFont(42,"xyz");
566 gStyle->SetTitleColor(1,"xyz");
567 //gStyle->SetPadTopMargin(0.1);
568 //gStyle->SetPadRightMargin(0.1);
569 gStyle->SetPadBottomMargin(0.2);
570 gStyle->SetPadLeftMargin(0.2);
573 const Int_t NCont=255;
574 const Int_t NRGBs = 5;
575 Double_t stops[NRGBs] = { 0.00, 0.34, 0.61, 0.84, 1.00 };
576 Double_t red[NRGBs] = { 0.00, 0.00, 0.87, 1.00, 0.51 };
577 Double_t green[NRGBs] = { 0.00, 0.81, 1.00, 0.20, 0.00 };
578 Double_t blue[NRGBs] = { 0.51, 1.00, 0.12, 0.00, 0.00 };
579 TColor::CreateGradientColorTable(NRGBs, stops, red, green, blue, NCont);
580 gStyle->SetNumberContours(NCont);