Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4ASCIITreeSceneHandler Class Reference

#include <G4ASCIITreeSceneHandler.hh>

+ Inheritance 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 G4NURBS &)
 
virtual void AddPrimitive (const G4Polymarker &)
 
virtual void AddPrimitive (const G4Scale &)
 
virtual void BeginModeling ()
 
virtual void EndModeling ()
 
- Public Member Functions inherited from G4VSceneHandler
 G4VSceneHandler (G4VGraphicsSystem &system, G4int id, const G4String &name="")
 
virtual ~G4VSceneHandler ()
 
virtual void PreAddSolid (const G4Transform3D &objectTransformation, const G4VisAttributes &)
 
virtual void PostAddSolid ()
 
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 BeginModeling ()
 
virtual void EndModeling ()
 
virtual void BeginPrimitives (const G4Transform3D &objectTransformation)
 
virtual void EndPrimitives ()
 
virtual void BeginPrimitives2D (const G4Transform3D &objectTransformation)
 
virtual void EndPrimitives2D ()
 
virtual void AddPrimitive (const G4Polyline &)=0
 
virtual void AddPrimitive (const G4Scale &)
 
virtual void AddPrimitive (const G4Text &)=0
 
virtual void AddPrimitive (const G4Circle &)=0
 
virtual void AddPrimitive (const G4Square &)=0
 
virtual void AddPrimitive (const G4Polymarker &)
 
virtual void AddPrimitive (const G4Polyhedron &)=0
 
virtual void AddPrimitive (const G4NURBS &)=0
 
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 ()
 
virtual void PreAddSolid (const G4Transform3D &objectTransformation, const G4VisAttributes &visAttribs)=0
 
virtual void PostAddSolid ()=0
 
virtual void AddSolid (const G4Box &)=0
 
virtual void AddSolid (const G4Cons &)=0
 
virtual void AddSolid (const G4Tubs &)=0
 
virtual void AddSolid (const G4Trd &)=0
 
virtual void AddSolid (const G4Trap &)=0
 
virtual void AddSolid (const G4Sphere &)=0
 
virtual void AddSolid (const G4Para &)=0
 
virtual void AddSolid (const G4Torus &)=0
 
virtual void AddSolid (const G4Polycone &)=0
 
virtual void AddSolid (const G4Polyhedra &)=0
 
virtual void AddSolid (const G4VSolid &)=0
 
virtual void AddCompound (const G4VTrajectory &)=0
 
virtual void AddCompound (const G4VHit &)=0
 
virtual void AddCompound (const G4VDigi &)=0
 
virtual void AddCompound (const G4THitsMap< G4double > &)=0
 
virtual void BeginPrimitives (const G4Transform3D &objectTransformation=G4Transform3D())=0
 
virtual void EndPrimitives ()=0
 
virtual void BeginPrimitives2D (const G4Transform3D &objectTransformation=G4Transform3D())=0
 
virtual void EndPrimitives2D ()=0
 
virtual void AddPrimitive (const G4Polyline &)=0
 
virtual void AddPrimitive (const G4Scale &)=0
 
virtual void AddPrimitive (const G4Text &)=0
 
virtual void AddPrimitive (const G4Circle &)=0
 
virtual void AddPrimitive (const G4Square &)=0
 
virtual void AddPrimitive (const G4Polymarker &)=0
 
virtual void AddPrimitive (const G4Polyhedron &)=0
 
virtual void AddPrimitive (const G4NURBS &)=0
 

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 void RequestPrimitives (const G4VSolid &solid)
 
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 
)

◆ ~G4ASCIITreeSceneHandler()

G4ASCIITreeSceneHandler::~G4ASCIITreeSceneHandler ( )
virtual

Definition at line 63 of file G4ASCIITreeSceneHandler.cc.

63{}

Member Function Documentation

◆ BeginModeling()

void G4ASCIITreeSceneHandler::BeginModeling ( )
virtual

Reimplemented from G4VTreeSceneHandler.

Definition at line 65 of file G4ASCIITreeSceneHandler.cc.

65 {
66
67 G4VTreeSceneHandler::BeginModeling (); // To re-use "culling off" code.
68
69 const G4ASCIITree* pSystem = (G4ASCIITree*)GetGraphicsSystem();
70 const G4String& outFileName = pSystem -> GetOutFileName();
71 if (outFileName == "G4cout") {
73 } else {
74 fOutFile.open (outFileName);
76 }
77
78 G4cout << "G4ASCIITreeSceneHandler::BeginModeling: writing to ";
79 if (outFileName == "G4cout") {
80 G4cout << "G4 standard output (G4cout)";
81 } else {
82 G4cout << "file \"" << outFileName << "\"";
83 }
84 G4cout << G4endl;
85
87 if (outFileName != "G4cout") {
88 WriteHeader (fOutFile); fOutFile << std::endl;
89 }
90}
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
void WriteHeader(std::ostream &)
G4VGraphicsSystem * GetGraphicsSystem() const
virtual void BeginModeling()

◆ EndModeling()

void G4ASCIITreeSceneHandler::EndModeling ( )
virtual

Reimplemented from G4VTreeSceneHandler.

Definition at line 113 of file G4ASCIITreeSceneHandler.cc.

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

◆ RequestPrimitives()

void G4ASCIITreeSceneHandler::RequestPrimitives ( const G4VSolid solid)
protectedvirtual

Reimplemented from G4VSceneHandler.

Definition at line 188 of file G4ASCIITreeSceneHandler.cc.

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

◆ WriteHeader()

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

Definition at line 92 of file G4ASCIITreeSceneHandler.cc.

93{
94 const G4ASCIITree* pSystem = (G4ASCIITree*)GetGraphicsSystem();
95 const G4int verbosity = pSystem->GetVerbosity();
96 const G4int detail = verbosity % 10;
97 os << "# Set verbosity with \"/vis/ASCIITree/verbose <verbosity>\":";
98 for (size_t i = 0;
101 }
102 os << "\n# Now printing with verbosity " << verbosity;
103 os << "\n# Format is: PV:n";
104 if (detail >= 1) os << " / LV (SD,RO)";
105 if (detail >= 2) os << " / Solid(type)";
106 if (detail >= 3) os << ", volume, density";
107 if (detail >= 5) os << ", daughter-subtracted volume and mass";
108 os <<
109 "\n# Abbreviations: PV = Physical Volume, LV = Logical Volume,"
110 "\n# SD = Sensitive Detector, RO = Read Out Geometry.";
111}
static std::vector< G4String > fVerbosityGuidance

Referenced by BeginModeling().

Member Data Documentation

◆ fLastCopyNo

G4int G4ASCIITreeSceneHandler::fLastCopyNo
protected

Definition at line 69 of file G4ASCIITreeSceneHandler.hh.

Referenced by EndModeling(), and RequestPrimitives().

◆ fLastNonSequentialCopyNo

G4int G4ASCIITreeSceneHandler::fLastNonSequentialCopyNo
protected

Definition at line 70 of file G4ASCIITreeSceneHandler.hh.

Referenced by EndModeling(), and RequestPrimitives().

◆ fLastPVName

G4String G4ASCIITreeSceneHandler::fLastPVName
protected

Definition at line 68 of file G4ASCIITreeSceneHandler.hh.

Referenced by EndModeling(), and RequestPrimitives().

◆ fLVSet

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

Definition at line 72 of file G4ASCIITreeSceneHandler.hh.

Referenced by EndModeling(), and RequestPrimitives().

◆ fOutFile

std::ofstream G4ASCIITreeSceneHandler::fOutFile
protected

Definition at line 63 of file G4ASCIITreeSceneHandler.hh.

Referenced by BeginModeling(), and EndModeling().

◆ fpLastPV

const G4VPhysicalVolume* G4ASCIITreeSceneHandler::fpLastPV
protected

Definition at line 67 of file G4ASCIITreeSceneHandler.hh.

Referenced by EndModeling(), and RequestPrimitives().

◆ fpOutFile

std::ostream* G4ASCIITreeSceneHandler::fpOutFile
protected

Definition at line 62 of file G4ASCIITreeSceneHandler.hh.

Referenced by BeginModeling(), EndModeling(), and RequestPrimitives().

◆ fReplicaSet

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

Definition at line 76 of file G4ASCIITreeSceneHandler.hh.

Referenced by EndModeling(), and RequestPrimitives().

◆ fRestOfLine

std::ostringstream G4ASCIITreeSceneHandler::fRestOfLine
protected

Definition at line 64 of file G4ASCIITreeSceneHandler.hh.

Referenced by EndModeling(), and RequestPrimitives().


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