FairRoot/PandaRoot
Public Types | Public Member Functions | Protected Member Functions | Private Attributes | Friends | List of all members
PndMQFileSamplerHits Class Reference

#include <PndMQFileSamplerHits.h>

Inheritance diagram for PndMQFileSamplerHits:

Public Types

enum  { InputFile = FairMQDevice::Last, Last }
 

Public Member Functions

 PndMQFileSamplerHits (std::string inputFileName="")
 
virtual ~PndMQFileSamplerHits ()
 
virtual void SetFileName (std::string fileName)
 
virtual void InitInputFile ()
 
template<class Archive >
void serialize (Archive &ar, const unsigned int version)
 
void SetProperty (const int key, const std::string &value)
 
std::string GetProperty (const int key, const std::string &default_)
 

Protected Member Functions

virtual void Run ()
 

Private Attributes

TFile * fInFile
 
TTree * fTree
 
TClonesArray * fInput
 
TList * fBranchNameList
 
std::string fInputFileName
 
vector< vector< PndSdsHit > > fHitVector
 
bool fHasBoostSerialization
 

Friends

class boost::serialization::access
 

Detailed Description

Definition at line 54 of file PndMQFileSamplerHits.h.

Member Enumeration Documentation

anonymous enum
Enumerator
InputFile 
Last 

Definition at line 59 of file PndMQFileSamplerHits.h.

Constructor & Destructor Documentation

PndMQFileSamplerHits::PndMQFileSamplerHits ( std::string  inputFileName = "")
inline

Definition at line 64 of file PndMQFileSamplerHits.h.

65  : fInFile(NULL)
66  , fTree(NULL)
67  , fInput(NULL)
68  , fHitVector()
69  , fHasBoostSerialization(false)
70  , fInputFileName(inputFileName)
71  , fBranchNameList(0)
72  {
73  gSystem->ResetSignal(kSigInterrupt);
74  gSystem->ResetSignal(kSigTermination);
75 
76  // Check if boost serialization is available if it is chosen
77  using namespace baseMQ::tools::resolve;
78  // coverity[pointless_expression]: suppress coverity warnings on apparant if(const).
79  if (is_same<boost::archive::binary_oarchive, boost::archive::binary_oarchive>::value)
80  {
81  if (has_BoostSerialization<PndSdsHit, void(boost::archive::binary_oarchive&, const unsigned int)>::value == 1)
82  {
84  }
85  }
86  }
vector< vector< PndSdsHit > > fHitVector
virtual PndMQFileSamplerHits::~PndMQFileSamplerHits ( )
inlinevirtual

Definition at line 88 of file PndMQFileSamplerHits.h.

89  {
90  //fBranchNameList->Write("BranchList", TObject::kSingleKey);
91  //fTree->Write();
92  fInFile->Close();
93  if (fHitVector.size() > 0)
94  {
95  fHitVector.clear();
96  }
97 
98  }
vector< vector< PndSdsHit > > fHitVector

Member Function Documentation

std::string PndMQFileSamplerHits::GetProperty ( const int  key,
const std::string &  default_ 
)
inline

Definition at line 132 of file PndMQFileSamplerHits.h.

133  {
134  switch (key)
135  {
136  case InputFile:
137  return fInputFileName;
138  default:
139  return FairMQDevice::GetProperty(key, default_);
140  }
141  }
virtual void PndMQFileSamplerHits::InitInputFile ( )
inlinevirtual

Definition at line 104 of file PndMQFileSamplerHits.h.

Referenced by Run().

105  {
106  fInput = new TClonesArray("PndSdsHit");
107 
108  fInFile = new TFile(fInputFileName.c_str(), "read");
109  fTree = (TTree*)fInFile->Get("pndsim");
110  fTree->SetBranchAddress("Output", &fInput);
111  }
void PndMQFileSamplerHits::Run ( )
protectedvirtual

Definition at line 13 of file PndMQFileSamplerHits.cxx.

References fHasBoostSerialization, fHitVector, fInput, fTree, hit(), i, InitInputFile(), nEvents, PndMQStatus::RUNNING, status, and PndMQStatus::STOP.

14 {
16  {
17  InitInputFile();
18  int sendMsgs = 0;
19 
20  // store the channel references to avoid traversing the map on every loop iteration
21  FairMQChannel& dataOutChannel = fChannels.at("data-out").at(0);
22 
23  int nEvents = fTree->GetEntries();
24  bool firstRun = true;
25  do
26  {
27  if (firstRun == true){
28  LOG(INFO) << "Number of events: " << nEvents;
29  for (int eventNr = 0; eventNr < nEvents; eventNr++){
30  std::vector<PndSdsHit> tempVector;
31  fTree->GetEntry(eventNr);
32  for (int i = 0; i < fInput->GetEntriesFast(); i++){
33  PndSdsHit* hit = static_cast<PndSdsHit*>(fInput->At(i));
34  if (!hit)
35  continue;
36  tempVector.push_back(*hit);
37  }
38 
39  fHitVector.push_back(tempVector);
40 
41  if (eventNr % 1000 == 0){
42 
43  unique_ptr<FairMQMessage> header(fTransportFactory->CreateMessage(sizeof(int)));
45  memcpy(header->GetData(), &status, sizeof(int));
46  dataOutChannel.SendPart(header);
47 
48  std::ostringstream obuffer;
49  boost::archive::binary_oarchive OutputArchive(obuffer);
50  OutputArchive << fHitVector;
51  int outputSize = obuffer.str().length();
52  unique_ptr<FairMQMessage> msg2(fTransportFactory->CreateMessage(outputSize));
53  memcpy(msg2->GetData(), obuffer.str().c_str(), outputSize);
54  dataOutChannel.Send(msg2);
55 
56 
57  if (fHitVector.size() > 0)
58  fHitVector.clear();
59 
60  if (!CheckCurrentState(RUNNING))
61  {
62  break;
63  }
64  }
65 
66  if (eventNr + 1 == nEvents){
67  unique_ptr<FairMQMessage> header(fTransportFactory->CreateMessage(sizeof(int)));
68  int status = PndMQStatus::RUNNING;
69  memcpy(header->GetData(), &status, sizeof(int));
70  dataOutChannel.SendPart(header);
71 
72  std::ostringstream obuffer;
73  boost::archive::binary_oarchive OutputArchive(obuffer);
74  OutputArchive << fHitVector;
75  int outputSize = obuffer.str().length();
76  unique_ptr<FairMQMessage> msg2(fTransportFactory->CreateMessage(outputSize));
77  memcpy(msg2->GetData(), obuffer.str().c_str(), outputSize);
78  dataOutChannel.Send(msg2);
79 
80 
81  if (fHitVector.size() > 0)
82  fHitVector.clear();
83 
84  if (!CheckCurrentState(RUNNING))
85  {
86  break;
87  }
88  }
89 
90  }
91  firstRun = false;
92  LOG(INFO) << "Finished reading data!";
93 
94  unique_ptr<FairMQMessage> header(fTransportFactory->CreateMessage(sizeof(int)));
95  int status = PndMQStatus::STOP;
96  memcpy(header->GetData(), &status, sizeof(int));
97  dataOutChannel.Send(header);
98  }
99  } while (CheckCurrentState(RUNNING));
100 
101  LOG(INFO) << "I've send " << sendMsgs << " messages!";
102  }
103  else
104  {
105  LOG(ERROR) << " Boost Serialization not ok";
106  }
107 }
Int_t i
Definition: run_full.C:25
vector< vector< PndSdsHit > > fHitVector
Int_t nEvents
Definition: hit_dirc.C:11
int hit(Int_t nEvents=0, TString inFile="sim.root", TString parFile="par.root", TString inDigi="digi.root", TString outFile="hit.root", Int_t timeBased=0)
Definition: hit.C:1
int status[10]
Definition: f_Init.h:28
template<class Archive >
void PndMQFileSamplerHits::serialize ( Archive &  ar,
const unsigned int  version 
)
inline

Definition at line 114 of file PndMQFileSamplerHits.h.

115  {
116  ar& fHitVector;
117  }
vector< vector< PndSdsHit > > fHitVector
virtual void PndMQFileSamplerHits::SetFileName ( std::string  fileName)
inlinevirtual

Definition at line 100 of file PndMQFileSamplerHits.h.

100  {
101  fInputFileName = fileName;
102  }
void PndMQFileSamplerHits::SetProperty ( const int  key,
const std::string &  value 
)
inline

Definition at line 119 of file PndMQFileSamplerHits.h.

Referenced by main().

120  {
121  switch (key)
122  {
123  case InputFile:
124  fInputFileName = value;
125  break;
126  default:
127  FairMQDevice::SetProperty(key, value);
128  break;
129  }
130  }

Friends And Related Function Documentation

friend class boost::serialization::access
friend

Definition at line 154 of file PndMQFileSamplerHits.h.

Member Data Documentation

TList* PndMQFileSamplerHits::fBranchNameList
private

Definition at line 150 of file PndMQFileSamplerHits.h.

bool PndMQFileSamplerHits::fHasBoostSerialization
private

Definition at line 157 of file PndMQFileSamplerHits.h.

Referenced by Run().

vector<vector<PndSdsHit> > PndMQFileSamplerHits::fHitVector
private

Definition at line 156 of file PndMQFileSamplerHits.h.

Referenced by Run().

TFile* PndMQFileSamplerHits::fInFile
private

Definition at line 147 of file PndMQFileSamplerHits.h.

TClonesArray* PndMQFileSamplerHits::fInput
private

Definition at line 149 of file PndMQFileSamplerHits.h.

Referenced by Run().

std::string PndMQFileSamplerHits::fInputFileName
private

Definition at line 151 of file PndMQFileSamplerHits.h.

TTree* PndMQFileSamplerHits::fTree
private

Definition at line 148 of file PndMQFileSamplerHits.h.

Referenced by Run().


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