CMS 3D CMS Logo

CSCGainAnalyzer Class Reference

#include <OnlineDB/CSCCondDB/interface/CSCGainAnalyzer.h>

Inheritance diagram for CSCGainAnalyzer:

edm::EDAnalyzer

List of all members.

Public Member Functions

virtual void analyze (edm::Event const &e, edm::EventSetup const &iSetup)
 CSCGainAnalyzer (edm::ParameterSet const &conf)
 Analyzer for reading CSC pedestals.
 ~CSCGainAnalyzer ()

Private Attributes

TH2F adcCharge_ch0
TH2F adcCharge_ch1
TH2F adcCharge_ch10
TH2F adcCharge_ch11
TH2F adcCharge_ch12
TH2F adcCharge_ch2
TH2F adcCharge_ch3
TH2F adcCharge_ch4
TH2F adcCharge_ch5
TH2F adcCharge_ch6
TH2F adcCharge_ch7
TH2F adcCharge_ch8
TH2F adcCharge_ch9
float adcMax [DDU_ga][CHAMBERS_ga][LAYERS_ga][STRIPS_ga]
float adcMean_max [DDU_ga][CHAMBERS_ga][LAYERS_ga][STRIPS_ga]
std::string chamber_id
int chamber_index
int chamber_num
std::string chamber_type
int chamberIndex
int counterzero
int crateID [CHAMBERS_ga]
bool debug
int dmbID [CHAMBERS_ga]
int eventNumber
int evt
int fff
std::ifstream filein
int first_strip_index
int flagGain
int flagIntercept
int flagRun
float gainIntercept
float gainSlope
int i_chamber
int i_layer
int length
int lines
float maxmodten [NUMMODTEN_ga][CHAMBERS_ga][LAYERS_ga][STRIPS_ga]
int misMatch
float myCharge [20]
int myIndex
int myNcham
std::string name
int NChambers
int Nddu
std::vector< intnewadc
float newChi2 [480]
float newGain [480]
float newIntercept [480]
std::string PSet
time_t rawtime
int record
int reportedChambers
int ret_code
int sector
int size [CHAMBERS_ga]
int strip
int strips_per_layer


Detailed Description

Definition at line 33 of file CSCGainAnalyzer.h.


Constructor & Destructor Documentation

CSCGainAnalyzer::CSCGainAnalyzer ( edm::ParameterSet const &  conf  )  [explicit]

Analyzer for reading CSC pedestals.

author S. Durkin, O.Boeriu 18/03/06 ripped from Jeremy's and Rick's analyzersto be used for old mapping for new mapping

Definition at line 39 of file CSCGainAnalyzer.cc.

References adcCharge_ch0, adcCharge_ch1, adcCharge_ch10, adcCharge_ch11, adcCharge_ch12, adcCharge_ch2, adcCharge_ch3, adcCharge_ch4, adcCharge_ch5, adcCharge_ch6, adcCharge_ch7, adcCharge_ch8, adcCharge_ch9, adcMax, adcMean_max, CHAMBERS_ga, counterzero, DDU_ga, debug, eventNumber, evt, flagGain, flagIntercept, flagRun, gainIntercept, gainSlope, edm::ParameterSet::getUntrackedParameter(), i, i_chamber, i_layer, j, k, edm::es::l(), LAYERS_ga, length, maxmodten, misMatch, myIndex, myNcham, Nddu, newChi2, newGain, newIntercept, NUMMODTEN_ga, reportedChambers, size, strip, and STRIPS_ga.

00039                                                             {
00040   debug = conf.getUntrackedParameter<bool>("debug",false);
00041   eventNumber=0,evt=0,counterzero=0,Nddu=0,flagGain=-9,flagIntercept=-9,flagRun=-9;
00042   strip=0,misMatch=0,myIndex=0,myNcham=-999;
00043   i_chamber=0,i_layer=0,reportedChambers=0;
00044   length=1,gainSlope=-999.0,gainIntercept=-999.0;
00045   
00046   adcCharge_ch0 = TH2F("adcCharge Cham 0","adcCharge", 100,0,300,100,0,3000);
00047   adcCharge_ch1 = TH2F("adcCharge Cham 1","adcCharge", 100,0,300,100,0,3000);
00048   adcCharge_ch2 = TH2F("adcCharge Cham 2","adcCharge", 100,0,300,100,0,3000);
00049   adcCharge_ch3 = TH2F("adcCharge Cham 3","adcCharge", 100,0,300,100,0,3000);
00050   adcCharge_ch4 = TH2F("adcCharge Cham 4","adcCharge", 100,0,300,100,0,3000);
00051   adcCharge_ch5 = TH2F("adcCharge Cham 5","adcCharge", 100,0,300,100,0,3000);
00052   adcCharge_ch6 = TH2F("adcCharge Cham 6","adcCharge", 100,0,300,100,0,3000);
00053   adcCharge_ch7 = TH2F("adcCharge Cham 7","adcCharge", 100,0,300,100,0,3000);
00054   adcCharge_ch8 = TH2F("adcCharge Cham 8","adcCharge", 100,0,300,100,0,3000);
00055   adcCharge_ch9 = TH2F("adcCharge Cham 9","adcCharge", 100,0,300,100,0,3000);
00056   adcCharge_ch10 = TH2F("adcCharge Cham 10","adcCharge", 100,0,300,100,0,3000);
00057   adcCharge_ch11 = TH2F("adcCharge Cham 11","adcCharge", 100,0,300,100,0,3000);
00058   adcCharge_ch12 = TH2F("adcCharge Cham 12","adcCharge", 100,0,300,100,0,3000);
00059 
00060   for (int i=0; i<NUMMODTEN_ga; i++){
00061     for (int j=0; j<CHAMBERS_ga; j++){
00062       for (int k=0; k<LAYERS_ga; k++){
00063         for (int l=0; l<STRIPS_ga; l++){
00064           maxmodten[i][j][k][l] = 0.0;
00065         }
00066       }
00067     }
00068   }
00069   
00070 
00071   for (int i=0; i<CHAMBERS_ga; i++){
00072     size[i]  = 0;
00073   }
00074 
00075   for (int iii=0;iii<DDU_ga;iii++){
00076     for (int i=0; i<CHAMBERS_ga; i++){
00077       for (int j=0; j<LAYERS_ga; j++){
00078         for (int k=0; k<STRIPS_ga; k++){
00079           adcMax[iii][i][j][k]            = -999.0;
00080           adcMean_max[iii][i][j][k]       = -999.0;
00081         }
00082       }
00083     }
00084   }
00085 
00086   for (int i=0; i<480; i++){
00087     newGain[i]     =0.0;
00088     newIntercept[i]=0.0;
00089     newChi2[i]     =0.0;
00090   }
00091 }

CSCGainAnalyzer::~CSCGainAnalyzer (  ) 

old DB map

new DB mapping

old mapping

new mapping

Definition at line 188 of file CSCGainAnalyzer.cc.

References adcCharge_ch0, adcCharge_ch1, adcCharge_ch10, adcCharge_ch11, adcCharge_ch12, adcCharge_ch2, adcCharge_ch3, adcCharge_ch4, adcCharge_ch5, adcCharge_ch6, adcCharge_ch7, adcCharge_ch8, adcCharge_ch9, condbon::cdbon_last_record(), condbon::cdbon_write(), TCalibGainEvt::cham, chamber_id, chamber_index, chamber_num, chamber_type, CSCMapItem::MapItem::chamberId, chamberIndex, CSCMapItem::MapItem::chamberLabel, TCalibGainEvt::chi2, counter(), GenMuonPlsPt100GeV_cfg::cout, cscmap1::cratedmb(), crateID, CSCMapItem::MapItem::cscIndex, debug, dmbID, lat::endl(), filein, first_strip_index, FITNUMBERS_ga, flagGain, TCalibGainEvt::flagGain, TCalibGainEvt::flagIntercept, flagIntercept, TCalibGainEvt::flagRun, flagRun, gainIntercept, gainSlope, TCalibGainEvt::id, TCalibGainEvt::intercept, edm::isnan(), j, k, TCalibGainEvt::layer, LAYERS_ga, lines, python::multivaluedict::map(), maxmodten, myCharge, myIndex, myNcham, name, Nddu, NULL, NUMBERPLOTTED_ga, CSCobject::obj, out, record, CSCMapItem::MapItem::sector, sector, size, TCalibGainEvt::slope, TCalibGainEvt::strip, CSCMapItem::MapItem::stripIndex, CSCMapItem::MapItem::strips, STRIPS_ga, and strips_per_layer.

00188                                  {
00189   //get time of Run file for DB transfer
00190   filein.open("../test/CSCgain.cfg");
00191   filein.ignore(1000,'\n');
00192   
00193   while(filein != NULL){
00194     lines++;
00195     getline(filein,PSet);
00196     
00197     if (lines==2){
00198       name=PSet;  
00199       std::cout<<name<<std::endl;
00200     }
00201   }
00202   std::string::size_type runNameStart = name.find("\"",0);
00203   std::string::size_type runNameEnd   = name.find("raw",0);
00204   std::string::size_type rootStart    = name.find("Gains",0);
00205   int nameSize = runNameEnd+2-runNameStart;
00206   int myRootSize = rootStart-runNameStart+8;
00207   std::string myname= name.substr(runNameStart+1,nameSize);
00208   std::string myRootName= name.substr(runNameStart+1,myRootSize);
00209   std::string myRootEnd = ".root";
00210   std::string myASCIIFileEnd = ".dat";
00211   std::string runFile= myRootName;
00212   std::string myRootFileName = runFile+myRootEnd;
00213   std::string myASCIIFileName= runFile+myASCIIFileEnd;
00214   const char *myNewName=myRootFileName.c_str();
00215   const char *myFileName=myASCIIFileName.c_str();
00216   
00217   struct tm* clock;                         
00218   struct stat attrib;                       
00219   stat(myname.c_str(), &attrib);          
00220   clock = localtime(&(attrib.st_mtime));  
00221   std::string myTime=asctime(clock);
00222   std::ofstream myGainsFile(myFileName,std::ios::out);
00223   std::ofstream badstripsFile("badstrips.dat",std::ios::out);
00224   
00226   //cscmap *map = new cscmap();
00228   CSCMapItem::MapItem mapitem;
00229   cscmap1 *map = new cscmap1();
00230 
00231   CSCobject *cn = new CSCobject();
00232   condbon *dbon = new condbon();
00233 
00234   //root ntuple information
00235   TCalibGainEvt calib_evt;
00236   TFile calibfile(myNewName, "RECREATE");
00237   TTree calibtree("Calibration","Gains");
00238   calibtree.Branch("EVENT", &calib_evt, "slope/F:intercept/F:chi2/F:strip/I:layer/I:cham/I:id/I:flagGain/I:flagIntercept/I:flagRun/I");
00239 
00240   for (int dduiter=0;dduiter<Nddu;dduiter++){
00241      for(int chamberiter=0; chamberiter<myNcham; chamberiter++){
00242        for (int cham=0;cham<myNcham;cham++){
00243          if (cham !=chamberiter) continue;
00244  
00245         //get chamber ID from DB mapping        
00246         int new_crateID = crateID[cham];
00247         int new_dmbID   = dmbID[cham];
00248         int counter=0;
00249         std::cout<<" Crate: "<<new_crateID<<" and DMB:  "<<new_dmbID<<std::endl;
00250         //      myIndex=0;
00251 
00253         //map->crate_chamber(new_crateID,new_dmbID,&chamber_id,&chamber_num,&sector,&first_strip_index,&strips_per_layer,&chamber_index);
00255         map->cratedmb(new_crateID,new_dmbID,&mapitem);
00256         chamber_num=mapitem.chamberId;
00257         sector= mapitem.sector;
00258         first_strip_index=mapitem.stripIndex;
00259         strips_per_layer=mapitem.strips;
00260         chamber_index=mapitem.chamberId;
00261         chamber_type = mapitem.chamberLabel;
00262         chamberIndex = mapitem.cscIndex;
00263 
00264         std::cout<<"Data is for chamber:: "<<chamber_type<<"  "<<chamber_id<<" in sector:  "<<sector<<" index "<<first_strip_index<<std::endl;
00265         
00266         calib_evt.id=chamber_num;
00267 
00268         for (int layeriter=0; layeriter<LAYERS_ga; layeriter++){
00269           for (int stripiter=0; stripiter<STRIPS_ga; stripiter++){
00270               
00271             for (int j=0; j<LAYERS_ga; j++){//layer
00272               if (j != layeriter) continue;
00273               
00274               int layer_id=chamber_num+j+1;
00275 
00276               //if(sector==-100)continue;
00277               cn->obj[layer_id].resize(size[cham]);
00278               for (int k=0; k<size[cham]; k++){//strip
00279                 if (k != stripiter) continue;
00280                 float sumOfX    = 0.0;
00281                 float sumOfY    = 0.0;
00282                 float sumOfXY   = 0.0;
00283                 float sumx2     = 0.0;
00284                 float gainSlope = 0.0;
00285                 float gainIntercept = 0.0;
00286                 float chi2      = 0.0;
00287                 int bad_chan=0;
00288                 int bad=0;
00289                 int maxBadChan=0;
00290                 std::vector<int> BadChan;
00291                 int flag1=0;
00292                 int flag2=0;
00293                 int flag3=0;
00294                 int pointer=0;
00295                 
00296                 //now start at 0.1V and step 0.25 inbetween
00297                 float charge[NUMBERPLOTTED_ga]={11.2, 39.2, 67.2, 95.2, 123.2, 151.2, 179.2, 207.2, 235.2, 263.2, 291.2, 319.2, 347.2, 375.2, 403.2, 431.2, 459.2, 487.2, 515.2, 543.2};
00298                 //float charge[NUMBERPLOTTED_ga]={11.2, 39.2, 67.2, 95.2, 123.2, 151.2, 179.2, 207.2, 235.2, 263.2};
00299                 for(int ii=0; ii<FITNUMBERS_ga; ii++){//numbers    
00300                   sumOfX  += charge[ii];
00301                   sumOfY  += maxmodten[ii][cham][j][k];
00302                   sumOfXY += (charge[ii]*maxmodten[ii][cham][j][k]);
00303                   sumx2   += (charge[ii]*charge[ii]);
00304                   myCharge[ii] = 11.2 +(28.0*ii);
00305                   if(cham==0) adcCharge_ch0.Fill(myCharge[ii],maxmodten[ii][cham][j][k]);
00306                   if(cham==1) adcCharge_ch1.Fill(myCharge[ii],maxmodten[ii][cham][j][k]);
00307                   if(cham==2) adcCharge_ch2.Fill(myCharge[ii],maxmodten[ii][cham][j][k]);
00308                   if(cham==3) adcCharge_ch3.Fill(myCharge[ii],maxmodten[ii][cham][j][k]);
00309                   if(cham==4) adcCharge_ch4.Fill(myCharge[ii],maxmodten[ii][cham][j][k]);
00310                   if(cham==5) adcCharge_ch5.Fill(myCharge[ii],maxmodten[ii][cham][j][k]);
00311                   if(cham==6) adcCharge_ch6.Fill(myCharge[ii],maxmodten[ii][cham][j][k]);
00312                   if(cham==7) adcCharge_ch7.Fill(myCharge[ii],maxmodten[ii][cham][j][k]);
00313                   if(cham==8) adcCharge_ch8.Fill(myCharge[ii],maxmodten[ii][cham][j][k]);
00314                   if(cham==9) adcCharge_ch9.Fill(myCharge[ii],maxmodten[ii][cham][j][k]);
00315                   if(cham==10) adcCharge_ch10.Fill(myCharge[ii],maxmodten[ii][cham][j][k]);
00316                   if(cham==11) adcCharge_ch11.Fill(myCharge[ii],maxmodten[ii][cham][j][k]);
00317                   if(cham==12) adcCharge_ch12.Fill(myCharge[ii],maxmodten[ii][cham][j][k]);
00318                 }
00319         
00320                 //Fit parameters for straight line
00321                 gainSlope     = ((FITNUMBERS_ga*sumOfXY) - (sumOfX * sumOfY))/((FITNUMBERS_ga*sumx2) - (sumOfX*sumOfX));//k
00322                 gainIntercept = ((sumOfY*sumx2)-(sumOfX*sumOfXY))/((FITNUMBERS_ga*sumx2)-(sumOfX*sumOfX));//m
00323                 
00324                 //if (gainSlope <3.0)  gainSlope = 0.0;
00325                 if (isnan(gainSlope)) gainSlope = 0.0;
00326                 
00327                 for(int ii=0; ii<FITNUMBERS_ga; ii++){
00328                   chi2  += (maxmodten[ii][cham][j][k]-(gainIntercept+(gainSlope*charge[ii])))*(maxmodten[ii][cham][j][k]-(gainIntercept+(gainSlope*charge[ii])))/(FITNUMBERS_ga*FITNUMBERS_ga);
00329                 }
00330                               
00331                 if (gainSlope>5.0 && gainSlope<10.0) flagGain=1; // ok
00332                 if (gainSlope<5.0)                   flagGain=2; // warning fit fails
00333                 if (gainSlope>10.0)                  flagGain=3; // warning fit fails
00334 
00335                 if (gainIntercept> -40. && gainIntercept<15.)  flagIntercept = 1 ;
00336                 if (gainIntercept< -40.)                       flagIntercept = 2 ;
00337                 if (gainIntercept> 15.)                        flagIntercept = 3 ;  
00338 
00339                 //dump values to ASCII file
00340                 counter++; 
00341                 myIndex = first_strip_index+counter-1;
00342                 if (size[cham] != strips_per_layer) flagRun = 1; //bad run
00343                 if (size[cham] == strips_per_layer) flagRun = 0; //good run 
00344                 if (counter>size[cham]*LAYERS_ga) counter=0;
00345 
00346                 //BadStripChannels!!!!
00347                 myGainsFile <<"  "<<myIndex-1<<"  "<<gainSlope<<"  "<<gainIntercept<<"  "<<chi2<<std::endl;
00348                 if(gainSlope<4.0){
00349                   gainSlope=0.0;
00350                   bad +=bad_chan++;
00351                   flag1=1;
00352                 }
00353                 maxBadChan=bad;
00354                 for(bad=0;bad<maxBadChan;bad++){
00355                   BadChan[bad]=bad_chan;
00356                   pointer = BadChan[bad];
00357                 }
00358 
00359                 if(gainSlope>9.0) flag2=2;
00360                 
00361                 if(gainSlope<4.0 || gainSlope>9.0){
00362                   badstripsFile<<"  "<<chamberIndex<<"  "<<pointer<<"  "<<bad_chan<<"  "<<j<<"  "<<k<<"  "<<flag1<<"  "<<flag2<<"  "<<flag3<<std::endl;
00363                     }
00364                 
00365                 calib_evt.slope     = gainSlope;
00366                 calib_evt.intercept = gainIntercept;
00367                 calib_evt.chi2      = chi2;
00368                 calib_evt.strip     = k;
00369                 calib_evt.layer     = j;
00370                 calib_evt.cham      = cham;
00371                 calib_evt.flagGain  = flagGain;
00372                 calib_evt.flagIntercept  = flagIntercept;
00373                 calib_evt.flagRun   = flagRun;
00374                 
00375                 calibtree.Fill();
00376                 
00377                 cn->obj[layer_id][k].resize(3);
00378                 cn->obj[layer_id][k][0] = gainSlope;
00379                 cn->obj[layer_id][k][1] = gainIntercept;
00380                 cn->obj[layer_id][k][2] = chi2;
00381               }//k loop
00382             }//j loop
00383           }//stripiter loop
00384         }//layiter loop
00385       }//cham 
00386     }//chamberiter
00387   }//dduiter
00388     
00389   //send data to DB
00390   dbon->cdbon_last_record("gains",&record);
00391   std::cout<<"Last gains record "<<record<<" for run file "<<myname<<" saved "<<myTime<<std::endl;
00392   if(debug) dbon->cdbon_write(cn,"gains",12,3498,myTime);
00393   adcCharge_ch0.Write();
00394   adcCharge_ch1.Write();
00395   adcCharge_ch2.Write();
00396   adcCharge_ch3.Write();
00397   adcCharge_ch4.Write();
00398   adcCharge_ch5.Write();
00399   adcCharge_ch6.Write();
00400   adcCharge_ch7.Write();
00401   adcCharge_ch8.Write();
00402   adcCharge_ch9.Write();
00403   adcCharge_ch10.Write();
00404   adcCharge_ch11.Write();
00405   adcCharge_ch12.Write();
00406 
00407   calibfile.Write();
00408   calibfile.Close();
00409   
00410 }


Member Function Documentation

void CSCGainAnalyzer::analyze ( edm::Event const &  e,
edm::EventSetup const &  iSetup 
) [virtual]

Take a reference to this FED's data

unpack data

get a pointer to data and pass it to constructor for unpacking

loop over DDUs

get a reference to chamber data

Implements edm::EDAnalyzer.

Definition at line 93 of file CSCGainAnalyzer.cc.

References ecalMGPA::adc(), adcMax, adcMean_max, CSCDMBHeader::cfebAvailable(), CHAMBERS_ga, counterzero, GenMuonPlsPt100GeV_cfg::cout, crateID, FEDRawData::data(), DDU_ga, dmbID, lat::endl(), eventNumber, evt, edm::Event::getByLabel(), edm::Event::getByType(), FEDNumbering::getCSCFEDIds(), i, i_chamber, i_layer, id, int, k, kk, LAYERS_ga, maxmodten, misMatch, myNcham, NChambers, Nddu, NUMBERPLOTTED_ga, NUMMODTEN_ga, PULSES_ga, reportedChambers, FEDRawData::size(), size, strip, and STRIPS_ga.

00093                                                                             {
00094   
00095   edm::Handle<CSCStripDigiCollection> strips;
00096   e.getByLabel("cscunpacker","MuonCSCStripDigi",strips);
00097   
00098   edm::Handle<FEDRawDataCollection> rawdata;
00099   e.getByType(rawdata);
00100   counterzero=counterzero+1;
00101   evt=(counterzero+1)/2;
00102 
00103   for (int id=FEDNumbering::getCSCFEDIds().first;
00104        id<=FEDNumbering::getCSCFEDIds().second; ++id){ //for each of our DCCs    
00105     
00107     const FEDRawData& fedData = rawdata->FEDData(id);
00108     if (fedData.size()){ 
00109       
00111       CSCDCCEventData dccData((short unsigned int *) fedData.data()); 
00112       
00113       const std::vector<CSCDDUEventData> & dduData = dccData.dduData(); 
00114       //evt++;      
00115       for (unsigned int iDDU=0; iDDU<dduData.size(); ++iDDU) {  
00116 
00118         const std::vector<CSCEventData> & cscData = dduData[iDDU].cscData();
00119         
00120         Nddu = dduData.size();
00121         reportedChambers += dduData[iDDU].header().ncsc();
00122         NChambers = cscData.size();
00123         int repChambers = dduData[iDDU].header().ncsc();
00124         std::cout << " Reported Chambers = " << repChambers <<"   "<<NChambers<< std::endl;
00125         if (NChambers!=repChambers) { std::cout<< "misMatched size!!!" << std::endl; misMatch++;}
00126         //if(NChambers==0) continue;
00127         if(NChambers > myNcham){
00128           myNcham=NChambers;
00129         }
00130 
00131         float ped=0.0;
00132 
00133         for (int i_chamber=0; i_chamber<NChambers; i_chamber++) {   
00134           
00135           for (int i_layer = 1; i_layer <=6; ++i_layer) {
00136             
00137             std::vector<CSCStripDigi> digis = cscData[i_chamber].stripDigis(i_layer) ;
00138             const CSCDMBHeader * thisDMBheader = cscData[i_chamber].dmbHeader();
00139             
00140             if (cscData[i_chamber].dmbHeader() && thisDMBheader->cfebAvailable()){
00141               dmbID[i_chamber] = cscData[i_chamber].dmbHeader()->dmbID();//get DMB ID
00142               crateID[i_chamber] = cscData[i_chamber].dmbHeader()->crateID();//get crate ID
00143               if(crateID[i_chamber] == 255) continue;
00144               
00145               for (unsigned int i=0; i<digis.size(); i++){
00146                 size[i_chamber] = digis.size();
00147                 std::vector<int> adc = digis[i].getADCCounts();
00148                 strip = digis[i].getStrip();
00149                 adcMax[iDDU][i_chamber][i_layer-1][strip-1]=-99.0; 
00150                 for(unsigned int k=0; k<adc.size(); k++){
00151                   ped=(adc[0]+adc[1])/2.;
00152                   if(adc[k]-ped > adcMax[iDDU][i_chamber][i_layer-1][strip-1]) {
00153                     adcMax[iDDU][i_chamber][i_layer-1][strip-1]=adc[k]-ped;
00154                   }
00155                 }
00156                 adcMean_max[iDDU][i_chamber][i_layer-1][strip-1] += adcMax[iDDU][i_chamber][i_layer-1][strip-1]/PULSES_ga;  
00157                 
00158                 //Every 25 event save
00159                 if (evt%PULSES_ga == 0 && (strip-1)%16 == (evt-1)/NUMMODTEN_ga){
00160                   int pulsesNr = int((evt-1)/PULSES_ga)%NUMBERPLOTTED_ga ;
00161                   maxmodten[pulsesNr][i_chamber][i_layer-1][strip-1] = adcMean_max[iDDU][i_chamber][i_layer-1][strip-1];
00162                 }
00163               }//end digis loop
00164             }//end cfeb.available loop
00165           }//end layer loop
00166         }//end chamber loop
00167         
00168         //every 25 events reset
00169         if((evt-1)%PULSES_ga==0){
00170           for(int iii=0;iii<DDU_ga;iii++){
00171             for(int ii=0; ii<CHAMBERS_ga; ii++){
00172               for(int jj=0; jj<LAYERS_ga; jj++){
00173                 for(int kk=0; kk<STRIPS_ga; kk++){
00174                   adcMean_max[iii][ii][jj][kk]=0.0;
00175                 }
00176               }
00177             }
00178           }
00179         }
00180         eventNumber++;
00181         edm::LogInfo ("CSCGainAnalyzer")  <<"end of event number "<<eventNumber<<" and non-zero event "<<evt;
00182       }
00183     }
00184   }
00185 }


Member Data Documentation

TH2F CSCGainAnalyzer::adcCharge_ch0 [private]

Definition at line 72 of file CSCGainAnalyzer.h.

Referenced by CSCGainAnalyzer(), and ~CSCGainAnalyzer().

TH2F CSCGainAnalyzer::adcCharge_ch1 [private]

Definition at line 73 of file CSCGainAnalyzer.h.

Referenced by CSCGainAnalyzer(), and ~CSCGainAnalyzer().

TH2F CSCGainAnalyzer::adcCharge_ch10 [private]

Definition at line 82 of file CSCGainAnalyzer.h.

Referenced by CSCGainAnalyzer(), and ~CSCGainAnalyzer().

TH2F CSCGainAnalyzer::adcCharge_ch11 [private]

Definition at line 83 of file CSCGainAnalyzer.h.

Referenced by CSCGainAnalyzer(), and ~CSCGainAnalyzer().

TH2F CSCGainAnalyzer::adcCharge_ch12 [private]

Definition at line 84 of file CSCGainAnalyzer.h.

Referenced by CSCGainAnalyzer(), and ~CSCGainAnalyzer().

TH2F CSCGainAnalyzer::adcCharge_ch2 [private]

Definition at line 74 of file CSCGainAnalyzer.h.

Referenced by CSCGainAnalyzer(), and ~CSCGainAnalyzer().

TH2F CSCGainAnalyzer::adcCharge_ch3 [private]

Definition at line 75 of file CSCGainAnalyzer.h.

Referenced by CSCGainAnalyzer(), and ~CSCGainAnalyzer().

TH2F CSCGainAnalyzer::adcCharge_ch4 [private]

Definition at line 76 of file CSCGainAnalyzer.h.

Referenced by CSCGainAnalyzer(), and ~CSCGainAnalyzer().

TH2F CSCGainAnalyzer::adcCharge_ch5 [private]

Definition at line 77 of file CSCGainAnalyzer.h.

Referenced by CSCGainAnalyzer(), and ~CSCGainAnalyzer().

TH2F CSCGainAnalyzer::adcCharge_ch6 [private]

Definition at line 78 of file CSCGainAnalyzer.h.

Referenced by CSCGainAnalyzer(), and ~CSCGainAnalyzer().

TH2F CSCGainAnalyzer::adcCharge_ch7 [private]

Definition at line 79 of file CSCGainAnalyzer.h.

Referenced by CSCGainAnalyzer(), and ~CSCGainAnalyzer().

TH2F CSCGainAnalyzer::adcCharge_ch8 [private]

Definition at line 80 of file CSCGainAnalyzer.h.

Referenced by CSCGainAnalyzer(), and ~CSCGainAnalyzer().

TH2F CSCGainAnalyzer::adcCharge_ch9 [private]

Definition at line 81 of file CSCGainAnalyzer.h.

Referenced by CSCGainAnalyzer(), and ~CSCGainAnalyzer().

float CSCGainAnalyzer::adcMax[DDU_ga][CHAMBERS_ga][LAYERS_ga][STRIPS_ga] [private]

Definition at line 61 of file CSCGainAnalyzer.h.

Referenced by analyze(), and CSCGainAnalyzer().

float CSCGainAnalyzer::adcMean_max[DDU_ga][CHAMBERS_ga][LAYERS_ga][STRIPS_ga] [private]

Definition at line 62 of file CSCGainAnalyzer.h.

Referenced by analyze(), and CSCGainAnalyzer().

std::string CSCGainAnalyzer::chamber_id [private]

Definition at line 54 of file CSCGainAnalyzer.h.

Referenced by ~CSCGainAnalyzer().

int CSCGainAnalyzer::chamber_index [private]

Definition at line 55 of file CSCGainAnalyzer.h.

Referenced by ~CSCGainAnalyzer().

int CSCGainAnalyzer::chamber_num [private]

Definition at line 55 of file CSCGainAnalyzer.h.

Referenced by ~CSCGainAnalyzer().

std::string CSCGainAnalyzer::chamber_type [private]

Definition at line 70 of file CSCGainAnalyzer.h.

Referenced by ~CSCGainAnalyzer().

int CSCGainAnalyzer::chamberIndex [private]

Definition at line 56 of file CSCGainAnalyzer.h.

Referenced by ~CSCGainAnalyzer().

int CSCGainAnalyzer::counterzero [private]

Definition at line 55 of file CSCGainAnalyzer.h.

Referenced by analyze(), and CSCGainAnalyzer().

int CSCGainAnalyzer::crateID[CHAMBERS_ga] [private]

Definition at line 59 of file CSCGainAnalyzer.h.

Referenced by analyze(), and ~CSCGainAnalyzer().

bool CSCGainAnalyzer::debug [private]

Definition at line 71 of file CSCGainAnalyzer.h.

Referenced by CSCGainAnalyzer(), and ~CSCGainAnalyzer().

int CSCGainAnalyzer::dmbID[CHAMBERS_ga] [private]

Definition at line 59 of file CSCGainAnalyzer.h.

Referenced by analyze(), and ~CSCGainAnalyzer().

int CSCGainAnalyzer::eventNumber [private]

Definition at line 55 of file CSCGainAnalyzer.h.

Referenced by analyze(), and CSCGainAnalyzer().

int CSCGainAnalyzer::evt [private]

Definition at line 55 of file CSCGainAnalyzer.h.

Referenced by analyze(), and CSCGainAnalyzer().

int CSCGainAnalyzer::fff [private]

Definition at line 57 of file CSCGainAnalyzer.h.

std::ifstream CSCGainAnalyzer::filein [private]

Definition at line 69 of file CSCGainAnalyzer.h.

Referenced by ~CSCGainAnalyzer().

int CSCGainAnalyzer::first_strip_index [private]

Definition at line 55 of file CSCGainAnalyzer.h.

Referenced by ~CSCGainAnalyzer().

int CSCGainAnalyzer::flagGain [private]

Definition at line 68 of file CSCGainAnalyzer.h.

Referenced by CSCGainAnalyzer(), and ~CSCGainAnalyzer().

int CSCGainAnalyzer::flagIntercept [private]

Definition at line 68 of file CSCGainAnalyzer.h.

Referenced by CSCGainAnalyzer(), and ~CSCGainAnalyzer().

int CSCGainAnalyzer::flagRun [private]

Definition at line 68 of file CSCGainAnalyzer.h.

Referenced by CSCGainAnalyzer(), and ~CSCGainAnalyzer().

float CSCGainAnalyzer::gainIntercept [private]

Definition at line 60 of file CSCGainAnalyzer.h.

Referenced by CSCGainAnalyzer(), and ~CSCGainAnalyzer().

float CSCGainAnalyzer::gainSlope [private]

Definition at line 60 of file CSCGainAnalyzer.h.

Referenced by CSCGainAnalyzer(), and ~CSCGainAnalyzer().

int CSCGainAnalyzer::i_chamber [private]

Definition at line 55 of file CSCGainAnalyzer.h.

Referenced by analyze(), and CSCGainAnalyzer().

int CSCGainAnalyzer::i_layer [private]

Definition at line 55 of file CSCGainAnalyzer.h.

Referenced by analyze(), and CSCGainAnalyzer().

int CSCGainAnalyzer::length [private]

Definition at line 57 of file CSCGainAnalyzer.h.

Referenced by CSCGainAnalyzer().

int CSCGainAnalyzer::lines [private]

Definition at line 68 of file CSCGainAnalyzer.h.

Referenced by ~CSCGainAnalyzer().

float CSCGainAnalyzer::maxmodten[NUMMODTEN_ga][CHAMBERS_ga][LAYERS_ga][STRIPS_ga] [private]

Definition at line 63 of file CSCGainAnalyzer.h.

Referenced by analyze(), CSCGainAnalyzer(), and ~CSCGainAnalyzer().

int CSCGainAnalyzer::misMatch [private]

Definition at line 57 of file CSCGainAnalyzer.h.

Referenced by analyze(), and CSCGainAnalyzer().

float CSCGainAnalyzer::myCharge[20] [private]

Definition at line 67 of file CSCGainAnalyzer.h.

Referenced by ~CSCGainAnalyzer().

int CSCGainAnalyzer::myIndex [private]

Definition at line 68 of file CSCGainAnalyzer.h.

Referenced by CSCGainAnalyzer(), and ~CSCGainAnalyzer().

int CSCGainAnalyzer::myNcham [private]

Definition at line 57 of file CSCGainAnalyzer.h.

Referenced by analyze(), CSCGainAnalyzer(), and ~CSCGainAnalyzer().

std::string CSCGainAnalyzer::name [private]

Definition at line 70 of file CSCGainAnalyzer.h.

Referenced by ~CSCGainAnalyzer().

int CSCGainAnalyzer::NChambers [private]

Definition at line 57 of file CSCGainAnalyzer.h.

Referenced by analyze().

int CSCGainAnalyzer::Nddu [private]

Definition at line 57 of file CSCGainAnalyzer.h.

Referenced by analyze(), CSCGainAnalyzer(), and ~CSCGainAnalyzer().

std::vector<int> CSCGainAnalyzer::newadc [private]

Definition at line 53 of file CSCGainAnalyzer.h.

float CSCGainAnalyzer::newChi2[480] [private]

Definition at line 66 of file CSCGainAnalyzer.h.

Referenced by CSCGainAnalyzer().

float CSCGainAnalyzer::newGain[480] [private]

Definition at line 64 of file CSCGainAnalyzer.h.

Referenced by CSCGainAnalyzer().

float CSCGainAnalyzer::newIntercept[480] [private]

Definition at line 65 of file CSCGainAnalyzer.h.

Referenced by CSCGainAnalyzer().

std::string CSCGainAnalyzer::PSet [private]

Definition at line 70 of file CSCGainAnalyzer.h.

time_t CSCGainAnalyzer::rawtime [private]

Definition at line 58 of file CSCGainAnalyzer.h.

int CSCGainAnalyzer::record [private]

Definition at line 57 of file CSCGainAnalyzer.h.

Referenced by ~CSCGainAnalyzer().

int CSCGainAnalyzer::reportedChambers [private]

Definition at line 55 of file CSCGainAnalyzer.h.

Referenced by analyze(), and CSCGainAnalyzer().

int CSCGainAnalyzer::ret_code [private]

Definition at line 57 of file CSCGainAnalyzer.h.

int CSCGainAnalyzer::sector [private]

Definition at line 55 of file CSCGainAnalyzer.h.

Referenced by ~CSCGainAnalyzer().

int CSCGainAnalyzer::size[CHAMBERS_ga] [private]

Definition at line 59 of file CSCGainAnalyzer.h.

Referenced by analyze(), CSCGainAnalyzer(), and ~CSCGainAnalyzer().

int CSCGainAnalyzer::strip [private]

Definition at line 57 of file CSCGainAnalyzer.h.

Referenced by analyze(), and CSCGainAnalyzer().

int CSCGainAnalyzer::strips_per_layer [private]

Definition at line 55 of file CSCGainAnalyzer.h.

Referenced by ~CSCGainAnalyzer().


The documentation for this class was generated from the following files:
Generated on Tue Jun 9 18:17:18 2009 for CMSSW by  doxygen 1.5.4