Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4VFieldModel Class Referenceabstract

#include <G4VFieldModel.hh>

+ Inheritance diagram for G4VFieldModel:

Public Types

enum  Representation { fullArrow , lightArrow }
 

Public Member Functions

 G4VFieldModel (const G4String &typeOfField, const G4String &symbol="", const G4VisExtent &extentForField=G4VisExtent(), const std::vector< G4PhysicalVolumesSearchScene::Findings > &pvFindings=std::vector< G4PhysicalVolumesSearchScene::Findings >(), G4int nDataPointsPerHalfScene=10, Representation representation=Representation::fullArrow, G4int arrow3DLineSegmentsPerCircle=6)
 
virtual ~G4VFieldModel ()
 
virtual void DescribeYourselfTo (G4VGraphicsScene &sceneHandler)
 
- Public Member Functions inherited from G4VModel
 G4VModel (const G4ModelingParameters *=0)
 
virtual ~G4VModel ()
 
virtual void DescribeYourselfTo (G4VGraphicsScene &)=0
 
const G4ModelingParametersGetModelingParameters () const
 
const G4StringGetType () const
 
virtual G4String GetCurrentDescription () const
 
virtual G4String GetCurrentTag () 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 Member Functions

virtual void GetFieldAtLocation (const G4Field *field, const G4Point3D &position, G4double time, G4Point3D &result) const =0
 

Additional Inherited Members

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

Detailed Description

Definition at line 46 of file G4VFieldModel.hh.

Member Enumeration Documentation

◆ Representation

Enumerator
fullArrow 
lightArrow 

Definition at line 50 of file G4VFieldModel.hh.

Constructor & Destructor Documentation

◆ G4VFieldModel()

G4VFieldModel::G4VFieldModel ( const G4String typeOfField,
const G4String symbol = "",
const G4VisExtent extentForField = G4VisExtent(),
const std::vector< G4PhysicalVolumesSearchScene::Findings > &  pvFindings = std::vector<G4PhysicalVolumesSearchScene::Findings>(),
G4int  nDataPointsPerHalfScene = 10,
Representation  representation = Representation::fullArrow,
G4int  arrow3DLineSegmentsPerCircle = 6 
)

Definition at line 61 of file G4VFieldModel.cc.

68: fExtentForField(extentForField)
69, fPVFindings(pvFindings)
70, fNDataPointsPerMaxHalfExtent(nDataPointsPerMaxHalfExtent)
71, fRepresentation(representation)
72, fArrow3DLineSegmentsPerCircle(arrow3DLineSegmentsPerCircle)
73, fTypeOfField(typeOfField)
74, fArrowPrefix(symbol)
75{
76 fType = "G4"+typeOfField+"FieldModel";
78
79 std::ostringstream oss;
80 oss << ':' << fNDataPointsPerMaxHalfExtent
81 << ':' << fArrow3DLineSegmentsPerCircle;
82 if (fExtentForField == G4VisExtent::GetNullExtent()) {
83 oss << " whole scene";
84 } else {
85 oss
86 << ':' << fExtentForField.GetXmin()
87 << ':' << fExtentForField.GetXmax()
88 << ':' << fExtentForField.GetYmin()
89 << ':' << fExtentForField.GetYmax()
90 << ':' << fExtentForField.GetZmin()
91 << ':' << fExtentForField.GetZmax();
92 }
93 for (const auto& findings: fPVFindings) {
94 oss
95 << ',' << findings.fpFoundPV->GetName()
96 << ':' << findings.fFoundPVCopyNo;
97 }
98 if (fRepresentation == Representation::fullArrow) {
99 oss << " full arrow";
100 } else if (fRepresentation == Representation::lightArrow) {
101 oss << " light arrow";
102 }
103
104 fGlobalDescription = fType + oss.str();
105}
G4String fGlobalDescription
Definition: G4VModel.hh:100
G4String fType
Definition: G4VModel.hh:98
G4String fGlobalTag
Definition: G4VModel.hh:99
static const G4VisExtent & GetNullExtent()
Definition: G4VisExtent.cc:60
G4double GetYmin() const
Definition: G4VisExtent.hh:101
G4double GetXmax() const
Definition: G4VisExtent.hh:100
G4double GetYmax() const
Definition: G4VisExtent.hh:102
G4double GetZmax() const
Definition: G4VisExtent.hh:104
G4double GetZmin() const
Definition: G4VisExtent.hh:103
G4double GetXmin() const
Definition: G4VisExtent.hh:99

◆ ~G4VFieldModel()

G4VFieldModel::~G4VFieldModel ( )
virtual

Definition at line 59 of file G4VFieldModel.cc.

59{;}

Member Function Documentation

◆ DescribeYourselfTo()

void G4VFieldModel::DescribeYourselfTo ( G4VGraphicsScene sceneHandler)
virtual

Implements G4VModel.

Definition at line 110 of file G4VFieldModel.cc.

110 {
111// G4cout << "G4VFieldModel::DescribeYourselfTo" << G4endl;
112
115 assert(tMgr);
116 G4Navigator* navigator = tMgr->GetNavigatorForTracking();
117 assert(navigator);
118
119 G4FieldManager* globalFieldMgr = tMgr->GetFieldManager();
120 const G4Field* globalField = 0;
121 const G4String intro = "G4VFieldModel::DescribeYourselfTo: ";
122 if (globalFieldMgr) {
123 if (globalFieldMgr->DoesFieldExist()) {
124 globalField = globalFieldMgr->GetDetectorField();
125 if (!globalField) {
126 static G4bool warned = false;
127 if (!warned) {
128 G4warn << intro << "Null global field pointer." << G4endl;
129 warned = true;
130 }
131 }
132 }
133 } else {
134 static G4bool warned = false;
135 if (!warned) {
136 G4warn << intro << "No global field manager." << G4endl;
137 warned = true;
138 }
139 }
140
141 G4VisExtent extent = sceneHandler.GetExtent();
142 if (fExtentForField == G4VisExtent::GetNullExtent()) {
143 extent = sceneHandler.GetExtent();
144 } else {
145 extent = fExtentForField;
146 }
147 const G4double& xMin = extent.GetXmin();
148 const G4double& yMin = extent.GetYmin();
149 const G4double& zMin = extent.GetZmin();
150 const G4double& xMax = extent.GetXmax();
151 const G4double& yMax = extent.GetYmax();
152 const G4double& zMax = extent.GetZmax();
153 const G4double xHalfScene = 0.5 * (xMax - xMin);
154 const G4double yHalfScene = 0.5 * (yMax - yMin);
155 const G4double zHalfScene = 0.5 * (zMax - zMin);
156 const G4double xSceneCentre = 0.5 * (xMax + xMin);
157 const G4double ySceneCentre = 0.5 * (yMax + yMin);
158 const G4double zSceneCentre = 0.5 * (zMax + zMin);
159 const G4double maxHalfScene =
160 std::max(xHalfScene,std::max(yHalfScene,zHalfScene));
161 if (maxHalfScene <= 0.) {
162 G4warn << "Scene extent non-positive." << G4endl;
163 return;
164 }
165
166 // Constants
167 const G4double interval = maxHalfScene / fNDataPointsPerMaxHalfExtent;
168 const G4int nDataPointsPerXHalfScene = G4int(xHalfScene / interval);
169 const G4int nDataPointsPerYHalfScene = G4int(yHalfScene / interval);
170 const G4int nDataPointsPerZHalfScene = G4int(zHalfScene / interval);
171 const G4int nXSamples = 2 * nDataPointsPerXHalfScene + 1;
172 const G4int nYSamples = 2 * nDataPointsPerYHalfScene + 1;
173 const G4int nZSamples = 2 * nDataPointsPerZHalfScene + 1;
174 const G4int nSamples = nXSamples * nYSamples * nZSamples;
175 const G4double arrowLengthMax = 0.8 * interval;
176
177 // Working vectors for field values, etc.
178 std::vector<G4Point3D> Field(nSamples); // Initialises to (0,0,0)
179 std::vector<G4Point3D> xyz(nSamples); // Initialises to (0,0,0)
180 G4double FieldMagnitudeMax = -std::numeric_limits<G4double>::max();
181
182 // Get field values and ascertain maximum field.
183 for (G4int i = 0; i < nXSamples; i++) {
184 G4double x = xSceneCentre + (i - nDataPointsPerXHalfScene) * interval;
185
186 for (G4int j = 0; j < nYSamples; j++) {
187 G4double y = ySceneCentre + (j - nDataPointsPerYHalfScene) * interval;
188
189 for (G4int k = 0; k < nZSamples; k++) {
190 G4double z = zSceneCentre + (k - nDataPointsPerZHalfScene) * interval;
191
192 // Calculate indices into working vectors
193 const G4int ijk = i * nYSamples * nZSamples + j * nZSamples + k;
194 xyz[ijk].set(x,y,z);
195
196 G4ThreeVector pos(x,y,z);
197
198 // Check if point is in findings
199 if (!fPVFindings.empty()) {
200 G4bool isInPV = false;
201 for (const auto& findings: fPVFindings) {
202 G4VPhysicalVolume* pv = findings.fpFoundPV;
203 G4int copyNo = findings.fFoundPVCopyNo;
204 G4VSolid* solid = pv->GetLogicalVolume()->GetSolid();
205 G4PVParameterised* pvParam = dynamic_cast<G4PVParameterised*>(pv);
206 if (pvParam) {
207 auto* param = pvParam->GetParameterisation();
208 solid = param->ComputeSolid(copyNo,pvParam);
209 solid->ComputeDimensions(param,copyNo,pvParam);
210 }
211 // Transform point to local coordinate system
212 const auto& transform = findings.fFoundObjectTransformation;
213 auto rotation = transform.getRotation();
214 auto translation = transform.getTranslation();
215 G4ThreeVector lPos = pos; lPos -= translation; lPos.transform(rotation.invert());
216 if (solid->Inside(lPos)==kInside) {
217 isInPV = true;
218 break;
219 }
220 }
221 if (!isInPV) continue;
222 }
223 // Point is in findings - or there were no findings
224
225 // Find volume and field at this location.
226 const G4VPhysicalVolume* pPV =
227 navigator->LocateGlobalPointAndSetup(pos,0,false,true);
228 const G4Field* field = globalField;
229 if (pPV) {
230 // Get logical volume.
231 const G4LogicalVolume* pLV = pPV->GetLogicalVolume();
232 if (pLV) {
233 // Value for Region, if any, overrides
234 G4Region* pRegion = pLV->GetRegion();
235 if (pRegion) {
236 G4FieldManager* pRegionFieldMgr = pRegion->GetFieldManager();
237 if (pRegionFieldMgr) {
238 field = pRegionFieldMgr->GetDetectorField();
239 // G4cout << "Region with field" << G4endl;
240 }
241 }
242 // 'Local' value from logical volume, if any, overrides
243 G4FieldManager* pLVFieldMgr = pLV->GetFieldManager();
244 if (pLVFieldMgr) {
245 field = pLVFieldMgr->GetDetectorField();
246 // G4cout << "Logical volume with field" << G4endl;
247 }
248 }
249 }
250
251 G4double time = 0.; // FIXME: Can we get event time in some way?
252
253 // Subclasses will have implemented this for their own field
254 GetFieldAtLocation(field, xyz[ijk], time, Field[ijk]);
255
256 G4double mag = Field[ijk].mag();
257 if (mag > FieldMagnitudeMax) FieldMagnitudeMax = mag;
258 } // for (k, z
259 } // for (j, y
260 } // for (i, x
261
262 if (FieldMagnitudeMax <= 0.) {
263 G4warn << "No " << fTypeOfField << " field in this extent." << G4endl;
264 return;
265 }
266
267 for (G4int i = 0; i < nSamples; i++) {
268 const G4double Fmag = Field[i].mag();
269 const G4double f = Fmag / FieldMagnitudeMax;
270 if (f <= 0.) continue; // Skip zero field locations
271
272 G4double red = 0., green = 0., blue = 0., alpha = 1.;
273 if (f < 0.5) { // Linear colour scale: 0->0.5->1 is red->green->blue.
274 green = 2. * f;
275 red = 2. * (0.5 - f);
276 } else {
277 blue = 2. * (f - 0.5);
278 green = 2. * (1.0 - f);
279 }
280 const G4Colour arrowColour(red,green,blue,alpha);
281
282 // Very small arrows are difficult to see. Better to draw a line.
283 G4bool drawAsLine = false;
284 switch (fRepresentation) {
286 if (f < 0.1) {
287 drawAsLine = true;
288 }
289 break;
291 drawAsLine = true;
292 break;
293 default:
294 break;
295 }
296
297 // Head of arrow depends on field direction and strength...
298 G4double arrowLength = arrowLengthMax * f;
299 // ...but limit the length so it's visible.
300 if (f < 0.01) arrowLength = arrowLengthMax * 0.01;
301 const G4Point3D head = xyz[i] + arrowLength*Field[i]/Fmag;
302
303 if (drawAsLine) {
304 G4Polyline FArrowLite;
305 G4VisAttributes va(arrowColour);
306 va.SetLineWidth(2.);
307 FArrowLite.SetVisAttributes(va);
308 FArrowLite.push_back(xyz[i]);
309 FArrowLite.push_back(head);
310 sceneHandler.BeginPrimitives();
311 sceneHandler.AddPrimitive(FArrowLite);
312 sceneHandler.EndPrimitives();
313 } else {
314 G4ArrowModel FArrow(xyz[i].x(), xyz[i].y(), xyz[i].z(),
315 head.x(), head.y(), head.z(),
316 arrowLength/5, arrowColour,
317 fArrowPrefix+"Field",
318 fArrow3DLineSegmentsPerCircle);
319 FArrow.DescribeYourselfTo(sceneHandler);
320 }
321 } // for (i, nSamples
322}
#define G4warn
Definition: G4Scene.cc:41
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
Hep3Vector & transform(const HepRotation &)
Definition: ThreeVectorR.cc:20
const G4Field * GetDetectorField() const
G4bool DoesFieldExist() const
G4VSolid * GetSolid() const
G4Region * GetRegion() const
G4FieldManager * GetFieldManager() const
virtual G4VPhysicalVolume * LocateGlobalPointAndSetup(const G4ThreeVector &point, const G4ThreeVector *direction=nullptr, const G4bool pRelativeSearch=true, const G4bool ignoreDirection=true)
Definition: G4Navigator.cc:132
G4VPVParameterisation * GetParameterisation() const
G4FieldManager * GetFieldManager() const
static G4TransportationManager * GetTransportationManager()
G4Navigator * GetNavigatorForTracking() const
G4FieldManager * GetFieldManager() const
virtual void GetFieldAtLocation(const G4Field *field, const G4Point3D &position, G4double time, G4Point3D &result) const =0
virtual void BeginPrimitives(const G4Transform3D &objectTransformation=G4Transform3D())=0
virtual void AddPrimitive(const G4Polyline &)=0
virtual const G4VisExtent & GetExtent() const
virtual void EndPrimitives()=0
G4LogicalVolume * GetLogicalVolume() const
virtual EInside Inside(const G4ThreeVector &p) const =0
virtual void ComputeDimensions(G4VPVParameterisation *p, const G4int n, const G4VPhysicalVolume *pRep)
Definition: G4VSolid.cc:137
void SetVisAttributes(const G4VisAttributes *)
Definition: G4Visible.cc:98
@ kInside
Definition: geomdefs.hh:70

◆ GetFieldAtLocation()

virtual void G4VFieldModel::GetFieldAtLocation ( const G4Field field,
const G4Point3D position,
G4double  time,
G4Point3D result 
) const
protectedpure virtual

Implemented in G4ElectricFieldModel, and G4MagneticFieldModel.

Referenced by DescribeYourselfTo().


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