1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Tests MatMatSolve() and MatMatTransposeSolve() for computing inv(A) with MUMPS.\n\ 3c4762a1bSJed Brown Example: mpiexec -n <np> ./ex214 -displ \n\n"; 4c4762a1bSJed Brown 5c4762a1bSJed Brown #include <petscmat.h> 6c4762a1bSJed Brown 7c4762a1bSJed Brown int main(int argc,char **args) 8c4762a1bSJed Brown { 9c4762a1bSJed Brown PetscMPIInt size,rank; 10c4762a1bSJed Brown #if defined(PETSC_HAVE_MUMPS) 11c4762a1bSJed Brown Mat A,RHS,C,F,X,AX,spRHST; 12c4762a1bSJed Brown PetscInt m,n,nrhs,M,N,i,Istart,Iend,Ii,j,J,test; 13c4762a1bSJed Brown PetscScalar v; 14c4762a1bSJed Brown PetscReal norm,tol=PETSC_SQRT_MACHINE_EPSILON; 15c4762a1bSJed Brown PetscRandom rand; 16c4762a1bSJed Brown PetscBool displ=PETSC_FALSE; 17c4762a1bSJed Brown char solver[256]; 18c4762a1bSJed Brown #endif 19c4762a1bSJed Brown 20*327415f7SBarry Smith PetscFunctionBeginUser; 219566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&args,(char*)0,help)); 229566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 239566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 24c4762a1bSJed Brown 25c4762a1bSJed Brown #if !defined(PETSC_HAVE_MUMPS) 269566063dSJacob Faibussowitsch if (rank == 0) PetscCall(PetscPrintf(PETSC_COMM_SELF,"This example requires MUMPS, exit...\n")); 279566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 28b122ec5aSJacob Faibussowitsch return 0; 29c4762a1bSJed Brown #else 30c4762a1bSJed Brown 319566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL,NULL,"-displ",&displ,NULL)); 32c4762a1bSJed Brown 33c4762a1bSJed Brown /* Create matrix A */ 34c4762a1bSJed Brown m = 4; n = 4; 359566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL)); 369566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL)); 37c4762a1bSJed Brown 389566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD,&A)); 399566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n)); 409566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A)); 419566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(A,5,NULL,5,NULL)); 429566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(A,5,NULL)); 43c4762a1bSJed Brown 449566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(A,&Istart,&Iend)); 45c4762a1bSJed Brown for (Ii=Istart; Ii<Iend; Ii++) { 46c4762a1bSJed Brown v = -1.0; i = Ii/n; j = Ii - i*n; 479566063dSJacob Faibussowitsch if (i>0) {J = Ii - n; PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,ADD_VALUES));} 489566063dSJacob Faibussowitsch if (i<m-1) {J = Ii + n; PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,ADD_VALUES));} 499566063dSJacob Faibussowitsch if (j>0) {J = Ii - 1; PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,ADD_VALUES));} 509566063dSJacob Faibussowitsch if (j<n-1) {J = Ii + 1; PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,ADD_VALUES));} 519566063dSJacob Faibussowitsch v = 4.0; PetscCall(MatSetValues(A,1,&Ii,1,&Ii,&v,ADD_VALUES)); 52c4762a1bSJed Brown } 539566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 549566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 55c4762a1bSJed Brown 569566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(A,&m,&n)); 579566063dSJacob Faibussowitsch PetscCall(MatGetSize(A,&M,&N)); 5808401ef6SPierre Jolivet PetscCheck(m == n,PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ, "This example is not intended for rectangular matrices (%" PetscInt_FMT ", %" PetscInt_FMT ")", m, n); 59c4762a1bSJed Brown 60a5b23f4aSJose E. Roman /* Create dense matrix C and X; C holds true solution with identical columns */ 61c4762a1bSJed Brown nrhs = N; 629566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-nrhs",&nrhs,NULL)); 639566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD,&C)); 649566063dSJacob Faibussowitsch PetscCall(MatSetSizes(C,m,PETSC_DECIDE,PETSC_DECIDE,nrhs)); 659566063dSJacob Faibussowitsch PetscCall(MatSetType(C,MATDENSE)); 669566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(C)); 679566063dSJacob Faibussowitsch PetscCall(MatSetUp(C)); 68c4762a1bSJed Brown 699566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_WORLD,&rand)); 709566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(rand)); 719566063dSJacob Faibussowitsch PetscCall(MatSetRandom(C,rand)); 729566063dSJacob Faibussowitsch PetscCall(MatDuplicate(C,MAT_DO_NOT_COPY_VALUES,&X)); 73c4762a1bSJed Brown 749566063dSJacob Faibussowitsch PetscCall(PetscStrcpy(solver,MATSOLVERMUMPS)); 759566063dSJacob Faibussowitsch if (rank == 0 && displ) PetscCall(PetscPrintf(PETSC_COMM_SELF,"Solving with %s: nrhs %" PetscInt_FMT ", size mat %" PetscInt_FMT " x %" PetscInt_FMT "\n",solver,nrhs,M,N)); 76c4762a1bSJed Brown 77c4762a1bSJed Brown for (test=0; test<2; test++) { 78c4762a1bSJed Brown if (test == 0) { 79c4762a1bSJed Brown /* Test LU Factorization */ 809566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"test LU factorization\n")); 819566063dSJacob Faibussowitsch PetscCall(MatGetFactor(A,solver,MAT_FACTOR_LU,&F)); 829566063dSJacob Faibussowitsch PetscCall(MatLUFactorSymbolic(F,A,NULL,NULL,NULL)); 839566063dSJacob Faibussowitsch PetscCall(MatLUFactorNumeric(F,A,NULL)); 84c4762a1bSJed Brown } else { 85c4762a1bSJed Brown /* Test Cholesky Factorization */ 86c4762a1bSJed Brown PetscBool flg; 879566063dSJacob Faibussowitsch PetscCall(MatIsSymmetric(A,0.0,&flg)); 8828b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"A must be symmetric for Cholesky factorization"); 89c4762a1bSJed Brown 909566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"test Cholesky factorization\n")); 919566063dSJacob Faibussowitsch PetscCall(MatGetFactor(A,solver,MAT_FACTOR_CHOLESKY,&F)); 929566063dSJacob Faibussowitsch PetscCall(MatCholeskyFactorSymbolic(F,A,NULL,NULL)); 939566063dSJacob Faibussowitsch PetscCall(MatCholeskyFactorNumeric(F,A,NULL)); 94c4762a1bSJed Brown } 95c4762a1bSJed Brown 96c4762a1bSJed Brown /* (1) Test MatMatSolve(): dense RHS = A*C, C: true solutions */ 97c4762a1bSJed Brown /* ---------------------------------------------------------- */ 989566063dSJacob Faibussowitsch PetscCall(MatMatMult(A,C,MAT_INITIAL_MATRIX,2.0,&RHS)); 999566063dSJacob Faibussowitsch PetscCall(MatMatSolve(F,RHS,X)); 100c4762a1bSJed Brown /* Check the error */ 1019566063dSJacob Faibussowitsch PetscCall(MatAXPY(X,-1.0,C,SAME_NONZERO_PATTERN)); 1029566063dSJacob Faibussowitsch PetscCall(MatNorm(X,NORM_FROBENIUS,&norm)); 103c4762a1bSJed Brown if (norm > tol) { 1049566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"(1) MatMatSolve: Norm of error %g\n",norm)); 105c4762a1bSJed Brown } 106c4762a1bSJed Brown 107c4762a1bSJed Brown /* Test X=RHS */ 1089566063dSJacob Faibussowitsch PetscCall(MatMatSolve(F,RHS,RHS)); 109c4762a1bSJed Brown /* Check the error */ 1109566063dSJacob Faibussowitsch PetscCall(MatAXPY(RHS,-1.0,C,SAME_NONZERO_PATTERN)); 1119566063dSJacob Faibussowitsch PetscCall(MatNorm(RHS,NORM_FROBENIUS,&norm)); 112c4762a1bSJed Brown if (norm > tol) { 1139566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"(1.1) MatMatSolve: Norm of error %g\n",norm)); 114c4762a1bSJed Brown } 115c4762a1bSJed Brown 116c4762a1bSJed Brown /* (2) Test MatMatSolve() for inv(A) with dense RHS: 117c4762a1bSJed Brown RHS = [e[0],...,e[nrhs-1]], dense X holds first nrhs columns of inv(A) */ 118c4762a1bSJed Brown /* -------------------------------------------------------------------- */ 1199566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(RHS)); 120c4762a1bSJed Brown for (i=0; i<nrhs; i++) { 121c4762a1bSJed Brown v = 1.0; 1229566063dSJacob Faibussowitsch PetscCall(MatSetValues(RHS,1,&i,1,&i,&v,INSERT_VALUES)); 123c4762a1bSJed Brown } 1249566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(RHS,MAT_FINAL_ASSEMBLY)); 1259566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(RHS,MAT_FINAL_ASSEMBLY)); 126c4762a1bSJed Brown 1279566063dSJacob Faibussowitsch PetscCall(MatMatSolve(F,RHS,X)); 128c4762a1bSJed Brown if (displ) { 1299566063dSJacob Faibussowitsch if (rank == 0) PetscCall(PetscPrintf(PETSC_COMM_SELF," \n(2) first %" PetscInt_FMT " columns of inv(A) with dense RHS:\n",nrhs)); 1309566063dSJacob Faibussowitsch PetscCall(MatView(X,PETSC_VIEWER_STDOUT_WORLD)); 131c4762a1bSJed Brown } 132c4762a1bSJed Brown 133c4762a1bSJed Brown /* Check the residual */ 1349566063dSJacob Faibussowitsch PetscCall(MatMatMult(A,X,MAT_INITIAL_MATRIX,2.0,&AX)); 1359566063dSJacob Faibussowitsch PetscCall(MatAXPY(AX,-1.0,RHS,SAME_NONZERO_PATTERN)); 1369566063dSJacob Faibussowitsch PetscCall(MatNorm(AX,NORM_INFINITY,&norm)); 137c4762a1bSJed Brown if (norm > tol) { 1389566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"(2) MatMatSolve: Norm of residual %g\n",norm)); 139c4762a1bSJed Brown } 1409566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(X)); 141c4762a1bSJed Brown 142c4762a1bSJed Brown /* (3) Test MatMatTransposeSolve() for inv(A) with sparse RHS stored in the host: 143c4762a1bSJed Brown spRHST = [e[0],...,e[nrhs-1]]^T, dense X holds first nrhs columns of inv(A) */ 144c4762a1bSJed Brown /* --------------------------------------------------------------------------- */ 145c4762a1bSJed Brown /* Create spRHST: PETSc does not support compressed column format which is required by MUMPS for sparse RHS matrix, 146c4762a1bSJed Brown thus user must create spRHST=spRHS^T and call MatMatTransposeSolve() */ 1479566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD,&spRHST)); 148dd400576SPatrick Sanan if (rank == 0) { 149c4762a1bSJed Brown /* MUMPS requirs RHS be centralized on the host! */ 1509566063dSJacob Faibussowitsch PetscCall(MatSetSizes(spRHST,nrhs,M,PETSC_DECIDE,PETSC_DECIDE)); 151c4762a1bSJed Brown } else { 1529566063dSJacob Faibussowitsch PetscCall(MatSetSizes(spRHST,0,0,PETSC_DECIDE,PETSC_DECIDE)); 153c4762a1bSJed Brown } 1549566063dSJacob Faibussowitsch PetscCall(MatSetType(spRHST,MATAIJ)); 1559566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(spRHST)); 1569566063dSJacob Faibussowitsch PetscCall(MatSetUp(spRHST)); 157dd400576SPatrick Sanan if (rank == 0) { 158c4762a1bSJed Brown v = 1.0; 159c4762a1bSJed Brown for (i=0; i<nrhs; i++) { 1609566063dSJacob Faibussowitsch PetscCall(MatSetValues(spRHST,1,&i,1,&i,&v,INSERT_VALUES)); 161c4762a1bSJed Brown } 162c4762a1bSJed Brown } 1639566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(spRHST,MAT_FINAL_ASSEMBLY)); 1649566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(spRHST,MAT_FINAL_ASSEMBLY)); 165c4762a1bSJed Brown 1669566063dSJacob Faibussowitsch PetscCall(MatMatTransposeSolve(F,spRHST,X)); 167c4762a1bSJed Brown 168c4762a1bSJed Brown if (displ) { 1699566063dSJacob Faibussowitsch if (rank == 0) PetscCall(PetscPrintf(PETSC_COMM_SELF," \n(3) first %" PetscInt_FMT " columns of inv(A) with sparse RHS:\n",nrhs)); 1709566063dSJacob Faibussowitsch PetscCall(MatView(X,PETSC_VIEWER_STDOUT_WORLD)); 171c4762a1bSJed Brown } 172c4762a1bSJed Brown 173c4762a1bSJed Brown /* Check the residual: R = A*X - RHS */ 1749566063dSJacob Faibussowitsch PetscCall(MatMatMult(A,X,MAT_REUSE_MATRIX,2.0,&AX)); 175c4762a1bSJed Brown 1769566063dSJacob Faibussowitsch PetscCall(MatAXPY(AX,-1.0,RHS,SAME_NONZERO_PATTERN)); 1779566063dSJacob Faibussowitsch PetscCall(MatNorm(AX,NORM_INFINITY,&norm)); 178c4762a1bSJed Brown if (norm > tol) { 1799566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"(3) MatMatSolve: Norm of residual %g\n",norm)); 180c4762a1bSJed Brown } 181c4762a1bSJed Brown 182c4762a1bSJed Brown /* (4) Test MatMatSolve() for inv(A) with selected entries: 183c4762a1bSJed Brown input: spRHS gives selected indices; output: spRHS holds selected entries of inv(A) */ 184c4762a1bSJed Brown /* --------------------------------------------------------------------------------- */ 185c4762a1bSJed Brown if (nrhs == N) { /* mumps requires nrhs = n */ 186c4762a1bSJed Brown /* Create spRHS on proc[0] */ 187c4762a1bSJed Brown Mat spRHS = NULL; 188c4762a1bSJed Brown 189c4762a1bSJed Brown /* Create spRHS = spRHST^T. Two matrices share internal matrix data structure */ 1909566063dSJacob Faibussowitsch PetscCall(MatCreateTranspose(spRHST,&spRHS)); 1919566063dSJacob Faibussowitsch PetscCall(MatMumpsGetInverse(F,spRHS)); 1929566063dSJacob Faibussowitsch PetscCall(MatDestroy(&spRHS)); 193c4762a1bSJed Brown 1949566063dSJacob Faibussowitsch PetscCall(MatMumpsGetInverseTranspose(F,spRHST)); 195c4762a1bSJed Brown if (displ) { 1969566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nSelected entries of inv(A^T):\n")); 1979566063dSJacob Faibussowitsch PetscCall(MatView(spRHST,PETSC_VIEWER_STDOUT_WORLD)); 198c4762a1bSJed Brown } 1999566063dSJacob Faibussowitsch PetscCall(MatDestroy(&spRHS)); 200c4762a1bSJed Brown } 2019566063dSJacob Faibussowitsch PetscCall(MatDestroy(&AX)); 2029566063dSJacob Faibussowitsch PetscCall(MatDestroy(&F)); 2039566063dSJacob Faibussowitsch PetscCall(MatDestroy(&RHS)); 2049566063dSJacob Faibussowitsch PetscCall(MatDestroy(&spRHST)); 205c4762a1bSJed Brown } 206c4762a1bSJed Brown 207c4762a1bSJed Brown /* Free data structures */ 2089566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 2099566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C)); 2109566063dSJacob Faibussowitsch PetscCall(MatDestroy(&X)); 2119566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rand)); 2129566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 213b122ec5aSJacob Faibussowitsch return 0; 214c4762a1bSJed Brown #endif 215c4762a1bSJed Brown } 216c4762a1bSJed Brown 217c4762a1bSJed Brown /*TEST 218c4762a1bSJed Brown 219c4762a1bSJed Brown test: 220c4762a1bSJed Brown requires: mumps double !complex 221c4762a1bSJed Brown 222c4762a1bSJed Brown test: 223c4762a1bSJed Brown suffix: 2 224c4762a1bSJed Brown requires: mumps double !complex 225c4762a1bSJed Brown nsize: 2 226c4762a1bSJed Brown 227c4762a1bSJed Brown TEST*/ 228