CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_1_8_patch13/src/FastSimulation/MaterialEffects/src/MaterialEffects.cc

Go to the documentation of this file.
00001 //Framework Headers
00002 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00003 
00004 //TrackingTools Headers
00005 
00006 // Famos Headers
00007 #include "FastSimulation/Event/interface/FSimEvent.h"
00008 #include "FastSimulation/ParticlePropagator/interface/ParticlePropagator.h"
00009 #include "FastSimulation/TrackerSetup/interface/TrackerLayer.h"
00010 #include "FastSimulation/MaterialEffects/interface/MaterialEffects.h"
00011 
00012 #include "FastSimulation/MaterialEffects/interface/PairProductionSimulator.h"
00013 #include "FastSimulation/MaterialEffects/interface/MultipleScatteringSimulator.h"
00014 #include "FastSimulation/MaterialEffects/interface/BremsstrahlungSimulator.h"
00015 #include "FastSimulation/MaterialEffects/interface/EnergyLossSimulator.h"
00016 #include "FastSimulation/MaterialEffects/interface/NuclearInteractionSimulator.h"
00017 
00018 #include <list>
00019 #include <map>
00020 #include <string>
00021 
00022 MaterialEffects::MaterialEffects(const edm::ParameterSet& matEff,
00023                                  const RandomEngine* engine)
00024   : PairProduction(0), Bremsstrahlung(0), 
00025     MultipleScattering(0), EnergyLoss(0), 
00026     NuclearInteraction(0),
00027     pTmin(999.), random(engine)
00028 {
00029   // Set the minimal photon energy for a Brem from e+/-
00030 
00031   bool doPairProduction     = matEff.getParameter<bool>("PairProduction");
00032   bool doBremsstrahlung     = matEff.getParameter<bool>("Bremsstrahlung");
00033   bool doEnergyLoss         = matEff.getParameter<bool>("EnergyLoss");
00034   bool doMultipleScattering = matEff.getParameter<bool>("MultipleScattering");
00035   bool doNuclearInteraction = matEff.getParameter<bool>("NuclearInteraction");
00036 
00037   double A = matEff.getParameter<double>("A");
00038   double Z = matEff.getParameter<double>("Z");
00039   double density = matEff.getParameter<double>("Density");
00040   double radLen = matEff.getParameter<double>("RadiationLength");
00041   
00042   // Set the minimal pT before giving up the dE/dx treatment
00043 
00044   if ( doPairProduction ) { 
00045 
00046     double photonEnergy = matEff.getParameter<double>("photonEnergy");
00047     PairProduction = new PairProductionSimulator(photonEnergy,
00048                                                  random);
00049 
00050   }
00051 
00052   if ( doBremsstrahlung ) { 
00053 
00054     double bremEnergy = matEff.getParameter<double>("bremEnergy");
00055     double bremEnergyFraction = matEff.getParameter<double>("bremEnergyFraction");
00056     Bremsstrahlung = new BremsstrahlungSimulator(bremEnergy,
00057                                                  bremEnergyFraction,
00058                                                  random);
00059 
00060   }
00061 
00062   if ( doEnergyLoss ) { 
00063 
00064     pTmin = matEff.getParameter<double>("pTmin");
00065     EnergyLoss = new EnergyLossSimulator(random,A,Z,density,radLen);
00066 
00067   }
00068 
00069   if ( doMultipleScattering ) { 
00070     
00071     MultipleScattering = new MultipleScatteringSimulator(random,A,Z,density,radLen);
00072 
00073   }
00074 
00075   if ( doNuclearInteraction ) { 
00076 
00077     // The energies simulated
00078     std::vector<double> pionEnergies 
00079       = matEff.getUntrackedParameter<std::vector<double> >("pionEnergies");
00080 
00081     // The particle types simulated
00082     std::vector<int> pionTypes 
00083       = matEff.getUntrackedParameter<std::vector<int> >("pionTypes");
00084 
00085     // The corresponding particle names
00086     std::vector<std::string> pionNames 
00087       = matEff.getUntrackedParameter<std::vector<std::string> >("pionNames");
00088 
00089     // The corresponding particle masses
00090     std::vector<double> pionMasses 
00091       = matEff.getUntrackedParameter<std::vector<double> >("pionMasses");
00092 
00093     // The smallest momentum for inelastic interactions
00094     std::vector<double> pionPMin 
00095       = matEff.getUntrackedParameter<std::vector<double> >("pionMinP");
00096 
00097     // The interaction length / radiation length ratio for each particle type
00098     std::vector<double> lengthRatio 
00099       = matEff.getParameter<std::vector<double> >("lengthRatio");
00100     //    std::map<int,double> lengthRatio;
00101     //    for ( unsigned i=0; i<theLengthRatio.size(); ++i )
00102     //      lengthRatio[ pionTypes[i] ] = theLengthRatio[i];
00103 
00104     // A global fudge factor for TEC layers (which apparently do not react to 
00105     // hadrons the same way as all other layers...
00106     theTECFudgeFactor = matEff.getParameter<double>("fudgeFactor");
00107 
00108     // The evolution of the interaction lengths with energy
00109     std::vector<double> theRatios  
00110       = matEff.getUntrackedParameter<std::vector<double> >("ratios");
00111     //std::map<int,std::vector<double> > ratios;
00112     //for ( unsigned i=0; i<pionTypes.size(); ++i ) { 
00113     //  for ( unsigned j=0; j<pionEnergies.size(); ++j ) { 
00114     //  ratios[ pionTypes[i] ].push_back(theRatios[ i*pionEnergies.size() + j ]);
00115     //  }
00116     //}
00117     std::vector< std::vector<double> > ratios;
00118     ratios.resize(pionTypes.size());
00119     for ( unsigned i=0; i<pionTypes.size(); ++i ) { 
00120       for ( unsigned j=0; j<pionEnergies.size(); ++j ) { 
00121         ratios[i].push_back(theRatios[ i*pionEnergies.size() + j ]);
00122       }
00123     }
00124 
00125     // The smallest momentum for elastic interactions
00126     double pionEnergy 
00127       = matEff.getParameter<double>("pionEnergy");
00128 
00129     // The algorithm to compute the distance between primary and secondaries
00130     // when a nuclear interaction occurs
00131     unsigned distAlgo 
00132       = matEff.getParameter<unsigned>("distAlgo");
00133     double distCut 
00134       = matEff.getParameter<double>("distCut");
00135 
00136     // The file to read the starting interaction in each files
00137     // (random reproducibility in case of a crash)
00138     std::string inputFile 
00139       = matEff.getUntrackedParameter<std::string>("inputFile");
00140 
00141     // Build the ID map (i.e., what is to be considered as a proton, etc...)
00142     std::map<int,int> idMap;
00143     // Protons
00144     std::vector<int> idProtons 
00145       = matEff.getUntrackedParameter<std::vector<int> >("protons");
00146     for ( unsigned i=0; i<idProtons.size(); ++i ) 
00147       idMap[idProtons[i]] = 2212;
00148     // Anti-Protons
00149     std::vector<int> idAntiProtons 
00150       = matEff.getUntrackedParameter<std::vector<int> >("antiprotons");
00151     for ( unsigned i=0; i<idAntiProtons.size(); ++i ) 
00152       idMap[idAntiProtons[i]] = -2212;
00153     // Neutrons
00154     std::vector<int> idNeutrons 
00155       = matEff.getUntrackedParameter<std::vector<int> >("neutrons");
00156     for ( unsigned i=0; i<idNeutrons.size(); ++i ) 
00157       idMap[idNeutrons[i]] = 2112;
00158     // Anti-Neutrons
00159     std::vector<int> idAntiNeutrons 
00160       = matEff.getUntrackedParameter<std::vector<int> >("antineutrons");
00161     for ( unsigned i=0; i<idAntiNeutrons.size(); ++i ) 
00162       idMap[idAntiNeutrons[i]] = -2112;
00163     // K0L's
00164     std::vector<int> idK0Ls 
00165       = matEff.getUntrackedParameter<std::vector<int> >("K0Ls");
00166     for ( unsigned i=0; i<idK0Ls.size(); ++i ) 
00167       idMap[idK0Ls[i]] = 130;
00168     // K+'s
00169     std::vector<int> idKplusses 
00170       = matEff.getUntrackedParameter<std::vector<int> >("Kplusses");
00171     for ( unsigned i=0; i<idKplusses.size(); ++i ) 
00172       idMap[idKplusses[i]] = 321;
00173     // K-'s
00174     std::vector<int> idKminusses 
00175       = matEff.getUntrackedParameter<std::vector<int> >("Kminusses");
00176     for ( unsigned i=0; i<idKminusses.size(); ++i ) 
00177       idMap[idKminusses[i]] = -321;
00178     // pi+'s
00179     std::vector<int> idPiplusses 
00180       = matEff.getUntrackedParameter<std::vector<int> >("Piplusses");
00181     for ( unsigned i=0; i<idPiplusses.size(); ++i ) 
00182       idMap[idPiplusses[i]] = 211;
00183     // pi-'s
00184     std::vector<int> idPiminusses 
00185       = matEff.getUntrackedParameter<std::vector<int> >("Piminusses");
00186     for ( unsigned i=0; i<idPiminusses.size(); ++i ) 
00187       idMap[idPiminusses[i]] = -211;
00188 
00189     // Construction
00190     NuclearInteraction = 
00191       new NuclearInteractionSimulator(pionEnergies, pionTypes, pionNames, 
00192                                       pionMasses, pionPMin, pionEnergy, 
00193                                       lengthRatio, ratios, idMap, 
00194                                       inputFile, distAlgo, distCut, random);
00195   }
00196 
00197 }
00198 
00199 
00200 MaterialEffects::~MaterialEffects() {
00201 
00202   if ( PairProduction ) delete PairProduction;
00203   if ( Bremsstrahlung ) delete Bremsstrahlung;
00204   if ( EnergyLoss ) delete EnergyLoss;
00205   if ( MultipleScattering ) delete MultipleScattering;
00206   if ( NuclearInteraction ) delete NuclearInteraction;
00207 
00208 }
00209 
00210 void MaterialEffects::interact(FSimEvent& mySimEvent, 
00211                                const TrackerLayer& layer,
00212                                ParticlePropagator& myTrack,
00213                                unsigned itrack) {
00214 
00215   MaterialEffectsSimulator::RHEP_const_iter DaughterIter;
00216   double radlen;
00217   theEnergyLoss = 0;
00218   theNormalVector = normalVector(layer,myTrack);
00219   radlen = radLengths(layer,myTrack);
00220 
00221 //-------------------
00222 //  Photon Conversion
00223 //-------------------
00224 
00225   if ( PairProduction && myTrack.pid()==22 ) {
00226     
00227     //
00228     PairProduction->updateState(myTrack,radlen);
00229 
00230     if ( PairProduction->nDaughters() ) {       
00231       //add a vertex to the mother particle
00232       int ivertex = mySimEvent.addSimVertex(myTrack.vertex(),itrack,
00233                                             FSimVertexType::PAIR_VERTEX);
00234       
00235       // This was a photon that converted
00236       for ( DaughterIter = PairProduction->beginDaughters();
00237             DaughterIter != PairProduction->endDaughters(); 
00238             ++DaughterIter) {
00239 
00240         mySimEvent.addSimTrack(&(*DaughterIter), ivertex);
00241 
00242       }
00243       // The photon converted. Return.
00244       return;
00245     }
00246   }
00247 
00248   if ( myTrack.pid() == 22 ) return;
00249 
00250 //------------------------
00251 //   Nuclear interactions
00252 //------------------------ 
00253 
00254   if ( NuclearInteraction && abs(myTrack.pid()) > 100 
00255                           && abs(myTrack.pid()) < 1000000) { 
00256 
00257     // Simulate a nuclear interaction
00258     double factor = 1.0;
00259     if (layer.layerNumber() >= 19 && layer.layerNumber() <= 27 ) 
00260       factor = theTECFudgeFactor;
00261     NuclearInteraction->updateState(myTrack,radlen*factor);
00262 
00263     if ( NuclearInteraction->nDaughters() ) { 
00264 
00265       //add a end vertex to the mother particle
00266       int ivertex = mySimEvent.addSimVertex(myTrack.vertex(),itrack,
00267                                             FSimVertexType::NUCL_VERTEX);
00268       
00269       // This was a hadron that interacted inelastically
00270       int idaugh = 0;
00271       for ( DaughterIter = NuclearInteraction->beginDaughters();
00272             DaughterIter != NuclearInteraction->endDaughters(); 
00273             ++DaughterIter) {
00274 
00275         // The daughter in the event
00276         int daughId = mySimEvent.addSimTrack(&(*DaughterIter), ivertex);
00277         
00278         // Store the closest daughter in the mother info (for later tracking purposes)
00279         if ( NuclearInteraction->closestDaughterId() == idaugh++ ) {
00280           if ( mySimEvent.track(itrack).vertex().position().Pt() < 4.0 ) 
00281             mySimEvent.track(itrack).setClosestDaughterId(daughId);
00282         }
00283 
00284       }
00285       // The hadron is destroyed. Return.
00286       return;
00287     }
00288     
00289   }
00290 
00291   if ( myTrack.charge() == 0 ) return;
00292 
00293   if ( !Bremsstrahlung && !EnergyLoss && !MultipleScattering ) return;
00294 
00295 //----------------
00296 //  Bremsstrahlung
00297 //----------------
00298 
00299   if ( Bremsstrahlung && abs(myTrack.pid())==11 ) {
00300         
00301     Bremsstrahlung->updateState(myTrack,radlen);
00302 
00303     if ( Bremsstrahlung->nDaughters() ) {
00304       
00305       // Add a vertex, but do not attach it to the electron, because it 
00306       // continues its way...
00307       int ivertex = mySimEvent.addSimVertex(myTrack.vertex(),itrack,
00308                                             FSimVertexType::BREM_VERTEX);
00309 
00310       for ( DaughterIter = Bremsstrahlung->beginDaughters();
00311             DaughterIter != Bremsstrahlung->endDaughters(); 
00312             ++DaughterIter) {
00313         mySimEvent.addSimTrack(&(*DaughterIter), ivertex);
00314 
00315       }
00316       
00317     }
00318     
00319   } 
00320 
00321 
00325 
00326   if ( EnergyLoss )
00327   {
00328     theEnergyLoss = myTrack.E();
00329     EnergyLoss->updateState(myTrack,radlen);
00330     theEnergyLoss -= myTrack.E();
00331   }
00332   
00333 
00337 
00338   if ( MultipleScattering && myTrack.Pt() > pTmin ) {
00339     //    MultipleScattering->setNormalVector(normalVector(layer,myTrack));
00340     MultipleScattering->setNormalVector(theNormalVector);
00341     MultipleScattering->updateState(myTrack,radlen);
00342   }
00343     
00344 }
00345 
00346 double
00347 MaterialEffects::radLengths(const TrackerLayer& layer,
00348                             ParticlePropagator& myTrack) {
00349 
00350   // Thickness of layer
00351   theThickness = layer.surface().mediumProperties()->radLen();
00352 
00353   GlobalVector P(myTrack.Px(),myTrack.Py(),myTrack.Pz());
00354   
00355   // Effective length of track inside layer (considering crossing angle)
00356   //  double radlen = theThickness / fabs(P.dot(theNormalVector)/(P.mag()*theNormalVector.mag()));
00357   double radlen = theThickness / fabs(P.dot(theNormalVector)) * P.mag();
00358 
00359   // This is a series of fudge factors (from the geometry description), 
00360   // to describe the layer inhomogeneities (services, cables, supports...)
00361   double rad = myTrack.R();
00362   double zed = fabs(myTrack.Z());
00363 
00364   double factor = 1;
00365 
00366   // Are there fudge factors for this layer
00367   if ( layer.fudgeNumber() ) 
00368 
00369     // If yes, loop on them
00370     for ( unsigned int iLayer=0; iLayer < layer.fudgeNumber(); ++iLayer ) { 
00371 
00372       // Apply to R if forward layer, to Z if barrel layer
00373       if ( (  layer.forward() && layer.fudgeMin(iLayer) < rad && rad < layer.fudgeMax(iLayer) )  ||
00374            ( !layer.forward() && layer.fudgeMin(iLayer) < zed && zed < layer.fudgeMax(iLayer) ) ) {
00375         factor = layer.fudgeFactor(iLayer);
00376         break;
00377       }
00378     
00379   }
00380 
00381   theThickness *= factor;
00382 
00383   return radlen * factor;
00384 
00385 }
00386 
00387 GlobalVector
00388 MaterialEffects::normalVector(const TrackerLayer& layer,
00389                               ParticlePropagator& myTrack ) const {
00390   return layer.forward() ?  
00391     layer.disk()->normalVector() :
00392     GlobalVector(myTrack.X(),myTrack.Y(),0.)/myTrack.R();
00393 }
00394 
00395 void 
00396 MaterialEffects::save() { 
00397 
00398   // Save current nuclear interactions in the event libraries.
00399   if ( NuclearInteraction ) NuclearInteraction->save();
00400 
00401 }