[petsc-users] TS for incompressible NS equations

Gianluca Meneghello gianmail at gmail.com
Wed Sep 28 10:37:30 CDT 2011


Jed,

thanks for your reply. I am not sure I understand everything...

I do have a constant Jacobian that never changes. What I understand is
the matrix I should pass to TSSetIJacobian looks like

J + a B

together with TSComputeIJacobianConstant. Something like


TSSetIFunction(ts,PETSC_NULL,TSComputeIFunctionLinear,PETSC_NULL);
TSSetIJacobian(ts,A,A,TSComputeIJacobianConstant,PETSC_NULL);

where A = J + a B

but then I don't understand how I can fix a before deciding the time
integration method or even the time step, so I guess I'm wrong.

Is there a way I can use TSComputeIJacobianConstant or am I wrong
about trying to use it?

Thanks

Gianluca







On 28 September 2011 15:25, Jed Brown <jedbrown at mcs.anl.gov> wrote:
> On Wed, Sep 28, 2011 at 08:02, Gianluca Meneghello <gianmail at gmail.com>
> wrote:
>>
>> I need to solve the linearized, incompressible Navier Stokes equation
>> in time. I have reformulated the problem as
>>
>> B du/dt = J u
>>
>> where J is the Jacobian of the Navier Stokes equations and u is a
>> vector containing {u,v,,w,p}.
>> B is a diagonal matrix which is 0 in the lines corresponding to the
>> pressure equation (enforcing the divergence free condition) and
>> Dirichlet BC's and 1 otherwise.
>>
>> Backward Euler would then solve it as
>>
>>  u(t+1) = (B-dt J)^-1 B u(t)
>>
>> as in page 115 of the user guide (is there a typing error there?
>> Should it be " (B - dt A) u(n+1) = B u(n) " ?)
>
> Yes, thank you. I should rewrite this section because it's confusing in its
> current form. It would be better to look at the later sections 6.1.2 Solving
> Differential Algebraic Equations and 6.1.3 on IMEX methods.
>
>>
>> Using PETSc 3.2, I pass J to the TS object as
>>
>>
>> TSSetRHSFunction(ts,PETSC_NULL,TSComputeRHSFunctionLinear,PETSC_NULL);CHKERRQ(ierr);
>>
>> TSSetRHSJacobian(ts,A,A,TSComputeRHSJacobianConstant,PETSC_NULL);CHKERRQ(ierr);
>
> This sets a constant Jacobian that never changes. I don't know if that's
> what you mean by linearized Navier-Stokes.
>
>>
>> but I don't understand how I can set B (or whether I am using a wrong
>> approach).
>
> If you have a B that is not the identity, then you should use
> TSSetIFunction() and/or TSSetIJacobian(). Note that incompressible flow in
> velocity-pressure form is an index 2 DAE, so not all coupled time
> integration methods perform well. There is very new code in petsc-dev that
> you may want to try for this problem. You have to enable these by
> configuring --with-rosw and then you can run with -ts_type rosw
> -ts_rosw_type ra34pw2. (This would use a four stage third order Rosenbrock-W
> scheme for stiff problems that performed well in the Navier-Stokes test
> problems of Rang & Angermann 2005, -ts_rosw_type ra3pw is a three-stage
> third order method that is also worth trying.)
> The need to configure --with-rosw will go away in a week or two, once that
> code stabilizes.



-- 
"[Je pense que] l'homme est un monde qui vaut des fois les mondes et
que les plus ardentes ambitions sont celles qui ont eu l'orgueil de
l'Anonymat" -- Non omnibus, sed mihi et tibi
Amedeo Modigliani


More information about the petsc-users mailing list