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 G4PVReplica Implementation 27 // 28 // 29.07.95, P.Kent - First non-stub version 29 // ---------------------------------------------------------------------- 30 31 #include "G4PVReplica.hh" 32 #include "G4LogicalVolume.hh" 33 34 // ---------------------------------------------------------------------- 35 G4PVRManager G4PVReplica::subInstanceManager; 36 // Helping in the use of the class G4PVRManager. 37 38 #define G4MT_copyNo ((subInstanceManager.offset[instanceID]).fcopyNo) 39 // This macro changes the references to fields that are now encapsulated 40 // in the class G4ReplicaData. 41 42 // ---------------------------------------------------------------------- 43 G4PVReplica::G4PVReplica( const G4String& pName, 44 G4LogicalVolume* pLogical, 45 G4VPhysicalVolume* pMother, 46 const EAxis pAxis, 47 const G4int nReplicas, 48 const G4double width, 49 const G4double offset ) 50 : G4VPhysicalVolume(nullptr, G4ThreeVector(), pName, pLogical, pMother) 51 { 52 53 instanceID = subInstanceManager.CreateSubInstance(); 54 55 if ((pMother == nullptr) || (pMother->GetLogicalVolume() == nullptr)) 56 { 57 std::ostringstream message; 58 message << "NULL pointer specified as mother volume." << G4endl 59 << "The world volume cannot be sliced or parameterised !"; 60 G4Exception("G4PVReplica::G4PVReplica()", "GeomVol0002", 61 FatalException, message); 62 return; 63 } 64 G4LogicalVolume* motherLogical = pMother->GetLogicalVolume(); 65 if (pLogical == motherLogical) 66 { 67 G4Exception("G4PVReplica::G4PVReplica()", "GeomVol0002", 68 FatalException, "Cannot place a volume inside itself!"); 69 return; 70 } 71 SetMotherLogical(motherLogical); 72 motherLogical->AddDaughter(this); 73 if (motherLogical->GetNoDaughters() != 1) 74 { 75 std::ostringstream message; 76 message << "Replica or parameterised volume must be the only daughter !" 77 << G4endl 78 << " Mother physical volume: " << pMother->GetName() << G4endl 79 << " Replicated volume: " << pName; 80 G4Exception("G4PVReplica::G4PVReplica()", "GeomVol0002", 81 FatalException, message); 82 return; 83 } 84 CheckAndSetParameters (pAxis, nReplicas, width, offset); 85 } 86 87 // ---------------------------------------------------------------------- 88 G4PVReplica::G4PVReplica( const G4String& pName, 89 G4LogicalVolume* pLogical, 90 G4LogicalVolume* pMotherLogical, 91 const EAxis pAxis, 92 const G4int nReplicas, 93 const G4double width, 94 const G4double offset ) 95 : G4VPhysicalVolume(nullptr, G4ThreeVector(), pName, pLogical, nullptr) 96 { 97 98 instanceID = subInstanceManager.CreateSubInstance(); 99 100 if (pMotherLogical == nullptr) 101 { 102 std::ostringstream message; 103 message << "NULL pointer specified as mother volume for " 104 << pName << "."; 105 G4Exception("G4PVReplica::G4PVReplica()", "GeomVol0002", 106 FatalException, message); 107 return; 108 } 109 if (pLogical == pMotherLogical) 110 { 111 G4Exception("G4PVReplica::G4PVReplica()", "GeomVol0002", 112 FatalException, "Cannot place a volume inside itself!"); 113 return; 114 } 115 116 pMotherLogical->AddDaughter(this); 117 SetMotherLogical(pMotherLogical); 118 if (pMotherLogical->GetNoDaughters() != 1) 119 { 120 std::ostringstream message; 121 message << "Replica or parameterised volume must be the only daughter !" 122 << G4endl 123 << " Mother logical volume: " << pMotherLogical->GetName() 124 << G4endl 125 << " Replicated volume: " << pName; 126 G4Exception("G4PVReplica::G4PVReplica()", "GeomVol0002", 127 FatalException, message); 128 return; 129 } 130 CheckAndSetParameters (pAxis, nReplicas, width, offset); 131 } 132 133 // ---------------------------------------------------------------------- 134 G4PVReplica::G4PVReplica( const G4String& pName, 135 G4int nReplicas, 136 EAxis pAxis, 137 G4LogicalVolume* pLogical, 138 G4LogicalVolume* pMotherLogical ) 139 : G4VPhysicalVolume(nullptr, G4ThreeVector(), pName, pLogical, nullptr) 140 { 141 // Constructor for derived type(s) 142 // Does not set mother volume or register this one in mother volume 143 // ( To allow the correct type to be found in mother->AddDaughter ) 144 145 instanceID = subInstanceManager.CreateSubInstance(); 146 147 if (pMotherLogical == nullptr) 148 { 149 std::ostringstream message; 150 message << "NULL pointer specified as mother volume for " 151 << pName << "."; 152 G4Exception("G4PVReplica::G4PVReplica()", "GeomVol0002", 153 FatalException, message); 154 return; 155 } 156 if (pLogical == pMotherLogical) 157 { 158 G4Exception("G4PVReplica::G4PVReplica()", "GeomVol0002", 159 FatalException, "Cannot place a volume inside itself!"); 160 return; 161 } 162 CheckOnlyDaughter(pMotherLogical); 163 /*** 164 if (pMotherLogical->GetNoDaughters() != 0) 165 { 166 std::ostringstream message; 167 message << "Replica or parameterised volume must be the only daughter !" 168 << G4endl 169 << " Mother logical volume: " << pMotherLogical->GetName() 170 << G4endl 171 << " Replicated volume: " << pName; 172 G4Exception("G4PVReplica::G4PVReplica()", "GeomVol0002", 173 FatalException, message); 174 return; 175 } 176 **/ 177 CheckAndSetParameters (pAxis, nReplicas, 0.0, 0.0); 178 } 179 180 // ---------------------------------------------------------------------- 181 void G4PVReplica::CheckOnlyDaughter(G4LogicalVolume* pMotherLogical) 182 { 183 if (pMotherLogical->GetNoDaughters() != 0) 184 { 185 std::ostringstream message; 186 message << "Replica or parameterised volume must be the only daughter !" 187 << G4endl 188 << " Mother logical volume: " << pMotherLogical->GetName() 189 << G4endl 190 << " Replicated volume: " << this->GetName() << G4endl 191 << " Existing 'sister': " << pMotherLogical->GetDaughter(0) 192 ->GetName(); 193 G4Exception("G4PVReplica::G4PVReplica()", "GeomVol0002", 194 FatalException, message); 195 return; 196 } 197 } 198 199 // ---------------------------------------------------------------------- 200 void G4PVReplica::CheckAndSetParameters( const EAxis pAxis, 201 const G4int nReplicas, 202 const G4double width, 203 const G4double offset) 204 { 205 if (nReplicas<1) 206 { 207 G4Exception("G4PVReplica::CheckAndSetParameters()", "GeomVol0002", 208 FatalException, "Illegal number of replicas."); 209 } 210 fnReplicas=nReplicas; 211 if (width<0) 212 { 213 G4Exception("G4PVReplica::CheckAndSetParameters()", "GeomVol0002", 214 FatalException, "Width must be positive."); 215 } 216 fwidth = width; 217 foffset = offset; 218 faxis = pAxis; 219 220 // Create rotation matrix for phi axis case & check axis is valid 221 // 222 G4RotationMatrix* pRMat = nullptr; 223 switch (faxis) 224 { 225 case kPhi: 226 pRMat = new G4RotationMatrix(); 227 if (pRMat == nullptr) 228 { 229 G4Exception("G4PVReplica::CheckAndSetParameters()", "GeomVol0003", 230 FatalException, "Rotation matrix allocation failed."); 231 } 232 SetRotation(pRMat); 233 break; 234 case kRho: 235 case kXAxis: 236 case kYAxis: 237 case kZAxis: 238 case kUndefined: 239 break; 240 default: 241 G4Exception("G4PVReplica::CheckAndSetParameters()", "GeomVol0002", 242 FatalException, "Unknown axis of replication."); 243 break; 244 } 245 } 246 247 // ---------------------------------------------------------------------- 248 G4PVReplica::G4PVReplica( __void__& a ) 249 : G4VPhysicalVolume(a), faxis(kZAxis), fnReplicas(0), fwidth(0.), foffset(0.) 250 { 251 instanceID = subInstanceManager.CreateSubInstance(); 252 } 253 254 // ---------------------------------------------------------------------- 255 G4PVReplica::~G4PVReplica() = default; 256 257 // ---------------------------------------------------------------------- 258 G4bool G4PVReplica::IsMany() const 259 { 260 return false; 261 } 262 263 // ---------------------------------------------------------------------- 264 G4int G4PVReplica::GetCopyNo() const 265 { 266 return G4MT_copyNo; 267 } 268 269 // ---------------------------------------------------------------------- 270 void G4PVReplica::SetCopyNo(G4int newCopyNo) 271 { 272 G4MT_copyNo = newCopyNo; 273 } 274 275 // ---------------------------------------------------------------------- 276 G4bool G4PVReplica::IsReplicated() const 277 { 278 return true; 279 } 280 281 // ---------------------------------------------------------------------- 282 G4bool G4PVReplica::IsParameterised() const 283 { 284 return false; 285 } 286 287 // ---------------------------------------------------------------------- 288 G4VPVParameterisation* G4PVReplica::GetParameterisation() const 289 { 290 return nullptr; 291 } 292 293 // ---------------------------------------------------------------------- 294 G4int G4PVReplica::GetMultiplicity() const 295 { 296 return fnReplicas; 297 } 298 299 // ---------------------------------------------------------------------- 300 EVolume G4PVReplica::VolumeType() const 301 { 302 return kReplica; 303 } 304 305 // ---------------------------------------------------------------------- 306 void G4PVReplica::GetReplicationData( EAxis& axis, 307 G4int& nReplicas, 308 G4double& width, 309 G4double& offset, 310 G4bool& consuming ) const 311 { 312 axis = faxis; 313 nReplicas = fnReplicas; 314 width = fwidth; 315 offset = foffset; 316 consuming = true; 317 } 318 319 // ---------------------------------------------------------------------- 320 G4bool G4PVReplica::IsRegularStructure() const 321 { 322 return (fRegularVolsId != 0); 323 } 324 325 // ---------------------------------------------------------------------- 326 G4int G4PVReplica::GetRegularStructureId() const 327 { 328 return fRegularVolsId; 329 } 330 331 // ---------------------------------------------------------------------- 332 void G4PVReplica::SetRegularStructureId( G4int code ) 333 { 334 fRegularVolsId = code; 335 } 336 337 // ---------------------------------------------------------------------- 338 // Returns the private data instance manager. 339 // 340 const G4PVRManager& G4PVReplica::GetSubInstanceManager() 341 { 342 return subInstanceManager; 343 } 344 345 // ---------------------------------------------------------------------- 346 // This method is similar to the constructor. It is used by each worker 347 // thread to achieve the same effect as that of the master thread exept 348 // to register the new created instance. This method is invoked explicitly. 349 // It does not create a new G4PVReplica instance. It only assigns the value 350 // for the fields encapsulated by the class G4ReplicaData. 351 // 352 void G4PVReplica::InitialiseWorker(G4PVReplica* pMasterObject) 353 { 354 355 G4VPhysicalVolume::InitialiseWorker( pMasterObject, nullptr, G4ThreeVector()); 356 subInstanceManager.SlaveCopySubInstanceArray(); 357 G4MT_copyNo = -1; 358 359 // This call causes "self-assignment" of the input parameters 360 // Issue reported by DRD since TerminateWorker() below can be called 361 // at the same time by another thread. 362 // What we need here is the split-class component of CheckAndSetParameters() 363 // funciton copied here. 364 365 // Create rotation matrix for phi axis case & check axis is valid 366 // 367 G4RotationMatrix* pRMat = nullptr; 368 switch (faxis) 369 { 370 case kPhi: 371 pRMat = new G4RotationMatrix(); 372 if (pRMat == nullptr) 373 { 374 G4Exception("G4PVReplica::InitialiseWorker(...)", "GeomVol0003", 375 FatalException, "Rotation matrix allocation failed."); 376 } 377 SetRotation(pRMat); 378 break; 379 case kRho: 380 case kXAxis: 381 case kYAxis: 382 case kZAxis: 383 case kUndefined: 384 break; 385 default: 386 G4Exception("G4PVReplica::InitialiseWorker(...)", "GeomVol0002", 387 FatalException, "Unknown axis of replication."); 388 break; 389 } 390 } 391 392 // ---------------------------------------------------------------------- 393 // This method is similar to the destructor. It is used by each worker 394 // thread to achieve the partial effect as that of the master thread. 395 // For G4PVReplica instances, it destroys the rotation matrix. 396 // 397 void G4PVReplica::TerminateWorker(G4PVReplica* /*pMasterObject*/) 398 { 399 if ( faxis==kPhi ) 400 { 401 delete GetRotation(); 402 } 403 } 404