source: JEM-EUSO/esaf_lal/tags/v1_r0/esaf/packages/simulation/detector/electronics/src/PTTTrigger.cc @ 117

Last change on this file since 117 was 117, checked in by moretto, 11 years ago

ESAF version compilable on mac OS

File size: 14.5 KB
Line 
1// $Id$
2// Author: fenu   2007/12/05
3
4/*****************************************************************************
5 * ESAF: Euso Simulation and Analysis Framework                              *
6 *                                                                           *
7 *  Id: BertainaTrigger                                                           *
8 *  Package: <packagename>                                                   *
9 *  Coordinator: <coordinator>                                               *
10 *                                                                           *
11 *****************************************************************************/
12
13//_____________________________________________________________________________
14//
15// PTTTrigger
16//
17// <extensive class description>
18//
19//   Config file parameters
20//   ======================
21//
22//   <parameter name>: <parameter description>
23//   -Valid options: <available options>
24/*
25 Algorithm description:
26 Looks for persistency of signal over time in some EC.
27
28 1) Checks all the pixels above a certain threshold. (PTTTrigger.PixelThr )
29 2) Calculates the excess on a 3X3 box around the pixel with respect to this threshold
30 3) Counts for how many consecutive GTUs the condition counts>PTTTrigger.PixelThr is
31 verified within one EC.
32 4) If the number of consecutive GTUs satisfying this condition is above a threshold
33 (PTTTrigger.Persistency) integration is started.
34 5) The integration has to be performed on the excesses from GTU(i-PTTTrigger.Persistency)
35 to GTU(i) on the entire EC for each pixel which was above PTTTrigger.PixelThr.
36 6) If the integration for at least one pixel is above a certain threshold
37 (PTTTrigger.Integration) a trigger
38 is issued.
39
40
41
42 Methods description
43 *  PTTTrigger: initializer
44 *  ~PTTTrigger destructor
45 *  Simulate the entire algorithm logic  is here
46 *  Transform performs all the necessary transformations
47 for properly building and translating the box
48 *  Debug table make a table for visualising the box
49 translation and the position of all the boxes
50
51 */
52// If you change something on the layout structure please make a test to
53// check if the box  is properly translated
54// Especially if you change the definition of the geometry
55// within one PMT this will for sure affect the PTT algorithm
56// put PTTTrigger.fDebug=1 and run just one event.
57// Look to Debug_Table1 and Debug_Table2 (end of this file) for understanding the meaning
58// of the debug table. Names used there are self explaining
59//
60// FOR TRIGGER DEFINITION: if Side cut is used put a limit of +-170km on the y axis.
61//      Otherwise in the efficiency plot a plateau of 80% will be reached nut never 100%
62//      This is because part of the FOV is not covered for sidecut.
63//      in Y direction FOV +-22deg while in X is +-30
64//
65//To read the .root file use
66//
67//   TFile  *f = new TFile("nameOffile.root");
68//   TTree *etree = (TTree*)f->Get("etree");
69//   EEvent *evs = new EEvent();
70//   evs->SetBranches(etree);
71//   etree->GetEntry(ii);//ii number of event to read
72// after that
73// evs->GetPTTTrigger()->GetPTTTrackNum()
74// will give the number of PDMs having at least one triggering track
75// The event will have triggered if evs->GetPTTTrigger()->GetPTTTrackNum()>0
76// For other useful informations on the trigger
77// look in the file packages/common/root/include/ELTTTriggerSegment.hh
78//
79
80//
81#include "PTTTrigger.hh"
82#include "MacroCellData.hh"
83#include "MacroCell.hh"
84#include "ERunParameters.hh"
85#include "NumbersFileParser.hh"
86#include "EPTTTriggerDataAdder.hh"
87#include "PTTTriggerSegment.hh"
88#include "EEvent.hh"
89#include "EHeader.hh"
90#include <cmath>
91
92ClassImp(PTTTrigger)
93//_____________________________________________________________________________
94PTTTrigger::PTTTrigger() :
95    TriggerEngine(string("PTTTrigger"), kPTTTrigger) {
96  //
97  // Constructor
98  //
99  fPixelThr = (Int_t) Conf()->GetNum("PTTTrigger.PixelThr");
100  fPersistency = (Int_t) Conf()->GetNum("PTTTrigger.Persistency");
101  fIntegration = (Int_t) Conf()->GetNum("PTTTrigger.Integration");
102  Flag_debug_table = (Int_t) Conf()->GetNum("PTTTrigger.fDebug");
103}
104
105//_____________________________________________________________________________
106PTTTrigger::~PTTTrigger() {
107  //
108  // Destructor
109  //
110
111  cout << "Destroy " << endl;
112  Clear();
113}
114
115//______________________________________________________________________________
116void PTTTrigger::Simulate(MacroCellData * pData) {
117  // Simulate PTT algorithm
118  Msg(EsafMsg::Info) << "Simulating PTT algorithm (Catalano) " << MsgDispatch;
119  map<Int_t, map<Int_t, ChipGtuData *>*>&m1 = pData->fChipData2;
120  map<Int_t, map<Int_t, ChipGtuData *>*>::const_iterator it1;
121  Int_t matsize = pData->Cell()->GetEC(0)->FrontEnd()->Channels();
122  side_EC = sqrt(matsize);
123  side_PMT = side_EC / 2;
124  Int_t maximum_int = 0;
125
126  Int_t flag2 = 0;
127  // iterate on CHIPS in this macrocell
128  for (it1 = m1.begin(); it1 != m1.end(); it1++) {
129
130    map<Int_t, ChipGtuData *>&m2 = *(it1->second);
131    map<Int_t, ChipGtuData *>::const_iterator it2;
132    // iterate on GTU for the current chip
133    Int_t pers = 0;
134    Int_t integr[side_EC][side_EC][fPersistency];
135    Int_t sum_int[side_EC][side_EC];
136
137    //initialize counters
138    for (Int_t i = 0; i < side_EC; i++)
139      for (Int_t ii = 0; ii < side_EC; ii++) {
140        sum_int[i][ii] = 0;
141        for (Int_t iii = 0; iii < fPersistency; iii++)
142          integr[i][ii][iii] = 0;
143
144      }
145    Int_t gtu_count = 0;
146    // iterate on GTUs in this chip
147    for (it2 = m2.begin(); it2 != m2.end(); it2++) {
148      if (gtu_count > 0 && gtu_count < m2.end()->first) {
149        Int_t integer_int, rest;
150        double integer;
151        rest = modf(it2->first / fPersistency, &integer);
152        integer_int = integer;
153        Int_t gtu = it2->first - integer_int * fPersistency;
154        for (Int_t i = 0; i < side_EC; i++)
155          for (Int_t ii = 0; ii < side_EC; ii++)
156            integr[i][ii][gtu] = 0;
157
158        //run over all the pixels in one GTU
159        Int_t flag = 0;
160        for (Int_t iteratore = 0; iteratore < matsize; iteratore++) {
161          Int_t nn_PMT = it2->second->FrontEnd()->PmtChannel(iteratore)
162              / (side_PMT);
163          Int_t mm_PMT = it2->second->FrontEnd()->PmtChannel(iteratore)
164              % (side_PMT);
165
166          if (Flag_debug_table)
167            Debug_Table1(
168                it2->second->FrontEnd()->Pmt(iteratore)->Geometry()->Position(
169                    nn_PMT, mm_PMT).X(),
170                it2->second->FrontEnd()->Pmt(iteratore)->Geometry()->Position(
171                    nn_PMT, mm_PMT).Y(),
172                it2->second->FrontEnd()->Pmt(iteratore)->Geometry()->Position(
173                    nn_PMT, mm_PMT).Z(), it2->first, it1->first,
174                it2->second->FrontEnd()->Row(iteratore),
175                it2->second->FrontEnd()->Column(iteratore));
176
177          Int_t nn_EC = it2->second->FrontEnd()->Row(iteratore);
178          Int_t mm_EC = it2->second->FrontEnd()->Column(iteratore);
179
180          Int_t center_row, center_col, center_row_aux, center_col_aux;
181          center_row = nn_EC;
182          center_col = mm_EC;
183
184          //transform
185          Transform(&center_row, &center_col);
186          center_row_aux = center_row;
187          center_col_aux = center_col;
188
189          //inverts
190          Transform(&center_row_aux, &center_col_aux);
191
192          Int_t new_row, new_col;
193          // since a 3X3 box is build you need to stay within 1 pixel border from the edge of the EC
194          if (center_row > 0 && center_row < side_EC - 1 && center_col > 0
195              && center_col < side_EC - 1) {
196            //build the 3X3 box
197            for (Int_t x = -1; x < 2; x++) {
198              for (Int_t y = -1; y < 2; y++) {
199                new_row = center_row + x;
200                new_col = center_col + y;
201
202                //Inverts transformation functions
203                Transform(&new_row, &new_col);
204                Int_t channel = it2->second->FrontEnd()->ChanRowCol(new_row,
205                    new_col);
206
207                nn_PMT = it2->second->FrontEnd()->PmtChannel(channel)
208                    / (side_PMT);
209                mm_PMT = it2->second->FrontEnd()->PmtChannel(channel)
210                    % (side_PMT);
211
212                if (Flag_debug_table)
213                  Debug_Table2(
214                      it2->second->FrontEnd()->Pmt(channel)->Geometry()->Position(
215                          nn_PMT, mm_PMT).X(),
216                      it2->second->FrontEnd()->Pmt(channel)->Geometry()->Position(
217                          nn_PMT, mm_PMT).Y(),
218                      it2->second->FrontEnd()->Pmt(channel)->Geometry()->Position(
219                          nn_PMT, mm_PMT).Z(), x, y, it2->first, it1->first,
220                      channel);
221
222                //excess is calculated
223                Int_t eccesso = it2->second->GetCounter(channel) - fPixelThr
224                    + 1;
225
226                //if above threshold
227                if (it2->second->GetCounter(channel) >= fPixelThr) {
228                  // calculate the excess of the 3X3 box
229                  //assigning it to the central pixel
230                  integr[center_row_aux][center_col_aux][gtu] += eccesso;
231                  // in this GTU something is above threshold
232                  flag = 1;
233                }
234
235              } //for(Int_t y=-1;y<2;y++){
236            } //     for(Int_t x=-1;x<2;x++)
237
238          } // if(center_row>0 && center_row<11 && center_col>0 && center_col<11)
239        } // for(Int_t iteratore=0;iteratore<matsize ;iteratore++) {
240
241        //if within this EC and GTU something is above threshold....
242        if (flag == 1) {
243          // update the persistency count
244          pers += 1;
245          // if now the persistency  is above a certain value...
246          if (pers >= fPersistency) {
247            for (Int_t i = 0; i < side_EC; i++) {
248              for (Int_t ii = 0; ii < side_EC; ii++) {
249                // calculates the sum of the excesses for each pixel in the
250                // fPersistency GTUs previous to the one under analysis
251                sum_int[i][ii] = 0;
252                for (Int_t iii = 0; iii < fPersistency; iii++)
253                  sum_int[i][ii] += integr[i][ii][iii];
254
255                // is there one pixel with integration above the threshold?
256                if (sum_int[i][ii] >= fIntegration) {
257                  flag2 = 1;
258                  //TRIGGER!!!!
259                  if (maximum_int < sum_int[i][ii]) {
260                    // is the trigger the maximum?
261                    maximum_int = sum_int[i][ii];
262                    row_up = i;
263                    col_up = ii;
264                    GTU_up = it2->first;
265                    chip_ID_up = it1->first;
266                    Int_t FEE_ch = it2->second->FrontEnd()->ChanRowCol(i, ii);
267                    PdmID = it2->second->FrontEnd()->Pmt(FEE_ch)->Cell()->Id();
268                  }
269                }
270              }
271            }
272          }
273        } else { //if in this GTU no pixel above threshold...
274          //reinitialize all the counters
275          pers = 0;
276          for (Int_t i = 0; i < side_EC; i++)
277            for (Int_t ii = 0; ii < side_EC; ii++)
278              for (Int_t iii = 0; iii < fPersistency; iii++)
279                integr[i][ii][iii] = 0;
280        }
281      }
282      gtu_count++;
283    } //  for(it2=m2.begin(); it2!=m2.end(); it2++) {
284
285  } // for(it1=m1.begin(); it1!=m1.end(); it1++) {
286
287  if (flag2 == 1) {
288    //FILLING the INFORMATION VECTOR
289    fTracks.push_back(PTTTriggerSegment());
290    PTTTriggerSegment & trkseg = fTracks.back();
291
292    trkseg.SetNumEvt(EEvent::GetCurrent()->GetHeader()->GetNum());
293    trkseg.SetMaxCount(maximum_int);
294    trkseg.SetMaxGtu(GTU_up);
295    trkseg.SetMaxRow(row_up);
296    trkseg.SetMaxCol(col_up);
297    trkseg.SetMaxChipID(chip_ID_up);
298    trkseg.SetPdmID(PdmID);
299    SetTrigger();
300
301    pData->SetPTTTrigger(EEvent::GetCurrent()->GetHeader()->GetNum());
302    pData->SetPTTTriggerMaxCounts(maximum_int);
303    pData->SetPTTTriggerMaxGtu(GTU_up);
304    pData->SetPTTTriggerMaxRow(row_up);
305    pData->SetPTTTriggerMaxCol(col_up);
306    pData->SetPTTTriggerChipID(chip_ID_up);
307    pData->SetPTTTriggerPdmID(PdmID);
308
309    FillEEvent(EEvent::GetCurrent()->GetHeader()->GetNum());
310
311    Msg(EsafMsg::Info) << "TRIGGERING... NO MORE THAN 1 TRIGGER X 1 Macro Cell "
312        << MsgDispatch;
313    Msg(EsafMsg::Info) << "Integrated " << maximum_int << " Chip ID "
314        << chip_ID_up << " GTU ID " << GTU_up << MsgDispatch;
315
316  } else {
317    pData->SetPTTTrigger(-100);
318    pData->SetPTTTriggerMaxCounts(-100);
319    pData->SetPTTTriggerMaxGtu(-100);
320    pData->SetPTTTriggerMaxRow(-100);
321    pData->SetPTTTriggerMaxCol(-100);
322    pData->SetPTTTriggerChipID(-100);
323    pData->SetPTTTriggerPdmID(-100);
324  }
325
326  Clear();
327
328}
329
330//______________________________________________________________________________
331void PTTTrigger::Clear() {
332/////////////////////////////////////////////////////////////////////////////////////
333//clear/////////////////////////////////////////////////////////////////////////////////////
334/////////////////////////////////////////////////////////////////////////////////////
335
336}
337
338//______________________________________________________________________________
339void PTTTrigger::Transform(Int_t * row, Int_t * col) {
340  //trasformation
341
342  //transformation
343  //first part
344  if (*row < side_PMT && *col < side_PMT) {
345    Int_t aux = *col;
346    *col = *row;
347    *row = aux;
348  }
349
350  if (*row > side_PMT - 1 && *col < side_PMT) {
351    Int_t aux = *col;
352    *col = *row - side_PMT;
353    *row = aux + side_PMT;
354  }
355
356  if (*row < side_PMT && *col > side_PMT - 1) {
357    Int_t aux = *col;
358    *col = *row + side_PMT;
359    *row = aux - side_PMT;
360  }
361
362  if (*row > side_PMT - 1 && *col > side_PMT - 1) {
363    Int_t aux = *col;
364    *col = *row;
365    *row = aux;
366  }
367
368  //second part
369  if (*row < side_PMT)
370    *row += side_PMT;
371  else
372    *row -= side_PMT;
373
374}
375
376void PTTTrigger::FillEEvent(Int_t Triggering) {
377  vector<PTTTriggerSegment>::iterator it;
378
379  for (it = fTracks.begin(); it != fTracks.end(); it++) {
380
381    if (EEvent::GetCurrent()) {
382
383      EPTTTriggerDataAdder a(&(*it), Triggering, 1);
384      EEvent::GetCurrent()->Fill(a);
385    }
386
387  }
388  fTracks.clear();
389
390}
391
392//PIXEL POSITION
393void PTTTrigger::Debug_Table1(Float_t x, Float_t y, Float_t z,
394    Int_t analyzed_gtu, Int_t chip_id, Int_t row, Int_t col) {
395  FILE *ff = fopen("Debug_Table1.txt", "a");
396
397  fprintf(ff, "%f  %f  %f  %d  %d  %d  %d\n", x, y, z, analyzed_gtu, chip_id,
398      row, col);
399
400  fclose(ff);
401}
402
403//INTEGRATION BOX POSITION
404void PTTTrigger::Debug_Table2(Float_t x, Float_t y, Float_t z, Int_t i,
405    Int_t dir, Int_t analyzed_gtu, Int_t chip_id, Int_t box) {
406  FILE *ff = fopen("Debug_Table2.txt", "a");
407
408  fprintf(ff, "%f  %f  %f  %d  %d  %d  %d  %d\n", x, y, z, i, dir, analyzed_gtu,
409      chip_id, box);
410
411  fclose(ff);
412}
Note: See TracBrowser for help on using the repository browser.