source: trunk/MagicSoft/Mars/mfileio/MReadMarsFile.cc@ 6230

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