Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
add additional necessary packages for pajengr installation
[simgrid.git] / docs / source / Tutorial_MPI_Applications.rst
1 .. _usecase_smpi:
2
3 Simulating MPI Applications
4 ===========================
5
6 Discover SMPI
7 -------------
8
9 SimGrid can not only :ref:`simulate algorithms <usecase_simalgo>`, but
10 it can also be used to execute real MPI applications on top of
11 virtual, simulated platforms with the SMPI module. Even complex
12 C/C++/F77/F90 applications should run out of the box in this
13 environment. Almost all proxy apps provided by the `ExaScale
14 Project <https://proxyapps.exascaleproject.org/>`_ only require minor
15 modifications to `run on top of SMPI
16 <https://framagit.org/simgrid/SMPI-proxy-apps>`_.
17
18 This setting permits one to debug your MPI applications in a perfectly
19 reproducible setup, with no Heisenbugs. Enjoy the full Clairvoyance
20 provided by the simulator while running what-if analyses on platforms
21 that are still to be built! Several `production-grade MPI applications
22 <https://framagit.org/simgrid/SMPI-proxy-apps#full-scale-applications>`_
23 use SimGrid for their integration and performance testing.
24
25 MPI 2.2 is already partially covered: over 160 primitives are
26 supported. Some parts of the standard are still missing: MPI-IO, MPI3
27 collectives, spawning ranks, inter-communicators, and some others. If
28 one of the functions you use is still missing, please drop us an
29 email. We may find the time to implement it for you.
30
31 Multi-threading support is very limited in SMPI. Only funneled
32 applications are supported: at most one thread per rank can issue any
33 MPI calls. For better timing predictions, your application should even
34 be completely mono-threaded. Using OpenMP (or pthreads directly) may
35 greatly decrease SimGrid predictive power. That may still be OK if you
36 only plan to debug your application in a reproducible setup, without
37 any performance-related analysis.
38
39 How does it work?
40 .................
41
42 In SMPI, communications are simulated while computations are
43 emulated. This means that while computations occur as they would in
44 the real systems, communication calls are intercepted and achieved by
45 the simulator.
46
47 To start using SMPI, you just need to compile your application with
48 ``smpicc`` instead of ``mpicc``, or with ``smpiff`` instead of
49 ``mpiff``, or with ``smpicxx`` instead of ``mpicxx``. Then, the only
50 difference between the classical ``mpirun`` and the new ``smpirun`` is
51 that it requires a new parameter ``-platform`` with a file describing
52 the simulated platform on which your application shall run.
53
54 Internally, all ranks of your application are executed as threads of a
55 single unix process. That's not a problem if your application has
56 global variables, because ``smpirun`` loads one application instance
57 per MPI rank as if it was another dynamic library. Then, MPI
58 communication calls are implemented using SimGrid: data is exchanged
59 through memory copy, while the simulator's performance models are used
60 to predict the time taken by each communication. Any computations
61 occurring between two MPI calls are benchmarked, and the corresponding
62 time is reported into the simulator.
63
64 .. image:: /tuto_smpi/img/big-picture.svg
65    :align: center
66
67 Describing Your Platform
68 ------------------------
69
70 As an SMPI user, you are supposed to provide a description of your
71 simulated platform, which is mostly a set of simulated hosts and network
72 links with some performance characteristics. SimGrid provides plenty
73 of :ref:`documentation <platform>` and examples (in the
74 `examples/platforms <https://framagit.org/simgrid/simgrid/tree/master/examples/platforms>`_
75 source directory), and this section only shows a small set of introductory
76 examples.
77
78 Feel free to skip this section if you want to jump right away to usage
79 examples.
80
81 Simple Example with 3 hosts
82 ...........................
83
84 Imagine you want to describe a little platform with three hosts,
85 interconnected as follows:
86
87 .. image:: /tuto_smpi/3hosts.png
88    :align: center
89
90 This can be done with the following platform file, which considers the
91 simulated platform as a graph of hosts and network links.
92
93 .. literalinclude:: /tuto_smpi/3hosts.xml
94    :language: xml
95
96 The elements basic elements (with :ref:`pf_tag_host` and
97 :ref:`pf_tag_link`) are described first, and then the routes between
98 any pair of hosts are explicitly given with :ref:`pf_tag_route`. 
99
100 Any host must be given a computational speed in flops while links must
101 be given a latency and a bandwidth. You can write 1Gf for
102 1,000,000,000 flops (full list of units in the reference guide of 
103 :ref:`pf_tag_host` and :ref:`pf_tag_link`). 
104
105 Routes defined with :ref:`pf_tag_route` are symmetrical by default,
106 meaning that the list of traversed links from A to B is the same as
107 from B to A. Explicitly define non-symmetrical routes if you prefer.
108
109 Cluster with a Crossbar
110 .......................
111
112 A very common parallel computing platform is a homogeneous cluster in
113 which hosts are interconnected via a crossbar switch with as many
114 ports as hosts so that any disjoint pairs of hosts can communicate
115 concurrently at full speed. For instance:
116
117 .. literalinclude:: ../../examples/platforms/cluster_crossbar.xml
118    :language: xml
119    :lines: 1-3,18-
120
121 One specifies a name prefix and suffix for each host and then gives an
122 integer range. In the example, the cluster contains 65535 hosts (!),
123 named ``node-0.simgrid.org`` to ``node-65534.simgrid.org``. All hosts
124 have the same power (1 Gflop/sec) and are connected to the switch via
125 links with the same bandwidth (125 MBytes/sec) and latency (50
126 microseconds).
127
128 .. todo::
129
130    Add the picture.
131
132 Cluster with a Shared Backbone
133 ..............................
134
135 Another popular model for a parallel platform is that of a set of
136 homogeneous hosts connected to a shared communication medium, a
137 backbone, with some finite bandwidth capacity, and on which
138 communicating host pairs can experience contention. For instance:
139
140
141 .. literalinclude:: ../../examples/platforms/cluster_backbone.xml
142    :language: xml
143    :lines: 1-3,18-
144
145 The only differences with the crossbar cluster above are the ``bb_bw``
146 and ``bb_lat`` attributes that specify the backbone characteristics
147 (here, a 500 microseconds latency and a 2.25 GByte/sec
148 bandwidth). This link is used for every communication within the
149 cluster. The route from ``node-0.simgrid.org`` to ``node-1.simgrid.org``
150 counts 3 links: the private link of ``node-0.simgrid.org``, the backbone,
151 and the private link of ``node-1.simgrid.org``.
152
153 .. todo::
154
155    Add the picture.
156
157 Torus Cluster
158 .............
159
160 Many HPC facilities use torus clusters to reduce sharing and
161 performance loss on concurrent internal communications. Modeling this
162 in SimGrid is very easy. Simply add a ``topology="TORUS"`` attribute
163 to your cluster. Configure it with the ``topo_parameters="X,Y,Z"``
164 attribute, where ``X``, ``Y``, and ``Z`` are the dimensions of your
165 torus.
166
167 .. image:: ../../examples/platforms/cluster_torus.svg
168    :align: center
169
170 .. literalinclude:: ../../examples/platforms/cluster_torus.xml
171    :language: xml
172
173 Note that in this example, we used ``loopback_bw`` and
174 ``loopback_lat`` to specify the characteristics of the loopback link
175 of each node (i.e., the link allowing each node to communicate with
176 itself). We could have done so in the previous example too. When no
177 loopback is given, the communication from a node to itself is handled
178 as if it were two distinct nodes: it goes twice through the private
179 link and through the backbone (if any).
180
181 Fat-Tree Cluster
182 ................
183
184 This topology was introduced to reduce the number of links in the
185 cluster (and thus reduce its price) while maintaining a high bisection
186 bandwidth and a relatively low diameter. To model this in SimGrid,
187 pass a ``topology="FAT_TREE"`` attribute to your cluster. The
188 ``topo_parameters=#levels;#downlinks;#uplinks;link count`` follows the
189 semantic introduced in `Figure 1B of this article
190 <http://webee.eedev.technion.ac.il/wp-content/uploads/2014/08/publication_574.pdf>`_.
191
192 Here is the meaning of this example: ``2 ; 4,4 ; 1,2 ; 1,2``
193
194 - That's a two-level cluster (thus the initial ``2``).
195 - Routers are connected to 4 elements below them, regardless of their
196   level. Thus the ``4,4`` component used as
197   ``#downlinks``. This means that the hosts are grouped by 4 on a
198   given router, and that there are 4 level-1 routers (in the middle of
199   the figure).
200 - Hosts are connected to only 1 router above them, while these routers
201   are connected to 2 routers above them (thus the ``1,2`` used as
202   ``#uplink``).
203 - Hosts have only one link to their router while every path between 
204   level-1 routers and level-2 routers uses 2 parallel links. Thus the
205   ``1,2`` used as ``link count``.
206
207 .. image:: ../../examples/platforms/cluster_fat_tree.svg
208    :align: center
209
210 .. literalinclude:: ../../examples/platforms/cluster_fat_tree.xml
211    :language: xml
212    :lines: 1-3,10-
213
214
215 Dragonfly Cluster
216 .................
217
218 This topology was introduced to further reduce the number of links
219 while maintaining a high bandwidth for local communications. To model
220 this in SimGrid, pass a ``topology="DRAGONFLY"`` attribute to your
221 cluster. It's based on the implementation of the topology used on
222 Cray XC systems, described in the paper
223 `Cray Cascade: A scalable HPC system based on a Dragonfly network <https://dl.acm.org/citation.cfm?id=2389136>`_.
224
225 System description follows the format ``topo_parameters=#groups;#chassis;#routers;#nodes``
226 For example, ``3,4 ; 3,2 ; 3,1 ; 2``:
227
228 - ``3,4``: There are 3 groups with 4 links between each (blue level).
229   Links to the nth group are attached to the nth router of the group
230   on our implementation.
231 - ``3,2``: In each group, there are 3 chassis with 2 links between each nth router
232   of each group (black level)
233 - ``3,1``: In each chassis, 3 routers are connected with a single link
234   (green level)
235 - ``2``: Each router has two nodes attached (single link)
236
237 .. image:: ../../examples/platforms/cluster_dragonfly.svg
238    :align: center
239
240 .. literalinclude:: ../../examples/platforms/cluster_dragonfly.xml
241    :language: xml
242
243 Final Word
244 ..........
245
246 We only glanced over the abilities offered by SimGrid to describe the
247 platform topology. Other networking zones model non-HPC platforms
248 (such as wide area networks, ISP networks comprising set-top boxes, or
249 even your own routing schema). You can interconnect several networking
250 zones in your platform to form a tree of zones, that is both a time-
251 and memory-efficient representation of distributed platforms. Please
252 head to the dedicated :ref:`documentation <platform>` for more
253 information.
254
255 Hands-on!
256 ---------
257
258 It is time to start using SMPI yourself. For that, you first need to
259 install it somehow, and then you will need an MPI application to play with.
260
261 Using Docker
262 ............
263
264 The easiest way to take the tutorial is to use the dedicated Docker
265 image. Once you `installed Docker itself
266 <https://docs.docker.com/install/>`_, simply do the following:
267
268 .. code-block:: shell
269
270    docker pull simgrid/tuto-smpi
271    docker run -it --rm --name simgrid --volume ~/smpi-tutorial:/source/tutorial simgrid/tuto-smpi bash
272
273 This will start a new container with all you need to take this
274 tutorial, and create a ``smpi-tutorial`` directory in your home on
275 your host machine that will be visible as ``/source/tutorial`` within the
276 container.  You can then edit the files you want with your favorite
277 editor in ``~/smpi-tutorial``, and compile them within the
278 container to enjoy the provided dependencies.
279
280 .. warning::
281
282    Any change to the container out of ``/source/tutorial`` will be lost
283    when you log out of the container, so don't edit the other files!
284
285 All needed dependencies are already installed in this container
286 (SimGrid, the C/C++/Fortran compilers, make, pajeng, R and pajengr). Vite being
287 only optional in this tutorial, it is not installed to reduce the
288 image size.
289
290 The container also includes the example platform files from the
291 previous section as well as the source code of the NAS Parallel
292 Benchmarks. These files are available under
293 ``/source/simgrid-template-smpi`` in the image. You should copy it to
294 your working directory when you first log in:
295
296 .. code-block:: shell
297
298    cp -r /source/simgrid-template-smpi/* /source/tutorial
299    cd /source/tutorial
300
301 Using your Computer Natively
302 ............................
303
304 To take the tutorial on your machine, you first need to :ref:`install
305 SimGrid <install>`, the C/C++/Fortran compilers, and also ``pajeng`` to
306 visualize the traces. You may want to install `Vite
307 <http://vite.gforge.inria.fr/>`_ to get a first glance at the
308 traces. The provided code template requires ``make`` to compile. On
309 Debian and Ubuntu, you can get them as follows:
310
311 .. code-block:: shell
312
313    sudo apt install simgrid pajeng make gcc g++ gfortran python3 vite
314
315 For R analysis of the produced traces, you may want to install R
316 and the `pajengr <https://github.com/schnorr/pajengr#installation/>`_ package.
317
318 .. code-block:: shell
319
320    # install R and necessary packages
321    sudo apt install r-base r-cran-devtools r-cran-tidyverse
322    # install pajengr dependencies
323    sudo apt install git cmake flex bison
324    # install the pajengr R package
325    Rscript -e "library(devtools); install_github('schnorr/pajengr');"
326
327 To take this tutorial, you will also need the platform files from the
328 previous section as well as the source code of the NAS Parallel
329 Benchmarks. Just  clone `this repository
330 <https://framagit.org/simgrid/simgrid-template-smpi>`_  to get them all:
331
332 .. code-block:: shell
333
334    git clone https://framagit.org/simgrid/simgrid-template-smpi.git
335    cd simgrid-template-smpi/
336
337 If you struggle with the compilation, then you should double-check
338 your :ref:`SimGrid installation <install>`.  On need, please refer to
339 the :ref:`Troubleshooting your Project Setup
340 <install_yours_troubleshooting>` section.
341
342 Lab 0: Hello World
343 ------------------
344
345 It is time to simulate your first MPI program. Use the simplistic
346 example `roundtrip.c
347 <https://framagit.org/simgrid/simgrid-template-smpi/raw/master/roundtrip.c?inline=false>`_
348 that comes with the template.
349
350 .. literalinclude:: /tuto_smpi/roundtrip.c
351    :language: c
352
353 Compiling and Executing
354 .......................
355
356 Compiling the program is straightforward (double-check your
357 :ref:`SimGrid installation <install>` if you get an error message):
358
359
360 .. code-block:: shell
361
362   $ smpicc -O3 roundtrip.c -o roundtrip
363
364
365 Once compiled, you can simulate the execution of this program on 16
366 nodes from the ``cluster_crossbar.xml`` platform as follows:
367
368 .. code-block:: shell
369
370    $ smpirun -np 16 -platform cluster_crossbar.xml -hostfile cluster_hostfile ./roundtrip
371
372 - The ``-np 16`` option, just like in regular MPI, specifies the
373   number of MPI processes to use.
374 - The ``-hostfile cluster_hostfile`` option, just like in regular
375   MPI, specifies the host file. If you omit this option, ``smpirun``
376   will deploy the application on the first machines of your platform.
377 - The ``-platform cluster_crossbar.xml`` option, **which doesn't exist
378   in regular MPI**, specifies the platform configuration to be
379   simulated.
380 - At the end of the line, one finds the executable name and
381   command-line arguments (if any -- roundtrip does not expect any arguments).
382
383 Feel free to tweak the content of the XML platform file and the
384 program to see the effect on the simulated execution time. It may be
385 easier to compare the executions with the extra option
386 ``--cfg=smpi/display-timing:yes``.  Note that the simulation accounts
387 for realistic network protocol effects and MPI implementation
388 effects. As a result, you may see "unexpected behavior" like in the
389 real world (e.g., sending a message 1 byte larger may lead to
390 significantly higher execution time).
391
392 Lab 1: Visualizing LU
393 ---------------------
394
395 We will now simulate a larger application: the LU benchmark of the NAS
396 suite. The version provided in the code template was modified to
397 compile with SMPI instead of the regular MPI. Compare the difference
398 between the original ``config/make.def.template`` and the
399 ``config/make.def`` that was adapted to SMPI. We use ``smpiff`` and
400 ``smpicc`` as compilers, and don't pass any additional library.
401
402 Now compile and execute the LU benchmark, class S (i.e., for `small
403 data size
404 <https://www.nas.nasa.gov/publications/npb_problem_sizes.html>`_) with
405 4 nodes.
406
407 .. code-block:: shell
408
409    $ make lu NPROCS=4 CLASS=S
410    (compilation logs)
411    $ smpirun -np 4 -platform ../cluster_backbone.xml bin/lu.S.4
412    (execution logs)
413
414 To get a better understanding of what is going on, activate the
415 visualization tracing, and convert the produced trace for later
416 use:
417
418 .. code-block:: shell
419
420    smpirun -np 4 -platform ../cluster_backbone.xml -trace --cfg=tracing/filename:lu.S.4.trace bin/lu.S.4
421
422 You can then produce a Gantt Chart with the following R chunk. You can
423 either copy/paste it in an R session, or `turn it into a Rscript executable
424 <https://swcarpentry.github.io/r-novice-inflammation/05-cmdline/>`_ to
425 run it again and again.
426
427 .. code-block:: R
428
429    # Read the data
430    library(tidyverse)
431    library(pajengr)
432    dta <- pajeng_read("lu.S.4.trace")
433
434    # Manipulate the data
435    dta$state %>%
436       # Remove some unnecessary columns for this example
437       select(-Type, -Imbrication) %>%
438       # Create the nice MPI rank and operations identifiers
439       mutate(Container = as.integer(gsub("rank-", "", Container)),
440              Value = gsub("^PMPI_", "MPI_", Value)) %>%
441       # Rename some columns so it can better fit MPI terminology
442       rename(Rank = Container,
443              Operation = Value) -> df.states
444
445    # Draw the Gantt Chart
446    df.states %>%
447       ggplot() +
448       # Each MPI operation is becoming a rectangle
449       geom_rect(aes(xmin=Start, xmax=End,
450                     ymin=Rank,  ymax=Rank + 1,
451                     fill=Operation)) +
452       # Cosmetics
453       xlab("Time [seconds]") +
454       ylab("Rank [count]") +
455       theme_bw(base_size=14) +
456       theme(
457         plot.margin = unit(c(0,0,0,0), "cm"),
458         legend.margin = margin(t = 0, unit='cm'),
459         panel.grid = element_blank(),
460         legend.position = "top",
461         legend.justification = "left",
462         legend.box.spacing = unit(0, "pt"),
463         legend.box.margin = margin(0,0,0,0),
464         legend.title = element_text(size=10)) -> plot
465
466     # Save the plot in a PNG file (dimensions in inches)
467     ggsave("smpi.png",
468            plot,
469            width = 10,
470            height = 3)
471
472 This produces a file called ``smpi.png`` with the following
473 content. You can find more visualization examples `online
474 <https://simgrid.org/contrib/R_visualization.html>`_.
475
476 .. image:: /tuto_smpi/img/lu.S.4.png
477    :align: center
478
479 Lab 2: Tracing and Replay of LU
480 -------------------------------
481
482 Now compile and execute the LU benchmark, class A, with 32 nodes.
483
484 .. code-block:: shell
485
486    $ make lu NPROCS=32 CLASS=A
487
488 This takes several minutes to simulate, because all code from all
489 processes has to be actually executed, and everything is serialized.
490
491 SMPI provides several methods to speed things up. One of them is to
492 capture a time-independent trace of the running application and
493 replay it on a different platform with the same amount of nodes. The
494 replay is much faster than live simulation, as the computations are
495 skipped (the application must be network-dependent for this to work).
496
497 You can even generate the trace during the live simulation as follows:
498
499 .. code-block:: shell
500
501    $ smpirun -trace-ti --cfg=tracing/filename:LU.A.32 -np 32 -platform ../cluster_backbone.xml bin/lu.A.32
502
503 The produced trace is composed of a file ``LU.A.32`` and a folder
504 ``LU.A.32_files``. You can replay this trace with SMPI thanks to ``smpirun``.
505 For example, the following command replays the trace on a different platform:
506
507 .. code-block:: shell
508
509    $ smpirun -np 32 -platform ../cluster_crossbar.xml -hostfile ../cluster_hostfile -replay LU.A.32
510
511 All the outputs are gone, as the application is not really simulated
512 here. Its trace is simply replayed. But if you visualize the live
513 simulation and the replay, you will see that the behavior is
514 unchanged. The simulation does not run much faster on this very
515 example, but this becomes very interesting when your application
516 is computationally hungry.
517
518 .. todo::
519
520     The commands should be separated and executed by some CI to make sure
521     the documentation is up-to-date.
522
523 Lab 3: Execution Sampling on Matrix Multiplication example
524 ----------------------------------------------------------
525
526 The second method to speed up simulations is to sample the computation
527 parts in the code.  This means that the person doing the simulation
528 needs to know the application and identify parts that are compute-intensive
529 and take time while being regular enough not to ruin
530 simulation accuracy. Furthermore, there should not be any MPI calls
531 inside such parts of the code.
532
533 Use for this part the `gemm_mpi.cpp
534 <https://gitlab.com/PRACE-4IP/CodeVault/raw/master/hpc_kernel_samples/dense_linear_algebra/gemm/mpi/src/gemm_mpi.cpp>`_
535 example, which is provided by the `PRACE Codevault repository
536 <http://www.prace-ri.eu/prace-codevault/>`_.
537
538 The computing part of this example is the matrix multiplication routine
539
540 .. literalinclude:: /tuto_smpi/gemm_mpi.cpp
541    :language: cpp
542    :lines: 4-19
543
544 .. code-block:: shell
545
546   $ smpicxx -O3 gemm_mpi.cpp -o gemm
547   $ time smpirun -np 16 -platform cluster_crossbar.xml -hostfile cluster_hostfile --cfg=smpi/display-timing:yes --cfg=smpi/running-power:1000000000 ./gemm
548
549 This should end quite quickly, as the size of each matrix is only 1000x1000. 
550 But what happens if we want to simulate larger runs?
551 Replace the size by 2000, 3000, and try again.
552
553 The simulation time increases a lot, while there are no more MPI calls performed, only computation.
554
555 The ``--cfg=smpi/display-timing`` option gives more details about execution 
556 and advises using sampling if the time spent in computing loops seems too high.
557
558 The ``--cfg=smpi/host-speed:1000000000`` option sets the speed of the processor used for 
559 running the simulation. Here we say that its speed is the same as one of the 
560 processors we are simulating (1Gf), so that 1 second of computation is injected 
561 as 1 second in the simulation.
562
563 .. code-block:: shell
564
565   [5.568556] [smpi_kernel/INFO] Simulated time: 5.56856 seconds.
566
567   The simulation took 24.9403 seconds (after parsing and platform setup)
568   24.0764 seconds were actual computation of the application
569   [5.568556] [smpi_kernel/INFO] More than 75% of the time was spent inside the application code.
570   You may want to use sampling functions or trace replay to reduce this.
571
572 So in our case (size 3000), the simulation ran for 25 seconds, and the simulated time was 5.57s at the end.
573 Computation by itself took 24 seconds, and can quickly grow with larger sizes 
574 (as computation is really performed, there will be variability between similar runs).
575
576 SMPI provides sampling macros to accelerate simulation by sampling iterations 
577 of large computation loops, and skip computation after a certain amount of iterations, 
578 or when the sampling is stable enough.
579
580 The two macros only slightly differ :
581
582 - ``SMPI_SAMPLE_GLOBAL`` : the specified number of samples is produced by all processors
583 - ``SMPI_SAMPLE_LOCAL`` : each process executes a specified number of iterations
584
585 So if the size of the computed part varies between processes (imbalance), 
586 it's safer to use the LOCAL one.
587
588 To use one of them, replacing the external for loop of the multiply routine:
589
590 .. code-block:: c
591
592   for (int i = istart; i <= iend; ++i)
593
594 by:
595
596 .. code-block:: c
597
598   SMPI_SAMPLE_GLOBAL(int i = istart, i <= iend, ++i, 10, 0.005)
599
600 The first three parameters are the ones from the loop, while the two last ones are for sampling.
601 They mean that at most 10 iterations will be performed and that the sampling phase can be exited 
602 earlier if the expected stability is reached after fewer samples.
603
604 Now run the code again with various sizes and parameters and check the time taken for the 
605 simulation, as well as the resulting simulated time.
606
607 .. code-block:: shell
608
609   [5.575691] [smpi_kernel/INFO] Simulated time: 5.57569 seconds.
610   The simulation took 1.23698 seconds (after parsing and platform setup)
611   0.0319454 seconds were actual computation of the application
612
613 In this case, the simulation only took 1.2 seconds, while the simulated time 
614 remained almost identical.
615
616 The computation results will obviously be altered since most computations are skipped.
617 These macros thus cannot be used when results are critical for the application behavior
618 (convergence estimation for instance will be wrong on some codes).
619
620
621 Lab 4: Memory folding on large allocations
622 ------------------------------------------
623
624 Another issue that can be encountered when simulation with SMPI is lack of memory.
625 Indeed we are executing all MPI processes on a single node, which can lead to crashes.
626 We will use the DT benchmark of the NAS suite to illustrate how to avoid such issues.
627
628 With 85 processes and class C, the DT simulated benchmark will try to allocate 35GB of memory,
629 which may not be available on the node you are using.
630
631 To avoid this we can simply replace the largest calls to malloc and free by calls 
632 to ``SMPI_SHARED_MALLOC`` and ``SMPI_SHARED_FREE``.
633 This means that all processes will share one single instance of this buffer.
634 As for sampling, results will be altered, and this should not be used for control structures.
635
636 For DT example, there are three different calls to malloc in the file, and one of them is for a needed structure.
637 Find it and replace the two other ones with ``SMPI_SHARED_MALLOC`` (there is only one free to replace for both of them).
638
639 Once done, you can now run
640
641 .. code-block:: shell
642
643    $ make dt NPROCS=85 CLASS=C
644    (compilation logs)
645    $ smpirun -np 85 -platform ../cluster_backbone.xml bin/dt.C.x BH
646    (execution logs)
647
648 And simulation should finish without swapping/crashing (Ignore the warning about the return value).
649
650 If control structures are also problematic, you can use ``SMPI_PARTIAL_SHARED_MALLOC(size, offsets, offsetscount)`` 
651 macro, which shares only specific parts of the structure between processes, 
652 and use specific memory for the important parts.
653 It can be freed afterward with SMPI_SHARED_FREE.
654
655 If allocations are performed with malloc or calloc, SMPI (from version 3.25) provides the option
656 ``--cfg=smpi/auto-shared-malloc-shared:n`` which will replace all allocations above size n bytes by
657 shared allocations. The value has to be carefully selected to avoid smaller control arrays, 
658 containing data necessary for the completion of the run.
659 Try to run the (non modified) DT example again, with values going from 10 to 100,000 to show that
660 too small values can cause crashes.
661
662 A useful option to identify the largest allocations in the code is ``--cfg=smpi/display-allocs:yes`` (from 3.27).
663 It will display at the end of a (successful) run the largest allocations and their locations, helping pinpoint the
664 targets for sharing, or setting the threshold for automatic ones.
665 For DT, the process would be to run a smaller class of problems, 
666
667 .. code-block:: shell
668
669    $ make dt NPROCS=21 CLASS=A
670    $ smpirun --cfg=smpi/display-allocs:yes -np 21 -platform ../cluster_backbone.xml bin/dt.A.x BH
671
672 Which should output:
673
674 .. code-block:: shell
675
676     [smpi_utils/INFO] Memory Usage: Simulated application allocated 198533192 bytes during its lifetime through malloc/calloc calls.
677     Largest allocation at once from a single process was 3553184 bytes, at dt.c:388. It was called 3 times during the whole simulation.
678     If this is too much, consider sharing allocations for computation buffers.
679     This can be done automatically by setting --cfg=smpi/auto-shared-malloc-thresh to the minimum size wanted size (this can alter execution if data content is necessary)
680
681 And from there we can identify dt.c:388 as the main allocation, and the best target to convert to
682 shared mallocs for larger simulations. Furthermore, with 21 processes, we see that this particular
683 allocation and size was only called 3 times, which means that other processes are likely to allocate
684 less memory here (imbalance). Using 3553184 as a threshold value might be unwise, as most processes
685 wouldn't share memory, so a lower threshold would be advisable.
686
687 Further Readings
688 ----------------
689
690 You may also be interested in the `SMPI reference article
691 <https://hal.inria.fr/hal-01415484>`_ or these `introductory slides
692 <http://simgrid.org/tutorials/simgrid-smpi-101.pdf>`_. The :ref:`SMPI
693 reference documentation <SMPI_doc>` covers much more content than
694 this short tutorial.
695
696 Finally, we regularly use SimGrid in our teachings on MPI. This way,
697 our student can experiment with platforms that they do not have access
698 to, and the associated visualization tools help them to understand
699 their work.  The whole material is available online, in a separate
700 project: the `SMPI CourseWare <https://simgrid.github.io/SMPI_CourseWare/>`_.
701
702 ..  LocalWords:  SimGrid