Geant4  10.02.p03
G4ASCIITreeSceneHandler Class Reference

#include <G4ASCIITreeSceneHandler.hh>

Inheritance diagram for G4ASCIITreeSceneHandler:
Collaboration diagram for G4ASCIITreeSceneHandler:

Public Member Functions

 G4ASCIITreeSceneHandler (G4VGraphicsSystem &system, const G4String &name)
 
virtual ~G4ASCIITreeSceneHandler ()
 
virtual void BeginModeling ()
 
virtual void EndModeling ()
 
- Public Member Functions inherited from G4VTreeSceneHandler
 G4VTreeSceneHandler (G4VGraphicsSystem &system, const G4String &name)
 
virtual ~G4VTreeSceneHandler ()
 
void PreAddSolid (const G4Transform3D &objectTransformation, const G4VisAttributes &)
 
void PostAddSolid ()
 
virtual void AddPrimitive (const G4Polyline &)
 
virtual void AddPrimitive (const G4Text &)
 
virtual void AddPrimitive (const G4Circle &)
 
virtual void AddPrimitive (const G4Square &)
 
virtual void AddPrimitive (const G4Polyhedron &)
 
virtual void AddPrimitive (const G4Polymarker &)
 
virtual void AddPrimitive (const G4Scale &)
 
- Public Member Functions inherited from G4VSceneHandler
 G4VSceneHandler (G4VGraphicsSystem &system, G4int id, const G4String &name="")
 
virtual ~G4VSceneHandler ()
 
virtual void AddSolid (const G4Box &)
 
virtual void AddSolid (const G4Cons &)
 
virtual void AddSolid (const G4Tubs &)
 
virtual void AddSolid (const G4Trd &)
 
virtual void AddSolid (const G4Trap &)
 
virtual void AddSolid (const G4Sphere &)
 
virtual void AddSolid (const G4Para &)
 
virtual void AddSolid (const G4Torus &)
 
virtual void AddSolid (const G4Polycone &)
 
virtual void AddSolid (const G4Polyhedra &)
 
virtual void AddSolid (const G4VSolid &)
 
virtual void AddCompound (const G4VTrajectory &)
 
virtual void AddCompound (const G4VHit &)
 
virtual void AddCompound (const G4VDigi &)
 
virtual void AddCompound (const G4THitsMap< G4double > &)
 
virtual void BeginPrimitives (const G4Transform3D &objectTransformation)
 
virtual void EndPrimitives ()
 
virtual void BeginPrimitives2D (const G4Transform3D &objectTransformation)
 
virtual void EndPrimitives2D ()
 
virtual const G4VisExtentGetExtent () const
 
const G4StringGetName () const
 
G4int GetSceneHandlerId () const
 
G4int GetViewCount () const
 
G4VGraphicsSystemGetGraphicsSystem () const
 
G4SceneGetScene () const
 
const G4ViewerListGetViewerList () const
 
G4VModelGetModel () const
 
G4VViewerGetCurrentViewer () const
 
G4bool GetMarkForClearingTransientStore () const
 
G4bool IsReadyForTransients () const
 
G4bool GetTransientsDrawnThisEvent () const
 
G4bool GetTransientsDrawnThisRun () const
 
const G4Transform3DGetObjectTransformation () const
 
void SetName (const G4String &)
 
void SetCurrentViewer (G4VViewer *)
 
virtual void SetScene (G4Scene *)
 
G4ViewerListSetViewerList ()
 
void SetModel (G4VModel *)
 
void SetMarkForClearingTransientStore (G4bool)
 
void SetTransientsDrawnThisEvent (G4bool)
 
void SetTransientsDrawnThisRun (G4bool)
 
void SetObjectTransformation (const G4Transform3D &)
 
const G4ColourGetColour (const G4Visible &)
 
const G4ColourGetColor (const G4Visible &)
 
const G4ColourGetTextColour (const G4Text &)
 
const G4ColourGetTextColor (const G4Text &)
 
G4double GetLineWidth (const G4VisAttributes *)
 
G4ViewParameters::DrawingStyle GetDrawingStyle (const G4VisAttributes *)
 
G4bool GetAuxEdgeVisible (const G4VisAttributes *)
 
G4int GetNoOfSides (const G4VisAttributes *)
 
G4double GetMarkerSize (const G4VMarker &, MarkerSizeType &)
 
G4double GetMarkerDiameter (const G4VMarker &, MarkerSizeType &)
 
G4double GetMarkerRadius (const G4VMarker &, MarkerSizeType &)
 
G4ModelingParametersCreateModelingParameters ()
 
void DrawEvent (const G4Event *)
 
void DrawEndOfRunModels ()
 
G4int IncrementViewCount ()
 
virtual void ClearStore ()
 
virtual void ClearTransientStore ()
 
void AddViewerToList (G4VViewer *pView)
 
void RemoveViewerFromList (G4VViewer *pView)
 
- Public Member Functions inherited from G4VGraphicsScene
 G4VGraphicsScene ()
 
virtual ~G4VGraphicsScene ()
 

Protected Types

typedef std::set< G4LogicalVolume * >::iterator LVSetIterator
 
typedef G4PhysicalVolumeModel::G4PhysicalVolumeNodeID PVNodeID
 
typedef std::vector< PVNodeIDPVPath
 
typedef std::set< PVPath >::iterator ReplicaSetIterator
 
typedef std::set< PVPath >::reverse_iterator ReplicaSetReverseIterator
 

Protected Member Functions

virtual void RequestPrimitives (const G4VSolid &)
 
void WriteHeader (std::ostream &)
 
- Protected Member Functions inherited from G4VSceneHandler
virtual void ProcessScene ()
 
virtual G4VSolidCreateSectionSolid ()
 
virtual G4VSolidCreateCutawaySolid ()
 
void LoadAtts (const G4Visible &, G4AttHolder *)
 

Protected Attributes

std::ostream * fpOutFile
 
std::ofstream fOutFile
 
std::ostringstream fRestOfLine
 
const G4VPhysicalVolumefpLastPV
 
G4String fLastPVName
 
G4int fLastCopyNo
 
G4int fLastNonSequentialCopyNo
 
std::set< G4LogicalVolume * > fLVSet
 
std::set< PVPathfReplicaSet
 
- Protected Attributes inherited from G4VTreeSceneHandler
const G4Transform3DfpCurrentObjectTransformation
 
std::set< G4LogicalVolume * > fDrawnLVStore
 
- Protected Attributes inherited from G4VSceneHandler
G4VGraphicsSystemfSystem
 
const G4int fSceneHandlerId
 
G4String fName
 
G4int fViewCount
 
G4ViewerList fViewerList
 
G4VViewerfpViewer
 
G4ScenefpScene
 
G4bool fMarkForClearingTransientStore
 
G4bool fReadyForTransients
 
G4bool fTransientsDrawnThisEvent
 
G4bool fTransientsDrawnThisRun
 
G4bool fProcessingSolid
 
G4bool fProcessing2D
 
G4VModelfpModel
 
G4Transform3D fObjectTransformation
 
G4int fNestingDepth
 
const G4VisAttributesfpVisAttribs
 
const G4Transform3D fIdentityTransformation
 

Additional Inherited Members

- Public Types inherited from G4VSceneHandler
enum  MarkerSizeType { world, screen }
 
- Static Protected Attributes inherited from G4VTreeSceneHandler
static G4int fSceneIdCount = 0
 

Detailed Description

Definition at line 48 of file G4ASCIITreeSceneHandler.hh.

Member Typedef Documentation

◆ LVSetIterator

typedef std::set<G4LogicalVolume*>::iterator G4ASCIITreeSceneHandler::LVSetIterator
protected

Definition at line 73 of file G4ASCIITreeSceneHandler.hh.

◆ PVNodeID

◆ PVPath

typedef std::vector<PVNodeID> G4ASCIITreeSceneHandler::PVPath
protected

Definition at line 75 of file G4ASCIITreeSceneHandler.hh.

◆ ReplicaSetIterator

typedef std::set<PVPath>::iterator G4ASCIITreeSceneHandler::ReplicaSetIterator
protected

Definition at line 77 of file G4ASCIITreeSceneHandler.hh.

◆ ReplicaSetReverseIterator

typedef std::set<PVPath>::reverse_iterator G4ASCIITreeSceneHandler::ReplicaSetReverseIterator
protected

Definition at line 78 of file G4ASCIITreeSceneHandler.hh.

Constructor & Destructor Documentation

◆ G4ASCIITreeSceneHandler()

G4ASCIITreeSceneHandler::G4ASCIITreeSceneHandler ( G4VGraphicsSystem system,
const G4String name 
)

Definition at line 56 of file G4ASCIITreeSceneHandler.cc.

57  :
58  G4VTreeSceneHandler(system, name),
59  fpOutFile(0),
60  fpLastPV(0),
61  fLastCopyNo(-99),
63 {}
G4VTreeSceneHandler(G4VGraphicsSystem &system, const G4String &name)
const G4VPhysicalVolume * fpLastPV

◆ ~G4ASCIITreeSceneHandler()

G4ASCIITreeSceneHandler::~G4ASCIITreeSceneHandler ( )
virtual

Definition at line 65 of file G4ASCIITreeSceneHandler.cc.

65 {}

Member Function Documentation

◆ BeginModeling()

void G4ASCIITreeSceneHandler::BeginModeling ( )
virtual

Reimplemented from G4VTreeSceneHandler.

Definition at line 67 of file G4ASCIITreeSceneHandler.cc.

67  {
68 
69  G4VTreeSceneHandler::BeginModeling (); // To re-use "culling off" code.
70 
71  const G4ASCIITree* pSystem = (G4ASCIITree*)GetGraphicsSystem();
72  const G4String& outFileName = pSystem -> GetOutFileName();
73  if (outFileName == "G4cout") {
74  fpOutFile = &G4cout;
75  } else {
76  fOutFile.open (outFileName);
78  }
79 
80  G4cout << "G4ASCIITreeSceneHandler::BeginModeling: writing to ";
81  if (outFileName == "G4cout") {
82  G4cout << "G4 standard output (G4cout)";
83  } else {
84  G4cout << "file \"" << outFileName << "\"";
85  }
86  G4cout << G4endl;
87 
89  if (outFileName != "G4cout") {
90  WriteHeader (fOutFile); fOutFile << std::endl;
91  }
92 }
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61
void WriteHeader(std::ostream &)
virtual void BeginModeling()
G4VGraphicsSystem * GetGraphicsSystem() const
Here is the call graph for this function:

◆ EndModeling()

void G4ASCIITreeSceneHandler::EndModeling ( )
virtual

Reimplemented from G4VTreeSceneHandler.

Definition at line 117 of file G4ASCIITreeSceneHandler.cc.

117  {
118  const G4ASCIITree* pSystem = (G4ASCIITree*) GetGraphicsSystem();
119  const G4int verbosity = pSystem->GetVerbosity();
120  const G4int detail = verbosity % 10;
121  const G4String& outFileName = pSystem -> GetOutFileName();
122 
123  // Output left over copy number, if any...
125  if (fLastCopyNo == fLastNonSequentialCopyNo + 1) *fpOutFile << ',';
126  else *fpOutFile << '-';
128  }
129  // Output outstanding rest of line, if any...
130  if (!fRestOfLine.str().empty()) *fpOutFile << fRestOfLine.str();
131  fRestOfLine.str("");
132  fpLastPV = 0;
133  fLastPVName.clear();
134  fLastCopyNo = -99;
136 
137  // This detail to G4cout regardless of outFileName...
138  if (detail >= 4) {
139  G4cout << "Calculating mass(es)..." << G4endl;
140  const std::vector<G4Scene::Model>& models = fpScene->GetRunDurationModelList();
141  std::vector<G4Scene::Model>::const_iterator i;
142  for (i = models.begin(); i != models.end(); ++i) {
143  G4PhysicalVolumeModel* pvModel =
144  dynamic_cast<G4PhysicalVolumeModel*>(i->fpModel);
145  if (pvModel) {
146  if (pvModel->GetTopPhysicalVolume() ==
149  const G4ModelingParameters* tempMP =
150  pvModel->GetModelingParameters();
151  G4ModelingParameters mp; // Default - no culling.
152  pvModel->SetModelingParameters (&mp);
153  G4PhysicalVolumeMassScene massScene(pvModel);
154  pvModel->DescribeYourselfTo (massScene);
155  G4double volume = massScene.GetVolume();
156  G4double mass = massScene.GetMass();
157 
158  G4cout << "Overall volume of \""
159  << pvModel->GetTopPhysicalVolume()->GetName()
160  << "\":"
161  << pvModel->GetTopPhysicalVolume()->GetCopyNo()
162  << ", is "
163  << G4BestUnit (volume, "Volume")
164  << " and the daughter-included mass";
165  G4int requestedDepth = pvModel->GetRequestedDepth();
166  if (requestedDepth == G4PhysicalVolumeModel::UNLIMITED) {
167  G4cout << " to unlimited depth";
168  } else {
169  G4cout << ", ignoring daughters at depth "
170  << requestedDepth
171  << " and below,";
172  }
173  G4cout << " is " << G4BestUnit (mass, "Mass")
174  << G4endl;
175 
176  pvModel->SetModelingParameters (tempMP);
177  }
178  }
179  }
180  }
181 
182  if (outFileName != "G4cout") {
183  fOutFile.close();
184  G4cout << "Output file \"" << outFileName << "\" closed." << G4endl;
185  }
186  fLVSet.clear();
187  fReplicaSet.clear();
188  G4cout << "G4ASCIITreeSceneHandler::EndModeling" << G4endl;
189  G4VTreeSceneHandler::EndModeling (); // To re-use "culling off" code.
190 }
G4VPhysicalVolume * GetTopPhysicalVolume() const
G4int GetVerbosity() const
Definition: G4ASCIITree.hh:47
const std::vector< Model > & GetRunDurationModelList() const
#define G4BestUnit(a, b)
#define G4_USE_G4BESTUNIT_FOR_VERBOSE 1
int G4int
Definition: G4Types.hh:78
void SetModelingParameters(const G4ModelingParameters *)
G4GLOB_DLL std::ostream G4cout
const G4ModelingParameters * GetModelingParameters() const
G4Navigator * GetNavigatorForTracking() const
static G4TransportationManager * GetTransportationManager()
const G4String & GetName() const
virtual void EndModeling()
virtual G4int GetCopyNo() const =0
const G4VPhysicalVolume * fpLastPV
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
std::set< G4LogicalVolume * > fLVSet
G4VPhysicalVolume * GetWorldVolume() const
void DescribeYourselfTo(G4VGraphicsScene &)
G4VGraphicsSystem * GetGraphicsSystem() const
Here is the call graph for this function:

◆ RequestPrimitives()

void G4ASCIITreeSceneHandler::RequestPrimitives ( const G4VSolid solid)
protectedvirtual

Reimplemented from G4VSceneHandler.

Definition at line 192 of file G4ASCIITreeSceneHandler.cc.

192  {
193 
194  G4PhysicalVolumeModel* pPVModel =
195  dynamic_cast<G4PhysicalVolumeModel*>(fpModel);
196  if (!pPVModel) return; // Not from a G4PhysicalVolumeModel.
197 
198  // This call comes from a G4PhysicalVolumeModel. drawnPVPath is
199  // the path of the current drawn (non-culled) volume in terms of
200  // drawn (non-culled) ancesters. Each node is identified by a
201  // PVNodeID object, which is a physical volume and copy number. It
202  // is a vector of PVNodeIDs corresponding to the geometry hierarchy
203  // actually selected, i.e., not culled.
204  // The following typedef's already set in header file...
205  //typedef G4PhysicalVolumeModel::G4PhysicalVolumeNodeID PVNodeID;
206  //typedef std::vector<PVNodeID> PVPath;
207  const PVPath& drawnPVPath = pPVModel->GetDrawnPVPath();
208  //G4int currentDepth = pPVModel->GetCurrentDepth();
209  G4VPhysicalVolume* pCurrentPV = pPVModel->GetCurrentPV();
210  G4LogicalVolume* pCurrentLV = pPVModel->GetCurrentLV();
211  G4Material* pCurrentMaterial = pPVModel->GetCurrentMaterial();
212 
214  G4int verbosity = pSystem->GetVerbosity();
215  G4int detail = verbosity % 10;
216 
217  if (verbosity < 10 && pCurrentPV->IsReplicated()) {
218  // See if this has been a replica found before with same mother LV...
219  PVPath::const_reverse_iterator thisID = drawnPVPath.rbegin();
220  PVPath::const_reverse_iterator motherID = ++drawnPVPath.rbegin();
221  G4bool ignore = false;
222  for (ReplicaSetIterator i = fReplicaSet.begin(); i != fReplicaSet.end();
223  ++i) {
224  if (i->back().GetPhysicalVolume()->GetLogicalVolume() ==
225  thisID->GetPhysicalVolume()->GetLogicalVolume()) {
226  // For each one previously found (if more than one, they must
227  // have different mothers)...
228  // To avoid compilation errors on VC++ .Net 7.1...
229  // Previously:
230  // PVNodeID previousMotherID = ++(i->rbegin());
231  // (Should that have been: PVNodeID::const_iterator previousMotherID?)
232  // Replace
233  // previousMotherID == i->rend()
234  // by
235  // i->size() <= 1
236  // Replace
237  // previousMotherID != i->rend()
238  // by
239  // i->size() > 1
240  // Replace
241  // previousMotherID->
242  // by
243  // (*i)[i->size() - 2].
244  if (motherID == drawnPVPath.rend() &&
245  i->size() <= 1)
246  ignore = true; // Both have no mother.
247  if (motherID != drawnPVPath.rend() &&
248  i->size() > 1 &&
249  motherID->GetPhysicalVolume()->GetLogicalVolume() ==
250  (*i)[i->size() - 2].GetPhysicalVolume()->GetLogicalVolume())
251  ignore = true; // Same mother LV...
252  }
253  }
254  if (ignore) {
255  pPVModel->CurtailDescent();
256  return;
257  }
258  }
259 
260  const G4String& currentPVName = pCurrentPV->GetName();
261  const G4int currentCopyNo = pCurrentPV->GetCopyNo();
262 
263  if (verbosity < 10 &&
264  currentPVName == fLastPVName &&
265  currentCopyNo != fLastCopyNo) {
266  // For low verbosity, trap series of G4PVPlacements with the same name but
267  // different copy number (replicas should have been taken out above)...
268  // Check...
269  if (pCurrentPV->IsReplicated()) {
270  G4Exception("G4ASCIITreeSceneHandler::RequestPrimitives",
271  "vistree0001",
272  JustWarning,
273  "Replica unexpected");
274  }
275  // Check mothers are identical...
276  else if (pCurrentLV == (fpLastPV? fpLastPV->GetLogicalVolume(): 0)) {
277  if (currentCopyNo != fLastCopyNo + 1) {
278  // Non-sequential copy number...
279  *fpOutFile << ',' << currentCopyNo;
280  fLastNonSequentialCopyNo = currentCopyNo;
281  }
282  fLastCopyNo = currentCopyNo;
283  pPVModel->CurtailDescent();
284  return;
285  }
286  }
287  fpLastPV = pCurrentPV;
288 
289  // High verbosity or a new G4VPhysicalVolume...
290  // Output copy number, if any, from previous invocation...
292  if (fLastCopyNo == fLastNonSequentialCopyNo + 1) *fpOutFile << ',';
293  else *fpOutFile << '-';
295  }
296  // Output rest of line, if any, from previous invocation...
297  if (!fRestOfLine.str().empty()) *fpOutFile << fRestOfLine.str();
298  fRestOfLine.str("");
299  fLastPVName = currentPVName;
300  fLastCopyNo = currentCopyNo;
301  fLastNonSequentialCopyNo = currentCopyNo;
302  // Indent according to drawn path depth...
303  for (size_t i = 0; i < drawnPVPath.size(); i++ ) *fpOutFile << " ";
304  // Start next line...
305  *fpOutFile << "\"" << currentPVName
306  << "\":" << currentCopyNo;
307 
308  if (pCurrentPV->IsReplicated()) {
309  if (verbosity < 10) {
310  // Add printing for replicas (when replicas are ignored)...
311  EAxis axis;
312  G4int nReplicas;
313  G4double width;
314  G4double offset;
315  G4bool consuming;
316  pCurrentPV->GetReplicationData(axis,nReplicas,width,offset,consuming);
317  G4VPVParameterisation* pP = pCurrentPV->GetParameterisation();
318  if (pP) {
319  if (detail < 3) {
320  fReplicaSet.insert(drawnPVPath);
321  if (nReplicas > 2) fRestOfLine << '-';
322  else fRestOfLine << ',';
323  fRestOfLine << nReplicas - 1
324  << " (" << nReplicas << " parametrised volumes)";
325  }
326  }
327  else {
328  fReplicaSet.insert(drawnPVPath);
329  if (nReplicas > 2) fRestOfLine << '-';
330  else fRestOfLine << ',';
331  fRestOfLine << nReplicas - 1
332  << " (" << nReplicas << " replicas)";
333  }
334  }
335  } else {
336  if (fLVSet.find(pCurrentLV) != fLVSet.end()) {
337  if (verbosity < 10) {
338  // Add printing for repeated LV (if it has daughters)...
339  if (pCurrentLV->GetNoDaughters()) fRestOfLine << " (repeated LV)";
340  // ...and curtail descent.
341  pPVModel->CurtailDescent();
342  }
343  }
344  }
345 
346  if (detail >= 1) {
347  fRestOfLine << " / \""
348  << pCurrentLV->GetName() << "\"";
349  G4VSensitiveDetector* sd = pCurrentLV->GetSensitiveDetector();
350  if (sd) {
351  fRestOfLine << " (SD=\"" << sd->GetName() << "\"";
352  G4VReadOutGeometry* roGeom = sd->GetROgeometry();
353  if (roGeom) {
354  fRestOfLine << ",RO=\"" << roGeom->GetName() << "\"";
355  }
356  fRestOfLine << ")";
357  }
358  }
359 
360  if (detail >= 2) {
361  fRestOfLine << " / \""
362  << solid.GetName()
363  << "\"("
364  << solid.GetEntityType() << ")";
365  }
366 
367  if (detail >= 3) {
368  fRestOfLine << ", "
369  << G4BestUnit(((G4VSolid&)solid).GetCubicVolume(),"Volume")
370  << ", ";
371  if (pCurrentMaterial) {
373  << G4BestUnit(pCurrentMaterial->GetDensity(), "Volumic Mass")
374  << " (" << pCurrentMaterial->GetName() << ")";
375  } else {
376  fRestOfLine << "(No material)";
377  }
378  }
379 
380  if (detail >= 5) {
381  if (pCurrentMaterial) {
382  G4Material* pMaterial = const_cast<G4Material*>(pCurrentMaterial);
383  // ...and find daughter-subtracted mass...
384  G4double daughter_subtracted_mass = pCurrentLV->GetMass
385  (pCurrentPV->IsParameterised(), // Force if parametrised.
386  false, // Do not propagate - calculate for this volume minus
387  // volume of daughters.
388  pMaterial);
389  G4double daughter_subtracted_volume =
390  daughter_subtracted_mass / pCurrentMaterial->GetDensity();
391  fRestOfLine << ", "
392  << G4BestUnit(daughter_subtracted_volume,"Volume")
393  << ", "
394  << G4BestUnit(daughter_subtracted_mass,"Mass");
395  }
396  }
397 
398  if (detail >= 6) {
399  std::vector<G4AttValue>* attValues = pPVModel->CreateCurrentAttValues();
400  const std::map<G4String,G4AttDef>* attDefs = pPVModel->GetAttDefs();
401  fRestOfLine << '\n' << G4AttCheck(attValues, attDefs);
402  delete attValues;
403  }
404 
405  if (detail >= 7) {
406  G4Polyhedron* polyhedron = solid.GetPolyhedron();
407  fRestOfLine << "\nLocal polyhedron coordinates:\n" << *polyhedron;
408  G4Transform3D* transform = pPVModel->GetCurrentTransform();
409  polyhedron->Transform(*transform);
410  fRestOfLine << "\nGlobal polyhedron coordinates:\n" << *polyhedron;
411  }
412 
413  if (fLVSet.find(pCurrentLV) == fLVSet.end()) {
414  fLVSet.insert(pCurrentLV); // Record new logical volume.
415  }
416 
417  fRestOfLine << std::endl;
418 
419  return;
420 }
const std::vector< G4PhysicalVolumeNodeID > & GetDrawnPVPath() const
G4int GetNoDaughters() const
virtual G4bool IsReplicated() const =0
const std::map< G4String, G4AttDef > * GetAttDefs() const
std::vector< G4AttValue > * CreateCurrentAttValues() const
G4int GetVerbosity() const
Definition: G4ASCIITree.hh:47
G4double GetDensity() const
Definition: G4Material.hh:180
#define width
#define G4BestUnit(a, b)
#define G4_USE_G4BESTUNIT_FOR_VERBOSE 1
virtual G4GeometryType GetEntityType() const =0
int G4int
Definition: G4Types.hh:78
G4String GetName() const
virtual G4Polyhedron * GetPolyhedron() const
Definition: G4VSolid.cc:644
std::vector< PVNodeID > PVPath
std::set< PVPath >::iterator ReplicaSetIterator
const G4String & GetName() const
bool G4bool
Definition: G4Types.hh:79
virtual G4VPVParameterisation * GetParameterisation() const =0
G4LogicalVolume * GetCurrentLV() const
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
virtual G4bool IsParameterised() const =0
const G4String & GetName() const
G4Transform3D * GetCurrentTransform() const
EAxis
Definition: geomdefs.hh:54
G4String GetName() const
virtual G4int GetCopyNo() const =0
G4VSensitiveDetector * GetSensitiveDetector() const
const G4VPhysicalVolume * fpLastPV
G4VPhysicalVolume * GetCurrentPV() const
G4Material * GetCurrentMaterial() const
virtual void GetReplicationData(EAxis &axis, G4int &nReplicas, G4double &width, G4double &offset, G4bool &consuming) const =0
G4double GetMass(G4bool forced=false, G4bool propagate=true, G4Material *parMaterial=0)
double G4double
Definition: G4Types.hh:76
G4VReadOutGeometry * GetROgeometry() const
std::set< G4LogicalVolume * > fLVSet
G4LogicalVolume * GetLogicalVolume() const
const G4String & GetName() const
Definition: G4Material.hh:178
G4VGraphicsSystem * GetGraphicsSystem() const
Here is the call graph for this function:

◆ WriteHeader()

void G4ASCIITreeSceneHandler::WriteHeader ( std::ostream &  os)
protected

Definition at line 94 of file G4ASCIITreeSceneHandler.cc.

95 {
96  const G4ASCIITree* pSystem = (G4ASCIITree*)GetGraphicsSystem();
97  const G4int verbosity = pSystem->GetVerbosity();
98  const G4int detail = verbosity % 10;
99  os << "# Set verbosity with \"/vis/ASCIITree/verbose <verbosity>\":";
100  for (size_t i = 0;
102  os << "\n# " << G4ASCIITreeMessenger::fVerbosityGuidance[i];
103  }
104  os << "\n# Now printing with verbosity " << verbosity;
105  os << "\n# Format is: PV:n";
106  if (detail >= 1) os << " / LV (SD,RO)";
107  if (detail >= 2) os << " / Solid(type)";
108  if (detail >= 3) os << ", volume, density";
109  if (detail >= 5) os << ", daughter-subtracted volume and mass";
110  if (detail >= 6) os << ", physical volume dump";
111  if (detail >= 7) os << ", polyhedron dump";
112  os <<
113  "\n# Abbreviations: PV = Physical Volume, LV = Logical Volume,"
114  "\n# SD = Sensitive Detector, RO = Read Out Geometry.";
115 }
G4int GetVerbosity() const
Definition: G4ASCIITree.hh:47
int G4int
Definition: G4Types.hh:78
static std::vector< G4String > fVerbosityGuidance
G4VGraphicsSystem * GetGraphicsSystem() const
Here is the call graph for this function:
Here is the caller graph for this function:

Member Data Documentation

◆ fLastCopyNo

G4int G4ASCIITreeSceneHandler::fLastCopyNo
protected

Definition at line 69 of file G4ASCIITreeSceneHandler.hh.

◆ fLastNonSequentialCopyNo

G4int G4ASCIITreeSceneHandler::fLastNonSequentialCopyNo
protected

Definition at line 70 of file G4ASCIITreeSceneHandler.hh.

◆ fLastPVName

G4String G4ASCIITreeSceneHandler::fLastPVName
protected

Definition at line 68 of file G4ASCIITreeSceneHandler.hh.

◆ fLVSet

std::set<G4LogicalVolume*> G4ASCIITreeSceneHandler::fLVSet
protected

Definition at line 72 of file G4ASCIITreeSceneHandler.hh.

◆ fOutFile

std::ofstream G4ASCIITreeSceneHandler::fOutFile
protected

Definition at line 63 of file G4ASCIITreeSceneHandler.hh.

◆ fpLastPV

const G4VPhysicalVolume* G4ASCIITreeSceneHandler::fpLastPV
protected

Definition at line 67 of file G4ASCIITreeSceneHandler.hh.

◆ fpOutFile

std::ostream* G4ASCIITreeSceneHandler::fpOutFile
protected

Definition at line 62 of file G4ASCIITreeSceneHandler.hh.

◆ fReplicaSet

std::set<PVPath> G4ASCIITreeSceneHandler::fReplicaSet
protected

Definition at line 76 of file G4ASCIITreeSceneHandler.hh.

◆ fRestOfLine

std::ostringstream G4ASCIITreeSceneHandler::fRestOfLine
protected

Definition at line 64 of file G4ASCIITreeSceneHandler.hh.


The documentation for this class was generated from the following files: