CMS 3D CMS Logo

TrackingElectronProducer Class Reference

Produces one composite TrackingParticle per electron, for electron tracking studies: More...

#include <SimGeneral/TrackingAnalysis/interface/TrackingElectronProducer.h>

Inheritance diagram for TrackingElectronProducer:

edm::EDProducer edm::ProducerBase edm::ProductRegistryHelper

List of all members.

Public Member Functions

 TrackingElectronProducer (const edm::ParameterSet &)

Private Member Functions

void addG4Track (TrackingParticle &, const TrackingParticle *) const
int layerFromDetid (const unsigned int &)
void listElectrons (const TrackingParticleCollection &tPC) const
void produce (edm::Event &, const edm::EventSetup &)

Private Attributes

edm::ParameterSet conf_


Detailed Description

Produces one composite TrackingParticle per electron, for electron tracking studies:

Definition at line 19 of file TrackingElectronProducer.h.


Constructor & Destructor Documentation

TrackingElectronProducer::TrackingElectronProducer ( const edm::ParameterSet conf  )  [explicit]

Definition at line 35 of file TrackingElectronProducer.cc.

References conf_.

00035                                                                               {
00036   produces<TrackingParticleCollection>("ElectronTrackTruth");
00037 
00038 //  std::cout << " TrackingElectronProducer CTOR " << std::endl;
00039 
00040 
00041   conf_ = conf;
00042 }


Member Function Documentation

void TrackingElectronProducer::addG4Track ( TrackingParticle e,
const TrackingParticle s 
) const [private]

Definition at line 171 of file TrackingElectronProducer.cc.

References TrackingParticle::addG4Track(), and TrackingParticle::addPSimHit().

Referenced by produce().

00173 {
00174 
00175   // add G4 tracks
00176   std::vector<SimTrack> g4Tracks = (*s).g4Tracks();
00177   for ( std::vector<SimTrack>::const_iterator ig4 = g4Tracks.begin();
00178         ig4 != g4Tracks.end(); ig4++ ) {
00179     e.addG4Track(*ig4);
00180   }
00181 
00182   // add hits
00183   // FIXME configurable for dropping delta-ray hits
00184   std::vector< PSimHit > hits = (*s).trackPSimHit();
00185   for ( std::vector<PSimHit>::const_iterator ih = hits.begin();
00186         ih != hits.end(); ih++ ) {
00187     e.addPSimHit(*ih);
00188   }
00189 }

int TrackingElectronProducer::layerFromDetid ( const unsigned int detid  )  [private]

Definition at line 192 of file TrackingElectronProducer.cc.

References detId, align::tpb::layerNumber(), PixelSubdetector::PixelBarrel, PixelSubdetector::PixelEndcap, DetId::rawId(), DetId::subdetId(), StripSubdetector::TEC, StripSubdetector::TIB, StripSubdetector::TID, and StripSubdetector::TOB.

00192                                                                        {
00193   DetId detId = DetId(detid);
00194   int layerNumber=0;
00195   unsigned int subdetId = static_cast<unsigned int>(detId.subdetId());
00196   if ( subdetId == StripSubdetector::TIB)
00197     {
00198       TIBDetId tibid(detId.rawId());
00199       layerNumber = tibid.layer();
00200     }
00201   else if ( subdetId ==  StripSubdetector::TOB )
00202     {
00203       TOBDetId tobid(detId.rawId());
00204       layerNumber = tobid.layer();
00205     }
00206   else if ( subdetId ==  StripSubdetector::TID)
00207     {
00208       TIDDetId tidid(detId.rawId());
00209       layerNumber = tidid.wheel();
00210     }
00211   else if ( subdetId ==  StripSubdetector::TEC )
00212     {
00213       TECDetId tecid(detId.rawId());
00214       layerNumber = tecid.wheel();
00215     }
00216   else if ( subdetId ==  PixelSubdetector::PixelBarrel )
00217     {
00218       PXBDetId pxbid(detId.rawId());
00219       layerNumber = pxbid.layer();
00220     }
00221   else if ( subdetId ==  PixelSubdetector::PixelEndcap )
00222     {
00223       PXFDetId pxfid(detId.rawId());
00224       layerNumber = pxfid.disk();
00225     }
00226   else
00227     edm::LogVerbatim("TrackingTruthProducer") << "Unknown subdetid: " <<  subdetId;
00228 
00229   return layerNumber;
00230 }

void TrackingElectronProducer::listElectrons ( const TrackingParticleCollection tPC  )  const [private]

Definition at line 140 of file TrackingElectronProducer.cc.

References funct::abs(), and it.

00142 {
00143 //  cout << "TrackingElectronProducer::printing electrons before assembly..."
00144 //       << endl;
00145   for (TrackingParticleCollection::const_iterator it = tPC.begin();
00146        it != tPC.end(); it++) {
00147     if (abs((*it).pdgId()) == 11) {
00148 //      cout << "Electron: sim tk " << (*it).g4Tracks().front().trackId()
00149 //         << endl;
00150 
00151       TrackingVertexRef parentV = (*it).parentVertex();
00152 //      if (parentV.isNull()) {
00153 //      cout << " No parent vertex" << endl;
00154 //      } else {
00155 //      cout << " Parent  vtx position " << parentV -> position() << endl;
00156 //      }
00157 
00158       TrackingVertexRefVector decayVertices = (*it).decayVertices();
00159 //      if ( decayVertices.empty() ) {
00160 //      cout << " No decay vertex" << endl;
00161 //      } else {
00162 //      cout << " Decay vtx position "
00163 //           << decayVertices.at(0) -> position() << endl;
00164 //      }
00165     }
00166   }
00167 }

void TrackingElectronProducer::produce ( edm::Event event,
const edm::EventSetup  
) [private, virtual]

Implements edm::EDProducer.

Definition at line 45 of file TrackingElectronProducer.cc.

References addG4Track(), TkNavigableSimElectronAssembler::assemble(), edm::RefVector< C, T, F >::at(), edm::RefVector< C, T, F >::empty(), edm::Ref< C, T, F >::isNonnull(), it, t, and trackingParticles_cfi::trackingParticles.

00045                                                                        {
00046 
00047   // Get information out of event record
00048 
00049   edm::Handle<TrackingParticleCollection>  TruthTrackContainer ;
00050   event.getByType(TruthTrackContainer );
00051 
00052   TkNavigableSimElectronAssembler assembler;
00053   std::vector<TrackingParticle*> particles;
00054 
00055   for (TrackingParticleCollection::const_iterator t = TruthTrackContainer->begin();
00056        t != TruthTrackContainer->end(); ++t) {
00057     particles.push_back(const_cast<TrackingParticle*>(&(*t)));
00058   }
00059 
00060   TkNavigableSimElectronAssembler::ElectronList electrons = assembler.assemble(particles);
00061 
00062   //
00063   // now add electron tracks and vertices to event
00064   //
00065   // prepare collections
00066   //
00067   std::auto_ptr<TrackingParticleCollection> trackingParticles(new TrackingParticleCollection);
00068 
00069   edm::RefProd<TrackingParticleCollection> refTPC =
00070     event.getRefBeforePut<TrackingParticleCollection>("ElectronTrackTruth");
00071 
00072   //
00073   // create TrackingParticles
00074   //
00075 
00076   int totsimhit = 0; //fixing the electron hits
00077  // loop over electrons
00078   for ( TkNavigableSimElectronAssembler::ElectronList::const_iterator ie
00079           = electrons.begin(); ie != electrons.end(); ie++ ) {
00080 
00081     // create TrackingParticle from first track segment
00082     TrackingParticle * tk = (*ie).front();
00083     LorentzVector hep4Pos = (*(*tk).parentVertex()).position();
00084     TrackingParticle::Point pos(hep4Pos.x(), hep4Pos.y(), hep4Pos.z());
00085     TrackingParticle tkp( (*tk).charge(), (*tk).p4(), pos, hep4Pos.t(),
00086                           (*tk).pdgId(), (*tk).status(), (*tk).eventId() );
00087 
00088     // add G4 tracks and hits of all segments
00089     int ngenp = 0;
00090     totsimhit = 0;//initialize the number of matchedHits for each track
00091     for (TkNavigableSimElectronAssembler::TrackList::const_iterator it = (*ie).begin();
00092          it != (*ie).end(); it++ ) {
00093 
00094       for (TrackingParticle::genp_iterator uz=(*it)->genParticle_begin();
00095            uz!=(*it)->genParticle_end();uz++) {
00096         tkp.addGenParticle(*uz);
00097         ngenp++;
00098       }
00099       addG4Track(tkp, *it);
00100       totsimhit +=(*tk).matchedHit();
00101       
00102       /*
00103         std::cout << "Electron list of tracks Original Segment  = " << (*tk) 
00104         << "\t Matched Hits = " << (*tk).matchedHit() 
00105         << "\t SimHits = " << (*tk).trackPSimHit().size() 
00106         << std::endl;
00107       */
00108       
00109     }
00110 
00111     //    int totsimhit = 20; // FIXME temp. hack
00112     tkp.setMatchedHit(totsimhit);
00113 
00114     //
00115     // add references to parent and decay vertices
00116     //
00117     TrackingVertexRef parentV = (*(*ie).front()).parentVertex();
00118     if (parentV.isNonnull()) {
00119       tkp.setParentVertex(parentV);
00120     }
00121 
00122     TrackingVertexRefVector decayVertices = (*(*ie).back()).decayVertices();
00123     if ( !decayVertices.empty() ) {
00124       // get first decay vertex
00125       TrackingVertexRef decayV = decayVertices.at(0);
00126       tkp.addDecayVertex(decayV);
00127     }
00128 
00129     //
00130     // put particle in transient store
00131     //
00132     (*trackingParticles).push_back(tkp);
00133   }
00134 
00135   event.put(trackingParticles,"ElectronTrackTruth");
00136 
00137 }


Member Data Documentation

edm::ParameterSet TrackingElectronProducer::conf_ [private]

Definition at line 36 of file TrackingElectronProducer.h.

Referenced by TrackingElectronProducer().


The documentation for this class was generated from the following files:
Generated on Tue Jun 9 18:34:12 2009 for CMSSW by  doxygen 1.5.4