1 static const char help[] = "Test parallel assembly of SBAIJ matrices\n\n"; 2 3 #include <petscmat.h> 4 5 PetscErrorCode Assemble(MPI_Comm comm, PetscInt n, MatType mtype) { 6 Mat A; 7 PetscInt first, last, i; 8 PetscMPIInt rank, size; 9 10 PetscFunctionBegin; 11 PetscCall(MatCreate(PETSC_COMM_WORLD, &A)); 12 PetscCall(MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, n, n)); 13 PetscCall(MatSetType(A, MATMPISBAIJ)); 14 PetscCall(MatSetFromOptions(A)); 15 PetscCallMPI(MPI_Comm_size(comm, &size)); 16 PetscCallMPI(MPI_Comm_rank(comm, &rank)); 17 if (rank < size - 1) { 18 PetscCall(MatMPISBAIJSetPreallocation(A, 1, 1, NULL, 1, NULL)); 19 } else { 20 PetscCall(MatMPISBAIJSetPreallocation(A, 1, 2, NULL, 0, NULL)); 21 } 22 PetscCall(MatGetOwnershipRange(A, &first, &last)); 23 PetscCall(MatSetOption(A, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_TRUE)); 24 last--; 25 for (i = first; i <= last; i++) { 26 PetscCall(MatSetValue(A, i, i, 2., INSERT_VALUES)); 27 if (i != n - 1) PetscCall(MatSetValue(A, i, n - 1, -1., INSERT_VALUES)); 28 } 29 PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 30 PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 31 PetscCall(MatDestroy(&A)); 32 PetscFunctionReturn(0); 33 } 34 35 int main(int argc, char *argv[]) { 36 MPI_Comm comm; 37 PetscInt n = 6; 38 39 PetscFunctionBeginUser; 40 PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 41 comm = PETSC_COMM_WORLD; 42 PetscCall(PetscOptionsGetInt(NULL, NULL, "-n", &n, NULL)); 43 PetscCall(Assemble(comm, n, MATMPISBAIJ)); 44 PetscCall(PetscFinalize()); 45 return 0; 46 } 47 48 /*TEST 49 50 test: 51 nsize: 4 52 args: -n 1000 -mat_view ascii::ascii_info_detail 53 requires: double !complex !defined(PETSC_USE_64BIT_INDICES) 54 55 TEST*/ 56