source: JEM-EUSO/esaf_cc_at_lal/packages/simulation/detector/electronics/src/PTTTrigger.cc @ 114

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

actual version of ESAF at CCin2p3

File size: 14.4 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
91ClassImp(PTTTrigger)
92//_____________________________________________________________________________
93PTTTrigger::PTTTrigger() :
94    TriggerEngine(string("PTTTrigger"), kPTTTrigger) {
95  //
96  // Constructor
97  //
98  fPixelThr = (Int_t) Conf()->GetNum("PTTTrigger.PixelThr");
99  fPersistency = (Int_t) Conf()->GetNum("PTTTrigger.Persistency");
100  fIntegration = (Int_t) Conf()->GetNum("PTTTrigger.Integration");
101  Flag_debug_table = (Int_t) Conf()->GetNum("PTTTrigger.fDebug");
102}
103
104//_____________________________________________________________________________
105PTTTrigger::~PTTTrigger() {
106  //
107  // Destructor
108  //
109
110  cout << "Destroy " << endl;
111  Clear();
112}
113
114//______________________________________________________________________________
115void PTTTrigger::Simulate(MacroCellData * pData) {
116  // Simulate PTT algorithm
117  Msg(EsafMsg::Info) << "Simulating PTT algorithm (Catalano) " << MsgDispatch;
118  map<Int_t, map<Int_t, ChipGtuData *>*>&m1 = pData->fChipData2;
119  map<Int_t, map<Int_t, ChipGtuData *>*>::const_iterator it1;
120  Int_t matsize = pData->Cell()->GetEC(0)->FrontEnd()->Channels();
121  side_EC = sqrt(matsize);
122  side_PMT = side_EC / 2;
123  Int_t maximum_int = 0;
124
125  Int_t flag2 = 0;
126  // iterate on CHIPS in this macrocell
127  for (it1 = m1.begin(); it1 != m1.end(); it1++) {
128
129    map<Int_t, ChipGtuData *>&m2 = *(it1->second);
130    map<Int_t, ChipGtuData *>::const_iterator it2;
131    // iterate on GTU for the current chip
132    Int_t pers = 0;
133    Int_t integr[side_EC][side_EC][fPersistency];
134    Int_t sum_int[side_EC][side_EC];
135
136    //initialize counters
137    for (Int_t i = 0; i < side_EC; i++)
138      for (Int_t ii = 0; ii < side_EC; ii++) {
139        sum_int[i][ii] = 0;
140        for (Int_t iii = 0; iii < fPersistency; iii++)
141          integr[i][ii][iii] = 0;
142
143      }
144    Int_t gtu_count = 0;
145    // iterate on GTUs in this chip
146    for (it2 = m2.begin(); it2 != m2.end(); it2++) {
147      if (gtu_count > 0 && gtu_count < m2.end()->first) {
148        Int_t integer_int, rest;
149        double integer;
150        rest = modf(it2->first / fPersistency, &integer);
151        integer_int = integer;
152        Int_t gtu = it2->first - integer_int * fPersistency;
153        for (Int_t i = 0; i < side_EC; i++)
154          for (Int_t ii = 0; ii < side_EC; ii++)
155            integr[i][ii][gtu] = 0;
156
157        //run over all the pixels in one GTU
158        Int_t flag = 0;
159        for (Int_t iteratore = 0; iteratore < matsize; iteratore++) {
160          Int_t nn_PMT = it2->second->FrontEnd()->PmtChannel(iteratore)
161              / (side_PMT);
162          Int_t mm_PMT = it2->second->FrontEnd()->PmtChannel(iteratore)
163              % (side_PMT);
164
165          if (Flag_debug_table)
166            Debug_Table1(
167                it2->second->FrontEnd()->Pmt(iteratore)->Geometry()->Position(
168                    nn_PMT, mm_PMT).X(),
169                it2->second->FrontEnd()->Pmt(iteratore)->Geometry()->Position(
170                    nn_PMT, mm_PMT).Y(),
171                it2->second->FrontEnd()->Pmt(iteratore)->Geometry()->Position(
172                    nn_PMT, mm_PMT).Z(), it2->first, it1->first,
173                it2->second->FrontEnd()->Row(iteratore),
174                it2->second->FrontEnd()->Column(iteratore));
175
176          Int_t nn_EC = it2->second->FrontEnd()->Row(iteratore);
177          Int_t mm_EC = it2->second->FrontEnd()->Column(iteratore);
178
179          Int_t center_row, center_col, center_row_aux, center_col_aux;
180          center_row = nn_EC;
181          center_col = mm_EC;
182
183          //transform
184          Transform(&center_row, &center_col);
185          center_row_aux = center_row;
186          center_col_aux = center_col;
187
188          //inverts
189          Transform(&center_row_aux, &center_col_aux);
190
191          Int_t new_row, new_col;
192          // since a 3X3 box is build you need to stay within 1 pixel border from the edge of the EC
193          if (center_row > 0 && center_row < side_EC - 1 && center_col > 0
194              && center_col < side_EC - 1) {
195            //build the 3X3 box
196            for (Int_t x = -1; x < 2; x++) {
197              for (Int_t y = -1; y < 2; y++) {
198                new_row = center_row + x;
199                new_col = center_col + y;
200
201                //Inverts transformation functions
202                Transform(&new_row, &new_col);
203                Int_t channel = it2->second->FrontEnd()->ChanRowCol(new_row,
204                    new_col);
205
206                nn_PMT = it2->second->FrontEnd()->PmtChannel(channel)
207                    / (side_PMT);
208                mm_PMT = it2->second->FrontEnd()->PmtChannel(channel)
209                    % (side_PMT);
210
211                if (Flag_debug_table)
212                  Debug_Table2(
213                      it2->second->FrontEnd()->Pmt(channel)->Geometry()->Position(
214                          nn_PMT, mm_PMT).X(),
215                      it2->second->FrontEnd()->Pmt(channel)->Geometry()->Position(
216                          nn_PMT, mm_PMT).Y(),
217                      it2->second->FrontEnd()->Pmt(channel)->Geometry()->Position(
218                          nn_PMT, mm_PMT).Z(), x, y, it2->first, it1->first,
219                      channel);
220
221                //excess is calculated
222                Int_t eccesso = it2->second->GetCounter(channel) - fPixelThr
223                    + 1;
224
225                //if above threshold
226                if (it2->second->GetCounter(channel) >= fPixelThr) {
227                  // calculate the excess of the 3X3 box
228                  //assigning it to the central pixel
229                  integr[center_row_aux][center_col_aux][gtu] += eccesso;
230                  // in this GTU something is above threshold
231                  flag = 1;
232                }
233
234              } //for(Int_t y=-1;y<2;y++){
235            } //     for(Int_t x=-1;x<2;x++)
236
237          } // if(center_row>0 && center_row<11 && center_col>0 && center_col<11)
238        } // for(Int_t iteratore=0;iteratore<matsize ;iteratore++) {
239
240        //if within this EC and GTU something is above threshold....
241        if (flag == 1) {
242          // update the persistency count
243          pers += 1;
244          // if now the persistency  is above a certain value...
245          if (pers >= fPersistency) {
246            for (Int_t i = 0; i < side_EC; i++) {
247              for (Int_t ii = 0; ii < side_EC; ii++) {
248                // calculates the sum of the excesses for each pixel in the
249                // fPersistency GTUs previous to the one under analysis
250                sum_int[i][ii] = 0;
251                for (Int_t iii = 0; iii < fPersistency; iii++)
252                  sum_int[i][ii] += integr[i][ii][iii];
253
254                // is there one pixel with integration above the threshold?
255                if (sum_int[i][ii] >= fIntegration) {
256                  flag2 = 1;
257                  //TRIGGER!!!!
258                  if (maximum_int < sum_int[i][ii]) {
259                    // is the trigger the maximum?
260                    maximum_int = sum_int[i][ii];
261                    row_up = i;
262                    col_up = ii;
263                    GTU_up = it2->first;
264                    chip_ID_up = it1->first;
265                    Int_t FEE_ch = it2->second->FrontEnd()->ChanRowCol(i, ii);
266                    PdmID = it2->second->FrontEnd()->Pmt(FEE_ch)->Cell()->Id();
267                  }
268                }
269              }
270            }
271          }
272        } else { //if in this GTU no pixel above threshold...
273          //reinitialize all the counters
274          pers = 0;
275          for (Int_t i = 0; i < side_EC; i++)
276            for (Int_t ii = 0; ii < side_EC; ii++)
277              for (Int_t iii = 0; iii < fPersistency; iii++)
278                integr[i][ii][iii] = 0;
279        }
280      }
281      gtu_count++;
282    } //  for(it2=m2.begin(); it2!=m2.end(); it2++) {
283
284  } // for(it1=m1.begin(); it1!=m1.end(); it1++) {
285
286  if (flag2 == 1) {
287    //FILLING the INFORMATION VECTOR
288    fTracks.push_back(PTTTriggerSegment());
289    PTTTriggerSegment & trkseg = fTracks.back();
290
291    trkseg.SetNumEvt(EEvent::GetCurrent()->GetHeader()->GetNum());
292    trkseg.SetMaxCount(maximum_int);
293    trkseg.SetMaxGtu(GTU_up);
294    trkseg.SetMaxRow(row_up);
295    trkseg.SetMaxCol(col_up);
296    trkseg.SetMaxChipID(chip_ID_up);
297    trkseg.SetPdmID(PdmID);
298    SetTrigger();
299
300    pData->SetPTTTrigger(EEvent::GetCurrent()->GetHeader()->GetNum());
301    pData->SetPTTTriggerMaxCounts(maximum_int);
302    pData->SetPTTTriggerMaxGtu(GTU_up);
303    pData->SetPTTTriggerMaxRow(row_up);
304    pData->SetPTTTriggerMaxCol(col_up);
305    pData->SetPTTTriggerChipID(chip_ID_up);
306    pData->SetPTTTriggerPdmID(PdmID);
307
308    FillEEvent(EEvent::GetCurrent()->GetHeader()->GetNum());
309
310    Msg(EsafMsg::Info) << "TRIGGERING... NO MORE THAN 1 TRIGGER X 1 Macro Cell "
311        << MsgDispatch;
312    Msg(EsafMsg::Info) << "Integrated " << maximum_int << " Chip ID "
313        << chip_ID_up << " GTU ID " << GTU_up << MsgDispatch;
314
315  } else {
316    pData->SetPTTTrigger(-100);
317    pData->SetPTTTriggerMaxCounts(-100);
318    pData->SetPTTTriggerMaxGtu(-100);
319    pData->SetPTTTriggerMaxRow(-100);
320    pData->SetPTTTriggerMaxCol(-100);
321    pData->SetPTTTriggerChipID(-100);
322    pData->SetPTTTriggerPdmID(-100);
323  }
324
325  Clear();
326
327}
328
329//______________________________________________________________________________
330void PTTTrigger::Clear() {
331/////////////////////////////////////////////////////////////////////////////////////
332//clear/////////////////////////////////////////////////////////////////////////////////////
333/////////////////////////////////////////////////////////////////////////////////////
334
335}
336
337//______________________________________________________________________________
338void PTTTrigger::Transform(Int_t * row, Int_t * col) {
339  //trasformation
340
341  //transformation
342  //first part
343  if (*row < side_PMT && *col < side_PMT) {
344    Int_t aux = *col;
345    *col = *row;
346    *row = aux;
347  }
348
349  if (*row > side_PMT - 1 && *col < side_PMT) {
350    Int_t aux = *col;
351    *col = *row - side_PMT;
352    *row = aux + side_PMT;
353  }
354
355  if (*row < side_PMT && *col > side_PMT - 1) {
356    Int_t aux = *col;
357    *col = *row + side_PMT;
358    *row = aux - side_PMT;
359  }
360
361  if (*row > side_PMT - 1 && *col > side_PMT - 1) {
362    Int_t aux = *col;
363    *col = *row;
364    *row = aux;
365  }
366
367  //second part
368  if (*row < side_PMT)
369    *row += side_PMT;
370  else
371    *row -= side_PMT;
372
373}
374
375void PTTTrigger::FillEEvent(Int_t Triggering) {
376  vector<PTTTriggerSegment>::iterator it;
377
378  for (it = fTracks.begin(); it != fTracks.end(); it++) {
379
380    if (EEvent::GetCurrent()) {
381
382      EPTTTriggerDataAdder a(&(*it), Triggering, 1);
383      EEvent::GetCurrent()->Fill(a);
384    }
385
386  }
387  fTracks.clear();
388
389}
390
391//PIXEL POSITION
392void PTTTrigger::Debug_Table1(Float_t x, Float_t y, Float_t z,
393    Int_t analyzed_gtu, Int_t chip_id, Int_t row, Int_t col) {
394  FILE *ff = fopen("Debug_Table1.txt", "a");
395
396  fprintf(ff, "%f  %f  %f  %d  %d  %d  %d\n", x, y, z, analyzed_gtu, chip_id,
397      row, col);
398
399  fclose(ff);
400}
401
402//INTEGRATION BOX POSITION
403void PTTTrigger::Debug_Table2(Float_t x, Float_t y, Float_t z, Int_t i,
404    Int_t dir, Int_t analyzed_gtu, Int_t chip_id, Int_t box) {
405  FILE *ff = fopen("Debug_Table2.txt", "a");
406
407  fprintf(ff, "%f  %f  %f  %d  %d  %d  %d  %d\n", x, y, z, i, dir, analyzed_gtu,
408      chip_id, box);
409
410  fclose(ff);
411}
Note: See TracBrowser for help on using the repository browser.