source: CIVL/examples/contracts/contractsMPI/allgather2bad.c@ 02876df

1.23 2.0 main test-branch
Last change on this file since 02876df was bd0cedf, checked in by Ziqing Luo <ziqing@…>, 10 years ago

add implementation of "assigns"

add transformation for "old" expression

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

  • Property mode set to 100644
File size: 4.2 KB
Line 
1// How to get rid of BOUNDS ?
2/************************** source code **************************/
3#include<mpi.h>
4#include<civl-mpi.cvh>
5#include<civlc.cvh>
6#include<string.h>
7
8#define BUFFER_BOUND 3
9
10/*@
11 @ \mpi_collective(comm, P2P):
12 @ requires 0 <= root && root < \mpi_comm_size;
13 @ requires \mpi_agree(root) && \mpi_agree(count * \mpi_extent(datatype));
14 @ requires 0 <= count && count * \mpi_extent(datatype) < 10;
15 @ requires \mpi_valid(buf, count, datatype);
16 @ behavior root:
17 @ assumes \mpi_comm_rank == root;
18 @ behavior others:
19 @ assumes \mpi_comm_rank != root;
20 @ assigns \mpi_region(buf, count, datatype);
21 @ ensures \mpi_equals(buf, count, datatype, \on(root, buf));
22 @ waitsfor (root .. root);
23 @*/
24int broadcast(void * buf, int count,
25 MPI_Datatype datatype, int root, MPI_Comm comm) {
26 int nprocs, rank;
27 int tag = 999;
28
29 MPI_Comm_size(comm, &nprocs);
30 MPI_Comm_rank(comm, &rank);
31 if (rank == root) {
32 for (int i = 0; i < nprocs; i++)
33 if (i != root)
34 MPI_Send(buf, count, datatype, i, tag, comm);
35 } else
36 MPI_Recv(buf, count, datatype, root, tag, comm,
37 MPI_STATUS_IGNORE);
38 return 0;
39}
40
41/*@
42 @ \mpi_collective(comm, P2P) :
43 @ requires \mpi_agree(root) && \mpi_agree(sendcount * \mpi_extent(sendtype));
44 @ requires sendcount >= 0 && sendcount * \mpi_extent(sendtype) < 10;
45 @ requires recvcount >= 0 && recvcount * \mpi_extent(recvtype) < 10;
46 @ requires 0 <= root && root < \mpi_comm_size;
47 @ requires \mpi_valid(sendbuf, sendcount, sendtype);
48 @ behavior imroot:
49 @ assumes \mpi_comm_rank == root;
50 @ assigns \mpi_region(recvbuf, recvcount * \mpi_comm_size, recvtype);
51 @ requires \mpi_valid(recvbuf, recvcount * \mpi_comm_size, recvtype);
52 @ requires recvcount * \mpi_extent(recvtype) ==
53 @ sendcount * \mpi_extent(sendtype);
54 @ ensures \mpi_equals(\mpi_offset(recvbuf, root * sendcount, sendtype), sendcount, sendtype, sendbuf);
55 @ waitsfor (0 .. \mpi_comm_size-1);
56 @ behavior imnroot:
57 @ assumes \mpi_comm_rank != root;
58 @ ensures \mpi_equals(sendbuf, sendcount, sendtype,
59 @ \mpi_offset(\on(root, recvbuf), \mpi_comm_rank * sendcount, sendtype));
60 @*/
61int gather(void* sendbuf, int sendcount, MPI_Datatype sendtype,
62 void* recvbuf, int recvcount, MPI_Datatype recvtype,
63 int root, MPI_Comm comm){
64 int rank, nprocs;
65 MPI_Status status;
66 int tag = 998;
67
68 MPI_Comm_rank(comm, &rank);
69 MPI_Comm_size(comm, &nprocs);
70 if(root == rank) {
71 void *ptr;
72
73 ptr = $mpi_pointer_add(recvbuf, root * recvcount, recvtype);
74 $elaborate(recvcount);
75 memcpy(ptr, sendbuf, recvcount * sizeofDatatype(recvtype));
76 }else
77 MPI_Send(sendbuf, sendcount, sendtype, root, tag, comm);
78 if(rank == root) {
79 int real_recvcount;
80 int offset;
81
82 for(int i=0; i<nprocs; i++) {
83 if(i != root) {
84 void * ptr;
85
86 offset = i * recvcount;
87 ptr = $mpi_pointer_add(recvbuf, offset, recvtype);
88 MPI_Recv(ptr, recvcount, recvtype, i, tag, comm,
89 &status);
90 }
91 }
92 }
93 return 0;
94}
95
96/*@ \mpi_collective(comm, P2P):
97 @ requires \mpi_agree(sendcount * \mpi_extent(sendtype));
98 @ requires sendcount >= 0 && sendcount * \mpi_extent(sendtype) < 10;
99 @ requires recvcount >= 0 && recvcount * \mpi_extent(recvtype) < 10;
100 @ requires \mpi_valid(sendbuf, sendcount, sendtype);
101 @ requires \mpi_valid(recvbuf, recvcount * \mpi_comm_size, recvtype);
102 @ requires \mpi_extent(recvtype) * recvcount == \mpi_extent(sendtype) * sendcount;
103 @ //assigns \mpi_region(recvbuf, recvcount * \mpi_comm_size, recvtype);
104 @ ensures \mpi_agree(\mpi_region(recvbuf, recvcount * \mpi_comm_size, recvtype));
105 @ ensures \mpi_equals(sendbuf, sendcount, sendtype,
106 @ \mpi_offset(recvbuf, \mpi_comm_rank * recvcount, recvtype));
107 @
108 */
109int allgather(void *sendbuf, int sendcount, MPI_Datatype sendtype,
110 void *recvbuf, int recvcount, MPI_Datatype recvtype,
111 MPI_Comm comm){
112 int place;
113 int nprocs;
114
115 MPI_Comm_rank(comm, &place);
116 MPI_Comm_size(comm, &nprocs);
117 gather(sendbuf, sendcount, sendtype,
118 recvbuf, recvcount, recvtype,
119 0, comm);
120 broadcast(recvbuf, recvcount*nprocs, recvtype, 0, comm);
121 return 0;
122}
Note: See TracBrowser for help on using the repository browser.