[petsc-users] PetscFV and TS implicit

Matthew Knepley knepley at gmail.com
Thu Aug 20 11:33:44 CDT 2020


On Thu, Aug 20, 2020 at 7:45 AM Thibault Bridel-Bertomeu <
thibault.bridelbertomeu at gmail.com> wrote:

> Dear all,
>
> I have a finite volume code inspired from TS example ex11.c (with a
> riemann solver, etc...).
> So far, I used only explicit time stepping through the TSSSP, and to set
> the RHS of my hyperbolic system I used :
>
> TSSetType(ts, TSSSP);
> DMTSSetRHSFunctionLocal(dm, DMPlexTSComputeRHSFunctionFVM <https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/TS/DMPlexTSComputeRHSFunctionFVM.html#DMPlexTSComputeRHSFunctionFVM>, &ctx);
>
> after setting the right Riemann solver in the TS associated to the DM.
> Now, in some cases where the physics is stationary, I would like to reach
> the steady state faster and use an implicit timestepping operator to do so.
> After looking through the examples, especially TS examples ex48.c and
> ex53.c, I simply tried to set
>
> TSSetType(ts, TSBEULER);
> DMTSSetIFunctionLocal(dm, DMPlexTSComputeIFunctionFEM <https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/TS/DMPlexTSComputeIFunctionFEM.html#DMPlexTSComputeIFunctionFEM>, &ctx);
> DMTSSetIJacobianLocal(dm, DMPlexTSComputeIJacobianFEM <https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/TS/DMPlexTSComputeIJacobianFEM.html#DMPlexTSComputeIJacobianFEM>, &ctx);
>
> instead of the previous calls. It compiles fine, and it runs. However,
> nothing happens : it behaves like there is no time evolution at all, the
> solution does not change from its initial state.
> From the source code, it is my understanding that the
> DMPlexTSComputeIFunctionFEM and DMPlexTSComputeIJacobianFEM methods, in
> spite of their names, call respectively DMPlexComputeResidual_Internal and
> DMPlexComputeJacobian_Internal that can handle a FVM discretization.
>
> What am I missing ? Are there other steps to take before I can simply try
> to run with a finite volume discretization and an implicit time stepping
> algorithm ?
>
> Thank you very much for your help in advance !
>

I could never get the FVM stuff to make sense to me for implicit methods.
Here is my problem understanding. If you have an FVM method, it decides
to move "stuff" from one cell to its neighboring cells depending on the
solution to the Riemann problem on each face, which computed the flux. This
is
fine unless the timestep is so big that material can flow through into the
cells beyond the neighbor. Then I should have considered the effect of the
Riemann problem for those interfaces. That would be in the Jacobian, but I
don't know how to compute that Jacobian. I guess you could do everything
matrix-free, but without a preconditioner it seems hard.

Operationally, I always mark FVM things explicit, so they are only handled
by RHSFunction/Jacobian, except for the u_t term which I stick in the
IJacobian.
I was never successful in running an implicit FVM and still do not
understand it. I could not find any literature that I understood either.

If there is a definite thing you want changed, I might be able to do it.

  Thanks,

     Matt


> Thibault Bridel-Bertomeu
>> Eng, MSc, PhD
> Research Engineer
> CEA/CESTA
> 33114 LE BARP
> Tel.: (+33)557046924
> Mob.: (+33)611025322
> Mail: thibault.bridelbertomeu at gmail.com
>


-- 
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/ <http://www.cse.buffalo.edu/~knepley/>
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20200820/d5dea3c4/attachment-0001.html>


More information about the petsc-users mailing list