Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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 99076 2016-09-01 12:40:47Z 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  static G4bool firstTime = true;
81  if (firstTime) {
82  firstTime = false;
83  G4cout << "G4ASCIITreeSceneHandler::BeginModeling: writing to ";
84  if (outFileName == "G4cout") {
85  G4cout << "G4 standard output (G4cout)";
86  } else {
87  G4cout << "file \"" << outFileName << "\"";
88  }
89  G4cout << G4endl;
90 
92  }
93 
94  if (outFileName != "G4cout") {
95  WriteHeader (fOutFile); fOutFile << std::endl;
96  }
97 }
98 
99 void G4ASCIITreeSceneHandler::WriteHeader (std::ostream& os)
100 {
101  const G4ASCIITree* pSystem = (G4ASCIITree*)GetGraphicsSystem();
102  const G4int verbosity = pSystem->GetVerbosity();
103  const G4int detail = verbosity % 10;
104  os << "# Set verbosity with \"/vis/ASCIITree/verbose <verbosity>\":";
105  for (size_t i = 0;
107  os << "\n# " << G4ASCIITreeMessenger::fVerbosityGuidance[i];
108  }
109  os << "\n# Now printing with verbosity " << verbosity;
110  os << "\n# Format is: PV:n";
111  if (detail >= 1) os << " / LV (SD,RO)";
112  if (detail >= 2) os << " / Solid(type)";
113  if (detail >= 3) os << ", volume, density";
114  if (detail >= 5) os << ", daughter-subtracted volume and mass";
115  if (detail >= 6) os << ", physical volume dump";
116  if (detail >= 7) os << ", polyhedron dump";
117  os <<
118  "\n# Abbreviations: PV = Physical Volume, LV = Logical Volume,"
119  "\n# SD = Sensitive Detector, RO = Read Out Geometry.";
120 }
121 
123  const G4ASCIITree* pSystem = (G4ASCIITree*) GetGraphicsSystem();
124  const G4int verbosity = pSystem->GetVerbosity();
125  const G4int detail = verbosity % 10;
126  const G4String& outFileName = pSystem -> GetOutFileName();
127 
128  // Output left over copy number, if any...
130  if (fLastCopyNo == fLastNonSequentialCopyNo + 1) *fpOutFile << ',';
131  else *fpOutFile << '-';
133  }
134  // Output outstanding rest of line, if any...
135  if (!fRestOfLine.str().empty()) *fpOutFile << fRestOfLine.str();
136  fRestOfLine.str("");
137  fpLastPV = 0;
138  fLastPVName.clear();
139  fLastCopyNo = -99;
141 
142  // This detail to G4cout regardless of outFileName...
143  if (detail >= 4) {
144  G4cout << "Calculating mass(es)..." << G4endl;
145  const std::vector<G4Scene::Model>& models = fpScene->GetRunDurationModelList();
146  std::vector<G4Scene::Model>::const_iterator i;
147  for (i = models.begin(); i != models.end(); ++i) {
148  G4PhysicalVolumeModel* pvModel =
149  dynamic_cast<G4PhysicalVolumeModel*>(i->fpModel);
150  if (pvModel) {
151  const G4ModelingParameters* tempMP =
152  pvModel->GetModelingParameters();
153  G4ModelingParameters mp; // Default - no culling.
154  pvModel->SetModelingParameters (&mp);
155  G4PhysicalVolumeMassScene massScene(pvModel);
156  pvModel->DescribeYourselfTo (massScene);
157  G4double volume = massScene.GetVolume();
158  G4double mass = massScene.GetMass();
159 
160  G4cout << "Overall volume of \""
161  << pvModel->GetTopPhysicalVolume()->GetName()
162  << "\":"
163  << pvModel->GetTopPhysicalVolume()->GetCopyNo()
164  << ", is "
165  << G4BestUnit (volume, "Volume")
166  << " and the daughter-included mass";
167  G4int requestedDepth = pvModel->GetRequestedDepth();
168  if (requestedDepth == G4PhysicalVolumeModel::UNLIMITED) {
169  G4cout << " to unlimited depth";
170  } else {
171  G4cout << ", ignoring daughters at depth "
172  << requestedDepth
173  << " and below,";
174  }
175  G4cout << " is " << G4BestUnit (mass, "Mass")
176  << G4endl;
177 
178  pvModel->SetModelingParameters (tempMP);
179  }
180  }
181  }
182 
183  if (outFileName != "G4cout") {
184  fOutFile.close();
185  G4cout << "Output file \"" << outFileName << "\" closed." << G4endl;
186  }
187  fLVSet.clear();
188  fReplicaSet.clear();
189  G4cout << "G4ASCIITreeSceneHandler::EndModeling" << G4endl;
190  G4VTreeSceneHandler::EndModeling (); // To re-use "culling off" code.
191 }
192 
194 
195  G4PhysicalVolumeModel* pPVModel =
196  dynamic_cast<G4PhysicalVolumeModel*>(fpModel);
197  if (!pPVModel) return; // Not from a G4PhysicalVolumeModel.
198 
199  // This call comes from a G4PhysicalVolumeModel. drawnPVPath is
200  // the path of the current drawn (non-culled) volume in terms of
201  // drawn (non-culled) ancesters. Each node is identified by a
202  // PVNodeID object, which is a physical volume and copy number. It
203  // is a vector of PVNodeIDs corresponding to the geometry hierarchy
204  // actually selected, i.e., not culled.
205  // The following typedef's already set in header file...
206  //typedef G4PhysicalVolumeModel::G4PhysicalVolumeNodeID PVNodeID;
207  //typedef std::vector<PVNodeID> PVPath;
208  const PVPath& drawnPVPath = pPVModel->GetDrawnPVPath();
209  //G4int currentDepth = pPVModel->GetCurrentDepth();
210  G4VPhysicalVolume* pCurrentPV = pPVModel->GetCurrentPV();
211  G4LogicalVolume* pCurrentLV = pPVModel->GetCurrentLV();
212  G4Material* pCurrentMaterial = pPVModel->GetCurrentMaterial();
213 
215  G4int verbosity = pSystem->GetVerbosity();
216  G4int detail = verbosity % 10;
217 
218  if (verbosity < 10 && pCurrentPV->IsReplicated()) {
219  // See if this has been a replica found before with same mother LV...
220  PVPath::const_reverse_iterator thisID = drawnPVPath.rbegin();
221  PVPath::const_reverse_iterator motherID = ++drawnPVPath.rbegin();
222  G4bool ignore = false;
223  for (ReplicaSetIterator i = fReplicaSet.begin(); i != fReplicaSet.end();
224  ++i) {
225  if (i->back().GetPhysicalVolume()->GetLogicalVolume() ==
226  thisID->GetPhysicalVolume()->GetLogicalVolume()) {
227  // For each one previously found (if more than one, they must
228  // have different mothers)...
229  // To avoid compilation errors on VC++ .Net 7.1...
230  // Previously:
231  // PVNodeID previousMotherID = ++(i->rbegin());
232  // (Should that have been: PVNodeID::const_iterator previousMotherID?)
233  // Replace
234  // previousMotherID == i->rend()
235  // by
236  // i->size() <= 1
237  // Replace
238  // previousMotherID != i->rend()
239  // by
240  // i->size() > 1
241  // Replace
242  // previousMotherID->
243  // by
244  // (*i)[i->size() - 2].
245  if (motherID == drawnPVPath.rend() &&
246  i->size() <= 1)
247  ignore = true; // Both have no mother.
248  if (motherID != drawnPVPath.rend() &&
249  i->size() > 1 &&
250  motherID->GetPhysicalVolume()->GetLogicalVolume() ==
251  (*i)[i->size() - 2].GetPhysicalVolume()->GetLogicalVolume())
252  ignore = true; // Same mother LV...
253  }
254  }
255  if (ignore) {
256  pPVModel->CurtailDescent();
257  return;
258  }
259  }
260 
261  const G4String& currentPVName = pCurrentPV->GetName();
262  const G4int currentCopyNo = pCurrentPV->GetCopyNo();
263 
264  if (verbosity < 10 &&
265  currentPVName == fLastPVName &&
266  currentCopyNo != fLastCopyNo) {
267  // For low verbosity, trap series of G4PVPlacements with the same name but
268  // different copy number (replicas should have been taken out above)...
269  // Check...
270  if (pCurrentPV->IsReplicated()) {
271  G4Exception("G4ASCIITreeSceneHandler::RequestPrimitives",
272  "vistree0001",
273  JustWarning,
274  "Replica unexpected");
275  }
276  // Check mothers are identical...
277  else if (pCurrentLV == (fpLastPV? fpLastPV->GetLogicalVolume(): 0)) {
278  if (currentCopyNo != fLastCopyNo + 1) {
279  // Non-sequential copy number...
280  *fpOutFile << ',' << currentCopyNo;
281  fLastNonSequentialCopyNo = currentCopyNo;
282  }
283  fLastCopyNo = currentCopyNo;
284  pPVModel->CurtailDescent();
285  return;
286  }
287  }
288  fpLastPV = pCurrentPV;
289 
290  // High verbosity or a new G4VPhysicalVolume...
291  // Output copy number, if any, from previous invocation...
293  if (fLastCopyNo == fLastNonSequentialCopyNo + 1) *fpOutFile << ',';
294  else *fpOutFile << '-';
296  }
297  // Output rest of line, if any, from previous invocation...
298  if (!fRestOfLine.str().empty()) *fpOutFile << fRestOfLine.str();
299  fRestOfLine.str("");
300  fLastPVName = currentPVName;
301  fLastCopyNo = currentCopyNo;
302  fLastNonSequentialCopyNo = currentCopyNo;
303  // Indent according to drawn path depth...
304  for (size_t i = 0; i < drawnPVPath.size(); i++ ) *fpOutFile << " ";
305  // Start next line...
306  *fpOutFile << "\"" << currentPVName
307  << "\":" << currentCopyNo;
308 
309  if (pCurrentPV->IsReplicated()) {
310  if (verbosity < 10) {
311  // Add printing for replicas (when replicas are ignored)...
312  EAxis axis;
313  G4int nReplicas;
314  G4double width;
315  G4double offset;
316  G4bool consuming;
317  pCurrentPV->GetReplicationData(axis,nReplicas,width,offset,consuming);
318  G4VPVParameterisation* pP = pCurrentPV->GetParameterisation();
319  if (pP) {
320  if (detail < 3) {
321  fReplicaSet.insert(drawnPVPath);
322  if (nReplicas > 2) fRestOfLine << '-';
323  else fRestOfLine << ',';
324  fRestOfLine << nReplicas - 1
325  << " (" << nReplicas << " parametrised volumes)";
326  }
327  }
328  else {
329  fReplicaSet.insert(drawnPVPath);
330  if (nReplicas > 2) fRestOfLine << '-';
331  else fRestOfLine << ',';
332  fRestOfLine << nReplicas - 1
333  << " (" << nReplicas << " replicas)";
334  }
335  }
336  } else {
337  if (fLVSet.find(pCurrentLV) != fLVSet.end()) {
338  if (verbosity < 10) {
339  // Add printing for repeated LV (if it has daughters)...
340  if (pCurrentLV->GetNoDaughters()) fRestOfLine << " (repeated LV)";
341  // ...and curtail descent.
342  pPVModel->CurtailDescent();
343  }
344  }
345  }
346 
347  if (detail >= 1) {
348  fRestOfLine << " / \""
349  << pCurrentLV->GetName() << "\"";
350  G4VSensitiveDetector* sd = pCurrentLV->GetSensitiveDetector();
351  if (sd) {
352  fRestOfLine << " (SD=\"" << sd->GetName() << "\"";
353  G4VReadOutGeometry* roGeom = sd->GetROgeometry();
354  if (roGeom) {
355  fRestOfLine << ",RO=\"" << roGeom->GetName() << "\"";
356  }
357  fRestOfLine << ")";
358  }
359  }
360 
361  if (detail >= 2) {
362  fRestOfLine << " / \""
363  << solid.GetName()
364  << "\"("
365  << solid.GetEntityType() << ")";
366  }
367 
368  if (detail >= 3) {
369  fRestOfLine << ", "
370  << G4BestUnit(((G4VSolid&)solid).GetCubicVolume(),"Volume")
371  << ", ";
372  if (pCurrentMaterial) {
374  << G4BestUnit(pCurrentMaterial->GetDensity(), "Volumic Mass")
375  << " (" << pCurrentMaterial->GetName() << ")";
376  } else {
377  fRestOfLine << "(No material)";
378  }
379  }
380 
381  if (detail >= 5) {
382  if (pCurrentMaterial) {
383  G4Material* pMaterial = const_cast<G4Material*>(pCurrentMaterial);
384  // ...and find daughter-subtracted mass...
385  G4double daughter_subtracted_mass = pCurrentLV->GetMass
386  (pCurrentPV->IsParameterised(), // Force if parametrised.
387  false, // Do not propagate - calculate for this volume minus
388  // volume of daughters.
389  pMaterial);
390  G4double daughter_subtracted_volume =
391  daughter_subtracted_mass / pCurrentMaterial->GetDensity();
392  fRestOfLine << ", "
393  << G4BestUnit(daughter_subtracted_volume,"Volume")
394  << ", "
395  << G4BestUnit(daughter_subtracted_mass,"Mass");
396  }
397  }
398 
399  if (detail >= 6) {
400  std::vector<G4AttValue>* attValues = pPVModel->CreateCurrentAttValues();
401  const std::map<G4String,G4AttDef>* attDefs = pPVModel->GetAttDefs();
402  fRestOfLine << '\n' << G4AttCheck(attValues, attDefs);
403  delete attValues;
404  }
405 
406  if (detail >= 7) {
407  G4Polyhedron* polyhedron = solid.GetPolyhedron();
408  fRestOfLine << "\nLocal polyhedron coordinates:\n" << *polyhedron;
409  G4Transform3D* transform = pPVModel->GetCurrentTransform();
410  polyhedron->Transform(*transform);
411  fRestOfLine << "\nGlobal polyhedron coordinates:\n" << *polyhedron;
412  }
413 
414  if (fLVSet.find(pCurrentLV) == fLVSet.end()) {
415  fLVSet.insert(pCurrentLV); // Record new logical volume.
416  }
417 
418  fRestOfLine << std::endl;
419 
420  return;
421 }
virtual G4Polyhedron * GetPolyhedron() const
Definition: G4VSolid.cc:665
G4String GetName() const
const XML_Char * name
Definition: expat.h:151
const std::map< G4String, G4AttDef > * GetAttDefs() const
virtual G4bool IsReplicated() const =0
const G4String & GetName() const
Definition: G4Material.hh:178
G4Material * GetCurrentMaterial() const
G4double GetDensity() const
Definition: G4Material.hh:180
const G4ModelingParameters * GetModelingParameters() const
#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
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
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
virtual void BeginModeling()
void DescribeYourselfTo(G4VGraphicsScene &)