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 Types | Private Attributes
EcalLaserDbService Class Reference

#include <EcalLaserDbService.h>

Public Member Functions

 EcalLaserDbService ()
 
 EcalLaserDbService (const edm::ParameterSet &)
 
const EcalLaserAlphasgetAlphas () const
 
const EcalLaserAPDPNRatiosgetAPDPNRatios () const
 
const EcalLaserAPDPNRatiosRefgetAPDPNRatiosRef () const
 
float getLaserCorrection (DetId const &xid, edm::Timestamp const &iTime) const
 
const EcalLinearCorrectionsgetLinearCorrections () const
 
void setAlphaData (const EcalLaserAlphas *fItem)
 
void setAPDPNData (const EcalLaserAPDPNRatios *fItem)
 
void setAPDPNRefData (const EcalLaserAPDPNRatiosRef *fItem)
 
void setLinearCorrectionsData (const EcalLinearCorrections *fItem)
 

Private Types

typedef
tbb::concurrent_unordered_set
< uint32_t > 
ErrorMapT
 

Private Attributes

ErrorMapT channelsWithInvalidCorrection_
 
const EcalLaserAlphasmAlphas_
 
const EcalLaserAPDPNRatiosmAPDPNRatios_
 
const EcalLaserAPDPNRatiosRefmAPDPNRatiosRef_
 
const EcalLinearCorrectionsmLinearCorrections_
 

Detailed Description

Definition at line 27 of file EcalLaserDbService.h.

Member Typedef Documentation

typedef tbb::concurrent_unordered_set<uint32_t> EcalLaserDbService::ErrorMapT
private

Definition at line 50 of file EcalLaserDbService.h.

Constructor & Destructor Documentation

EcalLaserDbService::EcalLaserDbService ( )

Definition at line 16 of file EcalLaserDbService.cc.

17  :
18  mAlphas_ (0),
19  mAPDPNRatiosRef_ (0),
20  mAPDPNRatios_ (0),
22  {}
const EcalLinearCorrections * mLinearCorrections_
const EcalLaserAPDPNRatiosRef * mAPDPNRatiosRef_
const EcalLaserAPDPNRatios * mAPDPNRatios_
const EcalLaserAlphas * mAlphas_
EcalLaserDbService::EcalLaserDbService ( const edm::ParameterSet )

Member Function Documentation

const EcalLaserAlphas * EcalLaserDbService::getAlphas ( ) const

Definition at line 26 of file EcalLaserDbService.cc.

References mAlphas_.

26  {
27  return mAlphas_;
28 }
const EcalLaserAlphas * mAlphas_
const EcalLaserAPDPNRatios * EcalLaserDbService::getAPDPNRatios ( ) const

Definition at line 34 of file EcalLaserDbService.cc.

References mAPDPNRatios_.

34  {
35  return mAPDPNRatios_;
36 }
const EcalLaserAPDPNRatios * mAPDPNRatios_
const EcalLaserAPDPNRatiosRef * EcalLaserDbService::getAPDPNRatiosRef ( ) const

Definition at line 30 of file EcalLaserDbService.cc.

References mAPDPNRatiosRef_.

30  {
31  return mAPDPNRatiosRef_;
32 }
const EcalLaserAPDPNRatiosRef * mAPDPNRatiosRef_
float EcalLaserDbService::getLaserCorrection ( DetId const &  xid,
edm::Timestamp const &  iTime 
) const

Definition at line 44 of file EcalLaserDbService.cc.

References alpha, channelsWithInvalidCorrection_, DetId::det(), DetId::Ecal, EcalBarrel, EcalEndcap, EcalCondObjectContainer< T >::end(), f, EcalCondObjectContainer< T >::find(), readConfiguration2004_v2_cff::getCond, EcalLaserAPDPNRatios::getLaserMap(), EcalCondObjectContainer< T >::getMap(), EcalTimeDependentCorrections::getTimeMap(), EcalLaserAPDPNRatios::getTimeMap(), EcalTimeDependentCorrections::getValueMap(), EBDetId::hashedIndex(), EEDetId::hashedIndex(), GeomDetEnumerators::isBarrel(), MEEBGeom::lmr(), MEEEGeom::lmr(), mAlphas_, mAPDPNRatios_, mAPDPNRatiosRef_, mLinearCorrections_, EcalTimeDependentCorrections::Values::p1, EcalLaserAPDPNRatios::EcalLaserAPDPNpair::p1, EcalTimeDependentCorrections::Values::p2, EcalLaserAPDPNRatios::EcalLaserAPDPNpair::p2, EcalLaserAPDPNRatios::EcalLaserAPDPNpair::p3, EcalTimeDependentCorrections::Values::p3, funct::pow(), DetId::rawId(), edm::second(), DetId::subdetId(), lumiQTWidget::t, EcalTimeDependentCorrections::Times::t1, EcalLaserAPDPNRatios::EcalLaserTimeStamp::t1, EcalLaserAPDPNRatios::EcalLaserTimeStamp::t2, EcalTimeDependentCorrections::Times::t2, EcalTimeDependentCorrections::Times::t3, EcalLaserAPDPNRatios::EcalLaserTimeStamp::t3, cond::timestamp, groupFilesInBlocks::tt, and edm::Timestamp::value().

Referenced by EcalHitResponse::findLaserConstant().

44  {
45 
46  float correctionFactor = 1.0;
47 
50  const EcalLaserAPDPNRatiosRefMap& laserRefMap = mAPDPNRatiosRef_->getMap();
51  const EcalLaserAlphaMap& laserAlphaMap = mAlphas_->getMap();
54 
57  EcalLaserAPDPNref apdpnref;
61 
62  if (xid.det()==DetId::Ecal) {
63 
64  } else {
65 
66  edm::LogError("EcalLaserDbService") << " DetId is NOT in ECAL" << endl;
67  return correctionFactor;
68  }
69 
70 
71  int iLM;
72  int xind;
73  bool isBarrel=true;
74  if (xid.subdetId()==EcalBarrel) {
75  EBDetId ebid( xid.rawId() );
76  xind = ebid.hashedIndex();
77  iLM = MEEBGeom::lmr(ebid.ieta(), ebid.iphi());
78  } else if (xid.subdetId()==EcalEndcap) {
79  isBarrel=false;
80  EEDetId eeid( xid.rawId() );
81  xind = eeid.hashedIndex();
82  // SuperCrystal coordinates
83  MEEEGeom::SuperCrysCoord iX = (eeid.ix()-1)/5 + 1;
84  MEEEGeom::SuperCrysCoord iY = (eeid.iy()-1)/5 + 1;
85  iLM = MEEEGeom::lmr(iX, iY, eeid.zside());
86  } else {
87  edm::LogError("EcalLaserDbService") << " DetId is NOT in ECAL Barrel or Endcap" << endl;
88  return correctionFactor;
89  }
90  // std::cout << " LM num ====> " << iLM << endl;
91 
92  // get alpha, apd/pn ref, apd/pn pairs and timestamps for interpolation
93 
94 
95 #ifdef VERIFY_LASER
97  if (itratio != laserRatiosMap.end()) {
98  apdpnpair = (*itratio);
99  } else {
100  edm::LogError("EcalLaserDbService") << "error with laserRatiosMap!" << endl;
101  return correctionFactor;
102  }
103 
104  if (iLM-1< (int)laserTimeMap.size()) {
105  timestamp = laserTimeMap[iLM-1];
106  } else {
107  edm::LogError("EcalLaserDbService") << "error with laserTimeMap!" << endl;
108  return correctionFactor;
109  }
110 
111  EcalLinearCorrections::EcalValueMap::const_iterator itlin = linearValueMap.find(xid);
112  if (itlin != linearValueMap.end()) {
113  linValues = (*itlin);
114  } else {
115  edm::LogError("EcalLaserDbService") << "error with linearValueMap!" << endl;
116  return correctionFactor;
117  }
118 
119  if (iLM-1< (int)linearTimeMap.size()) {
120  linTimes = linearTimeMap[iLM-1];
121  } else {
122  edm::LogError("EcalLaserDbService") << "error with laserTimeMap!" << endl;
123  return correctionFactor;
124  }
125 
126  EcalLaserAPDPNRatiosRefMap::const_iterator itref = laserRefMap.find(xid);
127  if ( itref != laserRefMap.end() ) {
128  apdpnref = (*itref);
129  } else {
130  edm::LogError("EcalLaserDbService") << "error with laserRefMap!" << endl;
131  return correctionFactor;
132  }
133 
134  EcalLaserAlphaMap::const_iterator italpha = laserAlphaMap.find(xid);
135  if ( italpha != laserAlphaMap.end() ) {
136  alpha = (*italpha);
137  } else {
138  edm::LogError("EcalLaserDbService") << "error with laserAlphaMap!" << endl;
139  return correctionFactor;
140  }
141 
142 #else
143 
144  // waiting for templated lambdas
145  auto getCond =[=](EcalFloatCondObjectContainer const & cond)->float {
146  return isBarrel ? cond.barrel(xind) : cond.endcap(xind);
147  };
148 
150  return isBarrel ? cond.barrel(xind) : cond.endcap(xind);
151  };
152 
153  auto getLinear =[=](EcalLinearCorrections::EcalValueMap const & cond)->EcalLinearCorrections::Values {
154  return isBarrel ? cond.barrel(xind) : cond.endcap(xind);
155  };
156 
157 
158  apdpnpair = getPair(laserRatiosMap);
159  linValues = getLinear(linearValueMap);
160  apdpnref = getCond(laserRefMap);
161  alpha = getCond(laserAlphaMap);
162 
163  if (iLM-1< (int)laserTimeMap.size()) {
164  timestamp = laserTimeMap[iLM-1];
165  } else {
166  edm::LogError("EcalLaserDbService") << "error with laserTimeMap!" << endl;
167  return correctionFactor;
168  }
169 
170  if (iLM-1< (int)linearTimeMap.size()) {
171  linTimes = linearTimeMap[iLM-1];
172  } else {
173  edm::LogError("EcalLaserDbService") << "error with laserTimeMap!" << endl;
174  return correctionFactor;
175  }
176 
177 
178 #endif
179 
180 
181 
182  // should implement some default in case of error...
183 
184  // should do some quality checks first
185  // ...
186 
187  // we will need to treat time differently...
188  // is time in DB same format as in MC? probably not...
189 
190  // interpolation
191 
192  edm::TimeValue_t t = iTime.value();
193  long long t_i = 0, t_f = 0;
194  float p_i = 0, p_f = 0;
195  long long lt_i = 0, lt_f = 0;
196  float lp_i = 0, lp_f = 0;
197 
198  if ( t >= timestamp.t1.value() && t < timestamp.t2.value() ) {
199  t_i = timestamp.t1.value();
200  t_f = timestamp.t2.value();
201  p_i = apdpnpair.p1;
202  p_f = apdpnpair.p2;
203  } else if ( t >= timestamp.t2.value() && t <= timestamp.t3.value() ) {
204  t_i = timestamp.t2.value();
205  t_f = timestamp.t3.value();
206  p_i = apdpnpair.p2;
207  p_f = apdpnpair.p3;
208  } else if ( t < timestamp.t1.value() ) {
209  t_i = timestamp.t1.value();
210  t_f = timestamp.t2.value();
211  p_i = apdpnpair.p1;
212  p_f = apdpnpair.p2;
213 
214  } else if ( t > timestamp.t3.value() ) {
215  t_i = timestamp.t2.value();
216  t_f = timestamp.t3.value();
217  p_i = apdpnpair.p2;
218  p_f = apdpnpair.p3;
219 
220  }
221 
222  if ( t >= linTimes.t1.value() && t < linTimes.t2.value() ) {
223  lt_i = linTimes.t1.value();
224  lt_f = linTimes.t2.value();
225  lp_i = linValues.p1;
226  lp_f = linValues.p2;
227  } else if ( t >= linTimes.t2.value() && t <= linTimes.t3.value() ) {
228  lt_i = linTimes.t2.value();
229  lt_f = linTimes.t3.value();
230  lp_i = linValues.p2;
231  lp_f = linValues.p3;
232  } else if ( t < linTimes.t1.value() ) {
233  lt_i = linTimes.t1.value();
234  lt_f = linTimes.t2.value();
235  lp_i = linValues.p1;
236  lp_f = linValues.p2;
237 
238  } else if ( t > linTimes.t3.value() ) {
239  lt_i = linTimes.t2.value();
240  lt_f = linTimes.t3.value();
241  lp_i = linValues.p2;
242  lp_f = linValues.p3;
243 
244  }
245 
246  if ( apdpnref != 0 && (t_i - t_f) != 0 && (lt_i - lt_f) != 0) {
247  long long tt = t; // never subtract two unsigned!
248  float interpolatedLaserResponse = p_i/apdpnref + float(tt-t_i)*(p_f-p_i)/(apdpnref*float(t_f-t_i));
249  float interpolatedLinearResponse = lp_i/apdpnref + float(tt-lt_i)*(lp_f-lp_i)/(apdpnref*float(lt_f-lt_i)); // FIXED BY FC
250 
251  if(interpolatedLinearResponse >2.f || interpolatedLinearResponse <0.1f)
252  interpolatedLinearResponse=1.f;
253  if ( interpolatedLaserResponse <= 0. ) {
254 
255  // print message only if it is the first time we see < 0
256  // on this detid
257  if(channelsWithInvalidCorrection_.insert(xid.rawId()).second) {
258  edm::LogError("EcalLaserDbService")
259  << "Interpolated Laser correction <0 for detid "<<xid.rawId();
260 
261  }
262 
263  return correctionFactor;
264 
265  } else {
266 
267  float interpolatedTransparencyResponse = interpolatedLaserResponse / interpolatedLinearResponse;
268 
269  correctionFactor = 1.f/( std::pow(interpolatedTransparencyResponse,alpha) *interpolatedLinearResponse );
270 
271  }
272 
273  } else {
274  edm::LogError("EcalLaserDbService")
275  << "apdpnref (" << apdpnref << ") "
276  << "or t_i-t_f (" << (t_i - t_f) << " is zero!";
277  return correctionFactor;
278  }
279 
280  return correctionFactor;
281 }
const EcalTimeMap & getTimeMap() const
float EcalLaserAlpha
float alpha
Definition: AMPTWrapper.h:95
int hashedIndex() const
get a compact index for arrays
Definition: EBDetId.h:86
const self & getMap() const
ErrorMapT channelsWithInvalidCorrection_
bool isBarrel(GeomDetEnumerators::SubDetector m)
static int lmr(EBGlobalCoord ieta, EBGlobalCoord iphi)
Definition: MEEBGeom.cc:121
float EcalLaserAPDPNref
U second(std::pair< T, U > const &p)
const EcalLinearCorrections * mLinearCorrections_
const EcalLaserAPDPNRatiosMap & getLaserMap() const
const EcalLaserAPDPNRatiosRef * mAPDPNRatiosRef_
const EcalValueMap & getValueMap() const
double f[11][100]
unsigned long long TimeValue_t
Definition: Timestamp.h:28
int SuperCrysCoord
Definition: MEEEGeom.h:22
int hashedIndex() const
Definition: EEDetId.h:182
const EcalLaserTimeStampMap & getTimeMap() const
std::vector< EcalLaserTimeStamp > EcalLaserTimeStampMap
const EcalLaserAPDPNRatios * mAPDPNRatios_
std::vector< Item >::const_iterator const_iterator
const EcalLaserAlphas * mAlphas_
const_iterator find(uint32_t rawId) const
const_iterator end() const
TimeValue_t value() const
Definition: Timestamp.h:56
static int lmr(SuperCrysCoord iX, SuperCrysCoord iY, int iz)
Definition: MEEEGeom.cc:250
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
const EcalLinearCorrections * EcalLaserDbService::getLinearCorrections ( ) const

Definition at line 38 of file EcalLaserDbService.cc.

References mLinearCorrections_.

38  {
39  return mLinearCorrections_;
40 }
const EcalLinearCorrections * mLinearCorrections_
void EcalLaserDbService::setAlphaData ( const EcalLaserAlphas fItem)
inline

Definition at line 38 of file EcalLaserDbService.h.

References mAlphas_.

38 {mAlphas_ = fItem;}
const EcalLaserAlphas * mAlphas_
void EcalLaserDbService::setAPDPNData ( const EcalLaserAPDPNRatios fItem)
inline

Definition at line 40 of file EcalLaserDbService.h.

References mAPDPNRatios_.

40 {mAPDPNRatios_ = fItem;}
const EcalLaserAPDPNRatios * mAPDPNRatios_
void EcalLaserDbService::setAPDPNRefData ( const EcalLaserAPDPNRatiosRef fItem)
inline

Definition at line 39 of file EcalLaserDbService.h.

References mAPDPNRatiosRef_.

39 {mAPDPNRatiosRef_ = fItem;}
const EcalLaserAPDPNRatiosRef * mAPDPNRatiosRef_
void EcalLaserDbService::setLinearCorrectionsData ( const EcalLinearCorrections fItem)
inline

Definition at line 41 of file EcalLaserDbService.h.

References mLinearCorrections_.

41 {mLinearCorrections_ = fItem;}
const EcalLinearCorrections * mLinearCorrections_

Member Data Documentation

ErrorMapT EcalLaserDbService::channelsWithInvalidCorrection_
mutableprivate

Definition at line 51 of file EcalLaserDbService.h.

Referenced by getLaserCorrection().

const EcalLaserAlphas* EcalLaserDbService::mAlphas_
private

Definition at line 45 of file EcalLaserDbService.h.

Referenced by getAlphas(), getLaserCorrection(), and setAlphaData().

const EcalLaserAPDPNRatios* EcalLaserDbService::mAPDPNRatios_
private

Definition at line 47 of file EcalLaserDbService.h.

Referenced by getAPDPNRatios(), getLaserCorrection(), and setAPDPNData().

const EcalLaserAPDPNRatiosRef* EcalLaserDbService::mAPDPNRatiosRef_
private

Definition at line 46 of file EcalLaserDbService.h.

Referenced by getAPDPNRatiosRef(), getLaserCorrection(), and setAPDPNRefData().

const EcalLinearCorrections* EcalLaserDbService::mLinearCorrections_
private