CMS 3D CMS Logo

PrintG4Solids.cc
Go to the documentation of this file.
4 
8 
9 #include "G4Box.hh"
10 #include "G4Cons.hh"
11 #include "G4ExtrudedSolid.hh"
12 #include "G4Polycone.hh"
13 #include "G4Polyhedra.hh"
14 #include "G4Trap.hh"
15 #include "G4Trd.hh"
16 #include "G4Tubs.hh"
17 #include "G4LogicalVolumeStore.hh"
18 #include "G4LogicalVolume.hh"
19 #include "G4VSolid.hh"
20 #include "G4NavigationHistory.hh"
21 #include "G4PhysicalVolumeStore.hh"
22 #include "G4TransportationManager.hh"
23 #include "G4VPhysicalVolume.hh"
24 
25 #include <algorithm>
26 #include <map>
27 #include <set>
28 #include <sstream>
29 #include <string>
30 #include <vector>
31 
33 
34 class PrintG4Solids : public SimWatcher, public Observer<const BeginOfRun *> {
35 public:
37  ~PrintG4Solids() override = default;
38 
39 private:
40  void update(const BeginOfRun *run) override;
41  void dumpSummary(std::ostream &out = G4cout);
42  G4VPhysicalVolume *getTopPV();
43 
44 private:
45  G4VPhysicalVolume *theTopPV_;
46 };
47 
49  G4cout << "PrintG4Solids:: initialised for printing information about G4VSolids" << G4endl;
50 }
51 
53  //Now take action
54  theTopPV_ = getTopPV();
55 
57 }
58 
59 void PrintG4Solids::dumpSummary(std::ostream &out) {
60  //---------- Dump number of objects of each class
61  out << " @@@@@@@@@@@@@@@@@@ Dumping G4 geometry objects Summary " << G4endl;
62  if (theTopPV_ == nullptr) {
63  out << " No volume created " << G4endl;
64  return;
65  }
66  out << " @@@ Geometry built inside world volume: " << theTopPV_->GetName() << G4endl;
67  // Get number of solids (< # LV if several LV share a solid)
68  const G4LogicalVolumeStore *lvs = G4LogicalVolumeStore::GetInstance();
69  std::vector<G4LogicalVolume *>::const_iterator lvcite;
70  std::set<G4VSolid *> theSolids;
71  for (lvcite = lvs->begin(); lvcite != lvs->end(); lvcite++)
72  theSolids.insert((*lvcite)->GetSolid());
73  out << " Number of G4VSolid's: " << theSolids.size() << G4endl;
74  std::set<G4VSolid *>::const_iterator solid;
75  for (solid = theSolids.begin(); solid != theSolids.end(); solid++) {
76  G4String type = (*solid)->GetEntityType();
77  out << (*solid)->GetName() << ":" << type << " Volume " << (*solid)->GetCubicVolume();
78  if (type == "G4Box") {
79  const G4Box *box = static_cast<const G4Box *>(*solid);
80  out << " dx:dy:dz " << box->GetXHalfLength() << ":" << box->GetYHalfLength() << ":" << box->GetZHalfLength();
81  } else if (type == "G4Tubs") {
82  const G4Tubs *tube = static_cast<const G4Tubs *>(*solid);
83  out << " rin:rout:dz:phistart:dphi " << tube->GetInnerRadius() << ":" << tube->GetOuterRadius() << ":"
84  << tube->GetZHalfLength() << ":" << convertRadToDeg(tube->GetStartPhiAngle()) << ":"
85  << convertRadToDeg(tube->GetDeltaPhiAngle());
86  } else if (type == "G4Cons") {
87  const G4Cons *cone = static_cast<const G4Cons *>(*solid);
88  out << " rinminus:routminus:rinplus:routplus:dz:phistart:dphi " << cone->GetInnerRadiusMinusZ() << ":"
89  << cone->GetOuterRadiusMinusZ() << ":" << cone->GetInnerRadiusPlusZ() << ":" << cone->GetOuterRadiusPlusZ()
90  << ":" << cone->GetZHalfLength() << ":" << convertRadToDeg(cone->GetStartPhiAngle()) << ":"
91  << convertRadToDeg(cone->GetDeltaPhiAngle());
92  } else if (type == "G4Trap") {
93  const G4Trap *trap = static_cast<const G4Trap *>(*solid);
94  out << " zhalf:yl1:xl11:xl12:tana1:yl2:xl21:xl22:tana2 " << trap->GetZHalfLength() << ":"
95  << trap->GetYHalfLength1() << ":" << trap->GetXHalfLength1() << ":" << trap->GetXHalfLength2() << ":"
96  << trap->GetTanAlpha1() << ":" << trap->GetYHalfLength2() << ":" << trap->GetXHalfLength3() << ":"
97  << trap->GetXHalfLength4() << ":" << trap->GetTanAlpha2();
98  } else if (type == "G4Trd") {
99  const G4Trd *trd = static_cast<const G4Trd *>(*solid);
100  out << " xl1:xl2:yl1:yl2:zhalf " << trd->GetXHalfLength1() << ":" << trd->GetXHalfLength2() << ":"
101  << trd->GetYHalfLength1() << ":" << trd->GetYHalfLength2() << ":" << trd->GetZHalfLength();
102  } else if (type == "G4Polycone") {
103  const G4Polycone *cone = static_cast<const G4Polycone *>(*solid);
104  const auto hist = cone->GetOriginalParameters();
105  int num = hist->Num_z_planes;
106  out << " angle " << convertRadToDeg(hist->Start_angle) << ":" << convertRadToDeg(hist->Opening_angle) << " with "
107  << num << " planes:";
108  for (int k = 0; k < num; ++k)
109  out << " [" << k << "] " << hist->Z_values[k] << ":" << hist->Rmin[k] << ":" << hist->Rmax[k];
110  } else if (type == "G4Polyhedra") {
111  const G4Polyhedra *pgon = static_cast<const G4Polyhedra *>(*solid);
112  const auto hist = pgon->GetOriginalParameters();
113  int num = hist->Num_z_planes;
114  out << " angle " << convertRadToDeg(hist->Start_angle) << ":" << convertRadToDeg(hist->Opening_angle) << " with "
115  << hist->numSide << " sides and " << num << " planes:";
116  for (int k = 0; k < num; ++k)
117  out << " [" << k << "] " << hist->Z_values[k] << ":" << hist->Rmin[k] << ":" << hist->Rmax[k];
118  } else if (type == "G4ExtrudedSolid") {
119  const G4ExtrudedSolid *pgon = static_cast<const G4ExtrudedSolid *>(*solid);
120  int vert = pgon->GetNofVertices();
121  int numz = pgon->GetNofZSections();
122  out << " " << vert << " vertices:";
123  for (int k = 0; k < vert; ++k)
124  out << " [" << k << "] " << pgon->GetVertex(k);
125  out << "; and " << numz << " z-sections:";
126  for (int k = 0; k < numz; ++k) {
127  const auto &zsec = pgon->GetZSection(k);
128  out << " [" << k << "] " << zsec.fZ << ":" << zsec.fScale << ":" << zsec.fOffset;
129  }
130  }
131  out << G4endl;
132  }
133 }
134 
135 G4VPhysicalVolume *PrintG4Solids::getTopPV() {
136  return G4TransportationManager::GetTransportationManager()->GetNavigatorForTracking()->GetWorldVolume();
137 }
138 
141 
#define DEFINE_SIMWATCHER(type)
void dumpSummary(std::ostream &out=G4cout)
constexpr NumType convertRadToDeg(NumType radians)
Definition: angle_units.h:21
G4VPhysicalVolume * theTopPV_
~PrintG4Solids() override=default
PrintG4Solids(edm::ParameterSet const &p)
void update(const BeginOfRun *run) override
This routine will be called when the appropriate signal arrives.
G4VPhysicalVolume * getTopPV()