FairRoot
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
PixelAltDigiWriteToRootVector.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  * PixelAltDigiWriteToRootVector.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 #include "PixelPayload.h"
24 
25 #include <TClonesArray.h>
26 #include <TFile.h>
27 #include <TTree.h>
28 
30  : PixelAltDigiWriteToRootVector("Pixel DigiWriter", 0)
31 {}
32 
34  : PixelAltDigiWriteToRootVector("Pixel DigiWriter", iVerbose)
35 {}
36 
38  : FairTask(name, iVerbose)
39  , fDigis(nullptr)
40  , fOutputRootFile()
41  , fOutputRootTree()
42  , fPixelEventHeader(0)
43  , fPixelDigiVector()
44  , fOutputFileName("test.root")
45  , fNofOutputFiles(0)
46  , fDivideLevel(0)
47  , fRunId(0)
48  , fMCEntryNo(0)
49  , fPartNo(0)
50 {
51  Reset();
52 }
53 
55 
56 void PixelAltDigiWriteToRootVector::Exec(Option_t* /*opt*/)
57 {
58  Reset();
59 
60  Int_t nofDigis = fDigis->GetEntriesFast();
61 
62  for (Int_t ifile = 0; ifile < fNofOutputFiles; ifile++) {
63  fPixelDigiVector.clear();
64  fPixelEventHeader->fRunId = fRunId;
65  fPixelEventHeader->fMCEntryNo = fMCEntryNo;
66  fPixelEventHeader->fPartNo = ifile;
67 
68  // std::cout << ifile << " > " << std::flush;
69  for (Int_t iDigi = 0; iDigi < nofDigis; iDigi++) {
70  PixelDigi* currentDigi = static_cast<PixelDigi*>(fDigis->At(iDigi));
71 
72  int detId = (int)currentDigi->GetDetectorID();
73  Int_t fileToSave = 0;
74  if (fDivideLevel == 1) {
75  fileToSave = detId / 256 - 1;
76  } else if (fDivideLevel == 2) {
77  fileToSave = (detId / 256 - 1) * 4 + (detId % 256 - 1);
78  }
79  // std::cout << fileToSave << " " << std::flush;
80  if (fileToSave != ifile)
81  continue;
82 
83  int feId = (int)currentDigi->GetFeID();
84  int col = (int)currentDigi->GetCol();
85  int row = (int)currentDigi->GetRow();
86  // double charge = (double)currentDigi->GetCharge();
87 
88  PixelPayload::Digi tempDigi; // = new PixelPayload::Digi();
89  tempDigi.fDetectorID = detId;
90  tempDigi.fFeID = feId;
91  tempDigi.fCharge = 0.;
92  tempDigi.fCol = col;
93  tempDigi.fRow = row;
94  fPixelDigiVector.push_back(tempDigi);
95  }
96  fPixelDigiVector.shrink_to_fit();
97  // std::cout << std::endl;
98  fOutputRootTree[ifile]->Fill();
99  }
100  fMCEntryNo++;
101 }
102 
103 void PixelAltDigiWriteToRootVector::SetParContainers() {}
104 
105 InitStatus PixelAltDigiWriteToRootVector::Init()
106 {
107 
108  // Get input array
110 
111  if (!ioman)
112  LOG(fatal) << "No FairRootManager";
113  fDigis = static_cast<TClonesArray*>(ioman->GetObject("PixelDigis"));
114 
115  if (!fDigis)
116  LOG(warn) << "PixelAltDigiWriteToRootVector::Init() No input PixelDigis array!";
117 
118  LOG(info) << "-I- " << fName.Data() << "::Init(). Initialization succesfull.";
119 
120  fRunId = ioman->GetRunId();
121 
122  if (fDivideLevel == 0) {
123  fNofOutputFiles = 1;
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) << "PixelAltDigiWriteToRootVector::Init(), fDivideLevel = " << fDivideLevel
131  << " unknown, it has to be in the range <0,2>";
132  return kFATAL;
133  }
134  }
135 
136  fPixelEventHeader = new PixelPayload::EventHeader();
137  for (Int_t ifile = 0; ifile < fNofOutputFiles; ifile++) {
138  TString fileName = fOutputFileName;
139  TString uniqFile =
140  (fDivideLevel == 0 ? "." : (fDivideLevel == 1 ? Form(".s%d.", ifile) : Form(".c%d.", ifile)));
141  fileName.Replace(fileName.Last('.'), 1, uniqFile.Data());
142 
143  fOutputRootFile[ifile] = new TFile(fileName.Data(), "RECREATE");
144  fOutputRootTree[ifile] = new TTree("fairdata", "example MQ 9 Pixel data");
145  LOG(debug) << "trying to create branch in file >>" << fOutputRootFile[ifile]->GetName() << "<<";
146  fOutputRootTree[ifile]->Branch("EventHeader.", "PixelPayload::EventHeader", &fPixelEventHeader, 0);
147  fOutputRootTree[ifile]->Branch("DigiVector.", &fPixelDigiVector, 32000, 0);
148  LOG(debug) << "branches created, fPixelDigiVector = >>" << &fPixelDigiVector << "<<";
149  }
150 
151  return kSUCCESS;
152 }
153 
154 InitStatus PixelAltDigiWriteToRootVector::ReInit()
155 {
157  InitStatus Status = kFATAL;
158  if (!ioman) {
159  LOG(fatal) << "No FairRootManager found.";
160  } else {
161  fRunId = ioman->GetRunId();
162  fMCEntryNo = 0;
163  Status = kSUCCESS;
164  }
165  return Status;
166 }
167 
168 void PixelAltDigiWriteToRootVector::Reset() {}
169 
170 void PixelAltDigiWriteToRootVector::Finish()
171 {
172  for (Int_t ifile = 0; ifile < fNofOutputFiles; ifile++) {
173  fOutputRootFile[ifile]->Write();
174  fOutputRootFile[ifile]->Close();
175  }
176 }
177 
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
Int_t GetFeID()
Definition: PixelDigi.h:46