10#include "BmnDetectorList.h"
11#include "BmnFHCalGeo.h"
12#include "BmnFHCalPoint.h"
13#include "TTreeFormula.h"
18 fDoMatchMC = (fIsExp) ?
false :
true;
19 fIsGlobal = isGlobalCoordinates;
20 LOG(detail) << Form(
"BmnFHCalReconstructor: isExp %i; isGlobal %i", fIsExp, fIsGlobal);
25 fRecoCutsFile = reco_cuts_file;
27 LOG(detail) << Form(
"BmnFHCalReconstructor: Reco cuts %s", fRecoCutsFile.Data());
33 delete fBmnFHCalEvent;
39 fpFairRootMgr = FairRootManager::Instance();
42 fArrayOfDigits = (TClonesArray*)fpFairRootMgr->GetObject(
"FHCalDigi");
43 if (!fArrayOfDigits) {
44 LOG(error) <<
"BmnFHCalReconstructor::Init(): branch with Digits not found! Task will be deactivated";
49 fArrayOfDigits = (TClonesArray*)fpFairRootMgr->GetObject(
"FHCalDigit");
50 if (!fArrayOfDigits) {
51 LOG(error) <<
"BmnFHCalReconstructor::Init(): branch with Digits not found! Task will be deactivated";
56 fArrayOfPoints = (TClonesArray*)fpFairRootMgr->GetObject(
"FHCalPoint");
57 if (!fArrayOfPoints) {
58 LOG(error) <<
"BmnFHCalReconstructor::Init(): branch with FHCalPoints not found! No MC match will be done";
64 fBmnFHCalEvent->
reset();
66 fpFairRootMgr->RegisterAny(
"FHCalEvent", fBmnFHCalEvent, kTRUE);
67 fpFairRootMgr->Register(
"FHCalEvent",
"FHCal", fBmnFHCalEvent, kFALSE);
69 LOG(debug) <<
"FHCal Reconstructor ready";
75 namespace po = boost::program_options;
77 TString dir = getenv(
"VMCWORKDIR");
78 TString path = dir +
"/parameters/fhcal/";
86 po::options_description desc(
"Options");
87 desc.add_options()(
"VERSION.id", po::value<float>(&version),
88 "version identificator")(
"COMMENT.str", po::value<std::string>(&comment),
"comment")(
89 "EXPERIMENT.cut", po::value<std::string>(&cutExp),
90 "cut for Experimental data")(
"SIMULATION.cut", po::value<std::string>(&cutSim),
"cut for Simulated data");
94 std::ifstream config_file((path + fRecoCutsFile).Data(), std::ifstream::in);
95 if (!config_file.is_open()) {
96 LOG(error) << Form(
"BmnFHCalReconstructor: Loading Cuts from file: %s - file open error!",
97 fRecoCutsFile.Data());
100 LOG(debug) << Form(
"BmnFHCalReconstructor: Loading Cuts from file: %s", fRecoCutsFile.Data());
101 po::store(po::parse_config_file(config_file, desc), vm);
106 fSelectionString = cutExp;
108 fSelectionString = cutSim;
116 LOG(error) << Form(
"BmnFHCalReconstructor::Check geometry file. Trying to set number of FHCal modules %d",
120 for (
auto it = positionMap.begin(); it != positionMap.end(); ++it) {
121 auto unique_address = it->first;
122 auto position = it->second.first;
126 auto this_mod = fBmnFHCalEvent->
GetModule(mod_id);
128 this_mod->SetPosition(position);
130 this_mod->SetNsections(7);
132 this_mod->SetNsections(10);
134 auto err = it->second.second;
135 this_mod->SetDxyz(err.X(), err.Y(), err.Z() * this_mod->GetNsections());
148 auto digiVector = GetSelectedDigiVector(fSelectionString);
150 for (
auto it : digiVector) {
154 if (mod_id <= 0 || (
int)mod_id > fBmnFHCalEvent->
GetTotalModules() || sec_id <= 0
158 LOG(warning) << Form(
"FHCal digi ignored. Mod %d Sec %d", mod_id, sec_id);
164 if (!fIsExp && fDoMatchMC)
169 fworkTime += sw.RealTime();
172std::vector<TObject*> BmnFHCalReconstructor::GetSelectedDigiVector(TString formulaString)
174 TTree* temp_tree =
new TTree(
"temp_tree",
"temp_tree");
175 TString temp_br_name = (fIsExp) ?
"FHCalDigi" :
"FHCalDigit";
176 fSelectedDigiVector.clear();
180 temp_tree->Branch(temp_br_name.Data(), &ptr);
181 for (
int i = 0;
i < fArrayOfDigits->GetEntriesFast();
i++) {
187 temp_tree->Branch(temp_br_name.Data(), &ptr);
188 for (
int i = 0;
i < fArrayOfDigits->GetEntriesFast();
i++) {
194 TTreeFormula* formula =
new TTreeFormula(
"TTreeformula", formulaString.Data(), temp_tree);
195 if (!formula || formula->GetNdim() == 0)
196 LOG(error) << Form(
"BmnFHCalReconstructor::GetSelectedDigiVector. Bad cuts. Check formula %s from file %s",
197 formulaString.Data(), fRecoCutsFile.Data());
199 TObject* ThisDigi =
nullptr;
200 temp_tree->SetBranchAddress(temp_br_name.Data(), &ThisDigi);
201 for (
int i = 0;
i < temp_tree->GetEntries();
i++) {
202 temp_tree->GetEntry(
i);
204 if (formula->EvalInstance() > 0)
205 fSelectedDigiVector.push_back(fArrayOfDigits->At(
i));
209 return fSelectedDigiVector;
212void BmnFHCalReconstructor::MatchMCTracks()
215 fpFairRootMgr->SetUseFairLinks(kTRUE);
218 for (uint imod = 1; imod < (uint)fBmnFHCalEvent->
GetTotalModules(); ++imod) {
219 auto* mod = fBmnFHCalEvent->
GetModule(imod);
223 for (
int iPt = 0; iPt < fArrayOfPoints->GetEntriesFast(); ++iPt) {
224 auto* point =
static_cast<BmnFHCalPoint*
>(fArrayOfPoints->At(iPt));
226 if (point->GetModuleId() != imod)
229 if (!point->GetPointerToLinks()) {
230 LOG(debug4) <<
"BmnFHCalReconstructor::MatchMCTracks. Link nullptr";
234 mod->AddLinks(*point->GetPointerToLinks());
235 LOG(debug4) << Form(
"BmnFHCalReconstructor::MatchMCTracks. Updated module with id "
236 "%d. Now it has %d links",
237 mod->GetModuleId(), mod->GetNLinks());
242 fpFairRootMgr->SetUseFairLinks(kFALSE);
250 resultTree->Branch(
"FHCalEvent", &fBmnFHCalEvent);
256 printf(
"Work time of BmnFHCalReconstructor: %4.2f sec.\n", fworkTime);
static uint32_t GetModuleId(uint32_t address)
Return Module id from address.
static uint32_t GetModuleType(uint32_t address)
Return Module id from address.
static uint32_t GetSectionId(uint32_t address)
Return Section id from address.
Data class for Bmn FHCal digital signal processing.
uint32_t GetSectionId() const
uint32_t GetModuleId() const
Class for Bmn FHCal data container in event.
BmnFHCalModule * GetModule(uint8_t mod_id)
Module info.
void SetTotalModules(int TotalModules)
Set TotalModules externally.
void reset()
Zero all fiels.
int GetTotalModules() const
Total number of modules.
int GetMaxModuleId() const
void SetSectionEnergy(int sec_id, float Energy)
void SetModuleId(int ModId)
int GetNsections() const
Sections number.
virtual InitStatus Init()
virtual void Exec(Option_t *opt)
BmnFHCalReconstructor(bool isExp, bool isGlobalCoordinates)
void SetRecoCutsFile(TString reco_cuts_file)
virtual void OnlineWrite(const std::unique_ptr< TTree > &resultTree)
Write task resul to tree.
const auto & GetPositionMap() const
BmnStatus ReadGeometryFromGeoManager(bool getGlobalPosition=true)