CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
EcalDeadChannelRecoveryNN.cc
Go to the documentation of this file.
3 
5 
6 #include <iostream>
7 #include <TMath.h>
8 
9 template <typename T>
11  for (int id = 0; id < 9; ++id) {
12  ctx_[id].mlp = NULL;
13  }
14 
15  this->load();
16 }
17 
18 template <typename T>
20  for (int id = 0; id < 9; ++id) {
21  if (ctx_[id].mlp) {
22  // @TODO segfaults for an uknown reason
23  // delete ctx[id].mlp;
24  // delete ctx[id].tree;
25  }
26  }
27 }
28 
29 template <>
31 {
33 }
34 
35 template <>
37 {
39 }
40 
41 template <typename T>
44 
45  TTree *t = new TTree("t", "dummy MLP tree");
46  t->SetDirectory(0);
47 
48  t->Branch("z1", &(ctx.tmp[0]), "z1/D");
49  t->Branch("z2", &(ctx.tmp[1]), "z2/D");
50  t->Branch("z3", &(ctx.tmp[2]), "z3/D");
51  t->Branch("z4", &(ctx.tmp[3]), "z4/D");
52  t->Branch("z5", &(ctx.tmp[4]), "z5/D");
53  t->Branch("z6", &(ctx.tmp[5]), "z6/D");
54  t->Branch("z7", &(ctx.tmp[6]), "z7/D");
55  t->Branch("z8", &(ctx.tmp[7]), "z8/D");
56  t->Branch("zf", &(ctx.tmp[8]), "zf/D");
57 
58  ctx.tree = t;
59  ctx.mlp =
60  new TMultiLayerPerceptron("@z1,@z2,@z3,@z4,@z5,@z6,@z7,@z8:10:5:zf", t);
61  ctx.mlp->LoadWeights(path.c_str());
62 }
63 
64 template <>
66  std::string p = "RecoLocalCalo/EcalDeadChannelRecoveryAlgos/data/NNWeights/";
67 
68  this->load_file(ctx_[CellID::CC], p + "EB_ccNNWeights.txt");
69  this->load_file(ctx_[CellID::RR], p + "EB_rrNNWeights.txt");
70  this->load_file(ctx_[CellID::LL], p + "EB_llNNWeights.txt");
71  this->load_file(ctx_[CellID::UU], p + "EB_uuNNWeights.txt");
72  this->load_file(ctx_[CellID::DD], p + "EB_ddNNWeights.txt");
73  this->load_file(ctx_[CellID::RU], p + "EB_ruNNWeights.txt");
74  this->load_file(ctx_[CellID::RD], p + "EB_rdNNWeights.txt");
75  this->load_file(ctx_[CellID::LU], p + "EB_luNNWeights.txt");
76  this->load_file(ctx_[CellID::LD], p + "EB_ldNNWeights.txt");
77 }
78 
79 template <>
81  std::string p = "RecoLocalCalo/EcalDeadChannelRecoveryAlgos/data/NNWeights/";
82 
83  this->load_file(ctx_[CellID::CC], p + "EE_ccNNWeights.txt");
84  this->load_file(ctx_[CellID::RR], p + "EE_rrNNWeights.txt");
85  this->load_file(ctx_[CellID::LL], p + "EE_llNNWeights.txt");
86  this->load_file(ctx_[CellID::UU], p + "EE_uuNNWeights.txt");
87  this->load_file(ctx_[CellID::DD], p + "EE_ddNNWeights.txt");
88  this->load_file(ctx_[CellID::RU], p + "EE_ruNNWeights.txt");
89  this->load_file(ctx_[CellID::RD], p + "EE_rdNNWeights.txt");
90  this->load_file(ctx_[CellID::LU], p + "EE_luNNWeights.txt");
91  this->load_file(ctx_[CellID::LD], p + "EE_ldNNWeights.txt");
92 }
93 
94 template <typename T>
95 double EcalDeadChannelRecoveryNN<T>::recover(const T id, const EcalRecHitCollection& hit_collection, double Sum8Cut, bool *AcceptFlag)
96 {
97  // use the correct probability matrix
98  typedef class EcalCrystalMatrixProbality<T> P;
99 
100  double NewEnergy = 0.0;
101 
102  double NewEnergy_RelMC = 0.0;
103  double NewEnergy_RelDC = 0.0;
104 
105  double MNxN_RelMC[9] = { 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0 } ;
106  double MNxN_RelDC[9] = { 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0 } ;
107 
108  double sum8 = 0.0;
109 
110  double sum8_RelMC = makeNxNMatrice_RelMC(id, hit_collection, MNxN_RelMC,AcceptFlag);
111  double sum8_RelDC = makeNxNMatrice_RelDC(id, hit_collection, MNxN_RelDC,AcceptFlag);
112 
113  // Only if "AcceptFlag" is true call the ANN
114  if ( *AcceptFlag ) {
115  if (sum8_RelDC > Sum8Cut && sum8_RelMC > Sum8Cut) {
116  NewEnergy_RelMC = estimateEnergy(MNxN_RelMC);
117  NewEnergy_RelDC = estimateEnergy(MNxN_RelDC);
118 
119  // Matrices "MNxN_RelMC" and "MNxN_RelDC" have now the full set of energies, the original ones plus
120  // whatever "estimates" by the ANN for the "dead" xtal. Use those full matrices and calculate probabilities.
121  //
122  double SumMNxN_RelMC = MNxN_RelMC[LU] + MNxN_RelMC[UU] + MNxN_RelMC[RU] +
123  MNxN_RelMC[LL] + MNxN_RelMC[CC] + MNxN_RelMC[RR] +
124  MNxN_RelMC[LD] + MNxN_RelMC[DD] + MNxN_RelMC[RD] ;
125 
126  double frMNxN_RelMC[9]; for (int i=0; i<9; i++) { frMNxN_RelMC[i] = MNxN_RelMC[i] / SumMNxN_RelMC ; }
127 
128  double prMNxN_RelMC = P::Diagonal( frMNxN_RelMC[LU] ) * P::UpDown( frMNxN_RelMC[UU] ) * P::Diagonal( frMNxN_RelMC[RU] ) *
129  P::ReftRight( frMNxN_RelMC[LL] ) * P::Central( frMNxN_RelMC[CC] ) * P::ReftRight( frMNxN_RelMC[RR] ) *
130  P::Diagonal( frMNxN_RelMC[LD] ) * P::UpDown( frMNxN_RelMC[DD] ) * P::Diagonal( frMNxN_RelMC[RD] ) ;
131 
132  double SumMNxN_RelDC = MNxN_RelDC[LU] + MNxN_RelDC[UU] + MNxN_RelDC[RU] +
133  MNxN_RelDC[LL] + MNxN_RelDC[CC] + MNxN_RelDC[RR] +
134  MNxN_RelDC[LD] + MNxN_RelDC[DD] + MNxN_RelDC[RD] ;
135 
136  double frMNxN_RelDC[9]; for (int i=0; i<9; i++) { frMNxN_RelDC[i] = MNxN_RelDC[i] / SumMNxN_RelDC ; }
137 
138  double prMNxN_RelDC = P::Diagonal( frMNxN_RelDC[LU] ) * P::UpDown( frMNxN_RelDC[UU] ) * P::Diagonal( frMNxN_RelDC[RU] ) *
139  P::ReftRight( frMNxN_RelDC[LL] ) * P::Central( frMNxN_RelDC[CC] ) * P::ReftRight( frMNxN_RelDC[RR] ) *
140  P::Diagonal( frMNxN_RelDC[LD] ) * P::UpDown( frMNxN_RelDC[DD] ) * P::Diagonal( frMNxN_RelDC[RD] ) ;
141 
142  if ( prMNxN_RelDC > prMNxN_RelMC ) { NewEnergy = NewEnergy_RelDC ; sum8 = sum8_RelDC ; }
143  if ( prMNxN_RelDC <= prMNxN_RelMC ) { NewEnergy = NewEnergy_RelMC ; sum8 = sum8_RelMC ; }
144 
145 
146  // If the return value of "CorrectDeadChannelsNN" is one of the followin negative values then
147  // it corresponds to an error condition. See "CorrectDeadChannelsNN.cc" for possible values.
148  if ( NewEnergy == -1000000.0 ||
149  NewEnergy == -1000001.0 ||
150  NewEnergy == -2000000.0 ||
151  NewEnergy == -3000000.0 ||
152  NewEnergy == -3000001.0 ) { *AcceptFlag=false ; NewEnergy = 0.0 ; }
153  }
154  }
155 
156  // Protect against non physical high values
157  // From the distribution of (max.cont.xtal / Sum8) we get as limit 5 (hard) and 10 (softer)
158  // Choose 10 as highest possible energy to be assigned to the dead channel under any scenario.
159  if ( NewEnergy > 10.0 * sum8 ) { *AcceptFlag=false ; NewEnergy = 0.0 ; }
160 
161  return NewEnergy;
162 }
163 
164 template <typename T>
165 double EcalDeadChannelRecoveryNN<T>::estimateEnergy(double *M3x3Input, double epsilon) {
166  int missing[9];
167  int missing_index = 0;
168 
169  for (int i = 0; i < 9; i++) {
170  if (TMath::Abs(M3x3Input[i]) < epsilon) {
171  missing[missing_index++] = i;
172  } else {
173  // Generally the "dead" cells are allowed to have negative energies (since they will be estimated by the ANN anyway).
174  // But all the remaining "live" ones must have positive values otherwise the logarithm fails.
175 
176  if (M3x3Input[i] < 0.0) { return -2000000.0; }
177  }
178  }
179 
180  // Currently EXACTLY ONE AND ONLY ONE dead cell is corrected. Return -1000000.0 if zero DC's detected and -101.0 if more than one DC's exist.
181  int idxDC = -1 ;
182  if (missing_index == 0) { return -1000000.0; } // Zero DC's were detected
183  if (missing_index > 1) { return -1000001.0; } // More than one DC's detected.
184  if (missing_index == 1) { idxDC = missing[0]; }
185 
186  // Arrange inputs into an array of 8, excluding the dead cell;
187  int input_order[9] = { CC, RR, LL, UU, DD, RU, RD, LU, LD };
188  int input_index = 0;
189  Double_t input[8];
190 
191  for (int id : input_order) {
192  if (id == idxDC)
193  continue;
194 
195  input[input_index++] = TMath::Log(M3x3Input[id]);
196  }
197 
198  // Select the case to apply the appropriate NN and return the result.
199  M3x3Input[idxDC] = TMath::Exp(ctx_[idxDC].mlp->Evaluate(0, input));
200  return M3x3Input[idxDC];
201 }
202 
203 template <typename DetIdT>
204 double EcalDeadChannelRecoveryNN<DetIdT>::makeNxNMatrice_RelMC(DetIdT itID, const EcalRecHitCollection& hit_collection, double *MNxN_RelMC, bool* AcceptFlag) {
205  // Since ANN corrects within a 3x3 window, the possible candidate 3x3 windows that contain
206  // the "dead" crystal form a 5x5 window around it (totaly eight 3x3 windows overlapping).
207  // Get this 5x5 and locate the Max.Contain.Crystal within.
208 
209  // Get the 5x5 window around the "dead" crystal -> vector "NxNaroundDC"
210  std::vector<DetId> NxNaroundDC = topology_->getWindow(itID, 5, 5);
211 
212  DetIdT CellMax ; // Create a null DetId
213  double EnergyMax = 0.0;
214 
215  // Loop over all cells in the vector "NxNaroundDC", and for each cell find it's energy
216  // (from the EcalRecHits collection). Use this energy to detect the Max.Cont.Crystal.
217  std::vector<DetId>::const_iterator theCells;
218 
219  for (theCells = NxNaroundDC.begin(); theCells != NxNaroundDC.end(); ++theCells) {
220  DetIdT cell = DetIdT(*theCells);
221 
222  if (! cell.null()) {
223  EcalRecHitCollection::const_iterator goS_it = hit_collection.find(cell);
224 
225  if ( goS_it != hit_collection.end() && goS_it->energy() >= EnergyMax ) {
226  EnergyMax = goS_it->energy();
227  CellMax = cell;
228  }
229  }
230  }
231 
232  // No Max.Cont.Crystal found, return back with no changes.
233  if ( CellMax.null() ) { *AcceptFlag=false ; return 0.0 ; }
234 
235 #if 1
236  // Not a smart check, because why not just get 4x4 matrix and have a guaranteed hit?
237 
238  // Check that the "dead" crystal belongs to the 3x3 around Max.Cont.Crystal
239  bool dcIn3x3 = false ;
240 
241  std::vector<DetId> NxNaroundMaxCont = topology_->getWindow(CellMax,3,3);
242  std::vector<DetId>::const_iterator testCell;
243  for (testCell = NxNaroundMaxCont.begin(); testCell != NxNaroundMaxCont.end(); ++testCell) {
244  if ( itID == DetIdT(*testCell) ) { dcIn3x3 = true ; }
245  }
246 
247  // If the "dead" crystal is outside the 3x3 then do nothing.
248  if (!dcIn3x3) { *AcceptFlag=false ; return 0.0 ; }
249 #endif
250 
251  return makeNxNMatrice_RelDC(CellMax, hit_collection, MNxN_RelMC, AcceptFlag);
252 }
253 
254 template <>
255 double EcalDeadChannelRecoveryNN<EBDetId>::makeNxNMatrice_RelDC(EBDetId id, const EcalRecHitCollection& hit_collection, double *MNxN, bool* AcceptFlag) {
256  // Make an ANN 3x3 energy matrix around the crystal.
257  // If the "dead" crystal is at the EB boundary (+/- 85) do nothing since we need a full 3x3 around it.
258  if ( id.ieta() == 85 || id.ieta() == -85 ) { *AcceptFlag=false ; return 0.0 ; }
259 
260  // Define the ieta and iphi steps (zero, plus, minus)
261  int ietaZ = id.ieta() ;
262  int ietaP = ( ietaZ == -1 ) ? 1 : ietaZ + 1 ;
263  int ietaN = ( ietaZ == 1 ) ? -1 : ietaZ - 1 ;
264 
265  int iphiZ = id.iphi() ;
266  int iphiP = ( iphiZ == 360 ) ? 1 : iphiZ + 1 ;
267  int iphiN = ( iphiZ == 1 ) ? 360 : iphiZ - 1 ;
268 
269  for (int i=0; i<9; i++) { MNxN[i] = 0.0 ; }
270 
271  // Loop over all cells in the vector "window", and fill the MNxN matrix
272  // to be passed to the ANN for prediction.
273  std::vector<DetId> window = topology_->getWindow(id, 3, 3);
274 
275  std::vector<DetId>::const_iterator itCells;
276  for (itCells = window.begin(); itCells != window.end(); ++itCells) {
277  EBDetId cell = EBDetId(*itCells);
278 
279  if (! cell.null()) {
280  EcalRecHitCollection::const_iterator goS_it = hit_collection.find(cell);
281 
282  if (goS_it != hit_collection.end()) {
283  double energy = goS_it->energy();
284 
285  if ( cell.ieta() == ietaP && cell.iphi() == iphiP ) { MNxN[RU] = energy; }
286  else if ( cell.ieta() == ietaP && cell.iphi() == iphiZ ) { MNxN[RR] = energy; }
287  else if ( cell.ieta() == ietaP && cell.iphi() == iphiN ) { MNxN[RD] = energy; }
288 
289  else if ( cell.ieta() == ietaZ && cell.iphi() == iphiP ) { MNxN[UU] = energy; }
290  else if ( cell.ieta() == ietaZ && cell.iphi() == iphiZ ) { MNxN[CC] = energy; }
291  else if ( cell.ieta() == ietaZ && cell.iphi() == iphiN ) { MNxN[DD] = energy; }
292 
293  else if ( cell.ieta() == ietaN && cell.iphi() == iphiP ) { MNxN[LU] = energy; }
294  else if ( cell.ieta() == ietaN && cell.iphi() == iphiZ ) { MNxN[LL] = energy; }
295  else if ( cell.ieta() == ietaN && cell.iphi() == iphiN ) { MNxN[LD] = energy; }
296 
297  else { *AcceptFlag=false ; return 0.0 ;}
298  }
299  }
300  }
301 
302  // Get the sum of 8
303  double ESUMis = 0.0 ;
304  for (int i=0; i<9; i++) { ESUMis = ESUMis + MNxN[i] ; }
305 
306  *AcceptFlag=true ;
307  return ESUMis;
308 }
309 
310 template <>
311 double EcalDeadChannelRecoveryNN<EEDetId>::makeNxNMatrice_RelDC(EEDetId itID,const EcalRecHitCollection& hit_collection, double *MNxN, bool* AcceptFlag) {
312  // Make an ANN 3x3 energy matrix around the crystal.
313  // If the "dead" crystal is at the EB boundary (+/- 85) do nothing since we need a full 3x3 around it.
314 
315  // If the "dead" crystal is at the EE boundary (inner or outer ring) do nothing since we need a full 3x3 around it.
316  if ( EEDetId::isNextToRingBoundary(itID) ) { *AcceptFlag=false ; return 0.0 ; }
317 
318  // Define the ix and iy steps (zero, plus, minus)
319  int ixZ = itID.ix() ;
320  int ixP = ixZ + 1 ;
321  int ixN = ixZ - 1 ;
322 
323  int iyZ = itID.iy() ;
324  int iyP = iyZ + 1 ;
325  int iyN = iyZ - 1 ;
326 
327  for (int i=0; i<9; i++) { MNxN[i] = 0.0 ; }
328 
329  // Take the "dead" crystal as reference and get the 3x3 around it.
330  std::vector<DetId> NxNaroundRefXtal = topology_->getWindow(itID,3,3);
331 
332  // Loop over all cells in the vector "NxNaroundRefXtal", and fill the MNxN matrix
333  // to be passed to the ANN for prediction.
334  std::vector<DetId>::const_iterator itCells;
335 
336  for (itCells = NxNaroundRefXtal.begin(); itCells != NxNaroundRefXtal.end(); ++itCells) {
337  EEDetId cell = EEDetId(*itCells);
338 
339  if (!cell.null()) {
340  EcalRecHitCollection::const_iterator goS_it = hit_collection.find(cell);
341 
342  if ( goS_it != hit_collection.end() ) {
343  double energy = goS_it->energy();
344 
345  if ( cell.ix() == ixP && cell.iy() == iyP ) { MNxN[RU] = energy; }
346  else if ( cell.ix() == ixP && cell.iy() == iyZ ) { MNxN[RR] = energy; }
347  else if ( cell.ix() == ixP && cell.iy() == iyN ) { MNxN[RD] = energy; }
348 
349  else if ( cell.ix() == ixZ && cell.iy() == iyP ) { MNxN[UU] = energy; }
350  else if ( cell.ix() == ixZ && cell.iy() == iyZ ) { MNxN[CC] = energy; }
351  else if ( cell.ix() == ixZ && cell.iy() == iyN ) { MNxN[DD] = energy; }
352 
353  else if ( cell.ix() == ixN && cell.iy() == iyP ) { MNxN[LU] = energy; }
354  else if ( cell.ix() == ixN && cell.iy() == iyZ ) { MNxN[LL] = energy; }
355  else if ( cell.ix() == ixN && cell.iy() == iyN ) { MNxN[LD] = energy; }
356 
357  else { *AcceptFlag=false ; return 0.0 ;}
358  }
359  }
360  }
361 
362  // Get the sum of 8
363  double ESUMis = 0.0 ;
364  for (int i=0; i<9; i++) { ESUMis = ESUMis + MNxN[i] ; }
365 
366  *AcceptFlag=true ;
367  return ESUMis;
368 }
369 
370 
dictionary missing
Definition: combine.py:4
int i
Definition: DBlmapReader.cc:9
def window
Definition: svgfig.py:642
int ix() const
Definition: EEDetId.h:76
std::vector< EcalRecHit >::const_iterator const_iterator
double makeNxNMatrice_RelMC(DetIdT itID, const EcalRecHitCollection &hit_collection, double *MNxN_RelMC, bool *AccFlag)
#define NULL
Definition: scimark2.h:8
double makeNxNMatrice_RelDC(DetIdT itID, const EcalRecHitCollection &hit_collection, double *MNxN_RelDC, bool *AccFlag)
static std::string const input
Definition: EdmProvDump.cc:43
int iphi() const
get the crystal iphi
Definition: EBDetId.h:53
def load
Definition: svgfig.py:546
T Abs(T a)
Definition: MathUtil.h:49
const CaloSubdetectorTopology * topology_
static bool isNextToRingBoundary(EEDetId id)
Definition: EEDetId.cc:375
double recover(const DetIdT id, const EcalRecHitCollection &hit_collection, double Sum8Cut, bool *AcceptFlag)
int iy() const
Definition: EEDetId.h:82
int ieta() const
get the crystal ieta
Definition: EBDetId.h:51
MultiLayerPerceptronContext ctx_[9]
const_iterator end() const
void setCaloTopology(const CaloTopology *topo)
double estimateEnergy(double *M3x3Input, double epsilon=0.0000001)
virtual std::vector< DetId > getWindow(const DetId &id, const int &northSouthSize, const int &eastWestSize) const
bool null() const
is this a null id ?
Definition: DetId.h:45
const CaloSubdetectorTopology * getSubdetectorTopology(const DetId &id) const
access the subdetector Topology for the given subdetector directly
Definition: CaloTopology.cc:26
iterator find(key_type k)
void load_file(MultiLayerPerceptronContext &ctx, std::string fn)
volatile std::atomic< bool > shutdown_flag false
std::string fullPath() const
Definition: FileInPath.cc:165
long double T