#ifdef _CIVL #include #endif /* FEVS: A Functional Equivalence Verification Suite for High-Performance * Scientific Computing * * Copyright (C) 2009-2010, Stephen F. Siegel, Timothy K. Zirkel, * University of Delaware * * This program is free software; you can redistribute it and/or * modify it under the terms of the GNU General Public License as * published by the Free Software Foundation; either version 3 of the * License, or (at your option) any later version. * * This program is distributed in the hope that it will be useful, but * WITHOUT ANY WARRANTY; without even the implied warranty of * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU * General Public License for more details. * * You should have received a copy of the GNU General Public License * along with this program; if not, write to the Free Software * Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA * 02110-1301 USA. */ /* * diffusion1d.c: parallel 1d-diffusion simulation. */ #include #include #include #include #include #include #define OWNER(index) ((nprocs*(index+1)-1)/nx) #pragma CIVL $input int NX; /* Upper bound for NX */ #pragma CIVL $input int NXB; #pragma CIVL $input double K; #pragma CIVL $input int NSTEPS; /* Upper bound for NSTEPS */ #pragma CIVL $input int NSTEPSB; #pragma CIVL $input int WSTEP; #pragma CIVL $output double u_output[NX]; #pragma CIVL $assume(2 < NX && NX <= NXB); #pragma CIVL $assume(0 < NSTEPS && NSTEPS <= NSTEPSB); /* Parameters: These are defined at the beginning of the input file: * * nx = number of points in x direction, including endpoints * k = D*dt/(dx*dx) * nsteps = number of time steps * wstep = write frame every this many steps * * Compiling with the flag -DDEBUG will also cause the data to be written * to a sequence of plain text files. */ /* Global variables */ int nx = -1; /* number of discrete points including endpoints */ double k = -1; /* D*dt/(dx*dx) */ int nsteps = -1; /* number of time steps */ int wstep = -1; /* write frame every this many time steps */ // Initialize pointers with NULL so that these // pointers can be checked if it's malloced already. double *u = NULL; /* temperature function */ double *u_new = NULL; /* temp. used to update u */ int nprocs; /* number of processes */ int rank; /* the rank of this process */ int left; /* rank of left neighbor */ int right; /* rank of right neighbor on torus */ int nxl; /* horizontal extent of one process */ int first; /* global index for local index 0 */ int start; /* first local index to update */ int stop; /* last local index to update */ double *buf; /* temp. buffer used on proc 0 only */ int firstForProc(int rank) { return (rank*nx)/nprocs; } int countForProc(int rank) { int a; int b; a = firstForProc(rank+1); b = firstForProc(rank); return a-b; } // Added a argument of File pointer so that file can be closed before the // program exit. void quit(FILE * file) { printf("Input file must have format:\n\n"); printf("nx = \n"); printf("k = \n"); printf("nsteps = \n"); printf("wstep = \n"); printf(" ...\n\n"); printf("where there are nx doubles at the end.\n"); fflush(stdout); // Missing free(u) and fclose(file) probably is a bug if(u != NULL) free(u); fclose(file); // Following statements should be added in $exit() or added by CIVL automatically exit(1); } void readint(FILE *file, char *keyword, int *ptr) { char buf[101]; int value; int returnval; returnval = fscanf(file, "%100s", buf); if (returnval != 1) quit(file); // fscanf can only make the output argument symbolic, so the // following comparison will always enable 2 transitions. // Trivially, for verification of the computing results of sequential program, // we can ignore the transition of executing "quit()". For concurrent programs, // ,especially for mpi programs, any one of processes terminated exceptionally // will make the some rest processes go into a dead lock (e.g.MPI_Recv forever). #pragma CIVL $assume(strcmp(keyword, buf) == 0); if (strcmp(keyword, buf) != 0) quit(file); returnval = fscanf(file, "%10s", buf); if (returnval != 1) quit(file); #pragma CIVL $assume(strcmp("=", buf) == 0); if (strcmp("=", buf) != 0) quit(file); returnval = fscanf(file, "%d", ptr); if (returnval != 1) quit(file); } void readdouble(FILE *file, char *keyword, double *ptr) { char buf[101]; int value; int returnval; returnval = fscanf(file, "%100s", buf); if (returnval != 1) quit(file); #pragma CIVL $assume(strcmp(keyword, buf) == 0); if (strcmp(keyword, buf) != 0) quit(file); returnval = fscanf(file, "%10s", buf); if (returnval != 1) quit(file); #pragma CIVL $assume(strcmp("=", buf) == 0); if (strcmp("=", buf) != 0) quit(file); returnval = fscanf(file, "%lf", ptr); if (returnval != 1) quit(file); } /* init: initializes global variables. */ void init(char* infilename) { char keyword[101]; FILE *infile = fopen(infilename, "r"); int i, j; MPI_Comm_size(MPI_COMM_WORLD, &nprocs); MPI_Comm_rank(MPI_COMM_WORLD, &rank); if (rank == 0) { assert(infile); readint(infile, "nx", &nx); readdouble(infile, "k", &k); readint(infile, "nsteps", &nsteps); readint(infile, "wstep", &wstep); #pragma CIVL $assume(nx == NX); #pragma CIVL $assume(nsteps == NSTEPS); #pragma CIVL $assume(wstep == WSTEP); #pragma CIVL $assume(0.0 < k && k < 0.5); printf("Diffusion1d with nx=%d, k=%f, nsteps=%d, wstep=%d nprocs=%d\n", nx, k, nsteps, wstep, nprocs); } MPI_Bcast(&nx, 1, MPI_INT, 0, MPI_COMM_WORLD); MPI_Bcast(&k, 1, MPI_DOUBLE, 0, MPI_COMM_WORLD); MPI_Bcast(&nsteps, 1, MPI_INT, 0, MPI_COMM_WORLD); MPI_Bcast(&wstep, 1, MPI_INT, 0, MPI_COMM_WORLD); assert(k>0 && k<.5); assert(nx>=2); assert(nsteps>=1); // nxl: number actual points (incl. end-points) // nxl+2: size of array (incl. ghost cells) first = firstForProc(rank); nxl = countForProc(rank); if (first == 0 || nxl == 0) left = MPI_PROC_NULL; else left = OWNER(first-1); if (first+nxl >= nx || nxl == 0) right = MPI_PROC_NULL; else right = OWNER(first+nxl); u = (double*)malloc((nxl+2)*sizeof(double)); assert(u!=NULL); u_new = (double*)malloc((nxl+2)*sizeof(double)); assert(u_new); if (rank == 0) { buf = (double*)malloc((1+nx/nprocs)*sizeof(double)); for (i=1; i <= nxl; i++){ if (fscanf(infile, "%lf", u+i) != 1) quit(infile); } for (i=1; i < nprocs; i++){ for (j=0; j