xref: /petsc/src/mat/tests/ex127.c (revision 3886731f808f6f9ca8348e6e791c818340d72f91)
1c4762a1bSJed Brown static char help[] = "Test MatMult() for Hermitian matrix.\n\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown #include <petscmat.h>
4c4762a1bSJed Brown 
5d71ae5a4SJacob Faibussowitsch int main(int argc, char **args)
6d71ae5a4SJacob Faibussowitsch {
7c4762a1bSJed Brown   Mat         A, As;
8c4762a1bSJed Brown   PetscBool   flg;
9c4762a1bSJed Brown   PetscMPIInt size;
10c4762a1bSJed Brown   PetscInt    i, j;
11c4762a1bSJed Brown   PetscScalar v, sigma2;
12c4762a1bSJed Brown   PetscReal   h2, sigma1 = 100.0;
13c4762a1bSJed Brown   PetscInt    dim, Ii, J, n = 3, rstart, rend;
14c4762a1bSJed Brown 
15327415f7SBarry Smith   PetscFunctionBeginUser;
16c8025a54SPierre Jolivet   PetscCall(PetscInitialize(&argc, &args, NULL, help));
179566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
189566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL, NULL, "-sigma1", &sigma1, NULL));
199566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-n", &n, NULL));
20c4762a1bSJed Brown   dim = n * n;
21c4762a1bSJed Brown 
229566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
239566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, dim, dim));
249566063dSJacob Faibussowitsch   PetscCall(MatSetType(A, MATAIJ));
259566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(A));
269566063dSJacob Faibussowitsch   PetscCall(MatSetUp(A));
27c4762a1bSJed Brown 
28c4762a1bSJed Brown   sigma2 = 10.0 * PETSC_i;
29c4762a1bSJed Brown   h2     = 1.0 / ((n + 1) * (n + 1));
30c4762a1bSJed Brown 
319566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(A, &rstart, &rend));
32c4762a1bSJed Brown   for (Ii = rstart; Ii < rend; Ii++) {
339371c9d4SSatish Balay     v = -1.0;
349371c9d4SSatish Balay     i = Ii / n;
359371c9d4SSatish Balay     j = Ii - i * n;
36c4762a1bSJed Brown     if (i > 0) {
379371c9d4SSatish Balay       J = Ii - n;
389371c9d4SSatish Balay       PetscCall(MatSetValues(A, 1, &Ii, 1, &J, &v, ADD_VALUES));
39c4762a1bSJed Brown     }
40c4762a1bSJed Brown     if (i < n - 1) {
419371c9d4SSatish Balay       J = Ii + n;
429371c9d4SSatish Balay       PetscCall(MatSetValues(A, 1, &Ii, 1, &J, &v, ADD_VALUES));
43c4762a1bSJed Brown     }
44c4762a1bSJed Brown     if (j > 0) {
459371c9d4SSatish Balay       J = Ii - 1;
469371c9d4SSatish Balay       PetscCall(MatSetValues(A, 1, &Ii, 1, &J, &v, ADD_VALUES));
47c4762a1bSJed Brown     }
48c4762a1bSJed Brown     if (j < n - 1) {
499371c9d4SSatish Balay       J = Ii + 1;
509371c9d4SSatish Balay       PetscCall(MatSetValues(A, 1, &Ii, 1, &J, &v, ADD_VALUES));
51c4762a1bSJed Brown     }
52c4762a1bSJed Brown     v = 4.0 - sigma1 * h2;
539566063dSJacob Faibussowitsch     PetscCall(MatSetValues(A, 1, &Ii, 1, &Ii, &v, ADD_VALUES));
54c4762a1bSJed Brown   }
559566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
569566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
57c4762a1bSJed Brown 
58c4762a1bSJed Brown   /* Check whether A is symmetric */
599566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-check_symmetric", &flg));
60c4762a1bSJed Brown   if (flg) {
619566063dSJacob Faibussowitsch     PetscCall(MatIsSymmetric(A, 0.0, &flg));
6228b400f6SJacob Faibussowitsch     PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_USER, "A is not symmetric");
63c4762a1bSJed Brown   }
649566063dSJacob Faibussowitsch   PetscCall(MatSetOption(A, MAT_SYMMETRIC, PETSC_TRUE));
65c4762a1bSJed Brown 
66c4762a1bSJed Brown   /* make A complex Hermitian */
679371c9d4SSatish Balay   Ii = 0;
689371c9d4SSatish Balay   J  = dim - 1;
69c4762a1bSJed Brown   if (Ii >= rstart && Ii < rend) {
70c4762a1bSJed Brown     v = sigma2 * h2; /* RealPart(v) = 0.0 */
719566063dSJacob Faibussowitsch     PetscCall(MatSetValues(A, 1, &Ii, 1, &J, &v, ADD_VALUES));
72c4762a1bSJed Brown     v = -sigma2 * h2;
739566063dSJacob Faibussowitsch     PetscCall(MatSetValues(A, 1, &J, 1, &Ii, &v, ADD_VALUES));
74c4762a1bSJed Brown   }
75c4762a1bSJed Brown 
769371c9d4SSatish Balay   Ii = dim - 2;
779371c9d4SSatish Balay   J  = dim - 1;
78c4762a1bSJed Brown   if (Ii >= rstart && Ii < rend) {
79c4762a1bSJed Brown     v = sigma2 * h2; /* RealPart(v) = 0.0 */
809566063dSJacob Faibussowitsch     PetscCall(MatSetValues(A, 1, &Ii, 1, &J, &v, ADD_VALUES));
81c4762a1bSJed Brown     v = -sigma2 * h2;
829566063dSJacob Faibussowitsch     PetscCall(MatSetValues(A, 1, &J, 1, &Ii, &v, ADD_VALUES));
83c4762a1bSJed Brown   }
84c4762a1bSJed Brown 
859566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
869566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
879566063dSJacob Faibussowitsch   PetscCall(MatViewFromOptions(A, NULL, "-disp_mat"));
88c4762a1bSJed Brown 
89c4762a1bSJed Brown   /* Check whether A is Hermitian, then set A->hermitian flag */
909566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-check_Hermitian", &flg));
91c4762a1bSJed Brown   if (flg && size == 1) {
929566063dSJacob Faibussowitsch     PetscCall(MatIsHermitian(A, 0.0, &flg));
9328b400f6SJacob Faibussowitsch     PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_USER, "A is not Hermitian");
94c4762a1bSJed Brown   }
959566063dSJacob Faibussowitsch   PetscCall(MatSetOption(A, MAT_HERMITIAN, PETSC_TRUE));
96c4762a1bSJed Brown 
97c4762a1bSJed Brown #if defined(PETSC_HAVE_SUPERLU_DIST)
98c4762a1bSJed Brown   /* Test Cholesky factorization */
999566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-test_choleskyfactor", &flg));
100c4762a1bSJed Brown   if (flg) {
101c4762a1bSJed Brown     Mat           F;
102c4762a1bSJed Brown     IS            perm, iperm;
103c4762a1bSJed Brown     MatFactorInfo info;
104c4762a1bSJed Brown     PetscInt      nneg, nzero, npos;
105c4762a1bSJed Brown 
1069566063dSJacob Faibussowitsch     PetscCall(MatGetFactor(A, MATSOLVERSUPERLU_DIST, MAT_FACTOR_CHOLESKY, &F));
1079566063dSJacob Faibussowitsch     PetscCall(MatGetOrdering(A, MATORDERINGND, &perm, &iperm));
1089566063dSJacob Faibussowitsch     PetscCall(MatCholeskyFactorSymbolic(F, A, perm, &info));
1099566063dSJacob Faibussowitsch     PetscCall(MatCholeskyFactorNumeric(F, A, &info));
110c4762a1bSJed Brown 
1119566063dSJacob Faibussowitsch     PetscCall(MatGetInertia(F, &nneg, &nzero, &npos));
1129566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, " MatInertia: nneg: %" PetscInt_FMT ", nzero: %" PetscInt_FMT ", npos: %" PetscInt_FMT "\n", nneg, nzero, npos));
1139566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&F));
1149566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&perm));
1159566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&iperm));
116c4762a1bSJed Brown   }
117c4762a1bSJed Brown #endif
118c4762a1bSJed Brown 
119c4762a1bSJed Brown   /* Create a Hermitian matrix As in sbaij format */
1209566063dSJacob Faibussowitsch   PetscCall(MatConvert(A, MATSBAIJ, MAT_INITIAL_MATRIX, &As));
1219566063dSJacob Faibussowitsch   PetscCall(MatViewFromOptions(As, NULL, "-disp_mat"));
122c4762a1bSJed Brown 
123c4762a1bSJed Brown   /* Test MatMult */
1249566063dSJacob Faibussowitsch   PetscCall(MatMultEqual(A, As, 10, &flg));
12528b400f6SJacob Faibussowitsch   PetscCheck(flg, PetscObjectComm((PetscObject)A), PETSC_ERR_PLIB, "MatMult not equal");
1269566063dSJacob Faibussowitsch   PetscCall(MatMultAddEqual(A, As, 10, &flg));
12728b400f6SJacob Faibussowitsch   PetscCheck(flg, PetscObjectComm((PetscObject)A), PETSC_ERR_PLIB, "MatMultAdd not equal");
128c4762a1bSJed Brown 
129c4762a1bSJed Brown   /* Free spaces */
1309566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
1319566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&As));
1329566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
133b122ec5aSJacob Faibussowitsch   return 0;
134c4762a1bSJed Brown }
135c4762a1bSJed Brown 
136c4762a1bSJed Brown /*TEST
137c4762a1bSJed Brown 
138c4762a1bSJed Brown    build:
139c4762a1bSJed Brown       requires: complex
140c4762a1bSJed Brown 
141c4762a1bSJed Brown    test:
142c4762a1bSJed Brown       args: -n 1000
143*3886731fSPierre Jolivet       output_file: output/empty.out
144c4762a1bSJed Brown 
145c4762a1bSJed Brown    test:
146c4762a1bSJed Brown       suffix: 2
147c4762a1bSJed Brown       nsize: 3
148c4762a1bSJed Brown       args: -n 1000
149*3886731fSPierre Jolivet       output_file: output/empty.out
150c4762a1bSJed Brown 
151c4762a1bSJed Brown    test:
152c4762a1bSJed Brown       suffix: superlu_dist
153c4762a1bSJed Brown       nsize: 3
154c4762a1bSJed Brown       requires: superlu_dist
155c4762a1bSJed Brown       args: -test_choleskyfactor -mat_superlu_dist_rowperm NOROWPERM
156c4762a1bSJed Brown       output_file: output/ex127_superlu_dist.out
157c4762a1bSJed Brown TEST*/
158