Geant4 Cross Reference |
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 // Implementation for G4UTrd wrapper class 27 // 28 // 13.09.13 G.Cosmo, CERN/PH 29 // -------------------------------------------------------------------- 30 31 #include "G4Trd.hh" 32 #include "G4UTrd.hh" 33 34 #if ( defined(G4GEOM_USE_USOLIDS) || defined(G4GEOM_USE_PARTIAL_USOLIDS) ) 35 36 #include "G4AffineTransform.hh" 37 #include "G4VPVParameterisation.hh" 38 #include "G4BoundingEnvelope.hh" 39 40 using namespace CLHEP; 41 42 ///////////////////////////////////////////////////////////////////////// 43 // 44 // Constructor - check & set half widths 45 // 46 G4UTrd::G4UTrd(const G4String& pName, 47 G4double pdx1, G4double pdx2, 48 G4double pdy1, G4double pdy2, 49 G4double pdz) 50 : Base_t(pName, pdx1, pdx2, pdy1, pdy2, pdz) 51 { 52 } 53 54 /////////////////////////////////////////////////////////////////////// 55 // 56 // Fake default constructor - sets only member data and allocates memory 57 // for usage restricted to object persistency. 58 // 59 G4UTrd::G4UTrd( __void__& a ) 60 : Base_t(a) 61 { 62 } 63 64 ////////////////////////////////////////////////////////////////////////// 65 // 66 // Destructor 67 // 68 G4UTrd::~G4UTrd() = default; 69 70 ////////////////////////////////////////////////////////////////////////// 71 // 72 // Copy constructor 73 // 74 G4UTrd::G4UTrd(const G4UTrd& rhs) 75 : Base_t(rhs) 76 { 77 } 78 79 ////////////////////////////////////////////////////////////////////////// 80 // 81 // Assignment operator 82 // 83 G4UTrd& G4UTrd::operator = (const G4UTrd& rhs) 84 { 85 // Check assignment to self 86 // 87 if (this == &rhs) { return *this; } 88 89 // Copy base class data 90 // 91 Base_t::operator=(rhs); 92 93 return *this; 94 } 95 96 ////////////////////////////////////////////////////////////////////////// 97 // 98 // Accessors & modifiers 99 100 G4double G4UTrd::GetXHalfLength1() const 101 { 102 return dx1(); 103 } 104 G4double G4UTrd::GetXHalfLength2() const 105 { 106 return dx2(); 107 } 108 G4double G4UTrd::GetYHalfLength1() const 109 { 110 return dy1(); 111 } 112 G4double G4UTrd::GetYHalfLength2() const 113 { 114 return dy2(); 115 } 116 G4double G4UTrd::GetZHalfLength() const 117 { 118 return dz(); 119 } 120 121 void G4UTrd::SetXHalfLength1(G4double val) 122 { 123 Base_t::SetXHalfLength1(val); 124 fRebuildPolyhedron = true; 125 } 126 void G4UTrd::SetXHalfLength2(G4double val) 127 { 128 Base_t::SetXHalfLength2(val); 129 fRebuildPolyhedron = true; 130 } 131 void G4UTrd::SetYHalfLength1(G4double val) 132 { 133 Base_t::SetYHalfLength1(val); 134 fRebuildPolyhedron = true; 135 } 136 void G4UTrd::SetYHalfLength2(G4double val) 137 { 138 Base_t::SetYHalfLength2(val); 139 fRebuildPolyhedron = true; 140 } 141 void G4UTrd::SetZHalfLength(G4double val) 142 { 143 Base_t::SetZHalfLength(val); 144 fRebuildPolyhedron = true; 145 } 146 void G4UTrd::SetAllParameters(G4double pdx1, G4double pdx2, 147 G4double pdy1, G4double pdy2, G4double pdz) 148 { 149 Base_t::SetAllParameters(pdx1, pdx2, pdy1, pdy2, pdz); 150 fRebuildPolyhedron = true; 151 } 152 153 ///////////////////////////////////////////////////////////////////////// 154 // 155 // Dispatch to parameterisation for replication mechanism dimension 156 // computation & modification. 157 // 158 void G4UTrd::ComputeDimensions( G4VPVParameterisation* p, 159 const G4int n, 160 const G4VPhysicalVolume* pRep) 161 { 162 p->ComputeDimensions(*(G4Trd*)this,n,pRep); 163 } 164 165 ////////////////////////////////////////////////////////////////////////// 166 // 167 // Make a clone of the object 168 169 G4VSolid* G4UTrd::Clone() const 170 { 171 return new G4UTrd(*this); 172 } 173 174 ////////////////////////////////////////////////////////////////////////// 175 // 176 // Get bounding box 177 178 void G4UTrd::BoundingLimits(G4ThreeVector& pMin, G4ThreeVector& pMax) const 179 { 180 static G4bool checkBBox = true; 181 182 G4double dx1 = GetXHalfLength1(); 183 G4double dx2 = GetXHalfLength2(); 184 G4double dy1 = GetYHalfLength1(); 185 G4double dy2 = GetYHalfLength2(); 186 G4double dz = GetZHalfLength(); 187 188 G4double xmax = std::max(dx1,dx2); 189 G4double ymax = std::max(dy1,dy2); 190 pMin.set(-xmax,-ymax,-dz); 191 pMax.set( xmax, ymax, dz); 192 193 // Check correctness of the bounding box 194 // 195 if (pMin.x() >= pMax.x() || pMin.y() >= pMax.y() || pMin.z() >= pMax.z()) 196 { 197 std::ostringstream message; 198 message << "Bad bounding box (min >= max) for solid: " 199 << GetName() << " !" 200 << "\npMin = " << pMin 201 << "\npMax = " << pMax; 202 G4Exception("G4UTrd::BoundingLimits()", "GeomMgt0001", 203 JustWarning, message); 204 StreamInfo(G4cout); 205 } 206 207 // Check consistency of bounding boxes 208 // 209 if (checkBBox) 210 { 211 U3Vector vmin, vmax; 212 Extent(vmin,vmax); 213 if (std::abs(pMin.x()-vmin.x()) > kCarTolerance || 214 std::abs(pMin.y()-vmin.y()) > kCarTolerance || 215 std::abs(pMin.z()-vmin.z()) > kCarTolerance || 216 std::abs(pMax.x()-vmax.x()) > kCarTolerance || 217 std::abs(pMax.y()-vmax.y()) > kCarTolerance || 218 std::abs(pMax.z()-vmax.z()) > kCarTolerance) 219 { 220 std::ostringstream message; 221 message << "Inconsistency in bounding boxes for solid: " 222 << GetName() << " !" 223 << "\nBBox min: wrapper = " << pMin << " solid = " << vmin 224 << "\nBBox max: wrapper = " << pMax << " solid = " << vmax; 225 G4Exception("G4UTrd::BoundingLimits()", "GeomMgt0001", 226 JustWarning, message); 227 checkBBox = false; 228 } 229 } 230 } 231 232 ////////////////////////////////////////////////////////////////////////// 233 // 234 // Calculate extent under transform and specified limit 235 236 G4bool 237 G4UTrd::CalculateExtent(const EAxis pAxis, 238 const G4VoxelLimits& pVoxelLimit, 239 const G4AffineTransform& pTransform, 240 G4double& pMin, G4double& pMax) const 241 { 242 G4ThreeVector bmin, bmax; 243 G4bool exist; 244 245 // Check bounding box (bbox) 246 // 247 BoundingLimits(bmin,bmax); 248 G4BoundingEnvelope bbox(bmin,bmax); 249 #ifdef G4BBOX_EXTENT 250 if (true) return bbox.CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax); 251 #endif 252 if (bbox.BoundingBoxVsVoxelLimits(pAxis,pVoxelLimit,pTransform,pMin,pMax)) 253 { 254 return exist = pMin < pMax; 255 } 256 257 // Set bounding envelope (benv) and calculate extent 258 // 259 G4double dx1 = GetXHalfLength1(); 260 G4double dx2 = GetXHalfLength2(); 261 G4double dy1 = GetYHalfLength1(); 262 G4double dy2 = GetYHalfLength2(); 263 G4double dz = GetZHalfLength(); 264 265 G4ThreeVectorList baseA(4), baseB(4); 266 baseA[0].set(-dx1,-dy1,-dz); 267 baseA[1].set( dx1,-dy1,-dz); 268 baseA[2].set( dx1, dy1,-dz); 269 baseA[3].set(-dx1, dy1,-dz); 270 baseB[0].set(-dx2,-dy2, dz); 271 baseB[1].set( dx2,-dy2, dz); 272 baseB[2].set( dx2, dy2, dz); 273 baseB[3].set(-dx2, dy2, dz); 274 275 std::vector<const G4ThreeVectorList *> polygons(2); 276 polygons[0] = &baseA; 277 polygons[1] = &baseB; 278 279 G4BoundingEnvelope benv(bmin,bmax,polygons); 280 exist = benv.CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax); 281 return exist; 282 } 283 284 ////////////////////////////////////////////////////////////////////////// 285 // 286 // Create polyhedron for visualization 287 // 288 G4Polyhedron* G4UTrd::CreatePolyhedron() const 289 { 290 return new G4PolyhedronTrd2(GetXHalfLength1(), 291 GetXHalfLength2(), 292 GetYHalfLength1(), 293 GetYHalfLength2(), 294 GetZHalfLength()); 295 } 296 297 #endif // G4GEOM_USE_USOLIDS 298