00001
00002
00003
00004
00005
00006
00007
00008
00009 #include <TFile.h>
00010 #include <TTree.h>
00011
00012 #include "EcalPerEvtLaserAnalyzer.h"
00013
00014 #include <sstream>
00015 #include <iostream>
00016 #include <iomanip>
00017 #include <ctime>
00018
00019 #include <FWCore/MessageLogger/interface/MessageLogger.h>
00020 #include <FWCore/Utilities/interface/Exception.h>
00021
00022 #include <CalibCalorimetry/EcalLaserAnalyzer/interface/MEEEGeom.h>
00023 #include <CalibCalorimetry/EcalLaserAnalyzer/interface/MEEBGeom.h>
00024
00025 #include <FWCore/Framework/interface/Event.h>
00026 #include <FWCore/Framework/interface/MakerMacros.h>
00027 #include <FWCore/ParameterSet/interface/ParameterSet.h>
00028 #include <FWCore/Framework/interface/EventSetup.h>
00029 #include <FWCore/Framework/interface/ESHandle.h>
00030
00031 #include <Geometry/EcalMapping/interface/EcalElectronicsMapping.h>
00032 #include <Geometry/EcalMapping/interface/EcalMappingRcd.h>
00033
00034 #include <DataFormats/EcalDigi/interface/EcalDigiCollections.h>
00035 #include <DataFormats/EcalDetId/interface/EcalElectronicsId.h>
00036 #include <DataFormats/EcalDetId/interface/EcalDetIdCollections.h>
00037 #include <DataFormats/EcalRawData/interface/EcalRawDataCollections.h>
00038
00039 #include <vector>
00040
00041 #include <CalibCalorimetry/EcalLaserAnalyzer/interface/TPNFit.h>
00042 #include <CalibCalorimetry/EcalLaserAnalyzer/interface/PulseFitWithFunction.h>
00043
00044 using namespace std;
00045
00046
00047
00048 EcalPerEvtLaserAnalyzer::EcalPerEvtLaserAnalyzer(const edm::ParameterSet& iConfig)
00049
00050 :
00051 iEvent(0),
00052
00053
00054 _nsamples( iConfig.getUntrackedParameter< unsigned int >( "nSamples", 10 ) ),
00055 _presample( iConfig.getUntrackedParameter< unsigned int >( "nPresamples", 3 ) ),
00056 _firstsample( iConfig.getUntrackedParameter< unsigned int >( "firstSample", 1 ) ),
00057 _lastsample( iConfig.getUntrackedParameter< unsigned int >( "lastSample", 2 ) ),
00058 _samplemin( iConfig.getUntrackedParameter< unsigned int >( "sampleMin", 3 ) ),
00059 _samplemax( iConfig.getUntrackedParameter< unsigned int >( "sampleMax", 9 ) ),
00060 _nsamplesPN( iConfig.getUntrackedParameter< unsigned int >( "nSamplesPN", 50 ) ),
00061 _presamplePN( iConfig.getUntrackedParameter< unsigned int >( "nPresamplesPN", 6 ) ),
00062 _firstsamplePN( iConfig.getUntrackedParameter< unsigned int >( "firstSamplePN", 7 ) ),
00063 _lastsamplePN( iConfig.getUntrackedParameter< unsigned int >( "lastSamplePN", 8 ) ),
00064 _timingcutlow( iConfig.getUntrackedParameter< unsigned int >( "timingCutLow", 3 ) ),
00065 _timingcuthigh( iConfig.getUntrackedParameter< unsigned int >( "timingCutHigh", 7 ) ),
00066 _niter( iConfig.getUntrackedParameter< unsigned int >( "nIter", 3 ) ),
00067 _fedid( iConfig.getUntrackedParameter< unsigned int >( "fedID", 0 ) ),
00068 _tower( iConfig.getUntrackedParameter< unsigned int >( "tower", 1 ) ),
00069 _channel( iConfig.getUntrackedParameter< unsigned int >( "channel", 1 ) ),
00070 _ecalPart( iConfig.getUntrackedParameter< std::string >( "ecalPart", "EB" ) ),
00071 nCrys( NCRYSEB),
00072 nTT( NTTEB),
00073 nSides( NSIDES),
00074 IsFileCreated(0), nCh(0),
00075 runType(-1), runNum(0), dccID(-1), lightside(2), doesRefFileExist(0),
00076 ttMat(-1), peakMat(-1), peak(-1), evtMat(-1), colMat(-1)
00077
00078
00079 {
00080
00081
00082
00083 resdir_ = iConfig.getUntrackedParameter<std::string>("resDir");
00084 refalphabeta_ = iConfig.getUntrackedParameter<std::string>("refAlphaBeta");
00085
00086 digiCollection_ = iConfig.getParameter<std::string>("digiCollection");
00087 digiPNCollection_ = iConfig.getParameter<std::string>("digiPNCollection");
00088 digiProducer_ = iConfig.getParameter<std::string>("digiProducer");
00089
00090 eventHeaderCollection_ = iConfig.getParameter<std::string>("eventHeaderCollection");
00091 eventHeaderProducer_ = iConfig.getParameter<std::string>("eventHeaderProducer");
00092
00093
00094
00095 if (_ecalPart == "EB") {
00096 nCrys = NCRYSEB;
00097 nTT = NTTEB;
00098 } else {
00099 nCrys = NCRYSEE;
00100 nTT = NTTEE;
00101 }
00102
00103 }
00104
00105
00106 EcalPerEvtLaserAnalyzer::~EcalPerEvtLaserAnalyzer(){
00107
00108
00109
00110
00111
00112
00113 }
00114
00115
00116
00117
00118 void EcalPerEvtLaserAnalyzer::beginJob() {
00119
00120
00121
00122
00123
00124 stringstream namefile1;
00125 namefile1 <<resdir_<< "/ADCSamples.root" ;
00126
00127 ADCfile=namefile1.str();
00128
00129
00130
00131 ADCFile = new TFile(ADCfile.c_str(),"RECREATE");
00132
00133
00134 stringstream name;
00135 name <<"ADCTree";
00136
00137 ADCtrees= new TTree(name.str().c_str(),name.str().c_str());
00138
00139 ADCtrees->Branch( "iphi", &phi, "phi/I" );
00140 ADCtrees->Branch( "ieta", &eta, "eta/I" );
00141 ADCtrees->Branch( "dccID", &dccID, "dccID/I" );
00142 ADCtrees->Branch( "event", &event, "event/I" );
00143 ADCtrees->Branch( "color", &color, "color/I" );
00144 ADCtrees->Branch( "adc", &adc , "adc[10]/D" );
00145 ADCtrees->Branch( "ttrigMatacq", &ttrig, "ttrig/D" );
00146 ADCtrees->Branch( "peakMatacq", &peak, "peak/D" );
00147 ADCtrees->Branch( "pn0", &pn0, "pn0/D" );
00148 ADCtrees->Branch( "pn1", &pn1, "pn1/D" );
00149
00150 ADCtrees->SetBranchAddress( "ieta", &eta );
00151 ADCtrees->SetBranchAddress( "iphi", &phi );
00152 ADCtrees->SetBranchAddress( "dccID", &dccID );
00153 ADCtrees->SetBranchAddress( "event", &event );
00154 ADCtrees->SetBranchAddress( "color", &color );
00155 ADCtrees->SetBranchAddress( "adc", adc );
00156 ADCtrees->SetBranchAddress( "ttrigMatacq", &ttrig );
00157 ADCtrees->SetBranchAddress( "peakMatacq", &peak );
00158 ADCtrees->SetBranchAddress( "pn0", &pn0 );
00159 ADCtrees->SetBranchAddress( "pn1", &pn1 );
00160
00161 IsFileCreated=0;
00162 nCh=nCrys;
00163
00164 }
00165
00166
00167
00168 void EcalPerEvtLaserAnalyzer:: analyze( const edm::Event & e, const edm::EventSetup& c){
00169
00170
00171 ++iEvent;
00172
00173
00174 edm::Handle<EcalRawDataCollection> pDCCHeader;
00175 const EcalRawDataCollection* DCCHeader=0;
00176 try {
00177 e.getByLabel(eventHeaderProducer_,eventHeaderCollection_, pDCCHeader);
00178 DCCHeader=pDCCHeader.product();
00179 }catch ( std::exception& ex ) {
00180 std::cerr << "Error! can't get the product retrieving DCC header" << eventHeaderCollection_.c_str() << std::endl;
00181 }
00182
00183
00184 edm::Handle<EBDigiCollection> pEBDigi;
00185 const EBDigiCollection* EBDigi=0;
00186 edm::Handle<EEDigiCollection> pEEDigi;
00187 const EEDigiCollection* EEDigi=0;
00188
00189 if (_ecalPart == "EB") {
00190 try {
00191 e.getByLabel(digiProducer_,digiCollection_, pEBDigi);
00192 EBDigi=pEBDigi.product();
00193 }catch ( std::exception& ex ) {
00194 std::cerr << "Error! can't get the product retrieving EB crystal data " << digiCollection_.c_str() << std::endl;
00195 }
00196 } else {
00197 try {
00198 e.getByLabel(digiProducer_,digiCollection_, pEEDigi);
00199 EEDigi=pEEDigi.product();
00200 }catch ( std::exception& ex ) {
00201 std::cerr << "Error! can't get the product retrieving EE crystal data " << digiCollection_.c_str() << std::endl;
00202 }
00203 }
00204
00205
00206
00207 edm::Handle<EcalPnDiodeDigiCollection> pPNDigi;
00208 const EcalPnDiodeDigiCollection* PNDigi=0;
00209 try {
00210 e.getByLabel(digiProducer_,digiPNCollection_, pPNDigi);
00211 PNDigi=pPNDigi.product();
00212 }catch ( std::exception& ex ) {
00213 std::cerr << "Error! can't get the product " << digiPNCollection_.c_str() << std::endl;
00214 }
00215
00216
00217 edm::ESHandle< EcalElectronicsMapping > ecalmapping;
00218 const EcalElectronicsMapping* TheMapping=0;
00219 try{
00220 c.get< EcalMappingRcd >().get(ecalmapping);
00221 TheMapping = ecalmapping.product();
00222 }catch ( std::exception& ex ) {
00223 std::cerr << "Error! can't get the product EcalMappingRcd"<< std::endl;
00224 }
00225
00226
00227
00228
00229
00230 for ( EcalRawDataCollection::const_iterator headerItr= DCCHeader->begin();headerItr != DCCHeader->end();
00231 ++headerItr ) {
00232
00233
00234
00235 int fed = headerItr->fedId();
00236
00237 if(fed!=_fedid && _fedid!=-999) continue;
00238
00239 runType=headerItr->getRunType();
00240 runNum=headerItr->getRunNumber();
00241 event=headerItr->getLV1();
00242 dccID=headerItr->getDccInTCCCommand();
00243 fedID=headerItr->fedId();
00244
00245
00246 if( 600+dccID != fedID ) continue;
00247
00248
00249 if(runType!=EcalDCCHeaderBlock::LASER_STD && runType!=EcalDCCHeaderBlock::LASER_GAP) return;
00250
00251
00252
00253 if (IsFileCreated==0){
00254
00255 stringstream namefile2;
00256 namefile2 << resdir_ <<"/APDAmpl_Run"<<runNum<<"_"<<_fedid <<"_"<<_tower<<"_"<<_channel<<".root";
00257 resfile=namefile2.str();
00258
00259
00260
00261 stringstream namefile;
00262 namefile << resdir_ <<"/Matacq-Run"<<runNum<<".root";
00263
00264 doesRefFileExist=0;
00265
00266 FILE *test;
00267 test = fopen(namefile.str().c_str(),"r");
00268 if (test) doesRefFileExist=1;
00269
00270 if(doesRefFileExist==1){
00271 matacqFile = new TFile((namefile.str().c_str()));
00272 matacqTree=(TTree*)matacqFile->Get("MatacqShape");
00273
00274 matacqTree->SetBranchAddress( "event", &evtMat );
00275 matacqTree->SetBranchAddress( "color", &colMat );
00276 matacqTree->SetBranchAddress( "peak", &peakMat );
00277 matacqTree->SetBranchAddress( "ttrig", &ttMat );
00278 }
00279
00280 IsFileCreated=1;
00281
00282 }
00283
00284
00285
00286 EcalDCCHeaderBlock::EcalDCCEventSettings settings = headerItr->getEventSettings();
00287 int color = settings.wavelength;
00288
00289 vector<int>::iterator iter = find( colors.begin(), colors.end(), color );
00290 if( iter==colors.end() ){
00291 colors.push_back( color );
00292 cout <<" new color found "<< color<<" "<< colors.size()<< endl;
00293 }
00294
00295 }
00296
00297
00298
00299 if(fedID!=_fedid && _fedid!=-999) return;
00300
00301
00302
00303
00304
00305
00306 TPNFit * pnfit = new TPNFit();
00307 pnfit -> init(_nsamplesPN,_firstsamplePN, _lastsamplePN);
00308
00309 double chi2pn=0;
00310 double ypnrange[50];
00311 double dsum=0.;
00312 double dsum1=0.;
00313 double bl=0.;
00314 double bl1=0.;
00315 double val_max=0.;
00316 unsigned int samplemax=0;
00317 unsigned int k;
00318
00319 std::vector<double> allPNAmpl;
00320
00321 for ( EcalPnDiodeDigiCollection::const_iterator pnItr = PNDigi->begin(); pnItr != PNDigi->end(); ++pnItr ) {
00322
00323
00324 for ( int samId=0; samId < (*pnItr).size() ; samId++ ) {
00325 pn[samId]=(*pnItr).sample(samId).adc();
00326 }
00327
00328
00329 for(dsum=0.,k=0;k<_presamplePN;k++) {
00330 dsum+=pn[k];
00331 }
00332 bl=dsum/((double) _presamplePN);
00333
00334 for(val_max=0.,k=0;k<_nsamplesPN;k++) {
00335 ypnrange[k]=pn[k] - bl;
00336
00337 if(ypnrange[k] > val_max) {
00338 val_max= ypnrange[k]; samplemax=k;
00339 }
00340 }
00341
00342 chi2pn = pnfit -> doFit(samplemax,&ypnrange[0]);
00343
00344 if(chi2pn == 101 || chi2pn == 102 || chi2pn == 103) pnAmpl=0.;
00345 else pnAmpl= pnfit -> getAmpl();
00346
00347 allPNAmpl.push_back(pnAmpl);
00348
00349 }
00350
00351
00352
00353
00354
00355 ttrig=-1.;
00356 peak=-1.;
00357
00358 if(doesRefFileExist==1){
00359
00360 if (color==0) matacqTree->GetEntry(event-1);
00361 else if(color==3) matacqTree->GetEntry(matacqTree->GetEntries("color==0")+event-1);
00362 ttrig=ttMat;
00363 peak=peakMat;
00364
00365 }
00366
00367
00368
00369
00370
00371
00372 double yrange[10];
00373 int adcGain=0;
00374 int side=0;
00375
00376 if (EBDigi){
00377
00378 for ( EBDigiCollection::const_iterator digiItr= EBDigi->begin(); digiItr != EBDigi->end(); ++digiItr ) {
00379
00380 EBDetId id_crystal(digiItr->id()) ;
00381 EBDataFrame df( *digiItr );
00382
00383 int etaG = id_crystal.ieta() ;
00384 int phiG = id_crystal.iphi() ;
00385
00386 int etaL ;
00387 int phiL ;
00388 std::pair<int, int> LocalCoord=MEEBGeom::localCoord( etaG , phiG );
00389
00390 etaL=LocalCoord.first ;
00391 phiL=LocalCoord.second ;
00392
00393 eta = etaG;
00394 phi = phiG;
00395
00396 side=MEEBGeom::side(etaG,phiG);
00397
00398 EcalElectronicsId elecid_crystal = TheMapping->getElectronicsId(id_crystal);
00399
00400 int towerID=elecid_crystal.towerId();
00401
00402 int strip=elecid_crystal.stripId();
00403 int xtal=elecid_crystal.xtalId();
00404 int channelID= 5*(strip-1) + xtal-1;
00405
00406 int module= MEEBGeom::lmmod(etaG, phiG);
00407
00408 std::pair<int,int> pnpair=MEEBGeom::pn(module);
00409 unsigned int MyPn0=pnpair.first;
00410 unsigned int MyPn1=pnpair.second;
00411
00412 unsigned int channel=MEEBGeom::electronic_channel( etaL, phiL );
00413 assert(channel < nCrys);
00414
00415 int iadcmax=0; double adcmax=0.0;
00416
00417 if(towerID!=int(_tower) || channelID!=int(_channel) || dccID!=int(_fedid-600)) continue;
00418 else channelNumber=channel;
00419
00420 for (unsigned int i=0; i< (*digiItr).size() ; ++i ) {
00421
00422 EcalMGPASample samp_crystal(df.sample(i));
00423 adc[i]=samp_crystal.adc() ;
00424 adcG[i]=samp_crystal.gainId();
00425
00426 if (i==0) adcGain=adcG[i];
00427 if (i>0) adcGain=TMath::Max(adcG[i],adcGain);
00428
00429 if (adc[i]>adcmax) {
00430 adcmax=adc[i];
00431 iadcmax=i;
00432 }
00433 }
00434
00435 for(dsum=0.,dsum1=0.,k=0;k<_presample;k++) {
00436 dsum+=adc[k];
00437 if(k<_presample-1) dsum1+=adc[k];
00438 }
00439
00440 bl=dsum/((double)_presample);
00441 bl1=dsum1/((double)_presample-1);
00442
00443 for(val_max=0.,k=0;k<_nsamples;k++) {
00444 yrange[k]=adc[k] - bl;
00445 if(yrange[k] > val_max) {
00446 val_max= yrange[k]; samplemax=k;
00447 }
00448 }
00449
00450 if(samplemax==4 || samplemax==5) {
00451 val_max=val_max+bl-bl1;
00452 for(k=0;k<_nsamples;k++) {
00453 yrange[k]=yrange[k]+bl-bl1;
00454 }
00455 }
00456
00457 for(unsigned int k=0;k<_nsamples;k++) {
00458 adc[k]=yrange[k];
00459 }
00460
00461 pn0=allPNAmpl[MyPn0];
00462 pn1=allPNAmpl[MyPn1];
00463
00464 if(samplemax>=_timingcutlow && samplemax<=_timingcuthigh && lightside==side) ADCtrees->Fill();
00465
00466 }
00467
00468 } else {
00469
00470 for ( EEDigiCollection::const_iterator digiItr= EEDigi->begin(); digiItr != EEDigi->end(); ++digiItr ) {
00471
00472 EEDetId id_crystal(digiItr->id()) ;
00473 EEDataFrame df( *digiItr );
00474
00475 phi = id_crystal.ix()-1 ;
00476 eta = id_crystal.iy()-1 ;
00477 side=0;
00478
00479
00480
00481 EcalElectronicsId elecid_crystal = TheMapping->getElectronicsId(id_crystal);
00482
00483 int towerID=elecid_crystal.towerId();
00484 int channelID=elecid_crystal.channelId()-1;
00485
00486 int module=MEEEGeom::lmmod( phi, eta );
00487
00488 std::pair<int,int> pnpair=MEEEGeom::pn( module, _fedid ) ;
00489 unsigned int MyPn0=pnpair.first;
00490 unsigned int MyPn1=pnpair.second;
00491
00492 unsigned int channel=MEEEGeom::crystal( phi, eta );
00493 assert(channel < nCrys);
00494
00495 int iadcmax=0; double adcmax=0.0;
00496
00497 if(towerID!=int(_tower) || channelID!=int(_channel) || dccID!=int(_fedid-600)) continue;
00498 else channelNumber=channel;
00499
00500 for (unsigned int i=0; i< (*digiItr).size() ; ++i ) {
00501
00502 EcalMGPASample samp_crystal(df.sample(i));
00503 adc[i]=samp_crystal.adc() ;
00504 adcG[i]=samp_crystal.gainId();
00505
00506 if (i==0) adcGain=adcG[i];
00507 if (i>0) adcGain=TMath::Max(adcG[i],adcGain);
00508
00509 if (adc[i]>adcmax) {
00510 adcmax=adc[i];
00511 iadcmax=i;
00512 }
00513 }
00514
00515 for(dsum=0.,dsum1=0.,k=0;k<_presample;k++) {
00516 dsum+=adc[k];
00517 if(k<_presample-1) dsum1+=adc[k];
00518 }
00519
00520 bl=dsum/((double)_presample);
00521 bl1=dsum1/((double)_presample-1);
00522
00523 for(val_max=0.,k=0;k<_nsamples;k++) {
00524 yrange[k]=adc[k] - bl;
00525 if(yrange[k] > val_max) {
00526 val_max= yrange[k]; samplemax=k;
00527 }
00528 }
00529
00530 if(samplemax==4 || samplemax==5) {
00531 val_max=val_max+bl-bl1;
00532 for(k=0;k<_nsamples;k++) {
00533 yrange[k]=yrange[k]+bl-bl1;
00534 }
00535 }
00536
00537 for(unsigned int k=0;k<_nsamples;k++) {
00538 adc[k]=yrange[k];
00539 }
00540
00541 pn0=allPNAmpl[MyPn0];
00542 pn1=allPNAmpl[MyPn1];
00543
00544 if(samplemax>=_timingcutlow && samplemax<=_timingcuthigh && lightside==side) ADCtrees->Fill();
00545
00546 }
00547 }
00548
00549 }
00550
00551
00552
00553 void EcalPerEvtLaserAnalyzer::endJob() {
00554
00555
00556
00557 assert( colors.size()<= nColor );
00558 unsigned int nCol=colors.size();
00559
00560
00561 ADCtrees->Write();
00562
00563
00564 cout << "\n\t+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+" << endl;
00565 cout << "\t+=+ Analyzing laser data: getting per event +=+" << endl;
00566 cout << "\t+=+ APD Amplitudes and ADC samples +=+"<< endl;
00567 cout << "\t+=+ for fed:" <<_fedid<<", tower:"<<_tower<<", and channel:" << _channel << endl;
00568
00569
00570
00571
00572 APDFile = new TFile(resfile.c_str(),"RECREATE");
00573
00574 int ieta, iphi, channelID, towerID, flag;
00575 double alpha, beta;
00576
00577
00578 colors.push_back( color );
00579
00580 for (unsigned int i=0;i<nCol;i++){
00581
00582 stringstream name1;
00583 name1<<"headerCol"<<colors[i];
00584
00585 header[i] = new TTree(name1.str().c_str(),name1.str().c_str());
00586
00587 header[i]->Branch( "alpha", &alpha, "alpha/D" );
00588 header[i]->Branch( "beta", &beta, "beta/D" );
00589 header[i]->Branch( "iphi", &iphi, "iphi/I" );
00590 header[i]->Branch( "ieta", &ieta, "ieta/I" );
00591 header[i]->Branch( "dccID", &dccID, "dccID/I" );
00592 header[i]->Branch( "towerID", &towerID, "towerID/I" );
00593 header[i]->Branch( "channelID", &channelID, "channelID/I" );
00594
00595 header[i]->SetBranchAddress( "alpha", &alpha );
00596 header[i]->SetBranchAddress( "beta", &beta );
00597 header[i]->SetBranchAddress( "iphi", &iphi );
00598 header[i]->SetBranchAddress( "ieta", &ieta );
00599 header[i]->SetBranchAddress( "dccID", &dccID );
00600 header[i]->SetBranchAddress( "towerID", &towerID );
00601 header[i]->SetBranchAddress( "channelID", &channelID );
00602
00603 }
00604
00605 stringstream name2;
00606 name2<<"APDTree";
00607 APDtrees= new TTree(name2.str().c_str(),name2.str().c_str());
00608
00609
00610
00611 APDtrees->Branch( "event", &event, "event/I" );
00612 APDtrees->Branch( "color", &color, "color/I" );
00613 APDtrees->Branch( "adc", &adc , "adc[10]/D" );
00614 APDtrees->Branch( "peakMatacq", &peak , "peak/D" );
00615 APDtrees->Branch( "ttrigMatacq", &ttrig , "ttrig/D" );
00616 APDtrees->Branch( "apdAmpl", &apdAmpl, "apdAmpl/D" );
00617 APDtrees->Branch( "apdTime", &apdTime, "apdTime/D" );
00618 APDtrees->Branch( "flag", &flag, "flag/I" );
00619 APDtrees->Branch( "pn0", &pn0, "pn0/D" );
00620 APDtrees->Branch( "pn1", &pn1, "pn1/D" );
00621
00622 APDtrees->SetBranchAddress( "event", &event );
00623 APDtrees->SetBranchAddress( "color", &color );
00624 APDtrees->SetBranchAddress( "adc", adc );
00625 APDtrees->SetBranchAddress( "peakMatacq", &peak );
00626 APDtrees->SetBranchAddress( "ttrigMatacq", &ttrig );
00627 APDtrees->SetBranchAddress( "apdAmpl", &apdAmpl );
00628 APDtrees->SetBranchAddress( "apdTime", &apdTime );
00629 APDtrees->SetBranchAddress( "flag", &flag );
00630 APDtrees->SetBranchAddress( "pn0", &pn0 );
00631 APDtrees->SetBranchAddress( "pn1", &pn1 );
00632
00633
00634
00635
00636 TFile *alphaFile = new TFile(refalphabeta_.c_str());
00637 TTree *alphaTree[2];
00638
00639 Double_t alphaRun, betaRun;
00640 int ietaRun, iphiRun, channelIDRun, towerIDRun, dccIDRun, flagRun;
00641
00642 for( unsigned int i=0;i<nCol;i++){
00643
00644 stringstream name3;
00645 name3<<"ABCol"<<i;
00646 alphaTree[i]=(TTree*)alphaFile->Get(name3.str().c_str());
00647 alphaTree[i]->SetBranchStatus( "*", 0 );
00648 alphaTree[i]->SetBranchStatus( "alpha", 1 );
00649 alphaTree[i]->SetBranchStatus( "beta", 1 );
00650 alphaTree[i]->SetBranchStatus( "iphi", 1 );
00651 alphaTree[i]->SetBranchStatus( "ieta", 1 );
00652 alphaTree[i]->SetBranchStatus( "dccID", 1 );
00653 alphaTree[i]->SetBranchStatus( "towerID", 1 );
00654 alphaTree[i]->SetBranchStatus( "channelID", 1 );
00655 alphaTree[i]->SetBranchStatus( "flag", 1 );
00656
00657 alphaTree[i]->SetBranchAddress( "alpha", &alphaRun );
00658 alphaTree[i]->SetBranchAddress( "beta", &betaRun );
00659 alphaTree[i]->SetBranchAddress( "iphi", &iphiRun );
00660 alphaTree[i]->SetBranchAddress( "ieta", &ietaRun );
00661 alphaTree[i]->SetBranchAddress( "dccID", &dccIDRun );
00662 alphaTree[i]->SetBranchAddress( "towerID", &towerIDRun );
00663 alphaTree[i]->SetBranchAddress( "channelID", &channelIDRun );
00664 alphaTree[i]->SetBranchAddress( "flag", &flagRun );
00665
00666 }
00667
00668 PulseFitWithFunction * pslsfit = new PulseFitWithFunction();
00669
00670 double chi2;
00671
00672 for (unsigned int icol=0;icol<nCol;icol++){
00673
00674 IsThereDataADC[icol]=1;
00675 stringstream cut;
00676 cut <<"color=="<<colors.at(icol);
00677 if(ADCtrees->GetEntries(cut.str().c_str())<10) IsThereDataADC[icol]=0;
00678 IsHeaderFilled[icol]=0;
00679
00680 }
00681
00682
00683
00684
00685 Long64_t nbytes = 0, nb = 0;
00686 for (Long64_t jentry=0; jentry< ADCtrees->GetEntriesFast();jentry++) {
00687 nb = ADCtrees->GetEntry(jentry); nbytes += nb;
00688
00689 int iCry=channelNumber;
00690
00691
00692
00693 unsigned int iCol=0;
00694 for(unsigned int i=0;i<nCol;i++){
00695 if(color==colors[i]) {
00696 iCol=i;
00697 i=colors.size();
00698 }
00699 }
00700
00701 alphaTree[iCol]->GetEntry(iCry);
00702
00703
00704 flag=flagRun;
00705 iphi=iphiRun;
00706 ieta=ietaRun;
00707 towerID=towerIDRun;
00708 channelID=channelIDRun;
00709 alpha=alphaRun;
00710 beta=betaRun;
00711
00712 if (IsHeaderFilled[iCol]==0){
00713 header[iCol]->Fill();
00714 IsHeaderFilled[iCol]=1;
00715
00716 }
00717
00718
00719 apdAmpl=0;
00720 apdTime=0;
00721
00722 pslsfit -> init(_nsamples,_firstsample,_lastsample,_niter, alphaRun, betaRun);
00723 chi2 = pslsfit -> doFit(&adc[0]);
00724
00725 if( chi2 < 0. || chi2 == 102 || chi2 == 101 ) {
00726
00727 apdAmpl=0;
00728 apdTime=0;
00729
00730 }else{
00731
00732 apdAmpl = pslsfit -> getAmpl();
00733 apdTime = pslsfit -> getTime();
00734
00735 }
00736
00737 APDtrees->Fill();
00738
00739 }
00740
00741 alphaFile->Close();
00742
00743 ADCFile->Close();
00744
00745 APDFile->Write();
00746 APDFile->Close();
00747
00748
00749
00750
00751
00752
00753 stringstream del;
00754 del << "rm " <<ADCfile;
00755 system(del.str().c_str());
00756
00757 cout << "\t+=+ .................................................. done +=+" << endl;
00758 cout << "\t+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+" << endl;
00759 }
00760
00761
00762
00763 DEFINE_FWK_MODULE(EcalPerEvtLaserAnalyzer);
00764