CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_9_patch3/src/FastSimulation/MuonSimHitProducer/src/MuonSimHitProducer.cc

Go to the documentation of this file.
00001 //
00002 // Package:    MuonSimHitProducer
00003 // Class:      MuonSimHitProducer
00004 // 
00014 //
00015 // Original Author:  Martijn Mulders/Matthew Jones
00016 //         Created:  Wed Jul 30 11:37:24 CET 2007
00017 //         Working:  Fri Nov  9 09:39:33 CST 2007
00018 //
00019 // $Id: MuonSimHitProducer.cc,v 1.36 2011/10/07 08:25:42 aperrott Exp $
00020 //
00021 //
00022 
00023 // CMSSW headers 
00024 #include "FWCore/ServiceRegistry/interface/Service.h"
00025 #include "FWCore/Utilities/interface/RandomNumberGenerator.h"
00026 #include "FWCore/Framework/interface/MakerMacros.h"
00027 
00028 // Fast Simulation headers
00029 #include "FastSimulation/Utilities/interface/RandomEngine.h"
00030 #include "FastSimulation/MuonSimHitProducer/interface/MuonSimHitProducer.h"
00031 #include "FastSimulation/MaterialEffects/interface/MaterialEffects.h"
00032 #include "FastSimulation/MaterialEffects/interface/MultipleScatteringSimulator.h"
00033 #include "FastSimulation/MaterialEffects/interface/EnergyLossSimulator.h"
00034 #include "FastSimulation/ParticlePropagator/interface/ParticlePropagator.h"
00035 #include "FastSimulation/MaterialEffects/interface/MuonBremsstrahlungSimulator.h"
00036 
00037 // SimTrack
00038 #include "SimDataFormats/Track/interface/SimTrack.h"
00039 #include "SimDataFormats/Vertex/interface/SimVertex.h"
00040 #include "SimDataFormats/TrackingHit/interface/PSimHitContainer.h"
00041 
00042 // STL headers 
00043 #include <vector>
00044 #include <iostream>
00045 
00046 // RecoMuon headers
00047 #include "RecoMuon/TrackingTools/interface/MuonServiceProxy.h"
00048 #include "RecoMuon/Navigation/interface/DirectMuonNavigation.h"
00049 #include "RecoMuon/MeasurementDet/interface/MuonDetLayerMeasurements.h"
00050 #include "RecoMuon/TrackingTools/interface/MuonPatternRecoDumper.h"
00051 
00052 // Tracking Tools
00053 #include "TrackingTools/GeomPropagators/interface/HelixArbitraryPlaneCrossing.h"
00054 #include "TrackPropagation/SteppingHelixPropagator/interface/SteppingHelixPropagator.h"
00055 
00056 // Data Formats
00057 #include "DataFormats/MuonDetId/interface/DTWireId.h"
00058 #include "DataFormats/GeometrySurface/interface/PlaneBuilder.h"
00059 #include "DataFormats/GeometrySurface/interface/TangentPlane.h"
00060 
00061 
00063 // Geometry, Magnetic Field
00064 #include "Geometry/DTGeometry/interface/DTGeometry.h"
00065 #include "Geometry/CSCGeometry/interface/CSCGeometry.h"
00066 #include "Geometry/RPCGeometry/interface/RPCGeometry.h"
00067 #include "Geometry/Records/interface/MuonGeometryRecord.h"
00068 #include "MagneticField/Records/interface/IdealMagneticFieldRecord.h"
00069 
00071 
00072 // #include "TrackingTools/TrackAssociator/interface/TrackDetectorAssociator.h"
00073 
00074 //for debug only 
00075 //#define FAMOS_DEBUG
00076 
00077 //
00078 // constructors and destructor
00079 //
00080 MuonSimHitProducer::MuonSimHitProducer(const edm::ParameterSet& iConfig):
00081   theEstimator(iConfig.getParameter<double>("Chi2EstimatorCut")),
00082   propagatorWithoutMaterial(0)
00083  {
00084 
00085   //
00086   //  Initialize the random number generator service
00087   //
00088   edm::Service<edm::RandomNumberGenerator> rng;
00089   if ( ! rng.isAvailable() ) {
00090     throw cms::Exception("Configuration") <<
00091       "MuonSimHitProducer requires the RandomGeneratorService \n"
00092       "which is not present in the configuration file. \n"
00093       "You must add the service in the configuration file\n"
00094       "or remove the module that requires it.";
00095   }
00096 
00097   random = new RandomEngine(&(*rng));
00098 
00099   // Read relevant parameters
00100   readParameters(iConfig.getParameter<edm::ParameterSet>("MUONS"),
00101                  iConfig.getParameter<edm::ParameterSet>("TRACKS"),
00102                  iConfig.getParameter<edm::ParameterSet>("MaterialEffectsForMuons"));
00103 
00104   //
00105   //  register your products ... need to declare at least one possible product...
00106   //
00107   produces<edm::PSimHitContainer>("MuonCSCHits");
00108   produces<edm::PSimHitContainer>("MuonDTHits");
00109   produces<edm::PSimHitContainer>("MuonRPCHits");
00110 
00111   edm::ParameterSet serviceParameters =
00112      iConfig.getParameter<edm::ParameterSet>("ServiceParameters");
00113   theService = new MuonServiceProxy(serviceParameters);
00114 }
00115 
00116 // ---- method called once each job just before starting event loop ----
00117 void 
00118 MuonSimHitProducer::beginRun(edm::Run & run, const edm::EventSetup & es) {
00119 
00120   //services
00121   edm::ESHandle<MagneticField>          magField;
00122   edm::ESHandle<DTGeometry>             dtGeometry;
00123   edm::ESHandle<CSCGeometry>            cscGeometry;
00124   edm::ESHandle<RPCGeometry>            rpcGeometry;
00125   edm::ESHandle<Propagator>             propagator;
00126 
00127   es.get<IdealMagneticFieldRecord>().get(magField);
00128   es.get<MuonGeometryRecord>().get("MisAligned",dtGeometry);
00129   es.get<MuonGeometryRecord>().get("MisAligned",cscGeometry);
00130   es.get<MuonGeometryRecord>().get(rpcGeometry);
00131 
00132   magfield = &(*magField);
00133   dtGeom = &(*dtGeometry);
00134   cscGeom = &(*cscGeometry);
00135   rpcGeom = &(*rpcGeometry);
00136 
00137   theService->update(es);
00138 
00139   // A few propagators
00140   propagatorWithMaterial = &(*(theService->propagator("SteppingHelixPropagatorAny")));
00141   propagatorWithoutMaterial = propagatorWithMaterial->clone();
00142   SteppingHelixPropagator* SHpropagator = dynamic_cast<SteppingHelixPropagator*> (propagatorWithoutMaterial); // Beuark!
00143   SHpropagator->setMaterialMode(true); // switches OFF material effects;
00144 
00145 }
00146   
00147 MuonSimHitProducer::~MuonSimHitProducer()
00148 {
00149   // do anything here that needs to be done at destruction time
00150   // (e.g. close files, deallocate resources etc.)
00151   
00152   if ( random ) { 
00153     delete random;
00154   }
00155 
00156   if ( theMaterialEffects ) delete theMaterialEffects;
00157   if ( propagatorWithoutMaterial) delete propagatorWithoutMaterial;
00158 }
00159 
00160 
00161 //
00162 // member functions
00163 //
00164 
00165 // ------------ method called to produce the data  ------------
00166 
00167 void 
00168 MuonSimHitProducer::produce(edm::Event& iEvent,const edm::EventSetup& iSetup) {
00169   // using namespace edm;
00170   // using namespace std;
00171 
00172   MuonPatternRecoDumper dumper;
00173 
00174   edm::Handle<std::vector<SimTrack> > simMuons;
00175   edm::Handle<std::vector<SimVertex> > simVertices;
00176   std::vector<PSimHit> theCSCHits;
00177   std::vector<PSimHit> theDTHits;
00178   std::vector<PSimHit> theRPCHits;
00179 
00180   DirectMuonNavigation navigation(theService->detLayerGeometry());
00181   iEvent.getByLabel(theSimModuleLabel_,theSimModuleProcess_,simMuons);
00182   iEvent.getByLabel(theSimModuleLabel_,simVertices);
00183 
00184   for ( unsigned int itrk=0; itrk<simMuons->size(); itrk++ ) {
00185     const SimTrack &mySimTrack = (*simMuons)[itrk];
00186     math::XYZTLorentzVector mySimP4(mySimTrack.momentum().x(),
00187                                     mySimTrack.momentum().y(),
00188                                     mySimTrack.momentum().z(),
00189                                     mySimTrack.momentum().t());
00190 
00191     // Decaying hadrons are now in the list, and so are their muon daughter
00192     // Ignore the hadrons here.
00193     int pid = mySimTrack.type(); 
00194     if ( abs(pid) != 13 ) continue;
00195 
00196     double t0 = 0;
00197     GlobalPoint initialPosition;
00198     int ivert = mySimTrack.vertIndex();
00199     if ( ivert >= 0 ) {
00200       t0 = (*simVertices)[ivert].position().t();
00201       GlobalPoint xyzzy((*simVertices)[ivert].position().x(),
00202                         (*simVertices)[ivert].position().y(),
00203                         (*simVertices)[ivert].position().z());
00204       initialPosition = xyzzy;
00205     }
00206 //
00207 //  Presumably t0 has dimensions of cm if not mm?
00208 //  Convert to ns for internal calculations.
00209 //  I wonder where we should get c from?
00210 //
00211     double tof = t0/29.98;
00212 
00213 #ifdef FAMOS_DEBUG
00214     std::cout << " ===> MuonSimHitProducer::reconstruct() found SIMTRACK - pid = "
00215               << pid ;
00216     std::cout << " : pT = " << mySimP4.Pt()
00217               << ", eta = " << mySimP4.Eta()
00218               << ", phi = " << mySimP4.Phi() << std::endl;
00219 #endif
00220 
00221 //
00222 //  Produce muons sim hits starting from undecayed simulated muons
00223 //
00224 
00225     GlobalPoint startingPosition(mySimTrack.trackerSurfacePosition().x(),
00226                                  mySimTrack.trackerSurfacePosition().y(),
00227                                  mySimTrack.trackerSurfacePosition().z());
00228     GlobalVector startingMomentum(mySimTrack.trackerSurfaceMomentum().x(),
00229                                   mySimTrack.trackerSurfaceMomentum().y(),
00230                                   mySimTrack.trackerSurfaceMomentum().z());
00231 //
00232 //  Crap... there's no time-of-flight to the trackerSurfacePosition()...
00233 //  So, this will be wrong when the curvature can't be neglected, but that
00234 //  will be rather seldom...  May as well ignore the mass too.
00235 //
00236     GlobalVector dtracker = startingPosition-initialPosition;
00237     tof += dtracker.mag()/29.98;
00238 
00239 #ifdef FAMOS_DEBUG
00240     std::cout << " the Muon START position " << startingPosition << std::endl;
00241     std::cout << " the Muon START momentum " << startingMomentum << std::endl;
00242 #endif
00243 
00244 // 
00245 //  Some magic to define a TrajectoryStateOnSurface
00246 //
00247     PlaneBuilder pb;
00248     GlobalVector zAxis = startingMomentum.unit();
00249     GlobalVector yAxis(zAxis.y(),-zAxis.x(),0); 
00250     GlobalVector xAxis = yAxis.cross(zAxis);
00251     Surface::RotationType rot = Surface::RotationType(xAxis,yAxis,zAxis);
00252     PlaneBuilder::ReturnType startingPlane = pb.plane(startingPosition,rot);
00253     GlobalTrajectoryParameters gtp(startingPosition,
00254                                    startingMomentum,
00255                                    (int)mySimTrack.charge(),
00256                                    magfield);
00257     TrajectoryStateOnSurface startingState(gtp,*startingPlane);
00258 
00259     std::vector<const DetLayer *> navLayers;
00260     if ( fabs(startingState.globalMomentum().eta()) > 4.5 ) {
00261       navLayers = navigation.compatibleEndcapLayers(*(startingState.freeState()),
00262                                                     alongMomentum);
00263     }
00264     else {
00265       navLayers = navigation.compatibleLayers(*(startingState.freeState()),
00266                                                alongMomentum);
00267     }
00268     /*
00269     edm::ESHandle<Propagator> propagator =
00270       theService->propagator("SteppingHelixPropagatorAny");
00271     */
00272 
00273     if ( navLayers.empty() ) continue;
00274 
00275 #ifdef FAMOS_DEBUG
00276       std::cout << "Found " << navLayers.size()
00277                 << " compatible DetLayers..." << std::endl;
00278 #endif
00279 
00280     TrajectoryStateOnSurface propagatedState = startingState;
00281     for ( unsigned int ilayer=0; ilayer<navLayers.size(); ilayer++ ) {
00282 
00283 #ifdef FAMOS_DEBUG
00284       std::cout << "Propagating to layer " << ilayer << " " 
00285                 << dumper.dumpLayer(navLayers[ilayer]) 
00286                 << std::endl;
00287 #endif
00288       
00289       std::vector<DetWithState> comps = 
00290         navLayers[ilayer]->compatibleDets(propagatedState,*propagatorWithMaterial,theEstimator);
00291       if ( comps.empty() ) continue;
00292 
00293 #ifdef FAMOS_DEBUG
00294       std::cout << "Propagating " << propagatedState << std::endl;
00295 #endif
00296 
00297       // Starting momentum
00298       double pi = propagatedState.globalMomentum().mag(); 
00299 
00300       // Propagate with material effects (dE/dx average only)
00301       SteppingHelixStateInfo shsStart(*(propagatedState.freeTrajectoryState()));
00302       const SteppingHelixStateInfo& shsDest = 
00303         ((const SteppingHelixPropagator*)propagatorWithMaterial)->propagate(shsStart,navLayers[ilayer]->surface());
00304       std::pair<TrajectoryStateOnSurface,double> next(shsDest.getStateOnSurface(navLayers[ilayer]->surface()),
00305                                                       shsDest.path());
00306       // No need to continue if there is no valid propagation available.
00307       // This happens rarely (~0.1% of ttbar events)
00308       if ( !next.first.isValid() ) continue; 
00309       // This is the estimate of the number of radiation lengths traversed, 
00310       // together with the total path length 
00311       double radPath = shsDest.radPath();
00312       double pathLength = next.second;
00313 
00314       // Now propagate without dE/dx (average) 
00315       // [To add the dE/dx fluctuations to the actual dE/dx]
00316       std::pair<TrajectoryStateOnSurface,double> nextNoMaterial = 
00317         propagatorWithoutMaterial->propagateWithPath(propagatedState,navLayers[ilayer]->surface());
00318 
00319       // Update the propagated state
00320       propagatedState = next.first;
00321       double pf = propagatedState.globalMomentum().mag();
00322 
00323       // Insert dE/dx fluctuations and multiple scattering
00324       // Skip this step if nextNoMaterial.first is not valid 
00325       // This happens rarely (~0.02% of ttbar events)
00326       if ( theMaterialEffects && nextNoMaterial.first.isValid() ) applyMaterialEffects(propagatedState,nextNoMaterial.first,radPath);
00327       // Check that the 'shaken' propagatedState is still valid, otherwise continue
00328       if ( !propagatedState.isValid() ) continue; 
00329       // (No evidence that this ever happens)
00330 //
00331 //  Consider this... 1 GeV muon has a velocity that is only 0.5% slower than c...
00332 //  We probably can safely ignore the mass for anything that makes it out to the
00333 //  muon chambers.
00334 //
00335       double pavg = 0.5*(pi+pf);
00336       double m = mySimP4.M();
00337       double rbeta = sqrt(1+m*m/(pavg*pavg))/29.98;
00338       double dtof = pathLength*rbeta;
00339 
00340 #ifdef FAMOS_DEBUG
00341       std::cout << "Propagated to next surface... path length = " << pathLength 
00342                 << " cm, dTOF = " << dtof << " ns" << std::endl;
00343 #endif
00344 
00345       tof += dtof;
00346 
00347       for ( unsigned int icomp=0; icomp<comps.size(); icomp++ )
00348         {
00349       const GeomDet *gd = comps[icomp].first;
00350       if ( gd->subDetector() == GeomDetEnumerators::DT ) {
00351         DTChamberId id(gd->geographicalId());
00352         const DTChamber *chamber = dtGeom->chamber(id);
00353         std::vector<const DTSuperLayer *> superlayer = chamber->superLayers();
00354         for ( unsigned int isl=0; isl<superlayer.size(); isl++ ) {
00355           std::vector<const DTLayer *> layer = superlayer[isl]->layers();
00356           for ( unsigned int ilayer=0; ilayer<layer.size(); ilayer++ ) {
00357             DTLayerId lid = layer[ilayer]->id();
00358 #ifdef FAMOS_DEBUG
00359             std::cout << "    Extrapolated to DT (" 
00360                       << lid.wheel() << "," 
00361                       << lid.station() << "," 
00362                       << lid.sector() << "," 
00363                       << lid.superlayer() << "," 
00364                       << lid.layer() << ")" << std::endl;
00365 #endif
00366 
00367             const GeomDetUnit *det = dtGeom->idToDetUnit(lid);
00368 
00369             HelixArbitraryPlaneCrossing crossing(propagatedState.globalPosition().basicVector(),
00370                                                  propagatedState.globalMomentum().basicVector(),
00371                                                  propagatedState.transverseCurvature(),
00372                                                  anyDirection);
00373             std::pair<bool,double> path = crossing.pathLength(det->surface());
00374             if ( ! path.first ) continue;
00375             LocalPoint lpos = det->toLocal(GlobalPoint(crossing.position(path.second)));
00376             if ( ! det->surface().bounds().inside(lpos) ) continue;
00377             const DTTopology& dtTopo = layer[ilayer]->specificTopology();
00378             int wire = dtTopo.channel(lpos);
00379             if (wire - dtTopo.firstChannel() < 0 || wire - dtTopo.lastChannel() > 0 ) continue;
00380             // no drift cell here (on the chamber edge or just outside)
00381             // this hit would otherwise be discarded downstream in the digitizer
00382 
00383             DTWireId wid(lid,wire);
00384             double thickness = det->surface().bounds().thickness();
00385             LocalVector lmom = det->toLocal(GlobalVector(crossing.direction(path.second)));
00386             lmom = lmom.unit()*propagatedState.localMomentum().mag();
00387 
00388             // Factor that takes into account the (rec)hits lost because of delta's, etc.:
00389             // (Not fully satisfactory patch, but it seems to work...)
00390             double pmu = lmom.mag();
00391             double theDTHitIneff = pmu>0? exp(kDT*log(pmu)+fDT):0.;
00392             if (random->flatShoot()<theDTHitIneff) continue;
00393 
00394             double eloss = 0;
00395             double pz = fabs(lmom.z());
00396             LocalPoint entry = lpos - 0.5*thickness*lmom/pz;
00397             LocalPoint exit = lpos + 0.5*thickness*lmom/pz;
00398             double dtof = path.second*rbeta;
00399             int trkid = mySimTrack.trackId();
00400             unsigned int id = wid.rawId();
00401             short unsigned int processType = 2;
00402             PSimHit hit(entry,exit,lmom.mag(),
00403                         tof+dtof,eloss,pid,id,trkid,lmom.theta(),lmom.phi(),processType);
00404             theDTHits.push_back(hit);
00405 
00406           }
00407         }
00408       }
00409       else if ( gd->subDetector() == GeomDetEnumerators::CSC ) {
00410         CSCDetId id(gd->geographicalId());
00411         const CSCChamber *chamber = cscGeom->chamber(id);
00412         std::vector<const CSCLayer *> layers = chamber->layers();
00413         for ( unsigned int ilayer=0; ilayer<layers.size(); ilayer++ ) {
00414           CSCDetId lid = layers[ilayer]->id();
00415 
00416 #ifdef FAMOS_DEBUG
00417             std::cout << "    Extrapolated to CSC (" 
00418                       << lid.endcap() << "," 
00419                       << lid.ring() << "," 
00420                       << lid.station() << "," 
00421                       << lid.layer() << ")" << std::endl;
00422 #endif
00423 
00424           const GeomDetUnit *det = cscGeom->idToDetUnit(lid);
00425           HelixArbitraryPlaneCrossing crossing(propagatedState.globalPosition().basicVector(),
00426                                                propagatedState.globalMomentum().basicVector(),
00427                                                propagatedState.transverseCurvature(),
00428                                                anyDirection);
00429           std::pair<bool,double> path = crossing.pathLength(det->surface());
00430           if ( ! path.first ) continue;
00431           LocalPoint lpos = det->toLocal(GlobalPoint(crossing.position(path.second)));
00432           // For CSCs the Bounds are for chamber frames not gas regions
00433           //      if ( ! det->surface().bounds().inside(lpos) ) continue;
00434           // New function knows where the 'active' volume is:
00435           const CSCLayerGeometry* laygeom = layers[ilayer]->geometry();
00436           if ( ! laygeom->inside( lpos ) ) continue; 
00437           //double thickness = laygeom->thickness(); gives number which is about 20 times too big
00438           double thickness = det->surface().bounds().thickness(); // this one works much better...
00439           LocalVector lmom = det->toLocal(GlobalVector(crossing.direction(path.second)));
00440           lmom = lmom.unit()*propagatedState.localMomentum().mag();
00441            
00442           // Factor that takes into account the (rec)hits lost because of delta's, etc.:
00443           // (Not fully satisfactory patch, but it seems to work...)
00444           double pmu = lmom.mag();
00445           double theCSCHitIneff = pmu>0? exp(kCSC*log(pmu)+fCSC):0.;
00446           // Take into account the different geometry in ME11:
00447           if (id.station()==1 && id.ring()==1)  theCSCHitIneff = theCSCHitIneff*0.442;
00448           if (random->flatShoot()<theCSCHitIneff) continue;
00449 
00450           double eloss = 0;
00451           double pz = fabs(lmom.z());
00452           LocalPoint entry = lpos - 0.5*thickness*lmom/pz;
00453           LocalPoint exit = lpos + 0.5*thickness*lmom/pz;
00454           double dtof = path.second*rbeta;
00455           int trkid = mySimTrack.trackId();
00456           unsigned int id = lid.rawId();
00457           short unsigned int processType = 2;
00458           PSimHit hit(entry,exit,lmom.mag(),
00459                       tof+dtof,eloss,pid,id,trkid,lmom.theta(),lmom.phi(),processType);
00460           theCSCHits.push_back(hit);
00461         }
00462       }
00463       else if ( gd->subDetector() == GeomDetEnumerators::RPCBarrel ||
00464                 gd->subDetector() == GeomDetEnumerators::RPCEndcap ) {
00465         RPCDetId id(gd->geographicalId());
00466         const RPCChamber *chamber = rpcGeom->chamber(id);
00467         std::vector<const RPCRoll *> roll = chamber->rolls();
00468         for ( unsigned int iroll=0; iroll<roll.size(); iroll++ ) {
00469           RPCDetId rid = roll[iroll]->id();
00470 
00471 #ifdef FAMOS_DEBUG
00472           std::cout << "    Extrapolated to RPC (" 
00473                     << rid.ring() << "," 
00474                     << rid.station() << ","
00475                     << rid.sector() << ","
00476                     << rid.subsector() << ","
00477                     << rid.layer() << ","
00478                     << rid.roll() << ")" << std::endl;
00479 #endif
00480 
00481           const GeomDetUnit *det = rpcGeom->idToDetUnit(rid);
00482           HelixArbitraryPlaneCrossing crossing(propagatedState.globalPosition().basicVector(),
00483                                                propagatedState.globalMomentum().basicVector(),
00484                                                propagatedState.transverseCurvature(),
00485                                                anyDirection);
00486           std::pair<bool,double> path = crossing.pathLength(det->surface());
00487           if ( ! path.first ) continue;
00488           LocalPoint lpos = det->toLocal(GlobalPoint(crossing.position(path.second)));
00489           if ( ! det->surface().bounds().inside(lpos) ) continue; 
00490           double thickness = det->surface().bounds().thickness();
00491           LocalVector lmom = det->toLocal(GlobalVector(crossing.direction(path.second)));
00492           lmom = lmom.unit()*propagatedState.localMomentum().mag();
00493           double eloss = 0;
00494           double pz = fabs(lmom.z());
00495           LocalPoint entry = lpos - 0.5*thickness*lmom/pz;
00496           LocalPoint exit = lpos + 0.5*thickness*lmom/pz;
00497           double dtof = path.second*rbeta;
00498           int trkid = mySimTrack.trackId();
00499           unsigned int id = rid.rawId();
00500           short unsigned int processType = 2;
00501           PSimHit hit(entry,exit,lmom.mag(),
00502                       tof+dtof,eloss,pid,id,trkid,lmom.theta(),lmom.phi(),processType);
00503           theRPCHits.push_back(hit);
00504         }
00505       }
00506       else {
00507         std::cout << "Extrapolated to unknown subdetector '" << gd->subDetector() << "'..." << std::endl;
00508       }
00509     }
00510   }
00511   }
00512 
00513   std::auto_ptr<edm::PSimHitContainer> pcsc(new edm::PSimHitContainer);
00514   int n = 0;
00515   for ( std::vector<PSimHit>::const_iterator i = theCSCHits.begin();
00516         i != theCSCHits.end(); i++ ) {
00517     pcsc->push_back(*i);
00518     n += 1;
00519   }
00520   iEvent.put(pcsc,"MuonCSCHits");
00521 
00522   std::auto_ptr<edm::PSimHitContainer> pdt(new edm::PSimHitContainer);
00523   n = 0;
00524   for ( std::vector<PSimHit>::const_iterator i = theDTHits.begin();
00525         i != theDTHits.end(); i++ ) {
00526     pdt->push_back(*i);
00527     n += 1;
00528   }
00529   iEvent.put(pdt,"MuonDTHits");
00530 
00531   std::auto_ptr<edm::PSimHitContainer> prpc(new edm::PSimHitContainer);
00532   n = 0;
00533   for ( std::vector<PSimHit>::const_iterator i = theRPCHits.begin();
00534         i != theRPCHits.end(); i++ ) {
00535     prpc->push_back(*i);
00536     n += 1;
00537   }
00538   iEvent.put(prpc,"MuonRPCHits");
00539 
00540 }
00541 
00542 
00543 // ------------ method called once each job just after ending the event loop  ------------
00544 void 
00545 MuonSimHitProducer::endJob() 
00546 {
00547 }
00548 
00549 
00550 void 
00551 MuonSimHitProducer::readParameters(const edm::ParameterSet& fastMuons, 
00552                                    const edm::ParameterSet& fastTracks,
00553                                    const edm::ParameterSet& matEff) {
00554   // Muons
00555   theSimModuleLabel_   = fastMuons.getParameter<std::string>("simModuleLabel");
00556   theSimModuleProcess_ = fastMuons.getParameter<std::string>("simModuleProcess");
00557   theTrkModuleLabel_   = fastMuons.getParameter<std::string>("trackModuleLabel");
00558   std::vector<double> simHitIneffDT  = fastMuons.getParameter<std::vector<double> >("simHitDTIneffParameters");
00559   std::vector<double> simHitIneffCSC = fastMuons.getParameter<std::vector<double> >("simHitCSCIneffParameters");
00560   kDT = simHitIneffDT[0];
00561   fDT = simHitIneffDT[1];
00562   kCSC = simHitIneffCSC[0];
00563   fCSC = simHitIneffCSC[1];  
00564 
00565   // Tracks
00566   fullPattern_  = fastTracks.getUntrackedParameter<bool>("FullPatternRecognition");
00567 
00568   // The following should be on LogInfo
00569   //  std::cout << " MUON SIM HITS: FastSimulation parameters " << std::endl;
00570   //  std::cout << " ============================================== " << std::endl;
00571   //  if ( fullPattern_ ) 
00572   //    std::cout << " The FULL pattern recognition option is turned ON" << std::endl;
00573   //  else
00574   //    std::cout << " The FAST tracking option is turned ON" << std::endl;
00575 
00576   // Material Effects
00577   theMaterialEffects = 0;
00578   if ( matEff.getParameter<bool>("PairProduction") || 
00579        matEff.getParameter<bool>("Bremsstrahlung") ||
00580        matEff.getParameter<bool>("MuonBremsstrahlung") ||
00581        matEff.getParameter<bool>("EnergyLoss") || 
00582        matEff.getParameter<bool>("MultipleScattering") )
00583     theMaterialEffects = new MaterialEffects(matEff,random);
00584 
00585 }
00586 
00587 void    
00588 MuonSimHitProducer::applyMaterialEffects(TrajectoryStateOnSurface& tsosWithdEdx,
00589                                          TrajectoryStateOnSurface& tsos,
00590                                          double radPath) { 
00591 
00592   // The energy loss simulator
00593   EnergyLossSimulator* energyLoss = theMaterialEffects->energyLossSimulator();
00594 
00595   // The multiple scattering simulator
00596   MultipleScatteringSimulator* multipleScattering = theMaterialEffects->multipleScatteringSimulator();
00597   
00598   // The Muon Bremsstrahlung simulator
00599   MuonBremsstrahlungSimulator* bremsstrahlung = theMaterialEffects->muonBremsstrahlungSimulator();
00600 
00601 
00602   // Initialize the Particle position, momentum and energy
00603   const Surface& nextSurface = tsos.surface();
00604   GlobalPoint gPos = energyLoss ? tsos.globalPosition() : tsosWithdEdx.globalPosition();
00605   GlobalVector gMom = energyLoss ? tsos.globalMomentum() : tsosWithdEdx.globalMomentum();  
00606   double mu = 0.1056583692;
00607   double en = std::sqrt(gMom.mag2()+mu*mu);
00608 
00609   // And now create the Particle
00610   XYZTLorentzVector position(gPos.x(),gPos.y(),gPos.z(),0.);
00611   XYZTLorentzVector momentum(gMom.x(),gMom.y(),gMom.z(),en);
00612   float charge = (float)(tsos.charge());
00613   ParticlePropagator theMuon(momentum,position,charge,0);
00614   theMuon.setID(-(int)charge*13);
00615 
00616   // Recompute the energy loss to get the fluctuations
00617   if ( energyLoss ) { 
00618     // Difference between with and without dE/dx (average only)
00619     // (for corrections once fluctuations are applied)
00620     GlobalPoint gPosWithdEdx = tsosWithdEdx.globalPosition();
00621     GlobalVector gMomWithdEdx = tsosWithdEdx.globalMomentum();
00622     double enWithdEdx = std::sqrt(gMomWithdEdx.mag2()+mu*mu);
00623     XYZTLorentzVector 
00624       deltaPos(gPosWithdEdx.x()-gPos.x(),gPosWithdEdx.y()-gPos.y(),
00625                gPosWithdEdx.z()-gPos.z(),0.);
00626     XYZTLorentzVector 
00627       deltaMom(gMomWithdEdx.x()-gMom.x(),gMomWithdEdx.y()-gMom.y(),
00628                gMomWithdEdx.z()-gMom.z(), enWithdEdx      -en);
00629 
00630     // Energy loss in iron (+ fluctuations)
00631     energyLoss->updateState(theMuon,radPath);
00632 
00633     // Correcting factors to account for slight differences in material descriptions
00634     // (Material description is more accurate in the stepping helix propagator)
00635     radPath *= -deltaMom.E()/energyLoss->mostLikelyLoss();
00636     double fac = energyLoss->deltaMom().E()/energyLoss->mostLikelyLoss();
00637 
00638     // Particle momentum & position after energy loss + fluctuation
00639     XYZTLorentzVector theNewMomentum = theMuon.momentum() + energyLoss->deltaMom() + fac * deltaMom;
00640     XYZTLorentzVector theNewPosition = theMuon.vertex() + fac * deltaPos;
00641     fac  = (theNewMomentum.E()*theNewMomentum.E()-mu*mu)/theNewMomentum.Vect().Mag2();
00642     fac  = fac>0.? std::sqrt(fac) : 1E-9;
00643     theMuon.SetXYZT(theNewMomentum.Px()*fac,theNewMomentum.Py()*fac,
00644                     theNewMomentum.Pz()*fac,theNewMomentum.E());    
00645     theMuon.setVertex(theNewPosition);
00646 
00647   }
00648 
00649   // Does the actual mutliple scattering
00650   if ( multipleScattering ) {
00651     // Pass the vector normal to the "next" surface 
00652     GlobalVector normal = nextSurface.tangentPlane(tsos.globalPosition())->normalVector();
00653     multipleScattering->setNormalVector(normal);
00654     // Compute the amount of multiple scattering after a given path length
00655     multipleScattering->updateState(theMuon,radPath);
00656   }
00657   
00658   // Muon Bremsstrahlung
00659   if ( bremsstrahlung ) {
00660     // Compute the amount of Muon Bremsstrahlung after given path length
00661     bremsstrahlung->updateState(theMuon,radPath);
00662   }
00663 
00664 
00665   // Fill the propagated state
00666   GlobalPoint propagatedPosition(theMuon.X(),theMuon.Y(),theMuon.Z());
00667   GlobalVector propagatedMomentum(theMuon.Px(),theMuon.Py(),theMuon.Pz());
00668   GlobalTrajectoryParameters propagatedGtp(propagatedPosition,propagatedMomentum,(int)charge,magfield);
00669   tsosWithdEdx = TrajectoryStateOnSurface(propagatedGtp,nextSurface);
00670 
00671 }
00672 
00673 
00674 
00675 //define this as a plug-in
00676 DEFINE_FWK_MODULE(MuonSimHitProducer);