CMS 3D CMS Logo

CTPPSPixelClusterProducer.cc
Go to the documentation of this file.
1 
4 
6  param_(conf) ,
7  clusterizer_(conf){
8 
9  src_ = conf.getParameter<std::string>("label");
10  verbosity_ = conf.getUntrackedParameter<int> ("RPixVerbosity");
11 
12  tokenCTPPSPixelDigi_ = consumes<edm::DetSetVector<CTPPSPixelDigi> >(edm::InputTag(src_));
13 
14  produces<edm::DetSetVector<CTPPSPixelCluster> > ();
15  }
16 
18 
19 }
20 
23  desc.addUntracked<int>("RPixVerbosity",0);
24  desc.add<std::string>("label", "ctppsPixelDigis");
25  desc.add<int>("SeedADCThreshold",2);
26  desc.add<int>("ADCThreshold",2);
27  desc.add<double>("ElectronADCGain",135.0);
28  desc.add<int>("VCaltoElectronGain",50);
29  desc.add<int>("VCaltoElectronOffset",-411);
30  desc.add<bool>("doSingleCalibration",false);
31  descriptions.add("ctppsPixelClusters", desc);
32 }
33 
35 
38  iEvent.getByToken(tokenCTPPSPixelDigi_, rpd);
39 
40 // get analysis mask to mask channels
42 
43  if(!rpd->empty())
44  iSetup.get<CTPPSPixelAnalysisMaskRcd>().get(aMask);
45 
47 
48 // run clusterisation
49  if (!rpd->empty()){
50 // get calibration DB
51  theGainCalibrationDB.getDB(iEvent,iSetup);
52  run(*rpd, output, aMask.product());
53  }
54 // write output
55  iEvent.put(std::make_unique<edm::DetSetVector<CTPPSPixelCluster> >(output));
56 
57 }
58 
61 
62  for (const auto &ds_digi : input)
63  {
64  edm::DetSet<CTPPSPixelCluster> &ds_cluster = output.find_or_insert(ds_digi.id);
65  clusterizer_.buildClusters(ds_digi.id, ds_digi.data, ds_cluster.data, theGainCalibrationDB.getCalibs(), mask);
66 
67  if(verbosity_){
68  unsigned int cluN=0;
69  for(std::vector<CTPPSPixelCluster>::iterator iit = ds_cluster.data.begin(); iit != ds_cluster.data.end(); iit++){
70  edm::LogInfo("CTPPSPixelClusterProducer") << "Cluster " << ++cluN <<" avg row "
71  << (*iit).avg_row()<< " avg col " << (*iit).avg_col()<<" ADC.size " << (*iit).size();
72  }
73  }
74  }
75 }
76 
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:136
ParameterDescriptionBase * addUntracked(U const &iLabel, T const &value)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:519
void run(const edm::DetSetVector< CTPPSPixelDigi > &input, edm::DetSetVector< CTPPSPixelCluster > &output, const CTPPSPixelAnalysisMask *mask)
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
edm::EDGetTokenT< edm::DetSetVector< CTPPSPixelDigi > > tokenCTPPSPixelDigi_
const CTPPSPixelGainCalibrations * getCalibs() const
static std::string const input
Definition: EdmProvDump.cc:44
reference find_or_insert(det_id_type id)
Definition: DetSetVector.h:254
int iEvent
Definition: GenABIO.cc:230
struct @608 param_
void get(HolderT &iHolder) const
Channel-mask mapping.
ParameterDescriptionBase * add(U const &iLabel, T const &value)
virtual void getDB(const edm::Event &e, const edm::EventSetup &c)
CTPPSPixelGainCalibrationDBService theGainCalibrationDB
void produce(edm::Event &, const edm::EventSetup &) override
const T & get() const
Definition: EventSetup.h:59
void add(std::string const &label, ParameterSetDescription const &psetDescription)
collection_type data
Definition: DetSet.h:78
void buildClusters(unsigned int detId, const std::vector< CTPPSPixelDigi > &digi, std::vector< CTPPSPixelCluster > &clusters, const CTPPSPixelGainCalibrations *pcalibration, const CTPPSPixelAnalysisMask *mask)
CTPPSPixelClusterProducer(const edm::ParameterSet &param)
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
T const * product() const
Definition: ESHandle.h:86