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); 30c4762a1bSJed Brown if (!flg) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_USER,"Must indicate binary file with the -f option"); 31c4762a1bSJed Brown 32c4762a1bSJed Brown /* Load matrix A */ 33c4762a1bSJed Brown ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&fd);CHKERRQ(ierr); 34c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr); 35c4762a1bSJed Brown ierr = MatSetFromOptions(A);CHKERRQ(ierr); 36c4762a1bSJed Brown ierr = MatLoad(A,fd);CHKERRQ(ierr); 37c4762a1bSJed Brown ierr = PetscViewerDestroy(&fd);CHKERRQ(ierr); 38c4762a1bSJed Brown ierr = MatGetLocalSize(A,&m,&n);CHKERRQ(ierr); 39c4762a1bSJed Brown if (m != n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ, "This example is not intended for rectangular matrices (%D, %D)", m, n); 40c4762a1bSJed Brown 41c4762a1bSJed Brown /* if A is symmetric, set its flag -- required by MatGetInertia() */ 42c4762a1bSJed Brown ierr = MatIsSymmetric(A,0.0,&flg);CHKERRQ(ierr); 43c4762a1bSJed Brown 44c4762a1bSJed Brown ierr = MatViewFromOptions(A,NULL,"-A_view");CHKERRQ(ierr); 45c4762a1bSJed Brown 46c4762a1bSJed Brown /* Create dense matrix C and X; C holds true solution with identical colums */ 47c4762a1bSJed Brown nrhs = 2; 48c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-nrhs",&nrhs,NULL);CHKERRQ(ierr); 49c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"ex125: nrhs %D\n",nrhs);CHKERRQ(ierr); 50c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD,&C);CHKERRQ(ierr); 51c4762a1bSJed Brown ierr = MatSetOptionsPrefix(C,"rhs_");CHKERRQ(ierr); 52c4762a1bSJed Brown ierr = MatSetSizes(C,m,PETSC_DECIDE,PETSC_DECIDE,nrhs);CHKERRQ(ierr); 53c4762a1bSJed Brown ierr = MatSetType(C,MATDENSE);CHKERRQ(ierr); 54c4762a1bSJed Brown ierr = MatSetFromOptions(C);CHKERRQ(ierr); 55c4762a1bSJed Brown ierr = MatSetUp(C);CHKERRQ(ierr); 56c4762a1bSJed Brown 57c4762a1bSJed Brown ierr = PetscOptionsGetBool(NULL,NULL,"-view_factor",&view,NULL);CHKERRQ(ierr); 58c4762a1bSJed Brown ierr = PetscOptionsGetBool(NULL,NULL,"-test_matmatsolve",&testMatMatSolve,NULL);CHKERRQ(ierr); 59c4762a1bSJed Brown ierr = PetscOptionsGetBool(NULL,NULL,"-cholesky",&chol,NULL);CHKERRQ(ierr); 60c4762a1bSJed Brown #if defined(PETSC_HAVE_MUMPS) 61c4762a1bSJed Brown ierr = PetscOptionsGetBool(NULL,NULL,"-test_mumps_opts",&test_mumps_opts,NULL);CHKERRQ(ierr); 62c4762a1bSJed Brown #endif 63c4762a1bSJed Brown 64c4762a1bSJed Brown ierr = PetscRandomCreate(PETSC_COMM_WORLD,&rand);CHKERRQ(ierr); 65c4762a1bSJed Brown ierr = PetscRandomSetFromOptions(rand);CHKERRQ(ierr); 66c4762a1bSJed Brown ierr = MatSetRandom(C,rand);CHKERRQ(ierr); 67c4762a1bSJed Brown ierr = MatDuplicate(C,MAT_DO_NOT_COPY_VALUES,&X);CHKERRQ(ierr); 68c4762a1bSJed Brown 69c4762a1bSJed Brown /* Create vectors */ 7038a8e8c1SStefano Zampini ierr = MatCreateVecs(A,&x,&b);CHKERRQ(ierr); 71c4762a1bSJed Brown ierr = VecDuplicate(x,&u);CHKERRQ(ierr); /* save the true solution */ 72c4762a1bSJed Brown 73c4762a1bSJed Brown /* Test Factorization */ 74c4762a1bSJed Brown ierr = MatGetOrdering(A,MATORDERINGND,&perm,&iperm);CHKERRQ(ierr); 75c4762a1bSJed Brown 76c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-mat_solver_type",&ipack,NULL);CHKERRQ(ierr); 77c4762a1bSJed Brown switch (ipack) { 78c4762a1bSJed Brown #if defined(PETSC_HAVE_SUPERLU) 79c4762a1bSJed Brown case 0: 80c4762a1bSJed Brown if (chol) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"SuperLU does not provide Cholesky!"); 81c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD," SUPERLU LU:\n");CHKERRQ(ierr); 82c4762a1bSJed Brown ierr = MatGetFactor(A,MATSOLVERSUPERLU,MAT_FACTOR_LU,&F);CHKERRQ(ierr); 83c4762a1bSJed Brown matsolvexx = PETSC_TRUE; 84c4762a1bSJed Brown break; 85c4762a1bSJed Brown #endif 86c4762a1bSJed Brown #if defined(PETSC_HAVE_SUPERLU_DIST) 87c4762a1bSJed Brown case 1: 88c4762a1bSJed Brown if (chol) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"SuperLU does not provide Cholesky!"); 89c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD," SUPERLU_DIST LU:\n");CHKERRQ(ierr); 90c4762a1bSJed Brown ierr = MatGetFactor(A,MATSOLVERSUPERLU_DIST,MAT_FACTOR_LU,&F);CHKERRQ(ierr); 91c4762a1bSJed Brown matsolvexx = PETSC_TRUE; 92c4762a1bSJed Brown break; 93c4762a1bSJed Brown #endif 94c4762a1bSJed Brown #if defined(PETSC_HAVE_MUMPS) 95c4762a1bSJed Brown case 2: 96c4762a1bSJed Brown if (chol) { 97c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD," MUMPS CHOLESKY:\n");CHKERRQ(ierr); 98c4762a1bSJed Brown ierr = MatGetFactor(A,MATSOLVERMUMPS,MAT_FACTOR_CHOLESKY,&F);CHKERRQ(ierr); 99c4762a1bSJed Brown } else { 100c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD," MUMPS LU:\n");CHKERRQ(ierr); 101c4762a1bSJed Brown ierr = MatGetFactor(A,MATSOLVERMUMPS,MAT_FACTOR_LU,&F);CHKERRQ(ierr); 102c4762a1bSJed Brown } 103c4762a1bSJed Brown matsolvexx = PETSC_TRUE; 104c4762a1bSJed Brown if (test_mumps_opts) { 105c4762a1bSJed Brown /* test mumps options */ 106c4762a1bSJed Brown PetscInt icntl; 107c4762a1bSJed Brown PetscReal cntl; 108c4762a1bSJed Brown 109c4762a1bSJed Brown icntl = 2; /* sequential matrix ordering */ 110c4762a1bSJed Brown ierr = MatMumpsSetIcntl(F,7,icntl);CHKERRQ(ierr); 111c4762a1bSJed Brown 112c4762a1bSJed Brown cntl = 1.e-6; /* threshold for row pivot detection */ 113c4762a1bSJed Brown ierr = MatMumpsSetIcntl(F,24,1);CHKERRQ(ierr); 114c4762a1bSJed Brown ierr = MatMumpsSetCntl(F,3,cntl);CHKERRQ(ierr); 115c4762a1bSJed Brown } 116c4762a1bSJed Brown break; 117c4762a1bSJed Brown #endif 118c4762a1bSJed Brown #if defined(PETSC_HAVE_MKL_PARDISO) 119c4762a1bSJed Brown case 3: 120c4762a1bSJed Brown if (chol) { 121c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD," MKL_PARDISO CHOLESKY:\n");CHKERRQ(ierr); 122c4762a1bSJed Brown ierr = MatGetFactor(A,MATSOLVERMKL_PARDISO,MAT_FACTOR_CHOLESKY,&F);CHKERRQ(ierr); 123c4762a1bSJed Brown } else { 124c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD," MKL_PARDISO LU:\n");CHKERRQ(ierr); 125c4762a1bSJed Brown ierr = MatGetFactor(A,MATSOLVERMKL_PARDISO,MAT_FACTOR_LU,&F);CHKERRQ(ierr); 126c4762a1bSJed Brown } 127c4762a1bSJed Brown break; 128c4762a1bSJed Brown #endif 12938a8e8c1SStefano Zampini #if defined(PETSC_HAVE_CUDA) 13038a8e8c1SStefano Zampini case 4: 13138a8e8c1SStefano Zampini if (chol) { 13238a8e8c1SStefano Zampini ierr = PetscPrintf(PETSC_COMM_WORLD," CUSPARSE CHOLESKY:\n");CHKERRQ(ierr); 13338a8e8c1SStefano Zampini ierr = MatGetFactor(A,MATSOLVERCUSPARSE,MAT_FACTOR_CHOLESKY,&F);CHKERRQ(ierr); 13438a8e8c1SStefano Zampini } else { 13538a8e8c1SStefano Zampini ierr = PetscPrintf(PETSC_COMM_WORLD," CUSPARSE LU:\n");CHKERRQ(ierr); 13638a8e8c1SStefano Zampini ierr = MatGetFactor(A,MATSOLVERCUSPARSE,MAT_FACTOR_LU,&F);CHKERRQ(ierr); 13738a8e8c1SStefano Zampini } 13838a8e8c1SStefano Zampini break; 13938a8e8c1SStefano Zampini #endif 140c4762a1bSJed Brown default: 141c4762a1bSJed Brown if (chol) { 142c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD," PETSC CHOLESKY:\n");CHKERRQ(ierr); 143c4762a1bSJed Brown ierr = MatGetFactor(A,MATSOLVERPETSC,MAT_FACTOR_CHOLESKY,&F);CHKERRQ(ierr); 144c4762a1bSJed Brown } else { 145c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD," PETSC LU:\n");CHKERRQ(ierr); 146c4762a1bSJed Brown ierr = MatGetFactor(A,MATSOLVERPETSC,MAT_FACTOR_LU,&F);CHKERRQ(ierr); 147c4762a1bSJed Brown } 148c4762a1bSJed Brown matsolvexx = PETSC_TRUE; 149c4762a1bSJed Brown } 150c4762a1bSJed Brown 151c4762a1bSJed Brown ierr = MatFactorInfoInitialize(&info);CHKERRQ(ierr); 152c4762a1bSJed Brown info.fill = 5.0; 153c4762a1bSJed Brown info.shifttype = (PetscReal) MAT_SHIFT_NONE; 154c4762a1bSJed Brown if (chol) { 155c4762a1bSJed Brown ierr = MatCholeskyFactorSymbolic(F,A,perm,&info);CHKERRQ(ierr); 156c4762a1bSJed Brown } else { 157c4762a1bSJed Brown ierr = MatLUFactorSymbolic(F,A,perm,iperm,&info);CHKERRQ(ierr); 158c4762a1bSJed Brown } 159c4762a1bSJed Brown 160c4762a1bSJed Brown for (nfact = 0; nfact < 2; nfact++) { 161c4762a1bSJed Brown if (chol) { 162c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD," %D-the CHOLESKY numfactorization \n",nfact);CHKERRQ(ierr); 163c4762a1bSJed Brown ierr = MatCholeskyFactorNumeric(F,A,&info);CHKERRQ(ierr); 164c4762a1bSJed Brown } else { 165c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD," %D-the LU numfactorization \n",nfact);CHKERRQ(ierr); 166c4762a1bSJed Brown ierr = MatLUFactorNumeric(F,A,&info);CHKERRQ(ierr); 167c4762a1bSJed Brown } 168c4762a1bSJed Brown if (view) { 169c4762a1bSJed Brown ierr = PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr); 170c4762a1bSJed Brown ierr = MatView(F,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 171c4762a1bSJed Brown ierr = PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 172c4762a1bSJed Brown view = PETSC_FALSE; 173c4762a1bSJed Brown } 174c4762a1bSJed Brown 175c4762a1bSJed Brown #if defined(PETSC_HAVE_SUPERLU_DIST) 176c4762a1bSJed Brown if (ipack == 1) { /* Test MatSuperluDistGetDiagU() 177c4762a1bSJed Brown -- input: matrix factor F; output: main diagonal of matrix U on all processes */ 178c4762a1bSJed Brown PetscInt M; 179c4762a1bSJed Brown PetscScalar *diag; 180c4762a1bSJed Brown #if !defined(PETSC_USE_COMPLEX) 181c4762a1bSJed Brown PetscInt nneg,nzero,npos; 182c4762a1bSJed Brown #endif 183c4762a1bSJed Brown 184c4762a1bSJed Brown ierr = MatGetSize(F,&M,NULL);CHKERRQ(ierr); 185c4762a1bSJed Brown ierr = PetscMalloc1(M,&diag);CHKERRQ(ierr); 186c4762a1bSJed Brown ierr = MatSuperluDistGetDiagU(F,diag);CHKERRQ(ierr); 187c4762a1bSJed Brown ierr = PetscFree(diag);CHKERRQ(ierr); 188c4762a1bSJed Brown 189c4762a1bSJed Brown #if !defined(PETSC_USE_COMPLEX) 190c4762a1bSJed Brown /* Test MatGetInertia() */ 191c4762a1bSJed Brown ierr = MatGetInertia(F,&nneg,&nzero,&npos);CHKERRQ(ierr); 192c4762a1bSJed Brown ierr = PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD," MatInertia: nneg: %D, nzero: %D, npos: %D\n",nneg,nzero,npos);CHKERRQ(ierr); 193c4762a1bSJed Brown #endif 194c4762a1bSJed Brown } 195c4762a1bSJed Brown #endif 196c4762a1bSJed Brown 197c4762a1bSJed Brown /* Test MatMatSolve() */ 198c4762a1bSJed Brown if (testMatMatSolve) { 199c4762a1bSJed Brown if (!nfact) { 200c4762a1bSJed Brown ierr = MatMatMult(A,C,MAT_INITIAL_MATRIX,2.0,&RHS);CHKERRQ(ierr); 201c4762a1bSJed Brown } else { 202c4762a1bSJed Brown ierr = MatMatMult(A,C,MAT_REUSE_MATRIX,2.0,&RHS);CHKERRQ(ierr); 203c4762a1bSJed Brown } 204c4762a1bSJed Brown for (nsolve = 0; nsolve < 2; nsolve++) { 205c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD," %D-the MatMatSolve \n",nsolve);CHKERRQ(ierr); 206c4762a1bSJed Brown ierr = MatMatSolve(F,RHS,X);CHKERRQ(ierr); 207c4762a1bSJed Brown 208c4762a1bSJed Brown /* Check the error */ 209c4762a1bSJed Brown ierr = MatAXPY(X,-1.0,C,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 210c4762a1bSJed Brown ierr = MatNorm(X,NORM_FROBENIUS,&norm);CHKERRQ(ierr); 211c4762a1bSJed Brown if (norm > tol) { 212c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"%D-the MatMatSolve: Norm of error %g, nsolve %D\n",nsolve,(double)norm,nsolve);CHKERRQ(ierr); 213c4762a1bSJed Brown } 214c4762a1bSJed Brown } 215c4762a1bSJed Brown if (matsolvexx) { 216c4762a1bSJed Brown /* Test MatMatSolve(F,RHS,RHS), RHS is a dense matrix */ 217c4762a1bSJed Brown ierr = MatCopy(RHS,X,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 218c4762a1bSJed Brown ierr = MatMatSolve(F,X,X);CHKERRQ(ierr); 219c4762a1bSJed Brown /* Check the error */ 220c4762a1bSJed Brown ierr = MatAXPY(X,-1.0,C,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 221c4762a1bSJed Brown ierr = MatNorm(X,NORM_FROBENIUS,&norm);CHKERRQ(ierr); 222c4762a1bSJed Brown if (norm > tol) { 223c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"MatMatSolve(F,RHS,RHS): Norm of error %g\n",(double)norm);CHKERRQ(ierr); 224c4762a1bSJed Brown } 225c4762a1bSJed Brown } 226c4762a1bSJed Brown 227c4762a1bSJed Brown if (ipack == 2 && size == 1) { 228c4762a1bSJed Brown Mat spRHS,spRHST,RHST; 229c4762a1bSJed Brown 230c4762a1bSJed Brown ierr = MatTranspose(RHS,MAT_INITIAL_MATRIX,&RHST);CHKERRQ(ierr); 231c4762a1bSJed Brown ierr = MatConvert(RHST,MATAIJ,MAT_INITIAL_MATRIX,&spRHST);CHKERRQ(ierr); 232c4762a1bSJed Brown ierr = MatCreateTranspose(spRHST,&spRHS);CHKERRQ(ierr); 233c4762a1bSJed Brown for (nsolve = 0; nsolve < 2; nsolve++) { 234c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD," %D-the sparse MatMatSolve \n",nsolve);CHKERRQ(ierr); 235c4762a1bSJed Brown ierr = MatMatSolve(F,spRHS,X);CHKERRQ(ierr); 236c4762a1bSJed Brown 237c4762a1bSJed Brown /* Check the error */ 238c4762a1bSJed Brown ierr = MatAXPY(X,-1.0,C,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 239c4762a1bSJed Brown ierr = MatNorm(X,NORM_FROBENIUS,&norm);CHKERRQ(ierr); 240c4762a1bSJed Brown if (norm > tol) { 241c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"%D-the sparse MatMatSolve: Norm of error %g, nsolve %D\n",nsolve,(double)norm,nsolve);CHKERRQ(ierr); 242c4762a1bSJed Brown } 243c4762a1bSJed Brown } 244c4762a1bSJed Brown ierr = MatDestroy(&spRHST);CHKERRQ(ierr); 245c4762a1bSJed Brown ierr = MatDestroy(&spRHS);CHKERRQ(ierr); 246c4762a1bSJed Brown ierr = MatDestroy(&RHST);CHKERRQ(ierr); 247c4762a1bSJed Brown } 248c4762a1bSJed Brown } 249c4762a1bSJed Brown 250c4762a1bSJed Brown /* Test MatSolve() */ 251c4762a1bSJed Brown if (testMatSolve) { 252c4762a1bSJed Brown for (nsolve = 0; nsolve < 2; nsolve++) { 253c4762a1bSJed Brown ierr = VecSetRandom(x,rand);CHKERRQ(ierr); 254c4762a1bSJed Brown ierr = VecCopy(x,u);CHKERRQ(ierr); 255c4762a1bSJed Brown ierr = MatMult(A,x,b);CHKERRQ(ierr); 256c4762a1bSJed Brown 257c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD," %D-the MatSolve \n",nsolve);CHKERRQ(ierr); 258c4762a1bSJed Brown ierr = MatSolve(F,b,x);CHKERRQ(ierr); 259c4762a1bSJed Brown 260c4762a1bSJed Brown /* Check the error */ 261c4762a1bSJed Brown ierr = VecAXPY(u,-1.0,x);CHKERRQ(ierr); /* u <- (-1.0)x + u */ 262c4762a1bSJed Brown ierr = VecNorm(u,NORM_2,&norm);CHKERRQ(ierr); 263c4762a1bSJed Brown if (norm > tol) { 264c4762a1bSJed Brown PetscReal resi; 265c4762a1bSJed Brown ierr = MatMult(A,x,u);CHKERRQ(ierr); /* u = A*x */ 266c4762a1bSJed Brown ierr = VecAXPY(u,-1.0,b);CHKERRQ(ierr); /* u <- (-1.0)b + u */ 267c4762a1bSJed Brown ierr = VecNorm(u,NORM_2,&resi);CHKERRQ(ierr); 268c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"MatSolve: Norm of error %g, resi %g, numfact %D\n",(double)norm,(double)resi,nfact);CHKERRQ(ierr); 269c4762a1bSJed Brown } 270c4762a1bSJed Brown } 271c4762a1bSJed Brown } 272c4762a1bSJed Brown } 273c4762a1bSJed Brown 274c4762a1bSJed Brown /* Free data structures */ 275c4762a1bSJed Brown ierr = MatDestroy(&A);CHKERRQ(ierr); 276c4762a1bSJed Brown ierr = MatDestroy(&C);CHKERRQ(ierr); 277c4762a1bSJed Brown ierr = MatDestroy(&F);CHKERRQ(ierr); 278c4762a1bSJed Brown ierr = MatDestroy(&X);CHKERRQ(ierr); 279c4762a1bSJed Brown if (testMatMatSolve) { 280c4762a1bSJed Brown ierr = MatDestroy(&RHS);CHKERRQ(ierr); 281c4762a1bSJed Brown } 282c4762a1bSJed Brown 283c4762a1bSJed Brown ierr = PetscRandomDestroy(&rand);CHKERRQ(ierr); 284c4762a1bSJed Brown ierr = ISDestroy(&perm);CHKERRQ(ierr); 285c4762a1bSJed Brown ierr = ISDestroy(&iperm);CHKERRQ(ierr); 286c4762a1bSJed Brown ierr = VecDestroy(&x);CHKERRQ(ierr); 287c4762a1bSJed Brown ierr = VecDestroy(&b);CHKERRQ(ierr); 288c4762a1bSJed Brown ierr = VecDestroy(&u);CHKERRQ(ierr); 289c4762a1bSJed Brown ierr = PetscFinalize(); 290c4762a1bSJed Brown return ierr; 291c4762a1bSJed Brown } 292c4762a1bSJed Brown 293c4762a1bSJed Brown /*TEST 294c4762a1bSJed Brown 295c4762a1bSJed Brown test: 296*dfd57a17SPierre Jolivet requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES) 297c4762a1bSJed Brown args: -f ${DATAFILESPATH}/matrices/small -mat_solver_type 10 298c4762a1bSJed Brown output_file: output/ex125.out 299c4762a1bSJed Brown 300c4762a1bSJed Brown test: 301c4762a1bSJed Brown suffix: mkl_pardiso 302*dfd57a17SPierre Jolivet requires: mkl_pardiso datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES) 303c4762a1bSJed Brown args: -f ${DATAFILESPATH}/matrices/small -mat_solver_type 3 304c4762a1bSJed Brown 305c4762a1bSJed Brown test: 306c4762a1bSJed Brown suffix: mumps 307*dfd57a17SPierre Jolivet requires: mumps datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES) 308c4762a1bSJed Brown args: -f ${DATAFILESPATH}/matrices/small -mat_solver_type 2 309c4762a1bSJed Brown output_file: output/ex125_mumps_seq.out 310c4762a1bSJed Brown 311c4762a1bSJed Brown test: 312c4762a1bSJed Brown suffix: mumps_2 313c4762a1bSJed Brown nsize: 3 314*dfd57a17SPierre Jolivet requires: mumps datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES) 315c4762a1bSJed Brown args: -f ${DATAFILESPATH}/matrices/small -mat_solver_type 2 316c4762a1bSJed Brown output_file: output/ex125_mumps_par.out 317c4762a1bSJed Brown 318c4762a1bSJed Brown test: 319c4762a1bSJed Brown suffix: superlu_dist 320*dfd57a17SPierre Jolivet requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES) superlu_dist 321c4762a1bSJed Brown args: -f ${DATAFILESPATH}/matrices/small -mat_solver_type 1 -mat_superlu_dist_rowperm NOROWPERM 322c4762a1bSJed Brown 323c4762a1bSJed Brown test: 324c4762a1bSJed Brown suffix: superlu_dist_2 325c4762a1bSJed Brown nsize: 3 326*dfd57a17SPierre Jolivet requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES) superlu_dist 327c4762a1bSJed Brown args: -f ${DATAFILESPATH}/matrices/small -mat_solver_type 1 -mat_superlu_dist_rowperm NOROWPERM 328c4762a1bSJed Brown output_file: output/ex125_superlu_dist.out 329c4762a1bSJed Brown 330c4762a1bSJed Brown test: 331c4762a1bSJed Brown suffix: superlu_dist_complex 332c4762a1bSJed Brown nsize: 3 333*dfd57a17SPierre Jolivet requires: datafilespath superlu_dist complex double !defined(PETSC_USE_64BIT_INDICES) 334c4762a1bSJed Brown args: -f ${DATAFILESPATH}/matrices/farzad_B_rhs -mat_solver_type 1 335c4762a1bSJed Brown output_file: output/ex125_superlu_dist_complex.out 336c4762a1bSJed Brown 33738a8e8c1SStefano Zampini test: 33838a8e8c1SStefano Zampini suffix: cusparse 339*dfd57a17SPierre Jolivet requires: cuda datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES) 34038a8e8c1SStefano Zampini args: -mat_type aijcusparse -f ${DATAFILESPATH}/matrices/small -mat_solver_type 4 -cholesky {{0 1}separate output} 34138a8e8c1SStefano Zampini 342c4762a1bSJed Brown TEST*/ 343