CMS 3D CMS Logo

Public Member Functions | Private Attributes

DTDDUUnpacker Class Reference

#include <DTDDUUnpacker.h>

Inheritance diagram for DTDDUUnpacker:
DTUnpacker

List of all members.

Public Member Functions

 DTDDUUnpacker (const edm::ParameterSet &ps)
 Constructor.
virtual void interpretRawData (const unsigned int *index, int datasize, int dduID, edm::ESHandle< DTReadOutMapping > &mapping, std::auto_ptr< DTDigiCollection > &product, std::auto_ptr< DTLocalTriggerCollection > &product2, uint16_t rosList=0)
virtual ~DTDDUUnpacker ()
 Destructor.

Private Attributes

DTDataMonitorInterfacedataMonitor
const edm::ParameterSet dduPSet
bool debug
bool localDAQ
 if data are read locally, status words are swapped
bool performDataIntegrityMonitor
 perform DQM for DDU
DTROS25Unpackerros25Unpacker

Detailed Description

The unpacker for DTs' FED.

Date:
2007/09/04 08:07:26
Revision:
1.4
Author:
M. Zanetti INFN Padova FRC 060906

Definition at line 19 of file DTDDUUnpacker.h.


Constructor & Destructor Documentation

DTDDUUnpacker::DTDDUUnpacker ( const edm::ParameterSet ps)

Constructor.

Definition at line 28 of file DTDDUUnpacker.cc.

References dataMonitor, dduPSet, debug, edm::ParameterSet::getParameter(), edm::ParameterSet::getUntrackedParameter(), edm::Service< T >::isAvailable(), localDAQ, cppFunctionSkipper::operator, performDataIntegrityMonitor, and ros25Unpacker.

                                                      : dduPSet(ps) { 
  
  // the ROS unpacker
  ros25Unpacker = new DTROS25Unpacker(dduPSet.getParameter<edm::ParameterSet>("rosParameters"));
  
  // parameters
  localDAQ = dduPSet.getUntrackedParameter<bool>("localDAQ",false);
  performDataIntegrityMonitor = dduPSet.getUntrackedParameter<bool>("performDataIntegrityMonitor",false);
  debug = dduPSet.getUntrackedParameter<bool>("debug",false);

  // enable DQM if Service is available
  if(performDataIntegrityMonitor) {
    if (edm::Service<DTDataMonitorInterface>().isAvailable()) {
      dataMonitor = edm::Service<DTDataMonitorInterface>().operator->(); 
    } else {
      LogWarning("DTRawToDigi|DTDDUUnpacker") << 
        "[DTDDUUnpacker] WARNING! Data Integrity Monitoring requested but no DTDataMonitorInterface Service available" << endl;
      performDataIntegrityMonitor = false;
    }
  }

}
DTDDUUnpacker::~DTDDUUnpacker ( ) [virtual]

Destructor.

Definition at line 52 of file DTDDUUnpacker.cc.

References ros25Unpacker.

                              {
  delete ros25Unpacker;
}

Member Function Documentation

void DTDDUUnpacker::interpretRawData ( const unsigned int *  index,
int  datasize,
int  dduID,
edm::ESHandle< DTReadOutMapping > &  mapping,
std::auto_ptr< DTDigiCollection > &  product,
std::auto_ptr< DTLocalTriggerCollection > &  product2,
uint16_t  rosList = 0 
) [virtual]

Unpacking method. index is the pointer to the beginning of the buffer. datasize is the size of the buffer in bytes

Implements DTUnpacker.

Definition at line 57 of file DTDDUUnpacker.cc.

References DTDDUData::addDDUStatusWord(), DTDDUData::addROSStatusWord(), FEDHeader::bxID(), FEDTrailer::check(), FEDHeader::check(), DTDDUData::checkCRCBit(), gather_cfg::cout, dataMonitor, debug, DTROS25Unpacker::getROSsControlData(), DTROS25Unpacker::interpretRawData(), FEDTrailer::lenght(), localDAQ, FEDHeader::lvl1ID(), performDataIntegrityMonitor, DTDataMonitorInterface::processFED(), ros25Unpacker, and DTDDUSecondStatusWord::rosList().

                                                       {

  // Definitions
  const int wordSize_32 = 4;
  const int wordSize_64 = 8;

  int numberOf32Words = datasize/wordSize_32;

  const unsigned char* index8 = reinterpret_cast<const unsigned char*>(index32);


  /*  D D U   d a t a */

  // DDU header
  FEDHeader dduHeader(index8);
  if (dduHeader.check()) {
    if(debug) cout << "[DTDDUUnpacker] FED Header. BXID: "<<dduHeader.bxID()
                   << " L1ID: "<<dduHeader.lvl1ID() <<endl;
  } else {
    LogWarning("DTRawToDigi|DTDDUUnpacker") << "[DTDDUUnpacker] WARNING!, this is not a DDU Header, FED ID: "
                                            << dduID << endl;
  }

  // DDU trailer
  // [BITS] stop before FED trailer := 8 bytes
  FEDTrailer dduTrailer(index8 + datasize - 1*wordSize_64); 

  if (dduTrailer.check()) {
    if(debug) cout << "[DTDDUUnpacker] FED Trailer. Lenght of the DT event: "
                   << dduTrailer.lenght() << endl;
  } else {
    LogWarning("DTRawToDigi|DTDDUUnpacker") << "[DTDDUUnpacker] WARNING!, this is not a DDU Trailer, FED ID: "
                                            << dduID << endl;
  }


  // Control DDU data
  DTDDUData controlData(dduHeader,dduTrailer);
  // check the CRC set in the FED trailer (FCRC errors)
  controlData.checkCRCBit(index8 + datasize - 1*wordSize_64);

  // Check Status Words 
  vector<DTDDUFirstStatusWord> rosStatusWords;
  // [BITS] 3 words of 8 bytes + "rosId" bytes
  // In the case we are reading from DMA, the status word are swapped as the ROS data
  if (localDAQ) {
    // DDU channels from 1 to 4
    for (int rosId = 0; rosId < 4; rosId++ ) {
      int wordIndex8 = numberOf32Words*wordSize_32 - 3*wordSize_64 + wordSize_32 + rosId; 
      controlData.addROSStatusWord(DTDDUFirstStatusWord(index8[wordIndex8]));
    }
    // DDU channels from 5 to 8
    for (int rosId = 0; rosId < 4; rosId++ ) {
      int wordIndex8 = numberOf32Words*wordSize_32 - 3*wordSize_64 + rosId; 
      controlData.addROSStatusWord(DTDDUFirstStatusWord(index8[wordIndex8]));
    }
    // DDU channels from 9 to 12
    for (int rosId = 0; rosId < 4; rosId++ ) {
      int wordIndex8 = numberOf32Words*wordSize_32 - 2*wordSize_64 + wordSize_32 + rosId; 
      controlData.addROSStatusWord(DTDDUFirstStatusWord(index8[wordIndex8]));
    }
  }
  else {
    for (int rosId = 0; rosId < 12; rosId++ ) {
      int wordIndex8 = numberOf32Words*wordSize_32 - 3*wordSize_64 + rosId; 
      controlData.addROSStatusWord(DTDDUFirstStatusWord(index8[wordIndex8]));
    }
  }

  int theROSList;
  // [BITS] 2 words of 8 bytes + 4 bytes (half 64 bit word)
  // In the case we are reading from DMA, the status word are swapped as the ROS data
  if (localDAQ) {
    DTDDUSecondStatusWord dduStatusWord(index32[numberOf32Words - 2*wordSize_64/wordSize_32]);
    controlData.addDDUStatusWord(dduStatusWord);
    theROSList =  dduStatusWord.rosList();
  }
  else {
    DTDDUSecondStatusWord dduStatusWord(index32[numberOf32Words - 2*wordSize_64/wordSize_32 + 1]);
    controlData.addDDUStatusWord(dduStatusWord);
    theROSList =  dduStatusWord.rosList();
  }


  /*  R O S   d a t a */

  // Set the index to start looping on ROS data
  // [BITS] one 8 bytes word
  index32 += (wordSize_64)/wordSize_32; 

  // Set the datasize to look only at ROS data 
  // [BITS] header, trailer, 2 status words
  datasize -= 4*wordSize_64; 

  // unpacking the ROS payload
  ros25Unpacker->interpretRawData(index32, datasize, dduID, mapping, detectorProduct, triggerProduct, theROSList);

  // Perform DQM if requested
  if (performDataIntegrityMonitor) 
    dataMonitor->processFED(controlData, ros25Unpacker->getROSsControlData(),dduID);  
  
}

Member Data Documentation

Definition at line 51 of file DTDDUUnpacker.h.

Referenced by DTDDUUnpacker(), and interpretRawData().

Definition at line 39 of file DTDDUUnpacker.h.

Referenced by DTDDUUnpacker().

bool DTDDUUnpacker::debug [private]

Definition at line 47 of file DTDDUUnpacker.h.

Referenced by DTDDUUnpacker(), and interpretRawData().

bool DTDDUUnpacker::localDAQ [private]

if data are read locally, status words are swapped

Definition at line 42 of file DTDDUUnpacker.h.

Referenced by DTDDUUnpacker(), and interpretRawData().

perform DQM for DDU

Definition at line 45 of file DTDDUUnpacker.h.

Referenced by DTDDUUnpacker(), and interpretRawData().

Definition at line 49 of file DTDDUUnpacker.h.

Referenced by DTDDUUnpacker(), interpretRawData(), and ~DTDDUUnpacker().