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

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