2 Dimensional Ising Model
Simulate the change in energy and magnetization in a ferro magnet
Loading...
Searching...
No Matches
phase_transition_mpi.cpp
Go to the documentation of this file.
1
16#include "data_type.hpp"
17#include "monte_carlo.hpp"
18#include "utils.hpp"
19
20#include <getopt.h>
21#include <mpi.h>
22#include <string>
23
25void usage(std::string filename)
26{
27 std::cout
28 << "Usage: " << filename
29 << " <start temperature> <end temperature> <lattice size> "
30 "<points> <cycles> <burn-in-time> <output file>\n"
31 << "This should be used with mpiexec or mpirun for maximum "
32 "performance\n\n"
33 << "\t[ -h | --help ]\n";
34 exit(-1);
35}
36
38int main(int argc, char **argv)
39{
40 // Command options
41 struct option long_options[] = {{"help", 0, 0, 0}, {NULL, 0, NULL, 0}};
42
43 int option_index = -1;
44 int c;
45
46 while (true) {
47 c = getopt_long(argc, argv, "h", long_options, &option_index);
48
49 if (c == -1)
50 break;
51
52 switch (c) {
53 case 0:
54 switch (option_index) {
55 case 0: // Not a mistake. This just goes to the default.
56 default:
57 usage(argv[0]);
58 }
59 break;
60 case 'h':
61 default:
62 usage(argv[0]);
63 }
64 }
65 // Check that the number of arguments is at least 8.
66 if (argc < 8) {
67 usage(argv[0]);
68 }
69
70 // Timing variables
71 double t0, t1;
72 t0 = MPI_Wtime();
73
74 // Define/initialize variables
75 double start = atof(argv[1]), end = atof(argv[2]);
76 int points = atoi(argv[3]), cycles = atoi(argv[5]), L = atoi(argv[4]),
77 burn_in_time = atoi(argv[6]), N = L * L;
78 double dt = (end - start) / points;
79 std::ofstream ofile;
80 std::string outfile = argv[7];
81 data_t data[points];
82
83 // MPI specific variables
84 int rank, cluster_size;
85
86 // Initialize MPI
87 MPI_Init(&argc, &argv);
88
89 // Get the cluster size and rank
90 MPI_Comm_size(MPI_COMM_WORLD, &cluster_size);
91 MPI_Comm_rank(MPI_COMM_WORLD, &rank);
92
93 int remainder = points % cluster_size;
94 double i_start; // What temperature to start from
95 int i_points; // How many points to simulate
96
97 // Distribute temperature points
98 if (rank < remainder) {
99 i_points = points / cluster_size + 1;
100 i_start = start + dt * i_points * rank;
101 }
102 else {
103 i_points = points / cluster_size;
104 i_start = start + dt * (i_points * rank + remainder);
105 }
106
107 // Initialize array to contains data for each temperature point
108 data_t i_data[i_points];
109
110 // Simulate and save data to array
111 for (size_t i = 0; i < i_points; i++) {
112 i_data[i] = montecarlo::mcmc_parallel(L, i_start + dt * i, cycles,
113 burn_in_time);
114 }
115
116 // Rank 0 collects all the data and copies it to the "master"
117 // data array.
118 if (rank == 0) {
119 // Copy its own i_data to the data array
120 std::copy_n(i_data, i_points, data);
121
122 // Collect i_data from other ranks in order and copy to data.
123 for (size_t i = 1; i < cluster_size; i++) {
124 if (rank < remainder) {
125 MPI_Recv((void *)i_data,
126 sizeof(data_t) * (points / cluster_size + 1), MPI_CHAR,
127 i, MPI_ANY_TAG, MPI_COMM_WORLD, MPI_STATUS_IGNORE);
128 std::copy_n(i_data, points / cluster_size + 1,
129 data + (points / cluster_size) * i);
130 }
131 else {
132 MPI_Recv((void *)i_data,
133 sizeof(data_t) * (points / cluster_size), MPI_CHAR, i,
134 MPI_ANY_TAG, MPI_COMM_WORLD, MPI_STATUS_IGNORE);
135 std::copy_n(i_data, points / cluster_size,
136 data + (points / cluster_size) * i + remainder);
137 }
138 }
139
140 // Write everything from data to file
141 utils::mkpath(utils::dirname(outfile));
142 ofile.open(outfile);
143
144 double temp, CV, X;
145
146 using utils::scientific_format;
147 for (size_t i = 0; i < points; i++) {
148 temp = start + dt * i;
149 CV = (data[i].E2 - data[i].E * data[i].E)
150 / ((double)N * temp * temp);
151 X = (data[i].M2 - data[i].M_abs * data[i].M_abs)
152 / ((double)N * temp);
153
154 ofile << scientific_format(temp) << ','
155 << scientific_format(data[i].E / N) << ','
156 << scientific_format(data[i].M_abs / N) << ','
157 << scientific_format(CV) << ',' << scientific_format(X)
158 << '\n';
159 }
160 ofile.close();
161 }
162 // For all other ranks, send the data to rank 0
163 else {
164 MPI_Send(i_data, i_points * sizeof(data_t), MPI_CHAR, 0, rank,
165 MPI_COMM_WORLD);
166 }
167
168 t1 = MPI_Wtime();
169
170 if (rank == 0) {
171 std::cout << "Time: " << t1 - t0 << " seconds\n";
172 }
173
174 MPI_Finalize();
175}
Type to use with the IsingModel class and montecarlo module.
Definition: data_type.hpp:19
double M_abs
Absolute Magnetization.
Definition: data_type.hpp:25
double E
Energy.
Definition: data_type.hpp:21
double M2
Magnetization squared.
Definition: data_type.hpp:24
double E2
Energy squared.
Definition: data_type.hpp:23
Header for the data_t type.
Functions for Monte Carlo simulations.
void usage(std::string filename)
A function that displays how to use the program and quits.
int main()
The main function.
Definition: test_suite.cpp:110
Function prototypes and macros that are useful.