/* Developer Notes: This should be implemented with a MatCreate_SchurComplement() as that is the standard design for new Mat classes. */ PetscErrorCode MatGetSchurComplement_Basic(Mat mat,IS isrow0,IS iscol0,IS isrow1,IS iscol1,MatReuse mreuse,Mat *newmat,MatReuse preuse,Mat *newpmat) { PetscErrorCode ierr; Mat A=0,Ap=0,B=0,C=0,D=0; PetscFunctionBegin; PetscValidHeaderSpecific(mat,MAT_CLASSID,1); PetscValidHeaderSpecific(isrow0,IS_CLASSID,2); PetscValidHeaderSpecific(iscol0,IS_CLASSID,3); PetscValidHeaderSpecific(isrow1,IS_CLASSID,4); PetscValidHeaderSpecific(iscol1,IS_CLASSID,5); if (mreuse == MAT_REUSE_MATRIX) PetscValidHeaderSpecific(*newmat,MAT_CLASSID,7); if (preuse == MAT_REUSE_MATRIX) PetscValidHeaderSpecific(*newpmat,MAT_CLASSID,9); PetscValidType(mat,1); if (mat->factortype) SETERRQ(((PetscObject)mat)->comm,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); if (mreuse != MAT_IGNORE_MATRIX) { /* Use MatSchurComplement */ if (mreuse == MAT_REUSE_MATRIX) { ierr = MatSchurComplementGetSubmatrices(*newmat,&A,&Ap,&B,&C,&D);CHKERRQ(ierr); if (!A || !Ap || !B || !C) SETERRQ(((PetscObject)mat)->comm,PETSC_ERR_ARG_WRONGSTATE,"Attempting to reuse matrix but Schur complement matrices unset"); if (A != Ap) SETERRQ(((PetscObject)mat)->comm,PETSC_ERR_ARG_WRONGSTATE,"Preconditioning matrix does not match operator"); ierr = MatDestroy(&Ap);CHKERRQ(ierr); /* get rid of extra reference */ } ierr = MatGetSubMatrix(mat,isrow0,iscol0,mreuse,&A);CHKERRQ(ierr); ierr = MatGetSubMatrix(mat,isrow0,iscol1,mreuse,&B);CHKERRQ(ierr); ierr = MatGetSubMatrix(mat,isrow1,iscol0,mreuse,&C);CHKERRQ(ierr); ierr = MatGetSubMatrix(mat,isrow1,iscol1,mreuse,&D);CHKERRQ(ierr); switch (mreuse) { case MAT_INITIAL_MATRIX: ierr = MatCreateSchurComplement(A,A,B,C,D,newmat);CHKERRQ(ierr); break; case MAT_REUSE_MATRIX: ierr = MatSchurComplementUpdate(*newmat,A,A,B,C,D,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); break; default: SETERRQ(((PetscObject)mat)->comm,PETSC_ERR_SUP,"Unrecognized value of mreuse"); } } if (preuse != MAT_IGNORE_MATRIX) { /* Use the diagonal part of A to form D - C inv(diag(A)) B */ Mat Ad,AdB,S; Vec diag; PetscInt i,m,n,mstart,mend; PetscScalar *x; /* We could compose these with newpmat so that the matrices can be reused. */ if (!A) {ierr = MatGetSubMatrix(mat,isrow0,iscol0,MAT_INITIAL_MATRIX,&A);CHKERRQ(ierr);} if (!B) {ierr = MatGetSubMatrix(mat,isrow0,iscol1,MAT_INITIAL_MATRIX,&B);CHKERRQ(ierr);} if (!C) {ierr = MatGetSubMatrix(mat,isrow1,iscol0,MAT_INITIAL_MATRIX,&C);CHKERRQ(ierr);} if (!D) {ierr = MatGetSubMatrix(mat,isrow1,iscol1,MAT_INITIAL_MATRIX,&D);CHKERRQ(ierr);} ierr = MatGetVecs(A,&diag,PETSC_NULL);CHKERRQ(ierr); ierr = MatGetDiagonal(A,diag);CHKERRQ(ierr); ierr = VecReciprocal(diag);CHKERRQ(ierr); ierr = MatGetLocalSize(A,&m,&n);CHKERRQ(ierr); /* We need to compute S = D - C inv(diag(A)) B. For row-oriented formats, it is easy to scale the rows of B and * for column-oriented formats the columns of C can be scaled. Would skip creating a silly diagonal matrix. */ ierr = MatCreate(((PetscObject)A)->comm,&Ad);CHKERRQ(ierr); ierr = MatSetSizes(Ad,m,n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); ierr = MatSetOptionsPrefix(Ad,((PetscObject)mat)->prefix);CHKERRQ(ierr); ierr = MatAppendOptionsPrefix(Ad,"diag_");CHKERRQ(ierr); ierr = MatSetFromOptions(Ad);CHKERRQ(ierr); ierr = MatSeqAIJSetPreallocation(Ad,1,PETSC_NULL);CHKERRQ(ierr); ierr = MatMPIAIJSetPreallocation(Ad,1,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr); ierr = MatGetOwnershipRange(Ad,&mstart,&mend);CHKERRQ(ierr); ierr = VecGetArray(diag,&x);CHKERRQ(ierr); for (i=mstart; i<mend; i++) { ierr = MatSetValue(Ad,i,i,x[i-mstart],INSERT_VALUES);CHKERRQ(ierr); } ierr = VecRestoreArray(diag,&x);CHKERRQ(ierr); ierr = MatAssemblyBegin(Ad,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); ierr = MatAssemblyEnd(Ad,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); ierr = VecDestroy(&diag);CHKERRQ(ierr); ierr = MatMatMult(Ad,B,MAT_INITIAL_MATRIX,1,&AdB);CHKERRQ(ierr); S = (preuse == MAT_REUSE_MATRIX) ? *newpmat : (Mat)0; ierr = MatMatMult(C,AdB,preuse,PETSC_DEFAULT,&S);CHKERRQ(ierr); ierr = MatAYPX(S,-1,D,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); *newpmat = S; ierr = MatDestroy(&Ad);CHKERRQ(ierr); ierr = MatDestroy(&AdB);CHKERRQ(ierr); } ierr = MatDestroy(&A);CHKERRQ(ierr); ierr = MatDestroy(&B);CHKERRQ(ierr); ierr = MatDestroy(&C);CHKERRQ(ierr); ierr = MatDestroy(&D);CHKERRQ(ierr); PetscFunctionReturn(0); }
static PetscErrorCode PCSetUp_PCD_Feelpp(PC pc) { PetscErrorCode ierr; PC_PCD_Feelpp *pcpcd = (PC_PCD_Feelpp*)pc->data; if ( !pcpcd->allocated ) { ierr = KSPCreate(PetscObjectComm((PetscObject)pc),&pcpcd->kspAp);CHKERRQ(ierr); ierr = PetscObjectIncrementTabLevel((PetscObject)pcpcd->kspAp,(PetscObject)pc,1);CHKERRQ(ierr); ierr = KSPSetType(pcpcd->kspAp,KSPPREONLY);CHKERRQ(ierr); ierr = KSPSetOptionsPrefix(pcpcd->kspAp,((PetscObject)pc)->prefix);CHKERRQ(ierr); ierr = KSPAppendOptionsPrefix(pcpcd->kspAp,"pcd_Ap");CHKERRQ(ierr); ierr = KSPSetFromOptions(pcpcd->kspAp);CHKERRQ(ierr); ierr = KSPCreate(PetscObjectComm((PetscObject)pc),&pcpcd->kspMp);CHKERRQ(ierr); ierr = PetscObjectIncrementTabLevel((PetscObject)pcpcd->kspMp,(PetscObject)pc,1);CHKERRQ(ierr); ierr = KSPSetType(pcpcd->kspMp,KSPPREONLY);CHKERRQ(ierr); ierr = KSPSetOptionsPrefix(pcpcd->kspMp,((PetscObject)pc)->prefix);CHKERRQ(ierr); ierr = KSPAppendOptionsPrefix(pcpcd->kspMp,"pcd_Mp");CHKERRQ(ierr); ierr = KSPSetFromOptions(pcpcd->kspMp);CHKERRQ(ierr); #if PETSC_VERSION_LESS_THAN(3,6,0) ierr = MatGetVecs(pcpcd->matFp,&pcpcd->x2,&pcpcd->x1);CHKERRQ(ierr); #else ierr = MatCreateVecs(pcpcd->matFp,&pcpcd->x2,&pcpcd->x1);CHKERRQ(ierr); #endif pcpcd->allocated = PETSC_TRUE; } PetscBool isBTBt; PetscStrcmp(pcpcd->pcdApType,"BTBt",&isBTBt); if ( isBTBt ) { if (!pcpcd->MvDiag ) { #if PETSC_VERSION_LESS_THAN(3,6,0) ierr = MatGetVecs(pcpcd->matMv,&pcpcd->MvDiag,NULL);CHKERRQ(ierr); #else ierr = MatCreateVecs(pcpcd->matMv,&pcpcd->MvDiag,NULL);CHKERRQ(ierr); #endif } Mat B, C; #if PETSC_VERSION_GREATER_OR_EQUAL_THAN( 3,5,0 ) ierr = MatSchurComplementGetSubMatrices(pc->mat,NULL,NULL,&B,&C,NULL);CHKERRQ(ierr); #else ierr = MatSchurComplementGetSubmatrices(pc->mat,NULL,NULL,&B,&C,NULL);CHKERRQ(ierr); #endif ierr = MatGetDiagonal(pcpcd->matMv,pcpcd->MvDiag);CHKERRQ(ierr); ierr = VecReciprocal(pcpcd->MvDiag);CHKERRQ(ierr); // diag(F)^-1 * B ierr = MatDiagonalScale( B, pcpcd->MvDiag ,NULL);CHKERRQ(ierr); // C* diag(F)^-1 * B if ( !pcpcd->matApBTBt ) ierr = MatMatMult(C,B,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&pcpcd->matApBTBt); else ierr = MatMatMult(C,B,MAT_REUSE_MATRIX,PETSC_DEFAULT,&pcpcd->matApBTBt); CHKERRQ(ierr); ierr = VecReciprocal(pcpcd->MvDiag);CHKERRQ(ierr); ierr = MatDiagonalScale( B, pcpcd->MvDiag ,NULL);CHKERRQ(ierr); pcpcd->matAp = pcpcd->matApBTBt; } else pcpcd->matAp = pcpcd->matApLaplacian; #if PETSC_VERSION_GREATER_OR_EQUAL_THAN( 3,5,0 ) ierr = KSPSetOperators(pcpcd->kspAp,pcpcd->matAp,pcpcd->matAp);CHKERRQ(ierr); ierr = KSPSetOperators(pcpcd->kspMp,pcpcd->matMp,pcpcd->matMp);CHKERRQ(ierr); #else ierr = KSPSetOperators(pcpcd->kspAp,pcpcd->matAp,pcpcd->matAp,SAME_PRECONDITIONER);CHKERRQ(ierr); ierr = KSPSetOperators(pcpcd->kspMp,pcpcd->matMp,pcpcd->matMp,SAME_PRECONDITIONER);CHKERRQ(ierr); #endif PetscFunctionReturn(0); }