Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4VoxelNavigation.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 G4VoxelNavigation Implementation
27//
28// Author: P.Kent, 1996
29//
30// --------------------------------------------------------------------
31#include <ostream>
32
33#include "G4VoxelNavigation.hh"
35#include "G4VoxelSafety.hh"
36
38
39// ********************************************************************
40// Constructor
41// ********************************************************************
42//
44 : fBList(),
45 fVoxelAxisStack(kNavigatorVoxelStackMax,kXAxis),
46 fVoxelNoSlicesStack(kNavigatorVoxelStackMax,0),
47 fVoxelSliceWidthStack(kNavigatorVoxelStackMax,0.),
48 fVoxelNodeNoStack(kNavigatorVoxelStackMax,0),
49 fVoxelHeaderStack(kNavigatorVoxelStackMax,(G4SmartVoxelHeader*)nullptr)
50{
51 fLogger= new G4NavigationLogger("G4VoxelNavigation");
54
55#ifdef G4DEBUG_NAVIGATION
56 SetVerboseLevel(5); // Reports most about daughter volumes
57#endif
58}
59
60// ********************************************************************
61// Destructor
62// ********************************************************************
63//
65{
66 delete fpVoxelSafety;
67 delete fLogger;
68}
69
70// --------------------------------------------------------------------------
71// Input:
72// exiting: : last step exited
73// blockedPhysical : phys volume last exited (if exiting)
74// blockedReplicaNo : copy/replica number of exited
75// Output:
76// entering : if true, found candidate volume to enter
77// blockedPhysical : candidate phys volume to enter - if entering
78// blockedReplicaNo : copy/replica number - if entering
79// exiting: : will exit current (mother) volume
80// In/Out
81// --------------------------------------------------------------------------
82
83// ********************************************************************
84// ComputeStep
85// ********************************************************************
86//
89 const G4ThreeVector& localDirection,
90 const G4double currentProposedStepLength,
91 G4double& newSafety,
92 /* const */ G4NavigationHistory& history,
93 G4bool& validExitNormal,
94 G4ThreeVector& exitNormal,
95 G4bool& exiting,
96 G4bool& entering,
97 G4VPhysicalVolume* (*pBlockedPhysical),
98 G4int& blockedReplicaNo )
99{
100 G4VPhysicalVolume *motherPhysical, *samplePhysical, *blockedExitedVol=nullptr;
101 G4LogicalVolume *motherLogical;
102 G4VSolid *motherSolid;
103 G4ThreeVector sampleDirection;
104 G4double ourStep=currentProposedStepLength, ourSafety;
105 G4double motherSafety, motherStep = DBL_MAX;
106 G4int localNoDaughters, sampleNo;
107
108 G4bool initialNode, noStep;
109 G4SmartVoxelNode *curVoxelNode;
110 G4long curNoVolumes, contentNo;
111 G4double voxelSafety;
112
113 motherPhysical = history.GetTopVolume();
114 motherLogical = motherPhysical->GetLogicalVolume();
115 motherSolid = motherLogical->GetSolid();
116
117 //
118 // Compute mother safety
119 //
120
121 motherSafety = motherSolid->DistanceToOut(localPoint);
122 ourSafety = motherSafety; // Working isotropic safety
123
124#ifdef G4VERBOSE
125 if ( fCheck )
126 {
127 fLogger->PreComputeStepLog (motherPhysical, motherSafety, localPoint);
128 }
129#endif
130
131 //
132 // Compute daughter safeties & intersections
133 //
134
135 // Exiting normal optimisation
136 //
137 if ( exiting && validExitNormal )
138 {
139 if ( localDirection.dot(exitNormal)>=kMinExitingNormalCosine )
140 {
141 // Block exited daughter volume
142 //
143 blockedExitedVol = *pBlockedPhysical;
144 ourSafety = 0;
145 }
146 }
147 exiting = false;
148 entering = false;
149
150 // For extra checking, get the distance to Mother early !!
151 G4bool motherValidExitNormal = false;
152 G4ThreeVector motherExitNormal(0.0, 0.0, 0.0);
153
154#ifdef G4VERBOSE
155 if ( fCheck )
156 {
157 // Compute early -- a) for validity
158 // b) to check against answer of daughters!
159 motherStep = motherSolid->DistanceToOut(localPoint,
160 localDirection,
161 true,
162 &motherValidExitNormal,
163 &motherExitNormal);
164 }
165#endif
166
167 localNoDaughters = (G4int)motherLogical->GetNoDaughters();
168
169 fBList.Enlarge(localNoDaughters);
170 fBList.Reset();
171
172 initialNode = true;
173 noStep = true;
174
175 while (noStep)
176 {
177 curVoxelNode = fVoxelNode;
178 curNoVolumes = curVoxelNode->GetNoContained();
179 for (contentNo=curNoVolumes-1; contentNo>=0; contentNo--)
180 {
181 sampleNo = curVoxelNode->GetVolume((G4int)contentNo);
182 if ( !fBList.IsBlocked(sampleNo) )
183 {
184 fBList.BlockVolume(sampleNo);
185 samplePhysical = motherLogical->GetDaughter(sampleNo);
186 if ( samplePhysical!=blockedExitedVol )
187 {
188 G4AffineTransform sampleTf(samplePhysical->GetRotation(),
189 samplePhysical->GetTranslation());
190 sampleTf.Invert();
191 const G4ThreeVector samplePoint =
192 sampleTf.TransformPoint(localPoint);
193 const G4VSolid *sampleSolid =
194 samplePhysical->GetLogicalVolume()->GetSolid();
195 const G4double sampleSafety =
196 sampleSolid->DistanceToIn(samplePoint);
197
198 if ( sampleSafety<ourSafety )
199 {
200 ourSafety = sampleSafety;
201 }
202 if ( sampleSafety<=ourStep )
203 {
204 sampleDirection = sampleTf.TransformAxis(localDirection);
205 G4double sampleStep =
206 sampleSolid->DistanceToIn(samplePoint, sampleDirection);
207#ifdef G4VERBOSE
208 if( fCheck )
209 {
210 fLogger->PrintDaughterLog(sampleSolid, samplePoint,
211 sampleSafety, true,
212 sampleDirection, sampleStep);
213 }
214#endif
215 if ( sampleStep<=ourStep )
216 {
217 ourStep = sampleStep;
218 entering = true;
219 exiting = false;
220 *pBlockedPhysical = samplePhysical;
221 blockedReplicaNo = -1;
222#ifdef G4VERBOSE
223 // Check to see that the resulting point is indeed in/on volume.
224 // This could be done only for successful candidate.
225 if ( fCheck )
226 {
227 fLogger->AlongComputeStepLog (sampleSolid, samplePoint,
228 sampleDirection, localDirection, sampleSafety, sampleStep);
229 }
230#endif
231 }
232#ifdef G4VERBOSE
233 if ( fCheck && ( sampleStep < kInfinity )
234 && ( sampleStep >= motherStep ) )
235 {
236 // The intersection point with the daughter is after the exit
237 // point from the mother volume. Double check this !!
238 fLogger->CheckDaughterEntryPoint(sampleSolid,
239 samplePoint, sampleDirection,
240 motherSolid,
241 localPoint, localDirection,
242 motherStep, sampleStep);
243 }
244#endif
245 }
246#ifdef G4VERBOSE
247 else // ie if sampleSafety > outStep
248 {
249 if( fCheck )
250 {
251 fLogger->PrintDaughterLog(sampleSolid, samplePoint,
252 sampleSafety, false,
253 G4ThreeVector(0.,0.,0.), -1.0 );
254 }
255 }
256#endif
257 }
258 }
259 }
260 if (initialNode)
261 {
262 initialNode = false;
263 voxelSafety = ComputeVoxelSafety(localPoint);
264 if ( voxelSafety<ourSafety )
265 {
266 ourSafety = voxelSafety;
267 }
268 if ( currentProposedStepLength<ourSafety )
269 {
270 // Guaranteed physics limited
271 //
272 noStep = false;
273 entering = false;
274 exiting = false;
275 *pBlockedPhysical = nullptr;
276 ourStep = kInfinity;
277 }
278 else
279 {
280 //
281 // Compute mother intersection if required
282 //
283 if ( motherSafety<=ourStep )
284 {
285 // In case of check mode this is a duplicate call -- acceptable
286 motherStep = motherSolid->DistanceToOut(localPoint, localDirection,
287 true, &motherValidExitNormal, &motherExitNormal);
288#ifdef G4VERBOSE
289 if ( fCheck )
290 {
291 fLogger->PostComputeStepLog(motherSolid, localPoint, localDirection,
292 motherStep, motherSafety);
293 if( motherValidExitNormal )
294 {
295 fLogger->CheckAndReportBadNormal(motherExitNormal,
296 localPoint, localDirection,
297 motherStep, motherSolid,
298 "From motherSolid::DistanceToOut" );
299 }
300 }
301#endif
302 if( (motherStep >= kInfinity) || (motherStep < 0.0) )
303 {
304#ifdef G4VERBOSE
305 if( fCheck ) // Error - indication of being outside solid !!
306 {
307 fLogger->ReportOutsideMother(localPoint, localDirection,
308 motherPhysical);
309 }
310#endif
311 motherStep = 0.0;
312 ourStep = 0.0;
313 exiting = true;
314 entering = false;
315
316 // validExitNormal= motherValidExitNormal;
317 // exitNormal= motherExitNormal;
318 // Useful only if the point is very close to surface
319 // => but it would need to be rotated to grand-mother ref frame !
320 validExitNormal= false;
321
322 *pBlockedPhysical = nullptr; // or motherPhysical ?
323 blockedReplicaNo = 0; // or motherReplicaNumber ?
324
325 newSafety = 0.0;
326 return ourStep;
327 }
328
329 if ( motherStep<=ourStep )
330 {
331 ourStep = motherStep;
332 exiting = true;
333 entering = false;
334
335 // Exit normal: Natural location to set these;confirmed short step
336 //
337 validExitNormal = motherValidExitNormal;
338 exitNormal = motherExitNormal;
339
340 if ( validExitNormal )
341 {
342 const G4RotationMatrix *rot = motherPhysical->GetRotation();
343 if (rot)
344 {
345 exitNormal *= rot->inverse();
346#ifdef G4VERBOSE
347 if( fCheck )
348 fLogger->CheckAndReportBadNormal(exitNormal, // rotated
349 motherExitNormal, // original
350 *rot,
351 "From RotationMatrix" );
352#endif
353 }
354 }
355 }
356 else
357 {
358 validExitNormal = false;
359 }
360 }
361 }
362 newSafety = ourSafety;
363 }
364 if (noStep)
365 {
366 noStep = LocateNextVoxel(localPoint, localDirection, ourStep);
367 }
368 } // end -while (noStep)- loop
369
370 return ourStep;
371}
372
373// ********************************************************************
374// ComputeVoxelSafety
375//
376// Computes safety from specified point to voxel boundaries
377// using already located point
378// o collected boundaries for most derived level
379// o adjacent boundaries for previous levels
380// ********************************************************************
381//
384{
385 G4SmartVoxelHeader *curHeader;
386 G4double voxelSafety, curNodeWidth;
387 G4double curNodeOffset, minCurCommonDelta, maxCurCommonDelta;
388 G4int minCurNodeNoDelta, maxCurNodeNoDelta;
389 G4int localVoxelDepth, curNodeNo;
390 EAxis curHeaderAxis;
391
392 localVoxelDepth = fVoxelDepth;
393
394 curHeader = fVoxelHeaderStack[localVoxelDepth];
395 curHeaderAxis = fVoxelAxisStack[localVoxelDepth];
396 curNodeNo = fVoxelNodeNoStack[localVoxelDepth];
397 curNodeWidth = fVoxelSliceWidthStack[localVoxelDepth];
398
399 // Compute linear intersection distance to boundaries of max/min
400 // to collected nodes at current level
401 //
402 curNodeOffset = curNodeNo*curNodeWidth;
403 maxCurNodeNoDelta = fVoxelNode->GetMaxEquivalentSliceNo()-curNodeNo;
404 minCurNodeNoDelta = curNodeNo-fVoxelNode->GetMinEquivalentSliceNo();
405 minCurCommonDelta = localPoint(curHeaderAxis)
406 - curHeader->GetMinExtent() - curNodeOffset;
407 maxCurCommonDelta = curNodeWidth-minCurCommonDelta;
408
409 if ( minCurNodeNoDelta<maxCurNodeNoDelta )
410 {
411 voxelSafety = minCurNodeNoDelta*curNodeWidth;
412 voxelSafety += minCurCommonDelta;
413 }
414 else if (maxCurNodeNoDelta < minCurNodeNoDelta)
415 {
416 voxelSafety = maxCurNodeNoDelta*curNodeWidth;
417 voxelSafety += maxCurCommonDelta;
418 }
419 else // (maxCurNodeNoDelta == minCurNodeNoDelta)
420 {
421 voxelSafety = minCurNodeNoDelta*curNodeWidth;
422 voxelSafety += std::min(minCurCommonDelta,maxCurCommonDelta);
423 }
424
425 // Compute isotropic safety to boundaries of previous levels
426 // [NOT to collected boundaries]
427
428 // Loop checking, 07.10.2016, JA
429 while ( (localVoxelDepth>0) && (voxelSafety>0) )
430 {
431 localVoxelDepth--;
432 curHeader = fVoxelHeaderStack[localVoxelDepth];
433 curHeaderAxis = fVoxelAxisStack[localVoxelDepth];
434 curNodeNo = fVoxelNodeNoStack[localVoxelDepth];
435 curNodeWidth = fVoxelSliceWidthStack[localVoxelDepth];
436 curNodeOffset = curNodeNo*curNodeWidth;
437 minCurCommonDelta = localPoint(curHeaderAxis)
438 - curHeader->GetMinExtent() - curNodeOffset;
439 maxCurCommonDelta = curNodeWidth-minCurCommonDelta;
440
441 if ( minCurCommonDelta<voxelSafety )
442 {
443 voxelSafety = minCurCommonDelta;
444 }
445 if ( maxCurCommonDelta<voxelSafety )
446 {
447 voxelSafety = maxCurCommonDelta;
448 }
449 }
450 if ( voxelSafety<0 )
451 {
452 voxelSafety = 0;
453 }
454
455 return voxelSafety;
456}
457
458// ********************************************************************
459// LocateNextVoxel
460//
461// Finds the next voxel from the current voxel and point
462// in the specified direction
463//
464// Returns false if all voxels considered
465// [current Step ends inside same voxel or leaves all voxels]
466// true otherwise
467// [the information on the next voxel is put into the set of
468// fVoxel* variables & "stacks"]
469// ********************************************************************
470//
471G4bool
473 const G4ThreeVector& localDirection,
474 const G4double currentStep)
475{
476 G4SmartVoxelHeader *workHeader=nullptr, *newHeader=nullptr;
477 G4SmartVoxelProxy *newProxy=nullptr;
478 G4SmartVoxelNode *newVoxelNode=nullptr;
479 G4ThreeVector targetPoint, voxelPoint;
480 G4double workNodeWidth, workMinExtent, workCoord;
481 G4double minVal, maxVal, newDistance=0.;
482 G4double newHeaderMin, newHeaderNodeWidth;
483 G4int depth=0, newDepth=0, workNodeNo=0, newNodeNo=0, newHeaderNoSlices=0;
484 EAxis workHeaderAxis, newHeaderAxis;
485 G4bool isNewVoxel = false;
486
487 G4double currentDistance = currentStep;
488
489 // Determine if end of Step within current voxel
490 //
491 for (depth=0; depth<fVoxelDepth; ++depth)
492 {
493 targetPoint = localPoint+localDirection*currentDistance;
494 newDistance = currentDistance;
495 workHeader = fVoxelHeaderStack[depth];
496 workHeaderAxis = fVoxelAxisStack[depth];
497 workNodeNo = fVoxelNodeNoStack[depth];
498 workNodeWidth = fVoxelSliceWidthStack[depth];
499 workMinExtent = workHeader->GetMinExtent();
500 workCoord = targetPoint(workHeaderAxis);
501 minVal = workMinExtent+workNodeNo*workNodeWidth;
502
503 if ( minVal<=workCoord+fHalfTolerance )
504 {
505 maxVal = minVal+workNodeWidth;
506 if ( maxVal<=workCoord-fHalfTolerance )
507 {
508 // Must consider next voxel
509 //
510 newNodeNo = workNodeNo+1;
511 newHeader = workHeader;
512 newDistance = (maxVal-localPoint(workHeaderAxis))
513 / localDirection(workHeaderAxis);
514 isNewVoxel = true;
515 newDepth = depth;
516 }
517 }
518 else
519 {
520 newNodeNo = workNodeNo-1;
521 newHeader = workHeader;
522 newDistance = (minVal-localPoint(workHeaderAxis))
523 / localDirection(workHeaderAxis);
524 isNewVoxel = true;
525 newDepth = depth;
526 }
527 currentDistance = newDistance;
528 }
529 targetPoint = localPoint+localDirection*currentDistance;
530
531 // Check if end of Step within collected boundaries of current voxel
532 //
533 depth = fVoxelDepth;
534 {
535 workHeader = fVoxelHeaderStack[depth];
536 workHeaderAxis = fVoxelAxisStack[depth];
537 workNodeNo = fVoxelNodeNoStack[depth];
538 workNodeWidth = fVoxelSliceWidthStack[depth];
539 workMinExtent = workHeader->GetMinExtent();
540 workCoord = targetPoint(workHeaderAxis);
541 minVal = workMinExtent+fVoxelNode->GetMinEquivalentSliceNo()*workNodeWidth;
542
543 if ( minVal<=workCoord+fHalfTolerance )
544 {
545 maxVal = workMinExtent+(fVoxelNode->GetMaxEquivalentSliceNo()+1)
546 *workNodeWidth;
547 if ( maxVal<=workCoord-fHalfTolerance )
548 {
549 newNodeNo = fVoxelNode->GetMaxEquivalentSliceNo()+1;
550 newHeader = workHeader;
551 newDistance = (maxVal-localPoint(workHeaderAxis))
552 / localDirection(workHeaderAxis);
553 isNewVoxel = true;
554 newDepth = depth;
555 }
556 }
557 else
558 {
559 newNodeNo = fVoxelNode->GetMinEquivalentSliceNo()-1;
560 newHeader = workHeader;
561 newDistance = (minVal-localPoint(workHeaderAxis))
562 / localDirection(workHeaderAxis);
563 isNewVoxel = true;
564 newDepth = depth;
565 }
566 currentDistance = newDistance;
567 }
568 if (isNewVoxel)
569 {
570 // Compute new voxel & adjust voxel stack
571 //
572 // newNodeNo=Candidate node no at
573 // newDepth =refinement depth of crossed voxel boundary
574 // newHeader=Header for crossed voxel
575 // newDistance=distance to crossed voxel boundary (along the track)
576 //
577 if ( (newNodeNo<0) || (newNodeNo>=G4int(newHeader->GetNoSlices())))
578 {
579 // Leaving mother volume
580 //
581 isNewVoxel = false;
582 }
583 else
584 {
585 // Compute intersection point on the least refined
586 // voxel boundary that is hit
587 //
588 voxelPoint = localPoint+localDirection*newDistance;
589 fVoxelNodeNoStack[newDepth] = newNodeNo;
590 fVoxelDepth = newDepth;
591 newVoxelNode = 0;
592 while ( !newVoxelNode )
593 {
594 newProxy = newHeader->GetSlice(newNodeNo);
595 if (newProxy->IsNode())
596 {
597 newVoxelNode = newProxy->GetNode();
598 }
599 else
600 {
601 ++fVoxelDepth;
602 newHeader = newProxy->GetHeader();
603 newHeaderAxis = newHeader->GetAxis();
604 newHeaderNoSlices = (G4int)newHeader->GetNoSlices();
605 newHeaderMin = newHeader->GetMinExtent();
606 newHeaderNodeWidth = (newHeader->GetMaxExtent()-newHeaderMin)
607 / newHeaderNoSlices;
608 newNodeNo = G4int( (voxelPoint(newHeaderAxis)-newHeaderMin)
609 / newHeaderNodeWidth );
610 // Rounding protection
611 //
612 if ( newNodeNo<0 )
613 {
614 newNodeNo=0;
615 }
616 else if ( newNodeNo>=newHeaderNoSlices )
617 {
618 newNodeNo = newHeaderNoSlices-1;
619 }
620 // Stack info for stepping
621 //
622 fVoxelAxisStack[fVoxelDepth] = newHeaderAxis;
623 fVoxelNoSlicesStack[fVoxelDepth] = newHeaderNoSlices;
624 fVoxelSliceWidthStack[fVoxelDepth] = newHeaderNodeWidth;
625 fVoxelNodeNoStack[fVoxelDepth] = newNodeNo;
626 fVoxelHeaderStack[fVoxelDepth] = newHeader;
627 }
628 }
629 fVoxelNode = newVoxelNode;
630 }
631 }
632 return isNewVoxel;
633}
634
635// ********************************************************************
636// ComputeSafety
637//
638// Calculates the isotropic distance to the nearest boundary from the
639// specified point in the local coordinate system.
640// The localpoint utilised must be within the current volume.
641// ********************************************************************
642//
645 const G4NavigationHistory& history,
646 const G4double maxLength)
647{
648 G4VPhysicalVolume *motherPhysical, *samplePhysical;
649 G4LogicalVolume *motherLogical;
650 G4VSolid *motherSolid;
651 G4double motherSafety, ourSafety;
652 G4int sampleNo;
653 G4SmartVoxelNode *curVoxelNode;
654 G4long curNoVolumes, contentNo;
655 G4double voxelSafety;
656
657 motherPhysical = history.GetTopVolume();
658 motherLogical = motherPhysical->GetLogicalVolume();
659 motherSolid = motherLogical->GetSolid();
660
661 if( fBestSafety )
662 {
663 return fpVoxelSafety->ComputeSafety( localPoint,*motherPhysical,maxLength );
664 }
665
666 //
667 // Compute mother safety
668 //
669
670 motherSafety = motherSolid->DistanceToOut(localPoint);
671 ourSafety = motherSafety; // Working isotropic safety
672
673 if( motherSafety == 0.0 )
674 {
675#ifdef G4DEBUG_NAVIGATION
676 // Check that point is inside mother volume
677 EInside insideMother = motherSolid->Inside(localPoint);
678
679 if( insideMother == kOutside )
680 {
682 message << "Safety method called for location outside current Volume." << G4endl
683 << "Location for safety is Outside this volume. " << G4endl
684 << "The approximate distance to the solid "
685 << "(safety from outside) is: "
686 << motherSolid->DistanceToIn( localPoint ) << G4endl;
687 message << " Problem occurred with physical volume: "
688 << " Name: " << motherPhysical->GetName()
689 << " Copy No: " << motherPhysical->GetCopyNo() << G4endl
690 << " Local Point = " << localPoint << G4endl;
691 message << " Description of solid: " << G4endl
692 << *motherSolid << G4endl;
693 G4Exception("G4VoxelNavigation::ComputeSafety()", "GeomNav0003",
694 JustWarning, message);
695 }
696
697 // Following check is NOT for an issue - it is only for information
698 // It is allowed that a solid gives approximate safety - even zero.
699 //
700 if( insideMother == kInside ) // && fVerbose )
701 {
702 G4ExceptionDescription messageIn;
703
704 messageIn << " Point is Inside, but safety is Zero ." << G4endl;
705 messageIn << " Inexact safety for volume " << motherPhysical->GetName() << G4endl
706 << " Solid: Name= " << motherSolid->GetName()
707 << " Type= " << motherSolid->GetEntityType() << G4endl;
708 messageIn << " Local point= " << localPoint << G4endl;
709 messageIn << " Solid parameters: " << G4endl << *motherSolid << G4endl;
710 G4Exception("G4VoxelNavigation::ComputeSafety()", "GeomNav0003",
711 JustWarning, messageIn);
712 }
713#endif
714 // if( insideMother != kInside )
715 return 0.0;
716 }
717
718#ifdef G4VERBOSE
719 if( fCheck )
720 {
721 fLogger->ComputeSafetyLog (motherSolid,localPoint,motherSafety,true,true);
722 }
723#endif
724 //
725 // Compute daughter safeties
726 //
727 // Look only inside the current Voxel only (in the first version).
728 //
729 curVoxelNode = fVoxelNode;
730 curNoVolumes = curVoxelNode->GetNoContained();
731
732 for ( contentNo=curNoVolumes-1; contentNo>=0; contentNo-- )
733 {
734 sampleNo = curVoxelNode->GetVolume((G4int)contentNo);
735 samplePhysical = motherLogical->GetDaughter(sampleNo);
736
737 G4AffineTransform sampleTf(samplePhysical->GetRotation(),
738 samplePhysical->GetTranslation());
739 sampleTf.Invert();
740 const G4ThreeVector samplePoint = sampleTf.TransformPoint(localPoint);
741 const G4VSolid* sampleSolid= samplePhysical->GetLogicalVolume()->GetSolid();
742 G4double sampleSafety = sampleSolid->DistanceToIn(samplePoint);
743 if ( sampleSafety<ourSafety )
744 {
745 ourSafety = sampleSafety;
746 }
747#ifdef G4VERBOSE
748 if( fCheck )
749 {
750 fLogger->ComputeSafetyLog(sampleSolid, samplePoint,
751 sampleSafety, false, false);
752 }
753#endif
754 }
755 voxelSafety = ComputeVoxelSafety(localPoint);
756 if ( voxelSafety<ourSafety )
757 {
758 ourSafety = voxelSafety;
759 }
760 return ourSafety;
761}
762
763// ********************************************************************
764// SetVerboseLevel
765// ********************************************************************
766//
768{
769 if( fLogger ) fLogger->SetVerboseLevel(level);
771}
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
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
#define G4endl
Definition: G4ios.hh:57
double dot(const Hep3Vector &) const
HepRotation inverse() const
G4AffineTransform & Invert()
G4ThreeVector TransformPoint(const G4ThreeVector &vec) const
G4ThreeVector TransformAxis(const G4ThreeVector &axis) const
void BlockVolume(const G4int v)
void Enlarge(const G4int nv)
G4bool IsBlocked(const G4int v) const
G4double GetSurfaceTolerance() const
static G4GeometryTolerance * GetInstance()
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
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
G4double GetMinExtent() const
EAxis GetAxis() const
G4int GetMaxEquivalentSliceNo() const
G4int GetVolume(G4int pVolumeNo) const
std::size_t GetNoContained() const
G4int GetMinEquivalentSliceNo() const
G4bool IsNode() const
G4SmartVoxelNode * GetNode() const
G4SmartVoxelHeader * GetHeader() const
const G4RotationMatrix * GetRotation() const
const G4ThreeVector GetTranslation() const
G4LogicalVolume * GetLogicalVolume() const
virtual G4int GetCopyNo() const =0
const G4String & GetName() const
G4String GetName() const
virtual EInside Inside(const G4ThreeVector &p) const =0
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
virtual G4GeometryType GetEntityType() const =0
G4NavigationLogger * fLogger
virtual G4double ComputeSafety(const G4ThreeVector &globalpoint, const G4NavigationHistory &history, const G4double pMaxLength=DBL_MAX)
virtual G4double ComputeStep(const G4ThreeVector &globalPoint, const G4ThreeVector &globalDirection, const G4double currentProposedStepLength, G4double &newSafety, G4NavigationHistory &history, G4bool &validExitNormal, G4ThreeVector &exitNormal, G4bool &exiting, G4bool &entering, G4VPhysicalVolume *(*pBlockedPhysical), G4int &blockedReplicaNo)
std::vector< G4double > fVoxelSliceWidthStack
std::vector< EAxis > fVoxelAxisStack
G4VoxelSafety * fpVoxelSafety
G4SmartVoxelNode * fVoxelNode
G4bool LocateNextVoxel(const G4ThreeVector &localPoint, const G4ThreeVector &localDirection, const G4double currentStep)
virtual ~G4VoxelNavigation()
std::vector< G4int > fVoxelNodeNoStack
G4BlockingList fBList
std::vector< G4SmartVoxelHeader * > fVoxelHeaderStack
std::vector< G4int > fVoxelNoSlicesStack
G4double ComputeVoxelSafety(const G4ThreeVector &localPoint) const
void SetVerboseLevel(G4int level)
void SetVerboseLevel(G4int level)
G4double ComputeSafety(const G4ThreeVector &localPoint, const G4VPhysicalVolume &currentPhysical, G4double maxLength=DBL_MAX)
EAxis
Definition: geomdefs.hh:54
@ kXAxis
Definition: geomdefs.hh:55
EInside
Definition: geomdefs.hh:67
@ kInside
Definition: geomdefs.hh:70
@ kOutside
Definition: geomdefs.hh:68
#define DBL_MAX
Definition: templates.hh:62