Geant4 11.1.1
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 G4Polymarker &)
 
virtual void BeginModeling ()
 
virtual void EndModeling ()
 
virtual void AddPrimitive (const G4Polyline &)=0
 
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 G4Plotter &)
 
- 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 G4Orb &)
 
virtual void AddSolid (const G4Para &)
 
virtual void AddSolid (const G4Sphere &)
 
virtual void AddSolid (const G4Torus &)
 
virtual void AddSolid (const G4Trap &)
 
virtual void AddSolid (const G4Trd &)
 
virtual void AddSolid (const G4Tubs &)
 
virtual void AddSolid (const G4Ellipsoid &)
 
virtual void AddSolid (const G4Polycone &)
 
virtual void AddSolid (const G4Polyhedra &)
 
virtual void AddSolid (const G4TessellatedSolid &)
 
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 AddCompound (const G4THitsMap< G4StatDouble > &)
 
virtual void AddCompound (const G4Mesh &)
 
virtual void BeginModeling ()
 
virtual void EndModeling ()
 
virtual void BeginPrimitives (const G4Transform3D &objectTransformation=G4Transform3D())
 
virtual void EndPrimitives ()
 
virtual void BeginPrimitives2D (const G4Transform3D &objectTransformation=G4Transform3D())
 
virtual void EndPrimitives2D ()
 
virtual void AddPrimitive (const G4Polyline &)=0
 
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 G4Plotter &)
 
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 G4ColourGetColor ()
 
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 *)
 
G4int GetNumberOfCloudPoints (const G4VisAttributes *) const
 
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 ()
 
template<class T >
void AddSolidT (const T &solid)
 
template<class T >
void AddSolidWithAuxiliaryEdges (const T &solid)
 
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 G4Orb &)=0
 
virtual void AddSolid (const G4Para &)=0
 
virtual void AddSolid (const G4Sphere &)=0
 
virtual void AddSolid (const G4Torus &)=0
 
virtual void AddSolid (const G4Trap &)=0
 
virtual void AddSolid (const G4Trd &)=0
 
virtual void AddSolid (const G4Tubs &)=0
 
virtual void AddSolid (const G4Ellipsoid &)=0
 
virtual void AddSolid (const G4Polycone &)=0
 
virtual void AddSolid (const G4Polyhedra &)=0
 
virtual void AddSolid (const G4TessellatedSolid &)=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 AddCompound (const G4THitsMap< G4StatDouble > &)=0
 
virtual void AddCompound (const G4Mesh &)=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 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 G4Plotter &)=0
 
virtual const G4VisExtentGetExtent () const
 

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 G4DisplacedSolidCreateSectionSolid ()
 
virtual G4DisplacedSolidCreateCutawaySolid ()
 
void LoadAtts (const G4Visible &, G4AttHolder *)
 
void StandardSpecialMeshRendering (const G4Mesh &)
 
void Draw3DRectMeshAsDots (const G4Mesh &)
 
void Draw3DRectMeshAsSurfaces (const G4Mesh &)
 
void DrawTetMeshAsDots (const G4Mesh &)
 
void DrawTetMeshAsSurfaces (const G4Mesh &)
 
G4ThreeVector GetPointInBox (const G4ThreeVector &pos, G4double halfX, G4double halfY, G4double halfZ) const
 
G4ThreeVector GetPointInTet (const std::vector< G4ThreeVector > &vertices) const
 

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 47 of file G4ASCIITreeSceneHandler.hh.

Member Typedef Documentation

◆ LVSetIterator

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

Definition at line 72 of file G4ASCIITreeSceneHandler.hh.

◆ PVNodeID

◆ PVPath

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

Definition at line 74 of file G4ASCIITreeSceneHandler.hh.

◆ ReplicaSetIterator

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

Definition at line 76 of file G4ASCIITreeSceneHandler.hh.

◆ ReplicaSetReverseIterator

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

Definition at line 77 of file G4ASCIITreeSceneHandler.hh.

Constructor & Destructor Documentation

◆ G4ASCIITreeSceneHandler()

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

◆ ~G4ASCIITreeSceneHandler()

G4ASCIITreeSceneHandler::~G4ASCIITreeSceneHandler ( )
virtual

Definition at line 64 of file G4ASCIITreeSceneHandler.cc.

64{}

Member Function Documentation

◆ BeginModeling()

void G4ASCIITreeSceneHandler::BeginModeling ( )
virtual

Reimplemented from G4VTreeSceneHandler.

Definition at line 66 of file G4ASCIITreeSceneHandler.cc.

66 {
67
68 G4VTreeSceneHandler::BeginModeling (); // To re-use "culling off" code.
69
70 const G4ASCIITree* pSystem = (G4ASCIITree*)GetGraphicsSystem();
71 const G4String& outFileName = pSystem -> GetOutFileName();
72 if (outFileName == "G4cout") {
74 } else {
75 fOutFile.open (outFileName);
77 }
78
79 static G4bool firstTime = true;
80 if (firstTime) {
81 firstTime = false;
82 G4cout << "G4ASCIITreeSceneHandler::BeginModeling: writing to ";
83 if (outFileName == "G4cout") {
84 G4cout << "G4 standard output (G4cout)";
85 } else {
86 G4cout << "file \"" << outFileName << "\"";
87 }
88 G4cout << G4endl;
89
91 }
92
93 if (outFileName != "G4cout") {
94 WriteHeader (fOutFile); fOutFile << std::endl;
95 }
96}
bool G4bool
Definition: G4Types.hh:86
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
void WriteHeader(std::ostream &)
G4VGraphicsSystem * GetGraphicsSystem() const
virtual void BeginModeling()

◆ EndModeling()

void G4ASCIITreeSceneHandler::EndModeling ( )
virtual

Reimplemented from G4VTreeSceneHandler.

Definition at line 121 of file G4ASCIITreeSceneHandler.cc.

121 {
122 const G4ASCIITree* pSystem = (G4ASCIITree*) GetGraphicsSystem();
123 const G4int verbosity = pSystem->GetVerbosity();
124 const G4int detail = verbosity % 10;
125 const G4String& outFileName = pSystem -> GetOutFileName();
126
127 // Output left over copy number, if any...
130 else *fpOutFile << '-';
132 }
133 // Output outstanding rest of line, if any...
134 if (!fRestOfLine.str().empty()) *fpOutFile << fRestOfLine.str();
135 fRestOfLine.str("");
136 fpLastPV = 0;
137 fLastPVName.clear();
138 fLastCopyNo = -99;
140
141 // This detail to G4cout regardless of outFileName...
142 if (detail >= 4) {
143 G4cout << "Calculating mass(es)..." << G4endl;
144 const std::vector<G4Scene::Model>& models = fpScene->GetRunDurationModelList();
145 std::vector<G4Scene::Model>::const_iterator i;
146 for (i = models.begin(); i != models.end(); ++i) {
147 G4PhysicalVolumeModel* pvModel =
148 dynamic_cast<G4PhysicalVolumeModel*>(i->fpModel);
149 if (pvModel) {
150 const G4ModelingParameters* tempMP =
151 pvModel->GetModelingParameters();
152 G4ModelingParameters mp; // Default - no culling.
153 pvModel->SetModelingParameters (&mp);
154 G4PhysicalVolumeMassScene massScene(pvModel);
155 pvModel->DescribeYourselfTo (massScene);
156 G4double volume = massScene.GetVolume();
157 G4double mass = massScene.GetMass();
158
159 G4cout << "Overall volume of \""
160 << pvModel->GetTopPhysicalVolume()->GetName()
161 << "\":"
162 << pvModel->GetTopPhysicalVolume()->GetCopyNo()
163 << ", is "
164 << G4BestUnit (volume, "Volume")
165 << " and the daughter-included mass";
166 G4int requestedDepth = pvModel->GetRequestedDepth();
167 if (requestedDepth == G4PhysicalVolumeModel::UNLIMITED) {
168 G4cout << " to unlimited depth";
169 } else {
170 G4cout << ", ignoring daughters at depth "
171 << requestedDepth
172 << " and below,";
173 }
174 G4cout << " is " << G4BestUnit (mass, "Mass")
175 << G4endl;
176
177 pvModel->SetModelingParameters (tempMP);
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}
#define G4BestUnit(a, b)
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
std::set< G4LogicalVolume * > fLVSet
G4int GetVerbosity() const
Definition: G4ASCIITree.hh:46
void DescribeYourselfTo(G4VGraphicsScene &)
G4VPhysicalVolume * GetTopPhysicalVolume() const
const std::vector< Model > & GetRunDurationModelList() 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 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) ancestors. 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 const G4String& currentPVName = pCurrentPV->GetName();
211 const G4int currentCopyNo = pCurrentPV->GetCopyNo();
212 G4LogicalVolume* pCurrentLV = pPVModel->GetCurrentLV();
213 G4Material* pCurrentMaterial = pPVModel->GetCurrentMaterial();
214
216 G4int verbosity = pSystem->GetVerbosity();
217 G4int detail = verbosity % 10;
218
219 // If verbosity < 10 suppress unnecessary repeated printing.
220 // Repeated simple replicas can always be suppressed.
221 // Paramaterisations can only be suppressed if verbosity < 3, since their
222 // size, density, etc., are in principle different.
223 const G4bool isParameterised = pCurrentPV->GetParameterisation();
224 const G4bool isSimpleReplica = pCurrentPV->IsReplicated() && !isParameterised;
225 const G4bool isAmenableToSupression =
226 (verbosity < 10 && isSimpleReplica) || (verbosity < 3 && isParameterised);
227 if (isAmenableToSupression) {
228 // See if this has been found before with same mother LV...
229 PVPath::const_reverse_iterator thisID = drawnPVPath.rbegin();
230 PVPath::const_reverse_iterator motherID = ++drawnPVPath.rbegin();
231 G4bool ignore = false;
232 for (ReplicaSetIterator i = fReplicaSet.begin(); i != fReplicaSet.end();
233 ++i) {
234 if (i->back().GetPhysicalVolume()->GetLogicalVolume() ==
235 thisID->GetPhysicalVolume()->GetLogicalVolume()) {
236 // For each one previously found (if more than one, they must
237 // have different mothers)...
238 // To avoid compilation errors on VC++ .Net 7.1...
239 // Previously:
240 // PVNodeID previousMotherID = ++(i->rbegin());
241 // (Should that have been: PVNodeID::const_iterator previousMotherID?)
242 // Replace
243 // previousMotherID == i->rend()
244 // by
245 // i->size() <= 1
246 // Replace
247 // previousMotherID != i->rend()
248 // by
249 // i->size() > 1
250 // Replace
251 // previousMotherID->
252 // by
253 // (*i)[i->size() - 2].
254 if (motherID == drawnPVPath.rend() &&
255 i->size() <= 1)
256 ignore = true; // Both have no mother.
257 if (motherID != drawnPVPath.rend() &&
258 i->size() > 1 &&
259 motherID->GetPhysicalVolume()->GetLogicalVolume() ==
260 (*i)[i->size() - 2].GetPhysicalVolume()->GetLogicalVolume())
261 ignore = true; // Same mother LV...
262 }
263 }
264 if (ignore) {
265 pPVModel->CurtailDescent();
266 return;
267 }
268 }
269
270 // Now suppress printing for volumes with the same name but different
271 // copy number - but not if they are parameterisations that have not been
272 // taken out above (those that are "amenable to suppression" will have been
273 // taken out).
274 if (verbosity < 10 && !isParameterised &&
275 currentPVName == fLastPVName &&
276 currentCopyNo != fLastCopyNo) {
277 // Check...
278 if (isAmenableToSupression) {
279 G4Exception("G4ASCIITreeSceneHandler::RequestPrimitives",
280 "vistree0001",
282 "Volume amenable to suppressed printing unexpected");
283 }
284 // Check mothers are identical...
285 else if (pCurrentLV == (fpLastPV? fpLastPV->GetLogicalVolume(): 0)) {
286 if (currentCopyNo != fLastCopyNo + 1) {
287 // Non-sequential copy number...
288 *fpOutFile << ',' << currentCopyNo;
289 fLastNonSequentialCopyNo = currentCopyNo;
290 }
291 fLastCopyNo = currentCopyNo;
292 pPVModel->CurtailDescent();
293 return;
294 }
295 }
296 fpLastPV = pCurrentPV;
297
298 // High verbosity or a new or otherwise non-amenable volume...
299 // Output copy number, if any, from previous invocation...
302 else *fpOutFile << '-';
304 }
305 // Output rest of line, if any, from previous invocation...
306 if (!fRestOfLine.str().empty()) *fpOutFile << fRestOfLine.str();
307 fRestOfLine.str("");
308 fLastPVName = currentPVName;
309 fLastCopyNo = currentCopyNo;
310 fLastNonSequentialCopyNo = currentCopyNo;
311 // Start next line...
312 // Indent according to drawn path depth...
313 for (size_t i = 0; i < drawnPVPath.size(); i++ ) *fpOutFile << " ";
314 *fpOutFile << "\"" << currentPVName
315 << "\":" << currentCopyNo;
316
317 if (pCurrentPV->IsReplicated()) {
318 if (verbosity < 10) {
319 // Add printing for replicas (when replicas are ignored)...
320 EAxis axis;
321 G4int nReplicas;
322 G4double width;
323 G4double offset;
324 G4bool consuming;
325 pCurrentPV->GetReplicationData(axis,nReplicas,width,offset,consuming);
326 G4VPVParameterisation* pP = pCurrentPV->GetParameterisation();
327 if (pP) {
328 if (detail < 3) {
329 fReplicaSet.insert(drawnPVPath);
330 if (nReplicas > 2) fRestOfLine << '-';
331 else fRestOfLine << ',';
332 fRestOfLine << nReplicas - 1
333 << " (" << nReplicas << " parametrised volumes)";
334 }
335 }
336 else {
337 fReplicaSet.insert(drawnPVPath);
338 if (nReplicas > 2) fRestOfLine << '-';
339 else fRestOfLine << ',';
340 fRestOfLine << nReplicas - 1
341 << " (" << nReplicas << " replicas)";
342 }
343 }
344 } else {
345 if (fLVSet.find(pCurrentLV) != fLVSet.end()) {
346 if (verbosity < 10) {
347 // Add printing for repeated LV (if it has daughters)...
348 if (pCurrentLV->GetNoDaughters()) fRestOfLine << " (repeated LV)";
349 // ...and curtail descent.
350 pPVModel->CurtailDescent();
351 }
352 }
353 }
354
355 if (detail >= 1) {
356 fRestOfLine << " / \""
357 << pCurrentLV->GetName() << "\"";
358 G4VSensitiveDetector* sd = pCurrentLV->GetSensitiveDetector();
359 if (sd) {
360 fRestOfLine << " (SD=\"" << sd->GetName() << "\"";
361 G4VReadOutGeometry* roGeom = sd->GetROgeometry();
362 if (roGeom) {
363 fRestOfLine << ",RO=\"" << roGeom->GetName() << "\"";
364 }
365 fRestOfLine << ")";
366 }
367 }
368
369 if (detail >= 2) {
370 fRestOfLine << " / \""
371 << solid.GetName()
372 << "\"("
373 << solid.GetEntityType() << ")";
374 }
375
376 if (detail >= 3) {
377 fRestOfLine << ", "
378 << G4BestUnit(((G4VSolid&)solid).GetCubicVolume(),"Volume")
379 << ", ";
380 if (pCurrentMaterial) {
382 << G4BestUnit(pCurrentMaterial->GetDensity(), "Volumic Mass")
383 << " (" << pCurrentMaterial->GetName() << ")";
384 } else {
385 fRestOfLine << "(No material)";
386 }
387 }
388
389 if (detail >= 5) {
390 if (pCurrentMaterial) {
391 G4Material* pMaterial = const_cast<G4Material*>(pCurrentMaterial);
392 // ...and find daughter-subtracted mass...
393 G4double daughter_subtracted_mass = pCurrentLV->GetMass
394 (pCurrentPV->IsParameterised(), // Force if parametrised.
395 false, // Do not propagate - calculate for this volume minus
396 // volume of daughters.
397 pMaterial);
398 G4double daughter_subtracted_volume =
399 daughter_subtracted_mass / pCurrentMaterial->GetDensity();
400 fRestOfLine << ", "
401 << G4BestUnit(daughter_subtracted_volume,"Volume")
402 << ", "
403 << G4BestUnit(daughter_subtracted_mass,"Mass");
404 }
405 }
406
407 if (detail >= 6) {
408 std::vector<G4AttValue>* attValues = pPVModel->CreateCurrentAttValues();
409 const std::map<G4String,G4AttDef>* attDefs = pPVModel->GetAttDefs();
410 fRestOfLine << '\n' << G4AttCheck(attValues, attDefs);
411 delete attValues;
412 }
413
414 if (detail >= 7) {
415 G4Polyhedron* polyhedron = solid.GetPolyhedron();
416 fRestOfLine << "\nLocal polyhedron coordinates:\n" << *polyhedron;
417 const G4Transform3D& transform = pPVModel->GetCurrentTransform();
418 polyhedron->Transform(transform);
419 fRestOfLine << "\nGlobal polyhedron coordinates:\n" << *polyhedron;
420 }
421
422 if (fLVSet.find(pCurrentLV) == fLVSet.end()) {
423 fLVSet.insert(pCurrentLV); // Record new logical volume.
424 }
425
426 fRestOfLine << std::endl;
427
428 return;
429}
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
std::vector< PVNodeID > PVPath
std::set< PVPath >::iterator ReplicaSetIterator
G4double GetMass(G4bool forced=false, G4bool propagate=true, G4Material *parMaterial=nullptr)
G4VSensitiveDetector * GetSensitiveDetector() const
std::size_t GetNoDaughters() const
const G4String & GetName() const
G4double GetDensity() const
Definition: G4Material.hh:175
const G4String & GetName() const
Definition: G4Material.hh:172
const std::vector< G4PhysicalVolumeNodeID > & GetDrawnPVPath() const
const G4Transform3D & GetCurrentTransform() const
G4VPhysicalVolume * GetCurrentPV() const
std::vector< G4AttValue > * CreateCurrentAttValues() const
G4LogicalVolume * GetCurrentLV() const
G4Material * GetCurrentMaterial() const
const std::map< G4String, G4AttDef > * GetAttDefs() 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 G4Polyhedron * GetPolyhedron() const
Definition: G4VSolid.cc:705
virtual G4GeometryType GetEntityType() const =0
EAxis
Definition: geomdefs.hh:54

◆ WriteHeader()

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

Definition at line 98 of file G4ASCIITreeSceneHandler.cc.

99{
100 const G4ASCIITree* pSystem = (G4ASCIITree*)GetGraphicsSystem();
101 const G4int verbosity = pSystem->GetVerbosity();
102 const G4int detail = verbosity % 10;
103 os << "# Set verbosity with \"/vis/ASCIITree/verbose <verbosity>\":";
104 for (size_t i = 0;
107 }
108 os << "\n# Now printing with verbosity " << verbosity;
109 os << "\n# Format is: PV:n";
110 if (detail >= 1) os << " / LV (SD,RO)";
111 if (detail >= 2) os << " / Solid(type)";
112 if (detail >= 3) os << ", volume, density";
113 if (detail >= 5) os << ", daughter-subtracted volume and mass";
114 if (detail >= 6) os << ", physical volume dump";
115 if (detail >= 7) os << ", polyhedron dump";
116 os <<
117 "\n# Abbreviations: PV = Physical Volume, LV = Logical Volume,"
118 "\n# SD = Sensitive Detector, RO = Read Out Geometry.";
119}
static std::vector< G4String > fVerbosityGuidance

Referenced by BeginModeling().

Member Data Documentation

◆ fLastCopyNo

G4int G4ASCIITreeSceneHandler::fLastCopyNo
protected

Definition at line 68 of file G4ASCIITreeSceneHandler.hh.

Referenced by EndModeling(), and RequestPrimitives().

◆ fLastNonSequentialCopyNo

G4int G4ASCIITreeSceneHandler::fLastNonSequentialCopyNo
protected

Definition at line 69 of file G4ASCIITreeSceneHandler.hh.

Referenced by EndModeling(), and RequestPrimitives().

◆ fLastPVName

G4String G4ASCIITreeSceneHandler::fLastPVName
protected

Definition at line 67 of file G4ASCIITreeSceneHandler.hh.

Referenced by EndModeling(), and RequestPrimitives().

◆ fLVSet

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

Definition at line 71 of file G4ASCIITreeSceneHandler.hh.

Referenced by EndModeling(), and RequestPrimitives().

◆ fOutFile

std::ofstream G4ASCIITreeSceneHandler::fOutFile
protected

Definition at line 62 of file G4ASCIITreeSceneHandler.hh.

Referenced by BeginModeling(), and EndModeling().

◆ fpLastPV

const G4VPhysicalVolume* G4ASCIITreeSceneHandler::fpLastPV
protected

Definition at line 66 of file G4ASCIITreeSceneHandler.hh.

Referenced by EndModeling(), and RequestPrimitives().

◆ fpOutFile

std::ostream* G4ASCIITreeSceneHandler::fpOutFile
protected

Definition at line 61 of file G4ASCIITreeSceneHandler.hh.

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

◆ fReplicaSet

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

Definition at line 75 of file G4ASCIITreeSceneHandler.hh.

Referenced by EndModeling(), and RequestPrimitives().

◆ fRestOfLine

std::ostringstream G4ASCIITreeSceneHandler::fRestOfLine
protected

Definition at line 63 of file G4ASCIITreeSceneHandler.hh.

Referenced by EndModeling(), and RequestPrimitives().


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