StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
StFtpcMixerMaker.cxx
1 //==================================================================================================
2 // StFtpcMixerMaker: Merges simulated and real data (embedding), needs two input chains
3 //
4 // Author: Frank Simon (fsimon@bnl.gov)
5 //==================================================================================================
6 
8 
21 #include <stdio.h>
22 #include <Stiostream.h>
23 #include "Stiostream.h"
24 #include <ctime> // time functions
25 
26 #include "StFtpcMixerMaker.h"
27 
28 // Message System
29 #include "StMessMgr.h"
30 
31 // SCL
32 #include "St_DataSetIter.h"
33 #include "St_ObjectSet.h"
34 
35 // DataBase
36 #include "StFtpcClusterMaker/StFtpcDbReader.hh"
37 
38 // outPut Data--decoder
39 //VP#include "StDaqLib/GENERIC/EventReader.hh"
40 #include "StSequence.hh"
41 
42 #include "StDAQMaker/StDAQReader.h"
43 
44 // Raw Output stuff
45 #include "StFtpcSequencer.hh"
46 
47 
48 // FTPC tables for SlowSimulator
49 
50 #include "tables/St_fcl_ftpcsqndx_Table.h"
51 #include "tables/St_fcl_ftpcadc_Table.h"
52 
53 
54 ClassImp(StFtpcMixerMaker)
55 
56 
57 StFtpcMixerMaker::StFtpcMixerMaker(const char *name, const char *kind1, const char *kind2):
59 StMaker(name),
60  mConfig1(kind1),
61 mConfig2(kind2)
62 {
63  LOG_INFO << "StFtpcMixerMaker constructed" <<endm;
64 }
65 
66 
69 
70 
71 
73 Int_t StFtpcMixerMaker::InitRun(int RunId) {
74 
75  // get Database
76  // as long as calibrations do not change from event to event InitRun is the right place
77 
78  St_DataSet *ftpc_geometry_db = GetDataBase("Geometry/ftpc");
79  if ( !ftpc_geometry_db ){
80  LOG_WARN << "StFtpcMixerMaker::Error Getting FTPC database: Geometry"<<endm;
81  return kStWarn;
82  }
83  St_DataSetIter dblocal_geometry(ftpc_geometry_db);
84 
85  m_dimensions = (St_ftpcDimensions *)dblocal_geometry("ftpcDimensions");
86  m_padrow_z = (St_ftpcPadrowZ *)dblocal_geometry("ftpcPadrowZ");
87  m_asicmap = (St_ftpcAsicMap *)dblocal_geometry("ftpcAsicMap");
88 
89 
90  St_DataSet *ftpc_calibrations_db = GetDataBase("Calibrations/ftpc");
91  if ( !ftpc_calibrations_db ){
92  LOG_WARN << "StFtpcMixerMaker::Error Getting FTPC database: Calibrations"<<endm;
93  return kStWarn;
94  }
95  St_DataSetIter dblocal_calibrations(ftpc_calibrations_db);
96 
97  m_efield = (St_ftpcEField *)dblocal_calibrations("ftpcEField" );
98  m_vdrift = (St_ftpcVDrift *)dblocal_calibrations("ftpcVDrift" );
99  m_deflection = (St_ftpcDeflection *)dblocal_calibrations("ftpcDeflection" );
100  m_dvdriftdp = (St_ftpcdVDriftdP *)dblocal_calibrations("ftpcdVDriftdP" );
101  m_ddeflectiondp = (St_ftpcdDeflectiondP *)dblocal_calibrations("ftpcdDeflectiondP" );
102  m_ampslope = (St_ftpcAmpSlope *)dblocal_calibrations("ftpcAmpSlope" );
103  m_ampoffset = (St_ftpcAmpOffset *)dblocal_calibrations("ftpcAmpOffset");
104  m_timeoffset = (St_ftpcTimeOffset *)dblocal_calibrations("ftpcTimeOffset");
105  m_driftfield = (St_ftpcDriftField *)dblocal_calibrations("ftpcDriftField");
106  m_gas = (St_ftpcGas *)dblocal_calibrations("ftpcGas");
107  m_electronics = (St_ftpcElectronics *)dblocal_calibrations("ftpcElectronics");
108 
109 
110  return kStOk;
111 }
112 
113 
115 
121 
122  // used to keep track of which FTPC readers are created here
123  Int_t reader1Created = 0;
124  Int_t reader2Created = 0;
125 
126  LOG_INFO << "StFtpcMixerMaker starting..." <<endm;
127 
128  // create FTPC data base reader
129  StFtpcDbReader *dbReader = new StFtpcDbReader(m_dimensions,
130  m_padrow_z);
131 
132 
133  // get the two datasets
134 
135  if(!strcmp(getConfig1(),"daq")) {
136  // DAQ
137  LOG_INFO << "StFtpcMixerMaker: Getting dataset 1 (DAQ)" << endm;
138  St_DataSet *dataset1;
139  dataset1=GetDataSet("Input1"); // get DataSet from chain1
140  if (!dataset1) {
141  LOG_WARN << "StFtpcMixerMaker: Missing Data! Bailing out..." <<endm;
142  assert(dataset1);
143  }
144  daqr1=(StDAQReader*)(dataset1->GetObject()); // get DAQReader
145  if (!daqr1) {
146  LOG_WARN << "StFtpcMixerMaker: Missing Data! Bailing out..." <<endm;
147  assert(daqr1);
148  }
149  ftpcr1=daqr1->getFTPCReader(); // get FTPCReader
150  if (!ftpcr1) {
151  LOG_WARN << "StFtpcMixerMaker: Missing Data! Bailing out..." <<endm;
152  assert(ftpcr1);
153  }
154 
155  if (!ftpcr1->checkForData()) {
156  LOG_WARN << "No FTPC DAQ data available!" << endm;
157  return kStWarn;
158  }
159  else {
160  LOG_INFO << "FTPC DAQ Dataset found!" <<endm;
161  }
162 
163  } else {
164 
165  // FTPC SlowSimulator:
166  LOG_INFO << "StFtpcMixerMaker: Getting dataset 1 (SlowSimulator)" << endm;
167  St_DataSet *dataset1;
168  dataset1=GetDataSet("Input1"); // get DataSet from chain1
169  if (!dataset1) {
170  LOG_WARN << "StFtpcMixerMaker: Missing Data! Bailing out..." <<endm;
171  assert(dataset1);
172  }
173  St_DataSet *raw = GetDataSet("ftpc_raw");
174  if (raw) {
175  LOG_INFO << "FTPC SlowSimulator Dataset found!"<<endm;
176  // FCL
177  St_DataSetIter get(raw);
178 
179  St_fcl_ftpcsqndx *fcl_ftpcsqndx = (St_fcl_ftpcsqndx*)get("fcl_ftpcsqndx");
180  St_fcl_ftpcadc *fcl_ftpcadc = (St_fcl_ftpcadc* )get("fcl_ftpcadc");
181 
182  if (fcl_ftpcsqndx&&fcl_ftpcadc) {
183 
184  ftpcr1=new StFTPCReader((short unsigned int *) fcl_ftpcsqndx->GetTable(),
185  fcl_ftpcsqndx->GetNRows(),
186  (char *) fcl_ftpcadc->GetTable(),
187  fcl_ftpcadc->GetNRows());
188 
189  reader1Created = 1;
190  LOG_INFO << "created StFTPCReader from tables" << endm;
191  }
192  else {
193 
194  LOG_INFO <<"StFtpcMixerMaker: Tables are not found:"
195  << " fcl_ftpcsqndx = " << fcl_ftpcsqndx
196  << " fcl_ftpcadc = " << fcl_ftpcadc << endm;
197  LOG_WARN << "StFtpcMixerMaker exiting... "<<endm;
198  return kStOK;
199  }
200  }
201  else {
202  LOG_WARN << "Error getting FTPC SlowSimulator Dataset!" << endm;
203  }
204 
205  }
206 
207  if(!strcmp(getConfig2(),"daq")) {
208  // DAQ
209  LOG_INFO << "StFtpcMixerMaker: Getting dataset 2 (DAQ)" << endm;
210  St_DataSet *dataset2;
211  dataset2=GetDataSet("Input2"); // get DataSet from chain1
212  if (!dataset2) {
213  LOG_WARN << "StFtpcMixerMaker: Missing Data! Bailing out..." <<endm;
214  assert(dataset2);
215  }
216  daqr2=(StDAQReader*)(dataset2->GetObject()); // get DAQReader
217  if (!daqr2) {
218  LOG_WARN << "StFtpcMixerMaker: Missing Data! Bailing out..." <<endm;
219  assert(daqr2);
220  }
221  ftpcr2=daqr2->getFTPCReader(); // get FTPCReader
222  if (!ftpcr2) {
223  LOG_WARN << "StFtpcMixerMaker: Missing Data! Bailing out..." <<endm;
224  assert(ftpcr2);
225  }
226 
227  if (!ftpcr2->checkForData()) {
228  LOG_WARN << "No FTPC DAQ data available!" << endm;
229  return kStWarn;
230  }
231  else {
232  LOG_INFO << "FTPC DAQ Dataset found!" <<endm;
233  }
234 
235  } else {
236  // FTPC SlowSimulator:
237  LOG_INFO << "StFtpcMixerMaker: Getting dataset 2 (SlowSimulator)" << endm;
238  St_DataSet *dataset2;
239  dataset2=GetDataSet("Input2"); // get DataSet from chain1
240  if (!dataset2) {
241  LOG_WARN << "StFtpcMixerMaker: Missing Data! Bailing out..." <<endm;
242  assert(dataset2);
243  }
244  St_DataSet *raw = GetDataSet("ftpc_raw");
245  if (raw) {
246  LOG_INFO << "FTPC SlowSimulator Dataset found!" << endm;
247  // FCL
248  St_DataSetIter get(raw);
249 
250  St_fcl_ftpcsqndx *fcl_ftpcsqndx = (St_fcl_ftpcsqndx*)get("fcl_ftpcsqndx");
251  St_fcl_ftpcadc *fcl_ftpcadc = (St_fcl_ftpcadc* )get("fcl_ftpcadc");
252 
253 
254 
255  if (fcl_ftpcsqndx&&fcl_ftpcadc) {
256 
257  ftpcr2=new StFTPCReader((short unsigned int *) fcl_ftpcsqndx->GetTable(),
258  fcl_ftpcsqndx->GetNRows(),
259  (char *) fcl_ftpcadc->GetTable(),
260  fcl_ftpcadc->GetNRows());
261 
262  reader2Created = 1;
263  LOG_INFO << "created StFTPCReader from tables" << endm;
264  }
265  else {
266 
267  LOG_INFO <<"StFtpcMixerMaker: Tables are not found:"
268  << " fcl_ftpcsqndx = " << fcl_ftpcsqndx
269  << " fcl_ftpcadc = " << fcl_ftpcadc << endm;
270  LOG_WARN << "StFtpcMixerMaker exiting... " << endm;
271  return kStOK;
272  }
273  }
274  else {
275  LOG_WARN << "Error getting FTPC SlowSimulator Dataset!" << endm;
276  }
277 
278  }
279 
280 
281  // Ok, all the readers are there, now let's take a look at the data of the FTPC
282 
283 
284  // Array to store all ADC values & to do the mixing
285  int *iADC = new int[dbReader->numberOfPadrows()
286  *dbReader->numberOfSectors()
287  *dbReader->numberOfPads()
288  *dbReader->numberOfTimebins()];
289 
290  for (Int_t i = 0; i < dbReader->numberOfPadrows() *dbReader->numberOfSectors() *dbReader->numberOfPads() *dbReader->numberOfTimebins(); i++)
291  iADC[i] = 0;
292 
293 
294  int iRow, iSec, iPad, iHardSec, iHardRow;
295  int firstPadrowToSearch = dbReader->firstPadrowToSearch() - 1;
296 
297  TPCSequence *CurrentSequence;
298 
299  // Let's loop, do the mixing
300 
301  for (iRow = firstPadrowToSearch; iRow < dbReader->numberOfPadrows(); iRow++) // loop over Rows
302  {
303  for (iSec = dbReader->firstSectorToSearch()-1; iSec < dbReader->lastSectorToSearch(); iSec++) // loop over Sectors
304  {
305  // get hardware (daq) numbering
306  iHardSec = dbReader->numberOfSectors() * (int)(iRow/2) + iSec + 1;
307  iHardRow = iRow%2 + 1;
308 
309  // get occupied pads for first reader
310  unsigned char *(padlist[2]);
311  int iOccPads = ftpcr1->getPadList(iHardSec, iHardRow, padlist[iHardRow-1]);
312 
313  int numberOfSequences;
314 
315 
316  for (int iPadCounter = 0; iPadCounter < iOccPads; iPadCounter++) // Loop over occupied pads
317  {
318  iPad = padlist[iHardRow-1][iPadCounter]; // carefull: This pad counts from 1 to 160!
319 
320  ftpcr1->getSequences(iHardSec, iHardRow, iPad, &numberOfSequences, CurrentSequence);
321 
322  iPad--; // iADC array counts pads from 0 to 159, important for StFtpcSequencer
323 
324  for (int iSeqCounter = 0; iSeqCounter < numberOfSequences; iSeqCounter++) // Loop over Sequences in Pad
325  {
326  for (int iEntry = 0; iEntry < CurrentSequence[iSeqCounter].Length; iEntry++) // Loop over entries in sequence
327  {
328  iADC[iRow*dbReader->numberOfSectors() *dbReader->numberOfPads() *dbReader->numberOfTimebins()
329  + iSec*dbReader->numberOfPads() *dbReader->numberOfTimebins()
330  + iPad * dbReader->numberOfTimebins() + CurrentSequence[iSeqCounter].startTimeBin + iEntry]
331  += (int)CurrentSequence[iSeqCounter].FirstAdc[iEntry];
332  }
333  } // finish Loop over Sequences in Pad
334  } // finish Loop over occupied pads
335 
336  iOccPads = ftpcr2->getPadList(iHardSec, iHardRow, padlist[iHardRow-1]);
337 
338  for (int iPadCounter = 0; iPadCounter < iOccPads; iPadCounter++) // Loop over occupied pads
339  {
340  iPad = padlist[iHardRow-1][iPadCounter]; // carefull: This pad counts from 1 to 160!
341 
342  ftpcr2->getSequences(iHardSec, iHardRow, iPad, &numberOfSequences, CurrentSequence);
343 
344  iPad--; // iADC array counts pads from 0 to 159, important for StFtpcSequencer
345 
346  for (int iSeqCounter = 0; iSeqCounter < numberOfSequences; iSeqCounter++) // Loop over Sequences in Pad
347  {
348  for (int iEntry = 0; iEntry < CurrentSequence[iSeqCounter].Length; iEntry++) // Loop over entries in sequence
349  {
350  iADC[iRow*dbReader->numberOfSectors() *dbReader->numberOfPads() *dbReader->numberOfTimebins()
351  + iSec*dbReader->numberOfPads() *dbReader->numberOfTimebins()
352  + iPad * dbReader->numberOfTimebins() + CurrentSequence[iSeqCounter].startTimeBin + iEntry]
353  += (int)CurrentSequence[iSeqCounter].FirstAdc[iEntry];
354 
355  }
356 
357  } // finish Loop over Sequences in Pad
358  } // finish Loop over occupied pads
359 
360  } // finish Loop over Sectors
361  } // finish Loop over Rows
362 
363 
364  // Embedding done, now write out the Sequences
365 
366  LOG_INFO << "FTPC Embedding done... "<< endm;
367 
368  St_DataSetIter local((GetData()));
369 
370  St_fcl_ftpcndx *fcl_ftpcndx_out = new St_fcl_ftpcndx("fcl_ftpcndx", 2);
371  local.Add(fcl_ftpcndx_out);
372  St_fcl_ftpcsqndx *fcl_ftpcsqndx_out = new St_fcl_ftpcsqndx("fcl_ftpcsqndx", 500000);
373  local.Add(fcl_ftpcsqndx_out);
374  St_fcl_ftpcadc *fcl_ftpcadc_out = new St_fcl_ftpcadc("fcl_ftpcadc", 2000000);
375  local.Add(fcl_ftpcadc_out);
376 
377  LOG_INFO << "Output sequences created..." << endm;
378 
379  // create StFtpcSequencer to write ADC data to sequences (zero suppressed)
380  StFtpcSequencer *ftpcSequencer = new StFtpcSequencer(fcl_ftpcndx_out, fcl_ftpcsqndx_out, fcl_ftpcadc_out);
381 
382  ftpcSequencer->writeArray(iADC, dbReader->numberOfPadrows(), dbReader->numberOfSectors(), dbReader->numberOfPads(), dbReader->numberOfTimebins());
383 
384 
385  delete ftpcSequencer;
386  delete dbReader;
387  if (reader1Created) delete ftpcr1;
388  if (reader2Created) delete ftpcr2;
389  delete[] iADC;
390 
391  LOG_INFO << "FtpcMixerMaker done..." << endm;
392 
393  return kStOK;
394 } // Make()
395 
396 
398 void StFtpcMixerMaker::Clear(Option_t *opt)
399 {
400  StMaker::Clear();
401 }
402 
403 
406 {
407  return kStOK;
408 }
409 
410  /***************************************************************************
411  *
412  * $Id: StFtpcMixerMaker.cxx,v 1.8 2017/04/26 19:48:41 perev Exp $
413  *
414  * $Log: StFtpcMixerMaker.cxx,v $
415  * Revision 1.8 2017/04/26 19:48:41 perev
416  * Hide m_DataSet
417  *
418  * Revision 1.7 2007/04/28 17:56:11 perev
419  * Redundant StChain.h removed
420  *
421  * Revision 1.6 2007/01/15 15:02:12 jcs
422  * replace printf, cout and gMesMgr with Logger
423  *
424  * Revision 1.5 2004/03/04 15:49:00 jcs
425  * remove unnecessary fcl_fppoint include
426  *
427  * Revision 1.4 2003/09/02 17:58:15 perev
428  * gcc 3.2 updates + WarnOff
429  *
430  * Revision 1.3 2003/07/18 18:31:47 perev
431  * test for nonexistance of XXXReader added
432  *
433  * Revision 1.2 2003/06/13 12:12:44 jcs
434  * use the same StFtpcDbReader constructor as used by Sti/StFtpcDetectorBuilder
435  *
436  * Revision 1.1 2003/02/14 18:11:25 fsimon
437  * Initial commit of FTPC embedding code
438  *
439  *
440  ***************************************************************************/
Int_t InitRun(int)
InitRun method, reads in the database (FTPC geometry and calibration)
~StFtpcMixerMaker()
Destructor.
int writeArray(const int *cArray, const int numberPadrows, const int numberSectors, const int numberPads, const int numberTimebins)
writeArray method, fills zero-suppressed sequences from ADC array
Int_t Make()
Make method, does all the work.
virtual void Clear(Option_t *option="")
User defined functions.
Definition: StMaker.cxx:634
FTPC Sequencer, creates Simulator sequences from ADC values.
FTPC Mixer Maker, main part of FTPC embedding Framework.
virtual TObject * GetObject() const
The depricated method (left here for the sake of the backward compatibility)
Definition: TDataSet.cxx:428
Definition: Stypes.h:42
void Clear(Option_t *option="")
Clear method.
Definition: Stypes.h:40
Int_t Finish()
Finish method.
Definition: Stypes.h:41