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

fix:Incorrect MET correction method CorrectedMETFactory.py #1066

Open
wants to merge 4 commits into
base: master
Choose a base branch
from
Open
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
15 changes: 7 additions & 8 deletions src/coffea/jetmet_tools/CorrectedMETFactory.py
Original file line number Diff line number Diff line change
Expand Up @@ -4,11 +4,11 @@


def corrected_polar_met(
met_pt, met_phi, jet_pt, jet_phi, jet_pt_orig, positive=None, dx=None, dy=None
met_pt, met_phi, jet_pt, jet_phi, positive=None, dx=None, dy=None
):
sj, cj = numpy.sin(jet_phi), numpy.cos(jet_phi)
x = met_pt * numpy.cos(met_phi) + awkward.sum((jet_pt - jet_pt_orig) * cj, axis=1)
y = met_pt * numpy.sin(met_phi) + awkward.sum((jet_pt - jet_pt_orig) * sj, axis=1)
x = met_pt * numpy.cos(met_phi) - awkward.sum(jet_pt * cj, axis=1)
y = met_pt * numpy.sin(met_phi) - awkward.sum(jet_pt * sj, axis=1)
if positive is not None and dx is not None and dy is not None:
x = x + dx if positive else x - dx
y = y + dy if positive else y - dy
Expand Down Expand Up @@ -36,7 +36,7 @@ def __init__(self, name_map):

self.name_map = name_map

def build(self, in_MET, in_corrected_jets):
def build(self, in_MET, type1_MET, in_corrected_jets):
Copy link
Collaborator

Choose a reason for hiding this comment

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

Couldn't you just set the Type1 MET as the MET that is input to the build function? This seems odd.

Copy link
Author

Choose a reason for hiding this comment

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

  • In order to obtain the corresponding uncertainty values of MET after correction for JEC and JER, we must use RawMET as input for correction.RawMET means uncorrected PF MET, while MET refers to Type-1 corrected PF MET. If using MET as input, it means that it has been corrected twice.

  • But if we want to obtain the UnclusterdEnergy uncertainty value of MET, you must use MET as input, as in NanoAOD, only MET has the branch for MetUnclusterEnUpDeltaX/Y, while RawMET does not. Therefore, it is necessary to use MET to obtain the UnclusteredEnergy uncertainty value.

So I have to input both in_MET and type1_MET meanwhile, they represent RawMET and MET respectively, for example:
met = self._jmeu.corrected_met(event.RawMET,event.MET, jets_col, event.fixedGridRhoFastjetAll, event.caches[0])

Copy link
Collaborator

Choose a reason for hiding this comment

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

Perhaps a better way to put this is I think this can be solved by changing the flexibility in the configuration as opposed to changing the meaning of the physics objects that people get out from correctors.

Copy link
Author

Choose a reason for hiding this comment

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

But it is obvious that the CorrectedMETFactory function here cannot implement the correction of RawMET as input. Because RawMET does not have the MetUnclusterEnUpDeltaX/Y branch, a bug will occur. How can I modify the configuration file to solve this problem? At the same time, I also need to obtain the uncertainty of MET_UnclusterdEnergy_up/down

Copy link
Collaborator

Choose a reason for hiding this comment

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

This would be better solved by making a new met object with the appropriate substitutions.

raw_met = events.RawMET
met_to_correct = events.MET
met_to_correct["pt"] = raw_met.pt
met_to_correct["phi"] = raw_met.phi

# then go on to correct met_to_correct, adjusting configuration as appropriate

Copy link
Collaborator

Choose a reason for hiding this comment

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

Though I'm not sure if this spoils the validity of the systematic variations.

Copy link
Author

Choose a reason for hiding this comment

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

I know what you mean! Great!
This will not spoils the validity of the systemic variations, the uncertainty of MET_UnclusteredEnergy is calculated using the corrected MET_pt value through corrected_polar_met, rather than RawMET.pt.
You are incredibly smart, thank you so much for your patient and thoughtful answers.

if not isinstance(
in_MET, (awkward.highlevel.Array, dask_awkward.Array)
) or not isinstance(
Expand All @@ -60,7 +60,6 @@ def switch_properties(raw_met, corrected_jets, dx, dy, positive, save_orig):
raw_met[self.name_map["METphi"]],
corrected_jets[self.name_map["JetPt"]],
corrected_jets[self.name_map["JetPhi"]],
corrected_jets[self.name_map["ptRaw"]],
positive=positive,
dx=dx,
dy=dy,
Expand Down Expand Up @@ -144,10 +143,10 @@ def create_variants(raw_met, corrected_jets_or_variants, dx, dy):

out_dict["MET_UnclusteredEnergy"] = dask_awkward.map_partitions(
create_variants,
MET,
type1_MET,
corrected_jets,
MET[self.name_map["UnClusteredEnergyDeltaX"]],
MET[self.name_map["UnClusteredEnergyDeltaY"]],
type1_MET[self.name_map["UnClusteredEnergyDeltaX"]],
type1_MET[self.name_map["UnClusteredEnergyDeltaY"]],
label="UnclusteredEnergy_met",
)

Expand Down
4 changes: 2 additions & 2 deletions tests/test_jetmet_tools.py
Original file line number Diff line number Diff line change
Expand Up @@ -919,12 +919,12 @@ def smear_factor(jetPt, pt_gen, jersf):
toc = time.time()

print("setup corrected MET time =", toc - tic)

rawmet = events.RawMET
met = events.MET
tic = time.time()
# prof = pyinstrument.Profiler()
# prof.start()
corrected_met = met_factory.build(met, corrected_jets)
corrected_met = met_factory.build(rawmet, met, corrected_jets)
# prof.stop()
toc = time.time()

Expand Down