CIVL Description and Demo

CIVL: Concurrency Intermediate Verification Language is a framework for the verification of multiple concurrent C dialects and libraries, such as OpenMP, MPI, Pthreads, CUDA, etc. It contains transformers to transform programs of such language or libraries into CIVL-C program, which is a dialect of C with concurrent extension, and uses symbolic execution to model check the programs. CIVL-C language is designed to capture the common features of concurrency so that could be used to describe a large number of concurrent dialects or libraries. CIVL is open source, and is designed to be flexible to be extended for verifying more concurrent languages and libraries.

Properties that CIVL checks include: deadlock, memory leak, array index-out-of bound, invalid pointer dereference, user-specified assertions, and functional equivalence.

CIVL is also able to verify hybrid concurrent programs, i.e., programs containing more than one dialects or librabries, such as MPI-OpenMP, MPI-Pthreads, etc.

The following demonstrate some of CIVL's capabilities. This is a small sample; see the CIVL manual for the full story.

Dining philosopher in CIVL-C

CIVL-C implementation of dining philosopher, which contains a deadlock. Suppose the file name is diningBad.cvl.

#include <civlc.cvh>

$input int B = 4; // upper bound on number of philosophers
$input int n;     // number of philosophers
$assume(2<=n && n<=B);

_Bool forks[n]; // Each fork will be on the table ($true) or in a hand ($false).

void dine(int id) {
  int left = id;
  int right = (id + 1) % n;

  while (1) {
    $when (forks[left]) forks[left] = $false;
    $when (forks[right]) forks[right] = $false;
    forks[right] = $true;
    forks[left] = $true;
  }
}

void main() { 
  $for(int i: 0 .. n-1)
    forks[i] = $true;
  $parfor(int i: 0 .. n-1)
    dine(i);
}

		

Command: civl verify -inputB=4 diningBad.cvl

The above command triggers CIVL to verify the CIVL-C program for dining philosophers, with the upper bound of number of philosophers being 4. As shown by the output below, CIVL detects a deadlock violation of the program at depth 79, i.e., after 78 trace steps.

Output:

CIVL v1.4 of 2015-08-11 -- http://vsl.cis.udel.edu/civl

Violation 0 encountered at depth 36:
CIVL execution violation (kind: DEADLOCK, certainty: PROVEABLE)
at diningBad.cvl:31.11-12 ";":
A deadlock is possible:
  Path condition: true
  Enabling predicate: false
process p0 (id=0): false
process p1 (id=1): false
process p2 (id=2): false

Call stacks:

process p0 (id=0):
  _CIVL_system at diningBad.cvl:31.11-12 ";"

process p1 (id=1):
  dine at diningBad.cvl:21.25-37 "forks[right]"

process p2 (id=2):
  dine at diningBad.cvl:21.25-37 "forks[right]"

Logging new entry 0, writing trace to CIVLREP/diningBad_0.trace
Terminating search after finding 1 violation.

=== Command ===
civl verify diningBad.cvl 

=== Stats ===
   time (s)            : 1.17
   memory (bytes)      : 128974848
   max process count   : 3
   states              : 70
   states saved        : 37
   state matches       : 1
   transitions         : 69
   trace steps         : 37
   valid calls         : 175
   provers             : cvc4, z3, cvc3
   prover calls        : 5

=== Result ===
The program MAY NOT be correct.  See CIVLREP/diningBad_log.txt
Now that we know there is a violation, we probably want to find a minimal one. We can use the option "-min" for the verification, as shown by the following command.

Command: civl verify -inputB=4 -min diningBad.cvl

As shown in the following output, CIVL keeps reducing the search depth until it finds the minimal counterexample at depth 23, which is much smaller and is the case when there are two philosophers.

Output:

CIVL v1.3+ of 2015-07-23 -- http://vsl.cis.udel.edu/civl

Violation 0 encountered at depth 79:

...

Logging new entry 0, writing trace to CIVLREP/diningBad_0.trace
Restricting search depth to 78

Violation 1 encountered at depth 74:

...

New log entry is equivalent to previously encountered entry 0
Length of new trace (74) is less than length of old (79): replacing old with new...
Restricting search depth to 73

...

Violation 12 encountered at depth 23:

...

New log entry is equivalent to previously encountered entry 0
Length of new trace (23) is less than length of old (28): replacing old with new...
Restricting search depth to 22

=== Command ===
civl verify -inputB=4 -min diningBad.cvl 

=== Stats ===
   time (s)            : 1.42
   memory (bytes)      : 128974848
   max process count   : 5
   states              : 243
   states saved        : 230
   state matches       : 155
   transitions         : 578
   trace steps         : 561
   valid calls         : 6033
   provers             : cvc4, z3, cvc3
   prover calls        : 9

=== Result ===
The program MAY NOT be correct.  See CIVLREP/diningBad_log.txt

Now CIVL finds the minimal counterexample and we can use CIVL to reproduce it via the "replay" command.

Command: civl replay -showTransitions diningBad.cvl

Output:

CIVL v1.3+ of 2015-07-23 -- http://vsl.cis.udel.edu/civl

Initial state:

State 0
...

Step 0: State 0, p0:
  0->1: B=4 at diningBad.cvl:9.0-16 "$input int B = 4"
  1->2: n=InitialValue(n) [n:=X0] at diningBad.cvl:10.0-12 "$input int n"
--> State 1

Step 1: State 1, p0:
  2->3: $assume((2<=X0)&&(X0<=4)) at diningBad.cvl:11.0-21 "$assume(2<=n && n ... )"
  3->4: forks=InitialValue(forks) [forks:=(lambda i : int . false)] at diningBad.cvl:13.0-14 "_Bool forks[n]"
--> State 2

Step 2: State 2, p0:
  4->11: $elaborate(X0) at diningBad.cvl:28.2-15 "$elaborate(n)"
--> State 3

Step 3: State 3, p0:
  11->12: i=0 at civlc.cvl:18.6-13 "int i=0"
--> State 4

Step 4: State 4, p0:
  12->13: LOOP_BODY_ENTER at civlc.cvl:18.15-18 "i<x"
--> State 5

Step 5: State 5, p0:
  13->12: i=0+1 [i:=1] at civlc.cvl:18.20-23 "i++"
--> State 6

Step 6: State 6, p0:
  12->13: LOOP_BODY_ENTER at civlc.cvl:18.15-18 "i<x"
--> State 7

Step 7: State 7, p0:
  13->12: i=1+1 [i:=2] at civlc.cvl:18.20-23 "i++"
--> State 8

Step 8: State 8, p0:
  12->14: LOOP_BODY_EXIT at civlc.cvl:18.15-18 "i<x"
  14->RET: $elaborate(...) return at civlc.cvl:18.24-25 ";"
--> State 10

Step 9: State 10, p0:
  5->6: LOOP_BODY_ENTER at diningBad.cvl:29.14-22 "0 .. n-1"
  6->7: $for((NULL) has next in ($domain){0, 1#1} at diningBad.cvl:29.2-6 "$for"
--> State 11

Step 10: State 11, p0:
  7->5: forks[0]=true at diningBad.cvl:30.4-20 "forks[i] = $true"
--> State 12

Step 11: State 12, p0:
  5->6: LOOP_BODY_ENTER at diningBad.cvl:29.14-22 "0 .. n-1"
  6->7: $for((0) has next in ($domain){0, 1#1} at diningBad.cvl:29.2-6 "$for"
--> State 13

Step 12: State 13, p0:
  7->5: forks[1]=true at diningBad.cvl:30.4-20 "forks[i] = $true"
--> State 14

Step 13: State 14, p0:
  5->8: LOOP_BODY_EXIT at diningBad.cvl:29.14-22 "0 .. n-1"
  8->9: $parfor(i0: ($domain){0, 1#1}) $spawn dine(i0) at diningBad.cvl:31.2-9 "$parfor"
--> State 15

Step 14: State 15, p1:
  15->16: left=0 at diningBad.cvl:16.2-15 "int left = id"
--> State 16

Step 15: State 16, p1:
  16->17: right=(0+1)%2 [right:=1] at diningBad.cvl:17.2-26 "int right = (id  ... n"
--> State 17

Step 16: State 17, p1:
  17->18: LOOP_BODY_ENTER at diningBad.cvl:19.9-10 "1"
--> State 18

Step 17: State 18, p2:
  15->16: left=1 at diningBad.cvl:16.2-15 "int left = id"
--> State 19

Step 18: State 19, p2:
  16->17: right=(1+1)%2 [right:=0] at diningBad.cvl:17.2-26 "int right = (id  ... n"
--> State 20

Step 19: State 20, p2:
  17->18: LOOP_BODY_ENTER at diningBad.cvl:19.9-10 "1"
--> State 21

Step 20: State 21, p1:
  18->19: forks[0]=false at diningBad.cvl:20.24-44 "forks[left] = $false"
--> State 22

Step 21: State 22, p2:
  18->19: forks[1]=false at diningBad.cvl:20.24-44 "forks[left] = $false"
--> State 23

State 23
| Path condition
| | true
| Dynamic scopes
| | dyscope d0 (id=0, parent=NULL, static=0)
| | | variables
| | | | B = 4
| | | | n = 2
| | | | forks = (boolean[2]){false, false}
| | | | _dom_size0 = 2
| | | | _par_procs0 = (process[2]){p1, p2}
| | | ...
| Process states
| | process p0(id=0)
| | | call stack
| | | | Frame[function=_CIVL_system, location=9, diningBad.cvl:32.11-12 ";", dyscope=d0]
| | process p1(id=1)
| | | call stack
| | | | Frame[function=dine, location=19, diningBad.cvl:21.25-37 "forks[right]", dyscope=d5]
| | process p2(id=2)
| | | call stack
| | | | Frame[function=dine, location=19, diningBad.cvl:21.25-37 "forks[right]", dyscope=d7]

Violation of Deadlock found in State 23:
A deadlock is possible:
  Path condition: true
  Enabling predicate: false
process p0 (id=0): false
process p1 (id=1): false
process p2 (id=2): false

Trace ends after 22 trace steps.
Violation(s) found.

=== Command ===
civl replay -showTransitions diningBad.cvl 

=== Stats ===
   time (s)            : 0.85
   memory (bytes)      : 128974848
   max process count   : 3
   states              : 29
   states Saved        : 24
   valid calls         : 159
   provers             : cvc4, z3, cvc3
   prover calls        : 5

1-dimensional wave simulator in MPI

MPI implementation of 1d-wave simulator, which has a bug. The file name of this program is wave1d.c. Some additional CIVL-C code has been inserted to specify the intended behavior of the program; in the full source, this extra code is within preprocessor directives #ifdef _CIVL ... #endif which cause the extra code to be ignored during normal compilation but used by CIVL. The program starts with the string in some arbitrary initial position specified by the input variable u_init . The variable oracle keeps the result of a sequential run of the simulator at every time step, and at every time step of the parallel run an assertion is checked to see if the parallel result agrees with the sequential one.

/* Input parameters */
#ifdef _CIVL
$input int NXB = 5;               /* upper bound on nx */
$input int nx;                    /* number of discrete points including endpoints */
$assume(2 <= nx && nx <= NXB);     /* setting bounds */
$input double c;                  /* physical constant to do with string */
$assume(c > 0.0);       
$input int NSTEPSB = 5;        
$input int nsteps;                /* number of time steps */
$assume(0 < nsteps && nsteps <= NSTEPSB);
$input int wstep = 1;             /* number of time steps between printing */
double oracle[nsteps + 1][nx + 2];/* result of sequential run at every time step */
$input double u_init[nx];         /* arbitrary initial position of string */
...
#else
...
#endif

/* Global varibales */
double *u_prev, *u_curr, *u_next;
double k;
int nprocs, nxl, rank;
int first;                       /* global index of first cell owned by this process */
int left, right;                 /* ranks of left neighbor and right neighbor */
...
/* Update cells owned by processes */
void update() {
  int i;
  double *tmp;

  for (i = 1; i < nxl + 1; i++){
    u_next[i] = 2.0*u_curr[i] - u_prev[i] +
      k*(u_curr[i+1] + u_curr[i-1] -2.0*u_curr[i]);
  }
  //cycle pointers:
  tmp = u_prev;
  u_prev = u_curr;
  u_curr = u_next;
  u_next = tmp;
}

/* Initialization function, initializes all parameters and data array.
 * In CIVL mode, process 0 runs the complete problem sequentially to
 * create the oracle which will be compared with the results of the
 * parallel run.
 */
void initialization() {
  ...
  nxl = countForProc(rank);
  first = firstForProc(rank);
  u_prev = (double *)malloc((nxl + 2) * sizeof(double));
  assert(u_prev);
  u_curr = (double *)malloc((nxl + 2) * sizeof(double));
  assert(u_curr);
  u_next = (double *)malloc((nxl + 2) * sizeof(double));
  assert(u_next);
  // sets constant boundaries
  u_prev[0] = 0;
  u_curr[0] = 0;
  //u_next[0] = 0;
  u_prev[nxl + 1] = 0;
  u_curr[nxl + 1] = 0;
  //u_next[nxl + 1] = 0;
  // skip processes which are allocated no cells
  if(first != 0 && nxl != 0)
    left = OWNER(first - 1);
  else 
    left = MPI_PROC_NULL;
  if(first + nxl < nx && nxl != 0)
    right = OWNER(first + nxl);
  else
    right = MPI_PROC_NULL;
}

...
/* Print out the value of data cells; 
   Do comparison in CIVL mode */
void printData (int time, int first, int length, double * buf) {
  for(int i = 0; i < length; i++){
    printf("u_curr[%d]=%8.8f   ", first + i, buf[i]);
#ifdef _CIVL

     $assert((oracle[time + 1][first + i + 1] == buf[i]), \
    "Error: disagreement at time %d position %d: saw %lf, expected %lf", \
    time, first + i, buf[i], oracle[time + 1][first + i + 1]);
#endif
  }
}
...

/* Exchanging ghost cells */
void exchange(){
  MPI_Sendrecv(&u_curr[1], 1, MPI_DOUBLE, left, FROMRIGHT, &u_curr[nxl+1], 1, MPI_DOUBLE, 
	       right, FROMRIGHT, MPI_COMM_WORLD, MPI_STATUS_IGNORE);
  MPI_Sendrecv(&u_curr[nxl], 1, MPI_DOUBLE, right, FROMLEFT, &u_curr[0], 1, 
	       MPI_DOUBLE, left, FROMLEFT, MPI_COMM_WORLD, MPI_STATUS_IGNORE);
}

int main(int argc, char * argv[]) {
  int iter;

  MPI_Init(&argc, &argv);
  MPI_Comm_rank(MPI_COMM_WORLD, &rank);
  MPI_Comm_size(MPI_COMM_WORLD, &nprocs); 
  initialization();
  ...
  for(iter = 0; iter < nsteps; iter++) {
    if(iter % wstep == 0)
      write_frame(iter);
    if(nxl != 0){
      exchange();
      update();
    }
  }
  free(u_curr);
  free(u_prev);
  free(u_next);
  MPI_Finalize(); 
  return 0;
}		
		
With the upper bound of the size of initial vector being 5, and the bound of number of steps of the computation being 5, and the number of processes is restricted to be in the range [1, 4], CIVL finds a violation of the assertion of the program in a few seconds. The violation is at depth 1082.

Command: civl verify -enablePrintf=false wave1d.c

Output:

CIVL v1.3+ of 2015-07-23 -- http://vsl.cis.udel.edu/civl
Error: disagreement at time 2 position 0: 
saw -4*(X5[1]*(X3^4))+X5[2]*X3^4+5*(X5[0]*(X3^4))+Hop1s2f2o0[0]*X3^2+3*(X5[1]*(X3^2))+-6*(X5[0]*(X3^2))+X5[0], 
expected -4*(X5[1]*(X3^4))+X5[2]*X3^4+5*(X5[0]*(X3^4))+3*(X5[1]*(X3^2))+-6*(X5[0]*(X3^2))+X5[0]

Violation 0 encountered at depth 1082:
CIVL execution violation in p1(id=1) (kind: ASSERTION_VIOLATION, certainty: PROVEABLE)
at wave1d.c:199.5-192 "$assert((oracle[time +  ... )":

Input variables:
...
NXB=5
nx=5
c=X3
NSTEPSB=5
nsteps=5
wstep=1
u_init=(real[5])X5
_NPROCS=4
_NPROCS_LOWER_BOUND=1
_NPROCS_UPPER_BOUND=4

Context: 0<X3 && 0<=SIZEOF_REAL+-1 && 0<=-1*X0+9 && 0<=X0+-1

Assertion: $assert((((oracle)[(time+1)])[((first+i)+1)]==*((buf+i))), &((_anon_11)[0]), time, (first+i), *((buf+i)), ((oracle)[(time+1)])[((first+i)+1)])
-> -4*(X5[1]*(X3^4))+X5[2]*X3^4+5*(X5[0]*(X3^4))+3*(X5[1]*(X3^2))+-6*(X5[0]*(X3^2))+X5[0]==*(&heap.malloc0[0][1]+0)
-> -4*(X5[1]*(X3^4))+X5[2]*X3^4+5*(X5[0]*(X3^4))+3*(X5[1]*(X3^2))+-6*(X5[0]*(X3^2))+X5[0]==-4*(X5[1]*(X3^4))+X5[2]*X3^4+5*(X5[0]*(X3^4))+Hop1s2f2o0[0]*X3^2+3*(X5[1]*(X3^2))+-6*(X5[0]*(X3^2))+X5[0]
-> 0==Hop1s2f2o0[0]*X3^2

Call stacks:

process p0 (id=0):
  _CIVL_system at MPItoCIVLTransformer:0.-1-20 "$parfor(int i: 0 ...." inserted by MPItoCIVLTransformer.$parfor MPI_Process before op.h:21.1-10 "Operation" included from civlc.cvh:11.9-15 "" included from civlc.cvl:8.9-20 ""

process p1 (id=1):
  printData at wave1d.c:199.5-12 "$assert" called from
  write_frame at wave1d.c:213.4-13 "printData" called from
  s_main at wave1d.c:266.6-17 "write_frame" called from
  MPI_Process at GeneralTransformer:0.-1-20 "CIVL_argc, &_argv[..." inserted by GeneralTransformer.new main function before wave1d.c:235.13-17 "argc"

process p2 (id=2):
  ...
  
process p3 (id=3):
  ...
  
process p4 (id=4):
  ...
  
Logging new entry 0, writing trace to CIVLREP/wave1d_0.trace
Terminating search after finding 1 violation.

=== Command ===
civl verify -enablePrintf=false wave1d.c 

=== Stats ===
   time (s)            : 6.95
   memory (bytes)      : 382730240
   max process count   : 5
   states              : 1692
   states saved        : 1153
   state matches       : 0
   transitions         : 1692
   trace steps         : 1081
   valid calls         : 5523
   provers             : cvc4, z3, cvc3
   prover calls        : 46

=== Result ===
The program MAY NOT be correct.  See CIVLREP/wave1d_log.txt
		
		
		
		
		
With "-min" option enabled, CIVL finds a minimal counterexample at depth 203 after 32 seconds. With the help of the counterexample, one can be guided to fix wave1d.c by adding statements to initialize u_next[0] and u_next[nxl+1] to zero, yielding wave1d_fix.c. Then CIVL is applied again to verify the fix. In fact, the erroneous wave1d.c had been used for two years in a parallel computing course taught by one of us (Dr. Siegel), but the defect was never noticed, apparently because the allocated memory was almost always initialized with 0s, something which is certainly not guaranteed by the C Standard.

Dot product in OpenMP

Now let's look at a simple OpenMP program dotProduct_critical.c that performs the dot product of two vectors. The program is originally used as an example program for a parallel programming course.

#include <omp.h>
...		
#define N 100

int main (int argc, char *argv[]) {
  double a[N], b[N];
  double localsum, sum = 0.0;
  int i, tid, nthreads;
  
#pragma omp parallel shared(a,b,sum) private(i, localsum)
  {
    /* Get thread number */
    tid = omp_get_thread_num();
    
    /* Only master thread does this */
#pragma omp master
    {
      nthreads = omp_get_num_threads();
      printf("Number of threads = %d\n", nthreads);
    }

    /* Initialization */
#pragma omp for
    for (i=0; i < N; i++)
      a[i] = b[i] = (double)i;

    localsum = 0.0;

    /* Compute the local sums of all products */
#pragma omp for
    for (i=0; i < N; i++)
      localsum = localsum + (a[i] * b[i]);

#pragma omp critical
    sum = sum+localsum;

  }  /* End of parallel region */

  printf("   Sum = %2.1f\n",sum);

#ifdef _CIVL
  assert(sum==328350);
#endif
  exit(0);
}
		

CIVL transforms an OpenMP program to be purely CIVL-C before it performs verification on it, during which an input variable THREAD_MAX is introduced as the upper bound of the number of threads. To verify an OpenMP program, one needs to specify the value of THREAD_MAX . Here, we use 3 as THREAD_MAX for verification.

Command: civl verify -min -inputTHREAD_MAX=3 -enablePrintf=false dotProduct_critical.c

Output:

CIVL v1.3+ of 2015-07-23 -- http://vsl.cis.udel.edu/civl
Thread 1 can't perform $omp_write because thread 0 has written to the same variable and hasn't flushed yet.


Violation 0 encountered at depth 1638:
CIVL execution violation in p3(id=2) (kind: ASSERTION_VIOLATION, certainty: PROVEABLE)
at civl-omp.cvl:239.2-240.26 "$assert(otherTid < 0,  ... )":

Input variables:
CIVL_argc=X0
CIVL_argv=(CHAR[10][])X1
THREAD_MAX=3

Context: 0<=SIZEOF(struct OMP_gteam)+-1 && 0<=SIZEOF(struct OMP_team)+-1 && 0<=SIZEOF(struct OMP_gshared)+-1 && 0<=SIZEOF(struct OMP_shared)+-1 && 0<=-1*X0+9 && 0<=X0+-1

Assertion: $assert((otherTid<0), &((_anon_4)[0]), (*(shared)).1, otherTid)
-> 0<0
-> false

Call stacks:

process p0 (id=0):
  ...
  
process p1 (id=1):
  ...
  
process p3 (id=2):
  $omp_write at civl-omp.cvl:239.2-9 "$assert" called from
  _par_proc0 at OpenMPtoCIVLTransformer:0.-1-20 "tid_shared, &tid_l..." inserted by OpenMPtoCIVLTransformer.tid_sharedWriteCall before op.h:21.1-10 "Operation" included from civlc.cvh:11.9-15 "" included from omp.cvl:1.9-20 ""
	
Logging new entry 0, writing trace to CIVLREP/dotProduct_critical_0.trace
Restricting search depth to 1637

...

Violation 202 encountered at depth 221:
CIVL execution violation in p3(id=2) (kind: ASSERTION_VIOLATION, certainty: PROVEABLE)

...

New log entry is equivalent to previously encountered entry 0
Length of new trace (221) is less than length of old (232): replacing old with new...
Restricting search depth to 220

=== Command ===
civl verify -min -inputTHREAD_MAX=3 -enablePrintf=false dotProduct_critical.c 

=== Stats ===
   time (s)            : 8.8
   memory (bytes)      : 648019968
   max process count   : 4
   states              : 27347
   states saved        : 12033
   state matches       : 3
   transitions         : 27755
   trace steps         : 12188
   valid calls         : 29899
   provers             : cvc4, z3, cvc3
   prover calls        : 7

=== Result ===
The program MAY NOT be correct.  See CIVLREP/dotProduct_critical_log.tx.
		
		
		
		
		
		
		

As shown by the above output, wihtin 10 seconds, CIVL dectecs a race condition and finds the minimal counterexample at depth 221 after detecting 202 violations, which is much smaller than the first violation which is at depth 1638. With the help of trace provided by CIVL replay for the minimal counterexample, one is guided to fix the error by adding the variable tid to the private list of the parallel construct, yielding dotProduct_critical_fix.c, which is verified by CIVL successfully.


VSL