2 /**************************************************************************
3 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
5 * Author: The ALICE Off-line Project. *
6 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
17 /////////////////////////////////////////////////////////////////////////////////
21 // Offline TRD calibration task
24 // R. Bailhache (R.Bailhache@gsi.de)
26 //////////////////////////////////////////////////////////////////////////////////////
30 #include "Riostream.h"
34 #include "TProfile2D.h"
41 #include "TObjArray.h"
47 #include "TIterator.h"
49 #include "AliAnalysisTask.h"
50 #include "AliAnalysisManager.h"
52 #include "AliESDVertex.h"
53 #include "AliESDEvent.h"
54 #include "AliESDfriend.h"
55 #include "AliESDInputHandler.h"
56 #include "AliESDtrack.h"
57 #include "AliESDfriendTrack.h"
58 #include "AliTRDtrackV1.h"
59 #include "AliTRDseedV1.h"
60 #include "AliTRDcluster.h"
61 #include "AliTRDgeometry.h"
62 #include "AliESDtrackCuts.h"
63 #include "AliESDVertex.h"
64 #include "AliTRDCalDet.h"
66 #include "AliTRDCalibraVector.h"
67 #include "AliTRDCalibraFillHisto.h"
68 #include "AliTRDCalibraVdriftLinearFit.h"
70 #include "AliTRDcalibDB.h"
75 #include "AliTRDCalibTask.h"
78 ClassImp(AliTRDCalibTask)
80 //________________________________________________________________________
81 AliTRDCalibTask::AliTRDCalibTask(const char *name)
82 : AliAnalysisTaskSE(name), fESD(0),
90 fTRDCalibraFillHisto(0),
93 fNbTRDTrackOffline(0),
94 fNbTRDTrackStandalone(0),
98 fNbTimeBinStandalone(0),
100 fNbClustersOffline(0),
101 fNbClustersStandalone(0),
103 fNbTrackletsOffline(0),
104 fNbTrackletsStandalone(0),
112 fVdriftLinear(kTRUE),
114 fSelectedTrigger(new TObjArray()),
117 fRequirePrimaryVertex(kFALSE),
120 fMinNbContributors(0),
121 fRangePrimaryVertexZ(9999999.0),
125 fNormalizeNbOfCluster(kFALSE),
129 fOfflineTracks(kFALSE),
130 fStandaloneTracks(kFALSE),
131 fVersionGainUsed(-1),
132 fSubVersionGainUsed(-1),
133 fVersionGainLocalUsed(-1),
134 fSubVersionGainLocalUsed(-1),
135 fVersionVdriftUsed(-1),
136 fSubVersionVdriftUsed(-1),
143 // Default constructor
154 // Define input and output slots here
155 // Input slot #0 works with a TChain
156 DefineInput(0, TChain::Class());
158 // Output slot #0 writes into a TList container
159 DefineOutput(1, TList::Class());
163 //____________________________________________________________________________________
164 AliTRDCalibTask::~AliTRDCalibTask()
167 // AliTRDCalibTask destructor
171 if(fNEvents) delete fNEvents;
172 if(fNbTRDTrack) delete fNbTRDTrack;
173 if(fNbTRDTrackOffline) delete fNbTRDTrackOffline;
174 if(fNbTRDTrackStandalone) delete fNbTRDTrackStandalone;
175 if(fNbTPCTRDtrack) delete fNbTPCTRDtrack;
176 if(fNbTimeBin) delete fNbTimeBin;
177 if(fNbTimeBinOffline) delete fNbTimeBinOffline;
178 if(fNbTimeBinStandalone) delete fNbTimeBinStandalone;
179 if(fNbClusters) delete fNbClusters;
180 if(fNbClustersOffline) delete fNbClustersOffline;
181 if(fNbClustersStandalone) delete fNbClustersStandalone;
182 if(fNbTracklets) delete fNbTracklets;
183 if(fNbTrackletsOffline) delete fNbTrackletsOffline;
184 if(fNbTrackletsStandalone) delete fNbTrackletsStandalone;
185 if(fAbsoluteGain) delete fAbsoluteGain;
186 if(fCH2dSum) delete fCH2dSum;
187 if(fPH2dSum) delete fPH2dSum;
188 if(fCH2dSM) delete fCH2dSM;
189 if(fPH2dSM) delete fPH2dSM;
190 if(fCalDetGain) delete fCalDetGain;
192 if(fSelectedTrigger) {
193 fSelectedTrigger->Delete();
194 delete fSelectedTrigger;
197 delete fEsdTrackCuts;
203 //________________________________________________________________________
204 void AliTRDCalibTask::ConnectInputData(Option_t *)
206 // Connect ESD or AOD here
207 // Called once per event
209 cout << "AliTRDCalibTask::ConnectInputData() IN" << endl;
212 // TTree* tree = dynamic_cast<TTree*> (GetInputData(0)); //pointer wird "umgecastet" auf anderen Variablentyp
214 //Printf("ERROR: Could not read chain from input slot 0");
217 AliESDInputHandler *esdH = dynamic_cast<AliESDInputHandler*> (AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler());
220 //Printf("ERROR: Could not get ESDInputHandler");
222 fESD = esdH->GetEvent();
223 // esdH->SetActiveBranches("ESDfriend*");
224 if ((esdH->GetTree())->GetBranch("ESDfriend.")) fESDfriend = esdH->GetESDfriend();
225 //else printf("No friend ESD\n");
226 //Printf("*** CONNECTED NEW EVENT ****");
231 //cout << "AliTRDCalibTask::ConnectInputData() OUT" << endl;
236 //________________________________________________________________________
237 void AliTRDCalibTask::UserCreateOutputObjects()
242 //cout << "AliTRDCalibTask::CreateOutputObjects() IN" << endl;
244 // Number of time bins
246 AliTRDcalibDB *cal = AliTRDcalibDB::Instance();
247 fNbTimeBins = cal->GetNumberOfTimeBinsDCS();
248 if(fNbTimeBins <= 0){
249 AliWarning(Form("No of TimeBins from DB [%d] use default [30]", fNbTimeBins));
254 // instance calibration
255 fTRDCalibraFillHisto = AliTRDCalibraFillHisto::Instance();
256 fTRDCalibraFillHisto->SetHisto2d(fHisto2d); // choose to use histograms
257 fTRDCalibraFillHisto->SetVector2d(fVector2d); // choose to use vectors
258 fTRDCalibraFillHisto->SetCH2dOn(); // choose to calibrate the gain
259 fTRDCalibraFillHisto->SetPH2dOn(); // choose to calibrate the drift velocity
260 fTRDCalibraFillHisto->SetPRF2dOn(); // choose to look at the PRF
261 fTRDCalibraFillHisto->SetLinearFitterOn(fVdriftLinear); // Other possibility vdrift VDRIFT
262 fTRDCalibraFillHisto->SetLinearFitterDebugOn(fVdriftLinear); // Other possibility vdrift
263 for(Int_t k = 0; k < 3; k++){
264 if(((fNz[k] != 10) && (fNrphi[k] != 10)) && ((fNz[k] != 100) && (fNrphi[k] != 100))) {
265 fTRDCalibraFillHisto->SetNz(k,fNz[k]); // Mode calibration
266 fTRDCalibraFillHisto->SetNrphi(k,fNrphi[k]); // Mode calibration
269 if((fNz[k] == 100) && (fNrphi[k] == 100)) {
270 if(fVector2d) AliInfo("The mode all together is not supported by the vector method");
271 fTRDCalibraFillHisto->SetAllTogether(k);
273 if((fNz[k] == 10) && (fNrphi[k] == 10)) {
274 if(fVector2d) AliInfo("The mode per supermodule is not supported by the vector method");
275 fTRDCalibraFillHisto->SetPerSuperModule(k);
279 // Variables for how to fill
280 fTRDCalibraFillHisto->SetFillWithZero(fFillZero);
281 fTRDCalibraFillHisto->SetNormalizeNbOfCluster(fNormalizeNbOfCluster);
282 fTRDCalibraFillHisto->SetMaxCluster(fMaxCluster);
283 fTRDCalibraFillHisto->SetNbMaxCluster(fNbMaxCluster);
285 // Init with 30 timebins
286 fTRDCalibraFillHisto->Init2Dhistos(fNbTimeBins); // initialise the histos
287 fTRDCalibraFillHisto->SetNumberClusters(fLow); // At least 11 clusters
288 fTRDCalibraFillHisto->SetNumberClustersf(fHigh); // At least 11 clusters
289 fRelativeScale = fTRDCalibraFillHisto->GetRelativeScale(); // Get the relative scale for the gain
292 if(fDebug > 2) fTRDCalibraFillHisto->SetDebugLevel(1); //debug stuff
295 fListHist = new TList();
297 fListHist->Add(fTRDCalibraFillHisto->GetCH2d());
298 fListHist->Add(fTRDCalibraFillHisto->GetPH2d());
299 fListHist->Add(fTRDCalibraFillHisto->GetPRF2d());
301 if(fVdriftLinear) fListHist->Add((TObject *) fTRDCalibraFillHisto->GetVdriftLinearFit());
302 if(fVector2d) fListHist->Add((TObject *) fTRDCalibraFillHisto->GetCalibraVector()); //calibra vector
303 fNEvents = new TH1I("NEvents","NEvents", 2, 0, 2);
304 fListHist->Add(fNEvents);
306 // absolute gain calibration even without AliESDfriend
308 Double_t minPt = 0.001;
309 Double_t maxPt = 10.0;
311 Double_t *binLimLogPt = new Double_t[nBinsPt+1];
312 Double_t *binLimPt = new Double_t[nBinsPt+1];
313 for(Int_t i=0; i<=nBinsPt; i++) binLimLogPt[i]=(Double_t)TMath::Log10(minPt) + (TMath::Log10(maxPt)-TMath::Log10(minPt))/nBinsPt*(Double_t)i ;
314 for(Int_t i=0; i<=nBinsPt; i++) binLimPt[i]=(Double_t)TMath::Power(10,binLimLogPt[i]);
316 fAbsoluteGain = new TH2F("AbsoluteGain","AbsoluteGain", 200, 0.0, 700.0, nBinsPt, binLimPt);
317 fAbsoluteGain->SetYTitle("Momentum at TRD");
318 fAbsoluteGain->SetXTitle("charge deposit [a.u]");
319 fAbsoluteGain->SetZTitle("counts");
320 fAbsoluteGain->SetStats(0);
321 fAbsoluteGain->Sumw2();
322 fListHist->Add(fAbsoluteGain);
326 /////////////////////////////////////////
328 ///////////////////////////////////////
331 // Standart with AliESDfriend
332 fPH2dSM = new TProfile2D("PH2dSM","Nz10Nrphi10"
333 ,fNbTimeBins,-0.05,(Double_t)((fNbTimeBins-0.5)/10.0)
335 fPH2dSM->SetYTitle("Det/pad groups");
336 fPH2dSM->SetXTitle("time [#mus]");
337 fPH2dSM->SetZTitle("<PH> [a.u.]");
338 fPH2dSM->SetStats(0);
340 fCH2dSM = new TH2I("CH2dSM","Nz10Nrphi10",50,0,300,18,0,18);
341 fCH2dSM->SetYTitle("Det/pad groups");
342 fCH2dSM->SetXTitle("charge deposit [a.u]");
343 fCH2dSM->SetZTitle("counts");
344 fCH2dSM->SetStats(0);
347 fPH2dSum = new TProfile2D("PH2dSum","Nz100Nrphi100"
348 ,fNbTimeBins,-0.05,(Double_t)((fNbTimeBins-0.5)/10.0)
350 fPH2dSum->SetYTitle("Det/pad groups");
351 fPH2dSum->SetXTitle("time [#mus]");
352 fPH2dSum->SetZTitle("<PH> [a.u.]");
353 fPH2dSum->SetStats(0);
355 fCH2dSum = new TH2I("CH2dSum","Nz100Nrphi100",50,0,300,1,0,1);
356 fCH2dSum->SetYTitle("Det/pad groups");
357 fCH2dSum->SetXTitle("charge deposit [a.u]");
358 fCH2dSum->SetZTitle("counts");
359 fCH2dSum->SetStats(0);
363 fListHist->Add(fPH2dSM);
364 fListHist->Add(fCH2dSM);
365 fListHist->Add(fPH2dSum);
366 fListHist->Add(fCH2dSum);
369 /////////////////////////////////////////
370 // Second debug level
371 ///////////////////////////////////////
374 fNbTRDTrack = new TH1F("TRDTrack","TRDTrack",50,0,50);
375 fNbTRDTrack->Sumw2();
376 fNbTRDTrackOffline = new TH1F("TRDTrackOffline","TRDTrackOffline",50,0,50);
377 fNbTRDTrackOffline->Sumw2();
378 fNbTRDTrackStandalone = new TH1F("TRDTrackStandalone","TRDTrackStandalone",50,0,50);
379 fNbTRDTrackStandalone->Sumw2();
380 fNbTPCTRDtrack = new TH2F("NbTPCTRDtrack","NbTPCTRDtrack",100,0,100,100,0,100);
381 fNbTPCTRDtrack->Sumw2();
383 fNbTimeBin = new TH1F("NbTimeBin","NbTimeBin",35,0,35);
385 fNbTimeBinOffline = new TH1F("NbTimeBinOffline","NbTimeBinOffline",35,0,35);
386 fNbTimeBinOffline->Sumw2();
387 fNbTimeBinStandalone = new TH1F("NbTimeBinStandalone","NbTimeBinStandalone",35,0,35);
388 fNbTimeBinStandalone->Sumw2();
390 fNbClusters = new TH1F("NbClusters","",35,0,35);
391 fNbClusters->Sumw2();
392 fNbClustersOffline = new TH1F("NbClustersOffline","",35,0,35);
393 fNbClustersOffline->Sumw2();
394 fNbClustersStandalone = new TH1F("NbClustersStandalone","",35,0,35);
395 fNbClustersStandalone->Sumw2();
397 fNbTracklets = new TH1F("NbTracklets","NbTracklets",540,0.,540.);
398 fNbTracklets->Sumw2();
399 fNbTrackletsOffline = new TH1F("NbTrackletsOffline","NbTrackletsOffline",540,0.,540.);
400 fNbTrackletsOffline->Sumw2();
401 fNbTrackletsStandalone = new TH1F("NbTrackletsStandalone","NbTrackletsStandalone",540,0.,540.);
402 fNbTrackletsStandalone->Sumw2();
404 fListHist->Add(fNbTRDTrack);
405 fListHist->Add(fNbTRDTrackOffline);
406 fListHist->Add(fNbTRDTrackStandalone);
407 fListHist->Add(fNbTPCTRDtrack);
409 fListHist->Add(fNbTimeBin);
410 fListHist->Add(fNbTimeBinOffline);
411 fListHist->Add(fNbTimeBinStandalone);
412 fListHist->Add(fNbClusters);
413 fListHist->Add(fNbClustersOffline);
414 fListHist->Add(fNbClustersStandalone);
415 fListHist->Add(fNbTracklets);
416 fListHist->Add(fNbTrackletsOffline);
417 fListHist->Add(fNbTrackletsStandalone);
421 delete [] binLimLogPt;
424 //cout << "AliTRDCalibTask::UserCreateOutputObjects() OUT" << endl;
428 //________________________________________________________________________
429 void AliTRDCalibTask::UserExec(Option_t *)
432 // Filling of the histos
434 //cout << "AliTRDCalibTask::Exec() IN" << endl;
436 // Init Versions and subversions used
437 if((fVersionGainUsed==-1) || (fSubVersionGainUsed==-1) || (fVersionGainLocalUsed==-1) || (fSubVersionGainLocalUsed==-1) || (fVersionVdriftUsed==-1) || (fSubVersionVdriftUsed==-1)) {
438 if(!SetVersionSubversion()) {
440 fSubVersionGainUsed=0;
441 fVersionGainLocalUsed=0;
442 fSubVersionGainLocalUsed=0;
443 fVersionVdriftUsed=0;
444 fSubVersionVdriftUsed=0;
448 fTRDCalibraFillHisto->SetVersionGainUsed(fVersionGainUsed); // Gain Used
449 fTRDCalibraFillHisto->SetSubVersionGainUsed(fSubVersionGainUsed); // Gain Used
450 fTRDCalibraFillHisto->SetVersionGainLocalUsed(fVersionGainLocalUsed); // Gain Used
451 fTRDCalibraFillHisto->SetSubVersionGainLocalUsed(fSubVersionGainLocalUsed); // Gain Used
452 fTRDCalibraFillHisto->SetVersionVdriftUsed(fVersionVdriftUsed); // Vdrift Used
453 fTRDCalibraFillHisto->SetSubVersionVdriftUsed(fSubVersionVdriftUsed); // Vdrift Used
454 fTRDCalibraFillHisto->InitCalDet();
457 // AliLog::SetGlobalLogLevel(AliLog::kError);
458 // cout << "AliTRDCalibTask::Exec() 1" << endl;
459 fESD = dynamic_cast<AliESDEvent*>(fInputEvent);
461 AliError("ESD Event missing");
462 PostData(1, fListHist);
466 //printf("Counter %d\n",fCounter);
469 //cout << "maxEvent = " << fMaxEvent << endl;
470 //if(fCounter%100==0) cout << "fCounter = " << fCounter << endl;
471 if((fMaxEvent != 0) && (fMaxEvent < fCounter)) return;
472 //if(fCounter%100==0) cout << "fCounter1 = " << fCounter << endl;
473 //cout << "event = " << fCounter << endl;
475 //printf("Counter %d\n",fCounter);
481 Int_t numberOfTriggerSelected = fSelectedTrigger->GetEntriesFast();
482 //printf("numberofTriggerSelected %d\n",numberOfTriggerSelected);
485 for(Int_t k = 0; k < numberOfTriggerSelected; k++){
486 const TObjString *const obString=(TObjString*)fSelectedTrigger->At(k);
487 const TString tString=obString->GetString();
488 if(fESD->IsTriggerClassFired((const char*)tString)) {
495 for(Int_t k = 0; k < numberOfTriggerSelected; k++){
496 const TObjString *const obString=(TObjString*)fSelectedTrigger->At(k);
497 const TString tString=obString->GetString();
498 if(fESD->IsTriggerClassFired((const char*)tString)) {
504 PostData(1, fListHist);
507 //printf("Class Fired %s\n",(const char*)fESD->GetFiredTriggerClasses());
508 //printf("Trigger passed\n");
510 ///////////////////////////////
511 // Require a primary vertex
512 //////////////////////////////
513 if(fRequirePrimaryVertex) {
514 const AliESDVertex* vtxESD = 0x0;
515 if (fVtxTPC) vtxESD = fESD->GetPrimaryVertexTPC() ;
516 else if (fVtxSPD) vtxESD = fESD->GetPrimaryVertexSPD() ;
517 else vtxESD = fESD->GetPrimaryVertexTracks() ;
519 PostData(1, fListHist);
522 Int_t nCtrb = vtxESD->GetNContributors();
523 if(nCtrb < fMinNbContributors) {
524 PostData(1, fListHist);
527 Double_t zPosition = vtxESD->GetZ();
528 if(TMath::Abs(zPosition) > fRangePrimaryVertexZ) {
529 PostData(1, fListHist);
535 //printf("Primary vertex passed\n");
540 Int_t nbTrdTracks = 0;
542 Int_t nbTrdTracksStandalone = 0;
544 Int_t nbTrdTracksOffline = 0;
546 Int_t nbtrackTPC = 0;
548 Double_t nbTracks = fESD->GetNumberOfTracks();
549 //printf("Number of tracks %f\n",nbTracks);
551 if (nbTracks <= 0.0) {
554 fNbTRDTrack->Fill(nbTrdTracks);
555 fNbTRDTrackStandalone->Fill(nbTrdTracksStandalone);
556 fNbTRDTrackOffline->Fill(nbTrdTracksOffline);
558 PostData(1, fListHist);
563 fESDfriend = dynamic_cast<AliESDfriend*> (fESD->FindListObject("AliESDfriend"));
565 AliError("fESDfriend not available");
566 PostData(1, fListHist);
570 //printf("has friends\n");
573 ////////////////////////////////////
574 // Check the number of TPC tracks
575 ///////////////////////////////////
576 //printf("Nb of tracks %f\n",nbTracks);
577 for(Int_t itrk = 0; itrk < nbTracks; itrk++){
579 fkEsdTrack = fESD->GetTrack(itrk);
580 ULong_t status = fkEsdTrack->GetStatus();
581 if(status&(AliESDtrack::kTPCout)) nbtrackTPC++;
582 if((status&(AliESDtrack::kTRDout)) && (!(status&(AliESDtrack::kTRDin)))) {
584 nbTrdTracksStandalone++;
586 if((status&(AliESDtrack::kTRDin))) {
588 nbTrdTracksOffline++;
592 if((nbtrackTPC>0) && (nbTrdTracks > (3.0*nbtrackTPC))) pass = kFALSE;
596 fNbTRDTrack->Fill(nbTrdTracks);
597 fNbTRDTrackStandalone->Fill(nbTrdTracksStandalone);
598 fNbTRDTrackOffline->Fill(nbTrdTracksOffline);
599 fNbTPCTRDtrack->Fill(nbTrdTracks,nbtrackTPC);
604 PostData(1, fListHist);
609 /////////////////////////////////////
610 // Loop on AliESDtrack
611 ////////////////////////////////////
612 //printf("Nb of tracks %f\n",nbTracks);
613 for(int itrk=0; itrk < nbTracks; ++itrk){
616 fkEsdTrack = fESD->GetTrack(itrk);
617 if(!fkEsdTrack) continue;
618 ULong_t status = fkEsdTrack->GetStatus();
619 if(status&(AliESDtrack::kTPCout)) ++nbtrackTPC;
621 // Quality cuts on the AliESDtrack
622 if((fEsdTrackCuts) && (!fEsdTrackCuts->IsSelected((AliVParticle *)fkEsdTrack))) {
623 //printf("Not a good track\n");
627 // First Absolute gain calibration
628 Int_t trdNTracklets = (Int_t) fkEsdTrack->GetTRDntracklets();
629 Int_t trdNTrackletsPID = (Int_t) fkEsdTrack->GetTRDntrackletsPID();
630 if((trdNTracklets > 0) && (trdNTrackletsPID > 0)) {
631 for(Int_t iPlane = 0; iPlane < 6; ++iPlane){
632 //Double_t slide = fkEsdTrack->GetTRDslice(iPlane);
633 //printf("Number of slide %d\n",fkEsdTrack->GetNumberOfTRDslices());
634 //Double_t momentum = fkEsdTrack->GetTRDmomentum(iPlane);
635 //printf("momentum %f, slide %f\n",momentum,slide);
636 if(fkEsdTrack->GetTRDslice(iPlane) > 0.0)
637 fAbsoluteGain->Fill(fkEsdTrack->GetTRDslice(iPlane)*8.0/100.0,
638 fkEsdTrack->GetTRDmomentum(iPlane));
644 Bool_t standalonetrack = kFALSE;
645 Bool_t offlinetrack = kFALSE;
646 //ULong_t status = fkEsdTrack->GetStatus();
648 fFriendTrack = fESDfriend->GetTrack(itrk);
650 //printf("No friend track %d\n",itrk);
653 //////////////////////////////////////
654 // Loop on calibration objects
655 //////////////////////////////////////
658 while((fCalibObject = (TObject *)(fFriendTrack->GetCalibObject(icalib++)))){
659 //printf("Name %s\n",fCalibObject->IsA()->GetName());
660 if(strcmp(fCalibObject->IsA()->GetName(), "AliTRDtrackV1") != 0) continue;
661 //printf("Find the calibration object\n");
664 if((status&(AliESDtrack::kTRDout)) && (!(status&(AliESDtrack::kTRDin)))) {
665 standalonetrack = kTRUE;
667 if((status&(AliESDtrack::kTRDin))) {
668 offlinetrack = kTRUE;
675 else if(fStandaloneTracks){
676 if(!standalonetrack){
681 fTrdTrack = (AliTRDtrackV1 *)fCalibObject;
683 //cout << "good" << endl;
684 fTRDCalibraFillHisto->UpdateHistogramsV1(fTrdTrack);
685 //printf("Fill fTRDCalibraFillHisto\n");
688 //////////////////////////////////
690 ////////////////////////////////
694 //printf("Enter debug\n");
696 Int_t nbtracklets = 0;
699 Bool_t standalonetracklet = kFALSE;
700 const AliTRDseedV1 *tracklet = 0x0;
701 //////////////////////////////////////
703 /////////////////////////////////////
705 Double_t phtb[AliTRDseedV1::kNtb];
706 memset(phtb, 0, AliTRDseedV1::kNtb*sizeof(Double_t));
708 Float_t normalisation = 6.67;
711 for(Int_t itr = 0; itr < 6; ++itr){
713 if(!(tracklet = fTrdTrack->GetTracklet(itr))) continue;
714 if(!tracklet->IsOK()) continue;
716 standalonetracklet = kFALSE;
717 if(tracklet->IsStandAlone()) standalonetracklet = kTRUE;
720 memset(phtb, 0, AliTRDseedV1::kNtb*sizeof(Double_t));
724 //Int_t crossrow = 0;
726 // Check no shared clusters
727 //for(int icc=AliTRDseedV1::kNtb; icc<AliTRDseedV1::kNclusters; icc++){
728 // if((fcl = tracklet->GetClusters(icc))) crossrow = 1;
735 for(int ic=0; ic<AliTRDseedV1::kNtb; ++ic){
737 if(!(fCl = tracklet->GetClusters(ic))) continue;
739 time = fCl->GetPadTime();
740 ch = tracklet->GetdQdl(ic);
741 qcl = TMath::Abs(fCl->GetQ());
742 detector = fCl->GetDetector();
743 // Add the charge if shared cluster
744 if((ic+AliTRDseedV1::kNtb) < AliTRDseedV1::kNclusters) {
745 if((fCl = tracklet->GetClusters(ic+AliTRDseedV1::kNtb))) {
746 qcl += TMath::Abs(fCl->GetQ());
747 //printf("Add the cluster charge\n");
750 if((time>-1) && (time<fNbTimeBins)) phtb[time]=qcl;
751 if((fCalDetGain) && (fCalDetGain->GetValue(detector) > 0.0)) sum += ch*fCalDetGain->GetValue(detector)/normalisation;
752 else sum += ch/normalisation;
755 fNbTimeBin->Fill(time);
756 if(tracklet->IsStandAlone()) fNbTimeBinStandalone->Fill(time);
757 else fNbTimeBinOffline->Fill(time);
760 sector = AliTRDgeometry::GetSector(detector);
763 fNbTracklets->Fill(detector);
764 if(tracklet->IsStandAlone()) fNbTrackletsStandalone->Fill(detector);
765 else fNbTrackletsOffline->Fill(detector);
767 fNbClusters->Fill(nbclusters);
768 if(tracklet->IsStandAlone()) fNbClustersStandalone->Fill(nbclusters);
769 else fNbClustersOffline->Fill(nbclusters);
773 if((nbclusters > fLow) && (nbclusters < fHigh)){
774 if(fRelativeScale > 0.0) sum = sum/fRelativeScale;
775 fCH2dSM->Fill(sum,sector+0.5);
776 fCH2dSum->Fill(sum,0.5);
777 Bool_t checknoise = kTRUE;
778 if(fMaxCluster > 0) {
779 if(phtb[0] > fMaxCluster) checknoise = kFALSE;
780 if(fNbTimeBins > fNbMaxCluster) {
781 for(Int_t k = (fNbTimeBins-fNbMaxCluster); k < fNbTimeBins; k++){
782 if(phtb[k] > fMaxCluster) checknoise = kFALSE;
787 for(int ic=0; ic<fNbTimeBins; ic++){
789 fPH2dSum->Fill((Double_t)(ic/10.0),0.5,(Double_t)phtb[ic]);
790 fPH2dSM->Fill((Double_t)(ic/10.0),sector+0.5,(Double_t)phtb[ic]);
794 fPH2dSum->Fill((Double_t)(ic/10.0),0.0,(Double_t)phtb[ic]);
795 fPH2dSM->Fill((Double_t)(ic/10.0),sector+0.5,(Double_t)phtb[ic]);
802 } // loop on tracklets
806 }// while calibration objects
807 if(nTRDtrackV1 > 0) {
809 if((status&(AliESDtrack::kTRDout)) && (!(status&(AliESDtrack::kTRDin)))) {
810 ++nbTrdTracksStandalone;
812 if((status&(AliESDtrack::kTRDin))) {
813 ++nbTrdTracksOffline;
816 //delete fFriendTrack;
820 fNbTRDTrack->Fill(nbTrdTracks);
821 fNbTRDTrackStandalone->Fill(nbTrdTracksStandalone);
822 fNbTRDTrackOffline->Fill(nbTrdTracksOffline);
823 fNbTPCTRDtrack->Fill(nbTrdTracks,nbtrackTPC);
827 PostData(1, fListHist);
828 //cout << "AliTRDCalibTask::Exec() OUT" << endl;
831 //________________________________________________________________________
832 void AliTRDCalibTask::Terminate(Option_t *)
838 if(fTRDCalibraFillHisto) fTRDCalibraFillHisto->DestroyDebugStreamer();
842 //_______________________________________________________
843 Bool_t AliTRDCalibTask::Load(const Char_t *filename)
846 // Generic container loader
849 if(!TFile::Open(filename)){
850 //AliWarning(Form("Couldn't open file %s.", filename));
854 if(!(o = (TList*)gFile->Get(GetName()))){
855 //AliWarning("Missing histogram container.");
858 fListHist = (TList*)o->Clone(GetName());
862 //_______________________________________________________
863 Bool_t AliTRDCalibTask::Load(TList *lister)
866 // Generic container loader
869 fListHist = (TList*)lister->Clone(GetName());
872 //________________________________________________________________________
873 void AliTRDCalibTask::Plot()
876 // Plot the histos stored in the TList
879 if(!fListHist) return;
881 /////////////////////////////////////
882 // Take the debug stuff
883 /////////////////////////////////////
885 TH1I *nEvents = (TH1I *) fListHist->FindObject("NEvents");
887 TH2F *absoluteGain = (TH2F *) fListHist->FindObject("AbsoluteGain");
889 TH1F *trdTrack = (TH1F *) fListHist->FindObject("TRDTrack");
890 TH1F *trdTrackOffline = (TH1F *) fListHist->FindObject("TRDTrackOffline");
891 TH1F *trdTrackStandalone = (TH1F *) fListHist->FindObject("TRDTrackStandalone");
893 TH2F *tpctrdTrack = (TH2F *) fListHist->FindObject("NbTPCTRDtrack");
895 TH1F *nbTimeBin = (TH1F *) fListHist->FindObject("NbTimeBin");
896 TH1F *nbTimeBinOffline = (TH1F *) fListHist->FindObject("NbTimeBinOffline");
897 TH1F *nbTimeBinStandalone = (TH1F *) fListHist->FindObject("NbTimeBinStandalone");
899 TH1F *nbClusters = (TH1F *) fListHist->FindObject("NbClusters");
900 TH1F *nbClustersOffline = (TH1F *) fListHist->FindObject("NbClustersOffline");
901 TH1F *nbClustersStandalone = (TH1F *) fListHist->FindObject("NbClustersStandalone");
903 TH1F *nbTracklets = (TH1F *) fListHist->FindObject("NbTracklets");
904 TH1F *nbTrackletsOffline = (TH1F *) fListHist->FindObject("NbTrackletsOffline");
905 TH1F *nbTrackletsStandalone = (TH1F *) fListHist->FindObject("NbTrackletsStandalone");
907 /////////////////////////////////////
908 // Take the calibration objects
909 /////////////////////////////////////
911 TH2I *ch2d = (TH2I *) fListHist->FindObject("CH2d");
912 TProfile2D *ph2d = (TProfile2D *) fListHist->FindObject("PH2d");
914 TH2I *ch2dSum = (TH2I *) fListHist->FindObject("CH2dSum");
915 TProfile2D *ph2dSum = (TProfile2D *) fListHist->FindObject("PH2dSum");
917 TH2I *ch2dSM = (TH2I *) fListHist->FindObject("CH2dSM");
918 TProfile2D *ph2dSM = (TProfile2D *) fListHist->FindObject("PH2dSM");
920 AliTRDCalibraVdriftLinearFit *linearfit = (AliTRDCalibraVdriftLinearFit *) fListHist->FindObject("AliTRDCalibraVdriftLinearFit");
922 ////////////////////////////////////////////////
923 // Add the AliTRDCalibraVdriftLinearFit
924 ///////////////////////////////////////////////
927 TH2S *histolinearfitsum = 0x0;
930 for(Int_t det = 0; det < 540; det++) {
931 if(linearfit->GetLinearFitterHisto(det)){
932 if(TMath::Abs(first)<0.0001){
933 histolinearfitsum = linearfit->GetLinearFitterHisto(det);
937 histolinearfitsum ->Add(linearfit->GetLinearFitterHisto(det));
943 ///////////////////////////
945 //////////////////////////
947 gStyle->SetPalette(1);
948 gStyle->SetOptStat(1111);
949 gStyle->SetOptFit(1111);
950 gStyle->SetPadBorderMode(0);
951 gStyle->SetCanvasColor(10);
952 gStyle->SetPadLeftMargin(0.13);
953 gStyle->SetPadRightMargin(0.13);
955 /////////////////////////
957 ////////////////////////
961 TCanvas *debugEvents = new TCanvas("cNEvents","cNEvents",10,10,510,510);
963 if(nEvents) nEvents->Draw();
969 TCanvas *debugAbsoluteGain = new TCanvas("cAbsoluteGain","cAbsoluteGain",10,10,510,510);
970 debugAbsoluteGain->cd(1);
971 if(absoluteGain) absoluteGain->Draw();
975 if(trdTrack || tpctrdTrack) {
977 TCanvas *debugtrdtpcTrack = new TCanvas("TRDtracktpctrdtrack","TRDtracktpctrdtrack",10,10,510,510);
978 debugtrdtpcTrack->Divide(2,1);
979 debugtrdtpcTrack->cd(1);
980 if(trdTrack) trdTrack->Draw();
981 if(trdTrackOffline) trdTrackOffline->Draw("same");
982 if(trdTrackStandalone) trdTrackStandalone->Draw("same");
983 TLegend *leg = new TLegend(0.4,0.6,0.89,0.89);
984 if(trdTrack) leg->AddEntry(trdTrack,"All","p");
985 if(trdTrackOffline) leg->AddEntry(trdTrackOffline,"Offline","p");
986 if(trdTrackStandalone) leg->AddEntry(trdTrackStandalone,"Standalone","p");
988 debugtrdtpcTrack->cd(2);
989 if(tpctrdTrack) tpctrdTrack->Draw();
990 TLine *line = new TLine(0.0,0.0,100.0,100.0);
995 if(nbTimeBin || nbTracklets || nbClusters) {
997 TCanvas *debugTracklets = new TCanvas("TRDtimebintrackletcluster","TRDtimebintrackletcluster",10,10,510,510);
998 debugTracklets->Divide(3,1);
999 debugTracklets->cd(1);
1000 if(nbTimeBin) nbTimeBin->Draw();
1001 if(nbTimeBinOffline) nbTimeBinOffline->Draw("same");
1002 if(nbTimeBinStandalone) nbTimeBinStandalone->Draw("same");
1003 TLegend *lega = new TLegend(0.4,0.6,0.89,0.89);
1004 if(nbTimeBin) lega->AddEntry(nbTimeBin,"All","p");
1005 if(nbTimeBinOffline) lega->AddEntry(nbTimeBinOffline,"Offline","p");
1006 if(nbTimeBinStandalone) lega->AddEntry(nbTimeBinStandalone,"Standalone","p");
1008 debugTracklets->cd(2);
1009 if(nbTracklets) nbTracklets->Draw();
1010 if(nbTrackletsOffline) nbTrackletsOffline->Draw("same");
1011 if(nbTrackletsStandalone) nbTrackletsStandalone->Draw("same");
1012 TLegend *legb = new TLegend(0.4,0.6,0.89,0.89);
1013 if(nbTracklets) legb->AddEntry(nbTracklets,"All","p");
1014 if(nbTrackletsOffline) legb->AddEntry(nbTrackletsOffline,"Offline","p");
1015 if(nbTrackletsStandalone) legb->AddEntry(nbTrackletsStandalone,"Standalone","p");
1017 debugTracklets->cd(3);
1018 if(nbClusters) nbClusters->Draw();
1019 if(nbClustersOffline) nbClustersOffline->Draw("same");
1020 if(nbClustersStandalone) nbClustersStandalone->Draw("same");
1021 TLegend *legc = new TLegend(0.4,0.6,0.89,0.89);
1022 if(nbClusters) legc->AddEntry(nbClusters,"All","p");
1023 if(nbClustersOffline) legc->AddEntry(nbClustersOffline,"Offline","p");
1024 if(nbClustersStandalone) legc->AddEntry(nbClustersStandalone,"Standalone","p");
1029 if(ch2dSum || ph2dSum || histolinearfitsum) {
1031 TCanvas *debugSum = new TCanvas("SumCalibrationObjects","SumCalibrationObjects",10,10,510,510);
1032 debugSum->Divide(3,1);
1034 if(ch2dSum) ch2dSum->Draw("lego");
1036 if(ph2dSum) ph2dSum->Draw("lego");
1038 if(histolinearfitsum) histolinearfitsum->Draw();
1042 if(ch2dSM || ph2dSM) {
1044 TCanvas *debugSM = new TCanvas("SMCalibrationObjects","SMCalibrationObjects",10,10,510,510);
1045 debugSM->Divide(2,1);
1047 if(ch2dSM) ch2dSM->Draw("lego");
1049 if(ph2dSM) ph2dSM->Draw("lego");
1055 TCanvas *debug = new TCanvas("CalibrationObjects","CalibrationObjects",10,10,510,510);
1058 if(ch2d) ch2d->Draw("lego");
1060 if(ph2d) ph2d->Draw("lego");
1065 //_______________________________________________________________________________________
1066 void AliTRDCalibTask::AddTask(const AliTRDCalibTask * calibTask) {
1072 TList *listcalibTask = calibTask->GetList();
1073 if(!listcalibTask) return;
1075 TH1I *nEvents = (TH1I *) listcalibTask->FindObject("NEvents");
1076 TH2F *absoluteGain = (TH2F *) listcalibTask->FindObject("AbsoluteGain");
1078 TH1F *trdTrack = (TH1F *) listcalibTask->FindObject("TRDTrack");
1079 TH1F *trdTrackOffline = (TH1F *) listcalibTask->FindObject("TRDTrackOffline");
1080 TH1F *trdTrackStandalone = (TH1F *) listcalibTask->FindObject("TRDTrackStandalone");
1082 TH2F *tpctrdTrack = (TH2F *) listcalibTask->FindObject("NbTPCTRDtrack");
1084 TH1F *nbTimeBin = (TH1F *) listcalibTask->FindObject("NbTimeBin");
1085 TH1F *nbTimeBinOffline = (TH1F *) listcalibTask->FindObject("NbTimeBinOffline");
1086 TH1F *nbTimeBinStandalone = (TH1F *) listcalibTask->FindObject("NbTimeBinStandalone");
1088 TH1F *nbClusters = (TH1F *) listcalibTask->FindObject("NbClusters");
1089 TH1F *nbClustersOffline = (TH1F *) listcalibTask->FindObject("NbClustersOffline");
1090 TH1F *nbClustersStandalone = (TH1F *) listcalibTask->FindObject("NbClustersStandalone");
1092 TH1F *nbTracklets = (TH1F *) listcalibTask->FindObject("NbTracklets");
1093 TH1F *nbTrackletsOffline = (TH1F *) listcalibTask->FindObject("NbTrackletsOffline");
1094 TH1F *nbTrackletsStandalone = (TH1F *) listcalibTask->FindObject("NbTrackletsStandalone");
1096 TH2I *ch2d = (TH2I *) listcalibTask->FindObject("CH2d");
1097 TProfile2D *ph2d = (TProfile2D *) listcalibTask->FindObject("PH2d");
1098 TProfile2D *prf2d = (TProfile2D *) listcalibTask->FindObject("PRF2d");
1100 TH2I *ch2dSum = (TH2I *) listcalibTask->FindObject("CH2dSum");
1101 TProfile2D *ph2dSum = (TProfile2D *) listcalibTask->FindObject("PH2dSum");
1103 TH2I *ch2dSM = (TH2I *) listcalibTask->FindObject("CH2dSM");
1104 TProfile2D *ph2dSM = (TProfile2D *) listcalibTask->FindObject("PH2dSM");
1106 AliTRDCalibraVdriftLinearFit *linearfit = (AliTRDCalibraVdriftLinearFit *) listcalibTask->FindObject("AliTRDCalibraVdriftLinearFit");
1107 AliTRDCalibraVector *calibraVector = (AliTRDCalibraVector *) listcalibTask->FindObject("AliTRDCalibraVector");
1111 TH1I *inEvents = (TH1I *) fListHist->FindObject("NEvents");
1112 TH2F *iabsoluteGain = (TH2F *) fListHist->FindObject("AbsoluteGain");
1114 TH1F *itrdTrack = (TH1F *) fListHist->FindObject("TRDTrack");
1115 TH1F *itrdTrackOffline = (TH1F *) fListHist->FindObject("TRDTrackOffline");
1116 TH1F *itrdTrackStandalone = (TH1F *) fListHist->FindObject("TRDTrackStandalone");
1118 TH2F *itpctrdTrack = (TH2F *) fListHist->FindObject("NbTPCTRDtrack");
1120 TH1F *inbTimeBin = (TH1F *) fListHist->FindObject("NbTimeBin");
1121 TH1F *inbTimeBinOffline = (TH1F *) fListHist->FindObject("NbTimeBinOffline");
1122 TH1F *inbTimeBinStandalone = (TH1F *) fListHist->FindObject("NbTimeBinStandalone");
1124 TH1F *inbClusters = (TH1F *) fListHist->FindObject("NbClusters");
1125 TH1F *inbClustersOffline = (TH1F *) fListHist->FindObject("NbClustersOffline");
1126 TH1F *inbClustersStandalone = (TH1F *) fListHist->FindObject("NbClustersStandalone");
1128 TH1F *inbTracklets = (TH1F *) fListHist->FindObject("NbTracklets");
1129 TH1F *inbTrackletsOffline = (TH1F *) fListHist->FindObject("NbTrackletsOffline");
1130 TH1F *inbTrackletsStandalone = (TH1F *) fListHist->FindObject("NbTrackletsStandalone");
1132 TH2I *ich2d = (TH2I *) fListHist->FindObject("CH2d");
1133 TProfile2D *iph2d = (TProfile2D *) fListHist->FindObject("PH2d");
1134 TProfile2D *iprf2d = (TProfile2D *) fListHist->FindObject("PRF2d");
1136 TH2I *ich2dSum = (TH2I *) fListHist->FindObject("CH2dSum");
1137 TProfile2D *iph2dSum = (TProfile2D *) fListHist->FindObject("PH2dSum");
1139 TH2I *ich2dSM = (TH2I *) fListHist->FindObject("CH2dSM");
1140 TProfile2D *iph2dSM = (TProfile2D *) fListHist->FindObject("PH2dSM");
1142 AliTRDCalibraVdriftLinearFit *ilinearfit = (AliTRDCalibraVdriftLinearFit *) fListHist->FindObject("AliTRDCalibraVdriftLinearFit");
1143 AliTRDCalibraVector *icalibraVector = (AliTRDCalibraVector *) fListHist->FindObject("AliTRDCalibraVector");
1150 inEvents->Add(nEvents);
1151 //printf("Add Events\n");
1154 //printf("Create new Events\n");
1155 inEvents = new TH1I(*nEvents);
1156 fListHist->Add(inEvents);
1161 if(iabsoluteGain) iabsoluteGain->Add(absoluteGain);
1163 iabsoluteGain = new TH2F(*absoluteGain);
1164 fListHist->Add(iabsoluteGain);
1169 if(itrdTrack) itrdTrack->Add(trdTrack);
1171 itrdTrack = new TH1F(*trdTrack);
1172 fListHist->Add(itrdTrack);
1176 if(trdTrackOffline) {
1177 if(itrdTrackOffline) itrdTrackOffline->Add(trdTrackOffline);
1179 itrdTrackOffline = new TH1F(*trdTrackOffline);
1180 fListHist->Add(itrdTrackOffline);
1184 if(trdTrackStandalone) {
1185 if(itrdTrackStandalone) itrdTrackStandalone->Add(trdTrackStandalone);
1187 itrdTrackStandalone = new TH1F(*trdTrackStandalone);
1188 fListHist->Add(itrdTrackStandalone);
1193 if(itpctrdTrack) itpctrdTrack->Add(tpctrdTrack);
1195 itpctrdTrack = new TH2F(*tpctrdTrack);
1196 fListHist->Add(itpctrdTrack);
1201 if(inbTimeBin) inbTimeBin->Add(nbTimeBin);
1203 inbTimeBin = new TH1F(*inbTimeBin);
1204 fListHist->Add(inbTimeBin);
1208 if(nbTimeBinOffline) {
1209 if(inbTimeBinOffline) inbTimeBinOffline->Add(nbTimeBinOffline);
1211 inbTimeBinOffline = new TH1F(*nbTimeBinOffline);
1212 fListHist->Add(inbTimeBinOffline);
1216 if(nbTimeBinStandalone) {
1217 if(inbTimeBinStandalone) inbTimeBinStandalone->Add(nbTimeBinStandalone);
1219 inbTimeBinStandalone = new TH1F(*nbTimeBinStandalone);
1220 fListHist->Add(inbTimeBinStandalone);
1225 if(inbClusters) inbClusters->Add(nbClusters);
1227 inbClusters = new TH1F(*nbClusters);
1228 fListHist->Add(inbClusters);
1232 if(nbClustersOffline) {
1233 if(inbClustersOffline) inbClustersOffline->Add(nbClustersOffline);
1235 inbClustersOffline = new TH1F(*nbClustersOffline);
1236 fListHist->Add(inbClustersOffline);
1240 if(nbClustersStandalone) {
1241 if(inbClustersStandalone) inbClustersStandalone->Add(nbClustersStandalone);
1243 inbClustersStandalone = new TH1F(*nbClustersStandalone);
1244 fListHist->Add(inbClustersStandalone);
1249 if(inbTracklets) inbTracklets->Add(nbTracklets);
1251 inbTracklets = new TH1F(*nbTracklets);
1252 fListHist->Add(inbTracklets);
1256 if(nbTrackletsOffline) {
1257 if(inbTrackletsOffline) inbTrackletsOffline->Add(nbTrackletsOffline);
1259 inbTrackletsOffline = new TH1F(*nbTrackletsOffline);
1260 fListHist->Add(inbTrackletsOffline);
1264 if(nbTrackletsStandalone) {
1265 if(inbTrackletsStandalone) inbTrackletsStandalone->Add(nbTrackletsStandalone);
1267 inbTrackletsStandalone = new TH1F(*nbTrackletsStandalone);
1268 fListHist->Add(inbTrackletsStandalone);
1273 if(ich2d) ich2d->Add(ch2d);
1275 ich2d = new TH2I(*ch2d);
1276 fListHist->Add(ich2d);
1281 if(iph2d) iph2d->Add(ph2d);
1283 iph2d = new TProfile2D(*ph2d);
1284 fListHist->Add(iph2d);
1289 if(iprf2d) iprf2d->Add(prf2d);
1291 iprf2d = new TProfile2D(*prf2d);
1292 fListHist->Add(iprf2d);
1297 if(ich2dSum) ich2dSum->Add(ch2dSum);
1299 ich2dSum = new TH2I(*ch2dSum);
1300 fListHist->Add(ich2dSum);
1305 if(iph2dSum) iph2dSum->Add(ph2dSum);
1307 iph2dSum = new TProfile2D(*ph2dSum);
1308 fListHist->Add(iph2dSum);
1313 if(ich2dSM) ich2dSM->Add(ch2dSM);
1315 ich2dSM = new TH2I(*ch2dSM);
1316 fListHist->Add(ich2dSM);
1321 if(iph2dSM) iph2dSM->Add(ph2dSM);
1323 iph2dSM = new TProfile2D(*ph2dSM);
1324 fListHist->Add(iph2dSM);
1329 if(ilinearfit) ilinearfit->Add(linearfit);
1331 ilinearfit = new AliTRDCalibraVdriftLinearFit(*linearfit);
1332 fListHist->Add(ilinearfit);
1337 if(icalibraVector) icalibraVector->Add(calibraVector);
1339 icalibraVector = new AliTRDCalibraVector(*calibraVector);
1340 fListHist->Add(icalibraVector);
1345 //________________________________________________________________________________
1346 Long64_t AliTRDCalibTask::Merge(TCollection *li) {
1352 TIterator* iter = li->MakeIterator();
1353 AliTRDCalibTask* cal = 0;
1355 while ((cal = (AliTRDCalibTask*)iter->Next())) {
1356 if (!cal->InheritsFrom(AliTRDCalibTask::Class())) {
1357 Error("Merge","Attempt to add object of class %s to a %s", cal->ClassName(), this->ClassName());
1361 // add histograms here...
1369 //_____________________________________________________
1370 Bool_t AliTRDCalibTask::SetVersionSubversion(){
1372 // Load Chamber Gain factors into the Tender supply
1375 printf("SetVersionSubversion\n");
1377 //find previous entry from the UserInfo
1378 TTree *tree=((TChain*)GetInputData(0))->GetTree();
1380 AliError("Tree not found in ESDhandler");
1384 TList *userInfo=(TList*)tree->GetUserInfo();
1386 AliError("No UserInfo found in tree");
1390 TList *cdbList=(TList*)userInfo->FindObject("cdbList");
1392 AliError("No cdbList found in UserInfo");
1393 if (AliLog::GetGlobalLogLevel()>=AliLog::kError) userInfo->Print();
1397 TIter nextCDB(cdbList);
1399 while ( (os=(TObjString*)nextCDB()) ){
1400 if(os->GetString().Contains("TRD/Calib/ChamberGainFactor")){
1401 // Get Old gain calibration
1402 AliCDBId *id=AliCDBId::MakeFromString(os->GetString());
1403 fVersionGainUsed = id->GetVersion();
1404 fSubVersionGainUsed = id->GetSubVersion();
1405 } else if(os->GetString().Contains("TRD/Calib/ChamberVdrift")){
1406 // Get Old drift velocity calibration
1407 AliCDBId *id=AliCDBId::MakeFromString(os->GetString());
1408 fVersionVdriftUsed = id->GetVersion();
1409 fSubVersionVdriftUsed = id->GetSubVersion();
1410 } else if(os->GetString().Contains("TRD/Calib/LocalGainFactor")){
1411 // Get Old drift velocity calibration
1412 AliCDBId *id=AliCDBId::MakeFromString(os->GetString());
1413 fVersionGainLocalUsed = id->GetVersion();
1414 fSubVersionGainLocalUsed = id->GetSubVersion();
1418 //printf("VersionGain %d, SubversionGain %d, VersionLocalGain %d, Subversionlocalgain %d, Versionvdrift %d, Subversionvdrift %d\n",fVersionGainUsed,fSubVersionGainUsed,fVersionGainLocalUsed,fSubVersionGainLocalUsed,fVersionVdriftUsed,fSubVersionVdriftUsed);
1421 if((fVersionGainUsed < 0) || (fVersionGainLocalUsed < 0) || (fSubVersionGainUsed < 0) || (fSubVersionGainLocalUsed < 0) || (fVersionVdriftUsed < 0) || (fSubVersionVdriftUsed < 0)) {
1422 AliError("No recent calibration found");