/* Small macro to generate and update OCDB entries for a given run: To be used together with pother macros to create OCDB entries // How to use it... // // 1. go to the directory with the "CalibObjectsTrain1.root" file and start your aliroot // // 2. copy and paste the following commands // gSystem->Load("libSTEER"); gSystem->Load("libANALYSIS"); gSystem->Load("libSTAT"); gSystem->Load("libTPCcalib"); gSystem->AddIncludePath("-I$ALICE_ROOT/STEER"); gSystem->AddIncludePath("-I$ALICE_ROOT/TPC"); .L $ALICE_ROOT/TPC/CalibMacros/CalibTimeVdrift.C+ // // 3. setup your OCDB output path e.g: // ocdbStorage="local://"+gSystem->GetFromPipe("pwd")+"/OCDB"; // // 4. simple execution // CalibTimeVdriftGlobal("CalibObjectsTrain1.root", 110000, AliCDBRunRange::Infinity()) // // 5. try to visualize new entry // TObjArray *picArray = new TObjArray; MakeDefaultPlots(vdriftArray,picArray); */ #if !defined(__CINT__) || defined(__MAKECINT__) #include "TMap.h" #include "TGraphErrors.h" #include "AliExternalTrackParam.h" #include "TFile.h" #include "TGraph.h" #include "TMultiGraph.h" #include "TCanvas.h" #include "THnSparse.h" #include "TLegend.h" #include "TPad.h" #include "TH2D.h" #include "AliESDfriend.h" #include "AliTPCcalibTime.h" #include "AliSplineFit.h" #include "AliCDBMetaData.h" #include "AliCDBId.h" #include "AliCDBManager.h" #include "AliCDBStorage.h" #include "AliTPCcalibBase.h" #include "AliTPCcalibDB.h" #include "AliTPCcalibDButil.h" #include "AliRelAlignerKalman.h" #include "AliTPCParamSR.h" #endif TString ocdbStorage="dummy"; const Int_t kMinEntries=100; // minimal number of entries Int_t startRun=0, endRun=0; Int_t startTime=0, endTime=0; TObjArray * vdriftArray = 0; AliTPCcalibTime * timeDrift=0; // // functions: // void CalibTimeVdriftGlobal(Char_t* file="CalibObjectsTrain1.root", Int_t ustartRun=0, Int_t uendRun=AliCDBRunRange::Infinity()); void UpdateOCDBDrift(Int_t ustartRun, Int_t uendRun); void AddHistoGraphs( TObjArray * vdriftArray, AliTPCcalibTime *timeDrift, Int_t minEntries); void AddAlignmentGraphs( TObjArray * vdriftArray, AliTPCcalibTime *timeDrift); void AddLaserGraphs( TObjArray * vdriftArray, AliTPCcalibTime *timeDrift); void SetDefaultGraphDrift(TGraph *graph, Int_t color, Int_t style); TGraphErrors* FilterGraphMedianAbs(TGraphErrors * graph, Float_t cut,Double_t &medianY); TGraphErrors* FilterGraphDrift(TGraphErrors * graph, Float_t errSigmaCut, Float_t medianCutAbs); TGraphErrors * MakeGraphFilter0(THnSparse *hisN, Int_t itime, Int_t ival, Int_t minEntries, Double_t offset=0); void GetRunRange(AliTPCcalibTime* timeDrift){ // // find the fist and last run // TObjArray *hisArray =timeDrift->GetHistoDrift(); {for (Int_t i=0; iGetEntriesFast(); i++){ THnSparse* addHist=(THnSparse*)hisArray->UncheckedAt(i); if (addHist->GetEntries()Projection(3); TH1D* histoTime=addHist->Projection(0); printf("%s\t%f\t%d\t%d\n",histo->GetName(), histo->GetEntries(),histo->FindFirstBinAbove(0),histo->FindLastBinAbove(0)); if (startRun<=0){ startRun=histo->FindFirstBinAbove(0); endRun =histo->FindLastBinAbove(0); }else{ startRun=TMath::Min(histo->FindFirstBinAbove(0),startRun); endRun =TMath::Max(histo->FindLastBinAbove(0),endRun); } if (startTime==0){ startTime=histoTime->FindFirstBinAbove(0); endTime =histoTime->FindLastBinAbove(0); }else{ startTime=TMath::Min(histoTime->FindFirstBinAbove(0),startTime); endTime =TMath::Max(histoTime->FindLastBinAbove(0),endTime); } delete histo; delete histoTime; }} if (startRun<0) startRun=0; if (endRun<0) endRun=100000000; printf("Run range :\t%d-%d\n", startRun, endRun); printf("Time range :\t%d-%d\n", startTime, endTime); } void CalibTimeVdriftGlobal(Char_t* file, Int_t ustartRun, Int_t uendRun){ // // // const Int_t kMinEntries=500; // minimal number of entries // // 1. Initialization and run range setting TFile fcalib(file); timeDrift=(AliTPCcalibTime*)fcalib.Get("calibTime"); startRun=ustartRun; endRun=ustartRun; TObjArray *hisArray =timeDrift->GetHistoDrift(); GetRunRange(timeDrift); for (Int_t i=0; iGetEntriesFast(); i++){ THnSparse* addHist=(THnSparse*)hisArray->At(i); if (!addHist) continue; if (startTimeGetAxis(0)->SetRange(startTime-1,endTime+1); if (startRunGetAxis(3)->SetRange(startRun-1,endRun+1); } // // // 2. extraction of the information // vdriftArray = new TObjArray(); AddAlignmentGraphs(vdriftArray,timeDrift); AddHistoGraphs(vdriftArray,timeDrift,kMinEntries); AddLaserGraphs(vdriftArray,timeDrift); // // // 3. update of OCDB // // UpdateOCDBDrift(ustartRun,uendRun); } void UpdateOCDBDrift( Int_t ustartRun, Int_t uendRun){ // // Update OCDB // AliCDBMetaData *metaData= new AliCDBMetaData(); metaData->SetObjectClassName("TObjArray"); metaData->SetResponsible("Marian Ivanov"); metaData->SetBeamPeriod(1); metaData->SetAliRootVersion("05-25-01"); //root version metaData->SetComment("Calibration of the time dependence of the drift velocity"); AliCDBId* id1=NULL; id1=new AliCDBId("TPC/Calib/TimeDrift", ustartRun, uendRun); AliCDBStorage* gStorage = AliCDBManager::Instance()->GetStorage(ocdbStorage); gStorage->Put(vdriftArray, (*id1), metaData); } void UpdateDriftParam(AliTPCParam *param, TObjArray *arr, Int_t startRun){ // // update the OCDB entry for the nominal time0 // // // AliTPCParam * param = AliTPCcalibDB::Instance()->GetParameters(); AliTPCParam *paramNew = (AliTPCParam *)param->Clone(); TGraphErrors *grT = (TGraphErrors *)arr->FindObject("ALIGN_ITSM_TPC_T0"); Double_t deltaTcm = TMath::Median(grT->GetN(),grT->GetY()); Double_t deltaT = deltaTcm/param->GetDriftV(); paramNew->SetL1Delay(param->GetL1Delay()-deltaT); paramNew->Update(); AliCDBMetaData *metaData= new AliCDBMetaData(); metaData->SetObjectClassName("TObjArray"); metaData->SetResponsible("Marian Ivanov"); metaData->SetBeamPeriod(1); metaData->SetAliRootVersion("05-25-02"); //root version metaData->SetComment("Updated calibration of nominal time 0"); AliCDBId* id1=NULL; id1=new AliCDBId("TPC/Calib/Parameters", startRun, AliCDBRunRange::Infinity()); AliCDBStorage* gStorage = AliCDBManager::Instance()->GetStorage(ocdbStorage); gStorage->Put(param, (*id1), metaData); } void PrintArray(TObjArray *array){ // // // Int_t entries = array->GetEntries(); for (Int_t i=0; iAt(i)->GetName()); } } TGraphErrors* FilterGraphDrift(TGraphErrors * graph, Float_t errSigmaCut, Float_t medianCutAbs){ // 2 filters: // 1. filter graph - error cut errSigmaCut // 2. filter graph - medianCutAbs around median // // errSigmaCut - cut on error // medianCutAbs - cut on value around median Double_t dummy=0; // // // 1. filter graph - error cut errSigmaCut // TGraphErrors *graphF; graphF = AliTPCcalibDButil::FilterGraphMedianErr(graph,errSigmaCut,dummy); delete graph; if (!graphF) return 0; graph = AliTPCcalibDButil::FilterGraphMedianErr(graphF,errSigmaCut,dummy); delete graphF; if (!graph) return 0; // // filter graph - kMedianCutAbs around median // graphF=FilterGraphMedianAbs(graph, medianCutAbs,dummy); delete graph; if (!graphF) return 0; graph=FilterGraphMedianAbs(graphF, medianCutAbs,dummy); delete graphF; if (!graph) return 0; return graph; } TGraphErrors* FilterGraphMedianAbs(TGraphErrors * graph, Float_t cut,Double_t &medianY){ // // filter outlyer measurement // Only points around median +- cut filtered // if (!graph) return 0; Int_t kMinPoints=2; Int_t npoints0 = graph->GetN(); Int_t npoints=0; Float_t rmsY=0; Double_t *outx=new Double_t[npoints0]; Double_t *outy=new Double_t[npoints0]; Double_t *errx=new Double_t[npoints0]; Double_t *erry=new Double_t[npoints0]; // // if (npoints0GetY()[ipoint]==0) continue; if (iter>0 &&TMath::Abs(graph->GetY()[ipoint]-medianY)>cut) continue; outx[npoints] = graph->GetX()[ipoint]; outy[npoints] = graph->GetY()[ipoint]; errx[npoints] = graph->GetErrorX(ipoint); erry[npoints] = graph->GetErrorY(ipoint); npoints++; } if (npoints<=1) break; medianY =TMath::Median(npoints,outy); rmsY =TMath::RMS(npoints,outy); } TGraphErrors *graphOut=0; if (npoints>1) graphOut= new TGraphErrors(npoints,outx,outy,errx,erry); return graphOut; } void AddHistoGraphs( TObjArray * vdriftArray, AliTPCcalibTime *timeDrift, Int_t minEntries){ // // Add graphs corresponding to the alignment // const Double_t kErrSigmaCut=5; // error sigma cut - for filtering const Double_t kMedianCutAbs=0.03; // error sigma cut - for filtering // TObjArray * array=timeDrift->GetHistoDrift(); if (array){ THnSparse* hist=NULL; // 2.a) cosmics with different triggers for (Int_t i=0; iGetEntriesFast();i++){ hist=(THnSparseF*)array->UncheckedAt(i); if(!hist) continue; if (hist->GetEntries()Print(); TString name=hist->GetName(); Int_t dim[4]={0,1,2,3}; THnSparse* newHist=hist->Projection(4,dim); newHist->SetName(name); TGraphErrors* graph=AliTPCcalibBase::FitSlices(newHist,2,0,400,100,0.05,0.95, kTRUE); printf("name=%s graph=%i, N=%i\n", name.Data(), graph==0, graph->GetN()); Int_t pos=name.Index("_"); name=name(pos,name.Capacity()-pos); TString graphName=graph->ClassName(); graphName+=name; graphName.ToUpper(); // graph = FilterGraphDrift(graph, kErrSigmaCut, kMedianCutAbs); if (!graph) { printf("Graph =%s filtered out\n", name.Data()); continue; } // graph->SetMarkerStyle(i%8+20); graph->SetMarkerColor(i%7); graph->GetXaxis()->SetTitle("Time"); graph->GetYaxis()->SetTitle("v_{dcor}"); graph->SetName(graphName); graph->SetTitle(graphName); printf("Graph %d\t=\t%s\n", i, graphName.Data()); vdriftArray->Add(graph); } } } void AddAlignmentGraphs( TObjArray * vdriftArray, AliTPCcalibTime *timeDrift){ // // Add graphs corresponding to alignment to the object array // TObjArray *arrayITS=0; TObjArray *arrayTOF=0; TObjArray *arrayTRD=0; TMatrixD *mstatITS=0; TMatrixD *mstatTOF=0; TMatrixD *mstatTRD=0; // arrayITS=timeDrift->GetAlignITSTPC(); arrayTRD=timeDrift->GetAlignTRDTPC(); arrayTOF=timeDrift->GetAlignTOFTPC(); if (arrayITS->GetEntries()>0) mstatITS= AliTPCcalibDButil::MakeStatRelKalman(arrayITS,0.9,50,0.025); if (arrayTOF->GetEntries()>0) mstatTOF= AliTPCcalibDButil::MakeStatRelKalman(arrayTOF,0.9,1000,0.025); if (arrayTRD->GetEntries()>0) mstatTRD= AliTPCcalibDButil::MakeStatRelKalman(arrayTRD,0.9,50,0.025); // TObjArray * arrayITSP= AliTPCcalibDButil::SmoothRelKalman(arrayITS,*mstatITS, 0, 5.); TObjArray * arrayITSM= AliTPCcalibDButil::SmoothRelKalman(arrayITS,*mstatITS, 1, 5.); TObjArray * arrayITSB= AliTPCcalibDButil::SmoothRelKalman(arrayITSP,arrayITSM); TObjArray * arrayTOFP= AliTPCcalibDButil::SmoothRelKalman(arrayTOF,*mstatTOF, 0, 5.); TObjArray * arrayTOFM= AliTPCcalibDButil::SmoothRelKalman(arrayTOF,*mstatTOF, 1, 5.); TObjArray * arrayTOFB= AliTPCcalibDButil::SmoothRelKalman(arrayTOFP,arrayTOFM); TObjArray * arrayTRDP= 0x0; TObjArray * arrayTRDM= 0x0; TObjArray * arrayTRDB= 0x0; arrayTRDP= AliTPCcalibDButil::SmoothRelKalman(arrayTRD,*mstatTRD, 0, 5.); arrayTRDM= AliTPCcalibDButil::SmoothRelKalman(arrayTRD,*mstatTRD, 1, 5.); arrayTRDB= AliTPCcalibDButil::SmoothRelKalman(arrayTRDP,arrayTRDM); // // Int_t entries=TMath::Max(arrayITS->GetEntriesFast(),arrayTOF->GetEntriesFast()); TObjArray *arrays[12]={arrayITS, arrayITSP, arrayITSM, arrayITSB, arrayTRD, arrayTRDP, arrayTRDM, arrayTRDB, arrayTOF, arrayTOFP, arrayTOFM, arrayTOFB}; TString grnames[12]={"ALIGN_ITS", "ALIGN_ITSP", "ALIGN_ITSM", "ALIGN_ITSB", "ALIGN_TRD", "ALIGN_TRDP", "ALIGN_TRDM","ALIGN_TRDB", "ALIGN_TOF", "ALIGN_TOFP", "ALIGN_TOFM","ALIGN_TOFB"}; TString grpar[9]={"DELTAPSI", "DELTATHETA", "DELTAPHI", "DELTAX", "DELTAY", "DELTAZ", "DRIFTVD", "T0", "VDGY"}; TVectorD vX(entries); TVectorD vY(entries); TVectorD vEx(entries); TVectorD vEy(entries); TObjArray *arr=0; for (Int_t iarray=0; iarray<12; iarray++){ arr = arrays[iarray]; if (arr==0) continue; for (Int_t ipar=0; ipar<9; ipar++){ Int_t counter=0; for (Int_t itime=0; itimeGetEntriesFast(); itime++){ AliRelAlignerKalman * kalman = (AliRelAlignerKalman *) arr->UncheckedAt(itime); if (!kalman) continue; vX[counter]=kalman->GetTimeStamp(); vY[counter]=(*(kalman->GetState()))[ipar]; if (ipar==6) vY[counter]=1./(*(kalman->GetState()))[ipar]-1; vEx[counter]=0; vEy[counter]=TMath::Sqrt((*(kalman->GetStateCov()))(ipar,ipar)); counter++; } TGraphErrors * graph=new TGraphErrors(counter, vX.GetMatrixArray(), vY.GetMatrixArray(), vEx.GetMatrixArray(), vEy.GetMatrixArray()); TString grName=grnames[iarray]; grName+="_TPC_"; grName+=grpar[ipar]; graph->SetName(grName.Data()); vdriftArray->AddLast(graph); } } } void AddLaserGraphs( TObjArray * vdriftArray, AliTPCcalibTime *timeDrift){ // // add graphs for laser // const Double_t delayL0L1 = 0.07; //this is hack for 1/2 weeks THnSparse *hisN=0; TGraphErrors *grLaser[6]={0,0,0,0,0,0}; hisN = timeDrift->GetHistVdriftLaserA(0); if (timeDrift->GetHistVdriftLaserA(0)){ grLaser[0]=MakeGraphFilter0(timeDrift->GetHistVdriftLaserA(0),0,2,5,delayL0L1); grLaser[0]->SetName("GRAPH_MEAN_DELAY_LASER_ALL_A"); vdriftArray->AddLast(grLaser[0]); } if (timeDrift->GetHistVdriftLaserA(1)){ grLaser[1]=MakeGraphFilter0(timeDrift->GetHistVdriftLaserA(1),0,2,5); grLaser[1]->SetName("GRAPH_MEAN_DRIFT_LASER_ALL_A"); vdriftArray->AddLast(grLaser[1]); } if (timeDrift->GetHistVdriftLaserA(2)){ grLaser[2]=MakeGraphFilter0(timeDrift->GetHistVdriftLaserA(2),0,2,5); grLaser[2]->SetName("GRAPH_MEAN_GLOBALYGRADIENT_LASER_ALL_A"); vdriftArray->AddLast(grLaser[2]); } if (timeDrift->GetHistVdriftLaserC(0)){ grLaser[3]=MakeGraphFilter0(timeDrift->GetHistVdriftLaserC(0),0,2,5,delayL0L1); grLaser[3]->SetName("GRAPH_MEAN_DELAY_LASER_ALL_C"); vdriftArray->AddLast(grLaser[3]); } if (timeDrift->GetHistVdriftLaserC(1)){ grLaser[4]=MakeGraphFilter0(timeDrift->GetHistVdriftLaserC(1),0,2,5); grLaser[4]->SetName("GRAPH_MEAN_DRIFT_LASER_ALL_C"); vdriftArray->AddLast(grLaser[4]); } if (timeDrift->GetHistVdriftLaserC(2)){ grLaser[5]=MakeGraphFilter0(timeDrift->GetHistVdriftLaserC(2),0,2,5); grLaser[5]->SetName("GRAPH_MEAN_GLOBALYGRADIENT_LASER_ALL_C"); vdriftArray->AddLast(grLaser[5]); } for (Int_t i=0; i<6;i++){ if (grLaser[i]) { SetDefaultGraphDrift(grLaser[i], 1,(i+20)); grLaser[i]->GetYaxis()->SetTitle("Laser Correction"); } } } TGraphErrors * MakeGraphFilter0(THnSparse *hisN, Int_t itime, Int_t ival, Int_t minEntries, Double_t offset){ // // Make graph with mean values and rms // hisN->GetAxis(itime)->SetRange(0,100000000); hisN->GetAxis(ival)->SetRange(0,100000000); TH1 * hisT = hisN->Projection(itime); TH1 * hisV = hisN->Projection(ival); // Int_t firstBinA = hisT->FindFirstBinAbove(2); Int_t lastBinA = hisT->FindLastBinAbove(2); Int_t firstBinV = hisV->FindFirstBinAbove(0); Int_t lastBinV = hisV->FindLastBinAbove(0); hisN->GetAxis(itime)->SetRange(firstBinA,lastBinA); hisN->GetAxis(ival)->SetRange(firstBinV,lastBinV); Int_t entries=0; for (Int_t ibin=firstBinA; ibinGetBinContent(ibin); if (contGetBinContent(ibin); if (contGetAxis(itime)->SetRange(ibin-1,ibin+1); Double_t time = hisT->GetBinCenter(ibin); TH1 * his = hisN->Projection(ival); Double_t nentries0= his->GetBinContent(his->FindBin(0)); if (cont-nentries0SetBinContent(his->FindBin(0),0); vecTime[entries]=time; vecMean0[entries]=his->GetMean()+offset; vecMean1[entries]=his->GetMeanError(); vecRMS0[entries] =his->GetRMS(); vecRMS1[entries] =his->GetRMSError(); delete his; entries++; }} delete hisT; delete hisV; TGraphErrors * graph = new TGraphErrors(entries,vecTime.GetMatrixArray(), vecMean0.GetMatrixArray(), 0, vecMean1.GetMatrixArray()); return graph; } void SetDefaultGraphDrift(TGraph *graph, Int_t color, Int_t style){ // // // graph->GetXaxis()->SetTimeDisplay(kTRUE); graph->GetXaxis()->SetTimeFormat("#splitline{%d/%m}{%H:%M}"); graph->SetMaximum( 0.025); graph->SetMinimum(-0.025); graph->GetXaxis()->SetTitle("Time"); graph->GetYaxis()->SetTitle("v_{dcorr}"); // graph->GetYaxis()->SetLabelSize(0.03); graph->GetXaxis()->SetLabelSize(0.03); // graph->GetXaxis()->SetNdivisions(10,5,0); graph->GetYaxis()->SetNdivisions(10,5,0); // graph->GetXaxis()->SetLabelOffset(0.02); graph->GetYaxis()->SetLabelOffset(0.005); // graph->GetXaxis()->SetTitleOffset(1.3); graph->GetYaxis()->SetTitleOffset(1.2); // graph->SetMarkerColor(color); graph->SetLineColor(color); graph->SetMarkerStyle(style); } void SetPadStyle(TPad *pad, Float_t mx0, Float_t mx1, Float_t my0, Float_t my1){ // // pad->SetTicks(1,1); pad->SetMargin(mx0,mx1,my0,my1); } void MakeDefaultPlots(TObjArray * arr, TObjArray *picArray){ // // // // margins Float_t mx0=0.12, mx1=0.1, my0=0.15, my1=0.1; // TGraphErrors* laserA =(TGraphErrors*)arr->FindObject("GRAPH_MEAN_DRIFT_LASER_ALL_A"); TGraphErrors* laserC =(TGraphErrors*)arr->FindObject("GRAPH_MEAN_DRIFT_LASER_ALL_C"); TGraphErrors* cosmic =(TGraphErrors*)arr->FindObject("TGRAPHERRORS_MEAN_VDRIFT_COSMICS_ALL"); TGraphErrors* cross =(TGraphErrors*)arr->FindObject("TGRAPHERRORS_VDRIFT_CROSS_ALL"); TGraphErrors* itstpcP =(TGraphErrors*)arr->FindObject("ALIGN_ITSP_TPC_DRIFTVD"); TGraphErrors* itstpcM =(TGraphErrors*)arr->FindObject("ALIGN_ITSM_TPC_DRIFTVD"); TGraphErrors* itstpcB =(TGraphErrors*)arr->FindObject("ALIGN_ITSB_TPC_DRIFTVD"); // if (laserA) SetDefaultGraphDrift(laserA,2,25); if (laserC) SetDefaultGraphDrift(laserC,4,26); if (cosmic) SetDefaultGraphDrift(cosmic,3,27); if (cross) SetDefaultGraphDrift(cross,4,28); if (itstpcP) SetDefaultGraphDrift(itstpcP,2,29); if (itstpcM) SetDefaultGraphDrift(itstpcM,4,30); if (itstpcB) SetDefaultGraphDrift(itstpcB,1,31); // // TPad *pad=0; // // Laser-Laser // if (laserA&&laserC){ pad = new TCanvas("TPCLaserVDrift","TPCLaserVDrift"); laserA->Draw("alp"); SetPadStyle(pad,mx0,mx1,my0,my1); laserA->Draw("apl"); laserC->Draw("p"); TLegend *legend = new TLegend(mx0+0.01,1-my1-0.2, 0.5, 1-my1-0.01, "Drift velocity correction"); legend->AddEntry(laserA,"Laser A side"); legend->AddEntry(laserC,"Laser C side"); legend->Draw(); picArray->AddLast(pad); } if (itstpcP&&itstpcM){ pad = new TCanvas("ITSTPC","ITSTPC"); itstpcP->Draw("alp"); SetPadStyle(pad,mx0,mx1,my0,my1); itstpcP->Draw("alp"); gPad->Clear(); itstpcM->Draw("apl"); itstpcP->Draw("p"); itstpcB->Draw("p"); TLegend *legend = new TLegend(mx0+0.01,1-my1-0.2, 0.5, 1-my1-0.01, "Drift velocity correction"); legend->AddEntry(itstpcP,"ITS-TPC smooth plus"); legend->AddEntry(itstpcM,"ITS-TPC smooth minus"); legend->AddEntry(itstpcB,"ITS-TPC smooth "); legend->Draw(); picArray->AddLast(pad); } if (itstpcB&&laserA){ pad = new TCanvas("ITSTPC_LASER","ITSTPC_LASER"); SetPadStyle(pad,mx0,mx1,my0,my1); laserA->Draw("alp"); itstpcP->Draw("p"); itstpcM->Draw("p"); itstpcB->Draw("p"); TLegend *legend = new TLegend(mx0+0.01,1-my1-0.2, 0.5, 1-my1-0.01, "Drift velocity correction"); legend->AddEntry(laserA,"TPC laser"); legend->AddEntry(itstpcP,"ITS-TPC smooth plus"); legend->AddEntry(itstpcM,"ITS-TPC smooth minus"); legend->AddEntry(itstpcB,"ITS-TPC smooth "); legend->Draw(); picArray->AddLast(pad); } if (itstpcP&&cross){ pad = new TCanvas("ITSTPC_CROSS","ITSTPC_CROSS"); SetPadStyle(pad,mx0,mx1,my0,my1); itstpcP->Draw("alp"); pad->Clear(); cross->Draw("ap"); itstpcP->Draw("p"); // TLegend *legend = new TLegend(mx0+0.01,1-my1-0.2, 0.5, 1-my1-0.01, "Drift velocity correction"); legend->AddEntry(cross,"TPC cross tracks"); legend->AddEntry(itstpcB,"ITS-TPC smooth"); legend->Draw(); picArray->AddLast(pad); } if (itstpcP&&cosmic){ pad = new TCanvas("ITSTPC_COSMIC","ITSTPC_COSMIC"); SetPadStyle(pad,mx0,mx1,my0,my1); itstpcP->Draw("alp"); pad->Clear(); cosmic->Draw("ap"); itstpcP->Draw("p"); // TLegend *legend = new TLegend(mx0+0.01,1-my1-0.2, 0.5, 1-my1-0.01, "Drift velocity correction"); legend->AddEntry(cosmic,"TPC cross tracks0 up-down"); legend->AddEntry(itstpcB,"ITS-TPC smooth"); legend->Draw(); picArray->AddLast(pad); } } TObjArray * SelectEntries(TObjArray* array, char * mask){ // // // TObjArray * arraySelect = new TObjArray; TObjArray * maskArray=0; if (mask){ TString grmask(mask); maskArray = grmask.Tokenize("*"); } for (Int_t i=0; iGetEntries();i++){ TObject *graph= array->At(i); if (!graph) continue; if (maskArray){ Bool_t isOK=kTRUE; TString str(graph->GetName()); for (Int_t imask=0; imaskGetEntries();imask++) if (str.Contains(maskArray->At(imask)->GetName())==0) isOK=kFALSE; if (!isOK) continue; if (isOK) { printf("%s\n",graph->GetName()); arraySelect->AddLast(graph->Clone()); } } } return arraySelect; }