Skip to content
Snippets Groups Projects
Commit 2538d8c8 authored by Markus Blatt's avatar Markus Blatt
Browse files

If we read a char from a file we should have written it as char (not

an enum), too

[[Imported from SVN: r1483]]
parent ca58516e
No related branches found
No related tags found
No related merge requests found
......@@ -945,7 +945,7 @@ namespace Dune
for(Iterator iter = comm.indexSet().begin();
iter != comm.indexSet().end(); ++iter) {
file << iter->global()<<" "<<(std::size_t)iter->local()<<" "
<<iter->local().attribute()<<" "<<(int)iter->local().isPublic()<<std::endl;
<<(char)iter->local().attribute()<<" "<<(int)iter->local().isPublic()<<std::endl;
}
// Store neighbour information for efficient remote indices setup.
file<<"neighbours:";
......
......@@ -4,8 +4,10 @@
#include "mpi.h"
#include <dune/istl/io.hh>
#include <dune/istl/bvector.hh>
#include <dune/istl/schwarz.hh>
#include <dune/common/fmatrix.hh>
#include <dune/common/fvector.hh>
#include <dune/common/float_cmp.hh>
#include <dune/istl/paamg/test/anisotropic.hh>
#include <dune/common/timer.hh>
#include <dune/istl/matrixmarket.hh>
......@@ -40,7 +42,7 @@ int main(int argc, char** argv)
BCRSMat mat = setupAnisotropic2d<BS,double>(N, comm.indexSet(), comm.communicator(), &n, .011);
BVector bv = BVector(mat.N());
BVector bv(mat.N()), cv(mat.N());
typedef BVector::iterator VIter;
int i=0;
......@@ -52,18 +54,42 @@ int main(int argc, char** argv)
comm.remoteIndices().rebuild<false>();
comm.copyOwnerToAll(bv,bv);
Dune::OverlappingSchwarzOperator<BCRSMat,BVector,BVector,Communication> op(mat, comm);
op.apply(bv, cv);
storeMatrixMarket(mat, std::string("testmat"), comm);
storeMatrixMarket(bv, std::string("testvec"), comm, false);
BCRSMat mat1;
BVector bv1;
BVector bv1,cv1;
Communication comm1(MPI_COMM_WORLD);
loadMatrixMarket(mat1, std::string("testmat"), comm1);
loadMatrixMarket(bv1, std::string("testvec"), comm1, false);
int ret=0;
if(mat.N()!=mat1.N() || mat.M()!=mat1.M())
{
++ret;
std::cerr<<"matrix sizes do not match"<<std::endl;
}
typedef BCRSMat::const_iterator RowIterator;
typedef BCRSMat::ConstColIterator ColIterator;
for(RowIterator row=mat.begin(), row1=mat1.begin(); row!=mat.end(); ++row, ++row1)
for(ColIterator col=row->begin(), col1=row1->begin(); col!= row->end(); ++col, ++col1)
{
if(col.index()!=col1.index()) {
std::cerr <<"Column indices do not match"<<std::endl;
++ret;
}
if(!Dune::FloatCmp::eq(*col, *col1)) {
std::cerr <<"Matrix entries do not match"<<std::endl;
++ret;
}
}
for(VIter entry=bv.begin(), entry1=bv1.begin(); bv.end() != entry; ++entry, ++entry1)
if(*entry!=*entry1)
{
......@@ -71,6 +97,15 @@ int main(int argc, char** argv)
++ret;
}
cv1.resize(mat1.M());
Dune::OverlappingSchwarzOperator<BCRSMat,BVector,BVector,Communication> op1(mat1, comm1);
op1.apply(bv1, cv1);
for(VIter entry=cv.begin(), entry1=cv1.begin(); cv.end() != entry; ++entry, ++entry1)
if(*entry!=*entry1)
{
std::cerr<<"computed vectors do not match"<<std::endl;
++ret;
}
if(comm1.indexSet()!=comm.indexSet())
{
std::cerr<<"written and read idxset do not match"<<std::endl;
......
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Finish editing this message first!
Please register or to comment