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