[petsc-users] PetscFE and TS
Matthew Knepley
knepley at gmail.com
Wed Nov 9 07:33:09 CST 2016
On Wed, Nov 9, 2016 at 12:32 AM, Julian Andrej <juan at tf.uni-kiel.de> wrote:
> I tested a few other things and it turned out that the function added
> as a dirichlet boundary condition via PetscDSAddBoundary doesn't
> receive the time values correctly (PetscReal time is always 0.0).
>
I will check out your example and fix this. However, I wonder if there is a
bug in your implementation:
void f0_temp(PetscInt dim, PetscInt Nf, PetscInt NfAux,
const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar
u[], const PetscScalar u_t[], const PetscScalar u_x[],
const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar
a[], const PetscScalar a_t[], const PetscScalar a_x[],
PetscReal t, const PetscReal x[], PetscScalar f0[])
{
f0[0] = u_t[0] - 2.0;
}
I get
du/dt - \Delta u = -2
so that if u = x^2 + y^2 + 2t, then
2 - (2 + 2) = -2
so you would want u_t + 2.
If you give the Laplacian the opposite sign, this is an unstable
formulation.
Thanks,
Matt
> I also see no effect adding DMTSSetIJacobianLocal vs. not setting the
> jacobian function. The analytic solution i used in my attached example
> is the same as in the src/ts/examples/tutorials/ex32.c.
>
> Am i doing anything obviously wrong here? My next step would be to try
> if assembling the matrices myself and writing a custom
> IFunction/IJacobian which assembles the different parts of the matrix
> like M and J for the stiff ODE with nontrivial mass matrix (see manual
> section of TS -> F = M u' - f) but i think this should be obsolete
> right?
>
> Appreciate your help.
>
> Regards
> Julian
>
> On Mon, Nov 7, 2016 at 10:49 AM, Julian Andrej <juan at tf.uni-kiel.de>
> wrote:
> > Hello,
> >
> > i'm using PetscFE in combination with SNES and a hand written backward
> > euler solver right now and like to hand that over to TS. I have been
> > fiddling quite a while getting TS to work with PetscFE and have
> > encountered a few unclarities.
> >
> > I have looked at examples which are outdated i guess
> > (src/ts/examples/tutorials/ex32.c) and i am confused on how i have to
> > formulate the discretization of the jacobian and the residual. It
> > seems obvious to use the DMTSSetIFunctionLocal and
> > DMTSSetIJacobianLocal functions because there are prepared wrappers
> > from DMPlexTSComputeIFunctionFEM and DMPlexTSComputeIJacobianFEM.
> >
> > I need a Mass matrix for my problem, is that formed from the function
> > space information or do i have to form it myself? Is there any working
> > example which uses PetscFE and TS to form the DMTSSetIFunctionLocal
> > and DMTSSetIJacobianLocal? (grepping the tree tells me "no"
> > unfortunately).
> >
> > Regards
> > Julian
>
--
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/20161109/648df980/attachment.html>
More information about the petsc-users
mailing list