pxar
 All Classes Namespaces Functions Variables Typedefs Friends
PixTestDacScanPh.cc
1 // -- author: Daniel Pitzl
2 // PH vs DAC
3 
4 #include <stdlib.h> // atof, atoi
5 #include <algorithm> // std::find
6 
7 #include "PixTestDacScanPh.hh"
8 #include "log.h"
9 
10 using namespace std;
11 using namespace pxar;
12 
13 ClassImp(PixTestDacScanPh)
14 
15 //------------------------------------------------------------------------------
16 PixTestDacScanPh::PixTestDacScanPh( PixSetup *a, std::string name )
17 : PixTest(a, name), fParNtrig(-1),
18  fParDAC("nada"), fParLoDAC(-1), fParHiDAC(-1),
19  fParCals(0)
20 {
21  PixTest::init();
22  init();
23 }
24 
25 //------------------------------------------------------------------------------
26 PixTestDacScanPh::PixTestDacScanPh() : PixTest()
27 {
28  // LOG(logDEBUG) << "PixTestDacScanPh ctor()";
29 }
30 
31 //------------------------------------------------------------------------------
32 bool PixTestDacScanPh::setParameter( string parName, string sval )
33 {
34  bool found(false);
35 
36  for( uint32_t i = 0; i < fParameters.size(); ++i ) {
37 
38  if( fParameters[i].first == parName ) {
39 
40  found = true;
41 
42  sval.erase(remove(sval.begin(), sval.end(), ' '), sval.end() );
43 
44  if( !parName.compare( "ntrig" ) ) {
45  fParNtrig = atoi( sval.c_str() );
46  LOG(logDEBUG) << " setting fParNtrig ->" << fParNtrig
47  << "<- from sval = " << sval;
48  }
49 
50  if( !parName.compare( "dac" ) ) {
51  fParDAC = sval;
52  LOG(logDEBUG) << " setting fParDAC ->" << fParDAC
53  << "<- from sval = " << sval;
54  }
55 
56  if( !parName.compare( "daclo" ) ) {
57  fParLoDAC = atoi(sval.c_str() );
58  LOG(logDEBUG) << " setting fParLoDAC ->" << fParLoDAC
59  << "<- from sval = " << sval;
60  }
61  if( !parName.compare( "dachi" ) ) {
62  fParHiDAC = atoi( sval.c_str() );
63  LOG(logDEBUG) << " setting fParHiDAC ->" << fParHiDAC
64  << "<- from sval = " << sval;
65  }
66  if( !parName.compare( "cals" ) ) {
67  fParCals = atoi( sval.c_str() );
68  LOG(logDEBUG) << " setting fParCals ->" << fParCals
69  << "<- from sval = " << sval;
70  }
71 
72  if( !parName.compare( "pix1" ) ) {
73  string::size_type s1 = sval.find( "," );
74  if( string::npos != s1 ) {
75  string str1 = sval.substr(0, s1);
76  int pixc = atoi( str1.c_str() );
77  string str2 = sval.substr(s1+1);
78  int pixr = atoi( str2.c_str() );
79  fPIX.push_back( make_pair( pixc, pixr ) );
80  LOG(logDEBUG) << "new coordinates: (" << fPIX[fPIX.size()-1].first
81  << " | " << fPIX[fPIX.size()-1].second << ") with currently "
82  << fPIX.size() << " entries in the Vector" << endl;
83  }
84  else {
85  fPIX.push_back( make_pair( -1, -1 ) );
86  }
87  }
88  // FIXME: remove/update from fPIX if the user removes via the GUI!
89  break;
90  }
91  }
92  return found;
93 }
94 
95 //------------------------------------------------------------------------------
96 void PixTestDacScanPh::init()
97 {
98  fDirectory = gFile->GetDirectory( fName.c_str() );
99  if( !fDirectory )
100  fDirectory = gFile->mkdir( fName.c_str() );
101  fDirectory->cd();
102 }
103 
104 // ----------------------------------------------------------------------
106 {
107  fTestTip = string( "measure pixel pulse height vs DAC");
108  fSummaryTip = string("summary plot to be implemented");
109 }
110 
111 //------------------------------------------------------------------------------
112 void PixTestDacScanPh::bookHist(string name) // general booking routine
113 {
114  LOG(logDEBUG) << "nothing done with " << name;
115 }
116 
117 //------------------------------------------------------------------------------
118 PixTestDacScanPh::~PixTestDacScanPh()
119 {
120  std::list<TH1*>::iterator il;
121  fDirectory->cd();
122  for( il = fHistList.begin(); il != fHistList.end(); ++il ) {
123  LOG(logINFO) << "Write out " << (*il)->GetName();
124  (*il)->SetDirectory( fDirectory );
125  (*il)->Write();
126  }
127 }
128 
129 //------------------------------------------------------------------------------
131 {
132  fDirectory->cd();
133  fHistList.clear();
134  PixTest::update();
135 
136  LOG(logINFO) << "PixTestCurrentVsDac::doTest() DAC = " << fParDAC;
137 
138  uint8_t maxDac = fApi->getDACRange( fParDAC ); // 15 or 255
139 
140  if( maxDac == 0 ) {
141  LOG(logINFO) << "ERROR: " << fParDAC << " is not a ROC register";
142  return;
143  }
144 
145  if( fParHiDAC > maxDac ) {
146  LOG(logINFO) << fParDAC << " range only " << maxDac;
147  fParHiDAC = maxDac;
148  }
149 
150  LOG(logINFO) << "PixTestDacScanPh::doTest() ntrig = " << fParNtrig;
151 
152  // activate one pixel per ROC:
153 
154  fApi->_dut->testAllPixels(false);
155 
156  //coordinates of the last pair = presetly set pixel
157  int32_t col = fPIX[fPIX.size()-1].first;
158  int32_t row = fPIX[fPIX.size()-1].second;
159 
160  if( col > -1 )
161  fApi->_dut->testPixel( col, row, true );
162 
163  // measure:
164 
165  uint16_t flags = 0;
166  if( fParCals ) flags = FLAG_CALS;
167  LOG(logINFO) << "flag " << flags;
168  uint8_t ctl = fApi->_dut->getDAC( 0, "CtrlReg" );
169  if( fParCals ) {
170  fApi->setDAC( "CtrlReg", 4 ); // all ROCs large Vcal
171  LOG(logINFO) << "CtrlReg 4 (large Vcal)";
172  }
173 
174  vector<pair<uint8_t, vector<pixel> > > // dac and pix
175  result = fApi->getPulseheightVsDAC( fParDAC, fParLoDAC, fParHiDAC, flags, fParNtrig );
176 
177  if( fParCals ) {
178  fApi->setDAC( "CtrlReg", ctl ); // restore
179  LOG(logINFO) << "back to CtrlReg " << ctl;
180  }
181 
182  // book:
183 
184  vector<TH1D*> hsts;
185  TH1D *h1(0);
186 
187  uint32_t nRocs = fPixSetup->getConfigParameters()->getNrocs();
188 
189  for( uint32_t roc = 0; roc < nRocs; ++roc ) {
190 
191  h1 = new TH1D( Form( "PH_vs_%s_c%02d_r%02d_C%02d",
192  fParDAC.c_str(), col, row, roc ),
193  Form( "PH vs %s c%02d r%02d C%02d",
194  fParDAC.c_str(), col, row, roc ),
195  256, -0.5, 255.5 );
196  h1->SetMinimum(0);
197  h1->SetMaximum(256);
198  setTitles( h1, Form( "%s [DAC]", fParDAC.c_str() ), "<PH> [ADC]" );
199  hsts.push_back(h1);
200  fHistList.push_back(h1);
201 
202  } // rocs
203 
204  // plot data:
205 
206  for( size_t i = 0; i < result.size(); ++i ) {
207 
208  int idac = result[i].first;
209 
210  vector<pixel> vpix = result[i].second;
211 
212  for( size_t ipx = 0; ipx < vpix.size(); ++ipx ) {
213 
214  uint8_t roc = vpix[ipx].roc();
215 
216  if( roc < nRocs &&
217  vpix[ipx].column() == col &&
218  vpix[ipx].row() == row ) {
219  h1 = hsts.at(roc);
220  h1->Fill( idac, vpix[ipx].value()); // already averaged
221  } // valid
222 
223  } // pix
224 
225  } // dac values
226 
227  for( size_t roc = 0; roc < nRocs; ++roc ) {
228  hsts[roc]->Draw();
229  PixTest::update();
230  }
231  fDisplayedHist = find( fHistList.begin(), fHistList.end(), hsts[nRocs-1] );
232 
233  // test:
234 
235  flags = FLAG_RISING_EDGE;
236  vector<pixel> vpix9 = fApi->getThresholdMap( "Vcal", flags, fParNtrig );
237  for( size_t ipx = 0; ipx < vpix9.size(); ++ipx )
238  LOG(logINFO)
239  << "roc " << setw(2) << (int) vpix9[ipx].roc()
240  << " pix " << setw(2) << (int) vpix9[ipx].column()
241  << " " << setw(2) << (int) vpix9[ipx].row()
242  << " thr " << setw(3) << vpix9[ipx].value();
243 
244  if( col > -1 )
245  fApi->_dut->testPixel( col, row, false );
246 
247 }
std::vector< std::pair< uint8_t, std::vector< pixel > > > getPulseheightVsDAC(std::string dacName, uint8_t dacMin, uint8_t dacMax, uint16_t flags, uint16_t nTriggers)
Definition: api.cc:675
bool setDAC(std::string dacName, uint8_t dacValue, uint8_t rocI2C)
Definition: api.cc:546
void doTest()
function connected to "DoTest" button of PixTab
PixSetup * fPixSetup
all necessary stuff in one place
Definition: PixTest.hh:290
uint8_t getDACRange(std::string dacName)
Definition: api.cc:619
void setToolTips()
implement this to provide updated tool tips if the user changes test parameters
void testPixel(uint8_t column, uint8_t row, bool enable)
Definition: dut.cc:432
std::vector< std::pair< std::string, std::string > > fParameters
the parameters of this test
Definition: PixTest.hh:302
TDirectory * fDirectory
where the root histograms will end up
Definition: PixTest.hh:306
dut * _dut
Definition: api.h:728
std::vector< pixel > getThresholdMap(std::string dacName, uint8_t dacStep, uint8_t dacMin, uint8_t dacMax, uint16_t flags, uint16_t nTriggers)
Definition: api.cc:1080
void testAllPixels(bool enable)
Definition: dut.cc:503
std::list< TH1 * >::iterator fDisplayedHist
pointer to the histogram currently displayed
Definition: PixTest.hh:309
std::list< TH1 * > fHistList
list of histograms available in PixTab::next and PixTab::previous
Definition: PixTest.hh:307
void setTitles(TH1 *h, const char *sx, const char *sy, float size=0.05, float xoff=1.1, float yoff=1.1, float lsize=0.05, int font=42)
utility to set histogram titles
Definition: PixTest.cc:649
void update()
signal to PixTab to update the canvas
Definition: PixTest.cc:569
uint8_t getDAC(size_t rocId, std::string dacName)
Definition: dut.cc:314
pxar::pxarCore * fApi
pointer to the API
Definition: PixTest.hh:289
void init()
sets all test parameters
Definition: PixTest.cc:62
virtual bool setParameter(std::string parName, std::string sval)
set the string value of a parameter