petsc-3.10.5 2019-03-28
Report Typos and Errors

PetscDSGetBdResidual

Get the pointwise boundary residual function for a given test field

Synopsis

#include "petscds.h" 
PetscErrorCode PetscDSGetBdResidual(PetscDS prob, PetscInt f,
                                    void (**f0)(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[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]),
                                    void (**f1)(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[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]))
Not collective

Input Parameters

prob - The PetscDS
f - The test field number

Output Parameters

f0 - boundary integrand for the test function term
f1 - boundary integrand for the test function gradient term

Note: We are using a first order FEM model for the weak form

\int_\Gamma \phi {\vec f}_0(u, u_t, \nabla u, x, t) \cdot \hat n + \nabla\phi \cdot {\overleftrightarrow f}_1(u, u_t, \nabla u, x, t) \cdot \hat n

The calling sequence for the callbacks f0 and f1 is given by

f0(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[], const PetscReal n[], PetscScalar f0[])

dim - the spatial dimension
Nf - the number of fields
uOff - the offset into u[] and u_t[] for each field
uOff_x - the offset into u_x[] for each field
u - each field evaluated at the current point
u_t - the time derivative of each field evaluated at the current point
u_x - the gradient of each field evaluated at the current point
aOff - the offset into a[] and a_t[] for each auxiliary field
aOff_x - the offset into a_x[] for each auxiliary field
a - each auxiliary field evaluated at the current point
a_t - the time derivative of each auxiliary field evaluated at the current point
a_x - the gradient of auxiliary each field evaluated at the current point
t - current time
x - coordinates of the current point
n - unit normal at the current point
numConstants - number of constant parameters
constants - constant parameters
f0 - output values at the current point

See Also

PetscDSSetBdResidual()

Level

intermediate

Location

src/dm/dt/interface/dtds.c
Index of all DT routines
Table of Contents for all manual pages
Index of all manual pages