CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_10_patch2/src/CalibCalorimetry/HcalStandardModules/src/HcalCholeskyDecomp.cc

Go to the documentation of this file.
00001 // Original Author:  Annabel Downs
00002 //         Created:  Wed Jul 29 10:54:11 CEST 2009
00003 #include <memory>
00004 #include "CalibCalorimetry/HcalStandardModules/interface/HcalCholeskyDecomp.h"
00005 
00006 HcalCholeskyDecomp::HcalCholeskyDecomp(const edm::ParameterSet& iConfig)
00007 {
00008   outfile = iConfig.getUntrackedParameter<std::string>("outFile","CholeskyMatrices.txt");
00009 }
00010 
00011 
00012 HcalCholeskyDecomp::~HcalCholeskyDecomp()
00013 {
00014 }
00015 
00016 void
00017 HcalCholeskyDecomp::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup)
00018 {
00019    using namespace edm;
00020    using namespace std;
00021 
00022    edm::ESHandle<HcalCovarianceMatrices> refCov;
00023    iSetup.get<HcalCovarianceMatricesRcd>().get("reference",refCov);
00024    const HcalCovarianceMatrices* myCov = refCov.product(); //Fill emap from database
00025 
00026    double sig[4][10][10];
00027    double c[4][10][10], cikik[4], cikjk[4];
00028 
00029    HcalCholeskyMatrices * outMatrices = new HcalCholeskyMatrices();
00030 
00031    std::vector<DetId> listChan = myCov->getAllChannels();
00032    std::vector<DetId>::iterator cell;
00033 
00034    int HBcount = 0;
00035    int HEcount = 0;
00036    int HFcount = 0;
00037    int HOcount = 0;
00038 
00039    double HBmatrix[4][10][10];
00040    double HEmatrix[4][10][10];
00041    double HFmatrix[4][10][10];
00042    double HOmatrix[4][10][10];
00043 
00044    double tempmatrix[4][10][10] ;
00045 
00046 
00047    for(int m= 0; m != 4; m++){
00048       for(int i = 0; i != 10; i++){
00049          for(int j = 0; j!=10; j++){
00050           HBmatrix[m][i][j] = 0;
00051           HEmatrix[m][i][j] = 0;
00052           HFmatrix[m][i][j] = 0;
00053           HOmatrix[m][i][j] = 0;
00054          }
00055       }
00056    }
00057 
00058    for (std::vector<DetId>::iterator it = listChan.begin(); it != listChan.end(); it++)
00059    {
00060    for(int m= 0; m != 4; m++){
00061       for(int i = 0; i != 10; i++){
00062          for(int j = 0; j!=10; j++){
00063             sig[m][i][j] = 0;
00064             c[m][i][j] = 0;
00065             cikik[m] =0;
00066             cikjk[m] = 0;
00067             tempmatrix[m][i][j] = 0;
00068          }
00069       }
00070    }
00071 
00072    const HcalCovarianceMatrix * CMatrix = myCov->getValues(*it);
00073    HcalDetId hcalid(it->rawId());
00074 
00075    for(int m = 0; m != 4; m++){
00076       for(int i = 0; i != 10; i++){
00077          for(int j = 0; j != 10; j++) {sig[m][i][j] = CMatrix->getValue(m,i,j);}
00078         }
00079    }
00080 
00081 //Method to check 4x10 storage
00082 //  for(int m = 0; m != 4; m++){
00083 //      for(int i = 0; i != 6; i++){
00084 //         for(int j=(i+4); j!=10; j++) {sig[m][i][j] = 0;}
00085 //        }
00086 //   }
00087         
00088 
00089    for(int m = 0; m != 4; m++){
00090       for(int i = 0; i != 10; i++){
00091          for(int j = i; j != 10; j++){ sig[m][j][i] = sig[m][i][j];}
00092       }
00093    }
00094    //.......................Cholesky Decomposition.............................
00095    //Step 1a
00096    for(int m = 0; m!=4; m++)
00097          {
00098         
00099                 c[m][0][0]=sqrt(sig[m][0][0]);
00100         }
00101   for(int m = 0; m!=4; m++)
00102    {
00103       for(int i = 1; i != 10; i++)
00104         {
00105          c[m][i][0] = (sig[m][0][i] / c[m][0][0]);
00106         }
00107    }
00108    
00109  //Chelesky Matrices
00110    for(int m = 0; m!=4; m++)
00111    {
00112       c[m][1][1] = sqrt(sig[m][1][1] - (c[m][1][0]*c[m][1][0])); // i) step 2a of chelesky decomp
00113       if (((sig[m][1][1]) - (c[m][1][0]*c[m][1][0]))<=0) continue;
00114       for(int i = 2; i !=10; i++)
00115         {
00116         cikik[m] = 0;
00117          for(int j=1; j!= i; j++)
00118                 {
00119                 cikjk[m] = 0;
00120                     for(int k=0; k != j; k++)
00121                         {
00122                         cikjk[m] += (c[m][i][k]*c[m][j][k]); // ii)  step 2a of chelesky decomp
00123                         }
00124                 c[m][i][j] = (sig[m][i][j] - cikjk[m])/c[m][j][j]; // step 3a of chelesky decomp
00125                 }//end of j loop
00126          
00127                 for( int k = 0; k != i; k++)
00128                 {
00129                          cikik[m] += (c[m][i][k]*c[m][i][k]);
00130                 }
00131                 double test = (sig[m][i][i] - cikik[m]);
00132                 if(test > 0 )
00133                 {
00134                         c[m][i][i] = sqrt(test);
00135                 }
00136                 else
00137                 {
00138                         c[m][i][i] = 1000;
00139                 }
00140         }//end of i loop
00141 /*
00142  //Cholesky Matrix for rechit (2-5 ts)
00143         for (int i = 0; i!=2; i++)
00144         {
00145           for (int j=0; j!=10; j++)
00146           {
00147                 c[m][i][j] = 0;
00148           }
00149         }
00150         for (int i = 6; i!=10; i++)
00151         {
00152           for (int j=0; j!=10; j++)
00153           {
00154                 c[m][i][j] = 0;
00155           }
00156         }
00157 */
00158 
00159 
00160    }//end of m loop
00161    
00162 
00163    HcalCholeskyMatrix * item = new HcalCholeskyMatrix(it->rawId());
00164    for(int m = 0; m != 4; m++){
00165       for(int i = 0; i != 10; i++){
00166          for(int j = i; j != 10; j++){ 
00167              item->setValue(m,j,i, c[m][j][i]);
00168              tempmatrix[m][j][i] = c[m][j][i];
00169          }//sig[m][j][i]
00170       }
00171    }
00172 
00173    if(hcalid.subdet()==1)
00174    {
00175       for(int m = 0; m != 4; m++)
00176          for(int i = 0; i != 10; i++)
00177             for(int j = 0; j != 10; j++)
00178                HBmatrix[m][i][j] += tempmatrix[m][i][j];
00179       HBcount++;
00180    }
00181    if(hcalid.subdet()==2)
00182    {
00183       for(int m = 0; m != 4; m++)
00184          for(int i = 0; i != 10; i++)
00185             for(int j = 0; j != 10; j++)
00186                HEmatrix[m][i][j] += tempmatrix[m][i][j];
00187       HEcount++;
00188    }
00189    if(hcalid.subdet()==3)
00190    {
00191       for(int m = 0; m != 4; m++)
00192          for(int i = 0; i != 10; i++)
00193             for(int j = 0; j != 10; j++)
00194                HOmatrix[m][i][j] += tempmatrix[m][i][j];
00195       HOcount++;
00196    }
00197    if(hcalid.subdet()==4)
00198    {
00199       for(int m = 0; m != 4; m++)
00200          for(int i = 0; i != 10; i++)
00201             for(int j = 0; j != 10; j++)
00202                HFmatrix[m][i][j] += tempmatrix[m][i][j];
00203       HFcount++;
00204    }
00205 
00206 
00207    
00208    outMatrices->addValues(*item);
00209 
00210    }
00211 
00212    for(int m = 0; m != 4; m++)
00213    {
00214       for(int i = 0; i != 10; i++)
00215       {
00216          for(int j = 0; j != 10; j++)
00217          {
00218             HBmatrix [m][i][j] /= HBcount;
00219             HEmatrix [m][i][j] /= HEcount;
00220             HFmatrix [m][i][j] /= HFcount;
00221             HOmatrix [m][i][j] /= HOcount;
00222          }
00223       }
00224    }
00225 
00226   std::vector<DetId> listResult = outMatrices->getAllChannels();
00227 
00228   int n_avg = 0;
00229 
00230   edm::ESHandle<HcalElectronicsMap> refEMap;
00231   iSetup.get<HcalElectronicsMapRcd>().get(refEMap);
00232 //  iSetup.get<HcalElectronicsMapRcd>().get("reference",refEMap);
00233   const HcalElectronicsMap* myRefEMap = refEMap.product();
00234   std::vector<HcalGenericDetId> listEMap = myRefEMap->allPrecisionId();
00235 
00236     for (std::vector<HcalGenericDetId>::const_iterator it = listEMap.begin(); it != listEMap.end(); it++)
00237       {
00238       DetId mydetid = DetId(it->rawId());
00239         if (std::find(listResult.begin(), listResult.end(), mydetid ) == listResult.end()  )
00240           {
00241 //              std::cout << "Cholesky Matrix not found for cell " <<  HcalGenericDetId(it->rawId());
00242               if(it->isHcalDetId())
00243               {
00244               std::cout << "Cholesky Matrix not found for cell " <<  HcalGenericDetId(it->rawId());
00245                  HcalDetId hcalid2(it->rawId());
00246                  HcalCholeskyMatrix * item = new HcalCholeskyMatrix(it->rawId());
00247                  for(int m = 0; m != 4; m++)
00248                  {
00249                     for(int i = 0; i != 10; i++)
00250                     {
00251                        for(int j = 0; j != 10; j++)
00252                        {
00253                           if(j <= i){
00254                           if(hcalid2.subdet()==1) item->setValue(m,i,j,HBmatrix [m][i][j]);
00255                           if(hcalid2.subdet()==2) item->setValue(m,i,j,HEmatrix [m][i][j]);
00256                           if(hcalid2.subdet()==3) item->setValue(m,i,j,HFmatrix [m][i][j]);
00257                           if(hcalid2.subdet()==4) item->setValue(m,i,j,HOmatrix [m][i][j]);
00258                           }else{
00259                           if(hcalid2.subdet()==1) item->setValue(m,i,j,0);
00260                           if(hcalid2.subdet()==2) item->setValue(m,i,j,0);
00261                           if(hcalid2.subdet()==3) item->setValue(m,i,j,0);
00262                           if(hcalid2.subdet()==4) item->setValue(m,i,j,0);
00263                           }
00264                        }
00265                     }
00266                  }
00267                  outMatrices->addValues(*item);
00268                  std::cout << "... Average Matrix used\n";
00269                  n_avg++;
00270               }
00271           }
00272       }
00273 
00274    std::vector<DetId> checkResult = outMatrices->getAllChannels();
00275    int count = 0;
00276    for(std::vector<DetId>::const_iterator it = checkResult.begin(); it != checkResult.end(); it++) count ++;
00277 
00278    std::cout << "There are " << count << " channels with Cholesky matrices.\n" << "Used  " << n_avg << " average values.\n";
00279 
00280    ofstream cmatout(outfile.c_str());
00281    HcalDbASCIIIO::dumpObject(cmatout, *outMatrices);
00282 }
00283 
00284 // ------------ meth(d called once each job just before starting event loop  ------------
00285 void 
00286 HcalCholeskyDecomp::beginJob()
00287 {
00288 }
00289 
00290 // ------------ method called once each job just after ending the event loop  ------------
00291 void 
00292 HcalCholeskyDecomp::endJob() {
00293 }
00294 
00295 //define this as a plug-in
00296 DEFINE_FWK_MODULE(HcalCholeskyDecomp);