CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_6_1_2_SLHC2/src/CalibFormats/HcalObjects/src/HcalCoderDb.cc

Go to the documentation of this file.
00001 
00009 #include "CondFormats/HcalObjects/interface/HcalQIECoder.h"
00010 #include "CalibFormats/HcalObjects/interface/HcalCoderDb.h"
00011 #include "DataFormats/HcalDigi/interface/HcalUpgradeQIESample.h"
00012 
00013 HcalCoderDb::HcalCoderDb (const HcalQIECoder& fCoder, const HcalQIEShape& fShape)
00014   : mCoder (&fCoder),
00015     mShape (&fShape)
00016 {}
00017 
00018 template <class Digi> void HcalCoderDb::adc2fC_ (const Digi& df, CaloSamples& clf) const {
00019   clf=CaloSamples(df.id(),df.size());
00020   for (int i=0; i<df.size(); i++) {
00021     clf[i]=mCoder->charge (*mShape, df[i].adc (), df[i].capid ());
00022   }
00023   clf.setPresamples(df.presamples());
00024 }
00025 
00026 template <class Digi> void HcalCoderDb::fC2adc_ (const CaloSamples& clf, Digi& df, int fCapIdOffset) const {
00027   df = Digi (clf.id ());
00028   df.setSize (clf.size ());
00029   df.setPresamples (clf.presamples ());
00030   for (int i=0; i<clf.size(); i++) {
00031     int capId = (fCapIdOffset + i) % 4;
00032     df.setSample(i, HcalQIESample(mCoder->adc(*mShape, clf[i], capId), capId, 0, 0));
00033   }
00034 }
00035 
00036 template <class Digi> void HcalCoderDb::fCUpgrade2adc_ (const CaloSamples& clf, Digi& df, int fCapIdOffset) const {
00037   df = HcalUpgradeDataFrame(clf.id(), fCapIdOffset, clf.size(), 
00038                             clf.presamples());
00039   for (int i=0; i<clf.size(); ++i) {
00040     df.setSample(i, mCoder->adc(*mShape, clf[i], df.capId(i)),
00041                  0, true);
00042   }
00043 }
00044 
00045 
00046 void HcalCoderDb::adc2fC(const HBHEDataFrame& df, CaloSamples& lf) const {adc2fC_ (df, lf);}
00047 void HcalCoderDb::adc2fC(const HODataFrame& df, CaloSamples& lf) const {adc2fC_ (df, lf);}
00048 void HcalCoderDb::adc2fC(const HFDataFrame& df, CaloSamples& lf) const {adc2fC_ (df, lf);}
00049 void HcalCoderDb::adc2fC(const ZDCDataFrame& df, CaloSamples& lf) const {adc2fC_ (df, lf);}
00050 void HcalCoderDb::adc2fC(const HcalCalibDataFrame& df, CaloSamples& lf) const {adc2fC_ (df, lf);}
00051 void HcalCoderDb::adc2fC(const HcalUpgradeDataFrame& df, CaloSamples& lf) const {adc2fC_ (df, lf);}
00052 
00053 void HcalCoderDb::fC2adc(const CaloSamples& clf, HBHEDataFrame& df, int fCapIdOffset) const {fC2adc_ (clf, df, fCapIdOffset);}
00054 void HcalCoderDb::fC2adc(const CaloSamples& clf, HFDataFrame& df, int fCapIdOffset) const {fC2adc_ (clf, df, fCapIdOffset);}
00055 void HcalCoderDb::fC2adc(const CaloSamples& clf, HODataFrame& df, int fCapIdOffset) const {fC2adc_ (clf, df, fCapIdOffset);}
00056 void HcalCoderDb::fC2adc(const CaloSamples& clf, ZDCDataFrame& df, int fCapIdOffset) const {fC2adc_ (clf, df, fCapIdOffset);}
00057 void HcalCoderDb::fC2adc(const CaloSamples& clf, HcalCalibDataFrame& df, int fCapIdOffset) const {fC2adc_ (clf, df, fCapIdOffset);}
00058 void HcalCoderDb::fC2adc(const CaloSamples& clf, HcalUpgradeDataFrame& df, int fCapIdOffset) const {fCUpgrade2adc_ (clf, df, fCapIdOffset);}
00059