Geant4  10.02.p03
CML2SDWithParticle Class Reference

#include <ML2SDWithParticle.hh>

Inheritance diagram for CML2SDWithParticle:
Collaboration diagram for CML2SDWithParticle:

Public Member Functions

 CML2SDWithParticle ()
 
 CML2SDWithParticle (G4int idType, G4int max_N_particles_in_PhSp_File, G4int seed, G4int nMaxParticlesInRamPhaseSpace, G4String name, G4String PhaseSpaceOutFile, G4bool bSavePhaseSpace, G4bool bStopAtVolatilePhaseSpace, SPrimaryParticle *primaryParticleData, G4double accTargetZPosition)
 
 ~CML2SDWithParticle (void)
 
G4bool ProcessHits (G4Step *aStep, G4TouchableHistory *ROHist)
 
G4int getTotalNumberOfParticles ()
 
CML2SDWithParticlegetCML2SensitiveDetectorParticle ()
 
Sparticle getParticle (int i)
 
void setActive (G4bool act)
 
void save ()
 
- Public Member Functions inherited from G4VSensitiveDetector
 G4VSensitiveDetector (G4String name)
 
 G4VSensitiveDetector (const G4VSensitiveDetector &right)
 
virtual ~G4VSensitiveDetector ()
 
G4VSensitiveDetectoroperator= (const G4VSensitiveDetector &right)
 
G4int operator== (const G4VSensitiveDetector &right) const
 
G4int operator!= (const G4VSensitiveDetector &right) const
 
virtual void Initialize (G4HCofThisEvent *)
 
virtual void EndOfEvent (G4HCofThisEvent *)
 
virtual void clear ()
 
virtual void DrawAll ()
 
virtual void PrintAll ()
 
G4bool Hit (G4Step *aStep)
 
void SetROgeometry (G4VReadOutGeometry *value)
 
void SetFilter (G4VSDFilter *value)
 
G4int GetNumberOfCollections () const
 
G4String GetCollectionName (G4int id) const
 
void SetVerboseLevel (G4int vl)
 
void Activate (G4bool activeFlag)
 
G4bool isActive () const
 
G4String GetName () const
 
G4String GetPathName () const
 
G4String GetFullPathName () const
 
G4VReadOutGeometryGetROgeometry () const
 
G4VSDFilterGetFilter () const
 
virtual G4VSensitiveDetectorClone () const
 

Private Member Functions

void saveDataParticles (G4int nParticle)
 
void saveHeaderParticles ()
 

Private Attributes

G4ThreeVector halfSize
 
G4ThreeVector pos
 
SPrimaryParticleprimaryParticleData
 
Sparticleparticles
 
G4String fullOutFileData
 
G4int nTotalParticles
 
G4int nParticle
 
G4int idType
 
G4int nMaxParticlesInRamPhaseSpace
 
G4int max_N_particles_in_PhSp_File
 
G4bool bStopAtPhaseSpace
 
G4bool bSavePhaseSpace
 
G4bool bActive
 
G4double accTargetZPosition
 

Additional Inherited Members

- Protected Member Functions inherited from G4VSensitiveDetector
virtual G4int GetCollectionID (G4int i)
 
- Protected Attributes inherited from G4VSensitiveDetector
G4CollectionNameVector collectionName
 
G4String SensitiveDetectorName
 
G4String thePathName
 
G4String fullPathName
 
G4int verboseLevel
 
G4bool active
 
G4VReadOutGeometryROgeometry
 
G4VSDFilterfilter
 

Detailed Description

Definition at line 55 of file ML2SDWithParticle.hh.

Constructor & Destructor Documentation

◆ CML2SDWithParticle() [1/2]

CML2SDWithParticle::CML2SDWithParticle ( )

Definition at line 48 of file ML2SDWithParticle.cc.

◆ CML2SDWithParticle() [2/2]

CML2SDWithParticle::CML2SDWithParticle ( G4int  idType,
G4int  max_N_particles_in_PhSp_File,
G4int  seed,
G4int  nMaxParticlesInRamPhaseSpace,
G4String  name,
G4String  PhaseSpaceOutFile,
G4bool  bSavePhaseSpace,
G4bool  bStopAtVolatilePhaseSpace,
SPrimaryParticle primaryParticleData,
G4double  accTargetZPosition 
)

Definition at line 57 of file ML2SDWithParticle.cc.

62 {
63  accTargetZPosition=ZPos;
64  max_N_particles_in_PhSp_File=maxPartFile;
66  idType=id;
67  primaryParticleData=pData;
68  bActive=true;
69  nParticle=0;
71  bSavePhaseSpace=bSave;
72  bStopAtPhaseSpace=bStop;
73  if (bSavePhaseSpace)
74  {particles=new Sparticle[nMaxPart];}
75  G4String seedName;
76  char a[10];
77  sprintf(a,"%d", seed);
78  seedName=(G4String)a;
79 
80  fullOutFileData=PhaseSpaceOutFile+"_"+seedName+".txt";
81  fullOutFileData=PhaseSpaceOutFile+"_"+seedName+".txt";
82 }
G4VSensitiveDetector(G4String name)

◆ ~CML2SDWithParticle()

CML2SDWithParticle::~CML2SDWithParticle ( void  )

Definition at line 83 of file ML2SDWithParticle.cc.

84 {
85  if (particles!=0)
86  {
87  delete [] particles;
88  }
89 }

Member Function Documentation

◆ getCML2SensitiveDetectorParticle()

CML2SDWithParticle* CML2SDWithParticle::getCML2SensitiveDetectorParticle ( )
inline

Definition at line 63 of file ML2SDWithParticle.hh.

63 {return this;}
Here is the caller graph for this function:

◆ getParticle()

Sparticle CML2SDWithParticle::getParticle ( int  i)
inline

Definition at line 64 of file ML2SDWithParticle.hh.

64 {return particles[i];}

◆ getTotalNumberOfParticles()

G4int CML2SDWithParticle::getTotalNumberOfParticles ( )
inline

Definition at line 62 of file ML2SDWithParticle.hh.

62 {return nTotalParticles;}
Here is the caller graph for this function:

◆ ProcessHits()

G4bool CML2SDWithParticle::ProcessHits ( G4Step *  aStep,
G4TouchableHistory ROHist 
)
virtual

Implements G4VSensitiveDetector.

Definition at line 121 of file ML2SDWithParticle.cc.

122 {
123  if (bActive && (CML2AcceleratorConstruction::GetInstance()->getPhysicalVolume()->GetRotation()->isIdentity()))
124  {
125  G4double energyKin= aStep->GetTrack()->GetKineticEnergy();
126  static bool bFirstTime=true;
128  {
129  nTotalParticles++;
130  aStep->GetTrack()->SetTrackStatus(fStopAndKill);
131  }
132  else
133  {
134  if (energyKin>0.)
135  {
137  particles[nParticle].pos=aStep->GetPreStepPoint()->GetPosition();
138  particles[nParticle].dir=aStep->GetPreStepPoint()->GetMomentumDirection();
139  particles[nParticle].kinEnergy=energyKin;
140  particles[nParticle].partPDGE=aStep->GetTrack()->GetDefinition()->GetPDGEncoding();
143  nParticle++;
144  nTotalParticles++;
146  {
147  if (bFirstTime)
148  {
149  bFirstTime=false;
151  }
153  nParticle=0;
154  bActive =false;// to stop the phase space creation
155  }
157  {
158  if (bFirstTime)
159  {
160  bFirstTime=false;
162  }
164  nParticle=0;
165  }
166 
167  Sparticle *particle=new Sparticle;
168  particle->dir=aStep->GetPreStepPoint()->GetMomentumDirection();
169  particle->pos=aStep->GetPreStepPoint()->GetPosition();
170  particle->kinEnergy=energyKin;
171  particle->nPrimaryPart=nTotalParticles; // to pass the id of this phase space particle
172  particle->partPDGE=aStep->GetTrack()->GetDefinition()->GetPDGEncoding();
174  particle->volumeId=-1;
175  particle->volumeName="-1";
176  }
177  if (bStopAtPhaseSpace)
178  {aStep->GetTrack()->SetTrackStatus(fStopAndKill);}
179  }
180  }
181  else
182  {
183  if (bStopAtPhaseSpace)
184  {aStep->GetTrack()->SetTrackStatus(fStopAndKill);}
185  }
186 
187  return true;
188 }
SPrimaryParticle * primaryParticleData
static CML2AcceleratorConstruction * GetInstance(void)
G4int partPDGE
G4ThreeVector dir
G4int nPrimaryPart
G4int primaryParticlePDGE
G4ThreeVector pos
G4double kinEnergy
G4String volumeName
G4int volumeId
double G4double
Definition: G4Types.hh:76
void saveDataParticles(G4int nParticle)
Here is the call graph for this function:

◆ save()

void CML2SDWithParticle::save ( )

Definition at line 190 of file ML2SDWithParticle.cc.

191 {
192  if ((bActive) && (nParticle>0) && (CML2AcceleratorConstruction::GetInstance()->getPhysicalVolume()->GetRotation()->isIdentity()))
194 }
static CML2AcceleratorConstruction * GetInstance(void)
void saveDataParticles(G4int nParticle)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ saveDataParticles()

void CML2SDWithParticle::saveDataParticles ( G4int  nParticle)
private

Definition at line 98 of file ML2SDWithParticle.cc.

99 {
100  std::ofstream out;
101  out.open(fullOutFileData, std::ios::app);
102  static G4int nTotParticles=0;
103  for (int i=0; i< nPart; i++)
104  {
105  out << nTotParticles++ << '\t';
106 // out << particles[i].volumeName << '\t';
107  out << particles[i].pos.getX()/mm << '\t';
108  out << particles[i].pos.getY()/mm << '\t';
109  out << (accTargetZPosition + particles[i].pos.getZ())/mm << '\t'; // it translates the current z value in global coordinates to the accelerator local coordinates (only z)
110  out << particles[i].dir.getX() << '\t';
111  out << particles[i].dir.getY() << '\t';
112  out << particles[i].dir.getZ() << '\t';
113  out << particles[i].kinEnergy/MeV << '\t';
114  out << particles[i].partPDGE<< '\t';
115  out << particles[i].primaryParticlePDGE<< '\t';
116  out << particles[i].nPrimaryPart<< G4endl;
117  }
118  out.close();
119 }
static const double MeV
Definition: G4SIunits.hh:211
G4int partPDGE
G4ThreeVector dir
int G4int
Definition: G4Types.hh:78
G4int nPrimaryPart
double getY() const
G4int primaryParticlePDGE
double getX() const
double getZ() const
G4ThreeVector pos
G4double kinEnergy
#define G4endl
Definition: G4ios.hh:61
static const double mm
Definition: G4SIunits.hh:114
Here is the call graph for this function:
Here is the caller graph for this function:

◆ saveHeaderParticles()

void CML2SDWithParticle::saveHeaderParticles ( )
private

Definition at line 90 of file ML2SDWithParticle.cc.

91 {
92  std::ofstream out;
93  out.open(fullOutFileData, std::ios::out);
94  out << "Sensitive Detector-Particles"<<G4endl;
95  out << "n Total Events,\t x [mm],\t y [mm],\t z [mm],\t dirX,\t dirY,\t dirZ,\t KinEnergy [MeV],\t part Type,\t primary part type,\t nPrimaryPart" << G4endl;
96  out.close();
97 }
#define G4endl
Definition: G4ios.hh:61
Here is the caller graph for this function:

◆ setActive()

void CML2SDWithParticle::setActive ( G4bool  act)
inline

Definition at line 65 of file ML2SDWithParticle.hh.

65 {bActive=act;}
Here is the call graph for this function:

Member Data Documentation

◆ accTargetZPosition

G4double CML2SDWithParticle::accTargetZPosition
private

Definition at line 80 of file ML2SDWithParticle.hh.

◆ bActive

G4bool CML2SDWithParticle::bActive
private

Definition at line 78 of file ML2SDWithParticle.hh.

◆ bSavePhaseSpace

G4bool CML2SDWithParticle::bSavePhaseSpace
private

Definition at line 78 of file ML2SDWithParticle.hh.

◆ bStopAtPhaseSpace

G4bool CML2SDWithParticle::bStopAtPhaseSpace
private

Definition at line 78 of file ML2SDWithParticle.hh.

◆ fullOutFileData

G4String CML2SDWithParticle::fullOutFileData
private

Definition at line 75 of file ML2SDWithParticle.hh.

◆ halfSize

G4ThreeVector CML2SDWithParticle::halfSize
private

Definition at line 71 of file ML2SDWithParticle.hh.

◆ idType

G4int CML2SDWithParticle::idType
private

Definition at line 77 of file ML2SDWithParticle.hh.

◆ max_N_particles_in_PhSp_File

G4int CML2SDWithParticle::max_N_particles_in_PhSp_File
private

Definition at line 77 of file ML2SDWithParticle.hh.

◆ nMaxParticlesInRamPhaseSpace

G4int CML2SDWithParticle::nMaxParticlesInRamPhaseSpace
private

Definition at line 77 of file ML2SDWithParticle.hh.

◆ nParticle

G4int CML2SDWithParticle::nParticle
private

Definition at line 76 of file ML2SDWithParticle.hh.

◆ nTotalParticles

G4int CML2SDWithParticle::nTotalParticles
private

Definition at line 76 of file ML2SDWithParticle.hh.

◆ particles

Sparticle* CML2SDWithParticle::particles
private

Definition at line 74 of file ML2SDWithParticle.hh.

◆ pos

G4ThreeVector CML2SDWithParticle::pos
private

Definition at line 72 of file ML2SDWithParticle.hh.

◆ primaryParticleData

SPrimaryParticle* CML2SDWithParticle::primaryParticleData
private

Definition at line 73 of file ML2SDWithParticle.hh.


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