Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Add FlatPtGeneration to CloseByParticleGun #44008

Merged
merged 2 commits into from
Feb 29, 2024
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
12 changes: 8 additions & 4 deletions Configuration/Generator/python/CE_E_Front_120um_cfi.py
Original file line number Diff line number Diff line change
Expand Up @@ -3,9 +3,10 @@
generator = cms.EDProducer("CloseByParticleGunProducer",
PGunParameters = cms.PSet(PartID = cms.vint32(22),
ControlledByEta = cms.bool(False),
EnMin = cms.double(25.),
EnMax = cms.double(200.),
MaxEnSpread = cms.bool(False),
VarMin = cms.double(25.),
VarMax = cms.double(200.),
MaxVarSpread = cms.bool(False),
FlatPtGeneration = cms.bool(False),
RMin = cms.double(54.99),
RMax = cms.double(55.01),
ZMin = cms.double(320.99),
Expand All @@ -19,7 +20,10 @@
MinEta = cms.double(1.7),
MaxPhi = cms.double(3.14159265359),
MinPhi = cms.double(-3.14159265359),

UseDeltaT = cms.bool(False),
TMin = cms.double(0.),
TMax = cms.double(0.05),
OffsetFirst = cms.double(0.)
),
Verbosity = cms.untracked.int32(0),

Expand Down
12 changes: 8 additions & 4 deletions Configuration/Generator/python/CE_E_Front_200um_cfi.py
Original file line number Diff line number Diff line change
Expand Up @@ -3,9 +3,10 @@
generator = cms.EDProducer("CloseByParticleGunProducer",
PGunParameters = cms.PSet(PartID = cms.vint32(22),
ControlledByEta = cms.bool(False),
EnMin = cms.double(25.),
EnMax = cms.double(200.),
MaxEnSpread = cms.bool(False),
VarMin = cms.double(25.),
VarMax = cms.double(200.),
MaxVarSpread = cms.bool(False),
FlatPtGeneration = cms.bool(False),
RMin = cms.double(89.99),
RMax = cms.double(90.01),
ZMin = cms.double(320.99),
Expand All @@ -19,7 +20,10 @@
MinEta = cms.double(1.7),
MaxPhi = cms.double(3.14159265359),
MinPhi = cms.double(-3.14159265359),

UseDeltaT = cms.bool(False),
TMin = cms.double(0.),
TMax = cms.double(0.05),
OffsetFirst = cms.double(0.)
),
Verbosity = cms.untracked.int32(0),

Expand Down
12 changes: 8 additions & 4 deletions Configuration/Generator/python/CE_E_Front_300um_cfi.py
Original file line number Diff line number Diff line change
Expand Up @@ -3,9 +3,10 @@
generator = cms.EDProducer("CloseByParticleGunProducer",
PGunParameters = cms.PSet(PartID = cms.vint32(22),
ControlledByEta = cms.bool(False),
EnMin = cms.double(25.),
EnMax = cms.double(200.),
MaxEnSpread = cms.bool(False),
VarMin = cms.double(25.),
VarMax = cms.double(200.),
MaxVarSpread = cms.bool(False),
FlatPtGeneration = cms.bool(False),
RMin = cms.double(134.99),
RMax = cms.double(135.01),
ZMin = cms.double(320.99),
Expand All @@ -19,7 +20,10 @@
MinEta = cms.double(1.7),
MaxPhi = cms.double(3.14159265359),
MinPhi = cms.double(-3.14159265359),

UseDeltaT = cms.bool(False),
TMin = cms.double(0.),
TMax = cms.double(0.05),
OffsetFirst = cms.double(0.)
),
Verbosity = cms.untracked.int32(0),

Expand Down
12 changes: 8 additions & 4 deletions Configuration/Generator/python/CE_H_Coarse_300um_cfi.py
Original file line number Diff line number Diff line change
Expand Up @@ -3,9 +3,10 @@
generator = cms.EDProducer("CloseByParticleGunProducer",
PGunParameters = cms.PSet(PartID = cms.vint32(22),
ControlledByEta = cms.bool(False),
EnMin = cms.double(25.),
EnMax = cms.double(200.),
MaxEnSpread = cms.bool(False),
VarMin = cms.double(25.),
VarMax = cms.double(200.),
MaxVarSpread = cms.bool(False),
FlatPtGeneration = cms.bool(False),
RMin = cms.double(79.99),
RMax = cms.double(80.01),
ZMin = cms.double(429.99),
Expand All @@ -19,7 +20,10 @@
MinEta = cms.double(1.7),
MaxPhi = cms.double(3.14159265359),
MinPhi = cms.double(-3.14159265359),

UseDeltaT = cms.bool(False),
TMin = cms.double(0.),
TMax = cms.double(0.05),
OffsetFirst = cms.double(0.)
),
Verbosity = cms.untracked.int32(0),

Expand Down
12 changes: 8 additions & 4 deletions Configuration/Generator/python/CE_H_Coarse_Scint_cfi.py
Original file line number Diff line number Diff line change
Expand Up @@ -3,9 +3,10 @@
generator = cms.EDProducer("CloseByParticleGunProducer",
PGunParameters = cms.PSet(PartID = cms.vint32(22),
ControlledByEta = cms.bool(False),
EnMin = cms.double(25.),
EnMax = cms.double(200.),
MaxEnSpread = cms.bool(False),
VarMin = cms.double(25.),
VarMax = cms.double(200.),
MaxVarSpread = cms.bool(False),
FlatPtGeneration = cms.bool(False),
RMin = cms.double(179.99),
RMax = cms.double(180.01),
ZMin = cms.double(429.99),
Expand All @@ -19,7 +20,10 @@
MinEta = cms.double(1.7),
MaxPhi = cms.double(3.14159265359),
MinPhi = cms.double(-3.14159265359),

UseDeltaT = cms.bool(False),
TMin = cms.double(0.),
TMax = cms.double(0.05),
OffsetFirst = cms.double(0.)
),
Verbosity = cms.untracked.int32(0),

Expand Down
12 changes: 8 additions & 4 deletions Configuration/Generator/python/CE_H_Fine_120um_cfi.py
Original file line number Diff line number Diff line change
Expand Up @@ -3,9 +3,10 @@
generator = cms.EDProducer("CloseByParticleGunProducer",
PGunParameters = cms.PSet(PartID = cms.vint32(22),
ControlledByEta = cms.bool(False),
EnMin = cms.double(25.),
EnMax = cms.double(200.),
MaxEnSpread = cms.bool(False),
VarMin = cms.double(25.),
VarMax = cms.double(200.),
MaxVarSpread = cms.bool(False),
FlatPtGeneration = cms.bool(False),
RMin = cms.double(49.99),
RMax = cms.double(50.01),
ZMin = cms.double(362.519),
Expand All @@ -19,7 +20,10 @@
MinEta = cms.double(1.7),
MaxPhi = cms.double(3.14159265359),
MinPhi = cms.double(-3.14159265359),

UseDeltaT = cms.bool(False),
TMin = cms.double(0.),
TMax = cms.double(0.05),
OffsetFirst = cms.double(0.)
),
Verbosity = cms.untracked.int32(0),

Expand Down
12 changes: 8 additions & 4 deletions Configuration/Generator/python/CE_H_Fine_200um_cfi.py
Original file line number Diff line number Diff line change
Expand Up @@ -3,9 +3,10 @@
generator = cms.EDProducer("CloseByParticleGunProducer",
PGunParameters = cms.PSet(PartID = cms.vint32(22),
ControlledByEta = cms.bool(False),
EnMin = cms.double(25.),
EnMax = cms.double(200.),
MaxEnSpread = cms.bool(False),
VarMin = cms.double(25.),
VarMax = cms.double(200.),
MaxVarSpread = cms.bool(False),
FlatPtGeneration = cms.bool(False),
RMin = cms.double(89.99),
RMax = cms.double(90.01),
ZMin = cms.double(362.519),
Expand All @@ -19,7 +20,10 @@
MinEta = cms.double(1.7),
MaxPhi = cms.double(3.14159265359),
MinPhi = cms.double(-3.14159265359),

UseDeltaT = cms.bool(False),
TMin = cms.double(0.),
TMax = cms.double(0.05),
OffsetFirst = cms.double(0.)
),
Verbosity = cms.untracked.int32(0),

Expand Down
12 changes: 8 additions & 4 deletions Configuration/Generator/python/CE_H_Fine_300um_cfi.py
Original file line number Diff line number Diff line change
Expand Up @@ -3,9 +3,10 @@
generator = cms.EDProducer("CloseByParticleGunProducer",
PGunParameters = cms.PSet(PartID = cms.vint32(22),
ControlledByEta = cms.bool(False),
EnMin = cms.double(25.),
EnMax = cms.double(200.),
MaxEnSpread = cms.bool(False),
VarMin = cms.double(25.),
VarMax = cms.double(200.),
MaxVarSpread = cms.bool(False),
FlatPtGeneration = cms.bool(False),
RMin = cms.double(134.99),
RMax = cms.double(135.01),
ZMin = cms.double(362.519),
Expand All @@ -19,7 +20,10 @@
MinEta = cms.double(1.7),
MaxPhi = cms.double(3.14159265359),
MinPhi = cms.double(-3.14159265359),

UseDeltaT = cms.bool(False),
TMin = cms.double(0.),
TMax = cms.double(0.05),
OffsetFirst = cms.double(0.)
),
Verbosity = cms.untracked.int32(0),

Expand Down
5 changes: 3 additions & 2 deletions IOMC/ParticleGuns/interface/CloseByParticleGunProducer.h
Original file line number Diff line number Diff line change
Expand Up @@ -26,10 +26,11 @@ namespace edm {
protected:
// data members
bool fControlledByEta;
double fEnMin, fEnMax, fEtaMin, fEtaMax, fRMin, fRMax, fZMin, fZMax, fDelta, fPhiMin, fPhiMax, fTMin, fTMax,
double fVarMin, fVarMax, fEtaMin, fEtaMax, fRMin, fRMax, fZMin, fZMax, fDelta, fPhiMin, fPhiMax, fTMin, fTMax,
fOffsetFirst;
int fNParticles;
bool fMaxEnSpread = false;
bool fMaxVarSpread = false;
bool fFlatPtGeneration = false;
bool fPointing = false;
bool fOverlapping = false;
bool fRandomShoot = false;
Expand Down
64 changes: 41 additions & 23 deletions IOMC/ParticleGuns/src/CloseByParticleGunProducer.cc
Original file line number Diff line number Diff line change
Expand Up @@ -26,13 +26,13 @@ CloseByParticleGunProducer::CloseByParticleGunProducer(const ParameterSet& pset)
: BaseFlatGunProducer(pset), m_fieldToken(esConsumes()) {
ParameterSet pgun_params = pset.getParameter<ParameterSet>("PGunParameters");
fControlledByEta = pgun_params.getParameter<bool>("ControlledByEta");
fEnMax = pgun_params.getParameter<double>("EnMax");
fEnMin = pgun_params.getParameter<double>("EnMin");
if (fEnMin < 1)
fVarMax = pgun_params.getParameter<double>("VarMax");
fVarMin = pgun_params.getParameter<double>("VarMin");
fMaxVarSpread = pgun_params.getParameter<bool>("MaxVarSpread");
fFlatPtGeneration = pgun_params.getParameter<bool>("FlatPtGeneration");
if (fVarMin < 1 && !fFlatPtGeneration)
LogError("CloseByParticleGunProducer") << " Please choose a minimum energy greater than 1 GeV, otherwise time "
"information may be invalid or not reliable";

fMaxEnSpread = pgun_params.getParameter<bool>("MaxEnSpread");
if (fControlledByEta) {
fEtaMax = pgun_params.getParameter<double>("MaxEta");
fEtaMin = pgun_params.getParameter<double>("MinEta");
Expand All @@ -51,6 +51,9 @@ CloseByParticleGunProducer::CloseByParticleGunProducer(const ParameterSet& pset)
fPhiMax = pgun_params.getParameter<double>("MaxPhi");
fPointing = pgun_params.getParameter<bool>("Pointing");
fOverlapping = pgun_params.getParameter<bool>("Overlapping");
if (fFlatPtGeneration && !fPointing)
LogError("CloseByParticleGunProducer")
<< " Can't generate non pointing FlatPt samples; please disable FlatPt generation or generate pointing sample";
fRandomShoot = pgun_params.getParameter<bool>("RandomShoot");
fNParticles = pgun_params.getParameter<int>("NParticles");
fPartIDs = pgun_params.getParameter<vector<int>>("PartID");
Expand Down Expand Up @@ -79,9 +82,10 @@ void CloseByParticleGunProducer::fillDescriptions(ConfigurationDescriptions& des
edm::ParameterSetDescription psd0;
psd0.add<bool>("ControlledByEta", false);
psd0.add<double>("Delta", 10);
psd0.add<double>("EnMax", 200.0);
psd0.add<double>("EnMin", 25.0);
psd0.add<bool>("MaxEnSpread", false);
psd0.add<double>("VarMax", 200.0);
psd0.add<double>("VarMin", 25.0);
psd0.add<bool>("MaxVarSpread", false);
psd0.add<bool>("FlatPtGeneration", false);
psd0.add<double>("MaxEta", 2.7);
psd0.add<double>("MaxPhi", 3.14159265359);
psd0.add<double>("MinEta", 1.7);
Expand Down Expand Up @@ -126,13 +130,14 @@ void CloseByParticleGunProducer::produce(Event& e, const EventSetup& es) {

double phi = CLHEP::RandFlat::shoot(engine, fPhiMin, fPhiMax);
double fZ = CLHEP::RandFlat::shoot(engine, fZMin, fZMax);
double fR;
double fR, fEta;
double fT;

if (!fControlledByEta) {
fR = CLHEP::RandFlat::shoot(engine, fRMin, fRMax);
fEta = asinh(fZ / fR);
} else {
double fEta = CLHEP::RandFlat::shoot(engine, fEtaMin, fEtaMax);
fEta = CLHEP::RandFlat::shoot(engine, fEtaMin, fEtaMax);
fR = (fZ / sinh(fEta));
}

Expand All @@ -153,26 +158,39 @@ void CloseByParticleGunProducer::produce(Event& e, const EventSetup& es) {
} else
phi += fDelta / fR;

double fEn;
if (numParticles > 1 && fMaxEnSpread)
fEn = fEnMin + ip * (fEnMax - fEnMin) / (numParticles - 1);
double fVar;
if (numParticles > 1 && fMaxVarSpread)
fVar = fVarMin + ip * (fVarMax - fVarMin) / (numParticles - 1);
else
fEn = CLHEP::RandFlat::shoot(engine, fEnMin, fEnMax);
fVar = CLHEP::RandFlat::shoot(engine, fVarMin, fVarMax);

int partIdx = CLHEP::RandFlat::shoot(engine, 0, fPartIDs.size());
int PartID = fPartIDs[partIdx];
const HepPDT::ParticleData* PData = fPDGTable->particle(HepPDT::ParticleID(abs(PartID)));
double mass = PData->mass().value();
double mom2 = fEn * fEn - mass * mass;
double mom = 0.;
if (mom2 > 0.) {
mom = sqrt(mom2);
}
double px = 0.;
double py = 0.;
double pz = mom;
double energy = fEn;

double mom, px, py, pz;
double energy;

if (!fFlatPtGeneration) {
double mom2 = fVar * fVar - mass * mass;
mom = 0.;
if (mom2 > 0.) {
mom = sqrt(mom2);
}
px = 0.;
py = 0.;
pz = mom;
energy = fVar;
} else {
double theta = 2. * atan(exp(-fEta));
mom = fVar / sin(theta);
px = fVar * cos(phi);
py = fVar * sin(phi);
pz = mom * cos(theta);
double energy2 = mom * mom + mass * mass;
energy = sqrt(energy2);
}
// Compute Vertex Position
double x = fR * cos(phi);
double y = fR * sin(phi);
Expand Down