Nondeterministic Finite automata

Types

All typedefs and enums used by algorithms operating on NFAs.

namespace mata

Main namespace including structs and algorithms for all automata.

In particular, this includes:

  1. Alphabets,

  2. Formula graphs and nodes,

  3. Mintermization,

  4. Closed sets.

namespace nfa

Nondeterministic Finite Automata including structures, transitions and algorithms.

In particular, this includes:

  1. Structures (Automaton, Transitions, Results, Delta),

  2. Algorithms (operations, checks, tests),

  3. Constructions.

Other algorithms are included in mata::nfa::Plumbing (simplified API for, e.g., binding) and mata::nfa::algorithms (concrete implementations of algorithms, such as for complement).

Typedefs

using State = unsigned long
using StateSet = mata::utils::OrdVector<State>
using StateRenaming = std::unordered_map<State, State>
using ParameterMap = std::unordered_map<std::string, std::string>

Map of additional parameter name and value pairs.

Used by certain functions for specifying some additional parameters in the following format:

ParameterMap {
    { "algorithm", "classical" },
    { "minimize", "true" }
}

Enums

enum class EpsilonClosureOpt : unsigned

Values:

enumerator NONE

No epsilon closure.

enumerator BEFORE

Epsilon closure before the transition.

enumerator AFTER

Epsilon closure after the transition.

enumerator BEFORE_AND_AFTER

Epsilon closure before and after the transition.

enum class ProductFinalStateCondition

Values:

enumerator AND

Both original states have to be final.

enumerator OR

At least one of the original states has to be final.

Variables

const std::string TYPE_NFA
Symbol EPSILON = Limits::max_symbol

A non-deterministic finite automaton.

An epsilon symbol which is now defined as the maximal value of data type used for symbols.

struct Limits
#include <types.hh>

Public Static Attributes

static const State min_state = std::numeric_limits<State>::min()
static const State max_state = std::numeric_limits<State>::max()
static const Symbol min_symbol = std::numeric_limits<Symbol>::min()
static const Symbol max_symbol = std::numeric_limits<Symbol>::max()
struct Run
#include <types.hh>

Public Members

Word word = {}

A finite-length word.

std::vector<State> path = {}

A finite-length path through automaton.

NFA

NFA class and its methods.

namespace mata

Main namespace including structs and algorithms for all automata.

In particular, this includes:

  1. Alphabets,

  2. Formula graphs and nodes,

  3. Mintermization,

  4. Closed sets.

namespace nfa

Nondeterministic Finite Automata including structures, transitions and algorithms.

In particular, this includes:

  1. Structures (Automaton, Transitions, Results, Delta),

  2. Algorithms (operations, checks, tests),

  3. Constructions.

Other algorithms are included in mata::nfa::Plumbing (simplified API for, e.g., binding) and mata::nfa::algorithms (concrete implementations of algorithms, such as for complement).

Typedefs

template<typename ...Ts>
using conjunction = std::is_same<bool_pack<true, Ts::value...>, bool_pack<Ts::value..., true>>

Check that for all values in a pack Ts are ‘true’.

template<typename T, typename ...Ts>
using AreAllOfType = typename conjunction<std::is_same<Ts, T>...>::type

Check that all types in a sequence of parameters Ts are of type T.

Functions

template<typename ...Nfas, typename = AreAllOfType<const Nfa&, Nfas...>>
inline OnTheFlyAlphabet create_alphabet(const Nfas&... nfas)

Create alphabet from variadic number of NFAs given as arguments.

in] Nfas Type Nfa.

Tparam :

Parameters:

nfas[in] NFAs to create alphabet from.

Returns:

Created alphabet.

OnTheFlyAlphabet create_alphabet(const std::vector<std::reference_wrapper<const Nfa>> &nfas)

Create alphabet from a vector of NFAs.

Parameters:

nfas[in] Vector of NFAs to create alphabet from.

Returns:

Created alphabet.

OnTheFlyAlphabet create_alphabet(const std::vector<std::reference_wrapper<Nfa>> &nfas)

Create alphabet from a vector of NFAs.

Parameters:

nfas[in] Vector of NFAs to create alphabet from.

Returns:

Created alphabet.

OnTheFlyAlphabet create_alphabet(const std::vector<Nfa*> &nfas)

Create alphabet from a vector of NFAs.

Parameters:

nfas[in] Vector of pointers to NFAs to create alphabet from.

Returns:

Created alphabet.

OnTheFlyAlphabet create_alphabet(const std::vector<const Nfa*> &nfas)

Create alphabet from a vector of NFAs.

Parameters:

nfas[in] Vector of pointers to NFAs to create alphabet from.

Returns:

Created alphabet.

Nfa union_nondet(const Nfa &lhs, const Nfa &rhs)

Compute non-deterministic union.

Does not add epsilon transitions, just unites initial and final states.

Returns:

Non-deterministic union of lhs and rhs.

Nfa union_det_complete(const Nfa &lhs, const Nfa &rhs)

Compute union of two complete deterministic NFAs.

Perserves determinism.

The union is computed by product construction with OR condition on the final states.

Parameters:
  • lhs – First complete deterministic automaton.

  • rhs – Second complete deterministic automaton.

Nfa product(const Nfa &lhs, const Nfa &rhs, ProductFinalStateCondition final_condition = ProductFinalStateCondition::AND, Symbol first_epsilon = EPSILON, std::unordered_map<std::pair<State, State>, State> *prod_map = nullptr)

Compute product of two NFAs with OR condition on the final states.

Automata must share alphabets. //TODO: this is not implemented yet.

Parameters:
  • lhs – First NFA.

  • rhs – Second NFA.

  • final_condition – Condition for a product state to be final.

    • AND: both original states have to be final.

    • OR: at least one of the original states has to be final.

  • first_epsilon – Smallest epsilon symbol. //TODO: this should eventually be taken from the alphabet as anything larger than the largest symbol?

  • prod_map – Mapping of pairs of the original states (lhs_state, rhs_state) to new product states (not used internally, allocated only when !=nullptr, expensive).

Nfa lang_difference(const Nfa &nfa_included, const Nfa &nfa_excluded, std::optional<std::function<bool(const Nfa&, const Nfa&, const StateSet&, const StateSet&, const State, const Nfa&)>> macrostate_discover = std::nullopt)

Compute a language difference as nfa_included \ nfa_excluded.

Computed as a lazy intersection of nfa_included and a complement of nfa_excluded. The NFAs are lazily determinized and the complement is constructed lazily as well, guided by nfa_included.

Todo:

: TODO: Add support for specifying first epsilon symbol and compute epsilon closure during determinization.

Parameters:
  • nfa_included[in] NFA to include in the difference.

  • nfa_excluded[in] NFA to exclude from the difference.

  • macrostate_discover[in] Callback event handler for discovering a new macrostate in the language difference automaton for the first time. Return true if the computation should continue, and false if the computation should stop and return only the NFA for the language difference constructed so far. The parameters are: const Nfa& nfa_included, const Nfa& nfa_excluded, const StateSet& macrostate_included_state_set, const StateSet& macrostate_excluded_state_set, const State macrostate, const Nfa& nfa_lang_difference.

Nfa intersection(const Nfa &lhs, const Nfa &rhs, const Symbol first_epsilon = EPSILON, std::unordered_map<std::pair<State, State>, State> *prod_map = nullptr)

Compute intersection of two NFAs.

Both automata can contain ε-transitions. The product preserves the ε-transitions, i.e., for each each product state (s, t) withs -ε-> p, (s, t) -ε-> (p, t) is created, and vice versa.

Automata must share alphabets. //TODO: this is not implemented yet.

Parameters:
  • lhs[in] First NFA to compute intersection for.

  • rhs[in] Second NFA to compute intersection for.

  • first_epsilon[in] smallest epsilon. //TODO: this should eventually be taken from the alphabet as anything larger than the largest symbol?

  • prod_map[out] Mapping of pairs of the original states (lhs_state, rhs_state) to new product states (not used internally, allocated only when !=nullptr, expensive).

Returns:

NFA as a product of NFAs lhs and rhs with ε-transitions preserved.

Nfa concatenate(const Nfa &lhs, const Nfa &rhs, bool use_epsilon = false, StateRenaming *lhs_state_renaming = nullptr, StateRenaming *rhs_state_renaming = nullptr)

Concatenate two NFAs.

Supports epsilon symbols when use_epsilon is set to true.

Parameters:
  • lhs[in] First automaton to concatenate.

  • rhs[in] Second automaton to concatenate.

  • use_epsilon[in] Whether to concatenate over epsilon symbol.

  • lhs_state_renaming[out] Map mapping lhs states to result states.

  • rhs_state_renaming[out] Map mapping rhs states to result states.

Returns:

Concatenated automaton.

Nfa complement(const Nfa &aut, const Alphabet &alphabet, const ParameterMap &params = {{"algorithm", "classical"}})

Compute automaton accepting complement of aut.

Parameters:
  • aut[in] Automaton whose complement to compute.

  • alphabet[in] Alphabet used for complementation.

  • params[in] Optional parameters to control the complementation algorithm:

    • ”algorithm”:

      • ”classical”: The classical algorithm determinizes the automaton, makes it complete, and swaps final and non-final states.

      • ”brzozowski”: The Brzozowski algorithm determinizes the automaton using Brzozowski minimization, makes it complete, and swaps final and non-final states.

Returns:

Complemented automaton.

Nfa complement(const Nfa &aut, const utils::OrdVector<Symbol> &symbols, const ParameterMap &params = {{"algorithm", "classical"}})

Compute automaton accepting complement of aut.

This overloaded version complements over an already created ordered set of symbols instead of an alphabet. This is a more efficient solution in case you already have symbols precomputed or want to complement multiple automata over the same set of symbols: the function does not need to compute the ordered set of symbols from the alphabet again (and for each automaton).

Parameters:
  • aut[in] Automaton whose complement to compute.

  • symbols[in] Symbols to complement over.

  • params[in] Optional parameters to control the complementation algorithm:

    • ”algorithm”:

      • ”classical”: The classical algorithm determinizes the automaton, makes it complete, and swaps final and non-final states.

      • ”brzozowski”: The Brzozowski algorithm determinizes the automaton using Brzozowski minimization, makes it complete, and swaps final and non-final states.

Returns:

Complemented automaton.

Nfa minimize(const Nfa &aut, const ParameterMap &params = {{"algorithm", "brzozowski"}})

Compute minimal deterministic automaton.

Parameters:
  • aut[in] Automaton whose minimal version to compute.

  • params[in] Optional parameters to control the minimization algorithm:

    • ”algorithm”: “brzozowski”

Returns:

Minimal deterministic automaton.

Nfa determinize(const Nfa &aut, std::unordered_map<StateSet, State> *subset_map = nullptr, std::optional<std::function<bool(const Nfa&, const State, const StateSet&)>> macrostate_discover = std::nullopt)

Determinize automaton.

Todo:

: TODO: Add support for specifying first epsilon symbol and compute epsilon closure during determinization.

Parameters:
  • aut[in] Automaton to determinize.

  • subset_map[out] Map that maps sets of states of input automaton to states of determinized automaton.

  • macrostate_discover[in] Callback event handler for discovering a new macrostate for the first time. The parameters are the determinized NFA constructed so far, the current macrostate, and the set of the original states corresponding to the macrostate. Return true if the determinization should continue, and false if the determinization should stop and return only the determinized NFA constructed so far.

Returns:

Determinized automaton.

Nfa reduce(const Nfa &aut, StateRenaming *state_renaming = nullptr, const ParameterMap &params = {{"algorithm", "simulation"}, {"type", "after"}, {"direction", "forward"}})

Reduce the size of the automaton.

Parameters:
  • aut[in] Automaton to reduce.

  • state_renaming[out] Mapping of original states to reduced states.

  • params[in] Optional parameters to control the reduction algorithm:

    • ”algorithm”: “simulation”, “residual”, and options to parametrize residual reduction, not utilized in simulation

    • ”type”: “after”, “with”,

    • ”direction”: “forward”, “backward”.

Returns:

Reduced automaton.

bool is_included(const Nfa &smaller, const Nfa &bigger, Run *cex, const Alphabet *alphabet = nullptr, const ParameterMap &params = {{"algorithm", "antichains"}})

Checks inclusion of languages of two NFAs: smaller and bigger (smaller <= bigger).

Parameters:
  • smaller[in] First automaton to concatenate.

  • bigger[in] Second automaton to concatenate.

  • cex[out] Counterexample for the inclusion.

  • alphabet[in] Alphabet of both NFAs to compute with.

  • params[in] Optional parameters to control the equivalence check algorithm:

    • ”algorithm”: “naive”, “antichains” (Default: “antichains”)

Returns:

True if smaller is included in bigger, false otherwise.

inline bool is_included(const Nfa &smaller, const Nfa &bigger, const Alphabet *const alphabet = nullptr, const ParameterMap &params = {{"algorithm", "antichains"}})

Checks inclusion of languages of two NFAs: smaller and bigger (smaller <= bigger).

Parameters:
  • smaller[in] First automaton to concatenate.

  • bigger[in] Second automaton to concatenate.

  • alphabet[in] Alphabet of both NFAs to compute with.

  • params[in] Optional parameters to control the equivalence check algorithm:

    • ”algorithm”: “naive”, “antichains” (Default: “antichains”)

Returns:

True if smaller is included in bigger, false otherwise.

bool are_equivalent(const Nfa &lhs, const Nfa &rhs, const Alphabet *alphabet, const ParameterMap &params = {{"algorithm", "antichains"}})

Perform equivalence check of two NFAs: lhs and rhs.

Parameters:
  • lhs[in] First automaton to concatenate.

  • rhs[in] Second automaton to concatenate.

  • alphabet[in] Alphabet of both NFAs to compute with.

  • params[[in] Optional parameters to control the equivalence check algorithm:

    • ”algorithm”: “naive”, “antichains” (Default: “antichains”)

Returns:

True if lhs and rhs are equivalent, false otherwise.

bool are_equivalent(const Nfa &lhs, const Nfa &rhs, const ParameterMap &params = {{"algorithm", "antichains"}})

Perform equivalence check of two NFAs: lhs and rhs.

The current implementation of Nfa does not accept input alphabet. For this reason, an alphabet has to be created from all transitions each time an operation on alphabet is called. When calling this function, the alphabet has to be computed first.

Hence, this function is less efficient than its alternative taking already defined alphabet as its parameter. That way, alphabet has to be computed only once, as opposed to the current ad-hoc construction of the alphabet. The use of the alternative with defined alphabet should be preferred.

Parameters:
  • lhs[in] First automaton to concatenate.

  • rhs[in] Second automaton to concatenate.

  • params[in] Optional parameters to control the equivalence check algorithm:

    • ”algorithm”: “naive”, “antichains” (Default: “antichains”)

Returns:

True if lhs and rhs are equivalent, false otherwise.

Nfa revert(const Nfa &aut)
Nfa fragile_revert(const Nfa &aut)
Nfa simple_revert(const Nfa &aut)
Nfa somewhat_simple_revert(const Nfa &aut)
Nfa remove_epsilon(const Nfa &aut, Symbol epsilon = EPSILON)
Run encode_word(const Alphabet *alphabet, const std::vector<std::string> &input)

Encodes a vector of strings (each corresponding to one symbol) into a Word instance.

utils::OrdVector<Symbol> get_symbols_to_work_with(const nfa::Nfa &nfa, const Alphabet *const shared_alphabet = nullptr)

Get the set of symbols to work with during operations.

Parameters:

shared_alphabet[in] Optional alphabet shared between NFAs passed as an argument to a function.

std::optional<Word> get_word_from_lang_difference(const Nfa &nfa_included, const Nfa &nfa_excluded)

Get any arbitrary accepted word in the language difference of nfa_included without nfa_excluded.

The language difference automaton is lazily constructed without computing the whole determinized automata and the complememt of nfa_excluded. The algorithm returns an arbitrary word from the language difference constructed until the first macrostate with a final state in the original states in nfa_included and without any corresponding final states in nfa_excluded is encountered.

Parameters:
  • nfa_included[in] NFA to include in the language difference.

  • nfa_excluded[in] NFA to exclude in the language difference. TODO: Support lazy epsilon closure?

Pre:

The automaton does not contain any epsilon transitions.

Returns:

An arbitrary word from the language difference, or std::nullopt if the language difference automaton is universal on the set of symbols from transitions of nfa_included.

template<bool...>
struct bool_pack
#include <nfa.hh>

Pack of bools for reasoning about a sequence of parameters.

class Nfa
#include <nfa.hh>

A class representing an NFA.

Subclassed by Nft

Public Functions

inline explicit Nfa(Delta delta = {}, utils::SparseSet<State> initial_states = {}, utils::SparseSet<State> final_states = {}, Alphabet *alphabet = nullptr)
inline explicit Nfa(const size_t num_of_states, utils::SparseSet<State> initial_states = {}, utils::SparseSet<State> final_states = {}, Alphabet *alphabet = nullptr)

Construct a new explicit NFA with num_of_states states and optionally set initial and final states.

Parameters:

num_of_states[in] Number of states for which to preallocate Delta.

Nfa(const Nfa &other) = default

Construct a new explicit NFA from other NFA.

inline Nfa(Nfa &&other) noexcept
Nfa &operator=(const Nfa &other) = default
Nfa &operator=(Nfa &&other) noexcept
State add_state()

Add a new (fresh) state to the automaton.

Returns:

The newly created state.

State add_state(State state)

Add state state to delta if state is not in delta yet.

Returns:

The requested state.

State insert_word(State source, const Word &word, State target)

Inserts a word into the NFA from a source state source to a target state target.

Creates new states along the path of the word.

Parameters:
  • source – The source state where the word begins. It must already be a part of the automaton.

  • word – The nonempty word to be inserted into the NFA.

  • target – The target state where the word ends.

Returns:

The state target where the inserted word ends.

State insert_word(State source, const Word &word)

Inserts a word into the NFA from a source state source to a new target state.

Creates new states along the path of the word.

Parameters:
  • source – The source state where the word begins. It must already be a part of the automaton.

  • word – The nonempty word to be inserted into the NFA.

Returns:

The newly created target state where the inserted word ends.

size_t num_of_states() const

Get the current number of states in the whole automaton.

This includes the initial and final states as well as states in the transition relation.

Returns:

The number of states.

Nfa &unify_initial(bool force_new_state = false)

Unify initial states into a single new initial state.

Parameters:

force_new_state[in] Whether to force creating a new state even when initial states are already unified.

Returns:

this after unification.

Nfa &unify_final(bool force_new_state = false)

Unify final states into a single new final state.

Parameters:

force_new_state[in] Whether to force creating a new state even when final states are already unified.

Returns:

this after unification.

inline Nfa &swap_final_nonfinal()

Swap final and non-final states in-place.

inline bool is_state(const State &state_to_check) const
void clear()

Clear the underlying NFA to a blank NFA.

The whole NFA is cleared, each member is set to its zero value.

bool is_identical(const Nfa &aut) const

Check if this is exactly identical to aut.

This is exact equality of automata, including state numbering (so even stronger than isomorphism), essentially only useful for testing purposes.

Returns:

True if automata are exactly identical, false otherwise.

StateSet get_reachable_states() const

Get set of reachable states.

Reachable states are states accessible from any initial state.

Returns:

Set of reachable states. TODO: with the new get_useful_states, it might be useless now.

StateSet get_terminating_states() const

Get set of terminating states.

Terminating states are states leading to any final state.

Returns:

Set of terminating states. TODO: with the new get_useful_states, it might be useless now.

BoolVector get_useful_states() const

Get the useful states using a modified Tarjan’s algorithm.

A state is useful if it is reachable from an initial state and can reach a final state.

Returns:

BoolVector Bool vector whose ith value is true iff the state i is useful.

void tarjan_scc_discover(const TarjanDiscoverCallback &callback) const

Tarjan’s SCC discover algorihm.

Parameters:

callback – Callbacks class to instantiate callbacks for the Tarjan’s algorithm.

Nfa &trim(StateRenaming *state_renaming = nullptr)

Remove inaccessible (unreachable) and not co-accessible (non-terminating) states in-place.

Remove states which are not accessible (unreachable; state is accessible when the state is the endpoint of a path starting from an initial state) or not co-accessible (non-terminating; state is co-accessible when the state is the starting point of a path ending in a final state).

Parameters:

state_renaming[out] Mapping of trimmed states to new states.

Returns:

this after trimming.

Nfa decode_utf8() const

Decodes automaton from UTF-8 encoding.

Method removes unreachable states from delta.

Returns:

Decoded automaton.

std::vector<State> distances_from_initial() const

Returns vector ret where ret[q] is the length of the shortest path from any initial state to q.

std::vector<State> distances_to_final() const

Returns vector ret where ret[q] is the length of the shortest path from q to any final state.

Run get_shortest_accepting_run_from_state(State q, const std::vector<State> &distances_to_final) const

Get some shortest accepting run from state q.

Assumes that q is a state of this automaton and that there is some accepting run from q

Parameters:

distances_to_final – Vector of the lengths of the shortest runs from states (can be computed using distances_to_final())

void remove_epsilon(Symbol epsilon = EPSILON)

Remove epsilon transitions from the automaton.

Nfa &concatenate(const Nfa &aut)

In-place concatenation.

Nfa &unite_nondet_with(const Nfa &aut)

In-place nondeterministic union with aut.

Does not add epsilon transitions, just unites initial and final states.

Nfa get_one_letter_aut(Symbol abstract_symbol = 'x') const

Unify transitions to create a directed graph with at most a single transition between two states.

Parameters:

abstract_symbol[in] Abstract symbol to use for transitions in digraph.

Returns:

An automaton representing a directed graph.

inline bool is_epsilon(Symbol symbol) const

Check whether symbol is epsilon symbol or not.

Parameters:

symbolSymbol to check.

Returns:

True if the passed symbol is epsilon, false otherwise.

void get_one_letter_aut(Nfa &result) const

Unify transitions to create a directed graph with at most a single transition between two states.

Parameters:

result[out] An automaton representing a directed graph.

std::string print_to_dot(bool decode_ascii_chars = false, bool use_intervals = false, int max_label_length = -1) const

Prints the automaton in DOT format.

Parameters:
  • decode_ascii_chars[in] Whether to use ASCII characters for the output.

  • use_intervals[in] Whether to use intervals (e.g. [1-3] instead of 1,2,3) for labels.

  • max_label_length[in] Maximum label length for the output (-1 means no limit, 0 means no labels). If the label is longer than max_label_length, it will be truncated, with full label displayed on hover.

Returns:

automaton in DOT format

void print_to_dot(std::ostream &output, bool decode_ascii_chars = false, bool use_intervals = false, int max_label_length = -1) const

Prints the automaton to the output stream in DOT format.

Parameters:
  • decode_ascii_chars[in] Whether to use ASCII characters for the output.

  • use_intervals[in] Whether to use intervals (e.g. [1-3] instead of 1,2,3) for labels.

  • max_label_length[in] Maximum label length for the output (-1 means no limit, 0 means no labels). If the label is longer than max_label_length, it will be truncated, with full label displayed on hover.

void print_to_dot(const std::string &filename, bool decode_ascii_chars = false, bool use_intervals = false, int max_label_length = -1) const

Prints the automaton to the file in DOT format.

Parameters:
  • filename – Name of the file to print the automaton to

  • decode_ascii_chars[in] Whether to use ASCII characters for the output.

  • use_intervals[in] Whether to use intervals (e.g. [1-3] instead of 1,2,3) for labels.

  • max_label_length[in] Maximum label length for the output (-1 means no limit, 0 means no labels). If the label is longer than max_label_length, it will be truncated, with full label displayed on hover.

std::string print_to_mata(const Alphabet *alphabet = nullptr) const

Prints the automaton in mata format.

If you need to parse the automaton again, use IntAlphabet in construct()

Parameters:

alphabet[in] If specified, translates the symbols to their symbol names in the alphabet.

Returns:

automaton in mata format

void print_to_mata(std::ostream &output, const Alphabet *alphabet = nullptr) const

Prints the automaton to the output stream in mata format.

If you need to parse the automaton again, use IntAlphabet in construct()

Parameters:

alphabet[in] If specified, translates the symbols to their symbol names in the alphabet.

void print_to_mata(const std::string &filename, const Alphabet *alphabet = nullptr) const

Prints the automaton to the file in mata format.

If you need to parse the automaton again, use IntAlphabet in construct()

Parameters:
  • alphabet[in] If specified, translates the symbols to their symbol names in the alphabet.

  • filename – Name of the file to print the automaton to

StateSet post(const StateSet &states, const Symbol symbol, EpsilonClosureOpt epsilon_closure_opt = EpsilonClosureOpt::NONE) const

Get the set of states reachable from the given set of states over the given symbol.

TODO: Relict from VATA. What to do with inclusion/ universality/ this post function? Revise all of them.

Parameters:
  • states – Set of states to compute the post set from.

  • symbolSymbol to compute the post set for.

  • epsilon_closure_opt – Epsilon closure option. Perform epsilon closure before and/or after the post operation.

Returns:

Set of states reachable from the given set of states over the given symbol.

inline StateSet post(const State state, const Symbol symbol, EpsilonClosureOpt epsilon_closure_opt) const

Get the set of states reachable from the given state over the given symbol.

Parameters:
  • state – A state to compute the post set from.

  • symbolSymbol to compute the post set for.

  • epsilon_closure_opt – Epsilon closure option. Perform epsilon closure before and/or after the post operation.

Returns:

Set of states reachable from the given state over the given symbol.

inline const StateSet &post(const State state, const Symbol symbol) const

Returns a reference to targets (states) reachable from the given state over the given symbol.

This is an optimized shortcut for post(state, symbol, EpsilonClosureOpt::NONE).

Parameters:
  • state – A state to compute the post set from.

  • symbolSymbol to compute the post set for.

Returns:

Set of states reachable from the given state over the given symbol.

bool is_lang_empty(Run *cex = nullptr) const

Check whether the language of NFA is empty.

Currently calls is_lang_empty_scc if cex is null

Parameters:

cex[out] Counter-example path for a case the language is not empty.

Returns:

True if the language is empty, false otherwise.

bool is_lang_empty_scc() const

Check if the language is empty using Tarjan’s SCC discover algorithm.

Returns:

Language empty <-> True

bool is_deterministic() const

Test whether an automaton is deterministic.

I.e., whether it has exactly one initial state and every state has at most one outgoing transition over every symbol. Checks the whole automaton, not only the reachable part

bool is_complete(Alphabet const *alphabet = nullptr) const

Test for automaton completeness with regard to an alphabet.

An automaton is complete if every reachable state has at least one outgoing transition over every symbol.

bool is_acyclic() const

Is the automaton graph acyclic?

Used for checking language finiteness.

Returns:

true <-> Automaton graph is acyclic.

bool is_flat() const

Is the automaton flat?

Flat automaton is an NFA whose every SCC is a simple loop. Basically each state in an SCC has at most one successor within this SCC.

Returns:

true <-> Automaton graph is flat.

void fill_alphabet(mata::OnTheFlyAlphabet &alphabet_to_fill) const

Fill alphabet_to_fill with symbols from nfa.

Parameters:
  • nfa[in] NFA with symbols to fill alphabet_to_fill with.

  • alphabet_to_fill[out] Alphabet to be filled with symbols from nfa.

bool is_universal(const Alphabet &alphabet, Run *cex = nullptr, const ParameterMap &params = {{"algorithm", "antichains"}}) const

Check whether the language of the automaton is universal.

Parameters:
  • alphabetAlphabet to use for checking the universality.

  • cex – Counterexample path for a case the language is not universal.

  • params – Optional parameters to control the universality check algorithm:

    • ”algorithm”:

      • ”antichains”: The algorithm uses antichains to check the universality.

      • ”naive”: The algorithm uses the naive approach to check the universality.

Returns:

True if the language of the automaton is universal, false otherwise.

bool is_universal(const Alphabet &alphabet, const ParameterMap &params) const

Check whether the language of the automaton is universal.

Parameters:
  • alphabetAlphabet to use for checking the universality.

  • params – Optional parameters to control the universality check algorithm:

    • ”algorithm”:

      • ”antichains”: The algorithm uses antichains to check the universality.

      • ”naive”: The algorithm uses the naive approach to check the universality.

Returns:

True if the language of the automaton is universal, false otherwise.

bool is_in_lang(const Run &word, bool use_epsilon = false, bool match_prefix = false) const

Check whether a run over the word (or its prefix) is in the language of an automaton.

Parameters:
  • word – The run to check.

  • use_epsilon – Whether the automaton uses epsilon transitions.

  • match_prefix – Whether to also match the prefix of the word.

Returns:

True if the run (or its prefix) is in the language of the automaton, false otherwise.

inline bool is_in_lang(const Word &word, const bool use_epsilon = false, const bool match_prefix = false)

Check whether a word (or its prefix) is in the language of an automaton.

Parameters:
  • word – The word to check.

  • use_epsilon – Whether the automaton uses epsilon transitions.

  • match_prefix – Whether to also match the prefix of the word.

Returns:

True if the word (or its prefix) is in the language of the automaton, false otherwise.

inline bool is_prefix_in_lang(const Run &word, const bool use_epsilon = false) const

Check whether a prefix of a run is in the language of an automaton.

Parameters:
  • word – The run to check.

  • use_epsilon – Whether the automaton uses epsilon transitions.

Returns:

True if the prefix of the run is in the language of the automaton, false otherwise.

inline bool is_prefix_in_lang(const Word &word, const bool use_epsilon = false) const

Check whether a prefix of a word is in the language of an automaton.

Parameters:
  • word – The word to check.

  • use_epsilon – Whether the automaton uses epsilon transitions.

Returns:

True if the prefix of the word is in the language of the automaton, false otherwise.

std::pair<Run, bool> get_word_for_path(const Run &run) const
std::set<Word> get_words(size_t max_length) const

Get the set of all words in the language of the automaton whose length is <= max_length.

If you have an automaton with finite language (can be checked using is_acyclic), you can get all words by calling get_words(aut.num_of_states())

std::optional<Word> get_word(const std::optional<Symbol> first_epsilon = EPSILON) const

Get any arbitrary accepted word in the language of the automaton.

The automaton is searched using DFS, returning a word for the first reached final state.

Parameters:

first_epsilon – If defined, all symbols >=first_epsilon are assumed to be epsilon and therefore are not in the returned word.

Returns:

std::optional<Word> Some word from the language. If the language is empty, returns std::nullopt.

std::optional<Word> get_word_from_complement(const Alphabet *alphabet = nullptr) const

Get any arbitrary accepted word in the language of the complement of the automaton.

The automaton is lazily determinized and made complete. The algorithm returns an arbitrary word from the complemented NFA constructed until the first macrostate without any final states in the original automaton is encountered.

Parameters:

alphabet[in] Alphabet to use for computing the complement. If nullptr, uses this->alphabet when defined, otherwise uses this->delta.get_used_symbols().

Pre:

The automaton does not contain any epsilon transitions. TODO: Support lazy epsilon closure?

Returns:

An arbitrary word from the complemented automaton, or std::nullopt if the automaton is universal on the chosen set of symbols for the complement.

bool make_complete(const Alphabet *alphabet = nullptr, std::optional<State> sink_state = std::nullopt)

Make NFA complete in place.

For each state 0,…,this->num_of_states()-1, add transitions with “missing” symbols from alphabet (symbols that do not occur on transitions from given state) to sink_state. If sink_state does not belong to the NFA, it is added to it, but only in the case that some transition to sink_state was added. In the case that NFA does not contain any states, this function does nothing.

Parameters:
  • alphabet[in] Alphabet to use for computing “missing” symbols. If nullptr, use this->alphabet when defined, otherwise use this->delta.get_used_symbols().

  • sink_state[in] The state into which new transitions are added. If std::nullopt, add a new sink state.

Returns:

true if a new transition was added to the NFA.

bool make_complete(const utils::OrdVector<Symbol> &symbols, std::optional<State> sink_state = std::nullopt)

Make NFA complete in place.

For each state 0,…,this->num_of_states()-1, add transitions with “missing” symbols from alphabet (symbols that do not occur on transitions from given state) to sink_state. If sink_state does not belong to the NFA, it is added to it, but only in the case that some transition to sink_state was added. In the case that NFA does not contain any states, this function does nothing.

This overloaded version is a more efficient version which does not need to compute the set of symbols to complete to from the alphabet. Prefer this version when you already have the set of symbols precomputed or plan to complete multiple automata over the same set of symbols.

Parameters:
  • symbols[in] Symbols to compute “missing” symbols from.

  • sink_state[in] The state into which new transitions are added. If std::nullopt, add a new sink state.

Returns:

true if a new transition was added to the NFA.

Nfa &complement_deterministic(const mata::utils::OrdVector<Symbol> &symbols, std::optional<State> sink_state = std::nullopt)

Complement deterministic automaton in-place by adding a sink state and swapping final and non-final states.

Parameters:
  • symbols[in] Symbols needed to make the automaton complete.

  • sink_state[in] State to be used as a sink state. Adds a new sink state when not specified.

Returns:

DFA complemented in-place.

Pre:

this is a deterministic automaton.

Public Members

Delta delta

For state q, delta[q] keeps the list of transitions ordered by symbols.

The set of states of this automaton are the numbers from 0 to the number of states minus one.

utils::SparseSet<State> initial = {}
utils::SparseSet<State> final = {}
Alphabet *alphabet = nullptr

The alphabet which can be shared between multiple automata. Key value store for additional attributes for the NFA. Keys are attribute names as strings and the value types are up to the user. For example, we can set up attributes such as “state_dict” for state dictionary attribute mapping states to their respective names, or “transition_dict” for transition dictionary adding a human-readable meaning to each transition.

std::unordered_map<std::string, void*> attributes = {}
struct TarjanDiscoverCallback
#include <nfa.hh>

Structure for storing callback functions (event handlers) utilizing Tarjan’s SCC discover algorithm.

Public Members

std::function<bool(State)> state_discover
std::function<bool(const std::vector<State>&, const std::vector<State>&)> scc_discover
std::function<void(State)> scc_state_discover
std::function<void(State, State)> succ_state_discover
namespace std

STL namespace.

Functions

std::ostream &operator<<(std::ostream &os, const mata::nfa::Transition &trans)
std::ostream &operator<<(std::ostream &os, const mata::nfa::Nfa &nfa)
template<>
struct hash<mata::nfa::Transition>
#include <nfa.hh>

Public Functions

inline size_t operator()(const mata::nfa::Transition &trans) const

Delta

The delta function of an NFA, which maps states and input symbols to sets of states.

namespace mata

Main namespace including structs and algorithms for all automata.

In particular, this includes:

  1. Alphabets,

  2. Formula graphs and nodes,

  3. Mintermization,

  4. Closed sets.

namespace nfa

Nondeterministic Finite Automata including structures, transitions and algorithms.

In particular, this includes:

  1. Structures (Automaton, Transitions, Results, Delta),

  2. Algorithms (operations, checks, tests),

  3. Constructions.

Other algorithms are included in mata::nfa::Plumbing (simplified API for, e.g., binding) and mata::nfa::algorithms (concrete implementations of algorithms, such as for complement).

class Delta
#include <delta.hh>

Delta is a data structure for representing transition relation.

Transition is represented as a triple Trans(source state, symbol, target state). Move is the part (symbol, target state), specified for a single source state. Its underlying data structure is vector of StatePost classes. Each index to the vector corresponds to one source state, that is, a number for a certain state is an index to the vector of state posts. Transition relation (delta) in Mata stores a set of transitions in a four-level hierarchical structure: Delta, StatePost, SymbolPost, and a set of target states. A vector of ‘StatePost’s indexed by a source states on top, where the StatePost for a state ‘q’ (whose number is ‘q’ and it is the index to the vector of ‘StatePost’s) stores a set of ‘Move’s from the source state ‘q’. Namely, ‘StatePost’ has a vector of ‘SymbolPost’s, where each ‘SymbolPost’ stores a symbol ‘a’ and a vector of target states of ‘a’-moves from state ‘q’. ‘SymbolPost’s are ordered by the symbol, target states are ordered by the state number.

Public Types

using const_iterator = std::vector<StatePost>::const_iterator

Public Functions

inline Delta()
Delta(const Delta &other) = default
Delta(Delta &&other) = default
inline explicit Delta(size_t n)
Delta &operator=(const Delta &other) = default
Delta &operator=(Delta &&other) = default
bool operator==(const Delta &other) const
inline void reserve(size_t n)
inline const StatePost &state_post(const State source) const

Get constant reference to the state post of source.

If we try to access a state post of a source which is present in the automaton as an initial/final state, yet does not have allocated space in Delta, an empty_post is returned. Hence, the function has no side effects (no allocation is performed; iterators remain valid).

Parameters:

source[in] – Source state of a state post to access.

Returns:

State post of source.

inline const StatePost &operator[](const State source) const

Get constant reference to the state post of source.

If we try to access a state post of a source which is present in the automaton as an initial/final state, yet does not have allocated space in Delta, an empty_post is returned. Hence, the function has no side effects (no allocation is performed; iterators remain valid).

Parameters:

source[in] – Source state of a state post to access.

Returns:

State post of source.

StatePost &mutable_state_post(State source)

Get mutable (non-constant) reference to the state post of source.

The function allows modifying the state post.

BEWARE, IT HAS A SIDE EFFECT.

If we try to access a state post of a source which is present in the automaton as an initial/final state, yet does not have allocated space in Delta, a new state post for source will be allocated along with all state posts for all previous states. This in turn may cause that the entire post data structure is re-allocated. Iterators to Delta will get invalidated. Use the constant ‘state_post()’ is possible. Or, to prevent the side effect from causing issues, one might want to make sure that posts of all states in the automaton are allocated, e.g., write an NFA method that allocate Delta for all states of the NFA.

Parameters:

source[in] – Source state of a state post to access.

Returns:

State post of source.

void defragment(const BoolVector &is_staying, const std::vector<State> &renaming)
template<typename ...Args>
inline StatePost &emplace_back(Args&&... args)
inline void clear()
inline void allocate(const size_t num_of_states)

Allocate state posts up to num_of_states states, creating empty StatePost for yet unallocated state posts.

Parameters:

num_of_states[in] Number of states in Delta to allocate state posts for. Have to be at least num_of_states() + 1.

inline size_t num_of_states() const
Returns:

Number of states in the whole Delta, including both source and target states.

inline bool uses_state(const State state) const

Check whether the state is used in Delta.

size_t num_of_transitions() const
Returns:

Number of transitions in Delta.

void add(State source, Symbol symbol, State target)
inline void add(const Transition &trans)
void remove(State source, Symbol symbol, State target)
inline void remove(const Transition &transition)
bool contains(State source, Symbol symbol, State target) const

Check whether Delta contains a passed transition.

bool contains(const Transition &transition) const

Check whether Delta contains a transition passed as a triple.

bool empty() const

Check whether automaton contains no transitions.

Returns:

True if there are no transitions in the automaton, false otherwise.

inline void append(const std::vector<StatePost> &post_vector)

Append post vector to the delta.

Parameters:

post_vector – Vector of posts to be appended.

std::vector<StatePost> renumber_targets(const std::function<State(State)> &target_renumberer) const

Copy posts of delta and apply a lambda update function on each state from targets.

IMPORTANT: In order to work properly, the lambda function needs to be monotonic, that is, the order of states in targets cannot change.

Parameters:

target_renumberer – Monotonic lambda function mapping states to different states.

Returns:

std::vector<Post> Copied posts.

void add(const State source, const Symbol symbol, const StateSet &targets)

Add transitions to multiple destinations.

Parameters:
  • source – From

  • symbolSymbol

  • targets – Set of states to

inline const_iterator cbegin() const
inline const_iterator cend() const
inline const_iterator begin() const
inline const_iterator end() const
Transitions transitions() const

Iterator over transitions represented as Transition instances.

std::vector<Transition> get_transitions_to(State state_to) const

Get transitions leading to state_to.

Operation is slow, traverses over all symbol posts.

Parameters:

state_to[in] – Target state for transitions to get.

Returns:

Transitions leading to state_to.

std::vector<Transition> get_transitions_between(State state_from, State state_to) const

Get transitions from state_from to state_to.

Operation is slow, traverses over all symbol posts.

Parameters:
  • state_from[in] – Source state.

  • state_from[in] – Target state.

Returns:

Transitions from source to state_to.

StateSet get_successors(State state) const

Get the set of states that are successors of the given state.

Parameters:

state[in] State from which successors are checked.

Returns:

Set of states that are successors of the given state.

const StateSet &get_successors(State state, Symbol symbol) const
StateSet get_successors(State state, Symbol symbol, EpsilonClosureOpt epsilon_closure_opt) const
StatePost::const_iterator epsilon_symbol_posts(State state, Symbol epsilon = EPSILON) const

Iterate over epsilon symbol posts under the given state.

Parameters:
  • state[in] State from which epsilon transitions are checked.

  • epsilon[in] User can define his favourite epsilon or used default.

Returns:

An iterator to SymbolPost with epsilon symbol. End iterator when there are no epsilon transitions.

void add_symbols_to(OnTheFlyAlphabet &target_alphabet) const

Expand target_alphabet by symbols from this delta.

The value of the already existing symbols will NOT be overwritten.

utils::OrdVector<Symbol> get_used_symbols() const

Get the set of symbols used on the transitions in the automaton.

Does not necessarily have to equal the set of symbols in the alphabet used by the automaton.

Returns:

Set of symbols used on the transitions.

utils::OrdVector<Symbol> get_used_symbols_vec() const
std::set<Symbol> get_used_symbols_set() const
utils::SparseSet<Symbol> get_used_symbols_sps() const
std::vector<bool> get_used_symbols_bv() const
BoolVector get_used_symbols_chv() const
Symbol get_max_symbol() const

Get the maximum non-epsilon used symbol.

Public Static Functions

static StatePost::const_iterator epsilon_symbol_posts(const StatePost &state_post, Symbol epsilon = EPSILON)

Iterate over epsilon symbol posts under the given state_post.

Parameters:
  • state_post[in] State post from which epsilon transitions are checked.

  • epsilon[in] User can define his favourite epsilon or used default.

Returns:

An iterator to SymbolPost with epsilon symbol. End iterator when there are no epsilon transitions.

Public Static Attributes

static const StatePost empty_state_post

Protected Attributes

std::vector<StatePost> state_posts_
class Transitions
#include <delta.hh>

Iterator over transitions represented as Transition instances.

It iterates over triples (State source, Symbol symbol, State target).

Public Functions

Transitions() = default
inline explicit Transitions(const Delta *delta)
Transitions(Transitions&&) = default
Transitions(Transitions&) = default
Transitions &operator=(Transitions&&) = default
Transitions &operator=(Transitions&) = default
const_iterator begin() const
const_iterator end() const

Private Members

const Delta *delta_
class const_iterator
#include <delta.hh>

Iterator over transitions.

Public Types

using iterator_category = std::forward_iterator_tag
using value_type = Transition
using difference_type = size_t
using pointer = Transition*
using reference = Transition&

Public Functions

inline const_iterator()
explicit const_iterator(const Delta &delta)
const_iterator(const Delta &delta, State current_state)
const_iterator(const const_iterator &other) noexcept = default
const_iterator(const_iterator&&) = default
inline const Transition &operator*() const
inline const Transition *operator->() const
const_iterator &operator++()
const const_iterator operator++(int)
const_iterator &operator=(const const_iterator &other) noexcept = default
const_iterator &operator=(const_iterator&&) = default
bool operator==(const const_iterator &other) const

Private Members

const Delta *delta_ = nullptr
size_t current_state_ = {}
StatePost::const_iterator state_post_it_ = {}
StateSet::const_iterator symbol_post_it_ = {}
bool is_end_ = {false}
Transition transition_ = {}
class Move
#include <delta.hh>

Move from a StatePost for a single source state, represented as a pair of symbol and target state target.

Public Functions

bool operator==(const Move&) const = default

Public Members

Symbol symbol
State target
class StatePost : private OrdVector<SymbolPost>
#include <delta.hh>

A data structure representing possible transitions over different symbols from a source state.

It is an ordered vector containing possible SymbolPost (i.e., pair of symbol and target states). SymbolPosts in the vector are ordered by symbols in SymbolPosts.

Public Types

using iterator = typename VectorType::iterator
using const_iterator = typename VectorType::const_iterator

Public Functions

StatePost(const StatePost&) = default
StatePost(StatePost&&) = default
StatePost &operator=(const StatePost&) = default
StatePost &operator=(StatePost&&) = default
bool operator==(const StatePost&) const = default
inline iterator find(const Symbol symbol)
inline const_iterator find(const Symbol symbol) const
const_iterator first_epsilon_it(Symbol first_epsilon) const

returns an iterator to the smallest epsilon, or end() if there is no epsilon

StateSet get_successors() const

Get the set of all target states in the StatePost.

Returns:

Set of all target states in the StatePost.

const StateSet &get_successors(Symbol symbol) const

Returns a reference to target states for a given symbol in the StatePost.

If there is no such symbol, a static empty set is returned.

Parameters:

symbolSymbol to get the successors for.

Returns:

Set of target states for the given symbol.

inline Moves moves() const

Iterator over all moves (over all labels) in StatePost represented as Move instances.

Moves moves(StatePost::const_iterator symbol_post_it, StatePost::const_iterator symbol_post_end) const

Iterator over specified moves in StatePost represented as Move instances.

Parameters:
  • symbol_post_it[in] First iterator over symbol posts to iterate over.

  • symbol_post_end[in] End iterator over symbol posts (which functions as an sentinel, is not iterated over).

Moves moves_epsilons(Symbol first_epsilon = EPSILON) const

Iterator over epsilon moves in StatePost represented as Move instances.

Moves moves_symbols(Symbol last_symbol = EPSILON - 1) const

Iterator over alphabet (normal) symbols (not over epsilons) in StatePost represented as Move instances.

size_t num_of_moves() const

Count the number of all moves in StatePost.

inline virtual const_iterator begin() const
inline virtual iterator begin()
inline virtual const_iterator end() const
inline virtual iterator end()
inline virtual const_iterator cbegin() const
inline virtual const_iterator cend() const
inline OrdVector()
inline explicit OrdVector(const VectorType &vec)
inline explicit OrdVector(const std::set<SymbolPost> &set)
inline explicit OrdVector(const T &set)
inline OrdVector(std::initializer_list<SymbolPost> list)
OrdVector(const OrdVector &rhs) = default
inline OrdVector(OrdVector &&other) noexcept
inline explicit OrdVector(const SymbolPost &key)
inline explicit OrdVector(InputIterator first, InputIterator last)
inline OrdVector &operator=(const OrdVector &other)
inline OrdVector &operator=(OrdVector &&other) noexcept
inline bool operator==(const OrdVector &rhs) const
inline void insert(iterator itr, const SymbolPost &x)
inline void insert(const SymbolPost &x)
inline virtual void insert(const OrdVector &vec)
inline virtual void reserve(size_t size)
inline virtual bool empty() const
inline virtual size_t size() const
inline const std::vector<SymbolPost> &to_vector() const
inline reference push_back(const SymbolPost &t)
inline reference push_back(SymbolPost &&t)
inline reference emplace_back(Args&&... args)
inline virtual const SymbolPost &front() const
inline virtual SymbolPost &front()
inline virtual const_reference back() const
inline virtual reference back()

Get reference to the last element in the vector.

Modifying the underlying value in the reference could break sortedness.

inline virtual void pop_back()
inline void filter(const Fun &&is_staying)
inline void clear()
inline virtual iterator erase(const_iterator pos)
inline virtual iterator erase(const_iterator first, const_iterator last)
inline size_t erase(const SymbolPost &k)

Remove k from sorted vector.

This function expects the vector to be sorted.

inline virtual const_iterator find(const SymbolPost &key) const
inline virtual iterator find(const SymbolPost &key)

Private Types

using super = utils::OrdVector<SymbolPost>
class Moves
#include <delta.hh>

Iterator over moves represented as Move instances.

It iterates over pairs (symbol, target) for the given StatePost.

Public Functions

Moves() = default
Moves(const StatePost &state_post, StatePost::const_iterator symbol_post_it, StatePost::const_iterator symbol_post_end)

construct moves iterating over a range symbol_post_it (including) to symbol_post_end (excluding).

Parameters:
  • state_post[in] State post to iterate over.

  • symbol_post_it[in] First iterator over symbol posts to iterate over.

  • symbol_post_end[in] End iterator over symbol posts (which functions as an sentinel; is not iterated over).

Moves(Moves&&) = default
Moves(Moves&) = default
Moves &operator=(Moves &&other) noexcept
Moves &operator=(const Moves &other) noexcept
const_iterator begin() const
const_iterator end() const

Private Members

const StatePost *state_post_ = {nullptr}
StatePost::const_iterator symbol_post_it_ = {}

Current symbol post iterator to iterate over. End symbol post iterator which is no longer iterated over (one after the last symbol post iterated over or end()).

StatePost::const_iterator symbol_post_end_ = {}
class const_iterator
#include <delta.hh>

Iterator over moves.

Public Types

using iterator_category = std::forward_iterator_tag
using value_type = Move
using difference_type = size_t
using pointer = Move*
using reference = Move&

Public Functions

inline const_iterator()

Construct end iterator.

const_iterator(const StatePost &state_post)

Const all moves iterator.

const_iterator(const StatePost &state_post, StatePost::const_iterator symbol_post_it, StatePost::const_iterator symbol_post_it_end)

Construct iterator from symbol_post_it (including) to symbol_post_it_end (excluding).

const_iterator(const const_iterator &other) noexcept = default
const_iterator(const_iterator&&) = default
inline const Move &operator*() const
inline const Move *operator->() const
const_iterator &operator++()
const const_iterator operator++(int)
const_iterator &operator=(const const_iterator &other) noexcept = default
const_iterator &operator=(const_iterator&&) = default
bool operator==(const const_iterator &other) const

Private Members

const StatePost *state_post_ = {nullptr}
StatePost::const_iterator symbol_post_it_ = {}
StateSet::const_iterator target_it_ = {}
StatePost::const_iterator symbol_post_end_ = {}
bool is_end_ = {false}
Move move_ = {}

Internal allocated instance of Move which is set for the move currently iterated over and returned as a reference with operator*().

class SymbolPost
#include <delta.hh>

Structure represents a post of a single symbol: a set of target states in transitions.

A set of SymbolPost, called StatePost, is describing the automata transitions from a single source state.

Public Functions

SymbolPost() = default
inline explicit SymbolPost(Symbol symbol)
inline SymbolPost(Symbol symbol, State state_to)
inline SymbolPost(Symbol symbol, StateSet states_to)
inline SymbolPost(SymbolPost &&rhs) noexcept
SymbolPost(const SymbolPost &rhs) = default
SymbolPost &operator=(SymbolPost &&rhs) noexcept
SymbolPost &operator=(const SymbolPost &rhs) = default
inline std::weak_ordering operator<=>(const SymbolPost &other) const
inline bool operator==(const SymbolPost &other) const
inline StateSet::iterator begin()
inline StateSet::iterator end()
inline StateSet::const_iterator cbegin() const
inline StateSet::const_iterator cend() const
inline size_t count(State s) const
inline bool empty() const
inline size_t num_of_targets() const
void insert(State s)
void insert(const StateSet &states)
inline void push_back(const State s)
template<typename ...Args>
inline StateSet &emplace_back(Args&&... args)
inline void erase(State s)
inline std::vector<State>::const_iterator find(State s) const
inline std::vector<State>::iterator find(State s)

Public Members

Symbol symbol = {}
StateSet targets = {}
class SynchronizedExistentialSymbolPostIterator : public SynchronizedExistentialIterator<utils::OrdVector<SymbolPost>::const_iterator>
#include <delta.hh>

Specialization of utils::SynchronizedExistentialIterator for iterating over SymbolPosts.

Public Functions

StateSet unify_targets() const

Get union of all targets.

bool synchronize_with(const SymbolPost &sync)

Synchronize with the given SymbolPost sync.

Alignes the synchronized iterator to the same symbol as sync.

Returns:

True iff the synchronized iterator points to the same symbol as sync.

bool synchronize_with(Symbol sync_symbol)

Synchronize with the given symbol sync_symbol.

Alignes the synchronized iterator to the same symbol as sync_symbol.

Returns:

True iff the synchronized iterator points to the same symbol as sync.

struct Transition
#include <delta.hh>

A single transition in Delta represented as a triple(source, symbol, target).

Public Functions

inline Transition()
Transition(const Transition&) = default
Transition(Transition&&) = default
Transition &operator=(const Transition&) = default
Transition &operator=(Transition&&) = default
inline Transition(const State source, const Symbol symbol, const State target)
auto operator<=>(const Transition&) const = default

Public Members

State source

Source state.

Symbol symbol

Transition symbol.

State target

Target state.

Builder

Function to build predefined types of NFAs, to create them from regular expressions, and to load them from files.

namespace mata

Main namespace including structs and algorithms for all automata.

In particular, this includes:

  1. Alphabets,

  2. Formula graphs and nodes,

  3. Mintermization,

  4. Closed sets.

namespace nfa

Nondeterministic Finite Automata including structures, transitions and algorithms.

In particular, this includes:

  1. Structures (Automaton, Transitions, Results, Delta),

  2. Algorithms (operations, checks, tests),

  3. Constructions.

Other algorithms are included in mata::nfa::Plumbing (simplified API for, e.g., binding) and mata::nfa::algorithms (concrete implementations of algorithms, such as for complement).

namespace builder

Namespace providing options to build NFAs.

Typedefs

using NameStateMap = std::unordered_map<std::string, State>

Functions

Nfa create_single_word_nfa(const std::vector<Symbol> &word)

Create an automaton accepting only a single word.

Nfa create_single_word_nfa(const std::vector<std::string> &word, Alphabet *alphabet = nullptr)

Create an automaton accepting only a single word.

Parameters:
  • wordWord to accept.

  • alphabetAlphabet to use in NFA for translating word into symbols. If specified, the alphabet has to contain translations for all of the word symbols. If left empty, a new alphabet with only the symbols of the word will be created.

Nfa create_empty_string_nfa()

Create automaton accepting only epsilon string.

Nfa create_sigma_star_nfa(Alphabet *alphabet = new OnTheFlyAlphabet{})

Create automaton accepting sigma star over the passed alphabet.

Parameters:

alphabet[in] Alphabet to construct sigma star automaton with. When alphabet is left empty, the default empty alphabet is used, creating an automaton accepting only the empty string.

Nfa create_random_nfa_tabakov_vardi(const size_t num_of_states, const size_t alphabet_size, const double states_trans_ratio_per_symbol, const double final_state_density)

Creates Tabakov-Vardi random NFA.

The implementation is based on the paper “Experimental Evaluation of Classical Automata Constructions” by Tabakov and Vardi.

Parameters:
  • num_of_states – Number of states in the automaton.

  • alphabet_size – Size of the alphabet.

  • states_transitions_ratio_per_symbol – Ratio between number of transitions and number of states for each symbol. The value must be in range [0, num_of_states]. A value of 1 means that there will be num_of_states transitions for each symbol. A value of num_of_states means that there will be a transition between every pair of states for each symbol.

  • final_state_density – Density of final states in the automaton. The value must be in range [0, 1]. The state 0 is always final. If the density is 1, every state will be final.

Nfa construct(const mata::parser::ParsedSection &parsec, Alphabet *alphabet, NameStateMap *state_map = nullptr)

Loads an automaton from Parsed object.

Nfa construct(const mata::IntermediateAut &inter_aut, Alphabet *alphabet, NameStateMap *state_map = nullptr)

Loads an automaton from Parsed object.

void construct(Nfa *result, const mata::IntermediateAut &inter_aut, Alphabet *alphabet, NameStateMap *state_map = nullptr)

Loads an automaton from Parsed object; version for python binding.

template<class ParsedObject>
Nfa construct(const ParsedObject &parsed, Alphabet *alphabet = nullptr, NameStateMap *state_map = nullptr)
Nfa parse_from_mata(std::istream &nfa_stream)

Parse NFA from the mata format in an input stream.

Parameters:

nfa_stream – Input stream containing NFA in mata format.

Throws:

std::runtime_error – Parsing of NFA fails.

Nfa parse_from_mata(const std::string &nfa_in_mata)

Parse NFA from the mata format in a string.

Parameters:

nfa_stream – String containing NFA in mata format.

Throws:

std::runtime_error – Parsing of NFA fails.

Nfa parse_from_mata(const std::filesystem::path &nfa_file)

Parse NFA from the mata format in a file.

Parameters:

nfa_stream – Path to the file containing NFA in mata format.

Throws:
  • std::runtime_errornfa_file does not exist.

  • std::runtime_error – Parsing of NFA fails.

Nfa create_from_regex(const std::string &regex)

Create NFA from regex.

The created NFA does not contain epsilons, is trimmed and reduced. It uses the parser from RE2, see mata/parser/re2praser.hh for more details and options.

At https://github.com/google/re2/wiki/Syntax, you can find the syntax of regex with following futher limitations: 1) The allowed characters are the first 256 characters of Unicode, i.e., Latin1 encoding (ASCII + 128 more characters). For the full Unicode, check mata/parser/re2praser.hh. 2) The created automaton represents the language of the regex and is not expected to be used in regex matching. Therefore, stuff like ^, $, , etc. are ignored in the regex. For example, the regular expressions a*b and ^a*b will both result in the same NFA accepting the language of multiple ‘a’s followed by one ‘b’. See also issue #494.

See also

mata::parser::create_nfa()

Parameters:

regex – regular expression

Algorithms

Functions to operate on NFAs, such as inclusion and universality checking, minimization, etc.

namespace mata

Main namespace including structs and algorithms for all automata.

In particular, this includes:

  1. Alphabets,

  2. Formula graphs and nodes,

  3. Mintermization,

  4. Closed sets.

namespace nfa

Nondeterministic Finite Automata including structures, transitions and algorithms.

In particular, this includes:

  1. Structures (Automaton, Transitions, Results, Delta),

  2. Algorithms (operations, checks, tests),

  3. Constructions.

Other algorithms are included in mata::nfa::Plumbing (simplified API for, e.g., binding) and mata::nfa::algorithms (concrete implementations of algorithms, such as for complement).

namespace algorithms

Concrete NFA implementations of algorithms, such as complement, inclusion, or universality checking.

This is a separation of the implementation from the interface defined in mata::nfa. Note, that in mata::nfa interface, there are particular dispatch functions calling these function according to parameters provided by a user. E.g. we can call the following function: is_universal(aut, alph, {{‘algorithm, ‘antichains’}})` to check for universality based on antichain-based algorithm.

Functions

Nfa minimize_brzozowski(const Nfa &aut)

Brzozowski minimization of automata (revert -> determinize -> revert -> determinize).

Parameters:

aut[in] Automaton to be minimized.

Returns:

Minimized automaton.

Nfa minimize_hopcroft(const Nfa &dfa_trimmed)

Hopcroft minimization of automata.

Based on the algorithm from the paper: “Efficient Minimization of DFAs With Partial Transition Functions” by Antti Valmari and Petri Lehtinen. The algorithm works in O(a*n*log(n)) time and O(m+n+a) space, where: n is the number of states, a is the size of the alphabet, and m is the number of transitions. [https://dl.acm.org/doi/10.1016/j.ipl.2011.12.004]

Parameters:

dfa_trimmed[in] Deterministic automaton without useless states. Perform trimming before calling this function.

Returns:

Minimized deterministic automaton.

Nfa complement_classical(const Nfa &aut, const mata::utils::OrdVector<Symbol> &symbols)

Complement implemented by determization, adding sink state and making automaton complete.

Then it adds final states which were non final in the original automaton.

Parameters:
  • aut[in] Automaton to be complemented.

  • symbols[in] Symbols needed to make the automaton complete.

  • minimize_during_determinization[in] Whether the determinized automaton is computed by (brzozowski) minimization.

Returns:

Complemented automaton.

Nfa complement_brzozowski(const Nfa &aut, const mata::utils::OrdVector<Symbol> &symbols)

Complement implemented by determization using Brzozowski minimization, adding a sink state and making the automaton complete.

Then it swaps final and non-final states.

Parameters:
  • aut[in] Automaton to be complemented.

  • symbols[in] Symbols needed to make the automaton complete.

Returns:

Complemented automaton.

bool is_included_naive(const Nfa &smaller, const Nfa &bigger, const Alphabet *alphabet = nullptr, Run *cex = nullptr)

Inclusion implemented by complementation of bigger automaton, intersecting it with smaller and then it checks emptiness of intersection.

Parameters:
  • smaller[in] Automaton which language should be included in the bigger one.

  • bigger[in] Automaton which language should include the smaller one.

  • alphabet[in] Alphabet of both automata (it is computed automatically, but it is more efficient to set it if you have it).

  • cex[out] A potential counterexample word which breaks inclusion

Returns:

True if smaller language is included, i.e., if the final intersection of smaller complement of bigger is empty.

bool is_included_antichains(const Nfa &smaller, const Nfa &bigger, const Alphabet *alphabet = nullptr, Run *cex = nullptr)

Inclusion implemented by antichain algorithms.

Parameters:
  • smaller[in] Automaton which language should be included in the bigger one

  • bigger[in] Automaton which language should include the smaller one

  • alphabet[in] Alphabet of both automata (not needed for antichain algorithm)

  • cex[out] A potential counterexample word which breaks inclusion

Returns:

True if smaller language is included, i.e., if the final intersection of smaller complement of bigger is empty.

bool is_universal_naive(const Nfa &aut, const Alphabet &alphabet, Run *cex)

Universality check implemented by checking emptiness of complemented automaton.

Parameters:
  • aut[in] Automaton which universality is checked

  • alphabet[in] Alphabet of the automaton

  • cex[out] Counterexample word which eventually breaks the universality

Returns:

True if the complemented automaton has non empty language, i.e., the original one is not universal

bool is_universal_antichains(const Nfa &aut, const Alphabet &alphabet, Run *cex)

Universality checking based on subset construction with antichain.

Parameters:
  • aut[in] Automaton which universality is checked

  • alphabet[in] Alphabet of the automaton

  • cex[out] Counterexample word which eventually breaks the universality

Returns:

True if the automaton is universal, otherwise false.

Simlib::Util::BinaryRelation compute_relation(const Nfa &aut, const ParameterMap &params = {{"relation", "simulation"}, {"direction", "forward"}})
Nfa product(const Nfa &lhs, const Nfa &rhs, const std::function<bool(State, State)> &&final_condition, const Symbol first_epsilon = EPSILON, std::unordered_map<std::pair<State, State>, State> *prod_map = nullptr)

Compute product of two NFAs, final condition is to be specified, with a possibility of using multiple epsilons.

Parameters:
  • lhs[in] First NFA to compute intersection for.

  • rhs[in] Second NFA to compute intersection for.

  • first_epsilons[in] The smallest epsilon.

  • final_condition[in] The predicate that tells whether a pair of states is final (conjunction for intersection).

  • prod_map[out] Can be used to get the mapping of the pairs of the original states to product states. Mostly useless, it is only filled in and returned if !=nullptr, but the algorithm internally uses another data structures, because this one is too slow.

Returns:

NFA as a product of NFAs lhs and rhs with ε-transitions preserved.

Nfa concatenate_eps(const Nfa &lhs, const Nfa &rhs, const Symbol &epsilon, bool use_epsilon = false, StateRenaming *lhs_state_renaming = nullptr, StateRenaming *rhs_state_renaming = nullptr)

Concatenate two NFAs.

Supports epsilon symbols when use_epsilon is set to true.

Parameters:
  • lhs[in] First automaton to concatenate.

  • rhs[in] Second automaton to concatenate.

  • epsilon[in] Epsilon to be used for concatenation (provided use_epsilon is true)

  • use_epsilon[in] Whether to concatenate over epsilon symbol.

  • lhs_state_renaming[out] Map mapping lhs states to result states.

  • rhs_state_renaming[out] Map mapping rhs states to result states.

Returns:

Concatenated automaton.

Nfa reduce_simulation(const Nfa &nfa, StateRenaming &state_renaming)

Reduce NFA using (forward) simulation.

Parameters:
  • nfa[in] NFA to reduce

  • state_renaming[out] Map mapping original states to the reduced states.

Nfa reduce_residual(const Nfa &nfa, StateRenaming &state_renaming, const std::string &type, const std::string &direction)

Reduce NFA using residual construction.

Parameters:
  • nfa[in] NFA to reduce.

  • state_renaming[out] Map mapping original states to the reduced states.

  • type[in] Type of the residual construction (values: “after”, “with”).

  • direction[in] Direction of the residual construction (values: “forward”, “backward”).

Nfa reduce_residual_with(const Nfa &nfa)

Reduce NFA using residual construction.

The residual construction of the residual automaton and the removal of the covering states is done during the last determinization.

Similar performance to reduce_residual_after(). The output is almost the same except the transitions: transitions may slightly differ, but the number of states is the same for both algorithm types.

Nfa reduce_residual_after(const Nfa &nfa)

Reduce NFA using residual construction.

The residual construction of the residual automaton and the removal of the covering states is done after the final determinization.

Similar performance to reduce_residual_with(). The output is almost the same except the transitions: transitions may slightly differ, but the number of states is the same for both algorithm types.

Strings

NFA algorithms and classes used for solving string constraints.

namespace mata

Main namespace including structs and algorithms for all automata.

In particular, this includes:

  1. Alphabets,

  2. Formula graphs and nodes,

  3. Mintermization,

  4. Closed sets.

namespace strings

NFA algorithms usable for solving string constraints.

Typedefs

using Nfa = nfa::Nfa
using State = nfa::State
using StateSet = nfa::StateSet
using Transition = nfa::Transition
using ParameterMap = nfa::ParameterMap
using SymbolPost = nfa::SymbolPost
using Nft = nft::Nft

Functions

std::set<Word> get_shortest_words(const Nfa &nfa)

Get shortest words (regarding their length) of the automaton using BFS.

Returns:

Set of shortest words.

std::set<Symbol> get_accepted_symbols(const Nfa &nfa)

Get all the one symbol words accepted by nfa.

std::set<std::pair<int, int>> get_word_lengths(const Nfa &aut)

Get the lengths of all words in the automaton aut.

The function returns a set of pairs <u,v> where for each such a pair there is a word with length u+k*v for all ks. The disjunction of such formulae of all pairs hence precisely describe lengths of all words in the automaton.

Parameters:

aut – Input automaton

Returns:

Set of pairs describing lengths

Nfa reluctant_nfa(Nfa nfa)

Modify nfa in-place to remove outgoing transitions from final states.

If nfa accepts empty string, returned NFA will accept only the empty string.

Parameters:

nfa – NFA to modify.

Returns:

The reluctant version of nfa.

bool is_lang_eps(const Nfa &nfa)

Checks if the automaton nfa accepts only a single word \eps.

Parameters:

nfa – Input automaton

Returns:

true iff L(nfa) = {\eps}

class ShortestWordsMap
#include <strings.hh>

Class mapping states to the shortest words accepted by languages of the states.

Public Functions

inline explicit ShortestWordsMap(const Nfa &aut)

Maps states in the automaton aut to shortest words accepted by languages of the states.

Parameters:

aut – Automaton to compute shortest words for.

std::set<Word> get_shortest_words_from(const StateSet &states) const

Gets shortest words for the given states.

Parameters:

states[in] States to map shortest words for.

Returns:

Set of shortest words.

std::set<Word> get_shortest_words_from(State state) const

Gets shortest words for the given state.

Parameters:

state[in] State to map shortest words for.

Returns:

Set of shortest words.

Private Types

using WordLength = int

A length of a word. Pair binding the length of all words in the word set and word set with words of the given length.

using LengthWordsPair = std::pair<WordLength, std::set<Word>>

Private Functions

void insert_initial_lengths()

Inserts initial lengths into the shortest words map.

Inserts initial length of length 0 for final state in the automaton (initial states in the reversed automaton).

void compute()

Computes shortest words for all states in the automaton.

void compute_for_state(State state)

Computes shortest words for the given state.

Parameters:

state[in] State to compute shortest words for.

inline LengthWordsPair map_default_shortest_words(const State state)

Creates default shortest words mapping for yet unprocessed state.

Parameters:

state[in] State to map default shortest words.

Returns:

Created default shortest words map element for the given state.

Private Members

std::unordered_map<State, LengthWordsPair> shortest_words_map = {}

Map mapping states to the shortest words accepted by the automaton from the mapped state.

std::set<State> processed = {}

Set of already processed states.

std::deque<State> fifo_queue = {}

FIFO queue for states to process.

const Nfa reversed_automaton

Reversed input automaton.

Private Static Functions

static void update_current_words(LengthWordsPair &act, const LengthWordsPair &dst, Symbol symbol)

Update words for the current state.

Parameters:
  • act[out] Current state shortest words and length.

  • dst[in] Transition target state shortest words and length.

  • symbol[in] Symbol to update with.

namespace seg_nfa

Segment Automata including structs and algorithms.

These are automata whose state space can be split into several segments connected by ε-transitions in a chain. No other ε-transitions are allowed. As a consequence, no ε-transitions can appear in a cycle. Segment automaton can have initial states only in the first segment and final states only in the last segment.

Typedefs

using SegNfa = Nfa
using VisitedEpsMap = std::map<State, std::map<Symbol, unsigned>>
using VisitedEpsilonsCounterMap = std::map<Symbol, unsigned>

Number of visited epsilons.

using VisitedEpsilonsCounterVector = std::vector<unsigned>

Projection of VisitedEpsilonsNumberMap to sorted keys (in descending order).

using Noodle = std::vector<std::shared_ptr<SegNfa>>

A noodle is represented as a sequence of segments (a copy of the segment automata) created as if there was exactly one ε-transition between each two consecutive segments.

using SegmentWithEpsilonsCounter = std::pair<std::shared_ptr<Nfa>, VisitedEpsilonsCounterVector>

Segment with a counter of visited epsilons.

using NoodleWithEpsilonsCounter = std::vector<SegmentWithEpsilonsCounter>

Noodles as segments enriched with EpsCntMap.

using TransducerNoodle = std::vector<TransducerNoodleElement>

Functions

void segs_one_initial_final(const std::vector<Nfa> &segments, bool include_empty, const State &unused_state, std::map<std::pair<State, State>, std::shared_ptr<Nfa>> &out)

segs_one_initial_final

segments_one_initial_final[init, final] is the pointer to automaton created from one of the segments such that init and final are one of the initial and final states of the segment and the created automaton takes this segment, sets initial={init}, final={final} and trims it; also segments_one_initial_final[unused_state, final] is used for the first segment (where we always want all initial states, only final state changes) and segments_one_initial_final[init, unused_state] is similarly for the last segment TODO: should we use unordered_map? then we need hash

std::vector<Noodle> noodlify(const SegNfa &aut, Symbol epsilon, bool include_empty = false)

Create noodles from segment automaton aut.

Segment automaton is a chain of finite automata (segments) connected via ε-transitions. A noodle is a vector of pointers to copy of the segments automata created as if there was exactly one ε-transition between each two consecutive segments.

Parameters:
  • automaton[in] Segment automaton to noodlify.

  • epsilon[in] Epsilon symbol to noodlify for.

  • include_empty[in] Whether to also include empty noodles.

Returns:

A list of all (non-empty) noodles.

std::vector<NoodleWithEpsilonsCounter> noodlify_mult_eps(const SegNfa &aut, const std::set<Symbol> &epsilons, bool include_empty = false)

Create noodles from segment automaton aut.

Segment automaton is a chain of finite automata (segments) connected via ε-transitions. A noodle is a vector of pointers to copy of the segments automata created as if there was exactly one ε-transition between each two consecutive segments.

Parameters:
  • automaton[in] Segment automaton to noodlify.

  • epsilons[in] Epsilon symbols to noodlify for.

  • include_empty[in] Whether to also include empty noodles.

Returns:

A list of all (non-empty) noodles.

std::vector<Noodle> noodlify_for_equation(const std::vector<std::reference_wrapper<Nfa>> &lhs_automata, const Nfa &rhs_automaton, bool include_empty = false, const ParameterMap &params = {{"reduce", "false"}})

Create noodles for left and right side of equation.

Segment automaton is a chain of finite automata (segments) connected via ε-transitions. A noodle is a copy of the segment automaton with exactly one ε-transition between each two consecutive segments.

Mata cannot work with equations, queries etc. Hence, we compute the noodles for the equation, but represent the equation in a way that libMata understands. The left side automata represent the left side of the equation and the right automaton represents the right side of the equation. To create noodles, we need a segment automaton representing the intersection. That can be achieved by computing a product of both sides. First, the left side has to be concatenated over an epsilon transitions into a single automaton to compute the intersection on, though.

Parameters:
  • lhs_automata[in] Sequence of segment automata for left side of an equation to noodlify.

  • rhs_automaton[in] Segment automaton for right side of an equation to noodlify.

  • include_empty[in] Whether to also include empty noodles.

  • params[in] Additional parameters for the noodlification:

    • ”reduce”: “false”, “forward”, “backward”, “bidirectional”; Execute forward, backward or bidirectional simulation minimization before noodlification.

Returns:

A list of all (non-empty) noodles.

std::vector<Noodle> noodlify_for_equation(const std::vector<Nfa*> &lhs_automata, const Nfa &rhs_automaton, bool include_empty = false, const ParameterMap &params = {{"reduce", "false"}})

Create noodles for left and right side of equation.

Segment automaton is a chain of finite automata (segments) connected via ε-transitions. A noodle is a copy of the segment automaton with exactly one ε-transition between each two consecutive segments.

Mata cannot work with equations, queries etc. Hence, we compute the noodles for the equation, but represent the equation in a way that libMata understands. The left side automata represent the left side of the equation and the right automaton represents the right side of the equation. To create noodles, we need a segment automaton representing the intersection. That can be achieved by computing a product of both sides. First, the left side has to be concatenated over an epsilon transitions into a single automaton to compute the intersection on, though.

Parameters:
  • lhs_automata[in] Sequence of pointers to segment automata for left side of an equation to noodlify.

  • rhs_automaton[in] Segment automaton for right side of an equation to noodlify.

  • include_empty[in] Whether to also include empty noodles.

  • params[in] Additional parameters for the noodlification:

    • ”reduce”: “false”, “forward”, “backward”, “bidirectional”; Execute forward, backward or bidirectional simulation minimization before noodlification.

Returns:

A list of all (non-empty) noodles.

std::vector<NoodleWithEpsilonsCounter> noodlify_for_equation(const std::vector<std::shared_ptr<Nfa>> &lhs_automata, const std::vector<std::shared_ptr<Nfa>> &rhs_automata, bool include_empty = false, const ParameterMap &params = {{"reduce", "false"}})

Create noodles for left and right side of equation (both sides are given as a sequence of automata).

Parameters:
  • lhs_automata[in] Sequence of pointers to segment automata for left side of an equation to noodlify.

  • rhs_automaton[in] Sequence of pointers to segment automata for right side of an equation to noodlify.

  • include_empty[in] Whether to also include empty noodles.

  • params[in] Additional parameters for the noodlification:

    • ”reduce”: “false”, “forward”, “backward”, “bidirectional”; Execute forward, backward or bidirectional simulation minimization before noodlification.

Returns:

A list of all (non-empty) noodles together with the positions reached from the beginning of left/right side.

std::vector<TransducerNoodle> noodlify_for_transducer(std::shared_ptr<Nft> nft, const std::vector<std::shared_ptr<Nfa>> &input_automata, const std::vector<std::shared_ptr<Nfa>> &output_automata, bool reduce_intersection = false)
VisitedEpsilonsCounterVector process_eps_map(const VisitedEpsilonsCounterMap &eps_cnt)

Process epsilon map to a sequence of values (sorted according to key desc)

Parameters:

eps_cnt – Epsilon count

Returns:

Vector of keys (count of epsilons)

class Segmentation
#include <strings.hh>

Class executing segmentation operations for a given segment automaton.

Works only with segment automata.

Public Types

using EpsilonDepth = size_t

Depth of ε-transitions. Dictionary of lists of ε-transitions grouped by their depth. For each depth ‘i’ we have ‘depths[i]’ which contains a list of ε-transitions of depth ‘i’.

using EpsilonDepthTransitions = std::unordered_map<EpsilonDepth, std::vector<Transition>>
using EpsilonDepthTransitionMap = std::unordered_map<EpsilonDepth, std::unordered_map<State, std::vector<Transition>>>

Public Functions

inline explicit Segmentation(const SegNfa &aut, const std::set<Symbol> &epsilons)

Prepare automaton aut for segmentation.

Parameters:
  • aut[in] Segment automaton to make segments for.

  • epsilon[in] Symbol to execute segmentation for.

inline const EpsilonDepthTransitions &get_epsilon_depths() const

Get segmentation depths for ε-transitions.

Returns:

Map of depths to lists of ε-transitions.

inline const EpsilonDepthTransitionMap &get_epsilon_depth_trans_map() const

Get the epsilon depth trans map object (mapping of depths and states to eps-successors)

Returns:

Map of depths to a map of states to transitions

const std::vector<Nfa> &get_segments()

Get segment automata.

Returns:

A vector of segments for the segment automaton in the order from the left (initial state in segment automaton) to the right (final states of segment automaton).

const std::vector<Nfa> &get_untrimmed_segments()

Get raw segment automata.

Returns:

A vector of segments for the segment automaton in the order from the left (initial state in segment automaton) to the right (final states of segment automaton) without trimming (the states are same as in the original automaton).

inline const VisitedEpsMap &get_visited_eps() const

Private Functions

void compute_epsilon_depths()

Compute epsilon depths with their transitions.

void split_aut_into_segments()

Split segment automaton into segments.

void update_next_segment(size_t current_depth, const Transition &transition)

Propagate changes to the current segment automaton to the next segment automaton.

Parameters:
  • current_depth[in] Current depth.

  • transition[in] Current epsilon transition.

void update_current_segment(size_t current_depth, const Transition &transition)

Update current segment automaton.

Parameters:
  • current_depth[in] Current depth.

  • transition[in] Current epsilon transition.

std::unordered_map<State, bool> initialize_visited_map() const

Initialize map of visited states.

Returns:

Map of visited states.

std::deque<StateDepthTuple> initialize_worklist() const

Initialize worklist of states with depths to process.

Returns:

Queue of state and its depth pairs.

void process_state_depth_pair(const StateDepthTuple &state_depth_pair, std::deque<StateDepthTuple> &worklist)

Process pair of state and its depth.

Parameters:
  • state_depth_pair[in] Current state depth pair.

  • worklist[out] Worklist of state and depth pairs to process.

void add_transitions_to_worklist(const StateDepthTuple &state_depth_pair, const SymbolPost &move, std::deque<StateDepthTuple> &worklist)

Add states with non-epsilon transitions to the worklist.

Parameters:
  • move[in] – Move from current state.

  • depth[in] – Current depth.

  • worklist[out] – Worklist of state and depth pairs to process.

void handle_epsilon_transitions(const StateDepthTuple &state_depth_pair, const SymbolPost &move, std::deque<StateDepthTuple> &worklist)

Process epsilon transitions for the current state.

Parameters:
  • state_depth_pair[in] Current state depth pair.

  • move[in] Move from current state.

  • worklist[out] Worklist of state and depth pairs to process.

void remove_inner_initial_and_final_states()

Remove inner initial and final states.

Remove all initial states for all segments but the first one and all final states for all segments but the last one.

Private Members

const std::set<Symbol> epsilons

Symbol for which to execute segmentation. Automaton to execute segmentation for. Must be a segment automaton (can be split into segments).

const SegNfa &automaton
EpsilonDepthTransitions epsilon_depth_transitions = {}

Epsilon depths.

EpsilonDepthTransitionMap eps_depth_trans_map = {}
std::vector<SegNfa> segments = {}

Epsilon depths with mapping of states to epsilon transitions.

Segments for automaton.

std::vector<SegNfa> segments_raw = {}

Raw segments for automaton.

VisitedEpsMap visited_eps = {}
struct StateDepthTuple

number of visited eps for each state

Pair of state and its depth.

Public Members

State state

State with a depth.

EpsilonDepth depth

Depth of a state.

VisitedEpsilonsCounterMap eps
struct TransducerNoodleElement
#include <strings.hh>

Public Functions

inline TransducerNoodleElement(std::shared_ptr<Nft> transducer, std::shared_ptr<Nfa> input_aut, unsigned input_index, std::shared_ptr<Nfa> output_aut, unsigned output_index)

Public Members

std::shared_ptr<Nft> transducer
std::shared_ptr<Nfa> input_aut
unsigned input_index
std::shared_ptr<Nfa> output_aut
unsigned output_index