Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 11:58:20

0001 #ifndef _CSCCHIPSPEEDCORRECTIONDBCONDITIONS_H
0002 #define _CSCCHIPSPEEDCORRECTIONDBCONDITIONS_H
0003 
0004 #include "FWCore/Framework/interface/ESHandle.h"
0005 #include "FWCore/Framework/interface/ESProducer.h"
0006 #include "FWCore/Framework/interface/Event.h"
0007 #include "FWCore/Framework/interface/EventSetup.h"
0008 #include "FWCore/Framework/interface/EventSetupRecordIntervalFinder.h"
0009 #include "FWCore/Framework/interface/Frameworkfwd.h"
0010 #include "FWCore/Framework/interface/MakerMacros.h"
0011 #include "FWCore/Framework/interface/SourceFactory.h"
0012 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0013 #include <cmath>
0014 #include <memory>
0015 
0016 #include "CondFormats/CSCObjects/interface/CSCDBChipSpeedCorrection.h"
0017 #include "CondFormats/DataRecord/interface/CSCDBChipSpeedCorrectionRcd.h"
0018 #include "DataFormats/MuonDetId/interface/CSCIndexer.h"
0019 #include <DataFormats/MuonDetId/interface/CSCDetId.h>
0020 
0021 class CSCChipSpeedCorrectionDBConditions : public edm::ESProducer, public edm::EventSetupRecordIntervalFinder {
0022 public:
0023   CSCChipSpeedCorrectionDBConditions(const edm::ParameterSet &);
0024   ~CSCChipSpeedCorrectionDBConditions() override;
0025 
0026   inline static CSCDBChipSpeedCorrection *prefillDBChipSpeedCorrection(bool isForMC,
0027                                                                        std::string dataCorrFileName,
0028                                                                        float dataOffse);
0029 
0030   typedef std::unique_ptr<CSCDBChipSpeedCorrection> ReturnType;
0031 
0032   ReturnType produceDBChipSpeedCorrection(const CSCDBChipSpeedCorrectionRcd &);
0033 
0034 private:
0035   // ----------member data ---------------------------
0036   void setIntervalFor(const edm::eventsetup::EventSetupRecordKey &,
0037                       const edm::IOVSyncValue &,
0038                       edm::ValidityInterval &) override;
0039   CSCDBChipSpeedCorrection *cndbChipCorr;
0040 
0041   // Flag for determining if this is for setting MC or data corrections
0042   bool isForMC;
0043   // File for reading <= 15768 data chip corrections.  MC will be fake (only one
0044   // value for every chip);
0045   std::string dataCorrFileName;
0046   float dataOffset;
0047 };
0048 
0049 #include "CondFormats/CSCObjects/interface/CSCDBChipSpeedCorrection.h"
0050 #include "CondFormats/DataRecord/interface/CSCDBChipSpeedCorrectionRcd.h"
0051 #include <DataFormats/MuonDetId/interface/CSCDetId.h>
0052 
0053 #include <fstream>
0054 #include <iostream>
0055 #include <vector>
0056 
0057 // to workaround plugin library
0058 inline CSCDBChipSpeedCorrection *CSCChipSpeedCorrectionDBConditions::prefillDBChipSpeedCorrection(bool isMC,
0059                                                                                                   std::string filename,
0060                                                                                                   float dataOffset) {
0061   if (isMC)
0062     printf("\n Generating fake DB constants for MC\n");
0063   else {
0064     printf("\n Reading chip corrections from file %s \n", filename.data());
0065     printf("my data offset value is %f \n", dataOffset);
0066   }
0067 
0068   CSCIndexer indexer;
0069 
0070   const int CHIP_FACTOR = 100;
0071   const int MAX_SIZE = 15768;
0072   const int MAX_SHORT = 32767;
0073   CSCDBChipSpeedCorrection *cndbChipCorr = new CSCDBChipSpeedCorrection();
0074 
0075   CSCDBChipSpeedCorrection::ChipSpeedContainer &itemvector = cndbChipCorr->chipSpeedCorr;
0076   itemvector.resize(MAX_SIZE);
0077   cndbChipCorr->factor_speedCorr = int(CHIP_FACTOR);
0078 
0079   // Fill chip corrections for MC is very simple
0080   if (isMC) {
0081     for (int i = 0; i < MAX_SIZE; i++) {
0082       itemvector[i].speedCorr = 0;
0083     }
0084     return cndbChipCorr;
0085   }
0086 
0087   // Filling for data takes a little more time
0088   FILE *fin = fopen(filename.data(), "r");
0089 
0090   int serialChamber, endcap, station, ring, chamber, chip;
0091   float t, dt;
0092   int nPulses;
0093 
0094   std::vector<int> new_index_id;
0095   std::vector<float> new_chipPulse;
0096   double runningTotal = 0;
0097   int numNonZero = 0;
0098   while (!feof(fin)) {
0099     // note space at end of format string to convert last \n
0100     // int check = fscanf(fin,"%d %d %f %f \n",&serialChamber,&chip,&t,&dt);
0101     int check = fscanf(fin,
0102                        "%d %d %d %d %d %d %f %f %d \n",
0103                        &serialChamber,
0104                        &endcap,
0105                        &station,
0106                        &ring,
0107                        &chamber,
0108                        &chip,
0109                        &t,
0110                        &dt,
0111                        &nPulses);
0112     if (check != 9) {
0113       printf("The input file format is not as expected\n");
0114       assert(0);
0115     }
0116 
0117     int serialChamber_safecopy = serialChamber;
0118     // Now to map from S. Durkin's serialChamber index convention
0119     // ME+1/1-ME+41 = 1  -234
0120     // ME+4/2       = 235-270*
0121     // ME-1/1-ME-41 = 271-504*
0122     // ME-4/2       = 505-540
0123     // To the convention used in /DataFormats/MuonDetId/interface/CSCIndexer.h
0124     // ME+1/1-ME+41 = 1  -234
0125     // ME+4/2       = 469-504*
0126     // ME-1/1-ME-41 = 235-468*
0127     // ME-4/2       = 505-540
0128     // Only the chambers marked * need to be remapped
0129 
0130     if (serialChamber >= 235 && serialChamber <= 270)
0131       serialChamber += 234;
0132     else {  // not in ME+4\2
0133       if (serialChamber >= 271 && serialChamber <= 504)
0134         serialChamber -= 36;
0135     }
0136 
0137     CSCDetId chamberId = indexer.detIdFromChamberIndex(serialChamber);
0138     // Now to map from S. Durkin's chip index convention 0-29 (with
0139     // 4,9,14,19,24,29 unused in ME1/3) To the convention used in
0140     // /DataFormats/MuonDetId/interface/CSCIndexer.h Layer 1-6 and chip 1-5 for
0141     // all chambers except ME1/3 which is chip 1-4
0142     int layer = (chip) / 5 + 1;
0143     CSCDetId cscId(chamberId.endcap(), chamberId.station(), chamberId.ring(), chamberId.chamber(), layer);
0144     // This should yield the same CSCDetId as decoding the chamber serial does
0145     // If not, some debugging is needed
0146     CSCDetId cscId_doubleCheck(endcap, station, ring, chamber, layer);
0147     if (cscId != cscId_doubleCheck) {
0148       printf("Why doesn't chamberSerial %d map to e: %d s: %d r: %d c: %d ? \n",
0149              serialChamber_safecopy,
0150              endcap,
0151              station,
0152              ring,
0153              chamber);
0154       assert(0);
0155     }
0156 
0157     chip = (chip % 5) + 1;
0158     // The file produced by Stan starts from the geometrical strip number stored
0159     // in the CSCStripDigi When the strip digis are built, the conversion from
0160     // electronics channel to geometrical strip (reversing ME+1/1a and, more
0161     // importantly, ME-1/1b) is done before the digi is filled
0162     // (EventFilter/CSCRawToDigi/src/CSCCFEBData.cc). Since we're filling an
0163     // electronics channel database, I'll flip again ME-1/1b
0164     if (endcap == 2 && station == 1 && ring == 1 && chip < 5) {
0165       chip = 5 - chip;  // change 1-4 to 4-1
0166     }
0167 
0168     new_index_id.push_back(indexer.chipIndex(cscId, chip));
0169     new_chipPulse.push_back(t);
0170     if (t != 0) {
0171       runningTotal += t;
0172       numNonZero++;
0173     }
0174   }
0175   fclose(fin);
0176 
0177   // Fill the chip corrections with zeros to start
0178   for (int i = 0; i < MAX_SIZE; i++) {
0179     itemvector[i].speedCorr = 0;
0180   }
0181 
0182   for (unsigned int i = 0; i < new_index_id.size(); i++) {
0183     if ((short int)(fabs((dataOffset - new_chipPulse[i]) * CHIP_FACTOR + 0.5)) < MAX_SHORT)
0184       itemvector[new_index_id[i] - 1].speedCorr =
0185           (short int)((dataOffset - new_chipPulse[i]) * CHIP_FACTOR + 0.5 * (dataOffset >= new_chipPulse[i]) -
0186                       0.5 * (dataOffset < new_chipPulse[i]));
0187     // printf("i= %d \t new index id = %d \t corr = %f \n",i,new_index_id[i],
0188     // new_chipPulse[i]);
0189   }
0190 
0191   // For now, calculate the mean chip correction and use it for all chambers
0192   // that don't have calibration pulse data (speedCorr ==0) or had values of
0193   // zero (speedCorr == dataOffset) This should be a temporary fix until all
0194   // chips that will read out in data have calibration information Since there
0195   // is only a handful out of 15K chips with values more than 3 ns away from the
0196   // average, this is probably very safe to first order
0197   float ave = runningTotal / numNonZero;
0198   for (int i = 0; i < MAX_SIZE; i++) {
0199     if (itemvector[i].speedCorr == 0 || itemvector[i].speedCorr == (short int)(dataOffset * CHIP_FACTOR + 0.5))
0200       itemvector[i].speedCorr =
0201           (short int)((dataOffset - ave) * CHIP_FACTOR + 0.5 * (dataOffset >= ave) - 0.5 * (dataOffset < ave));
0202   }
0203 
0204   return cndbChipCorr;
0205 }
0206 
0207 #endif