CMS 3D CMS Logo

EcalSimpleProducer.cc
Go to the documentation of this file.
9 #include "TFormula.h"
10 #include <iostream>
11 #include <string>
12 
13 using namespace std;
14 using namespace edm;
15 
17  const int ievt = evt.id().event();
18  if (formula_.get() != nullptr) {
19  unique_ptr<EBDigiCollection> digis(new EBDigiCollection);
20 
21  digis->reserve(170 * 360);
22 
23  const int nSamples = digis->stride();
24  for (int iEta0 = 0; iEta0 < 170; ++iEta0) {
25  for (int iPhi0 = 0; iPhi0 < 360; ++iPhi0) {
26  int iEta1 = cIndex2iEta(iEta0);
27  int iPhi = cIndex2iPhi(iPhi0);
28  if (verbose_)
29  cout << "(" << iEta0 << "," << iPhi0 << "): ";
30  digis->push_back(EBDetId(iEta1, iPhi));
31  DataFrame dframe(digis->back());
32 
33  for (int t = 0; t < nSamples; ++t) {
34  uint16_t encodedAdc = (uint16_t)formula_->Eval(iEta0, iPhi0, ievt - 1, t);
35  if (verbose_)
36  cout << encodedAdc << ((t < (nSamples - 1)) ? "\t" : "\n");
37  dframe[t] = encodedAdc;
38  }
39  }
40  }
41  evt.put(std::move(digis));
42  // puts an empty digi collecion for endcap:
43  evt.put(unique_ptr<EEDigiCollection>(new EEDigiCollection()));
44  }
45  if (tpFormula_.get() != nullptr) {
46  unique_ptr<EcalTrigPrimDigiCollection> tps = unique_ptr<EcalTrigPrimDigiCollection>(new EcalTrigPrimDigiCollection);
47  tps->reserve(56 * 72);
48  const int nSamples = 5;
49  for (int iTtEta0 = 0; iTtEta0 < 56; ++iTtEta0) {
50  for (int iTtPhi0 = 0; iTtPhi0 < 72; ++iTtPhi0) {
51  int iTtEta1 = cIndex2iTtEta(iTtEta0);
52  int iTtPhi = cIndex2iTtPhi(iTtPhi0);
53 
54  if (verbose_)
55  cout << "(" << iTtEta0 << "," << iTtPhi0 << "): ";
56  int zside = iTtEta1 < 0 ? -1 : 1;
57  EcalTriggerPrimitiveDigi tpframe(EcalTrigTowerDetId(zside, EcalTriggerTower, abs(iTtEta1), iTtPhi));
58 
59  tpframe.setSize(nSamples);
60 
61  if (verbose_)
62  cout << "TP: ";
63  for (int t = 0; t < nSamples; ++t) {
64  uint16_t encodedTp = (uint16_t)tpFormula_->Eval(iTtEta0, iTtPhi0, ievt - 1, t);
65 
66  if (verbose_)
67  cout << "TP(" << iTtEta0 << "," << iTtPhi0 << ") = " << encodedTp << ((t < (nSamples - 1)) ? "\t" : "\n");
68  tpframe.setSample(t, EcalTriggerPrimitiveSample(encodedTp));
69  }
70  tps->push_back(tpframe);
71  }
72  }
73  evt.put(std::move(tps));
74  }
75  if (simHitFormula_.get() != nullptr) { // generation of barrel sim hits
76  unique_ptr<PCaloHitContainer> hits = unique_ptr<PCaloHitContainer>(new PCaloHitContainer);
77  for (int iEta0 = 0; iEta0 < 170; ++iEta0) {
78  for (int iPhi0 = 0; iPhi0 < 360; ++iPhi0) {
79  int iEta1 = cIndex2iEta(iEta0);
80  int iPhi = cIndex2iPhi(iPhi0);
81  if (verbose_)
82  cout << "(" << iEta0 << "," << iPhi0 << "): ";
83 
84  double em = simHitFormula_->Eval(iEta0, iPhi0, ievt - 1);
85  double eh = 0.;
86  double t = 0.;
87  const PCaloHit hit(EBDetId(iEta1, iPhi).rawId(), em, eh, t, 0);
88  hits->push_back(hit);
89  }
90  }
91  evt.put(std::move(hits), "EcalHitsEB");
92  // puts an empty digi collecion for endcap:
93  evt.put(unique_ptr<PCaloHitContainer>(new PCaloHitContainer()), "EcalHitsEE");
94  }
95 }
96 
98  string formula = pset.getParameter<string>("formula");
99  string tpFormula = pset.getParameter<string>("tpFormula");
100  string simHitFormula = pset.getParameter<string>("simHitFormula");
101 
102  verbose_ = pset.getUntrackedParameter<bool>("verbose", false);
103  // replaceAll(formula, "itt0",
104  // "((((ieta0<85)*(84-ieta0)+(ieta0>=85)*(ieta0-85))/5-18)*4+((iphi0/5+2)%4))");
105  replaceAll(formula, "ebm", "(ieta0<85)");
106  replaceAll(formula, "ebp", "(ieta0>84)");
107  replaceAll(formula, "ieta0", "x");
108  replaceAll(formula, "iphi0", "y");
109  replaceAll(formula, "ievt0", "z");
110  replaceAll(formula, "isample0", "t");
111  // cout << "----------> " << formula << endl;
112 
113  replaceAll(tpFormula, "itt0", "((ieta0<28)*(27-ieta0)+(ieta0>=28)*(ieta0-28))*4+(iphi0+2)%4");
114  replaceAll(tpFormula, "eb", "(ieta0>10 && ieta0<45)");
115  replaceAll(tpFormula, "ebm", "(ieta0>10 && ieta0<28)");
116  replaceAll(tpFormula, "ebp", "(ieta0>27 && ieta0<45)");
117  replaceAll(tpFormula, "ee", "(ieta0<11 || ieta0>44)");
118  replaceAll(tpFormula, "eem", "(ieta0<11)");
119  replaceAll(tpFormula, "eep", "(ieta0>44)");
120  replaceAll(tpFormula, "ieta0", "x");
121  replaceAll(tpFormula, "iphi0", "y");
122  replaceAll(tpFormula, "ievt0", "z");
123  replaceAll(tpFormula, "isample0", "t");
124  // cout << "----------> " << tpFormula << endl;
125 
126  // replaceAll(simHitormula, "itt0",
127  // "((((ieta0<85)*(84-ieta0)+(ieta0>=85)*(ieta0-85))/5-18)*4+((iphi0/5+2)%4))");
128  replaceAll(simHitFormula, "ebm", "(ieta0<85)");
129  replaceAll(simHitFormula, "ebp", "(ieta0>84)");
130  replaceAll(simHitFormula, "ieta0", "x");
131  replaceAll(simHitFormula, "iphi0", "y");
132  replaceAll(simHitFormula, "ievt0", "z");
133 
134  if (!formula.empty()) {
135  formula_ = unique_ptr<TFormula>(new TFormula("f", formula.c_str()));
136  Int_t err = formula_->Compile();
137  if (err != 0) {
138  throw cms::Exception("Error in EcalSimpleProducer 'formula' config.");
139  }
140  produces<EBDigiCollection>();
141  produces<EEDigiCollection>();
142  }
143  if (!tpFormula.empty()) {
144  tpFormula_ = unique_ptr<TFormula>(new TFormula("f", tpFormula.c_str()));
145  Int_t err = tpFormula_->Compile();
146  if (err != 0) {
147  throw cms::Exception("Error in EcalSimpleProducer 'tpFormula' config.");
148  }
149  produces<EcalTrigPrimDigiCollection>();
150  }
151  if (!simHitFormula.empty()) {
152  simHitFormula_ = unique_ptr<TFormula>(new TFormula("f", simHitFormula.c_str()));
153  Int_t err = simHitFormula_->Compile();
154  if (err != 0) {
155  throw cms::Exception(
156  "Error in EcalSimpleProducer "
157  "'simHitFormula' config.");
158  }
159  produces<edm::PCaloHitContainer>("EcalHitsEB");
160  produces<edm::PCaloHitContainer>("EcalHitsEE");
161  }
162 }
163 
166  // cout << "replaceAll(" << s << "," << from << "," << to << ")\n";
167  while ((pos = s.find(from, pos)) != string::npos) {
168  // cout << "replace(" << pos << "," << from.size() << "," << to << ")\n";
169  s.replace(pos, from.size(), to);
170  // cout << " -> " << s << "\n";
171  }
172 }
T getParameter(std::string const &) const
EventNumber_t event() const
Definition: EventID.h:41
T getUntrackedParameter(std::string const &, T const &) const
std::vector< PCaloHit > PCaloHitContainer
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:125
std::unique_ptr< TFormula > simHitFormula_
std::unique_ptr< TFormula > tpFormula_
EcalSimpleProducer(const edm::ParameterSet &pset)
int zside(DetId const &)
uint16_t size_type
void produce(edm::Event &evt, const edm::EventSetup &) override
void setSample(int i, const EcalTriggerPrimitiveSample &sam)
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
std::unique_ptr< TFormula > formula_
constexpr size_t nSamples
edm::SortedCollection< EcalTriggerPrimitiveDigi > EcalTrigPrimDigiCollection
edm::EventID id() const
Definition: EventBase.h:59
HLT enums.
void replaceAll(std::string &s, const std::string &from, const std::string &to) const
def move(src, dest)
Definition: eostools.py:511