00001
00013
00014 #include "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->specificType().specificTopology());
00148 const RectangularPixelTopology *rectPixelTopology = static_cast<const RectangularPixelTopology*>(theSpecificTopology);
00149
00150 const int nrows = theSpecificTopology->nrows();
00151 const int ncolumns = theSpecificTopology->ncolumns();
00152
00153 const Local3DPoint lp = simHit.localPosition();
00154
00155 const MeasurementPoint mp = rectPixelTopology->measurementPosition( lp );
00156 float mpy = mp.y();
00157 float mpx = mp.x();
00158
00159 float pixelCenterY = 0.5 + (int)mpy;
00160 float pixelCenterX = 0.5 + (int)mpx;
00161 #ifdef FAMOS_DEBUG
00162 cout<<"Struck pixel center at pitch units x: "<<pixelCenterX<<" y: "<<pixelCenterY<<endl;
00163 #endif
00164
00165 const MeasurementPoint mpCenter(pixelCenterX, pixelCenterY);
00166
00167 const Local3DPoint lpCenter = rectPixelTopology->localPosition( mpCenter );
00168 #ifdef FAMOS_DEBUG
00169 cout<<"Struck point at cm x: "<<lp.x()<<" y: "<<lp.y()<<endl;
00170 cout<<"Struck pixel center at cm x: "<<lpCenter.x()<<" y: "<<lpCenter.y()<<endl;
00171 cout<<"The boundX is "<<boundX<<" boundY is "<<boundY<<endl;
00172 #endif
00173
00174
00175 float xtrk = lp.x() - lpCenter.x();
00176 float ytrk = lp.y() - lpCenter.y();
00177
00178 const float ysize={0.015}, xsize={0.01};
00179
00180 float yhit = 20. + 8.*(ytrk/ysize);
00181 float xhit = 20. + 8.*(xtrk/xsize);
00182 int ybin = (int)yhit;
00183 int xbin = (int)xhit;
00184 float yfrac= yhit - (float)ybin;
00185 float xfrac= xhit - (float)xbin;
00186
00187 if( ybin < 0 ) ybin = 0;
00188 if( ybin > 39 ) ybin = 39;
00189 if( xbin < 0 ) xbin = 0;
00190 if( xbin > 39 ) xbin = 39;
00191
00192
00193
00194 float qbin_frac[4];
00195
00196 float ny1_frac, ny2_frac, nx1_frac, nx2_frac;
00197 bool singlex = false, singley = false;
00198 templ.interpolate(tempId, cotalpha, cotbeta);
00199 templ.qbin_dist(tempId, cotalpha, cotbeta, qbin_frac, ny1_frac, ny2_frac, nx1_frac, nx2_frac );
00200 int nqbin;
00201
00202 double xsizeProbability = random->flatShoot();
00203 double ysizeProbability = random->flatShoot();
00204 bool hitbigx = rectPixelTopology->isItBigPixelInX( (int)mpx );
00205 bool hitbigy = rectPixelTopology->isItBigPixelInY( (int)mpy );
00206
00207 if( hitbigx )
00208 if( xsizeProbability < nx2_frac ) singlex = true;
00209 else singlex = false;
00210 else
00211 if( xsizeProbability < nx1_frac ) singlex = true;
00212 else singlex = false;
00213
00214 if( hitbigy )
00215 if( ysizeProbability < ny2_frac ) singley = true;
00216 else singley = false;
00217 else
00218 if( ysizeProbability < ny1_frac ) singley = true;
00219 else singley = false;
00220
00221
00222
00223
00224 double qbinProbability = random->flatShoot();
00225 for(int i = 0; i<4; ++i) {
00226 nqbin = i;
00227 if(qbinProbability < qbin_frac[i]) break;
00228 }
00229
00230
00231
00232 float ytempl[41][BYSIZE] = {{0}}, xtempl[41][BXSIZE] = {{0}} ;
00233 templ.ytemp(0, 40, ytempl);
00234 templ.xtemp(0, 40, xtempl);
00235
00236 std::vector<double> ytemp(BYSIZE);
00237 for( int i=0; i<BYSIZE; ++i) {
00238 ytemp[i]=(1.-yfrac)*ytempl[ybin][i]+yfrac*ytempl[ybin+1][i];
00239 }
00240
00241 std::vector<double> xtemp(BXSIZE);
00242 for(int i=0; i<BXSIZE; ++i) {
00243 xtemp[i]=(1.-xfrac)*xtempl[xbin][i]+xfrac*xtempl[xbin+1][i];
00244 }
00245
00246
00247 const float qThreshold = templ.s50()*2.0;
00248
00249
00250
00251 int offsetX1=0, offsetX2=0, offsetY1=0, offsetY2=0;
00252 int firstY, lastY, firstX, lastX;
00253 for( firstY = 0; firstY < BYSIZE; ++firstY ) {
00254 bool yCluster = ytemp[firstY] > qThreshold ;
00255 if( yCluster )
00256 {
00257 offsetY1 = BHY -firstY;
00258 break;
00259 }
00260 }
00261 for( lastY = firstY; lastY < BYSIZE; ++lastY )
00262 {
00263 bool yCluster = ytemp[lastY] > qThreshold ;
00264 if( !yCluster )
00265 {
00266 lastY = lastY - 1;
00267 offsetY2 = lastY - BHY;
00268 break;
00269 }
00270 }
00271
00272 for( firstX = 0; firstX < BXSIZE; ++firstX ) {
00273 bool xCluster = xtemp[firstX] > qThreshold ;
00274 if( xCluster ) {
00275 offsetX1 = BHX - firstX;
00276 break;
00277 }
00278 }
00279 for( lastX = firstX; lastX < BXSIZE; ++ lastX ) {
00280 bool xCluster = xtemp[lastX] > qThreshold ;
00281 if( !xCluster ) {
00282 lastX = lastX - 1;
00283 offsetX2 = lastX - BHX;
00284 break;
00285 }
00286 }
00287
00288 bool edge, edgex, edgey;
00289
00290 unsigned int clslenx = offsetX1 + offsetX2 + 1;
00291 unsigned int clsleny = offsetY1 + offsetY2 + 1;
00292
00293 theClslenx = clslenx;
00294 theClsleny = clsleny;
00295
00296 int firstPixelInX = (int)mpx - offsetX1 ;
00297 int firstPixelInY = (int)mpy - offsetY1 ;
00298 int lastPixelInX = (int)mpx + offsetX2 ;
00299 int lastPixelInY = (int)mpy + offsetY2 ;
00300 firstPixelInX = (firstPixelInX >= 0) ? firstPixelInX : 0 ;
00301 firstPixelInY = (firstPixelInY >= 0) ? firstPixelInY : 0 ;
00302 lastPixelInX = (lastPixelInX < nrows ) ? lastPixelInX : nrows-1 ;
00303 lastPixelInY = (lastPixelInY < ncolumns ) ? lastPixelInY : ncolumns-1;
00304
00305 edgex = rectPixelTopology->isItEdgePixelInX( firstPixelInX ) || rectPixelTopology->isItEdgePixelInX( lastPixelInX );
00306 edgey = rectPixelTopology->isItEdgePixelInY( firstPixelInY ) || rectPixelTopology->isItEdgePixelInY( lastPixelInY );
00307 edge = edgex || edgey;
00308
00309
00310
00311 bool hasBigPixelInX = rectPixelTopology->containsBigPixelInX( firstPixelInX, lastPixelInX );
00312 bool hasBigPixelInY = rectPixelTopology->containsBigPixelInY( firstPixelInY, lastPixelInY );
00313
00314
00315 float sigmay, sigmax, sy1, sy2, sx1, sx2;
00316 templ.temperrors(tempId, cotalpha, cotbeta, nqbin, sigmay, sigmax, sy1, sy2, sx1, sx2 );
00317
00318
00319 if( edge ) {
00320 if( edgex && !edgey ) {
00321 theErrorX = 23.0*microntocm;
00322 theErrorY = 39.0*microntocm;
00323 }
00324 else if( !edgex && edgey ) {
00325 theErrorX = 24.0*microntocm;
00326 theErrorY = 96.0*microntocm;
00327 }
00328 else
00329 {
00330 theErrorX = 31.0*microntocm;
00331 theErrorY = 90.0*microntocm;
00332 }
00333
00334 }
00335 else {
00336 if( singlex )
00337 if ( hitbigx )
00338 theErrorX = sx2*microntocm;
00339 else
00340 theErrorX = sx1*microntocm;
00341 else theErrorX = sigmax*microntocm;
00342 if( singley )
00343 if( hitbigy )
00344 theErrorY = sy2*microntocm;
00345 else
00346 theErrorY = sy1*microntocm;
00347 else theErrorY = sigmay*microntocm;
00348 }
00349 theErrorZ = 1e-8;
00350 theError = LocalError( theErrorX*theErrorX, 0., theErrorY*theErrorY);
00351
00352
00353
00354 #ifdef FAMOS_DEBUG
00355 std::cout << " Pixel Errors "
00356 << "\talpha(x) = " << theErrorX
00357 << "\tbeta(y) = " << theErrorY
00358 << std::endl;
00359 #endif
00360
00361
00362 int cotalphaHistBin = (int)( ( cotalpha - rescotAlpha_binMin ) / rescotAlpha_binWidth + 1 );
00363 int cotbetaHistBin = (int)( ( cotbeta - rescotBeta_binMin ) / rescotBeta_binWidth + 1 );
00364
00365 if( cotalphaHistBin < 1 ) cotalphaHistBin = 1;
00366 if( cotbetaHistBin < 1 ) cotbetaHistBin = 1;
00367 if( cotalphaHistBin > (int)rescotAlpha_binN ) cotalphaHistBin = (int)rescotAlpha_binN;
00368 if( cotbetaHistBin > (int)rescotBeta_binN ) cotbetaHistBin = (int)rescotBeta_binN;
00369
00370 unsigned int theXHistN;
00371 unsigned int theYHistN;
00372 if( !isForward ) {
00373 if(edge)
00374 {
00375 theXHistN = cotalphaHistBin * 1000 + cotbetaHistBin * 10 + (nqbin+1);
00376 theYHistN = theXHistN;
00377 }
00378 else
00379 {
00380 if(singlex)
00381 {
00382 if(hitbigx) theXHistN = 1 * 100000 + cotalphaHistBin * 100 + cotbetaHistBin ;
00383 else theXHistN = 1 * 100000 + 1 * 1000 + cotalphaHistBin * 100 + cotbetaHistBin ;
00384 }
00385 else
00386 {
00387 if(hasBigPixelInX) theXHistN = 1 * 1000000 + 1 * 100000 + cotalphaHistBin * 1000 + cotbetaHistBin * 10 + (nqbin+1);
00388 else theXHistN = 1 * 1000000 + 1 * 100000 + 1 * 10000 + cotalphaHistBin * 1000 + cotbetaHistBin * 10 + (nqbin+1);
00389 }
00390 if(singley)
00391 {
00392 if(hitbigy) theYHistN = 1 * 100000 + cotalphaHistBin * 100 + cotbetaHistBin ;
00393 else theYHistN = 1 * 100000 + 1 * 1000 + cotalphaHistBin * 100 + cotbetaHistBin ;
00394 }
00395 else
00396 {
00397 if(hasBigPixelInY) theYHistN = 1 * 1000000 + 1 * 100000 + cotalphaHistBin * 1000 + cotbetaHistBin * 10 + (nqbin+1);
00398 else theYHistN = 1 * 1000000 + 1 * 100000 + 1 * 10000 + cotalphaHistBin * 1000 + cotbetaHistBin * 10 + (nqbin+1);
00399 }
00400 }
00401 }
00402 else
00403 {
00404 if(edge)
00405 {
00406 theXHistN = cotalphaHistBin * 1000 + cotbetaHistBin * 10 + (nqbin+1);
00407 theYHistN = theXHistN;
00408 }
00409 else
00410 {
00411 if( singlex )
00412 if( hitbigx )
00413 theXHistN = 100000 + cotalphaHistBin * 100 + cotbetaHistBin;
00414 else
00415 theXHistN = 100000 + 1000 + cotalphaHistBin * 100 + cotbetaHistBin;
00416 else
00417 theXHistN = 100000 + 10000 + cotalphaHistBin * 1000 + cotbetaHistBin * 10 + (nqbin+1);
00418 if( singley )
00419 if( hitbigy )
00420 theYHistN = 100000 + cotalphaHistBin * 100 + cotbetaHistBin;
00421 else
00422 theYHistN = 100000 + 1000 + cotalphaHistBin * 100 + cotbetaHistBin;
00423 else
00424 theYHistN = 100000 + 10000 + cotalphaHistBin * 1000 + cotbetaHistBin * 10 + (nqbin+1);
00425 }
00426 }
00427 unsigned int counter = 0;
00428 do {
00429
00430
00431 thePositionX = theXHistos[theXHistN]->generate();
00432 thePositionY = theYHistos[theYHistN]->generate();
00433 if( isForward ) thePositionY *= sign;
00434 thePositionZ = 0.0;
00435
00436
00437
00438 thePosition =
00439 Local3DPoint(simHit.localPosition().x() + thePositionX ,
00440 simHit.localPosition().y() + thePositionY ,
00441 simHit.localPosition().z() + thePositionZ );
00442 #ifdef FAMOS_DEBUG
00443 std::cout << " Detector bounds: "
00444 << "\t\tx = " << boundX
00445 << "\ty = " << boundY
00446 << std::endl;
00447 std::cout << " Generated local position "
00448 << "\tx = " << thePosition.x()
00449 << "\ty = " << thePosition.y()
00450 << std::endl;
00451 #endif
00452 counter++;
00453 if(counter > 20) {
00454 thePosition = Local3DPoint(simHit.localPosition().x(),
00455 simHit.localPosition().y(),
00456 simHit.localPosition().z());
00457 break;
00458 }
00459 } while(fabs(thePosition.x()) > boundX || fabs(thePosition.y()) > boundY);
00460
00461 }
00462
00463
00464
00465
00466
00467
00468
00469
00470
00471
00472
00473
00474
00475
00476 bool SiPixelGaussianSmearingRecHitConverterAlgorithm::isFlipped(const PixelGeomDetUnit* theDet) const {
00477
00478 float tmp1 = theDet->surface().toGlobal(Local3DPoint(0.,0.,0.)).perp();
00479 float tmp2 = theDet->surface().toGlobal(Local3DPoint(0.,0.,1.)).perp();
00480
00481 if ( tmp2<tmp1 ) return true;
00482 else return false;
00483 }
00484
00485 void SiPixelGaussianSmearingRecHitConverterAlgorithm::initializeBarrel()
00486 {
00487
00488 rescotAlpha_binMin = -0.2;
00489 rescotAlpha_binWidth = 0.08 ;
00490 rescotAlpha_binN = 5;
00491 rescotBeta_binMin = -5.5;
00492 rescotBeta_binWidth = 1.0;
00493 rescotBeta_binN = 11;
00494 resqbin_binMin = 0;
00495 resqbin_binWidth = 1;
00496 resqbin_binN = 4;
00497
00498
00499 for ( unsigned cotalphaHistBin=1; cotalphaHistBin<=rescotAlpha_binN; ++cotalphaHistBin )
00500 for ( unsigned cotbetaHistBin=1; cotbetaHistBin<=rescotBeta_binN; ++cotbetaHistBin ) {
00501 unsigned int singleBigPixelHistN = 1 * 100000
00502 + cotalphaHistBin * 100
00503 + cotbetaHistBin ;
00504 theXHistos[singleBigPixelHistN] = new SimpleHistogramGenerator(
00505 (TH1F*) thePixelResolutionFile1->Get( Form( "DQMData/clustBPIX/hx%u" , singleBigPixelHistN ) ),
00506 random);
00507 theYHistos[singleBigPixelHistN] = new SimpleHistogramGenerator(
00508 (TH1F*) thePixelResolutionFile1->Get( Form( "DQMData/clustBPIX/hy%u" , singleBigPixelHistN ) ),
00509 random);
00510 unsigned int singlePixelHistN = 1 * 100000 + 1 * 1000
00511 + cotalphaHistBin * 100
00512 + cotbetaHistBin ;
00513 theXHistos[singlePixelHistN] = new SimpleHistogramGenerator(
00514 (TH1F*) thePixelResolutionFile1->Get( Form( "DQMData/clustBPIX/hx%u" , singlePixelHistN ) ),
00515 random);
00516 theYHistos[singlePixelHistN] = new SimpleHistogramGenerator(
00517 (TH1F*) thePixelResolutionFile1->Get( Form( "DQMData/clustBPIX/hy%u" , singlePixelHistN ) ),
00518 random);
00519 for( unsigned qbinBin=1; qbinBin<=resqbin_binN; ++qbinBin ) {
00520 unsigned int edgePixelHistN = cotalphaHistBin * 1000
00521 + cotbetaHistBin * 10
00522 + qbinBin;
00523 theXHistos[edgePixelHistN] = new SimpleHistogramGenerator(
00524 (TH1F*) thePixelResolutionFile2->Get( Form( "DQMData/clustBPIX/hx0%u" ,edgePixelHistN ) ), random );
00525 theYHistos[edgePixelHistN] = new SimpleHistogramGenerator(
00526 (TH1F*) thePixelResolutionFile2->Get( Form( "DQMData/clustBPIX/hy0%u" ,edgePixelHistN ) ), random );
00527 unsigned int multiPixelBigHistN = 1 * 1000000 + 1 * 100000
00528 + cotalphaHistBin * 1000
00529 + cotbetaHistBin * 10
00530 + qbinBin;
00531 theXHistos[multiPixelBigHistN] = new SimpleHistogramGenerator(
00532 (TH1F*) thePixelResolutionFile1->Get( Form( "DQMData/clustBPIX/hx%u" ,multiPixelBigHistN ) ), random );
00533 theYHistos[multiPixelBigHistN] = new SimpleHistogramGenerator(
00534 (TH1F*) thePixelResolutionFile1->Get( Form( "DQMData/clustBPIX/hy%u" ,multiPixelBigHistN ) ), random );
00535 unsigned int multiPixelHistN = 1 * 1000000 + 1 * 100000 + 1 * 10000
00536 + cotalphaHistBin * 1000
00537 + cotbetaHistBin * 10
00538 + qbinBin;
00539 theXHistos[multiPixelHistN] = new SimpleHistogramGenerator(
00540 (TH1F*) thePixelResolutionFile1->Get( Form( "DQMData/clustBPIX/hx%u" , multiPixelHistN ) ),
00541 random);
00542 theYHistos[multiPixelHistN] = new SimpleHistogramGenerator(
00543 (TH1F*) thePixelResolutionFile1->Get( Form( "DQMData/clustBPIX/hy%u" , multiPixelHistN ) ),
00544 random);
00545 }
00546 }
00547 }
00548
00549 void SiPixelGaussianSmearingRecHitConverterAlgorithm::initializeForward()
00550 {
00551
00552 rescotAlpha_binMin = 0.1;
00553 rescotAlpha_binWidth = 0.1 ;
00554 rescotAlpha_binN = 4;
00555 rescotBeta_binMin = 0.;
00556 rescotBeta_binWidth = 0.15;
00557 rescotBeta_binN = 4;
00558 resqbin_binMin = 0;
00559 resqbin_binWidth = 1;
00560 resqbin_binN = 4;
00561
00562
00563 for ( unsigned cotalphaHistBin=1; cotalphaHistBin<=rescotAlpha_binN; ++cotalphaHistBin )
00564 for ( unsigned cotbetaHistBin=1; cotbetaHistBin<=rescotBeta_binN; ++cotbetaHistBin )
00565 for( unsigned qbinBin=1; qbinBin<=resqbin_binN; ++qbinBin ) {
00566 unsigned int edgePixelHistN = cotalphaHistBin * 1000 + cotbetaHistBin * 10 + qbinBin;
00567 theXHistos[edgePixelHistN] = new SimpleHistogramGenerator(
00568 (TH1F*) thePixelResolutionFile1->Get( Form( "DQMData/clustFPIX/fhx0%u" ,edgePixelHistN ) ), random );
00569 theYHistos[edgePixelHistN] = new SimpleHistogramGenerator(
00570 (TH1F*) thePixelResolutionFile1->Get( Form( "DQMData/clustFPIX/fhy0%u" ,edgePixelHistN ) ), random );
00571 unsigned int PixelHistN = 100000 + 10000 + cotalphaHistBin * 1000 + cotbetaHistBin * 10 + qbinBin;
00572 theXHistos[PixelHistN] = new SimpleHistogramGenerator(
00573 (TH1F*) thePixelResolutionFile1->Get( Form( "DQMData/clustFPIX/fhx%u" ,PixelHistN ) ), random );
00574 theYHistos[PixelHistN] = new SimpleHistogramGenerator(
00575 (TH1F*) thePixelResolutionFile1->Get( Form( "DQMData/clustFPIX/fhy%u" ,PixelHistN ) ), random );
00576 }
00577
00578 for ( unsigned cotalphaHistBin=1; cotalphaHistBin<=rescotAlpha_binN; ++cotalphaHistBin )
00579 for ( unsigned cotbetaHistBin=1; cotbetaHistBin<=rescotBeta_binN; ++cotbetaHistBin )
00580 {
00581 unsigned int SingleBigPixelHistN = 100000 + cotalphaHistBin * 100 + cotbetaHistBin;
00582 theXHistos[SingleBigPixelHistN] = new SimpleHistogramGenerator(
00583 (TH1F*) thePixelResolutionFile1->Get( Form( "DQMData/clustFPIX/fhx%u" ,SingleBigPixelHistN ) ), random );
00584 theYHistos[SingleBigPixelHistN] = new SimpleHistogramGenerator(
00585 (TH1F*) thePixelResolutionFile1->Get( Form( "DQMData/clustFPIX/fhy%u" ,SingleBigPixelHistN ) ), random );
00586 unsigned int SinglePixelHistN = 100000 + 1000 + cotalphaHistBin * 100 + cotbetaHistBin;
00587 theXHistos[SinglePixelHistN] = new SimpleHistogramGenerator(
00588 (TH1F*) thePixelResolutionFile1->Get( Form( "DQMData/clustFPIX/fhx%u" ,SinglePixelHistN ) ), random );
00589 theYHistos[SinglePixelHistN] = new SimpleHistogramGenerator(
00590 (TH1F*) thePixelResolutionFile1->Get( Form( "DQMData/clustFPIX/fhy%u" ,SinglePixelHistN ) ), random );
00591
00592 }
00593 }