00001
00002
00003
00004
00005
00006
00007
00008
00009
00011
00012 #include "Stiostream.h"
00013 #include "StVertexSeedMaker.h"
00014 #include "StDbLib/StDbManager.hh"
00015 #include "StDbLib/StDbConfigNode.hh"
00016 #include "StDbLib/StDataBaseI.hh"
00017 #include "StDetectorDbMaker/StDetectorDbTriggerID.h"
00018 #include "StMessMgr.h"
00019 #include "StVertexId.h"
00020 #include "tables/St_vertexSeed_Table.h"
00021 #include "tables/St_vertexSeedTriggers_Table.h"
00022 #include "tables/St_beamInfo_Table.h"
00023 #include "TSystem.h"
00024 #include "TFile.h"
00025 #include "TVirtualFitter.h"
00026 #include "TNtuple.h"
00027 #include "TEventList.h"
00028 #include "TArrayF.h"
00029 #include "StTree.h"
00030 #include "TMath.h"
00031
00032
00033
00034 static TArrayF xVert, yVert, zVert, multA, exVert, eyVert;
00035 int nverts,nsize;
00036 double beamWidth;
00037 Double_t funcX(float z,Double_t *par) {
00038 Double_t x = par[0] + par[1]*z;
00039 return x;
00040 }
00041 Double_t funcY(float z,Double_t *par) {
00042 Double_t y = par[2] + par[3]*z;
00043 return y;
00044 }
00045 void fnch(Int_t &npar, Double_t *gin, Double_t &f, Double_t *par, Int_t iflag) {
00046
00047 double chisq = 0;
00048 for (int i=0;i<nverts; i++) {
00049 if (exVert[i]==0) {
00050
00051 double delta_sq, error_sq;
00052
00053
00054
00055
00056
00057
00058
00059
00060 error_sq = 0.0016 + (2.45 / multA[i]);
00061 delta_sq = TMath::Power(xVert[i]-funcX(zVert[i],par),2) +
00062 TMath::Power(yVert[i]-funcY(zVert[i],par),2);
00063 chisq += (delta_sq/error_sq);
00064
00065 } else {
00066
00067
00068
00069
00070 double errx2 = exVert[i]*exVert[i] + beamWidth*beamWidth;
00071 double erry2 = eyVert[i]*eyVert[i] + beamWidth*beamWidth;
00072 chisq += TMath::Power(xVert[i]-funcX(zVert[i],par),2)/errx2 +
00073 TMath::Power(yVert[i]-funcY(zVert[i],par),2)/erry2;
00074
00075 }
00076 }
00077 f = chisq;
00078 }
00079 void setArraySize(int n) {
00080 xVert.Set(n);
00081 yVert.Set(n);
00082 zVert.Set(n);
00083 multA.Set(n);
00084 exVert.Set(n);
00085 eyVert.Set(n);
00086 nsize = n;
00087 }
00088 void addVert(float x, float y, float z, float m, float ex, float ey) {
00089 if (nverts >= nsize) setArraySize(nsize*2);
00090 xVert[nverts] = x;
00091 yVert[nverts] = y;
00092 zVert[nverts] = z;
00093 multA[nverts] = m;
00094 exVert[nverts] = ex;
00095 eyVert[nverts] = ey;
00096 nverts++;
00097 }
00098
00099
00100
00101 ClassImp(StVertexSeedMaker)
00102
00103 StVertexSeedMaker::StVertexSeedMaker(const char *name,
00104 const char* defaultDir):StMaker(name){
00105 xdist = new TH1F("xdist","xdist",1000,HIST_MIN,HIST_MAX);
00106 ydist = new TH1F("ydist","ydist",1000,HIST_MIN,HIST_MAX);
00107 xerr = new TH1F("xerr","x measured - x guess",1000,HIST_MIN,HIST_MAX);
00108 yerr = new TH1F("yerr","y measured - y guess",1000,HIST_MIN,HIST_MAX);
00109 defDir = defaultDir;
00110 mTempOut = 0;
00111 resNtuple = 0;
00112 nsize = 0;
00113 setArraySize(512);
00114 UseEventDateTime();
00115 useAllTriggers = kFALSE;
00116 dbTriggersTable = 0;
00117 Reset();
00118 }
00119
00120 void StVertexSeedMaker::Reset() {
00121 minEntries = 100;
00122 maxX0Err = 0.05;
00123 maxY0Err = 0.05;
00124 mHistOut=kTRUE;
00125 zVertexMax = 100.0;
00126 zVertexMin = -100.0;
00127 r2VertexMax = 15.0;
00128 nverts = 0;
00129 Clear("");
00130 xguess = 0;
00131 yguess = 0;
00132 HIST_MIN = -3.0;
00133 HIST_MAX = 3.0;
00134 xdist->Reset();
00135 ydist->Reset();
00136 xerr->Reset();
00137 yerr->Reset();
00138
00139 if (mTempOut) delete mTempOut;
00140 else if (resNtuple) delete resNtuple;
00141 mTempOut = new TFile(Form("%s/vertexseedhist.%d.root",
00142 gSystem->TempDirectory(),
00143 gSystem->GetPid()),"RECREATE");
00144 resNtuple = new TNtuple("resNtuple","resNtuple","event:x:y:z:mult:trig:run:fill:zdc:rank:itpc:otpc:detmap:ex:ey");
00145 LOG_INFO << "Opening new temp file at " << mTempOut->GetName() << endm;
00146
00147 date = 0;
00148 time = 0;
00149 fill = -1;
00150 run = -1;
00151 zdc = -1;
00152 rank = 0;
00153 itpc = 0;
00154 otpc = 0;
00155 detmap = 0;
00156 a[0] = -888.0;
00157
00158 a[1] = 0.0;
00159 a[2] = 0.0;
00160 a[3] = 0.0;
00161 chi = 0.0;
00162 }
00163
00164 StVertexSeedMaker::~StVertexSeedMaker(){
00165 }
00166
00167 Int_t StVertexSeedMaker::Init(){
00168 AddHist(xdist);
00169 AddHist(ydist);
00170 AddHist(xerr);
00171 AddHist(yerr);
00172 return StMaker::Init();
00173 }
00174
00175 void StVertexSeedMaker::Clear(Option_t *option){
00176 xguess = 0;
00177 yguess = 0;
00178 zvertex = -999.0;
00179 yvertex = -999.0;
00180 xvertex = -999.0;
00181 exvertex = 0;
00182 eyvertex = 0;
00183 }
00184
00185 Int_t StVertexSeedMaker::Make(){
00186 if (date==0) FillDateTime();
00187
00188
00189 if (a[0] == -888) {
00190 LOG_INFO << "Reading db values at the start..." << endm;
00191 Int_t status = FillAssumed();
00192 if (status != kStOk) return status;
00193 status = GetVertexSeedTriggers();
00194 if (status != kStOk) return status;
00195 }
00196
00197 if (CheckTriggers()) {
00198 LOG_INFO << "event does not satisfy triggers" << endm;
00199 return kStOk;
00200 }
00201
00202 Int_t eventResult = GetEventData();
00203 if (eventResult != kStOk) return eventResult;
00204
00205 if (zvertex<-998) {
00206 LOG_INFO << "No primary vertex" << endm;
00207 return kStOk;
00208 }
00209
00210
00211 xguess = a[0] + (a[1] * zvertex);
00212 yguess = a[2] + (a[3] * zvertex);
00213 LOG_INFO << "x guess = " << xguess << endm;
00214 LOG_INFO << "y guess = " << yguess << endm;
00215
00216
00217 float r2vertex = xvertex*xvertex + yvertex*yvertex;
00218 if ((zvertex > zVertexMin) && (zvertex < zVertexMax) &&
00219 (mult >= 5) && (r2vertex < r2VertexMax)){
00220 xdist->Fill(xvertex);
00221 xerr ->Fill(xvertex-xguess);
00222 ydist->Fill(yvertex);
00223 yerr ->Fill(yvertex-yguess);
00224 eventNumber = (float)GetEventNumber();
00225 resNtuple->Fill(eventNumber,xvertex,yvertex,zvertex,mult,trig,
00226 (float) run,(float) fill,zdc,rank,
00227 (float) itpc,(float) otpc,(float) detmap,exvertex,eyvertex);
00228 addVert(xvertex,yvertex,zvertex,mult,exvertex,eyvertex);
00229 }
00230
00231 return kStOk;
00232 }
00233
00234 Int_t StVertexSeedMaker::Finish() {
00235 FindResult();
00236 return StMaker::Finish();
00237 }
00238
00239 Bool_t StVertexSeedMaker::ValidTrigger(unsigned int tid) {
00240
00241 if (!dbTriggersTable) return kTRUE;
00242 vertexSeedTriggers_st* trigsTable = dbTriggersTable->GetTable();
00243 Int_t nTrigs = (Int_t) dbTriggersTable->GetNRows();
00244 for (Int_t i = 0; i < nTrigs; i++, trigsTable++) {
00245 unsigned int dbTid = trigsTable->trigid;
00246 if (useAllTriggers || dbTid == 9999999 || (dbTid > 0 && tid == dbTid)) {
00247 trig = (float) tid;
00248 return kTRUE;
00249 }
00250 }
00251 return kFALSE;
00252 }
00253
00254 void StVertexSeedMaker::FindResult(Bool_t checkDb) {
00255 Bool_t writeIt = kFALSE;
00256 if (nverts >= minEntries){
00257 FitData();
00258 if (ep[0] > maxX0Err){
00259 LOG_ERROR << "x unstable. x0 error = " << ep[0] << " cm." << endm;
00260 }
00261 if (ep[2] > maxY0Err){
00262 LOG_ERROR << "y unstable. y0 error = " << ep[2] << " cm." << endm;
00263 }
00264 if ((ep[0] <= maxX0Err) && (ep[2] <= maxY0Err)) {
00265 if (checkDb) {
00266
00267
00268 LOG_INFO << "Reading db values at the end..." << endm;
00269 Int_t status = FillAssumed();
00270 if (status == kStOk) {
00271 if (ChangedValues() || BetterErrors()) writeIt = kTRUE;
00272 else { LOG_INFO << "Values have not changed/improved." << endm; }
00273 } else {
00274 LOG_WARN << "Could not get db values." << endm;
00275 writeIt = kTRUE;
00276 }
00277 } else {
00278 writeIt = kTRUE;
00279 }
00280 }
00281 } else {
00282 LOG_ERROR << "Insufficient stats for " <<
00283 "mean vertex determination.\n Only " << nverts << " entries." << endm;
00284 }
00285
00286 if (writeIt) WriteTableToFile();
00287 else { LOG_WARN << "Not writing table!!!!!" << endm; }
00288
00289 if (mHistOut) WriteHistFile();
00290 }
00291
00292 void StVertexSeedMaker::PrintInfo() {
00293 LOG_INFO << "\n**************************************************************"
00294 << "\n* $Id: StVertexSeedMaker.cxx,v 1.47 2012/02/29 02:00:04 genevb Exp $"
00295 << "\n**************************************************************" << endm;
00296
00297 if (Debug()) StMaker::PrintInfo();
00298 }
00299
00300 void StVertexSeedMaker::WriteTableToFile(){
00301 char filename[80];
00302 if (defDir.Length()>0 && !defDir.EndsWith("/")) defDir.Append("/");
00303 sprintf(filename,"%svertexSeed.%08d.%06d.C",defDir.Data(),date,time);
00304 LOG_INFO << "Writing new table to:\n "
00305 << filename << endm;
00306 TString dirname = gSystem->DirName(filename);
00307 if (gSystem->OpenDirectory(dirname.Data())==0) {
00308 if (gSystem->mkdir(dirname.Data())) {
00309 LOG_WARN << "Directory creation failed for:\n " << dirname
00310 << "\n Putting table file in current directory" << endm;
00311 for (int i=0;i<80;i++){filename[i]=0;}
00312 sprintf(filename,"vertexSeed.%08d.%06d.C",date,time);
00313 }
00314 }
00315 ofstream *out = new ofstream(filename);
00316 VertexSeedTable()->SavePrimitive(*out,"");
00317 return;
00318 }
00319
00320 St_vertexSeed* StVertexSeedMaker::VertexSeedTable(){
00321 St_vertexSeed* table = new St_vertexSeed("vertexSeed",1);
00322 vertexSeed_st* row = table->GetTable();
00323 row->x0 = p[0];
00324 row->dxdz = p[1];
00325 row->y0 = p[2];
00326 row->dydz = p[3];
00327 row->err_x0 = ep[0];
00328 row->err_dxdz = ep[1];
00329 row->err_y0 = ep[2];
00330 row->err_dydz = ep[3];
00331 row->chisq_dof = chi;
00332 row->weight = 100.0;
00333 row->stats = (float) nverts;
00334 table->SetNRows(1);
00335 return table;
00336 }
00337
00338 void StVertexSeedMaker::WriteHistFile(){
00339 if (resNtuple->GetEntries() == 0) {
00340 LOG_INFO << "Not writing histograms - no entries!!!" << endm;
00341 return;
00342 }
00343
00344 if (defDir.Length()>0 && !defDir.EndsWith("/")) defDir.Append("/");
00345 TString fileNameBase = Form("vertexseedhist.%08d.%06d.ROOT",date,time);
00346 TString fileName = defDir;
00347 fileName += fileNameBase;
00348 LOG_INFO << "Writing new histograms to:\n "
00349 << fileName << endm;
00350 TString dirname = gSystem->DirName(fileName.Data());
00351 if (gSystem->OpenDirectory(dirname.Data())==0) {
00352 if (gSystem->mkdir(dirname.Data())) {
00353 LOG_WARN << "Directory creation failed for:\n " << dirname
00354 << "\n Putting histogram file in current directory" << endm;
00355 fileName = fileNameBase;
00356 }
00357 }
00358 if (mTempOut) {
00359 mTempOut->Write();
00360 mTempOut->Close();
00361 if (gSystem->CopyFile(mTempOut->GetName(),fileName.Data()) ||
00362 gSystem->Unlink(mTempOut->GetName())) {
00363 LOG_ERROR << "Could not copy and/or delete temp vertexseedhist file!" << endm;
00364 }
00365 resNtuple = 0;
00366 }
00367 TFile out(fileName.Data(),"UPDATE");
00368 GetHistList()->Write();
00369 out.Close();
00370 }
00371
00372 Int_t StVertexSeedMaker::FillAssumed(){
00373 TDataSet* dbDataSet = GetDataBase("Calibrations/rhic");
00374 if (!dbDataSet) {
00375 LOG_ERROR << "Could not find Calibrations/rhic database" << endm;
00376 return kStErr;
00377 }
00378 St_vertexSeed* dbTableC =
00379 (St_vertexSeed*) (dbDataSet->FindObject("vertexSeed"));
00380 if (!dbTableC) {
00381 LOG_ERROR << "Could not find vertexSeed in database" << endm;
00382 return kStErr;
00383 }
00384 vertexSeed_st* dbTable = dbTableC->GetTable();
00385 a[0] = dbTable->x0;
00386 a[1] = dbTable->dxdz;
00387 a[2] = dbTable->y0;
00388 a[3] = dbTable->dydz;
00389 ea[0] = dbTable->err_x0;
00390 ea[1] = dbTable->err_dxdz;
00391 ea[2] = dbTable->err_y0;
00392 ea[3] = dbTable->err_dydz;
00393 LOG_INFO << "Assumed values:"
00394 << "\n x0 assumed = " << a[0] << " +/- " << ea[0]
00395 << "\n dxdz assumed = " << a[1] << " +/- " << ea[1]
00396 << "\n y0 assumed = " << a[2] << " +/- " << ea[2]
00397 << "\n dydz assumed = " << a[3] << " +/- " << ea[3]
00398 << endm;
00399 return kStOk;
00400 }
00401
00402 Int_t StVertexSeedMaker::GetVertexSeedTriggers(){
00403 TDataSet* dbDataSet = GetDataBase("Calibrations/rhic");
00404 if (!dbDataSet) {
00405 LOG_ERROR << "Could not find Calibrations/rhic database" << endm;
00406 return kStErr;
00407 }
00408 dbTriggersTable =
00409 (St_vertexSeedTriggers*) (dbDataSet->FindObject("vertexSeedTriggers"));
00410 if (!dbTriggersTable) {
00411 LOG_ERROR << "Could not find vertexSeedTriggers in database" << endm;
00412 return kStErr;
00413 }
00414 return kStOk;
00415 }
00416
00417 Bool_t StVertexSeedMaker::BetterErrors(){
00418 Bool_t better = kFALSE;
00419 if ((ep[0] < ea[0]) || (ep[1] < ea[1]) ||
00420 (ep[2] < ea[2]) || (ep[3] < ea[3])) better = kTRUE;
00421 if (better) { LOG_INFO << "Values have improved" << endm; }
00422 return better;
00423 }
00424
00425 Bool_t StVertexSeedMaker::ChangedValues(){
00426 Bool_t changed = kFALSE;
00427 for (int i = 0; i<4; i++) {
00428 double diff = TMath::Abs(p[i] - a[i]);
00429 if ((diff >= ep[i]) || (diff >= ea[i])) changed = kTRUE;
00430 }
00431 if (changed) { LOG_INFO << "Values have changed" << endm; }
00432 return changed;
00433 }
00434
00435 void StVertexSeedMaker::FillDateTime() {
00436 date = GetDate();
00437 time = GetTime();
00438 LOG_INFO << "event date = " << date << endm;
00439 LOG_INFO << "event time = " << time << endm;
00440 if (!useEventDateTime) GetFillDateTime();
00441 }
00442
00443 void StVertexSeedMaker::GetFillDateTime() {
00444
00445 StDbManager* mgr=StDbManager::Instance();
00446 StDbConfigNode* node=mgr->initConfig("RunLog_onl");
00447 StDbTable* tab=node->addDbTable("beamInfo");
00448 StDataBaseI* db=mgr->findDb("RunLog_onl");
00449 unsigned int ts;
00450 char queryStr[128];
00451
00452
00453 sprintf(queryStr,"%08d %06d",date,time);
00454 TString tdStr = queryStr;
00455 tdStr.Insert(4,'-').Insert(7,'-').Insert(13,':').Insert(16,':');
00456 const char* tdstr = tdStr.Data();
00457 sprintf(queryStr,
00458 " where beginTime<='%s' and deactive=0 order by beginTime desc limit 1",
00459 tdstr);
00460 ts = db->QueryDb(tab,queryStr);
00461
00462 if (ts) {
00463
00464 run = *(int*) (tab->getDataValue("runNumber",0));
00465 LOG_INFO << tdstr << " is from run " << run << endm;
00466 float thisFill = *(float*) (tab->getDataValue("blueFillNumber",0));
00467 sprintf(queryStr,
00468 " where blueFillNumber=%f and deactive=0 order by beginTime asc limit 1",
00469 thisFill);
00470 ts = db->QueryDb(tab,queryStr);
00471
00472
00473 char* start = tab->getBeginDateTime();
00474 fill = (int) thisFill;
00475 time = atoi(&(start[8]));
00476 start[8] = 0;
00477 date = atoi(start);
00478
00479 LOG_INFO << "Using fill no. = " << fill << endm;
00480 LOG_INFO << "Using fill date = " << date << endm;
00481 LOG_INFO << "Using fill time = " << time << endm;
00482 } else {
00483 LOG_WARN << "Could not find beamInfo in database\n" <<
00484 " Using event date/time." << endm;
00485 UseEventDateTime();
00486 }
00487 }
00488
00489 void StVertexSeedMaker::FitData() {
00490 LOG_INFO << "Now fitting the data..." <<
00491 "\n *****************************************************" << endm;
00492 TVirtualFitter *minuit = TVirtualFitter::Fitter(0,26);
00493
00494 minuit->SetFCN(fnch);
00495
00496
00497 static Double_t pstart[4] = {0., 0., 0., 0.};
00498 static Double_t pstep[4] = {0.0001, 0.0000001, 0.0001, 0.0000001};
00499 static Double_t plow[4] = {-3., -0.05, -3., -0.05};
00500 static Double_t phigh[4] = { 3., 0.05, 3., 0.05};
00501 minuit->SetParameter(0, "x0" , pstart[0], pstep[0], plow[0], phigh[0]);
00502 minuit->SetParameter(1, "dxdz", pstart[1], pstep[1], plow[1], phigh[1]);
00503 minuit->SetParameter(2, "y0" , pstart[2], pstep[2], plow[2], phigh[2]);
00504 minuit->SetParameter(3, "dydz", pstart[3], pstep[3], plow[3], phigh[3]);
00505
00506
00507 double arglist[10];
00508 arglist[0] = 500;
00509
00510 beamWidth = 0;
00511 do {
00512
00513 int status = minuit->ExecuteCommand("MIGRAD", arglist ,1);
00514 if (status) {
00515 LOG_ERROR << "StVertexMaker: error on migrad call, err = "
00516 << status << endm;
00517 return;
00518 }
00519
00520 double amin,edm,errdef;
00521 int nvpar,nparx;
00522 minuit->GetStats(amin,edm,errdef,nvpar,nparx);
00523 chi = amin/((double) (nverts-4));
00524 LOG_INFO << "beamWidth = " << beamWidth << ", chisq = " << amin << ", chisq/dof = " << chi <<
00525 "\n *****************************************************" << endm;
00526
00527 beamWidth += 0.005;
00528 } while (chi>1.1 && beamWidth<=0.15);
00529
00530 char pname[10];
00531 for (int i=0; i<4; i++)
00532 minuit->GetParameter(i, pname, p[i], ep[i], plow[i], phigh[i]);
00533 }
00534
00535 Int_t StVertexSeedMaker::Aggregate(Char_t* dir, const Char_t* cuts) {
00536
00537
00538
00539
00540 const char* defaultDir = "./";
00541 TString dirStr = dir;
00542 if (!dir) dirStr = defaultDir;
00543 if (dirStr.EndsWith("/")) dirStr.Append("vertexseedhist.*.root");
00544 StFile allFiles;
00545 allFiles.AddFile(dirStr.Data());
00546
00547 TList fileList;
00548 const char* fileName;
00549 for (int fileb = 0; fileb < allFiles.GetNBundles(); fileb++) {
00550 allFiles.GetNextBundle();
00551 int filen=0;
00552 while ((fileName = allFiles.GetFileName(filen++))) {
00553 fileList.Add(new TNamed(fileName,fileName));
00554 }
00555 }
00556 fileList.Sort();
00557
00558 TFile* currentFile=0;
00559 float* vals=0;
00560 int nfiles = fileList.GetSize();
00561 for (int filen = 0; filen < nfiles; filen++) {
00562 int fillf = fill;
00563 int datef = date;
00564 int timef = time;
00565 fileName = fileList.At(filen)->GetName();
00566 TString dateTime = fileName;
00567 dateTime.Remove(0,dateTime.Last('/') + 1);
00568 dateTime.Remove(0,dateTime.First('.') + 1).Remove(15);
00569 TString dateStr = dateTime;
00570 date = atoi(dateStr.Remove(8).Data());
00571 time = atoi(dateTime.Remove(0,9).Remove(6).Data());
00572 GetFillDateTime();
00573 if ((currentFile) && (fill != fillf)) {
00574 LOG_INFO << "Fill number has changed\n"
00575 << " Processing data from previous fill before continuing" << endm;
00576 int fillp = fill;
00577 int datep = date;
00578 int timep = time;
00579 fill = fillf;
00580 date = datef;
00581 time = timef;
00582 currentFile->Close();
00583 currentFile = 0;
00584 FindResult(kFALSE);
00585 Reset();
00586 fill = fillp;
00587 date = datep;
00588 time = timep;
00589 }
00590
00591 LOG_INFO << "Now opening file:\n " << fileName << endm;
00592 if (currentFile) currentFile->Close();
00593 currentFile = new TFile(fileName);
00594 TNtuple* curNtuple = (TNtuple*) currentFile->Get("resNtuple");
00595 if (!curNtuple) {
00596 LOG_ERROR << "No resNtuple found in " << fileName << endm;
00597 continue;
00598 }
00599 curNtuple->Draw(">>elist",cuts);
00600 TEventList* elist = (TEventList*) gDirectory->Get("elist");
00601 Int_t nentries = (Int_t) elist->GetN();
00602 for (Int_t entryn = 0; entryn < nentries; entryn++) {
00603 curNtuple->GetEntry(elist->GetEntry(entryn));
00604 vals = curNtuple->GetArgs();
00605 unsigned int tid = (unsigned int) vals[5];
00606 if (ValidTrigger(tid)) {
00607 resNtuple->Fill(vals);
00608 if (curNtuple->GetNvar()>13)
00609 addVert(vals[1],vals[2],vals[3],vals[4],vals[13],vals[14]);
00610 else
00611 addVert(vals[1],vals[2],vals[3],vals[4],0.,0.);
00612 } else {
00613 LOG_INFO << "Invalid trigger: " << tid << endm;
00614 }
00615 }
00616 LOG_INFO << "Current statistics: " << nverts << endm;
00617 }
00618 if (currentFile) {
00619 currentFile->Close();
00620 currentFile = 0;
00621 FindResult(kFALSE);
00622 }
00623 LOG_INFO << "Examined " << nfiles << " files" << endm;
00624 return nfiles;
00625 }
00626
00627
00628
00629
00630
00631
00632
00633
00634
00635
00636
00637
00638
00639
00640
00641
00642
00643
00644
00645
00646
00647
00648
00649
00650
00651
00652
00653
00654
00655
00656
00657
00658
00659
00660
00661
00662
00663
00664
00665
00666
00667
00668
00669
00670
00671
00672
00673
00674
00675
00676
00677
00678
00679
00680
00681
00682
00683
00684
00685
00686
00687
00688
00689
00690
00691
00692
00693
00694
00695
00696
00697
00698
00699
00700
00701
00702
00703
00704
00705
00706
00707
00708
00709
00710
00711
00712
00713
00714
00715
00716
00717
00718
00719
00720
00721
00722
00723
00724
00725
00726
00727
00728
00729
00730
00731
00732
00733
00734
00735
00736
00737
00738
00739
00740
00741
00742
00743
00744
00745
00746
00747
00748
00749
00750
00751
00752
00753
00754
00755
00756
00757
00758
00759
00760
00761
00762
00763
00764
00765
00766
00767
00768
00769
00770
00771
00772