3 ===============================
4 SMPI: Simulate MPI Applications
5 ===============================
9 <object id="TOC" data="graphical-toc.svg" width="100%" type="image/svg+xml"></object>
11 window.onload=function() { // Wait for the SVG to be loaded before changing it
12 var elem=document.querySelector("#TOC").contentDocument.getElementById("SMPIBox")
13 elem.style="opacity:0.93999999;fill:#ff0000;fill-opacity:0.1";
19 SMPI enables the study of MPI application by emulating them on top of
20 the SimGrid simulator. This is particularly interesting to study
21 existing MPI applications within the comfort of the simulator.
23 To get started with SMPI, you should head to `the SMPI tutorial
24 <usecase_smpi>`_. You may also want to read the `SMPI reference
25 article <https://hal.inria.fr/hal-01415484>`_ or these `introductory
26 slides <http://simgrid.org/tutorials/simgrid-smpi-101.pdf>`_. If you
27 are new to MPI, you should first take our online `SMPI CourseWare
28 <https://simgrid.github.io/SMPI_CourseWare/>`_. It consists in several
29 projects that progressively introduce the MPI concepts. It proposes to
30 use SimGrid and SMPI to run the experiments, but the learning
31 objectives are centered on MPI itself.
33 Our goal is to enable the study of **unmodified MPI applications**.
34 Some constructs and features are still missing, but we can probably
35 add them on demand. If you already used MPI before, SMPI should sound
36 very familiar to you: Use smpicc instead of mpicc, and smpirun instead
37 of mpirun. The main difference is that smpirun takes a :ref:`simulated
38 platform <platform>` as an extra parameter.
40 For **further scalability**, you may modify your code to speed up your
41 studies or save memory space. Maximal **simulation accuracy**
42 requires some specific care from you.
50 In this mode, your application is actually executed. Every computation
51 occurs for real while every communication is simulated. In addition,
52 the executions are automatically benchmarked so that their timings can
53 be applied within the simulator.
55 SMPI can also go offline by replaying a trace. :ref:`Trace replay
56 <SMPI_offline>` is usually ways faster than online simulation (because
57 the computation are skipped), but it can only applied to applications
58 with constant execution and communication patterns (for the exact same
65 If your application is in C, then simply use ``smpicc`` as a
66 compiler just like you use mpicc with other MPI implementations. This
67 script still calls your default compiler (gcc, clang, ...) and adds
68 the right compilation flags along the way. If your application is in
69 C++, Fortran 77 or Fortran 90, use respectively ``smpicxx``,
70 ``smpiff`` or ``smpif90``.
76 Use the ``smpirun`` script as follows:
80 smpirun -hostfile my_hostfile.txt -platform my_platform.xml ./program -blah
82 - ``my_hostfile.txt`` is a classical MPI hostfile (that is, this file
83 lists the machines on which the processes must be dispatched, one
85 - ``my_platform.xml`` is a classical SimGrid platform file. Of course,
86 the hosts of the hostfile must exist in the provided platform.
87 - ``./program`` is the MPI program to simulate, that you compiled with ``smpicc``
88 - ``-blah`` is a command-line parameter passed to this program.
90 ``smpirun`` accepts other parameters, such as ``-np`` if you don't
91 want to use all the hosts defined in the hostfile, ``-map`` to display
92 on which host each rank gets mapped of ``-trace`` to activate the
93 tracing during the simulation. You can get the full list by running
96 ...............................
97 Debugging your Code within SMPI
98 ...............................
100 If you want to explore the automatic platform and deployment files
101 that are generated by ``smpirun``, add ``-keep-temps`` to the command
104 You can also run your simulation within valgrind or gdb using the
105 following commands. Once in GDB, each MPI ranks will be represented as
106 a regular thread, and you can explore the state of each of them as
109 .. code-block:: shell
111 smpirun -wrapper valgrind ...other args...
112 smpirun -wrapper "gdb --args" --cfg=contexts/factory:thread ...other args...
116 ................................
117 Simulating Collective Operations
118 ................................
120 MPI collective operations are crucial to the performance of MPI
121 applications and must be carefully optimized according to many
122 parameters. Every existing implementation provides several algorithms
123 for each collective operation, and selects by default the best suited
124 one, depending on the sizes sent, the number of nodes, the
125 communicator, or the communication library being used. These
126 decisions are based on empirical results and theoretical complexity
127 estimation, and are very different between MPI implementations. In
128 most cases, the users can also manually tune the algorithm used for
129 each collective operation.
131 SMPI can simulate the behavior of several MPI implementations:
132 OpenMPI, MPICH, `STAR-MPI <http://star-mpi.sourceforge.net/>`_, and
133 MVAPICH2. For that, it provides 115 collective algorithms and several
134 selector algorithms, that were collected directly in the source code
135 of the targeted MPI implementations.
137 You can switch the automatic selector through the
138 ``smpi/coll-selector`` configuration item. Possible values:
140 - **ompi:** default selection logic of OpenMPI (version 3.1.2)
141 - **mpich**: default selection logic of MPICH (version 3.3b)
142 - **mvapich2**: selection logic of MVAPICH2 (version 1.9) tuned
143 on the Stampede cluster
144 - **impi**: preliminary version of an Intel MPI selector (version
145 4.1.3, also tuned for the Stampede cluster). Due the closed source
146 nature of Intel MPI, some of the algorithms described in the
147 documentation are not available, and are replaced by mvapich ones.
148 - **default**: legacy algorithms used in the earlier days of
149 SimGrid. Do not use for serious perform performance studies.
151 .. todo:: default should not even exist.
157 You can also pick the algorithm used for each collective with the
158 corresponding configuration item. For example, to use the pairwise
159 alltoall algorithm, one should add ``--cfg=smpi/alltoall:pair`` to the
160 line. This will override the selector (if any) for this algorithm. It
161 means that the selected algorithm will be used
163 .. Warning:: Some collective may require specific conditions to be
164 executed correctly (for instance having a communicator with a power
165 of two number of nodes only), which are currently not enforced by
166 Simgrid. Some crashes can be expected while trying these algorithms
167 with unusual sizes/parameters
172 Most of these are best described in `STAR-MPI <http://www.cs.arizona.edu/~dkl/research/papers/ics06.pdf>`_.
174 - default: naive one, by default
175 - ompi: use openmpi selector for the alltoall operations
176 - mpich: use mpich selector for the alltoall operations
177 - mvapich2: use mvapich2 selector for the alltoall operations
178 - impi: use intel mpi selector for the alltoall operations
179 - automatic (experimental): use an automatic self-benchmarking algorithm
180 - bruck: Described by Bruck et.al. in <a href="http://ieeexplore.ieee.org/xpl/articleDetails.jsp?arnumber=642949">this paper</a>
181 - 2dmesh: organizes the nodes as a two dimensional mesh, and perform allgather
183 - 3dmesh: adds a third dimension to the previous algorithm
184 - rdb: recursive doubling: extends the mesh to a nth dimension, each one
186 - pair: pairwise exchange, only works for power of 2 procs, size-1 steps,
187 each process sends and receives from the same process at each step
188 - pair_light_barrier: same, with small barriers between steps to avoid
190 - pair_mpi_barrier: same, with MPI_Barrier used
191 - pair_one_barrier: only one barrier at the beginning
192 - ring: size-1 steps, at each step a process send to process (n+i)%size, and receives from (n-i)%size
193 - ring_light_barrier: same, with small barriers between some phases to avoid contention
194 - ring_mpi_barrier: same, with MPI_Barrier used
195 - ring_one_barrier: only one barrier at the beginning
196 - basic_linear: posts all receives and all sends,
197 starts the communications, and waits for all communication to finish
198 - mvapich2_scatter_dest: isend/irecv with scattered destinations, posting only a few messages at the same time
202 - default: naive one, by default
203 - ompi: use openmpi selector for the alltoallv operations
204 - mpich: use mpich selector for the alltoallv operations
205 - mvapich2: use mvapich2 selector for the alltoallv operations
206 - impi: use intel mpi selector for the alltoallv operations
207 - automatic (experimental): use an automatic self-benchmarking algorithm
208 - bruck: same as alltoall
209 - pair: same as alltoall
210 - pair_light_barrier: same as alltoall
211 - pair_mpi_barrier: same as alltoall
212 - pair_one_barrier: same as alltoall
213 - ring: same as alltoall
214 - ring_light_barrier: same as alltoall
215 - ring_mpi_barrier: same as alltoall
216 - ring_one_barrier: same as alltoall
217 - ompi_basic_linear: same as alltoall
222 - default: naive one, by default
223 - ompi: use openmpi selector for the gather operations
224 - mpich: use mpich selector for the gather operations
225 - mvapich2: use mvapich2 selector for the gather operations
226 - impi: use intel mpi selector for the gather operations
227 - automatic (experimental): use an automatic self-benchmarking algorithm which will iterate over all implemented versions and output the best
228 - ompi_basic_linear: basic linear algorithm from openmpi, each process sends to the root
229 - ompi_binomial: binomial tree algorithm
230 - ompi_linear_sync: same as basic linear, but with a synchronization at the
231 beginning and message cut into two segments.
232 - mvapich2_two_level: SMP-aware version from MVAPICH. Gather first intra-node (defaults to mpich's gather), and then exchange with only one process/node. Use mvapich2 selector to change these to tuned algorithms for Stampede cluster.
237 - default: naive one, by default
238 - ompi: use openmpi selector for the barrier operations
239 - mpich: use mpich selector for the barrier operations
240 - mvapich2: use mvapich2 selector for the barrier operations
241 - impi: use intel mpi selector for the barrier operations
242 - automatic (experimental): use an automatic self-benchmarking algorithm
243 - ompi_basic_linear: all processes send to root
244 - ompi_two_procs: special case for two processes
245 - ompi_bruck: nsteps = sqrt(size), at each step, exchange data with rank-2^k and rank+2^k
246 - ompi_recursivedoubling: recursive doubling algorithm
247 - ompi_tree: recursive doubling type algorithm, with tree structure
248 - ompi_doublering: double ring algorithm
249 - mvapich2_pair: pairwise algorithm
250 - mpich_smp: barrier intra-node, then inter-node
255 - default: naive one, by default
256 - ompi: use openmpi selector for the scatter operations
257 - mpich: use mpich selector for the scatter operations
258 - mvapich2: use mvapich2 selector for the scatter operations
259 - impi: use intel mpi selector for the scatter operations
260 - automatic (experimental): use an automatic self-benchmarking algorithm
261 - ompi_basic_linear: basic linear scatter
262 - ompi_binomial: binomial tree scatter
263 - mvapich2_two_level_direct: SMP aware algorithm, with an intra-node stage (default set to mpich selector), and then a basic linear inter node stage. Use mvapich2 selector to change these to tuned algorithms for Stampede cluster.
264 - mvapich2_two_level_binomial: SMP aware algorithm, with an intra-node stage (default set to mpich selector), and then a binomial phase. Use mvapich2 selector to change these to tuned algorithms for Stampede cluster.
269 - default: naive one, by default
270 - ompi: use openmpi selector for the reduce operations
271 - mpich: use mpich selector for the reduce operations
272 - mvapich2: use mvapich2 selector for the reduce operations
273 - impi: use intel mpi selector for the reduce operations
274 - automatic (experimental): use an automatic self-benchmarking algorithm
275 - arrival_pattern_aware: root exchanges with the first process to arrive
276 - binomial: uses a binomial tree
277 - flat_tree: uses a flat tree
278 - NTSL: Non-topology-specific pipelined linear-bcast function
279 0->1, 1->2 ,2->3, ....., ->last node: in a pipeline fashion, with segments
281 - scatter_gather: scatter then gather
282 - ompi_chain: openmpi reduce algorithms are built on the same basis, but the
283 topology is generated differently for each flavor
284 chain = chain with spacing of size/2, and segment size of 64KB
285 - ompi_pipeline: same with pipeline (chain with spacing of 1), segment size
286 depends on the communicator size and the message size
287 - ompi_binary: same with binary tree, segment size of 32KB
288 - ompi_in_order_binary: same with binary tree, enforcing order on the
290 - ompi_binomial: same with binomial algo (redundant with default binomial
292 - ompi_basic_linear: basic algorithm, each process sends to root
293 - mvapich2_knomial: k-nomial algorithm. Default factor is 4 (mvapich2 selector adapts it through tuning)
294 - mvapich2_two_level: SMP-aware reduce, with default set to mpich both for intra and inter communicators. Use mvapich2 selector to change these to tuned algorithms for Stampede cluster.
295 - rab: `Rabenseifner <https://fs.hlrs.de/projects/par/mpi//myreduce.html>`_'s reduce algorithm
300 - default: naive one, by default
301 - ompi: use openmpi selector for the allreduce operations
302 - mpich: use mpich selector for the allreduce operations
303 - mvapich2: use mvapich2 selector for the allreduce operations
304 - impi: use intel mpi selector for the allreduce operations
305 - automatic (experimental): use an automatic self-benchmarking algorithm
306 - lr: logical ring reduce-scatter then logical ring allgather
307 - rab1: variations of the <a href="https://fs.hlrs.de/projects/par/mpi//myreduce.html">Rabenseifner</a> algorithm: reduce_scatter then allgather
308 - rab2: variations of the <a href="https://fs.hlrs.de/projects/par/mpi//myreduce.html">Rabenseifner</a> algorithm: alltoall then allgather
309 - rab_rsag: variation of the <a href="https://fs.hlrs.de/projects/par/mpi//myreduce.html">Rabenseifner</a> algorithm: recursive doubling
310 reduce_scatter then recursive doubling allgather
311 - rdb: recursive doubling
312 - smp_binomial: binomial tree with smp: binomial intra
313 SMP reduce, inter reduce, inter broadcast then intra broadcast
314 - smp_binomial_pipeline: same with segment size = 4096 bytes
315 - smp_rdb: intra: binomial allreduce, inter: Recursive
316 doubling allreduce, intra: binomial broadcast
317 - smp_rsag: intra: binomial allreduce, inter: reduce-scatter,
318 inter:allgather, intra: binomial broadcast
319 - smp_rsag_lr: intra: binomial allreduce, inter: logical ring
320 reduce-scatter, logical ring inter:allgather, intra: binomial broadcast
321 - smp_rsag_rab: intra: binomial allreduce, inter: rab
322 reduce-scatter, rab inter:allgather, intra: binomial broadcast
323 - redbcast: reduce then broadcast, using default or tuned algorithms if specified
324 - ompi_ring_segmented: ring algorithm used by OpenMPI
325 - mvapich2_rs: rdb for small messages, reduce-scatter then allgather else
326 - mvapich2_two_level: SMP-aware algorithm, with mpich as intra algoritm, and rdb as inter (Change this behavior by using mvapich2 selector to use tuned values)
327 - rab: default `Rabenseifner <https://fs.hlrs.de/projects/par/mpi//myreduce.html>`_ implementation
332 - default: naive one, by default
333 - ompi: use openmpi selector for the reduce_scatter operations
334 - mpich: use mpich selector for the reduce_scatter operations
335 - mvapich2: use mvapich2 selector for the reduce_scatter operations
336 - impi: use intel mpi selector for the reduce_scatter operations
337 - automatic (experimental): use an automatic self-benchmarking algorithm
338 - ompi_basic_recursivehalving: recursive halving version from OpenMPI
339 - ompi_ring: ring version from OpenMPI
340 - mpich_pair: pairwise exchange version from MPICH
341 - mpich_rdb: recursive doubling version from MPICH
342 - mpich_noncomm: only works for power of 2 procs, recursive doubling for noncommutative ops
348 - default: naive one, by default
349 - ompi: use openmpi selector for the allgather operations
350 - mpich: use mpich selector for the allgather operations
351 - mvapich2: use mvapich2 selector for the allgather operations
352 - impi: use intel mpi selector for the allgather operations
353 - automatic (experimental): use an automatic self-benchmarking algorithm
354 - 2dmesh: see alltoall
355 - 3dmesh: see alltoall
356 - bruck: Described by Bruck et.al. in <a href="http://ieeexplore.ieee.org/xpl/articleDetails.jsp?arnumber=642949">
357 Efficient algorithms for all-to-all communications in multiport message-passing systems</a>
358 - GB: Gather - Broadcast (uses tuned version if specified)
359 - loosely_lr: Logical Ring with grouping by core (hardcoded, default
361 - NTSLR: Non Topology Specific Logical Ring
362 - NTSLR_NB: Non Topology Specific Logical Ring, Non Blocking operations
365 - rhv: only power of 2 number of processes
367 - SMP_NTS: gather to root of each SMP, then every root of each SMP node
368 post INTER-SMP Sendrecv, then do INTRA-SMP Bcast for each receiving message,
369 using logical ring algorithm (hardcoded, default processes/SMP: 8)
370 - smp_simple: gather to root of each SMP, then every root of each SMP node
371 post INTER-SMP Sendrecv, then do INTRA-SMP Bcast for each receiving message,
372 using simple algorithm (hardcoded, default processes/SMP: 8)
373 - spreading_simple: from node i, order of communications is i -> i + 1, i ->
374 i + 2, ..., i -> (i + p -1) % P
375 - ompi_neighborexchange: Neighbor Exchange algorithm for allgather.
376 Described by Chen et.al. in `Performance Evaluation of Allgather
377 Algorithms on Terascale Linux Cluster with Fast Ethernet <http://ieeexplore.ieee.org/xpl/articleDetails.jsp?tp=&arnumber=1592302>`_
378 - mvapich2_smp: SMP aware algorithm, performing intra-node gather, inter-node allgather with one process/node, and bcast intra-node
383 - default: naive one, by default
384 - ompi: use openmpi selector for the allgatherv operations
385 - mpich: use mpich selector for the allgatherv operations
386 - mvapich2: use mvapich2 selector for the allgatherv operations
387 - impi: use intel mpi selector for the allgatherv operations
388 - automatic (experimental): use an automatic self-benchmarking algorithm
389 - GB: Gatherv - Broadcast (uses tuned version if specified, but only for Bcast, gatherv is not tuned)
392 - ompi_neighborexchange: see allgather
393 - ompi_bruck: see allgather
394 - mpich_rdb: recursive doubling algorithm from MPICH
395 - mpich_ring: ring algorithm from MPICh - performs differently from the one from STAR-MPI
400 - default: naive one, by default
401 - ompi: use openmpi selector for the bcast operations
402 - mpich: use mpich selector for the bcast operations
403 - mvapich2: use mvapich2 selector for the bcast operations
404 - impi: use intel mpi selector for the bcast operations
405 - automatic (experimental): use an automatic self-benchmarking algorithm
406 - arrival_pattern_aware: root exchanges with the first process to arrive
407 - arrival_pattern_aware_wait: same with slight variation
408 - binomial_tree: binomial tree exchange
409 - flattree: flat tree exchange
410 - flattree_pipeline: flat tree exchange, message split into 8192 bytes pieces
411 - NTSB: Non-topology-specific pipelined binary tree with 8192 bytes pieces
412 - NTSL: Non-topology-specific pipelined linear with 8192 bytes pieces
413 - NTSL_Isend: Non-topology-specific pipelined linear with 8192 bytes pieces, asynchronous communications
414 - scatter_LR_allgather: scatter followed by logical ring allgather
415 - scatter_rdb_allgather: scatter followed by recursive doubling allgather
416 - arrival_scatter: arrival pattern aware scatter-allgather
417 - SMP_binary: binary tree algorithm with 8 cores/SMP
418 - SMP_binomial: binomial tree algorithm with 8 cores/SMP
419 - SMP_linear: linear algorithm with 8 cores/SMP
420 - ompi_split_bintree: binary tree algorithm from OpenMPI, with message split in 8192 bytes pieces
421 - ompi_pipeline: pipeline algorithm from OpenMPI, with message split in 128KB pieces
422 - mvapich2_inter_node: Inter node default mvapich worker
423 - mvapich2_intra_node: Intra node default mvapich worker
424 - mvapich2_knomial_intra_node: k-nomial intra node default mvapich worker. default factor is 4.
429 .. warning:: This is still very experimental.
431 An automatic version is available for each collective (or even as a selector). This specific
432 version will loop over all other implemented algorithm for this particular collective, and apply
433 them while benchmarking the time taken for each process. It will then output the quickest for
434 each process, and the global quickest. This is still unstable, and a few algorithms which need
435 specific number of nodes may crash.
440 To add a new algorithm, one should check in the src/smpi/colls folder
441 how other algorithms are coded. Using plain MPI code inside Simgrid
442 can't be done, so algorithms have to be changed to use smpi version of
443 the calls instead (MPI_Send will become smpi_mpi_send). Some functions
444 may have different signatures than their MPI counterpart, please check
445 the other algorithms or contact us using the `>SimGrid
446 developers mailing list <http://lists.gforge.inria.fr/mailman/listinfo/simgrid-devel>`_.
448 Example: adding a "pair" version of the Alltoall collective.
450 - Implement it in a file called alltoall-pair.c in the src/smpi/colls folder. This file should include colls_private.hpp.
452 - The name of the new algorithm function should be smpi_coll_tuned_alltoall_pair, with the same signature as MPI_Alltoall.
454 - Once the adaptation to SMPI code is done, add a reference to the file ("src/smpi/colls/alltoall-pair.c") in the SMPI_SRC part of the DefinePackages.cmake file inside buildtools/cmake, to allow the file to be built and distributed.
456 - To register the new version of the algorithm, simply add a line to the corresponding macro in src/smpi/colls/cools.h ( add a "COLL_APPLY(action, COLL_ALLTOALL_SIG, pair)" to the COLL_ALLTOALLS macro ). The algorithm should now be compiled and be selected when using --cfg=smpi/alltoall:pair at runtime.
458 - To add a test for the algorithm inside Simgrid's test suite, juste add the new algorithm name in the ALLTOALL_COLL list found inside teshsuite/smpi/CMakeLists.txt . When running ctest, a test for the new algorithm should be generated and executed. If it does not pass, please check your code or contact us.
460 - Please submit your patch for inclusion in SMPI, for example through a pull request on GitHub or directly per email.
463 Tracing of Internal Communications
464 ^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
466 By default, the collective operations are traced as a unique operation
467 because tracing all point-to-point communications composing them could
468 result in overloaded, hard to interpret traces. If you want to debug
469 and compare collective algorithms, you should set the
470 ``tracing/smpi/internals`` configuration item to 1 instead of 0.
472 Here are examples of two alltoall collective algorithms runs on 16 nodes,
473 the first one with a ring algorithm, the second with a pairwise one.
475 .. image:: /img/smpi_simgrid_alltoall_ring_16.png
478 Alltoall on 16 Nodes with the Ring Algorithm.
480 .. image:: /img/smpi_simgrid_alltoall_pair_16.png
483 Alltoall on 16 Nodes with the Pairwise Algorithm.
485 -------------------------
486 What can run within SMPI?
487 -------------------------
489 You can run unmodified MPI applications (both C/C++ and Fortran) within
490 SMPI, provided that you only use MPI calls that we implemented. Global
491 variables should be handled correctly on Linux systems.
497 Our coverage of the interface is very decent, but still incomplete;
498 Given the size of the MPI standard, we may well never manage to
499 implement absolutely all existing primitives. Currently, we have
500 almost no support for I/O primitives, but we still pass a very large
501 amount of the MPICH coverage tests.
503 The full list of not yet implemented functions is documented in the
504 file `include/smpi/smpi.h
505 <https://framagit.org/simgrid/simgrid/tree/master/include/smpi/smpi.h>`_
506 in your version of SimGrid, between two lines containing the ``FIXME``
507 marker. If you really miss a feature, please get in touch with us: we
508 can guide you though the SimGrid code to help you implementing it, and
509 we'd be glad to integrate your contribution to the main project.
511 .. _SMPI_what_globals:
513 .................................
514 Privatization of global variables
515 .................................
517 Concerning the globals, the problem comes from the fact that usually,
518 MPI processes run as real UNIX processes while they are all folded
519 into threads of a unique system process in SMPI. Global variables are
520 usually private to each MPI process while they become shared between
521 the processes in SMPI. The problem and some potential solutions are
522 discussed in this article: `Automatic Handling of Global Variables for
523 Multi-threaded MPI Programs
524 <http://charm.cs.illinois.edu/newPapers/11-23/paper.pdf>` (note that
525 this article does not deal with SMPI but with a competing solution
526 called AMPI that suffers of the same issue). This point used to be
527 problematic in SimGrid, but the problem should now be handled
528 automatically on Linux.
530 Older versions of SimGrid came with a script that automatically
531 privatized the globals through static analysis of the source code. But
532 our implementation was not robust enough to be used in production, so
533 it was removed at some point. Currently, SMPI comes with two
534 privatization mechanisms that you can :ref:`select at runtime
535 <options_smpi_privatization>`_. The dlopen approach is used by
536 default as it is much faster and still very robust. The mmap approach
537 is an older approach that proves to be slower.
539 With the **mmap approach**, SMPI duplicates and dynamically switch the
540 ``.data`` and ``.bss`` segments of the ELF process when switching the
541 MPI ranks. This allows each ranks to have its own copy of the global
542 variables. No copy actually occures as this mechanism uses ``mmap()``
543 for efficiency. This mechanism is considered to be very robust on all
544 systems supporting ``mmap()`` (Linux and most BSDs). Its performance
545 is questionable since each context switch between MPI ranks induces
546 several syscalls to change the ``mmap`` that redirects the ``.data``
547 and ``.bss`` segments to the copies of the new rank. The code will
548 also be copied several times in memory, inducing a slight increase of
551 Another limitation is that SMPI only accounts for global variables
552 defined in the executable. If the processes use external global
553 variables from dynamic libraries, they won't be switched
554 correctly. The easiest way to solve this is to statically link against
555 the library with these globals. This way, each MPI rank will get its
556 own copy of these libraries. Of course you should never statically
557 link against the SimGrid library itself.
559 With the **dlopen approach**, SMPI loads several copies of the same
560 executable in memory as if it were a library, so that the global
561 variables get naturally dupplicated. It first requires the executable
562 to be compiled as a relocatable binary, which is less common for
563 programs than for libraries. But most distributions are now compiled
564 this way for security reason as it allows to randomize the address
565 space layout. It should thus be safe to compile most (any?) program
566 this way. The second trick is that the dynamic linker refuses to link
567 the exact same file several times, be it a library or a relocatable
568 executable. It makes perfectly sense in the general case, but we need
569 to circumvent this rule of thumb in our case. To that extend, the
570 binary is copied in a temporary file before being re-linked against.
571 ``dlmopen()`` cannot be used as it only allows 256 contextes, and as it
572 would also dupplicate simgrid itself.
574 This approach greatly speeds up the context switching, down to about
575 40 CPU cycles with our raw contextes, instead of requesting several
576 syscalls with the ``mmap()`` approach. Another advantage is that it
577 permits to run the SMPI contexts in parallel, which is obviously not
578 possible with the ``mmap()`` approach. It was tricky to implement, but
579 we are not aware of any flaws, so smpirun activates it by default.
581 In the future, it may be possible to further reduce the memory and
582 disk consumption. It seems that we could `punch holes
583 <https://lwn.net/Articles/415889/>`_ in the files before dl-loading
584 them to remove the code and constants, and mmap these area onto a
585 unique copy. If done correctly, this would reduce the disk- and
586 memory- usage to the bare minimum, and would also reduce the pressure
587 on the CPU instruction cache. See the `relevant bug
588 <https://github.com/simgrid/simgrid/issues/137>`_ on github for
589 implementation leads.\n
591 Also, currently, only the binary is copied and dlopen-ed for each MPI
592 rank. We could probably extend this to external dependencies, but for
593 now, any external dependencies must be statically linked into your
594 application. As usual, simgrid itself shall never be statically linked
595 in your app. You don't want to give a copy of SimGrid to each MPI rank:
596 that's ways too much for them to deal with.
598 .. todo: speak of smpi/privatize-libs here
600 ----------------------------------------------
601 Adapting your MPI code for further scalability
602 ----------------------------------------------
604 As detailed in the `reference article
605 <http://hal.inria.fr/hal-01415484>`_, you may want to adapt your code
606 to improve the simulation performance. But these tricks may seriously
607 hinder the result quality (or even prevent the app to run) if used
608 wrongly. We assume that if you want to simulate an HPC application,
609 you know what you are doing. Don't prove us wrong!
611 ..............................
612 Reducing your memory footprint
613 ..............................
615 If you get short on memory (the whole app is executed on a single node when
616 simulated), you should have a look at the SMPI_SHARED_MALLOC and
617 SMPI_SHARED_FREE macros. It allows to share memory areas between processes: The
618 purpose of these macro is that the same line malloc on each process will point
619 to the exact same memory area. So if you have a malloc of 2M and you have 16
620 processes, this macro will change your memory consumption from 2M*16 to 2M
621 only. Only one block for all processes.
623 If your program is ok with a block containing garbage value because all
624 processes write and read to the same place without any kind of coordination,
625 then this macro can dramatically shrink your memory consumption. For example,
626 that will be very beneficial to a matrix multiplication code, as all blocks will
627 be stored on the same area. Of course, the resulting computations will useless,
628 but you can still study the application behavior this way.
630 Naturally, this won't work if your code is data-dependent. For example, a Jacobi
631 iterative computation depends on the result computed by the code to detect
632 convergence conditions, so turning them into garbage by sharing the same memory
633 area between processes does not seem very wise. You cannot use the
634 SMPI_SHARED_MALLOC macro in this case, sorry.
636 This feature is demoed by the example file
637 `examples/smpi/NAS/dt.c <https://framagit.org/simgrid/simgrid/tree/master/examples/smpi/NAS/dt.c>`_
639 .........................
640 Toward Faster Simulations
641 .........................
643 If your application is too slow, try using SMPI_SAMPLE_LOCAL,
644 SMPI_SAMPLE_GLOBAL and friends to indicate which computation loops can
645 be sampled. Some of the loop iterations will be executed to measure
646 their duration, and this duration will be used for the subsequent
647 iterations. These samples are done per processor with
648 SMPI_SAMPLE_LOCAL, and shared between all processors with
649 SMPI_SAMPLE_GLOBAL. Of course, none of this will work if the execution
650 time of your loop iteration are not stable.
652 This feature is demoed by the example file
653 `examples/smpi/NAS/ep.c <https://framagit.org/simgrid/simgrid/tree/master/examples/smpi/NAS/ep.c>`_
655 .............................
656 Ensuring Accurate Simulations
657 .............................
659 Out of the box, SimGrid may give you fairly accurate results, but
660 there is a plenty of factors that could go wrong and make your results
661 inaccurate or even plainly wrong. Actually, you can only get accurate
662 results of a nicely built model, including both the system hardware
663 and your application. Such models are hard to pass over and reuse in
664 other settings, because elements that are not relevant to an
665 application (say, the latency of point-to-point communications,
666 collective operation implementation details or CPU-network
667 interaction) may be irrelevant to another application. The dream of
668 the perfect model, encompassing every aspects is only a chimera, as
669 the only perfect model of the reality is the reality. If you go for
670 simulation, then you have to ignore some irrelevant aspects of the
671 reality, but which aspects are irrelevant is actually
672 application-dependent...
674 The only way to assess whether your settings provide accurate results
675 is to double-check these results. If possible, you should first run
676 the same experiment in simulation and in real life, gathering as much
677 information as you can. Try to understand the discrepancies in the
678 results that you observe between both settings (visualization can be
679 precious for that). Then, try to modify your model (of the platform,
680 of the collective operations) to reduce the most preeminent differences.
682 If the discrepancies come from the computing time, try adapting the
683 ``smpi/host-speed``: reduce it if your simulation runs faster than in
684 reality. If the error come from the communication, then you need to
685 fiddle with your platform file.
687 Be inventive in your modeling. Don't be afraid if the names given by
688 SimGrid does not match the real names: we got very good results by
689 modeling multicore/GPU machines with a set of separate hosts
690 interconnected with very fast networks (but don't trust your model
691 because it has the right names in the right place either).
693 Finally, you may want to check `this article
694 <https://hal.inria.fr/hal-00907887>`_ on the classical pitfalls in
695 modeling distributed systems.
697 -------------------------
698 Troubleshooting with SMPI
699 -------------------------
701 .................................
702 ./configure refuses to use smpicc
703 .................................
705 If your ``./configure`` reports that the compiler is not
706 functional or that you are cross-compiling, try to define the
707 ``SMPI_PRETEND_CC`` environment variable before running the
710 .. code-block:: shell
712 SMPI_PRETEND_CC=1 ./configure # here come the configure parameters
715 Indeed, the programs compiled with ``smpicc`` cannot be executed
716 without ``smpirun`` (they are shared libraries and do weird things on
717 startup), while configure wants to test them directly. With
718 ``SMPI_PRETEND_CC`` smpicc does not compile as shared, and the SMPI
719 initialization stops and returns 0 before doing anything that would
720 fail without ``smpirun``.
724 Make sure that SMPI_PRETEND_CC is only set when calling ./configure,
725 not during the actual execution, or any program compiled with smpicc
726 will stop before starting.
728 ..............................................
729 ./configure does not pick smpicc as a compiler
730 ..............................................
732 In addition to the previous answers, some projects also need to be
733 explicitely told what compiler to use, as follows:
735 .. code-block:: shell
737 SMPI_PRETEND_CC=1 ./configure CC=smpicc # here come the other configure parameters
740 Maybe your configure is using another variable, such as ``cc`` (in
741 lower case) or similar. Just check the logs.
743 .....................................
744 error: unknown type name 'useconds_t'
745 .....................................
747 Try to add ``-D_GNU_SOURCE`` to your compilation line to get ride
750 The reason is that SMPI provides its own version of ``usleep(3)``
751 to override it and to block in the simulation world, not in the real
752 one. It needs the ``useconds_t`` type for that, which is declared
753 only if you declare ``_GNU_SOURCE`` before including
754 ``unistd.h``. If your project includes that header file before
755 SMPI, then you need to ensure that you pass the right configuration
756 defines as advised above.
762 -----------------------------
763 Trace Replay and Offline SMPI
764 -----------------------------
766 Although SMPI is often used for :ref:`online simulation
767 <SMPI_online>`, where the application is executed for real, you can
768 also go for offline simulation through trace replay.
770 SimGrid uses time-independent traces, in which each actor is given a
771 script of the actions to do sequentially. These trace files can
772 actually be captured with the online version of SMPI, as follows:
774 .. code-block:: shell
776 $ smpirun -trace-ti --cfg=tracing/filename:LU.A.32 -np 32 -platform ../cluster_backbone.xml bin/lu.A.32
778 The produced trace is composed of a file ``LU.A.32`` and a folder
779 ``LU.A.32_files``. The file names don't match with the MPI ranks, but
782 To replay this with SMPI, you need to first compile the provided
783 ``smpi_replay.cpp`` file, that comes from
784 `simgrid/examples/smpi/replay
785 <https://framagit.org/simgrid/simgrid/tree/master/examples/smpi/replay>`_.
787 .. code-block:: shell
789 $ smpicxx ../replay.cpp -O3 -o ../smpi_replay
791 Afterward, you can replay your trace in SMPI as follows:
793 $ smpirun -np 32 -platform ../cluster_torus.xml -ext smpi_replay ../smpi_replay LU.A.32
795 All the outputs are gone, as the application is not really simulated
796 here. Its trace is simply replayed. But if you visualize the live
797 simulation and the replay, you will see that the behavior is
798 unchanged. The simulation does not run much faster on this very
799 example, but this becomes very interesting when your application
800 is computationally hungry.