xref: /petsc/src/mat/tests/ex127.c (revision 8bb6a241cc2b8f6ed3ced698b43bc00095eff8a0)
1c4762a1bSJed Brown static char help[] = "Test MatMult() for Hermitian matrix.\n\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown #include <petscmat.h>
4c4762a1bSJed Brown 
main(int argc,char ** args)5d71ae5a4SJacob Faibussowitsch int main(int argc, char **args)
6d71ae5a4SJacob Faibussowitsch {
7*e9b284a9SPierre Jolivet   Mat         A, As, B;
8*e9b284a9SPierre Jolivet   Vec         l, r;
9c4762a1bSJed Brown   PetscBool   flg;
10c4762a1bSJed Brown   PetscMPIInt size;
11c4762a1bSJed Brown   PetscInt    i, j;
12c4762a1bSJed Brown   PetscScalar v, sigma2;
13c4762a1bSJed Brown   PetscReal   h2, sigma1 = 100.0;
14c4762a1bSJed Brown   PetscInt    dim, Ii, J, n = 3, rstart, rend;
15c4762a1bSJed Brown 
16327415f7SBarry Smith   PetscFunctionBeginUser;
17c8025a54SPierre Jolivet   PetscCall(PetscInitialize(&argc, &args, NULL, help));
189566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
199566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL, NULL, "-sigma1", &sigma1, NULL));
209566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-n", &n, NULL));
21c4762a1bSJed Brown   dim = n * n;
22c4762a1bSJed Brown 
239566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
249566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, dim, dim));
259566063dSJacob Faibussowitsch   PetscCall(MatSetType(A, MATAIJ));
269566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(A));
279566063dSJacob Faibussowitsch   PetscCall(MatSetUp(A));
28c4762a1bSJed Brown 
29c4762a1bSJed Brown   sigma2 = 10.0 * PETSC_i;
30c4762a1bSJed Brown   h2     = 1.0 / ((n + 1) * (n + 1));
31c4762a1bSJed Brown 
329566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(A, &rstart, &rend));
33c4762a1bSJed Brown   for (Ii = rstart; Ii < rend; Ii++) {
349371c9d4SSatish Balay     v = -1.0;
359371c9d4SSatish Balay     i = Ii / n;
369371c9d4SSatish Balay     j = Ii - i * n;
37c4762a1bSJed Brown     if (i > 0) {
389371c9d4SSatish Balay       J = Ii - n;
399371c9d4SSatish Balay       PetscCall(MatSetValues(A, 1, &Ii, 1, &J, &v, ADD_VALUES));
40c4762a1bSJed Brown     }
41c4762a1bSJed Brown     if (i < n - 1) {
429371c9d4SSatish Balay       J = Ii + n;
439371c9d4SSatish Balay       PetscCall(MatSetValues(A, 1, &Ii, 1, &J, &v, ADD_VALUES));
44c4762a1bSJed Brown     }
45c4762a1bSJed Brown     if (j > 0) {
469371c9d4SSatish Balay       J = Ii - 1;
479371c9d4SSatish Balay       PetscCall(MatSetValues(A, 1, &Ii, 1, &J, &v, ADD_VALUES));
48c4762a1bSJed Brown     }
49c4762a1bSJed Brown     if (j < n - 1) {
509371c9d4SSatish Balay       J = Ii + 1;
519371c9d4SSatish Balay       PetscCall(MatSetValues(A, 1, &Ii, 1, &J, &v, ADD_VALUES));
52c4762a1bSJed Brown     }
53c4762a1bSJed Brown     v = 4.0 - sigma1 * h2;
549566063dSJacob Faibussowitsch     PetscCall(MatSetValues(A, 1, &Ii, 1, &Ii, &v, ADD_VALUES));
55c4762a1bSJed Brown   }
569566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
579566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
58c4762a1bSJed Brown 
59c4762a1bSJed Brown   /* Check whether A is symmetric */
609566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-check_symmetric", &flg));
61c4762a1bSJed Brown   if (flg) {
629566063dSJacob Faibussowitsch     PetscCall(MatIsSymmetric(A, 0.0, &flg));
63*e9b284a9SPierre Jolivet     PetscCheck(flg, PETSC_COMM_WORLD, PETSC_ERR_USER, "A is not symmetric");
64c4762a1bSJed Brown   }
659566063dSJacob Faibussowitsch   PetscCall(MatSetOption(A, MAT_SYMMETRIC, PETSC_TRUE));
66c4762a1bSJed Brown 
67c4762a1bSJed Brown   /* make A complex Hermitian */
689371c9d4SSatish Balay   Ii = 0;
699371c9d4SSatish Balay   J  = dim - 1;
70c4762a1bSJed Brown   if (Ii >= rstart && Ii < rend) {
71c4762a1bSJed Brown     v = sigma2 * h2; /* RealPart(v) = 0.0 */
729566063dSJacob Faibussowitsch     PetscCall(MatSetValues(A, 1, &Ii, 1, &J, &v, ADD_VALUES));
73c4762a1bSJed Brown     v = -sigma2 * h2;
749566063dSJacob Faibussowitsch     PetscCall(MatSetValues(A, 1, &J, 1, &Ii, &v, ADD_VALUES));
75c4762a1bSJed Brown   }
76c4762a1bSJed Brown 
779371c9d4SSatish Balay   Ii = dim - 2;
789371c9d4SSatish Balay   J  = dim - 1;
79c4762a1bSJed Brown   if (Ii >= rstart && Ii < rend) {
80c4762a1bSJed Brown     v = sigma2 * h2; /* RealPart(v) = 0.0 */
819566063dSJacob Faibussowitsch     PetscCall(MatSetValues(A, 1, &Ii, 1, &J, &v, ADD_VALUES));
82c4762a1bSJed Brown     v = -sigma2 * h2;
839566063dSJacob Faibussowitsch     PetscCall(MatSetValues(A, 1, &J, 1, &Ii, &v, ADD_VALUES));
84c4762a1bSJed Brown   }
85c4762a1bSJed Brown 
869566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
879566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
889566063dSJacob Faibussowitsch   PetscCall(MatViewFromOptions(A, NULL, "-disp_mat"));
89c4762a1bSJed Brown 
90*e9b284a9SPierre Jolivet   if (flg) {
91*e9b284a9SPierre Jolivet     PetscCall(MatIsSymmetric(A, 0.0, &flg));
92*e9b284a9SPierre Jolivet     PetscCheck(!flg, PETSC_COMM_WORLD, PETSC_ERR_USER, "A is symmetric");
93*e9b284a9SPierre Jolivet   }
94*e9b284a9SPierre Jolivet   PetscCall(MatSetOption(A, MAT_SYMMETRIC, PETSC_FALSE));
95c4762a1bSJed Brown   /* Check whether A is Hermitian, then set A->hermitian flag */
969566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-check_Hermitian", &flg));
97c4762a1bSJed Brown   if (flg && size == 1) {
989566063dSJacob Faibussowitsch     PetscCall(MatIsHermitian(A, 0.0, &flg));
99*e9b284a9SPierre Jolivet     PetscCheck(flg, PETSC_COMM_WORLD, PETSC_ERR_USER, "A is not Hermitian");
100c4762a1bSJed Brown   }
1019566063dSJacob Faibussowitsch   PetscCall(MatSetOption(A, MAT_HERMITIAN, PETSC_TRUE));
102c4762a1bSJed Brown 
103c4762a1bSJed Brown #if defined(PETSC_HAVE_SUPERLU_DIST)
104c4762a1bSJed Brown   /* Test Cholesky factorization */
1059566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-test_choleskyfactor", &flg));
106c4762a1bSJed Brown   if (flg) {
107c4762a1bSJed Brown     Mat           F;
108c4762a1bSJed Brown     IS            perm, iperm;
109c4762a1bSJed Brown     MatFactorInfo info;
110c4762a1bSJed Brown     PetscInt      nneg, nzero, npos;
111c4762a1bSJed Brown 
1129566063dSJacob Faibussowitsch     PetscCall(MatGetFactor(A, MATSOLVERSUPERLU_DIST, MAT_FACTOR_CHOLESKY, &F));
1139566063dSJacob Faibussowitsch     PetscCall(MatGetOrdering(A, MATORDERINGND, &perm, &iperm));
1149566063dSJacob Faibussowitsch     PetscCall(MatCholeskyFactorSymbolic(F, A, perm, &info));
1159566063dSJacob Faibussowitsch     PetscCall(MatCholeskyFactorNumeric(F, A, &info));
116c4762a1bSJed Brown 
1179566063dSJacob Faibussowitsch     PetscCall(MatGetInertia(F, &nneg, &nzero, &npos));
1189566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, " MatInertia: nneg: %" PetscInt_FMT ", nzero: %" PetscInt_FMT ", npos: %" PetscInt_FMT "\n", nneg, nzero, npos));
1199566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&F));
1209566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&perm));
1219566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&iperm));
122c4762a1bSJed Brown   }
123c4762a1bSJed Brown #endif
124c4762a1bSJed Brown 
125c4762a1bSJed Brown   /* Create a Hermitian matrix As in sbaij format */
1269566063dSJacob Faibussowitsch   PetscCall(MatConvert(A, MATSBAIJ, MAT_INITIAL_MATRIX, &As));
1279566063dSJacob Faibussowitsch   PetscCall(MatViewFromOptions(As, NULL, "-disp_mat"));
128c4762a1bSJed Brown 
129c4762a1bSJed Brown   /* Test MatMult */
1309566063dSJacob Faibussowitsch   PetscCall(MatMultEqual(A, As, 10, &flg));
13128b400f6SJacob Faibussowitsch   PetscCheck(flg, PetscObjectComm((PetscObject)A), PETSC_ERR_PLIB, "MatMult not equal");
1329566063dSJacob Faibussowitsch   PetscCall(MatMultAddEqual(A, As, 10, &flg));
13328b400f6SJacob Faibussowitsch   PetscCheck(flg, PetscObjectComm((PetscObject)A), PETSC_ERR_PLIB, "MatMultAdd not equal");
134c4762a1bSJed Brown 
135*e9b284a9SPierre Jolivet   PetscCall(MatDuplicate(A, MAT_COPY_VALUES, &B));
136*e9b284a9SPierre Jolivet   PetscCall(MatRealPart(B));
137*e9b284a9SPierre Jolivet   PetscCall(MatSetOption(B, MAT_SYMMETRIC, PETSC_TRUE));
138*e9b284a9SPierre Jolivet 
139*e9b284a9SPierre Jolivet   PetscCall(MatCreateVecs(A, &r, &l));
140*e9b284a9SPierre Jolivet   PetscCall(VecSetRandom(r, NULL));
141*e9b284a9SPierre Jolivet   PetscCall(VecCopy(r, l));
142*e9b284a9SPierre Jolivet   PetscCall(MatDiagonalScale(A, r, l));
143*e9b284a9SPierre Jolivet   PetscCall(PetscOptionsHasName(NULL, NULL, "-check_Hermitian", &flg));
144*e9b284a9SPierre Jolivet   if (flg && size == 1) {
145*e9b284a9SPierre Jolivet     PetscCall(MatIsHermitian(A, 0.0, &flg));
146*e9b284a9SPierre Jolivet     PetscCheck(!flg, PETSC_COMM_WORLD, PETSC_ERR_USER, "A is Hermitian");
147*e9b284a9SPierre Jolivet   }
148*e9b284a9SPierre Jolivet   PetscCall(MatSetOption(A, MAT_HERMITIAN, PETSC_FALSE));
149*e9b284a9SPierre Jolivet 
150c4762a1bSJed Brown   /* Free spaces */
1519566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
1529566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&As));
153*e9b284a9SPierre Jolivet   PetscCall(MatDestroy(&B));
154*e9b284a9SPierre Jolivet   PetscCall(VecDestroy(&r));
155*e9b284a9SPierre Jolivet   PetscCall(VecDestroy(&l));
1569566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
157b122ec5aSJacob Faibussowitsch   return 0;
158c4762a1bSJed Brown }
159c4762a1bSJed Brown 
160c4762a1bSJed Brown /*TEST
161c4762a1bSJed Brown 
162c4762a1bSJed Brown    build:
163c4762a1bSJed Brown       requires: complex
164c4762a1bSJed Brown 
165c4762a1bSJed Brown    test:
166*e9b284a9SPierre Jolivet       args: -n 1000 -check_symmetric -check_Hermitian
1673886731fSPierre Jolivet       output_file: output/empty.out
168c4762a1bSJed Brown 
169c4762a1bSJed Brown    test:
170c4762a1bSJed Brown       suffix: 2
171c4762a1bSJed Brown       nsize: 3
172c4762a1bSJed Brown       args: -n 1000
1733886731fSPierre Jolivet       output_file: output/empty.out
174c4762a1bSJed Brown 
175c4762a1bSJed Brown    test:
176c4762a1bSJed Brown       suffix: superlu_dist
177c4762a1bSJed Brown       nsize: 3
178c4762a1bSJed Brown       requires: superlu_dist
179c4762a1bSJed Brown       args: -test_choleskyfactor -mat_superlu_dist_rowperm NOROWPERM
180c4762a1bSJed Brown       output_file: output/ex127_superlu_dist.out
181c4762a1bSJed Brown TEST*/
182