Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4Voxelizer.cc
Go to the documentation of this file.
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// G4Voxelizer implementation
27//
28// 19.10.12 Marek Gayer, created
29// --------------------------------------------------------------------
30
31#include <iostream>
32#include <iomanip>
33#include <sstream>
34#include <algorithm>
35#include <set>
36
37#include "G4VSolid.hh"
38
39#include "G4Orb.hh"
40#include "G4Voxelizer.hh"
41#include "G4SolidStore.hh"
42#include "Randomize.hh"
45#include "G4CSGSolid.hh"
46#include "G4Orb.hh"
47#include "G4Types.hh"
48#include "geomdefs.hh"
49
50using namespace std;
51
52G4ThreadLocal G4int G4Voxelizer::fDefaultVoxelsCount = -1;
53
54//______________________________________________________________________________
56 : fBoundingBox("VoxBBox", 1, 1, 1)
57{
58 fCountOfVoxels = fNPerSlice = fTotalCandidates = 0;
59
61
62 SetMaxVoxels(fDefaultVoxelsCount);
63
64 G4SolidStore::GetInstance()->DeRegister(&fBoundingBox);
65}
66
67//______________________________________________________________________________
69{
70}
71
72//______________________________________________________________________________
73void G4Voxelizer::BuildEmpty()
74{
75 // by reserving the size of candidates, we would avoid reallocation of
76 // the vector which could cause fragmentation
77 //
78 std::vector<G4int> xyz(3), max(3), candidates(fTotalCandidates);
79 const std::vector<G4int> empty(0);
80
81 for (auto i = 0; i <= 2; ++i) max[i] = (G4int)fBoundaries[i].size();
82 unsigned int size = max[0] * max[1] * max[2];
83
84 fEmpty.Clear();
85 fEmpty.ResetBitNumber(size-1);
86 fEmpty.ResetAllBits(true);
87
88 for (xyz[2] = 0; xyz[2] < max[2]; ++xyz[2])
89 {
90 for (xyz[1] = 0; xyz[1] < max[1]; ++xyz[1])
91 {
92 for (xyz[0] = 0; xyz[0] < max[0]; ++xyz[0])
93 {
94 if (GetCandidatesVoxelArray(xyz, candidates))
95 {
96 G4int index = GetVoxelsIndex(xyz);
97 fEmpty.SetBitNumber(index, false);
98
99 // rather than assigning directly with:
100 // "fCandidates[index] = candidates;", in an effort to ensure that
101 // capacity would be just exact, we rather use following 3 lines
102 //
103 std::vector<G4int> &c = (fCandidates[index] = empty);
104 c.reserve(candidates.size());
105 c.assign(candidates.begin(), candidates.end());
106 }
107 }
108 }
109 }
110#ifdef G4SPECSDEBUG
111 G4cout << "Non-empty voxels count: " << fCandidates.size() << G4endl;
112#endif
113}
114
115//______________________________________________________________________________
116void G4Voxelizer::BuildVoxelLimits(std::vector<G4VSolid*>& solids,
117 std::vector<G4Transform3D>& transforms)
118{
119 // "BuildVoxelLimits"'s aim is to store the coordinates of the origin as
120 // well as the half lengths related to the bounding box of each node.
121 // These quantities are stored in the array "fBoxes" (6 different values per
122 // node
123 //
124 if (std::size_t numNodes = solids.size()) // Number of nodes in "multiUnion"
125 {
126 fBoxes.resize(numNodes); // Array which will store the half lengths
127 fNPerSlice = G4int(1 + (fBoxes.size() - 1) / (8 * sizeof(unsigned int)));
128
129 // related to a particular node, but also
130 // the coordinates of its origin
131
132 G4ThreeVector toleranceVector(fTolerance,fTolerance,fTolerance);
133
134 for (std::size_t i = 0; i < numNodes; ++i)
135 {
136 G4VSolid& solid = *solids[i];
137 G4Transform3D transform = transforms[i];
139
140 solid.BoundingLimits(min, max);
141 if (solid.GetEntityType() == "G4Orb")
142 {
143 G4Orb& orb = *(G4Orb*) &solid;
144 G4ThreeVector orbToleranceVector;
145 G4double tolerance = orb.GetRadialTolerance() / 2.0;
146 orbToleranceVector.set(tolerance,tolerance,tolerance);
147 min -= orbToleranceVector;
148 max += orbToleranceVector;
149 }
150 else
151 {
152 min -= toleranceVector;
153 max += toleranceVector;
154 }
155 TransformLimits(min, max, transform);
156 fBoxes[i].hlen = (max - min) / 2.;
157 fBoxes[i].pos = (max + min) / 2.;
158 }
159 fTotalCandidates = (G4int)fBoxes.size();
160 }
161}
162
163//______________________________________________________________________________
164void G4Voxelizer::BuildVoxelLimits(std::vector<G4VFacet*>& facets)
165{
166 // "BuildVoxelLimits"'s aim is to store the coordinates of the origin as well
167 // as the half lengths related to the bounding box of each node.
168 // These quantities are stored in the array "fBoxes" (6 different values per
169 // node.
170
171 if (std::size_t numNodes = facets.size()) // Number of nodes
172 {
173 fBoxes.resize(numNodes); // Array which will store the half lengths
174 fNPerSlice = G4int(1+(fBoxes.size()-1)/(8*sizeof(unsigned int)));
175
176 G4ThreeVector toleranceVector(10*fTolerance, 10*fTolerance, 10*fTolerance);
177
178 for (std::size_t i = 0; i < numNodes; ++i)
179 {
180 G4VFacet &facet = *facets[i];
182 G4ThreeVector x(1,0,0), y(0,1,0), z(0,0,1);
183 G4ThreeVector extent;
184 max.set (facet.Extent(x), facet.Extent(y), facet.Extent(z));
185 min.set (-facet.Extent(-x), -facet.Extent(-y), -facet.Extent(-z));
186 min -= toleranceVector;
187 max += toleranceVector;
188 G4ThreeVector hlen = (max - min) / 2;
189 fBoxes[i].hlen = hlen;
190 fBoxes[i].pos = min + hlen;
191 }
192 fTotalCandidates = (G4int)fBoxes.size();
193 }
194}
195
196//______________________________________________________________________________
198{
199 // "DisplayVoxelLimits" displays the dX, dY, dZ, pX, pY and pZ for each node
200
201 std::size_t numNodes = fBoxes.size();
202 G4long oldprec = G4cout.precision(16);
203 for(std::size_t i = 0; i < numNodes; ++i)
204 {
205 G4cout << setw(10) << setiosflags(ios::fixed) <<
206 " -> Node " << i+1 << ":\n" <<
207 "\t * [x,y,z] = " << fBoxes[i].hlen <<
208 "\t * [x,y,z] = " << fBoxes[i].pos << "\n";
209 }
210 G4cout.precision(oldprec);
211}
212
213//______________________________________________________________________________
214void G4Voxelizer::CreateSortedBoundary(std::vector<G4double>& boundary,
215 G4int axis)
216{
217 // "CreateBoundaries"'s aim is to determine the slices induced by the
218 // bounding fBoxes, along each axis. The created boundaries are stored
219 // in the array "boundariesRaw"
220
221 std::size_t numNodes = fBoxes.size(); // Number of nodes in structure
222
223 // Determination of the boundaries along x, y and z axis
224 //
225 for(std::size_t i = 0 ; i < numNodes; ++i)
226 {
227 // For each node, the boundaries are created by using the array "fBoxes"
228 // built in method "BuildVoxelLimits"
229 //
230 G4double p = fBoxes[i].pos[axis], d = fBoxes[i].hlen[axis];
231
232 // x boundaries
233 //
234#ifdef G4SPECSDEBUG
235 G4cout << "Boundary " << p - d << " - " << p + d << G4endl;
236#endif
237 boundary[2*i] = p - d;
238 boundary[2*i+1] = p + d;
239 }
240 std::sort(boundary.begin(), boundary.end());
241}
242
243//______________________________________________________________________________
244void G4Voxelizer::BuildBoundaries()
245{
246 // "SortBoundaries" orders the boundaries along each axis (increasing order)
247 // and also does not take into account redundant boundaries, i.e. if two
248 // boundaries are separated by a distance strictly inferior to "tolerance".
249 // The sorted boundaries are respectively stored in:
250 // * boundaries[0..2]
251
252 // In addition, the number of elements contained in the three latter arrays
253 // are precise thanks to variables: boundariesCountX, boundariesCountY and
254 // boundariesCountZ.
255
256 if (std::size_t numNodes = fBoxes.size())
257 {
258 const G4double tolerance = fTolerance / 100.0;
259 // Minimal distance to discriminate two boundaries.
260
261 std::vector<G4double> sortedBoundary(2*numNodes);
262
263 for (auto j = 0; j <= 2; ++j)
264 {
265 CreateSortedBoundary(sortedBoundary, j);
266 std::vector<G4double> &boundary = fBoundaries[j];
267 boundary.clear();
268
269 for(std::size_t i = 0 ; i < 2*numNodes; ++i)
270 {
271 G4double newBoundary = sortedBoundary[i];
272#ifdef G4SPECSDEBUG
273 if (j == 0) G4cout << "Examining " << newBoundary << "..." << G4endl;
274#endif
275 G4int size = (G4int)boundary.size();
276 if(!size || std::abs(boundary[size-1] - newBoundary) > tolerance)
277 {
278 {
279#ifdef G4SPECSDEBUG
280 if (j == 0) G4cout << "Adding boundary " << newBoundary << "..."
281 << G4endl;
282#endif
283 boundary.push_back(newBoundary);
284 continue;
285 }
286 }
287 // If two successive boundaries are too close from each other,
288 // only the first one is considered
289 }
290
291 G4int n = (G4int)boundary.size();
292 G4int max = 100000;
293 if (n > max/2)
294 {
295 G4int skip = n / (max /2); // n has to be 2x bigger then 50.000.
296 // therefore only from 100.000 reduced
297 std::vector<G4double> reduced;
298 for (G4int i = 0; i < n; ++i)
299 {
300 // 50 ok for 2k, 1000, 2000
301 G4int size = (G4int)boundary.size();
302 if (i % skip == 0 || i == 0 || i == size - 1)
303 {
304 // this condition of merging boundaries was wrong,
305 // it did not count with right part, which can be
306 // completely ommited and not included in final consideration.
307 // Now should be OK
308 //
309 reduced.push_back(boundary[i]);
310 }
311 }
312 boundary = reduced;
313 }
314 }
315 }
316}
317
318//______________________________________________________________________________
320{
321 char axis[3] = {'X', 'Y', 'Z'};
322 for (auto i = 0; i <= 2; ++i)
323 {
324 G4cout << " * " << axis[i] << " axis:" << G4endl << " | ";
325 DisplayBoundaries(fBoundaries[i]);
326 }
327}
328
329//______________________________________________________________________________
330void G4Voxelizer::DisplayBoundaries(std::vector<G4double> &boundaries)
331{
332 // Prints the positions of the boundaries of the slices on the three axes
333
334 std::size_t count = boundaries.size();
335 G4long oldprec = G4cout.precision(16);
336 for(std::size_t i = 0; i < count; ++i)
337 {
338 G4cout << setw(10) << setiosflags(ios::fixed) << boundaries[i];
339 if(i != count-1) G4cout << "-> ";
340 }
341 G4cout << "|" << G4endl << "Number of boundaries: " << count << G4endl;
342 G4cout.precision(oldprec);
343}
344
345//______________________________________________________________________________
346void G4Voxelizer::BuildBitmasks(std::vector<G4double> boundaries[],
347 G4SurfBits bitmasks[], G4bool countsOnly)
348{
349 // "BuildListNodes" stores in the bitmasks solids present in each slice
350 // along an axis.
351
352 std::size_t numNodes = fBoxes.size();
353 G4int bitsPerSlice = GetBitsPerSlice();
354
355 for (auto k = 0; k < 3; ++k)
356 {
357 std::vector<G4double>& boundary = boundaries[k];
358 G4int voxelsCount = (G4int)boundary.size() - 1;
359 G4SurfBits& bitmask = bitmasks[k];
360
361 if (!countsOnly)
362 {
363 bitmask.Clear();
364#ifdef G4SPECSDEBUG
365 G4cout << "Allocating bitmask..." << G4endl;
366#endif
367 bitmask.SetBitNumber(voxelsCount*bitsPerSlice-1, false);
368 // it is here so we can set the maximum number of bits. this line
369 // will rellocate the memory and set all to zero
370 }
371 std::vector<G4int>& candidatesCount = fCandidatesCounts[k];
372 candidatesCount.resize(voxelsCount);
373
374 for(G4int i = 0 ; i < voxelsCount; ++i) { candidatesCount[i] = 0; }
375
376 // Loop on the nodes, number of slices per axis
377 //
378 for(std::size_t j = 0 ; j < numNodes; ++j)
379 {
380 // Determination of the minimum and maximum position along x
381 // of the bounding boxe of each node
382 //
383 G4double p = fBoxes[j].pos[k], d = fBoxes[j].hlen[k];
384
385 G4double min = p - d; // - localTolerance;
386 G4double max = p + d; // + localTolerance;
387
388 G4int i = BinarySearch(boundary, min);
389 if (i < 0) { i = 0; }
390
391 do // Loop checking, 13.08.2015, G.Cosmo
392 {
393 if (!countsOnly)
394 {
395 bitmask.SetBitNumber(i*bitsPerSlice+(G4int)j);
396 }
397 candidatesCount[i]++;
398 ++i;
399 }
400 while (max > boundary[i] && i < voxelsCount);
401 }
402 }
403#ifdef G4SPECSDEBUG
404 G4cout << "Build list nodes completed." << G4endl;
405#endif
406}
407
408//______________________________________________________________________________
409G4String G4Voxelizer::GetCandidatesAsString(const G4SurfBits& bits) const
410{
411 // Decodes the candidates in mask as G4String.
412
413 stringstream ss;
414 G4int numNodes = (G4int)fBoxes.size();
415
416 for(auto i=0; i<numNodes; ++i)
417 {
418 if (bits.TestBitNumber(i)) { ss << i+1 << " "; }
419 }
420 return ss.str();
421}
422
423//______________________________________________________________________________
425{
426 // Prints which solids are present in the slices previously elaborated.
427
428 char axis[3] = {'X', 'Y', 'Z'};
429 G4int size=8*sizeof(G4int)*fNPerSlice;
430 G4SurfBits bits(size);
431
432 for (auto j = 0; j <= 2; ++j)
433 {
434 G4cout << " * " << axis[j] << " axis:" << G4endl;
435 G4int count = (G4int)fBoundaries[j].size();
436 for(G4int i=0; i < count-1; ++i)
437 {
438 G4cout << " Slice #" << i+1 << ": [" << fBoundaries[j][i]
439 << " ; " << fBoundaries[j][i+1] << "] -> ";
440 bits.set(size,(const char *)fBitmasks[j].fAllBits+i
441 *fNPerSlice*sizeof(G4int));
442 G4String result = GetCandidatesAsString(bits);
443 G4cout << "[ " << result.c_str() << "] " << G4endl;
444 }
445 }
446}
447
448//______________________________________________________________________________
449void G4Voxelizer::BuildBoundingBox()
450{
451 G4ThreeVector min(fBoundaries[0].front(),
452 fBoundaries[1].front(),
453 fBoundaries[2].front());
454 G4ThreeVector max(fBoundaries[0].back(),
455 fBoundaries[1].back(),
456 fBoundaries[2].back());
457 BuildBoundingBox(min, max);
458}
459
460//______________________________________________________________________________
461void G4Voxelizer::BuildBoundingBox(G4ThreeVector& amin,
462 G4ThreeVector& amax,
463 G4double tolerance)
464{
465 for (auto i = 0; i <= 2; ++i)
466 {
467 G4double min = amin[i];
468 G4double max = amax[i];
469 fBoundingBoxSize[i] = (max - min) / 2 + tolerance * 0.5;
470 fBoundingBoxCenter[i] = min + fBoundingBoxSize[i];
471 }
472 fBoundingBox.SetXHalfLength(fBoundingBoxSize.x());
473 fBoundingBox.SetYHalfLength(fBoundingBoxSize.y());
474 fBoundingBox.SetZHalfLength(fBoundingBoxSize.z());
475}
476
477// algorithm -
478
479// in order to get balanced voxels, merge should always unite those regions,
480// where the number of voxels is least the number.
481// We will keep sorted list (std::set) with all voxels. there will be
482// comparator function between two voxels, which will tell if voxel is less
483// by looking at his right neighbor.
484// First, we will add all the voxels into the tree.
485// We will be pick the first item in the tree, merging it, adding the right
486// merged voxel into the a list for future reduction (fBitmasks will be
487// rebuilded later, therefore they need not to be updated).
488// The merged voxel need to be added to the tree again, so it's position
489// would be updated.
490
491//______________________________________________________________________________
492void G4Voxelizer::SetReductionRatio(G4int maxVoxels,
493 G4ThreeVector& reductionRatio)
494{
495 G4double maxTotal = (G4double) fCandidatesCounts[0].size()
496 * fCandidatesCounts[1].size() * fCandidatesCounts[2].size();
497
498 if (maxVoxels > 0 && maxVoxels < maxTotal)
499 {
500 G4double ratio = (G4double) maxVoxels / maxTotal;
501 ratio = std::pow(ratio, 1./3.);
502 if (ratio > 1) { ratio = 1; }
503 reductionRatio.set(ratio,ratio,ratio);
504 }
505}
506
507//______________________________________________________________________________
508void G4Voxelizer::BuildReduceVoxels(std::vector<G4double> boundaries[],
509 G4ThreeVector reductionRatio)
510{
511 for (auto k = 0; k <= 2; ++k)
512 {
513 std::vector<G4int> &candidatesCount = fCandidatesCounts[k];
514 G4int max = (G4int)candidatesCount.size();
515 std::vector<G4VoxelInfo> voxels(max);
516 G4VoxelComparator comp(voxels);
517 std::set<G4int, G4VoxelComparator> voxelSet(comp);
518 std::vector<G4int> mergings;
519
520 for (G4int j = 0; j < max; ++j)
521 {
522 G4VoxelInfo &voxel = voxels[j];
523 voxel.count = candidatesCount[j];
524 voxel.previous = j - 1;
525 voxel.next = j + 1;
526 voxels[j] = voxel;
527 }
528
529 for (G4int j = 0; j < max - 1; ++j) { voxelSet.insert(j); }
530 // we go to size-1 to make sure we will not merge the last element
531
532 G4double reduction = reductionRatio[k];
533 if (reduction != 0)
534 {
535 G4int count = 0, currentCount;
536 while ((currentCount = (G4int)voxelSet.size()) > 2)
537 {
538 G4double currentRatio = 1 - (G4double) count / max;
539 if ((currentRatio <= reduction) && (currentCount <= 1000))
540 break;
541 const G4int pos = *voxelSet.begin();
542 mergings.push_back(pos + 1);
543
544 G4VoxelInfo& voxel = voxels[pos];
545 G4VoxelInfo& nextVoxel = voxels[voxel.next];
546
547 if (voxelSet.erase(pos) != 1)
548 {
549 ;// k = k;
550 }
551 if (voxel.next != max - 1)
552 if (voxelSet.erase(voxel.next) != 1)
553 {
554 ;// k = k;
555 }
556 if (voxel.previous != -1)
557 if (voxelSet.erase(voxel.previous) != 1)
558 {
559 ;// k = k;
560 }
561 nextVoxel.count += voxel.count;
562 voxel.count = 0;
563 nextVoxel.previous = voxel.previous;
564
565 if (voxel.next != max - 1)
566 voxelSet.insert(voxel.next);
567
568 if (voxel.previous != -1)
569 {
570 voxels[voxel.previous].next = voxel.next;
571 voxelSet.insert(voxel.previous);
572 }
573 ++count;
574 } // Loop checking, 13.08.2015, G.Cosmo
575 }
576
577 if (mergings.size())
578 {
579 std::sort(mergings.begin(), mergings.end());
580
581 const std::vector<G4double>& boundary = boundaries[k];
582 G4int mergingsSize = (G4int)mergings.size();
583 vector<G4double> reducedBoundary;
584 G4int skip = mergings[0], i = 0;
585 max = (G4int)boundary.size();
586 for (G4int j = 0; j < max; ++j)
587 {
588 if (j != skip)
589 {
590 reducedBoundary.push_back(boundary[j]);
591 }
592 else if (++i < mergingsSize)
593 {
594 skip = mergings[i];
595 }
596 }
597 boundaries[k] = reducedBoundary;
598 }
599/*
600 G4int count = 0;
601 while (true) // Loop checking, 13.08.2015, G.Cosmo
602 {
603 G4double reduction = reductionRatio[k];
604 if (reduction == 0)
605 break;
606 G4int currentCount = voxelSet.size();
607 if (currentCount <= 2)
608 break;
609 G4double currentRatio = 1 - (G4double) count / max;
610 if (currentRatio <= reduction && currentCount <= 1000)
611 break;
612 const G4int pos = *voxelSet.begin();
613 mergings.push_back(pos);
614
615 G4VoxelInfo &voxel = voxels[pos];
616 G4VoxelInfo &nextVoxel = voxels[voxel.next];
617
618 voxelSet.erase(pos);
619 if (voxel.next != max - 1) { voxelSet.erase(voxel.next); }
620 if (voxel.previous != -1) { voxelSet.erase(voxel.previous); }
621
622 nextVoxel.count += voxel.count;
623 voxel.count = 0;
624 nextVoxel.previous = voxel.previous;
625
626 if (voxel.next != max - 1)
627 voxelSet.insert(voxel.next);
628
629 if (voxel.previous != -1)
630 {
631 voxels[voxel.previous].next = voxel.next;
632 voxelSet.insert(voxel.previous);
633 }
634 ++count;
635 }
636
637 if (mergings.size())
638 {
639 std::sort(mergings.begin(), mergings.end());
640
641 std::vector<G4double> &boundary = boundaries[k];
642 std::vector<G4double> reducedBoundary(boundary.size() - mergings.size());
643 G4int skip = mergings[0] + 1, cur = 0, i = 0;
644 max = boundary.size();
645 for (G4int j = 0; j < max; ++j)
646 {
647 if (j != skip)
648 {
649 reducedBoundary[cur++] = boundary[j];
650 }
651 else
652 {
653 if (++i < (G4int)mergings.size()) { skip = mergings[i] + 1; }
654 }
655 }
656 boundaries[k] = reducedBoundary;
657 }
658*/
659 }
660}
661
662//______________________________________________________________________________
663void G4Voxelizer::BuildReduceVoxels2(std::vector<G4double> boundaries[],
664 G4ThreeVector reductionRatio)
665{
666 for (auto k = 0; k <= 2; ++k)
667 {
668 std::vector<G4int> &candidatesCount = fCandidatesCounts[k];
669 G4int max = (G4int)candidatesCount.size();
670 G4int total = 0;
671 for (G4int i = 0; i < max; ++i) total += candidatesCount[i];
672
673 G4double reduction = reductionRatio[k];
674 if (reduction == 0)
675 break;
676
677 G4int destination = (G4int) (reduction * max) + 1;
678 if (destination > 1000) destination = 1000;
679 if (destination < 2) destination = 2;
680 G4double average = ((G4double)total / max) / reduction;
681
682 std::vector<G4int> mergings;
683
684 std::vector<G4double> &boundary = boundaries[k];
685 std::vector<G4double> reducedBoundary(destination);
686
687 G4int sum = 0, cur = 0;
688 for (G4int i = 0; i < max; ++i)
689 {
690 sum += candidatesCount[i];
691 if (sum > average * (cur + 1) || i == 0)
692 {
693 G4double val = boundary[i];
694 reducedBoundary[cur] = val;
695 ++cur;
696 if (cur == destination)
697 break;
698 }
699 }
700 reducedBoundary[destination-1] = boundary[max];
701 boundaries[k] = reducedBoundary;
702 }
703}
704
705//______________________________________________________________________________
706void G4Voxelizer::Voxelize(std::vector<G4VSolid*>& solids,
707 std::vector<G4Transform3D>& transforms)
708{
709 BuildVoxelLimits(solids, transforms);
710 BuildBoundaries();
711 BuildBitmasks(fBoundaries, fBitmasks);
712 BuildBoundingBox();
713 BuildEmpty(); // this does not work well for multi-union,
714 // actually only makes performance slower,
715 // these are only pre-calculated but not used by multi-union
716
717 for (auto i = 0; i < 3; ++i)
718 {
719 fCandidatesCounts[i].resize(0);
720 }
721}
722
723//______________________________________________________________________________
724void G4Voxelizer::CreateMiniVoxels(std::vector<G4double> boundaries[],
725 G4SurfBits bitmasks[])
726{
727 std::vector<G4int> voxel(3), maxVoxels(3);
728 for (auto i = 0; i <= 2; ++i) maxVoxels[i] = (G4int)boundaries[i].size();
729
730 G4ThreeVector point;
731 for (voxel[2] = 0; voxel[2] < maxVoxels[2] - 1; ++voxel[2])
732 {
733 for (voxel[1] = 0; voxel[1] < maxVoxels[1] - 1; ++voxel[1])
734 {
735 for (voxel[0] = 0; voxel[0] < maxVoxels[0] - 1; ++voxel[0])
736 {
737 std::vector<G4int> candidates;
738 if (GetCandidatesVoxelArray(voxel, bitmasks, candidates, 0))
739 {
740 // find a box for corresponding non-empty voxel
741 G4VoxelBox box;
742 for (auto i = 0; i <= 2; ++i)
743 {
744 G4int index = voxel[i];
745 const std::vector<G4double> &boundary = boundaries[i];
746 G4double hlen = 0.5 * (boundary[index+1] - boundary[index]);
747 box.hlen[i] = hlen;
748 box.pos[i] = boundary[index] + hlen;
749 }
750 fVoxelBoxes.push_back(box);
751 std::vector<G4int>(candidates).swap(candidates);
752 fVoxelBoxesCandidates.push_back(candidates);
753 }
754 }
755 }
756 }
757}
758
759//______________________________________________________________________________
760void G4Voxelizer::Voxelize(std::vector<G4VFacet*>& facets)
761{
762 G4int maxVoxels = fMaxVoxels;
763 G4ThreeVector reductionRatio = fReductionRatio;
764
765 std::size_t size = facets.size();
766 if (size < 10)
767 {
768 for (std::size_t i = 0; i < facets.size(); ++i)
769 {
770 if (facets[i]->GetNumberOfVertices() > 3) ++size;
771 }
772 }
773
774 if ((size >= 10 || maxVoxels > 0) && maxVoxels != 0 && maxVoxels != 1)
775 {
776#ifdef G4SPECSDEBUG
777 G4cout << "Building voxel limits..." << G4endl;
778#endif
779
780 BuildVoxelLimits(facets);
781
782#ifdef G4SPECSDEBUG
783 G4cout << "Building boundaries..." << G4endl;
784#endif
785
786 BuildBoundaries();
787
788#ifdef G4SPECSDEBUG
789 G4cout << "Building bitmasks..." << G4endl;
790#endif
791
792 BuildBitmasks(fBoundaries, 0, true);
793
794 if (maxVoxels < 0 && reductionRatio == G4ThreeVector())
795 {
796 maxVoxels = fTotalCandidates;
797 if (fTotalCandidates > 1000000) maxVoxels = 1000000;
798 }
799
800 SetReductionRatio(maxVoxels, reductionRatio);
801
802 fCountOfVoxels = CountVoxels(fBoundaries);
803
804#ifdef G4SPECSDEBUG
805 G4cout << "Total number of voxels: " << fCountOfVoxels << G4endl;
806#endif
807
808 BuildReduceVoxels2(fBoundaries, reductionRatio);
809
810 fCountOfVoxels = CountVoxels(fBoundaries);
811
812#ifdef G4SPECSDEBUG
813 G4cout << "Total number of voxels after reduction: "
814 << fCountOfVoxels << G4endl;
815#endif
816
817#ifdef G4SPECSDEBUG
818 G4cout << "Building bitmasks..." << G4endl;
819#endif
820
821 BuildBitmasks(fBoundaries, fBitmasks);
822
823 G4ThreeVector reductionRatioMini;
824
825 G4SurfBits bitmasksMini[3];
826
827 // section for building mini voxels
828
829 std::vector<G4double> miniBoundaries[3];
830
831 for (auto i = 0; i <= 2; ++i) { miniBoundaries[i] = fBoundaries[i]; }
832
833 G4int voxelsCountMini = (fCountOfVoxels >= 1000)
834 ? 100 : G4int(fCountOfVoxels / 10);
835
836 SetReductionRatio(voxelsCountMini, reductionRatioMini);
837
838#ifdef G4SPECSDEBUG
839 G4cout << "Building reduced voxels..." << G4endl;
840#endif
841
842 BuildReduceVoxels(miniBoundaries, reductionRatioMini);
843
844#ifdef G4SPECSDEBUG
845 G4int total = CountVoxels(miniBoundaries);
846 G4cout << "Total number of mini voxels: " << total << G4endl;
847#endif
848
849#ifdef G4SPECSDEBUG
850 G4cout << "Building mini bitmasks..." << G4endl;
851#endif
852
853 BuildBitmasks(miniBoundaries, bitmasksMini);
854
855#ifdef G4SPECSDEBUG
856 G4cout << "Creating Mini Voxels..." << G4endl;
857#endif
858
859 CreateMiniVoxels(miniBoundaries, bitmasksMini);
860
861#ifdef G4SPECSDEBUG
862 G4cout << "Building bounding box..." << G4endl;
863#endif
864
865 BuildBoundingBox();
866
867#ifdef G4SPECSDEBUG
868 G4cout << "Building empty..." << G4endl;
869#endif
870
871 BuildEmpty();
872
873#ifdef G4SPECSDEBUG
874 G4cout << "Deallocating unnecessary fields during runtime..." << G4endl;
875#endif
876 // deallocate fields unnecessary during runtime
877 //
878 fBoxes.resize(0);
879 for (auto i = 0; i < 3; ++i)
880 {
881 fCandidatesCounts[i].resize(0);
882 fBitmasks[i].Clear();
883 }
884 }
885}
886
887
888//______________________________________________________________________________
889void G4Voxelizer::GetCandidatesVoxel(std::vector<G4int>& voxels)
890{
891 // "GetCandidates" should compute which solids are possibly contained in
892 // the voxel defined by the three slices characterized by the passed indexes.
893
894 G4cout << " Candidates in voxel [" << voxels[0] << " ; " << voxels[1]
895 << " ; " << voxels[2] << "]: ";
896 std::vector<G4int> candidates;
897 G4int count = GetCandidatesVoxelArray(voxels, candidates);
898 G4cout << "[ ";
899 for (G4int i = 0; i < count; ++i) G4cout << candidates[i];
900 G4cout << "] " << G4endl;
901}
902
903//______________________________________________________________________________
904void G4Voxelizer::FindComponentsFastest(unsigned int mask,
905 std::vector<G4int>& list, G4int i)
906{
907 for (G4int byte = 0; byte < (G4int) (sizeof(unsigned int)); ++byte)
908 {
909 if (G4int maskByte = mask & 0xFF)
910 {
911 for (G4int bit = 0; bit < 8; ++bit)
912 {
913 if (maskByte & 1)
914 { list.push_back(8*(sizeof(unsigned int)*i+ byte) + bit); }
915 if (!(maskByte >>= 1)) break;
916 }
917 }
918 mask >>= 8;
919 }
920}
921
922//______________________________________________________________________________
923void G4Voxelizer::TransformLimits(G4ThreeVector& min, G4ThreeVector& max,
924 const G4Transform3D& transformation) const
925{
926 // The goal of this method is to convert the quantities min and max
927 // (representing the bounding box of a given solid in its local frame)
928 // to the main frame, using "transformation"
929
930 G4ThreeVector vertices[8] = // Detemination of the vertices thanks to
931 { // the extension of each solid:
932 G4ThreeVector(min.x(), min.y(), min.z()), // 1st vertice:
933 G4ThreeVector(min.x(), max.y(), min.z()), // 2nd vertice:
934 G4ThreeVector(max.x(), max.y(), min.z()),
935 G4ThreeVector(max.x(), min.y(), min.z()),
936 G4ThreeVector(min.x(), min.y(), max.z()),
937 G4ThreeVector(min.x(), max.y(), max.z()),
938 G4ThreeVector(max.x(), max.y(), max.z()),
939 G4ThreeVector(max.x(), min.y(), max.z())
940 };
941
942 min.set(kInfinity,kInfinity,kInfinity);
943 max.set(-kInfinity,-kInfinity,-kInfinity);
944
945 // Loop on th vertices
946 G4int limit = sizeof(vertices) / sizeof(G4ThreeVector);
947 for (G4int i = 0 ; i < limit; ++i)
948 {
949 // From local frame to the global one:
950 // Current positions on the three axis:
951 G4ThreeVector current = GetGlobalPoint(transformation, vertices[i]);
952
953 // If need be, replacement of the min & max values:
954 if (current.x() > max.x()) max.setX(current.x());
955 if (current.x() < min.x()) min.setX(current.x());
956
957 if (current.y() > max.y()) max.setY(current.y());
958 if (current.y() < min.y()) min.setY(current.y());
959
960 if (current.z() > max.z()) max.setZ(current.z());
961 if (current.z() < min.z()) min.setZ(current.z());
962 }
963}
964
965//______________________________________________________________________________
967 std::vector<G4int> &list, G4SurfBits *crossed) const
968{
969 // Method returning the candidates corresponding to the passed point
970
971 list.clear();
972
973 for (auto i = 0; i <= 2; ++i)
974 {
975 if(point[i] < fBoundaries[i].front() || point[i] >= fBoundaries[i].back())
976 return 0;
977 }
978
979 if (fTotalCandidates == 1)
980 {
981 list.push_back(0);
982 return 1;
983 }
984 else
985 {
986 if (fNPerSlice == 1)
987 {
988 unsigned int mask = 0xFFffFFff;
989 G4int slice;
990 if (fBoundaries[0].size() > 2)
991 {
992 slice = BinarySearch(fBoundaries[0], point.x());
993 if (!(mask = ((unsigned int*) fBitmasks[0].fAllBits)[slice]))
994 return 0;
995 }
996 if (fBoundaries[1].size() > 2)
997 {
998 slice = BinarySearch(fBoundaries[1], point.y());
999 if (!(mask &= ((unsigned int*) fBitmasks[1].fAllBits)[slice]))
1000 return 0;
1001 }
1002 if (fBoundaries[2].size() > 2)
1003 {
1004 slice = BinarySearch(fBoundaries[2], point.z());
1005 if (!(mask &= ((unsigned int*) fBitmasks[2].fAllBits)[slice]))
1006 return 0;
1007 }
1008 if (crossed && (!(mask &= ~((unsigned int*)crossed->fAllBits)[0])))
1009 return 0;
1010
1011 FindComponentsFastest(mask, list, 0);
1012 }
1013 else
1014 {
1015 unsigned int* masks[3], mask; // masks for X,Y,Z axis
1016 for (auto i = 0; i <= 2; ++i)
1017 {
1018 G4int slice = BinarySearch(fBoundaries[i], point[i]);
1019 masks[i] = ((unsigned int*) fBitmasks[i].fAllBits)
1020 + slice * fNPerSlice;
1021 }
1022 unsigned int* maskCrossed = crossed
1023 ? (unsigned int*)crossed->fAllBits : 0;
1024
1025 for (G4int i = 0 ; i < fNPerSlice; ++i)
1026 {
1027 // Logic "and" of the masks along the 3 axes x, y, z:
1028 // removing "if (!" and ") continue" => slightly slower
1029 //
1030 if (!(mask = masks[0][i])) continue;
1031 if (!(mask &= masks[1][i])) continue;
1032 if (!(mask &= masks[2][i])) continue;
1033 if (maskCrossed && !(mask &= ~maskCrossed[i])) continue;
1034
1035 FindComponentsFastest(mask, list, i);
1036 }
1037 }
1038/*
1039 if (fNPerSlice == 1)
1040 {
1041 unsigned int mask;
1042 G4int slice = BinarySearch(fBoundaries[0], point.x());
1043 if (!(mask = ((unsigned int *) fBitmasks[0].fAllBits)[slice]
1044 )) return 0;
1045 slice = BinarySearch(fBoundaries[1], point.y());
1046 if (!(mask &= ((unsigned int *) fBitmasks[1].fAllBits)[slice]
1047 )) return 0;
1048 slice = BinarySearch(fBoundaries[2], point.z());
1049 if (!(mask &= ((unsigned int *) fBitmasks[2].fAllBits)[slice]
1050 )) return 0;
1051 if (crossed && (!(mask &= ~((unsigned int *)crossed->fAllBits)[0])))
1052 return 0;
1053
1054 FindComponentsFastest(mask, list, 0);
1055 }
1056 else
1057 {
1058 unsigned int *masks[3], mask; // masks for X,Y,Z axis
1059 for (auto i = 0; i <= 2; ++i)
1060 {
1061 G4int slice = BinarySearch(fBoundaries[i], point[i]);
1062 masks[i] = ((unsigned int *) fBitmasks[i].fAllBits) + slice*fNPerSlice;
1063 }
1064 unsigned int *maskCrossed =
1065 crossed ? (unsigned int *)crossed->fAllBits : 0;
1066
1067 for (G4int i = 0 ; i < fNPerSlice; ++i)
1068 {
1069 // Logic "and" of the masks along the 3 axes x, y, z:
1070 // removing "if (!" and ") continue" => slightly slower
1071 //
1072 if (!(mask = masks[0][i])) continue;
1073 if (!(mask &= masks[1][i])) continue;
1074 if (!(mask &= masks[2][i])) continue;
1075 if (maskCrossed && !(mask &= ~maskCrossed[i])) continue;
1076
1077 FindComponentsFastest(mask, list, i);
1078 }
1079 }
1080*/
1081 }
1082 return (G4int)list.size();
1083}
1084
1085//______________________________________________________________________________
1086G4int
1087G4Voxelizer::GetCandidatesVoxelArray(const std::vector<G4int>& voxels,
1088 const G4SurfBits bitmasks[],
1089 std::vector<G4int>& list,
1090 G4SurfBits* crossed) const
1091{
1092 list.clear();
1093
1094 if (fTotalCandidates == 1)
1095 {
1096 list.push_back(0);
1097 return 1;
1098 }
1099 else
1100 {
1101 if (fNPerSlice == 1)
1102 {
1103 unsigned int mask;
1104 if (!(mask = ((unsigned int *) bitmasks[0].fAllBits)[voxels[0]]))
1105 return 0;
1106 if (!(mask &= ((unsigned int *) bitmasks[1].fAllBits)[voxels[1]]))
1107 return 0;
1108 if (!(mask &= ((unsigned int *) bitmasks[2].fAllBits)[voxels[2]]))
1109 return 0;
1110 if (crossed && (!(mask &= ~((unsigned int *)crossed->fAllBits)[0])))
1111 return 0;
1112
1113 FindComponentsFastest(mask, list, 0);
1114 }
1115 else
1116 {
1117 unsigned int *masks[3], mask; // masks for X,Y,Z axis
1118 for (auto i = 0; i <= 2; ++i)
1119 {
1120 masks[i] = ((unsigned int *) bitmasks[i].fAllBits)
1121 + voxels[i]*fNPerSlice;
1122 }
1123 unsigned int *maskCrossed = crossed != nullptr
1124 ? (unsigned int *)crossed->fAllBits : 0;
1125
1126 for (G4int i = 0 ; i < fNPerSlice; ++i)
1127 {
1128 // Logic "and" of the masks along the 3 axes x, y, z:
1129 // removing "if (!" and ") continue" => slightly slower
1130 //
1131 if (!(mask = masks[0][i])) continue;
1132 if (!(mask &= masks[1][i])) continue;
1133 if (!(mask &= masks[2][i])) continue;
1134 if (maskCrossed && !(mask &= ~maskCrossed[i])) continue;
1135
1136 FindComponentsFastest(mask, list, i);
1137 }
1138 }
1139 }
1140 return (G4int)list.size();
1141}
1142
1143//______________________________________________________________________________
1144G4int
1145G4Voxelizer::GetCandidatesVoxelArray(const std::vector<G4int>& voxels,
1146 std::vector<G4int>& list, G4SurfBits* crossed) const
1147{
1148 // Method returning the candidates corresponding to the passed point
1149
1150 return GetCandidatesVoxelArray(voxels, fBitmasks, list, crossed);
1151}
1152
1153//______________________________________________________________________________
1155{
1156 for (auto i = 0; i < 3; ++i)
1157 {
1158 if (point[i] < fBoundaries[i].front() || point[i] > fBoundaries[i].back())
1159 return false;
1160 }
1161 return true;
1162}
1163
1164//______________________________________________________________________________
1167 const G4ThreeVector& direction) const
1168{
1169 G4ThreeVector pointShifted = point - fBoundingBoxCenter;
1170 G4double shift = fBoundingBox.DistanceToIn(pointShifted, direction);
1171 return shift;
1172}
1173
1174//______________________________________________________________________________
1177{
1178 G4ThreeVector pointShifted = point - fBoundingBoxCenter;
1179 G4double shift = MinDistanceToBox(pointShifted, fBoundingBoxSize);
1180 return shift;
1181}
1182
1183//______________________________________________________________________________
1186 const G4ThreeVector& f)
1187{
1188 // Estimates the isotropic safety from a point outside the current solid to
1189 // any of its surfaces. The algorithm may be accurate or should provide a
1190 // fast underestimate.
1191
1192 G4double safe, safx, safy, safz;
1193 safe = safx = -f.x() + std::abs(aPoint.x());
1194 safy = -f.y() + std::abs(aPoint.y());
1195 if ( safy > safe ) safe = safy;
1196 safz = -f.z() + std::abs(aPoint.z());
1197 if ( safz > safe ) safe = safz;
1198 if (safe < 0.0) return 0.0; // point is inside
1199
1200 G4double safsq = 0.0;
1201 G4int count = 0;
1202 if ( safx > 0 ) { safsq += safx*safx; ++count; }
1203 if ( safy > 0 ) { safsq += safy*safy; ++count; }
1204 if ( safz > 0 ) { safsq += safz*safz; ++count; }
1205 if (count == 1) return safe;
1206 return std::sqrt(safsq);
1207}
1208
1209//______________________________________________________________________________
1212 const G4ThreeVector& direction,
1213 std::vector<G4int>& curVoxel) const
1214{
1215 G4double shift = kInfinity;
1216
1217 G4int cur = 0; // the smallest index, which would be than increased
1218 for (G4int i = 0; i <= 2; ++i)
1219 {
1220 // Looking for the next voxels on the considered direction X,Y,Z axis
1221 //
1222 const std::vector<G4double>& boundary = fBoundaries[i];
1223 G4int index = curVoxel[i];
1224 if (direction[i] >= 1e-10)
1225 {
1226 ++index;
1227 }
1228 else
1229 {
1230 if (direction[i] > -1e-10)
1231 continue;
1232 }
1233 G4double dif = boundary[index] - point[i];
1234 G4double distance = dif / direction[i];
1235
1236 if (shift > distance)
1237 {
1238 shift = distance;
1239 cur = i;
1240 }
1241 }
1242
1243 if (shift != kInfinity)
1244 {
1245 // updating current voxel using the index corresponding
1246 // to the closest voxel boundary on the ray
1247
1248 if (direction[cur] > 0)
1249 {
1250 if (++curVoxel[cur] >= (G4int) fBoundaries[cur].size() - 1)
1251 shift = kInfinity;
1252 }
1253 else
1254 {
1255 if (--curVoxel[cur] < 0)
1256 shift = kInfinity;
1257 }
1258 }
1259
1260/*
1261 for (auto i = 0; i <= 2; ++i)
1262 {
1263 // Looking for the next voxels on the considered direction X,Y,Z axis
1264 //
1265 const std::vector<G4double> &boundary = fBoundaries[i];
1266 G4int cur = curVoxel[i];
1267 if(direction[i] >= 1e-10)
1268 {
1269 if (boundary[++cur] - point[i] < fTolerance) // make sure shift would
1270 if (++cur >= (G4int) boundary.size()) // be non-zero
1271 continue;
1272 }
1273 else
1274 {
1275 if(direction[i] <= -1e-10)
1276 {
1277 if (point[i] - boundary[cur] < fTolerance) // make sure shift would
1278 if (--cur < 0) // be non-zero
1279 continue;
1280 }
1281 else
1282 continue;
1283 }
1284 G4double dif = boundary[cur] - point[i];
1285 G4double distance = dif / direction[i];
1286
1287 if (shift > distance)
1288 shift = distance;
1289 }
1290*/
1291 return shift;
1292}
1293
1294//______________________________________________________________________________
1295G4bool
1297 const G4ThreeVector& direction,
1298 std::vector<G4int>& curVoxel) const
1299{
1300 for (auto i = 0; i <= 2; ++i)
1301 {
1302 G4int index = curVoxel[i];
1303 const std::vector<G4double> &boundary = fBoundaries[i];
1304
1305 if (direction[i] > 0)
1306 {
1307 if (point[i] >= boundary[++index])
1308 if (++curVoxel[i] >= (G4int) boundary.size() - 1)
1309 return false;
1310 }
1311 else
1312 {
1313 if (point[i] < boundary[index])
1314 if (--curVoxel[i] < 0)
1315 return false;
1316 }
1317#ifdef G4SPECSDEBUG
1318 G4int indexOK = BinarySearch(boundary, point[i]);
1319 if (curVoxel[i] != indexOK)
1320 curVoxel[i] = indexOK; // put breakpoint here
1321#endif
1322 }
1323 return true;
1324}
1325
1326//______________________________________________________________________________
1328{
1329 fMaxVoxels = max;
1330 fReductionRatio.set(0,0,0);
1331}
1332
1333//______________________________________________________________________________
1334void G4Voxelizer::SetMaxVoxels(const G4ThreeVector& ratioOfReduction)
1335{
1336 fMaxVoxels = -1;
1337 fReductionRatio = ratioOfReduction;
1338}
1339
1340//______________________________________________________________________________
1342{
1343 fDefaultVoxelsCount = count;
1344}
1345
1346//______________________________________________________________________________
1348{
1349 return fDefaultVoxelsCount;
1350}
1351
1352//______________________________________________________________________________
1354{
1355 G4int size = fEmpty.GetNbytes();
1356 size += fBoxes.capacity() * sizeof(G4VoxelBox);
1357 size += sizeof(G4double) * (fBoundaries[0].capacity()
1358 + fBoundaries[1].capacity() + fBoundaries[2].capacity());
1359 size += sizeof(G4int) * (fCandidatesCounts[0].capacity()
1360 + fCandidatesCounts[1].capacity() + fCandidatesCounts[2].capacity());
1361 size += fBitmasks[0].GetNbytes() + fBitmasks[1].GetNbytes()
1362 + fBitmasks[2].GetNbytes();
1363
1364 G4int csize = (G4int)fCandidates.size();
1365 for (G4int i = 0; i < csize; ++i)
1366 {
1367 size += sizeof(vector<G4int>) + fCandidates[i].capacity() * sizeof(G4int);
1368 }
1369
1370 return size;
1371}
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition: G4Types.hh:83
long G4long
Definition: G4Types.hh:87
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
double z() const
double x() const
double y() const
void set(double x, double y, double z)
void SetZHalfLength(G4double dz)
Definition: G4Box.cc:172
void SetYHalfLength(G4double dy)
Definition: G4Box.cc:149
void SetXHalfLength(G4double dx)
Definition: G4Box.cc:125
G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const
Definition: G4Box.cc:325
G4double GetSurfaceTolerance() const
static G4GeometryTolerance * GetInstance()
Definition: G4Orb.hh:56
G4double GetRadialTolerance() const
static void DeRegister(G4VSolid *pSolid)
static G4SolidStore * GetInstance()
unsigned int GetNbytes() const
Definition: G4SurfBits.hh:93
unsigned char * fAllBits
Definition: G4SurfBits.hh:104
void ResetAllBits(G4bool value=false)
Definition: G4SurfBits.cc:155
void Clear()
Definition: G4SurfBits.cc:89
void ResetBitNumber(unsigned int bitnumber)
Definition: G4SurfBits.hh:154
void SetBitNumber(unsigned int bitnumber, G4bool value=true)
Definition: G4SurfBits.hh:114
void set(unsigned int nbits, const char *array)
Definition: G4SurfBits.cc:176
G4bool TestBitNumber(unsigned int bitnumber) const
Definition: G4SurfBits.hh:141
virtual G4double Extent(const G4ThreeVector)=0
virtual void BoundingLimits(G4ThreeVector &pMin, G4ThreeVector &pMax) const
Definition: G4VSolid.cc:665
virtual G4GeometryType GetEntityType() const =0
long long CountVoxels(std::vector< G4double > boundaries[]) const
G4double DistanceToBoundingBox(const G4ThreeVector &point) const
G4bool UpdateCurrentVoxel(const G4ThreeVector &point, const G4ThreeVector &direction, std::vector< G4int > &curVoxel) const
void DisplayListNodes() const
Definition: G4Voxelizer.cc:424
void DisplayVoxelLimits() const
Definition: G4Voxelizer.cc:197
static void SetDefaultVoxelsCount(G4int count)
G4int GetBitsPerSlice() const
G4double DistanceToFirst(const G4ThreeVector &point, const G4ThreeVector &direction) const
G4int GetCandidatesVoxelArray(const G4ThreeVector &point, std::vector< G4int > &list, G4SurfBits *crossed=nullptr) const
Definition: G4Voxelizer.cc:966
static G4double MinDistanceToBox(const G4ThreeVector &aPoint, const G4ThreeVector &f)
void SetMaxVoxels(G4int max)
void GetCandidatesVoxel(std::vector< G4int > &voxels)
Definition: G4Voxelizer.cc:889
G4int AllocatedMemory()
static G4int GetDefaultVoxelsCount()
void DisplayBoundaries()
Definition: G4Voxelizer.cc:319
static G4int BinarySearch(const std::vector< T > &vec, T value)
G4double DistanceToNext(const G4ThreeVector &point, const G4ThreeVector &direction, std::vector< G4int > &curVoxel) const
void Voxelize(std::vector< G4VSolid * > &solids, std::vector< G4Transform3D > &transforms)
Definition: G4Voxelizer.cc:706
G4int GetVoxelsIndex(G4int x, G4int y, G4int z) const
G4bool Contains(const G4ThreeVector &point) const
G4double total(Particle const *const p1, Particle const *const p2)
T max(const T t1, const T t2)
brief Return the largest of the two arguments
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
G4ThreeVector hlen
Definition: G4Voxelizer.hh:51
G4ThreeVector pos
Definition: G4Voxelizer.hh:52
G4int previous
Definition: G4Voxelizer.hh:58
#define G4ThreadLocal
Definition: tls.hh:77