CMS 3D CMS Logo

Public Member Functions | Private Member Functions | Private Attributes

EcalCleaningAlgo Class Reference

#include <EcalCleaningAlgo.h>

List of all members.

Public Member Functions

EcalRecHit::Flags checkTopology (const DetId &id, const EcalRecHitCollection &rhs)
 EcalCleaningAlgo (const edm::ParameterSet &p)
void setFlags (EcalRecHitCollection &rhs)

Private Member Functions

float e4e1 (const DetId &id, const EcalRecHitCollection &rhs)
 yet another function to calculate swiss cross
float e6e2 (const DetId &id, const EcalRecHitCollection &rhs)
bool isNearCrack (const DetId &detid)
 in EB, check if we are near a crack
const std::vector< DetIdneighbours (const DetId &id)
 return the id of the 4 neighbours in the swiss cross
float recHitE (const DetId id, const EcalRecHitCollection &recHits, bool useTimingInfo)

Private Attributes

float cThreshold_barrel_
float cThreshold_double_
float cThreshold_endcap_
float e4e1_a_barrel_
float e4e1_a_endcap_
float e4e1_b_barrel_
float e4e1_b_endcap_
float e4e1Treshold_barrel_
float e4e1Treshold_endcap_
float e6e2thresh_
float ignoreOutOfTimeThresh_
 ignore kOutOfTime above threshold when calculating e4e1
float tightenCrack_e1_double_
float tightenCrack_e1_single_
float tightenCrack_e4e1_single_
float tightenCrack_e6e2_double_

Detailed Description

Definition at line 22 of file EcalCleaningAlgo.h.


Constructor & Destructor Documentation

EcalCleaningAlgo::EcalCleaningAlgo ( const edm::ParameterSet p)

Definition at line 15 of file EcalCleaningAlgo.cc.

References cThreshold_barrel_, cThreshold_double_, cThreshold_endcap_, e4e1_a_barrel_, e4e1_a_endcap_, e4e1_b_barrel_, e4e1_b_endcap_, e4e1Treshold_barrel_, e4e1Treshold_endcap_, e6e2thresh_, edm::ParameterSet::getParameter(), ignoreOutOfTimeThresh_, tightenCrack_e1_double_, tightenCrack_e1_single_, tightenCrack_e4e1_single_, and tightenCrack_e6e2_double_.

                                                           {
  
  cThreshold_barrel_   = p.getParameter<double>("cThreshold_barrel");;
  cThreshold_endcap_   = p.getParameter<double>("cThreshold_endcap");;  
  e4e1_a_barrel_       = p.getParameter<double>("e4e1_a_barrel");
  e4e1_b_barrel_       = p.getParameter<double>("e4e1_b_barrel");
  e4e1_a_endcap_       = p.getParameter<double>("e4e1_a_endcap");
  e4e1_b_endcap_       = p.getParameter<double>("e4e1_b_endcap");
  e4e1Treshold_barrel_ = p.getParameter<double>("e4e1Threshold_barrel");
  e4e1Treshold_endcap_ = p.getParameter<double>("e4e1Threshold_endcap");

  cThreshold_double_   = p.getParameter<double>("cThreshold_double"); 

  ignoreOutOfTimeThresh_   =p.getParameter<double>("ignoreOutOfTimeThresh");  
  tightenCrack_e1_single_  =p.getParameter<double>("tightenCrack_e1_single");
  tightenCrack_e4e1_single_=p.getParameter<double>("tightenCrack_e4e1_single");
  tightenCrack_e1_double_  =p.getParameter<double>("tightenCrack_e1_double");
  tightenCrack_e6e2_double_=p.getParameter<double>("tightenCrack_e6e2_double");
  e6e2thresh_=              p.getParameter<double>("e6e2thresh");


  
}

Member Function Documentation

EcalRecHit::Flags EcalCleaningAlgo::checkTopology ( const DetId id,
const EcalRecHitCollection rhs 
)

check topology, return : kGood : not anomalous kWeird : spike kDiWeird : dispike

Flag spikey channels

Mark single spikes. Spike definition:

Barrel: e> cThreshold_barrel_ && e4e1 > e4e1_a_barrel_ * log10(e) + e4e1_b_barrel_

Near cracks: energy threshold is multiplied by tightenCrack_e1_single e4e1 threshold is divided by tightenCrack_e4e1_single

Endcap : e> cThreshold_endcap_ && e4e1> e4e1_a_endcap_ * log10(e) + e4e1_b_endcap_

Mark double spikes (barrel only) e> cThreshold_double_ && e6e2 > e6e2thresh_;

Near cracks: energy threshold multiplied by tightenCrack_e1_double e6e2 threshold divided by tightenCrack_e6e2_double

Out of time hits above e4e1_IgnoreOutOfTimeThresh_ are ignored in topological quantities

Definition at line 70 of file EcalCleaningAlgo.cc.

References a, b, cThreshold_barrel_, cThreshold_double_, cThreshold_endcap_, e4e1(), e4e1_a_barrel_, e4e1_a_endcap_, e4e1_b_barrel_, e4e1_b_endcap_, e6e2(), e6e2thresh_, EcalBarrel, EcalEndcap, relval_parameters_module::energy, isNearCrack(), EcalRecHit::kDiWeird, EcalSeverityLevel::kGood, EcalSeverityLevel::kWeird, recHitE(), tightenCrack_e1_double_, tightenCrack_e1_single_, tightenCrack_e4e1_single_, and tightenCrack_e6e2_double_.

Referenced by setFlags().

                                                                {


  float a=0,b=0,e4e1thresh=0,ethresh=0;


  if( id.subdetId() == EcalBarrel) {
    a= e4e1_a_barrel_;
    b= e4e1_b_barrel_; 
    ethresh=cThreshold_barrel_;
   
  }
  else if( id.subdetId() == EcalEndcap){
    a= e4e1_a_endcap_;
    b= e4e1_b_endcap_;
    ethresh=cThreshold_endcap_;
   }
  
  

  // for energies below threshold, we don't apply e4e1 cut
  float energy = recHitE(id,rhs,false);
 
  if (energy< ethresh) return EcalRecHit::kGood;
  if (isNearCrack(id) && energy < ethresh*tightenCrack_e1_single_) 
    return EcalRecHit::kGood;


  float e4e1value = e4e1(id,rhs);
  e4e1thresh = a* log10(energy) + b;

  // near cracks the cut is tighter by a factor 
  if (isNearCrack(id)) {
    e4e1thresh/=tightenCrack_e4e1_single_;
  }

  // identify spike
  if (e4e1value < e4e1thresh) return EcalRecHit::kWeird; 
  


  // now for double spikes
 
  // no checking for double spikes in EE
  if( id.subdetId() == EcalEndcap) return EcalRecHit::kGood;

  float e6e2value = e6e2(id,rhs);
  float e6e2thresh = e6e2thresh_ ;
  if (isNearCrack(id) && energy < cThreshold_double_ *tightenCrack_e1_double_ )
    return EcalRecHit::kGood;

  if  (energy <  cThreshold_double_) return EcalRecHit::kGood;
  
  // near cracks the cut is tighter by a factor 
  if (id.subdetId() == EcalBarrel && isNearCrack(id)) 
    e6e2thresh/=tightenCrack_e6e2_double_;

  // identify double spike
  if (e6e2value < e6e2thresh) return EcalRecHit::kDiWeird; 

  return EcalRecHit::kGood;

}
float EcalCleaningAlgo::e4e1 ( const DetId id,
const EcalRecHitCollection rhs 
) [private]

yet another function to calculate swiss cross

Definition at line 140 of file EcalCleaningAlgo.cc.

References i, neighbours(), and recHitE().

Referenced by checkTopology(), and e6e2().

                                                             {

 
  float s4 = 0;
  float e1 = recHitE( id, rhs, false );
  
  
  if ( e1 == 0 ) return 0;
  const std::vector<DetId>& neighs =  neighbours(id);
  for (size_t i=0; i<neighs.size(); ++i)
    // avoid hits out of time when making s4
    s4+=recHitE(neighs[i],rhs, true);
  
  return s4 / e1;
  
 
}
float EcalCleaningAlgo::e6e2 ( const DetId id,
const EcalRecHitCollection rhs 
) [private]

Compute e6 over e2 around xtal 1, where 2 is the most energetic in the swiss cross around 1

| | | +-+-+-+-+ | |1|2| | +-+-+-+-+ | | |

Definition at line 162 of file EcalCleaningAlgo.cc.

References e4e1(), i, neighbours(), and recHitE().

Referenced by checkTopology().

                                                             {

    float s4_1 = 0;
    float s4_2 = 0;
    float e1 = recHitE( id, rhs , false );


    float maxene=0;
    DetId maxid;

    if ( e1 == 0 ) return 0;

    const std::vector<DetId>& neighs =  neighbours(id);

    // find the most energetic neighbour ignoring time info
    for (size_t i=0; i<neighs.size(); ++i){
      float ene = recHitE(neighs[i],rhs,false);
      if (ene>maxene)  {
        maxene=ene;
        maxid = neighs[i];
      }
    }

    float e2=maxene;

    s4_1 = e4e1(id,rhs)* e1;
    s4_2 = e4e1(maxid,rhs)* e2;

    return (s4_1 + s4_2) / (e1+e2) -1. ;

}
bool EcalCleaningAlgo::isNearCrack ( const DetId detid) [private]

in EB, check if we are near a crack

Definition at line 256 of file EcalCleaningAlgo.cc.

References EcalEndcap, isNextToBoundary(), and EEDetId::isNextToRingBoundary().

Referenced by checkTopology().

                                                 {

  if (id.subdetId() == EcalEndcap) { 
    return EEDetId::isNextToRingBoundary(id);
  } else {
    return EBDetId::isNextToBoundary(id);
  }
}
const std::vector< DetId > EcalCleaningAlgo::neighbours ( const DetId id) [private]

return the id of the 4 neighbours in the swiss cross

four neighbours in the swiss cross around id

Definition at line 230 of file EcalCleaningAlgo.cc.

References EcalBarrel, EcalEndcap, EEDetId::offsetBy(), EBDetId::offsetBy(), and run_regression::ret.

Referenced by e4e1(), and e6e2().

                                                                  {
 
  std::vector<DetId> ret;

  if ( id.subdetId() == EcalBarrel) {

    ret.push_back( EBDetId::offsetBy( id,  1, 0 ));
    ret.push_back( EBDetId::offsetBy( id, -1, 0 ));
    ret.push_back( EBDetId::offsetBy( id,  0, 1 ));
    ret.push_back( EBDetId::offsetBy( id,  0,-1 ));
  }
  // nobody understands what polymorphism is for, sgrunt !
  else  if (id.subdetId() == EcalEndcap) {
    ret.push_back( EEDetId::offsetBy( id,  1, 0 ));
    ret.push_back( EEDetId::offsetBy( id, -1, 0 ));
    ret.push_back( EEDetId::offsetBy( id,  0, 1 ));
    ret.push_back( EEDetId::offsetBy( id,  0,-1 ));

  }


  return ret;

} 
float EcalCleaningAlgo::recHitE ( const DetId  id,
const EcalRecHitCollection recHits,
bool  useTimingInfo 
) [private]

Definition at line 198 of file EcalCleaningAlgo.cc.

References e4e1Treshold_barrel_, e4e1Treshold_endcap_, EcalBarrel, EcalEndcap, edm::SortedCollection< T, SORT >::end(), edm::SortedCollection< T, SORT >::find(), ignoreOutOfTimeThresh_, EcalRecHit::kOutOfTime, and dtDQMClient_cfg::threshold.

Referenced by checkTopology(), e4e1(), and e6e2().

{
  if ( id.rawId() == 0 ) return 0;
  

  float threshold = e4e1Treshold_barrel_;
  if ( id.subdetId() == EcalEndcap) threshold = e4e1Treshold_endcap_; 

  EcalRecHitCollection::const_iterator it = recHits.find( id );
  if ( it != recHits.end() ){
    float ene= (*it).energy();

    // ignore out of time in EB when making e4e1 if so configured
    if (useTimingInfo){

      if (id.subdetId()==EcalBarrel &&
          it->checkFlag(EcalRecHit::kOutOfTime) 
          && ene>ignoreOutOfTimeThresh_) return 0;
    }

    // ignore hits below threshold
    if (ene < threshold) return 0;

    // else return the energy of this hit
    return ene;
  }
  return 0;
}
void EcalCleaningAlgo::setFlags ( EcalRecHitCollection rhs)

Definition at line 267 of file EcalCleaningAlgo.cc.

References edm::SortedCollection< T, SORT >::begin(), checkTopology(), edm::SortedCollection< T, SORT >::end(), and EcalRecHit::kGood.

Referenced by EcalRecHitProducer::produce(), and EcalRawToRecHitProducer::produce().

                                                        {
  EcalRecHitCollection::iterator rh;
  //changing the collection on place
  for (rh=rhs.begin(); rh!=rhs.end(); ++rh){
    EcalRecHit::Flags state=checkTopology(rh->id(),rhs);
    if (state!=EcalRecHit::kGood) { 
      rh->unsetFlag(EcalRecHit::kGood);
      rh->setFlag(state);
    }
  }
}

Member Data Documentation

Definition at line 68 of file EcalCleaningAlgo.h.

Referenced by checkTopology(), and EcalCleaningAlgo().

Definition at line 79 of file EcalCleaningAlgo.h.

Referenced by checkTopology(), and EcalCleaningAlgo().

Definition at line 69 of file EcalCleaningAlgo.h.

Referenced by checkTopology(), and EcalCleaningAlgo().

Definition at line 70 of file EcalCleaningAlgo.h.

Referenced by checkTopology(), and EcalCleaningAlgo().

Definition at line 72 of file EcalCleaningAlgo.h.

Referenced by checkTopology(), and EcalCleaningAlgo().

Definition at line 71 of file EcalCleaningAlgo.h.

Referenced by checkTopology(), and EcalCleaningAlgo().

Definition at line 73 of file EcalCleaningAlgo.h.

Referenced by checkTopology(), and EcalCleaningAlgo().

Definition at line 75 of file EcalCleaningAlgo.h.

Referenced by EcalCleaningAlgo(), and recHitE().

Definition at line 76 of file EcalCleaningAlgo.h.

Referenced by EcalCleaningAlgo(), and recHitE().

Definition at line 82 of file EcalCleaningAlgo.h.

Referenced by checkTopology(), and EcalCleaningAlgo().

ignore kOutOfTime above threshold when calculating e4e1

Definition at line 64 of file EcalCleaningAlgo.h.

Referenced by EcalCleaningAlgo(), and recHitE().

Definition at line 80 of file EcalCleaningAlgo.h.

Referenced by checkTopology(), and EcalCleaningAlgo().

Definition at line 77 of file EcalCleaningAlgo.h.

Referenced by checkTopology(), and EcalCleaningAlgo().

Definition at line 78 of file EcalCleaningAlgo.h.

Referenced by checkTopology(), and EcalCleaningAlgo().

Definition at line 81 of file EcalCleaningAlgo.h.

Referenced by checkTopology(), and EcalCleaningAlgo().