source: trunk/MagicSoft/Simulation/Detector/include-MTrigger/MTrigger.cxx@ 18477

Last change on this file since 18477 was 5319, checked in by moralejo, 20 years ago
*** empty log message ***
File size: 43.2 KB
Line 
1/////////////////////////////////////////////////////////////////
2//
3// MTrigger
4//
5//
6#include "MTrigger.hxx"
7
8#include "TROOT.h"
9#include "TFile.h"
10#include "TH1.h"
11#include "TObjArray.h"
12#include "MGTriggerSignal.hxx"
13#include "MGeomCam.h"
14#include "MGeomPix.h"
15
16using namespace std;
17
18MTrigger::MTrigger(int pix) {
19 // ============================================================
20 //
21 // default constructor
22 //
23 // The procedure is the following:
24 //
25 // 1. Allocation of some memory needed
26 // 2. some parameters of the trigger are set to default.
27 // 3. if a File MTrigger.card exists in the current directory,
28 // this parameters of the trigger may be changed
29 // 4. Then the all signals are set to zero
30
31 FILE *unit_mtrig ;
32 Int_t endflag = 1 ;
33 Int_t bthresholdpixel = FALSE;
34 char datac[256] ;
35 char dummy[50] ;
36 char input_thres[50];
37 Int_t i, ii ;
38
39 Float_t threshold ;
40
41 // Number of pixels in the trigger region
42 pixnum=pix;
43
44 //
45 // allocate the memory for the 2dim arrays (a_sig, d_sig )
46 //
47
48 used = new Bool_t[pix];
49 nphotshow = new Int_t[pix];
50 nphotnsb = new Int_t[pix];
51 nphotstar = new Int_t[pix];
52 a_sig = new Float_t * [pix];
53 d_sig = new Float_t * [pix];
54 baseline = new Float_t[pix];
55 dknt = new Bool_t[pix];
56 noise = new Float_t[TRIGGER_TIME_SLICES*1001];
57 chan_thres = new Float_t[pix];
58 for(Int_t j=0;j<6;j++)
59 NN[j] = new Int_t[pix];
60 for(Int_t j=0;j<TRIGGER_CELLS;j++)
61 TC[j] = new Int_t[pix];
62
63 for( Int_t j=0; j<pix; j++ ) {
64
65 a_sig[j] = new Float_t[TRIGGER_TIME_SLICES] ;
66
67 d_sig[j] = new Float_t[TRIGGER_TIME_SLICES] ;
68 }
69
70
71
72 //
73 // set the values for the standard response pulse
74 //
75
76 fwhm_resp = RESPONSE_FWHM;
77 ampl_resp = RESPONSE_AMPLITUDE;
78
79 overlaping_time = TRIGGER_OVERLAPING;
80
81 threshold = CHANNEL_THRESHOLD;
82
83
84 gate_leng = TRIGGER_GATE;
85 trigger_multi = TRIGGER_MULTI;
86 trigger_geometry = TRIGGER_GEOM;
87
88 //
89 // check if the file MTrigger.card exists
90 //
91
92 if ( (unit_mtrig = fopen ("MTrigger.card", "r")) != 0 ) {
93 cout << "[MTrigger] use the values from MTrigger.card "<< endl ;
94
95 while ( endflag == 1 ) {
96 //
97 //
98 fgets (datac, 255, unit_mtrig) ;
99 // printf ("--> %s <--", datac ) ;
100
101 //
102 // now compare the line with controlcard words
103 //
104
105 if ( strncmp (datac, "channel_threshold", 17 ) == 0 ) {
106 sscanf (datac, "%s %f", dummy, &threshold ) ;
107 }
108 else if ( strncmp (datac, "gate_length", 11 ) == 0 ) {
109 sscanf (datac, "%s %f", dummy, &gate_leng ) ;
110 }
111 else if ( strncmp (datac, "response_fwhm", 13 ) == 0 ) {
112 sscanf (datac, "%s %f", dummy, &fwhm_resp ) ;
113 }
114 else if ( strncmp (datac, "response_ampl", 13 ) == 0 ) {
115 sscanf (datac, "%s %f", dummy, &ampl_resp ) ;
116 }
117 else if ( strncmp (datac, "overlaping", 10 ) == 0 ) {
118 sscanf (datac, "%s %f", dummy, &overlaping_time ) ;
119 }
120 else if ( strncmp (datac, "multiplicity", 12 ) == 0 ) {
121 sscanf (datac, "%s %f", dummy, &trigger_multi ) ;
122 }
123 else if ( strncmp (datac, "topology", 8 ) == 0 ) {
124 sscanf (datac, "%s %i", dummy, &trigger_geometry ) ;
125 }
126 else if ( strncmp (datac, "threshold_file", 14 ) == 0 ) {
127 sscanf (datac, "%s %s", dummy, input_thres ) ;
128 bthresholdpixel=kTRUE;
129 }
130
131 if ( feof(unit_mtrig) != 0 ) {
132 endflag = 0 ;
133 }
134
135 }
136
137 fclose ( unit_mtrig ) ;
138 }
139 else {
140 cout << "[MTrigger] use the standard values for MTrigger "<< endl ;
141 }
142
143 cout << endl
144 << "[MTrigger] Setting up the MTrigger with this values "<< endl ;
145 if(bthresholdpixel){
146 cout<<endl
147 << "[MTrigger] ChannelThreshold from file: "<<input_thres
148 <<endl;
149 }
150 else{
151 cout << endl
152 << "[MTrigger] ChannelThreshold: " << threshold << " mV"
153 << endl ;
154 }
155 cout << "[MTrigger] Gate Length: " << gate_leng << " ns"
156 << endl ;
157 cout << "[MTrigger] Overlaping time: " << overlaping_time << " ns"
158 << endl ;
159 cout << "[MTrigger] Response FWHM: " << fwhm_resp << " ns"
160 << endl ;
161 cout << "[MTrigger] Response Amplitude: " << ampl_resp << " mV"
162 << endl ;
163 cout << "[MTrigger] Trigger Multiplicity: " << trigger_multi << " pixels"
164 << endl ;
165 cout << "[MTrigger] Trigger Topology: " << trigger_geometry
166 << endl ;
167
168 cout << endl ;
169
170
171 //
172 // we have introduced individual thresholds for all pixels
173 //
174 FILE *unit_thres;
175
176 if (bthresholdpixel == kTRUE) {
177 if ((unit_thres=fopen(input_thres, "r"))==0){
178 cout<<"WARNING: not able to read ..."<<input_thres<<endl;
179 cout<<"Threshold will be set to "<<threshold<<" for all pixels"<<endl;
180 for (Int_t k=0; k<pix; k++ ) {
181 chan_thres[k] = threshold ;
182 }
183 }
184 else {
185 for (i=0;i<pix;i++){
186 fscanf(unit_thres, "%f",&chan_thres[i]);
187 }
188 fclose (unit_thres);
189 }
190 }
191 else {
192 for (Int_t k=0; k<pix; k++ ) {
193 chan_thres[k] = threshold ;
194 }
195 }
196
197
198 //
199 // set up the response shape
200 //
201
202 Float_t sigma ;
203 Float_t x, x0 ;
204
205 sigma = fwhm_resp / 2.35 ;
206 x0 = 3*sigma ;
207
208 for (i=0; i< RESPONSE_SLICES_TRIG ; i++ ) {
209
210 x = i * (1./((Float_t)TRIG_SLICES_PER_NSEC))
211 + (1./( 2 * (Float_t)TRIG_SLICES_PER_NSEC )) ;
212
213 sing_resp[i] = (Float_t)
214 ampl_resp * expf(-0.5 * (x-x0)*(x-x0) / (sigma*sigma) ) ;
215
216 }
217
218 //
219 // look for the time between start of response function and the
220 // maximum value of the response function. This is needed by the
221 // member functions FillNSB() and FillStar()
222 //
223
224 Int_t imax = 0 ;
225 Float_t max = 0. ;
226 for (i=0; i< RESPONSE_SLICES_TRIG ; i++ ) {
227 if ( sing_resp[i] > max ) {
228 imax = i ;
229 max = sing_resp[i] ;
230 }
231 }
232
233 peak_time = ( (Float_t) imax ) / ( (Float_t) TRIG_SLICES_PER_NSEC ) ;
234
235
236 // Initialize histogram histPmt which represents the distribution of amplitudes
237 // of the PMT response to a single photpelectron.
238
239 InitGainFluctuations();
240
241
242 //
243 // create the random generator for the Electronic Noise
244 //
245
246 GenElec = new TRandom() ;
247
248 GenElec->SetSeed(0);
249
250 //
251 // Read in the lookup table for NN trigger
252 //
253
254 FILE *unit ;
255 int id ;
256
257 i = 0 ;
258
259 if ( (unit = fopen("../include-MTrigger/TABLE_NEXT_NEIGHBOUR", "r" )) == 0 ) {
260 cout << "ERROR: not able to read ../include-MTrigger/TABLE_NEXT_NEIGHBOUR"
261 << endl ;
262 exit(123) ;
263 }
264 else {
265 while ( i < pix )
266 {
267 fscanf ( unit, " %d", &id ) ;
268
269 for ( Int_t k=0; k<6; k++ ) {
270 fscanf ( unit, "%d ", &NN[k][i] ) ;
271 }
272 i++ ;
273 }
274
275 fclose (unit) ;
276 }
277
278
279 //
280 // Read in the lookup table for trigger cells
281 //
282
283 i = 0 ;
284
285 if ( (unit = fopen("../include-MTrigger/TABLE_PIXELS_IN_CELLS", "r" )) == 0 ) {
286 cout << "ERROR: not able to read ../include-MTrigger/TABLE_PIXELS_IN_CELLS"
287 << endl ;
288 exit(123) ;
289 }
290 else {
291 while ( i < pix )
292 {
293 for ( Int_t k=0; k<TRIGGER_CELLS; k++ ) {
294 TC[k][i]=kFALSE;
295 }
296 i++ ;
297 }
298 while ( feof(unit) == 0 ) {
299 for ( Int_t k=0; k<TRIGGER_CELLS; k++ ) {
300 fscanf ( unit, "%d ", &i ) ;
301 if ((i-1)<pix)
302 TC[k][i-1]=kTRUE;
303 }
304 }
305 fclose (unit) ;
306 }
307
308
309 //
310 //
311 // set all the booleans used to FALSE, indicating that the pixel is not
312 // used in this event.
313 //
314
315 for ( i =0 ; i <pix ; i++ ) {
316 used [i] = kFALSE ;
317 dknt [i] = kFALSE ;
318
319 nphotshow[i] = 0 ;
320 nphotnsb [i] = 0 ;
321 nphotstar[i] = 0 ;
322
323 baseline[i] = 0 ;
324 }
325
326 for ( ii=0 ; ii<TRIGGER_TIME_SLICES; ii++ ) {
327 sum_d_sig[ii] = 0. ;
328 }
329
330 //
331 // set the information about the Different Level Triggers to zero
332 //
333
334 nZero = nFirst = nSecond = 0 ;
335
336 for (ii=0 ; ii<TRIGGER_TIME_SLICES; ii++ ) {
337 SlicesZero[ii] = kFALSE;
338 }
339
340 for ( i = 0 ; i < 5 ; i++) {
341 SlicesFirst[i] = 0;
342 SlicesSecond[i] = 0;
343 PixelsFirst[i] = -1;
344 PixelsSecond[i] = -1;
345 }
346
347 cout << " end of MTrigger::MTrigger()" << endl ;
348}
349
350MTrigger::MTrigger(Int_t pix, MGeomCam *camgeom,
351 float gate, float overt, float ampl, float fwhm, int ct_id) {
352 // ============================================================
353 //
354 // constructor
355 //
356 // The procedure is the following:
357 //
358 // 1. Allocation of some memory needed
359 // 2. some parameters of the trigger are set.
360 // 3. Then the all signals are set to zero
361
362 Int_t i, ii ;
363
364 Float_t threshold ;
365
366 // Number of pixels in the trigger region
367 pixnum=pix;
368
369 //
370 // allocate the memory for the 2dim arrays (a_sig, d_sig )
371 //
372
373 used = new Bool_t[pix];
374 nphotshow = new Int_t[pix];
375 nphotnsb = new Int_t[pix];
376 nphotstar = new Int_t[pix];
377 a_sig = new Float_t * [pix];
378 d_sig = new Float_t * [pix];
379 baseline = new Float_t[pix];
380 dknt = new Bool_t[pix];
381 noise = new Float_t[TRIGGER_TIME_SLICES*1001];
382 chan_thres = new Float_t[pix];
383 for(Int_t j=0;j<6;j++)
384 NN[j] = new Int_t[pix];
385 for(Int_t j=0;j<TRIGGER_CELLS;j++)
386 TC[j] = new Int_t[pix];
387
388 for( Int_t j=0; j<pix; j++ ) {
389
390 a_sig[j] = new Float_t[TRIGGER_TIME_SLICES] ;
391
392 d_sig[j] = new Float_t[TRIGGER_TIME_SLICES] ;
393 }
394
395 //
396 // set the values for the standard response pulse
397 //
398
399 fwhm_resp = fwhm ;
400 ampl_resp = ampl ;
401
402 overlaping_time = overt;
403
404
405 threshold = CHANNEL_THRESHOLD ;
406
407
408 gate_leng = gate ;
409 trigger_multi = TRIGGER_MULTI ;
410 trigger_geometry = TRIGGER_GEOM ;
411
412 cout << endl
413 << "[MTrigger] Setting up the MTrigger with this values "<< endl ;
414 cout << "[MTrigger] Gate Length: " << gate_leng << " ns"
415 << endl ;
416 cout << "[MTrigger] Overlaping time: " << overlaping_time << " ns"
417 << endl ;
418 cout << "[MTrigger] Response FWHM: " << fwhm_resp << " ns"
419 << endl ;
420 cout << "[MTrigger] Response Amplitude: " << ampl_resp << " mV"
421 << endl ;
422 cout << endl ;
423
424 for (Int_t k=0; k<pixnum; k++ ) {
425 chan_thres[k] = threshold ;
426 }
427
428 //
429 // set up the response shape
430 //
431
432 Float_t sigma ;
433 Float_t x, x0 ;
434
435 sigma = fwhm_resp / 2.35 ;
436 x0 = 3*sigma ;
437
438 for (i=0; i< RESPONSE_SLICES_TRIG ; i++ ) {
439
440 x = i * (1./((Float_t)TRIG_SLICES_PER_NSEC))
441 + (1./( 2 * (Float_t)TRIG_SLICES_PER_NSEC )) ;
442
443 sing_resp[i] = (Float_t)
444 ampl_resp * expf(-0.5 * (x-x0)*(x-x0) / (sigma*sigma) ) ;
445
446 }
447
448 //
449 // look for the time between start of response function and the
450 // maximum value of the response function. This is needed by the
451 // member functions FillNSB() and FillStar()
452 //
453
454 Int_t imax = 0 ;
455 Float_t max = 0. ;
456 for (i=0; i< RESPONSE_SLICES_TRIG ; i++ ) {
457 if ( sing_resp[i] > max ) {
458 imax = i ;
459 max = sing_resp[i] ;
460 }
461 }
462
463 peak_time = ( (Float_t) imax ) / ( (Float_t) TRIG_SLICES_PER_NSEC ) ;
464
465
466 // Initialize histogram histPmt which represents the distribution of amplitudes
467 // of the PMT response to a single photpelectron.
468
469 InitGainFluctuations();
470
471 //
472 // Change name of histogram (tag with CT id) to avoid confusion in the case of
473 // multi-telescope simulations.
474 //
475 char histname[32];
476 sprintf(histname, "histPmt_%d", ct_id);
477 histPmt->SetName(histname);
478
479
480 //
481 // create the random generator for the Electronic Noise
482 //
483
484 GenElec = new TRandom() ;
485
486 //
487 // Read in the lookup table for NN trigger
488 //
489
490 for(i=0; i < pixnum;i++ )
491 {
492 MGeomPix &pixel = (*camgeom)[i];
493 for ( Int_t k=0; k<6; k++ ) {
494 NN[k][i]=pixel.GetNeighbor(k);
495 }
496 }
497
498 //
499 // Read in the lookup table for trigger cells
500 //
501
502 FILE *unit;
503
504 i = 0 ;
505
506 if ( (unit = fopen("../include-MTrigger/TABLE_PIXELS_IN_CELLS", "r" )) == 0 ) {
507 cout << "ERROR: not able to read ../include-MTrigger/TABLE_PIXELS_IN_CELLS"
508 << endl ;
509 exit(123) ;
510 }
511 else {
512 while ( i < pixnum )
513 {
514 for ( Int_t k=0; k<TRIGGER_CELLS; k++ ) {
515 TC[k][i]=kFALSE;
516 }
517 i++ ;
518 }
519 while ( feof(unit) == 0 ) {
520 for ( Int_t k=0; k<TRIGGER_CELLS; k++ ) {
521 fscanf ( unit, "%d ", &i ) ;
522 if((i-1)<pixnum)
523 TC[k][i-1]=kTRUE;
524 }
525 }
526 fclose (unit) ;
527 }
528
529 //
530 //
531 // set all the booleans used to FALSE, indicating that the pixel is not
532 // used in this event.
533 //
534
535 for ( i =0 ; i <pixnum ; i++ ) {
536 used [i] = kFALSE ;
537 dknt [i] = kFALSE ;
538
539 nphotshow[i] = 0 ;
540 nphotnsb [i] = 0 ;
541 nphotstar[i] = 0 ;
542
543 baseline[i] = 0 ;
544 }
545
546 for ( ii=0 ; ii<TRIGGER_TIME_SLICES; ii++ ) {
547 sum_d_sig[ii] = 0. ;
548 }
549
550 //
551 // set the information about the Different Level Triggers to zero
552 //
553
554 nZero = nFirst = nSecond = 0 ;
555
556 for (ii=0 ; ii<TRIGGER_TIME_SLICES; ii++ ) {
557 SlicesZero[ii] = kFALSE;
558 }
559
560 for ( i = 0 ; i < 5 ; i++) {
561 SlicesFirst[i] = 0 ;
562 SlicesSecond[i] = 0 ;
563 PixelsFirst[i] = -1;
564 PixelsSecond[i] = -1;
565 }
566 cout << " end of MTrigger::MTrigger()" << endl ;
567}
568MTrigger::MTrigger(Int_t pix,
569 float gate, float overt, float ampl, float fwhm) {
570 // ============================================================
571 //
572 // constructor
573 //
574 // The procedure is the following:
575 //
576 // 1. Allocation of some memory needed
577 // 2. some parameters of the trigger are set.
578 // 3. Then the all signals are set to zero
579
580 Int_t i, ii ;
581
582 Float_t threshold ;
583
584 // Number of pixels in the trigger region
585 pixnum=pix;
586
587 //
588 // allocate the memory for the 2dim arrays (a_sig, d_sig )
589 //
590
591 used = new Bool_t[pix];
592 nphotshow = new Int_t[pix];
593 nphotnsb = new Int_t[pix];
594 nphotstar = new Int_t[pix];
595 a_sig = new Float_t * [pix];
596 d_sig = new Float_t * [pix];
597 baseline = new Float_t[pix];
598 dknt = new Bool_t[pix];
599 noise = new Float_t[TRIGGER_TIME_SLICES*1001];
600 chan_thres = new Float_t[pix];
601 for(Int_t j=0;j<6;j++)
602 NN[j] = new Int_t[pix];
603 for(Int_t j=0;j<TRIGGER_CELLS;j++)
604 TC[j] = new Int_t[pix];
605
606 for( Int_t j=0; j<pix; j++ ) {
607
608 a_sig[j] = new Float_t[TRIGGER_TIME_SLICES] ;
609
610 d_sig[j] = new Float_t[TRIGGER_TIME_SLICES] ;
611 }
612
613 //
614 // set the values for the standard response pulse
615 //
616
617 fwhm_resp = fwhm;
618 ampl_resp = ampl;
619
620 overlaping_time = overt;
621
622
623 threshold = CHANNEL_THRESHOLD;
624
625
626 gate_leng = gate;
627 trigger_multi = TRIGGER_MULTI;
628 trigger_geometry = TRIGGER_GEOM;
629
630 cout << endl
631 << "[MTrigger] Setting up the MTrigger with this values "<< endl ;
632 cout << "[MTrigger] Gate Length: " << gate_leng << " ns"
633 << endl ;
634 cout << "[MTrigger] Overlaping time: " << overlaping_time << " ns"
635 << endl ;
636 cout << "[MTrigger] Response FWHM: " << fwhm_resp << " ns"
637 << endl ;
638 cout << "[MTrigger] Response Amplitude: " << ampl_resp << " mV"
639 << endl ;
640 cout << endl ;
641
642 for (Int_t k=0; k<pixnum; k++ ) {
643 chan_thres[k] = threshold ;
644 }
645
646 //
647 // set up the response shape
648 //
649
650 Float_t sigma ;
651 Float_t x, x0 ;
652
653 sigma = fwhm_resp / 2.35 ;
654 x0 = 3*sigma ;
655
656 for (i=0; i< RESPONSE_SLICES_TRIG ; i++ ) {
657
658 x = i * (1./((Float_t)TRIG_SLICES_PER_NSEC))
659 + (1./( 2 * (Float_t)TRIG_SLICES_PER_NSEC )) ;
660
661 sing_resp[i] = (Float_t)
662 ampl_resp * expf(-0.5 * (x-x0)*(x-x0) / (sigma*sigma) ) ;
663
664 }
665
666 //
667 // look for the time between start of response function and the
668 // maximum value of the response function. This is needed by the
669 // member functions FillNSB() and FillStar()
670 //
671
672 Int_t imax = 0 ;
673 Float_t max = 0. ;
674 for (i=0; i< RESPONSE_SLICES_TRIG ; i++ ) {
675 if ( sing_resp[i] > max ) {
676 imax = i ;
677 max = sing_resp[i] ;
678 }
679 }
680
681 peak_time = ( (Float_t) imax ) / ( (Float_t) TRIG_SLICES_PER_NSEC ) ;
682
683
684 // Initialize histogram histPmt which represents the distribution of amplitudes
685 // of the PMT response to a single photpelectron.
686
687 InitGainFluctuations();
688
689
690 //
691 // create the random generator for the Electronic Noise
692 //
693
694 GenElec = new TRandom() ;
695
696 //
697 //
698 // set all the booleans used to FALSE, indicating that the pixel is not
699 // used in this event.
700 //
701
702 for ( i =0 ; i <pixnum ; i++ ) {
703 used [i] = kFALSE ;
704 dknt [i] = kFALSE ;
705
706 nphotshow[i] = 0 ;
707 nphotnsb [i] = 0 ;
708 nphotstar[i] = 0 ;
709
710 baseline[i] = 0 ;
711 }
712
713 for ( ii=0 ; ii<TRIGGER_TIME_SLICES; ii++ ) {
714 sum_d_sig[ii] = 0. ;
715 }
716
717 //
718 // set the information about the Different Level Triggers to zero
719 //
720
721 nZero = nFirst = nSecond = 0 ;
722
723 for (ii=0 ; ii<TRIGGER_TIME_SLICES; ii++ ) {
724 SlicesZero[ii] = kFALSE;
725 }
726
727 for ( i = 0 ; i < 5 ; i++) {
728 SlicesFirst[i] = 0 ;
729 SlicesSecond[i] = 0 ;
730 PixelsFirst[i] = -1;
731 PixelsSecond[i] = -1;
732 }
733 cout << " end of MTrigger::MTrigger()" << endl ;
734}
735
736MTrigger::~MTrigger() {
737 // ============================================================//
738 // destructor
739 //
740 int i;
741
742 delete histPmt ;
743
744 for(i=0;i<pixnum;i++){
745 //delete [] a_sig[i];
746 //delete [] d_sig[i];
747 }
748
749 delete GenElec;
750}
751
752
753void MTrigger::Reset() {
754 // ============================================================
755 //
756 // reset all values of the signals to zero
757 //
758 Int_t i, ii ;
759
760 for ( i =0 ; i <pixnum ; i++ ) {
761 used [i] = kFALSE ;
762 dknt [i] = kFALSE ;
763
764 nphotshow[i] = 0 ;
765 nphotnsb [i] = 0 ;
766 nphotstar[i] = 0 ;
767 }
768
769 for ( ii=0 ; ii<TRIGGER_TIME_SLICES; ii++ ) {
770 sum_d_sig[ii] = 0. ;
771 }
772}
773
774void MTrigger::ClearZero() {
775 //
776 // set the information about the Zero Level Trigger to zero
777 //
778
779 Int_t i;
780
781 nZero = 0 ;
782
783 for (i=0 ; i<TRIGGER_TIME_SLICES; i++ ) {
784 SlicesZero[i] = kFALSE;
785 }
786
787}
788
789void MTrigger::ClearFirst() {
790 //
791 // set the information about the First Level Trigger to zero
792 //
793
794 Int_t i;
795
796 nFirst = 0 ;
797
798 for ( i = 0 ; i < 5 ; i++) {
799 SlicesFirst[i] = 0 ;
800 PixelsFirst[i] = -1;
801 }
802}
803
804Float_t MTrigger::FillShow(Int_t iPix, Float_t time) {
805 // ============================================================
806 //
807 // Fills the information of one single Phe electron that
808 // comes from the shower
809 //
810
811 //
812 // First check the time
813 //
814
815 if ( time < 0. || time > TOTAL_TRIGGER_TIME ) {
816 cout << " WARNING: time of phe out of time range: " << time << endl;
817 return 0. ;
818 }
819 else {
820 return ( Fill( iPix, time, CASE_SHOW ) ) ;
821 }
822}
823
824Float_t MTrigger::FillNSB(Int_t iPix, Float_t time) {
825 // ============================================================
826 //
827 // Fills the information of one single Phe electron that
828 // comes from the NSB
829 //
830
831 //
832 // First check the time
833 //
834
835 if ( time < 0. || time > TOTAL_TRIGGER_TIME ) {
836 cout << " WARNING: time of phe out of time range: " << time << endl;
837 return 0. ;
838 }
839 else {
840 return ( Fill( iPix, time - peak_time, CASE_NSB ) ) ;
841 }
842}
843
844Float_t MTrigger::FillStar(Int_t iPix, Float_t time) {
845 // ============================================================
846 //
847 // Fills the information of one single Phe electron that
848 // comes from a star
849
850 //
851 // First check the time
852 //
853
854 if ( time < 0. || time > TOTAL_TRIGGER_TIME ) {
855 cout << " WARNING: time of phe out of time range: " << time << endl;
856 return 0. ;
857 }
858 else {
859 return ( Fill( iPix, time - peak_time, CASE_STAR ) ) ;
860 }
861}
862
863Float_t MTrigger::Fill( Int_t iPix, Float_t time, Int_t fall ) {
864 // ============================================================
865 //
866 // Fills the information in the array for the analog signal
867 //
868
869 Float_t PmtAmp = 0 ; // Amplitude of the PMT signal (results from noise)
870
871 if ( iPix < 0 ) {
872 cout << " ERROR: in MTrigger::Fill() " << endl ;
873 cout << " ERROR: Pixel Id < 0 ---> Exit " << endl ;
874 exit (1) ;
875 }
876 else if ( iPix >= CAMERA_PIXELS ) {
877 cout << " ERROR: in MTrigger::Fill() " << endl ;
878 cout << " ERROR: Pixel Id > CAMERA_PIXELS ---> Exit " << endl ;
879 exit (1) ;
880 }
881 else if ( iPix >= pixnum ) {
882 //
883 // We do not have to fill information in the trigger part,
884 // but we must create the height of the pulse going into
885 // the FADC simulation
886 //
887 if (fGainFluctuations)
888 PmtAmp = (histPmt->GetRandom()/histMean) ;
889 else
890 PmtAmp = 1.0;
891
892 // AM April 2004: removed updating of counters nphotshow, nphotnsb,
893 // nphotstar for outer pixels... these arrays (see constructors) are
894 // initialized only with as many elements as trigger pixels!! This
895 // made the camera crash at random in some ocassions.
896 }
897 else {
898 //
899 // we have a trigger pixel and we fill it
900 //
901 Int_t i ;
902
903 //
904 // but at the beginning we must check if this pixel is
905 // hitted the first time
906 //
907
908 if ( used[iPix] == kFALSE ) {
909 used [iPix] = kTRUE ;
910 // baseline[iPix] = 0. ;
911
912 for (i=0; i < TRIGGER_TIME_SLICES; i++ ) {
913 a_sig[iPix][i] = 0. ;
914 d_sig[iPix][i] = 0. ;
915 }
916 }
917
918 //
919 // get the randomized amplitude (unless fGainFluctuations == kFALSE)
920 //
921 if (fGainFluctuations)
922 PmtAmp = (histPmt->GetRandom()/histMean) ;
923 else
924 PmtAmp = 1.0;
925
926 //
927 // select the first slice to fill
928 //
929
930 Int_t ichan = (Int_t) ( time * ((Float_t) TRIG_SLICES_PER_NSEC) ) ;
931
932 //
933 // look over the response signal and put it in the signal line
934 //
935
936 for ( i = 0 ; i<RESPONSE_SLICES_TRIG; i++ ) {
937
938 if ( (ichan+i) >= 0 &&
939 (ichan+i) < TRIGGER_TIME_SLICES ) {
940 a_sig[iPix][ichan+i] += PmtAmp * sing_resp[i] ;
941 }
942 }
943
944 //
945 // we fill the information in the counters of phe's
946 //
947
948 if ( fall == CASE_SHOW )
949 nphotshow[iPix]++ ;
950 else if ( fall == CASE_NSB )
951 nphotnsb[iPix]++ ;
952 else if ( fall == CASE_STAR )
953 nphotstar[iPix]++ ;
954
955 //
956 //
957 return PmtAmp ;
958 }
959 return PmtAmp ;
960}
961
962
963void MTrigger::AddNSB( Int_t iPix, Float_t resp[TRIGGER_TIME_SLICES]){
964 // ================================================================
965 //
966 // Sets the information in the array for the analog signal
967 // from a given array
968 //
969
970
971 if ( iPix < 0 ) {
972 cout << " ERROR: in MTrigger::SetNSB() " << endl ;
973 cout << " ERROR: Pixel Id < 0 ---> Exit " << endl ;
974 exit (1) ;
975 }
976 else if ( iPix >= CAMERA_PIXELS ) {
977 cout << " ERROR: in MTrigger::SetNSB() " << endl ;
978 cout << " ERROR: Pixel Id > CAMERA_PIXELS ---> Exit " << endl ;
979 exit (1) ;
980 }
981 else if ( iPix >= pixnum ) {
982 //
983 // We have not to fill information in the trigger part.
984 //
985 }
986 else {
987 //
988 // we have a trigger pixel and we fill it
989 //
990 Int_t i ;
991
992 //
993 // but at the beginning we must check if this pixel is
994 // hitted the first time
995 //
996
997 if ( used[iPix] == kFALSE ) {
998 used [iPix] = kTRUE ;
999
1000 for (i=0; i < TRIGGER_TIME_SLICES; i++ ) {
1001 a_sig[iPix][i] = 0. ;
1002 d_sig[iPix][i] = 0. ;
1003 }
1004 }
1005
1006 //
1007 // look over the response signal and put it in the signal line
1008 //
1009
1010 for ( i = 0 ; i<TRIGGER_TIME_SLICES; i++ ) {
1011
1012 a_sig[iPix][i] += resp[i];
1013 }
1014
1015 }
1016}
1017
1018void MTrigger::SetElecNoise(Float_t factor){
1019
1020 UInt_t i;
1021 Float_t rausch ;
1022
1023 rausch = ampl_resp * factor ;
1024
1025 cout<<"MTrigger::SetElecNoise ... generating database for electronic noise."
1026 <<endl;
1027
1028 for (i=0;i<TRIGGER_TIME_SLICES*1001;i++){
1029 noise[i]=GenElec->Gaus(0., rausch );
1030 }
1031
1032 cout<<"MTrigger::SetElecNoise ... done"<<endl;
1033
1034}
1035
1036void MTrigger::ElecNoise(Float_t factor) {
1037 // ============================================================
1038 //
1039 // adds the noise due to optronic and electronic
1040 // to the signal
1041 //
1042
1043 UInt_t startslice;
1044
1045 for ( Int_t i=0 ; i < pixnum; i++ ) {
1046 //
1047 // but at the beginning we must check if this pixel is
1048 // hitted the first time
1049 //
1050 startslice=GenElec->Integer(TRIGGER_TIME_SLICES*1000);
1051
1052 if ( used[i] == kFALSE ) {
1053 used [i] = kTRUE ;
1054
1055 memcpy( (Float_t*)a_sig[i],
1056 (Float_t*)&noise[startslice],
1057 TRIGGER_TIME_SLICES*sizeof(Float_t));
1058 memset( (Float_t*)d_sig[i],
1059 0,
1060 TRIGGER_TIME_SLICES*sizeof(Float_t));
1061
1062 }
1063 //
1064 // Then the noise is introduced for each time slice
1065 //
1066 else
1067 for ( Int_t ii=0 ; ii<TRIGGER_TIME_SLICES; ii++ ) {
1068
1069 a_sig [i][ii] += noise[startslice+ii] ;
1070
1071 }
1072 }
1073}
1074
1075void MTrigger::SetFwhm(Float_t fwhm){
1076 //===========================================================
1077 //
1078 // It sets the fwhm for the single phe signal and
1079 // updates the sing_resp for it
1080
1081 Float_t sigma ;
1082 Float_t x, x0 ;
1083 Int_t i;
1084
1085 fwhm_resp = fwhm;
1086
1087 sigma = fwhm_resp / 2.35 ;
1088 x0 = 3*sigma ;
1089
1090 for (i=0; i< RESPONSE_SLICES_TRIG ; i++ ) {
1091
1092 x = i * (1./((Float_t)TRIG_SLICES_PER_NSEC))
1093 + (1./( 2 * (Float_t)TRIG_SLICES_PER_NSEC )) ;
1094
1095 sing_resp[i] = (Float_t)
1096 ampl_resp * expf(-0.5 * (x-x0)*(x-x0) / (sigma*sigma) ) ;
1097
1098 }
1099
1100
1101}
1102
1103void MTrigger::SetMultiplicity(Int_t multi){
1104 //=============================================================
1105 //
1106 // It sets the private member trigger_multi
1107
1108 trigger_multi=multi;
1109}
1110
1111void MTrigger::SetTopology(Int_t topo){
1112 //=============================================================
1113 //
1114 // It sets the private member trigger_geometry
1115
1116 trigger_geometry=topo;
1117}
1118
1119void MTrigger::SetThreshold(Float_t thres[]){
1120 //=============================================================
1121 //
1122 // It sets the private member chan_thres[pixnum]
1123
1124 Int_t i;
1125
1126 for(i=0;i<pixnum;i++){
1127 chan_thres[i]=thres[i];
1128 }
1129}
1130
1131
1132void MTrigger::CheckThreshold(float *thres, int cells){
1133 //=============================================================
1134 //
1135 // Set Right Discriminator threshold, taking into account trigger pixels
1136
1137 FILE *unit;
1138
1139 float thres_aux[CAMERA_PIXELS];
1140 int id;
1141
1142 for (int i=0;i<CAMERA_PIXELS;i++){
1143 if(i<pixnum){
1144 thres_aux[i]=999999.99;
1145 thres[i]=thres[i];
1146 }
1147 else{
1148 thres_aux[i]=-10.0;
1149 thres[i]=-10.0;
1150 }
1151 }
1152
1153 if (cells==1){
1154 if((unit =fopen("../include-MTrigger/TABLE_PIXELS_IN_CELLS", "r" )) == 0 ){
1155 cout << "ERROR: not able to read ../include-MTrigger/TABLE_PIXELS_IN_CELLS"
1156 << endl ;
1157 exit(123) ;
1158 }
1159 else {
1160 while ( feof(unit) == 0 ) {
1161 for ( Int_t k=0; k<TRIGGER_CELLS; k++ ) {
1162 fscanf ( unit, "%d ", &id ) ;
1163 if ((id-1)<pixnum)
1164 thres_aux[id-1]=thres[id-1];
1165 }
1166 }
1167 }
1168 fclose (unit) ;
1169
1170 for (int i=0;i<CAMERA_PIXELS;i++){
1171 thres[i]=thres_aux[i];
1172 }
1173 }
1174
1175}
1176
1177void MTrigger::ReadThreshold(char name[]){
1178 //=============================================================
1179 //
1180 // It reads values for threshold of each pixel from file name
1181
1182 FILE *unit;
1183 Int_t i=0;
1184
1185 if ((unit=fopen(name, "r"))==0){
1186 cout<<"WARNING: not able to read ..."<<name<<endl;
1187 }
1188 else {
1189 while (i<pixnum){
1190 fscanf(unit, "%f",&chan_thres[i++]);
1191 }
1192 fclose (unit);
1193 }
1194
1195}
1196
1197void MTrigger::GetResponse(Float_t *resp) {
1198 // ============================================================
1199 //
1200 // puts the standard response function into the array resp
1201
1202 for ( Int_t i=0; i< RESPONSE_SLICES_TRIG; i++ ) {
1203
1204 resp[i] = sing_resp[i] ;
1205 }
1206
1207}
1208
1209void MTrigger::GetMapDiskriminator(Byte_t *map){
1210 //=============================================================
1211 //
1212 // Gives a map of the fired pixels (Bool_t dknt [pixnum])
1213 // in an array of Byte_t (each byte has the information for 8 pixels)
1214 //
1215
1216 Int_t i,ii;
1217
1218 for(i=0;i<pixnum/8+1;i++){
1219 map[i]=0;
1220 }
1221
1222 for(i=0;i<pixnum;i++){
1223 ii=(Int_t)i/8;
1224 if (dknt[i]==kTRUE){
1225 map[ii]=map[ii]+(Int_t)pow((double)2, (double)i-ii*8);
1226 }
1227 }
1228}
1229
1230
1231void MTrigger::Diskriminate() {
1232 // ============================================================
1233 //
1234 // Diskriminates the analog signal
1235 //
1236 // one very important part is the calucaltion of the baseline
1237 // shift. Because of the AC coupling of the PMT, only the
1238 // fluctuations are interesting. If there are a lot of phe,
1239 // a so-called shift of the baseline occurs.
1240 //
1241
1242 Int_t iM = 0 ;
1243 Int_t i, ii ;
1244
1245
1246 Int_t jmax = (Int_t) (gate_leng * TRIG_SLICES_PER_NSEC ) ;
1247
1248 //
1249 // first of all determine the integral of all signals to get
1250 // the baseline shift.
1251 //
1252
1253 for ( i=0 ; i < pixnum ; i++ ) {
1254 if ( used[i] == kTRUE ) {
1255 baseline[i] = 0. ;
1256
1257 for ( ii = 0 ; ii < TRIGGER_TIME_SLICES ; ii++ ) {
1258 baseline[i] += a_sig[i][ii] ;
1259 }
1260
1261 baseline[i] = baseline[i] / ( (Float_t ) TRIGGER_TIME_SLICES) ;
1262
1263 //
1264 // now correct the baseline shift in the analog signal!!
1265 //
1266 for ( ii = 0 ; ii < TRIGGER_TIME_SLICES ; ii++ ) {
1267 a_sig[i][ii] = a_sig[i][ii] - baseline[i] ;
1268 }
1269 }
1270 }
1271
1272 //
1273 // now the diskrimination is coming
1274 //
1275 // take only that pixel which are used
1276 //
1277
1278 for ( i=0 ; i < pixnum; i++ ) {
1279 if ( used [i] == kTRUE ) {
1280
1281 for ( ii=1 ; ii<TRIGGER_TIME_SLICES; ii++ ) {
1282 //
1283 // first check if the signal is crossing the CHANNEL_THRESHOLD
1284 // form low to big signals
1285 //
1286
1287 if ( a_sig[i][ii-1] < chan_thres[i] &&
1288 a_sig[i][ii] >= chan_thres[i] ) {
1289 {
1290 if ( dknt[i] == kFALSE ) {
1291 dknt [i] = kTRUE ;
1292 iM++ ;
1293 }
1294 // cout << " disk " << ii ;
1295 //
1296 // put the standard diskriminator signal in
1297 // the diskriminated signal
1298 //
1299 for ( Int_t j=0 ; j < jmax ; j++ ) {
1300
1301 if ( ii+j < TRIGGER_TIME_SLICES ) {
1302 d_sig [i][ii+j] = 1. ;
1303 }
1304 }
1305 ii = ii + jmax ;
1306 }
1307 }
1308 else d_sig[i][ii]=0.;
1309 }
1310 }
1311 }
1312}
1313
1314
1315void MTrigger::ShowSignal (MMcEvt *McEvt) {
1316 // ============================================================
1317 //
1318 // This method is used to book the histogramm to show the signal in
1319 // a special gui frame (class MGTriggerSignal). After the look onto the
1320 // signals for a better understanding of the things we will expect
1321 // the gui frame and all histogramms will be destroyed.
1322 //
1323
1324 //
1325 // first of all create a list of the histograms to show
1326 //
1327 // take only that one with a entry
1328
1329 TH1F *hist ;
1330 TH1F *dhist ;
1331 Char_t dumm[10];
1332 Char_t name[256];
1333
1334 TObjArray *AList ;
1335 AList = new TObjArray(10) ;
1336
1337 TObjArray *DList ;
1338 DList = new TObjArray(10) ;
1339
1340 // the list of analog signal histograms
1341 // at the beginning we initalise 10 elements
1342 // but this array expand automaticly if neccessay
1343
1344 Int_t ic = 0 ;
1345 for ( Int_t i=0 ; i < pixnum; i++ ) {
1346 if ( used [i] == kTRUE ) {
1347
1348 sprintf (dumm, "A_%d", i ) ;
1349 sprintf (name, "analog %d", i ) ;
1350
1351 hist = new TH1F(dumm, name, TRIGGER_TIME_SLICES, 0., TOTAL_TRIGGER_TIME);
1352 //
1353 // fill the histogram
1354 //
1355
1356 for (Int_t ibin=1; ibin <=TRIGGER_TIME_SLICES; ibin++) {
1357 hist->SetBinContent (ibin, a_sig[i][ibin-1]) ;
1358 }
1359 hist->SetMaximum(8.);
1360 hist->SetMinimum(-8.);
1361 hist->SetStats(kFALSE);
1362
1363 AList->Add(hist) ;
1364
1365 sprintf (dumm, "D_%d", i ) ;
1366 sprintf (name, "digital %d", i ) ;
1367
1368 dhist = new TH1F(dumm, name, TRIGGER_TIME_SLICES, 0., TOTAL_TRIGGER_TIME);
1369 if ( dknt[i] == kTRUE ) {
1370 //
1371 // fill the histogram of digital signal
1372 //
1373 for (Int_t ibin=1; ibin <=TRIGGER_TIME_SLICES; ibin++) {
1374 dhist->SetBinContent (ibin, d_sig[i][ibin-1]) ;
1375 dhist->SetStats(kFALSE);
1376 }
1377 }
1378 dhist->SetMaximum(1.5);
1379
1380 DList->Add(dhist);
1381
1382 ic++ ;
1383
1384 }
1385 }
1386
1387 //
1388 // create the Gui Tool
1389 //
1390 //
1391
1392 new MGTriggerSignal(McEvt,
1393 AList,
1394 DList,
1395 gClient->GetRoot(),
1396 gClient->GetRoot(),
1397 400, 400 ) ;
1398
1399 //
1400 // delete the List of histogramms
1401 //
1402
1403 AList->Delete() ;
1404 DList->Delete() ;
1405
1406 delete AList ;
1407 delete DList ;
1408}
1409
1410
1411Int_t MTrigger::ZeroLevel() {
1412 // ============================================================
1413 //
1414 // This is a level introduced just to speed up the program.
1415 // It makes sense to look for next neighbours only if there
1416 // are at least trigger_multi pixels with a diskriminator
1417 // signal.
1418 //
1419
1420 //
1421 // first count the pixels with a diskriminator signal
1422 //
1423 Int_t iMul = 0 ;
1424 for ( Int_t iP =0 ; iP < pixnum; iP++ ) {
1425 //
1426 //
1427
1428 if ( dknt[iP] == kTRUE ) {
1429 iMul++ ;
1430 }
1431 }
1432
1433
1434 //
1435 // only if there are at least more pixels than requested
1436 // it make sense to look into details
1437 if ( iMul >= trigger_multi ) {
1438 //
1439 // fill the sum signal of all diskriminator signals
1440 //
1441 for ( Int_t iP =0 ; iP < pixnum; iP++ ) {
1442 //
1443 //
1444 if ( dknt[iP] == kTRUE ) {
1445 //
1446 // sum it up
1447 //
1448 for (Int_t iS=0; iS< TRIGGER_TIME_SLICES; iS++ ) {
1449 //
1450 //
1451 sum_d_sig [iS] += d_sig[iP][iS] ;
1452 }
1453 }
1454 }
1455 //
1456 // run over the sum_d_sig and check each time slice
1457 //
1458 Int_t iReturn = 0 ;
1459
1460 for (Int_t iS=0; iS< TRIGGER_TIME_SLICES; iS++ ) {
1461
1462 if ( sum_d_sig[iS] >= trigger_multi ) {
1463 iReturn++ ;
1464 nZero++;
1465 SlicesZero[iS] = kTRUE ;
1466
1467 }
1468 else SlicesZero[iS] = kFALSE;
1469 }
1470 return ( iReturn ) ;
1471 }
1472 else {
1473 return 0 ;
1474 }
1475}
1476
1477Int_t MTrigger::FirstLevel() {
1478 //=================================================
1479 //
1480 // This is a level trigger which can look for several
1481 // multiplicities (trigger_multi)
1482 // and topologies (trigger_geometry)
1483 //
1484
1485 Int_t iReturn = 0 ; // Return value for this function
1486
1487 // Definition of needed variables
1488 Bool_t Muster[pixnum] ;
1489 Bool_t Neighb[pixnum] ;
1490 Int_t iMulti = 0 ;
1491
1492 // We put several wrong topologies which we already know that they
1493 // are not possible. It can save time.
1494
1495 if (trigger_geometry==0 && trigger_multi>7) {
1496 cout <<"You are looking for a topology that needs more than six neighbours of the same pixel"<<endl;
1497 cout <<" Topology "<<trigger_geometry<<" Multiplicity "<<trigger_multi<<endl;;
1498 return (kFALSE);
1499 }
1500
1501 if (trigger_geometry==2 && trigger_multi<3) {
1502 cout<<"Closed pack geometry with multiplicity "<<trigger_multi<<" does not make sense, I'll check simple neihgbour condition"<<endl;
1503 trigger_geometry=1;
1504 }
1505 if (trigger_geometry>2) {
1506 cout << "This trigger topology is not implemented"<<endl;
1507 return (kFALSE);
1508 }
1509
1510 //
1511 // loop over all ZeroLevel Trigger
1512 //
1513 // it is only neccessary to look after a ZeroLevel Trigger for
1514 // a FirstLevel (NextNeighbour) trigger.
1515 //
1516
1517 if (nZero) {
1518
1519 //
1520 // Then run over all slices
1521 //
1522
1523 for ( Int_t iSli = 0;
1524 iSli < TRIGGER_TIME_SLICES; iSli++ ) {
1525
1526 // Check if this time slice has more fired pixels than trigger_multi
1527
1528 if (SlicesZero[iSli]){
1529 //
1530 // Loop over trigger cells. It is topology analysis,
1531 // therefore it is kept here after multiplicity and
1532 // threshold checks.
1533 //
1534
1535 for(Int_t iCell=0; iCell<TRIGGER_CELLS; iCell++){
1536 //
1537 // then look in all pixel of that cell if the
1538 // diskriminated signal is 1
1539 //
1540 for ( Int_t iPix = 0 ; iPix < pixnum; iPix++ ) {
1541 Muster[iPix] = kFALSE ;
1542 Neighb[iPix] = kFALSE ;
1543 // Select pixels which are used and it the current cell
1544 if ( used [iPix] == kTRUE && TC[iCell][iPix]==kTRUE) {
1545 //
1546 // now check the diskriminated signal
1547 //
1548 if ( d_sig [iPix][iSli] > 0. ) {
1549 Muster[iPix] = kTRUE ;
1550 }
1551 }
1552 } // end of loop over the pixels
1553
1554 //
1555 // Here we check which of the "muster" pixels will be fired for
1556 // the minimum required overlaping time
1557 //
1558
1559 OverlapingTime(Muster, &Muster[0],iSli);
1560
1561 //
1562 // here we have to look for the topologies
1563 //
1564
1565 switch(trigger_geometry){
1566 case 0:{
1567
1568 // It looks for a pixel above threshold which has
1569 // trigger_multi-1 neighbour pixels above threshold
1570
1571 Bool_t Dummy[pixnum] ;
1572
1573 // Loop over all pixels
1574 for (int j=0;j<pixnum;j++){
1575
1576 for (int k=0; k<pixnum; k++){
1577 Neighb[k]=kFALSE;
1578
1579 Dummy[k] = Muster[k] ;
1580 }
1581 if(Muster[j]){
1582 // If pixel is fired, it checks how many fired neighbours it has
1583 for (iMulti=1;iMulti<trigger_multi; iMulti++) {
1584 Neighb[j] = kTRUE ;
1585 Dummy[j] = kTRUE ;
1586 if (!PassNextNeighbour(Dummy, &Neighb[0])){
1587 break;
1588 }
1589 for (int k=0; k<pixnum; k++){
1590 if (Neighb[k]){
1591 Dummy[k]=kFALSE;
1592 Neighb[k]=kFALSE;
1593 }
1594 }
1595 }
1596 if (iMulti==trigger_multi ) {
1597 //
1598 // A NN-Trigger is detected at time Slice
1599 //
1600 PixelsFirst[nFirst] = j; // We save pixel that triggers
1601 SlicesFirst[nFirst++] = iSli ; // We save time when it triggers
1602 iReturn++ ;
1603 iSli+=(LEVEL1_DEAD_TIME*TRIG_SLICES_PER_NSEC); // We skip the following 50 ns (dead time)
1604 iCell=TRIGGER_CELLS; // We skip the remaining trigger cells
1605 break ;
1606 }
1607 }
1608 }
1609 break;
1610 };
1611
1612 case 1:{
1613
1614 // It looks for trigger_multi neighbour pixels above the
1615 // threshold.
1616
1617 for (int j=0;j<pixnum;j++){
1618 if(Muster[j]){
1619 // It checks if you can find
1620 // trigger_multi fired neighbour pixels
1621 Neighb[j] = kTRUE ;
1622 for (iMulti=1;iMulti<trigger_multi; iMulti++) {
1623 if (!PassNextNeighbour(Muster, &Neighb[0]))
1624 break;
1625 }
1626 if (iMulti==trigger_multi ) {
1627 //
1628 // A NN-Trigger is detected at time Slice
1629 //
1630 PixelsFirst[nFirst] = j; // We save pixel that triggers
1631 SlicesFirst[nFirst++] = iSli ; // We save when it triggers
1632 iReturn++ ;
1633 iSli+=(LEVEL1_DEAD_TIME*TRIG_SLICES_PER_NSEC); // We skip the following 50 ns (dead time)
1634 iCell=TRIGGER_CELLS; // We skip the remaining trigger cells
1635 break ;
1636 }
1637 else {
1638 // We put Neighb to kFALSE to check an other pixel
1639 for (int k=0; k<pixnum; k++){
1640 if (Neighb[k]){
1641 Neighb[k]=kFALSE;
1642 }
1643 }
1644 }
1645 }
1646 }
1647 break;
1648 };
1649 case 2:{
1650
1651 // It looks for trigger_multi closed pack neighbours
1652 // above threshold
1653 // Closed pack means that you can take out any pixel
1654 // and you will still get a trigger for trigger_multi -1
1655 // The algorithm is not perfect, there still somes cases
1656 // that are not really well treated
1657
1658 Int_t closed_pack = 1;
1659
1660 for (int j=0;j<pixnum;j++){
1661 if(Muster[j]){
1662 // It checks if there are trigger_multi
1663 // neighbours above threshold
1664
1665 Neighb[j] = kTRUE ;
1666 iMulti=1;
1667
1668 //while(PassNextNeighbour(Muster, &Neighb[0])) iMulti++;
1669 for (iMulti=1;iMulti<trigger_multi;iMulti++){
1670 if (!PassNextNeighbour(Muster, &Neighb[0]))
1671 break;
1672 }
1673
1674 if (iMulti==trigger_multi ) {
1675 //
1676 // A NN-Trigger is detected at time Slice
1677 //
1678
1679 // Check if there is closed pack topology
1680
1681 Bool_t Aux1[pixnum];
1682 Bool_t Aux2[pixnum];
1683 for (int jj=0;jj<pixnum;jj++)
1684 Aux2[jj]=kFALSE;
1685
1686 for (int i=0;i<pixnum;i++){
1687 if (Neighb[i]) {
1688 // Loop over pixels that achive neighbouring condition
1689
1690 for (int jj=0;jj<pixnum;jj++) {
1691
1692 Aux1[jj] = Neighb[jj] ; // huschel
1693 Aux2[jj]=kFALSE;
1694 }
1695
1696 // It checks if taking out any of the pixels we lose
1697 // neighbouring condition for trigger_multi -1
1698
1699 Aux1[i]=kFALSE;
1700 closed_pack=0;
1701 for (int jj=0;jj<pixnum;jj++) {
1702 if (Aux1[jj]==kTRUE){
1703 Aux2[jj]=kTRUE;
1704 for (iMulti=1;iMulti<(trigger_multi-1);iMulti++){
1705 if (!PassNextNeighbour(Aux1, &Aux2[0]))
1706 break;
1707 }
1708 if (iMulti==(trigger_multi-1)){
1709 // We found a NN trigger for trigger_multi -1
1710 // taking out pixel jj
1711 closed_pack=1;
1712 break;
1713 }
1714 Aux2[jj]=kFALSE;
1715 }
1716 }
1717 if (!closed_pack) break;
1718 // For some pixell we did not found NN condition
1719 // for trigger_multi -1
1720 }
1721 }
1722 if (closed_pack){
1723 PixelsFirst[nFirst] = j; // We save pixel that triggers
1724 SlicesFirst[nFirst++] = iSli ; // We save time when it triggers
1725 iReturn++ ;
1726 iSli+=(LEVEL1_DEAD_TIME*TRIG_SLICES_PER_NSEC); // We skip the following 50 ns (dead time)
1727 iCell=TRIGGER_CELLS; // We skip the remaining trigger cells
1728 break ;
1729 }
1730 else {
1731 for (int k=0; k<pixnum; k++){
1732 if (Neighb[k]){
1733 Neighb[k]=kFALSE;
1734 }
1735 }
1736 }
1737 } // end if trigger multiplicity achived
1738 else{
1739 for (int k=0; k<pixnum; k++)
1740 Neighb[k]=kFALSE;
1741 }
1742 } // end if pixel fired
1743 } // end loop trigger pixels
1744 break;
1745 }; // end case 2:
1746 default:{
1747 cout << "This topology is not implemented yet"<<endl;
1748 break;
1749 }
1750 }
1751 } //end loop over trigger cells.
1752 }
1753 } // end of loop over the slices
1754 } // end of conditional for a trigger Zero
1755
1756 //
1757 // return the Number of FirstLevel Triggers
1758 //
1759 return iReturn ;
1760}
1761
1762
1763Bool_t MTrigger::PassNextNeighbour ( Bool_t m[], Bool_t *n) {
1764 //
1765 // This function is looking for a next neighbour of pixels in n[]
1766 // above triggers using a NNlookup table.
1767 // This table is builded by the default constructor
1768 //
1769
1770 //
1771 // loop over all trigger pixels
1772 //
1773
1774 Bool_t return_val = kFALSE;
1775
1776 for ( Int_t i=0; i<pixnum; i++) {
1777 //
1778 // check if this pixel has a diskrminator signal
1779 // (this is inside n[] )
1780 //
1781
1782 if ( n[i] && !return_val) {
1783
1784 //
1785 // look in the next neighbours from the lookuptable
1786 //
1787
1788 for ( Int_t kk=0; kk<6; kk++ ) {
1789 //
1790 // if the nextneighbour is outside the triggerarea do nothing
1791 //
1792 if (!return_val){
1793 if (NN[kk][i] >= pixnum ) {
1794
1795 }
1796 // the nextneighbour is not inside the pixnum
1797 else {
1798 //
1799 // look if the boolean of nn pixels is true
1800 //
1801
1802 if ( m[ NN[kk][i] ] && !n[NN[kk][i]] ) {
1803 n[NN[kk][i]]=kTRUE ;
1804 return_val =kTRUE;
1805 }
1806 }
1807 }
1808 else break;
1809 }
1810 }
1811 }
1812 return(return_val);
1813}
1814
1815Float_t MTrigger::GetFirstLevelTime( Int_t il ){
1816
1817 //=============================================================
1818 //
1819 // It gives the time for the il trigger at first level
1820
1821 return((Float_t) ((Float_t) SlicesFirst[il]/((Float_t) TRIG_SLICES_PER_NSEC)));
1822}
1823
1824Int_t MTrigger::GetFirstLevelPixel( Int_t il ){
1825
1826 //=============================================================
1827 //
1828 // It gives the pixel that triggers for the il trigger at first level
1829 return(PixelsFirst[il]);
1830}
1831
1832void MTrigger::OverlapingTime ( Bool_t m[], Bool_t *n, Int_t ifSli){
1833
1834 //============================================================
1835 //
1836 // It returns in n the pixels of m that are fired during the
1837 // required overlaping time for trigger after ifSli
1838
1839 int i,j;
1840 int iNumSli;
1841
1842 // Translation from ns to slices
1843 iNumSli=(int) (overlaping_time*TRIG_SLICES_PER_NSEC);
1844 if (iNumSli<1) iNumSli=1;
1845
1846 // Put pixels that fulfill the requirement in n
1847 for (i=0;i<pixnum;i++){
1848 if (m[i]==kTRUE){
1849 for(j=ifSli;j<ifSli+iNumSli;j++){
1850 if(!d_sig[i][j]){
1851 n[i]=kFALSE;
1852 break;
1853 }
1854 }
1855 }
1856 }
1857
1858}
1859
1860//=================================================================
1861
1862void MTrigger::InitGainFluctuations()
1863{
1864 //
1865 // The amplitude of one single photo electron is not a constant.
1866 // There exists a measured distribution from Razmik. This distribution
1867 // is used to simulate the noise of the amplitude.
1868 // For this a histogramm (histPmt) is created and filled with the
1869 // values.
1870 //
1871
1872 histPmt = new TH1F ("histPmt","Noise of PMT", 40, 0., 40.) ;
1873
1874 Stat_t ValRazmik[40] = { 2.14, 2.06, 2.05, 2.05, 2.06, 2.07, 2.08, 2.15,
1875 2.27, 2.40, 2.48, 2.55, 2.50, 2.35, 2.20, 2.10,
1876 1.90, 1.65, 1.40, 1.25, 1.00, 0.80, 0.65, 0.50,
1877 0.35, 0.27, 0.20, 0.18, 0.16, 0.14, 0.12, 0.10,
1878 0.08, 0.06, 0.04, 0.02, 0.01, 0.005,0.003, 0.001} ;
1879
1880 for (Int_t i = 0; i < 40; i++)
1881 histPmt->SetBinContent(i+1 ,ValRazmik[i]);
1882
1883 histMean = histPmt->GetMean() ;
1884
1885 fGainFluctuations = kTRUE; // Default value
1886
1887 return;
1888}
Note: See TracBrowser for help on using the repository browser.