Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:11:19

0001 #ifndef HDShower_H
0002 #define HDShower_H
0003 
0004 //FastSimulation Headers
0005 #include "FastSimulation/ShowerDevelopment/interface/HDShowerParametrization.h"
0006 
0007 #include "DataFormats/Math/interface/Vector3D.h"
0008 
0009 #include <vector>
0010 
0011 /** 
0012  * \author Salavat Abdullin
0013  * \date: 21-Oct-2004
0014  * \parameterized hadronic shower simulation 
0015  * \based on paper G.Gridhammer, M.Rudowicz, S.Peters, SLAC-PUB-5072
0016  */
0017 
0018 class EcalHitMaker;
0019 class HcalHitMaker;
0020 class RandomEngineAndDistribution;
0021 
0022 class HDShower {
0023 public:
0024   typedef math::XYZVector XYZPoint;
0025 
0026   typedef std::pair<XYZPoint, double> Spot;
0027   typedef std::pair<unsigned int, double> Step;
0028   typedef std::vector<Step> Steps;
0029   typedef Steps::const_iterator step_iterator;
0030 
0031   HDShower(const RandomEngineAndDistribution* engine,
0032            HDShowerParametrization* myParam,
0033            EcalHitMaker* myGrid,
0034            HcalHitMaker* myHcalHitMaker,
0035            int onECAL,
0036            double epart,
0037            double pmip);
0038 
0039   int getmip() { return mip; }
0040 
0041   virtual ~HDShower() { ; }
0042 
0043   /// Compute the shower longitudinal and lateral development
0044   bool compute();
0045 
0046 private:
0047   // The longitudinal development ersatzt.
0048   double gam(double x, double a) const { return pow(x, a - 1.) * exp(-x); }
0049 
0050   // Transverse integrated probability function (for 4R max size)
0051   // integral of the transverse ansatz f(r) =  2rR/(r**2 + R**2)**2 ->
0052   // 1/R - R/(R**2+r**2) | at limits 4R - 0
0053   double transProb(double factor, double R, double r) {
0054     double fsq = factor * factor;
0055     return ((fsq + 1.) / fsq) * r * r / (r * r + R * R);
0056   }
0057   // Compute longE[i] and transR[i] for all nsteps
0058   void makeSteps(int nsteps);
0059 
0060   int indexFinder(double x, const std::vector<double>& Fhist);
0061 
0062   // The parametrization
0063   HDShowerParametrization* theParam;
0064 
0065   // The Calorimeter properties
0066   const ECALProperties* theECALproperties;
0067   const HCALProperties* theHCALproperties;
0068 
0069   // Basic parameters of the simulation
0070   double theR1, theR2, theR3;
0071   double alpEM, betEM, alpHD, betHD, part, tgamEM, tgamHD;
0072 
0073   // The basic quantities for the shower development.
0074   double lambdaEM, lambdaHD, x0EM, x0HD;
0075   double depthStart;
0076   double aloge;
0077 
0078   std::vector<int> detector, nspots;
0079   std::vector<double> eStep, rlamStep;
0080   std::vector<double> x0curr, x0depth;
0081   std::vector<double> lamstep, lamcurr, lamdepth, lamtotal;
0082 
0083   int infinity;  // big number of cycles if exit is on a special condition
0084 
0085   // The crystal grid
0086   EcalHitMaker* theGrid;
0087 
0088   // The HCAL
0089   HcalHitMaker* theHcalHitMaker;
0090 
0091   // OnECAL flag as an input parameter ...
0092   int onEcal;
0093 
0094   // MIP in ECAL map flag
0095   int mip;
0096 
0097   // Input energy to distribute
0098   double e;
0099 
0100   // HCAL losses option (0-off, 1-on)
0101   int lossesOpt;
0102   // Number of longitudinal steps in HCAL
0103   int nDepthSteps;
0104   // Number of bins in the transverse probability histogram
0105   int nTRsteps;
0106   // Transverse size tuning factor
0107   double transParam;
0108   // Transverse normalization : 1 for HB/HE, 0.5 for HF (narrow showers)
0109   double transFactor;
0110   // HCAL energy spot size
0111   double eSpotSize;
0112   // Longitudinal step size (lambda units)
0113   double depthStep;
0114   // Energy threshold (one depth step -> nDepthSteps);
0115   double criticalEnergy;
0116   // Transverse size cut (in character transverse size units)
0117   double maxTRfactor;
0118   // Balance between ECAL and HCAL "visible" energy (default = 1.)
0119   double balanceEH;
0120   // Regulator of HCAL depth of the shower (to adjust/shrink it to CMS depth)
0121   double hcalDepthFactor;
0122 
0123   // Famos Random Engine
0124   const RandomEngineAndDistribution* random;
0125 
0126   //calorimeter depths
0127   double depthECAL, depthGAP, depthGAPx0, depthHCAL, depthToHCAL;
0128 };
0129 
0130 #endif