DMProjectFieldLabelLocal#

This projects the given function of the input fields into the function space provided, putting the coefficients in a local vector, calculating only over the portion of the domain specified by the label.

Synopsis#

#include "petscdm.h"          
#include "petscdmlabel.h"     
#include "petscds.h"     
PetscErrorCode DMProjectFieldLabelLocal(DM dm, PetscReal time, DMLabel label, PetscInt numIds, const PetscInt ids[], PetscInt Nc, const PetscInt comps[], Vec localU, void (**funcs)(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[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f[]), InsertMode mode, Vec localX)

Not Collective

Input Parameters#

  • dm - The DM

  • time - The time

  • label - The DMLabel marking the portion of the domain to output

  • numIds - The number of label ids to use

  • ids - The label ids to use for marking

  • Nc - The number of components to set in the output, or PETSC_DETERMINE for all components

  • comps - The components to set in the output, or NULL for all components

  • localU - The input field vector

  • funcs - The functions to evaluate, one per field

  • mode - The insertion mode for values

Output Parameter#

  • localX - The output vector

Calling sequence of funcs#

  • dim - The spatial dimension

  • Nf - The number of input fields

  • NfAux - The number of input auxiliary fields

  • uOff - The offset of each field in u[]

  • uOff_x - The offset of each field in u_x[]

  • u - The field values at this point in space

  • u_t - The field time derivative at this point in space (or NULL)

  • u_x - The field derivatives at this point in space

  • aOff - The offset of each auxiliary field in u[]

  • aOff_x - The offset of each auxiliary field in u_x[]

  • a - The auxiliary field values at this point in space

  • a_t - The auxiliary field time derivative at this point in space (or NULL)

  • a_x - The auxiliary field derivatives at this point in space

  • t - The current time

  • x - The coordinates of this point

  • numConstants - The number of constants

  • constants - The value of each constant

  • f - The value of the function at this point in space

Note#

There are three different DMs that potentially interact in this function. The output DM, dm, specifies the layout of the values calculates by funcs. The input DM, attached to localU, may be different. For example, you can input the solution over the full domain, but output over a piece of the boundary, or a subdomain. You can also output a different number of fields than the input, with different discretizations. Last the auxiliary DM, attached to the auxiliary field vector, which is attached to dm, can also be different. It can have a different topology, number of fields, and discretizations.

Developer Notes#

This API is specific to only particular usage of DM

The notes need to provide some information about what has to be provided to the DM to be able to perform the computation.

See Also#

DM Basics, DM, DMProjectField(), DMProjectFieldLabel(), DMProjectFunction(), DMComputeL2Diff()

Level#

intermediate

Location#

src/dm/interface/dm.c

Implementations#

DMProjectFieldLabelLocal_Plex(DM dm, PetscReal time, DMLabel label, PetscInt numIds, const PetscInt ids[], PetscInt Ncc, const PetscInt comps[], Vec localU, void (**funcs)() in src/dm/impls/plex/plexproject.c


Index of all DM routines
Table of Contents for all manual pages
Index of all manual pages