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 PetscMPIInt size; 13c4762a1bSJed Brown PetscInt isolver=0,size_schur,m,n,nfact,nsolve,nrhs; 14c4762a1bSJed Brown PetscReal norm,tol=PETSC_SQRT_MACHINE_EPSILON; 15c4762a1bSJed Brown PetscRandom rand; 16c4762a1bSJed Brown PetscBool data_provided,herm,symm,use_lu,cuda = PETSC_FALSE; 17c4762a1bSJed Brown PetscReal sratio = 5.1/12.; 18c4762a1bSJed Brown PetscViewer fd; /* viewer */ 19c4762a1bSJed Brown char solver[256]; 20c4762a1bSJed Brown char file[PETSC_MAX_PATH_LEN]; /* input file name */ 21c4762a1bSJed Brown 22*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc,&args,(char*)0,help)); 235f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 242c71b3e2SJacob Faibussowitsch PetscCheckFalse(size > 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor test"); 25c4762a1bSJed Brown /* Determine which type of solver we want to test for */ 26c4762a1bSJed Brown herm = PETSC_FALSE; 27c4762a1bSJed Brown symm = PETSC_FALSE; 285f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-symmetric_solve",&symm,NULL)); 295f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-hermitian_solve",&herm,NULL)); 30c4762a1bSJed Brown if (herm) symm = PETSC_TRUE; 315f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-cuda_solve",&cuda,NULL)); 325f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-tol",&tol,NULL)); 33c4762a1bSJed Brown 34c4762a1bSJed Brown /* Determine file from which we read the matrix A */ 355f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetString(NULL,NULL,"-f",file,sizeof(file),&data_provided)); 36c4762a1bSJed Brown if (!data_provided) { /* get matrices from PETSc distribution */ 375f80ce2aSJacob Faibussowitsch CHKERRQ(PetscStrncpy(file,"${PETSC_DIR}/share/petsc/datafiles/matrices/",sizeof(file))); 38c4762a1bSJed Brown if (symm) { 39c4762a1bSJed Brown #if defined (PETSC_USE_COMPLEX) 405f80ce2aSJacob Faibussowitsch CHKERRQ(PetscStrlcat(file,"hpd-complex-",sizeof(file))); 41c4762a1bSJed Brown #else 425f80ce2aSJacob Faibussowitsch CHKERRQ(PetscStrlcat(file,"spd-real-",sizeof(file))); 43c4762a1bSJed Brown #endif 44c4762a1bSJed Brown } else { 45c4762a1bSJed Brown #if defined (PETSC_USE_COMPLEX) 465f80ce2aSJacob Faibussowitsch CHKERRQ(PetscStrlcat(file,"nh-complex-",sizeof(file))); 47c4762a1bSJed Brown #else 485f80ce2aSJacob Faibussowitsch CHKERRQ(PetscStrlcat(file,"ns-real-",sizeof(file))); 49c4762a1bSJed Brown #endif 50c4762a1bSJed Brown } 51c4762a1bSJed Brown #if defined(PETSC_USE_64BIT_INDICES) 525f80ce2aSJacob Faibussowitsch CHKERRQ(PetscStrlcat(file,"int64-",sizeof(file))); 53c4762a1bSJed Brown #else 545f80ce2aSJacob Faibussowitsch CHKERRQ(PetscStrlcat(file,"int32-",sizeof(file))); 55c4762a1bSJed Brown #endif 56c4762a1bSJed Brown #if defined (PETSC_USE_REAL_SINGLE) 575f80ce2aSJacob Faibussowitsch CHKERRQ(PetscStrlcat(file,"float32",sizeof(file))); 58c4762a1bSJed Brown #else 595f80ce2aSJacob Faibussowitsch CHKERRQ(PetscStrlcat(file,"float64",sizeof(file))); 60c4762a1bSJed Brown #endif 61c4762a1bSJed Brown } 62c4762a1bSJed Brown /* Load matrix A */ 635f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&fd)); 645f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A)); 655f80ce2aSJacob Faibussowitsch CHKERRQ(MatLoad(A,fd)); 665f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDestroy(&fd)); 675f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetSize(A,&m,&n)); 682c71b3e2SJacob Faibussowitsch PetscCheckFalse(m != n,PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ, "This example is not intended for rectangular matrices (%" PetscInt_FMT ", %" PetscInt_FMT ")", m, n); 69c4762a1bSJed Brown 70a5b23f4aSJose E. Roman /* Create dense matrix C and X; C holds true solution with identical columns */ 71c4762a1bSJed Brown nrhs = 2; 725f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-nrhs",&nrhs,NULL)); 735f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD,&C)); 745f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(C,m,PETSC_DECIDE,PETSC_DECIDE,nrhs)); 755f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetType(C,MATDENSE)); 765f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(C)); 775f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(C)); 78c4762a1bSJed Brown 795f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomCreate(PETSC_COMM_WORLD,&rand)); 805f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomSetFromOptions(rand)); 815f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetRandom(C,rand)); 825f80ce2aSJacob Faibussowitsch CHKERRQ(MatDuplicate(C,MAT_DO_NOT_COPY_VALUES,&X)); 83c4762a1bSJed Brown 84c4762a1bSJed Brown /* Create vectors */ 855f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreate(PETSC_COMM_WORLD,&x)); 865f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetSizes(x,n,PETSC_DECIDE)); 875f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetFromOptions(x)); 885f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(x,&b)); 895f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(x,&u)); /* save the true solution */ 90c4762a1bSJed Brown 915f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-solver",&isolver,NULL)); 92c4762a1bSJed Brown switch (isolver) { 93c4762a1bSJed Brown #if defined(PETSC_HAVE_MUMPS) 94c4762a1bSJed Brown case 0: 955f80ce2aSJacob Faibussowitsch CHKERRQ(PetscStrcpy(solver,MATSOLVERMUMPS)); 96c4762a1bSJed Brown break; 97c4762a1bSJed Brown #endif 98c4762a1bSJed Brown #if defined(PETSC_HAVE_MKL_PARDISO) 99c4762a1bSJed Brown case 1: 1005f80ce2aSJacob Faibussowitsch CHKERRQ(PetscStrcpy(solver,MATSOLVERMKL_PARDISO)); 101c4762a1bSJed Brown break; 102c4762a1bSJed Brown #endif 103c4762a1bSJed Brown default: 1045f80ce2aSJacob Faibussowitsch CHKERRQ(PetscStrcpy(solver,MATSOLVERPETSC)); 105c4762a1bSJed Brown break; 106c4762a1bSJed Brown } 107c4762a1bSJed Brown 108c4762a1bSJed Brown #if defined (PETSC_USE_COMPLEX) 109c4762a1bSJed Brown if (isolver == 0 && symm && !data_provided) { /* MUMPS (5.0.0) does not have support for hermitian matrices, so make them symmetric */ 110c4762a1bSJed Brown PetscScalar im = PetscSqrtScalar((PetscScalar)-1.); 111c4762a1bSJed Brown PetscScalar val = -1.0; 112c4762a1bSJed Brown val = val + im; 1135f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValue(A,1,0,val,INSERT_VALUES)); 1145f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 1155f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 116c4762a1bSJed Brown } 117c4762a1bSJed Brown #endif 118c4762a1bSJed Brown 1195f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-schur_ratio",&sratio,NULL)); 1202c71b3e2SJacob Faibussowitsch PetscCheckFalse(sratio < 0. || sratio > 1.,PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ, "Invalid ratio for schur degrees of freedom %g", (double)sratio); 121c4762a1bSJed Brown size_schur = (PetscInt)(sratio*m); 122c4762a1bSJed Brown 1235f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"Solving with %s: nrhs %" PetscInt_FMT ", sym %d, herm %d, size schur %" PetscInt_FMT ", size mat %" PetscInt_FMT "\n",solver,nrhs,symm,herm,size_schur,m)); 124c4762a1bSJed Brown 125c4762a1bSJed Brown /* Test LU/Cholesky Factorization */ 126c4762a1bSJed Brown use_lu = PETSC_FALSE; 127c4762a1bSJed Brown if (!symm) use_lu = PETSC_TRUE; 128c4762a1bSJed Brown #if defined (PETSC_USE_COMPLEX) 129c4762a1bSJed Brown if (isolver == 1) use_lu = PETSC_TRUE; 130c4762a1bSJed Brown #endif 131c4762a1bSJed Brown if (cuda && symm && !herm) use_lu = PETSC_TRUE; 132c4762a1bSJed Brown 133c4762a1bSJed Brown if (herm && !use_lu) { /* test also conversion routines inside the solver packages */ 1345f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetOption(A,MAT_SYMMETRIC,PETSC_TRUE)); 1355f80ce2aSJacob Faibussowitsch CHKERRQ(MatConvert(A,MATSEQSBAIJ,MAT_INPLACE_MATRIX,&A)); 136c4762a1bSJed Brown } 137c4762a1bSJed Brown 138c4762a1bSJed Brown if (use_lu) { 1395f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetFactor(A,solver,MAT_FACTOR_LU,&F)); 140c4762a1bSJed Brown } else { 141c4762a1bSJed Brown if (herm) { 1425f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetOption(A,MAT_SYMMETRIC,PETSC_TRUE)); 1435f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetOption(A,MAT_SPD,PETSC_TRUE)); 144c4762a1bSJed Brown } else { 1455f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetOption(A,MAT_SYMMETRIC,PETSC_TRUE)); 1465f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetOption(A,MAT_SPD,PETSC_FALSE)); 147c4762a1bSJed Brown } 1485f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetFactor(A,solver,MAT_FACTOR_CHOLESKY,&F)); 149c4762a1bSJed Brown } 1505f80ce2aSJacob Faibussowitsch CHKERRQ(ISCreateStride(PETSC_COMM_SELF,size_schur,m-size_schur,1,&is_schur)); 1515f80ce2aSJacob Faibussowitsch CHKERRQ(MatFactorSetSchurIS(F,is_schur)); 152c4762a1bSJed Brown 1535f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&is_schur)); 154c4762a1bSJed Brown if (use_lu) { 1555f80ce2aSJacob Faibussowitsch CHKERRQ(MatLUFactorSymbolic(F,A,NULL,NULL,NULL)); 156c4762a1bSJed Brown } else { 1575f80ce2aSJacob Faibussowitsch CHKERRQ(MatCholeskyFactorSymbolic(F,A,NULL,NULL)); 158c4762a1bSJed Brown } 159c4762a1bSJed Brown 160c4762a1bSJed Brown for (nfact = 0; nfact < 3; nfact++) { 161c4762a1bSJed Brown Mat AD; 162c4762a1bSJed Brown 163c4762a1bSJed Brown if (!nfact) { 1645f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetRandom(x,rand)); 165c4762a1bSJed Brown if (symm && herm) { 1665f80ce2aSJacob Faibussowitsch CHKERRQ(VecAbs(x)); 167c4762a1bSJed Brown } 1685f80ce2aSJacob Faibussowitsch CHKERRQ(MatDiagonalSet(A,x,ADD_VALUES)); 169c4762a1bSJed Brown } 170c4762a1bSJed Brown if (use_lu) { 1715f80ce2aSJacob Faibussowitsch CHKERRQ(MatLUFactorNumeric(F,A,NULL)); 172c4762a1bSJed Brown } else { 1735f80ce2aSJacob Faibussowitsch CHKERRQ(MatCholeskyFactorNumeric(F,A,NULL)); 174c4762a1bSJed Brown } 175c4762a1bSJed Brown if (cuda) { 1765f80ce2aSJacob Faibussowitsch CHKERRQ(MatFactorGetSchurComplement(F,&S,NULL)); 1775f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetType(S,MATSEQDENSECUDA)); 1785f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateVecs(S,&xschur,&bschur)); 1795f80ce2aSJacob Faibussowitsch CHKERRQ(MatFactorRestoreSchurComplement(F,&S,MAT_FACTOR_SCHUR_UNFACTORED)); 180c4762a1bSJed Brown } 1815f80ce2aSJacob Faibussowitsch CHKERRQ(MatFactorCreateSchurComplement(F,&S,NULL)); 182c4762a1bSJed Brown if (!cuda) { 1835f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateVecs(S,&xschur,&bschur)); 184c4762a1bSJed Brown } 1855f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(xschur,&uschur)); 186c4762a1bSJed Brown if (nfact == 1 && (!cuda || (herm && symm))) { 1875f80ce2aSJacob Faibussowitsch CHKERRQ(MatFactorInvertSchurComplement(F)); 188c4762a1bSJed Brown } 189c4762a1bSJed Brown for (nsolve = 0; nsolve < 2; nsolve++) { 1905f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetRandom(x,rand)); 1915f80ce2aSJacob Faibussowitsch CHKERRQ(VecCopy(x,u)); 192c4762a1bSJed Brown 193c4762a1bSJed Brown if (nsolve) { 1945f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(A,x,b)); 1955f80ce2aSJacob Faibussowitsch CHKERRQ(MatSolve(F,b,x)); 196c4762a1bSJed Brown } else { 1975f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultTranspose(A,x,b)); 1985f80ce2aSJacob Faibussowitsch CHKERRQ(MatSolveTranspose(F,b,x)); 199c4762a1bSJed Brown } 200c4762a1bSJed Brown /* Check the error */ 2015f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(u,-1.0,x)); /* u <- (-1.0)x + u */ 2025f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(u,NORM_2,&norm)); 203c4762a1bSJed Brown if (norm > tol) { 204c4762a1bSJed Brown PetscReal resi; 205c4762a1bSJed Brown if (nsolve) { 2065f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(A,x,u)); /* u = A*x */ 207c4762a1bSJed Brown } else { 2085f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultTranspose(A,x,u)); /* u = A*x */ 209c4762a1bSJed Brown } 2105f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(u,-1.0,b)); /* u <- (-1.0)b + u */ 2115f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(u,NORM_2,&resi)); 212c4762a1bSJed Brown if (nsolve) { 2135f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"(f %" PetscInt_FMT ", s %" PetscInt_FMT ") MatSolve error: Norm of error %g, residual %g\n",nfact,nsolve,(double)norm,(double)resi)); 214c4762a1bSJed Brown } else { 2155f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"(f %" PetscInt_FMT ", s %" PetscInt_FMT ") MatSolveTranspose error: Norm of error %g, residual %f\n",nfact,nsolve,(double)norm,(double)resi)); 216c4762a1bSJed Brown } 217c4762a1bSJed Brown } 2185f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetRandom(xschur,rand)); 2195f80ce2aSJacob Faibussowitsch CHKERRQ(VecCopy(xschur,uschur)); 220c4762a1bSJed Brown if (nsolve) { 2215f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(S,xschur,bschur)); 2225f80ce2aSJacob Faibussowitsch CHKERRQ(MatFactorSolveSchurComplement(F,bschur,xschur)); 223c4762a1bSJed Brown } else { 2245f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultTranspose(S,xschur,bschur)); 2255f80ce2aSJacob Faibussowitsch CHKERRQ(MatFactorSolveSchurComplementTranspose(F,bschur,xschur)); 226c4762a1bSJed Brown } 227c4762a1bSJed Brown /* Check the error */ 2285f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(uschur,-1.0,xschur)); /* u <- (-1.0)x + u */ 2295f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(uschur,NORM_2,&norm)); 230c4762a1bSJed Brown if (norm > tol) { 231c4762a1bSJed Brown PetscReal resi; 232c4762a1bSJed Brown if (nsolve) { 2335f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(S,xschur,uschur)); /* u = A*x */ 234c4762a1bSJed Brown } else { 2355f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultTranspose(S,xschur,uschur)); /* u = A*x */ 236c4762a1bSJed Brown } 2375f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(uschur,-1.0,bschur)); /* u <- (-1.0)b + u */ 2385f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(uschur,NORM_2,&resi)); 239c4762a1bSJed Brown if (nsolve) { 2405f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"(f %" PetscInt_FMT ", s %" PetscInt_FMT ") MatFactorSolveSchurComplement error: Norm of error %g, residual %g\n",nfact,nsolve,(double)norm,(double)resi)); 241c4762a1bSJed Brown } else { 2425f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"(f %" PetscInt_FMT ", s %" PetscInt_FMT ") MatFactorSolveSchurComplementTranspose error: Norm of error %g, residual %f\n",nfact,nsolve,(double)norm,(double)resi)); 243c4762a1bSJed Brown } 244c4762a1bSJed Brown } 245c4762a1bSJed Brown } 2465f80ce2aSJacob Faibussowitsch CHKERRQ(MatConvert(A,MATSEQAIJ,MAT_INITIAL_MATRIX,&AD)); 247c4762a1bSJed Brown if (!nfact) { 2485f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatMult(AD,C,MAT_INITIAL_MATRIX,2.0,&RHS)); 249c4762a1bSJed Brown } else { 2505f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatMult(AD,C,MAT_REUSE_MATRIX,2.0,&RHS)); 251c4762a1bSJed Brown } 2525f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&AD)); 253c4762a1bSJed Brown for (nsolve = 0; nsolve < 2; nsolve++) { 2545f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatSolve(F,RHS,X)); 255c4762a1bSJed Brown 256c4762a1bSJed Brown /* Check the error */ 2575f80ce2aSJacob Faibussowitsch CHKERRQ(MatAXPY(X,-1.0,C,SAME_NONZERO_PATTERN)); 2585f80ce2aSJacob Faibussowitsch CHKERRQ(MatNorm(X,NORM_FROBENIUS,&norm)); 259c4762a1bSJed Brown if (norm > tol) { 2605f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"(f %" PetscInt_FMT ", s %" PetscInt_FMT ") MatMatSolve: Norm of error %g\n",nfact,nsolve,(double)norm)); 261c4762a1bSJed Brown } 262c4762a1bSJed Brown } 263c4762a1bSJed Brown if (isolver == 0) { 264c4762a1bSJed Brown Mat spRHS,spRHST,RHST; 265c4762a1bSJed Brown 2665f80ce2aSJacob Faibussowitsch CHKERRQ(MatTranspose(RHS,MAT_INITIAL_MATRIX,&RHST)); 2675f80ce2aSJacob Faibussowitsch CHKERRQ(MatConvert(RHST,MATSEQAIJ,MAT_INITIAL_MATRIX,&spRHST)); 2685f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateTranspose(spRHST,&spRHS)); 269c4762a1bSJed Brown for (nsolve = 0; nsolve < 2; nsolve++) { 2705f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatSolve(F,spRHS,X)); 271c4762a1bSJed Brown 272c4762a1bSJed Brown /* Check the error */ 2735f80ce2aSJacob Faibussowitsch CHKERRQ(MatAXPY(X,-1.0,C,SAME_NONZERO_PATTERN)); 2745f80ce2aSJacob Faibussowitsch CHKERRQ(MatNorm(X,NORM_FROBENIUS,&norm)); 275c4762a1bSJed Brown if (norm > tol) { 2765f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"(f %" PetscInt_FMT ", s %" PetscInt_FMT ") sparse MatMatSolve: Norm of error %g\n",nfact,nsolve,(double)norm)); 277c4762a1bSJed Brown } 278c4762a1bSJed Brown } 2795f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&spRHST)); 2805f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&spRHS)); 2815f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&RHST)); 282c4762a1bSJed Brown } 2835f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&S)); 2845f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&xschur)); 2855f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&bschur)); 2865f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&uschur)); 287c4762a1bSJed Brown } 288c4762a1bSJed Brown /* Free data structures */ 2895f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 2905f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&C)); 2915f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&F)); 2925f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&X)); 2935f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&RHS)); 2945f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomDestroy(&rand)); 2955f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&x)); 2965f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&b)); 2975f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&u)); 298*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 299*b122ec5aSJacob Faibussowitsch return 0; 300c4762a1bSJed Brown } 301c4762a1bSJed Brown 302c4762a1bSJed Brown /*TEST 303c4762a1bSJed Brown 304c4762a1bSJed Brown testset: 30538cf7977SStefano Zampini requires: mkl_pardiso double !complex 306c4762a1bSJed Brown args: -solver 1 307c4762a1bSJed Brown 308c4762a1bSJed Brown test: 309c4762a1bSJed Brown suffix: mkl_pardiso 310c4762a1bSJed Brown test: 311c4762a1bSJed Brown requires: cuda 312c4762a1bSJed Brown suffix: mkl_pardiso_cuda 313c4762a1bSJed Brown args: -cuda_solve 314c4762a1bSJed Brown output_file: output/ex192_mkl_pardiso.out 315c4762a1bSJed Brown test: 316c4762a1bSJed Brown suffix: mkl_pardiso_1 317c4762a1bSJed Brown args: -symmetric_solve 318c4762a1bSJed Brown output_file: output/ex192_mkl_pardiso_1.out 319c4762a1bSJed Brown test: 320c4762a1bSJed Brown requires: cuda 321c4762a1bSJed Brown suffix: mkl_pardiso_cuda_1 322c4762a1bSJed Brown args: -symmetric_solve -cuda_solve 323c4762a1bSJed Brown output_file: output/ex192_mkl_pardiso_1.out 324c4762a1bSJed Brown test: 325c4762a1bSJed Brown suffix: mkl_pardiso_3 326c4762a1bSJed Brown args: -symmetric_solve -hermitian_solve 327c4762a1bSJed Brown output_file: output/ex192_mkl_pardiso_3.out 328c4762a1bSJed Brown test: 329dfd57a17SPierre Jolivet requires: cuda defined(PETSC_HAVE_CUSOLVERDNDPOTRI) 330c4762a1bSJed Brown suffix: mkl_pardiso_cuda_3 331c4762a1bSJed Brown args: -symmetric_solve -hermitian_solve -cuda_solve 332c4762a1bSJed Brown output_file: output/ex192_mkl_pardiso_3.out 333c4762a1bSJed Brown 334c4762a1bSJed Brown testset: 335c4762a1bSJed Brown requires: mumps double !complex 336c4762a1bSJed Brown args: -solver 0 337c4762a1bSJed Brown 338c4762a1bSJed Brown test: 339c4762a1bSJed Brown suffix: mumps 340c4762a1bSJed Brown test: 341c4762a1bSJed Brown requires: cuda 342c4762a1bSJed Brown suffix: mumps_cuda 343c4762a1bSJed Brown args: -cuda_solve 344c4762a1bSJed Brown output_file: output/ex192_mumps.out 345c4762a1bSJed Brown test: 346c4762a1bSJed Brown suffix: mumps_2 347c4762a1bSJed Brown args: -symmetric_solve 348c4762a1bSJed Brown output_file: output/ex192_mumps_2.out 349c4762a1bSJed Brown test: 350c4762a1bSJed Brown requires: cuda 351c4762a1bSJed Brown suffix: mumps_cuda_2 352c4762a1bSJed Brown args: -symmetric_solve -cuda_solve 353c4762a1bSJed Brown output_file: output/ex192_mumps_2.out 354c4762a1bSJed Brown test: 355c4762a1bSJed Brown suffix: mumps_3 356c4762a1bSJed Brown args: -symmetric_solve -hermitian_solve 357c4762a1bSJed Brown output_file: output/ex192_mumps_3.out 358c4762a1bSJed Brown test: 359dfd57a17SPierre Jolivet requires: cuda defined(PETSC_HAVE_CUSOLVERDNDPOTRI) 360c4762a1bSJed Brown suffix: mumps_cuda_3 361c4762a1bSJed Brown args: -symmetric_solve -hermitian_solve -cuda_solve 362c4762a1bSJed Brown output_file: output/ex192_mumps_3.out 363c4762a1bSJed Brown 364c4762a1bSJed Brown TEST*/ 365