Skip to content
Snippets Groups Projects
Commit 531c2ae1 authored by Max Kahnt's avatar Max Kahnt
Browse files

Import generic vector truncation from dune-solvers.

Use sparse range for iteration to obtain usage on static
data structures for free in addition to former functionality.
parent 2b4816bf
No related branches found
No related tags found
No related merge requests found
Pipeline #
...@@ -24,6 +24,11 @@ void writeBinary(std::ostream& s, const Vector& v) { ...@@ -24,6 +24,11 @@ void writeBinary(std::ostream& s, const Vector& v) {
Helper<Vector>::writeBinary(s, v); Helper<Vector>::writeBinary(s, v);
} }
//! Truncate vector by given bit set
template <class Vector, class BitVector>
void truncate(Vector& v, const BitVector& tr) {
Helper<Vector>::truncate(v, tr);
}
template <class Vector> template <class Vector>
struct Helper { struct Helper {
...@@ -32,6 +37,12 @@ struct Helper { ...@@ -32,6 +37,12 @@ struct Helper {
Generic::writeBinary(s, vi); Generic::writeBinary(s, vi);
} }
template <class BitVector>
static void truncate(Vector& v, const BitVector& tr) {
sparseRangeFor(v, [&tr](auto&& vi, auto&& i) {
Generic::truncate(vi, tr[i]);
});
}
}; };
template <class Field, int n> template <class Field, int n>
...@@ -43,10 +54,16 @@ struct Helper<FieldVector<Field, n>> { ...@@ -43,10 +54,16 @@ struct Helper<FieldVector<Field, n>> {
s.write(reinterpret_cast<const char*>(&vi), sizeof(Field)); s.write(reinterpret_cast<const char*>(&vi), sizeof(Field));
} }
template<class BitVector>
static void truncate(Vector& v, const BitVector& tr) {
for (auto it = v.begin(), end = v.end(); it != end; ++it)
if(tr[it.index()])
*it = 0;
}
};
} // end namespace Generic } // end namespace Generic
} // end namespace MatrixVector } // end namespace MatrixVector
} // end namespace Dune } // end namespace Dune
#endif // DUNE_MATRIX_VECTOR_GENERICVECTORTOOLS_HH #endif // DUNE_MATRIX_VECTOR_GENERICVECTORTOOLS_HH
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Please register or to comment