ROOT logo

From $ROOTSYS/tutorials/net/hclient.C

void hclient(Bool_t evol=kFALSE) 
{
   // Client program which creates and fills a histogram. Every 1000 fills
   // the histogram is send to the server which displays the histogram.
   //
   // To run this demo do the following:
   //   - Open three windows
   //   - Start ROOT in all three windows
   //   - Execute in the first window: .x hserv.C (or hserv2.C)
   //   - Execute in the second and third windows: .x hclient.C
   // If you want to run the hserv.C on a different host, just change
   // "localhost" in the TSocket ctor below to the desired hostname.
   //
   // The script argument "evol" can be used when using a modified version
   // of the script where the clients and server are on systems with
   // different versions of ROOT. When evol is set to kTRUE the socket will
   // support automatic schema evolution between the client and the server.
   //
   //Author: Fons Rademakers
   
   gBenchmark->Start("hclient");

   // Open connection to server
   TSocket *sock = new TSocket("localhost", 9090);

   // Wait till we get the start message
   char str[32];
   sock->Recv(str, 32);

   // server tells us who we are
   int idx = !strcmp(str, "go 0") ? 0 : 1;

   Float_t messlen  = 0;
   Float_t cmesslen = 0;
   if (idx == 1)
      sock->SetCompressionLevel(1);

   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);
   }

   TMessage::EnableSchemaEvolutionForAll(evol);
   TMessage mess(kMESS_OBJECT);
   //TMessage mess(kMESS_OBJECT | kMESS_ACK);

   // Fill histogram randomly
   gRandom->SetSeed();
   Float_t px, py;
   const int kUPDATE = 1000;
   for (int i = 0; i < 25000; i++) {
      gRandom->Rannor(px,py);
      if (idx == 0)
         hpx->Fill(px);
      else
         hpx->Fill(px,py);
      if (i && (i%kUPDATE) == 0) {
         mess.Reset();              // re-use TMessage object
         mess.WriteObject(hpx);     // write object in message buffer
         sock->Send(mess);          // send message
         messlen  += mess.Length();
         cmesslen += mess.CompLength();
      }
   }
   sock->Send("Finished");          // tell server we are finished

   if (cmesslen > 0)
      printf("Average compression ratio: %g\n", messlen/cmesslen);

   gBenchmark->Show("hclient");

   // Close the socket
   sock->Close();
}
 hclient.C:1
 hclient.C:2
 hclient.C:3
 hclient.C:4
 hclient.C:5
 hclient.C:6
 hclient.C:7
 hclient.C:8
 hclient.C:9
 hclient.C:10
 hclient.C:11
 hclient.C:12
 hclient.C:13
 hclient.C:14
 hclient.C:15
 hclient.C:16
 hclient.C:17
 hclient.C:18
 hclient.C:19
 hclient.C:20
 hclient.C:21
 hclient.C:22
 hclient.C:23
 hclient.C:24
 hclient.C:25
 hclient.C:26
 hclient.C:27
 hclient.C:28
 hclient.C:29
 hclient.C:30
 hclient.C:31
 hclient.C:32
 hclient.C:33
 hclient.C:34
 hclient.C:35
 hclient.C:36
 hclient.C:37
 hclient.C:38
 hclient.C:39
 hclient.C:40
 hclient.C:41
 hclient.C:42
 hclient.C:43
 hclient.C:44
 hclient.C:45
 hclient.C:46
 hclient.C:47
 hclient.C:48
 hclient.C:49
 hclient.C:50
 hclient.C:51
 hclient.C:52
 hclient.C:53
 hclient.C:54
 hclient.C:55
 hclient.C:56
 hclient.C:57
 hclient.C:58
 hclient.C:59
 hclient.C:60
 hclient.C:61
 hclient.C:62
 hclient.C:63
 hclient.C:64
 hclient.C:65
 hclient.C:66
 hclient.C:67
 hclient.C:68
 hclient.C:69
 hclient.C:70
 hclient.C:71
 hclient.C:72
 hclient.C:73
 hclient.C:74
 hclient.C:75
 hclient.C:76
 hclient.C:77
 hclient.C:78
 hclient.C:79