Logo ROOT   6.08/07
Reference Guide
TProofLite.cxx
Go to the documentation of this file.
1 // @(#)root/proof:$Id: 7735e42a1b96a9f40ae76bd884acac883a178dee $
2 // Author: G. Ganis March 2008
3 
4 /*************************************************************************
5  * Copyright (C) 1995-2000, Rene Brun and Fons Rademakers. *
6  * All rights reserved. *
7  * *
8  * For the licensing terms see $ROOTSYS/LICENSE. *
9  * For the list of contributors see $ROOTSYS/README/CREDITS. *
10  *************************************************************************/
11 
12 /** \class TProofLite
13 \ingroup proofkernel
14 
15 This class starts a PROOF session on the local machine: no daemons,
16 client and master merged, communications via UNIX-like sockets.
17 By default the number of workers started is NumberOfCores+1; a
18 different number can be forced on construction.
19 
20 */
21 
22 #include "TProofLite.h"
23 
24 #ifdef WIN32
25 # include <io.h>
26 # include "snprintf.h"
27 #endif
28 #include "RConfigure.h"
29 #include "TDSet.h"
30 #include "TEnv.h"
31 #include "TError.h"
32 #include "TFile.h"
33 #include "TFileCollection.h"
34 #include "TFileInfo.h"
35 #include "THashList.h"
36 #include "TMessage.h"
37 #include "TMonitor.h"
38 #include "TObjString.h"
39 #include "TPluginManager.h"
40 #include "TDataSetManager.h"
41 #include "TDataSetManagerFile.h"
42 #include "TParameter.h"
43 #include "TPRegexp.h"
44 #include "TProofQueryResult.h"
45 #include "TProofServ.h"
46 #include "TQueryResultManager.h"
47 #include "TROOT.h"
48 #include "TServerSocket.h"
49 #include "TSlave.h"
50 #include "TSortedList.h"
51 #include "TTree.h"
52 #include "TVirtualProofPlayer.h"
53 #include "TSelector.h"
54 #include "TPackMgr.h"
55 
57 
58 Int_t TProofLite::fgWrksMax = -2; // Unitialized max number of workers
59 
60 ////////////////////////////////////////////////////////////////////////////////
61 /// Create a PROOF environment. Starting PROOF involves either connecting
62 /// to a master server, which in turn will start a set of slave servers, or
63 /// directly starting as master server (if master = ""). Masterurl is of
64 /// the form: [proof[s]://]host[:port]. Conffile is the name of the config
65 /// file describing the remote PROOF cluster (this argument alows you to
66 /// describe different cluster configurations).
67 /// The default is proof.conf. Confdir is the directory where the config
68 /// file and other PROOF related files are (like motd and noproof files).
69 /// Loglevel is the log level (default = 1). User specified custom config
70 /// files will be first looked for in $HOME/.conffile.
71 
72 TProofLite::TProofLite(const char *url, const char *conffile, const char *confdir,
73  Int_t loglevel, const char *alias, TProofMgr *mgr)
74 {
75  fUrl.SetUrl(url);
76 
77  // Default initializations
78  fServSock = 0;
79  fCacheLock = 0;
80  fQueryLock = 0;
81  fQMgr = 0;
82  fDataSetManager = 0;
83  fDataSetStgRepo = 0;
84  fReInvalid = new TPMERegexp("[^A-Za-z0-9._-]");
85  InitMembers();
86 
87  // This may be needed during init
88  fManager = mgr;
89 
90  // Default server type
92 
93  // Default query mode
94  fQueryMode = kSync;
95 
96  // Client and master are merged
100 
101  // Flag that we are a client
102  if (!gSystem->Getenv("ROOTPROOFCLIENT")) gSystem->Setenv("ROOTPROOFCLIENT","");
103 
104  // Protocol and Host
105  fUrl.SetProtocol("proof");
106  fUrl.SetHost("__lite__");
107  fUrl.SetPort(1093);
108 
109  // User
110  if (strlen(fUrl.GetUser()) <= 0) {
111  // Get user logon name
113  if (pw) {
114  fUrl.SetUser(pw->fUser);
115  delete pw;
116  }
117  }
118  fMaster = gSystem->HostName();
119 
120  // Analysise the conffile field
121  ParseConfigField(conffile);
122 
123  // Determine the number of workers giving priority to users request.
124  // Otherwise use the system information, if available, or just start
125  // the minimal number, i.e. 2 .
126  if ((fNWorkers = GetNumberOfWorkers(url)) > 0) {
127 
128  TString stup;
129  if (gProofServ) {
130  Int_t port = gEnv->GetValue("ProofServ.XpdPort", 1093);
131  stup.Form("%s @ %s:%d ", gProofServ->GetOrdinal(), gSystem->HostName(), port);
132  }
133  Printf(" +++ Starting PROOF-Lite %swith %d workers +++", stup.Data(), fNWorkers);
134  // Init the session now
135  Init(url, conffile, confdir, loglevel, alias);
136  }
137 
138  // For final cleanup
139  if (!gROOT->GetListOfProofs()->FindObject(this))
140  gROOT->GetListOfProofs()->Add(this);
141 
142  // Still needed by the packetizers: needs to be changed
143  gProof = this;
144 }
145 
146 ////////////////////////////////////////////////////////////////////////////////
147 /// Start the PROOF environment. Starting PROOF involves either connecting
148 /// to a master server, which in turn will start a set of slave servers, or
149 /// directly starting as master server (if master = ""). For a description
150 /// of the arguments see the TProof ctor. Returns the number of started
151 /// master or slave servers, returns 0 in case of error, in which case
152 /// fValid remains false.
153 
154 Int_t TProofLite::Init(const char *, const char *conffile,
155  const char *confdir, Int_t loglevel, const char *)
156 {
158 
159  fValid = kFALSE;
160 
161  // Connected to terminal?
162  fTty = (isatty(0) == 0 || isatty(1) == 0) ? kFALSE : kTRUE;
163 
164  if (TestBit(TProof::kIsMaster)) {
165  // Fill default conf file and conf dir
166  if (!conffile || !conffile[0])
168  if (!confdir || !confdir[0])
170  } else {
171  fConfDir = confdir;
172  fConfFile = conffile;
173  }
174 
175  // The sandbox for this session
176  if (CreateSandbox() != 0) {
177  Error("Init", "could not create/assert sandbox for this session");
178  return 0;
179  }
180 
181  // UNIX path for communication with workers
182  TString sockpathdir = gEnv->GetValue("ProofLite.SockPathDir", gSystem->TempDirectory());
183  if (sockpathdir.IsNull()) sockpathdir = gSystem->TempDirectory();
184  if (sockpathdir(sockpathdir.Length()-1) == '/') sockpathdir.Remove(sockpathdir.Length()-1);
185  fSockPath.Form("%s/plite-%d", sockpathdir.Data(), gSystem->GetPid());
186  if (fSockPath.Length() > 104) {
187  // Sort of hardcoded limit length for Unix systems
188  Error("Init", "Unix socket path '%s' is too long (%d bytes):",
190  Error("Init", "use 'ProofLite.SockPathDir' to create it under a directory different"
191  " from '%s'", sockpathdir.Data());
192  return 0;
193  }
194 
195  fLogLevel = loglevel;
198  fImage = "<local>";
199  fIntHandler = 0;
200  fStatus = 0;
201  fRecvMessages = new TList;
203  fSlaveInfo = 0;
204  fChains = new TList;
205  fAvailablePackages = 0;
206  fEnabledPackages = 0;
208  fInputData = 0;
210 
213 
214  // Timeout for some collect actions
215  fCollectTimeout = gEnv->GetValue("Proof.CollectTimeout", -1);
216 
217  // Should the workers be started dynamically; default: no
219  fDynamicStartupStep = -1;
220  fDynamicStartupNMax = -1;
221  TString dynconf = gEnv->GetValue("Proof.SimulateDynamicStartup", "");
222  if (dynconf.Length() > 0) {
224  fLastPollWorkers_s = time(0);
225  // Extract parameters
226  Int_t from = 0;
227  TString p;
228  if (dynconf.Tokenize(p, from, ":"))
229  if (p.IsDigit()) fDynamicStartupStep = p.Atoi();
230  if (dynconf.Tokenize(p, from, ":"))
231  if (p.IsDigit()) fDynamicStartupNMax = p.Atoi();
232  }
233 
234 
235  fProgressDialog = 0;
237 
238  // Client logging of messages from the workers
239  fRedirLog = kFALSE;
240  if (TestBit(TProof::kIsClient)) {
241  fLogFileName = Form("%s/session-%s.log", fWorkDir.Data(), GetName());
242  if ((fLogFileW = fopen(fLogFileName.Data(), "w")) == 0)
243  Error("Init", "could not create temporary logfile %s", fLogFileName.Data());
244  if ((fLogFileR = fopen(fLogFileName.Data(), "r")) == 0)
245  Error("Init", "could not open logfile %s for reading", fLogFileName.Data());
246  }
248 
251  TString(fCacheDir).ReplaceAll("/","%").Data()));
252 
253  // Create 'queries' locker instance and lock it
256  TString(fQueryDir).ReplaceAll("/","%").Data()));
257  fQueryLock->Lock();
258  // Create the query manager
261 
262  // Apply quotas, if any
263  Int_t maxq = gEnv->GetValue("ProofLite.MaxQueriesSaved", 10);
264  if (fQMgr && fQMgr->ApplyMaxQueries(maxq) != 0)
265  Warning("Init", "problems applying fMaxQueries");
266 
267  if (InitDataSetManager() != 0)
268  Warning("Init", "problems initializing the dataset manager");
269 
270  // Status of cluster
271  fNotIdle = 0;
272 
273  // Query type
274  fSync = kTRUE;
275 
276  // List of queries
277  fQueries = 0;
278  fOtherQueries = 0;
279  fDrawQueries = 0;
280  fMaxDrawQueries = 1;
281  fSeqNum = 0;
282 
283  // Remote ID of the session
284  fSessionID = -1;
285 
286  // Part of active query
287  fWaitingSlaves = 0;
288 
289  // Make remote PROOF player
290  fPlayer = 0;
291  MakePlayer("lite");
292 
293  fFeedback = new TList;
294  fFeedback->SetOwner();
295  fFeedback->SetName("FeedbackList");
297 
298  // Sort workers by descending performance index
300  fActiveSlaves = new TList;
301  fInactiveSlaves = new TList;
302  fUniqueSlaves = new TList;
303  fAllUniqueSlaves = new TList;
304  fNonUniqueMasters = new TList;
305  fBadSlaves = new TList;
306  fAllMonitor = new TMonitor;
307  fActiveMonitor = new TMonitor;
308  fUniqueMonitor = new TMonitor;
310  fCurrentMonitor = 0;
311  fServSock = 0;
312 
315 
316  // Control how to start the workers; copy-on-write (fork) is *very*
317  // experimental and available on Unix only.
319  if (gEnv->GetValue("ProofLite.ForkStartup", 0) != 0) {
320 #ifndef WIN32
322 #else
323  Warning("Init", "fork-based workers startup is not available on Windows - ignoring");
324 #endif
325  }
326 
327  fLoadedMacros = 0;
328  if (TestBit(TProof::kIsClient)) {
329 
330  // List of directories where to look for global packages
331  TString globpack = gEnv->GetValue("Proof.GlobalPackageDirs","");
332  TProofServ::ResolveKeywords(globpack);
333  Int_t nglb = TPackMgr::RegisterGlobalPath(globpack);
334  if (gDebug > 0)
335  Info("Init", " %d global package directories registered", nglb);
336  }
337 
338  // Start workers
339  if (SetupWorkers(0) != 0) {
340  Error("Init", "problems setting up workers");
341  return 0;
342  }
343 
344  // we are now properly initialized
345  fValid = kTRUE;
346 
347  // De-activate monitor (will be activated in Collect)
349 
350  // By default go into parallel mode
351  GoParallel(-1, kFALSE);
352 
353  // Send relevant initial state to slaves
355 
356  SetActive(kFALSE);
357 
358  if (IsValid()) {
359  // Activate input handler
361  // Set PROOF to running state
363  }
364  // We register the session as a socket so that cleanup is done properly
366  gROOT->GetListOfSockets()->Add(this);
367 
368  AskParallel();
369 
370  return fActiveSlaves->GetSize();
371 }
372 ////////////////////////////////////////////////////////////////////////////////
373 /// Destructor
374 
376 {
377  // Shutdown the workers
378  RemoveWorkers(0);
379 
380  if (!(fQMgr && fQMgr->Queries() && fQMgr->Queries()->GetSize())) {
381  // needed in case fQueryDir is on NFS ?!
382  gSystem->MakeDirectory(fQueryDir+"/.delete");
383  gSystem->Exec(Form("%s %s", kRM, fQueryDir.Data()));
384  }
385 
386  // Remove lock file
387  if (fQueryLock) {
389  fQueryLock->Unlock();
390  }
391 
395 
396  // Cleanup the socket
399 }
400 
401 ////////////////////////////////////////////////////////////////////////////////
402 /// Static method to determine the number of workers giving priority to users request.
403 /// Otherwise use the system information, if available, or just start
404 /// the minimal number, i.e. 2 .
405 
407 {
408  Bool_t notify = kFALSE;
409  if (fgWrksMax == -2) {
410  // Find the max number of workers, if any
411  TString sysname = "system.rootrc";
412 #ifdef ROOTETCDIR
413  char *s = gSystem->ConcatFileName(ROOTETCDIR, sysname);
414 #else
415  TString etc = gRootDir;
416 #ifdef WIN32
417  etc += "\\etc";
418 #else
419  etc += "/etc";
420 #endif
421  char *s = gSystem->ConcatFileName(etc, sysname);
422 #endif
423  TEnv sysenv(0);
424  sysenv.ReadFile(s, kEnvGlobal);
425  fgWrksMax = sysenv.GetValue("ProofLite.MaxWorkers", -1);
426  // Notify once the user if its will is changed
427  notify = kTRUE;
428  if (s) delete[] s;
429  }
430  if (fgWrksMax == 0) {
431  ::Error("TProofLite::GetNumberOfWorkers",
432  "PROOF-Lite disabled by the system administrator: sorry!");
433  return 0;
434  }
435 
436  TString nw;
437  Int_t nWorkers = -1;
438  Bool_t urlSetting = kFALSE;
439  if (url && strlen(url)) {
440  nw = url;
441  Int_t in = nw.Index("workers=");
442  if (in != kNPOS) {
443  nw.Remove(0, in + strlen("workers="));
444  while (!nw.IsDigit())
445  nw.Remove(nw.Length()-1);
446  if (!nw.IsNull()) {
447  if ((nWorkers = nw.Atoi()) <= 0) {
448  ::Warning("TProofLite::GetNumberOfWorkers",
449  "number of workers specified by 'workers='"
450  " is non-positive: using default");
451  } else {
452  urlSetting = kFALSE;
453  }
454  }
455  }
456  }
457  if (!urlSetting && fgProofEnvList) {
458  // Check PROOF_NWORKERS
459  TNamed *nm = (TNamed *) fgProofEnvList->FindObject("PROOF_NWORKERS");
460  if (nm) {
461  nw = nm->GetTitle();
462  if (nw.IsDigit()) {
463  if ((nWorkers = nw.Atoi()) == 0) {
464  ::Warning("TProofLite::GetNumberOfWorkers",
465  "number of workers specified by 'workers='"
466  " is non-positive: using default");
467  }
468  }
469  }
470  }
471  if (nWorkers <= 0) {
472  nWorkers = gEnv->GetValue("ProofLite.Workers", -1);
473  if (nWorkers <= 0) {
474  SysInfo_t si;
475  if (gSystem->GetSysInfo(&si) == 0 && si.fCpus > 2) {
476  nWorkers = si.fCpus;
477  } else {
478  // Two workers by default
479  nWorkers = 2;
480  }
481  if (notify) notify = kFALSE;
482  }
483  }
484  // Apply the max, if any
485  if (fgWrksMax > 0 && fgWrksMax < nWorkers) {
486  if (notify)
487  ::Warning("TProofLite::GetNumberOfWorkers", "number of PROOF-Lite workers limited by"
488  " the system administrator to %d", fgWrksMax);
489  nWorkers = fgWrksMax;
490  }
491 
492  // Done
493  return nWorkers;
494 }
495 
496 ////////////////////////////////////////////////////////////////////////////////
497 /// Start up PROOF workers.
498 
500 {
501  // Create server socket on the assigned UNIX sock path
502  if (!fServSock) {
503  if ((fServSock = new TServerSocket(fSockPath))) {
505  // Remove from the list so that cleanup can be done in the correct order
506  gROOT->GetListOfSockets()->Remove(fServSock);
507  }
508  }
509  if (!fServSock || !fServSock->IsValid()) {
510  Error("SetupWorkers",
511  "unable to create server socket for internal communications");
513  return -1;
514  }
515 
516  // Create a monitor and add the socket to it
517  TMonitor *mon = new TMonitor;
518  mon->Add(fServSock);
519 
520  TList started;
521  TSlave *wrk = 0;
522  Int_t nWrksDone = 0, nWrksTot = -1;
523  TString fullord;
524 
525  if (opt == 0) {
526  nWrksTot = fForkStartup ? 1 : fNWorkers;
527  // Now we create the worker applications which will call us back to finalize
528  // the setup
529  Int_t ord = 0;
530  for (; ord < nWrksTot; ord++) {
531 
532  // Ordinal for this worker server
533  const char *o = (gProofServ) ? gProofServ->GetOrdinal() : "0";
534  fullord.Form("%s.%d", o, ord);
535 
536  // Create environment files
537  SetProofServEnv(fullord);
538 
539  // Create worker server and add to the list
540  if ((wrk = CreateSlave("lite", fullord, 100, fImage, fWorkDir)))
541  started.Add(wrk);
542 
543  // Notify
544  NotifyStartUp("Opening connections to workers", ++nWrksDone, nWrksTot);
545 
546  } //end of worker loop
547  } else {
548  if (!fForkStartup) {
549  Warning("SetupWorkers", "standard startup: workers already started");
550  return -1;
551  }
552  nWrksTot = fNWorkers - 1;
553  // Now we create the worker applications which will call us back to finalize
554  // the setup
555  TString clones;
556  Int_t ord = 0;
557  for (; ord < nWrksTot; ord++) {
558 
559  // Ordinal for this worker server
560  const char *o = (gProofServ) ? gProofServ->GetOrdinal() : "0";
561  fullord.Form("%s.%d", o, ord + 1);
562  if (!clones.IsNull()) clones += " ";
563  clones += fullord;
564 
565  // Create worker server and add to the list
566  if ((wrk = CreateSlave("lite", fullord, -1, fImage, fWorkDir)))
567  started.Add(wrk);
568 
569  // Notify
570  NotifyStartUp("Opening connections to workers", ++nWrksDone, nWrksTot);
571 
572  } //end of worker loop
573 
574  // Send the request
576  m << clones;
577  Broadcast(m, kActive);
578  }
579 
580  // Wait for call backs
581  nWrksDone = 0;
582  nWrksTot = started.GetSize();
583  Int_t nSelects = 0;
584  Int_t to = gEnv->GetValue("ProofLite.StartupTimeOut", 5) * 1000;
585  while (started.GetSize() > 0 && nSelects < nWrksTot) {
586 
587  // Wait for activity on the socket for max 5 secs
588  TSocket *xs = mon->Select(to);
589 
590  // Count attempts and check
591  nSelects++;
592  if (xs == (TSocket *) -1) continue;
593 
594  // Get the connection
595  TSocket *s = fServSock->Accept();
596  if (s && s->IsValid()) {
597  // Receive ordinal
598  TMessage *msg = 0;
599  if (s->Recv(msg) < 0) {
600  Warning("SetupWorkers", "problems receiving message from accepted socket!");
601  } else {
602  if (msg) {
603  TString ord;
604  *msg >> ord;
605  // Find who is calling back
606  if ((wrk = (TSlave *) started.FindObject(ord))) {
607  // Remove it from the started list
608  started.Remove(wrk);
609 
610  // Assign tis socket the selected worker
611  wrk->SetSocket(s);
612  // Remove socket from global TROOT socket list. Only the TProof object,
613  // representing all worker sockets, will be added to this list. This will
614  // ensure the correct termination of all proof servers in case the
615  // root session terminates.
617  gROOT->GetListOfSockets()->Remove(s);
618  }
619  if (wrk->IsValid()) {
620  // Set the input handler
621  wrk->SetInputHandler(new TProofInputHandler(this, wrk->GetSocket()));
622  // Set fParallel to 1 for workers since they do not
623  // report their fParallel with a LOG_DONE message
624  wrk->fParallel = 1;
625  // Finalize setup of the server
626  wrk->SetupServ(TSlave::kSlave, 0);
627  }
628 
629  // Monitor good workers
630  fSlaves->Add(wrk);
631  if (wrk->IsValid()) {
632  if (opt == 1) fActiveSlaves->Add(wrk);
633  fAllMonitor->Add(wrk->GetSocket());
634  // Record also in the list for termination
635  if (startedWorkers) startedWorkers->Add(wrk);
636  // Notify startup operations
637  NotifyStartUp("Setting up worker servers", ++nWrksDone, nWrksTot);
638  } else {
639  // Flag as bad
640  fBadSlaves->Add(wrk);
641  }
642  }
643  } else {
644  Warning("SetupWorkers", "received empty message from accepted socket!");
645  }
646  }
647  }
648  }
649 
650  // Cleanup the monitor and the server socket
651  mon->DeActivateAll();
652  delete mon;
653 
654  // Create Progress dialog, if needed
655  if (!gROOT->IsBatch() && !fProgressDialog) {
656  if ((fProgressDialog =
657  gROOT->GetPluginManager()->FindHandler("TProofProgressDialog")))
658  if (fProgressDialog->LoadPlugin() == -1)
659  fProgressDialog = 0;
660  }
661 
662  if (opt == 1) {
663  // Collect replies
664  Collect(kActive);
665  // Update group view
666  SendGroupView();
667  // By default go into parallel mode
668  SetParallel(-1, 0);
669  }
670  // Done
671  return 0;
672 }
673 
674 ////////////////////////////////////////////////////////////////////////////////
675 /// Notify setting-up operation message
676 
677 void TProofLite::NotifyStartUp(const char *action, Int_t done, Int_t tot)
678 {
679  Int_t frac = (Int_t) (done*100.)/tot;
680  char msg[512] = {0};
681  if (frac >= 100) {
682  snprintf(msg, 512, "%s: OK (%d workers) \n",
683  action, tot);
684  } else {
685  snprintf(msg, 512, "%s: %d out of %d (%d %%)\r",
686  action, done, tot, frac);
687  }
688  fprintf(stderr,"%s", msg);
689 }
690 
691 ////////////////////////////////////////////////////////////////////////////////
692 /// Create environment files for worker 'ord'
693 
695 {
696  // Check input
697  if (!ord || strlen(ord) <= 0) {
698  Error("SetProofServEnv", "ordinal string undefined");
699  return -1;
700  }
701 
702  // ROOT env file
703  TString rcfile(Form("%s/worker-%s.rootrc", fWorkDir.Data(), ord));
704  FILE *frc = fopen(rcfile.Data(), "w");
705  if (!frc) {
706  Error("SetProofServEnv", "cannot open rc file %s", rcfile.Data());
707  return -1;
708  }
709 
710  // The session working dir depends on the role
711  fprintf(frc,"# The session working dir\n");
712  fprintf(frc,"ProofServ.SessionDir: %s/worker-%s\n", fWorkDir.Data(), ord);
713 
714  // The session unique tag
715  fprintf(frc,"# Session tag\n");
716  fprintf(frc,"ProofServ.SessionTag: %s\n", GetName());
717 
718  // Log / Debug level
719  fprintf(frc,"# Proof Log/Debug level\n");
720  fprintf(frc,"Proof.DebugLevel: %d\n", gDebug);
721 
722  // Ordinal number
723  fprintf(frc,"# Ordinal number\n");
724  fprintf(frc,"ProofServ.Ordinal: %s\n", ord);
725 
726  // ROOT Version tag
727  fprintf(frc,"# ROOT Version tag\n");
728  fprintf(frc,"ProofServ.RootVersionTag: %s\n", gROOT->GetVersion());
729 
730  // Work dir
731  TString sandbox = fSandbox;
732  if (GetSandbox(sandbox, kFALSE, "ProofServ.Sandbox") != 0)
733  Warning("SetProofServEnv", "problems getting sandbox string for worker");
734  fprintf(frc,"# Users sandbox\n");
735  fprintf(frc, "ProofServ.Sandbox: %s\n", sandbox.Data());
736 
737  // Cache dir
738  fprintf(frc,"# Users cache\n");
739  fprintf(frc, "ProofServ.CacheDir: %s\n", fCacheDir.Data());
740 
741  // Package dir
742  fprintf(frc,"# Users packages\n");
743  fprintf(frc, "ProofServ.PackageDir: %s\n", fPackMgr->GetDir());
744 
745  // Image
746  fprintf(frc,"# Server image\n");
747  fprintf(frc, "ProofServ.Image: %s\n", fImage.Data());
748 
749  // Set Open socket
750  fprintf(frc,"# Open socket\n");
751  fprintf(frc, "ProofServ.OpenSock: %s\n", fSockPath.Data());
752 
753  // Client Protocol
754  fprintf(frc,"# Client Protocol\n");
755  fprintf(frc, "ProofServ.ClientVersion: %d\n", kPROOF_Protocol);
756 
757  // ROOT env file created
758  fclose(frc);
759 
760  // System env file
761  TString envfile(Form("%s/worker-%s.env", fWorkDir.Data(), ord));
762  FILE *fenv = fopen(envfile.Data(), "w");
763  if (!fenv) {
764  Error("SetProofServEnv", "cannot open env file %s", envfile.Data());
765  return -1;
766  }
767  // ROOTSYS
768 #ifdef R__HAVE_CONFIG
769  fprintf(fenv, "export ROOTSYS=%s\n", ROOTPREFIX);
770 #else
771  fprintf(fenv, "export ROOTSYS=%s\n", gSystem->Getenv("ROOTSYS"));
772 #endif
773  // Conf dir
774 #ifdef R__HAVE_CONFIG
775  fprintf(fenv, "export ROOTCONFDIR=%s\n", ROOTETCDIR);
776 #else
777  fprintf(fenv, "export ROOTCONFDIR=%s\n", gSystem->Getenv("ROOTSYS"));
778 #endif
779  // TMPDIR
780  fprintf(fenv, "export TMPDIR=%s\n", gSystem->TempDirectory());
781  // Log file in the log dir
782  TString logfile(Form("%s/worker-%s.log", fWorkDir.Data(), ord));
783  fprintf(fenv, "export ROOTPROOFLOGFILE=%s\n", logfile.Data());
784  // RC file
785  fprintf(fenv, "export ROOTRCFILE=%s\n", rcfile.Data());
786  // ROOT version tag (needed in building packages)
787  fprintf(fenv, "export ROOTVERSIONTAG=%s\n", gROOT->GetVersion());
788  // This flag can be used to identify the type of worker; for example, in BUILD.sh or SETUP.C ...
789  fprintf(fenv, "export ROOTPROOFLITE=%d\n", fNWorkers);
790  // Local files are on the local file system
791  fprintf(fenv, "export LOCALDATASERVER=\"file://\"\n");
792  // Set the user envs
793  if (fgProofEnvList) {
794  TString namelist;
795  TIter nxenv(fgProofEnvList);
796  TNamed *env = 0;
797  while ((env = (TNamed *)nxenv())) {
798  TString senv(env->GetTitle());
799  ResolveKeywords(senv, ord, logfile.Data());
800  fprintf(fenv, "export %s=%s\n", env->GetName(), senv.Data());
801  if (namelist.Length() > 0)
802  namelist += ',';
803  namelist += env->GetName();
804  }
805  fprintf(fenv, "export PROOF_ALLVARS=%s\n", namelist.Data());
806  }
807 
808  // System env file created
809  fclose(fenv);
810 
811  // Done
812  return 0;
813 }
814 
815 ////////////////////////////////////////////////////////////////////////////////
816 /// Resolve some keywords in 's'
817 /// <logfilewrk>, <user>, <rootsys>, <cpupin>
818 
819 void TProofLite::ResolveKeywords(TString &s, const char *ord,
820  const char *logfile)
821 {
822  if (!logfile) return;
823 
824  // Log file
825  if (s.Contains("<logfilewrk>") && logfile) {
826  TString lfr(logfile);
827  if (lfr.EndsWith(".log")) lfr.Remove(lfr.Last('.'));
828  s.ReplaceAll("<logfilewrk>", lfr.Data());
829  }
830 
831  // user
832  if (gSystem->Getenv("USER") && s.Contains("<user>")) {
833  s.ReplaceAll("<user>", gSystem->Getenv("USER"));
834  }
835 
836  // rootsys
837  if (gSystem->Getenv("ROOTSYS") && s.Contains("<rootsys>")) {
838  s.ReplaceAll("<rootsys>", gSystem->Getenv("ROOTSYS"));
839  }
840 
841  // cpupin: pin to this CPU num (from 0 to ncpus-1)
842  if (s.Contains("<cpupin>")) {
843  TString o = ord;
844  Int_t n = o.Index('.');
845  if (n != kNPOS) {
846 
847  o.Remove(0, n+1);
848  n = o.Atoi(); // n is ord
849 
850  TString cpuPinList;
851  {
852  const TList *envVars = GetEnvVars();
853  TNamed *var;
854  if (envVars) {
855  var = dynamic_cast<TNamed *>(envVars->FindObject("PROOF_SLAVE_CPUPIN_ORDER"));
856  if (var) cpuPinList = var->GetTitle();
857  }
858  }
859 
860  UInt_t nCpus = 1;
861  {
862  SysInfo_t si;
863  if (gSystem->GetSysInfo(&si) == 0 && (si.fCpus > 0))
864  nCpus = si.fCpus;
865  else nCpus = 1; // fallback
866  }
867 
868  if (cpuPinList.IsNull() || (cpuPinList == "*")) {
869  // Use processors in order
870  n = n % nCpus;
871  }
872  else {
873  // Use processors in user's order
874  // n is now the ordinal, converting to idx
875  n = n % (cpuPinList.CountChar('+')+1);
876  TString tok;
877  Ssiz_t from = 0;
878  for (Int_t i=0; cpuPinList.Tokenize(tok, from, "\\+"); i++) {
879  if (i == n) {
880  n = (tok.Atoi() % nCpus);
881  break;
882  }
883  }
884  }
885 
886  o.Form("%d", n);
887  }
888  else {
889  o = "0"; // should not happen
890  }
891  s.ReplaceAll("<cpupin>", o);
892  }
893 }
894 
895 ////////////////////////////////////////////////////////////////////////////////
896 /// Create the sandbox for this session
897 
899 {
900  // Make sure the sandbox area exist and is writable
901  if (GetSandbox(fSandbox, kTRUE, "ProofLite.Sandbox") != 0) return -1;
902 
903  // Package Manager
904  TString packdir = gEnv->GetValue("Proof.PackageDir", "");
905  if (packdir.IsNull())
906  packdir.Form("%s/%s", fSandbox.Data(), kPROOF_PackDir);
907  if (AssertPath(packdir, kTRUE) != 0) return -1;
908  fPackMgr = new TPackMgr(packdir);
909 
910  // Cache Dir
911  fCacheDir = gEnv->GetValue("Proof.CacheDir", "");
912  if (fCacheDir.IsNull())
914  if (AssertPath(fCacheDir, kTRUE) != 0) return -1;
915 
916  // Data Set Dir
917  fDataSetDir = gEnv->GetValue("Proof.DataSetDir", "");
918  if (fDataSetDir.IsNull())
920  if (AssertPath(fDataSetDir, kTRUE) != 0) return -1;
921 
922  // Session unique tag (name of this TProof instance)
923  TString stag;
924  stag.Form("%s-%d-%d", gSystem->HostName(), (int)time(0), gSystem->GetPid());
925  SetName(stag.Data());
926 
927  Int_t subpath = gEnv->GetValue("ProofLite.SubPath", 1);
928  // Subpath for this session in the fSandbox (<sandbox>/path-to-working-dir)
929  TString sessdir;
930  if (subpath != 0) {
931  sessdir = gSystem->WorkingDirectory();
932  sessdir.ReplaceAll(gSystem->HomeDirectory(),"");
933  sessdir.ReplaceAll("/","-");
934  sessdir.Replace(0,1,"/",1);
935  sessdir.Insert(0, fSandbox.Data());
936  } else {
937  // USe the sandbox
938  sessdir = fSandbox;
939  }
940 
941  // Session working and queries dir
942  fWorkDir.Form("%s/session-%s", sessdir.Data(), stag.Data());
943  if (AssertPath(fWorkDir, kTRUE) != 0) return -1;
944 
945  // Create symlink to the last session
946  TString lastsess;
947  lastsess.Form("%s/last-lite-session", sessdir.Data());
948  gSystem->Unlink(lastsess);
949  gSystem->Symlink(fWorkDir, lastsess);
950 
951  // Queries Dir: local to the working dir, unless required differently
952  fQueryDir = gEnv->GetValue("Proof.QueryDir", "");
953  if (fQueryDir.IsNull())
954  fQueryDir.Form("%s/%s", sessdir.Data(), kPROOF_QueryDir);
955  if (AssertPath(fQueryDir, kTRUE) != 0) return -1;
956 
957  // Cleanup old sessions dirs
958  CleanupSandbox();
959 
960  // Done
961  return 0;
962 }
963 
964 ////////////////////////////////////////////////////////////////////////////////
965 /// Print status of PROOF-Lite cluster.
966 
967 void TProofLite::Print(Option_t *option) const
968 {
969  TString ord;
970  if (gProofServ) ord.Form("%s ", gProofServ->GetOrdinal());
971  if (IsParallel())
972  Printf("*** PROOF-Lite cluster %s(parallel mode, %d workers):", ord.Data(), GetParallel());
973  else
974  Printf("*** PROOF-Lite cluster %s(sequential mode)", ord.Data());
975 
976  if (gProofServ) {
977  TString url(gSystem->HostName());
978  // Add port to URL, if defined
979  Int_t port = gEnv->GetValue("ProofServ.XpdPort", 1093);
980  if (port > -1) url.Form("%s:%d",gSystem->HostName(), port);
981  Printf("URL: %s", url.Data());
982  } else {
983  Printf("Host name: %s", gSystem->HostName());
984  }
985  Printf("User: %s", GetUser());
986  TString ver(gROOT->GetVersion());
987  ver += TString::Format("|%s", gROOT->GetGitCommit());
988  if (gSystem->Getenv("ROOTVERSIONTAG"))
989  ver += TString::Format("|%s", gSystem->Getenv("ROOTVERSIONTAG"));
990  Printf("ROOT version|rev|tag: %s", ver.Data());
991  Printf("Architecture-Compiler: %s-%s", gSystem->GetBuildArch(),
993  Printf("Protocol version: %d", GetClientProtocol());
994  Printf("Working directory: %s", gSystem->WorkingDirectory());
995  Printf("Communication path: %s", fSockPath.Data());
996  Printf("Log level: %d", GetLogLevel());
997  Printf("Number of workers: %d", GetNumberOfSlaves());
998  Printf("Number of active workers: %d", GetNumberOfActiveSlaves());
999  Printf("Number of unique workers: %d", GetNumberOfUniqueSlaves());
1000  Printf("Number of inactive workers: %d", GetNumberOfInactiveSlaves());
1001  Printf("Number of bad workers: %d", GetNumberOfBadSlaves());
1002  Printf("Total MB's processed: %.2f", float(GetBytesRead())/(1024*1024));
1003  Printf("Total real time used (s): %.3f", GetRealTime());
1004  Printf("Total CPU time used (s): %.3f", GetCpuTime());
1005  if (TString(option).Contains("a", TString::kIgnoreCase) && GetNumberOfSlaves()) {
1006  Printf("List of workers:");
1007  TIter nextslave(fSlaves);
1008  while (TSlave* sl = dynamic_cast<TSlave*>(nextslave())) {
1009  if (sl->IsValid())
1010  sl->Print(option);
1011  }
1012  }
1013 }
1014 
1015 ////////////////////////////////////////////////////////////////////////////////
1016 /// Create a TProofQueryResult instance for this query.
1017 
1019  Long64_t fst, TDSet *dset,
1020  const char *selec)
1021 {
1022  // Increment sequential number
1023  Int_t seqnum = -1;
1024  if (fQMgr) {
1026  seqnum = fQMgr->SeqNum();
1027  }
1028 
1029  // Create the instance and add it to the list
1030  TProofQueryResult *pqr = new TProofQueryResult(seqnum, opt,
1031  fPlayer->GetInputList(), nent,
1032  fst, dset, selec,
1033  (dset ? dset->GetEntryList() : 0));
1034  // Title is the session identifier
1035  pqr->SetTitle(GetName());
1036 
1037  return pqr;
1038 }
1039 
1040 ////////////////////////////////////////////////////////////////////////////////
1041 /// Set query in running state.
1042 
1044 {
1045  // Record current position in the log file at start
1046  fflush(fLogFileW);
1047  Int_t startlog = lseek(fileno(fLogFileW), (off_t) 0, SEEK_END);
1048 
1049  // Add some header to logs
1050  Printf(" ");
1051  Info("SetQueryRunning", "starting query: %d", pq->GetSeqNum());
1052 
1053  // Build the list of loaded PAR packages
1054  TString parlist = "";
1055  fPackMgr->GetEnabledPackages(parlist);
1056 
1057  // Set in running state
1058  pq->SetRunning(startlog, parlist, GetParallel());
1059 
1060  // Bytes and CPU at start (we will calculate the differential at end)
1061  AskStatistics();
1063 }
1064 
1065 ////////////////////////////////////////////////////////////////////////////////
1066 /// Execute the specified drawing action on a data set (TDSet).
1067 /// Event- or Entry-lists should be set in the data set object using
1068 /// TDSet::SetEntryList.
1069 /// Returns -1 in case of error or number of selected events otherwise.
1070 
1071 Long64_t TProofLite::DrawSelect(TDSet *dset, const char *varexp,
1072  const char *selection, Option_t *option,
1074 {
1075  if (!IsValid()) return -1;
1076 
1077  // Make sure that asynchronous processing is not active
1078  if (!IsIdle()) {
1079  Info("DrawSelect","not idle, asynchronous Draw not supported");
1080  return -1;
1081  }
1082  TString opt(option);
1083  Int_t idx = opt.Index("ASYN", 0, TString::kIgnoreCase);
1084  if (idx != kNPOS)
1085  opt.Replace(idx,4,"");
1086 
1087  // Fill the internal variables
1088  fVarExp = varexp;
1089  fSelection = selection;
1090 
1091  return Process(dset, "draw:", opt, nentries, first);
1092 }
1093 
1094 ////////////////////////////////////////////////////////////////////////////////
1095 /// Process a data set (TDSet) using the specified selector (.C) file.
1096 /// Entry- or event-lists should be set in the data set object using
1097 /// TDSet::SetEntryList.
1098 /// The return value is -1 in case of error and TSelector::GetStatus() in
1099 /// in case of success.
1100 
1101 Long64_t TProofLite::Process(TDSet *dset, const char *selector, Option_t *option,
1103 {
1104  // For the time being cannot accept other queries if not idle, even if in async
1105  // mode; needs to set up an event handler to manage that
1106 
1107  TString opt(option), optfb, outfile;
1108  // Enable feedback, if required
1109  if (opt.Contains("fb=") || opt.Contains("feedback=")) SetFeedback(opt, optfb, 0);
1110  // Define output file, either from 'opt' or the default one
1111  if (HandleOutputOptions(opt, outfile, 0) != 0) return -1;
1112 
1113  // Resolve query mode
1114  fSync = (GetQueryMode(opt) == kSync);
1115  if (!fSync) {
1116  Info("Process","asynchronous mode not yet supported in PROOF-Lite");
1117  return -1;
1118  }
1119 
1120  if (!IsIdle()) {
1121  // Notify submission
1122  Info("Process", "not idle: cannot accept queries");
1123  return -1;
1124  }
1125 
1126  // Cleanup old temporary datasets
1127  if (IsIdle() && fRunningDSets && fRunningDSets->GetSize() > 0) {
1129  fRunningDSets->Delete();
1130  }
1131 
1132  if (!IsValid() || !fQMgr || !fPlayer) {
1133  Error("Process", "invalid sesion or query-result manager undefined!");
1134  return -1;
1135  }
1136 
1137  // Make sure that all enabled workers get some work, unless stated
1138  // differently
1139  if (!fPlayer->GetInputList()->FindObject("PROOF_MaxSlavesPerNode"))
1140  SetParameter("PROOF_MaxSlavesPerNode", (Long_t)0);
1141 
1142  Bool_t hasNoData = (!dset || (dset && dset->TestBit(TDSet::kEmpty))) ? kTRUE : kFALSE;
1143 
1144  // If just a name was given to identify the dataset, retrieve it from the
1145  // local files
1146  // Make sure the dataset contains the information needed
1147  TString emsg;
1148  if ((!hasNoData) && dset->GetListOfElements()->GetSize() == 0) {
1149  if (TProof::AssertDataSet(dset, fPlayer->GetInputList(), fDataSetManager, emsg) != 0) {
1150  Error("Process", "from AssertDataSet: %s", emsg.Data());
1151  return -1;
1152  }
1153  if (dset->GetListOfElements()->GetSize() == 0) {
1154  Error("Process", "no files to process!");
1155  return -1;
1156  }
1157  } else if (hasNoData) {
1158  // Check if we are required to process with TPacketizerFile a registered dataset
1159  TNamed *ftp = dynamic_cast<TNamed *>(fPlayer->GetInputList()->FindObject("PROOF_FilesToProcess"));
1160  if (ftp) {
1161  TString dsn(ftp->GetTitle());
1162  if (!dsn.Contains(":") || dsn.BeginsWith("dataset:")) {
1163  dsn.ReplaceAll("dataset:", "");
1164  // Make sure we have something in input and a dataset manager
1165  if (!fDataSetManager) {
1166  emsg.Form("dataset manager not initialized!");
1167  } else {
1168  TFileCollection *fc = 0;
1169  // Get the dataset
1170  if (!(fc = fDataSetManager->GetDataSet(dsn))) {
1171  emsg.Form("requested dataset '%s' does not exists", dsn.Data());
1172  } else {
1173  TMap *fcmap = TProofServ::GetDataSetNodeMap(fc, emsg);
1174  if (fcmap) {
1175  fPlayer->GetInputList()->Remove(ftp);
1176  delete ftp;
1177  fcmap->SetOwner(kTRUE);
1178  fcmap->SetName("PROOF_FilesToProcess");
1179  fPlayer->GetInputList()->Add(fcmap);
1180  }
1181  }
1182  }
1183  if (!emsg.IsNull()) {
1184  Error("HandleProcess", "%s", emsg.Data());
1185  return -1;
1186  }
1187  }
1188  }
1189  }
1190 
1191  TString selec(selector), varexp, selection, objname;
1192  // If a draw query, extract the relevant info
1193  if (selec.BeginsWith("draw:")) {
1194  varexp = fVarExp;
1195  selection = fSelection;
1196  // Decode now the expression
1197  if (fPlayer->GetDrawArgs(varexp, selection, opt, selec, objname) != 0) {
1198  Error("Process", "draw query: error parsing arguments '%s', '%s', '%s'",
1199  varexp.Data(), selection.Data(), opt.Data());
1200  return -1;
1201  }
1202  }
1203 
1204  // Create instance of query results (the data set is added after Process)
1205  TProofQueryResult *pq = MakeQueryResult(nentries, opt, first, 0, selec);
1206 
1207  // Check if queries must be saved into files
1208  // Automatic saving is controlled by ProofLite.AutoSaveQueries
1209  Bool_t savequeries =
1210  (!strcmp(gEnv->GetValue("ProofLite.AutoSaveQueries", "off"), "on")) ? kTRUE : kFALSE;
1211 
1212  // Keep queries in memory and how many (-1 = all, 0 = none, ...)
1213  Int_t memqueries = gEnv->GetValue("ProofLite.MaxQueriesMemory", 1);
1214 
1215  // If not a draw action add the query to the main list
1216  if (!(pq->IsDraw())) {
1217  if (fQMgr->Queries()) {
1218  if (memqueries != 0) fQMgr->Queries()->Add(pq);
1219  if (memqueries >= 0 && fQMgr->Queries()->GetSize() > memqueries) {
1220  // Remove oldest
1221  TObject *qfst = fQMgr->Queries()->First();
1222  fQMgr->Queries()->Remove(qfst);
1223  delete qfst;
1224  }
1225  }
1226  // Also save it to queries dir
1227  if (savequeries) fQMgr->SaveQuery(pq);
1228  }
1229 
1230  // Set the query number
1231  fSeqNum = pq->GetSeqNum();
1232 
1233  // Set in running state
1234  SetQueryRunning(pq);
1235 
1236  // Save to queries dir, if not standard draw
1237  if (!(pq->IsDraw())) {
1238  if (savequeries) fQMgr->SaveQuery(pq);
1239  } else {
1241  }
1242 
1243  // Start or reset the progress dialog
1244  if (!gROOT->IsBatch()) {
1245  Int_t dsz = (dset && dset->GetListOfElements()) ? dset->GetListOfElements()->GetSize() : -1;
1246  if (fProgressDialog &&
1248  if (!fProgressDialogStarted) {
1249  fProgressDialog->ExecPlugin(5, this, selec.Data(), dsz,
1250  first, nentries);
1252  } else {
1253  ResetProgressDialog(selec.Data(), dsz, first, nentries);
1254  }
1255  }
1257  }
1258 
1259  // Add query results to the player lists
1260  if (!(pq->IsDraw()))
1261  fPlayer->AddQueryResult(pq);
1262 
1263  // Set query currently processed
1264  fPlayer->SetCurrentQuery(pq);
1265 
1266  // Make sure the unique query tag is available as TNamed object in the
1267  // input list so that it can be used in TSelectors for monitoring
1268  TNamed *qtag = (TNamed *) fPlayer->GetInputList()->FindObject("PROOF_QueryTag");
1269  if (qtag) {
1270  qtag->SetTitle(Form("%s:%s",pq->GetTitle(),pq->GetName()));
1271  } else {
1272  TObject *o = fPlayer->GetInputList()->FindObject("PROOF_QueryTag");
1273  if (o) fPlayer->GetInputList()->Remove(o);
1274  fPlayer->AddInput(new TNamed("PROOF_QueryTag",
1275  Form("%s:%s",pq->GetTitle(),pq->GetName())));
1276  }
1277 
1278  // Set PROOF to running state
1280 
1281  // deactivate the default application interrupt handler
1282  // ctrl-c's will be forwarded to PROOF to stop the processing
1283  TSignalHandler *sh = 0;
1284  if (fSync) {
1285  if (gApplication)
1287  }
1288 
1289  // Make sure we get a fresh result
1290  fOutputList.Clear();
1291 
1292  // Start the additional workers now if using fork-based startup
1293  TList *startedWorkers = 0;
1294  if (fForkStartup) {
1295  startedWorkers = new TList;
1296  startedWorkers->SetOwner(kFALSE);
1297  SetupWorkers(1, startedWorkers);
1298  }
1299 
1300  // This is the end of preparation
1301  fQuerySTW.Reset();
1302 
1303  Long64_t rv = 0;
1304  if (!(pq->IsDraw())) {
1305  if (selector && strlen(selector)) {
1306  rv = fPlayer->Process(dset, selec, opt, nentries, first);
1307  } else {
1308  rv = fPlayer->Process(dset, fSelector, opt, nentries, first);
1309  }
1310  } else {
1311  rv = fPlayer->DrawSelect(dset, varexp, selection, opt, nentries, first);
1312  }
1313 
1314  // This is the end of merging
1315  fQuerySTW.Stop();
1316  Float_t rt = fQuerySTW.RealTime();
1317  // Update the query content
1318  TQueryResult *qr = GetQueryResult();
1319  if (qr) {
1320  qr->SetTermTime(rt);
1321  // Preparation time is always null in PROOF-Lite
1322  }
1323 
1324  // Disable feedback, if required
1325  if (!optfb.IsNull()) SetFeedback(opt, optfb, 1);
1326 
1327  if (fSync) {
1328 
1329  // Terminate additional workers if using fork-based startup
1330  if (fForkStartup && startedWorkers) {
1331  RemoveWorkers(startedWorkers);
1332  SafeDelete(startedWorkers);
1333  }
1334 
1335  // reactivate the default application interrupt handler
1336  if (sh)
1338 
1339  // Return number of events processed
1342  ? kTRUE : kFALSE;
1343  if (abort) fPlayer->StopProcess(kTRUE);
1344  Emit("StopProcess(Bool_t)", abort);
1345  }
1346 
1347  // In PROOFLite this has to be done once only in TProofLite::Process
1349  // If the last object, notify the GUI that the result arrived
1350  QueryResultReady(Form("%s:%s", pq->GetTitle(), pq->GetName()));
1351  // Processing is over
1352  UpdateDialog();
1353 
1354  // Save the data set into the TQueryResult (should be done after Process to avoid
1355  // improper deletion during collection)
1356  if (rv == 0 && dset && !dset->TestBit(TDSet::kEmpty) && pq->GetInputList()) {
1357  pq->GetInputList()->Add(dset);
1358  if (dset->GetEntryList())
1359  pq->GetInputList()->Add(dset->GetEntryList());
1360  }
1361 
1362  // Register any dataset produced during this processing, if required
1364  TNamed *psr = (TNamed *) fPlayer->GetOutputList()->FindObject("PROOFSERV_RegisterDataSet");
1365  if (psr) {
1366  TString err;
1368  fPlayer->GetOutputList(), fDataSetManager, err) != 0)
1369  Warning("ProcessNext", "problems registering produced datasets: %s", err.Data());
1370  fPlayer->GetOutputList()->Remove(psr);
1371  delete psr;
1372  }
1373  }
1374 
1375  // Complete filling of the TQueryResult instance
1376  AskStatistics();
1377  if (!(pq->IsDraw())) {
1378  if (fQMgr->FinalizeQuery(pq, this, fPlayer)) {
1379  if (savequeries) fQMgr->SaveQuery(pq, -1);
1380  }
1381  }
1382 
1383  // Remove aborted queries from the list
1386  if (fQMgr) fQMgr->RemoveQuery(pq);
1387  } else {
1388  // If the last object, notify the GUI that the result arrived
1389  QueryResultReady(Form("%s:%s", pq->GetTitle(), pq->GetName()));
1390  // Keep in memory only light info about a query
1391  if (!(pq->IsDraw()) && memqueries >= 0) {
1392  if (fQMgr && fQMgr->Queries()) {
1393  TQueryResult *pqr = pq->CloneInfo();
1394  if (pqr) fQMgr->Queries()->Add(pqr);
1395  // Remove from the fQueries list
1396  fQMgr->Queries()->Remove(pq);
1397  }
1398  }
1399  // To get the prompt back
1400  TString msg;
1401  msg.Form("Lite-0: all output objects have been merged ");
1402  fprintf(stderr, "%s\n", msg.Data());
1403  }
1404  // Save the performance info, if required
1405  if (!fPerfTree.IsNull()) {
1406  if (SavePerfTree() != 0) Error("Process", "saving performance info ...");
1407  // Must be re-enabled each time
1408  SetPerfTree(0);
1409  }
1410  }
1411  // Finalise output file settings (opt is ignored in here)
1412  if (HandleOutputOptions(opt, outfile, 1) != 0) return -1;
1413 
1414  // Retrieve status from the output list
1415  if (rv >= 0) {
1416  TParameter<Long64_t> *sst =
1417  (TParameter<Long64_t> *) fOutputList.FindObject("PROOF_SelectorStatus");
1418  if (sst) rv = sst->GetVal();
1419  }
1420 
1421 
1422  // Done
1423  return rv;
1424 }
1425 
1426 ////////////////////////////////////////////////////////////////////////////////
1427 /// Initialize the dataset manager from directives or from defaults
1428 /// Return 0 on success, -1 on failure
1429 
1431 {
1432  fDataSetManager = 0;
1433 
1434  // Default user and group
1435  TString user("???"), group("default");
1436  UserGroup_t *pw = gSystem->GetUserInfo();
1437  if (pw) {
1438  user = pw->fUser;
1439  delete pw;
1440  }
1441 
1442  // Dataset manager instance via plug-in
1443  TPluginHandler *h = 0;
1444  TString dsm = gEnv->GetValue("Proof.DataSetManager", "");
1445  if (!dsm.IsNull()) {
1446  // Get plugin manager to load the appropriate TDataSetManager
1447  if (gROOT->GetPluginManager()) {
1448  // Find the appropriate handler
1449  h = gROOT->GetPluginManager()->FindHandler("TDataSetManager", dsm);
1450  if (h && h->LoadPlugin() != -1) {
1451  // make instance of the dataset manager
1452  fDataSetManager =
1453  reinterpret_cast<TDataSetManager*>(h->ExecPlugin(3, group.Data(),
1454  user.Data(), dsm.Data()));
1455  }
1456  }
1457  }
1459  Warning("InitDataSetManager", "dataset manager plug-in initialization failed");
1461  }
1462 
1463  // If no valid dataset manager has been created we instantiate the default one
1464  if (!fDataSetManager) {
1465  TString opts("Av:");
1466  TString dsetdir = gEnv->GetValue("ProofServ.DataSetDir", "");
1467  if (dsetdir.IsNull()) {
1468  // Use the default in the sandbox
1469  dsetdir = fDataSetDir;
1470  opts += "Sb:";
1471  }
1472  // Find the appropriate handler
1473  if (!h) {
1474  h = gROOT->GetPluginManager()->FindHandler("TDataSetManager", "file");
1475  if (h && h->LoadPlugin() == -1) h = 0;
1476  }
1477  if (h) {
1478  // make instance of the dataset manager
1479  fDataSetManager = reinterpret_cast<TDataSetManager*>(h->ExecPlugin(3,
1480  group.Data(), user.Data(),
1481  Form("dir:%s opt:%s", dsetdir.Data(), opts.Data())));
1482  }
1484  Warning("InitDataSetManager", "default dataset manager plug-in initialization failed");
1486  }
1487  }
1488 
1489  if (gDebug > 0 && fDataSetManager) {
1490  Info("InitDataSetManager", "datasetmgr Cq: %d, Ar: %d, Av: %d, Ti: %d, Sb: %d",
1496  }
1497 
1498  // Dataset manager for staging requests
1499  TString dsReqCfg = gEnv->GetValue("Proof.DataSetStagingRequests", "");
1500  if (!dsReqCfg.IsNull()) {
1501  TPMERegexp reReqDir("(^| )(dir:)?([^ ]+)( |$)");
1502 
1503  if (reReqDir.Match(dsReqCfg) == 5) {
1504  TString dsDirFmt;
1505  dsDirFmt.Form("dir:%s perms:open", reReqDir[3].Data());
1506  fDataSetStgRepo = new TDataSetManagerFile("_stage_", "_stage_", dsDirFmt);
1508  Warning("InitDataSetManager", "failed init of dataset staging requests repository");
1510  }
1511  } else {
1512  Warning("InitDataSetManager", "specify, with [dir:]<path>, a valid path for staging requests");
1513  }
1514  } else if (gDebug > 0) {
1515  Warning("InitDataSetManager", "no repository for staging requests available");
1516  }
1517 
1518  // Done
1519  return (fDataSetManager ? 0 : -1);
1520 }
1521 
1522 ////////////////////////////////////////////////////////////////////////////////
1523 /// List contents of file cache. If all is true show all caches also on
1524 /// slaves. If everything is ok all caches are to be the same.
1525 
1527 {
1528  if (!IsValid()) return;
1529 
1530  Printf("*** Local file cache %s ***", fCacheDir.Data());
1531  gSystem->Exec(Form("%s %s", kLS, fCacheDir.Data()));
1532 }
1533 
1534 ////////////////////////////////////////////////////////////////////////////////
1535 /// Remove files from all file caches.
1536 
1537 void TProofLite::ClearCache(const char *file)
1538 {
1539  if (!IsValid()) return;
1540 
1541  fCacheLock->Lock();
1542  if (!file || strlen(file) <= 0) {
1543  gSystem->Exec(Form("%s %s/*", kRM, fCacheDir.Data()));
1544  } else {
1545  gSystem->Exec(Form("%s %s/%s", kRM, fCacheDir.Data(), file));
1546  }
1547  fCacheLock->Unlock();
1548 }
1549 
1550 ////////////////////////////////////////////////////////////////////////////////
1551 /// Copy the specified macro in the cache directory. The macro file is
1552 /// uploaded if new or updated. If existing, the corresponding header
1553 /// basename(macro).h or .hh, is also uploaded. For the other arguments
1554 /// see TProof::Load().
1555 /// Returns 0 in case of success and -1 in case of error.
1556 
1557 Int_t TProofLite::Load(const char *macro, Bool_t notOnClient, Bool_t uniqueOnly,
1558  TList *wrks)
1559 {
1560  if (!IsValid()) return -1;
1561 
1562  if (!macro || !macro[0]) {
1563  Error("Load", "need to specify a macro name");
1564  return -1;
1565  }
1566 
1567  TString macs(macro), mac;
1568  Int_t from = 0;
1569  while (macs.Tokenize(mac, from, ",")) {
1570  if (IsIdle()) {
1571  if (CopyMacroToCache(mac) < 0) return -1;
1572  } else {
1573  // The name
1574  TString macn = gSystem->BaseName(mac);
1575  macn.Remove(macn.Last('.'));
1576  // Relevant pointers
1577  TList cachedFiles;
1578  TString cacheDir = fCacheDir;
1579  gSystem->ExpandPathName(cacheDir);
1580  void * dirp = gSystem->OpenDirectory(cacheDir);
1581  if (dirp) {
1582  const char *e = 0;
1583  while ((e = gSystem->GetDirEntry(dirp))) {
1584  if (!strncmp(e, macn.Data(), macn.Length())) {
1585  TString fncache = Form("%s/%s", cacheDir.Data(), e);
1586  cachedFiles.Add(new TObjString(fncache.Data()));
1587  }
1588  }
1589  gSystem->FreeDirectory(dirp);
1590  }
1591  }
1592  }
1593 
1594  return TProof::Load(macro, notOnClient, uniqueOnly, wrks);
1595 }
1596 
1597 ////////////////////////////////////////////////////////////////////////////////
1598 /// Copy a macro, and its possible associated .h[h] file,
1599 /// to the cache directory, from where the workers can get the file.
1600 /// If headerRequired is 1, return -1 in case the header is not found.
1601 /// If headerRequired is 0, try to copy header too.
1602 /// If headerRequired is -1, don't look for header, only copy macro.
1603 /// If the selector pionter is not 0, consider the macro to be a selector
1604 /// and try to load the selector and set it to the pointer.
1605 /// The mask 'opt' is an or of ESendFileOpt:
1606 /// kCpBin (0x8) Retrieve from the cache the binaries associated
1607 /// with the file
1608 /// kCp (0x10) Retrieve the files from the cache
1609 /// Return -1 in case of error, 0 otherwise.
1610 
1611 Int_t TProofLite::CopyMacroToCache(const char *macro, Int_t headerRequired,
1612  TSelector **selector, Int_t opt, TList *)
1613 {
1614  // Relevant pointers
1615  TString cacheDir = fCacheDir;
1616  gSystem->ExpandPathName(cacheDir);
1617  TProofLockPath *cacheLock = fCacheLock;
1618 
1619  // Split out the aclic mode, if any
1620  TString name = macro;
1621  TString acmode, args, io;
1622  name = gSystem->SplitAclicMode(name, acmode, args, io);
1623 
1624  PDB(kGlobal,1)
1625  Info("CopyMacroToCache", "enter: names: %s, %s", macro, name.Data());
1626 
1627  // Make sure that the file exists
1628  if (gSystem->AccessPathName(name, kReadPermission)) {
1629  Error("CopyMacroToCache", "file %s not found or not readable", name.Data());
1630  return -1;
1631  }
1632 
1633  // Update the macro path
1635  TString np(gSystem->DirName(name));
1636  if (!np.IsNull()) {
1637  np += ":";
1638  if (!mp.BeginsWith(np) && !mp.Contains(":"+np)) {
1639  Int_t ip = (mp.BeginsWith(".:")) ? 2 : 0;
1640  mp.Insert(ip, np);
1641  TROOT::SetMacroPath(mp);
1642  PDB(kGlobal,1)
1643  Info("CopyMacroToCache", "macro path set to '%s'", TROOT::GetMacroPath());
1644  }
1645  }
1646 
1647  // Check the header file
1648  Int_t dot = name.Last('.');
1649  const char *hext[] = { ".h", ".hh", "" };
1650  TString hname, checkedext;
1651  Int_t i = 0;
1652  while (strlen(hext[i]) > 0) {
1653  hname = name(0, dot);
1654  hname += hext[i];
1655  if (!gSystem->AccessPathName(hname, kReadPermission))
1656  break;
1657  if (!checkedext.IsNull()) checkedext += ",";
1658  checkedext += hext[i];
1659  hname = "";
1660  i++;
1661  }
1662  if (hname.IsNull() && headerRequired == 1) {
1663  Error("CopyMacroToCache", "header file for %s not found or not readable "
1664  "(checked extensions: %s)", name.Data(), checkedext.Data());
1665  return -1;
1666  }
1667  if (headerRequired < 0)
1668  hname = "";
1669 
1670  cacheLock->Lock();
1671 
1672  // Check these files with those in the cache (if any)
1673  Bool_t useCacheBinaries = kFALSE;
1674  TString cachedname = Form("%s/%s", cacheDir.Data(), gSystem->BaseName(name));
1675  TString cachedhname;
1676  if (!hname.IsNull())
1677  cachedhname = Form("%s/%s", cacheDir.Data(), gSystem->BaseName(hname));
1678  if (!gSystem->AccessPathName(cachedname, kReadPermission)) {
1679  TMD5 *md5 = TMD5::FileChecksum(name);
1680  TMD5 *md5cache = TMD5::FileChecksum(cachedname);
1681  if (md5 && md5cache && (*md5 == *md5cache))
1682  useCacheBinaries = kTRUE;
1683  if (!hname.IsNull()) {
1684  if (!gSystem->AccessPathName(cachedhname, kReadPermission)) {
1685  TMD5 *md5h = TMD5::FileChecksum(hname);
1686  TMD5 *md5hcache = TMD5::FileChecksum(cachedhname);
1687  if (md5h && md5hcache && (*md5h != *md5hcache))
1688  useCacheBinaries = kFALSE;
1689  SafeDelete(md5h);
1690  SafeDelete(md5hcache);
1691  }
1692  }
1693  SafeDelete(md5);
1694  SafeDelete(md5cache);
1695  }
1696 
1697  // Create version file name template
1698  TString vername(Form(".%s", name.Data()));
1699  dot = vername.Last('.');
1700  if (dot != kNPOS)
1701  vername.Remove(dot);
1702  vername += ".binversion";
1703  Bool_t savever = kFALSE;
1704 
1705  // Check binary version
1706  if (useCacheBinaries) {
1707  TString v, r;
1708  FILE *f = fopen(Form("%s/%s", cacheDir.Data(), vername.Data()), "r");
1709  if (f) {
1710  v.Gets(f);
1711  r.Gets(f);
1712  fclose(f);
1713  }
1714  if (!f || v != gROOT->GetVersion() || r != gROOT->GetGitCommit())
1715  useCacheBinaries = kFALSE;
1716  }
1717 
1718  // Create binary name template
1719  TString binname = gSystem->BaseName(name);
1720  dot = binname.Last('.');
1721  if (dot != kNPOS)
1722  binname.Replace(dot,1,"_");
1723  TString pcmname = TString::Format("%s_ACLiC_dict_rdict.pcm", binname.Data());
1724  binname += ".";
1725 
1726  FileStat_t stlocal, stcache;
1727  void *dirp = 0;
1728  if (useCacheBinaries) {
1729  // Loop over binaries in the cache and copy them locally if newer then the local
1730  // versions or there is no local version
1731  dirp = gSystem->OpenDirectory(cacheDir);
1732  if (dirp) {
1733  const char *e = 0;
1734  while ((e = gSystem->GetDirEntry(dirp))) {
1735  if (!strncmp(e, binname.Data(), binname.Length()) ||
1736  !strncmp(e, pcmname.Data(), pcmname.Length())) {
1737  TString fncache = Form("%s/%s", cacheDir.Data(), e);
1738  Bool_t docp = kTRUE;
1739  if (!gSystem->GetPathInfo(fncache, stcache)) {
1740  Int_t rc = gSystem->GetPathInfo(e, stlocal);
1741  if (rc == 0 && (stlocal.fMtime >= stcache.fMtime))
1742  docp = kFALSE;
1743  // Copy the file, if needed
1744  if (docp) {
1745  gSystem->Exec(Form("%s %s", kRM, e));
1746  PDB(kGlobal,2)
1747  Info("CopyMacroToCache",
1748  "retrieving %s from cache", fncache.Data());
1749  gSystem->Exec(Form("%s %s %s", kCP, fncache.Data(), e));
1750  }
1751  }
1752  }
1753  }
1754  gSystem->FreeDirectory(dirp);
1755  }
1756  }
1757  cacheLock->Unlock();
1758 
1759  if (selector) {
1760  // Now init the selector in optimized way
1761  if (!(*selector = TSelector::GetSelector(macro))) {
1762  Error("CopyMacroToCache", "could not create a selector from %s", macro);
1763  return -1;
1764  }
1765  }
1766 
1767  cacheLock->Lock();
1768 
1769  TList *cachedFiles = new TList;
1770  // Save information in the cache now for later usage
1771  dirp = gSystem->OpenDirectory(".");
1772  if (dirp) {
1773  const char *e = 0;
1774  while ((e = gSystem->GetDirEntry(dirp))) {
1775  if (!strncmp(e, binname.Data(), binname.Length()) ||
1776  !strncmp(e, pcmname.Data(), pcmname.Length())) {
1777  Bool_t docp = kTRUE;
1778  if (!gSystem->GetPathInfo(e, stlocal)) {
1779  TString fncache = Form("%s/%s", cacheDir.Data(), e);
1780  Int_t rc = gSystem->GetPathInfo(fncache, stcache);
1781  if (rc == 0 && (stlocal.fMtime <= stcache.fMtime))
1782  docp = kFALSE;
1783  // Copy the file, if needed
1784  if (docp) {
1785  gSystem->Exec(Form("%s %s", kRM, fncache.Data()));
1786  PDB(kGlobal,2)
1787  Info("CopyMacroToCache","caching %s ...", e);
1788  gSystem->Exec(Form("%s %s %s", kCP, e, fncache.Data()));
1789  savever = kTRUE;
1790  }
1791  if (opt & kCpBin)
1792  cachedFiles->Add(new TObjString(fncache.Data()));
1793  }
1794  }
1795  }
1796  gSystem->FreeDirectory(dirp);
1797  }
1798 
1799  // Save binary version if requested
1800  if (savever) {
1801  FILE *f = fopen(Form("%s/%s", cacheDir.Data(), vername.Data()), "w");
1802  if (f) {
1803  fputs(gROOT->GetVersion(), f);
1804  fputs(Form("\n%s", gROOT->GetGitCommit()), f);
1805  fclose(f);
1806  }
1807  }
1808 
1809  // Save also the selector info, if needed
1810  if (!useCacheBinaries) {
1811  gSystem->Exec(Form("%s %s", kRM, cachedname.Data()));
1812  PDB(kGlobal,2)
1813  Info("CopyMacroToCache","caching %s ...", name.Data());
1814  gSystem->Exec(Form("%s %s %s", kCP, name.Data(), cachedname.Data()));
1815  if (!hname.IsNull()) {
1816  gSystem->Exec(Form("%s %s", kRM, cachedhname.Data()));
1817  PDB(kGlobal,2)
1818  Info("CopyMacroToCache","caching %s ...", hname.Data());
1819  gSystem->Exec(Form("%s %s %s", kCP, hname.Data(), cachedhname.Data()));
1820  }
1821  }
1822  if (opt & kCp) {
1823  cachedFiles->Add(new TObjString(cachedname.Data()));
1824  if (!hname.IsNull())
1825  cachedFiles->Add(new TObjString(cachedhname.Data()));
1826  }
1827 
1828  cacheLock->Unlock();
1829 
1830  cachedFiles->SetOwner();
1831  delete cachedFiles;
1832 
1833  return 0;
1834 }
1835 
1836 ////////////////////////////////////////////////////////////////////////////////
1837 /// Remove old sessions dirs keep at most 'Proof.MaxOldSessions' (default 10)
1838 
1840 {
1841  Int_t maxold = gEnv->GetValue("Proof.MaxOldSessions", 1);
1842 
1843  if (maxold < 0) return 0;
1844 
1845  TSortedList *olddirs = new TSortedList(kFALSE);
1846 
1847  TString sandbox = gSystem->DirName(fWorkDir.Data());
1848 
1849  void *dirp = gSystem->OpenDirectory(sandbox);
1850  if (dirp) {
1851  const char *e = 0;
1852  while ((e = gSystem->GetDirEntry(dirp))) {
1853  if (!strncmp(e, "session-", 8) && !strstr(e, GetName())) {
1854  TString d(e);
1855  Int_t i = d.Last('-');
1856  if (i != kNPOS) d.Remove(i);
1857  i = d.Last('-');
1858  if (i != kNPOS) d.Remove(0,i+1);
1859  TString path = Form("%s/%s", sandbox.Data(), e);
1860  olddirs->Add(new TNamed(d, path));
1861  }
1862  }
1863  gSystem->FreeDirectory(dirp);
1864  }
1865 
1866  // Clean it up, if required
1867  Bool_t notify = kTRUE;
1868  while (olddirs->GetSize() > maxold) {
1869  if (notify && gDebug > 0)
1870  Printf("Cleaning sandbox at: %s", sandbox.Data());
1871  notify = kFALSE;
1872  TNamed *n = (TNamed *) olddirs->Last();
1873  if (n) {
1874  gSystem->Exec(Form("%s %s", kRM, n->GetTitle()));
1875  olddirs->Remove(n);
1876  delete n;
1877  }
1878  }
1879 
1880  // Cleanup
1881  olddirs->SetOwner();
1882  delete olddirs;
1883 
1884  // Done
1885  return 0;
1886 }
1887 
1888 ////////////////////////////////////////////////////////////////////////////////
1889 /// Get the list of queries.
1890 
1892 {
1893  Bool_t all = ((strchr(opt,'A') || strchr(opt,'a'))) ? kTRUE : kFALSE;
1894 
1895  TList *ql = new TList;
1896  Int_t ntot = 0, npre = 0, ndraw= 0;
1897  if (fQMgr) {
1898  if (all) {
1899  // Rescan
1900  TString qdir = fQueryDir;
1901  Int_t idx = qdir.Index("session-");
1902  if (idx != kNPOS)
1903  qdir.Remove(idx);
1904  fQMgr->ScanPreviousQueries(qdir);
1905  // Gather also information about previous queries, if any
1906  if (fQMgr->PreviousQueries()) {
1907  TIter nxq(fQMgr->PreviousQueries());
1908  TProofQueryResult *pqr = 0;
1909  while ((pqr = (TProofQueryResult *)nxq())) {
1910  ntot++;
1911  pqr->fSeqNum = ntot;
1912  ql->Add(pqr);
1913  }
1914  }
1915  }
1916 
1917  npre = ntot;
1918  if (fQMgr->Queries()) {
1919  // Add info about queries in this session
1920  TIter nxq(fQMgr->Queries());
1921  TProofQueryResult *pqr = 0;
1922  TQueryResult *pqm = 0;
1923  while ((pqr = (TProofQueryResult *)nxq())) {
1924  ntot++;
1925  if ((pqm = pqr->CloneInfo())) {
1926  pqm->fSeqNum = ntot;
1927  ql->Add(pqm);
1928  } else {
1929  Warning("GetListOfQueries", "unable to clone TProofQueryResult '%s:%s'",
1930  pqr->GetName(), pqr->GetTitle());
1931  }
1932  }
1933  }
1934  // Number of draw queries
1935  ndraw = fQMgr->DrawQueries();
1936  }
1937 
1938  fOtherQueries = npre;
1939  fDrawQueries = ndraw;
1940  if (fQueries) {
1941  fQueries->Delete();
1942  delete fQueries;
1943  fQueries = 0;
1944  }
1945  fQueries = ql;
1946 
1947  // This should have been filled by now
1948  return fQueries;
1949 }
1950 
1951 ////////////////////////////////////////////////////////////////////////////////
1952 /// Register the 'dataSet' on the cluster under the current
1953 /// user, group and the given 'dataSetName'.
1954 /// Fails if a dataset named 'dataSetName' already exists, unless 'optStr'
1955 /// contains 'O', in which case the old dataset is overwritten.
1956 /// If 'optStr' contains 'V' the dataset files are verified (default no
1957 /// verification).
1958 /// Returns kTRUE on success.
1959 
1961  TFileCollection *dataSet, const char* optStr)
1962 {
1963  if (!fDataSetManager) {
1964  Info("RegisterDataSet", "dataset manager not available");
1965  return kFALSE;
1966  }
1967 
1968  if (!uri || strlen(uri) <= 0) {
1969  Info("RegisterDataSet", "specifying a dataset name is mandatory");
1970  return kFALSE;
1971  }
1972 
1973  Bool_t parallelverify = kFALSE;
1974  TString sopt(optStr);
1975  if (sopt.Contains("V") && !sopt.Contains("S")) {
1976  // We do verification in parallel later on; just register for now
1977  parallelverify = kTRUE;
1978  sopt.ReplaceAll("V", "");
1979  }
1980  // This would screw up things remotely, make sure is not there
1981  sopt.ReplaceAll("S", "");
1982 
1983  Bool_t result = kTRUE;
1985  // Check the list
1986  if (!dataSet || dataSet->GetList()->GetSize() == 0) {
1987  Error("RegisterDataSet", "can not save an empty list.");
1988  result = kFALSE;
1989  }
1990  // Register the dataset (quota checks are done inside here)
1991  result = (fDataSetManager->RegisterDataSet(uri, dataSet, sopt) == 0)
1992  ? kTRUE : kFALSE;
1993  } else {
1994  Info("RegisterDataSet", "dataset registration not allowed");
1995  result = kFALSE;
1996  }
1997 
1998  if (!result)
1999  Error("RegisterDataSet", "dataset was not saved");
2000 
2001  // If old server or not verifying in parallel we are done
2002  if (!parallelverify) return result;
2003 
2004  // If we are here it means that we will verify in parallel
2005  sopt += "V";
2006  if (VerifyDataSet(uri, sopt) < 0){
2007  Error("RegisterDataSet", "problems verifying dataset '%s'", uri);
2008  return kFALSE;
2009  }
2010 
2011  // Done
2012  return kTRUE;
2013 }
2014 
2015 ////////////////////////////////////////////////////////////////////////////////
2016 /// Set/Change the name of the default tree. The tree name may contain
2017 /// subdir specification in the form "subdir/name".
2018 /// Returns 0 on success, -1 otherwise.
2019 
2020 Int_t TProofLite::SetDataSetTreeName(const char *dataset, const char *treename)
2021 {
2022  if (!fDataSetManager) {
2023  Info("ExistsDataSet", "dataset manager not available");
2024  return kFALSE;
2025  }
2026 
2027  if (!dataset || strlen(dataset) <= 0) {
2028  Info("SetDataSetTreeName", "specifying a dataset name is mandatory");
2029  return -1;
2030  }
2031 
2032  if (!treename || strlen(treename) <= 0) {
2033  Info("SetDataSetTreeName", "specifying a tree name is mandatory");
2034  return -1;
2035  }
2036 
2037  TUri uri(dataset);
2038  TString fragment(treename);
2039  if (!fragment.BeginsWith("/")) fragment.Insert(0, "/");
2040  uri.SetFragment(fragment);
2041 
2042  return fDataSetManager->ScanDataSet(uri.GetUri().Data(),
2044 }
2045 
2046 ////////////////////////////////////////////////////////////////////////////////
2047 /// Returns kTRUE if 'dataset' described by 'uri' exists, kFALSE otherwise
2048 
2050 {
2051  if (!fDataSetManager) {
2052  Info("ExistsDataSet", "dataset manager not available");
2053  return kFALSE;
2054  }
2055 
2056  if (!uri || strlen(uri) <= 0) {
2057  Error("ExistsDataSet", "dataset name missing");
2058  return kFALSE;
2059  }
2060 
2061  // Check if the dataset exists
2062  return fDataSetManager->ExistsDataSet(uri);
2063 }
2064 
2065 ////////////////////////////////////////////////////////////////////////////////
2066 /// lists all datasets that match given uri
2067 
2068 TMap *TProofLite::GetDataSets(const char *uri, const char *srvex)
2069 {
2070  if (!fDataSetManager) {
2071  Info("GetDataSets", "dataset manager not available");
2072  return (TMap *)0;
2073  }
2074 
2075  // Get the datasets and return the map
2076  if (srvex && strlen(srvex) > 0) {
2077  return fDataSetManager->GetSubDataSets(uri, srvex);
2078  } else {
2080  return fDataSetManager->GetDataSets(uri, opt);
2081  }
2082 }
2083 
2084 ////////////////////////////////////////////////////////////////////////////////
2085 /// Shows datasets in locations that match the uri
2086 /// By default shows the user's datasets and global ones
2087 
2088 void TProofLite::ShowDataSets(const char *uri, const char *opt)
2089 {
2090  if (!fDataSetManager) {
2091  Info("GetDataSet", "dataset manager not available");
2092  return;
2093  }
2094 
2095  fDataSetManager->ShowDataSets(uri, opt);
2096 }
2097 
2098 ////////////////////////////////////////////////////////////////////////////////
2099 /// Get a list of TFileInfo objects describing the files of the specified
2100 /// dataset.
2101 
2102 TFileCollection *TProofLite::GetDataSet(const char *uri, const char *)
2103 {
2104  if (!fDataSetManager) {
2105  Info("GetDataSet", "dataset manager not available");
2106  return (TFileCollection *)0;
2107  }
2108 
2109  if (!uri || strlen(uri) <= 0) {
2110  Info("GetDataSet", "specifying a dataset name is mandatory");
2111  return 0;
2112  }
2113 
2114  // Return the list
2115  return fDataSetManager->GetDataSet(uri);
2116 }
2117 
2118 ////////////////////////////////////////////////////////////////////////////////
2119 /// Remove the specified dataset from the PROOF cluster.
2120 /// Files are not deleted.
2121 
2122 Int_t TProofLite::RemoveDataSet(const char *uri, const char *)
2123 {
2124  if (!fDataSetManager) {
2125  Info("RemoveDataSet", "dataset manager not available");
2126  return -1;
2127  }
2128 
2130  if (!fDataSetManager->RemoveDataSet(uri)) {
2131  // Failure
2132  return -1;
2133  }
2134  } else {
2135  Info("RemoveDataSet", "dataset creation / removal not allowed");
2136  return -1;
2137  }
2138 
2139  // Done
2140  return 0;
2141 }
2142 
2143 ////////////////////////////////////////////////////////////////////////////////
2144 /// Allows users to request staging of a particular dataset. Requests are
2145 /// saved in a special dataset repository and must be honored by the endpoint.
2146 /// This is the special PROOF-Lite re-implementation of the TProof function
2147 /// and includes code originally implemented in TProofServ.
2148 
2150 {
2151  if (!dataset) {
2152  Error("RequestStagingDataSet", "invalid dataset specified");
2153  return kFALSE;
2154  }
2155 
2156  if (!fDataSetStgRepo) {
2157  Error("RequestStagingDataSet", "no dataset staging request repository available");
2158  return kFALSE;
2159  }
2160 
2161  TString dsUser, dsGroup, dsName, dsTree;
2162 
2163  // Transform input URI in a valid dataset name
2164  TString validUri = dataset;
2165  while (fReInvalid->Substitute(validUri, "_")) {}
2166 
2167  // Check if dataset exists beforehand: if it does, staging has already been requested
2168  if (fDataSetStgRepo->ExistsDataSet(validUri.Data())) {
2169  Warning("RequestStagingDataSet", "staging of %s already requested", dataset);
2170  return kFALSE;
2171  }
2172 
2173  // Try to get dataset from current manager
2175  if (!fc || (fc->GetNFiles() == 0)) {
2176  Error("RequestStagingDataSet", "empty dataset or no dataset returned");
2177  if (fc) delete fc;
2178  return kFALSE;
2179  }
2180 
2181  // Reset all staged bits and remove unnecessary URLs (all but last)
2182  TIter it(fc->GetList());
2183  TFileInfo *fi;
2184  while ((fi = dynamic_cast<TFileInfo *>(it.Next()))) {
2186  Int_t nToErase = fi->GetNUrls() - 1;
2187  for (Int_t i=0; i<nToErase; i++)
2188  fi->RemoveUrlAt(0);
2189  }
2190 
2191  fc->Update(); // absolutely necessary
2192 
2193  // Save request
2194  fDataSetStgRepo->ParseUri(validUri, &dsGroup, &dsUser, &dsName);
2195  if (fDataSetStgRepo->WriteDataSet(dsGroup, dsUser, dsName, fc) == 0) {
2196  // Error, can't save dataset
2197  Error("RequestStagingDataSet", "can't register staging request for %s", dataset);
2198  delete fc;
2199  return kFALSE;
2200  }
2201 
2202  Info("RequestStagingDataSet", "Staging request registered for %s", dataset);
2203  delete fc;
2204 
2205  return kTRUE;
2206 }
2207 
2208 ////////////////////////////////////////////////////////////////////////////////
2209 /// Cancels a dataset staging request. Returns kTRUE on success, kFALSE on
2210 /// failure. Dataset not found equals to a failure. PROOF-Lite
2211 /// re-implementation of the equivalent function in TProofServ.
2212 
2214 {
2215  if (!dataset) {
2216  Error("CancelStagingDataSet", "invalid dataset specified");
2217  return kFALSE;
2218  }
2219 
2220  if (!fDataSetStgRepo) {
2221  Error("CancelStagingDataSet", "no dataset staging request repository available");
2222  return kFALSE;
2223  }
2224 
2225  // Transform URI in a valid dataset name
2226  TString validUri = dataset;
2227  while (fReInvalid->Substitute(validUri, "_")) {}
2228 
2229  if (!fDataSetStgRepo->RemoveDataSet(validUri.Data()))
2230  return kFALSE;
2231 
2232  return kTRUE;
2233 }
2234 
2235 ////////////////////////////////////////////////////////////////////////////////
2236 /// Obtains a TFileCollection showing the staging status of the specified
2237 /// dataset. A valid dataset manager and dataset staging requests repository
2238 /// must be present on the endpoint. PROOF-Lite version of the equivalent
2239 /// function from TProofServ.
2240 
2242 {
2243  if (!dataset) {
2244  Error("GetStagingStatusDataSet", "invalid dataset specified");
2245  return 0;
2246  }
2247 
2248  if (!fDataSetStgRepo) {
2249  Error("GetStagingStatusDataSet", "no dataset staging request repository available");
2250  return 0;
2251  }
2252 
2253  // Transform URI in a valid dataset name
2254  TString validUri = dataset;
2255  while (fReInvalid->Substitute(validUri, "_")) {}
2256 
2257  // Get the list
2259  if (!fc) {
2260  // No such dataset (not an error)
2261  Info("GetStagingStatusDataSet", "no pending staging request for %s", dataset);
2262  return 0;
2263  }
2264 
2265  // Dataset found: return it (must be cleaned by caller)
2266  return fc;
2267 }
2268 
2269 ////////////////////////////////////////////////////////////////////////////////
2270 /// Verify if all files in the specified dataset are available.
2271 /// Print a list and return the number of missing files.
2272 
2273 Int_t TProofLite::VerifyDataSet(const char *uri, const char *optStr)
2274 {
2275  if (!fDataSetManager) {
2276  Info("VerifyDataSet", "dataset manager not available");
2277  return -1;
2278  }
2279 
2280  Int_t rc = -1;
2281  TString sopt(optStr);
2282  if (sopt.Contains("S")) {
2283 
2285  rc = fDataSetManager->ScanDataSet(uri);
2286  } else {
2287  Info("VerifyDataSet", "dataset verification not allowed");
2288  rc = -1;
2289  }
2290  return rc;
2291  }
2292 
2293  // Done
2294  return VerifyDataSetParallel(uri, optStr);
2295 }
2296 
2297 ////////////////////////////////////////////////////////////////////////////////
2298 /// Clear the content of the dataset cache, if any (matching 'dataset', if defined).
2299 
2300 void TProofLite::ClearDataSetCache(const char *dataset)
2301 {
2303  // Done
2304  return;
2305 }
2306 
2307 ////////////////////////////////////////////////////////////////////////////////
2308 /// Display the content of the dataset cache, if any (matching 'dataset', if defined).
2309 
2310 void TProofLite::ShowDataSetCache(const char *dataset)
2311 {
2312  // For PROOF-Lite act locally
2313  if (fDataSetManager) fDataSetManager->ShowCache(dataset);
2314  // Done
2315  return;
2316 }
2317 
2318 ////////////////////////////////////////////////////////////////////////////////
2319 /// Make sure that the input data objects are available to the workers in a
2320 /// dedicated file in the cache; the objects are taken from the dedicated list
2321 /// and / or the specified file.
2322 /// If the fInputData is empty the specified file is sent over.
2323 /// If there is no specified file, a file named "inputdata.root" is created locally
2324 /// with the content of fInputData and sent over to the master.
2325 /// If both fInputData and the specified file are not empty, a copy of the file
2326 /// is made locally and augmented with the content of fInputData.
2327 
2329 {
2330  // Prepare the file
2331  TString dataFile;
2332  PrepareInputDataFile(dataFile);
2333 
2334  // Make sure it is in the cache, if not empty
2335  if (dataFile.Length() > 0) {
2336 
2337  if (!dataFile.BeginsWith(fCacheDir)) {
2338  // Destination
2339  TString dst;
2340  dst.Form("%s/%s", fCacheDir.Data(), gSystem->BaseName(dataFile));
2341  // Remove it first if it exists
2342  if (!gSystem->AccessPathName(dst))
2343  gSystem->Unlink(dst);
2344  // Copy the file
2345  if (gSystem->CopyFile(dataFile, dst) != 0)
2346  Warning("SendInputDataFile", "problems copying '%s' to '%s'",
2347  dataFile.Data(), dst.Data());
2348  }
2349 
2350  // Set the name in the input list so that the workers can find it
2351  AddInput(new TNamed("PROOF_InputDataFile", Form("%s", gSystem->BaseName(dataFile))));
2352  }
2353 }
2354 
2355 ////////////////////////////////////////////////////////////////////////////////
2356 /// Handle remove request.
2357 
2358 Int_t TProofLite::Remove(const char *ref, Bool_t all)
2359 {
2360  PDB(kGlobal, 1)
2361  Info("Remove", "Enter: %s, %d", ref, all);
2362 
2363  if (all) {
2364  // Remove also local copies, if any
2365  if (fPlayer)
2366  fPlayer->RemoveQueryResult(ref);
2367  }
2368 
2369  TString queryref(ref);
2370 
2371  if (queryref == "cleanupdir") {
2372 
2373  // Cleanup previous sessions results
2374  Int_t nd = (fQMgr) ? fQMgr->CleanupQueriesDir() : -1;
2375 
2376  // Notify
2377  Info("Remove", "%d directories removed", nd);
2378  // We are done
2379  return 0;
2380  }
2381 
2382 
2383  if (fQMgr) {
2384  TProofLockPath *lck = 0;
2385  if (fQMgr->LockSession(queryref, &lck) == 0) {
2386 
2387  // Remove query
2388  fQMgr->RemoveQuery(queryref, 0);
2389 
2390  // Unlock and remove the lock file
2391  if (lck) {
2392  gSystem->Unlink(lck->GetName());
2393  SafeDelete(lck);
2394  }
2395 
2396  // We are done
2397  return 0;
2398  }
2399  } else {
2400  Warning("Remove", "query result manager undefined!");
2401  }
2402 
2403  // Notify failure
2404  Info("Remove",
2405  "query %s could not be removed (unable to lock session)", queryref.Data());
2406 
2407  // Done
2408  return -1;
2409 }
2410 
2411 ////////////////////////////////////////////////////////////////////////////////
2412 /// Creates a tree header (a tree with nonexisting files) object for
2413 /// the DataSet.
2414 
2416 {
2417  TTree *t = 0;
2418  if (!dset) {
2419  Error("GetTreeHeader", "undefined TDSet");
2420  return t;
2421  }
2422 
2423  dset->Reset();
2424  TDSetElement *e = dset->Next();
2425  Long64_t entries = 0;
2426  TFile *f = 0;
2427  if (!e) {
2428  PDB(kGlobal, 1) Info("GetTreeHeader", "empty TDSet");
2429  } else {
2430  f = TFile::Open(e->GetFileName());
2431  t = 0;
2432  if (f) {
2433  t = (TTree*) f->Get(e->GetObjName());
2434  if (t) {
2435  t->SetMaxVirtualSize(0);
2436  t->DropBaskets();
2437  entries = t->GetEntries();
2438 
2439  // compute #entries in all the files
2440  while ((e = dset->Next()) != 0) {
2441  TFile *f1 = TFile::Open(e->GetFileName());
2442  if (f1) {
2443  TTree *t1 = (TTree*) f1->Get(e->GetObjName());
2444  if (t1) {
2445  entries += t1->GetEntries();
2446  delete t1;
2447  }
2448  delete f1;
2449  }
2450  }
2451  t->SetMaxEntryLoop(entries); // this field will hold the total number of entries ;)
2452  }
2453  }
2454  }
2455  // Done
2456  return t;
2457 }
2458 
2459 ////////////////////////////////////////////////////////////////////////////////
2460 /// Add to the fUniqueSlave list the active slaves that have a unique
2461 /// (user) file system image. This information is used to transfer files
2462 /// only once to nodes that share a file system (an image). Submasters
2463 /// which are not in fUniqueSlaves are put in the fNonUniqueMasters
2464 /// list. That list is used to trigger the transferring of files to
2465 /// the submaster's unique slaves without the need to transfer the file
2466 /// to the submaster.
2467 
2469 {
2470  fUniqueSlaves->Clear();
2475 
2476  if (fActiveSlaves->GetSize() <= 0) return;
2477 
2478  TSlave *wrk = dynamic_cast<TSlave*>(fActiveSlaves->First());
2479  if (!wrk) {
2480  Error("FindUniqueSlaves", "first object in fActiveSlaves not a TSlave: embarrasing!");
2481  return;
2482  }
2483  fUniqueSlaves->Add(wrk);
2484  fAllUniqueSlaves->Add(wrk);
2485  fUniqueMonitor->Add(wrk->GetSocket());
2486  fAllUniqueMonitor->Add(wrk->GetSocket());
2487 
2488  // will be actiavted in Collect()
2491 }
2492 
2493 ////////////////////////////////////////////////////////////////////////////////
2494 /// List contents of the data directory in the sandbox.
2495 /// This is the place where files produced by the client queries are kept
2496 
2498 {
2499  if (!IsValid()) return;
2500 
2501  // Get worker infos
2502  TList *wrki = GetListOfSlaveInfos();
2503  TSlaveInfo *wi = 0;
2504  TIter nxwi(wrki);
2505  while ((wi = (TSlaveInfo *) nxwi())) {
2506  ShowDataDir(wi->GetDataDir());
2507  }
2508 }
2509 
2510 ////////////////////////////////////////////////////////////////////////////////
2511 /// List contents of the data directory 'dirname'
2512 
2513 void TProofLite::ShowDataDir(const char *dirname)
2514 {
2515  if (!dirname) return;
2516 
2517  FileStat_t dirst;
2518  if (gSystem->GetPathInfo(dirname, dirst) != 0) return;
2519  if (!R_ISDIR(dirst.fMode)) return;
2520 
2521  void *dirp = gSystem->OpenDirectory(dirname);
2522  TString fn;
2523  const char *ent = 0;
2524  while ((ent = gSystem->GetDirEntry(dirp))) {
2525  fn.Form("%s/%s", dirname, ent);
2526  FileStat_t st;
2527  if (gSystem->GetPathInfo(fn.Data(), st) == 0) {
2528  if (R_ISREG(st.fMode)) {
2529  Printf("lite:0| %s", fn.Data());
2530  } else if (R_ISREG(st.fMode)) {
2531  ShowDataDir(fn.Data());
2532  }
2533  }
2534  }
2535  // Done
2536  return;
2537 }
2538 
2539 ////////////////////////////////////////////////////////////////////////////////
2540 /// Simulate dynamic addition, for test purposes.
2541 /// Here we decide how many workers to add, we create them and set the
2542 /// environment.
2543 /// This call is called regularly by Collect if the opton is enabled.
2544 /// Returns the number of new workers added, or <0 on errors.
2545 
2547 {
2548  // Max workers
2549  if (fDynamicStartupNMax <= 0) {
2550  SysInfo_t si;
2551  if (gSystem->GetSysInfo(&si) == 0 && si.fCpus > 2) {
2553  } else {
2554  fDynamicStartupNMax = 2;
2555  }
2556  }
2557  if (fNWorkers >= fDynamicStartupNMax) {
2558  // Max reached: disable
2559  Info("PollForNewWorkers", "max reached: %d workers started", fNWorkers);
2561  return 0;
2562  }
2563 
2564  // Number of new workers
2565  Int_t nAdd = (fDynamicStartupStep > 0) ? fDynamicStartupStep : 1;
2566 
2567  // Create a monitor and add the socket to it
2568  TMonitor *mon = new TMonitor;
2569  mon->Add(fServSock);
2570 
2571  TList started;
2572  TSlave *wrk = 0;
2573  Int_t nWrksDone = 0, nWrksTot = -1;
2574  TString fullord;
2575 
2576  nWrksTot = fNWorkers + nAdd;
2577  // Now we create the worker applications which will call us back to finalize
2578  // the setup
2579  Int_t ord = fNWorkers;
2580  for (; ord < nWrksTot; ord++) {
2581 
2582  // Ordinal for this worker server
2583  fullord = Form("0.%d", ord);
2584 
2585  // Create environment files
2586  SetProofServEnv(fullord);
2587 
2588  // Create worker server and add to the list
2589  if ((wrk = CreateSlave("lite", fullord, 100, fImage, fWorkDir)))
2590  started.Add(wrk);
2591 
2592  PDB(kGlobal, 3)
2593  Info("PollForNewWorkers", "additional worker '%s' started", fullord.Data());
2594 
2595  // Notify
2596  NotifyStartUp("Opening connections to workers", ++nWrksDone, nWrksTot);
2597 
2598  } //end of worker loop
2599  fNWorkers = nWrksTot;
2600 
2601  // A list of TSlave objects for workers that are being added
2602  TList *addedWorkers = new TList();
2603  addedWorkers->SetOwner(kFALSE);
2604 
2605  // Wait for call backs
2606  nWrksDone = 0;
2607  nWrksTot = started.GetSize();
2608  Int_t nSelects = 0;
2609  Int_t to = gEnv->GetValue("ProofLite.StartupTimeOut", 5) * 1000;
2610  while (started.GetSize() > 0 && nSelects < nWrksTot) {
2611 
2612  // Wait for activity on the socket for max 5 secs
2613  TSocket *xs = mon->Select(to);
2614 
2615  // Count attempts and check
2616  nSelects++;
2617  if (xs == (TSocket *) -1) continue;
2618 
2619  // Get the connection
2620  TSocket *s = fServSock->Accept();
2621  if (s && s->IsValid()) {
2622  // Receive ordinal
2623  TMessage *msg = 0;
2624  if (s->Recv(msg) < 0) {
2625  Warning("PollForNewWorkers", "problems receiving message from accepted socket!");
2626  } else {
2627  if (msg) {
2628  *msg >> fullord;
2629  // Find who is calling back
2630  if ((wrk = (TSlave *) started.FindObject(fullord))) {
2631  // Remove it from the started list
2632  started.Remove(wrk);
2633 
2634  // Assign tis socket the selected worker
2635  wrk->SetSocket(s);
2636  // Remove socket from global TROOT socket list. Only the TProof object,
2637  // representing all worker sockets, will be added to this list. This will
2638  // ensure the correct termination of all proof servers in case the
2639  // root session terminates.
2641  gROOT->GetListOfSockets()->Remove(s);
2642  }
2643  if (wrk->IsValid()) {
2644  // Set the input handler
2645  wrk->SetInputHandler(new TProofInputHandler(this, wrk->GetSocket()));
2646  // Set fParallel to 1 for workers since they do not
2647  // report their fParallel with a LOG_DONE message
2648  wrk->fParallel = 1;
2649  // Finalize setup of the server
2650  wrk->SetupServ(TSlave::kSlave, 0);
2651  }
2652 
2653  // Monitor good workers
2654  fSlaves->Add(wrk);
2655  if (wrk->IsValid()) {
2656  fActiveSlaves->Add(wrk); // Is this required? Check!
2657  fAllMonitor->Add(wrk->GetSocket());
2658  // Record also in the list for termination
2659  if (addedWorkers) addedWorkers->Add(wrk);
2660  // Notify startup operations
2661  NotifyStartUp("Setting up added worker servers", ++nWrksDone, nWrksTot);
2662  } else {
2663  // Flag as bad
2664  fBadSlaves->Add(wrk);
2665  }
2666  }
2667  } else {
2668  Warning("PollForNewWorkers", "received empty message from accepted socket!");
2669  }
2670  }
2671  }
2672  }
2673 
2674  // Cleanup the monitor and the server socket
2675  mon->DeActivateAll();
2676  delete mon;
2677 
2678  Broadcast(kPROOF_GETSTATS, addedWorkers);
2679  Collect(addedWorkers, fCollectTimeout);
2680 
2681  // Update group view
2682  // SendGroupView();
2683 
2684  // By default go into parallel mode
2685  // SetParallel(-1, 0);
2686  SendCurrentState(addedWorkers);
2687 
2688  // Set worker processing environment
2689  SetupWorkersEnv(addedWorkers, kTRUE);
2690 
2691  // We are adding workers dynamically to an existing process, we
2692  // should invoke a special player's Process() to set only added workers
2693  // to the proper state
2694  if (fPlayer) {
2695  PDB(kGlobal, 3)
2696  Info("PollForNewWorkers", "Will send the PROCESS message to selected workers");
2697  fPlayer->JoinProcess(addedWorkers);
2698  }
2699 
2700  // Cleanup fwhat remained from startup
2701  Collect(addedWorkers);
2702 
2703  // Activate
2704  TIter naw(addedWorkers);
2705  while ((wrk = (TSlave *)naw())) {
2706  fActiveMonitor->Add(wrk->GetSocket());
2707  }
2708  // Cleanup
2709  delete addedWorkers;
2710 
2711  // Done
2712  return nWrksDone;
2713 }
Int_t SetProofServEnv(const char *ord)
Create environment files for worker &#39;ord&#39;.
Definition: TProofLite.cxx:694
void SetQueryRunning(TProofQueryResult *pq)
Set query in running state.
virtual const char * BaseName(const char *pathname)
Base name of a file name. Base name of /user/root is root.
Definition: TSystem.cxx:929
virtual Int_t GetDrawArgs(const char *var, const char *sel, Option_t *opt, TString &selector, TString &objname)=0
virtual const char * GetName() const
Returns name of object.
Definition: TNamed.h:51
Long64_t GetEntries() const
Definition: TQueryResult.h:130
virtual Bool_t AccessPathName(const char *path, EAccessMode mode=kFileExists)
Returns FALSE if one can access a file using the specified access mode.
Definition: TSystem.cxx:1266
Int_t VerifyDataSet(const char *uri, const char *=0)
Verify if all files in the specified dataset are available.
This class starts a PROOF session on the local machine: no daemons, client and master merged...
Definition: TProofLite.h:42
Bool_t RequestStagingDataSet(const char *dataset)
Allows users to request staging of a particular dataset.
virtual TList * GetInputList() const =0
The PROOF package manager contains tools to manage packages.
Definition: TPackMgr.h:47
virtual Bool_t IsValid() const
Definition: TSocket.h:162
void GetEnabledPackages(TString &packlist)
Method to get a semi-colon separated list with the names of the enabled packages. ...
Definition: TPackMgr.cxx:515
TQueryResultManager * fQMgr
Definition: TProofLite.h:64
virtual int GetPid()
Get process id.
Definition: TSystem.cxx:712
TString fPerfTree
Definition: TProof.h:587
virtual void Delete(Option_t *option="")
Remove all objects from the list AND delete all heap based objects.
Definition: TList.cxx:405
void ActivateAsyncInput()
Activate the a-sync input handler.
Definition: TProof.cxx:4382
virtual void Info(const char *method, const char *msgfmt,...) const
Issue info message.
Definition: TObject.cxx:899
Int_t GetNumberOfActiveSlaves() const
Return number of active slaves, i.e.
Definition: TProof.cxx:1965
void AskParallel()
Ask the for the number of parallel slaves.
Definition: TProof.cxx:2055
void SetPort(Int_t port)
Definition: TUrl.h:97
TPMERegexp * fReInvalid
Definition: TProofLite.h:69
Long64_t GetBytesRead() const
Definition: TProof.h:959
Double_t RealTime()
Stop the stopwatch (if it is running) and return the realtime (in seconds) passed between the start a...
Definition: TStopwatch.cxx:110
void ShowDataSetCache(const char *dataset=0)
Display the content of the dataset cache, if any (matching &#39;dataset&#39;, if defined).
virtual const char * GetBuildCompilerVersion() const
Return the build compiler version.
Definition: TSystem.cxx:3769
TMonitor * fAllUniqueMonitor
Definition: TProof.h:516
virtual Int_t ClearCache(const char *uri)
Clear cached information matching uri.
virtual void AddInput(TObject *inp)=0
long long Long64_t
Definition: RtypesCore.h:69
TFileCollection * GetDataSet(const char *uri, const char *=0)
Get a list of TFileInfo objects describing the files of the specified dataset.
virtual TDSetElement * Next(Long64_t totalEntries=-1)
Returns next TDSetElement.
Definition: TDSet.cxx:394
void PrepareInputDataFile(TString &dataFile)
Prepare the file with the input data objects to be sent the master; the objects are taken from the de...
Definition: TProof.cxx:9611
Int_t Load(const char *macro, Bool_t notOnClient=kFALSE, Bool_t uniqueOnly=kTRUE, TList *wrks=0)
Copy the specified macro in the cache directory.
void SetPerfTree(const char *pf="perftree.root", Bool_t withWrks=kFALSE)
Enable/Disable saving of the performance tree.
Definition: TProof.cxx:12597
virtual const char * WorkingDirectory()
Return working directory.
Definition: TSystem.cxx:866
static TMD5 * FileChecksum(const char *file)
Returns checksum of specified file.
Definition: TMD5.cxx:474
void SetProtocol(const char *proto, Bool_t setDefaultPort=kFALSE)
Set protocol and, optionally, change the port accordingly.
Definition: TUrl.cxx:520
TString fConfDir
Definition: TProof.h:599
Bool_t ExistsDataSet(const char *uri)
Returns kTRUE if &#39;dataset&#39; described by &#39;uri&#39; exists, kFALSE otherwise.
const char *const kCP
Definition: TProof.h:171
const char *const kLS
Definition: TProof.h:173
Int_t fOtherQueries
Definition: TProof.h:553
Collectable string class.
Definition: TObjString.h:32
float Float_t
Definition: RtypesCore.h:53
virtual TVirtualProofPlayer * MakePlayer(const char *player=0, TSocket *s=0)
Construct a TProofPlayer object.
Definition: TProof.cxx:10183
const char Option_t
Definition: RtypesCore.h:62
virtual ~TProofLite()
Destructor.
Definition: TProofLite.cxx:375
Bool_t fSync
Definition: TProof.h:536
virtual TString SplitAclicMode(const char *filename, TString &mode, TString &args, TString &io) const
This method split a filename of the form: ~~~ {.cpp} [path/]macro.C[+|++[k|f|g|O|c|s|d|v|-]][(args)]...
Definition: TSystem.cxx:4109
Bool_t RegisterDataSet(const char *dsName, TFileCollection *ds, const char *opt="")
Register the &#39;dataSet&#39; on the cluster under the current user, group and the given &#39;dataSetName&#39;...
TFileCollection * GetStagingStatusDataSet(const char *dataset)
Obtains a TFileCollection showing the staging status of the specified dataset.
TString & ReplaceAll(const TString &s1, const TString &s2)
Definition: TString.h:635
void SetupWorkersEnv(TList *wrks, Bool_t increasingpool=kFALSE)
Set up packages, loaded macros, include and lib paths ...
Definition: TProof.cxx:1512
void SendInputDataFile()
Make sure that the input data objects are available to the workers in a dedicated file in the cache; ...
int GetPathInfo(const char *path, Long_t *id, Long_t *size, Long_t *flags, Long_t *modtime)
Get info about a file: id, size, flags, modification time.
Definition: TSystem.cxx:1364
virtual Long64_t DrawSelect(TDSet *set, const char *varexp, const char *selection, Option_t *option="", Long64_t nentries=-1, Long64_t firstentry=0)=0
Bool_t TestBit(UInt_t f) const
Definition: TObject.h:157
virtual Int_t Recv(TMessage *&mess)
Receive a TMessage object.
Definition: TSocket.cxx:818
This class implements a data set to be used for PROOF processing.
Definition: TDSet.h:153
virtual void SetName(const char *name)
Set the name of the TNamed.
Definition: TNamed.cxx:131
virtual void SetOwner(Bool_t enable=kTRUE)
Set whether this collection is the owner (enable==true) of its content.
TH1 * h
Definition: legend2.C:5
void SetParameter(const char *par, const char *value)
Set input list parameter.
Definition: TProof.cxx:9794
virtual Bool_t RemoveDataSet(const char *uri)
Removes the indicated dataset.
The PROOF manager interacts with the PROOF server coordinator to create or destroy a PROOF session...
Definition: TProofMgr.h:53
void SetUrl(const char *url, Bool_t defaultIsFile=kFALSE)
Parse url character string and split in its different subcomponents.
Definition: TUrl.cxx:110
static const TList * GetEnvVars()
Get environemnt variables.
Definition: TProof.cxx:11723
Int_t LockSession(const char *sessiontag, TProofLockPath **lck)
Try locking query area of session tagged sessiontag.
A ROOT file is a suite of consecutive data records (TKey instances) with a well defined format...
Definition: TFile.h:50
virtual EExitStatus GetExitStatus() const =0
virtual TObject * Last() const
Return the last object in the list. Returns 0 when list is empty.
Definition: TList.cxx:581
TList * fAllUniqueSlaves
Definition: TProof.h:512
virtual int MakeDirectory(const char *name)
Make a directory.
Definition: TSystem.cxx:822
virtual void AddSignalHandler(TSignalHandler *sh)
Add a signal handler to list of system signal handlers.
Definition: TSystem.cxx:537
virtual const char * HomeDirectory(const char *userName=0)
Return the user&#39;s home directory.
Definition: TSystem.cxx:882
TString fSelection
Definition: TProofLite.h:60
TString fImage
Definition: TProof.h:600
virtual TFileCollection * GetDataSet(const char *uri, const char *server=0)
Utility function used in various methods for user dataset upload.
Int_t PollForNewWorkers()
Simulate dynamic addition, for test purposes.
virtual TObject * Get(const char *namecycle)
Return pointer to object identified by namecycle.
void SetSocket(TSocket *s)
Definition: TSlave.h:116
#define R__ASSERT(e)
Definition: TError.h:98
#define gROOT
Definition: TROOT.h:364
TString fLogFileName
Definition: TProof.h:541
Ssiz_t Index(const char *pat, Ssiz_t i=0, ECaseCompare cmp=kExact) const
Definition: TString.h:582
TProofLockPath * fCacheLock
Definition: TProofLite.h:62
virtual void Add(TSocket *sock, Int_t interest=kRead)
Add socket to the monitor&#39;s active list.
Definition: TMonitor.cxx:168
Int_t LoadPlugin()
Load the plugin library for this handler.
virtual void SetCurrentQuery(TQueryResult *q)=0
TList * Queries() const
TList * fQueries
Definition: TProof.h:552
The TEnv class reads config files, by default named .rootrc.
Definition: TEnv.h:128
Basic string class.
Definition: TString.h:137
void ClearDataSetCache(const char *dataset=0)
Clear the content of the dataset cache, if any (matching &#39;dataset&#39;, if defined).
Int_t GetNumberOfSlaves() const
Return number of slaves as described in the config file.
Definition: TProof.cxx:1956
Int_t GetClientProtocol() const
Definition: TProof.h:944
Int_t SetDataSetTreeName(const char *dataset, const char *treename)
Set/Change the name of the default tree.
int Int_t
Definition: RtypesCore.h:41
virtual const char * DirName(const char *pathname)
Return the directory name in pathname.
Definition: TSystem.cxx:997
bool Bool_t
Definition: RtypesCore.h:59
void NotifyStartUp(const char *action, Int_t done, Int_t tot)
Notify setting-up operation message.
Definition: TProofLite.cxx:677
TQueryResult * CloneInfo()
Return an instance of TQueryResult containing only the local info fields, i.e.
TList * fUniqueSlaves
Definition: TProof.h:511
virtual Bool_t JoinProcess(TList *workers)=0
R__EXTERN TVirtualMutex * gROOTMutex
Definition: TROOT.h:63
const Bool_t kFALSE
Definition: Rtypes.h:92
void SetUser(const char *user)
Definition: TUrl.h:91
TList * fWaitingSlaves
Definition: TProof.h:551
Int_t Broadcast(const TMessage &mess, TList *slaves)
Broadcast a message to all slaves in the specified list.
Definition: TProof.cxx:2453
const char *const kRM
Definition: TProof.h:172
virtual void RemoveAll()
Remove all sockets from the monitor.
Definition: TMonitor.cxx:241
virtual void ShowDataSets(const char *uri="*", const char *opt="")
Prints formatted information about the dataset &#39;uri&#39;.
This class represents a RFC 3986 compatible URI.
Definition: TUri.h:39
Long_t fMtime
Definition: TSystem.h:142
static Int_t fgWrksMax
Definition: TProofLite.h:71
Int_t GetSeqNum() const
Definition: TQueryResult.h:123
R__EXTERN TApplication * gApplication
Definition: TApplication.h:171
virtual void DeActivateAll()
De-activate all activated sockets.
Definition: TMonitor.cxx:302
void ShowData()
List contents of the data directory in the sandbox.
void ShowDataDir(const char *dirname)
List contents of the data directory &#39;dirname&#39;.
void ScanPreviousQueries(const char *dir)
Scan the queries directory for the results of previous queries.
TLatex * t1
Definition: textangle.C:20
static void ResolveKeywords(TString &fname, const char *path=0)
Replace <ord>, <user>, <u>, <group>, <stag>, <qnum>, <file>, <rver> and <build> placeholders in fname...
TString & Insert(Ssiz_t pos, const char *s)
Definition: TString.h:592
TList * fInputData
Definition: TProof.h:565
Int_t SendCurrentState(ESlaves list=kActive)
Transfer the current state of the master to the active slave servers.
Definition: TProof.cxx:6730
TList * GetListOfElements() const
Definition: TDSet.h:231
TVirtualProofPlayer * fPlayer
Definition: TProof.h:524
Bool_t R_ISREG(Int_t mode)
Definition: TSystem.h:129
void SetBit(UInt_t f, Bool_t set)
Set or unset the user status bits as specified in f.
Definition: TObject.cxx:739
void ResolveKeywords(TString &s, const char *ord, const char *logfile)
Resolve some keywords in &#39;s&#39; <logfilewrk>, <user>, <rootsys>, <cpupin>
Definition: TProofLite.cxx:819
virtual TObject * FindObject(const char *name) const
Find an object in this list using its name.
Definition: TList.cxx:497
TString & Replace(Ssiz_t pos, Ssiz_t n, const char *s)
Definition: TString.h:625
static const char * GetMacroPath()
Get macro search path. Static utility function.
Definition: TROOT.cxx:2563
static TFile * Open(const char *name, Option_t *option="", const char *ftitle="", Int_t compress=1, Int_t netopt=0)
Create / open a file.
Definition: TFile.cxx:3907
Int_t fMode
Definition: TSystem.h:138
Bool_t fForkStartup
Definition: TProofLite.h:54
Bool_t fMasterServ
Definition: TProof.h:596
virtual Long64_t Process(TDSet *set, const char *selector, Option_t *option="", Long64_t nentries=-1, Long64_t firstentry=0)=0
Int_t SavePerfTree(const char *pf=0, const char *qref=0)
Save performance information from TPerfStats to file &#39;pf&#39;.
Definition: TProof.cxx:12619
TDataSetManagerFile * fDataSetStgRepo
Definition: TProofLite.h:67
TSignalHandler * fIntHandler
Definition: TProof.h:521
Manages an element of a TDSet.
Definition: TDSet.h:68
Int_t fDrawQueries
Definition: TProof.h:554
TList * fChains
Definition: TProof.h:526
static struct mg_connection * fc(struct mg_context *ctx)
Definition: civetweb.c:1956
Int_t Update(Long64_t avgsize=-1)
Update accumulated information about the elements of the collection (e.g.
virtual const char * GetDirEntry(void *dirp)
Get a directory entry. Returns 0 if no more entries.
Definition: TSystem.cxx:848
#define SafeDelete(p)
Definition: RConfig.h:507
TPluginHandler * fProgressDialog
Definition: TProof.h:522
TList * fBadSlaves
Definition: TProof.h:604
virtual TList * GetOutputList() const =0
virtual int Unlink(const char *name)
Unlink, i.e. remove, a file.
Definition: TSystem.cxx:1347
TString fDataSetDir
Definition: TProofLite.h:51
Bool_t fSendGroupView
list returned by kPROOF_GETSLAVEINFO
Definition: TProof.h:504
void Stop()
Stop the stopwatch.
Definition: TStopwatch.cxx:77
virtual Bool_t IsValid() const
Definition: TSlave.h:154
TObject * GetEntryList() const
Definition: TDSet.h:251
static TString Format(const char *fmt,...)
Static method which formats a string using a printf style format descriptor and return a TString...
Definition: TString.cxx:2335
virtual int GetSysInfo(SysInfo_t *info) const
Returns static system info, like OS type, CPU type, number of CPUs RAM size, etc into the SysInfo_t s...
Definition: TSystem.cxx:2440
#define PDB(mask, level)
Definition: TProofDebug.h:58
void UpdateDialog()
Final update of the progress dialog.
Definition: TProof.cxx:4325
TDataSetManager * fDataSetManager
Definition: TProofLite.h:66
TSlave * CreateSlave(const char *url, const char *ord, Int_t perf, const char *image, const char *workdir)
Create a new TSlave of type TSlave::kSlave.
Definition: TProof.cxx:1831
TSocket * GetSocket() const
Definition: TSlave.h:138
This code implements the MD5 message-digest algorithm.
Definition: TMD5.h:46
TString fCacheDir
Definition: TProofLite.h:49
The TNamed class is the base class for all named ROOT classes.
Definition: TNamed.h:33
Long64_t Process(TDSet *dset, const char *sel, Option_t *o="", Long64_t nent=-1, Long64_t fst=0)
Process a data set (TDSet) using the specified selector (.C) file.
TList * PreviousQueries() const
Int_t GetLogLevel() const
Definition: TProof.h:946
const char *const kPROOF_QueryDir
Definition: TProof.h:158
Int_t Init(const char *masterurl, const char *conffile, const char *confdir, Int_t loglevel, const char *alias=0)
Start the PROOF environment.
Definition: TProofLite.cxx:154
virtual Int_t RegisterDataSet(const char *uri, TFileCollection *dataSet, const char *opt)
Register a dataset, perfoming quota checkings, if needed.
Int_t fParallel
Definition: TSlave.h:101
Bool_t CancelStagingDataSet(const char *dataset)
Cancels a dataset staging request.
virtual const char * Getenv(const char *env)
Get environment variable.
Definition: TSystem.cxx:1628
Int_t Collect(const TSlave *sl, Long_t timeout=-1, Int_t endtype=-1, Bool_t deactonfail=kFALSE)
Collect responses from slave sl.
Definition: TProof.cxx:2647
static Int_t RegisterDataSets(TList *in, TList *out, TDataSetManager *dsm, TString &e)
Register TFileCollections in &#39;out&#39; as datasets according to the rules in &#39;in&#39;.
Int_t ApplyMaxQueries(Int_t mxq)
Scan the queries directory and remove the oldest ones (and relative dirs, if empty) in such a way onl...
A sorted doubly linked list.
Definition: TSortedList.h:30
TString fConfFile
Definition: TProof.h:598
std::vector< std::vector< double > > Data
Bool_t EndsWith(const char *pat, ECaseCompare cmp=kExact) const
Return true if string ends with the specified string.
Definition: TString.cxx:2221
virtual UserGroup_t * GetUserInfo(Int_t uid)
Returns all user info in the UserGroup_t structure.
Definition: TSystem.cxx:1564
TQueryResult * GetQueryResult(const char *ref=0)
Return pointer to the full TQueryResult instance owned by the player and referenced by &#39;ref&#39;...
Definition: TProof.cxx:2126
TList * fSlaves
Definition: TProof.h:602
Int_t fNotIdle
Definition: TProof.h:535
virtual void SetOutputList(TList *out, Bool_t adopt=kTRUE)
Set / change the output list.
virtual TSocket * Accept(UChar_t Opt=0)
Accept a connection on a server socket.
virtual const char * TempDirectory() const
Return a user configured or systemwide directory to create temporary files in.
Definition: TSystem.cxx:1448
const Bool_t kSortDescending
Definition: TList.h:41
friend class TProofInputHandler
Definition: TProof.h:354
const char * GetDataDir() const
Definition: TProof.h:259
TList * fInactiveSlaves
Definition: TProof.h:510
A container class for query results.
Definition: TQueryResult.h:44
Int_t fCollectTimeout
Definition: TProof.h:613
TList * GetListOfSlaveInfos()
Returns list of TSlaveInfo&#39;s. In case of error return 0.
Definition: TProof.cxx:2299
Int_t fDynamicStartupNMax
Definition: TProofLite.h:57
static TList * fgProofEnvList
Definition: TProof.h:575
Float_t GetRealTime() const
Definition: TProof.h:960
TSocket * Select()
Return pointer to socket for which an event is waiting.
Definition: TMonitor.cxx:322
virtual void SetProcessInfo(Long64_t ent, Float_t cpu=0., Long64_t siz=-1, Float_t inittime=0., Float_t proctime=0.)
Set processing info.
void Emit(const char *signal)
Acitvate signal without args.
Definition: TQObject.cxx:561
TString fQueryDir
Definition: TProofLite.h:50
void SetTermTime(Float_t termtime)
Definition: TQueryResult.h:108
void ParseConfigField(const char *config)
The config file field may contain special instructions which need to be parsed at the beginning...
Definition: TProof.cxx:1031
TString fSandbox
Definition: TProofLite.h:48
Int_t GetNumberOfInactiveSlaves() const
Return number of inactive slaves, i.e.
Definition: TProof.cxx:1974
A doubly linked list.
Definition: TList.h:47
const char * GetObjName() const
Definition: TDSet.h:122
const char * GetUser() const
Definition: TUrl.h:74
Bool_t fRedirLog
Definition: TProof.h:540
TMonitor * fUniqueMonitor
Definition: TProof.h:515
const char *const kPROOF_ConfFile
Definition: TProof.h:152
Bool_t RemoveDataSet(const char *group, const char *user, const char *dsName)
Removes the indicated dataset.
Int_t RemoveDataSet(const char *uri, const char *=0)
Remove the specified dataset from the PROOF cluster.
virtual Int_t ReadFile(const char *fname, EEnvLevel level)
Read and parse the resource file for a certain level.
Definition: TEnv.cxx:597
const char *const kPROOF_QueryLockFile
Definition: TProof.h:163
TMonitor * fAllMonitor
Definition: TProof.h:605
virtual void AddQueryResult(TQueryResult *q)=0
void InitMembers()
Default initializations.
Definition: TProof.cxx:524
void ShowDataSets(const char *uri="", const char *=0)
Shows datasets in locations that match the uri By default shows the user&#39;s datasets and global ones...
const char * GetDir() const
Definition: TPackMgr.h:87
virtual TMap * GetSubDataSets(const char *uri, const char *excludeservers)
Partition dataset &#39;ds&#39; accordingly to the servers.
TString fUser
Definition: TSystem.h:152
virtual const char * GetBuildArch() const
Return the build architecture.
Definition: TSystem.cxx:3753
Int_t fLogLevel
Definition: TProof.h:499
TProofMgr * fManager
Definition: TProof.h:617
TList * fActiveSlaves
Definition: TProof.h:507
Long64_t DrawSelect(TDSet *dset, const char *varexp, const char *selection="", Option_t *option="", Long64_t nentries=-1, Long64_t firstentry=0)
Execute the specified drawing action on a data set (TDSet).
void QueryResultReady(const char *ref)
Notify availability of a query result.
Definition: TProof.cxx:9341
const TString GetUri() const
Returns the whole URI - an implementation of chapter 5.3 component recomposition. ...
Definition: TUri.cxx:140
Int_t fNWorkers
Definition: TProofLite.h:47
Bool_t fValid
Definition: TProof.h:494
virtual void Setenv(const char *name, const char *value)
Set environment variable.
Definition: TSystem.cxx:1612
Int_t fCpus
Definition: TSystem.h:165
virtual TObject * First() const
Return the first object in the list. Returns 0 when list is empty.
Definition: TList.cxx:557
TRandom2 r(17)
Class managing the query-result area.
R__EXTERN TSystem * gSystem
Definition: TSystem.h:549
TString fWorkDir
Definition: TProof.h:497
SVector< double, 2 > v
Definition: Dict.h:5
if object ctor succeeded but object should not be used
Definition: TObject.h:63
const char *const kPROOF_DataSetDir
Definition: TProof.h:159
TMonitor * fCurrentMonitor
Definition: TProof.h:517
THashList * GetList()
Long_t ExecPlugin(int nargs, const T &... params)
Int_t SetupWorkers(Int_t opt=0, TList *wrks=0)
Start up PROOF workers.
Definition: TProofLite.cxx:499
Long64_t GetNFiles() const
virtual Int_t GetValue(const char *name, Int_t dflt)
Returns the integer value for a resource.
Definition: TEnv.cxx:496
virtual TObject * Remove(TObject *obj)
Remove object from the list.
Definition: TList.cxx:675
const char *const kPROOF_ConfDir
Definition: TProof.h:153
Bool_t Gets(FILE *fp, Bool_t chop=kTRUE)
Read one line from the stream, including the , or until EOF.
Definition: Stringio.cxx:198
virtual Bool_t ExistsDataSet(const char *uri)
Checks if the indicated dataset exits.
Bool_t IsDraw() const
Definition: TQueryResult.h:152
Bool_t ParseUri(const char *uri, TString *dsGroup=0, TString *dsUser=0, TString *dsName=0, TString *dsTree=0, Bool_t onlyCurrent=kFALSE, Bool_t wildcards=kFALSE)
Parses a (relative) URI that describes a DataSet on the cluster.
Bool_t BeginsWith(const char *s, ECaseCompare cmp=kExact) const
Definition: TString.h:558
Int_t fSeqNum
Definition: TProof.h:556
void ClearCache(const char *file=0)
Remove files from all file caches.
Int_t fDynamicStartupStep
Definition: TProofLite.h:56
virtual Int_t ShowCache(const char *uri)
Show cached information matching uri.
void SetActive(Bool_t=kTRUE)
Definition: TProof.h:1018
Int_t WriteDataSet(const char *group, const char *user, const char *dsName, TFileCollection *dataset, UInt_t option=0, TMD5 *checksum=0)
Writes indicated dataset.
void Form(const char *fmt,...)
Formats a string using a printf style format descriptor.
Definition: TString.cxx:2322
TList * fSlaveInfo
Definition: TProof.h:503
unsigned int UInt_t
Definition: RtypesCore.h:42
TList * GetInputList()
Definition: TQueryResult.h:128
TMarker * m
Definition: textangle.C:8
const char *const kPROOF_CacheLockFile
Definition: TProof.h:161
Int_t ScanDataSet(const char *uri, const char *opt)
Scans the dataset indicated by &#39;uri&#39; following the &#39;opts&#39; directives.
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition: TObject.cxx:925
char * Form(const char *fmt,...)
Ssiz_t Length() const
Definition: TString.h:390
void SetRunning(Int_t startlog, const char *par, Int_t nwrks)
Call when running starts.
Int_t fSeqNum
query unique sequential number
Definition: TQueryResult.h:60
TServerSocket * fServSock
Definition: TProofLite.h:53
Bool_t SetFragment(const TString &fragment)
Set fragment component of URI: fragment = *( pchar / "/" / "?" ).
Definition: TUri.cxx:498
TList * GetListOfQueries(Option_t *opt="")
Get the list of queries.
Bool_t fProgressDialogStarted
Definition: TProof.h:523
virtual Int_t Exec(const char *shellcmd)
Execute a command.
Definition: TSystem.cxx:658
const Int_t kPROOF_Protocol
Definition: TProof.h:150
Int_t SendGroupView()
Send to all active slaves servers the current slave group size and their unique id.
Definition: TProof.cxx:6432
TStopwatch fQuerySTW
Definition: TProof.h:623
Bool_t FinalizeQuery(TProofQueryResult *pq, TProof *proof, TVirtualProofPlayer *player)
Final steps after Process() to complete the TQueryResult instance.
Int_t RemoveWorkers(TList *wrks)
Used for shuting down the workres after a query is finished.
Definition: TProof.cxx:1575
FILE * fLogFileR
Definition: TProof.h:543
void SetName(const char *name)
Definition: TCollection.h:116
Int_t fProtocol
Definition: TProof.h:601
static Int_t GetNumberOfWorkers(const char *url=0)
Static method to determine the number of workers giving priority to users request.
Definition: TProofLite.cxx:406
Bool_t fLogToWindowOnly
Definition: TProof.h:544
TList * fFeedback
Definition: TProof.h:525
Bool_t IsIdle() const
Definition: TProof.h:970
virtual void FreeDirectory(void *dirp)
Free a directory.
Definition: TSystem.cxx:840
void AddInput(TObject *obj)
Add objects that might be needed during the processing of the selector (see Process()).
Definition: TProof.cxx:9706
#define Printf
Definition: TGeoToOCC.h:18
TMap * GetDataSets(const char *uri="", const char *=0)
lists all datasets that match given uri
TList * fRunningDSets
Definition: TProof.h:611
#define R__LOCKGUARD2(mutex)
Int_t VerifyDataSetParallel(const char *uri, const char *optStr)
Internal function for parallel dataset verification used TProof::VerifyDataSet and TProofLite::Verify...
Definition: TProof.cxx:11153
Bool_t fDynamicStartup
Definition: TProof.h:619
virtual void RemoveQueryResult(const char *ref)=0
static TSelector * GetSelector(const char *filename)
The code in filename is loaded (interpreted or compiled, see below), filename must contain a valid cl...
Definition: TSelector.cxx:142
Int_t AssertPath(const char *path, Bool_t writable)
Make sure that &#39;path&#39; exists; if &#39;writable&#39; is kTRUE, make also sure that the path is writable...
Definition: TProof.cxx:1253
void SetHost(const char *host)
Definition: TUrl.h:93
Int_t CreateSandbox()
Create the sandbox for this session.
Definition: TProofLite.cxx:898
TString & Remove(Ssiz_t pos)
Definition: TString.h:616
long Long_t
Definition: RtypesCore.h:50
int Ssiz_t
Definition: RtypesCore.h:63
TString fSockPath
Definition: TProofLite.h:52
Int_t HandleOutputOptions(TString &opt, TString &target, Int_t action)
Extract from opt information about output handling settings.
Definition: TProof.cxx:4909
R__EXTERN TProof * gProof
Definition: TProof.h:1107
virtual TSignalHandler * RemoveSignalHandler(TSignalHandler *sh)
Remove a signal handler from list of signal handlers.
Definition: TSystem.cxx:547
Bool_t fEndMaster
Definition: TProof.h:560
TObjArray * Tokenize(const TString &delim) const
This function is used to isolate sequential tokens in a TString.
Definition: TString.cxx:2241
void SetFeedback(TString &opt, TString &optfb, Int_t action)
Extract from opt in optfb information about wanted feedback settings.
Definition: TProof.cxx:5204
#define ClassImp(name)
Definition: Rtypes.h:279
double f(double x)
Int_t GoParallel(Int_t nodes, Bool_t accept=kFALSE, Bool_t random=kFALSE)
Go in parallel mode with at most "nodes" slaves.
Definition: TProof.cxx:7245
TSignalHandler * GetSignalHandler() const
Definition: TApplication.h:112
virtual const char * HostName()
Return the system&#39;s host name.
Definition: TSystem.cxx:308
virtual int Symlink(const char *from, const char *to)
Create a symbolic link from file1 to file2.
Definition: TSystem.cxx:1338
Ssiz_t Last(char c) const
Find last occurrence of a character c.
Definition: TString.cxx:865
Int_t Lock()
Locks the directory.
TList * fEnabledPackages
Definition: TProof.h:610
EQueryMode GetQueryMode(Option_t *mode=0) const
Find out the query mode based on the current setting and &#39;mode&#39;.
Definition: TProof.cxx:6091
Int_t GetParallel() const
Returns number of slaves active in parallel mode.
Definition: TProof.cxx:2282
TList * fNonUniqueMasters
Definition: TProof.h:513
TMap implements an associative array of (key,value) pairs using a THashTable for efficient retrieval ...
Definition: TMap.h:44
R__EXTERN TEnv * gEnv
Definition: TEnv.h:174
Int_t CleanupSandbox()
Remove old sessions dirs keep at most &#39;Proof.MaxOldSessions&#39; (default 10)
TList * fRecvMessages
Definition: TProof.h:502
TNamed()
Definition: TNamed.h:40
Int_t GetNumberOfBadSlaves() const
Return number of bad slaves.
Definition: TProof.cxx:1992
TList * fEnabledPackagesOnCluster
Definition: TProof.h:563
int nentries
Definition: THbookFile.cxx:89
Bool_t Contains(const char *pat, ECaseCompare cmp=kExact) const
Definition: TString.h:567
you should not use this method at all Int_t Int_t Double_t Double_t Double_t e
Definition: TRolke.cxx:630
virtual void Reset()
Reset or initialize access to the elements.
Definition: TDSet.cxx:1350
void SetRunStatus(ERunStatus rst)
Definition: TProof.h:702
const char * GetFileName() const
Definition: TDSet.h:113
static Int_t RegisterGlobalPath(const char *paths)
Parse one or more paths as possible sources of packages Returns number of paths added; or -1 in case ...
Definition: TPackMgr.cxx:882
void RemoveQuery(TQueryResult *qr, Bool_t soft=kFALSE)
Remove everything about query qr.
Int_t Unlock()
Unlock the directory.
virtual void StopProcess(Bool_t abort, Int_t timeout=-1)=0
virtual void Clear(Option_t *option="")
Remove all objects from the list.
Definition: TList.cxx:349
Int_t GetSandbox(TString &sb, Bool_t assert=kFALSE, const char *rc=0)
Set the sandbox path from &#39; Proof.Sandbox&#39; or the alternative var &#39;rc&#39;.
Definition: TProof.cxx:1004
TQueryResult version adapted to PROOF neeeds.
EQueryMode fQueryMode
Definition: TProof.h:618
static TMap * GetDataSetNodeMap(TFileCollection *fc, TString &emsg)
Get a map {server-name, list-of-files} for collection &#39;fc&#39; to be used in TPacketizerFile.
virtual Long64_t GetEntries() const
Definition: TTree.h:393
Bool_t IsNull() const
Definition: TString.h:387
void SaveQuery(TProofQueryResult *qr, const char *fout=0)
Save current status of query &#39;qr&#39; to file name fout.
Int_t Match(const TString &s, UInt_t start=0)
Runs a match on s against the regex &#39;this&#39; was created with.
Definition: TPRegexp.cxx:705
FILE * fLogFileW
Definition: TProof.h:542
Mother of all ROOT objects.
Definition: TObject.h:37
TPackMgr * fPackMgr
Definition: TProof.h:562
TSelector * fSelector
Definition: TProof.h:621
TString fVarExp
Definition: TProofLite.h:59
Bool_t ExistsDataSet(const char *group, const char *user, const char *dsName)
Checks if the indicated dataset exits.
Float_t GetCpuTime() const
Definition: TProof.h:961
Int_t InitDataSetManager()
Initialize the dataset manager from directives or from defaults Return 0 on success, -1 on failure.
Bool_t R_ISDIR(Int_t mode)
Definition: TSystem.h:126
static void SetMacroPath(const char *newpath)
Set or extend the macro search path.
Definition: TROOT.cxx:2597
R__EXTERN TProofServ * gProofServ
Definition: TProofServ.h:361
virtual void Add(TObject *obj)
Definition: TList.h:81
const Ssiz_t kNPOS
Definition: Rtypes.h:115
Wrapper for PCRE library (Perl Compatible Regular Expressions).
Definition: TPRegexp.h:103
Class that contains a list of TFileInfo&#39;s and accumulated meta data information about its entries...
Definition: file.py:1
TProofOutputList fOutputList
Definition: TProof.h:568
R__EXTERN const char * gRootDir
Definition: TSystem.h:233
TProofLockPath * fQueryLock
Definition: TProofLite.h:63
TFileCollection * GetDataSet(const char *uri, const char *srv=0)
Utility function used in various methods for user dataset upload.
Int_t CleanupQueriesDir()
Remove all queries results referring to previous sessions.
Int_t fSessionID
Definition: TProof.h:558
void AskStatistics()
Ask the for the statistics of the slaves.
Definition: TProof.cxx:2000
TF1 * f1
Definition: legend1.C:11
#define snprintf
Definition: civetweb.c:822
TTree * GetTreeHeader(TDSet *tdset)
Creates a tree header (a tree with nonexisting files) object for the DataSet.
Int_t fStatus
Definition: TProof.h:500
Int_t SetParallel(Int_t nodes=-1, Bool_t random=kFALSE)
Tell PROOF how many slaves to use in parallel.
Definition: TProof.cxx:7112
virtual int CopyFile(const char *from, const char *to, Bool_t overwrite=kFALSE)
Copy a file.
Definition: TSystem.cxx:1311
virtual void * OpenDirectory(const char *name)
Open a directory. Returns 0 if directory does not exist.
Definition: TSystem.cxx:831
Int_t CountChar(Int_t c) const
Return number of times character c occurs in the string.
Definition: TString.cxx:444
R__EXTERN Int_t gDebug
Definition: Rtypes.h:128
Int_t Atoi() const
Return integer value of string.
Definition: TString.cxx:1965
void Reset()
Definition: TStopwatch.h:54
Bool_t IsParallel() const
Definition: TProof.h:969
const char * GetUser() const
Definition: TProof.h:936
Bool_t IsDigit() const
Returns true if all characters in string are digits (0-9) or white spaces, i.e.
Definition: TString.cxx:1807
TList * fAvailablePackages
Definition: TProof.h:609
A TTree object has a header with a name and a title.
Definition: TTree.h:98
const AParamType & GetVal() const
Definition: TParameter.h:77
double result[121]
Class describing a generic file including meta information.
Definition: TFileInfo.h:50
Int_t SendInitialState()
Transfer the initial (i.e.
Definition: TProof.cxx:6746
void ResetBit(UInt_t f)
Definition: TObject.h:156
TProofMgr::EServType fServType
Definition: TProof.h:616
Bool_t IsValid() const
Definition: TProof.h:967
virtual Bool_t ExpandPathName(TString &path)
Expand a pathname getting rid of special shell characters like ~.
Definition: TSystem.cxx:1244
TList * fTerminatedSlaveInfos
Definition: TProof.h:603
Definition: first.py:1
void Add(TObject *obj)
Add object in sorted list.
Definition: TSortedList.cxx:27
TProofQueryResult * MakeQueryResult(Long64_t nent, const char *opt, Long64_t fst, TDSet *dset, const char *selec)
Create a TProofQueryResult instance for this query.
TList * fLoadedMacros
Definition: TProof.h:574
Int_t Remove(const char *ref, Bool_t all)
Handle remove request.
virtual Int_t Load(const char *macro, Bool_t notOnClient=kFALSE, Bool_t uniqueOnly=kTRUE, TList *wrks=0)
Load the specified macro on master, workers and, if notOnClient is kFALSE, on the client...
Definition: TProof.cxx:8600
virtual Int_t GetSize() const
Definition: TCollection.h:95
Int_t Substitute(TString &s, const TString &r, Bool_t doDollarSubst=kTRUE)
Substitute matching part of s with r, dollar back-ref substitution is performed if doDollarSubst is t...
Definition: TPRegexp.cxx:871
Class describing a PROOF worker server.
Definition: TSlave.h:50
void FindUniqueSlaves()
Add to the fUniqueSlave list the active slaves that have a unique (user) file system image...
A TSelector object is used by the TTree::Draw, TTree::Scan, TTree::Process to navigate in a TTree and...
Definition: TSelector.h:39
const Bool_t kTRUE
Definition: Rtypes.h:91
void ResetProgressDialog(const char *sel, Int_t sz, Long64_t fst, Long64_t ent)
Reset progress dialog.
Definition: TProof.cxx:9271
virtual void SetTitle(const char *title="")
Set the title of the TNamed.
Definition: TNamed.cxx:155
TString fMaster
Definition: TProof.h:496
Bool_t fTty
Definition: TProof.h:495
TUrl fUrl
Definition: TProof.h:597
virtual TMap * GetDataSets(const char *uri, UInt_t=TDataSetManager::kExport)
Returns all datasets for the <group> and <user> specified by <uri>.
void ShowCache(Bool_t all=kFALSE)
List contents of file cache.
virtual char * ConcatFileName(const char *dir, const char *name)
Concatenate a directory and a file name. User must delete returned string.
Definition: TSystem.cxx:1045
void Print(Option_t *option="") const
Print status of PROOF-Lite cluster.
Definition: TProofLite.cxx:967
const char * GetOrdinal() const
Definition: TProofServ.h:267
const Int_t n
Definition: legend1.C:16
virtual Int_t SetupServ(Int_t stype, const char *conffile)
Init a PROOF slave object.
Definition: TSlave.cxx:179
TMonitor * fActiveMonitor
Definition: TProof.h:514
virtual TList * GetListOfResults() const =0
Int_t GetNumberOfUniqueSlaves() const
Return number of unique slaves, i.e.
Definition: TProof.cxx:1983
char name[80]
Definition: TGX11.cxx:109
virtual void Warning(const char *method, const char *msgfmt,...) const
Issue warning message.
Definition: TObject.cxx:911
Long64_t fLastPollWorkers_s
Definition: TProof.h:506
Int_t fMaxDrawQueries
Definition: TProof.h:555
virtual const char * GetTitle() const
Returns title of object.
Definition: TNamed.h:52
Int_t CopyMacroToCache(const char *macro, Int_t headerRequired=0, TSelector **selector=0, Int_t opt=0, TList *wrks=0)
Copy a macro, and its possible associated .h[h] file, to the cache directory, from where the workers ...
const char *const kPROOF_CacheDir
Definition: TProof.h:155
const char * Data() const
Definition: TString.h:349
const char *const kPROOF_PackDir
Definition: TProof.h:156
void SetInputHandler(TFileHandler *ih)
Adopt and register input handler for this slave.
Definition: TSlave.cxx:394
static Int_t AssertDataSet(TDSet *dset, TList *input, TDataSetManager *mgr, TString &emsg)
Make sure that dataset is in the form to be processed.
Definition: TProof.cxx:11987