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