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 PetscErrorCode ierr; 10c4762a1bSJed Brown PetscMPIInt size,rank; 11c4762a1bSJed Brown #if defined(PETSC_HAVE_MUMPS) 12c4762a1bSJed Brown Mat A,RHS,C,F,X,AX,spRHST; 13c4762a1bSJed Brown PetscInt m,n,nrhs,M,N,i,Istart,Iend,Ii,j,J,test; 14c4762a1bSJed Brown PetscScalar v; 15c4762a1bSJed Brown PetscReal norm,tol=PETSC_SQRT_MACHINE_EPSILON; 16c4762a1bSJed Brown PetscRandom rand; 17c4762a1bSJed Brown PetscBool displ=PETSC_FALSE; 18c4762a1bSJed Brown char solver[256]; 19c4762a1bSJed Brown #endif 20c4762a1bSJed Brown 21c4762a1bSJed Brown ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr; 22ffc4695bSBarry Smith ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr); 23ffc4695bSBarry Smith ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRMPI(ierr); 24c4762a1bSJed Brown 25c4762a1bSJed Brown #if !defined(PETSC_HAVE_MUMPS) 26dd400576SPatrick Sanan if (rank == 0) {ierr = PetscPrintf(PETSC_COMM_SELF,"This example requires MUMPS, exit...\n");CHKERRQ(ierr);} 27c4762a1bSJed Brown ierr = PetscFinalize(); 28c4762a1bSJed Brown return ierr; 29c4762a1bSJed Brown #else 30c4762a1bSJed Brown 31c4762a1bSJed Brown ierr = PetscOptionsGetBool(NULL,NULL,"-displ",&displ,NULL);CHKERRQ(ierr); 32c4762a1bSJed Brown 33c4762a1bSJed Brown /* Create matrix A */ 34c4762a1bSJed Brown m = 4; n = 4; 35c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL);CHKERRQ(ierr); 36c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);CHKERRQ(ierr); 37c4762a1bSJed Brown 38c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr); 39c4762a1bSJed Brown ierr = MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n);CHKERRQ(ierr); 40c4762a1bSJed Brown ierr = MatSetFromOptions(A);CHKERRQ(ierr); 41c4762a1bSJed Brown ierr = MatMPIAIJSetPreallocation(A,5,NULL,5,NULL);CHKERRQ(ierr); 42c4762a1bSJed Brown ierr = MatSeqAIJSetPreallocation(A,5,NULL);CHKERRQ(ierr); 43c4762a1bSJed Brown 44c4762a1bSJed Brown ierr = MatGetOwnershipRange(A,&Istart,&Iend);CHKERRQ(ierr); 45c4762a1bSJed Brown for (Ii=Istart; Ii<Iend; Ii++) { 46c4762a1bSJed Brown v = -1.0; i = Ii/n; j = Ii - i*n; 47c4762a1bSJed Brown if (i>0) {J = Ii - n; ierr = MatSetValues(A,1,&Ii,1,&J,&v,ADD_VALUES);CHKERRQ(ierr);} 48c4762a1bSJed Brown if (i<m-1) {J = Ii + n; ierr = MatSetValues(A,1,&Ii,1,&J,&v,ADD_VALUES);CHKERRQ(ierr);} 49c4762a1bSJed Brown if (j>0) {J = Ii - 1; ierr = MatSetValues(A,1,&Ii,1,&J,&v,ADD_VALUES);CHKERRQ(ierr);} 50c4762a1bSJed Brown if (j<n-1) {J = Ii + 1; ierr = MatSetValues(A,1,&Ii,1,&J,&v,ADD_VALUES);CHKERRQ(ierr);} 51c4762a1bSJed Brown v = 4.0; ierr = MatSetValues(A,1,&Ii,1,&Ii,&v,ADD_VALUES);CHKERRQ(ierr); 52c4762a1bSJed Brown } 53c4762a1bSJed Brown ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 54c4762a1bSJed Brown ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 55c4762a1bSJed Brown 56c4762a1bSJed Brown ierr = MatGetLocalSize(A,&m,&n);CHKERRQ(ierr); 57c4762a1bSJed Brown ierr = MatGetSize(A,&M,&N);CHKERRQ(ierr); 58*8fffc762SJacob Faibussowitsch if (m != n) SETERRQ2(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; 62c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-nrhs",&nrhs,NULL);CHKERRQ(ierr); 63c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD,&C);CHKERRQ(ierr); 64c4762a1bSJed Brown ierr = MatSetSizes(C,m,PETSC_DECIDE,PETSC_DECIDE,nrhs);CHKERRQ(ierr); 65c4762a1bSJed Brown ierr = MatSetType(C,MATDENSE);CHKERRQ(ierr); 66c4762a1bSJed Brown ierr = MatSetFromOptions(C);CHKERRQ(ierr); 67c4762a1bSJed Brown ierr = MatSetUp(C);CHKERRQ(ierr); 68c4762a1bSJed Brown 69c4762a1bSJed Brown ierr = PetscRandomCreate(PETSC_COMM_WORLD,&rand);CHKERRQ(ierr); 70c4762a1bSJed Brown ierr = PetscRandomSetFromOptions(rand);CHKERRQ(ierr); 71c4762a1bSJed Brown ierr = MatSetRandom(C,rand);CHKERRQ(ierr); 72c4762a1bSJed Brown ierr = MatDuplicate(C,MAT_DO_NOT_COPY_VALUES,&X);CHKERRQ(ierr); 73c4762a1bSJed Brown 74c4762a1bSJed Brown ierr = PetscStrcpy(solver,MATSOLVERMUMPS);CHKERRQ(ierr); 75*8fffc762SJacob Faibussowitsch if (rank == 0 && displ) {ierr = PetscPrintf(PETSC_COMM_SELF,"Solving with %s: nrhs %" PetscInt_FMT ", size mat %" PetscInt_FMT " x %" PetscInt_FMT "\n",solver,nrhs,M,N);CHKERRQ(ierr);} 76c4762a1bSJed Brown 77c4762a1bSJed Brown for (test=0; test<2; test++) { 78c4762a1bSJed Brown if (test == 0) { 79c4762a1bSJed Brown /* Test LU Factorization */ 80c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"test LU factorization\n");CHKERRQ(ierr); 81c4762a1bSJed Brown ierr = MatGetFactor(A,solver,MAT_FACTOR_LU,&F);CHKERRQ(ierr); 82c4762a1bSJed Brown ierr = MatLUFactorSymbolic(F,A,NULL,NULL,NULL);CHKERRQ(ierr); 83c4762a1bSJed Brown ierr = MatLUFactorNumeric(F,A,NULL);CHKERRQ(ierr); 84c4762a1bSJed Brown } else { 85c4762a1bSJed Brown /* Test Cholesky Factorization */ 86c4762a1bSJed Brown PetscBool flg; 87c4762a1bSJed Brown ierr = MatIsSymmetric(A,0.0,&flg);CHKERRQ(ierr); 88c4762a1bSJed Brown if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"A must be symmetric for Cholesky factorization"); 89c4762a1bSJed Brown 90c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"test Cholesky factorization\n");CHKERRQ(ierr); 91c4762a1bSJed Brown ierr = MatGetFactor(A,solver,MAT_FACTOR_CHOLESKY,&F);CHKERRQ(ierr); 92c4762a1bSJed Brown ierr = MatCholeskyFactorSymbolic(F,A,NULL,NULL);CHKERRQ(ierr); 93c4762a1bSJed Brown ierr = MatCholeskyFactorNumeric(F,A,NULL);CHKERRQ(ierr); 94c4762a1bSJed Brown } 95c4762a1bSJed Brown 96c4762a1bSJed Brown /* (1) Test MatMatSolve(): dense RHS = A*C, C: true solutions */ 97c4762a1bSJed Brown /* ---------------------------------------------------------- */ 98c4762a1bSJed Brown ierr = MatMatMult(A,C,MAT_INITIAL_MATRIX,2.0,&RHS);CHKERRQ(ierr); 99c4762a1bSJed Brown ierr = MatMatSolve(F,RHS,X);CHKERRQ(ierr); 100c4762a1bSJed Brown /* Check the error */ 101c4762a1bSJed Brown ierr = MatAXPY(X,-1.0,C,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 102c4762a1bSJed Brown ierr = MatNorm(X,NORM_FROBENIUS,&norm);CHKERRQ(ierr); 103c4762a1bSJed Brown if (norm > tol) { 104c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"(1) MatMatSolve: Norm of error %g\n",norm);CHKERRQ(ierr); 105c4762a1bSJed Brown } 106c4762a1bSJed Brown 107c4762a1bSJed Brown /* Test X=RHS */ 108c4762a1bSJed Brown ierr = MatMatSolve(F,RHS,RHS);CHKERRQ(ierr); 109c4762a1bSJed Brown /* Check the error */ 110c4762a1bSJed Brown ierr = MatAXPY(RHS,-1.0,C,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 111c4762a1bSJed Brown ierr = MatNorm(RHS,NORM_FROBENIUS,&norm);CHKERRQ(ierr); 112c4762a1bSJed Brown if (norm > tol) { 113c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"(1.1) MatMatSolve: Norm of error %g\n",norm);CHKERRQ(ierr); 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 /* -------------------------------------------------------------------- */ 119c4762a1bSJed Brown ierr = MatZeroEntries(RHS);CHKERRQ(ierr); 120c4762a1bSJed Brown for (i=0; i<nrhs; i++) { 121c4762a1bSJed Brown v = 1.0; 122c4762a1bSJed Brown ierr = MatSetValues(RHS,1,&i,1,&i,&v,INSERT_VALUES);CHKERRQ(ierr); 123c4762a1bSJed Brown } 124c4762a1bSJed Brown ierr = MatAssemblyBegin(RHS,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 125c4762a1bSJed Brown ierr = MatAssemblyEnd(RHS,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 126c4762a1bSJed Brown 127c4762a1bSJed Brown ierr = MatMatSolve(F,RHS,X);CHKERRQ(ierr); 128c4762a1bSJed Brown if (displ) { 129*8fffc762SJacob Faibussowitsch if (rank == 0) {ierr = PetscPrintf(PETSC_COMM_SELF," \n(2) first %" PetscInt_FMT " columns of inv(A) with dense RHS:\n",nrhs);CHKERRQ(ierr);} 130c4762a1bSJed Brown ierr = MatView(X,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 131c4762a1bSJed Brown } 132c4762a1bSJed Brown 133c4762a1bSJed Brown /* Check the residual */ 134c4762a1bSJed Brown ierr = MatMatMult(A,X,MAT_INITIAL_MATRIX,2.0,&AX);CHKERRQ(ierr); 135c4762a1bSJed Brown ierr = MatAXPY(AX,-1.0,RHS,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 136c4762a1bSJed Brown ierr = MatNorm(AX,NORM_INFINITY,&norm);CHKERRQ(ierr); 137c4762a1bSJed Brown if (norm > tol) { 138c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"(2) MatMatSolve: Norm of residual %g\n",norm);CHKERRQ(ierr); 139c4762a1bSJed Brown } 140c4762a1bSJed Brown ierr = MatZeroEntries(X);CHKERRQ(ierr); 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() */ 147c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD,&spRHST);CHKERRQ(ierr); 148dd400576SPatrick Sanan if (rank == 0) { 149c4762a1bSJed Brown /* MUMPS requirs RHS be centralized on the host! */ 150c4762a1bSJed Brown ierr = MatSetSizes(spRHST,nrhs,M,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr); 151c4762a1bSJed Brown } else { 152c4762a1bSJed Brown ierr = MatSetSizes(spRHST,0,0,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr); 153c4762a1bSJed Brown } 154c4762a1bSJed Brown ierr = MatSetType(spRHST,MATAIJ);CHKERRQ(ierr); 155c4762a1bSJed Brown ierr = MatSetFromOptions(spRHST);CHKERRQ(ierr); 156c4762a1bSJed Brown ierr = MatSetUp(spRHST);CHKERRQ(ierr); 157dd400576SPatrick Sanan if (rank == 0) { 158c4762a1bSJed Brown v = 1.0; 159c4762a1bSJed Brown for (i=0; i<nrhs; i++) { 160c4762a1bSJed Brown ierr = MatSetValues(spRHST,1,&i,1,&i,&v,INSERT_VALUES);CHKERRQ(ierr); 161c4762a1bSJed Brown } 162c4762a1bSJed Brown } 163c4762a1bSJed Brown ierr = MatAssemblyBegin(spRHST,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 164c4762a1bSJed Brown ierr = MatAssemblyEnd(spRHST,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 165c4762a1bSJed Brown 166c4762a1bSJed Brown ierr = MatMatTransposeSolve(F,spRHST,X);CHKERRQ(ierr); 167c4762a1bSJed Brown 168c4762a1bSJed Brown if (displ) { 169*8fffc762SJacob Faibussowitsch if (rank == 0) {ierr = PetscPrintf(PETSC_COMM_SELF," \n(3) first %" PetscInt_FMT " columns of inv(A) with sparse RHS:\n",nrhs);CHKERRQ(ierr);} 170c4762a1bSJed Brown ierr = MatView(X,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 171c4762a1bSJed Brown } 172c4762a1bSJed Brown 173c4762a1bSJed Brown /* Check the residual: R = A*X - RHS */ 174c4762a1bSJed Brown ierr = MatMatMult(A,X,MAT_REUSE_MATRIX,2.0,&AX);CHKERRQ(ierr); 175c4762a1bSJed Brown 176c4762a1bSJed Brown ierr = MatAXPY(AX,-1.0,RHS,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 177c4762a1bSJed Brown ierr = MatNorm(AX,NORM_INFINITY,&norm);CHKERRQ(ierr); 178c4762a1bSJed Brown if (norm > tol) { 179c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"(3) MatMatSolve: Norm of residual %g\n",norm);CHKERRQ(ierr); 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 */ 190c4762a1bSJed Brown ierr = MatCreateTranspose(spRHST,&spRHS);CHKERRQ(ierr); 191c4762a1bSJed Brown ierr = MatMumpsGetInverse(F,spRHS);CHKERRQ(ierr); 192c4762a1bSJed Brown ierr = MatDestroy(&spRHS);CHKERRQ(ierr); 193c4762a1bSJed Brown 194c4762a1bSJed Brown ierr = MatMumpsGetInverseTranspose(F,spRHST);CHKERRQ(ierr); 195c4762a1bSJed Brown if (displ) { 196c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"\nSelected entries of inv(A^T):\n");CHKERRQ(ierr); 197c4762a1bSJed Brown ierr = MatView(spRHST,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 198c4762a1bSJed Brown } 199c4762a1bSJed Brown ierr = MatDestroy(&spRHS);CHKERRQ(ierr); 200c4762a1bSJed Brown } 201c4762a1bSJed Brown ierr = MatDestroy(&AX);CHKERRQ(ierr); 202c4762a1bSJed Brown ierr = MatDestroy(&F);CHKERRQ(ierr); 203c4762a1bSJed Brown ierr = MatDestroy(&RHS);CHKERRQ(ierr); 204c4762a1bSJed Brown ierr = MatDestroy(&spRHST);CHKERRQ(ierr); 205c4762a1bSJed Brown } 206c4762a1bSJed Brown 207c4762a1bSJed Brown /* Free data structures */ 208c4762a1bSJed Brown ierr = MatDestroy(&A);CHKERRQ(ierr); 209c4762a1bSJed Brown ierr = MatDestroy(&C);CHKERRQ(ierr); 210c4762a1bSJed Brown ierr = MatDestroy(&X);CHKERRQ(ierr); 211c4762a1bSJed Brown ierr = PetscRandomDestroy(&rand);CHKERRQ(ierr); 212c4762a1bSJed Brown ierr = PetscFinalize(); 213c4762a1bSJed Brown return ierr; 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