[petsc-dev] getting ready for PETSc release

Jed Brown jed at jedbrown.org
Sun Jun 8 17:00:26 CDT 2014


Matthew Knepley <knepley at gmail.com> writes:
>> What only works on Cartesian grids?  Wave decomposition and such works
>>
>
> My understanding was that the Clawpack algorithms worked only on Cartesian
> grids. Are these "mode decomposed"?

The precise formulation used in Clawpack only works on structured
meshes, but non-conservative hyperbolic systems can be solved on
unstructured meshes and would be implemented in terms of wave
decomposition.

>> on unstructured grids.  PetscProblem is the special-purpose thing that
>> excludes large classes of good methods.  I think it is inappropriate to
>> use the PetscProblem namespace for a specific class of methods.
>
>
> I think you have an idiosyncratic and expansive understanding of the
> word "Problem".  It can also mean a "difficulty", "unwelcome or
> harmful situation", or "complication". If its really down to the word,
> I could see changing it to reduce the number of ultimately unhelpful
> emails.

Programmers, computers, and software exist to solve "problems".  That
word is so phenomenally inspecific as to be useless in a general-purpose
library like PETSc, unless it encompasses the entire scope of problems
for which PETSc is applicable.  Clearly it does not even come close.

>> Why can't you develop this in a separate package?
>
> Too much overhead for too little benefit.

Can you state necessary and sufficient criteria for deciding when
something should be a separate package?

>> >> From the two examples that use PetscProblem, I think the issue is that
>> >> you want it to work with DMDA and DMPlex because the formulation does
>> >> not depend on grid topology, but the naming convention for
>> >> implementation-specific functions uses different namespaces.  But we can
>> >> have DM functions that are implemented by multiple implementations.  So
>> >> maybe this maps to
>> >> DMTSSetIFunctionMattsCoolDiscretizationVolumetricQuadrature().
>> >
>> >
>> > 1) Volumetric here is inappropriate because it supports boundary
>> integrals
>>
>> Through a different function, PetscProblemSetBdResidual().  (I will also
>> complain that this does not support different boundary types.)
>
>
> Criticize away. I did not want to give you the idea that criticism is
> unwelcome.
>
>
>> > 2) It is supposed to support both FEM and FVM
>>
>> But you need new functions for FV because PetscProblemSetResidual
>> definitely does not support it.
>
>
> Yep.

So DMPlexTSSetRHSFunctionLocal is implicitly for upwind-based
finite-volume discretizations of conservative hyperbolic systems and
PetscProblemSetResidual is implicitly for the volumetric term in H^1
finite-element discretizations.  New functions will be added as needed
for solving more general problems or using more general methods?  Why
stake out this canonical namespace for something you know a priori is a
tiny subset of the problems types that should eventually be supported?

>> PetscLayout contains nontrivial data size and can be shared between
>> multiple objects (and we should do a better job of that sharing because
>> PetscLayout storage is not trivial with a million MPI ranks if each
>> object has its own).
>>
>
> PetscProblem can be shared between multiple objects, and has
> non-trivial sizes.  There is no point to the above paragraph and it
> dodges the main issue of the question.

Where is the problem-sized data?

struct _p_PetscProblem {
  PETSCHEADER(struct _PetscProblemOps);
  void        *data;      /* Implementation object */
  PetscBool    setup;     /* Flag for setup */
  PetscInt     Nf;        /* The number of solution fields */
  PetscObject *disc;      /* The discretization for each solution field (PetscFE, PetscFV, etc.) */
  PetscObject *discBd;    /* The boundary discretization for each solution field (PetscFE, PetscFV, etc.) */
  PointFunc   *obj;       /* Scalar integral (like an objective function) */
  PointFunc   *f,   *g;   /* Weak form integrands f_0, f_1, g_0, g_1, g_2, g_3 */
  BdPointFunc *fBd, *gBd; /* Weak form boundary integrands f_0, f_1, g_0, g_1, g_2, g_3 */
  PetscInt     dim;       /* The spatial dimension */
  /* Computed sizes */
  PetscInt     totDim, totDimBd;       /* Total problem dimension */
  PetscInt     totComp;                /* Total field components */
  /* Work space */
  PetscReal  **basis,    **basisBd;    /* Default basis tabulation for each field */
  PetscReal  **basisDer, **basisDerBd; /* Default basis derivative tabulation for each field */
  PetscScalar *u;                      /* Field evaluation */
  PetscScalar *u_t;                    /* Field time derivative evaluation */
  PetscScalar *u_x;                    /* Field gradient evaluation */
  PetscScalar *refSpaceDer;            /* Workspace for computing derivative in the reference coordinates */
  PetscReal   *x;                      /* Workspace for computing real coordinates */
  PetscScalar *f0, *f1;                /* Point evaluations of weak form residual integrands */
  PetscScalar *g0, *g1, *g2, *g3;      /* Point evaluations of weak form Jacobian integrands */
};


I don't think disc and discBd makes sense since combinations may be
either non-unique or incompatible.

Note that this formulation also cannot represent second order elliptic
problems like |D^2 u|^2 = 0.
-------------- next part --------------
A non-text attachment was scrubbed...
Name: not available
Type: application/pgp-signature
Size: 818 bytes
Desc: not available
URL: <http://lists.mcs.anl.gov/pipermail/petsc-dev/attachments/20140609/2ede0e09/attachment.sig>


More information about the petsc-dev mailing list