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 209566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&args,(char*)0,help)); 219566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 229566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 23c4762a1bSJed Brown 24c4762a1bSJed Brown #if !defined(PETSC_HAVE_MUMPS) 259566063dSJacob Faibussowitsch if (rank == 0) PetscCall(PetscPrintf(PETSC_COMM_SELF,"This example requires MUMPS, exit...\n")); 269566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 27b122ec5aSJacob Faibussowitsch return 0; 28c4762a1bSJed Brown #else 29c4762a1bSJed Brown 309566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL,NULL,"-displ",&displ,NULL)); 31c4762a1bSJed Brown 32c4762a1bSJed Brown /* Create matrix A */ 33c4762a1bSJed Brown m = 4; n = 4; 349566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL)); 359566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL)); 36c4762a1bSJed Brown 379566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD,&A)); 389566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n)); 399566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A)); 409566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(A,5,NULL,5,NULL)); 419566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(A,5,NULL)); 42c4762a1bSJed Brown 439566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(A,&Istart,&Iend)); 44c4762a1bSJed Brown for (Ii=Istart; Ii<Iend; Ii++) { 45c4762a1bSJed Brown v = -1.0; i = Ii/n; j = Ii - i*n; 469566063dSJacob Faibussowitsch if (i>0) {J = Ii - n; PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,ADD_VALUES));} 479566063dSJacob Faibussowitsch if (i<m-1) {J = Ii + n; PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,ADD_VALUES));} 489566063dSJacob Faibussowitsch if (j>0) {J = Ii - 1; PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,ADD_VALUES));} 499566063dSJacob Faibussowitsch if (j<n-1) {J = Ii + 1; PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,ADD_VALUES));} 509566063dSJacob Faibussowitsch v = 4.0; PetscCall(MatSetValues(A,1,&Ii,1,&Ii,&v,ADD_VALUES)); 51c4762a1bSJed Brown } 529566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 539566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 54c4762a1bSJed Brown 559566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(A,&m,&n)); 569566063dSJacob Faibussowitsch PetscCall(MatGetSize(A,&M,&N)); 57*08401ef6SPierre 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); 58c4762a1bSJed Brown 59a5b23f4aSJose E. Roman /* Create dense matrix C and X; C holds true solution with identical columns */ 60c4762a1bSJed Brown nrhs = N; 619566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-nrhs",&nrhs,NULL)); 629566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD,&C)); 639566063dSJacob Faibussowitsch PetscCall(MatSetSizes(C,m,PETSC_DECIDE,PETSC_DECIDE,nrhs)); 649566063dSJacob Faibussowitsch PetscCall(MatSetType(C,MATDENSE)); 659566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(C)); 669566063dSJacob Faibussowitsch PetscCall(MatSetUp(C)); 67c4762a1bSJed Brown 689566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_WORLD,&rand)); 699566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(rand)); 709566063dSJacob Faibussowitsch PetscCall(MatSetRandom(C,rand)); 719566063dSJacob Faibussowitsch PetscCall(MatDuplicate(C,MAT_DO_NOT_COPY_VALUES,&X)); 72c4762a1bSJed Brown 739566063dSJacob Faibussowitsch PetscCall(PetscStrcpy(solver,MATSOLVERMUMPS)); 749566063dSJacob 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)); 75c4762a1bSJed Brown 76c4762a1bSJed Brown for (test=0; test<2; test++) { 77c4762a1bSJed Brown if (test == 0) { 78c4762a1bSJed Brown /* Test LU Factorization */ 799566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"test LU factorization\n")); 809566063dSJacob Faibussowitsch PetscCall(MatGetFactor(A,solver,MAT_FACTOR_LU,&F)); 819566063dSJacob Faibussowitsch PetscCall(MatLUFactorSymbolic(F,A,NULL,NULL,NULL)); 829566063dSJacob Faibussowitsch PetscCall(MatLUFactorNumeric(F,A,NULL)); 83c4762a1bSJed Brown } else { 84c4762a1bSJed Brown /* Test Cholesky Factorization */ 85c4762a1bSJed Brown PetscBool flg; 869566063dSJacob Faibussowitsch PetscCall(MatIsSymmetric(A,0.0,&flg)); 8728b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"A must be symmetric for Cholesky factorization"); 88c4762a1bSJed Brown 899566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"test Cholesky factorization\n")); 909566063dSJacob Faibussowitsch PetscCall(MatGetFactor(A,solver,MAT_FACTOR_CHOLESKY,&F)); 919566063dSJacob Faibussowitsch PetscCall(MatCholeskyFactorSymbolic(F,A,NULL,NULL)); 929566063dSJacob Faibussowitsch PetscCall(MatCholeskyFactorNumeric(F,A,NULL)); 93c4762a1bSJed Brown } 94c4762a1bSJed Brown 95c4762a1bSJed Brown /* (1) Test MatMatSolve(): dense RHS = A*C, C: true solutions */ 96c4762a1bSJed Brown /* ---------------------------------------------------------- */ 979566063dSJacob Faibussowitsch PetscCall(MatMatMult(A,C,MAT_INITIAL_MATRIX,2.0,&RHS)); 989566063dSJacob Faibussowitsch PetscCall(MatMatSolve(F,RHS,X)); 99c4762a1bSJed Brown /* Check the error */ 1009566063dSJacob Faibussowitsch PetscCall(MatAXPY(X,-1.0,C,SAME_NONZERO_PATTERN)); 1019566063dSJacob Faibussowitsch PetscCall(MatNorm(X,NORM_FROBENIUS,&norm)); 102c4762a1bSJed Brown if (norm > tol) { 1039566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"(1) MatMatSolve: Norm of error %g\n",norm)); 104c4762a1bSJed Brown } 105c4762a1bSJed Brown 106c4762a1bSJed Brown /* Test X=RHS */ 1079566063dSJacob Faibussowitsch PetscCall(MatMatSolve(F,RHS,RHS)); 108c4762a1bSJed Brown /* Check the error */ 1099566063dSJacob Faibussowitsch PetscCall(MatAXPY(RHS,-1.0,C,SAME_NONZERO_PATTERN)); 1109566063dSJacob Faibussowitsch PetscCall(MatNorm(RHS,NORM_FROBENIUS,&norm)); 111c4762a1bSJed Brown if (norm > tol) { 1129566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"(1.1) MatMatSolve: Norm of error %g\n",norm)); 113c4762a1bSJed Brown } 114c4762a1bSJed Brown 115c4762a1bSJed Brown /* (2) Test MatMatSolve() for inv(A) with dense RHS: 116c4762a1bSJed Brown RHS = [e[0],...,e[nrhs-1]], dense X holds first nrhs columns of inv(A) */ 117c4762a1bSJed Brown /* -------------------------------------------------------------------- */ 1189566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(RHS)); 119c4762a1bSJed Brown for (i=0; i<nrhs; i++) { 120c4762a1bSJed Brown v = 1.0; 1219566063dSJacob Faibussowitsch PetscCall(MatSetValues(RHS,1,&i,1,&i,&v,INSERT_VALUES)); 122c4762a1bSJed Brown } 1239566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(RHS,MAT_FINAL_ASSEMBLY)); 1249566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(RHS,MAT_FINAL_ASSEMBLY)); 125c4762a1bSJed Brown 1269566063dSJacob Faibussowitsch PetscCall(MatMatSolve(F,RHS,X)); 127c4762a1bSJed Brown if (displ) { 1289566063dSJacob Faibussowitsch if (rank == 0) PetscCall(PetscPrintf(PETSC_COMM_SELF," \n(2) first %" PetscInt_FMT " columns of inv(A) with dense RHS:\n",nrhs)); 1299566063dSJacob Faibussowitsch PetscCall(MatView(X,PETSC_VIEWER_STDOUT_WORLD)); 130c4762a1bSJed Brown } 131c4762a1bSJed Brown 132c4762a1bSJed Brown /* Check the residual */ 1339566063dSJacob Faibussowitsch PetscCall(MatMatMult(A,X,MAT_INITIAL_MATRIX,2.0,&AX)); 1349566063dSJacob Faibussowitsch PetscCall(MatAXPY(AX,-1.0,RHS,SAME_NONZERO_PATTERN)); 1359566063dSJacob Faibussowitsch PetscCall(MatNorm(AX,NORM_INFINITY,&norm)); 136c4762a1bSJed Brown if (norm > tol) { 1379566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"(2) MatMatSolve: Norm of residual %g\n",norm)); 138c4762a1bSJed Brown } 1399566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(X)); 140c4762a1bSJed Brown 141c4762a1bSJed Brown /* (3) Test MatMatTransposeSolve() for inv(A) with sparse RHS stored in the host: 142c4762a1bSJed Brown spRHST = [e[0],...,e[nrhs-1]]^T, dense X holds first nrhs columns of inv(A) */ 143c4762a1bSJed Brown /* --------------------------------------------------------------------------- */ 144c4762a1bSJed Brown /* Create spRHST: PETSc does not support compressed column format which is required by MUMPS for sparse RHS matrix, 145c4762a1bSJed Brown thus user must create spRHST=spRHS^T and call MatMatTransposeSolve() */ 1469566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD,&spRHST)); 147dd400576SPatrick Sanan if (rank == 0) { 148c4762a1bSJed Brown /* MUMPS requirs RHS be centralized on the host! */ 1499566063dSJacob Faibussowitsch PetscCall(MatSetSizes(spRHST,nrhs,M,PETSC_DECIDE,PETSC_DECIDE)); 150c4762a1bSJed Brown } else { 1519566063dSJacob Faibussowitsch PetscCall(MatSetSizes(spRHST,0,0,PETSC_DECIDE,PETSC_DECIDE)); 152c4762a1bSJed Brown } 1539566063dSJacob Faibussowitsch PetscCall(MatSetType(spRHST,MATAIJ)); 1549566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(spRHST)); 1559566063dSJacob Faibussowitsch PetscCall(MatSetUp(spRHST)); 156dd400576SPatrick Sanan if (rank == 0) { 157c4762a1bSJed Brown v = 1.0; 158c4762a1bSJed Brown for (i=0; i<nrhs; i++) { 1599566063dSJacob Faibussowitsch PetscCall(MatSetValues(spRHST,1,&i,1,&i,&v,INSERT_VALUES)); 160c4762a1bSJed Brown } 161c4762a1bSJed Brown } 1629566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(spRHST,MAT_FINAL_ASSEMBLY)); 1639566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(spRHST,MAT_FINAL_ASSEMBLY)); 164c4762a1bSJed Brown 1659566063dSJacob Faibussowitsch PetscCall(MatMatTransposeSolve(F,spRHST,X)); 166c4762a1bSJed Brown 167c4762a1bSJed Brown if (displ) { 1689566063dSJacob Faibussowitsch if (rank == 0) PetscCall(PetscPrintf(PETSC_COMM_SELF," \n(3) first %" PetscInt_FMT " columns of inv(A) with sparse RHS:\n",nrhs)); 1699566063dSJacob Faibussowitsch PetscCall(MatView(X,PETSC_VIEWER_STDOUT_WORLD)); 170c4762a1bSJed Brown } 171c4762a1bSJed Brown 172c4762a1bSJed Brown /* Check the residual: R = A*X - RHS */ 1739566063dSJacob Faibussowitsch PetscCall(MatMatMult(A,X,MAT_REUSE_MATRIX,2.0,&AX)); 174c4762a1bSJed Brown 1759566063dSJacob Faibussowitsch PetscCall(MatAXPY(AX,-1.0,RHS,SAME_NONZERO_PATTERN)); 1769566063dSJacob Faibussowitsch PetscCall(MatNorm(AX,NORM_INFINITY,&norm)); 177c4762a1bSJed Brown if (norm > tol) { 1789566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"(3) MatMatSolve: Norm of residual %g\n",norm)); 179c4762a1bSJed Brown } 180c4762a1bSJed Brown 181c4762a1bSJed Brown /* (4) Test MatMatSolve() for inv(A) with selected entries: 182c4762a1bSJed Brown input: spRHS gives selected indices; output: spRHS holds selected entries of inv(A) */ 183c4762a1bSJed Brown /* --------------------------------------------------------------------------------- */ 184c4762a1bSJed Brown if (nrhs == N) { /* mumps requires nrhs = n */ 185c4762a1bSJed Brown /* Create spRHS on proc[0] */ 186c4762a1bSJed Brown Mat spRHS = NULL; 187c4762a1bSJed Brown 188c4762a1bSJed Brown /* Create spRHS = spRHST^T. Two matrices share internal matrix data structure */ 1899566063dSJacob Faibussowitsch PetscCall(MatCreateTranspose(spRHST,&spRHS)); 1909566063dSJacob Faibussowitsch PetscCall(MatMumpsGetInverse(F,spRHS)); 1919566063dSJacob Faibussowitsch PetscCall(MatDestroy(&spRHS)); 192c4762a1bSJed Brown 1939566063dSJacob Faibussowitsch PetscCall(MatMumpsGetInverseTranspose(F,spRHST)); 194c4762a1bSJed Brown if (displ) { 1959566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nSelected entries of inv(A^T):\n")); 1969566063dSJacob Faibussowitsch PetscCall(MatView(spRHST,PETSC_VIEWER_STDOUT_WORLD)); 197c4762a1bSJed Brown } 1989566063dSJacob Faibussowitsch PetscCall(MatDestroy(&spRHS)); 199c4762a1bSJed Brown } 2009566063dSJacob Faibussowitsch PetscCall(MatDestroy(&AX)); 2019566063dSJacob Faibussowitsch PetscCall(MatDestroy(&F)); 2029566063dSJacob Faibussowitsch PetscCall(MatDestroy(&RHS)); 2039566063dSJacob Faibussowitsch PetscCall(MatDestroy(&spRHST)); 204c4762a1bSJed Brown } 205c4762a1bSJed Brown 206c4762a1bSJed Brown /* Free data structures */ 2079566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 2089566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C)); 2099566063dSJacob Faibussowitsch PetscCall(MatDestroy(&X)); 2109566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rand)); 2119566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 212b122ec5aSJacob Faibussowitsch return 0; 213c4762a1bSJed Brown #endif 214c4762a1bSJed Brown } 215c4762a1bSJed Brown 216c4762a1bSJed Brown /*TEST 217c4762a1bSJed Brown 218c4762a1bSJed Brown test: 219c4762a1bSJed Brown requires: mumps double !complex 220c4762a1bSJed Brown 221c4762a1bSJed Brown test: 222c4762a1bSJed Brown suffix: 2 223c4762a1bSJed Brown requires: mumps double !complex 224c4762a1bSJed Brown nsize: 2 225c4762a1bSJed Brown 226c4762a1bSJed Brown TEST*/ 227