Ignore:
Timestamp:
02/10/04 16:10:52 (21 years ago)
Author:
gaug
Message:
*** empty log message ***
File:
1 edited

Legend:

Unmodified
Added
Removed
  • trunk/MagicSoft/Mars/mcalib/MHCalibrationBlindPixel.cc

    r3075 r3081  
    8888const Double_t MHCalibrationBlindPixel::fgBlindPixelElectronicAmpError = 0.002;
    8989 
    90 const Axis_t  MHCalibrationBlindPixel::fNyquistFreq       = 1.0;
    91 const Axis_t  MHCalibrationBlindPixel::fMinFreq           = 0.;
    9290const Int_t   MHCalibrationBlindPixel::fPSDNbins          = 30;
    9391 
     
    266264  fPSDHiGain = fourier.PowerSpectrumDensity(fHiGains);
    267265
    268   Int_t entries;
    269   TArrayF *array;
    270  
    271   fHPSD = new TH1F("HBlindPixelPSD",
    272                    "Power Spectrum Density Projection Blind Pixel",
    273                    fPSDNbins,fMinFreq,fNyquistFreq);
    274 
    275   array   = fPSDHiGain;
    276   entries = array->GetSize();
    277  
    278   for (Int_t i=0;i<entries;i++)
    279     fHPSD->Fill(array->At(i));
    280 
     266  fHPSD = MH::ProjectArray(fPSDHiGain, fPSDNbins,
     267                           "PSDProjectionBlindPixel",
     268                           "Power Spectrum Density Projection HiGain Blind Pixel");
    281269  //
    282270  // First guesses for the fit (should be as close to reality as possible,
    283271  //
    284   const Double_t area_guess  = entries*10.;
    285 
    286   fPSDExpFit = new TF1("PSDExpFit","[0]*exp(-[1]*x)",0.,1.);
    287 
    288   fPSDExpFit->SetParameters(entries,10.);
    289   fPSDExpFit->SetParNames("Area","slope");
    290   fPSDExpFit->SetParLimits(0,0.,3.*area_guess);
    291   fPSDExpFit->SetParLimits(1,0.,20.);
    292   fPSDExpFit->SetRange(fMinFreq,fNyquistFreq);
     272  const Double_t xmax = fHPSD->GetXaxis()->GetXmax();
     273
     274  fPSDExpFit = new TF1("PSDExpFit","exp([0]-[1]*x)",0.,xmax);
     275
     276  const Double_t slope_guess = (TMath::Log(fHPSD->GetEntries())+1.)/xmax;
     277  const Double_t offset_guess = slope_guess*xmax;
     278
     279  fPSDExpFit->SetParameters(offset_guess, slope_guess);
     280  fPSDExpFit->SetParNames("Offset","Slope");
     281  fPSDExpFit->SetParLimits(0,offset_guess/2.,2.*offset_guess);
     282  fPSDExpFit->SetParLimits(1,slope_guess/2.,2.*slope_guess);
     283  fPSDExpFit->SetRange(0.,xmax);
    293284
    294285  fHPSD->Fit(fPSDExpFit,"RQL0");
Note: See TracChangeset for help on using the changeset viewer.