CMS 3D CMS Logo

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