Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2021-02-14 12:54:17

0001 #include "DataFormats/ParticleFlowReco/interface/PFBlockElement.h"
0002 #include "DataFormats/ParticleFlowReco/interface/PFBlockElementTrack.h"
0003 #include "DataFormats/ParticleFlowReco/interface/PFBlockElementGsfTrack.h"
0004 #include "DataFormats/ParticleFlowReco/interface/PFBlockElementBrem.h"
0005 #include "DataFormats/ParticleFlowReco/interface/PFBlockElementCluster.h"
0006 #include "DataFormats/ParticleFlowReco/interface/PFBlockElementSuperCluster.h"
0007 
0008 const reco::TrackRef reco::PFBlockElement::nullTrack_ = reco::TrackRef();
0009 const reco::PFRecTrackRef reco::PFBlockElement::nullPFRecTrack_ = reco::PFRecTrackRef();
0010 const reco::PFClusterRef reco::PFBlockElement::nullPFCluster_ = reco::PFClusterRef();
0011 const reco::PFDisplacedTrackerVertexRef reco::PFBlockElement::nullPFDispVertex_ = reco::PFDisplacedTrackerVertexRef();
0012 const reco::ConversionRefVector reco::PFBlockElement::nullConv_ = reco::ConversionRefVector();
0013 const reco::MuonRef reco::PFBlockElement::nullMuon_ = reco::MuonRef();
0014 const reco::VertexCompositeCandidateRef reco::PFBlockElement::nullVertex_ = reco::VertexCompositeCandidateRef();
0015 
0016 using namespace reco;
0017 
0018 // int PFBlockElement::instanceCounter_ = 0;
0019 
0020 // int PFBlockElement::instanceCounter() {
0021 //   return instanceCounter_;
0022 // }
0023 
0024 void PFBlockElement::Dump(std::ostream& out, const char* pad) const {
0025   if (!out)
0026     return;
0027   out << pad << "base element";
0028 }
0029 
0030 std::ostream& reco::operator<<(std::ostream& out, const PFBlockElement& element) {
0031   if (!out)
0032     return out;
0033 
0034   out << "element " << element.index() << "- type " << element.type() << " ";
0035 
0036   try {
0037     switch (element.type()) {
0038       case PFBlockElement::TRACK: {
0039         const reco::PFBlockElementTrack& et = dynamic_cast<const reco::PFBlockElementTrack&>(element);
0040         et.Dump(out);
0041         if (et.trackType(PFBlockElement::T_FROM_DISP))
0042           out << " from displaced;";
0043         if (et.trackType(PFBlockElement::T_TO_DISP))
0044           out << " to displaced;";
0045         if (et.trackType(PFBlockElement::T_FROM_GAMMACONV))
0046           out << " from gammaconv;";
0047         if (et.trackType(PFBlockElement::T_FROM_V0))
0048           out << " from v0 decay;";
0049         break;
0050       }
0051       case PFBlockElement::ECAL:
0052       case PFBlockElement::HCAL:
0053       case PFBlockElement::HGCAL:
0054       case PFBlockElement::HO:
0055       case PFBlockElement::HFEM:
0056       case PFBlockElement::HFHAD:
0057       case PFBlockElement::PS1:
0058       case PFBlockElement::PS2: {
0059         const reco::PFBlockElementCluster& ec = dynamic_cast<const reco::PFBlockElementCluster&>(element);
0060         ec.Dump(out);
0061         break;
0062       }
0063       case PFBlockElement::GSF: {
0064         const reco::PFBlockElementGsfTrack& eg = dynamic_cast<const reco::PFBlockElementGsfTrack&>(element);
0065         eg.Dump(out);
0066         out << " from gsf;";
0067         break;
0068       }
0069       case PFBlockElement::BREM: {
0070         const reco::PFBlockElementBrem& em = dynamic_cast<const reco::PFBlockElementBrem&>(element);
0071         em.Dump(out);
0072         out << " from brem;";
0073         break;
0074       }
0075       case PFBlockElement::SC: {
0076         const reco::PFBlockElementSuperCluster& sc = dynamic_cast<const reco::PFBlockElementSuperCluster&>(element);
0077         sc.Dump(out);
0078         out << " from SuperCluster;";
0079         break;
0080       }
0081       default:
0082         out << " unknown type" << std::endl;
0083         break;
0084     }
0085   } catch (std::exception& err) {
0086     out << err.what() << std::endl;
0087   }
0088 
0089   return out;
0090 }