diff --git a/macro/reco/reco_ts.C b/macro/reco/reco_ts.C index e520b169d94ecab63c2005abbd702c95d689118a..ea0cebc743deb64122323809161f4842fb639518 100644 --- a/macro/reco/reco_ts.C +++ b/macro/reco/reco_ts.C @@ -65,6 +65,7 @@ void reco_ts(TString tsaFile = "", TString outFile = "") // TODO: Would need a small up-to-date default input file; the one distributed with // the code is outdated. + // ----- Default file names --------------------------------------------- if (tsaFile.IsNull()) tsaFile = srcDir + "/input/mcbm_run399_first20Ts"; TString inFile = tsaFile + ".tsa"; @@ -121,6 +122,13 @@ void reco_ts(TString tsaFile = "", TString outFile = "") // ------------------------------------------------------------------------ + // ----- Event QA ----------------------------------------------------- + auto eventQa = std::make_unique<CbmTaskDigiEventQa>(); + LOG(info) << myName << ": Added task " << eventQa->GetName(); + run->AddTask(eventQa.release()); + // ------------------------------------------------------------------------ + + // ----- Logger settings ---------------------------------------------- FairLogger::GetLogger()->SetLogScreenLevel(logLevel.Data()); FairLogger::GetLogger()->SetLogVerbosityLevel(logVerbosity.Data()); diff --git a/reco/tasks/CMakeLists.txt b/reco/tasks/CMakeLists.txt index dd87b3a51de4bf5eee4813a1e30b2dc18cd07808..7d98de2467e7d82c137bfb90a86f20ceb3bf3e22 100644 --- a/reco/tasks/CMakeLists.txt +++ b/reco/tasks/CMakeLists.txt @@ -13,6 +13,7 @@ set(SRCS CbmReco.cxx CbmSourceTs.cxx CbmTaskBuildEvents.cxx +CbmTaskDigiEventQa.cxx CbmTaskMakeRecoEvents.cxx CbmTaskTriggerDigi.cxx CbmTaskUnpack.cxx diff --git a/reco/tasks/CbmRecoTasksLinkDef.h b/reco/tasks/CbmRecoTasksLinkDef.h index 9d2dbacc5b1bc63f4f357949980a9a098d3508d1..ac0ef949083e0ac7e4329d9aca8f1a5c0ea99b88 100644 --- a/reco/tasks/CbmRecoTasksLinkDef.h +++ b/reco/tasks/CbmRecoTasksLinkDef.h @@ -15,6 +15,7 @@ #pragma link C++ class CbmRecoConfig + ; #pragma link C++ class CbmSourceTs + ; #pragma link C++ class CbmTaskBuildEvents + ; +#pragma link C++ class CbmTaskDigiEventQa + ; #pragma link C++ class CbmTaskMakeRecoEvents + ; #pragma link C++ class CbmTaskTriggerDigi + ; #pragma link C++ class CbmTaskUnpack + ; diff --git a/reco/tasks/CbmTaskDigiEventQa.cxx b/reco/tasks/CbmTaskDigiEventQa.cxx new file mode 100644 index 0000000000000000000000000000000000000000..c775574f50cf303637df45dfbe138bb06e387604 --- /dev/null +++ b/reco/tasks/CbmTaskDigiEventQa.cxx @@ -0,0 +1,144 @@ +/* Copyright (C) 2022 GSI Helmholtzzentrum fuer Schwerionenforschung, Darmstadt + SPDX-License-Identifier: GPL-3.0-only + Authors: Volker Friese [committer] */ + +#include "CbmTaskDigiEventQa.h" + +#include "CbmDigiBranchBase.h" +#include "CbmDigiManager.h" +#include "CbmDigiTimeslice.h" +#include "CbmModuleList.h" + +#include <FairLogger.h> +#include <FairRunOnline.h> + +#include <TH1F.h> +#include <THttpServer.h> +#include <TStopwatch.h> + +#include <cassert> +#include <iomanip> +#include <iostream> + +using namespace std; + +// ----- Constructor ----------------------------------------------------- +CbmTaskDigiEventQa::CbmTaskDigiEventQa() : FairTask("DigiEventQa") {} +// --------------------------------------------------------------------------- + + +// ----- Destructor ------------------------------------------------------ +CbmTaskDigiEventQa::~CbmTaskDigiEventQa() {} +// --------------------------------------------------------------------------- + + +// ----- Execution ------------------------------------------------------- +void CbmTaskDigiEventQa::Exec(Option_t*) +{ + + // --- Timer and counters + TStopwatch timer; + timer.Start(); + double sumT = 0.; + double sumTsq = 0.; + double sumNumDigis = 0.; + double sumNumDigisSq = 0.; + size_t numEvents = 0; + size_t numDigis = 0; + + // --- Event loop + for (const auto& event : *fEvents) { + numEvents++; + double numDigisEvt = double(event.fData.fSts.fDigis.size()); + sumNumDigis += numDigisEvt; + sumNumDigisSq += numDigisEvt * numDigisEvt; + for (const auto& digi : event.fData.fSts.fDigis) { + numDigis++; + const double tDigi = digi.GetTime() - event.fTime; + fHistDigiTime->Fill(tDigi); + sumT += tDigi; + sumTsq += tDigi * tDigi; + } + } + + // --- First and second moments of digi times and digi numbers + double tMean = sumT / double(numDigis); + double tSqMean = sumTsq / double(numDigis); + double tRms = sqrt(tSqMean - tMean * tMean); + double numDigisMean = sumNumDigis / double(numEvents); + double numDigisSqMean = sumNumDigisSq / double(numEvents); + double numDigisRms = sqrt(numDigisSqMean - numDigisMean * numDigisMean); + + // --- Timeslice log + timer.Stop(); + fExecTime += timer.RealTime(); + stringstream logOut; + logOut << setw(15) << left << GetName() << " ["; + logOut << fixed << setw(8) << setprecision(1) << right << timer.RealTime() * 1000. << " ms] "; + logOut << "TS " << fNumTs << ", events " << numEvents << ", digis " << numDigis << ", digis per event (" + << numDigisMean << " +- " << numDigisRms << ")" + << ", digi time (" << tMean << " +- " << tRms << ") ns"; + LOG(info) << logOut.str(); + + // --- Run statistics + fNumTs++; + fNumEvents += numEvents; + fNumDigis += numDigis; + + + // See: macro/run/run_unpack_online.C + // See : recp/detectors/sts/unpack/CbmStsUnpackMonitor +} +// ---------------------------------------------------------------------------- + + +// ----- End-of-timeslice action ------------------------------------------ +void CbmTaskDigiEventQa::Finish() +{ + std::cout << std::endl; + LOG(info) << "====================================="; + LOG(info) << GetName() << ": Run summary"; + LOG(info) << "Timeslices : " << fNumTs; + LOG(info) << "Events : " << fNumEvents; + LOG(info) << "Digis : " << fNumDigis; + LOG(info) << "Exec time : " << fixed << setprecision(2) << 1000. * fExecTime / double(fNumTs) << " ms / TS"; + LOG(info) << "Digi times : entries " << fHistDigiTime->GetEntries() << ", mean " << fHistDigiTime->GetMean() + << ", stddev " << fHistDigiTime->GetStdDev(); + LOG(info) << "====================================="; + fHistDigiTime->Write(); +} +// ---------------------------------------------------------------------------- + + +// ----- Initialisation --------------------------------------------------- +InitStatus CbmTaskDigiEventQa::Init() +{ + // --- Get FairRootManager instance + FairRootManager* ioman = FairRootManager::Instance(); + assert(ioman); + + std::cout << std::endl; + 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"; + + // --- Histograms + fHistDigiTime = new TH1F("hDigiTime", "Digi time in event", 100, -50., 50.); + THttpServer* server = FairRunOnline::Instance()->GetHttpServer(); + if (server) server->Register("DigiTime", fHistDigiTime); + + LOG(info) << "=================================================="; + std::cout << std::endl; + + return kSUCCESS; +} +// ---------------------------------------------------------------------------- + +ClassImp(CbmTaskDigiEventQa) diff --git a/reco/tasks/CbmTaskDigiEventQa.h b/reco/tasks/CbmTaskDigiEventQa.h new file mode 100644 index 0000000000000000000000000000000000000000..6a5d60ebf8128e3a334625391a96f7c85cb748d8 --- /dev/null +++ b/reco/tasks/CbmTaskDigiEventQa.h @@ -0,0 +1,69 @@ +/* Copyright (C) 2022 GSI Helmholtzzentrum fuer Schwerionenforschung, Darmstadt + SPDX-License-Identifier: GPL-3.0-only + Authors: Volker Friese [committer] */ + +#ifndef CBMTASKDIGIEVENTQA_H +#define CBMTASKDIGIEVENTQA_H 1 + +#include "CbmDefs.h" +#include "CbmDigiEvent.h" + +#include <FairTask.h> + +#include <vector> + +#include "EventBuilder.h" + +class CbmDigiManager; +class TH1F; + +/** @class CbmTaskDigiEventQa + ** @brief QA task class for digi events produced by the event builder + ** @author Volker Friese <v.friese@gsi.de> + ** @since 15.03.2022 + **/ +class CbmTaskDigiEventQa : public FairTask { + +public: + /** @brief Constructor **/ + CbmTaskDigiEventQa(); + + + /** @brief Copy constructor (disabled) **/ + CbmTaskDigiEventQa(const CbmTaskDigiEventQa&) = delete; + + + /** @brief Destructor **/ + virtual ~CbmTaskDigiEventQa(); + + + /** @brief Task execution **/ + virtual void Exec(Option_t* opt); + + + /** @brief Finish timeslice **/ + virtual void Finish(); + + + /** @brief Assignment operator (disabled) **/ + CbmTaskDigiEventQa& operator=(const CbmTaskDigiEventQa&) = 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 + size_t fNumDigis = 0; ///< Number of analysed digis + double fExecTime = 0.; ///< Execution time [s] + TH1F* fHistDigiTime = nullptr; //! + + + ClassDef(CbmTaskDigiEventQa, 1); +}; + +#endif /* CBMTASKDIGIEVENTQA_H */