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
12#include "data_type.hpp"
13#include "monte_carlo.hpp"
14#include "utils.hpp"
15
16#include <algorithm>
17#include <fstream>
18#include <iostream>
19#include <iterator>
20#include <mpi.h>
21#include <sstream>
22
24int main(int argc, char **argv)
25{
26 if (argc < 5) {
27 std::cout << "You need at least 4 arguments" << std::endl;
28 abort();
29 }
30 double t0, t1;
31 t0 = MPI_Wtime();
32 double start = atof(argv[1]), end = atof(argv[2]);
33 int points = atoi(argv[3]), N;
34 int lattice_sizes[] = {20, 40, 60, 80, 100};
35 double dt = (end - start) / points;
36 int cycles = atoi(argv[4]);
37 std::ofstream ofile;
38
39 data_t data[points];
40
41 // MPI stuff
42 int rank, cluster_size;
43
44 // Initialize MPI
45 MPI_Init(&argc, &argv);
46
47 // Get the cluster size and rank
48 MPI_Comm_size(MPI_COMM_WORLD, &cluster_size);
49 MPI_Comm_rank(MPI_COMM_WORLD, &rank);
50
51 int remainder = points % cluster_size;
52 double i_start;
53 int i_points;
54 // Distribute temperature points
55 if (rank < remainder) {
56 i_points = points / cluster_size + 1;
57 i_start = start + dt * i_points * rank;
58 }
59 else {
60 i_points = points / cluster_size;
61 i_start = start + dt * (i_points * rank + remainder);
62 }
63
64 data_t i_data[i_points];
65 std::cout << "Rank " << rank << ": " << i_points << ',' << i_start << '\n';
66
67 for (int L : lattice_sizes) {
68 N = L * L;
69 for (size_t i = 0; i < i_points; i++) {
70 i_data[i] = monte_carlo_parallel(L, i_start + dt * i, cycles);
71 }
72
73 if (rank == 0) {
74 std::copy_n(i_data, i_points, data);
75 for (size_t i = 1; i < cluster_size; i++) {
76 if (rank < remainder) {
77 MPI_Recv((void *)i_data,
78 sizeof(data_t) * (points / cluster_size + 1),
79 MPI_CHAR, i, MPI_ANY_TAG, MPI_COMM_WORLD,
80 MPI_STATUS_IGNORE);
81 std::copy_n(i_data, points / cluster_size + 1,
82 data + (points / cluster_size) * i);
83 }
84 else {
85 MPI_Recv((void *)i_data,
86 sizeof(data_t) * (points / cluster_size), MPI_CHAR,
87 i, MPI_ANY_TAG, MPI_COMM_WORLD, MPI_STATUS_IGNORE);
88 std::copy_n(i_data, points / cluster_size,
89 data + (points / cluster_size) * i + remainder);
90 }
91 }
92 std::stringstream outfile;
93 outfile << "output/phase_transition/size_" << L << ".txt";
94 utils::mkpath(utils::dirname(outfile.str()));
95 ofile.open(outfile.str());
96
97 double temp, CV, X;
98
99 using utils::scientific_format;
100 for (size_t i = 0; i < points; i++) {
101 temp = start + dt * i;
102 CV = (data[i].E2 - data[i].E * data[i].E)
103 / ((double)N * temp * temp);
104 X = (data[i].M2 - data[i].M_abs * data[i].M_abs)
105 / ((double)N * temp);
106
107 ofile << scientific_format(temp) << ','
108 << scientific_format(data[i].E / N) << ','
109 << scientific_format(data[i].M_abs / N) << ','
110 << scientific_format(CV) << ',' << scientific_format(X)
111 << '\n';
112 }
113 ofile.close();
114 }
115 else {
116 MPI_Send(i_data, i_points * sizeof(data_t), MPI_CHAR, 0, rank,
117 MPI_COMM_WORLD);
118 }
119 }
120
121 t1 = MPI_Wtime();
122
123 if (rank == 0) {
124 std::cout << "Time: " << t1 - t0 << " seconds\n";
125 }
126
127 MPI_Finalize();
128}
Header for the data_t type.
Functions for monte carlo simulations.
data_t monte_carlo_parallel(int L, double T, int cycles)
Execute the Metropolis algorithm for a certain amount of Monte Carlo cycles in parallel.
int main()
The main function.
Definition: test_suite.cpp:148
Function prototypes and macros that are useful.