ViewVC Help
View File | Revision Log | Show Annotations | Root Listing
root/cvsroot/UserCode/ShallowTools/plugins/ShallowDigisProducer.cc
Revision: 1.2
Committed: Tue Sep 29 10:16:31 2009 UTC (15 years, 7 months ago) by bbetchar
Content type: text/plain
Branch: MAIN
CVS Tags: HEAD
Changes since 1.1: +3 -0 lines
Log Message:
updates for CMSSW_3_3_X

File Contents

# User Rev Content
1 bbetchar 1.1 #include "UserCode/ShallowTools/interface/ShallowDigisProducer.h"
2    
3 bbetchar 1.2 #include "FWCore/Framework/interface/Event.h"
4     #include "FWCore/Framework/interface/EventSetup.h"
5     #include "FWCore/MessageLogger/interface/MessageLogger.h"
6 bbetchar 1.1 #include "DataFormats/Common/interface/DetSetVector.h"
7     #include "DataFormats/Common/interface/DetSetVectorNew.h"
8     #include "DataFormats/SiStripDigi/interface/SiStripDigi.h"
9     #include "CondFormats/DataRecord/interface/SiStripNoisesRcd.h"
10     #include "CondFormats/SiStripObjects/interface/SiStripNoises.h"
11     #include "boost/foreach.hpp"
12    
13     ShallowDigisProducer::ShallowDigisProducer(const edm::ParameterSet& conf)
14     : inputTags(conf.getParameter<std::vector<edm::InputTag> >("DigiProducersList"))
15     {
16     produces<std::vector<unsigned> >("id");
17     produces<std::vector<unsigned> >("subdet");
18     produces<std::vector<unsigned> >("strip");
19     produces<std::vector<unsigned> >("adc");
20     produces<std::vector<float> >("noise");
21     }
22    
23     void ShallowDigisProducer::
24     insert(products& p, edm::Event& e) {
25     e.put(p.id, "id");
26     e.put(p.subdet, "subdet");
27     e.put(p.strip, "strip");
28     e.put(p.adc, "adc");
29     e.put(p.noise, "noise");
30     }
31    
32     template<class T>
33     inline
34     void ShallowDigisProducer::
35     recordDigis(const T& digiCollection, products& p) {
36     BOOST_FOREACH(const typename T::value_type set, digiCollection) {
37     SiStripNoises::Range detNoiseRange = noiseHandle->getRange(set.detId());
38     BOOST_FOREACH(const SiStripDigi digi, set) {
39     p.id->push_back(set.detId());
40     p.subdet->push_back((set.detId()>>25)&0x7);
41     p.strip->push_back(digi.strip());
42     p.adc->push_back(digi.adc());
43     p.noise->push_back(noiseHandle->getNoise( digi.strip(), detNoiseRange));
44     }
45     }
46     }
47    
48     void ShallowDigisProducer::
49     produce(edm::Event& e, const edm::EventSetup& es) {
50     products p;
51     edm::Handle< edm::DetSetVector<SiStripDigi> > inputOld;
52     edm::Handle< edmNew::DetSetVector<SiStripDigi> > inputNew;
53     es.get<SiStripNoisesRcd>().get(noiseHandle);
54     if( findInput(inputOld, e) ) recordDigis(*inputOld, p); else
55     if( findInput(inputNew, e) ) recordDigis(*inputNew, p); else
56     edm::LogWarning("Input Not Found");
57     insert(p,e);
58     }
59    
60     template<class T>
61     inline
62     bool ShallowDigisProducer::
63     findInput(edm::Handle<T>& handle, const edm::Event& e) {
64     BOOST_FOREACH( const edm::InputTag inputTag, inputTags) {
65     e.getByLabel(inputTag, handle);
66     if( handle.isValid() && !handle->empty() ) {
67     LogDebug("Input") << inputTag;
68     return true;
69     }
70     }
71     return false;
72     }