Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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"
37#include "G4LogicalVolume.hh"
38#include "G4VSolid.hh"
39#include "G4SubtractionSolid.hh"
41#include "G4Material.hh"
42#include "G4VisAttributes.hh"
46#include "G4Polyhedron.hh"
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
57namespace {
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;
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
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
191 (fpTopPV,
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;
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
231(G4VPhysicalVolume* pVPV,
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;
261 G4double offset;
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);
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.
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);
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
379(G4VPhysicalVolume* pVPV,
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) {
413 tempVisAtts = new G4VisAttributes(*fpMP->GetDefaultVisAttributes());
414 } else {
415 tempVisAtts = new G4VisAttributes;
416 }
417 // The user may request /vis/viewer/set/colourByDensity.
418 if (fpMP->GetCBDAlgorithmNumber() == 1) {
419 // Algorithm 1: 3 parameters: Simple rainbow mapping.
420 if (fpMP->GetCBDParameters().size() != 3) {
421 G4Exception("G4PhysicalVolumeModelTouchable::DescribeAndDescend",
422 "modeling0014",
424 "Algorithm-parameter mismatch for Colour By Density");
425 } else {
426 const G4double d = pMaterial? pMaterial->GetDensity(): 0.;
427 const G4double d0 = fpMP->GetCBDParameters()[0]; // Invisible d < d0.
428 const G4double d1 = fpMP->GetCBDParameters()[1]; // Rainbow d0->d1->d2.
429 const G4double d2 = fpMP->GetCBDParameters()[2]; // Blue d > d2.
430 if (d < d0) { // Density < d0 is invisible.
431 tempVisAtts->SetVisibility(false);
432 } else { // Intermediate densities are on a spectrum.
433 G4double red, green, blue;
434 if (d < d1) {
435 red = (d1-d)/(d1-d0); green = (d-d0)/(d1-d0); blue = 0.;
436 } else if (d < d2) {
437 red = 0.; green = (d2-d)/(d2-d1); blue = (d-d1)/(d2-d1);
438 } else { // Density >= d2 is blue.
439 red = 0.; green = 0.; blue = 1.;
440 }
441 tempVisAtts->SetColour(G4Colour(red,green,blue));
442 }
443 }
444 } else if (fpMP->GetCBDAlgorithmNumber() == 2) {
445 // Algorithm 2
446 // ...etc.
447 }
448 pVisAttribs = tempVisAtts;
449 }
450 // From here, can assume pVisAttribs is a valid pointer. This is necessary
451 // because PreAddSolid needs a vis attributes object.
452
453 // Make decision to draw...
454 G4bool thisToBeDrawn = true;
455
456 // Update full path of physical volumes...
457 G4int copyNo = fpCurrentPV->GetCopyNo();
458 fFullPVPath.push_back
461
462 // Check if vis attributes are to be modified by a /vis/touchable/set/ command.
463 const auto& vams = fpMP->GetVisAttributesModifiers();
464 if (vams.size()) {
465 // OK, we have some VAMs (Vis Attributes Modifiers).
466 for (const auto& vam: vams) {
467 const auto& vamPath = vam.GetPVNameCopyNoPath();
468 if (vamPath.size() == fFullPVPath.size()) {
469 // OK, we have a size match.
470 // Check the volume name/copy number path.
471 auto iVAMNameCopyNo = vamPath.begin();
472 auto iPVNodeId = fFullPVPath.begin();
473 for (; iVAMNameCopyNo != vamPath.end(); ++iVAMNameCopyNo, ++iPVNodeId) {
474 if (!(
475 iVAMNameCopyNo->GetName() ==
476 iPVNodeId->GetPhysicalVolume()->GetName() &&
477 iVAMNameCopyNo->GetCopyNo() ==
478 iPVNodeId->GetPhysicalVolume()->GetCopyNo()
479 )) {
480 // This path element does NOT match.
481 break;
482 }
483 }
484 if (iVAMNameCopyNo == vamPath.end()) {
485 // OK, the paths match (the above loop terminated normally).
486 // Create a vis atts object for the modified vis atts.
487 // It is static so that we may return a reliable pointer to it.
488 static G4VisAttributes modifiedVisAtts;
489 // Initialise it with the current vis atts and reset the pointer.
490 modifiedVisAtts = *pVisAttribs;
491 pVisAttribs = &modifiedVisAtts;
492 const G4VisAttributes& transVisAtts = vam.GetVisAttributes();
493 switch (vam.GetVisAttributesSignifier()) {
495 modifiedVisAtts.SetVisibility(transVisAtts.IsVisible());
496 break;
498 modifiedVisAtts.SetDaughtersInvisible
499 (transVisAtts.IsDaughtersInvisible());
500 break;
502 modifiedVisAtts.SetColour(transVisAtts.GetColour());
503 break;
505 modifiedVisAtts.SetLineStyle(transVisAtts.GetLineStyle());
506 break;
508 modifiedVisAtts.SetLineWidth(transVisAtts.GetLineWidth());
509 break;
511 if (transVisAtts.IsForceDrawingStyle()) {
512 if (transVisAtts.GetForcedDrawingStyle() ==
514 modifiedVisAtts.SetForceWireframe(true);
515 }
516 }
517 break;
519 if (transVisAtts.IsForceDrawingStyle()) {
520 if (transVisAtts.GetForcedDrawingStyle() ==
522 modifiedVisAtts.SetForceSolid(true);
523 }
524 }
525 break;
527 if (transVisAtts.IsForceDrawingStyle()) {
528 if (transVisAtts.GetForcedDrawingStyle() ==
530 modifiedVisAtts.SetForceCloud(true);
531 }
532 }
533 break;
535 modifiedVisAtts.SetForceNumberOfCloudPoints
536 (transVisAtts.GetForcedNumberOfCloudPoints());
537 break;
539 if (transVisAtts.IsForceAuxEdgeVisible()) {
540 modifiedVisAtts.SetForceAuxEdgeVisible
541 (transVisAtts.IsForcedAuxEdgeVisible());
542 }
543 break;
545 modifiedVisAtts.SetForceLineSegmentsPerCircle
546 (transVisAtts.GetForcedLineSegmentsPerCircle());
547 break;
548 }
549 }
550 }
551 }
552 }
553
554 // There are various reasons why this volume
555 // might not be drawn...
556 G4bool culling = fpMP->IsCulling();
557 G4bool cullingInvisible = fpMP->IsCullingInvisible();
558 G4bool markedVisible = pVisAttribs->IsVisible();
559 G4bool cullingLowDensity = fpMP->IsDensityCulling();
560 G4double density = pMaterial? pMaterial->GetDensity(): 0;
561 G4double densityCut = fpMP -> GetVisibleDensity ();
562
563 // 1) Global culling is on....
564 if (culling) {
565 // 2) Culling of invisible volumes is on...
566 if (cullingInvisible) {
567 // 3) ...and the volume is marked not visible...
568 if (!markedVisible) thisToBeDrawn = false;
569 }
570 // 4) Or culling of low density volumes is on...
571 if (cullingLowDensity) {
572 // 5) ...and density is less than cut value...
573 if (density < densityCut) thisToBeDrawn = false;
574 }
575 }
576 // 6) The user has asked for all further traversing to be aborted...
577 if (fAbort) thisToBeDrawn = false;
578
579 // Record thisToBeDrawn in path...
580 fFullPVPath.back().SetDrawn(thisToBeDrawn);
581
582 if (thisToBeDrawn) {
583
584 // Update path of drawn physical volumes...
585 fDrawnPVPath.push_back
587 (fpCurrentPV,copyNo,fCurrentDepth,*fpCurrentTransform,thisToBeDrawn));
588
589 if (fpMP->IsExplode() && fDrawnPVPath.size() == 1) {
590 // For top-level drawn volumes, explode along radius...
592 G4Transform3D centred = centering.inverse() * theNewAT;
593 G4Scale3D oldScale;
594 G4Rotate3D oldRotation;
595 G4Translate3D oldTranslation;
596 centred.getDecomposition(oldScale, oldRotation, oldTranslation);
597 G4double explodeFactor = fpMP->GetExplodeFactor();
598 G4Translate3D newTranslation =
599 G4Translate3D(explodeFactor * oldTranslation.dx(),
600 explodeFactor * oldTranslation.dy(),
601 explodeFactor * oldTranslation.dz());
602 theNewAT = centering * newTranslation * oldRotation * oldScale;
603 }
604
605 volumeCount++;
606 DescribeSolid (theNewAT, pSol, pVisAttribs, sceneHandler);
607
608 }
609
610 // Make decision to draw daughters, if any. There are various
611 // reasons why daughters might not be drawn...
612
613 // First, reasons that do not depend on culling policy...
614 G4int nDaughters = pLV->GetNoDaughters();
615 G4bool daughtersToBeDrawn = true;
616 // 1) There are no daughters...
617 if (!nDaughters) daughtersToBeDrawn = false;
618 // 2) We are at the limit if requested depth...
619 else if (requestedDepth == 0) daughtersToBeDrawn = false;
620 // 3) The user has asked for all further traversing to be aborted...
621 else if (fAbort) daughtersToBeDrawn = false;
622 // 4) The user has asked that the descent be curtailed...
623 else if (fCurtailDescent) daughtersToBeDrawn = false;
624
625 // Now, reasons that depend on culling policy...
626 else {
627 G4bool daughtersInvisible = pVisAttribs->IsDaughtersInvisible();
628 // Culling of covered daughters request. This is computed in
629 // G4VSceneHandler::CreateModelingParameters() depending on view
630 // parameters...
631 G4bool cullingCovered = fpMP->IsCullingCovered();
632 G4bool surfaceDrawing =
635 if (pVisAttribs->IsForceDrawingStyle()) {
636 switch (pVisAttribs->GetForcedDrawingStyle()) {
637 default:
638 case G4VisAttributes::wireframe: surfaceDrawing = false; break;
639 case G4VisAttributes::solid: surfaceDrawing = true; break;
640 }
641 }
642 G4bool opaque = pVisAttribs->GetColour().GetAlpha() >= 1.;
643 // 5) Global culling is on....
644 if (culling) {
645 // 6) ..and culling of invisible volumes is on...
646 if (cullingInvisible) {
647 // 7) ...and the mother requests daughters invisible
648 if (daughtersInvisible) daughtersToBeDrawn = false;
649 }
650 // 8) Or culling of covered daughters is requested...
651 if (cullingCovered) {
652 // 9) ...and surface drawing is operating...
653 if (surfaceDrawing) {
654 // 10) ...but only if mother is visible...
655 if (thisToBeDrawn) {
656 // 11) ...and opaque...
657 if (opaque) daughtersToBeDrawn = false;
658 }
659 }
660 }
661 }
662 }
663
664 if (daughtersToBeDrawn) {
665 for (G4int iDaughter = 0; iDaughter < nDaughters; iDaughter++) {
666 // Store daughter pVPV in local variable ready for recursion...
667 G4VPhysicalVolume* pDaughterVPV = pLV -> GetDaughter (iDaughter);
668 // Descend the geometry structure recursively...
671 (pDaughterVPV, requestedDepth - 1, theNewAT, sceneHandler);
673 }
674 }
675
676 delete tempVisAtts;
677
678 // Reset for normal descending of next volume at this level...
679 fCurtailDescent = false;
680
681 // Pop item from paths physical volumes...
682 fFullPVPath.pop_back();
683 if (thisToBeDrawn) {
684 fDrawnPVPath.pop_back();
685 }
686}
687
689(const G4Transform3D& theAT,
690 G4VSolid* pSol,
691 const G4VisAttributes* pVisAttribs,
692 G4VGraphicsScene& sceneHandler)
693{
694 G4DisplacedSolid* pSectionSolid = fpMP->GetSectionSolid();
695 G4DisplacedSolid* pCutawaySolid = fpMP->GetCutawaySolid();
696
697 if (!fpClippingSolid && !pSectionSolid && !pCutawaySolid) {
698
699 sceneHandler.PreAddSolid (theAT, *pVisAttribs);
700 pSol -> DescribeYourselfTo (sceneHandler); // Standard treatment.
701 sceneHandler.PostAddSolid ();
702
703 } else {
704
705 // Clipping, etc., performed by Boolean operations.
706
707 // First, get polyhedron for current solid...
708 if (pVisAttribs->IsForceLineSegmentsPerCircle())
710 (pVisAttribs->GetForcedLineSegmentsPerCircle());
711 else
713 const G4Polyhedron* pOriginalPolyhedron = pSol->GetPolyhedron();
715
716 if (!pOriginalPolyhedron) {
717
718 if (fpMP->IsWarning())
719 G4cout <<
720 "WARNING: G4PhysicalVolumeModel::DescribeSolid: solid\n \""
721 << pSol->GetName() <<
722 "\" has no polyhedron. Cannot by clipped."
723 << G4endl;
724 pSol -> DescribeYourselfTo (sceneHandler); // Standard treatment.
725
726 } else {
727
728 G4VSolid* pResultantSolid = 0;
729
730 if (fpClippingSolid) {
731 switch (fClippingMode) {
732 default:
733 case subtraction:
734 pResultantSolid = new G4SubtractionSolid
735 ("subtracted_clipped_solid", pSol, fpClippingSolid, theAT.inverse());
736 break;
737 case intersection:
738 pResultantSolid = new G4IntersectionSolid
739 ("intersected_clipped_solid", pSol, fpClippingSolid, theAT.inverse());
740 break;
741 }
742 }
743
744 if (pSectionSolid) {
745 pResultantSolid = new G4IntersectionSolid
746 ("sectioned_solid", pSol, pSectionSolid, theAT.inverse());
747 }
748
749 if (pCutawaySolid) {
750 // Follow above...
751 pResultantSolid = new G4SubtractionSolid
752 ("cutaway_solid", pSol, pCutawaySolid, theAT.inverse());
753 }
754
755 const G4Polyhedron* pResultantPolyhedron = pResultantSolid->GetPolyhedron();
756 if (!pResultantPolyhedron) {
757 if (fpMP->IsWarning())
758 G4cout <<
759 "WARNING: G4PhysicalVolumeModel::DescribeSolid: resultant polyhedron for"
760 "\n solid \"" << pSol->GetName() <<
761 "\" not defined due to error during Boolean processing."
762 << G4endl;
763 } else {
764 // It seems that if the sectioning solid does not intersect the
765 // original solid the Boolean Processor returns the original
766 // polyhedron, or a copy thereof. We do not want it.
767 // Check the number of facets, etc. If same, ignore.
768 // What we need from the Boolean Processor is a null pointer or a
769 // null polyhedron. It seems to return the original or a copy of it.
770 if (pResultantPolyhedron->GetNoFacets() == pOriginalPolyhedron->GetNoFacets())
771 // This works in most cases but I still get a box in test202 with
772 // /vis/viewer/set/sectionPlane on 0 0 0 m 0.1 0.1 1
773 {
774 pResultantPolyhedron = nullptr;
775 }
776 }
777
778 if (pResultantPolyhedron) {
779 // Finally, draw polyhedron...
780 sceneHandler.BeginPrimitives(theAT);
781 sceneHandler.AddPrimitive(*pResultantPolyhedron);
782 sceneHandler.EndPrimitives();
783 }
784
785 delete pResultantSolid;
786 }
787 }
788}
789
791{
792 G4TransportationManager* transportationManager =
794
795 size_t nWorlds = transportationManager->GetNoWorlds();
796
797 G4bool found = false;
798
799 std::vector<G4VPhysicalVolume*>::iterator iterWorld =
800 transportationManager->GetWorldsIterator();
801 for (size_t i = 0; i < nWorlds; ++i, ++iterWorld) {
802 G4VPhysicalVolume* world = (*iterWorld);
803 if (!world) break; // This can happen if geometry has been cleared/destroyed.
804 // The idea now is to seek a PV with the same name and copy no
805 // in the hope it's the same one!!
806 G4PhysicalVolumeModel searchModel (world);
807 G4int verbosity = 0; // Suppress messages from G4PhysicalVolumeSearchScene.
809 (&searchModel, fTopPVName, fTopPVCopyNo, verbosity);
810 G4ModelingParameters mp; // Default modeling parameters for this search.
812 searchModel.SetModelingParameters (&mp);
813 searchModel.DescribeYourselfTo (searchScene);
814 G4VPhysicalVolume* foundVolume = searchScene.GetFoundVolume ();
815 if (foundVolume) {
816 if (foundVolume != fpTopPV && warn) {
817 G4cout <<
818 "G4PhysicalVolumeModel::Validate(): A volume of the same name and"
819 "\n copy number (\""
820 << fTopPVName << "\", copy " << fTopPVCopyNo
821 << ") still exists and is being used."
822 "\n But it is not the same volume you originally specified"
823 "\n in /vis/scene/add/."
824 << G4endl;
825 }
826 fpTopPV = foundVolume;
828 found = true;
829 }
830 }
831 if (found) return true;
832 else {
833 if (warn) {
834 G4cout <<
835 "G4PhysicalVolumeModel::Validate(): No volume of name and"
836 "\n copy number (\""
837 << fTopPVName << "\", copy " << fTopPVCopyNo
838 << ") exists."
839 << G4endl;
840 }
841 return false;
842 }
843}
844
845const std::map<G4String,G4AttDef>* G4PhysicalVolumeModel::GetAttDefs() const
846{
847 G4bool isNew;
848 std::map<G4String,G4AttDef>* store
849 = G4AttDefStore::GetInstance("G4PhysicalVolumeModel", isNew);
850 if (isNew) {
851 (*store)["PVPath"] =
852 G4AttDef("PVPath","Physical Volume Path","Physics","","G4String");
853 (*store)["BasePVPath"] =
854 G4AttDef("BasePVPath","Base Physical Volume Path","Physics","","G4String");
855 (*store)["LVol"] =
856 G4AttDef("LVol","Logical Volume","Physics","","G4String");
857 (*store)["Solid"] =
858 G4AttDef("Solid","Solid Name","Physics","","G4String");
859 (*store)["EType"] =
860 G4AttDef("EType","Entity Type","Physics","","G4String");
861 (*store)["DmpSol"] =
862 G4AttDef("DmpSol","Dump of Solid properties","Physics","","G4String");
863 (*store)["LocalTrans"] =
864 G4AttDef("LocalTrans","Local transformation of volume","Physics","","G4String");
865 (*store)["GlobalTrans"] =
866 G4AttDef("GlobalTrans","Global transformation of volume","Physics","","G4String");
867 (*store)["Material"] =
868 G4AttDef("Material","Material Name","Physics","","G4String");
869 (*store)["Density"] =
870 G4AttDef("Density","Material Density","Physics","G4BestUnit","G4double");
871 (*store)["State"] =
872 G4AttDef("State","Material State (enum undefined,solid,liquid,gas)","Physics","","G4String");
873 (*store)["Radlen"] =
874 G4AttDef("Radlen","Material Radiation Length","Physics","G4BestUnit","G4double");
875 (*store)["Region"] =
876 G4AttDef("Region","Cuts Region","Physics","","G4String");
877 (*store)["RootRegion"] =
878 G4AttDef("RootRegion","Root Region (0/1 = false/true)","Physics","","G4bool");
879 }
880 return store;
881}
882
883static std::ostream& operator<< (std::ostream& o, const G4Transform3D t)
884{
885 using namespace std;
886
887 G4Scale3D sc;
888 G4Rotate3D r;
889 G4Translate3D tl;
890 t.getDecomposition(sc, r, tl);
891
892 const int w = 10;
893
894 // Transformation itself
895 o << setw(w) << t.xx() << setw(w) << t.xy() << setw(w) << t.xz() << setw(w) << t.dx() << endl;
896 o << setw(w) << t.yx() << setw(w) << t.yy() << setw(w) << t.yz() << setw(w) << t.dy() << endl;
897 o << setw(w) << t.zx() << setw(w) << t.zy() << setw(w) << t.zz() << setw(w) << t.dz() << endl;
898
899 // Translation
900 o << "= translation:" << endl;
901 o << setw(w) << tl.dx() << setw(w) << tl.dy() << setw(w) << tl.dz() << endl;
902
903 // Rotation
904 o << "* rotation:" << endl;
905 o << setw(w) << r.xx() << setw(w) << r.xy() << setw(w) << r.xz() << endl;
906 o << setw(w) << r.yx() << setw(w) << r.yy() << setw(w) << r.yz() << endl;
907 o << setw(w) << r.zx() << setw(w) << r.zy() << setw(w) << r.zz() << endl;
908
909 // Scale
910 o << "* scale:" << endl;
911 o << setw(w) << sc.xx() << setw(w) << sc.yy() << setw(w) << sc.zz() << endl;
912
913 // Transformed axes
914 o << "Transformed axes:" << endl;
915 o << "x': " << r * G4Vector3D(1., 0., 0.) << endl;
916 o << "y': " << r * G4Vector3D(0., 1., 0.) << endl;
917 o << "z': " << r * G4Vector3D(0., 0., 1.) << endl;
918
919 return o;
920}
921
922std::vector<G4AttValue>* G4PhysicalVolumeModel::CreateCurrentAttValues() const
923{
924 std::vector<G4AttValue>* values = new std::vector<G4AttValue>;
925
926 if (!fpCurrentLV) {
928 ("G4PhysicalVolumeModel::CreateCurrentAttValues",
929 "modeling0004",
931 "Current logical volume not defined.");
932 return values;
933 }
934
935 std::ostringstream oss; oss << fFullPVPath;
936 values->push_back(G4AttValue("PVPath", oss.str(),""));
937 oss.str(""); oss << fBaseFullPVPath;
938 values->push_back(G4AttValue("BasePVPath", oss.str(),""));
939 values->push_back(G4AttValue("LVol", fpCurrentLV->GetName(),""));
940 G4VSolid* pSol = fpCurrentLV->GetSolid();
941 values->push_back(G4AttValue("Solid", pSol->GetName(),""));
942 values->push_back(G4AttValue("EType", pSol->GetEntityType(),""));
943 oss.str(""); oss << '\n' << *pSol;
944 values->push_back(G4AttValue("DmpSol", oss.str(),""));
945 const G4RotationMatrix localRotation = fpCurrentPV->GetObjectRotationValue();
946 const G4ThreeVector& localTranslation = fpCurrentPV->GetTranslation();
947 oss.str(""); oss << '\n' << G4Transform3D(localRotation,localTranslation);
948 values->push_back(G4AttValue("LocalTrans", oss.str(),""));
949 oss.str(""); oss << '\n' << *fpCurrentTransform;
950 values->push_back(G4AttValue("GlobalTrans", oss.str(),""));
951 G4String matName = fpCurrentMaterial? fpCurrentMaterial->GetName(): G4String("No material");
952 values->push_back(G4AttValue("Material", matName,""));
954 values->push_back(G4AttValue("Density", G4BestUnit(matDensity,"Volumic Mass"),""));
956 oss.str(""); oss << matState;
957 values->push_back(G4AttValue("State", oss.str(),""));
959 values->push_back(G4AttValue("Radlen", G4BestUnit(matRadlen,"Length"),""));
960 G4Region* region = fpCurrentLV->GetRegion();
961 G4String regionName = region? region->GetName(): G4String("No region");
962 values->push_back(G4AttValue("Region", regionName,""));
963 oss.str(""); oss << fpCurrentLV->IsRootRegion();
964 values->push_back(G4AttValue("RootRegion", oss.str(),""));
965 return values;
966}
967
968G4bool G4PhysicalVolumeModel::G4PhysicalVolumeNodeID::operator<
970{
971 if (fpPV < right.fpPV) return true;
972 if (fpPV == right.fpPV) {
973 if (fCopyNo < right.fCopyNo) return true;
974 if (fCopyNo == right.fCopyNo)
975 return fNonCulledDepth < right.fNonCulledDepth;
976 }
977 return false;
978}
979
980G4bool G4PhysicalVolumeModel::G4PhysicalVolumeNodeID::operator!=
982{
983 if (fpPV != right.fpPV ||
984 fCopyNo != right.fCopyNo ||
985 fNonCulledDepth != right.fNonCulledDepth ||
986 fTransform != right.fTransform ||
987 fDrawn != right.fDrawn) return true;
988 return false;
989}
990
991std::ostream& operator<<
992 (std::ostream& os, const G4PhysicalVolumeModel::G4PhysicalVolumeNodeID& node)
993{
994 G4VPhysicalVolume* pPV = node.GetPhysicalVolume();
995 if (pPV) {
996 os << pPV->GetName()
997 << ' ' << node.GetCopyNo()
998// << '[' << node.GetNonCulledDepth() << ']'
999// << ':' << node.GetTransform()
1000 ;
1001// os << " (";
1002// if (!node.GetDrawn()) os << "not ";
1003// os << "drawn)";
1004 } else {
1005 os << " (Null PV node)";
1006 }
1007 return os;
1008}
1009
1010std::ostream& operator<<
1011(std::ostream& os, const std::vector<G4PhysicalVolumeModel::G4PhysicalVolumeNodeID>& path)
1012{
1013 if (path.empty()) {
1014 os << " TOP";
1015 } else {
1016 for (const auto& nodeID: path) {
1017 os << ' ' << nodeID;
1018 }
1019 }
1020 return os;
1021}
1022
1024(const std::vector<G4PhysicalVolumeNodeID>& fullPVPath):
1025 fFullPVPath(fullPVPath) {}
1026
1028{
1029 size_t i = fFullPVPath.size() - depth - 1;
1030 if (i >= fFullPVPath.size()) {
1031 G4Exception("G4PhysicalVolumeModelTouchable::GetTranslation",
1032 "modeling0005",
1034 "Index out of range. Asking for non-existent depth");
1035 }
1036 static G4ThreeVector tempTranslation;
1037 tempTranslation = fFullPVPath[i].GetTransform().getTranslation();
1038 return tempTranslation;
1039}
1040
1042{
1043 size_t i = fFullPVPath.size() - depth - 1;
1044 if (i >= fFullPVPath.size()) {
1045 G4Exception("G4PhysicalVolumeModelTouchable::GetRotation",
1046 "modeling0006",
1048 "Index out of range. Asking for non-existent depth");
1049 }
1050 static G4RotationMatrix tempRotation;
1051 tempRotation = fFullPVPath[i].GetTransform().getRotation();
1052 return &tempRotation;
1053}
1054
1056{
1057 size_t i = fFullPVPath.size() - depth - 1;
1058 if (i >= fFullPVPath.size()) {
1059 G4Exception("G4PhysicalVolumeModelTouchable::GetVolume",
1060 "modeling0007",
1062 "Index out of range. Asking for non-existent depth");
1063 }
1064 return fFullPVPath[i].GetPhysicalVolume();
1065}
1066
1068{
1069 size_t i = fFullPVPath.size() - depth - 1;
1070 if (i >= fFullPVPath.size()) {
1071 G4Exception("G4PhysicalVolumeModelTouchable::GetSolid",
1072 "modeling0008",
1074 "Index out of range. Asking for non-existent depth");
1075 }
1076 return fFullPVPath[i].GetPhysicalVolume()->GetLogicalVolume()->GetSolid();
1077}
1078
1080{
1081 size_t i = fFullPVPath.size() - depth - 1;
1082 if (i >= fFullPVPath.size()) {
1083 G4Exception("G4PhysicalVolumeModelTouchable::GetReplicaNumber",
1084 "modeling0009",
1086 "Index out of range. Asking for non-existent depth");
1087 }
1088 return fFullPVPath[i].GetCopyNo();
1089}
@ JustWarning
@ FatalException
@ FatalErrorInArgument
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
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:153
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:178
G4State GetState() const
Definition: G4Material.hh:179
G4double GetRadlen() const
Definition: G4Material.hh:218
const G4String & GetName() const
Definition: G4Material.hh:175
G4bool IsWarning() const
const G4VisAttributes * GetDefaultVisAttributes() const
const std::vector< VisAttributesModifier > & GetVisAttributesModifiers() const
const G4Point3D & GetExplodeCentre() const
std::vector< PVNameCopyNo > PVNameCopyNoPath
G4bool IsCullingInvisible() const
G4bool IsExplode() const
G4bool IsCulling() const
const std::vector< G4double > & GetCBDParameters() const
G4int GetNoOfSides() const
G4bool IsDensityCulling() const
void SetDefaultVisAttributes(const G4VisAttributes *pDefaultVisAttributes)
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 >())
G4Transform3D * fpCurrentTransform
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 > &)
G4VPhysicalVolume * GetFoundVolume() const
const G4String & GetName() const
static G4TransportationManager * GetTransportationManager()
std::vector< G4VPhysicalVolume * >::iterator GetWorldsIterator()
size_t GetNoWorlds() const
Definition: G4Tubs.hh:75
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:113
void SetModelingParameters(const G4ModelingParameters *)
G4String fGlobalDescription
Definition: G4VModel.hh:112
const G4VisExtent & GetExtent() const
G4String fType
Definition: G4VModel.hh:110
const G4ModelingParameters * fpMP
Definition: G4VModel.hh:115
G4Transform3D fTransform
Definition: G4VModel.hh:114
G4String fGlobalTag
Definition: G4VModel.hh:111
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:693
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
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)