Browse Source

some optimizations in explicit model building

main
dehnert 7 years ago
parent
commit
1f9e2967c8
  1. 54
      src/storm/builder/jit/Distribution.cpp
  2. 10
      src/storm/builder/jit/Distribution.h
  3. 33
      src/storm/builder/jit/DistributionEntry.cpp
  4. 8
      src/storm/builder/jit/DistributionEntry.h
  5. 2
      src/storm/builder/jit/ModelComponentsBuilder.cpp
  6. 2
      src/storm/generator/NextStateGenerator.cpp
  7. 47
      src/storm/generator/PrismNextStateGenerator.cpp
  8. 14
      src/storm/storage/BitVector.cpp
  9. 5
      src/storm/storage/BitVector.h
  10. 95
      src/storm/storage/BitVectorHashMap.cpp
  11. 23
      src/storm/storage/BitVectorHashMap.h

54
src/storm/builder/jit/Distribution.cpp

@ -2,6 +2,8 @@
#include "storm/adapters/RationalFunctionAdapter.h"
#include "storm/storage/BitVector.h"
namespace storm {
namespace builder {
namespace jit {
@ -11,16 +13,48 @@ namespace storm {
// Intentionally left empty.
}
template <typename IndexType, typename ValueType>
Distribution<IndexType, ValueType>::Distribution(Distribution<IndexType, ValueType> const& other) {
this->storage = other.storage;
this->compressed = other.compressed;
}
template <typename IndexType, typename ValueType>
Distribution<IndexType, ValueType>::Distribution(Distribution<IndexType, ValueType>&& other) {
this->storage = std::move(other.storage);
this->compressed = other.compressed;
other.compressed = true;
}
template <typename IndexType, typename ValueType>
Distribution<IndexType, ValueType>& Distribution<IndexType, ValueType>::operator=(Distribution<IndexType, ValueType> const& other) {
if (this != &other) {
this->storage = other.storage;
this->compressed = other.compressed;
}
return *this;
}
template <typename IndexType, typename ValueType>
Distribution<IndexType, ValueType>& Distribution<IndexType, ValueType>::operator=(Distribution<IndexType, ValueType>&& other) {
if (this != &other) {
this->storage = std::move(other.storage);
this->compressed = other.compressed;
other.compressed = true;
}
return *this;
}
template <typename IndexType, typename ValueType>
void Distribution<IndexType, ValueType>::add(DistributionEntry<IndexType, ValueType> const& entry) {
storage.push_back(entry);
compressed &= storage.back().getIndex() < entry.getIndex();
compressed &= storage.back().getState() < entry.getState();
}
template <typename IndexType, typename ValueType>
void Distribution<IndexType, ValueType>::add(IndexType const& index, ValueType const& value) {
storage.emplace_back(index, value);
compressed &= storage.back().getIndex() < index;
compressed &= storage.back().getState() < index;
}
template <typename IndexType, typename ValueType>
@ -34,7 +68,7 @@ namespace storm {
if (!compressed) {
std::sort(storage.begin(), storage.end(),
[] (DistributionEntry<IndexType, ValueType> const& a, DistributionEntry<IndexType, ValueType> const& b) {
return a.getIndex() < b.getIndex();
return a.getState() < b.getState();
}
);
@ -45,7 +79,7 @@ namespace storm {
if (first != last) {
auto result = first;
while (++first != last) {
if (!(result->getIndex() == first->getIndex())) {
if (!(result->getState() == first->getState())) {
if (++result != first) {
*result = std::move(*first);
}
@ -68,6 +102,12 @@ namespace storm {
}
}
template <typename IndexType, typename ValueType>
void Distribution<IndexType, ValueType>::clear() {
this->storage.clear();
this->compressed = true;
}
template <typename IndexType, typename ValueType>
typename Distribution<IndexType, ValueType>::ContainerType::iterator Distribution<IndexType, ValueType>::begin() {
return storage.begin();
@ -91,7 +131,11 @@ namespace storm {
template class Distribution<uint32_t, double>;
template class Distribution<uint32_t, storm::RationalNumber>;
template class Distribution<uint32_t, storm::RationalFunction>;
template class Distribution<storm::storage::BitVector, double>;
template class Distribution<storm::storage::BitVector, storm::RationalNumber>;
template class Distribution<storm::storage::BitVector, storm::RationalFunction>;
}
}
}

10
src/storm/builder/jit/Distribution.h

@ -15,6 +15,11 @@ namespace storm {
Distribution();
Distribution(Distribution const&);
Distribution(Distribution&&);
Distribution& operator=(Distribution const&);
Distribution& operator=(Distribution&&);
/*!
* Adds the given entry to the distribution.
*/
@ -41,6 +46,11 @@ namespace storm {
*/
void divide(ValueType const& value);
/*!
* Clears this distribution.
*/
void clear();
/*!
* Access to iterators over the entries of the distribution. Note that there may be multiple entries for
* the same index. Also, no order is guaranteed. After a call to compress, the order is guaranteed to be

33
src/storm/builder/jit/DistributionEntry.cpp

@ -2,37 +2,39 @@
#include "storm/adapters/RationalFunctionAdapter.h"
#include "storm/storage/BitVector.h"
namespace storm {
namespace builder {
namespace jit {
template <typename IndexType, typename ValueType>
DistributionEntry<IndexType, ValueType>::DistributionEntry() : index(0), value(0) {
template <typename StateType, typename ValueType>
DistributionEntry<StateType, ValueType>::DistributionEntry() : state(0), value(0) {
// Intentionally left empty.
}
template <typename IndexType, typename ValueType>
DistributionEntry<IndexType, ValueType>::DistributionEntry(IndexType const& index, ValueType const& value) : index(index), value(value) {
template <typename StateType, typename ValueType>
DistributionEntry<StateType, ValueType>::DistributionEntry(StateType const& state, ValueType const& value) : state(state), value(value) {
// Intentionally left empty.
}
template <typename IndexType, typename ValueType>
IndexType const& DistributionEntry<IndexType, ValueType>::getIndex() const {
return index;
template <typename StateType, typename ValueType>
StateType const& DistributionEntry<StateType, ValueType>::getState() const {
return state;
}
template <typename IndexType, typename ValueType>
ValueType const& DistributionEntry<IndexType, ValueType>::getValue() const {
template <typename StateType, typename ValueType>
ValueType const& DistributionEntry<StateType, ValueType>::getValue() const {
return value;
}
template <typename IndexType, typename ValueType>
void DistributionEntry<IndexType, ValueType>::addToValue(ValueType const& value) {
template <typename StateType, typename ValueType>
void DistributionEntry<StateType, ValueType>::addToValue(ValueType const& value) {
this->value += value;
}
template <typename IndexType, typename ValueType>
void DistributionEntry<IndexType, ValueType>::divide(ValueType const& value) {
template <typename StateType, typename ValueType>
void DistributionEntry<StateType, ValueType>::divide(ValueType const& value) {
this->value /= value;
}
@ -40,6 +42,11 @@ namespace storm {
template class DistributionEntry<uint32_t, storm::RationalNumber>;
template class DistributionEntry<uint32_t, storm::RationalFunction>;
template class DistributionEntry<storm::storage::BitVector, double>;
template class DistributionEntry<storm::storage::BitVector, storm::RationalNumber>;
template class DistributionEntry<storm::storage::BitVector, storm::RationalFunction>;
}
}
}

8
src/storm/builder/jit/DistributionEntry.h

@ -4,20 +4,20 @@ namespace storm {
namespace builder {
namespace jit {
template <typename IndexType, typename ValueType>
template <typename StateType, typename ValueType>
class DistributionEntry {
public:
DistributionEntry();
DistributionEntry(IndexType const& index, ValueType const& value);
DistributionEntry(StateType const& state, ValueType const& value);
IndexType const& getIndex() const;
StateType const& getState() const;
ValueType const& getValue() const;
void addToValue(ValueType const& value);
void divide(ValueType const& value);
private:
IndexType index;
StateType state;
ValueType value;
};

2
src/storm/builder/jit/ModelComponentsBuilder.cpp

@ -58,7 +58,7 @@ namespace storm {
for (auto const& choice : behaviour.getChoices()) {
// Add the elements to the transition matrix.
for (auto const& element : choice.getDistribution()) {
transitionMatrixBuilder->addNextValue(currentRow, element.getIndex(), element.getValue());
transitionMatrixBuilder->addNextValue(currentRow, element.getState(), element.getValue());
}
// Add state-action reward entries.

2
src/storm/generator/NextStateGenerator.cpp

@ -74,7 +74,7 @@ namespace storm {
result.addLabel(label.first);
}
storm::storage::BitVectorHashMap<StateType> const& states = stateStorage.stateToId;
auto const& states = stateStorage.stateToId;
for (auto const& stateIndexPair : states) {
unpackStateIntoEvaluator(stateIndexPair.first, variableInformation, *this->evaluator);

47
src/storm/generator/PrismNextStateGenerator.cpp

@ -8,6 +8,8 @@
#include "storm/storage/expressions/SimpleValuation.h"
#include "storm/storage/sparse/PrismChoiceOrigins.h"
#include "storm/builder/jit/Distribution.h"
#include "storm/solver/SmtSolver.h"
#include "storm/utility/constants.h"
@ -449,6 +451,9 @@ namespace storm {
std::vector<Choice<ValueType>> PrismNextStateGenerator<ValueType, StateType>::getLabeledChoices(CompressedState const& state, StateToIdCallback stateToIdCallback) {
std::vector<Choice<ValueType>> result;
storm::builder::jit::Distribution<CompressedState, ValueType> currentDistribution;
storm::builder::jit::Distribution<CompressedState, ValueType> nextDistribution;
for (uint_fast64_t actionIndex : program.getSynchronizingActionIndices()) {
boost::optional<std::vector<std::vector<std::reference_wrapper<storm::prism::Command const>>>> optionalActiveCommandLists = getActiveCommandsByActionIndex(actionIndex);
@ -465,40 +470,32 @@ namespace storm {
// As long as there is one feasible combination of commands, keep on expanding it.
bool done = false;
while (!done) {
boost::container::flat_map<CompressedState, ValueType>* currentTargetStates = new boost::container::flat_map<CompressedState, ValueType>();
boost::container::flat_map<CompressedState, ValueType>* newTargetStates = new boost::container::flat_map<CompressedState, ValueType>();
currentTargetStates->emplace(state, storm::utility::one<ValueType>());
currentDistribution.clear();
nextDistribution.clear();
currentDistribution.add(state, storm::utility::one<ValueType>());
for (uint_fast64_t i = 0; i < iteratorList.size(); ++i) {
storm::prism::Command const& command = *iteratorList[i];
for (uint_fast64_t j = 0; j < command.getNumberOfUpdates(); ++j) {
storm::prism::Update const& update = command.getUpdate(j);
for (auto const& stateProbabilityPair : *currentTargetStates) {
ValueType probability = stateProbabilityPair.second * this->evaluator->asRational(update.getLikelihoodExpression());
for (auto const& stateProbability : currentDistribution) {
ValueType probability = stateProbability.getValue() * this->evaluator->asRational(update.getLikelihoodExpression());
if (!storm::utility::isZero<ValueType>(probability)) {
// Compute the new state under the current update and add it to the set of new target states.
CompressedState newTargetState = applyUpdate(stateProbabilityPair.first, update);
// If the new state was already found as a successor state, update the probability
// and otherwise insert it.
auto targetStateIt = newTargetStates->find(newTargetState);
if (targetStateIt != newTargetStates->end()) {
targetStateIt->second += probability;
} else {
newTargetStates->emplace(newTargetState, probability);
}
CompressedState newTargetState = applyUpdate(stateProbability.getState(), update);
nextDistribution.add(newTargetState, probability);
}
}
}
nextDistribution.compress();
// If there is one more command to come, shift the target states one time step back.
if (i < iteratorList.size() - 1) {
delete currentTargetStates;
currentTargetStates = newTargetStates;
newTargetStates = new boost::container::flat_map<CompressedState, ValueType>();
currentDistribution = std::move(nextDistribution);
}
}
@ -524,11 +521,11 @@ namespace storm {
// Add the probabilities/rates to the newly created choice.
ValueType probabilitySum = storm::utility::zero<ValueType>();
for (auto const& stateProbabilityPair : *newTargetStates) {
StateType actualIndex = stateToIdCallback(stateProbabilityPair.first);
choice.addProbability(actualIndex, stateProbabilityPair.second);
for (auto const& stateProbability : nextDistribution) {
StateType actualIndex = stateToIdCallback(stateProbability.getState());
choice.addProbability(actualIndex, stateProbability.getValue());
if (this->options.isExplorationChecksSet()) {
probabilitySum += stateProbabilityPair.second;
probabilitySum += stateProbability.getValue();
}
}
@ -550,10 +547,6 @@ namespace storm {
choice.addReward(stateActionRewardValue);
}
// Dispose of the temporary maps.
delete currentTargetStates;
delete newTargetStates;
// Now, check whether there is one more command combination to consider.
bool movedIterator = false;
for (int_fast64_t j = iteratorList.size() - 1; !movedIterator && j >= 0; --j) {

14
src/storm/storage/BitVector.cpp

@ -999,6 +999,20 @@ namespace storm {
out << std::endl;
}
std::size_t FNV1aBitVectorHash::operator()(storm::storage::BitVector const& bv) const {
std::size_t seed = 14695981039346656037ull;
uint64_t prime = 1099511628211ull;
unsigned char const* it = reinterpret_cast<unsigned char const*>(bv.buckets);
unsigned char const* ite = it + 8 * bv.bucketCount();
for (; it != ite; ++it) {
seed = (*it ^ seed) * prime;
}
return seed;
}
// All necessary explicit template instantiations.
template BitVector::BitVector(uint_fast64_t length, std::vector<uint_fast64_t>::iterator begin, std::vector<uint_fast64_t>::iterator end);
template BitVector::BitVector(uint_fast64_t length, std::vector<uint_fast64_t>::const_iterator begin, std::vector<uint_fast64_t>::const_iterator end);

5
src/storm/storage/BitVector.h

@ -501,6 +501,7 @@ namespace storm {
friend std::ostream& operator<<(std::ostream& out, BitVector const& bitVector);
friend struct std::hash<storm::storage::BitVector>;
friend struct FNV1aBitVectorHash;
private:
/*!
@ -570,6 +571,10 @@ namespace storm {
static const uint_fast64_t mod64mask = (1 << 6) - 1;
};
struct FNV1aBitVectorHash {
std::size_t operator()(storm::storage::BitVector const& bv) const;
};
} // namespace storage
} // namespace storm

95
src/storm/storage/BitVectorHashMap.cpp

@ -9,9 +9,6 @@
namespace storm {
namespace storage {
template<class ValueType, class Hash>
const std::vector<std::size_t> BitVectorHashMap<ValueType, Hash>::sizes = {5, 13, 31, 79, 163, 277, 499, 1021, 2029, 3989, 8059, 16001, 32099, 64301, 127921, 256499, 511111, 1024901, 2048003, 4096891, 8192411, 15485863, 32142191, 64285127, 128572517, 257148523, 514299959, 1028599919, 2057199839, 4114399697, 8228799419, 16457598791, 32915197603, 65830395223};
template<class ValueType, class Hash>
BitVectorHashMap<ValueType, Hash>::BitVectorHashMapIterator::BitVectorHashMapIterator(BitVectorHashMap const& map, BitVector::const_iterator indexIt) : map(map), indexIt(indexIt) {
// Intentionally left empty.
@ -45,19 +42,20 @@ namespace storm {
}
template<class ValueType, class Hash>
BitVectorHashMap<ValueType, Hash>::BitVectorHashMap(uint64_t bucketSize, uint64_t initialSize, double loadFactor) : loadFactor(loadFactor), bucketSize(bucketSize), numberOfElements(0) {
BitVectorHashMap<ValueType, Hash>::BitVectorHashMap(uint64_t bucketSize, uint64_t initialSize, double loadFactor) : loadFactor(loadFactor), bucketSize(bucketSize), currentSize(1), numberOfElements(0) {
STORM_LOG_ASSERT(bucketSize % 64 == 0, "Bucket size must be a multiple of 64.");
currentSizeIterator = std::find_if(sizes.begin(), sizes.end(), [=] (uint64_t value) { return value > initialSize; } );
while (initialSize > 0) {
++currentSize;
initialSize >>= 1;
}
// Create the underlying containers.
buckets = storm::storage::BitVector(bucketSize * *currentSizeIterator);
occupied = storm::storage::BitVector(*currentSizeIterator);
values = std::vector<ValueType>(*currentSizeIterator);
buckets = storm::storage::BitVector(bucketSize * (1ull << currentSize));
occupied = storm::storage::BitVector(1ull << currentSize);
values = std::vector<ValueType>(1ull << currentSize);
#ifndef NDEBUG
for (uint64_t i = 0; i < sizes.size() - 1; ++i) {
STORM_LOG_ASSERT(sizes[i] < sizes[i + 1], "Expected stricly increasing sizes.");
}
numberOfInsertions = 0;
numberOfInsertionProbingSteps = 0;
numberOfFinds = 0;
@ -77,60 +75,32 @@ namespace storm {
template<class ValueType, class Hash>
std::size_t BitVectorHashMap<ValueType, Hash>::capacity() const {
return *currentSizeIterator;
return 1ull << currentSize;
}
template<class ValueType, class Hash>
void BitVectorHashMap<ValueType, Hash>::increaseSize() {
++currentSizeIterator;
STORM_LOG_ASSERT(currentSizeIterator != sizes.end(), "Hash map became to big.");
++currentSize;
#ifndef NDEBUG
STORM_LOG_TRACE("Increasing size of hash map from " << *(currentSizeIterator - 1) << " to " << *currentSizeIterator << ". Stats: " << numberOfFinds << " finds (avg. " << (numberOfFindProbingSteps / static_cast<double>(numberOfFinds)) << " probing steps), " << numberOfInsertions << " insertions (avg. " << (numberOfInsertionProbingSteps / static_cast<double>(numberOfInsertions)) << " probing steps).");
STORM_LOG_TRACE("Increasing size of hash map from " << (1ull << (currentSize - 1)) << " to " << (1ull << currentSize) << ". Stats: " << numberOfFinds << " finds (avg. " << (numberOfFindProbingSteps / static_cast<double>(numberOfFinds)) << " probing steps), " << numberOfInsertions << " insertions (avg. " << (numberOfInsertionProbingSteps / static_cast<double>(numberOfInsertions)) << " probing steps).");
#else
STORM_LOG_TRACE("Increasing size of hash map from " << *(currentSizeIterator - 1) << " to " << *currentSizeIterator << ".");
STORM_LOG_TRACE("Increasing size of hash map from " << (1ull << (currentSize - 1)) << " to " << (1ull << currentSize) << ".");
#endif
// Create new containers and swap them with the old ones.
numberOfElements = 0;
storm::storage::BitVector oldBuckets(bucketSize * *currentSizeIterator);
storm::storage::BitVector oldBuckets(bucketSize * (1ull << currentSize));
std::swap(oldBuckets, buckets);
storm::storage::BitVector oldOccupied = storm::storage::BitVector(*currentSizeIterator);
storm::storage::BitVector oldOccupied = storm::storage::BitVector(1ull << currentSize);
std::swap(oldOccupied, occupied);
std::vector<ValueType> oldValues = std::vector<ValueType>(*currentSizeIterator);
std::vector<ValueType> oldValues = std::vector<ValueType>(1ull << currentSize);
std::swap(oldValues, values);
// Now iterate through the elements and reinsert them in the new storage.
bool fail = false;
for (auto bucketIndex : oldOccupied) {
fail = !this->insertWithoutIncreasingSize(oldBuckets.get(bucketIndex * bucketSize, bucketSize), oldValues[bucketIndex]);
// If we failed to insert just one element, we have to redo the procedure with a larger container size.
if (fail) {
break;
}
}
uint_fast64_t failCount = 0;
while (fail) {
++failCount;
STORM_LOG_THROW(failCount < 2, storm::exceptions::InternalException, "Increasing size failed too often.");
++currentSizeIterator;
STORM_LOG_THROW(currentSizeIterator != sizes.end(), storm::exceptions::InternalException, "Hash map became to big.");
numberOfElements = 0;
buckets = storm::storage::BitVector(bucketSize * *currentSizeIterator);
occupied = storm::storage::BitVector(*currentSizeIterator);
values = std::vector<ValueType>(*currentSizeIterator);
for (auto bucketIndex : oldOccupied) {
fail = !this->insertWithoutIncreasingSize(oldBuckets.get(bucketIndex * bucketSize, bucketSize), oldValues[bucketIndex]);
// If we failed to insert just one element, we have to redo the procedure with a larger container size.
if (fail) {
break;
}
}
STORM_LOG_ASSERT(!fail, "Expected to be able to insert all elements.");
}
}
@ -166,7 +136,7 @@ namespace storm {
template<class ValueType, class Hash>
std::pair<ValueType, std::size_t> BitVectorHashMap<ValueType, Hash>::findOrAddAndGetBucket(storm::storage::BitVector const& key, ValueType const& value) {
// If the load of the map is too high, we increase the size.
if (numberOfElements >= loadFactor * *currentSizeIterator) {
if (numberOfElements >= loadFactor * (1ull << currentSize)) {
this->increaseSize();
}
@ -187,7 +157,7 @@ namespace storm {
template<class ValueType, class Hash>
std::size_t BitVectorHashMap<ValueType, Hash>::setOrAddAndGetBucket(storm::storage::BitVector const& key, ValueType const& value) {
// If the load of the map is too high, we increase the size.
if (numberOfElements >= loadFactor * *currentSizeIterator) {
if (numberOfElements >= loadFactor * (1ull << currentSize)) {
this->increaseSize();
}
@ -230,17 +200,12 @@ namespace storm {
return const_iterator(*this, occupied.end());
}
template<class ValueType, class Hash>
uint_fast64_t BitVectorHashMap<ValueType, Hash>::getNextBucketInProbingSequence(uint_fast64_t, uint_fast64_t currentValue, uint_fast64_t step) const {
return (currentValue + step + step*step) % *currentSizeIterator;
}
template<class ValueType, class Hash>
std::pair<bool, std::size_t> BitVectorHashMap<ValueType, Hash>::findBucket(storm::storage::BitVector const& key) const {
#ifndef NDEBUG
++numberOfFinds;
#endif
uint_fast64_t initialHash = hasher(key) % *currentSizeIterator;
uint_fast64_t initialHash = hasher(key) % (1ull << currentSize);
uint_fast64_t bucket = initialHash;
uint_fast64_t i = 0;
@ -252,8 +217,11 @@ namespace storm {
if (buckets.matches(bucket * bucketSize, key)) {
return std::make_pair(true, bucket);
}
bucket = getNextBucketInProbingSequence(initialHash, bucket, i);
bucket += 1;
if (bucket == (1ull << currentSize)) {
bucket = 0;
}
if (bucket == initialHash) {
return std::make_pair(false, bucket);
}
@ -268,7 +236,7 @@ namespace storm {
#ifndef NDEBUG
++numberOfInsertions;
#endif
uint_fast64_t initialHash = hasher(key) % *currentSizeIterator;
uint_fast64_t initialHash = hasher(key) % (1ull << currentSize);
uint_fast64_t bucket = initialHash;
uint64_t i = 0;
@ -280,12 +248,15 @@ namespace storm {
if (buckets.matches(bucket * bucketSize, key)) {
return std::make_tuple(true, bucket, false);
}
bucket = getNextBucketInProbingSequence(initialHash, bucket, i);
bucket += 1;
if (bucket == (1ull << currentSize)) {
bucket = 0;
}
if (bucket == initialHash) {
if (increaseStorage) {
this->increaseSize();
bucket = initialHash = hasher(key) % *currentSizeIterator;
bucket = initialHash = hasher(key) % (1ull << currentSize);
} else {
return std::make_tuple(false, bucket, true);
}

23
src/storm/storage/BitVectorHashMap.h

@ -14,7 +14,8 @@ namespace storm {
* queries and insertions are supported. Also, the keys must be bit vectors with a length that is a multiple of
* 64.
*/
template<typename ValueType, typename Hash = std::hash<storm::storage::BitVector>>
// template<typename ValueType, typename Hash = std::hash<storm::storage::BitVector>>
template<typename ValueType, typename Hash = FNV1aBitVectorHash>
class BitVectorHashMap {
public:
class BitVectorHashMapIterator {
@ -57,6 +58,11 @@ namespace storm {
*/
BitVectorHashMap(uint64_t bucketSize = 64, uint64_t initialSize = 1000, double loadFactor = 0.75);
BitVectorHashMap(BitVectorHashMap const&) = default;
BitVectorHashMap(BitVectorHashMap&&) = default;
BitVectorHashMap& operator=(BitVectorHashMap const&) = default;
BitVectorHashMap& operator=(BitVectorHashMap&&) = default;
/*!
* Searches for the given key in the map. If it is found, the mapped-to value is returned. Otherwise, the
* key is inserted with the given value.
@ -211,19 +217,14 @@ namespace storm {
*/
void increaseSize();
/*!
* Computes the next bucket in the probing sequence.
*/
uint_fast64_t getNextBucketInProbingSequence(uint_fast64_t initialValue, uint_fast64_t currentValue, uint_fast64_t step) const;
// The load factor determining when the size of the map is increased.
double loadFactor;
// The size of one bucket.
uint64_t bucketSize;
// The number of buckets.
std::size_t numberOfBuckets;
// The number of buckets is 2^currentSize.
uint64_t currentSize;
// The buckets that hold the elements of the map.
storm::storage::BitVector buckets;
@ -237,15 +238,9 @@ namespace storm {
// The number of elements in this map.
std::size_t numberOfElements;
// An iterator to a value in the static sizes table.
std::vector<std::size_t>::const_iterator currentSizeIterator;
// Functor object that are used to perform the actual hashing.
Hash hasher;
// A static table that produces the next possible size of the hash table.
static const std::vector<std::size_t> sizes;
#ifndef NDEBUG
// Some performance metrics.
mutable uint64_t numberOfInsertions;

Loading…
Cancel
Save