Main Page | Namespace List | Class Hierarchy | Alphabetical List | Class List | Directories | File List | Namespace Members | Class Members | File Members

TGo4FitPeakFinder.cxx

Go to the documentation of this file.
00001 //-------------------------------------------------------------
00002 //        Go4 Release Package v3.04-01 (build 30401)
00003 //                      28-November-2008
00004 //---------------------------------------------------------------
00005 //   The GSI Online Offline Object Oriented (Go4) Project
00006 //   Experiment Data Processing at EE department, GSI
00007 //---------------------------------------------------------------
00008 //
00009 //Copyright (C) 2000- Gesellschaft f. Schwerionenforschung, GSI
00010 //                    Planckstr. 1, 64291 Darmstadt, Germany
00011 //Contact:            http://go4.gsi.de
00012 //----------------------------------------------------------------
00013 //This software can be used under the license agreements as stated
00014 //in Go4License.txt file which is part of the distribution.
00015 //----------------------------------------------------------------
00016 #include "TGo4FitPeakFinder.h"
00017 
00018 #include "TMath.h"
00019 #include "TH1D.h"
00020 #include "TMatrixD.h"
00021 #include "TVectorD.h"
00022 #include "TArrayD.h"
00023 #include "TArrayI.h"
00024 #include "TSpectrum.h"
00025 
00026 #include "TGo4Fitter.h"
00027 #include "TGo4FitData.h"
00028 #include "TGo4FitModelPolynom.h"
00029 
00030 TGo4FitPeakFinder::TGo4FitPeakFinder() : TGo4FitterAction() {
00031 
00032 }
00033 
00034 TGo4FitPeakFinder::TGo4FitPeakFinder(const char* Name,  const char* DataName, Bool_t ClearModels, Int_t PolOrder) :
00035   TGo4FitterAction(Name, "Peak finder action"),
00036   fiPeakFinderType(0), fxDataName(DataName), fbClearModels(ClearModels),
00037   fbUsePolynom(PolOrder>=0), fiPolynomOrder(PolOrder>=0 ? PolOrder : 0),
00038   fd0MinWidth(5.), fd0MaxWidth(30.), fd0MaxAmplFactor(0.2),
00039   fd1LineWidth(5.),
00040   fd2NoiseFactor(2.), fd2NoiseMinimum(5.), fi2ChannelSum(2) {
00041 }
00042 
00043 TGo4FitPeakFinder::~TGo4FitPeakFinder() {
00044 }
00045 
00046 void TGo4FitPeakFinder::SetupPolynomialBackground(Int_t PolynomOrder) {
00047    if (PolynomOrder<0) SetUsePolynom(kFALSE);
00048    else {
00049      SetUsePolynom(kTRUE);
00050      SetPolynomOrder(PolynomOrder);
00051    }
00052 }
00053 
00054 void TGo4FitPeakFinder::SetupForFirst(Double_t MaxAmplFactor, Double_t MinWidth, Double_t MaxWidth) {
00055   SetPeakFinderType(0);
00056   Set0MaxAmplFactor(MaxAmplFactor);
00057   Set0MinWidth(MinWidth);
00058   Set0MaxWidth(MaxWidth);
00059 }
00060 
00061 void TGo4FitPeakFinder::SetupForSecond(Double_t LineWidth) {
00062   SetPeakFinderType(1);
00063   Set1LineWidth(LineWidth);
00064 }
00065 
00066 void TGo4FitPeakFinder::SetupForThird(Double_t NoiseFactor, Double_t NoiseMinimum, Int_t ChannelSum) {
00067   SetPeakFinderType(2);
00068   Set2NoiseFactor(NoiseFactor);
00069   Set2NoiseMinimum(NoiseMinimum);
00070   Set2ChannelSum(ChannelSum);
00071 }
00072 
00073 void TGo4FitPeakFinder::DoAction(TGo4FitterAbstract* Fitter) {
00074   TGo4Fitter* fitter = dynamic_cast<TGo4Fitter*> (Fitter);
00075   if (fitter==0) return;
00076 
00077   TGo4FitData* data = fitter->FindData(GetDataName());
00078   if (data==0) return;
00079 
00080   if (GetClearModels())
00081     fitter->DeleteModelsAssosiatedTo(data->GetName());
00082 
00083   if (GetPeakFinderType()==0)
00084      SergeyLinevPeakFinder(fitter,
00085                            data,
00086                            GetUsePolynom() ? GetPolynomOrder() : -1,
00087                            Get0MaxAmplFactor(),
00088                            Get0MinWidth(),
00089                            Get0MaxWidth());
00090 
00091    if (GetPeakFinderType()==1)
00092      ROOTPeakFinder(fitter,
00093                     data,
00094                     GetUsePolynom() ? GetPolynomOrder() : -1,
00095                     Get1LineWidth());
00096 
00097    if (GetPeakFinderType()==2)
00098      HansEsselPeakFinder(fitter,
00099                          data,
00100                          500,
00101                          Get2ChannelSum(),
00102                          Get2NoiseFactor(),
00103                          Get2NoiseMinimum(),
00104                          GetUsePolynom() ? GetPolynomOrder() : -1);
00105 }
00106 
00107 void TGo4FitPeakFinder::Print(Option_t* option) const {
00108   TGo4FitterAction::Print(option);
00109 }
00110 /*
00111 void TGo4FitPeakFinder::ROOTPeakFinder(TGo4Fitter* fitter, TGo4FitData* data, Int_t PolynomOrder, Double_t Sigma) {
00112    if ((fitter==0) || (data==0)) return;
00113 
00114    TGo4FitDataIter* iter = data->MakeIter();
00115    if (iter==0) return;
00116 
00117    Int_t size = iter->CountPoints();
00118 
00119    if ((size<10) || (iter->ScalesSize()!=1))  { delete iter; return; }
00120 
00121    TArrayF Bins(size), Scales(size);
00122 
00123    Int_t nbin = 0;
00124    if (iter->Reset()) do {
00125       Bins[nbin] = data->GetAmplValue() * iter->Value();
00126       Scales[nbin] = iter->x();
00127       nbin++;
00128    } while (iter->Next());
00129 
00130    delete iter;
00131 
00132    if (PolynomOrder>=0)
00133      fitter->AddPolynomX(data->GetName(), "Pol", PolynomOrder);
00134 
00135    TSpectrum sp(100);
00136    sp.Search1(Bins.GetArray(), size, Sigma);
00137 
00138    for(Int_t n=0;n<sp.GetNPeaks();n++) {
00139       Double_t dindx = sp.GetPositionX()[n];
00140       Int_t left = TMath::Nint(dindx-Sigma/2.);
00141       Int_t right = TMath::Nint(dindx+Sigma/2.);
00142 
00143       Double_t pos = Scales[TMath::Nint(dindx)];
00144       Double_t width = TMath::Abs(Scales[right]-Scales[left]);
00145 
00146       fitter->AddGauss1(data->GetName(), fitter->FindNextName("Gauss",n), pos, width);
00147     }
00148 }
00149 */
00150 
00151 void TGo4FitPeakFinder::ROOTPeakFinder(TGo4Fitter* fitter, TGo4FitData* data, Int_t PolynomOrder, Double_t Sigma) {
00152    if ((fitter==0) || (data==0)) return;
00153 
00154    TGo4FitDataIter* iter = data->MakeIter();
00155    if (iter==0) return;
00156 
00157    Int_t size = iter->CountPoints();
00158 
00159    if ((size<10) || (iter->ScalesSize()!=1))  { delete iter; return; }
00160 
00161    TArrayD Bins(size), Scales(size), HScales(size+1);
00162 
00163    Int_t nbin = 0;
00164    if (iter->Reset()) do {
00165       Bins[nbin] = data->GetAmplValue() * iter->Value();
00166       Double_t pos = iter->x();
00167       Double_t width = iter->HasWidths() ? iter->Widths()[0] : 1.;
00168       Scales[nbin] = pos;
00169       HScales[nbin] = pos - width/2.;
00170       HScales[nbin+1] = pos + width/2;
00171       nbin++;
00172    } while (iter->Next());
00173 
00174    delete iter;
00175 
00176    if (PolynomOrder>=0)
00177      fitter->AddPolynomX(data->GetName(), "Pol", PolynomOrder);
00178 
00179    TH1D histo("ROOTPeakFinder", "ROOTPeakFinder spectrum", size, HScales.GetArray());
00180    histo.SetDirectory(0);
00181    for (Int_t i=0;i<size;i++)
00182       histo.SetBinContent(i+1, Bins[i]);
00183 
00184    TSpectrum sp(100);
00185    sp.Search(&histo, Sigma);
00186 //   sp.Search1(Bins.GetArray(), size, Sigma);
00187 
00188    for(Int_t n=0;n<sp.GetNPeaks();n++) {
00189       Double_t pos = sp.GetPositionX()[n];
00190       Double_t width = Sigma;
00191 
00192 
00193 //   this is old code for ROOT 3.*, may be usefull for somebody
00194 /*
00195        Double_t dindx = sp.GetPositionX()[n];
00196       Int_t left = TMath::Nint(dindx-Sigma/2.);
00197       Int_t right = TMath::Nint(dindx+Sigma/2.);
00198       if (left<0) left = 0;
00199       if (right>=size) right = size-1;
00200 
00201       Double_t pos = Scales[TMath::Nint(dindx)];
00202       Double_t width = TMath::Abs(Scales[right]-Scales[left]);
00203 */
00204       fitter->AddGauss1(data->GetName(), fitter->FindNextName("Gauss",n), pos, width);
00205     }
00206 }
00207 
00208 
00209 Double_t TGo4FitPeakFinder::CalcPolynom(const TArrayD& Coef, Double_t x) {
00210    Double_t power = 1., sum = 0.;
00211    for(Int_t n=0;n<Coef.GetSize();n++)
00212      { sum+=Coef[n]*power; power*=x; }
00213    return sum;
00214 }
00215 
00216 void TGo4FitPeakFinder::DefinePolynom(Int_t size, Double_t* bins, Double_t* scales, TArrayD& Coef,
00217                    Double_t* weight, Double_t* backgr, Char_t* use) {
00218    Int_t maxorder = Coef.GetSize()-1;
00219 
00220    TMatrixD matr(maxorder+1,maxorder+1);
00221    TVectorD vector(maxorder+1);
00222 
00223    for(Int_t order1=0;order1<=maxorder;order1++)
00224      for(Int_t order2=order1;order2<=maxorder;order2++) {
00225         Double_t sum = 0;
00226         for(Int_t i=0;i<size;i++)
00227            if((use==0) || use[i]) {
00228               Double_t zn = TMath::Power(scales[i],order1+order2);
00229               if (weight) zn*=weight[i];
00230               sum+=zn;
00231            }
00232         matr(order1,order2) = sum;
00233         matr(order2,order1) = sum;
00234      }
00235 
00236    for(Int_t order1=0;order1<=maxorder;order1++)  {
00237       Double_t sum = 0;
00238       for(Int_t i=0;i<size;i++)
00239         if((use==0) || use[i]) {
00240           Double_t zn = bins[i];
00241           if (backgr) zn-=backgr[i];
00242           zn*=TMath::Power(scales[i],order1);
00243           if (weight) zn*=weight[i];
00244           sum+=zn;
00245         }
00246       vector(order1) = sum;
00247    }
00248 
00249    matr.Invert();
00250    vector*=matr;
00251 
00252    for(Int_t order=0;order<=maxorder;order++)
00253      Coef[order]=vector(order);
00254 }
00255 
00256 void TGo4FitPeakFinder::DefinePolynomEx(Int_t size, Double_t* bins, Double_t* scales, Double_t* weight, Double_t* backgr,
00257                                         Int_t lbound, Int_t rbound, TArrayD& Coef) {
00258     TArrayC use(size); use.Reset(0);
00259     if (lbound<0) lbound = 0;
00260     if (rbound>=size) rbound=size-1;
00261     for(Int_t i=lbound;i<=rbound;i++) use[i]=1;
00262     DefinePolynom(size, bins, scales, Coef, weight, backgr, use.GetArray());
00263 }
00264 
00265 Bool_t FindValue(Double_t* bins, Double_t* backgr, Char_t* use, Int_t size, Int_t& p, Int_t dir, Double_t value) {
00266    p += dir;
00267    while ((p>=0) && (p<size) && use[p]) {
00268       Double_t zn = bins[p]-backgr[p];
00269       if (zn<value) return kTRUE;
00270       p+=dir;
00271    }
00272    return kFALSE;
00273 }
00274 
00275 void TGo4FitPeakFinder::SergeyLinevPeakFinder(TGo4Fitter* fitter,
00276                                               TGo4FitData* data,
00277                                               Int_t PolOrder,
00278                                               Double_t AmplThreshold,
00279                                               Double_t MinWidth,
00280                                               Double_t MaxWidth) {
00281    if ((fitter==0) || (data==0)) return;
00282 
00283    TGo4FitDataIter* iter = data->MakeIter();
00284    if (iter==0) return;
00285 
00286    Int_t size = iter->CountPoints();
00287 
00288    if ((size<10) || (iter->ScalesSize()!=1))  { delete iter; return; }
00289 
00290    TArrayD Bins(size), Scales(size), Weights(size), Background(size);
00291 
00292    Weights.Reset(1.); Background.Reset(0.);
00293 
00294    Double_t* bins = Bins.GetArray();
00295    Double_t* scales = Scales.GetArray();
00296    Double_t* weights = Weights.GetArray();
00297    Double_t* backgr = Background.GetArray();
00298 
00299    Int_t nbin = 0;
00300    if (iter->Reset()) do {
00301       bins[nbin] = data->GetAmplValue() * iter->Value();
00302       if (iter->StandardDeviation()>0) weights[nbin] = 1./iter->StandardDeviation();
00303       scales[nbin] = iter->x();
00304       nbin++;
00305    } while (iter->Next());
00306 
00307    delete iter;
00308 
00309    if ((PolOrder>=0) && (PolOrder<(size+1)*5)) {
00310       Int_t nsegments = (PolOrder+1)*2;
00311       TArrayC usage(size); usage.Reset(0);
00312       for (Int_t nseg=0; nseg<nsegments; nseg++) {
00313          Int_t lbound = Int_t((nseg+0.)/nsegments*size);
00314          Int_t rbound = Int_t((nseg+1.)/nsegments*size-1.);
00315          Int_t nselect = (rbound-lbound)/10;
00316          if (nselect<1) nselect = 1;
00317          Int_t nsel = 0;
00318          while (nsel++<nselect) {
00319             Int_t pmin=-1;
00320             for(Int_t i=lbound;i<=rbound;i++)
00321               if(!usage[i])
00322                 if ((pmin<0) || (bins[i]<bins[pmin])) pmin = i;
00323             usage[pmin] = 1;
00324          }
00325       }
00326 
00327       TArrayD Coef(PolOrder+1);
00328       DefinePolynom(size, bins, scales, Coef, weights, 0, usage.GetArray());
00329 
00330       for(Int_t i=0;i<size;i++)
00331          backgr[i] = CalcPolynom(Coef, scales[i]);
00332 
00333       fitter->AddPolynomX(data->GetName(), "Pol", Coef);
00334    }
00335 
00336    Double_t max=0.;
00337    for(Int_t i=0;i<size;i++) {
00338       Double_t zn = bins[i]-backgr[i];
00339       if(zn>max) max = zn;
00340    }
00341    AmplThreshold = AmplThreshold*max;
00342 
00343    TArrayC usage(size); usage.Reset(1);
00344    Char_t* use = usage.GetArray();
00345 
00346    Bool_t validline = kFALSE;
00347 
00348    do {
00349 
00350    validline = kFALSE;
00351 
00352    Int_t pmax=-1; Double_t max=0.;
00353    for(Int_t i=0;i<size;i++)
00354      if (use[i]) {
00355        Double_t zn = bins[i]-backgr[i];
00356        if((pmax<0) || (zn>max)) { pmax = i; max = zn; }
00357      }
00358 
00359    Double_t sigma = TMath::Sqrt(1./weights[pmax]);
00360 
00361    if ((max<2.*sigma) || (max<AmplThreshold)) break;
00362 
00363    Int_t lbound = pmax, rbound = pmax;
00364 
00365    Bool_t lok = FindValue(bins, backgr, use, size, lbound, -1, max*0.5);
00366    Bool_t rok = FindValue(bins, backgr, use, size, rbound, +1, max*0.5);
00367 
00368    if(lok || rok) {
00369 
00370      validline = kTRUE;
00371 
00372      Int_t w = (lok && rok) ? (rbound-lbound)/4 : ( lok ? (pmax-lbound)/2 : (rbound-pmax)/2 ) ;
00373      if (w<2) w=2;
00374 
00375      TArrayD Coef1(3), Coef2(2);
00376      DefinePolynomEx(size, bins, scales, weights, backgr, pmax-w, pmax+w, Coef1);
00377      Double_t middle = -0.5*Coef1[1]/Coef1[2];
00378      Double_t amplitude = CalcPolynom(Coef1, middle);
00379 
00380      Double_t left = 0, right = 0;
00381      if (lok) {
00382        DefinePolynomEx(size, bins, scales, weights, backgr, lbound-w, lbound+w, Coef2);
00383        left = (amplitude*.5-Coef2[0])/Coef2[1];
00384      }
00385 
00386      if (rok) {
00387        DefinePolynomEx(size, bins, scales, weights, backgr, rbound-w, rbound+w, Coef2);
00388        right = (amplitude*.5-Coef2[0])/Coef2[1];
00389      }
00390 
00391      if (!lok) left = 2*middle-right; else
00392      if (!rok) right = 2*middle-left;
00393 
00394      Double_t width = (right-left)*0.5;
00395 
00396      lbound = pmax-4*w; if (lbound<0) lbound=0;
00397      rbound = pmax+4*w; if (rbound>=size) rbound=size-1;
00398      for(Int_t i=lbound;i<=rbound;i++) use[i]=0;
00399 
00400      if(lok && rok && (width>=MinWidth) && (width<=MaxWidth) && (amplitude>AmplThreshold) && (amplitude>2*sigma))
00401         fitter->AddGauss1(data->GetName(), fitter->FindNextName("Gauss",0), middle, width, amplitude);
00402    }
00403 
00404    } while(validline);
00405 }
00406 
00407 #include "f_find_peaks.c"
00408 
00409 void TGo4FitPeakFinder::HansEsselPeakFinder(TGo4Fitter* fitter,
00410                                             TGo4FitData* data,
00411                                             Int_t MaxNumPeaks,
00412                                             Int_t ChannelSum,
00413                                             Double_t NoiseFactor,
00414                                             Double_t NoiseMinimum,
00415                                             Int_t MinimasOrder) {
00416    if ((fitter==0) || (data==0) || (MaxNumPeaks<1)) return;
00417 
00418    TGo4FitDataIter* iter = data->MakeIter();
00419    if (iter==0) return;
00420 
00421    Int_t size = iter->CountPoints();
00422 
00423    if ((size<10) || (iter->ScalesSize()!=1))  { delete iter; return; }
00424 
00425    TArrayD Bins(size), Scales(size);
00426    Double_t* bins = Bins.GetArray();
00427    Double_t* scales = Scales.GetArray();
00428 
00429    Int_t nbin = 0;
00430    if (iter->Reset()) do {
00431       bins[nbin] = data->GetAmplValue() * iter->Value();
00432       scales[nbin] = iter->x();
00433       nbin++;
00434    } while (iter->Next());
00435 
00436    delete iter;
00437 
00438 
00439    int NumPeaks = 0;
00440    TArrayI Minimas(MaxNumPeaks), Maximas(MaxNumPeaks), MaximaWidths(MaxNumPeaks);
00441    TArrayI PeaksLeft(MaxNumPeaks), PeaksRight(MaxNumPeaks), PeaksMaximum(MaxNumPeaks);
00442    TArrayD MinimaValues(MaxNumPeaks), MinimaScales(MaxNumPeaks), MaximaIntegrals(MaxNumPeaks);
00443 
00444    go4fit_find_peaks(bins,
00445                 TYPE__DOUBLE,
00446                 0,
00447                 size,
00448                 ChannelSum,
00449                 NoiseFactor,
00450                 NoiseMinimum,
00451                 0,
00452                 MaxNumPeaks,
00453                 &NumPeaks,
00454                 Minimas.GetArray(),
00455                 MinimaValues.GetArray(),
00456                 Maximas.GetArray(),
00457                 MaximaWidths.GetArray(),
00458                 MaximaIntegrals.GetArray(),
00459                 PeaksLeft.GetArray(),
00460                 PeaksMaximum.GetArray(),
00461                 PeaksRight.GetArray(),
00462                 0,
00463                 0);
00464 
00465    if (NumPeaks<=0) return;
00466 
00467    Int_t NumMinimas = (NumPeaks<MaxNumPeaks) ? NumPeaks+1 : MaxNumPeaks;
00468    for(int n=0; n<NumMinimas; n++) {
00469       MinimaScales[n] = scales[Minimas[n]];
00470    }
00471 
00472    if ((MinimasOrder>=0) && (NumMinimas>1)) {
00473      Int_t order = MinimasOrder;
00474      if (order>=NumMinimas) order = NumMinimas-1;
00475      TArrayD Coef(order+1);
00476      DefinePolynom(NumMinimas, MinimaValues.GetArray(), MinimaScales.GetArray(), Coef);
00477      fitter->AddPolynomX(data->GetName(), "Pol", Coef);
00478    }
00479 
00480    for(int n=0;n<NumPeaks;n++) {
00481       Double_t pos = scales[PeaksMaximum[n]];
00482       Double_t width = 0.;
00483       if (PeaksRight[n]-PeaksLeft[n]<1)
00484         width = scales[PeaksMaximum[n]] - scales[PeaksMaximum[n]-1];
00485       else
00486         width = TMath::Abs(scales[PeaksRight[n]] - scales[PeaksLeft[n]])/2.;
00487       Double_t ampl = bins[PeaksMaximum[n]];
00488 
00489       fitter->AddGauss1(data->GetName(), fitter->FindNextName("Gauss",0), pos, width, ampl);
00490    }
00491 }
00492 
00493 
00494 
00495 
00496 //----------------------------END OF GO4 SOURCE FILE ---------------------

Generated on Fri Nov 28 12:59:13 2008 for Go4-v3.04-1 by  doxygen 1.4.2