StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
tpc_FCFReader.cxx
1 #include <string.h>
2 #include <stdio.h>
3 
4 #include <daqFormats.h>
5 #include <rtsSystems.h>
6 #include <rtsLog.h>
7 
8 #include <fcfClass.hh>
9 
10 //#include <evpSupport.h>
11 //#include <tpcReader.h>
12 #include <TPC/padfinder.h>
13 #include <TPC/rowlen.h>
14 #include <adcLogTable.h>
15 
16 #include "daq_tpc.h"
17 
18 /************************************************/
19 // This is not a READER!
20 //
21 // It calculates clusters from the filled tpc
22 // structure...
23 //
24 /************************************************/
25 
26 #define FCF_MAX_CLUSTERS 6000
27 
28 // Assume that tpc raw structure is alread filled
29 //
30 // t0corr --> int t0[46][242], for this sector
31 // gainCorr --> u_int gain[46][242], for this sector
32 //
33 int daq_tpc::fcfReader(int sector, int *t0c, u_int *gainc, tpc_t *tpc) // from zero
34 {
35  fcfClass *fcf = NULL;
36  fcfAfterburner *fcf_after = NULL;
37 
38  int row, p, t;
39  u_short tb;
40 
41  u_int adcOff[183] ;
42  u_short cppOff[183] ;
43  u_short startFlags[183];
44  u_int output[2*FCF_MAX_CLUSTERS];
45  u_int hits = 0;
46  u_int *fcf_ptrs[3];
47  int t0[242];
48  u_int gain[242];
49 
50  memset(tpc->cl_counts, 0, sizeof(tpc->cl_counts));
51  memset(tpc->cl, 0, sizeof(tpc->cl));
52  tpc->has_clusters = 1;
53 
54  if(tpc->mode != 0) return -1; // don't do for pedestals...
55 
56  fcf = new fcfClass(TPC_ID, NULL);
57  fcf_after = new fcfAfterburner();
58 
59  for(int r=0;r<45;r++) {
60  int have_data = 0;
61  row = r+1;
62 
63  // Padrow adc data stored in "adc"
64  u_char adc[182][512];
65 
66  memset(adc, 0, sizeof(adc));
67  int adcs=0;
68  for(p=0;p<182;p++) {
69  for(t=0;t<tpc->counts[r][p];t++) {
70  u_char val;
71  tb = tpc->timebin[r][p][t];
72  val = tpc->adc[r][p][t];
73 
74  adc[p][tb] = val;
75 
76 // if((sector == 1)) {
77 // printf("cl_raw %d %d %d %d %d\n",sector,r+1,p+1,tb,val) ;
78 // }
79 
80  if(val) {
81  have_data = 1;
82  adcs++;
83  }
84  }
85  }
86 
87  // Padrow cpp stored in "cpp"
88  u_short cpp[182][32*2];
89 
90  int seqs = 0;
91 
92  memset(cpp, 0xff, sizeof(cpp));
93  for(p=0;p<182;p++) {
94  int started, cou, max_seq;
95 
96  started = cou = max_seq = 0;
97 
98  for(tb=0;tb<512;tb++) {
99 
100  if(adc[p][tb] > 0) {
101  if(started) continue;
102 
103  // start a new cluster pointer
104  started = 1;
105 
106  if(cou <= 60) // more than 31 clusters?
107  cpp[p][cou++] = tb;
108  else
109  LOG(NOTE, "Cou is begger %d pad %d tb %d",cou,p,tb);
110  }
111 
112  if(started && (adc[p][tb] == 0)) { // finish cluster
113 
114  seqs++;
115 
116  if(cou <= 61)
117  cpp[p][cou++] = tb-1;
118  else
119  LOG(NOTE, "Cou is bigger %d pad %d tb %d",cou,p,tb);
120 
121 
122  //test
123 // if(r==0) {
124 // LOG(NOTE, "r=%d p=%d t1=%d t2=%d:",r,p,cpp[p][cou-2],cpp[p][cou-1]);
125 // for(int ii=cpp[p][cou-2]; ii <= cpp[p][cou-1];ii++) {
126 // LOG(NOTE, " adc[%d][%d]=%d",p,ii,adc[p][ii]);
127 // }
128 // }
129  //
130 
131  started = 0;
132  }
133  }
134 
135  // done last tb, stop if neccessary
136  if(started) {
137  seqs++;
138  cpp[p][cou++] = tb-1;
139  }
140 
141  // LOG(NOTE, "r=%d p=%d seq=%d",r,p,cou);
142  started = 0;
143  }
144 
145 
146 
147  // LOG(NOTE, "adcs=%d seqs=%d",adcs, seqs);
148 
149  // Point to correct locations...
150  for(int i=0;i<182;i++) {
151  adcOff[i+1] = (char *)(&adc[i][0]) - (char *)(&adc[0][0]);
152  cppOff[i+1] = (char *)(&cpp[i][0]) - (char *)(&cpp[0][0]) ;
153  }
154 
155 // for(int i=0;i<182;i++) {
156 // u_short *pp = (u_short *)(((u_char *)cpp) + cppOff[i+1]);
157 // if((sector == 1)) {
158 // while(*pp != 0xffff) {
159 // printf("cpp %d %d %d %d %d\n",sector,r+1,i+1,*pp,1) ;
160 // pp++;
161 // printf("cpp %d %d %d %d %d\n",sector,r+1,i+1,*pp,1) ;
162 // pp++;
163 // }
164 // }
165 // }
166 
167 
168  /*
169  if(r == 0) {
170  for(p=0;p<182;p++) {
171  u_char *x = ((u_char *)cpp) + cppOff[p+1];
172 
173  for(tb=0;tb<512;tb++) {
174 
175  LOG(NOTE, "cpp[%d][%d] = %d = cppOff = %d",
176  p,tb,cpp[p][tb], *(x+tb));
177 
178  if(cpp[p][tb] == 0xff) break;
179  }
180  }
181  }
182  */
183 
184 
185  fcf->adcOff = adcOff;
186  fcf->cppOff = cppOff;
187 
188  fcf->timebinLo = 0;
189  fcf->timebinHi = 400;
190  fcf->chargeMin = 40;
191  fcf->set8to10(log8to10_table);
192 
193 
194  if(t0c) {
195  for(int i=0;i<182;i++) t0[i+1] = t0c[46*r + i];
196  }
197  else {
198  for(int i=0;i<182;i++) t0[i+1] = 0;
199  }
200 
201  if(gainc) {
202  for(int i=0;i<182;i++) gain[i+1] = gainc[46*r + i];
203  }
204  else {
205  for(int i=0;i<182;i++) gain[i+1] = 64;
206  }
207 
208  fcf->t0Corr = t0;
209  fcf->gainCorr = gain;
210 
211  fcf->maxClusters = FCF_MAX_CLUSTERS;
212 
213  fcf->simIn = NULL;
214  fcf->simOut = NULL;
215 
216  fcf_after->setVerbose(false);
217 
218  fcf->startFlags = startFlags;
219 
220  memset(startFlags, 0, sizeof(startFlags));
221 
222  memset(output, 0, sizeof(output));
223  u_int *out = output;
224 
225  memset(fcf_ptrs, 0, sizeof(fcf_ptrs));
226 
227  fcf->row = r+1;
228 
229  for(int mz=0;mz<3;mz++) {
230 
231  if(padfinder[r][mz].rdo == 0) continue;
232 
233  fcf->padStart = padfinder[r][mz].minpad;
234  fcf->padStop = padfinder[r][mz].maxpad;
235 
236 
237  if(fcf->padStart == 1) startFlags[fcf->padStart] |= FCF_ROW_EDGE;
238  else(startFlags[fcf->padStart]) |= FCF_BROKEN_EDGE;
239 
240  if(tpc_rowlen[r+1] == fcf->padStop)
241  startFlags[fcf->padStop] |= FCF_ROW_EDGE;
242  else
243  startFlags[fcf->padStop] |= FCF_BROKEN_EDGE;
244 
245 
246 
247  u_int words = fcf->finder((u_char *)adc,
248  (u_short *)cpp,
249  (u_int *)out);
250 
251  // LOG(NOTE, "c(%d)%d-%d padrow %d words=%d", mz, fcf->padStart, fcf->padStop, r,words);
252 
253 
254  fcf_ptrs[mz] = out;
255 
256  if(words != (*(out+1))*2 + 2) {
257  LOG(ERR, "Words should match %d vs %d (%d)",
258  words, (*(out+1))*2 + 2, *(out+1));
259  }
260 
261  out += words;
262 
263  } // mz loop
264 
265  // Can burn this row...
266 
267 // for(int ii=0;ii<3;ii++) {
268 // uint *x = fcf_ptrs[ii];
269 // if(!x) continue;
270 // int row = *x;
271 // x++;
272 // int ncl = *x;
273 // x++;
274 
275 
276 // while(ncl) {
277 // mzCentroid *cent=(mzCentroid *)x;
278 // if(sector == 1)
279 // {
280 // printf("nrcl: %d %d %d %d %d %d %d %d\n",
281 // sector,r+1,row,ii,cent->x,cent->t,cent->flags,cent->charge);
282 // }
283 // ncl--;
284 // x++;
285 // x++;
286 // }
287 // }
288 
289  fcf_after->burn(fcf_ptrs);
290 
291  int ccou=0;
292  fcfHit h;
293  tpc_cl *cld;
294 
295  while(fcf_after->next(&h)) {
296  cld = &tpc->cl[r][ccou];
297  ccou++;
298 
299 
300  cld->p = (double)h.pad/64.0 + 0.5;
301  cld->t = (double)h.tm/64.0 + 0.5;
302 
303  // LOG(WARN," pad=%lf %lf",cld->p, cld->t);
304  cld->charge = h.c;
305  cld->flags = h.f;
306  cld->t1 = h.t1;
307  cld->t2 = h.t2;
308  cld->p1 = h.p1;
309  cld->p2 = h.p2;
310 
311  }
312 
313  hits += ccou;
314  tpc->cl_counts[r] = ccou;
315  } // row loop
316 
317  delete fcf_after;
318  delete fcf;
319 
320  return hits;
321 }
Definition: daq_tpc.h:22
Definition: daq_tpc.h:11