|
||||
File indexing completed on 2025-01-18 09:59:19
0001 // 0002 // ******************************************************************** 0003 // * License and Disclaimer * 0004 // * * 0005 // * The Geant4 software is copyright of the Copyright Holders of * 0006 // * the Geant4 Collaboration. It is provided under the terms and * 0007 // * conditions of the Geant4 Software License, included in the file * 0008 // * LICENSE and available at http://cern.ch/geant4/license . These * 0009 // * include a list of copyright holders. * 0010 // * * 0011 // * Neither the authors of this software system, nor their employing * 0012 // * institutes,nor the agencies providing financial support for this * 0013 // * work make any representation or warranty, express or implied, * 0014 // * regarding this software system or assume any liability for its * 0015 // * use. Please see the license in the file LICENSE and URL above * 0016 // * for the full disclaimer and the limitation of liability. * 0017 // * * 0018 // * This code implementation is the result of the scientific and * 0019 // * technical work of the GEANT4 collaboration. * 0020 // * By using, copying, modifying or distributing the software (or * 0021 // * any work based on the software) you agree to acknowledge its * 0022 // * use in resulting scientific publications, and indicate your * 0023 // * acceptance of all terms of the Geant4 Software license. * 0024 // ******************************************************************** 0025 // 0026 // G4VCSGface 0027 // 0028 // Class description: 0029 // 0030 // Definition of the virtual base class G4VCSGface, one side (or face) 0031 // of a CSG-like solid. It should be possible to build a CSG entirely 0032 // out of connecting CSG faces. 0033 // 0034 // Each face has an inside and outside surface, the former represents 0035 // the inside of the volume, the latter, the outside. 0036 // 0037 // Virtual members: 0038 // 0039 // ------------------------------------------------------------------- 0040 // Intersect( const G4ThreeVector& p, const G4ThreeVector& v, 0041 // G4bool outGoing, G4double surfTolerance, 0042 // G4double& distance, G4double& distFromSurface, 0043 // G4ThreeVector& normal, G4bool& allBehind ); 0044 // 0045 // p - (in) position 0046 // v - (in) direction (assumed to be a unit vector) 0047 // outgoing - (in) true, to consider only inside surfaces 0048 // false, to consider only outside surfaces 0049 // distance - (out) distance to intersection 0050 // distFromSurface - (out) distance from surface (along surface normal), 0051 // < 0 if the point is in front of the surface 0052 // normal - (out) normal of surface at intersection point 0053 // allBehind - (out) true, if entire surface is behind normal 0054 // 0055 // return value = true if there is an intersection, 0056 // false if there is no intersection 0057 // (all output arguments undefined) 0058 // 0059 // Determine the distance along a line to the face. 0060 // 0061 // ------------------------------------------------------------------- 0062 // Distance( const G4ThreeVector& p, const G4bool outgoing ); 0063 // 0064 // p - (in) position 0065 // outgoing - (in) true, to consider only inside surfaces 0066 // false, to consider only outside surfaces 0067 // 0068 // return value = distance to closest surface satisifying requirements 0069 // or kInfinity if no such surface exists 0070 // 0071 // Determine the distance of a point from either the inside or outside 0072 // surfaces of the face. 0073 // 0074 // ------------------------------------------------------------------- 0075 // Inside( const G4ThreeVector& p, const G4double tolerance, 0076 // G4double* bestDistance ); 0077 // 0078 // p - (in) position 0079 // tolerance - (in) tolerance defining the bounds of the "kSurface", 0080 // nominally equal to kCarTolerance/2 0081 // bestDistance - (out) distance to closest surface (in or out) 0082 // 0083 // return value = kInside if the point is closest to the inside surface 0084 // kOutside if the point is closest to the outside surface 0085 // kSurface if the point is withing tolerance of the surface 0086 // 0087 // Determine whether a point is inside, outside, or on the surface of 0088 // the face. 0089 // 0090 // ------------------------------------------------------------------- 0091 // Normal( const G4ThreeVector& p, G4double* bestDistance ); 0092 // 0093 // p - (in) position 0094 // bestDistance - (out) distance to closest surface (in or out) 0095 // 0096 // return value = the normal of the surface nearest the point 0097 // 0098 // Return normal of surface closest to the point. 0099 // 0100 // ------------------------------------------------------------------- 0101 // Extent( const G4ThreeVector axis ); 0102 // 0103 // axis - (in) unit vector defining direction 0104 // 0105 // return value = the largest point along the given axis of the 0106 // the face's extent. 0107 // 0108 // ------------------------------------------------------------------- 0109 // CalculateExtent( const EAxis pAxis, 0110 // const G4VoxelLimit& pVoxelLimit, 0111 // const G4AffineTransform& pTransform, 0112 // G4double& min, G4double& max ) 0113 // 0114 // pAxis - (in) The x,y, or z axis in which to check 0115 // the shapes 3D extent against 0116 // pVoxelLimit - (in) Limits along x, y, and/or z axes 0117 // pTransform - (in) A coordinate transformation on which 0118 // to apply to the shape before testing 0119 // min - (out) If the face has any point on its 0120 // surface after tranformation and limits 0121 // along pAxis that is smaller than the value 0122 // of min, than it is used to replace min. 0123 // Undefined if the return value is false. 0124 // max - (out) Same as min, except for the largest 0125 // point. 0126 // Undefined if the return value is false. 0127 // 0128 // return value = true if anything remains of the face 0129 // 0130 // Calculate the extent of the face for the voxel navigator. 0131 // In analogy with CalculateExtent for G4VCSGfaceted, this is 0132 // done in the following steps: 0133 // 0134 // 1. Transform the face using pTranform, an arbitrary 3D 0135 // rotation/offset/reflection 0136 // 2. Clip the face to those boundaries as specified in 0137 // pVoxelLimit. This may include limits in any number 0138 // of x, y, or z axes. 0139 // 3. For each part of the face that remains (there could 0140 // be many separate pieces in general): 0141 // 4. Check to see if the piece overlaps the currently 0142 // existing limits along axis pAxis. For 0143 // pVoxelLimit.IsLimited(pAxis) = false, there are 0144 // no limits. 0145 // 5. For a piece that does overlap, update min/max 0146 // accordingly (within confines of pre-existing 0147 // limits) along the direction pAxis. 0148 // 6. If min/max were updated, return true 0149 // 0150 // ------------------------------------------------------------------- 0151 // G3VCSGface *Clone() 0152 // 0153 // This method is invoked by G4CSGfaceted during the copy constructor 0154 // or the assignment operator. Its purpose is to return a pointer 0155 // (of type G4VCSGface) to a duplicate copy of the face. 0156 // The implementation is straight forward for inherited classes. Example: 0157 // 0158 // G4VCSGface G4PolySideFace::Clone() { return new G4PolySideFace(*this); } 0159 // 0160 // Of course, this assumes the copy constructor of G4PolySideFace is 0161 // correctly implemented. 0162 // 0163 // Implementation notes: 0164 // * distance. 0165 // The meaning of distance includes the boundaries of the face. 0166 // For example, for a rectangular, planer face: 0167 // 0168 // A | B | C 0169 // | | 0170 // -------+--------------+----- 0171 // D | I | E 0172 // | | 0173 // -------+--------------+----- 0174 // F | G | H 0175 // | | 0176 // 0177 // A, C, F, and H: closest distance is the distance to 0178 // the adjacent corner. 0179 // 0180 // B, D, E, and G: closest distance is the distance to 0181 // the adjacent line. 0182 // 0183 // I: normal distance to plane 0184 // 0185 // For non-planer faces, one can use the normal to decide when 0186 // a point falls off the edge and then act accordingly. 0187 // 0188 // 0189 // Usage: 0190 // 0191 // A CSG shape can be defined by putting together any number of generic 0192 // faces, as long as the faces cover the entire surface of the shape 0193 // without overlapping. 0194 // 0195 // G4VSolid::CalculateExtent 0196 // 0197 // Define unit vectors along the specified transform axis. 0198 // Use the inverse of the specified coordinate transformation to rotate 0199 // these unit vectors. Loop over each face, call face->Extent, and save 0200 // the maximum value. 0201 // 0202 // G4VSolid::Inside 0203 // 0204 // To decide if a point is inside, outside, or on the surface of the shape, 0205 // loop through all faces, and find the answer from face->Inside which gives 0206 // a value of "bestDistance" smaller than any other. While looping, if any 0207 // face->Inside returns kSurface, this value can be returned immediately. 0208 // 0209 // EInside answer; 0210 // G4VCSGface *face = faces; 0211 // G4double best = kInfinity; 0212 // do { 0213 // G4double distance; 0214 // EInside result = (*face)->Inside( p, kCarTolerance/2, distance ); 0215 // if (result == kSurface) return kSurface; 0216 // if (distance < best) { 0217 // best = distance; 0218 // answer = result; 0219 // } 0220 // } while( ++face < faces + numFaces ); 0221 // 0222 // return(answer); 0223 // 0224 // G4VSolid::SurfaceNormal 0225 // 0226 // Loop over all faces, call face->Normal, and return the normal to the face 0227 // that is closest to the point. 0228 // 0229 // G4VSolid::DistanceToIn(p) 0230 // 0231 // Loop over all faces, invoking face->Distance with outgoing = false, 0232 // and save the answer that is smallest. 0233 // 0234 // G4VSolid::DistanceToIn(p,v) 0235 // 0236 // Loop over all faces, invoking face->Intersect with outgoing = false, 0237 // and save the answer that is smallest. 0238 // 0239 // G4VSolid::DistanceToOut(p) 0240 // 0241 // Loop over all faces, invoking face->Distance with outgoing = true, 0242 // and save the answer that is smallest. 0243 // 0244 // G4VSolid::DistanceToOut(p,v) 0245 // 0246 // Loop over all faces, invoking face->Intersect with outgoing = true, 0247 // and save the answer that is smallest. If there is more than one answer, 0248 // or if allBehind is false for the one answer, return validNorm as false. 0249 0250 // Author: David C. Williams (davidw@scipp.ucsc.edu) 0251 // -------------------------------------------------------------------- 0252 #ifndef G4VCSGFACE_HH 0253 #define G4VCSGFACE_HH 0254 0255 #include "G4Types.hh" 0256 #include "G4ThreeVector.hh" 0257 #include "geomdefs.hh" 0258 #include "G4VSolid.hh" 0259 0260 class G4VoxelLimits; 0261 class G4AffineTransform; 0262 class G4SolidExtentList; 0263 0264 class G4VCSGface 0265 { 0266 public: 0267 0268 G4VCSGface() = default; 0269 virtual ~G4VCSGface() = default; 0270 0271 virtual G4bool Intersect( const G4ThreeVector& p, const G4ThreeVector& v, 0272 G4bool outgoing, G4double surfTolerance, 0273 G4double& distance, G4double& distFromSurface, 0274 G4ThreeVector& normal, G4bool& allBehind ) = 0; 0275 0276 virtual G4double Distance( const G4ThreeVector& p, G4bool outgoing ) = 0; 0277 0278 virtual EInside Inside( const G4ThreeVector& p, G4double tolerance, 0279 G4double* bestDistance ) = 0; 0280 0281 virtual G4ThreeVector Normal( const G4ThreeVector& p, 0282 G4double* bestDistance ) = 0; 0283 0284 virtual G4double Extent( const G4ThreeVector axis ) = 0; 0285 0286 virtual void CalculateExtent( const EAxis axis, 0287 const G4VoxelLimits& voxelLimit, 0288 const G4AffineTransform& tranform, 0289 G4SolidExtentList& extentList ) = 0; 0290 0291 virtual G4VCSGface* Clone() = 0; 0292 0293 virtual G4double SurfaceArea() = 0; 0294 virtual G4ThreeVector GetPointOnFace() = 0; 0295 }; 0296 0297 #endif
[ Source navigation ] | [ Diff markup ] | [ Identifier search ] | [ general search ] |
This page was automatically generated by the 2.3.7 LXR engine. The LXR team |