xref: /petsc/src/mat/tests/ex192.c (revision 5f80ce2ab25dff0f4601e710601cbbcecf323266)
1c4762a1bSJed Brown static char help[] = "Tests MatSolve() and MatMatSolve() with MUMPS or MKL_PARDISO sequential solvers in Schur complement mode.\n\
2c4762a1bSJed Brown Example: mpiexec -n 1 ./ex192 -f <matrix binary file> -nrhs 4 -symmetric_solve -hermitian_solve -schur_ratio 0.3\n\n";
3c4762a1bSJed Brown 
4c4762a1bSJed Brown #include <petscmat.h>
5c4762a1bSJed Brown 
6c4762a1bSJed Brown int main(int argc,char **args)
7c4762a1bSJed Brown {
8c4762a1bSJed Brown   Mat            A,RHS,C,F,X,S;
9c4762a1bSJed Brown   Vec            u,x,b;
10c4762a1bSJed Brown   Vec            xschur,bschur,uschur;
11c4762a1bSJed Brown   IS             is_schur;
12c4762a1bSJed Brown   PetscErrorCode ierr;
13c4762a1bSJed Brown   PetscMPIInt    size;
14c4762a1bSJed Brown   PetscInt       isolver=0,size_schur,m,n,nfact,nsolve,nrhs;
15c4762a1bSJed Brown   PetscReal      norm,tol=PETSC_SQRT_MACHINE_EPSILON;
16c4762a1bSJed Brown   PetscRandom    rand;
17c4762a1bSJed Brown   PetscBool      data_provided,herm,symm,use_lu,cuda = PETSC_FALSE;
18c4762a1bSJed Brown   PetscReal      sratio = 5.1/12.;
19c4762a1bSJed Brown   PetscViewer    fd;              /* viewer */
20c4762a1bSJed Brown   char           solver[256];
21c4762a1bSJed Brown   char           file[PETSC_MAX_PATH_LEN]; /* input file name */
22c4762a1bSJed Brown 
23c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
24*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
252c71b3e2SJacob Faibussowitsch   PetscCheckFalse(size > 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor test");
26c4762a1bSJed Brown   /* Determine which type of solver we want to test for */
27c4762a1bSJed Brown   herm = PETSC_FALSE;
28c4762a1bSJed Brown   symm = PETSC_FALSE;
29*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-symmetric_solve",&symm,NULL));
30*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-hermitian_solve",&herm,NULL));
31c4762a1bSJed Brown   if (herm) symm = PETSC_TRUE;
32*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-cuda_solve",&cuda,NULL));
33*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-tol",&tol,NULL));
34c4762a1bSJed Brown 
35c4762a1bSJed Brown   /* Determine file from which we read the matrix A */
36*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetString(NULL,NULL,"-f",file,sizeof(file),&data_provided));
37c4762a1bSJed Brown   if (!data_provided) { /* get matrices from PETSc distribution */
38*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscStrncpy(file,"${PETSC_DIR}/share/petsc/datafiles/matrices/",sizeof(file)));
39c4762a1bSJed Brown     if (symm) {
40c4762a1bSJed Brown #if defined (PETSC_USE_COMPLEX)
41*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscStrlcat(file,"hpd-complex-",sizeof(file)));
42c4762a1bSJed Brown #else
43*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscStrlcat(file,"spd-real-",sizeof(file)));
44c4762a1bSJed Brown #endif
45c4762a1bSJed Brown     } else {
46c4762a1bSJed Brown #if defined (PETSC_USE_COMPLEX)
47*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscStrlcat(file,"nh-complex-",sizeof(file)));
48c4762a1bSJed Brown #else
49*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscStrlcat(file,"ns-real-",sizeof(file)));
50c4762a1bSJed Brown #endif
51c4762a1bSJed Brown     }
52c4762a1bSJed Brown #if defined(PETSC_USE_64BIT_INDICES)
53*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscStrlcat(file,"int64-",sizeof(file)));
54c4762a1bSJed Brown #else
55*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscStrlcat(file,"int32-",sizeof(file)));
56c4762a1bSJed Brown #endif
57c4762a1bSJed Brown #if defined (PETSC_USE_REAL_SINGLE)
58*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscStrlcat(file,"float32",sizeof(file)));
59c4762a1bSJed Brown #else
60*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscStrlcat(file,"float64",sizeof(file)));
61c4762a1bSJed Brown #endif
62c4762a1bSJed Brown   }
63c4762a1bSJed Brown   /* Load matrix A */
64*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&fd));
65*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A));
66*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatLoad(A,fd));
67*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerDestroy(&fd));
68*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetSize(A,&m,&n));
692c71b3e2SJacob 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);
70c4762a1bSJed Brown 
71a5b23f4aSJose E. Roman   /* Create dense matrix C and X; C holds true solution with identical columns */
72c4762a1bSJed Brown   nrhs = 2;
73*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-nrhs",&nrhs,NULL));
74*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&C));
75*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(C,m,PETSC_DECIDE,PETSC_DECIDE,nrhs));
76*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetType(C,MATDENSE));
77*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetFromOptions(C));
78*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetUp(C));
79c4762a1bSJed Brown 
80*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomCreate(PETSC_COMM_WORLD,&rand));
81*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomSetFromOptions(rand));
82*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetRandom(C,rand));
83*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDuplicate(C,MAT_DO_NOT_COPY_VALUES,&X));
84c4762a1bSJed Brown 
85c4762a1bSJed Brown   /* Create vectors */
86*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreate(PETSC_COMM_WORLD,&x));
87*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetSizes(x,n,PETSC_DECIDE));
88*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetFromOptions(x));
89*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(x,&b));
90*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(x,&u)); /* save the true solution */
91c4762a1bSJed Brown 
92*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-solver",&isolver,NULL));
93c4762a1bSJed Brown   switch (isolver) {
94c4762a1bSJed Brown #if defined(PETSC_HAVE_MUMPS)
95c4762a1bSJed Brown     case 0:
96*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscStrcpy(solver,MATSOLVERMUMPS));
97c4762a1bSJed Brown       break;
98c4762a1bSJed Brown #endif
99c4762a1bSJed Brown #if defined(PETSC_HAVE_MKL_PARDISO)
100c4762a1bSJed Brown     case 1:
101*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscStrcpy(solver,MATSOLVERMKL_PARDISO));
102c4762a1bSJed Brown       break;
103c4762a1bSJed Brown #endif
104c4762a1bSJed Brown     default:
105*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscStrcpy(solver,MATSOLVERPETSC));
106c4762a1bSJed Brown       break;
107c4762a1bSJed Brown   }
108c4762a1bSJed Brown 
109c4762a1bSJed Brown #if defined (PETSC_USE_COMPLEX)
110c4762a1bSJed Brown   if (isolver == 0 && symm && !data_provided) { /* MUMPS (5.0.0) does not have support for hermitian matrices, so make them symmetric */
111c4762a1bSJed Brown     PetscScalar im = PetscSqrtScalar((PetscScalar)-1.);
112c4762a1bSJed Brown     PetscScalar val = -1.0;
113c4762a1bSJed Brown     val = val + im;
114*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetValue(A,1,0,val,INSERT_VALUES));
115*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
116*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
117c4762a1bSJed Brown   }
118c4762a1bSJed Brown #endif
119c4762a1bSJed Brown 
120*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-schur_ratio",&sratio,NULL));
1212c71b3e2SJacob Faibussowitsch   PetscCheckFalse(sratio < 0. || sratio > 1.,PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ, "Invalid ratio for schur degrees of freedom %g", (double)sratio);
122c4762a1bSJed Brown   size_schur = (PetscInt)(sratio*m);
123c4762a1bSJed Brown 
124*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"Solving with %s: nrhs %" PetscInt_FMT ", sym %d, herm %d, size schur %" PetscInt_FMT ", size mat %" PetscInt_FMT "\n",solver,nrhs,symm,herm,size_schur,m));
125c4762a1bSJed Brown 
126c4762a1bSJed Brown   /* Test LU/Cholesky Factorization */
127c4762a1bSJed Brown   use_lu = PETSC_FALSE;
128c4762a1bSJed Brown   if (!symm) use_lu = PETSC_TRUE;
129c4762a1bSJed Brown #if defined (PETSC_USE_COMPLEX)
130c4762a1bSJed Brown   if (isolver == 1) use_lu = PETSC_TRUE;
131c4762a1bSJed Brown #endif
132c4762a1bSJed Brown   if (cuda && symm && !herm) use_lu = PETSC_TRUE;
133c4762a1bSJed Brown 
134c4762a1bSJed Brown   if (herm && !use_lu) { /* test also conversion routines inside the solver packages */
135*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetOption(A,MAT_SYMMETRIC,PETSC_TRUE));
136*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatConvert(A,MATSEQSBAIJ,MAT_INPLACE_MATRIX,&A));
137c4762a1bSJed Brown   }
138c4762a1bSJed Brown 
139c4762a1bSJed Brown   if (use_lu) {
140*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatGetFactor(A,solver,MAT_FACTOR_LU,&F));
141c4762a1bSJed Brown   } else {
142c4762a1bSJed Brown     if (herm) {
143*5f80ce2aSJacob Faibussowitsch       CHKERRQ(MatSetOption(A,MAT_SYMMETRIC,PETSC_TRUE));
144*5f80ce2aSJacob Faibussowitsch       CHKERRQ(MatSetOption(A,MAT_SPD,PETSC_TRUE));
145c4762a1bSJed Brown     } else {
146*5f80ce2aSJacob Faibussowitsch       CHKERRQ(MatSetOption(A,MAT_SYMMETRIC,PETSC_TRUE));
147*5f80ce2aSJacob Faibussowitsch       CHKERRQ(MatSetOption(A,MAT_SPD,PETSC_FALSE));
148c4762a1bSJed Brown     }
149*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatGetFactor(A,solver,MAT_FACTOR_CHOLESKY,&F));
150c4762a1bSJed Brown   }
151*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISCreateStride(PETSC_COMM_SELF,size_schur,m-size_schur,1,&is_schur));
152*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatFactorSetSchurIS(F,is_schur));
153c4762a1bSJed Brown 
154*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&is_schur));
155c4762a1bSJed Brown   if (use_lu) {
156*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatLUFactorSymbolic(F,A,NULL,NULL,NULL));
157c4762a1bSJed Brown   } else {
158*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatCholeskyFactorSymbolic(F,A,NULL,NULL));
159c4762a1bSJed Brown   }
160c4762a1bSJed Brown 
161c4762a1bSJed Brown   for (nfact = 0; nfact < 3; nfact++) {
162c4762a1bSJed Brown     Mat AD;
163c4762a1bSJed Brown 
164c4762a1bSJed Brown     if (!nfact) {
165*5f80ce2aSJacob Faibussowitsch       CHKERRQ(VecSetRandom(x,rand));
166c4762a1bSJed Brown       if (symm && herm) {
167*5f80ce2aSJacob Faibussowitsch         CHKERRQ(VecAbs(x));
168c4762a1bSJed Brown       }
169*5f80ce2aSJacob Faibussowitsch       CHKERRQ(MatDiagonalSet(A,x,ADD_VALUES));
170c4762a1bSJed Brown     }
171c4762a1bSJed Brown     if (use_lu) {
172*5f80ce2aSJacob Faibussowitsch       CHKERRQ(MatLUFactorNumeric(F,A,NULL));
173c4762a1bSJed Brown     } else {
174*5f80ce2aSJacob Faibussowitsch       CHKERRQ(MatCholeskyFactorNumeric(F,A,NULL));
175c4762a1bSJed Brown     }
176c4762a1bSJed Brown     if (cuda) {
177*5f80ce2aSJacob Faibussowitsch       CHKERRQ(MatFactorGetSchurComplement(F,&S,NULL));
178*5f80ce2aSJacob Faibussowitsch       CHKERRQ(MatSetType(S,MATSEQDENSECUDA));
179*5f80ce2aSJacob Faibussowitsch       CHKERRQ(MatCreateVecs(S,&xschur,&bschur));
180*5f80ce2aSJacob Faibussowitsch       CHKERRQ(MatFactorRestoreSchurComplement(F,&S,MAT_FACTOR_SCHUR_UNFACTORED));
181c4762a1bSJed Brown     }
182*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatFactorCreateSchurComplement(F,&S,NULL));
183c4762a1bSJed Brown     if (!cuda) {
184*5f80ce2aSJacob Faibussowitsch       CHKERRQ(MatCreateVecs(S,&xschur,&bschur));
185c4762a1bSJed Brown     }
186*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDuplicate(xschur,&uschur));
187c4762a1bSJed Brown     if (nfact == 1 && (!cuda || (herm && symm))) {
188*5f80ce2aSJacob Faibussowitsch       CHKERRQ(MatFactorInvertSchurComplement(F));
189c4762a1bSJed Brown     }
190c4762a1bSJed Brown     for (nsolve = 0; nsolve < 2; nsolve++) {
191*5f80ce2aSJacob Faibussowitsch       CHKERRQ(VecSetRandom(x,rand));
192*5f80ce2aSJacob Faibussowitsch       CHKERRQ(VecCopy(x,u));
193c4762a1bSJed Brown 
194c4762a1bSJed Brown       if (nsolve) {
195*5f80ce2aSJacob Faibussowitsch         CHKERRQ(MatMult(A,x,b));
196*5f80ce2aSJacob Faibussowitsch         CHKERRQ(MatSolve(F,b,x));
197c4762a1bSJed Brown       } else {
198*5f80ce2aSJacob Faibussowitsch         CHKERRQ(MatMultTranspose(A,x,b));
199*5f80ce2aSJacob Faibussowitsch         CHKERRQ(MatSolveTranspose(F,b,x));
200c4762a1bSJed Brown       }
201c4762a1bSJed Brown       /* Check the error */
202*5f80ce2aSJacob Faibussowitsch       CHKERRQ(VecAXPY(u,-1.0,x));  /* u <- (-1.0)x + u */
203*5f80ce2aSJacob Faibussowitsch       CHKERRQ(VecNorm(u,NORM_2,&norm));
204c4762a1bSJed Brown       if (norm > tol) {
205c4762a1bSJed Brown         PetscReal resi;
206c4762a1bSJed Brown         if (nsolve) {
207*5f80ce2aSJacob Faibussowitsch           CHKERRQ(MatMult(A,x,u)); /* u = A*x */
208c4762a1bSJed Brown         } else {
209*5f80ce2aSJacob Faibussowitsch           CHKERRQ(MatMultTranspose(A,x,u)); /* u = A*x */
210c4762a1bSJed Brown         }
211*5f80ce2aSJacob Faibussowitsch         CHKERRQ(VecAXPY(u,-1.0,b));  /* u <- (-1.0)b + u */
212*5f80ce2aSJacob Faibussowitsch         CHKERRQ(VecNorm(u,NORM_2,&resi));
213c4762a1bSJed Brown         if (nsolve) {
214*5f80ce2aSJacob Faibussowitsch           CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"(f %" PetscInt_FMT ", s %" PetscInt_FMT ") MatSolve error: Norm of error %g, residual %g\n",nfact,nsolve,(double)norm,(double)resi));
215c4762a1bSJed Brown         } else {
216*5f80ce2aSJacob Faibussowitsch           CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"(f %" PetscInt_FMT ", s %" PetscInt_FMT ") MatSolveTranspose error: Norm of error %g, residual %f\n",nfact,nsolve,(double)norm,(double)resi));
217c4762a1bSJed Brown         }
218c4762a1bSJed Brown       }
219*5f80ce2aSJacob Faibussowitsch       CHKERRQ(VecSetRandom(xschur,rand));
220*5f80ce2aSJacob Faibussowitsch       CHKERRQ(VecCopy(xschur,uschur));
221c4762a1bSJed Brown       if (nsolve) {
222*5f80ce2aSJacob Faibussowitsch         CHKERRQ(MatMult(S,xschur,bschur));
223*5f80ce2aSJacob Faibussowitsch         CHKERRQ(MatFactorSolveSchurComplement(F,bschur,xschur));
224c4762a1bSJed Brown       } else {
225*5f80ce2aSJacob Faibussowitsch         CHKERRQ(MatMultTranspose(S,xschur,bschur));
226*5f80ce2aSJacob Faibussowitsch         CHKERRQ(MatFactorSolveSchurComplementTranspose(F,bschur,xschur));
227c4762a1bSJed Brown       }
228c4762a1bSJed Brown       /* Check the error */
229*5f80ce2aSJacob Faibussowitsch       CHKERRQ(VecAXPY(uschur,-1.0,xschur));  /* u <- (-1.0)x + u */
230*5f80ce2aSJacob Faibussowitsch       CHKERRQ(VecNorm(uschur,NORM_2,&norm));
231c4762a1bSJed Brown       if (norm > tol) {
232c4762a1bSJed Brown         PetscReal resi;
233c4762a1bSJed Brown         if (nsolve) {
234*5f80ce2aSJacob Faibussowitsch           CHKERRQ(MatMult(S,xschur,uschur)); /* u = A*x */
235c4762a1bSJed Brown         } else {
236*5f80ce2aSJacob Faibussowitsch           CHKERRQ(MatMultTranspose(S,xschur,uschur)); /* u = A*x */
237c4762a1bSJed Brown         }
238*5f80ce2aSJacob Faibussowitsch         CHKERRQ(VecAXPY(uschur,-1.0,bschur));  /* u <- (-1.0)b + u */
239*5f80ce2aSJacob Faibussowitsch         CHKERRQ(VecNorm(uschur,NORM_2,&resi));
240c4762a1bSJed Brown         if (nsolve) {
241*5f80ce2aSJacob Faibussowitsch           CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"(f %" PetscInt_FMT ", s %" PetscInt_FMT ") MatFactorSolveSchurComplement error: Norm of error %g, residual %g\n",nfact,nsolve,(double)norm,(double)resi));
242c4762a1bSJed Brown         } else {
243*5f80ce2aSJacob Faibussowitsch           CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"(f %" PetscInt_FMT ", s %" PetscInt_FMT ") MatFactorSolveSchurComplementTranspose error: Norm of error %g, residual %f\n",nfact,nsolve,(double)norm,(double)resi));
244c4762a1bSJed Brown         }
245c4762a1bSJed Brown       }
246c4762a1bSJed Brown     }
247*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatConvert(A,MATSEQAIJ,MAT_INITIAL_MATRIX,&AD));
248c4762a1bSJed Brown     if (!nfact) {
249*5f80ce2aSJacob Faibussowitsch       CHKERRQ(MatMatMult(AD,C,MAT_INITIAL_MATRIX,2.0,&RHS));
250c4762a1bSJed Brown     } else {
251*5f80ce2aSJacob Faibussowitsch       CHKERRQ(MatMatMult(AD,C,MAT_REUSE_MATRIX,2.0,&RHS));
252c4762a1bSJed Brown     }
253*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDestroy(&AD));
254c4762a1bSJed Brown     for (nsolve = 0; nsolve < 2; nsolve++) {
255*5f80ce2aSJacob Faibussowitsch       CHKERRQ(MatMatSolve(F,RHS,X));
256c4762a1bSJed Brown 
257c4762a1bSJed Brown       /* Check the error */
258*5f80ce2aSJacob Faibussowitsch       CHKERRQ(MatAXPY(X,-1.0,C,SAME_NONZERO_PATTERN));
259*5f80ce2aSJacob Faibussowitsch       CHKERRQ(MatNorm(X,NORM_FROBENIUS,&norm));
260c4762a1bSJed Brown       if (norm > tol) {
261*5f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"(f %" PetscInt_FMT ", s %" PetscInt_FMT ") MatMatSolve: Norm of error %g\n",nfact,nsolve,(double)norm));
262c4762a1bSJed Brown       }
263c4762a1bSJed Brown     }
264c4762a1bSJed Brown     if (isolver == 0) {
265c4762a1bSJed Brown       Mat spRHS,spRHST,RHST;
266c4762a1bSJed Brown 
267*5f80ce2aSJacob Faibussowitsch       CHKERRQ(MatTranspose(RHS,MAT_INITIAL_MATRIX,&RHST));
268*5f80ce2aSJacob Faibussowitsch       CHKERRQ(MatConvert(RHST,MATSEQAIJ,MAT_INITIAL_MATRIX,&spRHST));
269*5f80ce2aSJacob Faibussowitsch       CHKERRQ(MatCreateTranspose(spRHST,&spRHS));
270c4762a1bSJed Brown       for (nsolve = 0; nsolve < 2; nsolve++) {
271*5f80ce2aSJacob Faibussowitsch         CHKERRQ(MatMatSolve(F,spRHS,X));
272c4762a1bSJed Brown 
273c4762a1bSJed Brown         /* Check the error */
274*5f80ce2aSJacob Faibussowitsch         CHKERRQ(MatAXPY(X,-1.0,C,SAME_NONZERO_PATTERN));
275*5f80ce2aSJacob Faibussowitsch         CHKERRQ(MatNorm(X,NORM_FROBENIUS,&norm));
276c4762a1bSJed Brown         if (norm > tol) {
277*5f80ce2aSJacob Faibussowitsch           CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"(f %" PetscInt_FMT ", s %" PetscInt_FMT ") sparse MatMatSolve: Norm of error %g\n",nfact,nsolve,(double)norm));
278c4762a1bSJed Brown         }
279c4762a1bSJed Brown       }
280*5f80ce2aSJacob Faibussowitsch       CHKERRQ(MatDestroy(&spRHST));
281*5f80ce2aSJacob Faibussowitsch       CHKERRQ(MatDestroy(&spRHS));
282*5f80ce2aSJacob Faibussowitsch       CHKERRQ(MatDestroy(&RHST));
283c4762a1bSJed Brown     }
284*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDestroy(&S));
285*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDestroy(&xschur));
286*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDestroy(&bschur));
287*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDestroy(&uschur));
288c4762a1bSJed Brown   }
289c4762a1bSJed Brown   /* Free data structures */
290*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&A));
291*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&C));
292*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&F));
293*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&X));
294*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&RHS));
295*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomDestroy(&rand));
296*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&x));
297*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&b));
298*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&u));
299c4762a1bSJed Brown   ierr = PetscFinalize();
300c4762a1bSJed Brown   return ierr;
301c4762a1bSJed Brown }
302c4762a1bSJed Brown 
303c4762a1bSJed Brown /*TEST
304c4762a1bSJed Brown 
305c4762a1bSJed Brown    testset:
30638cf7977SStefano Zampini      requires: mkl_pardiso double !complex
307c4762a1bSJed Brown      args: -solver 1
308c4762a1bSJed Brown 
309c4762a1bSJed Brown      test:
310c4762a1bSJed Brown        suffix: mkl_pardiso
311c4762a1bSJed Brown      test:
312c4762a1bSJed Brown        requires: cuda
313c4762a1bSJed Brown        suffix: mkl_pardiso_cuda
314c4762a1bSJed Brown        args: -cuda_solve
315c4762a1bSJed Brown        output_file: output/ex192_mkl_pardiso.out
316c4762a1bSJed Brown      test:
317c4762a1bSJed Brown        suffix: mkl_pardiso_1
318c4762a1bSJed Brown        args: -symmetric_solve
319c4762a1bSJed Brown        output_file: output/ex192_mkl_pardiso_1.out
320c4762a1bSJed Brown      test:
321c4762a1bSJed Brown        requires: cuda
322c4762a1bSJed Brown        suffix: mkl_pardiso_cuda_1
323c4762a1bSJed Brown        args: -symmetric_solve -cuda_solve
324c4762a1bSJed Brown        output_file: output/ex192_mkl_pardiso_1.out
325c4762a1bSJed Brown      test:
326c4762a1bSJed Brown        suffix: mkl_pardiso_3
327c4762a1bSJed Brown        args: -symmetric_solve -hermitian_solve
328c4762a1bSJed Brown        output_file: output/ex192_mkl_pardiso_3.out
329c4762a1bSJed Brown      test:
330dfd57a17SPierre Jolivet        requires: cuda defined(PETSC_HAVE_CUSOLVERDNDPOTRI)
331c4762a1bSJed Brown        suffix: mkl_pardiso_cuda_3
332c4762a1bSJed Brown        args: -symmetric_solve -hermitian_solve -cuda_solve
333c4762a1bSJed Brown        output_file: output/ex192_mkl_pardiso_3.out
334c4762a1bSJed Brown 
335c4762a1bSJed Brown    testset:
336c4762a1bSJed Brown      requires: mumps double !complex
337c4762a1bSJed Brown      args: -solver 0
338c4762a1bSJed Brown 
339c4762a1bSJed Brown      test:
340c4762a1bSJed Brown        suffix: mumps
341c4762a1bSJed Brown      test:
342c4762a1bSJed Brown        requires: cuda
343c4762a1bSJed Brown        suffix: mumps_cuda
344c4762a1bSJed Brown        args: -cuda_solve
345c4762a1bSJed Brown        output_file: output/ex192_mumps.out
346c4762a1bSJed Brown      test:
347c4762a1bSJed Brown        suffix: mumps_2
348c4762a1bSJed Brown        args: -symmetric_solve
349c4762a1bSJed Brown        output_file: output/ex192_mumps_2.out
350c4762a1bSJed Brown      test:
351c4762a1bSJed Brown        requires: cuda
352c4762a1bSJed Brown        suffix: mumps_cuda_2
353c4762a1bSJed Brown        args: -symmetric_solve -cuda_solve
354c4762a1bSJed Brown        output_file: output/ex192_mumps_2.out
355c4762a1bSJed Brown      test:
356c4762a1bSJed Brown        suffix: mumps_3
357c4762a1bSJed Brown        args: -symmetric_solve -hermitian_solve
358c4762a1bSJed Brown        output_file: output/ex192_mumps_3.out
359c4762a1bSJed Brown      test:
360dfd57a17SPierre Jolivet        requires: cuda defined(PETSC_HAVE_CUSOLVERDNDPOTRI)
361c4762a1bSJed Brown        suffix: mumps_cuda_3
362c4762a1bSJed Brown        args: -symmetric_solve -hermitian_solve -cuda_solve
363c4762a1bSJed Brown        output_file: output/ex192_mumps_3.out
364c4762a1bSJed Brown 
365c4762a1bSJed Brown TEST*/
366