[petsc-users] Guidelines for Nested fieldsplit for domains with a hybrid mesh
Barry Smith
bsmith at petsc.dev
Sun Oct 18 16:32:57 CDT 2020
From src/ksp/pc/impls/fieldsplit/fieldsplit.c this is how the saddle point is detected and set into the PC
if (jac->detect) {
IS zerodiags,rest;
PetscInt nmin,nmax;
ierr = MatGetOwnershipRange(pc->mat,&nmin,&nmax);CHKERRQ(ierr);
if (jac->diag_use_amat) {
ierr = MatFindZeroDiagonals(pc->mat,&zerodiags);CHKERRQ(ierr);
} else {
ierr = MatFindZeroDiagonals(pc->pmat,&zerodiags);CHKERRQ(ierr);
}
ierr = ISComplement(zerodiags,nmin,nmax,&rest);CHKERRQ(ierr);
ierr = PCFieldSplitSetIS(pc,"0",rest);CHKERRQ(ierr);
ierr = PCFieldSplitSetIS(pc,"1",zerodiags);CHKERRQ(ierr);
ierr = ISDestroy(&zerodiags);CHKERRQ(ierr);
ierr = ISDestroy(&rest);CHKERRQ(ierr);
In addition these two options are set
PetscErrorCode PCFieldSplitSetDetectSaddlePoint(PC pc,PetscBool flg)
{
PC_FieldSplit *jac = (PC_FieldSplit*)pc->data;
PetscErrorCode ierr;
PetscFunctionBegin;
jac->detect = flg;
if (jac->detect) {
ierr = PCFieldSplitSetType(pc,PC_COMPOSITE_SCHUR);CHKERRQ(ierr);
ierr = PCFieldSplitSetSchurPre(pc,PC_FIELDSPLIT_SCHUR_PRE_SELF,NULL);CHKERRQ(ierr);
}
PetscFunctionReturn(0);
}
You can use these routines to directly manage the IS yourself in any manner you choice.
Good luck
Barry
> On Oct 18, 2020, at 2:55 PM, Namala, Solomon <namala2 at illinois.edu> wrote:
>
> Hello,
>
> I am working to solve Stokes problem on a domain that is discretized using two different types of mesh. A part of the mesh uses fem formulation and the rest uses nodal integral method (NIM) formulation (the details of which I will skip). However, the key takeaway is that NIM formulation of stokes uses pressure Poisson formulation instead of the continuity equation while FEM formulation uses the continuity equation. They are coupled at the interface. Right now, I am building a single matrix for the entire domain and solving it using fieldsplit option in a nested fashion. The matrix structure and the unknown vector are shown below.
>
> My questions are:
>
> Are there any basic guidelines to solve these kind of problems.
> As I have mentioned I am currently using nested fieldsplit. The first split is using indices and the other split is done using detect saddle point option. is there a way to avoid using that option and doing it by combining set of indices or fields.
> The matrix structure is
> [Au_fem Bp_fem 0 0]
> [Cu_fem. 0 0 0]
> [0 0 Du_nim Ep_nim]
> [0 0 0 Fp_nim]
>
> the unknown vector is given by
>
> [ufem pfem unim pnim]
>
> Let me know if any additional information is needed.
>
> Thanks,
> Solomon.
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20201018/a5408a51/attachment.html>
More information about the petsc-users
mailing list