00001 #include "RecoLocalTracker/SiPixelRecHits/interface/PixelCPEGeneric.h"
00002
00003 #include "Geometry/TrackerGeometryBuilder/interface/PixelGeomDetUnit.h"
00004 #include "Geometry/TrackerGeometryBuilder/interface/RectangularPixelTopology.h"
00005
00006
00007 #include "RecoLocalTracker/SiPixelRecHits/interface/SiPixelTemplate.h"
00008 #include "DataFormats/SiPixelDetId/interface/PXBDetId.h"
00009 #include "DataFormats/DetId/interface/DetId.h"
00010
00011
00012
00013 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00014 #include "MagneticField/Engine/interface/MagneticField.h"
00015
00016 #include "boost/multi_array.hpp"
00017
00018 #include <iostream>
00019 using namespace std;
00020
00021 const double HALF_PI = 1.57079632679489656;
00022
00023
00025
00026 PixelCPEGeneric::PixelCPEGeneric(edm::ParameterSet const & conf,
00027 const MagneticField * mag, const SiPixelLorentzAngle * lorentzAngle, const SiPixelCPEGenericErrorParm * genErrorParm, const SiPixelTemplateDBObject * templateDBobject)
00028 : PixelCPEBase(conf, mag, lorentzAngle, genErrorParm, templateDBobject)
00029 {
00030
00031 if (theVerboseLevel > 0)
00032 LogDebug("PixelCPEGeneric")
00033 << " constructing a generic algorithm for ideal pixel detector.\n"
00034 << " CPEGeneric:: VerboseLevel = " << theVerboseLevel;
00035
00036
00037 the_eff_charge_cut_lowX = conf.getParameter<double>("eff_charge_cut_lowX");
00038 the_eff_charge_cut_lowY = conf.getParameter<double>("eff_charge_cut_lowY");
00039 the_eff_charge_cut_highX = conf.getParameter<double>("eff_charge_cut_highX");
00040 the_eff_charge_cut_highY = conf.getParameter<double>("eff_charge_cut_highY");
00041 the_size_cutX = conf.getParameter<double>("size_cutX");
00042 the_size_cutY = conf.getParameter<double>("size_cutY");
00043
00044 EdgeClusterErrorX_ = conf.getParameter<double>("EdgeClusterErrorX");
00045 EdgeClusterErrorY_ = conf.getParameter<double>("EdgeClusterErrorY");
00046
00047
00048 inflate_errors = conf.getParameter<bool>("inflate_errors");
00049 inflate_all_errors_no_trk_angle = conf.getParameter<bool>("inflate_all_errors_no_trk_angle");
00050
00051 UseErrorsFromTemplates_ = conf.getParameter<bool>("UseErrorsFromTemplates");
00052 TruncatePixelCharge_ = conf.getParameter<bool>("TruncatePixelCharge");
00053 IrradiationBiasCorrection_ = conf.getParameter<bool>("IrradiationBiasCorrection");
00054 DoCosmics_ = conf.getParameter<bool>("DoCosmics");
00055 LoadTemplatesFromDB_ = conf.getParameter<bool>("LoadTemplatesFromDB");
00056
00057 if ( !UseErrorsFromTemplates_ && ( TruncatePixelCharge_ ||
00058 IrradiationBiasCorrection_ ||
00059 DoCosmics_ ||
00060 LoadTemplatesFromDB_ ) )
00061 {
00062 throw cms::Exception("PixelCPEGeneric::PixelCPEGeneric: ")
00063 << "\nERROR: UseErrorsFromTemplates_ is set to False in PixelCPEGeneric_cfi.py. "
00064 << " In this case it does not make sense to set any of the following to True: "
00065 << " TruncatePixelCharge_, IrradiationBiasCorrection_, DoCosmics_, LoadTemplatesFromDB_ !!!"
00066 << "\n\n";
00067 }
00068
00069 if ( UseErrorsFromTemplates_ )
00070 {
00071 templID_ = -999;
00072 if ( LoadTemplatesFromDB_ )
00073 {
00074
00075 if ( !templ_.pushfile( *templateDBobject_) )
00076 throw cms::Exception("InvalidCalibrationLoaded")
00077 << "ERROR: Templates not filled correctly. Check the sqlite file. Using SiPixelTemplateDBObject version "
00078 << ( *templateDBobject_ ).version() << ". Template ID is " << templID_;
00079 }
00080 else
00081 {
00082 if ( !templ_.pushfile( templID_ ) )
00083 throw cms::Exception("InvalidCalibrationLoaded")
00084 << "ERROR: Templates not loaded correctly from text file. Reconstruction will fail." << " Template ID is " << templID_;
00085 }
00086
00087 }
00088
00089
00090
00091
00092
00093
00094
00095
00096
00097
00098
00099
00100 xerr_barrel_l1_= {0.00115, 0.00120, 0.00088};
00101 xerr_barrel_l1_def_=0.01030;
00102 yerr_barrel_l1_= {0.00375,0.00230,0.00250,0.00250,0.00230,0.00230,0.00210,0.00210,0.00240};
00103 yerr_barrel_l1_def_=0.00210;
00104 xerr_barrel_ln_= {0.00115, 0.00120, 0.00088};
00105 xerr_barrel_ln_def_=0.01030;
00106 yerr_barrel_ln_= {0.00375,0.00230,0.00250,0.00250,0.00230,0.00230,0.00210,0.00210,0.00240};
00107 yerr_barrel_ln_def_=0.00210;
00108 xerr_endcap_= {0.0020, 0.0020};
00109 xerr_endcap_def_=0.0020;
00110 yerr_endcap_= {0.00210};
00111 yerr_endcap_def_=0.00075;
00112
00113 bool isUpgrade=false;
00114 if ( conf.exists("Upgrade") && conf.getParameter<bool>("Upgrade")) {
00115 isUpgrade=true;
00116 xerr_barrel_ln_= {0.00114,0.00104,0.00214};
00117 xerr_barrel_ln_def_=0.00425;
00118 yerr_barrel_ln_= {0.00299,0.00203,0.0023,0.00237,0.00233,0.00243,0.00232,0.00259,0.00176};
00119 yerr_barrel_ln_def_=0.00245;
00120 xerr_endcap_= {0.00151,0.000813,0.00221};
00121 xerr_endcap_def_=0.00218;
00122 yerr_endcap_= {0.00261,0.00107,0.00264};
00123 yerr_endcap_def_=0.00357;
00124
00125 if ( conf.exists("SmallPitch") && conf.getParameter<bool>("SmallPitch")) {
00126 xerr_barrel_l1_= {0.00104, 0.000691, 0.00122};
00127 xerr_barrel_l1_def_=0.00321;
00128 yerr_barrel_l1_= {0.00199,0.00136,0.0015,0.00153,0.00152,0.00171,0.00154,0.00157,0.00154};
00129 yerr_barrel_l1_def_=0.00164;
00130 }
00131 else{
00132 xerr_barrel_l1_= {0.00114,0.00104,0.00214};
00133 xerr_barrel_l1_def_=0.00425;
00134 yerr_barrel_l1_= {0.00299,0.00203,0.0023,0.00237,0.00233,0.00243,0.00232,0.00259,0.00176};
00135 yerr_barrel_l1_def_=0.00245;
00136 }
00137 }
00138
00139 isUpgrade_=isUpgrade;
00140
00141 }
00142
00143
00144
00145
00149
00150 LocalPoint
00151 PixelCPEGeneric::localPosition(const SiPixelCluster& cluster,
00152 const GeomDetUnit & det) const
00153 {
00154 setTheDet( det, cluster );
00155 computeLorentzShifts();
00156 if ( UseErrorsFromTemplates_ )
00157 {
00158 templID_ = templateDBobject_->getTemplateID(theDet->geographicalId().rawId());
00159
00160
00161
00162
00163
00164
00165 float qclus = cluster.charge();
00166
00167
00168 float locBz = (*theParam).bz;
00169
00170
00171 pixmx = -999.9;
00172 sigmay = -999.9;
00173 deltay = -999.9;
00174 sigmax = -999.9;
00175 deltax = -999.9;
00176 sy1 = -999.9;
00177 dy1 = -999.9;
00178 sy2 = -999.9;
00179 dy2 = -999.9;
00180 sx1 = -999.9;
00181 dx1 = -999.9;
00182 sx2 = -999.9;
00183 dx2 = -999.9;
00184
00185 qBin_ = templ_.qbin( templID_, cotalpha_, cotbeta_, locBz, qclus,
00186 pixmx,
00187 sigmay, deltay, sigmax, deltax,
00188 sy1, dy1, sy2, dy2, sx1, dx1, sx2, dx2 );
00189
00190
00191 const float micronsToCm = 1.0e-4;
00192
00193 deltax = deltax * micronsToCm;
00194 dx1 = dx1 * micronsToCm;
00195 dx2 = dx2 * micronsToCm;
00196
00197 deltay = deltay * micronsToCm;
00198 dy1 = dy1 * micronsToCm;
00199 dy2 = dy2 * micronsToCm;
00200
00201 sigmax = sigmax * micronsToCm;
00202 sx1 = sx1 * micronsToCm;
00203 sx2 = sx2 * micronsToCm;
00204
00205 sigmay = sigmay * micronsToCm;
00206 sy1 = sy1 * micronsToCm;
00207 sy2 = sy2 * micronsToCm;
00208
00209 }
00210
00211
00212 float Q_f_X = 0.0;
00213 float Q_l_X = 0.0;
00214 float Q_m_X = 0.0;
00215 float Q_f_Y = 0.0;
00216 float Q_l_Y = 0.0;
00217 float Q_m_Y = 0.0;
00218 collect_edge_charges( cluster,
00219 Q_f_X, Q_l_X, Q_m_X,
00220 Q_f_Y, Q_l_Y, Q_m_Y );
00221
00222
00223
00224
00225
00226
00227
00228
00229
00230 MeasurementPoint meas_URcorn_LLpix( cluster.minPixelRow()+1.0,
00231 cluster.minPixelCol()+1.0 );
00232
00233
00234 MeasurementPoint meas_LLcorn_URpix( cluster.maxPixelRow(),
00235 cluster.maxPixelCol() );
00236
00237
00238
00239 LocalPoint local_URcorn_LLpix;
00240 LocalPoint local_LLcorn_URpix;
00241
00242
00243
00244
00245 if ( with_track_angle )
00246 {
00247 local_URcorn_LLpix = theTopol->localPosition(meas_URcorn_LLpix, loc_trk_pred_);
00248 local_LLcorn_URpix = theTopol->localPosition(meas_LLcorn_URpix, loc_trk_pred_);
00249 }
00250 else
00251 {
00252 local_URcorn_LLpix = theTopol->localPosition(meas_URcorn_LLpix);
00253 local_LLcorn_URpix = theTopol->localPosition(meas_LLcorn_URpix);
00254 }
00255
00256 if (theVerboseLevel > 20) {
00257 cout
00258 << "\n\t >>> cluster.x = " << cluster.x()
00259 << "\n\t >>> cluster.y = " << cluster.y()
00260 << "\n\t >>> cluster: minRow = " << cluster.minPixelRow()
00261 << " minCol = " << cluster.minPixelCol()
00262 << "\n\t >>> cluster: maxRow = " << cluster.maxPixelRow()
00263 << " maxCol = " << cluster.maxPixelCol()
00264 << "\n\t >>> meas: inner lower left = " << meas_URcorn_LLpix.x()
00265 << "," << meas_URcorn_LLpix.y()
00266 << "\n\t >>> meas: inner upper right = " << meas_LLcorn_URpix.x()
00267 << "," << meas_LLcorn_URpix.y()
00268 << endl;
00269 }
00270
00271
00272
00273
00274
00275
00276
00277 if (theVerboseLevel > 20)
00278 cout << "\t >>> Generic:: processing X" << endl;
00279
00280 float xPos =
00281 generic_position_formula( cluster.sizeX(),
00282 Q_f_X, Q_l_X,
00283 local_URcorn_LLpix.x(), local_LLcorn_URpix.x(),
00284 0.5*lorentzShiftInCmX_,
00285 cotalpha_,
00286 thePitchX,
00287 theRecTopol->isItBigPixelInX( cluster.minPixelRow() ),
00288 theRecTopol->isItBigPixelInX( cluster.maxPixelRow() ),
00289 the_eff_charge_cut_lowX,
00290 the_eff_charge_cut_highX,
00291 the_size_cutX);
00292
00293
00294 if (theVerboseLevel > 20)
00295 cout << "\t >>> Generic:: processing Y" << endl;
00296 float yPos =
00297 generic_position_formula( cluster.sizeY(),
00298 Q_f_Y, Q_l_Y,
00299 local_URcorn_LLpix.y(), local_LLcorn_URpix.y(),
00300 0.5*lorentzShiftInCmY_,
00301 cotbeta_,
00302 thePitchY,
00303 theRecTopol->isItBigPixelInY( cluster.minPixelCol() ),
00304 theRecTopol->isItBigPixelInY( cluster.maxPixelCol() ),
00305 the_eff_charge_cut_lowY,
00306 the_eff_charge_cut_highY,
00307 the_size_cutY);
00308
00309
00310 if ( IrradiationBiasCorrection_ )
00311 {
00312 if ( cluster.sizeX() == 1 )
00313 {
00314
00315
00316 bool bigInX = theRecTopol->isItBigPixelInX( cluster.maxPixelRow() );
00317
00318 if ( !bigInX )
00319 {
00320
00321 xPos -= dx1;
00322 }
00323 else
00324 {
00325
00326 xPos -= dx2;
00327 }
00328 }
00329 else
00330 {
00331
00332 xPos -= deltax;
00333 }
00334
00335 if ( cluster.sizeY() == 1 )
00336 {
00337
00338
00339 bool bigInY = theRecTopol->isItBigPixelInY( cluster.maxPixelCol() );
00340
00341 if ( !bigInY )
00342 {
00343
00344 yPos -= dy1;
00345 }
00346 else
00347 {
00348
00349 yPos -= dy2;
00350 }
00351 }
00352 else
00353 {
00354
00355 yPos -= deltay;
00356 }
00357
00358 }
00359
00360
00361 LocalPoint pos_in_local( xPos, yPos );
00362 return pos_in_local;
00363 }
00364
00365
00366
00367
00372
00373 double
00374 PixelCPEGeneric::
00375 generic_position_formula( int size,
00376 double Q_f,
00377 double Q_l,
00378 double upper_edge_first_pix,
00379 double lower_edge_last_pix,
00380 double half_lorentz_shift,
00381 double cot_angle,
00382 double pitch,
00383 bool first_is_big,
00384 bool last_is_big,
00385 double eff_charge_cut_low,
00386 double eff_charge_cut_high,
00387 double size_cut
00388 ) const
00389 {
00390 double geom_center = 0.5 * ( upper_edge_first_pix + lower_edge_last_pix );
00391
00392
00393
00394
00395 if ( size == 1 )
00396 {
00397
00398 if ( IrradiationBiasCorrection_ )
00399 return geom_center;
00400 else
00401 return geom_center + half_lorentz_shift;
00402 }
00403
00404
00405
00406
00407 double W_inner = lower_edge_last_pix - upper_edge_first_pix;
00408
00409
00410
00411 double W_pred =
00412 theThickness * cot_angle
00413 - 2 * half_lorentz_shift;
00414
00415
00416
00417 double sum_of_edge = 0.0;
00418 if (first_is_big) sum_of_edge += 2.0;
00419 else sum_of_edge += 1.0;
00420
00421 if (last_is_big) sum_of_edge += 2.0;
00422 else sum_of_edge += 1.0;
00423
00424
00425
00426 double W_eff = fabs( W_pred ) - W_inner;
00427
00428
00429
00430
00431
00432
00433
00434 bool usedEdgeAlgo = false;
00435 if (( W_eff/pitch < eff_charge_cut_low ) ||
00436 ( W_eff/pitch > eff_charge_cut_high ) || (size >= size_cut))
00437 {
00438 W_eff = pitch * 0.5 * sum_of_edge;
00439 usedEdgeAlgo = true;
00440 nRecHitsUsedEdge_++;
00441 }
00442
00443
00444
00445 double Qdiff = Q_l - Q_f;
00446 double Qsum = Q_l + Q_f;
00447
00448
00449 if(Qsum==0) Qsum=1.0;
00450 double hit_pos = geom_center + 0.5*(Qdiff/Qsum) * W_eff + half_lorentz_shift;
00451
00452
00453 if (theVerboseLevel > 20) {
00454 if ( thePart == GeomDetEnumerators::PixelBarrel ) {
00455 cout << "\t >>> We are in the Barrel." ;
00456 } else {
00457 cout << "\t >>> We are in the Forward." ;
00458 }
00459 cout
00460 << "\n\t >>> cot(angle) = " << cot_angle << " pitch = " << pitch << " size = " << size
00461 << "\n\t >>> upper_edge_first_pix = " << upper_edge_first_pix
00462 << "\n\t >>> lower_edge_last_pix = " << lower_edge_last_pix
00463 << "\n\t >>> geom_center = " << geom_center
00464 << "\n\t >>> half_lorentz_shift = " << half_lorentz_shift
00465 << "\n\t >>> W_inner = " << W_inner
00466 << "\n\t >>> W_pred = " << W_pred
00467 << "\n\t >>> W_eff(orig) = " << fabs( W_pred ) - W_inner
00468 << "\n\t >>> W_eff(used) = " << W_eff
00469 << "\n\t >>> sum_of_edge = " << sum_of_edge
00470 << "\n\t >>> Qdiff = " << Qdiff << " Qsum = " << Qsum
00471 << "\n\t >>> hit_pos = " << hit_pos
00472 << "\n\t >>> RecHits: total = " << nRecHitsTotal_
00473 << " used edge = " << nRecHitsUsedEdge_
00474 << endl;
00475 if (usedEdgeAlgo)
00476 cout << "\n\t >>> Used Edge algorithm." ;
00477 else
00478 cout << "\n\t >>> Used angle information." ;
00479 cout << endl;
00480 }
00481
00482
00483 return hit_pos;
00484 }
00485
00486
00487
00488
00489
00490
00494
00495 void
00496 PixelCPEGeneric::
00497 collect_edge_charges(const SiPixelCluster& cluster,
00498 float & Q_f_X,
00499 float & Q_l_X,
00500 float & Q_m_X,
00501 float & Q_f_Y,
00502 float & Q_l_Y,
00503 float & Q_m_Y
00504 ) const
00505 {
00506
00507 Q_f_X = Q_l_X = Q_m_X = 0.0;
00508 Q_f_Y = Q_l_Y = Q_m_Y = 0.0;
00509
00510
00511
00512 int xmin = cluster.minPixelRow();
00513 int xmax = cluster.maxPixelRow();
00514 int ymin = cluster.minPixelCol();
00515 int ymax = cluster.maxPixelCol();
00516
00517
00518
00519
00520 int isize = cluster.size();
00521
00522 for (int i = 0; i != isize; ++i)
00523 {
00524 auto const & pixel = cluster.pixel(i);
00525
00526 float pix_adc = pixel.adc;
00527 if ( UseErrorsFromTemplates_ && TruncatePixelCharge_ )
00528 pix_adc = std::min(pix_adc, pixmx );
00529
00530
00531
00532 if ( pixel.x == xmin )
00533 Q_f_X += pix_adc;
00534 else if ( pixel.x == xmax )
00535 Q_l_X += pix_adc;
00536 else
00537 Q_m_X += pix_adc;
00538
00539
00540 if ( pixel.y == ymin )
00541 Q_f_Y += pix_adc;
00542 else if ( pixel.y == ymax )
00543 Q_l_Y += pix_adc;
00544 else
00545 Q_m_Y += pix_adc;
00546 }
00547
00548 return;
00549 }
00550
00551
00552
00553
00554
00555
00556
00557 LocalError
00558 PixelCPEGeneric::localError( const SiPixelCluster& cluster,
00559 const GeomDetUnit & det) const
00560 {
00561 setTheDet( det, cluster );
00562
00563
00564 float xerr_sq = -99999.9f;
00565 float yerr_sq = -99999.9f;
00566
00567
00568
00569
00570
00571
00572
00573
00574
00575
00576
00577
00578
00579
00580
00581
00582
00583
00584
00585 const float micronsToCm = 1.0e-4f;
00586 float xerr = EdgeClusterErrorX_ * micronsToCm;
00587 float yerr = EdgeClusterErrorY_ * micronsToCm;
00588
00589
00590 int maxPixelCol = cluster.maxPixelCol();
00591 int maxPixelRow = cluster.maxPixelRow();
00592 int minPixelCol = cluster.minPixelCol();
00593 int minPixelRow = cluster.minPixelRow();
00594 unsigned int sizex = maxPixelRow - minPixelRow+1;
00595 unsigned int sizey = maxPixelCol - minPixelCol+1;
00596
00597 bool edgex = ( theRecTopol->isItEdgePixelInX( minPixelRow ) ) || ( theRecTopol->isItEdgePixelInX( maxPixelRow ) );
00598 bool edgey = ( theRecTopol->isItEdgePixelInY( minPixelCol ) ) || ( theRecTopol->isItEdgePixelInY( maxPixelCol ) );
00599
00600
00601 bool bigInX = theRecTopol->containsBigPixelInX( minPixelRow, maxPixelRow );
00602 bool bigInY = theRecTopol->containsBigPixelInY( minPixelCol, maxPixelCol );
00603 if ( isUpgrade_ ||(!with_track_angle && DoCosmics_) )
00604 {
00605
00606
00607
00608
00609 if ( thePart == GeomDetEnumerators::PixelBarrel )
00610 {
00611 DetId id = (det.geographicalId());
00612 int layer=PXBDetId(id).layer();
00613 if ( layer==1 ) {
00614 if ( !edgex )
00615 {
00616 if ( sizex<=xerr_barrel_l1_.size() ) xerr=xerr_barrel_l1_[sizex-1];
00617 else xerr=xerr_barrel_l1_def_;
00618 }
00619
00620 if ( !edgey )
00621 {
00622 if ( sizey<=yerr_barrel_l1_.size() ) yerr=yerr_barrel_l1_[sizey-1];
00623 else yerr=yerr_barrel_l1_def_;
00624 }
00625 }
00626 else{
00627 if ( !edgex )
00628 {
00629 if ( sizex<=xerr_barrel_ln_.size() ) xerr=xerr_barrel_ln_[sizex-1];
00630 else xerr=xerr_barrel_ln_def_;
00631 }
00632
00633 if ( !edgey )
00634 {
00635 if ( sizey<=yerr_barrel_ln_.size() ) yerr=yerr_barrel_ln_[sizey-1];
00636 else yerr=yerr_barrel_ln_def_;
00637 }
00638 }
00639 }
00640 else
00641 {
00642 if ( !edgex )
00643 {
00644 if ( sizex<=xerr_endcap_.size() ) xerr=xerr_endcap_[sizex-1];
00645 else xerr=xerr_endcap_def_;
00646 }
00647
00648 if ( !edgey )
00649 {
00650 if ( sizey<=yerr_endcap_.size() ) yerr=yerr_endcap_[sizey-1];
00651 else yerr=yerr_endcap_def_;
00652 }
00653 }
00654
00655 }
00656 else
00657 {
00658
00659
00660 if ( UseErrorsFromTemplates_ )
00661 {
00662 if (qBin_ == 0 && inflate_errors )
00663 {
00664 int n_bigx = 0;
00665 int n_bigy = 0;
00666
00667 int row_offset = minPixelRow;
00668 int col_offset = minPixelCol;
00669
00670 for (int irow = 0; irow < 7; ++irow)
00671 {
00672 if ( theRecTopol->isItBigPixelInX( irow+row_offset ) )
00673 ++n_bigx;
00674 }
00675
00676 for (int icol = 0; icol < 21; ++icol)
00677 {
00678 if ( theRecTopol->isItBigPixelInY( icol+col_offset ) )
00679 ++n_bigy;
00680 }
00681
00682 xerr = (float)(sizex + n_bigx) * thePitchX / sqrt( 12.0f );
00683 yerr = (float)(sizey + n_bigy) * thePitchY / sqrt( 12.0f );
00684
00685 }
00686 else
00687 {
00688
00689
00690 if ( !edgex )
00691 {
00692 if ( sizex == 1 )
00693 {
00694 if ( !bigInX )
00695 xerr = sx1;
00696 else
00697 xerr = sx2;
00698 }
00699 else
00700 xerr = sigmax;
00701
00702 }
00703
00704 if ( !edgey )
00705 {
00706 if ( sizey == 1 )
00707 {
00708 if ( !bigInY )
00709 yerr = sy1;
00710 else
00711 yerr = sy2;
00712 }
00713 else
00714 yerr = sigmay;
00715
00716 }
00717 }
00718
00719 }
00720 else
00721 {
00722
00723
00724 if ( edgex && edgey )
00725 {
00726
00727
00728 }
00729 else
00730 {
00731 pair<float,float> errPair =
00732 genErrorsFromDB_->getError( genErrorParm_, thePart, sizex, sizey,
00733 alpha_, beta_, bigInX, bigInY );
00734 if ( !edgex )
00735 xerr = errPair.first;
00736 if ( !edgey )
00737 yerr = errPair.second;
00738 }
00739
00740 if (theVerboseLevel > 9)
00741 {
00742 LogDebug("PixelCPEGeneric") <<
00743 " Sizex = " << cluster.sizeX() << " Sizey = " << cluster.sizeY() <<
00744 " Edgex = " << edgex << " Edgey = " << edgey <<
00745 " ErrX = " << xerr << " ErrY = " << yerr;
00746 }
00747
00748 }
00749
00750 }
00751
00752 if ( !(xerr > 0.0) )
00753 throw cms::Exception("PixelCPEGeneric::localError")
00754 << "\nERROR: Negative pixel error xerr = " << xerr << "\n\n";
00755
00756 if ( !(yerr > 0.0) )
00757 throw cms::Exception("PixelCPEGeneric::localError")
00758 << "\nERROR: Negative pixel error yerr = " << yerr << "\n\n";
00759
00760 xerr_sq = xerr*xerr;
00761 yerr_sq = yerr*yerr;
00762
00763 return LocalError( xerr_sq, 0, yerr_sq );
00764
00765 }
00766
00767
00768
00769
00770
00771