diff --git a/dune/istl/io.hh b/dune/istl/io.hh index 61bcf8d58e2deecdcf20c8e270aa3465dd3178b5..37842597732ea0ff0c14915e3137725c0deb676e 100644 --- a/dune/istl/io.hh +++ b/dune/istl/io.hh @@ -570,6 +570,73 @@ namespace Dune { outStream.precision(oldPrecision); } + // Recursively print all the blocks + template<class V> + void writeVectorToMatlabHelper (const V& v, std::ostream& stream) + { + for (const auto& entry : v) + writeVectorToMatlabHelper(entry, stream); + } + + // Recursively print all the blocks -- specialization for FieldVector + template<class K, int n> + void writeVectorToMatlabHelper (const FieldVector<K,n>& v, std::ostream& s) + { + for (const auto& entry : v) + { + s << entry << std::endl; + } + } + + // Recursively print all the blocks -- specialization for std::vector + template<class K> + void writeVectorToMatlabHelper (const std::vector<K>& v, std::ostream& s) + { + for (const auto& entry : v) + { + s << entry << std::endl; + } + } + + // Recursively print all the blocks -- specialization for std::array + template<class K, unsigned long n> + void writeVectorToMatlabHelper (const std::array<K,n>& v, std::ostream& s) + { + for (const auto& entry : v) + { + s << entry << std::endl; + } + } + + /** + * \brief Writes vectors in a Matlab-readable format + * + * \code + * #include <dune/istl/io.hh> + * \endcode + * This routine writes the argument block vector to a file with the name given + * by the filename argument. The file format is ASCII, with no header, and + * a single data column. Such a file can be read from Matlab using the + * command + * \code + * new_vec = load('filename'); + * \endcode + * \param vector reference to vector to be printed to output file + * \param filename filename of output file + * \param outputPrecision (number of digits) which is used to write the output file + */ + template <class VectorType> + void writeVectorToMatlab(const VectorType& vector, + const std::string& filename, int outputPrecision = 18) + { + std::ofstream outStream(filename.c_str()); + int oldPrecision = outStream.precision(); + outStream.precision(outputPrecision); + + writeVectorToMatlabHelper(vector, outStream); + outStream.precision(oldPrecision); + } + /** @} end documentation */ } // namespace Dune diff --git a/dune/istl/test/iotest.cc b/dune/istl/test/iotest.cc index d2d63601399ad227a2f40456c7175dc75d4a6179..7d81c558dc9f009616fc9a891a536f9b3d7d86db 100644 --- a/dune/istl/test/iotest.cc +++ b/dune/istl/test/iotest.cc @@ -23,6 +23,19 @@ void testWriteMatrixToMatlab(BlockType b=BlockType(0.0)) writeMatrixToMatlabHelper(A, 0, 0, std::cout); } +/* uses the writeVectorToMatlab method, filled with dummy data */ +template <class VectorType> +void testWriteVectorToMatlab() +{ + VectorType v; + for (unsigned int i = 0; i < v.size(); ++i) + { + v[i] = i; + } + + Dune::writeVectorToMatlabHelper(v, std::cout); +} + int main(int argc, char** argv) { /* testing the writeMatrixToMatlabHelper method for BlockType=FieldMatrix with different field_types */ @@ -50,4 +63,21 @@ int main(int argc, char** argv) testWriteMatrixToMatlab<Dune::ScaledIdentityMatrix<std::complex<double>,1> >(Dune::ScaledIdentityMatrix<std::complex<double>,1>(std::complex<double>(0,1))); testWriteMatrixToMatlab<Dune::DiagonalMatrix<std::complex<double>,2> >(Dune::DiagonalMatrix<std::complex<double>,2>(std::complex<double>(0,1))); testWriteMatrixToMatlab<Dune::ScaledIdentityMatrix<std::complex<double>,2> >(Dune::ScaledIdentityMatrix<std::complex<double>,2>(std::complex<double>(0,1))); + + + /* testing the test writeMatrixToMatlabHelper for FieldVector */ + testWriteVectorToMatlab<Dune::FieldVector<double,1> >(); + testWriteVectorToMatlab<Dune::FieldVector<float,5> >(); + testWriteVectorToMatlab<Dune::FieldVector<std::complex<double>,5> >(); + /* testing the test writeMatrixToMatlabHelper for BlockVector */ + Dune::BlockVector<Dune::FieldVector<double,3> > v1 = {{1.0, 2.0, 3.0}}; + Dune::writeVectorToMatlabHelper(v1, std::cout); + Dune::BlockVector<Dune::BlockVector<Dune::FieldVector<double,1> > > v2 = {{1.0, 2.0}, {3.0, 4.0, 5.0}, {6.0}}; + Dune::writeVectorToMatlabHelper(v2, std::cout); + /* testing the test writeMatrixToMatlabHelper for STL containers */ + testWriteVectorToMatlab<std::array<double,5> >(); + std::vector<double> v3; + v3.push_back(1); + v3.push_back(2); + Dune::writeVectorToMatlabHelper(v3, std::cout); }