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

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