1 static char help[] = "Test MatAXPY \n\n"; 2 3 #include <petscmat.h> 4 int main(int argc, char **args) 5 { 6 Mat A = NULL, B = NULL; 7 PetscInt k; 8 const PetscInt M = 18, N = 18; 9 PetscMPIInt rank; 10 11 /* A, B are 18 x 18 nonsymmetric matrices. Initially, B's nonzero pattern is a subset of A. 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, 4, 4, 5, 5, 5, 6, 7, 7, 8, 8, 9, 9, 10, 10, 12, 12, 13, 13, 14, 14, 15, 15, 16, 16, 17, 17}; 17 PetscInt Bj[] = {0, 1, 2, 7, 3, 8, 5, 8, 2, 6, 11, 0, 1, 6, 2, 4, 10, 16, 11, 15, 12, 13, 14, 17, 11, 13, 3, 16, 9, 15, 11, 13}; 18 19 PetscInt Annz = PETSC_STATIC_ARRAY_LENGTH(Ai); 20 PetscInt Bnnz = PETSC_STATIC_ARRAY_LENGTH(Bi); 21 22 PetscFunctionBeginUser; 23 PetscCall(PetscInitialize(&argc, &args, (char *)0, help)); 24 PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 25 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)); 32 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 } 36 37 PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 38 PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 39 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)); 46 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)); 52 53 PetscCall(MatAXPY(A, 1.0, B, SUBSET_NONZERO_PATTERN)); /* A += B */ 54 PetscCall(MatView(A, PETSC_VIEWER_STDOUT_WORLD)); 55 56 /* Blow up B so it will have the same nonzero pattern as A */ 57 PetscCall(MatAXPY(B, 5, A, DIFFERENT_NONZERO_PATTERN)); /* B += 5A */ 58 PetscCall(MatView(B, PETSC_VIEWER_STDOUT_WORLD)); 59 60 PetscCall(MatAXPY(A, 3, B, SAME_NONZERO_PATTERN)); /* A += 7B */ 61 PetscCall(MatView(A, PETSC_VIEWER_STDOUT_WORLD)); 62 63 PetscCall(MatDestroy(&A)); 64 PetscCall(MatDestroy(&B)); 65 66 PetscCall(PetscFinalize()); 67 return 0; 68 } 69 70 /*TEST 71 testset: 72 nsize: {{1 2}} 73 filter: grep -ve type -ve "Mat Object" 74 output_file: output/ex251_1.out 75 76 test: 77 suffix: 1 78 args: -mat_type aij 79 80 test: 81 suffix: cuda 82 requires: cuda 83 args: -mat_type aijcusparse 84 85 test: 86 suffix: kok 87 requires: kokkos_kernels 88 args: -mat_type aijkokkos 89 90 TEST*/ 91