CMS 3D CMS Logo

CheckOverlap.cc
Go to the documentation of this file.
2 
4 
9 
10 #include "G4Run.hh"
11 #include "G4LogicalVolumeStore.hh"
12 #include "G4PVPlacement.hh"
13 #include "G4PVParameterised.hh"
14 #include "G4LogicalVolume.hh"
15 #include "G4Material.hh"
16 #include "G4TransportationManager.hh"
17 
18 #include <set>
19 
21  std::vector<std::string> defNames;
22  nodeNames = p.getUntrackedParameter<std::vector<std::string> >("NodeNames", defNames);
23  nPoints = p.getUntrackedParameter<int>("Resolution", 1000);
24  G4cout << "CheckOverlap:: initialised with " << nodeNames.size() << " Node Names and Resolution " << nPoints
25  << " the names are:" << G4endl;
26  for (unsigned int ii = 0; ii < nodeNames.size(); ii++)
27  G4cout << "CheckOverlap:: Node[" << ii << "] : " << nodeNames[ii] << G4endl;
28 }
29 
31 
33  if (!nodeNames.empty()) {
34  const G4LogicalVolumeStore* lvs = G4LogicalVolumeStore::GetInstance();
35  G4cout << "CheckOverlap::update nLV= " << lvs->size() << G4endl;
36  std::vector<G4LogicalVolume*>::const_iterator lvcite;
37  int i = 0;
38  for (lvcite = lvs->begin(); lvcite != lvs->end(); lvcite++) {
39  for (unsigned int ii = 0; ii < nodeNames.size(); ii++) {
40  if ((*lvcite)->GetName() == (G4String)(nodeNames[ii])) {
41  topLV.push_back((*lvcite));
42  break;
43  }
44  }
45  G4cout << "Name of node " << (++i) << " : " << (*lvcite)->GetName() << G4endl;
46  if (topLV.size() == nodeNames.size())
47  break;
48  }
49  } else {
50  G4VPhysicalVolume* theTopPV = getTopPV();
51  topLV.push_back(theTopPV->GetLogicalVolume());
52  }
53 
54  if (topLV.empty()) {
55  G4cout << "No Top LV Found" << G4endl;
56  } else {
57  for (unsigned int ii = 0; ii < topLV.size(); ++ii) {
58  G4cout << "Top LV Name " << topLV[ii]->GetName() << G4endl;
60  }
61  }
62 }
63 
64 void CheckOverlap::checkHierarchyLeafPVLV(G4LogicalVolume* lv, unsigned int leafDepth) {
65  //----- Get LV daughters from list of PV daughters
66  mmlvpv lvpvDaughters;
67  std::set<G4LogicalVolume*> lvDaughters;
68  int NoDaughters = lv->GetNoDaughters();
69  while ((NoDaughters--) > 0) {
70  G4VPhysicalVolume* pvD = lv->GetDaughter(NoDaughters);
71  lvpvDaughters.insert(mmlvpv::value_type(pvD->GetLogicalVolume(), pvD));
72  lvDaughters.insert(pvD->GetLogicalVolume());
73  }
74 
75  std::set<G4LogicalVolume*>::const_iterator scite;
76  mmlvpv::const_iterator mmcite;
77 
78  //----- Check daughters of LV
79  for (scite = lvDaughters.begin(); scite != lvDaughters.end(); scite++) {
80  std::pair<mmlvpv::iterator, mmlvpv::iterator> mmER = lvpvDaughters.equal_range(*scite);
81  //----- Check daughters PV of this LV
82  for (mmcite = mmER.first; mmcite != mmER.second; mmcite++)
83  checkPV((*mmcite).second, leafDepth + 1);
84  //----- Check daughters LV
85  checkHierarchyLeafPVLV(*scite, leafDepth + 1);
86  }
87 }
88 
89 void CheckOverlap::checkPV(G4VPhysicalVolume* pv, unsigned int leafDepth) {
90  //----- PV info
91 #ifndef G4V7
92  std::string mother = "DDDWorld";
93  if (pv->GetMotherLogical())
94  mother = pv->GetMotherLogical()->GetName();
95  if (!pv->IsReplicated()) {
96  G4PVPlacement* pvplace = dynamic_cast<G4PVPlacement*>(pv);
97  G4bool ok = pvplace->CheckOverlaps(nPoints);
98  G4cout << "Placed PV " << pvplace->GetName() << " Number " << pvplace->GetCopyNo() << " in mother " << mother
99  << " at depth " << leafDepth << " Status " << ok << G4endl;
100  if (ok) {
101  if (pv->GetRotation() == nullptr) {
102  G4cout << "Translation " << pv->GetTranslation() << " and with no rotation" << G4endl;
103  } else {
104  G4cout << "Translation " << pv->GetTranslation() << " and with rotation " << *(pv->GetRotation()) << G4endl;
105  }
106  G4LogicalVolume* lv = pv->GetLogicalVolume();
107  dumpLV(lv, "Self");
108  if (pv->GetMotherLogical()) {
109  lv = pv->GetMotherLogical();
110  dumpLV(lv, "Mother");
111  }
112  }
113  } else {
114  if (pv->GetParameterisation() != nullptr) {
115  G4PVParameterised* pvparam = dynamic_cast<G4PVParameterised*>(pv);
116  G4bool ok = pvparam->CheckOverlaps(nPoints);
117  G4cout << "Parametrized PV " << pvparam->GetName() << " in mother " << mother << " at depth " << leafDepth
118  << " Status " << ok << G4endl;
119  }
120  }
121 #endif
122 }
123 
124 G4VPhysicalVolume* CheckOverlap::getTopPV() {
125  return G4TransportationManager::GetTransportationManager()->GetNavigatorForTracking()->GetWorldVolume();
126 }
127 
128 void CheckOverlap::dumpLV(G4LogicalVolume* lv, std::string str) {
129  G4cout << "Dump of " << str << " Logical Volume " << lv->GetName() << " Solid: " << lv->GetSolid()->GetName()
130  << " Material: " << lv->GetMaterial()->GetName() << G4endl;
131  G4cout << *(lv->GetSolid()) << G4endl;
132 }
void checkPV(G4VPhysicalVolume *pv, unsigned int leafDepth)
Definition: CheckOverlap.cc:89
std::multimap< G4LogicalVolume *, G4VPhysicalVolume *, std::less< G4LogicalVolume * > > mmlvpv
void checkHierarchyLeafPVLV(G4LogicalVolume *lv, unsigned int leafDepth)
Definition: CheckOverlap.cc:64
~CheckOverlap() override
Definition: CheckOverlap.cc:30
void dumpLV(G4LogicalVolume *lv, std::string str)
def pv(vc)
Definition: MetAnalyzer.py:7
void update(const BeginOfRun *run) override
This routine will be called when the appropriate signal arrives.
Definition: CheckOverlap.cc:32
std::vector< std::string > nodeNames
Definition: CheckOverlap.h:31
ii
Definition: cuy.py:589
CheckOverlap(edm::ParameterSet const &p)
Definition: CheckOverlap.cc:20
G4VPhysicalVolume * getTopPV()
std::vector< G4LogicalVolume * > topLV
Definition: CheckOverlap.h:33
#define str(s)