File: onecpu.c

package info (click to toggle)
hpcc 1.5.0-4
  • links: PTS, VCS
  • area: main
  • in suites: forky, sid, trixie
  • size: 4,752 kB
  • sloc: ansic: 27,044; makefile: 50; sh: 24
file content (126 lines) | stat: -rw-r--r-- 5,854 bytes parent folder | download | duplicates (3)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
/* -*- mode: C; tab-width: 2; indent-tabs-mode: nil; fill-column: 79; coding: iso-latin-1-unix -*- */


#include <hpcc.h>

int
HPCC_StarStream(HPCC_Params *params) {
  int myRank, commSize;
  int rv, errCount, failure = 0, failureAll = 0;
  double copyLocalGBs, copyMinGBs, copyMaxGBs, copyAvgGBs;
  double scaleLocalGBs, scaleMinGBs, scaleMaxGBs, scaleAvgGBs;
  double addLocalGBs, addMinGBs, addMaxGBs, addAvgGBs;
  double triadLocalGBs, triadMinGBs, triadMaxGBs, triadAvgGBs;
  FILE *outputFile;

  copyLocalGBs = copyMinGBs = copyMaxGBs = copyAvgGBs =
  scaleLocalGBs = scaleMinGBs = scaleMaxGBs = scaleAvgGBs =
  addLocalGBs = addMinGBs = addMaxGBs = addAvgGBs =
  triadLocalGBs = triadMinGBs = triadMaxGBs = triadAvgGBs = 0.0;

  MPI_Comm_size( MPI_COMM_WORLD, &commSize );
  MPI_Comm_rank( MPI_COMM_WORLD, &myRank );

  rv = HPCC_Stream( params, 0 == myRank, MPI_COMM_WORLD, myRank,
      &copyLocalGBs, &scaleLocalGBs, &addLocalGBs, &triadLocalGBs, &failure );
  MPI_Reduce( &rv, &errCount, 1, MPI_INT, MPI_SUM, 0, MPI_COMM_WORLD );
  MPI_Allreduce( &failure, &failureAll, 1, MPI_INT, MPI_MAX, MPI_COMM_WORLD );
  if (failureAll) params->Failure = 1;

  MPI_Reduce( &copyLocalGBs, &copyMinGBs, 1, MPI_DOUBLE, MPI_MIN, 0, MPI_COMM_WORLD );
  MPI_Reduce( &copyLocalGBs, &copyAvgGBs, 1, MPI_DOUBLE, MPI_SUM, 0, MPI_COMM_WORLD );
  MPI_Reduce( &copyLocalGBs, &copyMaxGBs, 1, MPI_DOUBLE, MPI_MAX, 0, MPI_COMM_WORLD );
  copyAvgGBs /= commSize;
  MPI_Reduce( &scaleLocalGBs, &scaleMinGBs, 1, MPI_DOUBLE, MPI_MIN, 0, MPI_COMM_WORLD );
  MPI_Reduce( &scaleLocalGBs, &scaleAvgGBs, 1, MPI_DOUBLE, MPI_SUM, 0, MPI_COMM_WORLD );
  MPI_Reduce( &scaleLocalGBs, &scaleMaxGBs, 1, MPI_DOUBLE, MPI_MAX, 0, MPI_COMM_WORLD );
  scaleAvgGBs /= commSize;
  MPI_Reduce( &addLocalGBs, &addMinGBs, 1, MPI_DOUBLE, MPI_MIN, 0, MPI_COMM_WORLD );
  MPI_Reduce( &addLocalGBs, &addAvgGBs, 1, MPI_DOUBLE, MPI_SUM, 0, MPI_COMM_WORLD );
  MPI_Reduce( &addLocalGBs, &addMaxGBs, 1, MPI_DOUBLE, MPI_MAX, 0, MPI_COMM_WORLD );
  addAvgGBs /= commSize;
  MPI_Reduce( &triadLocalGBs, &triadMinGBs, 1, MPI_DOUBLE, MPI_MIN, 0, MPI_COMM_WORLD );
  MPI_Reduce( &triadLocalGBs, &triadAvgGBs, 1, MPI_DOUBLE, MPI_SUM, 0, MPI_COMM_WORLD );
  MPI_Reduce( &triadLocalGBs, &triadMaxGBs, 1, MPI_DOUBLE, MPI_MAX, 0, MPI_COMM_WORLD );
  triadAvgGBs /= commSize;

  MPI_Bcast( &copyAvgGBs, 1, MPI_DOUBLE, 0, MPI_COMM_WORLD ); params->StarStreamCopyGBs = copyAvgGBs;
  MPI_Bcast( &scaleAvgGBs, 1, MPI_DOUBLE, 0, MPI_COMM_WORLD ); params->StarStreamScaleGBs = scaleAvgGBs;
  MPI_Bcast( &addAvgGBs, 1, MPI_DOUBLE, 0, MPI_COMM_WORLD ); params->StarStreamAddGBs = addAvgGBs;
  MPI_Bcast( &triadAvgGBs, 1, MPI_DOUBLE, 0, MPI_COMM_WORLD ); params->StarStreamTriadGBs = triadAvgGBs;

  BEGIN_IO( myRank, params->outFname, outputFile);
  fprintf( outputFile, "Node(s) with error %d\n", errCount );
  fprintf( outputFile, "Minimum Copy GB/s %.6f\n", copyMinGBs );
  fprintf( outputFile, "Average Copy GB/s %.6f\n", copyAvgGBs );
  fprintf( outputFile, "Maximum Copy GB/s %.6f\n", copyMaxGBs );
  fprintf( outputFile, "Minimum Scale GB/s %.6f\n", scaleMinGBs );
  fprintf( outputFile, "Average Scale GB/s %.6f\n", scaleAvgGBs );
  fprintf( outputFile, "Maximum Scale GB/s %.6f\n", scaleMaxGBs );
  fprintf( outputFile, "Minimum Add GB/s %.6f\n", addMinGBs );
  fprintf( outputFile, "Average Add GB/s %.6f\n", addAvgGBs );
  fprintf( outputFile, "Maximum Add GB/s %.6f\n", addMaxGBs );
  fprintf( outputFile, "Minimum Triad GB/s %.6f\n", triadMinGBs );
  fprintf( outputFile, "Average Triad GB/s %.6f\n", triadAvgGBs );
  fprintf( outputFile, "Maximum Triad GB/s %.6f\n", triadMaxGBs );
  END_IO( myRank, outputFile );

  return 0;
}

int
HPCC_SingleStream(HPCC_Params *params) {
  int myRank, commSize;
  int rv, errCount, rank, failure = 0;
  double copyLocalGBs, scaleLocalGBs, addLocalGBs, triadLocalGBs;
  double scl = 1.0 / RAND_MAX;
  FILE *outputFile;

  copyLocalGBs = scaleLocalGBs = addLocalGBs = triadLocalGBs = 0.0;

  MPI_Comm_size( MPI_COMM_WORLD, &commSize );
  MPI_Comm_rank( MPI_COMM_WORLD, &myRank );

  srand(time(NULL));
  scl *= commSize;

  /* select a node at random, but not node 0 (unless there is just one node) */
  if (1 == commSize)
    rank = 0;
  else
    for (rank = 0; ; rank = (int)(scl * rand())) {
      if (rank > 0 && rank < commSize) break;
    }

  MPI_Bcast( &rank, 1, MPI_INT, 0, MPI_COMM_WORLD ); /* broadcast the rank selected on node 0 */

  if (myRank == rank) /* if this node has been selected */
    rv = HPCC_Stream( params, 0 == myRank, MPI_COMM_SELF, myRank,
        &copyLocalGBs, &scaleLocalGBs, &addLocalGBs, &triadLocalGBs, &failure );

  MPI_Bcast( &rv, 1, MPI_INT, rank, MPI_COMM_WORLD ); /* broadcast error code */
  MPI_Bcast( &failure, 1, MPI_INT, rank, MPI_COMM_WORLD ); /* broadcast failure indication */
  if (failure) params->Failure = 1;

  /* broadcast results */
  MPI_Bcast( &copyLocalGBs,  1, MPI_DOUBLE, rank, MPI_COMM_WORLD );
  MPI_Bcast( &scaleLocalGBs, 1, MPI_DOUBLE, rank, MPI_COMM_WORLD );
  MPI_Bcast( &addLocalGBs,   1, MPI_DOUBLE, rank, MPI_COMM_WORLD );
  MPI_Bcast( &triadLocalGBs, 1, MPI_DOUBLE, rank, MPI_COMM_WORLD );
  errCount = rv;
  params->SingleStreamCopyGBs = copyLocalGBs;
  params->SingleStreamScaleGBs = scaleLocalGBs;
  params->SingleStreamAddGBs = addLocalGBs;
  params->SingleStreamTriadGBs = triadLocalGBs;

  BEGIN_IO( myRank, params->outFname, outputFile);
  fprintf( outputFile, "Node(s) with error %d\n", errCount );
  fprintf( outputFile, "Node selected %d\n", rank );
  fprintf( outputFile, "Single STREAM Copy GB/s %.6f\n", copyLocalGBs );
  fprintf( outputFile, "Single STREAM Scale GB/s %.6f\n", scaleLocalGBs );
  fprintf( outputFile, "Single STREAM Add GB/s %.6f\n", addLocalGBs );
  fprintf( outputFile, "Single STREAM Triad GB/s %.6f\n", triadLocalGBs );
  END_IO( myRank, outputFile );

  return 0;
}