source: trunk/MagicSoft/Simulation/Detector/include-MFadc/MFadc.cxx

Last change on this file was 6714, checked in by moralejo, 20 years ago
*** empty log message ***
File size: 30.6 KB
Line 
1////////////////////////////////////////////////////////////////
2//
3// MFadc
4//
5//
6#include "MFadc.hxx"
7
8#include "MMcEvt.hxx"
9
10#include "TROOT.h"
11#include <TApplication.h>
12#include <TVirtualX.h>
13#include <TGClient.h>
14
15#include "TH1.h"
16#include "TObjArray.h"
17
18#include "MGFadcSignal.hxx"
19
20using namespace std;
21
22MFadc::MFadc(Int_t pix, Int_t shape, Float_t integral, Float_t fwhm,
23 Int_t shapeout, Float_t integralout, Float_t fwhmout,
24 Float_t trigger_delay, Float_t fadc_slices_per_ns,
25 Int_t fadc_slices_written, Int_t gainswitchamp,
26 Int_t shiftfromswitch2lowgain, Float_t hi2logainpeak) {
27 //
28 // Constructor overloaded II
29 //
30 // Input variables:
31 // 1. integral(out) = integration of the single phe response for inner
32 // (outer) pixels.
33 // 2. fwhm(out) = width at half high of the single phe response for
34 // inner (outer) pixels.
35 //
36 // trigger_delay: shift of signals towards later times in FADC, in order
37 // to center the signals in a good range. It acts as a sort of delay of
38 // the signals (before being sent to the FADC) with respect to the trigger.
39 //
40 // The procedure is the following:
41 // 1. some parameters of the trigger are set to default.
42 // this parameters of the trigger may be changed
43 // 3. Then the all signals are set to zero
44
45 numpix=pix;
46
47 fwhm_resp = fwhm;
48 integ_resp = integral;
49 fwhm_resp_outer = fwhmout;
50 integ_resp_outer = integralout;
51 shape_resp = shape;
52 shape_resp_outer = shapeout;
53 fFadcSlicesPerNanosec = fadc_slices_per_ns;
54 fFadcSlices = fadc_slices_written;
55 fGainSwitchAmp = gainswitchamp;
56 fShiftFromSwitch2LowGain = shiftfromswitch2lowgain;
57 fHi2LoGainPeak = hi2logainpeak;
58
59 fSlices_mFadc = (Int_t)(TOTAL_TRIGGER_TIME*fFadcSlicesPerNanosec);
60
61 for (Int_t i = 0; i < CAMERA_PIXELS; i++)
62 {
63 sig[i] = new Float_t[fSlices_mFadc];
64 sig_LG[i] = new Float_t[fSlices_mFadc];
65 }
66
67 noise = new Float_t[fSlices_mFadc*1001];
68 noise_outer = new Float_t[fSlices_mFadc*1001];
69 digital_noise = new Float_t[fSlices_mFadc*1001];
70
71 for (Int_t i = 0; i < CAMERA_PIXELS; i++)
72 {
73 output[i] = new Float_t[fFadcSlices];
74 output_lowgain[i] = new Float_t[fFadcSlices];
75 }
76
77 cout<< "[MFadc] Setting up the MFadc with this values "<< endl ;
78 cout<< "[MFadc] FADC sampling frequency: " << fFadcSlicesPerNanosec << " GHz" << endl ;
79 cout<< "[MFadc] - Inner pixels : "<< endl ;
80
81 switch(shape_resp){
82 case 0:
83 cout<< "[MFadc] Pulse shape : Gaussian ("<<shape_resp<<")"<< endl ;
84 cout<< "[MFadc] Response Area : "<<integ_resp<<" adc counts"<< endl ;
85 cout<< "[MFadc] Response FWHM : "<<fwhm_resp<<" ns"<< endl ;
86 break;
87 case 1:
88 cout<< "[MFadc] Pulse shape : From Pulpo ("<<shape_resp<<")"<< endl ;
89 cout<< "[MFadc] Response Area : "<<integ_resp<<" adc counts"<< endl ;
90 break;
91 default:
92 cout<< "[MFadc] Pulse shape unknown"<<endl;
93 }
94 cout<< "[MFadc] - Outer pixels : "<< endl ;
95 switch(shape_resp_outer){
96 case 0:
97 cout<< "[MFadc] Pulse shape : Gaussian ("<<shape_resp_outer<<")"<<endl;
98 cout<< "[MFadc] Response Area : "<<integ_resp_outer<<" adc counts"<<endl;
99 cout<< "[MFadc] Response FWHM : "<<fwhm_resp_outer<<" ns"<< endl ;
100 break;
101 case 1:
102 cout<< "[MFadc] Pulse shape : From Pulpo ("<<shape_resp_outer<<")"<<endl;
103 cout<< "[MFadc] Response Area : "<<integ_resp_outer<<" adc counts"<< endl ;
104 break;
105 default:
106 cout<< "[MFadc] Pulse shape unknown ("<<shape_resp_outer<<")"<<endl;
107 }
108
109
110 //
111 // set up the response shape
112 //
113
114 //
115 // First select number of bins for the histogram which will contain the single
116 // photoelectron response of the FADC. The width of these bins is smaller than that
117 // of the real FADC slices by a factor SUBBINS (see MFadcDefine.h):
118 //
119 if (shape_resp == 1)
120 fResponseSlicesFadc = (Int_t)(50.*fFadcSlicesPerNanosec*SUBBINS);
121 // 50 ns range
122
123 else
124 fResponseSlicesFadc = (Int_t)(7*fwhm_resp/2.35*fFadcSlicesPerNanosec*SUBBINS);
125 // 7 sigma range
126
127 sing_resp = new Float_t[fResponseSlicesFadc];
128 sing_resp_outer = new Float_t[fResponseSlicesFadc];
129
130 sing_resp_lowgain = new Float_t[fResponseSlicesFadc];
131 sing_resp_outer_lowgain = new Float_t[fResponseSlicesFadc];
132
133
134 // Parameters for real pulse shaped as measured with the Pulpo device:
135 fPulseParameters[0] = 2.066;
136 fPulseParameters[1] = 1.568;
137 fPulseParameters[2] = 3; // This will set the peak of the pulse at x ~ 3*3.3 = 10 ns
138 // It is just a safe value so that the pulse is well contained.
139 fPulseParameters[3] = 0.00282;
140 fPulseParameters[4] = 0.04093;
141 fPulseParameters[5] = 0.2411;
142 fPulseParameters[6] = -0.009442;
143
144 // Now for the low gain:
145 fPulseParametersLG[0] = 6.88031e-01;
146 fPulseParametersLG[1] = 1.19492e+00;
147 fPulseParametersLG[2] = 3.; // Same as above
148 fPulseParametersLG[3] = 8.10590e-03;
149 fPulseParametersLG[4] = 1.92851e-01;
150 fPulseParametersLG[5] = 3.30315e-01;
151 fPulseParametersLG[6] = 1.12980e-03;
152
153 Int_t i ;
154
155 Float_t sigma ;
156 Float_t x, x0 ;
157 Float_t dX, dX2 ;
158
159 Float_t response_sum_inner, response_sum_outer;
160 Float_t response_sum_inner_LG, response_sum_outer_LG;
161
162 response_sum_inner = 0.;
163 response_sum_outer = 0.;
164 response_sum_inner_LG = 0.;
165 response_sum_outer_LG = 0.;
166
167 dX = 1. / fFadcSlicesPerNanosec / SUBBINS ; // Units: ns
168 dX2 = dX/2. ;
169
170 switch(shape_resp){
171
172 case 0:
173 sigma = fwhm_resp / 2.35 ;
174 x0 = 3*sigma;
175 fadc_time_offset = trigger_delay-x0; // ns
176
177 for (i = 0; i < fResponseSlicesFadc ; i++ )
178 {
179 x = i * dX + dX2 ;
180
181 sing_resp[i] = (Float_t)(expf(-0.5*(x-x0)*(x-x0)/(sigma*sigma)));
182 sing_resp_lowgain[i] = sing_resp[i];
183
184 // Use the same pulse for high and low gain in case of gaussian shape
185 // FIXME: it would be nice to be able to choose a different shape for
186 // high and low gain:
187 response_sum_inner += sing_resp[i];
188 response_sum_inner_LG += sing_resp_lowgain[i];
189 }
190
191 break;
192 case 1:
193 float p0,p1,p2,p3,p4,p5,p6;
194 float d;
195 float zed_slices;
196 // Parameters values extracted from fitting a real FADC response
197 // gaussian electronic pulse passed through the whole chain from
198 // transmitter boards to FADC.
199 p0 = fPulseParameters[0];
200 p1 = fPulseParameters[1];
201 p2 = fPulseParameters[2];
202 p3 = fPulseParameters[3];
203 p4 = fPulseParameters[4];
204 p5 = fPulseParameters[5];
205 p6 = fPulseParameters[6];
206
207 float p0_LG, p1_LG, p2_LG, p3_LG, p4_LG, p5_LG, p6_LG;
208 p0_LG = fPulseParametersLG[0];
209 p1_LG = fPulseParametersLG[1];
210 p2_LG = fPulseParametersLG[2];
211 p3_LG = fPulseParametersLG[3];
212 p4_LG = fPulseParametersLG[4];
213 p5_LG = fPulseParametersLG[5];
214 p6_LG = fPulseParametersLG[6];
215
216 p2_LG += (fHi2LoGainPeak - fFadcSlices); // Adjust distance between high and low gain peaks
217
218
219 // Now define the time before trigger to read FADC signal when it
220 // has to be written. Here FADC_SLICES_PER_NSEC (=0.3) is the value
221 // for the 300 MHz MAGIC FADCs and must NOT be changed, even if you
222 // use a faster sampling in the simulation (through the input card
223 // command "fadc_GHz"), because this is just a conversion of units. The
224 // parameters of the "pulpo" pulse shape were obtained with the 300 MHz
225 // FADC and so we convert the time parameter to units of 3.3 ns slices
226 // just to use the provided parametrization, and no matter what sampling
227 // frequency we are simulating!
228
229 fadc_time_offset = trigger_delay - p2 / FADC_SLICES_PER_NSEC; // ns
230
231
232 for (i=0; i< fResponseSlicesFadc ; i++ )
233 {
234 x = i * dX + dX2;
235
236 // x has to be converted from ns to units FADC slices of the default
237 // FADC of 300 MHz (these are just units, and must be these even if you
238 // are using another sampling frequency!):
239 //
240 zed_slices = x * FADC_SLICES_PER_NSEC - p2;
241 d = (zed_slices>0)? 0.5 : -0.5;
242
243 sing_resp[i] = (Float_t) (p0*exp(-p1*(exp(-p1*zed_slices)+zed_slices))+
244 p3+p4*exp(-p1*(exp(-p1*zed_slices)+
245 p5*zed_slices))+p6*d);
246 response_sum_inner += sing_resp[i];
247
248
249 // Now the low gain:
250
251 zed_slices = x * FADC_SLICES_PER_NSEC - p2_LG;
252 d = (zed_slices>0)? 0.5 : -0.5;
253
254 sing_resp_lowgain[i] = (Float_t) (p0_LG*exp(-p1_LG*(exp(-p1_LG*zed_slices)+zed_slices))+
255 p3_LG+p4_LG*exp(-p1_LG*(exp(-p1_LG*zed_slices)+
256 p5_LG*zed_slices))+p6_LG*d);
257 response_sum_inner_LG += sing_resp_lowgain[i];
258
259 }
260
261 break;
262 default:
263 cout<<"[MFadc] MFadc::MFadc : Shape of FADC pulse for inner pixel unknown."
264 <<endl;
265 cout<<"[MFadc] MFadc::MFadc : Exiting Camera ..."
266 <<endl;
267 exit(1);
268 }
269
270 // Response for outer pixels
271
272 switch(shape_resp_outer){
273
274 case 0:
275 sigma = fwhm_resp_outer / 2.35 ;
276 x0 = 3*sigma ;
277
278 for (i = 0; i < fResponseSlicesFadc ; i++ )
279 {
280 x = i * dX + dX2 ;
281
282 //
283 // the value 1/sqrt(2*Pi*sigma^2) was introduced to normalize
284 // the area at the input value After this, the integral
285 // of the response will be integ_resp.
286 //
287 sing_resp_outer[i] = (Float_t) (expf(-0.5 * (x-x0)*(x-x0) /
288 (sigma*sigma)) ) ;
289 response_sum_outer += sing_resp_outer[i];
290 }
291 break;
292 case 1:
293 float p0,p1,p2,p3,p4,p5,p6;
294 float d;
295 float zed_slices;
296 // Parameters values extracted from fitting a real FADC response
297 // gaussian electronic pulse passed through the whole chain from
298 // transmitter boards to FADC.
299
300 p0 = fPulseParameters[0];
301 p1 = fPulseParameters[1];
302 p2 = fPulseParameters[2];
303 p3 = fPulseParameters[3];
304 p4 = fPulseParameters[4];
305 p5 = fPulseParameters[5];
306 p6 = fPulseParameters[6];
307
308 float p0_LG, p1_LG, p2_LG, p3_LG, p4_LG, p5_LG, p6_LG;
309 p0_LG = fPulseParametersLG[0];
310 p1_LG = fPulseParametersLG[1];
311 p2_LG = fPulseParametersLG[2];
312 p3_LG = fPulseParametersLG[3];
313 p4_LG = fPulseParametersLG[4];
314 p5_LG = fPulseParametersLG[5];
315 p6_LG = fPulseParametersLG[6];
316
317 p2_LG += (fHi2LoGainPeak - fFadcSlices); // Adjust distance between high and low gain peaks
318
319 for (i=0; i< fResponseSlicesFadc ; i++ )
320 {
321 x = i * dX + dX2;
322
323 // x has to be converted from ns to units FADC slices of the default
324 // FADC of 300 MHz (these are just units, and must be these even if you
325 // are using another sampling frequency!):
326 //
327 zed_slices = x * FADC_SLICES_PER_NSEC - p2;
328 d = (zed_slices>0)? 0.5 : -0.5;
329
330 sing_resp_outer[i] = (Float_t) (p0*exp(-p1*(exp(-p1*zed_slices)+
331 zed_slices))+p3+
332 p4*exp(-p1*(exp(-p1*zed_slices)+
333 p5*zed_slices))+p6*d);
334 response_sum_outer += sing_resp_outer[i];
335
336
337 // Now the low gain:
338
339 zed_slices = x * FADC_SLICES_PER_NSEC - p2_LG;
340 d = (zed_slices>0)? 0.5 : -0.5;
341
342 sing_resp_outer_lowgain[i] = (Float_t) (p0_LG*exp(-p1_LG*(exp(-p1_LG*zed_slices)+zed_slices))+
343 p3_LG+p4_LG*exp(-p1_LG*(exp(-p1_LG*zed_slices)+
344 p5_LG*zed_slices))+p6_LG*d);
345 response_sum_outer_LG += sing_resp_outer_lowgain[i];
346 }
347 break;
348 default:
349 cout<<"[MFadc] MFadc::MFadc : Shape of FADC pulse for inner pixel unknown."
350 <<endl;
351 cout<<"[MFadc] MFadc::MFadc : Exiting Camera ..."
352 <<endl;
353 exit(1);
354 }
355
356 //
357 // Normalize responses to values set through input card: (= set gain of electronic chain)
358 // Take into account that only 1 of every SUBBINS bins of sing_resp[] will be "sampled" by
359 // the FADC, so we have to correct for this to get the right "FADC integral" (=integ_resp)
360 // per photoelectron:
361 //
362
363 for (i = 0; i < fResponseSlicesFadc ; i++ )
364 {
365 sing_resp[i] *= integ_resp / response_sum_inner * SUBBINS;
366 sing_resp_outer[i] *= integ_resp_outer / response_sum_outer * SUBBINS;
367
368 // The low gain will be further scaled down later; for now we normalize it with to same
369 // integral as the high gain.
370 sing_resp_lowgain[i] *= integ_resp / response_sum_inner_LG * SUBBINS;
371 sing_resp_outer_lowgain[i] *= integ_resp_outer / response_sum_outer_LG * SUBBINS;
372 }
373
374 //
375 // init the Random Generator for Electronic Noise
376 //
377
378 GenElec = new TRandom () ;
379
380 Reset();
381
382 //
383 // set all pedestals to 0
384 //
385
386 for ( i =0 ; i <CAMERA_PIXELS ; i++ ) {
387 pedestal[i] = 0.0 ;
388 }
389
390 cout<<" end of MFadc::MFadc()"<<endl;
391}
392
393void MFadc::Reset() {
394 //
395 // set all values of the signals to zero
396 // set the values of FADC slices that would be read after trigger to zero
397 //
398 memset(used, 0, CAMERA_PIXELS*sizeof(Bool_t));
399
400 for (Int_t i = 0; i < CAMERA_PIXELS; i++)
401 {
402 memset(output[i], 0, fFadcSlices*sizeof(Float_t));
403 memset(output_lowgain[i], 0, fFadcSlices*sizeof(Float_t));
404 }
405
406 //
407
408 //
409 for (Int_t i = 0; i < CAMERA_PIXELS; i++)
410 {
411 memset(sig[i], 0, (Int_t)(fSlices_mFadc*sizeof(Float_t))); // Added 15 01 2004, AM
412 memset(sig_LG[i], 0, (Int_t)(fSlices_mFadc*sizeof(Float_t))); // Added 01 03 2005, AM
413 }
414}
415void MFadc::Fill( Int_t iPix, Float_t time,
416 Float_t amplitude, Int_t isinner ) {
417
418 // AM, Jan 2004 : added delay to shift the signal peak to the desired
419 // range in the FADC window (indicated through the trigger_delay command
420 // in the camera input card.
421
422 time += fadc_time_offset;
423
424 if(isinner)
425 Fill(iPix, time, amplitude);
426 else
427 FillOuter(iPix, time, amplitude);
428
429}
430void MFadc::Fill( Int_t iPix, Float_t time, Float_t amplitude ) {
431
432 //
433 // fills the information about one single Phe in the Trigger class
434 //
435 // Parameters are the number of the pixel and the time-difference to the
436 // first photon.
437 //
438 //
439 // AM, Jan 2004: Replaced former FADC simulation (integration of signal)
440 // with a more realistic one (measuring signal height at discrete points).
441 //
442
443
444 Int_t i, ichan, ichanfadc ;
445
446 //
447 // first we have to check if the pixel iPix is used or not until now
448 // if this is the first use, reset all signal for that pixel
449 //
450 if ( iPix > numpix )
451 {
452 cout << " WARNING: MFadc::Fill() : iPix greater than Pixels in Camera = "
453 << numpix
454 << endl;
455 exit(987);
456 }
457
458 if ( used[iPix] == FALSE )
459 {
460 used [iPix] = TRUE;
461
462 for (i=0; i < (Int_t) fSlices_mFadc; i++ )
463 {
464 sig[iPix][i] = 0.;
465 sig_LG[iPix][i] = 0.;
466 }
467 }
468
469 //
470 // then select the time slice to use (ichan)
471 //
472
473 if ( time < TOTAL_TRIGGER_TIME+fadc_time_offset ) {
474 //
475 // Convert time into units of the width of the analog
476 // signal histogram, sing_resp:
477 //
478 ichan = (Int_t) ( time * fFadcSlicesPerNanosec * SUBBINS);
479
480 //
481 // putting the response slices in the right sig slices.
482 // Be careful, because both slices have different widths.
483 //
484
485 // We want to put the single phe response given by sing_resp into the
486 // array sig[][], but only one of each SUBBINS bins, since the binning
487 // of sing_resp is finer than that of sig[][]. We want that the start of
488 // sing_resp coincides with the time "time" with respect to the begining
489 // of sig[][]
490
491 // We take the pulse height in the middle of FADC slices, we start in the
492 // first such point after the time "time" (=ichan in response bins). Each
493 // FADC slice corresponds to SUBBINS response bins.
494
495 Int_t first_i = Int_t(SUBBINS/2) - ichan%(Int_t)SUBBINS;
496 first_i = first_i < 0 ? (Int_t)SUBBINS+first_i : first_i; //
497 //
498 // first_i is the first bin of sing_resp which matches the center of one
499 // bin of sig[][]
500 //
501
502 for ( i = first_i ; i < (Int_t)fResponseSlicesFadc; i += (Int_t)SUBBINS)
503 {
504 ichanfadc = (Int_t) ((ichan+i)/SUBBINS) ;
505 if ( ichanfadc < 0 )
506 continue;
507
508 //
509 // fSlices_mFadc is by default 48. sig[][] is not the final FADC output; that
510 // will be later filled (from sig[][]) in MFadc::TriggeredFadc()
511 //
512 if ( (ichanfadc) < (Int_t) fSlices_mFadc )
513 {
514 sig[iPix][ichanfadc] += (amplitude * sing_resp[i] );
515 sig_LG[iPix][ichanfadc] += (amplitude * sing_resp_lowgain[i] );
516 }
517 }
518 }
519 else
520 cout << " WARNING! Fadc::Fill " << time << " out of TriggerTimeRange "
521 << TOTAL_TRIGGER_TIME+fadc_time_offset << endl ;
522
523}
524
525void MFadc::FillOuter( Int_t iPix, Float_t time, Float_t amplitude ) {
526
527 //
528 // fills the information about one single Phe in the Trigger class
529 // for an outer pixel
530 //
531 // See explanations of the code in function Fill() above
532 //
533
534 Int_t i, ichan, ichanfadc ;
535
536 if ( iPix > numpix )
537 {
538 cout << " WARNING: MFadc::FillOuter() : iPix greater than CAMERA_PIXELS"
539 << endl ;
540 exit(987) ;
541 }
542
543 if ( used[iPix] == FALSE )
544 {
545 used [iPix] = TRUE ;
546
547 for (i=0; i < (Int_t) fSlices_mFadc; i++)
548 {
549 sig[iPix][i] = 0.;
550 sig_LG[iPix][i] = 0.;
551 }
552 }
553
554
555 if ( time < TOTAL_TRIGGER_TIME+fadc_time_offset ) {
556
557 ichan = (Int_t) ( time * fFadcSlicesPerNanosec * SUBBINS);
558
559 Int_t first_i = Int_t(SUBBINS/2) - ichan%(Int_t)SUBBINS;
560 first_i = first_i < 0 ? (Int_t)SUBBINS+first_i : first_i;
561
562 for ( i = first_i ; i < (Int_t)fResponseSlicesFadc; i += (Int_t)SUBBINS)
563 {
564 ichanfadc = (Int_t) ((ichan+i)/SUBBINS);
565
566 if ( ichanfadc < 0 )
567 continue;
568
569 if ( (ichanfadc) < (Int_t)fSlices_mFadc )
570 {
571 sig[iPix][ichanfadc] += (amplitude * sing_resp_outer[i] );
572 sig_LG[iPix][ichanfadc] += (amplitude * sing_resp_outer_lowgain[i] );
573 }
574
575 }
576
577 }
578 else {
579 cout << " WARNING! Fadc::FillOuter " << time << " out of TriggerTimeRange "
580 << TOTAL_TRIGGER_TIME+fadc_time_offset << endl ;
581 }
582
583}
584
585void MFadc::AddSignal( Int_t iPix, Float_t *resp) {
586
587 //
588 // Adds signals to the fadc reponse from a given array
589 // Parameters are the number of the pixel and the values to be added
590 // With add the signal equally to the high and low gain branch. The low
591 // gain branch is not yet scaled down!!
592 //
593
594 Int_t i ;
595
596 //
597 // first we have to check if the pixel iPix is used or not until now
598 // if this is the first use, reset all signal for that pixels
599 //
600 if ( iPix > numpix )
601 {
602 cout << " WARNING: MFadc::Fill() : iPix greater than CAMERA_PIXELS"
603 << endl ;
604 exit(987) ;
605 }
606
607 if ( used[iPix] == FALSE )
608 {
609 used [iPix] = TRUE ;
610
611 for (i=0; i < (Int_t)fSlices_mFadc; i++ )
612 {
613 sig[iPix][i] = 0. ;
614 sig_LG[iPix][i] = 0. ;
615 }
616 }
617 for ( i = 0 ; i<(Int_t)fSlices_mFadc; i++ )
618 {
619 sig[iPix][i] += resp[i] ;
620 sig_LG[iPix][i] += resp[i] ;
621 }
622
623}
624
625void MFadc::SetPedestals( Int_t ped) {
626 // It sets pedestal for each pixel flat randomly dstributed between 0 and ped
627 // It uses the instance of TRandom GenElec.
628
629 Int_t i;
630
631 for(i=0;i<numpix;i++){
632 pedestal[i]= (Float_t)(ped* GenElec->Rndm());
633 }
634}
635
636void MFadc::SetPedestals( Float_t *ped) {
637 // It sets pedestal for each pixel from ped array
638
639 Int_t i;
640
641 for(i=0;i<numpix;i++){
642 pedestal[i]= ped[i];
643 }
644}
645
646
647void MFadc::Pedestals(){
648 //
649 // It shifts the FADC contents their pedestal values
650 // It shifts the values in the analog signal,
651 // therefore it has to be done before getting FADC output
652 //
653
654 Int_t i, j;
655
656 for(i=0;i<numpix;i++)
657 for(j=0;j<(Int_t)fSlices_mFadc;j++)
658 {
659 sig[i][j] += pedestal[i];
660 sig_LG[i][j] += pedestal[i];
661 }
662
663 //
664 // AM 15 01 2003: Formerly the above operation was performed only
665 // for pixels in which used[] was true. But to run camera with no noise
666 // and get the right baseline on the pixels with no C-photons, we have
667 // to do it for all pixels.
668 //
669}
670
671
672void MFadc::SetElecNoise(Float_t value1, Float_t value2, UInt_t n_in_pix){
673
674 UInt_t i;
675
676 fInnerPixelsNum = n_in_pix;
677
678 cout<<"MFadc::SetElecNoise ... generating database for electronic noise."
679 <<endl;
680
681 for (i=0;i<(UInt_t (fSlices_mFadc))*1001;i++){
682 noise[i]=GenElec->Gaus(0., value1 );
683 noise_outer[i]=GenElec->Gaus(0., value2 );
684 }
685
686 cout<<"MFadc::SetElecNoise ... done"<<endl;
687
688}
689
690void MFadc::ElecNoise() {
691 // ============================================================
692 //
693 // adds the noise due to optronics and electronics
694 // to the signal. This is noise which comes before the FADC,
695 // so it will be later scaled down in the low gain branch, if
696 // the switch to low gain occurs.
697 //
698 UInt_t startslice;
699
700 for ( Int_t i = 0 ; i < numpix; i++) {
701 //
702 // but at the beginning we must check if this pixel is
703 // hitted the first time
704 //
705
706 startslice = GenElec->Integer(((Int_t)fSlices_mFadc)*1000);
707
708 if ( used[i] == FALSE )
709 {
710 used [i] = TRUE ;
711 if (i < fInnerPixelsNum)
712 {
713 memcpy( (Float_t*)&sig[i][0],
714 (Float_t*)&noise[startslice],
715 ((Int_t) fSlices_mFadc)*sizeof(Float_t));
716 memcpy( (Float_t*)&sig_LG[i][0],
717 (Float_t*)&noise[startslice],
718 ((Int_t) fSlices_mFadc)*sizeof(Float_t));
719 }
720
721 else
722 {
723 memcpy( (Float_t*)&sig[i][0],
724 (Float_t*)&noise_outer[startslice],
725 ((Int_t) fSlices_mFadc)*sizeof(Float_t));
726 memcpy( (Float_t*)&sig_LG[i][0],
727 (Float_t*)&noise_outer[startslice],
728 ((Int_t) fSlices_mFadc)*sizeof(Float_t));
729 }
730 }
731
732 //
733 // If pixel already in use, the noise is added each time slice
734 //
735 else
736 {
737 if (i < fInnerPixelsNum)
738 for ( Int_t is = 0 ; is < (Int_t)fSlices_mFadc ; is++ )
739 {
740 sig[i][is] += noise[startslice+is];
741 sig_LG[i][is] += noise[startslice+is];
742 }
743 else
744 for ( Int_t is = 0 ; is < (Int_t)fSlices_mFadc ; is++ )
745 {
746 sig[i][is] += noise_outer[startslice+is];
747 sig_LG[i][is] += noise_outer[startslice+is];
748 }
749 }
750 }
751}
752
753void MFadc::SetDigitalNoise(Float_t value){
754
755 UInt_t i;
756 Float_t xrdm;
757
758 cout<<"MFadc::SetDigitalNoise ... generating database for electronic noise."
759 <<endl;
760
761 for (i=0;i<UInt_t(fSlices_mFadc*1001);i++){
762 xrdm=GenElec->Gaus(0., value);
763 digital_noise[i]=(xrdm>0?Int_t(xrdm+0.5):Int_t(xrdm-0.5));
764 }
765
766 cout<<"MFadc::SetDigitalNoise ... done"<<endl;
767
768}
769
770void MFadc::DigitalNoise() {
771 // ============================================================
772 //
773 // adds the noise due to FADC electronics to the signal. This
774 // noise affects equally the high and low gain branches, that is,
775 // it is not scaled down in the low gain branch.
776 //
777 UInt_t startslice;
778
779 for ( Int_t i = 0 ; i < numpix; i++)
780 {
781 if ( used[i] == FALSE )
782 continue;
783
784 startslice=GenElec->Integer((Int_t) fSlices_mFadc*999);
785 //
786 // Then the noise is introduced for each time slice
787 //
788 for ( Int_t is = 0 ; is < fFadcSlices; is++ )
789 {
790 output[i][is] += digital_noise[startslice+is];
791 output_lowgain[i][is] += digital_noise[startslice+fFadcSlices+is];
792 }
793 }
794}
795
796void MFadc::Scan() {
797
798
799 for ( Int_t ip=0; ip<numpix; ip++ ) {
800
801 if ( used[ip] == kTRUE ) {
802
803 printf ("Pid %3d", ip ) ;
804
805 for ( Int_t is=0 ; is < (Int_t)fSlices_mFadc; is++ ) {
806
807 if ( sig[ip][is] > 0. ) {
808 printf (" %4.1f/", sig[ip][is] ) ;
809 }
810 else {
811 printf ("----/" ) ;
812 }
813 }
814
815 printf ("\n");
816
817 }
818 }
819
820}
821
822void MFadc::Scan(Float_t time) {
823
824 //
825 // first of all we subtract from the time a offset (8 ns)
826 //
827
828 Float_t t ;
829
830 (0 > time - TIME_BEFORE_TRIGGER)? t=fadc_time_offset: t=(time-TIME_BEFORE_TRIGGER+fadc_time_offset) ; // to show also the start of the pulse before the trigger time
831
832 if ( t < 0. ) {
833 cout << " WARNING!! FROM MFADC::SCAN(t) " << endl ;
834 exit (776) ;
835 }
836
837 //
838 // calculate the first slice to write out
839 //
840
841 Int_t iFirstSlice ;
842
843 iFirstSlice = (Int_t) ( t * fFadcSlicesPerNanosec ) ;
844
845 for ( Int_t ip=0; ip<numpix; ip++ ) {
846
847 if ( used[ip] == kTRUE ) {
848
849 printf ("Pid %3d", ip ) ;
850
851 for ( Int_t is=iFirstSlice ; is < (iFirstSlice+15); is++ ) {
852 printf (" %5.2f /", sig[ip][is] ) ;
853 }
854
855 printf ("\n");
856
857 }
858 }
859}
860
861void MFadc::GetResponse( Float_t *resp ) {
862 // ============================================================
863 //
864 // puts the standard response function into the array resp
865
866 for ( Int_t i=0; i< fResponseSlicesFadc; i++ )
867 resp[i] = sing_resp[i];
868
869}
870
871void MFadc::GetPedestals( Float_t *offset) {
872 // ============================================================
873 //
874 // puts the pedestal values into the array offset
875
876 for ( Int_t i=0; i< numpix; i++ ) {
877
878 offset[i] = pedestal[i] ;
879 }
880}
881
882//===========================================================================
883//
884// Next function adds up the noise in pixel "pix", scaling down the part
885// of it which comes from before the receivers in the case we are dealing with
886// low gain (ishigh=0). The output is the sum of the readouts of a number
887// n_slices of FADC slices. For the case of low gain, the FADC contents we add
888// are not what we would have in a real pedestal event, but nevertheless this
889// is useful in the camera simulation to obtain what the pedestal fluctuations
890// are for the low gain. This will be written to the camera output, in the
891// MMcFadcHeader.
892//
893Float_t MFadc::AddNoiseInSlices( Int_t pix, Int_t ishigh, Int_t n_slices) {
894
895 Float_t sum = 0;
896 Float_t fvalue = 0.;
897 UChar_t value = 0;
898 UInt_t startslice;
899
900 //
901 // Get at random a point in the FADC presimulated digital noise:
902 //
903 startslice=GenElec->Integer((Int_t) fSlices_mFadc*999);
904
905 for ( Int_t is=0; is < n_slices ; is++ )
906 {
907 if (ishigh)
908 fvalue = sig[pix][is];
909 else
910 // If we deal with low gain, we have to scale the values in sig_LG[][] by
911 // the gain ratio (high2low_gain), since "sig_LG" contains at this point the
912 // noise produced before the receiver boards (for instance NSB noise).
913 //
914 fvalue = pedestal[pix]+(sig_LG[pix][is]-pedestal[pix]) / high2low_gain;
915
916 fvalue += digital_noise[startslice+is]; // We add the noise intrinsic to FADC
917
918 fvalue = fvalue < 0? fvalue-0.5 : fvalue+0.5;
919
920 value = fvalue < 0.? (UChar_t) 0 :
921 (fvalue > 255.? 255 : (UChar_t) fvalue);
922
923 // Add up slices:
924 sum += value - pedestal[pix];
925 }
926
927 return sum;
928}
929
930//=======================================================================
931
932void MFadc::TriggeredFadc(Float_t time) {
933
934 //
935 // Here the slices to write out are calculated. Warning: the digitalization
936 // is NOT done here (it is already done in MFadc::Fill). This procedure only
937 // selects which FADC slices to write out, out of those contained in the sig[][]
938 // and sig_LG[][] arrays.
939 //
940
941 //
942 // calculate the first slice to write out, according to trigger time:
943 //
944
945 Int_t iFirstSlice ;
946 Int_t i;
947
948 //
949 // We had 0.5 for the correct rounding:
950 //
951 iFirstSlice = (Int_t) ( 0.5 + time * fFadcSlicesPerNanosec ) ;
952
953 for ( Int_t ip = 0; ip < numpix; ip++ )
954 {
955
956 if ( used[ip] == kFALSE)
957 // Pixels with no C-photons, in the case that camera is being run with
958 // no noise (nor NSB neither electronic). We then set the mean pedestal as
959 // signal, since when analyzing the camera output file, MARS will subtract
960 // it anyway!
961 {
962 for ( Int_t i=0 ; i < fFadcSlices ; i++ )
963 {
964 output[ip][i] = pedestal[ip];
965 output_lowgain[ip][i] = pedestal[ip];
966 }
967 continue;
968 }
969
970
971 // First put the high gain in the output slices:
972 i = 0;
973 Int_t switch_i = 0;
974 for ( Int_t is = iFirstSlice; is < (iFirstSlice+fFadcSlices); is++ )
975 {
976 if (is < (Int_t)fSlices_mFadc)
977 {
978 output[ip][i] = sig[ip][is];
979
980 if (switch_i == 0) // Hi gain limit not yet surpassed before.
981 {
982 if (output[ip][i] > fGainSwitchAmp)
983 switch_i = i + fShiftFromSwitch2LowGain;
984 }
985 }
986
987 else // We are beyond the simulated signal history in sig[][]! Put just mean pedestal!
988 output[ip][i] = pedestal[ip];
989
990 i++;
991 }
992
993 // Now put the low gain:
994 i=0;
995 for ( Int_t is = iFirstSlice; is < (iFirstSlice+fFadcSlices); is++ )
996 {
997 if (is < (Int_t)fSlices_mFadc)
998 {
999 if (switch_i > 0 && (i+fFadcSlices) >= switch_i)
1000 output_lowgain[ip][i] = pedestal[ip] +
1001 (sig_LG[ip][is]-pedestal[ip])/high2low_gain;
1002 // Once the shift occurs, low gain is filled with the high
1003 // gain signal scaled down by the factor high2low_gain
1004
1005 else
1006 output_lowgain[ip][i] = sig[ip][is+fFadcSlices];
1007 // Write out high gain into low gain slices if there was no
1008 // switch, or before the switch occurs.
1009 }
1010
1011 else // We are beyond the simulated signal history in sig[][]! Put just mean pedestal!
1012 {
1013 output_lowgain[ip][i] = pedestal[ip];
1014 }
1015 i++;
1016 }
1017
1018 }
1019}
1020
1021
1022void MFadc::ShowSignal (MMcEvt *McEvt, Float_t trigTime) {
1023 // ============================================================
1024 //
1025 // This method is used to book the histogram to show the signal in
1026 // a special gui frame (class MGTriggerSignal). After the look onto the
1027 // signals for a better understanding of the things we will expect
1028 // the gui frame and all histogramms will be destroyed.
1029 //
1030
1031 //
1032 // first of all create a list of the histograms to show
1033 //
1034 // take only that one with a entry
1035
1036 TH1F *hist ;
1037 Char_t dumm[10];
1038 Char_t name[256];
1039
1040 TObjArray *AList ;
1041 AList = new TObjArray(10) ;
1042
1043 // the list of analog signal histograms
1044 // at the beginning we initalise 10 elements
1045 // but this array expand automatically if neccessay
1046
1047 Int_t ic = 0 ;
1048 for ( Int_t i=0 ; i < numpix; i++ ) {
1049 if ( used [i] == TRUE ) {
1050
1051 sprintf (dumm, "FADC_%d", i ) ;
1052 sprintf (name, "fadc signal %d", i ) ;
1053
1054 hist = new TH1F(dumm, name, (Int_t)fSlices_mFadc, fadc_time_offset, TOTAL_TRIGGER_TIME+fadc_time_offset);
1055 //
1056 // fill the histogram
1057 //
1058
1059 for (Int_t ibin = 1; ibin <= (Int_t)fSlices_mFadc; ibin++)
1060 hist->SetBinContent (ibin, sig[i][ibin-1]);
1061
1062
1063 // hist->SetMaximum( 5.);
1064 // hist->SetMinimum(-10.);
1065 hist->SetStats(kFALSE);
1066
1067 // hist->SetAxisRange(0., 80. ) ;
1068
1069 AList->Add(hist) ;
1070
1071 ic++ ;
1072 }
1073 }
1074
1075 //
1076 // create the Gui Tool
1077 //
1078 //
1079
1080 new MGFadcSignal(McEvt,
1081 AList,
1082 trigTime,
1083 gClient->GetRoot(),
1084 gClient->GetRoot(),
1085 400, 400 ) ;
1086
1087 //
1088 // delete the List of histogramms
1089 //
1090 AList->Delete() ;
1091
1092 delete AList ;
1093}
1094
1095UChar_t MFadc::GetFadcSignal(Int_t pixel, Int_t slice){
1096
1097 // It returns the analog signal for a given pixel and a given FADC
1098 // time slice which would be read.
1099
1100 // Since May 1 2004, we do the rounding and the truncation to the range
1101 // 0-255 counts here. (A. Moralejo)
1102
1103 Float_t out = output[pixel][slice] > 0. ?
1104 output[pixel][slice]+0.5 : output[pixel][slice]-0.5;
1105 // (add or subtract 0.5 for correct rounding)
1106
1107 return (out < 0.? (UChar_t) 0 :
1108 (out > 255.? (UChar_t) 255 :
1109 (UChar_t) out));
1110}
1111
1112
1113UChar_t MFadc::GetFadcLowGainSignal(Int_t pixel, Int_t slice){
1114
1115 // It returns the analog signal for a given pixel and a given FADC
1116 // time slice which would be read. Same comment as above.
1117
1118 Float_t outlow = output_lowgain[pixel][slice] > 0. ?
1119 output_lowgain[pixel][slice]+0.5 :
1120 output_lowgain[pixel][slice]-0.5;
1121 // (add or subtract 0.5 for correct rounding)
1122
1123 return (outlow < 0.? (UChar_t) 0 :
1124 (outlow > 255.? (UChar_t) 255 :
1125 (UChar_t) outlow));
1126}
1127
1128
1129
Note: See TracBrowser for help on using the repository browser.