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