xref: /petsc/src/mat/tests/ex215.c (revision 327415f76d85372a4417cf1aaa14db707d4d6c04)
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   PetscMPIInt    size;
10c4762a1bSJed Brown   PetscInt       m,n,nsolve,nrhs;
11c4762a1bSJed Brown   PetscReal      norm,tol=PETSC_SQRT_MACHINE_EPSILON;
12c4762a1bSJed Brown   PetscRandom    rand;
13c4762a1bSJed Brown   PetscBool      data_provided,herm,symm,hpd;
14c4762a1bSJed Brown   MatFactorType  ftyp;
15c4762a1bSJed Brown   PetscViewer    fd;
16c4762a1bSJed Brown   char           file[PETSC_MAX_PATH_LEN];
17c4762a1bSJed Brown 
18*327415f7SBarry Smith   PetscFunctionBeginUser;
199566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&args,(char*)0,help));
209566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
21be096a46SBarry Smith   PetscCheck(size == 1,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;
269566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL,NULL,"-symmetric_solve",&symm,NULL));
279566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL,NULL,"-hermitian_solve",&herm,NULL));
289566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL,NULL,"-hpd_solve",&hpd,NULL));
29c4762a1bSJed Brown 
30c4762a1bSJed Brown   /* Determine file from which we read the matrix A */
31c4762a1bSJed Brown   ftyp = MAT_FACTOR_LU;
329566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetString(NULL,NULL,"-f",file,sizeof(file),&data_provided));
33c4762a1bSJed Brown   if (!data_provided) { /* get matrices from PETSc distribution */
349566063dSJacob Faibussowitsch     PetscCall(PetscStrcpy(file,"${PETSC_DIR}/share/petsc/datafiles/matrices/"));
35c4762a1bSJed Brown     if (hpd) {
36c4762a1bSJed Brown #if defined(PETSC_USE_COMPLEX)
379566063dSJacob Faibussowitsch       PetscCall(PetscStrcat(file,"hpd-complex-"));
38c4762a1bSJed Brown #else
399566063dSJacob Faibussowitsch       PetscCall(PetscStrcat(file,"spd-real-"));
40c4762a1bSJed Brown #endif
41c4762a1bSJed Brown       ftyp = MAT_FACTOR_CHOLESKY;
42c4762a1bSJed Brown     } else {
43c4762a1bSJed Brown #if defined(PETSC_USE_COMPLEX)
449566063dSJacob Faibussowitsch       PetscCall(PetscStrcat(file,"nh-complex-"));
45c4762a1bSJed Brown #else
469566063dSJacob Faibussowitsch       PetscCall(PetscStrcat(file,"ns-real-"));
47c4762a1bSJed Brown #endif
48c4762a1bSJed Brown     }
49c4762a1bSJed Brown #if defined(PETSC_USE_64BIT_INDICES)
509566063dSJacob Faibussowitsch     PetscCall(PetscStrcat(file,"int64-"));
51c4762a1bSJed Brown #else
529566063dSJacob Faibussowitsch     PetscCall(PetscStrcat(file,"int32-"));
53c4762a1bSJed Brown #endif
54c4762a1bSJed Brown #if defined(PETSC_USE_REAL_SINGLE)
559566063dSJacob Faibussowitsch     PetscCall(PetscStrcat(file,"float32"));
56c4762a1bSJed Brown #else
579566063dSJacob Faibussowitsch     PetscCall(PetscStrcat(file,"float64"));
58c4762a1bSJed Brown #endif
59c4762a1bSJed Brown   }
60c4762a1bSJed Brown 
61c4762a1bSJed Brown   /* Load matrix A */
62c4762a1bSJed Brown #if defined(PETSC_USE_REAL___FLOAT128)
639566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInsertString(NULL,"-binary_read_double"));
64c4762a1bSJed Brown #endif
659566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&fd));
669566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
679566063dSJacob Faibussowitsch   PetscCall(MatLoad(A,fd));
689566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&fd));
699566063dSJacob Faibussowitsch   PetscCall(MatConvert(A,MATSEQDENSE,MAT_INPLACE_MATRIX,&A));
709566063dSJacob Faibussowitsch   PetscCall(MatGetSize(A,&m,&n));
7108401ef6SPierre Jolivet   PetscCheck(m == n,PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ, "This example is not intended for rectangular matrices (%" PetscInt_FMT ", %" PetscInt_FMT ")", m, n);
72c4762a1bSJed Brown 
73a5b23f4aSJose E. Roman   /* Create dense matrix C and X; C holds true solution with identical columns */
74c4762a1bSJed Brown   nrhs = 2;
759566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-nrhs",&nrhs,NULL));
769566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD,&C));
779566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(C,m,PETSC_DECIDE,PETSC_DECIDE,nrhs));
789566063dSJacob Faibussowitsch   PetscCall(MatSetType(C,MATDENSE));
799566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(C));
809566063dSJacob Faibussowitsch   PetscCall(MatSetUp(C));
81c4762a1bSJed Brown 
829566063dSJacob Faibussowitsch   PetscCall(PetscRandomCreate(PETSC_COMM_WORLD,&rand));
839566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetFromOptions(rand));
849566063dSJacob Faibussowitsch   PetscCall(MatSetRandom(C,rand));
859566063dSJacob Faibussowitsch   PetscCall(MatDuplicate(C,MAT_DO_NOT_COPY_VALUES,&X));
869566063dSJacob Faibussowitsch   PetscCall(MatDuplicate(C,MAT_DO_NOT_COPY_VALUES,&RHS));
87c4762a1bSJed Brown 
88c4762a1bSJed Brown   /* Create vectors */
899566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD,&x));
909566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(x,n,PETSC_DECIDE));
919566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(x));
929566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x,&b));
939566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x,&u)); /* save the true solution */
94c4762a1bSJed Brown 
95c4762a1bSJed Brown   /* make a symmetric matrix */
96c4762a1bSJed Brown   if (symm) {
97c4762a1bSJed Brown     Mat AT;
98c4762a1bSJed Brown 
999566063dSJacob Faibussowitsch     PetscCall(MatTranspose(A,MAT_INITIAL_MATRIX,&AT));
1009566063dSJacob Faibussowitsch     PetscCall(MatAXPY(A,1.0,AT,SAME_NONZERO_PATTERN));
1019566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&AT));
102c4762a1bSJed Brown     ftyp = MAT_FACTOR_CHOLESKY;
103c4762a1bSJed Brown   }
104c4762a1bSJed Brown   /* make an hermitian matrix */
105c4762a1bSJed Brown   if (herm) {
106c4762a1bSJed Brown     Mat AH;
107c4762a1bSJed Brown 
1089566063dSJacob Faibussowitsch     PetscCall(MatHermitianTranspose(A,MAT_INITIAL_MATRIX,&AH));
1099566063dSJacob Faibussowitsch     PetscCall(MatAXPY(A,1.0,AH,SAME_NONZERO_PATTERN));
1109566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&AH));
111c4762a1bSJed Brown     ftyp = MAT_FACTOR_CHOLESKY;
112c4762a1bSJed Brown   }
1139566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)A,"A"));
1149566063dSJacob Faibussowitsch   PetscCall(MatViewFromOptions(A,NULL,"-amat_view"));
115c4762a1bSJed Brown 
1169566063dSJacob Faibussowitsch   PetscCall(MatDuplicate(A,MAT_COPY_VALUES,&F));
1179566063dSJacob Faibussowitsch   PetscCall(MatSetOption(F,MAT_SYMMETRIC,symm));
118c4762a1bSJed Brown   /* it seems that the SPD concept in PETSc extends naturally to Hermitian Positive definitess */
1199566063dSJacob Faibussowitsch   PetscCall(MatSetOption(F,MAT_HERMITIAN,(PetscBool)(hpd || herm)));
1209566063dSJacob Faibussowitsch   PetscCall(MatSetOption(F,MAT_SPD,hpd));
1214396437dSToby Isaac   {
1224396437dSToby Isaac     PetscInt iftyp = ftyp;
1239566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetEList(NULL,NULL,"-ftype",MatFactorTypes,MAT_FACTOR_NUM_TYPES,&iftyp,NULL));
1244396437dSToby Isaac     ftyp = (MatFactorType) iftyp;
1254396437dSToby Isaac   }
126c4762a1bSJed Brown   if (ftyp == MAT_FACTOR_LU) {
1279566063dSJacob Faibussowitsch     PetscCall(MatLUFactor(F,NULL,NULL,NULL));
1284396437dSToby Isaac   } else if (ftyp == MAT_FACTOR_CHOLESKY) {
1299566063dSJacob Faibussowitsch     PetscCall(MatCholeskyFactor(F,NULL,NULL));
1304396437dSToby Isaac   } else if (ftyp == MAT_FACTOR_QR) {
1319566063dSJacob Faibussowitsch     PetscCall(MatQRFactor(F,NULL,NULL));
13298921bdaSJacob Faibussowitsch   } else SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "Factorization %s not supported in this example", MatFactorTypes[ftyp]);
133c4762a1bSJed Brown 
134c4762a1bSJed Brown   for (nsolve = 0; nsolve < 2; nsolve++) {
1359566063dSJacob Faibussowitsch     PetscCall(VecSetRandom(x,rand));
1369566063dSJacob Faibussowitsch     PetscCall(VecCopy(x,u));
137c4762a1bSJed Brown     if (nsolve) {
1389566063dSJacob Faibussowitsch       PetscCall(MatMult(A,x,b));
1399566063dSJacob Faibussowitsch       PetscCall(MatSolve(F,b,x));
140c4762a1bSJed Brown     } else {
1419566063dSJacob Faibussowitsch       PetscCall(MatMultTranspose(A,x,b));
1429566063dSJacob Faibussowitsch       PetscCall(MatSolveTranspose(F,b,x));
143c4762a1bSJed Brown     }
144c4762a1bSJed Brown     /* Check the error */
1459566063dSJacob Faibussowitsch     PetscCall(VecAXPY(u,-1.0,x));  /* u <- (-1.0)x + u */
1469566063dSJacob Faibussowitsch     PetscCall(VecNorm(u,NORM_2,&norm));
147c4762a1bSJed Brown     if (norm > tol) {
148c4762a1bSJed Brown       PetscReal resi;
149c4762a1bSJed Brown       if (nsolve) {
1509566063dSJacob Faibussowitsch         PetscCall(MatMult(A,x,u)); /* u = A*x */
151c4762a1bSJed Brown       } else {
1529566063dSJacob Faibussowitsch         PetscCall(MatMultTranspose(A,x,u)); /* u = A*x */
153c4762a1bSJed Brown       }
1549566063dSJacob Faibussowitsch       PetscCall(VecAXPY(u,-1.0,b));  /* u <- (-1.0)b + u */
1559566063dSJacob Faibussowitsch       PetscCall(VecNorm(u,NORM_2,&resi));
156c4762a1bSJed Brown       if (nsolve) {
1579566063dSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_SELF,"MatSolve error: Norm of error %g, residual %g\n",(double)norm,(double)resi));
158c4762a1bSJed Brown       } else {
1599566063dSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_SELF,"MatSolveTranspose error: Norm of error %g, residual %g\n",(double)norm,(double)resi));
160c4762a1bSJed Brown       }
161c4762a1bSJed Brown     }
162c4762a1bSJed Brown   }
1639566063dSJacob Faibussowitsch   PetscCall(MatMatMult(A,C,MAT_REUSE_MATRIX,2.0,&RHS));
1649566063dSJacob Faibussowitsch   PetscCall(MatMatSolve(F,RHS,X));
165c4762a1bSJed Brown 
166c4762a1bSJed Brown   /* Check the error */
1679566063dSJacob Faibussowitsch   PetscCall(MatAXPY(X,-1.0,C,SAME_NONZERO_PATTERN));
1689566063dSJacob Faibussowitsch   PetscCall(MatNorm(X,NORM_FROBENIUS,&norm));
169c4762a1bSJed Brown   if (norm > tol) {
1709566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_SELF,"MatMatSolve: Norm of error %g\n",(double)norm));
171c4762a1bSJed Brown   }
172c4762a1bSJed Brown 
173c4762a1bSJed Brown   /* Free data structures */
1749566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
1759566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&C));
1769566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&F));
1779566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&X));
1789566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&RHS));
1799566063dSJacob Faibussowitsch   PetscCall(PetscRandomDestroy(&rand));
1809566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
1819566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&b));
1829566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&u));
1839566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
184b122ec5aSJacob Faibussowitsch   return 0;
185c4762a1bSJed Brown }
186c4762a1bSJed Brown 
187c4762a1bSJed Brown /*TEST
188c4762a1bSJed Brown 
189c4762a1bSJed Brown   testset:
190c4762a1bSJed Brown     output_file: output/ex215.out
191c4762a1bSJed Brown     test:
192c4762a1bSJed Brown       suffix: ns
193c4762a1bSJed Brown     test:
194c4762a1bSJed Brown       suffix: sym
195c4762a1bSJed Brown       args: -symmetric_solve
196c4762a1bSJed Brown     test:
197c4762a1bSJed Brown       suffix: herm
198c4762a1bSJed Brown       args: -hermitian_solve
199c4762a1bSJed Brown     test:
200c4762a1bSJed Brown       suffix: hpd
201c4762a1bSJed Brown       args: -hpd_solve
2024396437dSToby Isaac     test:
2034396437dSToby Isaac       suffix: qr
2044396437dSToby Isaac       args: -ftype qr
205c4762a1bSJed Brown 
206c4762a1bSJed Brown TEST*/
207