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

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