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