CMS 3D CMS Logo

Public Member Functions | Private Attributes

SaveSimTrack Class Reference

#include <SaveSimTrack.h>

Inheritance diagram for SaveSimTrack:
SimWatcher Observer< const BeginOfTrack * >

List of all members.

Public Member Functions

 SaveSimTrack (edm::ParameterSet const &p)
void update (const BeginOfTrack *trk)
 This routine will be called when the appropriate signal arrives.
 ~SaveSimTrack ()

Private Attributes

int pdgMax
int pdgMin

Detailed Description

Definition at line 9 of file SaveSimTrack.h.


Constructor & Destructor Documentation

SaveSimTrack::SaveSimTrack ( edm::ParameterSet const &  p)

Definition at line 10 of file SaveSimTrack.cc.

References edm::ParameterSet::getUntrackedParameter(), pdgMax, and pdgMin.

                                                    {

  pdgMin     = p.getUntrackedParameter<int>("MinimumPDGCode", 1000000);
  pdgMax     = p.getUntrackedParameter<int>("MaximumPDGCode", 2000000);

  edm::LogInfo("SaveSimTrack") << "SaveSimTrack:: Save Sim Track if PDG code "
                               << "lies between "  << pdgMin << " and " 
                               << pdgMax;
}
SaveSimTrack::~SaveSimTrack ( )

Definition at line 20 of file SaveSimTrack.cc.

{}

Member Function Documentation

void SaveSimTrack::update ( const BeginOfTrack ) [virtual]

This routine will be called when the appropriate signal arrives.

Implements Observer< const BeginOfTrack * >.

Definition at line 22 of file SaveSimTrack.cc.

References abs, LogDebug, pdgMax, pdgMin, and TrackInformation::storeTrack().

                                                  {

  G4Track* theTrack = (G4Track*)((*trk)());
  TrackInformation * trkInfo = (TrackInformation *)(theTrack->GetUserInformation());
  if (trkInfo) {
    int pdg = std::abs(theTrack->GetDefinition()->GetPDGEncoding());
    if (pdg >= pdgMin && pdg <= pdgMax) {
      trkInfo->storeTrack(true);
      LogDebug("SaveSimTrack") << "Save SimTrack the Track " 
                               << theTrack->GetTrackID() << " Type " 
                               << theTrack->GetDefinition()->GetParticleName()
                               << " Momentum " << theTrack->GetMomentum()/MeV 
                               << " MeV/c";
    }
  }
}

Member Data Documentation

int SaveSimTrack::pdgMax [private]

Definition at line 18 of file SaveSimTrack.h.

Referenced by SaveSimTrack(), and update().

int SaveSimTrack::pdgMin [private]

Definition at line 18 of file SaveSimTrack.h.

Referenced by SaveSimTrack(), and update().