source: trunk/source/geometry/navigation/include/G4VoxelSafety.hh@ 1355

Last change on this file since 1355 was 1350, checked in by garnier, 15 years ago

update to last version 4.9.4

File size: 4.6 KB
RevLine 
[1350]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// $Id: G4VoxelSafety.hh,v 1.2 2010/11/04 17:32:46 japost Exp $
27// GEANT4 tag $Name: geant4-09-04-ref-00 $
28//
29// class G4VoxelSafety
30//
31// Class description:
32//
33// Utility for isotropic safety in volumes containing only G4PVPlacement
34// daughter volumes for which voxels have been constructed.
35
36// History:
37// - Created. John Apostolakis, 30 April 2010
38// --------------------------------------------------------------------
39#ifndef G4VOXELSAFETY_HH
40#define G4VOXELSAFETY_HH
41
42#include "geomdefs.hh"
43#include "G4NavigationHistory.hh"
44#include "G4AffineTransform.hh"
45#include "G4VPhysicalVolume.hh"
46#include "G4LogicalVolume.hh"
47#include "G4VSolid.hh"
48#include "G4ThreeVector.hh"
49
50#include "G4BlockingList.hh"
51
52// #include "G4AuxiliaryNavServices.hh"
53
54// Required for voxel handling & voxel stack
55//
56#include <vector>
57
58class G4SmartVoxelNode;
59class G4SmartVoxelHeader;
60
61class G4VoxelSafety
62{
63 public: // with description
64
65 G4VoxelSafety();
66 ~G4VoxelSafety();
67
68 G4SmartVoxelNode* VoxelLocate( G4SmartVoxelHeader* pHead,
69 const G4ThreeVector& localPoint );
70
71 G4double ComputeSafety( const G4ThreeVector& localPoint,
72 const G4VPhysicalVolume& currentPhysical,
73 G4double maxLength=DBL_MAX );
74
75 inline G4int GetVerboseLevel() const { return fVerbose; }
76 inline void SetVerboseLevel(G4int level) { fVerbose= level; }
77 // Get/Set Verbose(ness) level.
78 // [if level>0 && G4VERBOSE, printout can occur]
79
80 // inline void CheckMode(G4bool mode);
81 // Make extra checks
82
83 protected:
84 G4double SafetyForVoxelHeader( G4SmartVoxelHeader* pHead,
85 const G4ThreeVector& localPoint );
86
87 G4double SafetyForVoxelNode( G4SmartVoxelNode *curVoxelNode,
88 const G4ThreeVector& localPoint );
89
90 G4SmartVoxelNode* VoxelLocateLight( G4SmartVoxelHeader* pHead,
91 const G4ThreeVector& localPoint ) const;
92 private:
93 // BEGIN State
94 // - values used during computation of Safety
95 //
96 G4BlockingList fBlockList;
97 // Blocked volumes
98 G4LogicalVolume* fpMotherLogical;
99
100 //
101 // BEGIN Voxel Stack information
102 //
103
104 G4int fVoxelDepth;
105 // Note: fVoxelDepth==0+ => fVoxelAxisStack(0+) contains axes of voxel
106 // fVoxelDepth==-1 -> not in voxel
107
108 std::vector<EAxis> fVoxelAxisStack;
109 // Voxel axes
110
111 std::vector<G4int> fVoxelNoSlicesStack;
112 // No slices per voxel at each level
113
114 std::vector<G4double> fVoxelSliceWidthStack;
115 // Width of voxels at each level
116
117 std::vector<G4int> fVoxelNodeNoStack;
118 // Node no point is inside at each level
119
120 std::vector<G4SmartVoxelHeader*> fVoxelHeaderStack;
121 // Voxel headers at each level
122
123 G4SmartVoxelNode* fVoxelNode;
124 // Node containing last located point
125
126 //
127 // END Voxel Stack information
128 //
129
130 G4bool fCheck;
131 G4int fVerbose;
132 G4double kCarTolerance;
133};
134
135// #include "G4VoxelSafety.icc"
136
137#endif
Note: See TracBrowser for help on using the repository browser.