1cff58d65SJunchao Zhang static char help[] = "Test MatDuplicate() with new nonzeros on the duplicate\n\n"; 2cff58d65SJunchao Zhang 3cff58d65SJunchao Zhang #include <petscmat.h> 4cff58d65SJunchao Zhang int main(int argc, char **args) 5cff58d65SJunchao Zhang { 6cff58d65SJunchao Zhang Mat A, B, C; 7cff58d65SJunchao Zhang PetscInt k; 8cff58d65SJunchao Zhang const PetscInt M = 18, N = 18; 9cff58d65SJunchao Zhang PetscBool equal; 10cff58d65SJunchao Zhang PetscScalar *vals; 11cff58d65SJunchao Zhang PetscMPIInt rank; 12cff58d65SJunchao Zhang 13cff58d65SJunchao Zhang // clang-format off 14cff58d65SJunchao Zhang // i0/j0[] has a dense diagonal 15cff58d65SJunchao Zhang PetscInt i0[] = {7, 7, 8, 8, 9, 16, 17, 9, 10, 1, 1, -2, 2, 3, 3, 14, 4, 5, 10, 13, 9, 9, 10, 1, 0, 0, 5, 5, 6, 6, 13, 13, 14, -14, 4, 4, 5, 11, 11, 12, 15, 15, 16}; 16cff58d65SJunchao Zhang PetscInt j0[] = {7, 6, 8, 4, 9, 16, 17, 16, 10, 2, 1, 3, 2, 4, 3, 14, 4, 5, 15, 13, 10, 16, 11, 2, 0, 1, 5, -11, 0, 6, 15, 17, 11, 13, 4, 8, 2, 11, 17, 12, 3, 15, 9}; 17cff58d65SJunchao Zhang 18cff58d65SJunchao Zhang // i0/j0[] miss some diagonals 19cff58d65SJunchao Zhang PetscInt i1[] = {8, 5, 15, 16, 6, 13, 4, 17, 8, 9, 9, 10, -6, 12, 7, 3, -4, 1, 1, 2, 5, 5, 6, 14, 17, 8, 9, 9, 10, 4, 5, 10, 11, 1, 2}; 20cff58d65SJunchao Zhang PetscInt j1[] = {2, 3, 16, 9, 5, 17, 1, 13, 4, 10, 16, 11, -5, 12, 1, 7, -1, 2, 7, 3, 6, 11, 0, 11, 13, 4, 10, 16, 11, 8, -2, 15, 12, 7, 3}; 21cff58d65SJunchao Zhang 22cff58d65SJunchao Zhang PetscInt i2[] = {3, 4, 1, 10, 0, 1, 1, 2, 1, 1, 2, 2, 3, 3, 4, 4, 1, 2, 5, 5, 6, 4, 17, 0, 1, 1, 8, 5, 5, 6, 4, 7, 8, 5}; 23cff58d65SJunchao Zhang PetscInt j2[] = {7, 1, 2, 11, 5, 2, 7, 3, 2, 7, 3, 8, 4, 9, 3, 5, 7, 3, 6, 11, 0, 1, 13, 5, 2, 7, 4, 6, 11, 0, 1, 3, 4, 2}; 24cff58d65SJunchao Zhang // clang-format on 25cff58d65SJunchao Zhang 26cff58d65SJunchao Zhang typedef struct { 27cff58d65SJunchao Zhang PetscInt *i, *j, n; 28cff58d65SJunchao Zhang } coo_data; 29cff58d65SJunchao Zhang 30cff58d65SJunchao Zhang coo_data coos[3] = { 31cff58d65SJunchao Zhang {i0, j0, PETSC_STATIC_ARRAY_LENGTH(i0)}, 32cff58d65SJunchao Zhang {i1, j1, PETSC_STATIC_ARRAY_LENGTH(i1)}, 33cff58d65SJunchao Zhang {i2, j2, PETSC_STATIC_ARRAY_LENGTH(i2)} 34cff58d65SJunchao Zhang }; 35cff58d65SJunchao Zhang coo_data mycoo; 36cff58d65SJunchao Zhang 37cff58d65SJunchao Zhang PetscFunctionBeginUser; 38*c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &args, NULL, help)); 39cff58d65SJunchao Zhang PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 40cff58d65SJunchao Zhang 41cff58d65SJunchao Zhang mycoo = coos[rank / 3]; 42cff58d65SJunchao Zhang 43cff58d65SJunchao Zhang PetscCall(MatCreate(PETSC_COMM_WORLD, &A)); 44cff58d65SJunchao Zhang PetscCall(MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, M, N)); 45cff58d65SJunchao Zhang PetscCall(MatSetFromOptions(A)); 46cff58d65SJunchao Zhang 47cff58d65SJunchao Zhang // Assemble matrix A with the full arrays 48cff58d65SJunchao Zhang PetscCall(PetscMalloc1(mycoo.n, &vals)); 49cff58d65SJunchao Zhang for (k = 0; k < mycoo.n; k++) { 50cff58d65SJunchao Zhang vals[k] = mycoo.j[k]; 51cff58d65SJunchao Zhang PetscCall(MatSetValue(A, mycoo.i[k], mycoo.j[k], vals[k], ADD_VALUES)); 52cff58d65SJunchao Zhang } 53cff58d65SJunchao Zhang PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 54cff58d65SJunchao Zhang PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 55cff58d65SJunchao Zhang 56cff58d65SJunchao Zhang // Assemble matrix B with the 1st half of the arrays 57cff58d65SJunchao Zhang PetscCall(MatCreate(PETSC_COMM_WORLD, &B)); 58cff58d65SJunchao Zhang PetscCall(MatSetSizes(B, PETSC_DECIDE, PETSC_DECIDE, M, N)); 59cff58d65SJunchao Zhang PetscCall(MatSetFromOptions(B)); 60cff58d65SJunchao Zhang for (k = 0; k < mycoo.n / 2; k++) PetscCall(MatSetValue(B, mycoo.i[k], mycoo.j[k], vals[k], ADD_VALUES)); 61cff58d65SJunchao Zhang PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY)); 62cff58d65SJunchao Zhang PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY)); 63cff58d65SJunchao Zhang 64cff58d65SJunchao Zhang // Duplicate B to C and continue adding nozeros to C with the 2nd half 65cff58d65SJunchao Zhang PetscCall(MatDuplicate(B, MAT_COPY_VALUES, &C)); 66cff58d65SJunchao Zhang for (k = mycoo.n / 2; k < mycoo.n; k++) PetscCall(MatSetValue(C, mycoo.i[k], mycoo.j[k], vals[k], ADD_VALUES)); 67cff58d65SJunchao Zhang PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY)); 68cff58d65SJunchao Zhang PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY)); 69cff58d65SJunchao Zhang 70cff58d65SJunchao Zhang // Test if A == C 71cff58d65SJunchao Zhang PetscCall(MatMultEqual(A, C, 10, &equal)); 72cff58d65SJunchao Zhang if (!equal) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "MatDuplicate() on regular matrices failed\n")); 73cff58d65SJunchao Zhang 74cff58d65SJunchao Zhang PetscCall(PetscFree(vals)); 75cff58d65SJunchao Zhang PetscCall(MatDestroy(&A)); 76cff58d65SJunchao Zhang PetscCall(MatDestroy(&B)); 77cff58d65SJunchao Zhang PetscCall(MatDestroy(&C)); 78cff58d65SJunchao Zhang PetscCall(PetscFinalize()); 79cff58d65SJunchao Zhang return 0; 80cff58d65SJunchao Zhang } 81cff58d65SJunchao Zhang 82cff58d65SJunchao Zhang /*TEST 83cff58d65SJunchao Zhang 84cff58d65SJunchao Zhang test: 85cff58d65SJunchao Zhang nsize: 3 86cff58d65SJunchao Zhang output_file: output/empty.out 87cff58d65SJunchao Zhang 88cff58d65SJunchao Zhang TEST*/ 89