18c78258cSHong Zhang 28c78258cSHong Zhang static char help[] = "Tests MatOption MAT_FORCE_DIAGONAL_ENTRIES.\n\n"; 38c78258cSHong Zhang 48c78258cSHong Zhang #include <petscmat.h> 58c78258cSHong Zhang 6*d71ae5a4SJacob Faibussowitsch int main(int argc, char **args) 7*d71ae5a4SJacob Faibussowitsch { 88c78258cSHong Zhang Mat A, B; 98c78258cSHong Zhang Vec diag; 108c78258cSHong Zhang PetscInt i, n = 10, col[3], test; 118c78258cSHong Zhang PetscScalar v[3]; 128c78258cSHong Zhang 13327415f7SBarry Smith PetscFunctionBeginUser; 149566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &args, (char *)0, help)); 159566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-n", &n, NULL)); 168c78258cSHong Zhang 178c78258cSHong Zhang /* Create A which has empty 0-th row and column */ 189566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &A)); 199566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, n, n)); 209566063dSJacob Faibussowitsch PetscCall(MatSetType(A, MATAIJ)); 219566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A)); 229566063dSJacob Faibussowitsch PetscCall(MatSetUp(A)); 238c78258cSHong Zhang 249371c9d4SSatish Balay v[0] = -1.; 259371c9d4SSatish Balay v[1] = 2.; 269371c9d4SSatish Balay v[2] = -1.; 278c78258cSHong Zhang for (i = 2; i < n - 1; i++) { 289371c9d4SSatish Balay col[0] = i - 1; 299371c9d4SSatish Balay col[1] = i; 309371c9d4SSatish Balay col[2] = i + 1; 319566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &i, 3, col, v, INSERT_VALUES)); 328c78258cSHong Zhang } 339371c9d4SSatish Balay i = 1; 349371c9d4SSatish Balay col[0] = 1; 359371c9d4SSatish Balay col[1] = 2; 369566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &i, 2, col, v + 1, INSERT_VALUES)); 379371c9d4SSatish Balay i = n - 1; 389371c9d4SSatish Balay col[0] = n - 2; 399371c9d4SSatish Balay col[1] = n - 1; 409566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &i, 2, col, v, INSERT_VALUES)); 419566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 429566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 438c78258cSHong Zhang 448c78258cSHong Zhang for (test = 0; test < 2; test++) { 459566063dSJacob Faibussowitsch PetscCall(MatProductCreate(A, A, NULL, &B)); 468c78258cSHong Zhang 478c78258cSHong Zhang if (test == 0) { 488c78258cSHong Zhang /* Compute B = A*A; B misses 0-th diagonal */ 499566063dSJacob Faibussowitsch PetscCall(MatProductSetType(B, MATPRODUCT_AB)); 509566063dSJacob Faibussowitsch PetscCall(MatSetOptionsPrefix(B, "AB_")); 518c78258cSHong Zhang } else { 528c78258cSHong Zhang /* Compute B = A^t*A; B misses 0-th diagonal */ 539566063dSJacob Faibussowitsch PetscCall(MatProductSetType(B, MATPRODUCT_AtB)); 549566063dSJacob Faibussowitsch PetscCall(MatSetOptionsPrefix(B, "AtB_")); 558c78258cSHong Zhang } 568c78258cSHong Zhang 578c78258cSHong Zhang /* Force allocate missing diagonal entries of B */ 589566063dSJacob Faibussowitsch PetscCall(MatSetOption(B, MAT_FORCE_DIAGONAL_ENTRIES, PETSC_TRUE)); 599566063dSJacob Faibussowitsch PetscCall(MatProductSetFromOptions(B)); 608c78258cSHong Zhang 619566063dSJacob Faibussowitsch PetscCall(MatProductSymbolic(B)); 629566063dSJacob Faibussowitsch PetscCall(MatProductNumeric(B)); 638c78258cSHong Zhang 649566063dSJacob Faibussowitsch PetscCall(MatSetOption(B, MAT_NEW_NONZERO_LOCATION_ERR, PETSC_TRUE)); 658c78258cSHong Zhang 668c78258cSHong Zhang /* Insert entries to diagonal of B */ 679566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(B, NULL, &diag)); 689566063dSJacob Faibussowitsch PetscCall(MatGetDiagonal(B, diag)); 699566063dSJacob Faibussowitsch PetscCall(VecSetValue(diag, 0, 100.0, INSERT_VALUES)); 709566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(diag)); 719566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(diag)); 728c78258cSHong Zhang 739566063dSJacob Faibussowitsch PetscCall(MatDiagonalSet(B, diag, INSERT_VALUES)); 7448a46eb9SPierre Jolivet if (test == 1) PetscCall(MatView(B, PETSC_VIEWER_STDOUT_WORLD)); 759566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 769566063dSJacob Faibussowitsch PetscCall(VecDestroy(&diag)); 778c78258cSHong Zhang } 788c78258cSHong Zhang 799566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 809566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 81b122ec5aSJacob Faibussowitsch return 0; 828c78258cSHong Zhang } 838c78258cSHong Zhang 848c78258cSHong Zhang /*TEST 858c78258cSHong Zhang 868c78258cSHong Zhang test: 878c78258cSHong Zhang output_file: output/ex81_1.out 888c78258cSHong Zhang 898c78258cSHong Zhang test: 908c78258cSHong Zhang suffix: 2 913e662e0bSHong Zhang args: -AtB_mat_product_algorithm at*b 928c78258cSHong Zhang output_file: output/ex81_1.out 938c78258cSHong Zhang 948c78258cSHong Zhang test: 958c78258cSHong Zhang suffix: 3 963e662e0bSHong Zhang args: -AtB_mat_product_algorithm outerproduct 978c78258cSHong Zhang output_file: output/ex81_1.out 988c78258cSHong Zhang 998c78258cSHong Zhang test: 1008c78258cSHong Zhang suffix: 4 1018c78258cSHong Zhang nsize: 3 1023e662e0bSHong Zhang args: -AtB_mat_product_algorithm nonscalable 1038c78258cSHong Zhang output_file: output/ex81_3.out 1048c78258cSHong Zhang 1058c78258cSHong Zhang test: 1068c78258cSHong Zhang suffix: 5 1078c78258cSHong Zhang nsize: 3 1083e662e0bSHong Zhang args: -AtB_mat_product_algorithm scalable 1098c78258cSHong Zhang output_file: output/ex81_3.out 1108c78258cSHong Zhang 1118c78258cSHong Zhang test: 1128c78258cSHong Zhang suffix: 6 1138c78258cSHong Zhang nsize: 3 1143e662e0bSHong Zhang args: -AtB_mat_product_algorithm at*b 1158c78258cSHong Zhang output_file: output/ex81_3.out 1168c78258cSHong Zhang 1178c78258cSHong Zhang TEST*/ 118