[petsc-users] MatRARt
John Fettig
john.fettig at gmail.com
Wed Jul 18 12:09:01 CDT 2012
Should MatRARt work for SeqAIJ matrices? I don't understand what is
wrong with this test code:
#include <petsc.h>
#undef __FUNCT__
#define __FUNCT__ "main"
int main(int argc,char **args)
{
Mat A,R,C;
PetscInt ai[5] = {0,2,4,6,8};
PetscInt aj[8] = {0,3,1,2,1,2,0,3};
PetscReal av[8] = {1,2,3,4,5,6,7,8};
PetscInt ri[3] = {0,1,2};
PetscInt rj[2] = {1,2};
PetscReal rv[2] = {1,1};
PetscErrorCode ierr;
PetscInitialize(&argc,&args,(char *)0,(char *)0);
ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_WORLD,4,4,ai,aj,av,&A);CHKERRQ(ierr);
ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_WORLD,2,4,ri,rj,rv,&R);CHKERRQ(ierr);
//this is ok:
//ierr = MatMatMult(R,A,MAT_INITIAL_MATRIX,1.0,&C);CHKERRQ(ierr);
//this is not:
ierr = MatRARt(A,R,MAT_INITIAL_MATRIX,1.0,&C);CHKERRQ(ierr);
ierr = PetscFinalize();
return 0;
}
When I run it, it says:
[0]PETSC ERROR: --------------------- Error Message
------------------------------------
[0]PETSC ERROR: No support for this operation for this object type!
[0]PETSC ERROR: MatGetColumnIJ() not supported for matrix type seqaij!
[0]PETSC ERROR:
------------------------------------------------------------------------
[0]PETSC ERROR: Petsc Release Version 3.3.0, Patch 2, Fri Jul 13
15:42:00 CDT 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: See docs/index.html for manual pages.
[0]PETSC ERROR:
------------------------------------------------------------------------
[0]PETSC ERROR: ./rart_simple on a linux-deb named lagrange.tomato by
jfe Wed Jul 18 13:07:33 2012
[0]PETSC ERROR: Libraries linked from
/home/jfe/local/petsc-3.3-p2/linux-debug/lib
[0]PETSC ERROR: Configure run at Wed Jul 18 13:00:59 2012
[0]PETSC ERROR: Configure options --with-x=0
--download-f-blas-lapack=1 --with-mpi=1 --with-mpi-shared=1
--with-mpi=1 --download-mpich=1 --with-debugging=1
--with-gnu-compilers=yes --with-shared-libraries=1 --with-c++-support
--with-clanguage=C
[0]PETSC ERROR:
------------------------------------------------------------------------
[0]PETSC ERROR: MatTransposeColoringCreate_SeqAIJ() line 1294 in
/home/jfe/local/petsc-3.3-p2/src/mat/impls/aij/seq/matmatmult.c
[0]PETSC ERROR: MatTransposeColoringCreate() line 9318 in
/home/jfe/local/petsc-3.3-p2/src/mat/interface/matrix.c
[0]PETSC ERROR: MatRARtSymbolic_SeqAIJ_SeqAIJ() line 342 in
/home/jfe/local/petsc-3.3-p2/src/mat/impls/aij/seq/matrart.c
[0]PETSC ERROR: MatRARt_SeqAIJ_SeqAIJ() line 541 in
/home/jfe/local/petsc-3.3-p2/src/mat/impls/aij/seq/matrart.c
[0]PETSC ERROR: MatRARt() line 8405 in
/home/jfe/local/petsc-3.3-p2/src/mat/interface/matrix.c
[0]PETSC ERROR: main() line 25 in src/ksp/ksp/examples/tutorials/rart_simple.c
What am I doing wrong?
Thanks,
John
More information about the petsc-users
mailing list