source: trunk/Mars/mfileio/MReadMarsFile.cc@ 9969

Last change on this file since 9969 was 9866, checked in by tbretz, 14 years ago
Evaluate return code of call to fRun->Process in MReadMarsFile::Notify()
File size: 8.2 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// MReadMarsFile
28//
29// This task works more or less like MReadTree, but in addition PreProcess
30// reads all the information from the 'RunHeader' tree.
31//
32// Warning: Until now this works only for 'one run header per file'
33//
34/////////////////////////////////////////////////////////////////////////////
35#include "MReadMarsFile.h"
36
37#include <fstream>
38
39#include <TFile.h>
40#include <TTree.h>
41
42#include "MLog.h"
43#include "MLogManip.h"
44
45#include "MParList.h"
46#include "MTaskList.h"
47
48#include "MRawRunHeader.h"
49
50#include "MMcRunHeader.hxx"
51
52ClassImp(MReadMarsFile);
53
54using namespace std;
55
56// --------------------------------------------------------------------------
57//
58// Default constructor. Don't use it.
59//
60MReadMarsFile::MReadMarsFile() : fRun(NULL)
61{
62 fName = "MRead";
63 fTitle = "Read tree and run headers from Mars file.";
64}
65
66// --------------------------------------------------------------------------
67//
68// Default constructor. It creates a MReadTree object to read the
69// RunHeaders and disables Auto Scheme for this tree.
70//
71MReadMarsFile::MReadMarsFile(const char *tname, const char *fname,
72 const char *name, const char *title)
73 : MReadTree(tname, fname)
74{
75 fName = name ? name : "MRead";
76 fTitle = title ? title : "Read tree and run headers from Mars file.";
77
78 //
79 // open the input stream
80 //
81 fRun = new MReadTree("RunHeaders", fname, "RunHeaders");
82
83 //
84 // This disables the auto scheme. because reading new runheader is done
85 // at a low frequency we don't loose time if we always read all
86 // runheaders
87 //
88 fRun->DisableAutoScheme();
89}
90
91// --------------------------------------------------------------------------
92//
93// Destructor. Deleted the MReadTree object for the RunHeaders
94//
95MReadMarsFile::~MReadMarsFile()
96{
97 delete fRun;
98}
99
100Byte_t MReadMarsFile::IsFileValid(const char *name)
101{
102 ifstream fin(name);
103 if (!fin)
104 return 0;
105
106 Char_t c[4];
107 fin.read(c, 4);
108 if (!fin)
109 return 0;
110
111 if (!(c[0]=='r'&& c[1]=='o' && c[2]=='o' && c[3]=='t'))
112 return 0;
113
114 TFile f(name, "READ");
115
116 TTree *t = (TTree*)f.Get("Events");
117 if (!t)
118 return 0;
119
120 // FIXME: Replace numbers by enum! Maybe use bits?
121 if (t->FindBranch("MRawEvtData."))
122 return t->FindBranch("MMcEvt.") ? 2 : 1;
123
124 if (t->FindBranch("MCerPhotEvt."))
125 return t->FindBranch("MMcEvt.") ? 4 : 3;
126
127 return 0;
128}
129
130// --------------------------------------------------------------------------
131//
132// see MReadTree::AddFile, too. The file is also added to the chain for
133// the run headers. If adding file gives a different result for both
134// chains -1 is returned, otherwise the number of files which were added.
135//
136Int_t MReadMarsFile::AddFile(const char *fname, Int_t entries)
137{
138 //
139 // FIXME! A check is missing whether the file already exists or not.
140 //
141 //
142 // returns the number of file which were added
143 //
144 if (!fRun)
145 return -1;
146
147 Int_t n1 = fRun->AddFile(fname);
148 Int_t n2 = MReadTree::AddFile(fname, entries);
149
150 return n1 != n2 ? -1 : n1;
151}
152
153// --------------------------------------------------------------------------
154//
155// Sort the files by their file-names
156//
157void MReadMarsFile::SortFiles()
158{
159 fRun->SortFiles();
160 MReadTree::SortFiles();
161}
162
163// --------------------------------------------------------------------------
164//
165// This overload MReadTree::Notify. Before the MReadTree Notify
166// TObjects are called the RunHeaders of the next files are read.
167//
168// WARNING: This doesn't work correctly yet, if the files are not read in
169// increasing order.
170//
171Bool_t MReadMarsFile::Notify()
172{
173 UInt_t runtype = 0xffff;
174
175 const MRawRunHeader *rawheader = (MRawRunHeader*)fParList->FindObject("MRawRunHeader");
176 if (rawheader)
177 runtype = rawheader->GetRunType();
178
179 //
180 // Try to read the new run headers. If reading the new run header
181 // was successfull call the ReInits
182 //
183 const Int_t idx = GetFileIndex();
184 fRun->SetEventNum(idx<0?0:idx); // Assumption: One Entry per File!
185
186 const Int_t rc = fRun->Process();
187 if (rc!=kTRUE)
188 {
189 *fLog << err << "ERROR - Cannot read new runheaders #" << idx;
190 *fLog << " after reading event #" << GetNumEntry();
191 *fLog << " (MReadTree::Process() returned rc=" << rc << ")" << endl;
192 return kFALSE;
193 }
194
195 if (!MReadTree::Notify())
196 return kFALSE;
197
198 if (rawheader && runtype!=0xffff && runtype != rawheader->GetRunType())
199 {
200 *fLog << warn << "WARNING - You are mixing files with different run types (";
201 *fLog << runtype << ", " << rawheader->GetRunType() << ")" << endl;
202 }
203
204 if (!rawheader->IsValid())
205 {
206 *fLog << warn << "WARNING - The run header (MRawRunHeader) read from the file returns not IsValid()." << endl;
207 *fLog << " Please check if the file contents are ok." << endl;
208 rawheader->Print("header");
209 }
210
211 const MMcRunHeader *mcheader = (MMcRunHeader*)fParList->FindObject("MMcRunHeader");
212 if (mcheader)
213 {
214 if (mcheader->GetCamVersion()==50)
215 {
216 *fLog << warn << "WARNING - You are using a file created with Camera 0.5." << endl;
217 *fLog << "In this camera version some events have undefined Impact-Values" << endl;
218 *fLog << "(MMcEvt::fImpact) Please don't use it for MC studies using the" << endl;
219 *fLog << "impact parameter." << endl;
220 }
221 }
222
223 MTaskList *tlist = (MTaskList*)fParList->FindObject("MTaskList");
224 if (!tlist)
225 {
226 *fLog << err << dbginf << "ERROR - MTaskList not found in Parameter List." << endl;
227 return kFALSE;
228 }
229
230 if (tlist->ReInit())
231 return kTRUE;
232
233 //MReadTree::Notify();
234
235 *fLog << err << "ERROR - ReInit of '" << tlist->GetDescriptor() << "' failed." << endl;
236 return kFALSE;
237}
238
239// --------------------------------------------------------------------------
240//
241// PreProcessed the MReadTree to read the run headers and its base class.
242// see MReadTree::PreProcess for more information
243//
244Int_t MReadMarsFile::PreProcess(MParList *pList)
245{
246 if (!fRun)
247 {
248 *fLog << err << "Error - Cannot use MReadMarsFile instantiated with default constructor... aborting." << endl;
249 return kFALSE;
250 }
251
252 fParList = pList;
253
254 if (!fRun->PreProcess(pList))
255 {
256 *fLog << err << "Error - PreProcessing MReadMarsFile::fRun... aborting." << endl;
257 return kFALSE;
258 }
259
260 // Check if really one file conatains one RunHeader
261 // (FIXME: Do the check more properly)
262 if ((Int_t)fRun->GetEntries() < fRun->GetNumFiles())
263 {
264 *fLog << err << "Error - Number of files exceeds number of RunHeaders... aborting." << endl;
265 return kFALSE;
266 }
267
268 if ((Int_t)fRun->GetEntries() > fRun->GetNumFiles())
269 {
270 *fLog << err << "Error - Number of RunHeaders exceeds number of files... aborting." << endl;
271 return kFALSE;
272 }
273
274 /*
275 const Int_t idx = GetFileIndex();
276 fRun->SetEventNum(idx<0?0:idx); // Assumption: One Entry per File!
277 if (!fRun->Process())
278 {
279 *fLog << err << "Error - Processing MReadMarsFile::fRun... aborting." << endl;
280 return kFALSE;
281 }
282 */
283 return MReadTree::PreProcess(pList);
284}
Note: See TracBrowser for help on using the repository browser.