CMS 3D CMS Logo

Public Member Functions | Private Attributes

CleanAndMergeProducer Class Reference

#include <CleanAndMergeProducer.h>

Inheritance diagram for CleanAndMergeProducer:
edm::EDProducer edm::ProducerBase edm::ProductRegistryHelper

List of all members.

Public Member Functions

 CleanAndMergeProducer (const edm::ParameterSet &ps)
virtual void produce (edm::Event &, const edm::EventSetup &)
 ~CleanAndMergeProducer ()

Private Attributes

std::string bcCollection_
edm::InputTag cleanScInputTag_
std::string hitcollection_
std::string hitproducer_
std::string refScCollection_
std::string scCollection_
edm::InputTag uncleanScInputTag_

Detailed Description

Definition at line 16 of file CleanAndMergeProducer.h.


Constructor & Destructor Documentation

CleanAndMergeProducer::CleanAndMergeProducer ( const edm::ParameterSet ps)

Definition at line 54 of file CleanAndMergeProducer.cc.

References bcCollection_, cleanScInputTag_, edm::ParameterSet::getParameter(), hitcollection_, hitproducer_, refScCollection_, scCollection_, and uncleanScInputTag_.

{

        // get the parameters
        // the cleaned collection:
        cleanScInputTag_   = ps.getParameter<edm::InputTag>("cleanScInputTag");
        uncleanScInputTag_ = ps.getParameter<edm::InputTag>("uncleanScInputTag");          

        // the names of the products to be produced:
        bcCollection_ = ps.getParameter<std::string>("bcCollection");
        scCollection_ = ps.getParameter<std::string>("scCollection");
        refScCollection_ = ps.getParameter<std::string>("refScCollection");

        std::map<std::string,double> providedParameters;  
        providedParameters.insert(std::make_pair("LogWeighted",ps.getParameter<bool>("posCalc_logweight")));
        providedParameters.insert(std::make_pair("T0_barl",ps.getParameter<double>("posCalc_t0")));
        providedParameters.insert(std::make_pair("W0",ps.getParameter<double>("posCalc_w0")));
        providedParameters.insert(std::make_pair("X0",ps.getParameter<double>("posCalc_x0")));

        hitproducer_ = ps.getParameter<std::string>("ecalhitproducer");
        hitcollection_ =ps.getParameter<std::string>("ecalhitcollection");

        // the products:
        produces< reco::BasicClusterCollection >(bcCollection_);
        produces< reco::SuperClusterCollection >(scCollection_);
        produces< reco::SuperClusterRefVector >(refScCollection_);

}
CleanAndMergeProducer::~CleanAndMergeProducer ( )

Definition at line 83 of file CleanAndMergeProducer.cc.

{;}

Member Function Documentation

void CleanAndMergeProducer::produce ( edm::Event evt,
const edm::EventSetup es 
) [virtual]

Implements edm::EDProducer.

Definition at line 85 of file CleanAndMergeProducer.cc.

References bcCollection_, cleanScInputTag_, relval_parameters_module::energy, first, edm::Event::getByLabel(), i, edm::HandleBase::isValid(), edm::OrphanHandleBase::isValid(), LogDebug, LogTrace, edm::RefVector< C, T, F >::push_back(), edm::PtrVector< T >::push_back(), edm::Event::put(), refScCollection_, scCollection_, edm::RefVector< C, T, F >::size(), and uncleanScInputTag_.

{
        // get the input collections
        // ______________________________________________________________________________________
        //
        // cluster collections:

        edm::Handle<reco::SuperClusterCollection> pCleanSC;
        edm::Handle<reco::SuperClusterCollection> pUncleanSC;

        evt.getByLabel(cleanScInputTag_, pCleanSC);
        if (!(pCleanSC.isValid())) 
        {
          edm::LogWarning("MissingInput")<< "could not get a handle on the clean Super Clusters!";
          return;
        }

        evt.getByLabel(uncleanScInputTag_, pUncleanSC);
        if (!(pUncleanSC.isValid())) 
        {
          
          edm::LogWarning("MissingInput")<< "could not get a handle on the unclean Super Clusters!";
          return;
        }

        //
        // collections are all taken now _____________________________________________________
        //
        //
        // the collections to be produced:
        reco::BasicClusterCollection basicClusters;
        reco::SuperClusterCollection superClusters;
        reco::SuperClusterRefVector * scRefs = new reco::SuperClusterRefVector;

        //
        // run over the uncleaned SC and check how many of them are matched to the cleaned ones
        // if you find a matched one, create a reference to the cleaned collection and store it
        // if you find an unmatched one, then keep all its basic clusters in the basic Clusters 
        // vector
        int uncleanSize =  pUncleanSC->size();
        int cleanSize =    pCleanSC->size();

        LogTrace("EcalCleaning") << "Size of Clean Collection: " << cleanSize 
                << ", uncleanSize: " << uncleanSize  << std::endl;
        //
        // keep whether the SC in unique in the uncleaned collection
        std::vector<int> isUncleanOnly;     // 1 if unique in the uncleaned
        std::vector<int> basicClusterOwner; // contains the index of the SC that owns that BS
        std::vector<int> isSeed; // if this basic cluster is a seed it is 1
        for (int isc =0; isc< uncleanSize; ++isc) {
          reco::SuperClusterRef unscRef( pUncleanSC, isc );    
          const std::vector< std::pair<DetId, float> > & uhits = unscRef->hitsAndFractions();
          int uhitsSize = uhits.size();
          bool foundTheSame = false;
          for (int jsc=0; jsc < cleanSize; ++jsc) { // loop over the cleaned SC
            reco::SuperClusterRef cscRef( pCleanSC, jsc );
            const std::vector< std::pair<DetId, float> > & chits = cscRef->hitsAndFractions();
            int chitsSize =  chits.size();
            foundTheSame = true;
            if (unscRef->seed()->seed() == cscRef->seed()->seed() && chitsSize == uhitsSize) { 
              // if the clusters are exactly the same then because the clustering
              // algorithm works in a deterministic way, the order of the rechits
              // will be the same
              for (int i=0; i< chitsSize; ++i) {
                if (uhits[i].first != chits[i].first ) { foundTheSame=false;  break;}
              }
              if (foundTheSame) { // ok you have found it!
                // make the reference
                //scRefs->push_back( edm::Ref<reco::SuperClusterCollection>(pCleanSC, jsc) );
                scRefs->push_back( cscRef );
                isUncleanOnly.push_back(0);
                break;
              }
            }
          }
          if (not foundTheSame) {
            // mark it as unique in the uncleaned
            isUncleanOnly.push_back(1);
            // keep all its basic clusters
            for (reco::CaloCluster_iterator bciter = unscRef->clustersBegin(); bciter != unscRef->clustersEnd(); ++bciter) {
              // the basic clusters
              basicClusters.push_back(**bciter);
              basicClusterOwner.push_back(isc);
            }
          }
        }
        int bcSize =  basicClusters.size();
        
        LogDebug("EcalCleaning") << "Found cleaned SC: " << cleanSize <<  " uncleaned SC: " 
                << uncleanSize << " from which " << scRefs->size() 
                << " will become refs to the cleaned collection" ;
     

        // now you have the collection of basic clusters of the SC to be remain in the
        // in the clean collection, export them to the event
        // you will need to reread them later in order to make correctly the refs to the SC
        std::auto_ptr< reco::BasicClusterCollection> basicClusters_p(new reco::BasicClusterCollection);
        basicClusters_p->assign(basicClusters.begin(), basicClusters.end());
        edm::OrphanHandle<reco::BasicClusterCollection> bccHandle =  
                evt.put(basicClusters_p, bcCollection_);
        if (!(bccHandle.isValid())) {
          edm::LogWarning("MissingInput")<<"could not get a handle on the BasicClusterCollection!" << std::endl;
          return;
        }

        LogDebug("EcalCleaning")<< "Got the BasicClusters from the event again";
        // now you have to create again your superclusters
        // you run over the uncleaned SC, but now you know which of them are
        // the ones that are needed and which are their basic clusters
        for (int isc=0; isc< uncleanSize; ++isc) {
          if (isUncleanOnly[isc]==1) { // look for sc that are unique in the unclean collection
            // make the basic cluster collection
            reco::CaloClusterPtrVector clusterPtrVector;
            reco::CaloClusterPtr seed; // the seed is the basic cluster with the highest energy
            double energy = -1;
            for (int jbc=0; jbc< bcSize; ++jbc) {
              if (basicClusterOwner[jbc]==isc) {
                reco::CaloClusterPtr currentClu = reco::CaloClusterPtr(bccHandle, jbc);
                clusterPtrVector.push_back(currentClu);
                if (energy< currentClu->energy()) {
                  energy = currentClu->energy(); seed = currentClu;
                }
              }
            }
            reco::SuperClusterRef unscRef( pUncleanSC, isc ); 
            reco::SuperCluster newSC(unscRef->energy(), unscRef->position(), seed, clusterPtrVector );
            superClusters.push_back(newSC);
          }
          
        }
        // export the collection of references to the clean collection
        std::auto_ptr< reco::SuperClusterRefVector >  scRefs_p( scRefs );
        evt.put(scRefs_p, refScCollection_);

        // the collection of basic clusters is already in the event
        // the collection of uncleaned SC
        std::auto_ptr< reco::SuperClusterCollection > superClusters_p(new reco::SuperClusterCollection);
        superClusters_p->assign(superClusters.begin(), superClusters.end());
        evt.put(superClusters_p, scCollection_);
        LogDebug("EcalCleaning")<< "Hybrid Clusters (Basic/Super) added to the Event! :-)";
}

Member Data Documentation

std::string CleanAndMergeProducer::bcCollection_ [private]

Definition at line 34 of file CleanAndMergeProducer.h.

Referenced by CleanAndMergeProducer(), and produce().

Definition at line 30 of file CleanAndMergeProducer.h.

Referenced by CleanAndMergeProducer(), and produce().

Definition at line 39 of file CleanAndMergeProducer.h.

Referenced by CleanAndMergeProducer().

std::string CleanAndMergeProducer::hitproducer_ [private]

Definition at line 38 of file CleanAndMergeProducer.h.

Referenced by CleanAndMergeProducer().

Definition at line 36 of file CleanAndMergeProducer.h.

Referenced by CleanAndMergeProducer(), and produce().

std::string CleanAndMergeProducer::scCollection_ [private]

Definition at line 35 of file CleanAndMergeProducer.h.

Referenced by CleanAndMergeProducer(), and produce().

Definition at line 31 of file CleanAndMergeProducer.h.

Referenced by CleanAndMergeProducer(), and produce().