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