[petsc-dev] The questions about MatSetValuesLocal

Matthew Knepley knepley at gmail.com
Fri Sep 27 17:48:48 CDT 2013


On Fri, Sep 27, 2013 at 11:00 AM, Lulu Liu <lulu.liu at kaust.edu.sa> wrote:

> Try to build analytical Jacobian with 4 blocks, I modify the code src/snes/examples/tutorials/ex28.c
> (1D problem, however, my problem is 2D). When building one block of
> Jacobian, I got the errors like this:
> 0]PETSC ERROR: --------------------- Error Message
> ------------------------------------
> [0]PETSC ERROR: Argument out of range!
> [0]PETSC ERROR: New nonzero at (28,58) caused a malloc!
> [0]PETSC ERROR: MatSetValues_SeqAIJ() line 353 in src/mat/impls/aij/seq/
> aij.c
> [0]PETSC ERROR: MatSetValues() line 1124 in src/mat/interface/matrix.c
> [0]PETSC ERROR: MatSetValuesLocal() line 2012 in src
> /mat/interface/matrix.c
> [0]PETSC ERROR: FormJacobianLocal_P() line 556 in /Users/liul
> /Desktop/FieldSplitcase2.c
> [0]PETSC ERROR: FormJacobian_All() line 847 in /Users/liul
> /Desktop/FieldSplitcase2.c
> [0]PETSC ERROR: SNESComputeJacobian() line 2227 in src/snes/interface/snes
> .c
> [0]PETSC ERROR: SNESSolve_NEWTONLS() line 233 in src/snes/impls/ls/ls.c
> [0]PETSC ERROR: SNESSolve() line 3789 in src/snes/interface/snes.c
> [0]PETSC ERROR: main() line 1062 in /Users/liul/Desktop/FieldSplitcase2.c
> application called MPI_Abort(MPI_COMM_WORLD, 63) - process 0
> [unset]: aborting job:
> application called MPI_Abort(MPI_COMM_WORLD, 63) - process 0
>

1) Obviously your index calculation is off somewhere

2) This is exactly why we introduced the FormFunctionLocal() where you get
multi-dimensional arrays
     and out of range access is obvious. Why are you not using that?

  Matt


> I think that I am not call MatSetValuesLocal correctly, right? How to fix
> it? Thanks!
> for (j=yints; j<yinte; j++) {
>     for (i=xints; i<xinte; i++) {
>  //        row.i = i; row.j = j;
>            row = i+j*info->xm;
>
>            if (i==0 || j==0 || i==mx-1 || j==my-1) {
>           num = 0; v_l=0.0; v_r=0.0; v_t=0.0; v_b=0.0;
>           if (j!=0) {
>             v=-1;v_b=1;
>             J=row-info->xm;
>             ierr = MatSetValuesLocal(Bpp,1,&row,1,&J,&v,INSERT_VALUES);
> CHKERRQ(ierr);
>           }
>           if (i!=0) {
>             v = -1; v_l=1;
>             J=row-1;
>             ierr = MatSetValuesLocal(Bpp,1,&row,1,&J,&v,INSERT_VALUES);
> CHKERRQ(ierr);
>
>           }
>           if (i!=mx-1) {
>             v = -1;v_r=1;
>             J=row+1;
>             ierr = MatSetValuesLocal(Bpp,1,&row,1,&J,&v,INSERT_VALUES);
> CHKERRQ(ierr);
>           }
>           if (j!=my-1) {
>             v = -1;v_t=1;
>              J=row+info->xm;
>             ierr = MatSetValuesLocal(Bpp,1,&row,1,&J,&v,INSERT_VALUES);
> CHKERRQ(ierr);
>           }
>           v   = v_l+v_r+v_t+v_b;
>
>          if(i==0&&j==0){ v+=2;}
>           J=row;
>          ierr = MatSetValuesLocal(Bpp,1,&row,1,&J,&v,INSERT_VALUES);
> CHKERRQ(ierr);
>          }else{
>
>           PetscInt cols[]={row-info->xm,row-1,row+1,row+info->xm,row};
>           PetscScalar vals[]= {-1,-1,-1,-1,4};
>           ierr = MatSetValuesLocal(Bpp,1,&row,5,cols,vals,INSERT_VALUES);
> CHKERRQ(ierr);
>
>          }
>
>     }
>   }
>
>
>
>
> ------------------------------
> This message and its contents, including attachments are intended solely
> for the original recipient. If you are not the intended recipient or have
> received this message in error, please notify me immediately and delete
> this message from your computer system. Any unauthorized use or
> distribution is prohibited. Please consider the environment before printing
> this email.




-- 
What most experimenters take for granted before they begin their
experiments is infinitely more interesting than any results to which their
experiments lead.
-- Norbert Wiener
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-dev/attachments/20130927/283af083/attachment.html>


More information about the petsc-dev mailing list