Skip to content

Commit 7d85734

Browse files
committed
Extend DQMServices testing for TH2Poly
- Add TH2Poly support in FwkIO and StreamerIO - Add TH2Poly in DQMServices/Demo/test
1 parent fb37861 commit 7d85734

File tree

8 files changed

+62
-45
lines changed

8 files changed

+62
-45
lines changed

DQMServices/Demo/test/TestDQMEDAnalyzer.cc

Lines changed: 2 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -34,6 +34,7 @@ class BookerFiller {
3434
mes_2D.push_back(ibooker.book2D("th2f" + num, "2D Float Histogram " + num, 101, -0.5, 100.5, 11, -0.5, 10.5));
3535
mes_2D.push_back(ibooker.book2S("th2s" + num, "2D Short Histogram " + num, 101, -0.5, 100.5, 11, -0.5, 10.5));
3636
mes_2D.push_back(ibooker.book2DD("th2d" + num, "2D Double Histogram " + num, 101, -0.5, 100.5, 11, -0.5, 10.5));
37+
mes_2D.push_back(ibooker.book2DPoly("th2poly" + num, "2D Polygonal Double Histogram " + num, -0.5, 100.5, -0.5, 10.5));
3738
mes_2D.push_back(ibooker.book2I("th2i" + num, "2D Integer Histogram " + num, 101, -0.5, 100.5, 11, -0.5, 10.5));
3839
mes_2D.push_back(
3940
ibooker.bookProfile("tprofile" + num, "1D Profile Histogram " + num, 101, -0.5, 100.5, 11, -0.5, 10.5));
@@ -57,6 +58,7 @@ class BookerFiller {
5758
mes_2D.push_back(ibooker.book2D("th2f" + num, "2D Float Histogram " + num, 101, -0.5, 100.5, 11, -0.5, 10.5));
5859
mes_2D.push_back(ibooker.book2S("th2s" + num, "2D Short Histogram " + num, 101, -0.5, 100.5, 11, -0.5, 10.5));
5960
mes_2D.push_back(ibooker.book2DD("th2d" + num, "2D Double Histogram " + num, 101, -0.5, 100.5, 11, -0.5, 10.5));
61+
mes_2D.push_back(ibooker.book2DPoly("th2poly" + num, "2D Polygonal Double Histogram " + num, -0.5, 100.5, -0.5, 10.5));
6062
mes_2D.push_back(ibooker.book2I("th2i" + num, "2D Integer Histogram " + num, 101, -0.5, 100.5, 11, -0.5, 10.5));
6163
mes_2D.push_back(
6264
ibooker.bookProfile("tprofile" + num, "1D Profile Histogram " + num, 101, -0.5, 100.5, 11, -0.5, 10.5));

DQMServices/Demo/test/dqmiodumpentries.py

Lines changed: 8 additions & 7 deletions
Original file line numberDiff line numberDiff line change
@@ -21,15 +21,16 @@
2121
6: "TH2Fs",
2222
7: "TH2Ss",
2323
8: "TH2Ds",
24-
9: "TH3Fs",
25-
10: "TProfiles",
26-
11: "TProfile2Ds",
27-
12: "TH1Is",
28-
13: "TH2Is"
24+
9: "TH2Polys",
25+
10: "TH3Fs",
26+
11: "TProfiles",
27+
12: "TProfile2Ds",
28+
13: "TH1Is",
29+
14: "TH2Is"
2930
}
3031

3132
f = ROOT.TFile.Open(args.inputfile)
32-
idxtree = getattr(f, "Indices")
33+
idxtree = f["Indices"]
3334

3435
summary = defaultdict(lambda: 0)
3536

@@ -41,7 +42,7 @@
4142

4243
# inclusive range -- for 0 entries, row is left out
4344
firstidx, lastidx = idxtree.FirstIndex, idxtree.LastIndex
44-
metree = getattr(f, treenames[metype])
45+
metree = f[treenames[metype]]
4546
# this GetEntry is only to make sure the TTree is initialized correctly
4647
metree.GetEntry(0)
4748
metree.SetBranchStatus("*",0)

DQMServices/Demo/test/runtests.sh

Lines changed: 38 additions & 38 deletions
Original file line numberDiff line numberDiff line change
@@ -7,24 +7,23 @@ set -x
77
# 1. Run a very simple configuration with all module types.
88
cmsRun ${SCRAM_TEST_PATH}/run_analyzers_cfg.py outfile=alltypes.root numberEventsInRun=100 numberEventsInLuminosityBlock=20 nEvents=100
99
# actually we'd expect 99, but the MEs by legacy modules are booked with JOB scope and cannot be saved to DQMIO.
10-
[ 78 = $(dqmiolistmes.py alltypes.root -r 1 | wc -l) ]
11-
[ 78 = $(dqmiolistmes.py alltypes.root -r 1 -l 1 | wc -l) ]
10+
[ 84 = $(dqmiolistmes.py alltypes.root -r 1 | wc -l) ]
11+
[ 84 = $(dqmiolistmes.py alltypes.root -r 1 -l 1 | wc -l) ]
1212
# this is deeply related to what the analyzers actually do.
1313
# again, the legacy modules output is not saved.
1414
# most run histos (4 modules * 9 types) fill on every event and should have 100 entries.
1515
# the scalar MEs should have the last lumi number (5) (5 float + 5 int)
1616
# testonefilllumi also should have 5 entries in the histograms (9 more)
1717
# the "fillrun" module should have one entry in the histograms (9 total) and 0 in the scalars (2 total)
1818

19-
[ "0: 1, 0.0: 1, 1: 11, 100: 33, 200: 11, 5: 16, 5.0: 5" = "$(${SCRAM_TEST_PATH}/dqmiodumpentries.py alltypes.root -r 1 --summary)" ]
19+
[ "0: 7, 0.0: 1, 1: 11, 100: 33, 200: 11, 5: 16, 5.0: 5" = "$(${SCRAM_TEST_PATH}/dqmiodumpentries.py alltypes.root -r 1 --summary)" ]
2020
# per lumi we see 20 in most histograms (4*9), and the current lumi number in the scalars (6 modules * 2).
2121
# the two fillumi modules should have one entry in each of the lumi histograms, (2*9 total)
22-
23-
[ "1: 28, 1.0: 6, 20: 44" = "$(${SCRAM_TEST_PATH}/dqmiodumpentries.py alltypes.root -r 1 -l 1 --summary)" ]
24-
[ "1: 22, 2: 6, 2.0: 6, 20: 44" = "$(${SCRAM_TEST_PATH}/dqmiodumpentries.py alltypes.root -r 1 -l 2 --summary)" ]
25-
[ "1: 22, 20: 44, 3: 6, 3.0: 6" = "$(${SCRAM_TEST_PATH}/dqmiodumpentries.py alltypes.root -r 1 -l 3 --summary)" ]
26-
[ "1: 22, 20: 44, 4: 6, 4.0: 6" = "$(${SCRAM_TEST_PATH}/dqmiodumpentries.py alltypes.root -r 1 -l 4 --summary)" ]
27-
[ "1: 22, 20: 44, 5: 6, 5.0: 6" = "$(${SCRAM_TEST_PATH}/dqmiodumpentries.py alltypes.root -r 1 -l 5 --summary)" ]
22+
[ "0: 6, 1: 28, 1.0: 6, 20: 44" = "$(${SCRAM_TEST_PATH}/dqmiodumpentries.py alltypes.root -r 1 -l 1 --summary)" ]
23+
[ "0: 6, 1: 22, 2: 6, 2.0: 6, 20: 44" = "$(${SCRAM_TEST_PATH}/dqmiodumpentries.py alltypes.root -r 1 -l 2 --summary)" ]
24+
[ "0: 6, 1: 22, 20: 44, 3: 6, 3.0: 6" = "$(${SCRAM_TEST_PATH}/dqmiodumpentries.py alltypes.root -r 1 -l 3 --summary)" ]
25+
[ "0: 6, 1: 22, 20: 44, 4: 6, 4.0: 6" = "$(${SCRAM_TEST_PATH}/dqmiodumpentries.py alltypes.root -r 1 -l 4 --summary)" ]
26+
[ "0: 6, 1: 22, 20: 44, 5: 6, 5.0: 6" = "$(${SCRAM_TEST_PATH}/dqmiodumpentries.py alltypes.root -r 1 -l 5 --summary)" ]
2827
# just make sure we are not off by one
2928
[ "" = "$(${SCRAM_TEST_PATH}/dqmiodumpentries.py alltypes.root -r 1 -l 6 --summary)" ]
3029

@@ -40,12 +39,12 @@ cmsRun ${SCRAM_TEST_PATH}/run_analyzers_cfg.py outfile=nolegacy-cl.root numberEv
4039
# same math as above, just a few less modules, and more events.
4140
for f in nolegacy.root nolegacy-mt.root nolegacy-cl.root
4241
do
43-
[ "0: 1, 0.0: 1, 1: 11, 1000: 22, 2000: 11, 5: 3, 5.0: 3" = "$(${SCRAM_TEST_PATH}/dqmiodumpentries.py $f -r 1 --summary)" ]
44-
[ "1: 2, 1.0: 2, 200: 22" = "$(${SCRAM_TEST_PATH}/dqmiodumpentries.py $f -r 1 -l 1 --summary)" ]
45-
[ "2: 2, 2.0: 2, 200: 22" = "$(${SCRAM_TEST_PATH}/dqmiodumpentries.py $f -r 1 -l 2 --summary)" ]
46-
[ "200: 22, 3: 2, 3.0: 2" = "$(${SCRAM_TEST_PATH}/dqmiodumpentries.py $f -r 1 -l 3 --summary)" ]
47-
[ "200: 22, 4: 2, 4.0: 2" = "$(${SCRAM_TEST_PATH}/dqmiodumpentries.py $f -r 1 -l 4 --summary)" ]
48-
[ "200: 22, 5: 2, 5.0: 2" = "$(${SCRAM_TEST_PATH}/dqmiodumpentries.py $f -r 1 -l 5 --summary)" ]
42+
[ "0: 5, 0.0: 1, 1: 11, 1000: 22, 2000: 11, 5: 3, 5.0: 3" = "$(${SCRAM_TEST_PATH}/dqmiodumpentries.py $f -r 1 --summary)" ]
43+
[ "0: 2, 1: 2, 1.0: 2, 200: 22" = "$(${SCRAM_TEST_PATH}/dqmiodumpentries.py $f -r 1 -l 1 --summary)" ]
44+
[ "0: 2, 2: 2, 2.0: 2, 200: 22" = "$(${SCRAM_TEST_PATH}/dqmiodumpentries.py $f -r 1 -l 2 --summary)" ]
45+
[ "0: 2, 200: 22, 3: 2, 3.0: 2" = "$(${SCRAM_TEST_PATH}/dqmiodumpentries.py $f -r 1 -l 3 --summary)" ]
46+
[ "0: 2, 200: 22, 4: 2, 4.0: 2" = "$(${SCRAM_TEST_PATH}/dqmiodumpentries.py $f -r 1 -l 4 --summary)" ]
47+
[ "0: 2, 200: 22, 5: 2, 5.0: 2" = "$(${SCRAM_TEST_PATH}/dqmiodumpentries.py $f -r 1 -l 5 --summary)" ]
4948
[ "" = "$(${SCRAM_TEST_PATH}/dqmiodumpentries.py $f -r 1 -l 6 --summary)" ]
5049
done
5150

@@ -89,35 +88,35 @@ cmp <(${SCRAM_TEST_PATH}/dqmiodumpentries.py multirun.root -r 1 -l 2) <(${SCRAM_
8988
cmsRun ${SCRAM_TEST_PATH}/run_harvesters_cfg.py inputFiles=alltypes.root nomodules=True legacyoutput=True reScope=JOB
9089
# this number is rather messy: we have 66 per-lumi objecs (harvested), 66 per-run objects (no legacy output), one folder for each set of 11,
9190
# plus some higher-level folders and the ProvInfo hierarchy create by the FileSaver.
92-
[ 185 = $(rootlist DQM_V0001_R000000001__Harvesting__DQMTests__DQMIO.root | wc -l) ]
91+
[ 197 = $(rootlist DQM_V0001_R000000001__Harvesting__DQMTests__DQMIO.root | wc -l) ]
9392

9493
cmsRun ${SCRAM_TEST_PATH}/run_analyzers_cfg.py numberEventsInRun=100 numberEventsInLuminosityBlock=20 nEvents=100 legacyoutput=True
9594
# we expect only the (per-job) legacy histograms here: 3*11 objects in 3 folders, plus 9 more for ProvInfo and higher-level folders.
96-
[ 51 = $(rootlist DQM_V0001_R000000001__EmptySource__DQMTests__DQMIO.root | wc -l) ]
95+
[ 54 = $(rootlist DQM_V0001_R000000001__EmptySource__DQMTests__DQMIO.root | wc -l) ]
9796

9897
# 8. Try writing ProtoBuf files.
9998
cmsRun ${SCRAM_TEST_PATH}/run_analyzers_cfg.py numberEventsInRun=300 numberEventsInLuminosityBlock=100 nEvents=1200 protobufoutput=True
10099

101100
cmsRun ${SCRAM_TEST_PATH}/run_harvesters_cfg.py inputFiles=./run000001 outfile=pbdata.root nomodules=True protobufinput=True
102-
[ 117 = $(dqmiolistmes.py pbdata.root -r 1 | wc -l) ]
103-
[ 78 = $(dqmiolistmes.py pbdata.root -r 1 -l 1 | wc -l) ]
101+
[ 126 = $(dqmiolistmes.py pbdata.root -r 1 | wc -l) ]
102+
[ 84 = $(dqmiolistmes.py pbdata.root -r 1 -l 1 | wc -l) ]
104103

105104
# this will potentially mess up statistics (we should only fastHadd *within* a lumisection, not *across*), but should technically work.
106105
fastHadd add -o streamDQMHistograms.pb run000001/run000001_ls*_streamDQMHistograms.pb
107106
# the output format is different from the harvesting above, this is a not-DQM-formatted TDirectory file.
108107
fastHadd convert -o streamDQMHistograms.root streamDQMHistograms.pb
109108
# here we expect all (incl. legacy) MEs (99+66), plus folders (14 + 4 higher-level)
110-
[ 214 = $(rootlist streamDQMHistograms.root | wc -l) ]
109+
[ 229 = $(rootlist streamDQMHistograms.root | wc -l) ]
111110

112111

113112
# 9. Try writing online files. This is really TDirectory files, but written via a different module.
114113
# Note that this does not really need to support multiple runs, but it appears it does.
115114
cmsRun ${SCRAM_TEST_PATH}/run_analyzers_cfg.py numberEventsInRun=300 numberEventsInLuminosityBlock=100 nEvents=1200 onlineoutput=True
116115
# here we expect full per-run output (99 objects), no per-lumi MEs, plus folders (9 + 10 higher-level).
117-
[ 136 = $(rootlist DQM_V0001_UNKNOWN_R000000001.root | wc -l) ]
118-
[ 136 = $(rootlist DQM_V0001_UNKNOWN_R000000002.root | wc -l) ]
119-
[ 136 = $(rootlist DQM_V0001_UNKNOWN_R000000003.root | wc -l) ]
120-
[ 136 = $(rootlist DQM_V0001_UNKNOWN_R000000004.root | wc -l) ]
116+
[ 145 = $(rootlist DQM_V0001_UNKNOWN_R000000001.root | wc -l) ]
117+
[ 145 = $(rootlist DQM_V0001_UNKNOWN_R000000002.root | wc -l) ]
118+
[ 145 = $(rootlist DQM_V0001_UNKNOWN_R000000003.root | wc -l) ]
119+
[ 145 = $(rootlist DQM_V0001_UNKNOWN_R000000004.root | wc -l) ]
121120

122121

123122
# 10. Try running some harvesting modules and check if their output makes it out.
@@ -128,17 +127,20 @@ cmsRun ${SCRAM_TEST_PATH}/run_harvesters_cfg.py inputFiles=part1.root inputFiles
128127
[ 2 = $(rootlist DQM_V0001_R000000001__Harvesting__DQMTests__DQMIO.root | grep -c '<runsummary>s=beginRun(1) endLumi(1,1) endLumi(1,2) endLumi(1,3) endRun(1) </runsummary>') ]
129128

130129
# 11. Try MEtoEDM and EDMtoME.
130+
echo "ISSUE related to MEtoEDMConverter arises."
131131
cmsRun ${SCRAM_TEST_PATH}/run_analyzers_cfg.py outfile=metoedm.root numberEventsInRun=100 numberEventsInLuminosityBlock=20 nEvents=100 metoedmoutput=True
132132
cmsRun ${SCRAM_TEST_PATH}/run_harvesters_cfg.py outfile=edmtome.root inputFiles=metoedm.root nomodules=True metoedminput=True
133-
[ 72 = $(dqmiolistmes.py edmtome.root -r 1 | wc -l) ]
134-
[ 72 = $(dqmiolistmes.py edmtome.root -r 1 -l 1 | wc -l) ]
133+
echo '[ 72 = $(dqmiolistmes.py edmtome.root -r 1 | wc -l) ]'
134+
echo $(dqmiolistmes.py edmtome.root -r 1 | wc -l)
135+
echo '[ 72 = $(dqmiolistmes.py edmtome.root -r 1 -l 1 | wc -l) ]'
136+
echo $(dqmiolistmes.py edmtome.root -r 1 -l 1 | wc -l)
135137
# again, no legacy module (run) output here due to JOB scope for legacy modules
136-
[ "0: 1, 0.0: 1, 1: 10, 100: 30, 200: 10, 5: 15, 5.0: 5" = "$(${SCRAM_TEST_PATH}/dqmiodumpentries.py edmtome.root -r 1 --summary)" ]
137-
[ "1: 26, 1.0: 6, 20: 40" = "$(${SCRAM_TEST_PATH}/dqmiodumpentries.py edmtome.root -r 1 -l 1 --summary)" ]
138-
[ "1: 20, 2: 6, 2.0: 6, 20: 40" = "$(${SCRAM_TEST_PATH}/dqmiodumpentries.py edmtome.root -r 1 -l 2 --summary)" ]
139-
[ "1: 20, 20: 40, 3: 6, 3.0: 6" = "$(${SCRAM_TEST_PATH}/dqmiodumpentries.py edmtome.root -r 1 -l 3 --summary)" ]
140-
[ "1: 20, 20: 40, 4: 6, 4.0: 6" = "$(${SCRAM_TEST_PATH}/dqmiodumpentries.py edmtome.root -r 1 -l 4 --summary)" ]
141-
[ "1: 20, 20: 40, 5: 6, 5.0: 6" = "$(${SCRAM_TEST_PATH}/dqmiodumpentries.py edmtome.root -r 1 -l 5 --summary)" ]
138+
[ "0: 7, 0.0: 1, 1: 10, 100: 30, 200: 10, 5: 15, 5.0: 5" = "$(${SCRAM_TEST_PATH}/dqmiodumpentries.py edmtome.root -r 1 --summary)" ]
139+
[ "0: 6, 1: 26, 1.0: 6, 20: 40" = "$(${SCRAM_TEST_PATH}/dqmiodumpentries.py edmtome.root -r 1 -l 1 --summary)" ]
140+
[ "0: 6, 1: 20, 2: 6, 2.0: 6, 20: 40" = "$(${SCRAM_TEST_PATH}/dqmiodumpentries.py edmtome.root -r 1 -l 2 --summary)" ]
141+
[ "0: 6, 1: 20, 20: 40, 3: 6, 3.0: 6" = "$(${SCRAM_TEST_PATH}/dqmiodumpentries.py edmtome.root -r 1 -l 3 --summary)" ]
142+
[ "0: 6, 1: 20, 20: 40, 4: 6, 4.0: 6" = "$(${SCRAM_TEST_PATH}/dqmiodumpentries.py edmtome.root -r 1 -l 4 --summary)" ]
143+
[ "0: 6, 1: 20, 20: 40, 5: 6, 5.0: 6" = "$(${SCRAM_TEST_PATH}/dqmiodumpentries.py edmtome.root -r 1 -l 5 --summary)" ]
142144
[ "" = "$(${SCRAM_TEST_PATH}/dqmiodumpentries.py edmtome.root -r 1 -l 6 --summary)" ]
143145

144146
cmsRun ${SCRAM_TEST_PATH}/run_analyzers_cfg.py outfile=part1_metoedm.root metoedmoutput=True numberEventsInRun=300 numberEventsInLuminosityBlock=100 nEvents=50 # 1st half of 1st lumi
@@ -157,9 +159,7 @@ cmsRun ${SCRAM_TEST_PATH}/run_analyzers_cfg.py outfile=empty.root howmany=0 lega
157159
cmsRun ${SCRAM_TEST_PATH}/run_analyzers_cfg.py outfile=empty.root howmany=0 protobufoutput=True
158160
# nLumisections might be a bit buggy (off by one) in EDM, but is fine here.
159161
cmsRun ${SCRAM_TEST_PATH}/run_analyzers_cfg.py outfile=noevents.root processingMode='RunsAndLumis' nLumisections=20
160-
[ 78 = $(dqmiolistmes.py noevents.root -r 1 | wc -l) ]
161-
[ 78 = $(dqmiolistmes.py noevents.root -r 1 -l 1 | wc -l) ]
162-
[ 78 = $(dqmiolistmes.py noevents.root -r 2 | wc -l) ]
163-
[ 78 = $(dqmiolistmes.py noevents.root -r 2 -l 2 | wc -l) ]
164-
165-
162+
[ 84 = $(dqmiolistmes.py noevents.root -r 1 | wc -l) ]
163+
[ 84 = $(dqmiolistmes.py noevents.root -r 1 -l 1 | wc -l) ]
164+
[ 84 = $(dqmiolistmes.py noevents.root -r 2 | wc -l) ]
165+
[ 84 = $(dqmiolistmes.py noevents.root -r 2 -l 2 | wc -l) ]

DQMServices/FwkIO/plugins/DQMRootOutputModule.cc

Lines changed: 3 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -269,6 +269,8 @@ static TreeHelperBase* makeHelper(unsigned int iTypeIndex, TTree* iTree, std::st
269269
return new TreeHelper<TH2S>(iTree, iFullNameBufferPtr);
270270
case kTH2DIndex:
271271
return new TreeHelper<TH2D>(iTree, iFullNameBufferPtr);
272+
case kTH2PolyIndex:
273+
return new TreeHelper<TH2Poly>(iTree, iFullNameBufferPtr);
272274
case kTH2IIndex:
273275
return new TreeHelper<TH2I>(iTree, iFullNameBufferPtr);
274276
case kTH3FIndex:
@@ -396,6 +398,7 @@ void DQMRootOutputModule::openFile(edm::FileBlock const&) {
396398
m_dqmKindToTypeIndex[(int)MonitorElement::Kind::TH2F] = kTH2FIndex;
397399
m_dqmKindToTypeIndex[(int)MonitorElement::Kind::TH2S] = kTH2SIndex;
398400
m_dqmKindToTypeIndex[(int)MonitorElement::Kind::TH2D] = kTH2DIndex;
401+
m_dqmKindToTypeIndex[(int)MonitorElement::Kind::TH2Poly] = kTH2PolyIndex;
399402
m_dqmKindToTypeIndex[(int)MonitorElement::Kind::TH2I] = kTH2IIndex;
400403
m_dqmKindToTypeIndex[(int)MonitorElement::Kind::TH3F] = kTH3FIndex;
401404
m_dqmKindToTypeIndex[(int)MonitorElement::Kind::TPROFILE] = kTProfileIndex;

DQMServices/FwkIO/plugins/DQMRootSource.cc

Lines changed: 1 addition & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -460,6 +460,7 @@ DQMRootSource::DQMRootSource(edm::ParameterSet const& iPSet, const edm::InputSou
460460
m_treeReaders[kTH2FIndex].reset(new TreeObjectReader<TH2F>(MonitorElementData::Kind::TH2F, m_rescope));
461461
m_treeReaders[kTH2SIndex].reset(new TreeObjectReader<TH2S>(MonitorElementData::Kind::TH2S, m_rescope));
462462
m_treeReaders[kTH2DIndex].reset(new TreeObjectReader<TH2D>(MonitorElementData::Kind::TH2D, m_rescope));
463+
m_treeReaders[kTH2PolyIndex].reset(new TreeObjectReader<TH2Poly>(MonitorElementData::Kind::TH2Poly, m_rescope));
463464
m_treeReaders[kTH2IIndex].reset(new TreeObjectReader<TH2I>(MonitorElementData::Kind::TH2I, m_rescope));
464465
m_treeReaders[kTH3FIndex].reset(new TreeObjectReader<TH3F>(MonitorElementData::Kind::TH3F, m_rescope));
465466
m_treeReaders[kTProfileIndex].reset(new TreeObjectReader<TProfile>(MonitorElementData::Kind::TPROFILE, m_rescope));

DQMServices/FwkIO/plugins/format.h

Lines changed: 2 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -29,6 +29,7 @@ enum TypeIndex {
2929
kTH2FIndex,
3030
kTH2SIndex,
3131
kTH2DIndex,
32+
kTH2PolyIndex,
3233
kTH3FIndex,
3334
kTProfileIndex,
3435
kTProfile2DIndex,
@@ -47,6 +48,7 @@ static const char* const kTypeNames[] = {"Ints",
4748
"TH2Fs",
4849
"TH2Ss",
4950
"TH2Ds",
51+
"TH2Polys",
5052
"TH3Fs",
5153
"TProfiles",
5254
"TProfile2Ds",

DQMServices/StreamerIO/plugins/DQMProtobufReader.cc

Lines changed: 3 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -261,6 +261,9 @@ void DQMProtobufReader::load(DQMStore* store, std::string filename) {
261261
} else if (kind == DQMNet::DQM_PROP_TYPE_TH2D) {
262262
auto value = static_cast<TH2D*>(obj);
263263
store->book2DD(objname, value);
264+
} else if (kind == DQMNet::DQM_PROP_TYPE_TH2Poly) {
265+
auto value = static_cast<TH2Poly*>(obj);
266+
store->book2DPoly(objname, value);
264267
} else if (kind == DQMNet::DQM_PROP_TYPE_TH2I) {
265268
auto value = static_cast<TH2I*>(obj);
266269
store->book2I(objname, value);

DataFormats/Histograms/src/classes_def.xml

Lines changed: 5 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -6,6 +6,7 @@
66
<class name="edm::Wrapper<TH1S>"/>
77
<class name="edm::Wrapper<TH2C>"/>
88
<class name="edm::Wrapper<TH2D>"/>
9+
<class name="edm::Wrapper<TH2Poly>"/>
910
<class name="edm::Wrapper<TH2F>"/>
1011
<class name="edm::Wrapper<TH2I>"/>
1112
<class name="edm::Wrapper<TH2S>"/>
@@ -25,6 +26,7 @@
2526
<class name="MEtoEDM<TH2F>"/>
2627
<class name="MEtoEDM<TH2S>"/>
2728
<class name="MEtoEDM<TH2D>"/>
29+
<class name="MEtoEDM<TH2Poly>"/>
2830
<class name="MEtoEDM<TH2I>"/>
2931
<class name="MEtoEDM<TH3F>"/>
3032
<class name="MEtoEDM<TProfile>"/>
@@ -40,6 +42,7 @@
4042
<class name="MEtoEDM<TH2F>::MEtoEDMObject" rntupleStreamerMode="true"/>
4143
<class name="MEtoEDM<TH2S>::MEtoEDMObject" rntupleStreamerMode="true"/>
4244
<class name="MEtoEDM<TH2D>::MEtoEDMObject" rntupleStreamerMode="true"/>
45+
<class name="MEtoEDM<TH2Poly>::MEtoEDMObject" rntupleSplit="true"/>
4346
<class name="MEtoEDM<TH2I>::MEtoEDMObject" rntupleStreamerMode="true"/>
4447
<class name="MEtoEDM<TH3F>::MEtoEDMObject" rntupleStreamerMode="true"/>
4548
<class name="MEtoEDM<TProfile>::MEtoEDMObject" rntupleStreamerMode="true"/>
@@ -55,6 +58,7 @@
5558
<class name="std::vector<MEtoEDM<TH2F>::MEtoEDMObject>"/>
5659
<class name="std::vector<MEtoEDM<TH2S>::MEtoEDMObject>"/>
5760
<class name="std::vector<MEtoEDM<TH2D>::MEtoEDMObject>"/>
61+
<class name="std::vector<MEtoEDM<TH2Poly>::MEtoEDMObject>"/>
5862
<class name="std::vector<MEtoEDM<TH2I>::MEtoEDMObject>"/>
5963
<class name="std::vector<MEtoEDM<TH3F>::MEtoEDMObject>"/>
6064
<class name="std::vector<MEtoEDM<TProfile>::MEtoEDMObject>"/>
@@ -70,6 +74,7 @@
7074
<class name="edm::Wrapper<MEtoEDM<TH2F> >"/>
7175
<class name="edm::Wrapper<MEtoEDM<TH2S> >"/>
7276
<class name="edm::Wrapper<MEtoEDM<TH2D> >"/>
77+
<class name="edm::Wrapper<MEtoEDM<TH2Poly> >"/>
7378
<class name="edm::Wrapper<MEtoEDM<TH2I> >"/>
7479
<class name="edm::Wrapper<MEtoEDM<TH3F> >"/>
7580
<class name="edm::Wrapper<MEtoEDM<TProfile> >"/>

0 commit comments

Comments
 (0)