File: demo_timing.py

package info (click to toggle)
dolfin 2019.2.0~git20201207.b495043-5
  • links: PTS, VCS
  • area: main
  • in suites: bullseye
  • size: 30,988 kB
  • sloc: xml: 104,040; cpp: 102,020; python: 24,139; makefile: 300; javascript: 226; sh: 185
file content (65 lines) | stat: -rw-r--r-- 2,198 bytes parent folder | download | duplicates (5)
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
# Copyright (C) 2015 Jan Blechta
#
# This file is part of DOLFIN.
#
# DOLFIN is free software: you can redistribute it and/or modify
# it under the terms of the GNU Lesser General Public License as published by
# the Free Software Foundation, either version 3 of the License, or
# (at your option) any later version.
#
# DOLFIN is distributed in the hope that it will be useful,
# but WITHOUT ANY WARRANTY; without even the implied warranty of
# MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
# GNU Lesser General Public License for more details.
#
# You should have received a copy of the GNU Lesser General Public License
# along with DOLFIN. If not, see <http://www.gnu.org/licenses/>.

from dolfin import *

# Let's solve some variational problem to get non-trivial timings
mesh = UnitSquareMesh(32, 32)
V = FunctionSpace(mesh, "Lagrange", 1)
bc = DirichletBC(V, 0.0, lambda x: near(x[0], 0.0) or near(x[0], 1.0))
u, v = TrialFunction(V), TestFunction(V)
f = Expression("10*exp(-(pow(x[0] - 0.5, 2) + pow(x[1] - 0.5, 2)) / 0.02)", degree=2)
g = Expression("sin(5*x[0])", degree=2)
a = inner(grad(u), grad(v))*dx
L = f*v*dx + g*v*ds
u = Function(V)
solve(a == L, u, bc)

# List timings; average across processes in parallel
list_timings(TimingClear.keep, [TimingType.wall, TimingType.system])

# Get Table object with timings
t = timings(TimingClear.keep,
            [TimingType.wall, TimingType.user, TimingType.system])

# Use different MPI reductions
t_sum = MPI.sum(MPI.comm_world, t)
t_min = MPI.min(MPI.comm_world, t)
t_max = MPI.max(MPI.comm_world, t)
t_avg = MPI.avg(MPI.comm_world, t)

# Print aggregate timings to screen
print('\n'+t_sum.str(True))
print('\n'+t_min.str(True))
print('\n'+t_max.str(True))
print('\n'+t_avg.str(True))

# Store to XML file on rank 0
if MPI.rank(MPI.comm_world) == 0:
    f = File(MPI.comm_self, "timings_aggregate.xml")
    f << t_sum
    f << t_min
    f << t_max
    f << t_avg

# Store timings of each rank separately
f = File(MPI.comm_self, "timings_rank_%d.xml"
         % MPI.rank(MPI.comm_world))
f << t

# Helper function for storing rank-wise average, min and max
dump_timings_to_xml("timings_avg_min_max.xml", TimingClear.clear)