1c4762a1bSJed Brown static char help[] = "Test MatMult() for Hermitian matrix.\n\n"; 2c4762a1bSJed Brown 3c4762a1bSJed Brown #include <petscmat.h> 4c4762a1bSJed Brown 5*9371c9d4SSatish Balay int main(int argc, char **args) { 6c4762a1bSJed Brown Mat A, As; 7c4762a1bSJed Brown PetscBool flg; 8c4762a1bSJed Brown PetscMPIInt size; 9c4762a1bSJed Brown PetscInt i, j; 10c4762a1bSJed Brown PetscScalar v, sigma2; 11c4762a1bSJed Brown PetscReal h2, sigma1 = 100.0; 12c4762a1bSJed Brown PetscInt dim, Ii, J, n = 3, rstart, rend; 13c4762a1bSJed Brown 14327415f7SBarry Smith PetscFunctionBeginUser; 159566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &args, (char *)0, help)); 169566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 179566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL, NULL, "-sigma1", &sigma1, NULL)); 189566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-n", &n, NULL)); 19c4762a1bSJed Brown dim = n * n; 20c4762a1bSJed Brown 219566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &A)); 229566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, dim, dim)); 239566063dSJacob Faibussowitsch PetscCall(MatSetType(A, MATAIJ)); 249566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A)); 259566063dSJacob Faibussowitsch PetscCall(MatSetUp(A)); 26c4762a1bSJed Brown 27c4762a1bSJed Brown sigma2 = 10.0 * PETSC_i; 28c4762a1bSJed Brown h2 = 1.0 / ((n + 1) * (n + 1)); 29c4762a1bSJed Brown 309566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(A, &rstart, &rend)); 31c4762a1bSJed Brown for (Ii = rstart; Ii < rend; Ii++) { 32*9371c9d4SSatish Balay v = -1.0; 33*9371c9d4SSatish Balay i = Ii / n; 34*9371c9d4SSatish Balay j = Ii - i * n; 35c4762a1bSJed Brown if (i > 0) { 36*9371c9d4SSatish Balay J = Ii - n; 37*9371c9d4SSatish Balay PetscCall(MatSetValues(A, 1, &Ii, 1, &J, &v, ADD_VALUES)); 38c4762a1bSJed Brown } 39c4762a1bSJed Brown if (i < n - 1) { 40*9371c9d4SSatish Balay J = Ii + n; 41*9371c9d4SSatish Balay PetscCall(MatSetValues(A, 1, &Ii, 1, &J, &v, ADD_VALUES)); 42c4762a1bSJed Brown } 43c4762a1bSJed Brown if (j > 0) { 44*9371c9d4SSatish Balay J = Ii - 1; 45*9371c9d4SSatish Balay PetscCall(MatSetValues(A, 1, &Ii, 1, &J, &v, ADD_VALUES)); 46c4762a1bSJed Brown } 47c4762a1bSJed Brown if (j < n - 1) { 48*9371c9d4SSatish Balay J = Ii + 1; 49*9371c9d4SSatish Balay PetscCall(MatSetValues(A, 1, &Ii, 1, &J, &v, ADD_VALUES)); 50c4762a1bSJed Brown } 51c4762a1bSJed Brown v = 4.0 - sigma1 * h2; 529566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &Ii, 1, &Ii, &v, ADD_VALUES)); 53c4762a1bSJed Brown } 549566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 559566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 56c4762a1bSJed Brown 57c4762a1bSJed Brown /* Check whether A is symmetric */ 589566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-check_symmetric", &flg)); 59c4762a1bSJed Brown if (flg) { 609566063dSJacob Faibussowitsch PetscCall(MatIsSymmetric(A, 0.0, &flg)); 6128b400f6SJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_USER, "A is not symmetric"); 62c4762a1bSJed Brown } 639566063dSJacob Faibussowitsch PetscCall(MatSetOption(A, MAT_SYMMETRIC, PETSC_TRUE)); 64c4762a1bSJed Brown 65c4762a1bSJed Brown /* make A complex Hermitian */ 66*9371c9d4SSatish Balay Ii = 0; 67*9371c9d4SSatish Balay J = dim - 1; 68c4762a1bSJed Brown if (Ii >= rstart && Ii < rend) { 69c4762a1bSJed Brown v = sigma2 * h2; /* RealPart(v) = 0.0 */ 709566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &Ii, 1, &J, &v, ADD_VALUES)); 71c4762a1bSJed Brown v = -sigma2 * h2; 729566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &J, 1, &Ii, &v, ADD_VALUES)); 73c4762a1bSJed Brown } 74c4762a1bSJed Brown 75*9371c9d4SSatish Balay Ii = dim - 2; 76*9371c9d4SSatish Balay J = dim - 1; 77c4762a1bSJed Brown if (Ii >= rstart && Ii < rend) { 78c4762a1bSJed Brown v = sigma2 * h2; /* RealPart(v) = 0.0 */ 799566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &Ii, 1, &J, &v, ADD_VALUES)); 80c4762a1bSJed Brown v = -sigma2 * h2; 819566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &J, 1, &Ii, &v, ADD_VALUES)); 82c4762a1bSJed Brown } 83c4762a1bSJed Brown 849566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 859566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 869566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(A, NULL, "-disp_mat")); 87c4762a1bSJed Brown 88c4762a1bSJed Brown /* Check whether A is Hermitian, then set A->hermitian flag */ 899566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-check_Hermitian", &flg)); 90c4762a1bSJed Brown if (flg && size == 1) { 919566063dSJacob Faibussowitsch PetscCall(MatIsHermitian(A, 0.0, &flg)); 9228b400f6SJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_USER, "A is not Hermitian"); 93c4762a1bSJed Brown } 949566063dSJacob Faibussowitsch PetscCall(MatSetOption(A, MAT_HERMITIAN, PETSC_TRUE)); 95c4762a1bSJed Brown 96c4762a1bSJed Brown #if defined(PETSC_HAVE_SUPERLU_DIST) 97c4762a1bSJed Brown /* Test Cholesky factorization */ 989566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-test_choleskyfactor", &flg)); 99c4762a1bSJed Brown if (flg) { 100c4762a1bSJed Brown Mat F; 101c4762a1bSJed Brown IS perm, iperm; 102c4762a1bSJed Brown MatFactorInfo info; 103c4762a1bSJed Brown PetscInt nneg, nzero, npos; 104c4762a1bSJed Brown 1059566063dSJacob Faibussowitsch PetscCall(MatGetFactor(A, MATSOLVERSUPERLU_DIST, MAT_FACTOR_CHOLESKY, &F)); 1069566063dSJacob Faibussowitsch PetscCall(MatGetOrdering(A, MATORDERINGND, &perm, &iperm)); 1079566063dSJacob Faibussowitsch PetscCall(MatCholeskyFactorSymbolic(F, A, perm, &info)); 1089566063dSJacob Faibussowitsch PetscCall(MatCholeskyFactorNumeric(F, A, &info)); 109c4762a1bSJed Brown 1109566063dSJacob Faibussowitsch PetscCall(MatGetInertia(F, &nneg, &nzero, &npos)); 1119566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, " MatInertia: nneg: %" PetscInt_FMT ", nzero: %" PetscInt_FMT ", npos: %" PetscInt_FMT "\n", nneg, nzero, npos)); 1129566063dSJacob Faibussowitsch PetscCall(MatDestroy(&F)); 1139566063dSJacob Faibussowitsch PetscCall(ISDestroy(&perm)); 1149566063dSJacob Faibussowitsch PetscCall(ISDestroy(&iperm)); 115c4762a1bSJed Brown } 116c4762a1bSJed Brown #endif 117c4762a1bSJed Brown 118c4762a1bSJed Brown /* Create a Hermitian matrix As in sbaij format */ 1199566063dSJacob Faibussowitsch PetscCall(MatConvert(A, MATSBAIJ, MAT_INITIAL_MATRIX, &As)); 1209566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(As, NULL, "-disp_mat")); 121c4762a1bSJed Brown 122c4762a1bSJed Brown /* Test MatMult */ 1239566063dSJacob Faibussowitsch PetscCall(MatMultEqual(A, As, 10, &flg)); 12428b400f6SJacob Faibussowitsch PetscCheck(flg, PetscObjectComm((PetscObject)A), PETSC_ERR_PLIB, "MatMult not equal"); 1259566063dSJacob Faibussowitsch PetscCall(MatMultAddEqual(A, As, 10, &flg)); 12628b400f6SJacob Faibussowitsch PetscCheck(flg, PetscObjectComm((PetscObject)A), PETSC_ERR_PLIB, "MatMultAdd not equal"); 127c4762a1bSJed Brown 128c4762a1bSJed Brown /* Free spaces */ 1299566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 1309566063dSJacob Faibussowitsch PetscCall(MatDestroy(&As)); 1319566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 132b122ec5aSJacob Faibussowitsch return 0; 133c4762a1bSJed Brown } 134c4762a1bSJed Brown 135c4762a1bSJed Brown /*TEST 136c4762a1bSJed Brown 137c4762a1bSJed Brown build: 138c4762a1bSJed Brown requires: complex 139c4762a1bSJed Brown 140c4762a1bSJed Brown test: 141c4762a1bSJed Brown args: -n 1000 142c4762a1bSJed Brown output_file: output/ex127.out 143c4762a1bSJed Brown 144c4762a1bSJed Brown test: 145c4762a1bSJed Brown suffix: 2 146c4762a1bSJed Brown nsize: 3 147c4762a1bSJed Brown args: -n 1000 148c4762a1bSJed Brown output_file: output/ex127.out 149c4762a1bSJed Brown 150c4762a1bSJed Brown test: 151c4762a1bSJed Brown suffix: superlu_dist 152c4762a1bSJed Brown nsize: 3 153c4762a1bSJed Brown requires: superlu_dist 154c4762a1bSJed Brown args: -test_choleskyfactor -mat_superlu_dist_rowperm NOROWPERM 155c4762a1bSJed Brown output_file: output/ex127_superlu_dist.out 156c4762a1bSJed Brown TEST*/ 157