VerifyThis 2019: GHC sort
solution to VerifyThis 2019 challenge 1
Auteurs: Quentin Garchery
Catégories: List Data Structure / Sorting Algorithms
Outils: Why3
Références: VerifyThis @ ETAPS 2019
see also the index (by topic, by tool, by reference, by year)
VerifyThis @ ETAPS 2019 competition Challenge 1: Monotonic Segments and GHC sort
Author: Quentin Garchery (LRI, Université Paris-Sud)
use int.Int use seq.Seq use seq.OfList use seq.FreeMonoid use array.Array use map.Occ use list.ListRich use bool.Bool use exn.Exn
PART A : Monotonic Segments
clone list.Sorted as StrictIncr with type t = int, predicate le = (<) clone list.RevSorted with type t = int, predicate le = (<=) use list.FoldRight let function eqb (b1 b2 : bool) : bool ensures { result <-> b1 = b2 } = andb (implb b1 b2) (implb b2 b1) type list_seq = { list : list int; ghost seq : seq int } invariant { seq = of_list (reverse list) } by { list = Nil; seq = Seq.empty }
Use a type invariant that ensures that the sequence represents the list. It is useful to get a nice specification when dealing with first and last elements of the list (begin-to-end property)
let constant nil : list_seq = { list = Nil; seq = Seq.empty } let extend a lseq ensures { result.list = Cons a lseq.list } = { list = Cons a lseq.list; seq = snoc lseq.seq a }
Compute the monotonic cutpoints of a sequence
let cutpoints (s : array int) requires { Array.length s > 0 } (* Begin-to-end property *) ensures { get result.seq 0 = 0 } ensures { get result.seq (Seq.length result.seq - 1) = Array.length s } (* Non-empty property *) ensures { length result.list >= 2 } (* Within bounds property *) ensures { forall z. mem z result.list -> 0 <= z <= Array.length s } (* Monotonic property *) ensures { forall k. 0 <= k < Seq.length result.seq - 1 -> let ck = get result.seq k in let ck1 = get result.seq (k+1) in (forall z1 z2. ck <= z1 < z2 < ck1 -> s[z1] < s[z2]) \/ (forall z1 z2. ck <= z1 < z2 < ck1 -> s[z1] >= s[z2]) } (* For the next part, we also need the cutpoints list to be decreasing *) ensures { forall i j. 0 <= i < j < Seq.length result.seq -> get result.seq i <= get result.seq j } = let n = s.length in let ref cut = extend 0 nil in let ref x = 0 in let ref y = 1 in let ref increasing = True in while y < n do variant { n - y } invariant { y = x + 1 } invariant { 0 < y <= n+1 } invariant { Seq.length cut.seq > 0 } invariant { get cut.seq 0 = 0 } invariant { get cut.seq (Seq.length cut.seq - 1) = x } invariant { forall z. mem z cut.list -> 0 <= z <= n } invariant { forall k. 0 <= k < Seq.length cut.seq - 1 -> let ck = get cut.seq k in let ck1 = get cut.seq (k+1) in (forall z1 z2. ck <= z1 < z2 < ck1 -> s[z1] < s[z2]) \/ (forall z1 z2. ck <= z1 < z2 < ck1 -> s[z1] >= s[z2])} invariant { forall i j. 0 <= i < j < Seq.length cut.seq -> get cut.seq i <= get cut.seq j } label StartLoop in increasing <- (s[x] < s[y]); while y < n && eqb (s[y-1] < s[y]) increasing do variant { n - y } invariant { y at StartLoop <= y <= n } invariant { (forall z1 z2. x <= z1 < z2 < y -> s[z1] < s[z2]) \/ (forall z1 z2. x <= z1 < z2 < y -> s[z1] >= s[z2]) } y <- y + 1 done; cut <- extend y cut; assert { get (cut.seq at StartLoop) (Seq.length cut.seq - 2) = x }; assert { forall k. 0 <= k < Seq.length cut.seq - 2 -> get cut.seq k at StartLoop = get cut.seq k /\ get cut.seq (k+1) at StartLoop = get cut.seq (k+1) }; x <- y; y <- x+1; done; label AfterLoop in if x < n then cut <- extend n cut; assert { get cut.seq (Seq.length cut.seq - 1) = n }; assert { forall k. 0 <= k < Seq.length cut.seq - 2 -> get cut.seq k at AfterLoop = get cut.seq k /\ get cut.seq (k+1) at AfterLoop = get cut.seq (k+1) }; cut
PART B : GHC Sort
lemma reverse_sorted_incr : forall l. Decr.sorted l -> Incr.sorted (reverse l) by Decr.sorted l /\ Incr.sorted Nil /\ compat l Nil let rec lemma lt_le_sorted (l : list int) variant { l } requires { StrictIncr.sorted l} ensures { Incr.sorted l } = match l with | Cons _ (Cons h2 t) -> lt_le_sorted (Cons h2 t) | _ -> () end let function order l requires { StrictIncr.sorted l \/ Decr.sorted l } ensures { Incr.sorted result } ensures { permut l result } = match l with | Nil -> l | Cons _ Nil -> l | Cons h1 (Cons h2 _) -> if h1 < h2 then (assert { Incr.sorted l by StrictIncr.sorted l }; l) else (assert { Decr.sorted l }; reverse l) end
Get an ordered list from a monotonic list
let rec list_from (a : array int) s e variant { e - s } requires { 0 <= s <= Array.length a } requires { 0 <= e <= Array.length a } requires { (forall z1 z2. s <= z1 < z2 < e -> a[z1] < a[z2]) \/ (forall z1 z2. s <= z1 < z2 < e -> a[z1] >= a[z2]) } ensures { forall x. num_occ x result = occ x a.elts s e } ensures { forall x. mem x result -> exists z. s <= z < e /\ a[z] = x } ensures { (forall z1 z2. s <= z1 < z2 < e -> a[z1] < a[z2]) -> StrictIncr.sorted result } ensures { (forall z1 z2. s <= z1 < z2 < e -> a[z1] >= a[z2]) -> Decr.sorted result } ensures { StrictIncr.sorted result \/ Decr.sorted result } = if s >= e then Nil else Cons a[s] (list_from a (s+1) e)
Get a monotonic list from two cutpoints in the array
let rec lemma occ_slice (a : array 'a) (c1 c2 c3 : int) variant { c3 - c2 } requires { 0 <= c1 <= c2 <= c3 <= Array.length a } ensures { forall x. occ x a.elts c1 c3 = occ x a.elts c1 c2 + occ x a.elts c2 c3 } = if c3 <> c2 then occ_slice a c1 c2 (c3-1) let rec sorted_lists (s: array int) (cutp : list_seq) requires { length cutp.list > 0 } (* This is where we need cutpoints to be sorted *) requires { forall x y. 0 <= x < y < Seq.length cutp.seq -> get cutp.seq x <= get cutp.seq y } requires { forall z. mem z cutp.list -> 0 <= z <= Array.length s } requires { forall k. 0 <= k < Seq.length cutp.seq - 1 -> let ck = get cutp.seq k in let ck1 = get cutp.seq (k+1) in (forall z1 z2. ck <= z1 < z2 < ck1 -> s[z1] < s[z2]) \/ (forall z1 z2. ck <= z1 < z2 < ck1 -> s[z1] >= s[z2]) } variant { cutp.list } ensures { forall l. mem l result -> Incr.sorted l } ensures { forall x. let first = get cutp.seq 0 in let last = get cutp.seq (Seq.length cutp.seq - 1) in num_occ x (fold_right (++) result Nil) = occ x s.elts first last } = let ls = cutp.list in let seq = cutp.seq in match ls with | Nil | Cons _ Nil -> Nil | Cons h1 (Cons h2 t) -> assert { let k = Seq.length cutp.seq - 2 in h2 = get cutp.seq k /\ h1 = get cutp.seq (k+1) }; let seqi = list_from s h2 h1 in let lseq = { list = Cons h2 t; seq = seq[0..(length seq - 1)] } in (* occ_slice s.elts (get seq 0) (get lseq.seq (Seq.length lseq.seq - 1)) *) (* (get seq (Seq.length seq - 1)); *) assert { Incr.sorted (order seqi) }; Cons (order seqi) (sorted_lists s lseq) end
Get sorted lists from the cutpoints
let rec merge l1 l2 variant { length l1 + length l2 } requires { Incr.sorted l1 } requires { Incr.sorted l2 } ensures { Incr.sorted result } ensures { permut result (l1 ++ l2) } = match l1, l2 with | Nil, l | l, Nil -> l | Cons h1 t1, Cons h2 t2 -> if h1 < h2 then (assert { forall x. mem x (t1 ++ l2) -> h1 <= x }; Cons h1 (merge t1 l2)) else (assert { forall x. mem x (l1 ++ t2) -> h2 <= x }; Cons h2 (merge l1 t2)) end
The merge of mergesort !
let rec merge_pair ls variant { length ls } requires { forall l. mem l ls -> Incr.sorted l } ensures { length result <= length ls } ensures { length ls > 1 -> 0 < length result < length ls } ensures { forall l. mem l result -> Incr.sorted l } ensures { permut (fold_right (++) result Nil) (fold_right (++) ls Nil) } = match ls with | Nil | Cons _ Nil -> ls | Cons l1 (Cons l2 r) -> Cons (merge l1 l2) (merge_pair r) end
Merge pair by pair for efficiency
let rec mergerec ls requires { forall l. mem l ls -> Incr.sorted l } variant { length ls } ensures { Incr.sorted result } ensures { permut result (fold_right (++) ls Nil) } = match ls with | Nil -> Nil | Cons l Nil -> l | Cons _ (Cons _ _) -> mergerec (merge_pair ls) end
Repeat previous merge
use seq.Occ as SO use seq.OfList
Show that the result of <mergerec> has the same length has the initial array
let rec find (seq : seq int) (v : int) (s e : int) variant { e - s } requires { 0 <= s < e <= Seq.length seq } requires { SO.occ v seq s e >= 1 } ensures { s <= result < e /\ get seq result = v } = if get seq s = v then s else find seq v (s+1) e
By induction, when increasing the size of the sub-array by one, we remove the new element in the corresponding array (use <find> to find the index to remove)
let rec lemma same_occs_same_lengths (a : array int) (seq : seq int) (s : int) variant { Array.length a - s } requires { 0 <= s <= Array.length a } requires { forall x. occ x a.elts s (Array.length a) = SO.occ x seq 0 (Seq.length seq) } ensures { Seq.length seq = Array.length a - s } = let na = Array.length a in let ns = Seq.length seq in if s = na then (assert { ns > 0 -> SO.occ (get seq 0) seq 0 ns > 0 }; ()) else (assert { ns = 0 -> false by SO.occ a[s] seq 0 ns = 0 }; let i = find seq a[s] 0 ns in let rem_as = seq[0..i] ++ seq[i+1..ns] in (* seq where a[s] is removed *) assert { forall x. (if x = a[s] then SO.occ x seq 0 ns = SO.occ x rem_as 0 (ns-1) + 1 else SO.occ x seq 0 ns = SO.occ x rem_as 0 (ns-1)) by seq == Seq.(++) seq[0..i] (Seq.(++) seq[i..i+1] seq[i+1..ns]) so SO.occ x seq 0 ns = SO.occ x seq[0..i] 0 i + SO.occ x seq[i..i+1] 0 1 + SO.occ x seq[i+1..ns] 0 (ns-i-1) so SO.occ x rem_as 0 (ns-1) = SO.occ x seq[0..i] 0 i + SO.occ x seq[i+1..ns] 0 (ns-i-1) }; same_occs_same_lengths a rem_as (s+1)) let rec lemma num_occ_seq_occ (l : list 'a) (x : 'a) variant { l } ensures { num_occ x l = SO.occ x (of_list l) 0 (length l) } = match l with | Nil -> () | Cons h t -> assert { of_list l = Seq.(++) (singleton h) (of_list t) /\ if x = h then SO.occ x (singleton h) 0 1 = 1 else SO.occ x (singleton h) 0 1 = 0 }; num_occ_seq_occ t x end let sort_to_list a = requires { Array.length a > 0 } ensures { Incr.sorted result } ensures { forall x. occ x a.elts 0 (Array.length a) = num_occ x result } ensures { length result = Array.length a } let res = mergerec (sorted_lists a (cutpoints a)) in same_occs_same_lengths a (of_list res) 0; res use array.IntArraySorted use array.ArrayPermut use option.Option use list.NthNoOpt let rec copy_list (l : list int) (a : array int) (s : int) : unit variant { l } requires { s >= 0 } requires { length l = Array.length a - s } ensures { forall x. s <= x < Array.length a -> a[x] = nth (x - s) l } ensures { forall x. 0 <= x < s -> a[x] = (old a)[x] } ensures { forall x. occ x a.elts s (Array.length a) = num_occ x l } = match l with | Nil -> () | Cons h t -> a[s] <- h; copy_list t a (s+1) end
Copy a list in an array, element by element, starting at a given index
let rec lemma mem_nth_in_bounds (l : list 'a) (j : int) requires { 0 <= j < length l } ensures { mem (nth j l) l } = match l with | Nil -> () | Cons _ t -> if j > 0 then mem_nth_in_bounds t (j-1) end let rec lemma sorted_list_nth (l : list int) (i j : int) variant { l } requires { Incr.sorted l } requires { 0 <= i <= j < length l } ensures { nth i l <= nth j l } = match l with | Nil -> () | Cons _ t -> if i > 0 then sorted_list_nth t (i-1) (j-1) end
Useful to deduce sorted on arrays from sorted on lists
let ghc_sort a ensures { sorted a } ensures { permut_all a (old a) } = if Array.length a = 0 then () else let l = sort_to_list a in assert { length l = Array.length a }; copy_list l a 0
download ZIP archive
Why3 Proof Results for Project "verifythis_2019_ghc_sort"
Theory "verifythis_2019_ghc_sort.Top": fully verified
Obligations | Alt-Ergo 2.2.0 | Alt-Ergo 2.3.3 | Alt-Ergo 2.4.1 | CVC3 2.4.1 | CVC5 1.0.5 | |||
StrictIncr.Transitive.Trans | --- | 0.01 | --- | --- | --- | |||
RevSorted.Transitive.Trans | --- | 0.01 | --- | --- | --- | |||
VC for eqb | --- | 0.03 | --- | --- | --- | |||
VC for list_seq | --- | 0.02 | --- | --- | --- | |||
VC for nil | --- | 0.02 | --- | --- | --- | |||
VC for extend | --- | 0.06 | --- | --- | --- | |||
VC for cutpoints | --- | --- | --- | --- | --- | |||
split_vc | ||||||||
loop invariant init | --- | 0.01 | --- | --- | --- | |||
loop invariant init | --- | 0.02 | --- | --- | --- | |||
loop invariant init | --- | 0.02 | --- | --- | --- | |||
loop invariant init | --- | 0.05 | --- | --- | --- | |||
loop invariant init | --- | 0.03 | --- | --- | --- | |||
loop invariant init | --- | 0.03 | --- | --- | --- | |||
loop invariant init | --- | 0.03 | --- | --- | --- | |||
loop invariant init | --- | 0.05 | --- | --- | --- | |||
index in array bounds | --- | 0.02 | --- | --- | --- | |||
index in array bounds | --- | 0.02 | --- | --- | --- | |||
loop invariant init | --- | 0.02 | --- | --- | --- | |||
loop invariant init | --- | 0.02 | --- | --- | --- | |||
index in array bounds | --- | 0.02 | --- | --- | --- | |||
index in array bounds | --- | 0.02 | --- | --- | --- | |||
loop variant decrease | --- | 0.02 | --- | --- | --- | |||
loop invariant preservation | --- | 0.02 | --- | --- | --- | |||
loop invariant preservation | --- | 0.12 | --- | --- | --- | |||
assertion | --- | 0.04 | --- | --- | --- | |||
assertion | --- | --- | --- | --- | 5.47 | |||
loop variant decrease | --- | 0.02 | --- | --- | --- | |||
loop invariant preservation | --- | 0.02 | --- | --- | --- | |||
loop invariant preservation | --- | 0.02 | --- | --- | --- | |||
loop invariant preservation | --- | 0.05 | --- | --- | --- | |||
loop invariant preservation | --- | --- | --- | --- | 0.42 | |||
loop invariant preservation | --- | 0.40 | --- | --- | --- | |||
loop invariant preservation | --- | 0.09 | --- | --- | --- | |||
loop invariant preservation | 2.72 | --- | --- | --- | --- | |||
remove zero,one,(-),(>),(<=),(>=),get,set,([]'),([<-]'),andb,orb,notb,xorb,implb,is_nil,mem,is_none,nth,hd,tl,(++),reverse,rev_append,num_occ,permut,fold_right,([]''),singleton,cons,snoc,(++'),of_list,point_wise,([]),sorted2,ge,sorted1,sorted,compat,eqb,n,Assoc1,Unit_def_l,Unit_def_r,Inv_def_l,Inv_def_r,Comm1,Assoc,Mul_distr_l,Mul_distr_r,Comm,Unitary,NonTrivialRing,Refl,Trans2,Antisymm,Total,ZeroLessOne,CompatOrderAdd,CompatOrderMult,occ_empty,occ_right_no_add,occ_right_add,occ_left_no_add,occ_left_add,occ_bounds,occ_append,occ_neq,occ_exists,occ_pos,occ_eq,occ_exchange,is_nil'spec,Length_nonnegative,Length_nil,is_none'spec,Nth_tl,Nth0_head,Append_assoc,Append_l_nil,Append_length,mem_append,mem_decomp,reverse_append,reverse_cons,cons_reverse,reverse_reverse,reverse_mem,rev_append_append_l,rev_append_length,rev_append_def,rev_append_append_r,Num_Occ_NonNeg,Mem_Num_Occ,Append_Num_Occ,reverse_num_occ,Permut_refl,Permut_sym,Permut_trans,Permut_cons,Permut_swap,Permut_cons_append,Permut_assoc,Permut_append,Permut_append_swap,Permut_mem,Permut_length,fold_right_append,length_nonnegative,(==)'spec,create'spec,empty'def,set'spec,set'def,([<-])'def,singleton'spec,cons'spec,snoc'spec,([..])'spec,([..])'def,([_..])'def,([.._])'def,(++)'spec,associative,left_neutral,right_neutral,cons_def,snoc_def,double_sub_sequence,cons_back,snoc_back,cat_back,cons_dec,snoc_dec,cat_dec,empty_dec,singleton_dec,to_list_empty,to_list_cons,to_list_length,to_list_nth,to_list_def_cons,elts_seq_of_list,is_of_list,of_list_app,of_list_app_length,of_list_snoc,convolution_to_of_list,array'invariant,([<-])'spec,make_spec,Trans1,sorted_mem2,sorted_append2,Trans,sorted_mem1,sorted_append1,sorted_mem,sorted_append,rev_append_sorted_incr,rev_append_sorted_decr,eqb'spec,nil'def,Requires,Ensures1,LoopInvariant13,H11,LoopInvariant12,LoopInvariant10,LoopInvariant9,LoopInvariant7,H9,H8,H7,H6,H5,H4,H2,LoopInvariant5,LoopInvariant4,LoopInvariant3,LoopInvariant | ||||||||
loop invariant preservation | 0.21 | --- | --- | --- | --- | |||
loop invariant preservation | --- | 1.58 | --- | --- | --- | |||
assertion | --- | 0.50 | --- | --- | --- | |||
assertion | --- | --- | --- | --- | 4.84 | |||
postcondition | --- | --- | --- | --- | 0.36 | |||
postcondition | --- | 0.02 | --- | --- | --- | |||
postcondition | --- | 0.06 | --- | --- | --- | |||
postcondition | --- | 0.07 | --- | --- | --- | |||
postcondition | --- | --- | 4.52 | --- | --- | |||
postcondition | --- | 0.31 | --- | --- | --- | |||
assertion | --- | 0.02 | --- | --- | --- | |||
assertion | --- | 0.02 | --- | --- | --- | |||
postcondition | --- | 0.02 | --- | --- | --- | |||
postcondition | --- | 0.02 | --- | --- | --- | |||
postcondition | --- | 0.04 | --- | --- | --- | |||
postcondition | --- | 0.03 | --- | --- | --- | |||
postcondition | --- | 0.08 | --- | --- | --- | |||
postcondition | --- | 0.03 | --- | --- | --- | |||
reverse_sorted_incr | --- | --- | --- | --- | --- | |||
split_vc | ||||||||
reverse_sorted_incr.0 | --- | 0.02 | --- | --- | --- | |||
reverse_sorted_incr.1 | --- | 0.02 | --- | --- | --- | |||
reverse_sorted_incr.2 | --- | 0.02 | --- | --- | --- | |||
reverse_sorted_incr.3 | --- | 0.03 | --- | --- | --- | |||
VC for lt_le_sorted | --- | 0.06 | --- | --- | --- | |||
VC for order | --- | 0.20 | --- | --- | --- | |||
VC for list_from | --- | --- | --- | --- | --- | |||
split_vc | ||||||||
variant decrease | --- | 0.02 | --- | --- | --- | |||
precondition | --- | 0.02 | --- | --- | --- | |||
precondition | --- | 0.02 | --- | --- | --- | |||
precondition | --- | 0.03 | --- | --- | --- | |||
index in array bounds | --- | 0.01 | --- | --- | --- | |||
postcondition | --- | 0.59 | --- | --- | --- | |||
postcondition | --- | 0.06 | --- | --- | --- | |||
postcondition | --- | 0.18 | --- | --- | --- | |||
postcondition | --- | 0.43 | --- | --- | --- | |||
postcondition | --- | 0.06 | --- | --- | --- | |||
VC for occ_slice | --- | --- | --- | --- | --- | |||
split_vc | ||||||||
variant decrease | --- | 0.02 | --- | --- | --- | |||
precondition | --- | 0.02 | --- | --- | --- | |||
postcondition | --- | 0.03 | --- | --- | --- | |||
VC for sorted_lists | --- | --- | --- | --- | --- | |||
split_vc | ||||||||
assertion | --- | 1.42 | --- | --- | --- | |||
precondition | --- | --- | --- | 0.17 | --- | |||
remove zero,one,(-),(>),(<=),(>=),get,set,([]'),([<-]'),andb,orb,notb,xorb,implb,is_nil,length1,is_none,nth,hd,tl,(++),reverse,rev_append,num_occ,permut,fold_right,([]''),singleton,cons,snoc,(++'),of_list,point_wise,([]),sorted2,ge,sorted1,sorted,compat,eqb,seq,Assoc1,Unit_def_l,Unit_def_r,Inv_def_l,Inv_def_r,Comm1,Assoc,Mul_distr_l,Mul_distr_r,Comm,Unitary,NonTrivialRing,Refl,Trans2,Antisymm,Total,ZeroLessOne,CompatOrderAdd,CompatOrderMult,occ_empty,occ_right_no_add,occ_right_add,occ_left_no_add,occ_left_add,occ_bounds,occ_append,occ_neq,occ_exists,occ_pos,occ_eq,occ_exchange,is_nil'spec,Length_nonnegative,Length_nil,is_none'spec,Nth_tl,Nth0_head,Append_assoc,Append_l_nil,Append_length,mem_append,mem_decomp,reverse_append,reverse_cons,cons_reverse,reverse_reverse,reverse_mem,Reverse_length,rev_append_append_l,rev_append_length,rev_append_def,rev_append_append_r,Num_Occ_NonNeg,Mem_Num_Occ,Append_Num_Occ,reverse_num_occ,Permut_refl,Permut_sym,Permut_trans,Permut_cons,Permut_swap,Permut_cons_append,Permut_assoc,Permut_append,Permut_append_swap,Permut_mem,Permut_length,fold_right_append,length_nonnegative,(==)'spec,create'spec,empty'def,set'spec,set'def,([<-])'def,singleton'spec,cons'spec,snoc'spec,([..])'spec,([..])'def,([_..])'def,([.._])'def,(++)'spec,associative,left_neutral,right_neutral,cons_def,snoc_def,double_sub_sequence,cons_back,snoc_back,cat_back,cons_dec,snoc_dec,cat_dec,empty_dec,singleton_dec,to_list_empty,to_list_cons,to_list_length,to_list_nth,to_list_def_cons,length_of_list,elts_seq_of_list,is_of_list,of_list_app,of_list_app_length,of_list_snoc,convolution_to_of_list,array'invariant,([<-])'spec,make_spec,Trans1,sorted_mem2,sorted_append2,Trans,sorted_mem1,sorted_append1,sorted_mem,sorted_append,rev_append_sorted_incr,rev_append_sorted_decr,eqb'spec,list_seq'invariant,nil'def,reverse_sorted_incr,lt_le_sorted,order'spec,order'def,occ_slice,Requires3,Requires2,Requires,Assert | ||||||||
precondition | --- | 0.02 | --- | 0.02 | --- | |||
precondition | --- | --- | --- | 0.16 | --- | |||
remove zero,one,(-),(>),(<=),(>=),get,set,([]'),([<-]'),andb,orb,notb,xorb,implb,is_nil,length1,mem,is_none,nth,tl,(++),reverse,rev_append,permut,fold_right,([]''),singleton,cons,snoc,(++'),of_list,point_wise,([]),sorted2,ge,sorted1,sorted,compat,eqb,seq,Assoc1,Unit_def_l,Unit_def_r,Inv_def_l,Inv_def_r,Comm1,Assoc,Mul_distr_l,Mul_distr_r,Comm,Unitary,NonTrivialRing,Refl,Trans2,Antisymm,Total,ZeroLessOne,CompatOrderAdd,CompatOrderMult,occ_empty,occ_right_no_add,occ_right_add,occ_left_no_add,occ_left_add,occ_bounds,occ_append,occ_neq,occ_exists,occ_pos,occ_eq,is_nil'spec,Length_nonnegative,Length_nil,is_none'spec,Nth_tl,Nth0_head,Append_assoc,Append_l_nil,Append_length,mem_append,mem_decomp,reverse_append,reverse_cons,cons_reverse,reverse_reverse,reverse_mem,Reverse_length,rev_append_append_l,rev_append_length,rev_append_def,rev_append_append_r,Append_Num_Occ,reverse_num_occ,Permut_refl,Permut_sym,Permut_trans,Permut_cons,Permut_swap,Permut_cons_append,Permut_assoc,Permut_append,Permut_append_swap,Permut_mem,Permut_length,fold_right_append,length_nonnegative,(==)'spec,create'spec,empty'def,set'spec,set'def,([<-])'def,singleton'spec,cons'spec,snoc'spec,([..])'spec,([..])'def,([_..])'def,([.._])'def,(++)'spec,associative,left_neutral,right_neutral,cons_def,snoc_def,double_sub_sequence,cons_back,snoc_back,cat_back,cons_dec,snoc_dec,cat_dec,empty_dec,singleton_dec,to_list_empty,to_list_cons,to_list_length,to_list_nth,to_list_def_cons,length_of_list,elts_seq_of_list,is_of_list,of_list_app,of_list_app_length,of_list_snoc,convolution_to_of_list,array'invariant,([<-])'spec,make_spec,Trans1,sorted_mem2,sorted_append2,Trans,sorted_mem1,sorted_append1,sorted_mem,sorted_append,rev_append_sorted_incr,rev_append_sorted_decr,eqb'spec,list_seq'invariant,nil'def,reverse_sorted_incr,lt_le_sorted,order'spec,order'def,occ_slice,Requires4,Requires3,Requires1,Assert,Requires | ||||||||
precondition | --- | 0.03 | --- | 0.07 | --- | |||
precondition | --- | 0.08 | --- | --- | --- | |||
precondition | --- | 0.04 | --- | --- | --- | |||
precondition | --- | 1.33 | --- | --- | --- | |||
assertion | --- | 0.03 | --- | --- | --- | |||
variant decrease | --- | 0.06 | --- | --- | --- | |||
precondition | --- | 0.04 | --- | --- | --- | |||
precondition | --- | 0.08 | --- | --- | --- | |||
precondition | --- | --- | --- | --- | 1.04 | |||
precondition | --- | 0.51 | --- | --- | --- | |||
precondition | --- | 0.02 | --- | --- | --- | |||
postcondition | --- | --- | --- | --- | --- | |||
split_vc | ||||||||
postcondition | --- | 0.02 | --- | --- | --- | |||
postcondition | --- | 0.03 | --- | --- | --- | |||
postcondition | --- | 0.02 | --- | --- | --- | |||
postcondition | --- | --- | --- | --- | --- | |||
split_vc | ||||||||
postcondition | --- | 0.03 | --- | --- | --- | |||
postcondition | --- | 0.05 | --- | --- | --- | |||
postcondition | --- | --- | --- | --- | --- | |||
assert (occ x (elts s) (get1 (seq1 cutp) 0) (get1 (seq1 cutp) (length2 (seq1 cutp) - 1)) = occ x (elts s) (get1 (seq1 cutp) 0) (get1 (seq1 lseq) (length2 (seq1 lseq) - 1)) + occ x (elts s) (get1 (seq1 lseq) (length2 (seq1 lseq) - 1)) (get1 (seq1 cutp) (length2 (seq1 cutp) - 1))) | ||||||||
asserted formula | --- | 0.18 | --- | --- | --- | |||
postcondition | --- | --- | 0.69 | --- | --- | |||
VC for merge | --- | --- | --- | --- | --- | |||
split_all_full | ||||||||
assertion | --- | 0.49 | --- | --- | --- | |||
variant decrease | --- | 0.03 | --- | --- | --- | |||
precondition | --- | 0.02 | --- | --- | --- | |||
precondition | --- | 0.02 | --- | --- | --- | |||
assertion | --- | 0.15 | --- | --- | --- | |||
variant decrease | --- | 0.03 | --- | --- | --- | |||
precondition | --- | 0.01 | --- | --- | --- | |||
precondition | --- | 0.02 | --- | --- | --- | |||
postcondition | --- | 0.02 | --- | --- | --- | |||
postcondition | --- | 0.03 | --- | --- | --- | |||
postcondition | --- | 0.01 | --- | --- | --- | |||
postcondition | --- | 0.03 | --- | --- | --- | |||
postcondition | --- | 0.02 | --- | --- | --- | |||
postcondition | --- | 0.02 | --- | --- | --- | |||
postcondition | --- | 0.08 | --- | --- | --- | |||
postcondition | --- | 0.06 | --- | --- | --- | |||
postcondition | --- | 0.10 | --- | --- | --- | |||
postcondition | --- | 0.17 | --- | --- | --- | |||
VC for merge_pair | --- | --- | --- | --- | --- | |||
split_all_full | ||||||||
variant decrease | --- | 0.02 | --- | --- | --- | |||
precondition | --- | 0.05 | --- | --- | --- | |||
precondition | --- | 0.04 | --- | --- | --- | |||
precondition | --- | 0.05 | --- | --- | --- | |||
postcondition | --- | 0.02 | --- | --- | --- | |||
postcondition | --- | 0.02 | --- | --- | --- | |||
postcondition | --- | 0.03 | --- | --- | --- | |||
postcondition | --- | 0.02 | --- | --- | --- | |||
postcondition | --- | 0.02 | --- | --- | --- | |||
postcondition | --- | 0.02 | --- | --- | --- | |||
postcondition | --- | 0.02 | --- | --- | --- | |||
postcondition | --- | 0.02 | --- | --- | --- | |||
postcondition | --- | 0.04 | --- | --- | --- | |||
postcondition | --- | 0.04 | --- | --- | --- | |||
postcondition | --- | 0.03 | --- | --- | --- | |||
postcondition | --- | 0.68 | --- | --- | --- | |||
VC for mergerec | --- | 0.16 | --- | --- | --- | |||
VC for find | --- | 0.10 | --- | --- | --- | |||
VC for same_occs_same_lengths | --- | --- | --- | --- | --- | |||
split_vc | ||||||||
assertion | --- | 0.11 | --- | --- | --- | |||
assertion | --- | --- | --- | --- | --- | |||
split_vc | ||||||||
assertion | --- | 0.04 | --- | --- | --- | |||
VC for same_occs_same_lengths | --- | 0.04 | --- | --- | --- | |||
index in array bounds | --- | 0.02 | --- | --- | --- | |||
precondition | --- | 0.02 | --- | --- | --- | |||
precondition | --- | 0.03 | --- | --- | --- | |||
precondition | --- | 0.02 | --- | --- | --- | |||
precondition | --- | 0.02 | --- | --- | --- | |||
assertion | --- | --- | --- | --- | --- | |||
split_vc | ||||||||
assertion | --- | 0.23 | --- | --- | --- | |||
VC for same_occs_same_lengths | --- | 0.29 | --- | --- | --- | |||
VC for same_occs_same_lengths | --- | 0.06 | --- | --- | --- | |||
VC for same_occs_same_lengths | --- | 0.09 | --- | --- | --- | |||
VC for same_occs_same_lengths | --- | 0.08 | --- | --- | --- | |||
variant decrease | --- | 0.02 | --- | --- | --- | |||
precondition | --- | 0.02 | --- | --- | --- | |||
precondition | --- | 0.09 | --- | --- | --- | |||
postcondition | --- | 0.04 | --- | --- | --- | |||
VC for num_occ_seq_occ | --- | --- | --- | --- | --- | |||
split_all_full | ||||||||
assertion | --- | 0.11 | --- | --- | --- | |||
variant decrease | --- | 0.04 | --- | --- | --- | |||
postcondition | --- | 0.04 | --- | --- | --- | |||
postcondition | --- | 0.40 | --- | --- | --- | |||
VC for sort_to_list | --- | --- | --- | --- | --- | |||
split_vc | ||||||||
precondition | --- | 0.02 | --- | --- | --- | |||
precondition | --- | 0.02 | --- | --- | --- | |||
precondition | --- | 0.04 | --- | --- | --- | |||
precondition | --- | 0.03 | --- | --- | --- | |||
precondition | --- | 0.07 | --- | --- | --- | |||
precondition | --- | 0.03 | --- | --- | --- | |||
precondition | --- | 0.02 | --- | --- | --- | |||
precondition | --- | 0.06 | --- | --- | --- | |||
postcondition | --- | 0.02 | --- | --- | --- | |||
postcondition | --- | 0.06 | --- | --- | --- | |||
postcondition | --- | 0.03 | --- | --- | --- | |||
VC for copy_list | --- | --- | --- | --- | --- | |||
split_vc | ||||||||
postcondition | --- | 0.03 | --- | --- | --- | |||
postcondition | --- | 0.02 | --- | --- | --- | |||
postcondition | --- | 0.03 | --- | --- | --- | |||
index in array bounds | --- | 0.03 | --- | --- | --- | |||
variant decrease | --- | 0.04 | --- | --- | --- | |||
precondition | --- | 0.02 | --- | --- | --- | |||
precondition | --- | 0.03 | --- | --- | --- | |||
postcondition | --- | 0.04 | --- | --- | --- | |||
postcondition | --- | 0.04 | --- | --- | --- | |||
postcondition | --- | 0.11 | --- | --- | --- | |||
VC for mem_nth_in_bounds | --- | 0.04 | --- | --- | --- | |||
VC for sorted_list_nth | --- | --- | --- | --- | --- | |||
split_vc | ||||||||
variant decrease | --- | 0.04 | --- | --- | --- | |||
precondition | --- | 0.03 | --- | --- | --- | |||
precondition | --- | 0.03 | --- | --- | --- | |||
postcondition | --- | 0.12 | --- | --- | --- | |||
VC for ghc_sort | --- | --- | --- | --- | --- | |||
split_vc | ||||||||
postcondition | --- | 0.03 | --- | --- | --- | |||
postcondition | --- | 0.03 | --- | --- | --- | |||
precondition | --- | 0.03 | --- | --- | --- | |||
assertion | --- | 0.02 | --- | --- | --- | |||
precondition | --- | 0.02 | --- | --- | --- | |||
precondition | --- | 0.02 | --- | --- | --- | |||
postcondition | --- | 0.12 | --- | --- | --- | |||
postcondition | --- | 0.04 | --- | --- | --- |