xref: /petsc/src/dm/dt/interface/dtds.c (revision 57fc01e9b7c9cd48aeb2bb3add3941f64bfb3b3b)
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 
82764a2aaSMatthew G. Knepley /*@C
92764a2aaSMatthew G. Knepley   PetscDSRegister - Adds a new PetscDS implementation
102764a2aaSMatthew G. Knepley 
112764a2aaSMatthew G. Knepley   Not Collective
122764a2aaSMatthew G. Knepley 
132764a2aaSMatthew G. Knepley   Input Parameters:
142764a2aaSMatthew G. Knepley + name        - The name of a new user-defined creation routine
152764a2aaSMatthew G. Knepley - create_func - The creation routine itself
162764a2aaSMatthew G. Knepley 
172764a2aaSMatthew G. Knepley   Notes:
182764a2aaSMatthew G. Knepley   PetscDSRegister() may be called multiple times to add several user-defined PetscDSs
192764a2aaSMatthew G. Knepley 
202764a2aaSMatthew G. Knepley   Sample usage:
212764a2aaSMatthew G. Knepley .vb
222764a2aaSMatthew G. Knepley     PetscDSRegister("my_ds", MyPetscDSCreate);
232764a2aaSMatthew G. Knepley .ve
242764a2aaSMatthew G. Knepley 
252764a2aaSMatthew G. Knepley   Then, your PetscDS type can be chosen with the procedural interface via
262764a2aaSMatthew G. Knepley .vb
272764a2aaSMatthew G. Knepley     PetscDSCreate(MPI_Comm, PetscDS *);
282764a2aaSMatthew G. Knepley     PetscDSSetType(PetscDS, "my_ds");
292764a2aaSMatthew G. Knepley .ve
302764a2aaSMatthew G. Knepley    or at runtime via the option
312764a2aaSMatthew G. Knepley .vb
322764a2aaSMatthew G. Knepley     -petscds_type my_ds
332764a2aaSMatthew G. Knepley .ve
342764a2aaSMatthew G. Knepley 
352764a2aaSMatthew G. Knepley   Level: advanced
362764a2aaSMatthew G. Knepley 
372764a2aaSMatthew G. Knepley .keywords: PetscDS, register
382764a2aaSMatthew G. Knepley .seealso: PetscDSRegisterAll(), PetscDSRegisterDestroy()
392764a2aaSMatthew G. Knepley 
402764a2aaSMatthew G. Knepley @*/
412764a2aaSMatthew G. Knepley PetscErrorCode PetscDSRegister(const char sname[], PetscErrorCode (*function)(PetscDS))
422764a2aaSMatthew G. Knepley {
432764a2aaSMatthew G. Knepley   PetscErrorCode ierr;
442764a2aaSMatthew G. Knepley 
452764a2aaSMatthew G. Knepley   PetscFunctionBegin;
462764a2aaSMatthew G. Knepley   ierr = PetscFunctionListAdd(&PetscDSList, sname, function);CHKERRQ(ierr);
472764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
482764a2aaSMatthew G. Knepley }
492764a2aaSMatthew G. Knepley 
502764a2aaSMatthew G. Knepley /*@C
512764a2aaSMatthew G. Knepley   PetscDSSetType - Builds a particular PetscDS
522764a2aaSMatthew G. Knepley 
532764a2aaSMatthew G. Knepley   Collective on PetscDS
542764a2aaSMatthew G. Knepley 
552764a2aaSMatthew G. Knepley   Input Parameters:
562764a2aaSMatthew G. Knepley + prob - The PetscDS object
572764a2aaSMatthew G. Knepley - name - The kind of system
582764a2aaSMatthew G. Knepley 
592764a2aaSMatthew G. Knepley   Options Database Key:
602764a2aaSMatthew G. Knepley . -petscds_type <type> - Sets the PetscDS type; use -help for a list of available types
612764a2aaSMatthew G. Knepley 
622764a2aaSMatthew G. Knepley   Level: intermediate
632764a2aaSMatthew G. Knepley 
642764a2aaSMatthew G. Knepley .keywords: PetscDS, set, type
652764a2aaSMatthew G. Knepley .seealso: PetscDSGetType(), PetscDSCreate()
662764a2aaSMatthew G. Knepley @*/
672764a2aaSMatthew G. Knepley PetscErrorCode PetscDSSetType(PetscDS prob, PetscDSType name)
682764a2aaSMatthew G. Knepley {
692764a2aaSMatthew G. Knepley   PetscErrorCode (*r)(PetscDS);
702764a2aaSMatthew G. Knepley   PetscBool      match;
712764a2aaSMatthew G. Knepley   PetscErrorCode ierr;
722764a2aaSMatthew G. Knepley 
732764a2aaSMatthew G. Knepley   PetscFunctionBegin;
742764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
752764a2aaSMatthew G. Knepley   ierr = PetscObjectTypeCompare((PetscObject) prob, name, &match);CHKERRQ(ierr);
762764a2aaSMatthew G. Knepley   if (match) PetscFunctionReturn(0);
772764a2aaSMatthew G. Knepley 
780f51fdf8SToby Isaac   ierr = PetscDSRegisterAll();CHKERRQ(ierr);
792764a2aaSMatthew G. Knepley   ierr = PetscFunctionListFind(PetscDSList, name, &r);CHKERRQ(ierr);
802764a2aaSMatthew G. Knepley   if (!r) SETERRQ1(PetscObjectComm((PetscObject) prob), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown PetscDS type: %s", name);
812764a2aaSMatthew G. Knepley 
822764a2aaSMatthew G. Knepley   if (prob->ops->destroy) {
832764a2aaSMatthew G. Knepley     ierr             = (*prob->ops->destroy)(prob);CHKERRQ(ierr);
842764a2aaSMatthew G. Knepley     prob->ops->destroy = NULL;
852764a2aaSMatthew G. Knepley   }
862764a2aaSMatthew G. Knepley   ierr = (*r)(prob);CHKERRQ(ierr);
872764a2aaSMatthew G. Knepley   ierr = PetscObjectChangeTypeName((PetscObject) prob, name);CHKERRQ(ierr);
882764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
892764a2aaSMatthew G. Knepley }
902764a2aaSMatthew G. Knepley 
912764a2aaSMatthew G. Knepley /*@C
922764a2aaSMatthew G. Knepley   PetscDSGetType - Gets the PetscDS type name (as a string) from the object.
932764a2aaSMatthew G. Knepley 
942764a2aaSMatthew G. Knepley   Not Collective
952764a2aaSMatthew G. Knepley 
962764a2aaSMatthew G. Knepley   Input Parameter:
972764a2aaSMatthew G. Knepley . prob  - The PetscDS
982764a2aaSMatthew G. Knepley 
992764a2aaSMatthew G. Knepley   Output Parameter:
1002764a2aaSMatthew G. Knepley . name - The PetscDS type name
1012764a2aaSMatthew G. Knepley 
1022764a2aaSMatthew G. Knepley   Level: intermediate
1032764a2aaSMatthew G. Knepley 
1042764a2aaSMatthew G. Knepley .keywords: PetscDS, get, type, name
1052764a2aaSMatthew G. Knepley .seealso: PetscDSSetType(), PetscDSCreate()
1062764a2aaSMatthew G. Knepley @*/
1072764a2aaSMatthew G. Knepley PetscErrorCode PetscDSGetType(PetscDS prob, PetscDSType *name)
1082764a2aaSMatthew G. Knepley {
1092764a2aaSMatthew G. Knepley   PetscErrorCode ierr;
1102764a2aaSMatthew G. Knepley 
1112764a2aaSMatthew G. Knepley   PetscFunctionBegin;
1122764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
113c959eef4SJed Brown   PetscValidPointer(name, 2);
1140f51fdf8SToby Isaac   ierr = PetscDSRegisterAll();CHKERRQ(ierr);
1152764a2aaSMatthew G. Knepley   *name = ((PetscObject) prob)->type_name;
1162764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
1172764a2aaSMatthew G. Knepley }
1182764a2aaSMatthew G. Knepley 
1197d8a60eaSMatthew G. Knepley static PetscErrorCode PetscDSView_Ascii(PetscDS prob, PetscViewer viewer)
1207d8a60eaSMatthew G. Knepley {
1217d8a60eaSMatthew G. Knepley   PetscViewerFormat  format;
12297b6e6e8SMatthew G. Knepley   const PetscScalar *constants;
12397b6e6e8SMatthew G. Knepley   PetscInt           numConstants, f;
1247d8a60eaSMatthew G. Knepley   PetscErrorCode     ierr;
1257d8a60eaSMatthew G. Knepley 
1267d8a60eaSMatthew G. Knepley   PetscFunctionBegin;
1277d8a60eaSMatthew G. Knepley   ierr = PetscViewerGetFormat(viewer, &format);CHKERRQ(ierr);
1287d8a60eaSMatthew G. Knepley   ierr = PetscViewerASCIIPrintf(viewer, "Discrete System with %d fields\n", prob->Nf);CHKERRQ(ierr);
1297d8a60eaSMatthew G. Knepley   ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
1307d8a60eaSMatthew G. Knepley   for (f = 0; f < prob->Nf; ++f) {
1317d8a60eaSMatthew G. Knepley     PetscObject  obj;
1327d8a60eaSMatthew G. Knepley     PetscClassId id;
1337d8a60eaSMatthew G. Knepley     const char  *name;
1347d8a60eaSMatthew G. Knepley     PetscInt     Nc;
1357d8a60eaSMatthew G. Knepley 
1367d8a60eaSMatthew G. Knepley     ierr = PetscDSGetDiscretization(prob, f, &obj);CHKERRQ(ierr);
1377d8a60eaSMatthew G. Knepley     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
1387d8a60eaSMatthew G. Knepley     ierr = PetscObjectGetName(obj, &name);CHKERRQ(ierr);
1397d8a60eaSMatthew G. Knepley     ierr = PetscViewerASCIIPrintf(viewer, "Field %s", name ? name : "<unknown>");CHKERRQ(ierr);
1407d8a60eaSMatthew G. Knepley     if (id == PETSCFE_CLASSID)      {
1417d8a60eaSMatthew G. Knepley       ierr = PetscFEGetNumComponents((PetscFE) obj, &Nc);CHKERRQ(ierr);
1427d8a60eaSMatthew G. Knepley       ierr = PetscViewerASCIIPrintf(viewer, " FEM");CHKERRQ(ierr);
1437d8a60eaSMatthew G. Knepley     } else if (id == PETSCFV_CLASSID) {
1447d8a60eaSMatthew G. Knepley       ierr = PetscFVGetNumComponents((PetscFV) obj, &Nc);CHKERRQ(ierr);
1457d8a60eaSMatthew G. Knepley       ierr = PetscViewerASCIIPrintf(viewer, " FVM");CHKERRQ(ierr);
1467d8a60eaSMatthew G. Knepley     }
14797b6e6e8SMatthew G. Knepley     else SETERRQ1(PetscObjectComm((PetscObject) prob), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %D", f);
14897b6e6e8SMatthew G. Knepley     if (Nc > 1) {ierr = PetscViewerASCIIPrintf(viewer, "%D components", Nc);CHKERRQ(ierr);}
14997b6e6e8SMatthew G. Knepley     else        {ierr = PetscViewerASCIIPrintf(viewer, "%D component ", Nc);CHKERRQ(ierr);}
150249df284SMatthew G. Knepley     if (prob->implicit[f]) {ierr = PetscViewerASCIIPrintf(viewer, " (implicit)");CHKERRQ(ierr);}
151249df284SMatthew G. Knepley     else                   {ierr = PetscViewerASCIIPrintf(viewer, " (explicit)");CHKERRQ(ierr);}
152a6cbbb48SMatthew G. Knepley     if (prob->adjacency[f*2+0]) {
153a6cbbb48SMatthew G. Knepley       if (prob->adjacency[f*2+1]) {ierr = PetscViewerASCIIPrintf(viewer, " (adj FVM++)");CHKERRQ(ierr);}
154a6cbbb48SMatthew G. Knepley       else                        {ierr = PetscViewerASCIIPrintf(viewer, " (adj FVM)");CHKERRQ(ierr);}
155a6cbbb48SMatthew G. Knepley     } else {
156a6cbbb48SMatthew G. Knepley       if (prob->adjacency[f*2+1]) {ierr = PetscViewerASCIIPrintf(viewer, " (adj FEM)");CHKERRQ(ierr);}
157a6cbbb48SMatthew G. Knepley       else                        {ierr = PetscViewerASCIIPrintf(viewer, " (adj FUNKY)");CHKERRQ(ierr);}
158a6cbbb48SMatthew G. Knepley     }
1597d8a60eaSMatthew G. Knepley     ierr = PetscViewerASCIIPrintf(viewer, "\n");CHKERRQ(ierr);
1607d8a60eaSMatthew G. Knepley     if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
1617d8a60eaSMatthew G. Knepley       if (id == PETSCFE_CLASSID)      {ierr = PetscFEView((PetscFE) obj, viewer);CHKERRQ(ierr);}
1627d8a60eaSMatthew G. Knepley       else if (id == PETSCFV_CLASSID) {ierr = PetscFVView((PetscFV) obj, viewer);CHKERRQ(ierr);}
1637d8a60eaSMatthew G. Knepley     }
1647d8a60eaSMatthew G. Knepley   }
16597b6e6e8SMatthew G. Knepley   ierr = PetscDSGetConstants(prob, &numConstants, &constants);CHKERRQ(ierr);
16697b6e6e8SMatthew G. Knepley   if (numConstants) {
16797b6e6e8SMatthew G. Knepley     ierr = PetscViewerASCIIPrintf(viewer, "%D constants\n", numConstants);CHKERRQ(ierr);
16897b6e6e8SMatthew G. Knepley     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
169*57fc01e9SMatthew G. Knepley     for (f = 0; f < numConstants; ++f) {ierr = PetscViewerASCIIPrintf(viewer, "%g\n", (double) PetscRealPart(constants[f]));CHKERRQ(ierr);}
17097b6e6e8SMatthew G. Knepley     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
17197b6e6e8SMatthew G. Knepley   }
1727d8a60eaSMatthew G. Knepley   ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1737d8a60eaSMatthew G. Knepley   PetscFunctionReturn(0);
1747d8a60eaSMatthew G. Knepley }
1757d8a60eaSMatthew G. Knepley 
1762764a2aaSMatthew G. Knepley /*@C
1772764a2aaSMatthew G. Knepley   PetscDSView - Views a PetscDS
1782764a2aaSMatthew G. Knepley 
1792764a2aaSMatthew G. Knepley   Collective on PetscDS
1802764a2aaSMatthew G. Knepley 
1812764a2aaSMatthew G. Knepley   Input Parameter:
1822764a2aaSMatthew G. Knepley + prob - the PetscDS object to view
1832764a2aaSMatthew G. Knepley - v  - the viewer
1842764a2aaSMatthew G. Knepley 
1852764a2aaSMatthew G. Knepley   Level: developer
1862764a2aaSMatthew G. Knepley 
1872764a2aaSMatthew G. Knepley .seealso PetscDSDestroy()
1882764a2aaSMatthew G. Knepley @*/
1892764a2aaSMatthew G. Knepley PetscErrorCode PetscDSView(PetscDS prob, PetscViewer v)
1902764a2aaSMatthew G. Knepley {
1917d8a60eaSMatthew G. Knepley   PetscBool      iascii;
1922764a2aaSMatthew G. Knepley   PetscErrorCode ierr;
1932764a2aaSMatthew G. Knepley 
1942764a2aaSMatthew G. Knepley   PetscFunctionBegin;
1952764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1962764a2aaSMatthew G. Knepley   if (!v) {ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject) prob), &v);CHKERRQ(ierr);}
1977d8a60eaSMatthew G. Knepley   else    {PetscValidHeaderSpecific(v, PETSC_VIEWER_CLASSID, 2);}
1987d8a60eaSMatthew G. Knepley   ierr = PetscObjectTypeCompare((PetscObject) v, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
1997d8a60eaSMatthew G. Knepley   if (iascii) {ierr = PetscDSView_Ascii(prob, v);CHKERRQ(ierr);}
2002764a2aaSMatthew G. Knepley   if (prob->ops->view) {ierr = (*prob->ops->view)(prob, v);CHKERRQ(ierr);}
2012764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
2022764a2aaSMatthew G. Knepley }
2032764a2aaSMatthew G. Knepley 
2042764a2aaSMatthew G. Knepley /*@
2052764a2aaSMatthew G. Knepley   PetscDSSetFromOptions - sets parameters in a PetscDS from the options database
2062764a2aaSMatthew G. Knepley 
2072764a2aaSMatthew G. Knepley   Collective on PetscDS
2082764a2aaSMatthew G. Knepley 
2092764a2aaSMatthew G. Knepley   Input Parameter:
2102764a2aaSMatthew G. Knepley . prob - the PetscDS object to set options for
2112764a2aaSMatthew G. Knepley 
2122764a2aaSMatthew G. Knepley   Options Database:
2132764a2aaSMatthew G. Knepley 
2142764a2aaSMatthew G. Knepley   Level: developer
2152764a2aaSMatthew G. Knepley 
2162764a2aaSMatthew G. Knepley .seealso PetscDSView()
2172764a2aaSMatthew G. Knepley @*/
2182764a2aaSMatthew G. Knepley PetscErrorCode PetscDSSetFromOptions(PetscDS prob)
2192764a2aaSMatthew G. Knepley {
220f1fd5e65SToby Isaac   DSBoundary     b;
2212764a2aaSMatthew G. Knepley   const char    *defaultType;
2222764a2aaSMatthew G. Knepley   char           name[256];
2232764a2aaSMatthew G. Knepley   PetscBool      flg;
2242764a2aaSMatthew G. Knepley   PetscErrorCode ierr;
2252764a2aaSMatthew G. Knepley 
2262764a2aaSMatthew G. Knepley   PetscFunctionBegin;
2272764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
2282764a2aaSMatthew G. Knepley   if (!((PetscObject) prob)->type_name) {
2292764a2aaSMatthew G. Knepley     defaultType = PETSCDSBASIC;
2302764a2aaSMatthew G. Knepley   } else {
2312764a2aaSMatthew G. Knepley     defaultType = ((PetscObject) prob)->type_name;
2322764a2aaSMatthew G. Knepley   }
2330f51fdf8SToby Isaac   ierr = PetscDSRegisterAll();CHKERRQ(ierr);
2342764a2aaSMatthew G. Knepley 
2352764a2aaSMatthew G. Knepley   ierr = PetscObjectOptionsBegin((PetscObject) prob);CHKERRQ(ierr);
236f1fd5e65SToby Isaac   for (b = prob->boundary; b; b = b->next) {
237f1fd5e65SToby Isaac     char       optname[1024];
238f1fd5e65SToby Isaac     PetscInt   ids[1024], len = 1024;
239f1fd5e65SToby Isaac     PetscBool  flg;
240f1fd5e65SToby Isaac 
241f1fd5e65SToby Isaac     ierr = PetscSNPrintf(optname, sizeof(optname), "-bc_%s", b->name);CHKERRQ(ierr);
242f1fd5e65SToby Isaac     ierr = PetscMemzero(ids, sizeof(ids));CHKERRQ(ierr);
243f1fd5e65SToby Isaac     ierr = PetscOptionsIntArray(optname, "List of boundary IDs", "", ids, &len, &flg);CHKERRQ(ierr);
244f1fd5e65SToby Isaac     if (flg) {
245f1fd5e65SToby Isaac       b->numids = len;
246f1fd5e65SToby Isaac       ierr = PetscFree(b->ids);CHKERRQ(ierr);
247f1fd5e65SToby Isaac       ierr = PetscMalloc1(len, &b->ids);CHKERRQ(ierr);
248f1fd5e65SToby Isaac       ierr = PetscMemcpy(b->ids, ids, len*sizeof(PetscInt));CHKERRQ(ierr);
249f1fd5e65SToby Isaac     }
250e7b0402cSSander Arens     len = 1024;
251f1fd5e65SToby Isaac     ierr = PetscSNPrintf(optname, sizeof(optname), "-bc_%s_comp", b->name);CHKERRQ(ierr);
252f1fd5e65SToby Isaac     ierr = PetscMemzero(ids, sizeof(ids));CHKERRQ(ierr);
253f1fd5e65SToby Isaac     ierr = PetscOptionsIntArray(optname, "List of boundary field components", "", ids, &len, &flg);CHKERRQ(ierr);
254f1fd5e65SToby Isaac     if (flg) {
255f1fd5e65SToby Isaac       b->numcomps = len;
256f1fd5e65SToby Isaac       ierr = PetscFree(b->comps);CHKERRQ(ierr);
257f1fd5e65SToby Isaac       ierr = PetscMalloc1(len, &b->comps);CHKERRQ(ierr);
258f1fd5e65SToby Isaac       ierr = PetscMemcpy(b->comps, ids, len*sizeof(PetscInt));CHKERRQ(ierr);
259f1fd5e65SToby Isaac     }
260f1fd5e65SToby Isaac   }
2612764a2aaSMatthew G. Knepley   ierr = PetscOptionsFList("-petscds_type", "Discrete System", "PetscDSSetType", PetscDSList, defaultType, name, 256, &flg);CHKERRQ(ierr);
2622764a2aaSMatthew G. Knepley   if (flg) {
2632764a2aaSMatthew G. Knepley     ierr = PetscDSSetType(prob, name);CHKERRQ(ierr);
2642764a2aaSMatthew G. Knepley   } else if (!((PetscObject) prob)->type_name) {
2652764a2aaSMatthew G. Knepley     ierr = PetscDSSetType(prob, defaultType);CHKERRQ(ierr);
2662764a2aaSMatthew G. Knepley   }
2672764a2aaSMatthew G. Knepley   if (prob->ops->setfromoptions) {ierr = (*prob->ops->setfromoptions)(prob);CHKERRQ(ierr);}
2682764a2aaSMatthew G. Knepley   /* process any options handlers added with PetscObjectAddOptionsHandler() */
2690633abcbSJed Brown   ierr = PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject) prob);CHKERRQ(ierr);
2702764a2aaSMatthew G. Knepley   ierr = PetscOptionsEnd();CHKERRQ(ierr);
2712764a2aaSMatthew G. Knepley   ierr = PetscDSViewFromOptions(prob, NULL, "-petscds_view");CHKERRQ(ierr);
2722764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
2732764a2aaSMatthew G. Knepley }
2742764a2aaSMatthew G. Knepley 
2752764a2aaSMatthew G. Knepley /*@C
2762764a2aaSMatthew G. Knepley   PetscDSSetUp - Construct data structures for the PetscDS
2772764a2aaSMatthew G. Knepley 
2782764a2aaSMatthew G. Knepley   Collective on PetscDS
2792764a2aaSMatthew G. Knepley 
2802764a2aaSMatthew G. Knepley   Input Parameter:
2812764a2aaSMatthew G. Knepley . prob - the PetscDS object to setup
2822764a2aaSMatthew G. Knepley 
2832764a2aaSMatthew G. Knepley   Level: developer
2842764a2aaSMatthew G. Knepley 
2852764a2aaSMatthew G. Knepley .seealso PetscDSView(), PetscDSDestroy()
2862764a2aaSMatthew G. Knepley @*/
2872764a2aaSMatthew G. Knepley PetscErrorCode PetscDSSetUp(PetscDS prob)
2882764a2aaSMatthew G. Knepley {
2892764a2aaSMatthew G. Knepley   const PetscInt Nf = prob->Nf;
2902764a2aaSMatthew G. Knepley   PetscInt       dim, work, NcMax = 0, NqMax = 0, f;
2912764a2aaSMatthew G. Knepley   PetscErrorCode ierr;
2922764a2aaSMatthew G. Knepley 
2932764a2aaSMatthew G. Knepley   PetscFunctionBegin;
2942764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
2952764a2aaSMatthew G. Knepley   if (prob->setup) PetscFunctionReturn(0);
2962764a2aaSMatthew G. Knepley   /* Calculate sizes */
2972764a2aaSMatthew G. Knepley   ierr = PetscDSGetSpatialDimension(prob, &dim);CHKERRQ(ierr);
298f744cafaSSander Arens   prob->totDim = prob->totComp = 0;
29947e57110SSander Arens   ierr = PetscMalloc2(Nf,&prob->Nc,Nf,&prob->Nb);CHKERRQ(ierr);
300f744cafaSSander Arens   ierr = PetscCalloc2(Nf+1,&prob->off,Nf+1,&prob->offDer);CHKERRQ(ierr);
301f744cafaSSander Arens   ierr = PetscMalloc4(Nf,&prob->basis,Nf,&prob->basisDer,Nf,&prob->basisFace,Nf,&prob->basisDerFace);CHKERRQ(ierr);
3022764a2aaSMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
3039de99aefSMatthew G. Knepley     PetscObject     obj;
3049de99aefSMatthew G. Knepley     PetscClassId    id;
3052764a2aaSMatthew G. Knepley     PetscQuadrature q;
3069de99aefSMatthew G. Knepley     PetscInt        Nq = 0, Nb, Nc;
3072764a2aaSMatthew G. Knepley 
3089de99aefSMatthew G. Knepley     ierr = PetscDSGetDiscretization(prob, f, &obj);CHKERRQ(ierr);
3099de99aefSMatthew G. Knepley     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
3109de99aefSMatthew G. Knepley     if (id == PETSCFE_CLASSID)      {
3119de99aefSMatthew G. Knepley       PetscFE fe = (PetscFE) obj;
3129de99aefSMatthew G. Knepley 
3132764a2aaSMatthew G. Knepley       ierr = PetscFEGetQuadrature(fe, &q);CHKERRQ(ierr);
3142764a2aaSMatthew G. Knepley       ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr);
3152764a2aaSMatthew G. Knepley       ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
3162764a2aaSMatthew G. Knepley       ierr = PetscFEGetDefaultTabulation(fe, &prob->basis[f], &prob->basisDer[f], NULL);CHKERRQ(ierr);
3174d0b9603SSander Arens       ierr = PetscFEGetFaceTabulation(fe, &prob->basisFace[f], &prob->basisDerFace[f], NULL);CHKERRQ(ierr);
3189de99aefSMatthew G. Knepley     } else if (id == PETSCFV_CLASSID) {
3199de99aefSMatthew G. Knepley       PetscFV fv = (PetscFV) obj;
3209de99aefSMatthew G. Knepley 
3219de99aefSMatthew G. Knepley       ierr = PetscFVGetQuadrature(fv, &q);CHKERRQ(ierr);
3229de99aefSMatthew G. Knepley       ierr = PetscFVGetNumComponents(fv, &Nc);CHKERRQ(ierr);
3239c3cf19fSMatthew G. Knepley       Nb   = Nc;
3246c1a3d01SMatthew G. Knepley       ierr = PetscFVGetDefaultTabulation(fv, &prob->basis[f], &prob->basisDer[f], NULL);CHKERRQ(ierr);
3254d0b9603SSander Arens       /* TODO: should PetscFV also have face tabulation? Otherwise there will be a null pointer in prob->basisFace */
326abac5ca0SMatthew G. Knepley     } else SETERRQ1(PetscObjectComm((PetscObject) prob), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", f);
32747e57110SSander Arens     prob->Nc[f]       = Nc;
32847e57110SSander Arens     prob->Nb[f]       = Nb;
329194d53e6SMatthew G. Knepley     prob->off[f+1]    = Nc     + prob->off[f];
330194d53e6SMatthew G. Knepley     prob->offDer[f+1] = Nc*dim + prob->offDer[f];
331a6b92713SMatthew G. Knepley     if (q) {ierr = PetscQuadratureGetData(q, NULL, NULL, &Nq, NULL, NULL);CHKERRQ(ierr);}
3322764a2aaSMatthew G. Knepley     NqMax          = PetscMax(NqMax, Nq);
3332764a2aaSMatthew G. Knepley     NcMax          = PetscMax(NcMax, Nc);
3349c3cf19fSMatthew G. Knepley     prob->totDim  += Nb;
3352764a2aaSMatthew G. Knepley     prob->totComp += Nc;
3362764a2aaSMatthew G. Knepley   }
3372764a2aaSMatthew G. Knepley   work = PetscMax(prob->totComp*dim, PetscSqr(NcMax*dim));
3382764a2aaSMatthew G. Knepley   /* Allocate works space */
3392764a2aaSMatthew G. Knepley   ierr = PetscMalloc5(prob->totComp,&prob->u,prob->totComp,&prob->u_t,prob->totComp*dim,&prob->u_x,dim,&prob->x,work,&prob->refSpaceDer);CHKERRQ(ierr);
3402764a2aaSMatthew G. Knepley   ierr = PetscMalloc6(NqMax*NcMax,&prob->f0,NqMax*NcMax*dim,&prob->f1,NqMax*NcMax*NcMax,&prob->g0,NqMax*NcMax*NcMax*dim,&prob->g1,NqMax*NcMax*NcMax*dim,&prob->g2,NqMax*NcMax*NcMax*dim*dim,&prob->g3);CHKERRQ(ierr);
3412764a2aaSMatthew G. Knepley   if (prob->ops->setup) {ierr = (*prob->ops->setup)(prob);CHKERRQ(ierr);}
3422764a2aaSMatthew G. Knepley   prob->setup = PETSC_TRUE;
3432764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
3442764a2aaSMatthew G. Knepley }
3452764a2aaSMatthew G. Knepley 
3462764a2aaSMatthew G. Knepley static PetscErrorCode PetscDSDestroyStructs_Static(PetscDS prob)
3472764a2aaSMatthew G. Knepley {
3482764a2aaSMatthew G. Knepley   PetscErrorCode ierr;
3492764a2aaSMatthew G. Knepley 
3502764a2aaSMatthew G. Knepley   PetscFunctionBegin;
35147e57110SSander Arens   ierr = PetscFree2(prob->Nc,prob->Nb);CHKERRQ(ierr);
352f744cafaSSander Arens   ierr = PetscFree2(prob->off,prob->offDer);CHKERRQ(ierr);
353f744cafaSSander Arens   ierr = PetscFree4(prob->basis,prob->basisDer,prob->basisFace,prob->basisDerFace);CHKERRQ(ierr);
3542764a2aaSMatthew G. Knepley   ierr = PetscFree5(prob->u,prob->u_t,prob->u_x,prob->x,prob->refSpaceDer);CHKERRQ(ierr);
3552764a2aaSMatthew G. Knepley   ierr = PetscFree6(prob->f0,prob->f1,prob->g0,prob->g1,prob->g2,prob->g3);CHKERRQ(ierr);
3562764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
3572764a2aaSMatthew G. Knepley }
3582764a2aaSMatthew G. Knepley 
3592764a2aaSMatthew G. Knepley static PetscErrorCode PetscDSEnlarge_Static(PetscDS prob, PetscInt NfNew)
3602764a2aaSMatthew G. Knepley {
361f744cafaSSander Arens   PetscObject      *tmpd;
362a6cbbb48SMatthew G. Knepley   PetscBool        *tmpi, *tmpa;
36332d2bbc9SMatthew G. Knepley   PetscPointFunc   *tmpobj, *tmpf, *tmpup;
364b7e05686SMatthew G. Knepley   PetscPointJac    *tmpg, *tmpgp, *tmpgt;
3652aa1fc23SMatthew G. Knepley   PetscBdPointFunc *tmpfbd;
3662aa1fc23SMatthew G. Knepley   PetscBdPointJac  *tmpgbd;
367194d53e6SMatthew G. Knepley   PetscRiemannFunc *tmpr;
3680c2f2876SMatthew G. Knepley   void            **tmpctx;
369a6cbbb48SMatthew G. Knepley   PetscInt          Nf = prob->Nf, f, i;
3702764a2aaSMatthew G. Knepley   PetscErrorCode    ierr;
3712764a2aaSMatthew G. Knepley 
3722764a2aaSMatthew G. Knepley   PetscFunctionBegin;
3732764a2aaSMatthew G. Knepley   if (Nf >= NfNew) PetscFunctionReturn(0);
3742764a2aaSMatthew G. Knepley   prob->setup = PETSC_FALSE;
3752764a2aaSMatthew G. Knepley   ierr = PetscDSDestroyStructs_Static(prob);CHKERRQ(ierr);
376f744cafaSSander Arens   ierr = PetscMalloc3(NfNew, &tmpd, NfNew, &tmpi, NfNew*2, &tmpa);CHKERRQ(ierr);
377f744cafaSSander Arens   for (f = 0; f < Nf; ++f) {tmpd[f] = prob->disc[f]; tmpi[f] = prob->implicit[f]; for (i = 0; i < 2; ++i) tmpa[f*2+i] = prob->adjacency[f*2+i];}
378f744cafaSSander Arens   for (f = Nf; f < NfNew; ++f) {tmpd[f] = NULL; tmpi[f] = PETSC_TRUE; tmpa[f*2+0] = PETSC_FALSE; tmpa[f*2+1] = PETSC_TRUE;}
379f744cafaSSander Arens   ierr = PetscFree3(prob->disc, prob->implicit, prob->adjacency);CHKERRQ(ierr);
3802764a2aaSMatthew G. Knepley   prob->Nf        = NfNew;
3812764a2aaSMatthew G. Knepley   prob->disc      = tmpd;
382249df284SMatthew G. Knepley   prob->implicit  = tmpi;
383a6cbbb48SMatthew G. Knepley   prob->adjacency = tmpa;
384b7e05686SMatthew G. Knepley   ierr = PetscCalloc7(NfNew, &tmpobj, NfNew*2, &tmpf, NfNew*NfNew*4, &tmpg, NfNew*NfNew*4, &tmpgp, NfNew*NfNew*4, &tmpgt, NfNew, &tmpr, NfNew, &tmpctx);CHKERRQ(ierr);
38532d2bbc9SMatthew G. Knepley   ierr = PetscCalloc1(NfNew, &tmpup);CHKERRQ(ierr);
3862764a2aaSMatthew G. Knepley   for (f = 0; f < Nf; ++f) tmpobj[f] = prob->obj[f];
3872764a2aaSMatthew G. Knepley   for (f = 0; f < Nf*2; ++f) tmpf[f] = prob->f[f];
3882764a2aaSMatthew G. Knepley   for (f = 0; f < Nf*Nf*4; ++f) tmpg[f] = prob->g[f];
389475e0ac9SMatthew G. Knepley   for (f = 0; f < Nf*Nf*4; ++f) tmpgp[f] = prob->gp[f];
3900c2f2876SMatthew G. Knepley   for (f = 0; f < Nf; ++f) tmpr[f] = prob->r[f];
39132d2bbc9SMatthew G. Knepley   for (f = 0; f < Nf; ++f) tmpup[f] = prob->update[f];
3920c2f2876SMatthew G. Knepley   for (f = 0; f < Nf; ++f) tmpctx[f] = prob->ctx[f];
3932764a2aaSMatthew G. Knepley   for (f = Nf; f < NfNew; ++f) tmpobj[f] = NULL;
3942764a2aaSMatthew G. Knepley   for (f = Nf*2; f < NfNew*2; ++f) tmpf[f] = NULL;
3952764a2aaSMatthew G. Knepley   for (f = Nf*Nf*4; f < NfNew*NfNew*4; ++f) tmpg[f] = NULL;
396475e0ac9SMatthew G. Knepley   for (f = Nf*Nf*4; f < NfNew*NfNew*4; ++f) tmpgp[f] = NULL;
397b7e05686SMatthew G. Knepley   for (f = Nf*Nf*4; f < NfNew*NfNew*4; ++f) tmpgt[f] = NULL;
3980c2f2876SMatthew G. Knepley   for (f = Nf; f < NfNew; ++f) tmpr[f] = NULL;
39932d2bbc9SMatthew G. Knepley   for (f = Nf; f < NfNew; ++f) tmpup[f] = NULL;
4000c2f2876SMatthew G. Knepley   for (f = Nf; f < NfNew; ++f) tmpctx[f] = NULL;
401b7e05686SMatthew G. Knepley   ierr = PetscFree7(prob->obj, prob->f, prob->g, prob->gp, prob->gt, prob->r, prob->ctx);CHKERRQ(ierr);
40232d2bbc9SMatthew G. Knepley   ierr = PetscFree(prob->update);CHKERRQ(ierr);
4032764a2aaSMatthew G. Knepley   prob->obj = tmpobj;
4042764a2aaSMatthew G. Knepley   prob->f   = tmpf;
4052764a2aaSMatthew G. Knepley   prob->g   = tmpg;
406475e0ac9SMatthew G. Knepley   prob->gp  = tmpgp;
407b7e05686SMatthew G. Knepley   prob->gt  = tmpgt;
4080c2f2876SMatthew G. Knepley   prob->r   = tmpr;
40932d2bbc9SMatthew G. Knepley   prob->update = tmpup;
4100c2f2876SMatthew G. Knepley   prob->ctx = tmpctx;
4112764a2aaSMatthew G. Knepley   ierr = PetscCalloc2(NfNew*2, &tmpfbd, NfNew*NfNew*4, &tmpgbd);CHKERRQ(ierr);
4122764a2aaSMatthew G. Knepley   for (f = 0; f < Nf*2; ++f) tmpfbd[f] = prob->fBd[f];
4132764a2aaSMatthew G. Knepley   for (f = 0; f < Nf*Nf*4; ++f) tmpgbd[f] = prob->gBd[f];
4142764a2aaSMatthew G. Knepley   for (f = Nf*2; f < NfNew*2; ++f) tmpfbd[f] = NULL;
4152764a2aaSMatthew G. Knepley   for (f = Nf*Nf*4; f < NfNew*NfNew*4; ++f) tmpgbd[f] = NULL;
4162764a2aaSMatthew G. Knepley   ierr = PetscFree2(prob->fBd, prob->gBd);CHKERRQ(ierr);
4172764a2aaSMatthew G. Knepley   prob->fBd = tmpfbd;
4182764a2aaSMatthew G. Knepley   prob->gBd = tmpgbd;
4192764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
4202764a2aaSMatthew G. Knepley }
4212764a2aaSMatthew G. Knepley 
4222764a2aaSMatthew G. Knepley /*@
4232764a2aaSMatthew G. Knepley   PetscDSDestroy - Destroys a PetscDS object
4242764a2aaSMatthew G. Knepley 
4252764a2aaSMatthew G. Knepley   Collective on PetscDS
4262764a2aaSMatthew G. Knepley 
4272764a2aaSMatthew G. Knepley   Input Parameter:
4282764a2aaSMatthew G. Knepley . prob - the PetscDS object to destroy
4292764a2aaSMatthew G. Knepley 
4302764a2aaSMatthew G. Knepley   Level: developer
4312764a2aaSMatthew G. Knepley 
4322764a2aaSMatthew G. Knepley .seealso PetscDSView()
4332764a2aaSMatthew G. Knepley @*/
4342764a2aaSMatthew G. Knepley PetscErrorCode PetscDSDestroy(PetscDS *prob)
4352764a2aaSMatthew G. Knepley {
4362764a2aaSMatthew G. Knepley   PetscInt       f;
43758ebd649SToby Isaac   DSBoundary     next;
4382764a2aaSMatthew G. Knepley   PetscErrorCode ierr;
4392764a2aaSMatthew G. Knepley 
4402764a2aaSMatthew G. Knepley   PetscFunctionBegin;
4412764a2aaSMatthew G. Knepley   if (!*prob) PetscFunctionReturn(0);
4422764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific((*prob), PETSCDS_CLASSID, 1);
4432764a2aaSMatthew G. Knepley 
4442764a2aaSMatthew G. Knepley   if (--((PetscObject)(*prob))->refct > 0) {*prob = 0; PetscFunctionReturn(0);}
4452764a2aaSMatthew G. Knepley   ((PetscObject) (*prob))->refct = 0;
4462764a2aaSMatthew G. Knepley   ierr = PetscDSDestroyStructs_Static(*prob);CHKERRQ(ierr);
4472764a2aaSMatthew G. Knepley   for (f = 0; f < (*prob)->Nf; ++f) {
4482764a2aaSMatthew G. Knepley     ierr = PetscObjectDereference((*prob)->disc[f]);CHKERRQ(ierr);
4492764a2aaSMatthew G. Knepley   }
450f744cafaSSander Arens   ierr = PetscFree3((*prob)->disc, (*prob)->implicit, (*prob)->adjacency);CHKERRQ(ierr);
451b7e05686SMatthew G. Knepley   ierr = PetscFree7((*prob)->obj,(*prob)->f,(*prob)->g,(*prob)->gp,(*prob)->gt,(*prob)->r,(*prob)->ctx);CHKERRQ(ierr);
45232d2bbc9SMatthew G. Knepley   ierr = PetscFree((*prob)->update);CHKERRQ(ierr);
4532764a2aaSMatthew G. Knepley   ierr = PetscFree2((*prob)->fBd,(*prob)->gBd);CHKERRQ(ierr);
4542764a2aaSMatthew G. Knepley   if ((*prob)->ops->destroy) {ierr = (*(*prob)->ops->destroy)(*prob);CHKERRQ(ierr);}
45558ebd649SToby Isaac   next = (*prob)->boundary;
45658ebd649SToby Isaac   while (next) {
45758ebd649SToby Isaac     DSBoundary b = next;
45858ebd649SToby Isaac 
45958ebd649SToby Isaac     next = b->next;
46058ebd649SToby Isaac     ierr = PetscFree(b->comps);CHKERRQ(ierr);
46158ebd649SToby Isaac     ierr = PetscFree(b->ids);CHKERRQ(ierr);
46258ebd649SToby Isaac     ierr = PetscFree(b->name);CHKERRQ(ierr);
46358ebd649SToby Isaac     ierr = PetscFree(b->labelname);CHKERRQ(ierr);
46458ebd649SToby Isaac     ierr = PetscFree(b);CHKERRQ(ierr);
46558ebd649SToby Isaac   }
46697b6e6e8SMatthew G. Knepley   ierr = PetscFree((*prob)->constants);CHKERRQ(ierr);
4672764a2aaSMatthew G. Knepley   ierr = PetscHeaderDestroy(prob);CHKERRQ(ierr);
4682764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
4692764a2aaSMatthew G. Knepley }
4702764a2aaSMatthew G. Knepley 
4712764a2aaSMatthew G. Knepley /*@
4722764a2aaSMatthew G. Knepley   PetscDSCreate - Creates an empty PetscDS object. The type can then be set with PetscDSSetType().
4732764a2aaSMatthew G. Knepley 
4742764a2aaSMatthew G. Knepley   Collective on MPI_Comm
4752764a2aaSMatthew G. Knepley 
4762764a2aaSMatthew G. Knepley   Input Parameter:
4772764a2aaSMatthew G. Knepley . comm - The communicator for the PetscDS object
4782764a2aaSMatthew G. Knepley 
4792764a2aaSMatthew G. Knepley   Output Parameter:
4802764a2aaSMatthew G. Knepley . prob - The PetscDS object
4812764a2aaSMatthew G. Knepley 
4822764a2aaSMatthew G. Knepley   Level: beginner
4832764a2aaSMatthew G. Knepley 
4842764a2aaSMatthew G. Knepley .seealso: PetscDSSetType(), PETSCDSBASIC
4852764a2aaSMatthew G. Knepley @*/
4862764a2aaSMatthew G. Knepley PetscErrorCode PetscDSCreate(MPI_Comm comm, PetscDS *prob)
4872764a2aaSMatthew G. Knepley {
4882764a2aaSMatthew G. Knepley   PetscDS   p;
4892764a2aaSMatthew G. Knepley   PetscErrorCode ierr;
4902764a2aaSMatthew G. Knepley 
4912764a2aaSMatthew G. Knepley   PetscFunctionBegin;
4922764a2aaSMatthew G. Knepley   PetscValidPointer(prob, 2);
4932764a2aaSMatthew G. Knepley   *prob  = NULL;
4942764a2aaSMatthew G. Knepley   ierr = PetscDSInitializePackage();CHKERRQ(ierr);
4952764a2aaSMatthew G. Knepley 
49673107ff1SLisandro Dalcin   ierr = PetscHeaderCreate(p, PETSCDS_CLASSID, "PetscDS", "Discrete System", "PetscDS", comm, PetscDSDestroy, PetscDSView);CHKERRQ(ierr);
4972764a2aaSMatthew G. Knepley 
4982764a2aaSMatthew G. Knepley   p->Nf    = 0;
4992764a2aaSMatthew G. Knepley   p->setup = PETSC_FALSE;
50097b6e6e8SMatthew G. Knepley   p->numConstants = 0;
50197b6e6e8SMatthew G. Knepley   p->constants    = NULL;
5022764a2aaSMatthew G. Knepley 
5032764a2aaSMatthew G. Knepley   *prob = p;
5042764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
5052764a2aaSMatthew G. Knepley }
5062764a2aaSMatthew G. Knepley 
507bc4ae4beSMatthew G. Knepley /*@
508bc4ae4beSMatthew G. Knepley   PetscDSGetNumFields - Returns the number of fields in the DS
509bc4ae4beSMatthew G. Knepley 
510bc4ae4beSMatthew G. Knepley   Not collective
511bc4ae4beSMatthew G. Knepley 
512bc4ae4beSMatthew G. Knepley   Input Parameter:
513bc4ae4beSMatthew G. Knepley . prob - The PetscDS object
514bc4ae4beSMatthew G. Knepley 
515bc4ae4beSMatthew G. Knepley   Output Parameter:
516bc4ae4beSMatthew G. Knepley . Nf - The number of fields
517bc4ae4beSMatthew G. Knepley 
518bc4ae4beSMatthew G. Knepley   Level: beginner
519bc4ae4beSMatthew G. Knepley 
520bc4ae4beSMatthew G. Knepley .seealso: PetscDSGetSpatialDimension(), PetscDSCreate()
521bc4ae4beSMatthew G. Knepley @*/
5222764a2aaSMatthew G. Knepley PetscErrorCode PetscDSGetNumFields(PetscDS prob, PetscInt *Nf)
5232764a2aaSMatthew G. Knepley {
5242764a2aaSMatthew G. Knepley   PetscFunctionBegin;
5252764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
5262764a2aaSMatthew G. Knepley   PetscValidPointer(Nf, 2);
5272764a2aaSMatthew G. Knepley   *Nf = prob->Nf;
5282764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
5292764a2aaSMatthew G. Knepley }
5302764a2aaSMatthew G. Knepley 
531bc4ae4beSMatthew G. Knepley /*@
532bc4ae4beSMatthew G. Knepley   PetscDSGetSpatialDimension - Returns the spatial dimension of the DS
533bc4ae4beSMatthew G. Knepley 
534bc4ae4beSMatthew G. Knepley   Not collective
535bc4ae4beSMatthew G. Knepley 
536bc4ae4beSMatthew G. Knepley   Input Parameter:
537bc4ae4beSMatthew G. Knepley . prob - The PetscDS object
538bc4ae4beSMatthew G. Knepley 
539bc4ae4beSMatthew G. Knepley   Output Parameter:
540bc4ae4beSMatthew G. Knepley . dim - The spatial dimension
541bc4ae4beSMatthew G. Knepley 
542bc4ae4beSMatthew G. Knepley   Level: beginner
543bc4ae4beSMatthew G. Knepley 
544bc4ae4beSMatthew G. Knepley .seealso: PetscDSGetNumFields(), PetscDSCreate()
545bc4ae4beSMatthew G. Knepley @*/
5462764a2aaSMatthew G. Knepley PetscErrorCode PetscDSGetSpatialDimension(PetscDS prob, PetscInt *dim)
5472764a2aaSMatthew G. Knepley {
5482764a2aaSMatthew G. Knepley   PetscErrorCode ierr;
5492764a2aaSMatthew G. Knepley 
5502764a2aaSMatthew G. Knepley   PetscFunctionBegin;
5512764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
5522764a2aaSMatthew G. Knepley   PetscValidPointer(dim, 2);
5532764a2aaSMatthew G. Knepley   *dim = 0;
5549de99aefSMatthew G. Knepley   if (prob->Nf) {
5559de99aefSMatthew G. Knepley     PetscObject  obj;
5569de99aefSMatthew G. Knepley     PetscClassId id;
5579de99aefSMatthew G. Knepley 
5589de99aefSMatthew G. Knepley     ierr = PetscDSGetDiscretization(prob, 0, &obj);CHKERRQ(ierr);
5599de99aefSMatthew G. Knepley     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
5609de99aefSMatthew G. Knepley     if (id == PETSCFE_CLASSID)      {ierr = PetscFEGetSpatialDimension((PetscFE) obj, dim);CHKERRQ(ierr);}
5619de99aefSMatthew G. Knepley     else if (id == PETSCFV_CLASSID) {ierr = PetscFVGetSpatialDimension((PetscFV) obj, dim);CHKERRQ(ierr);}
5629de99aefSMatthew G. Knepley     else SETERRQ1(PetscObjectComm((PetscObject) prob), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", 0);
5639de99aefSMatthew G. Knepley   }
5642764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
5652764a2aaSMatthew G. Knepley }
5662764a2aaSMatthew G. Knepley 
567bc4ae4beSMatthew G. Knepley /*@
568bc4ae4beSMatthew G. Knepley   PetscDSGetTotalDimension - Returns the total size of the approximation space for this system
569bc4ae4beSMatthew G. Knepley 
570bc4ae4beSMatthew G. Knepley   Not collective
571bc4ae4beSMatthew G. Knepley 
572bc4ae4beSMatthew G. Knepley   Input Parameter:
573bc4ae4beSMatthew G. Knepley . prob - The PetscDS object
574bc4ae4beSMatthew G. Knepley 
575bc4ae4beSMatthew G. Knepley   Output Parameter:
576bc4ae4beSMatthew G. Knepley . dim - The total problem dimension
577bc4ae4beSMatthew G. Knepley 
578bc4ae4beSMatthew G. Knepley   Level: beginner
579bc4ae4beSMatthew G. Knepley 
580bc4ae4beSMatthew G. Knepley .seealso: PetscDSGetNumFields(), PetscDSCreate()
581bc4ae4beSMatthew G. Knepley @*/
5822764a2aaSMatthew G. Knepley PetscErrorCode PetscDSGetTotalDimension(PetscDS prob, PetscInt *dim)
5832764a2aaSMatthew G. Knepley {
5842764a2aaSMatthew G. Knepley   PetscErrorCode ierr;
5852764a2aaSMatthew G. Knepley 
5862764a2aaSMatthew G. Knepley   PetscFunctionBegin;
5872764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
5882764a2aaSMatthew G. Knepley   ierr = PetscDSSetUp(prob);CHKERRQ(ierr);
5892764a2aaSMatthew G. Knepley   PetscValidPointer(dim, 2);
5902764a2aaSMatthew G. Knepley   *dim = prob->totDim;
5912764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
5922764a2aaSMatthew G. Knepley }
5932764a2aaSMatthew G. Knepley 
594bc4ae4beSMatthew G. Knepley /*@
595bc4ae4beSMatthew G. Knepley   PetscDSGetTotalComponents - Returns the total number of components in this system
596bc4ae4beSMatthew G. Knepley 
597bc4ae4beSMatthew G. Knepley   Not collective
598bc4ae4beSMatthew G. Knepley 
599bc4ae4beSMatthew G. Knepley   Input Parameter:
600bc4ae4beSMatthew G. Knepley . prob - The PetscDS object
601bc4ae4beSMatthew G. Knepley 
602bc4ae4beSMatthew G. Knepley   Output Parameter:
603bc4ae4beSMatthew G. Knepley . dim - The total number of components
604bc4ae4beSMatthew G. Knepley 
605bc4ae4beSMatthew G. Knepley   Level: beginner
606bc4ae4beSMatthew G. Knepley 
607bc4ae4beSMatthew G. Knepley .seealso: PetscDSGetNumFields(), PetscDSCreate()
608bc4ae4beSMatthew G. Knepley @*/
6092764a2aaSMatthew G. Knepley PetscErrorCode PetscDSGetTotalComponents(PetscDS prob, PetscInt *Nc)
6102764a2aaSMatthew G. Knepley {
6112764a2aaSMatthew G. Knepley   PetscErrorCode ierr;
6122764a2aaSMatthew G. Knepley 
6132764a2aaSMatthew G. Knepley   PetscFunctionBegin;
6142764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
6152764a2aaSMatthew G. Knepley   ierr = PetscDSSetUp(prob);CHKERRQ(ierr);
6162764a2aaSMatthew G. Knepley   PetscValidPointer(Nc, 2);
6172764a2aaSMatthew G. Knepley   *Nc = prob->totComp;
6182764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
6192764a2aaSMatthew G. Knepley }
6202764a2aaSMatthew G. Knepley 
621bc4ae4beSMatthew G. Knepley /*@
622bc4ae4beSMatthew G. Knepley   PetscDSGetDiscretization - Returns the discretization object for the given field
623bc4ae4beSMatthew G. Knepley 
624bc4ae4beSMatthew G. Knepley   Not collective
625bc4ae4beSMatthew G. Knepley 
626bc4ae4beSMatthew G. Knepley   Input Parameters:
627bc4ae4beSMatthew G. Knepley + prob - The PetscDS object
628bc4ae4beSMatthew G. Knepley - f - The field number
629bc4ae4beSMatthew G. Knepley 
630bc4ae4beSMatthew G. Knepley   Output Parameter:
631bc4ae4beSMatthew G. Knepley . disc - The discretization object
632bc4ae4beSMatthew G. Knepley 
633bc4ae4beSMatthew G. Knepley   Level: beginner
634bc4ae4beSMatthew G. Knepley 
635f744cafaSSander Arens .seealso: PetscDSSetDiscretization(), PetscDSAddDiscretization(), PetscDSGetNumFields(), PetscDSCreate()
636bc4ae4beSMatthew G. Knepley @*/
6372764a2aaSMatthew G. Knepley PetscErrorCode PetscDSGetDiscretization(PetscDS prob, PetscInt f, PetscObject *disc)
6382764a2aaSMatthew G. Knepley {
6392764a2aaSMatthew G. Knepley   PetscFunctionBegin;
6402764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
6412764a2aaSMatthew G. Knepley   PetscValidPointer(disc, 3);
6422764a2aaSMatthew G. Knepley   if ((f < 0) || (f >= prob->Nf)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be in [0, %d)", f, prob->Nf);
6432764a2aaSMatthew G. Knepley   *disc = prob->disc[f];
6442764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
6452764a2aaSMatthew G. Knepley }
6462764a2aaSMatthew G. Knepley 
647bc4ae4beSMatthew G. Knepley /*@
648bc4ae4beSMatthew G. Knepley   PetscDSSetDiscretization - Sets the discretization object for the given field
649bc4ae4beSMatthew G. Knepley 
650bc4ae4beSMatthew G. Knepley   Not collective
651bc4ae4beSMatthew G. Knepley 
652bc4ae4beSMatthew G. Knepley   Input Parameters:
653bc4ae4beSMatthew G. Knepley + prob - The PetscDS object
654bc4ae4beSMatthew G. Knepley . f - The field number
655bc4ae4beSMatthew G. Knepley - disc - The discretization object
656bc4ae4beSMatthew G. Knepley 
657bc4ae4beSMatthew G. Knepley   Level: beginner
658bc4ae4beSMatthew G. Knepley 
659bc4ae4beSMatthew G. Knepley .seealso: PetscDSGetDiscretization(), PetscDSAddDiscretization(), PetscDSGetNumFields(), PetscDSCreate()
660bc4ae4beSMatthew G. Knepley @*/
6612764a2aaSMatthew G. Knepley PetscErrorCode PetscDSSetDiscretization(PetscDS prob, PetscInt f, PetscObject disc)
6622764a2aaSMatthew G. Knepley {
6632764a2aaSMatthew G. Knepley   PetscErrorCode ierr;
6642764a2aaSMatthew G. Knepley 
6652764a2aaSMatthew G. Knepley   PetscFunctionBegin;
6662764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
6672764a2aaSMatthew G. Knepley   PetscValidPointer(disc, 3);
6682764a2aaSMatthew G. Knepley   if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f);
6692764a2aaSMatthew G. Knepley   ierr = PetscDSEnlarge_Static(prob, f+1);CHKERRQ(ierr);
6702764a2aaSMatthew G. Knepley   if (prob->disc[f]) {ierr = PetscObjectDereference(prob->disc[f]);CHKERRQ(ierr);}
6712764a2aaSMatthew G. Knepley   prob->disc[f] = disc;
6722764a2aaSMatthew G. Knepley   ierr = PetscObjectReference(disc);CHKERRQ(ierr);
673249df284SMatthew G. Knepley   {
674249df284SMatthew G. Knepley     PetscClassId id;
675249df284SMatthew G. Knepley 
676249df284SMatthew G. Knepley     ierr = PetscObjectGetClassId(disc, &id);CHKERRQ(ierr);
677a6cbbb48SMatthew G. Knepley     if (id == PETSCFV_CLASSID) {
678a6cbbb48SMatthew G. Knepley       prob->implicit[f]      = PETSC_FALSE;
679a6cbbb48SMatthew G. Knepley       prob->adjacency[f*2+0] = PETSC_TRUE;
680a6cbbb48SMatthew G. Knepley       prob->adjacency[f*2+1] = PETSC_FALSE;
681a6cbbb48SMatthew G. Knepley     }
682249df284SMatthew G. Knepley   }
6832764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
6842764a2aaSMatthew G. Knepley }
6852764a2aaSMatthew G. Knepley 
686bc4ae4beSMatthew G. Knepley /*@
687bc4ae4beSMatthew G. Knepley   PetscDSAddDiscretization - Adds a discretization object
688bc4ae4beSMatthew G. Knepley 
689bc4ae4beSMatthew G. Knepley   Not collective
690bc4ae4beSMatthew G. Knepley 
691bc4ae4beSMatthew G. Knepley   Input Parameters:
692bc4ae4beSMatthew G. Knepley + prob - The PetscDS object
693bc4ae4beSMatthew G. Knepley - disc - The boundary discretization object
694bc4ae4beSMatthew G. Knepley 
695bc4ae4beSMatthew G. Knepley   Level: beginner
696bc4ae4beSMatthew G. Knepley 
697bc4ae4beSMatthew G. Knepley .seealso: PetscDSGetDiscretization(), PetscDSSetDiscretization(), PetscDSGetNumFields(), PetscDSCreate()
698bc4ae4beSMatthew G. Knepley @*/
6992764a2aaSMatthew G. Knepley PetscErrorCode PetscDSAddDiscretization(PetscDS prob, PetscObject disc)
7002764a2aaSMatthew G. Knepley {
7012764a2aaSMatthew G. Knepley   PetscErrorCode ierr;
7022764a2aaSMatthew G. Knepley 
7032764a2aaSMatthew G. Knepley   PetscFunctionBegin;
7042764a2aaSMatthew G. Knepley   ierr = PetscDSSetDiscretization(prob, prob->Nf, disc);CHKERRQ(ierr);
7052764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
7062764a2aaSMatthew G. Knepley }
7072764a2aaSMatthew G. Knepley 
708249df284SMatthew G. Knepley /*@
709249df284SMatthew G. Knepley   PetscDSGetImplicit - Returns the flag for implicit solve for this field. This is just a guide for IMEX
710249df284SMatthew G. Knepley 
711249df284SMatthew G. Knepley   Not collective
712249df284SMatthew G. Knepley 
713249df284SMatthew G. Knepley   Input Parameters:
714249df284SMatthew G. Knepley + prob - The PetscDS object
715249df284SMatthew G. Knepley - f - The field number
716249df284SMatthew G. Knepley 
717249df284SMatthew G. Knepley   Output Parameter:
718249df284SMatthew G. Knepley . implicit - The flag indicating what kind of solve to use for this field
719249df284SMatthew G. Knepley 
720249df284SMatthew G. Knepley   Level: developer
721249df284SMatthew G. Knepley 
722f744cafaSSander Arens .seealso: PetscDSSetImplicit(), PetscDSSetDiscretization(), PetscDSAddDiscretization(), PetscDSGetNumFields(), PetscDSCreate()
723249df284SMatthew G. Knepley @*/
724249df284SMatthew G. Knepley PetscErrorCode PetscDSGetImplicit(PetscDS prob, PetscInt f, PetscBool *implicit)
725249df284SMatthew G. Knepley {
726249df284SMatthew G. Knepley   PetscFunctionBegin;
727249df284SMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
728249df284SMatthew G. Knepley   PetscValidPointer(implicit, 3);
729249df284SMatthew G. Knepley   if ((f < 0) || (f >= prob->Nf)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be in [0, %d)", f, prob->Nf);
730249df284SMatthew G. Knepley   *implicit = prob->implicit[f];
731249df284SMatthew G. Knepley   PetscFunctionReturn(0);
732249df284SMatthew G. Knepley }
733249df284SMatthew G. Knepley 
734249df284SMatthew G. Knepley /*@
735249df284SMatthew G. Knepley   PetscDSSetImplicit - Set the flag for implicit solve for this field. This is just a guide for IMEX
736249df284SMatthew G. Knepley 
737249df284SMatthew G. Knepley   Not collective
738249df284SMatthew G. Knepley 
739249df284SMatthew G. Knepley   Input Parameters:
740249df284SMatthew G. Knepley + prob - The PetscDS object
741249df284SMatthew G. Knepley . f - The field number
742249df284SMatthew G. Knepley - implicit - The flag indicating what kind of solve to use for this field
743249df284SMatthew G. Knepley 
744249df284SMatthew G. Knepley   Level: developer
745249df284SMatthew G. Knepley 
746f744cafaSSander Arens .seealso: PetscDSGetImplicit(), PetscDSSetDiscretization(), PetscDSAddDiscretization(), PetscDSGetNumFields(), PetscDSCreate()
747249df284SMatthew G. Knepley @*/
748249df284SMatthew G. Knepley PetscErrorCode PetscDSSetImplicit(PetscDS prob, PetscInt f, PetscBool implicit)
749249df284SMatthew G. Knepley {
750249df284SMatthew G. Knepley   PetscFunctionBegin;
751249df284SMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
752249df284SMatthew G. Knepley   if ((f < 0) || (f >= prob->Nf)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be in [0, %d)", f, prob->Nf);
753249df284SMatthew G. Knepley   prob->implicit[f] = implicit;
754249df284SMatthew G. Knepley   PetscFunctionReturn(0);
755249df284SMatthew G. Knepley }
756249df284SMatthew G. Knepley 
757a6cbbb48SMatthew G. Knepley /*@
758a6cbbb48SMatthew G. Knepley   PetscDSGetAdjacency - Returns the flags for determining variable influence
759a6cbbb48SMatthew G. Knepley 
760a6cbbb48SMatthew G. Knepley   Not collective
761a6cbbb48SMatthew G. Knepley 
762a6cbbb48SMatthew G. Knepley   Input Parameters:
763a6cbbb48SMatthew G. Knepley + prob - The PetscDS object
764a6cbbb48SMatthew G. Knepley - f - The field number
765a6cbbb48SMatthew G. Knepley 
766a6cbbb48SMatthew G. Knepley   Output Parameter:
767a6cbbb48SMatthew G. Knepley + useCone    - Flag for variable influence starting with the cone operation
768a6cbbb48SMatthew G. Knepley - useClosure - Flag for variable influence using transitive closure
769a6cbbb48SMatthew G. Knepley 
770a6cbbb48SMatthew G. Knepley   Note: See the discussion in DMPlexGetAdjacencyUseCone() and DMPlexGetAdjacencyUseClosure()
771a6cbbb48SMatthew G. Knepley 
772a6cbbb48SMatthew G. Knepley   Level: developer
773a6cbbb48SMatthew G. Knepley 
774f744cafaSSander Arens .seealso: PetscDSSetAdjacency(), DMPlexGetAdjacencyUseCone(), DMPlexGetAdjacencyUseClosure(), PetscDSSetDiscretization(), PetscDSAddDiscretization(), PetscDSGetNumFields(), PetscDSCreate()
775a6cbbb48SMatthew G. Knepley @*/
776a6cbbb48SMatthew G. Knepley PetscErrorCode PetscDSGetAdjacency(PetscDS prob, PetscInt f, PetscBool *useCone, PetscBool *useClosure)
777a6cbbb48SMatthew G. Knepley {
778a6cbbb48SMatthew G. Knepley   PetscFunctionBegin;
779a6cbbb48SMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
780a6cbbb48SMatthew G. Knepley   PetscValidPointer(useCone, 3);
781a6cbbb48SMatthew G. Knepley   PetscValidPointer(useClosure, 4);
782a6cbbb48SMatthew G. Knepley   if ((f < 0) || (f >= prob->Nf)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be in [0, %d)", f, prob->Nf);
783a6cbbb48SMatthew G. Knepley   *useCone    = prob->adjacency[f*2+0];
784a6cbbb48SMatthew G. Knepley   *useClosure = prob->adjacency[f*2+1];
785a6cbbb48SMatthew G. Knepley   PetscFunctionReturn(0);
786a6cbbb48SMatthew G. Knepley }
787a6cbbb48SMatthew G. Knepley 
788a6cbbb48SMatthew G. Knepley /*@
789a6cbbb48SMatthew G. Knepley   PetscDSSetAdjacency - Set the flags for determining variable influence
790a6cbbb48SMatthew G. Knepley 
791a6cbbb48SMatthew G. Knepley   Not collective
792a6cbbb48SMatthew G. Knepley 
793a6cbbb48SMatthew G. Knepley   Input Parameters:
794a6cbbb48SMatthew G. Knepley + prob - The PetscDS object
795a6cbbb48SMatthew G. Knepley . f - The field number
796a6cbbb48SMatthew G. Knepley . useCone    - Flag for variable influence starting with the cone operation
797a6cbbb48SMatthew G. Knepley - useClosure - Flag for variable influence using transitive closure
798a6cbbb48SMatthew G. Knepley 
799a6cbbb48SMatthew G. Knepley   Note: See the discussion in DMPlexGetAdjacencyUseCone() and DMPlexGetAdjacencyUseClosure()
800a6cbbb48SMatthew G. Knepley 
801a6cbbb48SMatthew G. Knepley   Level: developer
802a6cbbb48SMatthew G. Knepley 
803f744cafaSSander Arens .seealso: PetscDSGetAdjacency(), DMPlexGetAdjacencyUseCone(), DMPlexGetAdjacencyUseClosure(), PetscDSSetDiscretization(), PetscDSAddDiscretization(), PetscDSGetNumFields(), PetscDSCreate()
804a6cbbb48SMatthew G. Knepley @*/
805a6cbbb48SMatthew G. Knepley PetscErrorCode PetscDSSetAdjacency(PetscDS prob, PetscInt f, PetscBool useCone, PetscBool useClosure)
806a6cbbb48SMatthew G. Knepley {
807a6cbbb48SMatthew G. Knepley   PetscFunctionBegin;
808a6cbbb48SMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
809a6cbbb48SMatthew G. Knepley   if ((f < 0) || (f >= prob->Nf)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be in [0, %d)", f, prob->Nf);
810a6cbbb48SMatthew G. Knepley   prob->adjacency[f*2+0] = useCone;
811a6cbbb48SMatthew G. Knepley   prob->adjacency[f*2+1] = useClosure;
812a6cbbb48SMatthew G. Knepley   PetscFunctionReturn(0);
813a6cbbb48SMatthew G. Knepley }
814a6cbbb48SMatthew G. Knepley 
8152764a2aaSMatthew G. Knepley PetscErrorCode PetscDSGetObjective(PetscDS prob, PetscInt f,
81630b9ff8bSMatthew G. Knepley                                    void (**obj)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
817194d53e6SMatthew G. Knepley                                                 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
818194d53e6SMatthew G. Knepley                                                 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
81997b6e6e8SMatthew G. Knepley                                                 PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar obj[]))
8202764a2aaSMatthew G. Knepley {
8212764a2aaSMatthew G. Knepley   PetscFunctionBegin;
8222764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
8232764a2aaSMatthew G. Knepley   PetscValidPointer(obj, 2);
8242764a2aaSMatthew G. Knepley   if ((f < 0) || (f >= prob->Nf)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be in [0, %d)", f, prob->Nf);
8252764a2aaSMatthew G. Knepley   *obj = prob->obj[f];
8262764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
8272764a2aaSMatthew G. Knepley }
8282764a2aaSMatthew G. Knepley 
8292764a2aaSMatthew G. Knepley PetscErrorCode PetscDSSetObjective(PetscDS prob, PetscInt f,
83030b9ff8bSMatthew G. Knepley                                    void (*obj)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
831194d53e6SMatthew G. Knepley                                                const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
832194d53e6SMatthew G. Knepley                                                const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
83397b6e6e8SMatthew G. Knepley                                                PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar obj[]))
8342764a2aaSMatthew G. Knepley {
8352764a2aaSMatthew G. Knepley   PetscErrorCode ierr;
8362764a2aaSMatthew G. Knepley 
8372764a2aaSMatthew G. Knepley   PetscFunctionBegin;
8382764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
839de716cbcSToby Isaac   if (obj) PetscValidFunction(obj, 2);
8402764a2aaSMatthew G. Knepley   if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f);
8412764a2aaSMatthew G. Knepley   ierr = PetscDSEnlarge_Static(prob, f+1);CHKERRQ(ierr);
8422764a2aaSMatthew G. Knepley   prob->obj[f] = obj;
8432764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
8442764a2aaSMatthew G. Knepley }
8452764a2aaSMatthew G. Knepley 
846194d53e6SMatthew G. Knepley /*@C
847194d53e6SMatthew G. Knepley   PetscDSGetResidual - Get the pointwise residual function for a given test field
848194d53e6SMatthew G. Knepley 
849194d53e6SMatthew G. Knepley   Not collective
850194d53e6SMatthew G. Knepley 
851194d53e6SMatthew G. Knepley   Input Parameters:
852194d53e6SMatthew G. Knepley + prob - The PetscDS
853194d53e6SMatthew G. Knepley - f    - The test field number
854194d53e6SMatthew G. Knepley 
855194d53e6SMatthew G. Knepley   Output Parameters:
856194d53e6SMatthew G. Knepley + f0 - integrand for the test function term
857194d53e6SMatthew G. Knepley - f1 - integrand for the test function gradient term
858194d53e6SMatthew G. Knepley 
859194d53e6SMatthew G. Knepley   Note: We are using a first order FEM model for the weak form:
860194d53e6SMatthew G. Knepley 
861194d53e6SMatthew G. Knepley   \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)
862194d53e6SMatthew G. Knepley 
863194d53e6SMatthew G. Knepley The calling sequence for the callbacks f0 and f1 is given by:
864194d53e6SMatthew G. Knepley 
86530b9ff8bSMatthew G. Knepley $ f0(PetscInt dim, PetscInt Nf, PetscInt NfAux,
866194d53e6SMatthew G. Knepley $    const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
867194d53e6SMatthew G. Knepley $    const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
86830b9ff8bSMatthew G. Knepley $    PetscReal t, const PetscReal x[], PetscScalar f0[])
869194d53e6SMatthew G. Knepley 
870194d53e6SMatthew G. Knepley + dim - the spatial dimension
871194d53e6SMatthew G. Knepley . Nf - the number of fields
872194d53e6SMatthew G. Knepley . uOff - the offset into u[] and u_t[] for each field
873194d53e6SMatthew G. Knepley . uOff_x - the offset into u_x[] for each field
874194d53e6SMatthew G. Knepley . u - each field evaluated at the current point
875194d53e6SMatthew G. Knepley . u_t - the time derivative of each field evaluated at the current point
876194d53e6SMatthew G. Knepley . u_x - the gradient of each field evaluated at the current point
877194d53e6SMatthew G. Knepley . aOff - the offset into a[] and a_t[] for each auxiliary field
878194d53e6SMatthew G. Knepley . aOff_x - the offset into a_x[] for each auxiliary field
879194d53e6SMatthew G. Knepley . a - each auxiliary field evaluated at the current point
880194d53e6SMatthew G. Knepley . a_t - the time derivative of each auxiliary field evaluated at the current point
881194d53e6SMatthew G. Knepley . a_x - the gradient of auxiliary each field evaluated at the current point
882194d53e6SMatthew G. Knepley . t - current time
883194d53e6SMatthew G. Knepley . x - coordinates of the current point
88497b6e6e8SMatthew G. Knepley . numConstants - number of constant parameters
88597b6e6e8SMatthew G. Knepley . constants - constant parameters
886194d53e6SMatthew G. Knepley - f0 - output values at the current point
887194d53e6SMatthew G. Knepley 
888194d53e6SMatthew G. Knepley   Level: intermediate
889194d53e6SMatthew G. Knepley 
890194d53e6SMatthew G. Knepley .seealso: PetscDSSetResidual()
891194d53e6SMatthew G. Knepley @*/
8922764a2aaSMatthew G. Knepley PetscErrorCode PetscDSGetResidual(PetscDS prob, PetscInt f,
89330b9ff8bSMatthew G. Knepley                                   void (**f0)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
894194d53e6SMatthew G. Knepley                                               const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
895194d53e6SMatthew G. Knepley                                               const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
89697b6e6e8SMatthew G. Knepley                                               PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]),
89730b9ff8bSMatthew G. Knepley                                   void (**f1)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
898194d53e6SMatthew G. Knepley                                               const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
899194d53e6SMatthew G. Knepley                                               const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
90097b6e6e8SMatthew G. Knepley                                               PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]))
9012764a2aaSMatthew G. Knepley {
9022764a2aaSMatthew G. Knepley   PetscFunctionBegin;
9032764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
9042764a2aaSMatthew G. Knepley   if ((f < 0) || (f >= prob->Nf)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be in [0, %d)", f, prob->Nf);
9052764a2aaSMatthew G. Knepley   if (f0) {PetscValidPointer(f0, 3); *f0 = prob->f[f*2+0];}
9062764a2aaSMatthew G. Knepley   if (f1) {PetscValidPointer(f1, 4); *f1 = prob->f[f*2+1];}
9072764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
9082764a2aaSMatthew G. Knepley }
9092764a2aaSMatthew G. Knepley 
910194d53e6SMatthew G. Knepley /*@C
911194d53e6SMatthew G. Knepley   PetscDSSetResidual - Set the pointwise residual function for a given test field
912194d53e6SMatthew G. Knepley 
913194d53e6SMatthew G. Knepley   Not collective
914194d53e6SMatthew G. Knepley 
915194d53e6SMatthew G. Knepley   Input Parameters:
916194d53e6SMatthew G. Knepley + prob - The PetscDS
917194d53e6SMatthew G. Knepley . f    - The test field number
918194d53e6SMatthew G. Knepley . f0 - integrand for the test function term
919194d53e6SMatthew G. Knepley - f1 - integrand for the test function gradient term
920194d53e6SMatthew G. Knepley 
921194d53e6SMatthew G. Knepley   Note: We are using a first order FEM model for the weak form:
922194d53e6SMatthew G. Knepley 
923194d53e6SMatthew G. Knepley   \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)
924194d53e6SMatthew G. Knepley 
925194d53e6SMatthew G. Knepley The calling sequence for the callbacks f0 and f1 is given by:
926194d53e6SMatthew G. Knepley 
92730b9ff8bSMatthew G. Knepley $ f0(PetscInt dim, PetscInt Nf, PetscInt NfAux,
928194d53e6SMatthew G. Knepley $    const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
929194d53e6SMatthew G. Knepley $    const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
93030b9ff8bSMatthew G. Knepley $    PetscReal t, const PetscReal x[], PetscScalar f0[])
931194d53e6SMatthew G. Knepley 
932194d53e6SMatthew G. Knepley + dim - the spatial dimension
933194d53e6SMatthew G. Knepley . Nf - the number of fields
934194d53e6SMatthew G. Knepley . uOff - the offset into u[] and u_t[] for each field
935194d53e6SMatthew G. Knepley . uOff_x - the offset into u_x[] for each field
936194d53e6SMatthew G. Knepley . u - each field evaluated at the current point
937194d53e6SMatthew G. Knepley . u_t - the time derivative of each field evaluated at the current point
938194d53e6SMatthew G. Knepley . u_x - the gradient of each field evaluated at the current point
939194d53e6SMatthew G. Knepley . aOff - the offset into a[] and a_t[] for each auxiliary field
940194d53e6SMatthew G. Knepley . aOff_x - the offset into a_x[] for each auxiliary field
941194d53e6SMatthew G. Knepley . a - each auxiliary field evaluated at the current point
942194d53e6SMatthew G. Knepley . a_t - the time derivative of each auxiliary field evaluated at the current point
943194d53e6SMatthew G. Knepley . a_x - the gradient of auxiliary each field evaluated at the current point
944194d53e6SMatthew G. Knepley . t - current time
945194d53e6SMatthew G. Knepley . x - coordinates of the current point
94697b6e6e8SMatthew G. Knepley . numConstants - number of constant parameters
94797b6e6e8SMatthew G. Knepley . constants - constant parameters
948194d53e6SMatthew G. Knepley - f0 - output values at the current point
949194d53e6SMatthew G. Knepley 
950194d53e6SMatthew G. Knepley   Level: intermediate
951194d53e6SMatthew G. Knepley 
952194d53e6SMatthew G. Knepley .seealso: PetscDSGetResidual()
953194d53e6SMatthew G. Knepley @*/
9542764a2aaSMatthew G. Knepley PetscErrorCode PetscDSSetResidual(PetscDS prob, PetscInt f,
95530b9ff8bSMatthew G. Knepley                                   void (*f0)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
956194d53e6SMatthew G. Knepley                                              const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
957194d53e6SMatthew G. Knepley                                              const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
95897b6e6e8SMatthew G. Knepley                                              PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]),
95930b9ff8bSMatthew G. Knepley                                   void (*f1)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
960194d53e6SMatthew G. Knepley                                              const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
961194d53e6SMatthew G. Knepley                                              const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
96297b6e6e8SMatthew G. Knepley                                              PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]))
9632764a2aaSMatthew G. Knepley {
9642764a2aaSMatthew G. Knepley   PetscErrorCode ierr;
9652764a2aaSMatthew G. Knepley 
9662764a2aaSMatthew G. Knepley   PetscFunctionBegin;
9672764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
968f866a1d0SMatthew G. Knepley   if (f0) PetscValidFunction(f0, 3);
969f866a1d0SMatthew G. Knepley   if (f1) PetscValidFunction(f1, 4);
9702764a2aaSMatthew G. Knepley   if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f);
9712764a2aaSMatthew G. Knepley   ierr = PetscDSEnlarge_Static(prob, f+1);CHKERRQ(ierr);
9722764a2aaSMatthew G. Knepley   prob->f[f*2+0] = f0;
9732764a2aaSMatthew G. Knepley   prob->f[f*2+1] = f1;
9742764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
9752764a2aaSMatthew G. Knepley }
9762764a2aaSMatthew G. Knepley 
9773e75805dSMatthew G. Knepley /*@C
9783e75805dSMatthew G. Knepley   PetscDSHasJacobian - Signals that Jacobian functions have been set
9793e75805dSMatthew G. Knepley 
9803e75805dSMatthew G. Knepley   Not collective
9813e75805dSMatthew G. Knepley 
9823e75805dSMatthew G. Knepley   Input Parameter:
9833e75805dSMatthew G. Knepley . prob - The PetscDS
9843e75805dSMatthew G. Knepley 
9853e75805dSMatthew G. Knepley   Output Parameter:
9863e75805dSMatthew G. Knepley . hasJac - flag that pointwise function for the Jacobian has been set
9873e75805dSMatthew G. Knepley 
9883e75805dSMatthew G. Knepley   Level: intermediate
9893e75805dSMatthew G. Knepley 
9903e75805dSMatthew G. Knepley .seealso: PetscDSGetJacobianPreconditioner(), PetscDSSetJacobianPreconditioner(), PetscDSGetJacobian()
9913e75805dSMatthew G. Knepley @*/
9923e75805dSMatthew G. Knepley PetscErrorCode PetscDSHasJacobian(PetscDS prob, PetscBool *hasJac)
9933e75805dSMatthew G. Knepley {
9943e75805dSMatthew G. Knepley   PetscInt f, g, h;
9953e75805dSMatthew G. Knepley 
9963e75805dSMatthew G. Knepley   PetscFunctionBegin;
9973e75805dSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
9983e75805dSMatthew G. Knepley   *hasJac = PETSC_FALSE;
9993e75805dSMatthew G. Knepley   for (f = 0; f < prob->Nf; ++f) {
10003e75805dSMatthew G. Knepley     for (g = 0; g < prob->Nf; ++g) {
10013e75805dSMatthew G. Knepley       for (h = 0; h < 4; ++h) {
10023e75805dSMatthew G. Knepley         if (prob->g[(f*prob->Nf + g)*4+h]) *hasJac = PETSC_TRUE;
10033e75805dSMatthew G. Knepley       }
10043e75805dSMatthew G. Knepley     }
10053e75805dSMatthew G. Knepley   }
10063e75805dSMatthew G. Knepley   PetscFunctionReturn(0);
10073e75805dSMatthew G. Knepley }
10083e75805dSMatthew G. Knepley 
1009194d53e6SMatthew G. Knepley /*@C
1010194d53e6SMatthew G. Knepley   PetscDSGetJacobian - Get the pointwise Jacobian function for given test and basis field
1011194d53e6SMatthew G. Knepley 
1012194d53e6SMatthew G. Knepley   Not collective
1013194d53e6SMatthew G. Knepley 
1014194d53e6SMatthew G. Knepley   Input Parameters:
1015194d53e6SMatthew G. Knepley + prob - The PetscDS
1016194d53e6SMatthew G. Knepley . f    - The test field number
1017194d53e6SMatthew G. Knepley - g    - The field number
1018194d53e6SMatthew G. Knepley 
1019194d53e6SMatthew G. Knepley   Output Parameters:
1020194d53e6SMatthew G. Knepley + g0 - integrand for the test and basis function term
1021194d53e6SMatthew G. Knepley . g1 - integrand for the test function and basis function gradient term
1022194d53e6SMatthew G. Knepley . g2 - integrand for the test function gradient and basis function term
1023194d53e6SMatthew G. Knepley - g3 - integrand for the test function gradient and basis function gradient term
1024194d53e6SMatthew G. Knepley 
1025194d53e6SMatthew G. Knepley   Note: We are using a first order FEM model for the weak form:
1026194d53e6SMatthew G. Knepley 
1027194d53e6SMatthew G. Knepley   \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
1028194d53e6SMatthew G. Knepley 
1029194d53e6SMatthew G. Knepley The calling sequence for the callbacks g0, g1, g2 and g3 is given by:
1030194d53e6SMatthew G. Knepley 
103130b9ff8bSMatthew G. Knepley $ g0(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1032194d53e6SMatthew G. Knepley $    const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1033194d53e6SMatthew G. Knepley $    const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
103430b9ff8bSMatthew G. Knepley $    PetscReal t, const PetscReal u_tShift, const PetscReal x[], PetscScalar g0[])
1035194d53e6SMatthew G. Knepley 
1036194d53e6SMatthew G. Knepley + dim - the spatial dimension
1037194d53e6SMatthew G. Knepley . Nf - the number of fields
1038194d53e6SMatthew G. Knepley . uOff - the offset into u[] and u_t[] for each field
1039194d53e6SMatthew G. Knepley . uOff_x - the offset into u_x[] for each field
1040194d53e6SMatthew G. Knepley . u - each field evaluated at the current point
1041194d53e6SMatthew G. Knepley . u_t - the time derivative of each field evaluated at the current point
1042194d53e6SMatthew G. Knepley . u_x - the gradient of each field evaluated at the current point
1043194d53e6SMatthew G. Knepley . aOff - the offset into a[] and a_t[] for each auxiliary field
1044194d53e6SMatthew G. Knepley . aOff_x - the offset into a_x[] for each auxiliary field
1045194d53e6SMatthew G. Knepley . a - each auxiliary field evaluated at the current point
1046194d53e6SMatthew G. Knepley . a_t - the time derivative of each auxiliary field evaluated at the current point
1047194d53e6SMatthew G. Knepley . a_x - the gradient of auxiliary each field evaluated at the current point
1048194d53e6SMatthew G. Knepley . t - current time
10492aa1fc23SMatthew G. Knepley . u_tShift - the multiplier a for dF/dU_t
1050194d53e6SMatthew G. Knepley . x - coordinates of the current point
105197b6e6e8SMatthew G. Knepley . numConstants - number of constant parameters
105297b6e6e8SMatthew G. Knepley . constants - constant parameters
1053194d53e6SMatthew G. Knepley - g0 - output values at the current point
1054194d53e6SMatthew G. Knepley 
1055194d53e6SMatthew G. Knepley   Level: intermediate
1056194d53e6SMatthew G. Knepley 
1057194d53e6SMatthew G. Knepley .seealso: PetscDSSetJacobian()
1058194d53e6SMatthew G. Knepley @*/
10592764a2aaSMatthew G. Knepley PetscErrorCode PetscDSGetJacobian(PetscDS prob, PetscInt f, PetscInt g,
106030b9ff8bSMatthew G. Knepley                                   void (**g0)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1061194d53e6SMatthew G. Knepley                                               const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1062194d53e6SMatthew G. Knepley                                               const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
106397b6e6e8SMatthew G. Knepley                                               PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]),
106430b9ff8bSMatthew G. Knepley                                   void (**g1)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1065194d53e6SMatthew G. Knepley                                               const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1066194d53e6SMatthew G. Knepley                                               const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
106797b6e6e8SMatthew G. Knepley                                               PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[]),
106830b9ff8bSMatthew G. Knepley                                   void (**g2)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1069194d53e6SMatthew G. Knepley                                               const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1070194d53e6SMatthew G. Knepley                                               const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
107197b6e6e8SMatthew G. Knepley                                               PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g2[]),
107230b9ff8bSMatthew G. Knepley                                   void (**g3)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1073194d53e6SMatthew G. Knepley                                               const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1074194d53e6SMatthew G. Knepley                                               const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
107597b6e6e8SMatthew G. Knepley                                               PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[]))
10762764a2aaSMatthew G. Knepley {
10772764a2aaSMatthew G. Knepley   PetscFunctionBegin;
10782764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
10792764a2aaSMatthew G. Knepley   if ((f < 0) || (f >= prob->Nf)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be in [0, %d)", f, prob->Nf);
10802764a2aaSMatthew G. Knepley   if ((g < 0) || (g >= prob->Nf)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be in [0, %d)", g, prob->Nf);
10812764a2aaSMatthew G. Knepley   if (g0) {PetscValidPointer(g0, 4); *g0 = prob->g[(f*prob->Nf + g)*4+0];}
10822764a2aaSMatthew G. Knepley   if (g1) {PetscValidPointer(g1, 5); *g1 = prob->g[(f*prob->Nf + g)*4+1];}
10832764a2aaSMatthew G. Knepley   if (g2) {PetscValidPointer(g2, 6); *g2 = prob->g[(f*prob->Nf + g)*4+2];}
10842764a2aaSMatthew G. Knepley   if (g3) {PetscValidPointer(g3, 7); *g3 = prob->g[(f*prob->Nf + g)*4+3];}
10852764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
10862764a2aaSMatthew G. Knepley }
10872764a2aaSMatthew G. Knepley 
1088194d53e6SMatthew G. Knepley /*@C
1089194d53e6SMatthew G. Knepley   PetscDSSetJacobian - Set the pointwise Jacobian function for given test and basis fields
1090194d53e6SMatthew G. Knepley 
1091194d53e6SMatthew G. Knepley   Not collective
1092194d53e6SMatthew G. Knepley 
1093194d53e6SMatthew G. Knepley   Input Parameters:
1094194d53e6SMatthew G. Knepley + prob - The PetscDS
1095194d53e6SMatthew G. Knepley . f    - The test field number
1096194d53e6SMatthew G. Knepley . g    - The field number
1097194d53e6SMatthew G. Knepley . g0 - integrand for the test and basis function term
1098194d53e6SMatthew G. Knepley . g1 - integrand for the test function and basis function gradient term
1099194d53e6SMatthew G. Knepley . g2 - integrand for the test function gradient and basis function term
1100194d53e6SMatthew G. Knepley - g3 - integrand for the test function gradient and basis function gradient term
1101194d53e6SMatthew G. Knepley 
1102194d53e6SMatthew G. Knepley   Note: We are using a first order FEM model for the weak form:
1103194d53e6SMatthew G. Knepley 
1104194d53e6SMatthew G. Knepley   \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
1105194d53e6SMatthew G. Knepley 
1106194d53e6SMatthew G. Knepley The calling sequence for the callbacks g0, g1, g2 and g3 is given by:
1107194d53e6SMatthew G. Knepley 
110830b9ff8bSMatthew G. Knepley $ g0(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1109194d53e6SMatthew G. Knepley $    const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1110194d53e6SMatthew G. Knepley $    const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
111130b9ff8bSMatthew G. Knepley $    PetscReal t, const PetscReal x[], PetscScalar g0[])
1112194d53e6SMatthew G. Knepley 
1113194d53e6SMatthew G. Knepley + dim - the spatial dimension
1114194d53e6SMatthew G. Knepley . Nf - the number of fields
1115194d53e6SMatthew G. Knepley . uOff - the offset into u[] and u_t[] for each field
1116194d53e6SMatthew G. Knepley . uOff_x - the offset into u_x[] for each field
1117194d53e6SMatthew G. Knepley . u - each field evaluated at the current point
1118194d53e6SMatthew G. Knepley . u_t - the time derivative of each field evaluated at the current point
1119194d53e6SMatthew G. Knepley . u_x - the gradient of each field evaluated at the current point
1120194d53e6SMatthew G. Knepley . aOff - the offset into a[] and a_t[] for each auxiliary field
1121194d53e6SMatthew G. Knepley . aOff_x - the offset into a_x[] for each auxiliary field
1122194d53e6SMatthew G. Knepley . a - each auxiliary field evaluated at the current point
1123194d53e6SMatthew G. Knepley . a_t - the time derivative of each auxiliary field evaluated at the current point
1124194d53e6SMatthew G. Knepley . a_x - the gradient of auxiliary each field evaluated at the current point
1125194d53e6SMatthew G. Knepley . t - current time
11262aa1fc23SMatthew G. Knepley . u_tShift - the multiplier a for dF/dU_t
1127194d53e6SMatthew G. Knepley . x - coordinates of the current point
112897b6e6e8SMatthew G. Knepley . numConstants - number of constant parameters
112997b6e6e8SMatthew G. Knepley . constants - constant parameters
1130194d53e6SMatthew G. Knepley - g0 - output values at the current point
1131194d53e6SMatthew G. Knepley 
1132194d53e6SMatthew G. Knepley   Level: intermediate
1133194d53e6SMatthew G. Knepley 
1134194d53e6SMatthew G. Knepley .seealso: PetscDSGetJacobian()
1135194d53e6SMatthew G. Knepley @*/
11362764a2aaSMatthew G. Knepley PetscErrorCode PetscDSSetJacobian(PetscDS prob, PetscInt f, PetscInt g,
113730b9ff8bSMatthew G. Knepley                                   void (*g0)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1138194d53e6SMatthew G. Knepley                                              const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1139194d53e6SMatthew G. Knepley                                              const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
114097b6e6e8SMatthew G. Knepley                                              PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]),
114130b9ff8bSMatthew G. Knepley                                   void (*g1)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1142194d53e6SMatthew G. Knepley                                              const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1143194d53e6SMatthew G. Knepley                                              const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
114497b6e6e8SMatthew G. Knepley                                              PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[]),
114530b9ff8bSMatthew G. Knepley                                   void (*g2)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1146194d53e6SMatthew G. Knepley                                              const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1147194d53e6SMatthew G. Knepley                                              const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
114897b6e6e8SMatthew G. Knepley                                              PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g2[]),
114930b9ff8bSMatthew G. Knepley                                   void (*g3)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1150194d53e6SMatthew G. Knepley                                              const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1151194d53e6SMatthew G. Knepley                                              const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
115297b6e6e8SMatthew G. Knepley                                              PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[]))
11532764a2aaSMatthew G. Knepley {
11542764a2aaSMatthew G. Knepley   PetscErrorCode ierr;
11552764a2aaSMatthew G. Knepley 
11562764a2aaSMatthew G. Knepley   PetscFunctionBegin;
11572764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
11582764a2aaSMatthew G. Knepley   if (g0) PetscValidFunction(g0, 4);
11592764a2aaSMatthew G. Knepley   if (g1) PetscValidFunction(g1, 5);
11602764a2aaSMatthew G. Knepley   if (g2) PetscValidFunction(g2, 6);
11612764a2aaSMatthew G. Knepley   if (g3) PetscValidFunction(g3, 7);
11622764a2aaSMatthew G. Knepley   if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f);
11632764a2aaSMatthew G. Knepley   if (g < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", g);
11642764a2aaSMatthew G. Knepley   ierr = PetscDSEnlarge_Static(prob, PetscMax(f, g)+1);CHKERRQ(ierr);
11652764a2aaSMatthew G. Knepley   prob->g[(f*prob->Nf + g)*4+0] = g0;
11662764a2aaSMatthew G. Knepley   prob->g[(f*prob->Nf + g)*4+1] = g1;
11672764a2aaSMatthew G. Knepley   prob->g[(f*prob->Nf + g)*4+2] = g2;
11682764a2aaSMatthew G. Knepley   prob->g[(f*prob->Nf + g)*4+3] = g3;
11692764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
11702764a2aaSMatthew G. Knepley }
11712764a2aaSMatthew G. Knepley 
1172475e0ac9SMatthew G. Knepley /*@C
1173475e0ac9SMatthew G. Knepley   PetscDSHasJacobianPreconditioner - Signals that a Jacobian preconditioner matrix has been set
1174475e0ac9SMatthew G. Knepley 
1175475e0ac9SMatthew G. Knepley   Not collective
1176475e0ac9SMatthew G. Knepley 
1177475e0ac9SMatthew G. Knepley   Input Parameter:
1178475e0ac9SMatthew G. Knepley . prob - The PetscDS
1179475e0ac9SMatthew G. Knepley 
1180475e0ac9SMatthew G. Knepley   Output Parameter:
1181475e0ac9SMatthew G. Knepley . hasJacPre - flag that pointwise function for Jacobian preconditioner matrix has been set
1182475e0ac9SMatthew G. Knepley 
1183475e0ac9SMatthew G. Knepley   Level: intermediate
1184475e0ac9SMatthew G. Knepley 
1185475e0ac9SMatthew G. Knepley .seealso: PetscDSGetJacobianPreconditioner(), PetscDSSetJacobianPreconditioner(), PetscDSGetJacobian()
1186475e0ac9SMatthew G. Knepley @*/
1187475e0ac9SMatthew G. Knepley PetscErrorCode PetscDSHasJacobianPreconditioner(PetscDS prob, PetscBool *hasJacPre)
1188475e0ac9SMatthew G. Knepley {
1189475e0ac9SMatthew G. Knepley   PetscInt f, g, h;
1190475e0ac9SMatthew G. Knepley 
1191475e0ac9SMatthew G. Knepley   PetscFunctionBegin;
1192475e0ac9SMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1193475e0ac9SMatthew G. Knepley   *hasJacPre = PETSC_FALSE;
1194475e0ac9SMatthew G. Knepley   for (f = 0; f < prob->Nf; ++f) {
1195475e0ac9SMatthew G. Knepley     for (g = 0; g < prob->Nf; ++g) {
1196475e0ac9SMatthew G. Knepley       for (h = 0; h < 4; ++h) {
1197475e0ac9SMatthew G. Knepley         if (prob->gp[(f*prob->Nf + g)*4+h]) *hasJacPre = PETSC_TRUE;
1198475e0ac9SMatthew G. Knepley       }
1199475e0ac9SMatthew G. Knepley     }
1200475e0ac9SMatthew G. Knepley   }
1201475e0ac9SMatthew G. Knepley   PetscFunctionReturn(0);
1202475e0ac9SMatthew G. Knepley }
1203475e0ac9SMatthew G. Knepley 
1204475e0ac9SMatthew G. Knepley /*@C
1205475e0ac9SMatthew G. Knepley   PetscDSGetJacobianPreconditioner - Get the pointwise Jacobian preconditioner function for given test and basis field. If this is missing, the system matrix is used to build the preconditioner.
1206475e0ac9SMatthew G. Knepley 
1207475e0ac9SMatthew G. Knepley   Not collective
1208475e0ac9SMatthew G. Knepley 
1209475e0ac9SMatthew G. Knepley   Input Parameters:
1210475e0ac9SMatthew G. Knepley + prob - The PetscDS
1211475e0ac9SMatthew G. Knepley . f    - The test field number
1212475e0ac9SMatthew G. Knepley - g    - The field number
1213475e0ac9SMatthew G. Knepley 
1214475e0ac9SMatthew G. Knepley   Output Parameters:
1215475e0ac9SMatthew G. Knepley + g0 - integrand for the test and basis function term
1216475e0ac9SMatthew G. Knepley . g1 - integrand for the test function and basis function gradient term
1217475e0ac9SMatthew G. Knepley . g2 - integrand for the test function gradient and basis function term
1218475e0ac9SMatthew G. Knepley - g3 - integrand for the test function gradient and basis function gradient term
1219475e0ac9SMatthew G. Knepley 
1220475e0ac9SMatthew G. Knepley   Note: We are using a first order FEM model for the weak form:
1221475e0ac9SMatthew G. Knepley 
1222475e0ac9SMatthew G. Knepley   \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
1223475e0ac9SMatthew G. Knepley 
1224475e0ac9SMatthew G. Knepley The calling sequence for the callbacks g0, g1, g2 and g3 is given by:
1225475e0ac9SMatthew G. Knepley 
1226475e0ac9SMatthew G. Knepley $ g0(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1227475e0ac9SMatthew G. Knepley $    const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1228475e0ac9SMatthew G. Knepley $    const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1229475e0ac9SMatthew G. Knepley $    PetscReal t, const PetscReal u_tShift, const PetscReal x[], PetscScalar g0[])
1230475e0ac9SMatthew G. Knepley 
1231475e0ac9SMatthew G. Knepley + dim - the spatial dimension
1232475e0ac9SMatthew G. Knepley . Nf - the number of fields
1233475e0ac9SMatthew G. Knepley . uOff - the offset into u[] and u_t[] for each field
1234475e0ac9SMatthew G. Knepley . uOff_x - the offset into u_x[] for each field
1235475e0ac9SMatthew G. Knepley . u - each field evaluated at the current point
1236475e0ac9SMatthew G. Knepley . u_t - the time derivative of each field evaluated at the current point
1237475e0ac9SMatthew G. Knepley . u_x - the gradient of each field evaluated at the current point
1238475e0ac9SMatthew G. Knepley . aOff - the offset into a[] and a_t[] for each auxiliary field
1239475e0ac9SMatthew G. Knepley . aOff_x - the offset into a_x[] for each auxiliary field
1240475e0ac9SMatthew G. Knepley . a - each auxiliary field evaluated at the current point
1241475e0ac9SMatthew G. Knepley . a_t - the time derivative of each auxiliary field evaluated at the current point
1242475e0ac9SMatthew G. Knepley . a_x - the gradient of auxiliary each field evaluated at the current point
1243475e0ac9SMatthew G. Knepley . t - current time
1244475e0ac9SMatthew G. Knepley . u_tShift - the multiplier a for dF/dU_t
1245475e0ac9SMatthew G. Knepley . x - coordinates of the current point
124697b6e6e8SMatthew G. Knepley . numConstants - number of constant parameters
124797b6e6e8SMatthew G. Knepley . constants - constant parameters
1248475e0ac9SMatthew G. Knepley - g0 - output values at the current point
1249475e0ac9SMatthew G. Knepley 
1250475e0ac9SMatthew G. Knepley   Level: intermediate
1251475e0ac9SMatthew G. Knepley 
1252475e0ac9SMatthew G. Knepley .seealso: PetscDSSetJacobianPreconditioner(), PetscDSGetJacobian()
1253475e0ac9SMatthew G. Knepley @*/
1254475e0ac9SMatthew G. Knepley PetscErrorCode PetscDSGetJacobianPreconditioner(PetscDS prob, PetscInt f, PetscInt g,
1255475e0ac9SMatthew G. Knepley                                   void (**g0)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1256475e0ac9SMatthew G. Knepley                                               const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1257475e0ac9SMatthew G. Knepley                                               const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
125897b6e6e8SMatthew G. Knepley                                               PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]),
1259475e0ac9SMatthew G. Knepley                                   void (**g1)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1260475e0ac9SMatthew G. Knepley                                               const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1261475e0ac9SMatthew G. Knepley                                               const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
126297b6e6e8SMatthew G. Knepley                                               PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[]),
1263475e0ac9SMatthew G. Knepley                                   void (**g2)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1264475e0ac9SMatthew G. Knepley                                               const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1265475e0ac9SMatthew G. Knepley                                               const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
126697b6e6e8SMatthew G. Knepley                                               PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g2[]),
1267475e0ac9SMatthew G. Knepley                                   void (**g3)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1268475e0ac9SMatthew G. Knepley                                               const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1269475e0ac9SMatthew G. Knepley                                               const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
127097b6e6e8SMatthew G. Knepley                                               PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[]))
1271475e0ac9SMatthew G. Knepley {
1272475e0ac9SMatthew G. Knepley   PetscFunctionBegin;
1273475e0ac9SMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1274475e0ac9SMatthew G. Knepley   if ((f < 0) || (f >= prob->Nf)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be in [0, %d)", f, prob->Nf);
1275475e0ac9SMatthew G. Knepley   if ((g < 0) || (g >= prob->Nf)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be in [0, %d)", g, prob->Nf);
1276475e0ac9SMatthew G. Knepley   if (g0) {PetscValidPointer(g0, 4); *g0 = prob->gp[(f*prob->Nf + g)*4+0];}
1277475e0ac9SMatthew G. Knepley   if (g1) {PetscValidPointer(g1, 5); *g1 = prob->gp[(f*prob->Nf + g)*4+1];}
1278475e0ac9SMatthew G. Knepley   if (g2) {PetscValidPointer(g2, 6); *g2 = prob->gp[(f*prob->Nf + g)*4+2];}
1279475e0ac9SMatthew G. Knepley   if (g3) {PetscValidPointer(g3, 7); *g3 = prob->gp[(f*prob->Nf + g)*4+3];}
1280475e0ac9SMatthew G. Knepley   PetscFunctionReturn(0);
1281475e0ac9SMatthew G. Knepley }
1282475e0ac9SMatthew G. Knepley 
1283475e0ac9SMatthew G. Knepley /*@C
1284475e0ac9SMatthew G. Knepley   PetscDSSetJacobianPreconditioner - Set the pointwise Jacobian preconditioner function for given test and basis fields. If this is missing, the system matrix is used to build the preconditioner.
1285475e0ac9SMatthew G. Knepley 
1286475e0ac9SMatthew G. Knepley   Not collective
1287475e0ac9SMatthew G. Knepley 
1288475e0ac9SMatthew G. Knepley   Input Parameters:
1289475e0ac9SMatthew G. Knepley + prob - The PetscDS
1290475e0ac9SMatthew G. Knepley . f    - The test field number
1291475e0ac9SMatthew G. Knepley . g    - The field number
1292475e0ac9SMatthew G. Knepley . g0 - integrand for the test and basis function term
1293475e0ac9SMatthew G. Knepley . g1 - integrand for the test function and basis function gradient term
1294475e0ac9SMatthew G. Knepley . g2 - integrand for the test function gradient and basis function term
1295475e0ac9SMatthew G. Knepley - g3 - integrand for the test function gradient and basis function gradient term
1296475e0ac9SMatthew G. Knepley 
1297475e0ac9SMatthew G. Knepley   Note: We are using a first order FEM model for the weak form:
1298475e0ac9SMatthew G. Knepley 
1299475e0ac9SMatthew G. Knepley   \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
1300475e0ac9SMatthew G. Knepley 
1301475e0ac9SMatthew G. Knepley The calling sequence for the callbacks g0, g1, g2 and g3 is given by:
1302475e0ac9SMatthew G. Knepley 
1303475e0ac9SMatthew G. Knepley $ g0(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1304475e0ac9SMatthew G. Knepley $    const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1305475e0ac9SMatthew G. Knepley $    const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1306475e0ac9SMatthew G. Knepley $    PetscReal t, const PetscReal x[], PetscScalar g0[])
1307475e0ac9SMatthew G. Knepley 
1308475e0ac9SMatthew G. Knepley + dim - the spatial dimension
1309475e0ac9SMatthew G. Knepley . Nf - the number of fields
1310475e0ac9SMatthew G. Knepley . uOff - the offset into u[] and u_t[] for each field
1311475e0ac9SMatthew G. Knepley . uOff_x - the offset into u_x[] for each field
1312475e0ac9SMatthew G. Knepley . u - each field evaluated at the current point
1313475e0ac9SMatthew G. Knepley . u_t - the time derivative of each field evaluated at the current point
1314475e0ac9SMatthew G. Knepley . u_x - the gradient of each field evaluated at the current point
1315475e0ac9SMatthew G. Knepley . aOff - the offset into a[] and a_t[] for each auxiliary field
1316475e0ac9SMatthew G. Knepley . aOff_x - the offset into a_x[] for each auxiliary field
1317475e0ac9SMatthew G. Knepley . a - each auxiliary field evaluated at the current point
1318475e0ac9SMatthew G. Knepley . a_t - the time derivative of each auxiliary field evaluated at the current point
1319475e0ac9SMatthew G. Knepley . a_x - the gradient of auxiliary each field evaluated at the current point
1320475e0ac9SMatthew G. Knepley . t - current time
1321475e0ac9SMatthew G. Knepley . u_tShift - the multiplier a for dF/dU_t
1322475e0ac9SMatthew G. Knepley . x - coordinates of the current point
132397b6e6e8SMatthew G. Knepley . numConstants - number of constant parameters
132497b6e6e8SMatthew G. Knepley . constants - constant parameters
1325475e0ac9SMatthew G. Knepley - g0 - output values at the current point
1326475e0ac9SMatthew G. Knepley 
1327475e0ac9SMatthew G. Knepley   Level: intermediate
1328475e0ac9SMatthew G. Knepley 
1329475e0ac9SMatthew G. Knepley .seealso: PetscDSGetJacobianPreconditioner(), PetscDSSetJacobian()
1330475e0ac9SMatthew G. Knepley @*/
1331475e0ac9SMatthew G. Knepley PetscErrorCode PetscDSSetJacobianPreconditioner(PetscDS prob, PetscInt f, PetscInt g,
1332475e0ac9SMatthew G. Knepley                                   void (*g0)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1333475e0ac9SMatthew G. Knepley                                              const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1334475e0ac9SMatthew G. Knepley                                              const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
133597b6e6e8SMatthew G. Knepley                                              PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]),
1336475e0ac9SMatthew G. Knepley                                   void (*g1)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1337475e0ac9SMatthew G. Knepley                                              const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1338475e0ac9SMatthew G. Knepley                                              const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
133997b6e6e8SMatthew G. Knepley                                              PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[]),
1340475e0ac9SMatthew G. Knepley                                   void (*g2)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1341475e0ac9SMatthew G. Knepley                                              const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1342475e0ac9SMatthew G. Knepley                                              const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
134397b6e6e8SMatthew G. Knepley                                              PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g2[]),
1344475e0ac9SMatthew G. Knepley                                   void (*g3)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1345475e0ac9SMatthew G. Knepley                                              const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1346475e0ac9SMatthew G. Knepley                                              const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
134797b6e6e8SMatthew G. Knepley                                              PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[]))
1348475e0ac9SMatthew G. Knepley {
1349475e0ac9SMatthew G. Knepley   PetscErrorCode ierr;
1350475e0ac9SMatthew G. Knepley 
1351475e0ac9SMatthew G. Knepley   PetscFunctionBegin;
1352475e0ac9SMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1353475e0ac9SMatthew G. Knepley   if (g0) PetscValidFunction(g0, 4);
1354475e0ac9SMatthew G. Knepley   if (g1) PetscValidFunction(g1, 5);
1355475e0ac9SMatthew G. Knepley   if (g2) PetscValidFunction(g2, 6);
1356475e0ac9SMatthew G. Knepley   if (g3) PetscValidFunction(g3, 7);
1357475e0ac9SMatthew G. Knepley   if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f);
1358475e0ac9SMatthew G. Knepley   if (g < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", g);
1359475e0ac9SMatthew G. Knepley   ierr = PetscDSEnlarge_Static(prob, PetscMax(f, g)+1);CHKERRQ(ierr);
1360475e0ac9SMatthew G. Knepley   prob->gp[(f*prob->Nf + g)*4+0] = g0;
1361475e0ac9SMatthew G. Knepley   prob->gp[(f*prob->Nf + g)*4+1] = g1;
1362475e0ac9SMatthew G. Knepley   prob->gp[(f*prob->Nf + g)*4+2] = g2;
1363475e0ac9SMatthew G. Knepley   prob->gp[(f*prob->Nf + g)*4+3] = g3;
1364475e0ac9SMatthew G. Knepley   PetscFunctionReturn(0);
1365475e0ac9SMatthew G. Knepley }
1366475e0ac9SMatthew G. Knepley 
1367b7e05686SMatthew G. Knepley /*@C
1368b7e05686SMatthew G. Knepley   PetscDSHasDynamicJacobian - Signals that a dynamic Jacobian, dF/du_t, has been set
1369b7e05686SMatthew G. Knepley 
1370b7e05686SMatthew G. Knepley   Not collective
1371b7e05686SMatthew G. Knepley 
1372b7e05686SMatthew G. Knepley   Input Parameter:
1373b7e05686SMatthew G. Knepley . prob - The PetscDS
1374b7e05686SMatthew G. Knepley 
1375b7e05686SMatthew G. Knepley   Output Parameter:
1376b7e05686SMatthew G. Knepley . hasDynJac - flag that pointwise function for dynamic Jacobian has been set
1377b7e05686SMatthew G. Knepley 
1378b7e05686SMatthew G. Knepley   Level: intermediate
1379b7e05686SMatthew G. Knepley 
1380b7e05686SMatthew G. Knepley .seealso: PetscDSGetDynamicJacobian(), PetscDSSetDynamicJacobian(), PetscDSGetJacobian()
1381b7e05686SMatthew G. Knepley @*/
1382b7e05686SMatthew G. Knepley PetscErrorCode PetscDSHasDynamicJacobian(PetscDS prob, PetscBool *hasDynJac)
1383b7e05686SMatthew G. Knepley {
1384b7e05686SMatthew G. Knepley   PetscInt f, g, h;
1385b7e05686SMatthew G. Knepley 
1386b7e05686SMatthew G. Knepley   PetscFunctionBegin;
1387b7e05686SMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1388b7e05686SMatthew G. Knepley   *hasDynJac = PETSC_FALSE;
1389b7e05686SMatthew G. Knepley   for (f = 0; f < prob->Nf; ++f) {
1390b7e05686SMatthew G. Knepley     for (g = 0; g < prob->Nf; ++g) {
1391b7e05686SMatthew G. Knepley       for (h = 0; h < 4; ++h) {
1392b7e05686SMatthew G. Knepley         if (prob->gt[(f*prob->Nf + g)*4+h]) *hasDynJac = PETSC_TRUE;
1393b7e05686SMatthew G. Knepley       }
1394b7e05686SMatthew G. Knepley     }
1395b7e05686SMatthew G. Knepley   }
1396b7e05686SMatthew G. Knepley   PetscFunctionReturn(0);
1397b7e05686SMatthew G. Knepley }
1398b7e05686SMatthew G. Knepley 
1399b7e05686SMatthew G. Knepley /*@C
1400b7e05686SMatthew G. Knepley   PetscDSGetDynamicJacobian - Get the pointwise dynamic Jacobian, dF/du_t, function for given test and basis field
1401b7e05686SMatthew G. Knepley 
1402b7e05686SMatthew G. Knepley   Not collective
1403b7e05686SMatthew G. Knepley 
1404b7e05686SMatthew G. Knepley   Input Parameters:
1405b7e05686SMatthew G. Knepley + prob - The PetscDS
1406b7e05686SMatthew G. Knepley . f    - The test field number
1407b7e05686SMatthew G. Knepley - g    - The field number
1408b7e05686SMatthew G. Knepley 
1409b7e05686SMatthew G. Knepley   Output Parameters:
1410b7e05686SMatthew G. Knepley + g0 - integrand for the test and basis function term
1411b7e05686SMatthew G. Knepley . g1 - integrand for the test function and basis function gradient term
1412b7e05686SMatthew G. Knepley . g2 - integrand for the test function gradient and basis function term
1413b7e05686SMatthew G. Knepley - g3 - integrand for the test function gradient and basis function gradient term
1414b7e05686SMatthew G. Knepley 
1415b7e05686SMatthew G. Knepley   Note: We are using a first order FEM model for the weak form:
1416b7e05686SMatthew G. Knepley 
1417b7e05686SMatthew G. Knepley   \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
1418b7e05686SMatthew G. Knepley 
1419b7e05686SMatthew G. Knepley The calling sequence for the callbacks g0, g1, g2 and g3 is given by:
1420b7e05686SMatthew G. Knepley 
1421b7e05686SMatthew G. Knepley $ g0(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1422b7e05686SMatthew G. Knepley $    const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1423b7e05686SMatthew G. Knepley $    const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1424b7e05686SMatthew G. Knepley $    PetscReal t, const PetscReal u_tShift, const PetscReal x[], PetscScalar g0[])
1425b7e05686SMatthew G. Knepley 
1426b7e05686SMatthew G. Knepley + dim - the spatial dimension
1427b7e05686SMatthew G. Knepley . Nf - the number of fields
1428b7e05686SMatthew G. Knepley . uOff - the offset into u[] and u_t[] for each field
1429b7e05686SMatthew G. Knepley . uOff_x - the offset into u_x[] for each field
1430b7e05686SMatthew G. Knepley . u - each field evaluated at the current point
1431b7e05686SMatthew G. Knepley . u_t - the time derivative of each field evaluated at the current point
1432b7e05686SMatthew G. Knepley . u_x - the gradient of each field evaluated at the current point
1433b7e05686SMatthew G. Knepley . aOff - the offset into a[] and a_t[] for each auxiliary field
1434b7e05686SMatthew G. Knepley . aOff_x - the offset into a_x[] for each auxiliary field
1435b7e05686SMatthew G. Knepley . a - each auxiliary field evaluated at the current point
1436b7e05686SMatthew G. Knepley . a_t - the time derivative of each auxiliary field evaluated at the current point
1437b7e05686SMatthew G. Knepley . a_x - the gradient of auxiliary each field evaluated at the current point
1438b7e05686SMatthew G. Knepley . t - current time
1439b7e05686SMatthew G. Knepley . u_tShift - the multiplier a for dF/dU_t
1440b7e05686SMatthew G. Knepley . x - coordinates of the current point
144197b6e6e8SMatthew G. Knepley . numConstants - number of constant parameters
144297b6e6e8SMatthew G. Knepley . constants - constant parameters
1443b7e05686SMatthew G. Knepley - g0 - output values at the current point
1444b7e05686SMatthew G. Knepley 
1445b7e05686SMatthew G. Knepley   Level: intermediate
1446b7e05686SMatthew G. Knepley 
1447b7e05686SMatthew G. Knepley .seealso: PetscDSSetJacobian()
1448b7e05686SMatthew G. Knepley @*/
1449b7e05686SMatthew G. Knepley PetscErrorCode PetscDSGetDynamicJacobian(PetscDS prob, PetscInt f, PetscInt g,
1450b7e05686SMatthew G. Knepley                                          void (**g0)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1451b7e05686SMatthew G. Knepley                                                      const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1452b7e05686SMatthew G. Knepley                                                      const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
145397b6e6e8SMatthew G. Knepley                                                      PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]),
1454b7e05686SMatthew G. Knepley                                          void (**g1)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1455b7e05686SMatthew G. Knepley                                                      const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1456b7e05686SMatthew G. Knepley                                                      const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
145797b6e6e8SMatthew G. Knepley                                                      PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[]),
1458b7e05686SMatthew G. Knepley                                          void (**g2)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1459b7e05686SMatthew G. Knepley                                                      const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1460b7e05686SMatthew G. Knepley                                                      const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
146197b6e6e8SMatthew G. Knepley                                                      PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g2[]),
1462b7e05686SMatthew G. Knepley                                          void (**g3)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1463b7e05686SMatthew G. Knepley                                                      const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1464b7e05686SMatthew G. Knepley                                                      const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
146597b6e6e8SMatthew G. Knepley                                                      PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[]))
1466b7e05686SMatthew G. Knepley {
1467b7e05686SMatthew G. Knepley   PetscFunctionBegin;
1468b7e05686SMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1469b7e05686SMatthew G. Knepley   if ((f < 0) || (f >= prob->Nf)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be in [0, %d)", f, prob->Nf);
1470b7e05686SMatthew G. Knepley   if ((g < 0) || (g >= prob->Nf)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be in [0, %d)", g, prob->Nf);
1471b7e05686SMatthew G. Knepley   if (g0) {PetscValidPointer(g0, 4); *g0 = prob->gt[(f*prob->Nf + g)*4+0];}
1472b7e05686SMatthew G. Knepley   if (g1) {PetscValidPointer(g1, 5); *g1 = prob->gt[(f*prob->Nf + g)*4+1];}
1473b7e05686SMatthew G. Knepley   if (g2) {PetscValidPointer(g2, 6); *g2 = prob->gt[(f*prob->Nf + g)*4+2];}
1474b7e05686SMatthew G. Knepley   if (g3) {PetscValidPointer(g3, 7); *g3 = prob->gt[(f*prob->Nf + g)*4+3];}
1475b7e05686SMatthew G. Knepley   PetscFunctionReturn(0);
1476b7e05686SMatthew G. Knepley }
1477b7e05686SMatthew G. Knepley 
1478b7e05686SMatthew G. Knepley /*@C
1479b7e05686SMatthew G. Knepley   PetscDSSetDynamicJacobian - Set the pointwise dynamic Jacobian, dF/du_t, function for given test and basis fields
1480b7e05686SMatthew G. Knepley 
1481b7e05686SMatthew G. Knepley   Not collective
1482b7e05686SMatthew G. Knepley 
1483b7e05686SMatthew G. Knepley   Input Parameters:
1484b7e05686SMatthew G. Knepley + prob - The PetscDS
1485b7e05686SMatthew G. Knepley . f    - The test field number
1486b7e05686SMatthew G. Knepley . g    - The field number
1487b7e05686SMatthew G. Knepley . g0 - integrand for the test and basis function term
1488b7e05686SMatthew G. Knepley . g1 - integrand for the test function and basis function gradient term
1489b7e05686SMatthew G. Knepley . g2 - integrand for the test function gradient and basis function term
1490b7e05686SMatthew G. Knepley - g3 - integrand for the test function gradient and basis function gradient term
1491b7e05686SMatthew G. Knepley 
1492b7e05686SMatthew G. Knepley   Note: We are using a first order FEM model for the weak form:
1493b7e05686SMatthew G. Knepley 
1494b7e05686SMatthew G. Knepley   \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
1495b7e05686SMatthew G. Knepley 
1496b7e05686SMatthew G. Knepley The calling sequence for the callbacks g0, g1, g2 and g3 is given by:
1497b7e05686SMatthew G. Knepley 
1498b7e05686SMatthew G. Knepley $ g0(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1499b7e05686SMatthew G. Knepley $    const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1500b7e05686SMatthew G. Knepley $    const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1501b7e05686SMatthew G. Knepley $    PetscReal t, const PetscReal x[], PetscScalar g0[])
1502b7e05686SMatthew G. Knepley 
1503b7e05686SMatthew G. Knepley + dim - the spatial dimension
1504b7e05686SMatthew G. Knepley . Nf - the number of fields
1505b7e05686SMatthew G. Knepley . uOff - the offset into u[] and u_t[] for each field
1506b7e05686SMatthew G. Knepley . uOff_x - the offset into u_x[] for each field
1507b7e05686SMatthew G. Knepley . u - each field evaluated at the current point
1508b7e05686SMatthew G. Knepley . u_t - the time derivative of each field evaluated at the current point
1509b7e05686SMatthew G. Knepley . u_x - the gradient of each field evaluated at the current point
1510b7e05686SMatthew G. Knepley . aOff - the offset into a[] and a_t[] for each auxiliary field
1511b7e05686SMatthew G. Knepley . aOff_x - the offset into a_x[] for each auxiliary field
1512b7e05686SMatthew G. Knepley . a - each auxiliary field evaluated at the current point
1513b7e05686SMatthew G. Knepley . a_t - the time derivative of each auxiliary field evaluated at the current point
1514b7e05686SMatthew G. Knepley . a_x - the gradient of auxiliary each field evaluated at the current point
1515b7e05686SMatthew G. Knepley . t - current time
1516b7e05686SMatthew G. Knepley . u_tShift - the multiplier a for dF/dU_t
1517b7e05686SMatthew G. Knepley . x - coordinates of the current point
151897b6e6e8SMatthew G. Knepley . numConstants - number of constant parameters
151997b6e6e8SMatthew G. Knepley . constants - constant parameters
1520b7e05686SMatthew G. Knepley - g0 - output values at the current point
1521b7e05686SMatthew G. Knepley 
1522b7e05686SMatthew G. Knepley   Level: intermediate
1523b7e05686SMatthew G. Knepley 
1524b7e05686SMatthew G. Knepley .seealso: PetscDSGetJacobian()
1525b7e05686SMatthew G. Knepley @*/
1526b7e05686SMatthew G. Knepley PetscErrorCode PetscDSSetDynamicJacobian(PetscDS prob, PetscInt f, PetscInt g,
1527b7e05686SMatthew G. Knepley                                          void (*g0)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1528b7e05686SMatthew G. Knepley                                                     const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1529b7e05686SMatthew G. Knepley                                                     const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
153097b6e6e8SMatthew G. Knepley                                                     PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]),
1531b7e05686SMatthew G. Knepley                                          void (*g1)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1532b7e05686SMatthew G. Knepley                                                     const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1533b7e05686SMatthew G. Knepley                                                     const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
153497b6e6e8SMatthew G. Knepley                                                     PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[]),
1535b7e05686SMatthew G. Knepley                                          void (*g2)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1536b7e05686SMatthew G. Knepley                                                     const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1537b7e05686SMatthew G. Knepley                                                     const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
153897b6e6e8SMatthew G. Knepley                                                     PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g2[]),
1539b7e05686SMatthew G. Knepley                                          void (*g3)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1540b7e05686SMatthew G. Knepley                                                     const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1541b7e05686SMatthew G. Knepley                                                     const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
154297b6e6e8SMatthew G. Knepley                                                     PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[]))
1543b7e05686SMatthew G. Knepley {
1544b7e05686SMatthew G. Knepley   PetscErrorCode ierr;
1545b7e05686SMatthew G. Knepley 
1546b7e05686SMatthew G. Knepley   PetscFunctionBegin;
1547b7e05686SMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1548b7e05686SMatthew G. Knepley   if (g0) PetscValidFunction(g0, 4);
1549b7e05686SMatthew G. Knepley   if (g1) PetscValidFunction(g1, 5);
1550b7e05686SMatthew G. Knepley   if (g2) PetscValidFunction(g2, 6);
1551b7e05686SMatthew G. Knepley   if (g3) PetscValidFunction(g3, 7);
1552b7e05686SMatthew G. Knepley   if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f);
1553b7e05686SMatthew G. Knepley   if (g < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", g);
1554b7e05686SMatthew G. Knepley   ierr = PetscDSEnlarge_Static(prob, PetscMax(f, g)+1);CHKERRQ(ierr);
1555b7e05686SMatthew G. Knepley   prob->gt[(f*prob->Nf + g)*4+0] = g0;
1556b7e05686SMatthew G. Knepley   prob->gt[(f*prob->Nf + g)*4+1] = g1;
1557b7e05686SMatthew G. Knepley   prob->gt[(f*prob->Nf + g)*4+2] = g2;
1558b7e05686SMatthew G. Knepley   prob->gt[(f*prob->Nf + g)*4+3] = g3;
1559b7e05686SMatthew G. Knepley   PetscFunctionReturn(0);
1560b7e05686SMatthew G. Knepley }
1561b7e05686SMatthew G. Knepley 
15620c2f2876SMatthew G. Knepley /*@C
15630c2f2876SMatthew G. Knepley   PetscDSGetRiemannSolver - Returns the Riemann solver for the given field
15640c2f2876SMatthew G. Knepley 
15650c2f2876SMatthew G. Knepley   Not collective
15660c2f2876SMatthew G. Knepley 
15670c2f2876SMatthew G. Knepley   Input Arguments:
15680c2f2876SMatthew G. Knepley + prob - The PetscDS object
15690c2f2876SMatthew G. Knepley - f    - The field number
15700c2f2876SMatthew G. Knepley 
15710c2f2876SMatthew G. Knepley   Output Argument:
15720c2f2876SMatthew G. Knepley . r    - Riemann solver
15730c2f2876SMatthew G. Knepley 
15740c2f2876SMatthew G. Knepley   Calling sequence for r:
15750c2f2876SMatthew G. Knepley 
15765db36cf9SMatthew G. Knepley $ r(PetscInt dim, PetscInt Nf, const PetscReal x[], const PetscReal n[], const PetscScalar uL[], const PetscScalar uR[], PetscScalar flux[], void *ctx)
15770c2f2876SMatthew G. Knepley 
15785db36cf9SMatthew G. Knepley + dim  - The spatial dimension
15795db36cf9SMatthew G. Knepley . Nf   - The number of fields
15805db36cf9SMatthew G. Knepley . x    - The coordinates at a point on the interface
15810c2f2876SMatthew G. Knepley . n    - The normal vector to the interface
15820c2f2876SMatthew G. Knepley . uL   - The state vector to the left of the interface
15830c2f2876SMatthew G. Knepley . uR   - The state vector to the right of the interface
15840c2f2876SMatthew G. Knepley . flux - output array of flux through the interface
158597b6e6e8SMatthew G. Knepley . numConstants - number of constant parameters
158697b6e6e8SMatthew G. Knepley . constants - constant parameters
15870c2f2876SMatthew G. Knepley - ctx  - optional user context
15880c2f2876SMatthew G. Knepley 
15890c2f2876SMatthew G. Knepley   Level: intermediate
15900c2f2876SMatthew G. Knepley 
15910c2f2876SMatthew G. Knepley .seealso: PetscDSSetRiemannSolver()
15920c2f2876SMatthew G. Knepley @*/
15930c2f2876SMatthew G. Knepley PetscErrorCode PetscDSGetRiemannSolver(PetscDS prob, PetscInt f,
159497b6e6e8SMatthew G. Knepley                                        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))
15950c2f2876SMatthew G. Knepley {
15960c2f2876SMatthew G. Knepley   PetscFunctionBegin;
15970c2f2876SMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
15980c2f2876SMatthew G. Knepley   if ((f < 0) || (f >= prob->Nf)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be in [0, %d)", f, prob->Nf);
15990c2f2876SMatthew G. Knepley   PetscValidPointer(r, 3);
16000c2f2876SMatthew G. Knepley   *r = prob->r[f];
16010c2f2876SMatthew G. Knepley   PetscFunctionReturn(0);
16020c2f2876SMatthew G. Knepley }
16030c2f2876SMatthew G. Knepley 
16040c2f2876SMatthew G. Knepley /*@C
16050c2f2876SMatthew G. Knepley   PetscDSSetRiemannSolver - Sets the Riemann solver for the given field
16060c2f2876SMatthew G. Knepley 
16070c2f2876SMatthew G. Knepley   Not collective
16080c2f2876SMatthew G. Knepley 
16090c2f2876SMatthew G. Knepley   Input Arguments:
16100c2f2876SMatthew G. Knepley + prob - The PetscDS object
16110c2f2876SMatthew G. Knepley . f    - The field number
16120c2f2876SMatthew G. Knepley - r    - Riemann solver
16130c2f2876SMatthew G. Knepley 
16140c2f2876SMatthew G. Knepley   Calling sequence for r:
16150c2f2876SMatthew G. Knepley 
16165db36cf9SMatthew G. Knepley $ r(PetscInt dim, PetscInt Nf, const PetscReal x[], const PetscReal n[], const PetscScalar uL[], const PetscScalar uR[], PetscScalar flux[], void *ctx)
16170c2f2876SMatthew G. Knepley 
16185db36cf9SMatthew G. Knepley + dim  - The spatial dimension
16195db36cf9SMatthew G. Knepley . Nf   - The number of fields
16205db36cf9SMatthew G. Knepley . x    - The coordinates at a point on the interface
16210c2f2876SMatthew G. Knepley . n    - The normal vector to the interface
16220c2f2876SMatthew G. Knepley . uL   - The state vector to the left of the interface
16230c2f2876SMatthew G. Knepley . uR   - The state vector to the right of the interface
16240c2f2876SMatthew G. Knepley . flux - output array of flux through the interface
162597b6e6e8SMatthew G. Knepley . numConstants - number of constant parameters
162697b6e6e8SMatthew G. Knepley . constants - constant parameters
16270c2f2876SMatthew G. Knepley - ctx  - optional user context
16280c2f2876SMatthew G. Knepley 
16290c2f2876SMatthew G. Knepley   Level: intermediate
16300c2f2876SMatthew G. Knepley 
16310c2f2876SMatthew G. Knepley .seealso: PetscDSGetRiemannSolver()
16320c2f2876SMatthew G. Knepley @*/
16330c2f2876SMatthew G. Knepley PetscErrorCode PetscDSSetRiemannSolver(PetscDS prob, PetscInt f,
163497b6e6e8SMatthew G. Knepley                                        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))
16350c2f2876SMatthew G. Knepley {
16360c2f2876SMatthew G. Knepley   PetscErrorCode ierr;
16370c2f2876SMatthew G. Knepley 
16380c2f2876SMatthew G. Knepley   PetscFunctionBegin;
16390c2f2876SMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1640de716cbcSToby Isaac   if (r) PetscValidFunction(r, 3);
16410c2f2876SMatthew G. Knepley   if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f);
16420c2f2876SMatthew G. Knepley   ierr = PetscDSEnlarge_Static(prob, f+1);CHKERRQ(ierr);
16430c2f2876SMatthew G. Knepley   prob->r[f] = r;
16440c2f2876SMatthew G. Knepley   PetscFunctionReturn(0);
16450c2f2876SMatthew G. Knepley }
16460c2f2876SMatthew G. Knepley 
164732d2bbc9SMatthew G. Knepley /*@C
164832d2bbc9SMatthew G. Knepley   PetscDSGetUpdate - Get the pointwise update function for a given field
164932d2bbc9SMatthew G. Knepley 
165032d2bbc9SMatthew G. Knepley   Not collective
165132d2bbc9SMatthew G. Knepley 
165232d2bbc9SMatthew G. Knepley   Input Parameters:
165332d2bbc9SMatthew G. Knepley + prob - The PetscDS
165432d2bbc9SMatthew G. Knepley - f    - The field number
165532d2bbc9SMatthew G. Knepley 
165632d2bbc9SMatthew G. Knepley   Output Parameters:
165732d2bbc9SMatthew G. Knepley + update - update function
165832d2bbc9SMatthew G. Knepley 
165932d2bbc9SMatthew G. Knepley   Note: The calling sequence for the callback update is given by:
166032d2bbc9SMatthew G. Knepley 
166132d2bbc9SMatthew G. Knepley $ update(PetscInt dim, PetscInt Nf, PetscInt NfAux,
166232d2bbc9SMatthew G. Knepley $        const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
166332d2bbc9SMatthew G. Knepley $        const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
166432d2bbc9SMatthew G. Knepley $        PetscReal t, const PetscReal x[], PetscScalar uNew[])
166532d2bbc9SMatthew G. Knepley 
166632d2bbc9SMatthew G. Knepley + dim - the spatial dimension
166732d2bbc9SMatthew G. Knepley . Nf - the number of fields
166832d2bbc9SMatthew G. Knepley . uOff - the offset into u[] and u_t[] for each field
166932d2bbc9SMatthew G. Knepley . uOff_x - the offset into u_x[] for each field
167032d2bbc9SMatthew G. Knepley . u - each field evaluated at the current point
167132d2bbc9SMatthew G. Knepley . u_t - the time derivative of each field evaluated at the current point
167232d2bbc9SMatthew G. Knepley . u_x - the gradient of each field evaluated at the current point
167332d2bbc9SMatthew G. Knepley . aOff - the offset into a[] and a_t[] for each auxiliary field
167432d2bbc9SMatthew G. Knepley . aOff_x - the offset into a_x[] for each auxiliary field
167532d2bbc9SMatthew G. Knepley . a - each auxiliary field evaluated at the current point
167632d2bbc9SMatthew G. Knepley . a_t - the time derivative of each auxiliary field evaluated at the current point
167732d2bbc9SMatthew G. Knepley . a_x - the gradient of auxiliary each field evaluated at the current point
167832d2bbc9SMatthew G. Knepley . t - current time
167932d2bbc9SMatthew G. Knepley . x - coordinates of the current point
168032d2bbc9SMatthew G. Knepley - uNew - new value for field at the current point
168132d2bbc9SMatthew G. Knepley 
168232d2bbc9SMatthew G. Knepley   Level: intermediate
168332d2bbc9SMatthew G. Knepley 
168432d2bbc9SMatthew G. Knepley .seealso: PetscDSSetUpdate(), PetscDSSetResidual()
168532d2bbc9SMatthew G. Knepley @*/
168632d2bbc9SMatthew G. Knepley PetscErrorCode PetscDSGetUpdate(PetscDS prob, PetscInt f,
168732d2bbc9SMatthew G. Knepley                                   void (**update)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
168832d2bbc9SMatthew G. Knepley                                                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
168932d2bbc9SMatthew G. Knepley                                                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
16903fa77dffSMatthew G. Knepley                                                   PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar uNew[]))
169132d2bbc9SMatthew G. Knepley {
169232d2bbc9SMatthew G. Knepley   PetscFunctionBegin;
169332d2bbc9SMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
169432d2bbc9SMatthew G. Knepley   if ((f < 0) || (f >= prob->Nf)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be in [0, %d)", f, prob->Nf);
169532d2bbc9SMatthew G. Knepley   if (update) {PetscValidPointer(update, 3); *update = prob->update[f];}
169632d2bbc9SMatthew G. Knepley   PetscFunctionReturn(0);
169732d2bbc9SMatthew G. Knepley }
169832d2bbc9SMatthew G. Knepley 
169932d2bbc9SMatthew G. Knepley /*@C
17003fa77dffSMatthew G. Knepley   PetscDSSetUpdate - Set the pointwise update function for a given field
170132d2bbc9SMatthew G. Knepley 
170232d2bbc9SMatthew G. Knepley   Not collective
170332d2bbc9SMatthew G. Knepley 
170432d2bbc9SMatthew G. Knepley   Input Parameters:
170532d2bbc9SMatthew G. Knepley + prob   - The PetscDS
170632d2bbc9SMatthew G. Knepley . f      - The field number
170732d2bbc9SMatthew G. Knepley - update - update function
170832d2bbc9SMatthew G. Knepley 
170932d2bbc9SMatthew G. Knepley   Note: The calling sequence for the callback update is given by:
171032d2bbc9SMatthew G. Knepley 
171132d2bbc9SMatthew G. Knepley $ update(PetscInt dim, PetscInt Nf, PetscInt NfAux,
171232d2bbc9SMatthew G. Knepley $        const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
171332d2bbc9SMatthew G. Knepley $        const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
171432d2bbc9SMatthew G. Knepley $        PetscReal t, const PetscReal x[], PetscScalar uNew[])
171532d2bbc9SMatthew G. Knepley 
171632d2bbc9SMatthew G. Knepley + dim - the spatial dimension
171732d2bbc9SMatthew G. Knepley . Nf - the number of fields
171832d2bbc9SMatthew G. Knepley . uOff - the offset into u[] and u_t[] for each field
171932d2bbc9SMatthew G. Knepley . uOff_x - the offset into u_x[] for each field
172032d2bbc9SMatthew G. Knepley . u - each field evaluated at the current point
172132d2bbc9SMatthew G. Knepley . u_t - the time derivative of each field evaluated at the current point
172232d2bbc9SMatthew G. Knepley . u_x - the gradient of each field evaluated at the current point
172332d2bbc9SMatthew G. Knepley . aOff - the offset into a[] and a_t[] for each auxiliary field
172432d2bbc9SMatthew G. Knepley . aOff_x - the offset into a_x[] for each auxiliary field
172532d2bbc9SMatthew G. Knepley . a - each auxiliary field evaluated at the current point
172632d2bbc9SMatthew G. Knepley . a_t - the time derivative of each auxiliary field evaluated at the current point
172732d2bbc9SMatthew G. Knepley . a_x - the gradient of auxiliary each field evaluated at the current point
172832d2bbc9SMatthew G. Knepley . t - current time
172932d2bbc9SMatthew G. Knepley . x - coordinates of the current point
173032d2bbc9SMatthew G. Knepley - uNew - new field values at the current point
173132d2bbc9SMatthew G. Knepley 
173232d2bbc9SMatthew G. Knepley   Level: intermediate
173332d2bbc9SMatthew G. Knepley 
173432d2bbc9SMatthew G. Knepley .seealso: PetscDSGetResidual()
173532d2bbc9SMatthew G. Knepley @*/
173632d2bbc9SMatthew G. Knepley PetscErrorCode PetscDSSetUpdate(PetscDS prob, PetscInt f,
173732d2bbc9SMatthew G. Knepley                                 void (*update)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
173832d2bbc9SMatthew G. Knepley                                                const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
173932d2bbc9SMatthew G. Knepley                                                const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
17403fa77dffSMatthew G. Knepley                                                PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar uNew[]))
174132d2bbc9SMatthew G. Knepley {
174232d2bbc9SMatthew G. Knepley   PetscErrorCode ierr;
174332d2bbc9SMatthew G. Knepley 
174432d2bbc9SMatthew G. Knepley   PetscFunctionBegin;
174532d2bbc9SMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
174632d2bbc9SMatthew G. Knepley   if (update) PetscValidFunction(update, 3);
174732d2bbc9SMatthew G. Knepley   if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f);
174832d2bbc9SMatthew G. Knepley   ierr = PetscDSEnlarge_Static(prob, f+1);CHKERRQ(ierr);
174932d2bbc9SMatthew G. Knepley   prob->update[f] = update;
175032d2bbc9SMatthew G. Knepley   PetscFunctionReturn(0);
175132d2bbc9SMatthew G. Knepley }
175232d2bbc9SMatthew G. Knepley 
17530c2f2876SMatthew G. Knepley PetscErrorCode PetscDSGetContext(PetscDS prob, PetscInt f, void **ctx)
17540c2f2876SMatthew G. Knepley {
17550c2f2876SMatthew G. Knepley   PetscFunctionBegin;
17560c2f2876SMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
17570c2f2876SMatthew G. Knepley   if ((f < 0) || (f >= prob->Nf)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be in [0, %d)", f, prob->Nf);
17580c2f2876SMatthew G. Knepley   PetscValidPointer(ctx, 3);
17590c2f2876SMatthew G. Knepley   *ctx = prob->ctx[f];
17600c2f2876SMatthew G. Knepley   PetscFunctionReturn(0);
17610c2f2876SMatthew G. Knepley }
17620c2f2876SMatthew G. Knepley 
17630c2f2876SMatthew G. Knepley PetscErrorCode PetscDSSetContext(PetscDS prob, PetscInt f, void *ctx)
17640c2f2876SMatthew G. Knepley {
17650c2f2876SMatthew G. Knepley   PetscErrorCode ierr;
17660c2f2876SMatthew G. Knepley 
17670c2f2876SMatthew G. Knepley   PetscFunctionBegin;
17680c2f2876SMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
17690c2f2876SMatthew G. Knepley   if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f);
17700c2f2876SMatthew G. Knepley   ierr = PetscDSEnlarge_Static(prob, f+1);CHKERRQ(ierr);
17710c2f2876SMatthew G. Knepley   prob->ctx[f] = ctx;
17720c2f2876SMatthew G. Knepley   PetscFunctionReturn(0);
17730c2f2876SMatthew G. Knepley }
17740c2f2876SMatthew G. Knepley 
1775194d53e6SMatthew G. Knepley /*@C
1776194d53e6SMatthew G. Knepley   PetscDSGetBdResidual - Get the pointwise boundary residual function for a given test field
1777194d53e6SMatthew G. Knepley 
1778194d53e6SMatthew G. Knepley   Not collective
1779194d53e6SMatthew G. Knepley 
1780194d53e6SMatthew G. Knepley   Input Parameters:
1781194d53e6SMatthew G. Knepley + prob - The PetscDS
1782194d53e6SMatthew G. Knepley - f    - The test field number
1783194d53e6SMatthew G. Knepley 
1784194d53e6SMatthew G. Knepley   Output Parameters:
1785194d53e6SMatthew G. Knepley + f0 - boundary integrand for the test function term
1786194d53e6SMatthew G. Knepley - f1 - boundary integrand for the test function gradient term
1787194d53e6SMatthew G. Knepley 
1788194d53e6SMatthew G. Knepley   Note: We are using a first order FEM model for the weak form:
1789194d53e6SMatthew G. Knepley 
1790194d53e6SMatthew G. Knepley   \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
1791194d53e6SMatthew G. Knepley 
1792194d53e6SMatthew G. Knepley The calling sequence for the callbacks f0 and f1 is given by:
1793194d53e6SMatthew G. Knepley 
179430b9ff8bSMatthew G. Knepley $ f0(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1795194d53e6SMatthew G. Knepley $    const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1796194d53e6SMatthew G. Knepley $    const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
179730b9ff8bSMatthew G. Knepley $    PetscReal t, const PetscReal x[], const PetscReal n[], PetscScalar f0[])
1798194d53e6SMatthew G. Knepley 
1799194d53e6SMatthew G. Knepley + dim - the spatial dimension
1800194d53e6SMatthew G. Knepley . Nf - the number of fields
1801194d53e6SMatthew G. Knepley . uOff - the offset into u[] and u_t[] for each field
1802194d53e6SMatthew G. Knepley . uOff_x - the offset into u_x[] for each field
1803194d53e6SMatthew G. Knepley . u - each field evaluated at the current point
1804194d53e6SMatthew G. Knepley . u_t - the time derivative of each field evaluated at the current point
1805194d53e6SMatthew G. Knepley . u_x - the gradient of each field evaluated at the current point
1806194d53e6SMatthew G. Knepley . aOff - the offset into a[] and a_t[] for each auxiliary field
1807194d53e6SMatthew G. Knepley . aOff_x - the offset into a_x[] for each auxiliary field
1808194d53e6SMatthew G. Knepley . a - each auxiliary field evaluated at the current point
1809194d53e6SMatthew G. Knepley . a_t - the time derivative of each auxiliary field evaluated at the current point
1810194d53e6SMatthew G. Knepley . a_x - the gradient of auxiliary each field evaluated at the current point
1811194d53e6SMatthew G. Knepley . t - current time
1812194d53e6SMatthew G. Knepley . x - coordinates of the current point
1813194d53e6SMatthew G. Knepley . n - unit normal at the current point
181497b6e6e8SMatthew G. Knepley . numConstants - number of constant parameters
181597b6e6e8SMatthew G. Knepley . constants - constant parameters
1816194d53e6SMatthew G. Knepley - f0 - output values at the current point
1817194d53e6SMatthew G. Knepley 
1818194d53e6SMatthew G. Knepley   Level: intermediate
1819194d53e6SMatthew G. Knepley 
1820194d53e6SMatthew G. Knepley .seealso: PetscDSSetBdResidual()
1821194d53e6SMatthew G. Knepley @*/
18222764a2aaSMatthew G. Knepley PetscErrorCode PetscDSGetBdResidual(PetscDS prob, PetscInt f,
182330b9ff8bSMatthew G. Knepley                                     void (**f0)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1824194d53e6SMatthew G. Knepley                                                 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1825194d53e6SMatthew G. Knepley                                                 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
182697b6e6e8SMatthew G. Knepley                                                 PetscReal t, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]),
182730b9ff8bSMatthew G. Knepley                                     void (**f1)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1828194d53e6SMatthew G. Knepley                                                 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1829194d53e6SMatthew G. Knepley                                                 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
183097b6e6e8SMatthew G. Knepley                                                 PetscReal t, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]))
18312764a2aaSMatthew G. Knepley {
18322764a2aaSMatthew G. Knepley   PetscFunctionBegin;
18332764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
18342764a2aaSMatthew G. Knepley   if ((f < 0) || (f >= prob->Nf)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be in [0, %d)", f, prob->Nf);
18352764a2aaSMatthew G. Knepley   if (f0) {PetscValidPointer(f0, 3); *f0 = prob->fBd[f*2+0];}
18362764a2aaSMatthew G. Knepley   if (f1) {PetscValidPointer(f1, 4); *f1 = prob->fBd[f*2+1];}
18372764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
18382764a2aaSMatthew G. Knepley }
18392764a2aaSMatthew G. Knepley 
1840194d53e6SMatthew G. Knepley /*@C
1841194d53e6SMatthew G. Knepley   PetscDSSetBdResidual - Get the pointwise boundary residual function for a given test field
1842194d53e6SMatthew G. Knepley 
1843194d53e6SMatthew G. Knepley   Not collective
1844194d53e6SMatthew G. Knepley 
1845194d53e6SMatthew G. Knepley   Input Parameters:
1846194d53e6SMatthew G. Knepley + prob - The PetscDS
1847194d53e6SMatthew G. Knepley . f    - The test field number
1848194d53e6SMatthew G. Knepley . f0 - boundary integrand for the test function term
1849194d53e6SMatthew G. Knepley - f1 - boundary integrand for the test function gradient term
1850194d53e6SMatthew G. Knepley 
1851194d53e6SMatthew G. Knepley   Note: We are using a first order FEM model for the weak form:
1852194d53e6SMatthew G. Knepley 
1853194d53e6SMatthew G. Knepley   \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
1854194d53e6SMatthew G. Knepley 
1855194d53e6SMatthew G. Knepley The calling sequence for the callbacks f0 and f1 is given by:
1856194d53e6SMatthew G. Knepley 
185730b9ff8bSMatthew G. Knepley $ f0(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1858194d53e6SMatthew G. Knepley $    const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1859194d53e6SMatthew G. Knepley $    const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
186030b9ff8bSMatthew G. Knepley $    PetscReal t, const PetscReal x[], const PetscReal n[], PetscScalar f0[])
1861194d53e6SMatthew G. Knepley 
1862194d53e6SMatthew G. Knepley + dim - the spatial dimension
1863194d53e6SMatthew G. Knepley . Nf - the number of fields
1864194d53e6SMatthew G. Knepley . uOff - the offset into u[] and u_t[] for each field
1865194d53e6SMatthew G. Knepley . uOff_x - the offset into u_x[] for each field
1866194d53e6SMatthew G. Knepley . u - each field evaluated at the current point
1867194d53e6SMatthew G. Knepley . u_t - the time derivative of each field evaluated at the current point
1868194d53e6SMatthew G. Knepley . u_x - the gradient of each field evaluated at the current point
1869194d53e6SMatthew G. Knepley . aOff - the offset into a[] and a_t[] for each auxiliary field
1870194d53e6SMatthew G. Knepley . aOff_x - the offset into a_x[] for each auxiliary field
1871194d53e6SMatthew G. Knepley . a - each auxiliary field evaluated at the current point
1872194d53e6SMatthew G. Knepley . a_t - the time derivative of each auxiliary field evaluated at the current point
1873194d53e6SMatthew G. Knepley . a_x - the gradient of auxiliary each field evaluated at the current point
1874194d53e6SMatthew G. Knepley . t - current time
1875194d53e6SMatthew G. Knepley . x - coordinates of the current point
1876194d53e6SMatthew G. Knepley . n - unit normal at the current point
187797b6e6e8SMatthew G. Knepley . numConstants - number of constant parameters
187897b6e6e8SMatthew G. Knepley . constants - constant parameters
1879194d53e6SMatthew G. Knepley - f0 - output values at the current point
1880194d53e6SMatthew G. Knepley 
1881194d53e6SMatthew G. Knepley   Level: intermediate
1882194d53e6SMatthew G. Knepley 
1883194d53e6SMatthew G. Knepley .seealso: PetscDSGetBdResidual()
1884194d53e6SMatthew G. Knepley @*/
18852764a2aaSMatthew G. Knepley PetscErrorCode PetscDSSetBdResidual(PetscDS prob, PetscInt f,
188630b9ff8bSMatthew G. Knepley                                     void (*f0)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1887194d53e6SMatthew G. Knepley                                                const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1888194d53e6SMatthew G. Knepley                                                const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
188997b6e6e8SMatthew G. Knepley                                                PetscReal t, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]),
189030b9ff8bSMatthew G. Knepley                                     void (*f1)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1891194d53e6SMatthew G. Knepley                                                const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1892194d53e6SMatthew G. Knepley                                                const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
189397b6e6e8SMatthew G. Knepley                                                PetscReal t, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]))
18942764a2aaSMatthew G. Knepley {
18952764a2aaSMatthew G. Knepley   PetscErrorCode ierr;
18962764a2aaSMatthew G. Knepley 
18972764a2aaSMatthew G. Knepley   PetscFunctionBegin;
18982764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
18992764a2aaSMatthew G. Knepley   if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f);
19002764a2aaSMatthew G. Knepley   ierr = PetscDSEnlarge_Static(prob, f+1);CHKERRQ(ierr);
19012764a2aaSMatthew G. Knepley   if (f0) {PetscValidFunction(f0, 3); prob->fBd[f*2+0] = f0;}
19022764a2aaSMatthew G. Knepley   if (f1) {PetscValidFunction(f1, 4); prob->fBd[f*2+1] = f1;}
19032764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
19042764a2aaSMatthew G. Knepley }
19052764a2aaSMatthew G. Knepley 
1906194d53e6SMatthew G. Knepley /*@C
1907194d53e6SMatthew G. Knepley   PetscDSGetBdJacobian - Get the pointwise boundary Jacobian function for given test and basis field
1908194d53e6SMatthew G. Knepley 
1909194d53e6SMatthew G. Knepley   Not collective
1910194d53e6SMatthew G. Knepley 
1911194d53e6SMatthew G. Knepley   Input Parameters:
1912194d53e6SMatthew G. Knepley + prob - The PetscDS
1913194d53e6SMatthew G. Knepley . f    - The test field number
1914194d53e6SMatthew G. Knepley - g    - The field number
1915194d53e6SMatthew G. Knepley 
1916194d53e6SMatthew G. Knepley   Output Parameters:
1917194d53e6SMatthew G. Knepley + g0 - integrand for the test and basis function term
1918194d53e6SMatthew G. Knepley . g1 - integrand for the test function and basis function gradient term
1919194d53e6SMatthew G. Knepley . g2 - integrand for the test function gradient and basis function term
1920194d53e6SMatthew G. Knepley - g3 - integrand for the test function gradient and basis function gradient term
1921194d53e6SMatthew G. Knepley 
1922194d53e6SMatthew G. Knepley   Note: We are using a first order FEM model for the weak form:
1923194d53e6SMatthew G. Knepley 
1924194d53e6SMatthew G. Knepley   \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
1925194d53e6SMatthew G. Knepley 
1926194d53e6SMatthew G. Knepley The calling sequence for the callbacks g0, g1, g2 and g3 is given by:
1927194d53e6SMatthew G. Knepley 
192830b9ff8bSMatthew G. Knepley $ g0(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1929194d53e6SMatthew G. Knepley $    const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1930194d53e6SMatthew G. Knepley $    const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
193130b9ff8bSMatthew G. Knepley $    PetscReal t, const PetscReal x[], const PetscReal n[], PetscScalar g0[])
1932194d53e6SMatthew G. Knepley 
1933194d53e6SMatthew G. Knepley + dim - the spatial dimension
1934194d53e6SMatthew G. Knepley . Nf - the number of fields
1935194d53e6SMatthew G. Knepley . uOff - the offset into u[] and u_t[] for each field
1936194d53e6SMatthew G. Knepley . uOff_x - the offset into u_x[] for each field
1937194d53e6SMatthew G. Knepley . u - each field evaluated at the current point
1938194d53e6SMatthew G. Knepley . u_t - the time derivative of each field evaluated at the current point
1939194d53e6SMatthew G. Knepley . u_x - the gradient of each field evaluated at the current point
1940194d53e6SMatthew G. Knepley . aOff - the offset into a[] and a_t[] for each auxiliary field
1941194d53e6SMatthew G. Knepley . aOff_x - the offset into a_x[] for each auxiliary field
1942194d53e6SMatthew G. Knepley . a - each auxiliary field evaluated at the current point
1943194d53e6SMatthew G. Knepley . a_t - the time derivative of each auxiliary field evaluated at the current point
1944194d53e6SMatthew G. Knepley . a_x - the gradient of auxiliary each field evaluated at the current point
1945194d53e6SMatthew G. Knepley . t - current time
19462aa1fc23SMatthew G. Knepley . u_tShift - the multiplier a for dF/dU_t
1947194d53e6SMatthew G. Knepley . x - coordinates of the current point
1948194d53e6SMatthew G. Knepley . n - normal at the current point
194997b6e6e8SMatthew G. Knepley . numConstants - number of constant parameters
195097b6e6e8SMatthew G. Knepley . constants - constant parameters
1951194d53e6SMatthew G. Knepley - g0 - output values at the current point
1952194d53e6SMatthew G. Knepley 
1953194d53e6SMatthew G. Knepley   Level: intermediate
1954194d53e6SMatthew G. Knepley 
1955194d53e6SMatthew G. Knepley .seealso: PetscDSSetBdJacobian()
1956194d53e6SMatthew G. Knepley @*/
19572764a2aaSMatthew G. Knepley PetscErrorCode PetscDSGetBdJacobian(PetscDS prob, PetscInt f, PetscInt g,
195830b9ff8bSMatthew G. Knepley                                     void (**g0)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1959194d53e6SMatthew G. Knepley                                                 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1960194d53e6SMatthew G. Knepley                                                 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
196197b6e6e8SMatthew G. Knepley                                                 PetscReal t, PetscReal u_tShift, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]),
196230b9ff8bSMatthew G. Knepley                                     void (**g1)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1963194d53e6SMatthew G. Knepley                                                 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1964194d53e6SMatthew G. Knepley                                                 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
196597b6e6e8SMatthew G. Knepley                                                 PetscReal t, PetscReal u_tShift, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[]),
196630b9ff8bSMatthew G. Knepley                                     void (**g2)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1967194d53e6SMatthew G. Knepley                                                 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1968194d53e6SMatthew G. Knepley                                                 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
196997b6e6e8SMatthew G. Knepley                                                 PetscReal t, PetscReal u_tShift, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g2[]),
197030b9ff8bSMatthew G. Knepley                                     void (**g3)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1971194d53e6SMatthew G. Knepley                                                 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1972194d53e6SMatthew G. Knepley                                                 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
197397b6e6e8SMatthew G. Knepley                                                 PetscReal t, PetscReal u_tShift, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[]))
19742764a2aaSMatthew G. Knepley {
19752764a2aaSMatthew G. Knepley   PetscFunctionBegin;
19762764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
19772764a2aaSMatthew G. Knepley   if ((f < 0) || (f >= prob->Nf)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be in [0, %d)", f, prob->Nf);
19782764a2aaSMatthew G. Knepley   if ((g < 0) || (g >= prob->Nf)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be in [0, %d)", g, prob->Nf);
19792764a2aaSMatthew G. Knepley   if (g0) {PetscValidPointer(g0, 4); *g0 = prob->gBd[(f*prob->Nf + g)*4+0];}
19802764a2aaSMatthew G. Knepley   if (g1) {PetscValidPointer(g1, 5); *g1 = prob->gBd[(f*prob->Nf + g)*4+1];}
19812764a2aaSMatthew G. Knepley   if (g2) {PetscValidPointer(g2, 6); *g2 = prob->gBd[(f*prob->Nf + g)*4+2];}
19822764a2aaSMatthew G. Knepley   if (g3) {PetscValidPointer(g3, 7); *g3 = prob->gBd[(f*prob->Nf + g)*4+3];}
19832764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
19842764a2aaSMatthew G. Knepley }
19852764a2aaSMatthew G. Knepley 
1986194d53e6SMatthew G. Knepley /*@C
1987194d53e6SMatthew G. Knepley   PetscDSSetBdJacobian - Set the pointwise boundary Jacobian function for given test and basis field
1988194d53e6SMatthew G. Knepley 
1989194d53e6SMatthew G. Knepley   Not collective
1990194d53e6SMatthew G. Knepley 
1991194d53e6SMatthew G. Knepley   Input Parameters:
1992194d53e6SMatthew G. Knepley + prob - The PetscDS
1993194d53e6SMatthew G. Knepley . f    - The test field number
1994194d53e6SMatthew G. Knepley . g    - The field number
1995194d53e6SMatthew G. Knepley . g0 - integrand for the test and basis function term
1996194d53e6SMatthew G. Knepley . g1 - integrand for the test function and basis function gradient term
1997194d53e6SMatthew G. Knepley . g2 - integrand for the test function gradient and basis function term
1998194d53e6SMatthew G. Knepley - g3 - integrand for the test function gradient and basis function gradient term
1999194d53e6SMatthew G. Knepley 
2000194d53e6SMatthew G. Knepley   Note: We are using a first order FEM model for the weak form:
2001194d53e6SMatthew G. Knepley 
2002194d53e6SMatthew G. Knepley   \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
2003194d53e6SMatthew G. Knepley 
2004194d53e6SMatthew G. Knepley The calling sequence for the callbacks g0, g1, g2 and g3 is given by:
2005194d53e6SMatthew G. Knepley 
200630b9ff8bSMatthew G. Knepley $ g0(PetscInt dim, PetscInt Nf, PetscInt NfAux,
2007194d53e6SMatthew G. Knepley $    const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
2008194d53e6SMatthew G. Knepley $    const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
200930b9ff8bSMatthew G. Knepley $    PetscReal t, const PetscReal x[], const PetscReal n[], PetscScalar g0[])
2010194d53e6SMatthew G. Knepley 
2011194d53e6SMatthew G. Knepley + dim - the spatial dimension
2012194d53e6SMatthew G. Knepley . Nf - the number of fields
2013194d53e6SMatthew G. Knepley . uOff - the offset into u[] and u_t[] for each field
2014194d53e6SMatthew G. Knepley . uOff_x - the offset into u_x[] for each field
2015194d53e6SMatthew G. Knepley . u - each field evaluated at the current point
2016194d53e6SMatthew G. Knepley . u_t - the time derivative of each field evaluated at the current point
2017194d53e6SMatthew G. Knepley . u_x - the gradient of each field evaluated at the current point
2018194d53e6SMatthew G. Knepley . aOff - the offset into a[] and a_t[] for each auxiliary field
2019194d53e6SMatthew G. Knepley . aOff_x - the offset into a_x[] for each auxiliary field
2020194d53e6SMatthew G. Knepley . a - each auxiliary field evaluated at the current point
2021194d53e6SMatthew G. Knepley . a_t - the time derivative of each auxiliary field evaluated at the current point
2022194d53e6SMatthew G. Knepley . a_x - the gradient of auxiliary each field evaluated at the current point
2023194d53e6SMatthew G. Knepley . t - current time
20242aa1fc23SMatthew G. Knepley . u_tShift - the multiplier a for dF/dU_t
2025194d53e6SMatthew G. Knepley . x - coordinates of the current point
2026194d53e6SMatthew G. Knepley . n - normal at the current point
202797b6e6e8SMatthew G. Knepley . numConstants - number of constant parameters
202897b6e6e8SMatthew G. Knepley . constants - constant parameters
2029194d53e6SMatthew G. Knepley - g0 - output values at the current point
2030194d53e6SMatthew G. Knepley 
2031194d53e6SMatthew G. Knepley   Level: intermediate
2032194d53e6SMatthew G. Knepley 
2033194d53e6SMatthew G. Knepley .seealso: PetscDSGetBdJacobian()
2034194d53e6SMatthew G. Knepley @*/
20352764a2aaSMatthew G. Knepley PetscErrorCode PetscDSSetBdJacobian(PetscDS prob, PetscInt f, PetscInt g,
203630b9ff8bSMatthew G. Knepley                                     void (*g0)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
2037194d53e6SMatthew G. Knepley                                                const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
2038194d53e6SMatthew G. Knepley                                                const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
203997b6e6e8SMatthew G. Knepley                                                PetscReal t, PetscReal u_tShift, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]),
204030b9ff8bSMatthew G. Knepley                                     void (*g1)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
2041194d53e6SMatthew G. Knepley                                                const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
2042194d53e6SMatthew G. Knepley                                                const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
204397b6e6e8SMatthew G. Knepley                                                PetscReal t, PetscReal u_tShift, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[]),
204430b9ff8bSMatthew G. Knepley                                     void (*g2)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
2045194d53e6SMatthew G. Knepley                                                const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
2046194d53e6SMatthew G. Knepley                                                const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
204797b6e6e8SMatthew G. Knepley                                                PetscReal t, PetscReal u_tShift, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g2[]),
204830b9ff8bSMatthew G. Knepley                                     void (*g3)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
2049194d53e6SMatthew G. Knepley                                                const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
2050194d53e6SMatthew G. Knepley                                                const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
205197b6e6e8SMatthew G. Knepley                                                PetscReal t, PetscReal u_tShift, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[]))
20522764a2aaSMatthew G. Knepley {
20532764a2aaSMatthew G. Knepley   PetscErrorCode ierr;
20542764a2aaSMatthew G. Knepley 
20552764a2aaSMatthew G. Knepley   PetscFunctionBegin;
20562764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
20572764a2aaSMatthew G. Knepley   if (g0) PetscValidFunction(g0, 4);
20582764a2aaSMatthew G. Knepley   if (g1) PetscValidFunction(g1, 5);
20592764a2aaSMatthew G. Knepley   if (g2) PetscValidFunction(g2, 6);
20602764a2aaSMatthew G. Knepley   if (g3) PetscValidFunction(g3, 7);
20612764a2aaSMatthew G. Knepley   if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f);
20622764a2aaSMatthew G. Knepley   if (g < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", g);
20632764a2aaSMatthew G. Knepley   ierr = PetscDSEnlarge_Static(prob, PetscMax(f, g)+1);CHKERRQ(ierr);
20642764a2aaSMatthew G. Knepley   prob->gBd[(f*prob->Nf + g)*4+0] = g0;
20652764a2aaSMatthew G. Knepley   prob->gBd[(f*prob->Nf + g)*4+1] = g1;
20662764a2aaSMatthew G. Knepley   prob->gBd[(f*prob->Nf + g)*4+2] = g2;
20672764a2aaSMatthew G. Knepley   prob->gBd[(f*prob->Nf + g)*4+3] = g3;
20682764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
20692764a2aaSMatthew G. Knepley }
20702764a2aaSMatthew G. Knepley 
20710d3e9b51SMatthew G. Knepley /*@C
207297b6e6e8SMatthew G. Knepley   PetscDSGetConstants - Returns the array of constants passed to point functions
207397b6e6e8SMatthew G. Knepley 
207497b6e6e8SMatthew G. Knepley   Not collective
207597b6e6e8SMatthew G. Knepley 
207697b6e6e8SMatthew G. Knepley   Input Parameter:
207797b6e6e8SMatthew G. Knepley . prob - The PetscDS object
207897b6e6e8SMatthew G. Knepley 
207997b6e6e8SMatthew G. Knepley   Output Parameters:
208097b6e6e8SMatthew G. Knepley + numConstants - The number of constants
208197b6e6e8SMatthew G. Knepley - constants    - The array of constants, NULL if there are none
208297b6e6e8SMatthew G. Knepley 
208397b6e6e8SMatthew G. Knepley   Level: intermediate
208497b6e6e8SMatthew G. Knepley 
208597b6e6e8SMatthew G. Knepley .seealso: PetscDSSetConstants(), PetscDSCreate()
208697b6e6e8SMatthew G. Knepley @*/
208797b6e6e8SMatthew G. Knepley PetscErrorCode PetscDSGetConstants(PetscDS prob, PetscInt *numConstants, const PetscScalar *constants[])
208897b6e6e8SMatthew G. Knepley {
208997b6e6e8SMatthew G. Knepley   PetscFunctionBegin;
209097b6e6e8SMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
209197b6e6e8SMatthew G. Knepley   if (numConstants) {PetscValidPointer(numConstants, 2); *numConstants = prob->numConstants;}
209297b6e6e8SMatthew G. Knepley   if (constants)    {PetscValidPointer(constants, 3);    *constants    = prob->constants;}
209397b6e6e8SMatthew G. Knepley   PetscFunctionReturn(0);
209497b6e6e8SMatthew G. Knepley }
209597b6e6e8SMatthew G. Knepley 
20960d3e9b51SMatthew G. Knepley /*@C
209797b6e6e8SMatthew G. Knepley   PetscDSSetConstants - Set the array of constants passed to point functions
209897b6e6e8SMatthew G. Knepley 
209997b6e6e8SMatthew G. Knepley   Not collective
210097b6e6e8SMatthew G. Knepley 
210197b6e6e8SMatthew G. Knepley   Input Parameters:
210297b6e6e8SMatthew G. Knepley + prob         - The PetscDS object
210397b6e6e8SMatthew G. Knepley . numConstants - The number of constants
210497b6e6e8SMatthew G. Knepley - constants    - The array of constants, NULL if there are none
210597b6e6e8SMatthew G. Knepley 
210697b6e6e8SMatthew G. Knepley   Level: intermediate
210797b6e6e8SMatthew G. Knepley 
210897b6e6e8SMatthew G. Knepley .seealso: PetscDSGetConstants(), PetscDSCreate()
210997b6e6e8SMatthew G. Knepley @*/
211097b6e6e8SMatthew G. Knepley PetscErrorCode PetscDSSetConstants(PetscDS prob, PetscInt numConstants, PetscScalar constants[])
211197b6e6e8SMatthew G. Knepley {
211297b6e6e8SMatthew G. Knepley   PetscErrorCode ierr;
211397b6e6e8SMatthew G. Knepley 
211497b6e6e8SMatthew G. Knepley   PetscFunctionBegin;
211597b6e6e8SMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
211697b6e6e8SMatthew G. Knepley   if (numConstants != prob->numConstants) {
211797b6e6e8SMatthew G. Knepley     ierr = PetscFree(prob->constants);CHKERRQ(ierr);
211897b6e6e8SMatthew G. Knepley     prob->constants = NULL;
211997b6e6e8SMatthew G. Knepley   }
212097b6e6e8SMatthew G. Knepley   prob->numConstants = numConstants;
212197b6e6e8SMatthew G. Knepley   if (prob->numConstants) {
212297b6e6e8SMatthew G. Knepley     PetscValidPointer(constants, 3);
212397b6e6e8SMatthew G. Knepley     ierr = PetscMalloc1(prob->numConstants, &prob->constants);CHKERRQ(ierr);
212497b6e6e8SMatthew G. Knepley     ierr = PetscMemcpy(prob->constants, constants, prob->numConstants * sizeof(PetscScalar));CHKERRQ(ierr);
212597b6e6e8SMatthew G. Knepley   }
212697b6e6e8SMatthew G. Knepley   PetscFunctionReturn(0);
212797b6e6e8SMatthew G. Knepley }
212897b6e6e8SMatthew G. Knepley 
21294cd1e086SMatthew G. Knepley /*@
21304cd1e086SMatthew G. Knepley   PetscDSGetFieldIndex - Returns the index of the given field
21314cd1e086SMatthew G. Knepley 
21324cd1e086SMatthew G. Knepley   Not collective
21334cd1e086SMatthew G. Knepley 
21344cd1e086SMatthew G. Knepley   Input Parameters:
21354cd1e086SMatthew G. Knepley + prob - The PetscDS object
21364cd1e086SMatthew G. Knepley - disc - The discretization object
21374cd1e086SMatthew G. Knepley 
21384cd1e086SMatthew G. Knepley   Output Parameter:
21394cd1e086SMatthew G. Knepley . f - The field number
21404cd1e086SMatthew G. Knepley 
21414cd1e086SMatthew G. Knepley   Level: beginner
21424cd1e086SMatthew G. Knepley 
2143f744cafaSSander Arens .seealso: PetscGetDiscretization(), PetscDSGetNumFields(), PetscDSCreate()
21444cd1e086SMatthew G. Knepley @*/
21454cd1e086SMatthew G. Knepley PetscErrorCode PetscDSGetFieldIndex(PetscDS prob, PetscObject disc, PetscInt *f)
21464cd1e086SMatthew G. Knepley {
21474cd1e086SMatthew G. Knepley   PetscInt g;
21484cd1e086SMatthew G. Knepley 
21494cd1e086SMatthew G. Knepley   PetscFunctionBegin;
21504cd1e086SMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
21514cd1e086SMatthew G. Knepley   PetscValidPointer(f, 3);
21524cd1e086SMatthew G. Knepley   *f = -1;
21534cd1e086SMatthew G. Knepley   for (g = 0; g < prob->Nf; ++g) {if (disc == prob->disc[g]) break;}
21544cd1e086SMatthew G. Knepley   if (g == prob->Nf) SETERRQ(PetscObjectComm((PetscObject) prob), PETSC_ERR_ARG_WRONG, "Field not found in PetscDS.");
21554cd1e086SMatthew G. Knepley   *f = g;
21564cd1e086SMatthew G. Knepley   PetscFunctionReturn(0);
21574cd1e086SMatthew G. Knepley }
21584cd1e086SMatthew G. Knepley 
21594cd1e086SMatthew G. Knepley /*@
21604cd1e086SMatthew G. Knepley   PetscDSGetFieldSize - Returns the size of the given field in the full space basis
21614cd1e086SMatthew G. Knepley 
21624cd1e086SMatthew G. Knepley   Not collective
21634cd1e086SMatthew G. Knepley 
21644cd1e086SMatthew G. Knepley   Input Parameters:
21654cd1e086SMatthew G. Knepley + prob - The PetscDS object
21664cd1e086SMatthew G. Knepley - f - The field number
21674cd1e086SMatthew G. Knepley 
21684cd1e086SMatthew G. Knepley   Output Parameter:
21694cd1e086SMatthew G. Knepley . size - The size
21704cd1e086SMatthew G. Knepley 
21714cd1e086SMatthew G. Knepley   Level: beginner
21724cd1e086SMatthew G. Knepley 
2173f744cafaSSander Arens .seealso: PetscDSGetFieldOffset(), PetscDSGetNumFields(), PetscDSCreate()
21744cd1e086SMatthew G. Knepley @*/
21754cd1e086SMatthew G. Knepley PetscErrorCode PetscDSGetFieldSize(PetscDS prob, PetscInt f, PetscInt *size)
21764cd1e086SMatthew G. Knepley {
21772166fd64SMatthew G. Knepley   PetscErrorCode ierr;
21782166fd64SMatthew G. Knepley 
21794cd1e086SMatthew G. Knepley   PetscFunctionBegin;
21804cd1e086SMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
21814cd1e086SMatthew G. Knepley   PetscValidPointer(size, 3);
21824cd1e086SMatthew G. Knepley   if ((f < 0) || (f >= prob->Nf)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be in [0, %d)", f, prob->Nf);
21832166fd64SMatthew G. Knepley   ierr = PetscDSSetUp(prob);CHKERRQ(ierr);
2184d4742ddaSMatthew G. Knepley   *size = prob->Nb[f];
21854cd1e086SMatthew G. Knepley   PetscFunctionReturn(0);
21864cd1e086SMatthew G. Knepley }
21874cd1e086SMatthew G. Knepley 
2188bc4ae4beSMatthew G. Knepley /*@
2189bc4ae4beSMatthew G. Knepley   PetscDSGetFieldOffset - Returns the offset of the given field in the full space basis
2190bc4ae4beSMatthew G. Knepley 
2191bc4ae4beSMatthew G. Knepley   Not collective
2192bc4ae4beSMatthew G. Knepley 
2193bc4ae4beSMatthew G. Knepley   Input Parameters:
2194bc4ae4beSMatthew G. Knepley + prob - The PetscDS object
2195bc4ae4beSMatthew G. Knepley - f - The field number
2196bc4ae4beSMatthew G. Knepley 
2197bc4ae4beSMatthew G. Knepley   Output Parameter:
2198bc4ae4beSMatthew G. Knepley . off - The offset
2199bc4ae4beSMatthew G. Knepley 
2200bc4ae4beSMatthew G. Knepley   Level: beginner
2201bc4ae4beSMatthew G. Knepley 
2202f744cafaSSander Arens .seealso: PetscDSGetFieldSize(), PetscDSGetNumFields(), PetscDSCreate()
2203bc4ae4beSMatthew G. Knepley @*/
22042764a2aaSMatthew G. Knepley PetscErrorCode PetscDSGetFieldOffset(PetscDS prob, PetscInt f, PetscInt *off)
22052764a2aaSMatthew G. Knepley {
22064cd1e086SMatthew G. Knepley   PetscInt       size, g;
22072764a2aaSMatthew G. Knepley   PetscErrorCode ierr;
22082764a2aaSMatthew G. Knepley 
22092764a2aaSMatthew G. Knepley   PetscFunctionBegin;
22102764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
22112764a2aaSMatthew G. Knepley   PetscValidPointer(off, 3);
22122764a2aaSMatthew G. Knepley   if ((f < 0) || (f >= prob->Nf)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be in [0, %d)", f, prob->Nf);
22132764a2aaSMatthew G. Knepley   *off = 0;
22142764a2aaSMatthew G. Knepley   for (g = 0; g < f; ++g) {
22154cd1e086SMatthew G. Knepley     ierr = PetscDSGetFieldSize(prob, g, &size);CHKERRQ(ierr);
22164cd1e086SMatthew G. Knepley     *off += size;
22172764a2aaSMatthew G. Knepley   }
22182764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
22192764a2aaSMatthew G. Knepley }
22202764a2aaSMatthew G. Knepley 
2221bc4ae4beSMatthew G. Knepley /*@
222247e57110SSander Arens   PetscDSGetDimensions - Returns the size of the approximation space for each field on an evaluation point
2223bc4ae4beSMatthew G. Knepley 
2224bc4ae4beSMatthew G. Knepley   Not collective
2225bc4ae4beSMatthew G. Knepley 
222647e57110SSander Arens   Input Parameter:
222747e57110SSander Arens . prob - The PetscDS object
2228bc4ae4beSMatthew G. Knepley 
2229bc4ae4beSMatthew G. Knepley   Output Parameter:
223047e57110SSander Arens . dimensions - The number of dimensions
2231bc4ae4beSMatthew G. Knepley 
2232bc4ae4beSMatthew G. Knepley   Level: beginner
2233bc4ae4beSMatthew G. Knepley 
223447e57110SSander Arens .seealso: PetscDSGetComponentOffsets(), PetscDSGetNumFields(), PetscDSCreate()
2235bc4ae4beSMatthew G. Knepley @*/
223647e57110SSander Arens PetscErrorCode PetscDSGetDimensions(PetscDS prob, PetscInt *dimensions[])
22372764a2aaSMatthew G. Knepley {
22382764a2aaSMatthew G. Knepley   PetscErrorCode ierr;
22392764a2aaSMatthew G. Knepley 
22402764a2aaSMatthew G. Knepley   PetscFunctionBegin;
22412764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
224247e57110SSander Arens   ierr = PetscDSSetUp(prob);CHKERRQ(ierr);
224347e57110SSander Arens   PetscValidPointer(dimensions, 2);
224447e57110SSander Arens   *dimensions = prob->Nb;
224547e57110SSander Arens   PetscFunctionReturn(0);
22466ce16762SMatthew G. Knepley }
224747e57110SSander Arens 
224847e57110SSander Arens /*@
224947e57110SSander Arens   PetscDSGetComponents - Returns the number of components for each field on an evaluation point
225047e57110SSander Arens 
225147e57110SSander Arens   Not collective
225247e57110SSander Arens 
225347e57110SSander Arens   Input Parameter:
225447e57110SSander Arens . prob - The PetscDS object
225547e57110SSander Arens 
225647e57110SSander Arens   Output Parameter:
225747e57110SSander Arens . components - The number of components
225847e57110SSander Arens 
225947e57110SSander Arens   Level: beginner
226047e57110SSander Arens 
226147e57110SSander Arens .seealso: PetscDSGetComponentOffsets(), PetscDSGetNumFields(), PetscDSCreate()
226247e57110SSander Arens @*/
226347e57110SSander Arens PetscErrorCode PetscDSGetComponents(PetscDS prob, PetscInt *components[])
226447e57110SSander Arens {
226547e57110SSander Arens   PetscErrorCode ierr;
226647e57110SSander Arens 
226747e57110SSander Arens   PetscFunctionBegin;
226847e57110SSander Arens   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
226947e57110SSander Arens   ierr = PetscDSSetUp(prob);CHKERRQ(ierr);
227047e57110SSander Arens   PetscValidPointer(components, 2);
227147e57110SSander Arens   *components = prob->Nc;
22726ce16762SMatthew G. Knepley   PetscFunctionReturn(0);
22736ce16762SMatthew G. Knepley }
22746ce16762SMatthew G. Knepley 
22756ce16762SMatthew G. Knepley /*@
22766ce16762SMatthew G. Knepley   PetscDSGetComponentOffset - Returns the offset of the given field on an evaluation point
22776ce16762SMatthew G. Knepley 
22786ce16762SMatthew G. Knepley   Not collective
22796ce16762SMatthew G. Knepley 
22806ce16762SMatthew G. Knepley   Input Parameters:
22816ce16762SMatthew G. Knepley + prob - The PetscDS object
22826ce16762SMatthew G. Knepley - f - The field number
22836ce16762SMatthew G. Knepley 
22846ce16762SMatthew G. Knepley   Output Parameter:
22856ce16762SMatthew G. Knepley . off - The offset
22866ce16762SMatthew G. Knepley 
22876ce16762SMatthew G. Knepley   Level: beginner
22886ce16762SMatthew G. Knepley 
2289f744cafaSSander Arens .seealso: PetscDSGetNumFields(), PetscDSCreate()
22906ce16762SMatthew G. Knepley @*/
22916ce16762SMatthew G. Knepley PetscErrorCode PetscDSGetComponentOffset(PetscDS prob, PetscInt f, PetscInt *off)
22926ce16762SMatthew G. Knepley {
22936ce16762SMatthew G. Knepley   PetscFunctionBegin;
22946ce16762SMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
22956ce16762SMatthew G. Knepley   PetscValidPointer(off, 3);
22966ce16762SMatthew G. Knepley   if ((f < 0) || (f >= prob->Nf)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be in [0, %d)", f, prob->Nf);
229747e57110SSander Arens   *off = prob->off[f];
22982764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
22992764a2aaSMatthew G. Knepley }
23002764a2aaSMatthew G. Knepley 
2301194d53e6SMatthew G. Knepley /*@
2302194d53e6SMatthew G. Knepley   PetscDSGetComponentOffsets - Returns the offset of each field on an evaluation point
2303194d53e6SMatthew G. Knepley 
2304194d53e6SMatthew G. Knepley   Not collective
2305194d53e6SMatthew G. Knepley 
2306194d53e6SMatthew G. Knepley   Input Parameter:
2307194d53e6SMatthew G. Knepley . prob - The PetscDS object
2308194d53e6SMatthew G. Knepley 
2309194d53e6SMatthew G. Knepley   Output Parameter:
2310194d53e6SMatthew G. Knepley . offsets - The offsets
2311194d53e6SMatthew G. Knepley 
2312194d53e6SMatthew G. Knepley   Level: beginner
2313194d53e6SMatthew G. Knepley 
2314f744cafaSSander Arens .seealso: PetscDSGetNumFields(), PetscDSCreate()
2315194d53e6SMatthew G. Knepley @*/
2316194d53e6SMatthew G. Knepley PetscErrorCode PetscDSGetComponentOffsets(PetscDS prob, PetscInt *offsets[])
2317194d53e6SMatthew G. Knepley {
2318194d53e6SMatthew G. Knepley   PetscFunctionBegin;
2319194d53e6SMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
2320194d53e6SMatthew G. Knepley   PetscValidPointer(offsets, 2);
2321194d53e6SMatthew G. Knepley   *offsets = prob->off;
2322194d53e6SMatthew G. Knepley   PetscFunctionReturn(0);
2323194d53e6SMatthew G. Knepley }
2324194d53e6SMatthew G. Knepley 
2325194d53e6SMatthew G. Knepley /*@
2326194d53e6SMatthew G. Knepley   PetscDSGetComponentDerivativeOffsets - Returns the offset of each field derivative on an evaluation point
2327194d53e6SMatthew G. Knepley 
2328194d53e6SMatthew G. Knepley   Not collective
2329194d53e6SMatthew G. Knepley 
2330194d53e6SMatthew G. Knepley   Input Parameter:
2331194d53e6SMatthew G. Knepley . prob - The PetscDS object
2332194d53e6SMatthew G. Knepley 
2333194d53e6SMatthew G. Knepley   Output Parameter:
2334194d53e6SMatthew G. Knepley . offsets - The offsets
2335194d53e6SMatthew G. Knepley 
2336194d53e6SMatthew G. Knepley   Level: beginner
2337194d53e6SMatthew G. Knepley 
2338f744cafaSSander Arens .seealso: PetscDSGetNumFields(), PetscDSCreate()
2339194d53e6SMatthew G. Knepley @*/
2340194d53e6SMatthew G. Knepley PetscErrorCode PetscDSGetComponentDerivativeOffsets(PetscDS prob, PetscInt *offsets[])
2341194d53e6SMatthew G. Knepley {
2342194d53e6SMatthew G. Knepley   PetscFunctionBegin;
2343194d53e6SMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
2344194d53e6SMatthew G. Knepley   PetscValidPointer(offsets, 2);
2345194d53e6SMatthew G. Knepley   *offsets = prob->offDer;
2346194d53e6SMatthew G. Knepley   PetscFunctionReturn(0);
2347194d53e6SMatthew G. Knepley }
2348194d53e6SMatthew G. Knepley 
234968c9edb9SMatthew G. Knepley /*@C
235068c9edb9SMatthew G. Knepley   PetscDSGetTabulation - Return the basis tabulation at quadrature points for the volume discretization
235168c9edb9SMatthew G. Knepley 
235268c9edb9SMatthew G. Knepley   Not collective
235368c9edb9SMatthew G. Knepley 
235468c9edb9SMatthew G. Knepley   Input Parameter:
235568c9edb9SMatthew G. Knepley . prob - The PetscDS object
235668c9edb9SMatthew G. Knepley 
235768c9edb9SMatthew G. Knepley   Output Parameters:
235868c9edb9SMatthew G. Knepley + basis - The basis function tabulation at quadrature points
235968c9edb9SMatthew G. Knepley - basisDer - The basis function derivative tabulation at quadrature points
236068c9edb9SMatthew G. Knepley 
236168c9edb9SMatthew G. Knepley   Level: intermediate
236268c9edb9SMatthew G. Knepley 
2363f744cafaSSander Arens .seealso: PetscDSCreate()
236468c9edb9SMatthew G. Knepley @*/
23652764a2aaSMatthew G. Knepley PetscErrorCode PetscDSGetTabulation(PetscDS prob, PetscReal ***basis, PetscReal ***basisDer)
23662764a2aaSMatthew G. Knepley {
23672764a2aaSMatthew G. Knepley   PetscErrorCode ierr;
23682764a2aaSMatthew G. Knepley 
23692764a2aaSMatthew G. Knepley   PetscFunctionBegin;
23702764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
23712764a2aaSMatthew G. Knepley   ierr = PetscDSSetUp(prob);CHKERRQ(ierr);
23722764a2aaSMatthew G. Knepley   if (basis)    {PetscValidPointer(basis, 2);    *basis    = prob->basis;}
23732764a2aaSMatthew G. Knepley   if (basisDer) {PetscValidPointer(basisDer, 3); *basisDer = prob->basisDer;}
23742764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
23752764a2aaSMatthew G. Knepley }
23762764a2aaSMatthew G. Knepley 
237768c9edb9SMatthew G. Knepley /*@C
23784d0b9603SSander Arens   PetscDSGetFaceTabulation - Return the basis tabulation at quadrature points on the faces
237968c9edb9SMatthew G. Knepley 
238068c9edb9SMatthew G. Knepley   Not collective
238168c9edb9SMatthew G. Knepley 
238268c9edb9SMatthew G. Knepley   Input Parameter:
238368c9edb9SMatthew G. Knepley . prob - The PetscDS object
238468c9edb9SMatthew G. Knepley 
238568c9edb9SMatthew G. Knepley   Output Parameters:
23864d0b9603SSander Arens + basisFace - The basis function tabulation at quadrature points
23874d0b9603SSander Arens - basisDerFace - The basis function derivative tabulation at quadrature points
238868c9edb9SMatthew G. Knepley 
238968c9edb9SMatthew G. Knepley   Level: intermediate
239068c9edb9SMatthew G. Knepley 
239168c9edb9SMatthew G. Knepley .seealso: PetscDSGetTabulation(), PetscDSCreate()
239268c9edb9SMatthew G. Knepley @*/
23934d0b9603SSander Arens PetscErrorCode PetscDSGetFaceTabulation(PetscDS prob, PetscReal ***basis, PetscReal ***basisDer)
23942764a2aaSMatthew G. Knepley {
23952764a2aaSMatthew G. Knepley   PetscErrorCode ierr;
23962764a2aaSMatthew G. Knepley 
23972764a2aaSMatthew G. Knepley   PetscFunctionBegin;
23982764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
23992764a2aaSMatthew G. Knepley   ierr = PetscDSSetUp(prob);CHKERRQ(ierr);
24004d0b9603SSander Arens   if (basis)    {PetscValidPointer(basis, 2);    *basis    = prob->basisFace;}
24014d0b9603SSander Arens   if (basisDer) {PetscValidPointer(basisDer, 3); *basisDer = prob->basisDerFace;}
24022764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
24032764a2aaSMatthew G. Knepley }
24042764a2aaSMatthew G. Knepley 
24052764a2aaSMatthew G. Knepley PetscErrorCode PetscDSGetEvaluationArrays(PetscDS prob, PetscScalar **u, PetscScalar **u_t, PetscScalar **u_x)
24062764a2aaSMatthew G. Knepley {
24072764a2aaSMatthew G. Knepley   PetscErrorCode ierr;
24082764a2aaSMatthew G. Knepley 
24092764a2aaSMatthew G. Knepley   PetscFunctionBegin;
24102764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
24112764a2aaSMatthew G. Knepley   ierr = PetscDSSetUp(prob);CHKERRQ(ierr);
24122764a2aaSMatthew G. Knepley   if (u)   {PetscValidPointer(u, 2);   *u   = prob->u;}
24132764a2aaSMatthew G. Knepley   if (u_t) {PetscValidPointer(u_t, 3); *u_t = prob->u_t;}
24142764a2aaSMatthew G. Knepley   if (u_x) {PetscValidPointer(u_x, 4); *u_x = prob->u_x;}
24152764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
24162764a2aaSMatthew G. Knepley }
24172764a2aaSMatthew G. Knepley 
24182764a2aaSMatthew G. Knepley PetscErrorCode PetscDSGetWeakFormArrays(PetscDS prob, PetscScalar **f0, PetscScalar **f1, PetscScalar **g0, PetscScalar **g1, PetscScalar **g2, PetscScalar **g3)
24192764a2aaSMatthew G. Knepley {
24202764a2aaSMatthew G. Knepley   PetscErrorCode ierr;
24212764a2aaSMatthew G. Knepley 
24222764a2aaSMatthew G. Knepley   PetscFunctionBegin;
24232764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
24242764a2aaSMatthew G. Knepley   ierr = PetscDSSetUp(prob);CHKERRQ(ierr);
24252764a2aaSMatthew G. Knepley   if (f0) {PetscValidPointer(f0, 2); *f0 = prob->f0;}
24262764a2aaSMatthew G. Knepley   if (f1) {PetscValidPointer(f1, 3); *f1 = prob->f1;}
24272764a2aaSMatthew G. Knepley   if (g0) {PetscValidPointer(g0, 4); *g0 = prob->g0;}
24282764a2aaSMatthew G. Knepley   if (g1) {PetscValidPointer(g1, 5); *g1 = prob->g1;}
24292764a2aaSMatthew G. Knepley   if (g2) {PetscValidPointer(g2, 6); *g2 = prob->g2;}
24302764a2aaSMatthew G. Knepley   if (g3) {PetscValidPointer(g3, 7); *g3 = prob->g3;}
24312764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
24322764a2aaSMatthew G. Knepley }
24332764a2aaSMatthew G. Knepley 
24342764a2aaSMatthew G. Knepley PetscErrorCode PetscDSGetRefCoordArrays(PetscDS prob, PetscReal **x, PetscScalar **refSpaceDer)
24352764a2aaSMatthew G. Knepley {
24362764a2aaSMatthew G. Knepley   PetscErrorCode ierr;
24372764a2aaSMatthew G. Knepley 
24382764a2aaSMatthew G. Knepley   PetscFunctionBegin;
24392764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
24402764a2aaSMatthew G. Knepley   ierr = PetscDSSetUp(prob);CHKERRQ(ierr);
24412764a2aaSMatthew G. Knepley   if (x)           {PetscValidPointer(x, 2);           *x           = prob->x;}
24422764a2aaSMatthew G. Knepley   if (refSpaceDer) {PetscValidPointer(refSpaceDer, 3); *refSpaceDer = prob->refSpaceDer;}
24432764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
24442764a2aaSMatthew G. Knepley }
24452764a2aaSMatthew G. Knepley 
244658ebd649SToby Isaac /*@C
244758ebd649SToby Isaac   PetscDSAddBoundary - Add a boundary condition to the model
244858ebd649SToby Isaac 
244958ebd649SToby Isaac   Input Parameters:
245058ebd649SToby Isaac + ds          - The PetscDS object
24512d47a189SJulian Andrej . type        - The type of condition, e.g. DM_BC_ESSENTIAL/DM_BC_ESSENTIAL_FIELD (Dirichlet), or DM_BC_NATURAL (Neumann)
245258ebd649SToby Isaac . name        - The BC name
245358ebd649SToby Isaac . labelname   - The label defining constrained points
245458ebd649SToby Isaac . field       - The field to constrain
245558ebd649SToby Isaac . numcomps    - The number of constrained field components
245658ebd649SToby Isaac . comps       - An array of constrained component numbers
245758ebd649SToby Isaac . bcFunc      - A pointwise function giving boundary values
245858ebd649SToby Isaac . numids      - The number of DMLabel ids for constrained points
245958ebd649SToby Isaac . ids         - An array of ids for constrained points
246058ebd649SToby Isaac - ctx         - An optional user context for bcFunc
246158ebd649SToby Isaac 
246258ebd649SToby Isaac   Options Database Keys:
246358ebd649SToby Isaac + -bc_<boundary name> <num> - Overrides the boundary ids
246458ebd649SToby Isaac - -bc_<boundary name>_comp <num> - Overrides the boundary components
246558ebd649SToby Isaac 
246658ebd649SToby Isaac   Level: developer
246758ebd649SToby Isaac 
246858ebd649SToby Isaac .seealso: PetscDSGetBoundary()
246958ebd649SToby Isaac @*/
2470a30ec4eaSSatish Balay PetscErrorCode PetscDSAddBoundary(PetscDS ds, DMBoundaryConditionType type, const char name[], const char labelname[], PetscInt field, PetscInt numcomps, const PetscInt *comps, void (*bcFunc)(void), PetscInt numids, const PetscInt *ids, void *ctx)
247158ebd649SToby Isaac {
247258ebd649SToby Isaac   DSBoundary     b;
247358ebd649SToby Isaac   PetscErrorCode ierr;
247458ebd649SToby Isaac 
247558ebd649SToby Isaac   PetscFunctionBegin;
247658ebd649SToby Isaac   PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1);
247758ebd649SToby Isaac   ierr = PetscNew(&b);CHKERRQ(ierr);
247858ebd649SToby Isaac   ierr = PetscStrallocpy(name, (char **) &b->name);CHKERRQ(ierr);
247958ebd649SToby Isaac   ierr = PetscStrallocpy(labelname, (char **) &b->labelname);CHKERRQ(ierr);
248058ebd649SToby Isaac   ierr = PetscMalloc1(numcomps, &b->comps);CHKERRQ(ierr);
248158ebd649SToby Isaac   if (numcomps) {ierr = PetscMemcpy(b->comps, comps, numcomps*sizeof(PetscInt));CHKERRQ(ierr);}
248258ebd649SToby Isaac   ierr = PetscMalloc1(numids, &b->ids);CHKERRQ(ierr);
248358ebd649SToby Isaac   if (numids) {ierr = PetscMemcpy(b->ids, ids, numids*sizeof(PetscInt));CHKERRQ(ierr);}
2484f971fd6bSMatthew G. Knepley   b->type            = type;
248558ebd649SToby Isaac   b->field           = field;
248658ebd649SToby Isaac   b->numcomps        = numcomps;
248758ebd649SToby Isaac   b->func            = bcFunc;
248858ebd649SToby Isaac   b->numids          = numids;
248958ebd649SToby Isaac   b->ctx             = ctx;
249058ebd649SToby Isaac   b->next            = ds->boundary;
249158ebd649SToby Isaac   ds->boundary       = b;
249258ebd649SToby Isaac   PetscFunctionReturn(0);
249358ebd649SToby Isaac }
249458ebd649SToby Isaac 
249558ebd649SToby Isaac /*@
249658ebd649SToby Isaac   PetscDSGetNumBoundary - Get the number of registered BC
249758ebd649SToby Isaac 
249858ebd649SToby Isaac   Input Parameters:
249958ebd649SToby Isaac . ds - The PetscDS object
250058ebd649SToby Isaac 
250158ebd649SToby Isaac   Output Parameters:
250258ebd649SToby Isaac . numBd - The number of BC
250358ebd649SToby Isaac 
250458ebd649SToby Isaac   Level: intermediate
250558ebd649SToby Isaac 
250658ebd649SToby Isaac .seealso: PetscDSAddBoundary(), PetscDSGetBoundary()
250758ebd649SToby Isaac @*/
250858ebd649SToby Isaac PetscErrorCode PetscDSGetNumBoundary(PetscDS ds, PetscInt *numBd)
250958ebd649SToby Isaac {
251058ebd649SToby Isaac   DSBoundary b = ds->boundary;
251158ebd649SToby Isaac 
251258ebd649SToby Isaac   PetscFunctionBegin;
251358ebd649SToby Isaac   PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1);
251458ebd649SToby Isaac   PetscValidPointer(numBd, 2);
251558ebd649SToby Isaac   *numBd = 0;
251658ebd649SToby Isaac   while (b) {++(*numBd); b = b->next;}
251758ebd649SToby Isaac   PetscFunctionReturn(0);
251858ebd649SToby Isaac }
251958ebd649SToby Isaac 
252058ebd649SToby Isaac /*@C
252158ebd649SToby Isaac   PetscDSGetBoundary - Add a boundary condition to the model
252258ebd649SToby Isaac 
252358ebd649SToby Isaac   Input Parameters:
252458ebd649SToby Isaac + ds          - The PetscDS object
252558ebd649SToby Isaac - bd          - The BC number
252658ebd649SToby Isaac 
252758ebd649SToby Isaac   Output Parameters:
25282d47a189SJulian Andrej + type        - The type of condition, e.g. DM_BC_ESSENTIAL/DM_BC_ESSENTIAL_FIELD (Dirichlet), or DM_BC_NATURAL (Neumann)
252958ebd649SToby Isaac . name        - The BC name
253058ebd649SToby Isaac . labelname   - The label defining constrained points
253158ebd649SToby Isaac . field       - The field to constrain
253258ebd649SToby Isaac . numcomps    - The number of constrained field components
253358ebd649SToby Isaac . comps       - An array of constrained component numbers
253458ebd649SToby Isaac . bcFunc      - A pointwise function giving boundary values
253558ebd649SToby Isaac . numids      - The number of DMLabel ids for constrained points
253658ebd649SToby Isaac . ids         - An array of ids for constrained points
253758ebd649SToby Isaac - ctx         - An optional user context for bcFunc
253858ebd649SToby Isaac 
253958ebd649SToby Isaac   Options Database Keys:
254058ebd649SToby Isaac + -bc_<boundary name> <num> - Overrides the boundary ids
254158ebd649SToby Isaac - -bc_<boundary name>_comp <num> - Overrides the boundary components
254258ebd649SToby Isaac 
254358ebd649SToby Isaac   Level: developer
254458ebd649SToby Isaac 
254558ebd649SToby Isaac .seealso: PetscDSAddBoundary()
254658ebd649SToby Isaac @*/
2547a30ec4eaSSatish Balay PetscErrorCode PetscDSGetBoundary(PetscDS ds, PetscInt bd, DMBoundaryConditionType *type, const char **name, const char **labelname, PetscInt *field, PetscInt *numcomps, const PetscInt **comps, void (**func)(void), PetscInt *numids, const PetscInt **ids, void **ctx)
254858ebd649SToby Isaac {
254958ebd649SToby Isaac   DSBoundary b    = ds->boundary;
255058ebd649SToby Isaac   PetscInt   n    = 0;
255158ebd649SToby Isaac 
255258ebd649SToby Isaac   PetscFunctionBegin;
255358ebd649SToby Isaac   PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1);
255458ebd649SToby Isaac   while (b) {
255558ebd649SToby Isaac     if (n == bd) break;
255658ebd649SToby Isaac     b = b->next;
255758ebd649SToby Isaac     ++n;
255858ebd649SToby Isaac   }
255958ebd649SToby Isaac   if (!b) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Boundary %d is not in [0, %d)", bd, n);
2560f971fd6bSMatthew G. Knepley   if (type) {
2561f971fd6bSMatthew G. Knepley     PetscValidPointer(type, 3);
2562f971fd6bSMatthew G. Knepley     *type = b->type;
256358ebd649SToby Isaac   }
256458ebd649SToby Isaac   if (name) {
256558ebd649SToby Isaac     PetscValidPointer(name, 4);
256658ebd649SToby Isaac     *name = b->name;
256758ebd649SToby Isaac   }
256858ebd649SToby Isaac   if (labelname) {
256958ebd649SToby Isaac     PetscValidPointer(labelname, 5);
257058ebd649SToby Isaac     *labelname = b->labelname;
257158ebd649SToby Isaac   }
257258ebd649SToby Isaac   if (field) {
257358ebd649SToby Isaac     PetscValidPointer(field, 6);
257458ebd649SToby Isaac     *field = b->field;
257558ebd649SToby Isaac   }
257658ebd649SToby Isaac   if (numcomps) {
257758ebd649SToby Isaac     PetscValidPointer(numcomps, 7);
257858ebd649SToby Isaac     *numcomps = b->numcomps;
257958ebd649SToby Isaac   }
258058ebd649SToby Isaac   if (comps) {
258158ebd649SToby Isaac     PetscValidPointer(comps, 8);
258258ebd649SToby Isaac     *comps = b->comps;
258358ebd649SToby Isaac   }
258458ebd649SToby Isaac   if (func) {
258558ebd649SToby Isaac     PetscValidPointer(func, 9);
258658ebd649SToby Isaac     *func = b->func;
258758ebd649SToby Isaac   }
258858ebd649SToby Isaac   if (numids) {
258958ebd649SToby Isaac     PetscValidPointer(numids, 10);
259058ebd649SToby Isaac     *numids = b->numids;
259158ebd649SToby Isaac   }
259258ebd649SToby Isaac   if (ids) {
259358ebd649SToby Isaac     PetscValidPointer(ids, 11);
259458ebd649SToby Isaac     *ids = b->ids;
259558ebd649SToby Isaac   }
259658ebd649SToby Isaac   if (ctx) {
259758ebd649SToby Isaac     PetscValidPointer(ctx, 12);
259858ebd649SToby Isaac     *ctx = b->ctx;
259958ebd649SToby Isaac   }
260058ebd649SToby Isaac   PetscFunctionReturn(0);
260158ebd649SToby Isaac }
260258ebd649SToby Isaac 
2603dff059c6SToby Isaac PetscErrorCode PetscDSCopyBoundary(PetscDS probA, PetscDS probB)
2604dff059c6SToby Isaac {
2605dff059c6SToby Isaac   DSBoundary     b, next, *lastnext;
2606dff059c6SToby Isaac   PetscErrorCode ierr;
2607dff059c6SToby Isaac 
2608dff059c6SToby Isaac   PetscFunctionBegin;
2609dff059c6SToby Isaac   PetscValidHeaderSpecific(probA, PETSCDS_CLASSID, 1);
2610dff059c6SToby Isaac   PetscValidHeaderSpecific(probB, PETSCDS_CLASSID, 2);
2611dff059c6SToby Isaac   if (probA == probB) PetscFunctionReturn(0);
2612dff059c6SToby Isaac   next = probB->boundary;
2613dff059c6SToby Isaac   while (next) {
2614dff059c6SToby Isaac     DSBoundary b = next;
2615dff059c6SToby Isaac 
2616dff059c6SToby Isaac     next = b->next;
2617dff059c6SToby Isaac     ierr = PetscFree(b->comps);CHKERRQ(ierr);
2618dff059c6SToby Isaac     ierr = PetscFree(b->ids);CHKERRQ(ierr);
2619dff059c6SToby Isaac     ierr = PetscFree(b->name);CHKERRQ(ierr);
2620dff059c6SToby Isaac     ierr = PetscFree(b->labelname);CHKERRQ(ierr);
2621dff059c6SToby Isaac     ierr = PetscFree(b);CHKERRQ(ierr);
2622dff059c6SToby Isaac   }
2623dff059c6SToby Isaac   lastnext = &(probB->boundary);
2624dff059c6SToby Isaac   for (b = probA->boundary; b; b = b->next) {
2625dff059c6SToby Isaac     DSBoundary bNew;
2626dff059c6SToby Isaac 
2627459726d8SSatish Balay     ierr = PetscNew(&bNew);CHKERRQ(ierr);
2628dff059c6SToby Isaac     bNew->numcomps = b->numcomps;
2629dff059c6SToby Isaac     ierr = PetscMalloc1(bNew->numcomps, &bNew->comps);CHKERRQ(ierr);
2630dff059c6SToby Isaac     ierr = PetscMemcpy(bNew->comps, b->comps, bNew->numcomps*sizeof(PetscInt));CHKERRQ(ierr);
2631dff059c6SToby Isaac     bNew->numids = b->numids;
2632dff059c6SToby Isaac     ierr = PetscMalloc1(bNew->numids, &bNew->ids);CHKERRQ(ierr);
2633dff059c6SToby Isaac     ierr = PetscMemcpy(bNew->ids, b->ids, bNew->numids*sizeof(PetscInt));CHKERRQ(ierr);
2634dff059c6SToby Isaac     ierr = PetscStrallocpy(b->labelname,(char **) &(bNew->labelname));CHKERRQ(ierr);
2635dff059c6SToby Isaac     ierr = PetscStrallocpy(b->name,(char **) &(bNew->name));CHKERRQ(ierr);
2636dff059c6SToby Isaac     bNew->ctx   = b->ctx;
2637f971fd6bSMatthew G. Knepley     bNew->type  = b->type;
2638dff059c6SToby Isaac     bNew->field = b->field;
2639dff059c6SToby Isaac     bNew->func  = b->func;
2640dff059c6SToby Isaac 
2641dff059c6SToby Isaac     *lastnext = bNew;
2642dff059c6SToby Isaac     lastnext = &(bNew->next);
2643dff059c6SToby Isaac   }
2644dff059c6SToby Isaac   PetscFunctionReturn(0);
2645dff059c6SToby Isaac }
2646dff059c6SToby Isaac 
2647da51fcedSMatthew G. Knepley /*@
2648da51fcedSMatthew G. Knepley   PetscDSCopyEquations - Copy all pointwise function pointers to the new problem
2649da51fcedSMatthew G. Knepley 
2650da51fcedSMatthew G. Knepley   Not collective
2651da51fcedSMatthew G. Knepley 
2652da51fcedSMatthew G. Knepley   Input Parameter:
2653da51fcedSMatthew G. Knepley . prob - The PetscDS object
2654da51fcedSMatthew G. Knepley 
2655da51fcedSMatthew G. Knepley   Output Parameter:
2656da51fcedSMatthew G. Knepley . newprob - The PetscDS copy
2657da51fcedSMatthew G. Knepley 
2658da51fcedSMatthew G. Knepley   Level: intermediate
2659da51fcedSMatthew G. Knepley 
2660da51fcedSMatthew G. Knepley .seealso: PetscDSSetResidual(), PetscDSSetJacobian(), PetscDSSetRiemannSolver(), PetscDSSetBdResidual(), PetscDSSetBdJacobian(), PetscDSCreate()
2661da51fcedSMatthew G. Knepley @*/
2662da51fcedSMatthew G. Knepley PetscErrorCode PetscDSCopyEquations(PetscDS prob, PetscDS newprob)
2663da51fcedSMatthew G. Knepley {
2664da51fcedSMatthew G. Knepley   PetscInt       Nf, Ng, f, g;
2665da51fcedSMatthew G. Knepley   PetscErrorCode ierr;
2666da51fcedSMatthew G. Knepley 
2667da51fcedSMatthew G. Knepley   PetscFunctionBegin;
2668da51fcedSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
2669da51fcedSMatthew G. Knepley   PetscValidHeaderSpecific(newprob, PETSCDS_CLASSID, 2);
2670da51fcedSMatthew G. Knepley   ierr = PetscDSGetNumFields(prob, &Nf);CHKERRQ(ierr);
2671da51fcedSMatthew G. Knepley   ierr = PetscDSGetNumFields(newprob, &Ng);CHKERRQ(ierr);
2672da51fcedSMatthew G. Knepley   if (Nf != Ng) SETERRQ2(PetscObjectComm((PetscObject) prob), PETSC_ERR_ARG_SIZ, "Number of fields must match %D != %D", Nf, Ng);CHKERRQ(ierr);
2673da51fcedSMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
2674da51fcedSMatthew G. Knepley     PetscPointFunc   obj;
2675da51fcedSMatthew G. Knepley     PetscPointFunc   f0, f1;
2676da51fcedSMatthew G. Knepley     PetscPointJac    g0, g1, g2, g3;
2677da51fcedSMatthew G. Knepley     PetscBdPointFunc f0Bd, f1Bd;
2678da51fcedSMatthew G. Knepley     PetscBdPointJac  g0Bd, g1Bd, g2Bd, g3Bd;
2679da51fcedSMatthew G. Knepley     PetscRiemannFunc r;
2680da51fcedSMatthew G. Knepley 
2681da51fcedSMatthew G. Knepley     ierr = PetscDSGetObjective(prob, f, &obj);CHKERRQ(ierr);
2682da51fcedSMatthew G. Knepley     ierr = PetscDSGetResidual(prob, f, &f0, &f1);CHKERRQ(ierr);
2683da51fcedSMatthew G. Knepley     ierr = PetscDSGetBdResidual(prob, f, &f0Bd, &f1Bd);CHKERRQ(ierr);
2684da51fcedSMatthew G. Knepley     ierr = PetscDSGetRiemannSolver(prob, f, &r);CHKERRQ(ierr);
2685da51fcedSMatthew G. Knepley     ierr = PetscDSSetObjective(newprob, f, obj);CHKERRQ(ierr);
2686da51fcedSMatthew G. Knepley     ierr = PetscDSSetResidual(newprob, f, f0, f1);CHKERRQ(ierr);
2687da51fcedSMatthew G. Knepley     ierr = PetscDSSetBdResidual(newprob, f, f0Bd, f1Bd);CHKERRQ(ierr);
2688da51fcedSMatthew G. Knepley     ierr = PetscDSSetRiemannSolver(newprob, f, r);CHKERRQ(ierr);
2689da51fcedSMatthew G. Knepley     for (g = 0; g < Nf; ++g) {
2690da51fcedSMatthew G. Knepley       ierr = PetscDSGetJacobian(prob, f, g, &g0, &g1, &g2, &g3);CHKERRQ(ierr);
2691da51fcedSMatthew G. Knepley       ierr = PetscDSGetBdJacobian(prob, f, g, &g0Bd, &g1Bd, &g2Bd, &g3Bd);CHKERRQ(ierr);
2692da51fcedSMatthew G. Knepley       ierr = PetscDSSetJacobian(newprob, f, g, g0, g1, g2, g3);CHKERRQ(ierr);
2693da51fcedSMatthew G. Knepley       ierr = PetscDSSetBdJacobian(newprob, f, g, g0Bd, g1Bd, g2Bd, g3Bd);CHKERRQ(ierr);
2694da51fcedSMatthew G. Knepley     }
2695da51fcedSMatthew G. Knepley   }
2696da51fcedSMatthew G. Knepley   PetscFunctionReturn(0);
2697da51fcedSMatthew G. Knepley }
2698da51fcedSMatthew G. Knepley 
2699bc4ae4beSMatthew G. Knepley static PetscErrorCode PetscDSDestroy_Basic(PetscDS prob)
27002764a2aaSMatthew G. Knepley {
2701931fb3b8SToby Isaac   PetscErrorCode      ierr;
2702931fb3b8SToby Isaac 
27032764a2aaSMatthew G. Knepley   PetscFunctionBegin;
2704931fb3b8SToby Isaac   ierr = PetscFree(prob->data);CHKERRQ(ierr);
27052764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
27062764a2aaSMatthew G. Knepley }
27072764a2aaSMatthew G. Knepley 
2708bc4ae4beSMatthew G. Knepley static PetscErrorCode PetscDSInitialize_Basic(PetscDS prob)
27092764a2aaSMatthew G. Knepley {
27102764a2aaSMatthew G. Knepley   PetscFunctionBegin;
27112764a2aaSMatthew G. Knepley   prob->ops->setfromoptions = NULL;
27122764a2aaSMatthew G. Knepley   prob->ops->setup          = NULL;
27132764a2aaSMatthew G. Knepley   prob->ops->view           = NULL;
27142764a2aaSMatthew G. Knepley   prob->ops->destroy        = PetscDSDestroy_Basic;
27152764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
27162764a2aaSMatthew G. Knepley }
27172764a2aaSMatthew G. Knepley 
27182764a2aaSMatthew G. Knepley /*MC
27192764a2aaSMatthew G. Knepley   PETSCDSBASIC = "basic" - A discrete system with pointwise residual and boundary residual functions
27202764a2aaSMatthew G. Knepley 
27212764a2aaSMatthew G. Knepley   Level: intermediate
27222764a2aaSMatthew G. Knepley 
27232764a2aaSMatthew G. Knepley .seealso: PetscDSType, PetscDSCreate(), PetscDSSetType()
27242764a2aaSMatthew G. Knepley M*/
27252764a2aaSMatthew G. Knepley 
27262764a2aaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDSCreate_Basic(PetscDS prob)
27272764a2aaSMatthew G. Knepley {
27282764a2aaSMatthew G. Knepley   PetscDS_Basic *b;
27292764a2aaSMatthew G. Knepley   PetscErrorCode      ierr;
27302764a2aaSMatthew G. Knepley 
27312764a2aaSMatthew G. Knepley   PetscFunctionBegin;
2732931fb3b8SToby Isaac   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
27332764a2aaSMatthew G. Knepley   ierr       = PetscNewLog(prob, &b);CHKERRQ(ierr);
27342764a2aaSMatthew G. Knepley   prob->data = b;
27352764a2aaSMatthew G. Knepley 
27362764a2aaSMatthew G. Knepley   ierr = PetscDSInitialize_Basic(prob);CHKERRQ(ierr);
27372764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
27382764a2aaSMatthew G. Knepley }
2739