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