Skip to content
Snippets Groups Projects
Commit bad3dd33 authored by Ansgar Burchardt's avatar Ansgar Burchardt
Browse files

do not catch polymorphic types by value

This addresses several warnings from GCC 8 of the form

    ../dune/istl/test/matrixredisttest.cc:102:10: warning:
    catching polymorphic type ‘class Dune::ISTLError’ by value [-Wcatch-value=]
parent ae1ec50a
No related branches found
No related tags found
No related merge requests found
...@@ -578,7 +578,7 @@ namespace Dune ...@@ -578,7 +578,7 @@ namespace Dune
cont.sparsity[column].insert(i); cont.sparsity[column].insert(i);
} }
} }
catch(Dune::RangeError er) { catch(const Dune::RangeError&) {
// Entry not present in the new index set. Ignore! // Entry not present in the new index set. Ignore!
#ifdef DEBUG_REPART #ifdef DEBUG_REPART
typedef typename Container::LookupIndexSet GlobalLookup; typedef typename Container::LookupIndexSet GlobalLookup;
...@@ -593,7 +593,7 @@ namespace Dune ...@@ -593,7 +593,7 @@ namespace Dune
MPI_Comm_rank(MPI_COMM_WORLD,&rank); MPI_Comm_rank(MPI_COMM_WORLD,&rank);
std::cout<<rank<<cont.aggidxset<<std::endl; std::cout<<rank<<cont.aggidxset<<std::endl;
std::cout<<rank<<": row "<<i<<" (global="<<gi <<") not in index set for owner index "<<pi->global()<<std::endl; std::cout<<rank<<": row "<<i<<" (global="<<gi <<") not in index set for owner index "<<pi->global()<<std::endl;
throw er; throw;
} }
#endif #endif
} }
...@@ -657,7 +657,7 @@ namespace Dune ...@@ -657,7 +657,7 @@ namespace Dune
cont.matrix[i][column]=data.second; cont.matrix[i][column]=data.second;
} }
} }
catch(Dune::RangeError er) { catch(const Dune::RangeError&) {
// This an overlap row and might therefore lack some entries! // This an overlap row and might therefore lack some entries!
} }
...@@ -752,7 +752,7 @@ namespace Dune ...@@ -752,7 +752,7 @@ namespace Dune
{ {
try{ try{
newMatrix[col.index()][row.index()]; newMatrix[col.index()][row.index()];
}catch(Dune::ISTLError e) { }catch(const Dune::ISTLError&) {
std::cerr<<newComm.communicator().rank()<<": entry (" std::cerr<<newComm.communicator().rank()<<": entry ("
<<col.index()<<","<<row.index()<<") missing! for symmetry!"<<std::endl; <<col.index()<<","<<row.index()<<") missing! for symmetry!"<<std::endl;
ret=1; ret=1;
......
...@@ -6,41 +6,35 @@ using namespace Dune; ...@@ -6,41 +6,35 @@ using namespace Dune;
int main (int argc, char** argv) int main (int argc, char** argv)
{ {
try typedef BCRSMatrix<FieldMatrix<double,2,2> > Mat;
{
typedef BCRSMatrix<FieldMatrix<double,2,2> > Mat;
Mat A(1,1, Mat::random); Mat A(1,1, Mat::random);
A.setrowsize(0,1); A.setrowsize(0,1);
A.endrowsizes(); A.endrowsizes();
A.addindex(0, 0); A.addindex(0, 0);
A.endindices(); A.endindices();
A = 0; A = 0;
Mat B(2,2, Mat::random); Mat B(2,2, Mat::random);
B.setrowsize(0,2); B.setrowsize(0,2);
B.setrowsize(1,1); B.setrowsize(1,1);
B.endrowsizes(); B.endrowsizes();
B.addindex(0, 0); B.addindex(0, 0);
B.addindex(0, 1); B.addindex(0, 1);
B.addindex(1, 1); B.addindex(1, 1);
B.endindices(); B.endindices();
B = 0; B = 0;
B = A; B = A;
} catch(Exception e){ return 0;
std::cout<<e<<std::endl;
return 1;
}
return 0;
} }
...@@ -157,20 +157,15 @@ void testDoubleSetSize() ...@@ -157,20 +157,15 @@ void testDoubleSetSize()
int main() int main()
{ {
try{ Dune::TestSuite testSuite;
Dune::TestSuite testSuite;
Builder<Dune::BCRSMatrix<Dune::FieldMatrix<double,1,1> > > builder; Builder<Dune::BCRSMatrix<Dune::FieldMatrix<double,1,1> > > builder;
testSuite.subTest(builder.randomBuild(5,4)); testSuite.subTest(builder.randomBuild(5,4));
testSuite.subTest(builder.rowWiseBuild(5,4,13)); testSuite.subTest(builder.rowWiseBuild(5,4,13));
testSuite.subTest(builder.rowWiseBuild(5,4)); testSuite.subTest(builder.rowWiseBuild(5,4));
testDoubleSetSize(); testDoubleSetSize();
return testSuite.exit(); return testSuite.exit();
}catch(Dune::Exception e) {
std::cerr << e<<std::endl;
return 1;
}
} }
...@@ -99,7 +99,7 @@ int testRepart(int N, int coarsenTarget) ...@@ -99,7 +99,7 @@ int testRepart(int N, int coarsenTarget)
if(col.index()<=row.index()) if(col.index()<=row.index())
try{ try{
newMat[col.index()][row.index()]; newMat[col.index()][row.index()];
}catch(Dune::ISTLError e) { }catch(const Dune::ISTLError&) {
std::cerr<<coarseComm->communicator().rank()<<": entry (" std::cerr<<coarseComm->communicator().rank()<<": entry ("
<<col.index()<<","<<row.index()<<") missing!"<<std::endl; <<col.index()<<","<<row.index()<<") missing!"<<std::endl;
ret=1; ret=1;
......
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