00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031
00032
00033
00034 #include "G4PVPlacement.hh"
00035 #include "G4AffineTransform.hh"
00036 #include "G4UnitsTable.hh"
00037 #include "G4LogicalVolume.hh"
00038 #include "G4VSolid.hh"
00039
00040
00041
00042
00043 G4PVPlacement::G4PVPlacement( G4RotationMatrix *pRot,
00044 const G4ThreeVector &tlate,
00045 const G4String& pName,
00046 G4LogicalVolume *pLogical,
00047 G4VPhysicalVolume *pMother,
00048 G4bool pMany,
00049 G4int pCopyNo,
00050 G4bool pSurfChk )
00051 : G4VPhysicalVolume(pRot,tlate,pName,pLogical,pMother),
00052 fmany(pMany), fallocatedRotM(false), fcopyNo(pCopyNo)
00053 {
00054 if (pMother)
00055 {
00056 G4LogicalVolume* motherLogical = pMother->GetLogicalVolume();
00057 if (pLogical == motherLogical)
00058 {
00059 G4Exception("G4PVPlacement::G4PVPlacement()", "GeomVol0002",
00060 FatalException, "Cannot place a volume inside itself!");
00061 }
00062 SetMotherLogical(motherLogical);
00063 motherLogical->AddDaughter(this);
00064 if (pSurfChk) { CheckOverlaps(); }
00065 }
00066 }
00067
00068
00069
00070
00071 G4PVPlacement::G4PVPlacement( const G4Transform3D &Transform3D,
00072 const G4String& pName,
00073 G4LogicalVolume *pLogical,
00074 G4VPhysicalVolume *pMother,
00075 G4bool pMany,
00076 G4int pCopyNo,
00077 G4bool pSurfChk )
00078 : G4VPhysicalVolume(NewPtrRotMatrix(Transform3D.getRotation().inverse()),
00079 Transform3D.getTranslation(),pName,pLogical,pMother),
00080 fmany(pMany), fcopyNo(pCopyNo)
00081 {
00082 fallocatedRotM = (GetRotation() != 0);
00083 if (pMother)
00084 {
00085 G4LogicalVolume* motherLogical = pMother->GetLogicalVolume();
00086 if (pLogical == motherLogical)
00087 G4Exception("G4PVPlacement::G4PVPlacement()", "GeomVol0002",
00088 FatalException, "Cannot place a volume inside itself!");
00089 SetMotherLogical(motherLogical);
00090 motherLogical->AddDaughter(this);
00091 if (pSurfChk) { CheckOverlaps(); }
00092 }
00093 }
00094
00095
00096
00097
00098
00099
00100 G4PVPlacement::G4PVPlacement( G4RotationMatrix *pRot,
00101 const G4ThreeVector &tlate,
00102 G4LogicalVolume *pCurrentLogical,
00103 const G4String& pName,
00104 G4LogicalVolume *pMotherLogical,
00105 G4bool pMany,
00106 G4int pCopyNo,
00107 G4bool pSurfChk )
00108 : G4VPhysicalVolume(pRot,tlate,pName,pCurrentLogical,0),
00109 fmany(pMany), fallocatedRotM(false), fcopyNo(pCopyNo)
00110 {
00111 if (pCurrentLogical == pMotherLogical)
00112 {
00113 G4Exception("G4PVPlacement::G4PVPlacement()", "GeomVol0002",
00114 FatalException, "Cannot place a volume inside itself!");
00115 }
00116 SetMotherLogical(pMotherLogical);
00117 if (pMotherLogical) { pMotherLogical->AddDaughter(this); }
00118 if ((pSurfChk) && (pMotherLogical)) { CheckOverlaps(); }
00119 }
00120
00121
00122
00123
00124
00125 G4PVPlacement::G4PVPlacement( const G4Transform3D &Transform3D,
00126 G4LogicalVolume *pCurrentLogical,
00127 const G4String& pName,
00128 G4LogicalVolume *pMotherLogical,
00129 G4bool pMany,
00130 G4int pCopyNo,
00131 G4bool pSurfChk )
00132 : G4VPhysicalVolume(0,Transform3D.getTranslation(),pName,pCurrentLogical,0),
00133 fmany(pMany), fcopyNo(pCopyNo)
00134 {
00135 if (pCurrentLogical == pMotherLogical)
00136 {
00137 G4Exception("G4PVPlacement::G4PVPlacement()", "GeomVol0002",
00138 FatalException, "Cannot place a volume inside itself!");
00139 }
00140 SetRotation( NewPtrRotMatrix(Transform3D.getRotation().inverse()) );
00141 fallocatedRotM = (GetRotation() != 0);
00142 SetMotherLogical(pMotherLogical);
00143 if (pMotherLogical) { pMotherLogical->AddDaughter(this); }
00144 if ((pSurfChk) && (pMotherLogical)) { CheckOverlaps(); }
00145 }
00146
00147
00148
00149
00150
00151 G4PVPlacement::G4PVPlacement( __void__& a )
00152 : G4VPhysicalVolume(a), fmany(false), fallocatedRotM(0), fcopyNo(0)
00153 {
00154 }
00155
00156
00157
00158
00159 G4PVPlacement::~G4PVPlacement()
00160 {
00161 if( fallocatedRotM ){ delete frot; }
00162 }
00163
00164
00165
00166
00167 G4bool G4PVPlacement::IsMany() const
00168 {
00169 return fmany;
00170 }
00171
00172
00173
00174
00175 G4int G4PVPlacement::GetCopyNo() const
00176 {
00177 return fcopyNo;
00178 }
00179
00180
00181
00182
00183 void G4PVPlacement::SetCopyNo(G4int newCopyNo)
00184 {
00185 fcopyNo= newCopyNo;
00186 }
00187
00188
00189
00190
00191 G4bool G4PVPlacement::IsReplicated() const
00192 {
00193 return false;
00194 }
00195
00196
00197
00198
00199 G4bool G4PVPlacement::IsParameterised() const
00200 {
00201 return false;
00202 }
00203
00204
00205
00206
00207 G4VPVParameterisation* G4PVPlacement::GetParameterisation() const
00208 {
00209 return 0;
00210 }
00211
00212
00213
00214
00215 void G4PVPlacement::
00216 GetReplicationData( EAxis&, G4int&, G4double&, G4double&, G4bool& ) const
00217 {
00218
00219 }
00220
00221
00222
00223
00224
00225
00226 G4bool G4PVPlacement::IsRegularStructure() const
00227 {
00228 return false;
00229 }
00230
00231
00232
00233
00234
00235
00236 G4int G4PVPlacement::GetRegularStructureId() const
00237 {
00238 return 0;
00239 }
00240
00241
00242
00243
00244 G4bool G4PVPlacement::CheckOverlaps(G4int res, G4double tol, G4bool verbose)
00245 {
00246 if (res<=0) { return false; }
00247
00248 G4VSolid* solid = GetLogicalVolume()->GetSolid();
00249 G4LogicalVolume* motherLog = GetMotherLogical();
00250 if (!motherLog) { return false; }
00251
00252 G4VSolid* motherSolid = motherLog->GetSolid();
00253
00254 if (verbose)
00255 {
00256 G4cout << "Checking overlaps for volume " << GetName() << " ... ";
00257 }
00258
00259
00260
00261 G4AffineTransform Tm( GetRotation(), GetTranslation() );
00262
00263 for (G4int n=0; n<res; n++)
00264 {
00265
00266
00267 G4ThreeVector point = solid->GetPointOnSurface();
00268
00269
00270
00271 G4ThreeVector mp = Tm.TransformPoint(point);
00272
00273
00274
00275 if (motherSolid->Inside(mp)==kOutside)
00276 {
00277 G4double distin = motherSolid->DistanceToIn(mp);
00278 if (distin > tol)
00279 {
00280 std::ostringstream message;
00281 message << "Overlap with mother volume !" << G4endl
00282 << " Overlap is detected for volume "
00283 << GetName() << G4endl
00284 << " with its mother volume "
00285 << motherLog->GetName() << G4endl
00286 << " at mother local point " << mp << ", "
00287 << "overlapping by at least: "
00288 << G4BestUnit(distin, "Length");
00289 G4Exception("G4PVPlacement::CheckOverlaps()",
00290 "GeomVol1002", JustWarning, message);
00291 return true;
00292 }
00293 }
00294
00295
00296
00297 for (G4int i=0; i<motherLog->GetNoDaughters(); i++)
00298 {
00299 G4VPhysicalVolume* daughter = motherLog->GetDaughter(i);
00300
00301 if (daughter == this) { continue; }
00302
00303
00304
00305 G4AffineTransform Td( daughter->GetRotation(),
00306 daughter->GetTranslation() );
00307 G4ThreeVector md = Td.Inverse().TransformPoint(mp);
00308
00309 G4VSolid* daughterSolid = daughter->GetLogicalVolume()->GetSolid();
00310 if (daughterSolid->Inside(md)==kInside)
00311 {
00312 G4double distout = daughterSolid->DistanceToOut(md);
00313 if (distout > tol)
00314 {
00315 std::ostringstream message;
00316 message << "Overlap with volume already placed !" << G4endl
00317 << " Overlap is detected for volume "
00318 << GetName() << G4endl
00319 << " with " << daughter->GetName() << " volume's"
00320 << G4endl
00321 << " local point " << md << ", "
00322 << "overlapping by at least: "
00323 << G4BestUnit(distout,"Length");
00324 G4Exception("G4PVPlacement::CheckOverlaps()",
00325 "GeomVol1002", JustWarning, message);
00326 return true;
00327 }
00328 }
00329
00330
00331
00332
00333 if (n==0)
00334 {
00335
00336
00337
00338 G4ThreeVector dPoint = daughterSolid->GetPointOnSurface();
00339
00340
00341
00342
00343 G4ThreeVector mp2 = Td.TransformPoint(dPoint);
00344 G4ThreeVector msi = Tm.Inverse().TransformPoint(mp2);
00345
00346 if (solid->Inside(msi)==kInside)
00347 {
00348 std::ostringstream message;
00349 message << "Overlap with volume already placed !" << G4endl
00350 << " Overlap is detected for volume "
00351 << GetName() << G4endl
00352 << " apparently fully encapsulating volume "
00353 << daughter->GetName() << G4endl
00354 << " at the same level !";
00355 G4Exception("G4PVPlacement::CheckOverlaps()",
00356 "GeomVol1002", JustWarning, message);
00357 return true;
00358 }
00359 }
00360 }
00361 }
00362
00363 if (verbose)
00364 {
00365 G4cout << "OK! " << G4endl;
00366 }
00367
00368 return false;
00369 }
00370
00371
00372
00373
00374
00375
00376
00377
00378
00379
00380
00381 G4RotationMatrix*
00382 G4PVPlacement::NewPtrRotMatrix(const G4RotationMatrix &RotMat)
00383 {
00384 G4RotationMatrix *pRotMatrix;
00385 if ( RotMat.isIdentity() )
00386 {
00387 pRotMatrix = 0;
00388 }
00389 else
00390 {
00391 pRotMatrix = new G4RotationMatrix(RotMat);
00392 }
00393
00394
00395 return pRotMatrix;
00396 }