source: CIVL/examples/contracts/contractsMPI/diffusion2d.c@ bb03188

main test-branch
Last change on this file since bb03188 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: 5.7 KB
Line 
1#include <mpi.h>
2
3#define FROMLEFT 0
4#define FROMRIGHT 1
5#define FROMBOTTOM 2
6#define FROMTOP 3
7#define comm MPI_COMM_WORLD
8
9#pragma CIVL ACSL
10
11double (*u)[];
12double (*u_new)[];
13double k;
14int nxl, nyl;
15int rank, nprocsx, nprocsy, left, right, top, bottom;
16
17/*@ requires \valid((double (*)[nxl+2])u + (0 .. (nyl + 1)));
18 @ requires \valid(((double (*)[nxl+2])u_new) + (0 .. (nyl + 1)));
19 @ requires k > 0;
20 @ requires 0 <= nxl && 0 <= nyl;
21 @ assigns u[1 .. nyl][1 .. nxl];
22 @ ensures \forall int i; 1 <= i && i <= nyl ==>
23 @ ( \forall int j; 1 <= j && j <= nxl ==>
24 @ u[i][j] == \old(u[i][j] +
25 @ k*(u[i+1][j] + u[i-1][j] +
26 @ u[i][j+1] + u[i][j-1] - 4*u[i][j]))
27 @ );
28 @*/
29void update() {
30 double (* tmp)[];
31
32 /*@ loop invariant 1 <= i && i <= nyl+1;
33 @ loop invariant \forall int t; 1 <= t && t < i ==>
34 @ ( \forall int j0; 1 <= j0 && j0 <= nxl ==>
35 @ u_new[t][j0] == u[t][j0] + k*(u[t-1][j0] + u[t+1][j0] + u[t][j0+1] + u[t][j0-1] - 4 * u[t][j0])
36 @ );
37 @ loop assigns u_new[1 .. nyl][1 .. nxl], i;
38 @*/
39 for (int i = 1; i < nyl + 1; i++) {
40 /*@ loop invariant 1 <= j && j <= nxl+1;
41 @ loop invariant \forall int t; 1 <= t && t < j
42 @ ==> u_new[i][t] == u[i][t] + k*(u[i-1][t] + u[i+1][t] + u[i][t+1] + u[i][t-1] - 4 * u[i][t]);
43 @ loop assigns u_new[i][1 .. nxl], j;
44 @*/
45 for (int j = 1; j < nxl + 1; j++) {
46 u_new[i][j] = u[i][j] +
47 k*(u[i+1][j] + u[i-1][j] +
48 u[i][j+1] + u[i][j-1] - 4*u[i][j]);
49 }
50 }
51 tmp = u;
52 u = u_new;
53 u_new = tmp;
54}
55
56
57/*@ requires \valid((double (*)[nxl+2])u + (0 .. (nyl + 1)));
58 @ requires \valid((double (*)[nxl+2])u_new + (0 .. (nyl + 1)));
59 @ requires 0 < nxl && 0 < nyl;
60 @ \mpi_collective(MPI_COMM_WORLD, P2P):
61 @ requires rank == \mpi_comm_rank;
62 @ requires left == MPI_PROC_NULL ||
63 @ (0 <= left && left < \mpi_comm_size && left != rank);
64 @ requires right == MPI_PROC_NULL ||
65 @ (0 <= right && right < \mpi_comm_size && right != rank);
66 @ requires top == MPI_PROC_NULL ||
67 @ (0 <= top && top < \mpi_comm_size && top != rank);
68 @ requires bottom == MPI_PROC_NULL ||
69 @ (0 <= bottom && bottom < \mpi_comm_size && bottom != rank);
70 @ behavior hasTop:
71 @ assumes top != MPI_PROC_NULL;
72 @ requires rank == \on(top, bottom) && top != left && top != right && top != bottom;
73 @ requires nxl == \on(top, nxl);
74 @ assigns u[0][1 .. nxl];
75 @ ensures \forall int i; 1 <= i && i <= nxl ==> u[0][i] == \on(top, u[nyl][i]);
76 @ behavior hasBottom:
77 @ assumes bottom != MPI_PROC_NULL;
78 @ requires rank == \on(bottom, top) && bottom != left && bottom != right && bottom != top;
79 @ requires nxl == \on(bottom, nxl);
80 @ assigns u[nyl+1][1 .. nxl];
81 @ ensures \forall int i; 1 <= i && i <= nxl ==> u[nyl+1][i] == \on(bottom, u[1][i]);
82 @ behavior hasLeft:
83 @ assumes left != MPI_PROC_NULL;
84 @ requires rank == \on(left, right) && left != right && left != bottom && left != top;
85 @ requires nyl == \on(left, nyl);
86 @ assigns u[0][1 .. nyl];
87 @ ensures \forall int i; 1 <= i && i <= nyl ==> u[i][0] == \on(left, u[i][nxl]);
88 @ behavior hasRight:
89 @ assumes right != MPI_PROC_NULL;
90 @ requires rank == \on(right, left) && right != left && right != top && right != bottom;
91 @ requires nyl == \on(right, nyl);
92 @ assigns u[nxl+1][1 .. nyl];
93 @ ensures \forall int i; 1 <= i && i <= nyl ==> u[i][nxl+1] == \on(right, u[i][1]);
94 @*/
95void exchange() {
96 // sends top border row, receives into bottom ghost cell row
97 MPI_Sendrecv(&u[1][1], nxl, MPI_DOUBLE, top, FROMBOTTOM, &u[nyl+1][1], nxl,
98 MPI_DOUBLE, bottom, FROMBOTTOM, comm, MPI_STATUS_IGNORE);
99 // sends bottom border row, receives into top ghost cell row
100 MPI_Sendrecv(&u[nyl][1], nxl, MPI_DOUBLE, bottom, FROMTOP, &u[0][1], nxl,
101 MPI_DOUBLE, top, FROMTOP, comm, MPI_STATUS_IGNORE);
102
103 if (nyl > 0) {
104 double sendbuf[nyl];
105 double recvbuf[nyl];
106
107 // sends left border column, receives into temporary buffer
108 /*@ loop invariant 0 <= i && i <= nyl;
109 @ loop invariant \forall int t; 0 <= t && t < i
110 @ ==> sendbuf[t] == u[t+1][1];
111 @ loop assigns sendbuf[0 .. nyl-1], i;
112 @*/
113 for (int i = 0; i < nyl; i++) sendbuf[i] = u[i+1][1];
114 MPI_Sendrecv(sendbuf, nyl, MPI_DOUBLE, left, FROMRIGHT, recvbuf, nyl,
115 MPI_DOUBLE, right, FROMRIGHT, comm, MPI_STATUS_IGNORE);
116 // copies temporary buffer into right ghost cell column
117 if (right != MPI_PROC_NULL) {
118 /*@ loop invariant 0 <= i && i <= nyl;
119 @ loop invariant \forall int t; 0 <= t && t < i ==>
120 @ recvbuf[t] == u[t+1][nxl+1];
121 @ loop assigns u[1 .. nyl][nxl + 1], i;
122 @*/
123 for (int i = 0; i < nyl; i++) u[i+1][nxl+1] = recvbuf[i];
124 }
125 // sends right border column, receives into temporary buffer
126 /*@ loop invariant 0 <= i && i <= nyl;
127 @ loop invariant \forall int t; 0 <= t && t < i ==>
128 @ sendbuf[t] == u[t+1][nxl];
129 @ loop assigns sendbuf[0 .. nyl-1], i;
130 @*/
131 for (int i = 0; i < nyl; i++) sendbuf[i] = u[i+1][nxl];
132 MPI_Sendrecv(sendbuf, nyl, MPI_DOUBLE, right, FROMLEFT, recvbuf, nyl,
133 MPI_DOUBLE, left, FROMLEFT, comm, MPI_STATUS_IGNORE);
134 // copies temporary buffer into left ghost cell column
135 if (left != MPI_PROC_NULL) {
136 /*@ loop invariant 0 <= i && i <= nyl;
137 @ loop invariant \forall int t; 0 <= t && t < i ==>
138 @ recvbuf[t] == u[t+1][0];
139 @ loop assigns u[1 .. nyl][0], i;
140 @*/
141 for (int i = 0; i < nyl; i++) u[i+1][0] = recvbuf[i];
142 }
143 }
144}
145
146
Note: See TracBrowser for help on using the repository browser.