Go to the documentation of this file.00001 #include "RecoLocalCalo/EcalRecAlgos/interface/ESRecHitSimAlgo.h"
00002 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00003
00004 #include <iostream>
00005 #include <math.h>
00006
00007 ESRecHitSimAlgo::ESRecHitSimAlgo() {
00008
00009 }
00010
00011 double* ESRecHitSimAlgo::EvalAmplitude(const ESDataFrame& digi, const double& ped, const double& w0, const double& w1, const double& w2) const {
00012
00013 double *results = new double[4];
00014 float energy = 0;
00015 double adc[3];
00016 float pw[3];
00017 pw[0] = w0;
00018 pw[1] = w1;
00019 pw[2] = w2;
00020
00021 for (int i=0; i<digi.size(); i++) {
00022 energy += pw[i]*(digi.sample(i).adc()-ped);
00023 LogDebug("ESRecHitSimAlgo") << "ESRecHitSimAlgo : Digi "<<i<<" ADC counts "<<digi.sample(i).adc()<<" Ped "<<ped;
00024
00025 adc[i] = digi.sample(i).adc() - ped;
00026 }
00027
00028 double status = 0;
00029 if (adc[0] > 20) status = 14;
00030 if (adc[1] <= 0 || adc[2] <= 0) status = 10;
00031 if (adc[0] > adc[1] && adc[0] > adc[2]) status = 8;
00032 if (adc[2] > adc[1] && adc[2] > adc[0]) status = 9;
00033 double r12 = (adc[1] != 0) ? adc[0]/adc[1] : 99;
00034 double r23 = (adc[2] != 0) ? adc[1]/adc[2] : 99;
00035 if (r12 > ratioCuts_->getR12High()) status = 5;
00036 if (r23 > ratioCuts_->getR23High()) status = 6;
00037 if (r23 < ratioCuts_->getR23Low()) status = 7;
00038
00039 double A1 = adc[1];
00040 double A2 = adc[2];
00041
00042
00043 double n = 1.798;
00044 double w = 0.07291;
00045 double DeltaT = 25.;
00046 double aaa = (A2 > 0 && A1 > 0) ? log(A2/A1)/n : 20.;
00047 double bbb = w/n*DeltaT;
00048 double ccc= exp(aaa+bbb);
00049
00050 double t0 = (2.-ccc)/(1.-ccc) * DeltaT - 5;
00051
00052
00053 double t1 = 20.;
00054 double A_1 = pow(w/n*(t1),n) * exp(n-w*(t1));
00055 double AA1 = (A_1 != 0.) ? A1 / A_1 : 0.;
00056
00057 if (adc[1] > 2800 && adc[2] > 2800) status = 11;
00058 else if (adc[1] > 2800) status = 12;
00059 else if (adc[2] > 2800) status = 13;
00060
00061 results[0] = energy;
00062 results[1] = t0;
00063 results[2] = status;
00064 results[3] = AA1;
00065
00066 return results;
00067 }
00068
00069 EcalRecHit ESRecHitSimAlgo::reconstruct(const ESDataFrame& digi) const {
00070
00071 ESPedestals::const_iterator it_ped = peds_->find(digi.id());
00072
00073 ESIntercalibConstantMap::const_iterator it_mip = mips_->getMap().find(digi.id());
00074 ESAngleCorrectionFactors::const_iterator it_ang = ang_->getMap().find(digi.id());
00075
00076 ESChannelStatusMap::const_iterator it_status = channelStatus_->getMap().find(digi.id());
00077
00078 double* results;
00079
00080 results = EvalAmplitude(digi, it_ped->getMean(), w0_, w1_, w2_);
00081
00082 double energy = results[0];
00083 double t0 = results[1];
00084 int status = (int) results[2];
00085 double otenergy = results[3] * 1000000.;
00086 delete[] results;
00087
00088 double mipCalib = (fabs(cos(*it_ang)) != 0.) ? (*it_mip)/fabs(cos(*it_ang)) : 0.;
00089 energy *= (mipCalib != 0.) ? MIPGeV_/mipCalib : 0.;
00090 otenergy *= (mipCalib != 0.) ? MIPGeV_/mipCalib : 0.;
00091
00092 LogDebug("ESRecHitSimAlgo") << "ESRecHitSimAlgo : reconstructed energy "<<energy;
00093
00094 EcalRecHit rechit(digi.id(), energy, t0);
00095 rechit.setOutOfTimeEnergy(otenergy);
00096
00097 if (it_status->getStatusCode() == 1) {
00098 rechit.setFlag(EcalRecHit::kESDead);
00099 } else {
00100 if (status == 0)
00101 rechit.setFlag(EcalRecHit::kESGood);
00102 else if (status == 5)
00103 rechit.setFlag(EcalRecHit::kESBadRatioFor12);
00104 else if (status == 6)
00105 rechit.setFlag(EcalRecHit::kESBadRatioFor23Upper);
00106 else if (status == 7)
00107 rechit.setFlag(EcalRecHit::kESBadRatioFor23Lower);
00108 else if (status == 8)
00109 rechit.setFlag(EcalRecHit::kESTS1Largest);
00110 else if (status == 9)
00111 rechit.setFlag(EcalRecHit::kESTS3Largest);
00112 else if (status == 10)
00113 rechit.setFlag(EcalRecHit::kESTS3Negative);
00114 else if (status == 11)
00115 rechit.setFlag(EcalRecHit::kESSaturated);
00116 else if (status == 12)
00117 rechit.setFlag(EcalRecHit::kESTS2Saturated);
00118 else if (status == 13)
00119 rechit.setFlag(EcalRecHit::kESTS3Saturated);
00120 else if (status == 14)
00121 rechit.setFlag(EcalRecHit::kESTS13Sigmas);
00122 }
00123
00124 return rechit;
00125 }
00126