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