1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
7 * Permission to use, copy, modify and distribute this software and its *
8 * documentation strictly for non-commercial purposes is hereby granted *
9 * without fee, provided that the above copyright notice appears in all *
10 * copies and that both the copyright notice and this permission notice *
11 * appear in the supporting documentation. The authors make no claims *
12 * about the suitability of this software for any purpose. It is *
13 * provided "as is" without express or implied warranty. *
14 **************************************************************************/
16 /* $Id: AliTRDCalibraExbAltFit.cxx 46327 2011-01-10 13:29:56Z cblume $ */
18 ////////////////////////////////////////////////////////////////////////////
20 // AliTRDCalibraExbAltFit //
22 // Does the ExB calibration by applying a quadratic fit //
25 // R. Bailhache (R.Bailhache@gsi.de) //
27 ////////////////////////////////////////////////////////////////////////////
30 #include <TObjArray.h>
35 #include <TLinearFitter.h>
37 #include <TDirectory.h>
38 #include <TTreeStream.h>
39 #include <TGraphErrors.h>
43 #include "AliTRDCalibraExbAltFit.h"
45 ClassImp(AliTRDCalibraExbAltFit) /*FOLD00*/
47 //_____________________________________________________________________
48 AliTRDCalibraExbAltFit::AliTRDCalibraExbAltFit() : /*FOLD00*/
51 fFitterHistoArray(540),
59 // default constructor
62 //_____________________________________________________________________
63 AliTRDCalibraExbAltFit::AliTRDCalibraExbAltFit(const AliTRDCalibraExbAltFit &ped) : /*FOLD00*/
65 fVersion(ped.fVersion),
66 fFitterHistoArray(540),
76 for (Int_t idet = 0; idet < 540; idet++){
78 const TVectorD *vectorE = (TVectorD*)ped.fFitterEArray.UncheckedAt(idet);
79 const TVectorD *vectorP = (TVectorD*)ped.fFitterPArray.UncheckedAt(idet);
80 const TH2S *hped = (TH2S*)ped.fFitterHistoArray.UncheckedAt(idet);
82 if ( vectorE != 0x0 ) fFitterEArray.AddAt(new TVectorD(*vectorE), idet);
83 if ( vectorP != 0x0 ) fFitterPArray.AddAt(new TVectorD(*vectorP), idet);
85 TH2S *hNew = (TH2S *)hped->Clone();
86 //hNew->SetDirectory(0);
87 fFitterHistoArray.AddAt(hNew,idet);
91 //_____________________________________________________________________
92 AliTRDCalibraExbAltFit::AliTRDCalibraExbAltFit(const TObjArray &obja) : /*FOLD00*/
95 fFitterHistoArray(540),
103 // constructor from a TObjArray
105 for (Int_t idet = 0; idet < 540; idet++){
106 const TH2S *hped = (TH2S*)obja.UncheckedAt(idet);
108 TH2S *hNew = (TH2S *)hped->Clone();
109 //hNew->SetDirectory(0);
110 fFitterHistoArray.AddAt(hNew,idet);
114 //_____________________________________________________________________
115 AliTRDCalibraExbAltFit& AliTRDCalibraExbAltFit::operator = (const AliTRDCalibraExbAltFit &source)
118 // assignment operator
120 if (&source == this) return *this;
121 new (this) AliTRDCalibraExbAltFit(source);
125 //_____________________________________________________________________
126 AliTRDCalibraExbAltFit::~AliTRDCalibraExbAltFit() /*FOLD00*/
131 fFitterHistoArray.SetOwner();
132 fFitterPArray.SetOwner();
133 fFitterEArray.SetOwner();
135 fFitterHistoArray.Delete();
136 fFitterPArray.Delete();
137 fFitterEArray.Delete();
139 if ( fDebugStreamer ) delete fDebugStreamer;
142 //_____________________________________________________________________________
143 void AliTRDCalibraExbAltFit::Copy(TObject &c) const
149 AliTRDCalibraExbAltFit& target = (AliTRDCalibraExbAltFit &) c;
151 // Copy only the histos
152 for (Int_t idet = 0; idet < 540; idet++){
153 if(fFitterHistoArray.UncheckedAt(idet)){
154 TH2S *hped1 = (TH2S *)target.GetFitterHisto(idet,kTRUE);
155 //hped1->SetDirectory(0);
156 hped1->Add((const TH2S *)fFitterHistoArray.UncheckedAt(idet));
163 //_____________________________________________________________________________
164 Long64_t AliTRDCalibraExbAltFit::Merge(const TCollection* list)
166 // Merge list of objects (needed by PROOF)
174 TIterator* iter = list->MakeIterator();
177 // collection of generated histograms
179 while((obj = iter->Next()) != 0)
181 AliTRDCalibraExbAltFit* entry = dynamic_cast<AliTRDCalibraExbAltFit*>(obj);
182 if (entry == 0) continue;
184 // Copy only the histos
185 for (Int_t idet = 0; idet < 540; idet++){
186 if(entry->GetFitterHisto(idet)){
187 TH2S *hped1 = (TH2S *)GetFitterHisto(idet,kTRUE);
188 Double_t entriesa = hped1->GetEntries();
189 Double_t entriesb = ((TH2S *)entry->GetFitterHisto(idet))->GetEntries();
190 if((entriesa + entriesb) < 5*32767) hped1->Add(entry->GetFitterHisto(idet));
199 //_____________________________________________________________________
200 void AliTRDCalibraExbAltFit::Add(const AliTRDCalibraExbAltFit *ped)
208 for (Int_t idet = 0; idet < 540; idet++){
209 const TH2S *hped = (TH2S*)ped->GetFitterHistoNoForce(idet);
210 //printf("idet %d\n",idet);
213 TH2S *hped1 = (TH2S *)GetFitterHisto(idet,kTRUE);
214 Double_t entriesa = hped1->GetEntries();
215 Double_t entriesb = hped->GetEntries();
216 if((entriesa + entriesb) < 5*32767) hped1->Add(hped);
220 //______________________________________________________________________________________
221 TH2S* AliTRDCalibraExbAltFit::GetFitterHisto(Int_t detector, Bool_t force)
224 // return pointer to TH2F histo
225 // if force is true create a new histo if it doesn't exist allready
227 if ( !force || fFitterHistoArray.UncheckedAt(detector) )
228 return (TH2S*)fFitterHistoArray.UncheckedAt(detector);
230 return GetFitterHistoForce(detector);
233 //______________________________________________________________________________________
234 TH2S* AliTRDCalibraExbAltFit::GetFitterHistoForce(Int_t detector)
237 // return pointer to TH2F histo
238 // if NULL create a new histo if it doesn't exist allready
240 if (fFitterHistoArray.UncheckedAt(detector))
241 return (TH2S*)fFitterHistoArray.UncheckedAt(detector);
243 // if we are forced and TLinearFitter doesn't yes exist create it
246 TString name("LFEXB");
251 TH2S *lfdv = new TH2S((const Char_t *)name,(const Char_t *) name,
252 30, -TMath::DegToRad()*45, TMath::DegToRad()*45,
254 lfdv->SetXTitle("tan(phi_{track})");
255 lfdv->SetYTitle("rms");
256 lfdv->SetZTitle("Number of tracklets");
258 lfdv->SetDirectory(0);
260 fFitterHistoArray.AddAt(lfdv,detector);
263 //______________________________________________________________________________________
264 Bool_t AliTRDCalibraExbAltFit::GetParam(Int_t detector, TVectorD *param)
267 // return param for this detector
269 if ( fFitterPArray.UncheckedAt(detector) ){
270 const TVectorD *vectorP = (TVectorD*)fFitterPArray.UncheckedAt(detector);
271 if(!param) param = new TVectorD(vectorP->GetNoElements());
272 for(Int_t k = 0; k < vectorP->GetNoElements(); k++){
273 (*param)[k] = (*vectorP)[k];
280 //______________________________________________________________________________________
281 Bool_t AliTRDCalibraExbAltFit::GetError(Int_t detector, TVectorD *error)
284 // return error for this detector
286 if ( fFitterEArray.UncheckedAt(detector) ){
287 const TVectorD *vectorE = (TVectorD*)fFitterEArray.UncheckedAt(detector);
288 if(!error) error = new TVectorD(vectorE->GetNoElements());
289 for(Int_t k = 0; k < vectorE->GetNoElements(); k++){
290 (*error)[k] = (*vectorE)[k];
297 //______________________________________________________________________________________
298 void AliTRDCalibraExbAltFit::Update(Int_t detector, Float_t tnp, Float_t pars1)
301 // Fill the 2D histos for debugging
304 TH2S *h = ((TH2S *) GetFitterHisto(detector,kTRUE));
305 Double_t nbentries = h->GetEntries();
306 if(nbentries < 5*32767) h->Fill(tnp,pars1);
309 //____________Functions fit Online CH2d________________________________________
310 void AliTRDCalibraExbAltFit::FillPEArray()
313 // Fill fFitterPArray and fFitterEArray from inside
317 Int_t *arrayI = new Int_t[540];
318 for(Int_t k = 0; k< 540; k++){
323 for(Int_t cb = 0; cb < 540; cb++){
324 const TH2S *fitterhisto = (TH2S*)fFitterHistoArray.UncheckedAt(cb);
325 //printf("Processing the detector cb %d we find %d\n",cb, (Bool_t) fitterhisto);
327 if ( fitterhisto != 0 ){
330 TAxis *xaxis = fitterhisto->GetXaxis();
331 TAxis *yaxis = fitterhisto->GetYaxis();
332 TLinearFitter fitter = TLinearFitter(3,"pol2");
334 Double_t integral = fitterhisto->Integral();
335 //printf("Integral is %f\n",integral);
336 Bool_t securitybreaking = kFALSE;
337 if(TMath::Abs(integral-1199) < 0.00001) securitybreaking = kTRUE;
338 for(Int_t ibinx = 0; ibinx < fitterhisto->GetNbinsX(); ibinx++){
339 for(Int_t ibiny = 0; ibiny < fitterhisto->GetNbinsY(); ibiny++){
340 if(fitterhisto->GetBinContent(ibinx+1,ibiny+1)>0){
341 Double_t x = xaxis->GetBinCenter(ibinx+1);
342 Double_t y = yaxis->GetBinCenter(ibiny+1);
344 for(Int_t k = 0; k < (Int_t)fitterhisto->GetBinContent(ibinx+1,ibiny+1); k++){
345 if(!securitybreaking){
346 fitter.AddPoint(&x,y);
350 if(arrayI[cb]< 1198){
351 fitter.AddPoint(&x,y);
361 //printf("Find %d entries for the detector %d\n",arrayI[cb],cb);
365 TVectorD *par = new TVectorD(3);
366 TVectorD pare = TVectorD(2);
367 TVectorD *parE = new TVectorD(3);
369 //if((fitter.EvalRobust(0.8)==0)) {
370 //if((fitter.EvalRobust()==0)) {
371 if(((fRobustFit) && (fitter.EvalRobust(0.8)==0)) || ((!fRobustFit) && (fitter.Eval()==0))) {
372 //if((fitter.Eval()==0)) {
373 //printf("Take the param\n");
374 fitter.GetParameters(*par);
376 //fitter.GetErrors(*pare);
377 //Float_t ppointError = TMath::Sqrt(TMath::Abs(fitter.GetChisquare())/arrayI[cb]);
378 //(*parE)[0] = pare[0]*ppointError;
379 //(*parE)[1] = pare[1]*ppointError;
382 (*parE)[2] = (Double_t) arrayI[cb];
383 fFitterPArray.AddAt(par,cb);
384 fFitterEArray.AddAt(parE,cb);
389 //printf("Finish\n");
392 //delete fitterhisto;
402 //____________Functions fit Online CH2d________________________________________
403 void AliTRDCalibraExbAltFit::FillPEArray2()
406 // Fill fFitterPArray and fFitterEArray from inside
411 for(Int_t cb = 0; cb < 540; cb++){
412 TH2S *fitterhisto = (TH2S*)fFitterHistoArray.UncheckedAt(cb);
413 //printf("Processing the detector cb %d we find %d\n",cb, (Bool_t) fitterhisto);
415 if ( fitterhisto != 0 ){
416 TF1 *f1 = new TF1("f1","[0]*(x-[1])**2+[2]",-1,1);
419 TGraphErrors *gg=DrawMS(fitterhisto,nEntries);
421 // printf("N: %i\n",gg->GetN());
422 // printf("entries: %i\n",nEntries);
424 if(gg->GetN() < 20) {
430 TVectorD *par = new TVectorD(3);
431 TVectorD *parE = new TVectorD(3);
434 (*parE)[2] = (Double_t) nEntries;
435 (*par)[0] = f1->GetParameter(0)*TMath::Power(f1->GetParameter(1),2)+f1->GetParameter(2);
436 (*par)[1] = -2*f1->GetParameter(0)*f1->GetParameter(1);
437 (*par)[2] = f1->GetParameter(0);
438 fFitterPArray.AddAt(par,cb);
439 fFitterEArray.AddAt(parE,cb);
441 // if ( !fDebugStreamer ) {
443 // TDirectory *backup = gDirectory;
444 // fDebugStreamer = new TTreeSRedirector("TRDdebugCalibraFill.root");
445 // if ( backup ) backup->cd(); //we don't want to be cd'd to the debug streamer
449 // fitterhisto->Draw();
453 // double exb=f1->GetParameter(1);
454 // (* fDebugStreamer) << "FindP1TrackPHtrackletV1"<<
455 // //"snpright="<<snpright<<
457 // "h2="<<(TObject*)fitterhisto<<
469 //if(fDebugStreamer) delete fDebugStreamer;
472 //_________Helper function__________________________________________________
473 TGraphErrors* AliTRDCalibraExbAltFit::DrawMS(const TH2 *const h2, Int_t &nEntries)
480 TF1 fg("fg", "gaus", -10., 30.);
481 TGraphErrors *gp = new TGraphErrors();
483 TAxis *ax(h2->GetXaxis());
484 TAxis *ay(h2->GetYaxis());
486 for(Int_t ipt(0), jpt(1), ig(0); ipt<ax->GetNbins(); ipt++, jpt++){
487 h1 = h2->ProjectionY("py", jpt, jpt);
488 fg.SetParameter(1, h1->GetMean());
489 //Float_t x(ax->GetBinCenter(jpt));
490 Int_t n=Int_t(h1->Integral(1, h1->GetNbinsX()));
493 //Warning("drawMS()", Form("reject x[%d]=%f on n=%d", jpt, x, n));
496 h1->Fit(&fg, "WWQ0");
498 //Warning("drawMS()", Form("reject x[%d]=%f on NDF=%d", jpt, x, fg.GetNDF()));
501 if(((fg.GetParameter(1)+fg.GetParameter(2)/2)>ay->GetXmax()) || ((fg.GetParameter(1)-fg.GetParameter(2)/2)<ay->GetXmin()) || (TMath::Abs(fg.GetParameter(0))< 0.00001)) continue;
502 gp->SetPoint(ig, ax->GetBinCenter(jpt), fg.GetParameter(1));
503 gp->SetPointError(ig, 0, TMath::Sqrt(pow(fg.GetParError(1),2) + (1/pow(fg.GetParameter(0),2))));