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