source: trunk/MagicSoft/Mars/mfileio/MCT1ReadPreProc.cc@ 2313

Last change on this file since 2313 was 2296, checked in by tbretz, 21 years ago
*** empty log message ***
File size: 37.1 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, 11/2002 <mailto:tbretz@astro.uni-wuerzburg.de>
19!
20! Copyright: MAGIC Software Development, 2000-2003
21!
22!
23\* ======================================================================== */
24
25/////////////////////////////////////////////////////////////////////////////
26//
27// MCT1ReadPreProc
28//
29// Reads a output file of the CT1 preproc.
30//
31// Implements usage of a selector (see MRead) Use such a filter to skip
32// events before reading! But never use a filter which needs read data
33// as input...
34//
35// Input Containers:
36// -/-
37//
38// Output Containers:
39// MCerPhotEvt the data container for all data.
40// MPedestalCam ct1 pedestals
41// MMcEvt monte carlo data container for MC files
42// MMcTrig mc data container for trigger information
43// MSrcPosCam source position in the camera
44// MBlindPixels Array holding blind pixels
45//
46/////////////////////////////////////////////////////////////////////////////
47#include "MCT1ReadPreProc.h"
48
49#include <fstream>
50
51#include <TList.h>
52#include <TSystem.h>
53
54#define LINUX
55#define HISTO void
56#define HBOOK_FILE int
57#include "defines.h"
58#include "structures.h"
59
60#include "MLog.h"
61#include "MLogManip.h"
62
63#include "MTime.h"
64#include "MFilter.h"
65
66#include "MParList.h"
67#include "MCerPhotEvt.h"
68
69#include "MPedestalPix.h"
70#include "MPedestalCam.h"
71
72#include "MGeomCam.h"
73#include "MSrcPosCam.h"
74#include "MBlindPixels.h"
75
76#include "MRawRunHeader.h"
77#include "MTaskList.h"
78
79#include "MMcEvt.hxx"
80#include "MMcTrig.hxx"
81#include "MBinning.h"
82
83#include "TRandom3.h"
84#include "MParameters.h"
85
86ClassImp(MCT1ReadPreProc);
87
88using namespace std;
89
90// --------------------------------------------------------------------------
91//
92// Default constructor. Creates an array which stores the file names of
93// the files which should be read. If a filename is given it is added
94// to the list.
95//
96MCT1ReadPreProc::MCT1ReadPreProc(const char *fname, const char *name,
97 const char *title) : fIn(NULL), fEntries(0)
98{
99 fName = name ? name : "MRead";
100 fTitle = title ? title : "Reads a CT1 preproc data file";
101
102 //
103 // remember file name for opening the file in the preprocessor
104 //
105 fFileNames = new TList;
106 fFileNames->SetOwner();
107
108 if (fname)
109 AddFile(fname);
110}
111
112// --------------------------------------------------------------------------
113//
114// Delete the filename list and the input stream if one exists.
115//
116MCT1ReadPreProc::~MCT1ReadPreProc()
117{
118 delete fFileNames;
119 if (fIn)
120 delete fIn;
121}
122
123// --------------------------------------------------------------------------
124//
125// Add this file as the last entry in the chain
126//
127Int_t MCT1ReadPreProc::AddFile(const char *txt, int)
128{
129 const char *name = gSystem->ExpandPathName(txt);
130
131 TString fname(name);
132 delete [] name;
133
134 if (!CheckHeader(fname))
135 {
136 *fLog << warn << "WARNING - Problem reading header... ignored." << endl;
137 return 0;
138 }
139
140 const Int_t n = GetNumEvents(fname);
141 if (n==0)
142 {
143 *fLog << warn << "WARNING - File contains no data... ignored." << endl;
144 return 0;
145 }
146
147 fEntries += n;
148
149 *fLog << inf << "File " << txt << " contains " << n << " events (Total=" << fEntries << ")" << endl;
150
151 fFileNames->AddLast(new TNamed(txt, ""));
152 return 1;
153}
154
155// --------------------------------------------------------------------------
156//
157// Print data from the header to the screen and analyse the header data,
158// means store and copy the needed data into Mars structures and
159// data members
160//
161void MCT1ReadPreProc::ProcessRunHeader(const struct outputpars &outpars)
162{
163 if (outpars.inumpixels != iMAXNUMPIX)
164 *fLog << warn << "WARNING! File doesn't contain " << iMAXNUMPIX << " Pixels... maybe corrupt." << endl;
165
166 fNumEventsInRun = 0;
167
168 //
169 // ------------------- Output some stuff -----------------------
170 //
171
172 // int itelescope; // number of the CT which took the data
173 *fLog << inf << "Telescope: CT" << outpars.itelescope;
174
175 // float flongitude_deg; // longitude (counted positive towards West) of CT position */
176 // float flatitude_deg; // latitude (counted positive towards North) of CT position */
177 *fLog << " located @ Longitude=" << outpars.flongitude_deg;
178 *fLog << "deg Latitude=" << outpars.flatitude_deg << "deg" << endl;
179
180 // int irunnum; // run number (from parameters file)
181 fRawRunHeader->SetRunNumber(outpars.irunnum);
182 fRawRunHeader->SetReadyToSave();
183
184 // enum onoroff {NEITHER_ON_NOR_OFF, OFF_SOURCE, ON_SOURCE} eruntype; // runtype
185 *fLog << "Run: #" << outpars.irunnum << " (";
186 switch (outpars.eruntype)
187 {
188 case NEITHER_ON_NOR_OFF: *fLog << "unknown"; break;
189 case OFF_SOURCE: *fLog << "off-source"; break;
190 case ON_SOURCE: *fLog << "on-source"; break;
191 default: *fLog << (int)outpars.eruntype; break;
192 }
193 *fLog << ", ";
194 switch (outpars.etrackmode)
195 {
196 case NORMAL: *fLog << "normal tracking"; break;
197 case REVERSE: *fLog << "reverse tracking"; break;
198 case DUNNO: *fLog << "unknown tracking"; break;
199 default: *fLog << (int)outpars.etrackmode; break;
200 }
201 *fLog << ")" << endl;
202
203 //double dsourcera_hours; // right ascension of observed source in hours
204 //double dsourcedec_deg; // declination of observed source in degrees
205 *fLog << "Source: RA=" << outpars.dsourcera_hours << "h DEC=";
206 *fLog << outpars.dsourcedec_deg << "deg" << endl;
207
208 //int inummuonpixels; // number of pixels in the muon shield
209 //int inumcointdcs; // number of coincidence tdcs recorded in the runfile
210 //float fpixdiameter_deg; // smallest pixel diameter (degrees) (from parameters file) */
211
212 // enum axes {RA, DEC, ALT, AZ} ese1_is; // name of the axis to which shaft encoder 1 is attached (implies the type of mount)
213 *fLog << "Shaftencoder 1 @ ";
214 switch (outpars.ese1_is)
215 {
216 case RA: *fLog << "RA"; break;
217 case DEC: *fLog << "DEC"; break;
218 case ALT: *fLog << "ALT"; break;
219 case AZ: *fLog << "AZ"; break;
220 default: *fLog << (int)outpars.ese1_is; break;
221 }
222 *fLog << endl;
223
224 // int isezeropos[2]; // zero position of shaftencoders 1 and 2 (from parameters file)
225 *fLog << "SE Zero: SE(1)=" << outpars.isezeropos[0] << " ";
226 *fLog << "SE(2)=" << outpars.isezeropos[1] << endl;
227
228 // int iaz_rev_track_corr; // correction for the azimuth shaft encoder (ALT/AZ mount only) in reverse tracking mode
229 // int ialt_rev_track_corr; // correction for the altitude shaft encoder (ALT/AZ mount only) in reverse tracking mode
230 *fLog << "Reverse tracking corrections: SE(az)=" << outpars.iaz_rev_track_corr;
231 *fLog << " SE(alt)=" << outpars.ialt_rev_track_corr << endl;
232
233 // float fbendingcorr; // bending correction factor (ALT/AZ mount only)
234 // float fextinction; // atmospheric extinction (typically taken from the Carlsberg Meridian Circle data)
235 *fLog << "Bending: Correction factor=" << outpars.fbendingcorr << " ";
236 *fLog << "Extinction=" << outpars.fextinction << endl;
237
238 // Boolean bdontusepix[iMAXNUMPIX]; // bdontusepix is set true if the pixel should not be used in image analysis, otherwise it is true;
239
240 fBlinds->Clear();
241 *fLog << "Don't use pixels: ";
242 for (int i=0; i<iMAXNUMPIX; i++)
243 if (outpars.bdontusepix[i])
244 {
245 *fLog << i << " ";
246 fBlinds->SetPixelBlind(i);
247 }
248 *fLog << endl;
249
250 *fLog << "Exclude pixels: ";
251 // Boolean bexcludepix[iMAXNUMPIX];
252 for (int i=0; i<iMAXNUMPIX; i++)
253 if (outpars.bexcludepix[i])
254 {
255 *fLog << i << " ";
256 fBlinds->SetPixelBlind(i);
257 }
258 *fLog << endl;
259
260 // save blind pixels for all events of this run
261 fBlnd.Set(iMAXNUMPIX);
262 for (int i=0; i<iMAXNUMPIX; i++)
263 fBlnd[i] = fBlinds->IsBlind(i) ? 1 : 0;
264
265 fBlinds->SetReadyToSave();
266
267 /* bexcludepix[] is set TRUE (== exclude from pedestal, Laser
268 * calibration and the further analysis) when the Mean value
269 * of a pixel inside a pedestal Run is larger than 50 or ( || )
270 * if the pedestal RMS value of this pixel is larger than 5.0
271 * This is reflected in the (new for versions >= 0.4)
272 * variable "pixonoff" in the ntuple written by preproc:
273 * preproc.nt.hbook
274 *
275 * When the pixel is excluded by the user it will get a -2 otherwise
276 * pixonoff = 0.
277 * Additive to this a -1 is added when preproc excludes the pixel
278 * for a given Run. So the actual value tells you whether you caught
279 * it already by hand or not.
280 *
281 * A plot of pixonoff may also be handy to tell you the status of your
282 * ADC equipment. */
283
284 // float fphotoel_per_adccnt[iMAXNUMPIX]; // conversion factors for the pixel signals */
285 /*
286 float padc = outpars.fphotoel_per_adccnt[0];
287 *fLog << "Phe/ADC (pixel 0): " << padc << endl;
288 for (int i=0; i<iMAXNUMPIX; i++)
289 *fLog << outpars.fphotoel_per_adccnt[i] << " ";
290 *fLog << endl;
291 */
292 /*
293 --- USEFULL? NEEDED? ---
294 int irubminusutc_usecs; // difference between rubidium clock and UTC in microseconds
295 int isum_thresh_phot; // threshold for the total sum of photoelectrons filter
296 int i2out_thresh_phot; // threshold for the two-pixels-out-of-all software
297 int imuoncut_thresh_adccnt[iMAXNUMMUONPIX]; // thresholds for the muon cut
298 Boolean bmuon_suppression; // "Dolby" noise reduction flag
299 float ftolerated_pointerror_deg; // maximum tolerated pointing error in the position
300 */
301
302 // float fxpointcorr_deg; // pointing correction (to be added along the camera x axis) e.g. from point run */
303 // float fypointcorr_deg; // pointing correction (to be added along the camera y axis) e.g. from point run */
304 *fLog << "Pointing correction: dx=" << outpars.fxpointcorr_deg << "deg ";
305 *fLog << "dy=" << outpars.fypointcorr_deg << "deg" << endl;
306
307 // FIXME? Is x-y echanged between Mars CT1 geometry and CT1 definition?
308 fSrcPos->SetXY(-outpars.fypointcorr_deg/fGeom->GetConvMm2Deg(),
309 -outpars.fxpointcorr_deg/fGeom->GetConvMm2Deg());
310 fSrcPos->SetReadyToSave();
311
312 /*
313 --- USEFULL? NEEDED? ---
314 float fcamera_align_angle_deg; // the angle between the camera y-axis and the meridian when a culminating object is observed (defined counter-clockwise looking at the sky)
315 int iratecalc_numevents_odd; // number of events used in the rate calculation (must be odd)
316 int inumpedfile; // number of the pedestal file used
317 int inumpedrun; // number of the pedestal run used in the file (starting at 0)
318 int inumcalfile; // number of the calibration file used
319 int inumlaserrun; // number of the laserrun used in the file (starting at 0)
320 int inumtellogfile; // number of the TelLog file to be used
321 int inumtellogrun; // number of the tellog entry (Runnumber) used from the log file
322 int imcparticle; // CORSIKA-coded Monte Carlo particle type.
323 */
324
325 // ----- preprocessing results -----
326
327 // int istart_mjdate_day; // MJD of run start (first event) */
328 // int iend_mjdate_day; // MJD of run end (last event) */
329 // int irunduration_secs; // difference between start and end time (secs) */
330 *fLog << "Run Time: From " << outpars.istart_mjdate_day << " to ";
331 *fLog << outpars.iend_mjdate_day << " (MJD), Duration=";
332 *fLog << outpars.irunduration_secs/3600 << "h";
333 *fLog << (outpars.irunduration_secs/60)%60 << "m";
334 *fLog << outpars.irunduration_secs%60 << "s" << endl;
335
336 /*
337 --- USEFULL? NEEDED? ---
338 int iproc_mjdate; // MJD of data processing (i.e. creation of this file)
339 */
340
341 // int iproc_evts; // number of events processed */
342 *fLog << "Number of processed events: " << outpars.iproc_evts << endl;
343
344 // --- USEFULL? NEEDED? ---
345 // double dactual_sourcera_hours; // for off runs: the false source (that should have been) observed */
346
347 // float frms_pedsig_phot[iMAXNUMPIX]; // standard deviation of the calibrated signals from the pedestal run */
348 fPedest->InitSize(iMAXNUMPIX);
349
350 fPedRMS.Set(iMAXNUMPIX);
351
352 *fLog << "PedestalRMS : ";
353 for (Int_t i=0; i<iMAXNUMPIX; i++)
354 {
355 (*fPedest)[i].SetPedestalRms(outpars.frms_pedsig_phot[i]);
356 *fLog << outpars.frms_pedsig_phot[i] << " ";
357 fPedRMS[i] = outpars.frms_pedsig_phot[i];
358 }
359 *fLog << endl;
360
361 fPedest->SetReadyToSave();
362
363 // Used to communicate the mean over all pixels
364 // pedestal RMS into the Runs NTuple, as it might
365 // be used for e.g. quality cuts.
366 // float fpedrms_mean;
367 *fLog << "Pedestal RMS: " << outpars.fpedrms_mean << endl;
368
369 // The average current over the active pixels
370 // for this run. This is done separately for
371 // ON and OF runs.
372 //float fcurrent_mean;
373
374 // enum eERRORTOLERANCE {CAUTIOUS=0, GOODPHYSICS, TANK} eerrortolerance;
375 /* 0 == "cautious", exits on any reason (but tells in
376 * the .err file,
377 * 1 == "goodphysics", exits when physics could be affected
378 * by the error,
379 * 2 == "tank", never exits except on coredumps and when
380 * all files have been processed. Do not use such files for
381 * physics analysis!
382 *
383 * NOTE: the capital letter words are the enums, the small letter
384 * words must be used inside the parameter file. */
385
386 // enum eMCTRIGGERFLAG {ALL=0, FLAG, NOFLAG} emctriggerflag;
387 /* all: all events which survive the filter are written to the
388 * events NTuple.
389 * flag: When Dorota's triggerflag is set to 1 for a particular
390 * event, it shall be written to the output. All others shall
391 * just be disregarded. (Default)
392 * noflag: Opposite of 'flag': only events with triggerflag = 0 shall
393 * be treated further. */
394
395 *fLog << "Particle Id #" << outpars.imcparticle << endl;
396 *fLog << "Right Ascension: " << outpars.dsourcera_hours << "h" << endl;
397 *fLog << "Declination: " << outpars.dsourcedec_deg << "deg" << endl;
398
399 // Next statement commented out because bmontecarlo was set wrongly
400 //fIsMcFile = outpars.bmontecarlo==TRUE;
401 fIsMcFile = (outpars.dsourcera_hours==0 && outpars.dsourcedec_deg==0 &&
402 outpars.imcparticle != 0);
403
404 if (fIsMcFile != (outpars.bmontecarlo==TRUE))
405 {
406 *fLog << "File tells you that it is a ";
407 *fLog << (outpars.bmontecarlo ? "Monte Carlo" : "Real Data");
408 *fLog << " file." << endl;
409 }
410
411 *fLog << "File detected as a ";
412 *fLog << (fIsMcFile ? "Monte Carlo" : "Real Data");
413 *fLog << " file." << endl;
414 *fLog << " " << endl;
415}
416
417// --------------------------------------------------------------------------
418//
419// Read CT1 PreProc File Header:
420//
421Int_t MCT1ReadPreProc::ReadRunHeader()
422{
423 char cheadertitle[iHEADERTITLELENGTH];
424 fIn->read(cheadertitle, iHEADERTITLELENGTH);
425
426 TString s = cheadertitle;
427 TString m = cTITLE_TEMPLATE;
428
429 if (!s.BeginsWith(m(0, m.First('%'))))
430 return kFALSE;
431
432 *fLog << cheadertitle << flush;
433
434 // cTITLE_TEMPLATE "PREPROC V%f/S%f CT %d RUN %d %d PROCMJD %d\n"
435 struct outputpars outpars;
436
437 int dummy;
438
439 Float_t fpreprocversion, structversion;
440 sscanf(cheadertitle, cTITLE_TEMPLATE,
441 &fpreprocversion, &structversion,
442 &outpars.itelescope, &outpars.irunnum,
443 &dummy/*&outpars.eruntype*/, &outpars.iproc_mjdate);
444
445 if (fpreprocversion<0.6)
446 {
447 *fLog << err << "Sorry, only files from PreProc V0.6 and newer are supported." << endl;
448 return kFALSE;
449 }
450
451 //
452 // This is a stupid way of getting rid of numerical uncertanties when
453 // comparing floating point numbers (Argh...)
454 //
455 TString s1 = Form("%.2f", structversion);
456 TString s2 = Form("%.2f", STRUCT_VERSION);
457
458 if (s1 != s2)
459 {
460 *fLog << warn << "WARNING: Version of C-structures of file (V";
461 *fLog << s1 << ") not identical with current structures (V";
462 *fLog << s2 << ")" << endl;
463 }
464
465 fIn->read((char*)&outpars, sizeof(struct outputpars));
466
467 ProcessRunHeader(outpars);
468
469 //rwagner: ReInit whenever new run commences
470 // rc==-1 means: ReInit didn't work out
471
472 MTaskList *tlist = (MTaskList*)fParList->FindCreateObj("MTaskList");
473 if (!tlist)
474 return -1;
475
476 if (!tlist->ReInit(fParList))
477 return -1;
478
479 return kTRUE;
480}
481
482Int_t MCT1ReadPreProc::ReadRunFooter()
483{
484 char cheadertitle[iHEADERTITLELENGTH];
485 fIn->read(cheadertitle, iHEADERTITLELENGTH);
486 /*
487 sscanf(cheadertitle, cEND_EVENTS_TEMPLATE,
488 &filterres.ifilter_passed_evts);
489 */
490
491 TString s = cheadertitle;
492 TString m = cEND_EVENTS_TEMPLATE;
493 Int_t p = m.First('%');
494
495
496 if (!s.BeginsWith(m(0,p)))
497 {
498 fIn->seekg(-iHEADERTITLELENGTH, ios::cur);
499 return 0;
500 }
501
502 *fLog << inf << cheadertitle << flush;
503
504 struct filterresults filterres;
505 fIn->read((char*)&filterres, sizeof(struct filterresults));
506 /*
507 int imax_alt_arcs; // maximum altitude reached during the run
508 int iaz_at_max_alt_arcs; // azimuth at the time the max. alt. was reached
509 int itimeaverage_alt_arcs; // altitude averaged over the runtime
510 int icoord_inconsist_evts; // number of events with time-coordinate inconsistency in this run
511 int ifilter_passed_evts; // number of events which passed the filter
512 int imuon_fail_evts; // number of events rejected as muons (other filters passed)
513 int i2out_fail_evts; // number of events which failed in the two out of all pixels software trigger
514 int imuon_and_2out_fail_evts; // number of events failing in both muon and 2out filter
515 int isum_fail_evts; // number of events which failed the sum-of-all-calibrated ADC counts filter
516 int isum_and_muon_fail_evts; // number of events which failed in both the sum and the muon filter
517 int isum_and_2out_fail_evts; // number of events which failed in both the sum and the 2out filter
518 int iall_filters_fail_evts; // number of events which failed in all filters
519 float favg_event_rate_hz; // average rate before filtering
520 float fstddev_event_rate_hz; // standard deviation of the rate before filtering
521 */
522
523 if (fNumEventsInRun!=(UInt_t)filterres.ifilter_passed_evts)
524 {
525 *fLog << err << "ERROR! Number of events in run (" << (UInt_t)filterres.ifilter_passed_evts;
526 *fLog << ") doesn't match number of read events (";
527 *fLog << fNumEventsInRun << ")" << endl;
528 *fLog << " File corrupted." << endl;
529 return -1;
530 }
531
532 fNumFilterEvts += fNumEventsInRun;
533 fNumRuns++;
534
535 *fLog << inf << "Read " << fNumEventsInRun << " events from run (Total=";
536 *fLog << fNumFilterEvts << "/" << fEntries << " [";
537 *fLog << 100*fNumFilterEvts/fEntries << "%], Runs=" << fNumRuns << ")";
538 *fLog << endl;
539
540 return 1;
541}
542
543// --------------------------------------------------------------------------
544//
545// This opens the next file in the list and deletes its name from the list.
546//
547Bool_t MCT1ReadPreProc::OpenNextFile()
548{
549 //
550 // open the input stream and check if it is really open (file exists?)
551 //
552 if (fIn)
553 delete fIn;
554 fIn = NULL;
555
556 //
557 // Check for the existence of a next file to read
558 //
559 if (fNumFile >= (UInt_t)fFileNames->GetSize())
560 return kFALSE;
561
562 TNamed *file = (TNamed*)fFileNames->At(fNumFile);
563
564 //TNamed *file = (TNamed*)fFileNames->GetFirst();
565 //if (!file)
566 // return kFALSE;
567
568 //
569 // open the file which is the first one in the chain
570 //
571 const TString name = file->GetName();
572
573 const char *expname = gSystem->ExpandPathName(name);
574 const TString fname(expname);
575 delete [] expname;
576
577 //
578 // Remove this file from the list of pending files
579 //
580 //fFileNames->Remove(file);
581
582 *fLog << inf << "Open file: '" << name << "'" << endl;
583
584 if (!CheckHeader(fname))
585 return kFALSE;
586
587 fNumFile++;
588
589 fIn = new ifstream(fname);
590
591 *fLog << inf << "-----------------------------------------------------------------------" << endl;
592
593 switch (ReadRunHeader())
594 {
595 case kFALSE:
596 *fLog << warn << "Unable to read first run header... skipping file." << endl;
597 return kFALSE;
598 case -1:
599 *fLog << warn << "ReInit of Tasklist didn't succeed." << endl;
600 return kFALSE;
601 default:
602 *fLog << "After opening next file: Number of Events #" << fNumEventsInRun << endl;
603 return kTRUE;
604 }
605}
606
607Bool_t MCT1ReadPreProc::CheckHeader(const TString fname) const
608{
609 ifstream fin(fname);
610 if (!fin)
611 {
612 *fLog << dbginf << err << "ERROR - Cannot open file '" << fname << "'" << endl;
613 return kFALSE;
614 }
615
616 char cheadertitle[iHEADERTITLELENGTH];
617 fin.read(cheadertitle, iHEADERTITLELENGTH);
618
619 Float_t fpreprocversion, structversion;
620 Int_t dummyi;
621
622 sscanf(cheadertitle, cTITLE_TEMPLATE,
623 &fpreprocversion, &structversion,
624 &dummyi, &dummyi, &dummyi, &dummyi);
625
626 if (fpreprocversion < 0.6)
627 {
628 *fLog << dbginf << err << "ERROR - You must use PreProc V0.6 or higher." << endl;
629 return kFALSE;
630 }
631
632 if (STRUCT_VERSION > structversion)
633 {
634 *fLog << warn << "WARNING: Version of C-structures of file (V";
635 *fLog << structversion << ") newer than current structures (V";
636 *fLog << STRUCT_VERSION << ")" << endl;
637 }
638
639 *fLog << "Current structures: " << STRUCT_VERSION << " ";
640 *fLog << "Structures in file: " << structversion << " ";
641 *fLog << "Used preproc version: " << fpreprocversion << endl;
642
643 return kTRUE;
644}
645
646
647Int_t MCT1ReadPreProc::GetNumEvents(const TString fname) const
648{
649 *fLog << inf << "Scanning file " << fname << " for size" << flush;
650
651 ifstream fin(fname);
652 if (!fin)
653 {
654 *fLog << dbginf << err << "ERROR - Opening file." << endl;
655 return 0;
656 }
657
658 const TString m(cEND_EVENTS_TEMPLATE);
659 const Int_t p = m.First('%');
660 const TString test = m(0, p);
661
662 Int_t nevts = 0;
663 Int_t nruns = 0;
664
665 while (!fin.eof() && fin.peek()!=EOF)
666 {
667 fin.seekg(iHEADERTITLELENGTH, ios::cur);
668 fin.seekg(sizeof(struct outputpars), ios::cur);
669
670 while (1)
671 {
672 if (fin.peek()==cEND_EVENTS_TEMPLATE[0])
673 {
674 char cheadertitle[iHEADERTITLELENGTH];
675 fin.read(cheadertitle, iHEADERTITLELENGTH);
676
677 const TString s = cheadertitle;
678 if (s.BeginsWith(test))
679 {
680 fin.seekg(sizeof(struct filterresults), ios::cur);
681 nruns++;
682 break;
683 }
684
685 fin.seekg(-iHEADERTITLELENGTH, ios::cur);
686 }
687
688 fin.seekg(sizeof(struct eventrecord), ios::cur);
689 if (fin.eof())
690 break;
691
692 nevts++;
693 }
694 *fLog << "." << flush;
695 }
696
697 *fLog << "done." << endl;
698 *fLog << "Found " << nevts << " events in " << nruns << " runs." << endl;
699
700 return nevts;
701}
702
703Bool_t MCT1ReadPreProc::Rewind()
704{
705 fNumFilterEvts = 0;
706 fNumEvents = 0;
707 fNumRuns = 0;
708 fNumFile = 0;
709 if (fIn)
710 delete fIn;
711
712 fIn=NULL;
713
714 return kTRUE;
715}
716
717// --------------------------------------------------------------------------
718//
719// Open the first file in the list. Check for the output containers or create
720// them if they don't exist.
721//
722// Initialize the size of the MPedestalCam container to 127 pixels (CT1 camera)
723//
724Int_t MCT1ReadPreProc::PreProcess(MParList *pList)
725{
726 fParList = pList;
727
728 //
729 // look for the HourAngle container in the plist
730 //
731 fHourAngle = (MParameterD*)pList->FindCreateObj("MParameterD", "HourAngle");
732 if (!fHourAngle)
733 return kFALSE;
734 fHourAngle->SetTitle("Store the CT1 hour angle [deg]");
735
736 //
737 // look for the ThetaOrig container in the plist
738 //
739 fThetaOrig = (MParameterD*)pList->FindCreateObj("MParameterD", "ThetaOrig");
740 if (!fThetaOrig)
741 return kFALSE;
742 fThetaOrig->SetTitle("Store the original CT1 zenith angle [rad]");
743
744 //
745 // look for the MCerPhotEvt class in the plist
746 //
747 fNphot = (MCerPhotEvt*)pList->FindCreateObj("MCerPhotEvt");
748 if (!fNphot)
749 return kFALSE;
750
751 //
752 // look for the pedestal class in the plist
753 //
754 fPedest = (MPedestalCam*)pList->FindCreateObj("MPedestalCam");
755 if (!fPedest)
756 return kFALSE;
757
758 //
759 // look for the time class in the plist
760 //
761 fTime = (MTime*)pList->FindCreateObj("MTime");
762 if (!fTime)
763 return kFALSE;
764
765 //
766 // look for the pedestal class in the plist
767 //
768 fBlinds = (MBlindPixels*)pList->FindCreateObj("MBlindPixels");
769 if (!fBlinds)
770 return kFALSE;
771
772 //
773 // look for the source position in the camera
774 //
775 fSrcPos = (MSrcPosCam*)pList->FindCreateObj("MSrcPosCam");
776 if (!fSrcPos)
777 return kFALSE;
778
779 //
780 // look for the camera geometry
781 //
782 fGeom = (MGeomCam*)pList->FindCreateObj("MGeomCamCT1", "MGeomCam");
783 if (!fGeom)
784 return kFALSE;
785
786 //
787 // look for the mc event class
788 //
789 fMcEvt = (MMcEvt*)pList->FindCreateObj("MMcEvt");
790 if (!fMcEvt)
791 return kFALSE;
792
793 //
794 // look for the mc trigger class
795 //
796 fMcTrig = (MMcTrig*)pList->FindCreateObj("MMcTrig");
797 if (!fMcTrig)
798 return kFALSE;
799
800 //
801 // look for the raw run header class
802 //
803 fRawRunHeader = (MRawRunHeader*)pList->FindCreateObj("MRawRunHeader");
804 if (!fRawRunHeader)
805 return kFALSE;
806
807 fBinningT = (MBinning*)pList->FindObject("BinningTheta");
808
809 Rewind();
810
811 fPedest->InitSize(iMAXNUMPIX);
812
813 return GetSelector() ? GetSelector()->CallPreProcess(pList) : kTRUE;
814}
815
816
817// --------------------------------------------------------------------------
818//
819// Smear Theta uniformly in a bin of Theta
820// theta [rad]
821// SmearTheta [rad]
822//
823Double_t MCT1ReadPreProc::SmearTheta(Double_t theta)
824{
825 if (!fBinningT)
826 return theta;
827
828 const Int_t bin = fBinningT->FindLoEdge(theta * 180/TMath::Pi());
829 if (bin<0)
830 return theta;
831
832 // smear Theta within the Theta bin
833 const Double_t low = fBinningT->GetEdges()[bin];
834 const Double_t up = fBinningT->GetEdges()[bin+1];
835
836 // "up-": Do not allow the upper edge
837 return (up - gRandom->Uniform() * (up-low)) * TMath::Pi()/180;
838}
839
840// --------------------------------------------------------------------------
841//
842// Discretize Theta according to the binning in Theta
843// theta [rad]
844// DiscreteTheta [rad] (= bin center)
845//
846Double_t MCT1ReadPreProc::DiscreteTheta(Double_t theta)
847{
848 if (!fBinningT)
849 return theta;
850
851 const Int_t bin = fBinningT->FindLoEdge(theta * 180/TMath::Pi());
852 if (bin<0)
853 return theta;
854
855 const Double_t low = fBinningT->GetEdges()[bin];
856 const Double_t up = fBinningT->GetEdges()[bin+1];
857
858 return 0.5*(up+low) * TMath::Pi()/180;
859}
860
861// --------------------------------------------------------------------------
862//
863// Analyse the event data, means store and copy the needed data into
864// Mars structures and data members
865//
866Bool_t MCT1ReadPreProc::ProcessEvent(const struct eventrecord &event)
867{
868 /*
869 if (fRawRunHeader->GetRunNumber() == 1)
870 {
871 *fLog << "eventrecord" << endl;
872 *fLog << "isecs_since_midday = " << event.isecs_since_midday << endl;
873 *fLog << "isecfrac_200ns = " << event.isecfrac_200ns << endl;
874 *fLog << "snot_ok_flags = " << event.snot_ok_flags << endl;
875 *fLog << "ialt_arcs = " << event.ialt_arcs << endl;
876 *fLog << "iaz_arcs = " << event.iaz_arcs << endl;
877 *fLog << "ipreproc_alt_arcs = " << event.ipreproc_alt_arcs << endl;
878 *fLog << "ipreproc_az_arcs = " << event.ipreproc_az_arcs << endl;
879 *fLog << "ifieldrot_arcs = " << event.ifieldrot_arcs << endl;
880
881 *fLog << "srate_millihz = " << event.srate_millihz << endl;
882 *fLog << "fhourangle = " << event.fhourangle << endl;
883 *fLog << "fmcenergy_tev = " << event.fmcenergy_tev << endl;
884 *fLog << "fmcsize_phel = " << event.fmcsize_phel << endl;
885 *fLog << "imcimpact_m = " << event.imcimpact_m << endl;
886 *fLog << "imcparticle = " << event.imcparticle << endl;
887 *fLog << "imctriggerflag = " << event.imctriggerflag << endl;
888 }
889 */
890
891 // reset blind pixels for this event
892 fBlinds->Clear();
893 for (int i=0; i<iMAXNUMPIX; i++)
894 if ( fBlnd[i]==1 )
895 {
896 fBlinds->SetPixelBlind(i);
897 }
898
899 // reset pedestal RMS for this event
900 for (Int_t i=0; i<iMAXNUMPIX; i++)
901 (*fPedest)[i].SetPedestalRms(fPedRMS[i]);
902
903 // int isecs_since_midday; // seconds passed since midday before sunset (JD of run start)
904 // int isecfrac_200ns; // fractional part of isecs_since_midday
905 fTime->SetTime(event.isecfrac_200ns, event.isecs_since_midday);
906 fTime->SetReadyToSave();
907
908 /*
909 --- USEFULL? NEEDED? ---
910 short snot_ok_flags; // the bits in these two bytes are flags for additional information on the event: Everything OK =: all Bits = 0
911
912 // for ALT-AZ mount telescopes: rotation angle of the field of
913 // view; this angle is defined mathematically positive looking
914 // towards the sky as the angle between the hour circle through
915 // the object being tracked and the line through pixel 1 and 2
916 int ifieldrot_arcs;
917
918 // event rate in milli Hertz before filtering calculated by
919 // iratecalc_numevents_odd/(time[i+iratecalc_numevents_odd/2] -
920 // time[i-iratecalc_numevents_odd/2])
921 // For the first and the last iratecalc_numevents_odd/2
922 // events the rate is assumed to be constant
923 unsigned short srate_millihz;
924
925 // This is the angle between the observation of this event and the
926 // culmination point. It is going to be written into the events NTuple.
927 float fhourangle;
928 */
929
930 //
931 // read in the number of cerenkov photons and add the 'new' pixel
932 // too the list with it's id, number of photons and error
933 // number of photoelectrons measured in each pixel only the
934 // actual number of pixels (outputpars.inumpixels) is written out
935 // short spixsig_10thphot[iMAXNUMPIX];
936 //
937 for (Int_t i=0; i<iMAXNUMPIX; i++)
938 {
939 //*fLog << event.spixsig_10thphot[i] << " ";
940
941 if (event.spixsig_10thphot[i]==0)
942 continue;
943
944 fNphot->AddPixel(i, 0.1*event.spixsig_10thphot[i],
945 (*fPedest)[i].GetPedestalRms());
946 }
947 fNphot->FixSize();
948 fNphot->SetReadyToSave();
949
950 // int ipreproc_alt_arcs; // "should be" alt according to preproc (arcseconds)
951 // int ipreproc_az_arcs; // "should be" az according to preproc (arcseconds)
952
953 // smear Theta in its Theta bin
954 const Double_t theta = TMath::Pi()*(0.5-1./180*event.ialt_arcs/3600);
955 fThetaOrig->SetVal(theta);
956
957 // store hour angle
958 fHourAngle->SetVal(event.fhourangle);
959
960 fMcEvt->Fill(event.isecs_since_midday, //0, /*fEvtNum*/
961 fIsMcFile ? event.imcparticle : 0, /*corsika particle type*/
962 fIsMcFile ? event.fmcenergy_tev*1000 : 0,
963 0, /* fThi0 */
964 0, /* fFirTar */
965 0, /* fzFirInt */
966 // 0, /* fThet*/
967 // rwagner: The following should be theta, right? Check with
968 // altitude fill some lines down...
969 0, // altitude (arcseconds)
970 0, /* fPhii */
971 0, /* fCorD */
972 0, /* fCorX */
973 0, /* fCorY */
974 fIsMcFile ? event.imcimpact_m*100 : 0,
975 TMath::Pi()/180*event.iaz_arcs/3600, // azimuth (arcseconds)
976 // fIsMcFile ? SmearTheta(theta) : theta,
977 DiscreteTheta(theta),
978 0, /* fTFirst */
979 0, /* fTLast */
980 0, /* fL_Nmax */
981 0, /* fL_t0 */
982 0, /* fL_tmax */
983 0, /* fL_a */
984 0, /* fL_b */
985 0, /* fL_c */
986 0, /* fL_chi2 */
987 0, /* uiPin */
988 0, /* uiPat */
989 0, /* uiPre */
990 0, /* uiPco */
991 0, /* uiPelS */
992 (int)(fIsMcFile ? event.fmcsize_phel : 0), /* uiPelC, Simulated SIZE */
993 0, /* elec */
994 0, /* muon */
995 0 /* other */
996 );
997
998 fMcTrig->SetFirstLevel(event.imctriggerflag); // MC data from Dorota get a triggerflag: 1 means triggered, 0 not. */
999
1000 fMcTrig->SetReadyToSave();
1001 fMcEvt->SetReadyToSave();
1002
1003 return kTRUE;
1004}
1005
1006// --------------------------------------------------------------------------
1007//
1008// Because of the file format of the preproc output we have to check at any
1009// event where in the file stream we are...
1010//
1011Bool_t MCT1ReadPreProc::CheckFilePosition()
1012{
1013 //
1014 // Means: If no file is open (first call) try to open the first file
1015 //
1016 if (!fIn)
1017 return kFALSE;
1018
1019 //
1020 // Because we can have 0-event runs in the file we loop as often
1021 // as we don't find a new footer-header combination.
1022 //
1023 while (1)
1024 {
1025 //
1026 // If the first character isn't the first of the footer it must be
1027 // an event
1028 //
1029 if (fIn->peek()!=cEND_EVENTS_TEMPLATE[0])
1030 return kTRUE;
1031
1032 //
1033 // Try reading the footer. If this isn't successful...
1034 // must be an event
1035 //
1036 switch (ReadRunFooter())
1037 {
1038 case -1:
1039 return kFALSE;
1040 case 0:
1041 return kTRUE;
1042 }
1043
1044 *fLog << inf << "Footer found." << endl;
1045
1046 const char c = fIn->peek();
1047
1048 //
1049 // No after reading the footer check if we reached the end of the file
1050 //
1051 if (fIn->eof() || c==EOF)
1052 {
1053 *fLog << "End of file." << endl;
1054 return kFALSE;
1055 }
1056
1057 //
1058 // If the eof isn't reached a new header must follow. Check for it.
1059 //
1060 if (c!=cTITLE_TEMPLATE[0])
1061 {
1062 *fLog << inf << "Error finding new run header in file (possible EOF)... skipping rest of file." << endl;
1063 return kFALSE;
1064 }
1065
1066 *fLog << "-----------------------------------------------------------------------" << endl;
1067
1068
1069 if (ReadRunHeader() < 0)
1070 {
1071 *fLog << warn << "ReInit of Tasklist didn't succeed." << endl;
1072 return kFALSE;
1073 }
1074 }
1075}
1076
1077// --------------------------------------------------------------------------
1078//
1079// Check for the event number and depending on this number decide if
1080// pedestals or event data has to be read.
1081//
1082// If the end of the file is reached try to open the next in the list. If
1083// there is now next file stop the eventloop.
1084//
1085Int_t MCT1ReadPreProc::Process()
1086{
1087 //
1088 // Check where in the file we are. If neither a new event, nor a
1089 // footer/header combination is detected go to the next file.
1090 //
1091 if (!CheckFilePosition())
1092 if (!OpenNextFile())
1093 return kFALSE;
1094
1095 //
1096 // Check for a selector. If one is given and returns kFALSE
1097 // skip this event.
1098 //
1099 if (GetSelector())
1100 {
1101 //
1102 // Make sure selector is processed
1103 //
1104 if (!GetSelector()->CallProcess())
1105 {
1106 *fLog << err << dbginf << "Processing Selector failed." << endl;
1107 return kFALSE;
1108 }
1109
1110 //
1111 // Skip Event
1112 //
1113 if (!GetSelector()->IsConditionTrue())
1114 {
1115 fIn->seekg(sizeof(struct eventrecord), ios::cur);
1116
1117 fNumEvents++;
1118 fNumEventsInRun++;
1119
1120 return kCONTINUE;
1121 }
1122 }
1123
1124 // event data to be read from the file
1125 struct eventrecord event;
1126
1127
1128
1129 // read the eventrecord from the file
1130 fIn->read((char*)&event, sizeof(struct eventrecord));
1131
1132 switch (ProcessEvent(event))
1133 {
1134 case kFALSE:
1135 return kFALSE;
1136 case kCONTINUE:
1137 return kCONTINUE;
1138 }
1139
1140 fNumEvents++;
1141 fNumEventsInRun++;
1142
1143 return kTRUE;
1144}
1145
1146Int_t MCT1ReadPreProc::PostProcess()
1147{
1148 *fLog << all;
1149 *fLog << "Number events passed the filter: " << fNumFilterEvts << endl;
1150 *fLog << "Number of Events read from file: " << fNumEvents << endl;
1151 *fLog << "Number of Runs read from file: " << fNumRuns << endl;
1152 *fLog << "Number of events detected first: " << fEntries << endl;
1153
1154 if (fNumEvents!=fNumFilterEvts)
1155 {
1156 *fLog << warn << "WARNING! Number of events in file doesn't match number of read events..." << endl;
1157 *fLog << " File might be corrupt." << endl;
1158 }
1159
1160 delete fIn;
1161 fIn = NULL;
1162
1163 return GetSelector() ? GetSelector()->CallPostProcess() : kTRUE;
1164}
Note: See TracBrowser for help on using the repository browser.