CMS 3D CMS Logo

SiPixelGaussianSmearingRecHitConverterAlgorithm Class Reference

#include <FastSimulation/TrackingRecHitProducer/interface/SiPixelGaussianSmearingRecHitConverterAlgorithm.h>

List of all members.

Public Member Functions

LocalError getError ()
double getErrorX ()
double getErrorY ()
double getErrorZ ()
unsigned int getPixelMultiplicityAlpha ()
unsigned int getPixelMultiplicityBeta ()
Local3DPoint getPosition ()
double getPositionX ()
double getPositionY ()
double getPositionZ ()
 SiPixelGaussianSmearingRecHitConverterAlgorithm (const edm::ParameterSet &pset, GeomDetType::SubDetector pixelPart, std::vector< TH1F * > &alphaMultiplicityCumulativeProbabilities, std::vector< TH1F * > &betaMultiplicityCumulativeProbabilities, TFile *pixelResolutionFile, const RandomEngine *engine)
void smearHit (const PSimHit &simHit, const PixelGeomDetUnit *detUnit, const double boundX, const double boundY)
virtual ~SiPixelGaussianSmearingRecHitConverterAlgorithm ()

Private Member Functions

bool isFlipped (const PixelGeomDetUnit *theDet) const

Private Attributes

bool negativeErrorProtection
PixelErrorParametrizationpixelError
edm::ParameterSet pset_
const RandomEnginerandom
double resAlpha_binMin
unsigned int resAlpha_binN
double resAlpha_binWidth
double resBeta_binMin
unsigned int resBeta_binN
double resBeta_binWidth
std::map< unsigned, const
HistogramGenerator * > 
theAlphaHistos
std::vector< TH1F * > theAlphaMultiplicityCumulativeProbabilities
std::map< unsigned, const
HistogramGenerator * > 
theBetaHistos
std::vector< TH1F * > theBetaMultiplicityCumulativeProbabilities
LocalError theError
double theErrorX
double theErrorY
double theErrorZ
unsigned int theLayer
unsigned int thePixelMultiplicityAlpha
unsigned int thePixelMultiplicityBeta
GeomDetType::SubDetector thePixelPart
TFile * thePixelResolutionFile
Local3DPoint thePosition
double thePositionX
double thePositionY
double thePositionZ
bool useCMSSWPixelParameterization


Detailed Description

Definition at line 36 of file SiPixelGaussianSmearingRecHitConverterAlgorithm.h.


Constructor & Destructor Documentation

SiPixelGaussianSmearingRecHitConverterAlgorithm::SiPixelGaussianSmearingRecHitConverterAlgorithm ( const edm::ParameterSet pset,
GeomDetType::SubDetector  pixelPart,
std::vector< TH1F * > &  alphaMultiplicityCumulativeProbabilities,
std::vector< TH1F * > &  betaMultiplicityCumulativeProbabilities,
TFile *  pixelResolutionFile,
const RandomEngine engine 
) [explicit]

Definition at line 35 of file SiPixelGaussianSmearingRecHitConverterAlgorithm.cc.

References edm::ParameterSet::getParameter(), python::StorageManager_cfg::maxSize, GeomDetEnumerators::PixelBarrel, GeomDetEnumerators::PixelEndcap, pixelError, HLT_VtxMuL3::PixelErrorParametrization, pset_, random, resAlpha_binMin, resAlpha_binN, resAlpha_binWidth, resBeta_binMin, resBeta_binN, resBeta_binWidth, theAlphaHistos, theAlphaMultiplicityCumulativeProbabilities, theBetaHistos, theBetaMultiplicityCumulativeProbabilities, thePixelPart, thePixelResolutionFile, and useCMSSWPixelParameterization.

00042 :
00043   pset_(pset),
00044   thePixelPart(pixelPart),
00045   theAlphaMultiplicityCumulativeProbabilities(alphaMultiplicityCumulativeProbabilities),
00046   theBetaMultiplicityCumulativeProbabilities(betaMultiplicityCumulativeProbabilities),
00047   thePixelResolutionFile(pixelResolutionFile),
00048   random(engine)
00049 {
00050   // Switch between old (ORCA) and new (CMSSW) pixel parameterization
00051   useCMSSWPixelParameterization = pset.getParameter<bool>("UseCMSSWPixelParametrization");
00052 
00053   if(useCMSSWPixelParameterization) {
00054     if( thePixelPart == GeomDetEnumerators::PixelBarrel ) {
00055       // Resolution Barrel    
00056       resAlpha_binMin   = 
00057         pset.getParameter<double>("AlphaBarrel_BinMinNew"  );
00058       resAlpha_binWidth = 
00059         pset.getParameter<double>("AlphaBarrel_BinWidthNew");
00060       resAlpha_binN     = 
00061         pset.getParameter<int>("AlphaBarrel_BinNNew"       );
00062       resBeta_binMin    = 
00063         pset.getParameter<double>("BetaBarrel_BinMinNew"   );
00064       resBeta_binWidth  = 
00065         pset.getParameter<double>("BetaBarrel_BinWidthNew" );
00066       resBeta_binN      = 
00067         pset.getParameter<int>(   "BetaBarrel_BinNNew"     );
00068       //
00069     } else if( thePixelPart == GeomDetEnumerators::PixelEndcap ) {
00070       // Resolution Forward
00071       resAlpha_binMin   = 
00072         pset.getParameter<double>("AlphaForward_BinMinNew"  );
00073       resAlpha_binWidth = 
00074         pset.getParameter<double>("AlphaForward_BinWidthNew");
00075       resAlpha_binN     = 
00076         pset.getParameter<int>("AlphaBarrel_BinNNew"        );
00077       resBeta_binMin    = 
00078         pset.getParameter<double>("BetaForward_BinMinNew"   );
00079       resBeta_binWidth  = 
00080         pset.getParameter<double>("BetaForward_BinWidthNew" );
00081       resBeta_binN      = 
00082         pset.getParameter<int>(   "BetaBarrel_BinNNew"      );
00083     }
00084   } else {
00085         if( thePixelPart == GeomDetEnumerators::PixelBarrel ) {
00086       // Resolution Barrel    
00087       resAlpha_binMin   = 
00088         pset.getParameter<double>("AlphaBarrel_BinMin"  );
00089       resAlpha_binWidth = 
00090         pset.getParameter<double>("AlphaBarrel_BinWidth");
00091       resAlpha_binN     = 
00092         pset.getParameter<int>("AlphaBarrel_BinN"       );
00093       resBeta_binMin    = 
00094         pset.getParameter<double>("BetaBarrel_BinMin"   );
00095       resBeta_binWidth  = 
00096         pset.getParameter<double>("BetaBarrel_BinWidth" );
00097       resBeta_binN      = 
00098         pset.getParameter<int>(   "BetaBarrel_BinN"     );
00099       //
00100     } else if( thePixelPart == GeomDetEnumerators::PixelEndcap ) {
00101       // Resolution Forward
00102       resAlpha_binMin   = 
00103         pset.getParameter<double>("AlphaForward_BinMin"  );
00104       resAlpha_binWidth = 
00105         pset.getParameter<double>("AlphaForward_BinWidth");
00106       resAlpha_binN     = 
00107         pset.getParameter<int>("AlphaBarrel_BinN"        );
00108       resBeta_binMin    = 
00109         pset.getParameter<double>("BetaForward_BinMin"   );
00110       resBeta_binWidth  = 
00111         pset.getParameter<double>("BetaForward_BinWidth" );
00112       resBeta_binN      = 
00113         pset.getParameter<int>(   "BetaBarrel_BinN"      );
00114     }
00115   }
00116   // Initialize PixelErrorParametrization (time consuming!)
00117   pixelError = new PixelErrorParametrization(pset_);
00118 
00119   // Initialize the histos once and for all, and prepare the random generation
00120   for ( unsigned alphaHistBin=1; alphaHistBin<=resAlpha_binN; ++alphaHistBin ) {
00121     unsigned int maxSize;
00122     if(useCMSSWPixelParameterization)
00123       maxSize = theAlphaMultiplicityCumulativeProbabilities.size() / 2;
00124     else
00125       maxSize = theAlphaMultiplicityCumulativeProbabilities.size();
00126     for ( unsigned alphaMultiplicity=1; 
00127           alphaMultiplicity<=maxSize;
00128           ++alphaMultiplicity ) {
00129       unsigned int alphaHistN = (resAlpha_binWidth != 0. ?
00130                                  100 * alphaHistBin
00131                                  + 10
00132                                  + alphaMultiplicity
00133                                  :
00134                                  1110
00135                                  + alphaMultiplicity);    //
00136       theAlphaHistos[alphaHistN] = new HistogramGenerator(
00137         (TH1F*) thePixelResolutionFile->Get(  Form( "h%u" , alphaHistN ) ),
00138         random);
00139       // Fill also big pixels if new parametrization is used. Their code is 10000 + histogram number
00140       if(useCMSSWPixelParameterization) {
00141         theAlphaHistos[alphaHistN+10000] = new HistogramGenerator(
00142           (TH1F*) thePixelResolutionFile->Get(  Form( "h%ub" , alphaHistN ) ),
00143           random);
00144       }
00145     }
00146   }
00147 
00148 
00149   //
00150   for ( unsigned betaHistBin=1; betaHistBin<=resBeta_binN; ++betaHistBin ) {
00151     unsigned int maxSize;
00152     if(useCMSSWPixelParameterization)
00153       maxSize = theBetaMultiplicityCumulativeProbabilities.size() / 2;
00154     else
00155       maxSize = theBetaMultiplicityCumulativeProbabilities.size();
00156     for ( unsigned betaMultiplicity=1; 
00157           betaMultiplicity<=maxSize;
00158           ++betaMultiplicity ) {
00159       unsigned int betaHistN = (resBeta_binWidth != 0. ?
00160                                 100 * betaHistBin
00161                                 + betaMultiplicity
00162                                 :
00163                                 1100 + betaMultiplicity);    //
00164       theBetaHistos[betaHistN] = new HistogramGenerator(
00165         (TH1F*) thePixelResolutionFile->Get(  Form( "h%u" , betaHistN  ) ),
00166         random);
00167       // Fill also big pixels if new parametrization is used. Their code is 10000 + histogram number
00168       if(useCMSSWPixelParameterization)
00169         theBetaHistos[betaHistN+10000] = new HistogramGenerator(
00170           (TH1F*) thePixelResolutionFile->Get(  Form( "h%ub" , betaHistN  ) ),
00171           random);
00172     }
00173   }
00174 }

SiPixelGaussianSmearingRecHitConverterAlgorithm::~SiPixelGaussianSmearingRecHitConverterAlgorithm (  )  [virtual]

Definition at line 176 of file SiPixelGaussianSmearingRecHitConverterAlgorithm.cc.

References it, pixelError, theAlphaHistos, and theBetaHistos.

00177 {
00178   
00179   // Some cleaning
00180   delete pixelError;
00181 
00182   std::map<unsigned,const HistogramGenerator*>::const_iterator it;
00183 
00184   for ( it=theAlphaHistos.begin(); it!=theAlphaHistos.end(); ++it ) 
00185     delete it->second;
00186 
00187   for ( it=theBetaHistos.begin(); it!=theBetaHistos.end(); ++it ) 
00188     delete it->second;
00189 
00190   theAlphaHistos.clear();
00191   theBetaHistos.clear();
00192 
00193 }


Member Function Documentation

LocalError SiPixelGaussianSmearingRecHitConverterAlgorithm::getError (  )  [inline]

Definition at line 55 of file SiPixelGaussianSmearingRecHitConverterAlgorithm.h.

References theError.

Referenced by SiTrackerGaussianSmearingRecHitConverter::gaussianSmearing().

00055 {return theError;}

double SiPixelGaussianSmearingRecHitConverterAlgorithm::getErrorX (  )  [inline]

Definition at line 56 of file SiPixelGaussianSmearingRecHitConverterAlgorithm.h.

References theErrorX.

00056 {return theErrorX;}

double SiPixelGaussianSmearingRecHitConverterAlgorithm::getErrorY (  )  [inline]

Definition at line 57 of file SiPixelGaussianSmearingRecHitConverterAlgorithm.h.

References theErrorY.

00057 {return theErrorY;}

double SiPixelGaussianSmearingRecHitConverterAlgorithm::getErrorZ (  )  [inline]

Definition at line 58 of file SiPixelGaussianSmearingRecHitConverterAlgorithm.h.

References theErrorZ.

00058 {return theErrorZ;}

unsigned int SiPixelGaussianSmearingRecHitConverterAlgorithm::getPixelMultiplicityAlpha (  )  [inline]

Definition at line 59 of file SiPixelGaussianSmearingRecHitConverterAlgorithm.h.

References thePixelMultiplicityAlpha.

Referenced by SiTrackerGaussianSmearingRecHitConverter::gaussianSmearing().

00059 {return thePixelMultiplicityAlpha;}

unsigned int SiPixelGaussianSmearingRecHitConverterAlgorithm::getPixelMultiplicityBeta (  )  [inline]

Definition at line 60 of file SiPixelGaussianSmearingRecHitConverterAlgorithm.h.

References thePixelMultiplicityBeta.

Referenced by SiTrackerGaussianSmearingRecHitConverter::gaussianSmearing().

00060 {return thePixelMultiplicityBeta;}

Local3DPoint SiPixelGaussianSmearingRecHitConverterAlgorithm::getPosition (  )  [inline]

Definition at line 51 of file SiPixelGaussianSmearingRecHitConverterAlgorithm.h.

References thePosition.

Referenced by SiTrackerGaussianSmearingRecHitConverter::gaussianSmearing().

00051 {return thePosition;}

double SiPixelGaussianSmearingRecHitConverterAlgorithm::getPositionX (  )  [inline]

Definition at line 52 of file SiPixelGaussianSmearingRecHitConverterAlgorithm.h.

References thePositionX.

00052 {return thePositionX;}

double SiPixelGaussianSmearingRecHitConverterAlgorithm::getPositionY (  )  [inline]

Definition at line 53 of file SiPixelGaussianSmearingRecHitConverterAlgorithm.h.

References thePositionY.

00053 {return thePositionY;}

double SiPixelGaussianSmearingRecHitConverterAlgorithm::getPositionZ (  )  [inline]

Definition at line 54 of file SiPixelGaussianSmearingRecHitConverterAlgorithm.h.

References thePositionZ.

00054 {return thePositionZ;}

bool SiPixelGaussianSmearingRecHitConverterAlgorithm::isFlipped ( const PixelGeomDetUnit theDet  )  const [private]

Definition at line 495 of file SiPixelGaussianSmearingRecHitConverterAlgorithm.cc.

References PV3DBase< T, PVType, FrameType >::perp(), GeomDet::surface(), tmp1, tmp2, and Surface::toGlobal().

Referenced by smearHit().

00495                                                                                                     {
00496   // Check the relative position of the local +/- z in global coordinates.
00497   float tmp1 = theDet->surface().toGlobal(Local3DPoint(0.,0.,0.)).perp();
00498   float tmp2 = theDet->surface().toGlobal(Local3DPoint(0.,0.,1.)).perp();
00499   //  std::cout << " 1: " << tmp1 << " 2: " << tmp2 << std::endl;
00500   if ( tmp2<tmp1 ) return true;
00501   else return false;    
00502 }

void SiPixelGaussianSmearingRecHitConverterAlgorithm::smearHit ( const PSimHit simHit,
const PixelGeomDetUnit detUnit,
const double  boundX,
const double  boundY 
)

Definition at line 195 of file SiPixelGaussianSmearingRecHitConverterAlgorithm.cc.

References DeDxTools::beta(), counter(), GenMuonPlsPt100GeV_cfg::cout, e, lat::endl(), PSimHit::entryPoint(), PSimHit::exitPoint(), first, RandomEngine::flatShoot(), PixelErrorParametrization::getError(), int, isFlipped(), PSimHit::localDirection(), PSimHit::localPosition(), PSimHit::momentumAtEntry(), PixelTopology::ncolumns(), PixelTopology::nrows(), PI, PixelTopology::pitch(), GeomDetEnumerators::PixelBarrel, pixelError, random, resAlpha_binMin, resAlpha_binN, resAlpha_binWidth, resBeta_binMin, resBeta_binN, resBeta_binWidth, edm::second(), PixelGeomDetUnit::specificTopology(), funct::sqrt(), theAlphaHistos, theAlphaMultiplicityCumulativeProbabilities, theBetaHistos, theBetaMultiplicityCumulativeProbabilities, theError, theErrorX, theErrorY, theErrorZ, thePixelMultiplicityAlpha, thePixelMultiplicityBeta, thePixelPart, thePosition, thePositionX, thePositionY, thePositionZ, Vector3DBase< T, FrameTag >::unit(), useCMSSWPixelParameterization, PV3DBase< T, PVType, FrameType >::x(), PV3DBase< T, PVType, FrameType >::y(), and PV3DBase< T, PVType, FrameType >::z().

Referenced by SiTrackerGaussianSmearingRecHitConverter::gaussianSmearing().

00200 {
00201 
00202 #ifdef FAMOS_DEBUG
00203   std::cout << " Pixel smearing in " << thePixelPart 
00204             << std::endl;
00205 #endif
00206   //
00207   // at the beginning the position is the Local Point in the local pixel module reference frame
00208   // same code as in PixelCPEBase
00209   LocalVector localDir = simHit.momentumAtEntry().unit();
00210   float locx = localDir.x();
00211   float locy = localDir.y();
00212   float locz = localDir.z();
00213 
00214   //
00215   bool hasBigPixelInX = false;
00216   bool hasBigPixelInY = false;
00217 
00218   if(useCMSSWPixelParameterization) {
00219     // If the sim track crosses a region in which there are big pixels,
00220     // then we set to true the variables above
00221 
00222     // Get the topology of the pixel module
00223     const PixelTopology* theSpecificTopology = &(detUnit->specificTopology());
00224     RectangularPixelTopology rectPixelTopology(theSpecificTopology->nrows(), 
00225                                                theSpecificTopology->ncolumns(), 
00226                                                theSpecificTopology->pitch().first, 
00227                                                theSpecificTopology->pitch().second);
00228     
00229     // Get the rows and columns of entry and exit points
00230     // FIXME - these are not guaranteed to be the same as the cluster limits (as they should be)
00231     const int firstPixelInX = int(rectPixelTopology.pixel(simHit.entryPoint()).first);
00232     const int firstPixelInY = int(rectPixelTopology.pixel(simHit.entryPoint()).second);
00233     const int lastPixelInX = int(rectPixelTopology.pixel(simHit.exitPoint()).first);
00234     const int lastPixelInY = int(rectPixelTopology.pixel(simHit.exitPoint()).second);
00235     
00236     // Check if there is a big pixel inside and set hasBigPixelInX and hasBigPixelInY accordingly
00237     // This function only works if first <= last
00238     if(rectPixelTopology.containsBigPixelInX(firstPixelInX < lastPixelInX ? firstPixelInX : lastPixelInX,
00239                                              firstPixelInX > lastPixelInX ? firstPixelInX : lastPixelInX))
00240       hasBigPixelInX = true;
00241     if(rectPixelTopology.containsBigPixelInY(firstPixelInY < lastPixelInY ? firstPixelInY : lastPixelInY,
00242                                              firstPixelInY > lastPixelInY ? firstPixelInY : lastPixelInY))
00243       hasBigPixelInY = true;
00244 #ifdef FAMOS_DEBUG
00245     std::cout << " Simhit first pixel in X = " << (firstPixelInX < lastPixelInX ? firstPixelInX : lastPixelInX)
00246               << " last pixel in X = " << (firstPixelInX > lastPixelInX ? firstPixelInX : lastPixelInX)
00247               << " has big pixel in X is " << (hasBigPixelInX ? " true" : " false")
00248               << std::endl;
00249     std::cout << " Simhit first pixel in Y = " << (firstPixelInY < lastPixelInY ? firstPixelInY : lastPixelInY)
00250               << " last pixel in Y = " << (firstPixelInY > lastPixelInY ? firstPixelInY : lastPixelInY)
00251               << " has big pixel in Y is " << (hasBigPixelInY ? " true" : " false")
00252               << std::endl;
00253 #endif
00254   }
00255 
00256   // alpha: angle with respect to local x axis in local (x,z) plane
00257   float alpha = std::acos(locx/std::sqrt(locx*locx+locz*locz));
00258   if ( isFlipped( detUnit ) ) { // &&& check for FPIX !!!
00259 #ifdef FAMOS_DEBUG
00260     std::cout << " isFlipped " << std::endl;
00261 #endif
00262     alpha = PI - alpha ;
00263   }
00264   // beta: angle with respect to local y axis in local (y,z) plane
00265   float beta = std::acos(locy/std::sqrt(locy*locy+locz*locz));
00266   
00267   // look old FAMOS: FamosGeneric/FamosTracker/src/FamosPixelErrorParametrization
00268   float alphaToBeUsedForRootFiles = alpha;
00269   float betaToBeUsedForRootFiles  = beta;
00270   if( thePixelPart == GeomDetEnumerators::PixelBarrel ) { // BARREL
00271     alphaToBeUsedForRootFiles = PI/2. - alpha;
00272     betaToBeUsedForRootFiles  = fabs( PI/2. - beta );
00273   } else { // FORWARD
00274     betaToBeUsedForRootFiles = fabs( PI/2. - beta );
00275     alphaToBeUsedForRootFiles  = fabs( PI/2. - alpha );    
00276   }
00277   //
00278 #ifdef FAMOS_DEBUG
00279   std::cout << " Local Direction " << simHit.localDirection()
00280             << " alpha(x) = " << alpha
00281             << " beta(y) = "  << beta
00282             << " alpha for root files = " << alphaToBeUsedForRootFiles
00283             << " beta for root files = "  << betaToBeUsedForRootFiles
00284             << std::endl;
00285 #endif
00286 
00287   // Generate alpha and beta multiplicity
00288   unsigned int alphaMultiplicity = 0;
00289   unsigned int betaMultiplicity  = 0;
00290   // random multiplicity for alpha and beta
00291   double alphaProbability = random->flatShoot();
00292   double betaProbability  = random->flatShoot();
00293 
00294   // search which multiplicity correspond
00295   int alphaBin = 
00296     theAlphaMultiplicityCumulativeProbabilities.front()->GetXaxis()->FindFixBin(alphaToBeUsedForRootFiles);
00297   int betaBin = 
00298     theBetaMultiplicityCumulativeProbabilities.front()->GetXaxis()->FindFixBin(betaToBeUsedForRootFiles);
00299 
00300   // protection against out-of-range (undeflows and overflows)
00301   if( alphaBin == 0 ) alphaBin = 1;
00302   if( alphaBin > theAlphaMultiplicityCumulativeProbabilities.front()->GetNbinsX() ) 
00303     alphaBin = theAlphaMultiplicityCumulativeProbabilities.front()->GetNbinsX();
00304   if( betaBin == 0 ) betaBin = 1;
00305   if( betaBin > theBetaMultiplicityCumulativeProbabilities.front()->GetNbinsX() )   
00306     betaBin = theBetaMultiplicityCumulativeProbabilities.front()->GetNbinsX();
00307   
00308   unsigned int iMult, multSize;
00309   const unsigned int maxMultX = theAlphaMultiplicityCumulativeProbabilities.size();
00310   const unsigned int maxMultY = theBetaMultiplicityCumulativeProbabilities.size();
00311   if(useCMSSWPixelParameterization) {
00312     if(hasBigPixelInX) {     // Big pixels: second half of histograms vector
00313       iMult = maxMultX / 2;
00314       multSize = maxMultX;
00315     } else {                 // Normal pixels: first half of histograms vector
00316       iMult = 0;
00317       multSize = maxMultX / 2;
00318     }
00319   } else {
00320     iMult = 0;
00321     multSize = maxMultX;
00322   }
00323   for(/* void */; iMult < multSize; iMult++) {
00324     if(alphaProbability < theAlphaMultiplicityCumulativeProbabilities[iMult]->GetBinContent(alphaBin) ) {
00325       alphaMultiplicity = iMult+1;
00326       break;
00327     }
00328   }
00329   
00330   if(useCMSSWPixelParameterization) {
00331     if(hasBigPixelInY) {     // Big pixels: second half of histograms vector
00332       iMult = maxMultY / 2;
00333       multSize = maxMultY;
00334     } else {                 // Normal pixels: first half of histograms vector
00335       iMult = 0;
00336       multSize = maxMultY / 2;
00337     }
00338   } else {
00339     iMult = 0;
00340     multSize = maxMultY;
00341   }
00342   for(/* void */; iMult < multSize; iMult++) {
00343     if(betaProbability < theBetaMultiplicityCumulativeProbabilities[iMult]->GetBinContent(betaBin) ) {
00344       betaMultiplicity = iMult+1;
00345       break;
00346     }
00347   }
00348 
00349   // Correct multiplicity for big pixels
00350   if(hasBigPixelInX)
00351     alphaMultiplicity -= maxMultX / 2;
00352   if(hasBigPixelInY)
00353     betaMultiplicity -= maxMultY / 2;
00354 
00355   // protection against 0 or max multiplicity
00356   if(useCMSSWPixelParameterization) {
00357     if( alphaMultiplicity == 0 || alphaMultiplicity > maxMultX / 2 ) 
00358       alphaMultiplicity = maxMultX / 2;
00359     if( betaMultiplicity  == 0 || betaMultiplicity  > maxMultY / 2  ) 
00360       betaMultiplicity  = maxMultY / 2;
00361   } else {
00362     if( alphaMultiplicity == 0 || alphaMultiplicity > maxMultX ) 
00363       alphaMultiplicity = maxMultX;
00364     if( betaMultiplicity  == 0 || betaMultiplicity  > maxMultY ) 
00365       betaMultiplicity  = maxMultY;
00366   }
00367   //
00368   
00369 //
00370 #ifdef FAMOS_DEBUG
00371   std::cout << " Multiplicity set to"
00372             << "\talpha = " << alphaMultiplicity
00373             << "\tbeta = "  << betaMultiplicity
00374             << "\n"
00375             << "  from random probability"
00376             << "\talpha = " << alphaProbability
00377             << "\tbeta = "  << betaProbability
00378             << "\n"
00379             << "  taken from bin         "
00380             << "\talpha = " << alphaBin
00381             << "\tbeta = "  << betaBin
00382             << std::endl;       
00383 #endif
00384   
00385   // Compute pixel errors
00386   std::pair<float,float> theErrors = pixelError->getError( thePixelPart ,
00387                                                           (int)alphaMultiplicity , (int)betaMultiplicity ,
00388                                                            alpha                 , beta                  ,
00389                                                            hasBigPixelInX        , hasBigPixelInY         );
00390   // define private mebers --> Errors
00391   theErrorX = theErrors.first;  // PixelErrorParametrization returns sigma, not sigma^2
00392   theErrorY = theErrors.second; // PixelErrorParametrization returns sigma, not sigma^2
00393   theErrorZ = 1e-8; // 1 um means zero
00394   theError = LocalError( theErrorX*theErrorX, 0., theErrorY*theErrorY);
00395   // Local Error is 2D: (xx,xy,yy), square of sigma in first an third position 
00396   // as for resolution matrix
00397   //
00398 #ifdef FAMOS_DEBUG
00399   std::cout << " Pixel Errors "
00400             << "\talpha(x) = " << theErrorX
00401             << "\tbeta(y) = "  << theErrorY
00402             << std::endl;       
00403 #endif
00404   
00405   // 
00406   // Generate position
00407   // get resolution histograms
00408   int alphaHistBin = (int)( ( alphaToBeUsedForRootFiles - resAlpha_binMin ) / resAlpha_binWidth + 1 );
00409   int betaHistBin  = (int)( ( betaToBeUsedForRootFiles  - resBeta_binMin )  / resBeta_binWidth + 1 );
00410   // protection against out-of-range (undeflows and overflows)
00411   if( alphaHistBin < 1 ) alphaHistBin = 1; 
00412   if( betaHistBin  < 1 ) betaHistBin  = 1; 
00413   if( alphaHistBin > (int)resAlpha_binN ) alphaHistBin = (int)resAlpha_binN; 
00414   if( betaHistBin  > (int)resBeta_binN  ) betaHistBin  = (int)resBeta_binN; 
00415   //  
00416   unsigned int alphaHistN = (resAlpha_binWidth != 0. ?
00417                              100 * alphaHistBin
00418                              + 10
00419                              + alphaMultiplicity
00420                              :
00421                              1110
00422                              + alphaMultiplicity);    //
00423   //
00424   unsigned int betaHistN = (resBeta_binWidth != 0. ?
00425                             100 * betaHistBin
00426                             + betaMultiplicity
00427                             :
00428                             1100 + betaMultiplicity);    //
00429   //
00430   if(hasBigPixelInX)
00431     alphaHistN += 10000;
00432   if(hasBigPixelInY)
00433     betaHistN += 10000;
00434 #ifdef FAMOS_DEBUG
00435   std::cout << " Resolution histograms chosen "
00436             << "\talpha = " << alphaHistN
00437             << "\tbeta = "  << betaHistN
00438             << std::endl;       
00439 #endif
00440 
00441   unsigned int counter = 0;
00442   do {
00443     //
00444     // Smear the hit Position
00445     thePositionX = theAlphaHistos[alphaHistN]->generate();
00446     thePositionY = theBetaHistos[betaHistN]->generate();
00447     //  thePositionX = theAlphaHistos[alphaHistN]->getHisto()->GetRandom();
00448     //  thePositionY = theBetaHistos[betaHistN]->getHisto()->GetRandom();
00449     thePositionZ = 0.0; // set at the centre of the active area
00450     thePosition = 
00451       Local3DPoint(simHit.localPosition().x() + thePositionX , 
00452                    simHit.localPosition().y() + thePositionY , 
00453                    simHit.localPosition().z() + thePositionZ );
00454 #ifdef FAMOS_DEBUG
00455     std::cout << " Detector bounds: "
00456               << "\t\tx = " << boundX
00457               << "\ty = " << boundY
00458               << std::endl;
00459     std::cout << " Generated local position "
00460               << "\tx = " << thePosition.x()
00461               << "\ty = " << thePosition.y()
00462               << std::endl;       
00463 #endif  
00464     counter++;
00465     if(counter > 20) {
00466       thePosition = Local3DPoint(simHit.localPosition().x(), 
00467                                  simHit.localPosition().y(), 
00468                                  simHit.localPosition().z());
00469       break;
00470     }
00471   } while(fabs(thePosition.x()) > boundX  || fabs(thePosition.y()) > boundY);
00472   
00473 
00474   
00475   // define private mebers --> Multiplicities
00476   thePixelMultiplicityAlpha = alphaMultiplicity;
00477   thePixelMultiplicityBeta  = betaMultiplicity;
00478   //
00479  
00480 }


Member Data Documentation

bool SiPixelGaussianSmearingRecHitConverterAlgorithm::negativeErrorProtection [private]

Definition at line 72 of file SiPixelGaussianSmearingRecHitConverterAlgorithm.h.

PixelErrorParametrization* SiPixelGaussianSmearingRecHitConverterAlgorithm::pixelError [private]

Definition at line 105 of file SiPixelGaussianSmearingRecHitConverterAlgorithm.h.

Referenced by SiPixelGaussianSmearingRecHitConverterAlgorithm(), smearHit(), and ~SiPixelGaussianSmearingRecHitConverterAlgorithm().

edm::ParameterSet SiPixelGaussianSmearingRecHitConverterAlgorithm::pset_ [private]

Definition at line 80 of file SiPixelGaussianSmearingRecHitConverterAlgorithm.h.

Referenced by SiPixelGaussianSmearingRecHitConverterAlgorithm().

const RandomEngine* SiPixelGaussianSmearingRecHitConverterAlgorithm::random [private]

Definition at line 108 of file SiPixelGaussianSmearingRecHitConverterAlgorithm.h.

Referenced by SiPixelGaussianSmearingRecHitConverterAlgorithm(), and smearHit().

double SiPixelGaussianSmearingRecHitConverterAlgorithm::resAlpha_binMin [private]

Definition at line 75 of file SiPixelGaussianSmearingRecHitConverterAlgorithm.h.

Referenced by SiPixelGaussianSmearingRecHitConverterAlgorithm(), and smearHit().

unsigned int SiPixelGaussianSmearingRecHitConverterAlgorithm::resAlpha_binN [private]

Definition at line 76 of file SiPixelGaussianSmearingRecHitConverterAlgorithm.h.

Referenced by SiPixelGaussianSmearingRecHitConverterAlgorithm(), and smearHit().

double SiPixelGaussianSmearingRecHitConverterAlgorithm::resAlpha_binWidth [private]

Definition at line 75 of file SiPixelGaussianSmearingRecHitConverterAlgorithm.h.

Referenced by SiPixelGaussianSmearingRecHitConverterAlgorithm(), and smearHit().

double SiPixelGaussianSmearingRecHitConverterAlgorithm::resBeta_binMin [private]

Definition at line 77 of file SiPixelGaussianSmearingRecHitConverterAlgorithm.h.

Referenced by SiPixelGaussianSmearingRecHitConverterAlgorithm(), and smearHit().

unsigned int SiPixelGaussianSmearingRecHitConverterAlgorithm::resBeta_binN [private]

Definition at line 78 of file SiPixelGaussianSmearingRecHitConverterAlgorithm.h.

Referenced by SiPixelGaussianSmearingRecHitConverterAlgorithm(), and smearHit().

double SiPixelGaussianSmearingRecHitConverterAlgorithm::resBeta_binWidth [private]

Definition at line 77 of file SiPixelGaussianSmearingRecHitConverterAlgorithm.h.

Referenced by SiPixelGaussianSmearingRecHitConverterAlgorithm(), and smearHit().

std::map<unsigned,const HistogramGenerator*> SiPixelGaussianSmearingRecHitConverterAlgorithm::theAlphaHistos [private]

Definition at line 87 of file SiPixelGaussianSmearingRecHitConverterAlgorithm.h.

Referenced by SiPixelGaussianSmearingRecHitConverterAlgorithm(), smearHit(), and ~SiPixelGaussianSmearingRecHitConverterAlgorithm().

std::vector<TH1F*> SiPixelGaussianSmearingRecHitConverterAlgorithm::theAlphaMultiplicityCumulativeProbabilities [private]

Definition at line 84 of file SiPixelGaussianSmearingRecHitConverterAlgorithm.h.

Referenced by SiPixelGaussianSmearingRecHitConverterAlgorithm(), and smearHit().

std::map<unsigned,const HistogramGenerator*> SiPixelGaussianSmearingRecHitConverterAlgorithm::theBetaHistos [private]

Definition at line 88 of file SiPixelGaussianSmearingRecHitConverterAlgorithm.h.

Referenced by SiPixelGaussianSmearingRecHitConverterAlgorithm(), smearHit(), and ~SiPixelGaussianSmearingRecHitConverterAlgorithm().

std::vector<TH1F*> SiPixelGaussianSmearingRecHitConverterAlgorithm::theBetaMultiplicityCumulativeProbabilities [private]

Definition at line 85 of file SiPixelGaussianSmearingRecHitConverterAlgorithm.h.

Referenced by SiPixelGaussianSmearingRecHitConverterAlgorithm(), and smearHit().

LocalError SiPixelGaussianSmearingRecHitConverterAlgorithm::theError [private]

Definition at line 98 of file SiPixelGaussianSmearingRecHitConverterAlgorithm.h.

Referenced by getError(), and smearHit().

double SiPixelGaussianSmearingRecHitConverterAlgorithm::theErrorX [private]

Definition at line 99 of file SiPixelGaussianSmearingRecHitConverterAlgorithm.h.

Referenced by getErrorX(), and smearHit().

double SiPixelGaussianSmearingRecHitConverterAlgorithm::theErrorY [private]

Definition at line 100 of file SiPixelGaussianSmearingRecHitConverterAlgorithm.h.

Referenced by getErrorY(), and smearHit().

double SiPixelGaussianSmearingRecHitConverterAlgorithm::theErrorZ [private]

Definition at line 101 of file SiPixelGaussianSmearingRecHitConverterAlgorithm.h.

Referenced by getErrorZ(), and smearHit().

unsigned int SiPixelGaussianSmearingRecHitConverterAlgorithm::theLayer [private]

Definition at line 92 of file SiPixelGaussianSmearingRecHitConverterAlgorithm.h.

unsigned int SiPixelGaussianSmearingRecHitConverterAlgorithm::thePixelMultiplicityAlpha [private]

Definition at line 102 of file SiPixelGaussianSmearingRecHitConverterAlgorithm.h.

Referenced by getPixelMultiplicityAlpha(), and smearHit().

unsigned int SiPixelGaussianSmearingRecHitConverterAlgorithm::thePixelMultiplicityBeta [private]

Definition at line 103 of file SiPixelGaussianSmearingRecHitConverterAlgorithm.h.

Referenced by getPixelMultiplicityBeta(), and smearHit().

GeomDetType::SubDetector SiPixelGaussianSmearingRecHitConverterAlgorithm::thePixelPart [private]

Definition at line 82 of file SiPixelGaussianSmearingRecHitConverterAlgorithm.h.

Referenced by SiPixelGaussianSmearingRecHitConverterAlgorithm(), and smearHit().

TFile* SiPixelGaussianSmearingRecHitConverterAlgorithm::thePixelResolutionFile [private]

Definition at line 90 of file SiPixelGaussianSmearingRecHitConverterAlgorithm.h.

Referenced by SiPixelGaussianSmearingRecHitConverterAlgorithm().

Local3DPoint SiPixelGaussianSmearingRecHitConverterAlgorithm::thePosition [private]

Definition at line 94 of file SiPixelGaussianSmearingRecHitConverterAlgorithm.h.

Referenced by getPosition(), and smearHit().

double SiPixelGaussianSmearingRecHitConverterAlgorithm::thePositionX [private]

Definition at line 95 of file SiPixelGaussianSmearingRecHitConverterAlgorithm.h.

Referenced by getPositionX(), and smearHit().

double SiPixelGaussianSmearingRecHitConverterAlgorithm::thePositionY [private]

Definition at line 96 of file SiPixelGaussianSmearingRecHitConverterAlgorithm.h.

Referenced by getPositionY(), and smearHit().

double SiPixelGaussianSmearingRecHitConverterAlgorithm::thePositionZ [private]

Definition at line 97 of file SiPixelGaussianSmearingRecHitConverterAlgorithm.h.

Referenced by getPositionZ(), and smearHit().

bool SiPixelGaussianSmearingRecHitConverterAlgorithm::useCMSSWPixelParameterization [private]

Definition at line 68 of file SiPixelGaussianSmearingRecHitConverterAlgorithm.h.

Referenced by SiPixelGaussianSmearingRecHitConverterAlgorithm(), and smearHit().


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