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