00001
00013
00014 #include "FastSimulation/TrackingRecHitProducer/interface/SiPixelGaussianSmearingRecHitConverterAlgorithm.h"
00015
00016
00017
00018 #include "Geometry/CommonDetUnit/interface/GeomDetUnit.h"
00019
00020 #include "Geometry/TrackerGeometryBuilder/interface/RectangularPixelTopology.h"
00021 #include "DataFormats/GeometryVector/interface/LocalPoint.h"
00022 #include "DataFormats/GeometryCommonDetAlgo/interface/MeasurementPoint.h"
00023
00024
00025 #include "FastSimulation/Utilities/interface/RandomEngine.h"
00026 #include "FastSimulation/Utilities/interface/SimpleHistogramGenerator.h"
00027
00028
00029
00030
00031 #include <TFile.h>
00032 #include <TH1F.h>
00033
00034
00035
00036
00037 const double PI = 3.14159265358979323;
00038 const double microntocm = 0.0001;
00039 using namespace std;
00040
00041 SiPixelGaussianSmearingRecHitConverterAlgorithm::SiPixelGaussianSmearingRecHitConverterAlgorithm(
00042 const edm::ParameterSet& pset,
00043 GeomDetType::SubDetector pixelPart,
00044 const RandomEngine* engine)
00045 :
00046 pset_(pset),
00047 thePixelPart(pixelPart),
00048 random(engine)
00049 {
00050
00051 useCMSSWPixelParameterization = pset.getParameter<bool>("UseCMSSWPixelParametrization");
00052
00053 if( thePixelPart == GeomDetEnumerators::PixelBarrel ) {
00054 isForward = false;
00055 thePixelResolutionFileName1 = pset_.getParameter<string>( "NewPixelBarrelResolutionFile1" );
00056 thePixelResolutionFile1 = new
00057 TFile( edm::FileInPath( thePixelResolutionFileName1 ).fullPath().c_str() ,"READ");
00058 thePixelResolutionFileName2 = pset_.getParameter<string>( "NewPixelBarrelResolutionFile2" );
00059 thePixelResolutionFile2 = new
00060 TFile( edm::FileInPath( thePixelResolutionFileName2 ).fullPath().c_str() ,"READ");
00061 initializeBarrel();
00062 tempId = pset_.getParameter<int> ( "templateIdBarrel" );
00063 if( ! templ.pushfile(tempId) )
00064 throw cms::Exception("SiPixelGaussianSmearingRecHitConverterAlgorithm:")
00065 <<"SiPixel Barrel Template Not Loaded Correctly!"<<endl;
00066 #ifdef FAMOS_DEBUG
00067 cout<<"The Barrel map size is "<<theXHistos.size()<<" and "<<theYHistos.size()<<endl;
00068 #endif
00069 }
00070 else
00071 if ( thePixelPart == GeomDetEnumerators::PixelEndcap ) {
00072 isForward = true;
00073 thePixelResolutionFileName1 = pset_.getParameter<string>( "NewPixelForwardResolutionFile" );
00074 thePixelResolutionFile1 = new
00075 TFile( edm::FileInPath( thePixelResolutionFileName1 ).fullPath().c_str() ,"READ");
00076 initializeForward();
00077 tempId = pset_.getParameter<int> ( "templateIdForward" );
00078 if( ! templ.pushfile(tempId) )
00079 throw cms::Exception("SiPixelGaussianSmearingRecHitConverterAlgorithm:")
00080 <<"SiPixel Forward Template Not Loaded Correctly!"<<endl;
00081 #ifdef FAMOS_DEBUG
00082 cout<<"The Forward map size is "<<theXHistos.size()<<" and "<<theYHistos.size()<<endl;
00083 #endif
00084 }
00085 else
00086 throw cms::Exception("SiPixelGaussianSmearingRecHitConverterAlgorithm :")
00087 <<"Not a pixel detector"<<endl;
00088 }
00089
00090 SiPixelGaussianSmearingRecHitConverterAlgorithm::~SiPixelGaussianSmearingRecHitConverterAlgorithm()
00091 {
00092
00093 std::map<unsigned,const SimpleHistogramGenerator*>::const_iterator it;
00094 for ( it=theXHistos.begin(); it!=theXHistos.end(); ++it )
00095 delete it->second;
00096 for ( it=theYHistos.begin(); it!=theYHistos.end(); ++it )
00097 delete it->second;
00098
00099 theXHistos.clear();
00100 theYHistos.clear();
00101
00102 }
00103
00104 void SiPixelGaussianSmearingRecHitConverterAlgorithm::smearHit(
00105 const PSimHit& simHit,
00106 const PixelGeomDetUnit* detUnit,
00107 const double boundX,
00108 const double boundY)
00109 {
00110
00111 #ifdef FAMOS_DEBUG
00112 std::cout << " Pixel smearing in " << thePixelPart
00113 << std::endl;
00114 #endif
00115
00116
00117
00118 LocalVector localDir = simHit.momentumAtEntry().unit();
00119 float locx = localDir.x();
00120 float locy = localDir.y();
00121 float locz = localDir.z();
00122
00123
00124 float cotalpha = locx/locz;
00125 if ( isFlipped( detUnit ) ) {
00126 #ifdef FAMOS_DEBUG
00127 std::cout << " isFlipped " << std::endl;
00128 #endif
00129 }
00130
00131 float cotbeta = locy/locz;
00132 float sign=1.;
00133 if( isForward ) {
00134 if( cotbeta < 0 ) sign=-1.;
00135 cotbeta = sign*cotbeta;
00136 }
00137
00138
00139
00140 #ifdef FAMOS_DEBUG
00141 std::cout << " Local Direction " << simHit.localDirection()
00142 << " cotalpha(x) = " << cotalpha
00143 << " cotbeta(y) = " << cotbeta
00144 << std::endl;
00145 #endif
00146
00147 const PixelTopology* theSpecificTopology = &(detUnit->specificTopology());
00148 RectangularPixelTopology rectPixelTopology(theSpecificTopology->nrows(),
00149 theSpecificTopology->ncolumns(),
00150 theSpecificTopology->pitch().first,
00151 theSpecificTopology->pitch().second);
00152 const int nrows = theSpecificTopology->nrows();
00153 const int ncolumns = theSpecificTopology->ncolumns();
00154
00155 const Local3DPoint lp = simHit.localPosition();
00156
00157 const MeasurementPoint mp = rectPixelTopology.measurementPosition( lp );
00158 float mpy = mp.y();
00159 float mpx = mp.x();
00160
00161 float pixelCenterY = 0.5 + (int)mpy;
00162 float pixelCenterX = 0.5 + (int)mpx;
00163 #ifdef FAMOS_DEBUG
00164 cout<<"Struck pixel center at pitch units x: "<<pixelCenterX<<" y: "<<pixelCenterY<<endl;
00165 #endif
00166
00167 const MeasurementPoint mpCenter(pixelCenterX, pixelCenterY);
00168
00169 const Local3DPoint lpCenter = rectPixelTopology.localPosition( mpCenter );
00170 #ifdef FAMOS_DEBUG
00171 cout<<"Struck point at cm x: "<<lp.x()<<" y: "<<lp.y()<<endl;
00172 cout<<"Struck pixel center at cm x: "<<lpCenter.x()<<" y: "<<lpCenter.y()<<endl;
00173 cout<<"The boundX is "<<boundX<<" boundY is "<<boundY<<endl;
00174 #endif
00175
00176
00177 float xtrk = lp.x() - lpCenter.x();
00178 float ytrk = lp.y() - lpCenter.y();
00179
00180 const float ysize={0.015}, xsize={0.01};
00181
00182 float yhit = 20. + 8.*(ytrk/ysize);
00183 float xhit = 20. + 8.*(xtrk/xsize);
00184 int ybin = (int)yhit;
00185 int xbin = (int)xhit;
00186 float yfrac= yhit - (float)ybin;
00187 float xfrac= xhit - (float)xbin;
00188
00189 if( ybin < 0 ) ybin = 0;
00190 if( ybin > 39 ) ybin = 39;
00191 if( xbin < 0 ) xbin = 0;
00192 if( xbin > 39 ) xbin = 39;
00193
00194
00195
00196 float qbin_frac[4];
00197
00198 float ny1_frac, ny2_frac, nx1_frac, nx2_frac;
00199 bool singlex = false, singley = false;
00200 templ.interpolate(tempId, cotalpha, cotbeta);
00201 templ.qbin_dist(tempId, cotalpha, cotbeta, qbin_frac, ny1_frac, ny2_frac, nx1_frac, nx2_frac );
00202 int nqbin;
00203
00204 double xsizeProbability = random->flatShoot();
00205 double ysizeProbability = random->flatShoot();
00206 bool hitbigx = rectPixelTopology.isItBigPixelInX( (int)mpx );
00207 bool hitbigy = rectPixelTopology.isItBigPixelInY( (int)mpy );
00208
00209 if( hitbigx )
00210 if( xsizeProbability < nx2_frac ) singlex = true;
00211 else singlex = false;
00212 else
00213 if( xsizeProbability < nx1_frac ) singlex = true;
00214 else singlex = false;
00215
00216 if( hitbigy )
00217 if( ysizeProbability < ny2_frac ) singley = true;
00218 else singley = false;
00219 else
00220 if( ysizeProbability < ny1_frac ) singley = true;
00221 else singley = false;
00222
00223
00224
00225
00226 double qbinProbability = random->flatShoot();
00227 for(int i = 0; i<4; ++i) {
00228 nqbin = i;
00229 if(qbinProbability < qbin_frac[i]) break;
00230 }
00231
00232
00233
00234 float ytempl[41][BYSIZE] = {{0}}, xtempl[41][BXSIZE] = {{0}} ;
00235 templ.ytemp(0, 40, ytempl);
00236 templ.xtemp(0, 40, xtempl);
00237
00238 std::vector<double> ytemp(BYSIZE);
00239 for( int i=0; i<BYSIZE; ++i) {
00240 ytemp[i]=(1.-yfrac)*ytempl[ybin][i]+yfrac*ytempl[ybin+1][i];
00241 }
00242
00243 std::vector<double> xtemp(BXSIZE);
00244 for(int i=0; i<BXSIZE; ++i) {
00245 xtemp[i]=(1.-xfrac)*xtempl[xbin][i]+xfrac*xtempl[xbin+1][i];
00246 }
00247
00248
00249 const float qThreshold = templ.s50()*2.0;
00250
00251
00252
00253 int offsetX1=0, offsetX2=0, offsetY1=0, offsetY2=0;
00254 int firstY, lastY, firstX, lastX;
00255 for( firstY = 0; firstY < BYSIZE; ++firstY ) {
00256 bool yCluster = ytemp[firstY] > qThreshold ;
00257 if( yCluster )
00258 {
00259 offsetY1 = BHY -firstY;
00260 break;
00261 }
00262 }
00263 for( lastY = firstY; lastY < BYSIZE; ++lastY )
00264 {
00265 bool yCluster = ytemp[lastY] > qThreshold ;
00266 if( !yCluster )
00267 {
00268 lastY = lastY - 1;
00269 offsetY2 = lastY - BHY;
00270 break;
00271 }
00272 }
00273
00274 for( firstX = 0; firstX < BXSIZE; ++firstX ) {
00275 bool xCluster = xtemp[firstX] > qThreshold ;
00276 if( xCluster ) {
00277 offsetX1 = BHX - firstX;
00278 break;
00279 }
00280 }
00281 for( lastX = firstX; lastX < BXSIZE; ++ lastX ) {
00282 bool xCluster = xtemp[lastX] > qThreshold ;
00283 if( !xCluster ) {
00284 lastX = lastX - 1;
00285 offsetX2 = lastX - BHX;
00286 break;
00287 }
00288 }
00289
00290 bool edge, edgex, edgey;
00291 bool bigx, bigy;
00292 unsigned int clslenx = offsetX1 + offsetX2 + 1;
00293 unsigned int clsleny = offsetY1 + offsetY2 + 1;
00294
00295 theClslenx = clslenx;
00296 theClsleny = clsleny;
00297
00298 int firstPixelInX = (int)mpx - offsetX1 ;
00299 int firstPixelInY = (int)mpy - offsetY1 ;
00300 int lastPixelInX = (int)mpx + offsetX2 ;
00301 int lastPixelInY = (int)mpy + offsetY2 ;
00302 firstPixelInX = (firstPixelInX >= 0) ? firstPixelInX : 0 ;
00303 firstPixelInY = (firstPixelInY >= 0) ? firstPixelInY : 0 ;
00304 lastPixelInX = (lastPixelInX < nrows ) ? lastPixelInX : nrows-1 ;
00305 lastPixelInY = (lastPixelInY < ncolumns ) ? lastPixelInY : ncolumns-1;
00306
00307 edgex = rectPixelTopology.isItEdgePixelInX( firstPixelInX ) || rectPixelTopology.isItEdgePixelInX( lastPixelInX );
00308 edgey = rectPixelTopology.isItEdgePixelInY( firstPixelInY ) || rectPixelTopology.isItEdgePixelInY( lastPixelInY );
00309 edge = edgex || edgey;
00310
00311 bigx = rectPixelTopology.isItBigPixelInX( firstPixelInX ) || rectPixelTopology.isItBigPixelInX( lastPixelInX );
00312 bigy = rectPixelTopology.isItBigPixelInY( firstPixelInY ) || rectPixelTopology.isItBigPixelInY( lastPixelInY );
00313 bool hasBigPixelInX = rectPixelTopology.containsBigPixelInX( firstPixelInX, lastPixelInX );
00314 bool hasBigPixelInY = rectPixelTopology.containsBigPixelInY( firstPixelInY, lastPixelInY );
00315
00316
00317 float sigmay, sigmax, sy1, sy2, sx1, sx2;
00318 templ.temperrors(tempId, cotalpha, cotbeta, nqbin, sigmay, sigmax, sy1, sy2, sx1, sx2 );
00319
00320
00321 if( edge ) {
00322 if( edgex && !edgey ) {
00323 theErrorX = 23.0*microntocm;
00324 theErrorY = 39.0*microntocm;
00325 }
00326 else if( !edgex && edgey ) {
00327 theErrorX = 24.0*microntocm;
00328 theErrorY = 96.0*microntocm;
00329 }
00330 else
00331 {
00332 theErrorX = 31.0*microntocm;
00333 theErrorY = 90.0*microntocm;
00334 }
00335
00336 }
00337 else {
00338 if( singlex )
00339 if ( hitbigx )
00340 theErrorX = sx2*microntocm;
00341 else
00342 theErrorX = sx1*microntocm;
00343 else theErrorX = sigmax*microntocm;
00344 if( singley )
00345 if( hitbigy )
00346 theErrorY = sy2*microntocm;
00347 else
00348 theErrorY = sy1*microntocm;
00349 else theErrorY = sigmay*microntocm;
00350 }
00351 theErrorZ = 1e-8;
00352 theError = LocalError( theErrorX*theErrorX, 0., theErrorY*theErrorY);
00353
00354
00355
00356 #ifdef FAMOS_DEBUG
00357 std::cout << " Pixel Errors "
00358 << "\talpha(x) = " << theErrorX
00359 << "\tbeta(y) = " << theErrorY
00360 << std::endl;
00361 #endif
00362
00363
00364 int cotalphaHistBin = (int)( ( cotalpha - rescotAlpha_binMin ) / rescotAlpha_binWidth + 1 );
00365 int cotbetaHistBin = (int)( ( cotbeta - rescotBeta_binMin ) / rescotBeta_binWidth + 1 );
00366
00367 if( cotalphaHistBin < 1 ) cotalphaHistBin = 1;
00368 if( cotbetaHistBin < 1 ) cotbetaHistBin = 1;
00369 if( cotalphaHistBin > (int)rescotAlpha_binN ) cotalphaHistBin = (int)rescotAlpha_binN;
00370 if( cotbetaHistBin > (int)rescotBeta_binN ) cotbetaHistBin = (int)rescotBeta_binN;
00371
00372 unsigned int theXHistN;
00373 unsigned int theYHistN;
00374 if( !isForward ) {
00375 if(edge)
00376 {
00377 theXHistN = cotalphaHistBin * 1000 + cotbetaHistBin * 10 + (nqbin+1);
00378 theYHistN = theXHistN;
00379 }
00380 else
00381 {
00382 if(singlex)
00383 {
00384 if(hitbigx) theXHistN = 1 * 100000 + cotalphaHistBin * 100 + cotbetaHistBin ;
00385 else theXHistN = 1 * 100000 + 1 * 1000 + cotalphaHistBin * 100 + cotbetaHistBin ;
00386 }
00387 else
00388 {
00389 if(hasBigPixelInX) theXHistN = 1 * 1000000 + 1 * 100000 + cotalphaHistBin * 1000 + cotbetaHistBin * 10 + (nqbin+1);
00390 else theXHistN = 1 * 1000000 + 1 * 100000 + 1 * 10000 + cotalphaHistBin * 1000 + cotbetaHistBin * 10 + (nqbin+1);
00391 }
00392 if(singley)
00393 {
00394 if(hitbigy) theYHistN = 1 * 100000 + cotalphaHistBin * 100 + cotbetaHistBin ;
00395 else theYHistN = 1 * 100000 + 1 * 1000 + cotalphaHistBin * 100 + cotbetaHistBin ;
00396 }
00397 else
00398 {
00399 if(hasBigPixelInY) theYHistN = 1 * 1000000 + 1 * 100000 + cotalphaHistBin * 1000 + cotbetaHistBin * 10 + (nqbin+1);
00400 else theYHistN = 1 * 1000000 + 1 * 100000 + 1 * 10000 + cotalphaHistBin * 1000 + cotbetaHistBin * 10 + (nqbin+1);
00401 }
00402 }
00403 }
00404 else
00405 {
00406 if(edge)
00407 {
00408 theXHistN = cotalphaHistBin * 1000 + cotbetaHistBin * 10 + (nqbin+1);
00409 theYHistN = theXHistN;
00410 }
00411 else
00412 {
00413 if( singlex )
00414 if( hitbigx )
00415 theXHistN = 100000 + cotalphaHistBin * 100 + cotbetaHistBin;
00416 else
00417 theXHistN = 100000 + 1000 + cotalphaHistBin * 100 + cotbetaHistBin;
00418 else
00419 theXHistN = 100000 + 10000 + cotalphaHistBin * 1000 + cotbetaHistBin * 10 + (nqbin+1);
00420 if( singley )
00421 if( hitbigy )
00422 theYHistN = 100000 + cotalphaHistBin * 100 + cotbetaHistBin;
00423 else
00424 theYHistN = 100000 + 1000 + cotalphaHistBin * 100 + cotbetaHistBin;
00425 else
00426 theYHistN = 100000 + 10000 + cotalphaHistBin * 1000 + cotbetaHistBin * 10 + (nqbin+1);
00427 }
00428 }
00429 unsigned int counter = 0;
00430 do {
00431
00432
00433 thePositionX = theXHistos[theXHistN]->generate();
00434 thePositionY = theYHistos[theYHistN]->generate();
00435 if( isForward ) thePositionY *= sign;
00436 thePositionZ = 0.0;
00437
00438
00439
00440 thePosition =
00441 Local3DPoint(simHit.localPosition().x() + thePositionX ,
00442 simHit.localPosition().y() + thePositionY ,
00443 simHit.localPosition().z() + thePositionZ );
00444 #ifdef FAMOS_DEBUG
00445 std::cout << " Detector bounds: "
00446 << "\t\tx = " << boundX
00447 << "\ty = " << boundY
00448 << std::endl;
00449 std::cout << " Generated local position "
00450 << "\tx = " << thePosition.x()
00451 << "\ty = " << thePosition.y()
00452 << std::endl;
00453 #endif
00454 counter++;
00455 if(counter > 20) {
00456 thePosition = Local3DPoint(simHit.localPosition().x(),
00457 simHit.localPosition().y(),
00458 simHit.localPosition().z());
00459 break;
00460 }
00461 } while(fabs(thePosition.x()) > boundX || fabs(thePosition.y()) > boundY);
00462
00463 }
00464
00465
00466
00467
00468
00469
00470
00471
00472
00473
00474
00475
00476
00477
00478 bool SiPixelGaussianSmearingRecHitConverterAlgorithm::isFlipped(const PixelGeomDetUnit* theDet) const {
00479
00480 float tmp1 = theDet->surface().toGlobal(Local3DPoint(0.,0.,0.)).perp();
00481 float tmp2 = theDet->surface().toGlobal(Local3DPoint(0.,0.,1.)).perp();
00482
00483 if ( tmp2<tmp1 ) return true;
00484 else return false;
00485 }
00486
00487 void SiPixelGaussianSmearingRecHitConverterAlgorithm::initializeBarrel()
00488 {
00489
00490 rescotAlpha_binMin = -0.2;
00491 rescotAlpha_binWidth = 0.08 ;
00492 rescotAlpha_binN = 5;
00493 rescotBeta_binMin = -5.5;
00494 rescotBeta_binWidth = 1.0;
00495 rescotBeta_binN = 11;
00496 resqbin_binMin = 0;
00497 resqbin_binWidth = 1;
00498 resqbin_binN = 4;
00499
00500
00501 for ( unsigned cotalphaHistBin=1; cotalphaHistBin<=rescotAlpha_binN; ++cotalphaHistBin )
00502 for ( unsigned cotbetaHistBin=1; cotbetaHistBin<=rescotBeta_binN; ++cotbetaHistBin ) {
00503 unsigned int singleBigPixelHistN = 1 * 100000
00504 + cotalphaHistBin * 100
00505 + cotbetaHistBin ;
00506 theXHistos[singleBigPixelHistN] = new SimpleHistogramGenerator(
00507 (TH1F*) thePixelResolutionFile1->Get( Form( "DQMData/clustBPIX/hx%u" , singleBigPixelHistN ) ),
00508 random);
00509 theYHistos[singleBigPixelHistN] = new SimpleHistogramGenerator(
00510 (TH1F*) thePixelResolutionFile1->Get( Form( "DQMData/clustBPIX/hy%u" , singleBigPixelHistN ) ),
00511 random);
00512 unsigned int singlePixelHistN = 1 * 100000 + 1 * 1000
00513 + cotalphaHistBin * 100
00514 + cotbetaHistBin ;
00515 theXHistos[singlePixelHistN] = new SimpleHistogramGenerator(
00516 (TH1F*) thePixelResolutionFile1->Get( Form( "DQMData/clustBPIX/hx%u" , singlePixelHistN ) ),
00517 random);
00518 theYHistos[singlePixelHistN] = new SimpleHistogramGenerator(
00519 (TH1F*) thePixelResolutionFile1->Get( Form( "DQMData/clustBPIX/hy%u" , singlePixelHistN ) ),
00520 random);
00521 for( unsigned qbinBin=1; qbinBin<=resqbin_binN; ++qbinBin ) {
00522 unsigned int edgePixelHistN = cotalphaHistBin * 1000
00523 + cotbetaHistBin * 10
00524 + qbinBin;
00525 theXHistos[edgePixelHistN] = new SimpleHistogramGenerator(
00526 (TH1F*) thePixelResolutionFile2->Get( Form( "DQMData/clustBPIX/hx0%u" ,edgePixelHistN ) ), random );
00527 theYHistos[edgePixelHistN] = new SimpleHistogramGenerator(
00528 (TH1F*) thePixelResolutionFile2->Get( Form( "DQMData/clustBPIX/hy0%u" ,edgePixelHistN ) ), random );
00529 unsigned int multiPixelBigHistN = 1 * 1000000 + 1 * 100000
00530 + cotalphaHistBin * 1000
00531 + cotbetaHistBin * 10
00532 + qbinBin;
00533 theXHistos[multiPixelBigHistN] = new SimpleHistogramGenerator(
00534 (TH1F*) thePixelResolutionFile1->Get( Form( "DQMData/clustBPIX/hx%u" ,multiPixelBigHistN ) ), random );
00535 theYHistos[multiPixelBigHistN] = new SimpleHistogramGenerator(
00536 (TH1F*) thePixelResolutionFile1->Get( Form( "DQMData/clustBPIX/hy%u" ,multiPixelBigHistN ) ), random );
00537 unsigned int multiPixelHistN = 1 * 1000000 + 1 * 100000 + 1 * 10000
00538 + cotalphaHistBin * 1000
00539 + cotbetaHistBin * 10
00540 + qbinBin;
00541 theXHistos[multiPixelHistN] = new SimpleHistogramGenerator(
00542 (TH1F*) thePixelResolutionFile1->Get( Form( "DQMData/clustBPIX/hx%u" , multiPixelHistN ) ),
00543 random);
00544 theYHistos[multiPixelHistN] = new SimpleHistogramGenerator(
00545 (TH1F*) thePixelResolutionFile1->Get( Form( "DQMData/clustBPIX/hy%u" , multiPixelHistN ) ),
00546 random);
00547 }
00548 }
00549 }
00550
00551 void SiPixelGaussianSmearingRecHitConverterAlgorithm::initializeForward()
00552 {
00553
00554 rescotAlpha_binMin = 0.1;
00555 rescotAlpha_binWidth = 0.1 ;
00556 rescotAlpha_binN = 4;
00557 rescotBeta_binMin = 0.;
00558 rescotBeta_binWidth = 0.15;
00559 rescotBeta_binN = 4;
00560 resqbin_binMin = 0;
00561 resqbin_binWidth = 1;
00562 resqbin_binN = 4;
00563
00564
00565 for ( unsigned cotalphaHistBin=1; cotalphaHistBin<=rescotAlpha_binN; ++cotalphaHistBin )
00566 for ( unsigned cotbetaHistBin=1; cotbetaHistBin<=rescotBeta_binN; ++cotbetaHistBin )
00567 for( unsigned qbinBin=1; qbinBin<=resqbin_binN; ++qbinBin ) {
00568 unsigned int edgePixelHistN = cotalphaHistBin * 1000 + cotbetaHistBin * 10 + qbinBin;
00569 theXHistos[edgePixelHistN] = new SimpleHistogramGenerator(
00570 (TH1F*) thePixelResolutionFile1->Get( Form( "DQMData/clustFPIX/fhx0%u" ,edgePixelHistN ) ), random );
00571 theYHistos[edgePixelHistN] = new SimpleHistogramGenerator(
00572 (TH1F*) thePixelResolutionFile1->Get( Form( "DQMData/clustFPIX/fhy0%u" ,edgePixelHistN ) ), random );
00573 unsigned int PixelHistN = 100000 + 10000 + cotalphaHistBin * 1000 + cotbetaHistBin * 10 + qbinBin;
00574 theXHistos[PixelHistN] = new SimpleHistogramGenerator(
00575 (TH1F*) thePixelResolutionFile1->Get( Form( "DQMData/clustFPIX/fhx%u" ,PixelHistN ) ), random );
00576 theYHistos[PixelHistN] = new SimpleHistogramGenerator(
00577 (TH1F*) thePixelResolutionFile1->Get( Form( "DQMData/clustFPIX/fhy%u" ,PixelHistN ) ), random );
00578 }
00579
00580 for ( unsigned cotalphaHistBin=1; cotalphaHistBin<=rescotAlpha_binN; ++cotalphaHistBin )
00581 for ( unsigned cotbetaHistBin=1; cotbetaHistBin<=rescotBeta_binN; ++cotbetaHistBin )
00582 {
00583 unsigned int SingleBigPixelHistN = 100000 + cotalphaHistBin * 100 + cotbetaHistBin;
00584 theXHistos[SingleBigPixelHistN] = new SimpleHistogramGenerator(
00585 (TH1F*) thePixelResolutionFile1->Get( Form( "DQMData/clustFPIX/fhx%u" ,SingleBigPixelHistN ) ), random );
00586 theYHistos[SingleBigPixelHistN] = new SimpleHistogramGenerator(
00587 (TH1F*) thePixelResolutionFile1->Get( Form( "DQMData/clustFPIX/fhy%u" ,SingleBigPixelHistN ) ), random );
00588 unsigned int SinglePixelHistN = 100000 + 1000 + cotalphaHistBin * 100 + cotbetaHistBin;
00589 theXHistos[SinglePixelHistN] = new SimpleHistogramGenerator(
00590 (TH1F*) thePixelResolutionFile1->Get( Form( "DQMData/clustFPIX/fhx%u" ,SinglePixelHistN ) ), random );
00591 theYHistos[SinglePixelHistN] = new SimpleHistogramGenerator(
00592 (TH1F*) thePixelResolutionFile1->Get( Form( "DQMData/clustFPIX/fhy%u" ,SinglePixelHistN ) ), random );
00593
00594 }
00595 }