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

add support for channel masking in reconstruction. Only available for

TRD2D
parent 70f419a3
No related branches found
No related tags found
No related merge requests found
...@@ -39,11 +39,11 @@ CbmTrdCluster::CbmTrdCluster(const std::vector<int32_t>& indices, int32_t addres ...@@ -39,11 +39,11 @@ CbmTrdCluster::CbmTrdCluster(const std::vector<int32_t>& indices, int32_t addres
} }
//____________________________________________________________________ //____________________________________________________________________
CbmTrdCluster::CbmTrdCluster(int32_t address, int32_t idx, int32_t ch, int32_t row, int32_t time) CbmTrdCluster::CbmTrdCluster(int32_t address, int32_t idx, uint16_t chT, uint16_t chR, int32_t row, int32_t time)
: CbmCluster() : CbmCluster()
{ {
ReInit(address, row, time); ReInit(address, row, time);
AddDigi(idx, ch); AddDigi(idx, chT, chR);
} }
//____________________________________________________________________ //____________________________________________________________________
...@@ -62,42 +62,67 @@ CbmTrdCluster& CbmTrdCluster::operator=(const CbmTrdCluster& ref) ...@@ -62,42 +62,67 @@ CbmTrdCluster& CbmTrdCluster::operator=(const CbmTrdCluster& ref)
} }
//____________________________________________________________________ //____________________________________________________________________
bool CbmTrdCluster::AddDigi(int32_t idx, int32_t channel, int32_t terminator, int32_t dt) bool CbmTrdCluster::AddChannel(bool r)
{ {
/** Extend basic functionality of CbmCluster::AddDigi(). if (!r) {
* If channel>=0 add this info to channel map. if (!fStartCh) return false;
fStartCh--;
}
fNCols++;
return true;
}
//____________________________________________________________________
bool CbmTrdCluster::AddDigi(int32_t idx, uint16_t chT, uint16_t chR, int32_t dt)
{
/** Extend basic functionality of CbmCluster::AddDigi() for the case of 2D.
* If chT < 0 use the basic functionality [default].
*
* For the 2D the parameters are intergpreted as follows
* chT : tilted paired channel [default 0x0fffffff]
* chR : rectangular paired channel
* dt : offset in clks of the prompt signal
*
* if chT and chR positive the (chT, chR) are interpreted as the 2 channels
* of the digi specific to the 2D version. The following specific cases
* can be distinguished :
* - ch == 0 : no data, cluster signal sequence terminator
* - ch == -ch : no data, channel masked in HW
*/ */
if (channel < 0) { // basic functionality for rectangular pads if (chT == 0xffff) { // basic functionality for rectangular pads
CbmCluster::AddDigi(idx); CbmCluster::AddDigi(idx);
return true; return true;
} }
if (channel >= 0xffff) LOG(warn) << GetName() << "::AddDigi: pad-channel truncated to 2bytes."; uint16_t chMin = (chT != 0 ? chT : chR), chMax = (chR != 0 ? chR : chT);
// assume triangular pads only // assume triangular pads only
if (!fNCols) { // first digi if (!fNCols) { // first digi
fStartCh = channel; fStartCh = chMin;
CbmCluster::AddDigi(idx); CbmCluster::AddDigi(idx);
if (terminator < 0) SetProfileStart();
else if (terminator > 0)
SetProfileStop();
} }
else if (channel > GetEndCh()) { // digi @ end else if (chMin > GetEndCh()) { // digi @ end
if (HasOpenStop()) return false; //if (HasStop()) return false;
CbmCluster::AddDigi(idx); CbmCluster::AddDigi(idx);
if (terminator > 0) SetProfileStop();
} }
else if (channel < fStartCh) { // digi @ beginning else if (chMax < fStartCh) { // digi @ beginning
if (HasOpenStart()) return false; //if (HasStart()) return false;
fStartCh = channel; fStartCh = chMin;
vector<int32_t> vec = GetDigis(); vector<int32_t> vec = GetDigis();
ClearDigis(); ClearDigis();
CbmCluster::AddDigi(idx); CbmCluster::AddDigi(idx);
AddDigis(vec); AddDigis(vec);
if (terminator < 0) SetProfileStart();
} }
fNCols++; int nch(0);
if (chT == 0) SetStart();
else
nch++;
if (chR == 0) SetStop();
else
nch++;
fNCols += nch;
if (dt > 0) fStartTime -= dt; if (dt > 0) fStartTime -= dt;
return true; return true;
...@@ -122,45 +147,27 @@ void CbmTrdCluster::ReInit(int32_t address, int32_t row, int32_t time) ...@@ -122,45 +147,27 @@ void CbmTrdCluster::ReInit(int32_t address, int32_t row, int32_t time)
// check truncation // check truncation
if (row >= 0x1f) LOG(warn) << GetName() << "::ReInit: pad-row truncated to 5bits."; if (row >= 0x1f) LOG(warn) << GetName() << "::ReInit: pad-row truncated to 5bits.";
SetNRows(row); SetNRows(row);
SetProfileStart(0); SetStart(false);
SetProfileStop(0); SetStop(false);
if (std::abs(time) >= 0x7fffffff) LOG(warn) << GetName() << "::ReInit: buffer time truncated to 4bytes."; if (std::abs(time) >= 0x7fffffff) LOG(warn) << GetName() << "::ReInit: buffer time truncated to 4bytes.";
fStartTime = time; fStartTime = time;
} }
//____________________________________________________________________ //____________________________________________________________________
int32_t CbmTrdCluster::IsChannelInRange(int32_t ch) const int32_t CbmTrdCluster::IsChannelInRange(uint16_t chT, uint16_t chR) const
{ {
if (!fNCols) return -2; if (!fNCols) return -2;
// if(IsTerminatedLeft() && fAddressCh[0]>ch) return -1; // if(IsTerminatedLeft() && fAddressCh[0]>ch) return -1;
// if(IsTerminatedRight() && fAddressCh[clSize-1]<ch) return 1; // if(IsTerminatedRight() && fAddressCh[clSize-1]<ch) return 1;
if (fStartCh > ch + 1) return -1; uint16_t chMin = (chT != 0 ? chT : chR), chMax = (chR != 0 ? chR : chT);
if (fStartCh + fNCols < ch) return 1; if (fStartCh > chMax + 1) return -1;
if (fStartCh + fNCols < chMin) return 1;
return 0; return 0;
} }
//____________________________________________________________________ //____________________________________________________________________
uint16_t CbmTrdCluster::GetSize() const bool CbmTrdCluster::Merge(CbmTrdCluster* second)
{
uint16_t size(0);
if (HasFaspDigis()) {
size = GetNCols() << 1;
if (HasOpenStart()) size--;
if (HasOpenStop()) size--;
if (size <= 0) {
LOG(warn) << GetName() << "::GetSize: Fasp cluster meta-info corrupt.";
std::cout << ToString();
return 0;
}
}
else
size = GetNCols();
return size;
}
//____________________________________________________________________
bool CbmTrdCluster::Merge(CbmTrdCluster* second, bool typ)
{ {
if (GetRow() != second->GetRow()) return false; if (GetRow() != second->GetRow()) return false;
// time difference condition // time difference condition
...@@ -170,8 +177,8 @@ bool CbmTrdCluster::Merge(CbmTrdCluster* second, bool typ) ...@@ -170,8 +177,8 @@ bool CbmTrdCluster::Merge(CbmTrdCluster* second, bool typ)
else if (abs(int32_t(second->fStartTime - fStartTime)) > 20) else if (abs(int32_t(second->fStartTime - fStartTime)) > 20)
return false; return false;
// look before current // look before current
if (second->fStartCh + second->fNCols == fStartCh && !second->HasOpenStop() && !HasOpenStart()) { if (second->fStartCh + second->fNCols == fStartCh) {
// cout<<"Merge before with "<<second->ToString(); // std::cout<<"Merge before with "<<second->ToString();
fStartCh = second->fStartCh; fStartCh = second->fStartCh;
fNCols += second->fNCols; fNCols += second->fNCols;
fStartTime = std::min(fStartTime, second->fStartTime); fStartTime = std::min(fStartTime, second->fStartTime);
...@@ -180,47 +187,19 @@ bool CbmTrdCluster::Merge(CbmTrdCluster* second, bool typ) ...@@ -180,47 +187,19 @@ bool CbmTrdCluster::Merge(CbmTrdCluster* second, bool typ)
ClearDigis(); ClearDigis();
AddDigis(second->GetDigis()); AddDigis(second->GetDigis());
AddDigis(vec); AddDigis(vec);
if (second->HasOpenStart()) SetProfileStart(); if (second->HasStart()) SetStart();
return true; return true;
} }
// special care for clusters which can be merged also on pairing neighboring on 2D read-out
if (typ) {
if ((second->fStartCh + second->fNCols - 1 == fStartCh) && second->HasOpenStop()
&& HasOpenStart()) { // need to merge digi
fStartCh = second->fStartCh;
fNCols += second->fNCols - 1;
fStartTime = std::min(fStartTime, second->fStartTime);
vector<int32_t> vec = GetDigis();
ClearDigis();
AddDigis(second->GetDigis());
AddDigis(vec);
SetProfileStart(false);
return true;
}
}
// look after current // look after current
if (fStartCh + fNCols == second->fStartCh && !HasOpenStop() && !second->HasOpenStart()) { if (fStartCh + fNCols == second->fStartCh) {
// cout<<"Merge after with "<<second->ToString(); // std::cout<<"Merge after with "<<second->ToString();
fNCols += second->fNCols; fNCols += second->fNCols;
fStartTime = std::min(fStartTime, second->fStartTime); fStartTime = std::min(fStartTime, second->fStartTime);
AddDigis(second->GetDigis()); AddDigis(second->GetDigis());
if (second->HasOpenStop()) SetProfileStop(); if (second->HasStop()) SetStop();
return true; return true;
} }
// special care for clusters which can be merged also on pairing neighboring on 2D read-out
if (typ) {
if ((fStartCh + fNCols - 1 == second->fStartCh) && HasOpenStop() && second->HasOpenStart()) { // need to merge digi
fNCols += second->fNCols - 1;
fStartTime = std::min(fStartTime, second->fStartTime);
AddDigis(second->GetDigis());
SetProfileStop(false);
return true;
}
}
return false; return false;
} }
...@@ -231,10 +210,10 @@ string CbmTrdCluster::ToString() const ...@@ -231,10 +210,10 @@ string CbmTrdCluster::ToString() const
ss << CbmCluster::ToString(); ss << CbmCluster::ToString();
ss << "CbmTrdCluster: mod=" << GetAddress() << " row=" << (int32_t) GetRow() << " " ss << "CbmTrdCluster: mod=" << GetAddress() << " row=" << (int32_t) GetRow() << " "
<< (HasFaspDigis() ? "Fasp_" : "Spadic_") << "Chs="; << (HasFaspDigis() ? "Fasp_" : "Spadic_") << "Chs=";
ss << (HasOpenStart() ? "/" : "|"); ss << (HasStart() ? "|" : "/");
for (int32_t i(0); i < fNCols; i++) for (int32_t i(0); i < fNCols; i++)
ss << fStartCh + i << " "; ss << fStartCh + i << " ";
ss << (HasOpenStop() ? "\\" : "|"); ss << (HasStop() ? "|" : "/");
ss << endl; ss << endl;
return ss.str(); return ss.str();
} }
......
...@@ -30,9 +30,9 @@ public: ...@@ -30,9 +30,9 @@ public:
{ {
kFasp = 5 ///< set type of FEE digis contained kFasp = 5 ///< set type of FEE digis contained
, ,
kProfileStart ///< only for triangular if no T in first col kStart ///< only for triangular if no T in first col
, ,
kProfileStop ///< only for triangular if no R in last col kStop ///< only for triangular if no R in last col
}; };
/** /**
* \brief Default constructor. * \brief Default constructor.
...@@ -41,14 +41,15 @@ public: ...@@ -41,14 +41,15 @@ public:
CbmTrdCluster(const CbmTrdCluster& ref); CbmTrdCluster(const CbmTrdCluster& ref);
CbmTrdCluster(const std::vector<int32_t>& indices, int32_t address); CbmTrdCluster(const std::vector<int32_t>& indices, int32_t address);
/** /**
* \brief Constructor starting from first digit. * \brief Constructor starting from first digit (FASP specific).
* \param[in] address global module address * \param[in] address global module address
* \param[in] idx global digi index in the TClonesArray * \param[in] idx global digi index in the TClonesArray
* \param[in] ch RO channel address within the module * \param[in] chT RO channel address within the module for tilt pairing
* \param[in] chR RO channel address within the module for rect pairing
* \param[in] r module row for the RO channel * \param[in] r module row for the RO channel
* \param[in] time relative buffer DAQ time * \param[in] time relative buffer DAQ time
*/ */
CbmTrdCluster(int32_t address, int32_t idx, int32_t ch, int32_t r, int32_t time); CbmTrdCluster(int32_t address, int32_t idx, uint16_t chT, uint16_t chR, int32_t r, int32_t time);
/** /**
* \brief Destructor. * \brief Destructor.
*/ */
...@@ -56,14 +57,19 @@ public: ...@@ -56,14 +57,19 @@ public:
CbmTrdCluster& operator=(const CbmTrdCluster& ref); CbmTrdCluster& operator=(const CbmTrdCluster& ref);
/** \brief Append a channel to cluster edge. The usage is to account for the masked channels.
* The mask status is assumed to be performed in the calling function.
* \param[in] r by default apply to the right edge. If false apply to left
*/
bool AddChannel(bool r = true);
/** \brief Append digi to cluster /** \brief Append digi to cluster
* \param[in] idx index of digi in TClonesArray * \param[in] idx index of digi in TClonesArray
* \param[in] channel RO channel for digi * \param[in] chT RO channel for digi (tilt pairing for FASP) default 0xffff (SPADIC)
* \param[in] terminator state of digi for triangular pads: 0 if complete, -/+ if T/R misses * \param[in] chR RO channel for rect pairing (only for FASP)
* \param[in] dt update start time of cluster if current digi is prompt * \param[in] dt update start time of cluster if current digi is prompt
* \return true if successful * \return true if successful
*/ */
bool AddDigi(int32_t idx, int32_t channel = -1, int32_t terminator = 0, int32_t dt = 0); bool AddDigi(int32_t idx, uint16_t chT = 0xffff, uint16_t chR = 0, int32_t dt = 0);
/** \brief reset cluster data*/ /** \brief reset cluster data*/
void Clear(Option_t*); void Clear(Option_t*);
/** Accessors **/ /** Accessors **/
...@@ -71,24 +77,25 @@ public: ...@@ -71,24 +77,25 @@ public:
uint16_t GetNRows() const { return fNRows & 0x1f; } uint16_t GetNRows() const { return fNRows & 0x1f; }
uint16_t GetEndCh() const { return fStartCh + fNCols - 1; } uint16_t GetEndCh() const { return fStartCh + fNCols - 1; }
uint16_t GetRow() const { return GetNRows(); } uint16_t GetRow() const { return GetNRows(); }
uint16_t GetSize() const; uint16_t GetSize() const { return GetNCols(); }
uint16_t GetStartCh() const { return fStartCh; } uint16_t GetStartCh() const { return fStartCh; }
uint32_t GetStartTime() const { return fStartTime; } uint32_t GetStartTime() const { return fStartTime; }
bool HasFaspDigis() const { return TESTBIT(fNRows, kFasp); } bool HasFaspDigis() const { return TESTBIT(fNRows, kFasp); }
bool HasOpenStart() const { return TESTBIT(fNRows, kProfileStart); } bool HasStart() const { return TESTBIT(fNRows, kStart); }
bool HasOpenStop() const { return TESTBIT(fNRows, kProfileStop); } bool HasStop() const { return TESTBIT(fNRows, kStop); }
/** \brief Query on RO channel list /** \brief Query on RO channels list
* \param[in] channel RO channel for digi * \param[in] chT tilted RO channel for digi
* \param[in] chR rectangular RO channel for digi
* \return -1 before range; 0 in range; 1 after range; -2 cluster empty of digits * \return -1 before range; 0 in range; 1 after range; -2 cluster empty of digits
*/ */
int32_t IsChannelInRange(int32_t ch) const; int32_t IsChannelInRange(uint16_t chT, uint16_t chR) const;
/** \brief Merge current cluster with info from second /** \brief Merge current cluster with info from second
* \param[in] second cluster to be added * \param[in] second cluster to be added
* \param[in] typ the type of pad-plane of the source chamber; true if Trd2d * \param[in] typ the type of pad-plane of the source chamber; true if Trd2d
* \return success or fail * \return success or fail
*/ */
bool Merge(CbmTrdCluster* second, bool typ = true); bool Merge(CbmTrdCluster* second);
/** \brief Initialize basic parameters of the cluster /** \brief Initialize basic parameters of the cluster
* \param[in] address global module address * \param[in] address global module address
* \param[in] row cluster row in the module * \param[in] row cluster row in the module
...@@ -103,8 +110,8 @@ public: ...@@ -103,8 +110,8 @@ public:
fNRows |= (nrows & 0x1f); fNRows |= (nrows & 0x1f);
} }
void SetFaspDigis(bool set = true) { set ? SETBIT(fNRows, kFasp) : CLRBIT(fNRows, kFasp); } void SetFaspDigis(bool set = true) { set ? SETBIT(fNRows, kFasp) : CLRBIT(fNRows, kFasp); }
void SetProfileStart(bool set = true) { set ? SETBIT(fNRows, kProfileStart) : CLRBIT(fNRows, kProfileStart); } void SetStart(bool set = true) { set ? SETBIT(fNRows, kStart) : CLRBIT(fNRows, kStart); }
void SetProfileStop(bool set = true) { set ? SETBIT(fNRows, kProfileStop) : CLRBIT(fNRows, kProfileStop); } void SetStop(bool set = true) { set ? SETBIT(fNRows, kStop) : CLRBIT(fNRows, kStop); }
/** \brief Extended functionality*/ /** \brief Extended functionality*/
virtual std::string ToString() const; virtual std::string ToString() const;
......
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