Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 11:59:28

0001 /** \class InvMatrixUtils
0002 
0003     \brief various utilities
0004 
0005 */
0006 
0007 #ifndef InvMatrixUtils_h
0008 #define InvMatrixUtils_h
0009 
0010 #include <string>
0011 #include <map>
0012 #include "TObject.h"
0013 #include "TF1.h"
0014 #include "TH2.h"
0015 #include "TProfile.h"
0016 #include "TCanvas.h"
0017 #include "CLHEP/Geometry/Point3D.h"
0018 #include "CLHEP/Matrix/GenMatrix.h"
0019 #include "CLHEP/Matrix/Matrix.h"
0020 #include "CLHEP/Matrix/Vector.h"
0021 #include "TFile.h"
0022 
0023 //#include "ConfigParser.h"
0024 #include "Calibration/Tools/interface/InvMatrixCommonDefs.h"
0025 
0026 /** set the style for the printout*/
0027 void setStyle();
0028 
0029 /** search for an existing canvas with the name 
0030 and returns the poiter to it */
0031 TCanvas* getGlobalCanvas(std::string name = "Inv MatrixCanvas");
0032 
0033 /** search for an existing TFile with the name 
0034 and returns the poiter to it */
0035 TFile* getGlobalTFile(std::string name = "Inv MatrixTFile.root");
0036 //TFile * getGlobalTFile (const char* name) ;
0037 
0038 /** search for an existing TFile with the name 
0039 and saves it to disk with his name */
0040 int saveGlobalTFile(std::string name = "Inv MatrixFile.root");
0041 
0042 /** search for an existing calib matrix saved with the name 
0043 and returns the poiter to it,
0044 the deletion is responsiblity of the user */
0045 CLHEP::HepMatrix* getSavedMatrix(const std::string& name);
0046 
0047 /** return the impact position of the electron over ECAL */
0048 HepGeom::Point3D<Float_t> TBposition(const Float_t amplit[7][7],
0049                                      const Float_t beamEne,
0050                                      const Float_t w0 = 4.0,
0051                                      const Float_t x0 = 8.9,  //mm
0052                                      const Float_t a0 = 6.2,
0053                                      const Float_t sideX = 24.06,   //mm
0054                                      const Float_t sideY = 22.02);  //mm
0055 
0056 /** get the energy in the 5x5 
0057 from the 7x7 array around the most energetic crystal*/
0058 double get5x5(const Float_t energy[7][7]);
0059 
0060 /** get the energy in the 3x3 
0061 from the 7x7 array around the most energetic crystal*/
0062 double get3x3(const Float_t energy[7][7]);
0063 
0064 /**to get the parameters from a congiguration file*/
0065 int parseConfigFile(const TString& config);
0066 
0067 /**to get the crystal number from eta and phi*/
0068 int xtalFromEtaPhi(const int& myEta, const int& myPhi);
0069 
0070 /**to get the crystal number from iEta and iPhi
0071 iEta runs from 1 to 85
0072 iPhi runs from 1 to 20
0073 */
0074 int xtalFromiEtaiPhi(const int& iEta, const int& iPhi);
0075 
0076 /** get the eta coord [0,84] */
0077 int etaFromXtal(const int& xtal);
0078 
0079 /** get the phi coord [0,19] */
0080 int phiFromXtal(const int& xtal);
0081 
0082 /** get the eta coord [1,85] */
0083 int ietaFromXtal(const int& xtal);
0084 
0085 /** get the phi coord [1,20] */
0086 int iphiFromXtal(const int& xtal);
0087 
0088 /** to read a file containing unserted integers
0089 while avoiding comment lines */
0090 int extract(std::vector<int>* output, const std::string& dati);
0091 
0092 /** to write the calibration constants file */
0093 int writeCalibTxt(const CLHEP::HepMatrix& AmplitudeMatrix,
0094                   const CLHEP::HepMatrix& SigmaMatrix,
0095                   const CLHEP::HepMatrix& StatisticMatrix,
0096                   std::string fileName = "calibOutput.txt");
0097 
0098 /** to write the file fpr the CMSSW in the DB compliant format (using Energy as reference)*/
0099 int writeCMSSWCoeff(const CLHEP::HepMatrix& amplMatrix,
0100                     double calibThres,
0101                     float ERef,
0102                     const CLHEP::HepMatrix& sigmaMatrix,
0103                     const CLHEP::HepMatrix& statisticMatrix,
0104                     std::string fileName = "calibOutput.txt",
0105                     std::string genTag = "CAL_GENTAG",
0106                     std::string method = "CAL_METHOD",
0107                     std::string version = "CAL_VERSION",
0108                     std::string type = "CAL_TYPE");
0109 
0110 /** to write the file fpr the CMSSW in the DB compliant format 
0111     (using Crystal as reference) */
0112 int writeCMSSWCoeff(const CLHEP::HepMatrix& amplMatrix,
0113                     double calibThres,
0114                     int etaRef,
0115                     int phiRef,
0116                     const CLHEP::HepMatrix& sigmaMatrix,
0117                     const CLHEP::HepMatrix& statisticMatrix,
0118                     std::string fileName = "calibOutput.txt",
0119                     std::string genTag = "CAL_GENTAG",
0120                     std::string method = "CAL_METHOD",
0121                     std::string version = "CAL_VERSION",
0122                     std::string type = "CAL_TYPE");
0123 
0124 /** translates the calib coefficients format,
0125     from the TB06Studies one to the CMSSSW one */
0126 int translateCoeff(const CLHEP::HepMatrix& calibcoeff,
0127                    const CLHEP::HepMatrix& sigmaMatrix,
0128                    const CLHEP::HepMatrix& statisticMatrix,
0129                    std::string SMnumber = "1",
0130                    double calibThres = 0.01,
0131                    std::string fileName = "calibOutput.txt",
0132                    std::string genTag = "CAL_GENTAG",
0133                    std::string method = "CAL_METHOD",
0134                    std::string version = "CAL_VERSION",
0135                    std::string type = "CAL_TYPE");
0136 
0137 /** translates the calib coefficients format,
0138     from the CMSSW one to the TB06Studies one */
0139 int readCMSSWcoeff(CLHEP::HepMatrix& calibcoeff, const std::string& inputFileName, double defaultVal = 1.);
0140 
0141 /** translates the calib coefficients format,
0142     from the CMSSW one to the TB06Studies one */
0143 int readCMSSWcoeffForComparison(CLHEP::HepMatrix& calibcoeff, const std::string& inputFileName);
0144 
0145 /** smart profiling by double averaging */
0146 TH1D* smartProfile(TH2F* strip, double width);
0147 
0148 /** smart profiling by fixing gaussian parameters and
0149     range from a first averaging */
0150 TH1D* smartGausProfile(TH2F* strip, double width);
0151 
0152 /**
0153    */
0154 TH1D* smartError(TH1D* strip);
0155 
0156 /**
0157 find the effective sigma as the half width of the sub-distribution
0158 containing 68.3% of the total distribution
0159 */
0160 double effectiveSigma(TH1F& histogram, int vSteps = 100);
0161 
0162 /**
0163 find the support of the histogram above a threshold
0164 return the min and max bins
0165 */
0166 std::pair<int, int> findSupport(TH1F& histogram, double thres = 0.);
0167 
0168 /** 
0169 transfers a CLHEP matrix into a double array
0170 with the size of a supermodule
0171 */
0172 void mtrTransfer(double output[SCMaxEta][SCMaxPhi], CLHEP::HepMatrix* input, double Default);
0173 
0174 /**
0175 reset the matrices f the size of a supermodule
0176 */
0177 template <class Type>
0178 void mtrReset(Type superModules[SCMaxEta][SCMaxPhi], const Type val) {
0179   for (int e = 0; e < SCMaxEta; ++e)
0180     for (int p = 0; p < SCMaxPhi; ++p) {
0181       superModules[e][p] = val;
0182     }
0183 }
0184 
0185 /**
0186 correction for eta containment for 3*3 cluster */
0187 double etaCorrE1E9(int eta);
0188 /**
0189 correction for eta containment for 7*7 cluster */
0190 double etaCorrE1E49(int eta);
0191 /**
0192 correction for eta containment for 5*5 cluster */
0193 double etaCorrE1E25(int eta);
0194 
0195 #endif