Rework memoMap to squash space leaks.

The main benefit of this reorganization is that it
notices when a memoized `SeqMap` has been forced
at all of its locations.  This allows us to
discard the underlying computation, which will
never need to be consulted again.  This, in turn,
allows the garbage collector to reclaim the associated
memory and help prevent certain classes of space leaks.
This commit is contained in:
Rob Dockins 2021-04-13 16:14:57 -07:00
parent 29c5bc9224
commit 3f66dbf713
7 changed files with 101 additions and 76 deletions

View File

@ -62,24 +62,43 @@ import Cryptol.Backend.Concrete (Concrete)
import Cryptol.Backend.Monad (Unsupported(..))
import Cryptol.TypeCheck.Solver.InfNat(Nat'(..))
import Cryptol.Utils.Panic
-- | A sequence map represents a mapping from nonnegative integer indices
-- to values. These are used to represent both finite and infinite sequences.
data SeqMap sym a
= IndexSeqMap !(Integer -> SEval sym a)
| UpdateSeqMap !(Map Integer (SEval sym a))
!(Integer -> SEval sym a)
!(SeqMap sym a)
| MemoSeqMap
!Nat'
!(IORef (Map Integer a))
!(IORef (Integer -> SEval sym a))
indexSeqMap :: (Integer -> SEval sym a) -> SeqMap sym a
indexSeqMap = IndexSeqMap
lookupSeqMap :: SeqMap sym a -> Integer -> SEval sym a
lookupSeqMap :: Backend sym => SeqMap sym a -> Integer -> SEval sym a
lookupSeqMap (IndexSeqMap f) i = f i
lookupSeqMap (UpdateSeqMap m f) i =
lookupSeqMap (UpdateSeqMap m xs) i =
case Map.lookup i m of
Just x -> x
Nothing -> f i
Nothing -> lookupSeqMap xs i
lookupSeqMap (MemoSeqMap sz cache eval) i =
do mz <- liftIO (Map.lookup i <$> readIORef cache)
case mz of
Just z -> return z
Nothing ->
do f <- liftIO (readIORef eval)
v <- f i
msz <- liftIO $ atomicModifyIORef' cache (\m ->
let m' = Map.insert i v m in (m', Map.size m'))
-- If we memoize the entire map, overwrite the evaluation closure to let
-- the garbage collector reap it
when (case sz of Inf -> False; Nat sz' -> toInteger msz >= sz')
(liftIO (writeIORef eval
(\j -> panic "lookupSeqMap" ["Messed up size accounting", show sz, show j])))
return v
instance Backend sym => Functor (SeqMap sym) where
fmap f xs = IndexSeqMap (\i -> f <$> lookupSeqMap xs i)
@ -89,36 +108,37 @@ finiteSeqMap :: Backend sym => sym -> [SEval sym a] -> SeqMap sym a
finiteSeqMap sym xs =
UpdateSeqMap
(Map.fromList (zip [0..] xs))
(\i -> invalidIndex sym i)
(IndexSeqMap (\i -> invalidIndex sym i))
-- | Generate an infinite sequence map from a stream of values
infiniteSeqMap :: Backend sym => sym -> [SEval sym a] -> SEval sym (SeqMap sym a)
infiniteSeqMap sym xs =
-- TODO: use an int-trie?
memoMap sym (IndexSeqMap $ \i -> genericIndex xs i)
memoMap sym Inf (IndexSeqMap $ \i -> genericIndex xs i)
-- | Create a finite list of length @n@ of the values from @[0..n-1]@ in
-- the given the sequence emap.
enumerateSeqMap :: (Integral n) => n -> SeqMap sym a -> [SEval sym a]
enumerateSeqMap :: (Backend sym, Integral n) => n -> SeqMap sym a -> [SEval sym a]
enumerateSeqMap n m = [ lookupSeqMap m i | i <- [0 .. (toInteger n)-1] ]
-- | Create an infinite stream of all the values in a sequence map
streamSeqMap :: SeqMap sym a -> [SEval sym a]
streamSeqMap :: Backend sym => SeqMap sym a -> [SEval sym a]
streamSeqMap m = [ lookupSeqMap m i | i <- [0..] ]
-- | Reverse the order of a finite sequence map
reverseSeqMap :: Integer -- ^ Size of the sequence map
-> SeqMap sym a
-> SeqMap sym a
reverseSeqMap :: Backend sym =>
Integer {- ^ Size of the sequence map -} ->
SeqMap sym a ->
SeqMap sym a
reverseSeqMap n vals = IndexSeqMap $ \i -> lookupSeqMap vals (n - 1 - i)
updateSeqMap :: SeqMap sym a -> Integer -> SEval sym a -> SeqMap sym a
updateSeqMap (UpdateSeqMap m sm) i x = UpdateSeqMap (Map.insert i x m) sm
updateSeqMap (IndexSeqMap f) i x = UpdateSeqMap (Map.singleton i x) f
updateSeqMap xs i x = UpdateSeqMap (Map.singleton i x) xs
-- | Concatenate the first @n@ values of the first sequence map onto the
-- beginning of the second sequence map.
concatSeqMap :: Integer -> SeqMap sym a -> SeqMap sym a -> SeqMap sym a
concatSeqMap :: Backend sym => Integer -> SeqMap sym a -> SeqMap sym a -> SeqMap sym a
concatSeqMap n x y =
IndexSeqMap $ \i ->
if i < n
@ -128,14 +148,14 @@ concatSeqMap n x y =
-- | Given a number @n@ and a sequence map, return two new sequence maps:
-- the first containing the values from @[0..n-1]@ and the next containing
-- the values from @n@ onward.
splitSeqMap :: Integer -> SeqMap sym a -> (SeqMap sym a, SeqMap sym a)
splitSeqMap :: Backend sym => Integer -> SeqMap sym a -> (SeqMap sym a, SeqMap sym a)
splitSeqMap n xs = (hd,tl)
where
hd = xs
tl = IndexSeqMap $ \i -> lookupSeqMap xs (i+n)
-- | Drop the first @n@ elements of the given 'SeqMap'.
dropSeqMap :: Integer -> SeqMap sym a -> SeqMap sym a
dropSeqMap :: Backend sym => Integer -> SeqMap sym a -> SeqMap sym a
dropSeqMap 0 xs = xs
dropSeqMap n xs = IndexSeqMap $ \i -> lookupSeqMap xs (i+n)
@ -146,23 +166,20 @@ delaySeqMap sym xs =
-- | Given a sequence map, return a new sequence map that is memoized using
-- a finite map memo table.
memoMap :: Backend sym => sym -> SeqMap sym a -> SEval sym (SeqMap sym a)
memoMap sym x = do
memoMap :: Backend sym => sym -> Nat' -> SeqMap sym a -> SEval sym (SeqMap sym a)
-- Sequence is alreay memoized, just return it
memoMap _sym _sz x@(MemoSeqMap{}) = pure x
memoMap sym sz x = do
stk <- sGetCallStack sym
cache <- liftIO $ newIORef $ Map.empty
return $ IndexSeqMap (memo cache stk)
evalRef <- liftIO $ newIORef $ eval stk
return (MemoSeqMap sz cache evalRef)
where
memo cache stk i = do
mz <- liftIO (Map.lookup i <$> readIORef cache)
case mz of
Just z -> return z
Nothing -> sWithCallStack sym stk (doEval cache i)
eval stk i = sWithCallStack sym stk (lookupSeqMap x i)
doEval cache i = do
v <- lookupSeqMap x i
liftIO $ atomicModifyIORef' cache (\m -> (Map.insert i v m, ()))
return v
-- | Apply the given evaluation function pointwise to the two given
-- sequence maps.
@ -170,20 +187,23 @@ zipSeqMap ::
Backend sym =>
sym ->
(a -> a -> SEval sym a) ->
Nat' ->
SeqMap sym a ->
SeqMap sym a ->
SEval sym (SeqMap sym a)
zipSeqMap sym f x y =
memoMap sym (IndexSeqMap $ \i -> join (f <$> lookupSeqMap x i <*> lookupSeqMap y i))
zipSeqMap sym f sz x y =
memoMap sym sz (IndexSeqMap $ \i -> join (f <$> lookupSeqMap x i <*> lookupSeqMap y i))
-- | Apply the given function to each value in the given sequence map
mapSeqMap ::
Backend sym =>
sym ->
(a -> SEval sym a) ->
SeqMap sym a -> SEval sym (SeqMap sym a)
mapSeqMap sym f x =
memoMap sym (IndexSeqMap $ \i -> f =<< lookupSeqMap x i)
Nat' ->
SeqMap sym a ->
SEval sym (SeqMap sym a)
mapSeqMap sym f sz x =
memoMap sym sz (IndexSeqMap $ \i -> f =<< lookupSeqMap x i)
{-# INLINE mergeSeqMap #-}
@ -213,7 +233,7 @@ shiftSeqByInteger sym merge reindex zro m xs idx
| Just j <- integerAsLit sym idx = shiftOp xs j
| otherwise =
do (n, idx_bits) <- enumerateIntBits sym m idx
barrelShifter sym merge shiftOp xs n (map BitIndexSegment idx_bits)
barrelShifter sym merge shiftOp m xs n (map BitIndexSegment idx_bits)
where
shiftOp vs shft =
pure $ indexSeqMap $ \i ->
@ -231,6 +251,7 @@ data IndexSegment sym
Concrete ->
(SBit Concrete -> a -> a -> SEval Concrete a) ->
(SeqMap Concrete a -> Integer -> SEval Concrete (SeqMap Concrete a)) ->
Nat' ->
SeqMap Concrete a ->
Integer ->
[IndexSegment Concrete] ->
@ -241,11 +262,12 @@ barrelShifter :: Backend sym =>
(SBit sym -> a -> a -> SEval sym a) ->
(SeqMap sym a -> Integer -> SEval sym (SeqMap sym a))
{- ^ concrete shifting operation -} ->
Nat' {- ^ Size of the map being shifted -} ->
SeqMap sym a {- ^ initial value -} ->
Integer {- Number of bits in shift amount -} ->
[IndexSegment sym] {- ^ segments of the shift amount, in big-endian order -} ->
SEval sym (SeqMap sym a)
barrelShifter sym mux shift_op x0 n0 bs0
barrelShifter sym mux shift_op sz x0 n0 bs0
| n0 >= toInteger (maxBound :: Int) =
liftIO (X.throw (UnsupportedSymbolicOp ("Barrel shifter with too many bits in shift amount: " ++ show n0)))
| otherwise = go x0 (fromInteger n0) bs0
@ -273,5 +295,5 @@ barrelShifter sym mux shift_op x0 n0 bs0
go x_shft n' bs
Nothing ->
do x_shft <- shift_op x (bit n')
x' <- memoMap sym (mergeSeqMap sym mux b x_shft x)
x' <- memoMap sym sz (mergeSeqMap sym mux b x_shft x)
go x' n' bs

View File

@ -66,7 +66,7 @@ import Cryptol.Backend.Concrete (Concrete(..))
import Cryptol.Backend.Monad (EvalError(..))
import Cryptol.Backend.SeqMap
import Cryptol.TypeCheck.Solver.InfNat(widthInteger)
import Cryptol.TypeCheck.Solver.InfNat(widthInteger, Nat'(..))
-- | Force the evaluation of a word value
forceWordValue :: Backend sym => WordValue sym -> SEval sym ()
@ -258,19 +258,19 @@ wordValLogicOp _sym _ wop (WordVal w1) (WordVal w2) = WordVal <$> wop w1 w2
wordValLogicOp sym bop wop (WordVal w1) (BitmapVal n2 packed2 bs2) =
isReady sym packed2 >>= \case
Just w2 -> WordVal <$> wop w1 w2
Nothing -> bitmapWordVal sym n2 =<< zipSeqMap sym bop (unpackBitmap sym w1) bs2
Nothing -> bitmapWordVal sym n2 =<< zipSeqMap sym bop (Nat n2) (unpackBitmap sym w1) bs2
wordValLogicOp sym bop wop (BitmapVal n1 packed1 bs1) (WordVal w2) =
isReady sym packed1 >>= \case
Just w1 -> WordVal <$> wop w1 w2
Nothing -> bitmapWordVal sym n1 =<< zipSeqMap sym bop bs1 (unpackBitmap sym w2)
Nothing -> bitmapWordVal sym n1 =<< zipSeqMap sym bop (Nat n1) bs1 (unpackBitmap sym w2)
wordValLogicOp sym bop wop (BitmapVal n1 packed1 bs1) (BitmapVal _n2 packed2 bs2) =
do r1 <- isReady sym packed1
r2 <- isReady sym packed2
case (r1,r2) of
(Just w1, Just w2) -> WordVal <$> wop w1 w2
_ -> bitmapWordVal sym n1 =<< zipSeqMap sym bop bs1 bs2
_ -> bitmapWordVal sym n1 =<< zipSeqMap sym bop (Nat n1) bs1 bs2
wordValLogicOp sym bop wop (ThunkWordVal _ m1) w2 =
do w1 <- m1
@ -293,7 +293,7 @@ wordValUnaryOp sym bop wop (ThunkWordVal _ m) = wordValUnaryOp sym bop wop =<< m
wordValUnaryOp sym bop wop (BitmapVal n packed xs) =
isReady sym packed >>= \case
Just w -> WordVal <$> wop w
Nothing -> bitmapWordVal sym n =<< mapSeqMap sym bop xs
Nothing -> bitmapWordVal sym n =<< mapSeqMap sym bop (Nat n) xs
{-# SPECIALIZE joinWords ::
Concrete ->
@ -537,7 +537,7 @@ shiftWordByInteger sym wop reindex x idx =
Nothing ->
do (numbits, idx_bits) <- enumerateIntBits' sym n idx
bitmapWordVal sym n =<<
barrelShifter sym (iteBit sym) shiftOp bs0 numbits (map BitIndexSegment idx_bits)
barrelShifter sym (iteBit sym) shiftOp (Nat n) bs0 numbits (map BitIndexSegment idx_bits)
where
shiftOp vs shft =
@ -577,7 +577,7 @@ shiftWordByWord sym wop reindex x idx =
bitmapWordVal sym n =<< shiftOp bs0 j
_ ->
do idx_segs <- enumerateIndexSegments sym idx
bitmapWordVal sym n =<< barrelShifter sym (iteBit sym) shiftOp bs0 n idx_segs
bitmapWordVal sym n =<< barrelShifter sym (iteBit sym) shiftOp (Nat n) bs0 n idx_segs
where
shiftOp vs shft =
@ -654,12 +654,13 @@ shiftSeqByWord ::
(SBit sym -> a -> a -> SEval sym a) ->
(Integer -> Integer -> Maybe Integer) ->
SEval sym a ->
Nat' ->
SeqMap sym a ->
WordValue sym ->
SEval sym (SeqMap sym a)
shiftSeqByWord sym merge reindex zro xs idx =
shiftSeqByWord sym merge reindex zro sz xs idx =
do idx_segs <- enumerateIndexSegments sym idx
barrelShifter sym merge shiftOp xs (wordValueSize sym idx) idx_segs
barrelShifter sym merge shiftOp sz xs (wordValueSize sym idx) idx_segs
where
shiftOp vs shft =
pure $ indexSeqMap $ \i ->
@ -754,7 +755,7 @@ mergeBitmaps :: Backend sym =>
SeqMap sym (SBit sym) ->
SEval sym (WordValue sym)
mergeBitmaps sym c sz bs1 bs2 =
do bs <- memoMap sym (mergeSeqMap sym (iteBit sym) c bs1 bs2)
do bs <- memoMap sym (Nat sz) (mergeSeqMap sym (iteBit sym) c bs1 bs2)
bitmapWordVal sym sz bs
{-# INLINE mergeWord' #-}

View File

@ -603,7 +603,7 @@ evalComp ::
SEval sym (GenValue sym)
evalComp sym env len elty body ms =
do lenv <- mconcat <$> mapM (branchEnvs sym (toListEnv env)) ms
mkSeq sym len elty =<< memoMap sym (indexSeqMap $ \i -> do
mkSeq sym len elty =<< memoMap sym len (indexSeqMap $ \i -> do
evalExpr sym (evalListEnv lenv i) body)
{-# SPECIALIZE branchEnvs ::
@ -621,24 +621,24 @@ branchEnvs ::
ListEnv sym ->
[Match] ->
SEval sym (ListEnv sym)
branchEnvs sym env matches = foldM (evalMatch sym) env matches
branchEnvs sym env matches = snd <$> foldM (evalMatch sym) (1, env) matches
{-# SPECIALIZE evalMatch ::
(?range :: Range, ConcPrims) =>
Concrete ->
ListEnv Concrete ->
(Integer, ListEnv Concrete) ->
Match ->
SEval Concrete (ListEnv Concrete)
SEval Concrete (Integer, ListEnv Concrete)
#-}
-- | Turn a match into the list of environments it represents.
evalMatch ::
(?range :: Range, EvalPrims sym) =>
sym ->
ListEnv sym ->
(Integer, ListEnv sym) ->
Match ->
SEval sym (ListEnv sym)
evalMatch sym lenv m = case m of
SEval sym (Integer, ListEnv sym)
evalMatch sym (lsz, lenv) m = seq lsz $ case m of
-- many envs
From n l _ty expr ->
@ -646,7 +646,7 @@ evalMatch sym lenv m = case m of
-- Select from a sequence of finite length. This causes us to 'stutter'
-- through our previous choices `nLen` times.
Nat nLen -> do
vss <- memoMap sym $ indexSeqMap $ \i -> evalExpr sym (evalListEnv lenv i) expr
vss <- memoMap sym (Nat lsz) $ indexSeqMap $ \i -> evalExpr sym (evalListEnv lenv i) expr
let stutter xs = \i -> xs (i `div` nLen)
let lenv' = lenv { leVars = fmap stutter (leVars lenv) }
let vs i = do let (q, r) = i `divMod` nLen
@ -655,7 +655,7 @@ evalMatch sym lenv m = case m of
VSeq _ xs' -> lookupSeqMap xs' r
VStream xs' -> lookupSeqMap xs' r
_ -> evalPanic "evalMatch" ["Not a list value"]
return $ bindVarList n vs lenv'
return (lsz * nLen, bindVarList n vs lenv')
-- Select from a sequence of infinite length. Note that this means we
-- will never need to backtrack into previous branches. Thus, we can convert
@ -673,7 +673,9 @@ evalMatch sym lenv m = case m of
VSeq _ xs' -> lookupSeqMap xs' i
VStream xs' -> lookupSeqMap xs' i
_ -> evalPanic "evalMatch" ["Not a list value"]
return $ bindVarList n vs lenv'
-- Selecting from an infinite list effectively resets the length of the
-- list environment, so return 1 as the length
return (1, bindVarList n vs lenv')
where
len = evalNumType (leTypes lenv) l
@ -681,7 +683,7 @@ evalMatch sym lenv m = case m of
-- XXX we don't currently evaluate these as though they could be recursive, as
-- they are typechecked that way; the read environment to evalExpr is the same
-- as the environment to bind a new name in.
Let d -> return $ bindVarList (dName d) (\i -> f (evalListEnv lenv i)) lenv
Let d -> return (lsz, bindVarList (dName d) (\i -> f (evalListEnv lenv i)) lenv)
where
f env =
case dDefinition d of

View File

@ -221,13 +221,13 @@ ringBinary sym opw opi opz opq opfp = loop
rw <- fromVWord sym "ringRight" r
stk <- sGetCallStack sym
VWord w . wordVal <$> (sWithCallStack sym stk (opw w lw rw))
| otherwise -> VSeq w <$> (join (zipSeqMap sym (loop a) <$>
| otherwise -> VSeq w <$> (join (zipSeqMap sym (loop a) (Nat w) <$>
(fromSeq "ringBinary left" l) <*>
(fromSeq "ringBinary right" r)))
TVStream a ->
-- streams
VStream <$> (join (zipSeqMap sym (loop a) <$>
VStream <$> (join (zipSeqMap sym (loop a) Inf <$>
(fromSeq "ringBinary left" l) <*>
(fromSeq "ringBinary right" r)))
@ -307,10 +307,10 @@ ringUnary sym opw opi opz opq opfp = loop
wx <- fromVWord sym "ringUnary" v
stk <- sGetCallStack sym
VWord w . wordVal <$> sWithCallStack sym stk (opw w wx)
| otherwise -> VSeq w <$> (mapSeqMap sym (loop a) =<< fromSeq "ringUnary" v)
| otherwise -> VSeq w <$> (mapSeqMap sym (loop a) (Nat w) =<< fromSeq "ringUnary" v)
TVStream a ->
VStream <$> (mapSeqMap sym (loop a) =<< fromSeq "ringUnary" v)
VStream <$> (mapSeqMap sym (loop a) Inf =<< fromSeq "ringUnary" v)
-- functions
TVFun _ ety ->
@ -1260,12 +1260,12 @@ logicBinary sym opb opw = loop
-- finite sequences
| otherwise -> VSeq w <$>
(join (zipSeqMap sym (loop aty) <$>
(join (zipSeqMap sym (loop aty) (Nat w) <$>
(fromSeq "logicBinary left" l)
<*> (fromSeq "logicBinary right" r)))
TVStream aty ->
VStream <$> (join (zipSeqMap sym (loop aty) <$>
VStream <$> (join (zipSeqMap sym (loop aty) Inf <$>
(fromSeq "logicBinary left" l) <*>
(fromSeq "logicBinary right" r)))
@ -1324,11 +1324,11 @@ logicUnary sym opb opw = loop
-- finite sequences
| otherwise
-> VSeq w <$> (mapSeqMap sym (loop ety) =<< fromSeq "logicUnary" val)
-> VSeq w <$> (mapSeqMap sym (loop ety) (Nat w) =<< fromSeq "logicUnary" val)
-- streams
TVStream ety ->
VStream <$> (mapSeqMap sym (loop ety) =<< fromSeq "logicUnary" val)
VStream <$> (mapSeqMap sym (loop ety) Inf =<< fromSeq "logicUnary" val)
TVTuple etys ->
do as <- mapM (sDelay sym) (fromVTuple val)
@ -1620,8 +1620,8 @@ wordShifter :: Backend sym =>
wordShifter sym nm wop reindex m a xs idx =
case xs of
VWord w x -> VWord w <$> shiftWordByWord sym wop (reindex m) x idx
VSeq w vs -> VSeq w <$> shiftSeqByWord sym (mergeValue sym) (reindex m) (zeroV sym a) vs idx
VStream vs -> VStream <$> shiftSeqByWord sym (mergeValue sym) (reindex m) (zeroV sym a) vs idx
VSeq w vs -> VSeq w <$> shiftSeqByWord sym (mergeValue sym) (reindex m) (zeroV sym a) (Nat w) vs idx
VStream vs -> VStream <$> shiftSeqByWord sym (mergeValue sym) (reindex m) (zeroV sym a) Inf vs idx
_ -> evalPanic "expected sequence value in shift operation" [nm]

View File

@ -113,8 +113,8 @@ indexFront_segs ::
indexFront_segs sym mblen a xs ix _idx_bits [WordIndexSegment w] =
indexFront sym mblen a xs ix w
indexFront_segs sym _mblen _a xs _ix idx_bits segs =
do xs' <- barrelShifter sym (mergeValue sym) shiftOp xs idx_bits segs
indexFront_segs sym mblen _a xs _ix idx_bits segs =
do xs' <- barrelShifter sym (mergeValue sym) shiftOp mblen xs idx_bits segs
lookupSeqMap xs' 0
where
shiftOp vs amt = pure (indexSeqMap (\i -> lookupSeqMap vs $! amt+i))

View File

@ -492,8 +492,8 @@ mergeValue sym c v1 v2 =
(VRational q1, VRational q2) -> VRational <$> iteRational sym c q1 q2
(VFloat f1 , VFloat f2) -> VFloat <$> iteFloat sym c f1 f2
(VWord n1 w1 , VWord n2 w2 ) | n1 == n2 -> VWord n1 <$> mergeWord sym c w1 w2
(VSeq n1 vs1 , VSeq n2 vs2 ) | n1 == n2 -> VSeq n1 <$> memoMap sym (mergeSeqMapVal sym c vs1 vs2)
(VStream vs1 , VStream vs2 ) -> VStream <$> memoMap sym (mergeSeqMapVal sym c vs1 vs2)
(VSeq n1 vs1 , VSeq n2 vs2 ) | n1 == n2 -> VSeq n1 <$> memoMap sym (Nat n1) (mergeSeqMapVal sym c vs1 vs2)
(VStream vs1 , VStream vs2 ) -> VStream <$> memoMap sym Inf (mergeSeqMapVal sym c vs1 vs2)
(f1@VFun{} , f2@VFun{} ) -> lam sym $ \x -> mergeValue' sym c (fromVFun sym f1 x) (fromVFun sym f2 x)
(f1@VPoly{} , f2@VPoly{} ) -> tlam sym $ \x -> mergeValue' sym c (fromVPoly sym f1 x) (fromVPoly sym f2 x)
(_ , _ ) -> panic "Cryptol.Eval.Value"

View File

@ -586,8 +586,8 @@ indexFront_segs sym mblen _a xs _ix _idx_bits [WordIndexSegment idx]
Just (lo, hi) -> [lo .. min hi maxIdx]
_ -> [0 .. maxIdx]
indexFront_segs sym _mblen _a xs _ix idx_bits segs =
do xs' <- barrelShifter sym (mergeValue sym) shiftOp xs idx_bits segs
indexFront_segs sym mblen _a xs _ix idx_bits segs =
do xs' <- barrelShifter sym (mergeValue sym) shiftOp mblen xs idx_bits segs
lookupSeqMap xs' 0
where
shiftOp vs amt = pure (indexSeqMap (\i -> lookupSeqMap vs $! amt+i))
@ -609,11 +609,11 @@ updateFrontSym sym _len _eltTy vs (Left idx) val =
do b <- intEq sym idx =<< integerLit sym i
iteValue sym b val (lookupSeqMap vs i)
updateFrontSym sym _len _eltTy vs (Right wv) val =
updateFrontSym sym len _eltTy vs (Right wv) val =
wordValAsLit sym wv >>= \case
Just j -> return $ updateSeqMap vs j val
Nothing ->
memoMap sym $ indexSeqMap $ \i ->
memoMap sym len $ indexSeqMap $ \i ->
do b <- wordValueEqualsInteger sym wv i
iteValue sym b val (lookupSeqMap vs i)
@ -640,7 +640,7 @@ updateBackSym sym (Nat n) _eltTy vs (Right wv) val =
Just j ->
return $ updateSeqMap vs (n - 1 - j) val
Nothing ->
memoMap sym $ indexSeqMap $ \i ->
memoMap sym (Nat n) $ indexSeqMap $ \i ->
do b <- wordValueEqualsInteger sym wv (n - 1 - i)
iteValue sym b val (lookupSeqMap vs i)