CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_6/src/SimTracker/SiPixelDigitizer/src/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 <vector>
00043 #include <iostream>
00044 
00045 #include "DataFormats/SiPixelDigi/interface/PixelDigi.h"
00046 #include "SimDataFormats/TrackingHit/interface/PSimHitContainer.h"
00047 #include "SimTracker/SiPixelDigitizer/interface/SiPixelDigitizerAlgorithm.h"
00048 
00049 #include <gsl/gsl_sf_erf.h>
00050 #include "FWCore/Utilities/interface/RandomNumberGenerator.h"
00051 #include "CLHEP/Random/RandGaussQ.h"
00052 #include "CLHEP/Random/RandFlat.h"
00053 
00054 //#include "SimTracker/SiPixelDigitizer/interface/PixelIndices.h"
00055 #include "DataFormats/SiPixelDetId/interface/PixelSubdetector.h"
00056 #include "DataFormats/SiPixelDetId/interface/PXBDetId.h"
00057 
00058 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00059 #include "FWCore/ServiceRegistry/interface/Service.h"
00060 #include "FWCore/Utilities/interface/Exception.h"
00061 #include "CalibTracker/SiPixelESProducers/interface/SiPixelGainCalibrationOfflineSimService.h"
00062 
00063 // Accessing dead pixel modules from the DB:
00064 #include "DataFormats/DetId/interface/DetId.h"
00065 
00066 #include "CondFormats/SiPixelObjects/interface/GlobalPixel.h"
00067 #include "FWCore/Framework/interface/ESHandle.h"
00068 #include "CondFormats/DataRecord/interface/SiPixelFedCablingMapRcd.h"
00069 #include "CondFormats/SiPixelObjects/interface/SiPixelFedCablingMap.h"
00070 #include "CondFormats/SiPixelObjects/interface/SiPixelFedCablingTree.h"
00071 #include "CondFormats/SiPixelObjects/interface/SiPixelFedCabling.h"
00072 #include "CondFormats/SiPixelObjects/interface/PixelROC.h"
00073 #include "CondFormats/SiPixelObjects/interface/LocalPixel.h"
00074 #include "CondFormats/SiPixelObjects/interface/CablingPathToDetUnit.h"
00075 
00076 #include "CondFormats/SiPixelObjects/interface/SiPixelFrameReverter.h"
00077 #include "CondFormats/SiPixelObjects/interface/PixelFEDCabling.h"
00078 #include "CondFormats/SiPixelObjects/interface/PixelFEDLink.h"
00079 #include "DataFormats/FEDRawData/interface/FEDNumbering.h"
00080 #include "DataFormats/SiPixelDetId/interface/PixelBarrelName.h"
00081 #include "DataFormats/SiPixelDetId/interface/PixelEndcapName.h"
00082 
00083 // Geometry
00084 #include "Geometry/TrackerGeometryBuilder/interface/TrackerGeometry.h"
00085 #include "Geometry/Records/interface/TrackerDigiGeometryRecord.h"
00086 #include "Geometry/TrackerGeometryBuilder/interface/PixelGeomDetUnit.h"
00087 #include "Geometry/CommonTopologies/interface/PixelTopology.h"
00088 
00089 #include "CondFormats/SiPixelObjects/interface/PixelROC.h"
00090 
00091 using namespace std;
00092 using namespace edm;
00093 using namespace sipixelobjects;
00094 
00095 #define TP_DEBUG // protect all LogDebug with ifdef. Takes too much CPU
00096 
00097 void SiPixelDigitizerAlgorithm::init(const edm::EventSetup& es){
00098   if(use_ineff_from_db_){// load gain calibration service fromdb...
00099     theSiPixelGainCalibrationService_= new SiPixelGainCalibrationOfflineSimService(conf_);
00100     theSiPixelGainCalibrationService_->setESObjects( es );
00101   }
00102 
00103   fillDeadModules(es); // gets the dead module from config file or DB.
00104   fillLorentzAngle(es); // gets the Lorentz angle from the config file or DB.
00105   fillMapandGeom(es);  //gets the map and geometry from the DB (to kill ROCs)
00106 }
00107 
00108 //=========================================================================
00109 
00110 void SiPixelDigitizerAlgorithm::fillDeadModules(const edm::EventSetup& es){
00111   if(!use_deadmodule_DB_){
00112     DeadModules = conf_.getParameter<Parameters>("DeadModules"); // get dead module from cfg file
00113   }
00114   else{  // Get dead module from DB record
00115     // ESHandle was defined in the header file   edm::ESHandle<SiPixelQuality> SiPixelBadModule_;
00116     es.get<SiPixelQualityRcd>().get(SiPixelBadModule_);
00117   }
00118 }
00119 
00120 void SiPixelDigitizerAlgorithm::fillLorentzAngle(const edm::EventSetup& es){
00121   if(!use_LorentzAngle_DB_){
00122     // Get the Lorentz angle from the cfg file:
00123     tanLorentzAnglePerTesla_FPix=conf_.getParameter<double>("TanLorentzAnglePerTesla_FPix");
00124     tanLorentzAnglePerTesla_BPix=conf_.getParameter<double>("TanLorentzAnglePerTesla_BPix");
00125   }
00126   else {
00127     // Get Lorentz angle from DB record
00128     // ESHandle was defined in the header file edm::ESHandle<SiPixelLorentzAngle> SiPixelLorentzAngle_;
00129     es.get<SiPixelLorentzAngleSimRcd>().get(SiPixelLorentzAngle_);
00130   }
00131 }
00132 
00133 // For ROC killing:
00134 void SiPixelDigitizerAlgorithm::fillMapandGeom(const edm::EventSetup& es){
00135    es.get<SiPixelFedCablingMapRcd>().get(map_);
00136    es.get<TrackerDigiGeometryRecord>().get(geom_);
00137 }
00138 
00139 //=========================================================================
00140 
00141 SiPixelDigitizerAlgorithm::SiPixelDigitizerAlgorithm(const edm::ParameterSet& conf, CLHEP::HepRandomEngine& eng) :
00142   conf_(conf) , fluctuate(0), theNoiser(0), pIndexConverter(0),
00143   use_ineff_from_db_(conf_.getParameter<bool>("useDB")),
00144   use_module_killing_(conf_.getParameter<bool>("killModules")), // boolean to kill or not modules
00145   use_deadmodule_DB_(conf_.getParameter<bool>("DeadModules_DB")), // boolean to access dead modules from DB
00146   use_LorentzAngle_DB_(conf_.getParameter<bool>("LorentzAngle_DB")), // boolean to access Lorentz angle from DB
00147   theSiPixelGainCalibrationService_(0),rndEngine(eng)
00148 {
00149   using std::cout;
00150   using std::endl;
00151 
00152   // Common pixel parameters
00153   // These are parameters which are not likely to be changed
00154   NumberOfSegments = 20; // Default number of track segment divisions
00155   ClusterWidth = 3.;     // Charge integration spread on the collection plane
00156   GeVperElectron = 3.61E-09; // 1 electron =3.61eV, 1keV=277e, mod 9/06 d.k.
00157   Sigma0 = 0.00037;           // Charge diffusion constant 7->3.7
00158   Dist300 = 0.0300;          //   normalized to 300micron Silicon
00159 
00160   alpha2Order = conf_.getParameter<bool>("Alpha2Order");   // switch on/off of E.B effect
00161 
00162   // get external parameters:
00163   // ADC calibration 1adc count = 135e.
00164   // Corresponds to 2adc/kev, 270[e/kev]/135[e/adc]=2[adc/kev]
00165   // Be carefull, this parameter is also used in SiPixelDet.cc to
00166   // calculate the noise in adc counts from noise in electrons.
00167   // Both defaults should be the same.
00168   theElectronPerADC=conf_.getParameter<double>("ElectronPerAdc");
00169 
00170   // ADC saturation value, 255=8bit adc.
00171   //theAdcFullScale=conf_.getUntrackedParameter<int>("AdcFullScale",255);
00172   theAdcFullScale=conf_.getParameter<int>("AdcFullScale");
00173 
00174   // Pixel threshold in units of noise:
00175   // thePixelThreshold=conf_.getParameter<double>("ThresholdInNoiseUnits");
00176   // Pixel threshold in electron units.
00177   theThresholdInE_FPix=conf_.getParameter<double>("ThresholdInElectrons_FPix");
00178   theThresholdInE_BPix=conf_.getParameter<double>("ThresholdInElectrons_BPix");
00179 
00180   // Add threshold gaussian smearing:
00181   addThresholdSmearing = conf_.getParameter<bool>("AddThresholdSmearing");
00182   theThresholdSmearing_FPix = conf_.getParameter<double>("ThresholdSmearing_FPix");
00183   theThresholdSmearing_BPix = conf_.getParameter<double>("ThresholdSmearing_BPix");
00184 
00185   // signal response new parameterization: split Fpix and BPix
00186   FPix_p0 = conf_.getParameter<double>("FPix_SignalResponse_p0");
00187   FPix_p1 = conf_.getParameter<double>("FPix_SignalResponse_p1");
00188   FPix_p2 = conf_.getParameter<double>("FPix_SignalResponse_p2");
00189   FPix_p3 = conf_.getParameter<double>("FPix_SignalResponse_p3");
00190 
00191   BPix_p0 = conf_.getParameter<double>("BPix_SignalResponse_p0");
00192   BPix_p1 = conf_.getParameter<double>("BPix_SignalResponse_p1");
00193   BPix_p2 = conf_.getParameter<double>("BPix_SignalResponse_p2");
00194   BPix_p3 = conf_.getParameter<double>("BPix_SignalResponse_p3");
00195 
00196   // electrons to VCAL conversion needed in misscalibrate()
00197   electronsPerVCAL = conf_.getParameter<double>("ElectronsPerVcal");
00198   electronsPerVCAL_Offset = conf_.getParameter<double>("ElectronsPerVcal_Offset");
00199 
00200   // Add noise
00201   addNoise=conf_.getParameter<bool>("AddNoise");
00202 
00203   // Smear the pixel charge with a gaussian which RMS is a function of the
00204   // pixel charge (Danek's study)
00205   addChargeVCALSmearing=conf_.getParameter<bool>("ChargeVCALSmearing");
00206 
00207   // Add noisy pixels
00208   addNoisyPixels=conf_.getParameter<bool>("AddNoisyPixels");
00209   // Noise in electrons:
00210   // Pixel cell noise, relevant for generating noisy pixels
00211   theNoiseInElectrons=conf_.getParameter<double>("NoiseInElectrons");
00212   // Fill readout noise, including all readout chain, relevant for smearing
00213   //theReadoutNoise=conf_.getUntrackedParameter<double>("ReadoutNoiseInElec",500.);
00214   theReadoutNoise=conf_.getParameter<double>("ReadoutNoiseInElec");
00215 
00216   //theTofCut 12.5, cut in particle TOD +/- 12.5ns
00217   //theTofCut=conf_.getUntrackedParameter<double>("TofCut",12.5);
00218   theTofLowerCut=conf_.getParameter<double>("TofLowerCut");
00219   theTofUpperCut=conf_.getParameter<double>("TofUpperCut");
00220 
00221   // Fluctuate charge in track subsegments
00222   fluctuateCharge=conf_.getUntrackedParameter<bool>("FluctuateCharge",true);
00223 
00224   // delta cutoff in MeV, has to be same as in OSCAR=0.030/cmsim=1.0 MeV
00225   //tMax = 0.030; // In MeV.
00226   //tMax =conf_.getUntrackedParameter<double>("DeltaProductionCut",0.030);
00227   tMax =conf_.getParameter<double>("DeltaProductionCut");
00228 
00229   // Control the pixel inefficiency
00230   thePixelLuminosity=conf_.getParameter<int>("AddPixelInefficiency");
00231 
00232   // Get the constants for the miss-calibration studies
00233   doMissCalibrate=conf_.getParameter<bool>("MissCalibrate"); // Enable miss-calibration
00234   theGainSmearing=conf_.getParameter<double>("GainSmearing"); // sigma of the gain smearing
00235   theOffsetSmearing=conf_.getParameter<double>("OffsetSmearing"); //sigma of the offset smearing
00236 
00237   //pixel inefficiency
00238  if (thePixelLuminosity==-20){
00239                      thePixelColEfficiency[0] = conf_.getParameter<double>("thePixelColEfficiency_BPix1");
00240                      thePixelColEfficiency[1] = conf_.getParameter<double>("thePixelColEfficiency_BPix2");
00241                      thePixelColEfficiency[2] = conf_.getParameter<double>("thePixelColEfficiency_BPix3");
00242                      thePixelColEfficiency[3] = conf_.getParameter<double>("thePixelColEfficiency_FPix1");
00243                      thePixelColEfficiency[4] = conf_.getParameter<double>("thePixelColEfficiency_FPix2"); // Not used, but leave it in in case we want use it to later
00244                      cout<<"\nReading in custom Pixel efficiencies "<<thePixelColEfficiency[0]<<" , "<<thePixelColEfficiency[1]<<" , "
00245                                    <<thePixelColEfficiency[2]<<" , "<<thePixelColEfficiency[3]<<" , "<<thePixelColEfficiency[4]<<"\n";
00246                      if (thePixelColEfficiency[0]<=0.5) {cout <<"\n\nDid you mean to set the Pixel efficiency at "<<thePixelColEfficiency[0]
00247                                                               <<", or did you mean for this to be the inefficiency?\n\n\n";}
00248                      }
00249   // the first 3 settings [0],[1],[2] are for the barrel pixels
00250   // the next  3 settings [3],[4],[5] are for the endcaps (undecided how)
00251 
00252   if(thePixelLuminosity==-1) {  // No inefficiency, all 100% efficient
00253     pixelInefficiency=false;
00254     for (int i=0; i<6;i++) {
00255       thePixelEfficiency[i]     = 1.;  // pixels = 100%
00256       thePixelColEfficiency[i]  = 1.;  // columns = 100%
00257       thePixelChipEfficiency[i] = 1.; // chips = 100%
00258     }
00259 
00260     // include only the static (non rate depedent) efficiency
00261     // Usefull for very low rates (luminosity)
00262   } else if(thePixelLuminosity==0) { // static effciency
00263     pixelInefficiency=true;
00264     // Default efficiencies
00265     for (int i=0; i<6;i++) {
00266       if(i<3) {  // For the barrel
00267         // Assume 1% inefficiency for single pixels,
00268         // this is given by faulty bump-bonding and seus.
00269         thePixelEfficiency[i]     = 1.-0.001;  // pixels = 99.9%
00270         // For columns make 0.1% default.
00271         thePixelColEfficiency[i]  = 1.-0.001;  // columns = 99.9%
00272         // A flat 0.1% inefficiency due to lost rocs
00273         thePixelChipEfficiency[i] = 1.-0.001; // chips = 99.9%
00274       } else { // For the endcaps
00275         // Assume 1% inefficiency for single pixels,
00276         // this is given by faulty bump-bonding and seus.
00277         thePixelEfficiency[i]     = 1.-0.001;  // pixels = 99.9%
00278         // For columns make 0.1% default.
00279         thePixelColEfficiency[i]  = 1.-0.001;  // columns = 99.9%
00280         // A flat 0.1% inefficiency due to lost rocs
00281         thePixelChipEfficiency[i] = 1.-0.001; // chips = 99.9%
00282       }
00283     }
00284 
00285     // Include also luminosity rate dependent inefficieny
00286   } else if(thePixelLuminosity>0) { // Include effciency
00287     pixelInefficiency=true;
00288     // Default efficiencies
00289     for (int i=0; i<6;i++) {
00290       if(i<3) { // For the barrel
00291         // Assume 1% inefficiency for single pixels,
00292         // this is given by faulty bump-bonding and seus.
00293         thePixelEfficiency[i]     = 1.-0.01;  // pixels = 99%
00294         // For columns make 1% default.
00295         thePixelColEfficiency[i]  = 1.-0.01;  // columns = 99%
00296         // A flat 0.25% inefficiency due to lost data packets from TBM
00297         thePixelChipEfficiency[i] = 1.-0.0025; // chips = 99.75%
00298       } else { // For the endcaps
00299         // Assume 1% inefficiency for single pixels,
00300         // this is given by faulty bump-bonding and seus.
00301         thePixelEfficiency[i]     = 1.-0.01;  // pixels = 99%
00302         // For columns make 1% default.
00303         thePixelColEfficiency[i]  = 1.-0.01;  // columns = 99%
00304         // A flat 0.25% inefficiency due to lost data packets from TBM
00305         thePixelChipEfficiency[i] = 1.-0.0025; // chips = 99.75%
00306       }
00307     }
00308 
00309     // Special cases ( High-lumi for 4cm layer) where the readout losses are higher
00310     if(thePixelLuminosity==10) { // For high luminosity, bar layer 1
00311       thePixelColEfficiency[0] = 1.-0.034; // 3.4% for r=4 only
00312       thePixelEfficiency[0]    = 1.-0.015; // 1.5% for r=4
00313     }
00314 
00315   } // end the pixel inefficiency part
00316 
00317   // Init the random number services
00318   if(addNoise || thePixelLuminosity || fluctuateCharge || addThresholdSmearing ) {
00319     gaussDistribution_ = new CLHEP::RandGaussQ(rndEngine, 0., theReadoutNoise);
00320     gaussDistributionVCALNoise_ = new CLHEP::RandGaussQ(rndEngine, 0., 1.);
00321     flatDistribution_ = new CLHEP::RandFlat(rndEngine, 0., 1.);
00322 
00323     if(addNoise) {
00324       theNoiser = new GaussianTailNoiseGenerator(rndEngine);
00325     }
00326 
00327     if(fluctuateCharge) {
00328       fluctuate = new SiG4UniversalFluctuation(rndEngine);
00329     }
00330 
00331     // Threshold smearing with gaussian distribution:
00332     if(addThresholdSmearing) {
00333       smearedThreshold_FPix_ = new CLHEP::RandGaussQ(rndEngine, theThresholdInE_FPix , theThresholdSmearing_FPix);
00334       smearedThreshold_BPix_ = new CLHEP::RandGaussQ(rndEngine, theThresholdInE_BPix , theThresholdSmearing_BPix);
00335     }
00336     } //end Init the random number services
00337 
00338 
00339   // Prepare for the analog amplitude miss-calibration
00340   if(doMissCalibrate) {
00341     LogDebug ("PixelDigitizer ")
00342       << " miss-calibrate the pixel amplitude ";
00343 
00344     const bool ReadCalParameters = false;
00345     if(ReadCalParameters) {   // Read the calibration files from file
00346       // read the calibration constants from a file (testing only)
00347       ifstream in_file;  // data file pointer
00348       char filename[80] = "phCalibrationFit_C0.dat";
00349 
00350       in_file.open(filename, ios::in ); // in C++
00351       if(in_file.bad()) {
00352         cout << " File not found " << endl;
00353         return; // signal error
00354       }
00355       cout << " file opened : " << filename << endl;
00356 
00357       char line[500];
00358       for (int i = 0; i < 3; i++) {
00359         in_file.getline(line, 500,'\n');
00360         cout<<line<<endl;
00361       }
00362 
00363       cout << " test map" << endl;
00364 
00365       float par0,par1,par2,par3;
00366       int colid,rowid;
00367       string name;
00368       // Read MC tracks
00369       for(int i=0;i<(52*80);i++)  { // loop over tracks
00370         in_file >> par0 >> par1 >> par2 >> par3 >> name >> colid
00371                 >> rowid;
00372         if(in_file.bad()) { // check for errors
00373           cerr << "Cannot read data file" << endl;
00374           return;
00375         }
00376         if( in_file.eof() != 0 ) {
00377           cerr << in_file.eof() << " " << in_file.gcount() << " "
00378                << in_file.fail() << " " << in_file.good() << " end of file "
00379                << endl;
00380           return;
00381         }
00382 
00383         //cout << " line " << i << " " <<par0<<" "<<par1<<" "<<par2<<" "<<par3<<" "
00384         //   <<colid<<" "<<rowid<<endl;
00385 
00386         CalParameters onePix;
00387         onePix.p0=par0;
00388         onePix.p1=par1;
00389         onePix.p2=par2;
00390         onePix.p3=par3;
00391 
00392         // Convert ROC pixel index to channel
00393         int chan = PixelIndices::pixelToChannelROC(rowid,colid);
00394         calmap.insert(pair<int,CalParameters>(chan,onePix));
00395 
00396         // Testing the index conversion, can be skipped
00397         pair<int,int> p = PixelIndices::channelToPixelROC(chan);
00398         if(rowid!=p.first) cout<<" wrong channel row "<<rowid<<" "<<p.first<<endl;
00399         if(colid!=p.second) cout<<" wrong channel col "<<colid<<" "<<p.second<<endl;
00400 
00401       } // pixel loop in a ROC
00402 
00403       cout << " map size  " << calmap.size() <<" max "<<calmap.max_size() << " "
00404            <<calmap.empty()<< endl;
00405 
00406 //     cout << " map size  " << calmap.size()  << endl;
00407 //     map<int,CalParameters,less<int> >::iterator ix,it;
00408 //     map<int,CalParameters,less<int> >::const_iterator ip;
00409 //     for (ix = calmap.begin(); ix != calmap.end(); ++ix) {
00410 //       int i = (*ix).first;
00411 //       pair<int,int> p = channelToPixelROC(i);
00412 //       it  = calmap.find(i);
00413 //       CalParameters y  = (*it).second;
00414 //       CalParameters z = (*ix).second;
00415 //       cout << i <<" "<<p.first<<" "<<p.second<<" "<<y.p0<<" "<<z.p0<<" "<<calmap[i].p0<<endl;
00416 
00417 //       //int dummy=0;
00418 //       //cin>>dummy;
00419 //     }
00420 
00421     } // end if readparameters
00422   } // end if missCalibration
00423 
00424   LogInfo ("PixelDigitizer ") <<"SiPixelDigitizerAlgorithm constructed"
00425                               <<"Configuration parameters:"
00426                               << "Threshold/Gain = "
00427                               << "threshold in electron FPix = "
00428                               << theThresholdInE_FPix
00429                               << "threshold in electron BPix = "
00430                               << theThresholdInE_BPix
00431                               <<" " << theElectronPerADC << " " << theAdcFullScale
00432                               << " The delta cut-off is set to " << tMax
00433                               << " pix-inefficiency "<<thePixelLuminosity;
00434 
00435 }
00436 //=========================================================================
00437 SiPixelDigitizerAlgorithm::~SiPixelDigitizerAlgorithm() {
00438 
00439   LogDebug ("PixelDigitizer")<<"SiPixelDigitizerAlgorithm deleted";
00440 
00441   // Destructor
00442   delete gaussDistribution_;
00443   delete gaussDistributionVCALNoise_;
00444   delete flatDistribution_;
00445 
00446   delete theSiPixelGainCalibrationService_;
00447 
00448   // Threshold gaussian smearing:
00449   if(addThresholdSmearing) {
00450     delete smearedThreshold_FPix_;
00451     delete smearedThreshold_BPix_;
00452   }
00453 
00454   if(addNoise) delete theNoiser;
00455   if(fluctuateCharge) delete fluctuate;
00456 
00457 }
00458 
00459 //=========================================================================
00460 edm::DetSet<PixelDigi>::collection_type
00461 SiPixelDigitizerAlgorithm::run(
00462                                const std::vector<PSimHit> &input,
00463                                PixelGeomDetUnit *pixdet,
00464                                GlobalVector bfield) {
00465 
00466   _detp = pixdet; //cache the PixGeomDetUnit
00467   _PixelHits=input; //cache the SimHit
00468   _bfield=bfield; //cache the drift direction
00469 
00470   // Pixel Efficiency moved from the constructor to the method run because
00471   // the information of the det are not available in the constructor
00472   // Effciency parameters. 0 - no inefficiency, 1-low lumi, 10-high lumi
00473 
00474   detID= _detp->geographicalId().rawId();
00475 
00476   _signal.clear();
00477 
00478   // initalization  of pixeldigisimlinks
00479   link_coll.clear();
00480 
00481   //Digitization of the SimHits of a given pixdet
00482   vector<PixelDigi> collector =digitize(pixdet);
00483 
00484   // edm::DetSet<PixelDigi> collector;
00485 
00486 #ifdef TP_DEBUG
00487   LogDebug ("PixelDigitizer") << "[SiPixelDigitizerAlgorithm] converted " << collector.size() << " PixelDigis in DetUnit" << detID;
00488 #endif
00489 
00490   return collector;
00491 }
00492 
00493 //============================================================================
00494 vector<PixelDigi> SiPixelDigitizerAlgorithm::digitize(PixelGeomDetUnit *det){
00495 
00496   if( _PixelHits.size() > 0 || addNoisyPixels) {
00497 
00498     topol=&det->specificTopology(); // cache topology
00499     numColumns = topol->ncolumns();  // det module number of cols&rows
00500     numRows = topol->nrows();
00501 
00502     // full detector thickness
00503     moduleThickness = det->specificSurface().bounds().thickness();
00504 
00505     // The index converter is only needed when inefficiencies or misscalibration
00506     // are simulated.
00507     if((pixelInefficiency>0) || doMissCalibrate ) {  // Init pixel indices
00508       pIndexConverter = new PixelIndices(numColumns,numRows);
00509     }
00510 
00511     // Noise already defined in electrons
00512     //thePixelThresholdInE = thePixelThreshold * theNoiseInElectrons ;
00513     // Find the threshold in noise units, needed for the noiser.
00514 
00515   unsigned int Sub_detid=DetId(detID).subdetId();
00516 
00517   if(theNoiseInElectrons>0.){
00518     if(Sub_detid == PixelSubdetector::PixelBarrel){ // Barrel modules
00519       if(addThresholdSmearing) {
00520         thePixelThresholdInE = smearedThreshold_BPix_->fire(); // gaussian smearing
00521       } else {
00522         thePixelThresholdInE = theThresholdInE_BPix; // no smearing
00523       }
00524 
00525       thePixelThreshold = thePixelThresholdInE/theNoiseInElectrons;
00526 
00527     } else { // Forward disks modules
00528       if(addThresholdSmearing) {
00529         thePixelThresholdInE = smearedThreshold_FPix_->fire(); // gaussian smearing
00530       } else {
00531         thePixelThresholdInE = theThresholdInE_FPix; // no smearing
00532       }
00533 
00534       thePixelThreshold = thePixelThresholdInE/theNoiseInElectrons;
00535 
00536     }
00537   } else {
00538     thePixelThreshold = 0.;
00539   }
00540 
00541 
00542 #ifdef TP_DEBUG
00543     LogDebug ("PixelDigitizer")
00544       << " PixelDigitizer "
00545       << numColumns << " " << numRows << " " << moduleThickness;
00546 #endif
00547 
00548     // produce SignalPoint's for all SimHit's in detector
00549     // Loop over hits
00550 
00551     vector<PSimHit>::const_iterator ssbegin;
00552     for (ssbegin= _PixelHits.begin();ssbegin !=_PixelHits.end(); ++ssbegin) {
00553 
00554 #ifdef TP_DEBUG
00555       LogDebug ("Pixel Digitizer")
00556         << (*ssbegin).particleType() << " " << (*ssbegin).pabs() << " "
00557         << (*ssbegin).energyLoss() << " " << (*ssbegin).tof() << " "
00558         << (*ssbegin).trackId() << " " << (*ssbegin).processType() << " "
00559         << (*ssbegin).detUnitId()
00560         << (*ssbegin).entryPoint() << " " << (*ssbegin).exitPoint() ;
00561 #endif
00562 
00563       _collection_points.clear();  // Clear the container
00564       // fill _collection_points for this SimHit, indpendent of topology
00565       // Check the TOF cut
00566       //if(std::abs( (*ssbegin).tof() )<theTofCut){ // old cut
00567       if( ((*ssbegin).tof() >= theTofLowerCut) && ((*ssbegin).tof() <= theTofUpperCut) ) {
00568         primary_ionization(*ssbegin); // fills _ionization_points
00569         drift(*ssbegin);  // transforms _ionization_points to _collection_points
00570         // compute induced signal on readout elements and add to _signal
00571         induce_signal(*ssbegin); // *ihit needed only for SimHit<-->Digi link
00572       } //  end if
00573     } // end for
00574 
00575     if(addNoise) add_noise();  // generate noise
00576     // Do only if needed
00577 
00578     if((pixelInefficiency>0) && (_signal.size()>0))
00579       pixel_inefficiency(); // Kill some pixels
00580 
00581     if(use_ineff_from_db_ && (_signal.size()>0))
00582       pixel_inefficiency_db();
00583 
00584     delete pIndexConverter;
00585 
00586     if(use_module_killing_ && use_deadmodule_DB_) // remove dead modules using DB
00587       module_killing_DB();
00588 
00589     if(use_module_killing_ && !use_deadmodule_DB_) // remove dead modules using the list in cfg file
00590       module_killing_conf();
00591 
00592   }
00593 
00594   make_digis();
00595 
00596   return internal_coll;
00597 }
00598 
00599 //***********************************************************************/
00600 // Generate primary ionization along the track segment.
00601 // Divide the track into small sub-segments
00602 void SiPixelDigitizerAlgorithm::primary_ionization(const PSimHit& hit) {
00603 
00604 // Straight line approximation for trajectory inside active media
00605 
00606   const float SegmentLength = 0.0010; //10microns in cm
00607   float energy;
00608 
00609   // Get the 3D segment direction vector
00610   LocalVector direction = hit.exitPoint() - hit.entryPoint();
00611 
00612   float eLoss = hit.energyLoss();  // Eloss in GeV
00613   float length = direction.mag();  // Track length in Silicon
00614 
00615   NumberOfSegments = int ( length / SegmentLength); // Number of segments
00616   if(NumberOfSegments < 1) NumberOfSegments = 1;
00617 
00618 #ifdef TP_DEBUG
00619   LogDebug ("Pixel Digitizer")
00620     << " enter primary_ionzation " << NumberOfSegments
00621     << " shift = "
00622     << (hit.exitPoint().x()-hit.entryPoint().x()) << " "
00623     << (hit.exitPoint().y()-hit.entryPoint().y()) << " "
00624     << (hit.exitPoint().z()-hit.entryPoint().z()) << " "
00625     << hit.particleType() <<" "<< hit.pabs() ;
00626 #endif
00627 
00628   float* elossVector = new float[NumberOfSegments];  // Eloss vector
00629 
00630   if( fluctuateCharge ) {
00631     //MP DA RIMUOVERE ASSOLUTAMENTE
00632     int pid = hit.particleType();
00633     //int pid=211;  // assume it is a pion
00634 
00635     float momentum = hit.pabs();
00636     // Generate fluctuated charge points
00637     fluctuateEloss(pid, momentum, eLoss, length, NumberOfSegments,
00638                    elossVector);
00639   }
00640 
00641   _ionization_points.resize( NumberOfSegments); // set size
00642 
00643   // loop over segments
00644   for ( int i = 0; i != NumberOfSegments; i++) {
00645     // Divide the segment into equal length subsegments
00646     Local3DPoint point = hit.entryPoint() +
00647       float((i+0.5)/NumberOfSegments) * direction;
00648 
00649     if( fluctuateCharge )
00650       energy = elossVector[i]/GeVperElectron; // Convert charge to elec.
00651     else
00652       energy = hit.energyLoss()/GeVperElectron/float(NumberOfSegments);
00653 
00654     EnergyDepositUnit edu( energy, point); //define position,energy point
00655     _ionization_points[i] = edu; // save
00656 
00657 #ifdef TP_DEBUG
00658     LogDebug ("Pixel Digitizer")
00659       << i << " " << _ionization_points[i].x() << " "
00660       << _ionization_points[i].y() << " "
00661       << _ionization_points[i].z() << " "
00662       << _ionization_points[i].energy();
00663 #endif
00664 
00665   }  // end for loop
00666 
00667   delete[] elossVector;
00668 
00669 }
00670 //******************************************************************************
00671 
00672 // Fluctuate the charge comming from a small (10um) track segment.
00673 // Use the G4 routine. For mip pions for the moment.
00674 void SiPixelDigitizerAlgorithm::fluctuateEloss(int pid, float particleMomentum,
00675                                       float eloss, float length,
00676                                       int NumberOfSegs,float elossVector[]) {
00677 
00678   // Get dedx for this track
00679   //float dedx;
00680   //if( length > 0.) dedx = eloss/length;
00681   //else dedx = eloss;
00682 
00683   double particleMass = 139.6; // Mass in MeV, Assume pion
00684   pid = abs(pid);
00685   if(pid!=211) {       // Mass in MeV
00686     if(pid==11)        particleMass = 0.511;
00687     else if(pid==13)   particleMass = 105.7;
00688     else if(pid==321)  particleMass = 493.7;
00689     else if(pid==2212) particleMass = 938.3;
00690   }
00691   // What is the track segment length.
00692   float segmentLength = length/NumberOfSegs;
00693 
00694   // Generate charge fluctuations.
00695   float de=0.;
00696   float sum=0.;
00697   double segmentEloss = (1000.*eloss)/NumberOfSegs; //eloss in MeV
00698   for (int i=0;i<NumberOfSegs;i++) {
00699     //       material,*,   momentum,energy,*, *,  mass
00700     //myglandz_(14.,segmentLength,2.,2.,dedx,de,0.14);
00701     // The G4 routine needs momentum in MeV, mass in Mev, delta-cut in MeV,
00702     // track segment length in mm, segment eloss in MeV
00703     // Returns fluctuated eloss in MeV
00704     double deltaCutoff = tMax; // the cutoff is sometimes redefined inside, so fix it.
00705     de = fluctuate->SampleFluctuations(double(particleMomentum*1000.),
00706                                       particleMass, deltaCutoff,
00707                                       double(segmentLength*10.),
00708                                       segmentEloss )/1000.; //convert to GeV
00709 
00710     elossVector[i]=de;
00711     sum +=de;
00712   }
00713 
00714   if(sum>0.) {  // If fluctuations give eloss>0.
00715     // Rescale to the same total eloss
00716     float ratio = eloss/sum;
00717 
00718     for (int ii=0;ii<NumberOfSegs;ii++) elossVector[ii]= ratio*elossVector[ii];
00719   } else {  // If fluctuations gives 0 eloss
00720     float averageEloss = eloss/NumberOfSegs;
00721     for (int ii=0;ii<NumberOfSegs;ii++) elossVector[ii]= averageEloss;
00722   }
00723   return;
00724 }
00725 
00726 //*******************************************************************************
00727 // Drift the charge segments to the sensor surface (collection plane)
00728 // Include the effect of E-field and B-field
00729 void SiPixelDigitizerAlgorithm::drift(const PSimHit& hit){
00730 
00731 
00732 #ifdef TP_DEBUG
00733   LogDebug ("Pixel Digitizer") << " enter drift " ;
00734 #endif
00735 
00736   _collection_points.resize( _ionization_points.size()); // set size
00737 
00738   if(use_LorentzAngle_DB_){
00739     //float LorentzAng = SiPixelLorentzAngle_->getLorentzAngle(hit.detUnitId());
00740     lorentzAngle = SiPixelLorentzAngle_->getLorentzAngle(hit.detUnitId());
00741   } 
00742 
00743   //LocalVector driftDir=DriftDirection(LorentzAng);  // get the charge drift direction
00744   LocalVector driftDir=DriftDirection();  // get the charge drift direction
00745   if(driftDir.z() ==0.) {
00746     LogWarning("Magnetic field") << " pxlx: drift in z is zero ";
00747     return;
00748   }
00749   
00750   // tangent of Lorentz angle
00751   //float TanLorenzAngleX = driftDir.x()/driftDir.z();
00752   //float TanLorenzAngleY = 0.; // force to 0, driftDir.y()/driftDir.z();
00753   
00754   float TanLorenzAngleX, TanLorenzAngleY,dir_z, CosLorenzAngleX,
00755     CosLorenzAngleY;
00756   if( alpha2Order) {
00757 
00758     TanLorenzAngleX = driftDir.x(); // tangen of Lorentz angle
00759     TanLorenzAngleY = driftDir.y();
00760     dir_z = driftDir.z(); // The z drift direction
00761     CosLorenzAngleX = 1./sqrt(1.+TanLorenzAngleX*TanLorenzAngleX); //cosine
00762     CosLorenzAngleY = 1./sqrt(1.+TanLorenzAngleY*TanLorenzAngleY); //cosine;
00763 
00764   } else{
00765 
00766     TanLorenzAngleX = driftDir.x();
00767     TanLorenzAngleY = 0.; // force to 0, driftDir.y()/driftDir.z();
00768     dir_z = driftDir.z(); // The z drift direction
00769     CosLorenzAngleX = 1./sqrt(1.+TanLorenzAngleX*TanLorenzAngleX); //cosine to estimate the path length
00770     CosLorenzAngleY = 1.;
00771   }
00772 
00773 
00774 #ifdef TP_DEBUG
00775   LogDebug ("Pixel Digitizer")
00776     << " Lorentz Tan " << TanLorenzAngleX << " " << TanLorenzAngleY <<" "
00777     << CosLorenzAngleX << " " << CosLorenzAngleY << " "
00778     << moduleThickness*TanLorenzAngleX << " " << driftDir;
00779 #endif
00780 
00781   float Sigma_x = 1.;  // Charge spread
00782   float Sigma_y = 1.;
00783   float DriftDistance; // Distance between charge generation and collection
00784   float DriftLength;   // Actual Drift Lentgh
00785   float Sigma;
00786 
00787 
00788   for (unsigned int i = 0; i != _ionization_points.size(); i++) {
00789 
00790     float SegX, SegY, SegZ; // position
00791     SegX = _ionization_points[i].x();
00792     SegY = _ionization_points[i].y();
00793     SegZ = _ionization_points[i].z();
00794 
00795     // Distance from the collection plane
00796     //DriftDistance = (moduleThickness/2. + SegZ); // Drift to -z
00797     // Include explixitely the E drift direction (for CMS dir_z=-1)
00798     DriftDistance = moduleThickness/2. - (dir_z * SegZ); // Drift to -z
00799 
00800     //if( DriftDistance <= 0.)
00801     //cout<<" <=0 "<<DriftDistance<<" "<<i<<" "<<SegZ<<" "<<dir_z<<" "
00802     //  <<SegX<<" "<<SegY<<" "<<(moduleThickness/2)<<" "
00803     //  <<_ionization_points[i].energy()<<" "
00804     //  <<hit.particleType()<<" "<<hit.pabs()<<" "<<hit.energyLoss()<<" "
00805     //  <<hit.entryPoint()<<" "<<hit.exitPoint()
00806     //  <<endl;
00807 
00808     if( DriftDistance < 0.) {
00809       DriftDistance = 0.;
00810     } else if( DriftDistance > moduleThickness )
00811       DriftDistance = moduleThickness;
00812 
00813     // Assume full depletion now, partial depletion will come later.
00814     float XDriftDueToMagField = DriftDistance * TanLorenzAngleX;
00815     float YDriftDueToMagField = DriftDistance * TanLorenzAngleY;
00816 
00817     // Shift cloud center
00818     float CloudCenterX = SegX + XDriftDueToMagField;
00819     float CloudCenterY = SegY + YDriftDueToMagField;
00820 
00821     // Calculate how long is the charge drift path
00822     DriftLength = sqrt( DriftDistance*DriftDistance +
00823                         XDriftDueToMagField*XDriftDueToMagField +
00824                         YDriftDueToMagField*YDriftDueToMagField );
00825 
00826     // What is the charge diffusion after this path
00827     Sigma = sqrt(DriftLength/Dist300) * Sigma0;
00828 
00829     // Project the diffusion sigma on the collection plane
00830     Sigma_x = Sigma / CosLorenzAngleX ;
00831     Sigma_y = Sigma / CosLorenzAngleY ;
00832 
00833     SignalPoint sp( CloudCenterX, CloudCenterY,
00834      Sigma_x, Sigma_y, hit.tof(), _ionization_points[i].energy() );
00835 
00836     // Load the Charge distribution parameters
00837     _collection_points[i] = (sp);
00838 
00839   } // loop over ionization points, i.
00840 
00841 } // end drift
00842 
00843 //*************************************************************************
00844 // Induce the signal on the collection plane of the active sensor area.
00845 void SiPixelDigitizerAlgorithm::induce_signal( const PSimHit& hit) {
00846 
00847   // X  - Rows, Left-Right, 160, (1.6cm)   for barrel
00848   // Y  - Columns, Down-Up, 416, (6.4cm)
00849 
00850 #ifdef TP_DEBUG
00851     LogDebug ("Pixel Digitizer")
00852       << " enter induce_signal, "
00853       << topol->pitch().first << " " << topol->pitch().second; //OK
00854 #endif
00855 
00856    // local map to store pixels hit by 1 Hit.
00857    typedef map< int, float, less<int> > hit_map_type;
00858    hit_map_type hit_signal;
00859 
00860    // map to store pixel integrals in the x and in the y directions
00861    map<int, float, less<int> > x,y;
00862 
00863    // Assign signals to readout channels and store sorted by channel number
00864 
00865    // Iterate over collection points on the collection plane
00866    for ( vector<SignalPoint>::const_iterator i=_collection_points.begin();
00867          i != _collection_points.end(); i++) {
00868 
00869      float CloudCenterX = i->position().x(); // Charge position in x
00870      float CloudCenterY = i->position().y(); //                 in y
00871      float SigmaX = i->sigma_x();            // Charge spread in x
00872      float SigmaY = i->sigma_y();            //               in y
00873      float Charge = i->amplitude();          // Charge amplitude
00874 
00875 
00876      //if(SigmaX==0 || SigmaY==0) {
00877      //cout<<SigmaX<<" "<<SigmaY
00878      //   << " cloud " << i->position().x() << " " << i->position().y() << " "
00879      //   << i->sigma_x() << " " << i->sigma_y() << " " << i->amplitude()<<endl;
00880      //}
00881 
00882 #ifdef TP_DEBUG
00883        LogDebug ("Pixel Digitizer")
00884          << " cloud " << i->position().x() << " " << i->position().y() << " "
00885          << i->sigma_x() << " " << i->sigma_y() << " " << i->amplitude();
00886 #endif
00887 
00888      // Find the maximum cloud spread in 2D plane , assume 3*sigma
00889      float CloudRight = CloudCenterX + ClusterWidth*SigmaX;
00890      float CloudLeft  = CloudCenterX - ClusterWidth*SigmaX;
00891      float CloudUp    = CloudCenterY + ClusterWidth*SigmaY;
00892      float CloudDown  = CloudCenterY - ClusterWidth*SigmaY;
00893 
00894      // Define 2D cloud limit points
00895      LocalPoint PointRightUp  = LocalPoint(CloudRight,CloudUp);
00896      LocalPoint PointLeftDown = LocalPoint(CloudLeft,CloudDown);
00897 
00898      // This points can be located outside the sensor area.
00899      // The conversion to measurement point does not check for that
00900      // so the returned pixel index might be wrong (outside range).
00901      // We rely on the limits check below to fix this.
00902      // But remember whatever we do here THE CHARGE OUTSIDE THE ACTIVE
00903      // PIXEL ARE IS LOST, it should not be collected.
00904 
00905      // Convert the 2D points to pixel indices
00906      MeasurementPoint mp = topol->measurementPosition(PointRightUp ); //OK
00907 
00908      int IPixRightUpX = int( floor( mp.x()));
00909      int IPixRightUpY = int( floor( mp.y()));
00910 
00911 #ifdef TP_DEBUG
00912      LogDebug ("Pixel Digitizer") << " right-up " << PointRightUp << " "
00913                                   << mp.x() << " " << mp.y() << " "
00914                                   << IPixRightUpX << " " << IPixRightUpY ;
00915 #endif
00916 
00917      mp = topol->measurementPosition(PointLeftDown ); //OK
00918 
00919      int IPixLeftDownX = int( floor( mp.x()));
00920      int IPixLeftDownY = int( floor( mp.y()));
00921 
00922 #ifdef TP_DEBUG
00923      LogDebug ("Pixel Digitizer") << " left-down " << PointLeftDown << " "
00924                                   << mp.x() << " " << mp.y() << " "
00925                                   << IPixLeftDownX << " " << IPixLeftDownY ;
00926 #endif
00927 
00928      // Check detector limits to correct for pixels outside range.
00929      IPixRightUpX = numRows>IPixRightUpX ? IPixRightUpX : numRows-1 ;
00930      IPixRightUpY = numColumns>IPixRightUpY ? IPixRightUpY : numColumns-1 ;
00931      IPixLeftDownX = 0<IPixLeftDownX ? IPixLeftDownX : 0 ;
00932      IPixLeftDownY = 0<IPixLeftDownY ? IPixLeftDownY : 0 ;
00933 
00934      x.clear(); // clear temporary integration array
00935      y.clear();
00936 
00937      // First integrate cahrge strips in x
00938      int ix; // TT for compatibility
00939      for (ix=IPixLeftDownX; ix<=IPixRightUpX; ix++) {  // loop over x index
00940        float xUB, xLB, UpperBound, LowerBound;
00941 
00942        // Why is set to 0 if ix=0, does it meen that we accept charge
00943        // outside the sensor? CHeck How it was done in ORCA?
00944        //if(ix == 0) LowerBound = 0.;
00945        if(ix == 0 || SigmaX==0. )  // skip for surface segemnts
00946          LowerBound = 0.;
00947        else {
00948          mp = MeasurementPoint( float(ix), 0.0);
00949          xLB = topol->localPosition(mp).x();
00950          gsl_sf_result result;
00951          int status = gsl_sf_erf_Q_e( (xLB-CloudCenterX)/SigmaX, &result);
00952          if(status != 0)
00953            LogWarning ("Integration")<<"could not compute gaussian probability";
00954          LowerBound = 1-result.val;
00955        }
00956 
00957        if(ix == numRows-1 || SigmaX==0. )
00958          UpperBound = 1.;
00959        else {
00960          mp = MeasurementPoint( float(ix+1), 0.0);
00961          xUB = topol->localPosition(mp).x();
00962          gsl_sf_result result;
00963          int status = gsl_sf_erf_Q_e( (xUB-CloudCenterX)/SigmaX, &result);
00964          if(status != 0)
00965            LogWarning ("Integration")<<"could not compute gaussian probability";
00966          UpperBound = 1. - result.val;
00967        }
00968 
00969        float   TotalIntegrationRange = UpperBound - LowerBound; // get strip
00970        x[ix] = TotalIntegrationRange; // save strip integral
00971        //if(SigmaX==0 || SigmaY==0)
00972        //cout<<TotalIntegrationRange<<" "<<ix<<endl;
00973 
00974      }
00975 
00976     // Now integarte strips in y
00977     int iy; // TT for compatibility
00978     for (iy=IPixLeftDownY; iy<=IPixRightUpY; iy++) { //loope over y ind
00979       float yUB, yLB, UpperBound, LowerBound;
00980 
00981       if(iy == 0 || SigmaY==0.)
00982         LowerBound = 0.;
00983       else {
00984         mp = MeasurementPoint( 0.0, float(iy) );
00985         yLB = topol->localPosition(mp).y();
00986         gsl_sf_result result;
00987         int status = gsl_sf_erf_Q_e( (yLB-CloudCenterY)/SigmaY, &result);
00988         if(status != 0)
00989           LogWarning ("Integration")<<"could not compute gaussian probability";
00990         LowerBound = 1. - result.val;
00991       }
00992 
00993       if(iy == numColumns-1 || SigmaY==0. )
00994         UpperBound = 1.;
00995       else {
00996         mp = MeasurementPoint( 0.0, float(iy+1) );
00997         yUB = topol->localPosition(mp).y();
00998         gsl_sf_result result;
00999         int status = gsl_sf_erf_Q_e( (yUB-CloudCenterY)/SigmaY, &result);
01000         if(status != 0)
01001           LogWarning ("Integration")<<"could not compute gaussian probability";
01002         UpperBound = 1. - result.val;
01003       }
01004 
01005       float   TotalIntegrationRange = UpperBound - LowerBound;
01006       y[iy] = TotalIntegrationRange; // save strip integral
01007       //if(SigmaX==0 || SigmaY==0)
01008       //cout<<TotalIntegrationRange<<" "<<iy<<endl;
01009     }
01010 
01011     // Get the 2D charge integrals by folding x and y strips
01012     int chan;
01013     for (ix=IPixLeftDownX; ix<=IPixRightUpX; ix++) {  // loop over x index
01014       for (iy=IPixLeftDownY; iy<=IPixRightUpY; iy++) { //loope over y ind
01015 
01016         float ChargeFraction = Charge*x[ix]*y[iy];
01017 
01018         if( ChargeFraction > 0. ) {
01019           chan = PixelDigi::pixelToChannel( ix, iy);  // Get index
01020           // Load the amplitude
01021           hit_signal[chan] += ChargeFraction;
01022         } // endif
01023 
01024 
01025         mp = MeasurementPoint( float(ix), float(iy) );
01026         LocalPoint lp = topol->localPosition(mp);
01027         chan = topol->channel(lp);
01028 
01029 #ifdef TP_DEBUG
01030         LogDebug ("Pixel Digitizer")
01031           << " pixel " << ix << " " << iy << " - "<<" "
01032           << chan << " " << ChargeFraction<<" "
01033           << mp.x() << " " << mp.y() <<" "
01034           << lp.x() << " " << lp.y() << " "  // givex edge position
01035           << chan; // edge belongs to previous ?
01036 #endif
01037 
01038       } // endfor iy
01039     } //endfor ix
01040 
01041 
01042     // Test conversions (THIS IS FOR TESTING ONLY) comment-out.
01043     //     mp = topol->measurementPosition( i->position() ); //OK
01044     //     LocalPoint lp = topol->localPosition(mp);     //OK
01045     //     pair<float,float> p = topol->pixel( i->position() );  //OK
01046     //     chan = PixelDigi::pixelToChannel( int(p.first), int(p.second));
01047     //     pair<int,int> ip = PixelDigi::channelToPixel(chan);
01048     //     MeasurementPoint mp1 = MeasurementPoint( float(ip.first),
01049     //                                       float(ip.second) );
01050     //     LogDebug ("Pixel Digitizer") << " Test "<< mp.x() << " " << mp.y()
01051     //                           << " "<< lp.x() << " " << lp.y() << " "<<" "
01052     //                           <<p.first <<" "<<p.second<<" "<<chan<< " "
01053     //                           <<" " << ip.first << " " << ip.second << " "
01054     //                           << mp1.x() << " " << mp1.y() << " " //OK
01055     //                           << topol->localPosition(mp1).x() << " "  //OK
01056     //                           << topol->localPosition(mp1).y() << " "
01057     //                           << topol->channel( i->position() ); //OK
01058 
01059 
01060   } // loop over charge distributions
01061 
01062   // Fill the global map with all hit pixels from this event
01063 
01064   for ( hit_map_type::const_iterator im = hit_signal.begin();
01065         im != hit_signal.end(); im++) {
01066     _signal[(*im).first] += Amplitude( (*im).second, &hit, (*im).second);
01067 
01068     int chan =  (*im).first;
01069     pair<int,int> ip = PixelDigi::channelToPixel(chan);
01070 
01071 #ifdef TP_DEBUG
01072     LogDebug ("Pixel Digitizer")
01073       << " pixel " << ip.first << " " << ip.second << " "
01074       << _signal[(*im).first];
01075 #endif
01076   }
01077 
01078 } // end induce_signal
01079 
01080 /***********************************************************************/
01081 
01082 // Build pixels, check threshold, add misscalibration, ...
01083 void SiPixelDigitizerAlgorithm::make_digis() {
01084   internal_coll.reserve(50); internal_coll.clear();
01085 
01086 #ifdef TP_DEBUG
01087   LogDebug ("Pixel Digitizer") << " make digis "<<" "
01088                                << " pixel threshold FPix" << theThresholdInE_FPix << " "
01089                                << " pixel threshold BPix" << theThresholdInE_BPix << " "
01090                                << " List pixels passing threshold ";
01091 #endif
01092 
01093   // Loop over hit pixels
01094 
01095   for ( signal_map_iterator i = _signal.begin(); i != _signal.end(); i++) {
01096 
01097         float signalInElectrons = (*i).second ;   // signal in electrons
01098 
01099     // Do the miss calibration for calibration studies only.
01100     //if(doMissCalibrate) signalInElectrons = missCalibrate(signalInElectrons)
01101 
01102     // Do only for pixels above threshold
01103 
01104        if( signalInElectrons >= thePixelThresholdInE) { // check threshold
01105 
01106       int chan =  (*i).first;  // channel number
01107       pair<int,int> ip = PixelDigi::channelToPixel(chan);
01108       int adc=0;  // ADC count as integer
01109 
01110       // Do the miss calibration for calibration studies only.
01111       if(doMissCalibrate) {
01112         int row = ip.first;  // X in row
01113         int col = ip.second; // Y is in col
01114         adc = int(missCalibrate(col,row,signalInElectrons)); //full misscalib.
01115       } else { // Just do a simple electron->adc conversion
01116         adc = int( signalInElectrons / theElectronPerADC ); // calibrate gain
01117       }
01118       adc = min(adc, theAdcFullScale); // Check maximum value
01119 
01120 #ifdef TP_DEBUG
01121       LogDebug ("Pixel Digitizer")
01122         << (*i).first << " " << (*i).second << " " << signalInElectrons
01123         << " " << adc << ip.first << " " << ip.second ;
01124 #endif
01125 
01126       // Load digis
01127       internal_coll.push_back( PixelDigi( ip.first, ip.second, adc));
01128 
01129       //digilink
01130       if((*i).second.hits().size()>0){
01131         simi.clear();
01132         unsigned int il=0;
01133         for( vector<const PSimHit*>::const_iterator ihit = (*i).second.hits().begin();
01134              ihit != (*i).second.hits().end(); ihit++) {
01135           simi[(**ihit).trackId()].push_back((*i).second.individualampl()[il]);
01136           il++;
01137         }
01138 
01139         //sum the contribution of the same trackid
01140         for( simlink_map::iterator simiiter=simi.begin();
01141              simiiter!=simi.end();
01142              simiiter++){
01143 
01144           float sum_samechannel=0;
01145           for (unsigned int iii=0;iii<(*simiiter).second.size();iii++){
01146             sum_samechannel+=(*simiiter).second[iii];
01147           }
01148           float fraction=sum_samechannel/(*i).second;
01149           if(fraction>1.) fraction=1.;
01150           link_coll.push_back(PixelDigiSimLink((*i).first,(*simiiter).first,((*i).second.hits().front())->eventId(),fraction));
01151         }
01152 
01153       }
01154     }
01155 
01156   }
01157 }
01158 
01159 /***********************************************************************/
01160 
01161 //  Add electronic noise to pixel charge
01162 void SiPixelDigitizerAlgorithm::add_noise() {
01163 
01164 #ifdef TP_DEBUG
01165   LogDebug ("Pixel Digitizer") << " enter add_noise " << theNoiseInElectrons;
01166 #endif
01167 
01168   // First add noise to hit pixels
01169 
01170   for ( signal_map_iterator i = _signal.begin(); i != _signal.end(); i++) {
01171 
01172          if(addChargeVCALSmearing)
01173       {
01174         if((*i).second < 3000)
01175           {
01176             theSmearedChargeRMS = 543.6 - (*i).second * 0.093;
01177           } else if((*i).second < 6000){
01178             theSmearedChargeRMS = 307.6 - (*i).second * 0.01;
01179           } else{
01180             theSmearedChargeRMS = -432.4 +(*i).second * 0.123;
01181         }
01182 
01183         // Noise from Vcal smearing:
01184         float noise_ChargeVCALSmearing = theSmearedChargeRMS * gaussDistributionVCALNoise_->fire() ;
01185         // Noise from full readout:
01186         float noise  = gaussDistribution_->fire() ;
01187 
01188                 if(((*i).second + Amplitude(noise+noise_ChargeVCALSmearing,0,-1.)) < 0. ) {
01189                   (*i).second.set(0);}
01190                 else{
01191         (*i).second +=Amplitude(noise+noise_ChargeVCALSmearing,0,-1.);
01192                 }
01193 
01194       } // End if addChargeVCalSmearing
01195          else
01196      {
01197         // Noise: ONLY full READOUT Noise.
01198         // Use here the FULL readout noise, including TBM,ALT,AOH,OPT-REC.
01199 
01200         float noise  = gaussDistribution_->fire() ;
01201                 if(((*i).second + Amplitude(noise,0,-1.)) < 0. ) {
01202                   (*i).second.set(0);}
01203                 else{
01204         (*i).second +=Amplitude(noise,0,-1.);
01205                 }
01206      } // end if only Noise from full readout
01207 
01208   }
01209 
01210   if(!addNoisyPixels)  // Option to skip noise in non-hit pixels
01211     return;
01212 
01213   // Add noise on non-hit pixels
01214   // Use here the pixel noise
01215   int numberOfPixels = (numRows * numColumns);
01216   map<int,float, less<int> > otherPixels;
01217   map<int,float, less<int> >::iterator mapI;
01218 
01219   theNoiser->generate(numberOfPixels,
01220                       thePixelThreshold, //thr. in un. of nois
01221                       theNoiseInElectrons, // noise in elec.
01222                       otherPixels );
01223 
01224 #ifdef TP_DEBUG
01225   LogDebug ("Pixel Digitizer")
01226     <<  " Add noisy pixels " << numRows << " "
01227     << numColumns << " " << theNoiseInElectrons << " "
01228     << theThresholdInE_FPix << theThresholdInE_BPix <<" "<< numberOfPixels<<" "
01229     << otherPixels.size() ;
01230 #endif
01231 
01232   // Add noisy pixels
01233   for (mapI = otherPixels.begin(); mapI!= otherPixels.end(); mapI++) {
01234     int iy = ((*mapI).first) / numRows;
01235     int ix = ((*mapI).first) - (iy*numRows);
01236 
01237     // Keep for a while for testing.
01238     if( iy < 0 || iy > (numColumns-1) )
01239       LogWarning ("Pixel Geometry") << " error in iy " << iy ;
01240     if( ix < 0 || ix > (numRows-1) )
01241       LogWarning ("Pixel Geometry")  << " error in ix " << ix ;
01242 
01243     int chan = PixelDigi::pixelToChannel(ix, iy);
01244 
01245 #ifdef TP_DEBUG
01246     LogDebug ("Pixel Digitizer")
01247       <<" Storing noise = " << (*mapI).first << " " << (*mapI).second
01248       << " " << ix << " " << iy << " " << chan ;
01249 #endif
01250 
01251     if(_signal[chan] == 0){
01252       //      float noise = float( (*mapI).second );
01253       int noise=int( (*mapI).second );
01254       _signal[chan] = Amplitude (noise, 0,-1.);
01255     }
01256   }
01257 
01258 }
01259 
01260 /***********************************************************************/
01261 
01262 // Simulate the readout inefficiencies.
01263 // Delete a selected number of single pixels, dcols and rocs.
01264 void SiPixelDigitizerAlgorithm::pixel_inefficiency() {
01265 
01266 
01267   // Predefined efficiencies
01268   float pixelEfficiency  = 1.0;
01269   float columnEfficiency = 1.0;
01270   float chipEfficiency   = 1.0;
01271 
01272   // setup the chip indices conversion
01273   unsigned int Subid=DetId(detID).subdetId();
01274   if    (Subid==  PixelSubdetector::PixelBarrel){// barrel layers
01275     int layerIndex=PXBDetId(detID).layer();
01276     pixelEfficiency  = thePixelEfficiency[layerIndex-1];
01277     columnEfficiency = thePixelColEfficiency[layerIndex-1];
01278     chipEfficiency   = thePixelChipEfficiency[layerIndex-1];
01279 
01280     // This should never happen
01281     if(numColumns>416)  LogWarning ("Pixel Geometry") <<" wrong columns in barrel "<<numColumns;
01282     if(numRows>160)  LogWarning ("Pixel Geometry") <<" wrong rows in barrel "<<numRows;
01283 
01284   } else {                // forward disks
01285 
01286     // For endcaps take same for each endcap
01287     pixelEfficiency  = thePixelEfficiency[3];
01288     columnEfficiency = thePixelColEfficiency[3];
01289     chipEfficiency   = thePixelChipEfficiency[3];
01290 
01291     // Sometimes the forward pixels have wrong size,
01292     // this crashes the index conversion, so exit.
01293     if(numColumns>260 || numRows>160) {
01294       if(numColumns>260)  LogWarning ("Pixel Geometry") <<" wrong columns in endcaps "<<numColumns;
01295       if(numRows>160)  LogWarning ("Pixel Geometry") <<" wrong rows in endcaps "<<numRows;
01296       return;
01297     }
01298   } // if barrel/forward
01299 
01300 #ifdef TP_DEBUG
01301   LogDebug ("Pixel Digitizer") << " enter pixel_inefficiency " << pixelEfficiency << " "
01302                                << columnEfficiency << " " << chipEfficiency;
01303 #endif
01304 
01305   // Initilize the index converter
01306   //PixelIndices indexConverter(numColumns,numRows);
01307   int chipIndex = 0;
01308   int rowROC = 0;
01309   int colROC = 0;
01310   map<int, int, less<int> >chips, columns;
01311   map<int, int, less<int> >::iterator iter;
01312 
01313   // Find out the number of columns and rocs hits
01314   // Loop over hit pixels, amplitude in electrons, channel = coded row,col
01315   for (signal_map_iterator i = _signal.begin();i != _signal.end();i++) {
01316 
01317     int chan = i->first;
01318     pair<int,int> ip = PixelDigi::channelToPixel(chan);
01319     int row = ip.first;  // X in row
01320     int col = ip.second; // Y is in col
01321     //transform to ROC index coordinates
01322     pIndexConverter->transformToROC(col,row,chipIndex,colROC,rowROC);
01323     int dColInChip = pIndexConverter->DColumn(colROC); // get ROC dcol from ROC col
01324     //dcol in mod
01325     int dColInDet = pIndexConverter->DColumnInModule(dColInChip,chipIndex);
01326 
01327     chips[chipIndex]++;
01328     columns[dColInDet]++;
01329   }
01330 
01331   // Delete some ROC hits.
01332   for ( iter = chips.begin(); iter != chips.end() ; iter++ ) {
01333     //float rand  = RandFlat::shoot();
01334     float rand  = flatDistribution_->fire();
01335     if( rand > chipEfficiency ) chips[iter->first]=0;
01336   }
01337 
01338   // Delete some Dcol hits.
01339   for ( iter = columns.begin(); iter != columns.end() ; iter++ ) {
01340     //float rand  = RandFlat::shoot();
01341     float rand  = flatDistribution_->fire();
01342     if( rand > columnEfficiency ) columns[iter->first]=0;
01343   }
01344 
01345   // Now loop again over pixels to kill some of them.
01346   // Loop over hit pixels, amplitude in electrons, channel = coded row,col
01347   for(signal_map_iterator i = _signal.begin();i != _signal.end(); i++) {
01348 
01349     //    int chan = i->first;
01350     pair<int,int> ip = PixelDigi::channelToPixel(i->first);//get pixel pos
01351     int row = ip.first;  // X in row
01352     int col = ip.second; // Y is in col
01353     //transform to ROC index coordinates
01354     pIndexConverter->transformToROC(col,row,chipIndex,colROC,rowROC);
01355     int dColInChip = pIndexConverter->DColumn(colROC); //get ROC dcol from ROC col
01356     //dcol in mod
01357     int dColInDet = pIndexConverter->DColumnInModule(dColInChip,chipIndex);
01358 
01359 
01360     //float rand  = RandFlat::shoot();
01361     float rand  = flatDistribution_->fire();
01362     if( chips[chipIndex]==0 || columns[dColInDet]==0
01363         || rand>pixelEfficiency ) {
01364       // make pixel amplitude =0, pixel will be lost at clusterization
01365       i->second.set(0.); // reset amplitude,
01366     } // end if
01367 
01368   } // end pixel loop
01369 
01370 } // end pixel_indefficiency
01371 
01372 //***********************************************************************
01373 
01374 // Fluctuate the gain and offset for the amplitude calibration
01375 // Use gaussian smearing.
01376 //float SiPixelDigitizerAlgorithm::missCalibrate(const float amp) const {
01377   //float gain  = RandGaussQ::shoot(1.,theGainSmearing);
01378   //float offset  = RandGaussQ::shoot(0.,theOffsetSmearing);
01379   //float newAmp = amp * gain + offset;
01380   // More complex misscalibration
01381 float SiPixelDigitizerAlgorithm::missCalibrate(int col,int row,
01382                                  const float signalInElectrons) const {
01383 
01384   // Central values
01385   //const float p0=0.00352, p1=0.868, p2=112., p3=113.; // pix(0,0,0)
01386   //  const float p0=0.00382, p1=0.886, p2=112.7, p3=113.0; // average roc=0
01387   //const float p0=0.00492, p1=1.998, p2=90.6, p3=134.1; // average roc=6
01388   // Smeared (rms)
01389   //const float s0=0.00020, s1=0.051, s2=5.4, s3=4.4; // average roc=0
01390   //const float s0=0.00015, s1=0.043, s2=3.2, s3=3.1; // col average roc=0
01391 
01392   // Make 2 sets of parameters for Fpix and BPIx:
01393 
01394   float p0=0.0;
01395   float p1=0.0;
01396   float p2=0.0;
01397   float p3=0.0;
01398 
01399   unsigned int Sub_detid=DetId(detID).subdetId();
01400 
01401     if(Sub_detid == PixelSubdetector::PixelBarrel){// barrel layers
01402       p0 = BPix_p0;
01403       p1 = BPix_p1;
01404       p2 = BPix_p2;
01405       p3 = BPix_p3;
01406     } else {// forward disks
01407       p0 = FPix_p0;
01408       p1 = FPix_p1;
01409       p2 = FPix_p2;
01410       p3 = FPix_p3;
01411     }
01412 
01413   //  const float electronsPerVCAL = 65.5; // our present VCAL calibration (feb 2009)
01414   //  const float electronsPerVCAL_Offset = -414.0; // our present VCAL calibration (feb 2009)
01415   float newAmp = 0.; //Modified signal
01416 
01417   // Convert electrons to VCAL units
01418   float signal = (signalInElectrons-electronsPerVCAL_Offset)/electronsPerVCAL;
01419 
01420   // Simulate the analog response with fixed parametrization
01421   newAmp = p3 + p2 * tanh(p0*signal - p1);
01422 
01423 
01424   // Use the pixel-by-pixel calibrations
01425   //transform to ROC index coordinates
01426   //int chipIndex=0, colROC=0, rowROC=0;
01427   //pIndexConverter->transformToROC(col,row,chipIndex,colROC,rowROC);
01428 
01429   // Use calibration from a file
01430   //int chanROC = PixelIndices::pixelToChannelROC(rowROC,colROC); // use ROC coordinates
01431   //float pp0=0, pp1=0,pp2=0,pp3=0;
01432   //map<int,CalParameters,less<int> >::const_iterator it=calmap.find(chanROC);
01433   //CalParameters y  = (*it).second;
01434   //pp0 = y.p0;
01435   //pp1 = y.p1;
01436   //pp2 = y.p2;
01437   //pp3 = y.p3;
01438 
01439   //
01440   // Use random smearing
01441   // Randomize the pixel response
01442   //float pp0  = RandGaussQ::shoot(p0,s0);
01443   //float pp1  = RandGaussQ::shoot(p1,s1);
01444   //float pp2  = RandGaussQ::shoot(p2,s2);
01445   //float pp3  = RandGaussQ::shoot(p3,s3);
01446 
01447   //newAmp = pp3 + pp2 * tanh(pp0*signal - pp1); // Final signal
01448 
01449   //cout<<" misscalibrate "<<col<<" "<<row<<" "<<chipIndex<<" "<<colROC<<" "
01450   //  <<rowROC<<" "<<signalInElectrons<<" "<<signal<<" "<<newAmp<<" "
01451   //  <<(signalInElectrons/theElectronPerADC)<<endl;
01452 
01453   return newAmp;
01454 }
01455 //******************************************************************************
01456 
01457 // Set the drift direction accoring to the Bfield in local det-unit frame
01458 // Works for both barrel and forward pixels.
01459 // Replace the sign convention to fit M.Swartz's formulaes.
01460 // Configurations for barrel and foward pixels possess different tanLorentzAngleperTesla
01461 // parameter value
01462 
01463 //LocalVector SiPixelDigitizerAlgorithm::DriftDirection(float LoAng){
01464 LocalVector SiPixelDigitizerAlgorithm::DriftDirection(){
01465   Frame detFrame(_detp->surface().position(),_detp->surface().rotation());
01466   LocalVector Bfield=detFrame.toLocal(_bfield);
01467   
01468   float alpha2_FPix;
01469   float alpha2_BPix;
01470   float alpha2;
01471   
01472   //float dir_x = -tanLorentzAnglePerTesla * Bfield.y();
01473   //float dir_y = +tanLorentzAnglePerTesla * Bfield.x();
01474   //float dir_z = -1.; // E field always in z direction, so electrons go to -z
01475   // The dir_z has to be +/- 1. !
01476   // LocalVector theDriftDirection = LocalVector(dir_x,dir_y,dir_z);
01477 
01478   float dir_x = 0.0;
01479   float dir_y = 0.0;
01480   float dir_z = 0.0;
01481   float scale = 0.0;
01482 
01483   unsigned int Sub_detid=DetId(detID).subdetId();
01484 
01485   // Read Lorentz angle from cfg file:**************************************************************
01486 
01487   if(!use_LorentzAngle_DB_){
01488     
01489     if( alpha2Order) {
01490       alpha2_FPix = tanLorentzAnglePerTesla_FPix*tanLorentzAnglePerTesla_FPix;
01491       alpha2_BPix = tanLorentzAnglePerTesla_BPix*tanLorentzAnglePerTesla_BPix;
01492     }else {
01493       alpha2_FPix = 0.0;
01494       alpha2_BPix = 0.0;
01495     }
01496     
01497     if(Sub_detid == PixelSubdetector::PixelBarrel){// barrel layers
01498       dir_x = -( tanLorentzAnglePerTesla_BPix * Bfield.y() + alpha2_BPix* Bfield.z()* Bfield.x() );
01499       dir_y = +( tanLorentzAnglePerTesla_BPix * Bfield.x() - alpha2_BPix* Bfield.z()* Bfield.y() );
01500       dir_z = -(1 + alpha2_BPix* Bfield.z()*Bfield.z() );
01501       scale = (1 + alpha2_BPix* Bfield.z()*Bfield.z() );
01502 
01503     } else {// forward disks
01504       dir_x = -( tanLorentzAnglePerTesla_FPix * Bfield.y() + alpha2_FPix* Bfield.z()* Bfield.x() );
01505       dir_y = +( tanLorentzAnglePerTesla_FPix * Bfield.x() - alpha2_FPix* Bfield.z()* Bfield.y() );
01506       dir_z = -(1 + alpha2_FPix* Bfield.z()*Bfield.z() );
01507       scale = (1 + alpha2_FPix* Bfield.z()*Bfield.z() );
01508     }
01509     } // end: Read LA from cfg file.
01510 
01511   //Read Lorentz angle from DB:********************************************************************
01512   if(use_LorentzAngle_DB_){
01513     alpha2 = lorentzAngle * lorentzAngle;
01514     //std::cout << "detID is: "<< it->first <<"The LA per tesla is: "<< it->second << std::endl;
01515     dir_x = -( lorentzAngle * Bfield.y() + alpha2 * Bfield.z()* Bfield.x() );
01516     dir_y = +( lorentzAngle * Bfield.x() - alpha2 * Bfield.z()* Bfield.y() );
01517     dir_z = -(1 + alpha2 * Bfield.z()*Bfield.z() );
01518     scale = (1 + alpha2 * Bfield.z()*Bfield.z() );
01519   }// end: Read LA from DataBase.
01520   
01521   LocalVector theDriftDirection = LocalVector(dir_x/scale, dir_y/scale, dir_z/scale );
01522   
01523 #ifdef TP_DEBUG
01524   LogDebug ("Pixel Digitizer") << " The drift direction in local coordinate is "
01525                                << theDriftDirection ;
01526 #endif
01527   
01528   return theDriftDirection;
01529 }
01530 
01531 //****************************************************************************************************
01532 
01533 void SiPixelDigitizerAlgorithm::pixel_inefficiency_db(void){
01534   if(!use_ineff_from_db_)
01535     return;
01536   
01537   // Loop over hit pixels, amplitude in electrons, channel = coded row,col
01538   for(signal_map_iterator i = _signal.begin();i != _signal.end(); i++) {
01539 
01540     //    int chan = i->first;
01541     pair<int,int> ip = PixelDigi::channelToPixel(i->first);//get pixel pos
01542     int row = ip.first;  // X in row
01543     int col = ip.second; // Y is in col
01544     uint32_t detid = detID;
01545     //transform to ROC index coordinates
01546     if(theSiPixelGainCalibrationService_->isDead(detid, col, row)){
01547       //      std::cout << "now in isdead check, row " << detid << " " << col << "," << row << std::endl;
01548       // make pixel amplitude =0, pixel will be lost at clusterization
01549       i->second.set(0.); // reset amplitude,
01550     } // end if
01551   } // end pixel loop
01552 } // end pixel_indefficiency
01553 
01554 
01555 //****************************************************************************************************
01556 
01557 void SiPixelDigitizerAlgorithm::module_killing_conf(void){
01558   if(!use_module_killing_)
01559     return;
01560   
01561   bool isbad=false;
01562   int detid = detID;
01563   
01564   Parameters::iterator itDeadModules=DeadModules.begin();
01565   
01566   for(; itDeadModules != DeadModules.end(); ++itDeadModules){
01567     int Dead_detID = itDeadModules->getParameter<int>("Dead_detID");
01568     if(detid==Dead_detID){
01569       isbad=true;
01570       break;
01571     }
01572   }
01573   
01574   if(!isbad)
01575     return;
01576   
01577   std::string Module = itDeadModules->getParameter<std::string>("Module");
01578   
01579   if(Module=="whole"){
01580     for(signal_map_iterator i = _signal.begin();i != _signal.end(); i++) {
01581       i->second.set(0.); // reset amplitude
01582     }
01583   }
01584   
01585   for(signal_map_iterator i = _signal.begin();i != _signal.end(); i++) {
01586     pair<int,int> ip = PixelDigi::channelToPixel(i->first);//get pixel pos
01587 
01588     if(Module=="tbmA" && ip.first>=80 && ip.first<=159){
01589       i->second.set(0.);
01590     }
01591 
01592     if( Module=="tbmB" && ip.first<=79){
01593       i->second.set(0.);
01594 
01595     }
01596   }
01597 }
01598 //****************************************************************************************************
01599 void SiPixelDigitizerAlgorithm::module_killing_DB(void){
01600   if(!use_module_killing_)
01601     return;
01602   
01603   bool isbad=false;
01604   uint32_t detid = detID;
01605   
01606   std::vector<SiPixelQuality::disabledModuleType>disabledModules = SiPixelBadModule_->getBadComponentList();
01607   
01608   SiPixelQuality::disabledModuleType badmodule;
01609   
01610   for (size_t id=0;id<disabledModules.size();id++)
01611     {
01612       if(detid==disabledModules[id].DetID){
01613         isbad=true;
01614         badmodule = disabledModules[id];
01615         break;
01616       }
01617     }
01618   
01619   if(!isbad)
01620     return;
01621 
01622   //cout<<"Hit in: "<< detid <<" errorType "<< badmodule.errorType<<" BadRocs="<<hex<<SiPixelBadModule_->getBadRocs(detid)<<dec<<" "<<endl;
01623   if(badmodule.errorType == 0){ // this is a whole dead module.
01624     
01625     for(signal_map_iterator i = _signal.begin();i != _signal.end(); i++) {
01626       i->second.set(0.); // reset amplitude
01627     }
01628   }
01629   else { // all other module types: half-modules and single ROCs.
01630     // Get Bad ROC position:
01631     //follow the example of getBadRocPositions in CondFormats/SiPixelObjects/src/SiPixelQuality.cc
01632     std::vector<GlobalPixel> badrocpositions (0);
01633     for(unsigned int j = 0; j < 16; j++){
01634       if(SiPixelBadModule_->IsRocBad(detid, j) == true){
01635         
01636         std::vector<CablingPathToDetUnit> path = map_.product()->pathToDetUnit(detid);
01637         typedef  std::vector<CablingPathToDetUnit>::const_iterator IT;
01638         for  (IT it = path.begin(); it != path.end(); ++it) {
01639           const PixelROC* myroc = map_.product()->findItem(*it);
01640           if( myroc->idInDetUnit() == j) {
01641             LocalPixel::RocRowCol  local = { 39, 25};   //corresponding to center of ROC row, col
01642             GlobalPixel global = myroc->toGlobal( LocalPixel(local) );
01643             badrocpositions.push_back(global);
01644             break;
01645           }
01646         }
01647       }
01648     }// end of getBadRocPositions
01649     
01650     
01651     for(signal_map_iterator i = _signal.begin();i != _signal.end(); i++) {
01652       pair<int,int> ip = PixelDigi::channelToPixel(i->first);//get pixel pos
01653       
01654       for(std::vector<GlobalPixel>::const_iterator it = badrocpositions.begin(); it != badrocpositions.end(); ++it){
01655         if(it->row >= 80 && ip.first >= 80 ){
01656           if((fabs(ip.second - it->col) < 26) ) {i->second.set(0.);}
01657           else if(it->row==120 && ip.second-it->col==26){i->second.set(0.);}
01658           else if(it->row==119 && it->col-ip.second==26){i->second.set(0.);}
01659         }
01660         else if(it->row < 80 && ip.first < 80 ){
01661           if((fabs(ip.second - it->col) < 26) ){i->second.set(0.);}
01662           else if(it->row==40 && ip.second-it->col==26){i->second.set(0.);}
01663           else if(it->row==39 && it->col-ip.second==26){i->second.set(0.);}
01664        }
01665       }
01666     }
01667   }
01668 }
01669