xref: /petsc/src/mat/tests/ex125.c (revision 8fffc7623b428232ccbb2ec708553adbd8d47aa1)
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 {
8c4762a1bSJed Brown   Mat            A,RHS,C,F,X;
9c4762a1bSJed Brown   Vec            u,x,b;
10c4762a1bSJed Brown   PetscErrorCode ierr;
11c4762a1bSJed Brown   PetscMPIInt    size;
12c4762a1bSJed Brown   PetscInt       m,n,nfact,nsolve,nrhs,ipack=0;
13c4762a1bSJed Brown   PetscReal      norm,tol=1.e-10;
14c4762a1bSJed Brown   IS             perm,iperm;
15c4762a1bSJed Brown   MatFactorInfo  info;
16c4762a1bSJed Brown   PetscRandom    rand;
17c4762a1bSJed Brown   PetscBool      flg,testMatSolve=PETSC_TRUE,testMatMatSolve=PETSC_TRUE;
18c4762a1bSJed Brown   PetscBool      chol=PETSC_FALSE,view=PETSC_FALSE,matsolvexx = PETSC_FALSE;
19c4762a1bSJed Brown #if defined(PETSC_HAVE_MUMPS)
20c4762a1bSJed Brown   PetscBool      test_mumps_opts=PETSC_FALSE;
21c4762a1bSJed Brown #endif
22c4762a1bSJed Brown   PetscViewer    fd;              /* viewer */
23c4762a1bSJed Brown   char           file[PETSC_MAX_PATH_LEN]; /* input file name */
24c4762a1bSJed Brown 
25c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
26ffc4695bSBarry Smith   ierr = MPI_Comm_size(PETSC_COMM_WORLD, &size);CHKERRMPI(ierr);
27c4762a1bSJed Brown 
28c4762a1bSJed Brown   /* Determine file from which we read the matrix A */
29589a23caSBarry Smith   ierr = PetscOptionsGetString(NULL,NULL,"-f",file,sizeof(file),&flg);CHKERRQ(ierr);
309a14fc28SStefano Zampini   if (flg) { /* Load matrix A */
31c4762a1bSJed Brown     ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&fd);CHKERRQ(ierr);
32c4762a1bSJed Brown     ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr);
33c4762a1bSJed Brown     ierr = MatSetFromOptions(A);CHKERRQ(ierr);
34c4762a1bSJed Brown     ierr = MatLoad(A,fd);CHKERRQ(ierr);
35c4762a1bSJed Brown     ierr = PetscViewerDestroy(&fd);CHKERRQ(ierr);
369a14fc28SStefano Zampini   } else {
379a14fc28SStefano Zampini     n = 13;
389a14fc28SStefano Zampini     ierr = PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);CHKERRQ(ierr);
399a14fc28SStefano Zampini     ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr);
409a14fc28SStefano Zampini     ierr = MatSetType(A,MATAIJ);CHKERRQ(ierr);
419a14fc28SStefano Zampini     ierr = MatSetFromOptions(A);CHKERRQ(ierr);
429a14fc28SStefano Zampini     ierr = MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n);CHKERRQ(ierr);
439a14fc28SStefano Zampini     ierr = MatSetUp(A);CHKERRQ(ierr);
449a14fc28SStefano Zampini     ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
459a14fc28SStefano Zampini     ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
469a14fc28SStefano Zampini     ierr = MatShift(A,1.0);CHKERRQ(ierr);
479a14fc28SStefano Zampini   }
48c4762a1bSJed Brown   ierr = MatGetLocalSize(A,&m,&n);CHKERRQ(ierr);
49*8fffc762SJacob Faibussowitsch   if (m != n) SETERRQ2(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() */
52c4762a1bSJed Brown   ierr = MatIsSymmetric(A,0.0,&flg);CHKERRQ(ierr);
53c4762a1bSJed Brown 
54c4762a1bSJed Brown   ierr = MatViewFromOptions(A,NULL,"-A_view");CHKERRQ(ierr);
55c4762a1bSJed Brown 
56a5b23f4aSJose E. Roman   /* Create dense matrix C and X; C holds true solution with identical columns */
57c4762a1bSJed Brown   nrhs = 2;
58c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL,"-nrhs",&nrhs,NULL);CHKERRQ(ierr);
59*8fffc762SJacob Faibussowitsch   ierr = PetscPrintf(PETSC_COMM_WORLD,"ex125: nrhs %" PetscInt_FMT "\n",nrhs);CHKERRQ(ierr);
60c4762a1bSJed Brown   ierr = MatCreate(PETSC_COMM_WORLD,&C);CHKERRQ(ierr);
61c4762a1bSJed Brown   ierr = MatSetOptionsPrefix(C,"rhs_");CHKERRQ(ierr);
62c4762a1bSJed Brown   ierr = MatSetSizes(C,m,PETSC_DECIDE,PETSC_DECIDE,nrhs);CHKERRQ(ierr);
63c4762a1bSJed Brown   ierr = MatSetType(C,MATDENSE);CHKERRQ(ierr);
64c4762a1bSJed Brown   ierr = MatSetFromOptions(C);CHKERRQ(ierr);
65c4762a1bSJed Brown   ierr = MatSetUp(C);CHKERRQ(ierr);
66c4762a1bSJed Brown 
67c4762a1bSJed Brown   ierr = PetscOptionsGetBool(NULL,NULL,"-view_factor",&view,NULL);CHKERRQ(ierr);
68c4762a1bSJed Brown   ierr = PetscOptionsGetBool(NULL,NULL,"-test_matmatsolve",&testMatMatSolve,NULL);CHKERRQ(ierr);
69c4762a1bSJed Brown   ierr = PetscOptionsGetBool(NULL,NULL,"-cholesky",&chol,NULL);CHKERRQ(ierr);
70c4762a1bSJed Brown #if defined(PETSC_HAVE_MUMPS)
71c4762a1bSJed Brown   ierr = PetscOptionsGetBool(NULL,NULL,"-test_mumps_opts",&test_mumps_opts,NULL);CHKERRQ(ierr);
72c4762a1bSJed Brown #endif
73c4762a1bSJed Brown 
74c4762a1bSJed Brown   ierr = PetscRandomCreate(PETSC_COMM_WORLD,&rand);CHKERRQ(ierr);
75c4762a1bSJed Brown   ierr = PetscRandomSetFromOptions(rand);CHKERRQ(ierr);
76c4762a1bSJed Brown   ierr = MatSetRandom(C,rand);CHKERRQ(ierr);
77c4762a1bSJed Brown   ierr = MatDuplicate(C,MAT_DO_NOT_COPY_VALUES,&X);CHKERRQ(ierr);
78c4762a1bSJed Brown 
79c4762a1bSJed Brown   /* Create vectors */
8038a8e8c1SStefano Zampini   ierr = MatCreateVecs(A,&x,&b);CHKERRQ(ierr);
81c4762a1bSJed Brown   ierr = VecDuplicate(x,&u);CHKERRQ(ierr); /* save the true solution */
82c4762a1bSJed Brown 
83c4762a1bSJed Brown   /* Test Factorization */
84c4762a1bSJed Brown   ierr = MatGetOrdering(A,MATORDERINGND,&perm,&iperm);CHKERRQ(ierr);
85c4762a1bSJed Brown 
86c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL,"-mat_solver_type",&ipack,NULL);CHKERRQ(ierr);
87c4762a1bSJed Brown   switch (ipack) {
88c4762a1bSJed Brown #if defined(PETSC_HAVE_SUPERLU)
89c4762a1bSJed Brown   case 0:
90c4762a1bSJed Brown     if (chol) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"SuperLU does not provide Cholesky!");
91c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD," SUPERLU LU:\n");CHKERRQ(ierr);
92c4762a1bSJed Brown     ierr = MatGetFactor(A,MATSOLVERSUPERLU,MAT_FACTOR_LU,&F);CHKERRQ(ierr);
93c4762a1bSJed Brown     matsolvexx = PETSC_TRUE;
94c4762a1bSJed Brown     break;
95c4762a1bSJed Brown #endif
96c4762a1bSJed Brown #if defined(PETSC_HAVE_SUPERLU_DIST)
97c4762a1bSJed Brown   case 1:
98c4762a1bSJed Brown     if (chol) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"SuperLU does not provide Cholesky!");
99c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD," SUPERLU_DIST LU:\n");CHKERRQ(ierr);
100c4762a1bSJed Brown     ierr = MatGetFactor(A,MATSOLVERSUPERLU_DIST,MAT_FACTOR_LU,&F);CHKERRQ(ierr);
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) {
107c4762a1bSJed Brown       ierr = PetscPrintf(PETSC_COMM_WORLD," MUMPS CHOLESKY:\n");CHKERRQ(ierr);
108c4762a1bSJed Brown       ierr = MatGetFactor(A,MATSOLVERMUMPS,MAT_FACTOR_CHOLESKY,&F);CHKERRQ(ierr);
109c4762a1bSJed Brown     } else {
110c4762a1bSJed Brown       ierr = PetscPrintf(PETSC_COMM_WORLD," MUMPS LU:\n");CHKERRQ(ierr);
111c4762a1bSJed Brown       ierr = MatGetFactor(A,MATSOLVERMUMPS,MAT_FACTOR_LU,&F);CHKERRQ(ierr);
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 */
120c4762a1bSJed Brown       ierr  = MatMumpsSetIcntl(F,7,icntl);CHKERRQ(ierr);
121c4762a1bSJed Brown 
122c4762a1bSJed Brown       cntl = 1.e-6; /* threshold for row pivot detection */
123c4762a1bSJed Brown       ierr = MatMumpsSetIcntl(F,24,1);CHKERRQ(ierr);
124c4762a1bSJed Brown       ierr = MatMumpsSetCntl(F,3,cntl);CHKERRQ(ierr);
125c4762a1bSJed Brown     }
126c4762a1bSJed Brown     break;
127c4762a1bSJed Brown #endif
128c4762a1bSJed Brown #if defined(PETSC_HAVE_MKL_PARDISO)
129c4762a1bSJed Brown   case 3:
130c4762a1bSJed Brown     if (chol) {
131c4762a1bSJed Brown       ierr = PetscPrintf(PETSC_COMM_WORLD," MKL_PARDISO CHOLESKY:\n");CHKERRQ(ierr);
132c4762a1bSJed Brown       ierr = MatGetFactor(A,MATSOLVERMKL_PARDISO,MAT_FACTOR_CHOLESKY,&F);CHKERRQ(ierr);
133c4762a1bSJed Brown     } else {
134c4762a1bSJed Brown       ierr = PetscPrintf(PETSC_COMM_WORLD," MKL_PARDISO LU:\n");CHKERRQ(ierr);
135c4762a1bSJed Brown       ierr = MatGetFactor(A,MATSOLVERMKL_PARDISO,MAT_FACTOR_LU,&F);CHKERRQ(ierr);
136c4762a1bSJed Brown     }
137c4762a1bSJed Brown     break;
138c4762a1bSJed Brown #endif
13938a8e8c1SStefano Zampini #if defined(PETSC_HAVE_CUDA)
14038a8e8c1SStefano Zampini   case 4:
14138a8e8c1SStefano Zampini     if (chol) {
14238a8e8c1SStefano Zampini       ierr = PetscPrintf(PETSC_COMM_WORLD," CUSPARSE CHOLESKY:\n");CHKERRQ(ierr);
14338a8e8c1SStefano Zampini       ierr = MatGetFactor(A,MATSOLVERCUSPARSE,MAT_FACTOR_CHOLESKY,&F);CHKERRQ(ierr);
14438a8e8c1SStefano Zampini     } else {
14538a8e8c1SStefano Zampini       ierr = PetscPrintf(PETSC_COMM_WORLD," CUSPARSE LU:\n");CHKERRQ(ierr);
14638a8e8c1SStefano Zampini       ierr = MatGetFactor(A,MATSOLVERCUSPARSE,MAT_FACTOR_LU,&F);CHKERRQ(ierr);
14738a8e8c1SStefano Zampini     }
14838a8e8c1SStefano Zampini     break;
14938a8e8c1SStefano Zampini #endif
150c4762a1bSJed Brown   default:
151c4762a1bSJed Brown     if (chol) {
152c4762a1bSJed Brown       ierr = PetscPrintf(PETSC_COMM_WORLD," PETSC CHOLESKY:\n");CHKERRQ(ierr);
153c4762a1bSJed Brown       ierr = MatGetFactor(A,MATSOLVERPETSC,MAT_FACTOR_CHOLESKY,&F);CHKERRQ(ierr);
154c4762a1bSJed Brown     } else {
155c4762a1bSJed Brown       ierr = PetscPrintf(PETSC_COMM_WORLD," PETSC LU:\n");CHKERRQ(ierr);
156c4762a1bSJed Brown       ierr = MatGetFactor(A,MATSOLVERPETSC,MAT_FACTOR_LU,&F);CHKERRQ(ierr);
157c4762a1bSJed Brown     }
158c4762a1bSJed Brown     matsolvexx = PETSC_TRUE;
159c4762a1bSJed Brown   }
160c4762a1bSJed Brown 
161c4762a1bSJed Brown   ierr           = MatFactorInfoInitialize(&info);CHKERRQ(ierr);
162c4762a1bSJed Brown   info.fill      = 5.0;
163c4762a1bSJed Brown   info.shifttype = (PetscReal) MAT_SHIFT_NONE;
164c4762a1bSJed Brown   if (chol) {
165c4762a1bSJed Brown     ierr = MatCholeskyFactorSymbolic(F,A,perm,&info);CHKERRQ(ierr);
166c4762a1bSJed Brown   } else {
167c4762a1bSJed Brown     ierr = MatLUFactorSymbolic(F,A,perm,iperm,&info);CHKERRQ(ierr);
168c4762a1bSJed Brown   }
169c4762a1bSJed Brown 
170c4762a1bSJed Brown   for (nfact = 0; nfact < 2; nfact++) {
171c4762a1bSJed Brown     if (chol) {
172*8fffc762SJacob Faibussowitsch       ierr = PetscPrintf(PETSC_COMM_WORLD," %" PetscInt_FMT "-the CHOLESKY numfactorization \n",nfact);CHKERRQ(ierr);
173c4762a1bSJed Brown       ierr = MatCholeskyFactorNumeric(F,A,&info);CHKERRQ(ierr);
174c4762a1bSJed Brown     } else {
175*8fffc762SJacob Faibussowitsch       ierr = PetscPrintf(PETSC_COMM_WORLD," %" PetscInt_FMT "-the LU numfactorization \n",nfact);CHKERRQ(ierr);
176c4762a1bSJed Brown       ierr = MatLUFactorNumeric(F,A,&info);CHKERRQ(ierr);
177c4762a1bSJed Brown     }
178c4762a1bSJed Brown     if (view) {
179c4762a1bSJed Brown       ierr = PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr);
180c4762a1bSJed Brown       ierr = MatView(F,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
181c4762a1bSJed Brown       ierr = PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
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 
194c4762a1bSJed Brown       ierr = MatGetSize(F,&M,NULL);CHKERRQ(ierr);
195c4762a1bSJed Brown       ierr = PetscMalloc1(M,&diag);CHKERRQ(ierr);
196c4762a1bSJed Brown       ierr = MatSuperluDistGetDiagU(F,diag);CHKERRQ(ierr);
197c4762a1bSJed Brown       ierr = PetscFree(diag);CHKERRQ(ierr);
198c4762a1bSJed Brown 
199c4762a1bSJed Brown #if !defined(PETSC_USE_COMPLEX)
200c4762a1bSJed Brown       /* Test MatGetInertia() */
201c4762a1bSJed Brown       ierr = MatGetInertia(F,&nneg,&nzero,&npos);CHKERRQ(ierr);
202*8fffc762SJacob Faibussowitsch       ierr = PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD," MatInertia: nneg: %" PetscInt_FMT ", nzero: %" PetscInt_FMT ", npos: %" PetscInt_FMT "\n",nneg,nzero,npos);CHKERRQ(ierr);
203c4762a1bSJed Brown #endif
204c4762a1bSJed Brown     }
205c4762a1bSJed Brown #endif
206c4762a1bSJed Brown 
207c4762a1bSJed Brown     /* Test MatMatSolve() */
208c4762a1bSJed Brown     if (testMatMatSolve) {
209c4762a1bSJed Brown       if (!nfact) {
210c4762a1bSJed Brown         ierr = MatMatMult(A,C,MAT_INITIAL_MATRIX,2.0,&RHS);CHKERRQ(ierr);
211c4762a1bSJed Brown       } else {
212c4762a1bSJed Brown         ierr = MatMatMult(A,C,MAT_REUSE_MATRIX,2.0,&RHS);CHKERRQ(ierr);
213c4762a1bSJed Brown       }
214c4762a1bSJed Brown       for (nsolve = 0; nsolve < 2; nsolve++) {
215*8fffc762SJacob Faibussowitsch         ierr = PetscPrintf(PETSC_COMM_WORLD,"   %" PetscInt_FMT "-the MatMatSolve \n",nsolve);CHKERRQ(ierr);
216c4762a1bSJed Brown         ierr = MatMatSolve(F,RHS,X);CHKERRQ(ierr);
217c4762a1bSJed Brown 
218c4762a1bSJed Brown         /* Check the error */
219c4762a1bSJed Brown         ierr = MatAXPY(X,-1.0,C,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
220c4762a1bSJed Brown         ierr = MatNorm(X,NORM_FROBENIUS,&norm);CHKERRQ(ierr);
221c4762a1bSJed Brown         if (norm > tol) {
222*8fffc762SJacob Faibussowitsch           ierr = PetscPrintf(PETSC_COMM_WORLD,"%" PetscInt_FMT "-the MatMatSolve: Norm of error %g, nsolve %" PetscInt_FMT "\n",nsolve,(double)norm,nsolve);CHKERRQ(ierr);
223c4762a1bSJed Brown         }
224c4762a1bSJed Brown       }
225c4762a1bSJed Brown       if (matsolvexx) {
226c4762a1bSJed Brown         /* Test MatMatSolve(F,RHS,RHS), RHS is a dense matrix */
227c4762a1bSJed Brown         ierr = MatCopy(RHS,X,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
228c4762a1bSJed Brown         ierr = MatMatSolve(F,X,X);CHKERRQ(ierr);
229c4762a1bSJed Brown         /* Check the error */
230c4762a1bSJed Brown         ierr = MatAXPY(X,-1.0,C,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
231c4762a1bSJed Brown         ierr = MatNorm(X,NORM_FROBENIUS,&norm);CHKERRQ(ierr);
232c4762a1bSJed Brown         if (norm > tol) {
233c4762a1bSJed Brown           ierr = PetscPrintf(PETSC_COMM_WORLD,"MatMatSolve(F,RHS,RHS): Norm of error %g\n",(double)norm);CHKERRQ(ierr);
234c4762a1bSJed Brown         }
235c4762a1bSJed Brown       }
236c4762a1bSJed Brown 
237c4762a1bSJed Brown       if (ipack == 2 && size == 1) {
238c4762a1bSJed Brown         Mat spRHS,spRHST,RHST;
239c4762a1bSJed Brown 
240c4762a1bSJed Brown         ierr = MatTranspose(RHS,MAT_INITIAL_MATRIX,&RHST);CHKERRQ(ierr);
241c4762a1bSJed Brown         ierr = MatConvert(RHST,MATAIJ,MAT_INITIAL_MATRIX,&spRHST);CHKERRQ(ierr);
242c4762a1bSJed Brown         ierr = MatCreateTranspose(spRHST,&spRHS);CHKERRQ(ierr);
243c4762a1bSJed Brown         for (nsolve = 0; nsolve < 2; nsolve++) {
244*8fffc762SJacob Faibussowitsch           ierr = PetscPrintf(PETSC_COMM_WORLD,"   %" PetscInt_FMT "-the sparse MatMatSolve \n",nsolve);CHKERRQ(ierr);
245c4762a1bSJed Brown           ierr = MatMatSolve(F,spRHS,X);CHKERRQ(ierr);
246c4762a1bSJed Brown 
247c4762a1bSJed Brown           /* Check the error */
248c4762a1bSJed Brown           ierr = MatAXPY(X,-1.0,C,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
249c4762a1bSJed Brown           ierr = MatNorm(X,NORM_FROBENIUS,&norm);CHKERRQ(ierr);
250c4762a1bSJed Brown           if (norm > tol) {
251*8fffc762SJacob Faibussowitsch             ierr = PetscPrintf(PETSC_COMM_WORLD,"%" PetscInt_FMT "-the sparse MatMatSolve: Norm of error %g, nsolve %" PetscInt_FMT "\n",nsolve,(double)norm,nsolve);CHKERRQ(ierr);
252c4762a1bSJed Brown           }
253c4762a1bSJed Brown         }
254c4762a1bSJed Brown         ierr = MatDestroy(&spRHST);CHKERRQ(ierr);
255c4762a1bSJed Brown         ierr = MatDestroy(&spRHS);CHKERRQ(ierr);
256c4762a1bSJed Brown         ierr = MatDestroy(&RHST);CHKERRQ(ierr);
257c4762a1bSJed Brown       }
258c4762a1bSJed Brown     }
259c4762a1bSJed Brown 
260c4762a1bSJed Brown     /* Test MatSolve() */
261c4762a1bSJed Brown     if (testMatSolve) {
262c4762a1bSJed Brown       for (nsolve = 0; nsolve < 2; nsolve++) {
263c4762a1bSJed Brown         ierr = VecSetRandom(x,rand);CHKERRQ(ierr);
264c4762a1bSJed Brown         ierr = VecCopy(x,u);CHKERRQ(ierr);
265c4762a1bSJed Brown         ierr = MatMult(A,x,b);CHKERRQ(ierr);
266c4762a1bSJed Brown 
267*8fffc762SJacob Faibussowitsch         ierr = PetscPrintf(PETSC_COMM_WORLD,"   %" PetscInt_FMT "-the MatSolve \n",nsolve);CHKERRQ(ierr);
268c4762a1bSJed Brown         ierr = MatSolve(F,b,x);CHKERRQ(ierr);
269c4762a1bSJed Brown 
270c4762a1bSJed Brown         /* Check the error */
271c4762a1bSJed Brown         ierr = VecAXPY(u,-1.0,x);CHKERRQ(ierr);  /* u <- (-1.0)x + u */
272c4762a1bSJed Brown         ierr = VecNorm(u,NORM_2,&norm);CHKERRQ(ierr);
273c4762a1bSJed Brown         if (norm > tol) {
274c4762a1bSJed Brown           PetscReal resi;
275c4762a1bSJed Brown           ierr = MatMult(A,x,u);CHKERRQ(ierr); /* u = A*x */
276c4762a1bSJed Brown           ierr = VecAXPY(u,-1.0,b);CHKERRQ(ierr);  /* u <- (-1.0)b + u */
277c4762a1bSJed Brown           ierr = VecNorm(u,NORM_2,&resi);CHKERRQ(ierr);
278*8fffc762SJacob Faibussowitsch           ierr = PetscPrintf(PETSC_COMM_WORLD,"MatSolve: Norm of error %g, resi %g, numfact %" PetscInt_FMT "\n",(double)norm,(double)resi,nfact);CHKERRQ(ierr);
279c4762a1bSJed Brown         }
280c4762a1bSJed Brown       }
281c4762a1bSJed Brown     }
282c4762a1bSJed Brown   }
283c4762a1bSJed Brown 
284c4762a1bSJed Brown   /* Free data structures */
285c4762a1bSJed Brown   ierr = MatDestroy(&A);CHKERRQ(ierr);
286c4762a1bSJed Brown   ierr = MatDestroy(&C);CHKERRQ(ierr);
287c4762a1bSJed Brown   ierr = MatDestroy(&F);CHKERRQ(ierr);
288c4762a1bSJed Brown   ierr = MatDestroy(&X);CHKERRQ(ierr);
289c4762a1bSJed Brown   if (testMatMatSolve) {
290c4762a1bSJed Brown     ierr = MatDestroy(&RHS);CHKERRQ(ierr);
291c4762a1bSJed Brown   }
292c4762a1bSJed Brown 
293c4762a1bSJed Brown   ierr = PetscRandomDestroy(&rand);CHKERRQ(ierr);
294c4762a1bSJed Brown   ierr = ISDestroy(&perm);CHKERRQ(ierr);
295c4762a1bSJed Brown   ierr = ISDestroy(&iperm);CHKERRQ(ierr);
296c4762a1bSJed Brown   ierr = VecDestroy(&x);CHKERRQ(ierr);
297c4762a1bSJed Brown   ierr = VecDestroy(&b);CHKERRQ(ierr);
298c4762a1bSJed Brown   ierr = VecDestroy(&u);CHKERRQ(ierr);
299c4762a1bSJed Brown   ierr = PetscFinalize();
300c4762a1bSJed Brown   return ierr;
301c4762a1bSJed Brown }
302c4762a1bSJed Brown 
303c4762a1bSJed Brown /*TEST
304c4762a1bSJed Brown 
305c4762a1bSJed Brown    test:
306dfd57a17SPierre Jolivet       requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
307c4762a1bSJed Brown       args: -f ${DATAFILESPATH}/matrices/small -mat_solver_type 10
308c4762a1bSJed Brown       output_file: output/ex125.out
309c4762a1bSJed Brown 
310c4762a1bSJed Brown    test:
3119a14fc28SStefano Zampini       suffix: 2
3129a14fc28SStefano Zampini       args: -mat_solver_type 10
3139a14fc28SStefano Zampini       output_file: output/ex125.out
3149a14fc28SStefano Zampini 
3159a14fc28SStefano Zampini    test:
316c4762a1bSJed Brown       suffix: mkl_pardiso
317dfd57a17SPierre Jolivet       requires: mkl_pardiso datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
318c4762a1bSJed Brown       args: -f ${DATAFILESPATH}/matrices/small -mat_solver_type 3
319c4762a1bSJed Brown 
320c4762a1bSJed Brown    test:
3219a14fc28SStefano Zampini       suffix: mkl_pardiso_2
3229a14fc28SStefano Zampini       requires: mkl_pardiso
3239a14fc28SStefano Zampini       args: -mat_solver_type 3
3249a14fc28SStefano Zampini       output_file: output/ex125_mkl_pardiso.out
3259a14fc28SStefano Zampini 
3269a14fc28SStefano Zampini    test:
327c4762a1bSJed Brown       suffix: mumps
328dfd57a17SPierre Jolivet       requires: mumps datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
329c4762a1bSJed Brown       args: -f ${DATAFILESPATH}/matrices/small -mat_solver_type 2
330c4762a1bSJed Brown       output_file: output/ex125_mumps_seq.out
331c4762a1bSJed Brown 
332c4762a1bSJed Brown    test:
333c4762a1bSJed Brown       suffix: mumps_2
334c4762a1bSJed Brown       nsize: 3
335dfd57a17SPierre Jolivet       requires: mumps datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
336c4762a1bSJed Brown       args: -f ${DATAFILESPATH}/matrices/small -mat_solver_type 2
337c4762a1bSJed Brown       output_file: output/ex125_mumps_par.out
338c4762a1bSJed Brown 
339c4762a1bSJed Brown    test:
3409a14fc28SStefano Zampini       suffix: mumps_3
3419a14fc28SStefano Zampini       requires: mumps
3429a14fc28SStefano Zampini       args: -mat_solver_type 2
3439a14fc28SStefano Zampini       output_file: output/ex125_mumps_seq.out
3449a14fc28SStefano Zampini 
3459a14fc28SStefano Zampini    test:
3469a14fc28SStefano Zampini       suffix: mumps_4
3479a14fc28SStefano Zampini       nsize: 3
3489a14fc28SStefano Zampini       requires: mumps
3499a14fc28SStefano Zampini       args: -mat_solver_type 2
3509a14fc28SStefano Zampini       output_file: output/ex125_mumps_par.out
3519a14fc28SStefano Zampini 
3529a14fc28SStefano Zampini    test:
353c4762a1bSJed Brown       suffix: superlu_dist
3549a14fc28SStefano Zampini       nsize: {{1 3}}
355dfd57a17SPierre Jolivet       requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES) superlu_dist
356c4762a1bSJed Brown       args: -f ${DATAFILESPATH}/matrices/small -mat_solver_type 1 -mat_superlu_dist_rowperm NOROWPERM
357c4762a1bSJed Brown 
358c4762a1bSJed Brown    test:
359c4762a1bSJed Brown       suffix: superlu_dist_2
3609a14fc28SStefano Zampini       nsize: {{1 3}}
3619a14fc28SStefano Zampini       requires: superlu_dist !complex
3629a14fc28SStefano Zampini       args: -n 36 -mat_solver_type 1 -mat_superlu_dist_rowperm NOROWPERM
363c4762a1bSJed Brown       output_file: output/ex125_superlu_dist.out
364c4762a1bSJed Brown 
365c4762a1bSJed Brown    test:
366c4762a1bSJed Brown       suffix: superlu_dist_complex
367c4762a1bSJed Brown       nsize: 3
368dfd57a17SPierre Jolivet       requires: datafilespath superlu_dist complex double !defined(PETSC_USE_64BIT_INDICES)
369c4762a1bSJed Brown       args: -f ${DATAFILESPATH}/matrices/farzad_B_rhs -mat_solver_type 1
370c4762a1bSJed Brown       output_file: output/ex125_superlu_dist_complex.out
371c4762a1bSJed Brown 
37238a8e8c1SStefano Zampini    test:
3739a14fc28SStefano Zampini       suffix: superlu_dist_complex_2
3749a14fc28SStefano Zampini       nsize: 3
3759a14fc28SStefano Zampini       requires: superlu_dist complex
3769a14fc28SStefano Zampini       args: -mat_solver_type 1
3779a14fc28SStefano Zampini       output_file: output/ex125_superlu_dist_complex.out
3789a14fc28SStefano Zampini 
3799a14fc28SStefano Zampini    test:
38038a8e8c1SStefano Zampini       suffix: cusparse
381dfd57a17SPierre Jolivet       requires: cuda datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
38238a8e8c1SStefano Zampini       args: -mat_type aijcusparse -f ${DATAFILESPATH}/matrices/small -mat_solver_type 4 -cholesky {{0 1}separate output}
38338a8e8c1SStefano Zampini 
3849a14fc28SStefano Zampini    test:
3859a14fc28SStefano Zampini       suffix: cusparse_2
3869a14fc28SStefano Zampini       requires: cuda
3879a14fc28SStefano Zampini       args: -mat_type aijcusparse -mat_solver_type 4 -cholesky {{0 1}separate output}
3889a14fc28SStefano Zampini 
389c4762a1bSJed Brown TEST*/
390