297 lines
35 KiB
HTML
297 lines
35 KiB
HTML
<!DOCTYPE html PUBLIC "-//W3C//DTD XHTML 1.0 Transitional//EN" "https://www.w3.org/TR/xhtml1/DTD/xhtml1-transitional.dtd">
|
|
<html xmlns="http://www.w3.org/1999/xhtml" lang="en-US">
|
|
<head>
|
|
<meta http-equiv="Content-Type" content="text/xhtml;charset=UTF-8"/>
|
|
<meta http-equiv="X-UA-Compatible" content="IE=11"/>
|
|
<meta name="generator" content="Doxygen 1.9.6"/>
|
|
<meta name="viewport" content="width=device-width, initial-scale=1"/>
|
|
<title>2 Dimensional Ising Model: src/phase_transition_mpi.cpp Source File</title>
|
|
<link href="tabs.css" rel="stylesheet" type="text/css"/>
|
|
<script type="text/javascript" src="jquery.js"></script>
|
|
<script type="text/javascript" src="dynsections.js"></script>
|
|
<link href="navtree.css" rel="stylesheet" type="text/css"/>
|
|
<script type="text/javascript" src="resize.js"></script>
|
|
<script type="text/javascript" src="navtreedata.js"></script>
|
|
<script type="text/javascript" src="navtree.js"></script>
|
|
<link href="search/search.css" rel="stylesheet" type="text/css"/>
|
|
<script type="text/javascript" src="search/searchdata.js"></script>
|
|
<script type="text/javascript" src="search/search.js"></script>
|
|
<script type="text/x-mathjax-config">
|
|
MathJax.Hub.Config({
|
|
extensions: ["tex2jax.js"],
|
|
jax: ["input/TeX","output/HTML-CSS"],
|
|
});
|
|
</script>
|
|
<script type="text/javascript" async="async" src="https://cdn.jsdelivr.net/npm/mathjax@2/MathJax.js"></script>
|
|
<link href="doxygen.css" rel="stylesheet" type="text/css" />
|
|
<link href="doxygen-awesome.css" rel="stylesheet" type="text/css"/>
|
|
</head>
|
|
<body>
|
|
<div id="top"><!-- do not remove this div, it is closed by doxygen! -->
|
|
<div id="titlearea">
|
|
<table cellspacing="0" cellpadding="0">
|
|
<tbody>
|
|
<tr id="projectrow">
|
|
<td id="projectalign">
|
|
<div id="projectname">2 Dimensional Ising Model
|
|
</div>
|
|
<div id="projectbrief">Simulate the change in energy and magnetization in a ferro magnet</div>
|
|
</td>
|
|
</tr>
|
|
</tbody>
|
|
</table>
|
|
</div>
|
|
<!-- end header part -->
|
|
<!-- Generated by Doxygen 1.9.6 -->
|
|
<script type="text/javascript">
|
|
/* @license magnet:?xt=urn:btih:d3d9a9a6595521f9666a5e94cc830dab83b65699&dn=expat.txt MIT */
|
|
var searchBox = new SearchBox("searchBox", "search/",'.html');
|
|
/* @license-end */
|
|
</script>
|
|
<script type="text/javascript" src="menudata.js"></script>
|
|
<script type="text/javascript" src="menu.js"></script>
|
|
<script type="text/javascript">
|
|
/* @license magnet:?xt=urn:btih:d3d9a9a6595521f9666a5e94cc830dab83b65699&dn=expat.txt MIT */
|
|
$(function() {
|
|
initMenu('',true,false,'search.php','Search');
|
|
$(document).ready(function() { init_search(); });
|
|
});
|
|
/* @license-end */
|
|
</script>
|
|
<div id="main-nav"></div>
|
|
</div><!-- top -->
|
|
<div id="side-nav" class="ui-resizable side-nav-resizable">
|
|
<div id="nav-tree">
|
|
<div id="nav-tree-contents">
|
|
<div id="nav-sync" class="sync"></div>
|
|
</div>
|
|
</div>
|
|
<div id="splitbar" style="-moz-user-select:none;"
|
|
class="ui-resizable-handle">
|
|
</div>
|
|
</div>
|
|
<script type="text/javascript">
|
|
/* @license magnet:?xt=urn:btih:d3d9a9a6595521f9666a5e94cc830dab83b65699&dn=expat.txt MIT */
|
|
$(document).ready(function(){initNavTree('phase__transition__mpi_8cpp_source.html',''); initResizable(); });
|
|
/* @license-end */
|
|
</script>
|
|
<div id="doc-content">
|
|
<!-- window showing the filter options -->
|
|
<div id="MSearchSelectWindow"
|
|
onmouseover="return searchBox.OnSearchSelectShow()"
|
|
onmouseout="return searchBox.OnSearchSelectHide()"
|
|
onkeydown="return searchBox.OnSearchSelectKey(event)">
|
|
</div>
|
|
|
|
<!-- iframe showing the search results (closed by default) -->
|
|
<div id="MSearchResultsWindow">
|
|
<div id="MSearchResults">
|
|
<div class="SRPage">
|
|
<div id="SRIndex">
|
|
<div id="SRResults"></div>
|
|
<div class="SRStatus" id="Loading">Loading...</div>
|
|
<div class="SRStatus" id="Searching">Searching...</div>
|
|
<div class="SRStatus" id="NoMatches">No Matches</div>
|
|
</div>
|
|
</div>
|
|
</div>
|
|
</div>
|
|
|
|
<div class="header">
|
|
<div class="headertitle"><div class="title">phase_transition_mpi.cpp</div></div>
|
|
</div><!--header-->
|
|
<div class="contents">
|
|
<a href="phase__transition__mpi_8cpp.html">Go to the documentation of this file.</a><div class="fragment"><div class="line"><a id="l00001" name="l00001"></a><span class="lineno"> 1</span><span class="comment">/** @file phase_transition_mpi.cpp</span></div>
|
|
<div class="line"><span class="lineno"> 2</span><span class="comment"> *</span></div>
|
|
<div class="line"><span class="lineno"> 3</span><span class="comment"> * @author Cory Alexander Balaton (coryab)</span></div>
|
|
<div class="line"><span class="lineno"> 4</span><span class="comment"> * @author Janita Ovidie Sandtrøen Willumsen (janitaws)</span></div>
|
|
<div class="line"><span class="lineno"> 5</span><span class="comment"> *</span></div>
|
|
<div class="line"><span class="lineno"> 6</span><span class="comment"> * @version 1.0</span></div>
|
|
<div class="line"><span class="lineno"> 7</span><span class="comment"> *</span></div>
|
|
<div class="line"><span class="lineno"> 8</span><span class="comment"> * @brief Sweep over different temperatures and generate data.</span></div>
|
|
<div class="line"><span class="lineno"> 9</span><span class="comment"> *</span></div>
|
|
<div class="line"><span class="lineno"> 10</span><span class="comment"> * @details This program takes in 4 arguments: the start temperature,</span></div>
|
|
<div class="line"><span class="lineno"> 11</span><span class="comment"> * the end temperature, the amount of temperature points to simulate, and</span></div>
|
|
<div class="line"><span class="lineno"> 12</span><span class="comment"> * the amount of monte carlo samples to collect, in that order.</span></div>
|
|
<div class="line"><span class="lineno"> 13</span><span class="comment"> *</span></div>
|
|
<div class="line"><span class="lineno"> 14</span><span class="comment"> * @bug No known bugs</span></div>
|
|
<div class="line"><span class="lineno"> 15</span><span class="comment"> * */</span></div>
|
|
<div class="line"><a id="l00016" name="l00016"></a><span class="lineno"> 16</span><span class="preprocessor">#</span><span class="preprocessor">include</span> <a class="code" href="data__type_8hpp.html" title="Header for the data_t type.">"data_type.hpp"</a></div>
|
|
<div class="line"><a id="l00017" name="l00017"></a><span class="lineno"> 17</span><span class="preprocessor">#</span><span class="preprocessor">include</span> <a class="code" href="monte__carlo_8hpp.html" title="Functions for Monte Carlo simulations.">"monte_carlo.hpp"</a></div>
|
|
<div class="line"><a id="l00018" name="l00018"></a><span class="lineno"> 18</span><span class="preprocessor">#</span><span class="preprocessor">include</span> <a class="code" href="utils_8hpp.html" title="Function prototypes and macros that are useful.">"utils.hpp"</a></div>
|
|
<div class="line"><a id="l00019" name="l00019"></a><span class="lineno"> 19</span> </div>
|
|
<div class="line"><a id="l00020" name="l00020"></a><span class="lineno"> 20</span><span class="preprocessor">#</span><span class="preprocessor">include</span> <span class="preprocessor"><</span><span class="preprocessor">getopt</span><span class="preprocessor">.</span><span class="preprocessor">h</span><span class="preprocessor">></span></div>
|
|
<div class="line"><a id="l00021" name="l00021"></a><span class="lineno"> 21</span><span class="preprocessor">#</span><span class="preprocessor">include</span> <span class="preprocessor"><</span><span class="preprocessor">mpi</span><span class="preprocessor">.</span><span class="preprocessor">h</span><span class="preprocessor">></span></div>
|
|
<div class="line"><a id="l00022" name="l00022"></a><span class="lineno"> 22</span><span class="preprocessor">#</span><span class="preprocessor">include</span> <span class="preprocessor"><</span><span class="preprocessor">string</span><span class="preprocessor">></span></div>
|
|
<div class="line"><a id="l00023" name="l00023"></a><span class="lineno"> 23</span> </div>
|
|
<div class="line"><a id="l00024" name="l00024"></a><span class="lineno"> 24</span><span class="comment">/** @brief A function that displays how to use the program and quits.*/</span></div>
|
|
<div class="line"><a id="l00025" name="l00025"></a><span class="lineno"><a class="line" href="phase__transition__mpi_8cpp.html#ac907e18135856c90366aaa599a9e10b1"> 25</a></span><span class="keywordtype">void</span> <a class="code hl_function" href="main_8cpp.html#ac907e18135856c90366aaa599a9e10b1">usage</a>(std::string filename)</div>
|
|
<div class="line"><a id="l00026" name="l00026"></a><span class="lineno"> 26</span>{</div>
|
|
<div class="line"><a id="l00027" name="l00027"></a><span class="lineno"> 27</span> std::cout</div>
|
|
<div class="line"><a id="l00028" name="l00028"></a><span class="lineno"> 28</span> << <span class="stringliteral">"Usage: "</span> << filename</div>
|
|
<div class="line"><a id="l00029" name="l00029"></a><span class="lineno"> 29</span> << <span class="stringliteral">" <start temperature> <end temperature> <lattice size> "</span></div>
|
|
<div class="line"><a id="l00030" name="l00030"></a><span class="lineno"> 30</span> <span class="stringliteral">"<points> <cycles> <burn-in-time> <output file>\n"</span></div>
|
|
<div class="line"><a id="l00031" name="l00031"></a><span class="lineno"> 31</span> << <span class="stringliteral">"This should be used with mpiexec or mpirun for maximum "</span></div>
|
|
<div class="line"><a id="l00032" name="l00032"></a><span class="lineno"> 32</span> <span class="stringliteral">"performance\n\n"</span></div>
|
|
<div class="line"><a id="l00033" name="l00033"></a><span class="lineno"> 33</span> << <span class="stringliteral">"\t[ -h | --help ]\n"</span>;</div>
|
|
<div class="line"><a id="l00034" name="l00034"></a><span class="lineno"> 34</span> exit(-1);</div>
|
|
<div class="line"><a id="l00035" name="l00035"></a><span class="lineno"> 35</span>}</div>
|
|
<div class="line"><a id="l00036" name="l00036"></a><span class="lineno"> 36</span> </div>
|
|
<div class="line"><a id="l00037" name="l00037"></a><span class="lineno"> 37</span><span class="comment">/** @brief The main function.*/</span></div>
|
|
<div class="line"><a id="l00038" name="l00038"></a><span class="lineno"><a class="line" href="phase__transition__mpi_8cpp.html#a3c04138a5bfe5d72780bb7e82a18e627"> 38</a></span><span class="keywordtype">int</span> <a class="code hl_function" href="main_8cpp.html#a3c04138a5bfe5d72780bb7e82a18e627">main</a>(<span class="keywordtype">int</span> argc, <span class="keywordtype">char</span> **argv)</div>
|
|
<div class="line"><a id="l00039" name="l00039"></a><span class="lineno"> 39</span>{</div>
|
|
<div class="line"><a id="l00040" name="l00040"></a><span class="lineno"> 40</span> <span class="comment">// Command options</span></div>
|
|
<div class="line"><a id="l00041" name="l00041"></a><span class="lineno"> 41</span> <span class="keyword">struct</span> option long_options[] = {{<span class="stringliteral">"help"</span>, 0, 0, 0}, {NULL, 0, NULL, 0}};</div>
|
|
<div class="line"><a id="l00042" name="l00042"></a><span class="lineno"> 42</span> </div>
|
|
<div class="line"><a id="l00043" name="l00043"></a><span class="lineno"> 43</span> <span class="keywordtype">int</span> option_index = -1;</div>
|
|
<div class="line"><a id="l00044" name="l00044"></a><span class="lineno"> 44</span> <span class="keywordtype">int</span> c;</div>
|
|
<div class="line"><a id="l00045" name="l00045"></a><span class="lineno"> 45</span> </div>
|
|
<div class="line"><a id="l00046" name="l00046"></a><span class="lineno"> 46</span> <span class="keywordflow">while</span> (<span class="keyword">true</span>) {</div>
|
|
<div class="line"><a id="l00047" name="l00047"></a><span class="lineno"> 47</span> c = getopt_long(argc, argv, <span class="stringliteral">"h"</span>, long_options, &option_index);</div>
|
|
<div class="line"><a id="l00048" name="l00048"></a><span class="lineno"> 48</span> </div>
|
|
<div class="line"><a id="l00049" name="l00049"></a><span class="lineno"> 49</span> <span class="keywordflow">if</span> (c == -1)</div>
|
|
<div class="line"><a id="l00050" name="l00050"></a><span class="lineno"> 50</span> <span class="keywordflow">break</span>;</div>
|
|
<div class="line"><a id="l00051" name="l00051"></a><span class="lineno"> 51</span> </div>
|
|
<div class="line"><a id="l00052" name="l00052"></a><span class="lineno"> 52</span> <span class="keywordflow">switch</span> (c) {</div>
|
|
<div class="line"><a id="l00053" name="l00053"></a><span class="lineno"> 53</span> <span class="keywordflow">case</span> 0:</div>
|
|
<div class="line"><a id="l00054" name="l00054"></a><span class="lineno"> 54</span> <span class="keywordflow">switch</span> (option_index) {</div>
|
|
<div class="line"><a id="l00055" name="l00055"></a><span class="lineno"> 55</span> <span class="keywordflow">case</span> 0: <span class="comment">// Not a mistake. This just goes to the default.</span></div>
|
|
<div class="line"><a id="l00056" name="l00056"></a><span class="lineno"> 56</span> <span class="keywordflow">default</span>:</div>
|
|
<div class="line"><a id="l00057" name="l00057"></a><span class="lineno"> 57</span> usage(argv[0]);</div>
|
|
<div class="line"><a id="l00058" name="l00058"></a><span class="lineno"> 58</span> }</div>
|
|
<div class="line"><a id="l00059" name="l00059"></a><span class="lineno"> 59</span> <span class="keywordflow">break</span>;</div>
|
|
<div class="line"><a id="l00060" name="l00060"></a><span class="lineno"> 60</span> <span class="keywordflow">case</span> <span class="stringliteral">'h'</span>:</div>
|
|
<div class="line"><a id="l00061" name="l00061"></a><span class="lineno"> 61</span> <span class="keywordflow">default</span>:</div>
|
|
<div class="line"><a id="l00062" name="l00062"></a><span class="lineno"> 62</span> usage(argv[0]);</div>
|
|
<div class="line"><a id="l00063" name="l00063"></a><span class="lineno"> 63</span> }</div>
|
|
<div class="line"><a id="l00064" name="l00064"></a><span class="lineno"> 64</span> }</div>
|
|
<div class="line"><a id="l00065" name="l00065"></a><span class="lineno"> 65</span> <span class="comment">// Check that the number of arguments is at least 8.</span></div>
|
|
<div class="line"><a id="l00066" name="l00066"></a><span class="lineno"> 66</span> <span class="keywordflow">if</span> (argc < 8) {</div>
|
|
<div class="line"><a id="l00067" name="l00067"></a><span class="lineno"> 67</span> usage(argv[0]);</div>
|
|
<div class="line"><a id="l00068" name="l00068"></a><span class="lineno"> 68</span> }</div>
|
|
<div class="line"><a id="l00069" name="l00069"></a><span class="lineno"> 69</span> </div>
|
|
<div class="line"><a id="l00070" name="l00070"></a><span class="lineno"> 70</span> <span class="comment">// Timing variables</span></div>
|
|
<div class="line"><a id="l00071" name="l00071"></a><span class="lineno"> 71</span> <span class="keywordtype">double</span> t0, t1;</div>
|
|
<div class="line"><a id="l00072" name="l00072"></a><span class="lineno"> 72</span> t0 = MPI_Wtime();</div>
|
|
<div class="line"><a id="l00073" name="l00073"></a><span class="lineno"> 73</span> </div>
|
|
<div class="line"><a id="l00074" name="l00074"></a><span class="lineno"> 74</span> <span class="comment">// Define/initialize variables</span></div>
|
|
<div class="line"><a id="l00075" name="l00075"></a><span class="lineno"> 75</span> <span class="keywordtype">double</span> start = atof(argv[1]), end = atof(argv[2]);</div>
|
|
<div class="line"><a id="l00076" name="l00076"></a><span class="lineno"> 76</span> <span class="keywordtype">int</span> points = atoi(argv[3]), cycles = atoi(argv[5]), L = atoi(argv[4]),</div>
|
|
<div class="line"><a id="l00077" name="l00077"></a><span class="lineno"> 77</span> burn_in_time = atoi(argv[6]), N = L * L;</div>
|
|
<div class="line"><a id="l00078" name="l00078"></a><span class="lineno"> 78</span> <span class="keywordtype">double</span> dt = (end - start) / points;</div>
|
|
<div class="line"><a id="l00079" name="l00079"></a><span class="lineno"> 79</span> std::ofstream ofile;</div>
|
|
<div class="line"><a id="l00080" name="l00080"></a><span class="lineno"> 80</span> std::string outfile = argv[7];</div>
|
|
<div class="line"><a id="l00081" name="l00081"></a><span class="lineno"> 81</span> <a class="code hl_class" href="classdata__t.html">data_t</a> data[points];</div>
|
|
<div class="line"><a id="l00082" name="l00082"></a><span class="lineno"> 82</span> </div>
|
|
<div class="line"><a id="l00083" name="l00083"></a><span class="lineno"> 83</span> <span class="comment">// MPI specific variables</span></div>
|
|
<div class="line"><a id="l00084" name="l00084"></a><span class="lineno"> 84</span> <span class="keywordtype">int</span> rank, cluster_size;</div>
|
|
<div class="line"><a id="l00085" name="l00085"></a><span class="lineno"> 85</span> </div>
|
|
<div class="line"><a id="l00086" name="l00086"></a><span class="lineno"> 86</span> <span class="comment">// Initialize MPI</span></div>
|
|
<div class="line"><a id="l00087" name="l00087"></a><span class="lineno"> 87</span> MPI_Init(&argc, &argv);</div>
|
|
<div class="line"><a id="l00088" name="l00088"></a><span class="lineno"> 88</span> </div>
|
|
<div class="line"><a id="l00089" name="l00089"></a><span class="lineno"> 89</span> <span class="comment">// Get the cluster size and rank</span></div>
|
|
<div class="line"><a id="l00090" name="l00090"></a><span class="lineno"> 90</span> MPI_Comm_size(MPI_COMM_WORLD, &cluster_size);</div>
|
|
<div class="line"><a id="l00091" name="l00091"></a><span class="lineno"> 91</span> MPI_Comm_rank(MPI_COMM_WORLD, &rank);</div>
|
|
<div class="line"><a id="l00092" name="l00092"></a><span class="lineno"> 92</span> </div>
|
|
<div class="line"><a id="l00093" name="l00093"></a><span class="lineno"> 93</span> <span class="keywordtype">int</span> remainder = points % cluster_size;</div>
|
|
<div class="line"><a id="l00094" name="l00094"></a><span class="lineno"> 94</span> <span class="keywordtype">double</span> i_start; <span class="comment">// What temperature to start from</span></div>
|
|
<div class="line"><a id="l00095" name="l00095"></a><span class="lineno"> 95</span> <span class="keywordtype">int</span> i_points; <span class="comment">// How many points to simulate</span></div>
|
|
<div class="line"><a id="l00096" name="l00096"></a><span class="lineno"> 96</span> </div>
|
|
<div class="line"><a id="l00097" name="l00097"></a><span class="lineno"> 97</span> <span class="comment">// Distribute temperature points</span></div>
|
|
<div class="line"><a id="l00098" name="l00098"></a><span class="lineno"> 98</span> <span class="keywordflow">if</span> (rank < remainder) {</div>
|
|
<div class="line"><a id="l00099" name="l00099"></a><span class="lineno"> 99</span> i_points = points / cluster_size + 1;</div>
|
|
<div class="line"><a id="l00100" name="l00100"></a><span class="lineno"> 100</span> i_start = start + dt * i_points * rank;</div>
|
|
<div class="line"><a id="l00101" name="l00101"></a><span class="lineno"> 101</span> }</div>
|
|
<div class="line"><a id="l00102" name="l00102"></a><span class="lineno"> 102</span> <span class="keywordflow">else</span> {</div>
|
|
<div class="line"><a id="l00103" name="l00103"></a><span class="lineno"> 103</span> i_points = points / cluster_size;</div>
|
|
<div class="line"><a id="l00104" name="l00104"></a><span class="lineno"> 104</span> i_start = start + dt * (i_points * rank + remainder);</div>
|
|
<div class="line"><a id="l00105" name="l00105"></a><span class="lineno"> 105</span> }</div>
|
|
<div class="line"><a id="l00106" name="l00106"></a><span class="lineno"> 106</span> </div>
|
|
<div class="line"><a id="l00107" name="l00107"></a><span class="lineno"> 107</span> <span class="comment">// Initialize array to contains data for each temperature point</span></div>
|
|
<div class="line"><a id="l00108" name="l00108"></a><span class="lineno"> 108</span> <a class="code hl_class" href="classdata__t.html">data_t</a> i_data[i_points];</div>
|
|
<div class="line"><a id="l00109" name="l00109"></a><span class="lineno"> 109</span> </div>
|
|
<div class="line"><a id="l00110" name="l00110"></a><span class="lineno"> 110</span> <span class="comment">// Simulate and save data to array</span></div>
|
|
<div class="line"><a id="l00111" name="l00111"></a><span class="lineno"> 111</span> <span class="keywordflow">for</span> (size_t i = 0; i < i_points; i++) {</div>
|
|
<div class="line"><a id="l00112" name="l00112"></a><span class="lineno"> 112</span> i_data[i] = montecarlo<a class="code hl_function" href="monte__carlo_8hpp.html#ae1e7f904ecfc3d8f3c4dd1ef155dd771">::</a><a class="code hl_function" href="monte__carlo_8hpp.html#ae1e7f904ecfc3d8f3c4dd1ef155dd771">mcmc_parallel</a><a class="code hl_function" href="monte__carlo_8hpp.html#ae1e7f904ecfc3d8f3c4dd1ef155dd771">(</a>L<a class="code hl_function" href="monte__carlo_8hpp.html#ae1e7f904ecfc3d8f3c4dd1ef155dd771">,</a> i_start + dt * i<a class="code hl_function" href="monte__carlo_8hpp.html#ae1e7f904ecfc3d8f3c4dd1ef155dd771">,</a> cycles<a class="code hl_function" href="monte__carlo_8hpp.html#ae1e7f904ecfc3d8f3c4dd1ef155dd771">,</a></div>
|
|
<div class="line"><a id="l00113" name="l00113"></a><span class="lineno"> 113</span> burn_in_time<a class="code hl_function" href="monte__carlo_8hpp.html#ae1e7f904ecfc3d8f3c4dd1ef155dd771">)</a>;</div>
|
|
<div class="line"><a id="l00114" name="l00114"></a><span class="lineno"> 114</span> }</div>
|
|
<div class="line"><a id="l00115" name="l00115"></a><span class="lineno"> 115</span> </div>
|
|
<div class="line"><a id="l00116" name="l00116"></a><span class="lineno"> 116</span> <span class="comment">// Rank 0 collects all the data and copies it to the "master"</span></div>
|
|
<div class="line"><a id="l00117" name="l00117"></a><span class="lineno"> 117</span> <span class="comment">// data array.</span></div>
|
|
<div class="line"><a id="l00118" name="l00118"></a><span class="lineno"> 118</span> <span class="keywordflow">if</span> (rank == 0) {</div>
|
|
<div class="line"><a id="l00119" name="l00119"></a><span class="lineno"> 119</span> <span class="comment">// Copy its own i_data to the data array</span></div>
|
|
<div class="line"><a id="l00120" name="l00120"></a><span class="lineno"> 120</span> std::copy_n(i_data, i_points, data);</div>
|
|
<div class="line"><a id="l00121" name="l00121"></a><span class="lineno"> 121</span> </div>
|
|
<div class="line"><a id="l00122" name="l00122"></a><span class="lineno"> 122</span> <span class="comment">// Collect i_data from other ranks in order and copy to data.</span></div>
|
|
<div class="line"><a id="l00123" name="l00123"></a><span class="lineno"> 123</span> <span class="keywordflow">for</span> (size_t i = 1; i < cluster_size; i++) {</div>
|
|
<div class="line"><a id="l00124" name="l00124"></a><span class="lineno"> 124</span> <span class="keywordflow">if</span> (rank < remainder) {</div>
|
|
<div class="line"><a id="l00125" name="l00125"></a><span class="lineno"> 125</span> MPI_Recv((<span class="keywordtype">void</span> *)i_data,</div>
|
|
<div class="line"><a id="l00126" name="l00126"></a><span class="lineno"> 126</span> <span class="keyword">sizeof</span>(data_t) * (points / cluster_size + 1), MPI_CHAR,</div>
|
|
<div class="line"><a id="l00127" name="l00127"></a><span class="lineno"> 127</span> i, MPI_ANY_TAG, MPI_COMM_WORLD, MPI_STATUS_IGNORE);</div>
|
|
<div class="line"><a id="l00128" name="l00128"></a><span class="lineno"> 128</span> std::copy_n(i_data, points / cluster_size + 1,</div>
|
|
<div class="line"><a id="l00129" name="l00129"></a><span class="lineno"> 129</span> data + (points / cluster_size) * i);</div>
|
|
<div class="line"><a id="l00130" name="l00130"></a><span class="lineno"> 130</span> }</div>
|
|
<div class="line"><a id="l00131" name="l00131"></a><span class="lineno"> 131</span> <span class="keywordflow">else</span> {</div>
|
|
<div class="line"><a id="l00132" name="l00132"></a><span class="lineno"> 132</span> MPI_Recv((<span class="keywordtype">void</span> *)i_data,</div>
|
|
<div class="line"><a id="l00133" name="l00133"></a><span class="lineno"> 133</span> <span class="keyword">sizeof</span>(data_t) * (points / cluster_size), MPI_CHAR, i,</div>
|
|
<div class="line"><a id="l00134" name="l00134"></a><span class="lineno"> 134</span> MPI_ANY_TAG, MPI_COMM_WORLD, MPI_STATUS_IGNORE);</div>
|
|
<div class="line"><a id="l00135" name="l00135"></a><span class="lineno"> 135</span> std::copy_n(i_data, points / cluster_size,</div>
|
|
<div class="line"><a id="l00136" name="l00136"></a><span class="lineno"> 136</span> data + (points / cluster_size) * i + remainder);</div>
|
|
<div class="line"><a id="l00137" name="l00137"></a><span class="lineno"> 137</span> }</div>
|
|
<div class="line"><a id="l00138" name="l00138"></a><span class="lineno"> 138</span> }</div>
|
|
<div class="line"><a id="l00139" name="l00139"></a><span class="lineno"> 139</span> </div>
|
|
<div class="line"><a id="l00140" name="l00140"></a><span class="lineno"> 140</span> <span class="comment">// Write everything from data to file</span></div>
|
|
<div class="line"><a id="l00141" name="l00141"></a><span class="lineno"> 141</span> utils<a class="code hl_function" href="utils_8hpp.html#a2b45adc86b70f42021582994e83fa00d">::</a><a class="code hl_function" href="utils_8hpp.html#a2b45adc86b70f42021582994e83fa00d">mkpath</a><a class="code hl_function" href="utils_8hpp.html#a2b45adc86b70f42021582994e83fa00d">(</a>utils<a class="code hl_function" href="utils_8hpp.html#aed026119193a9bbe076671809ff0f430">::</a><a class="code hl_function" href="utils_8hpp.html#aed026119193a9bbe076671809ff0f430">dirname</a><a class="code hl_function" href="utils_8hpp.html#aed026119193a9bbe076671809ff0f430">(</a>outfile<a class="code hl_function" href="utils_8hpp.html#aed026119193a9bbe076671809ff0f430">)</a><a class="code hl_function" href="utils_8hpp.html#a2b45adc86b70f42021582994e83fa00d">)</a>;</div>
|
|
<div class="line"><a id="l00142" name="l00142"></a><span class="lineno"> 142</span> ofile.open(outfile);</div>
|
|
<div class="line"><a id="l00143" name="l00143"></a><span class="lineno"> 143</span> </div>
|
|
<div class="line"><a id="l00144" name="l00144"></a><span class="lineno"> 144</span> <span class="keywordtype">double</span> temp, CV, X;</div>
|
|
<div class="line"><a id="l00145" name="l00145"></a><span class="lineno"> 145</span> </div>
|
|
<div class="line"><a id="l00146" name="l00146"></a><span class="lineno"> 146</span> <span class="keyword">using</span> utils::scientific_format;</div>
|
|
<div class="line"><a id="l00147" name="l00147"></a><span class="lineno"> 147</span> <span class="keywordflow">for</span> (size_t i = 0; i < points; i++) {</div>
|
|
<div class="line"><a id="l00148" name="l00148"></a><span class="lineno"> 148</span> temp = start + dt * i;</div>
|
|
<div class="line"><a id="l00149" name="l00149"></a><span class="lineno"> 149</span> CV = (data[i].E2 - data[i].E * data[i].E)</div>
|
|
<div class="line"><a id="l00150" name="l00150"></a><span class="lineno"> 150</span> / ((<span class="keywordtype">double</span>)N * temp * temp);</div>
|
|
<div class="line"><a id="l00151" name="l00151"></a><span class="lineno"> 151</span> X = (data[i].M2 - data[i].M_abs * data[i].M_abs)</div>
|
|
<div class="line"><a id="l00152" name="l00152"></a><span class="lineno"> 152</span> / ((<span class="keywordtype">double</span>)N * temp);</div>
|
|
<div class="line"><a id="l00153" name="l00153"></a><span class="lineno"> 153</span> </div>
|
|
<div class="line"><a id="l00154" name="l00154"></a><span class="lineno"> 154</span> ofile << <a class="code hl_function" href="utils_8hpp.html#a3529a74fd2a25d24de73d9d4e1c90835">scientific_format</a><a class="code hl_function" href="utils_8hpp.html#a3529a74fd2a25d24de73d9d4e1c90835">(</a>temp<a class="code hl_function" href="utils_8hpp.html#a3529a74fd2a25d24de73d9d4e1c90835">)</a> << <span class="stringliteral">','</span></div>
|
|
<div class="line"><a id="l00155" name="l00155"></a><span class="lineno"> 155</span> << scientific_format(data[i].E / N) << <span class="stringliteral">','</span></div>
|
|
<div class="line"><a id="l00156" name="l00156"></a><span class="lineno"> 156</span> << scientific_format(data[i].M_abs / N) << <span class="stringliteral">','</span></div>
|
|
<div class="line"><a id="l00157" name="l00157"></a><span class="lineno"> 157</span> << <a class="code hl_function" href="utils_8hpp.html#a3529a74fd2a25d24de73d9d4e1c90835">scientific_format</a><a class="code hl_function" href="utils_8hpp.html#a3529a74fd2a25d24de73d9d4e1c90835">(</a>CV<a class="code hl_function" href="utils_8hpp.html#a3529a74fd2a25d24de73d9d4e1c90835">)</a> << <span class="stringliteral">','</span> << <a class="code hl_function" href="utils_8hpp.html#a3529a74fd2a25d24de73d9d4e1c90835">scientific_format</a><a class="code hl_function" href="utils_8hpp.html#a3529a74fd2a25d24de73d9d4e1c90835">(</a>X<a class="code hl_function" href="utils_8hpp.html#a3529a74fd2a25d24de73d9d4e1c90835">)</a></div>
|
|
<div class="line"><a id="l00158" name="l00158"></a><span class="lineno"> 158</span> << <span class="stringliteral">'\n'</span>;</div>
|
|
<div class="line"><a id="l00159" name="l00159"></a><span class="lineno"> 159</span> }</div>
|
|
<div class="line"><a id="l00160" name="l00160"></a><span class="lineno"> 160</span> ofile.close();</div>
|
|
<div class="line"><a id="l00161" name="l00161"></a><span class="lineno"> 161</span> }</div>
|
|
<div class="line"><a id="l00162" name="l00162"></a><span class="lineno"> 162</span> <span class="comment">// For all other ranks, send the data to rank 0</span></div>
|
|
<div class="line"><a id="l00163" name="l00163"></a><span class="lineno"> 163</span> <span class="keywordflow">else</span> {</div>
|
|
<div class="line"><a id="l00164" name="l00164"></a><span class="lineno"> 164</span> MPI_Send(i_data, i_points * <span class="keyword">sizeof</span>(data_t), MPI_CHAR, 0, rank,</div>
|
|
<div class="line"><a id="l00165" name="l00165"></a><span class="lineno"> 165</span> MPI_COMM_WORLD);</div>
|
|
<div class="line"><a id="l00166" name="l00166"></a><span class="lineno"> 166</span> }</div>
|
|
<div class="line"><a id="l00167" name="l00167"></a><span class="lineno"> 167</span> </div>
|
|
<div class="line"><a id="l00168" name="l00168"></a><span class="lineno"> 168</span> t1 = MPI_Wtime();</div>
|
|
<div class="line"><a id="l00169" name="l00169"></a><span class="lineno"> 169</span> </div>
|
|
<div class="line"><a id="l00170" name="l00170"></a><span class="lineno"> 170</span> <span class="keywordflow">if</span> (rank == 0) {</div>
|
|
<div class="line"><a id="l00171" name="l00171"></a><span class="lineno"> 171</span> std::cout << <span class="stringliteral">"Time: "</span> << t1 - t0 << <span class="stringliteral">" seconds\n"</span>;</div>
|
|
<div class="line"><a id="l00172" name="l00172"></a><span class="lineno"> 172</span> }</div>
|
|
<div class="line"><a id="l00173" name="l00173"></a><span class="lineno"> 173</span> </div>
|
|
<div class="line"><a id="l00174" name="l00174"></a><span class="lineno"> 174</span> MPI_Finalize();</div>
|
|
<div class="line"><a id="l00175" name="l00175"></a><span class="lineno"> 175</span>}</div>
|
|
<div class="ttc" id="aclassdata__t_html"><div class="ttname"><a href="classdata__t.html">data_t</a></div><div class="ttdoc">Type to use with the IsingModel class and montecarlo module.</div><div class="ttdef"><b>Definition:</b> <a href="data__type_8hpp_source.html#l00019">data_type.hpp:19</a></div></div>
|
|
<div class="ttc" id="amain_8cpp_html_a3c04138a5bfe5d72780bb7e82a18e627"><div class="ttname"><a href="main_8cpp.html#a3c04138a5bfe5d72780bb7e82a18e627">main</a></div><div class="ttdeci">int main(int argc, char **argv)</div><div class="ttdoc">The main function.</div><div class="ttdef"><b>Definition:</b> <a href="main_8cpp_source.html#l00125">main.cpp:125</a></div></div>
|
|
<div class="ttc" id="amain_8cpp_html_ac907e18135856c90366aaa599a9e10b1"><div class="ttname"><a href="main_8cpp.html#ac907e18135856c90366aaa599a9e10b1">usage</a></div><div class="ttdeci">void usage(std::string filename)</div><div class="ttdoc">A function that displays how to use the program and quits.</div><div class="ttdef"><b>Definition:</b> <a href="main_8cpp_source.html#l00110">main.cpp:110</a></div></div>
|
|
<div class="ttc" id="amonte__carlo_8hpp_html_ae1e7f904ecfc3d8f3c4dd1ef155dd771"><div class="ttname"><a href="monte__carlo_8hpp.html#ae1e7f904ecfc3d8f3c4dd1ef155dd771">montecarlo::mcmc_parallel</a></div><div class="ttdeci">data_t mcmc_parallel(int L, double T, int cycles, int burn_in_time=BURN_IN_TIME)</div><div class="ttdoc">Execute the Metropolis algorithm for a certain amount of Monte Carlo cycles in parallel.</div><div class="ttdef"><b>Definition:</b> <a href="monte__carlo_8cpp_source.html#l00127">monte_carlo.cpp:127</a></div></div>
|
|
<div class="ttc" id="autils_8hpp_html_a2b45adc86b70f42021582994e83fa00d"><div class="ttname"><a href="utils_8hpp.html#a2b45adc86b70f42021582994e83fa00d">utils::mkpath</a></div><div class="ttdeci">bool mkpath(std::string path, int mode=0777)</div><div class="ttdoc">Make path given.</div><div class="ttdef"><b>Definition:</b> <a href="utils_8cpp_source.html#l00032">utils.cpp:32</a></div></div>
|
|
<div class="ttc" id="autils_8hpp_html_a3529a74fd2a25d24de73d9d4e1c90835"><div class="ttname"><a href="utils_8hpp.html#a3529a74fd2a25d24de73d9d4e1c90835">utils::scientific_format</a></div><div class="ttdeci">std::string scientific_format(double d, int width=20, int prec=10)</div><div class="ttdoc">Turns a double into a string written in scientific format.</div><div class="ttdef"><b>Definition:</b> <a href="utils_8cpp_source.html#l00016">utils.cpp:16</a></div></div>
|
|
<div class="ttc" id="autils_8hpp_html_aed026119193a9bbe076671809ff0f430"><div class="ttname"><a href="utils_8hpp.html#aed026119193a9bbe076671809ff0f430">utils::dirname</a></div><div class="ttdeci">std::string dirname(const std::string &path)</div><div class="ttdoc">Get the directory name of the path.</div><div class="ttdef"><b>Definition:</b> <a href="utils_8cpp_source.html#l00058">utils.cpp:58</a></div></div>
|
|
</div><!-- fragment --></div><!-- contents -->
|
|
</div><!-- doc-content -->
|
|
<!-- start footer part -->
|
|
<div id="nav-path" class="navpath"><!-- id is needed for treeview function! -->
|
|
<ul>
|
|
<li class="navelem"><a class="el" href="dir_68267d1309a1af8e8297ef4c3efbcdba.html">src</a></li><li class="navelem"><a class="el" href="phase__transition__mpi_8cpp.html">phase_transition_mpi.cpp</a></li>
|
|
<li class="footer">Generated by <a href="https://www.doxygen.org/index.html"><img class="footer" src="doxygen.svg" width="104" height="31" alt="doxygen"/></a> 1.9.6 </li>
|
|
</ul>
|
|
</div>
|
|
</body>
|
|
</html>
|