CMS 3D CMS Logo

Public Member Functions | Private Types | Private Member Functions | Private Attributes

HcalTriggerPrimitiveAlgo Class Reference

#include <HcalTriggerPrimitiveAlgo.h>

List of all members.

Public Member Functions

 HcalTriggerPrimitiveAlgo (bool pf, const std::vector< double > &w, int latency, uint32_t FG_threshold, uint32_t ZS_threshold, int numberOfSamples, int numberOfPresamples, uint32_t minSignalThreshold=0, uint32_t PMT_NoiseThreshold=0)
void run (const HcalTPGCoder *incoder, const HcalTPGCompressor *outcoder, const HBHEDigiCollection &hbheDigis, const HFDigiCollection &hfDigis, HcalTrigPrimDigiCollection &result, float rctlsb)
void runFEFormatError (const FEDRawDataCollection *rawraw, const HcalElectronicsMap *emap, HcalTrigPrimDigiCollection &result)
void runZS (HcalTrigPrimDigiCollection &tp)
void setPeakFinderAlgorithm (int algo)
 ~HcalTriggerPrimitiveAlgo ()

Private Types

typedef std::map
< HcalTrigTowerDetId,
std::vector< bool > > 
FGbitMap
typedef std::vector
< IntegerCaloSamples
SumFGContainer
typedef std::map
< HcalTrigTowerDetId,
IntegerCaloSamples
SumMap
typedef std::map
< HcalTrigTowerDetId,
SumFGContainer
TowerMapFGSum
typedef std::map< uint32_t,
std::vector< bool > > 
TowerMapVeto

Private Member Functions

void addFG (const HcalTrigTowerDetId &id, std::vector< bool > &msb)
void addSignal (const HBHEDataFrame &frame)
 adds the signal to the map
void addSignal (const HFDataFrame &frame)
void addSignal (const IntegerCaloSamples &samples)
void analyze (IntegerCaloSamples &samples, HcalTriggerPrimitiveDigi &result)
 adds the actual RecHits
void analyzeHF (IntegerCaloSamples &samples, HcalTriggerPrimitiveDigi &result, float rctlsb)

Private Attributes

uint32_t FG_threshold_
FGbitMap fgMap_
TowerMapVeto HF_Veto
const HcaluLUTTPGCoderincoder_
int latency_
uint32_t minSignalThreshold_
int numberOfPresamples_
int numberOfSamples_
const HcalTPGCompressoroutcoder_
int peak_finder_algorithm_
bool peakfind_
uint32_t PMT_NoiseThreshold_
SumMap theSumMap
double theThreshold
TowerMapFGSum theTowerMapFGSum
HcalTrigTowerGeometry theTrigTowerGeometry
std::vector< double > weights_
uint32_t ZS_threshold_
int ZS_threshold_I_

Detailed Description

Definition at line 20 of file HcalTriggerPrimitiveAlgo.h.


Member Typedef Documentation

typedef std::map<HcalTrigTowerDetId, std::vector<bool> > HcalTriggerPrimitiveAlgo::FGbitMap [private]

Definition at line 98 of file HcalTriggerPrimitiveAlgo.h.

Definition at line 82 of file HcalTriggerPrimitiveAlgo.h.

Definition at line 79 of file HcalTriggerPrimitiveAlgo.h.

Definition at line 83 of file HcalTriggerPrimitiveAlgo.h.

typedef std::map<uint32_t, std::vector<bool> > HcalTriggerPrimitiveAlgo::TowerMapVeto [private]

Definition at line 95 of file HcalTriggerPrimitiveAlgo.h.


Constructor & Destructor Documentation

HcalTriggerPrimitiveAlgo::HcalTriggerPrimitiveAlgo ( bool  pf,
const std::vector< double > &  w,
int  latency,
uint32_t  FG_threshold,
uint32_t  ZS_threshold,
int  numberOfSamples,
int  numberOfPresamples,
uint32_t  minSignalThreshold = 0,
uint32_t  PMT_NoiseThreshold = 0 
)

Definition at line 15 of file HcalTriggerPrimitiveAlgo.cc.

References numberOfPresamples_, numberOfSamples_, peakfind_, ZS_threshold_, and ZS_threshold_I_.

                                                   : incoder_(0), outcoder_(0),
                                                   theThreshold(0), peakfind_(pf), weights_(w), latency_(latency),
                                                   FG_threshold_(FG_threshold), ZS_threshold_(ZS_threshold),
                                                   numberOfSamples_(numberOfSamples),
                                                   numberOfPresamples_(numberOfPresamples),
                                                   minSignalThreshold_(minSignalThreshold),
                                                   PMT_NoiseThreshold_(PMT_NoiseThreshold),
                                                   peak_finder_algorithm_(2)
{
   //No peak finding setting (for Fastsim)
   if (!peakfind_){
      numberOfSamples_ = 1; 
      numberOfPresamples_ = 0;
   }
   // Switch to integer for comparisons - remove compiler warning
   ZS_threshold_I_ = ZS_threshold_;
}
HcalTriggerPrimitiveAlgo::~HcalTriggerPrimitiveAlgo ( )

Definition at line 38 of file HcalTriggerPrimitiveAlgo.cc.

                                                    {
}

Member Function Documentation

void HcalTriggerPrimitiveAlgo::addFG ( const HcalTrigTowerDetId id,
std::vector< bool > &  msb 
) [private]

Definition at line 353 of file HcalTriggerPrimitiveAlgo.cc.

References fgMap_, and i.

Referenced by addSignal().

                                                                                      {
   FGbitMap::iterator itr = fgMap_.find(id);
   if (itr != fgMap_.end()){
      std::vector<bool>& _msb = itr->second;
      for (size_t i=0; i<msb.size(); ++i)
         _msb[i] = _msb[i] || msb[i];
   }
   else fgMap_[id] = msb;
}
void HcalTriggerPrimitiveAlgo::addSignal ( const HBHEDataFrame frame) [private]

adds the signal to the map

Definition at line 81 of file HcalTriggerPrimitiveAlgo.cc.

References HcaluLUTTPGCoder::adc2Linear(), addFG(), HcalDetId::depth(), i, HBHEDataFrame::id(), incoder_, HcaluLUTTPGCoder::lookupMSB(), HBHEDataFrame::presamples(), HBHEDataFrame::size(), theTrigTowerGeometry, and HcalTrigTowerGeometry::towerIds().

Referenced by addSignal(), and run().

                                                                    {
   //Hack for 300_pre10, should be removed.
   if (frame.id().depth()==5) return;

   std::vector<HcalTrigTowerDetId> ids = theTrigTowerGeometry.towerIds(frame.id());
   assert(ids.size() == 1 || ids.size() == 2);
   IntegerCaloSamples samples1(ids[0], int(frame.size()));

   samples1.setPresamples(frame.presamples());
   incoder_->adc2Linear(frame, samples1);

   std::vector<bool> msb;
   incoder_->lookupMSB(frame, msb);

   if(ids.size() == 2) {
      // make a second trigprim for the other one, and split the energy
      IntegerCaloSamples samples2(ids[1], samples1.size());
      for(int i = 0; i < samples1.size(); ++i) {
         samples1[i] = uint32_t(samples1[i]*0.5);
         samples2[i] = samples1[i];
      }
      samples2.setPresamples(frame.presamples());
      addSignal(samples2);
      addFG(ids[1], msb);
   }
   addSignal(samples1);
   addFG(ids[0], msb);
}
void HcalTriggerPrimitiveAlgo::addSignal ( const HFDataFrame frame) [private]

Definition at line 111 of file HcalTriggerPrimitiveAlgo.cc.

References HcaluLUTTPGCoder::adc2Linear(), addSignal(), HcalDetId::depth(), HF_Veto, i, HFDataFrame::id(), incoder_, minSignalThreshold_, HFDataFrame::presamples(), DetId::rawId(), IntegerCaloSamples::setPresamples(), HFDataFrame::size(), theTowerMapFGSum, theTrigTowerGeometry, and HcalTrigTowerGeometry::towerIds().

                                                                  {

   if(frame.id().depth() == 1 || frame.id().depth() == 2) {
      std::vector<HcalTrigTowerDetId> ids = theTrigTowerGeometry.towerIds(frame.id());
      assert(ids.size() == 1);
      IntegerCaloSamples samples(ids[0], frame.size());
      samples.setPresamples(frame.presamples());
      incoder_->adc2Linear(frame, samples);

      // Don't add to final collection yet
      // HF PMT veto sum is calculated in analyzerHF()
      IntegerCaloSamples zero_samples(ids[0], frame.size());
      zero_samples.setPresamples(frame.presamples());
      addSignal(zero_samples);

      // Mask off depths: fgid is the same for both depths
      uint32_t fgid = (frame.id().rawId() | 0x1c000) ;

      if ( theTowerMapFGSum.find(ids[0]) == theTowerMapFGSum.end() ) {
         SumFGContainer sumFG;
         theTowerMapFGSum.insert(std::pair<HcalTrigTowerDetId, SumFGContainer >(ids[0], sumFG));
      }

      SumFGContainer& sumFG = theTowerMapFGSum[ids[0]];
      SumFGContainer::iterator sumFGItr;
      for ( sumFGItr = sumFG.begin(); sumFGItr != sumFG.end(); ++sumFGItr) {
         if (sumFGItr->id() == fgid) break;
      }
      // If find
      if (sumFGItr != sumFG.end()) {
         for (int i=0; i<samples.size(); ++i) (*sumFGItr)[i] += samples[i];
      }
      else {
         //Copy samples (change to fgid)
         IntegerCaloSamples sumFGSamples(DetId(fgid), samples.size());
         sumFGSamples.setPresamples(samples.presamples());
         for (int i=0; i<samples.size(); ++i) sumFGSamples[i] = samples[i];
         sumFG.push_back(sumFGSamples);
      }

      // set veto to true if Long or Short less than threshold
      if (HF_Veto.find(fgid) == HF_Veto.end()) {
         vector<bool> vetoBits(samples.size(), false);
         HF_Veto[fgid] = vetoBits;
      }
      for (int i=0; i<samples.size(); ++i)
         if (samples[i] < minSignalThreshold_)
            HF_Veto[fgid][i] = true;
   }
}
void HcalTriggerPrimitiveAlgo::addSignal ( const IntegerCaloSamples samples) [private]

Definition at line 163 of file HcalTriggerPrimitiveAlgo.cc.

References i, IntegerCaloSamples::id(), IntegerCaloSamples::size(), and theSumMap.

                                                                           {
   HcalTrigTowerDetId id(samples.id());
   SumMap::iterator itr = theSumMap.find(id);
   if(itr == theSumMap.end()) {
      theSumMap.insert(std::make_pair(id, samples));
   }
   else {
      // wish CaloSamples had a +=
      for(int i = 0; i < samples.size(); ++i) {
         (itr->second)[i] += samples[i];
      }
   }
}
void HcalTriggerPrimitiveAlgo::analyze ( IntegerCaloSamples samples,
HcalTriggerPrimitiveDigi result 
) [private]

adds the actual RecHits

Definition at line 178 of file HcalTriggerPrimitiveAlgo.cc.

References HcalTPGCompressor::compress(), fgMap_, i, IntegerCaloSamples::id(), UserOptions_cff::idx, numberOfPresamples_, numberOfSamples_, outcoder_, convertSQLitetoXML_cfg::output, peak_finder_algorithm_, peakfind_, IntegerCaloSamples::presamples(), edm::shift, IntegerCaloSamples::size(), theThreshold, and weights_.

Referenced by run().

                                                                                                      {
   int shrink = weights_.size() - 1;
   std::vector<bool> finegrain(numberOfSamples_,false);
   std::vector<bool>& msb = fgMap_[samples.id()];
   IntegerCaloSamples sum(samples.id(), samples.size());

   //slide algo window
   for(int ibin = 0; ibin < int(samples.size())- shrink; ++ibin) {
      int algosumvalue = 0;
      for(unsigned int i = 0; i < weights_.size(); i++) {
         //add up value * scale factor
         algosumvalue += int(samples[ibin+i] * weights_[i]);
      }
      if (algosumvalue<0) sum[ibin]=0;            // low-side
                                                  //high-side
      //else if (algosumvalue>0x3FF) sum[ibin]=0x3FF;
      else sum[ibin] = algosumvalue;              //assign value to sum[]
   }

   IntegerCaloSamples output(samples.id(), numberOfSamples_);
   output.setPresamples(numberOfPresamples_);

   // Align digis and TP
   int shift = samples.presamples() - numberOfPresamples_;
   if (peakfind_) {
      assert (shift >= (peakfind_ ? shrink : 0));
      assert(shift + numberOfSamples_ + shrink <= samples.size() - (peak_finder_algorithm_ - 1));
   }


   for (int ibin = 0; ibin < numberOfSamples_; ++ibin) {
      // ibin - index for output TP
      // idx - index for samples + shift
      int idx = ibin + shift;

      //Peak finding
      if (peakfind_) {
         bool isPeak = false;
         switch (peak_finder_algorithm_) {
            case 1 :
               isPeak = (samples[idx] > samples[idx-1] && samples[idx] >= samples[idx+1] && samples[idx] > theThreshold);
               break;
            case 2:
               isPeak = (sum[idx] > sum[idx-1] && sum[idx] >= sum[idx+1] && sum[idx] > theThreshold);
               break;
            default:
               break;
         }

         if (isPeak){
            output[ibin] = std::min<unsigned int>(sum[idx],0x3FF);
            finegrain[ibin] = msb[idx];
         }
         // Not a peak
         else output[ibin] = 0;
      }
      else { // No peak finding, just output running sum
         output[ibin] = std::min<unsigned int>(sum[idx],0x3FF);
         finegrain[ibin] = msb[idx];
      }

      // Only Pegged for 1-TS algo.
      if (peak_finder_algorithm_ == 1) {
         if (samples[idx] >= 0x3FF)
            output[ibin] = 0x3FF;
      }
      outcoder_->compress(output, finegrain, result);
   }
}
void HcalTriggerPrimitiveAlgo::analyzeHF ( IntegerCaloSamples samples,
HcalTriggerPrimitiveDigi result,
float  rctlsb 
) [private]

Definition at line 249 of file HcalTriggerPrimitiveAlgo.cc.

References HcalTPGCompressor::compress(), FG_threshold_, HF_Veto, IntegerCaloSamples::id(), UserOptions_cff::idx, numberOfPresamples_, numberOfSamples_, outcoder_, convertSQLitetoXML_cfg::output, PMT_NoiseThreshold_, IntegerCaloSamples::presamples(), edm::shift, IntegerCaloSamples::size(), and theTowerMapFGSum.

Referenced by run().

                                                                                                                      {
   std::vector<bool> finegrain(numberOfSamples_, false);
   HcalTrigTowerDetId detId(samples.id());

   // Align digis and TP
   int shift = samples.presamples() - numberOfPresamples_;
   assert(shift >=  0);
   assert((shift + numberOfSamples_) <=  samples.size());

   TowerMapFGSum::const_iterator tower2fg = theTowerMapFGSum.find(detId);
   assert(tower2fg != theTowerMapFGSum.end());

   const SumFGContainer& sumFG = tower2fg->second;
   // Loop over all L+S pairs that mapped from samples.id()
   // Note: 1 samples.id() = 6 x (L+S) without noZS
   for (SumFGContainer::const_iterator sumFGItr = sumFG.begin(); sumFGItr != sumFG.end(); ++sumFGItr) {
      const std::vector<bool>& veto = HF_Veto[sumFGItr->id().rawId()];
      for (int ibin = 0; ibin < numberOfSamples_; ++ibin) {
         int idx = ibin + shift;
         // if not vetod, add L+S to total sum and calculate FG
         if (!(veto[idx] && (*sumFGItr)[idx] > PMT_NoiseThreshold_)) {
            samples[idx] += (*sumFGItr)[idx];
            finegrain[ibin] = (finegrain[ibin] || (*sumFGItr)[idx] >= FG_threshold_);
         }
      }
   }

   IntegerCaloSamples output(samples.id(), numberOfSamples_);
   output.setPresamples(numberOfPresamples_);

   for (int ibin = 0; ibin < numberOfSamples_; ++ibin) {
      int idx = ibin + shift;
      output[ibin] = samples[idx] / (rctlsb == 0.25 ? 4 : 8);
      if (output[ibin] > 0x3FF) output[ibin] = 0x3FF;
   }
   outcoder_->compress(output, finegrain, result);
}
void HcalTriggerPrimitiveAlgo::run ( const HcalTPGCoder incoder,
const HcalTPGCompressor outcoder,
const HBHEDigiCollection hbheDigis,
const HFDigiCollection hfDigis,
HcalTrigPrimDigiCollection result,
float  rctlsb 
)

Definition at line 42 of file HcalTriggerPrimitiveAlgo.cc.

References addSignal(), analyze(), analyzeHF(), edm::SortedCollection< T, SORT >::back(), edm::SortedCollection< T, SORT >::begin(), edm::SortedCollection< T, SORT >::end(), fgMap_, HcalTrigTowerGeometry::firstHFTower(), HF_Veto, incoder_, outcoder_, edm::SortedCollection< T, SORT >::push_back(), theSumMap, theTowerMapFGSum, and theTrigTowerGeometry.

Referenced by HcalTrigPrimDigiProducer::produce().

                                                 {

   incoder_=dynamic_cast<const HcaluLUTTPGCoder*>(incoder);
   outcoder_=outcoder;

   theSumMap.clear();
   theTowerMapFGSum.clear();
   HF_Veto.clear();
   fgMap_.clear();

   // do the HB/HE digis
   for(HBHEDigiCollection::const_iterator hbheItr = hbheDigis.begin();
   hbheItr != hbheDigis.end(); ++hbheItr) {
      addSignal(*hbheItr);
   }

   // and the HF digis
   for(HFDigiCollection::const_iterator hfItr = hfDigis.begin();
   hfItr != hfDigis.end(); ++hfItr) {
      addSignal(*hfItr);

   }

   for(SumMap::iterator mapItr = theSumMap.begin(); mapItr != theSumMap.end(); ++mapItr) {
      result.push_back(HcalTriggerPrimitiveDigi(mapItr->first));
      HcalTrigTowerDetId detId(mapItr->second.id());
      if(detId.ietaAbs() >= theTrigTowerGeometry.firstHFTower())
         { analyzeHF(mapItr->second, result.back(), rctlsb);}
         else{analyze(mapItr->second, result.back());}
   }
   return;
}
void HcalTriggerPrimitiveAlgo::runFEFormatError ( const FEDRawDataCollection rawraw,
const HcalElectronicsMap emap,
HcalTrigPrimDigiCollection result 
)

Definition at line 301 of file HcalTriggerPrimitiveAlgo.cc.

References edm::SortedCollection< T, SORT >::begin(), HcalHTRData::check(), FEDRawData::data(), DetId::det(), edm::SortedCollection< T, SORT >::end(), FEDRawDataCollection::FEDData(), HcalHTRData::getErrorsWord(), HcalDCCHeader::getSourceId(), HcalDCCHeader::getSpigotData(), HcalDCCHeader::getSpigotPresent(), HcalBarrel, HcalEndcap, HcalForward, HcalHTRData::htrSlot(), HcalHTRData::htrTopBottom(), i, HcalHTRData::isHistogramEvent(), HcalElectronicsMap::lookup(), FEDNumbering::MAXHCALFEDID, FEDNumbering::MINHCALFEDID, DetId::null(), HcalHTRData::readoutVMECrateId(), HcalElectronicsId::setHTR(), FEDRawData::size(), HcalDCCHeader::SPIGOT_COUNT, DetId::subdetId(), theTrigTowerGeometry, HcalTrigTowerGeometry::towerIds(), and TrackValidation_HighPurity_cff::valid.

Referenced by HcalTrigPrimDigiProducer::produce().

                                                 {
  std::set<uint32_t> FrontEndErrors;

  for(int i=FEDNumbering::MINHCALFEDID; i<=FEDNumbering::MAXHCALFEDID; ++i) {
    const FEDRawData& raw = rawraw->FEDData(i);
    if (raw.size()<12) continue;
    const HcalDCCHeader* dccHeader=(const HcalDCCHeader*)(raw.data());
    if(!dccHeader) continue;
    HcalHTRData htr;
    for (int spigot=0; spigot<HcalDCCHeader::SPIGOT_COUNT; spigot++) {
      if (!dccHeader->getSpigotPresent(spigot)) continue;
      dccHeader->getSpigotData(spigot,htr,raw.size());
      int dccid = dccHeader->getSourceId();
      int errWord = htr.getErrorsWord() & 0x1FFFF;
      bool HTRError = (!htr.check() || htr.isHistogramEvent() || (errWord & 0x800)!=0);

      if(HTRError) {
        bool valid =false;
        for(int fchan=0; fchan<3 && !valid; fchan++) {
          for(int fib=0; fib<9 && !valid; fib++) {
            HcalElectronicsId eid(fchan,fib,spigot,dccid-FEDNumbering::MINHCALFEDID);
            eid.setHTR(htr.readoutVMECrateId(),htr.htrSlot(),htr.htrTopBottom());
            DetId detId = emap->lookup(eid);
            if(detId.null()) continue;
            HcalSubdetector subdet=(HcalSubdetector(detId.subdetId()));
            if (detId.det()!=4||
              (subdet!=HcalBarrel && subdet!=HcalEndcap &&
              subdet!=HcalForward )) continue;
            std::vector<HcalTrigTowerDetId> ids = theTrigTowerGeometry.towerIds(detId);
            for (std::vector<HcalTrigTowerDetId>::const_iterator triggerId=ids.begin(); triggerId != ids.end(); ++triggerId) {
              FrontEndErrors.insert(triggerId->rawId());
            }
            //valid = true;
          }
        }
      }
    }
  }

  // Loop over TP collection
  // Set TP to zero if there is FE Format Error
  HcalTriggerPrimitiveSample zeroSample(0);
  for (HcalTrigPrimDigiCollection::iterator tp = result.begin(); tp != result.end(); ++tp){
    if (FrontEndErrors.find(tp->id().rawId()) != FrontEndErrors.end()) {
      for (int i=0; i<tp->size(); ++i) tp->setSample(i, zeroSample);
    }
  }
}
void HcalTriggerPrimitiveAlgo::runZS ( HcalTrigPrimDigiCollection tp)

Definition at line 287 of file HcalTriggerPrimitiveAlgo.cc.

References edm::SortedCollection< T, SORT >::begin(), edm::SortedCollection< T, SORT >::end(), i, and ZS_threshold_I_.

Referenced by HcalTrigPrimDigiProducer::produce().

                                                                       {
   for (HcalTrigPrimDigiCollection::iterator tp = result.begin(); tp != result.end(); ++tp){
      bool ZS = true;
      for (int i=0; i<tp->size(); ++i) {
         if (tp->sample(i).compressedEt()  > ZS_threshold_I_) {
            ZS=false;
            break;
         }
      }
      if (ZS) tp->setZSInfo(false,true);
      else tp->setZSInfo(true,false);
   }
}
void HcalTriggerPrimitiveAlgo::setPeakFinderAlgorithm ( int  algo)

Definition at line 363 of file HcalTriggerPrimitiveAlgo.cc.

References Exception, and peak_finder_algorithm_.

Referenced by HcalTrigPrimDigiProducer::HcalTrigPrimDigiProducer().

                                                             {
   if (algo <=0 && algo>2)
      throw cms::Exception("ERROR: Only algo 1 & 2 are supported.") << std::endl;
   peak_finder_algorithm_ = algo;
}

Member Data Documentation

Definition at line 60 of file HcalTriggerPrimitiveAlgo.h.

Referenced by analyzeHF().

Definition at line 99 of file HcalTriggerPrimitiveAlgo.h.

Referenced by addFG(), analyze(), and run().

Definition at line 96 of file HcalTriggerPrimitiveAlgo.h.

Referenced by addSignal(), analyzeHF(), and run().

Definition at line 54 of file HcalTriggerPrimitiveAlgo.h.

Referenced by addSignal(), and run().

Definition at line 59 of file HcalTriggerPrimitiveAlgo.h.

Definition at line 65 of file HcalTriggerPrimitiveAlgo.h.

Referenced by addSignal().

Definition at line 64 of file HcalTriggerPrimitiveAlgo.h.

Referenced by analyze(), analyzeHF(), and HcalTriggerPrimitiveAlgo().

Definition at line 63 of file HcalTriggerPrimitiveAlgo.h.

Referenced by analyze(), analyzeHF(), and HcalTriggerPrimitiveAlgo().

Definition at line 55 of file HcalTriggerPrimitiveAlgo.h.

Referenced by analyze(), analyzeHF(), and run().

Definition at line 72 of file HcalTriggerPrimitiveAlgo.h.

Referenced by analyze(), and setPeakFinderAlgorithm().

Definition at line 57 of file HcalTriggerPrimitiveAlgo.h.

Referenced by analyze(), and HcalTriggerPrimitiveAlgo().

Definition at line 66 of file HcalTriggerPrimitiveAlgo.h.

Referenced by analyzeHF().

Definition at line 80 of file HcalTriggerPrimitiveAlgo.h.

Referenced by addSignal(), and run().

Definition at line 56 of file HcalTriggerPrimitiveAlgo.h.

Referenced by analyze().

Definition at line 84 of file HcalTriggerPrimitiveAlgo.h.

Referenced by addSignal(), analyzeHF(), and run().

Definition at line 77 of file HcalTriggerPrimitiveAlgo.h.

Referenced by addSignal(), run(), and runFEFormatError().

std::vector<double> HcalTriggerPrimitiveAlgo::weights_ [private]

Definition at line 58 of file HcalTriggerPrimitiveAlgo.h.

Referenced by analyze().

Definition at line 61 of file HcalTriggerPrimitiveAlgo.h.

Referenced by HcalTriggerPrimitiveAlgo().

Definition at line 62 of file HcalTriggerPrimitiveAlgo.h.

Referenced by HcalTriggerPrimitiveAlgo(), and runZS().