CMS 3D CMS Logo

List of all members | Public Member Functions | Private Member Functions | Private Attributes
PrintSensitive Class Reference
Inheritance diagram for PrintSensitive:
SimWatcher Observer< const BeginOfRun *>

Public Member Functions

 PrintSensitive (edm::ParameterSet const &p)
 
 ~PrintSensitive () override
 
- Public Member Functions inherited from SimWatcher
virtual void beginRun (edm::EventSetup const &)
 
bool isMT () const
 
const SimWatcheroperator= (const SimWatcher &)=delete
 
virtual void registerConsumes (edm::ConsumesCollector)
 
 SimWatcher ()
 
 SimWatcher (const SimWatcher &)=delete
 
virtual ~SimWatcher ()
 
- Public Member Functions inherited from Observer< const BeginOfRun *>
 Observer ()
 
void slotForUpdate (const BeginOfRun * iT)
 
virtual ~Observer ()
 

Private Member Functions

int dumpTouch (G4VPhysicalVolume *pv, unsigned int leafDepth, bool printIt, int ns, std::ostream &out=G4cout)
 
G4VPhysicalVolume * getTopPV ()
 
void update (const BeginOfRun *run) override
 This routine will be called when the appropriate signal arrives. More...
 

Private Attributes

G4NavigationHistory fHistory
 
std::string name_
 
int nchar_
 

Additional Inherited Members

- Protected Member Functions inherited from SimWatcher
void setMT (bool val)
 

Detailed Description

Definition at line 17 of file PrintSensitive.cc.

Constructor & Destructor Documentation

◆ PrintSensitive()

PrintSensitive::PrintSensitive ( edm::ParameterSet const &  p)

Definition at line 33 of file PrintSensitive.cc.

References ecalTB2006H4_GenSimDigiReco_cfg::G4cout, name_, nchar_, AlCaHLTBitMon_ParallelJobs::p, and AlCaHLTBitMon_QueryRunRegistry::string.

33  {
34  name_ = p.getUntrackedParameter<std::string>("Name", "*");
35  nchar_ = name_.find('*');
36  name_.assign(name_, 0, nchar_);
37  G4cout << "PrintSensitive:: Print position of all Sensitive Touchables: "
38  << " for names (0-" << nchar_ << ") = " << name_ << G4endl;
39 }
std::string name_

◆ ~PrintSensitive()

PrintSensitive::~PrintSensitive ( )
override

Definition at line 41 of file PrintSensitive.cc.

41 {}

Member Function Documentation

◆ dumpTouch()

int PrintSensitive::dumpTouch ( G4VPhysicalVolume *  pv,
unsigned int  leafDepth,
bool  printIt,
int  ns,
std::ostream &  out = G4cout 
)
private

Definition at line 49 of file PrintSensitive.cc.

References fHistory, name_, nchar_, MillePedeFileConverter_cfg::out, riemannFit::printIt(), MetAnalyzer::pv(), and AlCaHLTBitMon_QueryRunRegistry::string.

Referenced by update().

49  {
50  if (leafDepth == 0)
51  fHistory.SetFirstEntry(pv);
52  else
53  fHistory.NewLevel(pv, kNormal, pv->GetCopyNo());
54 
55  int nsens(ns);
56  G4ThreeVector globalpoint = fHistory.GetTopTransform().Inverse().TransformPoint(G4ThreeVector(0, 0, 0));
57  G4LogicalVolume *lv = pv->GetLogicalVolume();
58 
59  std::string mother = "World";
60  if (pv->GetMotherLogical())
61  mother = pv->GetMotherLogical()->GetName();
62  std::string lvname = lv->GetName();
63  lvname.assign(lvname, 0, nchar_);
64  if (lvname == name_)
65  printIt = true;
66 
67  if (lv->GetSensitiveDetector() && printIt) {
68  ++nsens;
69  out << nsens << ":" << leafDepth << " ### VOLUME = " << lv->GetName() << " Copy No " << pv->GetCopyNo() << " in "
70  << mother << " global position of centre " << globalpoint << " (r=" << globalpoint.perp()
71  << ", phi=" << globalpoint.phi() / CLHEP::deg << ")\n";
72  }
73 
74  int NoDaughters = lv->GetNoDaughters();
75  while ((NoDaughters--) > 0) {
76  G4VPhysicalVolume *pvD = lv->GetDaughter(NoDaughters);
77  if (!pvD->IsReplicated())
78  nsens = dumpTouch(pvD, leafDepth + 1, printIt, nsens, out);
79  }
80 
81  if (leafDepth > 0)
82  fHistory.BackLevel();
83  return nsens;
84 }
G4NavigationHistory fHistory
def pv(vc)
Definition: MetAnalyzer.py:7
std::string name_
int dumpTouch(G4VPhysicalVolume *pv, unsigned int leafDepth, bool printIt, int ns, std::ostream &out=G4cout)
__host__ __device__ void printIt(C *m, const char *prefix="")
Definition: FitUtils.h:53

◆ getTopPV()

G4VPhysicalVolume * PrintSensitive::getTopPV ( )
private

Definition at line 86 of file PrintSensitive.cc.

Referenced by update().

86  {
87  return G4TransportationManager::GetTransportationManager()->GetNavigatorForTracking()->GetWorldVolume();
88 }

◆ update()

void PrintSensitive::update ( const BeginOfRun )
overrideprivatevirtual

This routine will be called when the appropriate signal arrives.

Implements Observer< const BeginOfRun *>.

Definition at line 43 of file PrintSensitive.cc.

References dumpTouch(), ecalTB2006H4_GenSimDigiReco_cfg::G4cout, getTopPV(), and name_.

Referenced by progressbar.ProgressBar::__next__(), MatrixUtil.Matrix::__setitem__(), MatrixUtil.Steps::__setitem__(), progressbar.ProgressBar::finish(), and MatrixUtil.Steps::overwrite().

43  {
44  G4VPhysicalVolume *theTopPV = getTopPV();
45  int nsens = dumpTouch(theTopPV, 0, false, 0, G4cout);
46  G4cout << "\nTotal number of sensitive detector volumes for " << name_ << " is " << nsens << G4endl;
47 }
G4VPhysicalVolume * getTopPV()
std::string name_
int dumpTouch(G4VPhysicalVolume *pv, unsigned int leafDepth, bool printIt, int ns, std::ostream &out=G4cout)

Member Data Documentation

◆ fHistory

G4NavigationHistory PrintSensitive::fHistory
private

Definition at line 30 of file PrintSensitive.cc.

Referenced by dumpTouch().

◆ name_

std::string PrintSensitive::name_
private

Definition at line 28 of file PrintSensitive.cc.

Referenced by dumpTouch(), PrintSensitive(), and update().

◆ nchar_

int PrintSensitive::nchar_
private

Definition at line 29 of file PrintSensitive.cc.

Referenced by dumpTouch(), and PrintSensitive().