UCFD_SPARSE  v1.0
Documentation
Loading...
Searching...
No Matches
mpi3d.py
Go to the documentation of this file.
1from ctypes import CDLL, c_double, c_void_p, c_int
2from mpi4py import MPI
3import numpy as np
4import sys
5
6from utils.pbutils import make_nei_ele, read_input
7
8"""
9doxygen?
10"""
11def run(params, comm, lib):
12 # Constants
13 nvars, nface, ndims = 5, 6, 3
14 gamma = 1.4
15 nx, ny, nz, npx, npy, npz, nstep = params
16 m, n, l = nx//npx, ny//npy, nz//npz
17 nele = m*n*l
18 rank = comm.Get_rank()
19
20 # Ctypes functions
21 update = lib.serial_update
22 pre_lusgs = lib.serial_pre_lusgs
23 lower_sweep = lib.ns_serial_lower_sweep
24 upper_sweep = lib.ns_serial_upper_sweep
25
26 update.argtypes = [c_int, c_int, c_void_p, c_void_p]
27 pre_lusgs.argtypes = [c_int, c_int, c_double, \
28 c_void_p, c_void_p, c_void_p, c_void_p]
29 lower_sweep.argtypes = [c_int, c_int, c_int, c_int, \
30 c_void_p, c_void_p, c_void_p, c_void_p, c_void_p, \
31 c_void_p, c_void_p, c_void_p, c_void_p, c_void_p]
32 upper_sweep.argtypes = [c_int, c_int, c_int, c_int, \
33 c_void_p, c_void_p, c_void_p, c_void_p, c_void_p, \
34 c_void_p, c_void_p, c_void_p, c_void_p, c_void_p]
35
36 if rank == 0:
37 print("[Run] Allocating arrays...")
38
39 # Array allocation
40 nei_ele = np.empty((nface, nele), dtype=np.int32)
41 mapping = np.arange(nele, dtype=np.int32)
42 unmapping = np.arange(nele, dtype=np.int32)
43
44 # Flow variables & arrays
45 rho, u, v, w, p = gamma, 1.5, 0, 0, 1
46 et = p / (gamma - 1) + 0.5*rho*u**2
47 upts = np.empty((nvars, nele), dtype=np.float64)
48 rhs = np.empty_like(upts)
49 dub = np.empty_like(upts)
50 diag = np.ones((nele,), dtype=np.float64)
51 fspr = np.ones((6, nele), dtype=np.float64)
52 dt = np.ones((nele,), dtype=np.float64)*0.1
53 fnorm_vol = np.ones((nface, nele), dtype=np.float64)
54 vfi = np.array([[0, -1, 0], [-1, 0, 0], [0, 0, 1], [1, 0, 0], [0, 0, -1], [0, 1, 0]], dtype=np.float64)
55 vec_fnorm = np.repeat(vfi[:, :, np.newaxis], nele, axis=2)
56
57 # Initial condition
58 if rank == 0:
59 print("[Run] Initializing arrays...")
60 upts[0] = rhs[0] = dub[0] = rho
61 upts[1] = rhs[1] = dub[1] = rho*u
62 upts[2] = rhs[2] = dub[2] = rho*v
63 upts[3] = rhs[3] = dub[3] = rho*w
64 upts[4] = rhs[4] = dub[4] = et
65
66 make_nei_ele(m, n, l, nei_ele)
67
68 # Iteration
69 if rank == 0:
70 print("[Run] Starting iteration...")
71 tarr = np.zeros(nstep, dtype=np.float64)
72 comm.Barrier()
73 for i in range(nstep):
74 start = MPI.Wtime()
75 pre_lusgs(nele, nface, 1.0, fnorm_vol.ctypes.data, dt.ctypes.data, \
76 diag.ctypes.data, fspr.ctypes.data)
77 lower_sweep(nele, nvars, nface, ndims, nei_ele.ctypes.data, mapping.ctypes.data, \
78 unmapping.ctypes.data, fnorm_vol.ctypes.data, vec_fnorm.ctypes.data, \
79 upts.ctypes.data, rhs.ctypes.data, dub.ctypes.data, \
80 diag.ctypes.data, fspr.ctypes.data)
81 upper_sweep(nele, nvars, nface, ndims, nei_ele.ctypes.data, mapping.ctypes.data, \
82 unmapping.ctypes.data, fnorm_vol.ctypes.data, vec_fnorm.ctypes.data, \
83 upts.ctypes.data, rhs.ctypes.data, dub.ctypes.data, \
84 diag.ctypes.data, fspr.ctypes.data)
85 update(nele, nvars, upts.ctypes.data, rhs.ctypes.data)
86 comm.Barrier()
87 tarr[i] = MPI.Wtime() - start
88
89 if nstep > 100:
90 return tarr[15:].mean()*1000
91 else:
92 return tarr.mean()*1000
93
94
95def write_output(nprocs, params, times):
96 m, n, l, dm, dn, dl, nstep = params
97 neles = m*n*l
98
99 f = open("MPIOUTPUT_py_{}_{}.txt".format(nprocs, neles), 'w')
100 f.write("========= LU-SGS example output =========\n\n")
101 f.write("*** Problem setup ***\n")
102 f.write("Number of cores: {}\n".format(nprocs))
103 f.write("Number of iteration step: {}\n".format(nstep))
104 f.write("Number of elements: {} = {} x {} x {}\n".format(neles, m, n, l))
105 f.write("Partitioning with {}: {} x {} x {}\n".format(dm*dn*dl, dm, dn, dl))
106 f.write("Elements per core: {} = {} x {} x {}\n\n".format(neles//nprocs, m//dm, n//dn, l//dl))
107 f.write("*** Average runtime for each processor [ms] ***\n")
108 for i in range(nprocs):
109 f.write("{:.6f}\n".format(times[i]))
110 f.write("\n*** Average runtime for entire processors [ms] ***\n")
111 f.write("{:.6f}\n".format(sum(times)/nprocs))
112 f.close()
113
114
115if __name__ == "__main__":
116 # MPI initialize
117 comm = MPI.COMM_WORLD
118 rank = comm.Get_rank()
119 nproc = comm.Get_size()
120
121 if rank == 0:
122 if (len(sys.argv)) < 3:
123 print("[Error] Input file, dynamic library argument is needed")
124 err = 1
125 else:
126 inputname = sys.argv[1]
127 params = read_input(inputname)
128 if params is None:
129 err = 1
130 else:
131 dm = params[3]
132 dn = params[4]
133 dl = params[5]
134
135 if nproc != dm*dn*dl:
136 print("[Error] Partitioning number mismatch")
137 print("{} processors, but total {} partitioning".format(nproc, dm*dn*dl))
138 err = 1
139 else:
140 print("[Main] 3D hexahedral example starts...")
141 err = 0
142 else:
143 err = 0
144 params = np.empty(7, dtype=np.int64)
145
146 err = comm.bcast(err, root=0)
147 if err == 1:
148 sys.exit()
149
150 params = comm.bcast(params, root=0)
151
152 # Dynamic Linking
153 try:
154 lib = CDLL(sys.argv[2])
155 except FileNotFoundError:
156 print("[Error] Dynamic library file not found")
157 sys.exit()
158
159 avtimei = run(params, comm, lib)
160 tbuf = comm.gather(avtimei)
161
162 if rank == 0:
163 print("[Main] Run finished. Writing result...")
164 write_output(nproc, params, tbuf)
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
def run(params, comm, lib)
Definition: mpi3d.py:11
def write_output(nprocs, params, times)
Definition: mpi3d.py:95