Storm
A Modern Probabilistic Model Checker
Loading...
Searching...
No Matches
MonotonicityHelper.cpp
Go to the documentation of this file.
2
5
9
12
15
16namespace storm {
17namespace analysis {
18/*** Constructor ***/
19template<typename ValueType, typename ConstantType>
21 std::vector<std::shared_ptr<logic::Formula const>> formulas,
22 std::vector<storage::ParameterRegion<ValueType>> regions, uint_fast64_t numberOfSamples,
23 double const& precision, bool dotOutput)
24 : assumptionMaker(model->getTransitionMatrix()) {
25 assert(model != nullptr);
26
27 this->model = model;
28 this->formulas = formulas;
29 this->precision = utility::convertNumber<ConstantType>(precision);
30 this->matrix = model->getTransitionMatrix();
31 this->dotOutput = dotOutput;
32
33 if (regions.size() == 1) {
34 this->region = *(regions.begin());
35 } else {
38 std::set<VariableType> vars;
40 for (auto var : vars) {
41 typename storage::ParameterRegion<ValueType>::CoefficientType lb = utility::convertNumber<CoefficientType>(0 + precision);
42 typename storage::ParameterRegion<ValueType>::CoefficientType ub = utility::convertNumber<CoefficientType>(1 - precision);
43 lowerBoundaries.insert(std::make_pair(var, lb));
44 upperBoundaries.insert(std::make_pair(var, ub));
45 }
46 this->region = storage::ParameterRegion<ValueType>(std::move(lowerBoundaries), std::move(upperBoundaries));
47 }
48
49 if (numberOfSamples > 2) {
50 // sampling
51 if (model->isOfType(models::ModelType::Dtmc)) {
52 checkMonotonicityOnSamples(model->template as<models::sparse::Dtmc<ValueType>>(), numberOfSamples);
53 } else if (model->isOfType(models::ModelType::Mdp)) {
54 checkMonotonicityOnSamples(model->template as<models::sparse::Mdp<ValueType>>(), numberOfSamples);
55 }
56 checkSamples = true;
57 } else {
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");
60 }
61 checkSamples = false;
62 }
63
64 this->extender = new analysis::OrderExtender<ValueType, ConstantType>(model, formulas[0]);
65
66 for (uint_fast64_t i = 0; i < matrix.getRowCount(); ++i) {
67 std::set<VariableType> occurringVariables;
68
69 for (auto& entry : matrix.getRow(i)) {
70 storm::utility::parametric::gatherOccurringVariables(entry.getValue(), occurringVariables);
71 }
72 for (auto& var : occurringVariables) {
73 occuringStatesAtVariable[var].push_back(i);
74 }
75 }
76}
77
78/*** Public methods ***/
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>>>>
82MonotonicityHelper<ValueType, ConstantType>::checkMonotonicityInBuild(std::ostream& outfile, bool usePLA, std::string dotOutfileName) {
83 if (usePLA) {
84 storm::utility::Stopwatch plaWatch(true);
85 this->extender->initializeMinMaxValues(region);
86 plaWatch.stop();
87 STORM_PRINT("\nTotal time for pla checking: " << plaWatch << ".\n\n");
88 }
89 createOrder();
90
91 // output of results
92 for (auto itr : monResults) {
93 if (itr.first != nullptr) {
94 std::cout << "Number of done states: " << itr.first->getNumberOfDoneStates() << '\n';
95 }
96 if (checkSamples) {
97 for (auto& entry : resultCheckOnSamples.getMonotonicityResult()) {
98 if (entry.second == Monotonicity::Not) {
99 itr.second.first->updateMonotonicityResult(entry.first, entry.second, true);
100 }
101 }
102 }
103 std::string temp = itr.second.first->toString();
104 bool first = true;
105 for (auto& assumption : itr.second.second) {
106 if (!first) {
107 outfile << " & ";
108 } else {
109 outfile << "Assumptions: \n"
110 << " ";
111 first = false;
112 }
113 outfile << *assumption;
114 }
115 if (!first) {
116 outfile << '\n';
117 } else {
118 outfile << "No Assumptions\n";
119 }
120 outfile << "Monotonicity Result: \n"
121 << " " << temp << "\n\n";
122 }
123
124 if (monResults.size() == 0) {
125 outfile << "No monotonicity found, as the order is insufficient\n";
126 if (checkSamples) {
127 outfile << "Monotonicity Result on samples: " << resultCheckOnSamples.toString() << '\n';
128 }
129 }
130
131 // dotoutput
132 if (dotOutput) {
133 STORM_LOG_WARN_COND(monResults.size() <= 10, "Too many Reachability Orders. Dot Output will only be created for 10.");
134 int i = 0;
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);
139 storm::io::openFile(name, dotOutfile);
140 dotOutfile << "Assumptions:\n";
141 auto assumptionItr = orderItr->second.second.begin();
142 while (assumptionItr != orderItr->second.second.end()) {
143 dotOutfile << *assumptionItr << '\n';
144 dotOutfile << '\n';
145 assumptionItr++;
146 }
147 dotOutfile << '\n';
148 orderItr->first->dotOutputToFile(dotOutfile);
149 storm::io::closeFile(dotOutfile);
150 i++;
151 orderItr++;
152 }
153 }
154 return monResults;
155}
156
157template<typename ValueType, typename ConstantType>
158std::shared_ptr<LocalMonotonicityResult<typename MonotonicityHelper<ValueType, ConstantType>::VariableType>>
160 LocalMonotonicityResult<VariableType> localMonRes(model->getNumberOfStates());
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));
164 }
165 }
166 localMonRes.setDone(order->getDoneBuilding());
167 return std::make_shared<LocalMonotonicityResult<VariableType>>(localMonRes);
168}
169
170/*** Private methods ***/
171template<typename ValueType, typename ConstantType>
173 // Transform to Orders
174 std::tuple<std::shared_ptr<Order>, uint_fast64_t, uint_fast64_t> criticalTuple;
175
176 // Create initial order
177 auto monRes = std::make_shared<MonotonicityResult<VariableType>>(MonotonicityResult<VariableType>());
178 criticalTuple = extender->toOrder(region, monRes);
179 // Continue based on not (yet) sorted states
180 std::map<std::shared_ptr<Order>, std::vector<std::shared_ptr<expressions::BinaryRelationExpression>>> result;
181
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;
186
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,
190 assumptions);
191 monResults.insert(
192 std::pair<std::shared_ptr<Order>,
193 std::pair<std::shared_ptr<MonotonicityResult<VariableType>>, std::vector<std::shared_ptr<expressions::BinaryRelationExpression>>>>(
194 std::get<0>(criticalTuple), resAssumptionPair));
195 } else if (val1 != numberOfStates && val2 != numberOfStates) {
196 extendOrderWithAssumptions(std::get<0>(criticalTuple), val1, val2, assumptions, monRes);
197 } else {
198 assert(false);
199 }
200}
201
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()) {
208 // We don't add anything as the order we created with assumptions turns out to be invalid
209 STORM_LOG_INFO(" The order was invalid, so we stop here");
210 return;
211 }
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,
218 assumptions);
219 monResults.insert(
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)));
223 } else {
224 // Make the three assumptions
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) {
235 break;
236 }
237 }
238 monRes->setDoneForVar(entry.first);
239 }
240 monResults.insert({order, {monRes, assumptions}});
241 STORM_LOG_INFO(" None of the assumptions were valid, we stop exploring the current order");
242 } else {
243 STORM_LOG_INFO(" Created " << newAssumptions.size() << " assumptions, we continue extending the current order");
244 }
245
246 while (itr != newAssumptions.end()) {
247 auto assumption = *itr;
248 ++itr;
249 if (assumption.second != AssumptionStatus::INVALID) {
250 if (itr != newAssumptions.end()) {
251 // We make a copy of the order and the assumptions
252 auto orderCopy = order->copy();
253 auto assumptionsCopy = std::vector<std::shared_ptr<expressions::BinaryRelationExpression>>(assumptions);
254 auto monResCopy = monRes->copy();
255
256 if (assumption.second == AssumptionStatus::UNKNOWN) {
257 // only add assumption to the set of assumptions if it is unknown whether it holds or not
258 assumptionsCopy.push_back(std::move(assumption.first));
259 }
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);
262 } else {
263 // It is the last one, so we don't need to create a copy.
264 if (assumption.second == AssumptionStatus::UNKNOWN) {
265 // only add assumption to the set of assumptions if it is unknown whether it holds or not
266 assumptions.push_back(std::move(assumption.first));
267 }
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);
270 }
271 }
272 }
273 }
274}
275
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);
280
281 auto instantiator = utility::ModelInstantiator<models::sparse::Dtmc<ValueType>, models::sparse::Dtmc<ConstantType>>(*model);
282 std::set<VariableType> variables = models::sparse::getProbabilityParameters(*model);
283 std::vector<std::vector<ConstantType>> samples;
284 // For each of the variables create a model in which we only change the value for this specific variable
285 for (auto itr = variables.begin(); itr != variables.end(); ++itr) {
286 ConstantType previous = -1;
287 bool monDecr = true;
288 bool monIncr = true;
289
290 // Check monotonicity in variable (*itr) by instantiating the model
291 // all other variables fixed on lb, only increasing (*itr)
292 for (uint_fast64_t i = 0; (monDecr || monIncr) && i < numberOfSamples; ++i) {
293 // Create valuation
294 auto valuation = utility::parametric::Valuation<ValueType>();
295 for (auto itr2 = variables.begin(); itr2 != variables.end(); ++itr2) {
296 // Only change value for current variable
297 if ((*itr) == (*itr2)) {
298 auto lb = region.getLowerBoundary(itr->name());
299 auto ub = region.getUpperBoundary(itr->name());
300 // Creates samples between lb and ub, that is: lb, lb + (ub-lb)/(#samples -1), lb + 2* (ub-lb)/(#samples -1), ..., ub
301 valuation[*itr2] = (lb + utility::convertNumber<CoefficientType>(i / (numberOfSamples - 1)) * (ub - lb));
302 } else {
303 auto lb = region.getLowerBoundary(itr2->name());
304 valuation[*itr2] = utility::convertNumber<typename utility::parametric::CoefficientType<ValueType>::type>(lb);
305 }
306 }
307
308 // Instantiate model and get result
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);
322 } else {
323 STORM_LOG_THROW(false, exceptions::NotSupportedException, "Expecting until or eventually formula");
324 }
325
326 auto quantitativeResult = checkResult->asExplicitQuantitativeCheckResult<ConstantType>();
327 std::vector<ConstantType> values = quantitativeResult.getValueVector();
328 auto initialStates = model->getInitialStates();
329 ConstantType initial = 0;
330 // Get total probability from initial states
331 for (auto j = initialStates.getNextSetIndex(0); j < model->getNumberOfStates(); j = initialStates.getNextSetIndex(j + 1)) {
332 initial += values[j];
333 }
334 // Calculate difference with result for previous valuation
335 assert(initial >= 0 - precision && initial <= 1 + precision);
336 ConstantType diff = previous - initial;
337 assert(previous == -1 || (diff >= -1 - precision && diff <= 1 + precision));
338
339 if (previous != -1 && (diff > precision || diff < -precision)) {
340 monDecr &= diff > precision; // then previous value is larger than the current value from the initial states
341 monIncr &= diff < -precision;
342 }
343 previous = initial;
344 samples.push_back(std::move(values));
345 }
346 auto res = (!monIncr && !monDecr) ? MonotonicityResult<VariableType>::Monotonicity::Not : MonotonicityResult<VariableType>::Monotonicity::Unknown;
347 resultCheckOnSamples.addMonotonicityResult(*itr, res);
348 }
349 assumptionMaker.setSampleValues(std::move(samples));
350}
351
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");
356}
357
358template class MonotonicityHelper<RationalFunction, double>;
359template class MonotonicityHelper<RationalFunction, RationalNumber>;
360} // namespace analysis
361} // namespace storm
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.
This class represents a discrete-time Markov chain.
Definition Dtmc.h:14
This class represents a (discrete-time) Markov decision process.
Definition Mdp.h:14
Base class for all sparse models.
Definition Model.h:33
storm::utility::parametric::CoefficientType< ParametricType >::type CoefficientType
storm::utility::parametric::Valuation< ParametricType > Valuation
A class that provides convenience operations to display run times.
Definition Stopwatch.h:14
void stop()
Stop stopwatch and add measured time to total time.
Definition Stopwatch.cpp:42
#define STORM_LOG_INFO(message)
Definition logging.h:29
#define STORM_LOG_WARN(message)
Definition logging.h:30
#define STORM_LOG_WARN_COND(cond, message)
Definition macros.h:38
#define STORM_LOG_THROW(cond, exception, message)
Definition macros.h:30
#define STORM_PRINT(message)
Define the macros that print information and optionally also log it.
Definition macros.h:62
typename utility::parametric::VariableType< FunctionType >::type VariableType
void closeFile(std::ofstream &stream)
Close the given file after writing.
Definition file.h:47
void openFile(std::string const &filepath, std::ofstream &filestream, bool append=false, bool silent=false)
Open the given file for writing.
Definition file.h:18
@ Unknown
the result is unknown
std::set< storm::RationalFunctionVariable > getProbabilityParameters(Model< storm::RationalFunction > const &model)
Get all probability parameters occurring on transitions.
Definition Model.cpp:703
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.
LabParser.cpp.
Definition cli.cpp:18