Skip to content

Commit 9dd2e0c

Browse files
committed
fixed ECAL GFlash
1 parent 4f16b1d commit 9dd2e0c

File tree

3 files changed

+50
-33
lines changed

3 files changed

+50
-33
lines changed

SimG4CMS/Forward/src/MtdSD.cc

Lines changed: 3 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -131,6 +131,9 @@ int MtdSD::getTrackID(const G4Track* aTrack) {
131131
edm::LogVerbatim("MtdSim") << "MtdSD: Track ID: " << aTrack->GetTrackID()
132132
<< " ETL Track ID: " << trkInfo->mcTruthID() << ":" << theID;
133133
#endif
134+
// In the case of ECAL GFlash fast spot may be inside MTD and should be ignored
135+
} else if (rname == "EcalRegion") {
136+
theID = -2;
134137
} else {
135138
throw cms::Exception("MtdSDError") << "MtdSD called in incorrect region " << rname;
136139
}

SimG4CMS/Forward/src/TimingSD.cc

Lines changed: 8 additions & 3 deletions
Original file line numberDiff line numberDiff line change
@@ -104,19 +104,25 @@ bool TimingSD::ProcessHits(G4Step* aStep, G4TouchableHistory*) {
104104
edeposit = aStep->GetTotalEnergyDeposit();
105105
if (edeposit > 0.f) {
106106
getStepInfo(aStep);
107-
if (!hitExists(aStep)) {
107+
// (primaryID = -2) means ECAL GFlash spot inside MTD
108+
if (-1 <= primaryID && !hitExists(aStep)) {
108109
createNewHit(aStep);
109110
}
110111
}
111112
return true;
112113
}
113114

114115
void TimingSD::getStepInfo(const G4Step* aStep) {
116+
const G4Track* newTrack = aStep->GetTrack();
117+
// exclude ECAL gflash spots inside MTD detectors,
118+
// which not possible correctly handle
119+
primaryID = getTrackID(newTrack);
120+
if (primaryID < -1) { return; }
121+
115122
preStepPoint = aStep->GetPreStepPoint();
116123
postStepPoint = aStep->GetPostStepPoint();
117124
hitPointExit = postStepPoint->GetPosition();
118125
setToLocal(preStepPoint, hitPointExit, hitPointLocalExit);
119-
const G4Track* newTrack = aStep->GetTrack();
120126

121127
// neutral particles deliver energy post step
122128
// charged particle start deliver energy pre step
@@ -185,7 +191,6 @@ void TimingSD::getStepInfo(const G4Step* aStep) {
185191

186192
setHitClassID(aStep);
187193
unitID = setDetUnitId(aStep);
188-
primaryID = getTrackID(theTrack);
189194
}
190195

191196
bool TimingSD::hitExists(const G4Step* aStep) {

SimG4Core/Application/src/LowEnergyFastSimModel.cc

Lines changed: 39 additions & 30 deletions
Original file line numberDiff line numberDiff line change
@@ -23,24 +23,24 @@ LowEnergyFastSimModel::LowEnergyFastSimModel(const G4String& name, G4Region* reg
2323
: G4VFastSimulationModel(name, region),
2424
fRegion(region),
2525
fTrackingAction(nullptr),
26-
fCheck(true),
26+
fCheck(false),
2727
fTailPos(0., 0., 0.) {
2828
fEmax = parSet.getParameter<double>("LowEnergyGflashEcalEmax") * CLHEP::GeV;
2929
fPositron = G4Positron::Positron();
3030
fMaterial = nullptr;
3131
auto table = G4Material::GetMaterialTable();
3232
for (auto const& mat : *table) {
33-
G4String nam = (G4String)(DD4hep2DDDName::noNameSpace(mat->GetName()));
34-
size_t n = nam.size();
33+
const G4String& nam = mat->GetName();
34+
std::size_t n = nam.size();
3535
if (n > 4) {
36-
G4String sn = nam.substr(n - 5, 5);
36+
const G4String& sn = nam.substr(n - 5, 5);
3737
if (sn == "PbWO4") {
3838
fMaterial = mat;
3939
break;
4040
}
4141
}
4242
}
43-
G4String nm = (nullptr == fMaterial) ? "not found!" : (G4String)(DD4hep2DDDName::noNameSpace(fMaterial->GetName()));
43+
const G4String& nm = (nullptr == fMaterial) ? "not found!" : fMaterial->GetName();
4444
edm::LogVerbatim("LowEnergyFastSimModel") << "LowEGFlash material: <" << nm << ">";
4545
}
4646

@@ -72,8 +72,6 @@ G4bool LowEnergyFastSimModel::ModelTrigger(const G4FastTrack& fastTrack) {
7272
}
7373

7474
void LowEnergyFastSimModel::DoIt(const G4FastTrack& fastTrack, G4FastStep& fastStep) {
75-
fastStep.KillPrimaryTrack();
76-
fastStep.ProposePrimaryTrackPathLength(0.0);
7775
auto track = fastTrack.GetPrimaryTrack();
7876
G4double energy = track->GetKineticEnergy() * scaleFactor;
7977

@@ -87,41 +85,52 @@ void LowEnergyFastSimModel::DoIt(const G4FastTrack& fastTrack, G4FastStep& fastS
8785

8886
const G4ThreeVector& momDir = track->GetMomentumDirection();
8987

90-
// in point energy deposition
91-
GFlashEnergySpot spot;
92-
spot.SetEnergy(inPointEnergy);
93-
spot.SetPosition(pos);
94-
fHitMaker.make(&spot, &fastTrack);
95-
9688
// Russian roulette
9789
double wt2 = track->GetWeight();
9890
if (wt2 <= 0.0) {
9991
wt2 = 1.0;
10092
}
101-
10293
// tail energy deposition
103-
const G4double etail = energy - inPointEnergy;
104-
const G4int nspots = etail;
105-
const G4double tailEnergy = etail * wt2 / (nspots + 1);
94+
G4double etail = energy - inPointEnergy;
95+
G4int nspots = G4lrint(etail);
96+
if (0 == nspots) {
97+
inPointEnergy = energy;
98+
etail = 0.0;
99+
}
100+
101+
// in point energy deposition
102+
GFlashEnergySpot spot;
103+
spot.SetEnergy(inPointEnergy*wt2);
104+
spot.SetPosition(pos);
105+
fHitMaker.make(&spot, &fastTrack);
106+
107+
G4double zz = 0.0;
108+
if (0 < nspots) {
109+
etail *= wt2 / (G4double)nspots;
106110
/*
107111
edm::LogVerbatim("LowEnergyFastSimModel") << track->GetDefinition()->GetParticleName()
108112
<< " Ekin(MeV)=" << energy << " material: <"
109113
<< track->GetMaterial()->GetName()
110114
<< "> Elocal=" << inPointEnergy
111115
<< " Etail=" << tailEnergy
112-
<< " Nspots=" << nspots+1;
116+
<< " Nspots=" << nspots;
113117
*/
114-
for (G4int i = 0; i <= nspots; ++i) {
115-
const G4double r = fParam.GetRadius(energy);
116-
const G4double z = fParam.GetZ();
117-
118-
const G4double phi = CLHEP::twopi * G4UniformRand();
119-
fTailPos.set(r * std::cos(phi), r * std::sin(phi), z);
120-
fTailPos.rotateUz(momDir);
121-
fTailPos += pos;
122-
123-
spot.SetEnergy(tailEnergy);
124-
spot.SetPosition(fTailPos);
125-
fHitMaker.make(&spot, &fastTrack);
118+
for (G4int i = 0; i < nspots; ++i) {
119+
const G4double r = fParam.GetRadius(energy);
120+
const G4double z = fParam.GetZ();
121+
zz += z;
122+
123+
const G4double phi = CLHEP::twopi * G4UniformRand();
124+
fTailPos.set(r * std::cos(phi), r * std::sin(phi), z);
125+
fTailPos.rotateUz(momDir);
126+
fTailPos += pos;
127+
128+
spot.SetEnergy(etail);
129+
spot.SetPosition(fTailPos);
130+
fHitMaker.make(&spot, &fastTrack);
131+
}
132+
zz /= (G4double)nspots;
126133
}
134+
fastStep.KillPrimaryTrack();
135+
fastStep.ProposePrimaryTrackPathLength(zz);
127136
}

0 commit comments

Comments
 (0)