[petsc-users] ERROR: Argument out of range!; Local index XX too large XX (max) at 0! with MatSetValuesStencil

Bishesh Khanal bisheshkh at gmail.com
Thu Jul 4 07:32:15 CDT 2013


Hi all,
I'm trying to use DMCreateMatrix and MatStencil to fill up a matrix that
results from a finite difference discretization of a PDE (objective is to
solve the resulting linear system in the form of Ax=b).

However there are some errors and I'm not sure about the issue!
A short demo:
A 2D mXn grid, with two variables at each node (dof=2),
so the resulting A would be 2mn X 2mn.
Let's say the two variables are vx and vy, and the two associated
equations discretized are x-eq and y-eq.

Here is the relevant part of the code: (mat1 variable for the A matrix):
    PetscInt m = 10, n=10;
     ierr = DMDACreate2d(PETSC_COMM_WORLD, DMDA_BOUNDARY_NONE,
DMDA_BOUNDARY_NONE,DMDA_STENCIL_BOX,m,n,

PETSC_DECIDE,PETSC_DECIDE,2,1,NULL,NULL,&da);CHKERRQ(ierr);

    Mat mat1;
    ierr = DMCreateMatrix(da,MATMPIAIJ,&mat1); CHKERRQ(ierr);
    MatStencil row, col[4]; //let's test 4 non-zeros in one row.
    PetscScalar val[4];
    PetscScalar coeff = 2.; //just a constant coefficient for testing.
    PetscInt i,j;
    DMDALocalInfo info;
    ierr = DMDAGetLocalInfo(da,&info); CHKERRQ(ierr);
//Now fill up the matrix:
for(i = info.ys; i < info.ys+info.ym; ++i){
        for(j = info.xs; j < info.xs+info.xm; ++j){
            row.i = i;  row.j = j;   //one node at a time
            if (j == 0 || j == info.mx-1){ //left and right borders
                //vx(i,j) = 0;
                row.c = 0;
                col[0].c = 0;   val[0] = 1;
                col[0].i = i;   col[0].j = j;
                ierr =
MatSetValuesStencil(mat1,1,&row,1,col,val,INSERT_VALUES);CHKERRQ(ierr);

                //vy: //vy(i,j) - c*vy(i,j+-1) = 0;
                row.c = 1;
                col[0].c = 1;   val[1] = coeff;
                col[1].c = 1;
                col[1].i = i;
                if(j == 0) //vy(i,j) - c*vy(i,j+1) = 0;
                    col[1].j = j+1;
                else //vy(i,j) - c*vy(i,j-1) = 0;
                    col[1].j = j-1;

                ierr =
MatSetValuesStencil(mat1,1,&row,2,col,val,INSERT_VALUES);CHKERRQ(ierr);
            }
            else if (i == 0 || i == info.my-1){  //top and bottom borders
                //vx: vx(i,j) - c* vx(i+-1,j) = 0;
                row.c = 0;
                col[0].c = 0;   val[0] = 1;
                col[0].i = i;   col[0].j = j;
                col[1].c = 0;   val[1] = coeff;
                if (i == 0) //vx(i,j) - c*vx(i+1,j) = 0;
                    col[1].i = i+1;
                else //vx(i,j) - c*vx(i-1,j) = 0;
                    col[1].i = i-1;
                col[1].j = j;
                ierr =
MatSetValuesStencil(mat1,1,&row,2,col,val,INSERT_VALUES);CHKERRQ(ierr);

                //vy(i,j) = 0;
                row.c = 1;
                col[0].c = 1;
                ierr =
MatSetValuesStencil(mat1,1,&row,1,col,val,INSERT_VALUES);CHKERRQ(ierr);
            }
            else { //Interior points:
                row.c = 0;//x-eq
                col[0].c = 0;   val[0] = 2*coeff;
                col[0].i = i;   col[0].j = j+1;
                col[1].c = 0;   val[1] = -val[0] - coeff;
                col[1].i = i;   col[1].j = j;
                col[2].c = 1;   val[2] = 4*coeff;
                col[2].i = i+1; col[2].j = j;
                col[3].c = 1;   val[3] = 4*coeff;
                col[3].i = i;   col[3].j = j;
                ierr =
MatSetValuesStencil(mat1,1,&row,4,col,val,INSERT_VALUES);CHKERRQ(ierr);

                row.c = 1; //y-eq
                col[0].c = 1;   val[0] = 2*coeff;
                col[0].i = i;   col[0].j = j+1;
                col[1].c = 1;   val[1] = -val[0] - coeff;
                col[1].i = i;   col[1].j = j;
                col[2].c = 0;   val[2] = 4*coeff;
                col[2].i = i+1; col[2].j = j;
                col[3].c = 0;   val[3] = 4*coeff;
                col[3].i = i;   col[3].j = j;
                ierr =
MatSetValuesStencil(mat1,1,&row,4,col,val,INSERT_VALUES);CHKERRQ(ierr);
            }
        }
    }

    MatAssemblyBegin(mat1,MAT_FINAL_ASSEMBLY);
    MatAssemblyEnd(mat1,MAT_FINAL_ASSEMBLY);

However I get the following error when run with 2 processors.
[0]PETSC ERROR: --------------------- Error Message
------------------------------------
[0]PETSC ERROR: Argument out of range!
[0]PETSC ERROR: Local index 120 too large 120 (max) at 0!
[0]PETSC ERROR:
------------------------------------------------------------------------
[0]PETSC ERROR: Petsc Release Version 3.4.1, Jun, 10, 2013
[0]PETSC ERROR: See docs/changes/index.html for recent updates.
[0]PETSC ERROR: See docs/faq.html for hints about trouble shooting.
[0]PETSC ERROR: See docs/index.html for manual pages.
[0]PETSC ERROR:
------------------------------------------------------------------------
[0]PETSC ERROR: examples/FDdmda_test on a arch-linux2-cxx-debug named
edwards by bkhanal Thu Jul  4 14:20:11 2013
[0]PETSC ERROR: Libraries linked from
/home/bkhanal/Documents/softwares/petsc-3.4.1/arch-linux2-cxx-debug/lib
[0]PETSC ERROR: Configure run at Wed Jun 19 11:04:51 2013
[0]PETSC ERROR: Configure options --with-cc=gcc --with-fc=g77
--with-cxx=g++ --download-f-blas-lapack=1 --download-mpich=1
-with-clanguage=cxx --download-hypre=1
[0]PETSC ERROR:
------------------------------------------------------------------------
[0]PETSC ERROR: ISLocalToGlobalMappingApply() line 444 in
/home/bkhanal/Documents/softwares/petsc-3.4.1/src/vec/is/utils/isltog.c
[0]PETSC ERROR: MatSetValuesLocal() line 1967 in
/home/bkhanal/Documents/softwares/petsc-3.4.1/src/mat/interface/matrix.c
[0]PETSC ERROR: MatSetValuesStencil() line 1339 in
/home/bkhanal/Documents/softwares/petsc-3.4.1/src/mat/interface/matrix.c
[0]PETSC ERROR: main() line 75 in
"unknowndirectory/"/user/bkhanal/home/works/cmake_tuts/petsc_test/examples/FDdmda_test.cxx
application called MPI_Abort(MPI_COMM_WORLD, 63) - process 0
[cli_0]: aborting job:
application called MPI_Abort(MPI_COMM_WORLD, 63) - process 0

===================================================================================
=   BAD TERMINATION OF ONE OF YOUR APPLICATION PROCESSES
=   EXIT CODE: 63
=   CLEANING UP REMAINING PROCESSES
=   YOU CAN IGNORE THE BELOW CLEANUP MESSAGES
===================================================================================

I do not understand how the index (i,j) are out of range when set to
corresponding fields in row and col variables. What could be the possible
problem ? And any suggestions on the way to debug this sort of issues ?

Thanks,
Bishesh
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20130704/22d5f3ca/attachment.html>


More information about the petsc-users mailing list