[petsc-users] Coloring of a finite volume unstructured mesh

Smith, Barry F. bsmith at mcs.anl.gov
Mon Nov 13 06:04:07 CST 2017



> On Nov 13, 2017, at 2:10 AM, SIERRA-AUSIN Javier <Javier.SIERRA-AUSIN at isae-supaero.fr> wrote:
> 
> Hi thanks for your answer,
> 
> I would like to precise that in my particular case I deal with an unstructured grid with an stencil that takes two distance neighbour (center of the cell). And indeed, you are right the coupling is between faces but since it is not order 1 (upwind), the reconstruction of the face state takes several (2 in my case) values at the cell of neighbors.
> So that in 2d if we perturbe a cell we would perturbe 13 around (distance 2 neighbours) and in 3d, 25 neighbors.
> My question is how I can provide PETSc that structure of the nonzero structure of the Jacobian? With AO?

  No, AO has nothing to do with it.

>  Would it be like src/snes/examples/tutorials/ex10d/ex10.c, but having an adjacency matrix with the stencil connectivity instead of the geometrical/parmetis connectivity?

   Yes, pretty much. Note that example uses MatSetValuesLocal() but most likely you will use MatSetValues()

    So basically you need to preallocate the matrix with enough space and then call MatSetValues() to indicate all possible nonzero entries in the Jacobian due to the stencil.




> If not I would really appreciate an example illustrating how I should do it.
> 
> Thanks in advance!!
> 
> Javier
> 
> 
> 
> On Friday, November 10, 2017 05:09 CET, "Smith, Barry F." <bsmith at mcs.anl.gov> wrote:
>  
>> 
>> To use the PETSc coloring based Jacobian computer (which uses finite differences) you absolutely have to be able to provide the nonzero structure of the Jacobian.
>> 
>> Now once you provide the nonzero structure of the Jacobian the PETSc MatColoring routines can actually compute the coloring for you.
>> 
>> So in other words you need not worry about the coloring at all, you just need to worry about providing the nonzero structure. Since you are using a finite volume method presumably all your coupling is between faces? In this case explicitly computing the nonzero structure of the Jacobian is probably pretty straightforward and you should just do it.
>> 
>> Barry
>> 
>> 
>> > On Nov 7, 2017, at 10:13 AM, Matthew Knepley <knepley at gmail.com> wrote:
>> >
>> > On Tue, Nov 7, 2017 at 10:18 AM, SIERRA-AUSIN Javier <Javier.SIERRA-AUSIN at isae-supaero.fr> wrote:
>> > Hi,
>> >
>> > I would like to ask you concerning the computation of the Jacobian matrix via finite difference and coloring of the connectivity graph.
>> > I wonder whether it is possible or not to color the Jacobian matrix of a given solver that evaluates the RHS with its associated connectivity in the global indeces of my solver (not PETSc).
>> > As well, if it is possible to do this from an already partioned domain in parallel.
>> > All of this is better explained in this post : https://scicomp.stackexchange.com/questions/28209/linking-petsc-with-an-already-parallel-in-house-finite-volume-solver
>> >
>> > The simplest thing you can do is to use the finite-difference Jacobian action (MatMFFD). This is setup automatically by SNES
>> > if you give a FormFunction pointer, but no FormJacobian routine. Just tell the PETSc Vecs to use your ParMetis layout (by
>> > setting the local sizes), and it should run fine in SNES.
>> >
>> > However, usually you need some kind of preconditioning. Thus you either have to form the Jacobian or some approximation. If
>> > you cannot form an approximation, then you can use coloring. Once option is to create a DMPlex with your mesh information.
>> > This can be done in parallel after you have already partitioned with ParMetis (as long as you know the "overlap" of vertices, or
>> > adjacency of cells). Then the coloring can be done automatically using that DM information. Otherwise, you will have to supply
>> > a coloring to the SNES.
>> >
>> > Thanks,
>> >
>> > Matt
>> >
>> > Thanks in advance,
>> >
>> > Javier.
>> >
>> >
>> >
>> >
>> > --
>> > 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
>> >
>> > https://www.cse.buffalo.edu/~knepley/
>>  
> 
> 
> 
>  



More information about the petsc-users mailing list