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