Storm
A Modern Probabilistic Model Checker
Loading...
Searching...
No Matches
QuickHull.cpp
Go to the documentation of this file.
2
3#include <algorithm>
5#include <unordered_map>
6
10
13
15
16namespace storm {
17namespace storage {
18namespace geometry {
19
20template<typename ValueType>
21void QuickHull<ValueType>::generateHalfspacesFromPoints(std::vector<EigenVector>& points, bool generateRelevantVerticesAndVertexSets) {
22 STORM_LOG_ASSERT(!points.empty(), "Invoked QuickHull with empty set of points.");
23 STORM_LOG_DEBUG("Invoked QuickHull on " << points.size() << " points");
24 const uint_fast64_t dimension = points.front().rows();
25
26 if (dimension == 1) {
27 handle1DPoints(points, generateRelevantVerticesAndVertexSets);
28 } else {
29 // Generate initial set of d+1 affine independent points (if such a set exists)
30 std::vector<uint_fast64_t> vertexIndices;
31 if (this->findInitialVertices(points, vertexIndices)) {
32 // compute point inside initial facet
33 EigenVector insidePoint(EigenVector::Zero(dimension));
34 for (uint_fast64_t vertexIndex : vertexIndices) {
35 insidePoint += points[vertexIndex];
36 }
37 insidePoint /= storm::utility::convertNumber<ValueType>(static_cast<uint_fast64_t>(vertexIndices.size()));
38
39 // Create the initial facets from the found vertices.
40 std::vector<Facet> facets = computeInitialFacets(points, vertexIndices, insidePoint);
41
42 // Enlarge the mesh
43 storm::storage::BitVector currentFacets(facets.size(), true);
44 this->extendMesh(points, facets, currentFacets, insidePoint);
45
46 // Finally retrieve the resulting constrants
47 this->getPolytopeFromMesh(points, facets, currentFacets, generateRelevantVerticesAndVertexSets);
48 STORM_LOG_DEBUG("QuickHull invokation yielded " << resultMatrix.rows() << " constraints");
49 } else {
50 // The points are affine dependent. This special case needs to be handled accordingly.
51 handleAffineDependentPoints(points, generateRelevantVerticesAndVertexSets);
52 }
53 }
54}
55
56template<typename ValueType>
58 return this->resultMatrix;
59}
60
61template<typename ValueType>
63 return this->resultVector;
64}
65
66template<typename ValueType>
67std::vector<typename QuickHull<ValueType>::EigenVector>& QuickHull<ValueType>::getRelevantVertices() {
68 return this->relevantVertices;
69}
70
71template<typename ValueType>
72std::vector<std::vector<uint_fast64_t>>& QuickHull<ValueType>::getVertexSets() {
73 return this->vertexSets;
74}
75
76template<typename ValueType>
77void QuickHull<ValueType>::handle1DPoints(std::vector<EigenVector>& points, bool generateRelevantVerticesAndVertexSets) {
78 ValueType minValue = points.front()(0);
79 ValueType maxValue = points.front()(0);
80 uint_fast64_t minIndex = 0;
81 uint_fast64_t maxIndex = 0;
82 for (uint_fast64_t pointIndex = 1; pointIndex < points.size(); ++pointIndex) {
83 ValueType const& pointValue = points[pointIndex](0);
84 if (pointValue < minValue) {
85 minValue = pointValue;
86 minIndex = pointIndex;
87 }
88 if (pointValue > maxValue) {
89 maxValue = pointValue;
90 maxIndex = pointIndex;
91 }
92 }
93 resultMatrix = EigenMatrix(2, 1);
94 resultMatrix(0, 0) = -storm::utility::one<ValueType>();
95 resultMatrix(1, 0) = storm::utility::one<ValueType>();
96 resultVector = EigenVector(2);
97 resultVector(0) = -minValue;
98 resultVector(1) = maxValue;
99 if (generateRelevantVerticesAndVertexSets) {
100 relevantVertices.push_back(points[minIndex]);
101 std::vector<uint_fast64_t> minVertexSet(1, 0);
102 vertexSets.push_back(std::move(minVertexSet));
103 if (minIndex != maxIndex && points[minIndex] != points[maxIndex]) {
104 relevantVertices.push_back(points[maxIndex]);
105 std::vector<uint_fast64_t> maxVertexSet(1, 1);
106 vertexSets.push_back(std::move(maxVertexSet));
107 } else {
108 vertexSets.push_back(vertexSets.front());
109 }
110 }
111}
112
113template<typename ValueType>
114bool QuickHull<ValueType>::affineFilter(std::vector<uint_fast64_t> const& subset, uint_fast64_t const& item, std::vector<EigenVector> const& points) {
115 EigenMatrix vectorMatrix(points[item].rows() + 1, subset.size() + 1);
116 for (uint_fast64_t i = 0; i < subset.size(); ++i) {
117 vectorMatrix.col(i) << points[subset[i]], storm::utility::one<ValueType>();
118 }
119 vectorMatrix.col(subset.size()) << points[item], storm::utility::one<ValueType>();
120 return (vectorMatrix.fullPivLu().rank() > (Eigen::Index)subset.size());
121}
122
123template<typename ValueType>
124void QuickHull<ValueType>::handleAffineDependentPoints(std::vector<EigenVector>& points, bool generateRelevantVerticesAndVertexSets) {
125 bool pointsAffineDependent = true;
126 std::vector<std::pair<EigenVector, ValueType>> additionalConstraints;
127 while (pointsAffineDependent) {
128 // get one hyperplane that holds all points
129 const uint_fast64_t dimension = points.front().rows();
130 EigenVector refPoint = points.front();
131 EigenVector normal;
132 if (points.size() == 1) {
133 normal.resize(dimension);
134 normal(0) = storm::utility::one<ValueType>();
135 } else {
136 EigenMatrix constraints(points.size() - 1, dimension);
137 for (unsigned row = 1; row < points.size(); ++row) {
138 constraints.row(row - 1) = points[row] - refPoint;
139 }
140 normal = constraints.fullPivLu().kernel().col(0);
141 }
142
143 // Eigen returns the column vector 0...0 if the kernel is empty (i.e., there is no such hyperplane)
144 if (normal.isZero()) {
145 pointsAffineDependent = false;
146 } else {
147 points.push_back(refPoint + normal);
148 ValueType offset = normal.dot(refPoint);
149 additionalConstraints.emplace_back(std::move(normal), std::move(offset));
150 }
151 }
152 // Now the points are affine independent so we can relaunch the procedure
153 generateHalfspacesFromPoints(points, generateRelevantVerticesAndVertexSets);
154
155 // Add the constraints obtained above
156 uint_fast64_t numOfAdditionalConstraints = additionalConstraints.size();
157 uint_fast64_t numOfRegularConstraints = resultMatrix.rows();
158 EigenMatrix newA(numOfRegularConstraints + numOfAdditionalConstraints, resultMatrix.cols());
159 EigenVector newb(static_cast<unsigned long int>(numOfRegularConstraints + numOfAdditionalConstraints));
160 uint_fast64_t row = 0;
161 for (; row < numOfRegularConstraints; ++row) {
162 newA.row(row) = resultMatrix.row(row);
163 newb(row) = resultVector(row);
164 }
165 for (; row < numOfRegularConstraints + numOfAdditionalConstraints; ++row) {
166 newA.row(row) = std::move(additionalConstraints[row - numOfRegularConstraints].first);
167 newb(row) = additionalConstraints[row - numOfRegularConstraints].second;
168 }
169 resultMatrix = std::move(newA);
170 resultVector = std::move(newb);
171
172 // clear the additionally added points. Note that the order of the points might have changed
173 storm::storage::BitVector keptPoints(points.size(), true);
174 for (uint_fast64_t pointIndex = 0; pointIndex < points.size(); ++pointIndex) {
175 for (Eigen::Index row = 0; row < resultMatrix.rows(); ++row) {
176 if ((resultMatrix.row(row) * points[pointIndex])(0) > resultVector(row)) {
177 keptPoints.set(pointIndex, false);
178 break;
179 }
180 }
181 }
183
184 if (generateRelevantVerticesAndVertexSets) {
185 storm::storage::BitVector keptVertices(relevantVertices.size(), true);
186 for (uint_fast64_t vertexIndex = 0; vertexIndex < relevantVertices.size(); ++vertexIndex) {
187 for (Eigen::Index row = 0; row < resultMatrix.rows(); ++row) {
188 if ((resultMatrix.row(row) * relevantVertices[vertexIndex])(0) > resultVector(row)) {
189 keptVertices.set(vertexIndex, false);
190 break;
191 }
192 }
193 }
194 storm::utility::vector::filterVectorInPlace(relevantVertices, keptVertices);
195
196 STORM_LOG_WARN("Can not retrieve vertex sets for degenerated polytope (not implemented)");
197 vertexSets.clear();
198 }
199}
200
201template<typename ValueType>
202bool QuickHull<ValueType>::findInitialVertices(std::vector<EigenVector>& points, std::vector<uint_fast64_t>& verticesOfInitialPolytope) const {
203 const uint_fast64_t dimension = points.front().rows();
204 if (points.size() < dimension + 1) {
205 // not enough points to obtain a (non-degenerated) polytope
206 return false;
207 }
208
209 // We first find some good candidates to get a (hopefully) large initial mesh.
210 storm::storage::BitVector notGoodCandidates(points.size(), true);
211 for (uint_fast64_t dim = 0; dim < dimension; ++dim) {
212 if (!notGoodCandidates.empty()) {
213 uint_fast64_t minIndex = *notGoodCandidates.begin();
214 uint_fast64_t maxIndex = *notGoodCandidates.begin();
215 for (uint_fast64_t pointIndex : notGoodCandidates) {
216 if (points[minIndex](dim) > points[pointIndex](dim)) {
217 minIndex = pointIndex;
218 }
219 if (points[maxIndex](dim) < points[pointIndex](dim)) {
220 maxIndex = pointIndex;
221 }
222 }
223 notGoodCandidates.set(minIndex, false);
224 notGoodCandidates.set(maxIndex, false);
225 }
226 }
227 storm::storage::BitVector goodCandidates = ~notGoodCandidates;
228
229 // Found candidates. Now swap them to the front.
230 const uint_fast64_t numOfGoodCandidates = goodCandidates.getNumberOfSetBits();
231 for (auto goodCandidate : goodCandidates) {
232 if (goodCandidate >= numOfGoodCandidates) {
233 uint_fast64_t notGoodCandidate = *notGoodCandidates.begin();
234 assert(notGoodCandidate < numOfGoodCandidates);
235 std::swap(points[notGoodCandidate], points[goodCandidate]);
236 notGoodCandidates.set(notGoodCandidate, false);
237 }
238 }
239
240 storm::storage::geometry::SubsetEnumerator<std::vector<EigenVector>> subsetEnum(points.size(), dimension + 1, points, affineFilter);
241 if (subsetEnum.setToFirstSubset()) {
242 verticesOfInitialPolytope = subsetEnum.getCurrentSubset();
243 return true;
244 } else {
245 return false;
246 }
247}
248
249template<typename ValueType>
250std::vector<typename QuickHull<ValueType>::Facet> QuickHull<ValueType>::computeInitialFacets(std::vector<EigenVector> const& points,
251 std::vector<uint_fast64_t> const& verticesOfInitialPolytope,
252 EigenVector const& insidePoint) const {
253 const uint_fast64_t dimension = points.front().rows();
254 assert(verticesOfInitialPolytope.size() == dimension + 1);
255 std::vector<Facet> result;
256 result.reserve(dimension + 1);
257 storm::storage::geometry::SubsetEnumerator<> subsetEnum(verticesOfInitialPolytope.size(), dimension);
258 if (!subsetEnum.setToFirstSubset()) {
259 STORM_LOG_THROW(false, storm::exceptions::UnexpectedException, "Could not find an initial subset.");
260 }
261 do {
262 Facet newFacet;
263 // set the points that lie on the new facet
264 std::vector<uint_fast64_t> const& subset(subsetEnum.getCurrentSubset());
265 newFacet.points.reserve(subset.size());
266 for (uint_fast64_t i : subset) {
267 newFacet.points.push_back(verticesOfInitialPolytope[i]);
268 }
269 // neighbors: these are always the remaining facets
270 newFacet.neighbors.reserve(dimension);
271 for (uint_fast64_t i = 0; i < dimension + 1; ++i) {
272 if (i != result.size()) { // initFacets.size() will be the index of this new facet!
273 newFacet.neighbors.push_back(i);
274 }
275 }
276 // normal and offset:
277 computeNormalAndOffsetOfFacet(points, insidePoint, newFacet);
278
279 result.push_back(std::move(newFacet));
280 } while (subsetEnum.incrementSubset());
281 assert(result.size() == dimension + 1);
282 return result;
283}
284
285template<typename ValueType>
286void QuickHull<ValueType>::computeNormalAndOffsetOfFacet(std::vector<EigenVector> const& points, EigenVector const& insidePoint, Facet& facet) const {
287 const uint_fast64_t dimension = points.front().rows();
288 assert(facet.points.size() == dimension);
289 EigenVector const& refPoint = points[facet.points.back()];
290 EigenMatrix constraints(dimension - 1, dimension);
291 for (unsigned row = 0; row < dimension - 1; ++row) {
292 constraints.row(row) = (points[facet.points[row]] - refPoint);
293 }
294 facet.normal = constraints.fullPivLu().kernel().col(0);
295 facet.offset = facet.normal.dot(refPoint);
296
297 // invert the plane if the insidePoint is not contained in it
298 if (facet.normal.dot(insidePoint) > facet.offset) {
299 facet.normal = -facet.normal;
300 facet.offset = -facet.offset;
301 }
302}
303
304template<typename ValueType>
305void QuickHull<ValueType>::extendMesh(std::vector<EigenVector>& points, std::vector<Facet>& facets, storm::storage::BitVector& currentFacets,
306 EigenVector& insidePoint) const {
307 storm::storage::BitVector currentOutsidePoints(points.size(), true);
308 // Compute initial outside Sets
309 for (uint_fast64_t facetIndex : currentFacets) {
310 computeOutsideSetOfFacet(facets[facetIndex], currentOutsidePoints, points);
311 }
312
313 for (uint_fast64_t facetCount = currentFacets.getNextSetIndex(0); facetCount != currentFacets.size();
314 facetCount = currentFacets.getNextSetIndex(facetCount + 1)) {
315 // set all points to false to get rid of points that lie within the polytope after each iteration
316 currentOutsidePoints.clear();
317 // Find a facet with a non-empty outside set
318 if (!facets[facetCount].outsideSet.empty()) {
319 uint_fast64_t numberOfNewFacets = 0;
320 // Now we compute the enlarged mesh
321 uint_fast64_t farAwayPointIndex = facets[facetCount].maxOutsidePointIndex;
322 // get Visible set from maxOutsidePoint of the current facet
323 std::set<uint_fast64_t> visibleSet = getVisibleSet(facets, facetCount, points[farAwayPointIndex]);
324 std::set<uint_fast64_t> invisibleSet = getInvisibleNeighbors(facets, visibleSet);
325 for (auto invisFacetIt = invisibleSet.begin(); invisFacetIt != invisibleSet.end(); ++invisFacetIt) {
326 for (auto visFacetIt = visibleSet.begin(); visFacetIt != visibleSet.end(); ++visFacetIt) {
327 if (std::find(facets[*invisFacetIt].neighbors.begin(), facets[*invisFacetIt].neighbors.end(), *visFacetIt) !=
328 facets[*invisFacetIt].neighbors.end()) {
329 Facet newFacet;
330 // Set points of Facet
331 newFacet.points = getCommonPoints(facets[*invisFacetIt], facets[*visFacetIt]);
332 newFacet.points.push_back(farAwayPointIndex);
333 // replace old facet index by new facet index in the current neighbor
334 replaceFacetNeighbor(facets, *visFacetIt, facets.size(), *invisFacetIt);
335 newFacet.neighbors.push_back(*invisFacetIt);
336 // Compute the normal and the offset
337 computeNormalAndOffsetOfFacet(points, insidePoint, newFacet);
338
339 // add new facet
340 facets.push_back(newFacet);
341 currentFacets.resize(currentFacets.size() + 1, true);
342 // increase Number Of new Facets
343 ++numberOfNewFacets;
344 }
345 }
346 }
347
348 for (auto visibleFacet : visibleSet) {
349 for (uint_fast64_t outsidePoint : facets[visibleFacet].outsideSet) {
350 currentOutsidePoints.set(outsidePoint, true);
351 }
352 currentFacets.set(visibleFacet, false);
353 }
354 // compute new outside sets
355 for (uint_fast64_t facetIndex : currentFacets) {
356 computeOutsideSetOfFacet(facets[facetIndex], currentOutsidePoints, points);
357 }
358
359 // find neighbors in new facets
360 setNeighborhoodOfNewFacets(facets, facets.size() - numberOfNewFacets, points.front().rows());
361 }
362 }
363}
364
365template<typename ValueType>
366void QuickHull<ValueType>::getPolytopeFromMesh(std::vector<EigenVector> const& points, std::vector<Facet> const& facets,
367 storm::storage::BitVector const& currentFacets, bool generateRelevantVerticesAndVertexSets) {
369 for (auto facet : currentFacets) {
370 hyperplaneCollector.insert(std::move(facets[facet].normal), std::move(facets[facet].offset),
371 generateRelevantVerticesAndVertexSets ? &facets[facet].points : nullptr);
372 }
373
374 if (generateRelevantVerticesAndVertexSets) {
375 // Get the mapping from a hyperplane to the set of vertices that lie on that plane, erase the duplicates, and count for each vertex the number of
376 // hyperplanes on which that vertex lies
377 vertexSets = hyperplaneCollector.getIndexLists();
378 std::vector<uint_fast64_t> hyperplanesOnVertexCounter(points.size(), 0);
379 for (auto& vertexVector : vertexSets) {
380 std::set<uint_fast64_t> vertexSet;
381 for (auto const& i : vertexVector) {
382 if (vertexSet.insert(i).second) {
383 ++hyperplanesOnVertexCounter[i];
384 }
385 }
386 vertexVector.assign(vertexSet.begin(), vertexSet.end());
387 }
388 // Now, we can erase all vertices which do not lie on at least dimension hyperplanes.
389 // Note that the indices of the HyperplaneToVerticesMapping needs to be refreshed according to the new set of vertices
390 // Therefore, we additionally store the old indices for every vertex to be able to translate from old to new indices
391 std::unordered_map<EigenVector, std::vector<uint_fast64_t>> relevantVerticesMap;
392 relevantVerticesMap.reserve(points.size());
393 for (uint_fast64_t vertexIndex = 0; vertexIndex < hyperplanesOnVertexCounter.size(); ++vertexIndex) {
394 if ((Eigen::Index)hyperplanesOnVertexCounter[vertexIndex] >= points.front().rows()) {
395 auto mapEntry = relevantVerticesMap
396 .insert(typename std::unordered_map<EigenVector, std::vector<uint_fast64_t>>::value_type(points[vertexIndex],
397 std::vector<uint_fast64_t>()))
398 .first;
399 mapEntry->second.push_back(vertexIndex);
400 }
401 }
402 // Fill in the relevant vertices and create a translation map from old to new indices
403 std::vector<uint_fast64_t> oldToNewIndexMapping(points.size(), points.size()); // Initialize with some illegal value
404 relevantVertices.clear();
405 relevantVertices.reserve(relevantVerticesMap.size());
406 for (auto const& mapEntry : relevantVerticesMap) {
407 for (auto const& oldIndex : mapEntry.second) {
408 oldToNewIndexMapping[oldIndex] = relevantVertices.size();
409 }
410 relevantVertices.push_back(mapEntry.first);
411 }
412 // Actually translate and erase duplicates
413 for (auto& vertexVector : vertexSets) {
414 std::set<uint_fast64_t> vertexSet;
415 for (auto const& oldIndex : vertexVector) {
416 if ((Eigen::Index)hyperplanesOnVertexCounter[oldIndex] >= points.front().rows()) {
417 vertexSet.insert(oldToNewIndexMapping[oldIndex]);
418 }
419 }
420 vertexVector.assign(vertexSet.begin(), vertexSet.end());
421 }
422 } else {
423 relevantVertices.clear();
424 vertexSets.clear();
425 }
426 auto matrixVector = hyperplaneCollector.getCollectedHyperplanesAsMatrixVector();
427 resultMatrix = std::move(matrixVector.first);
428 resultVector = std::move(matrixVector.second);
429}
430
431template<typename ValueType>
432std::set<uint_fast64_t> QuickHull<ValueType>::getVisibleSet(std::vector<Facet> const& facets, uint_fast64_t const& startIndex, EigenVector const& point) const {
433 std::set<uint_fast64_t> facetsToCheck;
434 std::set<uint_fast64_t> facetsChecked;
435 std::set<uint_fast64_t> visibleSet;
436 facetsChecked.insert(startIndex);
437 visibleSet.insert(startIndex);
438 for (uint_fast64_t i = 0; i < facets[startIndex].neighbors.size(); ++i) {
439 facetsToCheck.insert(facets[startIndex].neighbors[i]);
440 }
441 while (!facetsToCheck.empty()) {
442 auto elementIt = facetsToCheck.begin();
443 if (point.dot(facets[*elementIt].normal) > facets[*elementIt].offset) {
444 visibleSet.insert(*elementIt);
445 for (uint_fast64_t i = 0; i < facets[*elementIt].neighbors.size(); ++i) {
446 if (facetsChecked.find(facets[*elementIt].neighbors[i]) == facetsChecked.end()) {
447 facetsToCheck.insert(facets[*elementIt].neighbors[i]);
448 }
449 }
450 }
451 facetsChecked.insert(*elementIt);
452 facetsToCheck.erase(elementIt);
453 }
454 return visibleSet;
455}
456
457template<typename ValueType>
458void QuickHull<ValueType>::setNeighborhoodOfNewFacets(std::vector<Facet>& facets, uint_fast64_t firstNewFacet, uint_fast64_t dimension) const {
459 for (uint_fast64_t currentFacet = firstNewFacet; currentFacet < facets.size(); ++currentFacet) {
460 for (uint_fast64_t otherFacet = currentFacet + 1; otherFacet < facets.size(); ++otherFacet) {
461 if (getCommonPoints(facets[currentFacet], facets[otherFacet]).size() >= dimension - 1) {
462 facets[currentFacet].neighbors.push_back(otherFacet);
463 facets[otherFacet].neighbors.push_back(currentFacet);
464 }
465 }
466 }
467}
468
469template<typename ValueType>
470void QuickHull<ValueType>::replaceFacetNeighbor(std::vector<Facet>& facets, uint_fast64_t oldFacetIndex, uint_fast64_t newFacetIndex,
471 uint_fast64_t neighborIndex) const {
472 auto facetInNeighborListIt = std::find(facets[neighborIndex].neighbors.begin(), facets[neighborIndex].neighbors.end(), oldFacetIndex);
473 if (facetInNeighborListIt != facets[neighborIndex].neighbors.end()) {
474 *facetInNeighborListIt = newFacetIndex;
475 }
476}
477
478template<typename ValueType>
479void QuickHull<ValueType>::computeOutsideSetOfFacet(Facet& facet, storm::storage::BitVector& currentOutsidePoints,
480 std::vector<EigenVector> const& points) const {
481 ValueType maxMultiplicationResult = facet.offset;
482 for (uint_fast64_t pointIndex : currentOutsidePoints) {
483 ValueType multiplicationResult = points[pointIndex].dot(facet.normal);
484 if (multiplicationResult > facet.offset) {
485 currentOutsidePoints.set(pointIndex, false); // we already know that the point lies outside so it can be ignored for future facets
486 facet.outsideSet.push_back(pointIndex);
487 if (multiplicationResult > maxMultiplicationResult) {
488 maxMultiplicationResult = multiplicationResult;
489 facet.maxOutsidePointIndex = pointIndex;
490 }
491 }
492 }
493}
494
495template<typename ValueType>
496std::vector<uint_fast64_t> QuickHull<ValueType>::getCommonPoints(Facet const& lhs, Facet const& rhs) const {
497 std::vector<uint_fast64_t> commonPoints;
498 for (uint_fast64_t lhsPoint : lhs.points) {
499 for (uint_fast64_t rhsPoint : rhs.points) {
500 if (lhsPoint == rhsPoint) {
501 commonPoints.push_back(lhsPoint);
502 }
503 }
504 }
505 return commonPoints;
506}
507
508template<typename ValueType>
509std::set<uint_fast64_t> QuickHull<ValueType>::getInvisibleNeighbors(std::vector<Facet>& facets, std::set<uint_fast64_t> const& visibleSet) const {
510 std::set<uint_fast64_t> invisibleNeighbors;
511 for (auto currentFacetIt = visibleSet.begin(); currentFacetIt != visibleSet.end(); ++currentFacetIt) {
512 for (uint_fast64_t currentNeighbor = 0; currentNeighbor < facets[*currentFacetIt].neighbors.size(); ++currentNeighbor) {
513 if (visibleSet.find(facets[*currentFacetIt].neighbors[currentNeighbor]) == visibleSet.end()) {
514 invisibleNeighbors.insert(facets[*currentFacetIt].neighbors[currentNeighbor]);
515 }
516 }
517 }
518 return invisibleNeighbors;
519}
520
521template class QuickHull<double>;
522template class QuickHull<storm::RationalNumber>;
523
524} // namespace geometry
525} // namespace storage
526} // namespace storm
A bit vector that is internally represented as a vector of 64-bit values.
Definition BitVector.h:18
uint_fast64_t getNextSetIndex(uint_fast64_t startingIndex) const
Retrieves the index of the bit that is the next bit set to true in the bit vector.
void resize(uint_fast64_t newLength, bool init=false)
Resizes the bit vector to hold the given new number of bits.
void set(uint_fast64_t index, bool value=true)
Sets the given truth value at the given index.
size_t size() const
Retrieves the number of bits this bit vector can store.
uint_fast64_t getNumberOfSetBits() const
Returns the number of bits that are set to true in this bit vector.
This class can be used to collect a set of hyperplanes (without duplicates).
bool insert(EigenVector const &normal, ValueType const &offset, std::vector< uint_fast64_t > const *indexList=nullptr)
std::vector< std::vector< uint_fast64_t > > getIndexLists() const
std::pair< EigenMatrix, EigenVector > getCollectedHyperplanesAsMatrixVector() const
Eigen::Matrix< ValueType, Eigen::Dynamic, Eigen::Dynamic > EigenMatrix
Definition QuickHull.h:16
std::vector< EigenVector > & getRelevantVertices()
Returns the set of vertices which are not redundant.
Definition QuickHull.cpp:67
Eigen::Matrix< ValueType, Eigen::Dynamic, 1 > EigenVector
Definition QuickHull.h:17
void generateHalfspacesFromPoints(std::vector< EigenVector > &points, bool generateRelevantVerticesAndVertexSets)
Definition QuickHull.cpp:21
std::vector< std::vector< std::uint_fast64_t > > & getVertexSets()
Returns for each hyperplane the set of vertices that lie on that hyperplane.
Definition QuickHull.cpp:72
This class can be used to enumerate all k-sized subsets of {0,...,n-1}.
#define STORM_LOG_WARN(message)
Definition logging.h:30
#define STORM_LOG_DEBUG(message)
Definition logging.h:23
#define STORM_LOG_ASSERT(cond, message)
Definition macros.h:11
#define STORM_LOG_THROW(cond, exception, message)
Definition macros.h:30
SFTBDDChecker::ValueType ValueType
void filterVectorInPlace(std::vector< Type > &v, storm::storage::BitVector const &filter)
Definition vector.h:1224
LabParser.cpp.
Definition cli.cpp:18