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