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 PetscInt ia[3]={0,2,4}; 8c4762a1bSJed Brown PetscInt ja[4]={0,1,0,1}; 9c4762a1bSJed Brown PetscScalar c[16]={0,1,2,3,4,5,6,7,8,9,10,11,12,13,14,15}; 10c4762a1bSJed Brown PetscMPIInt size; 11c4762a1bSJed Brown Mat ssbaij,msbaij; 12c4762a1bSJed Brown Vec x,y; 13c4762a1bSJed Brown 14*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc,&argv,NULL,help)); 155f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 162c71b3e2SJacob Faibussowitsch PetscCheckFalse(size != 2,PETSC_COMM_WORLD,PETSC_ERR_SUP,"This is an example with two processors only!"); 175f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_SELF,&ssbaij)); 185f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetType(ssbaij,MATSEQSBAIJ)); 195f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetBlockSize(ssbaij,2)); 205f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(ssbaij,4,8,4,8)); 215f80ce2aSJacob Faibussowitsch CHKERRQ(MatSeqSBAIJSetPreallocationCSR(ssbaij,2,ia,ja,c)); 225f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateMPIMatConcatenateSeqMat(PETSC_COMM_WORLD,ssbaij,PETSC_DECIDE,MAT_INITIAL_MATRIX,&msbaij)); 235f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(msbaij,PETSC_VIEWER_STDOUT_WORLD)); 245f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&msbaij)); 255f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateMPIMatConcatenateSeqMat(PETSC_COMM_WORLD,ssbaij,4,MAT_INITIAL_MATRIX,&msbaij)); 265f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(msbaij,PETSC_VIEWER_STDOUT_WORLD)); 275f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateVecs(msbaij,&x,&y)); 285f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(x,1)); 295f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(msbaij,x,y)); 305f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(y,PETSC_VIEWER_STDOUT_WORLD)); 315f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultAdd(msbaij,x,x,y)); 325f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(y,PETSC_VIEWER_STDOUT_WORLD)); 335f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetDiagonal(msbaij,y)); 345f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(y,PETSC_VIEWER_STDOUT_WORLD)); 355f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&x)); 365f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&y)); 375f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&msbaij)); 385f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&ssbaij)); 39*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 40*b122ec5aSJacob Faibussowitsch return 0; 41c4762a1bSJed Brown } 42c4762a1bSJed Brown 43c4762a1bSJed Brown /*TEST 44c4762a1bSJed Brown 45c4762a1bSJed Brown test: 46c4762a1bSJed Brown nsize: 2 47c4762a1bSJed Brown filter: sed "s?\.??g" 48c4762a1bSJed Brown 49c4762a1bSJed Brown TEST*/ 50