source: CIVL/examples/mpi/mpi_pi.c

main
Last change on this file was ea777aa, checked in by Alex Wilton <awilton@…>, 3 years ago

Moved examples, include, build_default.properties, common.xml, and README out from dev.civl.com into the root of the repo.

git-svn-id: svn://vsl.cis.udel.edu/civl/trunk@5704 fb995dde-84ed-4084-dfe6-e5aef3e2452c

  • Property mode set to 100644
File size: 6.1 KB
RevLine 
[3ff27cf]1#ifdef _CIVL
2#include <civlc.cvh>
3#endif
[a96092d]4/* FILENAME: mpi_pi.c */
5/* Calculating value of pi using a "dartboard" algorithm. In CIVL
6 * mode, the program will first let process of rank 0 do a sequential
7 * run and the results of it will be used to compare with the result
8 * of each round during parallel run.
9 * To execute: mpicc mpi_pi.c && mpiexec -n 4 ./a.out
10 * Here '4' can be replaced by any number.
11 * To verify: civl verify mpi_pi.c
12 *
13 * Modified from the original program: mpi_pi_send.c
14 * Source: https://computing.llnl.gov/tutorials/mpi/samples/C/mpi_pi_send.c */
15#include <mpi.h>
16#include <stdio.h>
17#include <stdlib.h>
18#include <assert.h>
19
20#define MASTER 0
21#define SQR(X) ((X)*(X))
22
23#ifdef _CIVL
[58207ec]24$input int DARTSB = 2; // upper bound of DARTS
[a96092d]25$input int DARTS; // number of darts will be throwed
[3ff27cf]26$assume(0 < DARTS && DARTS <= DARTSB);
[58207ec]27$input int ROUNDSB = 2; // upper bound of ROUNDS
[a96092d]28$input int ROUNDS; // number of rounds of throwing darts
[3ff27cf]29$assume(0 < ROUNDS && ROUNDS <= ROUNDSB);
[db66e7f]30$input int _mpi_nprocs = 2;
[a96092d]31$input int N; // length of random data array
[58207ec]32$assume(N > 0);
[a96092d]33$input double RANDOM[N]; // random data array
34double oracle[ROUNDS]; // array of results of sequential run
35
36/* CIVL mode rand() function: returns unconstrained value as a random
37 number for a specific process throwing a specific dart right at the
38 given round. So that the sequential version and parallel version
39 are guaranteed to compute with a same set of ordered random
40 data. */
41double civlrand(int dart, int taskid, int curr_round) {
42 return RANDOM[taskid * (ROUNDS * DARTS) + dart * ROUNDS + curr_round];
43}
44#else
45int DARTS, ROUNDS;
46#endif
47int curr_round; // current round
48int taskid; // the identity of a task which is
49 // a process in parallel run.
50int tasks; // the number of tasks whichis the number
51 // of processes in parallel run
52
53/***********************************************************************
54 * Throw darts at board. Done by generating random numbers
55 * between 0 and 1 and converting them to values for x and y
56 * coordinates and then testing to see if they "land" in
57 * the circle." If so, score is incremented. After throwing the
58 * specified number of darts, pi is calculated. The computed value
59 * of pi is returned as the value of this function, dboard.
60 ************************************************************************/
61double dboard(int darts, int taskid) {
62 double x_coord, // x coordinates, between -1 and 1
63 y_coord, // y coordinates, between -1 and 1
64 pi, // pi
65 rx, // random number for x coordinate
66 ry; // random number for y coordinate
67 int score, // number of darts hit the circle
68 n;
69
70 score = 0;
71 for (n = 0; n < darts; n++) {
72 /* generate random numbers for x and y coordinates */
73#ifdef _CIVL
74 rx = civlrand(n, taskid, curr_round) / (double)RAND_MAX;
[db66e7f]75 ry = civlrand(n, taskid, curr_round+DARTSB * _mpi_nprocs * ROUNDSB) \
[a96092d]76 / (double)RAND_MAX;
77#else
78 rx = (double)rand()/(double)RAND_MAX;
79 ry = (double)rand()/(double)RAND_MAX;
80#endif
81 x_coord = (2.0 * rx) - 1.0;
82 y_coord = (2.0 * ry) - 1.0;
83
84 /* if dart lands in circle, increment score */
85 if ((SQR(x_coord) + SQR(y_coord)) <= 1.0)
86 score++;
87 }
88 /* calculate pi */
89 pi = 4.0 * (double)score/(double)darts;
90 return(pi);
91}
92
93/* In CIVL mode, process of rank 0 will do a sequential and store the
94 results in "oracle". MPI program will initialize global
95 variables */
96void initialization() {
97#ifdef _CIVL
98 double pi = 0.0;
99 double avepi = 0.0;
100
[db66e7f]101 $assume(N == DARTSB * _mpi_nprocs * ROUNDSB * 2);
[a96092d]102 if(taskid == 0) {
103 for(curr_round=0; curr_round < ROUNDS; curr_round++) {
104 for(int j=0; j < tasks; j++)
105 pi += dboard(DARTS, j);
106 pi = pi / tasks;
107 avepi = ((avepi * curr_round) + pi)/(curr_round + 1);
108 oracle[curr_round] = avepi;
109 pi = 0.0;
110 }
111 }
112#else
113 DARTS = 100;
114 ROUNDS = 500;
115#endif
116}
117
118int main(int argc, char * argv[]) {
119 double homepi; // local calculated pi
120 double avepi = 0; // average value of pi of all tasks
121 //and executed rounds
122 /* Each value of a round will be used as a message tag to
123 communicate among processes at that round */
124 int mtype;
125 int rc; // returned signal of MPI routines
126
127 MPI_Init(&argc, &argv);
128 MPI_Comm_size(MPI_COMM_WORLD, &tasks);
129 MPI_Comm_rank(MPI_COMM_WORLD, &taskid);
130 initialization();
131 for(curr_round=0; curr_round<ROUNDS;curr_round++) {
132 homepi = dboard(DARTS, taskid);
133 if(taskid != MASTER){
134 mtype = curr_round;
135 rc = MPI_Send(&homepi, 1, MPI_DOUBLE,
136 MASTER, mtype, MPI_COMM_WORLD);
137 if (rc != MPI_SUCCESS)
138 printf("%d: Send failure on round %d\n", taskid, mtype);
139 } else {
140 double pi; // value of pi in a single round
141 double pisum; // sum of the values of pi calculated by all tasks
142 double pirecv; // MPI receive buffer
143
144 mtype = curr_round;
145 pisum = 0;
146 for (int n = 1; n < tasks; n++) {
147 rc = MPI_Recv(&pirecv, 1, MPI_DOUBLE, MPI_ANY_SOURCE,
148 mtype, MPI_COMM_WORLD, MPI_STATUS_IGNORE);
149 if (rc != MPI_SUCCESS)
150 printf("%d: Receive failure on round %d\n", taskid, mtype);
151 /* keep running total of pi */
152 pisum = pisum + pirecv;
153 }
154 /* Master calculates the average value of pi for this iteration */
155 pi = (pisum + homepi)/tasks;
156 /* Master calculates the average value of pi over all iterations */
157 avepi = ((avepi * curr_round) + pi)/(curr_round + 1);
158 printf(" After %8d throws, average value of pi = %10.8f\n",
159 (DARTS * (curr_round + 1)),avepi);
160#ifdef _CIVL
[3ff27cf]161 $assert((avepi == oracle[curr_round]), "avepi is %f but oracle[%d]"
162 " is %f\n", avepi, curr_round, oracle[curr_round]);
[a96092d]163#endif
164 }
165 }
166 if (taskid == MASTER)
167 printf ("\nReal value of PI: 3.1415926535897 \n");
168 MPI_Finalize();
169 return 0;
170}
Note: See TracBrowser for help on using the repository browser.