00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031
00032
00033
00034
00035
00036
00037
00038
00039
00040
00041
00042 #include <iostream>
00043
00044 #include "SimGeneral/NoiseGenerators/interface/GaussianTailNoiseGenerator.h"
00045
00046 #include "DataFormats/SiPixelDigi/interface/PixelDigi.h"
00047 #include "SimDataFormats/TrackerDigiSimLink/interface/PixelDigiSimLink.h"
00048 #include "SimDataFormats/TrackingHit/interface/PSimHitContainer.h"
00049 #include "SimTracker/Common/interface/SiG4UniversalFluctuation.h"
00050 #include "SiPixelDigitizerAlgorithm.h"
00051
00052 #include <gsl/gsl_sf_erf.h>
00053 #include "FWCore/Utilities/interface/RandomNumberGenerator.h"
00054 #include "CLHEP/Random/RandGaussQ.h"
00055 #include "CLHEP/Random/RandFlat.h"
00056
00057
00058 #include "DataFormats/SiPixelDetId/interface/PixelSubdetector.h"
00059 #include "DataFormats/TrackerCommon/interface/TrackerTopology.h"
00060 #include "Geometry/Records/interface/IdealGeometryRecord.h"
00061
00062 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00063 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00064 #include "FWCore/ServiceRegistry/interface/Service.h"
00065 #include "FWCore/Utilities/interface/Exception.h"
00066 #include "CalibTracker/SiPixelESProducers/interface/SiPixelGainCalibrationOfflineSimService.h"
00067
00068
00069 #include "DataFormats/DetId/interface/DetId.h"
00070
00071 #include "CondFormats/SiPixelObjects/interface/GlobalPixel.h"
00072
00073 #include "CondFormats/DataRecord/interface/SiPixelQualityRcd.h"
00074 #include "CondFormats/DataRecord/interface/SiPixelFedCablingMapRcd.h"
00075 #include "CondFormats/DataRecord/interface/SiPixelLorentzAngleSimRcd.h"
00076 #include "CondFormats/SiPixelObjects/interface/SiPixelFedCablingMap.h"
00077 #include "CondFormats/SiPixelObjects/interface/SiPixelFedCablingTree.h"
00078 #include "CondFormats/SiPixelObjects/interface/SiPixelFedCabling.h"
00079 #include "CondFormats/SiPixelObjects/interface/PixelIndices.h"
00080 #include "CondFormats/SiPixelObjects/interface/SiPixelLorentzAngle.h"
00081 #include "CondFormats/SiPixelObjects/interface/SiPixelQuality.h"
00082 #include "CondFormats/SiPixelObjects/interface/PixelROC.h"
00083 #include "CondFormats/SiPixelObjects/interface/LocalPixel.h"
00084 #include "CondFormats/SiPixelObjects/interface/CablingPathToDetUnit.h"
00085
00086 #include "CondFormats/SiPixelObjects/interface/SiPixelFrameReverter.h"
00087 #include "CondFormats/SiPixelObjects/interface/PixelFEDCabling.h"
00088 #include "CondFormats/SiPixelObjects/interface/PixelFEDLink.h"
00089 #include "DataFormats/FEDRawData/interface/FEDNumbering.h"
00090 #include "DataFormats/SiPixelDetId/interface/PixelBarrelName.h"
00091 #include "DataFormats/SiPixelDetId/interface/PixelEndcapName.h"
00092
00093
00094 #include "Geometry/TrackerGeometryBuilder/interface/TrackerGeometry.h"
00095 #include "Geometry/Records/interface/TrackerDigiGeometryRecord.h"
00096 #include "Geometry/TrackerGeometryBuilder/interface/PixelGeomDetUnit.h"
00097 #include "Geometry/CommonTopologies/interface/PixelTopology.h"
00098
00099 #include "CondFormats/SiPixelObjects/interface/PixelROC.h"
00100
00101 using namespace edm;
00102 using namespace sipixelobjects;
00103
00104 #define TP_DEBUG // protect all LogDebug with ifdef. Takes too much CPU
00105
00106
00107 void SiPixelDigitizerAlgorithm::init(const edm::EventSetup& es) {
00108 if(use_ineff_from_db_){
00109 theSiPixelGainCalibrationService_->setESObjects( es );
00110 }
00111 if(use_deadmodule_DB_) {
00112 es.get<SiPixelQualityRcd>().get(SiPixelBadModule_);
00113 }
00114 if(use_LorentzAngle_DB_) {
00115
00116 es.get<SiPixelLorentzAngleSimRcd>().get(SiPixelLorentzAngle_);
00117 }
00118
00119 es.get<SiPixelFedCablingMapRcd>().get(map_);
00120 es.get<TrackerDigiGeometryRecord>().get(geom_);
00121 }
00122
00123
00124
00125 SiPixelDigitizerAlgorithm::SiPixelDigitizerAlgorithm(const edm::ParameterSet& conf, CLHEP::HepRandomEngine& eng) :
00126
00127 _signal(),
00128 makeDigiSimLinks_(conf.getUntrackedParameter<bool>("makeDigiSimLinks", true)),
00129 use_ineff_from_db_(conf.getParameter<bool>("useDB")),
00130 use_module_killing_(conf.getParameter<bool>("killModules")),
00131 use_deadmodule_DB_(conf.getParameter<bool>("DeadModules_DB")),
00132 use_LorentzAngle_DB_(conf.getParameter<bool>("LorentzAngle_DB")),
00133
00134 DeadModules(use_deadmodule_DB_ ? Parameters() : conf.getParameter<Parameters>("DeadModules")),
00135
00136
00137
00138 GeVperElectron(3.61E-09),
00139 Sigma0(0.00037),
00140 Dist300(0.0300),
00141 alpha2Order(conf.getParameter<bool>("Alpha2Order")),
00142 ClusterWidth(3.),
00143
00144
00145
00146
00147 NumberOfBarrelLayers(conf.exists("NumPixelBarrel")?conf.getParameter<int>("NumPixelBarrel"):3),
00148 NumberOfEndcapDisks(conf.exists("NumPixelEndcap")?conf.getParameter<int>("NumPixelEndcap"):2),
00149
00150
00151
00152
00153
00154
00155 theElectronPerADC(conf.getParameter<double>("ElectronPerAdc")),
00156
00157
00158
00159 theAdcFullScale(conf.getParameter<int>("AdcFullScale")),
00160 theAdcFullScaleStack(conf.exists("AdcFullScaleStack")?conf.getParameter<int>("AdcFullScaleStack"):255),
00161 theFirstStackLayer(conf.exists("FirstStackLayer")?conf.getParameter<int>("FirstStackLayer"):5),
00162
00163
00164
00165 theNoiseInElectrons(conf.getParameter<double>("NoiseInElectrons")),
00166
00167
00168
00169 theReadoutNoise(conf.getParameter<double>("ReadoutNoiseInElec")),
00170
00171
00172
00173
00174 theThresholdInE_FPix(conf.getParameter<double>("ThresholdInElectrons_FPix")),
00175 theThresholdInE_BPix(conf.getParameter<double>("ThresholdInElectrons_BPix")),
00176 theThresholdInE_BPix_L1(conf.exists("ThresholdInElectrons_BPix_L1")?conf.getParameter<double>("ThresholdInElectrons_BPix_L1"):theThresholdInE_BPix),
00177
00178
00179 theThresholdSmearing_FPix(conf.getParameter<double>("ThresholdSmearing_FPix")),
00180 theThresholdSmearing_BPix(conf.getParameter<double>("ThresholdSmearing_BPix")),
00181 theThresholdSmearing_BPix_L1(conf.exists("ThresholdSmearing_BPix_L1")?conf.getParameter<double>("ThresholdSmearing_BPix_L1"):theThresholdSmearing_BPix),
00182
00183
00184 electronsPerVCAL(conf.getParameter<double>("ElectronsPerVcal")),
00185 electronsPerVCAL_Offset(conf.getParameter<double>("ElectronsPerVcal_Offset")),
00186
00187
00188
00189 theTofLowerCut(conf.getParameter<double>("TofLowerCut")),
00190 theTofUpperCut(conf.getParameter<double>("TofUpperCut")),
00191
00192
00193 tanLorentzAnglePerTesla_FPix(use_LorentzAngle_DB_ ? 0.0 : conf.getParameter<double>("TanLorentzAnglePerTesla_FPix")),
00194 tanLorentzAnglePerTesla_BPix(use_LorentzAngle_DB_ ? 0.0 : conf.getParameter<double>("TanLorentzAnglePerTesla_BPix")),
00195
00196
00197 FPix_p0(conf.getParameter<double>("FPix_SignalResponse_p0")),
00198 FPix_p1(conf.getParameter<double>("FPix_SignalResponse_p1")),
00199 FPix_p2(conf.getParameter<double>("FPix_SignalResponse_p2")),
00200 FPix_p3(conf.getParameter<double>("FPix_SignalResponse_p3")),
00201
00202 BPix_p0(conf.getParameter<double>("BPix_SignalResponse_p0")),
00203 BPix_p1(conf.getParameter<double>("BPix_SignalResponse_p1")),
00204 BPix_p2(conf.getParameter<double>("BPix_SignalResponse_p2")),
00205 BPix_p3(conf.getParameter<double>("BPix_SignalResponse_p3")),
00206
00207
00208 addNoise(conf.getParameter<bool>("AddNoise")),
00209
00210
00211
00212 addChargeVCALSmearing(conf.getParameter<bool>("ChargeVCALSmearing")),
00213
00214
00215 addNoisyPixels(conf.getParameter<bool>("AddNoisyPixels")),
00216
00217
00218 fluctuateCharge(conf.getUntrackedParameter<bool>("FluctuateCharge",true)),
00219
00220
00221 AddPixelInefficiency(conf.getParameter<bool>("AddPixelInefficiencyFromPython")),
00222
00223
00224 addThresholdSmearing(conf.getParameter<bool>("AddThresholdSmearing")),
00225
00226
00227 doMissCalibrate(conf.getParameter<bool>("MissCalibrate")),
00228 theGainSmearing(conf.getParameter<double>("GainSmearing")),
00229 theOffsetSmearing(conf.getParameter<double>("OffsetSmearing")),
00230
00231
00232 pseudoRadDamage(conf.exists("PseudoRadDamage")?conf.getParameter<double>("PseudoRadDamage"):double(0.0)),
00233 pseudoRadDamageRadius(conf.exists("PseudoRadDamageRadius")?conf.getParameter<double>("PseudoRadDamageRadius"):double(0.0)),
00234
00235
00236
00237
00238 tMax(conf.getParameter<double>("DeltaProductionCut")),
00239
00240 fluctuate(fluctuateCharge ? new SiG4UniversalFluctuation(eng) : 0),
00241 theNoiser(addNoise ? new GaussianTailNoiseGenerator(eng) : 0),
00242 calmap(doMissCalibrate ? initCal() : std::map<int,CalParameters,std::less<int> >()),
00243 theSiPixelGainCalibrationService_(use_ineff_from_db_ ? new SiPixelGainCalibrationOfflineSimService(conf) : 0),
00244 pixelEfficiencies_(conf, AddPixelInefficiency,NumberOfBarrelLayers,NumberOfEndcapDisks),
00245 flatDistribution_((addNoise || AddPixelInefficiency || fluctuateCharge || addThresholdSmearing) ? new CLHEP::RandFlat(eng, 0., 1.) : 0),
00246 gaussDistribution_((addNoise || AddPixelInefficiency || fluctuateCharge || addThresholdSmearing) ? new CLHEP::RandGaussQ(eng, 0., theReadoutNoise) : 0),
00247 gaussDistributionVCALNoise_((addNoise || AddPixelInefficiency || fluctuateCharge || addThresholdSmearing) ? new CLHEP::RandGaussQ(eng, 0., 1.) : 0),
00248
00249 smearedThreshold_FPix_(addThresholdSmearing ? new CLHEP::RandGaussQ(eng, theThresholdInE_FPix , theThresholdSmearing_FPix) : 0),
00250 smearedThreshold_BPix_(addThresholdSmearing ? new CLHEP::RandGaussQ(eng, theThresholdInE_BPix , theThresholdSmearing_BPix) : 0),
00251 smearedThreshold_BPix_L1_(addThresholdSmearing ? new CLHEP::RandGaussQ(eng, theThresholdInE_BPix_L1 , theThresholdSmearing_BPix_L1) : 0)
00252 {
00253 LogInfo ("PixelDigitizer ") <<"SiPixelDigitizerAlgorithm constructed"
00254 <<"Configuration parameters:"
00255 << "Threshold/Gain = "
00256 << "threshold in electron FPix = "
00257 << theThresholdInE_FPix
00258 << "threshold in electron BPix = "
00259 << theThresholdInE_BPix
00260 << "threshold in electron BPix Layer1 = "
00261 << theThresholdInE_BPix_L1
00262 <<" " << theElectronPerADC << " " << theAdcFullScale
00263 << " The delta cut-off is set to " << tMax
00264 << " pix-inefficiency "<<AddPixelInefficiency;
00265
00266 }
00267
00268 std::map<int, SiPixelDigitizerAlgorithm::CalParameters, std::less<int> >
00269 SiPixelDigitizerAlgorithm::initCal() const {
00270
00271 using std::cerr;
00272 using std::cout;
00273 using std::endl;
00274
00275 std::map<int, SiPixelDigitizerAlgorithm::CalParameters, std::less<int> > calmap;
00276
00277 LogDebug ("PixelDigitizer ")
00278 << " miss-calibrate the pixel amplitude ";
00279
00280 const bool ReadCalParameters = false;
00281 if(ReadCalParameters) {
00282
00283 std::ifstream in_file;
00284 char filename[80] = "phCalibrationFit_C0.dat";
00285
00286 in_file.open(filename, std::ios::in );
00287 if(in_file.bad()) {
00288 cout << " File not found " << endl;
00289 return calmap;
00290 }
00291 cout << " file opened : " << filename << endl;
00292
00293 char line[500];
00294 for (int i = 0; i < 3; i++) {
00295 in_file.getline(line, 500,'\n');
00296 cout<<line<<endl;
00297 }
00298
00299 cout << " test map" << endl;
00300
00301 float par0,par1,par2,par3;
00302 int colid,rowid;
00303 std::string name;
00304
00305 for(int i=0;i<(52*80);i++) {
00306 in_file >> par0 >> par1 >> par2 >> par3 >> name >> colid >> rowid;
00307 if(in_file.bad()) {
00308 cerr << "Cannot read data file" << endl;
00309 return calmap;
00310 }
00311 if( in_file.eof() != 0 ) {
00312 cerr << in_file.eof() << " " << in_file.gcount() << " "
00313 << in_file.fail() << " " << in_file.good() << " end of file "
00314 << endl;
00315 return calmap;
00316 }
00317
00318
00319
00320
00321 CalParameters onePix;
00322 onePix.p0=par0;
00323 onePix.p1=par1;
00324 onePix.p2=par2;
00325 onePix.p3=par3;
00326
00327
00328 int chan = PixelIndices::pixelToChannelROC(rowid,colid);
00329 calmap.insert(std::pair<int,CalParameters>(chan,onePix));
00330
00331
00332 std::pair<int,int> p = PixelIndices::channelToPixelROC(chan);
00333 if(rowid!=p.first) cout<<" wrong channel row "<<rowid<<" "<<p.first<<endl;
00334 if(colid!=p.second) cout<<" wrong channel col "<<colid<<" "<<p.second<<endl;
00335
00336 }
00337
00338 cout << " map size " << calmap.size() <<" max "<<calmap.max_size() << " "
00339 <<calmap.empty()<< endl;
00340
00341
00342
00343
00344
00345
00346
00347
00348
00349
00350
00351
00352
00353
00354
00355
00356 }
00357 return calmap;
00358 }
00359
00360
00361 SiPixelDigitizerAlgorithm::~SiPixelDigitizerAlgorithm() {
00362 LogDebug ("PixelDigitizer")<<"SiPixelDigitizerAlgorithm deleted";
00363 }
00364
00365 SiPixelDigitizerAlgorithm::PixelEfficiencies::PixelEfficiencies(const edm::ParameterSet& conf, bool AddPixelInefficiency, int NumberOfBarrelLayers, int NumberOfEndcapDisks) {
00366
00367
00368 int NumberOfTotLayers = NumberOfBarrelLayers + NumberOfEndcapDisks;
00369 FPixIndex=NumberOfBarrelLayers;
00370 if (AddPixelInefficiency){
00371 int i=0;
00372 thePixelColEfficiency[i++] = conf.getParameter<double>("thePixelColEfficiency_BPix1");
00373 thePixelColEfficiency[i++] = conf.getParameter<double>("thePixelColEfficiency_BPix2");
00374 thePixelColEfficiency[i++] = conf.getParameter<double>("thePixelColEfficiency_BPix3");
00375 if (NumberOfBarrelLayers>=4){thePixelColEfficiency[i++] = conf.getParameter<double>("thePixelColEfficiency_BPix4");}
00376
00377 i=0;
00378 thePixelEfficiency[i++] = conf.getParameter<double>("thePixelEfficiency_BPix1");
00379 thePixelEfficiency[i++] = conf.getParameter<double>("thePixelEfficiency_BPix2");
00380 thePixelEfficiency[i++] = conf.getParameter<double>("thePixelEfficiency_BPix3");
00381 if (NumberOfBarrelLayers>=4){thePixelEfficiency[i++] = conf.getParameter<double>("thePixelEfficiency_BPix4");}
00382
00383 i=0;
00384 thePixelChipEfficiency[i++] = conf.getParameter<double>("thePixelChipEfficiency_BPix1");
00385 thePixelChipEfficiency[i++] = conf.getParameter<double>("thePixelChipEfficiency_BPix2");
00386 thePixelChipEfficiency[i++] = conf.getParameter<double>("thePixelChipEfficiency_BPix3");
00387 if (NumberOfBarrelLayers>=4){thePixelChipEfficiency[i++] = conf.getParameter<double>("thePixelChipEfficiency_BPix4");}
00388
00389 if (NumberOfBarrelLayers>=5){
00390 if (NumberOfTotLayers>20){throw cms::Exception("Configuration") <<"SiPixelDigitizer was given more layers than it can handle";}
00391
00392 for (int j=5 ; j<=NumberOfBarrelLayers ; j++){
00393 thePixelColEfficiency[j-1]=0.999;
00394 thePixelEfficiency[j-1]=0.999;
00395 thePixelChipEfficiency[j-1]=0.999;
00396 }
00397 }
00398
00399 i=FPixIndex;
00400 thePixelColEfficiency[i++] = conf.getParameter<double>("thePixelColEfficiency_FPix1");
00401 thePixelColEfficiency[i++] = conf.getParameter<double>("thePixelColEfficiency_FPix2");
00402 if (NumberOfEndcapDisks>=3){thePixelColEfficiency[i++] = conf.getParameter<double>("thePixelColEfficiency_FPix3");}
00403 i=FPixIndex;
00404 thePixelEfficiency[i++] = conf.getParameter<double>("thePixelEfficiency_FPix1");
00405 thePixelEfficiency[i++] = conf.getParameter<double>("thePixelEfficiency_FPix2");
00406 if (NumberOfEndcapDisks>=3){thePixelEfficiency[i++] = conf.getParameter<double>("thePixelEfficiency_FPix3");}
00407 i=FPixIndex;
00408 thePixelChipEfficiency[i++] = conf.getParameter<double>("thePixelChipEfficiency_FPix1");
00409 thePixelChipEfficiency[i++] = conf.getParameter<double>("thePixelChipEfficiency_FPix2");
00410 if (NumberOfEndcapDisks>=3){thePixelChipEfficiency[i++] = conf.getParameter<double>("thePixelChipEfficiency_FPix3");}
00411
00412 if (NumberOfEndcapDisks>=4){
00413 if (NumberOfTotLayers>20){throw cms::Exception("Configuration") <<"SiPixelDigitizer was given more layers than it can handle";}
00414
00415 for (int j=4+FPixIndex ; j<=NumberOfEndcapDisks+NumberOfBarrelLayers ; j++){
00416 thePixelColEfficiency[j-1]=0.999;
00417 thePixelEfficiency[j-1]=0.999;
00418 thePixelChipEfficiency[j-1]=0.999;
00419 }
00420 }
00421 }
00422
00423
00424 if(!AddPixelInefficiency) {
00425 for (int i=0; i<NumberOfTotLayers;i++) {
00426 thePixelEfficiency[i] = 1.;
00427 thePixelColEfficiency[i] = 1.;
00428 thePixelChipEfficiency[i] = 1.;
00429 }
00430 }
00431 }
00432
00433
00434 void SiPixelDigitizerAlgorithm::accumulateSimHits(std::vector<PSimHit>::const_iterator inputBegin,
00435 std::vector<PSimHit>::const_iterator inputEnd,
00436 const PixelGeomDetUnit* pixdet,
00437 const GlobalVector& bfield) {
00438
00439
00440
00441 uint32_t detId = pixdet->geographicalId().rawId();
00442 for (std::vector<PSimHit>::const_iterator ssbegin = inputBegin; ssbegin != inputEnd; ++ssbegin) {
00443
00444 if((*ssbegin).detUnitId() != detId) {
00445 continue;
00446 }
00447
00448 #ifdef TP_DEBUG
00449 LogDebug ("Pixel Digitizer")
00450 << (*ssbegin).particleType() << " " << (*ssbegin).pabs() << " "
00451 << (*ssbegin).energyLoss() << " " << (*ssbegin).tof() << " "
00452 << (*ssbegin).trackId() << " " << (*ssbegin).processType() << " "
00453 << (*ssbegin).detUnitId()
00454 << (*ssbegin).entryPoint() << " " << (*ssbegin).exitPoint() ;
00455 #endif
00456
00457
00458 std::vector<EnergyDepositUnit> ionization_points;
00459 std::vector<SignalPoint> collection_points;
00460
00461
00462
00463 if ( ((*ssbegin).tof() - pixdet->surface().toGlobal((*ssbegin).localPosition()).mag()/30.)>= theTofLowerCut &&
00464 ((*ssbegin).tof()- pixdet->surface().toGlobal((*ssbegin).localPosition()).mag()/30.) <= theTofUpperCut ) {
00465 primary_ionization(*ssbegin, ionization_points);
00466 drift(*ssbegin, pixdet, bfield, ionization_points, collection_points);
00467
00468 induce_signal(*ssbegin, pixdet, collection_points);
00469 }
00470 }
00471
00472 }
00473
00474
00475 void SiPixelDigitizerAlgorithm::digitize(const PixelGeomDetUnit* pixdet,
00476 std::vector<PixelDigi>& digis,
00477 std::vector<PixelDigiSimLink>& simlinks, const TrackerTopology *tTopo) {
00478
00479
00480
00481
00482
00483 uint32_t detID = pixdet->geographicalId().rawId();
00484 const signal_map_type& theSignal = _signal[detID];
00485
00486 const PixelTopology* topol=&pixdet->specificTopology();
00487 int numColumns = topol->ncolumns();
00488 int numRows = topol->nrows();
00489
00490
00491
00492
00493
00494 unsigned int Sub_detid=DetId(detID).subdetId();
00495
00496 float thePixelThresholdInE = 0.;
00497
00498 if(theNoiseInElectrons>0.){
00499 if(Sub_detid == PixelSubdetector::PixelBarrel){
00500 int lay = tTopo->pxbLayer(detID);
00501 if(addThresholdSmearing) {
00502 if(lay==1) {
00503 thePixelThresholdInE = smearedThreshold_BPix_L1_->fire();
00504 } else {
00505 thePixelThresholdInE = smearedThreshold_BPix_->fire();
00506 }
00507 } else {
00508 if(lay==1) {
00509 thePixelThresholdInE = theThresholdInE_BPix_L1;
00510 } else {
00511 thePixelThresholdInE = theThresholdInE_BPix;
00512 }
00513 }
00514 } else {
00515 if(addThresholdSmearing) {
00516 thePixelThresholdInE = smearedThreshold_FPix_->fire();
00517 } else {
00518 thePixelThresholdInE = theThresholdInE_FPix;
00519 }
00520 }
00521 }
00522
00523
00524 #ifdef TP_DEBUG
00525
00526 float moduleThickness = pixdet->specificSurface().bounds().thickness();
00527 LogDebug ("PixelDigitizer")
00528 << " PixelDigitizer "
00529 << numColumns << " " << numRows << " " << moduleThickness;
00530 #endif
00531
00532 if(addNoise) add_noise(pixdet, thePixelThresholdInE/theNoiseInElectrons);
00533
00534
00535
00536 if((AddPixelInefficiency) && (theSignal.size()>0))
00537 pixel_inefficiency(pixelEfficiencies_, pixdet, tTopo);
00538
00539 if(use_ineff_from_db_ && (theSignal.size()>0))
00540 pixel_inefficiency_db(detID);
00541
00542 if(use_module_killing_) {
00543 if (use_deadmodule_DB_) {
00544 module_killing_DB(detID);
00545 } else {
00546 module_killing_conf(detID);
00547 }
00548 }
00549
00550 make_digis(thePixelThresholdInE, detID, digis, simlinks, tTopo);
00551
00552 #ifdef TP_DEBUG
00553 LogDebug ("PixelDigitizer") << "[SiPixelDigitizerAlgorithm] converted " << digis.size() << " PixelDigis in DetUnit" << detID;
00554 #endif
00555 }
00556
00557
00558
00559
00560 void SiPixelDigitizerAlgorithm::primary_ionization(const PSimHit& hit, std::vector<EnergyDepositUnit>& ionization_points) const {
00561
00562
00563
00564 const float SegmentLength = 0.0010;
00565 float energy;
00566
00567
00568 LocalVector direction = hit.exitPoint() - hit.entryPoint();
00569
00570 float eLoss = hit.energyLoss();
00571 float length = direction.mag();
00572
00573 int NumberOfSegments = int ( length / SegmentLength);
00574 if(NumberOfSegments < 1) NumberOfSegments = 1;
00575
00576 #ifdef TP_DEBUG
00577 LogDebug ("Pixel Digitizer")
00578 << " enter primary_ionzation " << NumberOfSegments
00579 << " shift = "
00580 << (hit.exitPoint().x()-hit.entryPoint().x()) << " "
00581 << (hit.exitPoint().y()-hit.entryPoint().y()) << " "
00582 << (hit.exitPoint().z()-hit.entryPoint().z()) << " "
00583 << hit.particleType() <<" "<< hit.pabs() ;
00584 #endif
00585
00586 float* elossVector = new float[NumberOfSegments];
00587
00588 if( fluctuateCharge ) {
00589
00590 int pid = hit.particleType();
00591
00592
00593 float momentum = hit.pabs();
00594
00595 fluctuateEloss(pid, momentum, eLoss, length, NumberOfSegments,
00596 elossVector);
00597 }
00598
00599 ionization_points.resize( NumberOfSegments);
00600
00601
00602 for ( int i = 0; i != NumberOfSegments; i++) {
00603
00604 Local3DPoint point = hit.entryPoint() +
00605 float((i+0.5)/NumberOfSegments) * direction;
00606
00607 if( fluctuateCharge )
00608 energy = elossVector[i]/GeVperElectron;
00609 else
00610 energy = hit.energyLoss()/GeVperElectron/float(NumberOfSegments);
00611
00612 EnergyDepositUnit edu( energy, point);
00613 ionization_points[i] = edu;
00614
00615 #ifdef TP_DEBUG
00616 LogDebug ("Pixel Digitizer")
00617 << i << " " << ionization_points[i].x() << " "
00618 << ionization_points[i].y() << " "
00619 << ionization_points[i].z() << " "
00620 << ionization_points[i].energy();
00621 #endif
00622
00623 }
00624
00625 delete[] elossVector;
00626
00627 }
00628
00629
00630
00631
00632 void SiPixelDigitizerAlgorithm::fluctuateEloss(int pid, float particleMomentum,
00633 float eloss, float length,
00634 int NumberOfSegs,float elossVector[]) const {
00635
00636
00637
00638
00639
00640
00641 double particleMass = 139.6;
00642 pid = std::abs(pid);
00643 if(pid!=211) {
00644 if(pid==11) particleMass = 0.511;
00645 else if(pid==13) particleMass = 105.7;
00646 else if(pid==321) particleMass = 493.7;
00647 else if(pid==2212) particleMass = 938.3;
00648 }
00649
00650 float segmentLength = length/NumberOfSegs;
00651
00652
00653 float de=0.;
00654 float sum=0.;
00655 double segmentEloss = (1000.*eloss)/NumberOfSegs;
00656 for (int i=0;i<NumberOfSegs;i++) {
00657
00658
00659
00660
00661
00662 double deltaCutoff = tMax;
00663 de = fluctuate->SampleFluctuations(double(particleMomentum*1000.),
00664 particleMass, deltaCutoff,
00665 double(segmentLength*10.),
00666 segmentEloss )/1000.;
00667
00668 elossVector[i]=de;
00669 sum +=de;
00670 }
00671
00672 if(sum>0.) {
00673
00674 float ratio = eloss/sum;
00675
00676 for (int ii=0;ii<NumberOfSegs;ii++) elossVector[ii]= ratio*elossVector[ii];
00677 } else {
00678 float averageEloss = eloss/NumberOfSegs;
00679 for (int ii=0;ii<NumberOfSegs;ii++) elossVector[ii]= averageEloss;
00680 }
00681 return;
00682 }
00683
00684
00685
00686
00687 void SiPixelDigitizerAlgorithm::drift(const PSimHit& hit,
00688 const PixelGeomDetUnit* pixdet,
00689 const GlobalVector& bfield,
00690 const std::vector<EnergyDepositUnit>& ionization_points,
00691 std::vector<SignalPoint>& collection_points) const {
00692
00693 #ifdef TP_DEBUG
00694 LogDebug ("Pixel Digitizer") << " enter drift " ;
00695 #endif
00696
00697 collection_points.resize(ionization_points.size());
00698
00699 LocalVector driftDir=DriftDirection(pixdet, bfield, hit.detUnitId());
00700 if(driftDir.z() ==0.) {
00701 LogWarning("Magnetic field") << " pxlx: drift in z is zero ";
00702 return;
00703 }
00704
00705
00706
00707
00708
00709 float TanLorenzAngleX, TanLorenzAngleY,dir_z, CosLorenzAngleX,
00710 CosLorenzAngleY;
00711 if( alpha2Order) {
00712
00713 TanLorenzAngleX = driftDir.x();
00714 TanLorenzAngleY = driftDir.y();
00715 dir_z = driftDir.z();
00716 CosLorenzAngleX = 1./sqrt(1.+TanLorenzAngleX*TanLorenzAngleX);
00717 CosLorenzAngleY = 1./sqrt(1.+TanLorenzAngleY*TanLorenzAngleY);
00718
00719 } else{
00720
00721 TanLorenzAngleX = driftDir.x();
00722 TanLorenzAngleY = 0.;
00723 dir_z = driftDir.z();
00724 CosLorenzAngleX = 1./sqrt(1.+TanLorenzAngleX*TanLorenzAngleX);
00725 CosLorenzAngleY = 1.;
00726 }
00727
00728 float moduleThickness = pixdet->specificSurface().bounds().thickness();
00729 #ifdef TP_DEBUG
00730 LogDebug ("Pixel Digitizer")
00731 << " Lorentz Tan " << TanLorenzAngleX << " " << TanLorenzAngleY <<" "
00732 << CosLorenzAngleX << " " << CosLorenzAngleY << " "
00733 << moduleThickness*TanLorenzAngleX << " " << driftDir;
00734 #endif
00735
00736 float Sigma_x = 1.;
00737 float Sigma_y = 1.;
00738 float DriftDistance;
00739 float DriftLength;
00740 float Sigma;
00741
00742
00743 for (unsigned int i = 0; i != ionization_points.size(); i++) {
00744
00745 float SegX, SegY, SegZ;
00746 SegX = ionization_points[i].x();
00747 SegY = ionization_points[i].y();
00748 SegZ = ionization_points[i].z();
00749
00750
00751
00752
00753 DriftDistance = moduleThickness/2. - (dir_z * SegZ);
00754
00755
00756
00757
00758
00759
00760
00761
00762
00763 if( DriftDistance < 0.) {
00764 DriftDistance = 0.;
00765 } else if( DriftDistance > moduleThickness )
00766 DriftDistance = moduleThickness;
00767
00768
00769 float XDriftDueToMagField = DriftDistance * TanLorenzAngleX;
00770 float YDriftDueToMagField = DriftDistance * TanLorenzAngleY;
00771
00772
00773 float CloudCenterX = SegX + XDriftDueToMagField;
00774 float CloudCenterY = SegY + YDriftDueToMagField;
00775
00776
00777 DriftLength = sqrt( DriftDistance*DriftDistance +
00778 XDriftDueToMagField*XDriftDueToMagField +
00779 YDriftDueToMagField*YDriftDueToMagField );
00780
00781
00782 Sigma = sqrt(DriftLength/Dist300) * Sigma0;
00783
00784
00785 Sigma_x = Sigma / CosLorenzAngleX ;
00786 Sigma_y = Sigma / CosLorenzAngleY ;
00787
00788
00789 float energyOnCollector = ionization_points[i].energy();
00790
00791 if (pseudoRadDamage>=0.001){
00792 float moduleRadius = pixdet->surface().position().perp();
00793 if (moduleRadius<=pseudoRadDamageRadius){
00794 float kValue = pseudoRadDamage /(moduleRadius*moduleRadius);
00795
00796
00797 energyOnCollector = energyOnCollector * exp( -1*kValue*DriftDistance/moduleThickness );
00798 }
00799 }
00800 #ifdef TP_DEBUG
00801 LogDebug ("Pixel Digitizer")
00802 <<" Dift DistanceZ= "<<DriftDistance<<" module thickness= "<<moduleThickness
00803 <<" Start Energy= "<<ionization_points[i].energy()<<" Energy after loss= "<<energyOnCollector;
00804 #endif
00805 SignalPoint sp( CloudCenterX, CloudCenterY,
00806 Sigma_x, Sigma_y, hit.tof(), energyOnCollector );
00807
00808
00809 collection_points[i] = (sp);
00810
00811 }
00812
00813 }
00814
00815
00816
00817 void SiPixelDigitizerAlgorithm::induce_signal(const PSimHit& hit,
00818 const PixelGeomDetUnit* pixdet,
00819 const std::vector<SignalPoint>& collection_points) {
00820
00821
00822
00823
00824 const PixelTopology* topol=&pixdet->specificTopology();
00825 uint32_t detID= pixdet->geographicalId().rawId();
00826 signal_map_type& theSignal = _signal[detID];
00827
00828 #ifdef TP_DEBUG
00829 LogDebug ("Pixel Digitizer")
00830 << " enter induce_signal, "
00831 << topol->pitch().first << " " << topol->pitch().second;
00832 #endif
00833
00834
00835 typedef std::map< int, float, std::less<int> > hit_map_type;
00836 hit_map_type hit_signal;
00837
00838
00839 std::map<int, float, std::less<int> > x,y;
00840
00841
00842
00843
00844 for ( std::vector<SignalPoint>::const_iterator i=collection_points.begin();
00845 i != collection_points.end(); ++i) {
00846
00847 float CloudCenterX = i->position().x();
00848 float CloudCenterY = i->position().y();
00849 float SigmaX = i->sigma_x();
00850 float SigmaY = i->sigma_y();
00851 float Charge = i->amplitude();
00852
00853
00854
00855
00856
00857
00858
00859
00860 #ifdef TP_DEBUG
00861 LogDebug ("Pixel Digitizer")
00862 << " cloud " << i->position().x() << " " << i->position().y() << " "
00863 << i->sigma_x() << " " << i->sigma_y() << " " << i->amplitude();
00864 #endif
00865
00866
00867 float CloudRight = CloudCenterX + ClusterWidth*SigmaX;
00868 float CloudLeft = CloudCenterX - ClusterWidth*SigmaX;
00869 float CloudUp = CloudCenterY + ClusterWidth*SigmaY;
00870 float CloudDown = CloudCenterY - ClusterWidth*SigmaY;
00871
00872
00873 LocalPoint PointRightUp = LocalPoint(CloudRight,CloudUp);
00874 LocalPoint PointLeftDown = LocalPoint(CloudLeft,CloudDown);
00875
00876
00877
00878
00879
00880
00881
00882
00883
00884 MeasurementPoint mp = topol->measurementPosition(PointRightUp );
00885
00886 int IPixRightUpX = int( floor( mp.x()));
00887 int IPixRightUpY = int( floor( mp.y()));
00888
00889 #ifdef TP_DEBUG
00890 LogDebug ("Pixel Digitizer") << " right-up " << PointRightUp << " "
00891 << mp.x() << " " << mp.y() << " "
00892 << IPixRightUpX << " " << IPixRightUpY ;
00893 #endif
00894
00895 mp = topol->measurementPosition(PointLeftDown );
00896
00897 int IPixLeftDownX = int( floor( mp.x()));
00898 int IPixLeftDownY = int( floor( mp.y()));
00899
00900 #ifdef TP_DEBUG
00901 LogDebug ("Pixel Digitizer") << " left-down " << PointLeftDown << " "
00902 << mp.x() << " " << mp.y() << " "
00903 << IPixLeftDownX << " " << IPixLeftDownY ;
00904 #endif
00905
00906
00907 int numColumns = topol->ncolumns();
00908 int numRows = topol->nrows();
00909
00910 IPixRightUpX = numRows>IPixRightUpX ? IPixRightUpX : numRows-1 ;
00911 IPixRightUpY = numColumns>IPixRightUpY ? IPixRightUpY : numColumns-1 ;
00912 IPixLeftDownX = 0<IPixLeftDownX ? IPixLeftDownX : 0 ;
00913 IPixLeftDownY = 0<IPixLeftDownY ? IPixLeftDownY : 0 ;
00914
00915 x.clear();
00916 y.clear();
00917
00918
00919 int ix;
00920 for (ix=IPixLeftDownX; ix<=IPixRightUpX; ix++) {
00921 float xUB, xLB, UpperBound, LowerBound;
00922
00923
00924
00925
00926 if(ix == 0 || SigmaX==0. )
00927 LowerBound = 0.;
00928 else {
00929 mp = MeasurementPoint( float(ix), 0.0);
00930 xLB = topol->localPosition(mp).x();
00931 LowerBound = 1-calcQ((xLB-CloudCenterX)/SigmaX);
00932 }
00933
00934 if(ix == numRows-1 || SigmaX==0. )
00935 UpperBound = 1.;
00936 else {
00937 mp = MeasurementPoint( float(ix+1), 0.0);
00938 xUB = topol->localPosition(mp).x();
00939 UpperBound = 1. - calcQ((xUB-CloudCenterX)/SigmaX);
00940 }
00941
00942 float TotalIntegrationRange = UpperBound - LowerBound;
00943 x[ix] = TotalIntegrationRange;
00944
00945
00946
00947 }
00948
00949
00950 int iy;
00951 for (iy=IPixLeftDownY; iy<=IPixRightUpY; iy++) {
00952 float yUB, yLB, UpperBound, LowerBound;
00953
00954 if(iy == 0 || SigmaY==0.)
00955 LowerBound = 0.;
00956 else {
00957 mp = MeasurementPoint( 0.0, float(iy) );
00958 yLB = topol->localPosition(mp).y();
00959 LowerBound = 1. - calcQ((yLB-CloudCenterY)/SigmaY);
00960 }
00961
00962 if(iy == numColumns-1 || SigmaY==0. )
00963 UpperBound = 1.;
00964 else {
00965 mp = MeasurementPoint( 0.0, float(iy+1) );
00966 yUB = topol->localPosition(mp).y();
00967 UpperBound = 1. - calcQ((yUB-CloudCenterY)/SigmaY);
00968 }
00969
00970 float TotalIntegrationRange = UpperBound - LowerBound;
00971 y[iy] = TotalIntegrationRange;
00972
00973
00974 }
00975
00976
00977 int chan;
00978 for (ix=IPixLeftDownX; ix<=IPixRightUpX; ix++) {
00979 for (iy=IPixLeftDownY; iy<=IPixRightUpY; iy++) {
00980
00981 float ChargeFraction = Charge*x[ix]*y[iy];
00982
00983 if( ChargeFraction > 0. ) {
00984 chan = PixelDigi::pixelToChannel( ix, iy);
00985
00986 hit_signal[chan] += ChargeFraction;
00987 }
00988
00989
00990 mp = MeasurementPoint( float(ix), float(iy) );
00991 LocalPoint lp = topol->localPosition(mp);
00992 chan = topol->channel(lp);
00993
00994 #ifdef TP_DEBUG
00995 LogDebug ("Pixel Digitizer")
00996 << " pixel " << ix << " " << iy << " - "<<" "
00997 << chan << " " << ChargeFraction<<" "
00998 << mp.x() << " " << mp.y() <<" "
00999 << lp.x() << " " << lp.y() << " "
01000 << chan;
01001 #endif
01002
01003 }
01004 }
01005
01006
01007
01008
01009
01010
01011
01012
01013
01014
01015
01016
01017
01018
01019
01020
01021
01022
01023
01024
01025 }
01026
01027
01028
01029 for ( hit_map_type::const_iterator im = hit_signal.begin();
01030 im != hit_signal.end(); ++im) {
01031 int chan = (*im).first;
01032 theSignal[chan] += (makeDigiSimLinks_ ? Amplitude( (*im).second, &hit, (*im).second) : Amplitude( (*im).second, (*im).second) ) ;
01033
01034 #ifdef TP_DEBUG
01035 std::pair<int,int> ip = PixelDigi::channelToPixel(chan);
01036 LogDebug ("Pixel Digitizer")
01037 << " pixel " << ip.first << " " << ip.second << " "
01038 << theSignal[chan];
01039 #endif
01040 }
01041
01042 }
01043
01044
01045
01046
01047 void SiPixelDigitizerAlgorithm::make_digis(float thePixelThresholdInE,
01048 uint32_t detID,
01049 std::vector<PixelDigi>& digis,
01050 std::vector<PixelDigiSimLink>& simlinks,
01051 const TrackerTopology *tTopo) const {
01052
01053 #ifdef TP_DEBUG
01054 LogDebug ("Pixel Digitizer") << " make digis "<<" "
01055 << " pixel threshold FPix" << theThresholdInE_FPix << " "
01056 << " pixel threshold BPix" << theThresholdInE_BPix << " "
01057 << " pixel threshold BPix Layer1" << theThresholdInE_BPix_L1 << " "
01058 << " List pixels passing threshold ";
01059 #endif
01060
01061
01062
01063 signalMaps::const_iterator it = _signal.find(detID);
01064 if (it == _signal.end()) {
01065 return;
01066 }
01067
01068 const signal_map_type& theSignal = (*it).second;
01069
01070 for (signal_map_const_iterator i = theSignal.begin(); i != theSignal.end(); ++i) {
01071
01072 float signalInElectrons = (*i).second ;
01073
01074
01075
01076
01077
01078
01079 if( signalInElectrons >= thePixelThresholdInE) {
01080
01081 int chan = (*i).first;
01082 std::pair<int,int> ip = PixelDigi::channelToPixel(chan);
01083 int adc=0;
01084
01085
01086 if(doMissCalibrate) {
01087 int row = ip.first;
01088 int col = ip.second;
01089 adc = int(missCalibrate(detID, col, row, signalInElectrons));
01090 } else {
01091 adc = int( signalInElectrons / theElectronPerADC );
01092 }
01093 adc = std::min(adc, theAdcFullScale);
01094
01095 if (theAdcFullScale!=theAdcFullScaleStack){
01096 unsigned int Sub_detid=DetId(detID).subdetId();
01097 if(Sub_detid == PixelSubdetector::PixelBarrel){
01098 int lay = tTopo->pxbLayer(detID);
01099 if (lay>=theFirstStackLayer) {
01100
01101 if (theAdcFullScaleStack==1) {adc=1;}
01102
01103 if (theAdcFullScaleStack!=1&&theAdcFullScaleStack!=theAdcFullScale) {adc = int (1 + adc * (theAdcFullScaleStack-1)/float(theAdcFullScale) );}
01104 }
01105 }
01106 }
01107 #ifdef TP_DEBUG
01108 LogDebug ("Pixel Digitizer")
01109 << (*i).first << " " << (*i).second << " " << signalInElectrons
01110 << " " << adc << ip.first << " " << ip.second ;
01111 #endif
01112
01113
01114 digis.emplace_back(ip.first, ip.second, adc);
01115
01116 if (makeDigiSimLinks_ && (*i).second.hitInfo()!=0) {
01117
01118 if((*i).second.trackIds().size()>0){
01119 simlink_map simi;
01120 unsigned int il=0;
01121 for( std::vector<unsigned int>::const_iterator itid = (*i).second.trackIds().begin();
01122 itid != (*i).second.trackIds().end(); ++itid) {
01123 simi[*itid].push_back((*i).second.individualampl()[il]);
01124 il++;
01125 }
01126
01127
01128 for( simlink_map::iterator simiiter=simi.begin();
01129 simiiter!=simi.end();
01130 simiiter++){
01131
01132 float sum_samechannel=0;
01133 for (unsigned int iii=0;iii<(*simiiter).second.size();iii++){
01134 sum_samechannel+=(*simiiter).second[iii];
01135 }
01136 float fraction=sum_samechannel/(*i).second;
01137 if(fraction>1.) fraction=1.;
01138 simlinks.emplace_back((*i).first, (*simiiter).first, (*i).second.eventId(), fraction);
01139 }
01140 }
01141 }
01142 }
01143 }
01144 }
01145
01146
01147
01148
01149 void SiPixelDigitizerAlgorithm::add_noise(const PixelGeomDetUnit* pixdet,
01150 float thePixelThreshold) {
01151
01152 #ifdef TP_DEBUG
01153 LogDebug ("Pixel Digitizer") << " enter add_noise " << theNoiseInElectrons;
01154 #endif
01155
01156 uint32_t detID= pixdet->geographicalId().rawId();
01157 signal_map_type& theSignal = _signal[detID];
01158
01159
01160
01161 float theSmearedChargeRMS = 0.0;
01162
01163 for ( signal_map_iterator i = theSignal.begin(); i != theSignal.end(); i++) {
01164
01165 if(addChargeVCALSmearing)
01166 {
01167 if((*i).second < 3000)
01168 {
01169 theSmearedChargeRMS = 543.6 - (*i).second * 0.093;
01170 } else if((*i).second < 6000){
01171 theSmearedChargeRMS = 307.6 - (*i).second * 0.01;
01172 } else{
01173 theSmearedChargeRMS = -432.4 +(*i).second * 0.123;
01174 }
01175
01176
01177 float noise_ChargeVCALSmearing = theSmearedChargeRMS * gaussDistributionVCALNoise_->fire() ;
01178
01179 float noise = gaussDistribution_->fire() ;
01180
01181 if(((*i).second + Amplitude(noise+noise_ChargeVCALSmearing, -1.)) < 0. ) {
01182 (*i).second.set(0);}
01183 else{
01184 (*i).second +=Amplitude(noise+noise_ChargeVCALSmearing, -1.);
01185 }
01186
01187 }
01188 else
01189 {
01190
01191
01192 float noise = gaussDistribution_->fire() ;
01193
01194 if(((*i).second + Amplitude(noise, -1.)) < 0. ) {
01195 (*i).second.set(0);}
01196 else{
01197 (*i).second +=Amplitude(noise, -1.);
01198 }
01199 }
01200
01201 }
01202
01203 if(!addNoisyPixels)
01204 return;
01205
01206 const PixelTopology* topol=&pixdet->specificTopology();
01207 int numColumns = topol->ncolumns();
01208 int numRows = topol->nrows();
01209
01210
01211
01212 int numberOfPixels = (numRows * numColumns);
01213 std::map<int,float, std::less<int> > otherPixels;
01214 std::map<int,float, std::less<int> >::iterator mapI;
01215
01216 theNoiser->generate(numberOfPixels,
01217 thePixelThreshold,
01218 theNoiseInElectrons,
01219 otherPixels );
01220
01221 #ifdef TP_DEBUG
01222 LogDebug ("Pixel Digitizer")
01223 << " Add noisy pixels " << numRows << " "
01224 << numColumns << " " << theNoiseInElectrons << " "
01225 << theThresholdInE_FPix << theThresholdInE_BPix <<" "<< numberOfPixels<<" "
01226 << otherPixels.size() ;
01227 #endif
01228
01229
01230 for (mapI = otherPixels.begin(); mapI!= otherPixels.end(); mapI++) {
01231 int iy = ((*mapI).first) / numRows;
01232 int ix = ((*mapI).first) - (iy*numRows);
01233
01234
01235 if( iy < 0 || iy > (numColumns-1) )
01236 LogWarning ("Pixel Geometry") << " error in iy " << iy ;
01237 if( ix < 0 || ix > (numRows-1) )
01238 LogWarning ("Pixel Geometry") << " error in ix " << ix ;
01239
01240 int chan = PixelDigi::pixelToChannel(ix, iy);
01241
01242 #ifdef TP_DEBUG
01243 LogDebug ("Pixel Digitizer")
01244 <<" Storing noise = " << (*mapI).first << " " << (*mapI).second
01245 << " " << ix << " " << iy << " " << chan ;
01246 #endif
01247
01248 if(theSignal[chan] == 0){
01249
01250 int noise=int( (*mapI).second );
01251 theSignal[chan] = Amplitude (noise, -1.);
01252 }
01253 }
01254 }
01255
01256
01257
01258
01259
01260 void SiPixelDigitizerAlgorithm::pixel_inefficiency(const PixelEfficiencies& eff,
01261 const PixelGeomDetUnit* pixdet,
01262 const TrackerTopology *tTopo) {
01263
01264 uint32_t detID= pixdet->geographicalId().rawId();
01265
01266 signal_map_type& theSignal = _signal[detID];
01267
01268 const PixelTopology* topol=&pixdet->specificTopology();
01269 int numColumns = topol->ncolumns();
01270 int numRows = topol->nrows();
01271
01272
01273 float pixelEfficiency = 1.0;
01274 float columnEfficiency = 1.0;
01275 float chipEfficiency = 1.0;
01276
01277
01278 unsigned int Subid=DetId(detID).subdetId();
01279 if (Subid== PixelSubdetector::PixelBarrel){
01280 int layerIndex=tTopo->pxbLayer(detID);
01281 pixelEfficiency = eff.thePixelEfficiency[layerIndex-1];
01282 columnEfficiency = eff.thePixelColEfficiency[layerIndex-1];
01283 chipEfficiency = eff.thePixelChipEfficiency[layerIndex-1];
01284
01285
01286 if (NumberOfBarrelLayers==3){
01287 if(numColumns>416) LogWarning ("Pixel Geometry") <<" wrong columns in barrel "<<numColumns;
01288 if(numRows>160) LogWarning ("Pixel Geometry") <<" wrong rows in barrel "<<numRows;
01289 }
01290 } else {
01291 unsigned int diskIndex=tTopo->pxfDisk(detID)+eff.FPixIndex;
01292
01293
01294 pixelEfficiency = eff.thePixelEfficiency[diskIndex-1];
01295 columnEfficiency = eff.thePixelColEfficiency[diskIndex-1];
01296 chipEfficiency = eff.thePixelChipEfficiency[diskIndex-1];
01297
01298
01299
01300 if (NumberOfBarrelLayers==3){
01301 if(numColumns>260 || numRows>160) {
01302 if(numColumns>260) LogWarning ("Pixel Geometry") <<" wrong columns in endcaps "<<numColumns;
01303 if(numRows>160) LogWarning ("Pixel Geometry") <<" wrong rows in endcaps "<<numRows;
01304 return;
01305 }
01306 }
01307 }
01308
01309 #ifdef TP_DEBUG
01310 LogDebug ("Pixel Digitizer") << " enter pixel_inefficiency " << pixelEfficiency << " "
01311 << columnEfficiency << " " << chipEfficiency;
01312 #endif
01313
01314
01315
01316 std::auto_ptr<PixelIndices> pIndexConverter(new PixelIndices(numColumns,numRows));
01317
01318 int chipIndex = 0;
01319 int rowROC = 0;
01320 int colROC = 0;
01321 std::map<int, int, std::less<int> >chips, columns;
01322 std::map<int, int, std::less<int> >::iterator iter;
01323
01324
01325
01326 for (signal_map_const_iterator i = theSignal.begin(); i != theSignal.end(); ++i) {
01327
01328 int chan = i->first;
01329 std::pair<int,int> ip = PixelDigi::channelToPixel(chan);
01330 int row = ip.first;
01331 int col = ip.second;
01332
01333 pIndexConverter->transformToROC(col,row,chipIndex,colROC,rowROC);
01334 int dColInChip = pIndexConverter->DColumn(colROC);
01335
01336 int dColInDet = pIndexConverter->DColumnInModule(dColInChip,chipIndex);
01337
01338 chips[chipIndex]++;
01339 columns[dColInDet]++;
01340 }
01341
01342
01343 for ( iter = chips.begin(); iter != chips.end() ; iter++ ) {
01344
01345 float rand = flatDistribution_->fire();
01346 if( rand > chipEfficiency ) chips[iter->first]=0;
01347 }
01348
01349
01350 for ( iter = columns.begin(); iter != columns.end() ; iter++ ) {
01351
01352 float rand = flatDistribution_->fire();
01353 if( rand > columnEfficiency ) columns[iter->first]=0;
01354 }
01355
01356
01357
01358 for(signal_map_iterator i = theSignal.begin();i != theSignal.end(); ++i) {
01359
01360
01361 std::pair<int,int> ip = PixelDigi::channelToPixel(i->first);
01362 int row = ip.first;
01363 int col = ip.second;
01364
01365 pIndexConverter->transformToROC(col,row,chipIndex,colROC,rowROC);
01366 int dColInChip = pIndexConverter->DColumn(colROC);
01367
01368 int dColInDet = pIndexConverter->DColumnInModule(dColInChip,chipIndex);
01369
01370
01371
01372 float rand = flatDistribution_->fire();
01373 if( chips[chipIndex]==0 || columns[dColInDet]==0
01374 || rand>pixelEfficiency ) {
01375
01376 i->second.set(0.);
01377 }
01378
01379 }
01380 }
01381
01382
01383
01384
01385
01386
01387
01388
01389
01390
01391 float SiPixelDigitizerAlgorithm::missCalibrate(uint32_t detID, int col,int row,
01392 const float signalInElectrons) const {
01393
01394
01395
01396
01397
01398
01399
01400
01401
01402
01403 float p0=0.0;
01404 float p1=0.0;
01405 float p2=0.0;
01406 float p3=0.0;
01407
01408 unsigned int Sub_detid=DetId(detID).subdetId();
01409
01410 if(Sub_detid == PixelSubdetector::PixelBarrel){
01411 p0 = BPix_p0;
01412 p1 = BPix_p1;
01413 p2 = BPix_p2;
01414 p3 = BPix_p3;
01415 } else {
01416 p0 = FPix_p0;
01417 p1 = FPix_p1;
01418 p2 = FPix_p2;
01419 p3 = FPix_p3;
01420 }
01421
01422
01423
01424 float newAmp = 0.;
01425
01426
01427 float signal = (signalInElectrons-electronsPerVCAL_Offset)/electronsPerVCAL;
01428
01429
01430 newAmp = p3 + p2 * tanh(p0*signal - p1);
01431
01432
01433
01434
01435
01436
01437
01438
01439
01440
01441
01442
01443
01444
01445
01446
01447
01448
01449
01450
01451
01452
01453
01454
01455
01456
01457
01458
01459
01460
01461
01462
01463 return newAmp;
01464 }
01465
01466
01467
01468
01469
01470
01471
01472
01473 LocalVector SiPixelDigitizerAlgorithm::DriftDirection(const PixelGeomDetUnit* pixdet,
01474 const GlobalVector& bfield,
01475 const DetId& detId) const {
01476 Frame detFrame(pixdet->surface().position(),pixdet->surface().rotation());
01477 LocalVector Bfield=detFrame.toLocal(bfield);
01478
01479 float alpha2_FPix;
01480 float alpha2_BPix;
01481 float alpha2;
01482
01483
01484
01485
01486
01487
01488
01489 float dir_x = 0.0;
01490 float dir_y = 0.0;
01491 float dir_z = 0.0;
01492 float scale = 0.0;
01493
01494 uint32_t detID= pixdet->geographicalId().rawId();
01495
01496 unsigned int Sub_detid=DetId(detID).subdetId();
01497
01498
01499
01500 if(!use_LorentzAngle_DB_){
01501
01502 if( alpha2Order) {
01503 alpha2_FPix = tanLorentzAnglePerTesla_FPix*tanLorentzAnglePerTesla_FPix;
01504 alpha2_BPix = tanLorentzAnglePerTesla_BPix*tanLorentzAnglePerTesla_BPix;
01505 }else {
01506 alpha2_FPix = 0.0;
01507 alpha2_BPix = 0.0;
01508 }
01509
01510 if(Sub_detid == PixelSubdetector::PixelBarrel){
01511 dir_x = -( tanLorentzAnglePerTesla_BPix * Bfield.y() + alpha2_BPix* Bfield.z()* Bfield.x() );
01512 dir_y = +( tanLorentzAnglePerTesla_BPix * Bfield.x() - alpha2_BPix* Bfield.z()* Bfield.y() );
01513 dir_z = -(1 + alpha2_BPix* Bfield.z()*Bfield.z() );
01514 scale = (1 + alpha2_BPix* Bfield.z()*Bfield.z() );
01515
01516 } else {
01517 dir_x = -( tanLorentzAnglePerTesla_FPix * Bfield.y() + alpha2_FPix* Bfield.z()* Bfield.x() );
01518 dir_y = +( tanLorentzAnglePerTesla_FPix * Bfield.x() - alpha2_FPix* Bfield.z()* Bfield.y() );
01519 dir_z = -(1 + alpha2_FPix* Bfield.z()*Bfield.z() );
01520 scale = (1 + alpha2_FPix* Bfield.z()*Bfield.z() );
01521 }
01522 }
01523
01524
01525 if(use_LorentzAngle_DB_){
01526 float lorentzAngle = SiPixelLorentzAngle_->getLorentzAngle(detId);
01527 alpha2 = lorentzAngle * lorentzAngle;
01528
01529 dir_x = -( lorentzAngle * Bfield.y() + alpha2 * Bfield.z()* Bfield.x() );
01530 dir_y = +( lorentzAngle * Bfield.x() - alpha2 * Bfield.z()* Bfield.y() );
01531 dir_z = -(1 + alpha2 * Bfield.z()*Bfield.z() );
01532 scale = (1 + alpha2 * Bfield.z()*Bfield.z() );
01533 }
01534
01535 LocalVector theDriftDirection = LocalVector(dir_x/scale, dir_y/scale, dir_z/scale );
01536
01537 #ifdef TP_DEBUG
01538 LogDebug ("Pixel Digitizer") << " The drift direction in local coordinate is "
01539 << theDriftDirection ;
01540 #endif
01541
01542 return theDriftDirection;
01543 }
01544
01545
01546
01547 void SiPixelDigitizerAlgorithm::pixel_inefficiency_db(uint32_t detID) {
01548 if(!use_ineff_from_db_)
01549 return;
01550
01551 signal_map_type& theSignal = _signal[detID];
01552
01553
01554 for(signal_map_iterator i = theSignal.begin();i != theSignal.end(); ++i) {
01555
01556
01557 std::pair<int,int> ip = PixelDigi::channelToPixel(i->first);
01558 int row = ip.first;
01559 int col = ip.second;
01560
01561 if(theSiPixelGainCalibrationService_->isDead(detID, col, row)){
01562
01563
01564 i->second.set(0.);
01565 }
01566 }
01567 }
01568
01569
01570
01571
01572 void SiPixelDigitizerAlgorithm::module_killing_conf(uint32_t detID) {
01573 if(!use_module_killing_)
01574 return;
01575
01576 bool isbad=false;
01577
01578 Parameters::const_iterator itDeadModules=DeadModules.begin();
01579
01580 int detid = detID;
01581 for(; itDeadModules != DeadModules.end(); ++itDeadModules){
01582 int Dead_detID = itDeadModules->getParameter<int>("Dead_detID");
01583 if(detid == Dead_detID){
01584 isbad=true;
01585 break;
01586 }
01587 }
01588
01589 if(!isbad)
01590 return;
01591
01592 signal_map_type& theSignal = _signal[detID];
01593
01594 std::string Module = itDeadModules->getParameter<std::string>("Module");
01595
01596 if(Module=="whole"){
01597 for(signal_map_iterator i = theSignal.begin();i != theSignal.end(); ++i) {
01598 i->second.set(0.);
01599 }
01600 }
01601
01602 for(signal_map_iterator i = theSignal.begin();i != theSignal.end(); ++i) {
01603 std::pair<int,int> ip = PixelDigi::channelToPixel(i->first);
01604
01605 if(Module=="tbmA" && ip.first>=80 && ip.first<=159){
01606 i->second.set(0.);
01607 }
01608
01609 if( Module=="tbmB" && ip.first<=79){
01610 i->second.set(0.);
01611 }
01612 }
01613 }
01614
01615 void SiPixelDigitizerAlgorithm::module_killing_DB(uint32_t detID) {
01616
01617 if(!use_module_killing_)
01618 return;
01619
01620 bool isbad=false;
01621
01622 std::vector<SiPixelQuality::disabledModuleType>disabledModules = SiPixelBadModule_->getBadComponentList();
01623
01624 SiPixelQuality::disabledModuleType badmodule;
01625
01626 for (size_t id=0;id<disabledModules.size();id++)
01627 {
01628 if(detID==disabledModules[id].DetID){
01629 isbad=true;
01630 badmodule = disabledModules[id];
01631 break;
01632 }
01633 }
01634
01635 if(!isbad)
01636 return;
01637
01638 signal_map_type& theSignal = _signal[detID];
01639
01640
01641 if(badmodule.errorType == 0){
01642
01643 for(signal_map_iterator i = theSignal.begin();i != theSignal.end(); ++i) {
01644 i->second.set(0.);
01645 }
01646 }
01647 else {
01648
01649
01650 std::vector<GlobalPixel> badrocpositions (0);
01651 for(unsigned int j = 0; j < 16; j++){
01652 if(SiPixelBadModule_->IsRocBad(detID, j) == true){
01653
01654 std::vector<CablingPathToDetUnit> path = map_.product()->pathToDetUnit(detID);
01655 typedef std::vector<CablingPathToDetUnit>::const_iterator IT;
01656 for (IT it = path.begin(); it != path.end(); ++it) {
01657 const PixelROC* myroc = map_.product()->findItem(*it);
01658 if( myroc->idInDetUnit() == j) {
01659 LocalPixel::RocRowCol local = { 39, 25};
01660 GlobalPixel global = myroc->toGlobal( LocalPixel(local) );
01661 badrocpositions.push_back(global);
01662 break;
01663 }
01664 }
01665 }
01666 }
01667
01668
01669 for(signal_map_iterator i = theSignal.begin();i != theSignal.end(); ++i) {
01670 std::pair<int,int> ip = PixelDigi::channelToPixel(i->first);
01671
01672 for(std::vector<GlobalPixel>::const_iterator it = badrocpositions.begin(); it != badrocpositions.end(); ++it){
01673 if(it->row >= 80 && ip.first >= 80 ){
01674 if((fabs(ip.second - it->col) < 26) ) {i->second.set(0.);}
01675 else if(it->row==120 && ip.second-it->col==26){i->second.set(0.);}
01676 else if(it->row==119 && it->col-ip.second==26){i->second.set(0.);}
01677 }
01678 else if(it->row < 80 && ip.first < 80 ){
01679 if((fabs(ip.second - it->col) < 26) ){i->second.set(0.);}
01680 else if(it->row==40 && ip.second-it->col==26){i->second.set(0.);}
01681 else if(it->row==39 && it->col-ip.second==26){i->second.set(0.);}
01682 }
01683 }
01684 }
01685 }
01686 }
01687