Browse Source

Bindings for dd and hybrid model checking

refactoring
Matthias Volk 7 years ago
parent
commit
62f3d3630e
  1. 82
      lib/stormpy/__init__.py
  2. 22
      src/core/modelchecking.cpp
  3. 28
      src/core/result.cpp
  4. 21
      tests/core/test_modelchecking.py
  5. 28
      tests/pars/test_parametric.py

82
lib/stormpy/__init__.py

@ -93,6 +93,28 @@ def build_model(symbolic_description, properties=None):
""" """
Build a model in sparse representation from a symbolic description. Build a model in sparse representation from a symbolic description.
:param symbolic_description: Symbolic model description to translate into a model.
:param List[Property] properties: List of properties that should be preserved during the translation. If None, then all properties are preserved.
:return: Model in sparse representation.
"""
return build_sparse_model(symbolic_description, properties=properties)
def build_parametric_model(symbolic_description, properties=None):
"""
Build a parametric model in sparse representation from a symbolic description.
:param symbolic_description: Symbolic model description to translate into a model.
:param List[Property] properties: List of properties that should be preserved during the translation. If None, then all properties are preserved.
:return: Parametric model in sparse representation.
"""
return build_sparse_parametric_model(symbolic_description, properties=properties)
def build_sparse_model(symbolic_description, properties=None):
"""
Build a model in sparse representation from a symbolic description.
:param symbolic_description: Symbolic model description to translate into a model. :param symbolic_description: Symbolic model description to translate into a model.
:param List[Property] properties: List of properties that should be preserved during the translation. If None, then all properties are preserved. :param List[Property] properties: List of properties that should be preserved during the translation. If None, then all properties are preserved.
:return: Model in sparse representation. :return: Model in sparse representation.
@ -108,7 +130,7 @@ def build_model(symbolic_description, properties=None):
return _convert_sparse_model(intermediate, parametric=False) return _convert_sparse_model(intermediate, parametric=False)
def build_parametric_model(symbolic_description, properties=None):
def build_sparse_parametric_model(symbolic_description, properties=None):
""" """
Build a parametric model in sparse representation from a symbolic description. Build a parametric model in sparse representation from a symbolic description.
@ -203,6 +225,20 @@ def perform_bisimulation(model, properties, bisimulation_type):
def model_checking(model, property, only_initial_states=False, extract_scheduler=False): def model_checking(model, property, only_initial_states=False, extract_scheduler=False):
"""
Perform model checking on model for property.
:param model: Model.
:param property: Property to check for.
:param only_initial_states: If True, only results for initial states are computed, otherwise for all states.
:param extract_scheduler: If True, try to extract a scheduler
:return: Model checking result.
:rtype: CheckResult
"""
return check_model_sparse(model, property, only_initial_states=only_initial_states,
extract_scheduler=extract_scheduler)
def check_model_sparse(model, property, only_initial_states=False, extract_scheduler=False):
""" """
Perform model checking on model for property. Perform model checking on model for property.
:param model: Model. :param model: Model.
@ -227,6 +263,50 @@ def model_checking(model, property, only_initial_states=False, extract_scheduler
return core._model_checking_sparse_engine(model, task) return core._model_checking_sparse_engine(model, task)
def check_model_dd(model, property, only_initial_states=False):
"""
Perform model checking using dd engine.
:param model: Model.
:param property: Property to check for.
:param only_initial_states: If True, only results for initial states are computed, otherwise for all states.
:return: Model checking result.
:rtype: CheckResult
"""
if isinstance(property, Property):
formula = property.raw_formula
else:
formula = property
if model.supports_parameters:
task = core.ParametricCheckTask(formula, only_initial_states)
return core._parametric_model_checking_dd_engine(model, task)
else:
task = core.CheckTask(formula, only_initial_states)
return core._model_checking_dd_engine(model, task)
def check_model_hybrid(model, property, only_initial_states=False):
"""
Perform model checking using hybrid engine.
:param model: Model.
:param property: Property to check for.
:param only_initial_states: If True, only results for initial states are computed, otherwise for all states.
:return: Model checking result.
:rtype: CheckResult
"""
if isinstance(property, Property):
formula = property.raw_formula
else:
formula = property
if model.supports_parameters:
task = core.ParametricCheckTask(formula, only_initial_states)
return core._parametric_model_checking_hybrid_engine(model, task)
else:
task = core.CheckTask(formula, only_initial_states)
return core._model_checking_hybrid_engine(model, task)
def prob01min_states(model, eventually_formula): def prob01min_states(model, eventually_formula):
assert type(eventually_formula) == logic.EventuallyFormula assert type(eventually_formula) == logic.EventuallyFormula
labelform = eventually_formula.subformula labelform = eventually_formula.subformula

22
src/core/modelchecking.cpp

@ -4,12 +4,24 @@
template<typename ValueType> template<typename ValueType>
using CheckTask = storm::modelchecker::CheckTask<storm::logic::Formula, ValueType>; using CheckTask = storm::modelchecker::CheckTask<storm::logic::Formula, ValueType>;
// Thin wrapper for model checking
// Thin wrapper for model checking using sparse engine
template<typename ValueType> template<typename ValueType>
std::shared_ptr<storm::modelchecker::CheckResult> modelCheckingSparseEngine(std::shared_ptr<storm::models::sparse::Model<ValueType>> model, CheckTask<ValueType> const& task) { std::shared_ptr<storm::modelchecker::CheckResult> modelCheckingSparseEngine(std::shared_ptr<storm::models::sparse::Model<ValueType>> model, CheckTask<ValueType> const& task) {
return storm::api::verifyWithSparseEngine<ValueType>(model, task); return storm::api::verifyWithSparseEngine<ValueType>(model, task);
} }
// Thin wrapper for model checking using dd engine
template<storm::dd::DdType DdType, typename ValueType>
std::shared_ptr<storm::modelchecker::CheckResult> modelCheckingDdEngine(std::shared_ptr<storm::models::symbolic::Model<DdType, ValueType>> model, CheckTask<ValueType> const& task) {
return storm::api::verifyWithDdEngine<DdType, ValueType>(model, task);
}
// Thin wrapper for model checking using hybrid engine
template<storm::dd::DdType DdType, typename ValueType>
std::shared_ptr<storm::modelchecker::CheckResult> modelCheckingHybridEngine(std::shared_ptr<storm::models::symbolic::Model<DdType, ValueType>> model, CheckTask<ValueType> const& task) {
return storm::api::verifyWithHybridEngine<DdType, ValueType>(model, task);
}
// Thin wrapper for computing prob01 states // Thin wrapper for computing prob01 states
template<typename ValueType> template<typename ValueType>
std::pair<storm::storage::BitVector, storm::storage::BitVector> computeProb01(storm::models::sparse::Dtmc<ValueType> const& model, storm::storage::BitVector const& phiStates, storm::storage::BitVector const& psiStates) { std::pair<storm::storage::BitVector, storm::storage::BitVector> computeProb01(storm::models::sparse::Dtmc<ValueType> const& model, storm::storage::BitVector const& phiStates, storm::storage::BitVector const& psiStates) {
@ -42,8 +54,12 @@ void define_modelchecking(py::module& m) {
; ;
// Model checking // Model checking
m.def("_model_checking_sparse_engine", &modelCheckingSparseEngine<double>, "Perform model checking", py::arg("model"), py::arg("task"));
m.def("_parametric_model_checking_sparse_engine", &modelCheckingSparseEngine<storm::RationalFunction>, "Perform parametric model checking", py::arg("model"), py::arg("task"));
m.def("_model_checking_sparse_engine", &modelCheckingSparseEngine<double>, "Perform model checking using the sparse engine", py::arg("model"), py::arg("task"));
m.def("_parametric_model_checking_sparse_engine", &modelCheckingSparseEngine<storm::RationalFunction>, "Perform parametric model checking using the sparse engine", py::arg("model"), py::arg("task"));
m.def("_model_checking_dd_engine", &modelCheckingDdEngine<storm::dd::DdType::Sylvan, double>, "Perform model checking using the dd engine", py::arg("model"), py::arg("task"));
m.def("_parametric_model_checking_dd_engine", &modelCheckingDdEngine<storm::dd::DdType::Sylvan, storm::RationalFunction>, "Perform parametric model checking using the dd engine", py::arg("model"), py::arg("task"));
m.def("_model_checking_hybrid_engine", &modelCheckingHybridEngine<storm::dd::DdType::Sylvan, double>, "Perform model checking using the hybrid engine", py::arg("model"), py::arg("task"));
m.def("_parametric_model_checking_hybrid_engine", &modelCheckingHybridEngine<storm::dd::DdType::Sylvan, storm::RationalFunction>, "Perform parametric model checking using the hybrid engine", py::arg("model"), py::arg("task"));
m.def("_compute_prob01states_double", &computeProb01<double>, "Compute prob-0-1 states", py::arg("model"), py::arg("phi_states"), py::arg("psi_states")); m.def("_compute_prob01states_double", &computeProb01<double>, "Compute prob-0-1 states", py::arg("model"), py::arg("phi_states"), py::arg("psi_states"));
m.def("_compute_prob01states_rationalfunc", &computeProb01<storm::RationalFunction>, "Compute prob-0-1 states", py::arg("model"), py::arg("phi_states"), py::arg("psi_states")); m.def("_compute_prob01states_rationalfunc", &computeProb01<storm::RationalFunction>, "Compute prob-0-1 states", py::arg("model"), py::arg("phi_states"), py::arg("psi_states"));
m.def("_compute_prob01states_min_double", &computeProb01min<double>, "Compute prob-0-1 states (min)", py::arg("model"), py::arg("phi_states"), py::arg("psi_states")); m.def("_compute_prob01states_min_double", &computeProb01min<double>, "Compute prob-0-1 states (min)", py::arg("model"), py::arg("phi_states"), py::arg("psi_states"));

28
src/core/result.cpp

@ -1,11 +1,9 @@
#include "result.h" #include "result.h"
#include "storm/analysis/GraphConditions.h" #include "storm/analysis/GraphConditions.h"
// Thin wrapper
template<typename ValueType>
std::vector<ValueType> getValues(storm::modelchecker::ExplicitQuantitativeCheckResult<ValueType> const& result) {
return result.getValueVector();
}
#include "storm/modelchecker/results/SymbolicQualitativeCheckResult.h"
#include "storm/modelchecker/results/SymbolicQuantitativeCheckResult.h"
#include "storm/modelchecker/results/HybridQuantitativeCheckResult.h"
#include "storm/models/symbolic/StandardRewardModel.h"
// Define python bindings // Define python bindings
void define_result(py::module& m) { void define_result(py::module& m) {
@ -49,6 +47,8 @@ void define_result(py::module& m) {
}, py::arg("state"), "Get result for given state") }, py::arg("state"), "Get result for given state")
.def("get_truth_values", &storm::modelchecker::ExplicitQualitativeCheckResult::getTruthValuesVector, "Get BitVector representing the truth values") .def("get_truth_values", &storm::modelchecker::ExplicitQualitativeCheckResult::getTruthValuesVector, "Get BitVector representing the truth values")
; ;
py::class_<storm::modelchecker::SymbolicQualitativeCheckResult<storm::dd::DdType::Sylvan>, std::shared_ptr<storm::modelchecker::SymbolicQualitativeCheckResult<storm::dd::DdType::Sylvan>>>(m, "SymbolicQualitativeCheckResult", "Symbolic qualitative model checking result", qualitativeCheckResult)
;
// QuantitativeCheckResult // QuantitativeCheckResult
py::class_<storm::modelchecker::QuantitativeCheckResult<double>, std::shared_ptr<storm::modelchecker::QuantitativeCheckResult<double>>> quantitativeCheckResult(m, "_QuantitativeCheckResult", "Abstract class for quantitative model checking results", checkResult); py::class_<storm::modelchecker::QuantitativeCheckResult<double>, std::shared_ptr<storm::modelchecker::QuantitativeCheckResult<double>>> quantitativeCheckResult(m, "_QuantitativeCheckResult", "Abstract class for quantitative model checking results", checkResult);
@ -56,15 +56,27 @@ void define_result(py::module& m) {
.def("at", [](storm::modelchecker::ExplicitQuantitativeCheckResult<double> const& result, storm::storage::sparse::state_type state) { .def("at", [](storm::modelchecker::ExplicitQuantitativeCheckResult<double> const& result, storm::storage::sparse::state_type state) {
return result[state]; return result[state];
}, py::arg("state"), "Get result for given state") }, py::arg("state"), "Get result for given state")
.def("get_values", &getValues<double>, "Get model checking result values for all states")
.def("get_values", [](storm::modelchecker::ExplicitQuantitativeCheckResult<double> const& res) {return res.getValueVector();}, "Get model checking result values for all states")
.def_property_readonly("scheduler", [](storm::modelchecker::ExplicitQuantitativeCheckResult<double> const& res) {return res.getScheduler();}, "get scheduler") .def_property_readonly("scheduler", [](storm::modelchecker::ExplicitQuantitativeCheckResult<double> const& res) {return res.getScheduler();}, "get scheduler")
; ;
py::class_<storm::modelchecker::SymbolicQuantitativeCheckResult<storm::dd::DdType::Sylvan, double>, std::shared_ptr<storm::modelchecker::SymbolicQuantitativeCheckResult<storm::dd::DdType::Sylvan, double>>>(m, "SymbolicQuantitativeCheckResult", "Symbolic quantitative model checking result", quantitativeCheckResult)
;
py::class_<storm::modelchecker::HybridQuantitativeCheckResult<storm::dd::DdType::Sylvan, double>, std::shared_ptr<storm::modelchecker::HybridQuantitativeCheckResult<storm::dd::DdType::Sylvan, double>>>(m, "HybridQuantitativeCheckResult", "Hybrid quantitative model checking result", quantitativeCheckResult)
.def("get_values", &storm::modelchecker::HybridQuantitativeCheckResult<storm::dd::DdType::Sylvan, double>::getExplicitValueVector, "Get model checking result values for all states")
;
py::class_<storm::modelchecker::QuantitativeCheckResult<storm::RationalFunction>, std::shared_ptr<storm::modelchecker::QuantitativeCheckResult<storm::RationalFunction>>> parametricQuantitativeCheckResult(m, "_ParametricQuantitativeCheckResult", "Abstract class for parametric quantitative model checking results", checkResult); py::class_<storm::modelchecker::QuantitativeCheckResult<storm::RationalFunction>, std::shared_ptr<storm::modelchecker::QuantitativeCheckResult<storm::RationalFunction>>> parametricQuantitativeCheckResult(m, "_ParametricQuantitativeCheckResult", "Abstract class for parametric quantitative model checking results", checkResult);
py::class_<storm::modelchecker::ExplicitQuantitativeCheckResult<storm::RationalFunction>, std::shared_ptr<storm::modelchecker::ExplicitQuantitativeCheckResult<storm::RationalFunction>>>(m, "ExplicitParametricQuantitativeCheckResult", "Explicit parametric quantitative model checking result", parametricQuantitativeCheckResult) py::class_<storm::modelchecker::ExplicitQuantitativeCheckResult<storm::RationalFunction>, std::shared_ptr<storm::modelchecker::ExplicitQuantitativeCheckResult<storm::RationalFunction>>>(m, "ExplicitParametricQuantitativeCheckResult", "Explicit parametric quantitative model checking result", parametricQuantitativeCheckResult)
.def("at", [](storm::modelchecker::ExplicitQuantitativeCheckResult<storm::RationalFunction> const& result, storm::storage::sparse::state_type state) { .def("at", [](storm::modelchecker::ExplicitQuantitativeCheckResult<storm::RationalFunction> const& result, storm::storage::sparse::state_type state) {
return result[state]; return result[state];
}, py::arg("state"), "Get result for given state") }, py::arg("state"), "Get result for given state")
.def("get_values", &getValues<storm::RationalFunction>, "Get model checking result values for all states")
.def("get_values", [](storm::modelchecker::ExplicitQuantitativeCheckResult<storm::RationalFunction> const& res) { return res.getValueVector();}, "Get model checking result values for all states")
;
py::class_<storm::modelchecker::SymbolicQuantitativeCheckResult<storm::dd::DdType::Sylvan, storm::RationalFunction>, std::shared_ptr<storm::modelchecker::SymbolicQuantitativeCheckResult<storm::dd::DdType::Sylvan, storm::RationalFunction>>>(m, "SymbolicParametricQuantitativeCheckResult", "Symbolic parametric quantitative model checking result", quantitativeCheckResult)
; ;
py::class_<storm::modelchecker::HybridQuantitativeCheckResult<storm::dd::DdType::Sylvan, storm::RationalFunction>, std::shared_ptr<storm::modelchecker::HybridQuantitativeCheckResult<storm::dd::DdType::Sylvan, storm::RationalFunction>>>(m, "HybridParametricQuantitativeCheckResult", "Symbolic parametric hybrid quantitative model checking result", quantitativeCheckResult)
.def("get_values", &storm::modelchecker::HybridQuantitativeCheckResult<storm::dd::DdType::Sylvan, storm::RationalFunction>::getExplicitValueVector, "Get model checking result values for all states")
;
} }

21
tests/core/test_modelchecking.py

@ -127,3 +127,24 @@ class TestModelChecking:
assert initial_state == 1 assert initial_state == 1
result = stormpy.model_checking(model, formulas[0]) result = stormpy.model_checking(model, formulas[0])
assert math.isclose(result.at(initial_state), 4.166666667) assert math.isclose(result.at(initial_state), 4.166666667)
def test_model_checking_prism_dd_dtmc(self):
program = stormpy.parse_prism_program(get_example_path("dtmc", "die.pm"))
formulas = stormpy.parse_properties_for_prism_program("P=? [ F \"one\" ]", program)
model = stormpy.build_symbolic_model(program, formulas)
assert model.nr_states == 13
assert model.nr_transitions == 20
result = stormpy.check_model_dd(model, formulas[0])
assert type(result) is stormpy.SymbolicQuantitativeCheckResult
def test_model_checking_prism_hybrid_dtmc(self):
program = stormpy.parse_prism_program(get_example_path("dtmc", "die.pm"))
formulas = stormpy.parse_properties_for_prism_program("P=? [ F \"one\" ]", program)
model = stormpy.build_symbolic_model(program, formulas)
assert model.nr_states == 13
assert model.nr_transitions == 20
result = stormpy.check_model_hybrid(model, formulas[0])
assert type(result) is stormpy.HybridQuantitativeCheckResult
values = result.get_values()
assert len(values) == 3
assert math.isclose(values[0], 0.16666666666666663)

28
tests/pars/test_parametric.py

@ -8,7 +8,7 @@ from configurations import pars
@pars @pars
class TestParametric: class TestParametric:
def test_parametric_state_elimination(self):
def test_parametric_model_checking_sparse(self):
program = stormpy.parse_prism_program(get_example_path("pdtmc", "brp16_2.pm")) program = stormpy.parse_prism_program(get_example_path("pdtmc", "brp16_2.pm"))
prop = "P=? [F s=5]" prop = "P=? [F s=5]"
formulas = stormpy.parse_properties_for_prism_program(prop, program) formulas = stormpy.parse_properties_for_prism_program(prop, program)
@ -24,6 +24,32 @@ class TestParametric:
one = stormpy.FactorizedPolynomial(stormpy.RationalRF(1)) one = stormpy.FactorizedPolynomial(stormpy.RationalRF(1))
assert func.denominator == one assert func.denominator == one
def test_parametric_model_checking_dd(self):
program = stormpy.parse_prism_program(get_example_path("pdtmc", "parametric_die.pm"))
prop = "P=? [F s=5]"
formulas = stormpy.parse_properties_for_prism_program(prop, program)
model = stormpy.build_symbolic_parametric_model(program, formulas)
assert model.nr_states == 11
assert model.nr_transitions == 17
assert model.model_type == stormpy.ModelType.DTMC
assert model.has_parameters
result = stormpy.check_model_dd(model, formulas[0])
assert type(result) is stormpy.SymbolicParametricQuantitativeCheckResult
def test_parametric_model_checking_hybrid(self):
program = stormpy.parse_prism_program(get_example_path("pdtmc", "parametric_die.pm"))
prop = "P=? [F s=5]"
formulas = stormpy.parse_properties_for_prism_program(prop, program)
model = stormpy.build_symbolic_parametric_model(program, formulas)
assert model.nr_states == 11
assert model.nr_transitions == 17
assert model.model_type == stormpy.ModelType.DTMC
assert model.has_parameters
result = stormpy.check_model_hybrid(model, formulas[0])
assert type(result) is stormpy.HybridParametricQuantitativeCheckResult
values = result.get_values()
assert len(values) == 3
def test_constraints_collector(self): def test_constraints_collector(self):
from pycarl.formula import FormulaType, Relation from pycarl.formula import FormulaType, Relation
if stormpy.info.storm_ratfunc_use_cln(): if stormpy.info.storm_ratfunc_use_cln():

Loading…
Cancel
Save