CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_5_3_14/src/FastSimulation/MaterialEffects/src/MaterialEffects.cc

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