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

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