FairRoot
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
PixelDigiBinSource.cxx
Go to the documentation of this file.
1 /********************************************************************************
2  * Copyright (C) 2014 GSI Helmholtzzentrum fuer Schwerionenforschung GmbH *
3  * *
4  * This software is distributed under the terms of the *
5  * GNU Lesser General Public Licence (LGPL) version 3, *
6  * copied verbatim in the file "LICENSE" *
7  ********************************************************************************/
8 
9 //
10 // PixelDigiBinSource.cxx
11 // FAIRROOT
12 //
13 // Created by Radoslaw Karabowicz on 19.02.2016.
14 //
15 //
16 #include "PixelDigiBinSource.h"
17 
18 #include "FairLogger.h"
19 #include "FairRootManager.h"
20 #include "PixelDigi.h"
21 #include "PixelEventHeader.h"
22 
23 #include <TClonesArray.h>
24 #include <TObject.h>
25 #include <TString.h>
26 #include <cstring>
27 #include <string>
28 
30  : FairSource()
31  , fEventHeader(nullptr)
32  , fDigis(nullptr)
33  , fNDigis(0)
34  , fTNofEvents(0)
35  , fTNofDigis(0)
36  , fInputFileName(inputFileName)
37  , fInputFile()
38  , fCurrentEntryNo(0)
39  , fRunId(0)
40  , fMCEntryNo(0)
41  , fPartNo(0)
42 {
43  LOG(debug) << "PixelDigiBinSource created------------";
44 }
45 
47 
49 {
50  // Get input array
52 
53  LOG(info) << "PixelDigiBinSource::Init";
54  if (!ioman)
55  LOG(fatal) << "No FairRootManager";
56 
57  // Register output array StsDigi
58  fDigis = new TClonesArray("PixelDigi", 10000);
59  ioman->Register("PixelDigis", "Pixel", fDigis, kFALSE);
60 
61  fEventHeader = new PixelEventHeader();
62  fEventHeader->SetName("EventHeader.");
63  ioman->Register("EventHeader.", "EvtHeader", fEventHeader, kFALSE);
64 
65  fInputFile.open(fInputFileName.Data(), std::fstream::in | std::fstream::binary);
66 
67  if (!fInputFile.is_open()) {
68  LOG(fatal) << "PixelDigiBinSource::Init() fInputFile \"" << fInputFileName.Data() << "\" could not be open!";
69  return kFALSE;
70  }
71 
72  return kTRUE;
73 }
74 
76 {
77  fDigis->Clear();
78  fNDigis = 0;
79 
80  if (!fInputFile) {
81  return 1;
82  }
83 
84  if (i == 0) {
85  fInputFile.clear();
86  fInputFile.seekg(0, std::ios::beg);
87  }
88 
89  fCurrentEntryNo = i;
90 
91  std::string buffer;
92  LOG(debug) << "PixelDigiBinSource::ReadEvent() Begin of (" << fDigis->GetEntries() << ")";
93 
94  Int_t head[4]; // runId, MCEntryNo, PartNo, NofDigis
95  fInputFile.read((char*)head, sizeof(head));
96 
97  if (fInputFile.eof()) {
98  LOG(info) << "End of file reached!";
99  return 1;
100  }
101 
102  Int_t dataSize = 4; // detId, feId, col, row
103 
104  const Int_t constNofData = head[3] * dataSize;
105  short int dataCont[constNofData];
106  fInputFile.read((char*)dataCont, sizeof(dataCont));
107 
108  fRunId = head[0];
109  fMCEntryNo = head[1];
110  fPartNo = head[2];
111  fEventHeader->SetRunId(fRunId);
112  fEventHeader->SetMCEntryNumber(fMCEntryNo);
113  fEventHeader->SetPartNo(fPartNo);
114 
115  for (Int_t idata = 0; idata < head[3]; idata++) {
116  LOG(debug) << " --/" << idata << "/--> " << dataCont[idata * dataSize + 0] << " / "
117  << dataCont[idata * dataSize + 1] << " / " << dataCont[idata * dataSize + 2] << " / "
118  << dataCont[idata * dataSize + 3] << " / "
119  << " 0.";
120  new ((*fDigis)[fNDigis]) PixelDigi(-1,
121  (Int_t)dataCont[idata * dataSize + 0],
122  (Int_t)dataCont[idata * dataSize + 1],
123  (Int_t)dataCont[idata * dataSize + 2],
124  (Int_t)dataCont[idata * dataSize + 3],
125  0.,
126  0.);
127  fNDigis++;
128  }
129  LOG(debug) << "PixelDigiBinSource::ReadEvent() End of";
130 
131  if (!fInputFile) {
132  return 1;
133  }
134 
135  return 0;
136 }
137 
138 Bool_t PixelDigiBinSource::ActivateObject(TObject** obj, const char* BrName)
139 {
140  if (strcmp(BrName, "PixelDigis") == 0)
141  *obj = (TObject*)fDigis;
142  else if (strcmp(BrName, "EventHeader.") == 0)
143  *obj = (TObject*)fEventHeader;
144  else
145  return kFALSE;
146 
147  return kTRUE;
148 }
149 
150 void PixelDigiBinSource::Close() { fInputFile.close(); }
151 
153 
154 Int_t PixelDigiBinSource::CheckMaxEventNo(Int_t /*EvtEnd*/) { return -1; }
155 
157 {
158  ((PixelEventHeader*)feh)->SetRunId(fRunId);
159  ((PixelEventHeader*)feh)->SetMCEntryNumber(fMCEntryNo);
160  ((PixelEventHeader*)feh)->SetPartNo(fPartNo);
161 }
162 
void SetMCEntryNumber(Int_t id)
static FairRootManager * Instance()
void SetPartNo(Int_t ipart)
ClassImp(FairEventBuilder)
virtual Int_t CheckMaxEventNo(Int_t EvtEnd=0)
Int_t ReadEvent(UInt_t i=0)
PixelDigiBinSource(TString inputFileName="test.dat")
void Register(const char *name, const char *Foldername, TNamed *obj, Bool_t toFile)
virtual Bool_t ActivateObject(TObject **obj, const char *BrName)
virtual void FillEventHeader(FairEventHeader *feh)
void SetRunId(UInt_t runid)