CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Attributes
EcalUncalibRecHitMultiFitAlgo Class Reference

#include <EcalUncalibRecHitMultiFitAlgo.h>

Public Member Functions

void disableErrorCalculation ()
 
 EcalUncalibRecHitMultiFitAlgo ()
 
EcalUncalibratedRecHit makeRecHit (const EcalDataFrame &dataFrame, const EcalPedestals::Item *aped, const EcalMGPAGainRatio *aGain, const SampleMatrix &noisecor, const FullSampleVector &fullpulse, const FullSampleMatrix &fullpulsecov, const BXVector &activeBX)
 compute rechits More...
 
void setDoPrefit (bool b)
 
void setPrefitMaxChiSq (double x)
 
 ~EcalUncalibRecHitMultiFitAlgo ()
 

Private Attributes

bool _computeErrors
 
bool _doPrefit
 
double _prefitMaxChiSq
 
PulseChiSqSNNLS _pulsefunc
 
PulseChiSqSNNLS _pulsefuncSingle
 
BXVector _singlebx
 

Detailed Description

Amplitude reconstucted by the multi-template fit

Author
J.Bendavid, E.Di Marco

Definition at line 21 of file EcalUncalibRecHitMultiFitAlgo.h.

Constructor & Destructor Documentation

EcalUncalibRecHitMultiFitAlgo::EcalUncalibRecHitMultiFitAlgo ( )
EcalUncalibRecHitMultiFitAlgo::~EcalUncalibRecHitMultiFitAlgo ( )
inline

Definition at line 27 of file EcalUncalibRecHitMultiFitAlgo.h.

27 { };

Member Function Documentation

void EcalUncalibRecHitMultiFitAlgo::disableErrorCalculation ( )
inline
EcalUncalibratedRecHit EcalUncalibRecHitMultiFitAlgo::makeRecHit ( const EcalDataFrame dataFrame,
const EcalPedestals::Item aped,
const EcalMGPAGainRatio aGain,
const SampleMatrix noisecor,
const FullSampleVector fullpulse,
const FullSampleMatrix fullpulsecov,
const BXVector activeBX 
)

compute rechits

Definition at line 23 of file EcalUncalibRecHitMultiFitAlgo.cc.

References _computeErrors, _doPrefit, _prefitMaxChiSq, _pulsefunc, _pulsefuncSingle, _singlebx, EcalMGPASample::adc(), PulseChiSqSNNLS::BXs(), PulseChiSqSNNLS::ChiSq(), PulseChiSqSNNLS::disableErrorCalculation(), PulseChiSqSNNLS::DoFit(), PulseChiSqSNNLS::Errors(), flags, EcalMGPAGainRatio::gain12Over6(), EcalMGPAGainRatio::gain6Over1(), ecalMGPA::gainId(), EcalMGPASample::gainId(), EcalDataFrame::id(), bookConverter::max, EcalDataFrame::MAXSAMPLES, EcalPedestal::mean_x1, EcalPedestal::mean_x12, EcalPedestal::mean_x6, EcalCondDBWriter_cfi::pedestal, EcalPedestal::rms_x1, EcalPedestal::rms_x12, EcalPedestal::rms_x6, compare_using_db::sample, EcalDataFrame::sample(), EcalUncalibratedRecHit::setAmplitudeError(), mps_update::status, and PulseChiSqSNNLS::X().

Referenced by EcalUncalibRecHitWorkerMultiFit::run().

23  {
24 
25  uint32_t flags = 0;
26 
27  const unsigned int nsample = EcalDataFrame::MAXSAMPLES;
28 
29  double maxamplitude = -std::numeric_limits<double>::max();
30 
31  double pedval = 0.;
32  double pedrms = 0.;
33 
34  SampleVector amplitudes;
35  for(unsigned int iSample = 0; iSample < nsample; iSample++) {
36 
37  const EcalMGPASample &sample = dataFrame.sample(iSample);
38 
39  double amplitude = 0.;
40  int gainId = sample.gainId();
41 
42  double pedestal = 0.;
43  double pederr = 0.;
44  double gainratio = 1.;
45 
46  if (gainId==0 || gainId==3) {
47  pedestal = aped->mean_x1;
48  pederr = aped->rms_x1;
49  gainratio = aGain->gain6Over1()*aGain->gain12Over6();
50  }
51  else if (gainId==1) {
52  pedestal = aped->mean_x12;
53  pederr = aped->rms_x12;
54  gainratio = 1.;
55  }
56  else if (gainId==2) {
57  pedestal = aped->mean_x6;
58  pederr = aped->rms_x6;
59  gainratio = aGain->gain12Over6();
60  }
61 
62  amplitude = ((double)(sample.adc()) - pedestal) * gainratio;
63 
64  if (gainId == 0) {
65  //saturation
66  amplitude = (4095. - pedestal) * gainratio;
67  }
68 
69  amplitudes[iSample] = amplitude;
70 
71  if (amplitude>maxamplitude) {
72  //if (iSample==5) {
73  maxamplitude = amplitude;
74  pedval = pedestal;
75  pedrms = pederr*gainratio;
76  }
77 
78  }
79 
80  double amplitude, amperr, chisq;
81  bool status = false;
82 
83  //optimized one-pulse fit for hlt
84  bool usePrefit = false;
85  if (_doPrefit) {
86  status = _pulsefuncSingle.DoFit(amplitudes,noisecor,pedrms,_singlebx,fullpulse,fullpulsecov);
87  amplitude = status ? _pulsefuncSingle.X()[0] : 0.;
88  amperr = status ? _pulsefuncSingle.Errors()[0] : 0.;
89  chisq = _pulsefuncSingle.ChiSq();
90 
91  if (chisq < _prefitMaxChiSq) {
92  usePrefit = true;
93  }
94  }
95 
96  if (!usePrefit) {
97 
99  status = _pulsefunc.DoFit(amplitudes,noisecor,pedrms,activeBX,fullpulse,fullpulsecov);
100  chisq = _pulsefunc.ChiSq();
101 
102  if (!status) {
103  edm::LogWarning("EcalUncalibRecHitMultiFitAlgo::makeRecHit") << "Failed Fit" << std::endl;
104  }
105 
106  unsigned int ipulseintime = 0;
107  for (unsigned int ipulse=0; ipulse<_pulsefunc.BXs().rows(); ++ipulse) {
108  if (_pulsefunc.BXs().coeff(ipulse)==0) {
109  ipulseintime = ipulse;
110  break;
111  }
112  }
113 
114  amplitude = status ? _pulsefunc.X()[ipulseintime] : 0.;
115  amperr = status ? _pulsefunc.Errors()[ipulseintime] : 0.;
116 
117  }
118 
119  double jitter = 0.;
120 
121  //printf("status = %i\n",int(status));
122  //printf("amplitude = %5f +- %5f, chisq = %5f\n",amplitude,amperr,chisq);
123 
124  EcalUncalibratedRecHit rh( dataFrame.id(), amplitude , pedval, jitter, chisq, flags );
125  rh.setAmplitudeError(amperr);
126 
127  if (!usePrefit) {
128  for (unsigned int ipulse=0; ipulse<_pulsefunc.BXs().rows(); ++ipulse) {
129  int bx = _pulsefunc.BXs().coeff(ipulse);
130  if (bx!=0) {
131  rh.setOutOfTimeAmplitude(bx+5, status ? _pulsefunc.X().coeff(ipulse) : 0.);
132  }
133  }
134  }
135 
136  return rh;
137 }
DetId id() const
Definition: EcalDataFrame.h:24
int gainId(sample_type sample)
get the gainId (2 bits)
void disableErrorCalculation()
const PulseVector & Errors() const
Eigen::Matrix< double, 10, 1 > SampleVector
EcalMGPASample sample(int i) const
Definition: EcalDataFrame.h:29
std::vector< Variable::Flags > flags
Definition: MVATrainer.cc:135
int gainId() const
get the gainId (2 bits)
const PulseVector & X() const
float gain6Over1() const
void setAmplitudeError(float amplitudeerror)
double ChiSq() const
float gain12Over6() const
const BXVector & BXs() const
static const int MAXSAMPLES
Definition: EcalDataFrame.h:48
bool DoFit(const SampleVector &samples, const SampleMatrix &samplecor, double pederr, const BXVector &bxs, const FullSampleVector &fullpulse, const FullSampleMatrix &fullpulsecov)
int adc() const
get the ADC sample (12 bits)
tuple status
Definition: mps_update.py:57
void EcalUncalibRecHitMultiFitAlgo::setDoPrefit ( bool  b)
inline

Definition at line 30 of file EcalUncalibRecHitMultiFitAlgo.h.

References _doPrefit, and b.

Referenced by EcalUncalibRecHitWorkerMultiFit::run().

void EcalUncalibRecHitMultiFitAlgo::setPrefitMaxChiSq ( double  x)
inline

Member Data Documentation

bool EcalUncalibRecHitMultiFitAlgo::_computeErrors
private

Definition at line 36 of file EcalUncalibRecHitMultiFitAlgo.h.

Referenced by disableErrorCalculation(), and makeRecHit().

bool EcalUncalibRecHitMultiFitAlgo::_doPrefit
private

Definition at line 37 of file EcalUncalibRecHitMultiFitAlgo.h.

Referenced by makeRecHit(), and setDoPrefit().

double EcalUncalibRecHitMultiFitAlgo::_prefitMaxChiSq
private

Definition at line 38 of file EcalUncalibRecHitMultiFitAlgo.h.

Referenced by makeRecHit(), and setPrefitMaxChiSq().

PulseChiSqSNNLS EcalUncalibRecHitMultiFitAlgo::_pulsefunc
private

Definition at line 34 of file EcalUncalibRecHitMultiFitAlgo.h.

Referenced by makeRecHit().

PulseChiSqSNNLS EcalUncalibRecHitMultiFitAlgo::_pulsefuncSingle
private

Definition at line 35 of file EcalUncalibRecHitMultiFitAlgo.h.

Referenced by EcalUncalibRecHitMultiFitAlgo(), and makeRecHit().

BXVector EcalUncalibRecHitMultiFitAlgo::_singlebx
private

Definition at line 39 of file EcalUncalibRecHitMultiFitAlgo.h.

Referenced by EcalUncalibRecHitMultiFitAlgo(), and makeRecHit().