[petsc-users] Trying to set up a field-split preconditioner

Matthew Knepley knepley at gmail.com
Fri Jul 19 23:09:10 CDT 2013


On Fri, Jul 19, 2013 at 9:59 PM, subramanya sadasiva <potaman at outlook.com>wrote:

> Hi Matt,
> The DM being created is here (this is from Libmesh code
> (petscdmlibmesh.C )
>
>    01047 {
>    01048 PetscErrorCode ierr;
>    01049 PetscFunctionBegin;
>    01050 ierr = DMCreate(comm, dm); CHKERRQ(ierr);
>    01051 ierr = DMSetType(*dm, DMLIBMESH); CHKERRQ(ierr);
>    01052 ierr = DMLibMeshSetSystem(*dm, sys); CHKERRQ(ierr);
>    01053 PetscFunctionReturn(0);
>    01054 }
>
>
> This file has methods to access the  variables assigned to the DM (this
> seems to be stored in a struct.)
> So , I guess one should be able to add a bit of code to create sections as
> you mentioned somewhere around here.
>

Okay, they have their own DM. It must implement one of the interfaces for
field specification. They could provide


http://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DM/DMCreateFieldDecomposition.html

or at a lower level


http://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DM/DMCreateSubDM.html

which in turn can be constructed by specifying a default PetscSection


http://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DM/DMSetDefaultSection.html

    Matt


> Thanks,
> Subramanya
>
>
>
>
>
> Date: Fri, 19 Jul 2013 21:33:11 -0500
> Subject: Re: [petsc-users] Trying to set up a field-split preconditioner
> From: knepley at gmail.com
> To: potaman at outlook.com
> CC: petsc-users at mcs.anl.gov
>
> On Fri, Jul 19, 2013 at 9:17 PM, subramanya sadasiva <potaman at outlook.com>
> wrote:
> Hi Matt,
> I am using Libmesh so the DM stuff is actually in the background, and
> unfortunately the matrix doesn't have a saddle point,
> I thought that
>
>   -ch_solve_fieldsplit_block_size 2
>   -ch_solve_fieldsplit_0_fields 0
>   -ch_solve_fieldsplit_1_fields 1
>
> The block_size argument presumes you are using a DA. Are you?
>
> The other two options just say select the first DM field as field 0 in
> this PC, and the same with the second field. The
> DM must inform the PC about the initial field decomposition.
>
> would inform the solver of the structure. If this doesn't work owing to
> the fact that the problem is only being solved on a section of the mesh
> (because of the reduced space method), I guess I will have to use the
> PetscSection. Does that sound right?
>
> First, I think the right people to do this are the Libmesh people (we will
> of course help them). Second, you have not said
> whether you are using a structured or unstructured mesh. What DM class
> does the solver actually see?
>
>   Thanks,
>
>       Matt
>
> Thanks,
> Subramanya
>
>
> Subject: Re: [petsc-users] Trying to set up a field-split preconditioner
> From: knepley at gmail.com
> To: potaman at outlook.com
> CC: petsc-users at mcs.anl.gov; libmesh-users at lists.sourceforge.net
>
> On Fri, Jul 19, 2013 at 7:56 PM, subramanya sadasiva <potaman at outlook.com>
> wrote:
> Hi,
> I am trying to set up a fieldsplit preconditioner for my Cahn Hilliard
> solver and I get the following error,
>
> You have to tell the PCFIELDSPLIT about the dofs in each field. So
>
> 1) You are probably not using a DA, since it would tell it automatically
>
> 2) If you have a saddle point, you can use
> -pc_fieldsplit_detect_saddle_point
>
> 3) If none of those apply, you can set a PetscSection describing your
> layout to the DM for the solver.
>      Since this is new, I suspect you will need help, so mail back.
>
>   Thanks,
>
>      Matt
>
>
> [0]PETSC ERROR: --------------------- Error Message
> ------------------------------------
> [0]PETSC ERROR: Petsc has generated inconsistent data!
> [0]PETSC ERROR: Unhandled case, must have at least two fields, not 0!
> [0]PETSC ERROR:
> ------------------------------------------------------------------------
>
> These are the options that I am using,
> -ch_solve is just a prefix.
>
>
>
> -ch_solve_pc_type fieldsplit
> -ch_solve_pc_fieldsplit_type schur
> -ch_solve_fieldsplit_block_size 2
> -ch_solve_fieldsplit_0_field 1
> -ch_solve_fieldsplit_1_field 0
> -ch_solve_fieldsplit_0_ksp_type cg
> -ch_solve_fieldsplit_0_pc_type hypre
> -ch_solve_fieldsplit_0_pc_type_hypre boomeramg
> -ch_solve_fieldsplit_1_ksp_type cg
> -ch_solve_fieldsplit_1_pc_type hypre
> -ch_solve_fieldsplit_1_pc_type_hypre boomeramg
>
> Any ideas?
>
> Thanks,
> Subramanya
>
>
>
>
> --
> What most experimenters take for granted before they begin their
> experiments is infinitely more interesting than any results to which their
> experiments lead.
> -- Norbert Wiener
>
>
>
> --
> What most experimenters take for granted before they begin their
> experiments is infinitely more interesting than any results to which their
> experiments lead.
> -- Norbert Wiener
>



-- 
What most experimenters take for granted before they begin their
experiments is infinitely more interesting than any results to which their
experiments lead.
-- Norbert Wiener
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20130719/29785c37/attachment-0001.html>


More information about the petsc-users mailing list