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,v 1.66 2009/10/23 08:08:19 allison Exp $ |
---|
28 | // GEANT4 tag $Name: $ |
---|
29 | // |
---|
30 | // |
---|
31 | // John Allison 31st December 1997. |
---|
32 | // Model for physical volumes. |
---|
33 | |
---|
34 | #include "G4PhysicalVolumeModel.hh" |
---|
35 | |
---|
36 | #include "G4ModelingParameters.hh" |
---|
37 | #include "G4VGraphicsScene.hh" |
---|
38 | #include "G4VPhysicalVolume.hh" |
---|
39 | #include "G4VPVParameterisation.hh" |
---|
40 | #include "G4LogicalVolume.hh" |
---|
41 | #include "G4VSolid.hh" |
---|
42 | #include "G4Material.hh" |
---|
43 | #include "G4VisAttributes.hh" |
---|
44 | #include "G4BoundingSphereScene.hh" |
---|
45 | #include "G4PhysicalVolumeSearchScene.hh" |
---|
46 | #include "G4TransportationManager.hh" |
---|
47 | #include "G4Polyhedron.hh" |
---|
48 | #include "HepPolyhedronProcessor.h" |
---|
49 | #include "G4AttDefStore.hh" |
---|
50 | #include "G4AttDef.hh" |
---|
51 | #include "G4AttValue.hh" |
---|
52 | #include "G4UnitsTable.hh" |
---|
53 | #include "G4Vector3D.hh" |
---|
54 | |
---|
55 | #include <sstream> |
---|
56 | |
---|
57 | G4PhysicalVolumeModel::G4PhysicalVolumeModel |
---|
58 | (G4VPhysicalVolume* pVPV, |
---|
59 | G4int requestedDepth, |
---|
60 | const G4Transform3D& modelTransformation, |
---|
61 | const G4ModelingParameters* pMP, |
---|
62 | G4bool useFullExtent): |
---|
63 | G4VModel (modelTransformation, pMP), |
---|
64 | fpTopPV (pVPV), |
---|
65 | fTopPVName (pVPV -> GetName ()), |
---|
66 | fTopPVCopyNo (pVPV -> GetCopyNo ()), |
---|
67 | fRequestedDepth (requestedDepth), |
---|
68 | fUseFullExtent (useFullExtent), |
---|
69 | fCurrentDepth (0), |
---|
70 | fpCurrentPV (0), |
---|
71 | fpCurrentLV (0), |
---|
72 | fpCurrentMaterial (0), |
---|
73 | fpCurrentTransform (0), |
---|
74 | fCurtailDescent (false), |
---|
75 | fpClippingPolyhedron (0), |
---|
76 | fClippingMode (subtraction) |
---|
77 | { |
---|
78 | std::ostringstream o; |
---|
79 | o << fpTopPV -> GetCopyNo (); |
---|
80 | fGlobalTag = fpTopPV -> GetName () + "." + o.str(); |
---|
81 | fGlobalDescription = "G4PhysicalVolumeModel " + fGlobalTag; |
---|
82 | |
---|
83 | CalculateExtent (); |
---|
84 | } |
---|
85 | |
---|
86 | G4PhysicalVolumeModel::~G4PhysicalVolumeModel () |
---|
87 | { |
---|
88 | delete fpClippingPolyhedron; |
---|
89 | } |
---|
90 | |
---|
91 | void G4PhysicalVolumeModel::CalculateExtent () |
---|
92 | { |
---|
93 | if (fUseFullExtent) { |
---|
94 | fExtent = fpTopPV -> GetLogicalVolume () -> GetSolid () -> GetExtent (); |
---|
95 | } |
---|
96 | else { |
---|
97 | G4BoundingSphereScene bsScene(this); |
---|
98 | const G4int tempRequestedDepth = fRequestedDepth; |
---|
99 | fRequestedDepth = -1; // Always search to all depths to define extent. |
---|
100 | const G4ModelingParameters* tempMP = fpMP; |
---|
101 | G4ModelingParameters mParams |
---|
102 | (0, // No default vis attributes needed. |
---|
103 | G4ModelingParameters::wf, // wireframe (not relevant for this). |
---|
104 | true, // Global culling. |
---|
105 | true, // Cull invisible volumes. |
---|
106 | false, // Density culling. |
---|
107 | 0., // Density (not relevant if density culling false). |
---|
108 | true, // Cull daughters of opaque mothers. |
---|
109 | 24); // No of sides (not relevant for this operation). |
---|
110 | fpMP = &mParams; |
---|
111 | DescribeYourselfTo (bsScene); |
---|
112 | G4double radius = bsScene.GetRadius(); |
---|
113 | if (radius < 0.) { // Nothing in the scene. |
---|
114 | fExtent = fpTopPV -> GetLogicalVolume () -> GetSolid () -> GetExtent (); |
---|
115 | } else { |
---|
116 | // Transform back to coordinates relative to the top |
---|
117 | // transformation, which is in G4VModel::fTransform. This makes |
---|
118 | // it conform to all models, which are defined by a |
---|
119 | // transformation and an extent relative to that |
---|
120 | // transformation... |
---|
121 | G4Point3D centre = bsScene.GetCentre(); |
---|
122 | centre.transform(fTransform.inverse()); |
---|
123 | fExtent = G4VisExtent(centre, radius); |
---|
124 | } |
---|
125 | fpMP = tempMP; |
---|
126 | fRequestedDepth = tempRequestedDepth; |
---|
127 | } |
---|
128 | } |
---|
129 | |
---|
130 | void G4PhysicalVolumeModel::DescribeYourselfTo |
---|
131 | (G4VGraphicsScene& sceneHandler) |
---|
132 | { |
---|
133 | if (!fpMP) G4Exception |
---|
134 | ("G4PhysicalVolumeModel::DescribeYourselfTo: No modeling parameters."); |
---|
135 | |
---|
136 | // For safety... |
---|
137 | fCurrentDepth = 0; |
---|
138 | |
---|
139 | G4Transform3D startingTransformation = fTransform; |
---|
140 | |
---|
141 | VisitGeometryAndGetVisReps |
---|
142 | (fpTopPV, |
---|
143 | fRequestedDepth, |
---|
144 | startingTransformation, |
---|
145 | sceneHandler); |
---|
146 | |
---|
147 | // Clear data... |
---|
148 | fCurrentDepth = 0; |
---|
149 | fpCurrentPV = 0; |
---|
150 | fpCurrentLV = 0; |
---|
151 | fpCurrentMaterial = 0; |
---|
152 | fFullPVPath.clear(); |
---|
153 | fDrawnPVPath.clear(); |
---|
154 | } |
---|
155 | |
---|
156 | G4String G4PhysicalVolumeModel::GetCurrentTag () const |
---|
157 | { |
---|
158 | if (fpCurrentPV) { |
---|
159 | std::ostringstream o; |
---|
160 | o << fpCurrentPV -> GetCopyNo (); |
---|
161 | return fpCurrentPV -> GetName () + "." + o.str(); |
---|
162 | } |
---|
163 | else { |
---|
164 | return "WARNING: NO CURRENT VOLUME - global tag is " + fGlobalTag; |
---|
165 | } |
---|
166 | } |
---|
167 | |
---|
168 | G4String G4PhysicalVolumeModel::GetCurrentDescription () const |
---|
169 | { |
---|
170 | return "G4PhysicalVolumeModel " + GetCurrentTag (); |
---|
171 | } |
---|
172 | |
---|
173 | void G4PhysicalVolumeModel::VisitGeometryAndGetVisReps |
---|
174 | (G4VPhysicalVolume* pVPV, |
---|
175 | G4int requestedDepth, |
---|
176 | const G4Transform3D& theAT, |
---|
177 | G4VGraphicsScene& sceneHandler) |
---|
178 | { |
---|
179 | // Visits geometry structure to a given depth (requestedDepth), starting |
---|
180 | // at given physical volume with given starting transformation and |
---|
181 | // describes volumes to the scene handler. |
---|
182 | // requestedDepth < 0 (default) implies full visit. |
---|
183 | // theAT is the Accumulated Transformation. |
---|
184 | |
---|
185 | // Find corresponding logical volume and (later) solid, storing in |
---|
186 | // local variables to preserve re-entrancy. |
---|
187 | G4LogicalVolume* pLV = pVPV -> GetLogicalVolume (); |
---|
188 | |
---|
189 | G4VSolid* pSol; |
---|
190 | G4Material* pMaterial; |
---|
191 | |
---|
192 | if (!(pVPV -> IsReplicated ())) { |
---|
193 | // Non-replicated physical volume. |
---|
194 | pSol = pLV -> GetSolid (); |
---|
195 | pMaterial = pLV -> GetMaterial (); |
---|
196 | DescribeAndDescend (pVPV, requestedDepth, pLV, pSol, pMaterial, |
---|
197 | theAT, sceneHandler); |
---|
198 | } |
---|
199 | else { |
---|
200 | // Replicated or parametrised physical volume. |
---|
201 | EAxis axis; |
---|
202 | G4int nReplicas; |
---|
203 | G4double width; |
---|
204 | G4double offset; |
---|
205 | G4bool consuming; |
---|
206 | pVPV -> GetReplicationData (axis, nReplicas, width, offset, consuming); |
---|
207 | G4VPVParameterisation* pP = pVPV -> GetParameterisation (); |
---|
208 | if (pP) { // Parametrised volume. |
---|
209 | if (fCurrentDepth == 0) nReplicas = 1; // Just draw first |
---|
210 | for (int n = 0; n < nReplicas; n++) { |
---|
211 | pSol = pP -> ComputeSolid (n, pVPV); |
---|
212 | pP -> ComputeTransformation (n, pVPV); |
---|
213 | pSol -> ComputeDimensions (pP, n, pVPV); |
---|
214 | pVPV -> SetCopyNo (n); |
---|
215 | // Create a touchable of current parent for ComputeMaterial. |
---|
216 | // fFullPVPath has not been updated yet so at this point it |
---|
217 | // corresponds to the parent. |
---|
218 | G4PhysicalVolumeModelTouchable parentTouchable(fFullPVPath); |
---|
219 | pMaterial = pP -> ComputeMaterial (n, pVPV, &parentTouchable); |
---|
220 | DescribeAndDescend (pVPV, requestedDepth, pLV, pSol, pMaterial, |
---|
221 | theAT, sceneHandler); |
---|
222 | } |
---|
223 | } |
---|
224 | else { // Plain replicated volume. From geometry_guide.txt... |
---|
225 | // The replica's positions are claculated by means of a linear formula. |
---|
226 | // Replication may occur along: |
---|
227 | // |
---|
228 | // o Cartesian axes (kXAxis,kYAxis,kZAxis) |
---|
229 | // |
---|
230 | // The replications, of specified width have coordinates of |
---|
231 | // form (-width*(nReplicas-1)*0.5+n*width,0,0) where n=0.. nReplicas-1 |
---|
232 | // for the case of kXAxis, and are unrotated. |
---|
233 | // |
---|
234 | // o Radial axis (cylindrical polar) (kRho) |
---|
235 | // |
---|
236 | // The replications are cons/tubs sections, centred on the origin |
---|
237 | // and are unrotated. |
---|
238 | // They have radii of width*n+offset to width*(n+1)+offset |
---|
239 | // where n=0..nReplicas-1 |
---|
240 | // |
---|
241 | // o Phi axis (cylindrical polar) (kPhi) |
---|
242 | // The replications are `phi sections' or wedges, and of cons/tubs form |
---|
243 | // They have phi of offset+n*width to offset+(n+1)*width where |
---|
244 | // n=0..nReplicas-1 |
---|
245 | // |
---|
246 | pSol = pLV -> GetSolid (); |
---|
247 | pMaterial = pLV -> GetMaterial (); |
---|
248 | G4ThreeVector originalTranslation = pVPV -> GetTranslation (); |
---|
249 | G4RotationMatrix* pOriginalRotation = pVPV -> GetRotation (); |
---|
250 | G4double originalRMin = 0., originalRMax = 0.; |
---|
251 | if (axis == kRho && pSol->GetEntityType() == "G4Tubs") { |
---|
252 | originalRMin = ((G4Tubs*)pSol)->GetInnerRadius(); |
---|
253 | originalRMax = ((G4Tubs*)pSol)->GetOuterRadius(); |
---|
254 | } |
---|
255 | G4bool visualisable = true; |
---|
256 | for (int n = 0; n < nReplicas; n++) { |
---|
257 | G4ThreeVector translation; // Null. |
---|
258 | G4RotationMatrix rotation; // Null - life long enough for visualizing. |
---|
259 | G4RotationMatrix* pRotation = 0; |
---|
260 | switch (axis) { |
---|
261 | default: |
---|
262 | case kXAxis: |
---|
263 | translation = G4ThreeVector (-width*(nReplicas-1)*0.5+n*width,0,0); |
---|
264 | break; |
---|
265 | case kYAxis: |
---|
266 | translation = G4ThreeVector (0,-width*(nReplicas-1)*0.5+n*width,0); |
---|
267 | break; |
---|
268 | case kZAxis: |
---|
269 | translation = G4ThreeVector (0,0,-width*(nReplicas-1)*0.5+n*width); |
---|
270 | break; |
---|
271 | case kRho: |
---|
272 | if (pSol->GetEntityType() == "G4Tubs") { |
---|
273 | ((G4Tubs*)pSol)->SetInnerRadius(width*n+offset); |
---|
274 | ((G4Tubs*)pSol)->SetOuterRadius(width*(n+1)+offset); |
---|
275 | } else { |
---|
276 | if (fpMP->IsWarning()) |
---|
277 | G4cout << |
---|
278 | "G4PhysicalVolumeModel::VisitGeometryAndGetVisReps: WARNING:" |
---|
279 | "\n built-in replicated volumes replicated in radius for " |
---|
280 | << pSol->GetEntityType() << |
---|
281 | "-type\n solids (your solid \"" |
---|
282 | << pSol->GetName() << |
---|
283 | "\") are not visualisable." |
---|
284 | << G4endl; |
---|
285 | visualisable = false; |
---|
286 | } |
---|
287 | break; |
---|
288 | case kPhi: |
---|
289 | rotation.rotateZ (-(offset+(n+0.5)*width)); |
---|
290 | // Minus Sign because for the physical volume we need the |
---|
291 | // coordinate system rotation. |
---|
292 | pRotation = &rotation; |
---|
293 | break; |
---|
294 | } |
---|
295 | pVPV -> SetTranslation (translation); |
---|
296 | pVPV -> SetRotation (pRotation); |
---|
297 | pVPV -> SetCopyNo (n); |
---|
298 | if (visualisable) { |
---|
299 | DescribeAndDescend (pVPV, requestedDepth, pLV, pSol, pMaterial, |
---|
300 | theAT, sceneHandler); |
---|
301 | } |
---|
302 | } |
---|
303 | // Restore originals... |
---|
304 | pVPV -> SetTranslation (originalTranslation); |
---|
305 | pVPV -> SetRotation (pOriginalRotation); |
---|
306 | if (axis == kRho && pSol->GetEntityType() == "G4Tubs") { |
---|
307 | ((G4Tubs*)pSol)->SetInnerRadius(originalRMin); |
---|
308 | ((G4Tubs*)pSol)->SetOuterRadius(originalRMax); |
---|
309 | } |
---|
310 | } |
---|
311 | } |
---|
312 | |
---|
313 | return; |
---|
314 | } |
---|
315 | |
---|
316 | void G4PhysicalVolumeModel::DescribeAndDescend |
---|
317 | (G4VPhysicalVolume* pVPV, |
---|
318 | G4int requestedDepth, |
---|
319 | G4LogicalVolume* pLV, |
---|
320 | G4VSolid* pSol, |
---|
321 | G4Material* pMaterial, |
---|
322 | const G4Transform3D& theAT, |
---|
323 | G4VGraphicsScene& sceneHandler) |
---|
324 | { |
---|
325 | // Maintain useful data members... |
---|
326 | fpCurrentPV = pVPV; |
---|
327 | fpCurrentLV = pLV; |
---|
328 | fpCurrentMaterial = pMaterial; |
---|
329 | |
---|
330 | const G4RotationMatrix objectRotation = pVPV -> GetObjectRotationValue (); |
---|
331 | const G4ThreeVector& translation = pVPV -> GetTranslation (); |
---|
332 | G4Transform3D theLT (G4Transform3D (objectRotation, translation)); |
---|
333 | |
---|
334 | // Compute the accumulated transformation... |
---|
335 | // Note that top volume's transformation relative to the world |
---|
336 | // coordinate system is specified in theAT == startingTransformation |
---|
337 | // = fTransform (see DescribeYourselfTo), so first time through the |
---|
338 | // volume's own transformation, which is only relative to its |
---|
339 | // mother, i.e., not relative to the world coordinate system, should |
---|
340 | // not be accumulated. |
---|
341 | G4Transform3D theNewAT (theAT); |
---|
342 | if (fCurrentDepth != 0) theNewAT = theAT * theLT; |
---|
343 | fpCurrentTransform = &theNewAT; |
---|
344 | |
---|
345 | // Make decision to draw... |
---|
346 | const G4VisAttributes* pVisAttribs = pLV->GetVisAttributes(); |
---|
347 | if (!pVisAttribs) pVisAttribs = fpMP->GetDefaultVisAttributes(); |
---|
348 | // Beware - pVisAttribs might still be zero - create a temporary default one... |
---|
349 | G4bool visAttsCreated = false; |
---|
350 | if (!pVisAttribs) { |
---|
351 | pVisAttribs = new G4VisAttributes; |
---|
352 | visAttsCreated = true; |
---|
353 | } |
---|
354 | |
---|
355 | // From here, can assume pVisAttribs is a valid pointer. |
---|
356 | |
---|
357 | G4bool thisToBeDrawn = true; |
---|
358 | |
---|
359 | // There are various reasons why this volume |
---|
360 | // might not be drawn... |
---|
361 | G4bool culling = fpMP->IsCulling(); |
---|
362 | G4bool cullingInvisible = fpMP->IsCullingInvisible(); |
---|
363 | G4bool markedVisible = pVisAttribs->IsVisible(); |
---|
364 | G4bool cullingLowDensity = fpMP->IsDensityCulling(); |
---|
365 | G4double density = pMaterial? pMaterial->GetDensity(): 0; |
---|
366 | G4double densityCut = fpMP -> GetVisibleDensity (); |
---|
367 | |
---|
368 | // 1) Global culling is on.... |
---|
369 | if (culling) { |
---|
370 | // 2) Culling of invisible volumes is on... |
---|
371 | if (cullingInvisible) { |
---|
372 | // 3) ...and the volume is marked not visible... |
---|
373 | if (!markedVisible) thisToBeDrawn = false; |
---|
374 | } |
---|
375 | // 4) Or culling of low density volumes is on... |
---|
376 | if (cullingLowDensity) { |
---|
377 | // 5) ...and density is less than cut value... |
---|
378 | if (density < densityCut) thisToBeDrawn = false; |
---|
379 | } |
---|
380 | } |
---|
381 | |
---|
382 | // Update full path of physical volumes... |
---|
383 | G4int copyNo = fpCurrentPV->GetCopyNo(); |
---|
384 | fFullPVPath.push_back |
---|
385 | (G4PhysicalVolumeNodeID |
---|
386 | (fpCurrentPV,copyNo,fCurrentDepth,*fpCurrentTransform)); |
---|
387 | |
---|
388 | if (thisToBeDrawn) { |
---|
389 | |
---|
390 | // Update path of drawn physical volumes... |
---|
391 | G4int copyNo = fpCurrentPV->GetCopyNo(); |
---|
392 | fDrawnPVPath.push_back |
---|
393 | (G4PhysicalVolumeNodeID |
---|
394 | (fpCurrentPV,copyNo,fCurrentDepth,*fpCurrentTransform)); |
---|
395 | |
---|
396 | if (fpMP->IsExplode() && fDrawnPVPath.size() == 1) { |
---|
397 | // For top-level drawn volumes, explode along radius... |
---|
398 | G4Transform3D centering = G4Translate3D(fpMP->GetExplodeCentre()); |
---|
399 | G4Transform3D centred = centering.inverse() * theNewAT; |
---|
400 | G4Scale3D scale; |
---|
401 | G4Rotate3D rotation; |
---|
402 | G4Translate3D translation; |
---|
403 | centred.getDecomposition(scale, rotation, translation); |
---|
404 | G4double explodeFactor = fpMP->GetExplodeFactor(); |
---|
405 | G4Translate3D newTranslation = |
---|
406 | G4Translate3D(explodeFactor * translation.dx(), |
---|
407 | explodeFactor * translation.dy(), |
---|
408 | explodeFactor * translation.dz()); |
---|
409 | theNewAT = centering * newTranslation * rotation * scale; |
---|
410 | } |
---|
411 | |
---|
412 | DescribeSolid (theNewAT, pSol, pVisAttribs, sceneHandler); |
---|
413 | |
---|
414 | } |
---|
415 | |
---|
416 | // Make decision to draw daughters, if any. There are various |
---|
417 | // reasons why daughters might not be drawn... |
---|
418 | |
---|
419 | // First, reasons that do not depend on culling policy... |
---|
420 | G4int nDaughters = pLV->GetNoDaughters(); |
---|
421 | G4bool daughtersToBeDrawn = true; |
---|
422 | // 1) There are no daughters... |
---|
423 | if (!nDaughters) daughtersToBeDrawn = false; |
---|
424 | // 2) We are at the limit if requested depth... |
---|
425 | else if (requestedDepth == 0) daughtersToBeDrawn = false; |
---|
426 | // 3) The user has asked that the descent be curtailed... |
---|
427 | else if (fCurtailDescent) daughtersToBeDrawn = false; |
---|
428 | |
---|
429 | // Now, reasons that depend on culling policy... |
---|
430 | else { |
---|
431 | G4bool culling = fpMP->IsCulling(); |
---|
432 | G4bool cullingInvisible = fpMP->IsCullingInvisible(); |
---|
433 | G4bool daughtersInvisible = pVisAttribs->IsDaughtersInvisible(); |
---|
434 | // Culling of covered daughters request. This is computed in |
---|
435 | // G4VSceneHandler::CreateModelingParameters() depending on view |
---|
436 | // parameters... |
---|
437 | G4bool cullingCovered = fpMP->IsCullingCovered(); |
---|
438 | G4bool surfaceDrawing = |
---|
439 | fpMP->GetDrawingStyle() == G4ModelingParameters::hsr || |
---|
440 | fpMP->GetDrawingStyle() == G4ModelingParameters::hlhsr; |
---|
441 | if (pVisAttribs->IsForceDrawingStyle()) { |
---|
442 | switch (pVisAttribs->GetForcedDrawingStyle()) { |
---|
443 | default: |
---|
444 | case G4VisAttributes::wireframe: surfaceDrawing = false; break; |
---|
445 | case G4VisAttributes::solid: surfaceDrawing = true; break; |
---|
446 | } |
---|
447 | } |
---|
448 | G4bool opaque = pVisAttribs->GetColour().GetAlpha() >= 1.; |
---|
449 | // 4) Global culling is on.... |
---|
450 | if (culling) { |
---|
451 | // 5) ..and culling of invisible volumes is on... |
---|
452 | if (cullingInvisible) { |
---|
453 | // 6) ...and the mother requests daughters invisible |
---|
454 | if (daughtersInvisible) daughtersToBeDrawn = false; |
---|
455 | } |
---|
456 | // 7) Or culling of covered daughters is requested... |
---|
457 | if (cullingCovered) { |
---|
458 | // 8) ...and surface drawing is operating... |
---|
459 | if (surfaceDrawing) { |
---|
460 | // 9) ...but only if mother is visible... |
---|
461 | if (thisToBeDrawn) { |
---|
462 | // 10) ...and opaque... |
---|
463 | if (opaque) daughtersToBeDrawn = false; |
---|
464 | } |
---|
465 | } |
---|
466 | } |
---|
467 | } |
---|
468 | } |
---|
469 | |
---|
470 | // Vis atts for this volume no longer needed if created... |
---|
471 | if (visAttsCreated) delete pVisAttribs; |
---|
472 | |
---|
473 | if (daughtersToBeDrawn) { |
---|
474 | for (G4int iDaughter = 0; iDaughter < nDaughters; iDaughter++) { |
---|
475 | G4VPhysicalVolume* pVPV = pLV -> GetDaughter (iDaughter); |
---|
476 | // Descend the geometry structure recursively... |
---|
477 | fCurrentDepth++; |
---|
478 | VisitGeometryAndGetVisReps |
---|
479 | (pVPV, requestedDepth - 1, theNewAT, sceneHandler); |
---|
480 | fCurrentDepth--; |
---|
481 | } |
---|
482 | } |
---|
483 | |
---|
484 | // Reset for normal descending of next volume at this level... |
---|
485 | fCurtailDescent = false; |
---|
486 | |
---|
487 | // Pop item from paths physical volumes... |
---|
488 | fFullPVPath.pop_back(); |
---|
489 | if (thisToBeDrawn) { |
---|
490 | fDrawnPVPath.pop_back(); |
---|
491 | } |
---|
492 | } |
---|
493 | |
---|
494 | void G4PhysicalVolumeModel::DescribeSolid |
---|
495 | (const G4Transform3D& theAT, |
---|
496 | G4VSolid* pSol, |
---|
497 | const G4VisAttributes* pVisAttribs, |
---|
498 | G4VGraphicsScene& sceneHandler) |
---|
499 | { |
---|
500 | sceneHandler.PreAddSolid (theAT, *pVisAttribs); |
---|
501 | |
---|
502 | const G4Polyhedron* pSectionPolyhedron = fpMP->GetSectionPolyhedron(); |
---|
503 | const G4Polyhedron* pCutawayPolyhedron = fpMP->GetCutawayPolyhedron(); |
---|
504 | |
---|
505 | if (!fpClippingPolyhedron && !pSectionPolyhedron && !pCutawayPolyhedron) { |
---|
506 | |
---|
507 | pSol -> DescribeYourselfTo (sceneHandler); // Standard treatment. |
---|
508 | |
---|
509 | } else { |
---|
510 | |
---|
511 | // Clipping, etc., performed by Boolean operations on polyhedron objects. |
---|
512 | |
---|
513 | // First, get polyhedron for current solid... |
---|
514 | if (pVisAttribs->IsForceLineSegmentsPerCircle()) |
---|
515 | G4Polyhedron::SetNumberOfRotationSteps |
---|
516 | (pVisAttribs->GetForcedLineSegmentsPerCircle()); |
---|
517 | else |
---|
518 | G4Polyhedron::SetNumberOfRotationSteps(fpMP->GetNoOfSides()); |
---|
519 | G4Polyhedron* pOriginal = pSol->GetPolyhedron(); |
---|
520 | G4Polyhedron::ResetNumberOfRotationSteps(); |
---|
521 | |
---|
522 | if (!pOriginal) { |
---|
523 | |
---|
524 | if (fpMP->IsWarning()) |
---|
525 | G4cout << |
---|
526 | "WARNING: G4PhysicalVolumeModel::DescribeSolid: solid\n \"" |
---|
527 | << pSol->GetName() << |
---|
528 | "\" has no polyhedron. Cannot by clipped." |
---|
529 | << G4endl; |
---|
530 | pSol -> DescribeYourselfTo (sceneHandler); // Standard treatment. |
---|
531 | |
---|
532 | } else { |
---|
533 | |
---|
534 | G4Polyhedron resultant = *pOriginal; |
---|
535 | G4VisAttributes resultantVisAttribs(*pVisAttribs); |
---|
536 | |
---|
537 | if (fpClippingPolyhedron) { |
---|
538 | G4Polyhedron clipper = *fpClippingPolyhedron; // Local copy. |
---|
539 | clipper.Transform(theAT.inverse()); |
---|
540 | HepPolyhedronProcessor processor; |
---|
541 | switch (fClippingMode) { |
---|
542 | default: |
---|
543 | case subtraction: processor.push_back(HepPolyhedronProcessor::SUBTRACTION, clipper); break; |
---|
544 | case intersection: processor.push_back(HepPolyhedronProcessor::INTERSECTION, clipper); break; |
---|
545 | } |
---|
546 | if (!processor.execute(resultant)) { |
---|
547 | if (fpMP->IsWarning()) |
---|
548 | G4cout << |
---|
549 | "WARNING: G4PhysicalVolumeModel::DescribeSolid: clipped polyhedron for" |
---|
550 | "\n solid \"" << pSol->GetName() << |
---|
551 | "\" not defined due to error during Boolean processing." |
---|
552 | "\n It will be drawn in red." |
---|
553 | << G4endl; |
---|
554 | // Nevertheless, keep resultant, but draw it in red |
---|
555 | resultantVisAttribs.SetColour(G4Colour::Red()); |
---|
556 | } |
---|
557 | } |
---|
558 | |
---|
559 | if (pSectionPolyhedron) { |
---|
560 | G4Polyhedron sectioner = *pSectionPolyhedron; // Local copy. |
---|
561 | sectioner.Transform(theAT.inverse()); |
---|
562 | HepPolyhedronProcessor processor; |
---|
563 | processor.push_back(HepPolyhedronProcessor::INTERSECTION, sectioner); |
---|
564 | if (!processor.execute(resultant)) { |
---|
565 | if (fpMP->IsWarning()) |
---|
566 | G4cout << |
---|
567 | "WARNING: G4PhysicalVolumeModel::DescribeSolid: sectioned polyhedron for" |
---|
568 | "\n solid \"" << pSol->GetName() << |
---|
569 | "\" not defined due to error during Boolean processing." |
---|
570 | "\n It will be drawn in red." |
---|
571 | << G4endl; |
---|
572 | // Nevertheless, keep resultant, but draw it in red |
---|
573 | resultantVisAttribs.SetColour(G4Colour::Red()); |
---|
574 | } |
---|
575 | } |
---|
576 | |
---|
577 | if (pCutawayPolyhedron) { |
---|
578 | G4Polyhedron cutter = *pCutawayPolyhedron; // Local copy. |
---|
579 | cutter.Transform(theAT.inverse()); |
---|
580 | HepPolyhedronProcessor processor; |
---|
581 | processor.push_back(HepPolyhedronProcessor::SUBTRACTION, cutter); |
---|
582 | if (!processor.execute(resultant)) { |
---|
583 | if (fpMP->IsWarning()) |
---|
584 | G4cout << |
---|
585 | "WARNING: G4PhysicalVolumeModel::DescribeSolid: cutaway polyhedron for" |
---|
586 | "\n solid \"" << pSol->GetName() << |
---|
587 | "\" not defined due to error during Boolean processing." |
---|
588 | "\n It will be drawn in red." |
---|
589 | << G4endl; |
---|
590 | // Nevertheless, keep resultant, but draw it in red |
---|
591 | resultantVisAttribs.SetColour(G4Colour::Red()); |
---|
592 | } |
---|
593 | } |
---|
594 | |
---|
595 | // Finally, force polyhedron drawing... |
---|
596 | resultant.SetVisAttributes(resultantVisAttribs); |
---|
597 | sceneHandler.BeginPrimitives(theAT); |
---|
598 | sceneHandler.AddPrimitive(resultant); |
---|
599 | sceneHandler.EndPrimitives(); |
---|
600 | } |
---|
601 | } |
---|
602 | sceneHandler.PostAddSolid (); |
---|
603 | } |
---|
604 | |
---|
605 | G4bool G4PhysicalVolumeModel::Validate (G4bool warn) |
---|
606 | { |
---|
607 | G4VPhysicalVolume* world = |
---|
608 | G4TransportationManager::GetTransportationManager () |
---|
609 | -> GetNavigatorForTracking () -> GetWorldVolume (); |
---|
610 | // The idea now is to seek a PV with the same name and copy no |
---|
611 | // in the hope it's the same one!! |
---|
612 | if (warn) { |
---|
613 | G4cout << "G4PhysicalVolumeModel::Validate() called." << G4endl; |
---|
614 | } |
---|
615 | G4PhysicalVolumeModel searchModel (world); |
---|
616 | G4PhysicalVolumeSearchScene searchScene |
---|
617 | (&searchModel, fTopPVName, fTopPVCopyNo); |
---|
618 | G4ModelingParameters mp; // Default modeling parameters for this search. |
---|
619 | mp.SetDefaultVisAttributes(fpMP? fpMP->GetDefaultVisAttributes(): 0); |
---|
620 | searchModel.SetModelingParameters (&mp); |
---|
621 | searchModel.DescribeYourselfTo (searchScene); |
---|
622 | G4VPhysicalVolume* foundVolume = searchScene.GetFoundVolume (); |
---|
623 | if (foundVolume) { |
---|
624 | if (warn) { |
---|
625 | G4cout << " Volume of the same name and copy number (\"" |
---|
626 | << fTopPVName << "\", copy " << fTopPVCopyNo |
---|
627 | << ") still exists and is being used." |
---|
628 | "\n WARNING: This does not necessarily guarantee it's the same" |
---|
629 | "\n volume you originally specified in /vis/scene/add/." |
---|
630 | << G4endl; |
---|
631 | } |
---|
632 | fpTopPV = foundVolume; |
---|
633 | CalculateExtent (); |
---|
634 | return true; |
---|
635 | } |
---|
636 | else { |
---|
637 | if (warn) { |
---|
638 | G4cout << " A volume of the same name and copy number (\"" |
---|
639 | << fTopPVName << "\", copy " << fTopPVCopyNo |
---|
640 | << ") no longer exists." |
---|
641 | << G4endl; |
---|
642 | } |
---|
643 | return false; |
---|
644 | } |
---|
645 | } |
---|
646 | |
---|
647 | const std::map<G4String,G4AttDef>* G4PhysicalVolumeModel::GetAttDefs() const |
---|
648 | { |
---|
649 | G4bool isNew; |
---|
650 | std::map<G4String,G4AttDef>* store |
---|
651 | = G4AttDefStore::GetInstance("G4PhysicalVolumeModel", isNew); |
---|
652 | if (isNew) { |
---|
653 | (*store)["PVPath"] = |
---|
654 | G4AttDef("PVPath","Physical Volume Path","Physics","","G4String"); |
---|
655 | (*store)["LVol"] = |
---|
656 | G4AttDef("LVol","Logical Volume","Physics","","G4String"); |
---|
657 | (*store)["Solid"] = |
---|
658 | G4AttDef("Solid","Solid Name","Physics","","G4String"); |
---|
659 | (*store)["EType"] = |
---|
660 | G4AttDef("EType","Entity Type","Physics","","G4String"); |
---|
661 | (*store)["DmpSol"] = |
---|
662 | G4AttDef("DmpSol","Dump of Solid properties","Physics","","G4String"); |
---|
663 | (*store)["Trans"] = |
---|
664 | G4AttDef("Trans","Transformation of volume","Physics","","G4String"); |
---|
665 | (*store)["Material"] = |
---|
666 | G4AttDef("Material","Material Name","Physics","","G4String"); |
---|
667 | (*store)["Density"] = |
---|
668 | G4AttDef("Density","Material Density","Physics","G4BestUnit","G4double"); |
---|
669 | (*store)["State"] = |
---|
670 | G4AttDef("State","Material State (enum undefined,solid,liquid,gas)","Physics","","G4String"); |
---|
671 | (*store)["Radlen"] = |
---|
672 | G4AttDef("Radlen","Material Radiation Length","Physics","G4BestUnit","G4double"); |
---|
673 | } |
---|
674 | (*store)["Region"] = |
---|
675 | G4AttDef("Region","Cuts Region","Physics","","G4String"); |
---|
676 | (*store)["RootRegion"] = |
---|
677 | G4AttDef("RootRegion","Root Region (0/1 = false/true)","Physics","","G4bool"); |
---|
678 | return store; |
---|
679 | } |
---|
680 | |
---|
681 | #include <iomanip> |
---|
682 | |
---|
683 | static std::ostream& operator<< (std::ostream& o, const G4Transform3D t) |
---|
684 | { |
---|
685 | using namespace std; |
---|
686 | |
---|
687 | G4Scale3D s; |
---|
688 | G4Rotate3D r; |
---|
689 | G4Translate3D tl; |
---|
690 | t.getDecomposition(s, r, tl); |
---|
691 | |
---|
692 | const int w = 10; |
---|
693 | |
---|
694 | // Transformation itself |
---|
695 | o << setw(w) << t.xx() << setw(w) << t.xy() << setw(w) << t.xz() << setw(w) << t.dx() << endl; |
---|
696 | o << setw(w) << t.yx() << setw(w) << t.yy() << setw(w) << t.yz() << setw(w) << t.dy() << endl; |
---|
697 | o << setw(w) << t.zx() << setw(w) << t.zy() << setw(w) << t.zz() << setw(w) << t.dz() << endl; |
---|
698 | |
---|
699 | // Translation |
---|
700 | o << "= translation:" << endl; |
---|
701 | o << setw(w) << tl.dx() << setw(w) << tl.dy() << setw(w) << tl.dz() << endl; |
---|
702 | |
---|
703 | // Rotation |
---|
704 | o << "* rotation:" << endl; |
---|
705 | o << setw(w) << r.xx() << setw(w) << r.xy() << setw(w) << r.xz() << endl; |
---|
706 | o << setw(w) << r.yx() << setw(w) << r.yy() << setw(w) << r.yz() << endl; |
---|
707 | o << setw(w) << r.zx() << setw(w) << r.zy() << setw(w) << r.zz() << endl; |
---|
708 | |
---|
709 | // Scale |
---|
710 | o << "* scale:" << endl; |
---|
711 | o << setw(w) << s.xx() << setw(w) << s.yy() << setw(w) << s.zz() << endl; |
---|
712 | |
---|
713 | // Transformed axes |
---|
714 | o << "Transformed axes:" << endl; |
---|
715 | o << "x': " << r * G4Vector3D(1., 0., 0.) << endl; |
---|
716 | o << "y': " << r * G4Vector3D(0., 1., 0.) << endl; |
---|
717 | o << "z': " << r * G4Vector3D(0., 0., 1.) << endl; |
---|
718 | |
---|
719 | return o; |
---|
720 | } |
---|
721 | |
---|
722 | std::vector<G4AttValue>* G4PhysicalVolumeModel::CreateCurrentAttValues() const |
---|
723 | { |
---|
724 | std::vector<G4AttValue>* values = new std::vector<G4AttValue>; |
---|
725 | std::ostringstream oss; |
---|
726 | for (size_t i = 0; i < fFullPVPath.size(); ++i) { |
---|
727 | oss << fFullPVPath[i].GetPhysicalVolume()->GetName() |
---|
728 | << ':' << fFullPVPath[i].GetCopyNo(); |
---|
729 | if (i != fFullPVPath.size() - 1) oss << '/'; |
---|
730 | } |
---|
731 | |
---|
732 | if (!fpCurrentLV) { |
---|
733 | G4Exception |
---|
734 | ("G4PhysicalVolumeModel::CreateCurrentAttValues", |
---|
735 | "", |
---|
736 | JustWarning, |
---|
737 | "Current logical volume not defined."); |
---|
738 | return values; |
---|
739 | } |
---|
740 | |
---|
741 | values->push_back(G4AttValue("PVPath", oss.str(),"")); |
---|
742 | values->push_back(G4AttValue("LVol", fpCurrentLV->GetName(),"")); |
---|
743 | G4VSolid* pSol = fpCurrentLV->GetSolid(); |
---|
744 | values->push_back(G4AttValue("Solid", pSol->GetName(),"")); |
---|
745 | values->push_back(G4AttValue("EType", pSol->GetEntityType(),"")); |
---|
746 | oss.str(""); oss << '\n' << *pSol; |
---|
747 | values->push_back(G4AttValue("DmpSol", oss.str(),"")); |
---|
748 | oss.str(""); oss << '\n' << *fpCurrentTransform; |
---|
749 | values->push_back(G4AttValue("Trans", oss.str(),"")); |
---|
750 | G4String matName = fpCurrentMaterial? fpCurrentMaterial->GetName(): G4String("No material"); |
---|
751 | values->push_back(G4AttValue("Material", matName,"")); |
---|
752 | G4double matDensity = fpCurrentMaterial? fpCurrentMaterial->GetDensity(): 0.; |
---|
753 | values->push_back(G4AttValue("Density", G4BestUnit(matDensity,"Volumic Mass"),"")); |
---|
754 | G4State matState = fpCurrentMaterial? fpCurrentMaterial->GetState(): kStateUndefined; |
---|
755 | oss.str(""); oss << matState; |
---|
756 | values->push_back(G4AttValue("State", oss.str(),"")); |
---|
757 | G4double matRadlen = fpCurrentMaterial? fpCurrentMaterial->GetRadlen(): 0.; |
---|
758 | values->push_back(G4AttValue("Radlen", G4BestUnit(matRadlen,"Length"),"")); |
---|
759 | G4Region* region = fpCurrentLV->GetRegion(); |
---|
760 | G4String regionName = region? region->GetName(): G4String("No region"); |
---|
761 | values->push_back(G4AttValue("Region", regionName,"")); |
---|
762 | oss.str(""); oss << fpCurrentLV->IsRootRegion(); |
---|
763 | values->push_back(G4AttValue("RootRegion", oss.str(),"")); |
---|
764 | return values; |
---|
765 | } |
---|
766 | |
---|
767 | G4bool G4PhysicalVolumeModel::G4PhysicalVolumeNodeID::operator< |
---|
768 | (const G4PhysicalVolumeModel::G4PhysicalVolumeNodeID& right) const |
---|
769 | { |
---|
770 | if (fpPV < right.fpPV) return true; |
---|
771 | if (fpPV == right.fpPV) { |
---|
772 | if (fCopyNo < right.fCopyNo) return true; |
---|
773 | if (fCopyNo == right.fCopyNo) |
---|
774 | return fNonCulledDepth < right.fNonCulledDepth; |
---|
775 | } |
---|
776 | return false; |
---|
777 | } |
---|
778 | |
---|
779 | std::ostream& operator<< |
---|
780 | (std::ostream& os, const G4PhysicalVolumeModel::G4PhysicalVolumeNodeID node) |
---|
781 | { |
---|
782 | G4VPhysicalVolume* pPV = node.GetPhysicalVolume(); |
---|
783 | if (pPV) { |
---|
784 | os << pPV->GetName() |
---|
785 | << ':' << node.GetCopyNo() |
---|
786 | << '[' << node.GetNonCulledDepth() << ']' |
---|
787 | << ':' << node.GetTransform(); |
---|
788 | } else { |
---|
789 | os << "Null node"; |
---|
790 | } |
---|
791 | return os; |
---|
792 | } |
---|
793 | |
---|
794 | G4PhysicalVolumeModel::G4PhysicalVolumeModelTouchable::G4PhysicalVolumeModelTouchable |
---|
795 | (const std::vector<G4PhysicalVolumeNodeID>& fullPVPath): |
---|
796 | fFullPVPath(fullPVPath) {} |
---|
797 | |
---|
798 | const G4ThreeVector& G4PhysicalVolumeModel::G4PhysicalVolumeModelTouchable::GetTranslation(G4int depth) const |
---|
799 | { |
---|
800 | size_t i = fFullPVPath.size() - depth - 1; |
---|
801 | if (i >= fFullPVPath.size()) { |
---|
802 | G4Exception("G4PhysicalVolumeModelTouchable::GetTranslation", |
---|
803 | "Index out of range", |
---|
804 | FatalErrorInArgument, |
---|
805 | "Asking for non-existent depth"); |
---|
806 | } |
---|
807 | static G4ThreeVector tempTranslation; |
---|
808 | tempTranslation = fFullPVPath[i].GetTransform().getTranslation(); |
---|
809 | return tempTranslation; |
---|
810 | } |
---|
811 | |
---|
812 | const G4RotationMatrix* G4PhysicalVolumeModel::G4PhysicalVolumeModelTouchable::GetRotation(G4int depth) const |
---|
813 | { |
---|
814 | size_t i = fFullPVPath.size() - depth - 1; |
---|
815 | if (i >= fFullPVPath.size()) { |
---|
816 | G4Exception("G4PhysicalVolumeModelTouchable::GetRotation", |
---|
817 | "Index out of range", |
---|
818 | FatalErrorInArgument, |
---|
819 | "Asking for non-existent depth"); |
---|
820 | } |
---|
821 | static G4RotationMatrix tempRotation; |
---|
822 | tempRotation = fFullPVPath[i].GetTransform().getRotation(); |
---|
823 | return &tempRotation; |
---|
824 | } |
---|
825 | |
---|
826 | G4VPhysicalVolume* G4PhysicalVolumeModel::G4PhysicalVolumeModelTouchable::GetVolume(G4int depth) const |
---|
827 | { |
---|
828 | size_t i = fFullPVPath.size() - depth - 1; |
---|
829 | if (i >= fFullPVPath.size()) { |
---|
830 | G4Exception("G4PhysicalVolumeModelTouchable::GetVolume", |
---|
831 | "Index out of range", |
---|
832 | FatalErrorInArgument, |
---|
833 | "Asking for non-existent depth"); |
---|
834 | } |
---|
835 | return fFullPVPath[i].GetPhysicalVolume(); |
---|
836 | } |
---|
837 | |
---|
838 | G4VSolid* G4PhysicalVolumeModel::G4PhysicalVolumeModelTouchable::GetSolid(G4int depth) const |
---|
839 | { |
---|
840 | size_t i = fFullPVPath.size() - depth - 1; |
---|
841 | if (i >= fFullPVPath.size()) { |
---|
842 | G4Exception("G4PhysicalVolumeModelTouchable::GetSolid", |
---|
843 | "Index out of range", |
---|
844 | FatalErrorInArgument, |
---|
845 | "Asking for non-existent depth"); |
---|
846 | } |
---|
847 | return fFullPVPath[i].GetPhysicalVolume()->GetLogicalVolume()->GetSolid(); |
---|
848 | } |
---|
849 | |
---|
850 | G4int G4PhysicalVolumeModel::G4PhysicalVolumeModelTouchable::GetReplicaNumber(G4int depth) const |
---|
851 | { |
---|
852 | size_t i = fFullPVPath.size() - depth - 1; |
---|
853 | if (i >= fFullPVPath.size()) { |
---|
854 | G4Exception("G4PhysicalVolumeModelTouchable::GetReplicaNumber", |
---|
855 | "Index out of range", |
---|
856 | FatalErrorInArgument, |
---|
857 | "Asking for non-existent depth"); |
---|
858 | } |
---|
859 | return fFullPVPath[i].GetCopyNo(); |
---|
860 | } |
---|