CMS 3D CMS Logo

List of all members | Public Member Functions | Private Member Functions | Private Attributes
CSCComparatorDigiFitter Class Reference

#include <CSCComparatorDigiFitter.h>

Public Member Functions

 CSCComparatorDigiFitter ()
 
void fit (const CSCDetId &ch_id, const CSCCorrelatedLCTDigi &, const CSCComparatorDigiCollection &, std::vector< float > &fit_phi_layers, std::vector< float > &fit_z_layers, float &keyRadius)
 
void setGeometry (const CSCGeometry *csc_g)
 
void setStripBits (int bits)
 
void useKeyRadius (bool useKeyRadius)
 
 ~CSCComparatorDigiFitter ()
 

Private Member Functions

void calculateSlopeIntercept (float &alpha, float &beta)
 
bool comparatorInLCTPattern (int keyStrip, int pattern, int layer, int halfStrip) const
 
float cscHalfStripWidth (const CSCDetId &id) const
 
void getComparatorDigiCoordinates (const CSCDetId &ch_id, const CSCCorrelatedLCTDigi &stub)
 
void matchingComparatorDigisLCT (const CSCDetId &ch_id, const CSCCorrelatedLCTDigi &, const CSCComparatorDigiCollection &)
 

Private Attributes

CSCComparatorDigiContainerIds compDigisIds_
 
const CSCGeometrycscGeometry_
 
const std::vector< float > degrees_ = {10.,10.,10.,10.,20.,10.,20.,10.,20.,10.}
 
std::vector< float > ephis_
 
std::vector< float > ezs_
 
int nStripBits_
 
std::vector< float > phis_
 
float radius_
 
const std::vector< int > strips_ = {48,64,80,64, 80,80,80,80,80,80}
 
bool useKeyRadius_
 
std::vector< float > zs_
 

Detailed Description

Definition at line 38 of file CSCComparatorDigiFitter.h.

Constructor & Destructor Documentation

CSCComparatorDigiFitter::CSCComparatorDigiFitter ( )
inline

Definition at line 42 of file CSCComparatorDigiFitter.h.

42 {}
CSCComparatorDigiFitter::~CSCComparatorDigiFitter ( )
inline

Definition at line 43 of file CSCComparatorDigiFitter.h.

43 {}

Member Function Documentation

void CSCComparatorDigiFitter::calculateSlopeIntercept ( float &  alpha,
float &  beta 
)
private

Definition at line 266 of file CSCComparatorDigiFitter.cc.

References delta, ephis_, mps_fire::i, phis_, and zs_.

Referenced by fit(), and useKeyRadius().

267 {
268  // if there are at least 3 hits in the chamber, do a linear fit to the
269  // comparator digi positions with the chi2 method
270  if (phis_.size()>=3) {
271 
272  float Sxx = 0, Sxy = 0, Sx = 0, Sy = 0, S = 0;
273  for (unsigned i = 0; i<phis_.size(); ++i){
274  float sigma2_inv = 1./ephis_[i]*ephis_[i];
275  Sxx += zs_[i]*zs_[i] * sigma2_inv;
276  Sxy += zs_[i]*phis_[i] * sigma2_inv;
277  Sx += zs_[i]*zs_[i] * sigma2_inv;
278  Sy += phis_[i] * sigma2_inv;
279  S += sigma2_inv;
280  }
281  float delta = S * Sxx - Sx * Sx;
282  alpha = (Sxx * Sy - Sx * Sxy) / delta;
283  beta = (S * Sxy - Sx * Sy) / delta;
284  }
285  else {
286  alpha = -99; beta= 0.0;
287  }
288 }
dbl * delta
Definition: mlp_gen.cc:36
float alpha
Definition: AMPTWrapper.h:95
bool CSCComparatorDigiFitter::comparatorInLCTPattern ( int  keyStrip,
int  pattern,
int  layer,
int  halfStrip 
) const
private

Definition at line 301 of file CSCComparatorDigiFitter.cc.

References spr::find(), or, and listBenchmarks::pattern.

Referenced by matchingComparatorDigisLCT(), and useKeyRadius().

302 {
303  // get the (sub)pattern
304  const std::vector<int>& subpat = patIndexToPatternDelta[pattern].at(layer-1);
305 
306  // due to comparator digi time extension in the CLCT processor we need to
307  // search a bigger region around the key HS. +/-1, 0 should be sufficient
308  const int halfStripDelta = halfStrip - keyStrip;
309  return ( std::find(subpat.begin(), subpat.end(), halfStripDelta+1) != subpat.end() or
310  std::find(subpat.begin(), subpat.end(), halfStripDelta) != subpat.end() or
311  std::find(subpat.begin(), subpat.end(), halfStripDelta-1) != subpat.end() );
312 }
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:20
The Signals That Services Can Subscribe To This is based on ActivityRegistry and is current per Services can connect to the signals distributed by the ActivityRegistry in order to monitor the activity of the application Each possible callback has some defined which we here list in angle e< void, edm::EventID const &, edm::Timestamp const & > We also list in braces which AR_WATCH_USING_METHOD_ is used for those or
Definition: Activities.doc:12
float CSCComparatorDigiFitter::cscHalfStripWidth ( const CSCDetId id) const
private

Definition at line 291 of file CSCComparatorDigiFitter.cc.

References degrees_, M_PI, and strips_.

Referenced by getComparatorDigiCoordinates(), and useKeyRadius().

292 {
293  // what is the chamber type?
294  int index = id.iChamberType()-1;
295 
296  // calculate the half strip width of this chamber
297  return degrees_[index] * M_PI/180. / (2. * strips_[index]);
298 }
const std::vector< int > strips_
const std::vector< float > degrees_
#define M_PI
void CSCComparatorDigiFitter::fit ( const CSCDetId ch_id,
const CSCCorrelatedLCTDigi stub,
const CSCComparatorDigiCollection hCSCComparators,
std::vector< float > &  fit_phi_layers,
std::vector< float > &  fit_z_layers,
float &  keyRadius 
)

Definition at line 204 of file CSCComparatorDigiFitter.cc.

References funct::abs(), alpha, pfBoostedDoubleSVAK8TagInfos_cfi::beta, calculateSlopeIntercept(), CSCDetId::chamber(), CSCGeometry::chamber(), cscGeometry_, CSCDetId::endcap(), CSCLayer::geometry(), getComparatorDigiCoordinates(), CSCCorrelatedLCTDigi::getFractionalStrip(), CSCCorrelatedLCTDigi::getKeyWG(), mps_fire::i, CSCGeometry::idToDet(), CSCConstants::KEY_CLCT_LAYER, CSCChamber::layer(), matchingComparatorDigisLCT(), CSCLayerGeometry::middleWireOfGroup(), normalizedPhi(), nStripBits_, CSCConstants::NUM_LAYERS, or, PV3DBase< T, PVType, FrameType >::perp(), PV3DBase< T, PVType, FrameType >::phi(), phis_, radius_, CSCDetId::ring(), CSCDetId::station(), GeomDet::surface(), Surface::toGlobal(), and useKeyRadius_.

Referenced by trackingPlots.Iteration::modules(), and useKeyRadius().

208 {
209  // clear fit results
210  fit_phi_layers.clear();
211  fit_z_layers.clear();
212  keyRadius = 0;
213 
214  // first, match the comparator digis to the LCT
215  matchingComparatorDigisLCT(ch_id, stub, hCSCComparators);
216 
217  // second, get the coordinates
218  getComparatorDigiCoordinates(ch_id, stub);
219 
220  // get radius of the stub from key layer
221  const CSCDetId key_id(ch_id.endcap(), ch_id.station(), ch_id.ring(), ch_id.chamber(), CSCConstants::KEY_CLCT_LAYER);
222  const float fractional_strip = stub.getFractionalStrip();
223 
224  const auto& cscChamber = cscGeometry_->chamber(ch_id);
225  const auto& layer_geo = cscChamber->layer(CSCConstants::KEY_CLCT_LAYER)->geometry();
226 
227  // LCT::getKeyWG() also starts from 0
228  const float wire = layer_geo->middleWireOfGroup(stub.getKeyWG() + 1);
229  const LocalPoint& csc_intersect = layer_geo->intersectionOfStripAndWire(fractional_strip, wire);
230  const GlobalPoint& csc_gp = cscGeometry_->idToDet(key_id)->surface().toGlobal(csc_intersect);
231 
232  // get radius from key layer
233  if (useKeyRadius_)
234  radius_ = radius_/phis_.size();
235  else
236  radius_ = csc_gp.perp();
237 
238  float alpha = -99., beta = 0.;
239  // do a fit to the comparator digis
241  if (phis_.size() <= 2 or std::abs(alpha)>=99){
242  alpha = csc_gp.phi();
243  beta = 0.0;
244  }
245 
246  // determine the pitch of the chamber
247  // option to discretize the pitch
248  float stripPhiPitch = layer_geo->stripPhiPitch();
249  if (nStripBits_)
250  stripPhiPitch = stripPhiPitch/nStripBits_;
251 
252  // get the fit results
253  // option to discretize the fitted phi value
254  keyRadius = radius_;
255 
256  for (int i=0; i<CSCConstants::NUM_LAYERS; i++) {
257  const float fit_z = cscChamber->layer(i+1)->centerOfStrip(20).z();
258  const float fit_phi = normalizedPhi(alpha + beta * fit_z);
259  fit_z_layers.push_back(fit_z);
260  fit_phi_layers.push_back(fit_phi);
261  if (nStripBits_)
262  fit_phi_layers.push_back((std::floor(fit_phi/stripPhiPitch) + 0.5)*stripPhiPitch);
263  }
264 }
int chamber() const
Definition: CSCDetId.h:68
GlobalPoint toGlobal(const Point2DBase< Scalar, LocalTag > lp) const
Definition: Surface.h:106
float alpha
Definition: AMPTWrapper.h:95
T perp() const
Definition: PV3DBase.h:72
void getComparatorDigiCoordinates(const CSCDetId &ch_id, const CSCCorrelatedLCTDigi &stub)
const GeomDet * idToDet(DetId) const override
Definition: CSCGeometry.cc:114
constexpr T normalizedPhi(T phi)
Definition: normalizedPhi.h:9
Geom::Phi< T > phi() const
Definition: PV3DBase.h:69
const Plane & surface() const
The nominal surface of the GeomDet.
Definition: GeomDet.h:42
int endcap() const
Definition: CSCDetId.h:93
float getFractionalStrip() const
return the fractional strip. counts from 0.25
The Signals That Services Can Subscribe To This is based on ActivityRegistry and is current per Services can connect to the signals distributed by the ActivityRegistry in order to monitor the activity of the application Each possible callback has some defined which we here list in angle e< void, edm::EventID const &, edm::Timestamp const & > We also list in braces which AR_WATCH_USING_METHOD_ is used for those or
Definition: Activities.doc:12
float middleWireOfGroup(int wireGroup) const
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
const CSCLayer * layer(CSCDetId id) const
Return the layer corresponding to the given id.
Definition: CSCChamber.cc:39
void matchingComparatorDigisLCT(const CSCDetId &ch_id, const CSCCorrelatedLCTDigi &, const CSCComparatorDigiCollection &)
int ring() const
Definition: CSCDetId.h:75
const CSCChamber * chamber(CSCDetId id) const
Return the chamber corresponding to given DetId.
Definition: CSCGeometry.cc:133
const CSCGeometry * cscGeometry_
int station() const
Definition: CSCDetId.h:86
void calculateSlopeIntercept(float &alpha, float &beta)
const CSCLayerGeometry * geometry() const
Definition: CSCLayer.h:47
int getKeyWG() const
return the key wire group. counts from 0.
void CSCComparatorDigiFitter::getComparatorDigiCoordinates ( const CSCDetId ch_id,
const CSCCorrelatedLCTDigi stub 
)
private

Definition at line 148 of file CSCComparatorDigiFitter.cc.

References CSCGeometry::chamber(), compDigisIds_, cscGeometry_, cscHalfStripWidth(), ephis_, ezs_, geometry, CSCCorrelatedLCTDigi::getKeyWG(), CSCGeometry::idToDet(), M_PI, AlCaHLTBitMon_ParallelJobs::p, PV3DBase< T, PVType, FrameType >::perp(), PV3DBase< T, PVType, FrameType >::phi(), phis_, radius_, findQualityFiles::size, mathSSE::sqrt(), GeomDet::surface(), Surface::toGlobal(), PV3DBase< T, PVType, FrameType >::z(), and zs_.

Referenced by fit(), and useKeyRadius().

149 {
150  const auto& cscChamber = cscGeometry_->chamber(ch_id);
151 
152  // get the z and phi positions of the comparator digis
153  float radius_ = 0.0;
154 
155  // loop on all matching digis
156  for (const auto& p: compDigisIds_) {
157  const auto& detId = p.first;
158 
159  float phi_tmp = 0.0;
160  float radius_tmp = 0.0;
161  float z_tmp = 0.0;
162 
163  // ignore layers with no digis
164  if (p.second.empty()) continue;
165 
166  // loop on all matching digis in this layer
167  for (const auto& hit: p.second) {
168  const float fractional_strip = hit.getFractionalStrip();
169  const auto& layer_geo = cscChamber->layer(detId.layer())->geometry();
170  const float wire = layer_geo->middleWireOfGroup(stub.getKeyWG() + 1);
171 
172  // get the phi of each comparator digi
173  const LocalPoint& csc_intersect = layer_geo->intersectionOfStripAndWire(fractional_strip, wire);
174  const GlobalPoint& csc_gp = cscGeometry_->idToDet(detId)->surface().toGlobal(csc_intersect);
175  const float gpphi = csc_gp.phi();
176 
177  // normalize phi values according to first one
178  if (!phis_.empty() and gpphi>0 and phis_[0]<0 and (gpphi-phis_[0])>M_PI)
179  phi_tmp += (gpphi-2*M_PI);
180  else if (!phis_.empty() and gpphi<0 and phis_[0]>0 and (gpphi-phis_[0])<-M_PI)
181  phi_tmp += (gpphi+2*M_PI);
182  else
183  phi_tmp += (csc_gp.phi());
184 
185  z_tmp = csc_gp.z();
186  radius_tmp += csc_gp.perp();
187  }
188 
189  //in case there are more than one comparator digis in one layer
190  radius_tmp = radius_tmp/(p.second).size();
191  radius_ += radius_tmp;
192 
193  zs_.push_back(z_tmp);
194  ezs_.push_back(0);
195 
196  phi_tmp = phi_tmp/(p.second).size();
197  phis_.push_back(phi_tmp);
198  // assume that the charge is flat distributed across the half-strip
199  // this is only approximately valid, but good enough for now
200  ephis_.push_back(cscHalfStripWidth(detId)/sqrt(12));
201  }
202 }
size
Write out results.
GlobalPoint toGlobal(const Point2DBase< Scalar, LocalTag > lp) const
Definition: Surface.h:106
T perp() const
Definition: PV3DBase.h:72
const GeomDet * idToDet(DetId) const override
Definition: CSCGeometry.cc:114
Geom::Phi< T > phi() const
Definition: PV3DBase.h:69
float cscHalfStripWidth(const CSCDetId &id) const
const Plane & surface() const
The nominal surface of the GeomDet.
Definition: GeomDet.h:42
T sqrt(T t)
Definition: SSEVec.h:18
T z() const
Definition: PV3DBase.h:64
CSCComparatorDigiContainerIds compDigisIds_
#define M_PI
const CSCChamber * chamber(CSCDetId id) const
Return the chamber corresponding to given DetId.
Definition: CSCGeometry.cc:133
ESHandle< TrackerGeometry > geometry
const CSCGeometry * cscGeometry_
int getKeyWG() const
return the key wire group. counts from 0.
void CSCComparatorDigiFitter::matchingComparatorDigisLCT ( const CSCDetId ch_id,
const CSCCorrelatedLCTDigi stub,
const CSCComparatorDigiCollection hCSCComparators 
)
private

Definition at line 120 of file CSCComparatorDigiFitter.cc.

References funct::abs(), CSCDetId::chamber(), comparatorInLCTPattern(), compDigisIds_, CSCDetId::endcap(), CSCCorrelatedLCTDigi::getPattern(), CSCCorrelatedLCTDigi::getStrip(), CSCConstants::NUM_LAYERS, CSCDetId::ring(), and CSCDetId::station().

Referenced by fit(), and useKeyRadius().

121 {
122  // fetch the CSC comparator digis in this chamber
123  for (int iLayer=1; iLayer<=CSCConstants::NUM_LAYERS; ++iLayer) {
124  const CSCDetId layerId(ch_id.endcap(), ch_id.station(), ch_id.ring(), ch_id.chamber(), iLayer);
125 
126  // get the digis per layer
127  const auto& compRange = hCSCComparators.get(layerId);
128  CSCComparatorDigiContainer compDigis;
129 
130  for (auto compDigiItr = compRange.first; compDigiItr != compRange.second; compDigiItr++) {
131  const auto& compDigi = *compDigiItr;
132 
133  //if (stub.getTimeBin() < 4 or stub.getTimeBin() > 8) continue;
134  const int stubHalfStrip(compDigi.getHalfStrip());
135 
136  // these comparator digis never fit the pattern anyway!
137  if (std::abs(stubHalfStrip-stub.getStrip())>5) continue;
138 
139  // check if this comparator digi fits the pattern
140  if (comparatorInLCTPattern(stub.getStrip(), stub.getPattern(), iLayer, stubHalfStrip)) {
141  compDigis.push_back(compDigi);
142  }
143  }
144  compDigisIds_.emplace_back(layerId, compDigis);
145  }
146 }
int chamber() const
Definition: CSCDetId.h:68
int getStrip() const
return the key halfstrip from 0,159
std::vector< CSCComparatorDigi > CSCComparatorDigiContainer
int endcap() const
Definition: CSCDetId.h:93
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
CSCComparatorDigiContainerIds compDigisIds_
int ring() const
Definition: CSCDetId.h:75
bool comparatorInLCTPattern(int keyStrip, int pattern, int layer, int halfStrip) const
int getPattern() const
return pattern
int station() const
Definition: CSCDetId.h:86
void CSCComparatorDigiFitter::setGeometry ( const CSCGeometry csc_g)
inline

Definition at line 46 of file CSCComparatorDigiFitter.h.

References cscGeometry_.

46 {cscGeometry_= csc_g;}
const CSCGeometry * cscGeometry_
void CSCComparatorDigiFitter::setStripBits ( int  bits)
inline

Definition at line 49 of file CSCComparatorDigiFitter.h.

References bits, and nStripBits_.

49 {nStripBits_ = bits;}
How EventSelector::AcceptEvent() decides whether to accept an event for output otherwise it is excluding the probing of A single or multiple positive and the trigger will pass if any such matching triggers are PASS or EXCEPTION[A criterion thatmatches no triggers at all is detected and causes a throw.] A single negative with an expectation of appropriate bit checking in the decision bits
void CSCComparatorDigiFitter::useKeyRadius ( bool  useKeyRadius)
inline

Member Data Documentation

CSCComparatorDigiContainerIds CSCComparatorDigiFitter::compDigisIds_
private
const CSCGeometry* CSCComparatorDigiFitter::cscGeometry_
private

Definition at line 82 of file CSCComparatorDigiFitter.h.

Referenced by fit(), getComparatorDigiCoordinates(), and setGeometry().

const std::vector<float> CSCComparatorDigiFitter::degrees_ = {10.,10.,10.,10.,20.,10.,20.,10.,20.,10.}
private

Definition at line 100 of file CSCComparatorDigiFitter.h.

Referenced by cscHalfStripWidth().

std::vector<float> CSCComparatorDigiFitter::ephis_
private
std::vector<float> CSCComparatorDigiFitter::ezs_
private

Definition at line 93 of file CSCComparatorDigiFitter.h.

Referenced by getComparatorDigiCoordinates().

int CSCComparatorDigiFitter::nStripBits_
private

Definition at line 85 of file CSCComparatorDigiFitter.h.

Referenced by fit(), and setStripBits().

std::vector<float> CSCComparatorDigiFitter::phis_
private
float CSCComparatorDigiFitter::radius_
private

Definition at line 94 of file CSCComparatorDigiFitter.h.

Referenced by fit(), and getComparatorDigiCoordinates().

const std::vector<int> CSCComparatorDigiFitter::strips_ = {48,64,80,64, 80,80,80,80,80,80}
private

Definition at line 99 of file CSCComparatorDigiFitter.h.

Referenced by cscHalfStripWidth().

bool CSCComparatorDigiFitter::useKeyRadius_
private

Definition at line 95 of file CSCComparatorDigiFitter.h.

Referenced by fit(), and useKeyRadius().

std::vector<float> CSCComparatorDigiFitter::zs_
private