1c4762a1bSJed Brown static char help[] = "Basic test of various routines with SBAIJ matrices\n\n"; 2c4762a1bSJed Brown 3c4762a1bSJed Brown #include <petscmat.h> 4c4762a1bSJed Brown 5*d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 6*d71ae5a4SJacob Faibussowitsch { 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 14327415f7SBarry Smith PetscFunctionBeginUser; 159566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 169566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 17be096a46SBarry Smith PetscCheck(size == 2, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is an example with two processors only!"); 189566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_SELF, &ssbaij)); 199566063dSJacob Faibussowitsch PetscCall(MatSetType(ssbaij, MATSEQSBAIJ)); 209566063dSJacob Faibussowitsch PetscCall(MatSetBlockSize(ssbaij, 2)); 219566063dSJacob Faibussowitsch PetscCall(MatSetSizes(ssbaij, 4, 8, 4, 8)); 229566063dSJacob Faibussowitsch PetscCall(MatSeqSBAIJSetPreallocationCSR(ssbaij, 2, ia, ja, c)); 239566063dSJacob Faibussowitsch PetscCall(MatCreateMPIMatConcatenateSeqMat(PETSC_COMM_WORLD, ssbaij, PETSC_DECIDE, MAT_INITIAL_MATRIX, &msbaij)); 249566063dSJacob Faibussowitsch PetscCall(MatView(msbaij, PETSC_VIEWER_STDOUT_WORLD)); 259566063dSJacob Faibussowitsch PetscCall(MatDestroy(&msbaij)); 269566063dSJacob Faibussowitsch PetscCall(MatCreateMPIMatConcatenateSeqMat(PETSC_COMM_WORLD, ssbaij, 4, MAT_INITIAL_MATRIX, &msbaij)); 279566063dSJacob Faibussowitsch PetscCall(MatView(msbaij, PETSC_VIEWER_STDOUT_WORLD)); 289566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(msbaij, &x, &y)); 299566063dSJacob Faibussowitsch PetscCall(VecSet(x, 1)); 309566063dSJacob Faibussowitsch PetscCall(MatMult(msbaij, x, y)); 319566063dSJacob Faibussowitsch PetscCall(VecView(y, PETSC_VIEWER_STDOUT_WORLD)); 329566063dSJacob Faibussowitsch PetscCall(MatMultAdd(msbaij, x, x, y)); 339566063dSJacob Faibussowitsch PetscCall(VecView(y, PETSC_VIEWER_STDOUT_WORLD)); 349566063dSJacob Faibussowitsch PetscCall(MatGetDiagonal(msbaij, y)); 359566063dSJacob Faibussowitsch PetscCall(VecView(y, PETSC_VIEWER_STDOUT_WORLD)); 369566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 379566063dSJacob Faibussowitsch PetscCall(VecDestroy(&y)); 389566063dSJacob Faibussowitsch PetscCall(MatDestroy(&msbaij)); 399566063dSJacob Faibussowitsch PetscCall(MatDestroy(&ssbaij)); 409566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 41b122ec5aSJacob Faibussowitsch return 0; 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