//!///////////////////////////////////////////////////////////////////// // // camera // // @file camera.cxx // @title Camera simulation // @subtitle Code for the simulation of the camera phase // @desc Code for the simulation of the camera of CT1 and MAGIC // @author J C Gonzalez // @email gonzalez@mppmu.mpg.de // @date Thu May 7 16:24:22 1998 // //---------------------------------------------------------------------- // // Created: Thu May 7 16:24:22 1998 // Author: Jose Carlos Gonzalez // Purpose: Program for reflector simulation // Notes: See files README for details // //---------------------------------------------------------------------- // // $RCSfile: camera.cxx,v $ // $Revision: 1.8 $ // $Author: blanch $ // $Date: 2000-05-11 13:57:27 $ // //////////////////////////////////////////////////////////////////////// // @tableofcontents @coverpage //=----------------------------------------------------------- //!@section Source code of |camera.cxx|. /*!@" In this section we show the (commented) code of the program for the read-out of the output files generated by the simulator of the reflector, |reflector 0.3|. @"*/ //=----------------------------------------------------------- //!@subsection Includes and Global variables definition. //!@{ // includes for ROOT // BEWARE: the order matters! #include "TROOT.h" #include "TApplication.h" #include "TFile.h" #include "TTree.h" #include "TBranch.h" #include "TCanvas.h" #include "MTrigger.hxx" #include "MFadc.hxx" #include "MRawEvt.h" #include "MMcEvt.h" #include "MMcTrig.hxx" /*!@" All the defines are located in the file |camera.h|. @"*/ #include "camera.h" //!@} /*!@" The following set of flags are used in time of compilation. They do not affect directly the behaviour of the program at run-time (though, of course, if you disconnected the option for implementation of the Trigger logic, you will not be able to use any trigger at all. The 'default' values mean default in the sense of what you got from the server when you obtained this program. @"*/ //!@{ // flag for debugging (default: OFF ) #define __DEBUG__ #undef __DEBUG__ //!@} //=----------------------------------------------------------- //!@subsection Definition of global variables. /*!@" Now we define some global variables with data about the telescope, such as "focal distance", number of pixels/mirrors, "size of the camera", and so on. @"*/ /*!@" Depending on the telescope we are using (CT1 or MAGIC), the information stored in the definition file is different. The variable |ct_Type| has the value 0 when we use CT1, and 1 when we use MAGIC. @"*/ //!@{ static int ct_Type; //@< Type of telescope: 0:CT1, 1:MAGIC //!@} /*!@" And this is the information about the whole telescope. @"*/ //!@{ // parameters of the CT (from the CT definition file) ////@: Focal distances [cm] //static float *ct_Focal; //@: Mean Focal distances [cm] static float ct_Focal_mean; //@: STDev. Focal distances [cm] static float ct_Focal_std; //@: Mean Point Spread function [cm] static float ct_PSpread_mean; //@: STDev. Point Spread function [cm] static float ct_PSpread_std; //@: STDev. Adjustmente deviation [cm] static float ct_Adjustment_std; //@: Radius of the Black Spot in mirror [cm] static float ct_BlackSpot_rad; //@: Radius of one mirror [cm] static float ct_RMirror; //@: Camera width [cm] static float ct_CameraWidth; //@: Pixel width [cm] static float ct_PixelWidth; //@: ct_PixelWidth_corner_2_corner = ct_PixelWidth / cos(60) static float ct_PixelWidth_corner_2_corner; //@: ct_PixelWidth_corner_2_corner / 2 static float ct_PixelWidth_corner_2_corner_half; //@: Number of mirrors static int ct_NMirrors = 0; //@: Number of pixels static int ct_NPixels; //@: Number of pixels static int ct_NCentralPixels; //@: Number of pixels static int ct_NGapPixels; //@: ct_Apot = ct_PixelWidth / 2 static float ct_Apot; //@: ct_2Apot = 2 * ct_Apot = ct_PixelWidth static float ct_2Apot; //@: name of the CT definition file to use static char ct_filename[256]; //@: list of showers to be skipped static int *Skip; //@: number of showers to be skipped static int nSkip=0; //@: flag: TRUE: data come from STDIN; FALSE: from file static int Data_From_STDIN = FALSE; //@: flag: TRUE: write all images to output; FALSE: only triggered showers static int Write_All_Images = FALSE; //@: flag: TRUE: write all data to output; FALSE: only triggered showers static int Write_All_Data = FALSE; static int Write_McEvt = TRUE; static int Write_McTrig = FALSE; static int Write_RawEvt = FALSE; //@: flag: TRUE: selection on the energy static int Select_Energy = TRUE; //@: Lower edge of the selected energy range (in GeV) static float Select_Energy_le = 0.0; //@: Upper edge of the selected energy range (in GeV) static float Select_Energy_ue = 100000.0; //@: flag: TRUE: show all fadc singnal in the screen; FALSE: don't static int FADC_Scan = FALSE; //@: flag: TRUE: show all trigger singnal in the screen; FALSE: don't static int Trigger_Scan = FALSE; //@: flag: TRUE: loop trigger analysis over several thresholds, multiplicities and topologies; FALSE: a single trigger configuration static int Trigger_Loop = FALSE; //@: Upper and lower edges of the trigger loop static int Trigger_loop_lthres = 0; static int Trigger_loop_uthres = 10; static int Trigger_loop_lmult = 2; static int Trigger_loop_umult = 10; static int Trigger_loop_ltop = 0; static int Trigger_loop_utop = 2; //!@} /*!@" The following double-pointer is a 2-dimensional table with information about each pixel. The routine read_pixels will generate the information for filling it using igen_pixel_coordinates(). @"*/ //!@{ // Pointer to a tables/Arrays with information about the pixels // and data stored on them with information about the pixels //@: coordinates x,y for each pixel static float **pixary; //@: indexes of pixels neighbours of a given one static int **pixneig; //@: number of neighbours a pixel have static int *npixneig; //@: contents of the pixels (ph.e.) static float *fnpix; //@: contents of the pixels (ph.e.) after cleanning static float *fnpixclean; //!@} /*!@" The following double-pointer is a 2-dimensional table with the Quantum Efficiency @$QE@$ of each pixel in the camera, as a function of the wavelength @$\lambda@$. The routine |read_pixels()| will read also this information from the file |qe.dat|. @"*/ //!@{ // Pointer to a table with QE, number of datapoints, and wavelengths //@: table of QE static float ***QE; //@: number of datapoints for the QE curve static int pointsQE; //@: table of QE static float *QElambda; //!@} /*!@" The following double-pointer is a 2-dimensional table with information about each mirror in the dish. The routine |read_ct_file()| will read this information from the CT definition file. @"*/ //!@{ // Pointer to a table with the following info.: static float **ct_data; /* * TYPE=0 (CT1) * i s rho theta x y z thetan phin xn yn zn * * i : number of the mirror * s : arc length [cm] * rho : polar rho of the position of the center of the mirror [cm] * theta : polar angle of the position of the center of the mirror [cm] * x : x coordinate of the center of the mirror [cm] * y : y coordinate of the center of the mirror [cm] * z : z coordinate of the center of the mirror [cm] * thetan : polar theta angle of the direction where the mirror points to * phin : polar phi angle of the direction where the mirror points to * xn : xn coordinate of the normal vector in the center (normalized) * yn : yn coordinate of the normal vector in the center (normalized) * zn : zn coordinate of the normal vector in the center (normalized) * * TYPE=1 (MAGIC) * i f sx sy x y z thetan phin * * i : number of the mirror * f : focal distance of that mirror * sx : curvilinear coordinate of mirror's center in X[cm] * sy : curvilinear coordinate of mirror's center in X[cm] * x : x coordinate of the center of the mirror [cm] * y : y coordinate of the center of the mirror [cm] * z : z coordinate of the center of the mirror [cm] * thetan : polar theta angle of the direction where the mirror points to * phin : polar phi angle of the direction where the mirror points to * xn : xn coordinate of the normal vector in the center (normalized) * yn : yn coordinate of the normal vector in the center (normalized) * zn : zn coordinate of the normal vector in the center (normalized) */ //!@} /*!@" We define a table into where random numbers will be stored. The routines used for random number generation are provided by |RANLIB| (taken from NETLIB, |www.netlib.org|), and by the routine |double drand48(void)| (prototype defined in |stdlib.h|) through the macro |RandomNumber| defined in |camera.h|. @"*/ //!@} extern char FileName[]; //=----------------------------------------------------------- // @subsection Main program. //!@{ //++++++++++++++++++++++++++++++++++++++++ // MAIN PROGRAM //---------------------------------------- int main(int argc, char **argv) { //!@' @#### Definition of variables. //@' char inname[256]; //@< input file name char starfieldname[256]; //@< starfield input file name char datname[256]; //@< data (ASCII) output file name char diagname[256]; //@< diagnistic output file (ROOT format) char rootname[256] ; //@< ROOT file name char parname[256]; //@< parameters file name char sign[20]; //@< initialize sign char flag[SIZE_OF_FLAGS + 1]; //@< flags in the .rfl file FILE *inputfile; //@< stream for the input file ofstream datafile; //@< stream for the data file MCEventHeader mcevth; //@< Event Header class (MC) MCCphoton cphoton; //@< Cherenkov Photon class (MC) int inumphe; //@< number of photoelectrons in an event float arrtmin_ns; //@ arrival time of the first photoelectron float arrtmax_ns; //@ arrival time of the last photoelectron float thetaCT, phiCT; //@< parameters of a given shower float thetashw, phishw; //@< parameters of a given shower float coreD, coreX, coreY; //@< core position and distance float impactD; //@< impact parameter float l1, m1, n1; //@< auxiliary variables float l2, m2, n2; //@< auxiliary variables float num, den; //@< auxiliary variables int nshow=0; //@< partial number of shower in a given run int ntshow=0; //@< total number of showers int ncph=0; //@< partial number of photons in a given run int ntcph=0; //@< total number of photons int i, j, k; //@< simple counters int simulateNSB; //@< Will we simulate NSB? float meanNSB; //@< diffuse NSB mean value (phe per ns per central pixel) float diffnsb_phepns[iMAXNUMPIX]; //@< diffuse NSB values for each pixel derived //@< from meanNSB float nsbrate_phepns[iMAXNUMPIX][iNUMWAVEBANDS]; //@< non-diffuse nsb //@< photoelectron rates float ext[iNUMWAVEBANDS] = { //@< average atmospheric extinction in each waveband EXTWAVEBAND1, EXTWAVEBAND2, EXTWAVEBAND3, EXTWAVEBAND4, EXTWAVEBAND5 }; float baseline_mv[iMAXNUMPIX]; //@< The baseline (mV) caused by the NSB; to be subtracted //@< in order to simulate the preamps' AC coupling float qThreshold; //@< Threshold value float qTailCut; //@< Tail Cut value int nIslandsCut; //@< Islands Cut value int countIslands; //@< Will we count the islands? int anaPixels; float fCorrection; //@< Factor to apply to pixel values (def. 1.) int trigger; //@< trigger flag int itrigger; //@< index of pixel fired int ntrigger = 0; //@< number of triggers in the whole file int ithrescount; //@< counter for loop over threshold trigger int imulticount; //@< counter for loop over multiplicity trigger int itopocount; //@< counter for loop over topology trigger float fpixelthres[TRIGGER_PIXELS]; float fadcValues[(Int_t) SLICES_MFADC]; //@< the analog Fadc siganl for pixels float plateScale_cm2deg; //@< plate scale (deg/cm) float degTriggerZone; //@< trigger area in the camera (radius, in deg.) float dtheta, dphi; //@< deviations of CT from shower axis int still_in_loop = FALSE; float *image_data; int nvar, hidt; struct camera cam; // structure holding the camera definition //!@' @#### Definition of variables for |getopt()|. //@' int ch, errflg = 0; //@< used by getopt /*!@' @#### Beginning of the program. We start with the main program. First we (could) make some presentation, and follows the reading of the parameters file (now from the |stdin|), the reading of the CT parameters file, and the creation of the output file, where the processed data will be stored. */ //++ // START //-- // make unbuffered output cout.setf ( ios::stdio ); // parse command line options (see reflector.h) parname[0] = '\0'; optarg = NULL; while ( !errflg && ((ch = getopt(argc, argv, COMMAND_LINE_OPTIONS)) != -1) ) switch (ch) { case 'f': strcpy(parname, optarg); break; case 'h': usage(); break; default : errflg++; } // show help if error if ( errflg>0 ) usage(); // make some sort of presentation present(); // read parameters file if ( strlen(parname) < 1 ) readparam(NULL); else readparam(parname); // read data from file or from STDIN? Data_From_STDIN = get_data_from_stdin(); // write all images, even those without trigger? Write_All_Images = get_write_all_images(); // write all data (i.e., ph.e.s in pixels) Write_All_Data = get_write_all_data(); Write_McEvt = get_write_McEvt() ; Write_McTrig = get_write_McTrig() ; Write_RawEvt = get_write_RawEvt() ; FADC_Scan = get_FADC_Scan(); Trigger_Scan = get_Trigger_Scan(); Trigger_Loop = get_Trigger_Loop(&Trigger_loop_lthres, &Trigger_loop_uthres, &Trigger_loop_lmult, &Trigger_loop_umult, &Trigger_loop_ltop, &Trigger_loop_utop); // get filenames strcpy( inname, get_input_filename() ); strcpy( starfieldname, get_starfield_filename() ); strcpy( datname, get_data_filename() ); strcpy( diagname, get_diag_filename() ); strcpy( rootname, get_root_filename() ); strcpy( ct_filename, get_ct_filename() ); // get different parameters of the simulation qThreshold = get_threshold(); qTailCut = get_tail_cut(); simulateNSB = get_nsb( &meanNSB ); countIslands = get_islands_cut( &nIslandsCut ); // get selections on the parameters Select_Energy = get_select_energy( &Select_Energy_le, &Select_Energy_ue); // log filenames information log(SIGNATURE, "%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", "Filenames", "In", inname, "Stars", starfieldname, "CT", ct_filename, "Data", datname, "Diag", diagname, "ROOT", rootname ); // log flags information log(SIGNATURE, "%s:\n\t%20s: %s\n\t%20s: %s\n\t%20s: %s\n", "Flags", "Data_From_STDIN", ONoff(Data_From_STDIN), "Write_All_Images", ONoff(Write_All_Images), "Write_All_Data", ONoff(Write_All_Data)); // log flags information log(SIGNATURE, "%s:\n\t%20s: %s\n\t%20s: %s\n\t%20s: %s\n", "Root output", "Write_McEvt", ONoff(Write_McEvt), "Write_McTrig", ONoff(Write_McTrig), "Write_RawEvt", ONoff(Write_RawEvt)); // log parameters information log(SIGNATURE, "%s:\n\t%20s: %f\n\t%20s: %f\n\t%20s: %f %s\n\t%20s: %f %s\n", "Parameters", "q0 (Threshold)", qThreshold, "t0 (Tail-cut)", qTailCut, "NSB (phes/pixel)", meanNSB, ONoff(simulateNSB), "i0 (Islands-cut)", nIslandsCut, ONoff(countIslands)); // log selections log(SIGNATURE, "%s:\n\t%20s: %s (%f:%f)\n", "Selections:", "Energy", ONoff(Select_Energy), Select_Energy_le, Select_Energy_ue); // Definition and initialization of array to save trigger statistics int ntriggerloop[(int) (Trigger_loop_uthres+1)][Trigger_loop_umult+1][Trigger_loop_utop+1]; for (ithrescount=Trigger_loop_lthres;ithrescount<=Trigger_loop_uthres;ithrescount++){ for (imulticount=Trigger_loop_lmult;imulticount<=Trigger_loop_umult;imulticount++){ for(itopocount=Trigger_loop_ltop;itopocount<=Trigger_loop_utop;itopocount++){ ntriggerloop[ithrescount][imulticount][itopocount]=0; } } } // set all random numbers seeds setall( get_seeds(0), get_seeds(1) ); // get list of showers to evt. skip nSkip = get_nskip_showers(); if (nSkip > 0) { Skip = new int[ nSkip ]; get_skip_showers( Skip ); log(SIGNATURE, "There are some showers to skip:\n"); for (i=0; iIsBatch()) { fprintf(stderr, "%s: cannot run in batch mode\n", argv[0]); // return 1; } if(FADC_Scan){ TApplication theAppFadc("App", &argc, argv); if (gROOT->IsBatch()) { fprintf(stderr, "%s: cannot run in batch mode\n", argv[0]); // return 1; } } // for safety and for dimensioning image_data: count the elements in the // diagnostic data branch i=0; i++; // "n" i++; // "primary" i++; // "energy" i++; // "cored" i++; // "impact" i++; // "xcore" i++; // "ycore" i++; // "theta" i++; // "phi" i++; // "deviations" i++; // "dtheta" i++; // "dphi" i++; // "trigger" i++; // "ncphs" i++; // "maxpassthr_phe" i++; // "nphes" i++; // "nphes2" i++; // "length" i++; // "width" i++; // "dist" i++; // "xdist" i++; // "azw" i++; // "miss" i++; // "alpha" i++; // "conc2" i++; // "conc3" i++; // "conc4" i++; // "conc5" i++; // "conc6" i++; // "conc7" i++; // "conc8" i++; // "conc9" i++; // "conc10" i++; // "asymx" i++; // "asymy" i++; // "phiasym" nvar = i; image_data = new float[nvar]; // set plate scale (deg/cm) and trigger area (deg) plateScale_cm2deg = ( ct_Type == 0 ) ? (0.244/2.1) : 0.030952381; if ( ! get_trigger_radius( °TriggerZone ) ) degTriggerZone = ( ct_Type == 0 ) ? (5.0) : (5.0); if ( ! get_correction( &fCorrection ) ) fCorrection = 1.0; // number of pixels for parameters anaPixels = get_ana_pixels(); anaPixels = (anaPixels == -1) ? ct_NPixels : anaPixels; // prepare the NSB simulation if( simulateNSB ){ // // Calculate the non-diffuse NSB photoelectron rates // k = produce_nsbrates( starfieldname, &cam, nsbrate_phepns ); if (k != 0){ cout << "Error when reading starfield... \nExiting.\n"; exit(1); } // calculate diffuse rate correcting for the pixel size for(i=0; i 0 ) // [ |l1 m1|2 |m1 n1|2 |n1 l1|2 ]1/2 // [ | | + | | + | | ] // [ |l2 m2| |m2 n2| |n2 l2| ] // // playing a little bit, we get this reduced for in our case: // // // dist = (- m2 n1 x + m1 n2 x + l2 n1 y - l1 n2 y - l2 m1 z + l1 m2 z) / // [(l2^2 (m1^2 + n1^2) + (m2 n1 - m1 n2)^2 - // 2 l1 l2 (m1 m2 + n1 n2) + l1^2 (m2^2 + n2^2) ] ^(1/2) // read the direction of the incoming shower thetashw = mcevth.get_theta(); phishw = mcevth.get_phi(); // calculate vector for shower l1 = sin(thetashw)*cos(phishw); m1 = sin(thetashw)*sin(phishw); n1 = cos(thetashw); // read the deviation of the telescope with respect to the shower mcevth.get_deviations ( &thetaCT, &phiCT ); if ( (thetaCT == 0.) && (phiCT == 0.) ) { // CT was looking to the source (both lines are parallel) // therefore, we calculate the impact parameter as the distance // between the CT axis and the core position impactD = dist_r_P( 0., 0., 0., l1, m1, n1, coreX, coreY, 0. ); } else { // the shower comes off-axis // obtain with this the final direction of the CT thetaCT += thetashw; phiCT += phishw; // calculate vector for telescope l2 = sin(thetaCT)*cos(phiCT); m2 = sin(thetaCT)*sin(phiCT); n2 = cos(thetaCT); num = (m1*n2*coreX - m2*n1*coreX + l2*n1*coreY - l1*n2*coreY); den = (SQR(l1*m2 - l2*m1) + SQR(m1*n2 - m2*n1) + SQR(n1*l2 - n2*l1)); den = sqrt(den); impactD = fabs(num)/den; // fprintf(stderr, "[%f %f,%f %f] (%f %f %f) (%f %f %f) %f/%f = ", // thetashw, phishw, thetaCT, phiCT, l1, m1, n1, l2, m2, n2, // num, den); } // energy cut if ( Select_Energy ) { if (( mcevth.get_energy() < Select_Energy_le ) || ( mcevth.get_energy() > Select_Energy_ue )) { log(SIGNATURE, "select_energy: shower rejected.\n"); continue; } } // Read first and last time and put inumphe to 0 mcevth.get_times(&arrtmin_ns,&arrtmax_ns); inumphe=0; // NSB simulation if(simulateNSB){ k = produce_nsb_phes( &arrtmin_ns, // will be changed by the function! &arrtmax_ns, // will be changed by the function! thetaCT, &cam, nsbrate_phepns, diffnsb_phepns, ext, fnpix, // will be changed by the function! &Trigger, // will be changed by the function! &fadc, // will be changed by the function! &inumphe, // important for later: the size of photoe[] baseline_mv // will be generated by the function ); if( k != 0 ){ // non-zero returnvalue means error cout << "Exiting.\n"; exit(1); } }// end if(simulateNSB) ... // read the photons and produce the photoelectrons k = produce_phes( inputfile, &cam, WAVEBANDBOUND1, WAVEBANDBOUND6, &Trigger, // will be changed by the function! &fadc, // will be changed by the function! &inumphe, // important for later: the size of photoe[] fnpix, // will be changed by the function! &ncph, // will be changed by the function! &arrtmin_ns, // will be changed by the function! &arrtmax_ns // will be changed by the function! ); if( k != 0 ){ // non-zero returnvalue means error cout << "Exiting.\n"; exit(1); } log(SIGNATURE, "End of this event: %d cphs(+%d). . .\n", ncph, ntcph); ntcph += ncph; // skip it ? for ( i=0; i0){ Lev1=Lev2=0; // loop over topologies for(itopocount=Trigger_loop_ltop;itopocount<=Trigger_loop_utop;itopocount++){ if(itopocount==0 && imulticount>7) continue; if(itopocount==2 && imulticount<3) continue; Trigger.SetTopology(itopocount); Trigger.ClearFirst(); // // Start the First Level Trigger simulation // McTrig->SetFirstLevel (Lev1=Trigger.FirstLevel()); if(Lev1>0) { ntriggerloop[ithrescount][imulticount][itopocount]++; McTrig->SetTopology(itopocount); McTrig->SetMultiplicity(imulticount); McTrig->SetThreshold(fpixelthres); } if(Lev1==0 && Write_All_Images){ McTrig->SetTopology(itopocount); McTrig->SetMultiplicity(imulticount); McTrig->SetThreshold(fpixelthres); Lev1=1; } for (Int_t ii=0;iiSetTime(Trigger.GetFirstLevelTime(ii)); McTrig->SetPixel(Trigger.GetFirstLevelPixel(ii)); fadc.TriggeredFadc(Trigger.GetFirstLevelTime(ii)); // // Fill the header of this event // Evt->FillHeader ( (UShort_t) (ntshow + nshow) , 20 ) ; // fill pixel information for(i=0;iFillPixel(i,fadcValues); } // // fill the MMcEvt with all information // McEvt->Fill( (UShort_t) mcevth.get_primary() , mcevth.get_energy(), mcevth.get_theta(), mcevth.get_phi(), mcevth.get_core(), mcevth.get_coreX(), mcevth.get_coreY(), impactD, ulli, ulli, (UShort_t) ncph, ulli, (UShort_t) ncph) ; // We don not count phtons out of the camera. // // write it out to the file outfile // EvtTree.Fill() ; // clear all Evt->Clear() ; McEvt->Clear() ; } McTrig->Clear() ; } } else break; } } // We study a single trigger condition else { Trigger.Diskriminate() ; // // look if in all the signals in the trigger signal branch // is a possible Trigger. Therefore we habe to diskriminate all // the simulated analog signals (Method Diskriminate in class // MTrigger). We look simultanously for the moments at which // there are more than TRIGGER_MULTI pixels above the // CHANNEL_THRESHOLD. // Lev0 = (Short_t) Trigger.ZeroLevel() ; Lev1 = Lev2 = 0 ; // // Start the First Level Trigger simulation // if ( Lev0 > 0 ) { McTrig->SetFirstLevel (Lev1 = Trigger.FirstLevel()); } if (Lev1>0){ ++ntrigger; } if (Lev1==0 && Write_All_Images){ Lev1=1; } McTrig->SetTopology(Trigger.GetTopology()); McTrig->SetMultiplicity(Trigger.GetMultiplicity()); for(i=0;iSetThreshold(fpixelthres); for(Int_t ii=0;iiSetTime(Trigger.GetFirstLevelTime(ii)); McTrig->SetPixel(Trigger.GetFirstLevelPixel(ii)); // // Fill the header of this event // Evt->FillHeader ( (UShort_t) (ntshow + nshow) , 20 ) ; // fill pixel information for(i=0;iFillPixel(i,fadcValues); } // // fill the MMcEvt with all information // McEvt->Fill( (UShort_t) mcevth.get_primary() , mcevth.get_energy(), mcevth.get_theta(), mcevth.get_phi(), mcevth.get_core(), mcevth.get_coreX(), mcevth.get_coreY(), impactD, ulli, ulli, (UShort_t) ncph, ulli, (UShort_t) ncph) ; // We don not count photons out of the camera. // // write it out to the file outfile // EvtTree.Fill() ; // // if a first level trigger occurred, then // 1. do some other stuff (not implemented) // 2. start the gui tool if(FADC_Scan){ if ( Lev0 > 0 ) { fadc.ShowSignal( McEvt, (Float_t) 60. ) ; } } if(Trigger_Scan){ if ( Lev0 > 0 ) { Trigger.ShowSignal(McEvt) ; } } // clear all Evt->Clear() ; McEvt->Clear() ; } McTrig->Clear() ; } #ifdef __DEBUG__ printf("\n"); for ( ici=0; ici -1 ) { if ( fnpix[(int)pixels[ici][icj][PIXNUM]] > 0. ) { printf ("@@ %4d %4d %10f %10f %4f (%4d %4d)\n", nshow, (int)pixels[ici][icj][PIXNUM], pixels[ici][icj][PIXX], pixels[ici][icj][PIXY], fnpix[(int)pixels[ici][icj][PIXNUM]], ici, icj); } } } } for (i=0; i (0:CT1 ¦ 1:MAGIC) // get focal distance sscanf(line, "%s %d", token, &ct_Type); log( "read_ct_file", ": %s\n", ((ct_Type==0) ? "CT1" : "MAGIC") ); break; case focal_distance: // [cm] // get focal distance sscanf(line, "%s %f", token, &ct_Focal_mean); log( "read_ct_file", ": %f cm\n", ct_Focal_mean ); break; case focal_std: // s(focal distance) [cm] // get focal distance sscanf(line, "%s %f", token, &ct_Focal_std); log( "read_ct_file", "s(Focal distance): %f cm\n", ct_Focal_std ); break; case point_spread: // [cm] // get point spread sscanf(line, "%s %f", token, &ct_PSpread_mean); log( "read_ct_file", ": %f cm\n", ct_PSpread_mean ); break; case point_std: // s(point spread) [cm] // get point spread sscanf(line, "%s %f", token, &ct_PSpread_std); log( "read_ct_file", "s(Point spread): %f cm\n", ct_PSpread_std ); break; case adjustment_dev: // s(adjustment_dev) [cm] // get point spread sscanf(line, "%s %f", token, &ct_Adjustment_std); log( "read_ct_file", "s(Adjustment): %f cm\n", ct_Adjustment_std ); break; case black_spot: // radius of the black spot in the center [cm] // get black spot radius sscanf(line, "%s %f", token, &ct_BlackSpot_rad); log( "read_ct_file", "Radius of the black spots: %f cm\n", ct_BlackSpot_rad); break; case r_mirror: // radius of the mirrors [cm] // get radius of mirror sscanf(line, "%s %f", token, &ct_RMirror); log( "read_ct_file", "Radii of the mirrors: %f cm\n", ct_RMirror ); break; case n_mirrors: // number of mirrors // get the name of the output_file from the line sscanf(line, "%s %d", token, &ct_NMirrors); log( "read_ct_file", "Number of mirrors: %d\n", ct_NMirrors ); break; case camera_width: // camera width [cm] // get the name of the ct_file from the line sscanf(line, "%s %f", token, &ct_CameraWidth); log( "read_ct_file", "Camera width: %f cm\n", ct_CameraWidth ); break; case n_pixels: // number of pixels // get the name of the output_file from the line sscanf(line, "%s %d", token, &ct_NPixels); log( "read_ct_file", "Number of pixels: %d\n", ct_NPixels ); break; case n_centralpixels: // number of central pixels // get the name of the output_file from the line sscanf(line, "%s %d", token, &ct_NCentralPixels); log( "read_ct_file", "Number of central pixels: %d\n", ct_NCentralPixels ); break; case n_gappixels: // number of gap pixels // get the name of the output_file from the line sscanf(line, "%s %d", token, &ct_NGapPixels); log( "read_ct_file", "Number of gap pixels: %d\n", ct_NGapPixels ); break; case pixel_width: // pixel width [cm] // get the name of the ct_file from the line sscanf(line, "%s %f", token, &ct_PixelWidth); ct_PixelWidth_corner_2_corner = ct_PixelWidth / cos(RAD(30.0)); ct_PixelWidth_corner_2_corner_half = ct_PixelWidth_corner_2_corner * 0.50; ct_Apot = ct_PixelWidth / 2; ct_2Apot = ct_Apot * 2.0; log( "read_ct_file", "Pixel width: %f cm\n", ct_PixelWidth ); break; case define_mirrors: // read table with the parameters of the mirrors log( "read_ct_file", "Table of mirrors data:\n" ); // check whether the number of mirrors was already set if ( ct_NMirrors == 0 ) error( "read_ct_file", "NMirrors was not set.\n" ); // allocate memory for paths list log( "read_ct_file", "Allocating memory for ct_data\n" ); ct_data = new float*[ct_NMirrors]; for (i=0; i> ct_data[i][j]; break; } // switch ( i ) } // end while // end log( "read_ct_file", "done.\n" ); return; } //!@} //!----------------------------------------------------------- // @name read_pixels // // @desc read pixels data // // @date Fri Mar 12 16:33:34 MET 1999 //------------------------------------------------------------ // @function //!@{ void read_pixels(struct camera *pcam) { ifstream qefile; char line[LINE_MAX_LENGTH]; int n, i, j, icount; float qe; //------------------------------------------------------------ // first, pixels' coordinates pcam->inumpixels = ct_NPixels; pcam->inumcentralpixels = ct_NCentralPixels; pcam->inumgappixels = ct_NGapPixels; pcam->inumbigpixels = ct_NPixels - ct_NCentralPixels - ct_NGapPixels; pcam->dpixdiameter_cm = ct_PixelWidth; // initialize pixel numbers pixary = new float* [2*ct_NCentralPixels]; for ( i=0; i<2*ct_NCentralPixels; ++i ) pixary[i] = new float[2]; pixneig = new int* [ct_NCentralPixels]; for ( i=0; i -1) && ((j-1) < pointsQE) && ((j-1) > -1) ) { QE[i-1][0][j-1] = QElambda[j-1]; QE[i-1][1][j-1] = qe; } if ( i > ct_NPixels) break; icount++; } if(icount/pointsQE < ct_NPixels){ error( "read_pixels", "The quantum efficiency file is faulty\n (found only %d pixels instead of %d).\n", icount/pointsQE, ct_NPixels ); } // close file qefile.close(); // test QE for(icount=0; icount< ct_NPixels; icount++){ for(i=0; i 1000. || QE[icount][1][i] < 0. || QE[icount][1][i] > 100.){ error( "read_pixels", "The quantum efficiency file is faulty\n pixel %d, point %d is % f, %f\n", icount, i, QE[icount][0][i], QE[icount][1][i] ); } } } // end log("read_pixels", "Done.\n"); } //!@} //!----------------------------------------------------------- // @name pixels_are_neig // // @desc check whether two pixels are neighbours // // @var pix1 Number of the first pixel // @var pix2 Number of the second pixel // @return TRUE: both pixels are neighbours; FALSE: oth. // // @date Wed Sep 9 17:58:37 MET DST 1998 //------------------------------------------------------------ // @function //!@{ int pixels_are_neig(int pix1, int pix2) { if ( sqrt(SQR( pixary[pix1][0] - pixary[pix2][0] ) + SQR( pixary[pix1][1] - pixary[pix2][1] ) ) > ct_PixelWidth_corner_2_corner ) return ( FALSE ); else return ( TRUE ); } //!@} //!----------------------------------------------------------- // @name igen_pixel_coordinates // // @desc generate the pixel center coordinates // // @var *pcam structure camera containing all the // camera information // @return total number of pixels // // DP // // @date Thu Oct 14 10:41:03 CEST 1999 //------------------------------------------------------------ // @function //!@{ /******** igen_pixel_coordinates() *********************************/ int igen_pixel_coordinates(struct camera *pcam) { /* generate pixel coordinates, return value is number of pixels */ int i, itot_inside_ring, iN, in, ipixno, iring_no, ipix_in_ring, isegment; float fsegment_fract; double dtsize; double dhsize; double dpsize; double dxfirst_pix; double dyfirst_pix; double ddxseg1, ddxseg2, ddxseg3, ddxseg4, ddxseg5, ddxseg6; double ddyseg1, ddyseg2, ddyseg3, ddyseg4, ddyseg5, ddyseg6; double dstartx, dstarty; /* for the gap pixels and outer pixels */ int j, nrow; dpsize = pcam->dpixdiameter_cm; dtsize = dpsize * sqrt(3.) / 2.; dhsize = dpsize / 2.; /* Loop over central pixels to generate co-ordinates */ for(ipixno=1; ipixno <= pcam->inumcentralpixels; ipixno++){ /* Initialise variables. The central pixel = ipixno 1 in ring iring_no 0 */ pcam->dpixsizefactor[ipixno-1] = 1.; in = 0; i = 0; itot_inside_ring = 0; iring_no = 0; /* Calculate the number of pixels out to and including the ring containing pixel number */ /* ipixno e.g. for pixel number 17 in ring number 2 itot_inside_ring = 19 */ while (itot_inside_ring == 0){ iN = 3*(i*(i+1)) + 1; if (ipixno <= iN){ iring_no = i; itot_inside_ring = iN; } i++; } /* Find the number of pixels which make up ring number iring_no e.g. ipix_in_ring = 6 for ring 1 */ ipix_in_ring = 0; for (i = 0; i < iring_no; ++i){ ipix_in_ring = ipix_in_ring + 6; } /* The camera is viewed as 6 radial segments ("pie slices"). Knowing the number of pixels in its */ /* ring calculate which segment the pixel ipixno is in. Then find how far across this segment it is */ /* as a fraction of the number of pixels in this sixth of the ring (ask SMB). */ isegment = 0; fsegment_fract = 0.; if (iring_no > 0) { isegment = (int)((ipixno - itot_inside_ring + ipix_in_ring - 0.5) / iring_no + 1); /* integer division ! numbering starts at 1 */ fsegment_fract = (ipixno - (itot_inside_ring - ipix_in_ring)) - ((isegment-1)*iring_no) - 1 ; } /* the first pixel in each ring lies on the positive x axis at a distance dxfirst_pix = iring_no * the */ /* pixel width (flat to flat) dpsize. */ dxfirst_pix = dpsize*iring_no; dyfirst_pix = 0.; /* the vector between the first and last pixels in a segment n is (ddxsegn, ddysegn) */ ddxseg1 = - dhsize*iring_no; ddyseg1 = dtsize*iring_no; ddxseg2 = -dpsize*iring_no; ddyseg2 = 0.; ddxseg3 = ddxseg1; ddyseg3 = -ddyseg1; ddxseg4 = -ddxseg1; ddyseg4 = -ddyseg1; ddxseg5 = -ddxseg2; ddyseg5 = 0.; ddxseg6 = -ddxseg1; ddyseg6 = ddyseg1; /* to find the position of pixel ipixno take the position of the first pixel in the ring and move */ /* anti-clockwise around the ring by adding the segment to segment vectors. */ switch (isegment) { case 0: pcam->dxc[ipixno-1] = 0.; pcam->dyc[ipixno-1] = 0.; case 1: pcam->dxc[ipixno-1] = dxfirst_pix - dhsize*fsegment_fract; pcam->dyc[ipixno-1] = dyfirst_pix + dtsize*fsegment_fract; break; case 2: pcam->dxc[ipixno-1] = dxfirst_pix + ddxseg1 - dpsize*fsegment_fract; pcam->dyc[ipixno-1] = dyfirst_pix + ddyseg1 + 0.; break; case 3: pcam->dxc[ipixno-1] = dxfirst_pix + ddxseg1 + ddxseg2 - dhsize*fsegment_fract; pcam->dyc[ipixno-1] = dyfirst_pix + ddyseg1 + ddyseg2 - dtsize*fsegment_fract; break; case 4: pcam->dxc[ipixno-1] = dxfirst_pix + ddxseg1 + ddxseg2 + ddxseg3 + dhsize*fsegment_fract; pcam->dyc[ipixno-1] = dyfirst_pix + ddyseg1 + ddyseg2 + ddyseg3 - dtsize*fsegment_fract; break; case 5: pcam->dxc[ipixno-1] = dxfirst_pix + ddxseg1 + ddxseg2 + ddxseg3 + ddxseg4 + dpsize*fsegment_fract; pcam->dyc[ipixno-1] = dyfirst_pix + ddyseg1 + ddyseg2 + ddyseg3 + ddyseg4 + 0.; break; case 6: pcam->dxc[ipixno-1] = dxfirst_pix + ddxseg1 + ddxseg2 + ddxseg3 + ddxseg4 + ddxseg5 + dhsize*fsegment_fract; pcam->dyc[ipixno-1] = dyfirst_pix + ddyseg1 + ddyseg2 + ddyseg3 + ddyseg4 + ddyseg5 + dtsize*fsegment_fract; break; default: fprintf(stderr, "ERROR: problem in coordinate generation for pixel %d\n", ipixno); return(0); } /* end switch */ } /* end for */ dstartx = pcam->dxc[pcam->inumcentralpixels - 1] + dhsize; dstarty = pcam->dyc[pcam->inumcentralpixels - 1] + dtsize; if(pcam->inumgappixels > 0){ /* generate the positions of the gap pixels */ j = pcam->inumcentralpixels; for(i=0; iinumgappixels; i=i+6){ pcam->dxc[j + i ] = dstartx + 2. * (i/6 + 1) * dpsize; pcam->dyc[j + i ] = dstarty; pcam->dpixsizefactor[j + i] = 1.; pcam->dxc[j + i + 1] = pcam->dxc[j + i ] / 2.; pcam->dyc[j + i + 1] = sqrt(3.) * pcam->dxc[j + i + 1]; pcam->dpixsizefactor[j + i + 1] = 1.; pcam->dxc[j + i + 2] = - pcam->dxc[j + i + 1]; pcam->dyc[j + i + 2] = pcam->dyc[j + i + 1]; pcam->dpixsizefactor[j + i+ 2] = 1.; pcam->dxc[j + i + 3] = - pcam->dxc[j + i]; pcam->dyc[j + i + 3] = dstarty; pcam->dpixsizefactor[j + i+ 3] = 1.; pcam->dxc[j + i + 4] = pcam->dxc[j + i + 2]; pcam->dyc[j + i + 4] = - pcam->dyc[j + i + 2]; pcam->dpixsizefactor[j + i+ 4] = 1.; pcam->dxc[j + i + 5] = pcam->dxc[j + i + 1]; pcam->dyc[j + i + 5] = - pcam->dyc[j + i + 1]; pcam->dpixsizefactor[j + i + 5] = 1.; } /* end for */ } /* end if */ /* generate positions of the outer pixels */ if( pcam->inumbigpixels > 0 ){ j = pcam->inumcentralpixels + pcam->inumgappixels; for(i=0; iinumbigpixels; i++){ pcam->dpixsizefactor[j + i] = 2.; } in = 0; nrow = (int) ceil(dstartx / 2. / dpsize); while(in < pcam->inumbigpixels){ pcam->dxc[j + in] = dstartx + dpsize; pcam->dyc[j + in] = dstarty + 2 * dpsize / sqrt(3.); pcam->dxc[j + in + nrow] = dstartx / 2. - dpsize / 2.; pcam->dyc[j + in + nrow] = sqrt(3.)/2. * dstartx + 2.5 * dpsize/sqrt(3.); pcam->dxc[j + in + 3 * nrow - 1] = - pcam->dxc[j + in]; pcam->dyc[j + in + 3 * nrow - 1] = pcam->dyc[j + in]; pcam->dxc[j + in + 3 * nrow] = - pcam->dxc[j + in]; pcam->dyc[j + in + 3 * nrow] = - pcam->dyc[j + in]; pcam->dxc[j + in + 5 * nrow - 1] = pcam->dxc[j + in + nrow]; pcam->dyc[j + in + 5 * nrow - 1] = - pcam->dyc[j + in + nrow]; pcam->dxc[j + in + 6 * nrow - 1] = pcam->dxc[j + in]; pcam->dyc[j + in + 6 * nrow - 1] = - pcam->dyc[j + in]; for(i=1; idxc[j + in + i] = pcam->dxc[j + in] - i * dpsize; pcam->dyc[j + in + i] = pcam->dyc[j + in] + i * dpsize * sqrt(3.); pcam->dxc[j + in + i + nrow] = pcam->dxc[j + in + nrow] - i * 2 * dpsize; pcam->dyc[j + in + i + nrow] = pcam->dyc[j + in + nrow]; pcam->dxc[j + in + 3 * nrow - 1 - i] = - pcam->dxc[j + in + i]; pcam->dyc[j + in + 3 * nrow - 1- i] = pcam->dyc[j + in + i]; pcam->dxc[j + in + i + 3 * nrow] = - pcam->dxc[j + in + i]; pcam->dyc[j + in + i + 3 * nrow] = - pcam->dyc[j + in + i]; pcam->dxc[j + in + 5 * nrow - 1 - i] = pcam->dxc[j + in + i + nrow]; pcam->dyc[j + in + 5 * nrow - 1 - i] = - pcam->dyc[j + in + i + nrow]; pcam->dxc[j + in + 6 * nrow - 1 - i] = pcam->dxc[j + in + i]; pcam->dyc[j + in + 6 * nrow - 1 - i] = - pcam->dyc[j + in + i]; } in = in + 6 * nrow; dstartx = dstartx + 2. * dpsize; nrow = nrow + 1; } /* end while */ } /* end if */ return(pcam->inumpixels); } //!@} //!----------------------------------------------------------- // @name bpoint_is_in_pix // // @desc check if a point (x,y) in camera coordinates is inside a given pixel // // @var *pcam structure camera containing all the // camera information // @var dx, dy point coordinates in centimeters // @var ipixnum pixel number (starting at 0) // @return TRUE if the point is inside the pixel, FALSE otherwise // // DP // // @date Thu Oct 14 16:59:04 CEST 1999 //------------------------------------------------------------ // @function //!@{ /******** bpoint_is_in_pix() ***************************************/ int bpoint_is_in_pix(double dx, double dy, int ipixnum, struct camera *pcam){ /* return TRUE if point (dx, dy) is in pixel number ipixnum, else return FALSE (use camera coordinate system) */ /* the pixel is assumed to be a "closed set" */ double a, b; /* a = length of one of the edges of one pixel, b = half the width of one pixel */ double c, xx, yy; /* auxiliary variable */ b = pcam->dpixdiameter_cm / 2. * pcam->dpixsizefactor[ipixnum]; a = pcam->dpixdiameter_cm / sqrt(3.) * pcam->dpixsizefactor[ipixnum]; c = 1. - 1./sqrt(3.); if((ipixnum < 0)||(ipixnum >= pcam->inumpixels)){ fprintf(stderr, "Error in bpoint_is_in_pix: invalid pixel number %d\n", ipixnum); fprintf(stderr, "Exiting.\n"); exit(203); } xx = dx - pcam->dxc[ipixnum]; yy = dy - pcam->dyc[ipixnum]; if(((-b <= xx) && (xx <= 0.) && ((-c * xx - a) <= yy) && (yy <= ( c * xx + a))) || ((0. < xx) && (xx <= b ) && (( c * xx - a) <= yy) && (yy <= (-c * xx + a))) ){ return(TRUE); /* point is inside */ } else{ return(FALSE); /* point is outside */ } } //!@} //------------------------------------------------------------ // @name dist_r_P // // @desc distance straight line r - point P // // @date Sat Jun 27 05:58:56 MET DST 1998 // @function @code //------------------------------------------------------------ // dist_r_P // // distance straight line r - point P //------------------------------------------------------------ float dist_r_P(float a, float b, float c, float l, float m, float n, float x, float y, float z) { return ( sqrt((SQR((a-x)*m-(b-y)*l) + SQR((b-y)*n-(c-z)*m) + SQR((c-z)*l-(a-x)*n))/ (SQR(l)+SQR(m)+SQR(n)) ) ); } //------------------------------------------------------------ // @name check_reflector_file // // @desc check if a given reflector file has the right signature // @desc return TRUE or FALSE // // @date Mon Feb 14 16:44:21 CET 2000 // @function @code //------------------------------------------------------------ int check_reflector_file(FILE *infile){ char Signature[20]; // auxiliary variable char sign[20]; // auxiliary variable strcpy(Signature, REFL_SIGNATURE); strcpy(sign, Signature); fread( (char *)sign, strlen(Signature), 1, infile); if (strcmp(sign, Signature) != 0) { cout << "ERROR: Signature of .rfl file is not correct\n"; cout << '"' << sign << '"' << '\n'; cout << "should be: " << Signature << '\n'; return(FALSE); } fread( (char *)sign, 1, 1, infile); return(TRUE); } //------------------------------------------------------------ // @name lin_interpol // // @desc interpolate linearly between two points returning the // @desc y-value of the result // // @date Thu Feb 17 11:31:32 CET 2000 // @function @code //------------------------------------------------------------ float lin_interpol( float x1, float y1, float x2, float y2, float x){ if( (x2 - x1)<=0. ){ // avoid division by zero, return average cout << "Warning: lin_interpol was asked to interpolate between two points with x1>=x2.\n"; return((y1+y2)/2.); } else{ // note: the check whether x1 < x < x2 is omitted for speed reasons return((float) (y1 + (y2-y1)/(x2-x1)*(x-x1)) ); } } //------------------------------------------------------------ // @name produce_phes // // @desc read the photons of an event, pixelize them and simulate QE // @desc return various statistics and the array of Photoelectrons // // @date Mon Feb 14 16:44:21 CET 2000 // @function @code //------------------------------------------------------------ int produce_phes( FILE *sp, // the input file struct camera *cam, // the camera layout float minwl_nm, // the minimum accepted wavelength float maxwl_nm, // the maximum accepted wavelength class MTrigger *trigger, // the generated phes class MFadc *fadc, int *itotnphe, // total number of produced photoelectrons float nphe[iMAXNUMPIX], // number of photoelectrons in each pixel int *incph, // total number of cph read float *tmin_ns, // minimum arrival time of all phes float *tmax_ns // maximum arrival time of all phes ){ static int i, k, ipixnum; static float cx, cy, wl, qe, t; static MCCphoton photon; static float **qept; static char flag[SIZE_OF_FLAGS + 1]; static float radius; // reset variables for ( i=0; iinumpixels; ++i ){ nphe[i] = 0.0; } *incph = 0; radius = cam->dxc[cam->inumpixels-1] + 1.5*cam->dpixdiameter_cm*cam->dpixsizefactor[cam->inumpixels-1]; //- - - - - - - - - - - - - - - - - - - - - - - - - // read photons and "map" them into the pixels //-------------------------------------------------- // initialize CPhoton photon.fill(0., 0., 0., 0., 0., 0., 0., 0.); // read the photons data fread ( flag, SIZE_OF_FLAGS, 1, sp ); // loop over the photons while ( !isA( flag, FLAG_END_OF_EVENT ) ) { memcpy( (char*)&photon, flag, SIZE_OF_FLAGS ); fread( ((char*)&photon)+SIZE_OF_FLAGS, photon.mysize()-SIZE_OF_FLAGS, 1, sp ); // increase number of photons (*incph)++; // Chceck if photon is inside trigger time range t = photon.get_t() ; if (t-*tmin_ns>TOTAL_TRIGGER_TIME) { // read next Photon fread( flag, SIZE_OF_FLAGS, 1, sp ); // go to beginning of loop, the photon is lost continue; } // // Pixelization // cx = photon.get_x(); cy = photon.get_y(); // get wavelength wl = photon.get_wl(); // cout << "wl " << wl << " radius " << sqrt(cx*cx + cy*cy) << "\n"; // check if photon has valid wavelength and is inside outermost camera radius if( (wl > maxwl_nm) || (wl < minwl_nm) || (sqrt(cx*cx + cy*cy) > radius ) ){ // cout << " lost first check\n"; // read next CPhoton fread ( flag, SIZE_OF_FLAGS, 1, sp ); // go to beginning of loop, the photon is lost continue; } ipixnum = -1; for(i=0; iinumpixels; i++){ if( bpoint_is_in_pix( cx, cy, i, cam) ){ ipixnum = i; break; } } if(ipixnum==-1){// the photon is in none of the pixels // cout << " lost pixlization\n"; // read next CPhoton fread ( flag, SIZE_OF_FLAGS, 1, sp ); // go to beginning of loop, the photon is lost continue; } //+++ // QE simulation //--- // set pointer to the QE table of the relevant pixel qept = (float **)QE[ipixnum]; // check if wl is inside table; outside the table, QE is assumed to be zero if((wl < qept[0][0]) || (wl > qept[0][pointsQE-1])){ // cout << " lost\n"; // read next Photon fread ( flag, SIZE_OF_FLAGS, 1, sp ); // go to beginning of loop continue; } // find data point in the QE table (-> k) k = 1; // start at 1 because the condition was already tested for 0 while (k < pointsQE-1 && qept[0][k] < wl){ k++; } // calculate the qe value between 0. and 1. qe = lin_interpol(qept[0][k-1], qept[1][k-1], qept[0][k], qept[1][k], wl) / 100.0; // if random > quantum efficiency, reject it if ( RandomNumber > qe ) { // cout << " lost\n"; // read next Photon fread ( flag, SIZE_OF_FLAGS, 1, sp ); // go to beginning of loop continue; } //+++ // The photon has produced a photo electron //--- // cout << " accepted\n"; // increment the number of photoelectrons in the relevant pixel nphe[ipixnum] += 1.0; // store the new photoelectron fadc->Fill(ipixnum,(t-*tmin_ns) , trigger->FillShow(ipixnum,t-*tmin_ns)); *itotnphe += 1; // read next Photon fread( flag, SIZE_OF_FLAGS, 1, sp ); } // end while, i.e. found end of event return(0); } //------------------------------------------------------------ // @name produce_nsbrates // // @desc read the starfield file, call produce_phes on it in, // @desc different wavebands, calculate the nsbrates // // @date Mon Feb 14 16:44:21 CET 2000 // @function @code //------------------------------------------------------------ int produce_nsbrates( char *iname, // the starfield input file name struct camera *cam, // camera layout float rate_phepns[iMAXNUMPIX][iNUMWAVEBANDS] // the product of this function: // the NSB rates in phe/ns for each pixel ){ int i, j, k, ii; // counters MTrigger trigger; MFadc flashadc; static float wl_nm[iNUMWAVEBANDS + 1] = { WAVEBANDBOUND1, WAVEBANDBOUND2, WAVEBANDBOUND3, WAVEBANDBOUND4, WAVEBANDBOUND5, WAVEBANDBOUND6 }; FILE *infile; // the input file fpos_t fileposition; // marker on the input file static char cflag[SIZE_OF_FLAGS + 1]; // auxiliary variable static MCEventHeader evth; // the event header static float nphe[iMAXNUMPIX]; // the number of photoelectrons in each pixel int itnphe; // total number of produced photoelectrons int itotnphe; // total number of produced photoelectrons after averaging int incph; // total number of cph read float tmin_ns; // minimum arrival time of all phes float tmax_ns; // maximum arrival time of all phes float integtime_ns; // integration time from the starfield generator // open input file log(SIGNATURE, "Opening starfield input \"rfl\" file %s\n", iname ); infile = fopen( iname, "r" ); if ( infile == NULL ) error( SIGNATURE, "Cannot open starfield input file: %s\n", iname ); // get signature, and check it if(check_reflector_file(infile)==FALSE){ exit(1); } // initialize flag strcpy( cflag, " \0" ); // get flag fread( cflag, SIZE_OF_FLAGS, 1, infile ); if ( ! feof(infile)){ // reading .rfl file if(!isA( cflag, FLAG_START_OF_RUN )){ error( SIGNATURE, "Expected start of run flag, but found: %s\n", cflag ); } else { // found start of run fread( cflag, SIZE_OF_FLAGS, 1, infile ); if( isA( cflag, FLAG_START_OF_EVENT )){ // there is a event // get MCEventHeader fread( (char*)&evth, evth.mysize(), 1, infile ); integtime_ns = evth.get_energy(); // memorize where we are in the file if (fgetpos( infile, &fileposition ) != 0){ error( SIGNATURE, "Cannot position in file ...\n"); } // loop over the wavebands for(i=0; iinumpixels; j++){ // loop over pixels rate_phepns[j][i] = 0.; } itotnphe = 0; // read the photons and produce the photoelectrons // - in order to average over the QE simulation, call the // production function iNUMNSBPRODCALLS times for(ii=0; iiinumpixels; j++){ // loop over pixels rate_phepns[j][i] += nphe[j]/integtime_ns/(float)iNUMNSBPRODCALLS; } itotnphe += itnphe; } // end for(ii=0 ... fprintf(stdout, "Starfield, %6f - %6f nm: %d photoelectrons for %6f ns integration time\n", wl_nm[i], wl_nm[i+1], itotnphe/iNUMNSBPRODCALLS, integtime_ns); } // end for(i=0 ... } else{ cout << "Starfield file contains no event.\nExiting.\n"; exit(1); } // end if( isA ... event } // end if ( isA ... run } else{ cout << "Starfield file contains no run.\nExiting.\n"; exit(1); } fclose( infile ); return(0); } //------------------------------------------------------------ // @name produce_nsb_phes // // @desc produce the photoelectrons from the nsbrates // // @date Thu Feb 17 17:10:40 CET 2000 // @function @code //------------------------------------------------------------ int produce_nsb_phes( float *atmin_ns, float *atmax_ns, float theta_rad, struct camera *cam, float nsbr_phepns[iMAXNUMPIX][iNUMWAVEBANDS], float difnsbr_phepns[iMAXNUMPIX], float extinction[iNUMWAVEBANDS], float fnpx[iMAXNUMPIX], MTrigger *trigger, MFadc *fadc, int *inphe, float base_mv[iMAXNUMPIX]){ float simtime_ns; // NSB simulation time int i, j, k, ii; float zenfactor; // correction factor calculated from the extinction int inumnsbphe; // number of photoelectrons caused by NSB float t; ii = *inphe; // avoid dereferencing // check if the arrival times are set; if not generate them if(*atmin_ns *atmax_ns){ *atmin_ns = 0.; *atmax_ns = simtime_ns = TOTAL_TRIGGER_TIME; } else{ // extend the event time window by the given offsets *atmin_ns = *atmin_ns - SIMTIMEOFFSET_NS; *atmax_ns = *atmax_ns + SIMTIMEOFFSET_NS; simtime_ns = *atmax_ns - *atmin_ns; // make sure the simulated time is long enough for the FADC // simulation and not too long if(simtime_ns< TOTAL_TRIGGER_TIME){ *atmin_ns = *atmin_ns -(TOTAL_TRIGGER_TIME-simtime_ns)/2; *atmax_ns = *atmin_ns + TOTAL_TRIGGER_TIME; simtime_ns = TOTAL_TRIGGER_TIME; } if(simtime_ns> TOTAL_TRIGGER_TIME){ *atmax_ns =*atmin_ns + TOTAL_TRIGGER_TIME; simtime_ns = TOTAL_TRIGGER_TIME; } } // initialize baselines for(i=0; iinumpixels; i++){ base_mv[i] = 0.; } // calculate baselines and generate phes for(i=0; iinumpixels; j++){ // loop over the pixels inumnsbphe = (int) ((zenfactor * nsbr_phepns[j][i] + difnsbr_phepns[j]/iNUMWAVEBANDS) * simtime_ns ); base_mv[j] += inumnsbphe; // randomize inumnsbphe = ignpoi( inumnsbphe ); // create the photoelectrons for(k=0; k < inumnsbphe; k++){ t=(RandomNumber * simtime_ns); (*fadc).Fill(j,t ,(*trigger).FillNSB(j,t)); ii++; // increment total number of photoelectons fnpx[j] += 1.; // increment number of photoelectrons in this pixel } } // end for(j=0 ... } // end for(i=0 ... // finish baseline calculation for(i=0; iinumpixels; i++){ base_mv[i] *= RESPONSE_FWHM * RESPONSE_AMPLITUDE / simtime_ns; } *inphe = ii; // update the pointer return(0); } // @endcode //=------------------------------------------------------------ //!@subsection Log of this file. //!@{ // // $Log: not supported by cvs2svn $ // Revision 1.7 2000/03/24 18:10:46 blanch // A first FADC simulation and a trigger simulation are already implemented. // The calculation of the Hillas Parameters have been removed, since it was decided that it should be in the analysis software. // 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. // // Revision 1.6 2000/03/20 18:35:11 blanch // 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. // // Revision 1.5 2000/02/18 17:40:35 petry // This version includes drastic changes compared to camera.cxx 1.4. // It is not yet finished and not immediately useful because the // trigger simulation is not yet re-implemented. I had to take it // out together with some other stuff in order to tidy the whole // program up. This is not meant as an insult to anyone. I needed // to do this in order to be able to work on it. // // This version has been put in the repository in order to be // able to share the further development with others. // // If you need something working, wait or take an earlier one. // See file README. // // Revision 1.4 2000/01/25 08:36:23 petry // The pixelization in previous versions was buggy. // This is the first version with a correct pixelization. // // Revision 1.3 2000/01/20 18:22:17 petry // Found little bug which makes camera crash if it finds a photon // of invalid wavelength. This bug is now fixed and the range // of valid wavelengths extended to 290 - 800 nm. // This is in preparation for the NSB simulation to come. // Dirk // // Revision 1.2 1999/11/19 08:40:42 harald // Now it is possible to compile the camera programm under osf1. // // Revision 1.1.1.1 1999/11/05 11:59:31 harald // This the starting point for CVS controlled further developments of the // camera program. The program was originally written by Jose Carlos. // But here you can find a "rootified" version to the program. This means // that there is no hbook stuff in it now. Also the output of the // program changed to the MagicRawDataFormat. // // The "rootification" was done by Dirk Petry and Harald Kornmayer. // // Revision 1.3 1999/10/22 15:01:28 petry // version sent to H.K. and N.M. on Fri Oct 22 1999 // // Revision 1.2 1999/10/22 09:44:23 petry // first synthesized version which compiles and runs without crashing; // // Revision 1.1.1.1 1999/10/21 16:35:10 petry // first synthesised version // // Revision 1.13 1999/03/15 14:59:05 gonzalez // camera-1_1 // // Revision 1.12 1999/03/02 09:56:10 gonzalez // *** empty log message *** // // //!@} //=EOF