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