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