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 new configuration options to CloseByParticleGun #45515

Merged
merged 3 commits into from
Jul 30, 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
4 changes: 3 additions & 1 deletion IOMC/ParticleGuns/interface/CloseByParticleGunProducer.h
Original file line number Diff line number Diff line change
Expand Up @@ -25,10 +25,12 @@ namespace edm {

protected:
// data members
bool fControlledByEta;
bool fControlledByEta, fControlledByREta;
double fVarMin, fVarMax, fEtaMin, fEtaMax, fRMin, fRMax, fZMin, fZMax, fDelta, fPhiMin, fPhiMax, fTMin, fTMax,
fOffsetFirst;
double log_fVarMin = 0., log_fVarMax = 0.;
int fNParticles;
bool fLogSpacedVar = false;
bool fMaxVarSpread = false;
bool fFlatPtGeneration = false;
bool fPointing = false;
Expand Down
62 changes: 48 additions & 14 deletions IOMC/ParticleGuns/src/CloseByParticleGunProducer.cc
Original file line number Diff line number Diff line change
Expand Up @@ -26,33 +26,52 @@ CloseByParticleGunProducer::CloseByParticleGunProducer(const ParameterSet& pset)
: BaseFlatGunProducer(pset), m_fieldToken(esConsumes()) {
ParameterSet pgun_params = pset.getParameter<ParameterSet>("PGunParameters");
fControlledByEta = pgun_params.getParameter<bool>("ControlledByEta");
fControlledByREta = pgun_params.getParameter<bool>("ControlledByREta");
if (fControlledByEta and fControlledByREta)
throw cms::Exception("CloseByParticleGunProducer")
<< " Conflicting configuration, cannot have both ControlledByEta and ControlledByREta ";

fVarMax = pgun_params.getParameter<double>("VarMax");
fVarMin = pgun_params.getParameter<double>("VarMin");
fMaxVarSpread = pgun_params.getParameter<bool>("MaxVarSpread");
fLogSpacedVar = pgun_params.getParameter<bool>("LogSpacedVar");
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";
if (fControlledByEta) {
throw cms::Exception("CloseByParticleGunProducer")
<< " Please choose a minimum energy greater than 1 GeV, otherwise time "
"information may be invalid or not reliable";
if (fVarMin < 0 && fLogSpacedVar)
throw cms::Exception("CloseByParticleGunProducer") << " Minimum energy must be greater than zero for log spacing";
else {
log_fVarMin = std::log(fVarMin);
log_fVarMax = std::log(fVarMax);
}

if (fControlledByEta || fControlledByREta) {
fEtaMax = pgun_params.getParameter<double>("MaxEta");
fEtaMin = pgun_params.getParameter<double>("MinEta");
if (fEtaMax <= fEtaMin)
LogError("CloseByParticleGunProducer") << " Please fix MinEta and MaxEta values in the configuration";
throw cms::Exception("CloseByParticleGunProducer") << " Please fix MinEta and MaxEta values in the configuration";
} else {
fRMax = pgun_params.getParameter<double>("RMax");
fRMin = pgun_params.getParameter<double>("RMin");
if (fRMax <= fRMin)
LogError("CloseByParticleGunProducer") << " Please fix RMin and RMax values in the configuration";
throw cms::Exception("CloseByParticleGunProducer") << " Please fix RMin and RMax values in the configuration";
}
if (!fControlledByREta) {
fZMax = pgun_params.getParameter<double>("ZMax");
fZMin = pgun_params.getParameter<double>("ZMin");

if (fZMax <= fZMin)
throw cms::Exception("CloseByParticleGunProducer") << " Please fix ZMin and ZMax values in the configuration";
}
fZMax = pgun_params.getParameter<double>("ZMax");
fZMin = pgun_params.getParameter<double>("ZMin");
fDelta = pgun_params.getParameter<double>("Delta");
fPhiMin = pgun_params.getParameter<double>("MinPhi");
fPhiMax = pgun_params.getParameter<double>("MaxPhi");
fPointing = pgun_params.getParameter<bool>("Pointing");
fOverlapping = pgun_params.getParameter<bool>("Overlapping");
if (fFlatPtGeneration && !fPointing)
LogError("CloseByParticleGunProducer")
throw cms::Exception("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");
Expand All @@ -63,7 +82,7 @@ CloseByParticleGunProducer::CloseByParticleGunProducer(const ParameterSet& pset)
fTMax = pgun_params.getParameter<double>("TMax");
fTMin = pgun_params.getParameter<double>("TMin");
if (fTMax <= fTMin)
LogError("CloseByParticleGunProducer") << " Please fix TMin and TMax values in the configuration";
throw cms::Exception("CloseByParticleGunProducer") << " Please fix TMin and TMax values in the configuration";
// set a fixed time offset for the particles
fOffsetFirst = pgun_params.getParameter<double>("OffsetFirst");

Expand All @@ -81,10 +100,12 @@ void CloseByParticleGunProducer::fillDescriptions(ConfigurationDescriptions& des
{
edm::ParameterSetDescription psd0;
psd0.add<bool>("ControlledByEta", false);
psd0.add<bool>("ControlledByREta", false);
psd0.add<double>("Delta", 10);
psd0.add<double>("VarMax", 200.0);
psd0.add<double>("VarMin", 25.0);
psd0.add<bool>("MaxVarSpread", false);
psd0.add<bool>("LogSpacedVar", false);
psd0.add<bool>("FlatPtGeneration", false);
psd0.add<double>("MaxEta", 2.7);
psd0.add<double>("MaxPhi", 3.14159265359);
Expand Down Expand Up @@ -129,16 +150,24 @@ void CloseByParticleGunProducer::produce(Event& e, const EventSetup& es) {
unsigned int numParticles = fRandomShoot ? CLHEP::RandFlat::shoot(engine, 1, fNParticles) : fNParticles;

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

if (!fControlledByEta) {
fR = CLHEP::RandFlat::shoot(engine, fRMin, fRMax);
fEta = asinh(fZ / fR);
if (!fControlledByREta) {
fZ = CLHEP::RandFlat::shoot(engine, fZMin, fZMax);

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

if (fUseDeltaT) {
Expand All @@ -161,6 +190,11 @@ void CloseByParticleGunProducer::produce(Event& e, const EventSetup& es) {
double fVar;
if (numParticles > 1 && fMaxVarSpread)
fVar = fVarMin + ip * (fVarMax - fVarMin) / (numParticles - 1);
else if (fLogSpacedVar) {
double fVar_log = CLHEP::RandFlat::shoot(engine, log_fVarMin, log_fVarMax);
fVar = std::exp(fVar_log);
}

else
fVar = CLHEP::RandFlat::shoot(engine, fVarMin, fVarMax);

Expand Down