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

Last change on this file since 2173 was 2173, checked in by tbretz, 21 years ago
*** empty log message ***
File size: 8.3 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@uni-sw.gwdg.de>
19!
20! Copyright: MAGIC Software Development, 2000-2001
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 fOut = new TFile(fname, opt, ftitle, comp);
76}
77
78MRawFileWrite::~MRawFileWrite()
79{
80 //
81 // delete instance, this also does a fOut->Close()
82 //
83 if (fOut->IsOpen())
84 fOut->Write();
85
86 delete fOut;
87
88 //
89 // Remark:
90 // - Trees are automatically deleted by the the file
91 // (unless file.SetDirectory(0) was called)
92 // - Branches are automatically deleted by the tree destructor
93 //
94}
95
96
97// --------------------------------------------------------------------------
98//
99// The PreProcess function checks for the following input containers:
100// - MRawEvtHeader
101// - MRawEvtData
102// - MRawCrateArray
103// - MRawEvtTime <MTime>
104// - MRawRunHeader
105// if a container isn't found the eventloop is stopped.
106//
107// The tree which should containe the run header is created. <RunHeaders>
108// The trees which contains the Events <Events>, <PedEvents>, <CalEvents>
109// are created.
110//
111Bool_t MRawFileWrite::PreProcess (MParList *pList)
112{
113 //
114 // test whether file is now open or not
115 //
116 if (!fOut->IsOpen())
117 {
118 *fLog << dbginf << "Error: Cannot open file '" << fOut->GetName() << "'" << endl;
119 return kFALSE;
120 }
121
122 //
123 // remember the pointer to the parameter list fur further usage
124 //
125 pParList = pList;
126
127 //
128 // check if MEvtHeader exists in the Parameter list already.
129 // if not create one and add them to the list
130 //
131 fRawEvtHeader = (MRawEvtHeader*)pList->FindObject("MRawEvtHeader");
132 if (!fRawEvtHeader)
133 {
134 *fLog << err << dbginf << "MRawEvtHeader not found... aborting." << endl;
135 return kFALSE;
136 }
137
138 fRawEvtData = (MRawEvtData*)pList->FindObject("MRawEvtData");
139 if (!fRawEvtData)
140 {
141 *fLog << err << dbginf << "MRawEvtData not found... aborting." << endl;
142 return kFALSE;
143 }
144
145 fRawCrateArray = (MRawCrateArray*)pList->FindObject("MRawCrateArray");
146 if (!fRawCrateArray)
147 {
148 *fLog << err << dbginf << "MRawCrateArray not found... aborting." << endl;
149 return kFALSE;
150 }
151
152 fRawEvtTime = (MTime*)pList->FindObject("MRawEvtTime");
153 if (!fRawEvtTime)
154 {
155 *fLog << err << dbginf << "MRawEvtTime not found... aborting." << endl;
156 return kFALSE;
157 }
158
159 fRawRunHeader = (MRawRunHeader*)pList->FindObject("MRawRunHeader");
160 if (!fRawRunHeader)
161 {
162 *fLog << err << dbginf << "MRawRunHeader not found... aborting." << endl;
163 return kFALSE;
164 }
165
166 //
167 // Remark:
168 // - Trees are automatically deleted by the the file
169 // (unless file.SetDirectory(0) was called)
170 // - Branches are automatically deleted by the tree destructor
171 //
172
173 //
174 // Write the run header information to the file
175 //
176 TTree *rh = new TTree("RunHeaders", "Run headers of all runs in this file");
177 rh->Branch("MRawRunHeader", "MRawRunHeader", &fRawRunHeader, 32000);
178 rh->Fill();
179 //rh->Write();
180
181 //
182 // create data trees for the three types of data
183 //
184 fTData = new TTree("Events", "Normal Triggered Events");
185 fTPedestal = new TTree("PedEvents", "Pedestal Triggered Events");
186 fTCalibration = new TTree("CalEvents", "Calibration Triggered Events");
187
188 //
189 // From the root dicumentation:
190 //
191 // Note that calling TTree::AutoSave too frequently (or similarly calling
192 // TTree::SetAutoSave with a small value) is an expensive operation.
193 // You should make tests for your own application to find a compromize
194 // between speed and the quantity of information you may loose in case of
195 // a job crash.
196 //
197 // In case your program crashes before closing the file holding this tree,
198 // the file will be automatically recovered when you will connect the file
199 // in UPDATE mode.
200 // The Tree will be recovered at the status corresponding to the last AutoSave.
201 //
202 fTData ->SetAutoSave(2000000000); // 2GB
203 fTPedestal ->SetAutoSave(2000000000); // 2GB
204 fTCalibration->SetAutoSave(2000000000); // 2GB
205
206 //
207 // create all branches which are necessary
208 //
209 // FIXME: Can we calculate a good buffer size out of the event size?
210 // using splitlevel=0 sppeds up writing by a factor of 5-10%
211 fTData ->Branch("MTime.", "MTime", &fRawEvtTime, 32000);
212 fTPedestal ->Branch("MTime.", "MTime", &fRawEvtTime, 32000);
213 fTCalibration->Branch("MTime.", "MTime", &fRawEvtTime, 32000);
214 fTData ->Branch("MRawEvtHeader.", "MRawEvtHeader", &fRawEvtHeader, 32000);
215 fTPedestal ->Branch("MRawEvtHeader.", "MRawEvtHeader", &fRawEvtHeader, 32000);
216 fTCalibration->Branch("MRawEvtHeader.", "MRawEvtHeader", &fRawEvtHeader, 32000);
217 fTData ->Branch("MRawEvtData.", "MRawEvtData", &fRawEvtData, 320000);
218 fTPedestal ->Branch("MRawEvtData.", "MRawEvtData", &fRawEvtData, 320000);
219 fTCalibration->Branch("MRawEvtData.", "MRawEvtData", &fRawEvtData, 320000);
220 //fTree->Branch("MRawCrateArray", fRawCrateArray->GetArray(), 32000, 1);
221 fTData ->Branch("MRawCrateArray.", "MRawCrateArray", &fRawCrateArray, 32000);
222 fTPedestal ->Branch("MRawCrateArray.", "MRawCrateArray", &fRawCrateArray, 32000);
223 fTCalibration->Branch("MRawCrateArray.", "MRawCrateArray", &fRawCrateArray, 32000);
224
225 return kTRUE;
226}
227
228// --------------------------------------------------------------------------
229//
230// Gets the trigger type from the run header to decide into which tree the
231// event should be filled in and fills it into this tree.
232//
233Bool_t MRawFileWrite::Process()
234{
235 //
236 // get the trigger type of the actual event
237 //
238 const UShort_t type = fRawEvtHeader->GetTrigType();
239
240 //
241 // writa data to the tree. the tree is choosen by the type of the event
242 //
243 switch (type)
244 {
245 case kTTEvent:
246 fTData->Fill();
247 return kTRUE;
248
249 case kTTPedestal:
250 fTPedestal->Fill();
251 return kTRUE;
252
253 case kTTCalibration:
254 fTCalibration->Fill();
255 return kTRUE;
256 }
257
258 *fLog << warn << dbginf << "Got wrong number for the trigger type: " << type;
259 *fLog << " - skipped" << endl;
260
261 return kCONTINUE;
262}
263
Note: See TracBrowser for help on using the repository browser.