Problem with SNESSolve()
Rafael Santos Coelho
rafaelsantoscoelho at gmail.com
Mon Mar 10 07:24:08 CDT 2008
Hello, guys! (:
I've been trying to code a program that solves a simple 3D
diffusion-convection problem, but there's something weird going on. Whenever
I run it, I get the following error message:
$ mpirun -np 1 ./ex8 -x 8 -y 8 -z 8
[0]PETSC ERROR: --------------------- Error Message
------------------------------------
[0]PETSC ERROR: Invalid argument!
[0]PETSC ERROR: Wrong type of object: Parameter # 2!
[0]PETSC ERROR:
------------------------------------------------------------------------
[0]PETSC ERROR: VecMaxPointwiseDivide() line 35 in
src/vec/vec/interface/rvector.c
[0]PETSC ERROR: SNESLineSearchCubic() line 544 in src/snes/impls/ls/ls.c
[0]PETSC ERROR: SNESSolve_LS() line 211 in src/snes/impls/ls/ls.c
[0]PETSC ERROR: SNESSolve() line 1871 in src/snes/interface/snes.c
[0]PETSC ERROR: main() line 135 in src/snes/examples/tutorials/ex8.c
Well, line 135 is:
ierr = SNESSolve(snes, PETSC_NULL, x); CHKERRQ(ierr);
Here's the main function: (I removed several comments of mine to make it
smaller)
typedef struct {
DA da;
} appContext;
PetscScalar Hx, Hy, Hz, Hx2, Hy2, Hz2;
PetscInt M = 9, N = 9, P = 9;
#undef __FUNCT__
#define __FUNCT__ "main"
int main(int argc, char **argv) {
SNES snes;
Mat J;
appContext ctx;
Vec x, r;
PetscErrorCode ierr;
PetscScalar hx, hy, hz;
PetscInitialize(&argc, &argv, (char*) 0, help);
ierr = PetscOptionsGetInt(PETSC_NULL, "-y", &M, PETSC_NULL);
CHKERRQ(ierr);
ierr = PetscOptionsGetInt(PETSC_NULL, "-x", &N, PETSC_NULL);
CHKERRQ(ierr);
ierr = PetscOptionsGetInt(PETSC_NULL, "-z", &P, PETSC_NULL);
CHKERRQ(ierr);
hx = 1.0 / (N+1); hy = 1.0 / (M+1); hz = 1.0 / (P+1);
// global variables
Hx2 = 1.0 / (2.0*hx*hx); Hy2 = 1.0 / (2.0*hy*hy); Hz2 = 1.0 / (2.0*hz*hz
);
Hx = 1.0 / (hx*hx); Hy = 1.0 / (hy*hy); Hz = 1.0 / (hz*hz);
ierr = SNESCreate(PETSC_COMM_WORLD, &snes); CHKERRQ(ierr);
ierr = DACreate3d(PETSC_COMM_WORLD, DA_NONPERIODIC, DA_STENCIL_STAR, M,
N, P, PETSC_DECIDE, PETSC_DECIDE, PETSC_DECIDE, 1, 1, PETSC_NULL,
PETSC_NULL, PETSC_NULL, &ctx.da);
ierr = DACreateGlobalVector(ctx.da, &x); CHKERRQ(ierr);
ierr = VecDuplicate(x, &r); CHKERRQ(ierr);
ierr = DAGetMatrix(ctx.da, MATAIJ, &J); CHKERRQ(ierr);
ierr = SNESSetFunction(snes, r, computeFVector, &ctx); CHKERRQ(ierr);
ierr = SNESSetJacobian(snes, J, J, computeJacobian, &ctx); CHKERRQ(ierr);
ierr = SNESSetFromOptions(snes); CHKERRQ(ierr);
ierr = formInitialGuess(x); CHKERRQ(ierr);
ierr = SNESSolve(snes, PETSC_NULL, x); CHKERRQ(ierr);
ierr = VecDestroy(x); CHKERRQ(ierr);
ierr = VecDestroy(r); CHKERRQ(ierr);
ierr = MatDestroy(J); CHKERRQ(ierr);
ierr = SNESDestroy(snes); CHKERRQ(ierr);
ierr = DADestroy(ctx.da); CHKERRQ(ierr);
ierr = PetscFinalize(); CHKERRQ(ierr);
PetscFunctionReturn(0);
}
What am I missing here?
Regards,
Rafael Santos Coelho
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20080310/6b39ca6a/attachment.htm>
More information about the petsc-users
mailing list