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

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