22 #define TSelEventGen_cxx 
   40 ClassImp(TSelEventGen);
 
   45 TSelEventGen::TSelEventGen()
 
   46              : fBaseDir(
""), fNEvents(100000), fNTracks(100), fNTracksMax(-1),
 
   47                fRegenerate(kFALSE), fTotalGen(0), fFilesGenerated(0),
 
   48                fGenerateFun(0), fChain(0)
 
   51       fBaseDir=gProofServ->GetDataDir();
 
   53       fBaseDir.Remove(fBaseDir.Last(
'/'));
 
   54       fBaseDir.Remove(fBaseDir.Last(
'/'));
 
   66 void TSelEventGen::Begin(TTree *)
 
   68    TString option = GetOption();
 
   70    TMap *filemap = 
dynamic_cast<TMap *
> 
   71                      (fInput->FindObject(
"PROOF_FilesToProcess"));
 
   76       if (fInput->FindObject(
"PROOF_FilesToProcess")) {
 
   77          Error(
"Begin", 
"object 'PROOF_FilesToProcess' found but not a map" 
   78               " (%s)", fInput->FindObject(
"PROOF_FilesToProcess")->ClassName());
 
   80          Error(
"Begin", 
"object 'PROOF_FilesToProcess' not found");
 
   90 void TSelEventGen::SlaveBegin(TTree *tree)
 
   94    TString option = GetOption();
 
   98    Bool_t found_basedir=kFALSE;
 
   99    Bool_t found_nevents=kFALSE;
 
  100    Bool_t found_ntracks=kFALSE;
 
  101    Bool_t found_ntrkmax=kFALSE;
 
  102    Bool_t found_regenerate=kFALSE;
 
  108    while ((obj = nxt())){
 
  110       sinput=obj->GetName();
 
  113       if (sinput.Contains(
"PROOF_BenchmarkBaseDir")){
 
  114          TNamed* a = 
dynamic_cast<TNamed*
>(obj);
 
  116             TString bdir = a->GetTitle();
 
  119                Bool_t isLocal = !strcmp(u.GetProtocol(), 
"file") ? kTRUE : kFALSE;
 
  120                if (isLocal && !gSystem->IsAbsoluteFileName(u.GetFile()))
 
  121                   u.SetFile(TString::Format(
"%s/%s", fBaseDir.Data(), u.GetFile()));
 
  123                   if ((gSystem->AccessPathName(u.GetFile()) &&
 
  124                      gSystem->mkdir(u.GetFile(), kTRUE) == 0) ||
 
  125                      !gSystem->AccessPathName(u.GetFile(), kWritePermission)) {
 
  127                      fBaseDir = u.GetFile();
 
  128                      Info(
"SlaveBegin", 
"using base directory \"%s\"", fBaseDir.Data());
 
  131                      Warning(
"SlaveBegin", 
"\"%s\" directory is not writable or not existing," 
  132                            " using default directory: %s",
 
  133                            bdir.Data(), fBaseDir.Data());
 
  138                   Info(
"SlaveBegin", 
"using non local base directory \"%s\"", fBaseDir.Data());
 
  141                Info(
"SlaveBegin", 
"using default directory: %s",
 
  147             Error(
"SlaveBegin", 
"PROOF_BenchmarkBaseDir not type TNamed");
 
  151       if (sinput.Contains(
"PROOF_BenchmarkNEvents")){
 
  152          TParameter<Long64_t>* a=
dynamic_cast<TParameter<Long64_t>*
>(obj);
 
  154             fNEvents= a->GetVal();
 
  158             Error(
"SlaveBegin", 
"PROOF_BenchmarkEvents not type TParameter" 
  163       if (sinput.Contains(
"PROOF_BenchmarkNTracks")){
 
  164          TParameter<Int_t>* a=
dynamic_cast<TParameter<Int_t>*
>(obj);
 
  166             fNTracks=a->GetVal();
 
  170             Error(
"SlaveBegin", 
"PROOF_BenchmarkNTracks not type TParameter" 
  175       if (sinput.Contains(
"PROOF_BenchmarkNTracksMax")){
 
  176          TParameter<Int_t>* a=
dynamic_cast<TParameter<Int_t>*
>(obj);
 
  178             fNTracksMax=a->GetVal();
 
  182             Error(
"SlaveBegin", 
"PROOF_BenchmarkNTracksMax not type TParameter" 
  187       if (sinput.Contains(
"PROOF_BenchmarkRegenerate")){
 
  188          TParameter<Int_t>* a=
dynamic_cast<TParameter<Int_t>*
>(obj);
 
  190             fRegenerate=a->GetVal();
 
  191             found_regenerate=kTRUE;
 
  194             Error(
"SlaveBegin", 
"PROOF_BenchmarkRegenerate not type TParameter" 
  199       if (sinput.Contains(
"PROOF_GenerateFun")){
 
  200          TNamed *a = 
dynamic_cast<TNamed*
>(obj);
 
  201          if (!(fGenerateFun = dynamic_cast<TMacro *>(fInput->FindObject(a->GetTitle())))) {
 
  202             Error(
"SlaveBegin", 
"PROOF_GenerateFun requires the TMacro object in the input list");
 
  209       Warning(
"SlaveBegin", 
"PROOF_BenchmarkBaseDir not found; using default:" 
  210                             " %s", fBaseDir.Data());
 
  213       Warning(
"SlaveBegin", 
"PROOF_BenchmarkNEvents not found; using default:" 
  217       Warning(
"SlaveBegin", 
"PROOF_BenchmarkNTracks not found; using default:" 
  221       Warning(
"SlaveBegin", 
"PROOF_BenchmarkNTracksMax not found; using default:" 
  223    } 
else if (fNTracksMax <= fNTracks) {
 
  224       Warning(
"SlaveBegin", 
"PROOF_BenchmarkNTracksMax must be larger then" 
  225                             " fNTracks=%d ; ignoring", fNTracks);
 
  227       found_ntrkmax = kFALSE;
 
  229    if (!found_regenerate){
 
  230       Warning(
"SlaveBegin", 
"PROOF_BenchmarkRegenerate not found; using" 
  231                             " default: %d", fRegenerate);
 
  234    fFilesGenerated = 
new TList();
 
  236    TString hostname(TUrl(gSystem->HostName()).GetHostFQDN());
 
  237    TString thisordinal = gProofServ ? gProofServ->GetOrdinal() : 
"n.d";
 
  238    TString sfilegenerated =
 
  239       TString::Format(
"PROOF_FilesGenerated_%s_%s", hostname.Data(), thisordinal.Data());
 
  240    fFilesGenerated->SetName(sfilegenerated);
 
  241    fFilesGenerated->SetOwner(kTRUE);
 
  258 Long64_t TSelEventGen::GenerateFiles(
const char *filename, Long64_t sizenevents)
 
  261    TDirectory* savedir = gDirectory;
 
  264    TFile *f = TFile::Open(filename, 
"RECREATE");
 
  268    if (!f || f->IsZombie()) 
return 0;
 
  270    Event *
event=
new Event();
 
  272    TTree* eventtree= 
new TTree(
"EventTree", 
"Event Tree");
 
  273    eventtree->SetDirectory(f);
 
  275    const Int_t buffersize=32000;
 
  276    eventtree->Branch(
"event", 
"Event", &ep, buffersize, 1);
 
  277    eventtree->AutoSave();
 
  280    Long64_t size_generated=0;
 
  283    Int_t ntrks = fNTracks;
 
  285    Info(
"GenerateFiles", 
"Generating %s", filename);
 
  286    while (sizenevents--){
 
  288       if (fNTracksMax > fNTracks) {
 
  290          ntrks = fNTracks + gRandom->Integer(fNTracksMax - fNTracks);
 
  292       event->Build(i++, ntrks, 0);
 
  293       size_generated+=eventtree->Fill();
 
  295    nentries=eventtree->GetEntries();
 
  296    Info(
"GenerateFiles", 
"%s generated with %lld entries", filename, nentries);
 
  297    savedir = gDirectory;
 
  299    f = eventtree->GetCurrentFile();
 
  302    eventtree->SetDirectory(0);
 
  327 Bool_t TSelEventGen::Process(Long64_t entry)
 
  335    TDSetElement *fCurrent = 0;
 
  337    if (fInput && (elemPair = dynamic_cast<TPair *>
 
  338                       (fInput->FindObject(
"PROOF_CurrentElement")))) {
 
  339       if ((fCurrent = dynamic_cast<TDSetElement *>(elemPair->Value()))) {
 
  340          Info(
"Process", 
"entry %lld: file: '%s'", entry, fCurrent->GetName());
 
  342          Error(
"Process", 
"entry %lld: no file specified!", entry);
 
  348    TString filename(fCurrent->GetName());
 
  349    if (!fBaseDir.IsNull()) {
 
  350       if (fBaseDir.Contains(
"<fn>")) {
 
  352          filename.ReplaceAll(
"<fn>", fCurrent->GetName());
 
  354          filename.Form(
"%s/%s", fBaseDir.Data(), fCurrent->GetName());
 
  357    TString fndset(filename);
 
  360    TString seed = TString::Format(
"%s/%s", gSystem->HostName(), filename.Data()), dsrv;
 
  361    TUrl basedirurl(filename, kTRUE);
 
  362    if (!strcmp(basedirurl.GetProtocol(), 
"file")) {
 
  363       TProofServ::GetLocalServer(dsrv);
 
  364       TProofServ::FilterLocalroot(fndset, dsrv);
 
  368    Long64_t neventstogenerate = fNEvents;
 
  370    Long64_t entries_file=0;
 
  372    Bool_t filefound=kFALSE;
 
  375    if (!fRegenerate && !gSystem->GetPathInfo(filename, filestat)) { 
 
  376       TFile *f = TFile::Open(filename);
 
  377       if (f && !f->IsZombie()){
 
  378          TTree* t = (TTree *) f->Get(
"EventTree");
 
  380             entries_file = t->GetEntries();
 
  381             if (entries_file == neventstogenerate) {
 
  383                Info(
"Process", 
"bench file (%s, entries=%lld) exists:" 
  384                                " skipping generation.", filename.Data(), entries_file);
 
  385                filesize = f->GetSize();
 
  396    TString bdir(fBaseDir);
 
  397    bdir.ReplaceAll(
"<fn>", 
"");
 
  398    if (!gSystem->AccessPathName(bdir)) {
 
  399       Long_t devid, devbsz, devbtot, devbfree;
 
  400       gSystem->GetFsInfo(bdir, &devid, &devbsz, &devbtot, &devbfree);
 
  402       Long_t szneed = 1024 * 1024 * 1024, tomb = 1024 * 1024;
 
  403       if (devbfree * devbsz < szneed || devbfree < 0.1 * devbtot) {
 
  404          Error(
"Process", 
"not enough free space on device (%ld MB < {%ld, %ld} MB):" 
  405                           " skipping generation of: %s",
 
  406                           (devbfree * devbsz) / tomb,
 
  407                           szneed / tomb, (Long_t) (0.1 * devbtot * devbsz / tomb),
 
  409          fStatus = TSelector::kAbortFile;
 
  414       gRandom->SetSeed(static_cast<UInt_t>(TMath::Hash(seed)));
 
  416          TString fargs = TString::Format(
"\"%s\",%lld", filename.Data(), neventstogenerate);
 
  417          entries_file = (Long64_t) fGenerateFun->Exec(fargs);
 
  419          entries_file = GenerateFiles(filename, neventstogenerate);
 
  422       TFile *f = TFile::Open(filename);
 
  423       if (f && !f->IsZombie()) {
 
  424          filesize = f->GetSize();
 
  428          Error(
"Process", 
"can not open generated file: %s", filename.Data());
 
  429          fStatus = TSelector::kAbortFile;
 
  437    TFileInfoMeta* fimeta = 
new TFileInfoMeta(
"/EventTree", 
"TTree", entries_file);
 
  439    if (!strcmp(TUrl(filename,kTRUE).GetProtocol(), 
"file"))
 
  440       md5 = TMD5::FileChecksum(filename);
 
  441    TString md5s = (md5) ? md5->AsString() : 
"";
 
  442    TFileInfo *fi = 
new TFileInfo(TString::Format(
"%s%s", dsrv.Data(), fndset.Data()),
 
  443                                  filesize, uuid.AsString(), md5s.Data(), fimeta);
 
  447    fi->SetBit(TFileInfo::kStaged);
 
  450    if (fFilesGenerated) fFilesGenerated->Add(fi);
 
  460 void TSelEventGen::SlaveTerminate()
 
  462    if (fFilesGenerated && fFilesGenerated->GetSize() > 0) {
 
  463       fOutput->Add(fFilesGenerated);
 
  464       Info(
"SlaveTerminate",
 
  465               "list '%s' of files generated by this worker added to the output list",
 
  466               fFilesGenerated->GetName());
 
  468       if (!fFilesGenerated) {
 
  469          Warning(
"SlaveTerminate", 
"no list of generated files defined!");
 
  471          Warning(
"SlaveTerminate", 
"list of generated files is empty!");
 
  481 void TSelEventGen::Terminate()
 
  487 void TSelEventGen::Print(Option_t *)
 const 
  489    Printf(
"fNEvents=%lld", fNEvents);
 
  490    Printf(
"fBaseDir=%s", fBaseDir.Data());
 
  491    Printf(
"fNTracks=%d", fNTracks);
 
  492    Printf(
"fRegenerate=%d", fRegenerate);