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

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