Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4SmartVoxelHeader.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// G4SmartVoxelHeader implementation
27//
28// Define G4GEOMETRY_VOXELDEBUG for debugging information on G4cout
29//
30// 13.07.95 Initial version, stubb definitions only - P.Kent.
31// 21.07.95 Full implementation, supporting non divided volumes - P.Kent.
32// 11.02.99 Voxels at lower levels now built for collapsed slices - S.Giani.
33// 12.02.99 Introduction of new quality/smartless - S.Giani.
34// 18.04.01 Migrated to STL vector - G.Cosmo.
35// 29.04.02 Use 3D voxelisation for non consuming replication - G.Cosmo.
36// --------------------------------------------------------------------
37
38#include "G4SmartVoxelHeader.hh"
39
40#include "G4ios.hh"
41
42#include "G4LogicalVolume.hh"
43#include "G4VPhysicalVolume.hh"
44#include "G4VoxelLimits.hh"
45
46#include "voxeldefs.hh"
47#include "G4AffineTransform.hh"
48#include "G4VSolid.hh"
50
51// ***************************************************************************
52// Constructor for topmost header, to begin voxel construction at a
53// given logical volume.
54// Constructs target list of volumes, calls "Build and refine" constructor.
55// Assumes all daughters represent single volumes (ie. no divisions
56// or parametric)
57// ***************************************************************************
58//
60 G4int pSlice)
61 : fminEquivalent(pSlice),
62 fmaxEquivalent(pSlice),
63 fparamAxis(kUndefined)
64{
65 std::size_t nDaughters = pVolume->GetNoDaughters();
66
67 // Determine whether daughter is replicated
68 //
69 if ((nDaughters!=1) || (!pVolume->GetDaughter(0)->IsReplicated()))
70 {
71 // Daughter not replicated => conventional voxel Build
72 // where each daughters extents are computed
73 //
74 BuildVoxels(pVolume);
75 }
76 else
77 {
78 // Single replicated daughter
79 //
80 BuildReplicaVoxels(pVolume);
81 }
82}
83
84// ***************************************************************************
85// Builds and refines voxels between specified limits, considering only
86// the physical volumes numbered `pCandidates'. `pSlice' is used to set max
87// and min equivalent slice nos for the header - they apply to the level
88// of the header, not its nodes.
89// ***************************************************************************
90//
92 const G4VoxelLimits& pLimits,
93 const G4VolumeNosVector* pCandidates,
94 G4int pSlice)
95 : fminEquivalent(pSlice),
96 fmaxEquivalent(pSlice),
97 fparamAxis(kUndefined)
98{
99#ifdef G4GEOMETRY_VOXELDEBUG
100 G4cout << "**** G4SmartVoxelHeader::G4SmartVoxelHeader" << G4endl
101 << " Limits " << pLimits << G4endl
102 << " Candidate #s = " ;
103 for (auto i=0; i<pCandidates->size(); ++i)
104 {
105 G4cout << (*pCandidates)[i] << " ";
106 }
107 G4cout << G4endl;
108#endif
109
110 BuildVoxelsWithinLimits(pVolume,pLimits,pCandidates);
111}
112
113// ***************************************************************************
114// Destructor: deletes all proxies and underlying objects.
115// ***************************************************************************
116//
118{
119 // Manually destroy underlying nodes/headers
120 // Delete collected headers and nodes once only
121 //
122 std::size_t node, proxy, maxNode=fslices.size();
123 G4SmartVoxelProxy* lastProxy = nullptr;
124 G4SmartVoxelNode *dyingNode, *lastNode=nullptr;
125 G4SmartVoxelHeader *dyingHeader, *lastHeader=nullptr;
126
127 for (node=0; node<maxNode; ++node)
128 {
129 if (fslices[node]->IsHeader())
130 {
131 dyingHeader = fslices[node]->GetHeader();
132 if (lastHeader != dyingHeader)
133 {
134 lastHeader = dyingHeader;
135 lastNode = nullptr;
136 delete dyingHeader;
137 }
138 }
139 else
140 {
141 dyingNode = fslices[node]->GetNode();
142 if (dyingNode != lastNode)
143 {
144 lastNode = dyingNode;
145 lastHeader = nullptr;
146 delete dyingNode;
147 }
148 }
149 }
150 // Delete proxies
151 //
152 for (proxy=0; proxy<maxNode; ++proxy)
153 {
154 if (fslices[proxy] != lastProxy)
155 {
156 lastProxy = fslices[proxy];
157 delete lastProxy;
158 }
159 }
160 // Don't need to clear slices
161 // fslices.clear();
162}
163
164// ***************************************************************************
165// Equality operator: returns true if contents are equivalent.
166// Implies a deep search through contained nodes/header.
167// Compares headers' axes,sizes,extents. Returns false if different.
168// For each contained proxy, determines whether node/header, compares and
169// returns if different. Compares and returns if proxied nodes/headers
170// are different.
171// ***************************************************************************
172//
174{
175 if ( (GetAxis() == pHead.GetAxis())
176 && (GetNoSlices() == pHead.GetNoSlices())
177 && (GetMinExtent() == pHead.GetMinExtent())
178 && (GetMaxExtent() == pHead.GetMaxExtent()) )
179 {
180 std::size_t node, maxNode;
181 G4SmartVoxelProxy *leftProxy, *rightProxy;
182 G4SmartVoxelHeader *leftHeader, *rightHeader;
183 G4SmartVoxelNode *leftNode, *rightNode;
184
185 maxNode = GetNoSlices();
186 for (node=0; node<maxNode; ++node)
187 {
188 leftProxy = GetSlice(node);
189 rightProxy = pHead.GetSlice(node);
190 if (leftProxy->IsHeader())
191 {
192 if (rightProxy->IsNode())
193 {
194 return false;
195 }
196 leftHeader = leftProxy->GetHeader();
197 rightHeader = rightProxy->GetHeader();
198 if (!(*leftHeader == *rightHeader))
199 {
200 return false;
201 }
202 }
203 else
204 {
205 if (rightProxy->IsHeader())
206 {
207 return false;
208 }
209 leftNode = leftProxy->GetNode();
210 rightNode = rightProxy->GetNode();
211 if (!(*leftNode == *rightNode))
212 {
213 return false;
214 }
215 }
216 }
217 return true;
218 }
219
220 return false;
221}
222
223// ***************************************************************************
224// Builds voxels for daughters specified volume, in NON-REPLICATED case
225// o Create List of target volume nos (all daughters; 0->noDaughters-1)
226// o BuildWithinLimits does Build & also determines mother dimensions.
227// ***************************************************************************
228//
229void G4SmartVoxelHeader::BuildVoxels(G4LogicalVolume* pVolume)
230{
231 G4VoxelLimits limits; // Create `unlimited' limits object
232 std::size_t nDaughters = pVolume->GetNoDaughters();
233
234 G4VolumeNosVector targetList;
235 targetList.reserve(nDaughters);
236 for (std::size_t i=0; i<nDaughters; ++i)
237 {
238 targetList.push_back((G4int)i);
239 }
240 BuildVoxelsWithinLimits(pVolume, limits, &targetList);
241}
242
243// ***************************************************************************
244// Builds voxels for specified volume containing a single replicated volume.
245// If axis is not specified (i.e. "kUndefined"), 3D voxelisation is applied,
246// and the best axis is determined according to heuristics as for placements.
247// ***************************************************************************
248//
249void G4SmartVoxelHeader::BuildReplicaVoxels(G4LogicalVolume* pVolume)
250{
251 G4VPhysicalVolume* pDaughter = nullptr;
252
253 // Replication data
254 //
255 EAxis axis;
256 G4int nReplicas;
257 G4double width,offset;
258 G4bool consuming;
259
260 // Consistency check: pVolume should contain single replicated volume
261 //
262 if ( (pVolume->GetNoDaughters()==1)
263 && (pVolume->GetDaughter(0)->IsReplicated()) )
264 {
265 // Obtain replication data
266 //
267 pDaughter = pVolume->GetDaughter(0);
268 pDaughter->GetReplicationData(axis,nReplicas,width,offset,consuming);
269 fparamAxis = axis;
270 if ( !consuming )
271 {
272 G4VoxelLimits limits; // Create `unlimited' limits object
273 G4VolumeNosVector targetList;
274 targetList.reserve(nReplicas);
275 for (auto i=0; i<nReplicas; ++i)
276 {
277 targetList.push_back(i);
278 }
279 if (axis != kUndefined)
280 {
281 // Apply voxelisation along the specified axis only
282
283 G4ProxyVector* pSlices=BuildNodes(pVolume,limits,&targetList,axis);
284 faxis = axis;
285 fslices = *pSlices;
286 delete pSlices;
287
288 // Calculate and set min and max extents given our axis
289 //
290 const G4AffineTransform origin;
291 pVolume->GetSolid()->CalculateExtent(faxis, limits, origin,
292 fminExtent, fmaxExtent);
293 // Calculate equivalent nos
294 //
295 BuildEquivalentSliceNos();
296 CollectEquivalentNodes(); // Collect common nodes
297 }
298 else
299 {
300 // Build voxels similarly as for normal placements considering
301 // all three cartesian axes.
302
303 BuildVoxelsWithinLimits(pVolume, limits, &targetList);
304 }
305 }
306 else
307 {
308 // Replication is consuming -> Build voxels directly
309 //
310 // o Cartesian axes - range is -width*nREplicas/2 to +width*nREplicas/2
311 // nReplicas replications result
312 // o Radial axis (rho) = range is 0 to width*nReplicas
313 // nReplicas replications result
314 // o Phi axi - range is offset to offset+width*nReplicas radians
315 //
316 // Equivalent slices no computation & collection not required - all
317 // slices are different
318 //
319 switch (axis)
320 {
321 case kXAxis:
322 case kYAxis:
323 case kZAxis:
324 fminExtent = -width*nReplicas*0.5;
325 fmaxExtent = width*nReplicas*0.5;
326 break;
327 case kRho:
328 fminExtent = offset;
329 fmaxExtent = width*nReplicas+offset;
330 break;
331 case kPhi:
332 fminExtent = offset;
333 fmaxExtent = offset+width*nReplicas;
334 break;
335 default:
336 G4Exception("G4SmartVoxelHeader::BuildReplicaVoxels()",
337 "GeomMgt0002", FatalException, "Illegal axis.");
338 break;
339 }
340 faxis = axis; // Set axis
341 BuildConsumedNodes(nReplicas);
342 if ( (axis==kXAxis) || (axis==kYAxis) || (axis==kZAxis) )
343 {
344 // Sanity check on extent
345 //
346 G4double emin = kInfinity, emax = -kInfinity;
347 G4VoxelLimits limits;
348 G4AffineTransform origin;
349 pVolume->GetSolid()->CalculateExtent(axis, limits, origin, emin, emax);
350 if ( (std::fabs((emin-fminExtent)/fminExtent) +
351 std::fabs((emax-fmaxExtent)/fmaxExtent)) > 0.05)
352 {
353 std::ostringstream message;
354 message << "Sanity check: wrong solid extent." << G4endl
355 << " Replicated geometry, logical volume: "
356 << pVolume->GetName();
357 G4Exception("G4SmartVoxelHeader::BuildReplicaVoxels",
358 "GeomMgt0002", FatalException, message);
359 }
360 }
361 }
362 }
363 else
364 {
365 G4Exception("G4SmartVoxelHeader::BuildReplicaVoxels", "GeomMgt0002",
366 FatalException, "Only one replicated daughter is allowed !");
367 }
368}
369
370// ***************************************************************************
371// Builds 'consumed nodes': nReplicas nodes each containing one replication,
372// numbered in sequence 0->nReplicas-1
373// o Modifies fslices 'in place'
374// o faxis,fminExtent,fmaxExtent NOT modified.
375// ***************************************************************************
376//
377void G4SmartVoxelHeader::BuildConsumedNodes(G4int nReplicas)
378{
379 G4int nNode, nVol;
380 G4SmartVoxelNode* pNode;
381 G4SmartVoxelProxy* pProxyNode;
382
383 // Create and fill nodes in temporary G4NodeVector (on stack)
384 //
385 G4NodeVector nodeList;
386 nodeList.reserve(nReplicas);
387 for (nNode=0; nNode<nReplicas; ++nNode)
388 {
389 pNode = new G4SmartVoxelNode(nNode);
390 if (pNode == nullptr)
391 {
392 G4Exception("G4SmartVoxelHeader::BuildConsumedNodes()", "GeomMgt0003",
393 FatalException, "Node allocation error.");
394 }
395 nodeList.push_back(pNode);
396 }
397 for (nVol=0; nVol<nReplicas; ++nVol)
398 {
399 nodeList[nVol]->Insert(nVol); // Insert replication of number
400 } // identical to voxel number
401
402 // Create & fill proxy List `in place' by modifying instance data fslices
403 //
404 fslices.clear();
405 for (nNode=0; nNode<nReplicas; ++nNode)
406 {
407 pProxyNode = new G4SmartVoxelProxy(nodeList[nNode]);
408 if (pProxyNode == nullptr)
409 {
410 G4Exception("G4SmartVoxelHeader::BuildConsumedNodes()", "GeomMgt0003",
411 FatalException, "Proxy node allocation error.");
412 }
413 fslices.push_back(pProxyNode);
414 }
415}
416
417// ***************************************************************************
418// Builds and refines voxels between specified limits, considering only
419// the physical volumes numbered 'pCandidates'.
420// o Chooses axis
421// o Determines min and max extents (of mother solid) within limits.
422// ***************************************************************************
423//
424void
425G4SmartVoxelHeader::BuildVoxelsWithinLimits(G4LogicalVolume* pVolume,
426 G4VoxelLimits pLimits,
427 const G4VolumeNosVector* pCandidates)
428{
429 // Choose best axis for slicing by:
430 // 1. Trying all unlimited cartesian axes
431 // 2. Select axis which gives greatest no slices
432
433 G4ProxyVector *pGoodSlices=nullptr, *pTestSlices;
434 G4double goodSliceScore=kInfinity, testSliceScore;
435 EAxis goodSliceAxis = kXAxis;
436 std::size_t node, maxNode;
437 G4VoxelLimits noLimits;
438
439 // Try all non-limited cartesian axes
440 //
441 for ( EAxis testAxis : { kXAxis, kYAxis, kZAxis } )
442 {
443 if (!pLimits.IsLimited(testAxis))
444 {
445 pTestSlices = BuildNodes(pVolume,pLimits,pCandidates,testAxis);
446 testSliceScore = CalculateQuality(pTestSlices);
447 if ( (pGoodSlices == nullptr) || (testSliceScore<goodSliceScore) )
448 {
449 goodSliceAxis = testAxis;
450 goodSliceScore = testSliceScore;
451 std::swap( pGoodSlices, pTestSlices);
452 }
453 if (pTestSlices != nullptr)
454 {
455 // Destroy pTestSlices and all its contents
456 //
457 maxNode=pTestSlices->size();
458 for (node=0; node<maxNode; ++node)
459 {
460 delete (*pTestSlices)[node]->GetNode();
461 }
462 G4SmartVoxelProxy* tmpProx;
463 while (!pTestSlices->empty()) // Loop checking, 06.08.2015, G.Cosmo
464 {
465 tmpProx = pTestSlices->back();
466 pTestSlices->pop_back();
467 for (auto i=pTestSlices->cbegin(); i!=pTestSlices->cend(); )
468 {
469 if (*i==tmpProx)
470 {
471 i = pTestSlices->erase(i);
472 }
473 else
474 {
475 ++i;
476 }
477 }
478 delete tmpProx;
479 }
480 delete pTestSlices;
481 }
482 }
483 }
484 // Check for error case.. when limits already 3d,
485 // so cannot select a new axis
486 //
487 if (pGoodSlices == nullptr)
488 {
489 G4Exception("G4SmartVoxelHeader::BuildVoxelsWithinLimits()",
490 "GeomMgt0002", FatalException,
491 "Cannot select more than 3 axis for optimisation.");
492 return;
493 }
494
495 //
496 // We have selected pGoodSlices, with a score testSliceScore
497 //
498
499 // Store chosen axis, slice ptr
500 //
501 fslices =* pGoodSlices; // Set slice information, copy ptrs in collection
502 delete pGoodSlices; // Destroy slices vector, but not contained
503 // proxies or nodes
504 faxis = goodSliceAxis;
505
506#ifdef G4GEOMETRY_VOXELDEBUG
507 G4cout << G4endl << " Volume = " << pVolume->GetName()
508 << G4endl << " Selected axis = " << faxis << G4endl;
509 for (auto islice=0; islice<fslices.size(); ++islice)
510 {
511 G4cout << " Node #" << islice << " = {";
512 for (auto j=0; j<fslices[islice]->GetNode()->GetNoContained(); ++j)
513 {
514 G4cout << " " << fslices[islice]->GetNode()->GetVolume(j);
515 }
516 G4cout << " }" << G4endl;
517 }
518 G4cout << G4endl;
519#endif
520
521 // Calculate and set min and max extents given our axis
522 //
523 G4VSolid* outerSolid = pVolume->GetSolid();
524 const G4AffineTransform origin;
525 if(!outerSolid->CalculateExtent(faxis,pLimits,origin,fminExtent,fmaxExtent))
526 {
527 outerSolid->CalculateExtent(faxis,noLimits,origin,fminExtent,fmaxExtent);
528 }
529
530 // Calculate equivalent nos
531 //
532 BuildEquivalentSliceNos();
533 CollectEquivalentNodes(); // Collect common nodes
534 RefineNodes(pVolume, pLimits); // Refine nodes creating headers
535
536 // No common headers can exist because collapsed by construction
537}
538
539// ***************************************************************************
540// Calculates and stores the minimum and maximum equivalent neighbour
541// values for all slices at our level.
542//
543// Precondition: all slices are nodes.
544// For each potential start of a group of equivalent nodes:
545// o searches forwards in fslices to find group end
546// o loops from start to end setting start and end slices.
547// ***************************************************************************
548//
549void G4SmartVoxelHeader::BuildEquivalentSliceNos()
550{
551 std::size_t sliceNo, minNo, maxNo, equivNo;
552 std::size_t maxNode = fslices.size();
553 G4SmartVoxelNode *startNode, *sampleNode;
554 for (sliceNo=0; sliceNo<maxNode; ++sliceNo)
555 {
556 minNo = sliceNo;
557
558 // Get first node (see preconditions - will throw exception if a header)
559 //
560 startNode = fslices[minNo]->GetNode();
561
562 // Find max equivalent
563 //
564 for (equivNo=minNo+1; equivNo<maxNode; ++equivNo)
565 {
566 sampleNode = fslices[equivNo]->GetNode();
567 if (!((*startNode) == (*sampleNode))) { break; }
568 }
569 maxNo = equivNo-1;
570 if (maxNo != minNo)
571 {
572 // Set min and max nos
573 //
574 for (equivNo=minNo; equivNo<=maxNo; ++equivNo)
575 {
576 sampleNode = fslices[equivNo]->GetNode();
577 sampleNode->SetMinEquivalentSliceNo((G4int)minNo);
578 sampleNode->SetMaxEquivalentSliceNo((G4int)maxNo);
579 }
580 // Advance outer loop to end of equivalent group
581 //
582 sliceNo = maxNo;
583 }
584 }
585}
586
587// ***************************************************************************
588// Collects common nodes at our level, deleting all but one to save
589// memory, and adjusting stored slice pointers appropriately.
590//
591// Preconditions:
592// o the slices have not previously be "collected"
593// o all of the slices are nodes.
594// ***************************************************************************
595//
596void G4SmartVoxelHeader::CollectEquivalentNodes()
597{
598 std::size_t sliceNo, maxNo, equivNo;
599 std::size_t maxNode=fslices.size();
600 G4SmartVoxelNode* equivNode;
601 G4SmartVoxelProxy* equivProxy;
602
603 for (sliceNo=0; sliceNo<maxNode; ++sliceNo)
604 {
605 equivProxy=fslices[sliceNo];
606
607 // Assumption (see preconditions): all slices are nodes
608 //
609 equivNode = equivProxy->GetNode();
610 maxNo = equivNode->GetMaxEquivalentSliceNo();
611 if (maxNo != sliceNo)
612 {
613#ifdef G4GEOMETRY_VOXELDEBUG
614 G4cout << "**** G4SmartVoxelHeader::CollectEquivalentNodes" << G4endl
615 << " Collecting Nodes = "
616 << sliceNo << " - " << maxNo << G4endl;
617#endif
618 // Do collection between sliceNo and maxNo inclusive
619 //
620 for (equivNo=sliceNo+1; equivNo<=maxNo; ++equivNo)
621 {
622 delete fslices[equivNo]->GetNode();
623 delete fslices[equivNo];
624 fslices[equivNo] = equivProxy;
625 }
626 sliceNo = maxNo;
627 }
628 }
629}
630
631// ***************************************************************************
632// Collects common headers at our level, deleting all but one to save
633// memory, and adjusting stored slice pointers appropriately.
634//
635// Preconditions:
636// o if a header forms part of a range of equivalent slices
637// (ie. GetMaxEquivalentSliceNo()>GetMinEquivalentSliceNo()),
638// it is assumed that all slices in the range are headers.
639// o this will be true if a constant Expression is used to evaluate
640// when to refine nodes.
641// ***************************************************************************
642//
643void G4SmartVoxelHeader::CollectEquivalentHeaders()
644{
645 std::size_t sliceNo, maxNo, equivNo;
646 std::size_t maxNode = fslices.size();
647 G4SmartVoxelHeader *equivHeader, *sampleHeader;
648 G4SmartVoxelProxy *equivProxy;
649
650 for (sliceNo=0; sliceNo<maxNode; ++sliceNo)
651 {
652 equivProxy = fslices[sliceNo];
653 if (equivProxy->IsHeader())
654 {
655 equivHeader = equivProxy->GetHeader();
656 maxNo = equivHeader->GetMaxEquivalentSliceNo();
657 if (maxNo != sliceNo)
658 {
659 // Attempt collection between sliceNo and maxNo inclusive:
660 // look for common headers. All slices between sliceNo and maxNo
661 // are guaranteed to be headers but may not have equal contents
662 //
663#ifdef G4GEOMETRY_VOXELDEBUG
664 G4cout << "**** G4SmartVoxelHeader::CollectEquivalentHeaders" << G4endl
665 << " Collecting Headers =";
666#endif
667 for (equivNo=sliceNo+1; equivNo<=maxNo; ++equivNo)
668 {
669 sampleHeader = fslices[equivNo]->GetHeader();
670 if ( (*sampleHeader) == (*equivHeader) )
671 {
672#ifdef G4GEOMETRY_VOXELDEBUG
673 G4cout << " " << equivNo;
674#endif
675 // Delete sampleHeader + proxy and replace with equivHeader/Proxy
676 //
677 delete sampleHeader;
678 delete fslices[equivNo];
679 fslices[equivNo] = equivProxy;
680 }
681 else
682 {
683 // Not equal. Set this header to be
684 // the current header for comparisons
685 //
686 equivProxy = fslices[equivNo];
687 equivHeader = equivProxy->GetHeader();
688 }
689
690 }
691#ifdef G4GEOMETRY_VOXELDEBUG
692 G4cout << G4endl;
693#endif
694 // Skip past examined slices
695 //
696 sliceNo = maxNo;
697 }
698 }
699 }
700}
701
702// ***************************************************************************
703// Builds the nodes corresponding to slices between the specified limits
704// and along the specified axis, using candidate volume no.s in the vector
705// pCandidates. If the 'daughters' are replicated volumes (ie. the logical
706// volume has a single replicated/parameterised volume for a daughter)
707// the candidate no.s are interpreted as PARAMETERISED volume no.s &
708// PARAMETERISATIONs are applied to compute transformations & solid
709// dimensions appropriately. The volume must be parameterised - ie. has a
710// parameterisation object & non-consuming) - in this case.
711//
712// Returns pointer to built node "structure" (guaranteed non NULL) consisting
713// of G4SmartVoxelNodeProxies referring to G4SmartVoxelNodes.
714// ***************************************************************************
715//
716G4ProxyVector* G4SmartVoxelHeader::BuildNodes(G4LogicalVolume* pVolume,
717 G4VoxelLimits pLimits,
718 const G4VolumeNosVector* pCandidates,
719 EAxis pAxis)
720{
721 G4double motherMinExtent= kInfinity, motherMaxExtent= -kInfinity,
722 targetMinExtent= kInfinity, targetMaxExtent= -kInfinity;
723 G4VPhysicalVolume* pDaughter = nullptr;
724 G4VPVParameterisation* pParam = nullptr;
725 G4VSolid *targetSolid;
726 G4AffineTransform targetTransform;
727 G4bool replicated;
728 std::size_t nCandidates = pCandidates->size();
729 std::size_t nVol, nNode, targetVolNo;
730 G4VoxelLimits noLimits;
731
732#ifdef G4GEOMETRY_VOXELDEBUG
733 G4cout << "**** G4SmartVoxelHeader::BuildNodes" << G4endl
734 << " Limits = " << pLimits << G4endl
735 << " Axis = " << pAxis << G4endl
736 << " Candidates = " << nCandidates << G4endl;
737#endif
738
739 // Compute extent of logical volume's solid along this axis
740 // NOTE: results stored locally and not preserved/reused
741 //
742 G4VSolid* outerSolid = pVolume->GetSolid();
743 const G4AffineTransform origin;
744 if( !outerSolid->CalculateExtent(pAxis, pLimits, origin,
745 motherMinExtent, motherMaxExtent) )
746 {
747 outerSolid->CalculateExtent(pAxis, noLimits, origin,
748 motherMinExtent, motherMaxExtent);
749 }
750 G4VolumeExtentVector minExtents(nCandidates,0.);
751 G4VolumeExtentVector maxExtents(nCandidates,0.);
752
753 if ( (pVolume->GetNoDaughters() == 1)
754 && (pVolume->GetDaughter(0)->IsReplicated()) )
755 {
756 // Replication data not required: only parameterisation object
757 // and volume no. List used
758 //
759 pDaughter = pVolume->GetDaughter(0);
760 pParam = pDaughter->GetParameterisation();
761 if (pParam == nullptr)
762 {
763 std::ostringstream message;
764 message << "PANIC! - Missing parameterisation." << G4endl
765 << " Replicated volume with no parameterisation object !";
766 G4Exception("G4SmartVoxelHeader::BuildNodes()", "GeomMgt0003",
767 FatalException, message);
768 return nullptr;
769 }
770
771 // Setup daughter's transformations
772 //
773 targetTransform = G4AffineTransform(pDaughter->GetRotation(),
774 pDaughter->GetTranslation());
775 replicated = true;
776 }
777 else
778 {
779 replicated = false;
780 }
781
782 // Compute extents
783 //
784 for (nVol=0; nVol<nCandidates; ++nVol)
785 {
786 targetVolNo = (*pCandidates)[nVol];
787 if (!replicated)
788 {
789 pDaughter = pVolume->GetDaughter(targetVolNo);
790
791 // Setup daughter's transformations
792 //
793 targetTransform = G4AffineTransform(pDaughter->GetRotation(),
794 pDaughter->GetTranslation());
795 // Get underlying (and setup) solid
796 //
797 targetSolid = pDaughter->GetLogicalVolume()->GetSolid();
798 }
799 else
800 {
801 // Find solid
802 //
803 targetSolid = pParam->ComputeSolid((G4int)targetVolNo,pDaughter);
804
805 // Setup solid
806 //
807 targetSolid->ComputeDimensions(pParam,(G4int)targetVolNo,pDaughter);
808
809 // Setup transform
810 //
811 pParam->ComputeTransformation((G4int)targetVolNo,pDaughter);
812 targetTransform = G4AffineTransform(pDaughter->GetRotation(),
813 pDaughter->GetTranslation());
814 }
815 // Calculate extents
816 //
817 if(!targetSolid->CalculateExtent(pAxis, pLimits, targetTransform,
818 targetMinExtent, targetMaxExtent))
819 {
820 targetSolid->CalculateExtent(pAxis, noLimits, targetTransform,
821 targetMinExtent,targetMaxExtent);
822 }
823 minExtents[nVol] = targetMinExtent;
824 maxExtents[nVol] = targetMaxExtent;
825
826#ifdef G4GEOMETRY_VOXELDEBUG
827 G4cout << "---------------------------------------------------" << G4endl
828 << " Volume = " << pDaughter->GetName() << G4endl
829 << " Min Extent = " << targetMinExtent << G4endl
830 << " Max Extent = " << targetMaxExtent << G4endl
831 << "---------------------------------------------------" << G4endl;
832#endif
833
834 // Check not entirely outside mother when processing toplevel nodes
835 //
836 if ( (!pLimits.IsLimited()) && ((targetMaxExtent<=motherMinExtent)
837 ||(targetMinExtent>=motherMaxExtent)) )
838 {
839 std::ostringstream message;
840 message << "PANIC! - Overlapping daughter with mother volume." << G4endl
841 << " Daughter physical volume "
842 << pDaughter->GetName() << G4endl
843 << " is entirely outside mother logical volume "
844 << pVolume->GetName() << " !!";
845 G4Exception("G4SmartVoxelHeader::BuildNodes()", "GeomMgt0002",
846 FatalException, message);
847 }
848
849#ifdef G4GEOMETRY_VOXELDEBUG
850 // Check for straddling volumes when debugging.
851 // If a volume is >kStraddlePercent percent over the mother
852 // boundary, print a warning.
853 //
854 if (!pLimits.IsLimited())
855 {
856 G4double width;
857 G4int kStraddlePercent = 5;
858 width = maxExtents[nVol]-minExtents[nVol];
859 if ( (((motherMinExtent-minExtents[nVol])*100/width) > kStraddlePercent)
860 ||(((maxExtents[nVol]-motherMaxExtent)*100/width) > kStraddlePercent) )
861 {
862 G4cout << "**** G4SmartVoxelHeader::BuildNodes" << G4endl
863 << " WARNING : Daughter # " << nVol
864 << " name = " << pDaughter->GetName() << G4endl
865 << " Crosses mother boundary of logical volume, name = "
866 << pVolume->GetName() << G4endl
867 << " by more than " << kStraddlePercent
868 << "%" << G4endl;
869 }
870 }
871#endif
872 }
873
874 // Extents of all daughters known
875
876 // Calculate minimum slice width, only including volumes inside the limits
877 //
878 G4double minWidth = kInfinity;
879 G4double currentWidth;
880 for (nVol=0; nVol<nCandidates; ++nVol)
881 {
882 // currentWidth should -always- be a positive value. Inaccurate computed extent
883 // from the solid or situations of malformed geometries (overlaps) may lead to
884 // negative values and therefore unpredictable crashes !
885 //
886 currentWidth = std::abs(maxExtents[nVol]-minExtents[nVol]);
887 if ( (currentWidth<minWidth)
888 && (maxExtents[nVol]>=pLimits.GetMinExtent(pAxis))
889 && (minExtents[nVol]<=pLimits.GetMaxExtent(pAxis)) )
890 {
891 minWidth = currentWidth;
892 }
893 }
894
895 // No. of Nodes formula - nearest integer to
896 // mother width/half min daughter width +1
897 //
898 G4double noNodesExactD = ((motherMaxExtent-motherMinExtent)*2.0/minWidth)+1.0;
899
900 // Compare with "smartless quality", i.e. the average number of slices
901 // used per contained volume.
902 //
903 G4double smartlessComputed = noNodesExactD / nCandidates;
904 G4double smartlessUser = pVolume->GetSmartless();
905 G4double smartless = (smartlessComputed <= smartlessUser)
906 ? smartlessComputed : smartlessUser;
907 G4double noNodesSmart = smartless*nCandidates;
908 auto noNodesExactI = G4int(noNodesSmart);
909 G4long noNodes = ((noNodesSmart-noNodesExactI)>=0.5)
910 ? noNodesExactI+1 : noNodesExactI;
911 if( noNodes == 0 ) { noNodes=1; }
912
913#ifdef G4GEOMETRY_VOXELDEBUG
914 G4cout << " Smartless computed = " << smartlessComputed << G4endl
915 << " Smartless volume = " << smartlessUser
916 << " => # Smartless = " << smartless << G4endl;
917 G4cout << " Min width = " << minWidth
918 << " => # Nodes = " << noNodes << G4endl;
919#endif
920
921 if (noNodes>kMaxVoxelNodes)
922 {
923 noNodes=kMaxVoxelNodes;
924#ifdef G4GEOMETRY_VOXELDEBUG
925 G4cout << " Nodes Clipped to = " << kMaxVoxelNodes << G4endl;
926#endif
927 }
928 G4double nodeWidth = (motherMaxExtent-motherMinExtent)/noNodes;
929
930 // Create G4VoxelNodes. Will Add proxies before setting fslices
931 //
932 auto* nodeList = new G4NodeVector();
933 if (nodeList == nullptr)
934 {
935 G4Exception("G4SmartVoxelHeader::BuildNodes()", "GeomMgt0003",
936 FatalException, "NodeList allocation error.");
937 return nullptr;
938 }
939 nodeList->reserve(noNodes);
940
941 for (nNode=0; G4long(nNode)<noNodes; ++nNode)
942 {
943 G4SmartVoxelNode *pNode;
944 pNode = new G4SmartVoxelNode((G4int)nNode);
945 if (pNode == nullptr)
946 {
947 G4Exception("G4SmartVoxelHeader::BuildNodes()", "GeomMgt0003",
948 FatalException, "Node allocation error.");
949 return nullptr;
950 }
951 nodeList->push_back(pNode);
952 }
953
954 // All nodes created (empty)
955
956 // Fill nodes: Step through extent lists
957 //
958 for (nVol=0; nVol<nCandidates; ++nVol)
959 {
960 G4long nodeNo, minContainingNode, maxContainingNode;
961 minContainingNode = (minExtents[nVol]-motherMinExtent)/nodeWidth;
962 maxContainingNode = (maxExtents[nVol]-motherMinExtent)/nodeWidth;
963
964 // Only add nodes that are inside the limits of the axis
965 //
966 if ( (maxContainingNode>=0) && (minContainingNode<noNodes) )
967 {
968 // If max extent is on max boundary => maxContainingNode=noNodes:
969 // should be one less as nodeList has noNodes entries
970 //
971 if (maxContainingNode>=noNodes)
972 {
973 maxContainingNode = noNodes-1;
974 }
975 //
976 // Protection against protruding volumes
977 //
978 if (minContainingNode<0)
979 {
980 minContainingNode = 0;
981 }
982 for (nodeNo=minContainingNode; nodeNo<=maxContainingNode; ++nodeNo)
983 {
984 (*nodeList)[nodeNo]->Insert((*pCandidates)[nVol]);
985 }
986 }
987 }
988
989 // All nodes filled
990
991 // Create proxy List : caller has deletion responsibility
992 // (but we must delete nodeList *itself* - not the contents)
993 //
994 auto* proxyList = new G4ProxyVector();
995 if (proxyList == nullptr)
996 {
997 G4Exception("G4SmartVoxelHeader::BuildNodes()", "GeomMgt0003",
998 FatalException, "Proxy list allocation error.");
999 return nullptr;
1000 }
1001 proxyList->reserve(noNodes);
1002
1003 //
1004 // Fill proxy List
1005 //
1006 for (nNode=0; G4long(nNode)<noNodes; ++nNode)
1007 {
1008 // Get rid of possible excess capacity in the internal node vector
1009 //
1010 ((*nodeList)[nNode])->Shrink();
1011 auto* pProxyNode = new G4SmartVoxelProxy((*nodeList)[nNode]);
1012 if (pProxyNode == nullptr)
1013 {
1014 G4Exception("G4SmartVoxelHeader::BuildNodes()", "GeomMgt0003",
1015 FatalException, "Proxy node allocation failed.");
1016 return nullptr;
1017 }
1018 proxyList->push_back(pProxyNode);
1019 }
1020 delete nodeList;
1021 return proxyList;
1022}
1023
1024// ***************************************************************************
1025// Calculate a "quality value" for the specified vector of voxels.
1026// The value returned should be >0 and such that the smaller the number
1027// the higher the quality of the slice.
1028//
1029// Preconditions: pSlice must consist of G4SmartVoxelNodeProxies only
1030// Process:
1031// o Examine each node in turn, summing:
1032// no. of non-empty nodes
1033// no. of volumes in each node
1034// o Calculate Quality=sigma(volumes in nod)/(no. of non-empty nodes)
1035// if all nodes empty, return kInfinity
1036// o Call G4Exception on finding a G4SmartVoxelHeaderProxy
1037// ***************************************************************************
1038//
1039G4double G4SmartVoxelHeader::CalculateQuality(G4ProxyVector *pSlice)
1040{
1041 G4double quality;
1042 std::size_t nNodes = pSlice->size();
1043 std::size_t noContained, maxContained=0, sumContained=0, sumNonEmptyNodes=0;
1044 G4SmartVoxelNode *node;
1045
1046 for (std::size_t i=0; i<nNodes; ++i)
1047 {
1048 if ((*pSlice)[i]->IsNode())
1049 {
1050 // Definitely a node. Add info to running totals
1051 //
1052 node = (*pSlice)[i]->GetNode();
1053 noContained = node->GetNoContained();
1054 if (noContained != 0)
1055 {
1056 ++sumNonEmptyNodes;
1057 sumContained += noContained;
1058 //
1059 // Calc maxContained for statistics
1060 //
1061 if (noContained>maxContained)
1062 {
1063 maxContained = noContained;
1064 }
1065 }
1066 }
1067 else
1068 {
1069 G4Exception("G4SmartVoxelHeader::CalculateQuality()", "GeomMgt0001",
1070 FatalException, "Not applicable to replicated volumes.");
1071 }
1072 }
1073
1074 // Calculate quality with protection against no non-empty nodes
1075 //
1076 if (sumNonEmptyNodes != 0)
1077 {
1078 quality = sumContained/sumNonEmptyNodes;
1079 }
1080 else
1081 {
1082 quality = kInfinity;
1083 }
1084
1085#ifdef G4GEOMETRY_VOXELDEBUG
1086 G4cout << "**** G4SmartVoxelHeader::CalculateQuality" << G4endl
1087 << " Quality = " << quality << G4endl
1088 << " Nodes = " << nNodes
1089 << " of which " << sumNonEmptyNodes << " non empty" << G4endl
1090 << " Max Contained = " << maxContained << G4endl;
1091#endif
1092
1093 return quality;
1094}
1095
1096// ***************************************************************************
1097// Examined each contained node, refines (creates a replacement additional
1098// dimension of voxels) when there is more than one voxel in the slice.
1099// Does not refine further if already limited in two dimensions (=> this
1100// is the third level of limits)
1101//
1102// Preconditions: slices (nodes) have been built.
1103// ***************************************************************************
1104//
1105void G4SmartVoxelHeader::RefineNodes(G4LogicalVolume* pVolume,
1106 G4VoxelLimits pLimits)
1107{
1108 std::size_t refinedDepth=0, minVolumes;
1109 std::size_t maxNode = fslices.size();
1110
1111 if (pLimits.IsXLimited())
1112 {
1113 ++refinedDepth;
1114 }
1115 if (pLimits.IsYLimited())
1116 {
1117 ++refinedDepth;
1118 }
1119 if (pLimits.IsZLimited())
1120 {
1121 ++refinedDepth;
1122 }
1123
1124 // Calculate minimum number of volumes necessary to refine
1125 //
1126 switch (refinedDepth)
1127 {
1128 case 0:
1129 minVolumes=kMinVoxelVolumesLevel2;
1130 break;
1131 case 1:
1132 minVolumes=kMinVoxelVolumesLevel3;
1133 break;
1134 default:
1135 minVolumes=10000; // catch refinedDepth=3 and errors
1136 break;
1137 }
1138
1139 if (refinedDepth<2)
1140 {
1141 std::size_t targetNo, noContainedDaughters, minNo, maxNo, replaceNo, i;
1142 G4double sliceWidth = (fmaxExtent-fminExtent)/maxNode;
1143 G4VoxelLimits newLimits;
1144 G4SmartVoxelNode* targetNode;
1145 G4SmartVoxelProxy* targetNodeProxy;
1146 G4SmartVoxelHeader* replaceHeader;
1147 G4SmartVoxelProxy* replaceHeaderProxy;
1148 G4VolumeNosVector* targetList;
1149 G4SmartVoxelProxy* lastProxy;
1150
1151 for (targetNo=0; targetNo<maxNode; ++targetNo)
1152 {
1153 // Assume all slices are nodes (see preconditions)
1154 //
1155 targetNodeProxy = fslices[targetNo];
1156 targetNode = targetNodeProxy->GetNode();
1157
1158 if (targetNode->GetNoContained() >= minVolumes)
1159 {
1160 noContainedDaughters = targetNode->GetNoContained();
1161 targetList = new G4VolumeNosVector();
1162 if (targetList == nullptr)
1163 {
1164 G4Exception("G4SmartVoxelHeader::RefineNodes()",
1165 "GeomMgt0003", FatalException,
1166 "Target volume node list allocation error.");
1167 return;
1168 }
1169 targetList->reserve(noContainedDaughters);
1170 for (i=0; i<noContainedDaughters; ++i)
1171 {
1172 targetList->push_back(targetNode->GetVolume((G4int)i));
1173 }
1174 minNo = targetNode->GetMinEquivalentSliceNo();
1175 maxNo = targetNode->GetMaxEquivalentSliceNo();
1176
1177#ifdef G4GEOMETRY_VOXELDEBUG
1178 G4cout << "**** G4SmartVoxelHeader::RefineNodes" << G4endl
1179 << " Refining nodes " << minNo
1180 << " - " << maxNo << " inclusive" << G4endl;
1181#endif
1182 if (minNo > maxNo) // Delete node and list to be replaced
1183 { // and avoid further action ...
1184 delete targetNode;
1185 delete targetList;
1186 return;
1187 }
1188
1189 // Delete node proxies at start of collected sets of nodes/headers
1190 //
1191 lastProxy=nullptr;
1192 for (replaceNo=minNo; replaceNo<=maxNo; ++replaceNo)
1193 {
1194 if (lastProxy != fslices[replaceNo])
1195 {
1196 lastProxy=fslices[replaceNo];
1197 delete lastProxy;
1198 }
1199 }
1200 // Delete node to be replaced
1201 //
1202 delete targetNode;
1203
1204 // Create new headers + proxies and replace in fslices
1205 //
1206 newLimits = pLimits;
1207 newLimits.AddLimit(faxis,fminExtent+sliceWidth*minNo,
1208 fminExtent+sliceWidth*(maxNo+1));
1209 replaceHeader = new G4SmartVoxelHeader(pVolume,newLimits,
1210 targetList,(G4int)replaceNo);
1211 if (replaceHeader == nullptr)
1212 {
1213 G4Exception("G4SmartVoxelHeader::RefineNodes()", "GeomMgt0003",
1214 FatalException, "Refined VoxelHeader allocation error.");
1215 return;
1216 }
1217 replaceHeader->SetMinEquivalentSliceNo((G4int)minNo);
1218 replaceHeader->SetMaxEquivalentSliceNo((G4int)maxNo);
1219 replaceHeaderProxy = new G4SmartVoxelProxy(replaceHeader);
1220 if (replaceHeaderProxy == nullptr)
1221 {
1222 G4Exception("G4SmartVoxelHeader::RefineNodes()", "GeomMgt0003",
1223 FatalException, "Refined VoxelProxy allocation error.");
1224 return;
1225 }
1226 for (replaceNo=minNo; replaceNo<=maxNo; ++replaceNo)
1227 {
1228 fslices[replaceNo] = replaceHeaderProxy;
1229 }
1230 // Finished replacing current `equivalent' group
1231 //
1232 delete targetList;
1233 targetNo=maxNo;
1234 }
1235 }
1236 }
1237}
1238
1239// ***************************************************************************
1240// Returns true if all slices have equal contents.
1241// Preconditions: all equal slices have been collected.
1242// Procedure:
1243// o checks all slice proxy pointers are equal
1244// o returns true if only one slice or all slice proxies pointers equal.
1245// ***************************************************************************
1246//
1248{
1249 std::size_t noSlices = fslices.size();
1250 G4SmartVoxelProxy* refProxy;
1251
1252 if (noSlices>1)
1253 {
1254 refProxy=fslices[0];
1255 for (std::size_t i=1; i<noSlices; ++i)
1256 {
1257 if (refProxy!=fslices[i])
1258 {
1259 return false;
1260 }
1261 }
1262 }
1263 return true;
1264}
1265
1266// ***************************************************************************
1267// Streaming operator for debugging.
1268// ***************************************************************************
1269//
1270std::ostream& operator << (std::ostream& os, const G4SmartVoxelHeader& h)
1271{
1272 os << "Axis = " << G4int(h.faxis) << G4endl;
1273 G4SmartVoxelProxy *collectNode=nullptr, *collectHead=nullptr;
1274 std::size_t collectNodeNo = 0;
1275 std::size_t collectHeadNo = 0;
1276 std::size_t i, j;
1277 G4bool haveHeaders = false;
1278
1279 for (i=0; i<h.fslices.size(); ++i)
1280 {
1281 os << "Slice #" << i << " = ";
1282 if (h.fslices[i]->IsNode())
1283 {
1284 if (h.fslices[i]!=collectNode)
1285 {
1286 os << "{";
1287 for (std::size_t k=0; k<h.fslices[i]->GetNode()->GetNoContained(); ++k)
1288 {
1289 os << " " << h.fslices[i]->GetNode()->GetVolume((G4int)k);
1290 }
1291 os << " }" << G4endl;
1292 collectNode = h.fslices[i];
1293 collectNodeNo = i;
1294 }
1295 else
1296 {
1297 os << "As slice #" << collectNodeNo << G4endl;
1298 }
1299 }
1300 else
1301 {
1302 haveHeaders=true;
1303 if (h.fslices[i] != collectHead)
1304 {
1305 os << "Header" << G4endl;
1306 collectHead = h.fslices[i];
1307 collectHeadNo = i;
1308 }
1309 else
1310 {
1311 os << "As slice #" << collectHeadNo << G4endl;
1312 }
1313 }
1314 }
1315
1316 if (haveHeaders)
1317 {
1318 collectHead=nullptr;
1319 for (j=0; j<h.fslices.size(); ++j)
1320 {
1321 if (h.fslices[j]->IsHeader())
1322 {
1323 os << "Header at Slice #" << j << " = ";
1324 if (h.fslices[j] != collectHead)
1325 {
1326 os << G4endl
1327 << (*(h.fslices[j]->GetHeader()));
1328 collectHead = h.fslices[j];
1329 collectHeadNo = j;
1330 }
1331 else
1332 {
1333 os << "As slice #" << collectHeadNo << G4endl;
1334 }
1335 }
1336 }
1337 }
1338 return os;
1339}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
G4ThreadLocal T * G4GeomSplitter< T >::offset
std::ostream & operator<<(std::ostream &os, const G4SmartVoxelHeader &h)
std::vector< G4SmartVoxelProxy * > G4ProxyVector
std::vector< G4SmartVoxelNode * > G4NodeVector
std::vector< G4double > G4VolumeExtentVector
std::vector< G4int > G4VolumeNosVector
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
G4LogicalVolume represents a leaf node or unpositioned subtree in the geometry hierarchy....
G4VSolid * GetSolid() const
std::size_t GetNoDaughters() const
G4double GetSmartless() const
G4VPhysicalVolume * GetDaughter(const std::size_t i) const
const G4String & GetName() const
G4int GetMaxEquivalentSliceNo() const
std::size_t GetNoSlices() const
G4double GetMaxExtent() const
void SetMinEquivalentSliceNo(G4int pMin)
G4double GetMinExtent() const
G4SmartVoxelProxy * GetSlice(std::size_t n) const
EAxis GetAxis() const
G4SmartVoxelHeader(G4LogicalVolume *pVolume, G4int pSlice=0)
G4bool operator==(const G4SmartVoxelHeader &pHead) const
G4bool AllSlicesEqual() const
void SetMaxEquivalentSliceNo(G4int pMax)
G4SmartVoxelNode defines a node in the smart voxel hierarchy, i.e. a 'slice' of space along a given a...
G4int GetMaxEquivalentSliceNo() const
G4int GetVolume(G4int pVolumeNo) const
void SetMaxEquivalentSliceNo(G4int pMax)
void SetMinEquivalentSliceNo(G4int pMin)
std::size_t GetNoContained() const
G4int GetMinEquivalentSliceNo() const
G4SmartVoxelProxy is a class for proxying smart voxels. The class represents either a header (in turn...
G4bool IsNode() const
G4SmartVoxelNode * GetNode() const
G4SmartVoxelHeader * GetHeader() const
G4bool IsHeader() const
virtual void ComputeTransformation(const G4int no, G4VPhysicalVolume *pv) const =0
virtual G4VSolid * ComputeSolid(const G4int no, G4VPhysicalVolume *pv)
virtual G4bool IsReplicated() const =0
const G4RotationMatrix * GetRotation() const
const G4ThreeVector GetTranslation() const
G4LogicalVolume * GetLogicalVolume() const
virtual void GetReplicationData(EAxis &axis, G4int &nReplicas, G4double &width, G4double &offset, G4bool &consuming) const =0
const G4String & GetName() const
virtual G4VPVParameterisation * GetParameterisation() const =0
virtual G4bool CalculateExtent(const EAxis pAxis, const G4VoxelLimits &pVoxelLimit, const G4AffineTransform &pTransform, G4double &pMin, G4double &pMax) const =0
virtual void ComputeDimensions(G4VPVParameterisation *p, const G4int n, const G4VPhysicalVolume *pRep)
Definition G4VSolid.cc:136
G4VoxelLimits represents limitation/restrictions of space, where restrictions are only made perpendic...
G4double GetMinExtent(const EAxis pAxis) const
G4bool IsYLimited() const
void AddLimit(const EAxis pAxis, const G4double pMin, const G4double pMax)
G4bool IsXLimited() const
G4double GetMaxExtent(const EAxis pAxis) const
G4bool IsZLimited() const
G4bool IsLimited() const
EAxis
Definition geomdefs.hh:54
@ kPhi
Definition geomdefs.hh:60
@ kYAxis
Definition geomdefs.hh:56
@ kXAxis
Definition geomdefs.hh:55
@ kZAxis
Definition geomdefs.hh:57
@ kUndefined
Definition geomdefs.hh:61
@ kRho
Definition geomdefs.hh:58
const axis_t axis_to_type< N >::axis
Definition pugixml.cc:9668
const G4int kMaxVoxelNodes
Definition voxeldefs.hh:36
const G4int kMinVoxelVolumesLevel3
Definition voxeldefs.hh:47
const G4int kMinVoxelVolumesLevel2
Definition voxeldefs.hh:43