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 // class G4PVParameterised implementation 27 // 28 // 29.07.95, P.Kent - first non-stub version 29 // ---------------------------------------------------------------------- 30 31 #include "G4PVParameterised.hh" 32 #include "G4VPVParameterisation.hh" 33 #include "G4AffineTransform.hh" 34 #include "G4UnitsTable.hh" 35 #include "G4VSolid.hh" 36 #include "G4LogicalVolume.hh" 37 38 // ---------------------------------------------------------------------- 39 // Constructor 40 // 41 G4PVParameterised::G4PVParameterised( const G4String& pName, 42 G4LogicalVolume* pLogical, 43 G4VPhysicalVolume* pMotherPhysical, 44 const EAxis pAxis, 45 const G4int nReplicas, 46 G4VPVParameterisation* pParam, 47 G4bool pSurfChk ) 48 : G4PVReplica(pName, nReplicas, pAxis, pLogical, 49 pMotherPhysical != nullptr ? pMotherPhysical->GetLogicalVolume() : nullptr ), 50 fparam(pParam) 51 { 52 G4LogicalVolume* motherLogical= pMotherPhysical != nullptr ? 53 pMotherPhysical->GetLogicalVolume() : nullptr; 54 55 SetMotherLogical( motherLogical ); 56 if( motherLogical != nullptr ) 57 { 58 // Registration moved here to ensure that the volume is recognised as Parameterised 59 motherLogical->AddDaughter(this); 60 } 61 62 #ifdef G4VERBOSE 63 if ((pMotherPhysical != nullptr) && (pMotherPhysical->IsParameterised())) 64 { 65 std::ostringstream message, hint; 66 message << "A parameterised volume is being placed" << G4endl 67 << "inside another parameterised volume !"; 68 hint << "To make sure that no overlaps are generated," << G4endl 69 << "you should verify the mother replicated shapes" << G4endl 70 << "are of the same type and dimensions." << G4endl 71 << " Mother physical volume: " << pMotherPhysical->GetName() << G4endl 72 << " Parameterised volume: " << pName << G4endl 73 << " (To switch this warning off, compile with G4_NO_VERBOSE)"; 74 G4Exception("G4PVParameterised::G4PVParameterised()", "GeomVol1002", 75 JustWarning, message, G4String(hint.str())); 76 } 77 #endif 78 if (pSurfChk) { CheckOverlaps(); } 79 } 80 81 // ---------------------------------------------------------------------- 82 // Constructor 83 // 84 G4PVParameterised::G4PVParameterised( const G4String& pName, 85 G4LogicalVolume* pLogical, 86 G4LogicalVolume* pMotherLogical, 87 const EAxis pAxis, 88 const G4int nReplicas, 89 G4VPVParameterisation* pParam, 90 G4bool pSurfChk ) 91 : G4PVReplica(pName, nReplicas, pAxis, pLogical, pMotherLogical ), 92 fparam(pParam) 93 { 94 SetMotherLogical( pMotherLogical ); 95 if( pMotherLogical != nullptr ) 96 { 97 // Registration moved here to ensure that the volume is recognised as Parameterised 98 pMotherLogical->AddDaughter(this); 99 } 100 if (pSurfChk) { CheckOverlaps(); } 101 } 102 103 // ---------------------------------------------------------------------- 104 // Fake default constructor - sets only member data and allocates memory 105 // for usage restricted to object persistency. 106 // 107 G4PVParameterised::G4PVParameterised( __void__& a ) 108 : G4PVReplica(a) 109 { 110 } 111 112 // ---------------------------------------------------------------------- 113 // Destructor 114 // 115 G4PVParameterised::~G4PVParameterised() = default; 116 117 // ---------------------------------------------------------------------- 118 // GetParameterisation 119 // 120 G4VPVParameterisation* G4PVParameterised::GetParameterisation() const 121 { 122 return fparam; 123 } 124 125 // ---------------------------------------------------------------------- 126 // IsParameterised 127 // 128 G4bool G4PVParameterised::IsParameterised() const 129 { 130 return true; 131 } 132 133 // ---------------------------------------------------------------------- 134 // VolumeType 135 // 136 EVolume G4PVParameterised::VolumeType() const 137 { 138 return kParameterised; 139 } 140 141 // ---------------------------------------------------------------------- 142 // GetReplicationData 143 // 144 void G4PVParameterised::GetReplicationData( EAxis& axis, 145 G4int& nReplicas, 146 G4double& width, 147 G4double& offset, 148 G4bool& consuming) const 149 { 150 axis = faxis; 151 nReplicas = fnReplicas; 152 width = fwidth; 153 offset = foffset; 154 consuming = false; 155 } 156 157 // ---------------------------------------------------------------------- 158 // SetRegularStructureId 159 // 160 void G4PVParameterised::SetRegularStructureId( G4int code ) 161 { 162 G4PVReplica::SetRegularStructureId( code ); 163 // To undertake additional preparation, a derived volume must 164 // redefine this method, while calling also the above method 165 } 166 167 168 // ---------------------------------------------------------------------- 169 // CheckOverlaps 170 // 171 G4bool 172 G4PVParameterised::CheckOverlaps(G4int res, G4double tol, 173 G4bool verbose, G4int maxErr) 174 { 175 if (res<=0) { return false; } 176 177 G4int trials = 0; 178 G4bool retval = false; 179 G4VSolid *solidA = nullptr, *solidB = nullptr; 180 G4LogicalVolume* motherLog = GetMotherLogical(); 181 G4VSolid *motherSolid = motherLog->GetSolid(); 182 std::vector<G4ThreeVector> points; 183 184 if (verbose) 185 { 186 G4cout << "Checking overlaps for parameterised volume " 187 << GetName() << " ... "; 188 } 189 190 for (auto i=0; i<GetMultiplicity(); ++i) 191 { 192 solidA = fparam->ComputeSolid(i, this); 193 solidA->ComputeDimensions(fparam, i, this); 194 fparam->ComputeTransformation(i, this); 195 196 // Create the transformation from daughter to mother 197 // 198 G4AffineTransform Tm( GetRotation(), GetTranslation() ); 199 200 // Generate random points on surface according to the given resolution, 201 // transform them to the mother's coordinate system and if no overlaps 202 // with the mother volume, cache them in a vector for later use with 203 // the daughters 204 // 205 for (auto n=0; n<res; ++n) 206 { 207 G4ThreeVector mp = Tm.TransformPoint(solidA->GetPointOnSurface()); 208 209 // Checking overlaps with the mother volume 210 // 211 if (motherSolid->Inside(mp)==kOutside) 212 { 213 G4double distin = motherSolid->DistanceToIn(mp); 214 if (distin > tol) 215 { 216 ++trials; retval = true; 217 std::ostringstream message; 218 message << "Overlap with mother volume !" << G4endl 219 << " Overlap is detected for volume " 220 << GetName() << ", parameterised instance: " << i << G4endl 221 << " with its mother volume " 222 << motherLog->GetName() << G4endl 223 << " at mother local point " << mp << ", " 224 << "overlapping by at least: " 225 << G4BestUnit(distin, "Length"); 226 if (trials>=maxErr) 227 { 228 message << G4endl 229 << "NOTE: Reached maximum fixed number -" << maxErr 230 << "- of overlaps reports for this volume !"; 231 } 232 G4Exception("G4PVParameterised::CheckOverlaps()", 233 "GeomVol1002", JustWarning, message); 234 if (trials>=maxErr) { return true; } 235 } 236 } 237 points.push_back(mp); 238 } 239 240 // Checking overlaps with each other parameterised instance 241 // 242 for (auto j=i+1; j<GetMultiplicity(); ++j) 243 { 244 solidB = fparam->ComputeSolid(j,this); 245 solidB->ComputeDimensions(fparam, j, this); 246 fparam->ComputeTransformation(j, this); 247 248 // Create the transformation for daughter volume 249 // 250 G4AffineTransform Td( GetRotation(), GetTranslation() ); 251 252 for (const auto & point : points) 253 { 254 // Transform each point according to daughter's frame 255 // 256 G4ThreeVector md = Td.InverseTransformPoint(point); 257 258 if (solidB->Inside(md)==kInside) 259 { 260 G4double distout = solidB->DistanceToOut(md); 261 if (distout > tol) 262 { 263 ++trials; retval = true; 264 std::ostringstream message; 265 message << "Overlap within parameterised volumes !" << G4endl 266 << " Overlap is detected for volume " 267 << GetName() << ", parameterised instance: " << i << G4endl 268 << " with parameterised volume instance: " << j 269 << G4endl 270 << " at local point " << md << ", " 271 << "overlapping by at least: " 272 << G4BestUnit(distout, "Length") 273 << ", related to volume instance: " << j << "."; 274 if (trials>=maxErr) 275 { 276 message << G4endl 277 << "NOTE: Reached maximum fixed number -" << maxErr 278 << "- of overlaps reports for this volume !"; 279 } 280 G4Exception("G4PVParameterised::CheckOverlaps()", 281 "GeomVol1002", JustWarning, message); 282 if (trials>=maxErr) { return true; } 283 } 284 } 285 } 286 } 287 } 288 if (verbose) 289 { 290 G4cout << "OK! " << G4endl; 291 } 292 293 return retval; 294 } 295