G4RichTrajectory Class Reference

#include <G4RichTrajectory.hh>

Inheritance diagram for G4RichTrajectory:

G4Trajectory G4VTrajectory

Public Member Functions

 G4RichTrajectory ()
 G4RichTrajectory (const G4Track *aTrack)
 G4RichTrajectory (G4RichTrajectory &)
virtual ~G4RichTrajectory ()
void * operator new (size_t)
void operator delete (void *)
int operator== (const G4RichTrajectory &right) const
void AppendStep (const G4Step *aStep)
void MergeTrajectory (G4VTrajectory *secondTrajectory)
int GetPointEntries () const
G4VTrajectoryPointGetPoint (G4int i) const
virtual const std::map< G4String,
G4AttDef > * 
GetAttDefs () const
virtual std::vector< G4AttValue > * CreateAttValues () const

Detailed Description

Definition at line 65 of file G4RichTrajectory.hh.


Constructor & Destructor Documentation

G4RichTrajectory::G4RichTrajectory (  ) 

Definition at line 62 of file G4RichTrajectory.cc.

00062                                   :
00063   fpRichPointsContainer(0),
00064   fpCreatorProcess(0),
00065   fpEndingProcess(0),
00066   fFinalKineticEnergy(0.)
00067 {}

G4RichTrajectory::G4RichTrajectory ( const G4Track aTrack  ) 

Definition at line 69 of file G4RichTrajectory.cc.

References G4Track::GetCreatorProcess(), G4Track::GetKineticEnergy(), G4Track::GetNextTouchableHandle(), and G4Track::GetTouchableHandle().

00069                                                        :
00070   G4Trajectory(aTrack)  // Note: this initialises the base class data
00071                         // members and, unfortunately but never mind,
00072                         // creates a G4TrajectoryPoint in
00073                         // TrajectoryPointContainer that we cannot
00074                         // access because it's private.  We store the
00075                         // same information (plus more) in a
00076                         // G4RichTrajectoryPoint in the
00077                         // RichTrajectoryPointsContainer
00078 {
00079   fpInitialVolume = aTrack->GetTouchableHandle();
00080   fpInitialNextVolume = aTrack->GetNextTouchableHandle();
00081   fpCreatorProcess = aTrack->GetCreatorProcess();
00082   // On construction, set final values to initial values.
00083   // Final values are updated at the addition of every step - see AppendStep.
00084   fpFinalVolume = aTrack->GetTouchableHandle();
00085   fpFinalNextVolume = aTrack->GetNextTouchableHandle();
00086   fpEndingProcess = aTrack->GetCreatorProcess();
00087   fFinalKineticEnergy = aTrack->GetKineticEnergy();
00088   // Insert the first rich trajectory point (see note above)...
00089   fpRichPointsContainer = new RichTrajectoryPointsContainer;
00090   fpRichPointsContainer->push_back(new G4RichTrajectoryPoint(aTrack));
00091 }

G4RichTrajectory::G4RichTrajectory ( G4RichTrajectory  ) 

Definition at line 93 of file G4RichTrajectory.cc.

References fFinalKineticEnergy, fpCreatorProcess, fpEndingProcess, fpFinalNextVolume, fpFinalVolume, fpInitialNextVolume, fpInitialVolume, and fpRichPointsContainer.

00093                                                           :
00094   G4Trajectory(right)
00095 {
00096   fpInitialVolume = right.fpInitialVolume;
00097   fpInitialNextVolume = right.fpInitialNextVolume;
00098   fpCreatorProcess = right.fpCreatorProcess;
00099   fpFinalVolume = right.fpFinalVolume;
00100   fpFinalNextVolume = right.fpFinalNextVolume;
00101   fpEndingProcess = right.fpEndingProcess;
00102   fFinalKineticEnergy = right.fFinalKineticEnergy;
00103   fpRichPointsContainer = new RichTrajectoryPointsContainer;
00104   for(size_t i=0;i<right.fpRichPointsContainer->size();i++)
00105   {
00106     G4RichTrajectoryPoint* rightPoint =
00107       (G4RichTrajectoryPoint*)((*(right.fpRichPointsContainer))[i]);
00108     fpRichPointsContainer->push_back(new G4RichTrajectoryPoint(*rightPoint));
00109   }
00110 }

G4RichTrajectory::~G4RichTrajectory (  )  [virtual]

Definition at line 112 of file G4RichTrajectory.cc.

00113 {
00114   if (fpRichPointsContainer) {
00115     //  fpRichPointsContainer->clearAndDestroy();
00116     size_t i;
00117     for(i=0;i<fpRichPointsContainer->size();i++){
00118       delete  (*fpRichPointsContainer)[i];
00119     }
00120     fpRichPointsContainer->clear();
00121     delete fpRichPointsContainer;
00122   }
00123 }


Member Function Documentation

void G4RichTrajectory::AppendStep ( const G4Step aStep  )  [virtual]

Reimplemented from G4Trajectory.

Definition at line 125 of file G4RichTrajectory.cc.

References G4Track::GetCurrentStepNumber(), G4Track::GetNextTouchableHandle(), G4StepPoint::GetProcessDefinedStep(), and G4Track::GetTouchableHandle().

00126 {
00127   fpRichPointsContainer->push_back(new G4RichTrajectoryPoint(aStep));
00128   // Except for first step, which is a sort of virtual step to start
00129   // the track, compute the final values...
00130   const G4Track* track = aStep->GetTrack();
00131   const G4StepPoint* postStepPoint = aStep->GetPostStepPoint();
00132   if (track->GetCurrentStepNumber() > 0) {
00133     fpFinalVolume = track->GetTouchableHandle();
00134     fpFinalNextVolume = track->GetNextTouchableHandle();
00135     fpEndingProcess = postStepPoint->GetProcessDefinedStep();
00136     fFinalKineticEnergy =
00137       aStep->GetPreStepPoint()->GetKineticEnergy() -
00138       aStep->GetTotalEnergyDeposit();
00139   }
00140 }

std::vector< G4AttValue > * G4RichTrajectory::CreateAttValues (  )  const [virtual]

Reimplemented from G4Trajectory.

Definition at line 222 of file G4RichTrajectory.cc.

References G4Trajectory::CreateAttValues(), G4BestUnit, G4cout, GetAttDefs(), G4VProcess::GetProcessName(), G4VProcess::GetProcessType(), and G4VProcess::GetProcessTypeName().

00223 {
00224   // Create base class att values...
00225   std::vector<G4AttValue>* values = G4Trajectory::CreateAttValues();
00226 
00227   if (fpInitialVolume && fpInitialVolume->GetVolume()) {
00228     values->push_back(G4AttValue("IVPath",Path(fpInitialVolume),""));
00229   } else {
00230     values->push_back(G4AttValue("IVPath","None",""));
00231   }
00232 
00233   if (fpInitialNextVolume && fpInitialNextVolume->GetVolume()) {
00234     values->push_back(G4AttValue("INVPath",Path(fpInitialNextVolume),""));
00235   } else {
00236     values->push_back(G4AttValue("INVPath","None",""));
00237   }
00238 
00239   if (fpCreatorProcess) {
00240     values->push_back(G4AttValue("CPN",fpCreatorProcess->GetProcessName(),""));
00241     G4ProcessType type = fpCreatorProcess->GetProcessType();
00242     values->push_back(G4AttValue("CPTN",G4VProcess::GetProcessTypeName(type),""));
00243   } else {
00244     values->push_back(G4AttValue("CPN","None",""));
00245     values->push_back(G4AttValue("CPTN","None",""));
00246   }
00247 
00248   if (fpFinalVolume && fpFinalVolume->GetVolume()) {
00249     values->push_back(G4AttValue("FVPath",Path(fpFinalVolume),""));
00250   } else {
00251     values->push_back(G4AttValue("FVPath","None",""));
00252   }
00253 
00254   if (fpFinalNextVolume && fpFinalNextVolume->GetVolume()) {
00255     values->push_back(G4AttValue("FNVPath",Path(fpFinalNextVolume),""));
00256   } else {
00257     values->push_back(G4AttValue("FNVPath","None",""));
00258   }
00259 
00260   if (fpEndingProcess) {
00261     values->push_back(G4AttValue("EPN",fpEndingProcess->GetProcessName(),""));
00262     G4ProcessType type = fpEndingProcess->GetProcessType();
00263     values->push_back(G4AttValue("EPTN",G4VProcess::GetProcessTypeName(type),""));
00264   } else {
00265     values->push_back(G4AttValue("EPN","None",""));
00266     values->push_back(G4AttValue("EPTN","None",""));
00267   }
00268 
00269   values->push_back
00270     (G4AttValue("FKE",G4BestUnit(fFinalKineticEnergy,"Energy"),""));
00271 
00272 #ifdef G4ATTDEBUG
00273   G4cout << G4AttCheck(values,GetAttDefs());
00274 #endif
00275 
00276   return values;
00277 }

const std::map< G4String, G4AttDef > * G4RichTrajectory::GetAttDefs (  )  const [virtual]

Reimplemented from G4Trajectory.

Definition at line 157 of file G4RichTrajectory.cc.

References G4Trajectory::GetAttDefs(), and G4AttDefStore::GetInstance().

Referenced by CreateAttValues().

00158 {
00159   G4bool isNew;
00160   std::map<G4String,G4AttDef>* store
00161     = G4AttDefStore::GetInstance("G4RichTrajectory",isNew);
00162   if (isNew) {
00163 
00164     // Get att defs from base class...
00165     *store = *(G4Trajectory::GetAttDefs());
00166 
00167     G4String ID;
00168 
00169     ID = "IVPath";
00170     (*store)[ID] = G4AttDef(ID,"Initial Volume Path",
00171                             "Physics","","G4String");
00172 
00173     ID = "INVPath";
00174     (*store)[ID] = G4AttDef(ID,"Initial Next Volume Path",
00175                             "Physics","","G4String");
00176 
00177     ID = "CPN";
00178     (*store)[ID] = G4AttDef(ID,"Creator Process Name",
00179                             "Physics","","G4String");
00180 
00181     ID = "CPTN";
00182     (*store)[ID] = G4AttDef(ID,"Creator Process Type Name",
00183                             "Physics","","G4String");
00184 
00185     ID = "FVPath";
00186     (*store)[ID] = G4AttDef(ID,"Final Volume Path",
00187                             "Physics","","G4String");
00188 
00189     ID = "FNVPath";
00190     (*store)[ID] = G4AttDef(ID,"Final Next Volume Path",
00191                             "Physics","","G4String");
00192 
00193     ID = "EPN";
00194     (*store)[ID] = G4AttDef(ID,"Ending Process Name",
00195                             "Physics","","G4String");
00196 
00197     ID = "EPTN";
00198     (*store)[ID] = G4AttDef(ID,"Ending Process Type Name",
00199                             "Physics","","G4String");
00200 
00201     ID = "FKE";
00202     (*store)[ID] = G4AttDef(ID,"Final kinetic energy",
00203                             "Physics","G4BestUnit","G4double");
00204 
00205   }
00206 
00207   return store;
00208 }

G4VTrajectoryPoint* G4RichTrajectory::GetPoint ( G4int  i  )  const [inline, virtual]

Reimplemented from G4Trajectory.

Definition at line 89 of file G4RichTrajectory.hh.

00090    { return (*fpRichPointsContainer)[i]; }

int G4RichTrajectory::GetPointEntries (  )  const [inline, virtual]

Reimplemented from G4Trajectory.

Definition at line 88 of file G4RichTrajectory.hh.

Referenced by MergeTrajectory().

00088 { return fpRichPointsContainer->size(); }

void G4RichTrajectory::MergeTrajectory ( G4VTrajectory secondTrajectory  )  [virtual]

Reimplemented from G4Trajectory.

Definition at line 142 of file G4RichTrajectory.cc.

References fpRichPointsContainer, and GetPointEntries().

00143 {
00144   if(!secondTrajectory) return;
00145 
00146   G4RichTrajectory* seco = (G4RichTrajectory*)secondTrajectory;
00147   G4int ent = seco->GetPointEntries();
00148   for(G4int i=1;i<ent;i++) {
00149     // initial point of the second trajectory should not be merged
00150     fpRichPointsContainer->push_back((*(seco->fpRichPointsContainer))[i]);
00151     //    fpRichPointsContainer->push_back(seco->fpRichPointsContainer->removeAt(1));
00152   }
00153   delete (*seco->fpRichPointsContainer)[0];
00154   seco->fpRichPointsContainer->clear();
00155 }

void G4RichTrajectory::operator delete ( void *   )  [inline]

Reimplemented from G4Trajectory.

Definition at line 125 of file G4RichTrajectory.hh.

References aRichTrajectoryAllocator.

00126 {
00127   aRichTrajectoryAllocator.FreeSingle
00128     ((G4RichTrajectory*)aRichTrajectory);
00129 }

void * G4RichTrajectory::operator new ( size_t   )  [inline]

Reimplemented from G4Trajectory.

Definition at line 118 of file G4RichTrajectory.hh.

References aRichTrajectoryAllocator.

00119 {
00120   void* aRichTrajectory;
00121   aRichTrajectory = (void*)aRichTrajectoryAllocator.MallocSingle();
00122   return aRichTrajectory;
00123 }

int G4RichTrajectory::operator== ( const G4RichTrajectory right  )  const [inline]

Definition at line 82 of file G4RichTrajectory.hh.

00083   {return (this==&right);} 


The documentation for this class was generated from the following files:
Generated on Mon May 27 17:53:17 2013 for Geant4 by  doxygen 1.4.7