source: CIVL/examples/compare/diffusion1d/diffusion1d_spec_revision.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.9 KB
Line 
1#ifdef _CIVL
2#include <civlc.cvh>
3#endif
4/* FEVS: A Functional Equivalence Verification Suite for High-Performance
5 * Scientific Computing
6 *
7 * Copyright (C) 2009-2010, Stephen F. Siegel, Timothy K. Zirkel,
8 * University of Delaware
9 *
10 * This program is free software; you can redistribute it and/or
11 * modify it under the terms of the GNU General Public License as
12 * published by the Free Software Foundation; either version 3 of the
13 * License, or (at your option) any later version.
14 *
15 * This program is distributed in the hope that it will be useful, but
16 * WITHOUT ANY WARRANTY; without even the implied warranty of
17 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
18 * General Public License for more details.
19 *
20 * You should have received a copy of the GNU General Public License
21 * along with this program; if not, write to the Free Software
22 * Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
23 * 02110-1301 USA.
24 */
25
26
27/* diffusion1d_seq.c: sequential version of 1d diffusion.
28 * The length of the rod is 1. The endpoints are frozen at the input
29 * temperature.
30 *
31 */
32#include <stdlib.h>
33#include <stdio.h>
34#include <assert.h>
35#include <math.h>
36#include <string.h>
37#pragma CIVL $input int NX;
38#pragma CIVL $input int NXB;
39#pragma CIVL $input double K;
40#pragma CIVL $input int NSTEPS;
41#pragma CIVL $input int NSTEPSB;
42#pragma CIVL $input int WSTEP;
43#pragma CIVL $output double u_output[NX];
44#pragma CIVL $assume(2 < NX && NX <= NXB);
45#pragma CIVL $assume(0 < NSTEPS && NSTEPS <= NSTEPSB);
46
47/* Parameters: These are defined at the beginning of the input file:
48 *
49 * nx = number of points in x direction, including endpoints
50 * k = D*dt/(dx*dx)
51 * nsteps = number of time steps
52 * wstep = write frame every this many steps
53 *
54 * Compiling with the flag -DDEBUG will also cause the data to be written
55 * to a sequence of plain text files.
56 */
57
58/* Global variables */
59int nx = -1; /* number of discrete points including endpoints */
60double k = -1; /* D*dt/(dx*dx) */
61int nsteps = -1; /* number of time steps */
62int wstep = -1; /* write frame every this many time steps */
63double *u = NULL; /* temperature function */
64
65void quit(FILE * file) {
66 printf("Input file must have format:\n\n");
67 printf("nx = <INTEGER>\n");
68 printf("k = <DOUBLE>\n");
69 printf("nsteps = <INTEGER>\n");
70 printf("wstep = <INTEGER>\n");
71 printf("<DOUBLE> <DOUBLE> ...\n\n");
72 printf("where there are nx doubles at the end.\n");
73 fflush(stdout);
74 // Missing free(u) and fclose(file) probably is a bug
75 if(u != NULL)
76 free(u);
77 fclose(file);
78 // Following statements should be added in $exit() or added by CIVL automatically
79 exit(0);
80}
81
82void readint(FILE *file, char *keyword, int *ptr) {
83 char buf[101];
84 int value;
85 int returnval;
86
87 returnval = fscanf(file, "%100s", buf);
88 if (returnval != 1) quit(file);
89 // fscanf can only make the output argument symbolic, so the
90 // following comparison will always enable 2 transitions.
91 // Trivially, for verification of the computing results of sequential program,
92 // we can ignore the transition of executing "quit()". For concurrent programs,
93 // ,especially for mpi programs, any one of processes terminated exceptionally
94 // will make the some rest processes go into a dead lock (e.g.MPI_Recv forever).
95#pragma CIVL $assume(strcmp(keyword, buf) == 0);
96 if (strcmp(keyword, buf) != 0) quit(file);
97 returnval = fscanf(file, "%10s", buf);
98 if (returnval != 1) quit(file);
99#pragma CIVL $assume(strcmp("=", buf) == 0);
100 if (strcmp("=", buf) != 0) quit(file);
101 returnval = fscanf(file, "%d", ptr);
102 if (returnval != 1) quit(file);
103}
104
105void readdouble(FILE *file, char *keyword, double *ptr) {
106 char buf[101];
107 int value;
108 int returnval;
109
110 returnval = fscanf(file, "%100s", buf);
111 if (returnval != 1) quit(file);
112#pragma CIVL $assume(strcmp(keyword, buf) == 0);
113 if (strcmp(keyword, buf) != 0) quit(file);
114 returnval = fscanf(file, "%10s", buf);
115 if (returnval != 1) quit(file);
116#pragma CIVL $assume(strcmp("=", buf) == 0);
117 if (strcmp("=", buf) != 0) quit(file);
118 returnval = fscanf(file, "%lf", ptr);
119 if (returnval != 1) quit(file);
120}
121
122
123/* init: initializes global variables. read nx, k, nsteps, u,
124 * from infile. Open GIF output file. */
125void init(char* infilename) {
126 char keyword[101];
127 FILE *infile = fopen(infilename, "r");
128 int i;
129
130 assert(infile);
131 readint(infile, "nx", &nx);
132 readdouble(infile, "k", &k);
133 readint(infile, "nsteps", &nsteps);
134 readint(infile, "wstep", &wstep);
135#pragma CIVL $assume(nx == NX);
136#pragma CIVL $assume(nsteps == NSTEPS);
137#pragma CIVL $assume(wstep == WSTEP);
138#pragma CIVL $assume(0.0 < k && k < 0.5);
139 printf("Diffusion1d with nx=%d, k=%f, nsteps=%d, wstep=%d\n",
140 nx, k, nsteps, wstep);
141 fflush(stdout);
142 assert(nx>=2);
143 assert(k>0 && k<.5);
144 assert(nsteps >= 1);
145 assert(wstep >= 1 && wstep <=nsteps);
146 u = (double*)malloc(nx*sizeof(double));
147 assert(u);
148 for (i = 0; i < nx; i++) {
149 if (fscanf(infile, "%lf", u+i) != 1) quit(infile);
150 }
151 // Missing fclose for FILE * infile
152 fclose(infile);
153}
154
155/* write_plain: write current data to plain text file and stdout */
156/* currently I ignore sprintf stuff and let write_palin only do assignment to ouput
157 arguments(which is the only meaningful thing this function does in a verification
158 view)*/
159void write_plain(int time) {
160 for (int i = 0; i < nx; i++)
161 u_output[i] = u[i];
162
163 return;
164}
165
166/* updates u for next time step. */
167void update() {
168 int i;
169 double u_new[nx];
170
171 for (i = 1; i < nx-1; i++)
172 u_new[i] = u[i] + k*(u[i+1] + u[i-1] -2*u[i]);
173 for (i = 1; i < nx-1; i++)
174 u[i] = u_new[i];
175}
176
177/* main: executes simulation, creates one output file for each time
178 * step */
179int main(int argc, char *argv[]) {
180 int iter;
181
182#pragma CIVL $assume((argc == 2));
183 assert(argc==2);
184 init(argv[1]);
185 write_plain(0);
186 for (iter = 1; iter <= nsteps; iter++) {
187 update();
188 if (iter%wstep==0) write_plain(iter);
189 }
190 free(u);
191 return 0;
192}
Note: See TracBrowser for help on using the repository browser.