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

Last change on this file since 9562 was 4738, checked in by tbretz, 20 years ago
*** empty log message ***
File size: 9.6 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
34// MRawEvtHeader
35// MRawEvtData
36// MRawEvtData2 [MRawEvtData]
37// MRawCrateArray
38// MTime
39//
40// Output Containers:
41// -/-
42//
43////////////////////////////////////////////////////////////////////////
44
45#include "MRawFileWrite.h"
46
47#include <TFile.h>
48#include <TTree.h>
49#include <TBranch.h>
50
51#include "MLog.h"
52#include "MLogManip.h"
53
54#include "MParList.h"
55#include "MRawRunHeader.h"
56#include "MRawEvtHeader.h"
57#include "MRawEvtData.h"
58#include "MRawCrateArray.h"
59
60ClassImp(MRawFileWrite);
61
62using namespace std;
63
64// --------------------------------------------------------------------------
65//
66// Default constructor. It opens the output file (root-file)
67//
68MRawFileWrite::MRawFileWrite(const char *fname, const Option_t *opt,
69 const char *ftitle, const Int_t comp,
70 const char *name, const char *title)
71 : fTData(0), fSplit(kFALSE)
72{
73 fName = name ? name : "MRawFileWrite";
74 fTitle = title ? title : "Write task to write DAQ root files";
75
76 //
77 // Open a rootfile
78 //
79 TString str(fname);
80 if (!str.EndsWith(".root", TString::kIgnoreCase))
81 str += ".root";
82
83 fOut = new TFile(str, opt, ftitle, comp);
84}
85
86MRawFileWrite::~MRawFileWrite()
87{
88 if (fTData)
89 {
90 if (fOut != fTData->GetCurrentFile())
91 {
92 *fLog << warn << endl;
93 *fLog << "WARNING - MWriteRootFile: Root's TTree/TFile has opened a new file" << endl;
94 *fLog << " automatically. You can change this behaviour using TTree::SetMaxTreeSize." << endl;
95 *fLog << " You won't be able to read splitted files correctly with MReadMarsFile if" << endl;
96 *fLog << " they have more than one entry in 'RunHeaders' or you try to read more than" << endl;
97 *fLog << " one of such sequences at once." << endl;
98 *fLog << endl;
99 }
100
101 //
102 // For more information see TTree:ChangeFile()
103 //
104 fOut = fTData->GetCurrentFile();
105 }
106
107 //
108 // delete instance, this also does a fOut->Close()
109 //
110 if (fOut->IsOpen())
111 fOut->Write();
112
113 delete fOut;
114
115 //
116 // Remark:
117 // - Trees are automatically deleted by the the file
118 // (unless file.SetDirectory(0) was called)
119 // - Branches are automatically deleted by the tree destructor
120 //
121}
122
123
124// --------------------------------------------------------------------------
125//
126// The PreProcess function checks for the following input containers:
127// - MRawEvtHeader
128// - MRawEvtData
129// - MRawCrateArray
130// - MTime
131// - MRawRunHeader
132// if a container isn't found the eventloop is stopped.
133//
134// The tree which should containe the run header is created. <RunHeaders>
135// The trees which contains the Events <Events>, <PedEvents>, <CalEvents>
136// are created.
137//
138Int_t MRawFileWrite::PreProcess (MParList *pList)
139{
140 //
141 // test whether file is now open or not
142 //
143 if (!fOut->IsOpen())
144 {
145 *fLog << dbginf << "Error: Cannot open file '" << fOut->GetName() << "'" << endl;
146 return kFALSE;
147 }
148
149 //
150 // check if MEvtHeader exists in the Parameter list already.
151 // if not create one and add them to the list
152 //
153 fRawEvtHeader = (MRawEvtHeader*)pList->FindObject("MRawEvtHeader");
154 if (!fRawEvtHeader)
155 {
156 *fLog << err << "MRawEvtHeader not found... aborting." << endl;
157 return kFALSE;
158 }
159
160 fRawEvtData1 = (MRawEvtData*)pList->FindObject("MRawEvtData");
161 if (!fRawEvtData1)
162 {
163 *fLog << err << "MRawEvtData not found... aborting." << endl;
164 return kFALSE;
165 }
166 fRawEvtData2 = (MRawEvtData*)pList->FindObject("MRawEvtData2", "MRawEvtData");
167 if (!fRawEvtData2)
168 {
169 *fLog << err << "MRawEvtData2 [MRawEvtData] not found... aborting." << endl;
170 return kFALSE;
171 }
172
173 fRawCrateArray = (MRawCrateArray*)pList->FindObject("MRawCrateArray");
174 if (!fRawCrateArray)
175 {
176 *fLog << err << "MRawCrateArray not found... aborting." << endl;
177 return kFALSE;
178 }
179
180 fTime = (MTime*)pList->FindObject("MTime");
181 if (!fTime)
182 {
183 *fLog << err << "MTime not found... aborting." << endl;
184 return kFALSE;
185 }
186
187 fRawRunHeader = (MRawRunHeader*)pList->FindObject("MRawRunHeader");
188 if (!fRawRunHeader)
189 {
190 *fLog << err << "MRawRunHeader not found... aborting." << endl;
191 return kFALSE;
192 }
193
194 //
195 // Remark:
196 // - Trees are automatically deleted by the the file
197 // (unless file.SetDirectory(0) was called)
198 // - Branches are automatically deleted by the tree destructor
199 //
200 fTRunHeader = new TTree("RunHeaders", "Run headers of all runs in this file");
201 fTRunHeader->Branch("MRawRunHeader.", "MRawRunHeader", &fRawRunHeader, 32000);
202
203 //
204 // create data trees for the three types of data
205 //
206 fTData = new TTree("Events", "Normal Triggered Events");
207 fTPedestal = fSplit ? new TTree("Pedestals", "Pedestal Triggered Events") : fTData;
208 fTCalibration = fSplit ? new TTree("Calibration", "Calibration Triggered Events") : fTData;
209
210 //
211 // From the root dicumentation:
212 //
213 // Note that calling TTree::AutoSave too frequently (or similarly calling
214 // TTree::SetAutoSave with a small value) is an expensive operation.
215 // You should make tests for your own application to find a compromize
216 // between speed and the quantity of information you may loose in case of
217 // a job crash.
218 //
219 // In case your program crashes before closing the file holding this tree,
220 // the file will be automatically recovered when you will connect the file
221 // in UPDATE mode.
222 // The Tree will be recovered at the status corresponding to the last AutoSave.
223 //
224 fTData ->SetAutoSave(2000000000); // 2GB
225 fTPedestal ->SetAutoSave(2000000000); // 2GB
226 fTCalibration->SetAutoSave(2000000000); // 2GB
227
228 //
229 // create all branches which are necessary
230 //
231 // FIXME: Can we calculate a good buffer size out of the event size?
232 // using splitlevel=0 sppeds up writing by a factor of 5-10%
233 fTData->Branch("MTime.", "MTime", &fTime, 32000);
234 fTData->Branch("MRawEvtHeader.", "MRawEvtHeader", &fRawEvtHeader, 32000);
235 fTData->Branch("MRawEvtData.", "MRawEvtData", &fRawEvtData1, 320000);
236 fTData->Branch("MRawEvtData2.", "MRawEvtData", &fRawEvtData2, 320000);
237 fTData->Branch("MRawCrateArray.", "MRawCrateArray", &fRawCrateArray, 32000);
238 if (!fSplit)
239 return kTRUE;
240
241 fTPedestal ->Branch("MTime.", "MTime", &fTime, 32000);
242 fTCalibration->Branch("MTime.", "MTime", &fTime, 32000);
243 fTPedestal ->Branch("MRawEvtHeader.", "MRawEvtHeader", &fRawEvtHeader, 32000);
244 fTCalibration->Branch("MRawEvtHeader.", "MRawEvtHeader", &fRawEvtHeader, 32000);
245 fTPedestal ->Branch("MRawEvtData.", "MRawEvtData", &fRawEvtData1, 320000);
246 fTCalibration->Branch("MRawEvtData.", "MRawEvtData", &fRawEvtData1, 320000);
247 fTPedestal ->Branch("MRawEvtData2.", "MRawEvtData", &fRawEvtData2, 320000);
248 fTCalibration->Branch("MRawEvtData2.", "MRawEvtData", &fRawEvtData2, 320000);
249 fTPedestal ->Branch("MRawCrateArray.", "MRawCrateArray", &fRawCrateArray, 32000);
250 fTCalibration->Branch("MRawCrateArray.", "MRawCrateArray", &fRawCrateArray, 32000);
251
252 return kTRUE;
253}
254
255// --------------------------------------------------------------------------
256//
257// Write the run header information to the file
258//
259Bool_t MRawFileWrite::ReInit(MParList *pList)
260{
261 fTRunHeader->Fill();
262 return kTRUE;
263}
264
265// --------------------------------------------------------------------------
266//
267// Gets the trigger type from the run header to decide into which tree the
268// event should be filled in and fills it into this tree.
269//
270Int_t MRawFileWrite::Process()
271{
272 //
273 // get the trigger type of the actual event
274 //
275 const UShort_t type = fRawEvtHeader->GetTrigType();
276
277 //
278 // writa data to the tree. the tree is choosen by the type of the event
279 //
280 switch (type)
281 {
282 case MRawEvtHeader::kTTEvent:
283 fTData->Fill();
284 return kTRUE;
285
286 case MRawEvtHeader::kTTPedestal:
287 fTPedestal->Fill();
288 return kTRUE;
289
290 case MRawEvtHeader::kTTCalibration:
291 fTCalibration->Fill();
292 return kTRUE;
293 }
294
295 *fLog << warn << dbginf << "Got wrong number for the trigger type: " << type;
296 *fLog << " - skipped" << endl;
297
298 return kCONTINUE;
299}
300
Note: See TracBrowser for help on using the repository browser.