Skip to content
Snippets Groups Projects
Commit b4b9bbd6 authored by Volker Friese's avatar Volker Friese Committed by Pierre-Alain Loizeau
Browse files

Implementation of CbmSourceDigiEvents, reading DigiEvents from archive into the ROOT tree

parent 5dfe6fd7
No related branches found
No related tags found
1 merge request!1717Implementation of CbmSourceDigiEvents, reading DigiEvents from archive into the ROOT tree
Pipeline #28032 passed
/* Copyright (C) 2024 GSI Helmholtzzentrum fuer Schwerionenforschung, Darmstadt
SPDX-License-Identifier: GPL-3.0-only
Authors: Volker Friese [committer] */
/** @file run_inspect_digievents.C
** @author Volker Friese <v.friese@gsi.de>
** @since 19 March 2024
**/
// --- Includes needed for IDE
#include <RtypesCore.h>
#if !defined(__CLING__)
#include "CbmSourceDigiEvents.h"
#include "CbmTaskInspectDigiEvents.h"
#include <FairRunAna.h>
#include <FairSystemInfo.h>
#include <TStopwatch.h>
#endif
/** @brief Macro for simple inspection of DigiEvents
** @author Volker Friese <v.friese@gsi.de>
** @since 219 March 2024
** @param input Name of input file
** @param nTimeSlices Number of time-slices to process
**
** This macro performs a simple inspection of DigiEvents produced an archived by online data processing.
** It prints statistics of events on the console. The purpose is to test CbmSourceDigiEvents and to
** demonstrate the event loop within a timeslice and the access to digi data in the event.
**/
void run_inspect_digievents(TString inputFileName, TString outputFileName, size_t numTimeslices = -1)
{
// ========================================================================
// Adjust this part according to your requirements
// --- Logger settings ----------------------------------------------------
TString logLevel = "INFO";
TString logVerbosity = "LOW";
// ------------------------------------------------------------------------
// ----- Environment --------------------------------------------------
TString myName = "run_inspect_digievents"; // this macro's name for screen output
TString srcDir = gSystem->Getenv("VMCWORKDIR"); // top source directory
// ------------------------------------------------------------------------
// ----- Timer --------------------------------------------------------
TStopwatch timer;
timer.Start();
// ------------------------------------------------------------------------
// ----- FairRunAna ---------------------------------------------------
FairRunAna* run = new FairRunAna();
FairSource* source = new CbmSourceDigiEvents(inputFileName);
run->SetSource(source);
run->SetOutputFile(outputFileName);
run->SetGenerateRunInfo(kTRUE);
// ------------------------------------------------------------------------
// ----- Logger settings ----------------------------------------------
FairLogger::GetLogger()->SetLogScreenLevel(logLevel.Data());
FairLogger::GetLogger()->SetLogVerbosityLevel(logVerbosity.Data());
// ------------------------------------------------------------------------
// ----- Event inspection ---------------------------------------------
//std::unique_ptr<CbmTaskInspectDigiEvents> inspect = std::make_unique<CbmTaskInspectDigiEvents>;
//LOG(info) << "-I- " << myName << ": Adding task " << inspect->GetName();
//run->AddTask(inspect.release());
FairTask* inspect = new CbmTaskInspectDigiEvents();
LOG(info) << "-I- " << myName << ": Adding task " << inspect->GetName();
run->AddTask(inspect);
// ------------------------------------------------------------------------
// ----- Run initialisation -------------------------------------------
std::cout << std::endl;
std::cout << "-I- " << myName << ": Initialise run" << std::endl;
run->Init();
// ------------------------------------------------------------------------
// ----- Start run ----------------------------------------------------
std::cout << std::endl << std::endl;
std::cout << "-I- " << myName << ": Starting run" << std::endl;
if (numTimeslices == -1)
run->Run();
else
run->Run(0, numTimeslices);
// ------------------------------------------------------------------------
// ----- Finish -------------------------------------------------------
timer.Stop();
FairMonitor::GetMonitor()->Print();
Double_t rtime = timer.RealTime();
Double_t ctime = timer.CpuTime();
std::cout << std::endl << std::endl;
std::cout << "Macro finished successfully." << std::endl;
std::cout << "Real time " << rtime << " s, CPU time " << ctime << " s" << std::endl;
FairSystemInfo sysInfo;
Float_t maxMemory = sysInfo.GetMaxMemory();
std::cout << "<DartMeasurement name=\"MaxMemory\" type=\"numeric/double\">";
std::cout << maxMemory;
std::cout << "</DartMeasurement>" << std::endl;
Float_t cpuUsage = ctime / rtime;
std::cout << "<DartMeasurement name=\"CpuLoad\" type=\"numeric/double\">";
std::cout << cpuUsage;
std::cout << "</DartMeasurement>" << std::endl;
// ------------------------------------------------------------------------
} // End of main macro function
......@@ -8,6 +8,7 @@ set(INCLUDE_DIRECTORIES
set(SRCS
CbmRecoUnpack.cxx
CbmSourceDigiEvents.cxx
CbmSourceTsArchive.cxx
)
......
......@@ -10,6 +10,7 @@
// --- Classes
#pragma link C++ class CbmRecoUnpack + ;
#pragma link C++ class CbmSourceDigiEvents + ;
#pragma link C++ class CbmSourceTsArchive + ;
#endif /* __CINT__ */
/* Copyright (C) 2024 GSI Helmholtzzentrum fuer Schwerionenforschung, Darmstadt
SPDX-License-Identifier: GPL-3.0-only
Authors: Volker Friese[committer] */
#include "CbmSourceDigiEvents.h"
#include <FairRootManager.h>
#include <Logger.h>
#include <utility>
// ----- Constructor ------------------------------------------------------
CbmSourceDigiEvents::CbmSourceDigiEvents(const char* fileName) : fInputFileName(fileName) {}
// ----------------------------------------------------------------------------
// ----- Close ------------------------------------------------------------
void CbmSourceDigiEvents::Close()
{
LOG(info) << "Source: Closing after " << fNumTs << "timeslices with " << fNumEvents << "events.";
}
// ----------------------------------------------------------------------------
// ----- Initialisation ---------------------------------------------------
Bool_t CbmSourceDigiEvents::Init()
{
// Create input archive
fArchive = std::make_unique<cbm::algo::RecoResultsInputArchive>(fInputFileName);
LOG(info) << "Source: Reading from input archive " << fInputFileName;
auto desc = fArchive->descriptor();
LOG(info) << " - Time created: " << desc.time_created();
LOG(info) << " - Host name : " << desc.hostname();
LOG(info) << " - User name : " << desc.username();
// Create and register the DigiEvent tree branch
FairRootManager* ioman = FairRootManager::Instance();
assert(ioman);
if (ioman->GetObject("DigiEvent")) {
LOG(fatal) << "Source: Branch DigiEvent already exists!";
return kFALSE;
}
fEvents = new std::vector<CbmDigiEvent>();
ioman->RegisterAny("DigiEvent", fEvents, true);
LOG(info) << "Source: Registered branch DigiEvent";
return kTRUE;
}
// ----------------------------------------------------------------------------
// ----- Read one time slice from archive ---------------------------------
Int_t CbmSourceDigiEvents::ReadEvent(UInt_t)
{
// Clear output event vector
fEvents->clear();
// Get next timeslice data from archive. Stop run if end of archive is reached.
auto results = fArchive->get();
if (fArchive->eos()) {
LOG(info) << "Source: End of input archive; terminating run";
return 1;
}
// Move event data from input archive to ROOT tree
if (results == nullptr) {
LOG(error) << "Source: Failed to read RecoResults from archive";
return 1;
}
size_t numEvents = results->DigiEvents().size();
LOG(info);
LOG(info) << "Source: Reading TS " << fNumTs << ", index " << results->TsIndex() << ", start "
<< results->TsStartTime() << ", events " << numEvents;
std::move(results->DigiEvents().begin(), results->DigiEvents().end(), std::back_inserter(*fEvents));
assert(fEvents->size() == numEvents);
// Update counters
fNumTs++;
fNumEvents += numEvents;
return 0;
}
// ----------------------------------------------------------------------------
ClassImp(CbmSourceDigiEvents)
/* Copyright (C) 2024 GSI Helmholtzzentrum fuer Schwerionenforschung, Darmstadt
SPDX-License-Identifier: GPL-3.0-only
Authors: Volker Friese[committer] */
#pragma once // include this header only once per compilation unit
#include "CbmDigiEvent.h"
#include "DigiData.h"
#include "RecoResultsInputArchive.h"
#include <FairSource.h>
#include <cstdint>
#include <memory>
#include <string>
#include <vector>
/** @class CbmSourceDigiEvents
** @brief Source class for reading from files resulting from online processing (containing DigiEvents)
** @author Volker Friese
** @since 18 March 2024
**
** The online process creates a std::vector of DigiEvents per timeslice. These are saved to file using the BOOST streamer.
** This class allows to read such files and get the DigiEvents into the ROOT tree for offline analysis. It creates a
** branch DigiEvents containing the DigiEvent vector; one tree entry corresponds to one timelice.
**/
class CbmSourceDigiEvents : public FairSource {
public:
/** @brief Constructor
** @param fileName Name of (single) input file.
**
** More input files can be added by the method AddInputFile.
*/
CbmSourceDigiEvents(const char* fileName = "");
/** @brief Copy constructor - not implemented **/
CbmSourceDigiEvents(const CbmSourceDigiEvents&) = delete;
/** @brief Assignment operator - not implemented **/
CbmSourceDigiEvents& operator=(const CbmSourceDigiEvents&) = delete;
/** @brief Destructor **/
virtual ~CbmSourceDigiEvents() {}
/** @brief Close source after end of run **/
virtual void Close();
/** @brief Source type
** @return FairSource::Source_Type
**/
virtual Source_Type GetSourceType() { return fSourceType; }
/** @brief Initialisation **/
virtual Bool_t Init();
/** @brief Initialise unpackers (forced by base class, not relevant here) **/
virtual Bool_t InitUnpackers() { return kTRUE; }
/** @brief Read one time slice from file **/
Int_t ReadEvent(UInt_t = 0);
/** @brief Re-Initialise unpackers (forced by base class, not relevant here) **/
virtual Bool_t ReInitUnpackers() { return kTRUE; }
/** @brief Reset (forced by base class, not relevant here) **/
virtual void Reset() {}
/** @brief Set unpacker parameters (forced by base class, not relevant here) **/
virtual void SetParUnpackers() {}
/** @brief Set the Source Type
** @param type Source type
**/
void SetSourceType(Source_Type type) { fSourceType = type; }
/** @brief Set Run ID (forced by base class, not relevant here) **/
Bool_t SpecifyRunId() { return kTRUE; }
private:
/** Input file name **/
std::string fInputFileName = {};
/** Input archive **/
std::unique_ptr<cbm::algo::RecoResultsInputArchive> fArchive = nullptr;
/** Branch vector of DigiEvents **/
std::vector<CbmDigiEvent>* fEvents = nullptr;
/** @brief type of source that is currently used @remark currently we use kONLINE as default, since, kFILE skipps the first TS probably due to obsolete reasons (to be checked PR072021) */
Source_Type fSourceType = Source_Type::kONLINE;
/** Time-slice counter **/
size_t fNumTs = 0;
/** Event counter **/
size_t fNumEvents = 0;
ClassDef(CbmSourceDigiEvents, 1)
};
......@@ -12,6 +12,7 @@ set(SRCS
CbmTaskBuildEvents.cxx
CbmTaskDigiEventQa.cxx
CbmTaskEventsCloneInToOut.cxx
CbmTaskInspectDigiEvents.cxx
CbmTaskMakeRecoEvents.cxx
CbmTaskTriggerDigi.cxx
CbmTaskTofHitFinder.cxx
......
......@@ -16,6 +16,7 @@
#pragma link C++ class CbmTaskBuildEvents + ;
#pragma link C++ class CbmTaskDigiEventQa + ;
#pragma link C++ class CbmTaskEventsCloneInToOut + ;
#pragma link C++ class CbmTaskInspectDigiEvents + ;
#pragma link C++ class CbmTaskMakeRecoEvents + ;
#pragma link C++ class CbmTaskTofHitFinder + ;
#pragma link C++ class CbmTaskTofClusterizer + ;
......
/* Copyright (C) 2022 GSI Helmholtzzentrum fuer Schwerionenforschung, Darmstadt
SPDX-License-Identifier: GPL-3.0-only
Authors: Volker Friese [committer] */
#include "CbmTaskInspectDigiEvents.h"
#include <FairRootManager.h>
#include <Logger.h>
// ----- Constructor -----------------------------------------------------
CbmTaskInspectDigiEvents::CbmTaskInspectDigiEvents() : FairTask("InspectDigiEvents") {}
// ---------------------------------------------------------------------------
// ----- Destructor ------------------------------------------------------
CbmTaskInspectDigiEvents::~CbmTaskInspectDigiEvents() {}
// ---------------------------------------------------------------------------
// ----- Execution -------------------------------------------------------
void CbmTaskInspectDigiEvents::Exec(Option_t*)
{
// --- Inspect event vector
LOG(info) << GetName() << ": timeslice " << fNumTs << " with " << fEvents->size() << " events"
<< (fEvents->size() > 10 ? " (showing the first 10 only)" : "");
size_t numEventsInTs = 0;
for (auto& event : *fEvents) {
size_t numBmon = event.fData.fBmon.Size();
size_t numSts = event.fData.fSts.Size();
size_t numMuch = event.fData.fMuch.Size();
size_t numRich = event.fData.fRich.Size();
size_t numTrd1d = event.fData.fTrd.Size();
size_t numTrd2d = event.fData.fTrd2d.Size();
size_t numTof = event.fData.fTof.Size();
size_t numFsd = event.fData.fFsd.Size();
LOG(info) << " Event " << event.fNumber << ", time " << event.fTime << ", digis: bmon " << numBmon << " sts "
<< numSts << " much " << numMuch << " rich " << numRich << " trd1d " << numTrd1d << " trd2d "
<< numTrd2d << " tof " << numTof << " fsd " << numFsd;
numEventsInTs++;
if (numEventsInTs > 9) break;
}
// --- Run statistics
fNumTs++;
fNumEvents += fEvents->size();
}
// ----------------------------------------------------------------------------
// ----- End-of-run action ------------------------------------------------
void CbmTaskInspectDigiEvents::Finish()
{
LOG(info) << "=====================================";
LOG(info) << GetName() << ": Run summary";
LOG(info) << "Timeslices : " << fNumTs;
LOG(info) << "Events : " << fNumEvents;
LOG(info) << "=====================================";
}
// ----------------------------------------------------------------------------
// ----- Initialisation ---------------------------------------------------
InitStatus CbmTaskInspectDigiEvents::Init()
{
// --- Get FairRootManager instance
FairRootManager* ioman = FairRootManager::Instance();
assert(ioman);
LOG(info) << "==================================================";
LOG(info) << GetName() << ": Initialising...";
// --- Input data
fEvents = ioman->InitObjectAs<const std::vector<CbmDigiEvent>*>("DigiEvent");
if (!fEvents) {
LOG(error) << GetName() << ": No input branch DigiEvent!";
return kFATAL;
}
LOG(info) << "--- Found branch DigiEvent";
LOG(info) << "==================================================";
return kSUCCESS;
}
// ----------------------------------------------------------------------------
ClassImp(CbmTaskInspectDigiEvents)
/* Copyright (C) 2024 GSI Helmholtzzentrum fuer Schwerionenforschung, Darmstadt
SPDX-License-Identifier: GPL-3.0-only
Authors: Volker Friese [committer] */
#pragma once
#include "CbmDefs.h"
#include "CbmDigiEvent.h"
#include <FairTask.h>
#include <vector>
/** @class CbmTaskInspectDigiEvents
** @brief Demonstrator class to look at digi events in the ROOT tree
** @author Volker Friese <v.friese@gsi.de>
** @since 19.03.2024
**
** This is just a demonstrator how to look at DigiEvents in the ROOT tree. It logs some information to the console.
**/
class CbmTaskInspectDigiEvents : public FairTask {
public:
/** @brief Constructor **/
CbmTaskInspectDigiEvents();
/** @brief Copy constructor (disabled) **/
CbmTaskInspectDigiEvents(const CbmTaskInspectDigiEvents&) = delete;
/** @brief Destructor **/
virtual ~CbmTaskInspectDigiEvents();
/** @brief Task execution **/
virtual void Exec(Option_t* opt);
/** @brief Finish timeslice **/
virtual void Finish();
/** @brief Assignment operator (disabled) **/
CbmTaskInspectDigiEvents& operator=(const CbmTaskInspectDigiEvents&) = delete;
private: // methods
/** @brief Task initialisation **/
virtual InitStatus Init();
private: // members
const std::vector<CbmDigiEvent>* fEvents = nullptr; //! Input data (events)
size_t fNumTs = 0; ///< Number of processed timeslices
size_t fNumEvents = 0; ///< Number of analysed events
ClassDef(CbmTaskInspectDigiEvents, 1);
};
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Finish editing this message first!
Please register or to comment