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

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