Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:27:27

0001 #ifndef KDTreeLinkerBase_h
0002 #define KDTreeLinkerBase_h
0003 
0004 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0005 #include "DataFormats/ParticleFlowReco/interface/PFRecHitFraction.h"
0006 #include "DataFormats/ParticleFlowReco/interface/PFBlockElement.h"
0007 #include "DataFormats/ParticleFlowReco/interface/PFRecHit.h"
0008 
0009 #include <vector>
0010 #include <map>
0011 #include <set>
0012 
0013 using BlockEltSet = std::set<reco::PFBlockElement *>;
0014 using RecHitSet = std::set<const reco::PFRecHit *>;
0015 
0016 using RecHit2BlockEltMap = std::map<const reco::PFRecHit *, BlockEltSet>;
0017 using BlockElt2BlockEltMap = std::map<reco::PFBlockElement *, BlockEltSet>;
0018 
0019 class KDTreeLinkerBase {
0020 public:
0021   KDTreeLinkerBase(const edm::ParameterSet &conf) {}
0022   virtual ~KDTreeLinkerBase() {}
0023 
0024   void setTargetType(const reco::PFBlockElement::Type &tgt) { _targetType = tgt; }
0025 
0026   void setFieldType(const reco::PFBlockElement::Type &fld) { _fieldType = fld; }
0027 
0028   const reco::PFBlockElement::Type &targetType() const { return _targetType; }
0029 
0030   const reco::PFBlockElement::Type &fieldType() const { return _fieldType; }
0031 
0032   // Get/Set of the maximal size of the cristal (ECAL, HCAL,...) in phi/eta and
0033   // X/Y. By default, thus value are set for the ECAL cristal.
0034 
0035   // Get/Set phi offset. See bellow in the description of phiOffset_ to understand
0036   // the application.
0037 
0038   // Debug flag.
0039   void setDebug(bool isDebug);
0040 
0041   // With this method, we create the list of elements that we want to link.
0042   virtual void insertTargetElt(reco::PFBlockElement *target) = 0;
0043 
0044   // Here, we create the list of cluster that we want to link. From cluster
0045   // and fraction, we will create a second list of rechits that will be used to
0046   // build the KDTree.
0047   virtual void insertFieldClusterElt(reco::PFBlockElement *cluster) = 0;
0048 
0049   // The KDTree building from rechits list.
0050   virtual void buildTree() = 0;
0051 
0052   // Here we will iterate over all target elements. For each one, we will search the closest
0053   // rechits in the KDTree, from rechits we will find the associated clusters and after that
0054   // we will check the links between the target and all closest clusters.
0055   virtual void searchLinks() = 0;
0056 
0057   // Here, we will store all target/cluster founded links in the PFBlockElement class
0058   // of each target in the PFmultilinks field.
0059   virtual void updatePFBlockEltWithLinks() = 0;
0060 
0061   // Here we free all allocated structures.
0062   virtual void clear() = 0;
0063 
0064   // This method calls is the good order buildTree(), searchLinks(),
0065   // updatePFBlockEltWithLinks() and clear()
0066   inline void process() {
0067     buildTree();
0068     searchLinks();
0069     updatePFBlockEltWithLinks();
0070     clear();
0071   }
0072 
0073 protected:
0074   // target and field
0075   reco::PFBlockElement::Type _targetType, _fieldType;
0076   // Cristal maximal size. By default, thus value are set for the ECAL cristal.
0077   float cristalPhiEtaMaxSize_ = 0.04;
0078   float cristalXYMaxSize_ = 3.;
0079 
0080   // Usually, phi is between -Pi and +Pi. But phi space is circular, that's why an element
0081   // with phi = 3.13 and another with phi = -3.14 are close. To solve this problem, during
0082   // the kdtree building step, we duplicate some elements close enough to +Pi (resp -Pi) by
0083   // substracting (adding) 2Pi. This field define the threshold of this operation.
0084   float phiOffset_ = 0.25;
0085 
0086   // rechit with fraction this value will be ignored in KDTreeLinker
0087   const float cutOffFrac = 1E-4;
0088 
0089   // Debug boolean. Not used until now.
0090   bool debug_ = false;
0091 
0092   // Sorted indexes
0093   template <typename T>
0094   static std::vector<size_t> sort_indexes(const std::vector<T> &v) {
0095     std::vector<size_t> idx(v.size());
0096     for (size_t i = 0; i != idx.size(); ++i)
0097       idx[i] = i;
0098     std::sort(idx.begin(), idx.end(), [&v](size_t i1, size_t i2) { return v[i1] < v[i2]; });
0099     return idx;
0100   }
0101 };
0102 
0103 #include "FWCore/PluginManager/interface/PluginFactory.h"
0104 typedef edmplugin::PluginFactory<KDTreeLinkerBase *(const edm::ParameterSet &)> KDTreeLinkerFactory;
0105 
0106 #endif /* !KDTreeLinkerBase_h */