CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_9_patch3/src/RecoLocalCalo/EcalRecAlgos/src/ESRecHitFitAlgo.cc

Go to the documentation of this file.
00001 #include "RecoLocalCalo/EcalRecAlgos/interface/ESRecHitFitAlgo.h"
00002 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00003 
00004 #include "TMath.h"
00005 #include "TGraph.h"
00006 
00007 #include <iostream>
00008 #include <cmath>
00009 
00010 // fit function
00011 Double_t fitf(Double_t *x, Double_t *par) {
00012 
00013   Double_t wc = 0.07291;
00014   Double_t n  = 1.798; // n-1 (in fact)
00015   Double_t v1 = pow(wc/n*(x[0]-par[1]), n);
00016   Double_t v2 = TMath::Exp(n-wc*(x[0]-par[1]));
00017   Double_t v  = par[0]*v1*v2;
00018 
00019   if (x[0] < par[1]) v = 0;
00020 
00021   return v;
00022 }
00023 
00024 ESRecHitFitAlgo::ESRecHitFitAlgo() {
00025 
00026   fit_ = new TF1("fit", fitf, -200, 200, 2);
00027   fit_->SetParameters(50, 10);
00028 
00029 }
00030 
00031 ESRecHitFitAlgo::~ESRecHitFitAlgo() {
00032   delete fit_;
00033 }
00034 
00035 double* ESRecHitFitAlgo::EvalAmplitude(const ESDataFrame& digi, double ped) const {
00036   
00037   double *fitresults = new double[3];
00038   double adc[3];  
00039   double tx[3];
00040   tx[0] = -5.; 
00041   tx[1] = 20.;
00042   tx[2] = 45.;
00043 
00044   for (int i=0; i<digi.size(); i++) 
00045     adc[i] = digi.sample(i).adc() - ped;
00046 
00047   double status = 0;
00048   if (adc[0] > 20) status = 14;
00049   if (adc[1] < 0 || adc[2] < 0) status = 10;
00050   if (adc[0] > adc[1] && adc[0] > adc[2]) status = 8;
00051   if (adc[2] > adc[1] && adc[2] > adc[0]) status = 9;  
00052   double r12 = (adc[1] != 0) ? adc[0]/adc[1] : 99999;
00053   double r23 = (adc[2] != 0) ? adc[1]/adc[2] : 99999;
00054   if (r12 > ratioCuts_->getR12High()) status = 5;
00055   if (r23 > ratioCuts_->getR23High()) status = 6;
00056   if (r23 < ratioCuts_->getR23Low()) status = 7;
00057 
00058   if (int(status) == 0) {
00059     double para[10];
00060     TGraph *gr = new TGraph(3, tx, adc);
00061     fit_->SetParameters(50, 10);
00062     gr->Fit("fit", "MQ");
00063     fit_->GetParameters(para); 
00064     fitresults[0] = para[0]; 
00065     fitresults[1] = para[1];  
00066     delete gr;
00067 
00068     if (adc[1] > 2800 && adc[2] > 2800) status = 11;
00069     else if (adc[1] > 2800) status = 12;
00070     else if (adc[2] > 2800) status = 13;
00071 
00072   } else {
00073     fitresults[0] = 0;
00074     fitresults[1] = -999;
00075   }
00076 
00077   fitresults[2] = status;
00078 
00079   return fitresults;
00080 }
00081 
00082 EcalRecHit ESRecHitFitAlgo::reconstruct(const ESDataFrame& digi) const {
00083   
00084   ESPedestals::const_iterator it_ped = peds_->find(digi.id());
00085   
00086   ESIntercalibConstantMap::const_iterator it_mip = mips_->getMap().find(digi.id());
00087   ESAngleCorrectionFactors::const_iterator it_ang = ang_->getMap().find(digi.id());
00088 
00089   ESChannelStatusMap::const_iterator it_status = channelStatus_->getMap().find(digi.id());
00090 
00091   double* results;
00092 
00093   results = EvalAmplitude(digi, it_ped->getMean());
00094 
00095   double energy = results[0];
00096   double t0 = results[1];
00097   int status = (int) results[2];
00098   delete[] results;
00099 
00100   double mipCalib = (std::fabs(cos(*it_ang)) != 0.) ? (*it_mip)/fabs(cos(*it_ang)) : 0.;
00101   energy *= (mipCalib != 0.) ? MIPGeV_/mipCalib : 0.;
00102 
00103   LogDebug("ESRecHitFitAlgo") << "ESRecHitFitAlgo : reconstructed energy "<<energy<<" T0 "<<t0;
00104 
00105   EcalRecHit rechit(digi.id(), energy, t0);
00106 
00107   if (it_status->getStatusCode() == 1) {
00108       rechit.setFlag(EcalRecHit::kESDead);
00109   } else {
00110     if (status == 0) 
00111       rechit.setFlag(EcalRecHit::kESGood);
00112     else if (status == 5) 
00113       rechit.setFlag(EcalRecHit::kESBadRatioFor12);
00114     else if (status == 6) 
00115       rechit.setFlag(EcalRecHit::kESBadRatioFor23Upper);
00116     else if (status == 7) 
00117       rechit.setFlag(EcalRecHit::kESBadRatioFor23Lower);
00118     else if (status == 8) 
00119       rechit.setFlag(EcalRecHit::kESTS1Largest);
00120     else if (status == 9) 
00121       rechit.setFlag(EcalRecHit::kESTS3Largest);
00122     else if (status == 10) 
00123       rechit.setFlag(EcalRecHit::kESTS3Negative);
00124     else if (status == 11) 
00125       rechit.setFlag(EcalRecHit::kESSaturated);
00126     else if (status == 12) 
00127       rechit.setFlag(EcalRecHit::kESTS2Saturated);
00128     else if (status == 13) 
00129       rechit.setFlag(EcalRecHit::kESTS3Saturated);
00130     else if (status == 14) 
00131       rechit.setFlag(EcalRecHit::kESTS13Sigmas);
00132   }
00133 
00134   return rechit;
00135 }
00136