Line Code
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399 400 401 402 403 404 405 406 407 408 409 410 411 412 413 414 415 416 417 418 419 420 421 422 423 424 425 426 427 428 429 430 431 432 433 434 435 436 437 438 439 440 441 442 443 444 445 446 447 448 449 450 451 452 453 454 455 456 457 458 459 460 461 462 463 464 465 466 467 468 469 470 471 472 473 474 475 476 477 478 479 480 481 482 483 484 485 486 487 488 489 490 491 492 493 494 495 496 497 498 499 500 501 502 503 504 505 506 507 508 509 510 511 512 513 514 515 516 517 518 519 520 521 522 523 524 525 526 527 528 529 530 531 532 533 534 535 536 537 538 539 540 541 542 543 544 545 546 547 548 549 550 551 552 553 554 555 556 557 558 559 560 561 562 563 564 565 566 567 568 569 570 571 572 573 574 575 576 577 578 579 580 581 582 583 584 585 586 587 588 589 590 591 592 593 594 595 596 597 598 599 600 601 602 603 604 605 606 607 608 609 610 611 612 613 614 615 616 617 618 619 620 621 622 623 624 625 626 627 628 629 630 631 632 633 634 635 636 637 638 639 640 641 642 643 644 645 646 647 648 649 650 651 652 653 654 655 656 657 658 659 660 661 662 663 664 665 666 667 668 669 670 671 672 673 674 675 676 677 678 679 680 681 682
/* 
 *  \class EcalMatacqAnalyzer
 *
 *  primary author: Gautier Hamel De Monchenault - CEA/Saclay
 *  author: Julie Malcles - CEA/Saclay
 */

#include <TFile.h>
#include <TTree.h>
#include <TProfile.h>
#include <TChain.h>
#include <vector>

#include <CalibCalorimetry/EcalLaserAnalyzer/plugins/EcalMatacqAnalyzer.h>

#include <sstream>
#include <iostream>
#include <iomanip>

#include <FWCore/MessageLogger/interface/MessageLogger.h>
#include <FWCore/Utilities/interface/Exception.h>

#include <DataFormats/EcalDetId/interface/EcalDetIdCollections.h>

#include <FWCore/Framework/interface/Event.h>
#include <FWCore/Framework/interface/MakerMacros.h>
#include <FWCore/ParameterSet/interface/ParameterSet.h>

#include <CalibCalorimetry/EcalLaserAnalyzer/interface/TMatacq.h>
#include <CalibCalorimetry/EcalLaserAnalyzer/interface/TMom.h>
#include <CalibCalorimetry/EcalLaserAnalyzer/interface/TMTQ.h>

//========================================================================
EcalMatacqAnalyzer::EcalMatacqAnalyzer(const edm::ParameterSet& iConfig)
    //========================================================================
    :

      iEvent(0),

      // framework parameters with default values

      _presample(iConfig.getUntrackedParameter<double>("nPresamples", 6.7)),
      _nsamplesaftmax(iConfig.getUntrackedParameter<unsigned int>("nSamplesAftMax", 80)),
      _noiseCut(iConfig.getUntrackedParameter<unsigned int>("noiseCut", 7)),
      _parabnbefmax(iConfig.getUntrackedParameter<unsigned int>("paraBeforeMax", 8)),
      _parabnaftmax(iConfig.getUntrackedParameter<unsigned int>("paraAfterMax", 7)),
      _thres(iConfig.getUntrackedParameter<unsigned int>("threshold", 10)),
      _lowlev(iConfig.getUntrackedParameter<unsigned int>("lowLevel", 20)),
      _highlev(iConfig.getUntrackedParameter<unsigned int>("highLevel", 80)),
      _nevlasers(iConfig.getUntrackedParameter<unsigned int>("nEventLaser", 600)),
      _timebefmax(iConfig.getUntrackedParameter<unsigned int>("timeBefMax", 100)),
      _timeaftmax(iConfig.getUntrackedParameter<unsigned int>("timeAftMax", 250)),
      _cutwindow(iConfig.getUntrackedParameter<double>("cutWindow", 0.1)),
      _nsamplesshape(iConfig.getUntrackedParameter<unsigned int>("nSamplesShape", 250)),
      _presampleshape(iConfig.getUntrackedParameter<unsigned int>("nPresamplesShape", 50)),
      _slide(iConfig.getUntrackedParameter<unsigned int>("nSlide", 100)),
      _fedid(iConfig.getUntrackedParameter<int>("fedID", -999)),
      _debug(iConfig.getUntrackedParameter<int>("debug", 0)),
      resdir_(iConfig.getUntrackedParameter<std::string>("resDir")),
      digiCollection_(iConfig.getParameter<std::string>("digiCollection")),
      digiProducer_(iConfig.getParameter<std::string>("digiProducer")),
      eventHeaderCollection_(iConfig.getParameter<std::string>("eventHeaderCollection")),
      eventHeaderProducer_(iConfig.getParameter<std::string>("eventHeaderProducer")),
      pmatToken_(consumes<EcalMatacqDigiCollection>(edm::InputTag(digiProducer_, digiCollection_))),
      dccToken_(consumes<EcalRawDataCollection>(edm::InputTag(eventHeaderProducer_, eventHeaderCollection_))),
      nSides(NSIDES),
      lightside(0),
      runType(-1),
      runNum(0),
      event(0),
      color(-1),
      maxsamp(0),
      nsamples(0),
      tt(0)

//========================================================================
{
  //now do what ever initialization is needed
}

//========================================================================
void EcalMatacqAnalyzer::beginJob() {
  //========================================================================

  // Define temporary file name

  sampfile = resdir_;
  sampfile += "/TmpTreeMatacqAnalyzer.root";

  sampFile = new TFile(sampfile.c_str(), "RECREATE");

  // declaration of the tree to fill

  tree = new TTree("MatacqTree", "MatacqTree");

  //List of branches

  tree->Branch("event", &event, "event/I");
  tree->Branch("color", &color, "color/I");
  tree->Branch("matacq", &matacq, "matacq[2560]/D");
  tree->Branch("nsamples", &nsamples, "nsamples/I");
  tree->Branch("maxsamp", &maxsamp, "maxsamp/I");
  tree->Branch("tt", &tt, "tt/D");
  tree->Branch("lightside", &lightside, "lightside/I");

  tree->SetBranchAddress("event", &event);
  tree->SetBranchAddress("color", &color);
  tree->SetBranchAddress("matacq", matacq);
  tree->SetBranchAddress("nsamples", &nsamples);
  tree->SetBranchAddress("maxsamp", &maxsamp);
  tree->SetBranchAddress("tt", &tt);
  tree->SetBranchAddress("lightside", &lightside);

  // Define output results files' names

  std::stringstream namefile;
  namefile << resdir_ << "/MATACQ.root";
  outfile = namefile.str();

  // Laser events counter
  laserEvents = 0;
  isThereMatacq = false;
}

//========================================================================
void EcalMatacqAnalyzer::analyze(const edm::Event& e, const edm::EventSetup& c) {
  //========================================================================

  ++iEvent;

  if (_debug == 1)
    edm::LogVerbatim("EcalMatacqAnalyzzer") << "-- debug test -- Entering Analyze -- event= " << iEvent;

  // retrieving MATACQ :
  const edm::Handle<EcalMatacqDigiCollection>& pmatacqDigi = e.getHandle(pmatToken_);
  const EcalMatacqDigiCollection* matacqDigi = (pmatacqDigi.isValid()) ? pmatacqDigi.product() : nullptr;
  if (pmatacqDigi.isValid()) {
    if (_debug == 1)
      edm::LogVerbatim("EcalMatacqAnalyzzer") << "-- debug test -- Matacq Digis Found -- ";
  } else {
    edm::LogError("EcalMatacqAnalyzzer") << "Error! can't get the product EcalMatacqDigi producer:"
                                         << digiProducer_.c_str() << " collection:" << digiCollection_.c_str();
    return;
  }

  // retrieving DCC header

  const edm::Handle<EcalRawDataCollection>& pDCCHeader = e.getHandle(dccToken_);
  const EcalRawDataCollection* DCCHeader = (pDCCHeader.isValid()) ? pDCCHeader.product() : nullptr;
  if (!pDCCHeader.isValid()) {
    edm::LogError("EcalMatacqAnalyzzer") << "Error! can't get the product EcalRawData producer:"
                                         << eventHeaderProducer_.c_str()
                                         << " collection:" << eventHeaderCollection_.c_str();
    return;
  }

  // ====================================
  // Decode Basic DCCHeader Information
  // ====================================

  if (_debug == 1)
    edm::LogVerbatim("EcalMatacqAnalyzzer") << "-- debug test -- Before header -- ";

  for (EcalRawDataCollection::const_iterator headerItr = DCCHeader->begin(); headerItr != DCCHeader->end();
       ++headerItr) {
    EcalDCCHeaderBlock::EcalDCCEventSettings settings = headerItr->getEventSettings();
    color = (int)settings.wavelength;
    if (color < 0)
      return;

    // Get run type and run number

    int fed = headerItr->fedId();

    if (fed != _fedid && _fedid != -999)
      continue;

    runType = headerItr->getRunType();
    runNum = headerItr->getRunNumber();
    event = headerItr->getLV1();

    if (_debug == 1)
      edm::LogVerbatim("EcalMatacqAnalyzzer")
          << "-- debug test -- runtype:" << runType << " event:" << event << " runNum:" << runNum;

    dccID = headerItr->getDccInTCCCommand();
    fedID = headerItr->fedId();
    lightside = headerItr->getRtHalf();

    //assert (lightside<2 && lightside>=0);

    if (lightside != 1 && lightside != 0) {
      edm::LogVerbatim("EcalMatacqAnalyzzer") << "Unexpected lightside: " << lightside << " for event " << iEvent;
      return;
    }
    if (_debug == 1) {
      edm::LogVerbatim("EcalMatacqAnalyzzer")
          << "-- debug test -- Inside header before fed cut -- color=" << color << ", dcc=" << dccID
          << ", fed=" << fedID << ",  lightside=" << lightside << ", runType=" << runType;
    }

    // take event only if the fed corresponds to the DCC in TCC
    if (600 + dccID != fedID)
      continue;

    if (_debug == 1) {
      edm::LogVerbatim("EcalMatacqAnalyzzer")
          << "-- debug test -- Inside header after fed cut -- color=" << color << ", dcc=" << dccID << ", fed=" << fedID
          << ",  lightside=" << lightside << ", runType=" << runType;
    }

    // Cut on runType

    if (runType != EcalDCCHeaderBlock::LASER_STD && runType != EcalDCCHeaderBlock::LASER_GAP &&
        runType != EcalDCCHeaderBlock::LASER_POWER_SCAN && runType != EcalDCCHeaderBlock::LASER_DELAY_SCAN)
      return;

    std::vector<int>::iterator iter = find(colors.begin(), colors.end(), color);
    if (iter == colors.end()) {
      colors.push_back(color);
      edm::LogVerbatim("EcalMatacqAnalyzzer") << " new color found " << color << " " << colors.size();
    }
  }

  if (_debug == 1)
    edm::LogVerbatim("EcalMatacqAnalyzzer") << "-- debug test -- Before digis -- Event:" << iEvent;

  // Count laser events
  laserEvents++;

  // ===========================
  // Decode Matacq Information
  // ===========================

  double max = 0;

  for (EcalMatacqDigiCollection::const_iterator it = matacqDigi->begin(); it != matacqDigi->end();
       ++it) {  // Loop on matacq channel

    //
    const EcalMatacqDigi& digis = *it;

    if (_debug == 1) {
      edm::LogVerbatim("EcalMatacqAnalyzzer") << "-- debug test -- Inside digis -- digi size=" << digis.size();
    }

    if (digis.size() == 0)
      continue;
    else
      isThereMatacq = true;

    max = 0;
    maxsamp = 0;
    nsamples = digis.size();
    tt = digis.tTrig();

    for (int i = 0; i < digis.size(); ++i) {  // Loop on matacq samples
      matacq[i] = -digis.adcCount(i);
      if (matacq[i] > max) {
        max = matacq[i];
        maxsamp = i;
      }
    }
    if (_debug == 1) {
      edm::LogVerbatim("EcalMatacqAnalyzzer")
          << "-- debug test -- Inside digis -- nsamples=" << nsamples << ", max=" << max;
    }
  }

  if (_debug == 1)
    edm::LogVerbatim("EcalMatacqAnalyzzer") << "-- debug test -- After digis -- Event: " << iEvent;
  tree->Fill();

}  // analyze

//========================================================================
void EcalMatacqAnalyzer::endJob() {
  // Don't do anything if there is no events
  if (!isThereMatacq) {
    edm::LogVerbatim("EcalMatacqAnalyzzer") << "\n\t+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+";
    edm::LogVerbatim("EcalMatacqAnalyzzer") << "\t+=+     WARNING! NO MATACQ        +=+";
    edm::LogVerbatim("EcalMatacqAnalyzzer") << "\t+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+";

    // Remove temporary file
    FILE* test;
    test = fopen(sampfile.c_str(), "r");
    if (test) {
      std::stringstream del2;
      del2 << "rm " << sampfile;
      system(del2.str().c_str());
    }
    return;
  }

  assert(colors.size() <= nColor);
  unsigned int nCol = colors.size();

  for (unsigned int iCol = 0; iCol < nCol; iCol++) {
    for (unsigned int iSide = 0; iSide < nSide; iSide++) {
      MTQ[iCol][iSide] = new TMTQ();
    }
  }

  outFile = new TFile(outfile.c_str(), "RECREATE");

  TProfile* shapeMat = new TProfile("shapeLaser", "shapeLaser", _nsamplesshape, -0.5, double(_nsamplesshape) - 0.5);
  TProfile* shapeMatTmp = new TProfile(
      "shapeLaserTmp", "shapeLaserTmp", _timeaftmax + _timebefmax, -0.5, double(_timeaftmax + _timebefmax) - 0.5);

  edm::LogVerbatim("EcalMatacqAnalyzzer") << "\n\t+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+";
  edm::LogVerbatim("EcalMatacqAnalyzzer") << "\t+=+     Analyzing MATACQ data     +=+";

  //
  // create output ntuple
  //

  mtqShape = new TTree("MatacqShape", "MatacqShape");

  // list of branches
  // keep Patrice's notations

  mtqShape->Branch("event", &event, "event/I");
  mtqShape->Branch("color", &color, "color/I");
  mtqShape->Branch("status", &status, "status/I");
  mtqShape->Branch("peak", &peak, "peak/D");
  mtqShape->Branch("sigma", &sigma, "sigma/D");
  mtqShape->Branch("fit", &fit, "fit/D");
  mtqShape->Branch("ampl", &ampl, "ampl/D");
  mtqShape->Branch("trise", &trise, "trise/D");
  mtqShape->Branch("fwhm", &fwhm, "fwhm/D");
  mtqShape->Branch("fw20", &fw20, "fw20/D");
  mtqShape->Branch("fw80", &fw80, "fw80/D");
  mtqShape->Branch("ped", &ped, "ped/D");
  mtqShape->Branch("pedsig", &pedsig, "pedsig/D");
  mtqShape->Branch("ttrig", &ttrig, "ttrig/D");
  mtqShape->Branch("sliding", &sliding, "sliding/D");

  mtqShape->SetBranchAddress("event", &event);
  mtqShape->SetBranchAddress("color", &color);
  mtqShape->SetBranchAddress("status", &status);
  mtqShape->SetBranchAddress("peak", &peak);
  mtqShape->SetBranchAddress("sigma", &sigma);
  mtqShape->SetBranchAddress("fit", &fit);
  mtqShape->SetBranchAddress("ampl", &ampl);
  mtqShape->SetBranchAddress("fwhm", &fwhm);
  mtqShape->SetBranchAddress("fw20", &fw20);
  mtqShape->SetBranchAddress("fw80", &fw80);
  mtqShape->SetBranchAddress("trise", &trise);
  mtqShape->SetBranchAddress("ped", &ped);
  mtqShape->SetBranchAddress("pedsig", &pedsig);
  mtqShape->SetBranchAddress("ttrig", &ttrig);
  mtqShape->SetBranchAddress("sliding", &sliding);

  unsigned int endsample;
  unsigned int presample;

  // loop over the entries of the tree
  TChain* fChain = (TChain*)tree;
  Long64_t nentries = fChain->GetEntriesFast();

  for (Long64_t jentry = 0; jentry < nentries; jentry++) {
    // load the event
    Long64_t ientry = fChain->LoadTree(jentry);
    if (ientry < 0)
      break;
    fChain->GetEntry(jentry);

    status = 0;
    peak = -1;
    sigma = 0;
    fit = -1;
    ampl = -1;
    trise = -1;
    ttrig = tt;
    fwhm = 0;
    fw20 = 0;
    fw80 = 0;
    ped = 0;
    pedsig = 0;
    sliding = 0;

    if (_debug == 1)
      edm::LogVerbatim("EcalMatacqAnalyzzer")
          << "-- debug test -- inside loop 1  -- jentry:" << jentry << " over nentries=" << nentries;

    // create the object for Matacq data analysis

    endsample = maxsamp + _nsamplesaftmax;
    presample = int(_presample * nsamples / 100.);
    TMatacq* mtq = new TMatacq(nsamples,
                               presample,
                               endsample,
                               _noiseCut,
                               _parabnbefmax,
                               _parabnaftmax,
                               _thres,
                               _lowlev,
                               _highlev,
                               _nevlasers,
                               _slide);

    if (_debug == 1)
      edm::LogVerbatim("EcalMatacqAnalyzzer") << "-- debug test -- inside loop 2  -- ";

    // analyse the Matacq data
    if (mtq->rawPulseAnalysis(nsamples, &matacq[0]) == 0) {
      status = 1;
      ped = mtq->getBaseLine();
      pedsig = mtq->getsigBaseLine();

      if (_debug == 1)
        edm::LogVerbatim("EcalMatacqAnalyzzer") << "-- debug test -- inside loop 3  -- ped:" << ped;
      if (mtq->findPeak() == 0) {
        peak = mtq->getTimpeak();
        sigma = mtq->getsigTimpeak();
      }
      if (_debug == 1)
        edm::LogVerbatim("EcalMatacqAnalyzzer") << "-- debug test -- inside loop 4  -- peak:" << peak;
      if (mtq->doFit() == 0) {
        fit = mtq->getTimax();
        ampl = mtq->getAmpl();
        fwhm = mtq->getFwhm();
        fw20 = mtq->getWidth20();
        fw80 = mtq->getWidth80();
        sliding = mtq->getSlide();
      }
      if (_debug == 1)
        edm::LogVerbatim("EcalMatacqAnalyzzer") << "-- debug test -- inside loop 4  -- ampl:" << ampl;
      if (mtq->compute_trise() == 0) {
        trise = mtq->getTrise();
      }
      if (_debug == 1)
        edm::LogVerbatim("EcalMatacqAnalyzzer") << "-- debug test -- inside loop 5  -- trise:" << trise;
    }

    if (_debug == 1)
      edm::LogVerbatim("EcalMatacqAnalyzzer") << "-- debug test -- inside loop 6  -- status:" << status;

    if (status == 1 && mtq->findPeak() == 0) {
      int firstS = int(peak - double(_timebefmax));
      int lastS = int(peak + double(_timeaftmax));

      // Fill histo if there are enough samples
      if (_debug == 1)
        edm::LogVerbatim("EcalMatacqAnalyzzer")
            << "-- debug test -- inside loop 7  -- firstS:" << firstS << ", nsamples:" << nsamples;

      if (firstS >= 0 && lastS <= nsamples) {
        for (int i = firstS; i < lastS; i++) {
          shapeMatTmp->Fill(double(i) - firstS, matacq[i]);
        }

      } else {  // else extrapolate

        int firstSBis;

        if (firstS < 0) {  // fill first bins with 0

          double thisped;
          thisped = (matacq[0] + matacq[1] + matacq[2] + matacq[4] + matacq[5]) / 5.0;

          for (int i = firstS; i < 0; i++) {
            shapeMatTmp->Fill(double(i) - firstS, thisped);
          }
          firstSBis = 0;

        } else {
          firstSBis = firstS;
        }

        if (lastS > nsamples) {
          for (int i = firstSBis; i < int(nsamples); i++) {
            shapeMatTmp->Fill(double(i) - firstS, matacq[i]);
          }

          //extrapolate with expo tail

          double expb = 0.998;
          double matacqval = expb * matacq[nsamples - 1];

          for (int i = nsamples; i < lastS; i++) {
            shapeMatTmp->Fill(double(i) - firstS, matacqval);
            matacqval *= expb;
          }

        } else {
          for (int i = firstSBis; i < lastS; i++) {
            shapeMatTmp->Fill(double(i) - firstS, matacq[i]);
          }
        }
      }
    }
    if (_debug == 1)
      edm::LogVerbatim("EcalMatacqAnalyzzer") << "-- debug test -- inside loop 8";

    // get back color

    int iCol = nCol;
    for (unsigned int i = 0; i < nCol; i++) {
      if (color == colors[i]) {
        iCol = i;
        i = nCol;
      }
    }
    if (_debug == 1)
      edm::LogVerbatim("EcalMatacqAnalyzzer")
          << "-- debug test -- inside loop 8bis color:" << color << " iCol:" << iCol << " nCol:" << nCol;

    // fill TMTQ

    if (status == 1 && mtq->findPeak() == 0 && mtq->doFit() == 0 && mtq->compute_trise() == 0)
      MTQ[iCol][lightside]->addEntry(peak, sigma, fit, ampl, trise, fwhm, fw20, fw80, ped, pedsig, sliding);

    // fill the output tree

    if (_debug == 1)
      edm::LogVerbatim("EcalMatacqAnalyzzer") << "-- debug test -- inside loop 9";
    mtqShape->Fill();

    // clean up
    delete mtq;
  }

  if (_debug == 1)
    edm::LogVerbatim("EcalMatacqAnalyzzer") << "-- debug test -- after loop ";
  sampFile->Close();

  double Peak[6], Sigma[6], Fit[6], Ampl[6], Trise[6], Fwhm[6], Fw20[6], Fw80[6], Ped[6], Pedsig[6], Sliding[6];
  int Side;

  for (unsigned int iColor = 0; iColor < nCol; iColor++) {
    std::stringstream nametree;
    nametree << "MatacqCol" << colors[iColor];
    meanTree[iColor] = new TTree(nametree.str().c_str(), nametree.str().c_str());
    meanTree[iColor]->Branch("side", &Side, "Side/I");
    meanTree[iColor]->Branch("peak", &Peak, "Peak[6]/D");
    meanTree[iColor]->Branch("sigma", &Sigma, "Sigma[6]/D");
    meanTree[iColor]->Branch("fit", &Fit, "Fit[6]/D");
    meanTree[iColor]->Branch("ampl", &Ampl, "Ampl[6]/D");
    meanTree[iColor]->Branch("trise", &Trise, "Trise[6]/D");
    meanTree[iColor]->Branch("fwhm", &Fwhm, "Fwhm[6]/D");
    meanTree[iColor]->Branch("fw20", &Fw20, "Fw20[6]/D");
    meanTree[iColor]->Branch("fw80", &Fw80, "Fw80[6]/D");
    meanTree[iColor]->Branch("ped", &Ped, "Ped[6]/D");
    meanTree[iColor]->Branch("pedsig", &Pedsig, "Pedsig[6]/D");
    meanTree[iColor]->Branch("sliding", &Sliding, "Sliding[6]/D");

    meanTree[iColor]->SetBranchAddress("side", &Side);
    meanTree[iColor]->SetBranchAddress("peak", Peak);
    meanTree[iColor]->SetBranchAddress("sigma", Sigma);
    meanTree[iColor]->SetBranchAddress("fit", Fit);
    meanTree[iColor]->SetBranchAddress("ampl", Ampl);
    meanTree[iColor]->SetBranchAddress("fwhm", Fwhm);
    meanTree[iColor]->SetBranchAddress("fw20", Fw20);
    meanTree[iColor]->SetBranchAddress("fw80", Fw80);
    meanTree[iColor]->SetBranchAddress("trise", Trise);
    meanTree[iColor]->SetBranchAddress("ped", Ped);
    meanTree[iColor]->SetBranchAddress("pedsig", Pedsig);
    meanTree[iColor]->SetBranchAddress("sliding", Sliding);
  }

  for (unsigned int iCol = 0; iCol < nCol; iCol++) {
    for (unsigned int iSide = 0; iSide < nSides; iSide++) {
      Side = iSide;
      std::vector<double> val[TMTQ::nOutVar];

      for (int iVar = 0; iVar < TMTQ::nOutVar; iVar++) {
        val[iVar] = MTQ[iCol][iSide]->get(iVar);

        for (unsigned int i = 0; i < val[iVar].size(); i++) {
          switch (iVar) {
            case TMTQ::iPeak:
              Peak[i] = val[iVar].at(i);
              break;
            case TMTQ::iSigma:
              Sigma[i] = val[iVar].at(i);
              break;
            case TMTQ::iFit:
              Fit[i] = val[iVar].at(i);
              break;
            case TMTQ::iAmpl:
              Ampl[i] = val[iVar].at(i);
              break;
            case TMTQ::iFwhm:
              Fwhm[i] = val[iVar].at(i);
              break;
            case TMTQ::iFw20:
              Fw20[i] = val[iVar].at(i);
              break;
            case TMTQ::iFw80:
              Fw80[i] = val[iVar].at(i);
              break;
            case TMTQ::iTrise:
              Trise[i] = val[iVar].at(i);
              break;
            case TMTQ::iPed:
              Ped[i] = val[iVar].at(i);
              break;
            case TMTQ::iPedsig:
              Pedsig[i] = val[iVar].at(i);
              break;
            case TMTQ::iSlide:
              Sliding[i] = val[iVar].at(i);
              break;
          }
        }
      }
      meanTree[iCol]->Fill();
      if (_debug == 1)
        edm::LogVerbatim("EcalMatacqAnalyzzer") << "-- debug test -- inside final loop  ";
    }
  }

  // Calculate maximum with pol 2

  int im = shapeMatTmp->GetMaximumBin();
  double q1 = shapeMatTmp->GetBinContent(im - 1);
  double q2 = shapeMatTmp->GetBinContent(im);
  double q3 = shapeMatTmp->GetBinContent(im + 1);

  double a2 = (q3 + q1) / 2.0 - q2;
  double a1 = q2 - q1 + a2 * (1 - 2 * im);
  double a0 = q2 - a1 * im - a2 * im * im;

  double am = a0 - a1 * a1 / (4 * a2);

  // Compute pedestal

  double bl = 0;
  for (unsigned int i = 1; i < _presampleshape + 1; i++) {
    bl += shapeMatTmp->GetBinContent(i);
  }
  bl /= _presampleshape;

  // Compute and save laser shape

  if (_debug == 1)
    edm::LogVerbatim("EcalMatacqAnalyzzer") << "-- debug test -- computing shape  ";

  int firstBin = 0;
  double height = 0.0;

  for (unsigned int i = _timebefmax; i > _presampleshape; i--) {
    height = shapeMatTmp->GetBinContent(i) - bl;

    if (height < (am - bl) * _cutwindow) {
      firstBin = i;
      i = _presampleshape;
    }
  }

  unsigned int lastBin = firstBin + _nsamplesshape;

  for (unsigned int i = firstBin; i < lastBin; i++) {
    shapeMat->Fill(i - firstBin, shapeMatTmp->GetBinContent(i) - bl);
  }

  mtqShape->Write();
  for (unsigned int iColor = 0; iColor < nCol; iColor++) {
    meanTree[iColor]->Write();
  }
  if (_debug == 1)
    edm::LogVerbatim("EcalMatacqAnalyzzer") << "-- debug test -- writing  ";
  shapeMat->Write();

  // close the output file
  outFile->Close();

  // Remove temporary file
  FILE* test;
  test = fopen(sampfile.c_str(), "r");
  if (test) {
    std::stringstream del2;
    del2 << "rm " << sampfile;
    system(del2.str().c_str());
  }

  edm::LogVerbatim("EcalMatacqAnalyzzer") << "\t+=+    .................... done  +=+";
  edm::LogVerbatim("EcalMatacqAnalyzzer") << "\t+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+";
}

DEFINE_FWK_MODULE(EcalMatacqAnalyzer);