| [1316] | 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 | //
|
|---|
| 27 | // SBTvoxel.cc
|
|---|
| 28 | //
|
|---|
| 29 | // Implementation of a batch based voxel test
|
|---|
| 30 | //
|
|---|
| 31 |
|
|---|
| 32 | #include "globals.hh"
|
|---|
| 33 | #include "Randomize.hh"
|
|---|
| 34 |
|
|---|
| 35 | #include "SBTvoxel.hh"
|
|---|
| 36 |
|
|---|
| 37 | #include "SBTVisManager.hh"
|
|---|
| 38 | #include "G4Polyline.hh"
|
|---|
| 39 | #include "G4Circle.hh"
|
|---|
| 40 | #include "G4Color.hh"
|
|---|
| 41 | #include "G4VisAttributes.hh"
|
|---|
| 42 |
|
|---|
| 43 | #include "G4VSolid.hh"
|
|---|
| 44 | #include "G4VoxelLimits.hh"
|
|---|
| 45 | #include "G4AffineTransform.hh"
|
|---|
| 46 | #include "G4GeometryTolerance.hh"
|
|---|
| 47 | #include <iomanip>
|
|---|
| 48 | #include <sstream>
|
|---|
| 49 |
|
|---|
| 50 | #include <time.h>
|
|---|
| 51 |
|
|---|
| 52 | //
|
|---|
| 53 | // Constructor
|
|---|
| 54 | //
|
|---|
| 55 | SBTvoxel::SBTvoxel()
|
|---|
| 56 | {
|
|---|
| 57 | SetDefaults();
|
|---|
| 58 | }
|
|---|
| 59 |
|
|---|
| 60 |
|
|---|
| 61 | //
|
|---|
| 62 | // Destructor
|
|---|
| 63 | //
|
|---|
| 64 | SBTvoxel::~SBTvoxel() {;}
|
|---|
| 65 |
|
|---|
| 66 |
|
|---|
| 67 | //
|
|---|
| 68 | // SetDefaults
|
|---|
| 69 | //
|
|---|
| 70 | // Set default values for test parameters
|
|---|
| 71 | //
|
|---|
| 72 | void SBTvoxel::SetDefaults()
|
|---|
| 73 | {
|
|---|
| 74 | target = G4ThreeVector( 0, 0, 0 );
|
|---|
| 75 | widths = G4ThreeVector( 1*m, 1*m, 1*m );
|
|---|
| 76 |
|
|---|
| 77 | maxVoxels = 100;
|
|---|
| 78 | maxErrors = 20;
|
|---|
| 79 | }
|
|---|
| 80 |
|
|---|
| 81 |
|
|---|
| 82 |
|
|---|
| 83 | //
|
|---|
| 84 | // Debug
|
|---|
| 85 | //
|
|---|
| 86 | // Invoke the CalculateExtent method of the target volume.
|
|---|
| 87 | // This can be particularly useful for debugging.
|
|---|
| 88 | //
|
|---|
| 89 | void SBTvoxel::Debug( const G4VSolid *testVolume, const EAxis axis,
|
|---|
| 90 | const G4VoxelLimits &voxel, const G4AffineTransform &transform,
|
|---|
| 91 | const G4ThreeVector *point ) const
|
|---|
| 92 | {
|
|---|
| 93 | G4double min, max;
|
|---|
| 94 |
|
|---|
| 95 | if (testVolume->CalculateExtent( axis, voxel, transform, min, max ))
|
|---|
| 96 | G4cout << "Solid is intersected with: " << min << " " << max << G4endl;
|
|---|
| 97 | else
|
|---|
| 98 | G4cout << "Voxel misses solid" << G4endl;
|
|---|
| 99 |
|
|---|
| 100 | if (point) testVolume->Inside( *point );
|
|---|
| 101 | }
|
|---|
| 102 |
|
|---|
| 103 |
|
|---|
| 104 | //
|
|---|
| 105 | // Draw
|
|---|
| 106 | //
|
|---|
| 107 | // Draw a test shape and a voxel with arbitrary limits, using
|
|---|
| 108 | // an arbitrary tranformation, and with any number of optional
|
|---|
| 109 | // markers
|
|---|
| 110 | //
|
|---|
| 111 | void SBTvoxel::Draw( const G4VSolid *testVolume,
|
|---|
| 112 | const G4VoxelLimits &voxel, const G4AffineTransform &transform,
|
|---|
| 113 | const G4ThreeVector *points, const G4int numPoints,
|
|---|
| 114 | const EAxis axis, const G4double limits[2],
|
|---|
| 115 | SBTVisManager *visManager ) const
|
|---|
| 116 | {
|
|---|
| 117 | //
|
|---|
| 118 | // Get inverse transform
|
|---|
| 119 | //
|
|---|
| 120 | G4AffineTransform inverseTransform = transform.Inverse();
|
|---|
| 121 |
|
|---|
| 122 | //
|
|---|
| 123 | // Prepare visualization
|
|---|
| 124 | //
|
|---|
| 125 | // visManager->ClearView();
|
|---|
| 126 |
|
|---|
| 127 | //
|
|---|
| 128 | // Draw voxel as a box.
|
|---|
| 129 | //
|
|---|
| 130 | G4VisAttributes blueStuff( G4Color(0,0,1) );
|
|---|
| 131 | G4VisAttributes yuckStuff( G4Color(0,1,0) );
|
|---|
| 132 |
|
|---|
| 133 | static const EAxis axes[3] = { kXAxis, kYAxis, kZAxis };
|
|---|
| 134 | static const G4ThreeVector axisVectors[3] = {G4ThreeVector(1,0,0),
|
|---|
| 135 | G4ThreeVector(0,1,0),
|
|---|
| 136 | G4ThreeVector(0,0,1) };
|
|---|
| 137 |
|
|---|
| 138 | G4bool drawLimits = limits[0] <= limits[1];
|
|---|
| 139 |
|
|---|
| 140 | G4ThreeVector dmin[3], dmax[3], lvec[2];
|
|---|
| 141 | G4int i;
|
|---|
| 142 | for( i=0; i<3; i++ ) {
|
|---|
| 143 | G4double min, max;
|
|---|
| 144 | if (voxel.IsLimited(axes[i])) {
|
|---|
| 145 | min = voxel.GetMinExtent(axes[i]);
|
|---|
| 146 | max = voxel.GetMaxExtent(axes[i]);
|
|---|
| 147 | }
|
|---|
| 148 | else {
|
|---|
| 149 | min = -4.0*m;
|
|---|
| 150 | max = +4.0*m;
|
|---|
| 151 | }
|
|---|
| 152 |
|
|---|
| 153 | dmin[i] = min*axisVectors[i];
|
|---|
| 154 | dmax[i] = max*axisVectors[i];
|
|---|
| 155 |
|
|---|
| 156 | if (drawLimits && axis==axes[i]) {
|
|---|
| 157 | lvec[0] = limits[0]*axisVectors[i];
|
|---|
| 158 | lvec[1] = limits[1]*axisVectors[i];
|
|---|
| 159 | }
|
|---|
| 160 | }
|
|---|
| 161 |
|
|---|
| 162 | G4Transform3D objectTransformation;
|
|---|
| 163 |
|
|---|
| 164 | for( i=0; i<3 ;i++ ) {
|
|---|
| 165 | G4int bitmask;
|
|---|
| 166 | for( bitmask=0; bitmask < 8; bitmask++ ) {
|
|---|
| 167 | if (bitmask&(1<<i)) continue;
|
|---|
| 168 |
|
|---|
| 169 | G4ThreeVector a = dmin[i], b = dmax[i];
|
|---|
| 170 | G4ThreeVector c = a, d = b;
|
|---|
| 171 | G4int j;
|
|---|
| 172 | for( j=0; j<3; j++ ) {
|
|---|
| 173 | if (i==j) continue;
|
|---|
| 174 | G4ThreeVector *use = bitmask&(1<<j) ? dmin+j : dmax+j;
|
|---|
| 175 | a += *use;
|
|---|
| 176 | b += *use;
|
|---|
| 177 | if (drawLimits) {
|
|---|
| 178 | if (axes[j]==axis) use = lvec + (bitmask&(1<<j) ? 0 : 1);
|
|---|
| 179 | c += *use;
|
|---|
| 180 | d += *use;
|
|---|
| 181 | }
|
|---|
| 182 | }
|
|---|
| 183 | G4Polyline polyline;
|
|---|
| 184 | polyline.SetVisAttributes( blueStuff );
|
|---|
| 185 | inverseTransform.ApplyPointTransform( a );
|
|---|
| 186 | inverseTransform.ApplyPointTransform( b );
|
|---|
| 187 | polyline.push_back( a );
|
|---|
| 188 | polyline.push_back( b );
|
|---|
| 189 | visManager->Draw( polyline, objectTransformation);
|
|---|
| 190 |
|
|---|
| 191 | if (drawLimits && axes[i]!=axis) {
|
|---|
| 192 | G4Polyline polyline;
|
|---|
| 193 | polyline.SetVisAttributes( yuckStuff );
|
|---|
| 194 | inverseTransform.ApplyPointTransform( c );
|
|---|
| 195 | inverseTransform.ApplyPointTransform( d );
|
|---|
| 196 | polyline.push_back( c );
|
|---|
| 197 | polyline.push_back( d );
|
|---|
| 198 | visManager->Draw( polyline, objectTransformation );
|
|---|
| 199 | }
|
|---|
| 200 | }
|
|---|
| 201 | }
|
|---|
| 202 |
|
|---|
| 203 |
|
|---|
| 204 | //
|
|---|
| 205 | // Draw points
|
|---|
| 206 | //
|
|---|
| 207 | G4VisAttributes whiteStuff( G4Color(1,1,1) );
|
|---|
| 208 |
|
|---|
| 209 | const G4ThreeVector *thisPoint;
|
|---|
| 210 | for (thisPoint = points; thisPoint < points+numPoints; thisPoint++ ) {
|
|---|
| 211 | G4Circle circle(*thisPoint);
|
|---|
| 212 | circle.SetWorldSize( 5*cm );
|
|---|
| 213 | circle.SetVisAttributes( whiteStuff );
|
|---|
| 214 | visManager->Draw( circle, objectTransformation );
|
|---|
| 215 | }
|
|---|
| 216 |
|
|---|
| 217 |
|
|---|
| 218 | //
|
|---|
| 219 | // This draws the target solid
|
|---|
| 220 | //
|
|---|
| 221 | G4VisAttributes redStuff( G4Color(1,0,0) );
|
|---|
| 222 | visManager->Draw( *testVolume, redStuff, objectTransformation );
|
|---|
| 223 |
|
|---|
| 224 | // visManager->Show();
|
|---|
| 225 | }
|
|---|
| 226 |
|
|---|
| 227 |
|
|---|
| 228 | //
|
|---|
| 229 | // RunTest
|
|---|
| 230 | //
|
|---|
| 231 | // Perform a test on the specified solid
|
|---|
| 232 | //
|
|---|
| 233 | void SBTvoxel::RunTest( const G4VSolid *testVolume, std::ostream &logger )
|
|---|
| 234 | {
|
|---|
| 235 | //
|
|---|
| 236 | // Output test parameters
|
|---|
| 237 | //
|
|---|
| 238 | time_t now = time(0);
|
|---|
| 239 | time(&now);
|
|---|
| 240 | G4String dateTime(ctime(&now));
|
|---|
| 241 |
|
|---|
| 242 | logger << "% SBT voxel logged output " << dateTime;
|
|---|
| 243 | logger << "% target = " << target << G4endl;
|
|---|
| 244 | logger << "% widths = " << widths << G4endl;
|
|---|
| 245 | logger << "% maxVoxels = " << maxVoxels << G4endl;
|
|---|
| 246 | logger << "% maxErrors = " << maxErrors << G4endl;
|
|---|
| 247 |
|
|---|
| 248 | G4int nVoxel = 0,
|
|---|
| 249 | nError = 0;
|
|---|
| 250 |
|
|---|
| 251 | //
|
|---|
| 252 | // Generate a list of 1000 random points inside the solid
|
|---|
| 253 | //
|
|---|
| 254 | G4ThreeVector inside[1000];
|
|---|
| 255 | G4int numInside;
|
|---|
| 256 |
|
|---|
| 257 | GetInsidePoints( testVolume, inside, &numInside, 1000, 100000 );
|
|---|
| 258 |
|
|---|
| 259 | G4RotationMatrix randomRotate1, randomRotate2;
|
|---|
| 260 |
|
|---|
| 261 | for(;;) {
|
|---|
| 262 | //
|
|---|
| 263 | // Generate a random voxel limit.
|
|---|
| 264 | // G4VoxelLimits has no "reset" method, so we need
|
|---|
| 265 | // to create a brand spanking new one each iteration
|
|---|
| 266 | //
|
|---|
| 267 | G4VoxelLimits *voxel = NewRandomVoxel( widths );
|
|---|
| 268 |
|
|---|
| 269 | //
|
|---|
| 270 | // Move it to random positions
|
|---|
| 271 | //
|
|---|
| 272 | G4int j;
|
|---|
| 273 | for( j=0; j < 20; j++ ) {
|
|---|
| 274 | G4ThreeVector offset = j > 0 ? GetRandomPoint() : G4ThreeVector(0,0,0);
|
|---|
| 275 |
|
|---|
| 276 | G4AffineTransform transform( offset );
|
|---|
| 277 |
|
|---|
| 278 | //
|
|---|
| 279 | // Test aligned
|
|---|
| 280 | //
|
|---|
| 281 | if (TestOneVoxel( testVolume, *voxel, transform, inside, numInside, logger )) {
|
|---|
| 282 | if (++nError >= maxErrors) break;
|
|---|
| 283 | }
|
|---|
| 284 |
|
|---|
| 285 | //
|
|---|
| 286 | // Generate a couple random orientations
|
|---|
| 287 | //
|
|---|
| 288 | randomRotate1.rotateZ( G4UniformRand() );
|
|---|
| 289 | transform.SetNetRotation( randomRotate1 );
|
|---|
| 290 |
|
|---|
| 291 | if (TestOneVoxel( testVolume, *voxel, transform, inside, numInside, logger )) {
|
|---|
| 292 | if (++nError >= maxErrors) break;
|
|---|
| 293 | }
|
|---|
| 294 |
|
|---|
| 295 | randomRotate2.rotateX( G4UniformRand() );
|
|---|
| 296 | transform.SetNetRotation( randomRotate2 );
|
|---|
| 297 |
|
|---|
| 298 | if (TestOneVoxel( testVolume, *voxel, transform, inside, numInside, logger )) {
|
|---|
| 299 | if (++nError >= maxErrors) break;
|
|---|
| 300 | }
|
|---|
| 301 |
|
|---|
| 302 | randomRotate2.rotateY( G4UniformRand() );
|
|---|
| 303 | transform.SetNetRotation( randomRotate2 );
|
|---|
| 304 |
|
|---|
| 305 | if (TestOneVoxel( testVolume, *voxel, transform, inside, numInside, logger )) {
|
|---|
| 306 | if (++nError >= maxErrors) break;
|
|---|
| 307 | }
|
|---|
| 308 |
|
|---|
| 309 | randomRotate2.rotateZ( G4UniformRand() );
|
|---|
| 310 | transform.SetNetRotation( randomRotate2 );
|
|---|
| 311 |
|
|---|
| 312 | if (TestOneVoxel( testVolume, *voxel, transform, inside, numInside, logger )) {
|
|---|
| 313 | if (++nError >= maxErrors) break;
|
|---|
| 314 | }
|
|---|
| 315 | }
|
|---|
| 316 | delete voxel;
|
|---|
| 317 |
|
|---|
| 318 | if (nError >= maxErrors) {
|
|---|
| 319 | logger << "% End of test (maximum number errors) ";
|
|---|
| 320 | break;
|
|---|
| 321 | }
|
|---|
| 322 | if (++nVoxel >= maxVoxels) {
|
|---|
| 323 | logger << "% End of test (maximum number voxels) ";
|
|---|
| 324 | break;
|
|---|
| 325 | }
|
|---|
| 326 | }
|
|---|
| 327 |
|
|---|
| 328 | now = time(0);
|
|---|
| 329 | G4String dateTime2(ctime(&now));
|
|---|
| 330 | logger << dateTime2;
|
|---|
| 331 |
|
|---|
| 332 | logger << "% Statistics: voxels=" << nVoxel << " errors=" << nError << G4endl;
|
|---|
| 333 |
|
|---|
| 334 | logger << "%(End of file)" << G4endl;
|
|---|
| 335 | }
|
|---|
| 336 |
|
|---|
| 337 |
|
|---|
| 338 | //
|
|---|
| 339 | // TestOneVoxel
|
|---|
| 340 | //
|
|---|
| 341 | G4bool SBTvoxel::TestOneVoxel( const G4VSolid *testVolume,
|
|---|
| 342 | const G4VoxelLimits &voxel,
|
|---|
| 343 | const G4AffineTransform &transform,
|
|---|
| 344 | const G4ThreeVector inside[], const G4int numInside,
|
|---|
| 345 | std::ostream &logger ) const
|
|---|
| 346 | {
|
|---|
| 347 | static const EAxis axes[3] = { kXAxis, kYAxis, kZAxis };
|
|---|
| 348 | G4int numError = 0;
|
|---|
| 349 | G4double kCarTolerance = G4GeometryTolerance::GetInstance()->GetSurfaceTolerance();
|
|---|
| 350 |
|
|---|
| 351 | //
|
|---|
| 352 | // Get inverse transform
|
|---|
| 353 | //
|
|---|
| 354 | G4AffineTransform inverseTransform = transform.Inverse();
|
|---|
| 355 |
|
|---|
| 356 | //
|
|---|
| 357 | // Loop over the points, collecting min/max for each axis
|
|---|
| 358 | //
|
|---|
| 359 | G4double pointMins[3] = {+kInfinity, +kInfinity, +kInfinity},
|
|---|
| 360 | pointMaxs[3] = {-kInfinity, -kInfinity, -kInfinity};
|
|---|
| 361 | G4int numPointInside = 0;
|
|---|
| 362 |
|
|---|
| 363 | G4int i = numInside;
|
|---|
| 364 | while( i-- > 0 ) {
|
|---|
| 365 | G4ThreeVector point = transform.TransformPoint( inside[i] );
|
|---|
| 366 | if (voxel.Inside(point)) {
|
|---|
| 367 | numPointInside++;
|
|---|
| 368 |
|
|---|
| 369 | G4int j;
|
|---|
| 370 | for (j=0; j<3; j++) {
|
|---|
| 371 | G4double pv = point(axes[j]);
|
|---|
| 372 | if (pv < pointMins[j]) pointMins[j] = pv;
|
|---|
| 373 | if (pv > pointMaxs[j]) pointMaxs[j] = pv;
|
|---|
| 374 | }
|
|---|
| 375 | }
|
|---|
| 376 | }
|
|---|
| 377 |
|
|---|
| 378 | //
|
|---|
| 379 | // Loop over axes
|
|---|
| 380 | //
|
|---|
| 381 | for( i=0; i<3; i++ ) {
|
|---|
| 382 | G4double min, max;
|
|---|
| 383 |
|
|---|
| 384 | //
|
|---|
| 385 | // Query the solid
|
|---|
| 386 | //
|
|---|
| 387 | if (testVolume->CalculateExtent( axes[i], voxel, transform, min, max )) {
|
|---|
| 388 | //
|
|---|
| 389 | // Compare min/max to the list of inside points
|
|---|
| 390 | //
|
|---|
| 391 | if (min > pointMins[i] || max < pointMaxs[i]) {
|
|---|
| 392 | numError++;
|
|---|
| 393 | logger << "ERROR: Voxel limits are incorrect, axis " << i << G4endl;
|
|---|
| 394 | logger << " reported limits were: " << min << " " << max << G4endl;
|
|---|
| 395 | logger << " but points were found at: " << pointMins[i] << " " << pointMaxs[i] << G4endl;
|
|---|
| 396 | }
|
|---|
| 397 |
|
|---|
| 398 | //
|
|---|
| 399 | // Min or max in reversee order?
|
|---|
| 400 | //
|
|---|
| 401 | if (min >= max) {
|
|---|
| 402 | numError++;
|
|---|
| 403 | logger << "ERROR: Voxel limits max <= min, axis " << i << G4endl;
|
|---|
| 404 | logger << " reported limits were: " << min << " " << max << G4endl;
|
|---|
| 405 | }
|
|---|
| 406 |
|
|---|
| 407 | //
|
|---|
| 408 | // Min or max outside limits?
|
|---|
| 409 | //
|
|---|
| 410 | // We give the solid an extra "kCarTolerance" space, since
|
|---|
| 411 | // some solids like to add this value to their return values
|
|---|
| 412 | //
|
|---|
| 413 |
|
|---|
| 414 | // logger << std::setw(20) << std::setprecision(20);
|
|---|
| 415 |
|
|---|
| 416 | if ( voxel.IsLimited(axes[i]) ) {
|
|---|
| 417 | if (min < voxel.GetMinExtent(axes[i])-1.1*kCarTolerance) {
|
|---|
| 418 | numError++;
|
|---|
| 419 | logger << "ERROR: Voxel min is below pre-existing limit, axis " << i << G4endl;
|
|---|
| 420 | logger << " reported limits were: " << min << " " << max << G4endl;
|
|---|
| 421 | }
|
|---|
| 422 | if (max > voxel.GetMaxExtent(axes[i])+1.1*kCarTolerance) {
|
|---|
| 423 | numError++;
|
|---|
| 424 | logger << "ERROR: Voxel max is above pre-existing limit, axis " << i << G4endl;
|
|---|
| 425 | logger << " reported limits were: " << min << " " << max << G4endl;
|
|---|
| 426 | }
|
|---|
| 427 | }
|
|---|
| 428 |
|
|---|
| 429 |
|
|---|
| 430 | if ( (!voxel.IsLimited(axes[i])) || max < voxel.GetMaxExtent(axes[i]) ) {
|
|---|
| 431 | //
|
|---|
| 432 | // Construct a set of points just outside the voxel limits
|
|---|
| 433 | // and make sure they are outside
|
|---|
| 434 | //
|
|---|
| 435 | G4ThreeVector testPoints[9];
|
|---|
| 436 | MakeVoxelTestPoints( voxel, axes[i], max, testPoints );
|
|---|
| 437 |
|
|---|
| 438 | G4ThreeVector *testPoint = testPoints;
|
|---|
| 439 | do {
|
|---|
| 440 | G4ThreeVector tp = inverseTransform.TransformPoint( *testPoint );
|
|---|
| 441 | if (testVolume->Inside( tp ) == kInside) {
|
|---|
| 442 | numError++;
|
|---|
| 443 | logger << "ERROR: Voxel MAX limit is too small, axis " << i << G4endl;
|
|---|
| 444 | logger << " reported limits were: " << min << " " << max << G4endl;
|
|---|
| 445 | logger << " test point " << *testPoint << " [" << tp << "] was inside" << G4endl;
|
|---|
| 446 |
|
|---|
| 447 | MakeVoxelTestPoints( voxel, axes[i], max+10.0*kCarTolerance, testPoints );
|
|---|
| 448 | tp = inverseTransform.TransformPoint( *testPoint );
|
|---|
| 449 | if (testVolume->Inside( tp ) == kInside)
|
|---|
| 450 | logger << " and so was a more tolerant test point" << G4endl;
|
|---|
| 451 | break;
|
|---|
| 452 | }
|
|---|
| 453 | } while( ++testPoint < testPoints+9 );
|
|---|
| 454 | }
|
|---|
| 455 |
|
|---|
| 456 | if ( (!voxel.IsLimited(axes[i])) || min > voxel.GetMinExtent(axes[i]) ) {
|
|---|
| 457 | G4ThreeVector testPoints[9];
|
|---|
| 458 | MakeVoxelTestPoints( voxel, axes[i], min, testPoints );
|
|---|
| 459 |
|
|---|
| 460 | G4ThreeVector *testPoint = testPoints;
|
|---|
| 461 | do {
|
|---|
| 462 | G4ThreeVector tp = inverseTransform.TransformPoint( *testPoint );
|
|---|
| 463 | if (testVolume->Inside( tp ) == kInside) {
|
|---|
| 464 | numError++;
|
|---|
| 465 | logger << "ERROR: Voxel MIN limit is too large, axis " << i << G4endl;
|
|---|
| 466 | logger << " reported limits were: " << min << " " << max << G4endl;
|
|---|
| 467 | logger << " test point " << *testPoint << " [" << tp << "] was inside" << G4endl;
|
|---|
| 468 |
|
|---|
| 469 | MakeVoxelTestPoints( voxel, axes[i], min-10.0*kCarTolerance, testPoints );
|
|---|
| 470 | tp = inverseTransform.TransformPoint( *testPoint );
|
|---|
| 471 | if (testVolume->Inside( tp ) == kInside)
|
|---|
| 472 | logger << " and so was a more tolerant test point" << G4endl;
|
|---|
| 473 | break;
|
|---|
| 474 | }
|
|---|
| 475 | } while( ++testPoint < testPoints+9 );
|
|---|
| 476 | }
|
|---|
| 477 |
|
|---|
| 478 | }
|
|---|
| 479 | else {
|
|---|
| 480 | //
|
|---|
| 481 | // The voxel does not intersect the solid
|
|---|
| 482 | // Make sure there are no inside points inside
|
|---|
| 483 | //
|
|---|
| 484 | if (numPointInside) {
|
|---|
| 485 | numError++;
|
|---|
| 486 | logger << "ERROR: Voxel isn't empty, axis " << i
|
|---|
| 487 | << ", number points inside: " << numPointInside << " out of " << numInside << G4endl;
|
|---|
| 488 | logger << " they have limits of: " << pointMins[i] << " " << pointMaxs[i] << G4endl;
|
|---|
| 489 | }
|
|---|
| 490 | }
|
|---|
| 491 | }
|
|---|
| 492 |
|
|---|
| 493 | if (numError) {
|
|---|
| 494 | DumpVoxel( voxel, logger );
|
|---|
| 495 | DumpTransform( transform, logger );
|
|---|
| 496 | return true;
|
|---|
| 497 | }
|
|---|
| 498 |
|
|---|
| 499 | return false;
|
|---|
| 500 | }
|
|---|
| 501 |
|
|---|
| 502 |
|
|---|
| 503 | //
|
|---|
| 504 | // DumpVoxel
|
|---|
| 505 | //
|
|---|
| 506 | void SBTvoxel::DumpVoxel( const G4VoxelLimits &voxel, std::ostream &logger ) const
|
|---|
| 507 | {
|
|---|
| 508 | logger << "VOXEL =";
|
|---|
| 509 |
|
|---|
| 510 | static const EAxis axes[3] = { kXAxis, kYAxis, kZAxis };
|
|---|
| 511 |
|
|---|
| 512 | const EAxis *axis = axes;
|
|---|
| 513 | do {
|
|---|
| 514 | logger << " (";
|
|---|
| 515 | if (voxel.IsLimited(*axis))
|
|---|
| 516 | logger << voxel.GetMinExtent(*axis) << " "
|
|---|
| 517 | << voxel.GetMaxExtent(*axis);
|
|---|
| 518 | else
|
|---|
| 519 | logger << "unlimited";
|
|---|
| 520 | logger << ")";
|
|---|
| 521 | } while( ++axis < axes+3 );
|
|---|
| 522 |
|
|---|
| 523 | logger << G4endl;
|
|---|
| 524 | }
|
|---|
| 525 |
|
|---|
| 526 |
|
|---|
| 527 | //
|
|---|
| 528 | // DumpTransform
|
|---|
| 529 | //
|
|---|
| 530 | void SBTvoxel::DumpTransform( const G4AffineTransform &transform, std::ostream &logger ) const
|
|---|
| 531 | {
|
|---|
| 532 | G4RotationMatrix rotate = transform.NetRotation();
|
|---|
| 533 |
|
|---|
| 534 | G4ThreeVector axis;
|
|---|
| 535 | G4double amount;
|
|---|
| 536 |
|
|---|
| 537 | rotate.getAngleAxis( amount, axis );
|
|---|
| 538 |
|
|---|
| 539 | logger << "TRANLATE = " << transform.NetTranslation() << G4endl;
|
|---|
| 540 | logger << "ROTATE = " << axis << " " << amount << G4endl;
|
|---|
| 541 | }
|
|---|
| 542 |
|
|---|
| 543 |
|
|---|
| 544 |
|
|---|
| 545 | //
|
|---|
| 546 | // GetInsidePoints
|
|---|
| 547 | //
|
|---|
| 548 | void SBTvoxel::GetInsidePoints( const G4VSolid *testVolume,
|
|---|
| 549 | G4ThreeVector inside[], G4int *numInside,
|
|---|
| 550 | const G4int numPoints,
|
|---|
| 551 | const G4int maxAttempts ) const
|
|---|
| 552 | {
|
|---|
| 553 | *numInside = 0;
|
|---|
| 554 |
|
|---|
| 555 | G4int i;
|
|---|
| 556 | for( i=0; i<maxAttempts; i++ ) {
|
|---|
| 557 | G4ThreeVector point = GetRandomPoint();
|
|---|
| 558 |
|
|---|
| 559 | if (testVolume->Inside( point ) == kInside) {
|
|---|
| 560 | inside[*numInside] = point;
|
|---|
| 561 | if (++(*numInside) == numPoints) return;
|
|---|
| 562 | }
|
|---|
| 563 | }
|
|---|
| 564 | }
|
|---|
| 565 |
|
|---|
| 566 |
|
|---|
| 567 | //
|
|---|
| 568 | // GetRandomPoint
|
|---|
| 569 | //
|
|---|
| 570 | // Return a random point in three dimensions using the current
|
|---|
| 571 | // specs
|
|---|
| 572 | //
|
|---|
| 573 | G4ThreeVector SBTvoxel::GetRandomPoint() const {
|
|---|
| 574 | G4double dx = widths.x()*GaussianRandom(10*m/widths.x()),
|
|---|
| 575 | dy = widths.y()*GaussianRandom(10*m/widths.y()),
|
|---|
| 576 | dz = widths.z()*GaussianRandom(10*m/widths.z());
|
|---|
| 577 |
|
|---|
| 578 |
|
|---|
| 579 | G4ThreeVector randvec( dx, dy, dz );
|
|---|
| 580 |
|
|---|
| 581 | return target + randvec;
|
|---|
| 582 | }
|
|---|
| 583 |
|
|---|
| 584 |
|
|---|
| 585 | //
|
|---|
| 586 | // GaussianRandom
|
|---|
| 587 | //
|
|---|
| 588 | // Return Gaussian random number of unit width
|
|---|
| 589 | //
|
|---|
| 590 | // A classic, slow, but remarkably effective algorithm. Certainly good
|
|---|
| 591 | // enough for our purposes.
|
|---|
| 592 | //
|
|---|
| 593 | G4double SBTvoxel::GaussianRandom(const G4double cutoff) const {
|
|---|
| 594 | if (cutoff <= 0) G4Exception( "Illegal cutoff" );
|
|---|
| 595 |
|
|---|
| 596 | G4double answer;
|
|---|
| 597 | do {
|
|---|
| 598 | answer = -3.0;
|
|---|
| 599 | for( G4int j = 0; j < 6; j++ ) answer += G4UniformRand();
|
|---|
| 600 | answer *= std::sqrt(2.0);
|
|---|
| 601 | } while( std::fabs(answer) > cutoff );
|
|---|
| 602 |
|
|---|
| 603 | return(answer);
|
|---|
| 604 | }
|
|---|
| 605 |
|
|---|
| 606 |
|
|---|
| 607 | //
|
|---|
| 608 | // GetRandomLimit
|
|---|
| 609 | //
|
|---|
| 610 | G4bool SBTvoxel::GetRandomLimit( G4double x[2], const G4double range ) const
|
|---|
| 611 | {
|
|---|
| 612 | //
|
|---|
| 613 | // Generate a random number
|
|---|
| 614 | //
|
|---|
| 615 | G4double rand = G4UniformRand();
|
|---|
| 616 |
|
|---|
| 617 | //
|
|---|
| 618 | // Let's say that, well, 20% of the time we are
|
|---|
| 619 | // not limited in this dimension
|
|---|
| 620 | //
|
|---|
| 621 | if (rand < 0.2) return false;
|
|---|
| 622 |
|
|---|
| 623 | //
|
|---|
| 624 | // Otherwise, construct limits
|
|---|
| 625 | //
|
|---|
| 626 | x[0] = range*(rand - 0.6)/0.4;
|
|---|
| 627 | x[1] = x[0] + range*G4UniformRand();
|
|---|
| 628 | return true;
|
|---|
| 629 | }
|
|---|
| 630 |
|
|---|
| 631 |
|
|---|
| 632 | //
|
|---|
| 633 | // GetRandomVoxel
|
|---|
| 634 | //
|
|---|
| 635 | // Make a random voxel
|
|---|
| 636 | //
|
|---|
| 637 | G4VoxelLimits *SBTvoxel::NewRandomVoxel( const G4ThreeVector & // theWidths
|
|---|
| 638 | ) const
|
|---|
| 639 | {
|
|---|
| 640 | G4double xlim[2];
|
|---|
| 641 |
|
|---|
| 642 | G4VoxelLimits *voxel = new G4VoxelLimits;
|
|---|
| 643 |
|
|---|
| 644 | if (GetRandomLimit(xlim,widths.x())) voxel->AddLimit( kXAxis, xlim[0], xlim[1] );
|
|---|
| 645 | if (GetRandomLimit(xlim,widths.y())) voxel->AddLimit( kYAxis, xlim[0], xlim[1] );
|
|---|
| 646 | if (GetRandomLimit(xlim,widths.z())) voxel->AddLimit( kZAxis, xlim[0], xlim[1] );
|
|---|
| 647 |
|
|---|
| 648 | return voxel;
|
|---|
| 649 | }
|
|---|
| 650 |
|
|---|
| 651 |
|
|---|
| 652 | //
|
|---|
| 653 | // MakeVoxelTestPoints
|
|---|
| 654 | //
|
|---|
| 655 | // Make a set of nine vectors that lay in a grid inside the
|
|---|
| 656 | // voxel at the specified location along the specified axis
|
|---|
| 657 | //
|
|---|
| 658 | void SBTvoxel::MakeVoxelTestPoints( const G4VoxelLimits &voxel,
|
|---|
| 659 | const EAxis valueAxis, const G4double value,
|
|---|
| 660 | G4ThreeVector testPoints[9] ) const
|
|---|
| 661 | {
|
|---|
| 662 | G4ThreeVector grid[6], offset;
|
|---|
| 663 |
|
|---|
| 664 | static const EAxis axes[3] = { kXAxis, kYAxis, kZAxis };
|
|---|
| 665 | static const G4ThreeVector axisVectors[3] = {G4ThreeVector(1,0,0),
|
|---|
| 666 | G4ThreeVector(0,1,0),
|
|---|
| 667 | G4ThreeVector(0,0,1) };
|
|---|
| 668 |
|
|---|
| 669 | const EAxis *axis = axes;
|
|---|
| 670 | G4ThreeVector *nextGrid = grid;
|
|---|
| 671 | const G4ThreeVector *axisVector = axisVectors;
|
|---|
| 672 | do {
|
|---|
| 673 | if (*axis == valueAxis) {
|
|---|
| 674 | offset = value*(*axisVector);
|
|---|
| 675 | }
|
|---|
| 676 | else if (voxel.IsLimited(*axis)) {
|
|---|
| 677 | G4double min = voxel.GetMinExtent(*axis);
|
|---|
| 678 | G4double max = voxel.GetMaxExtent(*axis);
|
|---|
| 679 |
|
|---|
| 680 | if (min <= -kInfinity)
|
|---|
| 681 | min = max - 10.0*m;
|
|---|
| 682 | else if (max >= kInfinity)
|
|---|
| 683 | max = min + 10.0*m;
|
|---|
| 684 |
|
|---|
| 685 | (*nextGrid++) = 0.5*(min+max)*(*axisVector);
|
|---|
| 686 | (*nextGrid++) = min*(*axisVector);
|
|---|
| 687 | (*nextGrid++) = max*(*axisVector);
|
|---|
| 688 | }
|
|---|
| 689 | else {
|
|---|
| 690 | nextGrid++; // zero vector
|
|---|
| 691 | (*nextGrid++) = +10.0*m*(*axisVector);
|
|---|
| 692 | (*nextGrid++) = -10.0*m*(*axisVector);
|
|---|
| 693 | }
|
|---|
| 694 | } while( ++axisVector, ++axis < axes+3 );
|
|---|
| 695 |
|
|---|
| 696 | testPoints[0] = offset + grid[0] + grid[3];
|
|---|
| 697 | testPoints[1] = offset + grid[1] + grid[3];
|
|---|
| 698 | testPoints[2] = offset + grid[2] + grid[3];
|
|---|
| 699 | testPoints[3] = offset + grid[0] + grid[4];
|
|---|
| 700 | testPoints[4] = offset + grid[1] + grid[4];
|
|---|
| 701 | testPoints[5] = offset + grid[2] + grid[4];
|
|---|
| 702 | testPoints[6] = offset + grid[0] + grid[5];
|
|---|
| 703 | testPoints[7] = offset + grid[1] + grid[5];
|
|---|
| 704 | testPoints[8] = offset + grid[2] + grid[5];
|
|---|
| 705 | }
|
|---|