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

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