CMS 3D CMS Logo

Public Member Functions | Private Member Functions | Private Attributes

PrintSensitive Class Reference

#include <PrintSensitive.h>

Inheritance diagram for PrintSensitive:
SimWatcher Observer< const BeginOfRun * >

List of all members.

Public Member Functions

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

Private Member Functions

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

Private Attributes

G4NavigationHistory fHistory
std::string name
int nchar

Detailed Description

Definition at line 16 of file PrintSensitive.h.


Constructor & Destructor Documentation

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

Definition at line 14 of file PrintSensitive.cc.

References gather_cfg::cout, edm::ParameterSet::getUntrackedParameter(), name, and nchar.

                                                       {
  name  = p.getUntrackedParameter<std::string>("Name","*");
  nchar = name.find("*");
  name.assign(name,0,nchar);
  std::cout << "PrintSensitive:: Print position of all Sensitive Touchables: "
            << " for names (0-" << nchar << ") = " << name << "\n";
}
PrintSensitive::~PrintSensitive ( )

Definition at line 22 of file PrintSensitive.cc.

{}

Member Function Documentation

void PrintSensitive::dumpTouch ( G4VPhysicalVolume *  pv,
unsigned int  leafDepth,
bool  printIt,
std::ostream &  out = std::cout 
) [private]

Definition at line 30 of file PrintSensitive.cc.

References fHistory, name, and nchar.

Referenced by update().

                                                               {

  if (leafDepth == 0) fHistory.SetFirstEntry(pv);
  else fHistory.NewLevel(pv, kNormal, pv->GetCopyNo());

  G4ThreeVector globalpoint = fHistory.GetTopTransform().Inverse().
    TransformPoint(G4ThreeVector(0,0,0));
  G4LogicalVolume * lv = pv->GetLogicalVolume();

  std::string mother = "World";
  if (pv->GetMotherLogical()) mother = pv->GetMotherLogical()->GetName();
  std::string lvname = lv->GetName();
  lvname.assign(lvname,0,nchar);
  if (lvname == name) printIt = true;

  if (lv->GetSensitiveDetector() && printIt) {
    out << leafDepth << " ### VOLUME = " << lv->GetName() 
        << " Copy No " << pv->GetCopyNo() << " in " << mother
        << " global position of centre " << globalpoint << " (r=" 
        <<  globalpoint.perp() << ", phi=" <<  globalpoint.phi()/deg
        << ")\n";
  }

  int NoDaughters = lv->GetNoDaughters();
  while ((NoDaughters--)>0)  {
    G4VPhysicalVolume * pvD = lv->GetDaughter(NoDaughters);
    if (!pvD->IsReplicated()) dumpTouch(pvD, leafDepth+1, printIt, out);
  }

  if (leafDepth > 0) fHistory.BackLevel();
}
G4VPhysicalVolume * PrintSensitive::getTopPV ( ) [private]

Definition at line 63 of file PrintSensitive.cc.

Referenced by update().

                                             {
  return G4TransportationManager::GetTransportationManager()
    ->GetNavigatorForTracking()->GetWorldVolume();
}
void PrintSensitive::update ( const BeginOfRun ) [private, virtual]

This routine will be called when the appropriate signal arrives.

Implements Observer< const BeginOfRun * >.

Definition at line 24 of file PrintSensitive.cc.

References gather_cfg::cout, dumpTouch(), and getTopPV().

                                                  {

  G4VPhysicalVolume * theTopPV = getTopPV();
  dumpTouch(theTopPV, 0, false, std::cout);
}

Member Data Documentation

G4NavigationHistory PrintSensitive::fHistory [private]

Definition at line 34 of file PrintSensitive.h.

Referenced by dumpTouch().

std::string PrintSensitive::name [private]

Definition at line 32 of file PrintSensitive.h.

Referenced by dumpTouch(), and PrintSensitive().

int PrintSensitive::nchar [private]

Definition at line 33 of file PrintSensitive.h.

Referenced by dumpTouch(), and PrintSensitive().