Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions
WLSTrajectory Class Reference

#include <WLSTrajectory.hh>

Inheritance diagram for WLSTrajectory:
G4VTrajectory

Public Member Functions

 WLSTrajectory ()
 
 WLSTrajectory (const G4Track *)
 
 WLSTrajectory (WLSTrajectory &)
 
virtual ~WLSTrajectory ()
 
voidoperator new (size_t)
 
void operator delete (void *)
 
int operator== (const WLSTrajectory &right) const
 
virtual G4int GetTrackID () const
 
virtual G4int GetParentID () const
 
virtual G4String GetParticleName () const
 
virtual G4double GetCharge () const
 
virtual G4int GetPDGEncoding () const
 
virtual G4ThreeVector GetInitialMomentum () const
 
virtual void ShowTrajectory (std::ostream &os=G4cout) const
 
virtual void DrawTrajectory () const
 
virtual void AppendStep (const G4Step *aStep)
 
virtual void MergeTrajectory (G4VTrajectory *secondTrajectory)
 
G4ParticleDefinitionGetParticleDefinition ()
 
virtual int GetPointEntries () const
 
virtual G4VTrajectoryPointGetPoint (G4int i) const
 
virtual const std::map
< G4String, G4AttDef > * 
GetAttDefs () const
 
virtual std::vector< G4AttValue > * CreateAttValues () const
 
void SetDrawTrajectory (G4bool b)
 
void WLS ()
 
void SetForceDrawTrajectory (G4bool b)
 
void SetForceNoDrawTrajectory (G4bool b)
 
- Public Member Functions inherited from G4VTrajectory
 G4VTrajectory ()
 
virtual ~G4VTrajectory ()
 
G4bool operator== (const G4VTrajectory &right) const
 

Detailed Description

Definition at line 53 of file WLSTrajectory.hh.

Constructor & Destructor Documentation

WLSTrajectory::WLSTrajectory ( )

Definition at line 59 of file WLSTrajectory.cc.

60  : fpPointsContainer(0), fTrackID(0), fParentID(0),
61  fPDGCharge(0.0), fPDGEncoding(0), fParticleName(""),
62  fInitialMomentum(G4ThreeVector())
63 {
64  fWLS = false;
65  fDrawIt = false;
66  fForceNoDraw = false;
67  fForceDraw = false;
68 
69  fParticleDefinition = NULL;
70 }
CLHEP::Hep3Vector G4ThreeVector
WLSTrajectory::WLSTrajectory ( const G4Track aTrack)

Definition at line 74 of file WLSTrajectory.cc.

References G4Track::GetDefinition(), G4Track::GetMomentum(), G4Track::GetParentID(), G4ParticleDefinition::GetParticleName(), G4ParticleDefinition::GetPDGCharge(), G4ParticleDefinition::GetPDGEncoding(), and G4Track::GetTrackID().

75 {
76  fParticleDefinition = aTrack->GetDefinition();
77  fParticleName = fParticleDefinition->GetParticleName();
78  fPDGCharge = fParticleDefinition->GetPDGCharge();
79  fPDGEncoding = fParticleDefinition->GetPDGEncoding();
80  fTrackID = aTrack->GetTrackID();
81  fParentID = aTrack->GetParentID();
82  fInitialMomentum = aTrack->GetMomentum();
83  fpPointsContainer = new WLSTrajectoryPointContainer();
84  // Following is for the first trajectory point
85  fpPointsContainer->push_back(new WLSTrajectoryPoint(aTrack));
86 
87  fWLS = false;
88  fDrawIt = false;
89  fForceNoDraw = false;
90  fForceDraw = false;
91 }
G4ParticleDefinition * GetDefinition() const
G4int GetParentID() const
const G4String & GetParticleName() const
G4int GetTrackID() const
G4ThreeVector GetMomentum() const
std::vector< G4VTrajectoryPoint * > WLSTrajectoryPointContainer
G4double GetPDGCharge() const
WLSTrajectory::WLSTrajectory ( WLSTrajectory right)

Definition at line 95 of file WLSTrajectory.cc.

95  : G4VTrajectory()
96 {
97  fParticleDefinition=right.fParticleDefinition;
98  fParticleName = right.fParticleName;
99  fPDGCharge = right.fPDGCharge;
100  fPDGEncoding = right.fPDGEncoding;
101  fTrackID = right.fTrackID;
102  fParentID = right.fParentID;
103  fInitialMomentum = right.fInitialMomentum;
104  fpPointsContainer = new WLSTrajectoryPointContainer();
105 
106  for(size_t i=0;i<right.fpPointsContainer->size();++i) {
107  WLSTrajectoryPoint* rightPoint
108  = (WLSTrajectoryPoint*)((*(right.fpPointsContainer))[i]);
109  fpPointsContainer->push_back(new WLSTrajectoryPoint(*rightPoint));
110  }
111 
112  fWLS = right.fWLS;
113  fDrawIt = right.fDrawIt;
114  fForceNoDraw = right.fForceNoDraw;
115  fForceDraw = right.fForceDraw;
116 }
std::vector< G4VTrajectoryPoint * > WLSTrajectoryPointContainer
WLSTrajectory::~WLSTrajectory ( )
virtual

Definition at line 120 of file WLSTrajectory.cc.

121 {
122  for(size_t i=0;i<fpPointsContainer->size();++i){
123  delete (*fpPointsContainer)[i];
124  }
125  fpPointsContainer->clear();
126 
127  delete fpPointsContainer;
128 }

Member Function Documentation

void WLSTrajectory::AppendStep ( const G4Step aStep)
virtual

Implements G4VTrajectory.

Definition at line 234 of file WLSTrajectory.cc.

235 {
236  fpPointsContainer->push_back(new WLSTrajectoryPoint(aStep));
237 }
std::vector< G4AttValue > * WLSTrajectory::CreateAttValues ( ) const
virtual

Reimplemented from G4VTrajectory.

Definition at line 306 of file WLSTrajectory.cc.

References G4UIcommand::ConvertToString(), G4BestUnit, G4cout, GetAttDefs(), GetPointEntries(), and CLHEP::Hep3Vector::mag().

307 {
308  std::vector<G4AttValue>* values = new std::vector<G4AttValue>;
309 
310  values->push_back
311  (G4AttValue("ID",G4UIcommand::ConvertToString(fTrackID),""));
312 
313  values->push_back
314  (G4AttValue("PID",G4UIcommand::ConvertToString(fParentID),""));
315 
316  values->push_back(G4AttValue("PN",fParticleName,""));
317 
318  values->push_back
319  (G4AttValue("Ch",G4UIcommand::ConvertToString(fPDGCharge),""));
320 
321  values->push_back
322  (G4AttValue("PDG",G4UIcommand::ConvertToString(fPDGEncoding),""));
323 
324  values->push_back
325  (G4AttValue("IMom",G4BestUnit(fInitialMomentum,"Energy"),""));
326 
327  values->push_back
328  (G4AttValue("IMag",G4BestUnit(fInitialMomentum.mag(),"Energy"),""));
329 
330  values->push_back
332 
333 #ifdef G4ATTDEBUG
334  G4cout << G4AttCheck(values,GetAttDefs());
335 #endif
336  return values;
337 }
virtual const std::map< G4String, G4AttDef > * GetAttDefs() const
static G4String ConvertToString(G4bool boolVal)
Definition: G4UIcommand.cc:357
#define G4BestUnit(a, b)
#define G4_USE_G4BESTUNIT_FOR_VERBOSE 1
virtual int GetPointEntries() const
G4GLOB_DLL std::ostream G4cout
double mag() const
void WLSTrajectory::DrawTrajectory ( ) const
virtual

Reimplemented from G4VTrajectory.

Definition at line 141 of file WLSTrajectory.cc.

References G4Polymarker::circles, G4VVisManager::Draw(), G4VMarker::filled, G4VTrajectoryPoint::GetAuxiliaryPoints(), G4VVisManager::GetConcreteInstance(), GetPoint(), GetPointEntries(), G4VTrajectoryPoint::GetPosition(), G4OpticalPhoton::OpticalPhotonDefinition(), G4VMarker::SetFillStyle(), G4Polymarker::SetMarkerType(), G4VMarker::SetScreenSize(), G4Visible::SetVisAttributes(), and G4Polymarker::squares.

Referenced by WLSEventAction::EndOfEventAction().

142 {
143  // i_mode is no longer available as an argument of G4VTrajectory.
144  // In this exampple it was always called with an argument of 50.
145  const G4int i_mode = 50;
146  // Consider using commands /vis/modeling/trajectories.
147 
148  // Invoke the default implementation in G4VTrajectory...
149  // G4VTrajectory::DrawTrajectory(i_mode);
150  // ... or override with your own code here.
151 
152  //Taken from G4VTrajectory and modified to select colours based on particle
153  //type and to selectively eliminate drawing of certain trajectories.
154 
155  if (!fForceDraw && (!fDrawIt || fForceNoDraw)) return;
156 
157  // If i_mode>=0, draws a trajectory as a polyline and, if i_mode!=0,
158  // adds markers - yellow circles for step points and magenta squares
159  // for auxiliary points, if any - whose screen size in pixels is
160  // given by std::abs(i_mode)/1000. E.g: i_mode = 5000 gives easily
161  // visible markers.
162 
164  if (!pVVisManager) return;
165 
166  const G4double markerSize = std::abs(i_mode)/1000;
167  G4bool lineRequired (i_mode >= 0);
168  G4bool markersRequired (markerSize > 0.);
169 
170  G4Polyline trajectoryLine;
171  G4Polymarker stepPoints;
172  G4Polymarker auxiliaryPoints;
173 
174  for (G4int i = 0; i < GetPointEntries() ; i++) {
175  G4VTrajectoryPoint* aTrajectoryPoint = GetPoint(i);
176  const std::vector<G4ThreeVector>* auxiliaries
177  = aTrajectoryPoint->GetAuxiliaryPoints();
178  if (auxiliaries) {
179  for (size_t iAux = 0; iAux < auxiliaries->size(); ++iAux) {
180  const G4ThreeVector pos((*auxiliaries)[iAux]);
181  if (lineRequired) {
182  trajectoryLine.push_back(pos);
183  }
184  if (markersRequired) {
185  auxiliaryPoints.push_back(pos);
186  }
187  }
188  }
189  const G4ThreeVector pos(aTrajectoryPoint->GetPosition());
190  if (lineRequired) {
191  trajectoryLine.push_back(pos);
192  }
193  if (markersRequired) {
194  stepPoints.push_back(pos);
195  }
196  }
197 
198  if (lineRequired) {
199  G4Colour colour;
200 
201  if(fParticleDefinition==G4OpticalPhoton::OpticalPhotonDefinition()){
202  if(fWLS) //WLS photons are red
203  colour = G4Colour(1.,0.,0.);
204  else{ //Scintillation and Cerenkov photons are green
205  colour = G4Colour(0.,1.,0.);
206  }
207  }
208  else //All other particles are blue
209  colour = G4Colour(0.,0.,1.);
210 
211  G4VisAttributes trajectoryLineAttribs(colour);
212  trajectoryLine.SetVisAttributes(&trajectoryLineAttribs);
213  pVVisManager->Draw(trajectoryLine);
214  }
215  if (markersRequired) {
216  auxiliaryPoints.SetMarkerType(G4Polymarker::squares);
217  auxiliaryPoints.SetScreenSize(markerSize);
218  auxiliaryPoints.SetFillStyle(G4VMarker::filled);
219  G4VisAttributes auxiliaryPointsAttribs(G4Colour(0.,1.,1.)); // Magenta
220  auxiliaryPoints.SetVisAttributes(&auxiliaryPointsAttribs);
221  pVVisManager->Draw(auxiliaryPoints);
222 
224  stepPoints.SetScreenSize(markerSize);
225  stepPoints.SetFillStyle(G4VMarker::filled);
226  G4VisAttributes stepPointsAttribs(G4Colour(1.,1.,0.)); // Yellow.
227  stepPoints.SetVisAttributes(&stepPointsAttribs);
228  pVVisManager->Draw(stepPoints);
229  }
230 }
virtual void Draw(const G4Circle &, const G4Transform3D &objectTransformation=G4Transform3D())=0
virtual const std::vector< G4ThreeVector > * GetAuxiliaryPoints() const
static G4VVisManager * GetConcreteInstance()
virtual G4VTrajectoryPoint * GetPoint(G4int i) const
void SetMarkerType(MarkerType)
void SetFillStyle(FillStyle)
int G4int
Definition: G4Types.hh:78
virtual int GetPointEntries() const
bool G4bool
Definition: G4Types.hh:79
void SetVisAttributes(const G4VisAttributes *)
Definition: G4Visible.cc:80
virtual const G4ThreeVector GetPosition() const =0
static G4OpticalPhoton * OpticalPhotonDefinition()
double G4double
Definition: G4Types.hh:76
void SetScreenSize(G4double)
const std::map< G4String, G4AttDef > * WLSTrajectory::GetAttDefs ( ) const
virtual

Reimplemented from G4VTrajectory.

Definition at line 264 of file WLSTrajectory.cc.

References G4AttDefStore::GetInstance().

Referenced by CreateAttValues().

265 {
266  G4bool isNew;
267  std::map<G4String,G4AttDef>* store
268  = G4AttDefStore::GetInstance("Trajectory",isNew);
269 
270  if (isNew) {
271 
272  G4String ID("ID");
273  (*store)[ID] = G4AttDef(ID,"Track ID","Bookkeeping","","G4int");
274 
275  G4String PID("PID");
276  (*store)[PID] = G4AttDef(PID,"Parent ID","Bookkeeping","","G4int");
277 
278  G4String PN("PN");
279  (*store)[PN] = G4AttDef(PN,"Particle Name","Physics","","G4String");
280 
281  G4String Ch("Ch");
282  (*store)[Ch] = G4AttDef(Ch,"Charge","Physics","e+","G4double");
283 
284  G4String PDG("PDG");
285  (*store)[PDG] = G4AttDef(PDG,"PDG Encoding","Physics","","G4int");
286 
287  G4String IMom("IMom");
288  (*store)[IMom] = G4AttDef(IMom,
289  "Momentum of track at start of trajectory",
290  "Physics","G4BestUnit","G4ThreeVector");
291 
292  G4String IMag("IMag");
293  (*store)[IMag] = G4AttDef(IMag,
294  "Magnitude of momentum of track at start of trajectory",
295  "Physics","G4BestUnit","G4double");
296 
297  G4String NTP("NTP");
298  (*store)[NTP] = G4AttDef(NTP,"No. of points","Bookkeeping","","G4int");
299 
300  }
301  return store;
302 }
bool G4bool
Definition: G4Types.hh:79
std::map< G4String, G4AttDef > * GetInstance(const G4String &storeKey, G4bool &isNew)
virtual G4double WLSTrajectory::GetCharge ( ) const
inlinevirtual

Implements G4VTrajectory.

Definition at line 79 of file WLSTrajectory.hh.

Referenced by WLSEventAction::EndOfEventAction().

79 { return fPDGCharge; }
virtual G4ThreeVector WLSTrajectory::GetInitialMomentum ( ) const
inlinevirtual

Implements G4VTrajectory.

Definition at line 81 of file WLSTrajectory.hh.

82  {return fInitialMomentum;}
virtual G4int WLSTrajectory::GetParentID ( ) const
inlinevirtual

Implements G4VTrajectory.

Definition at line 77 of file WLSTrajectory.hh.

77 { return fParentID; }
G4ParticleDefinition * WLSTrajectory::GetParticleDefinition ( )

Definition at line 241 of file WLSTrajectory.cc.

References G4ParticleTable::GetParticleTable().

242 {
243  return (G4ParticleTable::GetParticleTable()->FindParticle(fParticleName));
244 }
static G4ParticleTable * GetParticleTable()
virtual G4String WLSTrajectory::GetParticleName ( ) const
inlinevirtual

Implements G4VTrajectory.

Definition at line 78 of file WLSTrajectory.hh.

Referenced by WLSEventAction::EndOfEventAction().

78 { return fParticleName; }
virtual G4int WLSTrajectory::GetPDGEncoding ( ) const
inlinevirtual

Implements G4VTrajectory.

Definition at line 80 of file WLSTrajectory.hh.

80 { return fPDGEncoding; }
virtual G4VTrajectoryPoint* WLSTrajectory::GetPoint ( G4int  i) const
inlinevirtual

Implements G4VTrajectory.

Definition at line 95 of file WLSTrajectory.hh.

Referenced by DrawTrajectory().

96  { return (*fpPointsContainer)[i]; }
virtual int WLSTrajectory::GetPointEntries ( ) const
inlinevirtual

Implements G4VTrajectory.

Definition at line 93 of file WLSTrajectory.hh.

Referenced by CreateAttValues(), DrawTrajectory(), and MergeTrajectory().

94  { return fpPointsContainer->size(); }
virtual G4int WLSTrajectory::GetTrackID ( ) const
inlinevirtual

Implements G4VTrajectory.

Definition at line 76 of file WLSTrajectory.hh.

76 { return fTrackID; }
void WLSTrajectory::MergeTrajectory ( G4VTrajectory secondTrajectory)
virtual

Implements G4VTrajectory.

Definition at line 248 of file WLSTrajectory.cc.

References GetPointEntries(), and python.hepunit::second.

249 {
250  if(!secondTrajectory) return;
251 
252  WLSTrajectory* second = (WLSTrajectory*)secondTrajectory;
253  G4int ent = second->GetPointEntries();
254  // initial point of the second trajectory should not be merged
255  for(G4int i=1; i<ent; ++i) {
256  fpPointsContainer->push_back((*(second->fpPointsContainer))[i]);
257  }
258  delete (*second->fpPointsContainer)[0];
259  second->fpPointsContainer->clear();
260 }
int G4int
Definition: G4Types.hh:78
virtual int GetPointEntries() const
void WLSTrajectory::operator delete ( void aTrajectory)
inline

Definition at line 138 of file WLSTrajectory.hh.

References WLSTrajectoryAllocator.

138  {
139  WLSTrajectoryAllocator->FreeSingle((WLSTrajectory*)aTrajectory);
140 }
G4ThreadLocal G4Allocator< WLSTrajectory > * WLSTrajectoryAllocator
void * WLSTrajectory::operator new ( size_t  )
inline

Definition at line 132 of file WLSTrajectory.hh.

References WLSTrajectoryAllocator.

132  {
135  return (void*) WLSTrajectoryAllocator->MallocSingle();
136 }
G4ThreadLocal G4Allocator< WLSTrajectory > * WLSTrajectoryAllocator
int WLSTrajectory::operator== ( const WLSTrajectory right) const
inline

Definition at line 71 of file WLSTrajectory.hh.

72  { return (this==&right); }
void WLSTrajectory::SetDrawTrajectory ( G4bool  b)
inline

Definition at line 101 of file WLSTrajectory.hh.

References test::b.

Referenced by WLSTrackingAction::PostUserTrackingAction().

101 {fDrawIt=b;}
void WLSTrajectory::SetForceDrawTrajectory ( G4bool  b)
inline

Definition at line 103 of file WLSTrajectory.hh.

References test::b.

Referenced by WLSEventAction::EndOfEventAction().

103 {fForceDraw=b;}
void WLSTrajectory::SetForceNoDrawTrajectory ( G4bool  b)
inline

Definition at line 104 of file WLSTrajectory.hh.

References test::b.

Referenced by WLSEventAction::EndOfEventAction().

104 {fForceNoDraw=b;}
void WLSTrajectory::ShowTrajectory ( std::ostream &  os = G4cout) const
virtual

Reimplemented from G4VTrajectory.

Definition at line 132 of file WLSTrajectory.cc.

References G4VTrajectory::ShowTrajectory().

133 {
134  // Invoke the default implementation in G4VTrajectory...
136  // ... or override with your own code here.
137 }
virtual void ShowTrajectory(std::ostream &os=G4cout) const
void WLSTrajectory::WLS ( )
inline

Definition at line 102 of file WLSTrajectory.hh.

Referenced by WLSTrackingAction::PostUserTrackingAction().

102 {fWLS=true;}

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