ECCE @ EIC Software
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4PhysicalVolumeModel.cc
Go to the documentation of this file. Or view the newest version in sPHENIX GitHub for file G4PhysicalVolumeModel.cc
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 //
28 //
29 // John Allison 31st December 1997.
30 // Model for physical volumes.
31 
32 #include "G4PhysicalVolumeModel.hh"
33 
34 #include "G4VGraphicsScene.hh"
35 #include "G4VPhysicalVolume.hh"
36 #include "G4VPVParameterisation.hh"
37 #include "G4LogicalVolume.hh"
38 #include "G4VSolid.hh"
39 #include "G4SubtractionSolid.hh"
40 #include "G4IntersectionSolid.hh"
41 #include "G4Material.hh"
42 #include "G4VisAttributes.hh"
43 #include "G4BoundingExtentScene.hh"
46 #include "G4Polyhedron.hh"
47 #include "HepPolyhedronProcessor.h"
48 #include "G4AttDefStore.hh"
49 #include "G4AttDef.hh"
50 #include "G4AttValue.hh"
51 #include "G4UnitsTable.hh"
52 #include "G4Vector3D.hh"
53 
54 #include <sstream>
55 #include <iomanip>
56 
57 namespace {
58  G4int volumeCount = 0;
59 }
60 
63  , G4int requestedDepth
64  , const G4Transform3D& modelTransformation
65  , const G4ModelingParameters* pMP
66  , G4bool useFullExtent
67  , const std::vector<G4PhysicalVolumeNodeID>& baseFullPVPath)
68 : G4VModel (modelTransformation,pMP)
69 , fpTopPV (pVPV)
70 , fTopPVCopyNo (pVPV? pVPV->GetCopyNo(): 0)
71 , fRequestedDepth (requestedDepth)
72 , fUseFullExtent (useFullExtent)
73 , fCurrentDepth (0)
74 , fpCurrentPV (fpTopPV)
75 , fCurrentPVCopyNo (fpTopPV? fpTopPV->GetCopyNo(): 0)
76 , fpCurrentLV (fpTopPV? fpTopPV->GetLogicalVolume(): 0)
77 , fpCurrentMaterial (fpCurrentLV? fpCurrentLV->GetMaterial(): 0)
78 , fpCurrentTransform (const_cast<G4Transform3D*>(&modelTransformation))
79 , fBaseFullPVPath (baseFullPVPath)
80 , fAbort (false)
81 , fCurtailDescent (false)
82 , fpClippingSolid (0)
83 , fClippingMode (subtraction)
84 {
85  fType = "G4PhysicalVolumeModel";
86 
87  if (!fpTopPV) {
88 
89  // In some circumstances creating an "empty" G4PhysicalVolumeModel is
90  // allowed, so I have supressed the G4Exception below. If it proves to
91  // be a problem we might have to re-instate it, but it is unlikley to
92  // be used except by visualisation experts. See, for example, /vis/list,
93  // where it is used simply to get a list of G4AttDefs.
94  // G4Exception
95  // ("G4PhysicalVolumeModel::G4PhysicalVolumeModel",
96  // "modeling0010", FatalException, "Null G4PhysicalVolumeModel pointer.");
97 
98  fTopPVName = "NULL";
99  fGlobalTag = "Empty";
100  fGlobalDescription = "G4PhysicalVolumeModel " + fGlobalTag;
101 
102  } else {
103 
104  fTopPVName = fpTopPV -> GetName ();
105  std::ostringstream oss;
106  oss << fpTopPV->GetName() << ':' << fpTopPV->GetCopyNo()
107  << " BasePath:" << fBaseFullPVPath;
108  fGlobalTag = oss.str();
109  fGlobalDescription = "G4PhysicalVolumeModel " + fGlobalTag;
110  CalculateExtent ();
111  }
112 }
113 
115 {
116  delete fpClippingSolid;
117 }
118 
120 (const std::vector<G4PhysicalVolumeNodeID>& path)
121 {
123  for (const auto& node: path) {
124  PVNameCopyNoPath.push_back
126  (node.GetPhysicalVolume()->GetName(),node.GetCopyNo()));
127  }
128  return PVNameCopyNoPath;
129 }
130 
132 {
133  // To handle paramaterisations, set copy number and compute dimensions
134  // to get extent right
135  G4VPVParameterisation* pP = fpTopPV -> GetParameterisation ();
136  if (pP) {
137  fpTopPV -> SetCopyNo (fTopPVCopyNo);
138  G4VSolid* solid = pP -> ComputeSolid (fTopPVCopyNo, fpTopPV);
139  solid -> ComputeDimensions (pP, fTopPVCopyNo, fpTopPV);
140  }
141  if (fUseFullExtent) {
142  fExtent = fpTopPV -> GetLogicalVolume () -> GetSolid () -> GetExtent ();
143  } else {
144  // Calculate extent of *drawn* volumes, i.e., ignoring culled, e.g.,
145  // invisible volumes, by traversing the whole geometry hierarchy below
146  // this physical volume.
147  G4BoundingExtentScene beScene(this);
148  const G4int tempRequestedDepth = fRequestedDepth;
149  const G4Transform3D tempTransform = fTransform;
150  const G4ModelingParameters* tempMP = fpMP;
151  fRequestedDepth = -1; // Always search to all depths to define extent.
152  fTransform = G4Transform3D(); // Extent is in local cooridinates
153  G4ModelingParameters mParams
154  (0, // No default vis attributes needed.
155  G4ModelingParameters::wf, // wireframe (not relevant for this).
156  true, // Global culling.
157  true, // Cull invisible volumes.
158  false, // Density culling.
159  0., // Density (not relevant if density culling false).
160  true, // Cull daughters of opaque mothers.
161  24); // No of sides (not relevant for this operation).
162  fpMP = &mParams;
163  DescribeYourselfTo (beScene);
164  fExtent = beScene.GetBoundingExtent();
165  fpMP = tempMP;
166  fTransform = tempTransform;
167  fRequestedDepth = tempRequestedDepth;
168  }
170  if (radius < 0.) { // Nothing in the scene - revert to top extent
171  fExtent = fpTopPV -> GetLogicalVolume () -> GetSolid () -> GetExtent ();
172  }
173 }
174 
176 (G4VGraphicsScene& sceneHandler)
177 {
178  if (!fpTopPV) G4Exception
179  ("G4PhysicalVolumeModel::DescribeYourselfTo",
180  "modeling0012", FatalException, "No model.");
181 
182  if (!fpMP) G4Exception
183  ("G4PhysicalVolumeModel::DescribeYourselfTo",
184  "modeling0003", FatalException, "No modeling parameters.");
185 
186  G4Transform3D startingTransformation = fTransform;
187 
188  volumeCount = 0;
189 
190  VisitGeometryAndGetVisReps
191  (fpTopPV,
192  fRequestedDepth,
193  startingTransformation,
194  sceneHandler);
195 
196 // G4cout
197 // << "G4PhysicalVolumeModel::DescribeYourselfTo: volume count: "
198 // << volumeCount
199 // << G4endl;
200 
201  // Reset or clear data...
202  fCurrentDepth = 0;
203  fpCurrentPV = fpTopPV;
204  fCurrentPVCopyNo = fpTopPV->GetCopyNo();
205  fpCurrentLV = fpTopPV->GetLogicalVolume();
206  fpCurrentMaterial = fpCurrentLV? fpCurrentLV->GetMaterial(): 0;
207  fFullPVPath = fBaseFullPVPath;
208  fDrawnPVPath.clear();
209  fAbort = false;
210  fCurtailDescent = false;
211 }
212 
214 {
215  if (fpCurrentPV) {
216  std::ostringstream o;
217  o << fpCurrentPV -> GetCopyNo ();
218  return fpCurrentPV -> GetName () + "." + o.str();
219  }
220  else {
221  return "WARNING: NO CURRENT VOLUME - global tag is " + fGlobalTag;
222  }
223 }
224 
226 {
227  return "G4PhysicalVolumeModel " + GetCurrentTag ();
228 }
229 
232  G4int requestedDepth,
233  const G4Transform3D& theAT,
234  G4VGraphicsScene& sceneHandler)
235 {
236  // Visits geometry structure to a given depth (requestedDepth), starting
237  // at given physical volume with given starting transformation and
238  // describes volumes to the scene handler.
239  // requestedDepth < 0 (default) implies full visit.
240  // theAT is the Accumulated Transformation.
241 
242  // Find corresponding logical volume and (later) solid, storing in
243  // local variables to preserve re-entrancy.
244  G4LogicalVolume* pLV = pVPV -> GetLogicalVolume ();
245 
246  G4VSolid* pSol = nullptr;
247  G4Material* pMaterial = nullptr;
248 
249  if (!(pVPV -> IsReplicated ())) {
250  // Non-replicated physical volume.
251  pSol = pLV -> GetSolid ();
252  pMaterial = pLV -> GetMaterial ();
253  DescribeAndDescend (pVPV, requestedDepth, pLV, pSol, pMaterial,
254  theAT, sceneHandler);
255  }
256  else {
257  // Replicated or parametrised physical volume.
258  EAxis axis;
259  G4int nReplicas;
260  G4double width;
262  G4bool consuming;
263  pVPV -> GetReplicationData (axis, nReplicas, width, offset, consuming);
264  G4int nBegin = 0;
265  G4int nEnd = nReplicas;
266  if (fCurrentDepth == 0) { // i.e., top volume
267  nBegin = fTopPVCopyNo; // Describe only one volume, namely the one
268  nEnd = nBegin + 1; // specified by the given copy number.
269  }
270  G4VPVParameterisation* pP = pVPV -> GetParameterisation ();
271  if (pP) { // Parametrised volume.
272  for (int n = nBegin; n < nEnd; n++) {
273  pSol = pP -> ComputeSolid (n, pVPV);
274  pP -> ComputeTransformation (n, pVPV);
275  pSol -> ComputeDimensions (pP, n, pVPV);
276  pVPV -> SetCopyNo (n);
277  fCurrentPVCopyNo = n;
278  // Create a touchable of current parent for ComputeMaterial.
279  // fFullPVPath has not been updated yet so at this point it
280  // corresponds to the parent.
281  G4PhysicalVolumeModelTouchable parentTouchable(fFullPVPath);
282  pMaterial = pP -> ComputeMaterial (n, pVPV, &parentTouchable);
283  DescribeAndDescend (pVPV, requestedDepth, pLV, pSol, pMaterial,
284  theAT, sceneHandler);
285  }
286  }
287  else { // Plain replicated volume. From geometry_guide.txt...
288  // The replica's positions are claculated by means of a linear formula.
289  // Replication may occur along:
290  //
291  // o Cartesian axes (kXAxis,kYAxis,kZAxis)
292  //
293  // The replications, of specified width have coordinates of
294  // form (-width*(nReplicas-1)*0.5+n*width,0,0) where n=0.. nReplicas-1
295  // for the case of kXAxis, and are unrotated.
296  //
297  // o Radial axis (cylindrical polar) (kRho)
298  //
299  // The replications are cons/tubs sections, centred on the origin
300  // and are unrotated.
301  // They have radii of width*n+offset to width*(n+1)+offset
302  // where n=0..nReplicas-1
303  //
304  // o Phi axis (cylindrical polar) (kPhi)
305  // The replications are `phi sections' or wedges, and of cons/tubs form
306  // They have phi of offset+n*width to offset+(n+1)*width where
307  // n=0..nReplicas-1
308  //
309  pSol = pLV -> GetSolid ();
310  pMaterial = pLV -> GetMaterial ();
311  G4ThreeVector originalTranslation = pVPV -> GetTranslation ();
312  G4RotationMatrix* pOriginalRotation = pVPV -> GetRotation ();
313  G4double originalRMin = 0., originalRMax = 0.;
314  if (axis == kRho && pSol->GetEntityType() == "G4Tubs") {
315  originalRMin = ((G4Tubs*)pSol)->GetInnerRadius();
316  originalRMax = ((G4Tubs*)pSol)->GetOuterRadius();
317  }
318  G4bool visualisable = true;
319  for (int n = nBegin; n < nEnd; n++) {
320  G4ThreeVector translation; // Identity.
321  G4RotationMatrix rotation; // Identity - life enough for visualizing.
322  G4RotationMatrix* pRotation = 0;
323  switch (axis) {
324  default:
325  case kXAxis:
326  translation = G4ThreeVector (-width*(nReplicas-1)*0.5+n*width,0,0);
327  break;
328  case kYAxis:
329  translation = G4ThreeVector (0,-width*(nReplicas-1)*0.5+n*width,0);
330  break;
331  case kZAxis:
332  translation = G4ThreeVector (0,0,-width*(nReplicas-1)*0.5+n*width);
333  break;
334  case kRho:
335  if (pSol->GetEntityType() == "G4Tubs") {
336  ((G4Tubs*)pSol)->SetInnerRadius(width*n+offset);
337  ((G4Tubs*)pSol)->SetOuterRadius(width*(n+1)+offset);
338  } else {
339  if (fpMP->IsWarning())
340  G4cout <<
341  "G4PhysicalVolumeModel::VisitGeometryAndGetVisReps: WARNING:"
342  "\n built-in replicated volumes replicated in radius for "
343  << pSol->GetEntityType() <<
344  "-type\n solids (your solid \""
345  << pSol->GetName() <<
346  "\") are not visualisable."
347  << G4endl;
348  visualisable = false;
349  }
350  break;
351  case kPhi:
352  rotation.rotateZ (-(offset+(n+0.5)*width));
353  // Minus Sign because for the physical volume we need the
354  // coordinate system rotation.
355  pRotation = &rotation;
356  break;
357  }
358  pVPV -> SetTranslation (translation);
359  pVPV -> SetRotation (pRotation);
360  pVPV -> SetCopyNo (n);
361  fCurrentPVCopyNo = n;
362  if (visualisable) {
363  DescribeAndDescend (pVPV, requestedDepth, pLV, pSol, pMaterial,
364  theAT, sceneHandler);
365  }
366  }
367  // Restore originals...
368  pVPV -> SetTranslation (originalTranslation);
369  pVPV -> SetRotation (pOriginalRotation);
370  if (axis == kRho && pSol->GetEntityType() == "G4Tubs") {
371  ((G4Tubs*)pSol)->SetInnerRadius(originalRMin);
372  ((G4Tubs*)pSol)->SetOuterRadius(originalRMax);
373  }
374  }
375  }
376 }
377 
380  G4int requestedDepth,
381  G4LogicalVolume* pLV,
382  G4VSolid* pSol,
383  G4Material* pMaterial,
384  const G4Transform3D& theAT,
385  G4VGraphicsScene& sceneHandler)
386 {
387  // Maintain useful data members...
388  fpCurrentPV = pVPV;
389  fCurrentPVCopyNo = pVPV->GetCopyNo();
390  fpCurrentLV = pLV;
391  fpCurrentMaterial = pMaterial;
392 
393  const G4RotationMatrix objectRotation = pVPV -> GetObjectRotationValue ();
394  const G4ThreeVector& translation = pVPV -> GetTranslation ();
395  G4Transform3D theLT (G4Transform3D (objectRotation, translation));
396 
397  // Compute the accumulated transformation...
398  // Note that top volume's transformation relative to the world
399  // coordinate system is specified in theAT == startingTransformation
400  // = fTransform (see DescribeYourselfTo), so first time through the
401  // volume's own transformation, which is only relative to its
402  // mother, i.e., not relative to the world coordinate system, should
403  // not be accumulated.
404  G4Transform3D theNewAT (theAT);
405  if (fCurrentDepth != 0) theNewAT = theAT * theLT;
406  fpCurrentTransform = &theNewAT;
407 
408  const G4VisAttributes* pVisAttribs = pLV->GetVisAttributes();
409  // If the volume does not have any vis attributes, create it.
410  G4VisAttributes* tempVisAtts = nullptr;
411  if (!pVisAttribs) {
412  tempVisAtts = new G4VisAttributes; // Default value.
413  // The user may request /vis/viewer/set/colourByDensity.
414  if (fpMP->GetCBDAlgorithmNumber() == 1) {
415  // Algorithm 1: 3 parameters: Simple rainbow mapping.
416  if (fpMP->GetCBDParameters().size() != 3) {
417  G4Exception("G4PhysicalVolumeModelTouchable::DescribeAndDescend",
418  "modeling0014",
420  "Algorithm-parameter mismatch for Colour By Density");
421  } else {
422  const G4double d = pMaterial? pMaterial->GetDensity(): 0.;
423  const G4double d0 = fpMP->GetCBDParameters()[0]; // Invisible d < d0.
424  const G4double d1 = fpMP->GetCBDParameters()[1]; // Rainbow d0->d1->d2.
425  const G4double d2 = fpMP->GetCBDParameters()[2]; // Blue d > d2.
426  if (d < d0) { // Density < d0 is invisible.
427  tempVisAtts->SetVisibility(false);
428  } else { // Intermediate densities are on a spectrum.
429  G4double red, green, blue;
430  if (d < d1) {
431  red = (d1-d)/(d1-d0); green = (d-d0)/(d1-d0); blue = 0.;
432  } else if (d < d2) {
433  red = 0.; green = (d2-d)/(d2-d1); blue = (d-d1)/(d2-d1);
434  } else { // Density >= d2 is blue.
435  red = 0.; green = 0.; blue = 1.;
436  }
437  tempVisAtts->SetColour(G4Colour(red,green,blue));
438  }
439  }
440  } else if (fpMP->GetCBDAlgorithmNumber() == 2) {
441  // Algorithm 2
442  // ...etc.
443  }
444  pVisAttribs = tempVisAtts;
445  }
446  // From here, can assume pVisAttribs is a valid pointer. This is necessary
447  // because PreAddSolid needs a vis attributes object.
448 
449  // Make decision to draw...
450  G4bool thisToBeDrawn = true;
451 
452  // Update full path of physical volumes...
453  G4int copyNo = fpCurrentPV->GetCopyNo();
454  fFullPVPath.push_back
456  (fpCurrentPV,copyNo,fCurrentDepth,*fpCurrentTransform));
457 
458  // Check if vis attributes are to be modified by a /vis/touchable/set/ command.
459  const auto& vams = fpMP->GetVisAttributesModifiers();
460  if (vams.size()) {
461  // OK, we have some VAMs (Vis Attributes Modifiers).
462  for (const auto& vam: vams) {
463  const auto& vamPath = vam.GetPVNameCopyNoPath();
464  if (vamPath.size() == fFullPVPath.size()) {
465  // OK, we have a size match.
466  // Check the volume name/copy number path.
467  auto iVAMNameCopyNo = vamPath.begin();
468  auto iPVNodeId = fFullPVPath.begin();
469  for (; iVAMNameCopyNo != vamPath.end(); ++iVAMNameCopyNo, ++iPVNodeId) {
470  if (!(
471  iVAMNameCopyNo->GetName() ==
472  iPVNodeId->GetPhysicalVolume()->GetName() &&
473  iVAMNameCopyNo->GetCopyNo() ==
474  iPVNodeId->GetPhysicalVolume()->GetCopyNo()
475  )) {
476  // This path element does NOT match.
477  break;
478  }
479  }
480  if (iVAMNameCopyNo == vamPath.end()) {
481  // OK, the paths match (the above loop terminated normally).
482  // Create a vis atts object for the modified vis atts.
483  // It is static so that we may return a reliable pointer to it.
484  static G4VisAttributes modifiedVisAtts;
485  // Initialise it with the current vis atts and reset the pointer.
486  modifiedVisAtts = *pVisAttribs;
487  pVisAttribs = &modifiedVisAtts;
488  const G4VisAttributes& transVisAtts = vam.GetVisAttributes();
489  switch (vam.GetVisAttributesSignifier()) {
491  modifiedVisAtts.SetVisibility(transVisAtts.IsVisible());
492  break;
494  modifiedVisAtts.SetDaughtersInvisible
495  (transVisAtts.IsDaughtersInvisible());
496  break;
498  modifiedVisAtts.SetColour(transVisAtts.GetColour());
499  break;
501  modifiedVisAtts.SetLineStyle(transVisAtts.GetLineStyle());
502  break;
504  modifiedVisAtts.SetLineWidth(transVisAtts.GetLineWidth());
505  break;
507  if (transVisAtts.IsForceDrawingStyle()) {
508  if (transVisAtts.GetForcedDrawingStyle() ==
510  modifiedVisAtts.SetForceWireframe(true);
511  }
512  }
513  break;
515  if (transVisAtts.IsForceDrawingStyle()) {
516  if (transVisAtts.GetForcedDrawingStyle() ==
518  modifiedVisAtts.SetForceSolid(true);
519  }
520  }
521  break;
523  if (transVisAtts.IsForceDrawingStyle()) {
524  if (transVisAtts.GetForcedDrawingStyle() ==
526  modifiedVisAtts.SetForceCloud(true);
527  }
528  }
529  break;
531  modifiedVisAtts.SetForceNumberOfCloudPoints
532  (transVisAtts.GetForcedNumberOfCloudPoints());
533  break;
535  if (transVisAtts.IsForceAuxEdgeVisible()) {
536  modifiedVisAtts.SetForceAuxEdgeVisible
537  (transVisAtts.IsForcedAuxEdgeVisible());
538  }
539  break;
541  modifiedVisAtts.SetForceLineSegmentsPerCircle
542  (transVisAtts.GetForcedLineSegmentsPerCircle());
543  break;
544  }
545  }
546  }
547  }
548  }
549 
550  // There are various reasons why this volume
551  // might not be drawn...
552  G4bool culling = fpMP->IsCulling();
553  G4bool cullingInvisible = fpMP->IsCullingInvisible();
554  G4bool markedVisible = pVisAttribs->IsVisible();
555  G4bool cullingLowDensity = fpMP->IsDensityCulling();
556  G4double density = pMaterial? pMaterial->GetDensity(): 0;
557  G4double densityCut = fpMP -> GetVisibleDensity ();
558 
559  // 1) Global culling is on....
560  if (culling) {
561  // 2) Culling of invisible volumes is on...
562  if (cullingInvisible) {
563  // 3) ...and the volume is marked not visible...
564  if (!markedVisible) thisToBeDrawn = false;
565  }
566  // 4) Or culling of low density volumes is on...
567  if (cullingLowDensity) {
568  // 5) ...and density is less than cut value...
569  if (density < densityCut) thisToBeDrawn = false;
570  }
571  }
572  // 6) The user has asked for all further traversing to be aborted...
573  if (fAbort) thisToBeDrawn = false;
574 
575  // Record thisToBeDrawn in path...
576  fFullPVPath.back().SetDrawn(thisToBeDrawn);
577 
578  if (thisToBeDrawn) {
579 
580  // Update path of drawn physical volumes...
581  fDrawnPVPath.push_back
583  (fpCurrentPV,copyNo,fCurrentDepth,*fpCurrentTransform,thisToBeDrawn));
584 
585  if (fpMP->IsExplode() && fDrawnPVPath.size() == 1) {
586  // For top-level drawn volumes, explode along radius...
587  G4Transform3D centering = G4Translate3D(fpMP->GetExplodeCentre());
588  G4Transform3D centred = centering.inverse() * theNewAT;
589  G4Scale3D oldScale;
590  G4Rotate3D oldRotation;
591  G4Translate3D oldTranslation;
592  centred.getDecomposition(oldScale, oldRotation, oldTranslation);
593  G4double explodeFactor = fpMP->GetExplodeFactor();
594  G4Translate3D newTranslation =
595  G4Translate3D(explodeFactor * oldTranslation.dx(),
596  explodeFactor * oldTranslation.dy(),
597  explodeFactor * oldTranslation.dz());
598  theNewAT = centering * newTranslation * oldRotation * oldScale;
599  }
600 
601  volumeCount++;
602  DescribeSolid (theNewAT, pSol, pVisAttribs, sceneHandler);
603 
604  }
605 
606  // Make decision to draw daughters, if any. There are various
607  // reasons why daughters might not be drawn...
608 
609  // First, reasons that do not depend on culling policy...
610  G4int nDaughters = pLV->GetNoDaughters();
611  G4bool daughtersToBeDrawn = true;
612  // 1) There are no daughters...
613  if (!nDaughters) daughtersToBeDrawn = false;
614  // 2) We are at the limit if requested depth...
615  else if (requestedDepth == 0) daughtersToBeDrawn = false;
616  // 3) The user has asked for all further traversing to be aborted...
617  else if (fAbort) daughtersToBeDrawn = false;
618  // 4) The user has asked that the descent be curtailed...
619  else if (fCurtailDescent) daughtersToBeDrawn = false;
620 
621  // Now, reasons that depend on culling policy...
622  else {
623  G4bool daughtersInvisible = pVisAttribs->IsDaughtersInvisible();
624  // Culling of covered daughters request. This is computed in
625  // G4VSceneHandler::CreateModelingParameters() depending on view
626  // parameters...
627  G4bool cullingCovered = fpMP->IsCullingCovered();
628  G4bool surfaceDrawing =
629  fpMP->GetDrawingStyle() == G4ModelingParameters::hsr ||
630  fpMP->GetDrawingStyle() == G4ModelingParameters::hlhsr;
631  if (pVisAttribs->IsForceDrawingStyle()) {
632  switch (pVisAttribs->GetForcedDrawingStyle()) {
633  default:
634  case G4VisAttributes::wireframe: surfaceDrawing = false; break;
635  case G4VisAttributes::solid: surfaceDrawing = true; break;
636  }
637  }
638  G4bool opaque = pVisAttribs->GetColour().GetAlpha() >= 1.;
639  // 5) Global culling is on....
640  if (culling) {
641  // 6) ..and culling of invisible volumes is on...
642  if (cullingInvisible) {
643  // 7) ...and the mother requests daughters invisible
644  if (daughtersInvisible) daughtersToBeDrawn = false;
645  }
646  // 8) Or culling of covered daughters is requested...
647  if (cullingCovered) {
648  // 9) ...and surface drawing is operating...
649  if (surfaceDrawing) {
650  // 10) ...but only if mother is visible...
651  if (thisToBeDrawn) {
652  // 11) ...and opaque...
653  if (opaque) daughtersToBeDrawn = false;
654  }
655  }
656  }
657  }
658  }
659 
660  if (daughtersToBeDrawn) {
661  for (G4int iDaughter = 0; iDaughter < nDaughters; iDaughter++) {
662  // Store daughter pVPV in local variable ready for recursion...
663  G4VPhysicalVolume* pDaughterVPV = pLV -> GetDaughter (iDaughter);
664  // Descend the geometry structure recursively...
665  fCurrentDepth++;
666  VisitGeometryAndGetVisReps
667  (pDaughterVPV, requestedDepth - 1, theNewAT, sceneHandler);
668  fCurrentDepth--;
669  }
670  }
671 
672  delete tempVisAtts;
673 
674  // Reset for normal descending of next volume at this level...
675  fCurtailDescent = false;
676 
677  // Pop item from paths physical volumes...
678  fFullPVPath.pop_back();
679  if (thisToBeDrawn) {
680  fDrawnPVPath.pop_back();
681  }
682 }
683 
685 (const G4Transform3D& theAT,
686  G4VSolid* pSol,
687  const G4VisAttributes* pVisAttribs,
688  G4VGraphicsScene& sceneHandler)
689 {
690  G4DisplacedSolid* pSectionSolid = fpMP->GetSectionSolid();
691  G4DisplacedSolid* pCutawaySolid = fpMP->GetCutawaySolid();
692 
693  if (!fpClippingSolid && !pSectionSolid && !pCutawaySolid) {
694 
695  sceneHandler.PreAddSolid (theAT, *pVisAttribs);
696  pSol -> DescribeYourselfTo (sceneHandler); // Standard treatment.
697  sceneHandler.PostAddSolid ();
698 
699  } else {
700 
701  // Clipping, etc., performed by Boolean operations.
702 
703  // First, get polyhedron for current solid...
704  if (pVisAttribs->IsForceLineSegmentsPerCircle())
706  (pVisAttribs->GetForcedLineSegmentsPerCircle());
707  else
708  G4Polyhedron::SetNumberOfRotationSteps(fpMP->GetNoOfSides());
709  const G4Polyhedron* pOriginalPolyhedron = pSol->GetPolyhedron();
711 
712  if (!pOriginalPolyhedron) {
713 
714  if (fpMP->IsWarning())
715  G4cout <<
716  "WARNING: G4PhysicalVolumeModel::DescribeSolid: solid\n \""
717  << pSol->GetName() <<
718  "\" has no polyhedron. Cannot by clipped."
719  << G4endl;
720  pSol -> DescribeYourselfTo (sceneHandler); // Standard treatment.
721 
722  } else {
723 
724  G4VSolid* pResultantSolid = 0;
725 
726  if (fpClippingSolid) {
727  switch (fClippingMode) {
728  default:
729  case subtraction:
730  pResultantSolid = new G4SubtractionSolid
731  ("subtracted_clipped_solid", pSol, fpClippingSolid, theAT.inverse());
732  break;
733  case intersection:
734  pResultantSolid = new G4IntersectionSolid
735  ("intersected_clipped_solid", pSol, fpClippingSolid, theAT.inverse());
736  break;
737  }
738  }
739 
740  if (pSectionSolid) {
741  pResultantSolid = new G4IntersectionSolid
742  ("sectioned_solid", pSol, pSectionSolid, theAT.inverse());
743  }
744 
745  if (pCutawaySolid) {
746  // Follow above...
747  pResultantSolid = new G4SubtractionSolid
748  ("cutaway_solid", pSol, pCutawaySolid, theAT.inverse());
749  }
750 
751  const G4Polyhedron* pResultantPolyhedron = pResultantSolid->GetPolyhedron();
752  if (!pResultantPolyhedron) {
753  if (fpMP->IsWarning())
754  G4cout <<
755  "WARNING: G4PhysicalVolumeModel::DescribeSolid: resultant polyhedron for"
756  "\n solid \"" << pSol->GetName() <<
757  "\" not defined due to error during Boolean processing."
758  << G4endl;
759  } else {
760  // It seems that if the sectioning solid does not intersect the
761  // original solid the Boolean Processor returns the original
762  // polyhedron, or a copy thereof. We do not want it.
763  // Check the number of facets, etc. If same, ignore.
764  // What we need from the Boolean Processor is a null pointer or a
765  // null polyhedron. It seems to return the original or a copy of it.
766  if (pResultantPolyhedron->GetNoFacets() == pOriginalPolyhedron->GetNoFacets())
767  // This works in most cases but I still get a box in test202 with
768  // /vis/viewer/set/sectionPlane on 0 0 0 m 0.1 0.1 1
769  {
770  pResultantPolyhedron = nullptr;
771  }
772  }
773 
774  if (pResultantPolyhedron) {
775  // Finally, draw polyhedron...
776  sceneHandler.BeginPrimitives(theAT);
777  sceneHandler.AddPrimitive(*pResultantPolyhedron);
778  sceneHandler.EndPrimitives();
779  }
780 
781  delete pResultantSolid;
782  }
783  }
784 }
785 
787 {
788  G4TransportationManager* transportationManager =
790 
791  size_t nWorlds = transportationManager->GetNoWorlds();
792 
793  G4bool found = false;
794 
795  std::vector<G4VPhysicalVolume*>::iterator iterWorld =
796  transportationManager->GetWorldsIterator();
797  for (size_t i = 0; i < nWorlds; ++i, ++iterWorld) {
798  G4VPhysicalVolume* world = (*iterWorld);
799  if (!world) break; // This can happen if geometry has been cleared/destroyed.
800  // The idea now is to seek a PV with the same name and copy no
801  // in the hope it's the same one!!
802  G4PhysicalVolumeModel searchModel (world);
803  G4int verbosity = 0; // Suppress messages from G4PhysicalVolumeSearchScene.
804  G4PhysicalVolumeSearchScene searchScene
805  (&searchModel, fTopPVName, fTopPVCopyNo, verbosity);
806  G4ModelingParameters mp; // Default modeling parameters for this search.
808  searchModel.SetModelingParameters (&mp);
809  searchModel.DescribeYourselfTo (searchScene);
810  G4VPhysicalVolume* foundVolume = searchScene.GetFoundVolume ();
811  if (foundVolume) {
812  if (foundVolume != fpTopPV && warn) {
813  G4cout <<
814  "G4PhysicalVolumeModel::Validate(): A volume of the same name and"
815  "\n copy number (\""
816  << fTopPVName << "\", copy " << fTopPVCopyNo
817  << ") still exists and is being used."
818  "\n But it is not the same volume you originally specified"
819  "\n in /vis/scene/add/."
820  << G4endl;
821  }
822  fpTopPV = foundVolume;
823  CalculateExtent ();
824  found = true;
825  }
826  }
827  if (found) return true;
828  else {
829  if (warn) {
830  G4cout <<
831  "G4PhysicalVolumeModel::Validate(): No volume of name and"
832  "\n copy number (\""
833  << fTopPVName << "\", copy " << fTopPVCopyNo
834  << ") exists."
835  << G4endl;
836  }
837  return false;
838  }
839 }
840 
841 const std::map<G4String,G4AttDef>* G4PhysicalVolumeModel::GetAttDefs() const
842 {
843  G4bool isNew;
844  std::map<G4String,G4AttDef>* store
845  = G4AttDefStore::GetInstance("G4PhysicalVolumeModel", isNew);
846  if (isNew) {
847  (*store)["PVPath"] =
848  G4AttDef("PVPath","Physical Volume Path","Physics","","G4String");
849  (*store)["BasePVPath"] =
850  G4AttDef("BasePVPath","Base Physical Volume Path","Physics","","G4String");
851  (*store)["LVol"] =
852  G4AttDef("LVol","Logical Volume","Physics","","G4String");
853  (*store)["Solid"] =
854  G4AttDef("Solid","Solid Name","Physics","","G4String");
855  (*store)["EType"] =
856  G4AttDef("EType","Entity Type","Physics","","G4String");
857  (*store)["DmpSol"] =
858  G4AttDef("DmpSol","Dump of Solid properties","Physics","","G4String");
859  (*store)["LocalTrans"] =
860  G4AttDef("LocalTrans","Local transformation of volume","Physics","","G4String");
861  (*store)["GlobalTrans"] =
862  G4AttDef("GlobalTrans","Global transformation of volume","Physics","","G4String");
863  (*store)["Material"] =
864  G4AttDef("Material","Material Name","Physics","","G4String");
865  (*store)["Density"] =
866  G4AttDef("Density","Material Density","Physics","G4BestUnit","G4double");
867  (*store)["State"] =
868  G4AttDef("State","Material State (enum undefined,solid,liquid,gas)","Physics","","G4String");
869  (*store)["Radlen"] =
870  G4AttDef("Radlen","Material Radiation Length","Physics","G4BestUnit","G4double");
871  (*store)["Region"] =
872  G4AttDef("Region","Cuts Region","Physics","","G4String");
873  (*store)["RootRegion"] =
874  G4AttDef("RootRegion","Root Region (0/1 = false/true)","Physics","","G4bool");
875  }
876  return store;
877 }
878 
879 static std::ostream& operator<< (std::ostream& o, const G4Transform3D t)
880 {
881  using namespace std;
882 
883  G4Scale3D sc;
884  G4Rotate3D r;
885  G4Translate3D tl;
886  t.getDecomposition(sc, r, tl);
887 
888  const int w = 10;
889 
890  // Transformation itself
891  o << setw(w) << t.xx() << setw(w) << t.xy() << setw(w) << t.xz() << setw(w) << t.dx() << endl;
892  o << setw(w) << t.yx() << setw(w) << t.yy() << setw(w) << t.yz() << setw(w) << t.dy() << endl;
893  o << setw(w) << t.zx() << setw(w) << t.zy() << setw(w) << t.zz() << setw(w) << t.dz() << endl;
894 
895  // Translation
896  o << "= translation:" << endl;
897  o << setw(w) << tl.dx() << setw(w) << tl.dy() << setw(w) << tl.dz() << endl;
898 
899  // Rotation
900  o << "* rotation:" << endl;
901  o << setw(w) << r.xx() << setw(w) << r.xy() << setw(w) << r.xz() << endl;
902  o << setw(w) << r.yx() << setw(w) << r.yy() << setw(w) << r.yz() << endl;
903  o << setw(w) << r.zx() << setw(w) << r.zy() << setw(w) << r.zz() << endl;
904 
905  // Scale
906  o << "* scale:" << endl;
907  o << setw(w) << sc.xx() << setw(w) << sc.yy() << setw(w) << sc.zz() << endl;
908 
909  // Transformed axes
910  o << "Transformed axes:" << endl;
911  o << "x': " << r * G4Vector3D(1., 0., 0.) << endl;
912  o << "y': " << r * G4Vector3D(0., 1., 0.) << endl;
913  o << "z': " << r * G4Vector3D(0., 0., 1.) << endl;
914 
915  return o;
916 }
917 
918 std::vector<G4AttValue>* G4PhysicalVolumeModel::CreateCurrentAttValues() const
919 {
920  std::vector<G4AttValue>* values = new std::vector<G4AttValue>;
921 
922  if (!fpCurrentLV) {
924  ("G4PhysicalVolumeModel::CreateCurrentAttValues",
925  "modeling0004",
926  JustWarning,
927  "Current logical volume not defined.");
928  return values;
929  }
930 
931  std::ostringstream oss; oss << fFullPVPath;
932  values->push_back(G4AttValue("PVPath", oss.str(),""));
933  oss.str(""); oss << fBaseFullPVPath;
934  values->push_back(G4AttValue("BasePVPath", oss.str(),""));
935  values->push_back(G4AttValue("LVol", fpCurrentLV->GetName(),""));
936  G4VSolid* pSol = fpCurrentLV->GetSolid();
937  values->push_back(G4AttValue("Solid", pSol->GetName(),""));
938  values->push_back(G4AttValue("EType", pSol->GetEntityType(),""));
939  oss.str(""); oss << '\n' << *pSol;
940  values->push_back(G4AttValue("DmpSol", oss.str(),""));
941  const G4RotationMatrix localRotation = fpCurrentPV->GetObjectRotationValue();
942  const G4ThreeVector& localTranslation = fpCurrentPV->GetTranslation();
943  oss.str(""); oss << '\n' << G4Transform3D(localRotation,localTranslation);
944  values->push_back(G4AttValue("LocalTrans", oss.str(),""));
945  oss.str(""); oss << '\n' << *fpCurrentTransform;
946  values->push_back(G4AttValue("GlobalTrans", oss.str(),""));
947  G4String matName = fpCurrentMaterial? fpCurrentMaterial->GetName(): G4String("No material");
948  values->push_back(G4AttValue("Material", matName,""));
950  values->push_back(G4AttValue("Density", G4BestUnit(matDensity,"Volumic Mass"),""));
952  oss.str(""); oss << matState;
953  values->push_back(G4AttValue("State", oss.str(),""));
955  values->push_back(G4AttValue("Radlen", G4BestUnit(matRadlen,"Length"),""));
956  G4Region* region = fpCurrentLV->GetRegion();
957  G4String regionName = region? region->GetName(): G4String("No region");
958  values->push_back(G4AttValue("Region", regionName,""));
959  oss.str(""); oss << fpCurrentLV->IsRootRegion();
960  values->push_back(G4AttValue("RootRegion", oss.str(),""));
961  return values;
962 }
963 
964 G4bool G4PhysicalVolumeModel::G4PhysicalVolumeNodeID::operator<
966 {
967  if (fpPV < right.fpPV) return true;
968  if (fpPV == right.fpPV) {
969  if (fCopyNo < right.fCopyNo) return true;
970  if (fCopyNo == right.fCopyNo)
971  return fNonCulledDepth < right.fNonCulledDepth;
972  }
973  return false;
974 }
975 
976 std::ostream& operator<<
977  (std::ostream& os, const G4PhysicalVolumeModel::G4PhysicalVolumeNodeID& node)
978 {
979  G4VPhysicalVolume* pPV = node.GetPhysicalVolume();
980  if (pPV) {
981  os << pPV->GetName()
982  << ' ' << node.GetCopyNo()
983 // << '[' << node.GetNonCulledDepth() << ']'
984 // << ':' << node.GetTransform()
985  ;
986 // os << " (";
987 // if (!node.GetDrawn()) os << "not ";
988 // os << "drawn)";
989  } else {
990  os << " (Null node)";
991  }
992  return os;
993 }
994 
995 std::ostream& operator<<
996 (std::ostream& os, const std::vector<G4PhysicalVolumeModel::G4PhysicalVolumeNodeID>& path)
997 {
998  if (path.empty()) {
999  os << " TOP";
1000  } else {
1001  for (const auto& nodeID: path) {
1002  os << ' ' << nodeID;
1003  }
1004  }
1005  return os;
1006 }
1007 
1009 (const std::vector<G4PhysicalVolumeNodeID>& fullPVPath):
1010  fFullPVPath(fullPVPath) {}
1011 
1013 {
1014  size_t i = fFullPVPath.size() - depth - 1;
1015  if (i >= fFullPVPath.size()) {
1016  G4Exception("G4PhysicalVolumeModelTouchable::GetTranslation",
1017  "modeling0005",
1019  "Index out of range. Asking for non-existent depth");
1020  }
1021  static G4ThreeVector tempTranslation;
1022  tempTranslation = fFullPVPath[i].GetTransform().getTranslation();
1023  return tempTranslation;
1024 }
1025 
1027 {
1028  size_t i = fFullPVPath.size() - depth - 1;
1029  if (i >= fFullPVPath.size()) {
1030  G4Exception("G4PhysicalVolumeModelTouchable::GetRotation",
1031  "modeling0006",
1033  "Index out of range. Asking for non-existent depth");
1034  }
1035  static G4RotationMatrix tempRotation;
1036  tempRotation = fFullPVPath[i].GetTransform().getRotation();
1037  return &tempRotation;
1038 }
1039 
1041 {
1042  size_t i = fFullPVPath.size() - depth - 1;
1043  if (i >= fFullPVPath.size()) {
1044  G4Exception("G4PhysicalVolumeModelTouchable::GetVolume",
1045  "modeling0007",
1047  "Index out of range. Asking for non-existent depth");
1048  }
1049  return fFullPVPath[i].GetPhysicalVolume();
1050 }
1051 
1053 {
1054  size_t i = fFullPVPath.size() - depth - 1;
1055  if (i >= fFullPVPath.size()) {
1056  G4Exception("G4PhysicalVolumeModelTouchable::GetSolid",
1057  "modeling0008",
1059  "Index out of range. Asking for non-existent depth");
1060  }
1061  return fFullPVPath[i].GetPhysicalVolume()->GetLogicalVolume()->GetSolid();
1062 }
1063 
1065 {
1066  size_t i = fFullPVPath.size() - depth - 1;
1067  if (i >= fFullPVPath.size()) {
1068  G4Exception("G4PhysicalVolumeModelTouchable::GetReplicaNumber",
1069  "modeling0009",
1071  "Index out of range. Asking for non-existent depth");
1072  }
1073  return fFullPVPath[i].GetCopyNo();
1074 }