Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4NormalNavigation.cc
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26// class G4NormalNavigation Implementation
27//
28// Author: P.Kent, 1996
29//
30// --------------------------------------------------------------------
31
32#include "G4NormalNavigation.hh"
33#include "G4NavigationLogger.hh"
34#include "G4AffineTransform.hh"
35
36// ********************************************************************
37// Constructor
38// ********************************************************************
39//
41{
42 fLogger = new G4NavigationLogger("G4NormalNavigation");
43}
44
45// ********************************************************************
46// Destructor
47// ********************************************************************
48//
50{
51 delete fLogger;
52}
53
54// ********************************************************************
55// ComputeStep
56// ********************************************************************
57//
58// On entry
59// exitNormal, validExitNormal: for previous exited volume (daughter)
60//
61// On exit
62// exitNormal, validExitNormal: for mother, if exiting it (else unchanged)
65 const G4ThreeVector& localDirection,
66 const G4double currentProposedStepLength,
67 G4double& newSafety,
68 G4NavigationHistory& history,
69 G4bool& validExitNormal,
70 G4ThreeVector& exitNormal,
71 G4bool& exiting,
72 G4bool& entering,
73 G4VPhysicalVolume* (*pBlockedPhysical),
74 G4int& blockedReplicaNo)
75{
76 G4VPhysicalVolume *motherPhysical, *samplePhysical,
77 *blockedExitedVol = nullptr;
78 G4LogicalVolume *motherLogical;
79 G4VSolid *motherSolid;
80 G4ThreeVector sampleDirection;
81 G4double ourStep = currentProposedStepLength, ourSafety;
82 G4double motherSafety, motherStep = DBL_MAX;
83 G4long localNoDaughters, sampleNo;
84 G4bool motherValidExitNormal = false;
85 G4ThreeVector motherExitNormal;
86
87 motherPhysical = history.GetTopVolume();
88 motherLogical = motherPhysical->GetLogicalVolume();
89 motherSolid = motherLogical->GetSolid();
90
91 // Compute mother safety
92 //
93 motherSafety = motherSolid->DistanceToOut(localPoint);
94 ourSafety = motherSafety; // Working isotropic safety
95
96 localNoDaughters = motherLogical->GetNoDaughters();
97
98#ifdef G4VERBOSE
99 if ( fCheck && ( (localNoDaughters>0) || (ourStep < motherSafety) ) )
100 {
101 fLogger->PreComputeStepLog(motherPhysical, motherSafety, localPoint);
102 }
103#endif
104 // Compute daughter safeties & intersections
105 //
106
107 // Exiting normal optimisation
108 //
109 if ( exiting && validExitNormal )
110 {
111 if ( localDirection.dot(exitNormal)>=kMinExitingNormalCosine )
112 {
113 // Block exited daughter volume
114 //
115 blockedExitedVol = (*pBlockedPhysical);
116 ourSafety = 0;
117 }
118 }
119 exiting = false;
120 entering = false;
121
122#ifdef G4VERBOSE
123 if ( fCheck )
124 {
125 // Compute early:
126 // a) to check whether point is (wrongly) outside
127 // (signaled if step < 0 or step == kInfinity )
128 // b) to check value against answer of daughters!
129
130 motherStep = motherSolid->DistanceToOut(localPoint,
131 localDirection,
132 true,
133 &motherValidExitNormal,
134 &motherExitNormal);
135
136 if( (motherStep >= kInfinity) || (motherStep < 0.0) )
137 {
138 // Error - indication of being outside solid !!
139 fLogger->ReportOutsideMother(localPoint, localDirection, motherPhysical);
140
141 ourStep = motherStep = 0.0;
142
143 exiting = true;
144 entering = false;
145
146 // If we are outside the solid does the normal make sense?
147 validExitNormal = motherValidExitNormal;
148 exitNormal = motherExitNormal;
149
150 *pBlockedPhysical = nullptr; // or motherPhysical ?
151 blockedReplicaNo = 0; // or motherReplicaNumber ?
152
153 newSafety = 0.0;
154 return ourStep;
155 }
156 }
157#endif
158
159 for ( sampleNo=localNoDaughters-1; sampleNo>=0; sampleNo--)
160 {
161 samplePhysical = motherLogical->GetDaughter(sampleNo);
162 if ( samplePhysical!=blockedExitedVol )
163 {
164 G4AffineTransform sampleTf(samplePhysical->GetRotation(),
165 samplePhysical->GetTranslation());
166 sampleTf.Invert();
167 const G4ThreeVector samplePoint = sampleTf.TransformPoint(localPoint);
168 const G4VSolid *sampleSolid =
169 samplePhysical->GetLogicalVolume()->GetSolid();
170 const G4double sampleSafety =
171 sampleSolid->DistanceToIn(samplePoint);
172
173 if ( sampleSafety<ourSafety )
174 {
175 ourSafety=sampleSafety;
176 }
177
178 if ( sampleSafety<=ourStep )
179 {
180 sampleDirection = sampleTf.TransformAxis(localDirection);
181 const G4double sampleStep =
182 sampleSolid->DistanceToIn(samplePoint,sampleDirection);
183#ifdef G4VERBOSE
184 if( fCheck )
185 {
186 fLogger->PrintDaughterLog(sampleSolid, samplePoint,
187 sampleSafety, true,
188 sampleDirection, sampleStep);
189 }
190#endif
191 if ( sampleStep<=ourStep )
192 {
193 ourStep = sampleStep;
194 entering = true;
195 exiting = false;
196 *pBlockedPhysical = samplePhysical;
197 blockedReplicaNo = -1;
198#ifdef G4VERBOSE
199 if( fCheck )
200 {
201 fLogger->AlongComputeStepLog(sampleSolid, samplePoint,
202 sampleDirection, localDirection,
203 sampleSafety, sampleStep);
204 }
205#endif
206 }
207
208#ifdef G4VERBOSE
209 if( fCheck && (sampleStep < kInfinity) && (sampleStep >= motherStep) )
210 {
211 // The intersection point with the daughter is at or after the exit
212 // point from the mother volume. Double check!
213 fLogger->CheckDaughterEntryPoint(sampleSolid,
214 samplePoint, sampleDirection,
215 motherSolid,
216 localPoint, localDirection,
217 motherStep, sampleStep);
218 }
219#endif
220 } // end of if ( sampleSafety <= ourStep )
221#ifdef G4VERBOSE
222 else if ( fCheck )
223 {
224 fLogger->PrintDaughterLog(sampleSolid, samplePoint,
225 sampleSafety, false,
226 G4ThreeVector(0.,0.,0.), -1.0 );
227 }
228#endif
229 }
230 }
231 if ( currentProposedStepLength<ourSafety )
232 {
233 // Guaranteed physics limited
234 //
235 entering = false;
236 exiting = false;
237 *pBlockedPhysical = nullptr;
238 ourStep = kInfinity;
239 }
240 else
241 {
242 // Consider intersection with mother solid
243 //
244 if ( motherSafety<=ourStep )
245 {
246 if ( !fCheck ) // The call is moved above when running in check_mode
247 {
248 motherStep = motherSolid->DistanceToOut(localPoint,
249 localDirection,
250 true,
251 &motherValidExitNormal,
252 &motherExitNormal);
253 }
254#ifdef G4VERBOSE
255 else // check_mode
256 {
257 fLogger->PostComputeStepLog(motherSolid, localPoint, localDirection,
258 motherStep, motherSafety);
259 if( motherValidExitNormal )
260 {
261 fLogger->CheckAndReportBadNormal(motherExitNormal,
262 localPoint,
263 localDirection,
264 motherStep,
265 motherSolid,
266 "From motherSolid::DistanceToOut" );
267 }
268 }
269#endif
270
271 if( (motherStep >= kInfinity) || (motherStep < 0.0) )
272 {
273#ifdef G4VERBOSE
274 if( fCheck ) // Clearly outside the mother solid!
275 {
276 fLogger->ReportOutsideMother(localPoint, localDirection,
277 motherPhysical);
278 }
279#endif
280 ourStep = motherStep = 0.0;
281 exiting = true;
282 entering = false;
283 // validExitNormal= motherValidExitNormal;
284 // exitNormal= motherExitNormal;
285 // The normal could be useful - but only if near the mother
286 // But it could be unreliable!
287 validExitNormal = false;
288 *pBlockedPhysical = nullptr; // or motherPhysical ?
289 blockedReplicaNo = 0; // or motherReplicaNumber ?
290 newSafety= 0.0;
291 return ourStep;
292 }
293
294 if ( motherStep<=ourStep )
295 {
296 ourStep = motherStep;
297 exiting = true;
298 entering = false;
299 validExitNormal = motherValidExitNormal;
300 exitNormal = motherExitNormal;
301
302 if ( motherValidExitNormal )
303 {
304 const G4RotationMatrix *rot = motherPhysical->GetRotation();
305 if (rot != nullptr)
306 {
307 exitNormal *= rot->inverse();
308#ifdef G4VERBOSE
309 if( fCheck )
310 {
311 fLogger->CheckAndReportBadNormal(exitNormal, // rotated
312 motherExitNormal, // original
313 *rot,
314 "From RotationMatrix" );
315 }
316#endif
317 }
318 }
319 }
320 else
321 {
322 validExitNormal = false;
323 }
324 }
325 }
326 newSafety = ourSafety;
327 return ourStep;
328}
329
330// ********************************************************************
331// ComputeSafety
332// ********************************************************************
333//
335 const G4NavigationHistory& history,
336 const G4double)
337{
338 G4VPhysicalVolume *motherPhysical, *samplePhysical;
339 G4LogicalVolume *motherLogical;
340 G4VSolid *motherSolid;
341 G4double motherSafety, ourSafety;
342 G4long localNoDaughters, sampleNo;
343
344 motherPhysical = history.GetTopVolume();
345 motherLogical = motherPhysical->GetLogicalVolume();
346 motherSolid = motherLogical->GetSolid();
347
348 // Compute mother safety
349 //
350 motherSafety = motherSolid->DistanceToOut(localPoint);
351 ourSafety = motherSafety; // Working isotropic safety
352
353#ifdef G4VERBOSE
354 if( fCheck )
355 {
356 fLogger->ComputeSafetyLog(motherSolid,localPoint,motherSafety,true,1);
357 }
358#endif
359
360 // Compute daughter safeties
361 //
362 localNoDaughters = motherLogical->GetNoDaughters();
363 for ( sampleNo=localNoDaughters-1; sampleNo>=0; sampleNo-- )
364 {
365 samplePhysical = motherLogical->GetDaughter(sampleNo);
366 G4AffineTransform sampleTf(samplePhysical->GetRotation(),
367 samplePhysical->GetTranslation());
368 sampleTf.Invert();
369 const G4ThreeVector samplePoint =
370 sampleTf.TransformPoint(localPoint);
371 const G4VSolid *sampleSolid =
372 samplePhysical->GetLogicalVolume()->GetSolid();
373 const G4double sampleSafety =
374 sampleSolid->DistanceToIn(samplePoint);
375 if ( sampleSafety<ourSafety )
376 {
377 ourSafety = sampleSafety;
378 }
379#ifdef G4VERBOSE
380 if(fCheck)
381 {
382 fLogger->ComputeSafetyLog(sampleSolid, samplePoint,
383 sampleSafety, false, 0);
384 // Not mother, no banner
385 }
386#endif
387 }
388 return ourSafety;
389}
390
391// The following methods have been imported to this source file
392// in order to avoid dependency of the header file on the
393// header implementation of G4NavigationLogger.
394
395// ********************************************************************
396// GetVerboseLevel
397// ********************************************************************
398//
400{
401 return fLogger->GetVerboseLevel();
402}
403
404// ********************************************************************
405// SetVerboseLevel
406// ********************************************************************
407//
409{
410 fLogger->SetVerboseLevel(level);
411}
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition G4Types.hh:83
long G4long
Definition G4Types.hh:87
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
double dot(const Hep3Vector &) const
HepRotation inverse() const
G4AffineTransform & Invert()
G4ThreeVector TransformPoint(const G4ThreeVector &vec) const
G4ThreeVector TransformAxis(const G4ThreeVector &axis) const
G4VSolid * GetSolid() const
std::size_t GetNoDaughters() const
G4VPhysicalVolume * GetDaughter(const std::size_t i) const
G4VPhysicalVolume * GetTopVolume() const
void SetVerboseLevel(G4int level)
void PreComputeStepLog(const G4VPhysicalVolume *motherPhysical, G4double motherSafety, const G4ThreeVector &localPoint) const
void CheckDaughterEntryPoint(const G4VSolid *sampleSolid, const G4ThreeVector &samplePoint, const G4ThreeVector &sampleDirection, const G4VSolid *motherSolid, const G4ThreeVector &localPoint, const G4ThreeVector &localDirection, G4double motherStep, G4double sampleStep) const
G4int GetVerboseLevel() const
void PrintDaughterLog(const G4VSolid *sampleSolid, const G4ThreeVector &samplePoint, G4double sampleSafety, G4bool onlySafety, const G4ThreeVector &sampleDirection, G4double sampleStep) const
G4bool CheckAndReportBadNormal(const G4ThreeVector &unitNormal, const G4ThreeVector &localPoint, const G4ThreeVector &localDirection, G4double step, const G4VSolid *solid, const char *msg) const
void ComputeSafetyLog(const G4VSolid *solid, const G4ThreeVector &point, G4double safety, G4bool isMotherVolume, G4int banner=-1) const
void PostComputeStepLog(const G4VSolid *motherSolid, const G4ThreeVector &localPoint, const G4ThreeVector &localDirection, G4double motherStep, G4double motherSafety) const
void ReportOutsideMother(const G4ThreeVector &localPoint, const G4ThreeVector &localDirection, const G4VPhysicalVolume *motherPV, G4double tDist=30.0 *CLHEP::cm) const
void AlongComputeStepLog(const G4VSolid *sampleSolid, const G4ThreeVector &samplePoint, const G4ThreeVector &sampleDirection, const G4ThreeVector &localDirection, G4double sampleSafety, G4double sampleStep) const
virtual void SetVerboseLevel(G4int level) final
G4double ComputeSafety(const G4ThreeVector &globalpoint, const G4NavigationHistory &history, const G4double pMaxLength=DBL_MAX) final
G4double ComputeStep(const G4ThreeVector &localPoint, const G4ThreeVector &localDirection, const G4double currentProposedStepLength, G4double &newSafety, G4NavigationHistory &history, G4bool &validExitNormal, G4ThreeVector &exitNormal, G4bool &exiting, G4bool &entering, G4VPhysicalVolume *(*pBlockedPhysical), G4int &blockedReplicaNo) final
virtual G4int GetVerboseLevel() const final
const G4RotationMatrix * GetRotation() const
const G4ThreeVector GetTranslation() const
G4LogicalVolume * GetLogicalVolume() const
virtual G4double DistanceToOut(const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcNorm=false, G4bool *validNorm=nullptr, G4ThreeVector *n=nullptr) const =0
virtual G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const =0
#define DBL_MAX
Definition templates.hh:62