StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
getCluAvgEff.C
1 
2 void getCluAvgEff(Char_t* signalFile="signalShapes.root", Int_t diskNr=3, Bool_t onlyQuadB=false)
3 {
4  gStyle->SetPalette(1);
5  gStyle->SetOptStat(0);
6  plFgtQuads_dbMap();
7 
8  TFile f(signalFile);
9  //TFile f("signalShapes.root");
10  char buffer[100];
11 
12  sprintf(buffer,"allClusterCountsDisk_%d",diskNr);//counting different for this histo
13  TH2D* hEff=(TH2D*)f.Get(buffer);
14 
15  sprintf(buffer,"radioDiskNonEff_%d",diskNr-1);
16  TH2D* hNonEff=(TH2D*)f.Get(buffer);
17  // hEff->Rebin2D(2,2); //merge wo bins...
18  // hNonEff->Rebin2D(2,2);
19 
20  Double_t max=hEff->GetXaxis()->GetXmax();
21  Double_t min=hEff->GetXaxis()->GetXmin();
22  Int_t minCounts=10;
23  // cout <<"max: " << max << " min: " << min << " numBins: "<< h->GetNbinsX() <<endl;
24 
25  TH2D OverallEff("overallEff","overallEff",hEff->GetNbinsX(),min,max,hEff->GetNbinsY(),min,max);
26  TH2D OverallCounts("overallCounts","overallCounts",hEff->GetNbinsX(),min,max,hEff->GetNbinsY(),min,max);
27  TH2D OverallFound("overallFound","overallFound",hEff->GetNbinsX(),min,max,hEff->GetNbinsY(),min,max);
28  OverallEff.GetXaxis()->SetTitle("x [cm]");
29  OverallEff.GetYaxis()->SetTitle("y [cm]");
30  OverallCounts.GetXaxis()->SetTitle("x [cm]");
31  OverallCounts.GetYaxis()->SetTitle("y [cm]");
32  OverallFound.GetXaxis()->SetTitle("x [cm]");
33  OverallFound.GetYaxis()->SetTitle("y [cm]");
34 
35 
36 
37 Double_t eff=0;
38 Int_t count=0;
39 
40  Int_t sumEff=0;
41  Int_t sumNonEff=0;
42 
43  Double_t overallErr=0;
44 for(Int_t i=1;i<hEff->GetNbinsX()+1;i++)
45  {
46  for(Int_t j=1;j<hEff->GetNbinsY()+1;j++)
47  {
49  if(onlyQuadB)
50  {
51  if(i<(hEff->GetNbinsX()-1)/2|| j>(hEff->GetNbinsY()-1)/2)
52  continue;
53  }
54  //Int_t gBin=h->GetBin(i,j);
55  // Double_t xpos=h->GetXaxis()->GetBinCenter(gBin);
56  // Double_t ypos=h->GetYaxis()->GetBinCenter(gBin);
57  // if(xpos<0 || ypos >0)
58  // continue;
59 
60 
61  Int_t numEff=hEff->GetBinContent(i,j);
62  Int_t numNonEff=hNonEff->GetBinContent(i,j);
63  Int_t numCounts=numEff+numNonEff;
64  // cout <<"numEff: " << numEff <<" nonEff: " << numNonEff << " counts: " << numCounts<<endl;
65  Double_t efficiency=0;
66  if(numCounts>minCounts)
67  {
68  efficiency=(Double_t)numEff/(Double_t)numCounts;
69  }
70  if(numEff>0)
71  {
72  sumEff+=numEff;
73  sumNonEff+=numNonEff;
74  }
75  OverallEff.SetBinContent(i,j,efficiency);
76  OverallCounts.SetBinContent(i,j,numCounts);
77  OverallFound.SetBinContent(i,j,numEff);
78 
79  if(numCounts>minCounts)
80  {
81  //
82  Double_t relErr=999;
83  Double_t err=999;
84  if(numEff>0)
85  {
86  relErr=sqrt((1/(Double_t)(numEff+numNonEff)+1/(Double_t)numEff));
87  err=relErr*efficiency;
88 
89  }
90  else
91  {
92  err=1/sqrt(numNonEff);
93  }
94  eff+=efficiency/(err*err);
95  overallErr+=(1/(err*err));
96 
97  count++;
98  }
99  }
100  }
101 
102  Double_t avgEff=sumEff/(Double_t)(sumEff+sumNonEff);
103  //binomial error, beware of 0
104  Double_t errOnEffNum=sqrt(avgEff*(1-avgEff)*(sumEff+sumNonEff));
105  Double_t altErr= ((Double_t)1/(Double_t)(sumEff+sumNonEff))*sqrt(sumEff*(Double_t)(1-sumEff/(sumEff+sumNonEff)));
106  sprintf(buffer,"Average Efficiency is %f +- %f",avgEff,altErr);
107  //15 degrees
108  Float_t rotationRadians=-0.261799388;
109  Float_t rotationRadians2=-1.83259571;
110  Float_t rotation=-15;
111  Float_t innerR=11.5;
112  Float_t outerR=38;
113  TLatex t1(-30,0,buffer);
114  TArc outerA(0,0,outerR,90+rotation,-90+rotation);
115  TArc innerA(0,0,innerR,90+rotation,-90+rotation);
116  // cout <<"cos rot: " << cos(rotation) <<" sin: " << sin(rotation) <<endl;
117  TLine l1(innerR*cos(rotationRadians),innerR*sin(rotationRadians),outerR*cos(rotationRadians),outerR*sin(rotationRadians));
118  TLine l2(innerR*cos(rotationRadians2),innerR*sin(rotationRadians2),outerR*cos(rotationRadians2),outerR*sin(rotationRadians2));
119  l1.SetLineWidth(3);
120  l2.SetLineWidth(3);
121  outerA.SetLineWidth(3);
122  innerA.SetLineWidth(3);
123 
124  // plGood();
125 
126  int apvMap[]={5,6,7,8,9,12,13,14,15,16};
127  char buffer2[200];
128  TCanvas* c2=new TCanvas("ceff","ceff",1,1,800,800);
129  for(int i=0;i<10;i++)
130  {
131  // sprintf(buffer2,"cD%d_apv%d",diskNr,i);
132  // c2=new TCanvas(buffer2,buffer2,1,1,800,800);
133  if(i==0)
134  {
135  sprintf(buffer2,"overallEff_D%d.pdf(",diskNr);
136  Char_t buffer3[100];
137  sprintf(buffer3,"overallEff_D%d.png",diskNr);
138  OverallEff.Draw("colz");
139  outerA.Draw();
140  innerA.Draw();
141  OverallEff.Draw("same colz");
142 
143  t1.Draw();
144  l1.Draw();
145  l2.Draw();
146  c2.SaveAs(buffer3);
147  }
148  else
149  {
150  if(i==9)
151  sprintf(buffer2,"overallEff_D%d.pdf)",diskNr);
152  else
153  sprintf(buffer2,"overallEff_D%d.pdf",diskNr);
154  }
155  OverallEff.Draw("colz");
156  outerA.Draw();
157  innerA.Draw();
158 
159  OverallEff.Draw("same colz");
160  t1.Draw();
161  l1.Draw();
162  l2.Draw();
163  sprintf(buffer,"overallEff_D%d_APV%d.pdf",diskNr,apvMap[i]);
164  plAPV(apvMap[i],1,'B');
165  c2->SaveAs(buffer);
166  cout <<"printing to " << buffer2 <<endl;
167  c2->Print(buffer2,"pdf");
168  }
169  TCanvas c("effs","effs",1,1,800,800);
170 
171 
172  //counts
173 
174  TCanvas cC("counts","counts",1,1,800,800);
175  OverallCounts.Draw("colz");
176  outerA.Draw();
177  innerA.Draw();
178 
179  OverallCounts.Draw("same colz");
180  l1.Draw();
181  l2.Draw();
182  cC.SaveAs("overallCounts.png");
183 
184 
185  TCanvas cf("found","found",1,1,800,800);
186  OverallFound.Draw("colz");
187  outerA.Draw();
188  innerA.Draw();
189  OverallFound.Draw("same colz");
190  l1.Draw();
191  l2.Draw();
192  cf.SaveAs("overallFound.png");
193 
194  sprintf(buffer,"overall found");
195  sprintf(buffer2,"cFoundD",diskNr);
196  c2=new TCanvas(buffer2,buffer2,1,1,800,800);
197  TLatex t2(-30,0,buffer);
198  for(int i=0;i<10;i++)
199  {
200  if(i==0)
201  {
202  sprintf(buffer2,"overallFound_D%d.pdf(",diskNr);
203  Char_t buffer3[100];
204  sprintf(buffer3,"overallFound_D%d.png",diskNr);
205  sprintf(buffer2,"overallFound_D%d.pdf(",diskNr);
206  OverallFound.Draw("colz");
207  outerA.Draw();
208  innerA.Draw();
209  OverallFound.Draw("same colz");
210  l1.Draw();
211  l2.Draw();
212  c2.SaveAs(buffer3);
213  }
214  else
215  {
216  if(i==9)
217  sprintf(buffer2,"overallFound_D%d.pdf)",diskNr);
218  else
219  sprintf(buffer2,"overallFound_D%d.pdf",diskNr);
220  }
221  OverallEff.Draw("colz");
222  outerA.Draw();
223  innerA.Draw();
224 
225  // OverallFound.Draw("same colz");
226  t2.Draw();
227  l1.Draw();
228  l2.Draw();
229  sprintf(buffer,"overallFound_D%d_APV%d.png",diskNr,apvMap[i]);
230  plAPV(apvMap[i],1,'B');
231  c2.SaveAs(buffer);
232  cout <<"printing to " << buffer2 <<endl;
233  c2->Print(buffer2,"pdf");
234  }
235 
236 
237 
238 
239 
240 
241  // c.SaveAs("overallEff.C");
242 //cout <<"avg eff: " << eff/count <<endl;
243  cout <<"Hits found: " << sumEff <<" Hits not found: " << sumNonEff<< " efficiency: " <<avgEff<<" +- " << altErr<<endl;
244 
245 }
246 
247 double pi=2.*acos(0.);
248 double rad2deg=pi/180.;
249 //==================
250 enum {mxFgtElect=30720, mxFgtApv=22, mxFgtPln=2, kFgtPlnP=0, kFgtPlnR=1 };
251 char *plnC[mxFgtPln]={"Phi","R"};
252 
253 // plane[R,P], r1(cm), r2(cm), phi1(rad), phi2(rad)
254 
256  int electId,geoId, rdo,arm,apv,chan, disc,strip;
257  float r1, r2, phi1, phi2; // cm, rad
258  char layer; // P,R
259  char quad; // A-D
260  char name[10];
261  int stat; // 0 is good
262  float ped, sigPed;
263 };
264 
265 FgtStripDbItem stripDb[mxFgtElect];
266 
267 
268 // mmmmmmmmmmmmmmmmmmmmmmmmmmmmmmmmmmmmm
269 // MAIN
270 // mmmmmmmmmmmmmmmmmmmmmmmmmmmmmmmmmmmmm
271 
272 void plFgtQuads_dbMap() { // ipl: 0=Phi, 1=R; apv<0=all
273  memset(stripDb,0, sizeof(stripDb));
274  readFgtStripMap_CSV("fgtMapDump-jan28.csv");
275  printf("reading of CSV done\n");
276 
277  for(int i=0;i<10;i++) {
278  if(stripDb[i].rdo<=0) continue;
279  printStripItem(stripDb+i);
280  }
281 
282  printf(" execute: plAPV() , plPlane() , plGood(1) , plRstrips(1); \n");
283  // plAPV(); // plAPV( int apv=0, int disc=1, char quad='A' )
284  //plPlane(1); // 0=P, 1=R
285 
286  //plGood(2);
287  //plRstrips(1); // ioct=0=L, ioct=1=S
288 }
289 
290 //--------------------------------------------------
291 //--------------------------------------------------
292 <//--------------------------------------------------
293 void plGood( int disc=1, int ipl=-1 ) {
294  double par_fgtQuadMaxR=39.7;//cm, limit for displaiyng
295  int par_globView=0;
296 
297  TString canN=Form("fgt_disc%d_good",disc);
298  if(ipl>=0) canN=Form("fgt_%d_%s",disc,plnC[ipl]);
299  TH2F *hQuad=new TH2F( "aa","aa",20,-par_fgtQuadMaxR,par_fgtQuadMaxR,20,-par_fgtQuadMaxR,par_fgtQuadMaxR);
300  gStyle->SetOptStat(0) ;
301  gStyle->SetFillColor(kWhite);
302  c=new TCanvas( canN,canN,700,720);
303  hQuad->Draw("SAME");
304  // gPad->SetGrid();
305  int nstr=0;
306 
307  nstr+=plotStrips( -1, disc, 'A', ipl, -15);
308  nstr+=plotStrips( -1, disc, 'B', ipl, -15-90);
309  nstr+=plotStrips( -1, disc, 'D', ipl, 75);
310  nstr+=plotStrips( -1, disc, 'C', ipl, 75+90);
311 
312  TString tit;
313  tit=Form("FGT: disc %d, good strips=%d ; STAR ref. X (cm); STAR ref. Y (cm)",disc, nstr);
314 
315  if(ipl>=0) tit=Form("FGT: disc %d plane=%s, GOOD strip=%d ; local X (cm); local Y (cm)",disc, plnC[ipl], nstr);
316 
317  hQuad->SetTitle(tit);
318  c->Print();
319 }
320 
321 
322 
323 //--------------------------------------------------
324 //--------------------------------------------------
325 //--------------------------------------------------
326 void plRstrips( int ioct ) {
327  int ipl=1;
328  int disc=1; char quad='A';
329  double par_fgtQuadMaxR=39.7;//cm, limit for displaiyng
330  int par_globView=0;
331 
332  TString canN=Form("fgt_%d%c_%s_str_oct%d",disc,quad,plnC[ipl],ioct);
333  TH2F *hQuad=new TH2F( "aa","aa",200,-0.5,par_fgtQuadMaxR,200,-0.5,par_fgtQuadMaxR);
334  gStyle->SetOptStat(0) ;
335  gStyle->SetFillColor(kWhite);
336  c=new TCanvas( canN,canN,700,720);
337  hQuad->Draw("SAME");
338  // gPad->SetGrid();
339  int nstr=0;
340 
341  nstr=plotStrips( 0+17*ioct, disc, quad, ipl,0,1);
342  nstr=plotStrips( 1+17*ioct, disc, quad, ipl,0,2);
343  nstr=plotStrips( 2+17*ioct, disc, quad, ipl,0,8);
344  nstr=plotStrips( 3+17*ioct, disc, quad, ipl,0,4);
345  nstr=plotStrips( 4+17*ioct, disc, quad, ipl,0,6);
346 
347  TString tit;
348  tit=Form("FGT: quad %d%c, plane=%s, strip=%d ; local X (cm); local Y (cm)",disc, quad, plnC[ipl], nstr);
349  hQuad->SetTitle(tit);
350  c->Print();
351 }
352 
353 //--------------------------------------------------
354 //--------------------------------------------------
355 //--------------------------------------------------
356 void plPlane( int ipl=0, int disc=1, char quad='A' ) {
357  double par_fgtQuadMaxR=39.7;//cm, limit for displaiyng
358  int par_globView=0;
359 
360  TString canN=Form("fgt_%d%c_%s",disc,quad,plnC[ipl]);
361  TH2F *hQuad=new TH2F( "aa","aa",200,-0.5,par_fgtQuadMaxR,200,-0.5,par_fgtQuadMaxR);
362  gStyle->SetOptStat(0) ;
363  gStyle->SetFillColor(kWhite);
364  c=new TCanvas( canN,canN,700,720);
365  hQuad->Draw("SAME");
366  // gPad->SetGrid();
367  int nstr=0;
368 
369  nstr=plotStrips( -1, disc, quad, ipl);
370 
371  TString tit;
372  tit=Form("FGT: quad %d%c, plane=%s, strip=%d ; local X (cm); local Y (cm)",disc, quad, plnC[ipl], nstr);
373  hQuad->SetTitle(tit);
374  c->Print();
375 }
376 
377 //--------------------------------------------------
378 //--------------------------------------------------
379 //--------------------------------------------------
380 void plAPV( int apv=0, int disc=1, char quad='A' ) {
381 
382  double par_fgtQuadMaxR=39.7;//cm, limit for displaiyng
383  int par_globView=0;
384 
385  TString canN=Form("fgt_%d%c_apv%d",disc,quad,apv);
386  TH2F *hQuad=new TH2F( "aa","aa",200,-0.5,par_fgtQuadMaxR,200,-0.5,par_fgtQuadMaxR);
387  gStyle->SetOptStat(0) ;
388  gStyle->SetFillColor(kWhite);
389  // c=new TCanvas( canN,canN,700,720);
390  hQuad->Draw("SAME");
391  // gPad->SetGrid();
392  int nstr=0;
393  TText *tx;
394  // if(disc==1 && quad=='A')
395  // tx=new TText(5,40.3,"same mapping for quads 1A,1C,2A,3B,4A,5B,6A");
396  // if(disc==1 && quad=='B')
397  // tx=new TText(5,40.3,"same mapping for quads 1B,1D,2B,3A,4B,5A,6B");
398  // tx->Draw(); tx->SetTextSize(0.03);
399 
400  nstr=plotStrips( apv, disc, quad,-1);
401 
402  TString tit;
403  if( apv>=0) tit=Form("FGT: quad %d%c, APV=%d, Nstrip=%d ; local X (cm); local Y (cm)",disc, quad, apv, nstr);
404  hQuad->SetTitle(tit);
405  // c->Print(Form("fgt_%d%c_apv%02d.ps",disc,quad,apv));
406 }
407 
408 //--------------------------------------------------
409 //--------------------------------------------------
410 //--------------------------------------------------
411 int plotStrips(int apv, int disc, char quad, int ipl, int phiDeg=0, int iCol0) {
412 
413  Float_t rotDeg=105;
414  Bool_t secondOct=false;
415  // if(apv>4)
416  // {
417  // secondOct=true;
418  // }
419 
420  // apv-=5;
421  int lineWidth=0;
422  int nstr=0;
423  for(int i=0;i<mxFgtElect;i++) {
424  FgtStripDbItem *S=stripDb+i;
425  if(S->geoId<0) continue;
426  if(S->disc!=disc) continue;
427  if(S->quad!=quad) continue;
428  if(ipl>=0)
429  if(S->layer!=plnC[ipl][0]) continue;
430  if(phiDeg)
431  if(S->stat) continue;
432 
433  //printStripItem(S);
434  if(nstr<5) printStripItem(S);
435 
436  int iCol=kBlue;
437  if( S->strip <=360) iCol=kRed;
438 
439  if(apv>=0) {
440  iCol=kBlue;
441  lineWidth=1;
442  if(S->apv!=apv) continue;
443  }
444 
445 
446  if(iCol0) iCol=iCol0;
447 
448  nstr++;
449  if(apv>=0) { //print strip ID on the right
450  tx=new TText(42-(nstr%2)*2., 42.-(nstr+S->chan/32)/3., S->name+2);
451  if(S->layer=='P') tx->SetTextColor(kBlue);
452  //B tx->Draw();
453 tx->SetTextSize(0.014);
454  if(nstr==1) { tx=new TText(39.7, 42.9,"by electID"); // tx->Draw();
455  tx->SetTextSize(0.02);}
456  }
457 
458  if( S->layer=='P') { //ppppppppppppppppppppp
459  int keep= ( S->strip%720)==0;
460  if(S->strip%7!=0 && !keep) continue;
461 
462 
463  // phiDeg+=0.94;
464  // if(phiDeg>2*pi)
465  // phiDeg-=(2*pi);
466  // cout <<"phiDeg: " << phiDeg;
467  // phiDeg=phiDeg-rotDeg;
468  // cout << " now : " << phiDeg;
469  phiDeg=-rotDeg;
470  double phi1=S->phi1+phiDeg*rad2deg, phi2=S->phi2+phiDeg*rad2deg;
471 
472  // phi1-=rotDeg*rad2deg;
473  // phi2-=rotDeg*rad2deg;
474  cout <<" phi1: " << phi1 <<" phi2: " << phi2 <<endl;
475  double x1=S->r1*cos(phi1), x2=S->r2*cos(phi2);
476  double y1=S->r1*sin(phi1), y2=S->r2*sin(phi2);
477 
478 
479  ln=new TLine(x1,y1,x2,y2);
480  ln->Draw(); ln->SetLineColor(iCol);
481  ln->SetLineWidth( lineWidth);
482 
483  if(S->strip%35==0 || keep) {
484  double x3=x2+1., y3=y2+0.6;
485  ln=new TLine(x2*1.005,y2*1.005,x3,y3);
486  ln->Draw(); ln->SetLineColor(iCol);
487  ln->SetLineWidth(0);
488  tx=new TText(x3,y3,Form("%c%03d",S->layer,S->strip));
489  tx->Draw(); tx->SetTextColor(iCol); tx->SetTextSize(0.02);
490  tx->SetTextAngle(phiDeg+10.);
491  }
492  } // end of P-plane
493 
494  if( S->layer=='R' ) { //RRRRRRRRRRRRRRRRRRRRRRR
495  int kk;
496  if( S->strip <400) { kk=S->strip;}
497  else { kk=S->strip-400;}
498  int keep= (kk %279)==0;
499  if(kk%7!=0 && !keep) continue;
500 
501  double phi1Deg=S->phi1/rad2deg +phiDeg;
502  double phi2Deg=S->phi2/rad2deg+phiDeg;
503 
504  // phi1Deg-=rotDeg;
505  // phi2Deg-=rotDeg;
506 
507  TArc * ar=new TArc(0.,0,S->r1, phi1Deg, phi2Deg);
508  ar->SetFillStyle(0); ar->SetLineColor(iCol);
509  ar->Draw("only");
510  ar->SetLineWidth( lineWidth);
511 
512  if(kk%35==0 || keep) {
513  double x2,y2,x3,y3,dx,dy;
514  double phiAvr=(S->phi1+S->phi2)/2.;
515  if(phiDeg) { // rotate P-strips
516  double delPhi=phiDeg*rad2deg;
517  phiAvr+=delPhi;
518  }
519 
520  x2=S->r1*cos(phiAvr); y2=S->r2*sin(phiAvr);
521  x3=x2-0.3; y3=y2+0.6;
522  dx=-1.; dy=0.2;
523  if(S->apv<5 ) {
524  x3=x2+0.2; y3=y2+0.7;
525  dy=.3; dx=-0.7;
526  }
527  ln=new TLine(x2,y2,x3,y3);
528  ln->Draw(); ln->SetLineColor(iCol);
529  ln->SetLineWidth(0);
530  tx=new TText(x3+dx,y3+dy,Form("%c%03d",S->layer,S->strip));
531  tx->Draw(); tx->SetTextSize(0.02);//tx->SetTextColor(iCol);
532  tx->SetTextAngle(phiDeg-10.);
533  }
534  } // end of R-plane
535 
536  } // end of strip loop
537  return nstr;
538 }
539 
540 //--------------------------------------------------
541 //--------------------------------------------------
542 //--------------------------------------------------
543 
544 void printStripItem( FgtStripDbItem *S) {
545  printf("%d%c_%c-plane strip=%d geo=%d elec=%d apv%d chan%d, R[%.1f, %.1f]cm Phi[%.2f, %.2f]deg %s stat=%d ped=%.1f\n",S->disc,S->quad,S->layer,S->strip, S->geoId,S->electId,S->apv, S->chan, S->r1, S->r2, S->phi1/rad2deg, S->phi2/rad2deg, S->name, S->stat, S->ped);
546 
547 }
548 
549 //--------------------------------------------------
550 //--------------------------------------------------
551 //--------------------------------------------------
552 void readFgtStripMap_CSV(char *fname){
553  FILE *fd=fopen(fname,"r"); assert(fd);
554 
555  const int mx=1000;
556  char buf[mx];
557  int k=0, nok=0;
558 
559  float ord,lowSpan,upSpan;
560 
561  FgtStripDbItem S;
562  while(true) {
563  char * ret=fgets(buf,mx,fd);
564  if(ret==0) break;
565  if(buf[0]=='#') continue;
566 
567  char *item=strtok(buf,","); // init 'strtok'
568  int i=0;
569 
570  do {
571  i++;
572  //printf("i=%d, item=%s=\n",i,item);
573  switch(i){
574  case 1: S.electId=atoi(item); break;
575  case 2: S.geoId=atoi(item); break;
576 
577  case 3: S.rdo=atoi(item); break;
578  case 4: S.arm=atoi(item); break;
579  case 5: S.apv=atoi(item); break;
580  case 6: S.chan=atoi(item); break;
581 
582  case 7: S.disc=atoi(item); break;
583  case 8: S.quad=item[0]; break;
584  case 9: S.layer=item[0]; break;
585  case 10: S.strip=atoi(item); break;
586 
587  case 11: sscanf(item,"%f",&ord); break;
588  case 12: sscanf(item,"%f",&lowSpan); break;
589  case 13: sscanf(item,"%f",&upSpan); break;
590 
591  case 14: sscanf(item,"%s",S.name); break;
592 
593  case 15: S.stat=atoi(item); break;
594  case 16: sscanf(item,"%f",&S.ped); break;
595  case 17: sscanf(item,"%f",&S.sigPed); break;
596  default:
597  }
598 
599  } while (item=strtok(0,",")); // advance by one name
600 
601  k++;
602  if(S.geoId<0) {
603  stripDb[S.electId].geoId=-1;
604  continue; // skip invalid records
605  }
606  nok++;
607 
608  if (S.layer=='P') {
609  S.phi1=S.phi2=ord;
610  S.r1=lowSpan;
611  S.r2=upSpan;
612  } else {
613  S.r1=S.r2=ord;
614  S.phi1=lowSpan;
615  S.phi2=upSpan;
616  }
617 
618  // printf("k=%d nok=%d %d-%c\n",k,nok,S.quad,S.layer ); printStripItem(&S);
619 
620  assert(S.electId>=0);
621  assert(S.electId<mxFgtElect);
622 
623  stripDb[S.electId]=S;
624 
625  //if(k>282) break;
626  }
627  printf ("total # of APV channels mapped %d, see records %d\n",nok,k);
628 }
629 
630 
631 //--------------------------------------------------
632 //--------------------------------------------------
633 //--------------------------------------------------
634 void doAll() {
635  // print APVs in disc 1
636  char quad='B';
637  // for(int apv=0; apv<5;apv++) plAPV(apv,1,quad);
638  // for(int apv=17; apv<22;apv++) plAPV(apv,1,quad);
639  //for(int apv=5; apv<10;apv++) plAPV(apv,1,quad);
640  //for(int apv=12; apv<17;apv++) plAPV(apv,1,quad);
641 
642  for(int disc=1;disc<=6;disc++) plGood(disc);
643 
644 }