FairRoot
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
PixelDigiSource.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 // PixelDigiSource.cxx
11 // FAIRROOT
12 //
13 // Created by Radoslaw Karabowicz on 19.02.2016.
14 //
15 //
16 #include "PixelDigiSource.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 <TString.h>
25 #include <cstdlib>
26 #include <cstring>
27 #include <string>
28 
29 PixelDigiSource::PixelDigiSource(TString inputFileName)
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) << "PixelDigiSource created------------";
44 }
45 
47 
49 {
50  // Get input array
52 
53  LOG(info) << "PixelDigiSource::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);
66 
67  if (!fInputFile.is_open()) {
68  LOG(fatal) << "PixelDigiSource::Init() fInputFile \"" << fInputFileName.Data() << "\" could not be open!";
69  return kFALSE;
70  }
71 
72  return kTRUE;
73 }
74 
75 int ReadIntFromString(const std::string& wholestr, const std::string& pattern)
76 {
77  std::string tempstr = wholestr;
78  tempstr.replace(0, tempstr.find(pattern) + pattern.length(), "");
79  tempstr.replace(0, tempstr.find('=') + 1, "");
80  return atoi(tempstr.c_str());
81 }
82 
84 {
85  fDigis->Delete();
86  fNDigis = 0;
87 
88  if (i == 0) {
89  fInputFile.clear();
90  fInputFile.seekg(0, std::ios::beg);
91  }
92 
93  fCurrentEntryNo = i;
94 
95  std::string buffer;
96  LOG(debug) << "PixelDigiSource::ReadEvent() Begin of (" << fDigis->GetEntries() << ")";
97  do {
98  getline(fInputFile, buffer);
99  LOG(debug) << "read from file: \"" << buffer << "\"";
100  if (buffer.find("EVENT BEGIN") == 0) {
101  fRunId = ReadIntFromString(buffer, "RUNID");
102  fMCEntryNo = ReadIntFromString(buffer, "MCENTRYNO");
103  fPartNo = ReadIntFromString(buffer, "PARTNO");
104  fEventHeader->SetRunId(fRunId);
105  fEventHeader->SetMCEntryNumber(fMCEntryNo);
106  fEventHeader->SetPartNo(fPartNo);
107 
108  LOG(debug) << "GOT NEW EVENT " << fMCEntryNo << " (part " << fPartNo << ") with run id = " << fRunId;
109  }
110  if (buffer.find("EVENT") == 0)
111  continue;
112  Int_t detId = atoi(buffer.c_str());
113  buffer.erase(0, buffer.find(' ') + 1);
114  Int_t feId = atoi(buffer.c_str());
115  buffer.erase(0, buffer.find(' ') + 1);
116  Int_t col = atoi(buffer.c_str());
117  buffer.erase(0, buffer.find(' ') + 1);
118  Int_t row = atoi(buffer.c_str());
119  buffer.erase(0, buffer.find(' ') + 1);
120  Double_t charge = atof(buffer.c_str());
121  LOG(debug) << " --/" << fNDigis << "/--> " << detId << " / " << feId << " / " << col << " / " << row
122  << " / " << charge;
123  new ((*fDigis)[fNDigis]) PixelDigi(-1, detId, feId, col, row, charge, 0.);
124  fNDigis++;
125  } while (fInputFile && buffer.compare("EVENT END"));
126  LOG(debug) << "PixelDigiSource::ReadEvent() End of";
127 
128  if (!fInputFile) {
129  return 1;
130  }
131 
132  return 0;
133 }
134 
135 Bool_t PixelDigiSource::ActivateObject(TObject** obj, const char* BrName)
136 {
137  if (strcmp(BrName, "PixelDigis"))
138  *obj = (TObject*)fDigis;
139  else if (strcmp(BrName, "EventHeader."))
140  *obj = (TObject*)fEventHeader;
141  else
142  return kFALSE;
143 
144  return kTRUE;
145 }
146 
147 void PixelDigiSource::Close() { fInputFile.close(); }
148 
150 
151 Int_t PixelDigiSource::CheckMaxEventNo(Int_t /*EvtEnd*/) { return -1; }
152 
154 {
155  ((PixelEventHeader*)feh)->SetRunId(fRunId);
156  ((PixelEventHeader*)feh)->SetMCEntryNumber(fMCEntryNo);
157  ((PixelEventHeader*)feh)->SetPartNo(fPartNo);
158 }
159 
void SetMCEntryNumber(Int_t id)
virtual ~PixelDigiSource()
static FairRootManager * Instance()
int ReadIntFromString(const std::string &wholestr, const std::string &pattern)
void SetPartNo(Int_t ipart)
ClassImp(FairEventBuilder)
virtual Int_t CheckMaxEventNo(Int_t EvtEnd=0)
virtual Bool_t ActivateObject(TObject **obj, const char *BrName)
Int_t ReadEvent(UInt_t i=0)
void Register(const char *name, const char *Foldername, TNamed *obj, Bool_t toFile)
PixelDigiSource(TString inputFileName="test.dat")
virtual void FillEventHeader(FairEventHeader *feh)
void SetRunId(UInt_t runid)