Skip to content
Snippets Groups Projects
Commit e408d322 authored by Alexandru Bercuci's avatar Alexandru Bercuci
Browse files

add support for channel mask and processing over full TS

parent 1aae7687
No related branches found
No related tags found
1 merge request!1179WIP : TRD2D fixes asked by the community
......@@ -37,7 +37,7 @@ CbmTrdModuleRec2D::CbmTrdModuleRec2D() : CbmTrdModuleRec() {}
CbmTrdModuleRec2D::CbmTrdModuleRec2D(Int_t mod, Int_t ly, Int_t rot) : CbmTrdModuleRec(mod, ly, rot)
{
SetNameTitle(Form("Trd2dReco%d", mod), "Reconstructor for triangular pads.");
// printf("%s (%s)\n", GetName(), GetTitle()); Config(1,0);
// printf("%s (%s)\n", GetName(), GetTitle()); Config(0, 1, 0);
}
//_______________________________________________________________________________
......@@ -51,17 +51,32 @@ Bool_t CbmTrdModuleRec2D::AddDigi(const CbmTrdDigi* d, Int_t id)
* chunk is to have abs(dt)<5 wrt cluster t0
*/
if (CWRITE()) cout << "\nadd @" << id << " " << d->ToString();
int pad = d->GetAddressChannel(), col, row = GetPadRowCol(pad, col), dtime;
uint16_t chT = pad << 1, chR = chT + 1;
int faspAddress = fAsicPar->GetAsicAddress(chT);
const CbmTrdParFasp* p = static_cast<const CbmTrdParFasp*>(fAsicPar->GetAsicPar(faspAddress));
if (!p) {
LOG(error) << GetName() << "::AddDigi : Could not find FASP params for address=" << faspAddress << " @ pad=" << pad;
return false;
}
const CbmTrdParFaspChannel* daqFaspChT = p->GetChannel(pad, 0);
const CbmTrdParFaspChannel* daqFaspChR = p->GetChannel(pad, 1);
if (CWRITE(0)) {
cout << "\nadd @" << id << " " << d->ToString();
daqFaspChT->Print();
daqFaspChR->Print();
}
Int_t ch = d->GetAddressChannel(), col, row = GetPadRowCol(ch, col), dtime;
Double_t t, r = d->GetCharge(t, dtime);
Int_t tm = d->GetTimeDAQ() - fT0, terminator(0);
Int_t tm = d->GetTimeDAQ() - fT0;
if (dtime < 0) tm += dtime; // correct for the difference between tilt and rect
if (r < 1) terminator = 1;
else if (t < 1)
terminator = -1;
if (r < 1 && !daqFaspChR->IsMasked()) chR = 0;
if (t < 1 && !daqFaspChT->IsMasked()) chT = 0;
if (CWRITE()) printf("row[%2d] col[%2d] tm[%2d] terminator[%d]\n", row, col, tm, terminator);
if (CWRITE(0))
printf("row[%2d] col[%2d] tm[%2d] chT[%4d] chR[%4d]\n", row, col, tm, chT * (daqFaspChT->IsMasked() ? -1 : 1),
chR * (daqFaspChR->IsMasked() ? -1 : 1));
CbmTrdCluster* cl(nullptr);
// get the link to saved clusters
......@@ -71,23 +86,23 @@ Bool_t CbmTrdModuleRec2D::AddDigi(const CbmTrdDigi* d, Int_t id)
Bool_t kINSERT(kFALSE);
std::list<CbmTrdCluster*>::iterator itc = fBuffer[row].begin();
for (; itc != fBuffer[row].end(); itc++) {
//if (CWRITE()) cout << (*itc)->ToString();
//if (CWRITE(0)) cout << (*itc)->ToString();
UShort_t tc = (*itc)->GetStartTime();
Int_t dt = tc - tm;
if (dt < -5) continue;
else if (dt < 5) {
if (CWRITE()) printf("match time dt=%d\n", dt);
if ((*itc)->IsChannelInRange(ch) == 0) {
if (!(*itc)->AddDigi(id, ch, terminator, dt)) break;
if (CWRITE(0)) printf("match time dt=%d\n", dt);
if ((*itc)->IsChannelInRange(chT, chR) == 0) {
if (!(*itc)->AddDigi(id, chT, chR, dt)) break;
kINSERT = kTRUE;
if (CWRITE()) cout << " => Cl (Ad) " << (*itc)->ToString();
if (CWRITE(0)) cout << " => Cl (Ad) " << (*itc)->ToString();
break;
}
}
else {
if (CWRITE()) printf("break for time dt=%d\n", dt);
if (CWRITE(0)) printf("break for time dt=%d\n", dt);
break;
}
}
......@@ -95,28 +110,23 @@ Bool_t CbmTrdModuleRec2D::AddDigi(const CbmTrdDigi* d, Int_t id)
if (!kINSERT) {
if (itc != fBuffer[row].end() && itc != fBuffer[row].begin()) {
itc--;
fBuffer[row].insert(itc, cl = new CbmTrdCluster(fModAddress, id, ch, row, tm));
if (CWRITE()) cout << " => Cl (I) ";
fBuffer[row].insert(itc, cl = new CbmTrdCluster(fModAddress, id, chT, chR, row, tm));
if (CWRITE(0)) cout << " => Cl (I) ";
}
else {
fBuffer[row].push_back(cl = new CbmTrdCluster(fModAddress, id, ch, row, tm));
if (CWRITE()) cout << " => Cl (Pb) ";
fBuffer[row].push_back(cl = new CbmTrdCluster(fModAddress, id, chT, chR, row, tm));
if (CWRITE(0)) cout << " => Cl (Pb) ";
}
cl->SetFaspDigis((d->GetType() == CbmTrdDigi::eCbmTrdAsicType::kFASP));
if (terminator < 0) cl->SetProfileStart();
else if (terminator > 0)
cl->SetProfileStop();
if (CWRITE()) cout << cl->ToString();
if (CWRITE(0)) cout << cl->ToString();
}
}
else {
fBuffer[row].push_back(cl = new CbmTrdCluster(fModAddress, id, ch, row, tm));
fBuffer[row].push_back(cl = new CbmTrdCluster(fModAddress, id, chT, chR, row, tm));
cl->SetFaspDigis((d->GetType() == CbmTrdDigi::eCbmTrdAsicType::kFASP));
if (terminator < 0) cl->SetProfileStart();
else if (terminator > 0)
cl->SetProfileStop();
if (CWRITE()) cout << " => Cl (Nw) " << cl->ToString();
if (CWRITE(0)) cout << " => Cl (Nw) " << cl->ToString();
}
fTimeLast = tm;
return kTRUE;
}
......@@ -133,46 +143,98 @@ Int_t CbmTrdModuleRec2D::GetOverThreshold() const
}
//_______________________________________________________________________________
Int_t CbmTrdModuleRec2D::FindClusters()
const CbmTrdParFaspChannel* CbmTrdModuleRec2D::GetFaspChCalibrator(uint16_t ch) const
{
CbmTrdCluster* cl(nullptr);
int faspAddress = fAsicPar->GetAsicAddress(ch);
const CbmTrdParFasp* p = static_cast<const CbmTrdParFasp*>(fAsicPar->GetAsicPar(faspAddress));
if (!p) {
LOG(error) << GetName() << "::GetFaspChCalibrator : Could not find FASP params for address=" << faspAddress
<< " @ ch=" << ch;
return nullptr;
}
return p->GetChannel(ch / 2, ch % 2);
}
// get the link to saved clusters
Int_t ncl(0);
std::list<CbmTrdCluster*>::iterator itc0, itc1;
for (std::map<Int_t, std::list<CbmTrdCluster*>>::iterator ir = fBuffer.begin(); ir != fBuffer.end(); ir++) {
itc0 = (*ir).second.begin();
while ((*ir).second.size() > 1 && itc0 != (*ir).second.end()) { // try merge clusters
itc1 = itc0;
itc1++;
//_______________________________________________________________________________
int CbmTrdModuleRec2D::AddClusterEdges(CbmTrdCluster* cl)
{
bool left = false, right = !left;
int nchRow = 2 * GetNcols(), nchAdd(0);
Bool_t kMERGE = kFALSE;
while (itc1 != (*ir).second.end()) {
if (CWRITE()) cout << " base cl[0] : " << (*itc0)->ToString() << " + cl[1] : " << (*itc1)->ToString();
// check cluster is not at chmb left edge.
if (cl->GetStartCh() > 0 && (cl->GetStartCh() % nchRow != 0)) {
const CbmTrdParFaspChannel* daqCh = GetFaspChCalibrator(cl->GetStartCh() - 1);
if (daqCh && daqCh->IsMasked()) {
cl->AddChannel(left);
nchAdd++;
}
}
// check cluster is not at chmb right edge.
if (cl->GetEndCh() < NFASPMOD * NFASPCH && (cl->GetEndCh() % nchRow != nchRow - 1)) {
const CbmTrdParFaspChannel* daqCh = GetFaspChCalibrator(cl->GetEndCh() + 1);
if (daqCh && daqCh->IsMasked()) {
cl->AddChannel(right);
nchAdd++;
}
}
return nchAdd;
}
//_______________________________________________________________________________
Int_t CbmTrdModuleRec2D::FindClusters(bool clr)
{
int ncl0(0), ncl(0), mcl(0);
std::list<CbmTrdCluster*>::iterator itc0, itc1, itc;
for (auto& clRow : fBuffer) {
if (CWRITE(0)) cout << "\nrow=" << clRow.first << " ncl=" << clRow.second.size() << endl;
// Phase 0 : try merge clusters if more than one on a row
if (clRow.second.size() > 1) {
itc0 = clRow.second.begin();
// TODO look left and right for masked channels. If they exists add them to cluster.
// if (AddClusterEdges(*itc0) && CWRITE(0)) cout << " edge cl[0] : " << (*itc0)->ToString();
itc = std::prev(clRow.second.end());
while (itc0 != itc) { // try merge clusters
if (CWRITE(0)) cout << "->BASE cl[0] : " << (*itc0)->ToString();
bool kMERGE(false);
itc1 = std::next(itc0);
while (itc1 != clRow.second.end()) {
if (CWRITE(0)) cout << " + cl[1] : " << (*itc1)->ToString();
if ((*itc0)->Merge((*itc1))) {
if (CWRITE()) cout << " SUM : " << (*itc0)->ToString();
kMERGE = kTRUE;
if (CWRITE(0)) cout << " SUM : " << (*itc0)->ToString();
kMERGE = true;
delete (*itc1);
itc1 = (*ir).second.erase(itc1);
if (itc1 == (*ir).second.end()) break;
itc1 = clRow.second.erase(itc1);
itc = std::prev(clRow.second.end());
}
else
itc1++;
}
if (!kMERGE) itc0++;
}
}
mcl += clRow.second.size();
for (itc0 = (*ir).second.begin(); itc0 != (*ir).second.end(); itc0++) {
cl = (*itc0);
cl = new ((*fClusters)[ncl++]) CbmTrdCluster(*cl);
cl->SetFaspDigis((*itc0)->HasFaspDigis());
delete (*itc0);
// Phase 1 : copy older clusters from the buffer to the module wise storage
CbmTrdCluster* clDet(nullptr);
for (itc = clRow.second.begin(); itc != clRow.second.end();) {
if (!clr && fTimeLast - (*itc)->GetStartTime() < fTimeWinKeep) {
itc++;
continue;
}
clDet = new ((*fClusters)[ncl++]) CbmTrdCluster(*(*itc));
clDet->SetFaspDigis((*itc)->HasFaspDigis());
delete (*itc);
itc = clRow.second.erase(itc);
}
fBuffer.clear();
}
if (clr) fBuffer.clear();
for (auto clRow : fBuffer)
ncl0 += clRow.second.size();
if (CWRITE(0)) printf("AB :: FindClusters() cls_found = %d cls_write = %d cls_keep = %d\n", mcl, ncl, ncl0);
// printf("fClusters[%p] nCl[%d]\n", (void*)fClusters, fClusters->GetEntriesFast());
// LOG(info) << GetName() << "::FindClusters : " << ncl;
return ncl;
}
......@@ -187,7 +249,7 @@ Bool_t CbmTrdModuleRec2D::PreProcessHits()
*/
Int_t nhits = fHits->GetEntriesFast();
if (CWRITE()) LOG(info) << GetName() << "::PreProcessHits(" << nhits << ")";
if (CWRITE(1)) LOG(info) << GetName() << "::PreProcessHits(" << nhits << ")";
CbmTrdHit* hit(nullptr);
for (Int_t ih(0); ih < nhits; ih++) {
......@@ -196,7 +258,7 @@ Bool_t CbmTrdModuleRec2D::PreProcessHits()
nhits += Deconvolute(hit);
}
nhits = fHits->GetEntriesFast();
if (CWRITE()) LOG(info) << GetName() << "::Deconvolute(" << nhits << ")";
if (CWRITE(1)) LOG(info) << GetName() << "::Deconvolute(" << nhits << ")";
return kTRUE;
}
......@@ -212,7 +274,7 @@ Bool_t CbmTrdModuleRec2D::PostProcessHits()
/** Steering routine for classifying hits and apply further analysis
* -> hit merging for row-cross (see RowCross)
*/
return true;
CbmTrdHit *h0(nullptr), *h1(nullptr);
Int_t a0, nhits = fHits->GetEntriesFast();
......@@ -241,7 +303,7 @@ Bool_t CbmTrdModuleRec2D::PostProcessHits()
// call the working algorithm
if (MergeHits(h0, a0)) h0->SetRowCross(h1);
if (CWRITE()) {
if (CWRITE(1)) {
cout << ih << " : " << h0->ToString();
cout << jh << " : " << h1->ToString();
cout << "\n" << endl;
......@@ -265,7 +327,7 @@ Bool_t CbmTrdModuleRec2D::PostProcessHits()
}
fDigis.clear();
if (CWRITE()) LOG(info) << GetName() << "::MergeHits(" << nhits << ")";
if (CWRITE(1)) LOG(info) << GetName() << "::MergeHits(" << nhits << ")";
return kTRUE;
}
......@@ -526,7 +588,7 @@ Bool_t CbmTrdModuleRec2D::MergeHits(CbmTrdHit* h, Int_t a0)
h->SetMaxType(IsMaxTilt());
h->SetOverFlow(HasOvf());
if (CWRITE()) {
if (CWRITE(1)) {
printf("-> loc[%6.2f %6.2f %6.2f] err[%6.2f %6.2f %6.2f]\n", local_pad[0], local_pad[1], local_pad[2],
local_pad_err[0], local_pad_err[1], local_pad_err[2]);
printf("REC col[%2d] row[%2d] dx[%7.3f(pw) %7.3f(cm)] x[%7.2f] y[%7.2f] dy[%5.2f] t0[%llu]\n", vcM, vrM,
......@@ -688,7 +750,7 @@ Bool_t CbmTrdModuleRec2D::BuildHit(CbmTrdHit* h)
fgEdep->SetPoint(ip, xc, 0);
fgEdep->SetPointError(ip, 0., 300);
}
//if(CWRITE()) fgEdep->Print();
//if(CWRITE(1)) fgEdep->Print();
}
Double_t e(0.), xlo(*vx.begin()), xhi(*vx.rbegin());
......@@ -738,7 +800,7 @@ Bool_t CbmTrdModuleRec2D::BuildHit(CbmTrdHit* h)
h->SetMaxType(IsMaxTilt());
h->SetOverFlow(HasOvf());
if (CWRITE()) {
if (CWRITE(1)) {
printf("-> loc[%6.2f %6.2f %6.2f] err[%6.2f %6.2f %6.2f]\n", local_pad[0], local_pad[1], local_pad[2],
local_pad_err[0], local_pad_err[1], local_pad_err[2]);
printf("REC col[%2d] row[%2d] x[%7.2f] dx[%5.2f] y[%7.2f] dy[%5.2f] t0[%llu]\n", vcM, vrM, global[0], dx, global[1],
......@@ -790,7 +852,7 @@ CbmTrdHit* CbmTrdModuleRec2D::MakeHit(Int_t ic, const CbmTrdCluster* cl, std::ve
hf->Draw("p");
}
if (CWRITE()) cout << cl->ToString();
if (CWRITE(1)) cout << cl->ToString();
if (!LoadDigis(digis, ic)) return nullptr;
if (!ProjectDigis(ic)) return nullptr;
Int_t nofHits = fHits->GetEntriesFast();
......@@ -799,7 +861,7 @@ CbmTrdHit* CbmTrdModuleRec2D::MakeHit(Int_t ic, const CbmTrdCluster* cl, std::ve
hit->SetRefId(ic);
//hit->SetMatch();
BuildHit(hit);
if (CWRITE()) cout << hit->ToString();
if (CWRITE(1)) cout << hit->ToString();
if (CDRAW()) DrawHit(hit);
return hit;
}
......@@ -851,7 +913,7 @@ Double_t CbmTrdModuleRec2D::GetXoffset(Int_t n0) const
dx += vs[ir] * vx[ir];
}
if (TMath::Abs(R) > 0) return dx / R;
LOG(warn) << GetName() << "::GetXoffset : Unexpected null sum.";
LOG(warn) << GetName() << "::GetXoffset : Null total charge for hit size " << n;
return 0.;
}
......@@ -867,7 +929,8 @@ Double_t CbmTrdModuleRec2D::GetYoffset(Int_t n0) const
}
}
if (TMath::Abs(T) > 0) return dy / T;
LOG(warn) << GetName() << "::GetYoffset : Unexpected null sum.";
LOG(warn) << GetName() << "::GetYoffset : Null total charge for hit size " << n;
//if (CWRITE(1))
return 0.;
}
......@@ -1089,7 +1152,7 @@ Int_t CbmTrdModuleRec2D::ProjectDigis(Int_t cid, Int_t cjd)
dg = (*i);
dg0 = nullptr;
dg1 = nullptr;
if (CWRITE()) cout << "dg0 :" << dg->ToString();
if (CWRITE(1)) cout << "dg0 :" << dg->ToString();
// initialize
if (col < 0) {
......@@ -1138,7 +1201,7 @@ Int_t CbmTrdModuleRec2D::ProjectDigis(Int_t cid, Int_t cjd)
}
if (dg0) i1++;
}
if (CWRITE()) {
if (CWRITE(1)) {
if (dg0) cout << "dgR :" << dg0->ToString();
if (dg1) cout << "dgT :" << dg1->ToString();
cout << "-------------------------------------" << endl;
......@@ -1280,7 +1343,7 @@ Int_t CbmTrdModuleRec2D::ProjectDigis(Int_t cid, Int_t cjd)
else
SetBiasY(0);
if (CWRITE()) {
if (CWRITE(1)) {
printf("t0[clk]=%llu col[%2d] row[%2d] sz[%d] OVF[%c] %c", vt0, vcM, vrM, Int_t(vs.size() - 2),
(ovf < 0 ? 'y' : 'n'), IsOpenLeft() ? '(' : '[');
if (IsSymmHit()) {
......@@ -1339,7 +1402,7 @@ Int_t CbmTrdModuleRec2D::LoadDigis(vector<const CbmTrdDigi*>* digis, vector<CbmT
dg = (*idgM);
idgM++;
}
if (CWRITE()) cout << dg->ToString();
if (CWRITE(1)) cout << dg->ToString();
r = dg->GetCharge(t, dt);
if (t0 == 0) t0 = dg->GetTimeDAQ(); // set arbitrary t0 to avoid double digis loop
if (col0 < 0) GetPadRowCol(dg->GetAddressChannel(), col0); // initialilze
......@@ -1497,7 +1560,7 @@ Int_t CbmTrdModuleRec2D::LoadDigisRC(vector<const CbmTrdDigi*>* digis, const Int
col = col0;
step = 1;
}
if (CWRITE()) printf("col0[%d] col1[%d] step[%2d]\n", col0, col1, step);
if (CWRITE(1)) printf("col0[%d] col1[%d] step[%2d]\n", col0, col1, step);
const CbmTrdDigi *dg0(nullptr), *dg1(nullptr), *dg10(nullptr);
// always loop on the largest cluster
......@@ -1505,7 +1568,7 @@ Int_t CbmTrdModuleRec2D::LoadDigisRC(vector<const CbmTrdDigi*>* digis, const Int
dg0 = (*i0);
dg1 = nullptr;
dg10 = nullptr;
if (CWRITE()) cout << "dg0 :" << dg0->ToString();
if (CWRITE(1)) cout << "dg0 :" << dg0->ToString();
r = dg0->GetCharge(t, dt);
if (t > 0) t -= CbmTrdFASP::GetBaselineCorr();
......@@ -1554,7 +1617,7 @@ Int_t CbmTrdModuleRec2D::LoadDigisRC(vector<const CbmTrdDigi*>* digis, const Int
}
if (dg1) i1++;
if (CWRITE()) {
if (CWRITE(1)) {
if (dg1) cout << "dgR :" << dg1->ToString();
if (dg10) cout << "dgT :" << dg10->ToString();
cout << "-------------------------------------" << endl;
......@@ -1669,7 +1732,7 @@ Bool_t CbmTrdModuleRec2D::MergeDigis(vector<const CbmTrdDigi*>* digis, vector<Cb
dgM->SetCharge(t, r, dt);
Int_t rtrg(dgR->GetTriggerType() & 2), ttrg(dgT->GetTriggerType() & 1);
dgM->SetTriggerType(rtrg | ttrg); //merge the triggers
if (CWRITE()) {
if (CWRITE(1)) {
cout << "MERGE" << endl;
cout << dgT->ToString();
cout << dgR->ToString();
......
......@@ -19,6 +19,7 @@ using std::map;
using std::vector;
class TGraphErrors;
class CbmTrdDigiRec;
class CbmTrdParFaspChannel;
class TF1;
/** @class CbmTrdModuleRec2D
** @brief Cluster finding and hit reconstruction algorithms for the TRD(2D) module.
......@@ -37,9 +38,10 @@ class CbmTrdModuleRec2D : public CbmTrdModuleRec {
public:
enum ECbmTrdModuleRec2D
{
kVerbose = 0, ///< steer verbosity on/off
kDraw = 1, ///< steer graphic representation on/off
kHelpers = 2 ///< use helper graph for time and energy estimation
kVerbCluster = 0, ///< steer clusterizer verbosity on/off
kVerbReco = 1, ///< steer reconstructor verbosity on/off
kDraw = 2, ///< steer graphic representation on/off
kHelpers = 3 ///< use helper graph for time and energy estimation
};
/** \brief Default constructor.*/
......@@ -51,11 +53,6 @@ public:
/** \brief Add digi to local module **/
virtual Bool_t AddDigi(const CbmTrdDigi* d, Int_t id);
/** \brief Config task with the following settings
* \param[in] v verbosity toggle
* \param[in] d drawing toggle
*/
virtual inline void Config(Bool_t vb, Bool_t dw);
virtual void DrawHit(CbmTrdHit*) const { ; }
/** \brief Count RO channels (R or T) with data**/
virtual Int_t GetOverThreshold() const;
......@@ -64,7 +61,7 @@ public:
/** \brief Finalize hits (merge RC hits, etc)**/
virtual Bool_t PostProcessHits();
/** \brief Finalize clusters **/
virtual Int_t FindClusters();
virtual Int_t FindClusters(bool clr);
/** \brief Steering routine for building hits **/
virtual Bool_t MakeHits();
/** \brief Steering routine for converting cluster to hit **/
......@@ -135,10 +132,25 @@ private:
CbmTrdModuleRec2D(const CbmTrdModuleRec2D& ref);
const CbmTrdModuleRec2D& operator=(const CbmTrdModuleRec2D& ref);
/** \brief Config task with the following settings
* \param[in] vcl toggle verbosity clusterizer
* \param[in] vrc toggle verbosity reconsructor
* \param[in] dw drawing toggle
*/
virtual inline void Config(Bool_t vcl, Bool_t vrc, Bool_t dw);
Bool_t CDRAW() const { return TESTBIT(fConfigMap, ECbmTrdModuleRec2D::kDraw); }
Bool_t CWRITE() const { return TESTBIT(fConfigMap, ECbmTrdModuleRec2D::kVerbose); }
Bool_t CWRITE(int level) const
{
if (level) return TESTBIT(fConfigMap, ECbmTrdModuleRec2D::kVerbReco);
else
return TESTBIT(fConfigMap, ECbmTrdModuleRec2D::kVerbCluster);
}
Bool_t CHELPERS() const { return TESTBIT(fConfigMap, ECbmTrdModuleRec2D::kHelpers); }
/** \brief Add left and right edge channels to the cluster in case this are masked channels
* \return no of edges added to cluster
*/
int AddClusterEdges(CbmTrdCluster* cl);
/** \brief Implement cuts for hit convolution definition
* \param[in] h hit to be analysed.
* \return TRUE if double cluster
......@@ -151,6 +163,10 @@ private:
Bool_t Deconvolute(CbmTrdHit* h);
Double_t GetXoffset(Int_t n0 = 0) const;
Double_t GetYoffset(Int_t n0 = 0) const;
/** \brief Retrive FASP ch calibrator by RO ch number in the module
* \param[in] ch channel id in the current module.
*/
const CbmTrdParFaspChannel* GetFaspChCalibrator(uint16_t ch) const;
/** \brief Load digis info into local data structures
* \param[in] digis initial digis list shrinked for incomplete digis.
* \param[in] vdgM list of merged digis
......@@ -201,6 +217,8 @@ private:
UChar_t fConfigMap = 0; //! task configuration settings
ULong64_t fT0 = 0; //! start time of event/time slice [clk]
uint fTimeLast = 0; //! time of last digi processed in module [clk]
uint fTimeWinKeep = 11; //! time interval to still keep clusters in buffer [clk]
std::map<Int_t, std::list<CbmTrdCluster*>> fBuffer; //row-wise organized clusters
std::map<Int_t, vector<CbmTrdDigiRec*>> fDigis; //!cluster-wise organized calibrated digi
// working representation of a hit on which the reconstruction is performed
......@@ -230,12 +248,16 @@ private:
2) // Triangular pad module; Cluster finding and hit reconstruction algorithms
};
void CbmTrdModuleRec2D::Config(Bool_t vb, Bool_t dw)
void CbmTrdModuleRec2D::Config(Bool_t vcl, Bool_t vrc, Bool_t dw)
{
if (vb) SETBIT(fConfigMap, ECbmTrdModuleRec2D::kVerbose);
if (vcl) SETBIT(fConfigMap, ECbmTrdModuleRec2D::kVerbCluster);
else
CLRBIT(fConfigMap, ECbmTrdModuleRec2D::kVerbCluster);
printf("CbmTrdModuleRec2D::kVerbCluster[%c]\n", CWRITE(0) ? 'y' : 'n');
if (vrc) SETBIT(fConfigMap, ECbmTrdModuleRec2D::kVerbReco);
else
CLRBIT(fConfigMap, ECbmTrdModuleRec2D::kVerbose);
printf("CbmTrdModuleRec2D::Verbose[%c]\n", CWRITE() ? 'y' : 'n');
CLRBIT(fConfigMap, ECbmTrdModuleRec2D::kVerbReco);
printf("CbmTrdModuleRec2D::kVerbReco[%c]\n", CWRITE(1) ? 'y' : 'n');
if (dw) SETBIT(fConfigMap, ECbmTrdModuleRec2D::kDraw);
else
CLRBIT(fConfigMap, ECbmTrdModuleRec2D::kDraw);
......
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Please register or to comment