CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_1_8_patch12/src/SimG4CMS/ShowerLibraryProducer/src/HFChamberSD.cc

Go to the documentation of this file.
00001 #include "SimG4CMS/ShowerLibraryProducer/interface/HFChamberSD.h"
00002 #include "DataFormats/Math/interface/Point3D.h"
00003 
00004 #include "G4VPhysicalVolume.hh"
00005 #include "G4PVPlacement.hh"
00006 #include "G4HCofThisEvent.hh"
00007 #include "G4TouchableHistory.hh"
00008 #include "G4Track.hh"
00009 #include "G4Step.hh"
00010 #include "G4VSolid.hh"
00011 #include "G4DynamicParticle.hh"
00012 #include "G4ParticleDefinition.hh"
00013 #include "G4SDManager.hh"
00014 #include "G4ios.hh"
00015 
00016 HFChamberSD::HFChamberSD(std::string name, const DDCompactView & cpv,
00017                  SensitiveDetectorCatalog & clg, edm::ParameterSet const & p,
00018                  const SimTrackManager* manager) :
00019   SensitiveCaloDetector(name, cpv, clg, p), theName(name),
00020   m_trackManager(manager), theHCID(-1), theHC(0), theNSteps(0) {
00021 
00022   collectionName.insert(name);
00023   LogDebug("FiberSim") << "***************************************************"
00024                        << "\n"
00025                        << "*                                                 *"
00026                        << "\n"
00027                        << "* Constructing a HFChamberSD  with name " << GetName()
00028                        << "\n"
00029                        << "*                                                 *"
00030                        << "\n"
00031                        << "***************************************************";
00032   //
00033   // Now attach the right detectors (LogicalVolumes) to me
00034   //
00035   std::vector<std::string> lvNames = clg.logicalNames(name);
00036   this->Register();
00037   for (std::vector<std::string>::iterator it=lvNames.begin();
00038        it !=lvNames.end(); it++){
00039     this->AssignSD(*it);
00040     LogDebug("FiberSim") << "HFChamberSD : Assigns SD to LV " << (*it);
00041   }
00042 }
00043 
00044 HFChamberSD::~HFChamberSD() {
00045   if (theHC)    delete theHC;
00046 }
00047 
00048 void HFChamberSD::Initialize(G4HCofThisEvent * HCE) {
00049 
00050   LogDebug("FiberSim") << "HFChamberSD : Initialize called for " << GetName();
00051   theHC = new HFShowerG4HitsCollection(GetName(), collectionName[0]);
00052   if (theHCID<0)
00053     theHCID = G4SDManager::GetSDMpointer()->GetCollectionID(collectionName[0]);
00054   HCE->AddHitsCollection(theHCID, theHC);
00055 
00056 }
00057 
00058 G4bool HFChamberSD::ProcessHits(G4Step * aStep, G4TouchableHistory*) {
00059 
00060   //do not process hits other than primary particle hits:
00061   double charge = aStep->GetTrack()->GetDefinition()->GetPDGCharge();
00062   int trackID = aStep->GetTrack()->GetTrackID();
00063   if(charge == 0. || trackID != 1 ||aStep->GetTrack()->GetParentID() != 0 || aStep->GetTrack()->GetCreatorProcess() != NULL) return false;
00064   ++theNSteps;
00065   //if(theNSteps>1)return false;
00066 
00067   G4StepPoint* preStepPoint = aStep->GetPreStepPoint();
00068   const G4VTouchable* touch = preStepPoint->GetTouchable();
00069   int detID   = setDetUnitId(aStep);
00070 
00071   double edep = aStep->GetTotalEnergyDeposit();
00072   double time = (preStepPoint->GetGlobalTime())/ns;
00073 
00074   G4ThreeVector globalPos = preStepPoint->GetPosition();
00075   G4ThreeVector localPos  = touch->GetHistory()->GetTopTransform().TransformPoint(globalPos);
00076   const G4DynamicParticle* particle =  aStep->GetTrack()->GetDynamicParticle();
00077   G4ThreeVector momDir   = particle->GetMomentumDirection();
00078 
00079   HFShowerG4Hit *aHit = new HFShowerG4Hit(detID, trackID, edep, time);
00080   aHit->setLocalPos(localPos);
00081   aHit->setGlobalPos(globalPos);
00082   aHit->setPrimMomDir(momDir);
00083 
00084   LogDebug("FiberSim") << "HFChamberSD: Hit created in ("
00085                        << touch->GetVolume(0)->GetLogicalVolume()->GetName()
00086                        << ") " << " ID " << detID << " Track " << trackID
00087                        << " Edep: " << edep/MeV << " MeV; Time: " << time
00088                        << " ns; Position (local) " << localPos << " (global ) "
00089                        << globalPos << " direction " << momDir;
00090 
00091   theHC->insert(aHit);
00092   return true;
00093 }
00094 
00095 void HFChamberSD::EndOfEvent(G4HCofThisEvent * HCE) {
00096 
00097   LogDebug("FiberSim") << "HFChamberSD: Sees" << theHC->entries() << " hits";
00098   clear();
00099 }
00100 
00101 void HFChamberSD::clear() {
00102   theNSteps = 0;
00103 }
00104 
00105 void HFChamberSD::DrawAll()  {}
00106 
00107 void HFChamberSD::PrintAll() {}
00108 
00109 void HFChamberSD::clearHits() {}
00110 
00111 uint32_t HFChamberSD::setDetUnitId(G4Step* aStep) {
00112   const G4VTouchable* touch = aStep->GetPreStepPoint()->GetTouchable();
00113   return (touch->GetReplicaNumber(0));
00114 }
00115 
00116 void HFChamberSD::fillHits(edm::PCaloHitContainer&, std::string) {}