CMS 3D CMS Logo

Public Member Functions | Protected Member Functions | Private Attributes

FiberSD Class Reference

#include <FiberSD.h>

Inheritance diagram for FiberSD:
SensitiveCaloDetector Observer< const BeginOfRun * > Observer< const BeginOfEvent * > Observer< const EndOfEvent * > SensitiveDetector

List of all members.

Public Member Functions

virtual void clear ()
virtual void DrawAll ()
virtual void EndOfEvent (G4HCofThisEvent *HCE)
 FiberSD (std::string, const DDCompactView &, SensitiveDetectorCatalog &, edm::ParameterSet const &, const SimTrackManager *)
virtual void Initialize (G4HCofThisEvent *HCE)
virtual void PrintAll ()
virtual G4bool ProcessHits (G4Step *aStep, G4TouchableHistory *ROhist)
virtual ~FiberSD ()

Protected Member Functions

virtual void clearHits ()
virtual void fillHits (edm::PCaloHitContainer &, std::string)
virtual uint32_t setDetUnitId (G4Step *)
virtual void update (const BeginOfRun *)
 This routine will be called when the appropriate signal arrives.
virtual void update (const ::EndOfEvent *)
virtual void update (const BeginOfEvent *)
 This routine will be called when the appropriate signal arrives.

Private Attributes

const SimTrackManagerm_trackManager
FiberG4HitsCollectiontheHC
G4int theHCID
std::string theName
HFShowertheShower

Detailed Description

Definition at line 25 of file FiberSD.h.


Constructor & Destructor Documentation

FiberSD::FiberSD ( std::string  name,
const DDCompactView cpv,
SensitiveDetectorCatalog clg,
edm::ParameterSet const &  p,
const SimTrackManager manager 
)

Definition at line 17 of file FiberSD.cc.

References SensitiveDetector::AssignSD(), ecaldqm::collectionName, LogDebug, SensitiveDetectorCatalog::logicalNames(), SensitiveDetector::Register(), and theShower.

                                                 :
  SensitiveCaloDetector(name, cpv, clg, p), theName(name),
  m_trackManager(manager), theHCID(-1), theHC(0) {

  collectionName.insert(name);
  LogDebug("FiberSim") << "***************************************************"
                       << "\n"
                       << "*                                                 *"
                       << "\n"
                       << "* Constructing a FiberSD  with name " << GetName()
                       << "\n"
                       << "*                                                 *"
                       << "\n"
                       << "***************************************************";
  theShower = new HFShower(name, cpv, p, 1);

  //
  // Now attach the right detectors (LogicalVolumes) to me
  //
  std::vector<std::string> lvNames = clg.logicalNames(name);
  this->Register();
  for (std::vector<std::string>::iterator it=lvNames.begin();
       it !=lvNames.end(); it++){
    this->AssignSD(*it);
    LogDebug("FiberSim") << "FiberSD : Assigns SD to LV " << (*it);
  }
}
FiberSD::~FiberSD ( ) [virtual]

Definition at line 47 of file FiberSD.cc.

References theHC, and theShower.

                  {
 
  if (theShower) delete theShower;
  if (theHC)    delete theHC;
}

Member Function Documentation

void FiberSD::clear ( void  ) [virtual]

Definition at line 112 of file FiberSD.cc.

Referenced by EndOfEvent().

{}
void FiberSD::clearHits ( ) [protected, virtual]

Implements SensitiveDetector.

Definition at line 124 of file FiberSD.cc.

{}
void FiberSD::DrawAll ( ) [virtual]

Definition at line 114 of file FiberSD.cc.

{}
void FiberSD::EndOfEvent ( G4HCofThisEvent *  HCE) [virtual]

Reimplemented from SensitiveDetector.

Definition at line 106 of file FiberSD.cc.

References clear(), LogDebug, and theHC.

                                              {
 
  LogDebug("FiberSim") << "FiberSD: Sees" << theHC->entries() << " hits";
  clear();
}
void FiberSD::fillHits ( edm::PCaloHitContainer ,
std::string   
) [protected, virtual]

Implements SensitiveCaloDetector.

Definition at line 134 of file FiberSD.cc.

{}
void FiberSD::Initialize ( G4HCofThisEvent *  HCE) [virtual]

Reimplemented from SensitiveDetector.

Definition at line 53 of file FiberSD.cc.

References ecaldqm::collectionName, LogDebug, theHC, and theHCID.

                                              {

  LogDebug("FiberSim") << "FiberSD : Initialize called for " << GetName();
  theHC = new FiberG4HitsCollection(GetName(), collectionName[0]);
  if (theHCID<0)
    theHCID = G4SDManager::GetSDMpointer()->GetCollectionID(collectionName[0]);
  HCE->AddHitsCollection(theHCID, theHC);
  
}
void FiberSD::PrintAll ( ) [virtual]

Definition at line 116 of file FiberSD.cc.

{}
G4bool FiberSD::ProcessHits ( G4Step *  aStep,
G4TouchableHistory *  ROhist 
) [virtual]

Implements SensitiveDetector.

Definition at line 63 of file FiberSD.cc.

References HFShower::getHits(), i, LogDebug, position, setDetUnitId(), theHC, and theShower.

                                                               {

  std::vector<HFShower::Hit> hits = theShower->getHits(aStep);

  if (hits.size() > 0) {
    std::vector<HFShowerPhoton> thePE;
    for (unsigned int i=0; i<hits.size(); i++) {
      HFShowerPhoton pe = HFShowerPhoton(hits[i].position.x(), 
                                         hits[i].position.y(), 
                                         hits[i].position.z(), 
                                         hits[i].wavelength, hits[i].time);
      thePE.push_back(pe);
    }
    int trackID = aStep->GetTrack()->GetTrackID();
    G4StepPoint* preStepPoint = aStep->GetPreStepPoint();
    const G4VTouchable* touch = preStepPoint->GetTouchable();
    G4LogicalVolume* lv = touch->GetVolume(0)->GetLogicalVolume();
    int depth = (touch->GetReplicaNumber(0))%10;
    int detID = setDetUnitId(aStep);
    math::XYZPoint theHitPos(preStepPoint->GetPosition().x(),
                             preStepPoint->GetPosition().y(),
                             preStepPoint->GetPosition().z());
    
    FiberG4Hit *aHit = new FiberG4Hit(lv, detID, depth, trackID);
    aHit->setNpe(hits.size());
    aHit->setPos(theHitPos);
    aHit->setTime(preStepPoint->GetGlobalTime());
    aHit->setPhoton(thePE);

    LogDebug("FiberSim") << "FiberSD: Hit created at " << lv->GetName()
                         << " DetID: " << aHit->towerId() << " Depth: " 
                         << aHit->depth() << " Track ID: " << aHit->trackId() 
                         << " Nb. of Cerenkov Photons: " << aHit->npe()
                         << " Time: " << aHit->time() << " at " 
                         << aHit->hitPos();
    for (unsigned int i=0; i<thePE.size(); i++) 
      LogDebug("FiberSim") << "FiberSD: PE[" << i << "] " << thePE[i];

    theHC->insert(aHit);
  }
  return true;
}
uint32_t FiberSD::setDetUnitId ( G4Step *  aStep) [protected, virtual]

Implements SensitiveDetector.

Definition at line 126 of file FiberSD.cc.

Referenced by ProcessHits().

                                            {
  const G4VTouchable* touch = aStep->GetPreStepPoint()->GetTouchable();
  int fibre = (touch->GetReplicaNumber(1))%10;
  int cell  = (touch->GetReplicaNumber(2));
  int tower = (touch->GetReplicaNumber(3));
  return ((tower*1000+cell)*10+fibre);
}
void FiberSD::update ( const BeginOfRun ) [protected, virtual]

This routine will be called when the appropriate signal arrives.

Implements Observer< const BeginOfRun * >.

Definition at line 118 of file FiberSD.cc.

{}
void FiberSD::update ( const ::EndOfEvent ) [protected, virtual]

Definition at line 122 of file FiberSD.cc.

{}
void FiberSD::update ( const BeginOfEvent ) [protected, virtual]

This routine will be called when the appropriate signal arrives.

Implements Observer< const BeginOfEvent * >.

Definition at line 120 of file FiberSD.cc.

{}

Member Data Documentation

Definition at line 56 of file FiberSD.h.

Definition at line 60 of file FiberSD.h.

Referenced by EndOfEvent(), Initialize(), ProcessHits(), and ~FiberSD().

G4int FiberSD::theHCID [private]

Definition at line 59 of file FiberSD.h.

Referenced by Initialize().

std::string FiberSD::theName [private]

Definition at line 55 of file FiberSD.h.

Definition at line 57 of file FiberSD.h.

Referenced by FiberSD(), ProcessHits(), and ~FiberSD().