17887b771SAlex Lindsay static char help[] = "Tests resetting preallocation after filling the full sparsity pattern"; 27887b771SAlex Lindsay 37887b771SAlex Lindsay #include <petscmat.h> 47887b771SAlex Lindsay 57887b771SAlex Lindsay PetscErrorCode Assemble(Mat mat) 67887b771SAlex Lindsay { 77887b771SAlex Lindsay PetscInt idx[4], i; 87887b771SAlex Lindsay PetscScalar vals[16]; 97887b771SAlex Lindsay int rank; 107887b771SAlex Lindsay 117887b771SAlex Lindsay PetscFunctionBegin; 12*458b0db5SMartin Diehl PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 137887b771SAlex Lindsay for (i = 0; i < 16; ++i) vals[i] = 1; 147887b771SAlex Lindsay if (rank == 0) { 157887b771SAlex Lindsay // element 0 167887b771SAlex Lindsay idx[0] = 0; 177887b771SAlex Lindsay idx[1] = 1; 187887b771SAlex Lindsay idx[2] = 2; 197887b771SAlex Lindsay idx[3] = 3; 207887b771SAlex Lindsay PetscCall(MatSetValues(mat, 4, idx, 4, idx, vals, ADD_VALUES)); 217887b771SAlex Lindsay // element 1 227887b771SAlex Lindsay idx[0] = 3; 237887b771SAlex Lindsay idx[1] = 2; 247887b771SAlex Lindsay idx[2] = 4; 257887b771SAlex Lindsay idx[3] = 5; 267887b771SAlex Lindsay PetscCall(MatSetValues(mat, 4, idx, 4, idx, vals, ADD_VALUES)); 277887b771SAlex Lindsay } else { 287887b771SAlex Lindsay // element 2 297887b771SAlex Lindsay idx[0] = 6; 307887b771SAlex Lindsay idx[1] = 0; 317887b771SAlex Lindsay idx[2] = 3; 327887b771SAlex Lindsay idx[3] = 7; 337887b771SAlex Lindsay PetscCall(MatSetValues(mat, 4, idx, 4, idx, vals, ADD_VALUES)); 347887b771SAlex Lindsay // element 3 357887b771SAlex Lindsay idx[0] = 7; 367887b771SAlex Lindsay idx[1] = 3; 377887b771SAlex Lindsay idx[2] = 5; 387887b771SAlex Lindsay idx[3] = 8; 397887b771SAlex Lindsay PetscCall(MatSetValues(mat, 4, idx, 4, idx, vals, ADD_VALUES)); 407887b771SAlex Lindsay } 417887b771SAlex Lindsay PetscCall(MatAssemblyBegin(mat, MAT_FINAL_ASSEMBLY)); 427887b771SAlex Lindsay PetscCall(MatAssemblyEnd(mat, MAT_FINAL_ASSEMBLY)); 437887b771SAlex Lindsay PetscFunctionReturn(PETSC_SUCCESS); 447887b771SAlex Lindsay } 457887b771SAlex Lindsay 467887b771SAlex Lindsay int main(int argc, char **argv) 477887b771SAlex Lindsay { 487887b771SAlex Lindsay Mat mat; 497887b771SAlex Lindsay int rank; 507887b771SAlex Lindsay 517887b771SAlex Lindsay PetscFunctionBeginUser; 527887b771SAlex Lindsay PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 53*458b0db5SMartin Diehl PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 547887b771SAlex Lindsay PetscCall(MatCreate(PETSC_COMM_WORLD, &mat)); 557887b771SAlex Lindsay if (rank == 0) PetscCall(MatSetSizes(mat, 6, 6, PETSC_DETERMINE, PETSC_DETERMINE)); 567887b771SAlex Lindsay else PetscCall(MatSetSizes(mat, 3, 3, PETSC_DETERMINE, PETSC_DETERMINE)); 577887b771SAlex Lindsay PetscCall(MatSetFromOptions(mat)); 587887b771SAlex Lindsay if (rank == 0) { 597887b771SAlex Lindsay PetscInt ndz[6], noz[6]; 607887b771SAlex Lindsay ndz[0] = 4; 617887b771SAlex Lindsay noz[0] = 2; 627887b771SAlex Lindsay ndz[1] = 4; 637887b771SAlex Lindsay noz[1] = 0; 647887b771SAlex Lindsay ndz[2] = 6; 657887b771SAlex Lindsay noz[2] = 0; 667887b771SAlex Lindsay ndz[3] = 6; 677887b771SAlex Lindsay noz[3] = 3; 687887b771SAlex Lindsay ndz[4] = 4; 697887b771SAlex Lindsay noz[4] = 0; 707887b771SAlex Lindsay ndz[5] = 4; 717887b771SAlex Lindsay noz[5] = 2; 727887b771SAlex Lindsay PetscCall(MatMPIAIJSetPreallocation(mat, 0, ndz, 0, noz)); 737887b771SAlex Lindsay } else { 747887b771SAlex Lindsay PetscInt ndz[3], noz[3]; 757887b771SAlex Lindsay ndz[0] = 2; 767887b771SAlex Lindsay noz[0] = 2; 777887b771SAlex Lindsay ndz[1] = 3; 787887b771SAlex Lindsay noz[1] = 3; 797887b771SAlex Lindsay ndz[2] = 2; 807887b771SAlex Lindsay noz[2] = 2; 817887b771SAlex Lindsay PetscCall(MatMPIAIJSetPreallocation(mat, 0, ndz, 0, noz)); 827887b771SAlex Lindsay } 837887b771SAlex Lindsay PetscCall(MatSetUp(mat)); 847887b771SAlex Lindsay PetscCall(Assemble(mat)); 857887b771SAlex Lindsay PetscCall(MatView(mat, NULL)); 867887b771SAlex Lindsay PetscCall(MatResetPreallocation(mat)); 877887b771SAlex Lindsay PetscCall(Assemble(mat)); 887887b771SAlex Lindsay PetscCall(MatView(mat, NULL)); 897887b771SAlex Lindsay PetscCall(MatDestroy(&mat)); 907887b771SAlex Lindsay PetscCall(PetscFinalize()); 917887b771SAlex Lindsay return 0; 927887b771SAlex Lindsay } 937887b771SAlex Lindsay 947887b771SAlex Lindsay /*TEST 957887b771SAlex Lindsay 967887b771SAlex Lindsay test: 977887b771SAlex Lindsay suffix: 1 987887b771SAlex Lindsay nsize: 2 997887b771SAlex Lindsay 1007887b771SAlex Lindsay TEST*/ 101