static char help[] = "Tests MatMatSolve() and MatMatTransposeSolve() for computing inv(A) with MUMPS.\n\ Example: mpiexec -n ./ex214 -displ \n\n"; #include int main(int argc,char **args) { PetscErrorCode ierr; PetscMPIInt size,rank; #if defined(PETSC_HAVE_MUMPS) Mat A,RHS,C,F,X,AX,spRHST; PetscInt m,n,nrhs,M,N,i,Istart,Iend,Ii,j,J,test; PetscScalar v; PetscReal norm,tol=PETSC_SQRT_MACHINE_EPSILON; PetscRandom rand; PetscBool displ=PETSC_FALSE; char solver[256]; #endif ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr; CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); #if !defined(PETSC_HAVE_MUMPS) if (rank == 0) CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"This example requires MUMPS, exit...\n")); ierr = PetscFinalize(); return ierr; #else CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-displ",&displ,NULL)); /* Create matrix A */ m = 4; n = 4; CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL)); CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL)); CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A)); CHKERRQ(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n)); CHKERRQ(MatSetFromOptions(A)); CHKERRQ(MatMPIAIJSetPreallocation(A,5,NULL,5,NULL)); CHKERRQ(MatSeqAIJSetPreallocation(A,5,NULL)); CHKERRQ(MatGetOwnershipRange(A,&Istart,&Iend)); for (Ii=Istart; Ii0) {J = Ii - n; CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&v,ADD_VALUES));} if (i0) {J = Ii - 1; CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&v,ADD_VALUES));} if (j tol) { CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"(1) MatMatSolve: Norm of error %g\n",norm)); } /* Test X=RHS */ CHKERRQ(MatMatSolve(F,RHS,RHS)); /* Check the error */ CHKERRQ(MatAXPY(RHS,-1.0,C,SAME_NONZERO_PATTERN)); CHKERRQ(MatNorm(RHS,NORM_FROBENIUS,&norm)); if (norm > tol) { CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"(1.1) MatMatSolve: Norm of error %g\n",norm)); } /* (2) Test MatMatSolve() for inv(A) with dense RHS: RHS = [e[0],...,e[nrhs-1]], dense X holds first nrhs columns of inv(A) */ /* -------------------------------------------------------------------- */ CHKERRQ(MatZeroEntries(RHS)); for (i=0; i tol) { CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"(2) MatMatSolve: Norm of residual %g\n",norm)); } CHKERRQ(MatZeroEntries(X)); /* (3) Test MatMatTransposeSolve() for inv(A) with sparse RHS stored in the host: spRHST = [e[0],...,e[nrhs-1]]^T, dense X holds first nrhs columns of inv(A) */ /* --------------------------------------------------------------------------- */ /* Create spRHST: PETSc does not support compressed column format which is required by MUMPS for sparse RHS matrix, thus user must create spRHST=spRHS^T and call MatMatTransposeSolve() */ CHKERRQ(MatCreate(PETSC_COMM_WORLD,&spRHST)); if (rank == 0) { /* MUMPS requirs RHS be centralized on the host! */ CHKERRQ(MatSetSizes(spRHST,nrhs,M,PETSC_DECIDE,PETSC_DECIDE)); } else { CHKERRQ(MatSetSizes(spRHST,0,0,PETSC_DECIDE,PETSC_DECIDE)); } CHKERRQ(MatSetType(spRHST,MATAIJ)); CHKERRQ(MatSetFromOptions(spRHST)); CHKERRQ(MatSetUp(spRHST)); if (rank == 0) { v = 1.0; for (i=0; i tol) { CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"(3) MatMatSolve: Norm of residual %g\n",norm)); } /* (4) Test MatMatSolve() for inv(A) with selected entries: input: spRHS gives selected indices; output: spRHS holds selected entries of inv(A) */ /* --------------------------------------------------------------------------------- */ if (nrhs == N) { /* mumps requires nrhs = n */ /* Create spRHS on proc[0] */ Mat spRHS = NULL; /* Create spRHS = spRHST^T. Two matrices share internal matrix data structure */ CHKERRQ(MatCreateTranspose(spRHST,&spRHS)); CHKERRQ(MatMumpsGetInverse(F,spRHS)); CHKERRQ(MatDestroy(&spRHS)); CHKERRQ(MatMumpsGetInverseTranspose(F,spRHST)); if (displ) { CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"\nSelected entries of inv(A^T):\n")); CHKERRQ(MatView(spRHST,PETSC_VIEWER_STDOUT_WORLD)); } CHKERRQ(MatDestroy(&spRHS)); } CHKERRQ(MatDestroy(&AX)); CHKERRQ(MatDestroy(&F)); CHKERRQ(MatDestroy(&RHS)); CHKERRQ(MatDestroy(&spRHST)); } /* Free data structures */ CHKERRQ(MatDestroy(&A)); CHKERRQ(MatDestroy(&C)); CHKERRQ(MatDestroy(&X)); CHKERRQ(PetscRandomDestroy(&rand)); ierr = PetscFinalize(); return ierr; #endif } /*TEST test: requires: mumps double !complex test: suffix: 2 requires: mumps double !complex nsize: 2 TEST*/