[petsc-dev] FAS with TS

Matthew Knepley knepley at gmail.com
Tue Feb 7 23:57:56 CST 2012


On Tue, Feb 7, 2012 at 11:50 PM, Jed Brown <jedbrown at mcs.anl.gov> wrote:

> On Wed, Feb 8, 2012 at 04:19, Peter Brune <prbrune at gmail.com> wrote:
>
>>
>>
>> On Tue, Feb 7, 2012 at 6:40 PM, Jed Brown <jedbrown at mcs.anl.gov> wrote:
>>
>>> Suppose I want to solve a finite time step problem using FAS. The time
>>> step will be defined something like
>>>
>>> static PetscErrorCode SNESTSFormFunction_ARKIMEX(SNES snes,Vec X,Vec
>>> F,TS ts)
>>> {
>>>   TS_ARKIMEX     *ark = (TS_ARKIMEX*)ts->data;
>>>   PetscErrorCode ierr;
>>>
>>>   PetscFunctionBegin;
>>>   ierr =
>>> VecAXPBYPCZ(ark->Ydot,-ark->shift,ark->shift,0,ark->Z,X);CHKERRQ(ierr); /*
>>> Ydot = shift*(X-Z) */
>>>   ierr =
>>> TSComputeIFunction(ts,ark->stage_time,X,ark->Ydot,F,ark->imex);CHKERRQ(ierr);
>>>   PetscFunctionReturn(0);
>>> }
>>>
>>>
>>> where ark->Z was determined from previous stages and ark->Ydot is a work
>>> vector. So FAS creates some hierarchy, but it isn't very forthcoming about
>>> showing the hierarchy to the user. But I need a way to
>>>
>>> 1. Restrict ark->Z to the coarse level and obtain the work vector
>>> ark->Ydot somewhere that I can reuse it. The user might also have to
>>> restrict "auxiliary" variables to coarse grids. It's conceivable that they
>>> would want to do this in a different way from the restriction/prolongation
>>> for the solution variables. For example, if you were using upwinded R and P
>>> for a hyperbolic problem, but you needed to restrict coefficients (like
>>> bathymetry). Or they might want to sample differently out of some file. So
>>> I think we need an auxiliary transfer callback that is at least called
>>> every top-level SNESSolve() and optionally also called in every transfer
>>> (because it could be homogenizing coefficient structure that is a nonlinear
>>> function of state).
>>>
>>>
>> The restriction (injection, etc) might be one step too hard at this
>> point.  We have to define the problem at each level in FAS and grid
>> sequencing so we might as well make it easy.  Is it possible to have the
>> callback for the SNES point to things that have to be prolonged or
>> restricted?  This could be called both before FAS and before grid
>> sequencing to get the coefficients up to the next level.  We could provide
>> some sort of default that goes through a number of user-provided vectors on
>> a level given vectors on another level, and inject them.  A user with
>> user-created multilevel coefficients could provide their own similar
>> callbacks.
>>
>
> So coefficients usually won't live on the same DM as the solution. This
> case with TS is easier in that regard, but it can't be the sole driver of
> design. There needs to be a general function to perform level setup.
>
^^^^^^^^^^^^^^^^^^^^^^^^

I just want to point out that Jed envisions that coefficients (and maybe
subproblems, etc) cannot be accommodated on the
same DM. I agree. However, this silly idea that we can make DMs all over
the place with no cost, like DAs, if they contain
all the mesh information, is just wrong. I think this is a good argument
for having both a topology object and a DM handling
layout/solver information. What is the counter-argument?

    Matt


>
> /* called by SNESSetUp_FAS() before restricting a nonlinear solution to a
> coarser level (but usually only used the first time unless this is doing
> solution-dependent homogenization) */
> typedef PetscErrorCode (*SNESFASRestrictHook)(SNES fine,Vec Xfine,SNES
> coarse);
> PetscErrorCode SNESFASSetRestrictHook(SNES fine,SNESFASRestrictHook
> hook,void *ctx);
>
> /* called in SNESSolve() each time a new solution is prolongated */
> typedef PetscErrorCode (*SNESGridSequenceInterpolateHook)(SNES coarse,Vec
> Xcoarse,SNES fine);
> PetscErrorCode SNESGridSequenceSetInterpolateHook(SNES
> coarse,SNESGridSequenceInterpolateHook hook,void *ctx);
>
>
> I think we also need to add state to DMCreateInterpolation() so that we
> can make it nonlinear when appropriate.
>
>
>>
>>> 2. Put the DM from this SNES into the TS that calls the user's
>>> IFunction. The easiest way is to have TSPushDM(TS,DM) that just jams it in
>>> there, caching the old DM (so the user's call to TSGetDM() works
>>> correctly), and TSPopDM(TS). This is ugly and is scary if the user does
>>> something weird like TSGetSolution() (returning the state at the beginning
>>> of the step, not the current value at the stage). This is something that
>>> doesn't make semantic sense except maybe if they have some weird
>>> diagnostics, but TSGetIJacobian() might be used for caching, and I'm scared
>>> of the semantics it would involve to support this sort of caching.
>>>
>>>
>> This option seems like it would be scary complex.  I don't like it.
>>
>
> I think it's likely simpler than the others, at least as a first
> implementation. I'm not wild about it either, but I'm not sure the
> alternatives are better.
>
>
>>
>>
>>> The alternative is to make different TSs for each level and somehow
>>> locate the correct TS (maybe cache it in the SNES hierarchy). I think this
>>> could be quite confusing to restrict all the members to the coarse grid.
>>> But having a multi-level TS object may eventually be sensible because
>>> multilevel time integration schemes exist. Spectral deferred correction is
>>> one example where coarse grids (in space and time) can be used to
>>> accelerate convergence to a high-order solution. But these somewhat
>>> naturally expose additional parallelism in time (usually more
>>> communication, but also more concurrency), so it likely needs a more
>>> complete solution that includes ways to instantiate entire models on
>>> sub-communicators (we'll also want an interface of this type for UQ).
>>>
>>>
>> I'm not sure this would be harmonious with FAS given the coarse problem
>> is redefined with a defect from the fine problem (making a FAS cycle on a
>> actual fine solution do nothing) at each stage.  As I understand it, we'd
>> need a "real" time derivative rather than a derivative of the defect
>> problem.  Therefore, we would be stuck solving at each level and using each
>> level for FAS independently.
>>
>
> Multilevel time integration would involve having a fully functional TS on
> each level, each of which may or may not use FAS. A conventional time
> integration scheme solved by FAS need these things with type TS to pass to
> the user, that implement TSGetDM() and perhaps some similar functionality.
> They also need the affine vectors to define the residual.
>
> Maybe it would be better to make a little container that was basically
> just a struct, but masqueraded as a TS.
>
>
>> However, having a multilevel time integration scheme outside of FAS
>> accessible to FAS to accelerate solution might be an interesting project;
>> are there references? Would this look like grid sequencing + FAS with the
>> grid sequence solutions moved forward in time?
>>
>
> Hmm, I haven't seen anything really like that, but note that the stability
> threshold for explicit methods is larger on coarse levels.
>



-- 
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-dev/attachments/20120207/d07eef61/attachment.html>


More information about the petsc-dev mailing list