CMS 3D CMS Logo

SiPixelDBErrorParametrization Class Reference

#include <CondTools/SiPixel/interface/SiPixelDBErrorParametrization.h>

List of all members.

Public Member Functions

std::pair< float, float > getError (GeomDetType::SubDetector pixelPart, int sizex, int sizey, float alpha, float beta, bool bigInX=false, bool bigInY=false)
std::pair< float, float > getError (const SiPixelCPEParmErrors *parmErrors, GeomDetType::SubDetector pixelPart, int sizex, int sizey, float alpha, float beta, bool bigInX=false, bool bigInY=false)
float index (int ind_subpart, int size, float alpha, float beta, bool big)
void setDBAccess (const edm::EventSetup &es)
 SiPixelDBErrorParametrization ()
 ~SiPixelDBErrorParametrization ()

Private Attributes

edm::ESHandle
< SiPixelCPEParmErrors
errorsH

Static Private Attributes

static const float * a_max [4] = {by_a_max, bx_a_max, fy_a_max, fx_a_max}
static const float * a_min [4] = {by_a_min, bx_a_min, fy_a_min, fx_a_min}
static const int alpha_bin_size [4] = {10, 1, 10, 1}
static const float * b_max [4] = {by_b_max, garbage, fy_b_max, fx_b_max}
static const float * b_min [4] = {by_b_min, garbage, fy_b_min, fx_b_min}
static const int beta_bin_size [4] = { 1, 10, 1, 10}
static const float bx_a_max [3] = {1.725, 1.675, 1.625}
static const float bx_a_min [3] = {1.525, 1.475, 1.425}
static const float by_a_max [6] = {1.67078, 1.67078, 1.67078, 1.67078, 1.67078, 1.67078}
static const float by_a_min [6] = {1.47078, 1.47078, 1.47078, 1.47078, 1.47078, 1.47078}
static const float by_b_max [6] = {0.50, 0.90, 1.05, 1.15, 1.20, 1.40}
static const float by_b_min [6] = {0.05, 0.15, 0.70, 0.95, 1.15, 1.20}
static const float errors_big_pix [4] = {0.0070, 0.0030, 0.0068, 0.0040}
static const float fx_a_max [2] = {0.285, 0.465}
static const float fx_a_min [2] = {0.165, 0.185}
static const float fx_b_max [2] = {999, 999}
static const float fx_b_min [2] = {998, 998}
static const float fy_a_max [2] = {999, 999}
static const float fy_a_min [2] = {998, 998}
static const float fy_b_max [2] = {0.39, 0.39}
static const float fy_b_min [2] = {0.31, 0.31}
static const int part_bin_size [4] = { 0, 240, 360, 400}
static const int size_bin_size [4] = {40, 40, 40, 40}
static const int size_max [4] = {5, 2, 1, 1}


Detailed Description

Definition at line 11 of file SiPixelDBErrorParametrization.h.


Constructor & Destructor Documentation

SiPixelDBErrorParametrization::SiPixelDBErrorParametrization (  ) 

Definition at line 46 of file SiPixelDBErrorParametrization.cc.

00046 {}

SiPixelDBErrorParametrization::~SiPixelDBErrorParametrization (  ) 

Definition at line 48 of file SiPixelDBErrorParametrization.cc.

00048 {}


Member Function Documentation

std::pair< float, float > SiPixelDBErrorParametrization::getError ( GeomDetType::SubDetector  pixelPart,
int  sizex,
int  sizey,
float  alpha,
float  beta,
bool  bigInX = false,
bool  bigInY = false 
)

Definition at line 91 of file SiPixelDBErrorParametrization.cc.

References element(), HLT_VtxMuL3::errors, errorsH, Exception, if(), index(), int, GeomDetEnumerators::PixelBarrel, and GeomDetEnumerators::PixelEndcap.

00095 {
00096         std::pair<float,float> element;
00097   std::pair<float,float> errors;
00098         
00099         switch (pixelPart)
00100         {
00101                 case GeomDetEnumerators::PixelBarrel:
00102                         element = std::pair<float,float>(index(1, sizex, alpha, beta, bigInX),  //1 -- Bx
00103                                                        index(0, sizey, alpha, beta, bigInY)); //0 -- By
00104                         break;
00105                 case GeomDetEnumerators::PixelEndcap:
00106                         element = std::pair<float,float>(index(3, sizex, alpha, beta, bigInX),  //3 -- Fx
00107                                                              index(2, sizey, alpha, beta, bigInY)); //2 -- Fy
00108                         break;
00109                 default:
00110                         throw cms::Exception("PixelDBErrorParametrization::getError")
00111                                 << "Non-pixel detector type !!!" ;
00112         }
00113         
00114         const SiPixelCPEParmErrors::DbVector & db_errors = errorsH->errors();
00115 
00116         if (bigInX && sizex == 1) errors.first  = element.first;
00117         else                      errors.first  = db_errors[(int)element.first].sigma;
00118         if (bigInY && sizey == 1) errors.second = element.second;
00119         else                      errors.second = db_errors[(int)element.second].sigma;
00120         
00121         return errors;
00122 }

std::pair< float, float > SiPixelDBErrorParametrization::getError ( const SiPixelCPEParmErrors parmErrors,
GeomDetType::SubDetector  pixelPart,
int  sizex,
int  sizey,
float  alpha,
float  beta,
bool  bigInX = false,
bool  bigInY = false 
)

Definition at line 56 of file SiPixelDBErrorParametrization.cc.

References element(), SiPixelCPEParmErrors::errors(), HLT_VtxMuL3::errors, Exception, if(), index(), int, GeomDetEnumerators::PixelBarrel, and GeomDetEnumerators::PixelEndcap.

Referenced by PxCPEdbReader::beginJob().

00061 {
00062         std::pair<float,float> element;
00063   std::pair<float,float> errors;
00064         
00065         switch (pixelPart)
00066         {
00067                 case GeomDetEnumerators::PixelBarrel:
00068                         element = std::pair<float,float>(index(1, sizex, alpha, beta, bigInX),  //1 -- Bx
00069                                                        index(0, sizey, alpha, beta, bigInY)); //0 -- By
00070                         break;
00071                 case GeomDetEnumerators::PixelEndcap:
00072                         element = std::pair<float,float>(index(3, sizex, alpha, beta, bigInX),  //3 -- Fx
00073                                                              index(2, sizey, alpha, beta, bigInY)); //2 -- Fy
00074                         break;
00075                 default:
00076                         throw cms::Exception("PixelDBErrorParametrization::getError")
00077                                 << "Non-pixel detector type !!!" ;
00078         }
00079         
00080         const SiPixelCPEParmErrors::DbVector & db_errors = parmErrors->errors();
00081 
00082         if (bigInX && sizex == 1) errors.first  = element.first;
00083         else                      errors.first  = db_errors[(int)element.first].sigma;
00084         if (bigInY && sizey == 1) errors.second = element.second;
00085         else                      errors.second = db_errors[(int)element.second].sigma;
00086         
00087         return errors;
00088 }

float SiPixelDBErrorParametrization::index ( int  ind_subpart,
int  size,
float  alpha,
float  beta,
bool  big 
)

Definition at line 126 of file SiPixelDBErrorParametrization.cc.

References a_max, a_min, alpha_bin_size, b_max, b_min, beta_bin_size, errors_big_pix, int, math_pi, min, part_bin_size, size_bin_size, and size_max.

Referenced by getError().

00127 {
00128         //This is a check for big pixels. If it passes, the code returns a given error and the function ends.
00129         if ( big && size == 1) return errors_big_pix[ind_subpart];
00130         
00131         int ind_size = std::min(size - 1, size_max[ind_subpart]);
00132 
00133         float alpha_rad = -999.9;
00134         float betap_rad = -999.9;
00135 
00136         int ind_alpha = -99;
00137         int ind_beta  = -99;
00138 
00139         float binw_a = -999.9;
00140         float binw_b = -999.9;
00141         int maxbin_a = -99;
00142         int maxbin_b = -99;
00143 
00144         betap_rad = fabs(math_pi/2.0 - beta);
00145         //We must take into account that Fx(subpart=3) has different alpha parametrization
00146         if(ind_subpart == 3) alpha_rad = fabs(math_pi/2.0 - alpha);
00147         else                 alpha_rad = fabs(alpha);
00148 
00149         //Sets the correct binning for alpha and beta based on whether in x or y
00150         if(ind_subpart == 0||ind_subpart == 2)
00151         {
00152                 binw_a = (a_max[ind_subpart][ind_size] - a_min[ind_subpart][ind_size])/2.0;
00153                 binw_b = (b_max[ind_subpart][ind_size] - b_min[ind_subpart][ind_size])/8.0;
00154                 maxbin_a = 3;
00155                 maxbin_b = 9;
00156         }
00157         else
00158         {
00159                 binw_a = (a_max[ind_subpart][ind_size] - a_min[ind_subpart][ind_size])/8.0;
00160                 binw_b = (b_max[ind_subpart][ind_size] - b_min[ind_subpart][ind_size])/2.0;
00161                 maxbin_a = 3;
00162                 maxbin_b = 9;
00163         }
00164 
00165         //Binning for alpha
00166         if      ( alpha_rad <  a_min[ind_subpart][ind_size]) ind_alpha = 0;
00167         else if ( alpha_rad >= a_max[ind_subpart][ind_size]) ind_alpha = maxbin_a;
00168         else      ind_alpha  = 1 + (int)((alpha_rad - a_min[ind_subpart][ind_size])/binw_a);
00169 
00170         //Binning for beta -- we need to account for Bx(subpart=1) having uneven binning
00171         if(ind_subpart == 1)
00172         {
00173                 if      (                     betap_rad <= 0.7 ) ind_beta = 0;
00174                 else if ( 0.7 <  betap_rad && betap_rad <= 1.0 ) ind_beta = 1;
00175                 else if ( 1.0 <  betap_rad && betap_rad <= 1.2 ) ind_beta = 2;
00176                 else if ( 1.2 <= betap_rad                     ) ind_beta = 3;
00177         }
00178         else if ( betap_rad <  b_min[ind_subpart][ind_size]) ind_beta = 0;
00179         else if ( betap_rad >= b_max[ind_subpart][ind_size]) ind_beta = maxbin_b;
00180         else      ind_beta   = 1 + (int)((betap_rad - b_min[ind_subpart][ind_size])/binw_b);
00181 
00182         //Index to be used to find error in database
00183         int index = part_bin_size[ind_subpart] + size_bin_size[ind_subpart] * ind_size + alpha_bin_size[ind_subpart] * ind_alpha + beta_bin_size[ind_subpart] * ind_beta;
00184         
00185         return index;
00186 }

void SiPixelDBErrorParametrization::setDBAccess ( const edm::EventSetup es  ) 

Definition at line 50 of file SiPixelDBErrorParametrization.cc.

References errorsH, and edm::EventSetup::get().

Referenced by PxCPEdbReader::beginJob().

00051 {
00052         es.get<SiPixelCPEParmErrorsRcd>().get(errorsH);
00053 }


Member Data Documentation

const float * SiPixelDBErrorParametrization::a_max = {by_a_max, bx_a_max, fy_a_max, fx_a_max} [static, private]

Definition at line 68 of file SiPixelDBErrorParametrization.h.

Referenced by index().

const float * SiPixelDBErrorParametrization::a_min = {by_a_min, bx_a_min, fy_a_min, fx_a_min} [static, private]

Definition at line 67 of file SiPixelDBErrorParametrization.h.

Referenced by index().

const int SiPixelDBErrorParametrization::alpha_bin_size = {10, 1, 10, 1} [static, private]

Definition at line 64 of file SiPixelDBErrorParametrization.h.

Referenced by index().

const float * SiPixelDBErrorParametrization::b_max = {by_b_max, garbage, fy_b_max, fx_b_max} [static, private]

Definition at line 70 of file SiPixelDBErrorParametrization.h.

Referenced by index().

const float * SiPixelDBErrorParametrization::b_min = {by_b_min, garbage, fy_b_min, fx_b_min} [static, private]

Definition at line 69 of file SiPixelDBErrorParametrization.h.

Referenced by index().

const int SiPixelDBErrorParametrization::beta_bin_size = { 1, 10, 1, 10} [static, private]

Definition at line 65 of file SiPixelDBErrorParametrization.h.

Referenced by index().

const float SiPixelDBErrorParametrization::bx_a_max = {1.725, 1.675, 1.625} [static, private]

Definition at line 42 of file SiPixelDBErrorParametrization.h.

const float SiPixelDBErrorParametrization::bx_a_min = {1.525, 1.475, 1.425} [static, private]

Definition at line 41 of file SiPixelDBErrorParametrization.h.

const float SiPixelDBErrorParametrization::by_a_max = {1.67078, 1.67078, 1.67078, 1.67078, 1.67078, 1.67078} [static, private]

Definition at line 50 of file SiPixelDBErrorParametrization.h.

const float SiPixelDBErrorParametrization::by_a_min = {1.47078, 1.47078, 1.47078, 1.47078, 1.47078, 1.47078} [static, private]

Definition at line 49 of file SiPixelDBErrorParametrization.h.

const float SiPixelDBErrorParametrization::by_b_max = {0.50, 0.90, 1.05, 1.15, 1.20, 1.40} [static, private]

Definition at line 52 of file SiPixelDBErrorParametrization.h.

const float SiPixelDBErrorParametrization::by_b_min = {0.05, 0.15, 0.70, 0.95, 1.15, 1.20} [static, private]

Definition at line 51 of file SiPixelDBErrorParametrization.h.

const float SiPixelDBErrorParametrization::errors_big_pix = {0.0070, 0.0030, 0.0068, 0.0040} [static, private]

Definition at line 59 of file SiPixelDBErrorParametrization.h.

Referenced by index().

edm::ESHandle<SiPixelCPEParmErrors> SiPixelDBErrorParametrization::errorsH [private]

Definition at line 39 of file SiPixelDBErrorParametrization.h.

Referenced by getError(), and setDBAccess().

const float SiPixelDBErrorParametrization::fx_a_max = {0.285, 0.465} [static, private]

Definition at line 45 of file SiPixelDBErrorParametrization.h.

const float SiPixelDBErrorParametrization::fx_a_min = {0.165, 0.185} [static, private]

Definition at line 44 of file SiPixelDBErrorParametrization.h.

const float SiPixelDBErrorParametrization::fx_b_max = {999, 999} [static, private]

Definition at line 47 of file SiPixelDBErrorParametrization.h.

const float SiPixelDBErrorParametrization::fx_b_min = {998, 998} [static, private]

Definition at line 46 of file SiPixelDBErrorParametrization.h.

const float SiPixelDBErrorParametrization::fy_a_max = {999, 999} [static, private]

Definition at line 55 of file SiPixelDBErrorParametrization.h.

const float SiPixelDBErrorParametrization::fy_a_min = {998, 998} [static, private]

Definition at line 54 of file SiPixelDBErrorParametrization.h.

const float SiPixelDBErrorParametrization::fy_b_max = {0.39, 0.39} [static, private]

Definition at line 57 of file SiPixelDBErrorParametrization.h.

const float SiPixelDBErrorParametrization::fy_b_min = {0.31, 0.31} [static, private]

Definition at line 56 of file SiPixelDBErrorParametrization.h.

const int SiPixelDBErrorParametrization::part_bin_size = { 0, 240, 360, 400} [static, private]

Definition at line 62 of file SiPixelDBErrorParametrization.h.

Referenced by index().

const int SiPixelDBErrorParametrization::size_bin_size = {40, 40, 40, 40} [static, private]

Definition at line 63 of file SiPixelDBErrorParametrization.h.

Referenced by index().

const int SiPixelDBErrorParametrization::size_max = {5, 2, 1, 1} [static, private]

Definition at line 60 of file SiPixelDBErrorParametrization.h.

Referenced by index().


The documentation for this class was generated from the following files:
Generated on Tue Jun 9 18:31:53 2009 for CMSSW by  doxygen 1.5.4