xref: /petsc/src/benchmarks/benchmarkExample.py (revision 9895aa37ac365bac650f6bd8bf977519f7222510)
1#!/usr/bin/env python
2import os,sys
3sys.path.append(os.path.join(os.environ['PETSC_DIR'], 'config'))
4from builder2 import buildExample
5from benchmarkBatch import generateBatchScript
6
7class PETSc(object):
8  def __init__(self):
9    return
10
11  def dir(self):
12    '''Return the root directory for the PETSc tree (usually $PETSC_DIR)'''
13    # This should search for a valid PETSc
14    return os.environ['PETSC_DIR']
15
16  def arch(self):
17    '''Return the PETSc build label (usually $PETSC_ARCH)'''
18    # This should be configurable
19    return os.environ['PETSC_ARCH']
20
21  def mpiexec(self):
22    '''Return the path for the mpi launch executable'''
23    mpiexec = os.path.join(self.dir(), self.arch(), 'bin', 'mpiexec')
24    if not os.path.isfile(mpiexec):
25      return None
26    return mpiexec
27
28  def example(self, num):
29    '''Return the path to the executable for a given example number'''
30    return os.path.join(self.dir(), self.arch(), 'lib', 'ex'+str(num)+'-obj', 'ex'+str(num))
31
32  def source(self, library, num):
33    '''Return the path to the sources for a given example number'''
34    d = os.path.join(self.dir(), 'src', library.lower(), 'examples', 'tutorials')
35    name = 'ex'+str(num)
36    sources = []
37    for f in os.listdir(d):
38      if f == name+'.c':
39        sources.insert(0, f)
40      elif f.startswith(name) and f.endswith('.cu'):
41        sources.append(f)
42    return map(lambda f: os.path.join(d, f), sources)
43
44class PETScExample(object):
45  def __init__(self, library, num, **defaultOptions):
46    self.petsc   = PETSc()
47    self.library = library
48    self.num     = num
49    self.opts    = defaultOptions
50    return
51
52  @staticmethod
53  def runShellCommand(command, cwd = None, log = True):
54    import subprocess
55
56    Popen = subprocess.Popen
57    PIPE  = subprocess.PIPE
58    if log: print 'Executing: %s\n' % (command,)
59    pipe = Popen(command, cwd=cwd, stdin=None, stdout=PIPE, stderr=PIPE, bufsize=-1, shell=True, universal_newlines=True)
60    (out, err) = pipe.communicate()
61    ret = pipe.returncode
62    return (out, err, ret)
63
64  def optionsToString(self, **opts):
65    '''Convert a dictionary of options to a command line argument string'''
66    a = []
67    for key,value in opts.iteritems():
68      if value is None:
69        a.append('-'+key)
70      else:
71        a.append('-'+key+' '+str(value))
72    return ' '.join(a)
73
74  def run(self, numProcs = 1, log = True, **opts):
75    if self.petsc.mpiexec() is None:
76      cmd = self.petsc.example(self.num)
77    else:
78      cmd = ' '.join([self.petsc.mpiexec(), '-n', str(numProcs), self.petsc.example(self.num)])
79    cmd += ' '+self.optionsToString(**self.opts)+' '+self.optionsToString(**opts)
80    if 'batch' in opts and opts['batch']:
81      del opts['batch']
82      filename = generateBatchScript(self.num, numProcs, 120, ' '+self.optionsToString(**self.opts)+' '+self.optionsToString(**opts))
83      # Submit job
84      out, err, ret = self.runShellCommand('qsub -q gpu '+filename, log = log)
85      if ret:
86        print err
87        print out
88    else:
89      out, err, ret = self.runShellCommand(cmd, log = log)
90      if ret:
91        print err
92        print out
93    return out
94
95def processSummary(moduleName, defaultStage, eventNames, times, events):
96  '''Process the Python log summary into plot data'''
97  m = __import__(moduleName)
98  reload(m)
99  # Total Time
100  times.append(m.Time[0])
101  # Particular events
102  for name in eventNames:
103    if name.find(':') >= 0:
104      stageName, name = name.split(':', 1)
105      stage = getattr(m, stageName)
106    else:
107      stage = getattr(m, defaultStage)
108    if name in stage.event:
109      if not name in events:
110        events[name] = []
111      try:
112        events[name].append((stage.event[name].Time[0], stage.event[name].Flops[0]/(stage.event[name].Time[0] * 1e6)))
113      except ZeroDivisionError:
114        events[name].append((stage.event[name].Time[0], 0))
115  return
116
117def plotTime(library, num, eventNames, sizes, times, events):
118  from pylab import legend, plot, show, title, xlabel, ylabel
119  import numpy as np
120
121  arches = sizes.keys()
122  data   = []
123  for arch in arches:
124    data.append(sizes[arch])
125    data.append(times[arch])
126  plot(*data)
127  title('Performance on '+library+' Example '+str(num))
128  xlabel('Number of Dof')
129  ylabel('Time (s)')
130  legend(arches, 'upper left', shadow = True)
131  show()
132  return
133
134def plotEventTime(library, num, eventNames, sizes, times, events, filename = None):
135  from pylab import close, legend, plot, savefig, show, title, xlabel, ylabel
136  import numpy as np
137
138  close()
139  arches = sizes.keys()
140  bs     = events[arches[0]].keys()[0]
141  data   = []
142  names  = []
143  for event, color in zip(eventNames, ['b', 'g', 'r', 'y']):
144    for arch, style in zip(arches, ['-', ':']):
145      if event in events[arch][bs]:
146        names.append(arch+'-'+str(bs)+' '+event)
147        data.append(sizes[arch][bs])
148        data.append(np.array(events[arch][bs][event])[:,0])
149        data.append(color+style)
150      else:
151        print 'Could not find %s in %s-%d events' % (event, arch, bs)
152  print data
153  plot(*data)
154  title('Performance on '+library+' Example '+str(num))
155  xlabel('Number of Dof')
156  ylabel('Time (s)')
157  legend(names, 'upper left', shadow = True)
158  if filename is None:
159    show()
160  else:
161    savefig(filename)
162  return
163
164def plotEventFlop(library, num, eventNames, sizes, times, events, filename = None):
165  from pylab import legend, plot, savefig, semilogy, show, title, xlabel, ylabel
166  import numpy as np
167
168  arches = sizes.keys()
169  bs     = events[arches[0]].keys()[0]
170  data   = []
171  names  = []
172  for event, color in zip(eventNames, ['b', 'g', 'r', 'y']):
173    for arch, style in zip(arches, ['-', ':']):
174      if event in events[arch][bs]:
175        names.append(arch+'-'+str(bs)+' '+event)
176        data.append(sizes[arch][bs])
177        data.append(1e-3*np.array(events[arch][bs][event])[:,1])
178        data.append(color+style)
179      else:
180        print 'Could not find %s in %s-%d events' % (event, arch, bs)
181  semilogy(*data)
182  title('Performance on '+library+' Example '+str(num))
183  xlabel('Number of Dof')
184  ylabel('Computation Rate (GF/s)')
185  legend(names, 'upper left', shadow = True)
186  if filename is None:
187    show()
188  else:
189    savefig(filename)
190  return
191
192def plotSummaryLine(library, num, eventNames, sizes, times, events):
193  from pylab import legend, plot, show, title, xlabel, ylabel
194  import numpy as np
195  showTime       = False
196  showEventTime  = True
197  showEventFlops = True
198  arches         = sizes.keys()
199  # Time
200  if showTime:
201    data = []
202    for arch in arches:
203      data.append(sizes[arch])
204      data.append(times[arch])
205    plot(*data)
206    title('Performance on '+library+' Example '+str(num))
207    xlabel('Number of Dof')
208    ylabel('Time (s)')
209    legend(arches, 'upper left', shadow = True)
210    show()
211  # Common event time
212  #   We could make a stacked plot like Rio uses here
213  if showEventTime:
214    bs    = events[arches[0]].keys()[0]
215    data  = []
216    names = []
217    for event, color in zip(eventNames, ['b', 'g', 'r', 'y']):
218      for arch, style in zip(arches, ['-', ':']):
219        if event in events[arch][bs]:
220          names.append(arch+'-'+str(bs)+' '+event)
221          data.append(sizes[arch][bs])
222          data.append(np.array(events[arch][bs][event])[:,0])
223          data.append(color+style)
224        else:
225          print 'Could not find %s in %s-%d events' % (event, arch, bs)
226    print data
227    plot(*data)
228    title('Performance on '+library+' Example '+str(num))
229    xlabel('Number of Dof')
230    ylabel('Time (s)')
231    legend(names, 'upper left', shadow = True)
232    show()
233  # Common event flops
234  #   We could make a stacked plot like Rio uses here
235  if showEventFlops:
236    bs    = events[arches[0]].keys()[0]
237    data  = []
238    names = []
239    for event, color in zip(eventNames, ['b', 'g', 'r', 'y']):
240      for arch, style in zip(arches, ['-', ':']):
241        if event in events[arch][bs]:
242          names.append(arch+'-'+str(bs)+' '+event)
243          data.append(sizes[arch][bs])
244          data.append(np.array(events[arch][bs][event])[:,1])
245          data.append(color+style)
246        else:
247          print 'Could not find %s in %s-%d events' % (event, arch, bs)
248    plot(*data)
249    title('Performance on '+library+' Example '+str(num))
250    xlabel('Number of Dof')
251    ylabel('Computation Rate (MF/s)')
252    legend(names, 'upper left', shadow = True)
253    show()
254  return
255
256def plotSummaryBar(library, num, eventNames, sizes, times, events):
257  import numpy as np
258  import matplotlib.pyplot as plt
259
260  eventColors = ['b', 'g', 'r', 'y']
261  arches = sizes.keys()
262  names  = []
263  N      = len(sizes[arches[0]])
264  width  = 0.2
265  ind    = np.arange(N) - 0.25
266  bars   = {}
267  for arch in arches:
268    bars[arch] = []
269    bottom = np.zeros(N)
270    for event, color in zip(eventNames, eventColors):
271      names.append(arch+' '+event)
272      times = np.array(events[arch][event])[:,0]
273      bars[arch].append(plt.bar(ind, times, width, color=color, bottom=bottom))
274      bottom += times
275    ind += 0.3
276
277  plt.xlabel('Number of Dof')
278  plt.ylabel('Time (s)')
279  plt.title('GPU vs. CPU Performance on '+library+' Example '+str(num))
280  plt.xticks(np.arange(N), map(str, sizes[arches[0]]))
281  #plt.yticks(np.arange(0,81,10))
282  #plt.legend( (p1[0], p2[0]), ('Men', 'Women') )
283  plt.legend([bar[0] for bar in bars[arches[0]]], eventNames, 'upper right', shadow = True)
284
285  plt.show()
286  return
287
288def getDMComplexSize(dim, out):
289  '''Retrieves the number of cells from '''
290  size = 0
291  for line in out.split('\n'):
292    if line.strip().startswith(str(dim)+'-cells: '):
293      size = int(line.strip()[9:])
294      break
295  return size
296
297def run_DMDA(ex, name, opts, args, sizes, times, events, log=True):
298  for n in map(int, args.size):
299    ex.run(log=log, da_grid_x=n, da_grid_y=n, **opts)
300    sizes[name].append(n*n * args.comp)
301    processSummary('summary', args.stage, args.events, times[name], events[name])
302  return
303
304def run_DMComplex(ex, name, opts, args, sizes, times, events, log=True):
305  # This should eventually be replaced by a direct FFC/Ignition interface
306  if args.operator == 'laplacian':
307    numComp  = 1
308  elif args.operator == 'elasticity':
309    numComp  = args.dim
310  else:
311    raise RuntimeError('Unknown operator: %s' % args.operator)
312
313  for numBlock in [2**i for i in map(int, args.blockExp)]:
314    opts['gpu_blocks'] = numBlock
315    # Generate new block size
316    cmd = './bin/pythonscripts/PetscGenerateFEMQuadrature.py %d %d %d %d %s %s.h' % (args.dim, args.order, numComp, numBlock, args.operator, os.path.splitext(source[0])[0])
317    print(cmd)
318    ret = os.system('python '+cmd)
319    args.files = ['['+','.join(source)+']']
320    buildExample(args)
321    sizes[name][numBlock]  = []
322    times[name][numBlock]  = []
323    events[name][numBlock] = {}
324    for r in map(float, args.refine):
325      out = ex.run(log=log, refinement_limit=r, **opts)
326      sizes[name][numBlock].append(getDMComplexSize(args.dim, out))
327      processSummary('summary', args.stage, args.events, times[name][numBlock], events[name][numBlock])
328  return
329
330def outputData(sizes, times, events, name = 'output.py'):
331  if os.path.exists(name):
332    base, ext = os.path.splitext(name)
333    num = 1
334    while os.path.exists(base+str(num)+ext):
335      num += 1
336    name = base+str(num)+ext
337  with file(name, 'w') as f:
338    f.write('#PETSC_ARCH='+os.environ['PETSC_ARCH']+' '+' '.join(sys.argv)+'\n')
339    f.write('sizes  = '+repr(sizes)+'\n')
340    f.write('times  = '+repr(times)+'\n')
341    f.write('events = '+repr(events)+'\n')
342  return
343
344if __name__ == '__main__':
345  import argparse
346
347  parser = argparse.ArgumentParser(description     = 'PETSc Benchmarking',
348                                   epilog          = 'This script runs src/<library>/examples/tutorials/ex<num>, For more information, visit http://www.mcs.anl.gov/petsc',
349                                   formatter_class = argparse.ArgumentDefaultsHelpFormatter)
350  parser.add_argument('--library', default='SNES',                     help='The PETSc library used in this example')
351  parser.add_argument('--num',     type = int, default='5',            help='The example number')
352  parser.add_argument('--module',  default='summary',                  help='The module for timing output')
353  parser.add_argument('--stage',   default='Main_Stage',               help='The default logging stage')
354  parser.add_argument('--events',  nargs='+',                          help='Events to process')
355  parser.add_argument('--batch',   action='store_true', default=False, help='Generate batch files for the runs instead')
356  parser.add_argument('--daemon',  action='store_true', default=False, help='Run as a daemon')
357  subparsers = parser.add_subparsers(help='DM types')
358
359  parser_dmda = subparsers.add_parser('DMDA', help='Use a DMDA for the problem geometry')
360  parser_dmda.add_argument('--size', nargs='+',  default=['10'], help='Grid size (implementation dependent)')
361  parser_dmda.add_argument('--comp', type = int, default='1',    help='Number of field components')
362  parser_dmda.add_argument('runs',   nargs='*',                  help='Run descriptions: <name>=<args>')
363
364  parser_dmmesh = subparsers.add_parser('DMComplex', help='Use a DMComplex for the problem geometry')
365  parser_dmmesh.add_argument('--dim',      type = int, default='2',        help='Spatial dimension')
366  parser_dmmesh.add_argument('--refine',   nargs='+',  default=['0.0'],    help='List of refinement limits')
367  parser_dmmesh.add_argument('--order',    type = int, default='1',        help='Order of the finite element')
368  parser_dmmesh.add_argument('--operator', default='laplacian',            help='The operator name')
369  parser_dmmesh.add_argument('--blockExp', nargs='+', default=range(0, 5), help='List of block exponents j, block size is 2^j')
370  parser_dmmesh.add_argument('runs',       nargs='*',                      help='Run descriptions: <name>=<args>')
371
372  args = parser.parse_args()
373  print(args)
374  if hasattr(args, 'comp'):
375    args.dmType = 'DMDA'
376  else:
377    args.dmType = 'DMComplex'
378
379  ex     = PETScExample(args.library, args.num, log_summary='summary.dat', log_summary_python = None if args.batch else args.module+'.py', preload='off')
380  source = ex.petsc.source(args.library, args.num)
381  sizes  = {}
382  times  = {}
383  events = {}
384  log    = not args.daemon
385
386  if args.daemon:
387    import daemon
388    print 'Starting daemon'
389    daemon.createDaemon('.')
390
391  for run in args.runs:
392    name, stropts = run.split('=', 1)
393    opts = dict([t if len(t) == 2 else (t[0], None) for t in [arg.split('=', 1) for arg in stropts.split(' ')]])
394    if args.dmType == 'DMDA':
395      sizes[name]  = []
396      times[name]  = []
397      events[name] = {}
398      run_DMDA(ex, name, opts, args, sizes, times, events, log=log)
399    elif args.dmType == 'DMComplex':
400      sizes[name]  = {}
401      times[name]  = {}
402      events[name] = {}
403      run_DMComplex(ex, name, opts, args, sizes, times, events, log=log)
404  outputData(sizes, times, events)
405  if not args.batch and log: plotSummaryLine(args.library, args.num, args.events, sizes, times, events)
406# Benchmark for ex50
407# ./src/benchmarks/benchmarkExample.py --events VecMDot VecMAXPY KSPGMRESOrthog MatMult VecCUSPCopyTo VecCUSPCopyFrom MatCUSPCopyTo --num 50 DMDA --size 10 20 50 100 --comp 4 CPU='pc_type=none mat_no_inode dm_vec_type=seq dm_mat_type=seqaij' GPU='pc_type=none mat_no_inode dm_vec_type=seqcusp dm_mat_type=seqaijcusp cusp_synchronize'
408# Benchmark for ex52
409# ./src/benchmarks/benchmarkExample.py --events IntegBatchCPU IntegBatchGPU IntegGPUOnly --num 52 DMComplex --refine 0.0625 0.00625 0.000625 0.0000625 --blockExp 4 --order=1 CPU='dm_view show_residual=0 compute_function batch' GPU='dm_view show_residual=0 compute_function batch gpu gpu_batches=8'
410# ./src/benchmarks/benchmarkExample.py --events IntegBatchCPU IntegBatchGPU IntegGPUOnly --num 52 DMComplex --refine 0.0625 0.00625 0.000625 0.0000625 --blockExp 4 --order=1 --operator=elasticity CPU='dm_view op_type=elasticity show_residual=0 compute_function batch' GPU='dm_view op_type=elasticity show_residual=0 compute_function batch gpu gpu_batches=8'
411# ./src/benchmarks/benchmarkExample.py --events IntegBatchCPU IntegBatchGPU IntegGPUOnly --num 52 DMComplex --dim=3 --refine 0.0625 0.00625 0.000625 0.0000625 --blockExp 4 --order=1 CPU='dim=3 dm_view show_residual=0 compute_function batch' GPU='dim=3 dm_view show_residual=0 compute_function batch gpu gpu_batches=8'
412