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