1 static char help[] = "Test MatMatMult(), MatTranspose(), MatTransposeMatMult() for Dense and Elemental matrices.\n\n"; 2 /* 3 Example: 4 mpiexec -n <np> ./ex104 -mat_type elemental 5 */ 6 7 #include <petscmat.h> 8 9 int main(int argc, char **argv) 10 { 11 Mat A, B, C, D; 12 PetscInt i, M = 10, N = 5, j, nrows, ncols, am, an, rstart, rend; 13 PetscRandom r; 14 PetscBool equal, Aiselemental; 15 PetscReal fill = 1.0; 16 IS isrows, iscols; 17 const PetscInt *rows, *cols; 18 PetscScalar *v, rval; 19 PetscBool Test_MatMatMult = PETSC_TRUE; 20 PetscMPIInt size; 21 22 PetscFunctionBeginUser; 23 PetscCall(PetscInitialize(&argc, &argv, (char *)0, help)); 24 PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 25 26 PetscCall(PetscOptionsGetInt(NULL, NULL, "-M", &M, NULL)); 27 PetscCall(PetscOptionsGetInt(NULL, NULL, "-N", &N, NULL)); 28 PetscCall(MatCreate(PETSC_COMM_WORLD, &A)); 29 PetscCall(MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, M, N)); 30 PetscCall(MatSetType(A, MATDENSE)); 31 PetscCall(MatSetFromOptions(A)); 32 PetscCall(MatSetUp(A)); 33 PetscCall(PetscRandomCreate(PETSC_COMM_WORLD, &r)); 34 PetscCall(PetscRandomSetFromOptions(r)); 35 36 /* Set local matrix entries */ 37 PetscCall(MatGetOwnershipIS(A, &isrows, &iscols)); 38 PetscCall(ISGetLocalSize(isrows, &nrows)); 39 PetscCall(ISGetIndices(isrows, &rows)); 40 PetscCall(ISGetLocalSize(iscols, &ncols)); 41 PetscCall(ISGetIndices(iscols, &cols)); 42 PetscCall(PetscMalloc1(nrows * ncols, &v)); 43 for (i = 0; i < nrows; i++) { 44 for (j = 0; j < ncols; j++) { 45 PetscCall(PetscRandomGetValue(r, &rval)); 46 v[i * ncols + j] = rval; 47 } 48 } 49 PetscCall(MatSetValues(A, nrows, rows, ncols, cols, v, INSERT_VALUES)); 50 PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 51 PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 52 PetscCall(ISRestoreIndices(isrows, &rows)); 53 PetscCall(ISRestoreIndices(iscols, &cols)); 54 PetscCall(ISDestroy(&isrows)); 55 PetscCall(ISDestroy(&iscols)); 56 PetscCall(PetscRandomDestroy(&r)); 57 58 PetscCall(PetscObjectTypeCompare((PetscObject)A, MATELEMENTAL, &Aiselemental)); 59 60 /* Test MatCreateTranspose() and MatTranspose() */ 61 PetscCall(MatCreateTranspose(A, &C)); 62 PetscCall(MatTranspose(A, MAT_INITIAL_MATRIX, &B)); /* B = A^T */ 63 PetscCall(MatMultEqual(C, B, 10, &equal)); 64 PetscCheck(equal, PETSC_COMM_SELF, PETSC_ERR_PLIB, "A^T*x != (x^T*A)^T"); 65 PetscCall(MatDestroy(&B)); 66 67 PetscCall(MatDuplicate(A, MAT_COPY_VALUES, &B)); 68 if (!Aiselemental) { 69 PetscCall(MatTranspose(B, MAT_INPLACE_MATRIX, &B)); 70 PetscCall(MatMultEqual(C, B, 10, &equal)); 71 PetscCheck(equal, PETSC_COMM_SELF, PETSC_ERR_PLIB, "C*x != B*x"); 72 } 73 PetscCall(MatDestroy(&B)); 74 75 /* Test B = C*A for matrix type transpose and seqdense */ 76 if (size == 1 && !Aiselemental) { 77 PetscCall(MatMatMult(C, A, MAT_INITIAL_MATRIX, fill, &B)); 78 PetscCall(MatMatMultEqual(C, A, B, 10, &equal)); 79 PetscCheck(equal, PETSC_COMM_SELF, PETSC_ERR_PLIB, "B != C*A for matrix type transpose and seqdense"); 80 PetscCall(MatDestroy(&B)); 81 } 82 PetscCall(MatDestroy(&C)); 83 84 /* Test MatMatMult() */ 85 if (Test_MatMatMult) { 86 PetscCall(MatTranspose(A, MAT_INITIAL_MATRIX, &B)); /* B = A^T */ 87 PetscCall(MatMatMult(B, A, MAT_INITIAL_MATRIX, fill, &C)); /* C = B*A = A^T*A */ 88 PetscCall(MatMatMult(B, A, MAT_REUSE_MATRIX, fill, &C)); 89 PetscCall(MatMatMultEqual(B, A, C, 10, &equal)); 90 PetscCheck(equal, PETSC_COMM_SELF, PETSC_ERR_PLIB, "B*A*x != C*x"); 91 92 /* Test MatDuplicate for matrix product */ 93 PetscCall(MatDuplicate(C, MAT_COPY_VALUES, &D)); 94 95 PetscCall(MatDestroy(&D)); 96 PetscCall(MatDestroy(&C)); 97 PetscCall(MatDestroy(&B)); 98 } 99 100 /* Test MatTransposeMatMult() */ 101 if (!Aiselemental) { 102 PetscCall(MatTransposeMatMult(A, A, MAT_INITIAL_MATRIX, fill, &D)); /* D = A^T*A */ 103 PetscCall(MatTransposeMatMult(A, A, MAT_REUSE_MATRIX, fill, &D)); 104 PetscCall(MatTransposeMatMultEqual(A, A, D, 10, &equal)); 105 PetscCheck(equal, PETSC_COMM_SELF, PETSC_ERR_PLIB, "D*x != A^T*A*x"); 106 107 /* Test MatDuplicate for matrix product */ 108 PetscCall(MatDuplicate(D, MAT_COPY_VALUES, &C)); 109 PetscCall(MatDestroy(&C)); 110 PetscCall(MatDestroy(&D)); 111 112 /* Test D*x = A^T*C*A*x, where C is in AIJ format */ 113 PetscCall(MatGetLocalSize(A, &am, &an)); 114 PetscCall(MatCreate(PETSC_COMM_WORLD, &C)); 115 if (size == 1) { 116 PetscCall(MatSetSizes(C, PETSC_DECIDE, PETSC_DECIDE, am, am)); 117 } else { 118 PetscCall(MatSetSizes(C, am, am, PETSC_DECIDE, PETSC_DECIDE)); 119 } 120 PetscCall(MatSetFromOptions(C)); 121 PetscCall(MatSetUp(C)); 122 PetscCall(MatGetOwnershipRange(C, &rstart, &rend)); 123 v[0] = 1.0; 124 for (i = rstart; i < rend; i++) PetscCall(MatSetValues(C, 1, &i, 1, &i, v, INSERT_VALUES)); 125 PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY)); 126 PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY)); 127 128 /* B = C*A, D = A^T*B */ 129 PetscCall(MatMatMult(C, A, MAT_INITIAL_MATRIX, 1.0, &B)); 130 PetscCall(MatTransposeMatMult(A, B, MAT_INITIAL_MATRIX, fill, &D)); 131 PetscCall(MatTransposeMatMultEqual(A, B, D, 10, &equal)); 132 PetscCheck(equal, PETSC_COMM_SELF, PETSC_ERR_PLIB, "D*x != A^T*B*x"); 133 134 PetscCall(MatDestroy(&D)); 135 PetscCall(MatDestroy(&C)); 136 PetscCall(MatDestroy(&B)); 137 } 138 139 /* Test MatMatTransposeMult() */ 140 if (!Aiselemental) { 141 PetscReal diff, scale; 142 PetscInt am, an, aM, aN; 143 144 PetscCall(MatGetLocalSize(A, &am, &an)); 145 PetscCall(MatGetSize(A, &aM, &aN)); 146 PetscCall(MatCreateDense(PetscObjectComm((PetscObject)A), PETSC_DECIDE, an, aM + 10, aN, NULL, &B)); 147 PetscCall(MatSetRandom(B, NULL)); 148 PetscCall(MatMatTransposeMult(A, B, MAT_INITIAL_MATRIX, fill, &D)); /* D = A*A^T */ 149 150 /* Test MatDuplicate for matrix product */ 151 PetscCall(MatDuplicate(D, MAT_COPY_VALUES, &C)); 152 153 PetscCall(MatMatTransposeMult(A, B, MAT_REUSE_MATRIX, fill, &D)); 154 PetscCall(MatAXPY(C, -1., D, SAME_NONZERO_PATTERN)); 155 156 PetscCall(MatNorm(C, NORM_FROBENIUS, &diff)); 157 PetscCall(MatNorm(D, NORM_FROBENIUS, &scale)); 158 PetscCheck(diff <= PETSC_SMALL * scale, PetscObjectComm((PetscObject)D), PETSC_ERR_PLIB, "MatMatTransposeMult() differs between MAT_INITIAL_MATRIX and MAT_REUSE_MATRIX"); 159 PetscCall(MatDestroy(&C)); 160 161 PetscCall(MatMatTransposeMultEqual(A, B, D, 10, &equal)); 162 PetscCheck(equal, PETSC_COMM_SELF, PETSC_ERR_PLIB, "D*x != A^T*A*x"); 163 PetscCall(MatDestroy(&D)); 164 PetscCall(MatDestroy(&B)); 165 } 166 167 PetscCall(MatDestroy(&A)); 168 PetscCall(PetscFree(v)); 169 PetscCall(PetscFinalize()); 170 return 0; 171 } 172 173 /*TEST 174 175 test: 176 output_file: output/ex104.out 177 178 test: 179 suffix: 2 180 nsize: 2 181 output_file: output/ex104.out 182 183 test: 184 suffix: 3 185 nsize: 4 186 output_file: output/ex104.out 187 args: -M 23 -N 31 188 189 test: 190 suffix: 4 191 nsize: 4 192 output_file: output/ex104.out 193 args: -M 23 -N 31 -matmattransmult_mpidense_mpidense_via cyclic 194 195 test: 196 suffix: 5 197 nsize: 4 198 output_file: output/ex104.out 199 args: -M 23 -N 31 -matmattransmult_mpidense_mpidense_via allgatherv 200 201 test: 202 suffix: 6 203 args: -mat_type elemental 204 requires: elemental 205 output_file: output/ex104.out 206 207 testset: 208 nsize: 2 209 output_file: output/ex104.out 210 requires: elemental 211 test: 212 suffix: 7_dense 213 args: -mat_type dense -mat_product_algorithm elemental 214 test: 215 suffix: 7_elemental 216 args: -mat_type elemental 217 218 TEST*/ 219