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 149566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&argv,NULL,help)); 159566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 16*be096a46SBarry Smith PetscCheck(size == 2,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is an example with two processors only!"); 179566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_SELF,&ssbaij)); 189566063dSJacob Faibussowitsch PetscCall(MatSetType(ssbaij,MATSEQSBAIJ)); 199566063dSJacob Faibussowitsch PetscCall(MatSetBlockSize(ssbaij,2)); 209566063dSJacob Faibussowitsch PetscCall(MatSetSizes(ssbaij,4,8,4,8)); 219566063dSJacob Faibussowitsch PetscCall(MatSeqSBAIJSetPreallocationCSR(ssbaij,2,ia,ja,c)); 229566063dSJacob Faibussowitsch PetscCall(MatCreateMPIMatConcatenateSeqMat(PETSC_COMM_WORLD,ssbaij,PETSC_DECIDE,MAT_INITIAL_MATRIX,&msbaij)); 239566063dSJacob Faibussowitsch PetscCall(MatView(msbaij,PETSC_VIEWER_STDOUT_WORLD)); 249566063dSJacob Faibussowitsch PetscCall(MatDestroy(&msbaij)); 259566063dSJacob Faibussowitsch PetscCall(MatCreateMPIMatConcatenateSeqMat(PETSC_COMM_WORLD,ssbaij,4,MAT_INITIAL_MATRIX,&msbaij)); 269566063dSJacob Faibussowitsch PetscCall(MatView(msbaij,PETSC_VIEWER_STDOUT_WORLD)); 279566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(msbaij,&x,&y)); 289566063dSJacob Faibussowitsch PetscCall(VecSet(x,1)); 299566063dSJacob Faibussowitsch PetscCall(MatMult(msbaij,x,y)); 309566063dSJacob Faibussowitsch PetscCall(VecView(y,PETSC_VIEWER_STDOUT_WORLD)); 319566063dSJacob Faibussowitsch PetscCall(MatMultAdd(msbaij,x,x,y)); 329566063dSJacob Faibussowitsch PetscCall(VecView(y,PETSC_VIEWER_STDOUT_WORLD)); 339566063dSJacob Faibussowitsch PetscCall(MatGetDiagonal(msbaij,y)); 349566063dSJacob Faibussowitsch PetscCall(VecView(y,PETSC_VIEWER_STDOUT_WORLD)); 359566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 369566063dSJacob Faibussowitsch PetscCall(VecDestroy(&y)); 379566063dSJacob Faibussowitsch PetscCall(MatDestroy(&msbaij)); 389566063dSJacob Faibussowitsch PetscCall(MatDestroy(&ssbaij)); 399566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 40b122ec5aSJacob 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