#include "TMessage.h"
#include "TObjString.h"
#include "TParameter.h"
#include "TPerfStats.h"
#include "TPacketizerProgressive.h"
#include "TProofDebug.h"
#include "TProofServ.h"
#include "TSocket.h"
#include "Riostream.h"
#include "TClass.h"
#include "TDSet.h"
#include "TError.h"
#include "TList.h"
#include "TMap.h"
#include "TSlave.h"
#include "TTimer.h"
#include "TUrl.h"
TPacketizerProgressive::TFileNode::TFileNode(const char *name)
   : fNodeName(name), fFiles(new TList), fUnAllocFileNext(0),fActFiles(new TList),
     fActFileNext(0), fMySlaveCnt(0), fSlaveCnt(0)
{
   
   fFiles->SetOwner();
   fActFiles->SetOwner(kFALSE);
}
TPacketizerProgressive::TFileNode::~TFileNode()
{
   
   SafeDelete(fFiles);
   SafeDelete(fActFiles);
}
void TPacketizerProgressive::TFileNode::Add(TDSetElement *elem)
{
   
   TFileStat *f = new TFileStat(this,elem);
   fFiles->Add(f);
   if (fUnAllocFileNext == 0) fUnAllocFileNext = fFiles->First();
}
void TPacketizerProgressive::TFileNode::DecSlaveCnt(const char *wrk)
{
   
   if (fNodeName != wrk)
      fSlaveCnt--;
   R__ASSERT(fSlaveCnt >= 0);
}
TPacketizerProgressive::TFileStat *TPacketizerProgressive::TFileNode::GetNextUnAlloc()
{
   
   TObject *next = fUnAllocFileNext;
   if (next != 0) {
      
      fActFiles->Add(next);
      if (fActFileNext == 0) fActFileNext = fActFiles->First();
      
      fUnAllocFileNext = fFiles->After(fUnAllocFileNext);
   }
   return (TFileStat *) next;
}
TPacketizerProgressive::TFileStat *TPacketizerProgressive::TFileNode::GetNextActive()
{
   
   TObject *next = fActFileNext;
   if (fActFileNext != 0) {
      fActFileNext = fActFiles->After(fActFileNext);
      if (fActFileNext == 0) fActFileNext = fActFiles->First();
   }
   return (TFileStat *) next;
}
void TPacketizerProgressive::TFileNode::RemoveActive(TFileStat *file)
{
   
   if (fActFileNext == file) fActFileNext = fActFiles->After(file);
   fActFiles->Remove(file);
   if (fActFileNext == 0) fActFileNext = fActFiles->First();
}
Int_t TPacketizerProgressive::TFileNode::Compare(const TObject *other) const
{
   
   
   const TFileNode *obj = dynamic_cast<const TFileNode*>(other);
   R__ASSERT(obj != 0);
   Int_t myVal = GetSlaveCnt();
   Int_t otherVal = obj->GetSlaveCnt();
   if (myVal < otherVal) {
      return -1;
   } else if (myVal > otherVal) {
      return 1;
   } else {
      return 0;
   }
}
void TPacketizerProgressive::TFileNode::Print(Option_t *) const
{
   
   cout << "OBJ: " << IsA()->GetName() << "\t" << fNodeName
        << "\tMySlaveCount " << fMySlaveCnt
        << "\tSlaveCount " << fSlaveCnt << endl;
}
void TPacketizerProgressive::TFileNode::Reset()
{
   
   fUnAllocFileNext = fFiles->First();
   fActFiles->Clear();
   fActFileNext = 0;
   fSlaveCnt = 0;
   fMySlaveCnt = 0;
}
Int_t TPacketizerProgressive::TFileNode::GetNumberOfActiveFiles() const
{
   
   return (fActFiles ? fActFiles->GetSize() : 0);
}
Bool_t TPacketizerProgressive::TFileNode::HasActiveFiles()
{
   
   if (fActFiles && fActFiles->GetSize() > 0)
      return kTRUE;
   return kFALSE;
}
TPacketizerProgressive::TFileStat::TFileStat(TFileNode *node, TDSetElement *elem)
   : fIsDone(kFALSE), fNode(node), fElement(elem), fNextEntry(elem->GetFirst())
{
   
}
TPacketizerProgressive::TSlaveStat::TSlaveStat(TSlave *slave)
   : fSlave(slave), fFileNode(0), fCurFile(0), fCurElem(0), fProcessed(0)
{
   
}
const char *TPacketizerProgressive::TSlaveStat::GetName() const
{
   
   return (fSlave ? fSlave->GetName() : "");
}
ClassImp(TPacketizerProgressive)
TPacketizerProgressive::TPacketizerProgressive(TDSet* dset,
                                               TList* slaves,
                                               Long64_t first,
                                               Long64_t num,
                                               TList* )
   : fDset(dset), fSlaves(slaves), fSlavesRemaining(0), fFirstEvent(first),
     fTotalEvents(num), fEntriesSeen(0), fFilesOpened(0), fEstTotalEntries(0)
{
   
   PDB(kPacketizer,1) Info("TPacketizerProgressive",
                           "enter (first %lld, num %lld)", first, num);
   fProcessed = 0;
   if (fTotalEvents != -1) {
      
      Error("TPacketizerProgressive",
            "this packetizer does not handle TDSet regions");
   }
   fSlavesRemaining = new TList;
   fSlaveStats = new TMap;
   fUnAllocSlaves = new TList;
   fUnAllocNonSlaves = new TList;
   fActiveSlaves = new TList;
   fActiveNonSlaves = new TList;
   fLastEntrySizes = new TList;
   fNewFileSlaves = new THashTable;
   if (fSlaves)
      fSlavesRemaining->AddAll(fSlaves);
   fValid = kTRUE;
   Init();
}
TPacketizerProgressive::~TPacketizerProgressive()
{
   
   fSlaveStats->DeleteValues();
   delete fSlavesRemaining;
   delete fSlaveStats;
   delete fUnAllocSlaves;
   delete fUnAllocNonSlaves;
   delete fActiveSlaves;
   delete fActiveNonSlaves;
   delete fLastEntrySizes;
   delete fNewFileSlaves;
   SafeDelete(fProgress);
}
void TPacketizerProgressive::Init()
{
   
   
   TMap host_map; 
   TIter i(fSlaves);
   i.Reset();
   TSlave* slave;
   
   
   
   while ((slave = (TSlave*) i.Next())) {
      PDB(kPacketizer, 3)
         Info("Init", "adding info for slave %s", slave->GetName());
      TSlaveStat* ss = new TSlaveStat(slave);
      fSlaveStats->Add(slave, ss);
      
      
      TObjString host(slave->GetName());
      TFileNode* fn = (TFileNode*) host_map.GetValue(&host);
      if (! fn) {
         fn = new TFileNode(slave->GetName());
         host_map.Add(new TObjString(slave->GetName()), fn);
      }
      fn->IncMySlaveCnt();
      ss->SetFileNode(fn);
   }
   
   fDset->Lookup();
   
   THashTable slaves_added; 
                            
   TMap nonslaves_added; 
   fDset->Reset();
   TDSetElement* e;
   while ((e = (TDSetElement*) fDset->Next())) {
      TUrl url = e->GetFileName();
      TObjString host = url.GetHost();
      PDB(kPacketizer, 3) Info("Init", "found TDSetElement on host %s",
                               host.GetString().Data());
      TFileNode* fn = (TFileNode*) host_map.GetValue(&host);
      if (! fn) {
         if (! (fn = (TFileNode*) nonslaves_added.GetValue(&host))) {
            PDB(kPacketizer, 3) Info("Init", "adding info for non-slave %s",
                                     host.GetString().Data());
            
            
            fn = new TFileNode(host.GetString().Data());
            fUnAllocNonSlaves->Add(fn);
            nonslaves_added.Add(&host, fn);
         }
      } else {
         if (! slaves_added.FindObject(fn)) {
            fUnAllocSlaves->Add(fn);
            slaves_added.Add(fn);
         }
      }
      fn->Add(e);
   }
   host_map.DeleteKeys();
   
   fProgress = new TTimer;
   fProgress->SetObject(this);
   fProgress->Start(500, kFALSE);
   PDB(kPacketizer, 2) Info("Init", "finished init");
}
Long64_t TPacketizerProgressive::GetEntriesProcessed(TSlave* s) const
{
   
   TSlaveStat* stat = (TSlaveStat*) fSlaveStats->GetValue(s);
   return stat->GetEntriesProcessed();
}
TDSetElement *TPacketizerProgressive::BuildPacket(TSlaveStat* stat,
                                                  Long64_t size)
{
   
   
   
   TFileStat* fs = stat->GetCurrentFile();
   if (! fs) {
      Error("BuildPacket", "no TFileStat assigned");
      return 0;
   }
   Long64_t elem_entries = stat->GetCurrentElement()->GetNum();
   Long64_t num = size;
   Long64_t entries_remaining = elem_entries - fs->GetNextEntry();
   
   
   if (elem_entries != -1 && (num > entries_remaining ||
                              (entries_remaining < (num * 2)))) {
      num = entries_remaining;
   }
   PDB(kPacketizer,3) Info("BuildPacket",
                           "packet of size %lld requested (assigning %lld).  pos: %lld, num entries: %lld",
                           size, num, fs->GetNextEntry(), elem_entries);
   TDSetElement* base = stat->GetCurrentElement();
   TDSetElement* packet = CreateNewPacket(base, fs->GetNextEntry(), num);
   fs->MoveNextEntry(num);
   stat->IncEntriesProcessed(num);
   fProcessed += num;
   
   if ( fs->GetNextEntry() == elem_entries) {
      fs->SetDone();
      fs->GetNode()->RemoveActive(fs);
      if (! stat->GetFileNode()->HasActiveFiles()) {
         fActiveSlaves->Remove(stat->GetFileNode());
         fActiveNonSlaves->Remove(stat->GetFileNode());
      }
   }
   return packet;
}
void TPacketizerProgressive::RecalculatePacketSize(Long64_t newCount)
{
   
   if (fLastEntrySizes->GetSize() >= kEntryListSize) {
      while (fLastEntrySizes->GetSize() >= (kEntryListSize - 1)) {
         TParameter<Long64_t>* first = (TParameter<Long64_t>*) fLastEntrySizes->First();
         fLastEntrySizes->Remove(first);
         delete first;
      }
   }
   
   fLastEntrySizes->AddLast(new TParameter<Long64_t>("", newCount));
   
   Long64_t total = 0, mean = 0, elems_remaining;
   elems_remaining = fDset->GetListOfElements()->GetSize() - fFilesOpened;
   PDB(kPacketizer, 4) Info("RecalculatePacketSize", "files opened: %lld, fdset size: %d, elems remaining: %lld",
                            fFilesOpened, fDset->GetListOfElements()->GetSize(),
                            elems_remaining);
   TParameter<Long64_t>* count;
   TIter i(fLastEntrySizes);
   while ((count = (TParameter<Long64_t>*) i.Next())) {
      total += count->GetVal();
   }
   mean = total / fLastEntrySizes->GetSize();
   fEstTotalEntries = fEntriesSeen + mean * elems_remaining;
   Int_t num_slaves = fSlavesRemaining->GetSize();
   fPacketSize = fEstTotalEntries / (num_slaves * 20);
   if (fPacketSize == 0) {
      fPacketSize = 1;
   }
   PDB(kPacketizer, 3) Info("RecalculatePacketSize", "estimated number of entries: %lld, new packet size: %lld",
                            fEstTotalEntries, fPacketSize);
}
TPacketizerProgressive::TFileStat *TPacketizerProgressive::GetNextActive(TSlaveStat* stat)
{
   
   fActiveSlaves->Sort();
   fActiveNonSlaves->Sort();
   TFileStat* file = 0;
   
   if (stat->GetFileNode() && stat->GetFileNode()->HasActiveFiles()) {
      PDB(kPacketizer, 3) Info("GetNextActive",
                               "Assigning slave %s TDSetElement from current data source",
                               stat->GetName());
      file = stat->GetFileNode()->GetNextActive();
      return file;
   }
   
   if (fActiveNonSlaves->GetSize() &&
       ((TFileNode*) fActiveNonSlaves->First())->GetSlaveCnt() < kNonSlaveHostConnLim) {
      PDB(kPacketizer, 3) Info("GetNextActive",
                               "Assigning slave %s TDSetElement from non-slave data source",
                               stat->GetName());
      file = ((TFileNode*) fActiveNonSlaves->First())->GetNextActive();
      return file;
   }
   
   if (fActiveSlaves->GetSize() &&
       ((TFileNode*) fActiveSlaves->First())->GetSlaveCnt() < kSlaveHostConnLim) {
      PDB(kPacketizer, 3) Info("GetNextActive",
                               "Assigning slave %s TDSetElement from peer data source",
                               stat->GetName());
      file = ((TFileNode*) fActiveSlaves->First())->GetNextActive();
      return file;
   }
   return 0;
}
TPacketizerProgressive::TFileStat *TPacketizerProgressive::GetNextUnAlloc(TSlaveStat* stat)
{
   
   fUnAllocSlaves->Sort();
   fUnAllocNonSlaves->Sort();
   TFileStat* file = 0;
   TFileNode* fn = 0;
   
   if (stat->GetFileNode() && stat->GetFileNode()->HasUnAllocFiles()) {
      PDB(kPacketizer, 3) Info("GetNextUnAlloc",
                               "Assigning slave %s TDSetElement from current data source",
                               stat->GetName());
      fn = stat->GetFileNode();
      file = fn->GetNextUnAlloc();
      
      if (fUnAllocNonSlaves->FindObject(fn) &&
          ! fActiveNonSlaves->FindObject(fn)) {
         fActiveNonSlaves->Add(fn);
      }
      if (fUnAllocSlaves->FindObject(fn) &&
          ! fActiveSlaves->FindObject(fn)) {
         fActiveSlaves->Add(fn);
      }
      
      
      if (!fn->HasUnAllocFiles()) {
         if (fUnAllocNonSlaves->FindObject(fn)) {
            fUnAllocNonSlaves->Remove(fn);
         }
         if (fUnAllocSlaves->FindObject(fn)) {
            fUnAllocSlaves->Remove(fn);
         }
      }
      return file;
   }
   
   if (fUnAllocNonSlaves->GetSize() &&
       ((TFileNode*)fUnAllocNonSlaves->First())->GetSlaveCnt() < kNonSlaveHostConnLim) {
      PDB(kPacketizer, 3) Info("GetNextUnAlloc",
                               "Assigning slave %s TDSetElement from non-slave data source",
                               stat->GetName());
      fn = (TFileNode*) fUnAllocNonSlaves->First();
      file = fn->GetNextUnAlloc();
      if (! fn->HasUnAllocFiles()) {
         fUnAllocNonSlaves->Remove(fn);
      }
      if (! fActiveNonSlaves->FindObject(fn)) {
         fActiveNonSlaves->Add(fn);
      }
      return file;
   }
   
   if (fUnAllocSlaves->GetSize() &&
       ((TFileNode*) fUnAllocSlaves->First())->GetSlaveCnt() < kSlaveHostConnLim) {
      PDB(kPacketizer, 3) Info("GetNextUnAlloc",
                               "Assigning slave %s TDSetElement from peer data source",
                               stat->GetName());
      fn = (TFileNode*) fUnAllocSlaves->First();
      file = fn->GetNextUnAlloc();
      if (! fn->HasUnAllocFiles()) {
         fUnAllocSlaves->Remove(fn);
      }
      if (! fActiveNonSlaves->FindObject(fn)) {
         fActiveSlaves->Add(fn);
      }
      return file;
   }
   return 0;
}
TDSetElement *TPacketizerProgressive::GetNextPacket(TSlave *s, TMessage *r)
{
   
   PDB(kPacketizer, 3) Info("GetNextPacket", "enter with slave %s", s->GetName());
   TSlaveStat* stat = (TSlaveStat*) fSlaveStats->GetValue(s);
   PDB(kPacketizer, 4) Info("GetNextPacket", "current file node (%s) has %d connections",
                            stat->GetFileNode()->GetName(), stat->GetFileNode()->GetSlaveCnt());
   PDB(kPacketizer, 4) Info("GetNextPacket", "unalloc'd slaves: %d, unalloc'd non-slaves: %d, active slaves: %d, active non-slaves: %d",
                            fUnAllocSlaves->GetSize(), fUnAllocNonSlaves->GetSize(), fActiveSlaves->GetSize(), fActiveNonSlaves->GetSize());
   Double_t latency, proctime, proccpu;
   Long64_t bytesread, numentries;
   
   if (stat->GetCurrentElement() && ! stat->GetCurrentFile()->IsDone()) {
      (*r) >> latency >> proctime >> proccpu >> bytesread >> numentries;
      if (gPerfStats != 0) {
         gPerfStats->PacketEvent(s->GetOrdinal(), s->GetName(), stat->GetCurrentElement()->GetFileName(),
                                 numentries, latency, proctime, proccpu, bytesread);
      }
   }
   
   if (stat->GetCurrentFile() && stat->GetCurrentFile()->IsDone()) {
      if (gPerfStats != 0) {
         TFileStat* file = stat->GetCurrentFile();
         gPerfStats->FileEvent(s->GetOrdinal(), s->GetName(), file->GetNode()->GetName(),
                               file->GetElement()->GetFileName(), kFALSE);
      }
      stat->SetCurrentElement(0);
      stat->SetCurrentFile(0);
   }
   
   if (fNewFileSlaves->FindObject(stat)) {
      RecalculatePacketSize(numentries);
      PDB(kPacketizer, 3) Info("GetNextPacket",
         "Newly opened file has %lld entries; updated packet size to %lld",
         numentries, fPacketSize);
      fEntriesSeen += numentries;
      fFilesOpened++;
      stat->GetCurrentElement()->SetNum(numentries);
      fNewFileSlaves->Remove(stat);
   }
   
   if (stat->GetCurrentFile() && ! stat->GetCurrentFile()->IsDone()) {
      if (stat->GetCurrentElement()->GetNum() == -1) {
         
         PDB(kPacketizer, 3) Info("GetNextPacket",
            "working on a packet that isn't fully opened, waiting");
         return (TDSetElement*) -1;
      }
      return BuildPacket(stat, fPacketSize);
   }
   
   Int_t foundUnallocatedFile = -1; 
   stat->GetFileNode()->DecSlaveCnt(s->GetName());
   TFileStat* fs = 0;
   
   if ((fs = GetNextUnAlloc(stat))) {
      PDB(kPacketizer, 3) Info("AssignElement", "giving slave %s unallocated file",
                               stat->GetName());
      stat->SetFileNode(fs->GetNode());
      stat->GetFileNode()->IncSlaveCnt(s->GetName());
      stat->SetCurrentFile(fs);
      stat->SetCurrentElement(fs->GetElement());
      fNewFileSlaves->Add(stat);
      foundUnallocatedFile = 1;
   
   } else if ((fs = GetNextActive(stat))) {
      PDB(kPacketizer, 3) Info("AssignElement", "giving slave %s active file",
                               stat->GetName());
      stat->SetFileNode(fs->GetNode());
      stat->GetFileNode()->IncSlaveCnt(s->GetName());
      stat->SetCurrentFile(fs);
      stat->SetCurrentElement(fs->GetElement());
      if (stat->GetCurrentElement()->GetNum() == -1) {
         
         PDB(kPacketizer, 3) Info("AssignElement", "grabbed a packet that isn't fully opened, waiting");
         return (TDSetElement*) -1;
      }
      foundUnallocatedFile = 0;
   }
   if (foundUnallocatedFile == -1) { 
      
      PDB(kPacketizer, 3) Info("GetNextPacket", "no more packets available");
      fSlavesRemaining->Remove(s);
      return 0;
   } else {
      if (gPerfStats != 0) {
         TFileStat* file = stat->GetCurrentFile();
         gPerfStats->FileEvent(s->GetOrdinal(), s->GetName(), file->GetNode()->GetName(),
                               file->GetElement()->GetFileName(), kTRUE);
      }
      if (foundUnallocatedFile == 1)
         return BuildPacket(stat, 1);
      else
         return BuildPacket(stat, fPacketSize);
   }
}
Bool_t TPacketizerProgressive::HandleTimer(TTimer *)
{
   
   PDB(kPacketizer, 4) Info("HandlerTimer", "estimated total entries: %lld, entries processed: %lld",
                            fEstTotalEntries, fProcessed);
   if (fProgress == 0) return kFALSE; 
   
   if (fEstTotalEntries <= 0) return kFALSE;
   TMessage m(kPROOF_PROGRESS);
   m << fEstTotalEntries << fProcessed;
   
   gProofServ->GetSocket()->Send(m);
   return kFALSE; 
}
This page has been automatically generated. If you have any comments or suggestions about the page layout send a mail to ROOT support, or contact the developers with any questions or problems regarding ROOT.