[petsc-users] SuperLU_4.0 with petsc-3.1-p3 - flipped output in	KSPSolve.
    Chetan Jhurani 
    chetan.jhurani at gmail.com
       
    Mon Dec  6 21:47:36 CST 2010
    
    
  
Hello,
I have a small code that solves a 3x3 system with and without
SuperLU.  If SuperLU is used, vectors x and b (in A x = b) come
out flipped after KSPSolve.  Am I doing something stupid?  The two
outputs, the code, and the configure command are pasted below.
Thanks,
Chetan
--------------------------------------------------------------------------
Run with SuperLU:
~> lu_test -pc_factor_mat_solver_package superlu
x before solve
-10
-10
-10
b before solve
200
200
200
x after solve
200
200
200
b after solve
0.2
0.2
0.2
--------------------------------------------------------------------------
Run without SuperLU (default KSP, PC options)
~> lu_test
x before solve
-10
-10
-10
b before solve
200
200
200
x after solve
0.2
0.2
0.2
b after solve
200
200
200
--------------------------------------------------------------------------
Program:
int main(int argc, char* argv[])
{
    PetscErrorCode ierr;
    ierr = PetscInitialize(&argc, &argv, PETSC_NULL, PETSC_NULL);
CHKERRQ(ierr);
    Mat A;
    ierr = MatCreateSeqAIJ(PETSC_COMM_SELF, 3, 3, 1, PETSC_NULL, &A);
CHKERRQ(ierr);
    ierr = MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY);
    ierr = MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY);
    ierr = MatShift(A, 1000); CHKERRQ(ierr);
    Vec x, b;
    ierr = VecCreateSeq(PETSC_COMM_SELF, 3, &x); CHKERRQ(ierr);
    ierr = VecShift(x, -10); CHKERRQ(ierr);
    ierr = VecCreateSeq(PETSC_COMM_SELF, 3, &b); CHKERRQ(ierr);
    ierr = VecShift(b, 200); CHKERRQ(ierr);
    KSP ksp;
    ierr = KSPCreate(PETSC_COMM_SELF, &ksp); CHKERRQ(ierr);
    ierr = KSPSetFromOptions(ksp); CHKERRQ(ierr);
    printf("x before solve\n");
    ierr = VecView(x, PETSC_VIEWER_STDOUT_SELF); CHKERRQ(ierr);
    printf("b before solve\n");
    ierr = VecView(b, PETSC_VIEWER_STDOUT_SELF); CHKERRQ(ierr);
    ierr = KSPSetOperators(ksp, A, A, SAME_NONZERO_PATTERN); CHKERRQ(ierr);
    ierr = KSPSolve(ksp, b, x); CHKERRQ(ierr);
    printf("x after solve\n");
    ierr = VecView(x, PETSC_VIEWER_STDOUT_SELF); CHKERRQ(ierr);
    printf("b after solve\n");
    ierr = VecView(b, PETSC_VIEWER_STDOUT_SELF); CHKERRQ(ierr);
    ierr = PetscFinalize(); CHKERRQ(ierr);
}
--------------------------------------------------------------------------
configure command:
python ./config/configure.py -with-clanguage=C++ -with-debugging=no
--with-gnu-compilers=1 --with-mpi=1 --with-umfpack=1 --with-superlu=1
--with-hypre=1 --download-umfpack=1 --download-superlu=1
--download-hypre=1 --with-mumps=1 --download-mumps=1 --with-parmetis=1
--download-parmetis=1 --with-scalapack=1 --download-scalapack=1
--with-blacs=1 --download-blacs=1 -with-c-support=1
--------------------------------------------------------------------------
    
    
More information about the petsc-users
mailing list