CMS 3D CMS Logo

LinkDataXMLWriter Class Reference

#include <IORawData/RPCFileReader/interface/LinkDataXMLWriter.h>

Inheritance diagram for LinkDataXMLWriter:

edm::EDAnalyzer

List of all members.

Public Member Functions

virtual void analyze (const edm::Event &, const edm::EventSetup &)
 LinkDataXMLWriter (const edm::ParameterSet &)
virtual ~LinkDataXMLWriter ()

Private Member Functions

void addLinkData (int triggerCrateNum, int triggerBoardNum, int opticalLinkNum, int lbNumber, int partitionNumber, int partitionData, int halfPart, int eod)
void clear ()
std::pair< int, intgetTCandTBNumbers (int dccInputChannelNum, int dccFactor)
std::string IntToString (int i, int opt=0)
void writeLinkData ()

Private Attributes

edm::InputTag digiLabel
DOMDocument * doc
DOMElement * event
std::vector< std::vector
< std::vector< std::vector
< RPCPacData > > > > 
linkData
std::string m_xmlDir
XMLFormatTarget * myFormTarget
int nEvents
int nTB
int nTC
DOMElement * rootElem
DOMWriter * theSerializer

Static Private Attributes

static int m_instanceCount = 0


Detailed Description

Definition at line 35 of file LinkDataXMLWriter.h.


Constructor & Destructor Documentation

LinkDataXMLWriter::LinkDataXMLWriter ( const edm::ParameterSet iConfig  )  [explicit]

Set root element.

Setup output

Definition at line 71 of file LinkDataXMLWriter.cc.

References digiLabel, doc, event, Exception, edm::ParameterSet::getParameter(), IntToString(), linkData, m_instanceCount, m_xmlDir, myFormTarget, nEvents, nTB, nTC, NULL, rootElem, theSerializer, transcode(), and X.

00071                                                                   {
00072 
00073 
00074    digiLabel = iConfig.getParameter<edm::InputTag>("digisSource");
00075 
00076    m_xmlDir = iConfig.getParameter<std::string>("xmlDir");
00077 
00078   if (m_instanceCount == 0) {
00079     try {
00080         XMLPlatformUtils::Initialize();
00081         ++m_instanceCount;
00082     }
00083     catch(const XMLException &toCatch)  {
00084         throw cms::Exception("xmlError") << ("Error during Xerces-c Initialization: "
00085            + std::string(XMLString::transcode(toCatch.getMessage())));
00086     }
00087   }
00088 
00090   DOMImplementation* impl =  DOMImplementationRegistry::getDOMImplementation(X("Core"));
00091 
00092   if (impl != NULL) {
00093 
00094     doc = impl->createDocument(0,                    // root element namespace URI.
00095                                X("rpctDataStream"),         // root element name
00096                                0);                   // document type object (DTD).
00097     
00098     rootElem = doc->getDocumentElement();
00099   }
00100     else {
00101       throw cms::Exception("xmlError") << "Could'n get DOMImplementation\n";
00102     }
00103 
00105   XMLCh tempStr[100];
00106   XMLString::transcode("LS", tempStr, 99);
00107   DOMImplementation *impl_1          = DOMImplementationRegistry::getDOMImplementation(tempStr);
00108   theSerializer = ((DOMImplementationLS*)impl_1)->createDOMWriter();
00109   
00110   // set user specified output encoding
00111   theSerializer->setEncoding(X("UTF-8"));
00112   
00113   if (theSerializer->canSetFeature(XMLUni::fgDOMWRTFormatPrettyPrint, true))
00114     theSerializer->setFeature(XMLUni::fgDOMWRTFormatPrettyPrint, true);
00115   
00116   std::string outFileName = m_xmlDir;
00117   myFormTarget = new LocalFileFormatTarget(X(outFileName.c_str()));
00118   
00119   DOMNode* xmlstylesheet  = doc->createProcessingInstruction(X("xml-stylesheet"),
00120                                                              X("type=\"text/xsl\"href=\"default.xsl\""));
00121 
00122   doc->insertBefore(xmlstylesheet, rootElem);
00124   event = doc->createElement(X("event"));
00125   event->setAttribute(X("bx"), X( IntToString(0).c_str()));
00126   event->setAttribute(X("num"), X( IntToString(0).c_str()));
00127   rootElem->appendChild(event);
00128   
00129 
00130   //set vector sizes
00131   nTC = 12;
00132   nTB = 9;
00133   for(int iTC=0;iTC<nTC;iTC++){
00134     std::vector< RPCPacData> rpdv(18,  RPCPacData());
00135     std::vector<std::vector< RPCPacData> > rpdvv(18,rpdv); 
00136     std::vector<std::vector<std::vector< RPCPacData> > > rpdvvv(nTB,rpdvv); 
00137     linkData.push_back(rpdvvv);
00138   }
00139 
00140   nEvents = 0; 
00141 }

LinkDataXMLWriter::~LinkDataXMLWriter (  )  [virtual]

Definition at line 147 of file LinkDataXMLWriter.cc.

References doc, myFormTarget, and theSerializer.

00147                                      {
00148 
00149   theSerializer->writeNode(myFormTarget, *doc);
00150   doc->release();
00151 
00152   delete theSerializer;
00153   delete myFormTarget;
00154   
00155 
00156 
00157 }


Member Function Documentation

void LinkDataXMLWriter::addLinkData ( int  triggerCrateNum,
int  triggerBoardNum,
int  opticalLinkNum,
int  lbNumber,
int  partitionNumber,
int  partitionData,
int  halfPart,
int  eod 
) [private]

Definition at line 257 of file LinkDataXMLWriter.cc.

References linkData.

Referenced by analyze().

00260                                                           {
00261 
00262   int partDelay = 0;
00263 
00264   RPCPacData myLinkData( partitionData,  partitionNumber,
00265                          partDelay, eod, halfPart, lbNumber);
00266   /*
00267   std::cout<<" opticalLinkNum: "<<opticalLinkNum
00268            <<" lbNumber: "<<lbNumber
00269            <<" rawData: "<<std::hex<<myLinkData.toRaw()<<std::dec<<std::endl;
00270   */
00271   int delay = 2;
00272   for(delay=0;delay<15;delay++) {
00273     /*
00274     std::cout<<"delay: "<<delay<<" raw: "
00275              <<std::hex
00276              <<linkData[triggerCrateNum][triggerBoardNum][delay][opticalLinkNum-1].toRaw()
00277              <<std::dec<<std::endl;
00278 */
00279  
00280     if(!linkData[triggerCrateNum][triggerBoardNum][delay][opticalLinkNum].toRaw()) break;
00281   }
00282 
00283   linkData[triggerCrateNum][triggerBoardNum][delay][opticalLinkNum] = myLinkData;
00284 }

void LinkDataXMLWriter::analyze ( const edm::Event ev,
const edm::EventSetup es 
) [virtual]

Get Data from all FEDs

Implements edm::EDAnalyzer.

Definition at line 162 of file LinkDataXMLWriter.cc.

References a, addLinkData(), RPCEMap::convert(), GenMuonPlsPt100GeV_cfg::cout, detId, reco_application_tbsim_simpleTBanalysis_cfg::digiCollection, digiLabel, lat::endl(), RPCPackingModule::eventRecords(), edm::EventSetup::get(), edm::Event::getByLabel(), FEDNumbering::getRPCFEDIds(), getTCandTBNumbers(), id, edm::ESHandle< T >::product(), edm::Handle< T >::product(), range, RPCEMap::theVersion, RPCReadOutMapping::version(), and writeLinkData().

00162                                                                           {
00163   static RPCReadOutMapping * cabling =  new RPCReadOutMapping("");
00164 
00166   Handle< RPCDigiCollection > digiCollection;
00167   ev.getByLabel(digiLabel, digiCollection);
00169   RPCDigiCollection::DigiRangeIterator rpcDigiCI;
00170   for(rpcDigiCI = digiCollection->begin();rpcDigiCI!=digiCollection->end();rpcDigiCI++){
00171     cout<<(*rpcDigiCI).first<<endl;
00172     RPCDetId detId=(*rpcDigiCI).first;
00173     uint32_t id=detId();
00174 
00175     const RPCDigiCollection::Range& range = (*rpcDigiCI).second;
00176     for (RPCDigiCollection::const_iterator digiIt = range.first;
00177          digiIt!=range.second;++digiIt) cout<<"Digi: "<<*digiIt<<endl;
00178   }
00179  
00180  
00181 //  ESHandle<RPCReadOutMapping> readoutMapping;
00182 //  es.get<RPCReadOutMappingRcd>().get(readoutMapping);
00183 
00184   ESHandle<RPCEMap> readoutMapping;
00185   es.get<RPCEMapRcd>().get(readoutMapping);
00186   const RPCEMap* eMap=readoutMapping.product();
00187 
00188   if (eMap->theVersion != cabling->version()) {
00189     delete cabling;
00190     cabling = eMap->convert();
00191   }
00192 
00193   int trigger_BX = 200;
00194   int dccFactor = 3;
00195   
00196    pair<int,int> rpcFEDS=FEDNumbering::getRPCFEDIds();
00197    for (int id= rpcFEDS.first; id<=rpcFEDS.second; ++id){
00198      dccFactor--;
00199 
00200 //    RPCRecordFormatter formatter(id, readoutMapping.product()) ;
00201     RPCRecordFormatter formatter(id, cabling) ;
00202     std::vector<rpcrawtodigi::EventRecords> myEventRecords = RPCPackingModule::eventRecords(id,
00203                                                                                             trigger_BX,  
00204                                                                                             digiCollection.product(),
00205                                                                                             formatter); 
00206     std::cout<<" FED id: "<< id;
00207     std::cout<<" myEventRecords.size(): "<< myEventRecords.size()<<std::endl;
00208     
00209     std::vector<rpcrawtodigi::EventRecords>::const_iterator CI =  myEventRecords.begin();
00210     for(;CI!= myEventRecords.end();CI++){ 
00211 
00212       int dccInputChannelNum =  CI->recordSLD().rmb();
00213       int opticalLinkNum =   CI->recordSLD().tbLinkInputNumber();
00214       int partitionData =   CI->recordCD().partitionData(); 
00215       int halfP = CI->recordCD().halfP(); 
00216       int eod = CI->recordCD().eod(); 
00217       int partitionNumber = CI->recordCD().partitionNumber();  
00218       int lbNumber = CI->recordCD().lbInLink();
00219 
00220       std::pair<int,int> aPair = getTCandTBNumbers(dccInputChannelNum,dccFactor);
00221       int triggerCrateNum = aPair.first;
00222       int triggerBoardNum = aPair.second;
00223       int a = aPair.second;
00224     
00225       addLinkData(triggerCrateNum, triggerBoardNum, opticalLinkNum, 
00226                   lbNumber, partitionNumber,  partitionData, halfP, eod);
00227 
00228     }
00229    }
00230     
00231 
00232    /*
00233     int triggerCrateNum = 9;
00234     int halfP = 0;
00235     int eod = 0;
00236     for(int iTB=0;iTB<9;iTB++){
00237       int iOL = nEvents%18;
00238       //for(int iOL=0;iOL<18;iOL++){
00239       int iPartNum = nEvents%12;
00240       //for(int iPartNum=0;iPartNum<8;iPartNum++){          
00241         for(int iLbNum=0;iLbNum<3;iLbNum++){
00242           int partitionData = iPartNum*2+1;  
00243           addLinkData(triggerCrateNum, iTB, iOL, iLbNum, 
00244                       iPartNum, partitionData, halfP, eod);       
00245         }
00246         //}
00247         //}
00248     }
00249    */
00250     writeLinkData();  
00251     
00252 }

void LinkDataXMLWriter::clear ( void   )  [private]

Definition at line 372 of file LinkDataXMLWriter.cc.

References linkData, nTB, and nTC.

Referenced by writeLinkData().

00372                               {
00373 
00374   int nTC = 12;
00375   int nTB = 9;
00376 
00377   for(int triggerCrateNum=0;triggerCrateNum<nTC;triggerCrateNum++){
00378     for(int triggerBoardNum=0;triggerBoardNum<nTB;triggerBoardNum++){
00379       for(int opticalLinkNum=0;opticalLinkNum<18;opticalLinkNum++){
00380         for(int iDelay=0;iDelay<18;iDelay++){
00381           linkData[triggerCrateNum][triggerBoardNum][iDelay][opticalLinkNum] = RPCPacData();
00382         }
00383       }
00384     }
00385   }
00386 
00387 }

std::pair< int, int > LinkDataXMLWriter::getTCandTBNumbers ( int  dccInputChannelNum,
int  dccFactor 
) [private]

cout<<"dcc: "<<dccInputChannelNum<<endl;

Definition at line 390 of file LinkDataXMLWriter.cc.

References i.

Referenced by analyze().

00390                                                                                           {
00391 
00392   int tcNumber = 0;
00393   int tbNumber = 0;
00394 
00395   for(int i=0;i<9;i++){
00396     if(dccInputChannelNum==i ||
00397        dccInputChannelNum==9+i ||
00398        dccInputChannelNum==18+i ||
00399        dccInputChannelNum==27+i) tbNumber = i;
00400   }
00404  for(int i=0;i<4;i++){
00405     if(dccInputChannelNum>=i*9 &&
00406        dccInputChannelNum<9+i*9) tcNumber = i; //Count TC from 0, not from 1.
00407 
00408   }
00409 
00410  tcNumber+=4*dccFactor;
00411 
00412   return std::pair<int,int>(tcNumber,tbNumber);
00413 
00414 }

std::string LinkDataXMLWriter::IntToString ( int  i,
int  opt = 0 
) [private]

Definition at line 362 of file LinkDataXMLWriter.cc.

References ss.

Referenced by LinkDataXMLWriter(), and writeLinkData().

00362                                                       {
00363 
00364  std::stringstream ss;
00365  if(opt==1) ss << std::hex << i << std::dec;
00366  else ss << i ;
00367 
00368  return ss.str();
00369 
00370 }

void LinkDataXMLWriter::writeLinkData (  )  [private]

Definition at line 289 of file LinkDataXMLWriter.cc.

References clear(), doc, RPCPacData::endOfData(), RPCPacData::halfPartition(), IntToString(), RPCPacData::lbNum(), linkData, nEvents, nTB, nTC, NULL, RPCPacData::partitionData(), RPCPacData::partitionDelay(), RPCPacData::partitionNum(), tc, RPCPacData::toRaw(), and X.

Referenced by analyze().

00289                                      {
00290 
00291   int bxNum = nEvents*10; //Leave room for data multiplexing.
00292   nEvents++;
00293 
00294   time_t timer;
00295   struct tm *tblock;
00296   timer = time(NULL);
00297   tblock = localtime(&timer);
00298 
00299   DOMElement* bx = doc->createElement(X("bxData"));
00300   bx->setAttribute(X("num"), X( IntToString(bxNum).c_str()));
00301   event->appendChild(bx);
00302 
00303   DOMElement*  tc = 0;
00304   DOMElement*  tb = 0;
00305   DOMElement*  ol = 0;
00306 
00307   for(int triggerCrateNum=0;triggerCrateNum<nTC;triggerCrateNum++){
00309     bool nonEmpty = false;
00310     for(int triggerBoardNum=0;triggerBoardNum<nTB;triggerBoardNum++){
00311       for(int opticalLinkNum=0;opticalLinkNum<18;opticalLinkNum++){
00312         for(int iDelay=0;iDelay<18;iDelay++){
00313           int rawData =  linkData[triggerCrateNum][triggerBoardNum][iDelay][opticalLinkNum].toRaw();
00314           if(rawData) nonEmpty = true;
00315         }
00316       }
00317     }
00318     if(!nonEmpty) continue;
00320     tc = doc->createElement(X("tc"));
00321     tc->setAttribute(X("num"), X( IntToString(triggerCrateNum).c_str()));
00322     for(int triggerBoardNum=0;triggerBoardNum<nTB;triggerBoardNum++){
00324       bool nonEmpty = false;
00325       for(int opticalLinkNum=0;opticalLinkNum<18;opticalLinkNum++){
00326         for(int iDelay=0;iDelay<18;iDelay++){
00327           int rawData =  linkData[triggerCrateNum][triggerBoardNum][iDelay][opticalLinkNum].toRaw();
00328           if(rawData) nonEmpty = true;
00329         }
00330       }
00331       if(!nonEmpty) continue;
00333       tb = doc->createElement(X("tb"));
00334       tb->setAttribute(X("num"), X( IntToString(triggerBoardNum).c_str()));
00335       for(int opticalLinkNum=0;opticalLinkNum<18;opticalLinkNum++){
00336         if(!linkData[triggerCrateNum][triggerBoardNum][0][opticalLinkNum].toRaw()) continue;
00337         ol = doc->createElement(X("ol"));
00338         ol->setAttribute(X("num"), X( IntToString(opticalLinkNum).c_str()));
00339         for(int iDelay=0;iDelay<18;iDelay++){
00340           int rawData =  linkData[triggerCrateNum][triggerBoardNum][iDelay][opticalLinkNum].toRaw();
00341           if(!rawData) continue;
00342           const RPCPacData & myLinkData = linkData[triggerCrateNum][triggerBoardNum][iDelay][opticalLinkNum];
00343           DOMElement*  lmd = doc->createElement(X("lmd"));
00344           lmd->setAttribute(X("lb"), X( IntToString(myLinkData.lbNum()).c_str()));
00345           lmd->setAttribute(X("par"), X( IntToString(myLinkData.partitionNum()).c_str()));
00346           lmd->setAttribute(X("dat"), X( IntToString(myLinkData.partitionData(),1).c_str()));
00347           lmd->setAttribute(X("del"), X( IntToString(myLinkData.partitionDelay()).c_str()));
00348           lmd->setAttribute(X("eod"), X( IntToString(myLinkData.endOfData()).c_str()));
00349           lmd->setAttribute(X("hp"), X( IntToString(myLinkData.halfPartition()).c_str()));
00350           lmd->setAttribute(X("raw"), X( IntToString(myLinkData.toRaw(),1).c_str()));
00351           ol->appendChild(lmd);
00352         }
00353         tb->appendChild(ol);
00354       }
00355       tc->appendChild(tb);      
00356     }
00357     bx->appendChild(tc);
00358   }
00359   clear();
00360 }


Member Data Documentation

edm::InputTag LinkDataXMLWriter::digiLabel [private]

Definition at line 60 of file LinkDataXMLWriter.h.

Referenced by analyze(), and LinkDataXMLWriter().

DOMDocument* LinkDataXMLWriter::doc [private]

Definition at line 70 of file LinkDataXMLWriter.h.

Referenced by LinkDataXMLWriter(), writeLinkData(), and ~LinkDataXMLWriter().

DOMElement* LinkDataXMLWriter::event [private]

Definition at line 72 of file LinkDataXMLWriter.h.

Referenced by LinkDataXMLWriter().

std::vector<std::vector<std::vector<std::vector< RPCPacData> > > > LinkDataXMLWriter::linkData [private]

Definition at line 75 of file LinkDataXMLWriter.h.

Referenced by addLinkData(), clear(), LinkDataXMLWriter(), and writeLinkData().

int LinkDataXMLWriter::m_instanceCount = 0 [static, private]

Definition at line 62 of file LinkDataXMLWriter.h.

Referenced by LinkDataXMLWriter().

std::string LinkDataXMLWriter::m_xmlDir [private]

Definition at line 59 of file LinkDataXMLWriter.h.

Referenced by LinkDataXMLWriter().

XMLFormatTarget* LinkDataXMLWriter::myFormTarget [private]

Definition at line 69 of file LinkDataXMLWriter.h.

Referenced by LinkDataXMLWriter(), and ~LinkDataXMLWriter().

int LinkDataXMLWriter::nEvents [private]

Definition at line 66 of file LinkDataXMLWriter.h.

Referenced by LinkDataXMLWriter(), and writeLinkData().

int LinkDataXMLWriter::nTB [private]

Definition at line 67 of file LinkDataXMLWriter.h.

Referenced by clear(), LinkDataXMLWriter(), and writeLinkData().

int LinkDataXMLWriter::nTC [private]

Definition at line 67 of file LinkDataXMLWriter.h.

Referenced by clear(), LinkDataXMLWriter(), and writeLinkData().

DOMElement* LinkDataXMLWriter::rootElem [private]

Definition at line 71 of file LinkDataXMLWriter.h.

Referenced by LinkDataXMLWriter().

DOMWriter* LinkDataXMLWriter::theSerializer [private]

Definition at line 68 of file LinkDataXMLWriter.h.

Referenced by LinkDataXMLWriter(), and ~LinkDataXMLWriter().


The documentation for this class was generated from the following files:
Generated on Tue Jun 9 18:27:45 2009 for CMSSW by  doxygen 1.5.4