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