Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4XXXSGSceneHandler.cc
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26//
27// $Id$
28//
29//
30// John Allison 10th March 2006
31// A template for a sophisticated graphics driver with a scene graph.
32//?? Lines beginning like this require specialisation for your driver.
33
35
36#include "G4XXXSGViewer.hh"
39#include "G4VPhysicalVolume.hh"
40#include "G4LogicalVolume.hh"
41#include "G4Box.hh"
42#include "G4Polyline.hh"
43#include "G4Text.hh"
44#include "G4Circle.hh"
45#include "G4Square.hh"
46#include "G4Polyhedron.hh"
47#include "G4UnitsTable.hh"
48
49#include <sstream>
50
52// Counter for XXX scene handlers.
53
55 const G4String& name):
56 G4VSceneHandler(system, fSceneIdCount++, name)
57{}
58
60
61#ifdef G4XXXSGDEBUG
62// Useful function...
63void G4XXXSGSceneHandler::PrintThings() {
64 G4cout <<
65 " with transformation "
66 << (void*)fpObjectTransformation;
67 if (fpModel) {
68 G4cout << " from " << fpModel->GetCurrentDescription()
69 << " (tag " << fpModel->GetCurrentTag()
70 << ')';
71 } else {
72 G4cout << "(not from a model)";
73 }
74 G4PhysicalVolumeModel* pPVModel =
75 dynamic_cast<G4PhysicalVolumeModel*>(fpModel);
76 if (pPVModel) {
77 G4cout <<
78 "\n current physical volume: "
79 << pPVModel->GetCurrentPV()->GetName() <<
80 "\n current logical volume: "
81// There might be a problem with the LV pointer if this is a G4LogicalVolumeModel
82 << pPVModel->GetCurrentLV()->GetName() <<
83 "\n current depth of geometry tree: "
84 << pPVModel->GetCurrentDepth();
85 }
86 G4cout << G4endl;
87}
88#endif
89
91 // Utility for PreAddSolid and BeginPrimitives.
92
93 G4PhysicalVolumeModel* pPVModel =
94 dynamic_cast<G4PhysicalVolumeModel*>(fpModel);
95 G4LogicalVolumeModel* pLVModel =
96 dynamic_cast<G4LogicalVolumeModel*>(pPVModel);
97 if (pPVModel && !pLVModel) {
98
99 // This call comes from a G4PhysicalVolumeModel. drawnPVPath is
100 // the path of the current drawn (non-culled) volume in terms of
101 // drawn (non-culled) ancesters. Each node is identified by a
102 // PVNodeID object, which is a physical volume and copy number. It
103 // is a vector of PVNodeIDs corresponding to the geometry hierarchy
104 // actually selected, i.e., not culled.
106 typedef std::vector<PVNodeID> PVPath;
107 const PVPath& drawnPVPath = pPVModel->GetDrawnPVPath();
108 //G4int currentDepth = pPVModel->GetCurrentDepth();
109 //G4VPhysicalVolume* pCurrentPV = pPVModel->GetCurrentPV();
110 //G4LogicalVolume* pCurrentLV = pPVModel->GetCurrentLV();
111 //G4Material* pCurrentMaterial = pPVModel->GetCurrentMaterial();
112 // Note: pCurrentMaterial may be zero (parallel world).
113
114 // The simplest algorithm, used by the Open Inventor Driver
115 // developers, is to rely on the fact the G4PhysicalVolumeModel
116 // traverses the geometry hierarchy in an orderly manner. The last
117 // mother, if any, will be the node to which the volume should be
118 // added. So it is enough to keep a map of scene graph nodes keyed
119 // on the volume path ID. Actually, it is enough to use the logical
120 // volume as the key. (An alternative would be to keep the PVNodeID
121 // in the tree and match the PVPath from the root down.)
122
123 // BUT IN OPENGL, IF THERE ARE TRANSPARENT OBJECTS, VOLUMES DO NOT
124 // ARRIVE IN THE ABOVE ORDER. (TRANSPARENT OBJECTS ARE DRWAN
125 // LAST.) SO WE MUST BE MORE SOPHISTICATED IN CONSTRUCTING A
126 // TREE.
127
128 /* Debug
129 for (size_t i = 0; i < drawnPVPath.size(); ++i) {
130 std::cout << drawnPVPath[i].GetPhysicalVolume()->GetName() << ":"
131 << drawnPVPath[i].GetCopyNo() << " ("
132 << currentPOListIndex << "), ";
133 }
134 std::cout << std::endl;
135 */
136
137 static G4int index = 0; // Some index for future reference
138 JA::Insert(&drawnPVPath[0],drawnPVPath.size(),index++,&fSceneGraph);
139 //JA::PrintTree(std::cout,&root);
140
141 /*** Old algorithm, left here for historical interest!!
142 // Find mother. ri points to drawn mother, if any.
143 PVPath::const_reverse_iterator ri = ++drawnPVPath.rbegin();
144 if (ri != drawnPVPath.rend()) {
145 // This volume has a mother.
146 G4LogicalVolume* drawnMotherLV =
147 ri->GetPhysicalVolume()->GetLogicalVolume();
148 LVMapIterator mother = fLVMap.find(drawnMotherLV);
149 if (mother != fLVMap.end()) {
150 // This adds a child in Troy's tree...
151 fCurrentItem = mother->second.push_back(header);
152 } else {
153 // Mother not previously encountered. Shouldn't happen, since
154 // G4PhysicalVolumeModel sends volumes as it encounters them,
155 // i.e., mothers before daughters, in its descent of the
156 // geometry tree. Error!
157 G4cout << "ERROR: G4XXXSGSceneHandler::PreAddSolid: Mother "
158 << ri->GetPhysicalVolume()->GetName()
159 << ':' << ri->GetCopyNo()
160 << " not previously encountered."
161 "\nShouldn't happen! Please report to visualization coordinator."
162 << G4endl;
163 // Continue anyway. Add to root of scene graph tree...
164 fCurrentItem = fPermanentsRoot.push_back(header);
165 }
166 } else {
167 // This volume has no mother. Must be a top level un-culled
168 // volume. Add to root of scene graph tree...
169 fCurrentItem = fPermanentsRoot.push_back(header);
170 }
171
172 std::ostringstream oss;
173 oss << "Path of drawn PVs: ";
174 for (PVPath::const_iterator i = drawnPVPath.begin();
175 i != drawnPVPath.end(); ++i) {
176 oss << '/' << i->GetPhysicalVolume()->GetName()
177 << ':' << i->GetCopyNo();
178 }
179 oss << std::endl;
180 *fCurrentItem += oss.str();
181
182 // Store for future searches. Overwrites previous entries for this
183 // LV, so entry is always the *last* LV.
184 fLVMap[pCurrentLV] = fCurrentItem;
185 ***/
186
187 } else { // Not from a G4PhysicalVolumeModel.
188
189 /***
190 // Create a place for current solid in root...
191 if (fReadyForTransients) {
192 fCurrentItem = fTransientsRoot.push_back(header);
193 } else {
194 fCurrentItem = fPermanentsRoot.push_back(header);
195 }
196 ***/
197 }
198}
199
201(const G4Transform3D& objectTransformation,
202 const G4VisAttributes& visAttribs)
203{
204 G4VSceneHandler::PreAddSolid(objectTransformation, visAttribs);
205 CreateCurrentItem(G4String("\nPreAddSolid:\n"));
206}
207
209{
211}
212
214(const G4Transform3D& objectTransformation)
215{
216 G4VSceneHandler::BeginPrimitives(objectTransformation);
217}
218
220{
222}
223
224// Note: This function overrides G4VSceneHandler::AddSolid(const
225// G4Box&). You may not want to do this, but this is how it's done if
226// you do. Certain other specific solids may be treated this way -
227// see G4VSceneHandler.hh. The simplest possible driver would *not*
228// implement these polymorphic functions, with the effect that the
229// default versions in G4VSceneHandler are used, which simply call
230// G4VSceneHandler::RequestPrimitives to turn the solid into a
231// G4Polyhedron usually.
232// Don't forget, solids can be transients too (e.g., representing a hit).
234#ifdef G4XXXSGDEBUG
235 G4cout <<
236 "G4XXXSGSceneHandler::AddSolid(const G4Box& box) called for "
237 << box.GetName()
238 << G4endl;
239#endif
240 //?? Process your box...
241 std::ostringstream oss;
242 oss << "G4Box(" <<
246 (box.GetXHalfLength(), box.GetYHalfLength(), box.GetZHalfLength()),
247 "Length")).strip() << ')' << std::endl;
248 //*fCurrentItem += oss.str();
249}
250
252#ifdef G4XXXSGDEBUG
253 G4cout <<
254 "G4XXXSGSceneHandler::AddPrimitive(const G4Polyline& polyline) called.\n"
255 << polyline
256 << G4endl;
257#endif
258 // Get vis attributes - pick up defaults if none.
259 //const G4VisAttributes* pVA =
260 // fpViewer -> GetApplicableVisAttributes (polyline.GetVisAttributes ());
261 //?? Process polyline.
262 std::ostringstream oss;
263 oss << polyline << std::endl;
264 //*fCurrentItem += oss.str();
265}
266
268#ifdef G4XXXSGDEBUG
269 G4cout <<
270 "G4XXXSGSceneHandler::AddPrimitive(const G4Text& text) called.\n"
271 << text
272 << G4endl;
273#endif
274 // Get text colour - special method since default text colour is
275 // determined by the default text vis attributes, which may be
276 // specified independent of default vis attributes of other types of
277 // visible objects.
278 //const G4Colour& c = GetTextColour (text); // Picks up default if none.
279 //?? Process text.
280 std::ostringstream oss;
281 oss << text << std::endl;
282 //*fCurrentItem += oss.str();
283}
284
286#ifdef G4XXXSGDEBUG
287 G4cout <<
288 "G4XXXSGSceneHandler::AddPrimitive(const G4Circle& circle) called.\n"
289 << circle
290 << G4endl;
291 MarkerSizeType sizeType;
292 G4double size = GetMarkerSize (circle, sizeType);
293 switch (sizeType) {
294 default:
295 case screen:
296 // Draw in screen coordinates.
297 G4cout << "screen";
298 break;
299 case world:
300 // Draw in world coordinates.
301 G4cout << "world";
302 break;
303 }
304 G4cout << " size: " << size << G4endl;
305#endif
306 // Get vis attributes - pick up defaults if none.
307 //const G4VisAttributes* pVA =
308 // fpViewer -> GetApplicableVisAttributes (circle.GetVisAttributes ());
309 //?? Process circle.
310 std::ostringstream oss;
311 oss << circle << std::endl;
312 //*fCurrentItem += oss.str();
313}
314
316#ifdef G4XXXSGDEBUG
317 G4cout <<
318 "G4XXXSGSceneHandler::AddPrimitive(const G4Square& square) called.\n"
319 << square
320 << G4endl;
321 MarkerSizeType sizeType;
322 G4double size = GetMarkerSize (square, sizeType);
323 switch (sizeType) {
324 default:
325 case screen:
326 // Draw in screen coordinates.
327 G4cout << "screen";
328 break;
329 case world:
330 // Draw in world coordinates.
331 G4cout << "world";
332 break;
333 }
334 G4cout << " size: " << size << G4endl;
335#endif
336 // Get vis attributes - pick up defaults if none.
337 //const G4VisAttributes* pVA =
338 // fpViewer -> GetApplicableVisAttributes (square.GetVisAttributes ());
339 //?? Process square.
340 std::ostringstream oss;
341 oss << square << std::endl;
342 //*fCurrentItem += oss.str();
343}
344
346#ifdef G4XXXSGDEBUG
347 G4cout <<
348 "G4XXXSGSceneHandler::AddPrimitive(const G4Polyhedron& polyhedron) called.\n"
349 << polyhedron
350 << G4endl;
351#endif
352 //?? Process polyhedron.
353 std::ostringstream oss;
354 oss << polyhedron;
355 //*fCurrentItem += oss.str();
356
357 //?? Or... here are some ideas for decomposing into polygons...
358 //Assume all facets are convex quadrilaterals.
359 //Draw each G4Facet individually
360
361 //Get colour, etc..
362 if (polyhedron.GetNoFacets() == 0) return;
363
364 // Get vis attributes - pick up defaults if none.
365 const G4VisAttributes* pVA =
366 fpViewer -> GetApplicableVisAttributes (polyhedron.GetVisAttributes ());
367
368 // Get view parameters that the user can force through the vis
369 // attributes, thereby over-riding the current view parameter.
371 //G4bool isAuxEdgeVisible = GetAuxEdgeVisible (pVA);
372
373 //Get colour, etc..
374 //const G4Colour& c = pVA -> GetColour ();
375
376 // Initial action depending on drawing style.
377 switch (drawing_style) {
379 {
380 break;
381 }
383 {
384 break;
385 }
387 {
388 break;
389 }
390 default:
391 {
392 break;
393 }
394 }
395
396 // Loop through all the facets...
397
398 // Look at G4OpenGLSceneHandler::AddPrimitive(const G4Polyhedron&)
399 // for an example of how to get facets out of a G4Polyhedron,
400 // including how to cope with triangles if that's a problem.
401}
402
404#ifdef G4XXXSGDEBUG
405 G4cout <<
406 "G4XXXSGSceneHandler::AddPrimitive(const G4NURBS& nurbs) called."
407 << G4endl;
408#endif
409 //?? Don't bother implementing this. NURBS are not functional.
410}
411
413{
415}
416
418{
420}
421
422namespace JA {
423// Ad hoc tree class and utilities.
424
425#include "G4VPhysicalVolume.hh"
426
427void Insert(const PVNodeID* pvPath, size_t pathLength,
428 G4int index, Node* node) {
429 // Path passed as a PVNodeID* to avoid copying.
430
431 /* Debug
432 for (size_t i = 0; i < pathLength; ++i) {
433 std::cout << pvPath[i].GetPhysicalVolume()->GetName() << ":"
434 << pvPath[i].GetCopyNo() << " ("
435 << index << "), ";
436 }
437 */
438
439 // See if node has been encountered before
440 G4bool found = false; size_t foundPosition = 0;
441 for (size_t i = 0; i < node->fDaughters.size(); ++i) {
442 PVNodeID& daughterPVNodeID = node->fDaughters[i]->fPVNodeID;
443 // It is enough to compare volume and copy number at a given position in the tree
444 if (daughterPVNodeID.GetPhysicalVolume() == pvPath[0].GetPhysicalVolume() &&
445 daughterPVNodeID.GetCopyNo() == pvPath[0].GetCopyNo()) {
446 found = true;
447 foundPosition = i;
448 break;
449 }
450 }
451
452 if (pathLength == 1) { // This is a leaf
453 if (found) { // Update index
454 node->fDaughters[foundPosition]->fIndex = index;
455 } else { // Make a new full entry
456 node->fDaughters.push_back(new Node(pvPath[0],index));
457 }
458 /* Debug
459 std::cout << std::endl;
460 */
461 } else { // Not a leaf - carry on with rest of path
462 if (found) { // Just carry on
463 Insert(pvPath+1,--pathLength,index,
464 node->fDaughters[foundPosition]);
465 } else { // Insert place holder, then carry on
466 node->fDaughters.push_back(new Node(pvPath[0]));
467 Insert(pvPath+1,--pathLength,index,
468 node->fDaughters[node->fDaughters.size()-1]);
469 }
470 }
471}
472
473void PrintTree(std::ostream& os, Node* node)
474{
475 static G4int depth = -1;
476 depth++;
477 PVNodeID& thisPVNodeID = node->fPVNodeID;
478 G4int& thisIndex = node->fIndex;
479 const size_t& nDaughters = node->fDaughters.size();
480 G4VPhysicalVolume* thisPhysicalVolume= thisPVNodeID.GetPhysicalVolume();
481 if (!thisPhysicalVolume) os << "Root" << std::endl;
482 else {
483 for (G4int i = 0; i < depth; ++i) os << "__";
484 os << thisPVNodeID.GetPhysicalVolume()->GetName() << ":"
485 << thisPVNodeID.GetCopyNo() << " ("
486 << thisIndex << ")" << std::endl;;
487 }
488 for (size_t i = 0; i < nDaughters; ++i) {
489 PrintTree(os, node->fDaughters[i]);
490 }
491 depth--;
492}
493
494void Clear(Node* node)
495{
496 const size_t& nDaughters = node->fDaughters.size();
497 for (size_t i = 0; i < nDaughters; ++i) {
498 Clear(node->fDaughters[i]);
499 delete node->fDaughters[i];
500 }
501}
502
503} // End namespace JA
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
Definition: G4Box.hh:55
G4double GetYHalfLength() const
G4double GetZHalfLength() const
G4double GetXHalfLength() const
G4String GetName() const
const std::vector< G4PhysicalVolumeNodeID > & GetDrawnPVPath() const
G4VPhysicalVolume * GetCurrentPV() const
G4LogicalVolume * GetCurrentLV() const
G4String strip(G4int strip_Type=trailing, char c=' ')
Definition: G4Text.hh:73
virtual G4String GetCurrentDescription() const
Definition: G4VModel.cc:54
virtual G4String GetCurrentTag() const
Definition: G4VModel.cc:49
const G4String & GetName() const
virtual void EndPrimitives()
G4double GetMarkerSize(const G4VMarker &, MarkerSizeType &)
virtual void PreAddSolid(const G4Transform3D &objectTransformation, const G4VisAttributes &)
G4VViewer * fpViewer
virtual void PostAddSolid()
virtual void BeginPrimitives(const G4Transform3D &objectTransformation)
G4ViewParameters::DrawingStyle GetDrawingStyle(const G4VisAttributes *)
G4String GetName() const
const G4VisAttributes * GetVisAttributes() const
void CreateCurrentItem(const G4String &)
G4XXXSGSceneHandler(G4VGraphicsSystem &system, const G4String &name)
void AddPrimitive(const G4Polyline &)
void PreAddSolid(const G4Transform3D &objectTransformation, const G4VisAttributes &)
void BeginPrimitives(const G4Transform3D &objectTransformation)
void AddSolid(const G4Box &)
G4int GetNoFacets() const
void PrintTree(std::ostream &, Node *)
void Insert(const PVNodeID *pvPath, size_t pathLength, G4int index, Node *node)
void Clear(Node *)
PVNodeID fPVNodeID
std::vector< Node * > fDaughters