Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4GPSModel Class Reference

#include <G4GPSModel.hh>

+ Inheritance diagram for G4GPSModel:

Public Member Functions

 G4GPSModel (const G4Colour &colour=G4Colour(1., 0., 0., 0.3))
 
virtual ~G4GPSModel ()
 
void DescribeYourselfTo (G4VGraphicsScene &)
 
G4String GetCurrentDescription () const
 
G4String GetCurrentTag () const
 
- Public Member Functions inherited from G4VModel
 G4VModel (const G4ModelingParameters *=0)
 
virtual ~G4VModel ()
 
const G4ModelingParametersGetModelingParameters () const
 
const G4StringGetType () const
 
const G4VisExtentGetExtent () const
 
const G4StringGetGlobalDescription () const
 
const G4StringGetGlobalTag () const
 
void SetModelingParameters (const G4ModelingParameters *)
 
void SetExtent (const G4VisExtent &)
 
void SetType (const G4String &)
 
void SetGlobalDescription (const G4String &)
 
void SetGlobalTag (const G4String &)
 
virtual G4bool Validate (G4bool warn=true)
 

Protected Attributes

G4Colour fColour
 
- Protected Attributes inherited from G4VModel
G4String fType
 
G4String fGlobalTag
 
G4String fGlobalDescription
 
G4VisExtent fExtent
 
const G4ModelingParametersfpMP
 

Detailed Description

Definition at line 45 of file G4GPSModel.hh.

Constructor & Destructor Documentation

◆ G4GPSModel()

G4GPSModel::G4GPSModel ( const G4Colour & colour = G4Colour(1.,0.,0.,0.3))

Definition at line 51 of file G4GPSModel.cc.

52: fColour(colour)
53{
54 fType = "G4GPSModel";
55 std::ostringstream oss;
56 oss << "G4GPSModel for General Particle Source " << fColour;
57 fGlobalTag = oss.str();
59}
G4Colour fColour
Definition G4GPSModel.hh:67
G4String fGlobalDescription
Definition G4VModel.hh:96
G4String fType
Definition G4VModel.hh:94
G4String fGlobalTag
Definition G4VModel.hh:95

◆ ~G4GPSModel()

G4GPSModel::~G4GPSModel ( )
virtual

Definition at line 61 of file G4GPSModel.cc.

61{}

Member Function Documentation

◆ DescribeYourselfTo()

void G4GPSModel::DescribeYourselfTo ( G4VGraphicsScene & sceneHandler)
virtual

Implements G4VModel.

Definition at line 73 of file G4GPSModel.cc.

77{
79 // Note: As far as I can see, if this is the first time Instance has been
80 // called, it will, nevertheless, instantiate a default source, Type:"Point",
81 // Shaep: "NULL", which will be drawn as a small circle at the origin of
82 // coordinates whether you have set up GPS or not. Sorry, can't think of a
83 // way to avoid that. Mostly, of course, you will only invoke this function,
84 // if you have - or are about to - set up GPS, in which case all will be well.
85 if (!pGPSData) return;
86
87 G4int nSources = pGPSData->GetSourceVectorSize();
88 for (G4int iSource = 0; iSource < nSources; ++iSource) {
89
90 const G4SingleParticleSource* pCurrentSingleSource = pGPSData->GetCurrentSource(iSource);
91 if (!pCurrentSingleSource) return;
92
93 const G4SPSPosDistribution* pSPSPosDistribution = pCurrentSingleSource->GetPosDist();
94 if (!pSPSPosDistribution) return;
95
96 G4String Type = pSPSPosDistribution->GetPosDisType();
97 G4String Shape = pSPSPosDistribution->GetPosDisShape();
98 // Type can be: Point, Plane, Surface or Volume
99 // Shape can be: Square, Circle, Ellipse, Rectangle,
100 // Sphere, Ellipsoid, Cylinder, Parallelepiped
101// G4cout
102// << "G4GPSModel::DescribeYourselfTo"
103// << ": PosDisType: " << Type
104// << ", Shape: " << Shape
105// << G4endl;
106
107 const G4double& halfx = pSPSPosDistribution->GetHalfX();
108 const G4double& halfy = pSPSPosDistribution->GetHalfY();
109 const G4double& halfz = pSPSPosDistribution->GetHalfZ();
110 const G4double& Radius = pSPSPosDistribution->GetRadius();
111 const G4double& Radius0 = pSPSPosDistribution->GetRadius0();
112 const G4double& ParAlpha = pSPSPosDistribution->GetParAlpha();
113 const G4double& ParTheta = pSPSPosDistribution->GetParTheta();
114 const G4double& ParPhi = pSPSPosDistribution->GetParPhi();
115
116 const G4ThreeVector& Rotx = pSPSPosDistribution->GetRotx();
117 const G4ThreeVector& Roty = pSPSPosDistribution->GetRoty();
118 const G4ThreeVector& Rotz = pSPSPosDistribution->GetRotz();
119
120 const G4ThreeVector& position = pSPSPosDistribution->GetCentreCoords();
122 transform = G4Translate3D(position) * transform;
123
125 G4double smallHalfThickness = 10.*surfaceTolerance;
126
127 G4VisAttributes gpsAtts;
128 gpsAtts.SetColour(fColour);
129 gpsAtts.SetForceSolid();
130
131 if (Type == "Point") {
132
133 G4Circle circle;
134 circle.SetPosition(position);
135 circle.SetScreenDiameter(10.);
136 circle.SetVisAttributes(gpsAtts);
137 sceneHandler.BeginPrimitives(transform);
138 sceneHandler.AddPrimitive(circle);
139 sceneHandler.EndPrimitives();
140
141 } else if (Type == "Plane") {
142
143 // Code based on G4SPSPosDistribution::GeneratePointsInPlane.
144 sceneHandler.PreAddSolid(transform,gpsAtts);
145 if (Shape == "Circle") {
146 sceneHandler.AddSolid
147 (G4Tubs("GPS_Circle",0.,Radius,smallHalfThickness,0.,twopi));
148 } else if (Shape == "Annulus") {
149 sceneHandler.AddSolid
150 (G4Tubs("GPS_Annulus",Radius0,Radius,smallHalfThickness,0.,twopi));
151 } else if (Shape == "Ellipse") {
152 sceneHandler.AddSolid
153 (G4EllipticalTube("GPS_Ellipse",halfx,halfy,smallHalfThickness));
154 } else if (Shape == "Square") {
155 sceneHandler.AddSolid
156 (G4Box("GPS_Ellipse",halfx,halfx,smallHalfThickness));
157 } else if (Shape == "Rectangle") {
158 sceneHandler.AddSolid
159 (G4Box("GPS_Rectangle",halfx,halfy,smallHalfThickness));
160 }
161 sceneHandler.PostAddSolid();
162
163 } else if (Type == "Surface" || Type == "Volume") {
164
165 // Code based on G4SPSPosDistribution::GeneratePointsOnSurface.
166 // and G4SPSPosDistribution::GeneratePointsInVolume.
167 sceneHandler.PreAddSolid(transform,gpsAtts);
168 if (Shape == "Sphere") {
169 sceneHandler.AddSolid
170 (G4Orb("GPS_Sphere",Radius));
171 } else if (Shape == "Ellipsoid") {
172 sceneHandler.AddSolid
173 (G4Ellipsoid("GPS_Ellipsoid",halfx,halfy,halfz));
174 } else if (Shape == "Cylinder") {
175 sceneHandler.AddSolid
176 (G4Tubs("GPS_Cylinder",0.,Radius, halfz, 0., twopi));
177 } else if (Shape == "Para") {
178 sceneHandler.AddSolid
179 (G4Para("GPS_Para",halfx,halfy,halfz,ParAlpha,ParTheta,ParPhi));
180 }
181 sceneHandler.PostAddSolid();
182
183 }
184 }
185}
HepGeom::Translate3D G4Translate3D
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
Definition G4Box.hh:56
G4SingleParticleSource * GetCurrentSource(G4int idx)
static G4GeneralParticleSourceData * Instance()
G4double GetSurfaceTolerance() const
static G4GeometryTolerance * GetInstance()
Definition G4Orb.hh:56
const G4ThreeVector & GetCentreCoords() const
const G4String & GetPosDisType() const
const G4ThreeVector & GetRotx() const
const G4ThreeVector & GetRotz() const
const G4String & GetPosDisShape() const
const G4ThreeVector & GetRoty() const
G4SPSPosDistribution * GetPosDist() const
virtual void BeginPrimitives(const G4Transform3D &objectTransformation=G4Transform3D())=0
virtual void PostAddSolid()=0
virtual void AddPrimitive(const G4Polyline &)=0
virtual void AddSolid(const G4Box &)=0
virtual void EndPrimitives()=0
virtual void PreAddSolid(const G4Transform3D &objectTransformation, const G4VisAttributes &visAttribs)=0
void SetPosition(const G4Point3D &)
void SetScreenDiameter(G4double)
void SetColour(const G4Colour &)
void SetForceSolid(G4bool=true)
void SetVisAttributes(const G4VisAttributes *)
Definition G4Visible.cc:98
DLL_API const Hep3Vector HepZHat
DLL_API const Hep3Vector HepXHat
DLL_API const Hep3Vector HepYHat

◆ GetCurrentDescription()

G4String G4GPSModel::GetCurrentDescription ( ) const
virtual

Reimplemented from G4VModel.

Definition at line 68 of file G4GPSModel.cc.

69{
70 return "G4GPSModel " + GetCurrentTag ();
71}
G4String GetCurrentTag() const
Definition G4GPSModel.cc:63

◆ GetCurrentTag()

G4String G4GPSModel::GetCurrentTag ( ) const
virtual

Reimplemented from G4VModel.

Definition at line 63 of file G4GPSModel.cc.

64{
65return "";
66}

Referenced by GetCurrentDescription().

Member Data Documentation

◆ fColour

G4Colour G4GPSModel::fColour
protected

Definition at line 67 of file G4GPSModel.hh.

Referenced by DescribeYourselfTo(), and G4GPSModel().


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