Public Member Functions | Private Attributes | List of all members
LBNETrajectory Class Reference

#include <LBNETrajectory.hh>

Inheritance diagram for LBNETrajectory:

Public Member Functions

 LBNETrajectory ()
 
 LBNETrajectory (const G4Track *aTrack)
 
 LBNETrajectory (LBNETrajectory &)
 
 LBNETrajectory (const LBNETrajectory *right)
 
virtual ~LBNETrajectory ()
 
void * operator new (size_t)
 
void operator delete (void *)
 
int operator== (const LBNETrajectory &right) const
 
G4int GetTrackID () const
 
G4int GetParentID () const
 
G4String GetParticleName () const
 
G4double GetCharge () const
 
G4double GetMass () const
 
G4int GetPDGEncoding () const
 
const G4ThreeVector & GetVertexPosition () const
 
virtual int GetPointEntries () const
 
virtual G4VTrajectoryPoint * GetPoint (G4int i) const
 
virtual G4ThreeVector GetMomentum (G4int i) const
 
virtual G4String GetPreStepVolumeName (G4int i) const
 
G4ThreeVector GetInitialMomentum () const
 
virtual G4int GetTgen () const
 
G4int GetDecayCode () const
 
virtual G4double GetNImpWt () const
 
virtual G4double GetStepLength (G4int i) const
 
G4int GetMaterialNumber1rst () const
 
void SetMaterialNumber1rst (G4double matNum)
 
G4int GetMaterialNumberLast () const
 
void SetMaterialNumberLast (G4double matNum)
 
G4String GetVolName1rst () const
 
G4double GetTimeStart () const
 
G4ThreeVector GetPolarization () const
 
G4int GetPDGNucleus () const
 
G4String GetProcessName () const
 
G4String GetMaterialName1rst () const
 
virtual G4String GetMaterialName (G4int i) const
 
const G4ThreeVector GetParentMomentumAtThisProduction () const
 
virtual void ShowTrajectory () const
 
virtual void ShowTrajectory (std::ostream &o) const
 
virtual void DrawTrajectory () const
 
virtual const std::map< G4String, G4AttDef > * GetAttDefs () const
 
virtual std::vector< G4AttValue > * CreateAttValues () const
 
virtual void AppendStep (const G4Step *aStep)
 
virtual void MergeTrajectory (G4VTrajectory *secondTrajectory)
 
G4ParticleDefinition * GetParticleDefinition ()
 

Private Attributes

LBNETrajectoryMomentumContainerfMomentumRecord
 
LBNETrajectoryPointContainerfPositionRecord
 
G4int fDecayCode
 
G4int fEventNo
 
G4int fTgen
 
G4double fNImpWt
 
G4int fTrackID
 
G4int fParentID
 
G4ParticleDefinition * fParticleDefinition
 
G4String fParticleName
 
G4double fPDGCharge
 
G4int fPDGEncoding
 
G4ThreeVector fMomentum
 
G4ThreeVector fVertexPosition
 
G4double fParticleMass
 
LBNETrajectoryVolumeNamefPreStepVolume
 
DVecfStepLength
 
G4int fMaterialNumber1rst
 
G4String fVolName1rst
 
G4int fMaterialNumberLast
 
G4double fTimeStart
 
G4ThreeVector fPolarization
 
G4int fPDGNucleus
 
G4String fProcessName
 
G4ThreeVector fParentMomentumAtThisProduction
 
LBNETrajectoryMaterialNamefMaterialName
 
G4String fMaterialName1rst
 

Detailed Description

Definition at line 32 of file LBNETrajectory.hh.

Constructor & Destructor Documentation

LBNETrajectory::LBNETrajectory ( )

Definition at line 30 of file LBNETrajectory.cc.

31 {
33  fParticleName = "";
34  fPDGCharge = 0;
35  fPDGEncoding = 0;
36  fTrackID = 0;
37  fParentID = 0;
38  fPositionRecord = 0;
39  fMomentum = G4ThreeVector(0.,0.,0.);
40  fMomentumRecord = 0;
41  fVertexPosition = G4ThreeVector(0.,0.,0.);
42  fParticleMass = 0.;
43  fDecayCode = 0;
44  fTgen = 0;
45  fNImpWt = 1.;
46  fPreStepVolume = 0;
47  fStepLength = 0;
51  fTimeStart = 0.;
52  fPolarization = G4ThreeVector(0.,0.,0.);
53  fParentMomentumAtThisProduction = G4ThreeVector(0.,0.,0.); //for ppfx....A. Bashyal
54  fPDGNucleus = 0;
57 
58 }
std::string string
Definition: nybbler.cc:12
LBNETrajectoryVolumeName * fPreStepVolume
G4ParticleDefinition * fParticleDefinition
LBNETrajectoryMomentumContainer * fMomentumRecord
G4String fMaterialName1rst
G4String fVolName1rst
G4double fParticleMass
G4ThreeVector fVertexPosition
G4ThreeVector fPolarization
LBNETrajectoryPointContainer * fPositionRecord
G4ThreeVector fMomentum
G4ThreeVector fParentMomentumAtThisProduction
G4String fParticleName
G4String fProcessName
LBNETrajectory::LBNETrajectory ( const G4Track *  aTrack)

Definition at line 60 of file LBNETrajectory.cc.

61 {
62  fParticleDefinition = aTrack->GetDefinition();
63  fParticleName = fParticleDefinition->GetParticleName();
64  fPDGCharge = fParticleDefinition->GetPDGCharge();
65  fPDGEncoding = fParticleDefinition->GetPDGEncoding();
66  fTrackID = aTrack->GetTrackID();
67  fParentID = aTrack->GetParentID();
69  fPositionRecord -> push_back(new G4TrajectoryPoint(aTrack->GetPosition()));
71  fMomentumRecord -> push_back(aTrack->GetMomentum());
73  fPreStepVolume -> push_back(aTrack->GetVolume()->GetName());
75  fMaterialName -> push_back(aTrack->GetVolume()->GetLogicalVolume()->GetMaterial()->GetName());
76  fStepLength = new DVec();
77  fStepLength -> push_back(aTrack->GetStepLength());
78  fMomentum = aTrack->GetMomentum();
79  fVertexPosition = aTrack->GetPosition();
80  fParticleMass = aTrack->GetDefinition()->GetPDGMass();
81  if (!aTrack->GetLogicalVolumeAtVertex()) {
82 // std::cerr << " LBNETrajectory::LBNETrajectory No logical volume at vertex... ??? " << std::endl;
84  fMaterialName1rst = std::string("Unknown");
85  fVolName1rst = std::string("Unknown");
86  } else {
87  fMaterialNumber1rst = aTrack->GetLogicalVolumeAtVertex()->GetMaterial()->GetIndex();
88  fMaterialName1rst = aTrack->GetLogicalVolumeAtVertex()->GetMaterial()->GetName();
89  fVolName1rst = aTrack->GetLogicalVolumeAtVertex()->GetName();
90  }
91  fTimeStart = aTrack->GetGlobalTime();
92  fPolarization = aTrack->GetPolarization();
93 
94  LBNETrackInformation* info = (LBNETrackInformation*)(aTrack->GetUserInformation());
95  if (info!=0)
96  {
97  fDecayCode = info->GetDecayCode();
98  fTgen = info->GetTgen();
99  fNImpWt = info->GetNImpWt();
101 // std::cerr << " Transfering ImpWeight = " << fNImpWt << " from trInfo, track id " << aTrack->GetTrackID() << std::endl;
102  }
103  else
104  {
105  fDecayCode = 0;
106  fTgen = 0;
107  fNImpWt = 1.;
108  fParentMomentumAtThisProduction = G4ThreeVector(0.,0.,0.);
109  }
110  // Now the tricky bit: Dk2Nu package wants the process name and the nucleus
111  if (aTrack->GetCreatorProcess() != 0)
112  fProcessName = aTrack->GetCreatorProcess()->GetProcessName();
113  else fProcessName = std::string ("BeamParticle");
114  // Well, that one was indeed easy. The PDG of the target nucleus now..
115  fPDGNucleus = 0;
116 // if ((aTrack->GetCreatorProcess() != 0) && (aTrack->GetCreatorProcess()->GetProcessType() == fHadronic)) {
117 // std::cerr << " Hadronic interaction, what do I do to get the nucleus that got hit?!!!" << std::endl;
118 // exit(2);
119 // The goal is dig through the creator process ... but the info is lost..
120 //
121 // const G4HadronicProcess *hProc = static_cast<const G4HadronicProcess* >(aTrack->GetCreatorProcess());
122 // std::cerr << " Current interaction length " << hProc->GetCurrentInteractionLength()
123 // << " name " << hProc->GetProcessName() << " Process sub type " << hProc->GetProcessSubType() << std::endl;
124 // G4IsoParticleChange *ic = G4HadronicProcess::GetIsotopeProductionInfo();
125 // if (ic !=0) {
126 // std::cerr << " .... No Isotope production info.. " << std::endl;
127 // } else {
128 // std::cerr << " .... Isotope is " << ic->GetIsotope() << std::endl;
129 // This access actually causes a crash.. So, the trail is cold regarding getting info about the isotope in questions..
130 // Paul Lebrun, April 2014.
131 // }
132 // }
133 
134 }
std::vector< G4double > DVec
std::string string
Definition: nybbler.cc:12
LBNETrajectoryVolumeName * fPreStepVolume
std::vector< G4VTrajectoryPoint * > LBNETrajectoryPointContainer
G4ParticleDefinition * fParticleDefinition
LBNETrajectoryMaterialName * fMaterialName
LBNETrajectoryMomentumContainer * fMomentumRecord
G4String fMaterialName1rst
G4String fVolName1rst
G4double fParticleMass
G4ThreeVector fVertexPosition
G4ThreeVector fPolarization
std::vector< G4String > LBNETrajectoryVolumeName
std::vector< G4String > LBNETrajectoryMaterialName
G4double GetNImpWt() const
LBNETrajectoryPointContainer * fPositionRecord
G4ThreeVector fMomentum
G4ThreeVector fParentMomentumAtThisProduction
G4String fParticleName
std::vector< G4ThreeVector > LBNETrajectoryMomentumContainer
G4String fProcessName
G4ThreeVector GetParentMomentumAtThisProduction() const
LBNETrajectory::LBNETrajectory ( LBNETrajectory right)

Definition at line 137 of file LBNETrajectory.cc.

138  : G4VTrajectory()
139 {
142  fPDGCharge = right.fPDGCharge;
143  fPDGEncoding = right.fPDGEncoding;
144  fTrackID = right.fTrackID;
145  fParentID = right.fParentID;
146  fMomentum = right.fMomentum;
149  fDecayCode = right.fDecayCode;
150  fTgen = right.fTgen;
151  fNImpWt = right.fNImpWt;
155  fStepLength = new DVec();
158 
159  for(size_t i=0;i<right.fPositionRecord->size();i++)
160  {
161  G4TrajectoryPoint* rightPoint = (G4TrajectoryPoint*)((*(right.fPositionRecord))[i]);
162  fPositionRecord->push_back(new G4TrajectoryPoint(*rightPoint));
163  }
164  for(size_t i=0;i<right.fMomentumRecord->size();i++)
165  {
166  G4ThreeVector rightMomentum = (G4ThreeVector)((*(right.fMomentumRecord))[i]);
167  fMomentumRecord->push_back(rightMomentum);
168  }
169  for(size_t i=0;i<right.fPreStepVolume->size();i++)
170  {
171  G4String rightPreStepVolume=(G4String)((*(right.fPreStepVolume))[i]);
172  fPreStepVolume->push_back(rightPreStepVolume);
173  }
174  for(size_t i=0;i<right.fStepLength->size();i++)
175  {
176  G4double rightsteplength =(G4double)((*(right.fStepLength))[i]);
177  fStepLength->push_back(rightsteplength);
178  }
179  for(size_t i=0;i<right.fMaterialName->size();i++)
180  {
181  G4String rightMaterialName = (G4String)((*(right.fMaterialName))[i]);
182  fMaterialName->push_back(rightMaterialName);
183  }
185  fVolName1rst = right.fVolName1rst;
187  fTimeStart = right.fTimeStart;
189  fPDGNucleus = right.fPDGNucleus;
190  fProcessName = right.fProcessName;
193 
194 
195 }
std::vector< G4double > DVec
LBNETrajectoryVolumeName * fPreStepVolume
std::vector< G4VTrajectoryPoint * > LBNETrajectoryPointContainer
G4ParticleDefinition * fParticleDefinition
LBNETrajectoryMaterialName * fMaterialName
LBNETrajectoryMomentumContainer * fMomentumRecord
G4String fMaterialName1rst
G4String fVolName1rst
G4double fParticleMass
G4ThreeVector fVertexPosition
G4ThreeVector fPolarization
std::vector< G4String > LBNETrajectoryVolumeName
std::vector< G4String > LBNETrajectoryMaterialName
LBNETrajectoryPointContainer * fPositionRecord
G4ThreeVector fMomentum
G4ThreeVector fParentMomentumAtThisProduction
G4String fParticleName
std::vector< G4ThreeVector > LBNETrajectoryMomentumContainer
G4String fProcessName
LBNETrajectory::LBNETrajectory ( const LBNETrajectory right)

Definition at line 198 of file LBNETrajectory.cc.

199  : G4VTrajectory()
200 {
201  fParticleName = right -> fParticleName;
203  fPDGCharge = right -> fPDGCharge;
204  fPDGEncoding = right -> fPDGEncoding;
205  fTrackID = right -> fTrackID;
206  fParentID = right -> fParentID;
207  fMomentum = right -> fMomentum;
208  fVertexPosition = right -> fVertexPosition;
209  fParticleMass = right -> fParticleMass;
210  fDecayCode = right -> fDecayCode;
211  fTgen = right -> fTgen;
212  fNImpWt = right -> fNImpWt;
215 
216 
218  for(size_t i=0; i < (right->fPositionRecord)->size(); ++i)
219  {
220  G4TrajectoryPoint* rightPoint = (G4TrajectoryPoint*)((*(right->fPositionRecord))[i]);
221  fPositionRecord->push_back(new G4TrajectoryPoint(*rightPoint));
222  }
223 
225  for(size_t i=0; i < (right->fMomentumRecord)->size(); ++i)
226  {
227  G4ThreeVector rightMomentum = (G4ThreeVector)((*(right->fMomentumRecord))[i]);
228  fMomentumRecord->push_back(rightMomentum);
229  }
230 
232  for(size_t i=0; i < (right->fPreStepVolume)->size(); ++i)
233  {
234  G4String rightPreStepVolume =(G4String)((*(right->fPreStepVolume))[i]);
235  fPreStepVolume->push_back(rightPreStepVolume);
236  }
237 
238  fStepLength = new DVec();
239  for(size_t i=0; i < (right->fStepLength)->size(); ++i)
240  {
241  G4double rightsteplength =(G4double)((*(right->fStepLength))[i]);
242  fStepLength->push_back(rightsteplength);
243  }
244  for(size_t i=0;i<right->fMaterialName->size();i++)
245  {
246  G4String rightMaterialName = (G4String)((*(right->fMaterialName))[i]);
247  fMaterialName->push_back(rightMaterialName);
248  }
249 
251  fVolName1rst = right->fVolName1rst;
253  fTimeStart = right->fTimeStart;
254  fPolarization = right->fPolarization;
255  fPDGNucleus = right->fPDGNucleus;
256  fProcessName = right->fProcessName;
258  fMaterialName = right->fMaterialName;
259 
260 
261 }
std::vector< G4double > DVec
LBNETrajectoryVolumeName * fPreStepVolume
std::vector< G4VTrajectoryPoint * > LBNETrajectoryPointContainer
G4ParticleDefinition * fParticleDefinition
LBNETrajectoryMaterialName * fMaterialName
LBNETrajectoryMomentumContainer * fMomentumRecord
G4String fMaterialName1rst
G4String fVolName1rst
G4double fParticleMass
G4ThreeVector fVertexPosition
G4ThreeVector fPolarization
std::vector< G4String > LBNETrajectoryVolumeName
std::vector< G4String > LBNETrajectoryMaterialName
LBNETrajectoryPointContainer * fPositionRecord
G4ThreeVector fMomentum
G4ThreeVector fParentMomentumAtThisProduction
G4String fParticleName
std::vector< G4ThreeVector > LBNETrajectoryMomentumContainer
G4String fProcessName
LBNETrajectory::~LBNETrajectory ( )
virtual

Definition at line 264 of file LBNETrajectory.cc.

265 {
266  size_t i;
267  for(i=0;i<fPositionRecord->size();i++)
268  {
269  delete (*fPositionRecord)[i];
270  }
271  fPositionRecord->clear();
272  delete fPositionRecord;
273 
274 
275  fMomentumRecord->clear();
276  delete fMomentumRecord;
277 
278  fPreStepVolume->clear();
279  delete fPreStepVolume;
280 
281  fStepLength->clear();
282  delete fStepLength;
283 
284  fMaterialName->clear();
285  delete fMaterialName;
286 
287 }
LBNETrajectoryVolumeName * fPreStepVolume
LBNETrajectoryMaterialName * fMaterialName
LBNETrajectoryMomentumContainer * fMomentumRecord
LBNETrajectoryPointContainer * fPositionRecord

Member Function Documentation

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

Definition at line 484 of file LBNETrajectory.cc.

485 {
487  ->push_back(new G4TrajectoryPoint(aStep->GetPostStepPoint()
488  ->GetPosition() ));
489  fMomentumRecord->push_back(aStep->GetPostStepPoint()->GetMomentum());
490 
491  G4Track* aTrack=aStep->GetTrack();
492  LBNETrackInformation* info=(LBNETrackInformation*)(aTrack->GetUserInformation());
493  if (info!=0) {
494  fDecayCode=info->GetDecayCode();
496  fTgen=info->GetTgen();
497  }
498  else {fDecayCode=-1;
499  fParentMomentumAtThisProduction = G4ThreeVector(0.,0.,0.);}
500  G4StepPoint * steppoint=aStep->GetPreStepPoint();
501  G4String PreVolumeName=steppoint->GetPhysicalVolume()->GetName();
502  G4String MaterialName =
503  steppoint->GetPhysicalVolume()->GetLogicalVolume()->GetMaterial()->GetName();
504  fPreStepVolume->push_back(PreVolumeName);
505  fStepLength->push_back(aStep->GetStepLength());
506  fMaterialName->push_back(MaterialName);
507 }
LBNETrajectoryVolumeName * fPreStepVolume
LBNETrajectoryMaterialName * fMaterialName
LBNETrajectoryMomentumContainer * fMomentumRecord
LBNETrajectoryPointContainer * fPositionRecord
G4ThreeVector fParentMomentumAtThisProduction
G4ThreeVector GetParentMomentumAtThisProduction() const
std::vector< G4AttValue > * LBNETrajectory::CreateAttValues ( ) const
virtual

Definition at line 452 of file LBNETrajectory.cc.

453 {
454  std::vector<G4AttValue>* values = new std::vector<G4AttValue>;
455 
456  values->push_back
457  (G4AttValue("TrkID",G4UIcommand::ConvertToString(fTrackID),""));
458 
459  values->push_back
460  (G4AttValue("ParTrkID",G4UIcommand::ConvertToString(fParentID),""));
461 
462  values->push_back(G4AttValue("Name",fParticleName,""));
463 
464  values->push_back
465  (G4AttValue("Q",G4UIcommand::ConvertToString(fPDGCharge),""));
466 
467  values->push_back
468  (G4AttValue("PDG",G4UIcommand::ConvertToString(fPDGEncoding),""));
469 
470  values->push_back
471  (G4AttValue("IMom",G4BestUnit(GetInitialMomentum(),"Energy"),""));
472 
473  values->push_back
474  (G4AttValue("IMomMag",G4BestUnit((GetInitialMomentum()).mag(),"Energy"),""));
475 
476  values->push_back
477  (G4AttValue("NPts",G4UIcommand::ConvertToString(GetPointEntries()),""));
478 
479  return values;
480 }
Q_UINT16 values[128]
G4ThreeVector GetInitialMomentum() const
virtual int GetPointEntries() const
G4String fParticleName
void LBNETrajectory::DrawTrajectory ( ) const
virtual

Definition at line 314 of file LBNETrajectory.cc.

315 {
316  G4VVisManager* pVVisManager = G4VVisManager::GetConcreteInstance();
317  G4ThreeVector pos;
318 
319  G4Polyline pPolyline;
320  for (size_t i = 0; i < fPositionRecord->size() ; i++) {
321  G4TrajectoryPoint* aTrajectoryPoint = (G4TrajectoryPoint*)((*fPositionRecord)[i]);
322  pos = aTrajectoryPoint->GetPosition();
323  pPolyline.push_back( pos );
324  }
325 
326  G4Colour color(1.0,1.0,1.0); //white
327 
328  if(fPDGCharge == 0.) color = G4Colour::White();
329  if(fPDGCharge == 1.0) color = G4Colour(1.0, 0.5, 0.0); //orange
330  if(fPDGCharge == -1.0) color = G4Colour(0.5, 0.0, 0.5); //purple
331 
332 
333  if (fParticleDefinition==G4Proton::ProtonDefinition()) color = G4Color::Green();
334  if (fParticleDefinition==G4Neutron::NeutronDefinition()) color = G4Color::Grey();
335  if (fParticleDefinition==G4PionMinus::PionMinusDefinition()) color = G4Color::Cyan();
336  if (fParticleDefinition==G4PionPlus::PionPlusDefinition()) color = G4Color::Blue();
337  if (fParticleDefinition==G4KaonMinus::KaonMinusDefinition()) color = G4Color::Yellow();
338  if (fParticleDefinition==G4KaonPlus::KaonPlusDefinition()) color = G4Color::Yellow();
339  if (fParticleDefinition==G4MuonMinus::MuonMinusDefinition()) color = G4Color::Red();
340  if (fParticleDefinition==G4MuonPlus::MuonPlusDefinition()) color = G4Color::Red();
341 
342  if ((fParticleDefinition == G4NeutrinoE::NeutrinoEDefinition()) ||
343  (fParticleDefinition == G4NeutrinoMu::NeutrinoMuDefinition()) ||
344  (fParticleDefinition == G4NeutrinoTau::NeutrinoTauDefinition()) ||
345  (fParticleDefinition == G4AntiNeutrinoE::AntiNeutrinoEDefinition()) ||
346  (fParticleDefinition == G4AntiNeutrinoMu::AntiNeutrinoMuDefinition()) ||
347  (fParticleDefinition == G4AntiNeutrinoTau::AntiNeutrinoTauDefinition()))
348  color = G4Color::Magenta();
349 
350 
351 
352 
353 
354 
355 
356 G4VisAttributes attribs(color);
357  pPolyline.SetVisAttributes(attribs);
358  if(pVVisManager) pVVisManager->Draw(pPolyline);
359 
360 /*
361  G4Colour colour(0.2,0.2,0.2);
362  if(fParticleDefinition==G4Gamma::GammaDefinition())
363  colour = G4Colour(0.,0.,1.);
364  else if(fParticleDefinition==G4Electron::ElectronDefinition()
365  ||fParticleDefinition==G4Positron::PositronDefinition())
366  colour = G4Colour(1.,1.,0.);
367  else if(fParticleDefinition==G4MuonMinus::MuonMinusDefinition()
368  ||fParticleDefinition==G4MuonPlus::MuonPlusDefinition())
369  colour = G4Colour(0.,1.,0.);
370  else if(fParticleDefinition->GetParticleType()=="meson")
371  {
372  if(fPDGCharge!=0.)
373  colour = G4Colour(1.,0.,0.);
374  else
375  colour = G4Colour(0.5,0.,0.);
376  }
377  else if(fParticleDefinition->GetParticleType()=="baryon")
378  {
379  if(fPDGCharge!=0.)
380  colour = G4Colour(0.,1.,1.);
381  else
382  colour = G4Colour(0.,0.5,0.5);
383  }
384 
385  //
386 
387  //draw only protons,pi+ and pi-
388  G4VisAttributes attribs;
389  if (fParticleDefinition==G4Proton::ProtonDefinition()) {
390  colour=G4Colour(0.,0.,1.);
391  attribs=G4VisAttributes(colour);
392  pPolyline.SetVisAttributes(attribs);
393  if(pVVisManager) pVVisManager->Draw(pPolyline);
394  }
395  if (fParticleDefinition==G4PionMinus::PionMinusDefinition()) {
396  colour=G4Colour(1.,0.,0.);
397  attribs=G4VisAttributes(colour);
398  pPolyline.SetVisAttributes(attribs);
399  if(pVVisManager) pVVisManager->Draw(pPolyline);
400  }
401  if (fParticleDefinition==G4PionPlus::PionPlusDefinition()) {
402  colour=G4Colour(0.,1.,0.);
403  attribs=G4VisAttributes(colour);
404  pPolyline.SetVisAttributes(attribs);
405  if(pVVisManager) pVVisManager->Draw(pPolyline);
406  }
407 
408 */
409 
410 }
G4ParticleDefinition * fParticleDefinition
std::size_t color(std::string const &procname)
LBNETrajectoryPointContainer * fPositionRecord
const std::map< G4String, G4AttDef > * LBNETrajectory::GetAttDefs ( ) const
virtual

Definition at line 413 of file LBNETrajectory.cc.

414 {
415  G4bool isNew;
416  std::map<G4String,G4AttDef>* store
417  = G4AttDefStore::GetInstance("LBNETrajectory",isNew);
418  if (isNew)
419  {
420 
421  G4String ID("TrkID");
422  (*store)[ID] = G4AttDef(ID,"Track ID","Physics","","G4int");
423 
424  G4String PID("ParTrkID");
425  (*store)[PID] = G4AttDef(PID,"Parent ID","Physics","","G4int");
426 
427  G4String PN("Name");
428  (*store)[PN] = G4AttDef(PN,"Particle Name","Physics","","G4String");
429 
430  G4String Ch("Q");
431  (*store)[Ch] = G4AttDef(Ch,"Charge","Physics","e+","G4double");
432 
433  G4String PDG("PDG");
434  (*store)[PDG] = G4AttDef(PDG,"PDG Encoding","Physics","","G4int");
435 
436  G4String IMom("IMom");
437  (*store)[IMom] = G4AttDef(IMom, "Momentum at start of trajectory",
438  "Physics","G4BestUnit","G4ThreeVector");
439 
440  G4String IMag("IMomMag");
441  (*store)[IMag] =
442  G4AttDef(IMag, "Magnitude of momentum at start of trajectory",
443  "Physics","G4BestUnit","G4double");
444 
445  G4String NTP("NPts");
446  (*store)[NTP] = G4AttDef(NTP,"No. of points","Physics","","G4int");
447 
448  }
449  return store;
450 }
unsigned int ID
const uint PDG
Definition: qregexp.cpp:140
G4double LBNETrajectory::GetCharge ( ) const
inline

Definition at line 52 of file LBNETrajectory.hh.

53  { return fPDGCharge; }
G4int LBNETrajectory::GetDecayCode ( ) const
inline

Definition at line 72 of file LBNETrajectory.hh.

73  { return fDecayCode;}
G4ThreeVector LBNETrajectory::GetInitialMomentum ( ) const
inline

Definition at line 68 of file LBNETrajectory.hh.

69  { return fMomentum; }
G4ThreeVector fMomentum
G4double LBNETrajectory::GetMass ( ) const
inline

Definition at line 54 of file LBNETrajectory.hh.

55  { return fParticleMass; }
G4double fParticleMass
virtual G4String LBNETrajectory::GetMaterialName ( G4int  i) const
inlinevirtual

Definition at line 106 of file LBNETrajectory.hh.

107  {return (*fMaterialName)[i];}
LBNETrajectoryMaterialName * fMaterialName
G4String LBNETrajectory::GetMaterialName1rst ( ) const
inline

Definition at line 104 of file LBNETrajectory.hh.

105  { return fMaterialName1rst; }
G4String fMaterialName1rst
G4int LBNETrajectory::GetMaterialNumber1rst ( ) const
inline

Definition at line 79 of file LBNETrajectory.hh.

80  {return fMaterialNumber1rst;}
G4int LBNETrajectory::GetMaterialNumberLast ( ) const
inline

Definition at line 84 of file LBNETrajectory.hh.

85  {return fMaterialNumberLast;}
virtual G4ThreeVector LBNETrajectory::GetMomentum ( G4int  i) const
inlinevirtual

Definition at line 64 of file LBNETrajectory.hh.

65  { return (*fMomentumRecord)[i]; }
LBNETrajectoryMomentumContainer * fMomentumRecord
virtual G4double LBNETrajectory::GetNImpWt ( ) const
inlinevirtual

Definition at line 74 of file LBNETrajectory.hh.

75  { return fNImpWt;}
G4int LBNETrajectory::GetParentID ( ) const
inline

Definition at line 48 of file LBNETrajectory.hh.

49  { return fParentID; }
const G4ThreeVector LBNETrajectory::GetParentMomentumAtThisProduction ( ) const
inline

Definition at line 108 of file LBNETrajectory.hh.

108  {
110  }
G4ThreeVector fParentMomentumAtThisProduction
G4ParticleDefinition * LBNETrajectory::GetParticleDefinition ( )

Definition at line 509 of file LBNETrajectory.cc.

510 {
511  return (G4ParticleTable::GetParticleTable()->FindParticle(fParticleName));
512 }
G4String fParticleName
G4String LBNETrajectory::GetParticleName ( ) const
inline

Definition at line 50 of file LBNETrajectory.hh.

51  { return fParticleName; }
G4String fParticleName
G4int LBNETrajectory::GetPDGEncoding ( ) const
inline

Definition at line 56 of file LBNETrajectory.hh.

57  { return fPDGEncoding; }
G4int LBNETrajectory::GetPDGNucleus ( ) const
inline

Definition at line 98 of file LBNETrajectory.hh.

98  {
99  return fPDGNucleus; }
virtual G4VTrajectoryPoint* LBNETrajectory::GetPoint ( G4int  i) const
inlinevirtual

Definition at line 62 of file LBNETrajectory.hh.

63  { return (*fPositionRecord)[i]; }
LBNETrajectoryPointContainer * fPositionRecord
virtual int LBNETrajectory::GetPointEntries ( ) const
inlinevirtual

Definition at line 60 of file LBNETrajectory.hh.

61  { return fPositionRecord->size(); }
LBNETrajectoryPointContainer * fPositionRecord
G4ThreeVector LBNETrajectory::GetPolarization ( ) const
inline

Definition at line 95 of file LBNETrajectory.hh.

96  { return fPolarization;}
G4ThreeVector fPolarization
virtual G4String LBNETrajectory::GetPreStepVolumeName ( G4int  i) const
inlinevirtual

Definition at line 66 of file LBNETrajectory.hh.

67  { return (*fPreStepVolume)[i]; }
LBNETrajectoryVolumeName * fPreStepVolume
G4String LBNETrajectory::GetProcessName ( ) const
inline

Definition at line 101 of file LBNETrajectory.hh.

102  { return fProcessName; }
G4String fProcessName
virtual G4double LBNETrajectory::GetStepLength ( G4int  i) const
inlinevirtual

Definition at line 76 of file LBNETrajectory.hh.

77  {return (*fStepLength)[i];}
virtual G4int LBNETrajectory::GetTgen ( ) const
inlinevirtual

Definition at line 70 of file LBNETrajectory.hh.

71  { return fTgen;}
G4double LBNETrajectory::GetTimeStart ( ) const
inline

Definition at line 92 of file LBNETrajectory.hh.

93  {return fTimeStart;}
G4int LBNETrajectory::GetTrackID ( ) const
inline

Definition at line 46 of file LBNETrajectory.hh.

47  { return fTrackID; }
const G4ThreeVector& LBNETrajectory::GetVertexPosition ( ) const
inline

Definition at line 58 of file LBNETrajectory.hh.

59  { return fVertexPosition; }
G4ThreeVector fVertexPosition
G4String LBNETrajectory::GetVolName1rst ( ) const
inline

Definition at line 89 of file LBNETrajectory.hh.

90  {return fVolName1rst;}
G4String fVolName1rst
void LBNETrajectory::MergeTrajectory ( G4VTrajectory *  secondTrajectory)
virtual

Definition at line 514 of file LBNETrajectory.cc.

515 {
516  if(!secondTrajectory) return;
517 
518  LBNETrajectory* seco = (LBNETrajectory*)secondTrajectory;
519  G4int ent = seco->GetPointEntries();
520  for(int i=1;i<ent;i++) // initial point of the second trajectory should not be merged
521  {
522  fPositionRecord->push_back((*(seco->fPositionRecord))[i]);
523  }
524  delete (*seco->fPositionRecord)[0];
525  seco->fPositionRecord->clear();
526 
527 }
LBNETrajectoryPointContainer * fPositionRecord
virtual int GetPointEntries() const
void LBNETrajectory::operator delete ( void *  aTrajectory)
inline

Definition at line 165 of file LBNETrajectory.hh.

166 {
167  myTrajectoryAllocator.FreeSingle((LBNETrajectory*)aTrajectory);
168 }
G4Allocator< LBNETrajectory > myTrajectoryAllocator
void * LBNETrajectory::operator new ( size_t  )
inline

Definition at line 158 of file LBNETrajectory.hh.

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

Definition at line 43 of file LBNETrajectory.hh.

44  {return (this==&right);}
void LBNETrajectory::SetMaterialNumber1rst ( G4double  matNum)
inline

Definition at line 81 of file LBNETrajectory.hh.

82  {fMaterialNumber1rst = matNum;}
void LBNETrajectory::SetMaterialNumberLast ( G4double  matNum)
inline

Definition at line 86 of file LBNETrajectory.hh.

87  {fMaterialNumberLast = matNum;}
void LBNETrajectory::ShowTrajectory ( ) const
virtual

Definition at line 290 of file LBNETrajectory.cc.

291 {
292  G4cout << G4endl << "TrackID =" << fTrackID
293  << " : ParentID=" << fParentID << G4endl;
294  G4cout << "Particle name : " << fParticleName
295  << " Charge : " << fPDGCharge << G4endl;
296  G4cout << "Original momentum : " <<
297 G4BestUnit(fMomentum,"Energy") << G4endl;
298  G4cout << "Vertex : " << G4BestUnit(fVertexPosition,"Length") << G4endl;
299  G4cout << " Current trajectory has " << fPositionRecord->size()
300  << " points." << G4endl;
301 
302  for( size_t i=0 ; i < fPositionRecord->size() ; i++){
303  G4TrajectoryPoint* aTrajectoryPoint = (G4TrajectoryPoint*)((*fPositionRecord)[i]);
304  G4cout << "Point[" << i << "]"
305  << " Position= " << aTrajectoryPoint->GetPosition() << G4endl;
306  }
307 }
G4ThreeVector fVertexPosition
LBNETrajectoryPointContainer * fPositionRecord
G4ThreeVector fMomentum
G4String fParticleName
void LBNETrajectory::ShowTrajectory ( std::ostream &  o) const
virtual

Definition at line 309 of file LBNETrajectory.cc.

310 {
311  G4VTrajectory::ShowTrajectory(o);
312 }

Member Data Documentation

G4int LBNETrajectory::fDecayCode
private

Definition at line 126 of file LBNETrajectory.hh.

G4int LBNETrajectory::fEventNo
private

Definition at line 127 of file LBNETrajectory.hh.

LBNETrajectoryMaterialName* LBNETrajectory::fMaterialName
private

Definition at line 151 of file LBNETrajectory.hh.

G4String LBNETrajectory::fMaterialName1rst
private

Definition at line 152 of file LBNETrajectory.hh.

G4int LBNETrajectory::fMaterialNumber1rst
private

Definition at line 142 of file LBNETrajectory.hh.

G4int LBNETrajectory::fMaterialNumberLast
private

Definition at line 144 of file LBNETrajectory.hh.

G4ThreeVector LBNETrajectory::fMomentum
private

Definition at line 136 of file LBNETrajectory.hh.

LBNETrajectoryMomentumContainer* LBNETrajectory::fMomentumRecord
private

Definition at line 124 of file LBNETrajectory.hh.

G4double LBNETrajectory::fNImpWt
private

Definition at line 129 of file LBNETrajectory.hh.

G4int LBNETrajectory::fParentID
private

Definition at line 131 of file LBNETrajectory.hh.

G4ThreeVector LBNETrajectory::fParentMomentumAtThisProduction
private

Definition at line 150 of file LBNETrajectory.hh.

G4ParticleDefinition* LBNETrajectory::fParticleDefinition
private

Definition at line 132 of file LBNETrajectory.hh.

G4double LBNETrajectory::fParticleMass
private

Definition at line 138 of file LBNETrajectory.hh.

G4String LBNETrajectory::fParticleName
private

Definition at line 133 of file LBNETrajectory.hh.

G4double LBNETrajectory::fPDGCharge
private

Definition at line 134 of file LBNETrajectory.hh.

G4int LBNETrajectory::fPDGEncoding
private

Definition at line 135 of file LBNETrajectory.hh.

G4int LBNETrajectory::fPDGNucleus
private

Definition at line 147 of file LBNETrajectory.hh.

G4ThreeVector LBNETrajectory::fPolarization
private

Definition at line 146 of file LBNETrajectory.hh.

LBNETrajectoryPointContainer* LBNETrajectory::fPositionRecord
private

Definition at line 125 of file LBNETrajectory.hh.

LBNETrajectoryVolumeName* LBNETrajectory::fPreStepVolume
private

Definition at line 139 of file LBNETrajectory.hh.

G4String LBNETrajectory::fProcessName
private

Definition at line 149 of file LBNETrajectory.hh.

DVec* LBNETrajectory::fStepLength
private

Definition at line 140 of file LBNETrajectory.hh.

G4int LBNETrajectory::fTgen
private

Definition at line 128 of file LBNETrajectory.hh.

G4double LBNETrajectory::fTimeStart
private

Definition at line 145 of file LBNETrajectory.hh.

G4int LBNETrajectory::fTrackID
private

Definition at line 130 of file LBNETrajectory.hh.

G4ThreeVector LBNETrajectory::fVertexPosition
private

Definition at line 137 of file LBNETrajectory.hh.

G4String LBNETrajectory::fVolName1rst
private

Definition at line 143 of file LBNETrajectory.hh.


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