xref: /petsc/src/mat/tests/ex125.c (revision b18964edf8ca6736747b953a7963a9bb21bc6aad)
1c4762a1bSJed Brown static char help[] = "Tests MatSolve() and MatMatSolve() (interface to superlu_dist, mumps and mkl_pardiso).\n\
2c4762a1bSJed Brown Example: mpiexec -n <np> ./ex125 -f <matrix binary file> -nrhs 4 \n\n";
3c4762a1bSJed Brown 
4c4762a1bSJed Brown #include <petscmat.h>
5c4762a1bSJed Brown 
6c4762a1bSJed Brown int main(int argc,char **args)
7c4762a1bSJed Brown {
8*b18964edSHong Zhang   Mat            A,RHS=NULL,RHS1=NULL,C,F,X;
9c4762a1bSJed Brown   Vec            u,x,b;
10c4762a1bSJed Brown   PetscMPIInt    size;
11c4762a1bSJed Brown   PetscInt       m,n,nfact,nsolve,nrhs,ipack=0;
12c4762a1bSJed Brown   PetscReal      norm,tol=1.e-10;
13c4762a1bSJed Brown   IS             perm,iperm;
14c4762a1bSJed Brown   MatFactorInfo  info;
15c4762a1bSJed Brown   PetscRandom    rand;
16*b18964edSHong Zhang   PetscBool      flg,testMatSolve=PETSC_TRUE,testMatMatSolve=PETSC_TRUE,testMatMatSolveTranspose=PETSC_TRUE;
17c4762a1bSJed Brown   PetscBool      chol=PETSC_FALSE,view=PETSC_FALSE,matsolvexx = PETSC_FALSE;
18c4762a1bSJed Brown #if defined(PETSC_HAVE_MUMPS)
19c4762a1bSJed Brown   PetscBool      test_mumps_opts=PETSC_FALSE;
20c4762a1bSJed Brown #endif
21c4762a1bSJed Brown   PetscViewer    fd;              /* viewer */
22c4762a1bSJed Brown   char           file[PETSC_MAX_PATH_LEN]; /* input file name */
23c4762a1bSJed Brown 
24327415f7SBarry Smith   PetscFunctionBeginUser;
259566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&args,(char*)0,help));
269566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
27c4762a1bSJed Brown 
28c4762a1bSJed Brown   /* Determine file from which we read the matrix A */
299566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetString(NULL,NULL,"-f",file,sizeof(file),&flg));
309a14fc28SStefano Zampini   if (flg) { /* Load matrix A */
319566063dSJacob Faibussowitsch     PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&fd));
329566063dSJacob Faibussowitsch     PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
339566063dSJacob Faibussowitsch     PetscCall(MatSetFromOptions(A));
349566063dSJacob Faibussowitsch     PetscCall(MatLoad(A,fd));
359566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&fd));
369a14fc28SStefano Zampini   } else {
379a14fc28SStefano Zampini     n = 13;
389566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
399566063dSJacob Faibussowitsch     PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
409566063dSJacob Faibussowitsch     PetscCall(MatSetType(A,MATAIJ));
419566063dSJacob Faibussowitsch     PetscCall(MatSetFromOptions(A));
429566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n));
439566063dSJacob Faibussowitsch     PetscCall(MatSetUp(A));
449566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
459566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
469566063dSJacob Faibussowitsch     PetscCall(MatShift(A,1.0));
479a14fc28SStefano Zampini   }
489566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(A,&m,&n));
49e00437b9SBarry Smith   PetscCheck(m == n,PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ, "This example is not intended for rectangular matrices (%" PetscInt_FMT ", %" PetscInt_FMT ")", m, n);
50c4762a1bSJed Brown 
51c4762a1bSJed Brown   /* if A is symmetric, set its flag -- required by MatGetInertia() */
529566063dSJacob Faibussowitsch   PetscCall(MatIsSymmetric(A,0.0,&flg));
53c4762a1bSJed Brown 
549566063dSJacob Faibussowitsch   PetscCall(MatViewFromOptions(A,NULL,"-A_view"));
55c4762a1bSJed Brown 
56a5b23f4aSJose E. Roman   /* Create dense matrix C and X; C holds true solution with identical columns */
57c4762a1bSJed Brown   nrhs = 2;
589566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-nrhs",&nrhs,NULL));
599566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"ex125: nrhs %" PetscInt_FMT "\n",nrhs));
609566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD,&C));
619566063dSJacob Faibussowitsch   PetscCall(MatSetOptionsPrefix(C,"rhs_"));
629566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(C,m,PETSC_DECIDE,PETSC_DECIDE,nrhs));
639566063dSJacob Faibussowitsch   PetscCall(MatSetType(C,MATDENSE));
649566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(C));
659566063dSJacob Faibussowitsch   PetscCall(MatSetUp(C));
66c4762a1bSJed Brown 
679566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL,NULL,"-view_factor",&view,NULL));
689566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL,NULL,"-test_matmatsolve",&testMatMatSolve,NULL));
699566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL,NULL,"-cholesky",&chol,NULL));
70c4762a1bSJed Brown #if defined(PETSC_HAVE_MUMPS)
719566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL,NULL,"-test_mumps_opts",&test_mumps_opts,NULL));
72c4762a1bSJed Brown #endif
73c4762a1bSJed Brown 
749566063dSJacob Faibussowitsch   PetscCall(PetscRandomCreate(PETSC_COMM_WORLD,&rand));
759566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetFromOptions(rand));
769566063dSJacob Faibussowitsch   PetscCall(MatSetRandom(C,rand));
779566063dSJacob Faibussowitsch   PetscCall(MatDuplicate(C,MAT_DO_NOT_COPY_VALUES,&X));
78c4762a1bSJed Brown 
79c4762a1bSJed Brown   /* Create vectors */
809566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(A,&x,&b));
819566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x,&u)); /* save the true solution */
82c4762a1bSJed Brown 
83c4762a1bSJed Brown   /* Test Factorization */
849566063dSJacob Faibussowitsch   PetscCall(MatGetOrdering(A,MATORDERINGND,&perm,&iperm));
85c4762a1bSJed Brown 
869566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-mat_solver_type",&ipack,NULL));
87c4762a1bSJed Brown   switch (ipack) {
88c4762a1bSJed Brown #if defined(PETSC_HAVE_SUPERLU)
89c4762a1bSJed Brown   case 0:
9028b400f6SJacob Faibussowitsch     PetscCheck(!chol,PETSC_COMM_WORLD,PETSC_ERR_SUP,"SuperLU does not provide Cholesky!");
919566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD," SUPERLU LU:\n"));
929566063dSJacob Faibussowitsch     PetscCall(MatGetFactor(A,MATSOLVERSUPERLU,MAT_FACTOR_LU,&F));
93c4762a1bSJed Brown     matsolvexx = PETSC_TRUE;
94c4762a1bSJed Brown     break;
95c4762a1bSJed Brown #endif
96c4762a1bSJed Brown #if defined(PETSC_HAVE_SUPERLU_DIST)
97c4762a1bSJed Brown   case 1:
9828b400f6SJacob Faibussowitsch     PetscCheck(!chol,PETSC_COMM_WORLD,PETSC_ERR_SUP,"SuperLU does not provide Cholesky!");
999566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD," SUPERLU_DIST LU:\n"));
1009566063dSJacob Faibussowitsch     PetscCall(MatGetFactor(A,MATSOLVERSUPERLU_DIST,MAT_FACTOR_LU,&F));
101c4762a1bSJed Brown     matsolvexx = PETSC_TRUE;
102c4762a1bSJed Brown     break;
103c4762a1bSJed Brown #endif
104c4762a1bSJed Brown #if defined(PETSC_HAVE_MUMPS)
105c4762a1bSJed Brown   case 2:
106c4762a1bSJed Brown     if (chol) {
1079566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD," MUMPS CHOLESKY:\n"));
1089566063dSJacob Faibussowitsch       PetscCall(MatGetFactor(A,MATSOLVERMUMPS,MAT_FACTOR_CHOLESKY,&F));
109c4762a1bSJed Brown     } else {
1109566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD," MUMPS LU:\n"));
1119566063dSJacob Faibussowitsch       PetscCall(MatGetFactor(A,MATSOLVERMUMPS,MAT_FACTOR_LU,&F));
112c4762a1bSJed Brown     }
113c4762a1bSJed Brown     matsolvexx = PETSC_TRUE;
114c4762a1bSJed Brown     if (test_mumps_opts) {
115c4762a1bSJed Brown       /* test mumps options */
116c4762a1bSJed Brown       PetscInt  icntl;
117c4762a1bSJed Brown       PetscReal cntl;
118c4762a1bSJed Brown 
119c4762a1bSJed Brown       icntl = 2;        /* sequential matrix ordering */
1209566063dSJacob Faibussowitsch       PetscCall(MatMumpsSetIcntl(F,7,icntl));
121c4762a1bSJed Brown 
122c4762a1bSJed Brown       cntl = 1.e-6; /* threshold for row pivot detection */
1239566063dSJacob Faibussowitsch       PetscCall(MatMumpsSetIcntl(F,24,1));
1249566063dSJacob Faibussowitsch       PetscCall(MatMumpsSetCntl(F,3,cntl));
125c4762a1bSJed Brown     }
126c4762a1bSJed Brown     break;
127c4762a1bSJed Brown #endif
128c4762a1bSJed Brown #if defined(PETSC_HAVE_MKL_PARDISO)
129c4762a1bSJed Brown   case 3:
130c4762a1bSJed Brown     if (chol) {
1319566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD," MKL_PARDISO CHOLESKY:\n"));
1329566063dSJacob Faibussowitsch       PetscCall(MatGetFactor(A,MATSOLVERMKL_PARDISO,MAT_FACTOR_CHOLESKY,&F));
133c4762a1bSJed Brown     } else {
1349566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD," MKL_PARDISO LU:\n"));
1359566063dSJacob Faibussowitsch       PetscCall(MatGetFactor(A,MATSOLVERMKL_PARDISO,MAT_FACTOR_LU,&F));
136c4762a1bSJed Brown     }
137c4762a1bSJed Brown     break;
138c4762a1bSJed Brown #endif
13938a8e8c1SStefano Zampini #if defined(PETSC_HAVE_CUDA)
14038a8e8c1SStefano Zampini   case 4:
14138a8e8c1SStefano Zampini     if (chol) {
1429566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD," CUSPARSE CHOLESKY:\n"));
1439566063dSJacob Faibussowitsch       PetscCall(MatGetFactor(A,MATSOLVERCUSPARSE,MAT_FACTOR_CHOLESKY,&F));
14438a8e8c1SStefano Zampini     } else {
1459566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD," CUSPARSE LU:\n"));
1469566063dSJacob Faibussowitsch       PetscCall(MatGetFactor(A,MATSOLVERCUSPARSE,MAT_FACTOR_LU,&F));
14738a8e8c1SStefano Zampini     }
14838a8e8c1SStefano Zampini     break;
14938a8e8c1SStefano Zampini #endif
150c4762a1bSJed Brown   default:
151c4762a1bSJed Brown     if (chol) {
1529566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD," PETSC CHOLESKY:\n"));
1539566063dSJacob Faibussowitsch       PetscCall(MatGetFactor(A,MATSOLVERPETSC,MAT_FACTOR_CHOLESKY,&F));
154c4762a1bSJed Brown     } else {
1559566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD," PETSC LU:\n"));
1569566063dSJacob Faibussowitsch       PetscCall(MatGetFactor(A,MATSOLVERPETSC,MAT_FACTOR_LU,&F));
157c4762a1bSJed Brown     }
158c4762a1bSJed Brown     matsolvexx = PETSC_TRUE;
159c4762a1bSJed Brown   }
160c4762a1bSJed Brown 
1619566063dSJacob Faibussowitsch   PetscCall(MatFactorInfoInitialize(&info));
162c4762a1bSJed Brown   info.fill      = 5.0;
163c4762a1bSJed Brown   info.shifttype = (PetscReal) MAT_SHIFT_NONE;
164c4762a1bSJed Brown   if (chol) {
1659566063dSJacob Faibussowitsch     PetscCall(MatCholeskyFactorSymbolic(F,A,perm,&info));
166c4762a1bSJed Brown   } else {
1679566063dSJacob Faibussowitsch     PetscCall(MatLUFactorSymbolic(F,A,perm,iperm,&info));
168c4762a1bSJed Brown   }
169c4762a1bSJed Brown 
170c4762a1bSJed Brown   for (nfact = 0; nfact < 2; nfact++) {
171c4762a1bSJed Brown     if (chol) {
1729566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD," %" PetscInt_FMT "-the CHOLESKY numfactorization \n",nfact));
1739566063dSJacob Faibussowitsch       PetscCall(MatCholeskyFactorNumeric(F,A,&info));
174c4762a1bSJed Brown     } else {
1759566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD," %" PetscInt_FMT "-the LU numfactorization \n",nfact));
1769566063dSJacob Faibussowitsch       PetscCall(MatLUFactorNumeric(F,A,&info));
177c4762a1bSJed Brown     }
178c4762a1bSJed Brown     if (view) {
1799566063dSJacob Faibussowitsch       PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO));
1809566063dSJacob Faibussowitsch       PetscCall(MatView(F,PETSC_VIEWER_STDOUT_WORLD));
1819566063dSJacob Faibussowitsch       PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD));
182c4762a1bSJed Brown       view = PETSC_FALSE;
183c4762a1bSJed Brown     }
184c4762a1bSJed Brown 
185c4762a1bSJed Brown #if defined(PETSC_HAVE_SUPERLU_DIST)
186c4762a1bSJed Brown     if (ipack == 1) { /* Test MatSuperluDistGetDiagU()
187c4762a1bSJed Brown        -- input: matrix factor F; output: main diagonal of matrix U on all processes */
188c4762a1bSJed Brown       PetscInt    M;
189c4762a1bSJed Brown       PetscScalar *diag;
190c4762a1bSJed Brown #if !defined(PETSC_USE_COMPLEX)
191c4762a1bSJed Brown       PetscInt nneg,nzero,npos;
192c4762a1bSJed Brown #endif
193c4762a1bSJed Brown 
1949566063dSJacob Faibussowitsch       PetscCall(MatGetSize(F,&M,NULL));
1959566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(M,&diag));
1969566063dSJacob Faibussowitsch       PetscCall(MatSuperluDistGetDiagU(F,diag));
1979566063dSJacob Faibussowitsch       PetscCall(PetscFree(diag));
198c4762a1bSJed Brown 
199c4762a1bSJed Brown #if !defined(PETSC_USE_COMPLEX)
200c4762a1bSJed Brown       /* Test MatGetInertia() */
2019566063dSJacob Faibussowitsch       PetscCall(MatGetInertia(F,&nneg,&nzero,&npos));
2029566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD," MatInertia: nneg: %" PetscInt_FMT ", nzero: %" PetscInt_FMT ", npos: %" PetscInt_FMT "\n",nneg,nzero,npos));
203c4762a1bSJed Brown #endif
204c4762a1bSJed Brown     }
205c4762a1bSJed Brown #endif
206c4762a1bSJed Brown 
207d47f36abSHong Zhang #if defined(PETSC_HAVE_MUMPS)
208d47f36abSHong Zhang     /* mumps interface allows repeated call of MatCholeskyFactorSymbolic(), while the succession calls do nothing */
209d47f36abSHong Zhang     if (ipack == 2) {
210d47f36abSHong Zhang       if (chol) {
2119566063dSJacob Faibussowitsch         PetscCall(MatCholeskyFactorSymbolic(F,A,perm,&info));
2129566063dSJacob Faibussowitsch         PetscCall(MatCholeskyFactorNumeric(F,A,&info));
213d47f36abSHong Zhang       } else {
2149566063dSJacob Faibussowitsch         PetscCall(MatLUFactorSymbolic(F,A,perm,iperm,&info));
2159566063dSJacob Faibussowitsch         PetscCall(MatLUFactorNumeric(F,A,&info));
216d47f36abSHong Zhang       }
217d47f36abSHong Zhang     }
218d47f36abSHong Zhang #endif
219d47f36abSHong Zhang 
220*b18964edSHong Zhang     /* Test MatMatSolve(), A X = B, where B can be dense or sparse */
221c4762a1bSJed Brown     if (testMatMatSolve) {
222c4762a1bSJed Brown       if (!nfact) {
2239566063dSJacob Faibussowitsch         PetscCall(MatMatMult(A,C,MAT_INITIAL_MATRIX,2.0,&RHS));
224c4762a1bSJed Brown       } else {
2259566063dSJacob Faibussowitsch         PetscCall(MatMatMult(A,C,MAT_REUSE_MATRIX,2.0,&RHS));
226c4762a1bSJed Brown       }
227c4762a1bSJed Brown       for (nsolve = 0; nsolve < 2; nsolve++) {
2289566063dSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_WORLD,"   %" PetscInt_FMT "-the MatMatSolve \n",nsolve));
2299566063dSJacob Faibussowitsch         PetscCall(MatMatSolve(F,RHS,X));
230c4762a1bSJed Brown 
231c4762a1bSJed Brown         /* Check the error */
2329566063dSJacob Faibussowitsch         PetscCall(MatAXPY(X,-1.0,C,SAME_NONZERO_PATTERN));
2339566063dSJacob Faibussowitsch         PetscCall(MatNorm(X,NORM_FROBENIUS,&norm));
234*b18964edSHong Zhang         if (norm > tol) PetscCall(PetscPrintf(PETSC_COMM_WORLD,"%" PetscInt_FMT "-the MatMatSolve: Norm of error %g, nsolve %" PetscInt_FMT "\n",nsolve,(double)norm,nsolve));
235c4762a1bSJed Brown       }
236*b18964edSHong Zhang 
237c4762a1bSJed Brown       if (matsolvexx) {
238c4762a1bSJed Brown         /* Test MatMatSolve(F,RHS,RHS), RHS is a dense matrix */
2399566063dSJacob Faibussowitsch         PetscCall(MatCopy(RHS,X,SAME_NONZERO_PATTERN));
2409566063dSJacob Faibussowitsch         PetscCall(MatMatSolve(F,X,X));
241c4762a1bSJed Brown         /* Check the error */
2429566063dSJacob Faibussowitsch         PetscCall(MatAXPY(X,-1.0,C,SAME_NONZERO_PATTERN));
2439566063dSJacob Faibussowitsch         PetscCall(MatNorm(X,NORM_FROBENIUS,&norm));
244*b18964edSHong Zhang         if (norm > tol) PetscCall(PetscPrintf(PETSC_COMM_WORLD,"MatMatSolve(F,RHS,RHS): Norm of error %g\n",(double)norm));
245c4762a1bSJed Brown       }
246c4762a1bSJed Brown 
247c4762a1bSJed Brown       if (ipack == 2 && size == 1) {
248c4762a1bSJed Brown         Mat spRHS,spRHST,RHST;
249c4762a1bSJed Brown 
2509566063dSJacob Faibussowitsch         PetscCall(MatTranspose(RHS,MAT_INITIAL_MATRIX,&RHST));
2519566063dSJacob Faibussowitsch         PetscCall(MatConvert(RHST,MATAIJ,MAT_INITIAL_MATRIX,&spRHST));
2529566063dSJacob Faibussowitsch         PetscCall(MatCreateTranspose(spRHST,&spRHS));
253c4762a1bSJed Brown         for (nsolve = 0; nsolve < 2; nsolve++) {
2549566063dSJacob Faibussowitsch           PetscCall(PetscPrintf(PETSC_COMM_WORLD,"   %" PetscInt_FMT "-the sparse MatMatSolve \n",nsolve));
2559566063dSJacob Faibussowitsch           PetscCall(MatMatSolve(F,spRHS,X));
256c4762a1bSJed Brown 
257c4762a1bSJed Brown           /* Check the error */
2589566063dSJacob Faibussowitsch           PetscCall(MatAXPY(X,-1.0,C,SAME_NONZERO_PATTERN));
2599566063dSJacob Faibussowitsch           PetscCall(MatNorm(X,NORM_FROBENIUS,&norm));
260*b18964edSHong Zhang           if (norm > tol) PetscCall(PetscPrintf(PETSC_COMM_WORLD,"%" PetscInt_FMT "-the sparse MatMatSolve: Norm of error %g, nsolve %" PetscInt_FMT "\n",nsolve,(double)norm,nsolve));
261c4762a1bSJed Brown         }
262*b18964edSHong Zhang         PetscCall(MatDestroy(&spRHST));
263*b18964edSHong Zhang         PetscCall(MatDestroy(&spRHS));
264*b18964edSHong Zhang         PetscCall(MatDestroy(&RHST));
265*b18964edSHong Zhang       }
266*b18964edSHong Zhang     }
267*b18964edSHong Zhang 
268*b18964edSHong Zhang     /* Test testMatMatSolveTranspose(), A^T X = B, where B can be dense or sparse */
269*b18964edSHong Zhang     if (testMatMatSolveTranspose) {
270*b18964edSHong Zhang       if (!nfact) {
271*b18964edSHong Zhang         PetscCall(MatTransposeMatMult(A,C,MAT_INITIAL_MATRIX,2.0,&RHS1));
272*b18964edSHong Zhang       } else {
273*b18964edSHong Zhang         PetscCall(MatTransposeMatMult(A,C,MAT_REUSE_MATRIX,2.0,&RHS1));
274*b18964edSHong Zhang       }
275*b18964edSHong Zhang 
276*b18964edSHong Zhang       for (nsolve = 0; nsolve < 2; nsolve++) {
277*b18964edSHong Zhang         PetscCall(MatMatSolveTranspose(F,RHS1,X));
278*b18964edSHong Zhang 
279*b18964edSHong Zhang         /* Check the error */
280*b18964edSHong Zhang         PetscCall(MatAXPY(X,-1.0,C,SAME_NONZERO_PATTERN));
281*b18964edSHong Zhang         PetscCall(MatNorm(X,NORM_FROBENIUS,&norm));
282*b18964edSHong Zhang         if (norm > tol) PetscCall(PetscPrintf(PETSC_COMM_WORLD,"%" PetscInt_FMT "-the MatMatSolveTranspose: Norm of error %g, nsolve %" PetscInt_FMT "\n",nsolve,(double)norm,nsolve));
283*b18964edSHong Zhang       }
284*b18964edSHong Zhang 
285*b18964edSHong Zhang       if (ipack == 2 && size == 1) {
286*b18964edSHong Zhang         Mat spRHS,spRHST,RHST;
287*b18964edSHong Zhang 
288*b18964edSHong Zhang         PetscCall(MatTranspose(RHS1,MAT_INITIAL_MATRIX,&RHST));
289*b18964edSHong Zhang         PetscCall(MatConvert(RHST,MATAIJ,MAT_INITIAL_MATRIX,&spRHST));
290*b18964edSHong Zhang         PetscCall(MatCreateTranspose(spRHST,&spRHS));
291*b18964edSHong Zhang         for (nsolve = 0; nsolve < 2; nsolve++) {
292*b18964edSHong Zhang           PetscCall(MatMatSolveTranspose(F,spRHS,X));
293*b18964edSHong Zhang 
294*b18964edSHong Zhang           /* Check the error */
295*b18964edSHong Zhang           PetscCall(MatAXPY(X,-1.0,C,SAME_NONZERO_PATTERN));
296*b18964edSHong Zhang           PetscCall(MatNorm(X,NORM_FROBENIUS,&norm));
297*b18964edSHong Zhang           if (norm > tol) PetscCall(PetscPrintf(PETSC_COMM_WORLD,"%" PetscInt_FMT "-the sparse MatMatSolveTranspose: Norm of error %g, nsolve %" PetscInt_FMT "\n",nsolve,(double)norm,nsolve));
298c4762a1bSJed Brown         }
2999566063dSJacob Faibussowitsch         PetscCall(MatDestroy(&spRHST));
3009566063dSJacob Faibussowitsch         PetscCall(MatDestroy(&spRHS));
3019566063dSJacob Faibussowitsch         PetscCall(MatDestroy(&RHST));
302c4762a1bSJed Brown       }
303c4762a1bSJed Brown     }
304c4762a1bSJed Brown 
305c4762a1bSJed Brown     /* Test MatSolve() */
306c4762a1bSJed Brown     if (testMatSolve) {
307c4762a1bSJed Brown       for (nsolve = 0; nsolve < 2; nsolve++) {
3089566063dSJacob Faibussowitsch         PetscCall(VecSetRandom(x,rand));
3099566063dSJacob Faibussowitsch         PetscCall(VecCopy(x,u));
3109566063dSJacob Faibussowitsch         PetscCall(MatMult(A,x,b));
311c4762a1bSJed Brown 
3129566063dSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_WORLD,"   %" PetscInt_FMT "-the MatSolve \n",nsolve));
3139566063dSJacob Faibussowitsch         PetscCall(MatSolve(F,b,x));
314c4762a1bSJed Brown 
315c4762a1bSJed Brown         /* Check the error */
3169566063dSJacob Faibussowitsch         PetscCall(VecAXPY(u,-1.0,x));  /* u <- (-1.0)x + u */
3179566063dSJacob Faibussowitsch         PetscCall(VecNorm(u,NORM_2,&norm));
318c4762a1bSJed Brown         if (norm > tol) {
319c4762a1bSJed Brown           PetscReal resi;
3209566063dSJacob Faibussowitsch           PetscCall(MatMult(A,x,u)); /* u = A*x */
3219566063dSJacob Faibussowitsch           PetscCall(VecAXPY(u,-1.0,b));  /* u <- (-1.0)b + u */
3229566063dSJacob Faibussowitsch           PetscCall(VecNorm(u,NORM_2,&resi));
3239566063dSJacob Faibussowitsch           PetscCall(PetscPrintf(PETSC_COMM_WORLD,"MatSolve: Norm of error %g, resi %g, numfact %" PetscInt_FMT "\n",(double)norm,(double)resi,nfact));
324c4762a1bSJed Brown         }
325c4762a1bSJed Brown       }
326c4762a1bSJed Brown     }
327c4762a1bSJed Brown   }
328c4762a1bSJed Brown 
329c4762a1bSJed Brown   /* Free data structures */
3309566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
3319566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&C));
3329566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&F));
3339566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&X));
3349566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&RHS));
335*b18964edSHong Zhang   PetscCall(MatDestroy(&RHS1));
336c4762a1bSJed Brown 
3379566063dSJacob Faibussowitsch   PetscCall(PetscRandomDestroy(&rand));
3389566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&perm));
3399566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&iperm));
3409566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
3419566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&b));
3429566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&u));
3439566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
344b122ec5aSJacob Faibussowitsch   return 0;
345c4762a1bSJed Brown }
346c4762a1bSJed Brown 
347c4762a1bSJed Brown /*TEST
348c4762a1bSJed Brown 
349c4762a1bSJed Brown    test:
350dfd57a17SPierre Jolivet       requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
351c4762a1bSJed Brown       args: -f ${DATAFILESPATH}/matrices/small -mat_solver_type 10
352c4762a1bSJed Brown       output_file: output/ex125.out
353c4762a1bSJed Brown 
354c4762a1bSJed Brown    test:
3559a14fc28SStefano Zampini       suffix: 2
3569a14fc28SStefano Zampini       args: -mat_solver_type 10
3579a14fc28SStefano Zampini       output_file: output/ex125.out
3589a14fc28SStefano Zampini 
3599a14fc28SStefano Zampini    test:
360c4762a1bSJed Brown       suffix: mkl_pardiso
361dfd57a17SPierre Jolivet       requires: mkl_pardiso datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
362c4762a1bSJed Brown       args: -f ${DATAFILESPATH}/matrices/small -mat_solver_type 3
363c4762a1bSJed Brown 
364c4762a1bSJed Brown    test:
3659a14fc28SStefano Zampini       suffix: mkl_pardiso_2
3669a14fc28SStefano Zampini       requires: mkl_pardiso
3679a14fc28SStefano Zampini       args: -mat_solver_type 3
3689a14fc28SStefano Zampini       output_file: output/ex125_mkl_pardiso.out
3699a14fc28SStefano Zampini 
3709a14fc28SStefano Zampini    test:
371c4762a1bSJed Brown       suffix: mumps
372dfd57a17SPierre Jolivet       requires: mumps datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
373c4762a1bSJed Brown       args: -f ${DATAFILESPATH}/matrices/small -mat_solver_type 2
374c4762a1bSJed Brown       output_file: output/ex125_mumps_seq.out
375c4762a1bSJed Brown 
376c4762a1bSJed Brown    test:
377c4762a1bSJed Brown       suffix: mumps_2
378c4762a1bSJed Brown       nsize: 3
379dfd57a17SPierre Jolivet       requires: mumps datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
380c4762a1bSJed Brown       args: -f ${DATAFILESPATH}/matrices/small -mat_solver_type 2
381c4762a1bSJed Brown       output_file: output/ex125_mumps_par.out
382c4762a1bSJed Brown 
383c4762a1bSJed Brown    test:
3849a14fc28SStefano Zampini       suffix: mumps_3
3859a14fc28SStefano Zampini       requires: mumps
3869a14fc28SStefano Zampini       args: -mat_solver_type 2
3879a14fc28SStefano Zampini       output_file: output/ex125_mumps_seq.out
3889a14fc28SStefano Zampini 
3899a14fc28SStefano Zampini    test:
3909a14fc28SStefano Zampini       suffix: mumps_4
3919a14fc28SStefano Zampini       nsize: 3
3929a14fc28SStefano Zampini       requires: mumps
3939a14fc28SStefano Zampini       args: -mat_solver_type 2
3949a14fc28SStefano Zampini       output_file: output/ex125_mumps_par.out
3959a14fc28SStefano Zampini 
3969a14fc28SStefano Zampini    test:
397d47f36abSHong Zhang       suffix: mumps_5
398d47f36abSHong Zhang       nsize: 3
399d47f36abSHong Zhang       requires: mumps
400d47f36abSHong Zhang       args: -mat_solver_type 2 -cholesky
401d47f36abSHong Zhang       output_file: output/ex125_mumps_par_cholesky.out
402d47f36abSHong Zhang 
403d47f36abSHong Zhang    test:
404c4762a1bSJed Brown       suffix: superlu_dist
4059a14fc28SStefano Zampini       nsize: {{1 3}}
406d4783600SBarry Smith       requires: datafilespath double !complex !defined(PETSC_USE_64BIT_INDICES) superlu_dist
407c4762a1bSJed Brown       args: -f ${DATAFILESPATH}/matrices/small -mat_solver_type 1 -mat_superlu_dist_rowperm NOROWPERM
408c4762a1bSJed Brown 
409c4762a1bSJed Brown    test:
410c4762a1bSJed Brown       suffix: superlu_dist_2
4119a14fc28SStefano Zampini       nsize: {{1 3}}
4129a14fc28SStefano Zampini       requires: superlu_dist !complex
4139a14fc28SStefano Zampini       args: -n 36 -mat_solver_type 1 -mat_superlu_dist_rowperm NOROWPERM
414c4762a1bSJed Brown       output_file: output/ex125_superlu_dist.out
415c4762a1bSJed Brown 
416c4762a1bSJed Brown    test:
417c4762a1bSJed Brown       suffix: superlu_dist_complex
418c4762a1bSJed Brown       nsize: 3
419d4783600SBarry Smith       requires: datafilespath double superlu_dist complex !defined(PETSC_USE_64BIT_INDICES)
420c4762a1bSJed Brown       args: -f ${DATAFILESPATH}/matrices/farzad_B_rhs -mat_solver_type 1
421c4762a1bSJed Brown       output_file: output/ex125_superlu_dist_complex.out
422c4762a1bSJed Brown 
42338a8e8c1SStefano Zampini    test:
4249a14fc28SStefano Zampini       suffix: superlu_dist_complex_2
4259a14fc28SStefano Zampini       nsize: 3
4269a14fc28SStefano Zampini       requires: superlu_dist complex
4279a14fc28SStefano Zampini       args: -mat_solver_type 1
4289a14fc28SStefano Zampini       output_file: output/ex125_superlu_dist_complex.out
4299a14fc28SStefano Zampini 
4309a14fc28SStefano Zampini    test:
43138a8e8c1SStefano Zampini       suffix: cusparse
432dfd57a17SPierre Jolivet       requires: cuda datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
43338a8e8c1SStefano Zampini       args: -mat_type aijcusparse -f ${DATAFILESPATH}/matrices/small -mat_solver_type 4 -cholesky {{0 1}separate output}
43438a8e8c1SStefano Zampini 
4359a14fc28SStefano Zampini    test:
4369a14fc28SStefano Zampini       suffix: cusparse_2
4379a14fc28SStefano Zampini       requires: cuda
4389a14fc28SStefano Zampini       args: -mat_type aijcusparse -mat_solver_type 4 -cholesky {{0 1}separate output}
4399a14fc28SStefano Zampini 
440c4762a1bSJed Brown TEST*/
441