![]() |
![]() |
#include <RecoLocalMuon/CSCRecHitD/src/CSCHitFromStripOnly.h>
Public Types | |
typedef std::vector< CSCStripData > | PulseHeightMap |
Public Member Functions | |
CSCHitFromStripOnly (const edm::ParameterSet &ps) | |
std::vector< CSCStripHit > | runStrip (const CSCDetId &id, const CSCLayer *layer, const CSCStripDigiCollection::Range &rstripd) |
void | setConditions (const CSCRecoConditions *reco) |
~CSCHitFromStripOnly () | |
Protected Member Functions | |
void | fillPulseHeights (const CSCStripDigiCollection::Range &rstripd) |
Go through strip in layers and build a table with. | |
float | findHitOnStripPosition (const std::vector< CSCStripHitData > &data, const int ¢erStrip) |
Find position of hit in strip cluster in terms of strip #. | |
void | findMaxima () |
Find local maxima. | |
float | makeCluster (int centerStrip) |
Make clusters using local maxima. | |
Protected Attributes | |
CSCDetId | id_ |
const CSCLayer * | layer_ |
const CSCLayerGeometry * | layergeom_ |
const CSCChamberSpecs * | specs_ |
std::vector< int > | theClosestMaximum |
std::vector< int > | theConsecutiveStrips |
std::vector< int > | theMaxima |
PulseHeightMap | thePulseHeightMap |
Private Member Functions | |
bool | isNearDeadStrip (const CSCDetId &id, int centralStrip) |
CSCStripHitData | makeStripData (int centerStrip, int offset) |
makeStripData | |
Private Attributes | |
int | clusterSize |
float | gainWeight [80] |
These are the gain correction weights and X-talks read in from database. | |
unsigned | Nstrips |
CSCPeakBinOfStripPulse * | pulseheightOnStripFinder_ |
const CSCRecoConditions * | recoConditions_ |
Hold pointer to current conditions data. | |
std::vector< float > | strips_adc |
std::vector< float > | strips_adcRaw |
std::vector< int > | theStrips |
float | theThresholdForAPeak |
float | theThresholdForCluster |
int | tmax_cluster |
int | TmaxOfCluster |
bool | useCalib |
Static Private Attributes | |
static const int | theClusterSize = 3 |
For each of these strips, build a cluster of strip of size theClusterSize (typically 5 strips). Finally, make a Strip Hit out of these clusters by finding the center-of-mass position of the hit The DetId, strip hit position, and peaking time are stored in a CSCStripHit collection.
Be careful with the ME_1/a chambers: the 48 strips are ganged into 16 channels, Each channel contains signals from three strips, each separated by 16 strips from the next. This is the same for real data and for the MC. (This ME1a info is from Tim Cox.)
Definition at line 38 of file CSCHitFromStripOnly.h.
typedef std::vector<CSCStripData> CSCHitFromStripOnly::PulseHeightMap |
Definition at line 43 of file CSCHitFromStripOnly.h.
CSCHitFromStripOnly::CSCHitFromStripOnly | ( | const edm::ParameterSet & | ps | ) | [explicit] |
Definition at line 28 of file CSCHitFromStripOnly.cc.
References edm::ParameterSet::getUntrackedParameter(), pulseheightOnStripFinder_, theThresholdForAPeak, theThresholdForCluster, and useCalib.
00028 : recoConditions_(0) { 00029 00030 useCalib = ps.getUntrackedParameter<bool>("CSCUseCalibrations"); 00031 //theClusterSize = ps.getUntrackedParameter<int>("CSCStripClusterSize"); 00032 theThresholdForAPeak = ps.getUntrackedParameter<double>("CSCStripPeakThreshold"); 00033 theThresholdForCluster = ps.getUntrackedParameter<double>("CSCStripClusterChargeCut"); 00034 00035 pulseheightOnStripFinder_ = new CSCPeakBinOfStripPulse( ps ); 00036 00037 }
CSCHitFromStripOnly::~CSCHitFromStripOnly | ( | ) |
Definition at line 40 of file CSCHitFromStripOnly.cc.
References pulseheightOnStripFinder_.
00040 { 00041 delete pulseheightOnStripFinder_; 00042 }
void CSCHitFromStripOnly::fillPulseHeights | ( | const CSCStripDigiCollection::Range & | rstripd | ) | [protected] |
Go through strip in layers and build a table with.
Definition at line 327 of file CSCHitFromStripOnly.cc.
References edm::pset::fill(), gainWeight, id_, it, j, CSCPeakBinOfStripPulse::peakAboveBaseline(), pulseheightOnStripFinder_, CSCDetId::ring(), CSCDetId::station(), thePulseHeightMap, tmax, and useCalib.
Referenced by runStrip().
00327 { 00328 00329 // Loop over strip digis and fill the pulseheight map 00330 00331 thePulseHeightMap.clear(); 00332 thePulseHeightMap.resize(100); 00333 00334 //float maxADC = 0; 00335 //float maxADCCluster = 0; 00336 //float adc_cluster[3]; 00337 //adc_cluster[0] = 0.; 00338 //adc_cluster[1] = 0.; 00339 //adc_cluster[2] = 0.; 00340 00341 for ( CSCStripDigiCollection::const_iterator it = rstripd.first; it != rstripd.second; ++it ) { 00342 bool fill = true; 00343 int thisChannel = (*it).getStrip(); 00344 std::vector<int> sca = (*it).getADCCounts(); 00345 float height[6]; 00346 float hmax = 0.; 00347 int tmax = -1; 00348 00349 fill = pulseheightOnStripFinder_->peakAboveBaseline( (*it), hmax, tmax, height ); 00350 //float maxCluster = 0; 00351 00352 //adc_cluster[2] = adc_cluster[1]; 00353 //adc_cluster[1] = adc_cluster[0]; 00354 //adc_cluster[0] = hmax; 00355 00356 //for (int j = 0; j < 3; ++j ) maxCluster += adc_cluster[j]; 00357 00358 //if (hmax > maxADC) maxADC = hmax; 00359 //if (maxCluster > maxADCCluster ) maxADCCluster = maxCluster; 00360 00361 // Don't forget that the ME_11/a strips are ganged !!! 00362 // Have to loop 2 more times to populate strips 17-48. 00363 00364 if ( id_.station() == 1 && id_.ring() == 4 ) { 00365 for ( int j = 0; j < 3; ++j ) { 00366 thePulseHeightMap[thisChannel-1+16*j] = CSCStripData( float(thisChannel+16*j), hmax, tmax, height[0], height[1], height[2], height[3], height[4], height[5], height[0], height[1], height[2], height[3], height[4], height[5]); 00367 if ( useCalib ){ 00368 thePulseHeightMap[thisChannel-1+16*j] *= gainWeight[thisChannel-1]; 00369 } 00370 } 00371 } else { 00372 thePulseHeightMap[thisChannel-1] = CSCStripData( float(thisChannel), hmax, tmax, height[0], height[1], height[2], height[3], height[4], height[5], height[0], height[1], height[2], height[3], height[4], height[5]); 00373 if ( useCalib ){ 00374 thePulseHeightMap[thisChannel-1] *= gainWeight[thisChannel-1]; 00375 } 00376 } 00377 } 00378 }
float CSCHitFromStripOnly::findHitOnStripPosition | ( | const std::vector< CSCStripHitData > & | data, | |
const int & | centerStrip | |||
) | [protected] |
Find position of hit in strip cluster in terms of strip #.
Definition at line 474 of file CSCHitFromStripOnly.cc.
References i, LogTrace, strips_adc, strips_adcRaw, sum(), w1, w2, w3, and x.
Referenced by makeCluster().
00474 { 00475 00476 float strippos = -1.; 00477 00478 int nstrips = data.size(); 00479 if ( nstrips < 1 ) return strippos; 00480 00481 // biggestStrip is strip with largest pulse height 00482 // Use pointer subtraction 00483 00484 int biggestStrip = max_element(data.begin(), data.end()) - data.begin(); 00485 strippos = data[biggestStrip].x() * 1.; 00486 00487 // If more than one strip: use centroid to find center of cluster 00488 // but only use time bin == tmax (otherwise, bias centroid). 00489 float sum = 0.; 00490 float sum_w= 0.; 00491 00492 float stripLeft = 0.; 00493 float stripRight = 0.; 00494 00495 for ( unsigned i = 0; i != data.size(); ++i ) { 00496 float w0 = data[i].y0(); 00497 float w1 = data[i].y1(); 00498 float w2 = data[i].y2(); 00499 float w3 = data[i].y3(); 00500 float w0Raw = data[i].y0Raw(); 00501 float w1Raw = data[i].y1Raw(); 00502 float w2Raw = data[i].y2Raw(); 00503 float w3Raw = data[i].y3Raw(); 00504 00505 // Require ADC to be > 0. 00506 if (w0 < 0.) w0 = 0.001; 00507 if (w1 < 0.) w1 = 0.001; 00508 if (w2 < 0.) w2 = 0.001; 00509 if (w3 < 0.) w3 = 0.001; 00510 00511 00512 if (i == data.size()/2 -1) stripLeft = w1; 00513 if (i == data.size()/2 +1) stripRight = w1; 00514 00515 00516 // Fill the adcs to the strip hit --> needed for Gatti fitter 00517 strips_adc.push_back( w0 ); 00518 strips_adc.push_back( w1 ); 00519 strips_adc.push_back( w2 ); 00520 strips_adc.push_back( w3 ); 00521 strips_adcRaw.push_back( w0Raw ); 00522 strips_adcRaw.push_back( w1Raw ); 00523 strips_adcRaw.push_back( w2Raw ); 00524 strips_adcRaw.push_back( w3Raw ); 00525 00526 if ( data[i].x() < 1 ){ 00527 LogTrace("CSCRecHit") << "problem in indexing of strip, strip id is: " << data[i].x() << "\n"; 00528 } 00529 sum_w += w1; 00530 sum += w1 * data[i].x(); 00531 } 00532 00533 if ( sum_w > 0.) strippos = sum / sum_w; 00534 00535 return strippos; 00536 }
void CSCHitFromStripOnly::findMaxima | ( | ) | [protected] |
Find local maxima.
Definition at line 387 of file CSCHitFromStripOnly.cc.
References GenMuonPlsPt100GeV_cfg::cout, lat::endl(), i, j, k, edm::es::l(), t, theClosestMaximum, theConsecutiveStrips, theMaxima, thePulseHeightMap, theThresholdForAPeak, theThresholdForCluster, and gen_jpsi2muons_cfg::ymax.
Referenced by runStrip().
00387 { 00388 00389 theMaxima.clear(); 00390 theConsecutiveStrips.clear(); 00391 theClosestMaximum.clear(); 00392 for ( size_t i = 0; i < thePulseHeightMap.size(); ++i ) { 00393 00394 float heightPeak = thePulseHeightMap[i].ymax(); 00395 00396 // sum 3 strips so that hits between strips are not suppressed 00397 float heightCluster; 00398 00399 bool maximumFound = false; 00400 // Left edge of chamber 00401 if ( i == 0 ) { 00402 heightCluster = thePulseHeightMap[i].ymax()+thePulseHeightMap[i+1].ymax(); 00403 // Have found a strip Hit if... 00404 if (( heightPeak > theThresholdForAPeak ) && 00405 ( heightCluster > theThresholdForCluster ) && 00406 ( thePulseHeightMap[i].ymax() >= thePulseHeightMap[i+1].ymax()) && 00407 ( thePulseHeightMap[i].t() > 2 ) && 00408 ( thePulseHeightMap[i].t() < 7 )) { 00409 theMaxima.push_back(i); 00410 maximumFound = true; 00411 } 00412 // Right edge of chamber 00413 } else if ( i == thePulseHeightMap.size()-1) { 00414 heightCluster = thePulseHeightMap[i-1].ymax()+thePulseHeightMap[i].ymax(); 00415 // Have found a strip Hit if... 00416 if (( heightPeak > theThresholdForAPeak ) && 00417 ( heightCluster > theThresholdForCluster ) && 00418 ( thePulseHeightMap[i].ymax() > thePulseHeightMap[i-1].ymax()) && 00419 ( thePulseHeightMap[i].t() > 2 ) && 00420 ( thePulseHeightMap[i].t() < 7 )) { 00421 theMaxima.push_back(i); 00422 maximumFound = true; 00423 } 00424 // Any other strips 00425 } else { 00426 heightCluster = thePulseHeightMap[i-1].ymax()+thePulseHeightMap[i].ymax()+thePulseHeightMap[i+1].ymax(); 00427 // Have found a strip Hit if... 00428 if (( heightPeak > theThresholdForAPeak ) && 00429 ( heightCluster > theThresholdForCluster ) && 00430 ( thePulseHeightMap[i].ymax() > thePulseHeightMap[i-1].ymax()) && 00431 ( thePulseHeightMap[i].ymax() >= thePulseHeightMap[i+1].ymax()) && 00432 ( thePulseHeightMap[i].t() > 2 ) && 00433 ( thePulseHeightMap[i].t() < 7 )) { 00434 theMaxima.push_back(i); 00435 maximumFound = true; 00436 } 00437 } 00438 //---- Consecutive strips with charge (real cluster); if too wide - measurement is not accurate 00439 if(maximumFound){ 00440 int numberOfConsecutiveStrips = 1; 00441 float testThreshold = 10.;//---- ADC counts; 00442 //---- this is not XTalk corrected so it is correct in first approximation only 00443 int j = 0; 00444 for(int l = 0;l<8 ;l++){ 00445 if(j<0){ 00446 std::cout<<" CSCRecHitD: This should never occur!!! Contact CSC expert!"<<std::endl; 00447 } 00448 j++; 00449 bool signalPresent = false; 00450 for(int k =0;k<2;k++){ 00451 j*= -1;//---- check from left and right 00452 int anotherConsecutiveStrip = i+j; 00453 if(anotherConsecutiveStrip>=0 && anotherConsecutiveStrip<int(thePulseHeightMap.size() )){ 00454 if(thePulseHeightMap[anotherConsecutiveStrip].ymax()>testThreshold){ 00455 numberOfConsecutiveStrips++; 00456 signalPresent = true; 00457 } 00458 } 00459 } 00460 if(!signalPresent){ 00461 break; 00462 } 00463 } 00464 theConsecutiveStrips.push_back(numberOfConsecutiveStrips); 00465 } 00466 } 00467 }
Definition at line 538 of file CSCHitFromStripOnly.cc.
References CSCRecoConditions::nearBadStrip(), and recoConditions_.
Referenced by runStrip().
00538 { 00539 00540 // const std::bitset<80> & deadStrips = recoConditions_->badStripWord( id ); 00541 // bool isDead = false; 00542 // if(centralStrip>0 && centralStrip<79){ 00543 // isDead = (deadStrips.test(centralStrip+1) || deadStrips.test(centralStrip-1)); 00544 // } 00545 // return isDead; 00546 00547 //@@ Tim says: not sure I understand this properly... 00548 return recoConditions_->nearBadStrip( id, centralStrip ); 00549 00550 }
float CSCHitFromStripOnly::makeCluster | ( | int | centerStrip | ) | [protected] |
Make clusters using local maxima.
Definition at line 123 of file CSCHitFromStripOnly.cc.
References clusterSize, data, findHitOnStripPosition(), i, makeStripData(), CSCChamberSpecs::nStrips(), specs_, theClusterSize, and theStrips.
Referenced by runStrip().
00123 { 00124 00125 float strippos = -1.; 00126 clusterSize = theClusterSize; 00127 std::vector<CSCStripHitData> stripDataV; 00128 00129 // We only want to use strip position in terms of strip # for the strip hit. 00130 00131 // If the cluster size is such that you go beyond the edge of detector, shrink cluster appropriatly 00132 for ( int i = 1; i < theClusterSize/2 + 1; ++i) { 00133 00134 if ( centerStrip - i < 1 || centerStrip + i > specs_->nStrips() ) { 00135 00136 // Shrink cluster size, but keep it an odd number of strips. 00137 clusterSize = 2*i - 1; 00138 } 00139 } 00140 00141 for ( int i = -clusterSize/2; i <= clusterSize/2; ++i ) { 00142 CSCStripHitData data = makeStripData(centerStrip, i); 00143 stripDataV.push_back( data ); 00144 theStrips.push_back( centerStrip + i ); 00145 } 00146 00147 strippos = findHitOnStripPosition( stripDataV, centerStrip ); 00148 00149 return strippos; 00150 }
CSCStripHitData CSCHitFromStripOnly::makeStripData | ( | int | centerStrip, | |
int | offset | |||
) | [private] |
makeStripData
Definition at line 156 of file CSCHitFromStripOnly.cc.
References ecalMGPA::adc(), clusterSize, find(), i, k, LogTrace, theMaxima, thePulseHeightMap, tmax, and TmaxOfCluster.
Referenced by makeCluster().
00156 { 00157 00158 CSCStripHitData prelimData(-1.,0.,0.,0.,0.,0.,0.,0.,0.,0); 00159 int thisStrip = centerStrip+offset; 00160 00161 int tmax = thePulseHeightMap[centerStrip-1].t(); 00162 TmaxOfCluster = tmax; 00163 00164 float adc[12]; 00165 float adcRaw[12]; 00166 00167 if ( tmax == 3 ) { 00168 adc[0] = thePulseHeightMap[thisStrip-1].y2(); 00169 adc[1] = thePulseHeightMap[thisStrip-1].y3(); 00170 adc[2] = thePulseHeightMap[thisStrip-1].y4(); 00171 adc[3] = thePulseHeightMap[thisStrip-1].y5(); 00172 adcRaw[0] = thePulseHeightMap[thisStrip-1].y2Raw(); 00173 adcRaw[1] = thePulseHeightMap[thisStrip-1].y3Raw(); 00174 adcRaw[2] = thePulseHeightMap[thisStrip-1].y4Raw(); 00175 adcRaw[3] = thePulseHeightMap[thisStrip-1].y5Raw(); 00176 } else if ( tmax == 4 ) { 00177 adc[0] = thePulseHeightMap[thisStrip-1].y3(); 00178 adc[1] = thePulseHeightMap[thisStrip-1].y4(); 00179 adc[2] = thePulseHeightMap[thisStrip-1].y5(); 00180 adc[3] = thePulseHeightMap[thisStrip-1].y6(); 00181 adcRaw[0] = thePulseHeightMap[thisStrip-1].y3Raw(); 00182 adcRaw[1] = thePulseHeightMap[thisStrip-1].y4Raw(); 00183 adcRaw[2] = thePulseHeightMap[thisStrip-1].y5Raw(); 00184 adcRaw[3] = thePulseHeightMap[thisStrip-1].y6Raw(); 00185 } else if ( tmax == 5 ) { 00186 adc[0] = thePulseHeightMap[thisStrip-1].y4(); 00187 adc[1] = thePulseHeightMap[thisStrip-1].y5(); 00188 adc[2] = thePulseHeightMap[thisStrip-1].y6(); 00189 adc[3] = thePulseHeightMap[thisStrip-1].y7(); 00190 adcRaw[0] = thePulseHeightMap[thisStrip-1].y4Raw(); 00191 adcRaw[1] = thePulseHeightMap[thisStrip-1].y5Raw(); 00192 adcRaw[2] = thePulseHeightMap[thisStrip-1].y6Raw(); 00193 adcRaw[3] = thePulseHeightMap[thisStrip-1].y7Raw(); 00194 } else if ( tmax == 6 ) { 00195 adc[0] = thePulseHeightMap[thisStrip-1].y5(); 00196 adc[1] = thePulseHeightMap[thisStrip-1].y6(); 00197 adc[2] = thePulseHeightMap[thisStrip-1].y7(); 00198 adc[3] = 0.1; 00199 adcRaw[0] = thePulseHeightMap[thisStrip-1].y5Raw(); 00200 adcRaw[1] = thePulseHeightMap[thisStrip-1].y6Raw(); 00201 adcRaw[2] = thePulseHeightMap[thisStrip-1].y7Raw(); 00202 adcRaw[3] = 0.1; 00203 } else { 00204 adc[0] = 0.1; 00205 adc[1] = 0.1; 00206 adc[2] = 0.1; 00207 adc[3] = 0.1; 00208 adcRaw[0] = 0.1; 00209 adcRaw[1] = 0.1; 00210 adcRaw[2] = 0.1; 00211 adcRaw[3] = 0.1; 00212 LogTrace("CSCRecHit") << "[CSCHitFromStripOnly::makeStripData] Tmax out of range: contact CSC expert!!!" << "\n"; 00213 } 00214 00215 if ( offset == 0 ) { 00216 prelimData = CSCStripHitData(thisStrip, adc[0], adc[1], adc[2], adc[3], 00217 adcRaw[0], adcRaw[1], adcRaw[2], adcRaw[3],TmaxOfCluster); 00218 } else { 00219 int sign = offset>0 ? 1 : -1; 00220 // If there's another maximum that would like to use part of this cluster, 00221 // it gets shared in proportion to the height of the maxima 00222 for ( int i = 1; i <= clusterSize/2; ++i ) { 00223 00224 // Find the direction of the offset 00225 int testStrip = thisStrip + sign*i; 00226 std::vector<int>::iterator otherMax = find(theMaxima.begin(), theMaxima.end(), testStrip-1); 00227 00228 // No other maxima found, so just store 00229 if ( otherMax == theMaxima.end() ) { 00230 prelimData = CSCStripHitData(thisStrip, adc[0], adc[1], adc[2], adc[3], adcRaw[0], adcRaw[1], adcRaw[2], adcRaw[3],TmaxOfCluster); 00231 } else { 00232 if ( tmax == 3 ) { 00233 adc[4] = thePulseHeightMap[testStrip-1].y2(); 00234 adc[5] = thePulseHeightMap[testStrip-1].y3(); 00235 adc[6] = thePulseHeightMap[testStrip-1].y4(); 00236 adc[7] = thePulseHeightMap[testStrip-1].y5(); 00237 adc[8] = thePulseHeightMap[centerStrip-1].y2(); 00238 adc[9] = thePulseHeightMap[centerStrip-1].y3(); 00239 adc[10] = thePulseHeightMap[centerStrip-1].y4(); 00240 adc[11] = thePulseHeightMap[centerStrip-1].y5(); 00241 adcRaw[4] = thePulseHeightMap[testStrip-1].y2Raw(); 00242 adcRaw[5] = thePulseHeightMap[testStrip-1].y3Raw(); 00243 adcRaw[6] = thePulseHeightMap[testStrip-1].y4Raw(); 00244 adcRaw[7] = thePulseHeightMap[testStrip-1].y5Raw(); 00245 adcRaw[8] = thePulseHeightMap[centerStrip-1].y2Raw(); 00246 adcRaw[9] = thePulseHeightMap[centerStrip-1].y3Raw(); 00247 adcRaw[10] = thePulseHeightMap[centerStrip-1].y4Raw(); 00248 adcRaw[11] = thePulseHeightMap[centerStrip-1].y5Raw(); 00249 } else if ( tmax == 4 ) { 00250 adc[4] = thePulseHeightMap[testStrip-1].y3(); 00251 adc[5] = thePulseHeightMap[testStrip-1].y4(); 00252 adc[6] = thePulseHeightMap[testStrip-1].y5(); 00253 adc[7] = thePulseHeightMap[testStrip-1].y6(); 00254 adc[8] = thePulseHeightMap[centerStrip-1].y3(); 00255 adc[9] = thePulseHeightMap[centerStrip-1].y4(); 00256 adc[10] = thePulseHeightMap[centerStrip-1].y5(); 00257 adc[11] = thePulseHeightMap[centerStrip-1].y6(); 00258 adcRaw[4] = thePulseHeightMap[testStrip-1].y3Raw(); 00259 adcRaw[5] = thePulseHeightMap[testStrip-1].y4Raw(); 00260 adcRaw[6] = thePulseHeightMap[testStrip-1].y5Raw(); 00261 adcRaw[7] = thePulseHeightMap[testStrip-1].y6Raw(); 00262 adcRaw[8] = thePulseHeightMap[centerStrip-1].y3Raw(); 00263 adcRaw[9] = thePulseHeightMap[centerStrip-1].y4Raw(); 00264 adcRaw[10] = thePulseHeightMap[centerStrip-1].y5Raw(); 00265 adcRaw[11] = thePulseHeightMap[centerStrip-1].y6Raw(); 00266 } else if ( tmax == 5 ) { 00267 adc[4] = thePulseHeightMap[testStrip-1].y4(); 00268 adc[5] = thePulseHeightMap[testStrip-1].y5(); 00269 adc[6] = thePulseHeightMap[testStrip-1].y6(); 00270 adc[7] = thePulseHeightMap[testStrip-1].y7(); 00271 adc[8] = thePulseHeightMap[centerStrip-1].y4(); 00272 adc[9] = thePulseHeightMap[centerStrip-1].y5(); 00273 adc[10] = thePulseHeightMap[centerStrip-1].y6(); 00274 adc[11] = thePulseHeightMap[centerStrip-1].y7(); 00275 adcRaw[4] = thePulseHeightMap[testStrip-1].y4Raw(); 00276 adcRaw[5] = thePulseHeightMap[testStrip-1].y5Raw(); 00277 adcRaw[6] = thePulseHeightMap[testStrip-1].y6Raw(); 00278 adcRaw[7] = thePulseHeightMap[testStrip-1].y7Raw(); 00279 adcRaw[8] = thePulseHeightMap[centerStrip-1].y4Raw(); 00280 adcRaw[9] = thePulseHeightMap[centerStrip-1].y5Raw(); 00281 adcRaw[10] = thePulseHeightMap[centerStrip-1].y6Raw(); 00282 adcRaw[11] = thePulseHeightMap[centerStrip-1].y7Raw(); 00283 } else if ( tmax == 6 ) { 00284 adc[4] = thePulseHeightMap[testStrip-1].y5(); 00285 adc[5] = thePulseHeightMap[testStrip-1].y6(); 00286 adc[6] = thePulseHeightMap[testStrip-1].y7(); 00287 adc[7] = 0.1; 00288 adc[8] = thePulseHeightMap[centerStrip-1].y5(); 00289 adc[9] = thePulseHeightMap[centerStrip-1].y6(); 00290 adc[10] = thePulseHeightMap[centerStrip-1].y7(); 00291 adc[11] = 0.1; 00292 adcRaw[4] = thePulseHeightMap[testStrip-1].y5Raw(); 00293 adcRaw[5] = thePulseHeightMap[testStrip-1].y6Raw(); 00294 adcRaw[6] = thePulseHeightMap[testStrip-1].y7Raw(); 00295 adcRaw[7] = 0.1; 00296 adcRaw[8] = thePulseHeightMap[centerStrip-1].y5Raw(); 00297 adcRaw[9] = thePulseHeightMap[centerStrip-1].y6Raw(); 00298 adcRaw[10] = thePulseHeightMap[centerStrip-1].y7Raw(); 00299 adcRaw[11] = 0.1; 00300 } else { 00301 for (int k = 4; k < 12; ++k ){ 00302 adc[k] = adcRaw[k] = 0.1; 00303 } 00304 } 00305 // Scale shared strip B by ratio of peak of ADC counts from central strip A 00306 // and neighbouring maxima C 00307 for (int k = 0; k < 4; ++k){ 00308 if(adc[k+4]>0 && adc[k+8]>0){ 00309 adc[k] = adc[k] * adc[k+8] / (adc[k+4]+adc[k+8]); 00310 } 00311 if(adcRaw[k+4]>0 && adcRaw[k+8]>0){ 00312 adcRaw[k] = adcRaw[k] * adcRaw[k+8] / (adcRaw[k+4]+adcRaw[k+8]); 00313 } 00314 } 00315 prelimData = CSCStripHitData(thisStrip, adc[0], adc[1], adc[2], adc[3], 00316 adcRaw[0], adcRaw[1], adcRaw[2], adcRaw[3], TmaxOfCluster); 00317 } 00318 } 00319 } 00320 return prelimData; 00321 }
std::vector< CSCStripHit > CSCHitFromStripOnly::runStrip | ( | const CSCDetId & | id, | |
const CSCLayer * | layer, | |||
const CSCStripDigiCollection::Range & | rstripd | |||
) |
Definition at line 51 of file CSCHitFromStripOnly.cc.
References CSCLayer::chamber(), clusterSize, fillPulseHeights(), findMaxima(), gainWeight, CSCLayer::geometry(), id_, isNearDeadStrip(), layer_, layergeom_, makeCluster(), Nstrips, CSCChamberSpecs::nStrips(), recoConditions_, CSCChamber::specs(), specs_, strips_adc, strips_adcRaw, CSCRecoConditions::stripWeights(), theClosestMaximum, theClusterSize, theConsecutiveStrips, theMaxima, theStrips, TmaxOfCluster, and useCalib.
Referenced by CSCRecHitDBuilder::build().
00053 { 00054 00055 std::vector<CSCStripHit> hitsInLayer; 00056 00057 // cache layer info for ease of access 00058 id_ = id; 00059 layer_ = layer; 00060 layergeom_ = layer_->geometry(); 00061 specs_ = layer->chamber()->specs(); 00062 Nstrips = specs_->nStrips(); 00063 00064 TmaxOfCluster = 5; 00065 00066 // Get gain correction weights for all strips in layer, and cache in gainWeight. 00067 // They're used in fillPulseHeights below. 00068 if ( useCalib ) { 00069 recoConditions_->stripWeights( id, gainWeight ); 00070 } 00071 00072 // Fill adc map and find maxima (potential hits) 00073 fillPulseHeights( rstripd ); 00074 findMaxima(); 00075 00076 // Make a Strip Hit out of each strip local maximum 00077 for ( unsigned imax = 0; imax < theMaxima.size(); ++imax ) { 00078 00079 // Initialize parameters entering the CSCStripHit 00080 clusterSize = theClusterSize; 00081 theStrips.clear(); 00082 strips_adc.clear(); 00083 strips_adcRaw.clear(); 00084 00085 // This is where centroid position is determined 00086 // The strips_adc vector is also filled here 00087 // Remember, the array starts at 0, but the stripId starts at 1... 00088 float strippos = makeCluster( theMaxima[imax]+1 ); 00089 00090 if ( strippos < 0 || TmaxOfCluster < 3 ) continue; 00091 00092 //---- If two maxima are too close the error assigned will be width/sqrt(12) - see CSCXonStrip_MatchGatti.cc 00093 int maximum_to_left = 99; //---- If there is one maximum - the distance is set to 99 (strips) 00094 int maximum_to_right = 99; 00095 if(imax<theMaxima.size()-1){ 00096 maximum_to_right = theMaxima.at(imax+1) - theMaxima.at(imax); 00097 } 00098 if(imax>0 ){ 00099 maximum_to_left = theMaxima.at(imax-1) - theMaxima.at(imax); 00100 } 00101 if(fabs(maximum_to_right) < fabs(maximum_to_left)){ 00102 theClosestMaximum.push_back(maximum_to_right); 00103 } 00104 else{ 00105 theClosestMaximum.push_back(maximum_to_left); 00106 } 00107 00108 //---- Check if a neighbouring strip is a dead strip 00109 bool deadStrip = isNearDeadStrip(id, theMaxima.at(imax)); 00110 00111 CSCStripHit striphit( id, strippos, TmaxOfCluster, theStrips, strips_adc, strips_adcRaw, 00112 theConsecutiveStrips.at(imax), theClosestMaximum.at(imax), deadStrip); 00113 hitsInLayer.push_back( striphit ); 00114 } 00115 00116 return hitsInLayer; 00117 }
void CSCHitFromStripOnly::setConditions | ( | const CSCRecoConditions * | reco | ) | [inline] |
Definition at line 51 of file CSCHitFromStripOnly.h.
References recoConditions_.
Referenced by CSCRecHitDBuilder::setConditions().
00051 { 00052 recoConditions_ = reco; 00053 }
int CSCHitFromStripOnly::clusterSize [private] |
Definition at line 87 of file CSCHitFromStripOnly.h.
Referenced by makeCluster(), makeStripData(), and runStrip().
float CSCHitFromStripOnly::gainWeight[80] [private] |
These are the gain correction weights and X-talks read in from database.
Definition at line 100 of file CSCHitFromStripOnly.h.
Referenced by fillPulseHeights(), and runStrip().
CSCDetId CSCHitFromStripOnly::id_ [protected] |
Definition at line 75 of file CSCHitFromStripOnly.h.
Referenced by fillPulseHeights(), and runStrip().
const CSCLayer* CSCHitFromStripOnly::layer_ [protected] |
const CSCLayerGeometry* CSCHitFromStripOnly::layergeom_ [protected] |
unsigned CSCHitFromStripOnly::Nstrips [private] |
Definition at line 110 of file CSCHitFromStripOnly.h.
Referenced by CSCHitFromStripOnly(), fillPulseHeights(), and ~CSCHitFromStripOnly().
const CSCRecoConditions* CSCHitFromStripOnly::recoConditions_ [private] |
Hold pointer to current conditions data.
Definition at line 108 of file CSCHitFromStripOnly.h.
Referenced by isNearDeadStrip(), runStrip(), and setConditions().
const CSCChamberSpecs* CSCHitFromStripOnly::specs_ [protected] |
std::vector<float> CSCHitFromStripOnly::strips_adc [private] |
Definition at line 88 of file CSCHitFromStripOnly.h.
Referenced by findHitOnStripPosition(), and runStrip().
std::vector<float> CSCHitFromStripOnly::strips_adcRaw [private] |
Definition at line 89 of file CSCHitFromStripOnly.h.
Referenced by findHitOnStripPosition(), and runStrip().
std::vector<int> CSCHitFromStripOnly::theClosestMaximum [protected] |
const int CSCHitFromStripOnly::theClusterSize = 3 [static, private] |
std::vector<int> CSCHitFromStripOnly::theConsecutiveStrips [protected] |
std::vector<int> CSCHitFromStripOnly::theMaxima [protected] |
Definition at line 66 of file CSCHitFromStripOnly.h.
Referenced by findMaxima(), makeStripData(), and runStrip().
PulseHeightMap CSCHitFromStripOnly::thePulseHeightMap [protected] |
Definition at line 70 of file CSCHitFromStripOnly.h.
Referenced by fillPulseHeights(), findMaxima(), and makeStripData().
std::vector<int> CSCHitFromStripOnly::theStrips [private] |
float CSCHitFromStripOnly::theThresholdForAPeak [private] |
Definition at line 95 of file CSCHitFromStripOnly.h.
Referenced by CSCHitFromStripOnly(), and findMaxima().
float CSCHitFromStripOnly::theThresholdForCluster [private] |
Definition at line 96 of file CSCHitFromStripOnly.h.
Referenced by CSCHitFromStripOnly(), and findMaxima().
int CSCHitFromStripOnly::tmax_cluster [private] |
Definition at line 86 of file CSCHitFromStripOnly.h.
int CSCHitFromStripOnly::TmaxOfCluster [private] |
Definition at line 103 of file CSCHitFromStripOnly.h.
Referenced by makeStripData(), and runStrip().
bool CSCHitFromStripOnly::useCalib [private] |
Definition at line 93 of file CSCHitFromStripOnly.h.
Referenced by CSCHitFromStripOnly(), fillPulseHeights(), and runStrip().