Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:01:49

0001 /*!
0002   \file SiStripLorentzAngle_PayloadInspector
0003   \Payload Inspector Plugin for SiStrip Lorentz angles
0004   \author M. Musich
0005   \version $Revision: 1.0 $
0006   \date $Date: 2017/09/21 10:59:56 $
0007 */
0008 
0009 #include "CalibTracker/StandaloneTrackerTopology/interface/StandaloneTrackerTopology.h"
0010 #include "CommonTools/TrackerMap/interface/TrackerMap.h"
0011 #include "CondCore/CondDB/interface/Time.h"
0012 #include "CondCore/SiStripPlugins/interface/SiStripPayloadInspectorHelper.h"
0013 #include "CondCore/Utilities/interface/PayloadInspector.h"
0014 #include "CondCore/Utilities/interface/PayloadInspectorModule.h"
0015 #include "CondFormats/SiStripObjects/interface/SiStripDetSummary.h"
0016 #include "CondFormats/SiStripObjects/interface/SiStripLorentzAngle.h"
0017 #include "DQM/TrackerRemapper/interface/SiStripTkMaps.h"
0018 #include "SiStripCondObjectRepresent.h"
0019 
0020 #include <memory>
0021 #include <sstream>
0022 
0023 // include ROOT
0024 #include "TH2F.h"
0025 #include "TLegend.h"
0026 #include "TCanvas.h"
0027 #include "TLine.h"
0028 #include "TStyle.h"
0029 #include "TLatex.h"
0030 #include "TPave.h"
0031 #include "TPaveStats.h"
0032 #include "TGaxis.h"
0033 
0034 namespace {
0035 
0036   using namespace cond::payloadInspector;
0037 
0038   class SiStripLorentzAngleContainer
0039       : public SiStripCondObjectRepresent::SiStripDataContainer<SiStripLorentzAngle, float> {
0040   public:
0041     SiStripLorentzAngleContainer(const std::shared_ptr<SiStripLorentzAngle> &payload,
0042                                  const SiStripPI::MetaData &metadata,
0043                                  const std::string &tagName)
0044         : SiStripCondObjectRepresent::SiStripDataContainer<SiStripLorentzAngle, float>(payload, metadata, tagName) {
0045       payloadType_ = "SiStripLorentzAngle";
0046       setGranularity(SiStripCondObjectRepresent::PERMODULE);
0047     }
0048 
0049     void storeAllValues() override {
0050       auto LAMap_ = payload_->getLorentzAngles();
0051       for (const auto &element : LAMap_) {
0052         SiStripCondData_.fillByPushBack(element.first, element.second);
0053       }
0054     }
0055   };
0056 
0057   /************************************************
0058     testing the machinery
0059   ************************************************/
0060   class SiStripLorentzAngleTest : public PlotImage<SiStripLorentzAngle, SINGLE_IOV> {
0061   public:
0062     SiStripLorentzAngleTest() : PlotImage<SiStripLorentzAngle, SINGLE_IOV>("SiStrip LorentzAngle values") {}
0063 
0064     bool fill() override {
0065       auto tag = PlotBase::getTag<0>();
0066       auto iov = tag.iovs.front();
0067       auto tagname = tag.name;
0068       std::shared_ptr<SiStripLorentzAngle> payload = fetchPayload(std::get<1>(iov));
0069       if (payload.get()) {
0070         SiStripLorentzAngleContainer *objContainer = new SiStripLorentzAngleContainer(payload, iov, tagname);
0071         //objContainer->printAll();
0072 
0073         TCanvas canvas("Partion summary", "partition summary", 1200, 1000);
0074         objContainer->fillSummary(canvas);
0075 
0076         std::string fileName(m_imageFileName);
0077         canvas.SaveAs(fileName.c_str());
0078 
0079       }  // payload
0080       return true;
0081     }  // fill
0082   };
0083 
0084   class SiStripLorentzAngleByPartition : public PlotImage<SiStripLorentzAngle, SINGLE_IOV> {
0085   public:
0086     SiStripLorentzAngleByPartition()
0087         : PlotImage<SiStripLorentzAngle, SINGLE_IOV>("SiStrip LorentzAngle By Partition") {}
0088 
0089     bool fill() override {
0090       auto tag = PlotBase::getTag<0>();
0091       auto iov = tag.iovs.front();
0092       auto tagname = tag.name;
0093       std::shared_ptr<SiStripLorentzAngle> payload = fetchPayload(std::get<1>(iov));
0094       if (payload.get()) {
0095         SiStripLorentzAngleContainer *objContainer = new SiStripLorentzAngleContainer(payload, iov, tagname);
0096         objContainer->printAll();
0097 
0098         TCanvas canvas("Partition summary", "partition summary", 1400, 1000);
0099         objContainer->fillByPartition(canvas, 100, 0., 0.05);
0100 
0101         std::string fileName(m_imageFileName);
0102         canvas.SaveAs(fileName.c_str());
0103       }  // payload
0104       return true;
0105     }  // fill
0106   };
0107 
0108   class SiStripLorentzAngleCompareByRegion : public PlotImage<SiStripLorentzAngle, MULTI_IOV, 2> {
0109   public:
0110     SiStripLorentzAngleCompareByRegion()
0111         : PlotImage<SiStripLorentzAngle, MULTI_IOV, 2>("SiStrip LorentzAngle By Partition") {}
0112 
0113     bool fill() override {
0114       // trick to deal with the multi-ioved tag and two tag case at the same time
0115       auto theIOVs = PlotBase::getTag<0>().iovs;
0116       auto tagname1 = PlotBase::getTag<0>().name;
0117       auto tag2iovs = PlotBase::getTag<1>().iovs;
0118       auto tagname2 = PlotBase::getTag<1>().name;
0119       SiStripPI::MetaData firstiov = theIOVs.front();
0120       SiStripPI::MetaData lastiov = tag2iovs.front();
0121 
0122       std::shared_ptr<SiStripLorentzAngle> last_payload = fetchPayload(std::get<1>(lastiov));
0123       std::shared_ptr<SiStripLorentzAngle> first_payload = fetchPayload(std::get<1>(firstiov));
0124 
0125       SiStripLorentzAngleContainer *l_objContainer = new SiStripLorentzAngleContainer(last_payload, lastiov, tagname1);
0126       SiStripLorentzAngleContainer *f_objContainer =
0127           new SiStripLorentzAngleContainer(first_payload, firstiov, tagname2);
0128 
0129       l_objContainer->compare(f_objContainer);
0130 
0131       //l_objContainer->printAll();
0132 
0133       TCanvas canvas("Partition summary", "partition summary", 1400, 1000);
0134       l_objContainer->fillSummary(canvas);
0135 
0136       std::string fileName(m_imageFileName);
0137       canvas.SaveAs(fileName.c_str());
0138 
0139       return true;
0140     }  // fill
0141   };
0142 
0143   /************************************************
0144     1d histogram of SiStripLorentzAngle of 1 IOV 
0145   *************************************************/
0146 
0147   // inherit from one of the predefined plot class: Histogram1D
0148   class SiStripLorentzAngleValue : public Histogram1D<SiStripLorentzAngle, SINGLE_IOV> {
0149   public:
0150     SiStripLorentzAngleValue()
0151         : Histogram1D<SiStripLorentzAngle, SINGLE_IOV>(
0152               "SiStrip LorentzAngle values", "SiStrip LorentzAngle values", 100, 0.0, 0.05) {}
0153 
0154     bool fill() override {
0155       auto tag = PlotBase::getTag<0>();
0156       for (auto const &iov : tag.iovs) {
0157         std::shared_ptr<SiStripLorentzAngle> payload = Base::fetchPayload(std::get<1>(iov));
0158         if (payload.get()) {
0159           std::map<uint32_t, float> LAMap_ = payload->getLorentzAngles();
0160 
0161           for (const auto &element : LAMap_) {
0162             fillWithValue(element.second);
0163           }
0164         }  // payload
0165       }    // iovs
0166       return true;
0167     }  // fill
0168   };
0169 
0170   /************************************************
0171     TrackerMap of SiStrip Lorentz Angle
0172   *************************************************/
0173   class SiStripLorentzAngle_TrackerMap : public PlotImage<SiStripLorentzAngle, SINGLE_IOV> {
0174   public:
0175     SiStripLorentzAngle_TrackerMap()
0176         : PlotImage<SiStripLorentzAngle, SINGLE_IOV>("Tracker Map SiStrip Lorentz Angle") {}
0177 
0178     bool fill() override {
0179       auto tag = PlotBase::getTag<0>();
0180       auto iov = tag.iovs.front();
0181       std::shared_ptr<SiStripLorentzAngle> payload = fetchPayload(std::get<1>(iov));
0182 
0183       std::unique_ptr<TrackerMap> tmap = std::make_unique<TrackerMap>("SiStripLorentzAngle");
0184       tmap->setPalette(1);
0185       std::string titleMap = "TrackerMap of SiStrip Lorentz Angle per module, payload : " + std::get<1>(iov);
0186       tmap->setTitle(titleMap);
0187 
0188       std::map<uint32_t, float> LAMap_ = payload->getLorentzAngles();
0189 
0190       for (const auto &element : LAMap_) {
0191         tmap->fill(element.first, element.second);
0192       }  // loop over the LA MAP
0193 
0194       std::pair<float, float> extrema = tmap->getAutomaticRange();
0195 
0196       std::string fileName(m_imageFileName);
0197 
0198       // protect against uniform values (LA values are defined positive)
0199       if (extrema.first != extrema.second) {
0200         tmap->save(true, 0, 0, fileName);
0201       } else {
0202         tmap->save(true, extrema.first * 0.95, extrema.first * 1.05, fileName);
0203       }
0204 
0205       return true;
0206     }
0207   };
0208 
0209   /************************************************
0210     SiStripTkMaps of SiStrip Lorentz Angle
0211   *************************************************/
0212   class SiStripLorentzAngleTH2PolyTkMap : public PlotImage<SiStripLorentzAngle, SINGLE_IOV> {
0213   public:
0214     SiStripLorentzAngleTH2PolyTkMap()
0215         : PlotImage<SiStripLorentzAngle, SINGLE_IOV>("Tracker Map SiStrip Lorentz Angle") {}
0216 
0217     bool fill() override {
0218       //SiStripPI::setPaletteStyle(SiStripPI::DEFAULT);
0219       gStyle->SetPalette(1);
0220 
0221       auto tag = PlotBase::getTag<0>();
0222       auto iov = tag.iovs.front();
0223       auto tagname = PlotBase::getTag<0>().name;
0224 
0225       std::shared_ptr<SiStripLorentzAngle> payload = fetchPayload(std::get<1>(iov));
0226 
0227       auto theIOVsince = std::to_string(std::get<0>(iov));
0228       std::string titleMap = "SiStrip Lorentz Angle Map, Run: " + theIOVsince + " (tag:#color[2]{" + tagname + "})";
0229 
0230       SiStripTkMaps myMap("COLZA L");
0231       myMap.bookMap(titleMap, "SiStrip #mu_{H}=(tan#theta_{L}/B) [1/T]");
0232 
0233       std::map<uint32_t, float> LAMap_ = payload->getLorentzAngles();
0234 
0235       for (const auto &element : LAMap_) {
0236         myMap.fill(element.first, element.second);
0237       }  // loop over the LA MAP
0238 
0239       std::string fileName(m_imageFileName);
0240       TCanvas canvas("LA map", "LA map");
0241       myMap.drawMap(canvas, "");
0242       canvas.SaveAs(fileName.c_str());
0243 
0244 #ifdef MMDEBUG
0245       canvas.SaveAs("test.root");
0246 #endif
0247       return true;
0248     }
0249   };
0250 
0251   /************************************************
0252     Plot Lorentz Angle averages by partition 
0253   *************************************************/
0254 
0255   class SiStripLorentzAngleByRegion : public PlotImage<SiStripLorentzAngle, SINGLE_IOV> {
0256   public:
0257     SiStripLorentzAngleByRegion()
0258         : PlotImage<SiStripLorentzAngle, SINGLE_IOV>("SiStripLorentzAngle By Region"),
0259           m_trackerTopo{StandaloneTrackerTopology::fromTrackerParametersXMLFile(
0260               edm::FileInPath("Geometry/TrackerCommonData/data/trackerParameters.xml").fullPath())} {}
0261 
0262     bool fill() override {
0263       auto tag = PlotBase::getTag<0>();
0264       auto iov = tag.iovs.front();
0265       std::shared_ptr<SiStripLorentzAngle> payload = fetchPayload(std::get<1>(iov));
0266 
0267       SiStripDetSummary summaryLA{&m_trackerTopo};
0268 
0269       std::map<uint32_t, float> LAMap_ = payload->getLorentzAngles();
0270 
0271       for (const auto &element : LAMap_) {
0272         summaryLA.add(element.first, element.second);
0273       }
0274 
0275       std::map<unsigned int, SiStripDetSummary::Values> map = summaryLA.getCounts();
0276       //=========================
0277 
0278       TCanvas canvas("Partion summary", "partition summary", 1200, 1000);
0279       canvas.cd();
0280       auto h1 = std::make_unique<TH1F>("byRegion",
0281                                        "SiStrip LA average by partition;; average SiStrip Lorentz Angle [rad]",
0282                                        map.size(),
0283                                        0.,
0284                                        map.size());
0285       h1->SetStats(false);
0286       canvas.SetBottomMargin(0.18);
0287       canvas.SetLeftMargin(0.17);
0288       canvas.SetRightMargin(0.05);
0289       canvas.Modified();
0290 
0291       std::vector<int> boundaries;
0292       unsigned int iBin = 0;
0293 
0294       std::string detector;
0295       std::string currentDetector;
0296 
0297       for (const auto &element : map) {
0298         iBin++;
0299         int count = element.second.count;
0300         double mean = (element.second.mean) / count;
0301 
0302         if (currentDetector.empty())
0303           currentDetector = "TIB";
0304 
0305         switch ((element.first) / 1000) {
0306           case 1:
0307             detector = "TIB";
0308             break;
0309           case 2:
0310             detector = "TOB";
0311             break;
0312           case 3:
0313             detector = "TEC";
0314             break;
0315           case 4:
0316             detector = "TID";
0317             break;
0318         }
0319 
0320         h1->SetBinContent(iBin, mean);
0321         h1->GetXaxis()->SetBinLabel(iBin, SiStripPI::regionType(element.first).second);
0322         h1->GetXaxis()->LabelsOption("v");
0323 
0324         if (detector != currentDetector) {
0325           boundaries.push_back(iBin);
0326           currentDetector = detector;
0327         }
0328       }
0329 
0330       h1->GetYaxis()->SetRangeUser(0., h1->GetMaximum() * 1.30);
0331       h1->SetMarkerStyle(20);
0332       h1->SetMarkerSize(1);
0333       h1->Draw("HIST");
0334       h1->Draw("Psame");
0335 
0336       canvas.Update();
0337 
0338       TLine l[boundaries.size()];
0339       unsigned int i = 0;
0340       for (const auto &line : boundaries) {
0341         l[i] = TLine(h1->GetBinLowEdge(line), canvas.GetUymin(), h1->GetBinLowEdge(line), canvas.GetUymax());
0342         l[i].SetLineWidth(1);
0343         l[i].SetLineStyle(9);
0344         l[i].SetLineColor(2);
0345         l[i].Draw("same");
0346         i++;
0347       }
0348 
0349       TLegend legend = TLegend(0.52, 0.82, 0.95, 0.9);
0350       legend.SetHeader((std::get<1>(iov)).c_str(), "C");  // option "C" allows to center the header
0351       legend.AddEntry(h1.get(), ("IOV: " + std::to_string(std::get<0>(iov))).c_str(), "PL");
0352       legend.SetTextSize(0.025);
0353       legend.Draw("same");
0354 
0355       std::string fileName(m_imageFileName);
0356       canvas.SaveAs(fileName.c_str());
0357 
0358       return true;
0359     }
0360 
0361   private:
0362     TrackerTopology m_trackerTopo;
0363   };
0364 
0365   /************************************************
0366     Plot SiStripLorentz Angle averages by partition comparison
0367   *************************************************/
0368 
0369   template <int ntags, IOVMultiplicity nIOVs>
0370   class SiStripLorentzAngleComparatorByRegionBase : public PlotImage<SiStripLorentzAngle, nIOVs, ntags> {
0371   public:
0372     SiStripLorentzAngleComparatorByRegionBase()
0373         : PlotImage<SiStripLorentzAngle, nIOVs, ntags>("SiStripLorentzAngle By Region Comparison"),
0374           m_trackerTopo{StandaloneTrackerTopology::fromTrackerParametersXMLFile(
0375               edm::FileInPath("Geometry/TrackerCommonData/data/trackerParameters.xml").fullPath())} {}
0376 
0377     bool fill() override {
0378       // trick to deal with the multi-ioved tag and two tag case at the same time
0379       auto theIOVs = PlotBase::getTag<0>().iovs;
0380       auto tagname1 = PlotBase::getTag<0>().name;
0381       std::string tagname2 = "";
0382       auto firstiov = theIOVs.front();
0383       SiStripPI::MetaData lastiov;
0384 
0385       // we don't support (yet) comparison with more than 2 tags
0386       assert(this->m_plotAnnotations.ntags < 3);
0387 
0388       if (this->m_plotAnnotations.ntags == 2) {
0389         auto tag2iovs = PlotBase::getTag<1>().iovs;
0390         tagname2 = PlotBase::getTag<1>().name;
0391         lastiov = tag2iovs.front();
0392       } else {
0393         lastiov = theIOVs.back();
0394       }
0395 
0396       std::shared_ptr<SiStripLorentzAngle> f_payload = this->fetchPayload(std::get<1>(firstiov));
0397       std::shared_ptr<SiStripLorentzAngle> l_payload = this->fetchPayload(std::get<1>(lastiov));
0398 
0399       std::string lastIOVsince = std::to_string(std::get<0>(lastiov));
0400       std::string firstIOVsince = std::to_string(std::get<0>(firstiov));
0401 
0402       //=========================
0403       TCanvas canvas("Partion summary", "partition summary", 1200, 1000);
0404       canvas.cd();
0405       canvas.SetBottomMargin(0.18);
0406       canvas.SetLeftMargin(0.13);
0407       canvas.SetRightMargin(0.03);
0408       canvas.SetTopMargin(0.05);
0409       canvas.Modified();
0410 
0411       std::vector<int> boundaries;
0412       std::shared_ptr<TH1F> h_first;
0413       std::shared_ptr<TH1F> h_last;
0414 
0415       fillTheHistogram(f_payload, h_first, boundaries, 0);
0416       fillTheHistogram(l_payload, h_last, boundaries, 1);
0417 
0418       canvas.cd();
0419       h_first->Draw("HIST");
0420       h_last->Draw("Psame");
0421 
0422       canvas.Update();
0423 
0424       TLine l[boundaries.size()];
0425       unsigned int i = 0;
0426       for (const auto &line : boundaries) {
0427         l[i] = TLine(h_first->GetBinLowEdge(line), canvas.GetUymin(), h_first->GetBinLowEdge(line), canvas.GetUymax());
0428         l[i].SetLineWidth(1);
0429         l[i].SetLineStyle(9);
0430         l[i].SetLineColor(2);
0431         l[i].Draw("same");
0432         i++;
0433       }
0434 
0435       auto ltx = TLatex();
0436       ltx.SetTextFont(62);
0437       ltx.SetTextSize(0.045);
0438       ltx.SetTextAlign(11);
0439 
0440       std::unique_ptr<TLegend> legend = std::make_unique<TLegend>(0.50, 0.25, 0.80, 0.35);
0441       if (this->m_plotAnnotations.ntags == 2) {
0442         legend->AddEntry(h_last.get(), ("#color[2]{" + tagname2 + "}").c_str(), "P");
0443         legend->AddEntry(h_first.get(), ("#color[4]{" + tagname1 + "}").c_str(), "L");
0444         legend->SetTextSize(0.024);
0445         ltx.DrawLatexNDC(gPad->GetLeftMargin(),
0446                          1 - gPad->GetTopMargin() + 0.01,
0447                          ("IOV : #color[4]{" + std::to_string(std::get<0>(firstiov)) + "} vs #color[2]{" +
0448                           std::to_string(std::get<0>(lastiov)) + "}")
0449                              .c_str());
0450       } else {
0451         legend->AddEntry(h_last.get(), ("IOV: #color[2]{" + lastIOVsince + "}").c_str(), "P");
0452         legend->AddEntry(h_first.get(), ("IOV: #color[4]{" + firstIOVsince + "}").c_str(), "L");
0453         legend->SetTextSize(0.040);
0454         ltx.DrawLatexNDC(gPad->GetLeftMargin(), 1 - gPad->GetTopMargin() + 0.01, ("Tag: " + tagname1).c_str());
0455 
0456         legend->SetLineColor(kBlack);
0457       }
0458       legend->Draw("same");
0459 
0460       std::string fileName(this->m_imageFileName);
0461       canvas.SaveAs(fileName.c_str());
0462 
0463       return true;
0464     }
0465 
0466   private:
0467     TrackerTopology m_trackerTopo;
0468 
0469     void fillTheHistogram(const std::shared_ptr<SiStripLorentzAngle> &payload,
0470                           std::shared_ptr<TH1F> &hist,
0471                           std::vector<int> &boundaries,
0472                           unsigned int index = 0) {
0473       SiStripDetSummary summaryLA{&m_trackerTopo};
0474       auto LAMap_ = payload->getLorentzAngles();
0475       for (const auto &element : LAMap_) {
0476         summaryLA.add(element.first, element.second);
0477       }
0478 
0479       auto map = summaryLA.getCounts();
0480       hist = std::make_shared<TH1F>(
0481           (Form("byRegion_%i", index)), ";; average SiStrip Lorentz Angle #mu_{H} [1/T]", map.size(), 0., map.size());
0482 
0483       hist->SetStats(false);
0484       if (index == 0) {
0485         hist->SetLineColor(kBlue);
0486         hist->SetMarkerColor(kBlue);
0487         hist->SetLineWidth(2);
0488         hist->SetMarkerStyle(kFourSquaresX);
0489       } else {
0490         hist->SetMarkerStyle(kFourSquaresX);
0491         hist->SetLineColor(kRed);
0492         hist->SetMarkerColor(kRed);
0493       }
0494       unsigned int iBin = 0;
0495 
0496       std::string detector;
0497       std::string currentDetector;
0498 
0499       for (const auto &element : map) {
0500         iBin++;
0501         int count = element.second.count;
0502         double mean = (element.second.mean) / count;
0503 
0504         if (currentDetector.empty())
0505           currentDetector = "TIB";
0506 
0507         switch ((element.first) / 1000) {
0508           case 1:
0509             detector = "TIB";
0510             break;
0511           case 2:
0512             detector = "TOB";
0513             break;
0514           case 3:
0515             detector = "TEC";
0516             break;
0517           case 4:
0518             detector = "TID";
0519             break;
0520         }
0521 
0522         hist->SetBinContent(iBin, mean);
0523         hist->GetXaxis()->SetBinLabel(iBin, SiStripPI::regionType(element.first).second);
0524         hist->GetXaxis()->LabelsOption("v");
0525 
0526         if (detector != currentDetector) {
0527           if (index == 0) {
0528             boundaries.push_back(iBin);
0529           }
0530           currentDetector = detector;
0531         }
0532       }
0533 
0534       hist->GetYaxis()->SetTitleSize(0.04);
0535       hist->GetYaxis()->SetTitleOffset(1.55);
0536       hist->GetYaxis()->CenterTitle(true);
0537       hist->GetYaxis()->SetRangeUser(0., hist->GetMaximum() * 1.30);
0538       hist->SetMarkerSize(2);
0539     }
0540   };
0541 
0542   using SiStripLorentzAngleByRegionCompareSingleTag = SiStripLorentzAngleComparatorByRegionBase<1, MULTI_IOV>;
0543   using SiStripLorentzAngleByRegionCompareTwoTags = SiStripLorentzAngleComparatorByRegionBase<2, SINGLE_IOV>;
0544 
0545 }  // namespace
0546 
0547 PAYLOAD_INSPECTOR_MODULE(SiStripLorentzAngle) {
0548   PAYLOAD_INSPECTOR_CLASS(SiStripLorentzAngleTest);
0549   PAYLOAD_INSPECTOR_CLASS(SiStripLorentzAngleByPartition);
0550   PAYLOAD_INSPECTOR_CLASS(SiStripLorentzAngleValue);
0551   PAYLOAD_INSPECTOR_CLASS(SiStripLorentzAngleTH2PolyTkMap);
0552   PAYLOAD_INSPECTOR_CLASS(SiStripLorentzAngle_TrackerMap);
0553   PAYLOAD_INSPECTOR_CLASS(SiStripLorentzAngleByRegion);
0554   PAYLOAD_INSPECTOR_CLASS(SiStripLorentzAngleCompareByRegion);
0555   PAYLOAD_INSPECTOR_CLASS(SiStripLorentzAngleByRegionCompareSingleTag);
0556   PAYLOAD_INSPECTOR_CLASS(SiStripLorentzAngleByRegionCompareTwoTags);
0557 }