xref: /petsc/src/mat/tests/ex215.c (revision 589a23caa660d2a5f330cc8d1ed213e9cfaf51a7)
1c4762a1bSJed Brown static char help[] = "Tests MatSolve(), MatSolveTranspose() and MatMatSolve() with SEQDENSE\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown #include <petscmat.h>
4c4762a1bSJed Brown 
5c4762a1bSJed Brown int main(int argc,char **args)
6c4762a1bSJed Brown {
7c4762a1bSJed Brown   Mat            A,RHS,C,F,X;
8c4762a1bSJed Brown   Vec            u,x,b;
9c4762a1bSJed Brown   PetscErrorCode ierr;
10c4762a1bSJed Brown   PetscMPIInt    size;
11c4762a1bSJed Brown   PetscInt       m,n,nsolve,nrhs;
12c4762a1bSJed Brown   PetscReal      norm,tol=PETSC_SQRT_MACHINE_EPSILON;
13c4762a1bSJed Brown   PetscRandom    rand;
14c4762a1bSJed Brown   PetscBool      data_provided,herm,symm,hpd;
15c4762a1bSJed Brown   MatFactorType  ftyp;
16c4762a1bSJed Brown   PetscViewer    fd;
17c4762a1bSJed Brown   char           file[PETSC_MAX_PATH_LEN];
18c4762a1bSJed Brown 
19c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
20c4762a1bSJed Brown   ierr = MPI_Comm_size(PETSC_COMM_WORLD, &size);CHKERRQ(ierr);
21c4762a1bSJed Brown   if (size > 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor test");
22c4762a1bSJed Brown   /* Determine which type of solver we want to test for */
23c4762a1bSJed Brown   herm = PETSC_FALSE;
24c4762a1bSJed Brown   symm = PETSC_FALSE;
25c4762a1bSJed Brown   hpd  = PETSC_FALSE;
26c4762a1bSJed Brown   ierr = PetscOptionsGetBool(NULL,NULL,"-symmetric_solve",&symm,NULL);CHKERRQ(ierr);
27c4762a1bSJed Brown   ierr = PetscOptionsGetBool(NULL,NULL,"-hermitian_solve",&herm,NULL);CHKERRQ(ierr);
28c4762a1bSJed Brown   ierr = PetscOptionsGetBool(NULL,NULL,"-hpd_solve",&hpd,NULL);CHKERRQ(ierr);
29c4762a1bSJed Brown 
30c4762a1bSJed Brown   /* Determine file from which we read the matrix A */
31c4762a1bSJed Brown   ftyp = MAT_FACTOR_LU;
32*589a23caSBarry Smith   ierr = PetscOptionsGetString(NULL,NULL,"-f",file,sizeof(file),&data_provided);CHKERRQ(ierr);
33c4762a1bSJed Brown   if (!data_provided) { /* get matrices from PETSc distribution */
34c4762a1bSJed Brown     ierr = PetscStrcpy(file,"${PETSC_DIR}/share/petsc/datafiles/matrices/");CHKERRQ(ierr);
35c4762a1bSJed Brown     if (hpd) {
36c4762a1bSJed Brown #if defined(PETSC_USE_COMPLEX)
37c4762a1bSJed Brown       ierr = PetscStrcat(file,"hpd-complex-");CHKERRQ(ierr);
38c4762a1bSJed Brown #else
39c4762a1bSJed Brown       ierr = PetscStrcat(file,"spd-real-");CHKERRQ(ierr);
40c4762a1bSJed Brown #endif
41c4762a1bSJed Brown       ftyp = MAT_FACTOR_CHOLESKY;
42c4762a1bSJed Brown     } else {
43c4762a1bSJed Brown #if defined(PETSC_USE_COMPLEX)
44c4762a1bSJed Brown       ierr = PetscStrcat(file,"nh-complex-");CHKERRQ(ierr);
45c4762a1bSJed Brown #else
46c4762a1bSJed Brown       ierr = PetscStrcat(file,"ns-real-");CHKERRQ(ierr);
47c4762a1bSJed Brown #endif
48c4762a1bSJed Brown     }
49c4762a1bSJed Brown #if defined(PETSC_USE_64BIT_INDICES)
50c4762a1bSJed Brown     ierr = PetscStrcat(file,"int64-");CHKERRQ(ierr);
51c4762a1bSJed Brown #else
52c4762a1bSJed Brown     ierr = PetscStrcat(file,"int32-");CHKERRQ(ierr);
53c4762a1bSJed Brown #endif
54c4762a1bSJed Brown #if defined(PETSC_USE_REAL_SINGLE)
55c4762a1bSJed Brown     ierr = PetscStrcat(file,"float32");CHKERRQ(ierr);
56c4762a1bSJed Brown #else
57c4762a1bSJed Brown     ierr = PetscStrcat(file,"float64");CHKERRQ(ierr);
58c4762a1bSJed Brown #endif
59c4762a1bSJed Brown   }
60c4762a1bSJed Brown 
61c4762a1bSJed Brown   /* Load matrix A */
62c4762a1bSJed Brown #if defined(PETSC_USE_REAL___FLOAT128)
63c4762a1bSJed Brown   ierr = PetscOptionsInsertString(NULL,"-binary_read_double");CHKERRQ(ierr);
64c4762a1bSJed Brown #endif
65c4762a1bSJed Brown   ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&fd);CHKERRQ(ierr);
66c4762a1bSJed Brown   ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr);
67c4762a1bSJed Brown   ierr = MatLoad(A,fd);CHKERRQ(ierr);
68c4762a1bSJed Brown   ierr = PetscViewerDestroy(&fd);CHKERRQ(ierr);
69c4762a1bSJed Brown   ierr = MatConvert(A,MATSEQDENSE,MAT_INPLACE_MATRIX,&A);CHKERRQ(ierr);
70c4762a1bSJed Brown   ierr = MatGetSize(A,&m,&n);CHKERRQ(ierr);
71c4762a1bSJed Brown   if (m != n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ, "This example is not intended for rectangular matrices (%d, %d)", m, n);
72c4762a1bSJed Brown 
73c4762a1bSJed Brown   /* Create dense matrix C and X; C holds true solution with identical colums */
74c4762a1bSJed Brown   nrhs = 2;
75c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL,"-nrhs",&nrhs,NULL);CHKERRQ(ierr);
76c4762a1bSJed Brown   ierr = MatCreate(PETSC_COMM_WORLD,&C);CHKERRQ(ierr);
77c4762a1bSJed Brown   ierr = MatSetSizes(C,m,PETSC_DECIDE,PETSC_DECIDE,nrhs);CHKERRQ(ierr);
78c4762a1bSJed Brown   ierr = MatSetType(C,MATDENSE);CHKERRQ(ierr);
79c4762a1bSJed Brown   ierr = MatSetFromOptions(C);CHKERRQ(ierr);
80c4762a1bSJed Brown   ierr = MatSetUp(C);CHKERRQ(ierr);
81c4762a1bSJed Brown 
82c4762a1bSJed Brown   ierr = PetscRandomCreate(PETSC_COMM_WORLD,&rand);CHKERRQ(ierr);
83c4762a1bSJed Brown   ierr = PetscRandomSetFromOptions(rand);CHKERRQ(ierr);
84c4762a1bSJed Brown   ierr = MatSetRandom(C,rand);CHKERRQ(ierr);
85c4762a1bSJed Brown   ierr = MatDuplicate(C,MAT_DO_NOT_COPY_VALUES,&X);CHKERRQ(ierr);
86c4762a1bSJed Brown   ierr = MatDuplicate(C,MAT_DO_NOT_COPY_VALUES,&RHS);CHKERRQ(ierr);
87c4762a1bSJed Brown 
88c4762a1bSJed Brown   /* Create vectors */
89c4762a1bSJed Brown   ierr = VecCreate(PETSC_COMM_WORLD,&x);CHKERRQ(ierr);
90c4762a1bSJed Brown   ierr = VecSetSizes(x,n,PETSC_DECIDE);CHKERRQ(ierr);
91c4762a1bSJed Brown   ierr = VecSetFromOptions(x);CHKERRQ(ierr);
92c4762a1bSJed Brown   ierr = VecDuplicate(x,&b);CHKERRQ(ierr);
93c4762a1bSJed Brown   ierr = VecDuplicate(x,&u);CHKERRQ(ierr); /* save the true solution */
94c4762a1bSJed Brown 
95c4762a1bSJed Brown   /* make a symmetric matrix */
96c4762a1bSJed Brown   if (symm) {
97c4762a1bSJed Brown     Mat AT;
98c4762a1bSJed Brown 
99c4762a1bSJed Brown     ierr = MatTranspose(A,MAT_INITIAL_MATRIX,&AT);CHKERRQ(ierr);
100c4762a1bSJed Brown     ierr = MatAXPY(A,1.0,AT,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
101c4762a1bSJed Brown     ierr = MatDestroy(&AT);CHKERRQ(ierr);
102c4762a1bSJed Brown     ftyp = MAT_FACTOR_CHOLESKY;
103c4762a1bSJed Brown   }
104c4762a1bSJed Brown   /* make an hermitian matrix */
105c4762a1bSJed Brown   if (herm) {
106c4762a1bSJed Brown     Mat AH;
107c4762a1bSJed Brown 
108c4762a1bSJed Brown     ierr = MatHermitianTranspose(A,MAT_INITIAL_MATRIX,&AH);CHKERRQ(ierr);
109c4762a1bSJed Brown     ierr = MatAXPY(A,1.0,AH,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
110c4762a1bSJed Brown     ierr = MatDestroy(&AH);CHKERRQ(ierr);
111c4762a1bSJed Brown     ftyp = MAT_FACTOR_CHOLESKY;
112c4762a1bSJed Brown   }
113c4762a1bSJed Brown   ierr = PetscObjectSetName((PetscObject)A,"A");CHKERRQ(ierr);
114c4762a1bSJed Brown   ierr = MatViewFromOptions(A,NULL,"-amat_view");CHKERRQ(ierr);
115c4762a1bSJed Brown 
116c4762a1bSJed Brown   ierr = MatDuplicate(A,MAT_COPY_VALUES,&F);CHKERRQ(ierr);
117c4762a1bSJed Brown   ierr = MatSetOption(F,MAT_SYMMETRIC,symm);CHKERRQ(ierr);
118c4762a1bSJed Brown   /* it seems that the SPD concept in PETSc extends naturally to Hermitian Positive definitess */
119c4762a1bSJed Brown   ierr = MatSetOption(F,MAT_HERMITIAN,(PetscBool)(hpd || herm));CHKERRQ(ierr);
120c4762a1bSJed Brown   ierr = MatSetOption(F,MAT_SPD,hpd);CHKERRQ(ierr);
121c4762a1bSJed Brown   if (ftyp == MAT_FACTOR_LU) {
122c4762a1bSJed Brown     ierr = MatLUFactor(F,NULL,NULL,NULL);CHKERRQ(ierr);
123c4762a1bSJed Brown   } else {
124c4762a1bSJed Brown     ierr = MatCholeskyFactor(F,NULL,NULL);CHKERRQ(ierr);
125c4762a1bSJed Brown   }
126c4762a1bSJed Brown 
127c4762a1bSJed Brown   for (nsolve = 0; nsolve < 2; nsolve++) {
128c4762a1bSJed Brown     ierr = VecSetRandom(x,rand);CHKERRQ(ierr);
129c4762a1bSJed Brown     ierr = VecCopy(x,u);CHKERRQ(ierr);
130c4762a1bSJed Brown     if (nsolve) {
131c4762a1bSJed Brown       ierr = MatMult(A,x,b);CHKERRQ(ierr);
132c4762a1bSJed Brown       ierr = MatSolve(F,b,x);CHKERRQ(ierr);
133c4762a1bSJed Brown     } else {
134c4762a1bSJed Brown       ierr = MatMultTranspose(A,x,b);CHKERRQ(ierr);
135c4762a1bSJed Brown       ierr = MatSolveTranspose(F,b,x);CHKERRQ(ierr);
136c4762a1bSJed Brown     }
137c4762a1bSJed Brown     /* Check the error */
138c4762a1bSJed Brown     ierr = VecAXPY(u,-1.0,x);CHKERRQ(ierr);  /* u <- (-1.0)x + u */
139c4762a1bSJed Brown     ierr = VecNorm(u,NORM_2,&norm);CHKERRQ(ierr);
140c4762a1bSJed Brown     if (norm > tol) {
141c4762a1bSJed Brown       PetscReal resi;
142c4762a1bSJed Brown       if (nsolve) {
143c4762a1bSJed Brown         ierr = MatMult(A,x,u);CHKERRQ(ierr); /* u = A*x */
144c4762a1bSJed Brown       } else {
145c4762a1bSJed Brown         ierr = MatMultTranspose(A,x,u);CHKERRQ(ierr); /* u = A*x */
146c4762a1bSJed Brown       }
147c4762a1bSJed Brown       ierr = VecAXPY(u,-1.0,b);CHKERRQ(ierr);  /* u <- (-1.0)b + u */
148c4762a1bSJed Brown       ierr = VecNorm(u,NORM_2,&resi);CHKERRQ(ierr);
149c4762a1bSJed Brown       if (nsolve) {
150c4762a1bSJed Brown         ierr = PetscPrintf(PETSC_COMM_SELF,"MatSolve error: Norm of error %g, residual %f\n",norm,resi);CHKERRQ(ierr);
151c4762a1bSJed Brown       } else {
152c4762a1bSJed Brown         ierr = PetscPrintf(PETSC_COMM_SELF,"MatSolveTranspose error: Norm of error %g, residual %f\n",norm,resi);CHKERRQ(ierr);
153c4762a1bSJed Brown       }
154c4762a1bSJed Brown     }
155c4762a1bSJed Brown   }
156c4762a1bSJed Brown   ierr = MatMatMult(A,C,MAT_REUSE_MATRIX,2.0,&RHS);CHKERRQ(ierr);
157c4762a1bSJed Brown   ierr = MatMatSolve(F,RHS,X);CHKERRQ(ierr);
158c4762a1bSJed Brown 
159c4762a1bSJed Brown   /* Check the error */
160c4762a1bSJed Brown   ierr = MatAXPY(X,-1.0,C,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
161c4762a1bSJed Brown   ierr = MatNorm(X,NORM_FROBENIUS,&norm);CHKERRQ(ierr);
162c4762a1bSJed Brown   if (norm > tol) {
163c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_SELF,"MatMatSolve: Norm of error %g\n",norm);CHKERRQ(ierr);
164c4762a1bSJed Brown   }
165c4762a1bSJed Brown 
166c4762a1bSJed Brown   /* Free data structures */
167c4762a1bSJed Brown   ierr = MatDestroy(&A);CHKERRQ(ierr);
168c4762a1bSJed Brown   ierr = MatDestroy(&C);CHKERRQ(ierr);
169c4762a1bSJed Brown   ierr = MatDestroy(&F);CHKERRQ(ierr);
170c4762a1bSJed Brown   ierr = MatDestroy(&X);CHKERRQ(ierr);
171c4762a1bSJed Brown   ierr = MatDestroy(&RHS);CHKERRQ(ierr);
172c4762a1bSJed Brown   ierr = PetscRandomDestroy(&rand);CHKERRQ(ierr);
173c4762a1bSJed Brown   ierr = VecDestroy(&x);CHKERRQ(ierr);
174c4762a1bSJed Brown   ierr = VecDestroy(&b);CHKERRQ(ierr);
175c4762a1bSJed Brown   ierr = VecDestroy(&u);CHKERRQ(ierr);
176c4762a1bSJed Brown   ierr = PetscFinalize();
177c4762a1bSJed Brown   return ierr;
178c4762a1bSJed Brown }
179c4762a1bSJed Brown 
180c4762a1bSJed Brown 
181c4762a1bSJed Brown /*TEST
182c4762a1bSJed Brown 
183c4762a1bSJed Brown   testset:
184c4762a1bSJed Brown     output_file: output/ex215.out
185c4762a1bSJed Brown     test:
186c4762a1bSJed Brown       suffix: ns
187c4762a1bSJed Brown     test:
188c4762a1bSJed Brown       suffix: sym
189c4762a1bSJed Brown       args: -symmetric_solve
190c4762a1bSJed Brown     test:
191c4762a1bSJed Brown       suffix: herm
192c4762a1bSJed Brown       args: -hermitian_solve
193c4762a1bSJed Brown     test:
194c4762a1bSJed Brown       suffix: hpd
195c4762a1bSJed Brown       args: -hpd_solve
196c4762a1bSJed Brown 
197c4762a1bSJed Brown TEST*/
198