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
00035
00036
00037 #include <vector>
00038 #include <set>
00039 #include <algorithm>
00040 #include <iomanip>
00041
00042 #include "G4GeomTestVolume.hh"
00043
00044 #include "G4PhysicalConstants.hh"
00045 #include "G4GeomTestLogger.hh"
00046 #include "G4GeomTestVolPoint.hh"
00047 #include "G4GeomTestSegment.hh"
00048
00049 #include "G4VPhysicalVolume.hh"
00050 #include "G4LogicalVolume.hh"
00051 #include "G4VSolid.hh"
00052
00053
00054
00055
00056 G4GeomTestVolume::G4GeomTestVolume( const G4VPhysicalVolume *theTarget,
00057 G4GeomTestLogger *theLogger,
00058 G4double theTolerance )
00059 : target(theTarget),
00060 logger(theLogger),
00061 tolerance(theTolerance),
00062 extent(theTarget->GetLogicalVolume()->GetSolid()->GetExtent())
00063 {;}
00064
00065
00066
00067
00068 G4GeomTestVolume::~G4GeomTestVolume() {;}
00069
00070
00071
00072
00073 G4double G4GeomTestVolume::GetTolerance() const
00074 {
00075 return tolerance;
00076 }
00077
00078
00079
00080
00081 void G4GeomTestVolume::SetTolerance(G4double val)
00082 {
00083 tolerance = val;
00084 }
00085
00086
00087
00088
00089 void G4GeomTestVolume::TestCartGridXYZ( G4int nx, G4int ny, G4int nz )
00090 {
00091 TestCartGridX( ny, nz );
00092 TestCartGridY( nz, nx );
00093 TestCartGridZ( nx, ny );
00094 }
00095
00096
00097
00098
00099 void G4GeomTestVolume::TestCartGridX( G4int ny, G4int nz )
00100 {
00101 TestCartGrid( G4ThreeVector(0,1,0), G4ThreeVector(0,0,1),
00102 G4ThreeVector(1,0,0), ny, nz );
00103 }
00104
00105
00106
00107
00108 void G4GeomTestVolume::TestCartGridY( G4int nz, G4int nx )
00109 {
00110 TestCartGrid( G4ThreeVector(0,0,1), G4ThreeVector(1,0,0),
00111 G4ThreeVector(0,1,0), nz, nx );
00112 }
00113
00114
00115
00116
00117 void G4GeomTestVolume::TestCartGridZ( G4int nx, G4int ny )
00118 {
00119 TestCartGrid( G4ThreeVector(1,0,0), G4ThreeVector(0,1,0),
00120 G4ThreeVector(0,0,1), nx, ny );
00121 }
00122
00123
00124
00125
00126 void G4GeomTestVolume::TestRecursiveCartGrid( G4int nx, G4int ny, G4int nz,
00127 G4int slevel, G4int depth )
00128 {
00129
00130
00131
00132
00133 if (depth == 0) return;
00134 if (depth != -1) depth--;
00135 if (slevel != 0) slevel--;
00136
00137
00138
00139
00140
00141 if ( (!target->IsReplicated()) && (slevel==0) )
00142 {
00143 TestCartGridXYZ( nx, ny, nz );
00144 ReportErrors();
00145 }
00146
00147
00148
00149
00150 std::set<const G4LogicalVolume *> tested;
00151
00152 const G4LogicalVolume *logical = target->GetLogicalVolume();
00153 G4int nDaughter = logical->GetNoDaughters();
00154 G4int iDaughter;
00155 for( iDaughter=0; iDaughter<nDaughter; ++iDaughter )
00156 {
00157 const G4VPhysicalVolume *daughter =
00158 logical->GetDaughter(iDaughter);
00159 const G4LogicalVolume *daughterLogical =
00160 daughter->GetLogicalVolume();
00161
00162
00163
00164
00165 if (daughterLogical->GetNoDaughters() == 0) continue;
00166
00167
00168
00169
00170 std::pair<std::set<const G4LogicalVolume *>::iterator,G4bool>
00171 there = tested.insert(daughterLogical);
00172 if (!there.second) continue;
00173
00174
00175
00176
00177 G4GeomTestVolume vTest( daughter, logger, tolerance );
00178 vTest.TestRecursiveCartGrid( nx,ny,nz,slevel,depth );
00179 }
00180 }
00181
00182
00183
00184
00185 void
00186 G4GeomTestVolume::TestRecursiveCylinder( G4int nPhi, G4int nZ, G4int nRho,
00187 G4double fracZ, G4double fracRho,
00188 G4bool usePhi,
00189 G4int slevel, G4int depth )
00190 {
00191
00192
00193
00194
00195 if (depth == 0) return;
00196 if (depth != -1) depth--;
00197 if (slevel != 0) slevel--;
00198
00199
00200
00201
00202
00203 if ( (!target->IsReplicated()) && (slevel==0) )
00204 {
00205 TestCylinder( nPhi, nZ, nRho, fracZ, fracRho, usePhi );
00206 ReportErrors();
00207 }
00208
00209
00210
00211
00212 std::set<const G4LogicalVolume *> tested;
00213
00214 const G4LogicalVolume *logical = target->GetLogicalVolume();
00215 G4int nDaughter = logical->GetNoDaughters();
00216 G4int iDaughter;
00217 for( iDaughter=0; iDaughter<nDaughter; ++iDaughter )
00218 {
00219 const G4VPhysicalVolume *daughter =
00220 logical->GetDaughter(iDaughter);
00221 const G4LogicalVolume *daughterLogical =
00222 daughter->GetLogicalVolume();
00223
00224
00225
00226
00227 if (daughterLogical->GetNoDaughters() == 0) continue;
00228
00229
00230
00231
00232 std::pair<std::set<const G4LogicalVolume *>::iterator,G4bool>
00233 there = tested.insert(daughterLogical);
00234 if (!there.second) continue;
00235
00236
00237
00238
00239 G4GeomTestVolume vTest( daughter, logger, tolerance );
00240 vTest.TestRecursiveCylinder( nPhi, nZ, nRho,
00241 fracZ, fracRho, usePhi,
00242 slevel, depth );
00243 }
00244 }
00245
00246
00247
00248
00249 void G4GeomTestVolume::TestCylinder( G4int nPhi, G4int nZ, G4int nRho,
00250 G4double fracZ, G4double fracRho,
00251 G4bool usePhi )
00252 {
00253
00254
00255
00256 G4double xMax = std::max(extent.GetXmax(),-extent.GetXmin());
00257 G4double yMax = std::max(extent.GetYmax(),-extent.GetYmin());
00258 G4double rhoMax = std::sqrt(xMax*xMax + yMax*yMax);
00259
00260 G4double zMax = extent.GetZmax();
00261 G4double zMin = extent.GetZmin();
00262 G4double zHalfLength = 0.5*(zMax-zMin);
00263 G4double z0 = 0.5*(zMax+zMin);
00264
00265
00266
00267
00268 G4double deltaPhi = 2*pi/G4double(nPhi);
00269
00270 G4double phi = 0;
00271 G4int iPhi = nPhi;
00272 if ((iPhi&1) == 0) iPhi++;
00273 do
00274 {
00275 G4double cosPhi = std::cos(phi);
00276 G4double sinPhi = std::sin(phi);
00277 G4ThreeVector vPhi1(sinPhi,-cosPhi,0);
00278
00279
00280
00281
00282 G4double rho = rhoMax;
00283 G4int iRho = nRho;
00284 do
00285 {
00286 G4ThreeVector p(rho*cosPhi,rho*sinPhi,0);
00287 static G4ThreeVector v(0,0,1);
00288
00289 TestOneLine( p, v );
00290
00291 if (usePhi)
00292 {
00293
00294
00295
00296 G4double zScale = 1.0;
00297 G4int iZ=nZ;
00298 do
00299 {
00300 p.setZ( z0 + zScale*zHalfLength );
00301 TestOneLine(p,vPhi1);
00302 p.setZ( z0 - zScale*zHalfLength );
00303 TestOneLine(p,vPhi1);
00304 } while( zScale *= fracZ, --iZ );
00305 }
00306 } while( rho *= fracRho, --iRho );
00307
00308
00309
00310
00311 G4ThreeVector p(0,0,0);
00312 G4ThreeVector vPhi2(cosPhi,sinPhi,0);
00313
00314 G4double zScale = 1.0;
00315 G4int iZ=nZ;
00316 do
00317 {
00318 p.setZ( z0 + zScale*zHalfLength );
00319
00320 TestOneLine(p,vPhi2);
00321
00322 p.setZ( z0 - zScale*zHalfLength );
00323
00324 TestOneLine(p,vPhi2);
00325 } while( zScale *= fracZ, --iZ );
00326
00327 } while( phi += deltaPhi, --iPhi );
00328 }
00329
00330
00331
00332
00333 void G4GeomTestVolume::TestCartGrid( const G4ThreeVector &theG1,
00334 const G4ThreeVector &theG2,
00335 const G4ThreeVector &theV,
00336 G4int n1, G4int n2 )
00337 {
00338 if (n1 <= 0 || n2 <= 0)
00339 G4Exception( "G4GeomTestVolume::TestCartGrid()", "GeomNav0002",
00340 FatalException, "Arguments n1 and n2 must be >= 1" );
00341
00342 G4ThreeVector xMin( extent.GetXmin(), extent.GetYmin(),
00343 extent.GetZmin() );
00344 G4ThreeVector xMax( extent.GetXmax(), extent.GetYmax(),
00345 extent.GetZmax() );
00346
00347 G4ThreeVector g1(theG1.unit());
00348 G4ThreeVector g2(theG2.unit());
00349 G4ThreeVector v(theV.unit());
00350
00351 G4double gMin1 = g1.dot(xMin);
00352 G4double gMax1 = g1.dot(xMax);
00353
00354 G4double gMin2 = g2.dot(xMin);
00355 G4double gMax2 = g2.dot(xMax);
00356
00357 G4double delta1 = (gMax1-gMin1)/G4double(n1);
00358 G4double delta2 = (gMax2-gMin2)/G4double(n2);
00359
00360 G4int i1, i2;
00361
00362 for(i1=0;i1<=n1;++i1)
00363 {
00364 G4ThreeVector p1 = (gMin1 + G4double(i1)*delta1)*g1;
00365
00366 for(i2=0;i2<=n2;++i2)
00367 {
00368 G4ThreeVector p2 = (gMin2 + G4double(i2)*delta2)*g2;
00369
00370 TestOneLine( p1+p2, v );
00371 }
00372 }
00373 }
00374
00375
00376
00377
00378 void
00379 G4GeomTestVolume::TestRecursiveLine( const G4ThreeVector& p,
00380 const G4ThreeVector& v,
00381 G4int slevel, G4int depth )
00382 {
00383
00384
00385
00386
00387 if (depth == 0) return;
00388 if (depth != -1) depth--;
00389 if (slevel != 0) slevel--;
00390
00391
00392
00393
00394
00395 if ( (!target->IsReplicated()) && (slevel==0) )
00396 {
00397 TestOneLine( p, v );
00398 ReportErrors();
00399 }
00400
00401
00402
00403
00404 std::set<const G4LogicalVolume *> tested;
00405
00406 const G4LogicalVolume *logical = target->GetLogicalVolume();
00407 G4int nDaughter = logical->GetNoDaughters();
00408 G4int iDaughter;
00409 for( iDaughter=0; iDaughter<nDaughter; ++iDaughter )
00410 {
00411 const G4VPhysicalVolume *daughter =
00412 logical->GetDaughter(iDaughter);
00413 const G4LogicalVolume *daughterLogical =
00414 daughter->GetLogicalVolume();
00415
00416
00417
00418
00419 if (daughterLogical->GetNoDaughters() == 0) continue;
00420
00421
00422
00423
00424 std::pair<std::set<const G4LogicalVolume *>::iterator,G4bool>
00425 there = tested.insert(daughterLogical);
00426 if (!there.second) continue;
00427
00428
00429
00430
00431 G4GeomTestVolume vTest( daughter, logger, tolerance );
00432 vTest.TestRecursiveLine( p, v, slevel, depth );
00433 }
00434 }
00435
00436
00437
00438
00439
00440
00441 void G4GeomTestVolume::TestOneLine( const G4ThreeVector &p,
00442 const G4ThreeVector &v )
00443 {
00444
00445
00446
00447 std::vector<G4GeomTestVolPoint> points;
00448
00449
00450
00451
00452 G4GeomTestSegment targetSegment( target->GetLogicalVolume()->GetSolid(),
00453 p, v, logger );
00454
00455
00456
00457
00458 G4int n = targetSegment.GetNumberPoints();
00459 G4int i;
00460 for(i=0;i<n;++i)
00461 {
00462 points.push_back( G4GeomTestVolPoint( targetSegment.GetPoint(i), -1 ) );
00463 }
00464
00465
00466
00467
00468 const G4LogicalVolume *logical = target->GetLogicalVolume();
00469 G4int nDaughter = logical->GetNoDaughters();
00470 G4int iDaughter;
00471 for( iDaughter=0; iDaughter<nDaughter; ++iDaughter)
00472 {
00473 const G4VPhysicalVolume *daughter =
00474 logical->GetDaughter(iDaughter);
00475
00476
00477
00478
00479 const G4RotationMatrix *rotation = daughter->GetFrameRotation();
00480 const G4ThreeVector &translation = daughter->GetFrameTranslation();
00481
00482 G4ThreeVector pLocal = translation + p;
00483 G4ThreeVector vLocal = v;
00484
00485 if (rotation)
00486 {
00487 pLocal = (*rotation)*pLocal;
00488 vLocal = (*rotation)*vLocal;
00489 }
00490
00491
00492
00493
00494 G4GeomTestSegment
00495 daughterSegment( daughter->GetLogicalVolume()->GetSolid(),
00496 pLocal, vLocal, logger );
00497
00498
00499
00500
00501 n = daughterSegment.GetNumberPoints();
00502 for(i=0;i<n;++i)
00503 {
00504 points.push_back( G4GeomTestVolPoint( daughterSegment.GetPoint(i),
00505 iDaughter, translation, rotation ) );
00506 }
00507 }
00508
00509
00510
00511
00512 std::sort( points.begin(), points.end() );
00513
00514
00515
00516
00517
00518
00519
00520
00521
00522
00523
00524
00525
00526 n = points.size();
00527
00528
00529
00530
00531 std::vector<G4bool> checked( n, false );
00532
00533 for(i=0;i<n;++i)
00534 {
00535 if (checked[i]) continue;
00536
00537 G4int iDaug = points[i].GetDaughterIndex();
00538 if (iDaug < 0) continue;
00539
00540
00541
00542
00543 G4double iS = points[i].GetDistance();
00544 G4int j = i;
00545 while(++j<n)
00546 {
00547 if (iDaug == points[j].GetDaughterIndex()) break;
00548 }
00549 if (j>=n)
00550 {
00551
00552
00553
00554 logger->SolidProblem( logical->GetDaughter(iDaug)->
00555 GetLogicalVolume()->GetSolid(),
00556 "Unmatched intersection point",
00557 points[i].GetPosition() );
00558 continue;
00559 }
00560
00561 checked[j] = true;
00562
00563 G4double jS = points[j].GetDistance();
00564
00565
00566
00567
00568 G4int k = i;
00569 while(++k<j)
00570 {
00571 if (checked[k]) continue;
00572
00573 G4bool kEntering = points[k].Entering();
00574 G4double kS = points[k].GetDistance();
00575
00576
00577
00578
00579
00580 G4int kDaug = points[k].GetDaughterIndex();
00581 if (kDaug < 0)
00582 {
00583
00584
00585
00586 if (kS-iS < tolerance && kEntering ) continue;
00587 if (jS-kS < tolerance && (!kEntering)) continue;
00588
00589
00590
00591
00592 std::map<G4long,G4GeomTestOvershootList>::iterator overshoot =
00593 overshoots.find(iDaug);
00594 if (overshoot == overshoots.end())
00595 {
00596 std::pair<std::map<G4long,G4GeomTestOvershootList>::iterator,G4bool>
00597 result =
00598 overshoots.insert( std::pair<const G4long,G4GeomTestOvershootList>
00599 (iDaug,G4GeomTestOvershootList(target,iDaug)) );
00600 overshoot = result.first;
00601 }
00602
00603 if (kEntering)
00604 (*overshoot).second.AddError( points[i].GetPosition(),
00605 points[k].GetPosition() );
00606 else
00607 (*overshoot).second.AddError( points[k].GetPosition(),
00608 points[j].GetPosition() );
00609 }
00610 else
00611 {
00612
00613
00614
00615 if (kS-iS < tolerance && (!kEntering)) continue;
00616 if (jS-kS < tolerance && kEntering ) continue;
00617
00618
00619
00620
00621 G4long key = iDaug < kDaug ?
00622 (iDaug*nDaughter + kDaug) : (kDaug*nDaughter + iDaug);
00623
00624 std::map<G4long,G4GeomTestOverlapList>::iterator overlap =
00625 overlaps.find(key);
00626 if (overlap == overlaps.end())
00627 {
00628 std::pair<std::map<G4long,G4GeomTestOverlapList>::iterator,G4bool>
00629 result =
00630 overlaps.insert( std::pair<const G4long,G4GeomTestOverlapList>
00631 (key,G4GeomTestOverlapList(target,iDaug,kDaug)) );
00632 overlap = result.first;
00633 }
00634
00635 if (points[k].Entering())
00636 (*overlap).second.AddError( points[k].GetPosition(),
00637 points[j].GetPosition() );
00638 else
00639 (*overlap).second.AddError( points[i].GetPosition(),
00640 points[k].GetPosition() );
00641 }
00642 }
00643 }
00644 }
00645
00646
00647
00648
00649 void G4GeomTestVolume::ReportErrors()
00650 {
00651
00652
00653
00654 if (overshoots.empty())
00655 logger->NoProblem("GeomTest: no daughter volume extending outside mother detected.");
00656 else
00657 {
00658 std::map<G4long,G4GeomTestOvershootList>::iterator overshoot =
00659 overshoots.begin();
00660 while( overshoot != overshoots.end() )
00661 {
00662 logger->OvershootingDaughter( &(*overshoot).second );
00663 ++overshoot;
00664 }
00665 }
00666
00667
00668
00669
00670 if (overlaps.empty())
00671 logger->NoProblem("GeomTest: no overlapping daughters detected.");
00672 else
00673 {
00674 std::map<G4long,G4GeomTestOverlapList>::iterator overlap =
00675 overlaps.begin();
00676 while( overlap != overlaps.end() )
00677 {
00678 logger->OverlappingDaughters( &(*overlap).second );
00679 ++overlap;
00680 }
00681 }
00682 }
00683
00684
00685
00686
00687 void G4GeomTestVolume::ClearErrors()
00688 {
00689 overshoots.clear();
00690 overlaps.clear();
00691 }