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