[Trilinos-Users] Received elements must be stored after all x local elements

Daniele Avitabile d.avitabile at gmail.com
Tue Sep 16 02:34:11 MDT 2008


Mike,

I believe Daniele's problem is somehow connected to another other problem I
have reported to you lately. In the CVS repository, in

nox/examples/epetra/LOCA_ContinuationManager/LinearSystem.C

I build a very simple Epetra_CrsMatrix, even though I span the rows in a
different way, but I receive an error message similar to Daniele's one.

Best.

Daniele Avitabile

On Tue, Sep 16, 2008 at 9:00 AM, Daniele Bettella <sunriis at gmail.com> wrote:

> Ok, this is the code that faults
>
>    int NumGlobalElements = A.nrows();
>
>    // creates map
>    Epetra_Map gmap(NumGlobalElements, 0, Comm);
>    int NumMyElements = gmap.NumMyElements();
>    int * MyGlobalElements = gmap.MyGlobalElements();
>    Epetra_CrsMatrix Atr(Copy, gmap, A.ncols());
>
>    //adjusting the iterators for each process
>    typename gmm::linalg_traits<MAT>::const_row_iterator itm =
> mat_row_const_begin(A);
>    for(int i=0; i<MyGlobalElements[0]; i++) itm++;
>    typename gmm::linalg_traits<MAT>::const_row_iterator itme = itm;
>    for(int i=0; i<NumMyElements; i++) itme++;
>
>    //we access the row of the matrix
>    typedef typename gmm::linalg_traits<MAT>::const_sub_row_type rowtype;
>    rowtype row;
>    int nnzRow;
>    double* Values;
>    int* Indices;
>    int icount=0;
>
>    //now we copy the matrix by row
>    while(itm != itme){
>      row = gmm::linalg_traits<MAT>::row(itm);
>      nnzRow = gmm::nnz(row);
>      Values = new double[nnzRow];
>      Indices = new int[nnzRow];
>      typename gmm::linalg_traits<rowtype>::const_iterator itv =
> vect_const_begin(row);
>      typename gmm::linalg_traits<rowtype>::const_iterator itve =
> vect_const_end(row);
>      int j = 0;
>      while(itv != itve){
>        Values[j] = *itv;
>        Indices[j] = itv.index();
>        ++itv;
>        ++j;
>      }
>      // and we insert values and indices in epetra matrix
>      Atr.InsertGlobalValues(MyGlobalElements[icount], nnzRow, Values,
> Indices);
>      ++itm;
>      ++icount;
>    }
>    Atr.FillComplete();
>
>
> This is the code that copies the matrix; as said before, the error
> comes out when initializing an aztecoo object:
>    Epetra_LinearProblem problem(&Atr, &Xtr, &Btr);
>    AztecOO solver(problem);
>
> where Xtr and Btr can be any vector.
> Again, the error doesn't appear if the matrix has a main diagonal of
> nonzeros
>
> This is what I could come out with, hope it helps!
>
> Daniele
>
>
> On Wed, Sep 10, 2008 at 9:20 PM, Heroux, Michael A <maherou at sandia.gov>
> wrote:
> > Daniele,
> >
> > Do you have a simple example program that you could send?
> >
> > Mike
> >
> >
> > On 9/10/08 2:23 AM, "Daniele Bettella" <sunriis at gmail.com> wrote:
> >
> > Hi, I'm having a problem with AztecOO:
> > when initializing the aztecoo object the error written in the subject
> pops
> > out.
> > I'm translating a sparse matrix from another format by hand; let's say I
> > have the following nonzeros in a 3x3 matrix:
> > A(0,0), A(0,2), A(1,0), A(2,1), A(2,2)
> >
> > what I do is create a map with the number of rows (3)
> >
> > Epetra_Map gmap(nrows, 0, comm)
> >
> > then I do
> >
> > Epetra_CrsMatrix newA(Copy, gmap, ncols)
> >
> > I then insert the elements by hand, cycling through the rows. Using MPI
> and
> > more then one process I divide the matrix to be translated in partitions
> of
> > rows, assigned to each process as per map definition.
> >
> > With the example matrix written above I get the error
> > AZ_extract_comm_info: Received elements must be stored after all 2 local
> > elements
> >
> >  but if I add an element in position A(1,1), obtaining a matrix with a
> main
> > diagonal of nonzeros, the error disappears.
> >
> > Any clue would be greatly appreciated, thanks in advance.
> >
> > Daniele
> >
> >
> >
>
> _______________________________________________
> Trilinos-Users mailing list
> Trilinos-Users at software.sandia.gov
> http://software.sandia.gov/mailman/listinfo/trilinos-users
>
>
-------------- next part --------------
An HTML attachment was scrubbed...
URL: http://software.sandia.gov/mailman/private/trilinos-users/attachments/20080916/0c3e9b35/attachment-0001.html 


More information about the Trilinos-Users mailing list