CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_4_5_patch3/src/RecoLocalTracker/SiPixelRecHits/plugins/PixelCPEGeneric.cc

Go to the documentation of this file.
00001 #include "PixelCPEGeneric.h"
00002 
00003 #include "Geometry/TrackerGeometryBuilder/interface/PixelGeomDetUnit.h"
00004 #include "Geometry/TrackerGeometryBuilder/interface/RectangularPixelTopology.h"
00005 
00006 // this is needed to get errors from templates
00007 #include "SiPixelTemplate.h"
00008 
00009 // Services
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   // Externally settable cuts  
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   // Externally settable flags to inflate errors
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                         // Initialize template store to the selected ID [Morris, 6/25/08]  
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         } // if ( UseErrorsFromTemplates_ )
00085   
00086   //cout << endl;
00087   //cout << "From PixelCPEGeneric::PixelCPEGeneric(...)" << endl;
00088   //cout << "(int)UseErrorsFromTemplates_ = " << (int)UseErrorsFromTemplates_    << endl;
00089   //cout << "TruncatePixelCharge_         = " << (int)TruncatePixelCharge_       << endl;      
00090   //cout << "IrradiationBiasCorrection_   = " << (int)IrradiationBiasCorrection_ << endl;
00091   //cout << "(int)DoCosmics_              = " << (int)DoCosmics_                 << endl;
00092   //cout << "(int)LoadTemplatesFromDB_    = " << (int)LoadTemplatesFromDB_       << endl;
00093   //cout << endl;
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   // ggiurgiu@jhu.edu 12/09/2010 : trk angles needed for bow/kink correction
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       bool fpix;  
00129       if ( thePart == GeomDetEnumerators::PixelBarrel )
00130         fpix = false;    // no, it's not forward -- it's barrel
00131       else
00132         fpix = true;     // yes, it's forward
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       //cout << "PixelCPEGeneric::localPosition(...) : locBz = " << locBz << endl;
00142       
00143       pixmx  = -999.9; // max pixel charge for truncation of 2-D cluster
00144       sigmay = -999.9; // CPE Generic y-error for multi-pixel cluster
00145       deltay = -999.9; // CPE Generic y-bias for multi-pixel cluster
00146       sigmax = -999.9; // CPE Generic x-error for multi-pixel cluster
00147       deltax = -999.9; // CPE Generic x-bias for multi-pixel cluster
00148       sy1    = -999.9; // CPE Generic y-error for single single-pixel
00149       dy1    = -999.9; // CPE Generic y-bias for single single-pixel cluster
00150       sy2    = -999.9; // CPE Generic y-error for single double-pixel cluster
00151       dy2    = -999.9; // CPE Generic y-bias for single double-pixel cluster
00152       sx1    = -999.9; // CPE Generic x-error for single single-pixel cluster
00153       dx1    = -999.9; // CPE Generic x-bias for single single-pixel cluster
00154       sx2    = -999.9; // CPE Generic x-error for single double-pixel cluster
00155       dx2    = -999.9; // CPE Generic x-bias for single double-pixel cluster
00156       
00157       qBin_ = templ_.qbin( templID_, cotalpha_, cotbeta_, locBz, qclus,  // inputs
00158                            pixmx,                                       // returned by reference
00159                            sigmay, deltay, sigmax, deltax,              // returned by reference
00160                            sy1, dy1, sy2, dy2, sx1, dx1, sx2, dx2 );    // returned by reference
00161       
00162       // These numbers come in microns from the qbin(...) call. Transform them to cm.
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     } // if ( UseErrorsFromTemplates_ )
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   //--- Find the inner widths along X and Y in one shot.  We
00195   //--- compute the upper right corner of the inner pixels
00196   //--- (== lower left corner of upper right pixel) and
00197   //--- the lower left corner of the inner pixels
00198   //--- (== upper right corner of lower left pixel), and then
00199   //--- subtract these two points in the formula.
00200 
00201   //--- Upper Right corner of Lower Left pixel -- in measurement frame
00202   MeasurementPoint meas_URcorn_LLpix( cluster.minPixelRow()+1.0,
00203                                       cluster.minPixelCol()+1.0 );
00204 
00205   //--- Lower Left corner of Upper Right pixel -- in measurement frame
00206   MeasurementPoint meas_LLcorn_URpix( cluster.maxPixelRow(),
00207                                       cluster.maxPixelCol() );
00208 
00209   //--- These two now converted into the local
00210   
00211   LocalPoint local_URcorn_LLpix;
00212   LocalPoint local_LLcorn_URpix;
00213 
00214 
00215   // PixelCPEGeneric can be used with or without track angles
00216   // If PixelCPEGeneric is called with track angles, use them to correct for bows/kinks:
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   //--- &&& Note that the cuts below should not be hardcoded (like in Orca and
00244   //--- &&& CPEFromDetPosition/PixelCPEInitial), but rather be
00245   //--- &&& externally settable (but tracked) parameters.  
00246 
00247 
00248   //--- Position, including the half lorentz shift
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_,   // 0.5 * lorentz shift in 
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);           // cut for eff charge width &&&
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_,   // 0.5 * lorentz shift in cm
00272                               cotbeta_,
00273                               thePitchY,   // 0.5 * lorentz shift (may be 0)
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);           // cut for eff charge width &&&
00279                              
00280   // Apply irradiation corrections.
00281   if ( IrradiationBiasCorrection_ )
00282     {
00283       if ( cluster.sizeX() == 1 )
00284         {
00285           // sanity chack
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           // Find if pixel is double (big). 
00291           bool bigInX = theTopol->isItBigPixelInX( cluster.maxPixelRow() );
00292           
00293           if ( !bigInX ) 
00294             {
00295               //cout << "Apply correction dx1 = " << dx1 << " to xPos = " << xPos << endl;
00296               xPos -= dx1;
00297             }
00298           else           
00299             {
00300               //cout << "Apply correction dx2 = " << dx2 << " to xPos = " << xPos << endl;
00301               xPos -= dx2;
00302             }
00303         }
00304       else if ( cluster.sizeX() > 1 )
00305         {
00306           //cout << "Apply correction correction_deltax = " << deltax << " to xPos = " << xPos << endl;
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       // sanity chack
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       // Find if pixel is double (big). 
00321       bool bigInY = theTopol->isItBigPixelInY( cluster.maxPixelCol() );
00322       
00323       if ( !bigInY ) 
00324         {
00325           //cout << "Apply correction dy1 = " << dy1 << " to yPos = " << yPos  << endl;
00326           yPos -= dy1;
00327         }
00328       else           
00329         {
00330           //cout << "Apply correction dy2 = " << dy2  << " to yPos = " << yPos << endl;
00331           yPos -= dy2;
00332         }
00333     }
00334   else if ( cluster.sizeY() > 1 )
00335     {
00336       //cout << "Apply correction deltay = " << deltay << " to yPos = " << yPos << endl;
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         } // if ( IrradiationBiasCorrection_ )
00344         
00345   //--- Now put the two together
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   //--- The case of only one pixel in this projection is separate.  Note that
00378   //--- here first_pix == last_pix, so the average of the two is still the
00379   //--- center of the pixel.
00380   if ( size == 1 ) 
00381     {
00382       // ggiurgiu@jhu.edu, 02/03/09 : for size = 1, the Lorentz shift is already accounted by the irradiation correction
00383       if ( IrradiationBiasCorrection_ ) 
00384         return geom_center;
00385       else
00386         return geom_center + half_lorentz_shift;
00387     }
00388 
00389 
00390   //--- Width of the clusters minus the edge (first and last) pixels.
00391   //--- In the note, they are denoted x_F and x_L (and y_F and y_L)
00392   double W_inner      = lower_edge_last_pix - upper_edge_first_pix;  // in cm
00393 
00394 
00395   //--- Predicted charge width from geometry
00396   double W_pred = 
00397     theThickness * cot_angle                     // geometric correction (in cm)
00398     - 2 * half_lorentz_shift;                    // (in cm) &&& check fpix!  
00399   
00400 
00401   //--- Total length of the two edge pixels (first+last)
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   //--- The `effective' charge width -- particle's path in first and last pixels only
00411   double W_eff = fabs( W_pred ) - W_inner;
00412 
00413 
00414   //--- If the observed charge width is inconsistent with the expectations
00415   //--- based on the track, do *not* use W_pred-W_innner.  Instead, replace
00416   //--- it with an *average* effective charge width, which is the average
00417   //--- length of the edge pixels.
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;  // ave. length of edge pixels (first+last) (cm)
00424       usedEdgeAlgo = true;
00425       nRecHitsUsedEdge_++;
00426     }
00427 
00428   
00429   //--- Finally, compute the position in this projection
00430   double Qdiff = Q_l - Q_f;
00431   double Qsum  = Q_l + Q_f;
00432 
00433         //--- Temporary fix for clusters with both first and last pixel with charge = 0
00434         if(Qsum==0) Qsum=1.0;
00435   double hit_pos = geom_center + 0.5*(Qdiff/Qsum) * W_eff + half_lorentz_shift;
00436 
00437   //--- Debugging output
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   // Initialize return variables.
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   // Fetch the pixels vector from the cluster.
00496   const vector<SiPixelCluster::Pixel>& pixelsVec = cluster.pixels();
00497 
00498   // Obtain boundaries in index units
00499   int xmin = cluster.minPixelRow();
00500   int xmax = cluster.maxPixelRow();
00501   int ymin = cluster.minPixelCol();
00502   int ymax = cluster.maxPixelCol();
00503 
00504 //   // Obtain the cluster boundaries (note: in measurement units!)
00505 //   float xmin = cluster.minPixelRow()+0.5;
00506 //   float xmax = cluster.maxPixelRow()+0.5;  
00507 //   float ymin = cluster.minPixelCol()+0.5;
00508 //   float ymax = cluster.maxPixelCol()+0.5;
00509   
00510   // Iterate over the pixels.
00511   int isize = pixelsVec.size();
00512   
00513   for (int i = 0;  i < isize; ++i) 
00514     {
00515       float pix_adc = -999.9;
00516 
00517       // ggiurgiu@fnal.gov: add pixel charge truncation
00518       if ( UseErrorsFromTemplates_ && TruncatePixelCharge_ ) 
00519         pix_adc = min( (float)(pixelsVec[i].adc), pixmx );
00520       else 
00521         pix_adc = pixelsVec[i].adc;
00522 
00523       //
00524       // X projection
00525       if      ( pixelsVec[i].x == xmin )       // need to match with tolerance!!! &&&
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       // Y projection
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 //==============  INFLATED ERROR AND ERRORS FROM DB BELOW  ================
00546 
00547 //-------------------------------------------------------------------------
00548 //  Hit error in the local frame
00549 //-------------------------------------------------------------------------
00550 LocalError 
00551 PixelCPEGeneric::localError( const SiPixelCluster& cluster, 
00552                              const GeomDetUnit & det) const 
00553 {
00554   setTheDet( det, cluster );
00555 
00556   // The squared errors
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   // Default errors are the maximum error used for edge clusters.
00564   /*
00565     int row_offset = cluster.minPixelRow();
00566     int col_offset = cluster.minPixelCol();
00567     int n_bigInX = 0;
00568     for (int irow = 0; irow < sizex; ++irow)
00569     if ( theTopol->isItBigPixelInX( irow+row_offset ) )
00570     ++n_bigInX;
00571     int n_bigInY = 0;
00572     for (int icol = 0; icol < sizey; ++icol) 
00573     if ( theTopol->isItBigPixelInY( icol+col_offset ) )
00574     ++n_bigInX;
00575     float xerr = (float)(sizex + n_bigInX) * thePitchX / sqrt(12.0);      
00576     float yerr = (float)(sizey + n_bigInY) * thePitchY / sqrt(12.0); 
00577   */
00578 
00579   // These are determined by looking at residuals for edge clusters
00580   const float micronsToCm = 1.0e-4;
00581   float xerr = EdgeClusterErrorX_ * micronsToCm;
00582   float yerr = EdgeClusterErrorY_ * micronsToCm;
00583   
00584   // Find if cluster is at the module edge. 
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   // Find if cluster contains double (big) pixels. 
00593   bool bigInX = theTopol->containsBigPixelInX( minPixelRow, maxPixelRow );       
00594   bool bigInY = theTopol->containsBigPixelInY( minPixelCol, maxPixelCol );
00595 
00596   if ( !with_track_angle && DoCosmics_ )
00597     {
00598       //cout << "Track angles are not known and we are processing cosmics." << endl; 
00599       //cout << "Default angle estimation which assumes track from PV (0,0,0) does not work." << endl;
00600       //cout << "Use an error parameterization which only depends on cluster size (by Vincenzo Chiochia)." << endl; 
00601       
00602       if ( thePart == GeomDetEnumerators::PixelBarrel ) 
00603         {
00604           if ( !edgex )
00605             {
00606               if      ( sizex == 1 ) xerr = 0.00115; // Size = 1 -> Sigma = 11.5 um 
00607               else if ( sizex == 2 ) xerr = 0.00120; // Size = 2 -> Sigma = 12.0 um      
00608               else if ( sizex == 3 ) xerr = 0.00088; // Size = 3 -> Sigma =  8.8 um
00609               else                   xerr = 0.01030;
00610             }
00611           
00612           if ( !edgey )
00613             {
00614               if      ( sizey ==  1 ) yerr = 0.00375; // 37.5 um 
00615               else if ( sizey ==  2 ) yerr = 0.00230; // 23 um      
00616               else if ( sizey ==  3 ) yerr = 0.00250; // 25 um
00617               else if ( sizey ==  4 ) yerr = 0.00250; // 25 um
00618               else if ( sizey ==  5 ) yerr = 0.00230; // 23 um
00619               else if ( sizey ==  6 ) yerr = 0.00230; // 23 um
00620               else if ( sizey ==  7 ) yerr = 0.00210; // 21 um
00621               else if ( sizey ==  8 ) yerr = 0.00210; // 21 um
00622               else if ( sizey ==  9 ) yerr = 0.00240; // 24 um
00623               else                    yerr = 0.00210; // 21um
00624             }
00625         } 
00626       else // EndCap
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; // 21 um
00638               else               yerr = 0.00075; // 7.5 um
00639             }
00640         }
00641 
00642     } // if ( !with_track_angle )
00643   else
00644     {
00645       //cout << "Track angles are known. We can use either errors from templates or the error parameterization from DB." << endl;
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             } // if ( qbin == 0 && inflate_errors )
00673           else
00674             {
00675               // Default errors
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             } // if ( qbin == 0 && inflate_errors ) else
00709 
00710         } //if ( UseErrorsFromTemplates_ )
00711       else 
00712         {
00713           //cout << endl << "Use errors from DB:" << endl;
00714           
00715           if ( edgex && edgey ) 
00716             {    
00717               //--- Both axes on the edge, no point in calling PixelErrorParameterization,       
00718               //--- just return the max errors on both.          
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         } //if ( UseErrorsFromTemplates_ ) else 
00740       
00741     } // if ( !with_track_angle ) else
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