ROOT logo

From $ROOTSYS/tutorials/net/parallelMergeClient.C

#include "TMessage.h"
#include "TBenchmark.h"
#include "TSocket.h"
#include "TH2.h"
#include "TTree.h"
#include "TParallelMergingFile.h"
#include "TRandom.h"
#include "TError.h"

void parallelMergeClient() 
{
   // Client program which creates and fills 2 histograms and a TTree. 
   // Every 1000000 fills the histograms and TTree is send to the server which displays the histogram.
   //
   // To run this demo do the following:
   //   - Open at least 2 windows
   //   - Start ROOT in the first windows
   //   - Execute in the first window: .x fastMergeServer.C
   //   - Execute in the other windows: root.exe -b -l -q .x treeClient.C
   //     (You can put it in the background if wanted).
   // If you want to run the hserv.C on a different host, just change
   // "localhost" in the TSocket ctor below to the desired hostname.
   //
   //Author: Fons Rademakers, Philippe Canal
   
   gBenchmark->Start("treeClient");

   TParallelMergingFile *file = (TParallelMergingFile*)TFile::Open("mergedClient.root?pmerge=localhost:1095","RECREATE");
   
   file->Write();
   file->UploadAndReset();       // We do this early to get assigned an index.
   UInt_t idx = file->fServerIdx; // This works on in ACLiC.

   TH1 *hpx;
   if (idx == 0) {
      // Create the histogram
      hpx = new TH1F("hpx","This is the px distribution",100,-4,4);
      hpx->SetFillColor(48);  // set nice fillcolor
   } else {
      hpx = new TH2F("hpxpy","py vs px",40,-4,4,40,-4,4);
   }
   Float_t px, py;
   TTree *tree = new TTree("tree","tree");
   tree->SetAutoFlush(4000000);
   tree->Branch("px",&px);
   tree->Branch("py",&py);
 
   // Fill histogram randomly
   gRandom->SetSeed();
   const int kUPDATE = 1000000;
   for (int i = 0; i < 25000000; ) {
      gRandom->Rannor(px,py);
      if (idx%2 == 0)
         hpx->Fill(px);
      else
         hpx->Fill(px,py);
      tree->Fill();
      ++i;
      if (i && (i%kUPDATE) == 0) {
         file->Write();
      }
   }
   file->Write();
   delete file;

   gBenchmark->Show("treeClient");
}
 parallelMergeClient.C:1
 parallelMergeClient.C:2
 parallelMergeClient.C:3
 parallelMergeClient.C:4
 parallelMergeClient.C:5
 parallelMergeClient.C:6
 parallelMergeClient.C:7
 parallelMergeClient.C:8
 parallelMergeClient.C:9
 parallelMergeClient.C:10
 parallelMergeClient.C:11
 parallelMergeClient.C:12
 parallelMergeClient.C:13
 parallelMergeClient.C:14
 parallelMergeClient.C:15
 parallelMergeClient.C:16
 parallelMergeClient.C:17
 parallelMergeClient.C:18
 parallelMergeClient.C:19
 parallelMergeClient.C:20
 parallelMergeClient.C:21
 parallelMergeClient.C:22
 parallelMergeClient.C:23
 parallelMergeClient.C:24
 parallelMergeClient.C:25
 parallelMergeClient.C:26
 parallelMergeClient.C:27
 parallelMergeClient.C:28
 parallelMergeClient.C:29
 parallelMergeClient.C:30
 parallelMergeClient.C:31
 parallelMergeClient.C:32
 parallelMergeClient.C:33
 parallelMergeClient.C:34
 parallelMergeClient.C:35
 parallelMergeClient.C:36
 parallelMergeClient.C:37
 parallelMergeClient.C:38
 parallelMergeClient.C:39
 parallelMergeClient.C:40
 parallelMergeClient.C:41
 parallelMergeClient.C:42
 parallelMergeClient.C:43
 parallelMergeClient.C:44
 parallelMergeClient.C:45
 parallelMergeClient.C:46
 parallelMergeClient.C:47
 parallelMergeClient.C:48
 parallelMergeClient.C:49
 parallelMergeClient.C:50
 parallelMergeClient.C:51
 parallelMergeClient.C:52
 parallelMergeClient.C:53
 parallelMergeClient.C:54
 parallelMergeClient.C:55
 parallelMergeClient.C:56
 parallelMergeClient.C:57
 parallelMergeClient.C:58
 parallelMergeClient.C:59
 parallelMergeClient.C:60
 parallelMergeClient.C:61
 parallelMergeClient.C:62
 parallelMergeClient.C:63
 parallelMergeClient.C:64
 parallelMergeClient.C:65
 parallelMergeClient.C:66
 parallelMergeClient.C:67
 parallelMergeClient.C:68