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