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