QwGeant4
QweakSimTrajectory Class Reference

Stores the information about the various tracks. More...

#include <QweakSimTrajectory.hh>

Inherits G4VTrajectory.

Public Member Functions

 QweakSimTrajectory ()
 
 QweakSimTrajectory (const G4Track *aTrack)
 
 QweakSimTrajectory (QweakSimTrajectory &)
 
virtual ~QweakSimTrajectory ()
 
void * operator new (size_t)
 
void operator delete (void *)
 
int operator== (const QweakSimTrajectory &right) const
 
G4int GetTrackID () const
 
G4int GetParentID () const
 
G4String GetParticleName () const
 
G4double GetCharge () const
 
G4int GetPDGEncoding () const
 
G4ThreeVector GetInitialMomentum () const
 
virtual int GetPointEntries () const
 
virtual G4VTrajectoryPoint * GetPoint (G4int i) const
 
virtual void ShowTrajectory (std::ostream &os=G4cout) const
 
virtual void DrawTrajectory (G4int i_mode=0) const
 
virtual void AppendStep (const G4Step *aStep)
 
virtual void MergeTrajectory (G4VTrajectory *secondTrajectory)
 
void SetDrawTrajectory (G4bool b)
 
void SetForceDrawTrajectory (G4bool b)
 
void SetForceNoDrawTrajectory (G4bool b)
 
G4ParticleDefinition * GetParticleDefinition ()
 

Private Attributes

QweakSimTrajectoryPointContainerpositionRecord
 
G4int fTrackID
 
G4int fParentID
 
G4ParticleDefinition * fpParticleDefinition
 
G4String ParticleName
 
G4double PDGCharge
 
G4int PDGEncoding
 
G4ThreeVector InitialMomentum
 
G4ThreeVector momentum
 
G4ThreeVector vertexPosition
 
G4double globalTime
 
G4bool drawit
 
G4bool forceNoDraw
 
G4bool forceDraw
 

Detailed Description

Stores the information about the various tracks.

It can be viewed as a perfect track reconstruction in a hypothetical perfect Qweak setup

Placeholder for a long explaination

Definition at line 71 of file QweakSimTrajectory.hh.

Constructor & Destructor Documentation

QweakSimTrajectory::QweakSimTrajectory ( )

Definition at line 60 of file QweakSimTrajectory.cc.

References drawit, forceDraw, forceNoDraw, fParentID, fpParticleDefinition, fTrackID, globalTime, momentum, ParticleName, PDGCharge, PDGEncoding, positionRecord, and vertexPosition.

61 {
62  G4cout << G4endl << "###### Calling QweakSimTrajectory::QweakSimTrajectory()" << G4endl << G4endl;
63 
64  drawit = FALSE;
65  forceNoDraw = FALSE;
66  forceDraw = FALSE;
67 
68 
70  ParticleName = "";
71  PDGCharge = 0;
72  PDGEncoding = 0;
73  fTrackID = 0;
74  fParentID = 0;
75  positionRecord = 0;
76 
77  momentum = G4ThreeVector(0.,0.,0.);
78  vertexPosition = G4ThreeVector(0.,0.,0.);
79  globalTime = 0.;
80 
81  G4cout << G4endl << "###### Leaving QweakSimTrajectory::QweakSimTrajectory()" << G4endl << G4endl;
82 }
G4ThreeVector vertexPosition
G4ParticleDefinition * fpParticleDefinition
QweakSimTrajectoryPointContainer * positionRecord
QweakSimTrajectory::QweakSimTrajectory ( const G4Track *  aTrack)

Definition at line 85 of file QweakSimTrajectory.cc.

References fParentID, fpParticleDefinition, fTrackID, globalTime, momentum, ParticleName, PDGCharge, PDGEncoding, positionRecord, and vertexPosition.

86 {
87  // G4cout << G4endl << "###### Calling QweakSimTrajectory::QweakSimTrajectory(aTrack)" << G4endl << G4endl;
88 
89  fpParticleDefinition = aTrack->GetDefinition();
90  ParticleName = fpParticleDefinition->GetParticleName();
91  PDGCharge = fpParticleDefinition->GetPDGCharge();
92  PDGEncoding = fpParticleDefinition->GetPDGEncoding();
93  fTrackID = aTrack->GetTrackID();
94  fParentID = aTrack->GetParentID();
95 
97  positionRecord -> push_back(new G4TrajectoryPoint(aTrack->GetPosition()));
98 
99  momentum = aTrack->GetMomentum();
100  vertexPosition = aTrack->GetPosition();
101  globalTime = aTrack->GetGlobalTime();
102 
103 
104  // G4cout << G4endl << "###### Leaving QweakSimTrajectory::QweakSimTrajectory(aTrack)" << G4endl << G4endl;
105 }
G4ThreeVector vertexPosition
G4ParticleDefinition * fpParticleDefinition
std::vector< G4VTrajectoryPoint * > QweakSimTrajectoryPointContainer
QweakSimTrajectoryPointContainer * positionRecord
QweakSimTrajectory::QweakSimTrajectory ( QweakSimTrajectory right)

Definition at line 108 of file QweakSimTrajectory.cc.

References fParentID, fpParticleDefinition, fTrackID, globalTime, momentum, ParticleName, PDGCharge, PDGEncoding, positionRecord, and vertexPosition.

109  : G4VTrajectory()
110 {
111  ParticleName = right.ParticleName;
113  PDGCharge = right.PDGCharge;
114  PDGEncoding = right.PDGEncoding;
115  fTrackID = right.fTrackID;
116  fParentID = right.fParentID;
117 
119 
120  for(int i=0;i<(int)right.positionRecord->size();i++)
121  {
122  G4TrajectoryPoint* rightPoint = (G4TrajectoryPoint*)((*(right.positionRecord))[i]);
123 
124  positionRecord->push_back(new G4TrajectoryPoint(*rightPoint));
125  }
126 
127  momentum = right.momentum;
129  globalTime = right.globalTime;
130 }
G4ThreeVector vertexPosition
G4ParticleDefinition * fpParticleDefinition
std::vector< G4VTrajectoryPoint * > QweakSimTrajectoryPointContainer
QweakSimTrajectoryPointContainer * positionRecord
QweakSimTrajectory::~QweakSimTrajectory ( )
virtual

Definition at line 133 of file QweakSimTrajectory.cc.

References positionRecord.

134 {
135  size_t i;
136  for(i=0;i<positionRecord->size();i++){
137  delete (*positionRecord)[i];
138  }
139  positionRecord->clear();
140 
141  delete positionRecord;
142 }
QweakSimTrajectoryPointContainer * positionRecord

Member Function Documentation

void QweakSimTrajectory::AppendStep ( const G4Step *  aStep)
virtual

Definition at line 288 of file QweakSimTrajectory.cc.

References positionRecord.

289 {
290  positionRecord->push_back( new G4TrajectoryPoint(aStep->GetPostStepPoint()->GetPosition() ));
291 }
QweakSimTrajectoryPointContainer * positionRecord
void QweakSimTrajectory::DrawTrajectory ( G4int  i_mode = 0) const
virtual

Definition at line 172 of file QweakSimTrajectory.cc.

Referenced by QweakSimEventAction::EndOfEventAction().

176 {
177  // G4cout << G4endl << "###### Calling QweakSimTrajectory::DrawTrajectory()" << G4endl << G4endl;
178 
179  //if(!forceDraw && (!drawit || forceNoDraw))
180 
181 // if( (forceNoDraw == true))
182 // return;
183 
184  // If i_mode>=0, draws a trajectory as a polyline (blue for
185  // positive, red for negative, green for neutral) and, if i_mode!=0,
186  // adds markers - yellow circles for step points and magenta squares
187  // for auxiliary points, if any - whose screen size in pixels is
188  // given by abs(i_mode)/1000. E.g: i_mode = 5000 gives easily
189  // visible markers.
190 
191  G4VVisManager* pVVisManager = G4VVisManager::GetConcreteInstance();
192 
193  G4ThreeVector pos;
194  G4Polyline pPolyline;
195 
196  for (int i = 0; i < (int) positionRecord->size() ; i++) {
197 
198  G4TrajectoryPoint* aTrajectoryPoint = (G4TrajectoryPoint*)((*positionRecord)[i]);
199  pos = aTrajectoryPoint->GetPosition();
200  pPolyline.push_back( pos );
201  }
202 
203 
204  G4Colour colour(0.75,0.75,0.75); // LightGray
205 
206  G4Colour red ( 255/255., 0/255., 0/255.);
207  G4Colour blue ( 0/255., 0/255., 255/255.);
208  G4Colour green ( 0/255., 255/255., 0/255.);
209  G4Colour yellow ( 255/255., 255/255., 0/255.);
210 
211  G4Colour white ( 255/255., 255/255., 255/255.);
212 
213  G4Colour orange ( 255/255., 127/255., 0/255.);
214  G4Colour magenta ( 237/255., 173/255., 255/255.);
215  G4Colour magenta1 ( 104/255., 49/255., 94/255.);
216 
217 
218 
219  // Default colors for particlea
220  if ( PDGCharge < 0.) {colour = red;}
221  else if ( PDGCharge > 0.) {colour = blue;}
222  else
223  {
224  //G4cout << "%%%%% We have something neutral here !!!" << G4endl;
225  colour = green;
226  }
227 
228  if (fpParticleDefinition == G4Gamma ::GammaDefinition() ) {
229  //G4cout << "%%%%% We have a green Gamma here !!!" << G4endl;
230  colour = green;}
231  if (fpParticleDefinition == G4Electron ::ElectronDefinition() )
232  {
233  //G4cout << "%%%%% We have a red Electron here !!!" << G4endl;
234  colour = red;
235  //colour = white;
236  }
237  if (fpParticleDefinition == G4Positron ::PositronDefinition() )
238  {
239  //G4cout << "%%%%% We have a blue Positron here !!!" << G4endl;
240  colour = blue;
241  }
242 
243 // else if (fpParticleDefinition == G4MuonMinus ::MuonMinusDefinition()
244 // ||fpParticleDefinition == G4MuonPlus ::MuonPlusDefinition()) ) colour = G4Colour(1.,0.,1.); // Magenta
245 //
246 // else if(fpParticleDefinition->GetParticleType()=="meson")
247 // {
248 // if(PDGCharge!=0.)
249 // colour = yellow;
250 // else
251 // colour = G4Colour(0.5,0.,0.); // HalfRed
252 // }
253 // else if(fpParticleDefinition->GetParticleType()=="baryon")
254 // {
255 // if(PDGCharge!=0.)
256 // colour = orange; // Orange
257 // else
258 // colour = G4Colour(0.5,0.39,0.);// HalfOrange
259 // }
260 
261  if( fpParticleDefinition == G4Neutron ::NeutronDefinition())
262  {
263  //G4cout << "%%%%% We have a white Neutron here !!!" << G4endl;
264  colour = white;
265  }
266  if( fpParticleDefinition == G4Proton ::ProtonDefinition())
267  {
268  //G4cout << "%%%%% We have an orange Proton here !!!" << G4endl;
269  colour = orange;
270  }
271  if( fpParticleDefinition == G4OpticalPhoton ::OpticalPhotonDefinition())
272  {
273  //G4cout << "%%%%% We have a magenta OpticalPhoton here !!!" << G4endl;
274  colour = magenta;
275  }
276 
277  G4VisAttributes attribs(colour);
278  pPolyline.SetVisAttributes(attribs);
279 
280  // if(pVVisManager)
281  pVVisManager->Draw(pPolyline);
282 
283 
284  // G4cout << G4endl << "###### Leaving QweakSimTrajectory::DrawTrajectory()" << G4endl << G4endl;
285 }
G4ParticleDefinition * fpParticleDefinition
QweakSimTrajectoryPointContainer * positionRecord

+ Here is the caller graph for this function:

G4double QweakSimTrajectory::GetCharge ( ) const
inline

Definition at line 99 of file QweakSimTrajectory.hh.

References PDGCharge.

Referenced by QweakSimCerenkov_PMTSD::GetParentTrajectory(), and QweakSimPMTOnly_PMTSD::GetParentTrajectory().

100  { return PDGCharge; }

+ Here is the caller graph for this function:

G4ThreeVector QweakSimTrajectory::GetInitialMomentum ( ) const
inline

Definition at line 103 of file QweakSimTrajectory.hh.

References InitialMomentum.

104  { return InitialMomentum; }
G4ThreeVector InitialMomentum
G4int QweakSimTrajectory::GetParentID ( ) const
inline

Definition at line 95 of file QweakSimTrajectory.hh.

References fParentID.

96  { return fParentID; }
G4ParticleDefinition * QweakSimTrajectory::GetParticleDefinition ( )

Definition at line 294 of file QweakSimTrajectory.cc.

References ParticleName.

295 {
296  return (G4ParticleTable::GetParticleTable()->FindParticle(ParticleName));
297 }
G4String QweakSimTrajectory::GetParticleName ( ) const
inline

Definition at line 97 of file QweakSimTrajectory.hh.

References ParticleName.

98  { return ParticleName; }
G4int QweakSimTrajectory::GetPDGEncoding ( ) const
inline

Definition at line 101 of file QweakSimTrajectory.hh.

References PDGEncoding.

102  { return PDGEncoding; }
virtual G4VTrajectoryPoint* QweakSimTrajectory::GetPoint ( G4int  i) const
inlinevirtual

Definition at line 108 of file QweakSimTrajectory.hh.

References positionRecord.

109  { return (*positionRecord)[i]; }
QweakSimTrajectoryPointContainer * positionRecord
virtual int QweakSimTrajectory::GetPointEntries ( ) const
inlinevirtual

Definition at line 106 of file QweakSimTrajectory.hh.

References positionRecord.

Referenced by MergeTrajectory().

107  { return positionRecord->size(); }
QweakSimTrajectoryPointContainer * positionRecord

+ Here is the caller graph for this function:

G4int QweakSimTrajectory::GetTrackID ( ) const
inline

Definition at line 93 of file QweakSimTrajectory.hh.

References fTrackID.

Referenced by QweakSimCerenkov_PMTSD::GetParentTrajectory(), and QweakSimPMTOnly_PMTSD::GetParentTrajectory().

94  { return fTrackID; }

+ Here is the caller graph for this function:

void QweakSimTrajectory::MergeTrajectory ( G4VTrajectory *  secondTrajectory)
virtual

Definition at line 300 of file QweakSimTrajectory.cc.

References GetPointEntries(), and positionRecord.

301 {
302  if(!secondTrajectory) return;
303 
304  QweakSimTrajectory* seco = (QweakSimTrajectory*)secondTrajectory;
305  G4int ent = seco->GetPointEntries();
306 
307  for(int i=1;i<ent;i++) // initial point of the second trajectory should not be merged
308  {
309  positionRecord->push_back((*(seco->positionRecord))[i]);
310  }
311 
312  delete (*seco->positionRecord)[0];
313 
314  seco->positionRecord->clear();
315 
316 }
Stores the information about the various tracks.
virtual int GetPointEntries() const
QweakSimTrajectoryPointContainer * positionRecord

+ Here is the call graph for this function:

void QweakSimTrajectory::operator delete ( void *  aTrajectory)
inline

Definition at line 166 of file QweakSimTrajectory.hh.

References myTrajectoryAllocator.

167 {
168  myTrajectoryAllocator.FreeSingle((QweakSimTrajectory*)aTrajectory);
169 }
Stores the information about the various tracks.
G4Allocator< QweakSimTrajectory > myTrajectoryAllocator
void * QweakSimTrajectory::operator new ( size_t  )
inline

Definition at line 158 of file QweakSimTrajectory.hh.

References myTrajectoryAllocator.

159 {
160  void* aTrajectory;
161  aTrajectory = (void*)myTrajectoryAllocator.MallocSingle();
162  return aTrajectory;
163 }
G4Allocator< QweakSimTrajectory > myTrajectoryAllocator
int QweakSimTrajectory::operator== ( const QweakSimTrajectory right) const
inline

Definition at line 89 of file QweakSimTrajectory.hh.

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

Definition at line 123 of file QweakSimTrajectory.hh.

References drawit.

123 {drawit=b;}
void QweakSimTrajectory::SetForceDrawTrajectory ( G4bool  b)
inline

Definition at line 124 of file QweakSimTrajectory.hh.

References forceDraw.

void QweakSimTrajectory::SetForceNoDrawTrajectory ( G4bool  b)
inline

Definition at line 125 of file QweakSimTrajectory.hh.

References forceNoDraw.

void QweakSimTrajectory::ShowTrajectory ( std::ostream &  os = G4cout) const
virtual

Definition at line 145 of file QweakSimTrajectory.cc.

References fParentID, fTrackID, globalTime, momentum, ParticleName, PDGCharge, PDGEncoding, positionRecord, and vertexPosition.

146 {
147 
148  os << G4endl << "TrackID =" << fTrackID
149  << " : ParentID=" << fParentID << G4endl;
150 
151  os << "Particle name : " << ParticleName << " PDG code : " << PDGEncoding
152  << " Charge : " << PDGCharge << G4endl;
153 
154  os << "Original momentum : " <<
155  G4BestUnit(momentum,"Energy") << G4endl;
156 
157  os << "Vertex : " << G4BestUnit(vertexPosition,"Length")
158  << " Global time : " << G4BestUnit(globalTime,"Time") << G4endl;
159 
160  os << " Current trajectory has " << positionRecord->size()
161  << " points." << G4endl;
162 
163  for( size_t i=0 ; i < positionRecord->size() ; i++){
164  G4TrajectoryPoint* aTrajectoryPoint = (G4TrajectoryPoint*)((*positionRecord)[i]);
165 
166  G4cout << "Point[" << i << "]" << " Position= " << aTrajectoryPoint->GetPosition() << G4endl;
167  }
168 }
G4ThreeVector vertexPosition
QweakSimTrajectoryPointContainer * positionRecord

Field Documentation

G4bool QweakSimTrajectory::drawit
private

Definition at line 148 of file QweakSimTrajectory.hh.

Referenced by QweakSimTrajectory(), and SetDrawTrajectory().

G4bool QweakSimTrajectory::forceDraw
private

Definition at line 150 of file QweakSimTrajectory.hh.

Referenced by QweakSimTrajectory(), and SetForceDrawTrajectory().

G4bool QweakSimTrajectory::forceNoDraw
private

Definition at line 149 of file QweakSimTrajectory.hh.

Referenced by QweakSimTrajectory(), and SetForceNoDrawTrajectory().

G4int QweakSimTrajectory::fParentID
private

Definition at line 135 of file QweakSimTrajectory.hh.

Referenced by GetParentID(), QweakSimTrajectory(), and ShowTrajectory().

G4ParticleDefinition* QweakSimTrajectory::fpParticleDefinition
private

Definition at line 136 of file QweakSimTrajectory.hh.

Referenced by QweakSimTrajectory().

G4int QweakSimTrajectory::fTrackID
private

Definition at line 134 of file QweakSimTrajectory.hh.

Referenced by GetTrackID(), QweakSimTrajectory(), and ShowTrajectory().

G4double QweakSimTrajectory::globalTime
private

Definition at line 145 of file QweakSimTrajectory.hh.

Referenced by QweakSimTrajectory(), and ShowTrajectory().

G4ThreeVector QweakSimTrajectory::InitialMomentum
private

Definition at line 141 of file QweakSimTrajectory.hh.

Referenced by GetInitialMomentum().

G4ThreeVector QweakSimTrajectory::momentum
private

Definition at line 143 of file QweakSimTrajectory.hh.

Referenced by QweakSimTrajectory(), and ShowTrajectory().

G4String QweakSimTrajectory::ParticleName
private
G4double QweakSimTrajectory::PDGCharge
private

Definition at line 138 of file QweakSimTrajectory.hh.

Referenced by GetCharge(), QweakSimTrajectory(), and ShowTrajectory().

G4int QweakSimTrajectory::PDGEncoding
private

Definition at line 139 of file QweakSimTrajectory.hh.

Referenced by GetPDGEncoding(), QweakSimTrajectory(), and ShowTrajectory().

G4ThreeVector QweakSimTrajectory::vertexPosition
private

Definition at line 144 of file QweakSimTrajectory.hh.

Referenced by QweakSimTrajectory(), and ShowTrajectory().


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