FairRoot/PandaRoot
Public Member Functions | Static Public Member Functions | Protected Member Functions | Protected Attributes | Private Attributes | List of all members
PndMQMvdPixelDigiProcessorBursts Class Reference

#include <PndMQMvdPixelDigiProcessorBursts.h>

Inheritance diagram for PndMQMvdPixelDigiProcessorBursts:
PndMQBurstProcessor

Public Member Functions

 PndMQMvdPixelDigiProcessorBursts ()
 
virtual ~PndMQMvdPixelDigiProcessorBursts ()
 
virtual void SetParameters ()
 
virtual void ProcessData ()
 
virtual void UpdateParameters ()
 
virtual FairParGenericSet * UpdateParameter (FairParGenericSet *thisPar)
 
template<class Archive >
void serialize (Archive &ar, const unsigned int version)
 

Static Public Member Functions

static void CustomCleanupParameters (void *data, void *hint)
 
static void free_string (void *data, void *hint)
 

Protected Member Functions

virtual void Run ()
 

Protected Attributes

BurstData fBurstDataIn
 
BurstData fBurstDataOut
 
int fCurrentRunId
 
int fNewRunId
 
TList * fParCList
 

Private Attributes

FairGeoParSet * fGeoPar
 
PndSdsPixelDigiParfDigiPar
 
PndSdsTotDigiParfTotPar
 
PndSensorNameParfSensorPar
 
PndMQGapEventBuilder fEventBuilder
 
PndMapSorter fMapSorter
 
PndMvdPixelClusterFinder fClusterFinder
 
PndMQMvdChargeWeightedPixelMappingfPixelMapping
 
PndGeoHandlingfGeoHandler
 

Detailed Description

Definition at line 59 of file PndMQMvdPixelDigiProcessorBursts.h.

Constructor & Destructor Documentation

PndMQMvdPixelDigiProcessorBursts::PndMQMvdPixelDigiProcessorBursts ( )
inline

Definition at line 62 of file PndMQMvdPixelDigiProcessorBursts.h.

63  {
64  fGeoPar = new FairGeoParSet("FairGeoParSet");
65  fParCList = new TList();
66  fParCList->Add(fGeoPar);
67 
68  fDigiPar = new PndSdsPixelDigiPar("MVDPixelDigiPar");
69  fParCList->Add(fDigiPar);
70 
71  fTotPar = new PndSdsTotDigiPar("MVDPixelTotDigiPar");
72  fParCList->Add(fTotPar);
73 
74  fSensorPar = new PndSensorNamePar("PndSensorNamePar");
75  fParCList->Add(fSensorPar);
76  }
PndMQMvdChargeWeightedPixelMapping * fPixelMapping
Charge Digitization Parameter Class for SDS.
Unique match between SensorID and path in TGeoManager.
Digitization Parameter Class for SDS-Pixel part.
virtual PndMQMvdPixelDigiProcessorBursts::~PndMQMvdPixelDigiProcessorBursts ( )
inlinevirtual

Definition at line 78 of file PndMQMvdPixelDigiProcessorBursts.h.

79  {
80  if (fGeoPar != 0)
81  delete fGeoPar;
82  if (fDigiPar != 0)
83  delete fDigiPar;
84  if (fTotPar != 0)
85  delete fTotPar;
86  if (fSensorPar != 0)
87  delete fSensorPar;
88  if (fPixelMapping != 0)
89  delete fPixelMapping;
90  if (fGeoHandler != 0)
91  delete fGeoHandler;
92  }
PndMQMvdChargeWeightedPixelMapping * fPixelMapping

Member Function Documentation

void PndMQBurstProcessor::CustomCleanupParameters ( void *  data,
void *  hint 
)
staticinherited

Definition at line 26 of file PndMQBurstProcessor.cxx.

Referenced by PndMQBurstProcessor::UpdateParameter().

27 {
28  delete (std::string*)hint;
29 }
void PndMQBurstProcessor::free_string ( void *  data,
void *  hint 
)
staticinherited

Definition at line 31 of file PndMQBurstProcessor.cxx.

Referenced by PndMQBurstProcessor::Run().

32 {
33  delete static_cast<std::string*>(hint);
34 }
void PndMQMvdPixelDigiProcessorBursts::ProcessData ( )
virtual

Implements PndMQBurstProcessor.

Definition at line 38 of file PndMQMvdPixelDigiProcessorBursts.cxx.

References PndMapSorter::AddElements(), BurstHeader::fBranchName, PndMQBurstProcessor::fBurstDataIn, PndMQBurstProcessor::fBurstDataOut, fClusterFinder, BurstData::fData, fEventBuilder, BurstData::fHeader, PndMQGapEventBuilder::FillData(), fMapSorter, fPixelMapping, PndMQSdsChargeWeightedPixelMapping::GetCluster(), PndMvdPixelClusterFinder::GetClusters(), PndMQGapEventBuilder::GetLastData(), PndMapSorter::GetOutputData(), PndMQGapEventBuilder::GetSeparatedData(), hit, hits, pos, BurstData::Reset(), and PndMapSorter::WriteOutAll().

39 {
40 
42  for (auto dataItr : fBurstDataIn.fData)
43  fMapSorter.AddElements(dataItr);
45  std::vector<FairTimeStamp*> sortedData = fMapSorter.GetOutputData();
46 
47 // for (auto dataItr: sortedData)
48 // LOG(INFO) << "DataAfterSorter " << dataItr->GetTimeStamp() ;
49 
50 
51  fEventBuilder.FillData(sortedData);
52  std::vector<std::vector<FairTimeStamp*> > gapData = fEventBuilder.GetSeparatedData();
53  std::vector<std::vector<FairTimeStamp*> > lastData = fEventBuilder.GetLastData();
54  gapData.insert(gapData.end(), lastData.begin(), lastData.end());
55 // for (auto eventItr : gapData)
56 // for (auto dataItr : eventItr)
57 // LOG(INFO) << "DataAfterEventBuilder: " << dataItr->GetTimeStamp();
58  if (fBurstDataIn.fHeader.fBranchName == "MVDPixelDigis"){
59  std::vector<std::vector<FairTimeStamp*> > hitsInGaps;
60  for (auto gapIter : gapData){
61  std::vector<PndSdsDigiPixel*> digiPixel;
62  std::vector<FairTimeStamp*> hits;
63  if (gapIter.size() > 0){
64  // LOG(INFO) << "GapSize > 0 " << gapIter.size();
65  for (auto dataIter : gapIter){ //iterate hits in gap
66  digiPixel.push_back((PndSdsDigiPixel*) dataIter);
67 // LOG(INFO) << "DataInGap: " << *(PndSdsDigiPixel*) dataIter;
68  }
69  std::vector<std::vector< Int_t> > clusterInts = fClusterFinder.GetClusters(digiPixel);
70 
71  for (auto clusterIter : clusterInts) //iterate all found cluster
72  {
73  //LOG(INFO) << "ClusterSize: " << clusterIter.size();
74  std::vector<PndSdsDigiPixel*> digiInCluster;
75  for (auto digiInClusterItr : clusterIter){ //iterate digisInt in one cluster
76  int pos = digiInClusterItr;
77  digiInCluster.push_back(digiPixel[pos]);
78  }
79 // LOG(INFO) << "DigisInCluster: ";
80  if (digiInCluster.size() > 0){
81 // for (auto digiIter : digiInCluster)
82 // LOG(INFO) << digiIter->GetTimeStamp();
83  PndSdsHit* hit = new PndSdsHit(fPixelMapping->GetCluster(digiInCluster));
84  hits.push_back(hit);
85 // LOG(INFO) << "Hit: " << hit->GetTimeStamp();
86  digiInCluster.clear();
87  // LOG(INFO) << "Hit: " << *(PndSdsHit*)(*hits.rbegin());
88  }
89  }
90  hitsInGaps.push_back(hits);
91  }
92  }
93  fBurstDataOut.fHeader.fBranchName = "MVDHitsPixel";
94  fBurstDataOut.fData = hitsInGaps;
95  }
96 }
TVector3 pos
PndMQMvdChargeWeightedPixelMapping * fPixelMapping
virtual void Reset()
PndSdsHit GetCluster(std::vector< PndSdsDigiPixel * > &pixelArray)
Main function of class to calculate the PndSdsHit out of the given PndSdsDigis.
std::vector< std::vector< PndSdsDigiTopix4 > > GetLastData()
std::vector< std::vector< Int_t > > GetClusters(std::vector< PndSdsDigiPixel * > &hits)
std::vector< std::vector< PndSdsDigiTopix4 > > GetSeparatedData()
virtual void AddElements(std::vector< FairTimeStamp * > dataArray)
std::vector< std::vector< FairTimeStamp * > > fData
virtual std::vector< FairTimeStamp * > GetOutputData()
Data class to store the digi output of a pixel module.
PndSdsMCPoint * hit
Definition: anasim.C:70
CbmHit * hits[nHits]
Definition: RiemannTest.C:19
void PndMQBurstProcessor::Run ( )
protectedvirtualinherited

Definition at line 36 of file PndMQBurstProcessor.cxx.

References PndMQBurstProcessor::fBurstDataIn, PndMQBurstProcessor::fBurstDataOut, BurstHeader::fBurstID, PndMQBurstProcessor::fCurrentRunId, BurstData::fData, PndMQBurstProcessor::fHasBoostSerialization, BurstData::fHeader, PndMQBurstProcessor::fNewRunId, PndMQBurstProcessor::free_string(), BurstHeader::fRunID, PndMQBurstProcessor::ProcessData(), PndMQBurstProcessor::SetParameters(), and PndMQBurstProcessor::UpdateParameters().

37 {
39  {
40  while (CheckCurrentState(RUNNING))
41  {
42  FairMQParts parts;
43  std::unique_ptr<FairMQMessage> msg(fTransportFactory->CreateMessage());
44  if (Receive(msg, "data-in") >= 0)
45  {
46  LOG(INFO) << "Received data! " << msg->GetSize() << std::endl;
47  string msgStr(static_cast<char*>(msg->GetData()), msg->GetSize());
48  istringstream ibuffer(msgStr);
49  if(!ibuffer.good())
50  LOG(INFO) << "IBUFFER IS BAD!";
51  try {
52  boost::archive::binary_iarchive InputArchive(ibuffer);
53 
54  InputArchive >> fBurstDataIn;
55  }
56  catch (boost::archive::archive_exception& e)
57  {
58  LOG(ERROR) << e.what();
59 // continue;
60  }
62  if (fNewRunId != fCurrentRunId) {
65  SetParameters();
66  }
67 
68  if (fBurstDataIn.fData.size() < 1) continue;
69 // LOG(INFO) << "BurstData size: " << fBurstDataIn.fData[0].size() << " BurstID: " << fBurstDataIn.fHeader.fBurstID;
70 
71  ProcessData();
72  }
73  if (fBurstDataOut.fData.size() > 0){
75  std::ostringstream obuffer;
76  boost::archive::binary_oarchive OutputArchive(obuffer);
77  OutputArchive << fBurstDataOut;
78  std::string* strMsg = new std::string(obuffer.str());
79  unique_ptr<FairMQMessage> msgOut(NewMessage(const_cast<char*>(strMsg->c_str()), strMsg->length(), free_string, strMsg));
80  LOG(INFO) << "Data sent: " << fBurstDataOut.fHeader.fBranchName << " BurstID: " << fBurstDataOut.fHeader.fBurstID << " size: " << msgOut->GetSize();
81  int event = 0;
82 // for (auto eventItr : fBurstDataOut->fData){
83 // for (auto dataItr : eventItr){
84 // LOG(INFO) << event << " : " << dataItr->GetTimeStamp();
85 // }
86 // event++;
87 // }
88  Send(msgOut, "data-out");
89  }
90  }
91  }
92  else
93  {
94  LOG(ERROR) << " Boost Serialization not ok";
95  }
96 }
static void free_string(void *data, void *hint)
virtual void UpdateParameters()
virtual void SetParameters()
virtual void ProcessData()=0
std::vector< std::vector< FairTimeStamp * > > fData
template<class Archive >
void PndMQBurstProcessor::serialize ( Archive &  ar,
const unsigned int  version 
)
inlineinherited

Definition at line 89 of file PndMQBurstProcessor.h.

90  {
91  ar& fBurstDataIn;
92  ar& fBurstDataOut;
93  }
void PndMQMvdPixelDigiProcessorBursts::SetParameters ( )
virtual

Reimplemented from PndMQBurstProcessor.

Definition at line 13 of file PndMQMvdPixelDigiProcessorBursts.cxx.

References fDigiPar, fGeoHandler, fGeoPar, PndSensorNamePar::FillMap(), PndMQBurstProcessor::fParCList, fPixelMapping, fSensorPar, fTotPar, PndMQMvdChargeWeightedPixelMapping::SetDigiPar(), PndMQSdsPixelBackMapping::SetGeoHandling(), PndGeoHandling::SetSensorNamePar(), and PndMQMvdChargeWeightedPixelMapping::SetTotPar().

14 {
15  fGeoPar->GetGeometry();
16  fDigiPar = (PndSdsPixelDigiPar*)fParCList->FindObject("MVDPixelDigiPar");
17  fTotPar = (PndSdsTotDigiPar*)fParCList->FindObject("MVDPixelTotDigiPar");
18  fSensorPar = (PndSensorNamePar*)fParCList->FindObject("PndSensorNamePar");
20  if (fGeoHandler == 0){
22  } else {
24  }
25 
26  LOG(INFO) << "Conversion Method: " << *fDigiPar;
27  LOG(INFO) << "Tot Info: " << *fTotPar;
28  if (fPixelMapping == 0){
30  } else {
34  }
35 }
void SetSensorNamePar(PndSensorNamePar *par)
virtual void SetGeoHandling(PndGeoHandling *geo)
PndMQMvdChargeWeightedPixelMapping * fPixelMapping
Charge Digitization Parameter Class for SDS.
virtual void SetTotPar(PndSdsTotDigiPar *par)
Class to access the naming information of the MVD.
PndSdsChargedWeightedPixelMapping: Gets a vector of DigiHits and calculates the cluster center weight...
virtual void SetDigiPar(PndSdsPixelDigiPar *par)
Unique match between SensorID and path in TGeoManager.
Digitization Parameter Class for SDS-Pixel part.
FairParGenericSet * PndMQBurstProcessor::UpdateParameter ( FairParGenericSet *  thisPar)
virtualinherited

Definition at line 106 of file PndMQBurstProcessor.cxx.

References PndMQBurstProcessor::CustomCleanupParameters(), and PndMQBurstProcessor::fCurrentRunId.

Referenced by PndMQBurstProcessor::UpdateParameters().

106  {
107  std::string paramName = thisPar->GetName();
108  // boost::this_thread::sleep(boost::posix_time::milliseconds(1000));
109  std::string* reqStr = new std::string(paramName + "," + std::to_string(fCurrentRunId));
110  LOG(WARN) << "Requesting parameter \"" << paramName << "\" for Run ID " << fCurrentRunId << " (" << thisPar << ")";
111  std::unique_ptr<FairMQMessage> req(NewMessage(const_cast<char*>(reqStr->c_str()), reqStr->length(), CustomCleanupParameters, reqStr));
112  std::unique_ptr<FairMQMessage> rep(NewMessage());
113 
114  if (Send(req,"param") > 0)
115  {
116  if (Receive(rep,"param") > 0)
117  {
118  TMessage2 tm(rep->GetData(), rep->GetSize());
119  thisPar = (FairParGenericSet*)tm.ReadObject(tm.GetClass());
120  LOG(WARN) << "Received parameter"<< paramName <<" from the server (" << thisPar << ")" << tm.GetClass()->GetName() << " DataSize: " << rep->GetSize();
121  thisPar->print();
122  return thisPar;
123  }
124  }
125  return NULL;
126 }
static void CustomCleanupParameters(void *data, void *hint)
void PndMQBurstProcessor::UpdateParameters ( )
virtualinherited

Definition at line 98 of file PndMQBurstProcessor.cxx.

References PndMQBurstProcessor::fParCList, and PndMQBurstProcessor::UpdateParameter().

Referenced by PndMQBurstProcessor::Run().

98  {
99  for ( int iparC = 0 ; iparC < fParCList->GetEntries() ; iparC++ ) {
100  FairParGenericSet* tempObj = (FairParGenericSet*)(fParCList->At(iparC));
101  fParCList->Remove(tempObj);
102  fParCList->AddAt(UpdateParameter(tempObj),iparC);
103  }
104 }
virtual FairParGenericSet * UpdateParameter(FairParGenericSet *thisPar)

Member Data Documentation

BurstData PndMQBurstProcessor::fBurstDataIn
protectedinherited

Definition at line 97 of file PndMQBurstProcessor.h.

Referenced by ProcessData(), and PndMQBurstProcessor::Run().

BurstData PndMQBurstProcessor::fBurstDataOut
protectedinherited

Definition at line 98 of file PndMQBurstProcessor.h.

Referenced by ProcessData(), and PndMQBurstProcessor::Run().

PndMvdPixelClusterFinder PndMQMvdPixelDigiProcessorBursts::fClusterFinder
private

Definition at line 108 of file PndMQMvdPixelDigiProcessorBursts.h.

Referenced by ProcessData().

int PndMQBurstProcessor::fCurrentRunId
protectedinherited
PndSdsPixelDigiPar* PndMQMvdPixelDigiProcessorBursts::fDigiPar
private

Definition at line 102 of file PndMQMvdPixelDigiProcessorBursts.h.

Referenced by SetParameters().

PndMQGapEventBuilder PndMQMvdPixelDigiProcessorBursts::fEventBuilder
private

Definition at line 106 of file PndMQMvdPixelDigiProcessorBursts.h.

Referenced by ProcessData().

PndGeoHandling* PndMQMvdPixelDigiProcessorBursts::fGeoHandler
private

Definition at line 111 of file PndMQMvdPixelDigiProcessorBursts.h.

Referenced by SetParameters().

FairGeoParSet* PndMQMvdPixelDigiProcessorBursts::fGeoPar
private

Definition at line 101 of file PndMQMvdPixelDigiProcessorBursts.h.

Referenced by SetParameters().

PndMapSorter PndMQMvdPixelDigiProcessorBursts::fMapSorter
private

Definition at line 107 of file PndMQMvdPixelDigiProcessorBursts.h.

Referenced by ProcessData().

int PndMQBurstProcessor::fNewRunId
protectedinherited

Definition at line 100 of file PndMQBurstProcessor.h.

Referenced by PndMQBurstProcessor::Run().

TList* PndMQBurstProcessor::fParCList
protectedinherited

Definition at line 101 of file PndMQBurstProcessor.h.

Referenced by SetParameters(), and PndMQBurstProcessor::UpdateParameters().

PndMQMvdChargeWeightedPixelMapping* PndMQMvdPixelDigiProcessorBursts::fPixelMapping
private

Definition at line 110 of file PndMQMvdPixelDigiProcessorBursts.h.

Referenced by ProcessData(), and SetParameters().

PndSensorNamePar* PndMQMvdPixelDigiProcessorBursts::fSensorPar
private

Definition at line 104 of file PndMQMvdPixelDigiProcessorBursts.h.

Referenced by SetParameters().

PndSdsTotDigiPar* PndMQMvdPixelDigiProcessorBursts::fTotPar
private

Definition at line 103 of file PndMQMvdPixelDigiProcessorBursts.h.

Referenced by SetParameters().


The documentation for this class was generated from the following files: