Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4ASCIITreeSceneHandler.cc
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26//
27// $Id$
28//
29//
30// John Allison 5th April 2001
31// A scene handler to dump geometry hierarchy.
32// Based on a provisional G4ASCIITreeGraphicsScene (was in modeling).
33
35
36#include "G4ASCIITree.hh"
38#include "G4VSolid.hh"
40#include "G4VPhysicalVolume.hh"
41#include "G4LogicalVolume.hh"
43#include "G4Polyhedron.hh"
44#include "G4UnitsTable.hh"
45#include "G4Material.hh"
46#include "G4Scene.hh"
50#include "G4VReadOutGeometry.hh"
52
54(G4VGraphicsSystem& system,
55 const G4String& name):
56 G4VTreeSceneHandler(system, name),
57 fpOutFile(0),
58 fpLastPV(0),
59 fLastCopyNo(-99),
60 fLastNonSequentialCopyNo(-99)
61{}
62
64
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}
91
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}
112
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}
187
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
#define G4BestUnit(a, b)
#define G4_USE_G4BESTUNIT_FOR_VERBOSE 1
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
static std::vector< G4String > fVerbosityGuidance
const G4VPhysicalVolume * fpLastPV
std::vector< PVNodeID > PVPath
G4ASCIITreeSceneHandler(G4VGraphicsSystem &system, const G4String &name)
void WriteHeader(std::ostream &)
std::set< PVPath >::iterator ReplicaSetIterator
virtual void RequestPrimitives(const G4VSolid &)
std::set< G4LogicalVolume * > fLVSet
G4int GetVerbosity() const
Definition: G4ASCIITree.hh:47
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
G4VPhysicalVolume * GetWorldVolume() const
const std::vector< G4PhysicalVolumeNodeID > & GetDrawnPVPath() const
G4VPhysicalVolume * GetCurrentPV() const
void DescribeYourselfTo(G4VGraphicsScene &)
G4LogicalVolume * GetCurrentLV() const
G4VPhysicalVolume * GetTopPhysicalVolume() const
G4Material * GetCurrentMaterial() const
const std::vector< Model > & GetRunDurationModelList() const
static G4TransportationManager * GetTransportationManager()
G4Navigator * GetNavigatorForTracking() const
const G4ModelingParameters * GetModelingParameters() const
void SetModelingParameters(const G4ModelingParameters *)
virtual G4bool IsReplicated() const =0
G4LogicalVolume * GetLogicalVolume() const
virtual void GetReplicationData(EAxis &axis, G4int &nReplicas, G4double &width, G4double &offset, G4bool &consuming) const =0
virtual G4int GetCopyNo() const =0
const G4String & GetName() const
virtual G4VPVParameterisation * GetParameterisation() const =0
virtual G4bool IsParameterised() const =0
G4String GetName() const
G4VGraphicsSystem * GetGraphicsSystem() const
G4VReadOutGeometry * GetROgeometry() const
G4String GetName() const
virtual G4GeometryType GetEntityType() const =0
virtual void EndModeling()
virtual void BeginModeling()
EAxis
Definition: geomdefs.hh:54
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41