13.9.7. Macro “reoptimizeZoneBiFunMpi.py”
13.9.7.1. Objective
The objective of the macro is to give another example of a two level parallelism program using MPI paradigm. In the former example MPI function call is implicit using Uranie facilities. In this one, explicit calls to MPI functions is done. It’s presented to illustrate the case when the user evaluation fonction is an MPI function.
It takes the former example of zoning problem and adapts it. The intermediate level does not use a
TLauncher2 to run all different arrangements, but encodes it. Each Mpi ressources evaluates
different possible arrangements keeping its best one, and Mpi reduce these results to the final
result.
Warning
To be run, the code need the mpi4py python module to be installed. This dependancy is not tested when Uranie is installed
13.9.7.2. Macro Uranie
The low level evaluation function is the same than in previous example and is not shown again.
#!/usr/bin/env python
import ROOT
import URANIE.DataServer as DataServer
import URANIE.Relauncher as Relauncher
import URANIE.Reoptimizer as Reoptimizer
import URANIE.MpiRelauncher as MpiRelauncher
from mpi4py import MPI
import ctypes
from reoptimizeZoneCore import lowfun
# mpi4py helpers
def getCalculMpiComm() :
""" convert C MPI_Comm to python """
typemap = {
ctypes.sizeof(ctypes.c_int) : ctypes.c_int, #mpich
ctypes.sizeof(ctypes.c_void_p) : ctypes.c_void_p, #openmpi
}
comm = MPI.Comm()
handle_t = typemap[MPI._sizeof(comm)]
handle_new = handle_t.from_address(MPI._addressof(comm))
handle_new.value = MpiRelauncher.TBiMpiRun.getCalculMpiCommPy()
return comm
# doe
doe = [
[0, 1, 2], [1, 2, 3], [2, 3, 4], [3, 4, 6], [4, 6, 7], [0, 6, 7], [0, 1, 7], # 7
[0, 1, 3], [1, 2, 4], [2, 3, 6], [3, 4, 7], [0, 4, 6], [1, 6, 7], [0, 2, 7], # 14
[0, 1, 4], [1, 2, 6], [2, 3, 7], [0, 3, 4], [1, 4, 6], [2, 6, 7], [0, 3, 7], # 21
[0, 1, 6], [1, 2, 7], [0, 2, 3], [1, 3, 4], [2, 4, 6], [3, 6, 7], [0, 4, 7], # 28
[0, 2, 4], [1, 3, 6], [2, 4, 7], [0, 3, 6], [1, 4, 7], [0, 2, 6], [1, 3, 7], # 35
]
# doe fun
def doefun(emz0, emz1, emz2, emz3) :
""" find the best zoning for a given emission values """
# mpi context
comm = getCalculMpiComm()
rank = comm.Get_rank()
size = comm.Get_size()
# search local best
ind = rank
idz1 = doe[ind][0]
idz2 = doe[ind][1]
idz3 = doe[ind][2]
memo = lowfun(emz0, emz1, emz2, emz3, idz1, idz2, idz3)
for next in range(rank, len(doe), size) :
idz1 = doe[next][0]
idz2 = doe[next][1]
idz3 = doe[next][2]
curr = lowfun(emz0, emz1, emz2, emz3, idz1, idz2, idz3)
if curr[0] < memo[0] :
memo = curr
# global best
# where is min
ret = (memo[0], rank)
(val, fid) = comm.allreduce(ret, op=MPI.MINLOC)
# get all min value
rval = memo
if fid != 0 :
if rank == fid :
comm.send(memo, dest=0, tag=0)
elif rank == 0 :
rval = comm.recv(source=fid, tag=0)
return rval
# top level control
def coeurMpi():
""" optimization on emission values """
# inputs
nami = ["zone1", "zone2", "zone3", "zone4"]
funi = [DataServer.TAttribute(n, 0., 1.) for n in nami]
# outputs
namo = ["diff", "v0", "v1", "v2", "v3", "v4", "v5", "v6", "v7", "v8", "v9"]
funo = [DataServer.TAttribute(n) for n in namo]
# top fun
fun = Relauncher.TPythonEval(doefun)
for inp in funi:
fun.addInput(inp)
for out in funo:
fun.addOutput(out)
fun.setMpi()
# runner
# run = Relauncher.TSequentialRun(fun)
run = MpiRelauncher.TBiMpiRun(fun, 3, False)
run.startSlave()
if run.onMaster():
# data
data = DataServer.TDataServer("tdsvzr", "tds4optim")
fun.addAllInputs(data)
# optim
gen = Reoptimizer.TVizirGenetic()
gen.setSize(300, 200000, 100)
viz = Reoptimizer.TVizir2(data, run, gen)
viz.addObjective(funo[0])
viz.solverLoop()
# results
data.exportData("__coeurMpy__.dat")
# clean
run.stopSlave()
#ROOT.SetOwnership(run, True)
#ROOT.EnableThreadSafety()
coeurMpi()
The top level function (coeurMpi) does not change from the previous example. The only difference
is in the creation of the TBiMpiRun instances, which is done with an extra parameter set to False.
This tell the constructor to avoid calling MPI_Init function. With mpi4py, this call is already
done when the module is imported
The evaluation MPI fonction (doefun) is totaly different. It uses the function getCalculMpiComm
to get the MPI Communicator object (MPI.Comm) dedicated to the calcul ressources. With it,
different call to its method can be done : Get_rank and Get_size to get the context ; allreduce,
send and recv to communicate between calculs ressources.
Function getCalculMpiComm is an helper function use to convert the MPI communicator get from Uranie
MpiRelauncher.TBiMpiRun.getCalculMpiCommPy class method to the mpi4py equivalent object
You may run this script the same way as the precedent example, with the same constraint on the
ressource number. Also note that this version is really faster than the previous one, avoiding
creation and manipulation of many TDataServer.