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

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