Actual source code: ex250.c

  1: static char help[] = "Test Mat products \n\n";

  3: #include <petscmat.h>
  4: int main(int argc, char **args)
  5: {
  6:   Mat            A = NULL, B = NULL, C = NULL, D = NULL, E = NULL;
  7:   PetscInt       k;
  8:   const PetscInt M = 18, N = 18;
  9:   PetscMPIInt    rank;

 11:   /* A, B are 18 x 18 nonsymmetric matrices and have the same sparsity pattern but different values.
 12:      Big enough to have complex communication patterns but still small enough for debugging.
 13:   */
 14:   PetscInt Ai[] = {0, 0, 1, 1, 2, 2, 3, 3, 4, 4, 5, 5, 5, 6, 6, 7, 7, 8, 8, 9, 9, 10, 10, 11, 11, 12, 12, 13, 13, 13, 14, 14, 15, 15, 16, 16, 17, 17};
 15:   PetscInt Aj[] = {0, 1, 2, 7, 3, 8, 4, 9, 5, 8, 2, 6, 11, 0, 7, 1, 6, 2, 4, 10, 16, 11, 15, 12, 17, 12, 13, 14, 15, 17, 11, 13, 3, 16, 9, 15, 11, 13};
 16:   PetscInt Bi[] = {0, 0, 1, 1, 2, 2, 3, 3, 4, 4, 5, 5, 5, 6, 6, 7, 7, 8, 8, 9, 9, 10, 10, 11, 11, 12, 12, 13, 13, 13, 14, 14, 15, 15, 16, 16, 17, 17};
 17:   PetscInt Bj[] = {0, 1, 2, 7, 3, 8, 4, 9, 5, 8, 2, 6, 11, 0, 7, 1, 6, 2, 4, 10, 16, 11, 15, 12, 17, 12, 13, 14, 15, 17, 11, 13, 3, 16, 9, 15, 11, 13};

 19:   PetscInt Annz = PETSC_STATIC_ARRAY_LENGTH(Ai);
 20:   PetscInt Bnnz = PETSC_STATIC_ARRAY_LENGTH(Bi);

 22:   PetscFunctionBeginUser;
 23:   PetscCall(PetscInitialize(&argc, &args, (char *)0, help));
 24:   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));

 26:   PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
 27:   PetscCall(MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, M, N));
 28:   PetscCall(MatSetFromOptions(A));
 29:   PetscCall(MatSeqAIJSetPreallocation(A, 2, NULL));
 30:   PetscCall(MatMPIAIJSetPreallocation(A, 2, NULL, 2, NULL));
 31:   PetscCall(MatSetOption(A, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE));

 33:   if (rank == 0) {
 34:     for (k = 0; k < Annz; k++) PetscCall(MatSetValue(A, Ai[k], Aj[k], Ai[k] + Aj[k] + 1.0, INSERT_VALUES));
 35:   }

 37:   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
 38:   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));

 40:   PetscCall(MatCreate(PETSC_COMM_WORLD, &B));
 41:   PetscCall(MatSetSizes(B, PETSC_DECIDE, PETSC_DECIDE, M, N));
 42:   PetscCall(MatSetFromOptions(B));
 43:   PetscCall(MatSeqAIJSetPreallocation(B, 2, NULL));
 44:   PetscCall(MatMPIAIJSetPreallocation(B, 2, NULL, 2, NULL));
 45:   PetscCall(MatSetOption(B, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE));

 47:   if (rank == 0) {
 48:     for (k = 0; k < Bnnz; k++) PetscCall(MatSetValue(B, Bi[k], Bj[k], Bi[k] + Bj[k] + 2.0, INSERT_VALUES));
 49:   }
 50:   PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
 51:   PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));

 53:   PetscCall(MatMatMult(A, B, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &C));
 54:   PetscCall(MatView(C, PETSC_VIEWER_STDOUT_WORLD));

 56:   /* B, A have the same nonzero pattern, so it is legitimate to do so */
 57:   PetscCall(MatMatMult(B, A, MAT_REUSE_MATRIX, PETSC_DEFAULT, &C));
 58:   PetscCall(MatView(C, PETSC_VIEWER_STDOUT_WORLD));

 60:   PetscCall(MatTransposeMatMult(A, B, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &D));
 61:   PetscCall(MatView(D, PETSC_VIEWER_STDOUT_WORLD));

 63:   PetscCall(MatPtAP(A, B, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &E));
 64:   PetscCall(MatView(E, PETSC_VIEWER_STDOUT_WORLD));

 66:   PetscCall(MatDestroy(&A));
 67:   PetscCall(MatDestroy(&B));
 68:   PetscCall(MatDestroy(&C));
 69:   PetscCall(MatDestroy(&D));
 70:   PetscCall(MatDestroy(&E));

 72:   PetscCall(PetscFinalize());
 73:   return 0;
 74: }

 76: /*TEST
 77:   testset:
 78:     filter: grep -ve type -ve "Mat Object"
 79:     output_file: output/ex250_1.out

 81:     test:
 82:       suffix: 1
 83:       nsize: {{1 3}}
 84:       args: -mat_type aij

 86:     test:
 87:       suffix: 2
 88:       nsize: {{3 4}}
 89:       args: -mat_type aij -matmatmult_via backend -matptap_via backend -mattransposematmult_via backend

 91:     test:
 92:       suffix: cuda
 93:       requires: cuda
 94:       nsize: {{1 3 4}}
 95:       args: -mat_type aijcusparse

 97:     test:
 98:       suffix: kok
 99:       requires: kokkos_kernels
100:       nsize: {{1 3 4}}
101:       args: -mat_type aijkokkos

103: TEST*/