[petsc-users] Bug when multipling a SBAIJ matrix with a vector

Andreas Hauffe andreas.hauffe at tu-dresden.de
Thu Oct 21 03:59:34 CDT 2010


Hi,

I think there is a bug when multiplying an SBAIJ matrix with a vector, if the 
matrix has a zero/missing row. I tried to write a small example:
|x x| |1| = |0|
|x 1| |x|    |0|
The result since 3.1 is:
|x x| |1| = |1|
|x 1| |x|    |0|

I add the fortran code:
program main

  implicit none
  
#include "finclude/petscsys.h"
#include "finclude/petscvec.h"
#include "finclude/petscvec.h90"
#include "finclude/petscmat.h"
#include "finclude/petscmat.h90"

  Mat            :: KaaS
  Vec            :: v0,y
  PetscInt       :: m
  PetscInt       :: bs
  
  PetscErrorCode :: ierr
  
  call petscInitialize(PETSC_NULL_CHARACTER,ierr); CHKERRQ(ierr)

  m  = 2
  bs = 1
  call MatCreate(PETSC_COMM_WORLD,KaaS,ierr); CHKERRQ(ierr)
  call MatSetSizes(KaaS,PETSC_DECIDE,PETSC_DECIDE,m,m,ierr); CHKERRQ(ierr)
  call MatSetType(KaaS,MATSEQSBAIJ,ierr); CHKERRQ(ierr)
!  call MatSetType(KaaS,MATAIJ,ierr); CHKERRQ(ierr)
  call MatSetFromOptions(KaaS,ierr); CHKERRQ(ierr)
  
!  call MatSetValue(KaaS, 0, 0, 0.d0, ADD_VALUES, ierr); CHKERRQ(ierr)
  call MatSetValue(KaaS, 1, 1, 1.d0, ADD_VALUES, ierr); CHKERRQ(ierr)
  
  call MatAssemblyBegin(KaaS,MAT_FINAL_ASSEMBLY,ierr); CHKERRQ(ierr)
  call MatAssemblyEnd  (KaaS,MAT_FINAL_ASSEMBLY,ierr); CHKERRQ(ierr)

  call MatGetVecs(KaaS,y,v0,ierr); CHKERRQ(ierr)
  call VecSetValue(v0,0,1.D0,INSERT_VALUES,ierr); CHKERRQ(ierr)
  call MatMult(KaaS,v0,y,ierr); CHKERRQ(ierr)
  call VecView(y,PETSC_NULL_OBJECT,ierr); CHKERRQ(ierr)

  call VecDestroy(y,ierr); CHKERRQ(ierr)
  call VecDestroy(v0,ierr); CHKERRQ(ierr)
  call MatDestroy(KaaS,ierr); CHKERRQ(ierr)

  call petscFinalize(ierr)

end program main


best
-- 
Andreas Hauffe

----------------------------------------------------------------------------------------------------
Technische Universität Dresden
Institut für Luft- und Raumfahrttechnik / Institute of Aerospace Engineering
Lehrstuhl für Luftfahrzeugtechnik / Chair of Aircraft Engineering

D-01062 Dresden
Germany

phone   : (++49)351 463 38496
fax     : (++49)351 463 37263
mail    : andreas.hauffe at tu-dresden.de
Website : http://tu-dresden.de/mw/ilr/lft
----------------------------------------------------------------------------------------------------


More information about the petsc-users mailing list