[petsc-users] KSPSetDM problem

Michele Rosso mrosso at uci.edu
Thu Nov 14 19:49:31 CST 2013


Hi,

I am  solving a Poisson equation (cell-centered discretization) with 
multigrid.
I am using a 2D version of ksp/ksp/examples/tutorials/ex34.c as a 
reference.
The following code performs the same calculations as ex34:

   call PetscInitialize(PETSC_NULL_CHARACTER,ierr)

     call DMDACreate2d(PETSC_COMM_WORLD, DMDA_BOUNDARY_NONE, 
DMDA_BOUNDARY_NONE,&
          & DMDA_STENCIL_STAR, 16, 16, PETSC_DECIDE, PETSC_DECIDE, &
          & 1, 1, PETSC_NULL_INTEGER, PETSC_NULL_INTEGER,da,ierr)
     call DMDASetInterpolationType(da, DMDA_Q0, ierr);
     call DMCreateGlobalVector(da,b,ierr)
     call create_rhs(da,b)
     call DMCreateMatrix(da,'aij',A,ierr)
     call create_matrix(da,A)
     call VecDuplicate(b,x,ierr)

     ! Create solver
     call KSPCreate(PETSC_COMM_WORLD,solver,ierr)
     call KSPSetDM(solver,da,ierr)
     call KSPSetDMActive(solver,PETSC_FALSE,ierr)
     call KSPSetOperators(solver,A,A,SAME_NONZERO_PATTERN,ierr)
     call KSPSetType(solver,'cg',ierr)
     call KSPSetNormType(solver,KSP_NORM_UNPRECONDITIONED,ierr)
     call KSPSetFromOptions(solver,ierr)
     call KSPSolve(solver,b,x,ierr)

     rerr = error(da,x)

     call PetscFinalize(ierr)


Both ex34 and my code produce exactly the same results except when I run 
with

-ksp_view -ksp_monitor_true_residual -ksp_converged_reason -ksp_type cg 
-pc_type mg -pc_mg_levels 4

In this case my code results in the following error, despite 
KSPSetDMActive is used after KSPSetDM:

[0]PETSC ERROR: --------------------- Error Message 
------------------------------------
[0]PETSC ERROR: Object is in wrong state!
[0]PETSC ERROR: You called KSPSetDM() but did not use 
DMKSPSetComputeOperators() or KSPSetDMActive(dm,PETSC_FALSE);!
[0]PETSC ERROR: 
------------------------------------------------------------------------
[0]PETSC ERROR: Petsc Release Version 3.4.3, Oct, 15, 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: ./test on a linux-gnu-dbg named enterprise-A by mic Thu 
Nov 14 17:38:23 2013
[0]PETSC ERROR: Libraries linked from 
/opt/petsc/petsc-3.4.3/linux-gnu-dbg/lib
[0]PETSC ERROR: Configure run at Tue Oct 29 22:43:45 2013
[0]PETSC ERROR: Configure options --known-mpi-shared="0 " 
--known-memcmp-ok  --with-debugging="1 " --with-fortran-datatypes 
--with-shared-libraries=0 --with-dynamic-loading=0 
--with-mpi-compilers="1 " --download-blacs="1 " 
--download-superlu_dist="1 " --download-metis="1 " 
--download-parmetis="1 " --download-ml=1 PETSC_ARCH=linux-gnu-dbg
[0]PETSC ERROR: 
------------------------------------------------------------------------
[0]PETSC ERROR: KSPSetUp() line 230 in 
/opt/petsc/petsc-3.4.3/src/ksp/ksp/interface/itfunc.c
[0]PETSC ERROR: PCSetUp_MG() line 730 in 
/opt/petsc/petsc-3.4.3/src/ksp/pc/impls/mg/mg.c
[0]PETSC ERROR: PCSetUp() line 890 in 
/opt/petsc/petsc-3.4.3/src/ksp/pc/interface/precon.c
[0]PETSC ERROR: KSPSetUp() line 278 in 
/opt/petsc/petsc-3.4.3/src/ksp/ksp/interface/itfunc.c
[0]PETSC ERROR: KSPSolve() line 399 in 
/opt/petsc/petsc-3.4.3/src/ksp/ksp/interface/itfunc.c


Am I missing anything?
Thank you,

Michele







-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20131114/35e95d64/attachment-0001.html>


More information about the petsc-users mailing list