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