Browse Source

Binding for CheckResult

refactoring
Matthias Volk 8 years ago
parent
commit
95f83a7224
  1. 12
      src/core/modelchecking.cpp
  2. 47
      src/core/result.cpp
  3. 23
      tests/core/test_modelchecking.py
  4. 9
      tests/storage/test_matrix.py

12
src/core/modelchecking.cpp

@ -2,15 +2,8 @@
#include "result.h" #include "result.h"
// Thin wrapper for model checking // Thin wrapper for model checking
double modelChecking(std::shared_ptr<storm::models::sparse::Model<double>> model, std::shared_ptr<storm::logic::Formula const> const& formula) {
std::unique_ptr<storm::modelchecker::CheckResult> checkResult = storm::verifySparseModel<double>(model, formula);
return checkResult->asExplicitQuantitativeCheckResult<double>()[*model->getInitialStates().begin()];
}
// Thin wrapper for model checking for all states
std::vector<double> modelCheckingAll(std::shared_ptr<storm::models::sparse::Model<double>> model, std::shared_ptr<storm::logic::Formula const> const& formula) {
std::unique_ptr<storm::modelchecker::CheckResult> checkResult = storm::verifySparseModel<double>(model, formula);
return checkResult->asExplicitQuantitativeCheckResult<double>().getValueVector();
std::shared_ptr<storm::modelchecker::CheckResult> modelChecking(std::shared_ptr<storm::models::sparse::Model<double>> model, std::shared_ptr<storm::logic::Formula const> const& formula) {
return storm::verifySparseModel<double>(model, formula);
} }
// Thin wrapper for parametric model checking // Thin wrapper for parametric model checking
@ -33,7 +26,6 @@ void define_modelchecking(py::module& m) {
// Model checking // Model checking
m.def("_model_checking", &modelChecking, "Perform model checking", py::arg("model"), py::arg("formula")); m.def("_model_checking", &modelChecking, "Perform model checking", py::arg("model"), py::arg("formula"));
m.def("model_checking_all", &modelCheckingAll, "Perform model checking for all states", py::arg("model"), py::arg("formula"));
m.def("_parametric_model_checking", &parametricModelChecking, "Perform parametric model checking", py::arg("model"), py::arg("formula")); m.def("_parametric_model_checking", &parametricModelChecking, "Perform parametric model checking", py::arg("model"), py::arg("formula"));
m.def("compute_prob01states", &computeProb01, "Compute prob-0-1 states", py::arg("model"), py::arg("phi_states"), py::arg("psi_states")); m.def("compute_prob01states", &computeProb01, "Compute prob-0-1 states", py::arg("model"), py::arg("phi_states"), py::arg("psi_states"));
} }

47
src/core/result.cpp

@ -11,4 +11,51 @@ void define_result(py::module& m) {
.def_property_readonly("constraints_graph_preserving", &PmcResult::getConstraintsGraphPreserving, "Constraints ensuring graph preservation") .def_property_readonly("constraints_graph_preserving", &PmcResult::getConstraintsGraphPreserving, "Constraints ensuring graph preservation")
; ;
// CheckResult
py::class_<storm::modelchecker::CheckResult, std::shared_ptr<storm::modelchecker::CheckResult>> checkResult(m, "_CheckResult", "Base class for all modelchecking results");
checkResult.def_property_readonly("_symbolic", &storm::modelchecker::CheckResult::isSymbolic, "Flag if result is symbolic")
.def_property_readonly("_hybrid", &storm::modelchecker::CheckResult::isHybrid, "Flag if result is hybrid")
.def_property_readonly("_quantitative", &storm::modelchecker::CheckResult::isQuantitative, "Flag if result is quantitative")
.def_property_readonly("_qualitative", &storm::modelchecker::CheckResult::isQualitative, "Flag if result is qualitative")
.def_property_readonly("_explicit_qualitative", &storm::modelchecker::CheckResult::isExplicitQualitativeCheckResult, "Flag if result is explicit qualitative")
.def_property_readonly("_explicit_quantitative", &storm::modelchecker::CheckResult::isExplicitQuantitativeCheckResult, "Flag if result is explicit quantitative")
.def_property_readonly("_symbolic_qualitative", &storm::modelchecker::CheckResult::isSymbolicQualitativeCheckResult, "Flag if result is symbolic qualitative")
.def_property_readonly("_symbolic_quantitative", &storm::modelchecker::CheckResult::isSymbolicQuantitativeCheckResult, "Flag if result is symbolic quantitative")
.def_property_readonly("_hybrid_quantitative", &storm::modelchecker::CheckResult::isHybridQuantitativeCheckResult, "Flag if result is hybrid quantitative")
.def_property_readonly("_pareto_curve", &storm::modelchecker::CheckResult::isParetoCurveCheckResult, "Flag if result is a pareto curve")
.def_property_readonly("result_for_all_states", &storm::modelchecker::CheckResult::isResultForAllStates, "Flag if result is for all states")
.def("as_qualitative", [](storm::modelchecker::CheckResult const& result) {
return result.asQualitativeCheckResult();
}, "Convert into qualitative result")
.def("as_quantitative", [](storm::modelchecker::CheckResult const& result) {
return result.asQuantitativeCheckResult<double>();
}, "Convert into quantitative result")
.def("as_explicit_qualitative", [](storm::modelchecker::CheckResult const& result) {
return result.asExplicitQualitativeCheckResult();
}, "Convert into explicit qualitative result")
.def("as_explicit_quantitative", [](storm::modelchecker::CheckResult const& result) {
return result.asExplicitQuantitativeCheckResult<double>();
}, "Convert into explicit quantitative result")
.def("__str__", [](storm::modelchecker::CheckResult const& result) {
std::stringstream stream;
result.writeToStream(stream);
return stream.str();
})
;
// QualitativeCheckResult
py::class_<storm::modelchecker::QualitativeCheckResult, std::shared_ptr<storm::modelchecker::QualitativeCheckResult>> qualitativeCheckResult(m, "_QualitativeCheckResult", "Abstract class for qualitative model checking results", checkResult);
py::class_<storm::modelchecker::ExplicitQualitativeCheckResult, std::shared_ptr<storm::modelchecker::ExplicitQualitativeCheckResult>>(m, "ExplicitQualitativeCheckResult", "Explicit qualitative model checking result", qualitativeCheckResult)
.def("get_truth_values", &storm::modelchecker::ExplicitQualitativeCheckResult::getTruthValuesVector, "Get BitVector representing the truth values")
;
// 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::ExplicitQuantitativeCheckResult<double>, std::shared_ptr<storm::modelchecker::ExplicitQuantitativeCheckResult<double>>>(m, "ExplicitQuantitativeCheckResult", "Explicit quantitative model checking result", quantitativeCheckResult)
.def("at", [](storm::modelchecker::ExplicitQuantitativeCheckResult<double> const& result, storm::storage::sparse::state_type state) {
return result[state];
}, py::arg("state"), "Get result for given state")
.def("get_values", &storm::modelchecker::ExplicitQuantitativeCheckResult<double>::getValueVector, "Get model checking result values for all states")
;
} }

23
tests/core/test_modelchecking.py

@ -11,8 +11,11 @@ class TestModelChecking:
model = stormpy.build_model(program, formulas[0]) model = stormpy.build_model(program, formulas[0])
assert model.nr_states == 13 assert model.nr_states == 13
assert model.nr_transitions == 20 assert model.nr_transitions == 20
assert len(model.initial_states) == 1
initial_state = model.initial_states[0]
assert initial_state == 0
result = stormpy.model_checking(model, formulas[0]) result = stormpy.model_checking(model, formulas[0])
assert math.isclose(result, 0.16666666666666663)
assert math.isclose(result.at(initial_state), 0.16666666666666663)
def test_model_checking_all_dtmc(self): def test_model_checking_all_dtmc(self):
program = stormpy.parse_prism_program(get_example_path("dtmc", "die.pm")) program = stormpy.parse_prism_program(get_example_path("dtmc", "die.pm"))
@ -20,9 +23,10 @@ class TestModelChecking:
model = stormpy.build_model(program, formulas[0]) model = stormpy.build_model(program, formulas[0])
assert model.nr_states == 13 assert model.nr_states == 13
assert model.nr_transitions == 20 assert model.nr_transitions == 20
results = stormpy.model_checking_all(model, formulas[0])
result = stormpy.model_checking(model, formulas[0])
assert result.result_for_all_states
reference = [0.16666666666666663, 0.3333333333333333, 0, 0.6666666666666666, 0, 0, 0, 1, 0, 0, 0, 0, 0] reference = [0.16666666666666663, 0.3333333333333333, 0, 0.6666666666666666, 0, 0, 0, 1, 0, 0, 0, 0, 0]
assert all(map(math.isclose, results, reference))
assert all(map(math.isclose, result.get_values(), reference))
def test_parametric_state_elimination(self): def test_parametric_state_elimination(self):
import pycarl import pycarl
@ -48,10 +52,15 @@ class TestModelChecking:
def test_model_checking_prob01(self): def test_model_checking_prob01(self):
program = stormpy.parse_prism_program(get_example_path("dtmc", "die.pm")) program = stormpy.parse_prism_program(get_example_path("dtmc", "die.pm"))
formulas = stormpy.parse_formulas_for_prism_program("P=? [ F \"one\" ]", program)
model = stormpy.build_model(program, formulas[0])
phiStates = stormpy.BitVector(model.nr_states, True)
psiStates = stormpy.BitVector(model.nr_states, [model.nr_states-1])
formulaPhi = stormpy.parse_formulas("true")[0]
formulaPsi = stormpy.parse_formulas("\"six\"")[0]
model = stormpy.build_model(program, formulaPsi)
phiResult = stormpy.model_checking(model, formulaPhi)
phiStates = phiResult.get_truth_values()
assert phiStates.number_of_set_bits() == model.nr_states
psiResult = stormpy.model_checking(model, formulaPsi)
psiStates = psiResult.get_truth_values()
assert psiStates.number_of_set_bits() == 1
(prob0, prob1) = stormpy.compute_prob01states(model, phiStates, psiStates) (prob0, prob1) = stormpy.compute_prob01states(model, phiStates, psiStates)
assert prob0.number_of_set_bits() == 9 assert prob0.number_of_set_bits() == 9
assert prob1.number_of_set_bits() == 1 assert prob1.number_of_set_bits() == 1

9
tests/storage/test_matrix.py

@ -38,7 +38,8 @@ class TestMatrix:
# First model checking # First model checking
formulas = stormpy.parse_formulas("P=? [ F \"one\" ]") formulas = stormpy.parse_formulas("P=? [ F \"one\" ]")
result = stormpy.model_checking(model, formulas[0]) result = stormpy.model_checking(model, formulas[0])
assert math.isclose(result, 0.16666666666666663)
resValue = result.at(model.initial_states[0])
assert math.isclose(resValue, 0.16666666666666663)
# Change probabilities # Change probabilities
i = 0 i = 0
@ -53,7 +54,8 @@ class TestMatrix:
assert e.value() == 0.3 or e.value() == 0.7 or e.value() == 1 or e.value() == 0 assert e.value() == 0.3 or e.value() == 0.7 or e.value() == 1 or e.value() == 0
# Second model checking # Second model checking
result = stormpy.model_checking(model, formulas[0]) result = stormpy.model_checking(model, formulas[0])
assert math.isclose(result, 0.06923076923076932)
resValue = result.at(model.initial_states[0])
assert math.isclose(resValue, 0.06923076923076932)
# Change probabilities again # Change probabilities again
for state in stormpy.state.State(0, model): for state in stormpy.state.State(0, model):
@ -65,7 +67,8 @@ class TestMatrix:
transition.set_value(0.2) transition.set_value(0.2)
# Third model checking # Third model checking
result = stormpy.model_checking(model, formulas[0]) result = stormpy.model_checking(model, formulas[0])
assert math.isclose(result, 0.3555555555555556)
resValue = result.at(model.initial_states[0])
assert math.isclose(resValue, 0.3555555555555556)
def test_change_parametric_sparse_matrix_modelchecking(self): def test_change_parametric_sparse_matrix_modelchecking(self):
import stormpy.logic import stormpy.logic

Loading…
Cancel
Save