xref: /petsc/src/mat/tests/ex215.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
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*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&args,(char*)0,help));
195f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
202c71b3e2SJacob Faibussowitsch   PetscCheckFalse(size > 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor test");
21c4762a1bSJed Brown   /* Determine which type of solver we want to test for */
22c4762a1bSJed Brown   herm = PETSC_FALSE;
23c4762a1bSJed Brown   symm = PETSC_FALSE;
24c4762a1bSJed Brown   hpd  = PETSC_FALSE;
255f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-symmetric_solve",&symm,NULL));
265f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-hermitian_solve",&herm,NULL));
275f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-hpd_solve",&hpd,NULL));
28c4762a1bSJed Brown 
29c4762a1bSJed Brown   /* Determine file from which we read the matrix A */
30c4762a1bSJed Brown   ftyp = MAT_FACTOR_LU;
315f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetString(NULL,NULL,"-f",file,sizeof(file),&data_provided));
32c4762a1bSJed Brown   if (!data_provided) { /* get matrices from PETSc distribution */
335f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscStrcpy(file,"${PETSC_DIR}/share/petsc/datafiles/matrices/"));
34c4762a1bSJed Brown     if (hpd) {
35c4762a1bSJed Brown #if defined(PETSC_USE_COMPLEX)
365f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscStrcat(file,"hpd-complex-"));
37c4762a1bSJed Brown #else
385f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscStrcat(file,"spd-real-"));
39c4762a1bSJed Brown #endif
40c4762a1bSJed Brown       ftyp = MAT_FACTOR_CHOLESKY;
41c4762a1bSJed Brown     } else {
42c4762a1bSJed Brown #if defined(PETSC_USE_COMPLEX)
435f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscStrcat(file,"nh-complex-"));
44c4762a1bSJed Brown #else
455f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscStrcat(file,"ns-real-"));
46c4762a1bSJed Brown #endif
47c4762a1bSJed Brown     }
48c4762a1bSJed Brown #if defined(PETSC_USE_64BIT_INDICES)
495f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscStrcat(file,"int64-"));
50c4762a1bSJed Brown #else
515f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscStrcat(file,"int32-"));
52c4762a1bSJed Brown #endif
53c4762a1bSJed Brown #if defined(PETSC_USE_REAL_SINGLE)
545f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscStrcat(file,"float32"));
55c4762a1bSJed Brown #else
565f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscStrcat(file,"float64"));
57c4762a1bSJed Brown #endif
58c4762a1bSJed Brown   }
59c4762a1bSJed Brown 
60c4762a1bSJed Brown   /* Load matrix A */
61c4762a1bSJed Brown #if defined(PETSC_USE_REAL___FLOAT128)
625f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsInsertString(NULL,"-binary_read_double"));
63c4762a1bSJed Brown #endif
645f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&fd));
655f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A));
665f80ce2aSJacob Faibussowitsch   CHKERRQ(MatLoad(A,fd));
675f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerDestroy(&fd));
685f80ce2aSJacob Faibussowitsch   CHKERRQ(MatConvert(A,MATSEQDENSE,MAT_INPLACE_MATRIX,&A));
695f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetSize(A,&m,&n));
702c71b3e2SJacob Faibussowitsch   PetscCheckFalse(m != n,PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ, "This example is not intended for rectangular matrices (%" PetscInt_FMT ", %" PetscInt_FMT ")", m, n);
71c4762a1bSJed Brown 
72a5b23f4aSJose E. Roman   /* Create dense matrix C and X; C holds true solution with identical columns */
73c4762a1bSJed Brown   nrhs = 2;
745f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-nrhs",&nrhs,NULL));
755f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&C));
765f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(C,m,PETSC_DECIDE,PETSC_DECIDE,nrhs));
775f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetType(C,MATDENSE));
785f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetFromOptions(C));
795f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetUp(C));
80c4762a1bSJed Brown 
815f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomCreate(PETSC_COMM_WORLD,&rand));
825f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomSetFromOptions(rand));
835f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetRandom(C,rand));
845f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDuplicate(C,MAT_DO_NOT_COPY_VALUES,&X));
855f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDuplicate(C,MAT_DO_NOT_COPY_VALUES,&RHS));
86c4762a1bSJed Brown 
87c4762a1bSJed Brown   /* Create vectors */
885f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreate(PETSC_COMM_WORLD,&x));
895f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetSizes(x,n,PETSC_DECIDE));
905f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetFromOptions(x));
915f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(x,&b));
925f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(x,&u)); /* save the true solution */
93c4762a1bSJed Brown 
94c4762a1bSJed Brown   /* make a symmetric matrix */
95c4762a1bSJed Brown   if (symm) {
96c4762a1bSJed Brown     Mat AT;
97c4762a1bSJed Brown 
985f80ce2aSJacob Faibussowitsch     CHKERRQ(MatTranspose(A,MAT_INITIAL_MATRIX,&AT));
995f80ce2aSJacob Faibussowitsch     CHKERRQ(MatAXPY(A,1.0,AT,SAME_NONZERO_PATTERN));
1005f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDestroy(&AT));
101c4762a1bSJed Brown     ftyp = MAT_FACTOR_CHOLESKY;
102c4762a1bSJed Brown   }
103c4762a1bSJed Brown   /* make an hermitian matrix */
104c4762a1bSJed Brown   if (herm) {
105c4762a1bSJed Brown     Mat AH;
106c4762a1bSJed Brown 
1075f80ce2aSJacob Faibussowitsch     CHKERRQ(MatHermitianTranspose(A,MAT_INITIAL_MATRIX,&AH));
1085f80ce2aSJacob Faibussowitsch     CHKERRQ(MatAXPY(A,1.0,AH,SAME_NONZERO_PATTERN));
1095f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDestroy(&AH));
110c4762a1bSJed Brown     ftyp = MAT_FACTOR_CHOLESKY;
111c4762a1bSJed Brown   }
1125f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectSetName((PetscObject)A,"A"));
1135f80ce2aSJacob Faibussowitsch   CHKERRQ(MatViewFromOptions(A,NULL,"-amat_view"));
114c4762a1bSJed Brown 
1155f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDuplicate(A,MAT_COPY_VALUES,&F));
1165f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetOption(F,MAT_SYMMETRIC,symm));
117c4762a1bSJed Brown   /* it seems that the SPD concept in PETSc extends naturally to Hermitian Positive definitess */
1185f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetOption(F,MAT_HERMITIAN,(PetscBool)(hpd || herm)));
1195f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetOption(F,MAT_SPD,hpd));
1204396437dSToby Isaac   {
1214396437dSToby Isaac     PetscInt iftyp = ftyp;
1225f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscOptionsGetEList(NULL,NULL,"-ftype",MatFactorTypes,MAT_FACTOR_NUM_TYPES,&iftyp,NULL));
1234396437dSToby Isaac     ftyp = (MatFactorType) iftyp;
1244396437dSToby Isaac   }
125c4762a1bSJed Brown   if (ftyp == MAT_FACTOR_LU) {
1265f80ce2aSJacob Faibussowitsch     CHKERRQ(MatLUFactor(F,NULL,NULL,NULL));
1274396437dSToby Isaac   } else if (ftyp == MAT_FACTOR_CHOLESKY) {
1285f80ce2aSJacob Faibussowitsch     CHKERRQ(MatCholeskyFactor(F,NULL,NULL));
1294396437dSToby Isaac   } else if (ftyp == MAT_FACTOR_QR) {
1305f80ce2aSJacob Faibussowitsch     CHKERRQ(MatQRFactor(F,NULL,NULL));
13198921bdaSJacob Faibussowitsch   } else SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "Factorization %s not supported in this example", MatFactorTypes[ftyp]);
132c4762a1bSJed Brown 
133c4762a1bSJed Brown   for (nsolve = 0; nsolve < 2; nsolve++) {
1345f80ce2aSJacob Faibussowitsch     CHKERRQ(VecSetRandom(x,rand));
1355f80ce2aSJacob Faibussowitsch     CHKERRQ(VecCopy(x,u));
136c4762a1bSJed Brown     if (nsolve) {
1375f80ce2aSJacob Faibussowitsch       CHKERRQ(MatMult(A,x,b));
1385f80ce2aSJacob Faibussowitsch       CHKERRQ(MatSolve(F,b,x));
139c4762a1bSJed Brown     } else {
1405f80ce2aSJacob Faibussowitsch       CHKERRQ(MatMultTranspose(A,x,b));
1415f80ce2aSJacob Faibussowitsch       CHKERRQ(MatSolveTranspose(F,b,x));
142c4762a1bSJed Brown     }
143c4762a1bSJed Brown     /* Check the error */
1445f80ce2aSJacob Faibussowitsch     CHKERRQ(VecAXPY(u,-1.0,x));  /* u <- (-1.0)x + u */
1455f80ce2aSJacob Faibussowitsch     CHKERRQ(VecNorm(u,NORM_2,&norm));
146c4762a1bSJed Brown     if (norm > tol) {
147c4762a1bSJed Brown       PetscReal resi;
148c4762a1bSJed Brown       if (nsolve) {
1495f80ce2aSJacob Faibussowitsch         CHKERRQ(MatMult(A,x,u)); /* u = A*x */
150c4762a1bSJed Brown       } else {
1515f80ce2aSJacob Faibussowitsch         CHKERRQ(MatMultTranspose(A,x,u)); /* u = A*x */
152c4762a1bSJed Brown       }
1535f80ce2aSJacob Faibussowitsch       CHKERRQ(VecAXPY(u,-1.0,b));  /* u <- (-1.0)b + u */
1545f80ce2aSJacob Faibussowitsch       CHKERRQ(VecNorm(u,NORM_2,&resi));
155c4762a1bSJed Brown       if (nsolve) {
1565f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"MatSolve error: Norm of error %g, residual %g\n",(double)norm,(double)resi));
157c4762a1bSJed Brown       } else {
1585f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"MatSolveTranspose error: Norm of error %g, residual %g\n",(double)norm,(double)resi));
159c4762a1bSJed Brown       }
160c4762a1bSJed Brown     }
161c4762a1bSJed Brown   }
1625f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMatMult(A,C,MAT_REUSE_MATRIX,2.0,&RHS));
1635f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMatSolve(F,RHS,X));
164c4762a1bSJed Brown 
165c4762a1bSJed Brown   /* Check the error */
1665f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAXPY(X,-1.0,C,SAME_NONZERO_PATTERN));
1675f80ce2aSJacob Faibussowitsch   CHKERRQ(MatNorm(X,NORM_FROBENIUS,&norm));
168c4762a1bSJed Brown   if (norm > tol) {
1695f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"MatMatSolve: Norm of error %g\n",(double)norm));
170c4762a1bSJed Brown   }
171c4762a1bSJed Brown 
172c4762a1bSJed Brown   /* Free data structures */
1735f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&A));
1745f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&C));
1755f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&F));
1765f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&X));
1775f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&RHS));
1785f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomDestroy(&rand));
1795f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&x));
1805f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&b));
1815f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&u));
182*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
183*b122ec5aSJacob Faibussowitsch   return 0;
184c4762a1bSJed Brown }
185c4762a1bSJed Brown 
186c4762a1bSJed Brown /*TEST
187c4762a1bSJed Brown 
188c4762a1bSJed Brown   testset:
189c4762a1bSJed Brown     output_file: output/ex215.out
190c4762a1bSJed Brown     test:
191c4762a1bSJed Brown       suffix: ns
192c4762a1bSJed Brown     test:
193c4762a1bSJed Brown       suffix: sym
194c4762a1bSJed Brown       args: -symmetric_solve
195c4762a1bSJed Brown     test:
196c4762a1bSJed Brown       suffix: herm
197c4762a1bSJed Brown       args: -hermitian_solve
198c4762a1bSJed Brown     test:
199c4762a1bSJed Brown       suffix: hpd
200c4762a1bSJed Brown       args: -hpd_solve
2014396437dSToby Isaac     test:
2024396437dSToby Isaac       suffix: qr
2034396437dSToby Isaac       args: -ftype qr
204c4762a1bSJed Brown 
205c4762a1bSJed Brown TEST*/
206