source: trunk/MagicSoft/Simulation/Detector/Camera/camera.cxx@ 859

Last change on this file since 859 was 781, checked in by blanch, 24 years ago
New numbering for branches from different trigger conditions has been implemented. Now, they are calles: ClassName;1., ClasNema;2., ... The MontaCarlo Headers (MMcTrigHeader and MMcFadcHeader) have been move to the RunHeaders tree. Also the MRawRunHeader is thera with some of its information already filled.
File size: 93.8 KB
Line 
1//!/////////////////////////////////////////////////////////////////////
2//
3// camera
4//
5// @file camera.cxx
6// @title Camera simulation
7// @subtitle Code for the simulation of the camera phase
8// @desc Code for the simulation of the camera of CT1 and MAGIC
9// @author J C Gonzalez
10// @email gonzalez@mppmu.mpg.de
11// @date Thu May 7 16:24:22 1998
12//
13//----------------------------------------------------------------------
14//
15// Created: Thu May 7 16:24:22 1998
16// Author: Jose Carlos Gonzalez
17// Purpose: Program for reflector simulation
18// Notes: See files README for details
19//
20//----------------------------------------------------------------------
21//
22// $RCSfile: camera.cxx,v $
23// $Revision: 1.25 $
24// $Author: blanch $
25// $Date: 2001-05-08 08:07:54 $
26//
27////////////////////////////////////////////////////////////////////////
28// @tableofcontents @coverpage
29
30//=-----------------------------------------------------------
31//!@section Source code of |camera.cxx|.
32
33/*!@"
34
35 In this section we show the (commented) code of the program for the
36 read-out of the output files generated by the simulator of the
37 reflector, |reflector 0.3|.
38
39 @"*/
40
41//=-----------------------------------------------------------
42//!@subsection Includes and Global variables definition.
43
44//!@{
45
46// includes for ROOT
47// BEWARE: the order matters!
48
49#include "TROOT.h"
50
51#include "TRandom.h"
52#include "TApplication.h"
53
54#include "TFile.h"
55#include "TTree.h"
56#include "TBranch.h"
57#include "TCanvas.h"
58
59#include "TArrayC.h"
60
61#include "MTrigger.hxx"
62#include "MFadc.hxx"
63#include "MLons.hxx"
64
65#include "MRawRunHeader.h"
66#include "MRawEvtData.h"
67#include "MRawEvtHeader.h"
68#include "MRawCrateData.h"
69#include "MMcEvt.hxx"
70#include "MMcTrig.hxx"
71#include "MMcTrigHeader.hxx"
72#include "MMcFadcHeader.hxx"
73
74/*!@"
75
76 All the defines are located in the file |camera.h|.
77
78 @"*/
79
80#include "camera.h"
81
82//!@}
83
84/*!@"
85
86 The following set of flags are used in time of compilation. They do
87 not affect directly the behaviour of the program at run-time
88 (though, of course, if you disconnected the option for
89 implementation of the Trigger logic, you will not be able to use any
90 trigger at all. The 'default' values mean default in the sense of
91 what you got from the server when you obtained this program.
92
93 @"*/
94
95//!@{
96
97// flag for debugging (default: OFF )
98#define __DEBUG__
99#undef __DEBUG__
100
101
102//!@}
103
104//=-----------------------------------------------------------
105//!@subsection Definition of global variables.
106
107/*!@"
108
109 Now we define some global variables with data about the telescope,
110 such as "focal distance", number of pixels/mirrors,
111 "size of the camera", and so on.
112
113 @"*/
114
115/*!@"
116
117 Depending on the telescope we are using (CT1 or MAGIC), the
118 information stored in the definition file is different.
119 The variable |ct_Type| has the value 0 when we use
120 CT1, and 1 when we use MAGIC.
121
122 @"*/
123
124//!@{
125static int ct_Type; //@< Type of telescope: 0:CT1, 1:MAGIC
126//!@}
127
128/*!@"
129
130 And this is the information about the whole telescope.
131
132 @"*/
133
134//!@{
135
136// parameters of the CT (from the CT definition file)
137
138////@: Focal distances [cm]
139//static float *ct_Focal;
140
141//@: Mean Focal distances [cm]
142static float ct_Focal_mean;
143
144//@: STDev. Focal distances [cm]
145static float ct_Focal_std;
146
147//@: Mean Point Spread function [cm]
148static float ct_PSpread_mean;
149
150//@: STDev. Point Spread function [cm]
151static float ct_PSpread_std;
152
153//@: STDev. Adjustmente deviation [cm]
154static float ct_Adjustment_std;
155
156//@: Radius of the Black Spot in mirror [cm]
157static float ct_BlackSpot_rad;
158
159//@: Radius of one mirror [cm]
160static float ct_RMirror;
161
162//@: Camera width [cm]
163static float ct_CameraWidth;
164
165//@: Pixel width [cm]
166static float ct_PixelWidth;
167
168//@: ct_PixelWidth_corner_2_corner = ct_PixelWidth / cos(60)
169static float ct_PixelWidth_corner_2_corner;
170
171//@: ct_PixelWidth_corner_2_corner / 2
172static float ct_PixelWidth_corner_2_corner_half;
173
174//@: Number of mirrors
175static int ct_NMirrors = 0;
176
177//@: Number of pixels
178static int ct_NPixels;
179
180//@: Number of pixels
181static int ct_NCentralPixels;
182
183//@: Number of pixels
184static int ct_NGapPixels;
185
186//@: ct_Apot = ct_PixelWidth / 2
187static float ct_Apot;
188
189//@: ct_2Apot = 2 * ct_Apot = ct_PixelWidth
190static float ct_2Apot;
191
192//@: name of the CT definition file to use
193static char ct_filename[256];
194
195//@: list of showers to be skipped
196static int *Skip;
197
198//@: number of showers to be skipped
199static int nSkip=0;
200
201//@: flag: TRUE: data come from STDIN; FALSE: from file
202static int Data_From_STDIN = FALSE;
203
204//@: flag: TRUE: write all images to output; FALSE: only triggered showers
205static int Write_All_Images = FALSE;
206
207//@: flag: TRUE: write all data to output; FALSE: only triggered showers
208static int Write_All_Data = FALSE;
209
210static int Write_McEvt = TRUE;
211static int Write_McTrig = TRUE;
212static int Write_McFADC = TRUE;
213static int Write_RawEvt = FALSE;
214
215//@: flag: TRUE: selection on the energy
216static int Select_Energy = TRUE;
217
218//@: Lower edge of the selected energy range (in GeV)
219static float Select_Energy_le = 0.0;
220
221//@: Upper edge of the selected energy range (in GeV)
222static float Select_Energy_ue = 100000.0;
223
224//@: flag: TRUE: show all fadc singnal in the screen; FALSE: don't
225static int FADC_Scan = FALSE;
226
227//@: flag: TRUE: show all trigger singnal in the screen; FALSE: don't
228static int Trigger_Scan = FALSE;
229
230//@: flag: TRUE: loop trigger analysis over several thresholds, multiplicities and topologies; FALSE: a single trigger configuration
231static int Trigger_Loop = FALSE;
232
233//@: Properties of the trigger
234static float Trigger_gate_length = 6.0;
235static float Trigger_response_ampl = 1.0;
236static float Trigger_response_fwhm = 2.0;
237static float Trigger_overlaping_time= 0.25;
238
239//@: Properties of the FADC
240static float FADC_response_ampl = MFADC_RESPONSE_AMPLITUDE;
241static float FADC_response_fwhm = MFADC_RESPONSE_FWHM;
242
243//@: Trigger conditions for a single trigger mode
244static float Trigger_threshold = 7.0;
245static int Trigger_multiplicity = 4;
246static int Trigger_topology = 2;
247
248//@: Upper and lower edges of the trigger loop
249static int Trigger_loop_lthres = 0;
250static int Trigger_loop_uthres = 10;
251static int Trigger_loop_lmult = 2;
252static int Trigger_loop_umult = 10;
253static int Trigger_loop_ltop = 0;
254static int Trigger_loop_utop = 2;
255
256//!@}
257
258/*!@"
259
260 The following double-pointer is a 2-dimensional table with information
261 about each pixel. The routine read_pixels will generate
262 the information for filling it using igen_pixel_coordinates().
263
264 @"*/
265
266//!@{
267// Pointer to a tables/Arrays with information about the pixels
268// and data stored on them with information about the pixels
269
270
271//@: coordinates x,y for each pixel
272static float **pixary;
273
274//@: indexes of pixels neighbours of a given one
275static int **pixneig;
276
277//@: number of neighbours a pixel have
278static int *npixneig;
279
280//@: contents of the pixels (ph.e.)
281static float *fnpix;
282
283
284//!@}
285
286/*!@"
287
288 The following double-pointer is a 2-dimensional table with the
289 Quantum Efficiency @$QE@$ of each pixel in the camera, as a function
290 of the wavelength @$\lambda@$. The routine |read_pixels()| will read
291 also this information from the file |qe.dat|.
292
293 @"*/
294
295//!@{
296// Pointer to a table with QE, number of datapoints, and wavelengths
297
298//@: table of QE
299static float ***QE;
300
301//@: number of datapoints for the QE curve
302static int pointsQE;
303
304//@: table of QE
305static float *QElambda;
306//!@}
307
308/*!@"
309
310 The following double-pointer is a 2-dimensional table with information
311 about each mirror in the dish. The routine |read_ct_file()| will read
312 this information from the CT definition file.
313
314 @"*/
315
316//!@{
317// Pointer to a table with the following info.:
318
319static float **ct_data;
320
321/*
322 * TYPE=0 (CT1)
323 * i s rho theta x y z thetan phin xn yn zn
324 *
325 * i : number of the mirror
326 * s : arc length [cm]
327 * rho : polar rho of the position of the center of the mirror [cm]
328 * theta : polar angle of the position of the center of the mirror [cm]
329 * x : x coordinate of the center of the mirror [cm]
330 * y : y coordinate of the center of the mirror [cm]
331 * z : z coordinate of the center of the mirror [cm]
332 * thetan : polar theta angle of the direction where the mirror points to
333 * phin : polar phi angle of the direction where the mirror points to
334 * xn : xn coordinate of the normal vector in the center (normalized)
335 * yn : yn coordinate of the normal vector in the center (normalized)
336 * zn : zn coordinate of the normal vector in the center (normalized)
337 *
338 * TYPE=1 (MAGIC)
339 * i f sx sy x y z thetan phin
340 *
341 * i : number of the mirror
342 * f : focal distance of that mirror
343 * sx : curvilinear coordinate of mirror's center in X[cm]
344 * sy : curvilinear coordinate of mirror's center in X[cm]
345 * x : x coordinate of the center of the mirror [cm]
346 * y : y coordinate of the center of the mirror [cm]
347 * z : z coordinate of the center of the mirror [cm]
348 * thetan : polar theta angle of the direction where the mirror points to
349 * phin : polar phi angle of the direction where the mirror points to
350 * xn : xn coordinate of the normal vector in the center (normalized)
351 * yn : yn coordinate of the normal vector in the center (normalized)
352 * zn : zn coordinate of the normal vector in the center (normalized)
353 */
354//!@}
355
356/*!@"
357
358 We define a table into where random numbers will be stored.
359 The routines used for random number generation are provided by
360 |RANLIB| (taken from NETLIB, |www.netlib.org|), and by
361 the routine |double drand48(void)| (prototype defined in
362 |stdlib.h|) through the macro |RandomNumber| defined in
363 |camera.h|.
364
365 @"*/
366
367
368//!@}
369
370extern char FileName[];
371
372
373//=-----------------------------------------------------------
374// @subsection Main program.
375
376//!@{
377
378//++++++++++++++++++++++++++++++++++++++++
379// MAIN PROGRAM
380//----------------------------------------
381
382int main(int argc, char **argv)
383{
384
385 //!@' @#### Definition of variables.
386 //@'
387
388 int ioscar=0;
389
390 char inname[256]; //@< input file name
391 char starfieldname[256]; //@< starfield input file name
392
393 char datname[256]; //@< data (ASCII) output file name
394
395 char rootname[256] ; //@< ROOT file name
396 char rootname_loop[256] ; //@< ROOT file name
397
398 char parname[256]; //@< parameters file name
399
400 char nsbpathname[256]; //@< directory with the dataabse for th NSB
401
402 char flag[SIZE_OF_FLAGS + 1]; //@< flags in the .rfl file
403
404 FILE *inputfile; //@< stream for the input file
405 ofstream datafile; //@< stream for the data file
406
407 MCEventHeader mcevth; //@< Event Header class (MC)
408 MCCphoton cphoton; //@< Cherenkov Photon class (MC)
409
410 int inumphe; //@< number of photoelectrons in an event from showers
411 float inumphensb; //@< number of photoelectrons in an event from nsb
412
413 float arrtmin_ns; //@ arrival time of the first photoelectron
414 float arrtmax_ns; //@ arrival time of the last photoelectron
415
416 float thetaCT, phiCT; //@< parameters of a given shower
417 float thetashw, phishw; //@< parameters of a given shower
418 float coreD, coreX, coreY; //@< core position
419 float impactD; //@< impact parameter
420 float l1, m1, n1; //@< auxiliary variables
421 float l2, m2, n2; //@< auxiliary variables
422 float num, den; //@< auxiliary variables
423
424 int nshow=0; //@< partial number of shower in a given run
425 int ntshow=0; //@< total number of showers
426 int ncph=0; //@< partial number of photons in a given run
427 int ntcph=0; //@< total number of photons
428
429 int i, j, k; //@< simple counters
430
431 int simulateNSB; //@< Will we simulate NSB?
432 int nphe2NSB; //@< From how many ohe will we simulate NSB?
433 float meanNSB; //@< diffuse NSB mean value (phe per ns per central pixel)
434 float diffnsb_phepns[iMAXNUMPIX]; //@< diffuse NSB values for each pixel
435 //@< derived from meanNSB
436 float nsbrate_phepns[iMAXNUMPIX][iNUMWAVEBANDS]; //@< non-diffuse nsb
437 //@< photoelectron rates
438 float nsb_phepns[iMAXNUMPIX]; //@< NSB values for each pixel
439
440 Float_t nsb_trigresp[TRIGGER_TIME_SLICES]; //@< array to write the trigger
441 //@< response from the database
442 Float_t nsb_fadcresp[(Int_t) SLICES_MFADC]; //@< array to write the fadc
443 //@< response from the database
444 Byte_t trigger_map[((Int_t)(TRIGGER_PIXELS/8))+1];//@< Pixels on when the
445 //@< camera triggers
446 Float_t fadc_elecnoise[CAMERA_PIXELS]; //@< Electronic niose for each pixel
447 Float_t fadc_pedestals[CAMERA_PIXELS]; //@< array for fadc pedestals values
448
449 float ext[iNUMWAVEBANDS] = { //@< average atmospheric extinction in each waveband
450 EXTWAVEBAND1,
451 EXTWAVEBAND2,
452 EXTWAVEBAND3,
453 EXTWAVEBAND4,
454 EXTWAVEBAND5
455 };
456
457 float qThreshold; //@< Threshold value
458 float qTailCut; //@< Tail Cut value
459 int nIslandsCut; //@< Islands Cut value
460 int countIslands; //@< Will we count the islands?
461 int anaPixels;
462
463 float fCorrection; //@< Factor to apply to pixel values (def. 1.)
464
465 int ntrigger = 0; //@< number of triggers in the whole file
466 int btrigger = 0; //@< trigger flag
467 int ithrescount; //@< counter for loop over threshold trigger
468 int imulticount; //@< counter for loop over multiplicity trigger
469 int itopocount; //@< counter for loop over topology trigger
470 int icontrigger; //@< number of trigger conditions to be analised
471 UShort_t numPix; //@< number of sets of fadc written counts
472
473 float fpixelthres[TRIGGER_PIXELS];
474
475 TArrayC *fadcValues; //@< the analog Fadc siganl for pixels
476
477 float plateScale_cm2deg; //@< plate scale (deg/cm)
478 float degTriggerZone; //@< trigger area in the camera (radius, in deg.)
479
480 int still_in_loop = FALSE;
481
482 struct camera cam; // structure holding the camera definition
483
484
485 //!@' @#### Definition of variables for |getopt()|.
486 //@'
487
488 int ch, errflg = 0; //@< used by getopt
489
490 //!@' @#### Initialisation of some variables
491 //@'
492
493 for(i=0;i<iMAXNUMPIX;i++)
494 for(j=0;j<iNUMWAVEBANDS;j++)
495 nsbrate_phepns[i][j]=0.0; //@< Starlight rates initialised at 0.0
496
497 /*!@'
498
499 @#### Beginning of the program.
500
501 We start with the main program. First we (could) make some
502 presentation, and follows the reading of the parameters file (now
503 from the |stdin|), the reading of the CT parameters file, and the
504 creation of the output file, where the processed data will be
505 stored.
506
507 */
508
509 //++
510 // START
511 //--
512
513 // make unbuffered output
514
515 cout.setf ( ios::stdio );
516
517 // parse command line options (see reflector.h)
518
519 parname[0] = '\0';
520
521 optarg = NULL;
522 while ( !errflg && ((ch = getopt(argc, argv, COMMAND_LINE_OPTIONS)) != -1) )
523 switch (ch) {
524 case 'f':
525 strcpy(parname, optarg);
526 break;
527 case 'h':
528 usage();
529 break;
530 default :
531 errflg++;
532 }
533
534 // show help if error
535
536 if ( errflg>0 )
537 usage();
538
539 // make some sort of presentation
540
541 present();
542
543 // read parameters file
544
545 if ( strlen(parname) < 1 )
546 readparam(NULL);
547 else
548 readparam(parname);
549
550 // read data from file or from STDIN?
551
552 Data_From_STDIN = get_data_from_stdin();
553
554 // write all images, even those without trigger?
555
556 Write_All_Images = get_write_all_images();
557
558 // write all data (i.e., ph.e.s in pixels)
559
560 Write_All_Data = get_write_all_data();
561
562 Write_McEvt = get_write_McEvt() ;
563 Write_McTrig = get_write_McTrig() ;
564 Write_McFADC = get_write_McFadc() ;
565 Write_RawEvt = get_write_RawEvt() ;
566
567 FADC_Scan = get_FADC_Scan();
568 Trigger_Scan = get_Trigger_Scan();
569
570 get_FADC_properties( &FADC_response_ampl, &FADC_response_fwhm);
571
572 get_Trigger_properties( &Trigger_gate_length, &Trigger_overlaping_time, &Trigger_response_ampl, &Trigger_response_fwhm);
573
574 Trigger_Loop = get_Trigger_Loop(&Trigger_loop_lthres, &Trigger_loop_uthres, &Trigger_loop_lmult, &Trigger_loop_umult, &Trigger_loop_ltop, &Trigger_loop_utop);
575
576 icontrigger =(Trigger_loop_uthres-Trigger_loop_lthres+1)*
577 (Trigger_loop_umult-Trigger_loop_lmult+1)*
578 (Trigger_loop_utop-Trigger_loop_ltop+1);
579
580 if (!Trigger_Loop){
581 get_Trigger_Single (&Trigger_threshold, &Trigger_multiplicity, &Trigger_topology);
582 icontrigger=1;
583 }
584
585 // get filenames
586
587 strcpy( inname, get_input_filename() );
588 strcpy( starfieldname, get_starfield_filename() );
589 strcpy( datname, get_data_filename() );
590 strcpy( rootname, get_root_filename() );
591 strcpy( rootname_loop, get_loop_filename() );
592 strcpy( ct_filename, get_ct_filename() );
593 strcpy( nsbpathname, get_nsb_directory() );
594
595 // get different parameters of the simulation
596
597 qThreshold = get_threshold();
598 qTailCut = get_tail_cut();
599 simulateNSB = get_nsb( &meanNSB, &nphe2NSB );
600 countIslands = get_islands_cut( &nIslandsCut );
601
602 // get selections on the parameters
603
604 Select_Energy = get_select_energy( &Select_Energy_le, &Select_Energy_ue);
605
606 // log filenames information
607
608 log(SIGNATURE,
609 "%s:\n\t%20s:\t%s\n\t%20s:\t%s\n\t%20s:\t%s\n\t%20s:\t%s\n\t%20s:\t%s\n\t%20s:\t%s\n",
610 "Filenames",
611 "In", inname,
612 "Stars", starfieldname,
613 "NSB database", nsbpathname,
614 "CT", ct_filename,
615 "Data", datname,
616 "ROOT", rootname
617 );
618
619 // log Trigger information
620
621 if (Trigger_Loop) {
622 log(SIGNATURE,
623 "%s:\n\t%20s: %i - %i\n\t%20s: %i - %i\n\t%20s: %i - %i\n\t%20s\n",
624 "Trigger Loop mode",
625 "Threshold",Trigger_loop_lthres,Trigger_loop_uthres,
626 "Multiplicity",Trigger_loop_lmult,Trigger_loop_umult,
627 "Topology",Trigger_loop_ltop,Trigger_loop_utop,
628 rootname_loop);
629 }
630 else{
631 log(SIGNATURE,
632 "%s:\n\t%20s: %f\n\t%20s: %i\n\t%20s: %i\n",
633 "Single Trigger mode",
634 "Threshold",Trigger_threshold,
635 "Multiplicity",Trigger_multiplicity,
636 "Topology",Trigger_topology);
637 }
638
639 // log flags information
640
641 log(SIGNATURE,
642 "%s:\n\t%20s: %s\n\t%20s: %s\n\t%20s: %s\n",
643 "Flags",
644 "Data_From_STDIN", ONoff(Data_From_STDIN),
645 "Write_All_Images", ONoff(Write_All_Images),
646 "Write_All_Data", ONoff(Write_All_Data));
647
648 // log flags information
649
650 log(SIGNATURE,
651 "%s:\n\t%20s: %s\n\t%20s: %s\n\t%20s: %s\n",
652 "Root output",
653 "Write_McEvt", ONoff(Write_McEvt),
654 "Write_McTrig", ONoff(Write_McTrig),
655 "Write_McFADC", ONoff(Write_McFADC),
656 "Write_RawEvt", ONoff(Write_RawEvt));
657
658 // log parameters information
659
660 log(SIGNATURE,
661 "%s:\n\t%20s: %f\n\t%20s: %f\n\t%20s: %f %s\n\t%20s: %f %s\n",
662 "Parameters",
663 "q0 (Threshold)", qThreshold,
664 "t0 (Tail-cut)", qTailCut,
665 "NSB (phes/pixel)", meanNSB, ONoff(simulateNSB),
666 "i0 (Islands-cut)", nIslandsCut, ONoff(countIslands));
667
668 // log selections
669
670 log(SIGNATURE,
671 "%s:\n\t%20s: %s (%f:%f)\n",
672 "Selections:",
673 "Energy", ONoff(Select_Energy), Select_Energy_le, Select_Energy_ue);
674
675 // Definition and initialization of array to save trigger statistics
676
677 int ***ntriggerloop;
678
679 ntriggerloop= new int ** [(int) (Trigger_loop_uthres+1-Trigger_loop_lthres)];
680 for (ithrescount=0;ithrescount<=Trigger_loop_uthres-Trigger_loop_lthres;ithrescount++){
681 ntriggerloop[ithrescount]= new int * [Trigger_loop_umult-Trigger_loop_lmult+1];
682 for (imulticount=0;imulticount<=Trigger_loop_umult-Trigger_loop_lmult;imulticount++){
683 ntriggerloop[ithrescount][imulticount]= new int [Trigger_loop_utop-Trigger_loop_ltop+1];
684 for(itopocount=0;itopocount<=Trigger_loop_utop-Trigger_loop_ltop;itopocount++){
685 ntriggerloop[ithrescount][imulticount][itopocount]=0;
686 }
687 }
688 }
689
690 // set all random numbers seeds
691
692 setall( get_seeds(0), get_seeds(1) );
693
694 // get list of showers to evt. skip
695
696 nSkip = get_nskip_showers();
697
698 if (nSkip > 0) {
699 Skip = new int[ nSkip ];
700 get_skip_showers( Skip );
701
702 log(SIGNATURE, "There are some showers to skip:\n");
703 for (i=0; i<nSkip; ++i)
704 log(SIGNATURE, "\tshower # %d\n", Skip[i]);
705 }
706
707
708 // read parameters from the ct.def file
709
710 read_ct_file();
711
712 // read camera setup
713
714 read_pixels(&cam);
715
716 Int_t Lev0, Lev1, Lev2 ;
717
718 fadcValues = new TArrayC(FADC_SLICES);
719
720 // number of pixels for parameters
721
722 anaPixels = get_ana_pixels();
723 anaPixels = (anaPixels == -1) ? ct_NPixels : anaPixels;
724
725 // initialise ROOT
726
727 TROOT simple("simple", "MAGIC Telescope Monte Carlo");
728
729 // initialise instance of Trigger and FADC classes
730
731 MTrigger Trigger(Trigger_gate_length, Trigger_overlaping_time,
732 Trigger_response_ampl, Trigger_response_fwhm); //@< A instance of the Class MTrigger
733
734 // Initialise McTrig information class if we want to save trigger informtion
735
736 MMcTrig **McTrig;
737 MMcTrigHeader **HeaderTrig;
738 MMcFadcHeader **HeaderFadc;
739
740 if (Write_McTrig){
741
742 McTrig = new MMcTrig * [icontrigger];
743
744 for (i=0;i<icontrigger;i++) {
745 McTrig[i] = new MMcTrig();
746 }
747
748 HeaderTrig = new MMcTrigHeader * [icontrigger];
749
750 for (i=0;i<icontrigger;i++) {
751 HeaderTrig[i] = new MMcTrigHeader();
752 }
753 }
754
755 if (Write_McFADC){
756
757 HeaderFadc = new MMcFadcHeader * [icontrigger];
758
759 for (i=0;i<icontrigger;i++) {
760 HeaderFadc[i] = new MMcFadcHeader();
761 }
762 }
763
764 MFadc fadc(FADC_response_ampl,FADC_response_fwhm) ; //@< A instance of the Class MFadc
765
766 //!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
767 //
768 // Set the FADC pedestals for that run
769 // Some modifications
770 // mut be done to simulate a more realistic distribution of the pedestals.
771 // This simualtion is done int the SetPedestals methode inside the
772 // class MFadc
773 // Currentlly a input_pedestal array is declared with the pedestals.
774 // Thy can also be set randomly following a flat distribution.
775 //
776 /////////////////////////////////////////////////////////////////////
777
778 Float_t input_pedestals[CAMERA_PIXELS];
779
780 for(i=0;i<CAMERA_PIXELS;i++)
781 input_pedestals[i]=10.5;
782
783 fadc.SetPedestals(input_pedestals);
784
785 // Prepare the raw data output
786
787 // Header Tree
788
789 MRawRunHeader **RunHeader;
790
791 RunHeader = new MRawRunHeader * [1];
792 RunHeader[0] = new MRawRunHeader();
793
794 // Header branch
795
796 MRawEvtHeader **EvtHeader;
797
798 if (Write_RawEvt) {
799 EvtHeader = new MRawEvtHeader * [icontrigger];
800
801 for (i=0;i<icontrigger;i++) {
802 EvtHeader[i] = new MRawEvtHeader();
803 }
804 }
805
806 // Data branch
807
808 MRawEvtData **EvtData; // Data branch
809
810 if (Write_RawEvt) {
811 EvtData = new MRawEvtData * [icontrigger];
812
813 for (i=0;i<icontrigger;i++) {
814 EvtData[i] = new MRawEvtData();
815 }
816 }
817
818 MMcEvt *McEvt = new MMcEvt ();
819
820 //
821 // initalize a temporal ROOT file
822 //
823
824 TFile outfile_temp ( rootname , "RECREATE" );
825
826 Int_t bsize=128000; Int_t split=1;
827
828 // create a Tree for the Header Event
829 TTree HeaderTree("RunHeaders","Headers of Run");
830
831 // define branches of Header Tree
832
833 char help[4];
834
835 HeaderTree.Branch("MRawRunHeader","MRawRunHeader",
836 &RunHeader[0], bsize, split);
837
838 if(!Trigger_Loop && Write_McTrig){
839
840 HeaderTree.Branch("MMcTrigHeader","MMcTrigHeader",
841 &HeaderTrig[0], bsize, split);
842 }
843 if (Trigger_Loop && Write_McTrig){
844 for(char branchname[20],i=0;i<icontrigger;i++){
845
846 sprintf(help,"%i",i+1);
847 strcpy (branchname, "MMcTrigHeader;");
848 strcat (branchname, & help[0]);
849 strcat (branchname, ".");
850 HeaderTree.Branch(branchname,"MMcTrigHeader",
851 &HeaderTrig[i], bsize, split);
852 }
853 }
854
855 if(!Trigger_Loop && Write_McFADC){
856
857 HeaderTree.Branch("MMcFadcHeader","MMcFadcHeader",
858 &HeaderFadc[0], bsize, split);
859 }
860 if (Trigger_Loop && Write_McFADC){
861 for(char branchname[20],i=0;i<icontrigger;i++){
862
863 sprintf(help,"%i",i+1);
864 strcpy (branchname, "MMcFadcHeader;");
865 strcat (branchname, & help[0]);
866 strcat (branchname, ".");
867 HeaderTree.Branch(branchname,"MMcFadcHeader",
868 &HeaderFadc[i], bsize, split);
869 }
870 }
871
872 // Fill branches for MRawRunHeader
873
874 RunHeader[0]->SetMagicNumber(kMagicNumber);
875 RunHeader[0]->SetFormatVersion(2);
876 RunHeader[0]->SetSoftVersion((UShort_t) (VERSION*10));
877 RunHeader[0]->SetRunType(256);
878 RunHeader[0]->SetRunNumber(0);
879 RunHeader[0]->SetNumSamples(0,FADC_SLICES);
880
881
882 // Fill branches for MMcTrigHeader
883
884 if(!Trigger_Loop && Write_McTrig){
885
886 HeaderTrig[0]->SetTopology((Short_t) Trigger_topology);
887 HeaderTrig[0]->SetMultiplicity((Short_t) Trigger_multiplicity);
888 for(i=0;i<TRIGGER_PIXELS;i++){
889 fpixelthres[i]=(Float_t)Trigger_threshold;
890 }
891 HeaderTrig[0]->SetThreshold( fpixelthres);
892 HeaderTrig[0]->SetAmplitud(Trigger_response_ampl);
893 HeaderTrig[0]->SetFwhm(Trigger_response_fwhm);
894 HeaderTrig[0]->SetOverlap(Trigger_overlaping_time);
895 HeaderTrig[0]->SetGate(Trigger_gate_length);
896 }
897 if(Trigger_Loop && Write_McTrig){
898
899 for (int iconcount=0,ithrescount=0;ithrescount<=Trigger_loop_uthres-Trigger_loop_lthres;ithrescount++){
900 for (imulticount=0;imulticount<=Trigger_loop_umult-Trigger_loop_lmult;imulticount++){
901 for(itopocount=0;itopocount<=Trigger_loop_utop-Trigger_loop_ltop;itopocount++){
902 HeaderTrig[iconcount]->SetTopology((Short_t) itopocount+Trigger_loop_ltop);
903 HeaderTrig[iconcount]->SetMultiplicity((Short_t) imulticount+Trigger_loop_lmult);
904 for(i=0;i<TRIGGER_PIXELS;i++){
905 fpixelthres[i]=(Float_t)(ithrescount+Trigger_loop_lthres);
906 }
907 HeaderTrig[iconcount]->SetThreshold( fpixelthres);
908 HeaderTrig[iconcount]->SetAmplitud(Trigger_response_ampl);
909 HeaderTrig[iconcount]->SetFwhm(Trigger_response_fwhm);
910 HeaderTrig[iconcount]->SetOverlap(Trigger_overlaping_time);
911 HeaderTrig[iconcount]->SetGate(Trigger_gate_length);
912 iconcount++;
913 }
914 }
915 }
916 }
917
918 // Fill branches for MMcFadcHeader
919
920 for(i=0;i<CAMERA_PIXELS;i++){
921 fadc_elecnoise[i]=2.0;
922 }
923 fadc.GetPedestals(&fadc_pedestals[0]);
924
925 if(!Trigger_Loop && Write_McFADC){
926
927 HeaderFadc[0]->SetShape(0.0);
928 HeaderFadc[0]->SetAmplitud(FADC_response_ampl);
929 HeaderFadc[0]->SetFwhm(FADC_response_fwhm);
930 HeaderFadc[0]->SetPedestal(&fadc_pedestals[0],anaPixels);
931 HeaderFadc[0]->SetElecNoise(&fadc_elecnoise[0],anaPixels);
932 }
933 if(Trigger_Loop && Write_McFADC){
934
935 for (int iconcount=0,ithrescount=0;ithrescount<=Trigger_loop_uthres-Trigger_loop_lthres;ithrescount++){
936 for (imulticount=0;imulticount<=Trigger_loop_umult-Trigger_loop_lmult;imulticount++){
937 for(itopocount=0;itopocount<=Trigger_loop_utop-Trigger_loop_ltop;itopocount++){
938 HeaderFadc[iconcount]->SetShape(0.0);
939 HeaderFadc[iconcount]->SetAmplitud(Trigger_response_ampl);
940 HeaderFadc[iconcount]->SetFwhm(Trigger_response_fwhm);
941 HeaderFadc[iconcount]->SetPedestal(&fadc_pedestals[0],anaPixels);
942 HeaderFadc[iconcount]->SetElecNoise(&fadc_elecnoise[0],anaPixels);
943 iconcount++;
944 }
945 }
946 }
947 }
948
949 // Fill the Header Tree with the current leaves of each branch
950 HeaderTree.Fill() ;
951
952
953 // create a Tree for the Event data stream
954 TTree EvtTree("Events","Normal Triggered Events");
955
956 if (Write_McEvt){
957
958 EvtTree.Branch("MMcEvt","MMcEvt",
959 &McEvt, bsize, split);
960 }
961
962 if(!Trigger_Loop){
963
964 if (Write_RawEvt){
965 EvtTree.Branch("MRawEvtHeader","MRawEvtHeader",
966 &EvtHeader[0], bsize, split);
967 EvtTree.Branch("MRawEvtData","MRawEvtData",
968 &EvtData[0], bsize, split);
969 }
970 if (Write_McTrig){
971 EvtTree.Branch("MMcTrig","MMcTrig",
972 &McTrig[0], bsize, split);
973 }
974 }
975 else{
976 if (Write_McTrig){
977 for(char branchname[10],i=0;i<icontrigger;i++){
978
979 sprintf(help,"%i",i+1);
980 strcpy (branchname, "MMcTrig;");
981 strcat (branchname, & help[0]);
982 strcat (branchname, ".");
983 EvtTree.Branch(branchname,"MMcTrig",
984 &McTrig[i], bsize, split);
985 }
986 }
987 }
988
989 if (Trigger_Loop && Write_RawEvt){
990 for(char branchname[15],i=0;i<icontrigger;i++){
991
992 sprintf(help,"%i",i+1);
993 strcpy (branchname, "MRawEvtHeader;");
994 strcat (branchname, & help[0]);
995 strcat (branchname, ".");
996 EvtTree.Branch(branchname,"MRawEvtHeader",
997 &EvtHeader[i], bsize, split);
998 }
999 for(char branchname[15],i=0;i<icontrigger;i++){
1000
1001 sprintf(help,"%i",i+1);
1002 strcpy (branchname, "MRawEvtData;");
1003 strcat (branchname, & help[0]);
1004 strcat (branchname, ".");
1005 EvtTree.Branch(branchname,"MRawEvtData",
1006 &EvtData[i], bsize, split);
1007 }
1008 }
1009
1010 unsigned short ulli = 0 ;
1011
1012 TApplication theAppTrigger("App", &argc, argv);
1013
1014 if (gROOT->IsBatch()) {
1015 fprintf(stderr, "%s: cannot run in batch mode\n", argv[0]);
1016 // return 1;
1017 }
1018
1019 if(FADC_Scan){
1020 TApplication theAppFadc("App", &argc, argv);
1021
1022 if (gROOT->IsBatch()) {
1023 fprintf(stderr, "%s: cannot run in batch mode\n", argv[0]);
1024 // return 1;
1025 }
1026 }
1027
1028 // set plate scale (deg/cm) and trigger area (deg)
1029
1030 plateScale_cm2deg = ( ct_Type == 0 ) ? (0.244/2.1) : 0.030952381;
1031
1032 if ( ! get_trigger_radius( &degTriggerZone ) )
1033 degTriggerZone = ( ct_Type == 0 ) ? (5.0) : (5.0);
1034
1035 if ( ! get_correction( &fCorrection ) )
1036 fCorrection = 1.0;
1037
1038 // prepare the NSB simulation
1039
1040 // Instance of the Mlons class
1041 MLons lons(Trigger_response_ampl, Trigger_response_fwhm,
1042 MFADC_RESPONSE_AMPLITUDE, MFADC_RESPONSE_FWHM);
1043
1044 lons.SetPath(nsbpathname);
1045
1046 if( simulateNSB ){
1047
1048 //
1049 // Calculate the non-diffuse NSB photoelectron rates
1050 //
1051
1052 k = produce_nsbrates( starfieldname,
1053 &cam,
1054 nsbrate_phepns );
1055
1056 if (k != 0){
1057 cout << "Error when reading starfield... \nExiting.\n";
1058 exit(1);
1059 }
1060
1061 // calculate diffuse rate correcting for the pixel size
1062
1063 for(i=0; i<cam.inumpixels; i++){
1064 diffnsb_phepns[i] = meanNSB *
1065 cam.dpixsizefactor[i] * cam.dpixsizefactor[i];
1066 }
1067
1068 // calculate nsb rate including diffuse and starlight
1069 for(i=0; i<cam.inumpixels;i++){
1070 nsb_phepns[i]= diffnsb_phepns[i];
1071 for(j=0;j<iNUMWAVEBANDS;j++)
1072 nsb_phepns[i]=nsb_phepns[i]+ nsbrate_phepns[i][j];
1073 }
1074 }
1075
1076 //
1077 // Read the reflector file with the Cherenkov data
1078 //
1079
1080 // select input file
1081
1082 if ( Data_From_STDIN ) {
1083
1084 inputfile = stdin;
1085
1086 }
1087 else{
1088
1089 log( SIGNATURE, "Opening input \"rfl\" file %s\n", inname );
1090 inputfile = fopen( inname, "r" );
1091 if ( inputfile == NULL )
1092 error( SIGNATURE, "Cannot open input file: %s\n", inname );
1093
1094 }
1095
1096 // get signature, and check it
1097
1098 if(check_reflector_file( inputfile )==FALSE){
1099 exit(1);
1100 }
1101
1102 // open data file
1103
1104 log( SIGNATURE, "Opening data \"dat\" file %s\n", datname );
1105 datafile.open( datname );
1106
1107 if ( datafile.bad() )
1108 error( SIGNATURE, "Cannot open data file: %s\n", datname );
1109
1110 // initializes flag
1111
1112 strcpy( flag, " \0" );
1113
1114 // allocate space for PMTs numbers of pixels
1115
1116 fnpix = new float [ ct_NPixels ];
1117
1118
1119 //!@' @#### Main loop.
1120 //@'
1121
1122 // get flag
1123
1124 fread( flag, SIZE_OF_FLAGS, 1, inputfile );
1125
1126 // loop over the file
1127
1128 still_in_loop = TRUE;
1129
1130 while (
1131 ((! Data_From_STDIN) && ( !feof(inputfile) ))
1132 ||
1133 (Data_From_STDIN && still_in_loop)
1134 ) {
1135
1136 // reading .rfl files
1137 if(!isA( flag, FLAG_START_OF_RUN )){
1138
1139 // We exit
1140 //error( SIGNATURE, "Expected start of run flag, but found: %s\n", flag );
1141 // We break the main loop
1142 cout<<"Warning: Expected start of run flag, but found:"<<flag<<endl;
1143 cout<<" We break the main loop"<<endl;
1144 break;
1145 }
1146 else { // found start of run
1147
1148 nshow=0;
1149
1150 fread( flag, SIZE_OF_FLAGS, 1, inputfile );
1151
1152 while( isA( flag, FLAG_START_OF_EVENT )){ // while there is a next event
1153 //
1154 // Clear Trigger and Fadc
1155 //
1156 Trigger.Reset() ;
1157 Trigger.ClearFirst();
1158 Trigger.ClearZero();
1159 fadc.Reset() ;
1160
1161 ++nshow;
1162 log(SIGNATURE, "Event %d(+%d)\n", nshow, ntshow);
1163
1164 // get MCEventHeader
1165 fread( (char*)&mcevth, mcevth.mysize(), 1, inputfile );
1166
1167 // calculate core distance and impact parameter
1168
1169 coreD = mcevth.get_core(&coreX, &coreY);
1170
1171 // calculate impact parameter (shortest distance betwee the original
1172 // trajectory of the primary (assumed shower-axis) and the
1173 // direction where the telescope points to
1174 //
1175 // we use the following equation, given that the shower core position
1176 // is (x1,y1,z1)=(x,y,0),the trajectory is given by (l1,m1,n1),
1177 // and the telescope position and orientation are (x2,y2,z2)=(0,0,0)
1178 // and (l2,m2,n2)
1179 //
1180 // | |
1181 // | x1-x2 y1-y2 z1-z2 |
1182 // | |
1183 // + | l1 m1 n1 |
1184 // - | |
1185 // | l2 m2 n2 |
1186 // | |
1187 // dist = ------------------------------------ ( > 0 )
1188 // [ |l1 m1|2 |m1 n1|2 |n1 l1|2 ]1/2
1189 // [ | | + | | + | | ]
1190 // [ |l2 m2| |m2 n2| |n2 l2| ]
1191 //
1192 // playing a little bit, we get this reduced for in our case:
1193 //
1194 //
1195 // dist = (- m2 n1 x + m1 n2 x + l2 n1 y - l1 n2 y - l2 m1 z + l1 m2 z) /
1196 // [(l2^2 (m1^2 + n1^2) + (m2 n1 - m1 n2)^2 -
1197 // 2 l1 l2 (m1 m2 + n1 n2) + l1^2 (m2^2 + n2^2) ] ^(1/2)
1198
1199 // read the direction of the incoming shower
1200
1201 thetashw = mcevth.get_theta();
1202 phishw = mcevth.get_phi();
1203
1204 // calculate vector for shower
1205
1206 l1 = sin(thetashw)*cos(phishw);
1207 m1 = sin(thetashw)*sin(phishw);
1208 n1 = cos(thetashw);
1209
1210 // read the deviation of the telescope with respect to the shower
1211
1212 mcevth.get_deviations ( &thetaCT, &phiCT );
1213
1214 if ( (thetaCT == 0.) && (phiCT == 0.) ) {
1215
1216 // CT was looking to the source (both lines are parallel)
1217 // therefore, we calculate the impact parameter as the distance
1218 // between the CT axis and the core position
1219
1220 impactD = dist_r_P( 0., 0., 0., l1, m1, n1, coreX, coreY, 0. );
1221
1222 } else {
1223
1224 // the shower comes off-axis
1225
1226 // obtain with this the final direction of the CT
1227
1228 thetaCT += thetashw;
1229 phiCT += phishw;
1230
1231 // calculate vector for telescope
1232
1233 l2 = sin(thetaCT)*cos(phiCT);
1234 m2 = sin(thetaCT)*sin(phiCT);
1235 n2 = cos(thetaCT);
1236
1237 num = (m1*n2*coreX - m2*n1*coreX + l2*n1*coreY - l1*n2*coreY);
1238 den = (SQR(l1*m2 - l2*m1) +
1239 SQR(m1*n2 - m2*n1) +
1240 SQR(n1*l2 - n2*l1));
1241 den = sqrt(den);
1242
1243 impactD = fabs(num)/den;
1244
1245 }
1246
1247 // energy cut
1248
1249 if ( Select_Energy ) {
1250 if (( mcevth.get_energy() < Select_Energy_le ) ||
1251 ( mcevth.get_energy() > Select_Energy_ue )) {
1252 log(SIGNATURE, "select_energy: shower rejected.\n");
1253 continue;
1254 }
1255 }
1256
1257 // Read first and last time and put inumphe to 0
1258
1259 mcevth.get_times(&arrtmin_ns,&arrtmax_ns);
1260 inumphe=0;
1261
1262 // read the photons and produce the photoelectrons
1263
1264 k = produce_phes( inputfile,
1265 &cam,
1266 WAVEBANDBOUND1,
1267 WAVEBANDBOUND6,
1268 &Trigger, // will be changed by the function!
1269 &fadc, // will be changed by the function!
1270 &inumphe, // important for later: the size of photoe[]
1271 fnpix, // will be changed by the function!
1272 &ncph, // will be changed by the function!
1273 &arrtmin_ns, // will be changed by the function!
1274 &arrtmax_ns // will be changed by the function!
1275 );
1276
1277 if( k != 0 ){ // non-zero returnvalue means error
1278 cout << "Exiting.\n";
1279 exit(1);
1280 }
1281
1282 // NSB simulation
1283
1284 if(simulateNSB && nphe2NSB<=inumphe){
1285
1286 // Fill trigger and fadc response in the trigger class from the database
1287 for(i=0;i<cam.inumpixels;i++)
1288 if(nsb_phepns[i]>0.0){
1289 k=lons.GetResponse(nsb_phepns[i],0.1,
1290 & nsb_trigresp[0],& nsb_fadcresp[0]);
1291 if(k==0){
1292 cout << "Exiting.\n";
1293 exit(1);
1294 }
1295 Trigger.AddNSB(i,nsb_trigresp);
1296 fadc.AddSignal(i,nsb_fadcresp);
1297 }
1298 }// end if(simulateNSB && nphe2NSB<=inumphe) ...
1299
1300 log(SIGNATURE, "End of this event: %d cphs(+%d). . .\n",
1301 ncph, ntcph);
1302
1303 ntcph += ncph;
1304
1305 // skip it ?
1306
1307 for ( i=0; i<nSkip; ++i ) {
1308 if (Skip[i] == (nshow+ntshow)) {
1309 i = -1;
1310 break;
1311 }
1312 }
1313
1314 // if after the previous loop, the exit value of i is -1
1315 // then the shower number is in the list of showers to be
1316 // skipped
1317
1318 if (i == -1) {
1319 log(SIGNATURE, "\t\tskipped!\n");
1320 continue;
1321 }
1322
1323 cout << "Total number of phes: " << inumphe<<" + ";
1324 inumphensb=0;
1325 for (i=0;i<cam.inumpixels;i++){
1326 inumphensb+=nsb_phepns[i]*TOTAL_TRIGGER_TIME;
1327 }
1328 cout<<inumphensb<<endl;
1329
1330 //++++++++++++++++++++++++++++++++++++++++++++++++++
1331 // at this point we have a camera full of
1332 // ph.e.s
1333 // we should first apply the trigger condition,
1334 // and if there's trigger, then clean the image,
1335 // calculate the islands statistics and the
1336 // other parameters of the image (Hillas' parameters
1337 // and so on).
1338 //--------------------------------------------------
1339
1340 // TRIGGER HERE
1341
1342 //
1343 // now the noise of the electronic
1344 // (preamps, optical transmission,..) is introduced.
1345 // This is done inside the class MTrigger by the method ElecNoise.
1346 //
1347 Trigger.ElecNoise() ;
1348
1349 fadc.ElecNoise() ;
1350
1351 // now a shift in the fadc signal due to the pedestlas is
1352 // intorduced
1353 // This is done inside the class MFadc by the method Pedestals
1354 fadc.Pedestals();
1355
1356 // We study several trigger conditons
1357 if(Trigger_Loop){
1358 // Loop over trigger threshold
1359 for (int iconcount=0,ithrescount=Trigger_loop_lthres;ithrescount<=Trigger_loop_uthres;ithrescount++){
1360 for (i=0;i<TRIGGER_PIXELS;i++)
1361 fpixelthres[i]=(float) ithrescount;
1362 Trigger.SetThreshold(fpixelthres);
1363
1364 Trigger.Diskriminate();
1365 //
1366 // look if in all the signals in the trigger signal branch
1367 // is a possible Trigger. Therefore we habe to diskriminate all
1368 // the simulated analog signals (Method Diskriminate in class
1369 // MTrigger). We look simultanously for the moments at which
1370 // there are more than TRIGGER_MULTI pixels above the
1371 // CHANNEL_THRESHOLD.
1372 //
1373
1374 // Set trigger flags to zero
1375 Lev1=Lev2=0;
1376 btrigger=0;
1377
1378 // loop over multiplicity of trigger configuration
1379 for (imulticount=Trigger_loop_lmult;imulticount<=Trigger_loop_umult;imulticount++){
1380 Trigger.SetMultiplicity(imulticount);
1381 Trigger.ClearZero();
1382
1383 Lev0=(Short_t) Trigger.ZeroLevel();
1384 if (Lev0>0 || Write_All_Images || btrigger){
1385
1386 // loop over topologies
1387 for(itopocount=Trigger_loop_ltop;itopocount<=Trigger_loop_utop;itopocount++){
1388 Lev1=Lev2=0;
1389
1390 if(itopocount==0 && imulticount>7) continue;
1391 if(itopocount==2 && imulticount<3) continue;
1392 Trigger.SetTopology(itopocount);
1393 Trigger.ClearFirst();
1394
1395 //
1396 // Start the First Level Trigger simulation
1397 //
1398 Lev1=Trigger.FirstLevel();
1399 if (Write_McTrig)
1400 McTrig[iconcount]->SetFirstLevel (Lev1);
1401 if(Lev1>0) {
1402 btrigger= 1;
1403 ntriggerloop[ithrescount-Trigger_loop_lthres][imulticount-Trigger_loop_lmult][itopocount-Trigger_loop_ltop]++;
1404 }
1405
1406 if(Lev1==0 && (Write_All_Images || btrigger)){
1407 btrigger= 1;
1408 Lev1=1;
1409 }
1410 numPix=0;
1411 for (Int_t ii=0;ii<Lev1;ii++){
1412 if (Write_McTrig)
1413 McTrig[iconcount]->SetTime(Trigger.GetFirstLevelTime(ii),ii+1);
1414 if (Write_McTrig){
1415 Trigger.GetMapDiskriminator(trigger_map);
1416 McTrig[iconcount]->SetMapPixels(trigger_map,ii);
1417 }
1418 //
1419 // fill inside the class fadc the member output
1420 //
1421
1422 fadc.TriggeredFadc(Trigger.GetFirstLevelTime(ii));
1423
1424 if( Write_RawEvt ){
1425 //
1426 // Fill the header of this event
1427 //
1428
1429 EvtHeader[iconcount]->FillHeader ( (UInt_t) (ntshow + nshow),0);
1430 // fill pixel information
1431 for(i=0;i<ct_NPixels;i++){
1432 for (j=0;j<FADC_SLICES;j++){
1433 fadcValues->AddAt(fadc.GetFadcSignal(i,j),j);
1434 }
1435 EvtData[iconcount]->AddPixel(i+1000*ii,fadcValues,0);
1436 }
1437 }
1438 }
1439 //
1440 // Increase counter of analised trigger conditions
1441 //
1442 iconcount++;
1443 }
1444 }
1445 else{
1446 break;
1447 }
1448 }
1449 if (!btrigger) break;
1450 }
1451 if (btrigger){
1452
1453 //
1454 // fill the MMcEvt with all information
1455 //
1456
1457 if (Write_McEvt) {
1458 McEvt->Fill( (UShort_t) mcevth.get_primary() ,
1459 mcevth.get_energy(),
1460 mcevth.get_theta(),
1461 mcevth.get_phi(),
1462 mcevth.get_core(),
1463 mcevth.get_coreX(),
1464 mcevth.get_coreY(),
1465 impactD,
1466 (UInt_t)mcevth.get_CORSIKA(),
1467 (UInt_t)mcevth.get_AtmAbs(),
1468 (UInt_t)mcevth.get_MirrAbs()+mcevth.get_OutOfMirr()+mcevth.get_BlackSpot(),
1469 (UInt_t) ncph,
1470 (UInt_t) inumphe,
1471 (UInt_t) inumphensb+inumphe);
1472 }
1473 // Fill the Tree with the current leaves of each branch
1474 i=EvtTree.Fill() ;
1475
1476 // Clear the branches
1477 if(Write_McTrig){
1478 for(i=0;i<icontrigger;i++){
1479 McTrig[i]->Clear() ;
1480 }
1481 }
1482 if( Write_RawEvt ){
1483 for(i=0;i<icontrigger;i++){
1484 EvtHeader[i]->Clear() ;
1485 EvtData[i]->DeletePixels();
1486 }
1487 }
1488 if (Write_McEvt)
1489 McEvt->Clear() ;
1490 }
1491 }
1492 // We study a single trigger condition
1493 else {
1494
1495 // Setting trigger conditions
1496 Trigger.SetMultiplicity(Trigger_multiplicity);
1497 Trigger.SetTopology(Trigger_topology);
1498 for (i=0;i<TRIGGER_PIXELS;i++)
1499 fpixelthres[i]=Trigger_threshold;
1500 Trigger.SetThreshold(fpixelthres);
1501
1502 Trigger.Diskriminate() ;
1503
1504 //
1505 // look if in all the signals in the trigger signal branch
1506 // is a possible Trigger. Therefore we habe to diskriminate all
1507 // the simulated analog signals (Method Diskriminate in class
1508 // MTrigger). We look simultanously for the moments at which
1509 // there are more than TRIGGER_MULTI pixels above the
1510 // CHANNEL_THRESHOLD.
1511 //
1512
1513 Lev0 = (Short_t) Trigger.ZeroLevel() ;
1514
1515 Lev1 = Lev2 = 0 ;
1516
1517 //
1518 // Start the First Level Trigger simulation
1519 //
1520
1521 if ( Lev0 > 0 || Write_All_Images) {
1522 Lev1= Trigger.FirstLevel();
1523 if (Write_McTrig)
1524 McTrig[0]->SetFirstLevel (Lev1);
1525 }
1526 if (Lev1>0){
1527 ++ntrigger;
1528 }
1529 if (Lev1==0 && Write_All_Images){
1530 Lev1=1;
1531 }
1532
1533 numPix=0;
1534 for(Int_t ii=0;ii<Lev1;ii++){
1535 // Loop over different level one triggers
1536
1537 //
1538 // fill inside class fadc the member output
1539 //
1540 fadc.TriggeredFadc(Trigger.GetFirstLevelTime(ii));
1541
1542 if (Write_McTrig)
1543 McTrig[0]->SetTime(Trigger.GetFirstLevelTime(ii),ii+1);
1544
1545 if (Write_McTrig){
1546 Trigger.GetMapDiskriminator(trigger_map);
1547 McTrig[0]->SetMapPixels(trigger_map,ii);
1548 }
1549
1550 // Fill Evt information
1551
1552 if (Write_RawEvt){
1553
1554 //
1555 // Fill the header of this event
1556 //
1557
1558 EvtHeader[0]->FillHeader ( (UShort_t) (ntshow + nshow) , 20 ) ;
1559
1560 // fill pixel information
1561
1562 for(i=0;i<ct_NPixels;i++){
1563 for (j=0;j<FADC_SLICES;j++){
1564 fadcValues->AddAt(fadc.GetFadcSignal(i,j),j);
1565 }
1566 EvtData[0]->AddPixel(i+ii*1000,fadcValues,0);
1567 }
1568 }
1569 //
1570 // fill the MMcEvt with all information
1571 //
1572
1573 if (Write_McEvt){
1574 McEvt->Fill( (UShort_t) mcevth.get_primary() ,
1575 mcevth.get_energy(),
1576 mcevth.get_theta(),
1577 mcevth.get_phi(),
1578 mcevth.get_core(),
1579 mcevth.get_coreX(),
1580 mcevth.get_coreY(),
1581 impactD,
1582 (UInt_t)mcevth.get_CORSIKA(),
1583 (UInt_t)mcevth.get_AtmAbs(),
1584 (UInt_t)mcevth.get_MirrAbs()+mcevth.get_OutOfMirr()+mcevth.get_BlackSpot(),
1585 (UInt_t) ncph,
1586 (UInt_t) inumphe,
1587 (UInt_t) inumphensb+inumphe);
1588 }
1589 // We don not count photons out of the camera.
1590
1591
1592 //
1593 // write it out to the file outfile
1594 //
1595
1596 EvtTree.Fill() ;
1597
1598 //
1599 // if a first level trigger occurred, then
1600 // 1. do some other stuff (not implemented)
1601 // 2. start the gui tool
1602
1603 if(FADC_Scan){
1604 if ( Lev0 > 0 ) {
1605 fadc.ShowSignal( McEvt, (Float_t) 60. ) ;
1606 }
1607 }
1608
1609 if(Trigger_Scan){
1610 if ( Lev0 > 0 ) {
1611 Trigger.ShowSignal(McEvt) ;
1612 }
1613 }
1614
1615 // clear all
1616 if (Write_RawEvt) EvtHeader[0]->Clear() ;
1617 if (Write_RawEvt) EvtData[0]->DeletePixels();
1618 if (Write_McEvt) McEvt->Clear() ;
1619 }
1620 if (Write_McTrig) McTrig[0]->Clear() ;
1621 }
1622
1623#ifdef __DEBUG__
1624 printf("\n");
1625
1626 for ( ici=0; ici<PIX_ARRAY_SIDE; ++ici ) {
1627
1628 for ( icj=0; icj<PIX_ARRAY_SIDE; ++icj ) {
1629
1630 if ( (int)pixels[ici][icj][PIXNUM] > -1 ) {
1631
1632 if ( fnpix[(int)pixels[ici][icj][PIXNUM]] > 0. ) {
1633
1634 printf ("@@ %4d %4d %10f %10f %4f (%4d %4d)\n", nshow,
1635 (int)pixels[ici][icj][PIXNUM],
1636 pixels[ici][icj][PIXX],
1637 pixels[ici][icj][PIXY],
1638 fnpix[(int)pixels[ici][icj][PIXNUM]], ici, icj);
1639
1640 }
1641 }
1642 }
1643 }
1644
1645 for (i=0; i<ct_NPixels; ++i) {
1646 printf("%d (%d): ", i, npixneig[i]);
1647 for (j=0; j<npixneig[i]; ++i)
1648 printf(" %d", pixneig[i][j]);
1649 printf("\n");
1650 }
1651
1652#endif // __DEBUG__
1653
1654
1655 // look for the next event
1656
1657 fread( flag, SIZE_OF_FLAGS, 1, inputfile );
1658
1659 } // end while there is a next event
1660
1661 if( !isA( flag, FLAG_END_OF_RUN )){
1662 error( SIGNATURE, "Expected end of run flag, but found: %s\n", flag );
1663 }
1664 else { // found end of run
1665 ntshow += nshow;
1666 log(SIGNATURE, "End of this run with %d events . . .\n", nshow);
1667
1668 fread( flag, SIZE_OF_FLAGS, 1, inputfile );
1669
1670 if( isA( flag, FLAG_END_OF_FILE ) ){ // end of file
1671 log(SIGNATURE, "End of file . . .\n");
1672 still_in_loop = FALSE;
1673
1674 if ((! Data_From_STDIN) && ( !feof(inputfile) )){
1675
1676 // we have concatenated input files.
1677 // get signature of the next part and check it.
1678
1679 if(check_reflector_file( inputfile )==FALSE){
1680 exit(1);
1681 }
1682
1683 }
1684
1685 fread( flag, SIZE_OF_FLAGS, 1, inputfile );
1686
1687 } // end if found end of file
1688 } // end if found end of run
1689
1690 } // end if else found start of run
1691 } // end big while loop
1692
1693 //++
1694 // put the Event to the root file
1695 //--
1696
1697 outfile_temp.Write() ;
1698 outfile_temp.Close() ;
1699
1700
1701 // close input file
1702
1703 log( SIGNATURE, "%d event(s), with a total of %d C.photons\n",
1704 ntshow, ntcph );
1705 datafile<<ntshow<<" event(s), with a total of "<<ntcph<<" C.photons"<<endl;
1706 if (Trigger_Loop){
1707 log( SIGNATURE, "Trigger Mode. \n");
1708 log( SIGNATURE, "Fraction of triggers: \n");
1709 datafile<<"Fraction of triggers: "<<endl;
1710 for (ithrescount=Trigger_loop_lthres;ithrescount<=Trigger_loop_uthres;ithrescount++){
1711 for (imulticount=Trigger_loop_lmult;imulticount<=Trigger_loop_umult;imulticount++){
1712 for(itopocount=Trigger_loop_ltop;itopocount<=Trigger_loop_utop;itopocount++){
1713 log( SIGNATURE, "Thres %d, Multi %d, Topo %d: %5.1f%% (%d out of %d)\n",
1714 ithrescount,imulticount,itopocount,((float)ntriggerloop[ithrescount-Trigger_loop_lthres][imulticount-Trigger_loop_lmult][itopocount-Trigger_loop_ltop] / ((float)ntshow) * 100.0), ntriggerloop[ithrescount-Trigger_loop_lthres][imulticount-Trigger_loop_lmult][itopocount-Trigger_loop_ltop], ntshow);
1715 datafile<<"Thres "<<ithrescount<<", Multi "<<imulticount<<", Topo"<<itopocount<<": ";
1716 datafile<<((float)ntriggerloop[ithrescount-Trigger_loop_lthres][imulticount-Trigger_loop_lmult][itopocount-Trigger_loop_ltop] / ((float)ntshow) * 100.0)<<"% ("<<ntriggerloop[ithrescount-Trigger_loop_lthres][imulticount-Trigger_loop_lmult][itopocount-Trigger_loop_ltop]<<" out of "<<ntshow<<")"<<endl;
1717 }
1718 }
1719 }
1720 }
1721 else{
1722 log( SIGNATURE, "Fraction of triggers: %5.1f%% (%d out of %d)\n",
1723 ((float)ntrigger) / ((float)ntshow) * 100.0, ntrigger, ntshow);
1724 datafile<<"Fraction of triggers: "<<((float)ntrigger) / ((float)ntshow) * 100.0<<" ("<<ntrigger<<" out of "<<ntshow<<" )"<<endl;
1725 }
1726
1727 // close files
1728
1729 log( SIGNATURE, "Closing files\n" );
1730
1731 if( ! Data_From_STDIN ){
1732 fclose( inputfile );
1733 }
1734 datafile.close();
1735
1736 // program finished
1737
1738 log( SIGNATURE, "Done.\n");
1739
1740 return( 0 );
1741}
1742//!@}
1743
1744// @T \newpage
1745
1746//!@subsection Functions definition.
1747
1748//!-----------------------------------------------------------
1749// @name present
1750//
1751// @desc Make some presentation
1752//
1753// @date Sat Jun 27 05:58:56 MET DST 1998
1754//------------------------------------------------------------
1755// @function
1756
1757//!@{
1758void
1759present(void)
1760{
1761 cout << "##################################################\n"
1762 << SIGNATURE << '\n' << '\n'
1763 << "Processor of the reflector output\n"
1764 << "J C Gonzalez, Jun 1998\n"
1765 << "##################################################\n\n"
1766 << flush ;
1767}
1768//!@}
1769
1770
1771//!-----------------------------------------------------------
1772// @name usage
1773//
1774// @desc show help
1775//
1776// @date Tue Dec 15 16:23:30 MET 1998
1777//------------------------------------------------------------
1778// @function
1779
1780//!@{
1781void
1782usage(void)
1783{
1784 present();
1785 cout << "\nusage ::\n\n"
1786 << "\t camera "
1787 << " [ -@ paramfile ] "
1788 << " [ -h ] "
1789 << "\n\n or \n\n"
1790 << "\t camera < paramfile"
1791 << "\n\n";
1792 exit(0);
1793}
1794//!@}
1795
1796
1797//!-----------------------------------------------------------
1798// @name log
1799//
1800// @desc function to send log information
1801//
1802// @var funct Name of the caller function
1803// @var fmt Format to be used (message)
1804// @var ... Other information to be shown
1805//
1806// @date Sat Jun 27 05:58:56 MET DST 1998
1807//------------------------------------------------------------
1808// @function
1809
1810//!@{
1811void
1812log(const char *funct, char *fmt, ...)
1813{
1814 va_list args;
1815
1816 // Display the name of the function that called error
1817 printf("[%s]: ", funct);
1818
1819 // Display the remainder of the message
1820 va_start(args, fmt);
1821 vprintf(fmt, args);
1822 va_end(args);
1823}
1824//!@}
1825
1826
1827//!-----------------------------------------------------------
1828// @name error
1829//
1830// @desc function to send an error message, and abort the program
1831//
1832// @var funct Name of the caller function
1833// @var fmt Format to be used (message)
1834// @var ... Other information to be shown
1835//
1836// @date Sat Jun 27 05:58:56 MET DST 1998
1837//------------------------------------------------------------
1838// @function
1839
1840//!@{
1841void
1842error(const char *funct, char *fmt, ...)
1843{
1844 va_list args;
1845
1846 // Display the name of the function that called error
1847 fprintf(stdout, "ERROR in %s: ", funct);
1848
1849 // Display the remainder of the message
1850 va_start(args, fmt);
1851 vfprintf(stdout, fmt, args);
1852 va_end(args);
1853
1854 perror(funct);
1855
1856 exit(1);
1857}
1858//!@}
1859
1860
1861//!-----------------------------------------------------------
1862// @name isA
1863//
1864// @desc returns TRUE(FALSE), if the flag is(is not) the given
1865//
1866// @var s1 String to be searched
1867// @var flag Flag to compare with string s1
1868// @return TRUE: both strings match; FALSE: oth.
1869//
1870// @date Wed Jul 8 15:25:39 MET DST 1998
1871//------------------------------------------------------------
1872// @function
1873
1874//!@{
1875int
1876isA( char * s1, const char * flag ) {
1877 return ( (strncmp((char *)s1, flag, SIZE_OF_FLAGS)==0) ? 1 : 0 );
1878}
1879//!@}
1880
1881
1882//!-----------------------------------------------------------
1883// @name read_ct_file
1884//
1885// @desc read CT definition file
1886//
1887// @date Sat Jun 27 05:58:56 MET DST 1998
1888//------------------------------------------------------------
1889// @function
1890
1891//!@{
1892void
1893read_ct_file(void)
1894{
1895 char line[LINE_MAX_LENGTH]; //@< line to get from the ctin
1896 char token[ITEM_MAX_LENGTH]; //@< a single token
1897 int i, j; //@< dummy counters
1898
1899 log( "read_ct_file", "start.\n" );
1900
1901 ifstream ctin ( ct_filename );
1902
1903 if ( ctin.bad() )
1904 error( "read_ct_file",
1905 "Cannot open CT def. file: %s\n", ct_filename );
1906
1907 // loop till the "end" directive is reached
1908
1909 while (!ctin.eof()) {
1910
1911 // get line from stdin
1912
1913 ctin.getline(line, LINE_MAX_LENGTH);
1914
1915 // look for each item at the beginning of the line
1916
1917 for (i=0; i<=define_mirrors; i++)
1918 if (strstr(line, CT_ITEM_NAMES[i]) == line)
1919 break;
1920
1921 // if it is not a valid line, just ignore it
1922
1923 if (i == define_mirrors+1)
1924 continue;
1925
1926 // case block for each directive
1927
1928 switch ( i ) {
1929
1930 case type: // <type of telescope> (0:CT1 ¦ 1:MAGIC)
1931
1932 // get focal distance
1933
1934 sscanf(line, "%s %d", token, &ct_Type);
1935
1936 log( "read_ct_file", "<Type of Telescope>: %s\n",
1937 ((ct_Type==0) ? "CT1" : "MAGIC") );
1938
1939 break;
1940
1941 case focal_distance: // <focal distance> [cm]
1942
1943 // get focal distance
1944
1945 sscanf(line, "%s %f", token, &ct_Focal_mean);
1946
1947 log( "read_ct_file", "<Focal distance>: %f cm\n", ct_Focal_mean );
1948
1949 break;
1950
1951 case focal_std: // s(focal distance) [cm]
1952
1953 // get focal distance
1954
1955 sscanf(line, "%s %f", token, &ct_Focal_std);
1956
1957 log( "read_ct_file", "s(Focal distance): %f cm\n", ct_Focal_std );
1958
1959 break;
1960
1961 case point_spread: // <point spread> [cm]
1962
1963 // get point spread
1964
1965 sscanf(line, "%s %f", token, &ct_PSpread_mean);
1966
1967 log( "read_ct_file", "<Point spread>: %f cm\n", ct_PSpread_mean );
1968
1969 break;
1970
1971 case point_std: // s(point spread) [cm]
1972
1973 // get point spread
1974
1975 sscanf(line, "%s %f", token, &ct_PSpread_std);
1976
1977 log( "read_ct_file", "s(Point spread): %f cm\n", ct_PSpread_std );
1978
1979 break;
1980
1981 case adjustment_dev: // s(adjustment_dev) [cm]
1982
1983 // get point spread
1984
1985 sscanf(line, "%s %f", token, &ct_Adjustment_std);
1986
1987 log( "read_ct_file", "s(Adjustment): %f cm\n", ct_Adjustment_std );
1988
1989 break;
1990
1991 case black_spot: // radius of the black spot in the center [cm]
1992
1993 // get black spot radius
1994
1995 sscanf(line, "%s %f", token, &ct_BlackSpot_rad);
1996
1997 log( "read_ct_file", "Radius of the black spots: %f cm\n",
1998 ct_BlackSpot_rad);
1999
2000 break;
2001
2002 case r_mirror: // radius of the mirrors [cm]
2003
2004 // get radius of mirror
2005
2006 sscanf(line, "%s %f", token, &ct_RMirror);
2007
2008 log( "read_ct_file", "Radii of the mirrors: %f cm\n", ct_RMirror );
2009
2010 break;
2011
2012 case n_mirrors: // number of mirrors
2013
2014 // get the name of the output_file from the line
2015
2016 sscanf(line, "%s %d", token, &ct_NMirrors);
2017
2018 log( "read_ct_file", "Number of mirrors: %d\n", ct_NMirrors );
2019
2020 break;
2021
2022 case camera_width: // camera width [cm]
2023
2024 // get the name of the ct_file from the line
2025
2026 sscanf(line, "%s %f", token, &ct_CameraWidth);
2027
2028 log( "read_ct_file", "Camera width: %f cm\n", ct_CameraWidth );
2029
2030 break;
2031
2032 case n_pixels: // number of pixels
2033
2034 // get the name of the output_file from the line
2035
2036 sscanf(line, "%s %d", token, &ct_NPixels);
2037
2038 log( "read_ct_file", "Number of pixels: %d\n", ct_NPixels );
2039
2040 break;
2041
2042 case n_centralpixels: // number of central pixels
2043
2044 // get the name of the output_file from the line
2045
2046 sscanf(line, "%s %d", token, &ct_NCentralPixels);
2047
2048 log( "read_ct_file", "Number of central pixels: %d\n", ct_NCentralPixels );
2049
2050 break;
2051
2052 case n_gappixels: // number of gap pixels
2053
2054 // get the name of the output_file from the line
2055
2056 sscanf(line, "%s %d", token, &ct_NGapPixels);
2057
2058 log( "read_ct_file", "Number of gap pixels: %d\n", ct_NGapPixels );
2059
2060 break;
2061
2062 case pixel_width: // pixel width [cm]
2063
2064 // get the name of the ct_file from the line
2065
2066 sscanf(line, "%s %f", token, &ct_PixelWidth);
2067
2068 ct_PixelWidth_corner_2_corner = ct_PixelWidth / cos(RAD(30.0));
2069 ct_PixelWidth_corner_2_corner_half =
2070 ct_PixelWidth_corner_2_corner * 0.50;
2071 ct_Apot = ct_PixelWidth / 2;
2072 ct_2Apot = ct_Apot * 2.0;
2073
2074 log( "read_ct_file", "Pixel width: %f cm\n", ct_PixelWidth );
2075
2076 break;
2077
2078 case define_mirrors: // read table with the parameters of the mirrors
2079
2080 log( "read_ct_file", "Table of mirrors data:\n" );
2081
2082 // check whether the number of mirrors was already set
2083
2084 if ( ct_NMirrors == 0 )
2085 error( "read_ct_file", "NMirrors was not set.\n" );
2086
2087 // allocate memory for paths list
2088
2089 log( "read_ct_file", "Allocating memory for ct_data\n" );
2090
2091 ct_data = new float*[ct_NMirrors];
2092
2093 for (i=0; i<ct_NMirrors; i++)
2094 ct_data[i] = new float[CT_NDATA];
2095
2096 // read data
2097
2098 log( "read_ct_file", "Reading mirrors data...\n" );
2099
2100 for (i=0; i<ct_NMirrors; i++)
2101 for (j=0; j<CT_NDATA; j++)
2102 ctin >> ct_data[i][j];
2103
2104 break;
2105
2106 } // switch ( i )
2107
2108 } // end while
2109
2110 // end
2111
2112 log( "read_ct_file", "done.\n" );
2113
2114 return;
2115}
2116//!@}
2117
2118
2119//!-----------------------------------------------------------
2120// @name read_pixels
2121//
2122// @desc read pixels data
2123//
2124// @date Fri Mar 12 16:33:34 MET 1999
2125//------------------------------------------------------------
2126// @function
2127
2128//!@{
2129void
2130read_pixels(struct camera *pcam)
2131{
2132 ifstream qefile;
2133 char line[LINE_MAX_LENGTH];
2134 int n, i, j, icount;
2135 float qe;
2136
2137 //------------------------------------------------------------
2138 // first, pixels' coordinates
2139
2140 pcam->inumpixels = ct_NPixels;
2141 pcam->inumcentralpixels = ct_NCentralPixels;
2142 pcam->inumgappixels = ct_NGapPixels;
2143 pcam->inumbigpixels = ct_NPixels - ct_NCentralPixels - ct_NGapPixels;
2144 pcam->dpixdiameter_cm = ct_PixelWidth;
2145
2146 // initialize pixel numbers
2147
2148 pixary = new float* [2*ct_NCentralPixels];
2149 for ( i=0; i<2*ct_NCentralPixels; ++i )
2150 pixary[i] = new float[2];
2151
2152 pixneig = new int* [ct_NCentralPixels];
2153 for ( i=0; i<ct_NCentralPixels; ++i ) {
2154 pixneig[i] = new int[6];
2155 for ( j=0; j<6; ++j )
2156 pixneig[i][j] = -1;
2157 }
2158
2159 npixneig = new int[ct_NCentralPixels];
2160 for ( i=0; i<ct_NCentralPixels; ++i )
2161 npixneig[i] = 0;
2162
2163 // generate all coordinates
2164
2165 igen_pixel_coordinates(pcam);
2166
2167
2168 // calculate tables of neighbours
2169
2170#ifdef __DEBUG__
2171 for ( n=0 ; n<ct_NPixels ; ++n ) {
2172 cout << "Para el pixel " << n << ": ";
2173 for ( i=n+1 ; (i<ct_NPixels)&&(npixneig[n]<6) ; ++i) {
2174 if ( pixels_are_neig(n,i) == TRUE ) {
2175 pixneig[n][npixneig[n]] = i;
2176 pixneig[i][npixneig[i]] = n;
2177 cout << i << ' ';
2178 ++npixneig[n];
2179 ++npixneig[i];
2180 }
2181 }
2182 cout << endl << flush;
2183 }
2184#else // ! __DEBUG__
2185 for ( n=0 ; n<ct_NCentralPixels ; ++n )
2186 for ( i=n+1 ; (i<ct_NCentralPixels)&&(npixneig[n]<6) ; ++i)
2187 if ( pixels_are_neig(n,i) == TRUE ) {
2188 pixneig[n][npixneig[n]] = i;
2189 pixneig[i][npixneig[i]] = n;
2190 ++npixneig[n];
2191 ++npixneig[i];
2192 }
2193#endif // ! __DEBUG__
2194
2195#ifdef __DEBUG__
2196 for ( n=0 ; n<ct_NPixels ; ++n ) {
2197 cout << n << ':';
2198 for ( j=0; j<npixneig[n]; ++j)
2199 cout << ' ' << pixneig[n][j];
2200 cout << endl << flush;
2201 }
2202#endif // __DEBUG__
2203
2204 //------------------------------------------------------------
2205 // second, pixels' QE
2206
2207 // try to open the file
2208
2209 log("read_pixels", "Opening the file \"%s\" . . .\n", QE_FILE);
2210
2211 qefile.open( QE_FILE );
2212
2213 // if it is wrong or does not exist, exit
2214
2215 if ( qefile.bad() )
2216 error( "read_pixels", "Cannot open \"%s\". Exiting.\n", QE_FILE );
2217
2218 // read file
2219
2220 log("read_pixels", "Reading data . . .\n");
2221
2222 i=-1;
2223 icount = 0;
2224
2225 while ( ! qefile.eof() ) {
2226
2227 // get line from the file
2228
2229 qefile.getline(line, LINE_MAX_LENGTH);
2230
2231 // skip if comment
2232
2233 if ( *line == '#' )
2234 continue;
2235
2236 // if it is the first valid value, it is the number of QE data points
2237
2238 if ( i < 0 ) {
2239
2240 // get the number of datapoints
2241
2242 sscanf(line, "%d", &pointsQE);
2243
2244 // allocate memory for the table of QEs
2245
2246 QE = new float ** [ct_NPixels];
2247
2248 for ( i=0; i<ct_NPixels; ++i ) {
2249 QE[i] = new float * [2];
2250 QE[i][0] = new float[pointsQE];
2251 QE[i][1] = new float[pointsQE];
2252 }
2253
2254 QElambda = new float [pointsQE];
2255
2256 for ( i=0; i<pointsQE; ++i ) {
2257 qefile.getline(line, LINE_MAX_LENGTH);
2258 sscanf(line, "%f", &QElambda[i]);
2259 }
2260
2261 i=0;
2262
2263 continue;
2264 }
2265
2266 // get the values (num-pixel, num-datapoint, QE-value)
2267
2268 if( sscanf(line, "%d %d %f", &i, &j, &qe) != 3 )
2269 break;
2270
2271 if ( ((i-1) < ct_NPixels) && ((i-1) > -1) &&
2272 ((j-1) < pointsQE) && ((j-1) > -1) ) {
2273 QE[i-1][0][j-1] = QElambda[j-1];
2274 QE[i-1][1][j-1] = qe;
2275 }
2276
2277 if ( i > ct_NPixels) break;
2278
2279 icount++;
2280
2281 }
2282
2283 if(icount/pointsQE < ct_NPixels){
2284 error( "read_pixels", "The quantum efficiency file is faulty\n (found only %d pixels instead of %d).\n",
2285 icount/pointsQE, ct_NPixels );
2286 }
2287
2288 // close file
2289
2290 qefile.close();
2291
2292 // test QE
2293
2294 for(icount=0; icount< ct_NPixels; icount++){
2295 for(i=0; i<pointsQE; i++){
2296 if( QE[icount][0][i] < 100. || QE[icount][0][i] > 1000. ||
2297 QE[icount][1][i] < 0. || QE[icount][1][i] > 100.){
2298 error( "read_pixels", "The quantum efficiency file is faulty\n pixel %d, point %d is % f, %f\n",
2299 icount, i, QE[icount][0][i], QE[icount][1][i] );
2300 }
2301 }
2302 }
2303
2304 // end
2305
2306 log("read_pixels", "Done.\n");
2307
2308}
2309//!@}
2310
2311
2312//!-----------------------------------------------------------
2313// @name pixels_are_neig
2314//
2315// @desc check whether two pixels are neighbours
2316//
2317// @var pix1 Number of the first pixel
2318// @var pix2 Number of the second pixel
2319// @return TRUE: both pixels are neighbours; FALSE: oth.
2320//
2321// @date Wed Sep 9 17:58:37 MET DST 1998
2322//------------------------------------------------------------
2323// @function
2324
2325//!@{
2326int
2327pixels_are_neig(int pix1, int pix2)
2328{
2329 if ( sqrt(SQR( pixary[pix1][0] - pixary[pix2][0] ) +
2330 SQR( pixary[pix1][1] - pixary[pix2][1] ) )
2331 > ct_PixelWidth_corner_2_corner )
2332 return ( FALSE );
2333 else
2334 return ( TRUE );
2335}
2336//!@}
2337
2338//!-----------------------------------------------------------
2339// @name igen_pixel_coordinates
2340//
2341// @desc generate the pixel center coordinates
2342//
2343// @var *pcam structure camera containing all the
2344// camera information
2345// @return total number of pixels
2346//
2347// DP
2348//
2349// @date Thu Oct 14 10:41:03 CEST 1999
2350//------------------------------------------------------------
2351// @function
2352
2353//!@{
2354/******** igen_pixel_coordinates() *********************************/
2355
2356int igen_pixel_coordinates(struct camera *pcam) {
2357 /* generate pixel coordinates, return value is number of pixels */
2358
2359 int i, itot_inside_ring, iN, in, ipixno, iring_no, ipix_in_ring, isegment;
2360 float fsegment_fract;
2361 double dtsize;
2362 double dhsize;
2363 double dpsize;
2364 double dxfirst_pix;
2365 double dyfirst_pix;
2366 double ddxseg1, ddxseg2, ddxseg3, ddxseg4, ddxseg5, ddxseg6;
2367 double ddyseg1, ddyseg2, ddyseg3, ddyseg4, ddyseg5, ddyseg6;
2368
2369
2370 double dstartx, dstarty; /* for the gap pixels and outer pixels */
2371 int j, nrow;
2372
2373 dpsize = pcam->dpixdiameter_cm;
2374 dtsize = dpsize * sqrt(3.) / 2.;
2375 dhsize = dpsize / 2.;
2376
2377 /* Loop over central pixels to generate co-ordinates */
2378
2379 for(ipixno=1; ipixno <= pcam->inumcentralpixels; ipixno++){
2380
2381 /* Initialise variables. The central pixel = ipixno 1 in ring iring_no 0 */
2382
2383 pcam->dpixsizefactor[ipixno-1] = 1.;
2384
2385 in = 0;
2386
2387 i = 0;
2388 itot_inside_ring = 0;
2389 iring_no = 0;
2390
2391 /* Calculate the number of pixels out to and including the ring containing pixel number */
2392 /* ipixno e.g. for pixel number 17 in ring number 2 itot_inside_ring = 19 */
2393
2394 while (itot_inside_ring == 0){
2395
2396 iN = 3*(i*(i+1)) + 1;
2397
2398 if (ipixno <= iN){
2399 iring_no = i;
2400 itot_inside_ring = iN;
2401 }
2402
2403 i++;
2404 }
2405
2406
2407 /* Find the number of pixels which make up ring number iring_no e.g. ipix_in_ring = 6 for ring 1 */
2408
2409 ipix_in_ring = 0;
2410 for (i = 0; i < iring_no; ++i){
2411
2412 ipix_in_ring = ipix_in_ring + 6;
2413 }
2414
2415 /* The camera is viewed as 6 radial segments ("pie slices"). Knowing the number of pixels in its */
2416 /* ring calculate which segment the pixel ipixno is in. Then find how far across this segment it is */
2417 /* as a fraction of the number of pixels in this sixth of the ring (ask SMB). */
2418
2419 isegment = 0;
2420 fsegment_fract = 0.;
2421 if (iring_no > 0) {
2422
2423 isegment = (int)((ipixno - itot_inside_ring + ipix_in_ring - 0.5) / iring_no + 1); /* integer division ! numbering starts at 1 */
2424
2425 fsegment_fract = (ipixno - (itot_inside_ring - ipix_in_ring)) - ((isegment-1)*iring_no) - 1 ;
2426
2427 }
2428
2429 /* the first pixel in each ring lies on the positive x axis at a distance dxfirst_pix = iring_no * the */
2430 /* pixel width (flat to flat) dpsize. */
2431
2432 dxfirst_pix = dpsize*iring_no;
2433 dyfirst_pix = 0.;
2434
2435 /* the vector between the first and last pixels in a segment n is (ddxsegn, ddysegn) */
2436
2437 ddxseg1 = - dhsize*iring_no;
2438 ddyseg1 = dtsize*iring_no;
2439 ddxseg2 = -dpsize*iring_no;
2440 ddyseg2 = 0.;
2441 ddxseg3 = ddxseg1;
2442 ddyseg3 = -ddyseg1;
2443 ddxseg4 = -ddxseg1;
2444 ddyseg4 = -ddyseg1;
2445 ddxseg5 = -ddxseg2;
2446 ddyseg5 = 0.;
2447 ddxseg6 = -ddxseg1;
2448 ddyseg6 = ddyseg1;
2449
2450 /* to find the position of pixel ipixno take the position of the first pixel in the ring and move */
2451 /* anti-clockwise around the ring by adding the segment to segment vectors. */
2452
2453 switch (isegment) {
2454
2455 case 0:
2456
2457 pcam->dxc[ipixno-1] = 0.;
2458 pcam->dyc[ipixno-1] = 0.;
2459
2460 case 1:
2461 pcam->dxc[ipixno-1] = dxfirst_pix - dhsize*fsegment_fract;
2462 pcam->dyc[ipixno-1] = dyfirst_pix + dtsize*fsegment_fract;
2463
2464 break;
2465
2466 case 2:
2467
2468 pcam->dxc[ipixno-1] = dxfirst_pix + ddxseg1 - dpsize*fsegment_fract;
2469 pcam->dyc[ipixno-1] = dyfirst_pix + ddyseg1 + 0.;
2470
2471 break;
2472
2473 case 3:
2474
2475 pcam->dxc[ipixno-1] = dxfirst_pix + ddxseg1 + ddxseg2 - dhsize*fsegment_fract;
2476 pcam->dyc[ipixno-1] = dyfirst_pix + ddyseg1 + ddyseg2 - dtsize*fsegment_fract;
2477
2478 break;
2479
2480 case 4:
2481
2482 pcam->dxc[ipixno-1] = dxfirst_pix + ddxseg1 + ddxseg2 + ddxseg3 + dhsize*fsegment_fract;
2483 pcam->dyc[ipixno-1] = dyfirst_pix + ddyseg1 + ddyseg2 + ddyseg3 - dtsize*fsegment_fract;
2484
2485 break;
2486
2487 case 5:
2488
2489 pcam->dxc[ipixno-1] = dxfirst_pix + ddxseg1 + ddxseg2 + ddxseg3 + ddxseg4 + dpsize*fsegment_fract;
2490 pcam->dyc[ipixno-1] = dyfirst_pix + ddyseg1 + ddyseg2 + ddyseg3 + ddyseg4 + 0.;
2491
2492 break;
2493
2494 case 6:
2495
2496 pcam->dxc[ipixno-1] = dxfirst_pix + ddxseg1 + ddxseg2 + ddxseg3 + ddxseg4 + ddxseg5 + dhsize*fsegment_fract;
2497 pcam->dyc[ipixno-1] = dyfirst_pix + ddyseg1 + ddyseg2 + ddyseg3 + ddyseg4 + ddyseg5 + dtsize*fsegment_fract;
2498
2499 break;
2500
2501 default:
2502
2503 fprintf(stderr, "ERROR: problem in coordinate generation for pixel %d\n", ipixno);
2504 return(0);
2505
2506 } /* end switch */
2507
2508 } /* end for */
2509
2510 dstartx = pcam->dxc[pcam->inumcentralpixels - 1] + dhsize;
2511 dstarty = pcam->dyc[pcam->inumcentralpixels - 1] + dtsize;
2512
2513 if(pcam->inumgappixels > 0){ /* generate the positions of the gap pixels */
2514
2515 j = pcam->inumcentralpixels;
2516
2517 for(i=0; i<pcam->inumgappixels; i=i+6){
2518 pcam->dxc[j + i ] = dstartx + 2. * (i/6 + 1) * dpsize;
2519 pcam->dyc[j + i ] = dstarty;
2520 pcam->dpixsizefactor[j + i] = 1.;
2521 pcam->dxc[j + i + 1] = pcam->dxc[j + i ] / 2.;
2522 pcam->dyc[j + i + 1] = sqrt(3.) * pcam->dxc[j + i + 1];
2523 pcam->dpixsizefactor[j + i + 1] = 1.;
2524 pcam->dxc[j + i + 2] = - pcam->dxc[j + i + 1];
2525 pcam->dyc[j + i + 2] = pcam->dyc[j + i + 1];
2526 pcam->dpixsizefactor[j + i+ 2] = 1.;
2527 pcam->dxc[j + i + 3] = - pcam->dxc[j + i];
2528 pcam->dyc[j + i + 3] = dstarty;
2529 pcam->dpixsizefactor[j + i+ 3] = 1.;
2530 pcam->dxc[j + i + 4] = pcam->dxc[j + i + 2];
2531 pcam->dyc[j + i + 4] = - pcam->dyc[j + i + 2];
2532 pcam->dpixsizefactor[j + i+ 4] = 1.;
2533 pcam->dxc[j + i + 5] = pcam->dxc[j + i + 1];
2534 pcam->dyc[j + i + 5] = - pcam->dyc[j + i + 1];
2535 pcam->dpixsizefactor[j + i + 5] = 1.;
2536 } /* end for */
2537 } /* end if */
2538
2539 /* generate positions of the outer pixels */
2540
2541 if( pcam->inumbigpixels > 0 ){
2542
2543 j = pcam->inumcentralpixels + pcam->inumgappixels;
2544
2545 for(i=0; i<pcam->inumbigpixels; i++){
2546 pcam->dpixsizefactor[j + i] = 2.;
2547 }
2548
2549 in = 0;
2550
2551 nrow = (int) ceil(dstartx / 2. / dpsize);
2552
2553 while(in < pcam->inumbigpixels){
2554
2555 pcam->dxc[j + in] = dstartx + dpsize;
2556 pcam->dyc[j + in] = dstarty + 2 * dpsize / sqrt(3.);
2557 pcam->dxc[j + in + nrow] = dstartx / 2. - dpsize / 2.;
2558 pcam->dyc[j + in + nrow] = sqrt(3.)/2. * dstartx + 2.5 * dpsize/sqrt(3.);
2559 pcam->dxc[j + in + 3 * nrow - 1] = - pcam->dxc[j + in];
2560 pcam->dyc[j + in + 3 * nrow - 1] = pcam->dyc[j + in];
2561 pcam->dxc[j + in + 3 * nrow] = - pcam->dxc[j + in];
2562 pcam->dyc[j + in + 3 * nrow] = - pcam->dyc[j + in];
2563 pcam->dxc[j + in + 5 * nrow - 1] = pcam->dxc[j + in + nrow];
2564 pcam->dyc[j + in + 5 * nrow - 1] = - pcam->dyc[j + in + nrow];
2565 pcam->dxc[j + in + 6 * nrow - 1] = pcam->dxc[j + in];
2566 pcam->dyc[j + in + 6 * nrow - 1] = - pcam->dyc[j + in];
2567 for(i=1; i<nrow; i++){
2568 pcam->dxc[j + in + i] = pcam->dxc[j + in] - i * dpsize;
2569 pcam->dyc[j + in + i] = pcam->dyc[j + in] + i * dpsize * sqrt(3.);
2570 pcam->dxc[j + in + i + nrow] = pcam->dxc[j + in + nrow] - i * 2 * dpsize;
2571 pcam->dyc[j + in + i + nrow] = pcam->dyc[j + in + nrow];
2572 pcam->dxc[j + in + 3 * nrow - 1 - i] = - pcam->dxc[j + in + i];
2573 pcam->dyc[j + in + 3 * nrow - 1- i] = pcam->dyc[j + in + i];
2574 pcam->dxc[j + in + i + 3 * nrow] = - pcam->dxc[j + in + i];
2575 pcam->dyc[j + in + i + 3 * nrow] = - pcam->dyc[j + in + i];
2576 pcam->dxc[j + in + 5 * nrow - 1 - i] = pcam->dxc[j + in + i + nrow];
2577 pcam->dyc[j + in + 5 * nrow - 1 - i] = - pcam->dyc[j + in + i + nrow];
2578 pcam->dxc[j + in + 6 * nrow - 1 - i] = pcam->dxc[j + in + i];
2579 pcam->dyc[j + in + 6 * nrow - 1 - i] = - pcam->dyc[j + in + i];
2580 }
2581 in = in + 6 * nrow;
2582 dstartx = dstartx + 2. * dpsize;
2583 nrow = nrow + 1;
2584 } /* end while */
2585
2586 } /* end if */
2587
2588 return(pcam->inumpixels);
2589
2590}
2591//!@}
2592
2593//!-----------------------------------------------------------
2594// @name bpoint_is_in_pix
2595//
2596// @desc check if a point (x,y) in camera coordinates is inside a given pixel
2597//
2598// @var *pcam structure camera containing all the
2599// camera information
2600// @var dx, dy point coordinates in centimeters
2601// @var ipixnum pixel number (starting at 0)
2602// @return TRUE if the point is inside the pixel, FALSE otherwise
2603//
2604// DP
2605//
2606// @date Thu Oct 14 16:59:04 CEST 1999
2607//------------------------------------------------------------
2608// @function
2609
2610//!@{
2611
2612/******** bpoint_is_in_pix() ***************************************/
2613
2614int bpoint_is_in_pix(double dx, double dy, int ipixnum, struct camera *pcam){
2615 /* return TRUE if point (dx, dy) is in pixel number ipixnum, else return FALSE (use camera coordinate system) */
2616 /* the pixel is assumed to be a "closed set" */
2617
2618 double a, b; /* a = length of one of the edges of one pixel, b = half the width of one pixel */
2619 double c, xx, yy; /* auxiliary variable */
2620
2621 b = pcam->dpixdiameter_cm / 2. * pcam->dpixsizefactor[ipixnum];
2622 a = pcam->dpixdiameter_cm / sqrt(3.) * pcam->dpixsizefactor[ipixnum];
2623 c = 1./sqrt(3.);
2624
2625 if((ipixnum < 0)||(ipixnum >= pcam->inumpixels)){
2626 fprintf(stderr, "Error in bpoint_is_in_pix: invalid pixel number %d\n", ipixnum);
2627 fprintf(stderr, "Exiting.\n");
2628 exit(203);
2629 }
2630 xx = dx - pcam->dxc[ipixnum];
2631 yy = dy - pcam->dyc[ipixnum];
2632
2633 if(((-b <= xx) && (xx <= 0.) && ((-c * xx - a) <= yy) && (yy <= ( c * xx + a))) ||
2634 ((0. < xx) && (xx <= b ) && (( c * xx - a) <= yy) && (yy <= (-c * xx + a))) ){
2635 return(TRUE); /* point is inside */
2636 }
2637 else{
2638 return(FALSE); /* point is outside */
2639 }
2640}
2641
2642//!@}
2643
2644//------------------------------------------------------------
2645// @name dist_r_P
2646//
2647// @desc distance straight line r - point P
2648//
2649// @date Sat Jun 27 05:58:56 MET DST 1998
2650// @function @code
2651//------------------------------------------------------------
2652// dist_r_P
2653//
2654// distance straight line r - point P
2655//------------------------------------------------------------
2656
2657float
2658dist_r_P(float a, float b, float c,
2659 float l, float m, float n,
2660 float x, float y, float z)
2661{
2662 return (
2663 sqrt((SQR((a-x)*m-(b-y)*l) +
2664 SQR((b-y)*n-(c-z)*m) +
2665 SQR((c-z)*l-(a-x)*n))/
2666 (SQR(l)+SQR(m)+SQR(n))
2667 )
2668 );
2669}
2670
2671//------------------------------------------------------------
2672// @name check_reflector_file
2673//
2674// @desc check if a given reflector file has the right signature
2675// @desc return TRUE or FALSE
2676//
2677// @date Mon Feb 14 16:44:21 CET 2000
2678// @function @code
2679//------------------------------------------------------------
2680
2681int check_reflector_file(FILE *infile){
2682
2683 char Signature[20]; // auxiliary variable
2684 char sign[20]; // auxiliary variable
2685
2686 strcpy(Signature, REFL_SIGNATURE);
2687
2688 strcpy(sign, Signature);
2689
2690 fread( (char *)sign, strlen(Signature), 1, infile);
2691
2692 if (strcmp(sign, Signature) != 0) {
2693 cout << "ERROR: Signature of .rfl file is not correct\n";
2694 cout << '"' << sign << '"' << '\n';
2695 cout << "should be: " << Signature << '\n';
2696 return(FALSE);
2697 }
2698
2699 fread( (char *)sign, 1, 1, infile);
2700
2701 return(TRUE);
2702
2703}
2704
2705//------------------------------------------------------------
2706// @name lin_interpol
2707//
2708// @desc interpolate linearly between two points returning the
2709// @desc y-value of the result
2710//
2711// @date Thu Feb 17 11:31:32 CET 2000
2712// @function @code
2713//------------------------------------------------------------
2714
2715float lin_interpol( float x1, float y1, float x2, float y2, float x){
2716
2717 if( (x2 - x1)<=0. ){ // avoid division by zero, return average
2718 cout << "Warning: lin_interpol was asked to interpolate between two points with x1>=x2.\n";
2719 return((y1+y2)/2.);
2720 }
2721 else{ // note: the check whether x1 < x < x2 is omitted for speed reasons
2722 return((float) (y1 + (y2-y1)/(x2-x1)*(x-x1)) );
2723 }
2724}
2725
2726
2727//------------------------------------------------------------
2728// @name produce_phes
2729//
2730// @desc read the photons of an event, pixelize them and simulate QE
2731// @desc return various statistics and the array of Photoelectrons
2732//
2733// @date Mon Feb 14 16:44:21 CET 2000
2734// @function @code
2735//------------------------------------------------------------
2736
2737int produce_phes( FILE *sp, // the input file
2738 struct camera *cam, // the camera layout
2739 float minwl_nm, // the minimum accepted wavelength
2740 float maxwl_nm, // the maximum accepted wavelength
2741 class MTrigger *trigger, // the generated phes
2742 class MFadc *fadc,
2743 int *itotnphe, // total number of produced photoelectrons
2744 float nphe[iMAXNUMPIX], // number of photoelectrons in each pixel
2745 int *incph, // total number of cph read
2746 float *tmin_ns, // minimum arrival time of all phes
2747 float *tmax_ns // maximum arrival time of all phes
2748 ){
2749
2750 static int i, k, ipixnum;
2751 static float cx, cy, wl, qe, t;
2752 static MCCphoton photon;
2753 static float **qept;
2754 static char flag[SIZE_OF_FLAGS + 1];
2755 static float radius;
2756
2757
2758 // reset variables
2759
2760 for ( i=0; i<cam->inumpixels; ++i ){
2761
2762 nphe[i] = 0.0;
2763
2764 }
2765
2766 *incph = 0;
2767
2768 radius = cam->dxc[cam->inumpixels-1]
2769 + 1.5*cam->dpixdiameter_cm*cam->dpixsizefactor[cam->inumpixels-1];
2770
2771 //- - - - - - - - - - - - - - - - - - - - - - - - -
2772 // read photons and "map" them into the pixels
2773 //--------------------------------------------------
2774
2775 // initialize CPhoton
2776
2777 photon.fill(0., 0., 0., 0., 0., 0., 0., 0.);
2778
2779 // read the photons data
2780
2781 fread ( flag, SIZE_OF_FLAGS, 1, sp );
2782
2783 // loop over the photons
2784
2785 while ( !isA( flag, FLAG_END_OF_EVENT ) ) {
2786
2787 memcpy( (char*)&photon, flag, SIZE_OF_FLAGS );
2788
2789 fread( ((char*)&photon)+SIZE_OF_FLAGS, photon.mysize()-SIZE_OF_FLAGS, 1, sp );
2790
2791 // increase number of photons
2792
2793 (*incph)++;
2794
2795 // Chceck if photon is inside trigger time range
2796
2797 t = photon.get_t() ;
2798
2799 if (t-*tmin_ns>TOTAL_TRIGGER_TIME) {
2800
2801 // read next Photon
2802
2803 fread( flag, SIZE_OF_FLAGS, 1, sp );
2804
2805 // go to beginning of loop, the photon is lost
2806 continue;
2807 }
2808
2809 //
2810 // Pixelization
2811 //
2812
2813 cx = photon.get_x();
2814 cy = photon.get_y();
2815
2816 // get wavelength
2817
2818 wl = photon.get_wl();
2819
2820 // cout << "wl " << wl << " radius " << sqrt(cx*cx + cy*cy) << "\n";
2821
2822 // check if photon has valid wavelength and is inside outermost camera radius
2823
2824 if( (wl > maxwl_nm) || (wl < minwl_nm) || (sqrt(cx*cx + cy*cy) > radius ) ){
2825
2826 // cout << " lost first check\n";
2827
2828 // read next CPhoton
2829
2830 fread ( flag, SIZE_OF_FLAGS, 1, sp );
2831
2832 // go to beginning of loop, the photon is lost
2833 continue;
2834
2835 }
2836
2837 ipixnum = -1;
2838
2839 for(i=0; i<cam->inumpixels; i++){
2840 if( bpoint_is_in_pix( cx, cy, i, cam) ){
2841 ipixnum = i;
2842 break;
2843 }
2844 }
2845
2846 if(ipixnum==-1){// the photon is in none of the pixels
2847
2848 // cout << " lost pixlization\n";
2849
2850 // read next CPhoton
2851
2852 fread ( flag, SIZE_OF_FLAGS, 1, sp );
2853
2854 // go to beginning of loop, the photon is lost
2855 continue;
2856 }
2857
2858 if(ipixnum==0) {// the phton is in the central pixel, which is not used for trigger
2859 // read next CPhoton
2860
2861 fread ( flag, SIZE_OF_FLAGS, 1, sp );
2862
2863 // go to beginning of loop, the photon is lost
2864 continue;
2865 }
2866 //+++
2867 // QE simulation
2868 //---
2869
2870 // set pointer to the QE table of the relevant pixel
2871
2872 qept = (float **)QE[ipixnum];
2873
2874 // check if wl is inside table; outside the table, QE is assumed to be zero
2875
2876 if((wl < qept[0][0]) || (wl > qept[0][pointsQE-1])){
2877
2878 // cout << " lost\n";
2879
2880 // read next Photon
2881
2882 fread ( flag, SIZE_OF_FLAGS, 1, sp );
2883
2884 // go to beginning of loop
2885 continue;
2886
2887 }
2888
2889 // find data point in the QE table (-> k)
2890
2891 k = 1; // start at 1 because the condition was already tested for 0
2892 while (k < pointsQE-1 && qept[0][k] < wl){
2893 k++;
2894 }
2895
2896 // calculate the qe value between 0. and 1.
2897
2898 qe = lin_interpol(qept[0][k-1], qept[1][k-1], qept[0][k], qept[1][k], wl) / 100.0;
2899
2900 // if random > quantum efficiency, reject it
2901
2902 if ( (RandomNumber) > qe ) {
2903
2904 // cout << " lost\n";
2905
2906 // read next Photon
2907
2908 fread ( flag, SIZE_OF_FLAGS, 1, sp );
2909
2910 // go to beginning of loop
2911 continue;
2912
2913 }
2914
2915 //+++
2916 // The photon has produced a photo electron
2917 //---
2918
2919 // cout << " accepted\n";
2920
2921 // increment the number of photoelectrons in the relevant pixel
2922
2923 nphe[ipixnum] += 1.0;
2924
2925 // store the new photoelectron
2926
2927 fadc->Fill(ipixnum,(t-*tmin_ns) , trigger->FillShow(ipixnum,t-*tmin_ns));
2928
2929 *itotnphe += 1;
2930
2931 // read next Photon
2932
2933 fread( flag, SIZE_OF_FLAGS, 1, sp );
2934
2935 } // end while, i.e. found end of event
2936
2937 return(0);
2938
2939}
2940
2941
2942//------------------------------------------------------------
2943// @name produce_nsbrates
2944//
2945// @desc read the starfield file, call produce_phes on it in,
2946// @desc different wavebands, calculate the nsbrates
2947//
2948// @date Mon Feb 14 16:44:21 CET 2000
2949// @function @code
2950//------------------------------------------------------------
2951
2952int produce_nsbrates( char *iname, // the starfield input file name
2953 struct camera *cam, // camera layout
2954 float rate_phepns[iMAXNUMPIX][iNUMWAVEBANDS] // the product of this function:
2955 // the NSB rates in phe/ns for each pixel
2956 ){
2957
2958 int i, j, k, ii; // counters
2959
2960 MTrigger trigger(Trigger_gate_length, Trigger_overlaping_time, Trigger_response_ampl, Trigger_response_fwhm);
2961 MFadc flashadc;
2962
2963 static float wl_nm[iNUMWAVEBANDS + 1] = { WAVEBANDBOUND1,
2964 WAVEBANDBOUND2,
2965 WAVEBANDBOUND3,
2966 WAVEBANDBOUND4,
2967 WAVEBANDBOUND5,
2968 WAVEBANDBOUND6 };
2969
2970 FILE *infile; // the input file
2971 fpos_t fileposition; // marker on the input file
2972 static char cflag[SIZE_OF_FLAGS + 1]; // auxiliary variable
2973 static MCEventHeader evth; // the event header
2974 static float nphe[iMAXNUMPIX]; // the number of photoelectrons in each pixel
2975 int itnphe; // total number of produced photoelectrons
2976 int itotnphe; // total number of produced photoelectrons after averaging
2977 int incph; // total number of cph read
2978 float tmin_ns; // minimum arrival time of all phes
2979 float tmax_ns; // maximum arrival time of all phes
2980 float integtime_ns; // integration time from the starfield generator
2981
2982 // open input file
2983
2984 log(SIGNATURE, "Opening starfield input \"rfl\" file %s\n", iname );
2985
2986 infile = fopen( iname, "r" );
2987
2988 // check if the satrfield input file exists
2989
2990 if ( infile == NULL ) {
2991
2992 log( SIGNATURE, "Cannot open starfield input file: %s\n", iname );
2993
2994 log( SIGNATURE, "There is not NSB from the Stars\n");
2995
2996 return (0);
2997 }
2998
2999 // get signature, and check it
3000
3001 if(check_reflector_file(infile)==FALSE){
3002 exit(1);
3003 }
3004
3005 // initialize flag
3006
3007 strcpy( cflag, " \0" );
3008
3009 // get flag
3010
3011 fread( cflag, SIZE_OF_FLAGS, 1, infile );
3012
3013 if ( ! feof(infile)){
3014
3015 // reading .rfl file
3016
3017 if(!isA( cflag, FLAG_START_OF_RUN )){
3018 error( SIGNATURE, "Expected start of run flag, but found: %s\n", cflag );
3019 }
3020 else { // found start of run
3021
3022 fread( cflag, SIZE_OF_FLAGS, 1, infile );
3023
3024 if( isA( cflag, FLAG_START_OF_EVENT )){ // there is a event
3025
3026 // get MCEventHeader
3027
3028 fread( (char*)&evth, evth.mysize(), 1, infile );
3029
3030 integtime_ns = evth.get_energy();
3031
3032 // memorize where we are in the file
3033
3034 if (fgetpos( infile, &fileposition ) != 0){
3035 error( SIGNATURE, "Cannot position in file ...\n");
3036 }
3037
3038 // loop over the wavebands
3039
3040 for(i=0; i<iNUMWAVEBANDS; i++){
3041
3042 // initialize the rate array
3043
3044 for(j = 0; j<cam->inumpixels; j++){ // loop over pixels
3045 rate_phepns[j][i] = 0.;
3046 }
3047
3048 itotnphe = 0;
3049
3050 // read the photons and produce the photoelectrons
3051 // - in order to average over the QE simulation, call the
3052 // production function iNUMNSBPRODCALLS times
3053
3054 for(ii=0; ii<iNUMNSBPRODCALLS; ii++){
3055
3056 // position the file pointer to the beginning of the photons
3057
3058 fsetpos( infile, &fileposition);
3059
3060 // produce photoelectrons
3061
3062 k = produce_phes( infile,
3063 cam,
3064 wl_nm[i],
3065 wl_nm[i+1],
3066 &trigger, // this is a dummy here
3067 &flashadc, // this is a dummy here
3068 &itnphe,
3069 nphe, // we want this!
3070 &incph,
3071 &tmin_ns,
3072 &tmax_ns );
3073
3074 if( k != 0 ){ // non-zero returnvalue means error
3075 cout << "Exiting.\n";
3076 exit(1);
3077 }
3078
3079 for(j = 0; j<cam->inumpixels; j++){ // loop over pixels
3080 rate_phepns[j][i] += nphe[j]/integtime_ns/(float)iNUMNSBPRODCALLS;
3081 }
3082
3083 itotnphe += itnphe;
3084
3085 } // end for(ii=0 ...
3086
3087 fprintf(stdout, "Starfield, %6f - %6f nm: %d photoelectrons for %6f ns integration time\n",
3088 wl_nm[i], wl_nm[i+1], itotnphe/iNUMNSBPRODCALLS, integtime_ns);
3089
3090 } // end for(i=0 ...
3091
3092 }
3093 else{
3094 cout << "Starfield file contains no event.\nExiting.\n";
3095 exit(1);
3096 } // end if( isA ... event
3097 } // end if ( isA ... run
3098 }
3099 else{
3100 cout << "Starfield file contains no run.\nExiting.\n";
3101 exit(1);
3102 }
3103
3104 fclose( infile );
3105
3106 return(0);
3107
3108}
3109
3110
3111//------------------------------------------------------------
3112// @name produce_nsb_phes
3113//
3114// @desc produce the photoelectrons from the nsbrates
3115//
3116// @date Thu Feb 17 17:10:40 CET 2000
3117// @function @code
3118//------------------------------------------------------------
3119
3120int produce_nsb_phes( float *atmin_ns,
3121 float *atmax_ns,
3122 float theta_rad,
3123 struct camera *cam,
3124 float nsbr_phepns[iMAXNUMPIX][iNUMWAVEBANDS],
3125 float difnsbr_phepns[iMAXNUMPIX],
3126 float extinction[iNUMWAVEBANDS],
3127 float fnpx[iMAXNUMPIX],
3128 MTrigger *trigger,
3129 MFadc *fadc,
3130 int *inphe,
3131 float base_mv[iMAXNUMPIX]){
3132
3133 float simtime_ns; // NSB simulation time
3134 int i, j, k, ii;
3135 float zenfactor; // correction factor calculated from the extinction
3136 int inumnsbphe; // number of photoelectrons caused by NSB
3137
3138 float t;
3139 TRandom random; // Random numbers generator
3140
3141 random.SetSeed ((UInt_t) (RandomNumber*1000));
3142
3143 ii = *inphe; // avoid dereferencing
3144
3145 // check if the arrival times are set; if not generate them
3146
3147 if(*atmin_ns <SIMTIMEOFFSET_NS || *atmin_ns > *atmax_ns){
3148 *atmin_ns = 0.;
3149 *atmax_ns = simtime_ns = TOTAL_TRIGGER_TIME;
3150
3151 }
3152 else{ // extend the event time window by the given offsets
3153
3154 *atmin_ns = *atmin_ns - SIMTIMEOFFSET_NS;
3155 *atmax_ns = *atmax_ns + SIMTIMEOFFSET_NS;
3156
3157 simtime_ns = *atmax_ns - *atmin_ns;
3158
3159 // make sure the simulated time is long enough for the FADC
3160 // simulation and not too long
3161
3162 if(simtime_ns< TOTAL_TRIGGER_TIME){
3163 *atmin_ns = *atmin_ns -(TOTAL_TRIGGER_TIME-simtime_ns)/2;
3164 *atmax_ns = *atmin_ns + TOTAL_TRIGGER_TIME;
3165 simtime_ns = TOTAL_TRIGGER_TIME;
3166 }
3167
3168 if(simtime_ns> TOTAL_TRIGGER_TIME){
3169 *atmax_ns =*atmin_ns + TOTAL_TRIGGER_TIME;
3170 simtime_ns = TOTAL_TRIGGER_TIME;
3171 }
3172
3173 }
3174
3175 // initialize baselines
3176
3177 for(i=0; i<cam->inumpixels; i++){
3178 base_mv[i] = 0.;
3179 }
3180
3181 // calculate baselines and generate phes
3182
3183 for(i=0; i<iNUMWAVEBANDS; i++){ // loop over the wavebands
3184
3185 // calculate the effect of the atmospheric extinction
3186
3187 zenfactor = pow(10., -0.4 * extinction[i]/cos(theta_rad) );
3188
3189 for(j=0; j<cam->inumpixels; j++){ // loop over the pixels
3190
3191 inumnsbphe = (int) ((zenfactor * nsbr_phepns[j][i] + difnsbr_phepns[j]/iNUMWAVEBANDS)
3192 * simtime_ns );
3193
3194 base_mv[j] += inumnsbphe;
3195
3196 // randomize
3197
3198 if (inumnsbphe>0.0){
3199 inumnsbphe = random.Poisson (inumnsbphe );
3200 }
3201
3202 // create the photoelectrons
3203
3204 for(k=0; k < inumnsbphe; k++){
3205
3206 t=(RandomNumber * simtime_ns);
3207
3208 (*fadc).Fill(j,t ,(*trigger).FillNSB(j,t));
3209
3210 ii++; // increment total number of photoelectons
3211
3212 fnpx[j] += 1.; // increment number of photoelectrons in this pixel
3213
3214 }
3215
3216 } // end for(j=0 ...
3217 } // end for(i=0 ...
3218
3219 // finish baseline calculation
3220
3221 for(i=0; i<cam->inumpixels; i++){
3222 base_mv[i] *= RESPONSE_FWHM * RESPONSE_AMPLITUDE / simtime_ns;
3223 }
3224
3225 *inphe = ii; // update the pointer
3226
3227 return(0);
3228}
3229
3230
3231// @endcode
3232
3233
3234//=------------------------------------------------------------
3235//!@subsection Log of this file.
3236
3237//!@{
3238//
3239// $Log: not supported by cvs2svn $
3240// Revision 1.24 2001/04/06 16:48:09 magicsol
3241// New camera version able to read the new format of the reflector output:
3242// reflector 0.4
3243//
3244// Revision 1.23 2001/03/28 16:13:41 blanch
3245// While feeling the MMcFadeHeader branch the boolean conditoin was wrong. It has
3246// been solved.
3247//
3248// Revision 1.22 2001/03/20 18:52:43 blanch
3249// *** empty log message ***
3250//
3251// Revision 1.21 2001/03/19 19:53:03 blanch
3252// Some print outs have been removed.
3253//
3254// Revision 1.20 2001/03/19 19:30:06 magicsol
3255// Minor changes have been done to improve the FADC pedestals treatment.
3256//
3257// Revision 1.19 2001/03/05 11:14:41 magicsol
3258// I changed the position of readinf a parameter. It is a minnor change.
3259//
3260// Revision 1.18 2001/03/05 10:36:52 blanch
3261// A branch with information about the FADC simulation (MMcFadcHeader) is writen
3262// in the McHeaders tree of the aoutput root file.
3263// The NSB contribution is only applied if the the number of phe form the shower
3264// are above a given number.
3265//
3266// Revision 1.17 2001/02/23 11:05:57 magicsol
3267// Small changes due to slightly different output format and the introduction of
3268// pedesals for teh FADC.
3269//
3270// Revision 1.16 2001/01/18 18:44:40 magicsol
3271// *** empty log message ***
3272//
3273// Revision 1.15 2001/01/17 09:32:27 harald
3274// The changes are neccessary to have the same name for trees in MC and in
3275// data. So now there should be no differences in MC data and real data from
3276// FADC system.
3277//
3278// Revision 1.14 2001/01/15 12:33:34 magicsol
3279// Some modifications have been done to use the new (Dec'2000) Raw data format.
3280// There are also some minnors modifications to adapt the improvements in the
3281// MTrigger class (overlaping time and trigger cells).
3282//
3283// Revision 1.13 2000/10/25 08:14:23 magicsol
3284// The routine that produce poisson random numbers to decide how many phe
3285// form NSB are emmited in each pixel has been replaced. Now a ROOT routine
3286// is used.
3287//
3288// Revision 1.12 2000/09/22 17:40:18 harald
3289// Added a lot of changes done by oscar.
3290//
3291// Revision 1.11 2000/09/21 11:47:33 harald
3292// Oscar found some smaller errors in the calculation of the pixel shape and
3293// corrected it.
3294//
3295// $Log: not supported by cvs2svn $
3296// Revision 1.24 2001/04/06 16:48:09 magicsol
3297// New camera version able to read the new format of the reflector output:
3298// reflector 0.4
3299//
3300// Revision 1.23 2001/03/28 16:13:41 blanch
3301// While feeling the MMcFadeHeader branch the boolean conditoin was wrong. It has
3302// been solved.
3303//
3304// Revision 1.22 2001/03/20 18:52:43 blanch
3305// *** empty log message ***
3306//
3307// Revision 1.21 2001/03/19 19:53:03 blanch
3308// Some print outs have been removed.
3309//
3310// Revision 1.20 2001/03/19 19:30:06 magicsol
3311// Minor changes have been done to improve the FADC pedestals treatment.
3312//
3313// Revision 1.19 2001/03/05 11:14:41 magicsol
3314// I changed the position of readinf a parameter. It is a minnor change.
3315//
3316// Revision 1.18 2001/03/05 10:36:52 blanch
3317// A branch with information about the FADC simulation (MMcFadcHeader) is writen
3318// in the McHeaders tree of the aoutput root file.
3319// The NSB contribution is only applied if the the number of phe form the shower
3320// are above a given number.
3321//
3322// Revision 1.17 2001/02/23 11:05:57 magicsol
3323// Small changes due to slightly different output format and the introduction of
3324// pedesals for teh FADC.
3325//
3326// Revision 1.16 2001/01/18 18:44:40 magicsol
3327// *** empty log message ***
3328//
3329// Revision 1.15 2001/01/17 09:32:27 harald
3330// The changes are neccessary to have the same name for trees in MC and in
3331// data. So now there should be no differences in MC data and real data from
3332// FADC system.
3333//
3334// Revision 1.14 2001/01/15 12:33:34 magicsol
3335// Some modifications have been done to use the new (Dec'2000) Raw data format.
3336// There are also some minnors modifications to adapt the improvements in the
3337// MTrigger class (overlaping time and trigger cells).
3338//
3339// Revision 1.13 2000/10/25 08:14:23 magicsol
3340// The routine that produce poisson random numbers to decide how many phe
3341// form NSB are emmited in each pixel has been replaced. Now a ROOT routine
3342// is used.
3343//
3344// Revision 1.10 2000/07/04 14:10:20 MagicSol
3345// Some changes have been done in the root output file. The RawEvt tree is only
3346// stored in the single trigger mode.
3347// The trigger input parameters are also given by the general input card.
3348// The diffuse NSB and the star NSB have been decoupled. Now the contribution of
3349// each one can be studied seperately.
3350//
3351// Revision 1.9 2000/06/13 13:25:24 blanch
3352// The multiple arrays have been replaced, since they do not work
3353// in alpha machines. Now we are using pointers and the command new
3354// to allocate memory.
3355//
3356// Revision 1.8 2000/05/11 13:57:27 blanch
3357// The option to loop over several trigger configurations has been included.
3358// Some non-sense with triggertime range has been solved.
3359// Montecarlo information and ADC counts are saved in a root file.
3360// There was a problem with the maximum number of phe that the program could analyse. Now there is not limit.
3361//
3362// Revision 1.7 2000/03/24 18:10:46 blanch
3363// A first FADC simulation and a trigger simulation are already implemented.
3364// The calculation of the Hillas Parameters have been removed, since it was decided that it should be in the analysis software.
3365// A loop over trigger threshold and some corretcions in the time range where it looks for a trigger will be implemented as soon as possible.
3366//
3367// Revision 1.6 2000/03/20 18:35:11 blanch
3368// The trigger is already implemented but it does not save the trigger information in any file as it is implemented in timecam. In the next days there will be a version which also creates the files with the trigger information. It is going to be a mixing of the current camera and timecam programs.
3369//
3370// Revision 1.5 2000/02/18 17:40:35 petry
3371// This version includes drastic changes compared to camera.cxx 1.4.
3372// It is not yet finished and not immediately useful because the
3373// trigger simulation is not yet re-implemented. I had to take it
3374// out together with some other stuff in order to tidy the whole
3375// program up. This is not meant as an insult to anyone. I needed
3376// to do this in order to be able to work on it.
3377//
3378// This version has been put in the repository in order to be
3379// able to share the further development with others.
3380//
3381// If you need something working, wait or take an earlier one.
3382// See file README.
3383//
3384// Revision 1.4 2000/01/25 08:36:23 petry
3385// The pixelization in previous versions was buggy.
3386// This is the first version with a correct pixelization.
3387//
3388// Revision 1.3 2000/01/20 18:22:17 petry
3389// Found little bug which makes camera crash if it finds a photon
3390// of invalid wavelength. This bug is now fixed and the range
3391// of valid wavelengths extended to 290 - 800 nm.
3392// This is in preparation for the NSB simulation to come.
3393// Dirk
3394//
3395// Revision 1.2 1999/11/19 08:40:42 harald
3396// Now it is possible to compile the camera programm under osf1.
3397//
3398// Revision 1.1.1.1 1999/11/05 11:59:31 harald
3399// This the starting point for CVS controlled further developments of the
3400// camera program. The program was originally written by Jose Carlos.
3401// But here you can find a "rootified" version to the program. This means
3402// that there is no hbook stuff in it now. Also the output of the
3403// program changed to the MagicRawDataFormat.
3404//
3405// The "rootification" was done by Dirk Petry and Harald Kornmayer.
3406//
3407// Revision 1.3 1999/10/22 15:01:28 petry
3408// version sent to H.K. and N.M. on Fri Oct 22 1999
3409//
3410// Revision 1.2 1999/10/22 09:44:23 petry
3411// first synthesized version which compiles and runs without crashing;
3412//
3413// Revision 1.1.1.1 1999/10/21 16:35:10 petry
3414// first synthesised version
3415//
3416// Revision 1.13 1999/03/15 14:59:05 gonzalez
3417// camera-1_1
3418//
3419// Revision 1.12 1999/03/02 09:56:10 gonzalez
3420// *** empty log message ***
3421//
3422//
3423//!@}
3424
3425//=EOF
Note: See TracBrowser for help on using the repository browser.