CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_5_3_13_patch3/src/SimTracker/SiPixelDigitizer/interface/SiPixelDigitizerAlgorithm.h

Go to the documentation of this file.
00001 #ifndef SiPixelDigitizerAlgorithm_h
00002 #define SiPixelDigitizerAlgorithm_h
00003 
00004 
00005 #include <string>
00006 #include <map>
00007 
00008 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00009 #include "DataFormats/SiPixelDigi/interface/PixelDigi.h"
00010 #include "DataFormats/SiPixelDigi/interface/PixelDigiCollection.h"
00011 #include "SimDataFormats/TrackingHit/interface/PSimHitContainer.h"
00012 #include "SimDataFormats/TrackingHit/interface/PSimHit.h"
00013 
00014 //#include "SimGeneral/HepPDT/interface/HepPDTable.h"
00015 //#include "SimGeneral/HepPDTRecord/interface/ParticleDataTable.h"
00016 
00017 #include "SimTracker/Common/interface/SiG4UniversalFluctuation.h"
00018 #include "Geometry/CommonDetUnit/interface/GeomDetUnit.h"
00019 #include "DataFormats/DetId/interface/DetId.h"
00020 #include "Geometry/TrackerGeometryBuilder/interface/PixelGeomDetUnit.h"
00021 #include "Geometry/CommonTopologies/interface/PixelTopology.h"
00022 #include "Geometry/CommonDetUnit/interface/GeomDetType.h"
00023 //#include "Geometry/TrackerGeometryBuilder/interface/TrackerGeomFromDetUnits.h"
00024 #include "SimGeneral/NoiseGenerators/interface/GaussianTailNoiseGenerator.h"
00025 //#include "DataFormats/GeometrySurface/interface/TkRotation.h"
00026 //#include "DataFormats/GeometrySurface/interface/GloballyPositioned.h"
00027 #include "SimDataFormats/TrackerDigiSimLink/interface/PixelDigiSimLink.h"
00028 //#include "SimDataFormats/TrackerDigiSimLink/interface/PixelDigiSimLinkCollection.h"
00029 #include "DataFormats/Common/interface/DetSetVector.h"
00030 
00031 #include "CondFormats/SiPixelObjects/interface/PixelIndices.h"
00032 
00033 // pixel gain payload access (offline version for Simulation)
00034 #include "CalibTracker/SiPixelESProducers/interface/SiPixelGainCalibrationOfflineSimService.h"
00035 
00036 // Accessing Pixel Lorentz Angle from the DB:
00037 #include "CondFormats/SiPixelObjects/interface/SiPixelLorentzAngle.h"
00038 #include "CondFormats/DataRecord/interface/SiPixelLorentzAngleSimRcd.h"
00039 
00040 // Accessing Pixel dead modules from the DB:
00041 #include "CondFormats/SiPixelObjects/interface/SiPixelQuality.h"
00042 #include "CondFormats/DataRecord/interface/SiPixelQualityRcd.h"
00043 
00044 #include "CondFormats/SiPixelObjects/interface/DetectorIndex.h"
00045 #include "CondFormats/SiPixelObjects/interface/SiPixelFedCabling.h"
00046 
00047 #include "CondFormats/DataRecord/interface/SiPixelFedCablingMapRcd.h"
00048 #include "Geometry/Records/interface/TrackerDigiGeometryRecord.h"
00049 #include "Geometry/TrackerGeometryBuilder/interface/TrackerGeometry.h"
00050 
00051 
00052 // For the random numbers
00053 namespace CLHEP {
00054   class HepRandomEngine;
00055   class RandGaussQ;
00056   class RandFlat;
00057 }
00058 
00059 class SiPixelDigitizerAlgorithm  {
00060  public:
00061   
00062   SiPixelDigitizerAlgorithm(const edm::ParameterSet& conf, CLHEP::HepRandomEngine&);
00063   ~SiPixelDigitizerAlgorithm();
00064   
00065   //run the algorithm to digitize a single det
00066   edm::DetSet<PixelDigi>::collection_type  
00067     run(const std::vector<PSimHit> &input,PixelGeomDetUnit *pixdet,GlobalVector);
00068 
00069    //
00070   std::vector<PixelDigiSimLink> make_link() {
00071     return link_coll; }
00072   void init(const edm::EventSetup& es);
00073   void fillDeadModules(const edm::EventSetup& es);
00074   void fillLorentzAngle(const edm::EventSetup& es);
00075   void fillMapandGeom(const edm::EventSetup& es);
00076 
00077 
00078  private:
00079   
00080   //Accessing Lorentz angle from DB:
00081   edm::ESHandle<SiPixelLorentzAngle> SiPixelLorentzAngle_;
00082 
00083   //Accessing Dead pixel modules from DB:
00084   edm::ESHandle<SiPixelQuality> SiPixelBadModule_;
00085 
00086   //Accessing Map and Geom:
00087   edm::ESHandle<SiPixelFedCablingMap> map_;
00088   edm::ESHandle<TrackerGeometry> geom_;
00089 
00090   typedef std::vector<edm::ParameterSet> Parameters;
00091   Parameters DeadModules;
00092 
00093   // Define internal classes
00094   
00095   // Define a class to hold the calibration parameters per pixel
00096   // Internal
00097   class CalParameters {
00098   public:
00099     float p0;
00100     float p1;
00101     float p2;
00102     float p3;
00103   };
00104   //
00105   // Define a class for 3D ionization points and energy
00106   //
00110   class EnergyDepositUnit{
00111   public:
00112     EnergyDepositUnit(): _energy(0),_position(0,0,0){}
00113     EnergyDepositUnit(float energy,float x, float y, float z):
00114     _energy(energy),_position(x,y,z){}
00115     EnergyDepositUnit(float energy, Local3DPoint position):
00116     _energy(energy),_position(position){}
00117     float x() const{return _position.x();}
00118     float y() const{return _position.y();}
00119     float z() const{return _position.z();}
00120     float energy() const { return _energy;}
00121   private:
00122     float _energy;
00123     Local3DPoint _position;
00124   };
00125 
00126   //
00127   // define class to store signals on the collection surface
00128   //
00133   class SignalPoint {
00134   public:
00135     SignalPoint() : _pos(0,0), _time(0), _amplitude(0), 
00136       _sigma_x(1.), _sigma_y(1.), _hitp(0) {}
00137     
00138     SignalPoint( float x, float y, float sigma_x, float sigma_y,
00139                  float t, float a=1.0) :
00140     _pos(x,y), _time(t), _amplitude(a), _sigma_x(sigma_x), 
00141       _sigma_y(sigma_y), _hitp(0) {}
00142     
00143     SignalPoint( float x, float y, float sigma_x, float sigma_y,
00144                  float t, const PSimHit& hit, float a=1.0) :
00145     _pos(x,y), _time(t), _amplitude(a), _sigma_x(sigma_x), 
00146       _sigma_y(sigma_y),_hitp(&hit) {}
00147     
00148     const LocalPoint& position() const { return _pos;}
00149     float x()         const { return _pos.x();}
00150     float y()         const { return _pos.y();}
00151     float sigma_x()   const { return _sigma_x;}
00152     float sigma_y()   const { return _sigma_y;}
00153     float time()      const { return _time;}
00154     float amplitude() const { return _amplitude;}
00155     const PSimHit& hit()           { return *_hitp;}
00156     SignalPoint& set_amplitude( float amp) { _amplitude = amp; return *this;}
00157     
00158   
00159 
00160   private:
00161     LocalPoint         _pos;
00162     float              _time;
00163     float              _amplitude;
00164     float              _sigma_x;   // gaussian sigma in the x direction (cm)
00165     float              _sigma_y;   //    "       "          y direction (cm) */
00166     const PSimHit*   _hitp;
00167   };
00168  
00169   //
00170   // definition class
00171   //
00176   class Amplitude {
00177   public:
00178     Amplitude() : _amp(0.0) { _hits.reserve(1);}
00179     Amplitude( float amp, const PSimHit* hitp, float frac) :
00180       _amp(amp), _hits(1, hitp), _frac(1,frac) {
00181 
00182     //in case of digi from noisypixels
00183       //the MC information are removed 
00184       if (_frac[0]<-0.5) {
00185         _frac.pop_back();
00186         _hits.pop_back();
00187      }
00188 
00189     }
00190 
00191     // can be used as a float by convers.
00192     operator float() const { return _amp;}
00193     float ampl() const {return _amp;}
00194     std::vector<float> individualampl() const {return _frac;}
00195     const std::vector<const PSimHit*>& hits() { return _hits;}
00196 
00197     void operator+=( const Amplitude& other) {
00198       _amp += other._amp;
00199       //in case of contribution of noise to the digi
00200       //the MC information are removed 
00201       if (other._frac[0]>-0.5){
00202         _hits.insert( _hits.end(), other._hits.begin(), other._hits.end());
00203         _frac.insert(_frac.end(), other._frac.begin(), other._frac.end());
00204       }
00205    }
00206 
00207     void operator+=( const float& amp) {
00208       _amp += amp;
00209     }
00210    
00211     void set (const float amplitude) {  // Used to reset the amplitude
00212       _amp = amplitude;
00213     }
00214 /*     void setind (const float indamplitude) {  // Used to reset the amplitude */
00215 /*       _frac = idamplitude; */
00216 /*     } */
00217   private:
00218     float _amp;
00219     std::vector<const PSimHit*> _hits;
00220     std::vector<float> _frac;
00221   };  // end class Amplitude
00222 
00223 
00224  private:
00225 
00226     // Internal typedefs
00227     typedef std::map< int, Amplitude, std::less<int> >   signal_map_type;  // from
00228     typedef signal_map_type::iterator          signal_map_iterator; //Digi.Skel.  
00229     typedef std::map<unsigned int, std::vector<float>,std::less<unsigned int> > 
00230       simlink_map;
00231     typedef GloballyPositioned<double>      Frame;
00232 
00233     // Variables 
00234     edm::ParameterSet conf_;
00235     //external parameters 
00236     //-- primary ionization
00237     int    NumberOfSegments; // =20 does not work ;
00238     // go from Geant energy GeV to number of electrons
00239     float GeVperElectron; // 3.7E-09 
00240     
00241     //-- drift
00242     float Sigma0; //=0.0007  // Charge diffusion in microns for 300 micron Si
00243     float Dist300;  //=0.0300  // Define 300microns for normalization 
00244     bool alpha2Order;          // Switch on/off of E.B effect 
00245 
00246 
00247  
00248     //-- induce_signal
00249     float ClusterWidth;       // Gaussian charge cutoff width in sigma units
00250     //-- make_digis 
00251     float theElectronPerADC;     // Gain, number of electrons per adc count.
00252     int theAdcFullScale;         // Saturation count, 255=8bit.
00253     float theNoiseInElectrons;   // Noise (RMS) in units of electrons.
00254     float theReadoutNoise;       // Noise of the readount chain in elec,
00255                                  //inludes DCOL-Amp,TBM-Amp, Alt, AOH,OptRec.
00256 
00257     float theSmearedChargeRMS;
00258 
00259     float thePixelThreshold;     // Pixel threshold in units of noise.
00260 
00261     float thePixelThresholdInE;  // Pixel noise in electrons.
00262 
00263     float theThresholdInE_FPix;  // Pixel threshold in electrons FPix.
00264     float theThresholdInE_BPix;  // Pixel threshold in electrons BPix.
00265 
00266     double theThresholdSmearing_FPix;
00267     double theThresholdSmearing_BPix;
00268 
00269     double electronsPerVCAL;          // for electrons - VCAL conversion
00270     double electronsPerVCAL_Offset;   // in misscalibrate()
00271 
00272     float theTofLowerCut;             // Cut on the particle TOF
00273     float theTofUpperCut;             // Cut on the particle TOF
00274     float tanLorentzAnglePerTesla_FPix;   //FPix Lorentz angle tangent per Tesla
00275     float tanLorentzAnglePerTesla_BPix;   //BPix Lorentz angle tangent per Tesla
00276     float lorentzAngle;    
00277 
00278     float FPix_p0;
00279     float FPix_p1;
00280     float FPix_p2;
00281     float FPix_p3;
00282     float BPix_p0;
00283     float BPix_p1;
00284     float BPix_p2;
00285     float BPix_p3;
00286 
00287 
00288     //-- add_noise
00289     bool addNoise;
00290     bool addChargeVCALSmearing;
00291     bool addNoisyPixels;
00292     bool fluctuateCharge;
00293     bool addPixelInefficiency;
00294     //-- pixel efficiency
00295     bool pixelInefficiency;      // Switch on pixel ineffciency
00296     int  thePixelLuminosity;        // luminosity for inefficiency, 0,1,10
00297 
00298     bool addThresholdSmearing;
00299         
00300     int theColsInChip;           // num of columns per ROC (for pix ineff.)
00301     int theRowsInChip;           // num of rows per ROC
00302     
00303     int numColumns; // number of pixel columns in a module (detUnit)
00304     int numRows;    // number          rows
00305     float moduleThickness; // sensor thickness 
00306     //  int digis; 
00307     const PixelGeomDetUnit* _detp;
00308     uint32_t detID;     // Det id
00309     
00310 
00311     std::vector<PSimHit> _PixelHits; //cache
00312     const PixelTopology* topol;
00313     
00314     std::vector<PixelDigi> internal_coll; //empty vector of PixelDigi used in digitize
00315 
00316     std::vector<PixelDigiSimLink> link_coll;
00317     GlobalVector _bfield;
00318     
00319     float PixelEff;
00320     float PixelColEff;
00321     float PixelChipEff;
00322     float PixelEfficiency;
00323     float PixelColEfficiency;
00324     float PixelChipEfficiency;
00325     float thePixelEfficiency[6];     // Single pixel effciency
00326     float thePixelColEfficiency[6];  // Column effciency
00327     float thePixelChipEfficiency[6]; // ROC efficiency
00328     
00329     //-- calibration smearing
00330     bool doMissCalibrate;         // Switch on the calibration smearing
00331     float theGainSmearing;        // The sigma of the gain fluctuation (around 1)
00332     float theOffsetSmearing;      // The sigma of the offset fluct. (around 0)
00333     
00334 
00335     // The PDTable
00336     //HepPDTable *particleTable;
00337     //ParticleDataTable *particleTable;
00338 
00339     //-- charge fluctuation
00340     double tMax;  // The delta production cut, should be as in OSCAR = 30keV
00341     //                                           cmsim = 100keV
00342     // The eloss fluctuation class from G4. Is the right place? 
00343     SiG4UniversalFluctuation * fluctuate;   // make a pointer 
00344     GaussianTailNoiseGenerator * theNoiser; //
00345 
00346 
00347 
00348 
00349     PixelIndices * pIndexConverter;         // Pointer to the index converter 
00350    
00351     std::vector<EnergyDepositUnit> _ionization_points;
00352     std::vector<SignalPoint> _collection_points;
00353     
00354     simlink_map simi;
00355     signal_map_type     _signal;       // from Digi.Skel.
00356 
00357     // To store calibration constants
00358     std::map<int,CalParameters,std::less<int> > calmap;
00359 
00360 
00361     //-- additional member functions    
00362     // Private methods
00363     void primary_ionization( const PSimHit& hit);
00364     std::vector<PixelDigi> digitize(PixelGeomDetUnit *det);
00365     void drift(const PSimHit& hit);
00366     void induce_signal( const PSimHit& hit);
00367     void fluctuateEloss(int particleId, float momentum, float eloss, 
00368                         float length, int NumberOfSegments,
00369                         float elossVector[]);
00370     void add_noise();
00371     void make_digis();
00372     void pixel_inefficiency();
00373     bool use_ineff_from_db_;
00374 
00375     bool use_module_killing_; // remove or not the dead pixel modules
00376     bool use_deadmodule_DB_; // if we want to get dead pixel modules from the DataBase.
00377     bool use_LorentzAngle_DB_; // if we want to get Lorentz angle from the DataBase.
00378 
00379     void pixel_inefficiency_db(); 
00380        // access to the gain calibration payloads in the db. Only gets initialized if check_dead_pixels_ is set to true.
00381     SiPixelGainCalibrationOfflineSimService * theSiPixelGainCalibrationService_;    
00382     float missCalibrate(int col, int row, float amp) const;  
00383     LocalVector DriftDirection();
00384 
00385     void module_killing_conf(); // remove dead modules using the list in the configuration file PixelDigi_cfi.py
00386     void module_killing_DB();  // remove dead modules uisng the list in the DB
00387 
00388    // For random numbers
00389     CLHEP::HepRandomEngine& rndEngine;
00390 
00391     CLHEP::RandFlat *flatDistribution_;
00392     CLHEP::RandGaussQ *gaussDistribution_;
00393     CLHEP::RandGaussQ *gaussDistributionVCALNoise_;
00394 
00395     
00396     // Threshold gaussian smearing:
00397     CLHEP::RandGaussQ *smearedThreshold_FPix_;
00398     CLHEP::RandGaussQ *smearedThreshold_BPix_;
00399     
00400     CLHEP::RandGaussQ *smearedChargeDistribution_ ;
00401     
00402     // the random generator
00403     CLHEP::RandGaussQ* theGaussianDistribution;
00404 
00405 
00406 };
00407 
00408 #endif