CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_1_8_patch12/src/Calibration/EcalCalibAlgos/src/IMACalibBlock.cc

Go to the documentation of this file.
00001 
00008 #include "Calibration/EcalCalibAlgos/interface/IMACalibBlock.h"
00009 #include "Calibration/EcalCalibAlgos/interface/BlockSolver.h"
00010 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00011 #include "TH1F.h"
00012 #include "TFile.h"
00013 #include <cassert>
00014 
00015 // -----------------------------------------------------
00016 
00017 
00018 IMACalibBlock::IMACalibBlock (const int numberOfElements):
00019   VEcalCalibBlock (numberOfElements) ,
00020   m_kaliVector (m_numberOfElements) , 
00021   m_kaliMatrix (evalX2Size ())
00022 {
00023   reset () ;
00024 }
00025 
00026 
00027 // -----------------------------------------------------
00028 
00029 
00030 IMACalibBlock::~IMACalibBlock()
00031 {
00032 }
00033 
00034 
00035 // -----------------------------------------------------
00036 
00037 
00038 void
00039 IMACalibBlock::Fill (std::map<int,double>::const_iterator MapBegin,
00040                      std::map<int,double>::const_iterator MapEnd ,
00041                      double pTk,
00042                      double pSubtract,
00043                      double sigma)
00044 {
00045 //  std::cerr<<"\n\nfaccio il fill!\n";
00046   double inverror = 1./sigma ; 
00047   //LP fist loop over energies
00048   for (std::map<int,double>::const_iterator itMap1 =MapBegin ; 
00049        itMap1 !=MapEnd ; 
00050        ++itMap1)
00051     {
00052 //      std::cerr<<"numero "<<itMap1->first<<" vale "<<itMap1->second<<"\t"; 
00053       for (std::map<int,double>::const_iterator itMap2 = itMap1 ; 
00054       itMap2 !=MapEnd ; 
00055       ++itMap2)  
00056         {
00057           //LP calculate the chi square value
00058           double dummy = itMap1->second * itMap2->second;
00059           dummy *= inverror ;
00060           //LP fill the calib matrix
00061           m_kaliMatrix.at(itMap1->first*m_numberOfElements +itMap2->first) += dummy ;            
00062         } //LP second loop over xtals
00063 
00064       //LP calculate the vector value
00065       double dummy = pTk ;
00066       dummy -= pSubtract ;                        //LP borders
00067       dummy *= itMap1->second ;
00068       dummy *= inverror ;
00069       //LP fill the calib vector
00070       m_kaliVector.at(itMap1->first) += dummy ;
00071     } //LP first loop over energies
00072   return ;
00073 }
00074 
00075 
00076 //------------------------------------------------------------
00077 
00078 
00079 void 
00080 IMACalibBlock::complete ()
00081 {
00082  int bef;
00083  int aft;
00084  for (unsigned int i=0; i<m_numberOfElements;++i)
00085   {
00086     for (unsigned int j=i+1; j< m_numberOfElements; ++j) 
00087       {
00088          bef = (i*m_numberOfElements+j);
00089          aft =  (j*m_numberOfElements +i);
00090             m_kaliMatrix.at(aft) = m_kaliMatrix.at(bef) ;
00091           } //LP second loop over xtals
00092       } //LP first loop over xtals
00093 
00094     return ;
00095   }
00096 
00097 
00098 // ------------------------------------------------------------
00099 
00100 
00101 int
00102 IMACalibBlock::solve (int usingBlockSolver, double min, double max)
00103 {
00104  complete () ;
00105 // TH1F vettore ("vettore","vettore",10,-0.1,9.9);
00106 // TH1F matrice ("matrice","matrice",100,-0.1,99.9);
00107  CLHEP::HepMatrix kaliMatrix (m_numberOfElements,m_numberOfElements) ;
00108 // for (std::vector<double>::iterator it = m_kaliVector.begin();
00109 //               it!= m_kaliVector.end();++it)
00110 //       vettore.Fill(it-m_kaliVector.begin(),*it);
00111  riempiMtr (m_kaliMatrix , kaliMatrix) ;
00112 // for (std::vector<double>::iterator it = m_kaliMatrix.begin();
00113 //               it!= m_kaliMatrix.end();++it)
00114 //       matrice.Fill(it-m_kaliMatrix.begin(),*it);
00115 //  TFile f ("fileInteressante.root","RECREATE");
00116 //  vettore.Write();
00117 // matrice.Write();
00118  CLHEP::HepVector kaliVector (m_numberOfElements) ;
00119  riempiVtr (m_kaliVector , kaliVector) ;
00120  //PG linear system solution
00121  CLHEP::HepVector result = CLHEP::solve (kaliMatrix,kaliVector) ;
00122   for (int i = 0 ; i < kaliVector.num_row () ; ++i)
00123  if (result.normsq () < min * kaliVector.num_row () ||
00124      result.normsq () > max * kaliVector.num_row ()) 
00125    {
00126    if (usingBlockSolver)  
00127      {
00128         edm::LogWarning ("IML") << "using  blocSlover " << std::endl ;
00129         BlockSolver() (kaliMatrix,kaliVector,result) ;
00130      }
00131    else 
00132      {
00133        edm::LogWarning ("IML") <<"coeff out of range " <<std::endl;
00134        for (int i = 0 ; i < kaliVector.num_row () ; ++i)
00135              result[i] = 1. ;
00136      }
00137    }
00138  fillMap(result);
00139  return 0;
00140 }
00141 
00142 
00143 //------------------------------------------------------------
00144 
00145 
00146 inline int 
00147 IMACalibBlock::evalX2Size ()
00148   {
00149 
00150     return m_numberOfElements* m_numberOfElements;
00151   }
00152 
00153 
00154 // ------------------------------------------------------------
00155 
00156 
00157 void
00158 IMACalibBlock::riempiMtr (const std::vector<double> & piena, 
00159                           CLHEP::HepMatrix & vuota) 
00160   {
00161     unsigned int max = m_numberOfElements ;
00162 
00163     assert (piena.size () == max * max) ; 
00164     assert (vuota.num_row () == int(max)) ;
00165     assert (vuota.num_col () == int(max)) ;
00166     for (unsigned int i = 0 ; i < max ; ++i)
00167      for (unsigned int j = 0 ; j < max ; ++j)
00168        if (std::isnan (piena[i*max + j])) vuota[i][j] = 0. ;
00169          else vuota[i][j] = piena[i*max + j] ; 
00170 
00171     return ;
00172   }
00173 
00174 
00175 // ------------------------------------------------------------
00176 
00177 void
00178 IMACalibBlock::riempiVtr (const std::vector<double> & pieno, 
00179                           CLHEP::HepVector & vuoto) 
00180   {
00181 
00182     int max = m_numberOfElements ;
00183     assert (vuoto.num_row () == max) ;
00184     for (int i = 0 ; i < max ; ++i)
00185       if (std::isnan (pieno[i])) vuoto[i] = 0. ;
00186       else vuoto[i] = pieno[i] ; 
00187 
00188     return ;
00189   }
00190 
00191 
00192 // ------------------------------------------------------------
00193 
00194 
00195 void 
00196 IMACalibBlock::reset () 
00197 {
00198 
00199   for (std::vector<double>::iterator vecIt = m_kaliVector.begin () ;
00200        vecIt != m_kaliVector.end ();
00201        ++vecIt)
00202     {
00203       *vecIt = 0. ;
00204     }  
00205 
00206   for (std::vector<double>::iterator vecIt = m_kaliMatrix.begin () ;
00207        vecIt != m_kaliMatrix.end ();
00208        ++vecIt)
00209     {
00210       *vecIt = 0. ;
00211     }  
00212 
00213 }
00214 
00215 
00216 // ------------------------------------------------------------
00217 //LP sta provando a fare i seguenti due metodi come locali. Speriamo di non fare stronzate.
00218 
00219 
00220 void 
00221 IMACalibBlock::fillMap (const CLHEP::HepVector & result) 
00222 {
00223 
00224    for (unsigned int i=0; i < m_numberOfElements; ++i)
00225       {
00226          m_coefficients[i] = result[i] ;
00227       } 
00228 
00229    return ;
00230 }
00231