CMS 3D CMS Logo

SiHitDigitizer.cc

Go to the documentation of this file.
00001 #include "SimTracker/SiStripDigitizer/interface/SiHitDigitizer.h"
00002 #include "SimDataFormats/TrackingHit/interface/PSimHit.h"
00003 #include "SimTracker/SiStripDigitizer/interface/SiLinearChargeCollectionDrifter.h"
00004 #include "SimTracker/SiStripDigitizer/interface/SiLinearChargeDivider.h"
00005 
00006 #include "DataFormats/GeometryVector/interface/LocalPoint.h"
00007 #include "SimTracker/SiStripDigitizer/interface/SiTrivialInduceChargeOnStrips.h"
00008 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00009 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00010 
00011 #define CBOLTZ (1.38E-23)
00012 #define e_SI (1.6E-19)
00013 
00014 SiHitDigitizer::SiHitDigitizer(const edm::ParameterSet& conf,CLHEP::HepRandomEngine& eng ):conf_(conf),rndEngine(eng){
00015   
00016   //
00017   // Construct default classes
00018   //
00019   depletionVoltage       = conf_.getParameter<double>("DepletionVoltage");
00020   appliedVoltage         = conf_.getParameter<double>("AppliedVoltage");
00021   chargeMobility         = conf_.getParameter<double>("ChargeMobility");
00022   temperature            = conf_.getParameter<double>("Temperature");
00023   gevperelectron         =conf_.getParameter<double>("GevPerElectron");
00024   chargeDistributionRMS  =conf_.getParameter<double>("ChargeDistributionRMS");
00025   noDiffusion            = conf_.getParameter<bool>("noDiffusion");
00026   double diffusionConstant = CBOLTZ/e_SI * chargeMobility * temperature;
00027   if (noDiffusion) diffusionConstant *= 1.0e-3;
00028   
00029   theSiChargeDivider = new SiLinearChargeDivider(conf_,rndEngine);
00030   
00031   theSiChargeCollectionDrifter = 
00032     new SiLinearChargeCollectionDrifter(diffusionConstant,
00033                                         chargeDistributionRMS,
00034                                         depletionVoltage,
00035                                         appliedVoltage);
00036 
00037   theSiInduceChargeOnStrips = new SiTrivialInduceChargeOnStrips(conf,gevperelectron);
00038 }
00039 
00040 
00041 SiHitDigitizer::~SiHitDigitizer(){
00042   delete theSiChargeDivider;
00043   delete theSiChargeCollectionDrifter;
00044   delete theSiInduceChargeOnStrips;
00045 }
00046 
00047 
00048 void 
00049 SiHitDigitizer::processHit(const PSimHit& hit, const StripGeomDetUnit& det, GlobalVector bfield,float langle,
00050                            std::vector<double>& locAmpl, unsigned int& firstChannelWithSignal, unsigned int& lastChannelWithSignal){
00051   
00052   //
00053   // Compute the drift direction for this det
00054   //
00055   
00056   double moduleThickness = det.specificSurface().bounds().thickness(); // full detector thicness
00057   double timeNormalisation = (moduleThickness*moduleThickness)/(2.*depletionVoltage*chargeMobility);
00058   
00059   LocalVector driftDir = DriftDirection(&det,bfield,langle);
00060   
00061   //
00062   // Fully process one SimHit
00063   //
00064   
00065   SiChargeCollectionDrifter::ionization_type ion = theSiChargeDivider->divide(hit, driftDir, moduleThickness, det);
00066   
00067   theSiInduceChargeOnStrips->induce(theSiChargeCollectionDrifter->drift(ion,driftDir,moduleThickness,timeNormalisation),det,
00068                                     locAmpl,firstChannelWithSignal,lastChannelWithSignal);
00069 }
00070 
00071 LocalVector SiHitDigitizer::DriftDirection(const StripGeomDetUnit* _detp,GlobalVector _bfield,float langle){
00072   // taken from ORCA/Tracker/SiStripDet/src/SiStripDet.cc
00073   Frame detFrame(_detp->surface().position(),_detp->surface().rotation());
00074   LocalVector Bfield=detFrame.toLocal(_bfield);
00075   
00076   float dir_x = -langle * Bfield.y();
00077   float dir_y = +langle * Bfield.x();
00078   float dir_z = 1.; // E field always in z direction
00079   LocalVector theDriftDirection = LocalVector(dir_x,dir_y,dir_z);
00080   if ( conf_.getUntrackedParameter<int>("VerbosityLevel") > 0 ) {
00081     edm::LogInfo("StripDigiInfo")<< " The drift direction in local coordinate is "<<theDriftDirection;
00082     if(langle==0.)
00083       edm::LogWarning("StripDigiInfo")<< "ERROR: Lorentz angle = 0 for module "<<_detp->geographicalId().rawId();
00084   }
00085   return theDriftDirection;
00086   
00087 }
00088 void SiHitDigitizer::setParticleDataTable(const ParticleDataTable * pdt)
00089 {
00090   theSiChargeDivider->setParticleDataTable(pdt);
00091 }

Generated on Tue Jun 9 17:47:53 2009 for CMSSW by  doxygen 1.5.4