xref: /petsc/src/dm/dt/interface/dtds.c (revision a4e35b1925eceef64945ea472b84f2bf06a67b5e)
1af0996ceSBarry Smith #include <petsc/private/petscdsimpl.h> /*I "petscds.h" I*/
22764a2aaSMatthew G. Knepley 
32764a2aaSMatthew G. Knepley PetscClassId PETSCDS_CLASSID = 0;
42764a2aaSMatthew G. Knepley 
52764a2aaSMatthew G. Knepley PetscFunctionList PetscDSList              = NULL;
62764a2aaSMatthew G. Knepley PetscBool         PetscDSRegisterAllCalled = PETSC_FALSE;
72764a2aaSMatthew G. Knepley 
894dcdc3fSMatthew G. Knepley /* A PetscDS (Discrete System) encodes a set of equations posed in a discrete space, which represents a set of
994dcdc3fSMatthew G. Knepley    nonlinear continuum equations. The equations can have multiple fields, each field having a different
1094dcdc3fSMatthew G. Knepley    discretization. In addition, different pieces of the domain can have different field combinations and equations.
1194dcdc3fSMatthew G. Knepley 
1294dcdc3fSMatthew G. Knepley    The DS provides the user a description of the approximation space on any given cell. It also gives pointwise
1394dcdc3fSMatthew G. Knepley    functions representing the equations.
1494dcdc3fSMatthew G. Knepley 
1594dcdc3fSMatthew G. Knepley    Each field is associated with a label, marking the cells on which it is supported. Note that a field can be
1694dcdc3fSMatthew G. Knepley    supported on the closure of a cell not in the label due to overlap of the boundary of neighboring cells. The DM
1794dcdc3fSMatthew G. Knepley    then creates a DS for each set of cells with identical approximation spaces. When assembling, the user asks for
1894dcdc3fSMatthew G. Knepley    the space associated with a given cell. DMPlex uses the labels associated with each DS in the default integration loop.
1994dcdc3fSMatthew G. Knepley */
2094dcdc3fSMatthew G. Knepley 
212764a2aaSMatthew G. Knepley /*@C
22dce8aebaSBarry Smith   PetscDSRegister - Adds a new `PetscDS` implementation
232764a2aaSMatthew G. Knepley 
2420f4b53cSBarry Smith   Not Collective; No Fortran Support
252764a2aaSMatthew G. Knepley 
262764a2aaSMatthew G. Knepley   Input Parameters:
2720f4b53cSBarry Smith + sname    - The name of a new user-defined creation routine
2820f4b53cSBarry Smith - function - The creation routine itself
292764a2aaSMatthew G. Knepley 
3060225df5SJacob Faibussowitsch   Example Usage:
312764a2aaSMatthew G. Knepley .vb
322764a2aaSMatthew G. Knepley     PetscDSRegister("my_ds", MyPetscDSCreate);
332764a2aaSMatthew G. Knepley .ve
342764a2aaSMatthew G. Knepley 
352764a2aaSMatthew G. Knepley   Then, your PetscDS type can be chosen with the procedural interface via
362764a2aaSMatthew G. Knepley .vb
372764a2aaSMatthew G. Knepley     PetscDSCreate(MPI_Comm, PetscDS *);
382764a2aaSMatthew G. Knepley     PetscDSSetType(PetscDS, "my_ds");
392764a2aaSMatthew G. Knepley .ve
402764a2aaSMatthew G. Knepley   or at runtime via the option
412764a2aaSMatthew G. Knepley .vb
422764a2aaSMatthew G. Knepley     -petscds_type my_ds
432764a2aaSMatthew G. Knepley .ve
442764a2aaSMatthew G. Knepley 
452764a2aaSMatthew G. Knepley   Level: advanced
462764a2aaSMatthew G. Knepley 
47dce8aebaSBarry Smith   Note:
48dce8aebaSBarry Smith   `PetscDSRegister()` may be called multiple times to add several user-defined `PetscDSs`
49dce8aebaSBarry Smith 
50dce8aebaSBarry Smith .seealso: `PetscDSType`, `PetscDS`, `PetscDSRegisterAll()`, `PetscDSRegisterDestroy()`
512764a2aaSMatthew G. Knepley @*/
52d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDSRegister(const char sname[], PetscErrorCode (*function)(PetscDS))
53d71ae5a4SJacob Faibussowitsch {
542764a2aaSMatthew G. Knepley   PetscFunctionBegin;
559566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&PetscDSList, sname, function));
563ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
572764a2aaSMatthew G. Knepley }
582764a2aaSMatthew G. Knepley 
592764a2aaSMatthew G. Knepley /*@C
60dce8aebaSBarry Smith   PetscDSSetType - Builds a particular `PetscDS`
612764a2aaSMatthew G. Knepley 
6220f4b53cSBarry Smith   Collective; No Fortran Support
632764a2aaSMatthew G. Knepley 
642764a2aaSMatthew G. Knepley   Input Parameters:
65dce8aebaSBarry Smith + prob - The `PetscDS` object
66dce8aebaSBarry Smith - name - The `PetscDSType`
672764a2aaSMatthew G. Knepley 
682764a2aaSMatthew G. Knepley   Options Database Key:
692764a2aaSMatthew G. Knepley . -petscds_type <type> - Sets the PetscDS type; use -help for a list of available types
702764a2aaSMatthew G. Knepley 
712764a2aaSMatthew G. Knepley   Level: intermediate
722764a2aaSMatthew G. Knepley 
73dce8aebaSBarry Smith .seealso: `PetscDSType`, `PetscDS`, `PetscDSGetType()`, `PetscDSCreate()`
742764a2aaSMatthew G. Knepley @*/
75d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDSSetType(PetscDS prob, PetscDSType name)
76d71ae5a4SJacob Faibussowitsch {
772764a2aaSMatthew G. Knepley   PetscErrorCode (*r)(PetscDS);
782764a2aaSMatthew G. Knepley   PetscBool match;
792764a2aaSMatthew G. Knepley 
802764a2aaSMatthew G. Knepley   PetscFunctionBegin;
812764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
829566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)prob, name, &match));
833ba16761SJacob Faibussowitsch   if (match) PetscFunctionReturn(PETSC_SUCCESS);
842764a2aaSMatthew G. Knepley 
859566063dSJacob Faibussowitsch   PetscCall(PetscDSRegisterAll());
869566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListFind(PetscDSList, name, &r));
8728b400f6SJacob Faibussowitsch   PetscCheck(r, PetscObjectComm((PetscObject)prob), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown PetscDS type: %s", name);
882764a2aaSMatthew G. Knepley 
89dbbe0bcdSBarry Smith   PetscTryTypeMethod(prob, destroy);
902764a2aaSMatthew G. Knepley   prob->ops->destroy = NULL;
91dbbe0bcdSBarry Smith 
929566063dSJacob Faibussowitsch   PetscCall((*r)(prob));
939566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)prob, name));
943ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
952764a2aaSMatthew G. Knepley }
962764a2aaSMatthew G. Knepley 
972764a2aaSMatthew G. Knepley /*@C
98dce8aebaSBarry Smith   PetscDSGetType - Gets the `PetscDSType` name (as a string) from the `PetscDS`
992764a2aaSMatthew G. Knepley 
10020f4b53cSBarry Smith   Not Collective; No Fortran Support
1012764a2aaSMatthew G. Knepley 
1022764a2aaSMatthew G. Knepley   Input Parameter:
103dce8aebaSBarry Smith . prob - The `PetscDS`
1042764a2aaSMatthew G. Knepley 
1052764a2aaSMatthew G. Knepley   Output Parameter:
106dce8aebaSBarry Smith . name - The `PetscDSType` name
1072764a2aaSMatthew G. Knepley 
1082764a2aaSMatthew G. Knepley   Level: intermediate
1092764a2aaSMatthew G. Knepley 
110dce8aebaSBarry Smith .seealso: `PetscDSType`, `PetscDS`, `PetscDSSetType()`, `PetscDSCreate()`
1112764a2aaSMatthew G. Knepley @*/
112d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDSGetType(PetscDS prob, PetscDSType *name)
113d71ae5a4SJacob Faibussowitsch {
1142764a2aaSMatthew G. Knepley   PetscFunctionBegin;
1152764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1164f572ea9SToby Isaac   PetscAssertPointer(name, 2);
1179566063dSJacob Faibussowitsch   PetscCall(PetscDSRegisterAll());
1182764a2aaSMatthew G. Knepley   *name = ((PetscObject)prob)->type_name;
1193ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1202764a2aaSMatthew G. Knepley }
1212764a2aaSMatthew G. Knepley 
122d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscDSView_Ascii(PetscDS ds, PetscViewer viewer)
123d71ae5a4SJacob Faibussowitsch {
1247d8a60eaSMatthew G. Knepley   PetscViewerFormat  format;
12597b6e6e8SMatthew G. Knepley   const PetscScalar *constants;
1265fedec97SMatthew G. Knepley   PetscInt           Nf, numConstants, f;
1277d8a60eaSMatthew G. Knepley 
1287d8a60eaSMatthew G. Knepley   PetscFunctionBegin;
1299566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(ds, &Nf));
1309566063dSJacob Faibussowitsch   PetscCall(PetscViewerGetFormat(viewer, &format));
13163a3b9bcSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "Discrete System with %" PetscInt_FMT " fields\n", Nf));
1329566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPushTab(viewer));
13363a3b9bcSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "  cell total dim %" PetscInt_FMT " total comp %" PetscInt_FMT "\n", ds->totDim, ds->totComp));
1349566063dSJacob Faibussowitsch   if (ds->isCohesive) PetscCall(PetscViewerASCIIPrintf(viewer, "  cohesive cell\n"));
1355fedec97SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
13640967b3bSMatthew G. Knepley     DSBoundary      b;
1377d8a60eaSMatthew G. Knepley     PetscObject     obj;
1387d8a60eaSMatthew G. Knepley     PetscClassId    id;
139f35450b9SMatthew G. Knepley     PetscQuadrature q;
1407d8a60eaSMatthew G. Knepley     const char     *name;
141f35450b9SMatthew G. Knepley     PetscInt        Nc, Nq, Nqc;
1427d8a60eaSMatthew G. Knepley 
1439566063dSJacob Faibussowitsch     PetscCall(PetscDSGetDiscretization(ds, f, &obj));
1449566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetClassId(obj, &id));
1459566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetName(obj, &name));
1469566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "Field %s", name ? name : "<unknown>"));
1479566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_FALSE));
1487d8a60eaSMatthew G. Knepley     if (id == PETSCFE_CLASSID) {
1499566063dSJacob Faibussowitsch       PetscCall(PetscFEGetNumComponents((PetscFE)obj, &Nc));
1509566063dSJacob Faibussowitsch       PetscCall(PetscFEGetQuadrature((PetscFE)obj, &q));
1519566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, " FEM"));
1527d8a60eaSMatthew G. Knepley     } else if (id == PETSCFV_CLASSID) {
1539566063dSJacob Faibussowitsch       PetscCall(PetscFVGetNumComponents((PetscFV)obj, &Nc));
1549566063dSJacob Faibussowitsch       PetscCall(PetscFVGetQuadrature((PetscFV)obj, &q));
1559566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, " FVM"));
1569371c9d4SSatish Balay     } else SETERRQ(PetscObjectComm((PetscObject)ds), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %" PetscInt_FMT, f);
15763a3b9bcSJacob Faibussowitsch     if (Nc > 1) PetscCall(PetscViewerASCIIPrintf(viewer, " %" PetscInt_FMT " components", Nc));
15863a3b9bcSJacob Faibussowitsch     else PetscCall(PetscViewerASCIIPrintf(viewer, " %" PetscInt_FMT " component ", Nc));
1599566063dSJacob Faibussowitsch     if (ds->implicit[f]) PetscCall(PetscViewerASCIIPrintf(viewer, " (implicit)"));
1609566063dSJacob Faibussowitsch     else PetscCall(PetscViewerASCIIPrintf(viewer, " (explicit)"));
1613e60c2a6SMatthew G. Knepley     if (q) {
1629566063dSJacob Faibussowitsch       PetscCall(PetscQuadratureGetData(q, NULL, &Nqc, &Nq, NULL, NULL));
16363a3b9bcSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, " (Nq %" PetscInt_FMT " Nqc %" PetscInt_FMT ")", Nq, Nqc));
1643e60c2a6SMatthew G. Knepley     }
16563a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, " %" PetscInt_FMT "-jet", ds->jetDegree[f]));
1669566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
1679566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_TRUE));
1689566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushTab(viewer));
1699566063dSJacob Faibussowitsch     if (id == PETSCFE_CLASSID) PetscCall(PetscFEView((PetscFE)obj, viewer));
1709566063dSJacob Faibussowitsch     else if (id == PETSCFV_CLASSID) PetscCall(PetscFVView((PetscFV)obj, viewer));
1719566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopTab(viewer));
17240967b3bSMatthew G. Knepley 
1735fedec97SMatthew G. Knepley     for (b = ds->boundary; b; b = b->next) {
17406ad1575SMatthew G. Knepley       char    *name;
17540967b3bSMatthew G. Knepley       PetscInt c, i;
17640967b3bSMatthew G. Knepley 
17740967b3bSMatthew G. Knepley       if (b->field != f) continue;
1789566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPushTab(viewer));
1799566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "Boundary %s (%s) %s\n", b->name, b->lname, DMBoundaryConditionTypes[b->type]));
18045480ffeSMatthew G. Knepley       if (!b->Nc) {
1819566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  all components\n"));
18240967b3bSMatthew G. Knepley       } else {
1839566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  components: "));
1849566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_FALSE));
18545480ffeSMatthew G. Knepley         for (c = 0; c < b->Nc; ++c) {
1869566063dSJacob Faibussowitsch           if (c > 0) PetscCall(PetscViewerASCIIPrintf(viewer, ", "));
18763a3b9bcSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT, b->comps[c]));
18840967b3bSMatthew G. Knepley         }
1899566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
1909566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_TRUE));
19140967b3bSMatthew G. Knepley       }
1929566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "  values: "));
1939566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_FALSE));
19445480ffeSMatthew G. Knepley       for (i = 0; i < b->Nv; ++i) {
1959566063dSJacob Faibussowitsch         if (i > 0) PetscCall(PetscViewerASCIIPrintf(viewer, ", "));
19663a3b9bcSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT, b->values[i]));
19740967b3bSMatthew G. Knepley       }
1989566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
1999566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_TRUE));
20015943bb8SPierre Jolivet #if defined(__clang__)
201a8f51744SPierre Jolivet       PETSC_PRAGMA_DIAGNOSTIC_IGNORED_BEGIN("-Wformat-pedantic")
20215943bb8SPierre Jolivet #elif defined(__GNUC__) || defined(__GNUG__)
203a8f51744SPierre Jolivet       PETSC_PRAGMA_DIAGNOSTIC_IGNORED_BEGIN("-Wformat")
20415943bb8SPierre Jolivet #endif
2058e0d8d9cSMatthew G. Knepley       if (b->func) {
2069566063dSJacob Faibussowitsch         PetscCall(PetscDLAddr(b->func, &name));
2079566063dSJacob Faibussowitsch         if (name) PetscCall(PetscViewerASCIIPrintf(viewer, "  func: %s\n", name));
2089566063dSJacob Faibussowitsch         else PetscCall(PetscViewerASCIIPrintf(viewer, "  func: %p\n", b->func));
2099566063dSJacob Faibussowitsch         PetscCall(PetscFree(name));
2108e0d8d9cSMatthew G. Knepley       }
2118e0d8d9cSMatthew G. Knepley       if (b->func_t) {
2129566063dSJacob Faibussowitsch         PetscCall(PetscDLAddr(b->func_t, &name));
2139566063dSJacob Faibussowitsch         if (name) PetscCall(PetscViewerASCIIPrintf(viewer, "  func_t: %s\n", name));
2149566063dSJacob Faibussowitsch         else PetscCall(PetscViewerASCIIPrintf(viewer, "  func_t: %p\n", b->func_t));
2159566063dSJacob Faibussowitsch         PetscCall(PetscFree(name));
2168e0d8d9cSMatthew G. Knepley       }
217a8f51744SPierre Jolivet       PETSC_PRAGMA_DIAGNOSTIC_IGNORED_END()
2189566063dSJacob Faibussowitsch       PetscCall(PetscWeakFormView(b->wf, viewer));
2199566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPopTab(viewer));
22040967b3bSMatthew G. Knepley     }
2217d8a60eaSMatthew G. Knepley   }
2229566063dSJacob Faibussowitsch   PetscCall(PetscDSGetConstants(ds, &numConstants, &constants));
22397b6e6e8SMatthew G. Knepley   if (numConstants) {
22463a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT " constants\n", numConstants));
2259566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushTab(viewer));
2269566063dSJacob Faibussowitsch     for (f = 0; f < numConstants; ++f) PetscCall(PetscViewerASCIIPrintf(viewer, "%g\n", (double)PetscRealPart(constants[f])));
2279566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopTab(viewer));
22897b6e6e8SMatthew G. Knepley   }
2299566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormView(ds->wf, viewer));
2309566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPopTab(viewer));
2313ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2327d8a60eaSMatthew G. Knepley }
2337d8a60eaSMatthew G. Knepley 
2342764a2aaSMatthew G. Knepley /*@C
235dce8aebaSBarry Smith   PetscDSViewFromOptions - View a `PetscDS` based on values in the options database
236fe2efc57SMark 
23720f4b53cSBarry Smith   Collective
238fe2efc57SMark 
239fe2efc57SMark   Input Parameters:
240dce8aebaSBarry Smith + A    - the `PetscDS` object
24120f4b53cSBarry Smith . obj  - Optional object that provides the options prefix used in the search
242736c3998SJose E. Roman - name - command line option
243fe2efc57SMark 
244fe2efc57SMark   Level: intermediate
245dce8aebaSBarry Smith 
246dce8aebaSBarry Smith .seealso: `PetscDSType`, `PetscDS`, `PetscDSView()`, `PetscObjectViewFromOptions()`, `PetscDSCreate()`
247fe2efc57SMark @*/
248d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDSViewFromOptions(PetscDS A, PetscObject obj, const char name[])
249d71ae5a4SJacob Faibussowitsch {
250fe2efc57SMark   PetscFunctionBegin;
251fe2efc57SMark   PetscValidHeaderSpecific(A, PETSCDS_CLASSID, 1);
2529566063dSJacob Faibussowitsch   PetscCall(PetscObjectViewFromOptions((PetscObject)A, obj, name));
2533ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
254fe2efc57SMark }
255fe2efc57SMark 
256fe2efc57SMark /*@C
257dce8aebaSBarry Smith   PetscDSView - Views a `PetscDS`
2582764a2aaSMatthew G. Knepley 
25920f4b53cSBarry Smith   Collective
2602764a2aaSMatthew G. Knepley 
261d8d19677SJose E. Roman   Input Parameters:
262dce8aebaSBarry Smith + prob - the `PetscDS` object to view
2632764a2aaSMatthew G. Knepley - v    - the viewer
2642764a2aaSMatthew G. Knepley 
2652764a2aaSMatthew G. Knepley   Level: developer
2662764a2aaSMatthew G. Knepley 
26720f4b53cSBarry Smith .seealso: `PetscDSType`, `PetscDS`, `PetscViewer`, `PetscDSDestroy()`, `PetscDSViewFromOptions()`
2682764a2aaSMatthew G. Knepley @*/
269d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDSView(PetscDS prob, PetscViewer v)
270d71ae5a4SJacob Faibussowitsch {
2717d8a60eaSMatthew G. Knepley   PetscBool iascii;
2722764a2aaSMatthew G. Knepley 
2732764a2aaSMatthew G. Knepley   PetscFunctionBegin;
2742764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
2759566063dSJacob Faibussowitsch   if (!v) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)prob), &v));
276ad540459SPierre Jolivet   else PetscValidHeaderSpecific(v, PETSC_VIEWER_CLASSID, 2);
2779566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)v, PETSCVIEWERASCII, &iascii));
2789566063dSJacob Faibussowitsch   if (iascii) PetscCall(PetscDSView_Ascii(prob, v));
279dbbe0bcdSBarry Smith   PetscTryTypeMethod(prob, view, v);
2803ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2812764a2aaSMatthew G. Knepley }
2822764a2aaSMatthew G. Knepley 
2832764a2aaSMatthew G. Knepley /*@
284dce8aebaSBarry Smith   PetscDSSetFromOptions - sets parameters in a `PetscDS` from the options database
2852764a2aaSMatthew G. Knepley 
28620f4b53cSBarry Smith   Collective
2872764a2aaSMatthew G. Knepley 
2882764a2aaSMatthew G. Knepley   Input Parameter:
289dce8aebaSBarry Smith . prob - the `PetscDS` object to set options for
2902764a2aaSMatthew G. Knepley 
291dce8aebaSBarry Smith   Options Database Keys:
292dce8aebaSBarry Smith + -petscds_type <type>     - Set the `PetscDS` type
293dce8aebaSBarry Smith . -petscds_view <view opt> - View the `PetscDS`
294147403d9SBarry Smith . -petscds_jac_pre         - Turn formation of a separate Jacobian preconditioner on or off
295147403d9SBarry Smith . -bc_<name> <ids>         - Specify a list of label ids for a boundary condition
296147403d9SBarry Smith - -bc_<name>_comp <comps>  - Specify a list of field components to constrain for a boundary condition
2972764a2aaSMatthew G. Knepley 
298dce8aebaSBarry Smith   Level: intermediate
2992764a2aaSMatthew G. Knepley 
300dce8aebaSBarry Smith .seealso: `PetscDS`, `PetscDSView()`
3012764a2aaSMatthew G. Knepley @*/
302d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDSSetFromOptions(PetscDS prob)
303d71ae5a4SJacob Faibussowitsch {
304f1fd5e65SToby Isaac   DSBoundary  b;
3052764a2aaSMatthew G. Knepley   const char *defaultType;
3062764a2aaSMatthew G. Knepley   char        name[256];
3072764a2aaSMatthew G. Knepley   PetscBool   flg;
3082764a2aaSMatthew G. Knepley 
3092764a2aaSMatthew G. Knepley   PetscFunctionBegin;
3102764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
3112764a2aaSMatthew G. Knepley   if (!((PetscObject)prob)->type_name) {
3122764a2aaSMatthew G. Knepley     defaultType = PETSCDSBASIC;
3132764a2aaSMatthew G. Knepley   } else {
3142764a2aaSMatthew G. Knepley     defaultType = ((PetscObject)prob)->type_name;
3152764a2aaSMatthew G. Knepley   }
3169566063dSJacob Faibussowitsch   PetscCall(PetscDSRegisterAll());
3172764a2aaSMatthew G. Knepley 
318d0609cedSBarry Smith   PetscObjectOptionsBegin((PetscObject)prob);
319f1fd5e65SToby Isaac   for (b = prob->boundary; b; b = b->next) {
320f1fd5e65SToby Isaac     char      optname[1024];
321f1fd5e65SToby Isaac     PetscInt  ids[1024], len = 1024;
322f1fd5e65SToby Isaac     PetscBool flg;
323f1fd5e65SToby Isaac 
3249566063dSJacob Faibussowitsch     PetscCall(PetscSNPrintf(optname, sizeof(optname), "-bc_%s", b->name));
3259566063dSJacob Faibussowitsch     PetscCall(PetscMemzero(ids, sizeof(ids)));
3269566063dSJacob Faibussowitsch     PetscCall(PetscOptionsIntArray(optname, "List of boundary IDs", "", ids, &len, &flg));
327f1fd5e65SToby Isaac     if (flg) {
32845480ffeSMatthew G. Knepley       b->Nv = len;
3299566063dSJacob Faibussowitsch       PetscCall(PetscFree(b->values));
3309566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(len, &b->values));
3319566063dSJacob Faibussowitsch       PetscCall(PetscArraycpy(b->values, ids, len));
3329566063dSJacob Faibussowitsch       PetscCall(PetscWeakFormRewriteKeys(b->wf, b->label, len, b->values));
333f1fd5e65SToby Isaac     }
334e7b0402cSSander Arens     len = 1024;
3359566063dSJacob Faibussowitsch     PetscCall(PetscSNPrintf(optname, sizeof(optname), "-bc_%s_comp", b->name));
3369566063dSJacob Faibussowitsch     PetscCall(PetscMemzero(ids, sizeof(ids)));
3379566063dSJacob Faibussowitsch     PetscCall(PetscOptionsIntArray(optname, "List of boundary field components", "", ids, &len, &flg));
338f1fd5e65SToby Isaac     if (flg) {
33945480ffeSMatthew G. Knepley       b->Nc = len;
3409566063dSJacob Faibussowitsch       PetscCall(PetscFree(b->comps));
3419566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(len, &b->comps));
3429566063dSJacob Faibussowitsch       PetscCall(PetscArraycpy(b->comps, ids, len));
343f1fd5e65SToby Isaac     }
344f1fd5e65SToby Isaac   }
3459566063dSJacob Faibussowitsch   PetscCall(PetscOptionsFList("-petscds_type", "Discrete System", "PetscDSSetType", PetscDSList, defaultType, name, 256, &flg));
3462764a2aaSMatthew G. Knepley   if (flg) {
3479566063dSJacob Faibussowitsch     PetscCall(PetscDSSetType(prob, name));
3482764a2aaSMatthew G. Knepley   } else if (!((PetscObject)prob)->type_name) {
3499566063dSJacob Faibussowitsch     PetscCall(PetscDSSetType(prob, defaultType));
3502764a2aaSMatthew G. Knepley   }
3519566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-petscds_jac_pre", "Discrete System", "PetscDSUseJacobianPreconditioner", prob->useJacPre, &prob->useJacPre, &flg));
35212fc5b22SMatthew G. Knepley   PetscCall(PetscOptionsBool("-petscds_force_quad", "Discrete System", "PetscDSSetForceQuad", prob->forceQuad, &prob->forceQuad, &flg));
353dbbe0bcdSBarry Smith   PetscTryTypeMethod(prob, setfromoptions);
3542764a2aaSMatthew G. Knepley   /* process any options handlers added with PetscObjectAddOptionsHandler() */
355dbbe0bcdSBarry Smith   PetscCall(PetscObjectProcessOptionsHandlers((PetscObject)prob, PetscOptionsObject));
356d0609cedSBarry Smith   PetscOptionsEnd();
3579566063dSJacob Faibussowitsch   if (prob->Nf) PetscCall(PetscDSViewFromOptions(prob, NULL, "-petscds_view"));
3583ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3592764a2aaSMatthew G. Knepley }
3602764a2aaSMatthew G. Knepley 
3612764a2aaSMatthew G. Knepley /*@C
362dce8aebaSBarry Smith   PetscDSSetUp - Construct data structures for the `PetscDS`
3632764a2aaSMatthew G. Knepley 
36420f4b53cSBarry Smith   Collective
3652764a2aaSMatthew G. Knepley 
3662764a2aaSMatthew G. Knepley   Input Parameter:
367dce8aebaSBarry Smith . prob - the `PetscDS` object to setup
3682764a2aaSMatthew G. Knepley 
3692764a2aaSMatthew G. Knepley   Level: developer
3702764a2aaSMatthew G. Knepley 
371dce8aebaSBarry Smith .seealso: `PetscDS`, `PetscDSView()`, `PetscDSDestroy()`
3722764a2aaSMatthew G. Knepley @*/
373d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDSSetUp(PetscDS prob)
374d71ae5a4SJacob Faibussowitsch {
3752764a2aaSMatthew G. Knepley   const PetscInt Nf          = prob->Nf;
376f9244615SMatthew G. Knepley   PetscBool      hasH        = PETSC_FALSE;
37707218a29SMatthew G. Knepley   PetscInt       maxOrder[4] = {-1, -1, -1, -1};
3784bee2e38SMatthew G. Knepley   PetscInt       dim, dimEmbed, NbMax = 0, NcMax = 0, NqMax = 0, NsMax = 1, f;
3792764a2aaSMatthew G. Knepley 
3802764a2aaSMatthew G. Knepley   PetscFunctionBegin;
3812764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
3823ba16761SJacob Faibussowitsch   if (prob->setup) PetscFunctionReturn(PETSC_SUCCESS);
3832764a2aaSMatthew G. Knepley   /* Calculate sizes */
3849566063dSJacob Faibussowitsch   PetscCall(PetscDSGetSpatialDimension(prob, &dim));
3859566063dSJacob Faibussowitsch   PetscCall(PetscDSGetCoordinateDimension(prob, &dimEmbed));
386f744cafaSSander Arens   prob->totDim = prob->totComp = 0;
3879566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(Nf, &prob->Nc, Nf, &prob->Nb));
3889566063dSJacob Faibussowitsch   PetscCall(PetscCalloc2(Nf + 1, &prob->off, Nf + 1, &prob->offDer));
3899566063dSJacob Faibussowitsch   PetscCall(PetscCalloc6(Nf + 1, &prob->offCohesive[0], Nf + 1, &prob->offCohesive[1], Nf + 1, &prob->offCohesive[2], Nf + 1, &prob->offDerCohesive[0], Nf + 1, &prob->offDerCohesive[1], Nf + 1, &prob->offDerCohesive[2]));
3909566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(Nf, &prob->T, Nf, &prob->Tf));
39112fc5b22SMatthew G. Knepley   if (prob->forceQuad) {
39207218a29SMatthew G. Knepley     // Note: This assumes we have one kind of cell at each dimension.
39307218a29SMatthew G. Knepley     //       We can fix this by having quadrature hold the celltype
39407218a29SMatthew G. Knepley     PetscQuadrature maxQuad[4] = {NULL, NULL, NULL, NULL};
39512fc5b22SMatthew G. Knepley 
39612fc5b22SMatthew G. Knepley     for (f = 0; f < Nf; ++f) {
39712fc5b22SMatthew G. Knepley       PetscObject     obj;
39812fc5b22SMatthew G. Knepley       PetscClassId    id;
39912fc5b22SMatthew G. Knepley       PetscQuadrature q = NULL, fq = NULL;
40007218a29SMatthew G. Knepley       PetscInt        dim = -1, order = -1, forder = -1;
40112fc5b22SMatthew G. Knepley 
40212fc5b22SMatthew G. Knepley       PetscCall(PetscDSGetDiscretization(prob, f, &obj));
40312fc5b22SMatthew G. Knepley       if (!obj) continue;
40412fc5b22SMatthew G. Knepley       PetscCall(PetscObjectGetClassId(obj, &id));
40512fc5b22SMatthew G. Knepley       if (id == PETSCFE_CLASSID) {
40612fc5b22SMatthew G. Knepley         PetscFE fe = (PetscFE)obj;
40712fc5b22SMatthew G. Knepley 
40812fc5b22SMatthew G. Knepley         PetscCall(PetscFEGetQuadrature(fe, &q));
40912fc5b22SMatthew G. Knepley         PetscCall(PetscFEGetFaceQuadrature(fe, &fq));
41012fc5b22SMatthew G. Knepley       } else if (id == PETSCFV_CLASSID) {
41112fc5b22SMatthew G. Knepley         PetscFV fv = (PetscFV)obj;
41212fc5b22SMatthew G. Knepley 
41312fc5b22SMatthew G. Knepley         PetscCall(PetscFVGetQuadrature(fv, &q));
41412fc5b22SMatthew G. Knepley       }
41507218a29SMatthew G. Knepley       if (q) {
41607218a29SMatthew G. Knepley         PetscCall(PetscQuadratureGetData(q, &dim, NULL, NULL, NULL, NULL));
41707218a29SMatthew G. Knepley         PetscCall(PetscQuadratureGetOrder(q, &order));
41807218a29SMatthew G. Knepley         if (order > maxOrder[dim]) {
41907218a29SMatthew G. Knepley           maxOrder[dim] = order;
42007218a29SMatthew G. Knepley           maxQuad[dim]  = q;
42112fc5b22SMatthew G. Knepley         }
42207218a29SMatthew G. Knepley       }
42307218a29SMatthew G. Knepley       if (fq) {
42407218a29SMatthew G. Knepley         PetscCall(PetscQuadratureGetData(fq, &dim, NULL, NULL, NULL, NULL));
42507218a29SMatthew G. Knepley         PetscCall(PetscQuadratureGetOrder(fq, &forder));
42607218a29SMatthew G. Knepley         if (forder > maxOrder[dim]) {
42707218a29SMatthew G. Knepley           maxOrder[dim] = forder;
42807218a29SMatthew G. Knepley           maxQuad[dim]  = fq;
42907218a29SMatthew G. Knepley         }
43012fc5b22SMatthew G. Knepley       }
43112fc5b22SMatthew G. Knepley     }
43212fc5b22SMatthew G. Knepley     for (f = 0; f < Nf; ++f) {
43312fc5b22SMatthew G. Knepley       PetscObject     obj;
43412fc5b22SMatthew G. Knepley       PetscClassId    id;
43507218a29SMatthew G. Knepley       PetscQuadrature q;
43607218a29SMatthew G. Knepley       PetscInt        dim;
43712fc5b22SMatthew G. Knepley 
43812fc5b22SMatthew G. Knepley       PetscCall(PetscDSGetDiscretization(prob, f, &obj));
43912fc5b22SMatthew G. Knepley       if (!obj) continue;
44012fc5b22SMatthew G. Knepley       PetscCall(PetscObjectGetClassId(obj, &id));
44112fc5b22SMatthew G. Knepley       if (id == PETSCFE_CLASSID) {
44212fc5b22SMatthew G. Knepley         PetscFE fe = (PetscFE)obj;
44312fc5b22SMatthew G. Knepley 
44407218a29SMatthew G. Knepley         PetscCall(PetscFEGetQuadrature(fe, &q));
44507218a29SMatthew G. Knepley         PetscCall(PetscQuadratureGetData(q, &dim, NULL, NULL, NULL, NULL));
44607218a29SMatthew G. Knepley         PetscCall(PetscFESetQuadrature(fe, maxQuad[dim]));
447aa9788aaSMatthew G. Knepley         PetscCall(PetscFESetFaceQuadrature(fe, dim ? maxQuad[dim - 1] : NULL));
44812fc5b22SMatthew G. Knepley       } else if (id == PETSCFV_CLASSID) {
44912fc5b22SMatthew G. Knepley         PetscFV fv = (PetscFV)obj;
45012fc5b22SMatthew G. Knepley 
45107218a29SMatthew G. Knepley         PetscCall(PetscFVGetQuadrature(fv, &q));
45207218a29SMatthew G. Knepley         PetscCall(PetscQuadratureGetData(q, &dim, NULL, NULL, NULL, NULL));
45307218a29SMatthew G. Knepley         PetscCall(PetscFVSetQuadrature(fv, maxQuad[dim]));
45412fc5b22SMatthew G. Knepley       }
45512fc5b22SMatthew G. Knepley     }
45612fc5b22SMatthew G. Knepley   }
4572764a2aaSMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
4589de99aefSMatthew G. Knepley     PetscObject     obj;
4599de99aefSMatthew G. Knepley     PetscClassId    id;
460665f567fSMatthew G. Knepley     PetscQuadrature q  = NULL;
4619de99aefSMatthew G. Knepley     PetscInt        Nq = 0, Nb, Nc;
4622764a2aaSMatthew G. Knepley 
4639566063dSJacob Faibussowitsch     PetscCall(PetscDSGetDiscretization(prob, f, &obj));
464f9244615SMatthew G. Knepley     if (prob->jetDegree[f] > 1) hasH = PETSC_TRUE;
465665f567fSMatthew G. Knepley     if (!obj) {
466665f567fSMatthew G. Knepley       /* Empty mesh */
467665f567fSMatthew G. Knepley       Nb = Nc    = 0;
468665f567fSMatthew G. Knepley       prob->T[f] = prob->Tf[f] = NULL;
469665f567fSMatthew G. Knepley     } else {
4709566063dSJacob Faibussowitsch       PetscCall(PetscObjectGetClassId(obj, &id));
4719de99aefSMatthew G. Knepley       if (id == PETSCFE_CLASSID) {
4729de99aefSMatthew G. Knepley         PetscFE fe = (PetscFE)obj;
4739de99aefSMatthew G. Knepley 
4749566063dSJacob Faibussowitsch         PetscCall(PetscFEGetQuadrature(fe, &q));
47549ae0b56SMatthew G. Knepley         {
47649ae0b56SMatthew G. Knepley           PetscQuadrature fq;
47707218a29SMatthew G. Knepley           PetscInt        dim, order;
47849ae0b56SMatthew G. Knepley 
47907218a29SMatthew G. Knepley           PetscCall(PetscQuadratureGetData(q, &dim, NULL, NULL, NULL, NULL));
48049ae0b56SMatthew G. Knepley           PetscCall(PetscQuadratureGetOrder(q, &order));
48107218a29SMatthew G. Knepley           if (maxOrder[dim] < 0) maxOrder[dim] = order;
48207218a29SMatthew G. Knepley           PetscCheck(order == maxOrder[dim], PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Field %" PetscInt_FMT " cell quadrature order %" PetscInt_FMT " != %" PetscInt_FMT " DS cell quadrature order", f, order, maxOrder[dim]);
48349ae0b56SMatthew G. Knepley           PetscCall(PetscFEGetFaceQuadrature(fe, &fq));
48449ae0b56SMatthew G. Knepley           if (fq) {
48507218a29SMatthew G. Knepley             PetscCall(PetscQuadratureGetData(fq, &dim, NULL, NULL, NULL, NULL));
48607218a29SMatthew G. Knepley             PetscCall(PetscQuadratureGetOrder(fq, &order));
48707218a29SMatthew G. Knepley             if (maxOrder[dim] < 0) maxOrder[dim] = order;
48807218a29SMatthew G. Knepley             PetscCheck(order == maxOrder[dim], PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Field %" PetscInt_FMT " face quadrature order %" PetscInt_FMT " != %" PetscInt_FMT " DS face quadrature order", f, order, maxOrder[dim]);
48949ae0b56SMatthew G. Knepley           }
49049ae0b56SMatthew G. Knepley         }
4919566063dSJacob Faibussowitsch         PetscCall(PetscFEGetDimension(fe, &Nb));
4929566063dSJacob Faibussowitsch         PetscCall(PetscFEGetNumComponents(fe, &Nc));
4939566063dSJacob Faibussowitsch         PetscCall(PetscFEGetCellTabulation(fe, prob->jetDegree[f], &prob->T[f]));
4949566063dSJacob Faibussowitsch         PetscCall(PetscFEGetFaceTabulation(fe, prob->jetDegree[f], &prob->Tf[f]));
4959de99aefSMatthew G. Knepley       } else if (id == PETSCFV_CLASSID) {
4969de99aefSMatthew G. Knepley         PetscFV fv = (PetscFV)obj;
4979de99aefSMatthew G. Knepley 
4989566063dSJacob Faibussowitsch         PetscCall(PetscFVGetQuadrature(fv, &q));
4999566063dSJacob Faibussowitsch         PetscCall(PetscFVGetNumComponents(fv, &Nc));
5009c3cf19fSMatthew G. Knepley         Nb = Nc;
5019566063dSJacob Faibussowitsch         PetscCall(PetscFVGetCellTabulation(fv, &prob->T[f]));
5024d0b9603SSander Arens         /* TODO: should PetscFV also have face tabulation? Otherwise there will be a null pointer in prob->basisFace */
50363a3b9bcSJacob Faibussowitsch       } else SETERRQ(PetscObjectComm((PetscObject)prob), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %" PetscInt_FMT, f);
504665f567fSMatthew G. Knepley     }
50547e57110SSander Arens     prob->Nc[f]                    = Nc;
50647e57110SSander Arens     prob->Nb[f]                    = Nb;
507194d53e6SMatthew G. Knepley     prob->off[f + 1]               = Nc + prob->off[f];
508194d53e6SMatthew G. Knepley     prob->offDer[f + 1]            = Nc * dim + prob->offDer[f];
5099ee2af8cSMatthew G. Knepley     prob->offCohesive[0][f + 1]    = (prob->cohesive[f] ? Nc : Nc * 2) + prob->offCohesive[0][f];
5109ee2af8cSMatthew G. Knepley     prob->offDerCohesive[0][f + 1] = (prob->cohesive[f] ? Nc : Nc * 2) * dimEmbed + prob->offDerCohesive[0][f];
5119ee2af8cSMatthew G. Knepley     prob->offCohesive[1][f]        = (prob->cohesive[f] ? 0 : Nc) + prob->offCohesive[0][f];
5129ee2af8cSMatthew G. Knepley     prob->offDerCohesive[1][f]     = (prob->cohesive[f] ? 0 : Nc) * dimEmbed + prob->offDerCohesive[0][f];
5139ee2af8cSMatthew G. Knepley     prob->offCohesive[2][f + 1]    = (prob->cohesive[f] ? Nc : Nc * 2) + prob->offCohesive[2][f];
5149ee2af8cSMatthew G. Knepley     prob->offDerCohesive[2][f + 1] = (prob->cohesive[f] ? Nc : Nc * 2) * dimEmbed + prob->offDerCohesive[2][f];
5159566063dSJacob Faibussowitsch     if (q) PetscCall(PetscQuadratureGetData(q, NULL, NULL, &Nq, NULL, NULL));
5162764a2aaSMatthew G. Knepley     NqMax = PetscMax(NqMax, Nq);
5174bee2e38SMatthew G. Knepley     NbMax = PetscMax(NbMax, Nb);
5182764a2aaSMatthew G. Knepley     NcMax = PetscMax(NcMax, Nc);
5199c3cf19fSMatthew G. Knepley     prob->totDim += Nb;
5202764a2aaSMatthew G. Knepley     prob->totComp += Nc;
5215fedec97SMatthew G. Knepley     /* There are two faces for all fields on a cohesive cell, except for cohesive fields */
5225fedec97SMatthew G. Knepley     if (prob->isCohesive && !prob->cohesive[f]) prob->totDim += Nb;
5232764a2aaSMatthew G. Knepley   }
5249ee2af8cSMatthew G. Knepley   prob->offCohesive[1][Nf]    = prob->offCohesive[0][Nf];
5259ee2af8cSMatthew G. Knepley   prob->offDerCohesive[1][Nf] = prob->offDerCohesive[0][Nf];
5262764a2aaSMatthew G. Knepley   /* Allocate works space */
5275fedec97SMatthew G. Knepley   NsMax = 2; /* A non-cohesive discretizations can be used on a cohesive cell, so we need this extra workspace for all DS */
5289566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3(NsMax * prob->totComp, &prob->u, NsMax * prob->totComp, &prob->u_t, NsMax * prob->totComp * dimEmbed + (hasH ? NsMax * prob->totComp * dimEmbed * dimEmbed : 0), &prob->u_x));
5299566063dSJacob Faibussowitsch   PetscCall(PetscMalloc5(dimEmbed, &prob->x, NbMax * NcMax, &prob->basisReal, NbMax * NcMax * dimEmbed, &prob->basisDerReal, NbMax * NcMax, &prob->testReal, NbMax * NcMax * dimEmbed, &prob->testDerReal));
5309371c9d4SSatish Balay   PetscCall(PetscMalloc6(NsMax * NqMax * NcMax, &prob->f0, NsMax * NqMax * NcMax * dimEmbed, &prob->f1, NsMax * NsMax * NqMax * NcMax * NcMax, &prob->g0, NsMax * NsMax * NqMax * NcMax * NcMax * dimEmbed, &prob->g1, NsMax * NsMax * NqMax * NcMax * NcMax * dimEmbed,
5319371c9d4SSatish Balay                          &prob->g2, NsMax * NsMax * NqMax * NcMax * NcMax * dimEmbed * dimEmbed, &prob->g3));
532dbbe0bcdSBarry Smith   PetscTryTypeMethod(prob, setup);
5332764a2aaSMatthew G. Knepley   prob->setup = PETSC_TRUE;
5343ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5352764a2aaSMatthew G. Knepley }
5362764a2aaSMatthew G. Knepley 
537d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscDSDestroyStructs_Static(PetscDS prob)
538d71ae5a4SJacob Faibussowitsch {
5392764a2aaSMatthew G. Knepley   PetscFunctionBegin;
5409566063dSJacob Faibussowitsch   PetscCall(PetscFree2(prob->Nc, prob->Nb));
5419566063dSJacob Faibussowitsch   PetscCall(PetscFree2(prob->off, prob->offDer));
5429566063dSJacob Faibussowitsch   PetscCall(PetscFree6(prob->offCohesive[0], prob->offCohesive[1], prob->offCohesive[2], prob->offDerCohesive[0], prob->offDerCohesive[1], prob->offDerCohesive[2]));
5439566063dSJacob Faibussowitsch   PetscCall(PetscFree2(prob->T, prob->Tf));
5449566063dSJacob Faibussowitsch   PetscCall(PetscFree3(prob->u, prob->u_t, prob->u_x));
5459566063dSJacob Faibussowitsch   PetscCall(PetscFree5(prob->x, prob->basisReal, prob->basisDerReal, prob->testReal, prob->testDerReal));
5469566063dSJacob Faibussowitsch   PetscCall(PetscFree6(prob->f0, prob->f1, prob->g0, prob->g1, prob->g2, prob->g3));
5473ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5482764a2aaSMatthew G. Knepley }
5492764a2aaSMatthew G. Knepley 
550d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscDSEnlarge_Static(PetscDS prob, PetscInt NfNew)
551d71ae5a4SJacob Faibussowitsch {
552f744cafaSSander Arens   PetscObject          *tmpd;
55334aa8a36SMatthew G. Knepley   PetscBool            *tmpi;
554f9244615SMatthew G. Knepley   PetscInt             *tmpk;
5555fedec97SMatthew G. Knepley   PetscBool            *tmpc;
5566528b96dSMatthew G. Knepley   PetscPointFunc       *tmpup;
557f2cacb80SMatthew G. Knepley   PetscSimplePointFunc *tmpexactSol, *tmpexactSol_t;
558f2cacb80SMatthew G. Knepley   void                **tmpexactCtx, **tmpexactCtx_t;
5590c2f2876SMatthew G. Knepley   void                **tmpctx;
56034aa8a36SMatthew G. Knepley   PetscInt              Nf = prob->Nf, f;
5612764a2aaSMatthew G. Knepley 
5622764a2aaSMatthew G. Knepley   PetscFunctionBegin;
5633ba16761SJacob Faibussowitsch   if (Nf >= NfNew) PetscFunctionReturn(PETSC_SUCCESS);
5642764a2aaSMatthew G. Knepley   prob->setup = PETSC_FALSE;
5659566063dSJacob Faibussowitsch   PetscCall(PetscDSDestroyStructs_Static(prob));
5669566063dSJacob Faibussowitsch   PetscCall(PetscMalloc4(NfNew, &tmpd, NfNew, &tmpi, NfNew, &tmpc, NfNew, &tmpk));
5679371c9d4SSatish Balay   for (f = 0; f < Nf; ++f) {
5689371c9d4SSatish Balay     tmpd[f] = prob->disc[f];
5699371c9d4SSatish Balay     tmpi[f] = prob->implicit[f];
5709371c9d4SSatish Balay     tmpc[f] = prob->cohesive[f];
5719371c9d4SSatish Balay     tmpk[f] = prob->jetDegree[f];
5729371c9d4SSatish Balay   }
5739371c9d4SSatish Balay   for (f = Nf; f < NfNew; ++f) {
5749371c9d4SSatish Balay     tmpd[f] = NULL;
5759371c9d4SSatish Balay     tmpi[f] = PETSC_TRUE, tmpc[f] = PETSC_FALSE;
5769371c9d4SSatish Balay     tmpk[f] = 1;
5779371c9d4SSatish Balay   }
5789566063dSJacob Faibussowitsch   PetscCall(PetscFree4(prob->disc, prob->implicit, prob->cohesive, prob->jetDegree));
5799566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormSetNumFields(prob->wf, NfNew));
5802764a2aaSMatthew G. Knepley   prob->Nf        = NfNew;
5812764a2aaSMatthew G. Knepley   prob->disc      = tmpd;
582249df284SMatthew G. Knepley   prob->implicit  = tmpi;
5835fedec97SMatthew G. Knepley   prob->cohesive  = tmpc;
584f9244615SMatthew G. Knepley   prob->jetDegree = tmpk;
5859566063dSJacob Faibussowitsch   PetscCall(PetscCalloc2(NfNew, &tmpup, NfNew, &tmpctx));
58632d2bbc9SMatthew G. Knepley   for (f = 0; f < Nf; ++f) tmpup[f] = prob->update[f];
5870c2f2876SMatthew G. Knepley   for (f = 0; f < Nf; ++f) tmpctx[f] = prob->ctx[f];
58832d2bbc9SMatthew G. Knepley   for (f = Nf; f < NfNew; ++f) tmpup[f] = NULL;
5890c2f2876SMatthew G. Knepley   for (f = Nf; f < NfNew; ++f) tmpctx[f] = NULL;
5909566063dSJacob Faibussowitsch   PetscCall(PetscFree2(prob->update, prob->ctx));
59132d2bbc9SMatthew G. Knepley   prob->update = tmpup;
5920c2f2876SMatthew G. Knepley   prob->ctx    = tmpctx;
5939566063dSJacob Faibussowitsch   PetscCall(PetscCalloc4(NfNew, &tmpexactSol, NfNew, &tmpexactCtx, NfNew, &tmpexactSol_t, NfNew, &tmpexactCtx_t));
594c371a6d1SMatthew G. Knepley   for (f = 0; f < Nf; ++f) tmpexactSol[f] = prob->exactSol[f];
59595cbbfd3SMatthew G. Knepley   for (f = 0; f < Nf; ++f) tmpexactCtx[f] = prob->exactCtx[f];
596f2cacb80SMatthew G. Knepley   for (f = 0; f < Nf; ++f) tmpexactSol_t[f] = prob->exactSol_t[f];
597f2cacb80SMatthew G. Knepley   for (f = 0; f < Nf; ++f) tmpexactCtx_t[f] = prob->exactCtx_t[f];
598c371a6d1SMatthew G. Knepley   for (f = Nf; f < NfNew; ++f) tmpexactSol[f] = NULL;
59995cbbfd3SMatthew G. Knepley   for (f = Nf; f < NfNew; ++f) tmpexactCtx[f] = NULL;
600f2cacb80SMatthew G. Knepley   for (f = Nf; f < NfNew; ++f) tmpexactSol_t[f] = NULL;
601f2cacb80SMatthew G. Knepley   for (f = Nf; f < NfNew; ++f) tmpexactCtx_t[f] = NULL;
6029566063dSJacob Faibussowitsch   PetscCall(PetscFree4(prob->exactSol, prob->exactCtx, prob->exactSol_t, prob->exactCtx_t));
603c371a6d1SMatthew G. Knepley   prob->exactSol   = tmpexactSol;
60495cbbfd3SMatthew G. Knepley   prob->exactCtx   = tmpexactCtx;
605f2cacb80SMatthew G. Knepley   prob->exactSol_t = tmpexactSol_t;
606f2cacb80SMatthew G. Knepley   prob->exactCtx_t = tmpexactCtx_t;
6073ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6082764a2aaSMatthew G. Knepley }
6092764a2aaSMatthew G. Knepley 
6102764a2aaSMatthew G. Knepley /*@
61120f4b53cSBarry Smith   PetscDSDestroy - Destroys a `PetscDS` object
6122764a2aaSMatthew G. Knepley 
61320f4b53cSBarry Smith   Collective
6142764a2aaSMatthew G. Knepley 
6152764a2aaSMatthew G. Knepley   Input Parameter:
61660225df5SJacob Faibussowitsch . ds - the `PetscDS` object to destroy
6172764a2aaSMatthew G. Knepley 
6182764a2aaSMatthew G. Knepley   Level: developer
6192764a2aaSMatthew G. Knepley 
620dce8aebaSBarry Smith .seealso: `PetscDSView()`
6212764a2aaSMatthew G. Knepley @*/
622d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDSDestroy(PetscDS *ds)
623d71ae5a4SJacob Faibussowitsch {
6242764a2aaSMatthew G. Knepley   PetscInt f;
6252764a2aaSMatthew G. Knepley 
6262764a2aaSMatthew G. Knepley   PetscFunctionBegin;
6273ba16761SJacob Faibussowitsch   if (!*ds) PetscFunctionReturn(PETSC_SUCCESS);
6286528b96dSMatthew G. Knepley   PetscValidHeaderSpecific((*ds), PETSCDS_CLASSID, 1);
6292764a2aaSMatthew G. Knepley 
6309371c9d4SSatish Balay   if (--((PetscObject)(*ds))->refct > 0) {
6319371c9d4SSatish Balay     *ds = NULL;
6323ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
6339371c9d4SSatish Balay   }
6346528b96dSMatthew G. Knepley   ((PetscObject)(*ds))->refct = 0;
6356528b96dSMatthew G. Knepley   if ((*ds)->subprobs) {
636df3a45bdSMatthew G. Knepley     PetscInt dim, d;
637df3a45bdSMatthew G. Knepley 
6389566063dSJacob Faibussowitsch     PetscCall(PetscDSGetSpatialDimension(*ds, &dim));
6399566063dSJacob Faibussowitsch     for (d = 0; d < dim; ++d) PetscCall(PetscDSDestroy(&(*ds)->subprobs[d]));
640df3a45bdSMatthew G. Knepley   }
6419566063dSJacob Faibussowitsch   PetscCall(PetscFree((*ds)->subprobs));
6429566063dSJacob Faibussowitsch   PetscCall(PetscDSDestroyStructs_Static(*ds));
64348a46eb9SPierre Jolivet   for (f = 0; f < (*ds)->Nf; ++f) PetscCall(PetscObjectDereference((*ds)->disc[f]));
6449566063dSJacob Faibussowitsch   PetscCall(PetscFree4((*ds)->disc, (*ds)->implicit, (*ds)->cohesive, (*ds)->jetDegree));
6459566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormDestroy(&(*ds)->wf));
6469566063dSJacob Faibussowitsch   PetscCall(PetscFree2((*ds)->update, (*ds)->ctx));
6479566063dSJacob Faibussowitsch   PetscCall(PetscFree4((*ds)->exactSol, (*ds)->exactCtx, (*ds)->exactSol_t, (*ds)->exactCtx_t));
648dbbe0bcdSBarry Smith   PetscTryTypeMethod((*ds), destroy);
6499566063dSJacob Faibussowitsch   PetscCall(PetscDSDestroyBoundary(*ds));
6509566063dSJacob Faibussowitsch   PetscCall(PetscFree((*ds)->constants));
6514366bac7SMatthew G. Knepley   for (PetscInt c = 0; c < DM_NUM_POLYTOPES; ++c) {
6524366bac7SMatthew G. Knepley     const PetscInt Na = DMPolytopeTypeGetNumArrangments((DMPolytopeType)c);
6534366bac7SMatthew G. Knepley     if ((*ds)->quadPerm[c])
6544366bac7SMatthew G. Knepley       for (PetscInt o = 0; o < Na; ++o) PetscCall(ISDestroy(&(*ds)->quadPerm[c][o]));
6554366bac7SMatthew G. Knepley     PetscCall(PetscFree((*ds)->quadPerm[c]));
6564366bac7SMatthew G. Knepley     (*ds)->quadPerm[c] = NULL;
6574366bac7SMatthew G. Knepley   }
6589566063dSJacob Faibussowitsch   PetscCall(PetscHeaderDestroy(ds));
6593ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6602764a2aaSMatthew G. Knepley }
6612764a2aaSMatthew G. Knepley 
6622764a2aaSMatthew G. Knepley /*@
663dce8aebaSBarry Smith   PetscDSCreate - Creates an empty `PetscDS` object. The type can then be set with `PetscDSSetType()`.
6642764a2aaSMatthew G. Knepley 
665d083f849SBarry Smith   Collective
6662764a2aaSMatthew G. Knepley 
6672764a2aaSMatthew G. Knepley   Input Parameter:
668dce8aebaSBarry Smith . comm - The communicator for the `PetscDS` object
6692764a2aaSMatthew G. Knepley 
6702764a2aaSMatthew G. Knepley   Output Parameter:
671dce8aebaSBarry Smith . ds - The `PetscDS` object
6722764a2aaSMatthew G. Knepley 
6732764a2aaSMatthew G. Knepley   Level: beginner
6742764a2aaSMatthew G. Knepley 
675dce8aebaSBarry Smith .seealso: `PetscDS`, `PetscDSSetType()`, `PETSCDSBASIC`, `PetscDSType`
6762764a2aaSMatthew G. Knepley @*/
677d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDSCreate(MPI_Comm comm, PetscDS *ds)
678d71ae5a4SJacob Faibussowitsch {
6792764a2aaSMatthew G. Knepley   PetscDS p;
6802764a2aaSMatthew G. Knepley 
6812764a2aaSMatthew G. Knepley   PetscFunctionBegin;
6824f572ea9SToby Isaac   PetscAssertPointer(ds, 2);
6836528b96dSMatthew G. Knepley   *ds = NULL;
6849566063dSJacob Faibussowitsch   PetscCall(PetscDSInitializePackage());
6852764a2aaSMatthew G. Knepley 
6869566063dSJacob Faibussowitsch   PetscCall(PetscHeaderCreate(p, PETSCDS_CLASSID, "PetscDS", "Discrete System", "PetscDS", comm, PetscDSDestroy, PetscDSView));
6872764a2aaSMatthew G. Knepley 
6882764a2aaSMatthew G. Knepley   p->Nf           = 0;
6892764a2aaSMatthew G. Knepley   p->setup        = PETSC_FALSE;
69097b6e6e8SMatthew G. Knepley   p->numConstants = 0;
69197b6e6e8SMatthew G. Knepley   p->constants    = NULL;
692a859676bSMatthew G. Knepley   p->dimEmbed     = -1;
69355c1f793SMatthew G. Knepley   p->useJacPre    = PETSC_TRUE;
69412fc5b22SMatthew G. Knepley   p->forceQuad    = PETSC_TRUE;
6959566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormCreate(comm, &p->wf));
6964366bac7SMatthew G. Knepley   PetscCall(PetscArrayzero(p->quadPerm, DM_NUM_POLYTOPES));
6972764a2aaSMatthew G. Knepley 
6986528b96dSMatthew G. Knepley   *ds = p;
6993ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7002764a2aaSMatthew G. Knepley }
7012764a2aaSMatthew G. Knepley 
702bc4ae4beSMatthew G. Knepley /*@
703dce8aebaSBarry Smith   PetscDSGetNumFields - Returns the number of fields in the `PetscDS`
704bc4ae4beSMatthew G. Knepley 
70520f4b53cSBarry Smith   Not Collective
706bc4ae4beSMatthew G. Knepley 
707bc4ae4beSMatthew G. Knepley   Input Parameter:
70820f4b53cSBarry Smith . prob - The `PetscDS` object
709bc4ae4beSMatthew G. Knepley 
710bc4ae4beSMatthew G. Knepley   Output Parameter:
711bc4ae4beSMatthew G. Knepley . Nf - The number of fields
712bc4ae4beSMatthew G. Knepley 
713bc4ae4beSMatthew G. Knepley   Level: beginner
714bc4ae4beSMatthew G. Knepley 
715dce8aebaSBarry Smith .seealso: `PetscDS`, `PetscDSGetSpatialDimension()`, `PetscDSCreate()`
716bc4ae4beSMatthew G. Knepley @*/
717d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDSGetNumFields(PetscDS prob, PetscInt *Nf)
718d71ae5a4SJacob Faibussowitsch {
7192764a2aaSMatthew G. Knepley   PetscFunctionBegin;
7202764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
7214f572ea9SToby Isaac   PetscAssertPointer(Nf, 2);
7222764a2aaSMatthew G. Knepley   *Nf = prob->Nf;
7233ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7242764a2aaSMatthew G. Knepley }
7252764a2aaSMatthew G. Knepley 
726bc4ae4beSMatthew G. Knepley /*@
727dce8aebaSBarry Smith   PetscDSGetSpatialDimension - Returns the spatial dimension of the `PetscDS`, meaning the topological dimension of the discretizations
728bc4ae4beSMatthew G. Knepley 
72920f4b53cSBarry Smith   Not Collective
730bc4ae4beSMatthew G. Knepley 
731bc4ae4beSMatthew G. Knepley   Input Parameter:
732dce8aebaSBarry Smith . prob - The `PetscDS` object
733bc4ae4beSMatthew G. Knepley 
734bc4ae4beSMatthew G. Knepley   Output Parameter:
735bc4ae4beSMatthew G. Knepley . dim - The spatial dimension
736bc4ae4beSMatthew G. Knepley 
737bc4ae4beSMatthew G. Knepley   Level: beginner
738bc4ae4beSMatthew G. Knepley 
739dce8aebaSBarry Smith .seealso: `PetscDS`, `PetscDSGetCoordinateDimension()`, `PetscDSGetNumFields()`, `PetscDSCreate()`
740bc4ae4beSMatthew G. Knepley @*/
741d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDSGetSpatialDimension(PetscDS prob, PetscInt *dim)
742d71ae5a4SJacob Faibussowitsch {
7432764a2aaSMatthew G. Knepley   PetscFunctionBegin;
7442764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
7454f572ea9SToby Isaac   PetscAssertPointer(dim, 2);
7462764a2aaSMatthew G. Knepley   *dim = 0;
7479de99aefSMatthew G. Knepley   if (prob->Nf) {
7489de99aefSMatthew G. Knepley     PetscObject  obj;
7499de99aefSMatthew G. Knepley     PetscClassId id;
7509de99aefSMatthew G. Knepley 
7519566063dSJacob Faibussowitsch     PetscCall(PetscDSGetDiscretization(prob, 0, &obj));
752665f567fSMatthew G. Knepley     if (obj) {
7539566063dSJacob Faibussowitsch       PetscCall(PetscObjectGetClassId(obj, &id));
7549566063dSJacob Faibussowitsch       if (id == PETSCFE_CLASSID) PetscCall(PetscFEGetSpatialDimension((PetscFE)obj, dim));
7559566063dSJacob Faibussowitsch       else if (id == PETSCFV_CLASSID) PetscCall(PetscFVGetSpatialDimension((PetscFV)obj, dim));
75698921bdaSJacob Faibussowitsch       else SETERRQ(PetscObjectComm((PetscObject)prob), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", 0);
7579de99aefSMatthew G. Knepley     }
758665f567fSMatthew G. Knepley   }
7593ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7602764a2aaSMatthew G. Knepley }
7612764a2aaSMatthew G. Knepley 
762bc4ae4beSMatthew G. Knepley /*@
763dce8aebaSBarry Smith   PetscDSGetCoordinateDimension - Returns the coordinate dimension of the `PetscDS`, meaning the dimension of the space into which the discretiaztions are embedded
764a859676bSMatthew G. Knepley 
76520f4b53cSBarry Smith   Not Collective
766a859676bSMatthew G. Knepley 
767a859676bSMatthew G. Knepley   Input Parameter:
768dce8aebaSBarry Smith . prob - The `PetscDS` object
769a859676bSMatthew G. Knepley 
770a859676bSMatthew G. Knepley   Output Parameter:
771a859676bSMatthew G. Knepley . dimEmbed - The coordinate dimension
772a859676bSMatthew G. Knepley 
773a859676bSMatthew G. Knepley   Level: beginner
774a859676bSMatthew G. Knepley 
775dce8aebaSBarry Smith .seealso: `PetscDS`, `PetscDSSetCoordinateDimension()`, `PetscDSGetSpatialDimension()`, `PetscDSGetNumFields()`, `PetscDSCreate()`
776a859676bSMatthew G. Knepley @*/
777d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDSGetCoordinateDimension(PetscDS prob, PetscInt *dimEmbed)
778d71ae5a4SJacob Faibussowitsch {
779a859676bSMatthew G. Knepley   PetscFunctionBegin;
780a859676bSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
7814f572ea9SToby Isaac   PetscAssertPointer(dimEmbed, 2);
78208401ef6SPierre Jolivet   PetscCheck(prob->dimEmbed >= 0, PetscObjectComm((PetscObject)prob), PETSC_ERR_ARG_WRONGSTATE, "No coordinate dimension set for this DS");
783a859676bSMatthew G. Knepley   *dimEmbed = prob->dimEmbed;
7843ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
785a859676bSMatthew G. Knepley }
786a859676bSMatthew G. Knepley 
787a859676bSMatthew G. Knepley /*@
788dce8aebaSBarry Smith   PetscDSSetCoordinateDimension - Set the coordinate dimension of the `PetscDS`, meaning the dimension of the space into which the discretiaztions are embedded
789a859676bSMatthew G. Knepley 
79020f4b53cSBarry Smith   Logically Collective
791a859676bSMatthew G. Knepley 
792a859676bSMatthew G. Knepley   Input Parameters:
793dce8aebaSBarry Smith + prob     - The `PetscDS` object
794a859676bSMatthew G. Knepley - dimEmbed - The coordinate dimension
795a859676bSMatthew G. Knepley 
796a859676bSMatthew G. Knepley   Level: beginner
797a859676bSMatthew G. Knepley 
798dce8aebaSBarry Smith .seealso: `PetscDS`, `PetscDSGetCoordinateDimension()`, `PetscDSGetSpatialDimension()`, `PetscDSGetNumFields()`, `PetscDSCreate()`
799a859676bSMatthew G. Knepley @*/
800d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDSSetCoordinateDimension(PetscDS prob, PetscInt dimEmbed)
801d71ae5a4SJacob Faibussowitsch {
802a859676bSMatthew G. Knepley   PetscFunctionBegin;
803a859676bSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
80463a3b9bcSJacob Faibussowitsch   PetscCheck(dimEmbed >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Coordinate dimension must be non-negative, not %" PetscInt_FMT, dimEmbed);
805a859676bSMatthew G. Knepley   prob->dimEmbed = dimEmbed;
8063ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
807a859676bSMatthew G. Knepley }
808a859676bSMatthew G. Knepley 
809a859676bSMatthew G. Knepley /*@
81012fc5b22SMatthew G. Knepley   PetscDSGetForceQuad - Returns the flag to force matching quadratures among the field discretizations
81112fc5b22SMatthew G. Knepley 
81212fc5b22SMatthew G. Knepley   Not collective
81312fc5b22SMatthew G. Knepley 
81412fc5b22SMatthew G. Knepley   Input Parameter:
81560225df5SJacob Faibussowitsch . ds - The `PetscDS` object
81612fc5b22SMatthew G. Knepley 
81712fc5b22SMatthew G. Knepley   Output Parameter:
81812fc5b22SMatthew G. Knepley . forceQuad - The flag
81912fc5b22SMatthew G. Knepley 
82012fc5b22SMatthew G. Knepley   Level: intermediate
82112fc5b22SMatthew G. Knepley 
82212fc5b22SMatthew G. Knepley .seealso: `PetscDS`, `PetscDSSetForceQuad()`, `PetscDSGetDiscretization()`, `PetscDSGetNumFields()`, `PetscDSCreate()`
82312fc5b22SMatthew G. Knepley @*/
82412fc5b22SMatthew G. Knepley PetscErrorCode PetscDSGetForceQuad(PetscDS ds, PetscBool *forceQuad)
82512fc5b22SMatthew G. Knepley {
82612fc5b22SMatthew G. Knepley   PetscFunctionBegin;
82712fc5b22SMatthew G. Knepley   PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1);
8284f572ea9SToby Isaac   PetscAssertPointer(forceQuad, 2);
82912fc5b22SMatthew G. Knepley   *forceQuad = ds->forceQuad;
83012fc5b22SMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
83112fc5b22SMatthew G. Knepley }
83212fc5b22SMatthew G. Knepley 
83312fc5b22SMatthew G. Knepley /*@
83412fc5b22SMatthew G. Knepley   PetscDSSetForceQuad - Set the flag to force matching quadratures among the field discretizations
83512fc5b22SMatthew G. Knepley 
83612fc5b22SMatthew G. Knepley   Logically collective on ds
83712fc5b22SMatthew G. Knepley 
83812fc5b22SMatthew G. Knepley   Input Parameters:
83912fc5b22SMatthew G. Knepley + ds        - The `PetscDS` object
84012fc5b22SMatthew G. Knepley - forceQuad - The flag
84112fc5b22SMatthew G. Knepley 
84212fc5b22SMatthew G. Knepley   Level: intermediate
84312fc5b22SMatthew G. Knepley 
84412fc5b22SMatthew G. Knepley .seealso: `PetscDS`, `PetscDSGetForceQuad()`, `PetscDSGetDiscretization()`, `PetscDSGetNumFields()`, `PetscDSCreate()`
84512fc5b22SMatthew G. Knepley @*/
84612fc5b22SMatthew G. Knepley PetscErrorCode PetscDSSetForceQuad(PetscDS ds, PetscBool forceQuad)
84712fc5b22SMatthew G. Knepley {
84812fc5b22SMatthew G. Knepley   PetscFunctionBegin;
84912fc5b22SMatthew G. Knepley   PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1);
85012fc5b22SMatthew G. Knepley   ds->forceQuad = forceQuad;
85112fc5b22SMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
85212fc5b22SMatthew G. Knepley }
85312fc5b22SMatthew G. Knepley 
85412fc5b22SMatthew G. Knepley /*@
855dce8aebaSBarry Smith   PetscDSIsCohesive - Returns the flag indicating that this `PetscDS` is for a cohesive cell
8568edf6225SMatthew G. Knepley 
85720f4b53cSBarry Smith   Not Collective
8588edf6225SMatthew G. Knepley 
8598edf6225SMatthew G. Knepley   Input Parameter:
860dce8aebaSBarry Smith . ds - The `PetscDS` object
8618edf6225SMatthew G. Knepley 
8628edf6225SMatthew G. Knepley   Output Parameter:
8635fedec97SMatthew G. Knepley . isCohesive - The flag
8648edf6225SMatthew G. Knepley 
8658edf6225SMatthew G. Knepley   Level: developer
8668edf6225SMatthew G. Knepley 
867dce8aebaSBarry Smith .seealso: `PetscDS`, `PetscDSGetNumCohesive()`, `PetscDSGetCohesive()`, `PetscDSSetCohesive()`, `PetscDSCreate()`
8688edf6225SMatthew G. Knepley @*/
869d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDSIsCohesive(PetscDS ds, PetscBool *isCohesive)
870d71ae5a4SJacob Faibussowitsch {
8718edf6225SMatthew G. Knepley   PetscFunctionBegin;
8725fedec97SMatthew G. Knepley   PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1);
8734f572ea9SToby Isaac   PetscAssertPointer(isCohesive, 2);
8745fedec97SMatthew G. Knepley   *isCohesive = ds->isCohesive;
8753ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8768edf6225SMatthew G. Knepley }
8778edf6225SMatthew G. Knepley 
8788edf6225SMatthew G. Knepley /*@
879be87f6c0SPierre Jolivet   PetscDSGetNumCohesive - Returns the number of cohesive fields, meaning those defined on the interior of a cohesive cell
8805fedec97SMatthew G. Knepley 
88120f4b53cSBarry Smith   Not Collective
8825fedec97SMatthew G. Knepley 
8835fedec97SMatthew G. Knepley   Input Parameter:
884dce8aebaSBarry Smith . ds - The `PetscDS` object
8855fedec97SMatthew G. Knepley 
8865fedec97SMatthew G. Knepley   Output Parameter:
8875fedec97SMatthew G. Knepley . numCohesive - The number of cohesive fields
8885fedec97SMatthew G. Knepley 
8895fedec97SMatthew G. Knepley   Level: developer
8905fedec97SMatthew G. Knepley 
891dce8aebaSBarry Smith .seealso: `PetscDS`, `PetscDSSetCohesive()`, `PetscDSCreate()`
8925fedec97SMatthew G. Knepley @*/
893d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDSGetNumCohesive(PetscDS ds, PetscInt *numCohesive)
894d71ae5a4SJacob Faibussowitsch {
8955fedec97SMatthew G. Knepley   PetscInt f;
8965fedec97SMatthew G. Knepley 
8975fedec97SMatthew G. Knepley   PetscFunctionBegin;
8985fedec97SMatthew G. Knepley   PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1);
8994f572ea9SToby Isaac   PetscAssertPointer(numCohesive, 2);
9005fedec97SMatthew G. Knepley   *numCohesive = 0;
9015fedec97SMatthew G. Knepley   for (f = 0; f < ds->Nf; ++f) *numCohesive += ds->cohesive[f] ? 1 : 0;
9023ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
9035fedec97SMatthew G. Knepley }
9045fedec97SMatthew G. Knepley 
9055fedec97SMatthew G. Knepley /*@
9065fedec97SMatthew G. Knepley   PetscDSGetCohesive - Returns the flag indicating that a field is cohesive, meaning it is defined on the interior of a cohesive cell
9075fedec97SMatthew G. Knepley 
90820f4b53cSBarry Smith   Not Collective
9095fedec97SMatthew G. Knepley 
910f1a722f8SMatthew G. Knepley   Input Parameters:
911dce8aebaSBarry Smith + ds - The `PetscDS` object
9125fedec97SMatthew G. Knepley - f  - The field index
9135fedec97SMatthew G. Knepley 
9145fedec97SMatthew G. Knepley   Output Parameter:
9155fedec97SMatthew G. Knepley . isCohesive - The flag
9165fedec97SMatthew G. Knepley 
9175fedec97SMatthew G. Knepley   Level: developer
9185fedec97SMatthew G. Knepley 
919dce8aebaSBarry Smith .seealso: `PetscDS`, `PetscDSSetCohesive()`, `PetscDSIsCohesive()`, `PetscDSCreate()`
9205fedec97SMatthew G. Knepley @*/
921d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDSGetCohesive(PetscDS ds, PetscInt f, PetscBool *isCohesive)
922d71ae5a4SJacob Faibussowitsch {
9235fedec97SMatthew G. Knepley   PetscFunctionBegin;
9245fedec97SMatthew G. Knepley   PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1);
9254f572ea9SToby Isaac   PetscAssertPointer(isCohesive, 3);
92663a3b9bcSJacob Faibussowitsch   PetscCheck(!(f < 0) && !(f >= ds->Nf), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %" PetscInt_FMT " must be in [0, %" PetscInt_FMT ")", f, ds->Nf);
9275fedec97SMatthew G. Knepley   *isCohesive = ds->cohesive[f];
9283ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
9295fedec97SMatthew G. Knepley }
9305fedec97SMatthew G. Knepley 
9315fedec97SMatthew G. Knepley /*@
9325fedec97SMatthew G. Knepley   PetscDSSetCohesive - Set the flag indicating that a field is cohesive, meaning it is defined on the interior of a cohesive cell
9338edf6225SMatthew G. Knepley 
93420f4b53cSBarry Smith   Not Collective
9358edf6225SMatthew G. Knepley 
9368edf6225SMatthew G. Knepley   Input Parameters:
937dce8aebaSBarry Smith + ds         - The `PetscDS` object
9385fedec97SMatthew G. Knepley . f          - The field index
9395fedec97SMatthew G. Knepley - isCohesive - The flag for a cohesive field
9408edf6225SMatthew G. Knepley 
9418edf6225SMatthew G. Knepley   Level: developer
9428edf6225SMatthew G. Knepley 
943dce8aebaSBarry Smith .seealso: `PetscDS`, `PetscDSGetCohesive()`, `PetscDSIsCohesive()`, `PetscDSCreate()`
9448edf6225SMatthew G. Knepley @*/
945d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDSSetCohesive(PetscDS ds, PetscInt f, PetscBool isCohesive)
946d71ae5a4SJacob Faibussowitsch {
9475fedec97SMatthew G. Knepley   PetscInt i;
9485fedec97SMatthew G. Knepley 
9498edf6225SMatthew G. Knepley   PetscFunctionBegin;
9505fedec97SMatthew G. Knepley   PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1);
95163a3b9bcSJacob Faibussowitsch   PetscCheck(!(f < 0) && !(f >= ds->Nf), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %" PetscInt_FMT " must be in [0, %" PetscInt_FMT ")", f, ds->Nf);
9525fedec97SMatthew G. Knepley   ds->cohesive[f] = isCohesive;
9535fedec97SMatthew G. Knepley   ds->isCohesive  = PETSC_FALSE;
9545fedec97SMatthew G. Knepley   for (i = 0; i < ds->Nf; ++i) ds->isCohesive = ds->isCohesive || ds->cohesive[f] ? PETSC_TRUE : PETSC_FALSE;
9553ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
9568edf6225SMatthew G. Knepley }
9578edf6225SMatthew G. Knepley 
9588edf6225SMatthew G. Knepley /*@
959bc4ae4beSMatthew G. Knepley   PetscDSGetTotalDimension - Returns the total size of the approximation space for this system
960bc4ae4beSMatthew G. Knepley 
96120f4b53cSBarry Smith   Not Collective
962bc4ae4beSMatthew G. Knepley 
963bc4ae4beSMatthew G. Knepley   Input Parameter:
964dce8aebaSBarry Smith . prob - The `PetscDS` object
965bc4ae4beSMatthew G. Knepley 
966bc4ae4beSMatthew G. Knepley   Output Parameter:
967bc4ae4beSMatthew G. Knepley . dim - The total problem dimension
968bc4ae4beSMatthew G. Knepley 
969bc4ae4beSMatthew G. Knepley   Level: beginner
970bc4ae4beSMatthew G. Knepley 
971dce8aebaSBarry Smith .seealso: `PetscDS`, `PetscDSGetNumFields()`, `PetscDSCreate()`
972bc4ae4beSMatthew G. Knepley @*/
973d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDSGetTotalDimension(PetscDS prob, PetscInt *dim)
974d71ae5a4SJacob Faibussowitsch {
9752764a2aaSMatthew G. Knepley   PetscFunctionBegin;
9762764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
9779566063dSJacob Faibussowitsch   PetscCall(PetscDSSetUp(prob));
9784f572ea9SToby Isaac   PetscAssertPointer(dim, 2);
9792764a2aaSMatthew G. Knepley   *dim = prob->totDim;
9803ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
9812764a2aaSMatthew G. Knepley }
9822764a2aaSMatthew G. Knepley 
983bc4ae4beSMatthew G. Knepley /*@
984bc4ae4beSMatthew G. Knepley   PetscDSGetTotalComponents - Returns the total number of components in this system
985bc4ae4beSMatthew G. Knepley 
98620f4b53cSBarry Smith   Not Collective
987bc4ae4beSMatthew G. Knepley 
988bc4ae4beSMatthew G. Knepley   Input Parameter:
989dce8aebaSBarry Smith . prob - The `PetscDS` object
990bc4ae4beSMatthew G. Knepley 
991bc4ae4beSMatthew G. Knepley   Output Parameter:
99260225df5SJacob Faibussowitsch . Nc - The total number of components
993bc4ae4beSMatthew G. Knepley 
994bc4ae4beSMatthew G. Knepley   Level: beginner
995bc4ae4beSMatthew G. Knepley 
996dce8aebaSBarry Smith .seealso: `PetscDS`, `PetscDSGetNumFields()`, `PetscDSCreate()`
997bc4ae4beSMatthew G. Knepley @*/
998d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDSGetTotalComponents(PetscDS prob, PetscInt *Nc)
999d71ae5a4SJacob Faibussowitsch {
10002764a2aaSMatthew G. Knepley   PetscFunctionBegin;
10012764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
10029566063dSJacob Faibussowitsch   PetscCall(PetscDSSetUp(prob));
10034f572ea9SToby Isaac   PetscAssertPointer(Nc, 2);
10042764a2aaSMatthew G. Knepley   *Nc = prob->totComp;
10053ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
10062764a2aaSMatthew G. Knepley }
10072764a2aaSMatthew G. Knepley 
1008bc4ae4beSMatthew G. Knepley /*@
1009bc4ae4beSMatthew G. Knepley   PetscDSGetDiscretization - Returns the discretization object for the given field
1010bc4ae4beSMatthew G. Knepley 
101120f4b53cSBarry Smith   Not Collective
1012bc4ae4beSMatthew G. Knepley 
1013bc4ae4beSMatthew G. Knepley   Input Parameters:
1014dce8aebaSBarry Smith + prob - The `PetscDS` object
1015bc4ae4beSMatthew G. Knepley - f    - The field number
1016bc4ae4beSMatthew G. Knepley 
1017bc4ae4beSMatthew G. Knepley   Output Parameter:
1018bc4ae4beSMatthew G. Knepley . disc - The discretization object
1019bc4ae4beSMatthew G. Knepley 
1020bc4ae4beSMatthew G. Knepley   Level: beginner
1021bc4ae4beSMatthew G. Knepley 
1022dce8aebaSBarry Smith .seealso: `PetscDS`, `PetscFE`, `PetscFV`, `PetscDSSetDiscretization()`, `PetscDSAddDiscretization()`, `PetscDSGetNumFields()`, `PetscDSCreate()`
1023bc4ae4beSMatthew G. Knepley @*/
1024d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDSGetDiscretization(PetscDS prob, PetscInt f, PetscObject *disc)
1025d71ae5a4SJacob Faibussowitsch {
10266528b96dSMatthew G. Knepley   PetscFunctionBeginHot;
10272764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
10284f572ea9SToby Isaac   PetscAssertPointer(disc, 3);
102963a3b9bcSJacob Faibussowitsch   PetscCheck(!(f < 0) && !(f >= prob->Nf), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %" PetscInt_FMT " must be in [0, %" PetscInt_FMT ")", f, prob->Nf);
10302764a2aaSMatthew G. Knepley   *disc = prob->disc[f];
10313ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
10322764a2aaSMatthew G. Knepley }
10332764a2aaSMatthew G. Knepley 
1034bc4ae4beSMatthew G. Knepley /*@
1035bc4ae4beSMatthew G. Knepley   PetscDSSetDiscretization - Sets the discretization object for the given field
1036bc4ae4beSMatthew G. Knepley 
103720f4b53cSBarry Smith   Not Collective
1038bc4ae4beSMatthew G. Knepley 
1039bc4ae4beSMatthew G. Knepley   Input Parameters:
1040dce8aebaSBarry Smith + prob - The `PetscDS` object
1041bc4ae4beSMatthew G. Knepley . f    - The field number
1042bc4ae4beSMatthew G. Knepley - disc - The discretization object
1043bc4ae4beSMatthew G. Knepley 
1044bc4ae4beSMatthew G. Knepley   Level: beginner
1045bc4ae4beSMatthew G. Knepley 
1046dce8aebaSBarry Smith .seealso: `PetscDS`, `PetscFE`, `PetscFV`, `PetscDSGetDiscretization()`, `PetscDSAddDiscretization()`, `PetscDSGetNumFields()`, `PetscDSCreate()`
1047bc4ae4beSMatthew G. Knepley @*/
1048d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDSSetDiscretization(PetscDS prob, PetscInt f, PetscObject disc)
1049d71ae5a4SJacob Faibussowitsch {
10502764a2aaSMatthew G. Knepley   PetscFunctionBegin;
10512764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
10524f572ea9SToby Isaac   if (disc) PetscAssertPointer(disc, 3);
105363a3b9bcSJacob Faibussowitsch   PetscCheck(f >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %" PetscInt_FMT " must be non-negative", f);
10549566063dSJacob Faibussowitsch   PetscCall(PetscDSEnlarge_Static(prob, f + 1));
10559566063dSJacob Faibussowitsch   PetscCall(PetscObjectDereference(prob->disc[f]));
10562764a2aaSMatthew G. Knepley   prob->disc[f] = disc;
10579566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference(disc));
1058665f567fSMatthew G. Knepley   if (disc) {
1059249df284SMatthew G. Knepley     PetscClassId id;
1060249df284SMatthew G. Knepley 
10619566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetClassId(disc, &id));
10621cf84007SMatthew G. Knepley     if (id == PETSCFE_CLASSID) {
10639566063dSJacob Faibussowitsch       PetscCall(PetscDSSetImplicit(prob, f, PETSC_TRUE));
10641cf84007SMatthew G. Knepley     } else if (id == PETSCFV_CLASSID) {
10659566063dSJacob Faibussowitsch       PetscCall(PetscDSSetImplicit(prob, f, PETSC_FALSE));
1066a6cbbb48SMatthew G. Knepley     }
10679566063dSJacob Faibussowitsch     PetscCall(PetscDSSetJetDegree(prob, f, 1));
1068249df284SMatthew G. Knepley   }
10693ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
10702764a2aaSMatthew G. Knepley }
10712764a2aaSMatthew G. Knepley 
1072bc4ae4beSMatthew G. Knepley /*@
10736528b96dSMatthew G. Knepley   PetscDSGetWeakForm - Returns the weak form object
10746528b96dSMatthew G. Knepley 
107520f4b53cSBarry Smith   Not Collective
10766528b96dSMatthew G. Knepley 
10776528b96dSMatthew G. Knepley   Input Parameter:
1078dce8aebaSBarry Smith . ds - The `PetscDS` object
10796528b96dSMatthew G. Knepley 
10806528b96dSMatthew G. Knepley   Output Parameter:
10816528b96dSMatthew G. Knepley . wf - The weak form object
10826528b96dSMatthew G. Knepley 
10836528b96dSMatthew G. Knepley   Level: beginner
10846528b96dSMatthew G. Knepley 
1085dce8aebaSBarry Smith .seealso: `PetscWeakForm`, `PetscDSSetWeakForm()`, `PetscDSGetNumFields()`, `PetscDSCreate()`
10866528b96dSMatthew G. Knepley @*/
1087d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDSGetWeakForm(PetscDS ds, PetscWeakForm *wf)
1088d71ae5a4SJacob Faibussowitsch {
10896528b96dSMatthew G. Knepley   PetscFunctionBegin;
10906528b96dSMatthew G. Knepley   PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1);
10914f572ea9SToby Isaac   PetscAssertPointer(wf, 2);
10926528b96dSMatthew G. Knepley   *wf = ds->wf;
10933ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
10946528b96dSMatthew G. Knepley }
10956528b96dSMatthew G. Knepley 
10966528b96dSMatthew G. Knepley /*@
10976528b96dSMatthew G. Knepley   PetscDSSetWeakForm - Sets the weak form object
10986528b96dSMatthew G. Knepley 
109920f4b53cSBarry Smith   Not Collective
11006528b96dSMatthew G. Knepley 
11016528b96dSMatthew G. Knepley   Input Parameters:
1102dce8aebaSBarry Smith + ds - The `PetscDS` object
11036528b96dSMatthew G. Knepley - wf - The weak form object
11046528b96dSMatthew G. Knepley 
11056528b96dSMatthew G. Knepley   Level: beginner
11066528b96dSMatthew G. Knepley 
1107dce8aebaSBarry Smith .seealso: `PetscWeakForm`, `PetscDSGetWeakForm()`, `PetscDSGetNumFields()`, `PetscDSCreate()`
11086528b96dSMatthew G. Knepley @*/
1109d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDSSetWeakForm(PetscDS ds, PetscWeakForm wf)
1110d71ae5a4SJacob Faibussowitsch {
11116528b96dSMatthew G. Knepley   PetscFunctionBegin;
11126528b96dSMatthew G. Knepley   PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1);
11136528b96dSMatthew G. Knepley   PetscValidHeaderSpecific(wf, PETSCWEAKFORM_CLASSID, 2);
11149566063dSJacob Faibussowitsch   PetscCall(PetscObjectDereference((PetscObject)ds->wf));
11156528b96dSMatthew G. Knepley   ds->wf = wf;
11169566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)wf));
11179566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormSetNumFields(wf, ds->Nf));
11183ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
11196528b96dSMatthew G. Knepley }
11206528b96dSMatthew G. Knepley 
11216528b96dSMatthew G. Knepley /*@
1122bc4ae4beSMatthew G. Knepley   PetscDSAddDiscretization - Adds a discretization object
1123bc4ae4beSMatthew G. Knepley 
112420f4b53cSBarry Smith   Not Collective
1125bc4ae4beSMatthew G. Knepley 
1126bc4ae4beSMatthew G. Knepley   Input Parameters:
1127dce8aebaSBarry Smith + prob - The `PetscDS` object
1128bc4ae4beSMatthew G. Knepley - disc - The boundary discretization object
1129bc4ae4beSMatthew G. Knepley 
1130bc4ae4beSMatthew G. Knepley   Level: beginner
1131bc4ae4beSMatthew G. Knepley 
1132dce8aebaSBarry Smith .seealso: `PetscWeakForm`, `PetscDSGetDiscretization()`, `PetscDSSetDiscretization()`, `PetscDSGetNumFields()`, `PetscDSCreate()`
1133bc4ae4beSMatthew G. Knepley @*/
1134d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDSAddDiscretization(PetscDS prob, PetscObject disc)
1135d71ae5a4SJacob Faibussowitsch {
11362764a2aaSMatthew G. Knepley   PetscFunctionBegin;
11379566063dSJacob Faibussowitsch   PetscCall(PetscDSSetDiscretization(prob, prob->Nf, disc));
11383ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
11392764a2aaSMatthew G. Knepley }
11402764a2aaSMatthew G. Knepley 
1141249df284SMatthew G. Knepley /*@
1142dce8aebaSBarry Smith   PetscDSGetQuadrature - Returns the quadrature, which must agree for all fields in the `PetscDS`
1143083401c6SMatthew G. Knepley 
114420f4b53cSBarry Smith   Not Collective
1145083401c6SMatthew G. Knepley 
1146083401c6SMatthew G. Knepley   Input Parameter:
1147dce8aebaSBarry Smith . prob - The `PetscDS` object
1148083401c6SMatthew G. Knepley 
1149083401c6SMatthew G. Knepley   Output Parameter:
1150083401c6SMatthew G. Knepley . q - The quadrature object
1151083401c6SMatthew G. Knepley 
1152083401c6SMatthew G. Knepley   Level: intermediate
1153083401c6SMatthew G. Knepley 
1154dce8aebaSBarry Smith .seealso: `PetscDS`, `PetscQuadrature`, `PetscDSSetImplicit()`, `PetscDSSetDiscretization()`, `PetscDSAddDiscretization()`, `PetscDSGetNumFields()`, `PetscDSCreate()`
1155083401c6SMatthew G. Knepley @*/
1156d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDSGetQuadrature(PetscDS prob, PetscQuadrature *q)
1157d71ae5a4SJacob Faibussowitsch {
1158083401c6SMatthew G. Knepley   PetscObject  obj;
1159083401c6SMatthew G. Knepley   PetscClassId id;
1160083401c6SMatthew G. Knepley 
1161083401c6SMatthew G. Knepley   PetscFunctionBegin;
1162083401c6SMatthew G. Knepley   *q = NULL;
11633ba16761SJacob Faibussowitsch   if (!prob->Nf) PetscFunctionReturn(PETSC_SUCCESS);
11649566063dSJacob Faibussowitsch   PetscCall(PetscDSGetDiscretization(prob, 0, &obj));
11659566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetClassId(obj, &id));
11669566063dSJacob Faibussowitsch   if (id == PETSCFE_CLASSID) PetscCall(PetscFEGetQuadrature((PetscFE)obj, q));
11679566063dSJacob Faibussowitsch   else if (id == PETSCFV_CLASSID) PetscCall(PetscFVGetQuadrature((PetscFV)obj, q));
116898921bdaSJacob Faibussowitsch   else SETERRQ(PetscObjectComm((PetscObject)prob), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", 0);
11693ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1170083401c6SMatthew G. Knepley }
1171083401c6SMatthew G. Knepley 
1172083401c6SMatthew G. Knepley /*@
1173dce8aebaSBarry Smith   PetscDSGetImplicit - Returns the flag for implicit solve for this field. This is just a guide for `TSIMEX`
1174249df284SMatthew G. Knepley 
117520f4b53cSBarry Smith   Not Collective
1176249df284SMatthew G. Knepley 
1177249df284SMatthew G. Knepley   Input Parameters:
1178dce8aebaSBarry Smith + prob - The `PetscDS` object
1179249df284SMatthew G. Knepley - f    - The field number
1180249df284SMatthew G. Knepley 
1181249df284SMatthew G. Knepley   Output Parameter:
1182249df284SMatthew G. Knepley . implicit - The flag indicating what kind of solve to use for this field
1183249df284SMatthew G. Knepley 
1184249df284SMatthew G. Knepley   Level: developer
1185249df284SMatthew G. Knepley 
1186dce8aebaSBarry Smith .seealso: `TSIMEX`, `PetscDS`, `PetscDSSetImplicit()`, `PetscDSSetDiscretization()`, `PetscDSAddDiscretization()`, `PetscDSGetNumFields()`, `PetscDSCreate()`
1187249df284SMatthew G. Knepley @*/
1188d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDSGetImplicit(PetscDS prob, PetscInt f, PetscBool *implicit)
1189d71ae5a4SJacob Faibussowitsch {
1190249df284SMatthew G. Knepley   PetscFunctionBegin;
1191249df284SMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
11924f572ea9SToby Isaac   PetscAssertPointer(implicit, 3);
119363a3b9bcSJacob Faibussowitsch   PetscCheck(!(f < 0) && !(f >= prob->Nf), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %" PetscInt_FMT " must be in [0, %" PetscInt_FMT ")", f, prob->Nf);
1194249df284SMatthew G. Knepley   *implicit = prob->implicit[f];
11953ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1196249df284SMatthew G. Knepley }
1197249df284SMatthew G. Knepley 
1198249df284SMatthew G. Knepley /*@
1199dce8aebaSBarry Smith   PetscDSSetImplicit - Set the flag for implicit solve for this field. This is just a guide for `TSIMEX`
1200249df284SMatthew G. Knepley 
120120f4b53cSBarry Smith   Not Collective
1202249df284SMatthew G. Knepley 
1203249df284SMatthew G. Knepley   Input Parameters:
1204dce8aebaSBarry Smith + prob     - The `PetscDS` object
1205249df284SMatthew G. Knepley . f        - The field number
1206249df284SMatthew G. Knepley - implicit - The flag indicating what kind of solve to use for this field
1207249df284SMatthew G. Knepley 
1208249df284SMatthew G. Knepley   Level: developer
1209249df284SMatthew G. Knepley 
1210dce8aebaSBarry Smith .seealso: `TSIMEX`, `PetscDSGetImplicit()`, `PetscDSSetDiscretization()`, `PetscDSAddDiscretization()`, `PetscDSGetNumFields()`, `PetscDSCreate()`
1211249df284SMatthew G. Knepley @*/
1212d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDSSetImplicit(PetscDS prob, PetscInt f, PetscBool implicit)
1213d71ae5a4SJacob Faibussowitsch {
1214249df284SMatthew G. Knepley   PetscFunctionBegin;
1215249df284SMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
121663a3b9bcSJacob Faibussowitsch   PetscCheck(!(f < 0) && !(f >= prob->Nf), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %" PetscInt_FMT " must be in [0, %" PetscInt_FMT ")", f, prob->Nf);
1217249df284SMatthew G. Knepley   prob->implicit[f] = implicit;
12183ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1219249df284SMatthew G. Knepley }
1220249df284SMatthew G. Knepley 
1221f9244615SMatthew G. Knepley /*@
1222f9244615SMatthew G. Knepley   PetscDSGetJetDegree - Returns the highest derivative for this field equation, or the k-jet that the discretization needs to tabulate.
1223f9244615SMatthew G. Knepley 
122420f4b53cSBarry Smith   Not Collective
1225f9244615SMatthew G. Knepley 
1226f9244615SMatthew G. Knepley   Input Parameters:
1227dce8aebaSBarry Smith + ds - The `PetscDS` object
1228f9244615SMatthew G. Knepley - f  - The field number
1229f9244615SMatthew G. Knepley 
1230f9244615SMatthew G. Knepley   Output Parameter:
1231f9244615SMatthew G. Knepley . k - The highest derivative we need to tabulate
1232f9244615SMatthew G. Knepley 
1233f9244615SMatthew G. Knepley   Level: developer
1234f9244615SMatthew G. Knepley 
1235dce8aebaSBarry Smith .seealso: `PetscDS`, `PetscDSSetJetDegree()`, `PetscDSSetDiscretization()`, `PetscDSAddDiscretization()`, `PetscDSGetNumFields()`, `PetscDSCreate()`
1236f9244615SMatthew G. Knepley @*/
1237d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDSGetJetDegree(PetscDS ds, PetscInt f, PetscInt *k)
1238d71ae5a4SJacob Faibussowitsch {
1239f9244615SMatthew G. Knepley   PetscFunctionBegin;
1240f9244615SMatthew G. Knepley   PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1);
12414f572ea9SToby Isaac   PetscAssertPointer(k, 3);
124263a3b9bcSJacob Faibussowitsch   PetscCheck(!(f < 0) && !(f >= ds->Nf), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %" PetscInt_FMT " must be in [0, %" PetscInt_FMT ")", f, ds->Nf);
1243f9244615SMatthew G. Knepley   *k = ds->jetDegree[f];
12443ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1245f9244615SMatthew G. Knepley }
1246f9244615SMatthew G. Knepley 
1247f9244615SMatthew G. Knepley /*@
1248f9244615SMatthew G. Knepley   PetscDSSetJetDegree - Set the highest derivative for this field equation, or the k-jet that the discretization needs to tabulate.
1249f9244615SMatthew G. Knepley 
125020f4b53cSBarry Smith   Not Collective
1251f9244615SMatthew G. Knepley 
1252f9244615SMatthew G. Knepley   Input Parameters:
1253dce8aebaSBarry Smith + ds - The `PetscDS` object
1254f9244615SMatthew G. Knepley . f  - The field number
1255f9244615SMatthew G. Knepley - k  - The highest derivative we need to tabulate
1256f9244615SMatthew G. Knepley 
1257f9244615SMatthew G. Knepley   Level: developer
1258f9244615SMatthew G. Knepley 
125960225df5SJacob Faibussowitsch .seealso: ``PetscDS`, `PetscDSGetJetDegree()`, `PetscDSSetDiscretization()`, `PetscDSAddDiscretization()`, `PetscDSGetNumFields()`, `PetscDSCreate()`
1260f9244615SMatthew G. Knepley @*/
1261d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDSSetJetDegree(PetscDS ds, PetscInt f, PetscInt k)
1262d71ae5a4SJacob Faibussowitsch {
1263f9244615SMatthew G. Knepley   PetscFunctionBegin;
1264f9244615SMatthew G. Knepley   PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1);
126563a3b9bcSJacob Faibussowitsch   PetscCheck(!(f < 0) && !(f >= ds->Nf), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %" PetscInt_FMT " must be in [0, %" PetscInt_FMT ")", f, ds->Nf);
1266f9244615SMatthew G. Knepley   ds->jetDegree[f] = k;
12673ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1268f9244615SMatthew G. Knepley }
1269f9244615SMatthew G. Knepley 
1270c8943706SMatthew G. Knepley /*@C
1271c8943706SMatthew G. Knepley   PetscDSGetObjective - Get the pointwise objective function for a given test field
1272c8943706SMatthew G. Knepley 
1273c8943706SMatthew G. Knepley   Not Collective
1274c8943706SMatthew G. Knepley 
1275c8943706SMatthew G. Knepley   Input Parameters:
1276c8943706SMatthew G. Knepley + ds - The `PetscDS`
1277c8943706SMatthew G. Knepley - f  - The test field number
1278c8943706SMatthew G. Knepley 
1279*a4e35b19SJacob Faibussowitsch   Output Parameter:
1280c8943706SMatthew G. Knepley . obj - integrand for the test function term
1281c8943706SMatthew G. Knepley 
1282c8943706SMatthew G. Knepley   Calling sequence of `obj`:
1283c8943706SMatthew G. Knepley + dim          - the spatial dimension
1284c8943706SMatthew G. Knepley . Nf           - the number of fields
1285*a4e35b19SJacob Faibussowitsch . NfAux        - the number of auxiliary fields
1286c8943706SMatthew G. Knepley . uOff         - the offset into u[] and u_t[] for each field
1287c8943706SMatthew G. Knepley . uOff_x       - the offset into u_x[] for each field
1288c8943706SMatthew G. Knepley . u            - each field evaluated at the current point
1289c8943706SMatthew G. Knepley . u_t          - the time derivative of each field evaluated at the current point
1290c8943706SMatthew G. Knepley . u_x          - the gradient of each field evaluated at the current point
1291c8943706SMatthew G. Knepley . aOff         - the offset into a[] and a_t[] for each auxiliary field
1292c8943706SMatthew G. Knepley . aOff_x       - the offset into a_x[] for each auxiliary field
1293c8943706SMatthew G. Knepley . a            - each auxiliary field evaluated at the current point
1294c8943706SMatthew G. Knepley . a_t          - the time derivative of each auxiliary field evaluated at the current point
1295c8943706SMatthew G. Knepley . a_x          - the gradient of auxiliary each field evaluated at the current point
1296c8943706SMatthew G. Knepley . t            - current time
1297c8943706SMatthew G. Knepley . x            - coordinates of the current point
1298c8943706SMatthew G. Knepley . numConstants - number of constant parameters
1299c8943706SMatthew G. Knepley . constants    - constant parameters
1300c8943706SMatthew G. Knepley - obj          - output values at the current point
1301c8943706SMatthew G. Knepley 
1302c8943706SMatthew G. Knepley   Level: intermediate
1303c8943706SMatthew G. Knepley 
1304c8943706SMatthew G. Knepley   Note:
1305c8943706SMatthew G. Knepley   We are using a first order FEM model for the weak form:  \int_\Omega \phi obj(u, u_t, \nabla u, x, t)
1306c8943706SMatthew G. Knepley 
1307c8943706SMatthew G. Knepley .seealso: `PetscDS`, `PetscDSSetObjective()`, `PetscDSGetResidual()`
1308c8943706SMatthew G. Knepley @*/
1309d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDSGetObjective(PetscDS ds, PetscInt f, void (**obj)(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 obj[]))
1310d71ae5a4SJacob Faibussowitsch {
13116528b96dSMatthew G. Knepley   PetscPointFunc *tmp;
13126528b96dSMatthew G. Knepley   PetscInt        n;
13136528b96dSMatthew G. Knepley 
13142764a2aaSMatthew G. Knepley   PetscFunctionBegin;
13156528b96dSMatthew G. Knepley   PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1);
13164f572ea9SToby Isaac   PetscAssertPointer(obj, 3);
131763a3b9bcSJacob Faibussowitsch   PetscCheck(!(f < 0) && !(f >= ds->Nf), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %" PetscInt_FMT " must be in [0, %" PetscInt_FMT ")", f, ds->Nf);
13189566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormGetObjective(ds->wf, NULL, 0, f, 0, &n, &tmp));
13196528b96dSMatthew G. Knepley   *obj = tmp ? tmp[0] : NULL;
13203ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
13212764a2aaSMatthew G. Knepley }
13222764a2aaSMatthew G. Knepley 
1323c8943706SMatthew G. Knepley /*@C
1324c8943706SMatthew G. Knepley   PetscDSSetObjective - Set the pointwise objective function for a given test field
1325c8943706SMatthew G. Knepley 
1326c8943706SMatthew G. Knepley   Not Collective
1327c8943706SMatthew G. Knepley 
1328c8943706SMatthew G. Knepley   Input Parameters:
1329c8943706SMatthew G. Knepley + ds  - The `PetscDS`
1330c8943706SMatthew G. Knepley . f   - The test field number
1331c8943706SMatthew G. Knepley - obj - integrand for the test function term
1332c8943706SMatthew G. Knepley 
1333c8943706SMatthew G. Knepley   Calling sequence of `obj`:
1334c8943706SMatthew G. Knepley + dim          - the spatial dimension
1335c8943706SMatthew G. Knepley . Nf           - the number of fields
1336*a4e35b19SJacob Faibussowitsch . NfAux        - the number of auxiliary fields
1337c8943706SMatthew G. Knepley . uOff         - the offset into u[] and u_t[] for each field
1338c8943706SMatthew G. Knepley . uOff_x       - the offset into u_x[] for each field
1339c8943706SMatthew G. Knepley . u            - each field evaluated at the current point
1340c8943706SMatthew G. Knepley . u_t          - the time derivative of each field evaluated at the current point
1341c8943706SMatthew G. Knepley . u_x          - the gradient of each field evaluated at the current point
1342c8943706SMatthew G. Knepley . aOff         - the offset into a[] and a_t[] for each auxiliary field
1343c8943706SMatthew G. Knepley . aOff_x       - the offset into a_x[] for each auxiliary field
1344c8943706SMatthew G. Knepley . a            - each auxiliary field evaluated at the current point
1345c8943706SMatthew G. Knepley . a_t          - the time derivative of each auxiliary field evaluated at the current point
1346c8943706SMatthew G. Knepley . a_x          - the gradient of auxiliary each field evaluated at the current point
1347c8943706SMatthew G. Knepley . t            - current time
1348c8943706SMatthew G. Knepley . x            - coordinates of the current point
1349c8943706SMatthew G. Knepley . numConstants - number of constant parameters
1350c8943706SMatthew G. Knepley . constants    - constant parameters
1351c8943706SMatthew G. Knepley - obj          - output values at the current point
1352c8943706SMatthew G. Knepley 
1353c8943706SMatthew G. Knepley   Level: intermediate
1354c8943706SMatthew G. Knepley 
1355c8943706SMatthew G. Knepley   Note:
1356c8943706SMatthew G. Knepley   We are using a first order FEM model for the weak form:  \int_\Omega \phi obj(u, u_t, \nabla u, x, t)
1357c8943706SMatthew G. Knepley 
1358c8943706SMatthew G. Knepley .seealso: `PetscDS`, `PetscDSGetObjective()`, `PetscDSSetResidual()`
1359c8943706SMatthew G. Knepley @*/
1360d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDSSetObjective(PetscDS ds, PetscInt f, void (*obj)(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 obj[]))
1361d71ae5a4SJacob Faibussowitsch {
13622764a2aaSMatthew G. Knepley   PetscFunctionBegin;
13636528b96dSMatthew G. Knepley   PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1);
13646528b96dSMatthew G. Knepley   if (obj) PetscValidFunction(obj, 3);
136563a3b9bcSJacob Faibussowitsch   PetscCheck(f >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %" PetscInt_FMT " must be non-negative", f);
13669566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormSetIndexObjective(ds->wf, NULL, 0, f, 0, 0, obj));
13673ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
13682764a2aaSMatthew G. Knepley }
13692764a2aaSMatthew G. Knepley 
1370194d53e6SMatthew G. Knepley /*@C
1371194d53e6SMatthew G. Knepley   PetscDSGetResidual - Get the pointwise residual function for a given test field
1372194d53e6SMatthew G. Knepley 
137320f4b53cSBarry Smith   Not Collective
1374194d53e6SMatthew G. Knepley 
1375194d53e6SMatthew G. Knepley   Input Parameters:
1376dce8aebaSBarry Smith + ds - The `PetscDS`
1377194d53e6SMatthew G. Knepley - f  - The test field number
1378194d53e6SMatthew G. Knepley 
1379194d53e6SMatthew G. Knepley   Output Parameters:
1380194d53e6SMatthew G. Knepley + f0 - integrand for the test function term
1381194d53e6SMatthew G. Knepley - f1 - integrand for the test function gradient term
1382194d53e6SMatthew G. Knepley 
1383*a4e35b19SJacob Faibussowitsch   Calling sequence of `f0`:
1384194d53e6SMatthew G. Knepley + dim          - the spatial dimension
1385194d53e6SMatthew G. Knepley . Nf           - the number of fields
1386*a4e35b19SJacob Faibussowitsch . NfAux        - the number of auxiliary fields
1387194d53e6SMatthew G. Knepley . uOff         - the offset into u[] and u_t[] for each field
1388194d53e6SMatthew G. Knepley . uOff_x       - the offset into u_x[] for each field
1389194d53e6SMatthew G. Knepley . u            - each field evaluated at the current point
1390194d53e6SMatthew G. Knepley . u_t          - the time derivative of each field evaluated at the current point
1391194d53e6SMatthew G. Knepley . u_x          - the gradient of each field evaluated at the current point
1392194d53e6SMatthew G. Knepley . aOff         - the offset into a[] and a_t[] for each auxiliary field
1393194d53e6SMatthew G. Knepley . aOff_x       - the offset into a_x[] for each auxiliary field
1394194d53e6SMatthew G. Knepley . a            - each auxiliary field evaluated at the current point
1395194d53e6SMatthew G. Knepley . a_t          - the time derivative of each auxiliary field evaluated at the current point
1396194d53e6SMatthew G. Knepley . a_x          - the gradient of auxiliary each field evaluated at the current point
1397194d53e6SMatthew G. Knepley . t            - current time
1398194d53e6SMatthew G. Knepley . x            - coordinates of the current point
139997b6e6e8SMatthew G. Knepley . numConstants - number of constant parameters
140097b6e6e8SMatthew G. Knepley . constants    - constant parameters
1401194d53e6SMatthew G. Knepley - f0           - output values at the current point
1402194d53e6SMatthew G. Knepley 
1403194d53e6SMatthew G. Knepley   Level: intermediate
1404194d53e6SMatthew G. Knepley 
1405dce8aebaSBarry Smith   Note:
1406*a4e35b19SJacob Faibussowitsch   `f1` has an identical form and is omitted for brevity.
1407*a4e35b19SJacob Faibussowitsch 
1408dce8aebaSBarry Smith   We are using a first order FEM model for the weak form:  \int_\Omega \phi f_0(u, u_t, \nabla u, x, t) + \nabla\phi \cdot {\vec f}_1(u, u_t, \nabla u, x, t)
1409dce8aebaSBarry Smith 
1410dce8aebaSBarry Smith .seealso: `PetscDS`, `PetscDSSetResidual()`
1411194d53e6SMatthew G. Knepley @*/
1412*a4e35b19SJacob Faibussowitsch PetscErrorCode PetscDSGetResidual(PetscDS ds, 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[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]), void (**f1)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
1413d71ae5a4SJacob Faibussowitsch {
14146528b96dSMatthew G. Knepley   PetscPointFunc *tmp0, *tmp1;
14156528b96dSMatthew G. Knepley   PetscInt        n0, n1;
14166528b96dSMatthew G. Knepley 
14172764a2aaSMatthew G. Knepley   PetscFunctionBegin;
14186528b96dSMatthew G. Knepley   PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1);
141963a3b9bcSJacob Faibussowitsch   PetscCheck(!(f < 0) && !(f >= ds->Nf), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %" PetscInt_FMT " must be in [0, %" PetscInt_FMT ")", f, ds->Nf);
14209566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormGetResidual(ds->wf, NULL, 0, f, 0, &n0, &tmp0, &n1, &tmp1));
14216528b96dSMatthew G. Knepley   *f0 = tmp0 ? tmp0[0] : NULL;
14226528b96dSMatthew G. Knepley   *f1 = tmp1 ? tmp1[0] : NULL;
14233ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
14242764a2aaSMatthew G. Knepley }
14252764a2aaSMatthew G. Knepley 
1426194d53e6SMatthew G. Knepley /*@C
1427194d53e6SMatthew G. Knepley   PetscDSSetResidual - Set the pointwise residual function for a given test field
1428194d53e6SMatthew G. Knepley 
142920f4b53cSBarry Smith   Not Collective
1430194d53e6SMatthew G. Knepley 
1431194d53e6SMatthew G. Knepley   Input Parameters:
1432dce8aebaSBarry Smith + ds - The `PetscDS`
1433194d53e6SMatthew G. Knepley . f  - The test field number
1434194d53e6SMatthew G. Knepley . f0 - integrand for the test function term
1435194d53e6SMatthew G. Knepley - f1 - integrand for the test function gradient term
1436194d53e6SMatthew G. Knepley 
1437*a4e35b19SJacob Faibussowitsch   Calling sequence of `f0`:
1438194d53e6SMatthew G. Knepley + dim          - the spatial dimension
1439194d53e6SMatthew G. Knepley . Nf           - the number of fields
1440*a4e35b19SJacob Faibussowitsch . NfAux        - the number of auxiliary fields
1441194d53e6SMatthew G. Knepley . uOff         - the offset into u[] and u_t[] for each field
1442194d53e6SMatthew G. Knepley . uOff_x       - the offset into u_x[] for each field
1443194d53e6SMatthew G. Knepley . u            - each field evaluated at the current point
1444194d53e6SMatthew G. Knepley . u_t          - the time derivative of each field evaluated at the current point
1445194d53e6SMatthew G. Knepley . u_x          - the gradient of each field evaluated at the current point
1446194d53e6SMatthew G. Knepley . aOff         - the offset into a[] and a_t[] for each auxiliary field
1447194d53e6SMatthew G. Knepley . aOff_x       - the offset into a_x[] for each auxiliary field
1448194d53e6SMatthew G. Knepley . a            - each auxiliary field evaluated at the current point
1449194d53e6SMatthew G. Knepley . a_t          - the time derivative of each auxiliary field evaluated at the current point
1450194d53e6SMatthew G. Knepley . a_x          - the gradient of auxiliary each field evaluated at the current point
1451194d53e6SMatthew G. Knepley . t            - current time
1452194d53e6SMatthew G. Knepley . x            - coordinates of the current point
145397b6e6e8SMatthew G. Knepley . numConstants - number of constant parameters
145497b6e6e8SMatthew G. Knepley . constants    - constant parameters
1455194d53e6SMatthew G. Knepley - f0           - output values at the current point
1456194d53e6SMatthew G. Knepley 
1457194d53e6SMatthew G. Knepley   Level: intermediate
1458194d53e6SMatthew G. Knepley 
1459dce8aebaSBarry Smith   Note:
1460*a4e35b19SJacob Faibussowitsch   `f1` has an identical form and is omitted for brevity.
1461*a4e35b19SJacob Faibussowitsch 
1462dce8aebaSBarry Smith   We are using a first order FEM model for the weak form:  \int_\Omega \phi f_0(u, u_t, \nabla u, x, t) + \nabla\phi \cdot {\vec f}_1(u, u_t, \nabla u, x, t)
1463dce8aebaSBarry Smith 
1464dce8aebaSBarry Smith .seealso: `PetscDS`, `PetscDSGetResidual()`
1465194d53e6SMatthew G. Knepley @*/
1466*a4e35b19SJacob Faibussowitsch PetscErrorCode PetscDSSetResidual(PetscDS ds, 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[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]), void (*f1)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
1467d71ae5a4SJacob Faibussowitsch {
14682764a2aaSMatthew G. Knepley   PetscFunctionBegin;
14696528b96dSMatthew G. Knepley   PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1);
1470f866a1d0SMatthew G. Knepley   if (f0) PetscValidFunction(f0, 3);
1471f866a1d0SMatthew G. Knepley   if (f1) PetscValidFunction(f1, 4);
147263a3b9bcSJacob Faibussowitsch   PetscCheck(f >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %" PetscInt_FMT " must be non-negative", f);
14739566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormSetIndexResidual(ds->wf, NULL, 0, f, 0, 0, f0, 0, f1));
14743ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
14752764a2aaSMatthew G. Knepley }
14762764a2aaSMatthew G. Knepley 
14773e75805dSMatthew G. Knepley /*@C
1478cb36c0f9SMatthew G. Knepley   PetscDSGetRHSResidual - Get the pointwise RHS residual function for explicit timestepping for a given test field
1479cb36c0f9SMatthew G. Knepley 
148020f4b53cSBarry Smith   Not Collective
1481cb36c0f9SMatthew G. Knepley 
1482cb36c0f9SMatthew G. Knepley   Input Parameters:
1483dce8aebaSBarry Smith + ds - The `PetscDS`
1484cb36c0f9SMatthew G. Knepley - f  - The test field number
1485cb36c0f9SMatthew G. Knepley 
1486cb36c0f9SMatthew G. Knepley   Output Parameters:
1487cb36c0f9SMatthew G. Knepley + f0 - integrand for the test function term
1488cb36c0f9SMatthew G. Knepley - f1 - integrand for the test function gradient term
1489cb36c0f9SMatthew G. Knepley 
1490*a4e35b19SJacob Faibussowitsch   Calling sequence of `f0`:
1491cb36c0f9SMatthew G. Knepley + dim          - the spatial dimension
1492cb36c0f9SMatthew G. Knepley . Nf           - the number of fields
1493*a4e35b19SJacob Faibussowitsch . NfAux        - the number of auxiliary fields
1494cb36c0f9SMatthew G. Knepley . uOff         - the offset into u[] and u_t[] for each field
1495cb36c0f9SMatthew G. Knepley . uOff_x       - the offset into u_x[] for each field
1496cb36c0f9SMatthew G. Knepley . u            - each field evaluated at the current point
1497cb36c0f9SMatthew G. Knepley . u_t          - the time derivative of each field evaluated at the current point
1498cb36c0f9SMatthew G. Knepley . u_x          - the gradient of each field evaluated at the current point
1499cb36c0f9SMatthew G. Knepley . aOff         - the offset into a[] and a_t[] for each auxiliary field
1500cb36c0f9SMatthew G. Knepley . aOff_x       - the offset into a_x[] for each auxiliary field
1501cb36c0f9SMatthew G. Knepley . a            - each auxiliary field evaluated at the current point
1502cb36c0f9SMatthew G. Knepley . a_t          - the time derivative of each auxiliary field evaluated at the current point
1503cb36c0f9SMatthew G. Knepley . a_x          - the gradient of auxiliary each field evaluated at the current point
1504cb36c0f9SMatthew G. Knepley . t            - current time
1505cb36c0f9SMatthew G. Knepley . x            - coordinates of the current point
1506cb36c0f9SMatthew G. Knepley . numConstants - number of constant parameters
1507cb36c0f9SMatthew G. Knepley . constants    - constant parameters
1508cb36c0f9SMatthew G. Knepley - f0           - output values at the current point
1509cb36c0f9SMatthew G. Knepley 
1510cb36c0f9SMatthew G. Knepley   Level: intermediate
1511cb36c0f9SMatthew G. Knepley 
1512dce8aebaSBarry Smith   Note:
1513*a4e35b19SJacob Faibussowitsch   `f1` has an identical form and is omitted for brevity.
1514*a4e35b19SJacob Faibussowitsch 
1515dce8aebaSBarry Smith   We are using a first order FEM model for the weak form: \int_\Omega \phi f_0(u, u_t, \nabla u, x, t) + \nabla\phi \cdot {\vec f}_1(u, u_t, \nabla u, x, t)
1516dce8aebaSBarry Smith 
1517dce8aebaSBarry Smith .seealso: `PetscDS`, `PetscDSSetRHSResidual()`
1518cb36c0f9SMatthew G. Knepley @*/
1519*a4e35b19SJacob Faibussowitsch PetscErrorCode PetscDSGetRHSResidual(PetscDS ds, 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[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]), void (**f1)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
1520d71ae5a4SJacob Faibussowitsch {
1521cb36c0f9SMatthew G. Knepley   PetscPointFunc *tmp0, *tmp1;
1522cb36c0f9SMatthew G. Knepley   PetscInt        n0, n1;
1523cb36c0f9SMatthew G. Knepley 
1524cb36c0f9SMatthew G. Knepley   PetscFunctionBegin;
1525cb36c0f9SMatthew G. Knepley   PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1);
152663a3b9bcSJacob Faibussowitsch   PetscCheck(!(f < 0) && !(f >= ds->Nf), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %" PetscInt_FMT " must be in [0, %" PetscInt_FMT ")", f, ds->Nf);
15279566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormGetResidual(ds->wf, NULL, 0, f, 100, &n0, &tmp0, &n1, &tmp1));
1528cb36c0f9SMatthew G. Knepley   *f0 = tmp0 ? tmp0[0] : NULL;
1529cb36c0f9SMatthew G. Knepley   *f1 = tmp1 ? tmp1[0] : NULL;
15303ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1531cb36c0f9SMatthew G. Knepley }
1532cb36c0f9SMatthew G. Knepley 
1533cb36c0f9SMatthew G. Knepley /*@C
1534cb36c0f9SMatthew G. Knepley   PetscDSSetRHSResidual - Set the pointwise residual function for explicit timestepping for a given test field
1535cb36c0f9SMatthew G. Knepley 
153620f4b53cSBarry Smith   Not Collective
1537cb36c0f9SMatthew G. Knepley 
1538cb36c0f9SMatthew G. Knepley   Input Parameters:
1539dce8aebaSBarry Smith + ds - The `PetscDS`
1540cb36c0f9SMatthew G. Knepley . f  - The test field number
1541cb36c0f9SMatthew G. Knepley . f0 - integrand for the test function term
1542cb36c0f9SMatthew G. Knepley - f1 - integrand for the test function gradient term
1543cb36c0f9SMatthew G. Knepley 
1544*a4e35b19SJacob Faibussowitsch   Calling sequence for the callbacks `f0`:
1545cb36c0f9SMatthew G. Knepley + dim          - the spatial dimension
1546cb36c0f9SMatthew G. Knepley . Nf           - the number of fields
1547*a4e35b19SJacob Faibussowitsch . NfAux        - the number of auxiliary fields
1548cb36c0f9SMatthew G. Knepley . uOff         - the offset into u[] and u_t[] for each field
1549cb36c0f9SMatthew G. Knepley . uOff_x       - the offset into u_x[] for each field
1550cb36c0f9SMatthew G. Knepley . u            - each field evaluated at the current point
1551cb36c0f9SMatthew G. Knepley . u_t          - the time derivative of each field evaluated at the current point
1552cb36c0f9SMatthew G. Knepley . u_x          - the gradient of each field evaluated at the current point
1553cb36c0f9SMatthew G. Knepley . aOff         - the offset into a[] and a_t[] for each auxiliary field
1554cb36c0f9SMatthew G. Knepley . aOff_x       - the offset into a_x[] for each auxiliary field
1555cb36c0f9SMatthew G. Knepley . a            - each auxiliary field evaluated at the current point
1556cb36c0f9SMatthew G. Knepley . a_t          - the time derivative of each auxiliary field evaluated at the current point
1557cb36c0f9SMatthew G. Knepley . a_x          - the gradient of auxiliary each field evaluated at the current point
1558cb36c0f9SMatthew G. Knepley . t            - current time
1559cb36c0f9SMatthew G. Knepley . x            - coordinates of the current point
1560cb36c0f9SMatthew G. Knepley . numConstants - number of constant parameters
1561cb36c0f9SMatthew G. Knepley . constants    - constant parameters
1562cb36c0f9SMatthew G. Knepley - f0           - output values at the current point
1563cb36c0f9SMatthew G. Knepley 
1564cb36c0f9SMatthew G. Knepley   Level: intermediate
1565cb36c0f9SMatthew G. Knepley 
1566dce8aebaSBarry Smith   Note:
1567*a4e35b19SJacob Faibussowitsch   `f1` has an identical form and is omitted for brevity.
1568*a4e35b19SJacob Faibussowitsch 
1569dce8aebaSBarry Smith   We are using a first order FEM model for the weak form: \int_\Omega \phi f_0(u, u_t, \nabla u, x, t) + \nabla\phi \cdot {\vec f}_1(u, u_t, \nabla u, x, t)
1570dce8aebaSBarry Smith 
1571dce8aebaSBarry Smith .seealso: `PetscDS`, `PetscDSGetResidual()`
1572cb36c0f9SMatthew G. Knepley @*/
1573*a4e35b19SJacob Faibussowitsch PetscErrorCode PetscDSSetRHSResidual(PetscDS ds, 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[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]), void (*f1)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
1574d71ae5a4SJacob Faibussowitsch {
1575cb36c0f9SMatthew G. Knepley   PetscFunctionBegin;
1576cb36c0f9SMatthew G. Knepley   PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1);
1577cb36c0f9SMatthew G. Knepley   if (f0) PetscValidFunction(f0, 3);
1578cb36c0f9SMatthew G. Knepley   if (f1) PetscValidFunction(f1, 4);
157963a3b9bcSJacob Faibussowitsch   PetscCheck(f >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %" PetscInt_FMT " must be non-negative", f);
15809566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormSetIndexResidual(ds->wf, NULL, 0, f, 100, 0, f0, 0, f1));
15813ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1582cb36c0f9SMatthew G. Knepley }
1583cb36c0f9SMatthew G. Knepley 
1584cb36c0f9SMatthew G. Knepley /*@C
1585dce8aebaSBarry Smith   PetscDSHasJacobian - Checks that the Jacobian functions have been set
15863e75805dSMatthew G. Knepley 
158720f4b53cSBarry Smith   Not Collective
15883e75805dSMatthew G. Knepley 
15893e75805dSMatthew G. Knepley   Input Parameter:
159060225df5SJacob Faibussowitsch . ds - The `PetscDS`
15913e75805dSMatthew G. Knepley 
15923e75805dSMatthew G. Knepley   Output Parameter:
15933e75805dSMatthew G. Knepley . hasJac - flag that pointwise function for the Jacobian has been set
15943e75805dSMatthew G. Knepley 
15953e75805dSMatthew G. Knepley   Level: intermediate
15963e75805dSMatthew G. Knepley 
1597dce8aebaSBarry Smith .seealso: `PetscDS`, `PetscDSGetJacobianPreconditioner()`, `PetscDSSetJacobianPreconditioner()`, `PetscDSGetJacobian()`
15983e75805dSMatthew G. Knepley @*/
1599d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDSHasJacobian(PetscDS ds, PetscBool *hasJac)
1600d71ae5a4SJacob Faibussowitsch {
16013e75805dSMatthew G. Knepley   PetscFunctionBegin;
16026528b96dSMatthew G. Knepley   PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1);
16039566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormHasJacobian(ds->wf, hasJac));
16043ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
16053e75805dSMatthew G. Knepley }
16063e75805dSMatthew G. Knepley 
1607194d53e6SMatthew G. Knepley /*@C
1608194d53e6SMatthew G. Knepley   PetscDSGetJacobian - Get the pointwise Jacobian function for given test and basis field
1609194d53e6SMatthew G. Knepley 
161020f4b53cSBarry Smith   Not Collective
1611194d53e6SMatthew G. Knepley 
1612194d53e6SMatthew G. Knepley   Input Parameters:
1613dce8aebaSBarry Smith + ds - The `PetscDS`
1614194d53e6SMatthew G. Knepley . f  - The test field number
1615194d53e6SMatthew G. Knepley - g  - The field number
1616194d53e6SMatthew G. Knepley 
1617194d53e6SMatthew G. Knepley   Output Parameters:
1618194d53e6SMatthew G. Knepley + g0 - integrand for the test and basis function term
1619194d53e6SMatthew G. Knepley . g1 - integrand for the test function and basis function gradient term
1620194d53e6SMatthew G. Knepley . g2 - integrand for the test function gradient and basis function term
1621194d53e6SMatthew G. Knepley - g3 - integrand for the test function gradient and basis function gradient term
1622194d53e6SMatthew G. Knepley 
1623*a4e35b19SJacob Faibussowitsch   Calling sequence of `g0`:
1624194d53e6SMatthew G. Knepley + dim          - the spatial dimension
1625194d53e6SMatthew G. Knepley . Nf           - the number of fields
1626*a4e35b19SJacob Faibussowitsch . NfAux        - the number of auxiliary fields
1627194d53e6SMatthew G. Knepley . uOff         - the offset into u[] and u_t[] for each field
1628194d53e6SMatthew G. Knepley . uOff_x       - the offset into u_x[] for each field
1629194d53e6SMatthew G. Knepley . u            - each field evaluated at the current point
1630194d53e6SMatthew G. Knepley . u_t          - the time derivative of each field evaluated at the current point
1631194d53e6SMatthew G. Knepley . u_x          - the gradient of each field evaluated at the current point
1632194d53e6SMatthew G. Knepley . aOff         - the offset into a[] and a_t[] for each auxiliary field
1633194d53e6SMatthew G. Knepley . aOff_x       - the offset into a_x[] for each auxiliary field
1634194d53e6SMatthew G. Knepley . a            - each auxiliary field evaluated at the current point
1635194d53e6SMatthew G. Knepley . a_t          - the time derivative of each auxiliary field evaluated at the current point
1636194d53e6SMatthew G. Knepley . a_x          - the gradient of auxiliary each field evaluated at the current point
1637194d53e6SMatthew G. Knepley . t            - current time
16382aa1fc23SMatthew G. Knepley . u_tShift     - the multiplier a for dF/dU_t
1639194d53e6SMatthew G. Knepley . x            - coordinates of the current point
164097b6e6e8SMatthew G. Knepley . numConstants - number of constant parameters
164197b6e6e8SMatthew G. Knepley . constants    - constant parameters
1642194d53e6SMatthew G. Knepley - g0           - output values at the current point
1643194d53e6SMatthew G. Knepley 
1644194d53e6SMatthew G. Knepley   Level: intermediate
1645194d53e6SMatthew G. Knepley 
1646dce8aebaSBarry Smith   Note:
1647*a4e35b19SJacob Faibussowitsch   `g1`, `g2`, and `g3` have identical calling sequences to `g0` and are omitted for brevity.
164860225df5SJacob Faibussowitsch 
1649*a4e35b19SJacob Faibussowitsch   We are using a first order FEM model for the weak form\:
1650dce8aebaSBarry Smith   \int_\Omega \phi g_0(u, u_t, \nabla u, x, t) \psi + \phi {\vec g}_1(u, u_t, \nabla u, x, t) \nabla \psi + \nabla\phi \cdot {\vec g}_2(u, u_t, \nabla u, x, t) \psi + \nabla\phi \cdot {\overleftrightarrow g}_3(u, u_t, \nabla u, x, t) \cdot \nabla \psi
1651dce8aebaSBarry Smith 
1652dce8aebaSBarry Smith .seealso: `PetscDS`, `PetscDSSetJacobian()`
1653194d53e6SMatthew G. Knepley @*/
1654*a4e35b19SJacob Faibussowitsch PetscErrorCode PetscDSGetJacobian(PetscDS ds, PetscInt f, PetscInt g, void (**g0)(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, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]), void (**g1)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), void (**g2)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), void (**g3)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
1655d71ae5a4SJacob Faibussowitsch {
16566528b96dSMatthew G. Knepley   PetscPointJac *tmp0, *tmp1, *tmp2, *tmp3;
16576528b96dSMatthew G. Knepley   PetscInt       n0, n1, n2, n3;
16586528b96dSMatthew G. Knepley 
16592764a2aaSMatthew G. Knepley   PetscFunctionBegin;
16606528b96dSMatthew G. Knepley   PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1);
166163a3b9bcSJacob Faibussowitsch   PetscCheck(!(f < 0) && !(f >= ds->Nf), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %" PetscInt_FMT " must be in [0, %" PetscInt_FMT ")", f, ds->Nf);
166263a3b9bcSJacob Faibussowitsch   PetscCheck(!(g < 0) && !(g >= ds->Nf), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %" PetscInt_FMT " must be in [0, %" PetscInt_FMT ")", g, ds->Nf);
16639566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormGetJacobian(ds->wf, NULL, 0, f, g, 0, &n0, &tmp0, &n1, &tmp1, &n2, &tmp2, &n3, &tmp3));
16646528b96dSMatthew G. Knepley   *g0 = tmp0 ? tmp0[0] : NULL;
16656528b96dSMatthew G. Knepley   *g1 = tmp1 ? tmp1[0] : NULL;
16666528b96dSMatthew G. Knepley   *g2 = tmp2 ? tmp2[0] : NULL;
16676528b96dSMatthew G. Knepley   *g3 = tmp3 ? tmp3[0] : NULL;
16683ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
16692764a2aaSMatthew G. Knepley }
16702764a2aaSMatthew G. Knepley 
1671194d53e6SMatthew G. Knepley /*@C
1672194d53e6SMatthew G. Knepley   PetscDSSetJacobian - Set the pointwise Jacobian function for given test and basis fields
1673194d53e6SMatthew G. Knepley 
167420f4b53cSBarry Smith   Not Collective
1675194d53e6SMatthew G. Knepley 
1676194d53e6SMatthew G. Knepley   Input Parameters:
1677dce8aebaSBarry Smith + ds - The `PetscDS`
1678194d53e6SMatthew G. Knepley . f  - The test field number
1679194d53e6SMatthew G. Knepley . g  - The field number
1680194d53e6SMatthew G. Knepley . g0 - integrand for the test and basis function term
1681194d53e6SMatthew G. Knepley . g1 - integrand for the test function and basis function gradient term
1682194d53e6SMatthew G. Knepley . g2 - integrand for the test function gradient and basis function term
1683194d53e6SMatthew G. Knepley - g3 - integrand for the test function gradient and basis function gradient term
1684194d53e6SMatthew G. Knepley 
1685*a4e35b19SJacob Faibussowitsch   Calling sequence of `g0`:
1686194d53e6SMatthew G. Knepley + dim          - the spatial dimension
1687194d53e6SMatthew G. Knepley . Nf           - the number of fields
1688*a4e35b19SJacob Faibussowitsch . NfAux        - the number of auxiliary fields
1689194d53e6SMatthew G. Knepley . uOff         - the offset into u[] and u_t[] for each field
1690194d53e6SMatthew G. Knepley . uOff_x       - the offset into u_x[] for each field
1691194d53e6SMatthew G. Knepley . u            - each field evaluated at the current point
1692194d53e6SMatthew G. Knepley . u_t          - the time derivative of each field evaluated at the current point
1693194d53e6SMatthew G. Knepley . u_x          - the gradient of each field evaluated at the current point
1694194d53e6SMatthew G. Knepley . aOff         - the offset into a[] and a_t[] for each auxiliary field
1695194d53e6SMatthew G. Knepley . aOff_x       - the offset into a_x[] for each auxiliary field
1696194d53e6SMatthew G. Knepley . a            - each auxiliary field evaluated at the current point
1697194d53e6SMatthew G. Knepley . a_t          - the time derivative of each auxiliary field evaluated at the current point
1698194d53e6SMatthew G. Knepley . a_x          - the gradient of auxiliary each field evaluated at the current point
1699194d53e6SMatthew G. Knepley . t            - current time
17002aa1fc23SMatthew G. Knepley . u_tShift     - the multiplier a for dF/dU_t
1701194d53e6SMatthew G. Knepley . x            - coordinates of the current point
170297b6e6e8SMatthew G. Knepley . numConstants - number of constant parameters
170397b6e6e8SMatthew G. Knepley . constants    - constant parameters
1704194d53e6SMatthew G. Knepley - g0           - output values at the current point
1705194d53e6SMatthew G. Knepley 
1706194d53e6SMatthew G. Knepley   Level: intermediate
1707194d53e6SMatthew G. Knepley 
1708dce8aebaSBarry Smith   Note:
1709*a4e35b19SJacob Faibussowitsch   `g1`, `g2`, and `g3` have identical calling sequences to `g0` and are omitted for brevity.
171060225df5SJacob Faibussowitsch 
1711*a4e35b19SJacob Faibussowitsch   We are using a first order FEM model for the weak form\:
1712dce8aebaSBarry Smith   \int_\Omega \phi g_0(u, u_t, \nabla u, x, t) \psi + \phi {\vec g}_1(u, u_t, \nabla u, x, t) \nabla \psi + \nabla\phi \cdot {\vec g}_2(u, u_t, \nabla u, x, t) \psi + \nabla\phi \cdot {\overleftrightarrow g}_3(u, u_t, \nabla u, x, t) \cdot \nabla \psi
1713dce8aebaSBarry Smith 
1714dce8aebaSBarry Smith .seealso: `PetscDS`, `PetscDSGetJacobian()`
1715194d53e6SMatthew G. Knepley @*/
1716*a4e35b19SJacob Faibussowitsch PetscErrorCode PetscDSSetJacobian(PetscDS ds, PetscInt f, PetscInt g, void (*g0)(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, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]), void (*g1)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), void (*g2)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), void (*g3)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
1717d71ae5a4SJacob Faibussowitsch {
17182764a2aaSMatthew G. Knepley   PetscFunctionBegin;
17196528b96dSMatthew G. Knepley   PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1);
17202764a2aaSMatthew G. Knepley   if (g0) PetscValidFunction(g0, 4);
17212764a2aaSMatthew G. Knepley   if (g1) PetscValidFunction(g1, 5);
17222764a2aaSMatthew G. Knepley   if (g2) PetscValidFunction(g2, 6);
17232764a2aaSMatthew G. Knepley   if (g3) PetscValidFunction(g3, 7);
172463a3b9bcSJacob Faibussowitsch   PetscCheck(f >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %" PetscInt_FMT " must be non-negative", f);
172563a3b9bcSJacob Faibussowitsch   PetscCheck(g >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %" PetscInt_FMT " must be non-negative", g);
17269566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormSetIndexJacobian(ds->wf, NULL, 0, f, g, 0, 0, g0, 0, g1, 0, g2, 0, g3));
17273ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
17282764a2aaSMatthew G. Knepley }
17292764a2aaSMatthew G. Knepley 
1730475e0ac9SMatthew G. Knepley /*@C
1731dce8aebaSBarry Smith   PetscDSUseJacobianPreconditioner - Set whether to construct a Jacobian preconditioner
173255c1f793SMatthew G. Knepley 
173320f4b53cSBarry Smith   Not Collective
173455c1f793SMatthew G. Knepley 
173555c1f793SMatthew G. Knepley   Input Parameters:
1736dce8aebaSBarry Smith + prob      - The `PetscDS`
173755c1f793SMatthew G. Knepley - useJacPre - flag that enables construction of a Jacobian preconditioner
173855c1f793SMatthew G. Knepley 
173955c1f793SMatthew G. Knepley   Level: intermediate
174055c1f793SMatthew G. Knepley 
1741dce8aebaSBarry Smith .seealso: `PetscDS`, `PetscDSGetJacobianPreconditioner()`, `PetscDSSetJacobianPreconditioner()`, `PetscDSGetJacobian()`
174255c1f793SMatthew G. Knepley @*/
1743d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDSUseJacobianPreconditioner(PetscDS prob, PetscBool useJacPre)
1744d71ae5a4SJacob Faibussowitsch {
174555c1f793SMatthew G. Knepley   PetscFunctionBegin;
174655c1f793SMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
174755c1f793SMatthew G. Knepley   prob->useJacPre = useJacPre;
17483ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
174955c1f793SMatthew G. Knepley }
175055c1f793SMatthew G. Knepley 
175155c1f793SMatthew G. Knepley /*@C
1752dce8aebaSBarry Smith   PetscDSHasJacobianPreconditioner - Checks if a Jacobian preconditioner matrix has been set
1753475e0ac9SMatthew G. Knepley 
175420f4b53cSBarry Smith   Not Collective
1755475e0ac9SMatthew G. Knepley 
1756475e0ac9SMatthew G. Knepley   Input Parameter:
175760225df5SJacob Faibussowitsch . ds - The `PetscDS`
1758475e0ac9SMatthew G. Knepley 
1759475e0ac9SMatthew G. Knepley   Output Parameter:
1760475e0ac9SMatthew G. Knepley . hasJacPre - flag that pointwise function for Jacobian preconditioner matrix has been set
1761475e0ac9SMatthew G. Knepley 
1762475e0ac9SMatthew G. Knepley   Level: intermediate
1763475e0ac9SMatthew G. Knepley 
1764dce8aebaSBarry Smith .seealso: `PetscDS`, `PetscDSGetJacobianPreconditioner()`, `PetscDSSetJacobianPreconditioner()`, `PetscDSGetJacobian()`
1765475e0ac9SMatthew G. Knepley @*/
1766d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDSHasJacobianPreconditioner(PetscDS ds, PetscBool *hasJacPre)
1767d71ae5a4SJacob Faibussowitsch {
1768475e0ac9SMatthew G. Knepley   PetscFunctionBegin;
17696528b96dSMatthew G. Knepley   PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1);
1770475e0ac9SMatthew G. Knepley   *hasJacPre = PETSC_FALSE;
17713ba16761SJacob Faibussowitsch   if (!ds->useJacPre) PetscFunctionReturn(PETSC_SUCCESS);
17729566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormHasJacobianPreconditioner(ds->wf, hasJacPre));
17733ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1774475e0ac9SMatthew G. Knepley }
1775475e0ac9SMatthew G. Knepley 
1776475e0ac9SMatthew G. Knepley /*@C
1777dce8aebaSBarry Smith   PetscDSGetJacobianPreconditioner - Get the pointwise Jacobian preconditioner function for given test and basis field. If this is missing,
1778dce8aebaSBarry Smith   the system matrix is used to build the preconditioner.
1779475e0ac9SMatthew G. Knepley 
178020f4b53cSBarry Smith   Not Collective
1781475e0ac9SMatthew G. Knepley 
1782475e0ac9SMatthew G. Knepley   Input Parameters:
1783dce8aebaSBarry Smith + ds - The `PetscDS`
1784475e0ac9SMatthew G. Knepley . f  - The test field number
1785475e0ac9SMatthew G. Knepley - g  - The field number
1786475e0ac9SMatthew G. Knepley 
1787475e0ac9SMatthew G. Knepley   Output Parameters:
1788475e0ac9SMatthew G. Knepley + g0 - integrand for the test and basis function term
1789475e0ac9SMatthew G. Knepley . g1 - integrand for the test function and basis function gradient term
1790475e0ac9SMatthew G. Knepley . g2 - integrand for the test function gradient and basis function term
1791475e0ac9SMatthew G. Knepley - g3 - integrand for the test function gradient and basis function gradient term
1792475e0ac9SMatthew G. Knepley 
1793*a4e35b19SJacob Faibussowitsch   Calling sequence of `g0`:
1794475e0ac9SMatthew G. Knepley + dim          - the spatial dimension
1795475e0ac9SMatthew G. Knepley . Nf           - the number of fields
1796*a4e35b19SJacob Faibussowitsch . NfAux        - the number of auxiliary fields
1797475e0ac9SMatthew G. Knepley . uOff         - the offset into u[] and u_t[] for each field
1798475e0ac9SMatthew G. Knepley . uOff_x       - the offset into u_x[] for each field
1799475e0ac9SMatthew G. Knepley . u            - each field evaluated at the current point
1800475e0ac9SMatthew G. Knepley . u_t          - the time derivative of each field evaluated at the current point
1801475e0ac9SMatthew G. Knepley . u_x          - the gradient of each field evaluated at the current point
1802475e0ac9SMatthew G. Knepley . aOff         - the offset into a[] and a_t[] for each auxiliary field
1803475e0ac9SMatthew G. Knepley . aOff_x       - the offset into a_x[] for each auxiliary field
1804475e0ac9SMatthew G. Knepley . a            - each auxiliary field evaluated at the current point
1805475e0ac9SMatthew G. Knepley . a_t          - the time derivative of each auxiliary field evaluated at the current point
1806475e0ac9SMatthew G. Knepley . a_x          - the gradient of auxiliary each field evaluated at the current point
1807475e0ac9SMatthew G. Knepley . t            - current time
1808475e0ac9SMatthew G. Knepley . u_tShift     - the multiplier a for dF/dU_t
1809475e0ac9SMatthew G. Knepley . x            - coordinates of the current point
181097b6e6e8SMatthew G. Knepley . numConstants - number of constant parameters
181197b6e6e8SMatthew G. Knepley . constants    - constant parameters
1812475e0ac9SMatthew G. Knepley - g0           - output values at the current point
1813475e0ac9SMatthew G. Knepley 
1814475e0ac9SMatthew G. Knepley   Level: intermediate
1815475e0ac9SMatthew G. Knepley 
1816dce8aebaSBarry Smith   Note:
1817*a4e35b19SJacob Faibussowitsch   `g1`, `g2`, and `g3` have identical calling sequences to `g0` and are omitted for brevity.
1818*a4e35b19SJacob Faibussowitsch   We are using a first order FEM model for the weak form\:
1819dce8aebaSBarry Smith   \int_\Omega \phi g_0(u, u_t, \nabla u, x, t) \psi + \phi {\vec g}_1(u, u_t, \nabla u, x, t) \nabla \psi + \nabla\phi \cdot {\vec g}_2(u, u_t, \nabla u, x, t) \psi + \nabla\phi \cdot {\overleftrightarrow g}_3(u, u_t, \nabla u, x, t) \cdot \nabla \psi
1820dce8aebaSBarry Smith 
1821dce8aebaSBarry Smith .seealso: `PetscDS`, `PetscDSSetJacobianPreconditioner()`, `PetscDSGetJacobian()`
1822475e0ac9SMatthew G. Knepley @*/
1823*a4e35b19SJacob Faibussowitsch PetscErrorCode PetscDSGetJacobianPreconditioner(PetscDS ds, PetscInt f, PetscInt g, void (**g0)(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, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]), void (**g1)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), void (**g2)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), void (**g3)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
1824d71ae5a4SJacob Faibussowitsch {
18256528b96dSMatthew G. Knepley   PetscPointJac *tmp0, *tmp1, *tmp2, *tmp3;
18266528b96dSMatthew G. Knepley   PetscInt       n0, n1, n2, n3;
18276528b96dSMatthew G. Knepley 
1828475e0ac9SMatthew G. Knepley   PetscFunctionBegin;
18296528b96dSMatthew G. Knepley   PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1);
183063a3b9bcSJacob Faibussowitsch   PetscCheck(!(f < 0) && !(f >= ds->Nf), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %" PetscInt_FMT " must be in [0, %" PetscInt_FMT ")", f, ds->Nf);
183163a3b9bcSJacob Faibussowitsch   PetscCheck(!(g < 0) && !(g >= ds->Nf), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %" PetscInt_FMT " must be in [0, %" PetscInt_FMT ")", g, ds->Nf);
18329566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormGetJacobianPreconditioner(ds->wf, NULL, 0, f, g, 0, &n0, &tmp0, &n1, &tmp1, &n2, &tmp2, &n3, &tmp3));
18336528b96dSMatthew G. Knepley   *g0 = tmp0 ? tmp0[0] : NULL;
18346528b96dSMatthew G. Knepley   *g1 = tmp1 ? tmp1[0] : NULL;
18356528b96dSMatthew G. Knepley   *g2 = tmp2 ? tmp2[0] : NULL;
18366528b96dSMatthew G. Knepley   *g3 = tmp3 ? tmp3[0] : NULL;
18373ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1838475e0ac9SMatthew G. Knepley }
1839475e0ac9SMatthew G. Knepley 
1840475e0ac9SMatthew G. Knepley /*@C
1841dce8aebaSBarry Smith   PetscDSSetJacobianPreconditioner - Set the pointwise Jacobian preconditioner function for given test and basis fields.
1842dce8aebaSBarry Smith   If this is missing, the system matrix is used to build the preconditioner.
1843475e0ac9SMatthew G. Knepley 
184420f4b53cSBarry Smith   Not Collective
1845475e0ac9SMatthew G. Knepley 
1846475e0ac9SMatthew G. Knepley   Input Parameters:
1847dce8aebaSBarry Smith + ds - The `PetscDS`
1848475e0ac9SMatthew G. Knepley . f  - The test field number
1849475e0ac9SMatthew G. Knepley . g  - The field number
1850475e0ac9SMatthew G. Knepley . g0 - integrand for the test and basis function term
1851475e0ac9SMatthew G. Knepley . g1 - integrand for the test function and basis function gradient term
1852475e0ac9SMatthew G. Knepley . g2 - integrand for the test function gradient and basis function term
1853475e0ac9SMatthew G. Knepley - g3 - integrand for the test function gradient and basis function gradient term
1854475e0ac9SMatthew G. Knepley 
1855*a4e35b19SJacob Faibussowitsch   Calling sequence of `g0`:
1856475e0ac9SMatthew G. Knepley + dim          - the spatial dimension
1857475e0ac9SMatthew G. Knepley . Nf           - the number of fields
1858*a4e35b19SJacob Faibussowitsch . NfAux        - the number of auxiliary fields
1859475e0ac9SMatthew G. Knepley . uOff         - the offset into u[] and u_t[] for each field
1860475e0ac9SMatthew G. Knepley . uOff_x       - the offset into u_x[] for each field
1861475e0ac9SMatthew G. Knepley . u            - each field evaluated at the current point
1862475e0ac9SMatthew G. Knepley . u_t          - the time derivative of each field evaluated at the current point
1863475e0ac9SMatthew G. Knepley . u_x          - the gradient of each field evaluated at the current point
1864475e0ac9SMatthew G. Knepley . aOff         - the offset into a[] and a_t[] for each auxiliary field
1865475e0ac9SMatthew G. Knepley . aOff_x       - the offset into a_x[] for each auxiliary field
1866475e0ac9SMatthew G. Knepley . a            - each auxiliary field evaluated at the current point
1867475e0ac9SMatthew G. Knepley . a_t          - the time derivative of each auxiliary field evaluated at the current point
1868475e0ac9SMatthew G. Knepley . a_x          - the gradient of auxiliary each field evaluated at the current point
1869475e0ac9SMatthew G. Knepley . t            - current time
1870475e0ac9SMatthew G. Knepley . u_tShift     - the multiplier a for dF/dU_t
1871475e0ac9SMatthew G. Knepley . x            - coordinates of the current point
187297b6e6e8SMatthew G. Knepley . numConstants - number of constant parameters
187397b6e6e8SMatthew G. Knepley . constants    - constant parameters
1874475e0ac9SMatthew G. Knepley - g0           - output values at the current point
1875475e0ac9SMatthew G. Knepley 
1876475e0ac9SMatthew G. Knepley   Level: intermediate
1877475e0ac9SMatthew G. Knepley 
1878dce8aebaSBarry Smith   Note:
1879*a4e35b19SJacob Faibussowitsch   `g1`, `g2`, and `g3` have identical calling sequences to `g0` and are omitted for brevity.
188060225df5SJacob Faibussowitsch 
1881*a4e35b19SJacob Faibussowitsch   We are using a first order FEM model for the weak form\:
1882dce8aebaSBarry Smith   \int_\Omega \phi g_0(u, u_t, \nabla u, x, t) \psi + \phi {\vec g}_1(u, u_t, \nabla u, x, t) \nabla \psi + \nabla\phi \cdot {\vec g}_2(u, u_t, \nabla u, x, t) \psi + \nabla\phi \cdot {\overleftrightarrow g}_3(u, u_t, \nabla u, x, t) \cdot \nabla \psi
1883dce8aebaSBarry Smith 
1884dce8aebaSBarry Smith .seealso: `PetscDS`, `PetscDSGetJacobianPreconditioner()`, `PetscDSSetJacobian()`
1885475e0ac9SMatthew G. Knepley @*/
1886*a4e35b19SJacob Faibussowitsch PetscErrorCode PetscDSSetJacobianPreconditioner(PetscDS ds, PetscInt f, PetscInt g, void (*g0)(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, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]), void (*g1)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), void (*g2)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), void (*g3)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
1887d71ae5a4SJacob Faibussowitsch {
1888475e0ac9SMatthew G. Knepley   PetscFunctionBegin;
18896528b96dSMatthew G. Knepley   PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1);
1890475e0ac9SMatthew G. Knepley   if (g0) PetscValidFunction(g0, 4);
1891475e0ac9SMatthew G. Knepley   if (g1) PetscValidFunction(g1, 5);
1892475e0ac9SMatthew G. Knepley   if (g2) PetscValidFunction(g2, 6);
1893475e0ac9SMatthew G. Knepley   if (g3) PetscValidFunction(g3, 7);
189463a3b9bcSJacob Faibussowitsch   PetscCheck(f >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %" PetscInt_FMT " must be non-negative", f);
189563a3b9bcSJacob Faibussowitsch   PetscCheck(g >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %" PetscInt_FMT " must be non-negative", g);
18969566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormSetIndexJacobianPreconditioner(ds->wf, NULL, 0, f, g, 0, 0, g0, 0, g1, 0, g2, 0, g3));
18973ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1898475e0ac9SMatthew G. Knepley }
1899475e0ac9SMatthew G. Knepley 
1900b7e05686SMatthew G. Knepley /*@C
1901b7e05686SMatthew G. Knepley   PetscDSHasDynamicJacobian - Signals that a dynamic Jacobian, dF/du_t, has been set
1902b7e05686SMatthew G. Knepley 
190320f4b53cSBarry Smith   Not Collective
1904b7e05686SMatthew G. Knepley 
1905b7e05686SMatthew G. Knepley   Input Parameter:
1906dce8aebaSBarry Smith . ds - The `PetscDS`
1907b7e05686SMatthew G. Knepley 
1908b7e05686SMatthew G. Knepley   Output Parameter:
1909b7e05686SMatthew G. Knepley . hasDynJac - flag that pointwise function for dynamic Jacobian has been set
1910b7e05686SMatthew G. Knepley 
1911b7e05686SMatthew G. Knepley   Level: intermediate
1912b7e05686SMatthew G. Knepley 
1913dce8aebaSBarry Smith .seealso: `PetscDS`, `PetscDSGetDynamicJacobian()`, `PetscDSSetDynamicJacobian()`, `PetscDSGetJacobian()`
1914b7e05686SMatthew G. Knepley @*/
1915d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDSHasDynamicJacobian(PetscDS ds, PetscBool *hasDynJac)
1916d71ae5a4SJacob Faibussowitsch {
1917b7e05686SMatthew G. Knepley   PetscFunctionBegin;
19186528b96dSMatthew G. Knepley   PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1);
19199566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormHasDynamicJacobian(ds->wf, hasDynJac));
19203ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1921b7e05686SMatthew G. Knepley }
1922b7e05686SMatthew G. Knepley 
1923b7e05686SMatthew G. Knepley /*@C
1924b7e05686SMatthew G. Knepley   PetscDSGetDynamicJacobian - Get the pointwise dynamic Jacobian, dF/du_t, function for given test and basis field
1925b7e05686SMatthew G. Knepley 
192620f4b53cSBarry Smith   Not Collective
1927b7e05686SMatthew G. Knepley 
1928b7e05686SMatthew G. Knepley   Input Parameters:
1929dce8aebaSBarry Smith + ds - The `PetscDS`
1930b7e05686SMatthew G. Knepley . f  - The test field number
1931b7e05686SMatthew G. Knepley - g  - The field number
1932b7e05686SMatthew G. Knepley 
1933b7e05686SMatthew G. Knepley   Output Parameters:
1934b7e05686SMatthew G. Knepley + g0 - integrand for the test and basis function term
1935b7e05686SMatthew G. Knepley . g1 - integrand for the test function and basis function gradient term
1936b7e05686SMatthew G. Knepley . g2 - integrand for the test function gradient and basis function term
1937b7e05686SMatthew G. Knepley - g3 - integrand for the test function gradient and basis function gradient term
1938b7e05686SMatthew G. Knepley 
1939*a4e35b19SJacob Faibussowitsch   Calling sequence of `g0`:
1940b7e05686SMatthew G. Knepley + dim          - the spatial dimension
1941b7e05686SMatthew G. Knepley . Nf           - the number of fields
1942*a4e35b19SJacob Faibussowitsch . NfAux        - the number of auxiliary fields
1943b7e05686SMatthew G. Knepley . uOff         - the offset into u[] and u_t[] for each field
1944b7e05686SMatthew G. Knepley . uOff_x       - the offset into u_x[] for each field
1945b7e05686SMatthew G. Knepley . u            - each field evaluated at the current point
1946b7e05686SMatthew G. Knepley . u_t          - the time derivative of each field evaluated at the current point
1947b7e05686SMatthew G. Knepley . u_x          - the gradient of each field evaluated at the current point
1948b7e05686SMatthew G. Knepley . aOff         - the offset into a[] and a_t[] for each auxiliary field
1949b7e05686SMatthew G. Knepley . aOff_x       - the offset into a_x[] for each auxiliary field
1950b7e05686SMatthew G. Knepley . a            - each auxiliary field evaluated at the current point
1951b7e05686SMatthew G. Knepley . a_t          - the time derivative of each auxiliary field evaluated at the current point
1952b7e05686SMatthew G. Knepley . a_x          - the gradient of auxiliary each field evaluated at the current point
1953b7e05686SMatthew G. Knepley . t            - current time
1954b7e05686SMatthew G. Knepley . u_tShift     - the multiplier a for dF/dU_t
1955b7e05686SMatthew G. Knepley . x            - coordinates of the current point
195697b6e6e8SMatthew G. Knepley . numConstants - number of constant parameters
195797b6e6e8SMatthew G. Knepley . constants    - constant parameters
1958b7e05686SMatthew G. Knepley - g0           - output values at the current point
1959b7e05686SMatthew G. Knepley 
1960b7e05686SMatthew G. Knepley   Level: intermediate
1961b7e05686SMatthew G. Knepley 
1962dce8aebaSBarry Smith   Note:
1963*a4e35b19SJacob Faibussowitsch   `g1`, `g2`, and `g3` have identical calling sequences to `g0` and are omitted for brevity.
196460225df5SJacob Faibussowitsch 
1965*a4e35b19SJacob Faibussowitsch   We are using a first order FEM model for the weak form\:
1966dce8aebaSBarry Smith   \int_\Omega \phi g_0(u, u_t, \nabla u, x, t) \psi + \phi {\vec g}_1(u, u_t, \nabla u, x, t) \nabla \psi + \nabla\phi \cdot {\vec g}_2(u, u_t, \nabla u, x, t) \psi + \nabla\phi \cdot {\overleftrightarrow g}_3(u, u_t, \nabla u, x, t) \cdot \nabla \psi
1967dce8aebaSBarry Smith 
1968dce8aebaSBarry Smith .seealso: `PetscDS`, `PetscDSSetJacobian()`
1969b7e05686SMatthew G. Knepley @*/
1970*a4e35b19SJacob Faibussowitsch PetscErrorCode PetscDSGetDynamicJacobian(PetscDS ds, PetscInt f, PetscInt g, void (**g0)(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, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]), void (**g1)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), void (**g2)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), void (**g3)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
1971d71ae5a4SJacob Faibussowitsch {
19726528b96dSMatthew G. Knepley   PetscPointJac *tmp0, *tmp1, *tmp2, *tmp3;
19736528b96dSMatthew G. Knepley   PetscInt       n0, n1, n2, n3;
19746528b96dSMatthew G. Knepley 
1975b7e05686SMatthew G. Knepley   PetscFunctionBegin;
19766528b96dSMatthew G. Knepley   PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1);
197763a3b9bcSJacob Faibussowitsch   PetscCheck(!(f < 0) && !(f >= ds->Nf), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %" PetscInt_FMT " must be in [0, %" PetscInt_FMT ")", f, ds->Nf);
197863a3b9bcSJacob Faibussowitsch   PetscCheck(!(g < 0) && !(g >= ds->Nf), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %" PetscInt_FMT " must be in [0, %" PetscInt_FMT ")", g, ds->Nf);
19799566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormGetDynamicJacobian(ds->wf, NULL, 0, f, g, 0, &n0, &tmp0, &n1, &tmp1, &n2, &tmp2, &n3, &tmp3));
19806528b96dSMatthew G. Knepley   *g0 = tmp0 ? tmp0[0] : NULL;
19816528b96dSMatthew G. Knepley   *g1 = tmp1 ? tmp1[0] : NULL;
19826528b96dSMatthew G. Knepley   *g2 = tmp2 ? tmp2[0] : NULL;
19836528b96dSMatthew G. Knepley   *g3 = tmp3 ? tmp3[0] : NULL;
19843ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1985b7e05686SMatthew G. Knepley }
1986b7e05686SMatthew G. Knepley 
1987b7e05686SMatthew G. Knepley /*@C
1988b7e05686SMatthew G. Knepley   PetscDSSetDynamicJacobian - Set the pointwise dynamic Jacobian, dF/du_t, function for given test and basis fields
1989b7e05686SMatthew G. Knepley 
199020f4b53cSBarry Smith   Not Collective
1991b7e05686SMatthew G. Knepley 
1992b7e05686SMatthew G. Knepley   Input Parameters:
1993dce8aebaSBarry Smith + ds - The `PetscDS`
1994b7e05686SMatthew G. Knepley . f  - The test field number
1995b7e05686SMatthew G. Knepley . g  - The field number
1996b7e05686SMatthew G. Knepley . g0 - integrand for the test and basis function term
1997b7e05686SMatthew G. Knepley . g1 - integrand for the test function and basis function gradient term
1998b7e05686SMatthew G. Knepley . g2 - integrand for the test function gradient and basis function term
1999b7e05686SMatthew G. Knepley - g3 - integrand for the test function gradient and basis function gradient term
2000b7e05686SMatthew G. Knepley 
2001*a4e35b19SJacob Faibussowitsch   Calling sequence of `g0`:
2002b7e05686SMatthew G. Knepley + dim          - the spatial dimension
2003b7e05686SMatthew G. Knepley . Nf           - the number of fields
2004*a4e35b19SJacob Faibussowitsch . NfAux        - the number of auxiliary fields
2005b7e05686SMatthew G. Knepley . uOff         - the offset into u[] and u_t[] for each field
2006b7e05686SMatthew G. Knepley . uOff_x       - the offset into u_x[] for each field
2007b7e05686SMatthew G. Knepley . u            - each field evaluated at the current point
2008b7e05686SMatthew G. Knepley . u_t          - the time derivative of each field evaluated at the current point
2009b7e05686SMatthew G. Knepley . u_x          - the gradient of each field evaluated at the current point
2010b7e05686SMatthew G. Knepley . aOff         - the offset into a[] and a_t[] for each auxiliary field
2011b7e05686SMatthew G. Knepley . aOff_x       - the offset into a_x[] for each auxiliary field
2012b7e05686SMatthew G. Knepley . a            - each auxiliary field evaluated at the current point
2013b7e05686SMatthew G. Knepley . a_t          - the time derivative of each auxiliary field evaluated at the current point
2014b7e05686SMatthew G. Knepley . a_x          - the gradient of auxiliary each field evaluated at the current point
2015b7e05686SMatthew G. Knepley . t            - current time
2016b7e05686SMatthew G. Knepley . u_tShift     - the multiplier a for dF/dU_t
2017b7e05686SMatthew G. Knepley . x            - coordinates of the current point
201897b6e6e8SMatthew G. Knepley . numConstants - number of constant parameters
201997b6e6e8SMatthew G. Knepley . constants    - constant parameters
2020b7e05686SMatthew G. Knepley - g0           - output values at the current point
2021b7e05686SMatthew G. Knepley 
2022b7e05686SMatthew G. Knepley   Level: intermediate
2023b7e05686SMatthew G. Knepley 
2024dce8aebaSBarry Smith   Note:
2025*a4e35b19SJacob Faibussowitsch   `g1`, `g2`, and `g3` have identical calling sequences to `g0` and are omitted for brevity.
202660225df5SJacob Faibussowitsch 
2027*a4e35b19SJacob Faibussowitsch   We are using a first order FEM model for the weak form\:
2028dce8aebaSBarry Smith   \int_\Omega \phi g_0(u, u_t, \nabla u, x, t) \psi + \phi {\vec g}_1(u, u_t, \nabla u, x, t) \nabla \psi + \nabla\phi \cdot {\vec g}_2(u, u_t, \nabla u, x, t) \psi + \nabla\phi \cdot {\overleftrightarrow g}_3(u, u_t, \nabla u, x, t) \cdot \nabla \psi
2029dce8aebaSBarry Smith 
2030dce8aebaSBarry Smith .seealso: `PetscDS`, `PetscDSGetJacobian()`
2031b7e05686SMatthew G. Knepley @*/
2032*a4e35b19SJacob Faibussowitsch PetscErrorCode PetscDSSetDynamicJacobian(PetscDS ds, PetscInt f, PetscInt g, void (*g0)(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, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]), void (*g1)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), void (*g2)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), void (*g3)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
2033d71ae5a4SJacob Faibussowitsch {
2034b7e05686SMatthew G. Knepley   PetscFunctionBegin;
20356528b96dSMatthew G. Knepley   PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1);
2036b7e05686SMatthew G. Knepley   if (g0) PetscValidFunction(g0, 4);
2037b7e05686SMatthew G. Knepley   if (g1) PetscValidFunction(g1, 5);
2038b7e05686SMatthew G. Knepley   if (g2) PetscValidFunction(g2, 6);
2039b7e05686SMatthew G. Knepley   if (g3) PetscValidFunction(g3, 7);
204063a3b9bcSJacob Faibussowitsch   PetscCheck(f >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %" PetscInt_FMT " must be non-negative", f);
204163a3b9bcSJacob Faibussowitsch   PetscCheck(g >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %" PetscInt_FMT " must be non-negative", g);
20429566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormSetIndexDynamicJacobian(ds->wf, NULL, 0, f, g, 0, 0, g0, 0, g1, 0, g2, 0, g3));
20433ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2044b7e05686SMatthew G. Knepley }
2045b7e05686SMatthew G. Knepley 
20460c2f2876SMatthew G. Knepley /*@C
20470c2f2876SMatthew G. Knepley   PetscDSGetRiemannSolver - Returns the Riemann solver for the given field
20480c2f2876SMatthew G. Knepley 
204920f4b53cSBarry Smith   Not Collective
20500c2f2876SMatthew G. Knepley 
20514165533cSJose E. Roman   Input Parameters:
2052dce8aebaSBarry Smith + ds - The `PetscDS` object
20530c2f2876SMatthew G. Knepley - f  - The field number
20540c2f2876SMatthew G. Knepley 
20554165533cSJose E. Roman   Output Parameter:
20560c2f2876SMatthew G. Knepley . r - Riemann solver
20570c2f2876SMatthew G. Knepley 
205820f4b53cSBarry Smith   Calling sequence of `r`:
20595db36cf9SMatthew G. Knepley + dim          - The spatial dimension
20605db36cf9SMatthew G. Knepley . Nf           - The number of fields
20615db36cf9SMatthew G. Knepley . x            - The coordinates at a point on the interface
20620c2f2876SMatthew G. Knepley . n            - The normal vector to the interface
20630c2f2876SMatthew G. Knepley . uL           - The state vector to the left of the interface
20640c2f2876SMatthew G. Knepley . uR           - The state vector to the right of the interface
20650c2f2876SMatthew G. Knepley . flux         - output array of flux through the interface
206697b6e6e8SMatthew G. Knepley . numConstants - number of constant parameters
206797b6e6e8SMatthew G. Knepley . constants    - constant parameters
20680c2f2876SMatthew G. Knepley - ctx          - optional user context
20690c2f2876SMatthew G. Knepley 
20700c2f2876SMatthew G. Knepley   Level: intermediate
20710c2f2876SMatthew G. Knepley 
2072dce8aebaSBarry Smith .seealso: `PetscDS`, `PetscDSSetRiemannSolver()`
20730c2f2876SMatthew G. Knepley @*/
2074d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDSGetRiemannSolver(PetscDS ds, PetscInt f, void (**r)(PetscInt dim, PetscInt Nf, const PetscReal x[], const PetscReal n[], const PetscScalar uL[], const PetscScalar uR[], PetscInt numConstants, const PetscScalar constants[], PetscScalar flux[], void *ctx))
2075d71ae5a4SJacob Faibussowitsch {
20766528b96dSMatthew G. Knepley   PetscRiemannFunc *tmp;
20776528b96dSMatthew G. Knepley   PetscInt          n;
20786528b96dSMatthew G. Knepley 
20790c2f2876SMatthew G. Knepley   PetscFunctionBegin;
20806528b96dSMatthew G. Knepley   PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1);
20814f572ea9SToby Isaac   PetscAssertPointer(r, 3);
208263a3b9bcSJacob Faibussowitsch   PetscCheck(!(f < 0) && !(f >= ds->Nf), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %" PetscInt_FMT " must be in [0, %" PetscInt_FMT ")", f, ds->Nf);
20839566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormGetRiemannSolver(ds->wf, NULL, 0, f, 0, &n, &tmp));
20846528b96dSMatthew G. Knepley   *r = tmp ? tmp[0] : NULL;
20853ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
20860c2f2876SMatthew G. Knepley }
20870c2f2876SMatthew G. Knepley 
20880c2f2876SMatthew G. Knepley /*@C
20890c2f2876SMatthew G. Knepley   PetscDSSetRiemannSolver - Sets the Riemann solver for the given field
20900c2f2876SMatthew G. Knepley 
209120f4b53cSBarry Smith   Not Collective
20920c2f2876SMatthew G. Knepley 
20934165533cSJose E. Roman   Input Parameters:
2094dce8aebaSBarry Smith + ds - The `PetscDS` object
20950c2f2876SMatthew G. Knepley . f  - The field number
20960c2f2876SMatthew G. Knepley - r  - Riemann solver
20970c2f2876SMatthew G. Knepley 
209820f4b53cSBarry Smith   Calling sequence of `r`:
20995db36cf9SMatthew G. Knepley + dim          - The spatial dimension
21005db36cf9SMatthew G. Knepley . Nf           - The number of fields
21015db36cf9SMatthew G. Knepley . x            - The coordinates at a point on the interface
21020c2f2876SMatthew G. Knepley . n            - The normal vector to the interface
21030c2f2876SMatthew G. Knepley . uL           - The state vector to the left of the interface
21040c2f2876SMatthew G. Knepley . uR           - The state vector to the right of the interface
21050c2f2876SMatthew G. Knepley . flux         - output array of flux through the interface
210697b6e6e8SMatthew G. Knepley . numConstants - number of constant parameters
210797b6e6e8SMatthew G. Knepley . constants    - constant parameters
21080c2f2876SMatthew G. Knepley - ctx          - optional user context
21090c2f2876SMatthew G. Knepley 
21100c2f2876SMatthew G. Knepley   Level: intermediate
21110c2f2876SMatthew G. Knepley 
2112dce8aebaSBarry Smith .seealso: `PetscDS`, `PetscDSGetRiemannSolver()`
21130c2f2876SMatthew G. Knepley @*/
2114d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDSSetRiemannSolver(PetscDS ds, PetscInt f, void (*r)(PetscInt dim, PetscInt Nf, const PetscReal x[], const PetscReal n[], const PetscScalar uL[], const PetscScalar uR[], PetscInt numConstants, const PetscScalar constants[], PetscScalar flux[], void *ctx))
2115d71ae5a4SJacob Faibussowitsch {
21160c2f2876SMatthew G. Knepley   PetscFunctionBegin;
21176528b96dSMatthew G. Knepley   PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1);
2118de716cbcSToby Isaac   if (r) PetscValidFunction(r, 3);
211963a3b9bcSJacob Faibussowitsch   PetscCheck(f >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %" PetscInt_FMT " must be non-negative", f);
21209566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormSetIndexRiemannSolver(ds->wf, NULL, 0, f, 0, 0, r));
21213ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
21220c2f2876SMatthew G. Knepley }
21230c2f2876SMatthew G. Knepley 
212432d2bbc9SMatthew G. Knepley /*@C
212532d2bbc9SMatthew G. Knepley   PetscDSGetUpdate - Get the pointwise update function for a given field
212632d2bbc9SMatthew G. Knepley 
212720f4b53cSBarry Smith   Not Collective
212832d2bbc9SMatthew G. Knepley 
212932d2bbc9SMatthew G. Knepley   Input Parameters:
2130dce8aebaSBarry Smith + ds - The `PetscDS`
213132d2bbc9SMatthew G. Knepley - f  - The field number
213232d2bbc9SMatthew G. Knepley 
2133f899ff85SJose E. Roman   Output Parameter:
2134a2b725a8SWilliam Gropp . update - update function
213532d2bbc9SMatthew G. Knepley 
213620f4b53cSBarry Smith   Calling sequence of `update`:
213732d2bbc9SMatthew G. Knepley + dim          - the spatial dimension
213832d2bbc9SMatthew G. Knepley . Nf           - the number of fields
2139*a4e35b19SJacob Faibussowitsch . NfAux        - the number of auxiliary fields
214032d2bbc9SMatthew G. Knepley . uOff         - the offset into u[] and u_t[] for each field
214132d2bbc9SMatthew G. Knepley . uOff_x       - the offset into u_x[] for each field
214232d2bbc9SMatthew G. Knepley . u            - each field evaluated at the current point
214332d2bbc9SMatthew G. Knepley . u_t          - the time derivative of each field evaluated at the current point
214432d2bbc9SMatthew G. Knepley . u_x          - the gradient of each field evaluated at the current point
214532d2bbc9SMatthew G. Knepley . aOff         - the offset into a[] and a_t[] for each auxiliary field
214632d2bbc9SMatthew G. Knepley . aOff_x       - the offset into a_x[] for each auxiliary field
214732d2bbc9SMatthew G. Knepley . a            - each auxiliary field evaluated at the current point
214832d2bbc9SMatthew G. Knepley . a_t          - the time derivative of each auxiliary field evaluated at the current point
214932d2bbc9SMatthew G. Knepley . a_x          - the gradient of auxiliary each field evaluated at the current point
215032d2bbc9SMatthew G. Knepley . t            - current time
215132d2bbc9SMatthew G. Knepley . x            - coordinates of the current point
2152*a4e35b19SJacob Faibussowitsch . numConstants - number of constant parameters
2153*a4e35b19SJacob Faibussowitsch . constants    - constant parameters
215432d2bbc9SMatthew G. Knepley - uNew         - new value for field at the current point
215532d2bbc9SMatthew G. Knepley 
215632d2bbc9SMatthew G. Knepley   Level: intermediate
215732d2bbc9SMatthew G. Knepley 
2158dce8aebaSBarry Smith .seealso: `PetscDS`, `PetscDSSetUpdate()`, `PetscDSSetResidual()`
215932d2bbc9SMatthew G. Knepley @*/
2160d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDSGetUpdate(PetscDS ds, PetscInt f, void (**update)(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 uNew[]))
2161d71ae5a4SJacob Faibussowitsch {
216232d2bbc9SMatthew G. Knepley   PetscFunctionBegin;
21636528b96dSMatthew G. Knepley   PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1);
216463a3b9bcSJacob Faibussowitsch   PetscCheck(!(f < 0) && !(f >= ds->Nf), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %" PetscInt_FMT " must be in [0, %" PetscInt_FMT ")", f, ds->Nf);
21659371c9d4SSatish Balay   if (update) {
21664f572ea9SToby Isaac     PetscAssertPointer(update, 3);
21679371c9d4SSatish Balay     *update = ds->update[f];
21689371c9d4SSatish Balay   }
21693ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
217032d2bbc9SMatthew G. Knepley }
217132d2bbc9SMatthew G. Knepley 
217232d2bbc9SMatthew G. Knepley /*@C
21733fa77dffSMatthew G. Knepley   PetscDSSetUpdate - Set the pointwise update function for a given field
217432d2bbc9SMatthew G. Knepley 
217520f4b53cSBarry Smith   Not Collective
217632d2bbc9SMatthew G. Knepley 
217732d2bbc9SMatthew G. Knepley   Input Parameters:
2178dce8aebaSBarry Smith + ds     - The `PetscDS`
217932d2bbc9SMatthew G. Knepley . f      - The field number
218032d2bbc9SMatthew G. Knepley - update - update function
218132d2bbc9SMatthew G. Knepley 
218220f4b53cSBarry Smith   Calling sequence of `update`:
218332d2bbc9SMatthew G. Knepley + dim          - the spatial dimension
218432d2bbc9SMatthew G. Knepley . Nf           - the number of fields
2185*a4e35b19SJacob Faibussowitsch . NfAux        - the number of auxiliary fields
218632d2bbc9SMatthew G. Knepley . uOff         - the offset into u[] and u_t[] for each field
218732d2bbc9SMatthew G. Knepley . uOff_x       - the offset into u_x[] for each field
218832d2bbc9SMatthew G. Knepley . u            - each field evaluated at the current point
218932d2bbc9SMatthew G. Knepley . u_t          - the time derivative of each field evaluated at the current point
219032d2bbc9SMatthew G. Knepley . u_x          - the gradient of each field evaluated at the current point
219132d2bbc9SMatthew G. Knepley . aOff         - the offset into a[] and a_t[] for each auxiliary field
219232d2bbc9SMatthew G. Knepley . aOff_x       - the offset into a_x[] for each auxiliary field
219332d2bbc9SMatthew G. Knepley . a            - each auxiliary field evaluated at the current point
219432d2bbc9SMatthew G. Knepley . a_t          - the time derivative of each auxiliary field evaluated at the current point
219532d2bbc9SMatthew G. Knepley . a_x          - the gradient of auxiliary each field evaluated at the current point
219632d2bbc9SMatthew G. Knepley . t            - current time
219732d2bbc9SMatthew G. Knepley . x            - coordinates of the current point
2198*a4e35b19SJacob Faibussowitsch . numConstants - number of constant parameters
2199*a4e35b19SJacob Faibussowitsch . constants    - constant parameters
220032d2bbc9SMatthew G. Knepley - uNew         - new field values at the current point
220132d2bbc9SMatthew G. Knepley 
220232d2bbc9SMatthew G. Knepley   Level: intermediate
220332d2bbc9SMatthew G. Knepley 
2204dce8aebaSBarry Smith .seealso: `PetscDS`, `PetscDSGetResidual()`
220532d2bbc9SMatthew G. Knepley @*/
2206d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDSSetUpdate(PetscDS ds, PetscInt f, void (*update)(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 uNew[]))
2207d71ae5a4SJacob Faibussowitsch {
220832d2bbc9SMatthew G. Knepley   PetscFunctionBegin;
22096528b96dSMatthew G. Knepley   PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1);
221032d2bbc9SMatthew G. Knepley   if (update) PetscValidFunction(update, 3);
221163a3b9bcSJacob Faibussowitsch   PetscCheck(f >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %" PetscInt_FMT " must be non-negative", f);
22129566063dSJacob Faibussowitsch   PetscCall(PetscDSEnlarge_Static(ds, f + 1));
22136528b96dSMatthew G. Knepley   ds->update[f] = update;
22143ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
221532d2bbc9SMatthew G. Knepley }
221632d2bbc9SMatthew G. Knepley 
2217d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDSGetContext(PetscDS ds, PetscInt f, void *ctx)
2218d71ae5a4SJacob Faibussowitsch {
22190c2f2876SMatthew G. Knepley   PetscFunctionBegin;
22206528b96dSMatthew G. Knepley   PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1);
222163a3b9bcSJacob Faibussowitsch   PetscCheck(!(f < 0) && !(f >= ds->Nf), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %" PetscInt_FMT " must be in [0, %" PetscInt_FMT ")", f, ds->Nf);
22224f572ea9SToby Isaac   PetscAssertPointer(ctx, 3);
22233ec1f749SStefano Zampini   *(void **)ctx = ds->ctx[f];
22243ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
22250c2f2876SMatthew G. Knepley }
22260c2f2876SMatthew G. Knepley 
2227d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDSSetContext(PetscDS ds, PetscInt f, void *ctx)
2228d71ae5a4SJacob Faibussowitsch {
22290c2f2876SMatthew G. Knepley   PetscFunctionBegin;
22306528b96dSMatthew G. Knepley   PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1);
223163a3b9bcSJacob Faibussowitsch   PetscCheck(f >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %" PetscInt_FMT " must be non-negative", f);
22329566063dSJacob Faibussowitsch   PetscCall(PetscDSEnlarge_Static(ds, f + 1));
22336528b96dSMatthew G. Knepley   ds->ctx[f] = ctx;
22343ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
22350c2f2876SMatthew G. Knepley }
22360c2f2876SMatthew G. Knepley 
2237194d53e6SMatthew G. Knepley /*@C
2238194d53e6SMatthew G. Knepley   PetscDSGetBdResidual - Get the pointwise boundary residual function for a given test field
2239194d53e6SMatthew G. Knepley 
224020f4b53cSBarry Smith   Not Collective
2241194d53e6SMatthew G. Knepley 
2242194d53e6SMatthew G. Knepley   Input Parameters:
22436528b96dSMatthew G. Knepley + ds - The PetscDS
2244194d53e6SMatthew G. Knepley - f  - The test field number
2245194d53e6SMatthew G. Knepley 
2246194d53e6SMatthew G. Knepley   Output Parameters:
2247194d53e6SMatthew G. Knepley + f0 - boundary integrand for the test function term
2248194d53e6SMatthew G. Knepley - f1 - boundary integrand for the test function gradient term
2249194d53e6SMatthew G. Knepley 
2250*a4e35b19SJacob Faibussowitsch   Calling sequence of `f0`:
2251194d53e6SMatthew G. Knepley + dim          - the spatial dimension
2252194d53e6SMatthew G. Knepley . Nf           - the number of fields
2253*a4e35b19SJacob Faibussowitsch . NfAux        - the number of auxiliary fields
2254194d53e6SMatthew G. Knepley . uOff         - the offset into u[] and u_t[] for each field
2255194d53e6SMatthew G. Knepley . uOff_x       - the offset into u_x[] for each field
2256194d53e6SMatthew G. Knepley . u            - each field evaluated at the current point
2257194d53e6SMatthew G. Knepley . u_t          - the time derivative of each field evaluated at the current point
2258194d53e6SMatthew G. Knepley . u_x          - the gradient of each field evaluated at the current point
2259194d53e6SMatthew G. Knepley . aOff         - the offset into a[] and a_t[] for each auxiliary field
2260194d53e6SMatthew G. Knepley . aOff_x       - the offset into a_x[] for each auxiliary field
2261194d53e6SMatthew G. Knepley . a            - each auxiliary field evaluated at the current point
2262194d53e6SMatthew G. Knepley . a_t          - the time derivative of each auxiliary field evaluated at the current point
2263194d53e6SMatthew G. Knepley . a_x          - the gradient of auxiliary each field evaluated at the current point
2264194d53e6SMatthew G. Knepley . t            - current time
2265194d53e6SMatthew G. Knepley . x            - coordinates of the current point
2266194d53e6SMatthew G. Knepley . n            - unit normal at the current point
226797b6e6e8SMatthew G. Knepley . numConstants - number of constant parameters
226897b6e6e8SMatthew G. Knepley . constants    - constant parameters
2269194d53e6SMatthew G. Knepley - f0           - output values at the current point
2270194d53e6SMatthew G. Knepley 
2271194d53e6SMatthew G. Knepley   Level: intermediate
2272194d53e6SMatthew G. Knepley 
2273dce8aebaSBarry Smith   Note:
2274*a4e35b19SJacob Faibussowitsch   The calling sequence of `f1` is identical, and therefore omitted for brevity.
227560225df5SJacob Faibussowitsch 
2276*a4e35b19SJacob Faibussowitsch   We are using a first order FEM model for the weak form\:
2277dce8aebaSBarry Smith   \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
2278dce8aebaSBarry Smith 
2279dce8aebaSBarry Smith .seealso: `PetscDS`, `PetscDSSetBdResidual()`
2280194d53e6SMatthew G. Knepley @*/
2281*a4e35b19SJacob Faibussowitsch PetscErrorCode PetscDSGetBdResidual(PetscDS ds, 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, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
2282d71ae5a4SJacob Faibussowitsch {
22836528b96dSMatthew G. Knepley   PetscBdPointFunc *tmp0, *tmp1;
22846528b96dSMatthew G. Knepley   PetscInt          n0, n1;
22856528b96dSMatthew G. Knepley 
22862764a2aaSMatthew G. Knepley   PetscFunctionBegin;
22876528b96dSMatthew G. Knepley   PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1);
228863a3b9bcSJacob Faibussowitsch   PetscCheck(!(f < 0) && !(f >= ds->Nf), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %" PetscInt_FMT " must be in [0, %" PetscInt_FMT ")", f, ds->Nf);
22899566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormGetBdResidual(ds->wf, NULL, 0, f, 0, &n0, &tmp0, &n1, &tmp1));
22906528b96dSMatthew G. Knepley   *f0 = tmp0 ? tmp0[0] : NULL;
22916528b96dSMatthew G. Knepley   *f1 = tmp1 ? tmp1[0] : NULL;
22923ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
22932764a2aaSMatthew G. Knepley }
22942764a2aaSMatthew G. Knepley 
2295194d53e6SMatthew G. Knepley /*@C
2296194d53e6SMatthew G. Knepley   PetscDSSetBdResidual - Get the pointwise boundary residual function for a given test field
2297194d53e6SMatthew G. Knepley 
229820f4b53cSBarry Smith   Not Collective
2299194d53e6SMatthew G. Knepley 
2300194d53e6SMatthew G. Knepley   Input Parameters:
2301dce8aebaSBarry Smith + ds - The `PetscDS`
2302194d53e6SMatthew G. Knepley . f  - The test field number
2303194d53e6SMatthew G. Knepley . f0 - boundary integrand for the test function term
2304194d53e6SMatthew G. Knepley - f1 - boundary integrand for the test function gradient term
2305194d53e6SMatthew G. Knepley 
2306*a4e35b19SJacob Faibussowitsch   Calling sequence of `f0`:
2307194d53e6SMatthew G. Knepley + dim          - the spatial dimension
2308194d53e6SMatthew G. Knepley . Nf           - the number of fields
2309*a4e35b19SJacob Faibussowitsch . NfAux        - the number of auxiliary fields
2310194d53e6SMatthew G. Knepley . uOff         - the offset into u[] and u_t[] for each field
2311194d53e6SMatthew G. Knepley . uOff_x       - the offset into u_x[] for each field
2312194d53e6SMatthew G. Knepley . u            - each field evaluated at the current point
2313194d53e6SMatthew G. Knepley . u_t          - the time derivative of each field evaluated at the current point
2314194d53e6SMatthew G. Knepley . u_x          - the gradient of each field evaluated at the current point
2315194d53e6SMatthew G. Knepley . aOff         - the offset into a[] and a_t[] for each auxiliary field
2316194d53e6SMatthew G. Knepley . aOff_x       - the offset into a_x[] for each auxiliary field
2317194d53e6SMatthew G. Knepley . a            - each auxiliary field evaluated at the current point
2318194d53e6SMatthew G. Knepley . a_t          - the time derivative of each auxiliary field evaluated at the current point
2319194d53e6SMatthew G. Knepley . a_x          - the gradient of auxiliary each field evaluated at the current point
2320194d53e6SMatthew G. Knepley . t            - current time
2321194d53e6SMatthew G. Knepley . x            - coordinates of the current point
2322194d53e6SMatthew G. Knepley . n            - unit normal at the current point
232397b6e6e8SMatthew G. Knepley . numConstants - number of constant parameters
232497b6e6e8SMatthew G. Knepley . constants    - constant parameters
2325194d53e6SMatthew G. Knepley - f0           - output values at the current point
2326194d53e6SMatthew G. Knepley 
2327194d53e6SMatthew G. Knepley   Level: intermediate
2328194d53e6SMatthew G. Knepley 
2329dce8aebaSBarry Smith   Note:
2330*a4e35b19SJacob Faibussowitsch   The calling sequence of `f1` is identical, and therefore omitted for brevity.
233160225df5SJacob Faibussowitsch 
2332*a4e35b19SJacob Faibussowitsch   We are using a first order FEM model for the weak form\:
2333dce8aebaSBarry Smith   \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
2334dce8aebaSBarry Smith 
2335dce8aebaSBarry Smith .seealso: `PetscDS`, `PetscDSGetBdResidual()`
2336194d53e6SMatthew G. Knepley @*/
2337*a4e35b19SJacob Faibussowitsch PetscErrorCode PetscDSSetBdResidual(PetscDS ds, 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, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
2338d71ae5a4SJacob Faibussowitsch {
23392764a2aaSMatthew G. Knepley   PetscFunctionBegin;
23406528b96dSMatthew G. Knepley   PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1);
234163a3b9bcSJacob Faibussowitsch   PetscCheck(f >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %" PetscInt_FMT " must be non-negative", f);
23429566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormSetIndexBdResidual(ds->wf, NULL, 0, f, 0, 0, f0, 0, f1));
23433ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
23442764a2aaSMatthew G. Knepley }
23452764a2aaSMatthew G. Knepley 
234627f02ce8SMatthew G. Knepley /*@
2347dce8aebaSBarry Smith   PetscDSHasBdJacobian - Indicates that boundary Jacobian functions have been set
234827f02ce8SMatthew G. Knepley 
234920f4b53cSBarry Smith   Not Collective
235027f02ce8SMatthew G. Knepley 
235127f02ce8SMatthew G. Knepley   Input Parameter:
2352dce8aebaSBarry Smith . ds - The `PetscDS`
235327f02ce8SMatthew G. Knepley 
235427f02ce8SMatthew G. Knepley   Output Parameter:
235527f02ce8SMatthew G. Knepley . hasBdJac - flag that pointwise function for the boundary Jacobian has been set
235627f02ce8SMatthew G. Knepley 
235727f02ce8SMatthew G. Knepley   Level: intermediate
235827f02ce8SMatthew G. Knepley 
2359dce8aebaSBarry Smith .seealso: `PetscDS`, `PetscDSHasJacobian()`, `PetscDSSetBdJacobian()`, `PetscDSGetBdJacobian()`
236027f02ce8SMatthew G. Knepley @*/
2361d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDSHasBdJacobian(PetscDS ds, PetscBool *hasBdJac)
2362d71ae5a4SJacob Faibussowitsch {
236327f02ce8SMatthew G. Knepley   PetscFunctionBegin;
23646528b96dSMatthew G. Knepley   PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1);
23654f572ea9SToby Isaac   PetscAssertPointer(hasBdJac, 2);
23669566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormHasBdJacobian(ds->wf, hasBdJac));
23673ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
236827f02ce8SMatthew G. Knepley }
236927f02ce8SMatthew G. Knepley 
2370194d53e6SMatthew G. Knepley /*@C
2371194d53e6SMatthew G. Knepley   PetscDSGetBdJacobian - Get the pointwise boundary Jacobian function for given test and basis field
2372194d53e6SMatthew G. Knepley 
237320f4b53cSBarry Smith   Not Collective
2374194d53e6SMatthew G. Knepley 
2375194d53e6SMatthew G. Knepley   Input Parameters:
2376dce8aebaSBarry Smith + ds - The `PetscDS`
2377194d53e6SMatthew G. Knepley . f  - The test field number
2378194d53e6SMatthew G. Knepley - g  - The field number
2379194d53e6SMatthew G. Knepley 
2380194d53e6SMatthew G. Knepley   Output Parameters:
2381194d53e6SMatthew G. Knepley + g0 - integrand for the test and basis function term
2382194d53e6SMatthew G. Knepley . g1 - integrand for the test function and basis function gradient term
2383194d53e6SMatthew G. Knepley . g2 - integrand for the test function gradient and basis function term
2384194d53e6SMatthew G. Knepley - g3 - integrand for the test function gradient and basis function gradient term
2385194d53e6SMatthew G. Knepley 
2386*a4e35b19SJacob Faibussowitsch   Calling sequence of `g0`:
2387194d53e6SMatthew G. Knepley + dim          - the spatial dimension
2388194d53e6SMatthew G. Knepley . Nf           - the number of fields
2389*a4e35b19SJacob Faibussowitsch . NfAux        - the number of auxiliary fields
2390194d53e6SMatthew G. Knepley . uOff         - the offset into u[] and u_t[] for each field
2391194d53e6SMatthew G. Knepley . uOff_x       - the offset into u_x[] for each field
2392194d53e6SMatthew G. Knepley . u            - each field evaluated at the current point
2393194d53e6SMatthew G. Knepley . u_t          - the time derivative of each field evaluated at the current point
2394194d53e6SMatthew G. Knepley . u_x          - the gradient of each field evaluated at the current point
2395194d53e6SMatthew G. Knepley . aOff         - the offset into a[] and a_t[] for each auxiliary field
2396194d53e6SMatthew G. Knepley . aOff_x       - the offset into a_x[] for each auxiliary field
2397194d53e6SMatthew G. Knepley . a            - each auxiliary field evaluated at the current point
2398194d53e6SMatthew G. Knepley . a_t          - the time derivative of each auxiliary field evaluated at the current point
2399194d53e6SMatthew G. Knepley . a_x          - the gradient of auxiliary each field evaluated at the current point
2400194d53e6SMatthew G. Knepley . t            - current time
24012aa1fc23SMatthew G. Knepley . u_tShift     - the multiplier a for dF/dU_t
2402194d53e6SMatthew G. Knepley . x            - coordinates of the current point
2403194d53e6SMatthew G. Knepley . n            - normal at the current point
240497b6e6e8SMatthew G. Knepley . numConstants - number of constant parameters
240597b6e6e8SMatthew G. Knepley . constants    - constant parameters
2406194d53e6SMatthew G. Knepley - g0           - output values at the current point
2407194d53e6SMatthew G. Knepley 
2408194d53e6SMatthew G. Knepley   Level: intermediate
2409194d53e6SMatthew G. Knepley 
2410dce8aebaSBarry Smith   Note:
2411*a4e35b19SJacob Faibussowitsch   `g1`, `g2`, and `g3` have identical calling sequences to `g0` and are omitted for brevity.
241260225df5SJacob Faibussowitsch 
2413*a4e35b19SJacob Faibussowitsch   We are using a first order FEM model for the weak form\:
2414dce8aebaSBarry Smith   \int_\Gamma \phi {\vec g}_0(u, u_t, \nabla u, x, t) \cdot \hat n \psi + \phi {\vec g}_1(u, u_t, \nabla u, x, t) \cdot \hat n \nabla \psi + \nabla\phi \cdot {\vec g}_2(u, u_t, \nabla u, x, t) \cdot \hat n \psi + \nabla\phi \cdot {\overleftrightarrow g}_3(u, u_t, \nabla u, x, t) \cdot \hat n \cdot \nabla \psi
2415dce8aebaSBarry Smith 
2416dce8aebaSBarry Smith .seealso: `PetscDS`, `PetscDSSetBdJacobian()`
2417194d53e6SMatthew G. Knepley @*/
2418*a4e35b19SJacob Faibussowitsch PetscErrorCode PetscDSGetBdJacobian(PetscDS ds, PetscInt f, PetscInt g, void (**g0)(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, PetscReal u_tShift, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]), void (**g1)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), void (**g2)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), void (**g3)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
2419d71ae5a4SJacob Faibussowitsch {
24206528b96dSMatthew G. Knepley   PetscBdPointJac *tmp0, *tmp1, *tmp2, *tmp3;
24216528b96dSMatthew G. Knepley   PetscInt         n0, n1, n2, n3;
24226528b96dSMatthew G. Knepley 
24232764a2aaSMatthew G. Knepley   PetscFunctionBegin;
24246528b96dSMatthew G. Knepley   PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1);
242563a3b9bcSJacob Faibussowitsch   PetscCheck(!(f < 0) && !(f >= ds->Nf), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %" PetscInt_FMT " must be in [0, %" PetscInt_FMT ")", f, ds->Nf);
242663a3b9bcSJacob Faibussowitsch   PetscCheck(!(g < 0) && !(g >= ds->Nf), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %" PetscInt_FMT " must be in [0, %" PetscInt_FMT ")", g, ds->Nf);
24279566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormGetBdJacobian(ds->wf, NULL, 0, f, g, 0, &n0, &tmp0, &n1, &tmp1, &n2, &tmp2, &n3, &tmp3));
24286528b96dSMatthew G. Knepley   *g0 = tmp0 ? tmp0[0] : NULL;
24296528b96dSMatthew G. Knepley   *g1 = tmp1 ? tmp1[0] : NULL;
24306528b96dSMatthew G. Knepley   *g2 = tmp2 ? tmp2[0] : NULL;
24316528b96dSMatthew G. Knepley   *g3 = tmp3 ? tmp3[0] : NULL;
24323ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
24332764a2aaSMatthew G. Knepley }
24342764a2aaSMatthew G. Knepley 
2435194d53e6SMatthew G. Knepley /*@C
2436194d53e6SMatthew G. Knepley   PetscDSSetBdJacobian - Set the pointwise boundary Jacobian function for given test and basis field
2437194d53e6SMatthew G. Knepley 
243820f4b53cSBarry Smith   Not Collective
2439194d53e6SMatthew G. Knepley 
2440194d53e6SMatthew G. Knepley   Input Parameters:
24416528b96dSMatthew G. Knepley + ds - The PetscDS
2442194d53e6SMatthew G. Knepley . f  - The test field number
2443194d53e6SMatthew G. Knepley . g  - The field number
2444194d53e6SMatthew G. Knepley . g0 - integrand for the test and basis function term
2445194d53e6SMatthew G. Knepley . g1 - integrand for the test function and basis function gradient term
2446194d53e6SMatthew G. Knepley . g2 - integrand for the test function gradient and basis function term
2447194d53e6SMatthew G. Knepley - g3 - integrand for the test function gradient and basis function gradient term
2448194d53e6SMatthew G. Knepley 
2449*a4e35b19SJacob Faibussowitsch   Calling sequence of `g0`:
2450194d53e6SMatthew G. Knepley + dim          - the spatial dimension
2451194d53e6SMatthew G. Knepley . Nf           - the number of fields
2452*a4e35b19SJacob Faibussowitsch . NfAux        - the number of auxiliary fields
2453194d53e6SMatthew G. Knepley . uOff         - the offset into u[] and u_t[] for each field
2454194d53e6SMatthew G. Knepley . uOff_x       - the offset into u_x[] for each field
2455194d53e6SMatthew G. Knepley . u            - each field evaluated at the current point
2456194d53e6SMatthew G. Knepley . u_t          - the time derivative of each field evaluated at the current point
2457194d53e6SMatthew G. Knepley . u_x          - the gradient of each field evaluated at the current point
2458194d53e6SMatthew G. Knepley . aOff         - the offset into a[] and a_t[] for each auxiliary field
2459194d53e6SMatthew G. Knepley . aOff_x       - the offset into a_x[] for each auxiliary field
2460194d53e6SMatthew G. Knepley . a            - each auxiliary field evaluated at the current point
2461194d53e6SMatthew G. Knepley . a_t          - the time derivative of each auxiliary field evaluated at the current point
2462194d53e6SMatthew G. Knepley . a_x          - the gradient of auxiliary each field evaluated at the current point
2463194d53e6SMatthew G. Knepley . t            - current time
24642aa1fc23SMatthew G. Knepley . u_tShift     - the multiplier a for dF/dU_t
2465194d53e6SMatthew G. Knepley . x            - coordinates of the current point
2466194d53e6SMatthew G. Knepley . n            - normal at the current point
246797b6e6e8SMatthew G. Knepley . numConstants - number of constant parameters
246897b6e6e8SMatthew G. Knepley . constants    - constant parameters
2469194d53e6SMatthew G. Knepley - g0           - output values at the current point
2470194d53e6SMatthew G. Knepley 
2471194d53e6SMatthew G. Knepley   Level: intermediate
2472194d53e6SMatthew G. Knepley 
2473dce8aebaSBarry Smith   Note:
2474*a4e35b19SJacob Faibussowitsch   `g1`, `g2`, and `g3` have identical calling sequences to `g0` and are omitted for brevity.
247560225df5SJacob Faibussowitsch 
2476*a4e35b19SJacob Faibussowitsch   We are using a first order FEM model for the weak form\:
2477dce8aebaSBarry Smith   \int_\Gamma \phi {\vec g}_0(u, u_t, \nabla u, x, t) \cdot \hat n \psi + \phi {\vec g}_1(u, u_t, \nabla u, x, t) \cdot \hat n \nabla \psi + \nabla\phi \cdot {\vec g}_2(u, u_t, \nabla u, x, t) \cdot \hat n \psi + \nabla\phi \cdot {\overleftrightarrow g}_3(u, u_t, \nabla u, x, t) \cdot \hat n \cdot \nabla \psi
2478dce8aebaSBarry Smith 
2479dce8aebaSBarry Smith .seealso: `PetscDS`, `PetscDSGetBdJacobian()`
2480194d53e6SMatthew G. Knepley @*/
2481*a4e35b19SJacob Faibussowitsch PetscErrorCode PetscDSSetBdJacobian(PetscDS ds, PetscInt f, PetscInt g, void (*g0)(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, PetscReal u_tShift, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]), void (*g1)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), void (*g2)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), void (*g3)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
2482d71ae5a4SJacob Faibussowitsch {
24832764a2aaSMatthew G. Knepley   PetscFunctionBegin;
24846528b96dSMatthew G. Knepley   PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1);
24852764a2aaSMatthew G. Knepley   if (g0) PetscValidFunction(g0, 4);
24862764a2aaSMatthew G. Knepley   if (g1) PetscValidFunction(g1, 5);
24872764a2aaSMatthew G. Knepley   if (g2) PetscValidFunction(g2, 6);
24882764a2aaSMatthew G. Knepley   if (g3) PetscValidFunction(g3, 7);
248963a3b9bcSJacob Faibussowitsch   PetscCheck(f >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %" PetscInt_FMT " must be non-negative", f);
249063a3b9bcSJacob Faibussowitsch   PetscCheck(g >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %" PetscInt_FMT " must be non-negative", g);
24919566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormSetIndexBdJacobian(ds->wf, NULL, 0, f, g, 0, 0, g0, 0, g1, 0, g2, 0, g3));
24923ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
24932764a2aaSMatthew G. Knepley }
24942764a2aaSMatthew G. Knepley 
249527f02ce8SMatthew G. Knepley /*@
249627f02ce8SMatthew G. Knepley   PetscDSHasBdJacobianPreconditioner - Signals that boundary Jacobian preconditioner functions have been set
249727f02ce8SMatthew G. Knepley 
249820f4b53cSBarry Smith   Not Collective
249927f02ce8SMatthew G. Knepley 
250027f02ce8SMatthew G. Knepley   Input Parameter:
2501dce8aebaSBarry Smith . ds - The `PetscDS`
250227f02ce8SMatthew G. Knepley 
250327f02ce8SMatthew G. Knepley   Output Parameter:
250460225df5SJacob Faibussowitsch . hasBdJacPre - flag that pointwise function for the boundary Jacobian preconditioner has been set
250527f02ce8SMatthew G. Knepley 
250627f02ce8SMatthew G. Knepley   Level: intermediate
250727f02ce8SMatthew G. Knepley 
2508dce8aebaSBarry Smith .seealso: `PetscDS`, `PetscDSHasJacobian()`, `PetscDSSetBdJacobian()`, `PetscDSGetBdJacobian()`
250927f02ce8SMatthew G. Knepley @*/
2510d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDSHasBdJacobianPreconditioner(PetscDS ds, PetscBool *hasBdJacPre)
2511d71ae5a4SJacob Faibussowitsch {
251227f02ce8SMatthew G. Knepley   PetscFunctionBegin;
25136528b96dSMatthew G. Knepley   PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1);
25144f572ea9SToby Isaac   PetscAssertPointer(hasBdJacPre, 2);
25159566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormHasBdJacobianPreconditioner(ds->wf, hasBdJacPre));
25163ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
251727f02ce8SMatthew G. Knepley }
251827f02ce8SMatthew G. Knepley 
251927f02ce8SMatthew G. Knepley /*@C
252027f02ce8SMatthew G. Knepley   PetscDSGetBdJacobianPreconditioner - Get the pointwise boundary Jacobian preconditioner function for given test and basis field
252127f02ce8SMatthew G. Knepley 
252220f4b53cSBarry Smith   Not Collective; No Fortran Support
252327f02ce8SMatthew G. Knepley 
252427f02ce8SMatthew G. Knepley   Input Parameters:
2525dce8aebaSBarry Smith + ds - The `PetscDS`
252627f02ce8SMatthew G. Knepley . f  - The test field number
252727f02ce8SMatthew G. Knepley - g  - The field number
252827f02ce8SMatthew G. Knepley 
252927f02ce8SMatthew G. Knepley   Output Parameters:
253027f02ce8SMatthew G. Knepley + g0 - integrand for the test and basis function term
253127f02ce8SMatthew G. Knepley . g1 - integrand for the test function and basis function gradient term
253227f02ce8SMatthew G. Knepley . g2 - integrand for the test function gradient and basis function term
253327f02ce8SMatthew G. Knepley - g3 - integrand for the test function gradient and basis function gradient term
253427f02ce8SMatthew G. Knepley 
2535*a4e35b19SJacob Faibussowitsch   Calling sequence of `g0`:
253627f02ce8SMatthew G. Knepley + dim          - the spatial dimension
253727f02ce8SMatthew G. Knepley . Nf           - the number of fields
253827f02ce8SMatthew G. Knepley . NfAux        - the number of auxiliary fields
253927f02ce8SMatthew G. Knepley . uOff         - the offset into u[] and u_t[] for each field
254027f02ce8SMatthew G. Knepley . uOff_x       - the offset into u_x[] for each field
254127f02ce8SMatthew G. Knepley . u            - each field evaluated at the current point
254227f02ce8SMatthew G. Knepley . u_t          - the time derivative of each field evaluated at the current point
254327f02ce8SMatthew G. Knepley . u_x          - the gradient of each field evaluated at the current point
254427f02ce8SMatthew G. Knepley . aOff         - the offset into a[] and a_t[] for each auxiliary field
254527f02ce8SMatthew G. Knepley . aOff_x       - the offset into a_x[] for each auxiliary field
254627f02ce8SMatthew G. Knepley . a            - each auxiliary field evaluated at the current point
254727f02ce8SMatthew G. Knepley . a_t          - the time derivative of each auxiliary field evaluated at the current point
254827f02ce8SMatthew G. Knepley . a_x          - the gradient of auxiliary each field evaluated at the current point
254927f02ce8SMatthew G. Knepley . t            - current time
255027f02ce8SMatthew G. Knepley . u_tShift     - the multiplier a for dF/dU_t
255127f02ce8SMatthew G. Knepley . x            - coordinates of the current point
255227f02ce8SMatthew G. Knepley . n            - normal at the current point
255327f02ce8SMatthew G. Knepley . numConstants - number of constant parameters
255427f02ce8SMatthew G. Knepley . constants    - constant parameters
255527f02ce8SMatthew G. Knepley - g0           - output values at the current point
255627f02ce8SMatthew G. Knepley 
255727f02ce8SMatthew G. Knepley   Level: intermediate
255827f02ce8SMatthew G. Knepley 
2559dce8aebaSBarry Smith   Note:
2560*a4e35b19SJacob Faibussowitsch   `g1`, `g2`, and `g3` have identical calling sequences to `g0` and are omitted for brevity.
256160225df5SJacob Faibussowitsch 
2562*a4e35b19SJacob Faibussowitsch   We are using a first order FEM model for the weak form\:
2563dce8aebaSBarry Smith   \int_\Gamma \phi {\vec g}_0(u, u_t, \nabla u, x, t) \cdot \hat n \psi + \phi {\vec g}_1(u, u_t, \nabla u, x, t) \cdot \hat n \nabla \psi + \nabla\phi \cdot {\vec g}_2(u, u_t, \nabla u, x, t) \cdot \hat n \psi + \nabla\phi \cdot {\overleftrightarrow g}_3(u, u_t, \nabla u, x, t) \cdot \hat n \cdot \nabla \psi
2564dce8aebaSBarry Smith 
2565dce8aebaSBarry Smith .seealso: `PetscDS`, `PetscDSSetBdJacobianPreconditioner()`
256627f02ce8SMatthew G. Knepley @*/
2567*a4e35b19SJacob Faibussowitsch PetscErrorCode PetscDSGetBdJacobianPreconditioner(PetscDS ds, PetscInt f, PetscInt g, void (**g0)(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, PetscReal u_tShift, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]), void (**g1)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), void (**g2)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), void (**g3)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
2568d71ae5a4SJacob Faibussowitsch {
25696528b96dSMatthew G. Knepley   PetscBdPointJac *tmp0, *tmp1, *tmp2, *tmp3;
25706528b96dSMatthew G. Knepley   PetscInt         n0, n1, n2, n3;
25716528b96dSMatthew G. Knepley 
257227f02ce8SMatthew G. Knepley   PetscFunctionBegin;
25736528b96dSMatthew G. Knepley   PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1);
257463a3b9bcSJacob Faibussowitsch   PetscCheck(!(f < 0) && !(f >= ds->Nf), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %" PetscInt_FMT " must be in [0, %" PetscInt_FMT ")", f, ds->Nf);
257563a3b9bcSJacob Faibussowitsch   PetscCheck(!(g < 0) && !(g >= ds->Nf), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %" PetscInt_FMT " must be in [0, %" PetscInt_FMT ")", g, ds->Nf);
25769566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormGetBdJacobianPreconditioner(ds->wf, NULL, 0, f, g, 0, &n0, &tmp0, &n1, &tmp1, &n2, &tmp2, &n3, &tmp3));
25776528b96dSMatthew G. Knepley   *g0 = tmp0 ? tmp0[0] : NULL;
25786528b96dSMatthew G. Knepley   *g1 = tmp1 ? tmp1[0] : NULL;
25796528b96dSMatthew G. Knepley   *g2 = tmp2 ? tmp2[0] : NULL;
25806528b96dSMatthew G. Knepley   *g3 = tmp3 ? tmp3[0] : NULL;
25813ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
258227f02ce8SMatthew G. Knepley }
258327f02ce8SMatthew G. Knepley 
258427f02ce8SMatthew G. Knepley /*@C
258527f02ce8SMatthew G. Knepley   PetscDSSetBdJacobianPreconditioner - Set the pointwise boundary Jacobian preconditioner function for given test and basis field
258627f02ce8SMatthew G. Knepley 
258720f4b53cSBarry Smith   Not Collective; No Fortran Support
258827f02ce8SMatthew G. Knepley 
258927f02ce8SMatthew G. Knepley   Input Parameters:
2590dce8aebaSBarry Smith + ds - The `PetscDS`
259127f02ce8SMatthew G. Knepley . f  - The test field number
259227f02ce8SMatthew G. Knepley . g  - The field number
259327f02ce8SMatthew G. Knepley . g0 - integrand for the test and basis function term
259427f02ce8SMatthew G. Knepley . g1 - integrand for the test function and basis function gradient term
259527f02ce8SMatthew G. Knepley . g2 - integrand for the test function gradient and basis function term
259627f02ce8SMatthew G. Knepley - g3 - integrand for the test function gradient and basis function gradient term
259727f02ce8SMatthew G. Knepley 
2598*a4e35b19SJacob Faibussowitsch   Calling sequence of `g0':
259927f02ce8SMatthew G. Knepley + dim          - the spatial dimension
260027f02ce8SMatthew G. Knepley . Nf           - the number of fields
260127f02ce8SMatthew G. Knepley . NfAux        - the number of auxiliary fields
260227f02ce8SMatthew G. Knepley . uOff         - the offset into u[] and u_t[] for each field
260327f02ce8SMatthew G. Knepley . uOff_x       - the offset into u_x[] for each field
260427f02ce8SMatthew G. Knepley . u            - each field evaluated at the current point
260527f02ce8SMatthew G. Knepley . u_t          - the time derivative of each field evaluated at the current point
260627f02ce8SMatthew G. Knepley . u_x          - the gradient of each field evaluated at the current point
260727f02ce8SMatthew G. Knepley . aOff         - the offset into a[] and a_t[] for each auxiliary field
260827f02ce8SMatthew G. Knepley . aOff_x       - the offset into a_x[] for each auxiliary field
260927f02ce8SMatthew G. Knepley . a            - each auxiliary field evaluated at the current point
261027f02ce8SMatthew G. Knepley . a_t          - the time derivative of each auxiliary field evaluated at the current point
261127f02ce8SMatthew G. Knepley . a_x          - the gradient of auxiliary each field evaluated at the current point
261227f02ce8SMatthew G. Knepley . t            - current time
261327f02ce8SMatthew G. Knepley . u_tShift     - the multiplier a for dF/dU_t
261427f02ce8SMatthew G. Knepley . x            - coordinates of the current point
261527f02ce8SMatthew G. Knepley . n            - normal at the current point
261627f02ce8SMatthew G. Knepley . numConstants - number of constant parameters
261727f02ce8SMatthew G. Knepley . constants    - constant parameters
261827f02ce8SMatthew G. Knepley - g0           - output values at the current point
261927f02ce8SMatthew G. Knepley 
262027f02ce8SMatthew G. Knepley   Level: intermediate
262127f02ce8SMatthew G. Knepley 
2622dce8aebaSBarry Smith   Note:
2623*a4e35b19SJacob Faibussowitsch   `g1`, `g2`, and `g3` have identical calling sequences to `g0` and are omitted for brevity.
262460225df5SJacob Faibussowitsch 
2625*a4e35b19SJacob Faibussowitsch   We are using a first order FEM model for the weak form\:
2626dce8aebaSBarry Smith   \int_\Gamma \phi {\vec g}_0(u, u_t, \nabla u, x, t) \cdot \hat n \psi + \phi {\vec g}_1(u, u_t, \nabla u, x, t) \cdot \hat n \nabla \psi + \nabla\phi \cdot {\vec g}_2(u, u_t, \nabla u, x, t) \cdot \hat n \psi + \nabla\phi \cdot {\overleftrightarrow g}_3(u, u_t, \nabla u, x, t) \cdot \hat n \cdot \nabla \psi
2627dce8aebaSBarry Smith 
2628dce8aebaSBarry Smith .seealso: `PetscDS`, `PetscDSGetBdJacobianPreconditioner()`
262927f02ce8SMatthew G. Knepley @*/
2630*a4e35b19SJacob Faibussowitsch PetscErrorCode PetscDSSetBdJacobianPreconditioner(PetscDS ds, PetscInt f, PetscInt g, void (*g0)(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, PetscReal u_tShift, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]), void (*g1)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), void (*g2)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), void (*g3)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
2631d71ae5a4SJacob Faibussowitsch {
263227f02ce8SMatthew G. Knepley   PetscFunctionBegin;
26336528b96dSMatthew G. Knepley   PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1);
263427f02ce8SMatthew G. Knepley   if (g0) PetscValidFunction(g0, 4);
263527f02ce8SMatthew G. Knepley   if (g1) PetscValidFunction(g1, 5);
263627f02ce8SMatthew G. Knepley   if (g2) PetscValidFunction(g2, 6);
263727f02ce8SMatthew G. Knepley   if (g3) PetscValidFunction(g3, 7);
263863a3b9bcSJacob Faibussowitsch   PetscCheck(f >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %" PetscInt_FMT " must be non-negative", f);
263963a3b9bcSJacob Faibussowitsch   PetscCheck(g >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %" PetscInt_FMT " must be non-negative", g);
26409566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormSetIndexBdJacobianPreconditioner(ds->wf, NULL, 0, f, g, 0, 0, g0, 0, g1, 0, g2, 0, g3));
26413ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
264227f02ce8SMatthew G. Knepley }
264327f02ce8SMatthew G. Knepley 
26440d3e9b51SMatthew G. Knepley /*@C
2645c371a6d1SMatthew G. Knepley   PetscDSGetExactSolution - Get the pointwise exact solution function for a given test field
2646c371a6d1SMatthew G. Knepley 
264720f4b53cSBarry Smith   Not Collective
2648c371a6d1SMatthew G. Knepley 
2649c371a6d1SMatthew G. Knepley   Input Parameters:
2650c371a6d1SMatthew G. Knepley + prob - The PetscDS
2651c371a6d1SMatthew G. Knepley - f    - The test field number
2652c371a6d1SMatthew G. Knepley 
2653d8d19677SJose E. Roman   Output Parameters:
2654*a4e35b19SJacob Faibussowitsch + sol - exact solution for the test field
2655*a4e35b19SJacob Faibussowitsch - ctx - exact solution context
2656c371a6d1SMatthew G. Knepley 
265720f4b53cSBarry Smith   Calling sequence of `exactSol`:
2658c371a6d1SMatthew G. Knepley + dim - the spatial dimension
2659c371a6d1SMatthew G. Knepley . t   - current time
2660c371a6d1SMatthew G. Knepley . x   - coordinates of the current point
2661c371a6d1SMatthew G. Knepley . Nc  - the number of field components
2662*a4e35b19SJacob Faibussowitsch . u   - the solution field evaluated at the current point
2663c371a6d1SMatthew G. Knepley - ctx - a user context
2664c371a6d1SMatthew G. Knepley 
2665c371a6d1SMatthew G. Knepley   Level: intermediate
2666c371a6d1SMatthew G. Knepley 
2667dce8aebaSBarry Smith .seealso: `PetscDS`, `PetscDSSetExactSolution()`, `PetscDSGetExactSolutionTimeDerivative()`
2668c371a6d1SMatthew G. Knepley @*/
2669d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDSGetExactSolution(PetscDS prob, PetscInt f, PetscErrorCode (**sol)(PetscInt dim, PetscReal t, const PetscReal x[], PetscInt Nc, PetscScalar u[], void *ctx), void **ctx)
2670d71ae5a4SJacob Faibussowitsch {
2671c371a6d1SMatthew G. Knepley   PetscFunctionBegin;
2672c371a6d1SMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
267363a3b9bcSJacob Faibussowitsch   PetscCheck(!(f < 0) && !(f >= prob->Nf), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %" PetscInt_FMT " must be in [0, %" PetscInt_FMT ")", f, prob->Nf);
26749371c9d4SSatish Balay   if (sol) {
26754f572ea9SToby Isaac     PetscAssertPointer(sol, 3);
26769371c9d4SSatish Balay     *sol = prob->exactSol[f];
26779371c9d4SSatish Balay   }
26789371c9d4SSatish Balay   if (ctx) {
26794f572ea9SToby Isaac     PetscAssertPointer(ctx, 4);
26809371c9d4SSatish Balay     *ctx = prob->exactCtx[f];
26819371c9d4SSatish Balay   }
26823ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2683c371a6d1SMatthew G. Knepley }
2684c371a6d1SMatthew G. Knepley 
2685c371a6d1SMatthew G. Knepley /*@C
2686578a5ef5SMatthew Knepley   PetscDSSetExactSolution - Set the pointwise exact solution function for a given test field
2687c371a6d1SMatthew G. Knepley 
268820f4b53cSBarry Smith   Not Collective
2689c371a6d1SMatthew G. Knepley 
2690c371a6d1SMatthew G. Knepley   Input Parameters:
2691dce8aebaSBarry Smith + prob - The `PetscDS`
2692c371a6d1SMatthew G. Knepley . f    - The test field number
269395cbbfd3SMatthew G. Knepley . sol  - solution function for the test fields
269420f4b53cSBarry Smith - ctx  - solution context or `NULL`
2695c371a6d1SMatthew G. Knepley 
269620f4b53cSBarry Smith   Calling sequence of `sol`:
2697dce8aebaSBarry Smith .vb
269820f4b53cSBarry Smith   PetscErrorCode sol(PetscInt dim, PetscReal t, const PetscReal x[], PetscInt Nc, PetscScalar u[], void *ctx)
2699dce8aebaSBarry Smith .ve
2700c371a6d1SMatthew G. Knepley + dim - the spatial dimension
2701c371a6d1SMatthew G. Knepley . t   - current time
2702c371a6d1SMatthew G. Knepley . x   - coordinates of the current point
2703c371a6d1SMatthew G. Knepley . Nc  - the number of field components
2704c371a6d1SMatthew G. Knepley . u   - the solution field evaluated at the current point
2705c371a6d1SMatthew G. Knepley - ctx - a user context
2706c371a6d1SMatthew G. Knepley 
2707c371a6d1SMatthew G. Knepley   Level: intermediate
2708c371a6d1SMatthew G. Knepley 
2709dce8aebaSBarry Smith .seealso: `PetscDS`, `PetscDSGetExactSolution()`
2710c371a6d1SMatthew G. Knepley @*/
2711d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDSSetExactSolution(PetscDS prob, PetscInt f, PetscErrorCode (*sol)(PetscInt dim, PetscReal t, const PetscReal x[], PetscInt Nc, PetscScalar u[], void *ctx), void *ctx)
2712d71ae5a4SJacob Faibussowitsch {
2713c371a6d1SMatthew G. Knepley   PetscFunctionBegin;
2714c371a6d1SMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
271563a3b9bcSJacob Faibussowitsch   PetscCheck(f >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %" PetscInt_FMT " must be non-negative", f);
27169566063dSJacob Faibussowitsch   PetscCall(PetscDSEnlarge_Static(prob, f + 1));
27179371c9d4SSatish Balay   if (sol) {
27189371c9d4SSatish Balay     PetscValidFunction(sol, 3);
27199371c9d4SSatish Balay     prob->exactSol[f] = sol;
27209371c9d4SSatish Balay   }
27219371c9d4SSatish Balay   if (ctx) {
27229371c9d4SSatish Balay     PetscValidFunction(ctx, 4);
27239371c9d4SSatish Balay     prob->exactCtx[f] = ctx;
27249371c9d4SSatish Balay   }
27253ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2726c371a6d1SMatthew G. Knepley }
2727c371a6d1SMatthew G. Knepley 
27285638fd0eSMatthew G. Knepley /*@C
2729f2cacb80SMatthew G. Knepley   PetscDSGetExactSolutionTimeDerivative - Get the pointwise time derivative of the exact solution function for a given test field
2730f2cacb80SMatthew G. Knepley 
273120f4b53cSBarry Smith   Not Collective
2732f2cacb80SMatthew G. Knepley 
2733f2cacb80SMatthew G. Knepley   Input Parameters:
2734dce8aebaSBarry Smith + prob - The `PetscDS`
2735f2cacb80SMatthew G. Knepley - f    - The test field number
2736f2cacb80SMatthew G. Knepley 
2737d8d19677SJose E. Roman   Output Parameters:
2738*a4e35b19SJacob Faibussowitsch + sol - time derivative of the exact solution for the test field
2739*a4e35b19SJacob Faibussowitsch - ctx - time derivative of the exact solution context
2740f2cacb80SMatthew G. Knepley 
274120f4b53cSBarry Smith   Calling sequence of `exactSol`:
2742f2cacb80SMatthew G. Knepley + dim - the spatial dimension
2743f2cacb80SMatthew G. Knepley . t   - current time
2744f2cacb80SMatthew G. Knepley . x   - coordinates of the current point
2745f2cacb80SMatthew G. Knepley . Nc  - the number of field components
2746*a4e35b19SJacob Faibussowitsch . u   - the solution field evaluated at the current point
2747f2cacb80SMatthew G. Knepley - ctx - a user context
2748f2cacb80SMatthew G. Knepley 
2749f2cacb80SMatthew G. Knepley   Level: intermediate
2750f2cacb80SMatthew G. Knepley 
2751dce8aebaSBarry Smith .seealso: `PetscDS`, `PetscDSSetExactSolutionTimeDerivative()`, `PetscDSGetExactSolution()`
2752f2cacb80SMatthew G. Knepley @*/
2753d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDSGetExactSolutionTimeDerivative(PetscDS prob, PetscInt f, PetscErrorCode (**sol)(PetscInt dim, PetscReal t, const PetscReal x[], PetscInt Nc, PetscScalar u[], void *ctx), void **ctx)
2754d71ae5a4SJacob Faibussowitsch {
2755f2cacb80SMatthew G. Knepley   PetscFunctionBegin;
2756f2cacb80SMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
275763a3b9bcSJacob Faibussowitsch   PetscCheck(!(f < 0) && !(f >= prob->Nf), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %" PetscInt_FMT " must be in [0, %" PetscInt_FMT ")", f, prob->Nf);
27589371c9d4SSatish Balay   if (sol) {
27594f572ea9SToby Isaac     PetscAssertPointer(sol, 3);
27609371c9d4SSatish Balay     *sol = prob->exactSol_t[f];
27619371c9d4SSatish Balay   }
27629371c9d4SSatish Balay   if (ctx) {
27634f572ea9SToby Isaac     PetscAssertPointer(ctx, 4);
27649371c9d4SSatish Balay     *ctx = prob->exactCtx_t[f];
27659371c9d4SSatish Balay   }
27663ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2767f2cacb80SMatthew G. Knepley }
2768f2cacb80SMatthew G. Knepley 
2769f2cacb80SMatthew G. Knepley /*@C
2770f2cacb80SMatthew G. Knepley   PetscDSSetExactSolutionTimeDerivative - Set the pointwise time derivative of the exact solution function for a given test field
2771f2cacb80SMatthew G. Knepley 
277220f4b53cSBarry Smith   Not Collective
2773f2cacb80SMatthew G. Knepley 
2774f2cacb80SMatthew G. Knepley   Input Parameters:
2775dce8aebaSBarry Smith + prob - The `PetscDS`
2776f2cacb80SMatthew G. Knepley . f    - The test field number
2777f2cacb80SMatthew G. Knepley . sol  - time derivative of the solution function for the test fields
277820f4b53cSBarry Smith - ctx  - time derivative of the solution context or `NULL`
2779f2cacb80SMatthew G. Knepley 
278020f4b53cSBarry Smith   Calling sequence of `sol`:
2781dce8aebaSBarry Smith .vb
278220f4b53cSBarry Smith   PetscErrorCode sol(PetscInt dim, PetscReal t, const PetscReal x[], PetscInt Nc, PetscScalar u[], void *ctx)
2783dce8aebaSBarry Smith .ve
2784f2cacb80SMatthew G. Knepley + dim - the spatial dimension
2785f2cacb80SMatthew G. Knepley . t   - current time
2786f2cacb80SMatthew G. Knepley . x   - coordinates of the current point
2787f2cacb80SMatthew G. Knepley . Nc  - the number of field components
2788f2cacb80SMatthew G. Knepley . u   - the solution field evaluated at the current point
2789f2cacb80SMatthew G. Knepley - ctx - a user context
2790f2cacb80SMatthew G. Knepley 
2791f2cacb80SMatthew G. Knepley   Level: intermediate
2792f2cacb80SMatthew G. Knepley 
2793dce8aebaSBarry Smith .seealso: `PetscDS`, `PetscDSGetExactSolutionTimeDerivative()`, `PetscDSSetExactSolution()`
2794f2cacb80SMatthew G. Knepley @*/
2795d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDSSetExactSolutionTimeDerivative(PetscDS prob, PetscInt f, PetscErrorCode (*sol)(PetscInt dim, PetscReal t, const PetscReal x[], PetscInt Nc, PetscScalar u[], void *ctx), void *ctx)
2796d71ae5a4SJacob Faibussowitsch {
2797f2cacb80SMatthew G. Knepley   PetscFunctionBegin;
2798f2cacb80SMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
279963a3b9bcSJacob Faibussowitsch   PetscCheck(f >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %" PetscInt_FMT " must be non-negative", f);
28009566063dSJacob Faibussowitsch   PetscCall(PetscDSEnlarge_Static(prob, f + 1));
28019371c9d4SSatish Balay   if (sol) {
28029371c9d4SSatish Balay     PetscValidFunction(sol, 3);
28039371c9d4SSatish Balay     prob->exactSol_t[f] = sol;
28049371c9d4SSatish Balay   }
28059371c9d4SSatish Balay   if (ctx) {
28069371c9d4SSatish Balay     PetscValidFunction(ctx, 4);
28079371c9d4SSatish Balay     prob->exactCtx_t[f] = ctx;
28089371c9d4SSatish Balay   }
28093ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2810f2cacb80SMatthew G. Knepley }
2811f2cacb80SMatthew G. Knepley 
2812f2cacb80SMatthew G. Knepley /*@C
281397b6e6e8SMatthew G. Knepley   PetscDSGetConstants - Returns the array of constants passed to point functions
281497b6e6e8SMatthew G. Knepley 
281520f4b53cSBarry Smith   Not Collective
281697b6e6e8SMatthew G. Knepley 
281797b6e6e8SMatthew G. Knepley   Input Parameter:
2818dce8aebaSBarry Smith . prob - The `PetscDS` object
281997b6e6e8SMatthew G. Knepley 
282097b6e6e8SMatthew G. Knepley   Output Parameters:
282197b6e6e8SMatthew G. Knepley + numConstants - The number of constants
282297b6e6e8SMatthew G. Knepley - constants    - The array of constants, NULL if there are none
282397b6e6e8SMatthew G. Knepley 
282497b6e6e8SMatthew G. Knepley   Level: intermediate
282597b6e6e8SMatthew G. Knepley 
2826dce8aebaSBarry Smith .seealso: `PetscDS`, `PetscDSSetConstants()`, `PetscDSCreate()`
282797b6e6e8SMatthew G. Knepley @*/
2828d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDSGetConstants(PetscDS prob, PetscInt *numConstants, const PetscScalar *constants[])
2829d71ae5a4SJacob Faibussowitsch {
283097b6e6e8SMatthew G. Knepley   PetscFunctionBegin;
283197b6e6e8SMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
28329371c9d4SSatish Balay   if (numConstants) {
28334f572ea9SToby Isaac     PetscAssertPointer(numConstants, 2);
28349371c9d4SSatish Balay     *numConstants = prob->numConstants;
28359371c9d4SSatish Balay   }
28369371c9d4SSatish Balay   if (constants) {
28374f572ea9SToby Isaac     PetscAssertPointer(constants, 3);
28389371c9d4SSatish Balay     *constants = prob->constants;
28399371c9d4SSatish Balay   }
28403ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
284197b6e6e8SMatthew G. Knepley }
284297b6e6e8SMatthew G. Knepley 
28430d3e9b51SMatthew G. Knepley /*@C
284497b6e6e8SMatthew G. Knepley   PetscDSSetConstants - Set the array of constants passed to point functions
284597b6e6e8SMatthew G. Knepley 
284620f4b53cSBarry Smith   Not Collective
284797b6e6e8SMatthew G. Knepley 
284897b6e6e8SMatthew G. Knepley   Input Parameters:
2849dce8aebaSBarry Smith + prob         - The `PetscDS` object
285097b6e6e8SMatthew G. Knepley . numConstants - The number of constants
285197b6e6e8SMatthew G. Knepley - constants    - The array of constants, NULL if there are none
285297b6e6e8SMatthew G. Knepley 
285397b6e6e8SMatthew G. Knepley   Level: intermediate
285497b6e6e8SMatthew G. Knepley 
2855dce8aebaSBarry Smith .seealso: `PetscDS`, `PetscDSGetConstants()`, `PetscDSCreate()`
285697b6e6e8SMatthew G. Knepley @*/
2857d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDSSetConstants(PetscDS prob, PetscInt numConstants, PetscScalar constants[])
2858d71ae5a4SJacob Faibussowitsch {
285997b6e6e8SMatthew G. Knepley   PetscFunctionBegin;
286097b6e6e8SMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
286197b6e6e8SMatthew G. Knepley   if (numConstants != prob->numConstants) {
28629566063dSJacob Faibussowitsch     PetscCall(PetscFree(prob->constants));
286397b6e6e8SMatthew G. Knepley     prob->numConstants = numConstants;
286497b6e6e8SMatthew G. Knepley     if (prob->numConstants) {
28659566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(prob->numConstants, &prob->constants));
286620be0f5bSMatthew G. Knepley     } else {
286720be0f5bSMatthew G. Knepley       prob->constants = NULL;
286820be0f5bSMatthew G. Knepley     }
286920be0f5bSMatthew G. Knepley   }
287020be0f5bSMatthew G. Knepley   if (prob->numConstants) {
28714f572ea9SToby Isaac     PetscAssertPointer(constants, 3);
28729566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(prob->constants, constants, prob->numConstants));
287397b6e6e8SMatthew G. Knepley   }
28743ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
287597b6e6e8SMatthew G. Knepley }
287697b6e6e8SMatthew G. Knepley 
28774cd1e086SMatthew G. Knepley /*@
28784cd1e086SMatthew G. Knepley   PetscDSGetFieldIndex - Returns the index of the given field
28794cd1e086SMatthew G. Knepley 
288020f4b53cSBarry Smith   Not Collective
28814cd1e086SMatthew G. Knepley 
28824cd1e086SMatthew G. Knepley   Input Parameters:
2883dce8aebaSBarry Smith + prob - The `PetscDS` object
28844cd1e086SMatthew G. Knepley - disc - The discretization object
28854cd1e086SMatthew G. Knepley 
28864cd1e086SMatthew G. Knepley   Output Parameter:
28874cd1e086SMatthew G. Knepley . f - The field number
28884cd1e086SMatthew G. Knepley 
28894cd1e086SMatthew G. Knepley   Level: beginner
28904cd1e086SMatthew G. Knepley 
2891dce8aebaSBarry Smith .seealso: `PetscDS`, `PetscGetDiscretization()`, `PetscDSGetNumFields()`, `PetscDSCreate()`
28924cd1e086SMatthew G. Knepley @*/
2893d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDSGetFieldIndex(PetscDS prob, PetscObject disc, PetscInt *f)
2894d71ae5a4SJacob Faibussowitsch {
28954cd1e086SMatthew G. Knepley   PetscInt g;
28964cd1e086SMatthew G. Knepley 
28974cd1e086SMatthew G. Knepley   PetscFunctionBegin;
28984cd1e086SMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
28994f572ea9SToby Isaac   PetscAssertPointer(f, 3);
29004cd1e086SMatthew G. Knepley   *f = -1;
29019371c9d4SSatish Balay   for (g = 0; g < prob->Nf; ++g) {
29029371c9d4SSatish Balay     if (disc == prob->disc[g]) break;
29039371c9d4SSatish Balay   }
290408401ef6SPierre Jolivet   PetscCheck(g != prob->Nf, PetscObjectComm((PetscObject)prob), PETSC_ERR_ARG_WRONG, "Field not found in PetscDS.");
29054cd1e086SMatthew G. Knepley   *f = g;
29063ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
29074cd1e086SMatthew G. Knepley }
29084cd1e086SMatthew G. Knepley 
29094cd1e086SMatthew G. Knepley /*@
29104cd1e086SMatthew G. Knepley   PetscDSGetFieldSize - Returns the size of the given field in the full space basis
29114cd1e086SMatthew G. Knepley 
291220f4b53cSBarry Smith   Not Collective
29134cd1e086SMatthew G. Knepley 
29144cd1e086SMatthew G. Knepley   Input Parameters:
2915dce8aebaSBarry Smith + prob - The `PetscDS` object
29164cd1e086SMatthew G. Knepley - f    - The field number
29174cd1e086SMatthew G. Knepley 
29184cd1e086SMatthew G. Knepley   Output Parameter:
29194cd1e086SMatthew G. Knepley . size - The size
29204cd1e086SMatthew G. Knepley 
29214cd1e086SMatthew G. Knepley   Level: beginner
29224cd1e086SMatthew G. Knepley 
2923dce8aebaSBarry Smith .seealso: `PetscDS`, `PetscDSGetFieldOffset()`, `PetscDSGetNumFields()`, `PetscDSCreate()`
29244cd1e086SMatthew G. Knepley @*/
2925d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDSGetFieldSize(PetscDS prob, PetscInt f, PetscInt *size)
2926d71ae5a4SJacob Faibussowitsch {
29274cd1e086SMatthew G. Knepley   PetscFunctionBegin;
29284cd1e086SMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
29294f572ea9SToby Isaac   PetscAssertPointer(size, 3);
293063a3b9bcSJacob Faibussowitsch   PetscCheck(!(f < 0) && !(f >= prob->Nf), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %" PetscInt_FMT " must be in [0, %" PetscInt_FMT ")", f, prob->Nf);
29319566063dSJacob Faibussowitsch   PetscCall(PetscDSSetUp(prob));
2932d4742ddaSMatthew G. Knepley   *size = prob->Nb[f];
29333ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
29344cd1e086SMatthew G. Knepley }
29354cd1e086SMatthew G. Knepley 
2936bc4ae4beSMatthew G. Knepley /*@
2937bc4ae4beSMatthew G. Knepley   PetscDSGetFieldOffset - Returns the offset of the given field in the full space basis
2938bc4ae4beSMatthew G. Knepley 
293920f4b53cSBarry Smith   Not Collective
2940bc4ae4beSMatthew G. Knepley 
2941bc4ae4beSMatthew G. Knepley   Input Parameters:
2942dce8aebaSBarry Smith + prob - The `PetscDS` object
2943bc4ae4beSMatthew G. Knepley - f    - The field number
2944bc4ae4beSMatthew G. Knepley 
2945bc4ae4beSMatthew G. Knepley   Output Parameter:
2946bc4ae4beSMatthew G. Knepley . off - The offset
2947bc4ae4beSMatthew G. Knepley 
2948bc4ae4beSMatthew G. Knepley   Level: beginner
2949bc4ae4beSMatthew G. Knepley 
2950dce8aebaSBarry Smith .seealso: `PetscDS`, `PetscDSGetFieldSize()`, `PetscDSGetNumFields()`, `PetscDSCreate()`
2951bc4ae4beSMatthew G. Knepley @*/
2952d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDSGetFieldOffset(PetscDS prob, PetscInt f, PetscInt *off)
2953d71ae5a4SJacob Faibussowitsch {
29544cd1e086SMatthew G. Knepley   PetscInt size, g;
29552764a2aaSMatthew G. Knepley 
29562764a2aaSMatthew G. Knepley   PetscFunctionBegin;
29572764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
29584f572ea9SToby Isaac   PetscAssertPointer(off, 3);
295963a3b9bcSJacob Faibussowitsch   PetscCheck(!(f < 0) && !(f >= prob->Nf), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %" PetscInt_FMT " must be in [0, %" PetscInt_FMT ")", f, prob->Nf);
29602764a2aaSMatthew G. Knepley   *off = 0;
29612764a2aaSMatthew G. Knepley   for (g = 0; g < f; ++g) {
29629566063dSJacob Faibussowitsch     PetscCall(PetscDSGetFieldSize(prob, g, &size));
29634cd1e086SMatthew G. Knepley     *off += size;
29642764a2aaSMatthew G. Knepley   }
29653ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
29662764a2aaSMatthew G. Knepley }
29672764a2aaSMatthew G. Knepley 
2968bc4ae4beSMatthew G. Knepley /*@
29695fedec97SMatthew G. Knepley   PetscDSGetFieldOffsetCohesive - Returns the offset of the given field in the full space basis on a cohesive cell
29705fedec97SMatthew G. Knepley 
297120f4b53cSBarry Smith   Not Collective
29725fedec97SMatthew G. Knepley 
29735fedec97SMatthew G. Knepley   Input Parameters:
297460225df5SJacob Faibussowitsch + ds - The `PetscDS` object
29755fedec97SMatthew G. Knepley - f  - The field number
29765fedec97SMatthew G. Knepley 
29775fedec97SMatthew G. Knepley   Output Parameter:
29785fedec97SMatthew G. Knepley . off - The offset
29795fedec97SMatthew G. Knepley 
29805fedec97SMatthew G. Knepley   Level: beginner
29815fedec97SMatthew G. Knepley 
2982dce8aebaSBarry Smith .seealso: `PetscDS`, `PetscDSGetFieldSize()`, `PetscDSGetNumFields()`, `PetscDSCreate()`
29835fedec97SMatthew G. Knepley @*/
2984d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDSGetFieldOffsetCohesive(PetscDS ds, PetscInt f, PetscInt *off)
2985d71ae5a4SJacob Faibussowitsch {
29865fedec97SMatthew G. Knepley   PetscInt size, g;
29875fedec97SMatthew G. Knepley 
29885fedec97SMatthew G. Knepley   PetscFunctionBegin;
29895fedec97SMatthew G. Knepley   PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1);
29904f572ea9SToby Isaac   PetscAssertPointer(off, 3);
299163a3b9bcSJacob Faibussowitsch   PetscCheck(!(f < 0) && !(f >= ds->Nf), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %" PetscInt_FMT " must be in [0, %" PetscInt_FMT ")", f, ds->Nf);
29925fedec97SMatthew G. Knepley   *off = 0;
29935fedec97SMatthew G. Knepley   for (g = 0; g < f; ++g) {
29945fedec97SMatthew G. Knepley     PetscBool cohesive;
29955fedec97SMatthew G. Knepley 
29969566063dSJacob Faibussowitsch     PetscCall(PetscDSGetCohesive(ds, g, &cohesive));
29979566063dSJacob Faibussowitsch     PetscCall(PetscDSGetFieldSize(ds, g, &size));
29985fedec97SMatthew G. Knepley     *off += cohesive ? size : size * 2;
29995fedec97SMatthew G. Knepley   }
30003ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
30015fedec97SMatthew G. Knepley }
30025fedec97SMatthew G. Knepley 
30035fedec97SMatthew G. Knepley /*@
300447e57110SSander Arens   PetscDSGetDimensions - Returns the size of the approximation space for each field on an evaluation point
3005bc4ae4beSMatthew G. Knepley 
300620f4b53cSBarry Smith   Not Collective
3007bc4ae4beSMatthew G. Knepley 
300847e57110SSander Arens   Input Parameter:
3009dce8aebaSBarry Smith . prob - The `PetscDS` object
3010bc4ae4beSMatthew G. Knepley 
3011bc4ae4beSMatthew G. Knepley   Output Parameter:
301247e57110SSander Arens . dimensions - The number of dimensions
3013bc4ae4beSMatthew G. Knepley 
3014bc4ae4beSMatthew G. Knepley   Level: beginner
3015bc4ae4beSMatthew G. Knepley 
3016dce8aebaSBarry Smith .seealso: `PetscDS`, `PetscDSGetComponentOffsets()`, `PetscDSGetNumFields()`, `PetscDSCreate()`
3017bc4ae4beSMatthew G. Knepley @*/
3018d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDSGetDimensions(PetscDS prob, PetscInt *dimensions[])
3019d71ae5a4SJacob Faibussowitsch {
30202764a2aaSMatthew G. Knepley   PetscFunctionBegin;
30212764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
30229566063dSJacob Faibussowitsch   PetscCall(PetscDSSetUp(prob));
30234f572ea9SToby Isaac   PetscAssertPointer(dimensions, 2);
302447e57110SSander Arens   *dimensions = prob->Nb;
30253ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
30266ce16762SMatthew G. Knepley }
302747e57110SSander Arens 
302847e57110SSander Arens /*@
302947e57110SSander Arens   PetscDSGetComponents - Returns the number of components for each field on an evaluation point
303047e57110SSander Arens 
303120f4b53cSBarry Smith   Not Collective
303247e57110SSander Arens 
303347e57110SSander Arens   Input Parameter:
3034dce8aebaSBarry Smith . prob - The `PetscDS` object
303547e57110SSander Arens 
303647e57110SSander Arens   Output Parameter:
303747e57110SSander Arens . components - The number of components
303847e57110SSander Arens 
303947e57110SSander Arens   Level: beginner
304047e57110SSander Arens 
3041dce8aebaSBarry Smith .seealso: `PetscDS`, `PetscDSGetComponentOffsets()`, `PetscDSGetNumFields()`, `PetscDSCreate()`
304247e57110SSander Arens @*/
3043d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDSGetComponents(PetscDS prob, PetscInt *components[])
3044d71ae5a4SJacob Faibussowitsch {
304547e57110SSander Arens   PetscFunctionBegin;
304647e57110SSander Arens   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
30479566063dSJacob Faibussowitsch   PetscCall(PetscDSSetUp(prob));
30484f572ea9SToby Isaac   PetscAssertPointer(components, 2);
304947e57110SSander Arens   *components = prob->Nc;
30503ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
30516ce16762SMatthew G. Knepley }
30526ce16762SMatthew G. Knepley 
30536ce16762SMatthew G. Knepley /*@
30546ce16762SMatthew G. Knepley   PetscDSGetComponentOffset - Returns the offset of the given field on an evaluation point
30556ce16762SMatthew G. Knepley 
305620f4b53cSBarry Smith   Not Collective
30576ce16762SMatthew G. Knepley 
30586ce16762SMatthew G. Knepley   Input Parameters:
3059dce8aebaSBarry Smith + prob - The `PetscDS` object
30606ce16762SMatthew G. Knepley - f    - The field number
30616ce16762SMatthew G. Knepley 
30626ce16762SMatthew G. Knepley   Output Parameter:
30636ce16762SMatthew G. Knepley . off - The offset
30646ce16762SMatthew G. Knepley 
30656ce16762SMatthew G. Knepley   Level: beginner
30666ce16762SMatthew G. Knepley 
3067dce8aebaSBarry Smith .seealso: `PetscDS`, `PetscDSGetNumFields()`, `PetscDSCreate()`
30686ce16762SMatthew G. Knepley @*/
3069d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDSGetComponentOffset(PetscDS prob, PetscInt f, PetscInt *off)
3070d71ae5a4SJacob Faibussowitsch {
30716ce16762SMatthew G. Knepley   PetscFunctionBegin;
30726ce16762SMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
30734f572ea9SToby Isaac   PetscAssertPointer(off, 3);
307463a3b9bcSJacob Faibussowitsch   PetscCheck(!(f < 0) && !(f >= prob->Nf), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %" PetscInt_FMT " must be in [0, %" PetscInt_FMT ")", f, prob->Nf);
30759566063dSJacob Faibussowitsch   PetscCall(PetscDSSetUp(prob));
307647e57110SSander Arens   *off = prob->off[f];
30773ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
30782764a2aaSMatthew G. Knepley }
30792764a2aaSMatthew G. Knepley 
3080194d53e6SMatthew G. Knepley /*@
3081194d53e6SMatthew G. Knepley   PetscDSGetComponentOffsets - Returns the offset of each field on an evaluation point
3082194d53e6SMatthew G. Knepley 
308320f4b53cSBarry Smith   Not Collective
3084194d53e6SMatthew G. Knepley 
3085194d53e6SMatthew G. Knepley   Input Parameter:
3086dce8aebaSBarry Smith . prob - The `PetscDS` object
3087194d53e6SMatthew G. Knepley 
3088194d53e6SMatthew G. Knepley   Output Parameter:
3089194d53e6SMatthew G. Knepley . offsets - The offsets
3090194d53e6SMatthew G. Knepley 
3091194d53e6SMatthew G. Knepley   Level: beginner
3092194d53e6SMatthew G. Knepley 
3093dce8aebaSBarry Smith .seealso: `PetscDS`, `PetscDSGetNumFields()`, `PetscDSCreate()`
3094194d53e6SMatthew G. Knepley @*/
3095d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDSGetComponentOffsets(PetscDS prob, PetscInt *offsets[])
3096d71ae5a4SJacob Faibussowitsch {
3097194d53e6SMatthew G. Knepley   PetscFunctionBegin;
3098194d53e6SMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
30994f572ea9SToby Isaac   PetscAssertPointer(offsets, 2);
31009566063dSJacob Faibussowitsch   PetscCall(PetscDSSetUp(prob));
3101194d53e6SMatthew G. Knepley   *offsets = prob->off;
31023ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3103194d53e6SMatthew G. Knepley }
3104194d53e6SMatthew G. Knepley 
3105194d53e6SMatthew G. Knepley /*@
3106194d53e6SMatthew G. Knepley   PetscDSGetComponentDerivativeOffsets - Returns the offset of each field derivative on an evaluation point
3107194d53e6SMatthew G. Knepley 
310820f4b53cSBarry Smith   Not Collective
3109194d53e6SMatthew G. Knepley 
3110194d53e6SMatthew G. Knepley   Input Parameter:
3111dce8aebaSBarry Smith . prob - The `PetscDS` object
3112194d53e6SMatthew G. Knepley 
3113194d53e6SMatthew G. Knepley   Output Parameter:
3114194d53e6SMatthew G. Knepley . offsets - The offsets
3115194d53e6SMatthew G. Knepley 
3116194d53e6SMatthew G. Knepley   Level: beginner
3117194d53e6SMatthew G. Knepley 
3118dce8aebaSBarry Smith .seealso: `PetscDS`, `PetscDSGetNumFields()`, `PetscDSCreate()`
3119194d53e6SMatthew G. Knepley @*/
3120d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDSGetComponentDerivativeOffsets(PetscDS prob, PetscInt *offsets[])
3121d71ae5a4SJacob Faibussowitsch {
3122194d53e6SMatthew G. Knepley   PetscFunctionBegin;
3123194d53e6SMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
31244f572ea9SToby Isaac   PetscAssertPointer(offsets, 2);
31259566063dSJacob Faibussowitsch   PetscCall(PetscDSSetUp(prob));
3126194d53e6SMatthew G. Knepley   *offsets = prob->offDer;
31273ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3128194d53e6SMatthew G. Knepley }
3129194d53e6SMatthew G. Knepley 
31309ee2af8cSMatthew G. Knepley /*@
31319ee2af8cSMatthew G. Knepley   PetscDSGetComponentOffsetsCohesive - Returns the offset of each field on an evaluation point
31329ee2af8cSMatthew G. Knepley 
313320f4b53cSBarry Smith   Not Collective
31349ee2af8cSMatthew G. Knepley 
31359ee2af8cSMatthew G. Knepley   Input Parameters:
3136dce8aebaSBarry Smith + ds - The `PetscDS` object
31379ee2af8cSMatthew G. Knepley - s  - The cohesive side, 0 for negative, 1 for positive, 2 for cohesive
31389ee2af8cSMatthew G. Knepley 
31399ee2af8cSMatthew G. Knepley   Output Parameter:
31409ee2af8cSMatthew G. Knepley . offsets - The offsets
31419ee2af8cSMatthew G. Knepley 
31429ee2af8cSMatthew G. Knepley   Level: beginner
31439ee2af8cSMatthew G. Knepley 
3144dce8aebaSBarry Smith .seealso: `PetscDS`, `PetscDSGetNumFields()`, `PetscDSCreate()`
31459ee2af8cSMatthew G. Knepley @*/
3146d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDSGetComponentOffsetsCohesive(PetscDS ds, PetscInt s, PetscInt *offsets[])
3147d71ae5a4SJacob Faibussowitsch {
31489ee2af8cSMatthew G. Knepley   PetscFunctionBegin;
31499ee2af8cSMatthew G. Knepley   PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1);
31504f572ea9SToby Isaac   PetscAssertPointer(offsets, 3);
315128b400f6SJacob Faibussowitsch   PetscCheck(ds->isCohesive, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Cohesive offsets are only valid for a cohesive DS");
315263a3b9bcSJacob Faibussowitsch   PetscCheck(!(s < 0) && !(s > 2), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cohesive side %" PetscInt_FMT " is not in [0, 2]", s);
31539566063dSJacob Faibussowitsch   PetscCall(PetscDSSetUp(ds));
31549ee2af8cSMatthew G. Knepley   *offsets = ds->offCohesive[s];
31553ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
31569ee2af8cSMatthew G. Knepley }
31579ee2af8cSMatthew G. Knepley 
31589ee2af8cSMatthew G. Knepley /*@
31599ee2af8cSMatthew G. Knepley   PetscDSGetComponentDerivativeOffsetsCohesive - Returns the offset of each field derivative on an evaluation point
31609ee2af8cSMatthew G. Knepley 
316120f4b53cSBarry Smith   Not Collective
31629ee2af8cSMatthew G. Knepley 
31639ee2af8cSMatthew G. Knepley   Input Parameters:
3164dce8aebaSBarry Smith + ds - The `PetscDS` object
31659ee2af8cSMatthew G. Knepley - s  - The cohesive side, 0 for negative, 1 for positive, 2 for cohesive
31669ee2af8cSMatthew G. Knepley 
31679ee2af8cSMatthew G. Knepley   Output Parameter:
31689ee2af8cSMatthew G. Knepley . offsets - The offsets
31699ee2af8cSMatthew G. Knepley 
31709ee2af8cSMatthew G. Knepley   Level: beginner
31719ee2af8cSMatthew G. Knepley 
3172dce8aebaSBarry Smith .seealso: `PetscDS`, `PetscDSGetNumFields()`, `PetscDSCreate()`
31739ee2af8cSMatthew G. Knepley @*/
3174d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDSGetComponentDerivativeOffsetsCohesive(PetscDS ds, PetscInt s, PetscInt *offsets[])
3175d71ae5a4SJacob Faibussowitsch {
31769ee2af8cSMatthew G. Knepley   PetscFunctionBegin;
31779ee2af8cSMatthew G. Knepley   PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1);
31784f572ea9SToby Isaac   PetscAssertPointer(offsets, 3);
317928b400f6SJacob Faibussowitsch   PetscCheck(ds->isCohesive, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Cohesive offsets are only valid for a cohesive DS");
318063a3b9bcSJacob Faibussowitsch   PetscCheck(!(s < 0) && !(s > 2), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cohesive side %" PetscInt_FMT " is not in [0, 2]", s);
31819566063dSJacob Faibussowitsch   PetscCall(PetscDSSetUp(ds));
31829ee2af8cSMatthew G. Knepley   *offsets = ds->offDerCohesive[s];
31833ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
31849ee2af8cSMatthew G. Knepley }
31859ee2af8cSMatthew G. Knepley 
318668c9edb9SMatthew G. Knepley /*@C
318768c9edb9SMatthew G. Knepley   PetscDSGetTabulation - Return the basis tabulation at quadrature points for the volume discretization
318868c9edb9SMatthew G. Knepley 
318920f4b53cSBarry Smith   Not Collective
319068c9edb9SMatthew G. Knepley 
319168c9edb9SMatthew G. Knepley   Input Parameter:
3192dce8aebaSBarry Smith . prob - The `PetscDS` object
319368c9edb9SMatthew G. Knepley 
3194ef0bb6c7SMatthew G. Knepley   Output Parameter:
3195ef0bb6c7SMatthew G. Knepley . T - The basis function and derivatives tabulation at quadrature points for each field
319668c9edb9SMatthew G. Knepley 
319768c9edb9SMatthew G. Knepley   Level: intermediate
319868c9edb9SMatthew G. Knepley 
3199dce8aebaSBarry Smith .seealso: `PetscDS`, `PetscTabulation`, `PetscDSCreate()`
320068c9edb9SMatthew G. Knepley @*/
3201d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDSGetTabulation(PetscDS prob, PetscTabulation *T[])
3202d71ae5a4SJacob Faibussowitsch {
32032764a2aaSMatthew G. Knepley   PetscFunctionBegin;
32042764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
32054f572ea9SToby Isaac   PetscAssertPointer(T, 2);
32069566063dSJacob Faibussowitsch   PetscCall(PetscDSSetUp(prob));
3207ef0bb6c7SMatthew G. Knepley   *T = prob->T;
32083ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
32092764a2aaSMatthew G. Knepley }
32102764a2aaSMatthew G. Knepley 
321168c9edb9SMatthew G. Knepley /*@C
32124d0b9603SSander Arens   PetscDSGetFaceTabulation - Return the basis tabulation at quadrature points on the faces
321368c9edb9SMatthew G. Knepley 
321420f4b53cSBarry Smith   Not Collective
321568c9edb9SMatthew G. Knepley 
321668c9edb9SMatthew G. Knepley   Input Parameter:
3217dce8aebaSBarry Smith . prob - The `PetscDS` object
321868c9edb9SMatthew G. Knepley 
3219ef0bb6c7SMatthew G. Knepley   Output Parameter:
3220a5b23f4aSJose E. Roman . Tf - The basis function and derivative tabulation on each local face at quadrature points for each and field
322168c9edb9SMatthew G. Knepley 
322268c9edb9SMatthew G. Knepley   Level: intermediate
322368c9edb9SMatthew G. Knepley 
3224dce8aebaSBarry Smith .seealso: `PetscTabulation`, `PetscDS`, `PetscDSGetTabulation()`, `PetscDSCreate()`
322568c9edb9SMatthew G. Knepley @*/
3226d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDSGetFaceTabulation(PetscDS prob, PetscTabulation *Tf[])
3227d71ae5a4SJacob Faibussowitsch {
32282764a2aaSMatthew G. Knepley   PetscFunctionBegin;
32292764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
32304f572ea9SToby Isaac   PetscAssertPointer(Tf, 2);
32319566063dSJacob Faibussowitsch   PetscCall(PetscDSSetUp(prob));
3232ef0bb6c7SMatthew G. Knepley   *Tf = prob->Tf;
32333ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
32342764a2aaSMatthew G. Knepley }
32352764a2aaSMatthew G. Knepley 
3236d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDSGetEvaluationArrays(PetscDS prob, PetscScalar **u, PetscScalar **u_t, PetscScalar **u_x)
3237d71ae5a4SJacob Faibussowitsch {
32382764a2aaSMatthew G. Knepley   PetscFunctionBegin;
32392764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
32409566063dSJacob Faibussowitsch   PetscCall(PetscDSSetUp(prob));
32419371c9d4SSatish Balay   if (u) {
32424f572ea9SToby Isaac     PetscAssertPointer(u, 2);
32439371c9d4SSatish Balay     *u = prob->u;
32449371c9d4SSatish Balay   }
32459371c9d4SSatish Balay   if (u_t) {
32464f572ea9SToby Isaac     PetscAssertPointer(u_t, 3);
32479371c9d4SSatish Balay     *u_t = prob->u_t;
32489371c9d4SSatish Balay   }
32499371c9d4SSatish Balay   if (u_x) {
32504f572ea9SToby Isaac     PetscAssertPointer(u_x, 4);
32519371c9d4SSatish Balay     *u_x = prob->u_x;
32529371c9d4SSatish Balay   }
32533ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
32542764a2aaSMatthew G. Knepley }
32552764a2aaSMatthew G. Knepley 
3256d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDSGetWeakFormArrays(PetscDS prob, PetscScalar **f0, PetscScalar **f1, PetscScalar **g0, PetscScalar **g1, PetscScalar **g2, PetscScalar **g3)
3257d71ae5a4SJacob Faibussowitsch {
32582764a2aaSMatthew G. Knepley   PetscFunctionBegin;
32592764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
32609566063dSJacob Faibussowitsch   PetscCall(PetscDSSetUp(prob));
32619371c9d4SSatish Balay   if (f0) {
32624f572ea9SToby Isaac     PetscAssertPointer(f0, 2);
32639371c9d4SSatish Balay     *f0 = prob->f0;
32649371c9d4SSatish Balay   }
32659371c9d4SSatish Balay   if (f1) {
32664f572ea9SToby Isaac     PetscAssertPointer(f1, 3);
32679371c9d4SSatish Balay     *f1 = prob->f1;
32689371c9d4SSatish Balay   }
32699371c9d4SSatish Balay   if (g0) {
32704f572ea9SToby Isaac     PetscAssertPointer(g0, 4);
32719371c9d4SSatish Balay     *g0 = prob->g0;
32729371c9d4SSatish Balay   }
32739371c9d4SSatish Balay   if (g1) {
32744f572ea9SToby Isaac     PetscAssertPointer(g1, 5);
32759371c9d4SSatish Balay     *g1 = prob->g1;
32769371c9d4SSatish Balay   }
32779371c9d4SSatish Balay   if (g2) {
32784f572ea9SToby Isaac     PetscAssertPointer(g2, 6);
32799371c9d4SSatish Balay     *g2 = prob->g2;
32809371c9d4SSatish Balay   }
32819371c9d4SSatish Balay   if (g3) {
32824f572ea9SToby Isaac     PetscAssertPointer(g3, 7);
32839371c9d4SSatish Balay     *g3 = prob->g3;
32849371c9d4SSatish Balay   }
32853ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
32862764a2aaSMatthew G. Knepley }
32872764a2aaSMatthew G. Knepley 
3288d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDSGetWorkspace(PetscDS prob, PetscReal **x, PetscScalar **basisReal, PetscScalar **basisDerReal, PetscScalar **testReal, PetscScalar **testDerReal)
3289d71ae5a4SJacob Faibussowitsch {
32902764a2aaSMatthew G. Knepley   PetscFunctionBegin;
32912764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
32929566063dSJacob Faibussowitsch   PetscCall(PetscDSSetUp(prob));
32939371c9d4SSatish Balay   if (x) {
32944f572ea9SToby Isaac     PetscAssertPointer(x, 2);
32959371c9d4SSatish Balay     *x = prob->x;
32969371c9d4SSatish Balay   }
32979371c9d4SSatish Balay   if (basisReal) {
32984f572ea9SToby Isaac     PetscAssertPointer(basisReal, 3);
32999371c9d4SSatish Balay     *basisReal = prob->basisReal;
33009371c9d4SSatish Balay   }
33019371c9d4SSatish Balay   if (basisDerReal) {
33024f572ea9SToby Isaac     PetscAssertPointer(basisDerReal, 4);
33039371c9d4SSatish Balay     *basisDerReal = prob->basisDerReal;
33049371c9d4SSatish Balay   }
33059371c9d4SSatish Balay   if (testReal) {
33064f572ea9SToby Isaac     PetscAssertPointer(testReal, 5);
33079371c9d4SSatish Balay     *testReal = prob->testReal;
33089371c9d4SSatish Balay   }
33099371c9d4SSatish Balay   if (testDerReal) {
33104f572ea9SToby Isaac     PetscAssertPointer(testDerReal, 6);
33119371c9d4SSatish Balay     *testDerReal = prob->testDerReal;
33129371c9d4SSatish Balay   }
33133ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
33142764a2aaSMatthew G. Knepley }
33152764a2aaSMatthew G. Knepley 
331658ebd649SToby Isaac /*@C
3317*a4e35b19SJacob Faibussowitsch   PetscDSAddBoundary - Add a boundary condition to the model.
331858ebd649SToby Isaac 
331920f4b53cSBarry Smith   Collective
3320783e2ec8SMatthew G. Knepley 
332158ebd649SToby Isaac   Input Parameters:
332258ebd649SToby Isaac + ds       - The PetscDS object
3323dce8aebaSBarry Smith . type     - The type of condition, e.g. `DM_BC_ESSENTIAL`/`DM_BC_ESSENTIAL_FIELD` (Dirichlet), or `DM_BC_NATURAL` (Neumann)
332458ebd649SToby Isaac . name     - The BC name
332545480ffeSMatthew G. Knepley . label    - The label defining constrained points
3326dce8aebaSBarry Smith . Nv       - The number of `DMLabel` values for constrained points
332745480ffeSMatthew G. Knepley . values   - An array of label values for constrained points
332858ebd649SToby Isaac . field    - The field to constrain
332945480ffeSMatthew G. Knepley . Nc       - The number of constrained field components (0 will constrain all fields)
333058ebd649SToby Isaac . comps    - An array of constrained component numbers
333158ebd649SToby Isaac . bcFunc   - A pointwise function giving boundary values
3332a5b23f4aSJose E. Roman . bcFunc_t - A pointwise function giving the time derivative of the boundary values, or NULL
333358ebd649SToby Isaac - ctx      - An optional user context for bcFunc
333458ebd649SToby Isaac 
33352fe279fdSBarry Smith   Output Parameter:
333660225df5SJacob Faibussowitsch . bd - The boundary number
333745480ffeSMatthew G. Knepley 
333858ebd649SToby Isaac   Options Database Keys:
333958ebd649SToby Isaac + -bc_<boundary name> <num>      - Overrides the boundary ids
334058ebd649SToby Isaac - -bc_<boundary name>_comp <num> - Overrides the boundary components
334158ebd649SToby Isaac 
3342dce8aebaSBarry Smith   Level: developer
3343dce8aebaSBarry Smith 
334456cf3b9cSMatthew G. Knepley   Note:
3345*a4e35b19SJacob Faibussowitsch   Both `bcFunc` and `bcFunc_t` will depend on the boundary condition type. If the type if `DM_BC_ESSENTIAL`, then the calling sequence is\:
334656cf3b9cSMatthew G. Knepley 
334720f4b53cSBarry Smith $ void bcFunc(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar bcval[])
334856cf3b9cSMatthew G. Knepley 
3349*a4e35b19SJacob Faibussowitsch   If the type is `DM_BC_ESSENTIAL_FIELD` or other _FIELD value, then the calling sequence is\:
3350dce8aebaSBarry Smith .vb
335120f4b53cSBarry Smith   void bcFunc(PetscInt dim, PetscInt Nf, PetscInt NfAux,
3352dce8aebaSBarry Smith               const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
3353dce8aebaSBarry Smith               const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
3354dce8aebaSBarry Smith               PetscReal time, const PetscReal x[], PetscScalar bcval[])
3355dce8aebaSBarry Smith .ve
335656cf3b9cSMatthew G. Knepley + dim - the spatial dimension
335756cf3b9cSMatthew G. Knepley . Nf - the number of fields
335856cf3b9cSMatthew G. Knepley . uOff - the offset into u[] and u_t[] for each field
335956cf3b9cSMatthew G. Knepley . uOff_x - the offset into u_x[] for each field
336056cf3b9cSMatthew G. Knepley . u - each field evaluated at the current point
336156cf3b9cSMatthew G. Knepley . u_t - the time derivative of each field evaluated at the current point
336256cf3b9cSMatthew G. Knepley . u_x - the gradient of each field evaluated at the current point
336356cf3b9cSMatthew G. Knepley . aOff - the offset into a[] and a_t[] for each auxiliary field
336456cf3b9cSMatthew G. Knepley . aOff_x - the offset into a_x[] for each auxiliary field
336556cf3b9cSMatthew G. Knepley . a - each auxiliary field evaluated at the current point
336656cf3b9cSMatthew G. Knepley . a_t - the time derivative of each auxiliary field evaluated at the current point
336756cf3b9cSMatthew G. Knepley . a_x - the gradient of auxiliary each field evaluated at the current point
336856cf3b9cSMatthew G. Knepley . t - current time
336956cf3b9cSMatthew G. Knepley . x - coordinates of the current point
337056cf3b9cSMatthew G. Knepley . numConstants - number of constant parameters
337156cf3b9cSMatthew G. Knepley . constants - constant parameters
337256cf3b9cSMatthew G. Knepley - bcval - output values at the current point
337356cf3b9cSMatthew G. Knepley 
3374*a4e35b19SJacob Faibussowitsch   Notes:
3375*a4e35b19SJacob Faibussowitsch   The pointwise functions are used to provide boundary values for essential boundary
3376*a4e35b19SJacob Faibussowitsch   conditions. In FEM, they are acting upon by dual basis functionals to generate FEM
3377*a4e35b19SJacob Faibussowitsch   coefficients which are fixed. Natural boundary conditions signal to PETSc that boundary
3378*a4e35b19SJacob Faibussowitsch   integrals should be performed, using the kernels from `PetscDSSetBdResidual()`.
3379*a4e35b19SJacob Faibussowitsch 
3380dce8aebaSBarry Smith .seealso: `PetscDS`, `PetscWeakForm`, `DMLabel`, `DMBoundaryConditionType`, `PetscDSAddBoundaryByName()`, `PetscDSGetBoundary()`, `PetscDSSetResidual()`, `PetscDSSetBdResidual()`
338158ebd649SToby Isaac @*/
3382d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDSAddBoundary(PetscDS ds, DMBoundaryConditionType type, const char name[], DMLabel label, PetscInt Nv, const PetscInt values[], PetscInt field, PetscInt Nc, const PetscInt comps[], void (*bcFunc)(void), void (*bcFunc_t)(void), void *ctx, PetscInt *bd)
3383d71ae5a4SJacob Faibussowitsch {
338445480ffeSMatthew G. Knepley   DSBoundary  head = ds->boundary, b;
338545480ffeSMatthew G. Knepley   PetscInt    n    = 0;
338645480ffeSMatthew G. Knepley   const char *lname;
338758ebd649SToby Isaac 
338858ebd649SToby Isaac   PetscFunctionBegin;
338958ebd649SToby Isaac   PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1);
3390783e2ec8SMatthew G. Knepley   PetscValidLogicalCollectiveEnum(ds, type, 2);
33914f572ea9SToby Isaac   PetscAssertPointer(name, 3);
339245480ffeSMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 4);
339345480ffeSMatthew G. Knepley   PetscValidLogicalCollectiveInt(ds, Nv, 5);
339445480ffeSMatthew G. Knepley   PetscValidLogicalCollectiveInt(ds, field, 7);
339545480ffeSMatthew G. Knepley   PetscValidLogicalCollectiveInt(ds, Nc, 8);
3396dce9da9cSMatthew G. Knepley   PetscCheck(field >= 0 && field < ds->Nf, PetscObjectComm((PetscObject)ds), PETSC_ERR_ARG_OUTOFRANGE, "Field %" PetscInt_FMT " is not in [0, %" PetscInt_FMT ")", field, ds->Nf);
3397d57bb9dbSMatthew G. Knepley   if (Nc > 0) {
3398d57bb9dbSMatthew G. Knepley     PetscInt *fcomps;
3399d57bb9dbSMatthew G. Knepley     PetscInt  c;
3400d57bb9dbSMatthew G. Knepley 
34019566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponents(ds, &fcomps));
340263a3b9bcSJacob Faibussowitsch     PetscCheck(Nc <= fcomps[field], PetscObjectComm((PetscObject)ds), PETSC_ERR_ARG_OUTOFRANGE, "Number of constrained components %" PetscInt_FMT " > %" PetscInt_FMT " components for field %" PetscInt_FMT, Nc, fcomps[field], field);
3403d57bb9dbSMatthew G. Knepley     for (c = 0; c < Nc; ++c) {
34041dca8a05SBarry Smith       PetscCheck(comps[c] >= 0 && comps[c] < fcomps[field], PetscObjectComm((PetscObject)ds), PETSC_ERR_ARG_OUTOFRANGE, "Constrained component[%" PetscInt_FMT "] %" PetscInt_FMT " not in [0, %" PetscInt_FMT ") components for field %" PetscInt_FMT, c, comps[c], fcomps[field], field);
3405d57bb9dbSMatthew G. Knepley     }
3406d57bb9dbSMatthew G. Knepley   }
34079566063dSJacob Faibussowitsch   PetscCall(PetscNew(&b));
34089566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(name, (char **)&b->name));
34099566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormCreate(PETSC_COMM_SELF, &b->wf));
34109566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormSetNumFields(b->wf, ds->Nf));
34119566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(Nv, &b->values));
34129566063dSJacob Faibussowitsch   if (Nv) PetscCall(PetscArraycpy(b->values, values, Nv));
34139566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(Nc, &b->comps));
34149566063dSJacob Faibussowitsch   if (Nc) PetscCall(PetscArraycpy(b->comps, comps, Nc));
34159566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetName((PetscObject)label, &lname));
34169566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(lname, (char **)&b->lname));
3417f971fd6bSMatthew G. Knepley   b->type   = type;
341845480ffeSMatthew G. Knepley   b->label  = label;
341945480ffeSMatthew G. Knepley   b->Nv     = Nv;
342058ebd649SToby Isaac   b->field  = field;
342145480ffeSMatthew G. Knepley   b->Nc     = Nc;
342258ebd649SToby Isaac   b->func   = bcFunc;
342356cf3b9cSMatthew G. Knepley   b->func_t = bcFunc_t;
342458ebd649SToby Isaac   b->ctx    = ctx;
342545480ffeSMatthew G. Knepley   b->next   = NULL;
342645480ffeSMatthew G. Knepley   /* Append to linked list so that we can preserve the order */
342745480ffeSMatthew G. Knepley   if (!head) ds->boundary = b;
342845480ffeSMatthew G. Knepley   while (head) {
342945480ffeSMatthew G. Knepley     if (!head->next) {
343045480ffeSMatthew G. Knepley       head->next = b;
343145480ffeSMatthew G. Knepley       head       = b;
343245480ffeSMatthew G. Knepley     }
343345480ffeSMatthew G. Knepley     head = head->next;
343445480ffeSMatthew G. Knepley     ++n;
343545480ffeSMatthew G. Knepley   }
34369371c9d4SSatish Balay   if (bd) {
34374f572ea9SToby Isaac     PetscAssertPointer(bd, 13);
34389371c9d4SSatish Balay     *bd = n;
34399371c9d4SSatish Balay   }
34403ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
344145480ffeSMatthew G. Knepley }
344245480ffeSMatthew G. Knepley 
3443*a4e35b19SJacob Faibussowitsch // PetscClangLinter pragma ignore: -fdoc-section-header-unknown
344445480ffeSMatthew G. Knepley /*@C
3445*a4e35b19SJacob Faibussowitsch   PetscDSAddBoundaryByName - Add a boundary condition to the model.
344645480ffeSMatthew G. Knepley 
344720f4b53cSBarry Smith   Collective
344845480ffeSMatthew G. Knepley 
344945480ffeSMatthew G. Knepley   Input Parameters:
3450dce8aebaSBarry Smith + ds       - The `PetscDS` object
3451dce8aebaSBarry Smith . type     - The type of condition, e.g. `DM_BC_ESSENTIAL`/`DM_BC_ESSENTIAL_FIELD` (Dirichlet), or `DM_BC_NATURAL` (Neumann)
345245480ffeSMatthew G. Knepley . name     - The BC name
345345480ffeSMatthew G. Knepley . lname    - The naem of the label defining constrained points
3454dce8aebaSBarry Smith . Nv       - The number of `DMLabel` values for constrained points
345545480ffeSMatthew G. Knepley . values   - An array of label values for constrained points
345645480ffeSMatthew G. Knepley . field    - The field to constrain
345745480ffeSMatthew G. Knepley . Nc       - The number of constrained field components (0 will constrain all fields)
345845480ffeSMatthew G. Knepley . comps    - An array of constrained component numbers
345945480ffeSMatthew G. Knepley . bcFunc   - A pointwise function giving boundary values
3460a5b23f4aSJose E. Roman . bcFunc_t - A pointwise function giving the time derivative of the boundary values, or NULL
346145480ffeSMatthew G. Knepley - ctx      - An optional user context for bcFunc
346245480ffeSMatthew G. Knepley 
34632fe279fdSBarry Smith   Output Parameter:
346460225df5SJacob Faibussowitsch . bd - The boundary number
346545480ffeSMatthew G. Knepley 
346645480ffeSMatthew G. Knepley   Options Database Keys:
346745480ffeSMatthew G. Knepley + -bc_<boundary name> <num>      - Overrides the boundary ids
346845480ffeSMatthew G. Knepley - -bc_<boundary name>_comp <num> - Overrides the boundary components
346945480ffeSMatthew G. Knepley 
347020f4b53cSBarry Smith   Calling Sequence of `bcFunc` and `bcFunc_t`:
3471dce8aebaSBarry Smith   If the type is `DM_BC_ESSENTIAL`
3472dce8aebaSBarry Smith .vb
347320f4b53cSBarry Smith   void bcFunc(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar bcval[])
3474dce8aebaSBarry Smith .ve
3475dce8aebaSBarry Smith   If the type is `DM_BC_ESSENTIAL_FIELD` or other _FIELD value,
3476dce8aebaSBarry Smith .vb
347720f4b53cSBarry Smith   void bcFunc(PetscInt dim, PetscInt Nf, PetscInt NfAux,
3478dce8aebaSBarry Smith               const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
3479dce8aebaSBarry Smith               const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
3480dce8aebaSBarry Smith               PetscReal time, const PetscReal x[], PetscScalar bcval[])
3481dce8aebaSBarry Smith .ve
348245480ffeSMatthew G. Knepley + dim - the spatial dimension
348345480ffeSMatthew G. Knepley . Nf - the number of fields
348445480ffeSMatthew G. Knepley . uOff - the offset into u[] and u_t[] for each field
348545480ffeSMatthew G. Knepley . uOff_x - the offset into u_x[] for each field
348645480ffeSMatthew G. Knepley . u - each field evaluated at the current point
348745480ffeSMatthew G. Knepley . u_t - the time derivative of each field evaluated at the current point
348845480ffeSMatthew G. Knepley . u_x - the gradient of each field evaluated at the current point
348945480ffeSMatthew G. Knepley . aOff - the offset into a[] and a_t[] for each auxiliary field
349045480ffeSMatthew G. Knepley . aOff_x - the offset into a_x[] for each auxiliary field
349145480ffeSMatthew G. Knepley . a - each auxiliary field evaluated at the current point
349245480ffeSMatthew G. Knepley . a_t - the time derivative of each auxiliary field evaluated at the current point
349345480ffeSMatthew G. Knepley . a_x - the gradient of auxiliary each field evaluated at the current point
349445480ffeSMatthew G. Knepley . t - current time
349545480ffeSMatthew G. Knepley . x - coordinates of the current point
349645480ffeSMatthew G. Knepley . numConstants - number of constant parameters
349745480ffeSMatthew G. Knepley . constants - constant parameters
349845480ffeSMatthew G. Knepley - bcval - output values at the current point
349945480ffeSMatthew G. Knepley 
350045480ffeSMatthew G. Knepley   Level: developer
350145480ffeSMatthew G. Knepley 
3502*a4e35b19SJacob Faibussowitsch   Notes:
3503*a4e35b19SJacob Faibussowitsch   The pointwise functions are used to provide boundary values for essential boundary
3504*a4e35b19SJacob Faibussowitsch   conditions. In FEM, they are acting upon by dual basis functionals to generate FEM
3505*a4e35b19SJacob Faibussowitsch   coefficients which are fixed. Natural boundary conditions signal to PETSc that boundary
3506*a4e35b19SJacob Faibussowitsch   integrals should be performed, using the kernels from `PetscDSSetBdResidual()`.
3507*a4e35b19SJacob Faibussowitsch 
3508dce8aebaSBarry Smith   This function should only be used with `DMFOREST` currently, since labels cannot be defined before the underlying `DMPLEX` is built.
3509dce8aebaSBarry Smith 
3510dce8aebaSBarry Smith .seealso: `PetscDS`, `PetscWeakForm`, `DMLabel`, `DMBoundaryConditionType`, `PetscDSAddBoundary()`, `PetscDSGetBoundary()`, `PetscDSSetResidual()`, `PetscDSSetBdResidual()`
351145480ffeSMatthew G. Knepley @*/
3512d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDSAddBoundaryByName(PetscDS ds, DMBoundaryConditionType type, const char name[], const char lname[], PetscInt Nv, const PetscInt values[], PetscInt field, PetscInt Nc, const PetscInt comps[], void (*bcFunc)(void), void (*bcFunc_t)(void), void *ctx, PetscInt *bd)
3513d71ae5a4SJacob Faibussowitsch {
351445480ffeSMatthew G. Knepley   DSBoundary head = ds->boundary, b;
351545480ffeSMatthew G. Knepley   PetscInt   n    = 0;
351645480ffeSMatthew G. Knepley 
351745480ffeSMatthew G. Knepley   PetscFunctionBegin;
351845480ffeSMatthew G. Knepley   PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1);
351945480ffeSMatthew G. Knepley   PetscValidLogicalCollectiveEnum(ds, type, 2);
35204f572ea9SToby Isaac   PetscAssertPointer(name, 3);
35214f572ea9SToby Isaac   PetscAssertPointer(lname, 4);
352245480ffeSMatthew G. Knepley   PetscValidLogicalCollectiveInt(ds, Nv, 5);
352345480ffeSMatthew G. Knepley   PetscValidLogicalCollectiveInt(ds, field, 7);
352445480ffeSMatthew G. Knepley   PetscValidLogicalCollectiveInt(ds, Nc, 8);
35259566063dSJacob Faibussowitsch   PetscCall(PetscNew(&b));
35269566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(name, (char **)&b->name));
35279566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormCreate(PETSC_COMM_SELF, &b->wf));
35289566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormSetNumFields(b->wf, ds->Nf));
35299566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(Nv, &b->values));
35309566063dSJacob Faibussowitsch   if (Nv) PetscCall(PetscArraycpy(b->values, values, Nv));
35319566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(Nc, &b->comps));
35329566063dSJacob Faibussowitsch   if (Nc) PetscCall(PetscArraycpy(b->comps, comps, Nc));
35339566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(lname, (char **)&b->lname));
353445480ffeSMatthew G. Knepley   b->type   = type;
353545480ffeSMatthew G. Knepley   b->label  = NULL;
353645480ffeSMatthew G. Knepley   b->Nv     = Nv;
353745480ffeSMatthew G. Knepley   b->field  = field;
353845480ffeSMatthew G. Knepley   b->Nc     = Nc;
353945480ffeSMatthew G. Knepley   b->func   = bcFunc;
354045480ffeSMatthew G. Knepley   b->func_t = bcFunc_t;
354145480ffeSMatthew G. Knepley   b->ctx    = ctx;
354245480ffeSMatthew G. Knepley   b->next   = NULL;
354345480ffeSMatthew G. Knepley   /* Append to linked list so that we can preserve the order */
354445480ffeSMatthew G. Knepley   if (!head) ds->boundary = b;
354545480ffeSMatthew G. Knepley   while (head) {
354645480ffeSMatthew G. Knepley     if (!head->next) {
354745480ffeSMatthew G. Knepley       head->next = b;
354845480ffeSMatthew G. Knepley       head       = b;
354945480ffeSMatthew G. Knepley     }
355045480ffeSMatthew G. Knepley     head = head->next;
355145480ffeSMatthew G. Knepley     ++n;
355245480ffeSMatthew G. Knepley   }
35539371c9d4SSatish Balay   if (bd) {
35544f572ea9SToby Isaac     PetscAssertPointer(bd, 13);
35559371c9d4SSatish Balay     *bd = n;
35569371c9d4SSatish Balay   }
35573ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
355858ebd649SToby Isaac }
355958ebd649SToby Isaac 
3560b67eacb3SMatthew G. Knepley /*@C
3561*a4e35b19SJacob Faibussowitsch   PetscDSUpdateBoundary - Change a boundary condition for the model.
3562b67eacb3SMatthew G. Knepley 
3563b67eacb3SMatthew G. Knepley   Input Parameters:
3564dce8aebaSBarry Smith + ds       - The `PetscDS` object
3565b67eacb3SMatthew G. Knepley . bd       - The boundary condition number
3566dce8aebaSBarry Smith . type     - The type of condition, e.g. `DM_BC_ESSENTIAL`/`DM_BC_ESSENTIAL_FIELD` (Dirichlet), or `DM_BC_NATURAL` (Neumann)
3567b67eacb3SMatthew G. Knepley . name     - The BC name
356845480ffeSMatthew G. Knepley . label    - The label defining constrained points
3569dce8aebaSBarry Smith . Nv       - The number of `DMLabel` ids for constrained points
357045480ffeSMatthew G. Knepley . values   - An array of ids for constrained points
3571b67eacb3SMatthew G. Knepley . field    - The field to constrain
357245480ffeSMatthew G. Knepley . Nc       - The number of constrained field components
3573b67eacb3SMatthew G. Knepley . comps    - An array of constrained component numbers
3574b67eacb3SMatthew G. Knepley . bcFunc   - A pointwise function giving boundary values
3575a5b23f4aSJose E. Roman . bcFunc_t - A pointwise function giving the time derivative of the boundary values, or NULL
3576b67eacb3SMatthew G. Knepley - ctx      - An optional user context for bcFunc
3577b67eacb3SMatthew G. Knepley 
3578b67eacb3SMatthew G. Knepley   Level: developer
3579b67eacb3SMatthew G. Knepley 
3580*a4e35b19SJacob Faibussowitsch   Notes:
3581*a4e35b19SJacob Faibussowitsch   The pointwise functions are used to provide boundary values for essential boundary
3582*a4e35b19SJacob Faibussowitsch   conditions. In FEM, they are acting upon by dual basis functionals to generate FEM
3583*a4e35b19SJacob Faibussowitsch   coefficients which are fixed. Natural boundary conditions signal to PETSc that boundary
3584*a4e35b19SJacob Faibussowitsch   integrals should be performed, using the kernels from `PetscDSSetBdResidual()`.
3585*a4e35b19SJacob Faibussowitsch 
3586dce8aebaSBarry Smith   The boundary condition number is the order in which it was registered. The user can get the number of boundary conditions from `PetscDSGetNumBoundary()`.
3587dce8aebaSBarry Smith   See `PetscDSAddBoundary()` for a description of the calling sequences for the callbacks.
3588dce8aebaSBarry Smith 
3589dce8aebaSBarry Smith .seealso: `PetscDS`, `PetscWeakForm`, `DMBoundaryConditionType`, `PetscDSAddBoundary()`, `PetscDSGetBoundary()`, `PetscDSGetNumBoundary()`, `DMLabel`
3590b67eacb3SMatthew G. Knepley @*/
3591d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDSUpdateBoundary(PetscDS ds, PetscInt bd, DMBoundaryConditionType type, const char name[], DMLabel label, PetscInt Nv, const PetscInt values[], PetscInt field, PetscInt Nc, const PetscInt comps[], void (*bcFunc)(void), void (*bcFunc_t)(void), void *ctx)
3592d71ae5a4SJacob Faibussowitsch {
3593b67eacb3SMatthew G. Knepley   DSBoundary b = ds->boundary;
3594b67eacb3SMatthew G. Knepley   PetscInt   n = 0;
3595b67eacb3SMatthew G. Knepley 
3596b67eacb3SMatthew G. Knepley   PetscFunctionBegin;
3597b67eacb3SMatthew G. Knepley   PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1);
3598b67eacb3SMatthew G. Knepley   while (b) {
3599b67eacb3SMatthew G. Knepley     if (n == bd) break;
3600b67eacb3SMatthew G. Knepley     b = b->next;
3601b67eacb3SMatthew G. Knepley     ++n;
3602b67eacb3SMatthew G. Knepley   }
360363a3b9bcSJacob Faibussowitsch   PetscCheck(b, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Boundary %" PetscInt_FMT " is not in [0, %" PetscInt_FMT ")", bd, n);
3604b67eacb3SMatthew G. Knepley   if (name) {
36059566063dSJacob Faibussowitsch     PetscCall(PetscFree(b->name));
36069566063dSJacob Faibussowitsch     PetscCall(PetscStrallocpy(name, (char **)&b->name));
3607b67eacb3SMatthew G. Knepley   }
3608b67eacb3SMatthew G. Knepley   b->type = type;
360945480ffeSMatthew G. Knepley   if (label) {
361045480ffeSMatthew G. Knepley     const char *name;
361145480ffeSMatthew G. Knepley 
361245480ffeSMatthew G. Knepley     b->label = label;
36139566063dSJacob Faibussowitsch     PetscCall(PetscFree(b->lname));
36149566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetName((PetscObject)label, &name));
36159566063dSJacob Faibussowitsch     PetscCall(PetscStrallocpy(name, (char **)&b->lname));
361645480ffeSMatthew G. Knepley   }
361745480ffeSMatthew G. Knepley   if (Nv >= 0) {
361845480ffeSMatthew G. Knepley     b->Nv = Nv;
36199566063dSJacob Faibussowitsch     PetscCall(PetscFree(b->values));
36209566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(Nv, &b->values));
36219566063dSJacob Faibussowitsch     if (Nv) PetscCall(PetscArraycpy(b->values, values, Nv));
362245480ffeSMatthew G. Knepley   }
362345480ffeSMatthew G. Knepley   if (field >= 0) b->field = field;
362445480ffeSMatthew G. Knepley   if (Nc >= 0) {
362545480ffeSMatthew G. Knepley     b->Nc = Nc;
36269566063dSJacob Faibussowitsch     PetscCall(PetscFree(b->comps));
36279566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(Nc, &b->comps));
36289566063dSJacob Faibussowitsch     if (Nc) PetscCall(PetscArraycpy(b->comps, comps, Nc));
362945480ffeSMatthew G. Knepley   }
363045480ffeSMatthew G. Knepley   if (bcFunc) b->func = bcFunc;
363145480ffeSMatthew G. Knepley   if (bcFunc_t) b->func_t = bcFunc_t;
363245480ffeSMatthew G. Knepley   if (ctx) b->ctx = ctx;
36333ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3634b67eacb3SMatthew G. Knepley }
3635b67eacb3SMatthew G. Knepley 
363658ebd649SToby Isaac /*@
363758ebd649SToby Isaac   PetscDSGetNumBoundary - Get the number of registered BC
363858ebd649SToby Isaac 
36392fe279fdSBarry Smith   Input Parameter:
3640dce8aebaSBarry Smith . ds - The `PetscDS` object
364158ebd649SToby Isaac 
36422fe279fdSBarry Smith   Output Parameter:
364358ebd649SToby Isaac . numBd - The number of BC
364458ebd649SToby Isaac 
364558ebd649SToby Isaac   Level: intermediate
364658ebd649SToby Isaac 
3647dce8aebaSBarry Smith .seealso: `PetscDS`, `PetscDSAddBoundary()`, `PetscDSGetBoundary()`
364858ebd649SToby Isaac @*/
3649d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDSGetNumBoundary(PetscDS ds, PetscInt *numBd)
3650d71ae5a4SJacob Faibussowitsch {
365158ebd649SToby Isaac   DSBoundary b = ds->boundary;
365258ebd649SToby Isaac 
365358ebd649SToby Isaac   PetscFunctionBegin;
365458ebd649SToby Isaac   PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1);
36554f572ea9SToby Isaac   PetscAssertPointer(numBd, 2);
365658ebd649SToby Isaac   *numBd = 0;
36579371c9d4SSatish Balay   while (b) {
36589371c9d4SSatish Balay     ++(*numBd);
36599371c9d4SSatish Balay     b = b->next;
36609371c9d4SSatish Balay   }
36613ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
366258ebd649SToby Isaac }
366358ebd649SToby Isaac 
366458ebd649SToby Isaac /*@C
36659a6efb6aSMatthew G. Knepley   PetscDSGetBoundary - Gets a boundary condition to the model
366658ebd649SToby Isaac 
366758ebd649SToby Isaac   Input Parameters:
3668dce8aebaSBarry Smith + ds - The `PetscDS` object
366958ebd649SToby Isaac - bd - The BC number
367058ebd649SToby Isaac 
367158ebd649SToby Isaac   Output Parameters:
3672dce8aebaSBarry Smith + wf     - The `PetscWeakForm` holding the pointwise functions
3673dce8aebaSBarry Smith . type   - The type of condition, e.g. `DM_BC_ESSENTIAL`/`DM_BC_ESSENTIAL_FIELD` (Dirichlet), or `DM_BC_NATURAL` (Neumann)
367458ebd649SToby Isaac . name   - The BC name
367545480ffeSMatthew G. Knepley . label  - The label defining constrained points
3676dce8aebaSBarry Smith . Nv     - The number of `DMLabel` ids for constrained points
367745480ffeSMatthew G. Knepley . values - An array of ids for constrained points
367858ebd649SToby Isaac . field  - The field to constrain
367945480ffeSMatthew G. Knepley . Nc     - The number of constrained field components
368058ebd649SToby Isaac . comps  - An array of constrained component numbers
368160225df5SJacob Faibussowitsch . func   - A pointwise function giving boundary values
368260225df5SJacob Faibussowitsch . func_t - A pointwise function giving the time derivative of the boundary values
368358ebd649SToby Isaac - ctx    - An optional user context for bcFunc
368458ebd649SToby Isaac 
368558ebd649SToby Isaac   Options Database Keys:
368658ebd649SToby Isaac + -bc_<boundary name> <num>      - Overrides the boundary ids
368758ebd649SToby Isaac - -bc_<boundary name>_comp <num> - Overrides the boundary components
368858ebd649SToby Isaac 
368958ebd649SToby Isaac   Level: developer
369058ebd649SToby Isaac 
3691dce8aebaSBarry Smith .seealso: `PetscDS`, `PetscWeakForm`, `DMBoundaryConditionType`, `PetscDSAddBoundary()`, `DMLabel`
369258ebd649SToby Isaac @*/
3693d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDSGetBoundary(PetscDS ds, PetscInt bd, PetscWeakForm *wf, DMBoundaryConditionType *type, const char *name[], DMLabel *label, PetscInt *Nv, const PetscInt *values[], PetscInt *field, PetscInt *Nc, const PetscInt *comps[], void (**func)(void), void (**func_t)(void), void **ctx)
3694d71ae5a4SJacob Faibussowitsch {
369558ebd649SToby Isaac   DSBoundary b = ds->boundary;
369658ebd649SToby Isaac   PetscInt   n = 0;
369758ebd649SToby Isaac 
369858ebd649SToby Isaac   PetscFunctionBegin;
369958ebd649SToby Isaac   PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1);
370058ebd649SToby Isaac   while (b) {
370158ebd649SToby Isaac     if (n == bd) break;
370258ebd649SToby Isaac     b = b->next;
370358ebd649SToby Isaac     ++n;
370458ebd649SToby Isaac   }
370563a3b9bcSJacob Faibussowitsch   PetscCheck(b, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Boundary %" PetscInt_FMT " is not in [0, %" PetscInt_FMT ")", bd, n);
370645480ffeSMatthew G. Knepley   if (wf) {
37074f572ea9SToby Isaac     PetscAssertPointer(wf, 3);
370845480ffeSMatthew G. Knepley     *wf = b->wf;
370945480ffeSMatthew G. Knepley   }
3710f971fd6bSMatthew G. Knepley   if (type) {
37114f572ea9SToby Isaac     PetscAssertPointer(type, 4);
3712f971fd6bSMatthew G. Knepley     *type = b->type;
371358ebd649SToby Isaac   }
371458ebd649SToby Isaac   if (name) {
37154f572ea9SToby Isaac     PetscAssertPointer(name, 5);
371658ebd649SToby Isaac     *name = b->name;
371758ebd649SToby Isaac   }
371845480ffeSMatthew G. Knepley   if (label) {
37194f572ea9SToby Isaac     PetscAssertPointer(label, 6);
372045480ffeSMatthew G. Knepley     *label = b->label;
372145480ffeSMatthew G. Knepley   }
372245480ffeSMatthew G. Knepley   if (Nv) {
37234f572ea9SToby Isaac     PetscAssertPointer(Nv, 7);
372445480ffeSMatthew G. Knepley     *Nv = b->Nv;
372545480ffeSMatthew G. Knepley   }
372645480ffeSMatthew G. Knepley   if (values) {
37274f572ea9SToby Isaac     PetscAssertPointer(values, 8);
372845480ffeSMatthew G. Knepley     *values = b->values;
372958ebd649SToby Isaac   }
373058ebd649SToby Isaac   if (field) {
37314f572ea9SToby Isaac     PetscAssertPointer(field, 9);
373258ebd649SToby Isaac     *field = b->field;
373358ebd649SToby Isaac   }
373445480ffeSMatthew G. Knepley   if (Nc) {
37354f572ea9SToby Isaac     PetscAssertPointer(Nc, 10);
373645480ffeSMatthew G. Knepley     *Nc = b->Nc;
373758ebd649SToby Isaac   }
373858ebd649SToby Isaac   if (comps) {
37394f572ea9SToby Isaac     PetscAssertPointer(comps, 11);
374058ebd649SToby Isaac     *comps = b->comps;
374158ebd649SToby Isaac   }
374258ebd649SToby Isaac   if (func) {
37434f572ea9SToby Isaac     PetscAssertPointer(func, 12);
374458ebd649SToby Isaac     *func = b->func;
374558ebd649SToby Isaac   }
374656cf3b9cSMatthew G. Knepley   if (func_t) {
37474f572ea9SToby Isaac     PetscAssertPointer(func_t, 13);
374856cf3b9cSMatthew G. Knepley     *func_t = b->func_t;
374956cf3b9cSMatthew G. Knepley   }
375058ebd649SToby Isaac   if (ctx) {
37514f572ea9SToby Isaac     PetscAssertPointer(ctx, 14);
375258ebd649SToby Isaac     *ctx = b->ctx;
375358ebd649SToby Isaac   }
37543ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
375558ebd649SToby Isaac }
375658ebd649SToby Isaac 
3757d71ae5a4SJacob Faibussowitsch static PetscErrorCode DSBoundaryDuplicate_Internal(DSBoundary b, DSBoundary *bNew)
3758d71ae5a4SJacob Faibussowitsch {
375945480ffeSMatthew G. Knepley   PetscFunctionBegin;
37609566063dSJacob Faibussowitsch   PetscCall(PetscNew(bNew));
37619566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormCreate(PETSC_COMM_SELF, &(*bNew)->wf));
37629566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormCopy(b->wf, (*bNew)->wf));
37639566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(b->name, (char **)&((*bNew)->name)));
37649566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(b->lname, (char **)&((*bNew)->lname)));
376545480ffeSMatthew G. Knepley   (*bNew)->type  = b->type;
376645480ffeSMatthew G. Knepley   (*bNew)->label = b->label;
376745480ffeSMatthew G. Knepley   (*bNew)->Nv    = b->Nv;
37689566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(b->Nv, &(*bNew)->values));
37699566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy((*bNew)->values, b->values, b->Nv));
377045480ffeSMatthew G. Knepley   (*bNew)->field = b->field;
377145480ffeSMatthew G. Knepley   (*bNew)->Nc    = b->Nc;
37729566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(b->Nc, &(*bNew)->comps));
37739566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy((*bNew)->comps, b->comps, b->Nc));
377445480ffeSMatthew G. Knepley   (*bNew)->func   = b->func;
377545480ffeSMatthew G. Knepley   (*bNew)->func_t = b->func_t;
377645480ffeSMatthew G. Knepley   (*bNew)->ctx    = b->ctx;
37773ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
377845480ffeSMatthew G. Knepley }
377945480ffeSMatthew G. Knepley 
37809252d075SMatthew G. Knepley /*@
37819252d075SMatthew G. Knepley   PetscDSCopyBoundary - Copy all boundary condition objects to the new problem
37829252d075SMatthew G. Knepley 
378320f4b53cSBarry Smith   Not Collective
37849252d075SMatthew G. Knepley 
378536951cb5SMatthew G. Knepley   Input Parameters:
3786dce8aebaSBarry Smith + ds        - The source `PetscDS` object
3787dce8aebaSBarry Smith . numFields - The number of selected fields, or `PETSC_DEFAULT` for all fields
378836951cb5SMatthew G. Knepley - fields    - The selected fields, or NULL for all fields
37899252d075SMatthew G. Knepley 
37909252d075SMatthew G. Knepley   Output Parameter:
3791dce8aebaSBarry Smith . newds - The target `PetscDS`, now with a copy of the boundary conditions
37929252d075SMatthew G. Knepley 
37939252d075SMatthew G. Knepley   Level: intermediate
37949252d075SMatthew G. Knepley 
3795dce8aebaSBarry Smith .seealso: `PetscDS`, `DMBoundary`, `PetscDSCopyEquations()`, `PetscDSSetResidual()`, `PetscDSSetJacobian()`, `PetscDSSetRiemannSolver()`, `PetscDSSetBdResidual()`, `PetscDSSetBdJacobian()`, `PetscDSCreate()`
37969252d075SMatthew G. Knepley @*/
3797d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDSCopyBoundary(PetscDS ds, PetscInt numFields, const PetscInt fields[], PetscDS newds)
3798d71ae5a4SJacob Faibussowitsch {
379945480ffeSMatthew G. Knepley   DSBoundary b, *lastnext;
3800dff059c6SToby Isaac 
3801dff059c6SToby Isaac   PetscFunctionBegin;
380236951cb5SMatthew G. Knepley   PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1);
380336951cb5SMatthew G. Knepley   PetscValidHeaderSpecific(newds, PETSCDS_CLASSID, 4);
38043ba16761SJacob Faibussowitsch   if (ds == newds) PetscFunctionReturn(PETSC_SUCCESS);
38059566063dSJacob Faibussowitsch   PetscCall(PetscDSDestroyBoundary(newds));
380636951cb5SMatthew G. Knepley   lastnext = &(newds->boundary);
380736951cb5SMatthew G. Knepley   for (b = ds->boundary; b; b = b->next) {
3808dff059c6SToby Isaac     DSBoundary bNew;
380936951cb5SMatthew G. Knepley     PetscInt   fieldNew = -1;
3810dff059c6SToby Isaac 
381136951cb5SMatthew G. Knepley     if (numFields > 0 && fields) {
381236951cb5SMatthew G. Knepley       PetscInt f;
381336951cb5SMatthew G. Knepley 
38149371c9d4SSatish Balay       for (f = 0; f < numFields; ++f)
38159371c9d4SSatish Balay         if (b->field == fields[f]) break;
381636951cb5SMatthew G. Knepley       if (f == numFields) continue;
381736951cb5SMatthew G. Knepley       fieldNew = f;
381836951cb5SMatthew G. Knepley     }
38199566063dSJacob Faibussowitsch     PetscCall(DSBoundaryDuplicate_Internal(b, &bNew));
382036951cb5SMatthew G. Knepley     bNew->field = fieldNew < 0 ? b->field : fieldNew;
3821dff059c6SToby Isaac     *lastnext   = bNew;
3822dff059c6SToby Isaac     lastnext    = &(bNew->next);
3823dff059c6SToby Isaac   }
38243ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3825dff059c6SToby Isaac }
3826dff059c6SToby Isaac 
38276c1eb96dSMatthew G. Knepley /*@
3828dce8aebaSBarry Smith   PetscDSDestroyBoundary - Remove all `DMBoundary` objects from the `PetscDS`
382945480ffeSMatthew G. Knepley 
383020f4b53cSBarry Smith   Not Collective
383145480ffeSMatthew G. Knepley 
383245480ffeSMatthew G. Knepley   Input Parameter:
3833dce8aebaSBarry Smith . ds - The `PetscDS` object
383445480ffeSMatthew G. Knepley 
383545480ffeSMatthew G. Knepley   Level: intermediate
383645480ffeSMatthew G. Knepley 
3837dce8aebaSBarry Smith .seealso: `PetscDS`, `DMBoundary`, `PetscDSCopyBoundary()`, `PetscDSCopyEquations()`
383845480ffeSMatthew G. Knepley @*/
3839d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDSDestroyBoundary(PetscDS ds)
3840d71ae5a4SJacob Faibussowitsch {
384145480ffeSMatthew G. Knepley   DSBoundary next = ds->boundary;
384245480ffeSMatthew G. Knepley 
384345480ffeSMatthew G. Knepley   PetscFunctionBegin;
384445480ffeSMatthew G. Knepley   while (next) {
384545480ffeSMatthew G. Knepley     DSBoundary b = next;
384645480ffeSMatthew G. Knepley 
384745480ffeSMatthew G. Knepley     next = b->next;
38489566063dSJacob Faibussowitsch     PetscCall(PetscWeakFormDestroy(&b->wf));
38499566063dSJacob Faibussowitsch     PetscCall(PetscFree(b->name));
38509566063dSJacob Faibussowitsch     PetscCall(PetscFree(b->lname));
38519566063dSJacob Faibussowitsch     PetscCall(PetscFree(b->values));
38529566063dSJacob Faibussowitsch     PetscCall(PetscFree(b->comps));
38539566063dSJacob Faibussowitsch     PetscCall(PetscFree(b));
385445480ffeSMatthew G. Knepley   }
38553ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
385645480ffeSMatthew G. Knepley }
385745480ffeSMatthew G. Knepley 
385845480ffeSMatthew G. Knepley /*@
38596c1eb96dSMatthew G. Knepley   PetscDSSelectDiscretizations - Copy discretizations to the new problem with different field layout
38606c1eb96dSMatthew G. Knepley 
386120f4b53cSBarry Smith   Not Collective
38626c1eb96dSMatthew G. Knepley 
3863d8d19677SJose E. Roman   Input Parameters:
3864dce8aebaSBarry Smith + prob      - The `PetscDS` object
38656c1eb96dSMatthew G. Knepley . numFields - Number of new fields
38666c1eb96dSMatthew G. Knepley - fields    - Old field number for each new field
38676c1eb96dSMatthew G. Knepley 
38686c1eb96dSMatthew G. Knepley   Output Parameter:
3869dce8aebaSBarry Smith . newprob - The `PetscDS` copy
38706c1eb96dSMatthew G. Knepley 
38716c1eb96dSMatthew G. Knepley   Level: intermediate
38726c1eb96dSMatthew G. Knepley 
3873dce8aebaSBarry Smith .seealso: `PetscDS`, `PetscDSSelectEquations()`, `PetscDSCopyBoundary()`, `PetscDSSetResidual()`, `PetscDSSetJacobian()`, `PetscDSSetRiemannSolver()`, `PetscDSSetBdResidual()`, `PetscDSSetBdJacobian()`, `PetscDSCreate()`
38746c1eb96dSMatthew G. Knepley @*/
3875d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDSSelectDiscretizations(PetscDS prob, PetscInt numFields, const PetscInt fields[], PetscDS newprob)
3876d71ae5a4SJacob Faibussowitsch {
38776c1eb96dSMatthew G. Knepley   PetscInt Nf, Nfn, fn;
38786c1eb96dSMatthew G. Knepley 
38796c1eb96dSMatthew G. Knepley   PetscFunctionBegin;
38806c1eb96dSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
38814f572ea9SToby Isaac   if (fields) PetscAssertPointer(fields, 3);
38826c1eb96dSMatthew G. Knepley   PetscValidHeaderSpecific(newprob, PETSCDS_CLASSID, 4);
38839566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(prob, &Nf));
38849566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(newprob, &Nfn));
388545480ffeSMatthew G. Knepley   numFields = numFields < 0 ? Nf : numFields;
38866c1eb96dSMatthew G. Knepley   for (fn = 0; fn < numFields; ++fn) {
38876c1eb96dSMatthew G. Knepley     const PetscInt f = fields ? fields[fn] : fn;
38886c1eb96dSMatthew G. Knepley     PetscObject    disc;
38896c1eb96dSMatthew G. Knepley 
38906c1eb96dSMatthew G. Knepley     if (f >= Nf) continue;
38919566063dSJacob Faibussowitsch     PetscCall(PetscDSGetDiscretization(prob, f, &disc));
38929566063dSJacob Faibussowitsch     PetscCall(PetscDSSetDiscretization(newprob, fn, disc));
38936c1eb96dSMatthew G. Knepley   }
38943ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
38956c1eb96dSMatthew G. Knepley }
38966c1eb96dSMatthew G. Knepley 
38976c1eb96dSMatthew G. Knepley /*@
38989252d075SMatthew G. Knepley   PetscDSSelectEquations - Copy pointwise function pointers to the new problem with different field layout
38999252d075SMatthew G. Knepley 
390020f4b53cSBarry Smith   Not Collective
39019252d075SMatthew G. Knepley 
3902d8d19677SJose E. Roman   Input Parameters:
3903dce8aebaSBarry Smith + prob      - The `PetscDS` object
39049252d075SMatthew G. Knepley . numFields - Number of new fields
39059252d075SMatthew G. Knepley - fields    - Old field number for each new field
39069252d075SMatthew G. Knepley 
39079252d075SMatthew G. Knepley   Output Parameter:
3908dce8aebaSBarry Smith . newprob - The `PetscDS` copy
39099252d075SMatthew G. Knepley 
39109252d075SMatthew G. Knepley   Level: intermediate
39119252d075SMatthew G. Knepley 
3912dce8aebaSBarry Smith .seealso: `PetscDS`, `PetscDSSelectDiscretizations()`, `PetscDSCopyBoundary()`, `PetscDSSetResidual()`, `PetscDSSetJacobian()`, `PetscDSSetRiemannSolver()`, `PetscDSSetBdResidual()`, `PetscDSSetBdJacobian()`, `PetscDSCreate()`
39139252d075SMatthew G. Knepley @*/
3914d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDSSelectEquations(PetscDS prob, PetscInt numFields, const PetscInt fields[], PetscDS newprob)
3915d71ae5a4SJacob Faibussowitsch {
39169252d075SMatthew G. Knepley   PetscInt Nf, Nfn, fn, gn;
39179252d075SMatthew G. Knepley 
39189252d075SMatthew G. Knepley   PetscFunctionBegin;
39199252d075SMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
39204f572ea9SToby Isaac   if (fields) PetscAssertPointer(fields, 3);
39219252d075SMatthew G. Knepley   PetscValidHeaderSpecific(newprob, PETSCDS_CLASSID, 4);
39229566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(prob, &Nf));
39239566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(newprob, &Nfn));
392463a3b9bcSJacob Faibussowitsch   PetscCheck(numFields <= Nfn, PetscObjectComm((PetscObject)prob), PETSC_ERR_ARG_SIZ, "Number of fields %" PetscInt_FMT " to transfer must not be greater then the total number of fields %" PetscInt_FMT, numFields, Nfn);
39259252d075SMatthew G. Knepley   for (fn = 0; fn < numFields; ++fn) {
39269252d075SMatthew G. Knepley     const PetscInt   f = fields ? fields[fn] : fn;
39279252d075SMatthew G. Knepley     PetscPointFunc   obj;
39289252d075SMatthew G. Knepley     PetscPointFunc   f0, f1;
39299252d075SMatthew G. Knepley     PetscBdPointFunc f0Bd, f1Bd;
39309252d075SMatthew G. Knepley     PetscRiemannFunc r;
39319252d075SMatthew G. Knepley 
3932c52f1e13SMatthew G. Knepley     if (f >= Nf) continue;
39339566063dSJacob Faibussowitsch     PetscCall(PetscDSGetObjective(prob, f, &obj));
39349566063dSJacob Faibussowitsch     PetscCall(PetscDSGetResidual(prob, f, &f0, &f1));
39359566063dSJacob Faibussowitsch     PetscCall(PetscDSGetBdResidual(prob, f, &f0Bd, &f1Bd));
39369566063dSJacob Faibussowitsch     PetscCall(PetscDSGetRiemannSolver(prob, f, &r));
39379566063dSJacob Faibussowitsch     PetscCall(PetscDSSetObjective(newprob, fn, obj));
39389566063dSJacob Faibussowitsch     PetscCall(PetscDSSetResidual(newprob, fn, f0, f1));
39399566063dSJacob Faibussowitsch     PetscCall(PetscDSSetBdResidual(newprob, fn, f0Bd, f1Bd));
39409566063dSJacob Faibussowitsch     PetscCall(PetscDSSetRiemannSolver(newprob, fn, r));
39419252d075SMatthew G. Knepley     for (gn = 0; gn < numFields; ++gn) {
39429252d075SMatthew G. Knepley       const PetscInt  g = fields ? fields[gn] : gn;
39439252d075SMatthew G. Knepley       PetscPointJac   g0, g1, g2, g3;
39449252d075SMatthew G. Knepley       PetscPointJac   g0p, g1p, g2p, g3p;
39459252d075SMatthew G. Knepley       PetscBdPointJac g0Bd, g1Bd, g2Bd, g3Bd;
39469252d075SMatthew G. Knepley 
3947c52f1e13SMatthew G. Knepley       if (g >= Nf) continue;
39489566063dSJacob Faibussowitsch       PetscCall(PetscDSGetJacobian(prob, f, g, &g0, &g1, &g2, &g3));
39499566063dSJacob Faibussowitsch       PetscCall(PetscDSGetJacobianPreconditioner(prob, f, g, &g0p, &g1p, &g2p, &g3p));
39509566063dSJacob Faibussowitsch       PetscCall(PetscDSGetBdJacobian(prob, f, g, &g0Bd, &g1Bd, &g2Bd, &g3Bd));
39519566063dSJacob Faibussowitsch       PetscCall(PetscDSSetJacobian(newprob, fn, gn, g0, g1, g2, g3));
39529566063dSJacob Faibussowitsch       PetscCall(PetscDSSetJacobianPreconditioner(newprob, fn, gn, g0p, g1p, g2p, g3p));
39539566063dSJacob Faibussowitsch       PetscCall(PetscDSSetBdJacobian(newprob, fn, gn, g0Bd, g1Bd, g2Bd, g3Bd));
39549252d075SMatthew G. Knepley     }
39559252d075SMatthew G. Knepley   }
39563ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
39579252d075SMatthew G. Knepley }
39589252d075SMatthew G. Knepley 
3959da51fcedSMatthew G. Knepley /*@
3960dce8aebaSBarry Smith   PetscDSCopyEquations - Copy all pointwise function pointers to another `PetscDS`
3961da51fcedSMatthew G. Knepley 
396220f4b53cSBarry Smith   Not Collective
3963da51fcedSMatthew G. Knepley 
3964da51fcedSMatthew G. Knepley   Input Parameter:
3965dce8aebaSBarry Smith . prob - The `PetscDS` object
3966da51fcedSMatthew G. Knepley 
3967da51fcedSMatthew G. Knepley   Output Parameter:
3968dce8aebaSBarry Smith . newprob - The `PetscDS` copy
3969da51fcedSMatthew G. Knepley 
3970da51fcedSMatthew G. Knepley   Level: intermediate
3971da51fcedSMatthew G. Knepley 
3972dce8aebaSBarry Smith .seealso: `PetscDS`, `PetscDSCopyBoundary()`, `PetscDSSetResidual()`, `PetscDSSetJacobian()`, `PetscDSSetRiemannSolver()`, `PetscDSSetBdResidual()`, `PetscDSSetBdJacobian()`, `PetscDSCreate()`
3973da51fcedSMatthew G. Knepley @*/
3974d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDSCopyEquations(PetscDS prob, PetscDS newprob)
3975d71ae5a4SJacob Faibussowitsch {
3976b8025e53SMatthew G. Knepley   PetscWeakForm wf, newwf;
39779252d075SMatthew G. Knepley   PetscInt      Nf, Ng;
3978da51fcedSMatthew G. Knepley 
3979da51fcedSMatthew G. Knepley   PetscFunctionBegin;
3980da51fcedSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
3981da51fcedSMatthew G. Knepley   PetscValidHeaderSpecific(newprob, PETSCDS_CLASSID, 2);
39829566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(prob, &Nf));
39839566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(newprob, &Ng));
398463a3b9bcSJacob Faibussowitsch   PetscCheck(Nf == Ng, PetscObjectComm((PetscObject)prob), PETSC_ERR_ARG_SIZ, "Number of fields must match %" PetscInt_FMT " != %" PetscInt_FMT, Nf, Ng);
39859566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWeakForm(prob, &wf));
39869566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWeakForm(newprob, &newwf));
39879566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormCopy(wf, newwf));
39883ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
39899252d075SMatthew G. Knepley }
399045480ffeSMatthew G. Knepley 
39919252d075SMatthew G. Knepley /*@
3992dce8aebaSBarry Smith   PetscDSCopyConstants - Copy all constants to another `PetscDS`
3993da51fcedSMatthew G. Knepley 
399420f4b53cSBarry Smith   Not Collective
39959252d075SMatthew G. Knepley 
39969252d075SMatthew G. Knepley   Input Parameter:
3997dce8aebaSBarry Smith . prob - The `PetscDS` object
39989252d075SMatthew G. Knepley 
39999252d075SMatthew G. Knepley   Output Parameter:
4000dce8aebaSBarry Smith . newprob - The `PetscDS` copy
40019252d075SMatthew G. Knepley 
40029252d075SMatthew G. Knepley   Level: intermediate
40039252d075SMatthew G. Knepley 
4004dce8aebaSBarry Smith .seealso: `PetscDS`, `PetscDSCopyBoundary()`, `PetscDSCopyEquations()`, `PetscDSSetResidual()`, `PetscDSSetJacobian()`, `PetscDSSetRiemannSolver()`, `PetscDSSetBdResidual()`, `PetscDSSetBdJacobian()`, `PetscDSCreate()`
40059252d075SMatthew G. Knepley @*/
4006d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDSCopyConstants(PetscDS prob, PetscDS newprob)
4007d71ae5a4SJacob Faibussowitsch {
40089252d075SMatthew G. Knepley   PetscInt           Nc;
40099252d075SMatthew G. Knepley   const PetscScalar *constants;
40109252d075SMatthew G. Knepley 
40119252d075SMatthew G. Knepley   PetscFunctionBegin;
40129252d075SMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
40139252d075SMatthew G. Knepley   PetscValidHeaderSpecific(newprob, PETSCDS_CLASSID, 2);
40149566063dSJacob Faibussowitsch   PetscCall(PetscDSGetConstants(prob, &Nc, &constants));
40159566063dSJacob Faibussowitsch   PetscCall(PetscDSSetConstants(newprob, Nc, (PetscScalar *)constants));
40163ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4017da51fcedSMatthew G. Knepley }
4018da51fcedSMatthew G. Knepley 
401945480ffeSMatthew G. Knepley /*@
4020dce8aebaSBarry Smith   PetscDSCopyExactSolutions - Copy all exact solutions to another `PetscDS`
402145480ffeSMatthew G. Knepley 
402220f4b53cSBarry Smith   Not Collective
402345480ffeSMatthew G. Knepley 
402445480ffeSMatthew G. Knepley   Input Parameter:
4025dce8aebaSBarry Smith . ds - The `PetscDS` object
402645480ffeSMatthew G. Knepley 
402745480ffeSMatthew G. Knepley   Output Parameter:
4028dce8aebaSBarry Smith . newds - The `PetscDS` copy
402945480ffeSMatthew G. Knepley 
403045480ffeSMatthew G. Knepley   Level: intermediate
403145480ffeSMatthew G. Knepley 
4032dce8aebaSBarry Smith .seealso: `PetscDS`, `PetscDSCopyBoundary()`, `PetscDSCopyEquations()`, `PetscDSSetResidual()`, `PetscDSSetJacobian()`, `PetscDSSetRiemannSolver()`, `PetscDSSetBdResidual()`, `PetscDSSetBdJacobian()`, `PetscDSCreate()`
403345480ffeSMatthew G. Knepley @*/
4034d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDSCopyExactSolutions(PetscDS ds, PetscDS newds)
4035d71ae5a4SJacob Faibussowitsch {
403645480ffeSMatthew G. Knepley   PetscSimplePointFunc sol;
403745480ffeSMatthew G. Knepley   void                *ctx;
403845480ffeSMatthew G. Knepley   PetscInt             Nf, f;
403945480ffeSMatthew G. Knepley 
404045480ffeSMatthew G. Knepley   PetscFunctionBegin;
404145480ffeSMatthew G. Knepley   PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1);
404245480ffeSMatthew G. Knepley   PetscValidHeaderSpecific(newds, PETSCDS_CLASSID, 2);
40439566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(ds, &Nf));
404445480ffeSMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
40459566063dSJacob Faibussowitsch     PetscCall(PetscDSGetExactSolution(ds, f, &sol, &ctx));
40469566063dSJacob Faibussowitsch     PetscCall(PetscDSSetExactSolution(newds, f, sol, ctx));
40479566063dSJacob Faibussowitsch     PetscCall(PetscDSGetExactSolutionTimeDerivative(ds, f, &sol, &ctx));
40489566063dSJacob Faibussowitsch     PetscCall(PetscDSSetExactSolutionTimeDerivative(newds, f, sol, ctx));
404945480ffeSMatthew G. Knepley   }
40503ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
405145480ffeSMatthew G. Knepley }
405245480ffeSMatthew G. Knepley 
405307218a29SMatthew G. Knepley PetscErrorCode PetscDSCopy(PetscDS ds, DM dmNew, PetscDS dsNew)
405407218a29SMatthew G. Knepley {
405507218a29SMatthew G. Knepley   DSBoundary b;
405607218a29SMatthew G. Knepley   PetscInt   cdim, Nf, f, d;
405707218a29SMatthew G. Knepley   PetscBool  isCohesive;
405807218a29SMatthew G. Knepley   void      *ctx;
405907218a29SMatthew G. Knepley 
406007218a29SMatthew G. Knepley   PetscFunctionBegin;
406107218a29SMatthew G. Knepley   PetscCall(PetscDSCopyConstants(ds, dsNew));
406207218a29SMatthew G. Knepley   PetscCall(PetscDSCopyExactSolutions(ds, dsNew));
406307218a29SMatthew G. Knepley   PetscCall(PetscDSSelectDiscretizations(ds, PETSC_DETERMINE, NULL, dsNew));
406407218a29SMatthew G. Knepley   PetscCall(PetscDSCopyEquations(ds, dsNew));
406507218a29SMatthew G. Knepley   PetscCall(PetscDSGetNumFields(ds, &Nf));
406607218a29SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
406707218a29SMatthew G. Knepley     PetscCall(PetscDSGetContext(ds, f, &ctx));
406807218a29SMatthew G. Knepley     PetscCall(PetscDSSetContext(dsNew, f, ctx));
406907218a29SMatthew G. Knepley     PetscCall(PetscDSGetCohesive(ds, f, &isCohesive));
407007218a29SMatthew G. Knepley     PetscCall(PetscDSSetCohesive(dsNew, f, isCohesive));
407107218a29SMatthew G. Knepley     PetscCall(PetscDSGetJetDegree(ds, f, &d));
407207218a29SMatthew G. Knepley     PetscCall(PetscDSSetJetDegree(dsNew, f, d));
407307218a29SMatthew G. Knepley   }
407407218a29SMatthew G. Knepley   if (Nf) {
407507218a29SMatthew G. Knepley     PetscCall(PetscDSGetCoordinateDimension(ds, &cdim));
407607218a29SMatthew G. Knepley     PetscCall(PetscDSSetCoordinateDimension(dsNew, cdim));
407707218a29SMatthew G. Knepley   }
407807218a29SMatthew G. Knepley   PetscCall(PetscDSCopyBoundary(ds, PETSC_DETERMINE, NULL, dsNew));
407907218a29SMatthew G. Knepley   for (b = dsNew->boundary; b; b = b->next) {
408007218a29SMatthew G. Knepley     PetscCall(DMGetLabel(dmNew, b->lname, &b->label));
408107218a29SMatthew G. Knepley     /* Do not check if label exists here, since p4est calls this for the reference tree which does not have the labels */
408207218a29SMatthew G. Knepley     //PetscCheck(b->label,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Label %s missing in new DM", name);
408307218a29SMatthew G. Knepley   }
408407218a29SMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
408507218a29SMatthew G. Knepley }
408607218a29SMatthew G. Knepley 
4087d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDSGetHeightSubspace(PetscDS prob, PetscInt height, PetscDS *subprob)
4088d71ae5a4SJacob Faibussowitsch {
4089df3a45bdSMatthew G. Knepley   PetscInt dim, Nf, f;
4090b1353e8eSMatthew G. Knepley 
4091b1353e8eSMatthew G. Knepley   PetscFunctionBegin;
4092b1353e8eSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
40934f572ea9SToby Isaac   PetscAssertPointer(subprob, 3);
40949371c9d4SSatish Balay   if (height == 0) {
40959371c9d4SSatish Balay     *subprob = prob;
40963ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
40979371c9d4SSatish Balay   }
40989566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(prob, &Nf));
40999566063dSJacob Faibussowitsch   PetscCall(PetscDSGetSpatialDimension(prob, &dim));
410063a3b9bcSJacob Faibussowitsch   PetscCheck(height <= dim, PetscObjectComm((PetscObject)prob), PETSC_ERR_ARG_OUTOFRANGE, "DS can only handle height in [0, %" PetscInt_FMT "], not %" PetscInt_FMT, dim, height);
41019566063dSJacob Faibussowitsch   if (!prob->subprobs) PetscCall(PetscCalloc1(dim, &prob->subprobs));
4102df3a45bdSMatthew G. Knepley   if (!prob->subprobs[height - 1]) {
4103b1353e8eSMatthew G. Knepley     PetscInt cdim;
4104b1353e8eSMatthew G. Knepley 
41059566063dSJacob Faibussowitsch     PetscCall(PetscDSCreate(PetscObjectComm((PetscObject)prob), &prob->subprobs[height - 1]));
41069566063dSJacob Faibussowitsch     PetscCall(PetscDSGetCoordinateDimension(prob, &cdim));
41079566063dSJacob Faibussowitsch     PetscCall(PetscDSSetCoordinateDimension(prob->subprobs[height - 1], cdim));
4108b1353e8eSMatthew G. Knepley     for (f = 0; f < Nf; ++f) {
4109b1353e8eSMatthew G. Knepley       PetscFE      subfe;
4110b1353e8eSMatthew G. Knepley       PetscObject  obj;
4111b1353e8eSMatthew G. Knepley       PetscClassId id;
4112b1353e8eSMatthew G. Knepley 
41139566063dSJacob Faibussowitsch       PetscCall(PetscDSGetDiscretization(prob, f, &obj));
41149566063dSJacob Faibussowitsch       PetscCall(PetscObjectGetClassId(obj, &id));
41159566063dSJacob Faibussowitsch       if (id == PETSCFE_CLASSID) PetscCall(PetscFEGetHeightSubspace((PetscFE)obj, height, &subfe));
411663a3b9bcSJacob Faibussowitsch       else SETERRQ(PetscObjectComm((PetscObject)prob), PETSC_ERR_ARG_WRONG, "Unsupported discretization type for field %" PetscInt_FMT, f);
41179566063dSJacob Faibussowitsch       PetscCall(PetscDSSetDiscretization(prob->subprobs[height - 1], f, (PetscObject)subfe));
4118b1353e8eSMatthew G. Knepley     }
4119b1353e8eSMatthew G. Knepley   }
4120df3a45bdSMatthew G. Knepley   *subprob = prob->subprobs[height - 1];
41213ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4122b1353e8eSMatthew G. Knepley }
4123b1353e8eSMatthew G. Knepley 
41244366bac7SMatthew G. Knepley PetscErrorCode PetscDSPermuteQuadPoint(PetscDS ds, PetscInt ornt, PetscInt field, PetscInt q, PetscInt *qperm)
41254366bac7SMatthew G. Knepley {
41264366bac7SMatthew G. Knepley   IS              permIS;
41274366bac7SMatthew G. Knepley   PetscQuadrature quad;
41284366bac7SMatthew G. Knepley   DMPolytopeType  ct;
41294366bac7SMatthew G. Knepley   const PetscInt *perm;
41304366bac7SMatthew G. Knepley   PetscInt        Na, Nq;
41314366bac7SMatthew G. Knepley 
41324366bac7SMatthew G. Knepley   PetscFunctionBeginHot;
41334366bac7SMatthew G. Knepley   PetscCall(PetscFEGetQuadrature((PetscFE)ds->disc[field], &quad));
41344366bac7SMatthew G. Knepley   PetscCall(PetscQuadratureGetData(quad, NULL, NULL, &Nq, NULL, NULL));
41354366bac7SMatthew G. Knepley   PetscCall(PetscQuadratureGetCellType(quad, &ct));
41364366bac7SMatthew G. Knepley   PetscCheck(q >= 0 && q < Nq, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Quadrature point %" PetscInt_FMT " is not in [0, %" PetscInt_FMT ")", q, Nq);
41374366bac7SMatthew G. Knepley   Na = DMPolytopeTypeGetNumArrangments(ct) / 2;
41384366bac7SMatthew G. Knepley   PetscCheck(ornt >= -Na && ornt < Na, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Orientation %" PetscInt_FMT " of %s is not in [%" PetscInt_FMT ", %" PetscInt_FMT ")", ornt, DMPolytopeTypes[ct], -Na, Na);
41394366bac7SMatthew G. Knepley   if (!ds->quadPerm[(PetscInt)ct]) PetscCall(PetscQuadratureComputePermutations(quad, NULL, &ds->quadPerm[(PetscInt)ct]));
41404366bac7SMatthew G. Knepley   permIS = ds->quadPerm[(PetscInt)ct][ornt + Na];
41414366bac7SMatthew G. Knepley   PetscCall(ISGetIndices(permIS, &perm));
41424366bac7SMatthew G. Knepley   *qperm = perm[q];
41434366bac7SMatthew G. Knepley   PetscCall(ISRestoreIndices(permIS, &perm));
41444366bac7SMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
41454366bac7SMatthew G. Knepley }
41464366bac7SMatthew G. Knepley 
4147d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDSGetDiscType_Internal(PetscDS ds, PetscInt f, PetscDiscType *disctype)
4148d71ae5a4SJacob Faibussowitsch {
4149c7bd5f0bSMatthew G. Knepley   PetscObject  obj;
4150c7bd5f0bSMatthew G. Knepley   PetscClassId id;
4151c7bd5f0bSMatthew G. Knepley   PetscInt     Nf;
4152c7bd5f0bSMatthew G. Knepley 
4153c7bd5f0bSMatthew G. Knepley   PetscFunctionBegin;
4154c7bd5f0bSMatthew G. Knepley   PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1);
41554f572ea9SToby Isaac   PetscAssertPointer(disctype, 3);
4156665f567fSMatthew G. Knepley   *disctype = PETSC_DISC_NONE;
41579566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(ds, &Nf));
415863a3b9bcSJacob Faibussowitsch   PetscCheck(f < Nf, PetscObjectComm((PetscObject)ds), PETSC_ERR_ARG_SIZ, "Field %" PetscInt_FMT " must be in [0, %" PetscInt_FMT ")", f, Nf);
41599566063dSJacob Faibussowitsch   PetscCall(PetscDSGetDiscretization(ds, f, &obj));
4160665f567fSMatthew G. Knepley   if (obj) {
41619566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetClassId(obj, &id));
4162665f567fSMatthew G. Knepley     if (id == PETSCFE_CLASSID) *disctype = PETSC_DISC_FE;
4163665f567fSMatthew G. Knepley     else *disctype = PETSC_DISC_FV;
4164665f567fSMatthew G. Knepley   }
41653ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4166c7bd5f0bSMatthew G. Knepley }
4167c7bd5f0bSMatthew G. Knepley 
4168d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscDSDestroy_Basic(PetscDS ds)
4169d71ae5a4SJacob Faibussowitsch {
41702764a2aaSMatthew G. Knepley   PetscFunctionBegin;
41719566063dSJacob Faibussowitsch   PetscCall(PetscFree(ds->data));
41723ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
41732764a2aaSMatthew G. Knepley }
41742764a2aaSMatthew G. Knepley 
4175d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscDSInitialize_Basic(PetscDS ds)
4176d71ae5a4SJacob Faibussowitsch {
41772764a2aaSMatthew G. Knepley   PetscFunctionBegin;
41786528b96dSMatthew G. Knepley   ds->ops->setfromoptions = NULL;
41796528b96dSMatthew G. Knepley   ds->ops->setup          = NULL;
41806528b96dSMatthew G. Knepley   ds->ops->view           = NULL;
41816528b96dSMatthew G. Knepley   ds->ops->destroy        = PetscDSDestroy_Basic;
41823ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
41832764a2aaSMatthew G. Knepley }
41842764a2aaSMatthew G. Knepley 
41852764a2aaSMatthew G. Knepley /*MC
41862764a2aaSMatthew G. Knepley   PETSCDSBASIC = "basic" - A discrete system with pointwise residual and boundary residual functions
41872764a2aaSMatthew G. Knepley 
41882764a2aaSMatthew G. Knepley   Level: intermediate
41892764a2aaSMatthew G. Knepley 
4190db781477SPatrick Sanan .seealso: `PetscDSType`, `PetscDSCreate()`, `PetscDSSetType()`
41912764a2aaSMatthew G. Knepley M*/
41922764a2aaSMatthew G. Knepley 
4193d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscDSCreate_Basic(PetscDS ds)
4194d71ae5a4SJacob Faibussowitsch {
41952764a2aaSMatthew G. Knepley   PetscDS_Basic *b;
41962764a2aaSMatthew G. Knepley 
41972764a2aaSMatthew G. Knepley   PetscFunctionBegin;
41986528b96dSMatthew G. Knepley   PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1);
41994dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&b));
42006528b96dSMatthew G. Knepley   ds->data = b;
42012764a2aaSMatthew G. Knepley 
42029566063dSJacob Faibussowitsch   PetscCall(PetscDSInitialize_Basic(ds));
42033ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
42042764a2aaSMatthew G. Knepley }
4205