00001
00002
00003 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00004 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00005
00006
00007
00008
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), use_hardcoded(1)
00031 {
00032
00033
00034 use_hardcoded = matEff.getParameter<bool >("use_hardcoded_geometry");
00035
00036 bool doPairProduction = matEff.getParameter<bool>("PairProduction");
00037 bool doBremsstrahlung = matEff.getParameter<bool>("Bremsstrahlung");
00038 bool doEnergyLoss = matEff.getParameter<bool>("EnergyLoss");
00039 bool doMultipleScattering = matEff.getParameter<bool>("MultipleScattering");
00040 bool doNuclearInteraction = matEff.getParameter<bool>("NuclearInteraction");
00041 bool doMuonBremsstrahlung = matEff.getParameter<bool>("MuonBremsstrahlung");
00042
00043 double A = matEff.getParameter<double>("A");
00044 double Z = matEff.getParameter<double>("Z");
00045 double density = matEff.getParameter<double>("Density");
00046 double radLen = matEff.getParameter<double>("RadiationLength");
00047
00048
00049
00050 if ( doPairProduction ) {
00051
00052 double photonEnergy = matEff.getParameter<double>("photonEnergy");
00053 PairProduction = new PairProductionSimulator(photonEnergy,
00054 random);
00055
00056 }
00057
00058 if ( doBremsstrahlung ) {
00059
00060 double bremEnergy = matEff.getParameter<double>("bremEnergy");
00061 double bremEnergyFraction = matEff.getParameter<double>("bremEnergyFraction");
00062 Bremsstrahlung = new BremsstrahlungSimulator(bremEnergy,
00063 bremEnergyFraction,
00064 random);
00065
00066 }
00067
00068 if ( doMuonBremsstrahlung ) {
00069
00070 double bremEnergy = matEff.getParameter<double>("bremEnergy");
00071 double bremEnergyFraction = matEff.getParameter<double>("bremEnergyFraction");
00072 MuonBremsstrahlung = new MuonBremsstrahlungSimulator(random,A,Z,density,radLen,bremEnergy,
00073 bremEnergyFraction);
00074
00075 }
00076
00077
00078
00079
00080 if ( doEnergyLoss ) {
00081
00082 pTmin = matEff.getParameter<double>("pTmin");
00083 EnergyLoss = new EnergyLossSimulator(random,A,Z,density,radLen);
00084
00085 }
00086
00087 if ( doMultipleScattering ) {
00088
00089 MultipleScattering = new MultipleScatteringSimulator(random,A,Z,density,radLen);
00090
00091 }
00092
00093 if ( doNuclearInteraction ) {
00094
00095
00096 std::vector<double> pionEnergies
00097 = matEff.getUntrackedParameter<std::vector<double> >("pionEnergies");
00098
00099
00100 std::vector<int> pionTypes
00101 = matEff.getUntrackedParameter<std::vector<int> >("pionTypes");
00102
00103
00104 std::vector<std::string> pionNames
00105 = matEff.getUntrackedParameter<std::vector<std::string> >("pionNames");
00106
00107
00108 std::vector<double> pionMasses
00109 = matEff.getUntrackedParameter<std::vector<double> >("pionMasses");
00110
00111
00112 std::vector<double> pionPMin
00113 = matEff.getUntrackedParameter<std::vector<double> >("pionMinP");
00114
00115
00116 std::vector<double> lengthRatio
00117 = matEff.getParameter<std::vector<double> >("lengthRatio");
00118
00119
00120
00121
00122
00123
00124 theTECFudgeFactor = matEff.getParameter<double>("fudgeFactor");
00125
00126
00127 std::vector<double> theRatios
00128 = matEff.getUntrackedParameter<std::vector<double> >("ratios");
00129
00130
00131
00132
00133
00134
00135 std::vector< std::vector<double> > ratios;
00136 ratios.resize(pionTypes.size());
00137 for ( unsigned i=0; i<pionTypes.size(); ++i ) {
00138 for ( unsigned j=0; j<pionEnergies.size(); ++j ) {
00139 ratios[i].push_back(theRatios[ i*pionEnergies.size() + j ]);
00140 }
00141 }
00142
00143
00144 double pionEnergy
00145 = matEff.getParameter<double>("pionEnergy");
00146
00147
00148
00149 unsigned distAlgo
00150 = matEff.getParameter<unsigned>("distAlgo");
00151 double distCut
00152 = matEff.getParameter<double>("distCut");
00153
00154
00155
00156 std::string inputFile
00157 = matEff.getUntrackedParameter<std::string>("inputFile");
00158
00159
00160 std::map<int,int> idMap;
00161
00162 std::vector<int> idProtons
00163 = matEff.getUntrackedParameter<std::vector<int> >("protons");
00164 for ( unsigned i=0; i<idProtons.size(); ++i )
00165 idMap[idProtons[i]] = 2212;
00166
00167 std::vector<int> idAntiProtons
00168 = matEff.getUntrackedParameter<std::vector<int> >("antiprotons");
00169 for ( unsigned i=0; i<idAntiProtons.size(); ++i )
00170 idMap[idAntiProtons[i]] = -2212;
00171
00172 std::vector<int> idNeutrons
00173 = matEff.getUntrackedParameter<std::vector<int> >("neutrons");
00174 for ( unsigned i=0; i<idNeutrons.size(); ++i )
00175 idMap[idNeutrons[i]] = 2112;
00176
00177 std::vector<int> idAntiNeutrons
00178 = matEff.getUntrackedParameter<std::vector<int> >("antineutrons");
00179 for ( unsigned i=0; i<idAntiNeutrons.size(); ++i )
00180 idMap[idAntiNeutrons[i]] = -2112;
00181
00182 std::vector<int> idK0Ls
00183 = matEff.getUntrackedParameter<std::vector<int> >("K0Ls");
00184 for ( unsigned i=0; i<idK0Ls.size(); ++i )
00185 idMap[idK0Ls[i]] = 130;
00186
00187 std::vector<int> idKplusses
00188 = matEff.getUntrackedParameter<std::vector<int> >("Kplusses");
00189 for ( unsigned i=0; i<idKplusses.size(); ++i )
00190 idMap[idKplusses[i]] = 321;
00191
00192 std::vector<int> idKminusses
00193 = matEff.getUntrackedParameter<std::vector<int> >("Kminusses");
00194 for ( unsigned i=0; i<idKminusses.size(); ++i )
00195 idMap[idKminusses[i]] = -321;
00196
00197 std::vector<int> idPiplusses
00198 = matEff.getUntrackedParameter<std::vector<int> >("Piplusses");
00199 for ( unsigned i=0; i<idPiplusses.size(); ++i )
00200 idMap[idPiplusses[i]] = 211;
00201
00202 std::vector<int> idPiminusses
00203 = matEff.getUntrackedParameter<std::vector<int> >("Piminusses");
00204 for ( unsigned i=0; i<idPiminusses.size(); ++i )
00205 idMap[idPiminusses[i]] = -211;
00206
00207
00208 NuclearInteraction =
00209 new NuclearInteractionSimulator(pionEnergies, pionTypes, pionNames,
00210 pionMasses, pionPMin, pionEnergy,
00211 lengthRatio, ratios, idMap,
00212 inputFile, distAlgo, distCut, random);
00213 }
00214
00215 }
00216
00217
00218 MaterialEffects::~MaterialEffects() {
00219
00220 if ( PairProduction ) delete PairProduction;
00221 if ( Bremsstrahlung ) delete Bremsstrahlung;
00222 if ( EnergyLoss ) delete EnergyLoss;
00223 if ( MultipleScattering ) delete MultipleScattering;
00224 if ( NuclearInteraction ) delete NuclearInteraction;
00225
00226 if ( MuonBremsstrahlung ) delete MuonBremsstrahlung;
00227 }
00228
00229 void MaterialEffects::interact(FSimEvent& mySimEvent,
00230 const TrackerLayer& layer,
00231 ParticlePropagator& myTrack,
00232 unsigned itrack) {
00233
00234 MaterialEffectsSimulator::RHEP_const_iter DaughterIter;
00235 double radlen;
00236 theEnergyLoss = 0;
00237 theNormalVector = normalVector(layer,myTrack);
00238 radlen = radLengths(layer,myTrack);
00239
00240
00241
00242
00243
00244 if ( PairProduction && myTrack.pid()==22 ) {
00245
00246
00247 PairProduction->updateState(myTrack,radlen);
00248
00249 if ( PairProduction->nDaughters() ) {
00250
00251 int ivertex = mySimEvent.addSimVertex(myTrack.vertex(),itrack,
00252 FSimVertexType::PAIR_VERTEX);
00253
00254
00255 if (ivertex>=0) {
00256
00257 for ( DaughterIter = PairProduction->beginDaughters();
00258 DaughterIter != PairProduction->endDaughters();
00259 ++DaughterIter) {
00260
00261 mySimEvent.addSimTrack(&(*DaughterIter), ivertex);
00262
00263 }
00264
00265 return;
00266 }
00267 else {
00268 edm::LogWarning("MaterialEffects") << " WARNING: A non valid vertex was found in photon conv. -> " << ivertex << std::endl;
00269 }
00270
00271 }
00272
00273 }
00274
00275 if ( myTrack.pid() == 22 ) return;
00276
00277
00278
00279
00280
00281 if ( NuclearInteraction && abs(myTrack.pid()) > 100
00282 && abs(myTrack.pid()) < 1000000) {
00283
00284
00285 double factor = 1.0;
00286 if(use_hardcoded){
00287 if (layer.layerNumber() >= 19 && layer.layerNumber() <= 27 )
00288 factor = theTECFudgeFactor;
00289 }
00290 NuclearInteraction->updateState(myTrack,radlen*factor);
00291
00292 if ( NuclearInteraction->nDaughters() ) {
00293
00294
00295 int ivertex = mySimEvent.addSimVertex(myTrack.vertex(),itrack,
00296 FSimVertexType::NUCL_VERTEX);
00297
00298
00299 if (ivertex>=0) {
00300
00301 int idaugh = 0;
00302 for ( DaughterIter = NuclearInteraction->beginDaughters();
00303 DaughterIter != NuclearInteraction->endDaughters();
00304 ++DaughterIter) {
00305
00306
00307 int daughId = mySimEvent.addSimTrack(&(*DaughterIter), ivertex);
00308
00309
00310 if ( NuclearInteraction->closestDaughterId() == idaugh++ ) {
00311 if ( mySimEvent.track(itrack).vertex().position().Pt() < 4.0 )
00312 mySimEvent.track(itrack).setClosestDaughterId(daughId);
00313 }
00314
00315 }
00316
00317 return;
00318 }
00319 else {
00320 edm::LogWarning("MaterialEffects") << " WARNING: A non valid vertex was found in nucl. int. -> " << ivertex << std::endl;
00321 }
00322
00323 }
00324
00325 }
00326
00327 if ( myTrack.charge() == 0 ) return;
00328
00329 if ( !MuonBremsstrahlung && !Bremsstrahlung && !EnergyLoss && !MultipleScattering ) return;
00330
00331
00332
00333
00334
00335 if ( Bremsstrahlung && abs(myTrack.pid())==11 ) {
00336
00337 Bremsstrahlung->updateState(myTrack,radlen);
00338
00339 if ( Bremsstrahlung->nDaughters() ) {
00340
00341
00342
00343 int ivertex = mySimEvent.addSimVertex(myTrack.vertex(),itrack,
00344 FSimVertexType::BREM_VERTEX);
00345
00346
00347 if (ivertex>=0) {
00348 for ( DaughterIter = Bremsstrahlung->beginDaughters();
00349 DaughterIter != Bremsstrahlung->endDaughters();
00350 ++DaughterIter) {
00351 mySimEvent.addSimTrack(&(*DaughterIter), ivertex);
00352 }
00353 }
00354 else {
00355 edm::LogWarning("MaterialEffects") << " WARNING: A non valid vertex was found in brem -> " << ivertex << std::endl;
00356 }
00357
00358 }
00359
00360 }
00361
00362
00363
00364
00365
00366 if ( MuonBremsstrahlung && abs(myTrack.pid())==13 ) {
00367
00368 MuonBremsstrahlung->updateState(myTrack,radlen);
00369
00370 if ( MuonBremsstrahlung->nDaughters() ) {
00371
00372
00373
00374 int ivertex = mySimEvent.addSimVertex(myTrack.vertex(),itrack,
00375 FSimVertexType::BREM_VERTEX);
00376
00377
00378 if (ivertex>=0) {
00379 for ( DaughterIter = MuonBremsstrahlung->beginDaughters();
00380 DaughterIter != MuonBremsstrahlung->endDaughters();
00381 ++DaughterIter) {
00382 mySimEvent.addSimTrack(&(*DaughterIter), ivertex);
00383 }
00384 }
00385 else {
00386 edm::LogWarning("MaterialEffects") << " WARNING: A non valid vertex was found in muon brem -> " << ivertex << std::endl;
00387 }
00388
00389 }
00390
00391 }
00392
00396
00397 if ( EnergyLoss )
00398 {
00399 theEnergyLoss = myTrack.E();
00400 EnergyLoss->updateState(myTrack,radlen);
00401 theEnergyLoss -= myTrack.E();
00402 }
00403
00404
00408
00409 if ( MultipleScattering && myTrack.Pt() > pTmin ) {
00410
00411 MultipleScattering->setNormalVector(theNormalVector);
00412 MultipleScattering->updateState(myTrack,radlen);
00413 }
00414
00415 }
00416
00417 double
00418 MaterialEffects::radLengths(const TrackerLayer& layer,
00419 ParticlePropagator& myTrack) {
00420
00421
00422 theThickness = layer.surface().mediumProperties().radLen();
00423
00424 GlobalVector P(myTrack.Px(),myTrack.Py(),myTrack.Pz());
00425
00426
00427
00428 double radlen = theThickness / fabs(P.dot(theNormalVector)) * P.mag();
00429
00430
00431
00432 double rad = myTrack.R();
00433 double zed = fabs(myTrack.Z());
00434
00435 double factor = 1;
00436
00437
00438 if ( layer.fudgeNumber() )
00439
00440
00441 for ( unsigned int iLayer=0; iLayer < layer.fudgeNumber(); ++iLayer ) {
00442
00443
00444 if ( ( layer.forward() && layer.fudgeMin(iLayer) < rad && rad < layer.fudgeMax(iLayer) ) ||
00445 ( !layer.forward() && layer.fudgeMin(iLayer) < zed && zed < layer.fudgeMax(iLayer) ) ) {
00446 factor = layer.fudgeFactor(iLayer);
00447 break;
00448 }
00449
00450 }
00451
00452 theThickness *= factor;
00453
00454 return radlen * factor;
00455
00456 }
00457
00458 GlobalVector
00459 MaterialEffects::normalVector(const TrackerLayer& layer,
00460 ParticlePropagator& myTrack ) const {
00461 return layer.forward() ?
00462 layer.disk()->normalVector() :
00463 GlobalVector(myTrack.X(),myTrack.Y(),0.)/myTrack.R();
00464 }
00465
00466 void
00467 MaterialEffects::save() {
00468
00469
00470 if ( NuclearInteraction ) NuclearInteraction->save();
00471
00472 }