1*c4762a1bSJed Brown 2*c4762a1bSJed Brown static char help[] = "Tests MatMatSolve() and MatMatTransposeSolve() for computing inv(A) with MUMPS.\n\ 3*c4762a1bSJed Brown Example: mpiexec -n <np> ./ex214 -displ \n\n"; 4*c4762a1bSJed Brown 5*c4762a1bSJed Brown #include <petscmat.h> 6*c4762a1bSJed Brown 7*c4762a1bSJed Brown int main(int argc,char **args) 8*c4762a1bSJed Brown { 9*c4762a1bSJed Brown PetscErrorCode ierr; 10*c4762a1bSJed Brown PetscMPIInt size,rank; 11*c4762a1bSJed Brown #if defined(PETSC_HAVE_MUMPS) 12*c4762a1bSJed Brown Mat A,RHS,C,F,X,AX,spRHST; 13*c4762a1bSJed Brown PetscInt m,n,nrhs,M,N,i,Istart,Iend,Ii,j,J,test; 14*c4762a1bSJed Brown PetscScalar v; 15*c4762a1bSJed Brown PetscReal norm,tol=PETSC_SQRT_MACHINE_EPSILON; 16*c4762a1bSJed Brown PetscRandom rand; 17*c4762a1bSJed Brown PetscBool displ=PETSC_FALSE; 18*c4762a1bSJed Brown char solver[256]; 19*c4762a1bSJed Brown #endif 20*c4762a1bSJed Brown 21*c4762a1bSJed Brown ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr; 22*c4762a1bSJed Brown ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr); 23*c4762a1bSJed Brown ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr); 24*c4762a1bSJed Brown 25*c4762a1bSJed Brown #if !defined(PETSC_HAVE_MUMPS) 26*c4762a1bSJed Brown if (!rank) {ierr = PetscPrintf(PETSC_COMM_SELF,"This example requires MUMPS, exit...\n");CHKERRQ(ierr);} 27*c4762a1bSJed Brown ierr = PetscFinalize(); 28*c4762a1bSJed Brown return ierr; 29*c4762a1bSJed Brown #else 30*c4762a1bSJed Brown 31*c4762a1bSJed Brown ierr = PetscOptionsGetBool(NULL,NULL,"-displ",&displ,NULL);CHKERRQ(ierr); 32*c4762a1bSJed Brown 33*c4762a1bSJed Brown /* Create matrix A */ 34*c4762a1bSJed Brown m = 4; n = 4; 35*c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL);CHKERRQ(ierr); 36*c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);CHKERRQ(ierr); 37*c4762a1bSJed Brown 38*c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr); 39*c4762a1bSJed Brown ierr = MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n);CHKERRQ(ierr); 40*c4762a1bSJed Brown ierr = MatSetFromOptions(A);CHKERRQ(ierr); 41*c4762a1bSJed Brown ierr = MatMPIAIJSetPreallocation(A,5,NULL,5,NULL);CHKERRQ(ierr); 42*c4762a1bSJed Brown ierr = MatSeqAIJSetPreallocation(A,5,NULL);CHKERRQ(ierr); 43*c4762a1bSJed Brown 44*c4762a1bSJed Brown ierr = MatGetOwnershipRange(A,&Istart,&Iend);CHKERRQ(ierr); 45*c4762a1bSJed Brown for (Ii=Istart; Ii<Iend; Ii++) { 46*c4762a1bSJed Brown v = -1.0; i = Ii/n; j = Ii - i*n; 47*c4762a1bSJed Brown if (i>0) {J = Ii - n; ierr = MatSetValues(A,1,&Ii,1,&J,&v,ADD_VALUES);CHKERRQ(ierr);} 48*c4762a1bSJed Brown if (i<m-1) {J = Ii + n; ierr = MatSetValues(A,1,&Ii,1,&J,&v,ADD_VALUES);CHKERRQ(ierr);} 49*c4762a1bSJed Brown if (j>0) {J = Ii - 1; ierr = MatSetValues(A,1,&Ii,1,&J,&v,ADD_VALUES);CHKERRQ(ierr);} 50*c4762a1bSJed Brown if (j<n-1) {J = Ii + 1; ierr = MatSetValues(A,1,&Ii,1,&J,&v,ADD_VALUES);CHKERRQ(ierr);} 51*c4762a1bSJed Brown v = 4.0; ierr = MatSetValues(A,1,&Ii,1,&Ii,&v,ADD_VALUES);CHKERRQ(ierr); 52*c4762a1bSJed Brown } 53*c4762a1bSJed Brown ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 54*c4762a1bSJed Brown ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 55*c4762a1bSJed Brown 56*c4762a1bSJed Brown ierr = MatGetLocalSize(A,&m,&n);CHKERRQ(ierr); 57*c4762a1bSJed Brown ierr = MatGetSize(A,&M,&N);CHKERRQ(ierr); 58*c4762a1bSJed Brown if (m != n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ, "This example is not intended for rectangular matrices (%d, %d)", m, n); 59*c4762a1bSJed Brown 60*c4762a1bSJed Brown /* Create dense matrix C and X; C holds true solution with identical colums */ 61*c4762a1bSJed Brown nrhs = N; 62*c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-nrhs",&nrhs,NULL);CHKERRQ(ierr); 63*c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD,&C);CHKERRQ(ierr); 64*c4762a1bSJed Brown ierr = MatSetSizes(C,m,PETSC_DECIDE,PETSC_DECIDE,nrhs);CHKERRQ(ierr); 65*c4762a1bSJed Brown ierr = MatSetType(C,MATDENSE);CHKERRQ(ierr); 66*c4762a1bSJed Brown ierr = MatSetFromOptions(C);CHKERRQ(ierr); 67*c4762a1bSJed Brown ierr = MatSetUp(C);CHKERRQ(ierr); 68*c4762a1bSJed Brown 69*c4762a1bSJed Brown ierr = PetscRandomCreate(PETSC_COMM_WORLD,&rand);CHKERRQ(ierr); 70*c4762a1bSJed Brown ierr = PetscRandomSetFromOptions(rand);CHKERRQ(ierr); 71*c4762a1bSJed Brown ierr = MatSetRandom(C,rand);CHKERRQ(ierr); 72*c4762a1bSJed Brown ierr = MatDuplicate(C,MAT_DO_NOT_COPY_VALUES,&X);CHKERRQ(ierr); 73*c4762a1bSJed Brown 74*c4762a1bSJed Brown ierr = PetscStrcpy(solver,MATSOLVERMUMPS);CHKERRQ(ierr); 75*c4762a1bSJed Brown if (!rank && displ) {ierr = PetscPrintf(PETSC_COMM_SELF,"Solving with %s: nrhs %D, size mat %D x %D\n",solver,nrhs,M,N);CHKERRQ(ierr);} 76*c4762a1bSJed Brown 77*c4762a1bSJed Brown for (test=0; test<2; test++) { 78*c4762a1bSJed Brown if (test == 0) { 79*c4762a1bSJed Brown /* Test LU Factorization */ 80*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"test LU factorization\n");CHKERRQ(ierr); 81*c4762a1bSJed Brown ierr = MatGetFactor(A,solver,MAT_FACTOR_LU,&F);CHKERRQ(ierr); 82*c4762a1bSJed Brown ierr = MatLUFactorSymbolic(F,A,NULL,NULL,NULL);CHKERRQ(ierr); 83*c4762a1bSJed Brown ierr = MatLUFactorNumeric(F,A,NULL);CHKERRQ(ierr); 84*c4762a1bSJed Brown } else { 85*c4762a1bSJed Brown /* Test Cholesky Factorization */ 86*c4762a1bSJed Brown PetscBool flg; 87*c4762a1bSJed Brown ierr = MatIsSymmetric(A,0.0,&flg);CHKERRQ(ierr); 88*c4762a1bSJed Brown if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"A must be symmetric for Cholesky factorization"); 89*c4762a1bSJed Brown 90*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"test Cholesky factorization\n");CHKERRQ(ierr); 91*c4762a1bSJed Brown ierr = MatGetFactor(A,solver,MAT_FACTOR_CHOLESKY,&F);CHKERRQ(ierr); 92*c4762a1bSJed Brown ierr = MatCholeskyFactorSymbolic(F,A,NULL,NULL);CHKERRQ(ierr); 93*c4762a1bSJed Brown ierr = MatCholeskyFactorNumeric(F,A,NULL);CHKERRQ(ierr); 94*c4762a1bSJed Brown } 95*c4762a1bSJed Brown 96*c4762a1bSJed Brown /* (1) Test MatMatSolve(): dense RHS = A*C, C: true solutions */ 97*c4762a1bSJed Brown /* ---------------------------------------------------------- */ 98*c4762a1bSJed Brown ierr = MatMatMult(A,C,MAT_INITIAL_MATRIX,2.0,&RHS);CHKERRQ(ierr); 99*c4762a1bSJed Brown ierr = MatMatSolve(F,RHS,X);CHKERRQ(ierr); 100*c4762a1bSJed Brown /* Check the error */ 101*c4762a1bSJed Brown ierr = MatAXPY(X,-1.0,C,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 102*c4762a1bSJed Brown ierr = MatNorm(X,NORM_FROBENIUS,&norm);CHKERRQ(ierr); 103*c4762a1bSJed Brown if (norm > tol) { 104*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"(1) MatMatSolve: Norm of error %g\n",norm);CHKERRQ(ierr); 105*c4762a1bSJed Brown } 106*c4762a1bSJed Brown 107*c4762a1bSJed Brown /* Test X=RHS */ 108*c4762a1bSJed Brown ierr = MatMatSolve(F,RHS,RHS);CHKERRQ(ierr); 109*c4762a1bSJed Brown /* Check the error */ 110*c4762a1bSJed Brown ierr = MatAXPY(RHS,-1.0,C,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 111*c4762a1bSJed Brown ierr = MatNorm(RHS,NORM_FROBENIUS,&norm);CHKERRQ(ierr); 112*c4762a1bSJed Brown if (norm > tol) { 113*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"(1.1) MatMatSolve: Norm of error %g\n",norm);CHKERRQ(ierr); 114*c4762a1bSJed Brown } 115*c4762a1bSJed Brown 116*c4762a1bSJed Brown /* (2) Test MatMatSolve() for inv(A) with dense RHS: 117*c4762a1bSJed Brown RHS = [e[0],...,e[nrhs-1]], dense X holds first nrhs columns of inv(A) */ 118*c4762a1bSJed Brown /* -------------------------------------------------------------------- */ 119*c4762a1bSJed Brown ierr = MatZeroEntries(RHS);CHKERRQ(ierr); 120*c4762a1bSJed Brown for (i=0; i<nrhs; i++) { 121*c4762a1bSJed Brown v = 1.0; 122*c4762a1bSJed Brown ierr = MatSetValues(RHS,1,&i,1,&i,&v,INSERT_VALUES);CHKERRQ(ierr); 123*c4762a1bSJed Brown } 124*c4762a1bSJed Brown ierr = MatAssemblyBegin(RHS,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 125*c4762a1bSJed Brown ierr = MatAssemblyEnd(RHS,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 126*c4762a1bSJed Brown 127*c4762a1bSJed Brown ierr = MatMatSolve(F,RHS,X);CHKERRQ(ierr); 128*c4762a1bSJed Brown if (displ) { 129*c4762a1bSJed Brown if (!rank) {ierr = PetscPrintf(PETSC_COMM_SELF," \n(2) first %D columns of inv(A) with dense RHS:\n",nrhs);} 130*c4762a1bSJed Brown ierr = MatView(X,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 131*c4762a1bSJed Brown } 132*c4762a1bSJed Brown 133*c4762a1bSJed Brown /* Check the residual */ 134*c4762a1bSJed Brown ierr = MatMatMult(A,X,MAT_INITIAL_MATRIX,2.0,&AX);CHKERRQ(ierr); 135*c4762a1bSJed Brown ierr = MatAXPY(AX,-1.0,RHS,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 136*c4762a1bSJed Brown ierr = MatNorm(AX,NORM_INFINITY,&norm);CHKERRQ(ierr); 137*c4762a1bSJed Brown if (norm > tol) { 138*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"(2) MatMatSolve: Norm of residual %g\n",norm);CHKERRQ(ierr); 139*c4762a1bSJed Brown } 140*c4762a1bSJed Brown ierr = MatZeroEntries(X);CHKERRQ(ierr); 141*c4762a1bSJed Brown 142*c4762a1bSJed Brown /* (3) Test MatMatTransposeSolve() for inv(A) with sparse RHS stored in the host: 143*c4762a1bSJed Brown spRHST = [e[0],...,e[nrhs-1]]^T, dense X holds first nrhs columns of inv(A) */ 144*c4762a1bSJed Brown /* --------------------------------------------------------------------------- */ 145*c4762a1bSJed Brown /* Create spRHST: PETSc does not support compressed column format which is required by MUMPS for sparse RHS matrix, 146*c4762a1bSJed Brown thus user must create spRHST=spRHS^T and call MatMatTransposeSolve() */ 147*c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD,&spRHST);CHKERRQ(ierr); 148*c4762a1bSJed Brown if (!rank) { 149*c4762a1bSJed Brown /* MUMPS requirs RHS be centralized on the host! */ 150*c4762a1bSJed Brown ierr = MatSetSizes(spRHST,nrhs,M,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr); 151*c4762a1bSJed Brown } else { 152*c4762a1bSJed Brown ierr = MatSetSizes(spRHST,0,0,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr); 153*c4762a1bSJed Brown } 154*c4762a1bSJed Brown ierr = MatSetType(spRHST,MATAIJ);CHKERRQ(ierr); 155*c4762a1bSJed Brown ierr = MatSetFromOptions(spRHST);CHKERRQ(ierr); 156*c4762a1bSJed Brown ierr = MatSetUp(spRHST);CHKERRQ(ierr); 157*c4762a1bSJed Brown if (!rank) { 158*c4762a1bSJed Brown v = 1.0; 159*c4762a1bSJed Brown for (i=0; i<nrhs; i++) { 160*c4762a1bSJed Brown ierr = MatSetValues(spRHST,1,&i,1,&i,&v,INSERT_VALUES);CHKERRQ(ierr); 161*c4762a1bSJed Brown } 162*c4762a1bSJed Brown } 163*c4762a1bSJed Brown ierr = MatAssemblyBegin(spRHST,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 164*c4762a1bSJed Brown ierr = MatAssemblyEnd(spRHST,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 165*c4762a1bSJed Brown 166*c4762a1bSJed Brown ierr = MatMatTransposeSolve(F,spRHST,X);CHKERRQ(ierr); 167*c4762a1bSJed Brown 168*c4762a1bSJed Brown if (displ) { 169*c4762a1bSJed Brown if (!rank) {ierr = PetscPrintf(PETSC_COMM_SELF," \n(3) first %D columns of inv(A) with sparse RHS:\n",nrhs);} 170*c4762a1bSJed Brown ierr = MatView(X,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 171*c4762a1bSJed Brown } 172*c4762a1bSJed Brown 173*c4762a1bSJed Brown /* Check the residual: R = A*X - RHS */ 174*c4762a1bSJed Brown ierr = MatMatMult(A,X,MAT_REUSE_MATRIX,2.0,&AX);CHKERRQ(ierr); 175*c4762a1bSJed Brown 176*c4762a1bSJed Brown ierr = MatAXPY(AX,-1.0,RHS,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 177*c4762a1bSJed Brown ierr = MatNorm(AX,NORM_INFINITY,&norm);CHKERRQ(ierr); 178*c4762a1bSJed Brown if (norm > tol) { 179*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"(3) MatMatSolve: Norm of residual %g\n",norm);CHKERRQ(ierr); 180*c4762a1bSJed Brown } 181*c4762a1bSJed Brown 182*c4762a1bSJed Brown /* (4) Test MatMatSolve() for inv(A) with selected entries: 183*c4762a1bSJed Brown input: spRHS gives selected indices; output: spRHS holds selected entries of inv(A) */ 184*c4762a1bSJed Brown /* --------------------------------------------------------------------------------- */ 185*c4762a1bSJed Brown if (nrhs == N) { /* mumps requires nrhs = n */ 186*c4762a1bSJed Brown /* Create spRHS on proc[0] */ 187*c4762a1bSJed Brown Mat spRHS = NULL; 188*c4762a1bSJed Brown 189*c4762a1bSJed Brown /* Create spRHS = spRHST^T. Two matrices share internal matrix data structure */ 190*c4762a1bSJed Brown ierr = MatCreateTranspose(spRHST,&spRHS);CHKERRQ(ierr); 191*c4762a1bSJed Brown ierr = MatMumpsGetInverse(F,spRHS);CHKERRQ(ierr); 192*c4762a1bSJed Brown ierr = MatDestroy(&spRHS);CHKERRQ(ierr); 193*c4762a1bSJed Brown 194*c4762a1bSJed Brown ierr = MatMumpsGetInverseTranspose(F,spRHST);CHKERRQ(ierr); 195*c4762a1bSJed Brown if (displ) { 196*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"\nSelected entries of inv(A^T):\n");CHKERRQ(ierr); 197*c4762a1bSJed Brown ierr = MatView(spRHST,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 198*c4762a1bSJed Brown } 199*c4762a1bSJed Brown ierr = MatDestroy(&spRHS);CHKERRQ(ierr); 200*c4762a1bSJed Brown } 201*c4762a1bSJed Brown ierr = MatDestroy(&AX);CHKERRQ(ierr); 202*c4762a1bSJed Brown ierr = MatDestroy(&F);CHKERRQ(ierr); 203*c4762a1bSJed Brown ierr = MatDestroy(&RHS);CHKERRQ(ierr); 204*c4762a1bSJed Brown ierr = MatDestroy(&spRHST);CHKERRQ(ierr); 205*c4762a1bSJed Brown } 206*c4762a1bSJed Brown 207*c4762a1bSJed Brown /* Free data structures */ 208*c4762a1bSJed Brown ierr = MatDestroy(&A);CHKERRQ(ierr); 209*c4762a1bSJed Brown ierr = MatDestroy(&C);CHKERRQ(ierr); 210*c4762a1bSJed Brown ierr = MatDestroy(&X);CHKERRQ(ierr); 211*c4762a1bSJed Brown ierr = PetscRandomDestroy(&rand);CHKERRQ(ierr); 212*c4762a1bSJed Brown ierr = PetscFinalize(); 213*c4762a1bSJed Brown return ierr; 214*c4762a1bSJed Brown #endif 215*c4762a1bSJed Brown } 216*c4762a1bSJed Brown 217*c4762a1bSJed Brown /*TEST 218*c4762a1bSJed Brown 219*c4762a1bSJed Brown 220*c4762a1bSJed Brown test: 221*c4762a1bSJed Brown requires: mumps double !complex 222*c4762a1bSJed Brown 223*c4762a1bSJed Brown test: 224*c4762a1bSJed Brown suffix: 2 225*c4762a1bSJed Brown requires: mumps double !complex 226*c4762a1bSJed Brown nsize: 2 227*c4762a1bSJed Brown 228*c4762a1bSJed Brown TEST*/ 229