19template<
typename ValueType,
typename ConstantType>
21 std::vector<std::shared_ptr<logic::Formula const>> formulas,
23 double const& precision,
bool dotOutput)
24 : assumptionMaker(model->getTransitionMatrix()) {
25 assert(model !=
nullptr);
28 this->formulas = formulas;
29 this->precision = utility::convertNumber<ConstantType>(precision);
30 this->matrix = model->getTransitionMatrix();
31 this->dotOutput = dotOutput;
33 if (regions.size() == 1) {
34 this->region = *(regions.begin());
38 std::set<VariableType> vars;
40 for (
auto var : vars) {
43 lowerBoundaries.insert(std::make_pair(var, lb));
44 upperBoundaries.insert(std::make_pair(var, ub));
49 if (numberOfSamples > 2) {
58 if (numberOfSamples > 0) {
59 STORM_LOG_WARN(
"At least 3 sample points are needed to check for monotonicity on samples, not using samples for now");
66 for (uint_fast64_t i = 0; i < matrix.getRowCount(); ++i) {
67 std::set<VariableType> occurringVariables;
69 for (
auto& entry : matrix.getRow(i)) {
72 for (
auto& var : occurringVariables) {
73 occuringStatesAtVariable[var].push_back(i);
79template<
typename ValueType,
typename ConstantType>
80std::map<std::shared_ptr<Order>, std::pair<std::shared_ptr<MonotonicityResult<typename MonotonicityHelper<ValueType, ConstantType>::VariableType>>,
81 std::vector<std::shared_ptr<expressions::BinaryRelationExpression>>>>
85 this->extender->initializeMinMaxValues(region);
87 STORM_PRINT(
"\nTotal time for pla checking: " << plaWatch <<
".\n\n");
92 for (
auto itr : monResults) {
93 if (itr.first !=
nullptr) {
94 std::cout <<
"Number of done states: " << itr.first->getNumberOfDoneStates() <<
'\n';
97 for (
auto& entry : resultCheckOnSamples.getMonotonicityResult()) {
98 if (entry.second == Monotonicity::Not) {
99 itr.second.first->updateMonotonicityResult(entry.first, entry.second,
true);
103 std::string temp = itr.second.first->toString();
105 for (
auto& assumption : itr.second.second) {
109 outfile <<
"Assumptions: \n"
113 outfile << *assumption;
118 outfile <<
"No Assumptions\n";
120 outfile <<
"Monotonicity Result: \n"
121 <<
" " << temp <<
"\n\n";
124 if (monResults.size() == 0) {
125 outfile <<
"No monotonicity found, as the order is insufficient\n";
127 outfile <<
"Monotonicity Result on samples: " << resultCheckOnSamples.toString() <<
'\n';
133 STORM_LOG_WARN_COND(monResults.size() <= 10,
"Too many Reachability Orders. Dot Output will only be created for 10.");
135 auto orderItr = monResults.begin();
136 while (i < 10 && orderItr != monResults.end()) {
137 std::ofstream dotOutfile;
138 std::string name = dotOutfileName + std::to_string(i);
140 dotOutfile <<
"Assumptions:\n";
141 auto assumptionItr = orderItr->second.second.begin();
142 while (assumptionItr != orderItr->second.second.end()) {
143 dotOutfile << *assumptionItr <<
'\n';
148 orderItr->first->dotOutputToFile(dotOutfile);
157template<
typename ValueType,
typename ConstantType>
158std::shared_ptr<LocalMonotonicityResult<typename MonotonicityHelper<ValueType, ConstantType>::VariableType>>
161 for (uint_fast64_t state = 0; state < model->getNumberOfStates(); ++state) {
162 for (
auto& var : extender->getVariablesOccuringAtState()[state]) {
163 localMonRes.
setMonotonicity(state, var, extender->getMonotoncityChecker().checkLocalMonotonicity(order, state, var, region));
166 localMonRes.
setDone(order->getDoneBuilding());
167 return std::make_shared<LocalMonotonicityResult<VariableType>>(localMonRes);
171template<
typename ValueType,
typename ConstantType>
174 std::tuple<std::shared_ptr<Order>, uint_fast64_t, uint_fast64_t> criticalTuple;
178 criticalTuple = extender->toOrder(region, monRes);
180 std::map<std::shared_ptr<Order>, std::vector<std::shared_ptr<expressions::BinaryRelationExpression>>> result;
182 auto val1 = std::get<1>(criticalTuple);
183 auto val2 = std::get<2>(criticalTuple);
184 auto numberOfStates = model->getNumberOfStates();
185 std::vector<std::shared_ptr<expressions::BinaryRelationExpression>> assumptions;
187 if (val1 == numberOfStates && val2 == numberOfStates) {
188 auto resAssumptionPair =
189 std::pair<std::shared_ptr<MonotonicityResult<VariableType>>, std::vector<std::shared_ptr<expressions::BinaryRelationExpression>>>(monRes,
192 std::pair<std::shared_ptr<Order>,
194 std::get<0>(criticalTuple), resAssumptionPair));
195 }
else if (val1 != numberOfStates && val2 != numberOfStates) {
196 extendOrderWithAssumptions(std::get<0>(criticalTuple), val1, val2, assumptions, monRes);
202template<
typename ValueType,
typename ConstantType>
203void MonotonicityHelper<ValueType, ConstantType>::extendOrderWithAssumptions(std::shared_ptr<Order> order, uint_fast64_t val1, uint_fast64_t val2,
204 std::vector<std::shared_ptr<expressions::BinaryRelationExpression>> assumptions,
205 std::shared_ptr<MonotonicityResult<VariableType>> monRes) {
206 std::map<std::shared_ptr<Order>, std::vector<std::shared_ptr<expressions::BinaryRelationExpression>>> result;
207 if (order->isInvalid()) {
212 auto numberOfStates = model->getNumberOfStates();
213 if (val1 == numberOfStates || val2 == numberOfStates) {
214 assert(val1 == val2);
215 assert(order->getNumberOfAddedStates() == order->getNumberOfStates());
216 auto resAssumptionPair =
217 std::pair<std::shared_ptr<MonotonicityResult<VariableType>>, std::vector<std::shared_ptr<expressions::BinaryRelationExpression>>>(monRes,
220 std::pair<std::shared_ptr<Order>,
221 std::pair<std::shared_ptr<MonotonicityResult<VariableType>>, std::vector<std::shared_ptr<expressions::BinaryRelationExpression>>>>(
222 std::move(order), std::move(resAssumptionPair)));
225 STORM_LOG_INFO(
"Creating assumptions for " << val1 <<
" and " << val2 <<
". ");
226 auto newAssumptions = assumptionMaker.createAndCheckAssumptions(val1, val2, order, region);
227 assert(newAssumptions.size() <= 3);
228 auto itr = newAssumptions.begin();
229 if (newAssumptions.size() == 0) {
230 monRes = std::make_shared<MonotonicityResult<VariableType>>(MonotonicityResult<VariableType>());
231 for (
auto& entry : occuringStatesAtVariable) {
232 for (
auto& state : entry.second) {
233 extender->checkParOnStateMonRes(state, order, entry.first, monRes);
234 if (monRes->getMonotonicity(entry.first) == Monotonicity::Unknown) {
238 monRes->setDoneForVar(entry.first);
240 monResults.insert({order, {monRes, assumptions}});
241 STORM_LOG_INFO(
" None of the assumptions were valid, we stop exploring the current order");
243 STORM_LOG_INFO(
" Created " << newAssumptions.size() <<
" assumptions, we continue extending the current order");
246 while (itr != newAssumptions.end()) {
247 auto assumption = *itr;
250 if (itr != newAssumptions.end()) {
252 auto orderCopy = order->copy();
253 auto assumptionsCopy = std::vector<std::shared_ptr<expressions::BinaryRelationExpression>>(assumptions);
254 auto monResCopy = monRes->copy();
258 assumptionsCopy.push_back(std::move(assumption.first));
260 auto criticalTuple = extender->extendOrder(orderCopy, region, monResCopy, assumption.first);
261 extendOrderWithAssumptions(std::get<0>(criticalTuple), std::get<1>(criticalTuple), std::get<2>(criticalTuple), assumptionsCopy, monResCopy);
266 assumptions.push_back(std::move(assumption.first));
268 auto criticalTuple = extender->extendOrder(order, region, monRes, assumption.first);
269 extendOrderWithAssumptions(std::get<0>(criticalTuple), std::get<1>(criticalTuple), std::get<2>(criticalTuple), assumptions, monRes);
276template<
typename ValueType,
typename ConstantType>
277void MonotonicityHelper<ValueType, ConstantType>::checkMonotonicityOnSamples(std::shared_ptr<models::sparse::Dtmc<ValueType>> model,
278 uint_fast64_t numberOfSamples) {
279 assert(numberOfSamples > 2);
281 auto instantiator = utility::ModelInstantiator<models::sparse::Dtmc<ValueType>, models::sparse::Dtmc<ConstantType>>(*model);
283 std::vector<std::vector<ConstantType>> samples;
285 for (
auto itr = variables.begin(); itr != variables.end(); ++itr) {
286 ConstantType previous = -1;
292 for (uint_fast64_t i = 0; (monDecr || monIncr) && i < numberOfSamples; ++
i) {
294 auto valuation = utility::parametric::Valuation<ValueType>();
295 for (
auto itr2 = variables.begin(); itr2 != variables.end(); ++itr2) {
297 if ((*itr) == (*itr2)) {
298 auto lb = region.getLowerBoundary(itr->name());
299 auto ub = region.getUpperBoundary(itr->name());
301 valuation[*itr2] = (lb + utility::convertNumber<CoefficientType>(i / (numberOfSamples - 1)) * (ub - lb));
303 auto lb = region.getLowerBoundary(itr2->name());
304 valuation[*itr2] = utility::convertNumber<typename utility::parametric::CoefficientType<ValueType>::type>(lb);
309 models::sparse::Dtmc<ConstantType> sampleModel = instantiator.instantiate(valuation);
310 auto checker = modelchecker::SparseDtmcPrctlModelChecker<models::sparse::Dtmc<ConstantType>>(sampleModel);
311 std::unique_ptr<modelchecker::CheckResult> checkResult;
312 auto formula = formulas[0];
313 if (formula->isProbabilityOperatorFormula() && formula->asProbabilityOperatorFormula().getSubformula().isUntilFormula()) {
314 const modelchecker::CheckTask<logic::UntilFormula, ConstantType> checkTask =
315 modelchecker::CheckTask<logic::UntilFormula, ConstantType>((*formula).asProbabilityOperatorFormula().getSubformula().asUntilFormula());
316 checkResult = checker.computeUntilProbabilities(Environment(), checkTask);
317 }
else if (formula->isProbabilityOperatorFormula() && formula->asProbabilityOperatorFormula().getSubformula().isEventuallyFormula()) {
318 const modelchecker::CheckTask<logic::EventuallyFormula, ConstantType> checkTask =
319 modelchecker::CheckTask<logic::EventuallyFormula, ConstantType>(
320 (*formula).asProbabilityOperatorFormula().getSubformula().asEventuallyFormula());
321 checkResult = checker.computeReachabilityProbabilities(Environment(), checkTask);
323 STORM_LOG_THROW(
false, exceptions::NotSupportedException,
"Expecting until or eventually formula");
326 auto quantitativeResult = checkResult->asExplicitQuantitativeCheckResult<ConstantType>();
327 std::vector<ConstantType> values = quantitativeResult.getValueVector();
328 auto initialStates = model->getInitialStates();
329 ConstantType initial = 0;
331 for (
auto j = initialStates.getNextSetIndex(0); j < model->getNumberOfStates(); j = initialStates.getNextSetIndex(j + 1)) {
332 initial += values[j];
335 assert(initial >= 0 - precision && initial <= 1 + precision);
336 ConstantType diff = previous - initial;
337 assert(previous == -1 || (diff >= -1 - precision && diff <= 1 + precision));
339 if (previous != -1 && (diff > precision || diff < -precision)) {
340 monDecr &= diff > precision;
341 monIncr &= diff < -precision;
344 samples.push_back(std::move(values));
347 resultCheckOnSamples.addMonotonicityResult(*itr, res);
349 assumptionMaker.setSampleValues(std::move(samples));
352template<
typename ValueType,
typename ConstantType>
353void MonotonicityHelper<ValueType, ConstantType>::checkMonotonicityOnSamples(std::shared_ptr<models::sparse::Mdp<ValueType>> model,
354 uint_fast64_t numberOfSamples) {
355 STORM_LOG_THROW(
false, storm::exceptions::NotImplementedException,
"Checking monotonicity on samples not implemented for mdps");
358template class MonotonicityHelper<RationalFunction, double>;
359template class MonotonicityHelper<RationalFunction, RationalNumber>;
void setDone(bool done=true)
void setMonotonicity(uint_fast64_t state, VariableType var, Monotonicity mon)
Sets the local Monotonicity of a parameter at a given state.
std::map< std::shared_ptr< Order >, std::pair< std::shared_ptr< MonotonicityResult< VariableType > >, std::vector< std::shared_ptr< expressions::BinaryRelationExpression > > > > checkMonotonicityInBuild(std::ostream &outfile, bool usePLA=false, std::string dotOutfileName="dotOutput")
Builds Reachability Orders for the given model and simultaneously uses them to check for Monotonicity...
std::shared_ptr< LocalMonotonicityResult< VariableType > > createLocalMonotonicityResult(std::shared_ptr< Order > order, storage::ParameterRegion< ValueType > region)
Builds Reachability Orders for the given model and simultaneously uses them to check for Monotonicity...
MonotonicityHelper(std::shared_ptr< models::sparse::Model< ValueType > > model, std::vector< std::shared_ptr< logic::Formula const > > formulas, std::vector< storage::ParameterRegion< ValueType > > regions, uint_fast64_t numberOfSamples=0, double const &precision=0.000001, bool dotOutput=false)
Constructor of MonotonicityHelper.
@ Not
the result is not monotonic
This class represents a discrete-time Markov chain.
This class represents a (discrete-time) Markov decision process.
Base class for all sparse models.
storm::utility::parametric::CoefficientType< ParametricType >::type CoefficientType
storm::utility::parametric::Valuation< ParametricType > Valuation
A class that provides convenience operations to display run times.
void stop()
Stop stopwatch and add measured time to total time.
#define STORM_LOG_INFO(message)
#define STORM_LOG_WARN(message)
#define STORM_LOG_WARN_COND(cond, message)
#define STORM_LOG_THROW(cond, exception, message)
#define STORM_PRINT(message)
Define the macros that print information and optionally also log it.
typename utility::parametric::VariableType< FunctionType >::type VariableType
void closeFile(std::ofstream &stream)
Close the given file after writing.
void openFile(std::string const &filepath, std::ofstream &filestream, bool append=false, bool silent=false)
Open the given file for writing.
@ Unknown
the result is unknown
std::set< storm::RationalFunctionVariable > getProbabilityParameters(Model< storm::RationalFunction > const &model)
Get all probability parameters occurring on transitions.
void gatherOccurringVariables(FunctionType const &function, std::set< typename VariableType< FunctionType >::type > &variableSet)
Add all variables that occur in the given function to the the given set.