System representation

IntervalMDP.initial_statesMethod
initial_states(mp::IntervalMarkovProcess)

Return the initial states. If the initial states are not specified, return nothing.

source
IntervalMDP.AllStatesType
AllStates

A type to represent all states in a Markov process. This type is used to specify all states as the initial states.

source
IntervalMDP.IntervalMarkovDecisionProcessType
IntervalMarkovDecisionProcess{
    P <: IntervalProbabilities,
    VT <: AbstractVector{Int32},
    VI <: Union{AllStates, AbstractVector}
}

A type representing (stationary) Interval Markov Decision Processes (IMDP), which are Markov Decision Processes with uncertainty in the form of intervals on the transition probabilities.

Formally, let $(S, S_0, A, \Gamma)$ be an interval Markov decision process, where

  • $S$ is the set of states,
  • $S_0 \subseteq S$ is the set of initial states,
  • $A$ is the set of actions, and
  • $\Gamma = \{\Gamma_{s,a}\}_{s \in S, a \in A}$ is a set of interval ambiguity sets on the transition probabilities, for each source-action pair.

Then the IntervalMarkovDecisionProcess type is defined as follows: indices 1:num_states are the states in $S$, transition_prob represents $\Gamma$, actions are implicitly defined by stateptr (e.g. if stateptr[3] == 4 and stateptr[4] == 7 then the actions available to state 3 are [1, 2, 3]), and initial_states is the set of initial states $S_0$. If no initial states are specified, then the initial states are assumed to be all states in $S$ represented by AllStates. See IntervalProbabilities and Theory for more information on the structure of the transition probability ambiguity sets.

Fields

  • transition_prob::P: interval on transition probabilities where columns represent source/action pairs and rows represent target states.
  • stateptr::VT: pointer to the start of each source state in transition_prob (i.e. transition_prob[:, stateptr[j]:stateptr[j + 1] - 1] is the transition probability matrix for source state j) in the style of colptr for sparse matrices in CSC format.
  • initial_states::VI: initial states.
  • num_states::Int32: number of states.

Examples

transition_probs = IntervalProbabilities(;
    lower = [
        0.0 0.5 0.1 0.2 0.0
        0.1 0.3 0.2 0.3 0.0
        0.2 0.1 0.3 0.4 1.0
    ],
    upper = [
        0.5 0.7 0.6 0.6 0.0
        0.6 0.5 0.5 0.5 0.0
        0.7 0.3 0.4 0.4 1.0
    ],
)

stateptr = [1, 3, 5, 6]
initial_states = [1]

mdp = IntervalMarkovDecisionProcess(transition_probs, stateptr, initial_states)

There is also a constructor for IntervalMarkovDecisionProcess where the transition probabilities are given as a list of transition probabilities for each source state.

prob1 = IntervalProbabilities(;
    lower = [
        0.0 0.5
        0.1 0.3
        0.2 0.1
    ],
    upper = [
        0.5 0.7
        0.6 0.5
        0.7 0.3
    ],
)

prob2 = IntervalProbabilities(;
    lower = [
        0.1 0.2
        0.2 0.3
        0.3 0.4
    ],
    upper = [
        0.6 0.6
        0.5 0.5
        0.4 0.4
    ],
)

prob3 = IntervalProbabilities(;
    lower = [0.0; 0.0; 1.0],
    upper = [0.0; 0.0; 1.0]
)

transition_probs = [prob1, prob2, prob3]
initial_states = [1]

mdp = IntervalMarkovDecisionProcess(transition_probs, initial_states)
source
IntervalMDP.IntervalMarkovChainFunction
IntervalMarkovChain(transition_prob::IntervalProbabilities, initial_states::InitialStates = AllStates())

Construct an Interval Markov Chain from a square matrix pair of interval transition probabilities. The initial states are optional and if not specified, all states are assumed to be initial states. The number of states is inferred from the size of the transition probability matrix.

The returned type is an IntervalMarkovDecisionProcess with only one action per state (i.e. stateptr[j + 1] - stateptr[j] == 1 for all j). This is done to unify the interface for value iteration.

source
IntervalMDP.stateptrMethod
stateptr(mdp::IntervalMarkovDecisionProcess)

Return the state pointer of the Interval Markov Decision Process. The state pointer is a vector of integers where the i-th element is the index of the first element of the i-th state in the transition probability matrix. I.e. transition_prob[:, stateptr[j]:stateptr[j + 1] - 1] is the transition probability matrix for source state j.

source
IntervalMDP.OrthogonalIntervalMarkovDecisionProcessType
OrthogonalIntervalMarkovDecisionProcess{
    P <: OrthogonalIntervalProbabilities,
    VT <: AbstractVector{Int32},
    VI <: Union{AllStates, AbstractVector}
}

A type representing (stationary) Orthogonal Interval Markov Decision Processes (OIMDP), which are IMDPs where the transition probabilities for each state can be represented as the product of the transition probabilities of individual processes.

Formally, let $(S, S_0, A, \Gamma)$ be an orthogonal interval Markov decision process [1], where

  • $S = S_1 \times \cdots \times S_n$ is the set of joint states with $S_i$ the set of states for the i-th marginal,
  • $S_0 \subseteq S$ is the set of initial states,
  • $A$ is the set of actions, and
  • $\Gamma = \{\Gamma_{s,a}\}_{s \in S, a \in A}$ is a set of interval ambiguity sets on the transition probabilities, for each source-action pair, with $\Gamma_{s,a} = \bigotimes_{i=1}^n \Gamma_{s,a}^i$ and $\Gamma_{s,a}^i$ is a marginal interval ambiguity sets on the $i$-th marginal.

Then the OrthogonalIntervalMarkovDecisionProcess type is defined as follows: indices 1:num_states are the states in $S$ and transition_prob represents $\Gamma$. Actions are implicitly defined by stateptr (e.g. if source_dims in transition_prob is (2, 3, 2), and stateptr[3] == 4 and stateptr[4] == 7 then the actions available to state CartesianIndex(1, 2, 1) are [1, 2, 3]), and initial_states is the set of initial states $S_0$. If no initial states are specified, then the initial states are assumed to be all states in $S$ represented by AllStates. See OrthogonalIntervalProbabilities and Theory for more information on the structure of the transition probability ambiguity sets.

Fields

  • transition_prob::P: interval on transition probabilities where columns represent source/action pairs and rows represent target states along each marginal.
  • stateptr::VT: pointer to the start of each source state in transition_prob (i.e. transition_prob[l][:, stateptr[j]:stateptr[j + 1] - 1] is the transition probability matrix for source state j for each axis l) in the style of colptr for sparse matrices in CSC format.
  • initial_states::VI: initial states.
  • num_states::Int32: number of states.

Examples

Assume that prob1, prob2, and prob3 are IntervalProbabilities for the first, second, and third axis, respectively, defined as the example in OrthogonalIntervalProbabilities. Then the following code constructs an OrthogonalIntervalMarkovDecisionProcess with three axes of three states each. The number of actions per state is one, i.e. the model is a Markov chain. Therefore, the stateptr is a unit range 1:num_states + 1 and we can call the convenience constructor OrthogonalIntervalMarkovChain instead.

prob = OrthogonalIntervalProbabilities((prob1, prob2, prob3), (Int32(3), Int32(3), Int32(3)))
mc = OrthogonalIntervalMarkovChain(prob)

[1] Mathiesen, F. B., Haesaert, S., & Laurenti, L. (2024). Scalable control synthesis for stochastic systems via structural IMDP abstractions. arXiv preprint arXiv:2411.11803.

source
IntervalMDP.OrthogonalIntervalMarkovChainFunction
OrthogonalIntervalMarkovChain(transition_prob::OrthogonalIntervalProbabilities, initial_states::InitialStates = AllStates())

Construct a Orthogonal Interval Markov Chain from orthogonal interval transition probabilities. The initial states are optional and if not specified, all states are assumed to be initial states. The number of states is inferred from the size of the transition probability matrix.

The returned type is an OrthogonalIntervalMarkovDecisionProcess with only one action per state (i.e. stateptr[j + 1] - stateptr[j] == 1 for all j). This is done to unify the interface for value iteration.

source
IntervalMDP.stateptrMethod
stateptr(mdp::OrthogonalIntervalMarkovDecisionProcess)

Return the state pointer of the Interval Markov Decision Process. The state pointer is a vector of integers where the i-th element is the index of the first element of the i-th state in the transition probability matrix. I.e. transition_prob[l][:, stateptr[j]:stateptr[j + 1] - 1] is the transition probability matrix for (flattened) source state j for axis l.

source
IntervalMDP.MixtureIntervalMarkovDecisionProcessType
MixtureIntervalMarkovDecisionProcess{
    P <: IntervalProbabilities,
    VT <: AbstractVector{Int32},
    VI <: Union{AllStates, AbstractVector}
}

A type representing (stationary) Mixture Interval Markov Decision Processes (OIMDP), which are IMDPs where the transition probabilities for each state can be represented as the product of the transition probabilities of individual processes.

Formally, let $(S, S_0, A, \Gamma, \Gamma_\alpha)$ be an interval Markov decision processes, where

  • $S = S_1 \times \cdots \times S_n$ is the set of joint states with $S_i$ the set of states for the i-th marginal,
  • $S_0 \subseteq S$ is the set of initial states,
  • $A$ is the set of actions,
  • $\Gamma = \{\Gamma_{s,a}\}_{s \in S, a \in A}$ is a set of interval ambiguity sets on the transition probabilities, for each source-action pair, with $\Gamma_{s,a} = \bigotimes_{i=1}^n \Gamma_{s,a}^i$ and $\Gamma_{s,a}^i$ is a marginal interval ambiguity sets on the $i$-th marginal, and
  • $\Gamma^\alpha = \{\Gamma^\alpha_{s,a}\}_{s \in S, a \in A}$ is the interval ambiguity set for the mixture.

Then the MixtureIntervalMarkovDecisionProcess type is defined as follows: indices 1:num_states are the states in $S$ and transition_prob represents $\Gamma$ and $\Gamma^\alpha$. Actions are implicitly defined by stateptr (e.g. if source_dims in transition_prob is (2, 3, 2), and stateptr[3] == 4 and stateptr[4] == 7 then the actions available to state CartesianIndex(1, 2, 1) are [1, 2, 3]), and initial_states is the set of initial states $S_0$. If no initial states are specified, then the initial states are assumed to be all states in $S$ represented by AllStates. See MixtureIntervalProbabilities and Theory for more information on the structure of the transition probability ambiguity sets.

Fields

  • transition_prob::P: ambiguity set on transition probabilities (see MixtureIntervalProbabilities for the structure).
  • stateptr::VT: pointer to the start of each source state in transition_prob (i.e. transition_prob[k][l][:, stateptr[j]:stateptr[j + 1] - 1] is the transition probability matrix for source state j for each model k and axis l) in the style of colptr for sparse matrices in CSC format.
  • initial_states::VI: initial states.
  • num_states::Int32: number of states.

Examples

The following example is a simple mixture of two OrthogonalIntervalProbabilities with one dimension and the same source/action pairs. The first state has two actions and the second state has one action. The weighting ambiguity set is also specified for the same three source-action pairs.

prob1 = OrthogonalIntervalProbabilities(
    (
        IntervalProbabilities(;
            lower = [
                0.0 0.5 0.1
                0.1 0.3 0.2
            ],
            upper = [
                0.5 0.7 0.6
                0.7 0.4 0.8
            ],
        ),
    ),
    (Int32(2),),
)
prob2 = OrthogonalIntervalProbabilities(
    (
        IntervalProbabilities(;
            lower = [
                0.1 0.4 0.2
                0.3 0.0 0.1
            ],
            upper = [
                0.4 0.6 0.5
                0.7 0.5 0.7
            ],
        ),
    ),
    (Int32(2),),
)
weighting_probs = IntervalProbabilities(; lower = [
    0.3 0.5 0.4
    0.4 0.3 0.2
], upper = [
    0.8 0.7 0.7
    0.7 0.5 0.4
])
mixture_prob = MixtureIntervalProbabilities((prob1, prob2), weighting_probs)

stateptr = Int32[1, 3, 4]
mdp = MixtureIntervalMarkovDecisionProcess(mixture_prob, stateptr)
source
IntervalMDP.MixtureIntervalMarkovChainFunction
MixtureIntervalMarkovChain(transition_prob::MixtureIntervalProbabilities, initial_states::InitialStates = AllStates())

Construct a Mixture Interval Markov Chain from mixture interval transition probabilities. The initial states are optional and if not specified, all states are assumed to be initial states. The number of states is inferred from the size of the transition probability matrix.

The returned type is an MixtureIntervalMarkovDecisionProcess with only one action per state (i.e. stateptr[j + 1] - stateptr[j] == 1 for all j). This is done to unify the interface for value iteration.

source
IntervalMDP.stateptrMethod
stateptr(mdp::MixtureIntervalMarkovDecisionProcess)

Return the state pointer of the Interval Markov Decision Process. The state pointer is a vector of integers where the i-th element is the index of the first element of the i-th state in the transition probability matrix. I.e. mixture_probs(transition_prob)[k][l][:, stateptr[j]:stateptr[j + 1] - 1] is the independent transition probability matrix for (flattened) source state j for axis l and model k, and mixture_probs(transition_prob)[:, stateptr[j]:stateptr[j + 1] - 1] is the weighting matrix for j.

source
IntervalMDP.DFAType
DFA{
    T <: TransitionFunction,  
    VT <: AbstractVector{Int32},
    DA <: AbstractDict{String, Int32}
}

A type representing Deterministic Finite Automaton (DFA) which are finite automata with deterministic transitions.

Formally, let $(Z, 2^{AP}, \tau, z_0, Z_{ac})$ be an DFA, where

  • $Z$ is the set of states,
  • $Z_{ac} \subseteq Z$ is the set of accepting states,
  • $z_0$ is the initial state,
  • $2^{AP}$ is the power set of automic propositions, and
  • $\tau : |Z| \times |2^{AP}| => |Z|$ is the deterministic transition function, for each state-input pair.

Then the DFA type is defined as follows: indices 1:num_states are the states in $Z$, transition represents $\tau$, inputs are implicitly defined by indices enumerating the alphabet, and initial_states is the set of initial states $z_0$. See TransitionFunction for more information on the structure of the transition function.

Fields

  • transition::T: transition function where columns represent inputs and rows source states.
  • initial_state::Int32: initial states.
  • accepting_states::VT: vector of accepting states
  • alphabetptr::DA: mapping from input to index.
source
IntervalMDP.alphabet2indexMethod
Given vector of alphabet 2^L, maps its index for lookup
Returns dictionary, assume same indices used in Transition function.
source
IntervalMDP.alphabetptrMethod
alphabetptr(dfa::DFA)

Return the alphabet ($2^{AP}$) index mapping of the Deterministic Finite Automaton.

source
IntervalMDP.ProductIntervalMarkovDecisionProcessDFAType
struct ProductIntervalMarkovDecisionProcessDFA{
    D <: DFA,
    M <: IntervalMarkovDecisionProcess,
    L <: AbstractLabelling,
}

A type representing the product between Interval Markov Decision Processes (IMDP) and Deterministic Finite Automata (DFA).

Formally, let $(S^{\prime}, A, \Gamma^{\prime}, S^{\prime}_{ac}, S^{\prime}_{0}, L)$ be an product interval Markov decision process DFA, where

  • $S^{\prime} = S \times Z$ is the set of product states q = (s, z)``,
  • $S^{\prime}_{0} = S_0 \times z_0 \subset S^{\prime}$ is the set of initial product states $q = (s, z_0)$,
  • $S^{\prime}_{ac} = S \times Z_{ac} \subseteq S^{\prime}$ is the set of accepting product states,
  • $A$ is the set of actions,
  • $\Gamma^{\prime} = \{\Gamma^{\prime}_{q,a}\}_{q \in S^{\prime}, a \in A}$ is a set of interval ambiguity sets on the transition probabilities, for each product source-action pair, and
  • $L : S => 2^{AP}$ is the labelling function that maps a state in the IMDP to an input to the DFA.

Then the ProductIntervalMarkovDecisionProcessDFA type is defined as follows: DFA part of the product as a DFA object, IMDP part of the product as a IMDP object and a Labelling function to specify the relationshup between the DFA and IMDP.

See IntervalMarkovDecisionProcess and DFA for more information on the structure, definition, and usage of the DFA and IMDP.

Fields

  • imdp::M: contains details for the IMDP
  • dfa::D: contains details for the DFA
  • labelling_func::L: the labelling function from IMDP states to DFA actions
source
IntervalMDP.imdpMethod
imdp(md::ProductIntervalMarkovDecisionProcessDFA)

Return the interval markov decision process of the product

source
IntervalMDP.automatonMethod
automaton(md::ProductIntervalMarkovDecisionProcessDFA)

Return the deterministic finite automaton of the product

source

Probability representation

Interval ambiguity sets

IntervalMDP.IntervalProbabilitiesType
IntervalProbabilities{R, VR <: AbstractVector{R}, MR <: AbstractMatrix{R}}

A matrix pair to represent the lower and upper bound transition probabilities from all source states or source/action pairs to all target states. The matrices can be Matrix{R} or SparseMatrixCSC{R}, or their CUDA equivalents. For memory efficiency, it is recommended to use sparse matrices.

The columns represent the source and the rows represent the target, as if the probability matrix was a linear transformation. Mathematically, let $P$ be the probability matrix. Then $P_{ij}$ represents the probability of transitioning from state $j$ (or with state/action pair $j$) to state $i$. Due to the column-major format of Julia, this is also a more efficient representation (in terms of cache locality).

The lower bound is explicitly stored, while the upper bound is computed from the lower bound and the gap. This choice is because it simplifies repeated probability assignment using O-maximization [1].

Fields

  • lower::MR: The lower bound transition probabilities from a source state or source/action pair to a target state.
  • gap::MR: The gap between upper and lower bound transition probabilities from a source state or source/action pair to a target state.
  • sum_lower::VR: The sum of lower bound transition probabilities from a source state or source/action pair to all target states.

Examples

dense_prob = IntervalProbabilities(;
    lower = [0.0 0.5; 0.1 0.3; 0.2 0.1],
    upper = [0.5 0.7; 0.6 0.5; 0.7 0.3],
)

sparse_prob = IntervalProbabilities(;
    lower = sparse_hcat(
        SparseVector(15, [4, 10], [0.1, 0.2]),
        SparseVector(15, [5, 6, 7], [0.5, 0.3, 0.1]),
    ),
    upper = sparse_hcat(
        SparseVector(15, [1, 4, 10], [0.5, 0.6, 0.7]),
        SparseVector(15, [5, 6, 7], [0.7, 0.5, 0.3]),
    ),
)

[1] M. Lahijanian, S. B. Andersson and C. Belta, "Formal Verification and Synthesis for Discrete-Time Stochastic Systems," in IEEE Transactions on Automatic Control, vol. 60, no. 8, pp. 2031-2045, Aug. 2015, doi: 10.1109/TAC.2015.2398883.

source
IntervalMDP.lowerMethod
lower(p::IntervalProbabilities)

Return the lower bound transition probabilities from a source state or source/action pair to a target state.

source
IntervalMDP.lowerMethod
lower(p::IntervalProbabilities, i, j)

Return the lower bound transition probabilities from a source state or source/action pair to a target state.

source
IntervalMDP.upperMethod
upper(p::IntervalProbabilities)

Return the upper bound transition probabilities from a source state or source/action pair to a target state.

Note

It is not recommended to use this function for the hot loop of O-maximization. Because the IntervalProbabilities stores the lower and gap transition probabilities, fetching the upper bound requires allocation and computation.

source
IntervalMDP.upperMethod
upper(p::IntervalProbabilities, i, j)

Return the upper bound transition probabilities from a source state or source/action pair to a target state.

Note

It is not recommended to use this function for the hot loop of O-maximization. Because the IntervalProbabilities stores the lower and gap transition probabilities, fetching the upper bound requires allocation and computation.

source
IntervalMDP.gapMethod
gap(p::IntervalProbabilities)

Return the gap between upper and lower bound transition probabilities from a source state or source/action pair to a target state.

source
IntervalMDP.gapMethod
gap(p::IntervalProbabilities, i, j)

Return the gap between upper and lower bound transition probabilities from a source state or source/action pair to a target state.

source
IntervalMDP.sum_lowerMethod
sum_lower(p::IntervalProbabilities)

Return the sum of lower bound transition probabilities from a source state or source/action pair to all target states. This is useful in efficiently implementing O-maximization, where we start with a lower bound probability assignment and iteratively, according to the ordering, adding the gap until the sum of probabilities is 1.

source
IntervalMDP.sum_lowerMethod
sum_lower(p::IntervalProbabilities, j)

Return the sum of lower bound transition probabilities from a source state or source/action pair to all target states. This is useful in efficiently implementing O-maximization, where we start with a lower bound probability assignment and iteratively, according to the ordering, adding the gap until the sum of probabilities is 1.

source
IntervalMDP.axes_sourceMethod
axes_source(p::IntervalProbabilities)

Return the valid range of indices for the source states or source/action pairs.

source

Marginal interval ambiguity sets

IntervalMDP.OrthogonalIntervalProbabilitiesType
OrthogonalIntervalProbabilities{N, P <: IntervalProbabilities}

A tuple of IntervalProbabilities for (marginal) transition probabilities from all source/action pairs to the target states along each axis, with target states/marginals on the rows and source states or source/action pairs on the columns. The source states are ordered in a column-major order, i.e., the first axis of source states is the fastest, similar to the ordering of a multi-dimensional array in Julia. E.g. for an OrthogonalIntervalProbabilities with source_dims == (3, 3, 3) and 2 actions for each source state $\{a_1, a_2\}$, the columns in order represent the collowing:

\[ ((1, 1, 1), a_1), ((1, 1, 1), a_2), (2, 1, 1), a_1), ((2, 1, 1), a_2), ..., ((3, 3, 3), a_1), ((3, 3, 3), a_2).\]

The number of target states correspond to the number of rows in the transition probabilities of each axis.

Fields

  • probs::NTuple{N, P}: A tuple of IntervalProbabilities for (marginal) transition probabilities along each axis.
  • source_dims::NTuple{N, Int32}: The dimensions of the orthogonal probabilities for the source axis. This is flattened to a single dimension for indexing.

Examples

An example of OrthogonalIntervalProbabilities with 3 axes and 3 states for each axis, only one action per state. Therefore, the source_dims is (3, 3, 3) and the number of columns of the transition probabilities is 27.

lower1 = [
    1/15 3/10 1/15 3/10 1/30 1/3 7/30 4/15 1/6 1/5 1/10 1/5 0 7/30 7/30 1/5 2/15 1/6 1/10 1/30 1/10 1/15 1/10 1/15 4/15 4/15 1/3
    1/5 4/15 1/10 1/5 3/10 3/10 1/10 1/15 3/10 3/10 7/30 1/5 1/10 1/5 1/5 1/30 1/5 3/10 1/5 1/5 1/10 1/30 4/15 1/10 1/5 1/6 7/30
    4/15 1/30 1/5 1/5 7/30 4/15 2/15 7/30 1/5 1/3 2/15 1/6 1/6 1/3 4/15 3/10 1/30 3/10 3/10 1/10 1/15 1/30 2/15 1/6 1/5 1/10 4/15
]
upper1 = [
    7/15 17/30 13/30 3/5 17/30 17/30 17/30 13/30 3/5 2/3 11/30 7/15 0 1/2 17/30 13/30 7/15 13/30 17/30 13/30 2/5 2/5 2/3 2/5 17/30 2/5 19/30
    8/15 1/2 3/5 7/15 8/15 17/30 2/3 17/30 11/30 7/15 19/30 19/30 13/15 1/2 17/30 13/30 3/5 11/30 8/15 7/15 7/15 13/30 8/15 2/5 8/15 17/30 3/5
    11/30 1/3 2/5 8/15 7/15 3/5 2/3 17/30 2/3 8/15 2/15 3/5 2/3 3/5 17/30 2/3 7/15 8/15 2/5 2/5 11/30 17/30 17/30 1/2 2/5 19/30 13/30
]
prob1 = IntervalProbabilities(; lower = lower1, upper = upper1)

lower2 = [
    1/10 1/15 3/10 0 1/6 1/15 1/15 1/6 1/6 1/30 1/10 1/10 1/3 2/15 3/10 4/15 2/15 2/15 1/6 7/30 1/15 2/15 1/10 1/3 7/30 1/30 7/30
    3/10 1/5 3/10 2/15 0 1/30 0 1/15 1/30 7/30 1/30 1/15 7/30 1/15 1/6 1/30 1/10 1/15 3/10 0 3/10 1/6 3/10 1/5 0 7/30 2/15
    3/10 4/15 1/10 3/10 2/15 1/3 3/10 1/10 1/6 3/10 7/30 1/6 1/15 1/15 1/10 1/5 1/5 4/15 1/15 1/3 2/15 1/15 1/5 1/5 1/15 7/30 1/15
]
upper2 = [
    2/5 17/30 3/5 11/30 3/5 7/15 19/30 2/5 3/5 2/3 2/3 8/15 8/15 19/30 8/15 8/15 13/30 13/30 13/30 17/30 17/30 13/30 11/30 19/30 8/15 2/5 8/15
    1/3 13/30 11/30 2/5 2/3 2/3 0 13/30 1/2 17/30 17/30 1/3 2/5 1/3 13/30 11/30 8/15 1/3 1/2 8/15 8/15 8/15 8/15 2/5 3/5 2/3 13/30
    17/30 3/5 8/15 1/2 7/15 1/2 2/3 17/30 11/30 2/5 1/2 7/15 2/5 17/30 11/30 2/5 11/30 2/3 1/3 2/3 17/30 8/15 17/30 3/5 2/5 19/30 11/30
]
prob2 = IntervalProbabilities(; lower = lower2, upper = upper2)

lower3 = [
    4/15 1/5 3/10 3/10 4/15 7/30 1/5 4/15 7/30 1/6 1/5 0 1/15 1/30 3/10 1/3 2/15 1/15 7/30 4/15 1/10 1/3 1/5 7/30 1/30 1/5 7/30
    2/15 4/15 1/10 1/30 7/30 2/15 1/15 1/30 3/10 1/3 1/5 1/10 2/15 1/30 2/15 4/15 0 4/15 1/5 4/15 1/10 1/10 1/3 7/30 3/10 1/3 3/10
    1/5 1/3 3/10 1/10 1/15 1/10 1/30 1/5 2/15 7/30 1/3 2/15 1/10 1/6 3/10 1/5 7/30 1/30 0 1/30 1/15 2/15 1/6 7/30 4/15 4/15 7/30
]
upper3 = [
    3/5 17/30 1/2 3/5 19/30 2/5 8/15 1/3 11/30 2/5 17/30 13/30 2/5 3/5 3/5 11/30 1/2 11/30 2/3 17/30 3/5 7/15 19/30 1/2 3/5 1/3 19/30
    3/5 2/3 13/30 19/30 1/3 2/5 17/30 7/15 11/30 3/5 19/30 7/15 2/5 8/15 17/30 11/30 19/30 13/30 2/3 17/30 8/15 13/30 13/30 3/5 1/2 8/15 8/15
    3/5 2/3 1/2 1/2 2/3 7/15 3/5 3/5 1/2 1/3 2/5 8/15 2/5 11/30 1/3 8/15 7/15 13/30 0 2/5 11/30 19/30 19/30 2/5 1/2 7/15 7/15
]
prob3 = IntervalProbabilities(; lower = lower3, upper = upper3)

prob = OrthogonalIntervalProbabilities((prob1, prob2, prob3), (Int32(3), Int32(3), Int32(3)))
source
IntervalMDP.lowerMethod
lower(p::OrthogonalIntervalProbabilities, l)

Return the lower bound transition probabilities from a source state or source/action pair to a target axis.

source
IntervalMDP.lowerMethod
lower(p::OrthogonalIntervalProbabilities, l, i, j)

Return the lower bound transition probabilities from a source state or source/action pair to a target state.

source
IntervalMDP.upperMethod
upper(p::OrthogonalIntervalProbabilities, l)

Return the upper bound transition probabilities from a source state or source/action pair to a target state.

Note

It is not recommended to use this function for the hot loop of O-maximization. Because the IntervalProbabilities stores the lower and gap transition probabilities, fetching the upper bound requires allocation and computation.

source
IntervalMDP.upperMethod
upper(p::OrthogonalIntervalProbabilities, l, i, j)

Return the upper bound transition probabilities from a source state or source/action pair to a target state.

Note

It is not recommended to use this function for the hot loop of O-maximization. Because the IntervalProbabilities stores the lower and gap transition probabilities, fetching the upper bound requires allocation and computation.

source
IntervalMDP.gapMethod
gap(p::OrthogonalIntervalProbabilities, l)

Return the gap between upper and lower bound transition probabilities from a source state or source/action pair to a target axis.

source
IntervalMDP.gapMethod
gap(p::OrthogonalIntervalProbabilities, l, i, j)

Return the gap between upper and lower bound transition probabilities from a source state or source/action pair to a target state.

source
IntervalMDP.sum_lowerMethod
sum_lower(p::OrthogonalIntervalProbabilities, l)

Return the sum of lower bound transition probabilities from a source state or source/action pair to all target states on one axis. This is useful in efficiently implementing O-maximization, where we start with a lower bound probability assignment and iteratively, according to the ordering, adding the gap until the sum of probabilities is 1.

source
IntervalMDP.sum_lowerMethod
sum_lower(p::OrthogonalIntervalProbabilities, l, j)

Return the sum of lower bound transition probabilities from a source state or source/action pair to all target states. This is useful in efficiently implementing O-maximization, where we start with a lower bound probability assignment and iteratively, according to the ordering, adding the gap until the sum of probabilities is 1.

source
IntervalMDP.num_sourceMethod
num_source(p::OrthogonalIntervalProbabilities)

Return the number of source states or source/action pairs.

source
IntervalMDP.num_targetMethod
num_target(p::OrthogonalIntervalProbabilities)

Return the number of target states along each marginal.

source
IntervalMDP.axes_sourceMethod
axes_source(p::OrthogonalIntervalProbabilities)

Return the valid range of indices for the source states or source/action pairs.

source

Mixtures of marginal interval ambiguity sets

IntervalMDP.MixtureIntervalProbabilitiesType
MixtureIntervalProbabilities{N, P <: OrthogonalIntervalProbabilities, Q <: IntervalProbabilities}

A tuple of OrthogonalIntervalProbabilities for independent transition probabilities in a mixture that all share the same source/action pairs, and target states. See OrthogonalIntervalProbabilities for more information on the structure of the transition probabilities for each model in the mixture. The mixture is weighted by an IntervalProbabilities ambiguity set, called weighting_probs.

Fields

  • mixture_probs::NTuple{N, P}: A tuple of OrthogonalIntervalProbabilities transition probabilities along each axis.
  • weighting_probs::Q: The weighting ambiguity set for the mixture.

Examples

Below is a simple example of a mixture of two OrthogonalIntervalProbabilities with one dimension and the same source/action pairs and target states, and a weighting ambiguity set.

prob1 = OrthogonalIntervalProbabilities(
    (
        IntervalProbabilities(;
            lower = [
                0.0 0.5
                0.1 0.3
                0.2 0.1
            ],
            upper = [
                0.5 0.7
                0.6 0.5
                0.7 0.3
            ],
        ),
    ),
    (Int32(2),),
)
prob2 = OrthogonalIntervalProbabilities(
    (
        IntervalProbabilities(;
            lower = [
                0.1 0.4
                0.2 0.2
                0.3 0.0
            ],
            upper = [
                0.4 0.6
                0.5 0.4
                0.6 0.2
            ],
        ),
    ),
    (Int32(2),),
)
weighting_probs = IntervalProbabilities(; lower = [
    0.3 0.5
    0.4 0.3
], upper = [
    0.8 0.7
    0.7 0.5
])
mixture_prob = MixtureIntervalProbabilities((prob1, prob2), weighting_probs)
source
IntervalMDP.num_sourceMethod
num_source(p::MixtureIntervalProbabilities)

Return the number of source states or source/action pairs.

source
IntervalMDP.axes_sourceMethod
axes_source(p::MixtureIntervalProbabilities)

Return the valid range of indices for the source states or source/action pairs.

source
IntervalMDP.mixture_probsFunction
mixture_probs(p::MixtureIntervalProbabilities)

Return the tuple of OrthogonalIntervalProbabilities transition probabilities.

source
mixture_probs(p::MixtureIntervalProbabilities, k)

Return $k$-th OrthogonalIntervalProbabilities transition probabilities.

source

Labelling of IMDP states to Automaton alphabet

IntervalMDP.LabellingFunctionType
struct LabellingFunction{
    T  <:Unsigned, 
    VT <: AbstractVector{T}
}

A type representing the labelling of IMDP states into DFA inputs.

Formally, let $L : S => 2^{AP}$ be a labelling function, where

  • $S$ is the set of IMDP states, and
  • $2^{AP}$ is the power set of atomic propositions

Then the LabellingFunction type is defined as vector which stores the mapping.

Fields

  • map::VT: mapping function where indices are IMDP states and stored values are DFA inputs.
  • num_inputs::Int32: number of IMDP states accounted for in mapping.
  • num_outputs::Int32: number of DFA inputs accounted for in mapping.
source
IntervalMDP.mappingMethod
mapping(labelling_func::LabellingFunction)

Return the mapping array of the labelling function.

source

Transition function for DFA

IntervalMDP.TransitionFunctionType
struct TransitionFunction{
    T<:Unsigned, 
    MT <: AbstractMatrix{T}
}

A type representing the determininistic transition function of a DFA.

Formally, let $T : |2^{AP}| \times |Z| => |Z|$ be a transition function, where

  • $Z$ is the set of DFA states, and
  • $2^{AP}$ is the power set of atomic propositions

Then the TransitionFunction type is defined as matrix which stores the mapping. The row indices are the alphabet indices and the column indices represent the states.

Fields

  • transition::MT: transition function.
source
IntervalMDP.transitionMethod
transition(transition_func::TransitionFunction)

Return the transition matrix of the transition function.

source