FairRoot
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
PixelDigiWriteToBinFile.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  * PixelDigiWriteToBinFile.h
10  *
11  * Created on: 19.02.2016
12  * Author: R. Karabowicz
13  */
14 
16 
17 // Includes from base
18 #include "FairLogger.h"
19 #include "FairRootManager.h"
20 
21 // Includes from ROOT
22 #include "PixelDigi.h"
23 
24 #include <TClonesArray.h>
25 #include <TString.h>
26 
28  : PixelDigiWriteToBinFile("Pixel DigiWriter", 0)
29 {}
30 
32  : PixelDigiWriteToBinFile("Pixel DigiWriter", iVerbose)
33 {}
34 
35 PixelDigiWriteToBinFile::PixelDigiWriteToBinFile(const char* name, Int_t iVerbose)
36  : FairTask(name, iVerbose)
37  , fDigis(nullptr)
38  , fOutputFileName("test.dat")
39  , fNofOutputFiles(0)
40  , fOutputFiles()
41  , fDivideLevel(0)
42  , fRunId(0)
43  , fMCEntryNo(0)
44  , fPartNo(0)
45 
46 {
47  Reset();
48 }
49 
51 
52 void PixelDigiWriteToBinFile::Exec(Option_t* /*opt*/)
53 {
54  Reset();
55 
56  Int_t nofDigis = fDigis->GetEntriesFast();
57 
58  Int_t digisPerFile[12] = {0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0};
59 
60  for (Int_t iDigi = 0; iDigi < nofDigis; iDigi++) {
61  PixelDigi* currentDigi = static_cast<PixelDigi*>(fDigis->At(iDigi));
62 
63  short int detId = (short int)currentDigi->GetDetectorID();
64 
65  Int_t fileToSave = 0;
66  if (fDivideLevel == 1) {
67  fileToSave = detId / 256 - 1;
68  } else if (fDivideLevel == 2) {
69  fileToSave = (detId / 256 - 1) * 4 + (detId % 256 - 1);
70  }
71  digisPerFile[fileToSave] += 1;
72  }
73 
74  for (Int_t ifile = 0; ifile < fNofOutputFiles; ifile++) {
75  fOutputFiles[ifile].write((char*)&fRunId, sizeof(fRunId));
76  fOutputFiles[ifile].write((char*)&fMCEntryNo, sizeof(fMCEntryNo));
77  fOutputFiles[ifile].write((char*)&ifile, sizeof(ifile));
78  fOutputFiles[ifile].write((char*)&digisPerFile[ifile], sizeof(digisPerFile[ifile]));
79  }
80 
81  for (Int_t iDigi = 0; iDigi < nofDigis; iDigi++) {
82  PixelDigi* currentDigi = static_cast<PixelDigi*>(fDigis->At(iDigi));
83 
84  short int detId = (short int)currentDigi->GetDetectorID();
85  short int feId = (short int)currentDigi->GetFeID();
86  short int col = (short int)currentDigi->GetCol();
87  short int row = (short int)currentDigi->GetRow();
88 
89  Int_t fileToSave = 0;
90  if (fDivideLevel == 1) {
91  fileToSave = detId / 256 - 1;
92  } else if (fDivideLevel == 2) {
93  fileToSave = (detId / 256 - 1) * 4 + (detId % 256 - 1);
94  }
95  fOutputFiles[fileToSave].write((char*)&detId, sizeof(detId));
96  fOutputFiles[fileToSave].write((char*)&feId, sizeof(feId));
97  fOutputFiles[fileToSave].write((char*)&col, sizeof(col));
98  fOutputFiles[fileToSave].write((char*)&row, sizeof(row));
99  }
100 
101  fMCEntryNo++;
102 }
103 
104 InitStatus PixelDigiWriteToBinFile::Init()
105 {
106 
107  // Get input array
109 
110  if (!ioman)
111  LOG(fatal) << "No FairRootManager";
112  fDigis = static_cast<TClonesArray*>(ioman->GetObject("PixelDigis"));
113 
114  if (!fDigis)
115  LOG(warn) << "PixelDigiWriteToBinFile::Init() No input PixelDigis array!";
116 
117  LOG(info) << "-I- " << fName.Data() << "::Init(). Initialization succesfull.";
118 
119  fRunId = ioman->GetRunId();
120 
121  if (fDivideLevel == 0) {
122  fNofOutputFiles = 1;
123  fOutputFiles[0].open(fOutputFileName.Data(), std::fstream::out | std::fstream::binary);
124  } else {
125  if (fDivideLevel == 1) {
126  fNofOutputFiles = 3; // 1 file per station (3 stations)
127  } else if (fDivideLevel == 2) {
128  fNofOutputFiles = 12; // 1 file per sensor (3 stations times 4 sensors)
129  } else {
130  LOG(fatal) << "PixelDigiWriteToBinFile::Init(), fDivideLevel = " << fDivideLevel
131  << " unknown, it has to be in the range <0,2>";
132  return kFATAL;
133  }
134  for (Int_t ifile = 0; ifile < fNofOutputFiles; ifile++) {
135  TString fileName = fOutputFileName;
136  TString uniqFile = Form(".p%d.", ifile);
137  fileName.Replace(fileName.Last('.'), 1, uniqFile.Data());
138  fOutputFiles[ifile].open(fileName.Data(), std::fstream::out);
139  }
140  }
141 
142  return kSUCCESS;
143 }
144 
145 InitStatus PixelDigiWriteToBinFile::ReInit()
146 {
148 
149  InitStatus Status = kFATAL;
150  if (!ioman) {
151  LOG(fatal) << "No FairRootManager found.";
152  } else {
153  fRunId = ioman->GetRunId();
154  fMCEntryNo = 0;
155  Status = kSUCCESS;
156  }
157  return Status;
158 }
159 
160 void PixelDigiWriteToBinFile::Finish()
161 {
162  for (Int_t ifile = 0; ifile < fNofOutputFiles; ifile++) {
163  fOutputFiles[ifile].close();
164  }
165 }
166 
InitStatus
Definition: FairTask.h:33
Int_t GetRow()
Definition: PixelDigi.h:49
static FairRootManager * Instance()
ClassImp(FairEventBuilder)
TObject * GetObject(const char *BrName)
Int_t GetDetectorID()
Definition: PixelDigi.h:45
Int_t GetCol()
Definition: PixelDigi.h:48
virtual void Exec(Option_t *opt)
Int_t GetFeID()
Definition: PixelDigi.h:46