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

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