Skip to content

#202 exists() function for BCRSMatrix causes segmentation fault

Metadata

Property Value
Reported by Patrick Leidenberger (mailings.pl@gmail.com)
Reported at Nov 21, 2006 12:16
Type Bug Report
Version Git (pre2.4) [autotools]
Operating System Unspecified / All
Closed by Markus Blatt (markus@dr-blatt.de)
Closed at Nov 21, 2006 15:37
Closed in version Unknown
Resolution Fixed
Comment Fixed in revision 686.

Added check whether row has nonzero size. |

Description

The exists () function for hierarchical nested BCRS matrices causes an segmentation fault if a test to a non existing entry is made.

Below the modified dune-istl/istl/tutorial/example.cc and it's result to demonstrate the problem.

Regards Patrick

// start with including some headers #include "config.h"

#include // for input/output to shell #include // for input/output to files #include // STL vector class #include

#include<math.h> // Yes, we do some math here #include<stdio.h> // There is nothing better than sprintf #include<sys/times.h> // for timing measurements

#include<dune/istl/istlexception.hh> #include<dune/istl/basearray.hh> #include<dune/common/fvector.hh> #include<dune/common/fmatrix.hh> #include<dune/istl/bvector.hh> #include<dune/istl/vbvector.hh> #include<dune/istl/bcrsmatrix.hh> #include<dune/istl/io.hh> #include<dune/istl/gsetc.hh> #include<dune/istl/ilu.hh> #include<dune/istl/operators.hh> #include<dune/istl/solvers.hh> #include<dune/istl/preconditioners.hh> #include<dune/istl/scalarproducts.hh>

// a simple stop watch class Timer { public: Timer () { struct tms buf; cstart = times(&buf); }

void start () { struct tms buf; cstart = times(&buf); }

double stop () { struct tms buf; cend = times(&buf); return ((double)(cend-cstart))/100.0; }

double gettime () { return ((double)(cend-cstart))/100.0; }

private: clock_t cstart,cend;
};

void test_BCRSMatrix () { const int N=5, M =4; typedef double ValueType; typedef Dune::BCRSMatrix<Dune::FieldMatrix<ValueType,1,1> > EmpMatrixTopoElemOrder; typedef Dune::BCRSMatrix EmpMatrixTopoElem;

Dune::FieldMatrix<ValueType,1,1> F; F = 42.0;

EmpMatrixTopoElemOrder E(N,N,N*N,EmpMatrixTopoElemOrder::row_wise); for (EmpMatrixTopoElemOrder::CreateIterator i=E.createbegin(); i!=E.createend(); ++i) for (int j=0; j<N; ++j) i.insert(j);

for (EmpMatrixTopoElemOrder::RowIterator i=E.begin(); i!=E.end(); ++i) for (EmpMatrixTopoElemOrder::ColIterator j=(*i).begin(); j!=(*i).end(); ++j) *j = F;

std::cout << "============================================." << std::endl; printmatrix(std::cout,E,"E","row",10,2); std::cout << "=============================================."<< std::endl;

for(unsigned int row = 0; row < N; ++row) for(unsigned int column = 0; column < N; ++column) if (E.exists(row,column)) std::cout << "exists: row = " << row << " column = " << column << std::endl;

EmpMatrixTopoElem D(M,M,2,EmpMatrixTopoElem::row_wise); for (EmpMatrixTopoElem::CreateIterator i=D.createbegin(); i!=D.createend(); ++i) for (int j=0; j<M; ++j) { if ((i.index() == 1) && (j == 1)) { i.insert(j); } }

for (EmpMatrixTopoElem::RowIterator i=D.begin(); i!=D.end(); ++i) for (EmpMatrixTopoElem::ColIterator j=(*i).begin(); j!=(*i).end(); ++j) *j = E;

std::cout << "============================================." << std::endl; printmatrix(std::cout,D,"D","row",10,2); std::cout << "=============================================."<< std::endl;

for(unsigned int row = 0; row < M; ++row) for(unsigned int column = 0; column < M; ++column) if (D.exists(row,column)) std::cout << "exists: row = " << row << " column = " << column << std::endl;

}

int main (int argc , char ** argv) { try { test_BCRSMatrix(); } catch (Dune::ISTLError& error) { std::cout << error << std::endl; } catch (Dune::Exception& error) { std::cout << error << std::endl; } catch (const std::bad_alloc& e) { std::cout << "memory exhausted" << std::endl; } catch (...) { std::cout << "unknown exception caught" << std::endl; }

return 0; }

============================================. E [n=5,m=5,rowdim=5,coldim=5] row 0 4.20e+01 4.20e+01 4.20e+01 4.20e+01 4.20e+01 row 1 4.20e+01 4.20e+01 4.20e+01 4.20e+01 4.20e+01 row 2 4.20e+01 4.20e+01 4.20e+01 4.20e+01 4.20e+01 row 3 4.20e+01 4.20e+01 4.20e+01 4.20e+01 4.20e+01 row 4 4.20e+01 4.20e+01 4.20e+01 4.20e+01 4.20e+01 =============================================. exists: row = 0 column = 0 exists: row = 0 column = 1 exists: row = 0 column = 2 exists: row = 0 column = 3 exists: row = 0 column = 4 exists: row = 1 column = 0 exists: row = 1 column = 1 exists: row = 1 column = 2 exists: row = 1 column = 3 exists: row = 1 column = 4 exists: row = 2 column = 0 exists: row = 2 column = 1 exists: row = 2 column = 2 exists: row = 2 column = 3 exists: row = 2 column = 4 exists: row = 3 column = 0 exists: row = 3 column = 1 exists: row = 3 column = 2 exists: row = 3 column = 3 exists: row = 3 column = 4 exists: row = 4 column = 0 exists: row = 4 column = 1 exists: row = 4 column = 2 exists: row = 4 column = 3 exists: row = 4 column = 4 ============================================. D [n=4,m=4,rowdim=5,coldim=5] row 0 4.20e+01 4.20e+01 4.20e+01 4.20e+01 4.20e+01 row 1 4.20e+01 4.20e+01 4.20e+01 4.20e+01 4.20e+01 row 2 4.20e+01 4.20e+01 4.20e+01 4.20e+01 4.20e+01 row 3 4.20e+01 4.20e+01 4.20e+01 4.20e+01 4.20e+01 row 4 4.20e+01 4.20e+01 4.20e+01 4.20e+01 4.20e+01 =============================================. Segmentation fault