[petsc-users] Problem of solving simple example in SNES in parallel

Gaurish Telang gaurish108 at gmail.com
Fri Mar 22 10:59:26 CDT 2013


Hello,

I am trying to learn SNES to solve some non-linear equations.

Consider the following non-linear equation which has solution (x,y,z) =
(2,3,4)

x^2 - y^2 + z^2 = 11
xy  + y^2 - 3z   = 3
x   - xz    + yz   = 6

For solving this, I am able to get the following simple code I wrote to
work on 1 processor, but it seems to be failing for 2 or more processors.

The errors I am getting when I run with two processors (mpirun  -n 2
./a.out -snes_monitor )
Can anybody suggestions to fix this issue? It looks like a very trivial
problem, but I am not sure what is wrong.

(The code I wrote is after the error)

%---------------------------------------------------
The error is
[0]PETSC ERROR: PetscMallocValidate: error detected at
SNESComputeFunction() line 1889 in
/home/gaurish108/Software/petsc-3.3-p5/src/snes/interface/snes.c
[0]PETSC ERROR: Memory [id=0(16)] at address 0xa1933f0 is corrupted
(probably write past end of array)
[0]PETSC ERROR: Memory originally allocated in VecCreate_MPI_Private() line
193 in /home/gaurish108/Software/petsc-3.3-p5/src/vec/vec/impls/mpi/pbvec.c
[1]PETSC ERROR: PetscMallocValidate: error detected at
SNESComputeFunction() line 1889 in
/home/gaurish108/Software/petsc-3.3-p5/src/snes/interface/snes.c
[1]PETSC ERROR: Memory [id=0(16)] at address 0xa04c350 is corrupted
(probably write past end of array)
[1]PETSC ERROR: Memory originally allocated in VecCreate_MPI_Private() line
193 in /home/gaurish108/Software/petsc-3.3-p5/src/vec/vec/impls/mpi/pbvec.c
[1]PETSC ERROR: --------------------- Error Message
------------------------------------
[1]PETSC ERROR: Memory corruption!
[1]PETSC ERROR:  !
[1]PETSC ERROR:
------------------------------------------------------------------------
[1]PETSC ERROR: Petsc Release Version 3.3.0, Patch 5, Sat Dec  1 15:10:41
CST 2012
[1]PETSC ERROR: See docs/changes/index.html for recent updates.
[1]PETSC ERROR: See docs/faq.html for hints about trouble shooting.
[0]PETSC ERROR: --------------------- Error Message
------------------------------------
[0]PETSC ERROR: Memory corruption!
[0]PETSC ERROR:  !
[0]PETSC ERROR:
------------------------------------------------------------------------
[0]PETSC ERROR: Petsc Release Version 3.3.0, Patch 5, Sat Dec  1 15:10:41
CST 2012
[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: [1]PETSC ERROR: See docs/index.html for manual pages.
[1]PETSC ERROR:
------------------------------------------------------------------------
[1]PETSC ERROR: ./gaurish on a linux-gnu named gaurish108 by gaurish108 Fri
Mar 22 11:58:01 2013
[1]PETSC ERROR: Libraries linked from
/home/gaurish108/Software/petsc-3.3-p5/linux-gnu-c-yes-debug-refblas/lib
[1]PETSC ERROR: Configure run at Tue Feb 19 20:22:20 2013
[1]PETSC ERROR: Configure options --download-f-blas-lapack=1
[1]PETSC ERROR:
------------------------------------------------------------------------
[1]PETSC ERROR: PetscMallocValidate() line 150 in
/home/gaurish108/Software/petsc-3.3-p5/src/sys/memory/mtr.c
[1]PETSC ERROR: SNESComputeFunction() line 1889 in
/home/gaurish108/Software/petsc-3.3-p5/src/snes/interface/snes.c
[1]PETSC ERROR: SNESSolve_LS() line 161 in
/home/gaurish108/Software/petsc-3.3-p5/src/snes/impls/ls/ls.c
[1]PETSC ERROR: SNESSolve() line 3536 in
/home/gaurish108/Software/petsc-3.3-p5/src/snes/interface/snes.c
[1]PETSC ERROR: main() line 69 in gaurish.c
See docs/index.html for manual pages.
[0]PETSC ERROR:
------------------------------------------------------------------------
[0]PETSC ERROR: ./gaurish on a linux-gnu named gaurish108 by gaurish108 Fri
Mar 22 11:58:01 2013
[0]PETSC ERROR: Libraries linked from
/home/gaurish108/Software/petsc-3.3-p5/linux-gnu-c-yes-debug-refblas/lib
[0]PETSC ERROR: Configure run at Tue Feb 19 20:22:20 2013
[0]PETSC ERROR: Configure options --download-f-blas-lapack=1
[0]PETSC ERROR:
------------------------------------------------------------------------
[0]PETSC ERROR: PetscMallocValidate() line 150 in
/home/gaurish108/Software/petsc-3.3-p5/src/sys/memory/mtr.c
[0]PETSC ERROR: SNESComputeFunction() line 1889 in
/home/gaurish108/Software/petsc-3.3-p5/src/snes/interface/snes.c
[0]PETSC ERROR: SNESSolve_LS() line 161 in
/home/gaurish108/Software/petsc-3.3-p5/src/snes/impls/ls/ls.c
[0]PETSC ERROR: SNESSolve() line 3536 in
/home/gaurish108/Software/petsc-3.3-p5/src/snes/interface/snes.c
[0]PETSC ERROR: main() line 69 in gaurish.c
--------------------------------------------------------------------------
MPI_ABORT was invoked on rank 1 in communicator MPI_COMM_WORLD
with errorcode 78.

NOTE: invoking MPI_ABORT causes Open MPI to kill all MPI processes.
You may or may not see output from other processes, depending on
exactly when Open MPI kills them.
--------------------------------------------------------------------------
--------------------------------------------------------------------------
mpirun has exited due to process rank 1 with PID 17933 on
node gaurish108 exiting without calling "finalize". This may
have caused other processes in the application to be
terminated by signals sent by mpirun (as reported here).
--------------------------------------------------------------------------
[gaurish108:17931] 1 more process has sent help message help-mpi-api.txt /
mpi-abort
[gaurish108:17931] Set MCA parameter "orte_base_help_aggregate" to 0 to see
all help / error messages



%------------------------------------------------------
The code is

static char help[] = "Shows solve a nonlinear equation using Newton
iteration method.\n\n";
#include <petscsnes.h>
#include <math.h>
#include <petscsys.h>
#include <stdio.h>
#include <petsctime.h>
/* User-defined routines*/
extern PetscErrorCode FormJacobian(SNES,Vec,Mat*,Mat*,MatStructure*,void*);
extern PetscErrorCode FormFunction(SNES,Vec,Vec,void*);
extern PetscErrorCode FormInitialGuess(Vec);
extern PetscErrorCode Monitor(SNES,PetscInt,PetscReal,void *);
#undef __FUNCT__
#define __FUNCT__ "main"
int main(int argc,char **argv)
{ SNES           snes;              /* SNES context */
  Vec            x,r,F;             /* vectors */
  Mat            J;                 /* Jacobian matrix */
  PetscErrorCode ierr;
  PetscInt       its,n = 3,i,maxit,maxf;
  PetscMPIInt    size;
  PetscScalar    iniRHS[3];

  PetscInitialize(&argc,&argv,(char *)0,help);
  ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr);
  ierr = PetscOptionsGetInt(PETSC_NULL,"-n",&n,PETSC_NULL);CHKERRQ(ierr);

//-------------------------------------------------------------------------------
  ierr = SNESCreate(PETSC_COMM_WORLD,&snes);CHKERRQ(ierr);

  ierr = VecCreate(PETSC_COMM_WORLD,&x);CHKERRQ(ierr);
  ierr = VecSetSizes(x,PETSC_DECIDE,n);CHKERRQ(ierr);
  ierr = VecSetFromOptions(x);CHKERRQ(ierr);
  ierr = VecDuplicate(x,&r);CHKERRQ(ierr);
  ierr = VecDuplicate(x,&F);CHKERRQ(ierr);

  ierr = SNESSetFunction(snes,r,FormFunction,(void*)F);CHKERRQ(ierr);

  ierr = MatCreate(PETSC_COMM_WORLD,&J);CHKERRQ(ierr);
  ierr = MatSetSizes(J,PETSC_DECIDE,PETSC_DECIDE,n,n);CHKERRQ(ierr);
  ierr = MatSetFromOptions(J);CHKERRQ(ierr);
  ierr = MatSetUp(J);CHKERRQ(ierr);
  ierr = SNESSetJacobian(snes,J,J,FormJacobian,PETSC_NULL);CHKERRQ(ierr);

  ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);

  // Set RHS of equation to F
  iniRHS[0]=11; iniRHS[1]=3; iniRHS[2]=6;
  for (i=0;i<n;i++)
   {
       ierr = VecSetValues(F,1,&i,&iniRHS[i],INSERT_VALUES);CHKERRQ(ierr);
   }
  // initial guess
  ierr = FormInitialGuess(x);CHKERRQ(ierr);

  ierr = SNESSolve(snes,PETSC_NULL,x);CHKERRQ(ierr); // The command to call
to solve the equation.

  /* Free work space.  All PETSc objects should be destroyed when they are
no longer needed*/
  ierr = VecDestroy(&x);CHKERRQ(ierr);  ierr = VecDestroy(&r);CHKERRQ(ierr);
  ierr = VecDestroy(&F);CHKERRQ(ierr);
  ierr = MatDestroy(&J);CHKERRQ(ierr);  ierr =
SNESDestroy(&snes);CHKERRQ(ierr);
  ierr = PetscFinalize();
  return 0;
}


#undef __FUNCT__
#define __FUNCT__ "FormInitialGuess"
/*
   FormInitialGuess - Computes initial guess.

   Input/Output Parameter:
.  x - the solution vector
*/
PetscErrorCode FormInitialGuess(Vec x)
{
   PetscErrorCode ierr;
   PetscScalar    pfive =10;
   ierr = VecSet(x,pfive);CHKERRQ(ierr);
   return 0;
}


/* ------------------------------------------------------------------- */
#undef __FUNCT__
#define __FUNCT__ "FormFunction"
/*
   FormFunction - Evaluates nonlinear function, F(x).

   Input Parameters:
.  snes - the SNES context
.  x - input vector
.  ctx - optional user-defined context, as set by SNESSetFunction()

   Output Parameter:
.  f - function vector

   Note:
   The user-defined context can contain any application-specific data
   needed for the function evaluation (such as various parameters, work
   vectors, and grid information).  In this program the context is just
   a vector containing the right-hand-side of the discretized PDE.
 */

PetscErrorCode FormFunction(SNES snes,Vec x,Vec f,void *ctx)
{
   Vec            g = (Vec)ctx;
   PetscScalar    *xx,*ff,*gg;
   PetscErrorCode ierr;
   PetscInt       n;

  /*
     Get pointers to vector data.
       - For default PETSc vectors, VecGetArray() returns a pointer to
         the data array.  Otherwise, the routine is implementation
dependent.
       - You MUST call VecRestoreArray() when you no longer need access to
         the array.
  */
   ierr = VecGetArray(x,&xx);CHKERRQ(ierr);
   ierr = VecGetArray(f,&ff);CHKERRQ(ierr);
   ierr = VecGetArray(g,&gg);CHKERRQ(ierr);

  /*
     Compute function
  */
   ierr = VecGetSize(x,&n);CHKERRQ(ierr);

   ff[0]=pow(xx[0],2)-pow(xx[1],2)+pow(xx[2],2)-gg[0];
   ff[1]=xx[0]*xx[1]+pow(xx[1],2)-3.0*xx[2]-gg[1];
   ff[2]=xx[0]-xx[0]*xx[2]+xx[1]*xx[2]-gg[2];

  /*
     Restore vectors
  */
  ierr = VecRestoreArray(x,&xx);CHKERRQ(ierr);
  ierr = VecRestoreArray(f,&ff);CHKERRQ(ierr);
  ierr = VecRestoreArray(g,&gg);CHKERRQ(ierr);
  return 0;
}



/* ------------------------------------------------------------------- */
#undef __FUNCT__
#define __FUNCT__ "FormJacobian"
/*
   FormJacobian - Evaluates Jacobian matrix.

   Input Parameters:
.  snes - the SNES context
.  x - input vector
.  dummy - optional user-defined context (not used here)

   Output Parameters:
.  jac - Jacobian matrix
.  B - optionally different preconditioning matrix
.  flag - flag indicating matrix structure
*/

PetscErrorCode FormJacobian(SNES snes,Vec x,Mat *jac,Mat
*B,MatStructure*flag,void *dummy)
{
  PetscScalar    *xx,A[3];
  PetscErrorCode ierr;
  PetscInt       i,n,j[3];

  /*
     Get pointer to vector data
  */
  ierr = VecGetArray(x,&xx);CHKERRQ(ierr);

  /*
     Compute Jacobian entries and insert into matrix.
      - Note that in this case we set all elements for a particular
        row at once.
  */
  ierr = VecGetSize(x,&n);CHKERRQ(ierr);
  j[0] = 0; j[1] = 1; j[2] = 2;
  /*set the first row of jacobian*/
  i=0;
  A[0]=2*xx[0]; A[1]=-2*xx[1]; A[2]=2*xx[2];
  ierr = MatSetValues(*jac,1,&i,3,j,A,INSERT_VALUES);CHKERRQ(ierr);
  /*set the second row of jacobian*/
  i=1;
  A[0]=xx[1]; A[1]=xx[0]+2*xx[1]; A[2]=-3.0;
  ierr = MatSetValues(*jac,1,&i,3,j,A,INSERT_VALUES);CHKERRQ(ierr);
  /*set the third row of jacobian*/
  i=2;
  A[0]=1-xx[2]; A[1]=xx[2]; A[2]=-xx[0]+xx[1];
  ierr = MatSetValues(*jac,1,&i,3,j,A,INSERT_VALUES);CHKERRQ(ierr);



  ierr = VecRestoreArray(x,&xx);CHKERRQ(ierr);

  /*
     Assemble matrix
  */
  ierr = MatAssemblyBegin(*jac,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
  ierr = MatAssemblyEnd(*jac,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);

  return 0;
}
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20130322/4792f3bb/attachment.html>


More information about the petsc-users mailing list