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

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