Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
start reformating the manual for the AS -> NetZone transition (ignorable changes)
[simgrid.git] / doc / doxygen / module-smpi.doc
1 /** 
2 @defgroup SMPI_API      SMPI: Simulate real MPI applications
3 @brief Programming environment for the simulation of MPI applications
4
5 @tableofcontents
6
7 [TOC]
8
9 This programming environment enables the study of MPI application by
10 emulating them on top of the SimGrid simulator. This is particularly
11 interesting to study existing MPI applications within the comfort of
12 the simulator. The motivation for this work is detailed in the
13 reference article (available at http://hal.inria.fr/inria-00527150).
14
15
16 Our goal is to enable the study of **unmodified MPI applications**,
17 and even if some constructs and features are still missing, we
18 consider SMPI to be stable and usable in production.  For **further
19 scalability**, you may modify your code to speed up your studies or
20 save memory space.  Improved **simulation accuracy** requires some
21 specific care from you.
22
23  - @ref SMPI_use
24    - @ref SMPI_use_compile
25    - @ref SMPI_use_exec
26    - @ref SMPI_use_colls
27      - @ref SMPI_use_colls_algos
28      - @ref SMPI_use_colls_tracing
29  - @ref SMPI_what
30    - @ref SMPI_what_coverage
31    - @ref SMPI_what_globals
32  - @ref SMPI_adapting
33    - @ref SMPI_adapting_size
34    - @ref SMPI_adapting_speed
35  - @ref SMPI_accuracy
36
37
38 @section SMPI_use Using SMPI
39
40 If you're absolutely new to MPI, you should first take our online
41 [SMPI CourseWare](https://simgrid.github.io/SMPI_CourseWare/), and/or
42 take a MPI course in your favorite university.  If you already know
43 MPI, SMPI should sound very familiar to you: Use smpicc instead of
44 mpicc, and smpirun instead of mpirun, and you're almost set.  Once you
45 get a virtual platform description (see @ref platform), you're good to
46 go.
47
48 @subsection SMPI_use_compile Compiling your code
49
50 For that, simply use <tt>smpicc</tt> as a compiler just
51 like you use mpicc with other MPI implementations. This script
52 still calls your default compiler (gcc, clang, ...) and adds the right
53 compilation flags along the way.
54
55 Alas, some building infrastructures cannot cope with that and your
56 <tt>./configure</tt> may fail, reporting that the compiler is not
57 functional. If this happens, define the <tt>SMPI_PRETEND_CC</tt>
58 environment variable before running the configuration. Do not define
59 it when using SMPI!
60
61 @verbatim
62 SMPI_PRETEND_CC=1 ./configure # here come the configure parameters
63 make
64 @endverbatim
65
66 \warn
67   Again, make sure that SMPI_PRETEND_CC is not set when you actually 
68   compile your application. It is just a work-around for some configure-scripts
69   and replaces some internals by "return 0;". Your simulation will not
70   work with this variable set!
71
72 @subsection SMPI_use_exec Executing your code on the simulator
73
74 Use the <tt>smpirun</tt> script as follows for that:
75
76 @verbatim
77 smpirun -hostfile my_hostfile.txt -platform my_platform.xml ./program -blah
78 @endverbatim
79
80  - <tt>my_hostfile.txt</tt> is a classical MPI hostfile (that is, this
81    file lists the machines on which the processes must be dispatched, one
82    per line)
83  - <tt>my_platform.xml</tt> is a classical SimGrid platform file. Of
84    course, the hosts of the hostfile must exist in the provided
85    platform.
86  - <tt>./program</tt> is the MPI program to simulate, that you
87    compiled with <tt>smpicc</tt>
88  - <tt>-blah</tt> is a command-line parameter passed to this program.
89
90 <tt>smpirun</tt> accepts other parameters, such as <tt>-np</tt> if you
91 don't want to use all the hosts defined in the hostfile, <tt>-map</tt>
92 to display on which host each rank gets mapped of <tt>-trace</tt> to
93 activate the tracing during the simulation. You can get the full list
94 by running
95
96 @verbatim
97 smpirun -help
98 @endverbatim
99
100 @subsection SMPI_use_colls Simulating collective operations
101
102 MPI collective operations are crucial to the performance of MPI
103 applications and must be carefully optimized according to many
104 parameters. Every existing implementation provides several algorithms
105 for each collective operation, and selects by default the best suited
106 one, depending on the sizes sent, the number of nodes, the
107 communicator, or the communication library being used.  These
108 decisions are based on empirical results and theoretical complexity
109 estimation, and are very different between MPI implementations. In
110 most cases, the users can also manually tune the algorithm used for
111 each collective operation.
112
113 SMPI can simulate the behavior of several MPI implementations:
114 OpenMPI, MPICH,
115 <a href="http://star-mpi.sourceforge.net/">STAR-MPI</a>, and
116 MVAPICH2. For that, it provides 115 collective algorithms and several
117 selector algorithms, that were collected directly in the source code
118 of the targeted MPI implementations.
119
120 You can switch the automatic selector through the
121 \c smpi/coll_selector configuration item. Possible values:
122
123  - <b>ompi</b>: default selection logic of OpenMPI (version 1.7)
124  - <b>mpich</b>: default selection logic of MPICH (version 3.0.4)
125  - <b>mvapich2</b>: selection logic of MVAPICH2 (version 1.9) tuned
126    on the Stampede cluster   
127  - <b>impi</b>: preliminary version of an Intel MPI selector (version
128    4.1.3, also tuned for the Stampede cluster). Due the closed source
129    nature of Intel MPI, some of the algorithms described in the
130    documentation are not available, and are replaced by mvapich ones.   
131  - <b>default</b>: legacy algorithms used in the earlier days of
132    SimGrid. Do not use for serious perform performance studies.
133
134
135 @subsubsection SMPI_use_colls_algos Available algorithms
136
137 You can also pick the algorithm used for each collective with the
138 corresponding configuration item. For example, to use the pairwise
139 alltoall algorithm, one should add \c --cfg=smpi/alltoall:pair to the
140 line. This will override the selector (if any) for this algorithm.
141 It means that the selected algorithm will be used 
142
143 Warning: Some collective may require specific conditions to be
144 executed correctly (for instance having a communicator with a power of
145 two number of nodes only), which are currently not enforced by
146 Simgrid. Some crashes can be expected while trying these algorithms
147 with unusual sizes/parameters
148
149 #### MPI_Alltoall
150
151 Most of these are best described in <a href="http://www.cs.arizona.edu/~dkl/research/papers/ics06.pdf">STAR-MPI</a>
152
153  - default: naive one, by default
154  - ompi: use openmpi selector for the alltoall operations
155  - mpich: use mpich selector for the alltoall operations
156  - mvapich2: use mvapich2 selector for the alltoall operations
157  - impi: use intel mpi selector for the alltoall operations
158  - automatic (experimental): use an automatic self-benchmarking algorithm 
159  - 2dmesh: organizes the nodes as a two dimensional mesh, and perform allgather 
160    along the dimensions
161  - 3dmesh: adds a third dimension to the previous algorithm
162  - rdb: recursive doubling : extends the mesh to a nth dimension, each one 
163    containing two nodes
164  - pair: pairwise exchange, only works for power of 2 procs, size-1 steps,
165    each process sends and receives from the same process at each step
166  - pair_light_barrier: same, with small barriers between steps to avoid
167    contention
168  - pair_mpi_barrier: same, with MPI_Barrier used
169  - pair_one_barrier: only one barrier at the beginning
170  - ring: size-1 steps, at each step a process send to process (n+i)%size, and receives from (n-i)%size
171  - ring_light_barrier: same, with small barriers between some phases to avoid contention
172  - ring_mpi_barrier: same, with MPI_Barrier used
173  - ring_one_barrier: only one barrier at the beginning
174  - basic_linear: posts all receives and all sends,
175 starts the communications, and waits for all communication to finish
176  - mvapich2_scatter_dest: isend/irecv with scattered destinations, posting only a few messages at the same time
177
178 #### MPI_Alltoallv
179
180  - default: naive one, by default
181  - ompi: use openmpi selector for the alltoallv operations
182  - mpich: use mpich selector for the alltoallv operations
183  - mvapich2: use mvapich2 selector for the alltoallv operations
184  - impi: use intel mpi selector for the alltoallv operations
185  - automatic (experimental): use an automatic self-benchmarking algorithm 
186  - bruck: same as alltoall
187  - pair: same as alltoall
188  - pair_light_barrier: same as alltoall
189  - pair_mpi_barrier: same as alltoall
190  - pair_one_barrier: same as alltoall
191  - ring: same as alltoall
192  - ring_light_barrier: same as alltoall
193  - ring_mpi_barrier: same as alltoall
194  - ring_one_barrier: same as alltoall
195  - ompi_basic_linear: same as alltoall
196
197 #### MPI_Gather
198
199  - default: naive one, by default
200  - ompi: use openmpi selector for the gather operations
201  - mpich: use mpich selector for the gather operations
202  - mvapich2: use mvapich2 selector for the gather operations
203  - impi: use intel mpi selector for the gather operations
204  - automatic (experimental): use an automatic self-benchmarking algorithm 
205 which will iterate over all implemented versions and output the best
206  - ompi_basic_linear: basic linear algorithm from openmpi, each process sends to the root
207  - ompi_binomial: binomial tree algorithm
208  - ompi_linear_sync: same as basic linear, but with a synchronization at the
209  beginning and message cut into two segments.
210  - 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.
211
212 #### MPI_Barrier
213
214  - default: naive one, by default
215  - ompi: use openmpi selector for the barrier operations
216  - mpich: use mpich selector for the barrier operations
217  - mvapich2: use mvapich2 selector for the barrier operations
218  - impi: use intel mpi selector for the barrier operations
219  - automatic (experimental): use an automatic self-benchmarking algorithm 
220  - ompi_basic_linear: all processes send to root
221  - ompi_two_procs: special case for two processes
222  - ompi_bruck: nsteps = sqrt(size), at each step, exchange data with rank-2^k and rank+2^k
223  - ompi_recursivedoubling: recursive doubling algorithm
224  - ompi_tree: recursive doubling type algorithm, with tree structure
225  - ompi_doublering: double ring algorithm
226  - mvapich2_pair: pairwise algorithm
227
228 #### MPI_Scatter
229
230  - default: naive one, by default
231  - ompi: use openmpi selector for the scatter operations
232  - mpich: use mpich selector for the scatter operations
233  - mvapich2: use mvapich2 selector for the scatter operations
234  - impi: use intel mpi selector for the scatter operations
235  - automatic (experimental): use an automatic self-benchmarking algorithm 
236  - ompi_basic_linear: basic linear scatter 
237  - ompi_binomial: binomial tree scatter
238  - 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. 
239  - 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.
240
241 #### MPI_Reduce
242
243  - default: naive one, by default
244  - ompi: use openmpi selector for the reduce operations
245  - mpich: use mpich selector for the reduce operations
246  - mvapich2: use mvapich2 selector for the reduce operations
247  - impi: use intel mpi selector for the reduce operations
248  - automatic (experimental): use an automatic self-benchmarking algorithm 
249  - arrival_pattern_aware: root exchanges with the first process to arrive
250  - binomial: uses a binomial tree
251  - flat_tree: uses a flat tree
252  - NTSL: Non-topology-specific pipelined linear-bcast function 
253    0->1, 1->2 ,2->3, ....., ->last node: in a pipeline fashion, with segments
254  of 8192 bytes
255  - scatter_gather: scatter then gather
256  - ompi_chain: openmpi reduce algorithms are built on the same basis, but the
257  topology is generated differently for each flavor
258 chain = chain with spacing of size/2, and segment size of 64KB 
259  - ompi_pipeline: same with pipeline (chain with spacing of 1), segment size 
260 depends on the communicator size and the message size
261  - ompi_binary: same with binary tree, segment size of 32KB
262  - ompi_in_order_binary: same with binary tree, enforcing order on the 
263 operations
264  - ompi_binomial: same with binomial algo (redundant with default binomial 
265 one in most cases)
266  - ompi_basic_linear: basic algorithm, each process sends to root
267  - mvapich2_knomial: k-nomial algorithm. Default factor is 4 (mvapich2 selector adapts it through tuning)
268  - 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.
269  - rab: <a href="https://fs.hlrs.de/projects/par/mpi//myreduce.html">Rabenseifner</a>'s reduce algorithm 
270
271 #### MPI_Allreduce
272
273  - default: naive one, by default
274  - ompi: use openmpi selector for the allreduce operations
275  - mpich: use mpich selector for the allreduce operations
276  - mvapich2: use mvapich2 selector for the allreduce operations
277  - impi: use intel mpi selector for the allreduce operations
278  - automatic (experimental): use an automatic self-benchmarking algorithm 
279  - lr: logical ring reduce-scatter then logical ring allgather
280  - rab1: variations of the  <a href="https://fs.hlrs.de/projects/par/mpi//myreduce.html">Rabenseifner</a> algorithm: reduce_scatter then allgather
281  - rab2: variations of the  <a href="https://fs.hlrs.de/projects/par/mpi//myreduce.html">Rabenseifner</a> algorithm: alltoall then allgather
282  - rab_rsag: variation of the  <a href="https://fs.hlrs.de/projects/par/mpi//myreduce.html">Rabenseifner</a> algorithm: recursive doubling 
283 reduce_scatter then recursive doubling allgather 
284  - rdb: recursive doubling
285  - smp_binomial: binomial tree with smp: binomial intra 
286 SMP reduce, inter reduce, inter broadcast then intra broadcast
287  - smp_binomial_pipeline: same with segment size = 4096 bytes
288  - smp_rdb: intra: binomial allreduce, inter: Recursive 
289 doubling allreduce, intra: binomial broadcast
290  - smp_rsag: intra: binomial allreduce, inter: reduce-scatter, 
291 inter:allgather, intra: binomial broadcast
292  - smp_rsag_lr: intra: binomial allreduce, inter: logical ring 
293 reduce-scatter, logical ring inter:allgather, intra: binomial broadcast
294  - smp_rsag_rab: intra: binomial allreduce, inter: rab
295 reduce-scatter, rab inter:allgather, intra: binomial broadcast
296  - redbcast: reduce then broadcast, using default or tuned algorithms if specified
297  - ompi_ring_segmented: ring algorithm used by OpenMPI
298  - mvapich2_rs: rdb for small messages, reduce-scatter then allgather else
299  - 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)
300  - rab: default <a href="https://fs.hlrs.de/projects/par/mpi//myreduce.html">Rabenseifner</a> implementation
301
302 #### MPI_Reduce_scatter
303
304  - default: naive one, by default
305  - ompi: use openmpi selector for the reduce_scatter operations
306  - mpich: use mpich selector for the reduce_scatter operations
307  - mvapich2: use mvapich2 selector for the reduce_scatter operations
308  - impi: use intel mpi selector for the reduce_scatter operations
309  - automatic (experimental): use an automatic self-benchmarking algorithm 
310  - ompi_basic_recursivehalving: recursive halving version from OpenMPI
311  - ompi_ring: ring version from OpenMPI
312  - mpich_pair: pairwise exchange version from MPICH
313  - mpich_rdb: recursive doubling version from MPICH
314  - mpich_noncomm: only works for power of 2 procs, recursive doubling for noncommutative ops
315
316
317 #### MPI_Allgather
318
319  - default: naive one, by default
320  - ompi: use openmpi selector for the allgather operations
321  - mpich: use mpich selector for the allgather operations
322  - mvapich2: use mvapich2 selector for the allgather operations
323  - impi: use intel mpi selector for the allgather operations
324  - automatic (experimental): use an automatic self-benchmarking algorithm 
325  - 2dmesh: see alltoall
326  - 3dmesh: see alltoall
327  - bruck: Described by Bruck et.al. in <a href="http://ieeexplore.ieee.org/xpl/articleDetails.jsp?arnumber=642949">
328 Efficient algorithms for all-to-all communications in multiport message-passing systems</a> 
329  - GB: Gather - Broadcast (uses tuned version if specified)
330  - loosely_lr: Logical Ring with grouping by core (hardcoded, default 
331 processes/node: 4)
332  - NTSLR: Non Topology Specific Logical Ring
333  - NTSLR_NB: Non Topology Specific Logical Ring, Non Blocking operations
334  - pair: see alltoall
335  - rdb: see alltoall
336  - rhv: only power of 2 number of processes
337  - ring: see alltoall
338  - SMP_NTS: gather to root of each SMP, then every root of each SMP node 
339 post INTER-SMP Sendrecv, then do INTRA-SMP Bcast for each receiving message, 
340 using logical ring algorithm (hardcoded, default processes/SMP: 8)
341  - smp_simple: gather to root of each SMP, then every root of each SMP node 
342 post INTER-SMP Sendrecv, then do INTRA-SMP Bcast for each receiving message, 
343 using simple algorithm (hardcoded, default processes/SMP: 8)
344  - spreading_simple: from node i, order of communications is i -> i + 1, i ->
345  i + 2, ..., i -> (i + p -1) % P
346  - ompi_neighborexchange: Neighbor Exchange algorithm for allgather. 
347 Described by Chen et.al. in  <a href="http://ieeexplore.ieee.org/xpl/articleDetails.jsp?tp=&arnumber=1592302">Performance Evaluation of Allgather Algorithms on Terascale Linux Cluster with Fast Ethernet</a>
348  - mvapich2_smp: SMP aware algorithm, performing intra-node gather, inter-node allgather with one process/node, and bcast intra-node
349
350
351 ####_Allgatherv
352
353  - default: naive one, by default
354  - ompi: use openmpi selector for the allgatherv operations
355  - mpich: use mpich selector for the allgatherv operations
356  - mvapich2: use mvapich2 selector for the allgatherv operations
357  - impi: use intel mpi selector for the allgatherv operations
358  - automatic (experimental): use an automatic self-benchmarking algorithm 
359  - GB: Gatherv - Broadcast (uses tuned version if specified, but only for 
360 Bcast, gatherv is not tuned)
361  - pair: see alltoall
362  - ring: see alltoall
363  - ompi_neighborexchange: see allgather
364  - ompi_bruck: see allgather
365  - mpich_rdb: recursive doubling algorithm from MPICH
366  - mpich_ring: ring algorithm from MPICh - performs differently from the  one from STAR-MPI
367
368 #### MPI_Bcast
369
370  - default: naive one, by default
371  - ompi: use openmpi selector for the bcast operations
372  - mpich: use mpich selector for the bcast operations
373  - mvapich2: use mvapich2 selector for the bcast operations
374  - impi: use intel mpi selector for the bcast operations
375  - automatic (experimental): use an automatic self-benchmarking algorithm 
376  - arrival_pattern_aware: root exchanges with the first process to arrive
377  - arrival_pattern_aware_wait: same with slight variation
378  - binomial_tree: binomial tree exchange
379  - flattree: flat tree exchange
380  - flattree_pipeline: flat tree exchange, message split into 8192 bytes pieces
381  - NTSB: Non-topology-specific pipelined binary tree with 8192 bytes pieces
382  - NTSL: Non-topology-specific pipelined linear with 8192 bytes pieces
383  - NTSL_Isend: Non-topology-specific pipelined linear with 8192 bytes pieces, asynchronous communications
384  - scatter_LR_allgather: scatter followed by logical ring allgather
385  - scatter_rdb_allgather: scatter followed by recursive doubling allgather
386  - arrival_scatter: arrival pattern aware scatter-allgather
387  - SMP_binary: binary tree algorithm with 8 cores/SMP
388  - SMP_binomial: binomial tree algorithm with 8 cores/SMP
389  - SMP_linear: linear algorithm with 8 cores/SMP
390  - ompi_split_bintree: binary tree algorithm from OpenMPI, with message split in 8192 bytes pieces
391  - ompi_pipeline: pipeline algorithm from OpenMPI, with message split in 128KB pieces
392  - mvapich2_inter_node: Inter node default mvapich worker 
393  - mvapich2_intra_node: Intra node default mvapich worker
394  - mvapich2_knomial_intra_node:  k-nomial intra node default mvapich worker. default factor is 4.
395
396 #### Automatic evaluation 
397
398 (Warning: This is experimental and may be removed or crash easily)
399
400 An automatic version is available for each collective (or even as a selector). This specific 
401 version will loop over all other implemented algorithm for this particular collective, and apply 
402 them while benchmarking the time taken for each process. It will then output the quickest for 
403 each process, and the global quickest. This is still unstable, and a few algorithms which need 
404 specific number of nodes may crash.
405
406 #### Adding an algorithm
407
408 To add a new algorithm, one should check in the src/smpi/colls folder how other algorithms 
409 are coded. Using plain MPI code inside Simgrid can't be done, so algorithms have to be 
410 changed to use smpi version of the calls instead (MPI_Send will become smpi_mpi_send). Some functions may have different signatures than their MPI counterpart, please check the other algorithms or contact us using <a href="http://lists.gforge.inria.fr/mailman/listinfo/simgrid-devel">SimGrid developers mailing list</a>.
411
412 Example: adding a "pair" version of the Alltoall collective.
413
414  - Implement it in a file called alltoall-pair.c in the src/smpi/colls folder. This file should include colls_private.h.
415
416  - The name of the new algorithm function should be smpi_coll_tuned_alltoall_pair, with the same signature as MPI_Alltoall.
417
418  - 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.
419
420  - 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.
421
422  - 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.
423
424  - Please submit your patch for inclusion in SMPI, for example through a pull request on GitHub or directly per email.
425
426 @subsubsection SMPI_use_colls_tracing Tracing of internal communications
427
428 By default, the collective operations are traced as a unique operation
429 because tracing all point-to-point communications composing them could
430 result in overloaded, hard to interpret traces. If you want to debug
431 and compare collective algorithms, you should set the
432 \c tracing/smpi/internals configuration item to 1 instead of 0.
433
434 Here are examples of two alltoall collective algorithms runs on 16 nodes, 
435 the first one with a ring algorithm, the second with a pairwise one:
436
437 @htmlonly
438 <a href="smpi_simgrid_alltoall_ring_16.png" border=0><img src="smpi_simgrid_alltoall_ring_16.png" width="30%" border=0 align="center"></a>
439 <a href="smpi_simgrid_alltoall_pair_16.png" border=0><img src="smpi_simgrid_alltoall_pair_16.png" width="30%" border=0 align="center"></a>
440 <br/>
441 @endhtmlonly
442
443 @section SMPI_what What can run within SMPI?
444
445 You can run unmodified MPI applications (both C and Fortran) within
446 SMPI, provided that you only use MPI calls that we implemented. Global
447 variables should be handled correctly on Linux systems.
448
449 @subsection SMPI_what_coverage MPI coverage of SMPI
450
451 Our coverage of the interface is very decent, but still incomplete;
452 Given the size of the MPI standard, we may well never manage to 
453 implement absolutely all existing primitives. Currently, we have a
454 very sparse support for one-sided communications, and almost none for
455 I/O primitives. But our coverage is still very decent: we pass a very
456 large amount of the MPICH coverage tests.
457
458 The full list of not yet implemented functions is documented in the
459 file @ref include/smpi/smpi.h, between two lines containing the
460 <tt>FIXME</tt> marker. If you really need a missing feature, please
461 get in touch with us: we can guide you though the SimGrid code to help
462 you implementing it, and we'd glad to integrate it in the main project
463 afterward if you contribute them back.
464
465 @subsection SMPI_what_globals Global variables
466
467 Concerning the globals, the problem comes from the fact that usually,
468 MPI processes run as real UNIX processes while they are all folded
469 into threads of a unique system process in SMPI. Global variables are
470 usually private to each MPI process while they become shared between
471 the processes in SMPI. This point is rather problematic, and currently
472 forces to modify your application to privatize the global variables.
473
474 We tried several techniques to work this around. We used to have a
475 script that privatized automatically the globals through static
476 analysis of the source code, but it was not robust enough to be used
477 in production. This issue, as well as several potential solutions, is
478 discussed in this article: "Automatic Handling of Global Variables for
479 Multi-threaded MPI Programs",
480 available at http://charm.cs.illinois.edu/newPapers/11-23/paper.pdf
481 (note that this article does not deal with SMPI but with a competing
482 solution called AMPI that suffers of the same issue). 
483
484 SimGrid can duplicate and dynamically switch the .data and .bss
485 segments of the ELF process when switching the MPI ranks, allowing
486 each ranks to have its own copy of the global variables. This feature
487 is expected to work correctly on Linux and BSD, so smpirun activates
488 it by default. As no copy is involved, performance should not be
489 altered (but memory occupation will be higher).
490
491 If you want to turn it off, pass \c -no-privatize to smpirun. This may
492 be necessary if your application uses dynamic libraries as the global
493 variables of these libraries will not be privatized. You can fix this
494 by linking statically with these libraries (but NOT with libsimgrid,
495 as we need SimGrid's own global variables).
496
497 @section SMPI_adapting Adapting your MPI code for further scalability
498
499 As detailed in the reference article (available at
500 http://hal.inria.fr/inria-00527150), you may want to adapt your code
501 to improve the simulation performance. But these tricks may seriously
502 hinder the result quality (or even prevent the app to run) if used
503 wrongly. We assume that if you want to simulate an HPC application,
504 you know what you are doing. Don't prove us wrong!
505
506 @subsection SMPI_adapting_size Reducing your memory footprint
507
508 If you get short on memory (the whole app is executed on a single node when
509 simulated), you should have a look at the SMPI_SHARED_MALLOC and
510 SMPI_SHARED_FREE macros. It allows to share memory areas between processes: The
511 purpose of these macro is that the same line malloc on each process will point
512 to the exact same memory area. So if you have a malloc of 2M and you have 16
513 processes, this macro will change your memory consumption from 2M*16 to 2M
514 only. Only one block for all processes.
515
516 If your program is ok with a block containing garbage value because all
517 processes write and read to the same place without any kind of coordination,
518 then this macro can dramatically shrink your memory consumption. For example,
519 that will be very beneficial to a matrix multiplication code, as all blocks will
520 be stored on the same area. Of course, the resulting computations will useless,
521 but you can still study the application behavior this way. 
522
523 Naturally, this won't work if your code is data-dependent. For example, a Jacobi
524 iterative computation depends on the result computed by the code to detect
525 convergence conditions, so turning them into garbage by sharing the same memory
526 area between processes does not seem very wise. You cannot use the
527 SMPI_SHARED_MALLOC macro in this case, sorry.
528
529 This feature is demoed by the example file
530 <tt>examples/smpi/NAS/DT-folding/dt.c</tt>
531
532 @subsection SMPI_adapting_speed Toward faster simulations
533
534 If your application is too slow, try using SMPI_SAMPLE_LOCAL,
535 SMPI_SAMPLE_GLOBAL and friends to indicate which computation loops can
536 be sampled. Some of the loop iterations will be executed to measure
537 their duration, and this duration will be used for the subsequent
538 iterations. These samples are done per processor with
539 SMPI_SAMPLE_LOCAL, and shared between all processors with
540 SMPI_SAMPLE_GLOBAL. Of course, none of this will work if the execution
541 time of your loop iteration are not stable.
542
543 This feature is demoed by the example file 
544 <tt>examples/smpi/NAS/EP-sampling/ep.c</tt>
545
546 @section SMPI_accuracy Ensuring accurate simulations
547
548 Out of the box, SimGrid may give you fairly accurate results, but
549 there is a plenty of factors that could go wrong and make your results
550 inaccurate or even plainly wrong. Actually, you can only get accurate
551 results of a nicely built model, including both the system hardware
552 and your application. Such models are hard to pass over and reuse in
553 other settings, because elements that are not relevant to an
554 application (say, the latency of point-to-point communications,
555 collective operation implementation details or CPU-network
556 interaction) may be irrelevant to another application. The dream of
557 the perfect model, encompassing every aspects is only a chimera, as
558 the only perfect model of the reality is the reality. If you go for
559 simulation, then you have to ignore some irrelevant aspects of the
560 reality, but which aspects are irrelevant is actually
561 application-dependent...
562
563 The only way to assess whether your settings provide accurate results
564 is to double-check these results. If possible, you should first run
565 the same experiment in simulation and in real life, gathering as much
566 information as you can. Try to understand the discrepancies in the
567 results that you observe between both settings (visualization can be
568 precious for that). Then, try to modify your model (of the platform,
569 of the collective operations) to reduce the most preeminent differences.
570
571 If the discrepancies come from the computing time, try adapting the \c
572 smpi/host-speed: reduce it if your simulation runs faster than in
573 reality. If the error come from the communication, then you need to
574 fiddle with your platform file.
575
576 Be inventive in your modeling. Don't be afraid if the names given by
577 SimGrid does not match the real names: we got very good results by
578 modeling multicore/GPU machines with a set of separate hosts
579 interconnected with very fast networks (but don't trust your model
580 because it has the right names in the right place either).
581
582 Finally, you may want to check [this
583 article](https://hal.inria.fr/hal-00907887) on the classical pitfalls
584 in modeling distributed systems.
585
586 */
587
588
589 /** @example include/smpi/smpi.h */