Skip to content
Snippets Groups Projects
Commit 685dcb76 authored by Volker Friese's avatar Volker Friese Committed by Jan de Cuveland
Browse files

QA for digi events. Checks the digi time distribution within the event.

parent e465a41f
No related branches found
No related tags found
1 merge request!784QA for digi events. Checks the digi time distribution within the event.
......@@ -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());
......
......@@ -13,6 +13,7 @@ set(SRCS
CbmReco.cxx
CbmSourceTs.cxx
CbmTaskBuildEvents.cxx
CbmTaskDigiEventQa.cxx
CbmTaskMakeRecoEvents.cxx
CbmTaskTriggerDigi.cxx
CbmTaskUnpack.cxx
......
......@@ -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 + ;
......
/* 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)
/* 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 */
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