xref: /petsc/src/mat/tests/ex127.c (revision 2c71b3e237ead271e4f3aa1505f92bf476e3413d)
1c4762a1bSJed Brown static char help[] = "Test MatMult() for Hermitian matrix.\n\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown #include <petscmat.h>
4c4762a1bSJed Brown 
5c4762a1bSJed Brown int main(int argc,char **args)
6c4762a1bSJed Brown {
7c4762a1bSJed Brown   Mat            A,As;
8c4762a1bSJed Brown   PetscBool      flg;
9c4762a1bSJed Brown   PetscErrorCode ierr;
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 
16c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
17ffc4695bSBarry Smith   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr);
18c4762a1bSJed Brown   ierr = PetscOptionsGetReal(NULL,NULL,"-sigma1",&sigma1,NULL);CHKERRQ(ierr);
19c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);CHKERRQ(ierr);
20c4762a1bSJed Brown   dim  = n*n;
21c4762a1bSJed Brown 
22c4762a1bSJed Brown   ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr);
23c4762a1bSJed Brown   ierr = MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,dim,dim);CHKERRQ(ierr);
24c4762a1bSJed Brown   ierr = MatSetType(A,MATAIJ);CHKERRQ(ierr);
25c4762a1bSJed Brown   ierr = MatSetFromOptions(A);CHKERRQ(ierr);
26c4762a1bSJed Brown   ierr = MatSetUp(A);CHKERRQ(ierr);
27c4762a1bSJed Brown 
28c4762a1bSJed Brown   sigma2 = 10.0*PETSC_i;
29c4762a1bSJed Brown   h2 = 1.0/((n+1)*(n+1));
30c4762a1bSJed Brown 
31c4762a1bSJed Brown   ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr);
32c4762a1bSJed Brown   for (Ii=rstart; Ii<rend; Ii++) {
33c4762a1bSJed Brown     v = -1.0; i = Ii/n; j = Ii - i*n;
34c4762a1bSJed Brown     if (i>0) {
35c4762a1bSJed Brown       J = Ii-n; ierr = MatSetValues(A,1,&Ii,1,&J,&v,ADD_VALUES);CHKERRQ(ierr);
36c4762a1bSJed Brown     }
37c4762a1bSJed Brown     if (i<n-1) {
38c4762a1bSJed Brown       J = Ii+n; ierr = MatSetValues(A,1,&Ii,1,&J,&v,ADD_VALUES);CHKERRQ(ierr);
39c4762a1bSJed Brown     }
40c4762a1bSJed Brown     if (j>0) {
41c4762a1bSJed Brown       J = Ii-1; ierr = MatSetValues(A,1,&Ii,1,&J,&v,ADD_VALUES);CHKERRQ(ierr);
42c4762a1bSJed Brown     }
43c4762a1bSJed Brown     if (j<n-1) {
44c4762a1bSJed Brown       J = Ii+1; ierr = MatSetValues(A,1,&Ii,1,&J,&v,ADD_VALUES);CHKERRQ(ierr);
45c4762a1bSJed Brown     }
46c4762a1bSJed Brown     v    = 4.0 - sigma1*h2;
47c4762a1bSJed Brown     ierr = MatSetValues(A,1,&Ii,1,&Ii,&v,ADD_VALUES);CHKERRQ(ierr);
48c4762a1bSJed Brown   }
49c4762a1bSJed Brown   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
50c4762a1bSJed Brown   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
51c4762a1bSJed Brown 
52c4762a1bSJed Brown   /* Check whether A is symmetric */
53c4762a1bSJed Brown   ierr = PetscOptionsHasName(NULL,NULL, "-check_symmetric", &flg);CHKERRQ(ierr);
54c4762a1bSJed Brown   if (flg) {
55c4762a1bSJed Brown     ierr = MatIsSymmetric(A,0.0,&flg);CHKERRQ(ierr);
56*2c71b3e2SJacob Faibussowitsch     PetscCheckFalse(!flg,PETSC_COMM_SELF,PETSC_ERR_USER,"A is not symmetric");
57c4762a1bSJed Brown   }
58c4762a1bSJed Brown   ierr = MatSetOption(A,MAT_SYMMETRIC,PETSC_TRUE);CHKERRQ(ierr);
59c4762a1bSJed Brown 
60c4762a1bSJed Brown   /* make A complex Hermitian */
61c4762a1bSJed Brown   Ii = 0; J = dim-1;
62c4762a1bSJed Brown   if (Ii >= rstart && Ii < rend) {
63c4762a1bSJed Brown     v    = sigma2*h2; /* RealPart(v) = 0.0 */
64c4762a1bSJed Brown     ierr = MatSetValues(A,1,&Ii,1,&J,&v,ADD_VALUES);CHKERRQ(ierr);
65c4762a1bSJed Brown     v    = -sigma2*h2;
66c4762a1bSJed Brown     ierr = MatSetValues(A,1,&J,1,&Ii,&v,ADD_VALUES);CHKERRQ(ierr);
67c4762a1bSJed Brown   }
68c4762a1bSJed Brown 
69c4762a1bSJed Brown   Ii = dim-2; J = dim-1;
70c4762a1bSJed Brown   if (Ii >= rstart && Ii < rend) {
71c4762a1bSJed Brown     v    = sigma2*h2; /* RealPart(v) = 0.0 */
72c4762a1bSJed Brown     ierr = MatSetValues(A,1,&Ii,1,&J,&v,ADD_VALUES);CHKERRQ(ierr);
73c4762a1bSJed Brown     v    = -sigma2*h2;
74c4762a1bSJed Brown     ierr = MatSetValues(A,1,&J,1,&Ii,&v,ADD_VALUES);CHKERRQ(ierr);
75c4762a1bSJed Brown   }
76c4762a1bSJed Brown 
77c4762a1bSJed Brown   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
78c4762a1bSJed Brown   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
79c4762a1bSJed Brown   ierr = MatViewFromOptions(A,NULL,"-disp_mat");CHKERRQ(ierr);
80c4762a1bSJed Brown 
81c4762a1bSJed Brown   /* Check whether A is Hermitian, then set A->hermitian flag */
82c4762a1bSJed Brown   ierr = PetscOptionsHasName(NULL,NULL, "-check_Hermitian", &flg);CHKERRQ(ierr);
83c4762a1bSJed Brown   if (flg && size == 1) {
84c4762a1bSJed Brown     ierr = MatIsHermitian(A,0.0,&flg);CHKERRQ(ierr);
85*2c71b3e2SJacob Faibussowitsch     PetscCheckFalse(!flg,PETSC_COMM_SELF,PETSC_ERR_USER,"A is not Hermitian");
86c4762a1bSJed Brown   }
87c4762a1bSJed Brown   ierr = MatSetOption(A,MAT_HERMITIAN,PETSC_TRUE);CHKERRQ(ierr);
88c4762a1bSJed Brown 
89c4762a1bSJed Brown #if defined(PETSC_HAVE_SUPERLU_DIST)
90c4762a1bSJed Brown   /* Test Cholesky factorization */
91c4762a1bSJed Brown   ierr = PetscOptionsHasName(NULL,NULL, "-test_choleskyfactor", &flg);CHKERRQ(ierr);
92c4762a1bSJed Brown   if (flg) {
93c4762a1bSJed Brown     Mat      F;
94c4762a1bSJed Brown     IS       perm,iperm;
95c4762a1bSJed Brown     MatFactorInfo info;
96c4762a1bSJed Brown     PetscInt nneg,nzero,npos;
97c4762a1bSJed Brown 
98c4762a1bSJed Brown     ierr = MatGetFactor(A,MATSOLVERSUPERLU_DIST,MAT_FACTOR_CHOLESKY,&F);CHKERRQ(ierr);
99c4762a1bSJed Brown     ierr = MatGetOrdering(A,MATORDERINGND,&perm,&iperm);CHKERRQ(ierr);
100c4762a1bSJed Brown     ierr = MatCholeskyFactorSymbolic(F,A,perm,&info);CHKERRQ(ierr);
101c4762a1bSJed Brown     ierr = MatCholeskyFactorNumeric(F,A,&info);CHKERRQ(ierr);
102c4762a1bSJed Brown 
103c4762a1bSJed Brown     ierr = MatGetInertia(F,&nneg,&nzero,&npos);CHKERRQ(ierr);
1048fffc762SJacob Faibussowitsch     ierr = PetscPrintf(PETSC_COMM_WORLD," MatInertia: nneg: %" PetscInt_FMT ", nzero: %" PetscInt_FMT ", npos: %" PetscInt_FMT "\n",nneg,nzero,npos);CHKERRQ(ierr);
105c4762a1bSJed Brown     ierr = MatDestroy(&F);CHKERRQ(ierr);
106c4762a1bSJed Brown     ierr = ISDestroy(&perm);CHKERRQ(ierr);
107c4762a1bSJed Brown     ierr = ISDestroy(&iperm);CHKERRQ(ierr);
108c4762a1bSJed Brown   }
109c4762a1bSJed Brown #endif
110c4762a1bSJed Brown 
111c4762a1bSJed Brown   /* Create a Hermitian matrix As in sbaij format */
112c4762a1bSJed Brown   ierr = MatConvert(A,MATSBAIJ,MAT_INITIAL_MATRIX,&As);CHKERRQ(ierr);
113c4762a1bSJed Brown   ierr = MatViewFromOptions(As,NULL,"-disp_mat");CHKERRQ(ierr);
114c4762a1bSJed Brown 
115c4762a1bSJed Brown   /* Test MatMult */
116c4762a1bSJed Brown   ierr = MatMultEqual(A,As,10,&flg);CHKERRQ(ierr);
117*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(!flg,PetscObjectComm((PetscObject)A),PETSC_ERR_PLIB,"MatMult not equal");
118c4762a1bSJed Brown   ierr = MatMultAddEqual(A,As,10,&flg);CHKERRQ(ierr);
119*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(!flg,PetscObjectComm((PetscObject)A),PETSC_ERR_PLIB,"MatMultAdd not equal");
120c4762a1bSJed Brown 
121c4762a1bSJed Brown   /* Free spaces */
122c4762a1bSJed Brown   ierr = MatDestroy(&A);CHKERRQ(ierr);
123c4762a1bSJed Brown   ierr = MatDestroy(&As);CHKERRQ(ierr);
124c4762a1bSJed Brown   ierr = PetscFinalize();
125c4762a1bSJed Brown   return ierr;
126c4762a1bSJed Brown }
127c4762a1bSJed Brown 
128c4762a1bSJed Brown /*TEST
129c4762a1bSJed Brown 
130c4762a1bSJed Brown    build:
131c4762a1bSJed Brown       requires: complex
132c4762a1bSJed Brown 
133c4762a1bSJed Brown    test:
134c4762a1bSJed Brown       args: -n 1000
135c4762a1bSJed Brown       output_file: output/ex127.out
136c4762a1bSJed Brown 
137c4762a1bSJed Brown    test:
138c4762a1bSJed Brown       suffix: 2
139c4762a1bSJed Brown       nsize: 3
140c4762a1bSJed Brown       args: -n 1000
141c4762a1bSJed Brown       output_file: output/ex127.out
142c4762a1bSJed Brown 
143c4762a1bSJed Brown    test:
144c4762a1bSJed Brown       suffix: superlu_dist
145c4762a1bSJed Brown       nsize: 3
146c4762a1bSJed Brown       requires: superlu_dist
147c4762a1bSJed Brown       args: -test_choleskyfactor -mat_superlu_dist_rowperm NOROWPERM
148c4762a1bSJed Brown       output_file: output/ex127_superlu_dist.out
149c4762a1bSJed Brown TEST*/
150