1c4762a1bSJed Brown static const char help[] = "Test parallel assembly of SBAIJ matrices\n\n"; 2c4762a1bSJed Brown 3c4762a1bSJed Brown #include <petscmat.h> 4c4762a1bSJed Brown 5d71ae5a4SJacob Faibussowitsch PetscErrorCode Assemble(MPI_Comm comm, PetscInt n, MatType mtype) 6d71ae5a4SJacob Faibussowitsch { 7c4762a1bSJed Brown Mat A; 8c4762a1bSJed Brown PetscInt first, last, i; 9c4762a1bSJed Brown PetscMPIInt rank, size; 10c4762a1bSJed Brown 11c4762a1bSJed Brown PetscFunctionBegin; 129566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &A)); 139566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, n, n)); 149566063dSJacob Faibussowitsch PetscCall(MatSetType(A, MATMPISBAIJ)); 159566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A)); 169566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 179566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 18c4762a1bSJed Brown if (rank < size - 1) { 199566063dSJacob Faibussowitsch PetscCall(MatMPISBAIJSetPreallocation(A, 1, 1, NULL, 1, NULL)); 20c4762a1bSJed Brown } else { 219566063dSJacob Faibussowitsch PetscCall(MatMPISBAIJSetPreallocation(A, 1, 2, NULL, 0, NULL)); 22c4762a1bSJed Brown } 239566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(A, &first, &last)); 249566063dSJacob Faibussowitsch PetscCall(MatSetOption(A, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_TRUE)); 25c4762a1bSJed Brown last--; 26c4762a1bSJed Brown for (i = first; i <= last; i++) { 279566063dSJacob Faibussowitsch PetscCall(MatSetValue(A, i, i, 2., INSERT_VALUES)); 289566063dSJacob Faibussowitsch if (i != n - 1) PetscCall(MatSetValue(A, i, n - 1, -1., INSERT_VALUES)); 29c4762a1bSJed Brown } 309566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 319566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 329566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 33*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 34c4762a1bSJed Brown } 35c4762a1bSJed Brown 36d71ae5a4SJacob Faibussowitsch int main(int argc, char *argv[]) 37d71ae5a4SJacob Faibussowitsch { 38c4762a1bSJed Brown MPI_Comm comm; 39c4762a1bSJed Brown PetscInt n = 6; 40c4762a1bSJed Brown 41327415f7SBarry Smith PetscFunctionBeginUser; 429566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 43c4762a1bSJed Brown comm = PETSC_COMM_WORLD; 449566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-n", &n, NULL)); 459566063dSJacob Faibussowitsch PetscCall(Assemble(comm, n, MATMPISBAIJ)); 469566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 47b122ec5aSJacob Faibussowitsch return 0; 48c4762a1bSJed Brown } 49c4762a1bSJed Brown 50c4762a1bSJed Brown /*TEST 51c4762a1bSJed Brown 52c4762a1bSJed Brown test: 53c4762a1bSJed Brown nsize: 4 5497929ea7SJunchao Zhang args: -n 1000 -mat_view ascii::ascii_info_detail 55dfd57a17SPierre Jolivet requires: double !complex !defined(PETSC_USE_64BIT_INDICES) 56c4762a1bSJed Brown 57c4762a1bSJed Brown TEST*/ 58