CMS 3D CMS Logo

CMSSQLoopProcess.cc
Go to the documentation of this file.
1 
3 #include "G4SystemOfUnits.hh"
4 #include "G4Step.hh"
5 #include "G4ParticleDefinition.hh"
6 #include "G4VParticleChange.hh"
7 
9 
10 CMSSQLoopProcess::CMSSQLoopProcess(const G4String& name, G4ProcessType type) : G4VContinuousProcess(name, type) {
11  fParticleChange = new G4ParticleChange();
12 }
13 
15 
16 G4VParticleChange* CMSSQLoopProcess::AlongStepDoIt(const G4Track& track, const G4Step& step) {
17  if (track.GetPosition() == posini)
18  edm::LogInfo("CMSSQLoopProcess::AlongStepDoIt")
19  << "CMSSQLoopProcess::AlongStepDoIt: CMSSQLoopProcess::AlongStepDoIt MomentumDirection "
20  << track.GetMomentumDirection().eta() << " track GetPostion " << track.GetPosition() / cm << " trackId "
21  << track.GetTrackID() << " parentId: " << track.GetParentID() << " GlobalTime " << track.GetGlobalTime() / ns
22  << " TotalEnergy: " << track.GetTotalEnergy() / GeV << " Velocity " << track.GetVelocity() / m / ns
23  << std::endl;
24 
25  fParticleChange->Clear();
26  fParticleChange->Initialize(track);
27  fParticleChange->ProposeWeight(track.GetWeight());
28  //Sbar not passing the following criteria are not of interest. They will not be reconstructable. A cut like this is required otherwise you will get Sbar infinitely looping.
29  if (fabs(track.GetMomentumDirection().eta()) > 999. || fabs(track.GetPosition().z()) > 160 * centimeter) {
30  edm::LogInfo("CMSSQLoopProcess::AlongStepDoIt") << "Particle getting killed because too large z" << std::endl;
31  fParticleChange->ProposeTrackStatus(fStopAndKill);
32  }
33 
34  return fParticleChange;
35 }
36 
38  G4double previousStepSize,
39  G4double currentMinimumStep,
40  G4double& proposedSafety,
41  G4GPILSelection* selection) {
42  return 1. * centimeter;
43 }
44 
45 G4double CMSSQLoopProcess::GetContinuousStepLimit(const G4Track& track, G4double, G4double, G4double&) {
46  return 1. * centimeter; // seems irrelevant
47 }
48 
49 void CMSSQLoopProcess::StartTracking(G4Track* aTrack) { posini = aTrack->GetPosition(); }
type
Definition: HCALResponse.h:21
const double GeV
Definition: MathUtil.h:16
virtual G4double GetContinuousStepLimit(const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety)
selection
main part
Definition: corrVsCorr.py:100
virtual G4double AlongStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)
CMSSQLoopProcess(const G4String &name="SQLooper", G4ProcessType type=fUserDefined)
virtual ~CMSSQLoopProcess()
virtual void StartTracking(G4Track *aTrack)
virtual G4VParticleChange * AlongStepDoIt(const G4Track &, const G4Step &)
G4ThreeVector posini
G4ParticleChange * fParticleChange
step
Definition: StallMonitor.cc:94