CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
IMACalibBlock.cc
Go to the documentation of this file.
1 
11 #include "TH1F.h"
12 #include "TFile.h"
13 #include <cassert>
14 
15 // -----------------------------------------------------
16 
17 
18 IMACalibBlock::IMACalibBlock (const int numberOfElements):
19  VEcalCalibBlock (numberOfElements) ,
20  m_kaliVector (m_numberOfElements) ,
21  m_kaliMatrix (evalX2Size ())
22 {
23  reset () ;
24 }
25 
26 
27 // -----------------------------------------------------
28 
29 
31 {
32 }
33 
34 
35 // -----------------------------------------------------
36 
37 
38 void
39 IMACalibBlock::Fill (std::map<int,double>::const_iterator MapBegin,
40  std::map<int,double>::const_iterator MapEnd ,
41  double pTk,
42  double pSubtract,
43  double sigma)
44 {
45 // std::cerr<<"\n\nfaccio il fill!\n";
46  double inverror = 1./sigma ;
47  //LP fist loop over energies
48  for (std::map<int,double>::const_iterator itMap1 =MapBegin ;
49  itMap1 !=MapEnd ;
50  ++itMap1)
51  {
52 // std::cerr<<"numero "<<itMap1->first<<" vale "<<itMap1->second<<"\t";
53  for (std::map<int,double>::const_iterator itMap2 = itMap1 ;
54  itMap2 !=MapEnd ;
55  ++itMap2)
56  {
57  //LP calculate the chi square value
58  double dummy = itMap1->second * itMap2->second;
59  dummy *= inverror ;
60  //LP fill the calib matrix
61  m_kaliMatrix.at(itMap1->first*m_numberOfElements +itMap2->first) += dummy ;
62  } //LP second loop over xtals
63 
64  //LP calculate the vector value
65  double dummy = pTk ;
66  dummy -= pSubtract ; //LP borders
67  dummy *= itMap1->second ;
68  dummy *= inverror ;
69  //LP fill the calib vector
70  m_kaliVector.at(itMap1->first) += dummy ;
71  } //LP first loop over energies
72  return ;
73 }
74 
75 
76 //------------------------------------------------------------
77 
78 
79 void
81 {
82  int bef;
83  int aft;
84  for (unsigned int i=0; i<m_numberOfElements;++i)
85  {
86  for (unsigned int j=i+1; j< m_numberOfElements; ++j)
87  {
88  bef = (i*m_numberOfElements+j);
89  aft = (j*m_numberOfElements +i);
90  m_kaliMatrix.at(aft) = m_kaliMatrix.at(bef) ;
91  } //LP second loop over xtals
92  } //LP first loop over xtals
93 
94  return ;
95  }
96 
97 
98 // ------------------------------------------------------------
99 
100 
101 int
102 IMACalibBlock::solve (int usingBlockSolver, double min, double max)
103 {
104  complete () ;
105 // TH1F vettore ("vettore","vettore",10,-0.1,9.9);
106 // TH1F matrice ("matrice","matrice",100,-0.1,99.9);
107  CLHEP::HepMatrix kaliMatrix (m_numberOfElements,m_numberOfElements) ;
108 // for (std::vector<double>::iterator it = m_kaliVector.begin();
109 // it!= m_kaliVector.end();++it)
110 // vettore.Fill(it-m_kaliVector.begin(),*it);
111  riempiMtr (m_kaliMatrix , kaliMatrix) ;
112 // for (std::vector<double>::iterator it = m_kaliMatrix.begin();
113 // it!= m_kaliMatrix.end();++it)
114 // matrice.Fill(it-m_kaliMatrix.begin(),*it);
115 // TFile f ("fileInteressante.root","RECREATE");
116 // vettore.Write();
117 // matrice.Write();
118  CLHEP::HepVector kaliVector (m_numberOfElements) ;
119  riempiVtr (m_kaliVector , kaliVector) ;
120  //PG linear system solution
121  CLHEP::HepVector result = CLHEP::solve (kaliMatrix,kaliVector) ;
122  for (int i = 0 ; i < kaliVector.num_row () ; ++i)
123  if (result.normsq () < min * kaliVector.num_row () ||
124  result.normsq () > max * kaliVector.num_row ())
125  {
126  if (usingBlockSolver)
127  {
128  edm::LogWarning ("IML") << "using blocSlover " << std::endl ;
129  BlockSolver() (kaliMatrix,kaliVector,result) ;
130  }
131  else
132  {
133  edm::LogWarning ("IML") <<"coeff out of range " <<std::endl;
134  for (int i = 0 ; i < kaliVector.num_row () ; ++i)
135  result[i] = 1. ;
136  }
137  }
138  fillMap(result);
139  return 0;
140 }
141 
142 
143 //------------------------------------------------------------
144 
145 
146 inline int
148  {
149 
151  }
152 
153 
154 // ------------------------------------------------------------
155 
156 
157 void
158 IMACalibBlock::riempiMtr (const std::vector<double> & piena,
159  CLHEP::HepMatrix & vuota)
160  {
161  unsigned int max = m_numberOfElements ;
162 
163  assert (piena.size () == max * max) ;
164  assert (vuota.num_row () == int(max)) ;
165  assert (vuota.num_col () == int(max)) ;
166  for (unsigned int i = 0 ; i < max ; ++i)
167  for (unsigned int j = 0 ; j < max ; ++j)
168  if (std::isnan (piena[i*max + j])) vuota[i][j] = 0. ;
169  else vuota[i][j] = piena[i*max + j] ;
170 
171  return ;
172  }
173 
174 
175 // ------------------------------------------------------------
176 
177 void
178 IMACalibBlock::riempiVtr (const std::vector<double> & pieno,
179  CLHEP::HepVector & vuoto)
180  {
181 
182  int max = m_numberOfElements ;
183  assert (vuoto.num_row () == max) ;
184  for (int i = 0 ; i < max ; ++i)
185  if (std::isnan (pieno[i])) vuoto[i] = 0. ;
186  else vuoto[i] = pieno[i] ;
187 
188  return ;
189  }
190 
191 
192 // ------------------------------------------------------------
193 
194 
195 void
197 {
198 
199  for (std::vector<double>::iterator vecIt = m_kaliVector.begin () ;
200  vecIt != m_kaliVector.end ();
201  ++vecIt)
202  {
203  *vecIt = 0. ;
204  }
205 
206  for (std::vector<double>::iterator vecIt = m_kaliMatrix.begin () ;
207  vecIt != m_kaliMatrix.end ();
208  ++vecIt)
209  {
210  *vecIt = 0. ;
211  }
212 
213 }
214 
215 
216 // ------------------------------------------------------------
217 //LP sta provando a fare i seguenti due metodi come locali. Speriamo di non fare stronzate.
218 
219 
220 void
221 IMACalibBlock::fillMap (const CLHEP::HepVector & result)
222 {
223 
224  for (unsigned int i=0; i < m_numberOfElements; ++i)
225  {
226  m_coefficients[i] = result[i] ;
227  }
228 
229  return ;
230 }
231 
int i
Definition: DBlmapReader.cc:9
IMACalibBlock(const int)
ctor
void Fill(std::map< int, double >::const_iterator, std::map< int, double >::const_iterator, double pTk, double pSubtract, double sigma=1.)
insert an entry
element for the single ECAL block intercalibration
#define min(a, b)
Definition: mlp_lapack.h:161
std::vector< double > m_kaliVector
vector for the chi2 inversion
Definition: IMACalibBlock.h:61
return((rh^lh)&mask)
void riempiMtr(const std::vector< double > &piena, CLHEP::HepMatrix &vuota)
copy a vector into a CLHEP object
const T & max(const T &a, const T &b)
bool isnan(float x)
Definition: math.h:13
tuple result
Definition: query.py:137
int j
Definition: DBlmapReader.cc:9
std::vector< double > m_kaliMatrix
matrix for the chi2 inversion
Definition: IMACalibBlock.h:63
void riempiVtr(const std::vector< double > &pieno, CLHEP::HepVector &vuoto)
copy a vector into a CLHEP object
~IMACalibBlock()
dtor
int evalX2Size()
give the size of a chi2 matrix
solves at best the matrix invertion for calibration
Definition: BlockSolver.h:25
int solve(int usingBlockSolver, double min, double max)
solve the chi2 linear system
std::map< unsigned int, float > m_coefficients
map of coefficients
void reset()
reset the chi2 matrices
unsigned int m_numberOfElements
The only parameter!
void fillMap(const CLHEP::HepVector &result)
fill the coefficients map from the CLHEP vector solution
void complete()
complete the triangolar chi2 matrix to a sym one