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

Use Lorentz Angle from data base for pixel template CPE's in case of error #33495

Merged
merged 4 commits into from Apr 28, 2021
Merged
Show file tree
Hide file tree
Changes from 3 commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Jump to
Jump to file
Failed to load files.
Diff view
Diff view
15 changes: 15 additions & 0 deletions HLTrigger/Configuration/python/customizeHLTforCMSSW.py
Expand Up @@ -129,11 +129,26 @@ def customiseFor2018Input(process):

return process

def customiseFor33495(process):
"""Customize HLT menu to remove deprecated parameters for pixel Generic and Template CPE's """
for producer in esproducers_by_type(process, "PixelCPEGenericESProducer"):
if hasattr(producer, "DoLorentz"):
del producer.DoLorentz
if hasattr(producer, "useLAAlignmentOffsets"):
del producer.useLAAlignmentOffsets

for producer in esproducers_by_type(process, "PixelCPETemplateRecoESProducer"):
if hasattr(producer, "DoLorentz"):
del producer.DoLorentz
return process



# CMSSW version specific customizations
def customizeHLTforCMSSW(process, menuType="GRun"):

# add call to action function in proper order: newest last!
# process = customiseFor12718(process)
process = customiseFor33495(process)

return process
4 changes: 3 additions & 1 deletion RecoLocalTracker/SiPixelRecHits/interface/PixelCPEBase.h
Expand Up @@ -238,7 +238,9 @@ class PixelCPEBase : public PixelClusterParameterEstimator {
const SiPixelTemplateDBObject* templateDBobject_;
bool alpha2Order; // switch on/off E.B effect.

bool DoLorentz_;
bool useLAFromDB_; //Use LA value from the database (used for generic CPE or in template CPE if an error)

bool doLorentzFromAlignment_;
bool LoadTemplatesFromDB_;

//errors for template reco for edge hits, based on observed residuals from
Expand Down
Expand Up @@ -34,16 +34,17 @@ class PixelCPEClusterRepairESProducer : public edm::ESProducer {
edm::ESGetToken<SiPixel2DTemplateDBObject, SiPixel2DTemplateDBObjectESProducerRcd> templateDBobject2DToken_;

edm::ParameterSet pset_;
bool DoLorentz_;
bool doLorentzFromAlignment_;
bool useLAFromDB_;
};

using namespace edm;

PixelCPEClusterRepairESProducer::PixelCPEClusterRepairESProducer(const edm::ParameterSet& p) {
std::string myname = p.getParameter<std::string>("ComponentName");

//DoLorentz_ = p.getParameter<bool>("DoLorentz"); // True when LA from alignment is used
DoLorentz_ = p.getParameter<bool>("DoLorentz");
useLAFromDB_ = p.getParameter<bool>("useLAFromDB");
doLorentzFromAlignment_ = p.getParameter<bool>("doLorentzFromAlignment");
Comment on lines -46 to +47
Copy link
Contributor

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

there is no matching update in the ::fillDescriptions method.
There is a comment that DoLorentz was specific to this ES producer; since it's no longer needed, please remove it or mark it as obsolete

Copy link
Contributor Author

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Ok I removed the the line from the fillDescriptions method


pset_ = p;
auto c = setWhatProduced(this, myname);
Expand All @@ -52,11 +53,15 @@ PixelCPEClusterRepairESProducer::PixelCPEClusterRepairESProducer(const edm::Para
hTTToken_ = c.consumes();
templateDBobjectToken_ = c.consumes();
templateDBobject2DToken_ = c.consumes();
if (DoLorentz_) {

char const* laLabel = ""; // standard LA, from calibration, label=""
lorentzAngleToken_ = c.consumes(edm::ESInputTag("", laLabel));

if (doLorentzFromAlignment_) {
lorentzAngleToken_ = c.consumes(edm::ESInputTag("", "fromAlignment"));
Copy link
Contributor

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Suggested change
char const* laLabel = ""; // standard LA, from calibration, label=""
lorentzAngleToken_ = c.consumes(edm::ESInputTag("", laLabel));
if (doLorentzFromAlignment_) {
lorentzAngleToken_ = c.consumes(edm::ESInputTag("", "fromAlignment"));
if (useLAFromDB_ || doLorentzFromAlignment_) {
char const* laLabel = doLorentzFromAlignment_ ? "fromAlignment" : "";
lorentzAngleToken_ = c.consumes(edm::ESInputTag("", laLabel));

this is matching the actual use case in ::produce (only one payload is actually consumed; the current implementation can declare two different consumes). However I'm not sure what your actual goal is for the case with both flags set to true.

Copy link
Contributor Author

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Ok I updated to your suggestion

}

//std::cout<<" from ES Producer Templates "<<myname<<" "<<DoLorentz_<<std::endl; //dk
//std::cout<<" from ES Producer Templates "<<myname<<" "<<DoLorentz_<<std::endl; //dk
}
}

PixelCPEClusterRepairESProducer::~PixelCPEClusterRepairESProducer() {}
Expand All @@ -79,10 +84,11 @@ void PixelCPEClusterRepairESProducer::fillDescriptions(edm::ConfigurationDescrip

std::unique_ptr<PixelClusterParameterEstimator> PixelCPEClusterRepairESProducer::produce(
const TkPixelCPERecord& iRecord) {
// Normal, default LA actually is NOT needed
// null is ok becuse LA is not use by templates in this mode
// Normal, default LA is used in case of template failure, load it unless
// turned off
// if turned off, null is ok, becomes zero
const SiPixelLorentzAngle* lorentzAngleProduct = nullptr;
if (DoLorentz_) { // LA correction from alignment
if (useLAFromDB_ || doLorentzFromAlignment_) {
lorentzAngleProduct = &iRecord.get(lorentzAngleToken_);
}

Expand Down
Expand Up @@ -89,8 +89,6 @@ void PixelCPEFastESProducer::fillDescriptions(edm::ConfigurationDescriptions& de
// specific to PixelCPEFastESProducer
desc.add<std::string>("ComponentName", "PixelCPEFast");
desc.add<edm::ESInputTag>("MagneticFieldRecord", edm::ESInputTag());
desc.add<bool>("useLAAlignmentOffsets", false);
desc.add<bool>("DoLorentz", false);

descriptions.add("PixelCPEFastESProducer", desc);
}
Expand Down
Expand Up @@ -45,9 +45,9 @@ PixelCPEGenericESProducer::PixelCPEGenericESProducer(const edm::ParameterSet& p)
// Use LA-width from DB. If both (upper and this) are false LA-width is calcuated from LA-offset
useLAWidthFromDB_ = p.getParameter<bool>("useLAWidthFromDB");
// Use Alignment LA-offset
const bool useLAAlignmentOffsets = p.getParameter<bool>("useLAAlignmentOffsets");
const bool doLorentzFromAlignment = p.getParameter<bool>("doLorentzFromAlignment");
Copy link
Contributor

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

there is a commented out cout on L70 which is using useLAAlignmentOffsets_ (already broken, but still close/recognizable); please update there as well or drop that cout

Copy link
Contributor Author

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Ok I removed the cout

char const* laLabel = ""; // standard LA, from calibration, label=""
if (useLAAlignmentOffsets) {
if (doLorentzFromAlignment) {
laLabel = "fromAlignment";
}

Expand Down Expand Up @@ -108,8 +108,6 @@ void PixelCPEGenericESProducer::fillDescriptions(edm::ConfigurationDescriptions&
// specific to PixelCPEGenericESProducer
desc.add<std::string>("ComponentName", "PixelCPEGeneric");
desc.add<edm::ESInputTag>("MagneticFieldRecord", edm::ESInputTag(""));
desc.add<bool>("useLAAlignmentOffsets", false);
desc.add<bool>("DoLorentz", false);
descriptions.add("_generic_default", desc);
}

Expand Down
Expand Up @@ -30,35 +30,41 @@ class PixelCPETemplateRecoESProducer : public edm::ESProducer {
edm::ESGetToken<SiPixelTemplateDBObject, SiPixelTemplateDBObjectESProducerRcd> templateDBobjectToken_;

edm::ParameterSet pset_;
bool DoLorentz_;
bool doLorentzFromAlignment_;
bool useLAFromDB_;
};

using namespace edm;

PixelCPETemplateRecoESProducer::PixelCPETemplateRecoESProducer(const edm::ParameterSet& p) {
std::string myname = p.getParameter<std::string>("ComponentName");

//DoLorentz_ = p.getParameter<bool>("DoLorentz"); // True when LA from alignment is used
DoLorentz_ = p.getParameter<bool>("DoLorentz");
useLAFromDB_ = p.getParameter<bool>("useLAFromDB");
doLorentzFromAlignment_ = p.getParameter<bool>("doLorentzFromAlignment");

pset_ = p;
auto c = setWhatProduced(this, myname);
magfieldToken_ = c.consumes();
pDDToken_ = c.consumes();
hTTToken_ = c.consumes();
templateDBobjectToken_ = c.consumes();
if (DoLorentz_) {

char const* laLabel = ""; // standard LA, from calibration, label=""
lorentzAngleToken_ = c.consumes(edm::ESInputTag("", laLabel));

if (doLorentzFromAlignment_) {
lorentzAngleToken_ = c.consumes(edm::ESInputTag("", "fromAlignment"));
Copy link
Contributor

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

(copied from PixelCPEClusterRepairESProducer.cc)

Suggested change
char const* laLabel = ""; // standard LA, from calibration, label=""
lorentzAngleToken_ = c.consumes(edm::ESInputTag("", laLabel));
if (doLorentzFromAlignment_) {
lorentzAngleToken_ = c.consumes(edm::ESInputTag("", "fromAlignment"));
if (useLAFromDB_ || doLorentzFromAlignment_) {
char const* laLabel = doLorentzFromAlignment_ ? "fromAlignment" : "";
lorentzAngleToken_ = c.consumes(edm::ESInputTag("", laLabel));

this is matching the actual use case in ::produce (only one payload is actually consumed; the current implementation can declare two different consumes). However I'm not sure what your actual goal is for the case with both flags set to true.

Copy link
Contributor Author

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Ok I updated to your suggestion

}
//std::cout<<" from ES Producer Templates "<<myname<<" "<<DoLorentz_<<std::endl; //dk
}

std::unique_ptr<PixelClusterParameterEstimator> PixelCPETemplateRecoESProducer::produce(
const TkPixelCPERecord& iRecord) {
// Normal, deafult LA actually is NOT needed
// null is ok becuse LA is not use by templates in this mode
// Normal, default LA is used in case of template failure, load it unless
// turned off
// if turned off, null is ok, becomes zero
const SiPixelLorentzAngle* lorentzAngleProduct = nullptr;
if (DoLorentz_) { // LA correction from alignment
if (useLAFromDB_ || doLorentzFromAlignment_) {
lorentzAngleProduct = &iRecord.get(lorentzAngleToken_);
}

Expand All @@ -81,7 +87,6 @@ void PixelCPETemplateRecoESProducer::fillDescriptions(edm::ConfigurationDescript

// specific to PixelCPETemplateRecoESProducer
desc.add<std::string>("ComponentName", "PixelCPETemplateReco");
desc.add<bool>("DoLorentz", true);
descriptions.add("_templates_default", desc);
}

Expand Down
8 changes: 6 additions & 2 deletions RecoLocalTracker/SiPixelRecHits/src/PixelCPEBase.cc
Expand Up @@ -106,7 +106,8 @@ PixelCPEBase::PixelCPEBase(edm::ParameterSet const& conf,

// For Templates only
// Compute the Lorentz shifts for this detector element for templates (from Alignment)
DoLorentz_ = conf.getParameter<bool>("DoLorentz");
doLorentzFromAlignment_ = conf.getParameter<bool>("doLorentzFromAlignment");
useLAFromDB_ = conf.getParameter<bool>("useLAFromDB");

LogDebug("PixelCPEBase") << " LA constants - " << lAOffset_ << " " << lAWidthBPix_ << " " << lAWidthFPix_
<< endl; //dk
Expand Down Expand Up @@ -194,7 +195,8 @@ void PixelCPEBase::fillDetParams() {
p.bx = Bfield.x();

//--- Compute the Lorentz shifts for this detector element
if ((theFlag_ == 0) || DoLorentz_) { // do always for generic and if(DOLorentz) for templates
if ((theFlag_ == 0) || useLAFromDB_ ||
doLorentzFromAlignment_) { // do always for generic and if(DOLorentz) for templates
Copy link
Contributor

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

the comment about if(DOLorentz) does not seem to match the new implementation; please update

p.driftDirection = driftDirection(p, Bfield);
computeLorentzShifts(p);
}
Expand Down Expand Up @@ -470,4 +472,6 @@ void PixelCPEBase::fillPSetDescription(edm::ParameterSetDescription& desc) {
desc.add<double>("lAOffset", 0.0);
desc.add<double>("lAWidthBPix", 0.0);
desc.add<double>("lAWidthFPix", 0.0);
desc.add<bool>("doLorentzFromAlignment", false);
desc.add<bool>("useLAFromDB", true);
}
59 changes: 30 additions & 29 deletions RecoLocalTracker/SiPixelRecHits/src/PixelCPEClusterRepair.cc
Expand Up @@ -324,6 +324,11 @@ void PixelCPEClusterRepair::callTempReco1D(DetParam const& theDetParam,
// We have a boolean denoting whether the reco failed or not
theClusterParam.hasFilledProb_ = false;

// In case of template reco failure, these are the lorentz drift corrections
// to be applied
float lorentzshiftX = 0.5f * theDetParam.lorentzShiftInCmX;
float lorentzshiftY = 0.5f * theDetParam.lorentzShiftInCmY;

// ******************************************************************
//--- Call normal TemplateReco
//
Expand Down Expand Up @@ -363,30 +368,24 @@ void PixelCPEClusterRepair::callTempReco1D(DetParam const& theDetParam,

theClusterParam.probabilityX_ = theClusterParam.probabilityY_ = theClusterParam.probabilityQ_ = 0.f;
theClusterParam.qBin_ = 0;
//
// Template reco has failed, compute position estimates based on cluster center of gravity + Lorentz drift
// Future improvement would be to call generic reco instead

// Gavril: what do we do in this case ? For now, just return the cluster center of gravity in microns
// In the x case, apply a rough Lorentz drift average correction
// To do: call PixelCPEGeneric whenever PixelTempReco1D fails
float lorentz_drift = -999.9;
if (!GeomDetEnumerators::isEndcap(theDetParam.thePart))
lorentz_drift = 60.0f; // in microns
else
lorentz_drift = 10.0f; // in microns
// GG: trk angles needed to correct for bows/kinks
if (theClusterParam.with_track_angle) {
theClusterParam.templXrec_ =
theDetParam.theTopol->localX(theClusterParam.theCluster->x(), theClusterParam.loc_trk_pred) -
lorentz_drift * micronsToCm; // rough Lorentz drift correction
theDetParam.theTopol->localX(theClusterParam.theCluster->x(), theClusterParam.loc_trk_pred) + lorentzshiftX;
theClusterParam.templYrec_ =
theDetParam.theTopol->localY(theClusterParam.theCluster->y(), theClusterParam.loc_trk_pred);
theDetParam.theTopol->localY(theClusterParam.theCluster->y(), theClusterParam.loc_trk_pred) + lorentzshiftY;
} else {
edm::LogError("PixelCPEClusterRepair") << "@SUB = PixelCPEClusterRepair::localPosition"
<< "Should never be here. PixelCPEClusterRepair should always be called "
"with track angles. This is a bad error !!! ";

theClusterParam.templXrec_ = theDetParam.theTopol->localX(theClusterParam.theCluster->x()) -
lorentz_drift * micronsToCm; // rough Lorentz drift correction
theClusterParam.templYrec_ = theDetParam.theTopol->localY(theClusterParam.theCluster->y());
theClusterParam.templXrec_ =
theDetParam.theTopol->localX(theClusterParam.theCluster->x(), theClusterParam.loc_trk_pred) + lorentzshiftX;
theClusterParam.templYrec_ =
theDetParam.theTopol->localY(theClusterParam.theCluster->y(), theClusterParam.loc_trk_pred) + lorentzshiftY;
}
Copy link
Contributor

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

the implementation seems identical except for the error message

Suggested change
if (theClusterParam.with_track_angle) {
theClusterParam.templXrec_ =
theDetParam.theTopol->localX(theClusterParam.theCluster->x(), theClusterParam.loc_trk_pred) -
lorentz_drift * micronsToCm; // rough Lorentz drift correction
theDetParam.theTopol->localX(theClusterParam.theCluster->x(), theClusterParam.loc_trk_pred) + lorentzshiftX;
theClusterParam.templYrec_ =
theDetParam.theTopol->localY(theClusterParam.theCluster->y(), theClusterParam.loc_trk_pred);
theDetParam.theTopol->localY(theClusterParam.theCluster->y(), theClusterParam.loc_trk_pred) + lorentzshiftY;
} else {
edm::LogError("PixelCPEClusterRepair") << "@SUB = PixelCPEClusterRepair::localPosition"
<< "Should never be here. PixelCPEClusterRepair should always be called "
"with track angles. This is a bad error !!! ";
theClusterParam.templXrec_ = theDetParam.theTopol->localX(theClusterParam.theCluster->x()) -
lorentz_drift * micronsToCm; // rough Lorentz drift correction
theClusterParam.templYrec_ = theDetParam.theTopol->localY(theClusterParam.theCluster->y());
theClusterParam.templXrec_ =
theDetParam.theTopol->localX(theClusterParam.theCluster->x(), theClusterParam.loc_trk_pred) + lorentzshiftX;
theClusterParam.templYrec_ =
theDetParam.theTopol->localY(theClusterParam.theCluster->y(), theClusterParam.loc_trk_pred) + lorentzshiftY;
}
if (!theClusterParam.with_track_angle) {
edm::LogError("PixelCPEClusterRepair") << "@SUB = PixelCPEClusterRepair::localPosition"
<< "Should never be here. PixelCPEClusterRepair should always be called "
"with track angles. This is a bad error !!! ";
}
theClusterParam.templXrec_ = theDetParam.theTopol->localX(theClusterParam.theCluster->x(), theClusterParam.loc_trk_pred) + lorentzshiftX;
theClusterParam.templYrec_ = theDetParam.theTopol->localY(theClusterParam.theCluster->y(), theClusterParam.loc_trk_pred) + lorentzshiftY;

Copy link
Contributor Author

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Looking into it further, the loc_trk_pred should not be used when with_track_angle is false (never used in standard reconstruction) so the two cases are different and I updated the PR to match this.

} else {
//--- Template Reco succeeded. The probabilities are filled.
Expand Down Expand Up @@ -423,6 +422,11 @@ void PixelCPEClusterRepair::callTempReco2D(DetParam const& theDetParam,
// We have a boolean denoting whether the reco failed or not
theClusterParam.hasFilledProb_ = false;

// In case of template reco failure, these are the lorentz drift corrections
// to be applied
float lorentzshiftX = 0.5f * theDetParam.lorentzShiftInCmX;
float lorentzshiftY = 0.5f * theDetParam.lorentzShiftInCmY;

// ******************************************************************
//--- Call 2D TemplateReco
//
Expand Down Expand Up @@ -482,29 +486,26 @@ void PixelCPEClusterRepair::callTempReco2D(DetParam const& theDetParam,

theClusterParam.probabilityX_ = theClusterParam.probabilityY_ = theClusterParam.probabilityQ_ = 0.f;
theClusterParam.qBin_ = 0;
// GG: what do we do in this case? For now, just return the cluster center of gravity in microns
// In the x case, apply a rough Lorentz drift average correction
float lorentz_drift = -999.9;
if (!GeomDetEnumerators::isEndcap(theDetParam.thePart))
lorentz_drift = 60.0f; // in microns // &&& replace with a constant (globally)
else
lorentz_drift = 10.0f; // in microns
// GG: trk angles needed to correct for bows/kinks

// 2D Template reco has failed, compute position estimates based on cluster center of gravity + Lorentz drift
// Future improvement would be to call generic reco instead

if (theClusterParam.with_track_angle) {
theClusterParam.templXrec_ =
theDetParam.theTopol->localX(theClusterParam.theCluster->x(), theClusterParam.loc_trk_pred) -
lorentz_drift * micronsToCm; // rough Lorentz drift correction
theDetParam.theTopol->localX(theClusterParam.theCluster->x(), theClusterParam.loc_trk_pred) + lorentzshiftX;
theClusterParam.templYrec_ =
theDetParam.theTopol->localY(theClusterParam.theCluster->y(), theClusterParam.loc_trk_pred);
theDetParam.theTopol->localY(theClusterParam.theCluster->y(), theClusterParam.loc_trk_pred) + lorentzshiftY;
} else {
edm::LogError("PixelCPEClusterRepair") << "@SUB = PixelCPEClusterRepair::localPosition"
<< "Should never be here. PixelCPEClusterRepair should always be called "
"with track angles. This is a bad error !!! ";

theClusterParam.templXrec_ = theDetParam.theTopol->localX(theClusterParam.theCluster->x()) -
lorentz_drift * micronsToCm; // rough Lorentz drift correction
theClusterParam.templYrec_ = theDetParam.theTopol->localY(theClusterParam.theCluster->y());
theClusterParam.templXrec_ =
theDetParam.theTopol->localX(theClusterParam.theCluster->x(), theClusterParam.loc_trk_pred) + lorentzshiftX;
theClusterParam.templYrec_ =
theDetParam.theTopol->localY(theClusterParam.theCluster->y(), theClusterParam.loc_trk_pred) + lorentzshiftY;
Copy link
Contributor

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

it looks like the same code is repeated here as well as in the callTempReco1D

Copy link
Contributor Author

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Same comment/ change as above

}

} else {
//--- Template Reco succeeded.
theClusterParam.hasFilledProb_ = true;
Expand Down