Skip to content
Snippets Groups Projects
Cbm2021EventBuilderAlgo.cxx 41.9 KiB
Newer Older
/********************************************************************************
 *    Copyright (C) 2014 GSI Helmholtzzentrum fuer Schwerionenforschung GmbH    *
 *                                                                              *
 *              This software is distributed under the terms of the             *
 *              GNU Lesser General Public Licence (LGPL) version 3,             *
 *                  copied verbatim in the file "LICENSE"                       *
 ********************************************************************************/
#include "Cbm2021EventBuilderAlgo.h"

/// CBM headers
#include "CbmEvent.h"
#include "CbmMuchBeamTimeDigi.h"
#include "CbmMuchDigi.h"
#include "CbmPsdDigi.h"
#include "CbmRichDigi.h"
#include "CbmStsDigi.h"
#include "CbmTofDigi.h"
#include "CbmTrdDigi.h"
#include "TimesliceMetaData.h"

/// FAIRROOT headers
#include "FairLogger.h"
#include "FairRootManager.h"
#include "FairRunOnline.h"

/// FAIRSOFT headers (geant, boost, ...)
#include "TCanvas.h"
#include "TClonesArray.h"
#include "TH1.h"
#include "TH2.h"
#include "THttpServer.h"

/// C/C++ headers

// ---- Default constructor --------------------------------------------
Cbm2021EventBuilderAlgo::Cbm2021EventBuilderAlgo() {}

// ---- Destructor -----------------------------------------------------
Cbm2021EventBuilderAlgo::~Cbm2021EventBuilderAlgo() {}

// ---- Init -----------------------------------------------------------
Bool_t Cbm2021EventBuilderAlgo::InitAlgo() {
  LOG(info)
    << "Cbm2021EventBuilderAlgo::InitAlgo => Starting sequence";

  // Get a handle from the IO manager
  FairRootManager* ioman = FairRootManager::Instance();

  /// Check if reference detector data are available
  if (kFALSE == CheckDataAvailable(fRefDet)) {
    LOG(fatal) << "No digi input for reference detector, stopping there!";
  }  // if( kFALSE == CheckDataAvailable( fRefDet ) )

  /// Check if data for detectors in selection list are available
  for (std::vector<EventBuilderDetector>::iterator det = fvDets.begin();
       det != fvDets.end();
       ++det) {
    if (kFALSE == CheckDataAvailable(*det)) {
      LOG(fatal)
        << "No digi input for one of selection detector, stopping there!";
    }  // if( kFALSE == CheckDataAvailable( *det ) )
  }  // for (std::vector<EventBuilderDetector>::iterator det = fvDets.begin(); det != fvDets.end(); ++det)

  /// Access the TS metadata to know TS start time if needed
  if (fdTsStartTime < 0 || fdTsLength < 0 || fdTsOverLength < 0) {
    fTimeSliceMetaDataArray =
      dynamic_cast<TClonesArray*>(ioman->GetObject("TimesliceMetaData"));
    if (!fTimeSliceMetaDataArray) {
      LOG(fatal)
        << "No TS metadata input found"
        << " => Please check in the unpacking macro if the following line was "
           "present!"
        << std::endl
        << "source->SetWriteOutputFlag(kTRUE);  // For writing TS metadata";
    }  // if (!fTimeSliceMetaDataArray)
  }    // if ( fdTsStartTime < 0 || fdTsLength < 0 || fdTsOverLength < 0 )

  if (fbFillHistos) { CreateHistograms(); }  // if( fbFillHistos )

  LOG(info) << "Cbm2021EventBuilderAlgo::InitAlgo => Done";

  return kTRUE;
}

// ---- ProcessTs ------------------------------------------------------
void Cbm2021EventBuilderAlgo::ProcessTs() {
  LOG_IF(info, fuNrTs % 1000 == 0) << "Begin of TS " << fuNrTs;

  InitTs();

  BuildEvents();

  /// Store last event with trigger if not done by other seed
  if (nullptr != fCurrentEvent) {
    /// TODO: store start time of current event ?
    //        fCurrentEvent->SetStartTime( fPrevTime ); // Replace Seed time with time of first digi in event?
    fCurrentEvent->SetEndTime(fdPrevEvtEndTime);
    fEventVector.push_back(fCurrentEvent);
    fuCurEv++;

    /// Prevent building over TS edge
    fCurrentEvent = nullptr;
  }  // if( nullptr != fCurrentEvent )

  LOG(debug) << "Found " << fEventVector.size() << " triggered events";

  if (fbFillHistos) { FillHistos(); }  // if( fbFillHistos )

  fuNrTs++;
}
void Cbm2021EventBuilderAlgo::ClearEventVector() {
  /// Need to delete the object the pointer points to first
  int counter = 0;
  for (CbmEvent* event : fEventVector) {
    LOG(debug) << "Event " << counter << " has " << event->GetNofData()
               << " digis";
    delete event;
    counter++;
  }  // for( CbmEvent* event: fEventVector)

  fEventVector.clear();
}
// ---- Finish ---------------------------------------------------------
void Cbm2021EventBuilderAlgo::Finish() {
  LOG(info) << "Total errors: " << fuErrors;
}

// ---------------------------------------------------------------------
Bool_t Cbm2021EventBuilderAlgo::CheckDataAvailable(
  EventBuilderDetector& det) {
  // Get a handle from the IO manager
  FairRootManager* ioman = FairRootManager::Instance();

  if (ECbmModuleId::kT0 == det.detId) {
    // T0 is not included in DigiManager
    fT0DigiVec = ioman->InitObjectAs<std::vector<CbmTofDigi> const*>("T0Digi");
    if (!fT0DigiVec) {
      LOG(info) << "No T0 digi input found.";
      return kFALSE;
    }    
  } 
  else if (ECbmModuleId::kSts == det.detId) {
    if (!fStsDigis) {
      LOG(info) << "No " << det.sName << " digi input found.";
      return kFALSE;
    }
  }
  else if (ECbmModuleId::kMuch == det.detId) {
    if (!fMuchDigis && !fMuchBeamTimeDigis ) {
      LOG(info) << "No " << det.sName << " digi input found.";
      return kFALSE;
    }
  }
  else if (ECbmModuleId::kTrd == det.detId) {
    if (!fTrdDigis) {
      LOG(info) << "No " << det.sName << " digi input found.";
      return kFALSE;
    }
  }
  else if (ECbmModuleId::kTof == det.detId) {
    if (!fTofDigis) {
      LOG(info) << "No " << det.sName << " digi input found.";
      return kFALSE;
    }
  }
  else if (ECbmModuleId::kRich == det.detId) {
    if (!fRichDigis) {
      LOG(info) << "No " << det.sName << " digi input found.";
      return kFALSE;
    }
  }
  else if (ECbmModuleId::kPsd == det.detId) {
    if (!fPsdDigis) {
      LOG(info) << "No " << det.sName << " digi input found.";
      return kFALSE;
    }
  }
  return kTRUE;
}
// ---------------------------------------------------------------------
void Cbm2021EventBuilderAlgo::InitTs() {
  /// Reset TS based variables (analysis per TS = no building over the border)
  /// Reference detector
  fRefDet.fuStartIndex = 0;
  fRefDet.fuEndIndex   = 0;
  /// Loop on detectors in selection list
  for (std::vector<EventBuilderDetector>::iterator det = fvDets.begin();
       det != fvDets.end();
       ++det) {
    (*det).fuStartIndex = 0;
    (*det).fuEndIndex   = 0;
  }  // for (std::vector<EventBuilderDetector>::iterator det = fvDets.begin(); det != fvDets.end(); ++det)
}

void Cbm2021EventBuilderAlgo::BuildEvents() {
  /// Call LoopOnSeed with proper template argument
  switch (fRefDet.detId) {
    case ECbmModuleId::kSts: {
      LoopOnSeeds<CbmStsDigi>();
      break;
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 683 684 685 686 687 688 689 690 691 692 693 694 695 696 697 698 699 700 701 702 703 704 705 706 707 708 709 710 711 712 713 714 715 716 717 718 719 720 721 722 723 724 725 726 727 728 729 730 731 732 733 734 735 736 737 738 739 740 741 742 743 744 745 746 747 748 749 750 751 752 753 754 755 756 757 758 759 760 761 762 763 764 765 766 767 768 769 770 771 772 773 774 775 776 777 778 779 780 781 782 783 784 785 786 787 788 789 790 791 792 793 794 795 796 797 798 799 800 801 802 803 804 805 806 807 808 809 810 811 812 813 814 815 816 817 818 819 820 821 822 823 824 825 826 827 828 829 830 831 832 833 834 835 836 837 838 839 840 841 842 843 844 845 846 847 848 849 850 851 852 853 854 855 856 857 858 859 860 861 862 863 864 865 866 867 868 869 870 871 872 873 874 875 876 877 878 879 880 881 882 883 884 885 886 887 888 889 890 891 892 893 894 895 896 897 898 899 900 901 902 903 904 905 906 907 908 909 910 911 912 913 914 915 916 917 918 919 920 921 922 923 924 925 926 927 928 929 930 931 932 933 934 935 936 937 938 939 940 941 942 943 944 945 946 947 948 949 950 951 952 953 954 955 956 957 958 959 960 961 962 963 964 965 966 967 968 969 970 971 972 973 974 975 976 977 978 979 980 981 982 983 984 985 986 987 988 989 990 991 992 993 994 995 996 997 998 999 1000 1001 1002 1003 1004 1005 1006 1007 1008 1009 1010 1011 1012 1013 1014 1015 1016 1017 1018 1019 1020 1021 1022 1023 1024 1025 1026 1027 1028 1029 1030 1031 1032 1033 1034 1035 1036 1037 1038 1039 1040 1041 1042 1043 1044 1045 1046 1047 1048 1049 1050 1051 1052 1053 1054 1055 1056 1057 1058 1059 1060 1061 1062 1063 1064 1065 1066 1067 1068 1069 1070 1071 1072 1073 1074 1075 1076 1077 1078 1079 1080 1081 1082 1083 1084 1085 1086 1087 1088 1089 1090 1091 1092 1093 1094 1095 1096 1097 1098 1099 1100 1101 1102 1103 1104 1105 1106 1107
    }  // case ECbmModuleId::kSts:
    case ECbmModuleId::kMuch: {
      if (fbUseMuchBeamtimeDigi) {
        LoopOnSeeds<CbmMuchBeamTimeDigi>();
      }  // if (fbUseMuchBeamtimeDigi)
      else {
        LoopOnSeeds<CbmMuchDigi>();
      }  // else of if (fbUseMuchBeamtimeDigi)
      break;
    }  // case ECbmModuleId::kMuch:
    case ECbmModuleId::kTrd: {
      LoopOnSeeds<CbmTrdDigi>();
      break;
    }  // case ECbmModuleId::kTrd:
    case ECbmModuleId::kTof: {
      LoopOnSeeds<CbmTofDigi>();
      break;
    }  // case ECbmModuleId::kTof:
    case ECbmModuleId::kRich: {
      LoopOnSeeds<CbmRichDigi>();
      break;
    }  // case ECbmModuleId::kRich:
    case ECbmModuleId::kPsd: {
      LoopOnSeeds<CbmPsdDigi>();
      break;
    }  // case ECbmModuleId::kPsd:
    case ECbmModuleId::kT0: {
      LoopOnSeeds<CbmTofDigi>();
      break;
    }  // case ECbmModuleId::kT0:
    default: {
      LOG(fatal) << "Cbm2021EventBuilderAlgo::BuildEvents => "
                 << "Trying to search event seeds with unsupported det: "
                 << fRefDet.sName;
      break;
    }  // default:
  }    // switch( *det )
}

UInt_t Cbm2021EventBuilderAlgo::GetNofDigis( ECbmModuleId detId ) {
  switch (detId) {
    case ECbmModuleId::kSts: {
	return fStsDigis->size();
    }  
    case ECbmModuleId::kMuch: {
      if (fbUseMuchBeamtimeDigi) {
	return fMuchBeamTimeDigis->size();
      }  
      else {
	return fMuchDigis->size();
      }  
    }  
    case ECbmModuleId::kTrd: {
      return fTrdDigis->size();
    }  
    case ECbmModuleId::kTof: {
      return fTofDigis->size();
    }  
    case ECbmModuleId::kRich: {
      return fRichDigis->size();
    }  
    case ECbmModuleId::kPsd: {
      return fPsdDigis->size();
    }  
    case ECbmModuleId::kT0: {
      return fT0DigiVec->size(); //what to do here? Not in digi manager.
    }  
    default: {
      LOG(fatal) << "Cbm2021EventBuilderAlgo::GetNofDigis => "
                 << "Trying to get digi number with unsupported detector.";
      return -1;
    }  
  }    
}

bool Cbm2021EventBuilderAlgo::DetIsPresent( ECbmModuleId detId ) {
  switch (detId) {
    case ECbmModuleId::kSts: {
	return fStsDigis!=nullptr;
    }  
    case ECbmModuleId::kMuch: {
      if (fbUseMuchBeamtimeDigi) {
	return fMuchBeamTimeDigis!=nullptr;
      }  
      else {
	return fMuchDigis!=nullptr;
      }  
    }  
    case ECbmModuleId::kTrd: {
      return fTrdDigis!=nullptr;
    }  
    case ECbmModuleId::kTof: {
      return fTofDigis!=nullptr;
    }  
    case ECbmModuleId::kRich: {
      return fRichDigis!=nullptr;
    }  
    case ECbmModuleId::kPsd: {
      return fPsdDigis!=nullptr;
    }  
    case ECbmModuleId::kT0: {
      return fT0DigiVec!=nullptr; //what to do here? Not in digi manager.
    }  
    default: {
      LOG(fatal) << "Cbm2021EventBuilderAlgo::GetNofDigis => "
                 << "Trying to get digi number with unsupported detector.";
      return -1;
    }  
  }    
}

template<> const CbmStsDigi* Cbm2021EventBuilderAlgo::GetDigi( UInt_t uDigi ) {
	return &((*fStsDigis)[uDigi]);
}
template<> const CbmMuchBeamTimeDigi* Cbm2021EventBuilderAlgo::GetDigi( UInt_t uDigi ) {
	return &((*fMuchBeamTimeDigis)[uDigi]);
}
template<> const CbmMuchDigi* Cbm2021EventBuilderAlgo::GetDigi( UInt_t uDigi ) {
	return &((*fMuchDigis)[uDigi]);
}
template<> const CbmTrdDigi* Cbm2021EventBuilderAlgo::GetDigi( UInt_t uDigi ) {
	return &((*fTrdDigis)[uDigi]);
}
template<> const CbmTofDigi* Cbm2021EventBuilderAlgo::GetDigi( UInt_t uDigi ) {
	return &((*fTofDigis)[uDigi]);
}
template<> const CbmRichDigi* Cbm2021EventBuilderAlgo::GetDigi( UInt_t uDigi ) {
	return &((*fRichDigis)[uDigi]);
}
template<> const CbmPsdDigi* Cbm2021EventBuilderAlgo::GetDigi( UInt_t uDigi ) {
	return &((*fPsdDigis)[uDigi]);
}

template<class DigiSeed>
void Cbm2021EventBuilderAlgo::LoopOnSeeds() {
  /// Access the TS metadata if needed to know TS start time and overlap size
  Double_t dTsStartTime  = fdTsStartTime;
  Double_t dOverlapStart = fdTsStartTime + fdTsLength;
  Double_t dOverlapSize  = fdTsOverLength;
  if (fdTsStartTime < 0 || fdTsLength < 0 || fdTsOverLength < 0) {
    pTsMetaData =
      dynamic_cast<TimesliceMetaData*>(fTimeSliceMetaDataArray->At(0));
    if (nullptr == pTsMetaData)
      LOG(fatal) << Form("Cbm2021EventBuilderAlgo::LoopOnSeeds => "
                         "No TS metadata found for TS %6u.",
                         fuNrTs);

    dTsStartTime  = pTsMetaData->GetStartTime();
    dOverlapStart = pTsMetaData->GetOverlapStartTime();
    dOverlapSize  = pTsMetaData->GetOverlapDuration();
  }  // if ( fdTsStartTime < 0 || fdTsLength < 0  || fdTsOverLength < 0 )

  /// Print warning in first TS if time window borders out of potential overlap
  if ((0.0 < fdEarliestTimeWinBeg && dOverlapSize < fdLatestTimeWinEnd)
      || (dOverlapSize < fdWidestTimeWinRange)) {
    LOG(warning) << "Cbm2021EventBuilderAlgo::LoopOnSeeds => "
                 << Form("Event window not fitting in TS overlap, risk of "
                         "incomplete events: %f %f %f %f",
                         fdEarliestTimeWinBeg,
                         fdLatestTimeWinEnd,
                         fdWidestTimeWinRange,
                         dOverlapSize);
  }  // if end of event window does not fit in overlap for a seed at edge of TS core

  /// Define an acceptance window for the seeds in order to use the overlap
  /// part of the TS to avoid incomplete events
  Double_t dSeedWindowBeg =
    dTsStartTime + (0.0 < fdEarliestTimeWinBeg ? 0.0 : -fdEarliestTimeWinBeg);
  Double_t dSeedWindowEnd =
    dOverlapStart + (0.0 < fdEarliestTimeWinBeg ? 0.0 : -fdEarliestTimeWinBeg);
  if (fbIgnoreTsOverlap) {
    dSeedWindowBeg = dTsStartTime;
    dSeedWindowEnd = dOverlapStart;
  }  // if( fbIgnoreTsOverlap )

  if (ECbmModuleId::kT0 == fRefDet.detId) {
    if (fT0DigiVec) {
      /// Loop on size of vector
      UInt_t uNbRefDigis = fT0DigiVec->size();
      /// Loop on size of vector
      for (UInt_t uDigi = 0; uDigi < uNbRefDigis; ++uDigi) {
        LOG(debug) << Form("Checking seed %6u / %6u", uDigi, uNbRefDigis);

        Double_t dTime = fT0DigiVec->at(uDigi).GetTime();

        /// Check Seed and build event if needed
        CheckSeed(dTime, uDigi);
      }  // for( UInt_t uDigi = 0; uDigi < uNbRefDigis; ++uDigi )
    }    // if ( fT0DigiVec )
    else
      LOG(fatal) << "Cbm2021EventBuilderAlgo::LoopOnSeeds => "
                 << "T0 as reference detector but vector not found!";
  }  // if (ECbmModuleId::kT0 == fRefDet.detId)
  else {
    UInt_t uNbRefDigis = (0 < GetNofDigis(fRefDet.detId)
                            ? GetNofDigis(fRefDet.detId)
                            : 0);
    /// Loop on size of vector
    for (UInt_t uDigi = 0; uDigi < uNbRefDigis; ++uDigi) {
      LOG(debug) << Form("Checking seed %6u / %6u", uDigi, uNbRefDigis);
      const DigiSeed* pDigi = GetDigi<DigiSeed>(uDigi);
      /// Check that _entry is not out of range
      if (nullptr != pDigi) {
        Double_t dTime = pDigi->GetTime();

        /// Check if seed in acceptance window
        if (dTime < dSeedWindowBeg) {
          continue;
        }  // if( dTime < dSeedWindowBeg )
        else if (dSeedWindowEnd < dTime) {
          break;
        }  // else if( dSeedWindowEnd < dTime )

        /// Check Seed and build event if needed
        CheckSeed(dTime, uDigi);
      }  // if( nullptr != pDigi )
    }    // for( UInt_t uDigi = 0; uDigi < uNbRefDigis; ++uDigi )
  }  // else of if (ECbmModuleId::kT0 == fRefDet.detId) => Digi containers controlled by DigiManager
}

void Cbm2021EventBuilderAlgo::CheckSeed(Double_t dSeedTime,
                                                   UInt_t uSeedDigiIdx) {
  /// If previous event valid and event overlap not allowed, check if we are in overlap
  /// and react accordingly
  if (nullptr != fCurrentEvent
      && (EOverlapMode::AllowOverlap != fOverMode
          || dSeedTime - fdPrevEvtTime < fRefDet.GetTimeWinRange())
      && dSeedTime - fdPrevEvtTime < fdWidestTimeWinRange) {
    /// Within overlap range
    switch (fOverMode) {
      case EOverlapMode::NoOverlap: {
        /// No overlap allowed => reject
        LOG(debug1) << "Reject seed due to overlap";
        return;
        break;
      }  // case EOverlapMode::NoOverlap:
      case EOverlapMode::MergeOverlap: {
        /// Merge overlap mode => do nothing and go on filling current event
        break;
      }  // case EOverlapMode::MergeOverlap:
      case EOverlapMode::AllowOverlap: {
        /// In allow overlap mode => reject only if reference det is in overlap
        /// to avoid cloning events due to single seed cluster
        LOG(debug1) << "Reject seed because part of cluster of previous one";
        return;
        break;
      }  // case EOverlapMode::AllowOverlap:
    }    // switch( fOverMode )
  }      // if( prev Event exists and mode forbiden overlap present )
  else {
    /// Out of overlap range or in overlap allowed mode
    /// => store previous event if not empty and create new one
    if (nullptr != fCurrentEvent) {
      /// TODO: store start time of current event ?
      //        fCurrentEvent->SetStartTime( fPrevTime ); // Replace Seed time with time of first digi in event?
      fCurrentEvent->SetEndTime(fdPrevEvtEndTime);
      fEventVector.push_back(fCurrentEvent);
      fuCurEv++;
    }  // if( nullptr != fCurrentEvent )
    fCurrentEvent = new CbmEvent(fuCurEv, dSeedTime, 0.);
  }  // else of if( prev Event exists and mode forbiden overlap present )

  /// If window open for reference detector, search for other reference Digis matching it
  /// Otherwise only add the current seed
  if (fRefDet.fdTimeWinBeg < fRefDet.fdTimeWinEnd) {
    switch (fRefDet.detId) {
      case ECbmModuleId::kSts: {
        SearchMatches<CbmStsDigi>(dSeedTime, fRefDet);
        break;
      }  // case ECbmModuleId::kSts:
      case ECbmModuleId::kMuch: {
        if (fbUseMuchBeamtimeDigi) {
          SearchMatches<CbmMuchBeamTimeDigi>(dSeedTime, fRefDet);
        }  // if (fbUseMuchBeamtimeDigi)
        else {
          SearchMatches<CbmMuchDigi>(dSeedTime, fRefDet);
        }  // else of if (fbUseMuchBeamtimeDigi)
        break;
      }  // case ECbmModuleId::kMuch:
      case ECbmModuleId::kTrd: {
        SearchMatches<CbmTrdDigi>(dSeedTime, fRefDet);
        break;
      }  // case ECbmModuleId::kTrd:
      case ECbmModuleId::kTof: {
        SearchMatches<CbmTofDigi>(dSeedTime, fRefDet);
        break;
      }  // case ECbmModuleId::kTof:
      case ECbmModuleId::kRich: {
        SearchMatches<CbmRichDigi>(dSeedTime, fRefDet);
        break;
      }  // case ECbmModuleId::kRich:
      case ECbmModuleId::kPsd: {
        SearchMatches<CbmPsdDigi>(dSeedTime, fRefDet);
        break;
      }  // case ECbmModuleId::kPsd:
      case ECbmModuleId::kT0: {
        SearchMatches<CbmTofDigi>(dSeedTime, fRefDet);
        break;
      }  // case ECbmModuleId::kT0:
      default: {
        LOG(fatal) << "Cbm2021EventBuilderAlgo::LoopOnSeeds => "
                   << "Trying to search matches with unsupported det: "
                   << fRefDet.sName << std::endl
                   << "You may want to add support for it in the method.";
        break;
      }  // default:
    }    // switch( fRefDet )

    /// Also add the seed if the window starts after the seed
    if (0 < fRefDet.fdTimeWinBeg) AddDigiToEvent(fRefDet, uSeedDigiIdx);
  }  // if( fdRefTimeWinBeg < fdRefTimeWinEnd )
  else {
    AddDigiToEvent(fRefDet, uSeedDigiIdx);
  }  // else of if( fdRefTimeWinBeg < fdRefTimeWinEnd )

  /// Search for matches for each detector in selection list
  for (std::vector<EventBuilderDetector>::iterator det = fvDets.begin();
       det != fvDets.end();
       ++det) {
    switch ((*det).detId) {
      case ECbmModuleId::kSts: {
        SearchMatches<CbmStsDigi>(dSeedTime, *det);
        break;
      }  // case ECbmModuleId::kSts:
      case ECbmModuleId::kMuch: {
        if (fbUseMuchBeamtimeDigi) {
          SearchMatches<CbmMuchBeamTimeDigi>(dSeedTime, *det);
        }  // if (fbUseMuchBeamtimeDigi)
        else {
          SearchMatches<CbmMuchDigi>(dSeedTime, *det);
        }  // else of if (fbUseMuchBeamtimeDigi)
        break;
      }  // case ECbmModuleId::kMuch:
      case ECbmModuleId::kTrd: {
        SearchMatches<CbmTrdDigi>(dSeedTime, *det);
        break;
      }  // case ECbmModuleId::kTrd:
      case ECbmModuleId::kTof: {
        SearchMatches<CbmTofDigi>(dSeedTime, *det);
        break;
      }  // case ECbmModuleId::kTof:
      case ECbmModuleId::kRich: {
        SearchMatches<CbmRichDigi>(dSeedTime, *det);
        break;
      }  // case ECbmModuleId::kRich:
      case ECbmModuleId::kPsd: {
        SearchMatches<CbmPsdDigi>(dSeedTime, *det);
        break;
      }  // case ECbmModuleId::kPsd:
      case ECbmModuleId::kT0: {
        SearchMatches<CbmTofDigi>(dSeedTime, *det);
        break;
      }  // case ECbmModuleId::kT0:
      default: {
        LOG(fatal) << "Cbm2021EventBuilderAlgo::LoopOnSeeds => "
                   << "Trying to search matches with unsupported det: "
                   << (*det).sName << std::endl
                   << "You may want to add support for it in the method.";
        break;
      }  // default:
    }    // switch( *det )
  }  // for( std::vector< EventBuilderDetector >::iterator det = fvDets.begin(); det != fvDets.end(); ++det )

  /// Check if event is filling trigger conditions and clear it if not
  if (HasTrigger(fCurrentEvent)) {
    fdPrevEvtTime = dSeedTime;

    /// In case of NoOverlap or MergeOverlap, we can and should start checking the next window
    /// from end of current window in order to save CPU and avoid duplicating
    if (EOverlapMode::NoOverlap == fOverMode
        || EOverlapMode::MergeOverlap == fOverMode) {

      /// Update reference detector
      fRefDet.fuStartIndex = fRefDet.fuEndIndex;

      /// Loop on selection detectors
      for (std::vector<EventBuilderDetector>::iterator det = fvDets.begin();
           det != fvDets.end();
           ++det) {
        (*det).fuStartIndex = (*det).fuEndIndex;
      }  // for( std::vector< EventBuilderDetector >::iterator det = fvDets.begin(); det != fvDets.end(); ++det )
    }    // If no overlap or merge overlap
  }      // if( !HasTrigger( fCurrentEvent ) )
  else {
    LOG(debug1) << "Reject seed due to Trigger requirements";
    delete fCurrentEvent;
    fCurrentEvent = nullptr;  /// delete does NOT set a pointer to nullptr...
  }                           // else of if( !HasTrigger( fCurrentEvent ) )
}

template<class DigiCheck>
void Cbm2021EventBuilderAlgo::SearchMatches(
  Double_t dSeedTime,
  EventBuilderDetector& detMatch) {
  /// This algo relies on time sorted vectors for the selected detectors
  UInt_t uLocalIndexStart = detMatch.fuStartIndex;
  UInt_t uLocalIndexEnd   = detMatch.fuStartIndex;

  /// Check the Digis until out of window
  if (ECbmModuleId::kT0 == detMatch.detId) {
    if (fT0DigiVec) {
      /// Loop on size of vector
      UInt_t uNbSelDigis = fT0DigiVec->size();
      /// Loop on size of vector
      for (UInt_t uDigi = detMatch.fuStartIndex; uDigi < uNbSelDigis; ++uDigi) {
        Double_t dTime = fT0DigiVec->at(uDigi).GetTime();

        Double_t dTimeDiff = dTime - dSeedTime;

        /// Check if within time window, update start/stop indices if needed
        if (dTimeDiff < detMatch.fdTimeWinBeg) {
          ++uLocalIndexStart;
          continue;
        }  // if( dTimeDiff < detMatch.fdTimeWinBeg )
        else if (detMatch.fdTimeWinEnd < dTimeDiff) {
          /// Store as end the first digi out of window to avoid double counting in case of
          /// merged overlap event mode
          uLocalIndexEnd = uDigi;
          break;
        }  // else if( detMatch.fdTimeWinEnd < dTimeDiff ) of if( dTimeDiff < detMatch.fdTimeWinBeg )

        AddDigiToEvent(detMatch, uDigi);

        if (fdPrevEvtEndTime < dTime) fdPrevEvtEndTime = dTime;
      }  // for( UInt_t uDigi = 0; uDigi < uNbSelDigis; ++uDigi )

      /// catch the case where we reach the end of the vector before being out of the time window
      if (uLocalIndexEnd < uLocalIndexStart) uLocalIndexEnd = uNbSelDigis;
    }  // if ( fT0DigiVec )
    else
      LOG(fatal) << "Cbm2021EventBuilderAlgo::SearchMatches => "
                 << "T0 as selection detector but vector not found!";
  }  // if( ECbmModuleId::kT0 == detMatch.detId )
  else {
    UInt_t uNbSelDigis = (0 < GetNofDigis(detMatch.detId)
                            ? GetNofDigis(detMatch.detId)
                            : 0);
    /// Loop on size of vector
    for (UInt_t uDigi = detMatch.fuStartIndex; uDigi < uNbSelDigis; ++uDigi) {
      const DigiCheck* pDigi = GetDigi<DigiCheck>(uDigi);
      /// Check that _entry is not out of range
      if (nullptr != pDigi) {
        Double_t dTime     = pDigi->GetTime();
        Double_t dTimeDiff = dTime - dSeedTime;

        LOG(debug4) << detMatch.sName
                    << Form(" => Checking match %6u / %6u, dt %f",
                            uDigi,
                            uNbSelDigis,
                            dTimeDiff);

        /// Check if within time window, update start/stop indices if needed
        if (dTimeDiff < detMatch.fdTimeWinBeg) {
          ++uLocalIndexStart;
          continue;
        }  // if( dTimeDiff < detMatch.fdTimeWinBeg )
        else if (detMatch.fdTimeWinEnd < dTimeDiff) {
          /// Store as end the first digi out of window to avoid double counting in case of
          /// merged overlap event mode
          uLocalIndexEnd = uDigi;
          break;
        }  // else if( detMatch.fdTimeWinEnd < dTimeDiff ) of if( dTimeDiff < detMatch.fdTimeWinBeg )

        AddDigiToEvent(detMatch, uDigi);

        if (fdPrevEvtEndTime < dTime) fdPrevEvtEndTime = dTime;
      }  // if( nullptr != pDigi )
    }    // for( UInt_t uDigi = 0; uDigi < uNbSelDigis; ++uDigi )

    /// catch the case where we reach the end of the vector before being out of the time window
    if (uLocalIndexEnd < uLocalIndexStart) uLocalIndexEnd = uNbSelDigis;
  }  // else of if( ECbmModuleId::kT0 == detMatch.detId ) => Digi containers controlled by DigiManager

  /// Update the StartIndex and EndIndex for the next event seed
  detMatch.fuStartIndex = uLocalIndexStart;
  detMatch.fuEndIndex   = uLocalIndexEnd;
}

void Cbm2021EventBuilderAlgo::AddDigiToEvent(
  EventBuilderDetector& det,
  Int_t _entry) {
  fCurrentEvent->AddData(det.dataType, _entry);
}

Bool_t Cbm2021EventBuilderAlgo::HasTrigger(CbmEvent* event) {
  /// Check first reference detector
  if (kFALSE == CheckTriggerConditions(event, fRefDet)) {
    return kFALSE;
  }  // if (kFALSE == CheckTriggerConditions(event, fRefDet) )

  /// Loop on selection detectors
  for (std::vector<EventBuilderDetector>::iterator det = fvDets.begin();
       det != fvDets.end();
       ++det) {
    if (kFALSE == CheckTriggerConditions(event, *det)) return kFALSE;
  }  // for( std::vector< EventBuilderDetector >::iterator det = fvDets.begin(); det != fvDets.end(); ++det )

  /// All Ok, trigger is there
  return kTRUE;
}

Bool_t Cbm2021EventBuilderAlgo::CheckTriggerConditions(
  CbmEvent* event,
  EventBuilderDetector& det) {
  /// Check if both Trigger conditions disabled for this detector
  if (0 == det.fuTriggerMinDigis && det.fiTriggerMaxDigis < 0) {
    return kTRUE;
  }  // if( 0 == det.fuTriggerMinDigis && det.fiTriggerMaxDigis < 0 )

  /// Check if detector present
  if (ECbmModuleId::kT0 == det.detId) {
    /// FIXME: special case to be removed once T0 supported by DigiManager
    if (!(fT0DigiVec)) {
      LOG(warning) << "Event does not have digis storage for T0"
                   << " while the following trigger minimum are defined: "
                   << det.fuTriggerMinDigis << " " << det.fiTriggerMaxDigis;
      return kFALSE;
    }  // if( !(fT0DigiVec) )
  }    // if( ECbmDataType::kT0Digi == det.detId )
  else {
    if (!DetIsPresent(det.detId)) {
      LOG(warning) << "Event does not have digis storage for " << det.sName
                   << " while the following trigger min/max are defined: "
                   << det.fuTriggerMinDigis << " " << det.fiTriggerMaxDigis;
      return kFALSE;
    }  // if( !fDigiMan->IsPresent( det ) )
  }    // else of if( ECbmDataType::kT0Digi == det )

  /// Check trigger rejection by minimal number or absence
  Int_t iNbDigis = event->GetNofData(det.dataType);
  if ((-1 == iNbDigis)
      || (static_cast<UInt_t>(iNbDigis) < det.fuTriggerMinDigis)) {
    LOG(debug2) << "Event does not have enough digis: " << iNbDigis << " vs "
                << det.fuTriggerMinDigis << " for " << det.sName;
    return kFALSE;
  }  // if((-1 == iNbDigis) || (static_cast<UInt_t>(iNbDigis) < det.fuTriggerMinDigis))
  /// Check trigger rejection by maximal number
  else if (0 < det.fiTriggerMaxDigis && det.fiTriggerMaxDigis < iNbDigis) {
    LOG(debug2) << "Event Has too many digis: " << iNbDigis << " vs "
                << det.fiTriggerMaxDigis << " for " << det.sName;
    return kFALSE;
  }  // else if( iNbDigis < det.fiTriggerMaxDigis )
  else {
    return kTRUE;
  }  // else of else if( iNbDigis < det.fiTriggerMaxDigis )
}
//----------------------------------------------------------------------
void Cbm2021EventBuilderAlgo::CreateHistograms() {
  /// FIXME: Disable clang formatting for histograms declaration for now
  /* clang-format off */
  fhEventTime = new TH1F("hEventTime",
                         "seed time of the events; Seed time [s]; Events",
                         60000, 0, 600);
  fhEventDt   = new TH1F( "fhEventDt",
                          "interval in seed time of consecutive events; Seed time [s]; Events",
                          2100, -100.5, 1999.5);
  fhEventSize =
    new TH1F("hEventSize",
             "nb of all  digis in the event; Nb Digis []; Events []",
             10000, 0, 10000);
  fhNbDigiPerEvtTime =
    new TH2I("hNbDigiPerEvtTime",
             "nb of all  digis per event vs seed time of the events; Seed time "
             "[s]; Nb Digis []; Events []",
              600, 0,   600,
             1000, 0, 10000);

  /// Loop on selection detectors
  for (std::vector< EventBuilderDetector >::iterator det = fvDets.begin(); det != fvDets.end(); ++det) {
    /// In case name not provided, do not create the histo to avoid name conflicts!
    if( "Invalid" == (*det).sName )
    {
      fvhNbDigiPerEvtTimeDet.push_back( nullptr );
      continue;
    } // if( "Invalid" == (*det).sName )

    fvhNbDigiPerEvtTimeDet.push_back(
      new TH2I( Form( "hNbDigiPerEvtTime%s", (*det).sName.data() ),
                Form( "nb of %s digis per event vs seed time of the events; Seed time "
                      "[s]; Nb Digis []; Events []",
                      (*det).sName.data() ),
                 600, 0,  600,
                4000, 0, 4000) );
  }  // for( std::vector< EventBuilderDetector >::iterator det = fvDets.begin(); det != fvDets.end(); ++det )

  AddHistoToVector(fhEventTime,            "evtbuild");
  AddHistoToVector(fhEventDt,              "evtbuild");
  AddHistoToVector(fhEventSize,            "evtbuild");
  AddHistoToVector(fhNbDigiPerEvtTime,     "evtbuild");
  for (std::vector<TH2*>::iterator itHist = fvhNbDigiPerEvtTimeDet.begin();
       itHist != fvhNbDigiPerEvtTimeDet.end();
       ++itHist) {
    if( nullptr != (*itHist) )
    {
      AddHistoToVector((*itHist),   "evtbuild");
    } // if( nullptr != (*itHist) )
  }  // for( std::vector<TH2*>::iterator itHist = fvhNbDigiPerEvtTimeDet.begin(); itHist != fvhNbDigiPerEvtTimeDet.end(); ++itHist )

  /// FIXME: Re-enable clang formatting after histograms declaration
  /* clang-format on */
}
void Cbm2021EventBuilderAlgo::FillHistos() {
  Double_t dPreEvtTime = -1.0;
  for (CbmEvent* evt : fEventVector) {
    fhEventTime->Fill(evt->GetStartTime() * 1e-9);
    if (0.0 <= dPreEvtTime) {
      fhEventDt->Fill(evt->GetStartTime() - dPreEvtTime);
    }  // if( 0.0 <= dPreEvtTime )
    fhEventSize->Fill(evt->GetNofData());
    fhNbDigiPerEvtTime->Fill(evt->GetStartTime() * 1e-9, evt->GetNofData());

    /// Loop on selection detectors
    for (UInt_t uDetIdx = 0; uDetIdx < fvDets.size(); ++uDetIdx) {
      if (nullptr == fvhNbDigiPerEvtTimeDet[uDetIdx]) continue;

      fvhNbDigiPerEvtTimeDet[uDetIdx]->Fill(
        evt->GetStartTime() * 1e-9, evt->GetNofData(fvDets[uDetIdx].dataType));
    }  // for( std::vector< EventBuilderDetector >::iterator det = fvDets.begin(); det != fvDets.end(); ++det )

    dPreEvtTime = evt->GetStartTime();
  }  // for( CbmEvent * evt: fEventVector )
}
void Cbm2021EventBuilderAlgo::ResetHistograms(
  Bool_t /*bResetTime*/) {
  fhEventTime->Reset();
  fhEventDt->Reset();
  fhEventSize->Reset();

  fhNbDigiPerEvtTime->Reset();
  /// Loop on histograms
  for (std::vector<TH2*>::iterator itHist = fvhNbDigiPerEvtTimeDet.begin();
       itHist != fvhNbDigiPerEvtTimeDet.end();
       ++itHist) {
    (*itHist)->Reset();
  }  // for( std::vector<TH2*>::iterator itHist = fvhNbDigiPerEvtTimeDet.begin(); itHist != fvhNbDigiPerEvtTimeDet.end(); ++itHist )

  /*
   if( kTRUE == bResetTime )
   {
      /// Also reset the Start time for the evolution plots!
      fdStartTime = -1.0;
   } // if( kTRUE == bResetTime )
*/
}
//----------------------------------------------------------------------
void Cbm2021EventBuilderAlgo::SetReferenceDetector(
  ECbmModuleId refDet,
  ECbmDataType dataTypeIn,
  std::string sNameIn,
  UInt_t uTriggerMinDigisIn,
  Int_t iTriggerMaxDigisIn,
  Double_t fdTimeWinBegIn,
  Double_t fdTimeWinEndIn) {

  /// FIXME: Deprecated method to be removed later. For now create temp object.
  SetReferenceDetector(EventBuilderDetector(refDet,
                                            dataTypeIn,
                                            sNameIn,
                                            uTriggerMinDigisIn,
                                            iTriggerMaxDigisIn,
                                            fdTimeWinBegIn,
                                            fdTimeWinEndIn));
}
void Cbm2021EventBuilderAlgo::AddDetector(ECbmModuleId selDet,
                                                     ECbmDataType dataTypeIn,
                                                     std::string sNameIn,
                                                     UInt_t uTriggerMinDigisIn,
                                                     Int_t iTriggerMaxDigisIn,
                                                     Double_t fdTimeWinBegIn,
                                                     Double_t fdTimeWinEndIn) {

  /// FIXME: Deprecated method to be removed later. For now create temp object.
  AddDetector(EventBuilderDetector(selDet,
                                   dataTypeIn,
                                   sNameIn,
                                   uTriggerMinDigisIn,
                                   iTriggerMaxDigisIn,
                                   fdTimeWinBegIn,
                                   fdTimeWinEndIn));
}
//----------------------------------------------------------------------
//----------------------------------------------------------------------
void Cbm2021EventBuilderAlgo::SetReferenceDetector(
  EventBuilderDetector refDetIn) {
  /// Loop on selection detectors
  for (std::vector<EventBuilderDetector>::iterator det = fvDets.begin();
       det != fvDets.end();
       ++det) {
    if ((*det) == refDetIn) {
      LOG(warning)
        << "Cbm2021EventBuilderAlgo::SetReferenceDetector => "
           "Reference detector already in selection detector list!"
        << refDetIn.sName;
      LOG(warning)
        << "                                                         => "
           "It will be automatically removed from selection detector list!";
      LOG(warning)
        << "                                                         => "
           "Please also remember to update the selection windows to store "
           "clusters!";
      RemoveDetector(refDetIn);
    }  // if( (*det)  == refDetIn )
  }  // for( std::vector< EventBuilderDetector >::iterator det = fvDets.begin(); det != fvDets.end(); ++det )

  if (fRefDet == refDetIn) {
    LOG(warning)
      << "Cbm2021EventBuilderAlgo::SetReferenceDetector => "
         "Doing nothing, identical reference detector already in use";
  }  // if( fRefDet == refDetIn )
  else {
    LOG(info) << "Cbm2021EventBuilderAlgo::SetReferenceDetector => "
              << "Replacing " << fRefDet.sName << " with " << refDetIn.sName
              << " as reference detector";
    LOG(warning)
      << "                                                         => "
         "You may want to use AddDetector after this command to add in "
         "selection "
      << refDetIn.sName;
    LOG(warning)
      << "                                                         => "
         "Please also remember to update the selection windows!";
  }  // else of if( fRefDet == refDetIn )
  fRefDet = refDetIn;

  /// Update the variables storing the earliest and latest time window boundaries
  UpdateTimeWinBoundariesExtrema();
  /// Update the variable storing the size if widest time window for overlap detection
  UpdateWidestTimeWinRange();
}
void Cbm2021EventBuilderAlgo::AddDetector(
  EventBuilderDetector selDet) {
  if (fRefDet == selDet) {
    LOG(fatal) << "Cbm2021EventBuilderAlgo::AddDetector => Cannot "
                  "add the reference detector as selection detector!"
               << std::endl
               << "=> Maybe first change the reference detector with "
                  "SetReferenceDetector?";
  }  // if( fRefDet == selDet )

  /// Loop on selection detectors
  for (std::vector<EventBuilderDetector>::iterator det = fvDets.begin();
       det != fvDets.end();
       ++det) {
    if ((*det) == selDet) {
      LOG(warning) << "Cbm2021EventBuilderAlgo::AddDetector => "
                      "Doing nothing, selection detector already in list!"
                   << selDet.sName;
      return;
    }  // if( (*det)  == selDet )
  }  // for( std::vector< EventBuilderDetector >::iterator det = fvDets.begin(); det != fvDets.end(); ++det )
  fvDets.push_back(selDet);

  /// Update the variables storing the earliest and latest time window boundaries
  UpdateTimeWinBoundariesExtrema();
  /// Update the variable storing the size if widest time window for overlap detection
  UpdateWidestTimeWinRange();
}
void Cbm2021EventBuilderAlgo::RemoveDetector(
  EventBuilderDetector selDet) {
  /// Loop on selection detectors
  for (std::vector<EventBuilderDetector>::iterator det = fvDets.begin();
       det != fvDets.end();
       ++det) {
    if ((*det) == selDet) {
      fvDets.erase(det);
      return;
    }  // if( (*det)  == selDet )
  }  // for( std::vector< EventBuilderDetector >::iterator det = fvDets.begin(); det != fvDets.end(); ++det )
  LOG(warning) << "Cbm2021EventBuilderAlgo::RemoveDetector => Doing "
                  "nothing, selection detector not in list!"
               << selDet.sName;
}
//----------------------------------------------------------------------
void Cbm2021EventBuilderAlgo::SetTriggerMinNumber(
  ECbmModuleId selDet,
  UInt_t uVal) {
  /// Check first if reference detector
  if (fRefDet.detId == selDet) {
    fRefDet.fuTriggerMinDigis = uVal;

    LOG(debug) << "Set Trigger min limit for " << fRefDet.sName << " to "
               << uVal;

    return;
  }  // if( fRefDet == selDet )

  /// Loop on selection detectors
  for (std::vector<EventBuilderDetector>::iterator det = fvDets.begin();
       det != fvDets.end();
       ++det) {
    if ((*det).detId == selDet) {
      (*det).fuTriggerMinDigis = uVal;

      LOG(debug) << "Set Trigger min limit for " << (*det).sName << " to "
                 << uVal;

      return;
    }  // if( (*det).detId  == selDet )
  }  // for( std::vector< EventBuilderDetector >::iterator det = fvDets.begin(); det != fvDets.end(); ++det )

  LOG(warning)
    << "Cbm2021EventBuilderAlgo::SetTriggerMinNumber => "
       "Doing nothing, detector neither reference nor in selection list!"
    << selDet;
}
void Cbm2021EventBuilderAlgo::SetTriggerMaxNumber(
  ECbmModuleId selDet,
  Int_t iVal) {
  /// Check first if reference detector
  if (fRefDet.detId == selDet) {
    fRefDet.fiTriggerMaxDigis = iVal;

    LOG(debug) << "Set Trigger min limit for " << fRefDet.sName << " to "
               << iVal;

    return;
  }  // if( fRefDet == selDet )

  /// Loop on selection detectors
  for (std::vector<EventBuilderDetector>::iterator det = fvDets.begin();
       det != fvDets.end();
       ++det) {
    if ((*det).detId == selDet) {
      (*det).fiTriggerMaxDigis = iVal;

      LOG(debug) << "Set Trigger min limit for " << (*det).sName << " to "
                 << iVal;

      return;
    }  // if( (*det).detId  == selDet )
  }  // for( std::vector< EventBuilderDetector >::iterator det = fvDets.begin(); det != fvDets.end(); ++det )

  LOG(warning)
    << "Cbm2021EventBuilderAlgo::SetTriggerMaxNumber => "
       "Doing nothing, detector neither reference nor in selection list!"
    << selDet;
}
void Cbm2021EventBuilderAlgo::SetTriggerWindow(ECbmModuleId selDet,
                                                          Double_t dWinBeg,
                                                          Double_t dWinEnd) {
  /// Check if valid time window: end strictly after beginning
  if (dWinEnd <= dWinBeg)
    LOG(fatal) << "Cbm2021EventBuilderAlgo::SetTriggerWindow => "
                  "Invalid time window: [ "
               << dWinBeg << ", " << dWinEnd << " ]";

  Bool_t bFound = kFALSE;
  /// Check first if reference detector
  if (fRefDet.detId == selDet) {
    fRefDet.fdTimeWinBeg = dWinBeg;
    fRefDet.fdTimeWinEnd = dWinEnd;

    bFound = kTRUE;
  }  // if( fRefDet == selDet )

  /// Loop on selection detectors
  for (std::vector<EventBuilderDetector>::iterator det = fvDets.begin();
       det != fvDets.end();
       ++det) {
    if ((*det).detId == selDet) {
      (*det).fdTimeWinBeg = dWinBeg;
      (*det).fdTimeWinEnd = dWinEnd;

      bFound = kTRUE;
    }  // if( (*det).detId  == selDet )
  }  // for( std::vector< EventBuilderDetector >::iterator det = fvDets.begin(); det != fvDets.end(); ++det )

  if (kFALSE == bFound) {
    LOG(warning)
      << "Cbm2021EventBuilderAlgo::SetTriggerWindow => "
         "Doing nothing, detector neither reference nor in selection list!"
      << selDet;
  }  // if( kFALSE == bFound )

  /// Update the variables storing the earliest and latest time window boundaries
  UpdateTimeWinBoundariesExtrema();
  /// Update the variable storing the size if widest time window for overlap detection
  UpdateWidestTimeWinRange();
}
void Cbm2021EventBuilderAlgo::UpdateTimeWinBoundariesExtrema() {
  /// Initialize with reference detector
  fdEarliestTimeWinBeg = fRefDet.fdTimeWinBeg;
  fdLatestTimeWinEnd   = fRefDet.fdTimeWinEnd;

  /// Loop on selection detectors
  for (std::vector<EventBuilderDetector>::iterator det = fvDets.begin();
       det != fvDets.end();
       ++det) {
    fdEarliestTimeWinBeg = std::min(fdEarliestTimeWinBeg, (*det).fdTimeWinBeg);
    fdLatestTimeWinEnd   = std::max(fdLatestTimeWinEnd, (*det).fdTimeWinEnd);
  }  // for( std::vector< EventBuilderDetector >::iterator det = fvDets.begin(); det != fvDets.end(); ++det )
}
void Cbm2021EventBuilderAlgo::UpdateWidestTimeWinRange() {
  /// Initialize with reference detector
  fdWidestTimeWinRange = fRefDet.fdTimeWinEnd - fRefDet.fdTimeWinBeg;

  /// Loop on selection detectors
  for (std::vector<EventBuilderDetector>::iterator det = fvDets.begin();
       det != fvDets.end();
       ++det) {
    fdWidestTimeWinRange =
      std::max(fdWidestTimeWinRange, (*det).fdTimeWinEnd - (*det).fdTimeWinBeg);
  }  // for( std::vector< EventBuilderDetector >::iterator det = fvDets.begin(); det != fvDets.end(); ++det )
}
//----------------------------------------------------------------------

ClassImp(Cbm2021EventBuilderAlgo)