1 static char help[] = "Basic test of various routines with SBAIJ matrices\n\n"; 2 3 #include <petscmat.h> 4 5 int main(int argc,char **argv) 6 { 7 PetscErrorCode ierr; 8 PetscInt ia[3]={0,2,4}; 9 PetscInt ja[4]={0,1,0,1}; 10 PetscScalar c[16]={0,1,2,3,4,5,6,7,8,9,10,11,12,13,14,15}; 11 PetscMPIInt size; 12 Mat ssbaij,msbaij; 13 Vec x,y; 14 15 ierr = PetscInitialize(&argc,&argv,NULL,help);if (ierr) return ierr; 16 CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 17 PetscCheckFalse(size != 2,PETSC_COMM_WORLD,PETSC_ERR_SUP,"This is an example with two processors only!"); 18 CHKERRQ(MatCreate(PETSC_COMM_SELF,&ssbaij)); 19 CHKERRQ(MatSetType(ssbaij,MATSEQSBAIJ)); 20 CHKERRQ(MatSetBlockSize(ssbaij,2)); 21 CHKERRQ(MatSetSizes(ssbaij,4,8,4,8)); 22 CHKERRQ(MatSeqSBAIJSetPreallocationCSR(ssbaij,2,ia,ja,c)); 23 CHKERRQ(MatCreateMPIMatConcatenateSeqMat(PETSC_COMM_WORLD,ssbaij,PETSC_DECIDE,MAT_INITIAL_MATRIX,&msbaij)); 24 CHKERRQ(MatView(msbaij,PETSC_VIEWER_STDOUT_WORLD)); 25 CHKERRQ(MatDestroy(&msbaij)); 26 CHKERRQ(MatCreateMPIMatConcatenateSeqMat(PETSC_COMM_WORLD,ssbaij,4,MAT_INITIAL_MATRIX,&msbaij)); 27 CHKERRQ(MatView(msbaij,PETSC_VIEWER_STDOUT_WORLD)); 28 CHKERRQ(MatCreateVecs(msbaij,&x,&y)); 29 CHKERRQ(VecSet(x,1)); 30 CHKERRQ(MatMult(msbaij,x,y)); 31 CHKERRQ(VecView(y,PETSC_VIEWER_STDOUT_WORLD)); 32 CHKERRQ(MatMultAdd(msbaij,x,x,y)); 33 CHKERRQ(VecView(y,PETSC_VIEWER_STDOUT_WORLD)); 34 CHKERRQ(MatGetDiagonal(msbaij,y)); 35 CHKERRQ(VecView(y,PETSC_VIEWER_STDOUT_WORLD)); 36 CHKERRQ(VecDestroy(&x)); 37 CHKERRQ(VecDestroy(&y)); 38 CHKERRQ(MatDestroy(&msbaij)); 39 CHKERRQ(MatDestroy(&ssbaij)); 40 ierr = PetscFinalize(); 41 return ierr; 42 } 43 44 /*TEST 45 46 test: 47 nsize: 2 48 filter: sed "s?\.??g" 49 50 TEST*/ 51