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

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