LArSoft  v09_90_00
Liquid Argon Software toolkit - https://larsoft.org/
LArStackingAction Class Reference

#include "LArStackingAction.h"

Inheritance diagram for LArStackingAction:

Public Member Functions

 LArStackingAction (int)
 
virtual ~LArStackingAction ()
 
virtual G4ClassificationOfNewTrack ClassifyNewTrack (const G4Track *aTrack)
 
virtual void NewStage ()
 
virtual void PrepareNewEvent ()
 
void SetNRequestMuon (G4int val)
 
G4int GetNRequestMuon () const
 
void SetNRequestIsoMuon (G4int val)
 
G4int GetNRequestIsoMuon () const
 
void SetNIsolation (G4int val)
 
G4int GetNIsolation () const
 
void SetRoIAngle (G4double val)
 
G4double GetRoIAngle () const
 

Public Attributes

G4int fStack
 

Private Member Functions

std::string InsideTPC (const G4Track *aTrack)
 

Private Attributes

G4int fstage
 
G4int freqMuon
 
G4int freqIsoMuon
 
G4int freqIso
 
G4double fangRoI
 

Detailed Description

Definition at line 41 of file LArStackingAction.h.

Constructor & Destructor Documentation

LArStackingAction::LArStackingAction ( int  )

Definition at line 50 of file LArStackingAction.cxx.

References fStack.

51  : fstage(0), freqMuon(2), freqIsoMuon(0), freqIso(10), fangRoI(30. * CLHEP::deg)
52 {
53  //theMessenger = new LArStackingActionMessenger(this);
54  fStack = dum;
55  // Positive values effect action in this routine. Negative values
56  // effect action in G4BadIdeaAction.
57 }
LArStackingAction::~LArStackingAction ( )
virtual

Definition at line 59 of file LArStackingAction.cxx.

60 { //delete theMessenger;
61 }

Member Function Documentation

G4ClassificationOfNewTrack LArStackingAction::ClassifyNewTrack ( const G4Track *  aTrack)
virtual

Definition at line 63 of file LArStackingAction.cxx.

References geo::GeometryCore::DetHalfHeight(), geo::GeometryCore::DetHalfWidth(), geo::GeometryCore::DetLength(), fStack, fstage, geo::GeometryCore::GetLArTPCVolumeName(), InsideTPC(), and x4.

64 {
65  G4ClassificationOfNewTrack classification = fWaiting;
67  TString volName(InsideTPC(aTrack));
68  Double_t buffer = 500; // Keep muNucl neutrals within 5m (for now) of larTPC.
69 
70  // These 3 for now for investigation in gdb.
71  //int pdg = aTrack->GetDefinition()->GetPDGEncoding();
72  int ppdg = aTrack->GetParentID();
73  TString process("NA");
74  if (ppdg) process = (TString)aTrack->GetCreatorProcess()->GetProcessName();
75 
76  switch (fstage) {
77  case 0: // Fstage 0 : Primary muons only
78  if (aTrack->GetParentID() == 0) {
79  G4ParticleDefinition* particleType = aTrack->GetDefinition();
80  if (((particleType == G4MuonPlus::MuonPlusDefinition()) ||
81  (particleType == G4MuonMinus::MuonMinusDefinition())) &&
82  !volName.Contains("unknown")) {
83  classification = fUrgent;
84  }
85  }
86  if (volName.Contains("unknown")) classification = fKill;
87  break;
88 
89  case 1: // Stage 1 : K0,Lambda,n's made urgent here.
90  // Suspended tracks will be sent to the waiting stack
91  if (aTrack->GetTrackStatus() == fSuspend) { break; }
92 
93  if ((aTrack->GetDefinition()->GetPDGEncoding() == 2112 ||
94  aTrack->GetDefinition()->GetPDGEncoding() == 130 ||
95  aTrack->GetDefinition()->GetPDGEncoding() == 310 ||
96  aTrack->GetDefinition()->GetPDGEncoding() == 311 ||
97  aTrack->GetDefinition()->GetPDGEncoding() == 3122) &&
98  (aTrack->GetParentID() == 1) && !volName.Contains("unknown")) {
99 
100  const G4ThreeVector tr4Pos = aTrack->GetPosition();
101  // G4 returns positions in mm, have to convert to cm for LArSoft coordinate systems
102  const TVector3 trPos(tr4Pos.x() / CLHEP::cm, tr4Pos.y() / CLHEP::cm, tr4Pos.z() / CLHEP::cm);
103  //double locNeut = trPos.Mag();
104  classification = fUrgent;
105  // std::cout << "LArStackingAction: DetHalfWidth, Height, FullLength: " << geom->DetHalfWidth() << ", " << geom->DetHalfHeight() << ", " << geom->DetLength() << std::endl;
106 
107  if (trPos.X() < (geom->DetHalfWidth() * 2.0 + buffer) && trPos.X() > (-buffer) &&
108  trPos.Y() < (geom->DetHalfHeight() * 2.0 + buffer) &&
109  trPos.Y() > (-geom->DetHalfHeight() * 2.0 - buffer) &&
110  trPos.Z() < (geom->DetLength() + buffer) && trPos.Z() > (-buffer))
111 
112  {
113  classification = fUrgent;
114  break;
115  }
116  // These tracks need to be "scored" cuz every now and then they
117  // might get to the LAr.
118  else {
119  classification = fKill;
120  break;
121  }
122  }
123 
124  // if(aTrack->GetDefinition()->GetPDGCharge()==0.) { break; }
125  if (volName.Contains("unknown")) classification = fKill;
126  break;
127 
128  default:
129  // Track all other Primaries. Accept all secondaries in TPC.
130  // Kill muon ionization electrons outside TPC
131  // ignore primaries since they have no creator process
132 
133  if (aTrack->GetParentID() == 0 && !volName.Contains("unknown")) {
134  classification = fUrgent;
135  break;
136  }
137 
138  if (volName.Contains(geom->GetLArTPCVolumeName()) && aTrack->GetParentID() != 0) {
139  classification = fUrgent;
140  if (fStack & 0x4 && aTrack->GetCreatorProcess()->GetProcessName().contains("muIoni")) {
141  classification = fKill;
142  }
143  break;
144  }
145  else if (volName.Contains("unknown")) {
146  classification = fKill;
147  break;
148  }
149  // Leave this here, even though I claim we've Killed these in stage 2.
150  if (aTrack->GetDefinition()->GetPDGEncoding() == 11 &&
151  aTrack->GetCreatorProcess()->GetProcessName().contains("muIoni")) {
152  classification = fKill;
153  break;
154  }
155  // For now, kill every other thing, no matter where it is.
156  classification = fKill;
157 
158  } // end switch
159 
160  return classification;
161 }
std::string GetLArTPCVolumeName(TPCID const &tpcid=tpc_zero) const
Return the name of specified LAr TPC volume.
Length_t DetHalfWidth(TPCID const &tpcid=tpc_zero) const
Returns the half width of the active volume of the specified TPC.
std::string InsideTPC(const G4Track *aTrack)
Length_t DetLength(TPCID const &tpcid=tpc_zero) const
Returns the length of the active volume of the specified TPC.
Float_t x4[n_points_geant4]
Length_t DetHalfHeight(TPCID const &tpcid=tpc_zero) const
Returns the half height of the active volume of the specified TPC.
G4int LArStackingAction::GetNIsolation ( ) const
inline

Definition at line 74 of file LArStackingAction.h.

References freqIso.

74 { return freqIso; }
G4int LArStackingAction::GetNRequestIsoMuon ( ) const
inline

Definition at line 72 of file LArStackingAction.h.

References freqIsoMuon.

72 { return freqIsoMuon; }
G4int LArStackingAction::GetNRequestMuon ( ) const
inline

Definition at line 70 of file LArStackingAction.h.

References freqMuon.

70 { return freqMuon; }
G4double LArStackingAction::GetRoIAngle ( ) const
inline

Definition at line 76 of file LArStackingAction.h.

References fangRoI.

76 { return fangRoI; }
std::string LArStackingAction::InsideTPC ( const G4Track *  aTrack)
private

Definition at line 163 of file LArStackingAction.cxx.

References geo::GeometryCore::VolumeName().

Referenced by ClassifyNewTrack().

164 {
166  const G4ThreeVector tr4Pos = aTrack->GetPosition();
167 
168  // G4 returns positions in mm, have to convert to cm for LArSoft coordinate systems
169  const geo::Point_t trPos(tr4Pos.x() / CLHEP::cm, tr4Pos.y() / CLHEP::cm, tr4Pos.z() / CLHEP::cm);
170 
171  return geom->VolumeName(trPos);
172 }
std::string VolumeName(Point_t const &point) const
Returns the name of the deepest volume containing specified point.
ROOT::Math::PositionVector3D< ROOT::Math::Cartesian3D< double >, ROOT::Math::GlobalCoordinateSystemTag > Point_t
Type for representation of position in physical 3D space.
Definition: geo_vectors.h:180
void LArStackingAction::NewStage ( )
virtual

Definition at line 174 of file LArStackingAction.cxx.

References fstage.

175 {
176 
177  // Here when Urgent stack is empty. Waiting stack about to be made Urgent,
178  // upon saying ReClassify().
179  fstage++;
180 
181  // I yanked the ExN04's use here of stackManager->clear(), which clears stack
182  // and prepares to end the event. Think I may wanna do something like this if
183  // muon has been tracked, its doca is large, there are no hit voxels in the TPC,
184  // and I'm in further need of optimization.
185 
186  if (fstage == 1) {
187  // Stage 0->1 : check if at least "reqMuon" hits on muon chamber
188  // otherwise abort current event
189 
190  stackManager->ReClassify();
191  return;
192  }
193 
194  else if (fstage == 2) {
195  // Stage 1->2 : check the isolation of muon tracks
196  // at least "reqIsoMuon" isolated muons
197  // otherwise abort current event.
198  // Isolation requires "reqIso" or less hits
199  // (including own hits) in the RoI region
200  // in the tracker layers.
201  stackManager->ReClassify();
202  return;
203  }
204 
205  else {
206  // Other stage change : just re-classify
207  stackManager->ReClassify();
208  }
209 }
void LArStackingAction::PrepareNewEvent ( )
virtual

Definition at line 211 of file LArStackingAction.cxx.

References fstage.

212 {
213  fstage = 0;
214  //trkHits = 0;
215  //muonHits = 0;
216 }
void LArStackingAction::SetNIsolation ( G4int  val)
inline

Definition at line 73 of file LArStackingAction.h.

73 { freqIso = val; }
void LArStackingAction::SetNRequestIsoMuon ( G4int  val)
inline

Definition at line 71 of file LArStackingAction.h.

71 { freqIsoMuon = val; }
void LArStackingAction::SetNRequestMuon ( G4int  val)
inline

Definition at line 69 of file LArStackingAction.h.

69 { freqMuon = val; }
void LArStackingAction::SetRoIAngle ( G4double  val)
inline

Definition at line 75 of file LArStackingAction.h.

75 { fangRoI = val; }

Member Data Documentation

G4double LArStackingAction::fangRoI
private

Definition at line 66 of file LArStackingAction.h.

Referenced by GetRoIAngle().

G4int LArStackingAction::freqIso
private

Definition at line 65 of file LArStackingAction.h.

Referenced by GetNIsolation().

G4int LArStackingAction::freqIsoMuon
private

Definition at line 64 of file LArStackingAction.h.

Referenced by GetNRequestIsoMuon().

G4int LArStackingAction::freqMuon
private

Definition at line 63 of file LArStackingAction.h.

Referenced by GetNRequestMuon().

G4int LArStackingAction::fStack

Definition at line 51 of file LArStackingAction.h.

Referenced by ClassifyNewTrack(), and LArStackingAction().

G4int LArStackingAction::fstage
private

Definition at line 62 of file LArStackingAction.h.

Referenced by ClassifyNewTrack(), NewStage(), and PrepareNewEvent().


The documentation for this class was generated from the following files: