source: trunk/MagicSoft/Mars/mraw/MRawFileWrite.cc@ 4308

Last change on this file since 4308 was 3387, checked in by tbretz, 21 years ago
*** empty log message ***
File size: 8.4 KB
Line 
1/* ======================================================================== *\
2!
3! *
4! * This file is part of MARS, the MAGIC Analysis and Reconstruction
5! * Software. It is distributed to you in the hope that it can be a useful
6! * and timesaving tool in analysing Data of imaging Cerenkov telescopes.
7! * It is distributed WITHOUT ANY WARRANTY.
8! *
9! * Permission to use, copy, modify and distribute this software and its
10! * documentation for any purpose is hereby granted without fee,
11! * provided that the above copyright notice appear in all copies and
12! * that both that copyright notice and this permission notice appear
13! * in supporting documentation. It is provided "as is" without express
14! * or implied warranty.
15! *
16!
17!
18! Author(s): Thomas Bretz 12/2000 <mailto:tbretz@astro.uni-wuerzburg.de>
19!
20! Copyright: MAGIC Software Development, 2000-2004
21!
22!
23\* ======================================================================== */
24
25////////////////////////////////////////////////////////////////////////
26//
27// MRawFileWrite
28//
29// Here we write the root containers which contains the data from a
30// root binary file to a root file. See also MRawFileRead
31//
32// Input Containers:
33// MRawRunHeader, MRawEvtHeader, MRawEvtData, MRawCrateArray, MRawEvtTime
34//
35// Output Containers:
36// -/-
37//
38////////////////////////////////////////////////////////////////////////
39
40#include "MRawFileWrite.h"
41
42#include <TFile.h>
43#include <TTree.h>
44#include <TBranch.h>
45
46#include "MLog.h"
47#include "MLogManip.h"
48
49#include "MParList.h"
50#include "MRawRunHeader.h"
51#include "MRawEvtHeader.h"
52#include "MRawEvtData.h"
53#include "MRawCrateArray.h"
54
55ClassImp(MRawFileWrite);
56
57using namespace std;
58
59// --------------------------------------------------------------------------
60//
61// Default constructor. It opens the output file (root-file)
62//
63MRawFileWrite::MRawFileWrite(const char *fname,
64 const Option_t *opt,
65 const char *ftitle,
66 const Int_t comp,
67 const char *name, const char *title)
68{
69 fName = name ? name : "MRawFileWrite";
70 fTitle = title ? title : "Write task to write DAQ root files";
71
72 //
73 // Open a rootfile
74 //
75 TString str(fname);
76 if (!str.EndsWith(".root", TString::kIgnoreCase))
77 str += ".root";
78
79 fOut = new TFile(str, opt, ftitle, comp);
80}
81
82MRawFileWrite::~MRawFileWrite()
83{
84 //
85 // delete instance, this also does a fOut->Close()
86 //
87 if (fOut->IsOpen())
88 fOut->Write();
89
90 delete fOut;
91
92 //
93 // Remark:
94 // - Trees are automatically deleted by the the file
95 // (unless file.SetDirectory(0) was called)
96 // - Branches are automatically deleted by the tree destructor
97 //
98}
99
100
101// --------------------------------------------------------------------------
102//
103// The PreProcess function checks for the following input containers:
104// - MRawEvtHeader
105// - MRawEvtData
106// - MRawCrateArray
107// - MTime
108// - MRawRunHeader
109// if a container isn't found the eventloop is stopped.
110//
111// The tree which should containe the run header is created. <RunHeaders>
112// The trees which contains the Events <Events>, <PedEvents>, <CalEvents>
113// are created.
114//
115Int_t MRawFileWrite::PreProcess (MParList *pList)
116{
117 //
118 // test whether file is now open or not
119 //
120 if (!fOut->IsOpen())
121 {
122 *fLog << dbginf << "Error: Cannot open file '" << fOut->GetName() << "'" << endl;
123 return kFALSE;
124 }
125
126 //
127 // remember the pointer to the parameter list fur further usage
128 //
129 pParList = pList;
130
131 //
132 // check if MEvtHeader exists in the Parameter list already.
133 // if not create one and add them to the list
134 //
135 fRawEvtHeader = (MRawEvtHeader*)pList->FindObject("MRawEvtHeader");
136 if (!fRawEvtHeader)
137 {
138 *fLog << err << dbginf << "MRawEvtHeader not found... aborting." << endl;
139 return kFALSE;
140 }
141
142 fRawEvtData = (MRawEvtData*)pList->FindObject("MRawEvtData");
143 if (!fRawEvtData)
144 {
145 *fLog << err << dbginf << "MRawEvtData not found... aborting." << endl;
146 return kFALSE;
147 }
148
149 fRawCrateArray = (MRawCrateArray*)pList->FindObject("MRawCrateArray");
150 if (!fRawCrateArray)
151 {
152 *fLog << err << dbginf << "MRawCrateArray not found... aborting." << endl;
153 return kFALSE;
154 }
155
156 fTime = (MTime*)pList->FindObject("MTime");
157 if (!fTime)
158 {
159 *fLog << err << dbginf << "MTime not found... aborting." << endl;
160 return kFALSE;
161 }
162
163 fRawRunHeader = (MRawRunHeader*)pList->FindObject("MRawRunHeader");
164 if (!fRawRunHeader)
165 {
166 *fLog << err << dbginf << "MRawRunHeader not found... aborting." << endl;
167 return kFALSE;
168 }
169
170 //
171 // Remark:
172 // - Trees are automatically deleted by the the file
173 // (unless file.SetDirectory(0) was called)
174 // - Branches are automatically deleted by the tree destructor
175 //
176
177 //
178 // Write the run header information to the file
179 //
180 TTree *rh = new TTree("RunHeaders", "Run headers of all runs in this file");
181 rh->Branch("MRawRunHeader.", "MRawRunHeader", &fRawRunHeader, 32000);
182 rh->Fill();
183 //rh->Write();
184
185 //
186 // create data trees for the three types of data
187 //
188 fTData = new TTree("Events", "Normal Triggered Events");
189 fTPedestal = new TTree("Pedestals", "Pedestal Triggered Events");
190 fTCalibration = new TTree("Calibration", "Calibration Triggered Events");
191
192 //
193 // From the root dicumentation:
194 //
195 // Note that calling TTree::AutoSave too frequently (or similarly calling
196 // TTree::SetAutoSave with a small value) is an expensive operation.
197 // You should make tests for your own application to find a compromize
198 // between speed and the quantity of information you may loose in case of
199 // a job crash.
200 //
201 // In case your program crashes before closing the file holding this tree,
202 // the file will be automatically recovered when you will connect the file
203 // in UPDATE mode.
204 // The Tree will be recovered at the status corresponding to the last AutoSave.
205 //
206 fTData ->SetAutoSave(2000000000); // 2GB
207 fTPedestal ->SetAutoSave(2000000000); // 2GB
208 fTCalibration->SetAutoSave(2000000000); // 2GB
209
210 //
211 // create all branches which are necessary
212 //
213 // FIXME: Can we calculate a good buffer size out of the event size?
214 // using splitlevel=0 sppeds up writing by a factor of 5-10%
215 fTData ->Branch("MTime.", "MTime", &fTime, 32000);
216 fTPedestal ->Branch("MTime.", "MTime", &fTime, 32000);
217 fTCalibration->Branch("MTime.", "MTime", &fTime, 32000);
218 fTData ->Branch("MRawEvtHeader.", "MRawEvtHeader", &fRawEvtHeader, 32000);
219 fTPedestal ->Branch("MRawEvtHeader.", "MRawEvtHeader", &fRawEvtHeader, 32000);
220 fTCalibration->Branch("MRawEvtHeader.", "MRawEvtHeader", &fRawEvtHeader, 32000);
221 fTData ->Branch("MRawEvtData.", "MRawEvtData", &fRawEvtData, 320000);
222 fTPedestal ->Branch("MRawEvtData.", "MRawEvtData", &fRawEvtData, 320000);
223 fTCalibration->Branch("MRawEvtData.", "MRawEvtData", &fRawEvtData, 320000);
224 //fTree->Branch("MRawCrateArray", fRawCrateArray->GetArray(), 32000, 1);
225 fTData ->Branch("MRawCrateArray.", "MRawCrateArray", &fRawCrateArray, 32000);
226 fTPedestal ->Branch("MRawCrateArray.", "MRawCrateArray", &fRawCrateArray, 32000);
227 fTCalibration->Branch("MRawCrateArray.", "MRawCrateArray", &fRawCrateArray, 32000);
228
229 return kTRUE;
230}
231
232// --------------------------------------------------------------------------
233//
234// Gets the trigger type from the run header to decide into which tree the
235// event should be filled in and fills it into this tree.
236//
237Int_t MRawFileWrite::Process()
238{
239 //
240 // get the trigger type of the actual event
241 //
242 const UShort_t type = fRawEvtHeader->GetTrigType();
243
244 //
245 // writa data to the tree. the tree is choosen by the type of the event
246 //
247 switch (type)
248 {
249 case MRawEvtHeader::kTTEvent:
250 fTData->Fill();
251 return kTRUE;
252
253 case MRawEvtHeader::kTTPedestal:
254 fTPedestal->Fill();
255 return kTRUE;
256
257 case MRawEvtHeader::kTTCalibration:
258 fTCalibration->Fill();
259 return kTRUE;
260 }
261
262 *fLog << warn << dbginf << "Got wrong number for the trigger type: " << type;
263 *fLog << " - skipped" << endl;
264
265 return kCONTINUE;
266}
267
Note: See TracBrowser for help on using the repository browser.