#include <SimG4CMS/Forward/interface/TotemSD.h>
Public Member Functions | |
virtual void | clear () |
virtual void | DrawAll () |
virtual void | EndOfEvent (G4HCofThisEvent *eventHC) |
void | fillHits (edm::PSimHitContainer &, std::string use) |
virtual void | Initialize (G4HCofThisEvent *HCE) |
virtual void | PrintAll () |
virtual bool | ProcessHits (G4Step *, G4TouchableHistory *) |
virtual uint32_t | setDetUnitId (G4Step *) |
TotemSD (std::string, const DDCompactView &, SensitiveDetectorCatalog &, edm::ParameterSet const &, const SimTrackManager *) | |
virtual | ~TotemSD () |
Private Member Functions | |
virtual void | clearHits () |
void | CreateNewHit () |
void | CreateNewHitEvo () |
void | GetStepInfo (G4Step *aStep) |
bool | HitExists () |
G4ThreeVector | PosizioEvo (G4ThreeVector, double, double, double, double, int &) |
void | ResetForNewPrimary () |
G4ThreeVector | SetToLocal (G4ThreeVector globalPoint) |
void | StoreHit (TotemG4Hit *) |
void | Summarize () |
void | update (const ::EndOfEvent *) |
void | update (const BeginOfEvent *) |
This routine will be called when the appropriate signal arrives. | |
void | UpdateHit () |
Private Attributes | |
TotemG4Hit * | currentHit |
G4VPhysicalVolume * | currentPV |
float | edeposit |
float | Eloss |
G4ThreeVector | entrancePoint |
int | eventno |
G4int | hcID |
G4ThreeVector | hitPoint |
float | incidentEnergy |
std::string | name |
TotemVDetectorOrganization * | numberingScheme |
float | Pabs |
int | ParentId |
short | ParticleType |
float | PhiAtEntry |
G4ThreeVector | Posizio |
G4StepPoint * | postStepPoint |
G4StepPoint * | preStepPoint |
uint32_t | previousUnitID |
int | primaryID |
G4int | primID |
TrackingSlaveSD * | slave |
TotemG4HitCollection * | theHC |
const SimTrackManager * | theManager |
float | ThetaAtEntry |
G4Track * | theTrack |
float | Tof |
int | tsID |
double | tSlice |
int | tSliceID |
uint32_t | unitID |
float | Vx |
float | Vy |
float | Vz |
Description: Stores hits of Totem in appropriate container
Usage: Used in sensitive detector builder
TotemSD::TotemSD | ( | std::string | name, |
const DDCompactView & | cpv, | ||
SensitiveDetectorCatalog & | clg, | ||
edm::ParameterSet const & | p, | ||
const SimTrackManager * | manager | ||
) |
Definition at line 45 of file TotemSD.cc.
References SensitiveDetector::AssignSD(), ecaldqm::collectionName, edm::ParameterSet::getParameter(), edm::ParameterSet::getUntrackedParameter(), LogDebug, SensitiveDetectorCatalog::logicalNames(), numberingScheme, SensitiveDetector::Register(), and slave.
: SensitiveTkDetector(name, cpv, clg, p), numberingScheme(0), name(name), hcID(-1), theHC(0), theManager(manager), currentHit(0), theTrack(0), currentPV(0), unitID(0), previousUnitID(0), preStepPoint(0), postStepPoint(0), eventno(0){ //Add Totem Sentitive Detector Names collectionName.insert(name); //Parameters edm::ParameterSet m_p = p.getParameter<edm::ParameterSet>("TotemSD"); int verbn = m_p.getUntrackedParameter<int>("Verbosity"); SetVerboseLevel(verbn); LogDebug("ForwardSim") << "*******************************************************\n" << "* *\n" << "* Constructing a TotemSD with name " << name << "\n" << "* *\n" << "*******************************************************"; slave = new TrackingSlaveSD(name); // // 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); edm::LogInfo("ForwardSim") << "TotemSD : Assigns SD to LV " << (*it); } if (name == "TotemHitsT1") { numberingScheme = dynamic_cast<TotemVDetectorOrganization*>(new TotemT1NumberingScheme(1)); } else if (name == "TotemHitsT2Si") { numberingScheme = dynamic_cast<TotemVDetectorOrganization*>(new TotemT2NumberingSchemeGem(3)); } else if (name == "TotemHitsT2Gem") { numberingScheme = dynamic_cast<TotemVDetectorOrganization*>(new TotemT2NumberingSchemeGem(4)); } else if (name == "TotemHitsRP") { numberingScheme = dynamic_cast<TotemVDetectorOrganization*>(new TotemRPNumberingScheme(3)); } else { edm::LogWarning("ForwardSim") << "TotemSD: ReadoutName not supported\n"; } edm::LogInfo("ForwardSim") << "TotemSD: Instantiation completed"; }
TotemSD::~TotemSD | ( | ) | [virtual] |
Definition at line 96 of file TotemSD.cc.
References numberingScheme, and slave.
{ if (slave) delete slave; if (numberingScheme) delete numberingScheme; }
void TotemSD::clear | ( | void | ) | [virtual] |
Definition at line 165 of file TotemSD.cc.
{ }
void TotemSD::clearHits | ( | ) | [private, virtual] |
Implements SensitiveDetector.
Definition at line 190 of file TotemSD.cc.
References TrackingSlaveSD::Initialize(), and slave.
Referenced by update().
{ slave->Initialize(); }
void TotemSD::CreateNewHit | ( | ) | [private] |
Definition at line 285 of file TotemSD.cc.
References gather_cfg::cout, currentHit, currentPV, Eloss, incidentEnergy, LogDebug, NULL, Pabs, ParentId, ParticleType, PhiAtEntry, Posizio, primaryID, TotemG4Hit::setEnergyLoss(), TotemG4Hit::setEntry(), TotemG4Hit::setIncidentEnergy(), TotemG4Hit::setPabs(), TotemG4Hit::setParentId(), TotemG4Hit::setParticleType(), TotemG4Hit::setPhiAtEntry(), TotemG4Hit::setThetaAtEntry(), TotemG4Hit::setTimeSlice(), TotemG4Hit::setTof(), TotemG4Hit::setTrackID(), TotemG4Hit::setUnitID(), TotemG4Hit::setVx(), TotemG4Hit::setVy(), TotemG4Hit::setVz(), StoreHit(), ThetaAtEntry, theTrack, Tof, tSlice, tSliceID, unitID, UpdateHit(), Vx, Vy, and Vz.
Referenced by ProcessHits().
{ #ifdef debug LogDebug("ForwardSim") << "TotemSD CreateNewHit for" << " PV " << currentPV->GetName() << " PVid = " << currentPV->GetCopyNo() << " MVid = " << currentPV->GetMother()->GetCopyNo() << " Unit " << unitID << "\n" << " primary " << primaryID << " time slice " << tSliceID << " For Track " << theTrack->GetTrackID() << " which is a " << theTrack->GetDefinition()->GetParticleName(); if (theTrack->GetTrackID()==1) { LogDebug("ForwardSim") << " of energy " << theTrack->GetTotalEnergy(); } else { LogDebug("ForwardSim") << " daughter of part. " << theTrack->GetParentID(); } cout << " and created by " ; if (theTrack->GetCreatorProcess()!=NULL) LogDebug("ForwardSim") << theTrack->GetCreatorProcess()->GetProcessName() ; else LogDebug("ForwardSim") << "NO process"; #endif currentHit = new TotemG4Hit; currentHit->setTrackID(primaryID); currentHit->setTimeSlice(tSlice); currentHit->setUnitID(unitID); currentHit->setIncidentEnergy(incidentEnergy); currentHit->setPabs(Pabs); currentHit->setTof(Tof); currentHit->setEnergyLoss(Eloss); currentHit->setParticleType(ParticleType); currentHit->setThetaAtEntry(ThetaAtEntry); currentHit->setPhiAtEntry(PhiAtEntry); currentHit->setEntry(Posizio.x(),Posizio.y(),Posizio.z()); currentHit->setParentId(ParentId); currentHit->setVx(Vx); currentHit->setVy(Vy); currentHit->setVz(Vz); UpdateHit(); StoreHit(currentHit); }
void TotemSD::CreateNewHitEvo | ( | ) | [private] |
Definition at line 338 of file TotemSD.cc.
References currentHit, Eloss, incidentEnergy, Pabs, ParentId, ParticleType, PhiAtEntry, Posizio, PosizioEvo(), primaryID, TotemG4Hit::setEnergyLoss(), TotemG4Hit::setEntry(), TotemG4Hit::setIncidentEnergy(), TotemG4Hit::setPabs(), TotemG4Hit::setParentId(), TotemG4Hit::setParticleType(), TotemG4Hit::setPhiAtEntry(), TotemG4Hit::setThetaAtEntry(), TotemG4Hit::setTimeSlice(), TotemG4Hit::setTof(), TotemG4Hit::setTrackID(), TotemG4Hit::setUnitID(), TotemG4Hit::setVx(), TotemG4Hit::setVy(), TotemG4Hit::setVz(), StoreHit(), ThetaAtEntry, Tof, tSlice, unitID, UpdateHit(), Vx, Vy, and Vz.
Referenced by ProcessHits().
{ // LogDebug("ForwardSim") << "INSIDE CREATE NEW HIT EVO "; currentHit = new TotemG4Hit; currentHit->setTrackID(primaryID); currentHit->setTimeSlice(tSlice); currentHit->setUnitID(unitID); currentHit->setIncidentEnergy(incidentEnergy); currentHit->setPabs(Pabs); currentHit->setTof(Tof); currentHit->setEnergyLoss(Eloss); currentHit->setParticleType(ParticleType); currentHit->setThetaAtEntry(ThetaAtEntry); currentHit->setPhiAtEntry(PhiAtEntry); // LogDebug("ForwardSim") << Posizio.x() << " " << Posizio.y() << " " << Posizio.z(); currentHit->setParentId(ParentId); currentHit->setVx(Vx); currentHit->setVy(Vy); currentHit->setVz(Vz); G4ThreeVector _PosizioEvo; int flagAcc=0; _PosizioEvo=PosizioEvo(Posizio,Vx,Vy,Vz,Pabs,flagAcc); if(flagAcc==1){ currentHit->setEntry(_PosizioEvo.x(),_PosizioEvo.y(),_PosizioEvo.z()); // if(flagAcc==1) UpdateHit(); StoreHit(currentHit); } // LogDebug("ForwardSim") << "STORED HIT IN: " << unitID; }
void TotemSD::DrawAll | ( | ) | [virtual] |
Definition at line 168 of file TotemSD.cc.
{ }
void TotemSD::EndOfEvent | ( | G4HCofThisEvent * | eventHC | ) | [virtual] |
Reimplemented from SensitiveDetector.
Definition at line 138 of file TotemSD.cc.
References TotemG4Hit::getEnergyLoss(), TotemG4Hit::getEntry(), TotemG4Hit::getPabs(), TotemG4Hit::getParticleType(), TotemG4Hit::getPhiAtEntry(), TotemG4Hit::getThetaAtEntry(), TotemG4Hit::getTof(), TotemG4Hit::getTrackID(), TotemG4Hit::getUnitID(), j, LogDebug, TrackingSlaveSD::processHits(), slave, Summarize(), and theHC.
{ // here we loop over transient hits and make them persistent for (int j=0; j<theHC->entries() && j<15000; j++) { TotemG4Hit* aHit = (*theHC)[j]; #ifdef ddebug LogDebug("ForwardSim") << "HIT NUMERO " << j << "unit ID = " << aHit->getUnitID() << "\n" << " " << "enrty z " << aHit->getEntry().z() << "\n" << " " << "theta " << aHit->getThetaAtEntry() << "\n"; #endif Local3DPoint theExitPoint(0,0,0); Local3DPoint Entrata(aHit->getEntry().x(), aHit->getEntry().y(), aHit->getEntry().z()); slave->processHits(PSimHit(Entrata,theExitPoint, aHit->getPabs(), aHit->getTof(), aHit->getEnergyLoss(), aHit->getParticleType(), aHit->getUnitID(), aHit->getTrackID(), aHit->getThetaAtEntry(),aHit->getPhiAtEntry())); } Summarize(); }
void TotemSD::fillHits | ( | edm::PSimHitContainer & | c, |
std::string | use | ||
) | [virtual] |
Implements SensitiveTkDetector.
Definition at line 176 of file TotemSD.cc.
References TrackingSlaveSD::hits(), n, TrackingSlaveSD::name(), and slave.
void TotemSD::GetStepInfo | ( | G4Step * | aStep | ) | [private] |
Definition at line 202 of file TotemSD.cc.
References currentPV, edeposit, Eloss, hitPoint, LogDebug, name, Pabs, ParentId, ParticleType, PhiAtEntry, Posizio, postStepPoint, preStepPoint, primaryID, setDetUnitId(), ThetaAtEntry, theTrack, Tof, tSlice, tSliceID, unitID, Vx, Vy, and Vz.
Referenced by ProcessHits().
{ preStepPoint = aStep->GetPreStepPoint(); postStepPoint= aStep->GetPostStepPoint(); theTrack = aStep->GetTrack(); //Local3DPoint theEntryPoint = SensitiveDetector::InitialStepPosition(aStep,LocalCoordinates); //Local3DPoint theExitPoint = SensitiveDetector::FinalStepPosition(aStep,LocalCoordinates); hitPoint = preStepPoint->GetPosition(); currentPV = preStepPoint->GetPhysicalVolume(); // double weight = 1; G4String name = currentPV->GetName(); name.assign(name,0,4); G4String particleType = theTrack->GetDefinition()->GetParticleName(); edeposit = aStep->GetTotalEnergyDeposit(); tSlice = (postStepPoint->GetGlobalTime() )/nanosecond; tSliceID = (int) tSlice; unitID = setDetUnitId(aStep); #ifdef debug LogDebug("ForwardSim") << "UNITa " << unitID; #endif primaryID = theTrack->GetTrackID(); Posizio = hitPoint; Pabs = aStep->GetPreStepPoint()->GetMomentum().mag()/GeV; Tof = aStep->GetPostStepPoint()->GetGlobalTime()/nanosecond; Eloss = aStep->GetTotalEnergyDeposit()/GeV; ParticleType = theTrack->GetDefinition()->GetPDGEncoding(); ThetaAtEntry = aStep->GetPreStepPoint()->GetPosition().theta()/deg; PhiAtEntry = aStep->GetPreStepPoint()->GetPosition().phi()/deg; ParentId = theTrack->GetParentID(); Vx = theTrack->GetVertexPosition().x(); Vy = theTrack->GetVertexPosition().y(); Vz = theTrack->GetVertexPosition().z(); }
bool TotemSD::HitExists | ( | ) | [private] |
Definition at line 243 of file TotemSD.cc.
References currentHit, newFWLiteAna::found, TotemG4Hit::getTimeSliceID(), TotemG4Hit::getTrackID(), TotemG4Hit::getUnitID(), j, previousUnitID, primaryID, primID, ResetForNewPrimary(), theHC, tsID, tSliceID, unitID, and UpdateHit().
Referenced by ProcessHits().
{ if (primaryID<1) { edm::LogWarning("ForwardSim") << "***** TotemSD error: primaryID = " << primaryID << " maybe detector name changed"; } // Update if in the same detector, time-slice and for same track // if (primaryID == primID && tSliceID == tsID && unitID==previousUnitID) { if (tSliceID == tsID && unitID==previousUnitID) { UpdateHit(); return true; } // Reset entry point for new primary if (primaryID != primID) ResetForNewPrimary(); //look in the HitContainer whether a hit with the same primID, unitID, //tSliceID already exists: bool found = false; for (int j=0; j<theHC->entries()&&!found; j++) { TotemG4Hit* aPreviousHit = (*theHC)[j]; if (aPreviousHit->getTrackID() == primaryID && aPreviousHit->getTimeSliceID() == tSliceID && aPreviousHit->getUnitID() == unitID ) { currentHit = aPreviousHit; found = true; } } if (found) { UpdateHit(); return true; } else { return false; } }
void TotemSD::Initialize | ( | G4HCofThisEvent * | HCE | ) | [virtual] |
Reimplemented from SensitiveDetector.
Definition at line 125 of file TotemSD.cc.
References ecaldqm::collectionName, hcID, LogDebug, name, primID, theHC, and tsID.
{ LogDebug("ForwardSim") << "TotemSD : Initialize called for " << name; theHC = new TotemG4HitCollection(name, collectionName[0]); if (hcID<0) hcID = G4SDManager::GetSDMpointer()->GetCollectionID(collectionName[0]); HCE->AddHitsCollection(hcID, theHC); tsID = -2; primID = -2; }
G4ThreeVector TotemSD::PosizioEvo | ( | G4ThreeVector | Pos, |
double | vx, | ||
double | vy, | ||
double | vz, | ||
double | pabs, | ||
int & | accettanza | ||
) | [private] |
Definition at line 377 of file TotemSD.cc.
References j, and mathSSE::sqrt().
Referenced by CreateNewHitEvo().
{ accettanza=0; //Pos.xyz() in mm G4ThreeVector PosEvo; double ThetaX=atan((Pos.x()-vx)/(Pos.z()-vz)); double ThetaY=atan((Pos.y()-vy)/(Pos.z()-vz)); double X_at_0 =(vx-((Pos.x()-vx)/(Pos.z()-vz))*vz)/1000.; double Y_at_0 =(vy-((Pos.y()-vy)/(Pos.z()-vz))*vz)/1000.; // double temp_evo_X; // double temp_evo_Y; // double temp_evo_Z; // temp_evo_Z = fabs(Pos.z())/Pos.z()*220000.; //csi=-dp/d double csi = fabs((7000.-pabs)/7000.); // all in m const int no_rp=4; double x_par[no_rp+1]; double y_par[no_rp+1]; //rp z position double rp[no_rp]={141.,149.,198.,220.}; //{lx0,mlx} for each rp; Lx=lx0+mlx*csi double leffx[][2]={{122.5429,-46.9312},{125.4194,-49.1849},{152.6,-81.157},{98.8914,-131.8390}}; //{ly0,mly} for each rp; Ly=ly0+mly*csi double leffy[][2]={{124.2314,-55.4852},{127.7825,-57.4503},{179.455,-76.274},{273.0931,-40.4626}}; //{vx0,mvx0} for each rp; vx=vx0+mvx*csi double avx[][2]={{0.515483,-1.0123},{0.494122,-1.0534},{0.2217,-1.483},{0.004633,-1.0719}}; //{vy0,mvy0} for each rp; vy=vy0+mvy*csi double avy[][2]={{0.371418,-1.6327},{0.349035,-1.6955},{0.0815,-2.59},{0.007592,-4.0841}}; //{D0,md,a,b} for each rp; D=D0+(md+a*thetax)*csi+b*thetax double ddx[][4]= {{-0.082336,-0.092513,112.3436,-82.5029},{-0.086927,-0.097670,114.9513,-82.9835}, {-0.092117,-0.0915,180.6236,-82.443},{-0.050470,0.058837,208.1106,20.8198}}; // {10sigma_x+0.5mm,10sigma_y+0.5mm} double detlim[][2]={{0,0},{0.0028,0.0021},{0,0},{0.0008,0.0013}}; //{rmax,dmax} double pipelim[][2]={{0.026,0.026},{0.04,0.04},{0.0226,0.0177},{0.04,0.04}}; for(int j=0; j<no_rp ; j++) { //y=Ly*thetay+vy*y0 //x=Lx*thetax+vx*x0-csi*D y_par[j]=ThetaY*(leffy[j][0]+leffy[j][1]*csi)+(avy[j][0]+avy[j][1]*csi)*Y_at_0; x_par[j]=ThetaX*(leffx[j][0]+leffx[j][1]*csi)+(avx[j][0]+avx[j][1]*csi)*X_at_0- csi*(ddx[j][0]+(ddx[j][1]+ddx[j][2]*ThetaX)*csi+ddx[j][3]*ThetaX); } //pass TAN@141 if (fabs(y_par[0])<pipelim[0][1] && sqrt((y_par[0]*y_par[0])+(x_par[0]*x_par[0]))<pipelim[0][0]) { //pass 149 if ((sqrt((y_par[1]*y_par[1])+(x_par[1]*x_par[1]))<pipelim[1][0]) && (fabs(y_par[1])>detlim[1][1] || x_par[1]>detlim[1][0])) { accettanza = 1; } } //pass TAN@141 if (fabs(y_par[0])<pipelim[0][1] && sqrt((y_par[0])*(y_par[0])+(x_par[0])*(x_par[0]))<pipelim[0][0]) { //pass Q5@198 if (fabs(y_par[2])<pipelim[2][1] && sqrt((y_par[2]*y_par[2])+(x_par[2]*x_par[2]))<pipelim[2][0]) { //pass 220 if ((sqrt((y_par[3]*y_par[3])+(x_par[3]*x_par[3]))<pipelim[3][0]) && (fabs(y_par[3])>detlim[3][1] || x_par[3]>detlim[3][0])) { accettanza = 1; PosEvo.setX(1000*x_par[3]); PosEvo.setY(1000*y_par[3]); PosEvo.setZ(1000*rp[3]); if(Pos.z()<vz)PosEvo.setZ(-1000*rp[3]); } } } /* LogDebug("ForwardSim") << "\n" << "ACCETTANZA: "<<accettanza << "\n" << "CSI: "<< csi << "\n" << "Theta_X: " << ThetaX << "\n" << "Theta_Y: " << ThetaY << "\n" << "X_at_0: "<< X_at_0 << "\n" << "Y_at_0: "<< Y_at_0 << "\n" << "x_par[3]: "<< x_par[3] << "\n" << "y_par[3]: "<< y_par[3] << "\n" << "pos " << Pos.x() << " " << Pos.y() << " " << Pos.z() << "\n" << "V "<< vx << " " << vy << " " << vz << "\n" */ // -------------- return PosEvo; }
void TotemSD::PrintAll | ( | ) | [virtual] |
bool TotemSD::ProcessHits | ( | G4Step * | aStep, |
G4TouchableHistory * | |||
) | [virtual] |
Implements SensitiveDetector.
Definition at line 101 of file TotemSD.cc.
References CreateNewHit(), CreateNewHitEvo(), edeposit, GetStepInfo(), HitExists(), NULL, ParentId, ParticleType, and unitID.
{ if (aStep == NULL) { return true; } else { GetStepInfo(aStep); if (HitExists() == false && edeposit>0.) { CreateNewHit(); return true; } if (HitExists() == false && (((unitID==1111 || unitID==2222) && ParentId==0 && ParticleType==2212))) { CreateNewHitEvo(); return true; } } return true; }
void TotemSD::ResetForNewPrimary | ( | ) | [private] |
Definition at line 516 of file TotemSD.cc.
References entrancePoint, hitPoint, incidentEnergy, preStepPoint, and SetToLocal().
Referenced by HitExists().
{ entrancePoint = SetToLocal(hitPoint); incidentEnergy = preStepPoint->GetKineticEnergy(); }
uint32_t TotemSD::setDetUnitId | ( | G4Step * | aStep | ) | [virtual] |
Implements SensitiveDetector.
Definition at line 120 of file TotemSD.cc.
References TotemVDetectorOrganization::GetUnitID(), and numberingScheme.
Referenced by GetStepInfo().
{ return (numberingScheme == 0 ? 0 : numberingScheme->GetUnitID(aStep)); }
G4ThreeVector TotemSD::SetToLocal | ( | G4ThreeVector | globalPoint | ) | [private] |
Definition at line 194 of file TotemSD.cc.
References preStepPoint.
Referenced by ResetForNewPrimary().
{ G4ThreeVector localPoint; const G4VTouchable* touch= preStepPoint->GetTouchable(); localPoint = touch->GetHistory()->GetTopTransform().TransformPoint(global); return localPoint; }
void TotemSD::StoreHit | ( | TotemG4Hit * | hit | ) | [private] |
Definition at line 505 of file TotemSD.cc.
Referenced by CreateNewHit(), and CreateNewHitEvo().
{ if (primID<0) return; if (hit == 0 ) { edm::LogWarning("ForwardSim") << "TotemSD: hit to be stored is NULL !!"; return; } theHC->insert( hit ); }
void TotemSD::Summarize | ( | ) | [private] |
void TotemSD::update | ( | const ::EndOfEvent * | ) | [private] |
Definition at line 187 of file TotemSD.cc.
{ }
void TotemSD::update | ( | const BeginOfEvent * | ) | [private, virtual] |
This routine will be called when the appropriate signal arrives.
Implements Observer< const BeginOfEvent * >.
Definition at line 180 of file TotemSD.cc.
References clearHits(), eventno, and LogDebug.
void TotemSD::UpdateHit | ( | ) | [private] |
Definition at line 473 of file TotemSD.cc.
References currentHit, Eloss, LogDebug, postStepPoint, previousUnitID, primaryID, primID, TotemG4Hit::setEnergyLoss(), tsID, tSliceID, and unitID.
Referenced by CreateNewHit(), CreateNewHitEvo(), and HitExists().
{ // if (Eloss > 0.) { // currentHit->addEnergyDeposit(edepositEM,edepositHAD); #ifdef debug LogDebug("ForwardSim") << "G4TotemT1SD updateHit: add eloss " << Eloss << "\nCurrentHit=" << currentHit << ", PostStepPoint=" << postStepPoint->GetPosition(); #endif currentHit->setEnergyLoss(Eloss); } // if(PostStepPoint->GetPhysicalVolume() != CurrentPV){ // currentHit->setExitPoint(SetToLocal(postStepPoint->GetPosition())); // Local3DPoint exit=currentHit->exitPoint(); /* #ifdef debug LogDebug("ForwardSim") << "G4TotemT1SD updateHit: exit point " << exit.x() << " " << exit.y() << " " << exit.z(); // LogDebug("ForwardSim") << "Energy deposit in Unit " << unitID << " em " << edepositEM/MeV // << " hadronic " << edepositHAD/MeV << " MeV"; #endif */ // buffer for next steps: tsID = tSliceID; primID = primaryID; previousUnitID = unitID; }
TotemG4Hit* TotemSD::currentHit [private] |
Definition at line 103 of file TotemSD.h.
Referenced by CreateNewHit(), CreateNewHitEvo(), HitExists(), and UpdateHit().
G4VPhysicalVolume* TotemSD::currentPV [private] |
Definition at line 105 of file TotemSD.h.
Referenced by CreateNewHit(), and GetStepInfo().
float TotemSD::edeposit [private] |
Definition at line 112 of file TotemSD.h.
Referenced by GetStepInfo(), and ProcessHits().
float TotemSD::Eloss [private] |
Definition at line 118 of file TotemSD.h.
Referenced by CreateNewHit(), CreateNewHitEvo(), GetStepInfo(), and UpdateHit().
G4ThreeVector TotemSD::entrancePoint [private] |
Definition at line 93 of file TotemSD.h.
Referenced by ResetForNewPrimary().
int TotemSD::eventno [private] |
G4int TotemSD::hcID [private] |
Definition at line 98 of file TotemSD.h.
Referenced by Initialize().
G4ThreeVector TotemSD::hitPoint [private] |
Definition at line 113 of file TotemSD.h.
Referenced by GetStepInfo(), and ResetForNewPrimary().
float TotemSD::incidentEnergy [private] |
Definition at line 94 of file TotemSD.h.
Referenced by CreateNewHit(), CreateNewHitEvo(), and ResetForNewPrimary().
std::string TotemSD::name [private] |
Reimplemented from SensitiveDetector.
Definition at line 97 of file TotemSD.h.
Referenced by GetStepInfo(), and Initialize().
Definition at line 86 of file TotemSD.h.
Referenced by setDetUnitId(), TotemSD(), and ~TotemSD().
float TotemSD::Pabs [private] |
Definition at line 116 of file TotemSD.h.
Referenced by CreateNewHit(), CreateNewHitEvo(), and GetStepInfo().
int TotemSD::ParentId [private] |
Definition at line 124 of file TotemSD.h.
Referenced by CreateNewHit(), CreateNewHitEvo(), GetStepInfo(), and ProcessHits().
short TotemSD::ParticleType [private] |
Definition at line 119 of file TotemSD.h.
Referenced by CreateNewHit(), CreateNewHitEvo(), GetStepInfo(), and ProcessHits().
float TotemSD::PhiAtEntry [private] |
Definition at line 122 of file TotemSD.h.
Referenced by CreateNewHit(), CreateNewHitEvo(), and GetStepInfo().
G4ThreeVector TotemSD::Posizio [private] |
Definition at line 115 of file TotemSD.h.
Referenced by CreateNewHit(), CreateNewHitEvo(), and GetStepInfo().
G4StepPoint* TotemSD::postStepPoint [private] |
Definition at line 111 of file TotemSD.h.
Referenced by GetStepInfo(), and UpdateHit().
G4StepPoint* TotemSD::preStepPoint [private] |
Definition at line 110 of file TotemSD.h.
Referenced by GetStepInfo(), ResetForNewPrimary(), and SetToLocal().
uint32_t TotemSD::previousUnitID [private] |
Definition at line 106 of file TotemSD.h.
Referenced by HitExists(), and UpdateHit().
int TotemSD::primaryID [private] |
Definition at line 107 of file TotemSD.h.
Referenced by CreateNewHit(), CreateNewHitEvo(), GetStepInfo(), HitExists(), and UpdateHit().
G4int TotemSD::primID [private] |
Definition at line 95 of file TotemSD.h.
Referenced by HitExists(), Initialize(), StoreHit(), and UpdateHit().
TrackingSlaveSD* TotemSD::slave [private] |
Definition at line 85 of file TotemSD.h.
Referenced by clearHits(), EndOfEvent(), fillHits(), TotemSD(), and ~TotemSD().
TotemG4HitCollection* TotemSD::theHC [private] |
Definition at line 99 of file TotemSD.h.
Referenced by EndOfEvent(), HitExists(), Initialize(), PrintAll(), and StoreHit().
const SimTrackManager* TotemSD::theManager [private] |
float TotemSD::ThetaAtEntry [private] |
Definition at line 121 of file TotemSD.h.
Referenced by CreateNewHit(), CreateNewHitEvo(), and GetStepInfo().
G4Track* TotemSD::theTrack [private] |
Definition at line 104 of file TotemSD.h.
Referenced by CreateNewHit(), and GetStepInfo().
float TotemSD::Tof [private] |
Definition at line 117 of file TotemSD.h.
Referenced by CreateNewHit(), CreateNewHitEvo(), and GetStepInfo().
int TotemSD::tsID [private] |
Definition at line 102 of file TotemSD.h.
Referenced by HitExists(), Initialize(), and UpdateHit().
double TotemSD::tSlice [private] |
Definition at line 108 of file TotemSD.h.
Referenced by CreateNewHit(), CreateNewHitEvo(), and GetStepInfo().
int TotemSD::tSliceID [private] |
Definition at line 107 of file TotemSD.h.
Referenced by CreateNewHit(), GetStepInfo(), HitExists(), and UpdateHit().
uint32_t TotemSD::unitID [private] |
Definition at line 106 of file TotemSD.h.
Referenced by CreateNewHit(), CreateNewHitEvo(), GetStepInfo(), HitExists(), ProcessHits(), and UpdateHit().
float TotemSD::Vx [private] |
Definition at line 125 of file TotemSD.h.
Referenced by CreateNewHit(), CreateNewHitEvo(), and GetStepInfo().
float TotemSD::Vy [private] |
Definition at line 125 of file TotemSD.h.
Referenced by CreateNewHit(), CreateNewHitEvo(), and GetStepInfo().
float TotemSD::Vz [private] |
Definition at line 125 of file TotemSD.h.
Referenced by CreateNewHit(), CreateNewHitEvo(), and GetStepInfo().