00001
00002
00003
00004
00005
00006
00008
00021 #include <stdio.h>
00022 #include <Stiostream.h>
00023 #include "Stiostream.h"
00024 #include <ctime>
00025
00026 #include "StFtpcMixerMaker.h"
00027
00028
00029 #include "StMessMgr.h"
00030
00031
00032 #include "St_DataSetIter.h"
00033 #include "St_ObjectSet.h"
00034
00035
00036 #include "StFtpcClusterMaker/StFtpcDbReader.hh"
00037
00038
00039
00040 #include "StSequence.hh"
00041
00042 #include "StDAQMaker/StDAQReader.h"
00043
00044
00045 #include "StFtpcSequencer.hh"
00046
00047
00048
00049
00050 #include "tables/St_fcl_ftpcsqndx_Table.h"
00051 #include "tables/St_fcl_ftpcadc_Table.h"
00052
00053
00054 ClassImp(StFtpcMixerMaker)
00055
00056
00057
00058 StFtpcMixerMaker::StFtpcMixerMaker(const char *name, const char *kind1, const char *kind2):
00059 StMaker(name),
00060 mConfig1(kind1),
00061 mConfig2(kind2)
00062 {
00063 LOG_INFO << "StFtpcMixerMaker constructed" <<endm;
00064 }
00065
00066
00068 StFtpcMixerMaker::~StFtpcMixerMaker() { }
00069
00070
00071
00073 Int_t StFtpcMixerMaker::InitRun(int RunId) {
00074
00075
00076
00077
00078 St_DataSet *ftpc_geometry_db = GetDataBase("Geometry/ftpc");
00079 if ( !ftpc_geometry_db ){
00080 LOG_WARN << "StFtpcMixerMaker::Error Getting FTPC database: Geometry"<<endm;
00081 return kStWarn;
00082 }
00083 St_DataSetIter dblocal_geometry(ftpc_geometry_db);
00084
00085 m_dimensions = (St_ftpcDimensions *)dblocal_geometry("ftpcDimensions");
00086 m_padrow_z = (St_ftpcPadrowZ *)dblocal_geometry("ftpcPadrowZ");
00087 m_asicmap = (St_ftpcAsicMap *)dblocal_geometry("ftpcAsicMap");
00088
00089
00090 St_DataSet *ftpc_calibrations_db = GetDataBase("Calibrations/ftpc");
00091 if ( !ftpc_calibrations_db ){
00092 LOG_WARN << "StFtpcMixerMaker::Error Getting FTPC database: Calibrations"<<endm;
00093 return kStWarn;
00094 }
00095 St_DataSetIter dblocal_calibrations(ftpc_calibrations_db);
00096
00097 m_efield = (St_ftpcEField *)dblocal_calibrations("ftpcEField" );
00098 m_vdrift = (St_ftpcVDrift *)dblocal_calibrations("ftpcVDrift" );
00099 m_deflection = (St_ftpcDeflection *)dblocal_calibrations("ftpcDeflection" );
00100 m_dvdriftdp = (St_ftpcdVDriftdP *)dblocal_calibrations("ftpcdVDriftdP" );
00101 m_ddeflectiondp = (St_ftpcdDeflectiondP *)dblocal_calibrations("ftpcdDeflectiondP" );
00102 m_ampslope = (St_ftpcAmpSlope *)dblocal_calibrations("ftpcAmpSlope" );
00103 m_ampoffset = (St_ftpcAmpOffset *)dblocal_calibrations("ftpcAmpOffset");
00104 m_timeoffset = (St_ftpcTimeOffset *)dblocal_calibrations("ftpcTimeOffset");
00105 m_driftfield = (St_ftpcDriftField *)dblocal_calibrations("ftpcDriftField");
00106 m_gas = (St_ftpcGas *)dblocal_calibrations("ftpcGas");
00107 m_electronics = (St_ftpcElectronics *)dblocal_calibrations("ftpcElectronics");
00108
00109
00110 return kStOk;
00111 }
00112
00113
00115
00120 Int_t StFtpcMixerMaker::Make() {
00121
00122
00123 Int_t reader1Created = 0;
00124 Int_t reader2Created = 0;
00125
00126 LOG_INFO << "StFtpcMixerMaker starting..." <<endm;
00127
00128
00129 StFtpcDbReader *dbReader = new StFtpcDbReader(m_dimensions,
00130 m_padrow_z);
00131
00132
00133
00134
00135 if(!strcmp(getConfig1(),"daq")) {
00136
00137 LOG_INFO << "StFtpcMixerMaker: Getting dataset 1 (DAQ)" << endm;
00138 St_DataSet *dataset1;
00139 dataset1=GetDataSet("Input1");
00140 if (!dataset1) {
00141 LOG_WARN << "StFtpcMixerMaker: Missing Data! Bailing out..." <<endm;
00142 assert(dataset1);
00143 }
00144 daqr1=(StDAQReader*)(dataset1->GetObject());
00145 if (!daqr1) {
00146 LOG_WARN << "StFtpcMixerMaker: Missing Data! Bailing out..." <<endm;
00147 assert(daqr1);
00148 }
00149 ftpcr1=daqr1->getFTPCReader();
00150 if (!ftpcr1) {
00151 LOG_WARN << "StFtpcMixerMaker: Missing Data! Bailing out..." <<endm;
00152 assert(ftpcr1);
00153 }
00154
00155 if (!ftpcr1->checkForData()) {
00156 LOG_WARN << "No FTPC DAQ data available!" << endm;
00157 return kStWarn;
00158 }
00159 else {
00160 LOG_INFO << "FTPC DAQ Dataset found!" <<endm;
00161 }
00162
00163 } else {
00164
00165
00166 LOG_INFO << "StFtpcMixerMaker: Getting dataset 1 (SlowSimulator)" << endm;
00167 St_DataSet *dataset1;
00168 dataset1=GetDataSet("Input1");
00169 if (!dataset1) {
00170 LOG_WARN << "StFtpcMixerMaker: Missing Data! Bailing out..." <<endm;
00171 assert(dataset1);
00172 }
00173 St_DataSet *raw = GetDataSet("ftpc_raw");
00174 if (raw) {
00175 LOG_INFO << "FTPC SlowSimulator Dataset found!"<<endm;
00176
00177 St_DataSetIter get(raw);
00178
00179 St_fcl_ftpcsqndx *fcl_ftpcsqndx = (St_fcl_ftpcsqndx*)get("fcl_ftpcsqndx");
00180 St_fcl_ftpcadc *fcl_ftpcadc = (St_fcl_ftpcadc* )get("fcl_ftpcadc");
00181
00182 if (fcl_ftpcsqndx&&fcl_ftpcadc) {
00183
00184 ftpcr1=new StFTPCReader((short unsigned int *) fcl_ftpcsqndx->GetTable(),
00185 fcl_ftpcsqndx->GetNRows(),
00186 (char *) fcl_ftpcadc->GetTable(),
00187 fcl_ftpcadc->GetNRows());
00188
00189 reader1Created = 1;
00190 LOG_INFO << "created StFTPCReader from tables" << endm;
00191 }
00192 else {
00193
00194 LOG_INFO <<"StFtpcMixerMaker: Tables are not found:"
00195 << " fcl_ftpcsqndx = " << fcl_ftpcsqndx
00196 << " fcl_ftpcadc = " << fcl_ftpcadc << endm;
00197 LOG_WARN << "StFtpcMixerMaker exiting... "<<endm;
00198 return kStOK;
00199 }
00200 }
00201 else {
00202 LOG_WARN << "Error getting FTPC SlowSimulator Dataset!" << endm;
00203 }
00204
00205 }
00206
00207 if(!strcmp(getConfig2(),"daq")) {
00208
00209 LOG_INFO << "StFtpcMixerMaker: Getting dataset 2 (DAQ)" << endm;
00210 St_DataSet *dataset2;
00211 dataset2=GetDataSet("Input2");
00212 if (!dataset2) {
00213 LOG_WARN << "StFtpcMixerMaker: Missing Data! Bailing out..." <<endm;
00214 assert(dataset2);
00215 }
00216 daqr2=(StDAQReader*)(dataset2->GetObject());
00217 if (!daqr2) {
00218 LOG_WARN << "StFtpcMixerMaker: Missing Data! Bailing out..." <<endm;
00219 assert(daqr2);
00220 }
00221 ftpcr2=daqr2->getFTPCReader();
00222 if (!ftpcr2) {
00223 LOG_WARN << "StFtpcMixerMaker: Missing Data! Bailing out..." <<endm;
00224 assert(ftpcr2);
00225 }
00226
00227 if (!ftpcr2->checkForData()) {
00228 LOG_WARN << "No FTPC DAQ data available!" << endm;
00229 return kStWarn;
00230 }
00231 else {
00232 LOG_INFO << "FTPC DAQ Dataset found!" <<endm;
00233 }
00234
00235 } else {
00236
00237 LOG_INFO << "StFtpcMixerMaker: Getting dataset 2 (SlowSimulator)" << endm;
00238 St_DataSet *dataset2;
00239 dataset2=GetDataSet("Input2");
00240 if (!dataset2) {
00241 LOG_WARN << "StFtpcMixerMaker: Missing Data! Bailing out..." <<endm;
00242 assert(dataset2);
00243 }
00244 St_DataSet *raw = GetDataSet("ftpc_raw");
00245 if (raw) {
00246 LOG_INFO << "FTPC SlowSimulator Dataset found!" << endm;
00247
00248 St_DataSetIter get(raw);
00249
00250 St_fcl_ftpcsqndx *fcl_ftpcsqndx = (St_fcl_ftpcsqndx*)get("fcl_ftpcsqndx");
00251 St_fcl_ftpcadc *fcl_ftpcadc = (St_fcl_ftpcadc* )get("fcl_ftpcadc");
00252
00253
00254
00255 if (fcl_ftpcsqndx&&fcl_ftpcadc) {
00256
00257 ftpcr2=new StFTPCReader((short unsigned int *) fcl_ftpcsqndx->GetTable(),
00258 fcl_ftpcsqndx->GetNRows(),
00259 (char *) fcl_ftpcadc->GetTable(),
00260 fcl_ftpcadc->GetNRows());
00261
00262 reader2Created = 1;
00263 LOG_INFO << "created StFTPCReader from tables" << endm;
00264 }
00265 else {
00266
00267 LOG_INFO <<"StFtpcMixerMaker: Tables are not found:"
00268 << " fcl_ftpcsqndx = " << fcl_ftpcsqndx
00269 << " fcl_ftpcadc = " << fcl_ftpcadc << endm;
00270 LOG_WARN << "StFtpcMixerMaker exiting... " << endm;
00271 return kStOK;
00272 }
00273 }
00274 else {
00275 LOG_WARN << "Error getting FTPC SlowSimulator Dataset!" << endm;
00276 }
00277
00278 }
00279
00280
00281
00282
00283
00284
00285 int *iADC = new int[dbReader->numberOfPadrows()
00286 *dbReader->numberOfSectors()
00287 *dbReader->numberOfPads()
00288 *dbReader->numberOfTimebins()];
00289
00290 for (Int_t i = 0; i < dbReader->numberOfPadrows() *dbReader->numberOfSectors() *dbReader->numberOfPads() *dbReader->numberOfTimebins(); i++)
00291 iADC[i] = 0;
00292
00293
00294 int iRow, iSec, iPad, iHardSec, iHardRow;
00295 int firstPadrowToSearch = dbReader->firstPadrowToSearch() - 1;
00296
00297 TPCSequence *CurrentSequence;
00298
00299
00300
00301 for (iRow = firstPadrowToSearch; iRow < dbReader->numberOfPadrows(); iRow++)
00302 {
00303 for (iSec = dbReader->firstSectorToSearch()-1; iSec < dbReader->lastSectorToSearch(); iSec++)
00304 {
00305
00306 iHardSec = dbReader->numberOfSectors() * (int)(iRow/2) + iSec + 1;
00307 iHardRow = iRow%2 + 1;
00308
00309
00310 unsigned char *(padlist[2]);
00311 int iOccPads = ftpcr1->getPadList(iHardSec, iHardRow, padlist[iHardRow-1]);
00312
00313 int numberOfSequences;
00314
00315
00316 for (int iPadCounter = 0; iPadCounter < iOccPads; iPadCounter++)
00317 {
00318 iPad = padlist[iHardRow-1][iPadCounter];
00319
00320 ftpcr1->getSequences(iHardSec, iHardRow, iPad, &numberOfSequences, CurrentSequence);
00321
00322 iPad--;
00323
00324 for (int iSeqCounter = 0; iSeqCounter < numberOfSequences; iSeqCounter++)
00325 {
00326 for (int iEntry = 0; iEntry < CurrentSequence[iSeqCounter].Length; iEntry++)
00327 {
00328 iADC[iRow*dbReader->numberOfSectors() *dbReader->numberOfPads() *dbReader->numberOfTimebins()
00329 + iSec*dbReader->numberOfPads() *dbReader->numberOfTimebins()
00330 + iPad * dbReader->numberOfTimebins() + CurrentSequence[iSeqCounter].startTimeBin + iEntry]
00331 += (int)CurrentSequence[iSeqCounter].FirstAdc[iEntry];
00332 }
00333 }
00334 }
00335
00336 iOccPads = ftpcr2->getPadList(iHardSec, iHardRow, padlist[iHardRow-1]);
00337
00338 for (int iPadCounter = 0; iPadCounter < iOccPads; iPadCounter++)
00339 {
00340 iPad = padlist[iHardRow-1][iPadCounter];
00341
00342 ftpcr2->getSequences(iHardSec, iHardRow, iPad, &numberOfSequences, CurrentSequence);
00343
00344 iPad--;
00345
00346 for (int iSeqCounter = 0; iSeqCounter < numberOfSequences; iSeqCounter++)
00347 {
00348 for (int iEntry = 0; iEntry < CurrentSequence[iSeqCounter].Length; iEntry++)
00349 {
00350 iADC[iRow*dbReader->numberOfSectors() *dbReader->numberOfPads() *dbReader->numberOfTimebins()
00351 + iSec*dbReader->numberOfPads() *dbReader->numberOfTimebins()
00352 + iPad * dbReader->numberOfTimebins() + CurrentSequence[iSeqCounter].startTimeBin + iEntry]
00353 += (int)CurrentSequence[iSeqCounter].FirstAdc[iEntry];
00354
00355 }
00356
00357 }
00358 }
00359
00360 }
00361 }
00362
00363
00364
00365
00366 LOG_INFO << "FTPC Embedding done... "<< endm;
00367
00368 St_DataSetIter local(m_DataSet);
00369
00370 St_fcl_ftpcndx *fcl_ftpcndx_out = new St_fcl_ftpcndx("fcl_ftpcndx", 2);
00371 local.Add(fcl_ftpcndx_out);
00372 St_fcl_ftpcsqndx *fcl_ftpcsqndx_out = new St_fcl_ftpcsqndx("fcl_ftpcsqndx", 500000);
00373 local.Add(fcl_ftpcsqndx_out);
00374 St_fcl_ftpcadc *fcl_ftpcadc_out = new St_fcl_ftpcadc("fcl_ftpcadc", 2000000);
00375 local.Add(fcl_ftpcadc_out);
00376
00377 LOG_INFO << "Output sequences created..." << endm;
00378
00379
00380 StFtpcSequencer *ftpcSequencer = new StFtpcSequencer(fcl_ftpcndx_out, fcl_ftpcsqndx_out, fcl_ftpcadc_out);
00381
00382 ftpcSequencer->writeArray(iADC, dbReader->numberOfPadrows(), dbReader->numberOfSectors(), dbReader->numberOfPads(), dbReader->numberOfTimebins());
00383
00384
00385 delete ftpcSequencer;
00386 delete dbReader;
00387 if (reader1Created) delete ftpcr1;
00388 if (reader2Created) delete ftpcr2;
00389 delete[] iADC;
00390
00391 LOG_INFO << "FtpcMixerMaker done..." << endm;
00392
00393 return kStOK;
00394 }
00395
00396
00398 void StFtpcMixerMaker::Clear(Option_t *opt)
00399 {
00400 StMaker::Clear();
00401 }
00402
00403
00405 Int_t StFtpcMixerMaker::Finish()
00406 {
00407 return kStOK;
00408 }
00409
00410
00411
00412
00413
00414
00415
00416
00417
00418
00419
00420
00421
00422
00423
00424
00425
00426
00427
00428
00429
00430
00431
00432
00433
00434
00435
00436
00437