test
CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
EcalSimpleProducer.cc
Go to the documentation of this file.
2 #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()!=0){
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_) cout << "(" << iEta0 << "," << iPhi0 << "): ";
29  digis->push_back(EBDetId(iEta1,iPhi));
30  DataFrame dframe(digis->back());
31 
32  for(int t = 0; t < nSamples; ++t){
33  uint16_t encodedAdc = (uint16_t)formula_->Eval(iEta0, iPhi0, ievt-1, t);
34  if(verbose_) cout << encodedAdc << ((t<(nSamples-1))?"\t":"\n");
35  dframe[t]=encodedAdc;
36  }
37  }
38  }
39  evt.put(std::move(digis));
40  //puts an empty digi collecion for endcap:
41  evt.put(unique_ptr<EEDigiCollection>(new EEDigiCollection()));
42  }
43  if(tpFormula_.get()!=0){
44  unique_ptr<EcalTrigPrimDigiCollection> tps
45  = unique_ptr<EcalTrigPrimDigiCollection>(new EcalTrigPrimDigiCollection);
46  tps->reserve(56*72);
47  const int nSamples = 5;
48  for(int iTtEta0=0; iTtEta0<56; ++iTtEta0){
49  for(int iTtPhi0=0; iTtPhi0<72; ++iTtPhi0){
50  int iTtEta1 = cIndex2iTtEta(iTtEta0);
51  int iTtPhi = cIndex2iTtPhi(iTtPhi0);
52 
53  if(verbose_) cout << "(" << iTtEta0 << "," << iTtPhi0 << "): ";
54  int zside = iTtEta1<0?-1:1;
56  tpframe(EcalTrigTowerDetId(zside, EcalTriggerTower,
57  abs(iTtEta1), iTtPhi));
58 
59  tpframe.setSize(nSamples);
60 
61  if(verbose_) cout << "TP: ";
62  for(int t = 0; t < nSamples; ++t){
63  uint16_t encodedTp = (uint16_t)tpFormula_->Eval(iTtEta0, iTtPhi0, ievt-1, t);
64 
65  if(verbose_) cout << "TP(" << iTtEta0 << "," << iTtPhi0 << ") = "
66  << encodedTp << ((t<(nSamples-1))?"\t":"\n");
67  tpframe.setSample(t, EcalTriggerPrimitiveSample(encodedTp));
68  }
69  tps->push_back(tpframe);
70  }
71  }
72  evt.put(std::move(tps));
73  }
74  if(simHitFormula_.get()!=0){//generation of barrel sim hits
75  unique_ptr<PCaloHitContainer> hits
76  = 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_) cout << "(" << iEta0 << "," << iPhi0 << "): ";
82 
83  double em = simHitFormula_->Eval(iEta0, iPhi0, ievt-1);
84  double eh = 0.;
85  double t = 0.;
86  const PCaloHit hit(EBDetId(iEta1,iPhi).rawId(), em, eh, t, 0);
87  hits->push_back(hit);
88  }
89  }
90  evt.put(std::move(hits), "EcalHitsEB");
91  //puts an empty digi collecion for endcap:
92  evt.put(unique_ptr<PCaloHitContainer>(new PCaloHitContainer()),
93  "EcalHitsEE");
94  }
95 }
96 
98  EDProducer(){
99  string formula = pset.getParameter<string>("formula");
100  string tpFormula = pset.getParameter<string>("tpFormula");
101  string simHitFormula = pset.getParameter<string>("simHitFormula");
102 
103  verbose_ = pset.getUntrackedParameter<bool>("verbose", false);
104  // replaceAll(formula, "itt0", "((((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 
127  //replaceAll(simHitormula, "itt0", "((((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.size()!=0){
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.size()!=0){
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.size()!=0){
153  = unique_ptr<TFormula>(new TFormula("f", simHitFormula.c_str()));
154  Int_t err = simHitFormula_->Compile();
155  if(err!=0){
156  throw cms::Exception("Error in EcalSimpleProducer "
157  "'simHitFormula' config.");
158  }
159  produces<edm::PCaloHitContainer>("EcalHitsEB");
160  produces<edm::PCaloHitContainer>("EcalHitsEE");
161  }
162 }
163 
165  const std::string& to) const{
166  string::size_type pos = 0;
167  // cout << "replaceAll(" << s << "," << from << "," << to << ")\n";
168  while((pos=s.find(from, pos))!=string::npos){
169  // cout << "replace(" << pos << "," << from.size() << "," << to << ")\n";
170  s.replace(pos, from.size(), to);
171  // cout << " -> " << s << "\n";
172  }
173 }
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:122
std::unique_ptr< TFormula > simHitFormula_
std::unique_ptr< TFormula > tpFormula_
virtual void produce(edm::Event &evt, const edm::EventSetup &)
EcalSimpleProducer(const edm::ParameterSet &pset)
int zside(DetId const &)
uint16_t size_type
void setSample(int i, const EcalTriggerPrimitiveSample &sam)
def move
Definition: eostools.py:510
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:58
tuple cout
Definition: gather_cfg.py:145
void replaceAll(std::string &s, const std::string &from, const std::string &to) const