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