Geant4  10.02
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: G4ASCIITreeSceneHandler.cc 88761 2015-03-09 12:23:46Z gcosmo $
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"
37 #include "G4ASCIITreeMessenger.hh"
38 #include "G4VSolid.hh"
39 #include "G4PhysicalVolumeModel.hh"
40 #include "G4VPhysicalVolume.hh"
41 #include "G4LogicalVolume.hh"
42 #include "G4VPVParameterisation.hh"
43 #include "G4Polyhedron.hh"
44 #include "G4UnitsTable.hh"
45 #include "G4Material.hh"
46 #include "G4Scene.hh"
47 #include "G4ModelingParameters.hh"
49 #include "G4VSensitiveDetector.hh"
50 #include "G4VReadOutGeometry.hh"
52 #include "G4AttCheck.hh"
53 #include "G4AttValue.hh"
54 
57  const G4String& name):
58  G4VTreeSceneHandler(system, name),
59  fpOutFile(0),
60  fpLastPV(0),
61  fLastCopyNo(-99),
62  fLastNonSequentialCopyNo(-99)
63 {}
64 
66 
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 }
93 
94 void G4ASCIITreeSceneHandler::WriteHeader (std::ostream& os)
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 }
116 
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 }
191 
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 }
virtual G4Polyhedron * GetPolyhedron() const
Definition: G4VSolid.cc:644
G4String GetName() const
const std::map< G4String, G4AttDef > * GetAttDefs() const
virtual G4bool IsReplicated() const =0
G4String name
Definition: TRTMaterials.hh:40
const G4String & GetName() const
Definition: G4Material.hh:178
G4Material * GetCurrentMaterial() const
G4double GetDensity() const
Definition: G4Material.hh:180
G4Navigator * GetNavigatorForTracking() const
const G4ModelingParameters * GetModelingParameters() const
#define width
#define G4BestUnit(a, b)
#define G4_USE_G4BESTUNIT_FOR_VERBOSE 1
G4ASCIITreeSceneHandler(G4VGraphicsSystem &system, const G4String &name)
virtual G4GeometryType GetEntityType() const =0
const std::vector< G4PhysicalVolumeNodeID > & GetDrawnPVPath() const
int G4int
Definition: G4Types.hh:78
void SetModelingParameters(const G4ModelingParameters *)
std::vector< PVNodeID > PVPath
static std::vector< G4String > fVerbosityGuidance
virtual void RequestPrimitives(const G4VSolid &)
std::set< PVPath >::iterator ReplicaSetIterator
G4GLOB_DLL std::ostream G4cout
G4Transform3D * GetCurrentTransform() const
const G4String & GetName() const
bool G4bool
Definition: G4Types.hh:79
const std::vector< Model > & GetRunDurationModelList() const
virtual G4VPVParameterisation * GetParameterisation() const =0
G4String GetName() const
G4VReadOutGeometry * GetROgeometry() const
G4VGraphicsSystem * GetGraphicsSystem() const
HepGeom::Transform3D G4Transform3D
G4int GetNoDaughters() const
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
G4VPhysicalVolume * GetTopPhysicalVolume() const
virtual G4bool IsParameterised() const =0
static G4TransportationManager * GetTransportationManager()
G4int GetVerbosity() const
Definition: G4ASCIITree.hh:47
G4LogicalVolume * GetLogicalVolume() const
std::vector< G4AttValue > * CreateCurrentAttValues() const
virtual void EndModeling()
EAxis
Definition: geomdefs.hh:54
virtual G4int GetCopyNo() const =0
const G4VPhysicalVolume * fpLastPV
virtual void GetReplicationData(EAxis &axis, G4int &nReplicas, G4double &width, G4double &offset, G4bool &consuming) const =0
#define G4endl
Definition: G4ios.hh:61
void WriteHeader(std::ostream &)
G4double GetMass(G4bool forced=false, G4bool propagate=true, G4Material *parMaterial=0)
const G4String & GetName() const
double G4double
Definition: G4Types.hh:76
G4LogicalVolume * GetCurrentLV() const
G4VPhysicalVolume * GetCurrentPV() const
std::set< G4LogicalVolume * > fLVSet
G4VSensitiveDetector * GetSensitiveDetector() const
G4VPhysicalVolume * GetWorldVolume() const
virtual void BeginModeling()
void DescribeYourselfTo(G4VGraphicsScene &)