CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_6_1_1/src/SimTracker/SiPixelDigitizer/plugins/SiPixelDigitizerAlgorithm.cc

Go to the documentation of this file.
00001 //class SiPixelDigitizerAlgorithm SimTracker/SiPixelDigitizer/src/SiPixelDigitizerAlgoithm.cc
00002 
00003 // Original Author Danek Kotlinski
00004 // Ported in CMSSW by  Michele Pioppi-INFN perugia
00005 // Added DB capabilities by F.Blekman, Cornell University
00006 //         Created:  Mon Sep 26 11:08:32 CEST 2005
00007 // Add tof, change AddNoise to tracked. 4/06
00008 // Change drift direction. 6/06 d.k.
00009 // Add the statuis (non-rate dependent) inefficiency.
00010 //     -1 - no ineffciency
00011 //      0 - static inefficency only
00012 //    1,2 - low-lumi rate dependent inefficency added
00013 //     10 - high-lumi inefficiency added
00014 // Adopt the correct drift sign convetion from Morris Swartz. d.k. 8/06
00015 // Add more complex misscalinbration, change kev/e to 3.61, diff=3.7,d.k.9/06
00016 // Add the readout channel electronic noise. d.k. 3/07
00017 // Lower the pixel noise from 500 to 175elec.
00018 // Change the input threshold from noise units to electrons.
00019 // Lower the amount of static dead pixels from 0.01 to 0.001.
00020 // Modify to the new random number services. d.k. 5/07
00021 // Protect against sigma=0 (delta tracks on the surface). d.k.5/07
00022 // Change the TOF cut to lower and upper limit. d.k. 7/07
00023 //
00024 // July 2008: Split Lorentz Angle configuration in BPix/FPix (V. Cuplov)
00025 // tanLorentzAngleperTesla_FPix=0.0912 and tanLorentzAngleperTesla_BPix=0.106
00026 // Sept. 2008: Disable Pixel modules which are declared dead in the configuration python file. (V. Cuplov)
00027 // Oct. 2008: Accessing/Reading the Lorentz angle from the DataBase instead of the cfg file. (V. Cuplov)
00028 // Accessing dead modules from the DB. Implementation done and tested on a test.db
00029 // Do not use this option for now. The PixelQuality Objects are not in the official DB yet.
00030 // Feb. 2009: Split Fpix and Bpix threshold and use official numbers (V. Cuplov)
00031 // ThresholdInElectrons_FPix = 2870 and ThresholdInElectrons_BPix = 3700
00032 // update the electron to VCAL conversion using: VCAL_electrons = VCAL * 65.5 - 414
00033 // Feb. 2009: Threshold gaussian smearing (V. Cuplov)
00034 // March, 2009: changed DB access to *SimRcd objects (to de-couple the DB objects from reco chain) (F. Blekman)
00035 // May, 2009: Pixel charge VCAL smearing. (V. Cuplov)
00036 // November, 2009: new parameterization of the pixel response. (V. Cuplov)
00037 // December, 2009: Fix issue with different compilers.
00038 // October, 2010: Improvement: Removing single dead ROC (V. Cuplov)
00039 // November, 2010: Bug fix in removing TBMB/A half-modules (V. Cuplov)
00040 // February, 2011: Time improvement in DriftDirection()  (J. Bashir Butt)
00041 // June, 2011: Bug Fix for pixels on ROC edges in module_killing_DB() (J. Bashir Butt)
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 //#include "PixelIndices.h"
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 // Accessing dead pixel modules from the DB:
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 // Geometry
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_){// load gain calibration service fromdb...
00108     theSiPixelGainCalibrationService_->setESObjects( es );
00109   }
00110   if(use_deadmodule_DB_) {
00111     es.get<SiPixelQualityRcd>().get(SiPixelBadModule_);
00112   }
00113   if(use_LorentzAngle_DB_) {
00114     // Get Lorentz angle from DB record
00115     es.get<SiPixelLorentzAngleSimRcd>().get(SiPixelLorentzAngle_);
00116   }
00117   //gets the map and geometry from the DB (to kill ROCs)
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")), // boolean to kill or not modules
00130   use_deadmodule_DB_(conf.getParameter<bool>("DeadModules_DB")), // boolean to access dead modules from DB
00131   use_LorentzAngle_DB_(conf.getParameter<bool>("LorentzAngle_DB")), // boolean to access Lorentz angle from DB
00132 
00133   DeadModules(use_deadmodule_DB_ ? Parameters() : conf.getParameter<Parameters>("DeadModules")), // get dead module from cfg file
00134 
00135   // Common pixel parameters
00136   // These are parameters which are not likely to be changed
00137   GeVperElectron(3.61E-09), // 1 electron(3.61eV, 1keV(277e, mod 9/06 d.k.
00138   Sigma0(0.00037),           // Charge diffusion constant 7->3.7
00139   Dist300(0.0300),          //   normalized to 300micron Silicon
00140   alpha2Order(conf.getParameter<bool>("Alpha2Order")),   // switch on/off of E.B effect
00141   ClusterWidth(3.),     // Charge integration spread on the collection plane
00142 
00143   // get external parameters:
00144   // To account for upgrade geometries do not assume the number 
00145   // of layers or disks.
00146   NumberOfBarrelLayers(conf.exists("NumPixelBarrel")?conf.getParameter<int>("NumPixelBarrel"):3),
00147   NumberOfEndcapDisks(conf.exists("NumPixelEndcap")?conf.getParameter<int>("NumPixelEndcap"):2),
00148 
00149   // ADC calibration 1adc count(135e.
00150   // Corresponds to 2adc/kev, 270[e/kev]/135[e/adc](2[adc/kev]
00151   // Be carefull, this parameter is also used in SiPixelDet.cc to
00152   // calculate the noise in adc counts from noise in electrons.
00153   // Both defaults should be the same.
00154   theElectronPerADC(conf.getParameter<double>("ElectronPerAdc")),
00155 
00156   // ADC saturation value, 255(8bit adc.
00157   //theAdcFullScale(conf.getUntrackedParameter<int>("AdcFullScale",255)),
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   // Noise in electrons:
00163   // Pixel cell noise, relevant for generating noisy pixels
00164   theNoiseInElectrons(conf.getParameter<double>("NoiseInElectrons")),
00165 
00166   // Fill readout noise, including all readout chain, relevant for smearing
00167   //theReadoutNoise(conf.getUntrackedParameter<double>("ReadoutNoiseInElec",500.)),
00168   theReadoutNoise(conf.getParameter<double>("ReadoutNoiseInElec")),
00169 
00170   // Pixel threshold in units of noise:
00171   // thePixelThreshold(conf.getParameter<double>("ThresholdInNoiseUnits")),
00172   // Pixel threshold in electron units.
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   // Add threshold gaussian smearing:
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   // electrons to VCAL conversion needed in misscalibrate()
00183   electronsPerVCAL(conf.getParameter<double>("ElectronsPerVcal")),
00184   electronsPerVCAL_Offset(conf.getParameter<double>("ElectronsPerVcal_Offset")),
00185 
00186   //theTofCut 12.5, cut in particle TOD +/- 12.5ns
00187   //theTofCut(conf.getUntrackedParameter<double>("TofCut",12.5)),
00188   theTofLowerCut(conf.getParameter<double>("TofLowerCut")),
00189   theTofUpperCut(conf.getParameter<double>("TofUpperCut")),
00190 
00191   // Get the Lorentz angle from the cfg file:
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   // signal response new parameterization: split Fpix and BPix
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   // Add noise
00207   addNoise(conf.getParameter<bool>("AddNoise")),
00208 
00209   // Smear the pixel charge with a gaussian which RMS is a function of the
00210   // pixel charge (Danek's study)
00211   addChargeVCALSmearing(conf.getParameter<bool>("ChargeVCALSmearing")),
00212 
00213   // Add noisy pixels
00214   addNoisyPixels(conf.getParameter<bool>("AddNoisyPixels")),
00215 
00216   // Fluctuate charge in track subsegments
00217   fluctuateCharge(conf.getUntrackedParameter<bool>("FluctuateCharge",true)),
00218 
00219   // Control the pixel inefficiency
00220   AddPixelInefficiency(conf.getParameter<bool>("AddPixelInefficiencyFromPython")),
00221 
00222   // Add threshold gaussian smearing:
00223   addThresholdSmearing(conf.getParameter<bool>("AddThresholdSmearing")),
00224 
00225   // Get the constants for the miss-calibration studies
00226   doMissCalibrate(conf.getParameter<bool>("MissCalibrate")), // Enable miss-calibration
00227   theGainSmearing(conf.getParameter<double>("GainSmearing")), // sigma of the gain smearing
00228   theOffsetSmearing(conf.getParameter<double>("OffsetSmearing")), //sigma of the offset smearing
00229 
00230   // delta cutoff in MeV, has to be same as in OSCAR(0.030/cmsim=1.0 MeV
00231   //tMax(0.030), // In MeV.
00232   //tMax(conf.getUntrackedParameter<double>("DeltaProductionCut",0.030)),
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   // Threshold smearing with gaussian distribution:
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   // Prepare for the analog amplitude miss-calibration
00272   LogDebug ("PixelDigitizer ")
00273     << " miss-calibrate the pixel amplitude ";
00274 
00275   const bool ReadCalParameters = false;
00276   if(ReadCalParameters) {   // Read the calibration files from file
00277     // read the calibration constants from a file (testing only)
00278     std::ifstream in_file;  // data file pointer
00279     char filename[80] = "phCalibrationFit_C0.dat";
00280 
00281     in_file.open(filename, std::ios::in ); // in C++
00282     if(in_file.bad()) {
00283       cout << " File not found " << endl;
00284       return calmap; // signal error
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     // Read MC tracks
00300     for(int i=0;i<(52*80);i++)  { // loop over tracks
00301       in_file >> par0 >> par1 >> par2 >> par3 >> name >> colid >> rowid;
00302       if(in_file.bad()) { // check for errors
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       //cout << " line " << i << " " <<par0<<" "<<par1<<" "<<par2<<" "<<par3<<" "
00314       //   <<colid<<" "<<rowid<<endl;
00315 
00316       CalParameters onePix;
00317       onePix.p0=par0;
00318       onePix.p1=par1;
00319       onePix.p2=par2;
00320       onePix.p3=par3;
00321 
00322       // Convert ROC pixel index to channel
00323       int chan = PixelIndices::pixelToChannelROC(rowid,colid);
00324       calmap.insert(std::pair<int,CalParameters>(chan,onePix));
00325 
00326       // Testing the index conversion, can be skipped
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     } // pixel loop in a ROC
00332 
00333     cout << " map size  " << calmap.size() <<" max "<<calmap.max_size() << " "
00334          <<calmap.empty()<< endl;
00335 
00336 //     cout << " map size  " << calmap.size()  << endl;
00337 //     map<int,CalParameters,std::less<int> >::iterator ix,it;
00338 //     map<int,CalParameters,std::less<int> >::const_iterator ip;
00339 //     for (ix = calmap.begin(); ix != calmap.end(); ++ix) {
00340 //       int i = (*ix).first;
00341 //       std::pair<int,int> p = channelToPixelROC(i);
00342 //       it  = calmap.find(i);
00343 //       CalParameters y  = (*it).second;
00344 //       CalParameters z = (*ix).second;
00345 //       cout << i <<" "<<p.first<<" "<<p.second<<" "<<y.p0<<" "<<z.p0<<" "<<calmap[i].p0<<endl;
00346 
00347 //       //int dummy=0;
00348 //       //cin>>dummy;
00349 //     }
00350 
00351   } // end if readparameters
00352   return calmap;
00353 } // end initCal()
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   // pixel inefficiency
00362   // Don't use Hard coded values, read inefficiencies in from python or don't use any
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                      // The next is needed for Phase2 Tracker studies
00384                      if (NumberOfBarrelLayers>=5){
00385                         if (NumberOfTotLayers>20){throw cms::Exception("Configuration") <<"SiPixelDigitizer was given more layers than it can handle";}
00386                         // For Phase2 tracker layers just set the outermost BPix inefficiency to 99.9%
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   // the first "NumberOfBarrelLayers" settings [0],[1], ... , [NumberOfBarrelLayers-1] are for the barrel pixels
00399   // the next  "NumberOfEndcapDisks"  settings [NumberOfBarrelLayers],[NumberOfBarrelLayers+1], ... are for the endcaps (undecided how)
00400   if(!AddPixelInefficiency) {  // No inefficiency, all 100% efficient
00401     for (int i=0; i<NumberOfTotLayers;i++) {
00402       thePixelEfficiency[i]     = 1.;  // pixels = 100%
00403       thePixelColEfficiency[i]  = 1.;  // columns = 100%
00404       thePixelChipEfficiency[i] = 1.;  // chips = 100%
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     // produce SignalPoint's for all SimHit's in detector
00415     // Loop over hits
00416 
00417     uint32_t detId = pixdet->geographicalId().rawId();
00418     for (std::vector<PSimHit>::const_iterator ssbegin = inputBegin; ssbegin != inputEnd; ++ssbegin) {
00419       // skip hits not in this detector.
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       // fill collection_points for this SimHit, indpendent of topology
00438       // Check the TOF cut
00439       //if(std::abs( (*ssbegin).tof() )<theTofCut) // old cut
00440       if( ((*ssbegin).tof() >= theTofLowerCut) && ((*ssbegin).tof() <= theTofUpperCut) ) {
00441         primary_ionization(*ssbegin, ionization_points); // fills _ionization_points
00442         drift(*ssbegin, pixdet, bfield, ionization_points, collection_points);  // transforms _ionization_points to collection_points
00443         // compute induced signal on readout elements and add to _signal
00444         induce_signal(*ssbegin, pixdet, collection_points); // *ihit needed only for SimHit<-->Digi link
00445       } //  end if
00446     } // end for
00447 
00448 }
00449 
00450 //============================================================================
00451 void SiPixelDigitizerAlgorithm::digitize(const PixelGeomDetUnit* pixdet,
00452                                          std::vector<PixelDigi>& digis,
00453                                          std::vector<PixelDigiSimLink>& simlinks) {
00454 
00455    // Pixel Efficiency moved from the constructor to this method because
00456    // the information of the det are not available in the constructor
00457    // Effciency parameters. 0 - no inefficiency, 1-low lumi, 10-high lumi
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();  // det module number of cols&rows
00464    int numRows = topol->nrows();
00465 
00466     // Noise already defined in electrons
00467     //thePixelThresholdInE = thePixelThreshold * theNoiseInElectrons ;
00468     // Find the threshold in noise units, needed for the noiser.
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){ // Barrel modules
00476       int lay = PXBDetId(detID).layer();
00477       if(addThresholdSmearing) {
00478         if(lay==1) {
00479           thePixelThresholdInE = smearedThreshold_BPix_L1_->fire(); // gaussian smearing
00480         } else {
00481           thePixelThresholdInE = smearedThreshold_BPix_->fire(); // gaussian smearing
00482         }
00483       } else {
00484         if(lay==1) {
00485           thePixelThresholdInE = theThresholdInE_BPix_L1;
00486         } else {
00487           thePixelThresholdInE = theThresholdInE_BPix; // no smearing
00488         }
00489       }
00490     } else { // Forward disks modules
00491       if(addThresholdSmearing) {
00492         thePixelThresholdInE = smearedThreshold_FPix_->fire(); // gaussian smearing
00493       } else {
00494         thePixelThresholdInE = theThresholdInE_FPix; // no smearing
00495       }
00496     }
00497   }
00498 
00499 
00500 #ifdef TP_DEBUG
00501     // full detector thickness
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);  // generate noise
00509 
00510     // Do only if needed
00511 
00512     if((AddPixelInefficiency) && (theSignal.size()>0))
00513       pixel_inefficiency(pixelEfficiencies_, pixdet); // Kill some pixels
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_) {  // remove dead modules using DB
00520         module_killing_DB(detID);
00521       } else { // remove dead modules using the list in cfg file
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 // Generate primary ionization along the track segment.
00535 // Divide the track into small sub-segments
00536 void SiPixelDigitizerAlgorithm::primary_ionization(const PSimHit& hit, std::vector<EnergyDepositUnit>& ionization_points) const {
00537 
00538 // Straight line approximation for trajectory inside active media
00539 
00540   const float SegmentLength = 0.0010; //10microns in cm
00541   float energy;
00542 
00543   // Get the 3D segment direction vector
00544   LocalVector direction = hit.exitPoint() - hit.entryPoint();
00545 
00546   float eLoss = hit.energyLoss();  // Eloss in GeV
00547   float length = direction.mag();  // Track length in Silicon
00548 
00549   int NumberOfSegments = int ( length / SegmentLength); // Number of segments
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];  // Eloss vector
00563 
00564   if( fluctuateCharge ) {
00565     //MP DA RIMUOVERE ASSOLUTAMENTE
00566     int pid = hit.particleType();
00567     //int pid=211;  // assume it is a pion
00568 
00569     float momentum = hit.pabs();
00570     // Generate fluctuated charge points
00571     fluctuateEloss(pid, momentum, eLoss, length, NumberOfSegments,
00572                    elossVector);
00573   }
00574 
00575   ionization_points.resize( NumberOfSegments); // set size
00576 
00577   // loop over segments
00578   for ( int i = 0; i != NumberOfSegments; i++) {
00579     // Divide the segment into equal length subsegments
00580     Local3DPoint point = hit.entryPoint() +
00581       float((i+0.5)/NumberOfSegments) * direction;
00582 
00583     if( fluctuateCharge )
00584       energy = elossVector[i]/GeVperElectron; // Convert charge to elec.
00585     else
00586       energy = hit.energyLoss()/GeVperElectron/float(NumberOfSegments);
00587 
00588     EnergyDepositUnit edu( energy, point); //define position,energy point
00589     ionization_points[i] = edu; // save
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   }  // end for loop
00600 
00601   delete[] elossVector;
00602 
00603 }
00604 //******************************************************************************
00605 
00606 // Fluctuate the charge comming from a small (10um) track segment.
00607 // Use the G4 routine. For mip pions for the moment.
00608 void SiPixelDigitizerAlgorithm::fluctuateEloss(int pid, float particleMomentum,
00609                                       float eloss, float length,
00610                                       int NumberOfSegs,float elossVector[]) const {
00611 
00612   // Get dedx for this track
00613   //float dedx;
00614   //if( length > 0.) dedx = eloss/length;
00615   //else dedx = eloss;
00616 
00617   double particleMass = 139.6; // Mass in MeV, Assume pion
00618   pid = std::abs(pid);
00619   if(pid!=211) {       // Mass in MeV
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   // What is the track segment length.
00626   float segmentLength = length/NumberOfSegs;
00627 
00628   // Generate charge fluctuations.
00629   float de=0.;
00630   float sum=0.;
00631   double segmentEloss = (1000.*eloss)/NumberOfSegs; //eloss in MeV
00632   for (int i=0;i<NumberOfSegs;i++) {
00633     //       material,*,   momentum,energy,*, *,  mass
00634     //myglandz_(14.,segmentLength,2.,2.,dedx,de,0.14);
00635     // The G4 routine needs momentum in MeV, mass in Mev, delta-cut in MeV,
00636     // track segment length in mm, segment eloss in MeV
00637     // Returns fluctuated eloss in MeV
00638     double deltaCutoff = tMax; // the cutoff is sometimes redefined inside, so fix it.
00639     de = fluctuate->SampleFluctuations(double(particleMomentum*1000.),
00640                                       particleMass, deltaCutoff,
00641                                       double(segmentLength*10.),
00642                                       segmentEloss )/1000.; //convert to GeV
00643 
00644     elossVector[i]=de;
00645     sum +=de;
00646   }
00647 
00648   if(sum>0.) {  // If fluctuations give eloss>0.
00649     // Rescale to the same total eloss
00650     float ratio = eloss/sum;
00651 
00652     for (int ii=0;ii<NumberOfSegs;ii++) elossVector[ii]= ratio*elossVector[ii];
00653   } else {  // If fluctuations gives 0 eloss
00654     float averageEloss = eloss/NumberOfSegs;
00655     for (int ii=0;ii<NumberOfSegs;ii++) elossVector[ii]= averageEloss;
00656   }
00657   return;
00658 }
00659 
00660 //*******************************************************************************
00661 // Drift the charge segments to the sensor surface (collection plane)
00662 // Include the effect of E-field and B-field
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()); // set size
00674 
00675   LocalVector driftDir=DriftDirection(pixdet, bfield, hit.detUnitId());  // get the charge drift direction
00676   if(driftDir.z() ==0.) {
00677     LogWarning("Magnetic field") << " pxlx: drift in z is zero ";
00678     return;
00679   }
00680   
00681   // tangent of Lorentz angle
00682   //float TanLorenzAngleX = driftDir.x()/driftDir.z();
00683   //float TanLorenzAngleY = 0.; // force to 0, driftDir.y()/driftDir.z();
00684   
00685   float TanLorenzAngleX, TanLorenzAngleY,dir_z, CosLorenzAngleX,
00686     CosLorenzAngleY;
00687   if( alpha2Order) {
00688 
00689     TanLorenzAngleX = driftDir.x(); // tangen of Lorentz angle
00690     TanLorenzAngleY = driftDir.y();
00691     dir_z = driftDir.z(); // The z drift direction
00692     CosLorenzAngleX = 1./sqrt(1.+TanLorenzAngleX*TanLorenzAngleX); //cosine
00693     CosLorenzAngleY = 1./sqrt(1.+TanLorenzAngleY*TanLorenzAngleY); //cosine;
00694 
00695   } else{
00696 
00697     TanLorenzAngleX = driftDir.x();
00698     TanLorenzAngleY = 0.; // force to 0, driftDir.y()/driftDir.z();
00699     dir_z = driftDir.z(); // The z drift direction
00700     CosLorenzAngleX = 1./sqrt(1.+TanLorenzAngleX*TanLorenzAngleX); //cosine to estimate the path length
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.;  // Charge spread
00713   float Sigma_y = 1.;
00714   float DriftDistance; // Distance between charge generation and collection
00715   float DriftLength;   // Actual Drift Lentgh
00716   float Sigma;
00717 
00718 
00719   for (unsigned int i = 0; i != ionization_points.size(); i++) {
00720 
00721     float SegX, SegY, SegZ; // position
00722     SegX = ionization_points[i].x();
00723     SegY = ionization_points[i].y();
00724     SegZ = ionization_points[i].z();
00725 
00726     // Distance from the collection plane
00727     //DriftDistance = (moduleThickness/2. + SegZ); // Drift to -z
00728     // Include explixitely the E drift direction (for CMS dir_z=-1)
00729     DriftDistance = moduleThickness/2. - (dir_z * SegZ); // Drift to -z
00730 
00731     //if( DriftDistance <= 0.)
00732     //cout<<" <=0 "<<DriftDistance<<" "<<i<<" "<<SegZ<<" "<<dir_z<<" "
00733     //  <<SegX<<" "<<SegY<<" "<<(moduleThickness/2)<<" "
00734     //  <<ionization_points[i].energy()<<" "
00735     //  <<hit.particleType()<<" "<<hit.pabs()<<" "<<hit.energyLoss()<<" "
00736     //  <<hit.entryPoint()<<" "<<hit.exitPoint()
00737     //  <<std::endl;
00738 
00739     if( DriftDistance < 0.) {
00740       DriftDistance = 0.;
00741     } else if( DriftDistance > moduleThickness )
00742       DriftDistance = moduleThickness;
00743 
00744     // Assume full depletion now, partial depletion will come later.
00745     float XDriftDueToMagField = DriftDistance * TanLorenzAngleX;
00746     float YDriftDueToMagField = DriftDistance * TanLorenzAngleY;
00747 
00748     // Shift cloud center
00749     float CloudCenterX = SegX + XDriftDueToMagField;
00750     float CloudCenterY = SegY + YDriftDueToMagField;
00751 
00752     // Calculate how long is the charge drift path
00753     DriftLength = sqrt( DriftDistance*DriftDistance +
00754                         XDriftDueToMagField*XDriftDueToMagField +
00755                         YDriftDueToMagField*YDriftDueToMagField );
00756 
00757     // What is the charge diffusion after this path
00758     Sigma = sqrt(DriftLength/Dist300) * Sigma0;
00759 
00760     // Project the diffusion sigma on the collection plane
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     // Load the Charge distribution parameters
00768     collection_points[i] = (sp);
00769 
00770   } // loop over ionization points, i.
00771 
00772 } // end drift
00773 
00774 //*************************************************************************
00775 // Induce the signal on the collection plane of the active sensor area.
00776 void SiPixelDigitizerAlgorithm::induce_signal(const PSimHit& hit,
00777                                               const PixelGeomDetUnit* pixdet,
00778                                               const std::vector<SignalPoint>& collection_points) {
00779 
00780   // X  - Rows, Left-Right, 160, (1.6cm)   for barrel
00781   // Y  - Columns, Down-Up, 416, (6.4cm)
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; //OK
00791 #endif
00792 
00793    // local map to store pixels hit by 1 Hit.
00794    typedef std::map< int, float, std::less<int> > hit_map_type;
00795    hit_map_type hit_signal;
00796 
00797    // map to store pixel integrals in the x and in the y directions
00798    std::map<int, float, std::less<int> > x,y;
00799 
00800    // Assign signals to readout channels and store sorted by channel number
00801 
00802    // Iterate over collection points on the collection plane
00803    for ( std::vector<SignalPoint>::const_iterator i=collection_points.begin();
00804          i != collection_points.end(); ++i) {
00805 
00806      float CloudCenterX = i->position().x(); // Charge position in x
00807      float CloudCenterY = i->position().y(); //                 in y
00808      float SigmaX = i->sigma_x();            // Charge spread in x
00809      float SigmaY = i->sigma_y();            //               in y
00810      float Charge = i->amplitude();          // Charge amplitude
00811 
00812 
00813      //if(SigmaX==0 || SigmaY==0) {
00814      //cout<<SigmaX<<" "<<SigmaY
00815      //   << " cloud " << i->position().x() << " " << i->position().y() << " "
00816      //   << i->sigma_x() << " " << i->sigma_y() << " " << i->amplitude()<<std::endl;
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      // Find the maximum cloud spread in 2D plane , assume 3*sigma
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      // Define 2D cloud limit points
00832      LocalPoint PointRightUp  = LocalPoint(CloudRight,CloudUp);
00833      LocalPoint PointLeftDown = LocalPoint(CloudLeft,CloudDown);
00834 
00835      // This points can be located outside the sensor area.
00836      // The conversion to measurement point does not check for that
00837      // so the returned pixel index might be wrong (outside range).
00838      // We rely on the limits check below to fix this.
00839      // But remember whatever we do here THE CHARGE OUTSIDE THE ACTIVE
00840      // PIXEL ARE IS LOST, it should not be collected.
00841 
00842      // Convert the 2D points to pixel indices
00843      MeasurementPoint mp = topol->measurementPosition(PointRightUp ); //OK
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 ); //OK
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      // Check detector limits to correct for pixels outside range.
00866      int numColumns = topol->ncolumns();  // det module number of cols&rows
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(); // clear temporary integration array
00875      y.clear();
00876 
00877      // First integrate cahrge strips in x
00878      int ix; // TT for compatibility
00879      for (ix=IPixLeftDownX; ix<=IPixRightUpX; ix++) {  // loop over x index
00880        float xUB, xLB, UpperBound, LowerBound;
00881 
00882        // Why is set to 0 if ix=0, does it meen that we accept charge
00883        // outside the sensor? CHeck How it was done in ORCA?
00884        //if(ix == 0) LowerBound = 0.;
00885        if(ix == 0 || SigmaX==0. )  // skip for surface segemnts
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; // get strip
00910        x[ix] = TotalIntegrationRange; // save strip integral
00911        //if(SigmaX==0 || SigmaY==0)
00912        //cout<<TotalIntegrationRange<<" "<<ix<<std::endl;
00913 
00914      }
00915 
00916     // Now integarte strips in y
00917     int iy; // TT for compatibility
00918     for (iy=IPixLeftDownY; iy<=IPixRightUpY; iy++) { //loope over y ind
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; // save strip integral
00947       //if(SigmaX==0 || SigmaY==0)
00948       //cout<<TotalIntegrationRange<<" "<<iy<<std::endl;
00949     }
00950 
00951     // Get the 2D charge integrals by folding x and y strips
00952     int chan;
00953     for (ix=IPixLeftDownX; ix<=IPixRightUpX; ix++) {  // loop over x index
00954       for (iy=IPixLeftDownY; iy<=IPixRightUpY; iy++) { //loope over y ind
00955 
00956         float ChargeFraction = Charge*x[ix]*y[iy];
00957 
00958         if( ChargeFraction > 0. ) {
00959           chan = PixelDigi::pixelToChannel( ix, iy);  // Get index
00960           // Load the amplitude
00961           hit_signal[chan] += ChargeFraction;
00962         } // endif
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() << " "  // givex edge position
00975           << chan; // edge belongs to previous ?
00976 #endif
00977 
00978       } // endfor iy
00979     } //endfor ix
00980 
00981 
00982     // Test conversions (THIS IS FOR TESTING ONLY) comment-out.
00983     //     mp = topol->measurementPosition( i->position() ); //OK
00984     //     LocalPoint lp = topol->localPosition(mp);     //OK
00985     //     std::pair<float,float> p = topol->pixel( i->position() );  //OK
00986     //     chan = PixelDigi::pixelToChannel( int(p.first), int(p.second));
00987     //     std::pair<int,int> ip = PixelDigi::channelToPixel(chan);
00988     //     MeasurementPoint mp1 = MeasurementPoint( float(ip.first),
00989     //                                       float(ip.second) );
00990     //     LogDebug ("Pixel Digitizer") << " Test "<< mp.x() << " " << mp.y()
00991     //                           << " "<< lp.x() << " " << lp.y() << " "<<" "
00992     //                           <<p.first <<" "<<p.second<<" "<<chan<< " "
00993     //                           <<" " << ip.first << " " << ip.second << " "
00994     //                           << mp1.x() << " " << mp1.y() << " " //OK
00995     //                           << topol->localPosition(mp1).x() << " "  //OK
00996     //                           << topol->localPosition(mp1).y() << " "
00997     //                           << topol->channel( i->position() ); //OK
00998 
00999 
01000   } // loop over charge distributions
01001 
01002   // Fill the global map with all hit pixels from this event
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 } // end induce_signal
01018 
01019 /***********************************************************************/
01020 
01021 // Build pixels, check threshold, add misscalibration, ...
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   // Loop over hit pixels
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 ;   // signal in electrons
01047 
01048     // Do the miss calibration for calibration studies only.
01049     //if(doMissCalibrate) signalInElectrons = missCalibrate(signalInElectrons)
01050 
01051     // Do only for pixels above threshold
01052 
01053     if( signalInElectrons >= thePixelThresholdInE) { // check threshold
01054 
01055       int chan =  (*i).first;  // channel number
01056       std::pair<int,int> ip = PixelDigi::channelToPixel(chan);
01057       int adc=0;  // ADC count as integer
01058 
01059       // Do the miss calibration for calibration studies only.
01060       if(doMissCalibrate) {
01061         int row = ip.first;  // X in row
01062         int col = ip.second; // Y is in col
01063         adc = int(missCalibrate(detID, col, row, signalInElectrons)); //full misscalib.
01064       } else { // Just do a simple electron->adc conversion
01065         adc = int( signalInElectrons / theElectronPerADC ); // calibrate gain
01066       }
01067       adc = std::min(adc, theAdcFullScale); // Check maximum value
01068 // Calculate layerIndex
01069      if (theAdcFullScale!=theAdcFullScaleStack){
01070         unsigned int Sub_detid=DetId(detID).subdetId();
01071         if(Sub_detid == PixelSubdetector::PixelBarrel){ // Barrel modules
01072           int lay = PXBDetId(detID).layer();
01073           if (lay>=theFirstStackLayer) {
01074             // Set to 1 if over the threshold
01075             if (theAdcFullScaleStack==1) {adc=1;}
01076             // Make it a linear fit to the full scale of the normal adc count.   Start new adc from 1 not zero.
01077             if (theAdcFullScaleStack!=1&&theAdcFullScaleStack!=theAdcFullScale) {adc = int (1 + adc * (theAdcFullScaleStack-1)/float(theAdcFullScale) );}
01078           }
01079         }
01080      } // Only enter this if the Adc changes for the outer layers
01081 #ifdef TP_DEBUG
01082       LogDebug ("Pixel Digitizer")
01083         << (*i).first << " " << (*i).second << " " << signalInElectrons
01084         << " " << adc << ip.first << " " << ip.second ;
01085 #endif
01086 
01087       // Load digis
01088       digis.emplace_back(ip.first, ip.second, adc);
01089 
01090       if (makeDigiSimLinks_ && (*i).second.hitInfo()!=0) {
01091         //digilink
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           //sum the contribution of the same trackid
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 //  Add electronic noise to pixel charge
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   // First add noise to hit pixels
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         // Noise from Vcal smearing:
01151         float noise_ChargeVCALSmearing = theSmearedChargeRMS * gaussDistributionVCALNoise_->fire() ;
01152         // Noise from full readout:
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       } // End if addChargeVCalSmearing
01162          else
01163      {
01164         // Noise: ONLY full READOUT Noise.
01165         // Use here the FULL readout noise, including TBM,ALT,AOH,OPT-REC.
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      } // end if only Noise from full readout
01174 
01175   }
01176 
01177   if(!addNoisyPixels)  // Option to skip noise in non-hit pixels
01178     return;
01179 
01180   const PixelTopology* topol=&pixdet->specificTopology();
01181   int numColumns = topol->ncolumns();  // det module number of cols&rows
01182   int numRows = topol->nrows();
01183 
01184   // Add noise on non-hit pixels
01185   // Use here the pixel noise
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, //thr. in un. of nois
01192                       theNoiseInElectrons, // noise in elec.
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   // Add noisy pixels
01204   for (mapI = otherPixels.begin(); mapI!= otherPixels.end(); mapI++) {
01205     int iy = ((*mapI).first) / numRows;
01206     int ix = ((*mapI).first) - (iy*numRows);
01207 
01208     // Keep for a while for testing.
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       //      float noise = float( (*mapI).second );
01224       int noise=int( (*mapI).second );
01225       theSignal[chan] = Amplitude (noise, -1.);
01226     }
01227   }
01228 }
01229 
01230 /***********************************************************************/
01231 
01232 // Simulate the readout inefficiencies.
01233 // Delete a selected number of single pixels, dcols and rocs.
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();  // det module number of cols&rows
01243   int numRows = topol->nrows();
01244 
01245   // Predefined efficiencies
01246   float pixelEfficiency  = 1.0;
01247   float columnEfficiency = 1.0;
01248   float chipEfficiency   = 1.0;
01249 
01250   // setup the chip indices conversion
01251   unsigned int Subid=DetId(detID).subdetId();
01252   if    (Subid==  PixelSubdetector::PixelBarrel){// barrel layers
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     // This should never happen, but only check if it is not an upgrade geometry
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 {                // forward disks
01264 
01265     // For endcaps take same for each endcap
01266     pixelEfficiency  = eff.thePixelEfficiency[eff.FPixIndex];
01267     columnEfficiency = eff.thePixelColEfficiency[eff.FPixIndex];
01268     chipEfficiency   = eff.thePixelChipEfficiency[eff.FPixIndex];
01269     // Sometimes the forward pixels have wrong size,
01270     // this crashes the index conversion, so exit, but only check if it is not an upgrade geometry
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   } // if barrel/forward
01279 
01280 #ifdef TP_DEBUG
01281   LogDebug ("Pixel Digitizer") << " enter pixel_inefficiency " << pixelEfficiency << " "
01282                                << columnEfficiency << " " << chipEfficiency;
01283 #endif
01284 
01285   // Initilize the index converter
01286   //PixelIndices indexConverter(numColumns,numRows);
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   // Find out the number of columns and rocs hits
01296   // Loop over hit pixels, amplitude in electrons, channel = coded row,col
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;  // X in row
01302     int col = ip.second; // Y is in col
01303     //transform to ROC index coordinates
01304     pIndexConverter->transformToROC(col,row,chipIndex,colROC,rowROC);
01305     int dColInChip = pIndexConverter->DColumn(colROC); // get ROC dcol from ROC col
01306     //dcol in mod
01307     int dColInDet = pIndexConverter->DColumnInModule(dColInChip,chipIndex);
01308 
01309     chips[chipIndex]++;
01310     columns[dColInDet]++;
01311   }
01312 
01313   // Delete some ROC hits.
01314   for ( iter = chips.begin(); iter != chips.end() ; iter++ ) {
01315     //float rand  = RandFlat::shoot();
01316     float rand  = flatDistribution_->fire();
01317     if( rand > chipEfficiency ) chips[iter->first]=0;
01318   }
01319 
01320   // Delete some Dcol hits.
01321   for ( iter = columns.begin(); iter != columns.end() ; iter++ ) {
01322     //float rand  = RandFlat::shoot();
01323     float rand  = flatDistribution_->fire();
01324     if( rand > columnEfficiency ) columns[iter->first]=0;
01325   }
01326 
01327   // Now loop again over pixels to kill some of them.
01328   // Loop over hit pixels, amplitude in electrons, channel = coded row,col
01329   for(signal_map_iterator i = theSignal.begin();i != theSignal.end(); ++i) {
01330 
01331     //    int chan = i->first;
01332     std::pair<int,int> ip = PixelDigi::channelToPixel(i->first);//get pixel pos
01333     int row = ip.first;  // X in row
01334     int col = ip.second; // Y is in col
01335     //transform to ROC index coordinates
01336     pIndexConverter->transformToROC(col,row,chipIndex,colROC,rowROC);
01337     int dColInChip = pIndexConverter->DColumn(colROC); //get ROC dcol from ROC col
01338     //dcol in mod
01339     int dColInDet = pIndexConverter->DColumnInModule(dColInChip,chipIndex);
01340 
01341 
01342     //float rand  = RandFlat::shoot();
01343     float rand  = flatDistribution_->fire();
01344     if( chips[chipIndex]==0 || columns[dColInDet]==0
01345         || rand>pixelEfficiency ) {
01346       // make pixel amplitude =0, pixel will be lost at clusterization
01347       i->second.set(0.); // reset amplitude,
01348     } // end if
01349 
01350   } // end pixel loop
01351 } // end pixel_indefficiency
01352 
01353 //***********************************************************************
01354 
01355 // Fluctuate the gain and offset for the amplitude calibration
01356 // Use gaussian smearing.
01357 //float SiPixelDigitizerAlgorithm::missCalibrate(const float amp) const {
01358   //float gain  = RandGaussQ::shoot(1.,theGainSmearing);
01359   //float offset  = RandGaussQ::shoot(0.,theOffsetSmearing);
01360   //float newAmp = amp * gain + offset;
01361   // More complex misscalibration
01362 float SiPixelDigitizerAlgorithm::missCalibrate(uint32_t detID, int col,int row,
01363                                  const float signalInElectrons) const {
01364   // Central values
01365   //const float p0=0.00352, p1=0.868, p2=112., p3=113.; // pix(0,0,0)
01366   //  const float p0=0.00382, p1=0.886, p2=112.7, p3=113.0; // average roc=0
01367   //const float p0=0.00492, p1=1.998, p2=90.6, p3=134.1; // average roc=6
01368   // Smeared (rms)
01369   //const float s0=0.00020, s1=0.051, s2=5.4, s3=4.4; // average roc=0
01370   //const float s0=0.00015, s1=0.043, s2=3.2, s3=3.1; // col average roc=0
01371 
01372   // Make 2 sets of parameters for Fpix and BPIx:
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){// barrel layers
01382       p0 = BPix_p0;
01383       p1 = BPix_p1;
01384       p2 = BPix_p2;
01385       p3 = BPix_p3;
01386     } else {// forward disks
01387       p0 = FPix_p0;
01388       p1 = FPix_p1;
01389       p2 = FPix_p2;
01390       p3 = FPix_p3;
01391     }
01392 
01393   //  const float electronsPerVCAL = 65.5; // our present VCAL calibration (feb 2009)
01394   //  const float electronsPerVCAL_Offset = -414.0; // our present VCAL calibration (feb 2009)
01395   float newAmp = 0.; //Modified signal
01396 
01397   // Convert electrons to VCAL units
01398   float signal = (signalInElectrons-electronsPerVCAL_Offset)/electronsPerVCAL;
01399 
01400   // Simulate the analog response with fixed parametrization
01401   newAmp = p3 + p2 * tanh(p0*signal - p1);
01402 
01403 
01404   // Use the pixel-by-pixel calibrations
01405   //transform to ROC index coordinates
01406   //int chipIndex=0, colROC=0, rowROC=0;
01407   //std::auto_ptr<PixelIndices> pIndexConverter(new PixelIndices(numColumns,numRows));
01408   //pIndexConverter->transformToROC(col,row,chipIndex,colROC,rowROC);
01409 
01410   // Use calibration from a file
01411   //int chanROC = PixelIndices::pixelToChannelROC(rowROC,colROC); // use ROC coordinates
01412   //float pp0=0, pp1=0,pp2=0,pp3=0;
01413   //map<int,CalParameters,std::less<int> >::const_iterator it=calmap.find(chanROC);
01414   //CalParameters y  = (*it).second;
01415   //pp0 = y.p0;
01416   //pp1 = y.p1;
01417   //pp2 = y.p2;
01418   //pp3 = y.p3;
01419 
01420   //
01421   // Use random smearing
01422   // Randomize the pixel response
01423   //float pp0  = RandGaussQ::shoot(p0,s0);
01424   //float pp1  = RandGaussQ::shoot(p1,s1);
01425   //float pp2  = RandGaussQ::shoot(p2,s2);
01426   //float pp3  = RandGaussQ::shoot(p3,s3);
01427 
01428   //newAmp = pp3 + pp2 * tanh(pp0*signal - pp1); // Final signal
01429 
01430   //cout<<" misscalibrate "<<col<<" "<<row<<" "<<chipIndex<<" "<<colROC<<" "
01431   //  <<rowROC<<" "<<signalInElectrons<<" "<<signal<<" "<<newAmp<<" "
01432   //  <<(signalInElectrons/theElectronPerADC)<<std::endl;
01433 
01434   return newAmp;
01435 }
01436 //******************************************************************************
01437 
01438 // Set the drift direction accoring to the Bfield in local det-unit frame
01439 // Works for both barrel and forward pixels.
01440 // Replace the sign convention to fit M.Swartz's formulaes.
01441 // Configurations for barrel and foward pixels possess different tanLorentzAngleperTesla
01442 // parameter value
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   //float dir_x = -tanLorentzAnglePerTesla * Bfield.y();
01455   //float dir_y = +tanLorentzAnglePerTesla * Bfield.x();
01456   //float dir_z = -1.; // E field always in z direction, so electrons go to -z
01457   // The dir_z has to be +/- 1. !
01458   // LocalVector theDriftDirection = LocalVector(dir_x,dir_y,dir_z);
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   // Read Lorentz angle from cfg file:**************************************************************
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){// barrel layers
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 {// forward disks
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   } // end: Read LA from cfg file.
01494 
01495   //Read Lorentz angle from DB:********************************************************************
01496   if(use_LorentzAngle_DB_){
01497     float lorentzAngle = SiPixelLorentzAngle_->getLorentzAngle(detId);
01498     alpha2 = lorentzAngle * lorentzAngle;
01499     //std::cout << "detID is: "<< it->first <<"The LA per tesla is: "<< it->second << std::std::endl;
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   }// end: Read LA from DataBase.
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   // Loop over hit pixels, amplitude in electrons, channel = coded row,col
01525   for(signal_map_iterator i = theSignal.begin();i != theSignal.end(); ++i) {
01526 
01527     //    int chan = i->first;
01528     std::pair<int,int> ip = PixelDigi::channelToPixel(i->first);//get pixel pos
01529     int row = ip.first;  // X in row
01530     int col = ip.second; // Y is in col
01531     //transform to ROC index coordinates
01532     if(theSiPixelGainCalibrationService_->isDead(detID, col, row)){
01533       //      std::cout << "now in isdead check, row " << detID << " " << col << "," << row << std::std::endl;
01534       // make pixel amplitude =0, pixel will be lost at clusterization
01535       i->second.set(0.); // reset amplitude,
01536     } // end if
01537   } // end pixel loop
01538 } // end pixel_indefficiency
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.); // reset amplitude
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);//get pixel pos
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 // Not SLHC safe for now
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   //std::cout<<"Hit in: "<< detID <<" errorType "<< badmodule.errorType<<" BadRocs="<<std::hex<<SiPixelBadModule_->getBadRocs(detID)<<dec<<" "<<std::endl;
01612   if(badmodule.errorType == 0){ // this is a whole dead module.
01613     
01614     for(signal_map_iterator i = theSignal.begin();i != theSignal.end(); ++i) {
01615       i->second.set(0.); // reset amplitude
01616     }
01617   }
01618   else { // all other module types: half-modules and single ROCs.
01619     // Get Bad ROC position:
01620     //follow the example of getBadRocPositions in CondFormats/SiPixelObjects/src/SiPixelQuality.cc
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};   //corresponding to center of ROC row, col
01631             GlobalPixel global = myroc->toGlobal( LocalPixel(local) );
01632             badrocpositions.push_back(global);
01633             break;
01634           }
01635         }
01636       }
01637     }// end of getBadRocPositions
01638     
01639     
01640     for(signal_map_iterator i = theSignal.begin();i != theSignal.end(); ++i) {
01641       std::pair<int,int> ip = PixelDigi::channelToPixel(i->first);//get pixel pos
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