Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:11:50

0001 // ROOT includes
0002 #include "TEveScene.h"
0003 #include "TEveManager.h"
0004 #include "TEveStraightLineSet.h"
0005 #include "TEveTrack.h"
0006 #include "TEveTrackPropagator.h"
0007 #include "TEveTrans.h"
0008 #include "TEveText.h"
0009 #include "TEveGeoShape.h"
0010 #include "TGLViewer.h"
0011 #include "TGLScenePad.h"
0012 #include "TCanvas.h"
0013 #include "TLatex.h"
0014 #include "TLegend.h"
0015 
0016 #include "TH2.h"
0017 
0018 #include "TAxis.h"
0019 #include "TGSlider.h"
0020 #include "TGButton.h"
0021 #include "TGLabel.h"
0022 #include "TGLCameraOverlay.h"
0023 
0024 #include "Fireworks/ParticleFlow/plugins/FWPFCandidateDetailView.h"
0025 #include "Fireworks/Core/interface/fw3dlego_xbins.h"
0026 
0027 #include "Fireworks/Core/interface/FWGUIManager.h"
0028 #include "Fireworks/Core/interface/FWColorManager.h"
0029 #include "Fireworks/Core/interface/FWGeometry.h"
0030 #include "Fireworks/Core/interface/FWEventItem.h"
0031 #include "Fireworks/Core/interface/CSGAction.h"
0032 #include "Fireworks/Core/interface/FWIntValueListener.h"
0033 #include "Fireworks/Core/interface/FWMagField.h"
0034 #include "Fireworks/Core/interface/fwLog.h"
0035 #include "Fireworks/ParticleFlow/interface/FWPFMaths.h"
0036 
0037 #include "FWCore/Utilities/interface/Exception.h"
0038 
0039 #include "FWCore/Framework/interface/Event.h"
0040 #include "DataFormats/Candidate/interface/Candidate.h"
0041 #include "DataFormats/ParticleFlowReco/interface/PFRecHit.h"
0042 #include "DataFormats/ParticleFlowReco/interface/PFCluster.h"
0043 #include "DataFormats/ParticleFlowReco/interface/PFRecTrack.h"
0044 #include "DataFormats/ParticleFlowReco/interface/PFRecHitFraction.h"
0045 
0046 #include <functional>
0047 
0048 FWPFCandidateDetailView::FWPFCandidateDetailView()
0049     : m_range(1),
0050       m_candidate(nullptr),
0051       m_legend(nullptr),
0052       m_slider(nullptr),
0053       m_sliderListener(),
0054       m_eventList(nullptr),
0055       m_plotEt(true),
0056       m_rnrHcal(true) {}
0057 
0058 FWPFCandidateDetailView::~FWPFCandidateDetailView() {}
0059 
0060 float FWPFCandidateDetailView::eta() { return m_candidate->eta(); }
0061 
0062 float FWPFCandidateDetailView::phi() { return m_candidate->phi(); }
0063 
0064 bool FWPFCandidateDetailView::isPntInRng(float x, float y) {
0065   float dx = m_candidate->eta() - x;
0066   float dy = m_candidate->phi() - y;
0067   float sd = TMath::Sqrt(dx * dx + dy * dy);
0068   return sd < m_range;
0069 }
0070 
0071 //______________________________________________________________________________
0072 
0073 void FWPFCandidateDetailView::makeLegend() {
0074   m_legend = new TLegend(0.01, 0.01, 0.99, 0.99, nullptr, "NDC");
0075   m_legend->SetFillColor(kWhite);
0076   m_legend->SetTextSize(0.07);
0077   m_legend->SetBorderSize(0);
0078   m_legend->SetMargin(0.15);
0079   m_legend->SetEntrySeparation(0.01);
0080 }
0081 
0082 //______________________________________________________________________________
0083 
0084 void FWPFCandidateDetailView::rangeChanged(int x) {
0085   static float kSliderRangeFactor = 0.2;
0086 
0087   m_range = x * kSliderRangeFactor;
0088 
0089   if (m_eventList)
0090     buildGLEventScene();
0091 
0092   gEve->Redraw3D();
0093 }
0094 
0095 //______________________________________________________________________________
0096 
0097 void FWPFCandidateDetailView::setTextInfo(const FWModelId& id, const reco::PFCandidate* track) {
0098   m_infoCanvas->cd();
0099 
0100   float_t x = 0.02;
0101   float y = 0.95;
0102 
0103   TLatex* latex = new TLatex(x, y, "");
0104   const double textsize(0.07);
0105   latex->SetTextSize(textsize);
0106 
0107   latex->DrawLatex(x, y, id.item()->modelName(id.index()).c_str());
0108   y -= latex->GetTextSize() * 0.6;
0109 
0110   latex->SetTextSize(textsize);
0111   float lineH = latex->GetTextSize() * 0.6;
0112 
0113   latex->DrawLatex(
0114       x, y, Form(" P_{T} = %.1f GeV, #eta = %0.2f, #varphi = %0.2f", track->pt(), track->eta(), track->phi()));
0115   y -= lineH;
0116 
0117   if (track->charge() > 0)
0118     latex->DrawLatex(x, y, " charge = +1");
0119   else
0120     latex->DrawLatex(x, y, " charge = -1");
0121   y -= lineH;
0122   y -= lineH;
0123 
0124   m_legend->SetY2(y);
0125   m_legend->Draw();
0126   m_legend = nullptr;  // Deleted together with TPad.
0127 }
0128 
0129 void FWPFCandidateDetailView::plotEtChanged() {
0130   printf("plotEt = %d \n", m_plotEt);
0131   m_plotEt = !m_plotEt;
0132   buildGLEventScene();
0133 }
0134 
0135 void FWPFCandidateDetailView::rnrHcalChanged() {
0136   printf("rnrHcal = %d \n", m_rnrHcal);
0137   m_rnrHcal = !m_rnrHcal;
0138   buildGLEventScene();
0139 }
0140 
0141 //______________________________________________________________________________
0142 
0143 void FWPFCandidateDetailView::build(const FWModelId& id, const reco::PFCandidate* candidate) {
0144   m_candidate = candidate;
0145 
0146   // ROOT GUI
0147   //
0148   {
0149     TGCompositeFrame* f = new TGVerticalFrame(m_guiFrame);
0150     m_guiFrame->AddFrame(f);
0151     f->AddFrame(new TGLabel(f, "Rng:"), new TGLayoutHints(kLHintsLeft, 2, 2, 0, 0));
0152     m_slider = new TGHSlider(f, 120, kSlider1 | kScaleNo);
0153     f->AddFrame(m_slider, new TGLayoutHints(kLHintsTop | kLHintsLeft, 2, 2, 1, 4));
0154     m_slider->SetRange(1, 50);
0155     m_slider->SetPosition(8);
0156 
0157     m_sliderListener = new FWIntValueListener();
0158     TQObject::Connect(
0159         m_slider, "PositionChanged(Int_t)", "FWIntValueListenerBase", m_sliderListener, "setValue(Int_t)");
0160     m_sliderListener->valueChanged_.connect(
0161         std::bind(&FWPFCandidateDetailView::rangeChanged, this, std::placeholders::_1));
0162     {
0163       CSGAction* action = new CSGAction(this, "Scale Et");
0164       TGCheckButton* b = new TGCheckButton(m_guiFrame, action->getName().c_str());
0165       b->SetState(kButtonDown, true);
0166       m_guiFrame->AddFrame(b, new TGLayoutHints(kLHintsNormal, 2, 3, 1, 4));
0167       TQObject::Connect(b, "Clicked()", "CSGAction", action, "activate()");
0168       action->activated.connect(sigc::mem_fun(*this, &FWPFCandidateDetailView::plotEtChanged));
0169     }
0170     {
0171       CSGAction* action = new CSGAction(this, "RnrHcal");
0172       TGCheckButton* b = new TGCheckButton(m_guiFrame, action->getName().c_str());
0173       b->SetState(kButtonDown, true);
0174       m_guiFrame->AddFrame(b, new TGLayoutHints(kLHintsNormal, 2, 3, 1, 4));
0175       TQObject::Connect(b, "Clicked()", "CSGAction", action, "activate()");
0176       action->activated.connect(sigc::mem_fun(*this, &FWPFCandidateDetailView::rnrHcalChanged));
0177     }
0178   }
0179   makeLegend();
0180   setTextInfo(id, candidate);
0181 
0182   TGCompositeFrame* p = (TGCompositeFrame*)m_guiFrame->GetParent();
0183   p->MapSubwindows();
0184   p->Layout();
0185 
0186   ///////////////
0187   // GL stuff
0188 
0189   m_candidate = candidate;
0190 
0191   try {
0192     const edm::EventBase* event = FWGUIManager::getGUIManager()->getCurrentEvent();
0193     edm::Handle<std::vector<reco::PFRecHit> > ecalH;
0194     event->getByLabel(edm::InputTag("particleFlowRecHitECAL"), ecalH);
0195     if (ecalH.product())
0196       voteMaxEtEVal(ecalH.product());
0197 
0198     edm::Handle<std::vector<reco::PFRecHit> > hcalH;
0199     event->getByLabel(edm::InputTag("particleFlowRecHitHBHEHO"), hcalH);
0200     if (hcalH.product())
0201       voteMaxEtEVal(hcalH.product());
0202   } catch (const cms::Exception& iE) {
0203     std::cerr << iE.what();
0204   }
0205 
0206   m_eveScene->GetGLScene()->SetSelectable(false);
0207   m_eventList = new TEveElementList("PFDetailView");
0208   m_eveScene->AddElement(m_eventList);
0209 
0210   viewerGL()->SetStyle(TGLRnrCtx::kOutline);
0211   viewerGL()->SetCurrentCamera(TGLViewer::kCameraOrthoXOY);
0212 
0213   TGLCameraOverlay* co = viewerGL()->GetCameraOverlay();
0214   co->SetShowOrthographic(kTRUE);
0215   co->SetOrthographicMode(TGLCameraOverlay::kAxis);
0216 
0217   viewerGL()->ResetCamerasAfterNextUpdate();
0218   try {
0219     buildGLEventScene();
0220   } catch (...) {
0221     printf("unknown exception \n");
0222   }
0223 
0224   viewerGL()->UpdateScene(kFALSE);
0225 
0226   gEve->Redraw3D();
0227 
0228   // gEve->AddToListTree(m_eventList, true);//debug, used with --eve option
0229 }
0230 
0231 //______________________________________________________________________________
0232 
0233 void FWPFCandidateDetailView::voteMaxEtEVal(const std::vector<reco::PFRecHit>* hits) {
0234   if (!hits)
0235     return;
0236 
0237   // FIXME: require access to geometry while reading from reco file
0238   if ((!hits->empty()) && hits->front().hasCaloCell())
0239     for (std::vector<reco::PFRecHit>::const_iterator it = hits->begin(); it != hits->end(); ++it) {
0240       TEveVector centre(it->position().x(), it->position().y(), it->position().z());
0241       float E = it->energy();
0242       float Et = FWPFMaths::calculateEt(centre, E);
0243       item()->context().voteMaxEtAndEnergy(Et, E);
0244     }
0245 }
0246 
0247 //______________________________________________________________________________
0248 
0249 void FWPFCandidateDetailView::addTracks(const std::vector<reco::PFRecTrack>* tracks) {
0250   for (std::vector<reco::PFRecTrack>::const_iterator it = tracks->begin(); it != tracks->end(); ++it) {
0251     /// AMT trackRef() is a collection !!!
0252     /*
0253       if (!isPntInRng(it->trackRef().innerMomentum().Eta(), it->position().Phi()))
0254          continue;
0255 
0256       TEveLine* line = new TEveLine("Track");
0257       line->SetMainColor(kYellow);
0258       int  N = it->nTrajectoryPoints();
0259 
0260       for (int p = 0 ; p<N; ++p) {
0261          pos = track.extrapolatedPoint(p).position();
0262 
0263          if( pos.Eta() !=0 and pos.Phi() !=0)
0264             line->SetNextPoint(pos.Eta(), pos.Phi(), 0);
0265       }
0266       m_eventList->AddElement(line);
0267       */
0268   }
0269 }
0270 
0271 //______________________________________________________________________________
0272 
0273 void FWPFCandidateDetailView::addClusters(const std::vector<reco::PFCluster>* cluster) {
0274   if (!cluster)
0275     return;
0276 
0277   Color_t col = kViolet + 9;
0278 
0279   TEveStraightLineSet* ls = new TEveStraightLineSet("cluster_ls");
0280   ls->SetMainColor(col);
0281   m_eventList->AddElement(ls);
0282 
0283   TEvePointSet* ps = new TEvePointSet("cluster_ps");
0284   ps->SetMainColor(col);
0285   ps->SetMarkerStyle(2);
0286   ps->SetMarkerSize(0.005);
0287   m_eventList->AddElement(ps);
0288 
0289   for (std::vector<reco::PFCluster>::const_iterator it = cluster->begin(); it != cluster->end(); ++it) {
0290     if (!isPntInRng(it->position().Eta(), it->position().Phi()))
0291       continue;
0292 
0293     ps->SetNextPoint(it->position().Eta(), it->position().Phi(), 0);
0294 
0295     /*
0296       const std::vector< reco::PFRecHitFraction >& fractions = it->recHitFractions();
0297       for (std::vector< reco::PFRecHitFraction >::const_iterator fi = fractions.begin(); fi != fractions.end(); ++fi)
0298       {
0299          // !!! AMT can't get  fi->recHitRef().position()
0300          //    ls->AddLine(it->position().Eta(), it->position().Phi(), 0,
0301          //      fi->recHitRef().position().Eta(), fi->recHitRef().position().Phi(), 0);
0302       }
0303       */
0304   }
0305 }
0306 namespace {
0307   void WrapTwoPi(std::vector<TEveVector>& hc, float y) {
0308     if (TMath::Abs(hc[0].fY) < 2)
0309       return;
0310 
0311     if (hc[0].fY > 0 && hc[1].fY > 0 && hc[2].fY > 0 && hc[3].fY > 0)
0312       return;
0313     if (hc[0].fY < 0 && hc[1].fY < 0 && hc[2].fY < 0 && hc[3].fY < 0)
0314       return;
0315 
0316     for (int i = 0; i < 4; ++i)
0317       if (y > 0 && hc[i].fY < 0)
0318         hc[i].fY += TMath::TwoPi();
0319       else if (y < 0 && hc[i].fY > 0)
0320         hc[i].fY -= TMath::TwoPi();
0321   }
0322 }  // namespace
0323 //______________________________________________________________________________
0324 namespace {
0325   TEveStraightLineSet::Line_t* AddLineToLineSet(TEveStraightLineSet* ls,
0326                                                 const std::vector<TEveVector>& pnts,
0327                                                 int i0,
0328                                                 int i1) {
0329     if (false) {
0330       printf("add line \n");
0331       pnts[i0].Dump();
0332       pnts[i1].Dump();
0333     }
0334     return ls->AddLine(pnts[i0], pnts[i1]);
0335     // return ls->AddLine(pnts[i0].Eta(),pnts[i0].Phi(), 0 , pnts[i1].Eta(),pnts[i1].Phi(), 0);
0336   }
0337 }  // namespace
0338 void FWPFCandidateDetailView::addHits(const std::vector<reco::PFRecHit>* hits) {
0339   TEveStraightLineSet* lsOutline = (TEveStraightLineSet*)m_eventList->FindChild("outlines");
0340 
0341   TEvePointSet* ps = new TEvePointSet("test");
0342   m_eventList->AddElement(ps);
0343   ps->SetMainColor(kOrange);
0344 
0345   // FIXME, requires access to geometry
0346   if ((!hits->empty()) && hits->front().hasCaloCell())
0347     for (std::vector<reco::PFRecHit>::const_iterator it = hits->begin(); it != hits->end(); ++it) {
0348       const auto& corners = it->getCornersXYZ();
0349       if (!isPntInRng(corners[0].eta(), corners[0].phi()))
0350         continue;
0351 
0352       std::vector<TEveVector> hc;
0353       hc.reserve(4);
0354       for (int k = 0; k < 4; ++k) {
0355         hc.push_back(TEveVector(corners[k].eta(), corners[k].phi(), 0));
0356         // ps->SetNextPoint(corners[k].eta(),corners[k].phi(),0 ); //debug
0357       }
0358 
0359       WrapTwoPi(hc, corners[0].phi());
0360 
0361       AddLineToLineSet(lsOutline, hc, 0, 1);
0362       AddLineToLineSet(lsOutline, hc, 1, 2);
0363       AddLineToLineSet(lsOutline, hc, 2, 3);
0364       AddLineToLineSet(lsOutline, hc, 3, 0);
0365 
0366       // get scaled corners
0367       TEveVector centerOfGravity = hc[0] + hc[1] + hc[2] + hc[3];
0368       centerOfGravity *= 0.25;
0369 
0370       std::vector<TEveVector> radialVectors;
0371       radialVectors.reserve(4);
0372       for (int k = 0; k < 4; ++k)
0373         radialVectors.push_back(TEveVector(hc[k] - centerOfGravity));
0374 
0375       float factor = 1;
0376       if (m_plotEt) {
0377         float Et = FWPFMaths::calculateEt(TEveVector(corners[0].x(), corners[0].y(), corners[0].z()), it->energy());
0378         factor = Et / context().getMaxEnergyInEvent(m_plotEt);
0379       } else
0380         factor = it->energy() / context().getMaxEnergyInEvent(false);
0381 
0382       std::vector<TEveVector> scaledCorners;
0383       for (int k = 0; k < 4; ++k) {
0384         radialVectors[k] *= factor;
0385         scaledCorners.push_back(TEveVector(radialVectors[k] + centerOfGravity));
0386       }
0387 
0388       TEveStraightLineSet* ls = (TEveStraightLineSet*)m_eventList->FindChild(Form("%d_rechit", it->depth()));
0389       AddLineToLineSet(ls, scaledCorners, 0, 1);
0390       AddLineToLineSet(ls, scaledCorners, 1, 2);
0391       AddLineToLineSet(ls, scaledCorners, 2, 3);
0392       // AddLineToLineSet(ls, scaledCorners, 3, 0);
0393       TEveStraightLineSet::Line_t* li = AddLineToLineSet(ls, scaledCorners, 3, 0);
0394       ls->AddMarker(centerOfGravity, li->fId);
0395     }
0396 }
0397 
0398 //______________________________________________________________________________
0399 
0400 void FWPFCandidateDetailView::buildGLEventScene() {
0401   if (m_eventList->HasChildren())
0402     m_eventList->DestroyElements();
0403 
0404   for (int depth = 0; depth < 6; ++depth) {
0405     TEveStraightLineSet* ls = new TEveStraightLineSet(Form("%d_rechit", depth));
0406 
0407     if (depth == 0)
0408       ls->SetLineColor(kGray);
0409     else if (depth == 1)
0410       ls->SetLineColor(kRed);
0411     else if (depth == 2)
0412       ls->SetLineColor(kGreen);
0413     else if (depth == 3)
0414       ls->SetLineColor(kMagenta);
0415     else if (depth == 4)
0416       ls->SetLineColor(kOrange);
0417     else if (depth == 5)
0418       ls->SetLineColor(kYellow);
0419 
0420     ls->SetMarkerStyle(1);
0421     m_eventList->AddElement(ls);
0422   }
0423 
0424   TEveStraightLineSet* ls = new TEveStraightLineSet("outlines");
0425   ls->SetLineColor(kGray);
0426   ls->SetMainTransparency(80);
0427   m_eventList->AddElement(ls);
0428 
0429   const edm::EventBase* event = FWGUIManager::getGUIManager()->getCurrentEvent();
0430 
0431   //
0432   // recHits
0433   //
0434   try {
0435     edm::Handle<std::vector<reco::PFRecHit> > ecalH;
0436     event->getByLabel(edm::InputTag("particleFlowRecHitECAL"), ecalH);
0437     addHits(ecalH.product());
0438   } catch (const cms::Exception& iE) {
0439     std::cerr << iE.what();
0440   }
0441 
0442   if (m_rnrHcal) {
0443     try {
0444       edm::Handle<std::vector<reco::PFRecHit> > ecalH;
0445       event->getByLabel(edm::InputTag("particleFlowRecHitHF"), ecalH);
0446       addHits(ecalH.product());
0447     } catch (const cms::Exception& iE) {
0448       std::cerr << iE.what();
0449     }
0450 
0451     try {
0452       edm::Handle<std::vector<reco::PFRecHit> > hcalH;
0453       event->getByLabel(edm::InputTag("particleFlowRecHitHBHEHO"), hcalH);
0454       addHits(hcalH.product());
0455     } catch (const cms::Exception& iE) {
0456       std::cerr << iE.what();
0457     }
0458   }
0459 
0460   //
0461   // clusters
0462   //
0463   try {
0464     edm::Handle<std::vector<reco::PFCluster> > ecalClustersH;
0465     event->getByLabel(edm::InputTag("particleFlowClusterECAL"), ecalClustersH);
0466     addClusters(ecalClustersH.product());
0467   } catch (const cms::Exception& iE) {
0468     std::cerr << iE.what();
0469   }
0470 
0471   //
0472   // tracks
0473   //
0474   try {
0475     edm::Handle<std::vector<reco::PFRecTrack> > trackH;
0476     event->getByLabel(edm::InputTag("pfTrack"), trackH);
0477     addTracks(trackH.product());
0478   } catch (const cms::Exception& iE) {
0479     std::cerr << iE.what();
0480   }
0481 }
0482 
0483 REGISTER_FWDETAILVIEW(FWPFCandidateDetailView,
0484                       PF Candidate,
0485                       particleFlowRecHitECAL& particleFlowRecHitHF& particleFlowClusterECAL);