CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
CSCGasGainCorrectionDBConditions.h
Go to the documentation of this file.
1 #ifndef _CSCGASGAINCORRECTIONDBCONDITIONS_H
2 #define _CSCGASGAINCORRECTIONDBCONDITIONS_H
3 
4 #include <memory>
5 #include <cmath>
15 
20 
22  public:
25 
27 
29 
31 
32  private:
33  // ----------member data ---------------------------
36 
37  //Flag for determining if this is for setting MC or data corrections
38  bool isForMC;
39  //File for reading 55944 gas gain corrections. MC will be fake;
40  std::string dataCorrFileName;
41 
42 };
43 
47 
48 #include<fstream>
49 #include<vector>
50 #include<iostream>
51 
52 
53 
54 // to workaround plugin library
56 {
57  if (isMC)
58  printf("\n Generating fake DB constants for MC\n");
59  else {
60  printf("\n Reading gas gain corrections from file %s \n",filename.data());
61  }
62 
63  CSCIndexer indexer;
64 
65  const int MAX_SIZE = 55944;
66 
68 
69  CSCDBGasGainCorrection::GasGainContainer & itemvector = cndbGasGainCorr->gasGainCorr;
70  itemvector.resize(MAX_SIZE);
71 
72  //Filling corrections for MC is very simple
73  if (isMC){
74  for (int i=0;i<MAX_SIZE;i++){
75  itemvector[i].gainCorr = 1.;
76  }
77  return cndbGasGainCorr;
78  }
79 
80  struct gain_info {
81  int gas_gain_index ;
82  int endcap ;
83  int station ;
84  int ring ;
85  int chamber ;
86  int layer ;
87  int hvsegment ;
88  int cfeb ;
89  int nentries ;
90  float mean ;
91  float truncated_mean ;
92  float gas_gain_correction;
93  } gains[MAX_SIZE];
94 
95 
96  for (int j=0; j<MAX_SIZE; j++) {
97  gains[j].gas_gain_index = -999;
98  gains[j].endcap = -999;
99  gains[j].station = -999;
100  gains[j].ring = -999;
101  gains[j].chamber = -999;
102  gains[j].layer = -999;
103  gains[j].hvsegment = -999;
104  gains[j].cfeb = -999;
105  gains[j].nentries = -999;
106  gains[j].mean = -999.;
107  gains[j].truncated_mean = -999.;
108  gains[j].gas_gain_correction= -999.;
109  }
110 
111  FILE *fin = fopen(filename.data(),"r");
112 
113  int linecounter = 0; // set the line counter to the first serial number in the file....
114 
115  while (!feof(fin)){
116  //note space at end of format string to convert last \n
117  int check = fscanf(fin,"%d %d %d %d %d %d %d %d %d %f %f %f \n",
118  &gains[linecounter].gas_gain_index ,
119  &gains[linecounter].endcap ,
120  &gains[linecounter].station ,
121  &gains[linecounter].ring ,
122  &gains[linecounter].chamber ,
123  &gains[linecounter].layer ,
124  &gains[linecounter].hvsegment ,
125  &gains[linecounter].cfeb ,
126  &gains[linecounter].nentries ,
127  &gains[linecounter].mean ,
128  &gains[linecounter].truncated_mean ,
129  &gains[linecounter].gas_gain_correction);
130 
131  if (check != 12){
132  printf("The input file format is not as expected\n");
133  assert(0);
134  }
135 
136  linecounter++;
137 
138  }
139 
140  fclose(fin);
141 
142  if (linecounter == MAX_SIZE) {
143  std::cout << "Total number of gas gains read in = " << linecounter << std::endl;
144  } else {
145  std::cout << "ERROR: Total number of gas-gains read in = " << linecounter
146  << " while total number expected = " << MAX_SIZE << std::endl;
147  }
148 
149  // Fill the chip corrections with values from the file
150  for (int i=0;i<MAX_SIZE;i++){
151 
152  itemvector[i].gainCorr = 0.;
153 
154  if (gains[i].gas_gain_correction > 0.) {
155  itemvector[i].gainCorr = gains[i].gas_gain_correction;
156  } else {
157  // if there is no value, this should be fixed...
158  std::cout << "ERROR: gas_gain_correction < 0 for index " << gains[i].gas_gain_index << std::endl;
159  }
160  }
161 
162  return cndbGasGainCorr;
163 }
164 
165 #endif
int i
Definition: DBlmapReader.cc:9
void setIntervalFor(const edm::eventsetup::EventSetupRecordKey &, const edm::IOVSyncValue &, edm::ValidityInterval &)
CSCGasGainCorrectionDBConditions(const edm::ParameterSet &)
bool check(const DataFrame &df, bool capcheck, bool dvercheck)
int j
Definition: DBlmapReader.cc:9
static CSCDBGasGainCorrection * prefillDBGasGainCorrection(bool isForMC, std::string dataCorrFileName)
ReturnType produceDBGasGainCorrection(const CSCDBGasGainCorrectionRcd &)
std::vector< Item > GasGainContainer
tuple filename
Definition: lut2db_cfg.py:20
tuple cout
Definition: gather_cfg.py:121