xref: /petsc/src/dm/dt/interface/dtds.c (revision 475e0ac9a142c3397a2a9b490533a78f024949c6)
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 #undef __FUNCT__
92764a2aaSMatthew G. Knepley #define __FUNCT__ "PetscDSRegister"
102764a2aaSMatthew G. Knepley /*@C
112764a2aaSMatthew G. Knepley   PetscDSRegister - Adds a new PetscDS implementation
122764a2aaSMatthew G. Knepley 
132764a2aaSMatthew G. Knepley   Not Collective
142764a2aaSMatthew G. Knepley 
152764a2aaSMatthew G. Knepley   Input Parameters:
162764a2aaSMatthew G. Knepley + name        - The name of a new user-defined creation routine
172764a2aaSMatthew G. Knepley - create_func - The creation routine itself
182764a2aaSMatthew G. Knepley 
192764a2aaSMatthew G. Knepley   Notes:
202764a2aaSMatthew G. Knepley   PetscDSRegister() may be called multiple times to add several user-defined PetscDSs
212764a2aaSMatthew G. Knepley 
222764a2aaSMatthew G. Knepley   Sample usage:
232764a2aaSMatthew G. Knepley .vb
242764a2aaSMatthew G. Knepley     PetscDSRegister("my_ds", MyPetscDSCreate);
252764a2aaSMatthew G. Knepley .ve
262764a2aaSMatthew G. Knepley 
272764a2aaSMatthew G. Knepley   Then, your PetscDS type can be chosen with the procedural interface via
282764a2aaSMatthew G. Knepley .vb
292764a2aaSMatthew G. Knepley     PetscDSCreate(MPI_Comm, PetscDS *);
302764a2aaSMatthew G. Knepley     PetscDSSetType(PetscDS, "my_ds");
312764a2aaSMatthew G. Knepley .ve
322764a2aaSMatthew G. Knepley    or at runtime via the option
332764a2aaSMatthew G. Knepley .vb
342764a2aaSMatthew G. Knepley     -petscds_type my_ds
352764a2aaSMatthew G. Knepley .ve
362764a2aaSMatthew G. Knepley 
372764a2aaSMatthew G. Knepley   Level: advanced
382764a2aaSMatthew G. Knepley 
392764a2aaSMatthew G. Knepley .keywords: PetscDS, register
402764a2aaSMatthew G. Knepley .seealso: PetscDSRegisterAll(), PetscDSRegisterDestroy()
412764a2aaSMatthew G. Knepley 
422764a2aaSMatthew G. Knepley @*/
432764a2aaSMatthew G. Knepley PetscErrorCode PetscDSRegister(const char sname[], PetscErrorCode (*function)(PetscDS))
442764a2aaSMatthew G. Knepley {
452764a2aaSMatthew G. Knepley   PetscErrorCode ierr;
462764a2aaSMatthew G. Knepley 
472764a2aaSMatthew G. Knepley   PetscFunctionBegin;
482764a2aaSMatthew G. Knepley   ierr = PetscFunctionListAdd(&PetscDSList, sname, function);CHKERRQ(ierr);
492764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
502764a2aaSMatthew G. Knepley }
512764a2aaSMatthew G. Knepley 
522764a2aaSMatthew G. Knepley #undef __FUNCT__
532764a2aaSMatthew G. Knepley #define __FUNCT__ "PetscDSSetType"
542764a2aaSMatthew G. Knepley /*@C
552764a2aaSMatthew G. Knepley   PetscDSSetType - Builds a particular PetscDS
562764a2aaSMatthew G. Knepley 
572764a2aaSMatthew G. Knepley   Collective on PetscDS
582764a2aaSMatthew G. Knepley 
592764a2aaSMatthew G. Knepley   Input Parameters:
602764a2aaSMatthew G. Knepley + prob - The PetscDS object
612764a2aaSMatthew G. Knepley - name - The kind of system
622764a2aaSMatthew G. Knepley 
632764a2aaSMatthew G. Knepley   Options Database Key:
642764a2aaSMatthew G. Knepley . -petscds_type <type> - Sets the PetscDS type; use -help for a list of available types
652764a2aaSMatthew G. Knepley 
662764a2aaSMatthew G. Knepley   Level: intermediate
672764a2aaSMatthew G. Knepley 
682764a2aaSMatthew G. Knepley .keywords: PetscDS, set, type
692764a2aaSMatthew G. Knepley .seealso: PetscDSGetType(), PetscDSCreate()
702764a2aaSMatthew G. Knepley @*/
712764a2aaSMatthew G. Knepley PetscErrorCode PetscDSSetType(PetscDS prob, PetscDSType name)
722764a2aaSMatthew G. Knepley {
732764a2aaSMatthew G. Knepley   PetscErrorCode (*r)(PetscDS);
742764a2aaSMatthew G. Knepley   PetscBool      match;
752764a2aaSMatthew G. Knepley   PetscErrorCode ierr;
762764a2aaSMatthew G. Knepley 
772764a2aaSMatthew G. Knepley   PetscFunctionBegin;
782764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
792764a2aaSMatthew G. Knepley   ierr = PetscObjectTypeCompare((PetscObject) prob, name, &match);CHKERRQ(ierr);
802764a2aaSMatthew G. Knepley   if (match) PetscFunctionReturn(0);
812764a2aaSMatthew G. Knepley 
820f51fdf8SToby Isaac   ierr = PetscDSRegisterAll();CHKERRQ(ierr);
832764a2aaSMatthew G. Knepley   ierr = PetscFunctionListFind(PetscDSList, name, &r);CHKERRQ(ierr);
842764a2aaSMatthew G. Knepley   if (!r) SETERRQ1(PetscObjectComm((PetscObject) prob), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown PetscDS type: %s", name);
852764a2aaSMatthew G. Knepley 
862764a2aaSMatthew G. Knepley   if (prob->ops->destroy) {
872764a2aaSMatthew G. Knepley     ierr             = (*prob->ops->destroy)(prob);CHKERRQ(ierr);
882764a2aaSMatthew G. Knepley     prob->ops->destroy = NULL;
892764a2aaSMatthew G. Knepley   }
902764a2aaSMatthew G. Knepley   ierr = (*r)(prob);CHKERRQ(ierr);
912764a2aaSMatthew G. Knepley   ierr = PetscObjectChangeTypeName((PetscObject) prob, name);CHKERRQ(ierr);
922764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
932764a2aaSMatthew G. Knepley }
942764a2aaSMatthew G. Knepley 
952764a2aaSMatthew G. Knepley #undef __FUNCT__
962764a2aaSMatthew G. Knepley #define __FUNCT__ "PetscDSGetType"
972764a2aaSMatthew G. Knepley /*@C
982764a2aaSMatthew G. Knepley   PetscDSGetType - Gets the PetscDS type name (as a string) from the object.
992764a2aaSMatthew G. Knepley 
1002764a2aaSMatthew G. Knepley   Not Collective
1012764a2aaSMatthew G. Knepley 
1022764a2aaSMatthew G. Knepley   Input Parameter:
1032764a2aaSMatthew G. Knepley . prob  - The PetscDS
1042764a2aaSMatthew G. Knepley 
1052764a2aaSMatthew G. Knepley   Output Parameter:
1062764a2aaSMatthew G. Knepley . name - The PetscDS type name
1072764a2aaSMatthew G. Knepley 
1082764a2aaSMatthew G. Knepley   Level: intermediate
1092764a2aaSMatthew G. Knepley 
1102764a2aaSMatthew G. Knepley .keywords: PetscDS, get, type, name
1112764a2aaSMatthew G. Knepley .seealso: PetscDSSetType(), PetscDSCreate()
1122764a2aaSMatthew G. Knepley @*/
1132764a2aaSMatthew G. Knepley PetscErrorCode PetscDSGetType(PetscDS prob, PetscDSType *name)
1142764a2aaSMatthew G. Knepley {
1152764a2aaSMatthew G. Knepley   PetscErrorCode ierr;
1162764a2aaSMatthew G. Knepley 
1172764a2aaSMatthew G. Knepley   PetscFunctionBegin;
1182764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
119c959eef4SJed Brown   PetscValidPointer(name, 2);
1200f51fdf8SToby Isaac   ierr = PetscDSRegisterAll();CHKERRQ(ierr);
1212764a2aaSMatthew G. Knepley   *name = ((PetscObject) prob)->type_name;
1222764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
1232764a2aaSMatthew G. Knepley }
1242764a2aaSMatthew G. Knepley 
1252764a2aaSMatthew G. Knepley #undef __FUNCT__
1267d8a60eaSMatthew G. Knepley #define __FUNCT__ "PetscDSView_Ascii"
1277d8a60eaSMatthew G. Knepley static PetscErrorCode PetscDSView_Ascii(PetscDS prob, PetscViewer viewer)
1287d8a60eaSMatthew G. Knepley {
1297d8a60eaSMatthew G. Knepley   PetscViewerFormat format;
1307d8a60eaSMatthew G. Knepley   PetscInt          f;
1317d8a60eaSMatthew G. Knepley   PetscErrorCode    ierr;
1327d8a60eaSMatthew G. Knepley 
1337d8a60eaSMatthew G. Knepley   PetscFunctionBegin;
1347d8a60eaSMatthew G. Knepley   ierr = PetscViewerGetFormat(viewer, &format);CHKERRQ(ierr);
1357d8a60eaSMatthew G. Knepley   ierr = PetscViewerASCIIPrintf(viewer, "Discrete System with %d fields\n", prob->Nf);CHKERRQ(ierr);
1367d8a60eaSMatthew G. Knepley   ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
1377d8a60eaSMatthew G. Knepley   for (f = 0; f < prob->Nf; ++f) {
1387d8a60eaSMatthew G. Knepley     PetscObject  obj;
1397d8a60eaSMatthew G. Knepley     PetscClassId id;
1407d8a60eaSMatthew G. Knepley     const char  *name;
1417d8a60eaSMatthew G. Knepley     PetscInt     Nc;
1427d8a60eaSMatthew G. Knepley 
1437d8a60eaSMatthew G. Knepley     ierr = PetscDSGetDiscretization(prob, f, &obj);CHKERRQ(ierr);
1447d8a60eaSMatthew G. Knepley     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
1457d8a60eaSMatthew G. Knepley     ierr = PetscObjectGetName(obj, &name);CHKERRQ(ierr);
1467d8a60eaSMatthew G. Knepley     ierr = PetscViewerASCIIPrintf(viewer, "Field %s", name ? name : "<unknown>");CHKERRQ(ierr);
1477d8a60eaSMatthew G. Knepley     if (id == PETSCFE_CLASSID)      {
1487d8a60eaSMatthew G. Knepley       ierr = PetscFEGetNumComponents((PetscFE) obj, &Nc);CHKERRQ(ierr);
1497d8a60eaSMatthew G. Knepley       ierr = PetscViewerASCIIPrintf(viewer, " FEM");CHKERRQ(ierr);
1507d8a60eaSMatthew G. Knepley     } else if (id == PETSCFV_CLASSID) {
1517d8a60eaSMatthew G. Knepley       ierr = PetscFVGetNumComponents((PetscFV) obj, &Nc);CHKERRQ(ierr);
1527d8a60eaSMatthew G. Knepley       ierr = PetscViewerASCIIPrintf(viewer, " FVM");CHKERRQ(ierr);
1537d8a60eaSMatthew G. Knepley     }
1547d8a60eaSMatthew G. Knepley     else SETERRQ1(PetscObjectComm((PetscObject) prob), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", f);
1557d8a60eaSMatthew G. Knepley     if (Nc > 1) {ierr = PetscViewerASCIIPrintf(viewer, "%d components", Nc);CHKERRQ(ierr);}
1567d8a60eaSMatthew G. Knepley     else        {ierr = PetscViewerASCIIPrintf(viewer, "%d component ", Nc);CHKERRQ(ierr);}
157249df284SMatthew G. Knepley     if (prob->implicit[f]) {ierr = PetscViewerASCIIPrintf(viewer, " (implicit)");CHKERRQ(ierr);}
158249df284SMatthew G. Knepley     else                   {ierr = PetscViewerASCIIPrintf(viewer, " (explicit)");CHKERRQ(ierr);}
159a6cbbb48SMatthew G. Knepley     if (prob->adjacency[f*2+0]) {
160a6cbbb48SMatthew G. Knepley       if (prob->adjacency[f*2+1]) {ierr = PetscViewerASCIIPrintf(viewer, " (adj FVM++)");CHKERRQ(ierr);}
161a6cbbb48SMatthew G. Knepley       else                        {ierr = PetscViewerASCIIPrintf(viewer, " (adj FVM)");CHKERRQ(ierr);}
162a6cbbb48SMatthew G. Knepley     } else {
163a6cbbb48SMatthew G. Knepley       if (prob->adjacency[f*2+1]) {ierr = PetscViewerASCIIPrintf(viewer, " (adj FEM)");CHKERRQ(ierr);}
164a6cbbb48SMatthew G. Knepley       else                        {ierr = PetscViewerASCIIPrintf(viewer, " (adj FUNKY)");CHKERRQ(ierr);}
165a6cbbb48SMatthew G. Knepley     }
1667d8a60eaSMatthew G. Knepley     ierr = PetscViewerASCIIPrintf(viewer, "\n");CHKERRQ(ierr);
1677d8a60eaSMatthew G. Knepley     if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
1687d8a60eaSMatthew G. Knepley       if (id == PETSCFE_CLASSID)      {ierr = PetscFEView((PetscFE) obj, viewer);CHKERRQ(ierr);}
1697d8a60eaSMatthew G. Knepley       else if (id == PETSCFV_CLASSID) {ierr = PetscFVView((PetscFV) obj, viewer);CHKERRQ(ierr);}
1707d8a60eaSMatthew G. Knepley     }
1717d8a60eaSMatthew G. Knepley   }
1727d8a60eaSMatthew G. Knepley   ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1737d8a60eaSMatthew G. Knepley   PetscFunctionReturn(0);
1747d8a60eaSMatthew G. Knepley }
1757d8a60eaSMatthew G. Knepley 
1767d8a60eaSMatthew G. Knepley #undef __FUNCT__
1772764a2aaSMatthew G. Knepley #define __FUNCT__ "PetscDSView"
1782764a2aaSMatthew G. Knepley /*@C
1792764a2aaSMatthew G. Knepley   PetscDSView - Views a PetscDS
1802764a2aaSMatthew G. Knepley 
1812764a2aaSMatthew G. Knepley   Collective on PetscDS
1822764a2aaSMatthew G. Knepley 
1832764a2aaSMatthew G. Knepley   Input Parameter:
1842764a2aaSMatthew G. Knepley + prob - the PetscDS object to view
1852764a2aaSMatthew G. Knepley - v  - the viewer
1862764a2aaSMatthew G. Knepley 
1872764a2aaSMatthew G. Knepley   Level: developer
1882764a2aaSMatthew G. Knepley 
1892764a2aaSMatthew G. Knepley .seealso PetscDSDestroy()
1902764a2aaSMatthew G. Knepley @*/
1912764a2aaSMatthew G. Knepley PetscErrorCode PetscDSView(PetscDS prob, PetscViewer v)
1922764a2aaSMatthew G. Knepley {
1937d8a60eaSMatthew G. Knepley   PetscBool      iascii;
1942764a2aaSMatthew G. Knepley   PetscErrorCode ierr;
1952764a2aaSMatthew G. Knepley 
1962764a2aaSMatthew G. Knepley   PetscFunctionBegin;
1972764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1982764a2aaSMatthew G. Knepley   if (!v) {ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject) prob), &v);CHKERRQ(ierr);}
1997d8a60eaSMatthew G. Knepley   else    {PetscValidHeaderSpecific(v, PETSC_VIEWER_CLASSID, 2);}
2007d8a60eaSMatthew G. Knepley   ierr = PetscObjectTypeCompare((PetscObject) v, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
2017d8a60eaSMatthew G. Knepley   if (iascii) {ierr = PetscDSView_Ascii(prob, v);CHKERRQ(ierr);}
2022764a2aaSMatthew G. Knepley   if (prob->ops->view) {ierr = (*prob->ops->view)(prob, v);CHKERRQ(ierr);}
2032764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
2042764a2aaSMatthew G. Knepley }
2052764a2aaSMatthew G. Knepley 
2062764a2aaSMatthew G. Knepley #undef __FUNCT__
2072764a2aaSMatthew G. Knepley #define __FUNCT__ "PetscDSSetFromOptions"
2082764a2aaSMatthew G. Knepley /*@
2092764a2aaSMatthew G. Knepley   PetscDSSetFromOptions - sets parameters in a PetscDS from the options database
2102764a2aaSMatthew G. Knepley 
2112764a2aaSMatthew G. Knepley   Collective on PetscDS
2122764a2aaSMatthew G. Knepley 
2132764a2aaSMatthew G. Knepley   Input Parameter:
2142764a2aaSMatthew G. Knepley . prob - the PetscDS object to set options for
2152764a2aaSMatthew G. Knepley 
2162764a2aaSMatthew G. Knepley   Options Database:
2172764a2aaSMatthew G. Knepley 
2182764a2aaSMatthew G. Knepley   Level: developer
2192764a2aaSMatthew G. Knepley 
2202764a2aaSMatthew G. Knepley .seealso PetscDSView()
2212764a2aaSMatthew G. Knepley @*/
2222764a2aaSMatthew G. Knepley PetscErrorCode PetscDSSetFromOptions(PetscDS prob)
2232764a2aaSMatthew G. Knepley {
2242764a2aaSMatthew G. Knepley   const char    *defaultType;
2252764a2aaSMatthew G. Knepley   char           name[256];
2262764a2aaSMatthew G. Knepley   PetscBool      flg;
2272764a2aaSMatthew G. Knepley   PetscErrorCode ierr;
2282764a2aaSMatthew G. Knepley 
2292764a2aaSMatthew G. Knepley   PetscFunctionBegin;
2302764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
2312764a2aaSMatthew G. Knepley   if (!((PetscObject) prob)->type_name) {
2322764a2aaSMatthew G. Knepley     defaultType = PETSCDSBASIC;
2332764a2aaSMatthew G. Knepley   } else {
2342764a2aaSMatthew G. Knepley     defaultType = ((PetscObject) prob)->type_name;
2352764a2aaSMatthew G. Knepley   }
2360f51fdf8SToby Isaac   ierr = PetscDSRegisterAll();CHKERRQ(ierr);
2372764a2aaSMatthew G. Knepley 
2382764a2aaSMatthew G. Knepley   ierr = PetscObjectOptionsBegin((PetscObject) prob);CHKERRQ(ierr);
2392764a2aaSMatthew G. Knepley   ierr = PetscOptionsFList("-petscds_type", "Discrete System", "PetscDSSetType", PetscDSList, defaultType, name, 256, &flg);CHKERRQ(ierr);
2402764a2aaSMatthew G. Knepley   if (flg) {
2412764a2aaSMatthew G. Knepley     ierr = PetscDSSetType(prob, name);CHKERRQ(ierr);
2422764a2aaSMatthew G. Knepley   } else if (!((PetscObject) prob)->type_name) {
2432764a2aaSMatthew G. Knepley     ierr = PetscDSSetType(prob, defaultType);CHKERRQ(ierr);
2442764a2aaSMatthew G. Knepley   }
2452764a2aaSMatthew G. Knepley   if (prob->ops->setfromoptions) {ierr = (*prob->ops->setfromoptions)(prob);CHKERRQ(ierr);}
2462764a2aaSMatthew G. Knepley   /* process any options handlers added with PetscObjectAddOptionsHandler() */
2470633abcbSJed Brown   ierr = PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject) prob);CHKERRQ(ierr);
2482764a2aaSMatthew G. Knepley   ierr = PetscOptionsEnd();CHKERRQ(ierr);
2492764a2aaSMatthew G. Knepley   ierr = PetscDSViewFromOptions(prob, NULL, "-petscds_view");CHKERRQ(ierr);
2502764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
2512764a2aaSMatthew G. Knepley }
2522764a2aaSMatthew G. Knepley 
2532764a2aaSMatthew G. Knepley #undef __FUNCT__
2542764a2aaSMatthew G. Knepley #define __FUNCT__ "PetscDSSetUp"
2552764a2aaSMatthew G. Knepley /*@C
2562764a2aaSMatthew G. Knepley   PetscDSSetUp - Construct data structures for the PetscDS
2572764a2aaSMatthew G. Knepley 
2582764a2aaSMatthew G. Knepley   Collective on PetscDS
2592764a2aaSMatthew G. Knepley 
2602764a2aaSMatthew G. Knepley   Input Parameter:
2612764a2aaSMatthew G. Knepley . prob - the PetscDS object to setup
2622764a2aaSMatthew G. Knepley 
2632764a2aaSMatthew G. Knepley   Level: developer
2642764a2aaSMatthew G. Knepley 
2652764a2aaSMatthew G. Knepley .seealso PetscDSView(), PetscDSDestroy()
2662764a2aaSMatthew G. Knepley @*/
2672764a2aaSMatthew G. Knepley PetscErrorCode PetscDSSetUp(PetscDS prob)
2682764a2aaSMatthew G. Knepley {
2692764a2aaSMatthew G. Knepley   const PetscInt Nf = prob->Nf;
2702764a2aaSMatthew G. Knepley   PetscInt       dim, work, NcMax = 0, NqMax = 0, f;
2712764a2aaSMatthew G. Knepley   PetscErrorCode ierr;
2722764a2aaSMatthew G. Knepley 
2732764a2aaSMatthew G. Knepley   PetscFunctionBegin;
2742764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
2752764a2aaSMatthew G. Knepley   if (prob->setup) PetscFunctionReturn(0);
2762764a2aaSMatthew G. Knepley   /* Calculate sizes */
2772764a2aaSMatthew G. Knepley   ierr = PetscDSGetSpatialDimension(prob, &dim);CHKERRQ(ierr);
2782764a2aaSMatthew G. Knepley   prob->totDim = prob->totDimBd = prob->totComp = 0;
279194d53e6SMatthew G. Knepley   ierr = PetscCalloc4(Nf+1,&prob->off,Nf+1,&prob->offDer,Nf+1,&prob->offBd,Nf+1,&prob->offDerBd);CHKERRQ(ierr);
2802764a2aaSMatthew G. Knepley   ierr = PetscMalloc4(Nf,&prob->basis,Nf,&prob->basisDer,Nf,&prob->basisBd,Nf,&prob->basisDerBd);CHKERRQ(ierr);
2812764a2aaSMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
2822764a2aaSMatthew G. Knepley     PetscFE         feBd = (PetscFE) prob->discBd[f];
2839de99aefSMatthew G. Knepley     PetscObject     obj;
2849de99aefSMatthew G. Knepley     PetscClassId    id;
2852764a2aaSMatthew G. Knepley     PetscQuadrature q;
2869de99aefSMatthew G. Knepley     PetscInt        Nq = 0, Nb, Nc;
2872764a2aaSMatthew G. Knepley 
2889de99aefSMatthew G. Knepley     ierr = PetscDSGetDiscretization(prob, f, &obj);CHKERRQ(ierr);
2899de99aefSMatthew G. Knepley     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
2909de99aefSMatthew G. Knepley     if (id == PETSCFE_CLASSID)      {
2919de99aefSMatthew G. Knepley       PetscFE fe = (PetscFE) obj;
2929de99aefSMatthew G. Knepley 
2932764a2aaSMatthew G. Knepley       ierr = PetscFEGetQuadrature(fe, &q);CHKERRQ(ierr);
2942764a2aaSMatthew G. Knepley       ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr);
2952764a2aaSMatthew G. Knepley       ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
2962764a2aaSMatthew G. Knepley       ierr = PetscFEGetDefaultTabulation(fe, &prob->basis[f], &prob->basisDer[f], NULL);CHKERRQ(ierr);
2979de99aefSMatthew G. Knepley     } else if (id == PETSCFV_CLASSID) {
2989de99aefSMatthew G. Knepley       PetscFV fv = (PetscFV) obj;
2999de99aefSMatthew G. Knepley 
3009de99aefSMatthew G. Knepley       ierr = PetscFVGetQuadrature(fv, &q);CHKERRQ(ierr);
3019de99aefSMatthew G. Knepley       Nb   = 1;
3029de99aefSMatthew G. Knepley       ierr = PetscFVGetNumComponents(fv, &Nc);CHKERRQ(ierr);
3036c1a3d01SMatthew G. Knepley       ierr = PetscFVGetDefaultTabulation(fv, &prob->basis[f], &prob->basisDer[f], NULL);CHKERRQ(ierr);
304abac5ca0SMatthew G. Knepley     } else SETERRQ1(PetscObjectComm((PetscObject) prob), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", f);
305194d53e6SMatthew G. Knepley     prob->off[f+1]    = Nc     + prob->off[f];
306194d53e6SMatthew G. Knepley     prob->offDer[f+1] = Nc*dim + prob->offDer[f];
3079de99aefSMatthew G. Knepley     if (q) {ierr = PetscQuadratureGetData(q, NULL, &Nq, NULL, NULL);CHKERRQ(ierr);}
3082764a2aaSMatthew G. Knepley     NqMax          = PetscMax(NqMax, Nq);
3092764a2aaSMatthew G. Knepley     NcMax          = PetscMax(NcMax, Nc);
3102764a2aaSMatthew G. Knepley     prob->totDim  += Nb*Nc;
3112764a2aaSMatthew G. Knepley     prob->totComp += Nc;
3122764a2aaSMatthew G. Knepley     if (feBd) {
3132764a2aaSMatthew G. Knepley       ierr = PetscFEGetDimension(feBd, &Nb);CHKERRQ(ierr);
3142764a2aaSMatthew G. Knepley       ierr = PetscFEGetNumComponents(feBd, &Nc);CHKERRQ(ierr);
3152764a2aaSMatthew G. Knepley       ierr = PetscFEGetDefaultTabulation(feBd, &prob->basisBd[f], &prob->basisDerBd[f], NULL);CHKERRQ(ierr);
3162764a2aaSMatthew G. Knepley       prob->totDimBd += Nb*Nc;
317194d53e6SMatthew G. Knepley       prob->offBd[f+1]    = Nc     + prob->offBd[f];
318194d53e6SMatthew G. Knepley       prob->offDerBd[f+1] = Nc*dim + prob->offDerBd[f];
3192764a2aaSMatthew G. Knepley     }
3202764a2aaSMatthew G. Knepley   }
3212764a2aaSMatthew G. Knepley   work = PetscMax(prob->totComp*dim, PetscSqr(NcMax*dim));
3222764a2aaSMatthew G. Knepley   /* Allocate works space */
3232764a2aaSMatthew 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);
3242764a2aaSMatthew 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);
3252764a2aaSMatthew G. Knepley   if (prob->ops->setup) {ierr = (*prob->ops->setup)(prob);CHKERRQ(ierr);}
3262764a2aaSMatthew G. Knepley   prob->setup = PETSC_TRUE;
3272764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
3282764a2aaSMatthew G. Knepley }
3292764a2aaSMatthew G. Knepley 
3302764a2aaSMatthew G. Knepley #undef __FUNCT__
3312764a2aaSMatthew G. Knepley #define __FUNCT__ "PetscDSDestroyStructs_Static"
3322764a2aaSMatthew G. Knepley static PetscErrorCode PetscDSDestroyStructs_Static(PetscDS prob)
3332764a2aaSMatthew G. Knepley {
3342764a2aaSMatthew G. Knepley   PetscErrorCode ierr;
3352764a2aaSMatthew G. Knepley 
3362764a2aaSMatthew G. Knepley   PetscFunctionBegin;
337194d53e6SMatthew G. Knepley   ierr = PetscFree4(prob->off,prob->offDer,prob->offBd,prob->offDerBd);CHKERRQ(ierr);
3382764a2aaSMatthew G. Knepley   ierr = PetscFree4(prob->basis,prob->basisDer,prob->basisBd,prob->basisDerBd);CHKERRQ(ierr);
3392764a2aaSMatthew G. Knepley   ierr = PetscFree5(prob->u,prob->u_t,prob->u_x,prob->x,prob->refSpaceDer);CHKERRQ(ierr);
3402764a2aaSMatthew G. Knepley   ierr = PetscFree6(prob->f0,prob->f1,prob->g0,prob->g1,prob->g2,prob->g3);CHKERRQ(ierr);
3412764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
3422764a2aaSMatthew G. Knepley }
3432764a2aaSMatthew G. Knepley 
3442764a2aaSMatthew G. Knepley #undef __FUNCT__
3452764a2aaSMatthew G. Knepley #define __FUNCT__ "PetscDSEnlarge_Static"
3462764a2aaSMatthew G. Knepley static PetscErrorCode PetscDSEnlarge_Static(PetscDS prob, PetscInt NfNew)
3472764a2aaSMatthew G. Knepley {
3482764a2aaSMatthew G. Knepley   PetscObject      *tmpd, *tmpdbd;
349a6cbbb48SMatthew G. Knepley   PetscBool        *tmpi, *tmpa;
3502aa1fc23SMatthew G. Knepley   PetscPointFunc   *tmpobj, *tmpf;
351*475e0ac9SMatthew G. Knepley   PetscPointJac    *tmpg, *tmpgp;
3522aa1fc23SMatthew G. Knepley   PetscBdPointFunc *tmpfbd;
3532aa1fc23SMatthew G. Knepley   PetscBdPointJac  *tmpgbd;
354194d53e6SMatthew G. Knepley   PetscRiemannFunc *tmpr;
3550c2f2876SMatthew G. Knepley   void            **tmpctx;
356a6cbbb48SMatthew G. Knepley   PetscInt          Nf = prob->Nf, f, i;
3572764a2aaSMatthew G. Knepley   PetscErrorCode    ierr;
3582764a2aaSMatthew G. Knepley 
3592764a2aaSMatthew G. Knepley   PetscFunctionBegin;
3602764a2aaSMatthew G. Knepley   if (Nf >= NfNew) PetscFunctionReturn(0);
3612764a2aaSMatthew G. Knepley   prob->setup = PETSC_FALSE;
3622764a2aaSMatthew G. Knepley   ierr = PetscDSDestroyStructs_Static(prob);CHKERRQ(ierr);
363a6cbbb48SMatthew G. Knepley   ierr = PetscMalloc4(NfNew, &tmpd, NfNew, &tmpdbd, NfNew, &tmpi, NfNew*2, &tmpa);CHKERRQ(ierr);
364a6cbbb48SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {tmpd[f] = prob->disc[f]; tmpdbd[f] = prob->discBd[f]; tmpi[f] = prob->implicit[f]; for (i = 0; i < 2; ++i) tmpa[f*2+i] = prob->adjacency[f*2+i];}
365a6cbbb48SMatthew G. Knepley   for (f = Nf; f < NfNew; ++f) {ierr = PetscContainerCreate(PetscObjectComm((PetscObject) prob), (PetscContainer *) &tmpd[f]);CHKERRQ(ierr);
366a6cbbb48SMatthew G. Knepley     tmpdbd[f] = NULL; tmpi[f] = PETSC_TRUE; tmpa[f*2+0] = PETSC_FALSE; tmpa[f*2+1] = PETSC_TRUE;}
367a6cbbb48SMatthew G. Knepley   ierr = PetscFree4(prob->disc, prob->discBd, prob->implicit, prob->adjacency);CHKERRQ(ierr);
3682764a2aaSMatthew G. Knepley   prob->Nf        = NfNew;
3692764a2aaSMatthew G. Knepley   prob->disc      = tmpd;
370a6cbbb48SMatthew G. Knepley   prob->discBd    = tmpdbd;
371249df284SMatthew G. Knepley   prob->implicit  = tmpi;
372a6cbbb48SMatthew G. Knepley   prob->adjacency = tmpa;
373*475e0ac9SMatthew G. Knepley   ierr = PetscCalloc6(NfNew, &tmpobj, NfNew*2, &tmpf, NfNew*NfNew*4, &tmpg, NfNew*NfNew*4, &tmpgp, NfNew, &tmpr, NfNew, &tmpctx);CHKERRQ(ierr);
3742764a2aaSMatthew G. Knepley   for (f = 0; f < Nf; ++f) tmpobj[f] = prob->obj[f];
3752764a2aaSMatthew G. Knepley   for (f = 0; f < Nf*2; ++f) tmpf[f] = prob->f[f];
3762764a2aaSMatthew G. Knepley   for (f = 0; f < Nf*Nf*4; ++f) tmpg[f] = prob->g[f];
377*475e0ac9SMatthew G. Knepley   for (f = 0; f < Nf*Nf*4; ++f) tmpgp[f] = prob->gp[f];
3780c2f2876SMatthew G. Knepley   for (f = 0; f < Nf; ++f) tmpr[f] = prob->r[f];
3790c2f2876SMatthew G. Knepley   for (f = 0; f < Nf; ++f) tmpctx[f] = prob->ctx[f];
3802764a2aaSMatthew G. Knepley   for (f = Nf; f < NfNew; ++f) tmpobj[f] = NULL;
3812764a2aaSMatthew G. Knepley   for (f = Nf*2; f < NfNew*2; ++f) tmpf[f] = NULL;
3822764a2aaSMatthew G. Knepley   for (f = Nf*Nf*4; f < NfNew*NfNew*4; ++f) tmpg[f] = NULL;
383*475e0ac9SMatthew G. Knepley   for (f = Nf*Nf*4; f < NfNew*NfNew*4; ++f) tmpgp[f] = NULL;
3840c2f2876SMatthew G. Knepley   for (f = Nf; f < NfNew; ++f) tmpr[f] = NULL;
3850c2f2876SMatthew G. Knepley   for (f = Nf; f < NfNew; ++f) tmpctx[f] = NULL;
386*475e0ac9SMatthew G. Knepley   ierr = PetscFree6(prob->obj, prob->f, prob->g, prob->gp, prob->r, prob->ctx);CHKERRQ(ierr);
3872764a2aaSMatthew G. Knepley   prob->obj = tmpobj;
3882764a2aaSMatthew G. Knepley   prob->f   = tmpf;
3892764a2aaSMatthew G. Knepley   prob->g   = tmpg;
390*475e0ac9SMatthew G. Knepley   prob->gp  = tmpgp;
3910c2f2876SMatthew G. Knepley   prob->r   = tmpr;
3920c2f2876SMatthew G. Knepley   prob->ctx = tmpctx;
3932764a2aaSMatthew G. Knepley   ierr = PetscCalloc2(NfNew*2, &tmpfbd, NfNew*NfNew*4, &tmpgbd);CHKERRQ(ierr);
3942764a2aaSMatthew G. Knepley   for (f = 0; f < Nf*2; ++f) tmpfbd[f] = prob->fBd[f];
3952764a2aaSMatthew G. Knepley   for (f = 0; f < Nf*Nf*4; ++f) tmpgbd[f] = prob->gBd[f];
3962764a2aaSMatthew G. Knepley   for (f = Nf*2; f < NfNew*2; ++f) tmpfbd[f] = NULL;
3972764a2aaSMatthew G. Knepley   for (f = Nf*Nf*4; f < NfNew*NfNew*4; ++f) tmpgbd[f] = NULL;
3982764a2aaSMatthew G. Knepley   ierr = PetscFree2(prob->fBd, prob->gBd);CHKERRQ(ierr);
3992764a2aaSMatthew G. Knepley   prob->fBd = tmpfbd;
4002764a2aaSMatthew G. Knepley   prob->gBd = tmpgbd;
4012764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
4022764a2aaSMatthew G. Knepley }
4032764a2aaSMatthew G. Knepley 
4042764a2aaSMatthew G. Knepley #undef __FUNCT__
4052764a2aaSMatthew G. Knepley #define __FUNCT__ "PetscDSDestroy"
4062764a2aaSMatthew G. Knepley /*@
4072764a2aaSMatthew G. Knepley   PetscDSDestroy - Destroys a PetscDS object
4082764a2aaSMatthew G. Knepley 
4092764a2aaSMatthew G. Knepley   Collective on PetscDS
4102764a2aaSMatthew G. Knepley 
4112764a2aaSMatthew G. Knepley   Input Parameter:
4122764a2aaSMatthew G. Knepley . prob - the PetscDS object to destroy
4132764a2aaSMatthew G. Knepley 
4142764a2aaSMatthew G. Knepley   Level: developer
4152764a2aaSMatthew G. Knepley 
4162764a2aaSMatthew G. Knepley .seealso PetscDSView()
4172764a2aaSMatthew G. Knepley @*/
4182764a2aaSMatthew G. Knepley PetscErrorCode PetscDSDestroy(PetscDS *prob)
4192764a2aaSMatthew G. Knepley {
4202764a2aaSMatthew G. Knepley   PetscInt       f;
4212764a2aaSMatthew G. Knepley   PetscErrorCode ierr;
4222764a2aaSMatthew G. Knepley 
4232764a2aaSMatthew G. Knepley   PetscFunctionBegin;
4242764a2aaSMatthew G. Knepley   if (!*prob) PetscFunctionReturn(0);
4252764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific((*prob), PETSCDS_CLASSID, 1);
4262764a2aaSMatthew G. Knepley 
4272764a2aaSMatthew G. Knepley   if (--((PetscObject)(*prob))->refct > 0) {*prob = 0; PetscFunctionReturn(0);}
4282764a2aaSMatthew G. Knepley   ((PetscObject) (*prob))->refct = 0;
4292764a2aaSMatthew G. Knepley   ierr = PetscDSDestroyStructs_Static(*prob);CHKERRQ(ierr);
4302764a2aaSMatthew G. Knepley   for (f = 0; f < (*prob)->Nf; ++f) {
4312764a2aaSMatthew G. Knepley     ierr = PetscObjectDereference((*prob)->disc[f]);CHKERRQ(ierr);
4322764a2aaSMatthew G. Knepley     ierr = PetscObjectDereference((*prob)->discBd[f]);CHKERRQ(ierr);
4332764a2aaSMatthew G. Knepley   }
434a6cbbb48SMatthew G. Knepley   ierr = PetscFree4((*prob)->disc, (*prob)->discBd, (*prob)->implicit, (*prob)->adjacency);CHKERRQ(ierr);
435*475e0ac9SMatthew G. Knepley   ierr = PetscFree6((*prob)->obj,(*prob)->f,(*prob)->g,(*prob)->gp,(*prob)->r,(*prob)->ctx);CHKERRQ(ierr);
4362764a2aaSMatthew G. Knepley   ierr = PetscFree2((*prob)->fBd,(*prob)->gBd);CHKERRQ(ierr);
4372764a2aaSMatthew G. Knepley   if ((*prob)->ops->destroy) {ierr = (*(*prob)->ops->destroy)(*prob);CHKERRQ(ierr);}
4382764a2aaSMatthew G. Knepley   ierr = PetscHeaderDestroy(prob);CHKERRQ(ierr);
4392764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
4402764a2aaSMatthew G. Knepley }
4412764a2aaSMatthew G. Knepley 
4422764a2aaSMatthew G. Knepley #undef __FUNCT__
4432764a2aaSMatthew G. Knepley #define __FUNCT__ "PetscDSCreate"
4442764a2aaSMatthew G. Knepley /*@
4452764a2aaSMatthew G. Knepley   PetscDSCreate - Creates an empty PetscDS object. The type can then be set with PetscDSSetType().
4462764a2aaSMatthew G. Knepley 
4472764a2aaSMatthew G. Knepley   Collective on MPI_Comm
4482764a2aaSMatthew G. Knepley 
4492764a2aaSMatthew G. Knepley   Input Parameter:
4502764a2aaSMatthew G. Knepley . comm - The communicator for the PetscDS object
4512764a2aaSMatthew G. Knepley 
4522764a2aaSMatthew G. Knepley   Output Parameter:
4532764a2aaSMatthew G. Knepley . prob - The PetscDS object
4542764a2aaSMatthew G. Knepley 
4552764a2aaSMatthew G. Knepley   Level: beginner
4562764a2aaSMatthew G. Knepley 
4572764a2aaSMatthew G. Knepley .seealso: PetscDSSetType(), PETSCDSBASIC
4582764a2aaSMatthew G. Knepley @*/
4592764a2aaSMatthew G. Knepley PetscErrorCode PetscDSCreate(MPI_Comm comm, PetscDS *prob)
4602764a2aaSMatthew G. Knepley {
4612764a2aaSMatthew G. Knepley   PetscDS   p;
4622764a2aaSMatthew G. Knepley   PetscErrorCode ierr;
4632764a2aaSMatthew G. Knepley 
4642764a2aaSMatthew G. Knepley   PetscFunctionBegin;
4652764a2aaSMatthew G. Knepley   PetscValidPointer(prob, 2);
4662764a2aaSMatthew G. Knepley   *prob  = NULL;
4672764a2aaSMatthew G. Knepley   ierr = PetscDSInitializePackage();CHKERRQ(ierr);
4682764a2aaSMatthew G. Knepley 
46973107ff1SLisandro Dalcin   ierr = PetscHeaderCreate(p, PETSCDS_CLASSID, "PetscDS", "Discrete System", "PetscDS", comm, PetscDSDestroy, PetscDSView);CHKERRQ(ierr);
4702764a2aaSMatthew G. Knepley 
4712764a2aaSMatthew G. Knepley   p->Nf    = 0;
4722764a2aaSMatthew G. Knepley   p->setup = PETSC_FALSE;
4732764a2aaSMatthew G. Knepley 
4742764a2aaSMatthew G. Knepley   *prob = p;
4752764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
4762764a2aaSMatthew G. Knepley }
4772764a2aaSMatthew G. Knepley 
4782764a2aaSMatthew G. Knepley #undef __FUNCT__
4792764a2aaSMatthew G. Knepley #define __FUNCT__ "PetscDSGetNumFields"
480bc4ae4beSMatthew G. Knepley /*@
481bc4ae4beSMatthew G. Knepley   PetscDSGetNumFields - Returns the number of fields in the DS
482bc4ae4beSMatthew G. Knepley 
483bc4ae4beSMatthew G. Knepley   Not collective
484bc4ae4beSMatthew G. Knepley 
485bc4ae4beSMatthew G. Knepley   Input Parameter:
486bc4ae4beSMatthew G. Knepley . prob - The PetscDS object
487bc4ae4beSMatthew G. Knepley 
488bc4ae4beSMatthew G. Knepley   Output Parameter:
489bc4ae4beSMatthew G. Knepley . Nf - The number of fields
490bc4ae4beSMatthew G. Knepley 
491bc4ae4beSMatthew G. Knepley   Level: beginner
492bc4ae4beSMatthew G. Knepley 
493bc4ae4beSMatthew G. Knepley .seealso: PetscDSGetSpatialDimension(), PetscDSCreate()
494bc4ae4beSMatthew G. Knepley @*/
4952764a2aaSMatthew G. Knepley PetscErrorCode PetscDSGetNumFields(PetscDS prob, PetscInt *Nf)
4962764a2aaSMatthew G. Knepley {
4972764a2aaSMatthew G. Knepley   PetscFunctionBegin;
4982764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
4992764a2aaSMatthew G. Knepley   PetscValidPointer(Nf, 2);
5002764a2aaSMatthew G. Knepley   *Nf = prob->Nf;
5012764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
5022764a2aaSMatthew G. Knepley }
5032764a2aaSMatthew G. Knepley 
5042764a2aaSMatthew G. Knepley #undef __FUNCT__
5052764a2aaSMatthew G. Knepley #define __FUNCT__ "PetscDSGetSpatialDimension"
506bc4ae4beSMatthew G. Knepley /*@
507bc4ae4beSMatthew G. Knepley   PetscDSGetSpatialDimension - Returns the spatial dimension of the DS
508bc4ae4beSMatthew G. Knepley 
509bc4ae4beSMatthew G. Knepley   Not collective
510bc4ae4beSMatthew G. Knepley 
511bc4ae4beSMatthew G. Knepley   Input Parameter:
512bc4ae4beSMatthew G. Knepley . prob - The PetscDS object
513bc4ae4beSMatthew G. Knepley 
514bc4ae4beSMatthew G. Knepley   Output Parameter:
515bc4ae4beSMatthew G. Knepley . dim - The spatial dimension
516bc4ae4beSMatthew G. Knepley 
517bc4ae4beSMatthew G. Knepley   Level: beginner
518bc4ae4beSMatthew G. Knepley 
519bc4ae4beSMatthew G. Knepley .seealso: PetscDSGetNumFields(), PetscDSCreate()
520bc4ae4beSMatthew G. Knepley @*/
5212764a2aaSMatthew G. Knepley PetscErrorCode PetscDSGetSpatialDimension(PetscDS prob, PetscInt *dim)
5222764a2aaSMatthew G. Knepley {
5232764a2aaSMatthew G. Knepley   PetscErrorCode ierr;
5242764a2aaSMatthew G. Knepley 
5252764a2aaSMatthew G. Knepley   PetscFunctionBegin;
5262764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
5272764a2aaSMatthew G. Knepley   PetscValidPointer(dim, 2);
5282764a2aaSMatthew G. Knepley   *dim = 0;
5299de99aefSMatthew G. Knepley   if (prob->Nf) {
5309de99aefSMatthew G. Knepley     PetscObject  obj;
5319de99aefSMatthew G. Knepley     PetscClassId id;
5329de99aefSMatthew G. Knepley 
5339de99aefSMatthew G. Knepley     ierr = PetscDSGetDiscretization(prob, 0, &obj);CHKERRQ(ierr);
5349de99aefSMatthew G. Knepley     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
5359de99aefSMatthew G. Knepley     if (id == PETSCFE_CLASSID)      {ierr = PetscFEGetSpatialDimension((PetscFE) obj, dim);CHKERRQ(ierr);}
5369de99aefSMatthew G. Knepley     else if (id == PETSCFV_CLASSID) {ierr = PetscFVGetSpatialDimension((PetscFV) obj, dim);CHKERRQ(ierr);}
5379de99aefSMatthew G. Knepley     else SETERRQ1(PetscObjectComm((PetscObject) prob), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", 0);
5389de99aefSMatthew G. Knepley   }
5392764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
5402764a2aaSMatthew G. Knepley }
5412764a2aaSMatthew G. Knepley 
5422764a2aaSMatthew G. Knepley #undef __FUNCT__
5432764a2aaSMatthew G. Knepley #define __FUNCT__ "PetscDSGetTotalDimension"
544bc4ae4beSMatthew G. Knepley /*@
545bc4ae4beSMatthew G. Knepley   PetscDSGetTotalDimension - Returns the total size of the approximation space for this system
546bc4ae4beSMatthew G. Knepley 
547bc4ae4beSMatthew G. Knepley   Not collective
548bc4ae4beSMatthew G. Knepley 
549bc4ae4beSMatthew G. Knepley   Input Parameter:
550bc4ae4beSMatthew G. Knepley . prob - The PetscDS object
551bc4ae4beSMatthew G. Knepley 
552bc4ae4beSMatthew G. Knepley   Output Parameter:
553bc4ae4beSMatthew G. Knepley . dim - The total problem dimension
554bc4ae4beSMatthew G. Knepley 
555bc4ae4beSMatthew G. Knepley   Level: beginner
556bc4ae4beSMatthew G. Knepley 
557bc4ae4beSMatthew G. Knepley .seealso: PetscDSGetNumFields(), PetscDSCreate()
558bc4ae4beSMatthew G. Knepley @*/
5592764a2aaSMatthew G. Knepley PetscErrorCode PetscDSGetTotalDimension(PetscDS prob, PetscInt *dim)
5602764a2aaSMatthew G. Knepley {
5612764a2aaSMatthew G. Knepley   PetscErrorCode ierr;
5622764a2aaSMatthew G. Knepley 
5632764a2aaSMatthew G. Knepley   PetscFunctionBegin;
5642764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
5652764a2aaSMatthew G. Knepley   ierr = PetscDSSetUp(prob);CHKERRQ(ierr);
5662764a2aaSMatthew G. Knepley   PetscValidPointer(dim, 2);
5672764a2aaSMatthew G. Knepley   *dim = prob->totDim;
5682764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
5692764a2aaSMatthew G. Knepley }
5702764a2aaSMatthew G. Knepley 
5712764a2aaSMatthew G. Knepley #undef __FUNCT__
5722764a2aaSMatthew G. Knepley #define __FUNCT__ "PetscDSGetTotalBdDimension"
573bc4ae4beSMatthew G. Knepley /*@
574c3ac4435SMatthew G. Knepley   PetscDSGetTotalBdDimension - Returns the total size of the boundary approximation space for this system
575bc4ae4beSMatthew G. Knepley 
576bc4ae4beSMatthew G. Knepley   Not collective
577bc4ae4beSMatthew G. Knepley 
578bc4ae4beSMatthew G. Knepley   Input Parameter:
579bc4ae4beSMatthew G. Knepley . prob - The PetscDS object
580bc4ae4beSMatthew G. Knepley 
581bc4ae4beSMatthew G. Knepley   Output Parameter:
582bc4ae4beSMatthew G. Knepley . dim - The total boundary problem dimension
583bc4ae4beSMatthew G. Knepley 
584bc4ae4beSMatthew G. Knepley   Level: beginner
585bc4ae4beSMatthew G. Knepley 
586bc4ae4beSMatthew G. Knepley .seealso: PetscDSGetNumFields(), PetscDSCreate()
587bc4ae4beSMatthew G. Knepley @*/
5882764a2aaSMatthew G. Knepley PetscErrorCode PetscDSGetTotalBdDimension(PetscDS prob, PetscInt *dim)
5892764a2aaSMatthew G. Knepley {
5902764a2aaSMatthew G. Knepley   PetscErrorCode ierr;
5912764a2aaSMatthew G. Knepley 
5922764a2aaSMatthew G. Knepley   PetscFunctionBegin;
5932764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
5942764a2aaSMatthew G. Knepley   ierr = PetscDSSetUp(prob);CHKERRQ(ierr);
5952764a2aaSMatthew G. Knepley   PetscValidPointer(dim, 2);
5962764a2aaSMatthew G. Knepley   *dim = prob->totDimBd;
5972764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
5982764a2aaSMatthew G. Knepley }
5992764a2aaSMatthew G. Knepley 
6002764a2aaSMatthew G. Knepley #undef __FUNCT__
6012764a2aaSMatthew G. Knepley #define __FUNCT__ "PetscDSGetTotalComponents"
602bc4ae4beSMatthew G. Knepley /*@
603bc4ae4beSMatthew G. Knepley   PetscDSGetTotalComponents - Returns the total number of components in this system
604bc4ae4beSMatthew G. Knepley 
605bc4ae4beSMatthew G. Knepley   Not collective
606bc4ae4beSMatthew G. Knepley 
607bc4ae4beSMatthew G. Knepley   Input Parameter:
608bc4ae4beSMatthew G. Knepley . prob - The PetscDS object
609bc4ae4beSMatthew G. Knepley 
610bc4ae4beSMatthew G. Knepley   Output Parameter:
611bc4ae4beSMatthew G. Knepley . dim - The total number of components
612bc4ae4beSMatthew G. Knepley 
613bc4ae4beSMatthew G. Knepley   Level: beginner
614bc4ae4beSMatthew G. Knepley 
615bc4ae4beSMatthew G. Knepley .seealso: PetscDSGetNumFields(), PetscDSCreate()
616bc4ae4beSMatthew G. Knepley @*/
6172764a2aaSMatthew G. Knepley PetscErrorCode PetscDSGetTotalComponents(PetscDS prob, PetscInt *Nc)
6182764a2aaSMatthew G. Knepley {
6192764a2aaSMatthew G. Knepley   PetscErrorCode ierr;
6202764a2aaSMatthew G. Knepley 
6212764a2aaSMatthew G. Knepley   PetscFunctionBegin;
6222764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
6232764a2aaSMatthew G. Knepley   ierr = PetscDSSetUp(prob);CHKERRQ(ierr);
6242764a2aaSMatthew G. Knepley   PetscValidPointer(Nc, 2);
6252764a2aaSMatthew G. Knepley   *Nc = prob->totComp;
6262764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
6272764a2aaSMatthew G. Knepley }
6282764a2aaSMatthew G. Knepley 
6292764a2aaSMatthew G. Knepley #undef __FUNCT__
6302764a2aaSMatthew G. Knepley #define __FUNCT__ "PetscDSGetDiscretization"
631bc4ae4beSMatthew G. Knepley /*@
632bc4ae4beSMatthew G. Knepley   PetscDSGetDiscretization - Returns the discretization object for the given field
633bc4ae4beSMatthew G. Knepley 
634bc4ae4beSMatthew G. Knepley   Not collective
635bc4ae4beSMatthew G. Knepley 
636bc4ae4beSMatthew G. Knepley   Input Parameters:
637bc4ae4beSMatthew G. Knepley + prob - The PetscDS object
638bc4ae4beSMatthew G. Knepley - f - The field number
639bc4ae4beSMatthew G. Knepley 
640bc4ae4beSMatthew G. Knepley   Output Parameter:
641bc4ae4beSMatthew G. Knepley . disc - The discretization object
642bc4ae4beSMatthew G. Knepley 
643bc4ae4beSMatthew G. Knepley   Level: beginner
644bc4ae4beSMatthew G. Knepley 
645bc4ae4beSMatthew G. Knepley .seealso: PetscDSSetDiscretization(), PetscDSAddDiscretization(), PetscDSGetBdDiscretization(), PetscDSGetNumFields(), PetscDSCreate()
646bc4ae4beSMatthew G. Knepley @*/
6472764a2aaSMatthew G. Knepley PetscErrorCode PetscDSGetDiscretization(PetscDS prob, PetscInt f, PetscObject *disc)
6482764a2aaSMatthew G. Knepley {
6492764a2aaSMatthew G. Knepley   PetscFunctionBegin;
6502764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
6512764a2aaSMatthew G. Knepley   PetscValidPointer(disc, 3);
6522764a2aaSMatthew 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);
6532764a2aaSMatthew G. Knepley   *disc = prob->disc[f];
6542764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
6552764a2aaSMatthew G. Knepley }
6562764a2aaSMatthew G. Knepley 
6572764a2aaSMatthew G. Knepley #undef __FUNCT__
6582764a2aaSMatthew G. Knepley #define __FUNCT__ "PetscDSGetBdDiscretization"
659bc4ae4beSMatthew G. Knepley /*@
660bc4ae4beSMatthew G. Knepley   PetscDSGetBdDiscretization - Returns the boundary discretization object for the given field
661bc4ae4beSMatthew G. Knepley 
662bc4ae4beSMatthew G. Knepley   Not collective
663bc4ae4beSMatthew G. Knepley 
664bc4ae4beSMatthew G. Knepley   Input Parameters:
665bc4ae4beSMatthew G. Knepley + prob - The PetscDS object
666bc4ae4beSMatthew G. Knepley - f - The field number
667bc4ae4beSMatthew G. Knepley 
668bc4ae4beSMatthew G. Knepley   Output Parameter:
669bc4ae4beSMatthew G. Knepley . disc - The boundary discretization object
670bc4ae4beSMatthew G. Knepley 
671bc4ae4beSMatthew G. Knepley   Level: beginner
672bc4ae4beSMatthew G. Knepley 
673bc4ae4beSMatthew G. Knepley .seealso: PetscDSSetBdDiscretization(), PetscDSAddBdDiscretization(), PetscDSGetDiscretization(), PetscDSGetNumFields(), PetscDSCreate()
674bc4ae4beSMatthew G. Knepley @*/
6752764a2aaSMatthew G. Knepley PetscErrorCode PetscDSGetBdDiscretization(PetscDS prob, PetscInt f, PetscObject *disc)
6762764a2aaSMatthew G. Knepley {
6772764a2aaSMatthew G. Knepley   PetscFunctionBegin;
6782764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
6792764a2aaSMatthew G. Knepley   PetscValidPointer(disc, 3);
6802764a2aaSMatthew 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);
6812764a2aaSMatthew G. Knepley   *disc = prob->discBd[f];
6822764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
6832764a2aaSMatthew G. Knepley }
6842764a2aaSMatthew G. Knepley 
6852764a2aaSMatthew G. Knepley #undef __FUNCT__
6862764a2aaSMatthew G. Knepley #define __FUNCT__ "PetscDSSetDiscretization"
687bc4ae4beSMatthew G. Knepley /*@
688bc4ae4beSMatthew G. Knepley   PetscDSSetDiscretization - Sets the discretization object for the given field
689bc4ae4beSMatthew G. Knepley 
690bc4ae4beSMatthew G. Knepley   Not collective
691bc4ae4beSMatthew G. Knepley 
692bc4ae4beSMatthew G. Knepley   Input Parameters:
693bc4ae4beSMatthew G. Knepley + prob - The PetscDS object
694bc4ae4beSMatthew G. Knepley . f - The field number
695bc4ae4beSMatthew G. Knepley - disc - The discretization object
696bc4ae4beSMatthew G. Knepley 
697bc4ae4beSMatthew G. Knepley   Level: beginner
698bc4ae4beSMatthew G. Knepley 
699bc4ae4beSMatthew G. Knepley .seealso: PetscDSGetDiscretization(), PetscDSAddDiscretization(), PetscDSGetNumFields(), PetscDSCreate()
700bc4ae4beSMatthew G. Knepley @*/
7012764a2aaSMatthew G. Knepley PetscErrorCode PetscDSSetDiscretization(PetscDS prob, PetscInt f, PetscObject disc)
7022764a2aaSMatthew G. Knepley {
7032764a2aaSMatthew G. Knepley   PetscErrorCode ierr;
7042764a2aaSMatthew G. Knepley 
7052764a2aaSMatthew G. Knepley   PetscFunctionBegin;
7062764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
7072764a2aaSMatthew G. Knepley   PetscValidPointer(disc, 3);
7082764a2aaSMatthew G. Knepley   if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f);
7092764a2aaSMatthew G. Knepley   ierr = PetscDSEnlarge_Static(prob, f+1);CHKERRQ(ierr);
7102764a2aaSMatthew G. Knepley   if (prob->disc[f]) {ierr = PetscObjectDereference(prob->disc[f]);CHKERRQ(ierr);}
7112764a2aaSMatthew G. Knepley   prob->disc[f] = disc;
7122764a2aaSMatthew G. Knepley   ierr = PetscObjectReference(disc);CHKERRQ(ierr);
713249df284SMatthew G. Knepley   {
714249df284SMatthew G. Knepley     PetscClassId id;
715249df284SMatthew G. Knepley 
716249df284SMatthew G. Knepley     ierr = PetscObjectGetClassId(disc, &id);CHKERRQ(ierr);
717a6cbbb48SMatthew G. Knepley     if (id == PETSCFV_CLASSID) {
718a6cbbb48SMatthew G. Knepley       prob->implicit[f]      = PETSC_FALSE;
719a6cbbb48SMatthew G. Knepley       prob->adjacency[f*2+0] = PETSC_TRUE;
720a6cbbb48SMatthew G. Knepley       prob->adjacency[f*2+1] = PETSC_FALSE;
721a6cbbb48SMatthew G. Knepley     }
722249df284SMatthew G. Knepley   }
7232764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
7242764a2aaSMatthew G. Knepley }
7252764a2aaSMatthew G. Knepley 
7262764a2aaSMatthew G. Knepley #undef __FUNCT__
7272764a2aaSMatthew G. Knepley #define __FUNCT__ "PetscDSSetBdDiscretization"
728bc4ae4beSMatthew G. Knepley /*@
729bc4ae4beSMatthew G. Knepley   PetscDSSetBdDiscretization - Sets the boundary discretization object for the given field
730bc4ae4beSMatthew G. Knepley 
731bc4ae4beSMatthew G. Knepley   Not collective
732bc4ae4beSMatthew G. Knepley 
733bc4ae4beSMatthew G. Knepley   Input Parameters:
734bc4ae4beSMatthew G. Knepley + prob - The PetscDS object
735bc4ae4beSMatthew G. Knepley . f - The field number
736bc4ae4beSMatthew G. Knepley - disc - The boundary discretization object
737bc4ae4beSMatthew G. Knepley 
738bc4ae4beSMatthew G. Knepley   Level: beginner
739bc4ae4beSMatthew G. Knepley 
740bc4ae4beSMatthew G. Knepley .seealso: PetscDSGetBdDiscretization(), PetscDSAddBdDiscretization(), PetscDSGetNumFields(), PetscDSCreate()
741bc4ae4beSMatthew G. Knepley @*/
7422764a2aaSMatthew G. Knepley PetscErrorCode PetscDSSetBdDiscretization(PetscDS prob, PetscInt f, PetscObject disc)
7432764a2aaSMatthew G. Knepley {
7442764a2aaSMatthew G. Knepley   PetscErrorCode ierr;
7452764a2aaSMatthew G. Knepley 
7462764a2aaSMatthew G. Knepley   PetscFunctionBegin;
7472764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
7482764a2aaSMatthew G. Knepley   if (disc) PetscValidPointer(disc, 3);
7492764a2aaSMatthew G. Knepley   if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f);
7502764a2aaSMatthew G. Knepley   ierr = PetscDSEnlarge_Static(prob, f+1);CHKERRQ(ierr);
7512764a2aaSMatthew G. Knepley   if (prob->discBd[f]) {ierr = PetscObjectDereference(prob->discBd[f]);CHKERRQ(ierr);}
7522764a2aaSMatthew G. Knepley   prob->discBd[f] = disc;
7532764a2aaSMatthew G. Knepley   ierr = PetscObjectReference(disc);CHKERRQ(ierr);
7542764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
7552764a2aaSMatthew G. Knepley }
7562764a2aaSMatthew G. Knepley 
7572764a2aaSMatthew G. Knepley #undef __FUNCT__
7582764a2aaSMatthew G. Knepley #define __FUNCT__ "PetscDSAddDiscretization"
759bc4ae4beSMatthew G. Knepley /*@
760bc4ae4beSMatthew G. Knepley   PetscDSAddDiscretization - Adds a discretization object
761bc4ae4beSMatthew G. Knepley 
762bc4ae4beSMatthew G. Knepley   Not collective
763bc4ae4beSMatthew G. Knepley 
764bc4ae4beSMatthew G. Knepley   Input Parameters:
765bc4ae4beSMatthew G. Knepley + prob - The PetscDS object
766bc4ae4beSMatthew G. Knepley - disc - The boundary discretization object
767bc4ae4beSMatthew G. Knepley 
768bc4ae4beSMatthew G. Knepley   Level: beginner
769bc4ae4beSMatthew G. Knepley 
770bc4ae4beSMatthew G. Knepley .seealso: PetscDSGetDiscretization(), PetscDSSetDiscretization(), PetscDSGetNumFields(), PetscDSCreate()
771bc4ae4beSMatthew G. Knepley @*/
7722764a2aaSMatthew G. Knepley PetscErrorCode PetscDSAddDiscretization(PetscDS prob, PetscObject disc)
7732764a2aaSMatthew G. Knepley {
7742764a2aaSMatthew G. Knepley   PetscErrorCode ierr;
7752764a2aaSMatthew G. Knepley 
7762764a2aaSMatthew G. Knepley   PetscFunctionBegin;
7772764a2aaSMatthew G. Knepley   ierr = PetscDSSetDiscretization(prob, prob->Nf, disc);CHKERRQ(ierr);
7782764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
7792764a2aaSMatthew G. Knepley }
7802764a2aaSMatthew G. Knepley 
7812764a2aaSMatthew G. Knepley #undef __FUNCT__
7822764a2aaSMatthew G. Knepley #define __FUNCT__ "PetscDSAddBdDiscretization"
783bc4ae4beSMatthew G. Knepley /*@
784bc4ae4beSMatthew G. Knepley   PetscDSAddBdDiscretization - Adds a boundary discretization object
785bc4ae4beSMatthew G. Knepley 
786bc4ae4beSMatthew G. Knepley   Not collective
787bc4ae4beSMatthew G. Knepley 
788bc4ae4beSMatthew G. Knepley   Input Parameters:
789bc4ae4beSMatthew G. Knepley + prob - The PetscDS object
790bc4ae4beSMatthew G. Knepley - disc - The boundary discretization object
791bc4ae4beSMatthew G. Knepley 
792bc4ae4beSMatthew G. Knepley   Level: beginner
793bc4ae4beSMatthew G. Knepley 
794bc4ae4beSMatthew G. Knepley .seealso: PetscDSGetBdDiscretization(), PetscDSSetBdDiscretization(), PetscDSGetNumFields(), PetscDSCreate()
795bc4ae4beSMatthew G. Knepley @*/
7962764a2aaSMatthew G. Knepley PetscErrorCode PetscDSAddBdDiscretization(PetscDS prob, PetscObject disc)
7972764a2aaSMatthew G. Knepley {
7982764a2aaSMatthew G. Knepley   PetscErrorCode ierr;
7992764a2aaSMatthew G. Knepley 
8002764a2aaSMatthew G. Knepley   PetscFunctionBegin;
8012764a2aaSMatthew G. Knepley   ierr = PetscDSSetBdDiscretization(prob, prob->Nf, disc);CHKERRQ(ierr);
8022764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
8032764a2aaSMatthew G. Knepley }
8042764a2aaSMatthew G. Knepley 
8052764a2aaSMatthew G. Knepley #undef __FUNCT__
806249df284SMatthew G. Knepley #define __FUNCT__ "PetscDSGetImplicit"
807249df284SMatthew G. Knepley /*@
808249df284SMatthew G. Knepley   PetscDSGetImplicit - Returns the flag for implicit solve for this field. This is just a guide for IMEX
809249df284SMatthew G. Knepley 
810249df284SMatthew G. Knepley   Not collective
811249df284SMatthew G. Knepley 
812249df284SMatthew G. Knepley   Input Parameters:
813249df284SMatthew G. Knepley + prob - The PetscDS object
814249df284SMatthew G. Knepley - f - The field number
815249df284SMatthew G. Knepley 
816249df284SMatthew G. Knepley   Output Parameter:
817249df284SMatthew G. Knepley . implicit - The flag indicating what kind of solve to use for this field
818249df284SMatthew G. Knepley 
819249df284SMatthew G. Knepley   Level: developer
820249df284SMatthew G. Knepley 
821249df284SMatthew G. Knepley .seealso: PetscDSSetImplicit(), PetscDSSetDiscretization(), PetscDSAddDiscretization(), PetscDSGetBdDiscretization(), PetscDSGetNumFields(), PetscDSCreate()
822249df284SMatthew G. Knepley @*/
823249df284SMatthew G. Knepley PetscErrorCode PetscDSGetImplicit(PetscDS prob, PetscInt f, PetscBool *implicit)
824249df284SMatthew G. Knepley {
825249df284SMatthew G. Knepley   PetscFunctionBegin;
826249df284SMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
827249df284SMatthew G. Knepley   PetscValidPointer(implicit, 3);
828249df284SMatthew 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);
829249df284SMatthew G. Knepley   *implicit = prob->implicit[f];
830249df284SMatthew G. Knepley   PetscFunctionReturn(0);
831249df284SMatthew G. Knepley }
832249df284SMatthew G. Knepley 
833249df284SMatthew G. Knepley #undef __FUNCT__
834249df284SMatthew G. Knepley #define __FUNCT__ "PetscDSSetImplicit"
835249df284SMatthew G. Knepley /*@
836249df284SMatthew G. Knepley   PetscDSSetImplicit - Set the flag for implicit solve for this field. This is just a guide for IMEX
837249df284SMatthew G. Knepley 
838249df284SMatthew G. Knepley   Not collective
839249df284SMatthew G. Knepley 
840249df284SMatthew G. Knepley   Input Parameters:
841249df284SMatthew G. Knepley + prob - The PetscDS object
842249df284SMatthew G. Knepley . f - The field number
843249df284SMatthew G. Knepley - implicit - The flag indicating what kind of solve to use for this field
844249df284SMatthew G. Knepley 
845249df284SMatthew G. Knepley   Level: developer
846249df284SMatthew G. Knepley 
847249df284SMatthew G. Knepley .seealso: PetscDSGetImplicit(), PetscDSSetDiscretization(), PetscDSAddDiscretization(), PetscDSGetBdDiscretization(), PetscDSGetNumFields(), PetscDSCreate()
848249df284SMatthew G. Knepley @*/
849249df284SMatthew G. Knepley PetscErrorCode PetscDSSetImplicit(PetscDS prob, PetscInt f, PetscBool implicit)
850249df284SMatthew G. Knepley {
851249df284SMatthew G. Knepley   PetscFunctionBegin;
852249df284SMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
853249df284SMatthew 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);
854249df284SMatthew G. Knepley   prob->implicit[f] = implicit;
855249df284SMatthew G. Knepley   PetscFunctionReturn(0);
856249df284SMatthew G. Knepley }
857249df284SMatthew G. Knepley 
858249df284SMatthew G. Knepley #undef __FUNCT__
859a6cbbb48SMatthew G. Knepley #define __FUNCT__ "PetscDSGetAdjacency"
860a6cbbb48SMatthew G. Knepley /*@
861a6cbbb48SMatthew G. Knepley   PetscDSGetAdjacency - Returns the flags for determining variable influence
862a6cbbb48SMatthew G. Knepley 
863a6cbbb48SMatthew G. Knepley   Not collective
864a6cbbb48SMatthew G. Knepley 
865a6cbbb48SMatthew G. Knepley   Input Parameters:
866a6cbbb48SMatthew G. Knepley + prob - The PetscDS object
867a6cbbb48SMatthew G. Knepley - f - The field number
868a6cbbb48SMatthew G. Knepley 
869a6cbbb48SMatthew G. Knepley   Output Parameter:
870a6cbbb48SMatthew G. Knepley + useCone    - Flag for variable influence starting with the cone operation
871a6cbbb48SMatthew G. Knepley - useClosure - Flag for variable influence using transitive closure
872a6cbbb48SMatthew G. Knepley 
873a6cbbb48SMatthew G. Knepley   Note: See the discussion in DMPlexGetAdjacencyUseCone() and DMPlexGetAdjacencyUseClosure()
874a6cbbb48SMatthew G. Knepley 
875a6cbbb48SMatthew G. Knepley   Level: developer
876a6cbbb48SMatthew G. Knepley 
877a6cbbb48SMatthew G. Knepley .seealso: PetscDSSetAdjacency(), DMPlexGetAdjacencyUseCone(), DMPlexGetAdjacencyUseClosure(), PetscDSSetDiscretization(), PetscDSAddDiscretization(), PetscDSGetBdDiscretization(), PetscDSGetNumFields(), PetscDSCreate()
878a6cbbb48SMatthew G. Knepley @*/
879a6cbbb48SMatthew G. Knepley PetscErrorCode PetscDSGetAdjacency(PetscDS prob, PetscInt f, PetscBool *useCone, PetscBool *useClosure)
880a6cbbb48SMatthew G. Knepley {
881a6cbbb48SMatthew G. Knepley   PetscFunctionBegin;
882a6cbbb48SMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
883a6cbbb48SMatthew G. Knepley   PetscValidPointer(useCone, 3);
884a6cbbb48SMatthew G. Knepley   PetscValidPointer(useClosure, 4);
885a6cbbb48SMatthew 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);
886a6cbbb48SMatthew G. Knepley   *useCone    = prob->adjacency[f*2+0];
887a6cbbb48SMatthew G. Knepley   *useClosure = prob->adjacency[f*2+1];
888a6cbbb48SMatthew G. Knepley   PetscFunctionReturn(0);
889a6cbbb48SMatthew G. Knepley }
890a6cbbb48SMatthew G. Knepley 
891a6cbbb48SMatthew G. Knepley #undef __FUNCT__
892a6cbbb48SMatthew G. Knepley #define __FUNCT__ "PetscDSSetAdjacency"
893a6cbbb48SMatthew G. Knepley /*@
894a6cbbb48SMatthew G. Knepley   PetscDSSetAdjacency - Set the flags for determining variable influence
895a6cbbb48SMatthew G. Knepley 
896a6cbbb48SMatthew G. Knepley   Not collective
897a6cbbb48SMatthew G. Knepley 
898a6cbbb48SMatthew G. Knepley   Input Parameters:
899a6cbbb48SMatthew G. Knepley + prob - The PetscDS object
900a6cbbb48SMatthew G. Knepley . f - The field number
901a6cbbb48SMatthew G. Knepley . useCone    - Flag for variable influence starting with the cone operation
902a6cbbb48SMatthew G. Knepley - useClosure - Flag for variable influence using transitive closure
903a6cbbb48SMatthew G. Knepley 
904a6cbbb48SMatthew G. Knepley   Note: See the discussion in DMPlexGetAdjacencyUseCone() and DMPlexGetAdjacencyUseClosure()
905a6cbbb48SMatthew G. Knepley 
906a6cbbb48SMatthew G. Knepley   Level: developer
907a6cbbb48SMatthew G. Knepley 
908a6cbbb48SMatthew G. Knepley .seealso: PetscDSGetAdjacency(), DMPlexGetAdjacencyUseCone(), DMPlexGetAdjacencyUseClosure(), PetscDSSetDiscretization(), PetscDSAddDiscretization(), PetscDSGetBdDiscretization(), PetscDSGetNumFields(), PetscDSCreate()
909a6cbbb48SMatthew G. Knepley @*/
910a6cbbb48SMatthew G. Knepley PetscErrorCode PetscDSSetAdjacency(PetscDS prob, PetscInt f, PetscBool useCone, PetscBool useClosure)
911a6cbbb48SMatthew G. Knepley {
912a6cbbb48SMatthew G. Knepley   PetscFunctionBegin;
913a6cbbb48SMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
914a6cbbb48SMatthew 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);
915a6cbbb48SMatthew G. Knepley   prob->adjacency[f*2+0] = useCone;
916a6cbbb48SMatthew G. Knepley   prob->adjacency[f*2+1] = useClosure;
917a6cbbb48SMatthew G. Knepley   PetscFunctionReturn(0);
918a6cbbb48SMatthew G. Knepley }
919a6cbbb48SMatthew G. Knepley 
920a6cbbb48SMatthew G. Knepley #undef __FUNCT__
9212764a2aaSMatthew G. Knepley #define __FUNCT__ "PetscDSGetObjective"
9222764a2aaSMatthew G. Knepley PetscErrorCode PetscDSGetObjective(PetscDS prob, PetscInt f,
92330b9ff8bSMatthew G. Knepley                                    void (**obj)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
924194d53e6SMatthew G. Knepley                                                 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
925194d53e6SMatthew G. Knepley                                                 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
92630b9ff8bSMatthew G. Knepley                                                 PetscReal t, const PetscReal x[], PetscScalar obj[]))
9272764a2aaSMatthew G. Knepley {
9282764a2aaSMatthew G. Knepley   PetscFunctionBegin;
9292764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
9302764a2aaSMatthew G. Knepley   PetscValidPointer(obj, 2);
9312764a2aaSMatthew 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);
9322764a2aaSMatthew G. Knepley   *obj = prob->obj[f];
9332764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
9342764a2aaSMatthew G. Knepley }
9352764a2aaSMatthew G. Knepley 
9362764a2aaSMatthew G. Knepley #undef __FUNCT__
9372764a2aaSMatthew G. Knepley #define __FUNCT__ "PetscDSSetObjective"
9382764a2aaSMatthew G. Knepley PetscErrorCode PetscDSSetObjective(PetscDS prob, PetscInt f,
93930b9ff8bSMatthew G. Knepley                                    void (*obj)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
940194d53e6SMatthew G. Knepley                                                const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
941194d53e6SMatthew G. Knepley                                                const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
94230b9ff8bSMatthew G. Knepley                                                PetscReal t, const PetscReal x[], PetscScalar obj[]))
9432764a2aaSMatthew G. Knepley {
9442764a2aaSMatthew G. Knepley   PetscErrorCode ierr;
9452764a2aaSMatthew G. Knepley 
9462764a2aaSMatthew G. Knepley   PetscFunctionBegin;
9472764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
9482764a2aaSMatthew G. Knepley   PetscValidFunction(obj, 2);
9492764a2aaSMatthew G. Knepley   if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f);
9502764a2aaSMatthew G. Knepley   ierr = PetscDSEnlarge_Static(prob, f+1);CHKERRQ(ierr);
9512764a2aaSMatthew G. Knepley   prob->obj[f] = obj;
9522764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
9532764a2aaSMatthew G. Knepley }
9542764a2aaSMatthew G. Knepley 
9552764a2aaSMatthew G. Knepley #undef __FUNCT__
9562764a2aaSMatthew G. Knepley #define __FUNCT__ "PetscDSGetResidual"
957194d53e6SMatthew G. Knepley /*@C
958194d53e6SMatthew G. Knepley   PetscDSGetResidual - Get the pointwise residual function for a given test field
959194d53e6SMatthew G. Knepley 
960194d53e6SMatthew G. Knepley   Not collective
961194d53e6SMatthew G. Knepley 
962194d53e6SMatthew G. Knepley   Input Parameters:
963194d53e6SMatthew G. Knepley + prob - The PetscDS
964194d53e6SMatthew G. Knepley - f    - The test field number
965194d53e6SMatthew G. Knepley 
966194d53e6SMatthew G. Knepley   Output Parameters:
967194d53e6SMatthew G. Knepley + f0 - integrand for the test function term
968194d53e6SMatthew G. Knepley - f1 - integrand for the test function gradient term
969194d53e6SMatthew G. Knepley 
970194d53e6SMatthew G. Knepley   Note: We are using a first order FEM model for the weak form:
971194d53e6SMatthew G. Knepley 
972194d53e6SMatthew 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)
973194d53e6SMatthew G. Knepley 
974194d53e6SMatthew G. Knepley The calling sequence for the callbacks f0 and f1 is given by:
975194d53e6SMatthew G. Knepley 
97630b9ff8bSMatthew G. Knepley $ f0(PetscInt dim, PetscInt Nf, PetscInt NfAux,
977194d53e6SMatthew G. Knepley $    const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
978194d53e6SMatthew G. Knepley $    const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
97930b9ff8bSMatthew G. Knepley $    PetscReal t, const PetscReal x[], PetscScalar f0[])
980194d53e6SMatthew G. Knepley 
981194d53e6SMatthew G. Knepley + dim - the spatial dimension
982194d53e6SMatthew G. Knepley . Nf - the number of fields
983194d53e6SMatthew G. Knepley . uOff - the offset into u[] and u_t[] for each field
984194d53e6SMatthew G. Knepley . uOff_x - the offset into u_x[] for each field
985194d53e6SMatthew G. Knepley . u - each field evaluated at the current point
986194d53e6SMatthew G. Knepley . u_t - the time derivative of each field evaluated at the current point
987194d53e6SMatthew G. Knepley . u_x - the gradient of each field evaluated at the current point
988194d53e6SMatthew G. Knepley . aOff - the offset into a[] and a_t[] for each auxiliary field
989194d53e6SMatthew G. Knepley . aOff_x - the offset into a_x[] for each auxiliary field
990194d53e6SMatthew G. Knepley . a - each auxiliary field evaluated at the current point
991194d53e6SMatthew G. Knepley . a_t - the time derivative of each auxiliary field evaluated at the current point
992194d53e6SMatthew G. Knepley . a_x - the gradient of auxiliary each field evaluated at the current point
993194d53e6SMatthew G. Knepley . t - current time
994194d53e6SMatthew G. Knepley . x - coordinates of the current point
995194d53e6SMatthew G. Knepley - f0 - output values at the current point
996194d53e6SMatthew G. Knepley 
997194d53e6SMatthew G. Knepley   Level: intermediate
998194d53e6SMatthew G. Knepley 
999194d53e6SMatthew G. Knepley .seealso: PetscDSSetResidual()
1000194d53e6SMatthew G. Knepley @*/
10012764a2aaSMatthew G. Knepley PetscErrorCode PetscDSGetResidual(PetscDS prob, PetscInt f,
100230b9ff8bSMatthew G. Knepley                                   void (**f0)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1003194d53e6SMatthew G. Knepley                                               const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1004194d53e6SMatthew G. Knepley                                               const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
100530b9ff8bSMatthew G. Knepley                                               PetscReal t, const PetscReal x[], PetscScalar f0[]),
100630b9ff8bSMatthew G. Knepley                                   void (**f1)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1007194d53e6SMatthew G. Knepley                                               const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1008194d53e6SMatthew G. Knepley                                               const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
100930b9ff8bSMatthew G. Knepley                                               PetscReal t, const PetscReal x[], PetscScalar f1[]))
10102764a2aaSMatthew G. Knepley {
10112764a2aaSMatthew G. Knepley   PetscFunctionBegin;
10122764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
10132764a2aaSMatthew 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);
10142764a2aaSMatthew G. Knepley   if (f0) {PetscValidPointer(f0, 3); *f0 = prob->f[f*2+0];}
10152764a2aaSMatthew G. Knepley   if (f1) {PetscValidPointer(f1, 4); *f1 = prob->f[f*2+1];}
10162764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
10172764a2aaSMatthew G. Knepley }
10182764a2aaSMatthew G. Knepley 
10192764a2aaSMatthew G. Knepley #undef __FUNCT__
10202764a2aaSMatthew G. Knepley #define __FUNCT__ "PetscDSSetResidual"
1021194d53e6SMatthew G. Knepley /*@C
1022194d53e6SMatthew G. Knepley   PetscDSSetResidual - Set the pointwise residual function for a given test field
1023194d53e6SMatthew G. Knepley 
1024194d53e6SMatthew G. Knepley   Not collective
1025194d53e6SMatthew G. Knepley 
1026194d53e6SMatthew G. Knepley   Input Parameters:
1027194d53e6SMatthew G. Knepley + prob - The PetscDS
1028194d53e6SMatthew G. Knepley . f    - The test field number
1029194d53e6SMatthew G. Knepley . f0 - integrand for the test function term
1030194d53e6SMatthew G. Knepley - f1 - integrand for the test function gradient term
1031194d53e6SMatthew G. Knepley 
1032194d53e6SMatthew G. Knepley   Note: We are using a first order FEM model for the weak form:
1033194d53e6SMatthew G. Knepley 
1034194d53e6SMatthew 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)
1035194d53e6SMatthew G. Knepley 
1036194d53e6SMatthew G. Knepley The calling sequence for the callbacks f0 and f1 is given by:
1037194d53e6SMatthew G. Knepley 
103830b9ff8bSMatthew G. Knepley $ f0(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1039194d53e6SMatthew G. Knepley $    const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1040194d53e6SMatthew G. Knepley $    const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
104130b9ff8bSMatthew G. Knepley $    PetscReal t, const PetscReal x[], PetscScalar f0[])
1042194d53e6SMatthew G. Knepley 
1043194d53e6SMatthew G. Knepley + dim - the spatial dimension
1044194d53e6SMatthew G. Knepley . Nf - the number of fields
1045194d53e6SMatthew G. Knepley . uOff - the offset into u[] and u_t[] for each field
1046194d53e6SMatthew G. Knepley . uOff_x - the offset into u_x[] for each field
1047194d53e6SMatthew G. Knepley . u - each field evaluated at the current point
1048194d53e6SMatthew G. Knepley . u_t - the time derivative of each field evaluated at the current point
1049194d53e6SMatthew G. Knepley . u_x - the gradient of each field evaluated at the current point
1050194d53e6SMatthew G. Knepley . aOff - the offset into a[] and a_t[] for each auxiliary field
1051194d53e6SMatthew G. Knepley . aOff_x - the offset into a_x[] for each auxiliary field
1052194d53e6SMatthew G. Knepley . a - each auxiliary field evaluated at the current point
1053194d53e6SMatthew G. Knepley . a_t - the time derivative of each auxiliary field evaluated at the current point
1054194d53e6SMatthew G. Knepley . a_x - the gradient of auxiliary each field evaluated at the current point
1055194d53e6SMatthew G. Knepley . t - current time
1056194d53e6SMatthew G. Knepley . x - coordinates of the current point
1057194d53e6SMatthew G. Knepley - f0 - output values at the current point
1058194d53e6SMatthew G. Knepley 
1059194d53e6SMatthew G. Knepley   Level: intermediate
1060194d53e6SMatthew G. Knepley 
1061194d53e6SMatthew G. Knepley .seealso: PetscDSGetResidual()
1062194d53e6SMatthew G. Knepley @*/
10632764a2aaSMatthew G. Knepley PetscErrorCode PetscDSSetResidual(PetscDS prob, PetscInt f,
106430b9ff8bSMatthew G. Knepley                                   void (*f0)(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[],
106730b9ff8bSMatthew G. Knepley                                              PetscReal t, const PetscReal x[], PetscScalar f0[]),
106830b9ff8bSMatthew G. Knepley                                   void (*f1)(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[],
107130b9ff8bSMatthew G. Knepley                                              PetscReal t, const PetscReal x[], PetscScalar f1[]))
10722764a2aaSMatthew G. Knepley {
10732764a2aaSMatthew G. Knepley   PetscErrorCode ierr;
10742764a2aaSMatthew G. Knepley 
10752764a2aaSMatthew G. Knepley   PetscFunctionBegin;
10762764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1077f866a1d0SMatthew G. Knepley   if (f0) PetscValidFunction(f0, 3);
1078f866a1d0SMatthew G. Knepley   if (f1) PetscValidFunction(f1, 4);
10792764a2aaSMatthew G. Knepley   if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f);
10802764a2aaSMatthew G. Knepley   ierr = PetscDSEnlarge_Static(prob, f+1);CHKERRQ(ierr);
10812764a2aaSMatthew G. Knepley   prob->f[f*2+0] = f0;
10822764a2aaSMatthew G. Knepley   prob->f[f*2+1] = f1;
10832764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
10842764a2aaSMatthew G. Knepley }
10852764a2aaSMatthew G. Knepley 
10862764a2aaSMatthew G. Knepley #undef __FUNCT__
10872764a2aaSMatthew G. Knepley #define __FUNCT__ "PetscDSGetJacobian"
1088194d53e6SMatthew G. Knepley /*@C
1089194d53e6SMatthew G. Knepley   PetscDSGetJacobian - Get the pointwise Jacobian function for given test and basis field
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 
1098194d53e6SMatthew G. Knepley   Output Parameters:
1099194d53e6SMatthew G. Knepley + g0 - integrand for the test and basis function term
1100194d53e6SMatthew G. Knepley . g1 - integrand for the test function and basis function gradient term
1101194d53e6SMatthew G. Knepley . g2 - integrand for the test function gradient and basis function term
1102194d53e6SMatthew G. Knepley - g3 - integrand for the test function gradient and basis function gradient term
1103194d53e6SMatthew G. Knepley 
1104194d53e6SMatthew G. Knepley   Note: We are using a first order FEM model for the weak form:
1105194d53e6SMatthew G. Knepley 
1106194d53e6SMatthew 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
1107194d53e6SMatthew G. Knepley 
1108194d53e6SMatthew G. Knepley The calling sequence for the callbacks g0, g1, g2 and g3 is given by:
1109194d53e6SMatthew G. Knepley 
111030b9ff8bSMatthew G. Knepley $ g0(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1111194d53e6SMatthew G. Knepley $    const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1112194d53e6SMatthew G. Knepley $    const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
111330b9ff8bSMatthew G. Knepley $    PetscReal t, const PetscReal u_tShift, const PetscReal x[], PetscScalar g0[])
1114194d53e6SMatthew G. Knepley 
1115194d53e6SMatthew G. Knepley + dim - the spatial dimension
1116194d53e6SMatthew G. Knepley . Nf - the number of fields
1117194d53e6SMatthew G. Knepley . uOff - the offset into u[] and u_t[] for each field
1118194d53e6SMatthew G. Knepley . uOff_x - the offset into u_x[] for each field
1119194d53e6SMatthew G. Knepley . u - each field evaluated at the current point
1120194d53e6SMatthew G. Knepley . u_t - the time derivative of each field evaluated at the current point
1121194d53e6SMatthew G. Knepley . u_x - the gradient of each field evaluated at the current point
1122194d53e6SMatthew G. Knepley . aOff - the offset into a[] and a_t[] for each auxiliary field
1123194d53e6SMatthew G. Knepley . aOff_x - the offset into a_x[] for each auxiliary field
1124194d53e6SMatthew G. Knepley . a - each auxiliary field evaluated at the current point
1125194d53e6SMatthew G. Knepley . a_t - the time derivative of each auxiliary field evaluated at the current point
1126194d53e6SMatthew G. Knepley . a_x - the gradient of auxiliary each field evaluated at the current point
1127194d53e6SMatthew G. Knepley . t - current time
11282aa1fc23SMatthew G. Knepley . u_tShift - the multiplier a for dF/dU_t
1129194d53e6SMatthew G. Knepley . x - coordinates of the current point
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: PetscDSSetJacobian()
1135194d53e6SMatthew G. Knepley @*/
11362764a2aaSMatthew G. Knepley PetscErrorCode PetscDSGetJacobian(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[],
11402aa1fc23SMatthew G. Knepley                                               PetscReal t, PetscReal u_tShift, const PetscReal x[], 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[],
11442aa1fc23SMatthew G. Knepley                                               PetscReal t, PetscReal u_tShift, const PetscReal x[], 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[],
11482aa1fc23SMatthew G. Knepley                                               PetscReal t, PetscReal u_tShift, const PetscReal x[], 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[],
11522aa1fc23SMatthew G. Knepley                                               PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscScalar g3[]))
11532764a2aaSMatthew G. Knepley {
11542764a2aaSMatthew G. Knepley   PetscFunctionBegin;
11552764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
11562764a2aaSMatthew 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);
11572764a2aaSMatthew 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);
11582764a2aaSMatthew G. Knepley   if (g0) {PetscValidPointer(g0, 4); *g0 = prob->g[(f*prob->Nf + g)*4+0];}
11592764a2aaSMatthew G. Knepley   if (g1) {PetscValidPointer(g1, 5); *g1 = prob->g[(f*prob->Nf + g)*4+1];}
11602764a2aaSMatthew G. Knepley   if (g2) {PetscValidPointer(g2, 6); *g2 = prob->g[(f*prob->Nf + g)*4+2];}
11612764a2aaSMatthew G. Knepley   if (g3) {PetscValidPointer(g3, 7); *g3 = prob->g[(f*prob->Nf + g)*4+3];}
11622764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
11632764a2aaSMatthew G. Knepley }
11642764a2aaSMatthew G. Knepley 
11652764a2aaSMatthew G. Knepley #undef __FUNCT__
11662764a2aaSMatthew G. Knepley #define __FUNCT__ "PetscDSSetJacobian"
1167194d53e6SMatthew G. Knepley /*@C
1168194d53e6SMatthew G. Knepley   PetscDSSetJacobian - Set the pointwise Jacobian function for given test and basis fields
1169194d53e6SMatthew G. Knepley 
1170194d53e6SMatthew G. Knepley   Not collective
1171194d53e6SMatthew G. Knepley 
1172194d53e6SMatthew G. Knepley   Input Parameters:
1173194d53e6SMatthew G. Knepley + prob - The PetscDS
1174194d53e6SMatthew G. Knepley . f    - The test field number
1175194d53e6SMatthew G. Knepley . g    - The field number
1176194d53e6SMatthew G. Knepley . g0 - integrand for the test and basis function term
1177194d53e6SMatthew G. Knepley . g1 - integrand for the test function and basis function gradient term
1178194d53e6SMatthew G. Knepley . g2 - integrand for the test function gradient and basis function term
1179194d53e6SMatthew G. Knepley - g3 - integrand for the test function gradient and basis function gradient term
1180194d53e6SMatthew G. Knepley 
1181194d53e6SMatthew G. Knepley   Note: We are using a first order FEM model for the weak form:
1182194d53e6SMatthew G. Knepley 
1183194d53e6SMatthew 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
1184194d53e6SMatthew G. Knepley 
1185194d53e6SMatthew G. Knepley The calling sequence for the callbacks g0, g1, g2 and g3 is given by:
1186194d53e6SMatthew G. Knepley 
118730b9ff8bSMatthew G. Knepley $ g0(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1188194d53e6SMatthew G. Knepley $    const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1189194d53e6SMatthew G. Knepley $    const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
119030b9ff8bSMatthew G. Knepley $    PetscReal t, const PetscReal x[], PetscScalar g0[])
1191194d53e6SMatthew G. Knepley 
1192194d53e6SMatthew G. Knepley + dim - the spatial dimension
1193194d53e6SMatthew G. Knepley . Nf - the number of fields
1194194d53e6SMatthew G. Knepley . uOff - the offset into u[] and u_t[] for each field
1195194d53e6SMatthew G. Knepley . uOff_x - the offset into u_x[] for each field
1196194d53e6SMatthew G. Knepley . u - each field evaluated at the current point
1197194d53e6SMatthew G. Knepley . u_t - the time derivative of each field evaluated at the current point
1198194d53e6SMatthew G. Knepley . u_x - the gradient of each field evaluated at the current point
1199194d53e6SMatthew G. Knepley . aOff - the offset into a[] and a_t[] for each auxiliary field
1200194d53e6SMatthew G. Knepley . aOff_x - the offset into a_x[] for each auxiliary field
1201194d53e6SMatthew G. Knepley . a - each auxiliary field evaluated at the current point
1202194d53e6SMatthew G. Knepley . a_t - the time derivative of each auxiliary field evaluated at the current point
1203194d53e6SMatthew G. Knepley . a_x - the gradient of auxiliary each field evaluated at the current point
1204194d53e6SMatthew G. Knepley . t - current time
12052aa1fc23SMatthew G. Knepley . u_tShift - the multiplier a for dF/dU_t
1206194d53e6SMatthew G. Knepley . x - coordinates of the current point
1207194d53e6SMatthew G. Knepley - g0 - output values at the current point
1208194d53e6SMatthew G. Knepley 
1209194d53e6SMatthew G. Knepley   Level: intermediate
1210194d53e6SMatthew G. Knepley 
1211194d53e6SMatthew G. Knepley .seealso: PetscDSGetJacobian()
1212194d53e6SMatthew G. Knepley @*/
12132764a2aaSMatthew G. Knepley PetscErrorCode PetscDSSetJacobian(PetscDS prob, PetscInt f, PetscInt g,
121430b9ff8bSMatthew G. Knepley                                   void (*g0)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1215194d53e6SMatthew G. Knepley                                              const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1216194d53e6SMatthew G. Knepley                                              const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
121730b9ff8bSMatthew G. Knepley                                              PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscScalar g0[]),
121830b9ff8bSMatthew G. Knepley                                   void (*g1)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1219194d53e6SMatthew G. Knepley                                              const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1220194d53e6SMatthew G. Knepley                                              const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
122130b9ff8bSMatthew G. Knepley                                              PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscScalar g1[]),
122230b9ff8bSMatthew G. Knepley                                   void (*g2)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1223194d53e6SMatthew G. Knepley                                              const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1224194d53e6SMatthew G. Knepley                                              const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
122530b9ff8bSMatthew G. Knepley                                              PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscScalar g2[]),
122630b9ff8bSMatthew G. Knepley                                   void (*g3)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1227194d53e6SMatthew G. Knepley                                              const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1228194d53e6SMatthew G. Knepley                                              const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
122930b9ff8bSMatthew G. Knepley                                              PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscScalar g3[]))
12302764a2aaSMatthew G. Knepley {
12312764a2aaSMatthew G. Knepley   PetscErrorCode ierr;
12322764a2aaSMatthew G. Knepley 
12332764a2aaSMatthew G. Knepley   PetscFunctionBegin;
12342764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
12352764a2aaSMatthew G. Knepley   if (g0) PetscValidFunction(g0, 4);
12362764a2aaSMatthew G. Knepley   if (g1) PetscValidFunction(g1, 5);
12372764a2aaSMatthew G. Knepley   if (g2) PetscValidFunction(g2, 6);
12382764a2aaSMatthew G. Knepley   if (g3) PetscValidFunction(g3, 7);
12392764a2aaSMatthew G. Knepley   if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f);
12402764a2aaSMatthew G. Knepley   if (g < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", g);
12412764a2aaSMatthew G. Knepley   ierr = PetscDSEnlarge_Static(prob, PetscMax(f, g)+1);CHKERRQ(ierr);
12422764a2aaSMatthew G. Knepley   prob->g[(f*prob->Nf + g)*4+0] = g0;
12432764a2aaSMatthew G. Knepley   prob->g[(f*prob->Nf + g)*4+1] = g1;
12442764a2aaSMatthew G. Knepley   prob->g[(f*prob->Nf + g)*4+2] = g2;
12452764a2aaSMatthew G. Knepley   prob->g[(f*prob->Nf + g)*4+3] = g3;
12462764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
12472764a2aaSMatthew G. Knepley }
12482764a2aaSMatthew G. Knepley 
12492764a2aaSMatthew G. Knepley #undef __FUNCT__
1250*475e0ac9SMatthew G. Knepley #define __FUNCT__ "PetscDSHasJacobianPreconditioner"
1251*475e0ac9SMatthew G. Knepley /*@C
1252*475e0ac9SMatthew G. Knepley   PetscDSHasJacobianPreconditioner - Signals that a Jacobian preconditioner matrix has been set
1253*475e0ac9SMatthew G. Knepley 
1254*475e0ac9SMatthew G. Knepley   Not collective
1255*475e0ac9SMatthew G. Knepley 
1256*475e0ac9SMatthew G. Knepley   Input Parameter:
1257*475e0ac9SMatthew G. Knepley . prob - The PetscDS
1258*475e0ac9SMatthew G. Knepley 
1259*475e0ac9SMatthew G. Knepley   Output Parameter:
1260*475e0ac9SMatthew G. Knepley . hasJacPre - flag that pointwise function for Jacobian preconditioner matrix has been set
1261*475e0ac9SMatthew G. Knepley 
1262*475e0ac9SMatthew G. Knepley   Level: intermediate
1263*475e0ac9SMatthew G. Knepley 
1264*475e0ac9SMatthew G. Knepley .seealso: PetscDSGetJacobianPreconditioner(), PetscDSSetJacobianPreconditioner(), PetscDSGetJacobian()
1265*475e0ac9SMatthew G. Knepley @*/
1266*475e0ac9SMatthew G. Knepley PetscErrorCode PetscDSHasJacobianPreconditioner(PetscDS prob, PetscBool *hasJacPre)
1267*475e0ac9SMatthew G. Knepley {
1268*475e0ac9SMatthew G. Knepley   PetscInt f, g, h;
1269*475e0ac9SMatthew G. Knepley 
1270*475e0ac9SMatthew G. Knepley   PetscFunctionBegin;
1271*475e0ac9SMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1272*475e0ac9SMatthew G. Knepley   *hasJacPre = PETSC_FALSE;
1273*475e0ac9SMatthew G. Knepley   for (f = 0; f < prob->Nf; ++f) {
1274*475e0ac9SMatthew G. Knepley     for (g = 0; g < prob->Nf; ++g) {
1275*475e0ac9SMatthew G. Knepley       for (h = 0; h < 4; ++h) {
1276*475e0ac9SMatthew G. Knepley         if (prob->gp[(f*prob->Nf + g)*4+h]) *hasJacPre = PETSC_TRUE;
1277*475e0ac9SMatthew G. Knepley       }
1278*475e0ac9SMatthew G. Knepley     }
1279*475e0ac9SMatthew G. Knepley   }
1280*475e0ac9SMatthew G. Knepley   PetscFunctionReturn(0);
1281*475e0ac9SMatthew G. Knepley }
1282*475e0ac9SMatthew G. Knepley 
1283*475e0ac9SMatthew G. Knepley #undef __FUNCT__
1284*475e0ac9SMatthew G. Knepley #define __FUNCT__ "PetscDSGetJacobianPreconditioner"
1285*475e0ac9SMatthew G. Knepley /*@C
1286*475e0ac9SMatthew 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.
1287*475e0ac9SMatthew G. Knepley 
1288*475e0ac9SMatthew G. Knepley   Not collective
1289*475e0ac9SMatthew G. Knepley 
1290*475e0ac9SMatthew G. Knepley   Input Parameters:
1291*475e0ac9SMatthew G. Knepley + prob - The PetscDS
1292*475e0ac9SMatthew G. Knepley . f    - The test field number
1293*475e0ac9SMatthew G. Knepley - g    - The field number
1294*475e0ac9SMatthew G. Knepley 
1295*475e0ac9SMatthew G. Knepley   Output Parameters:
1296*475e0ac9SMatthew G. Knepley + g0 - integrand for the test and basis function term
1297*475e0ac9SMatthew G. Knepley . g1 - integrand for the test function and basis function gradient term
1298*475e0ac9SMatthew G. Knepley . g2 - integrand for the test function gradient and basis function term
1299*475e0ac9SMatthew G. Knepley - g3 - integrand for the test function gradient and basis function gradient term
1300*475e0ac9SMatthew G. Knepley 
1301*475e0ac9SMatthew G. Knepley   Note: We are using a first order FEM model for the weak form:
1302*475e0ac9SMatthew G. Knepley 
1303*475e0ac9SMatthew 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
1304*475e0ac9SMatthew G. Knepley 
1305*475e0ac9SMatthew G. Knepley The calling sequence for the callbacks g0, g1, g2 and g3 is given by:
1306*475e0ac9SMatthew G. Knepley 
1307*475e0ac9SMatthew G. Knepley $ g0(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1308*475e0ac9SMatthew G. Knepley $    const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1309*475e0ac9SMatthew G. Knepley $    const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1310*475e0ac9SMatthew G. Knepley $    PetscReal t, const PetscReal u_tShift, const PetscReal x[], PetscScalar g0[])
1311*475e0ac9SMatthew G. Knepley 
1312*475e0ac9SMatthew G. Knepley + dim - the spatial dimension
1313*475e0ac9SMatthew G. Knepley . Nf - the number of fields
1314*475e0ac9SMatthew G. Knepley . uOff - the offset into u[] and u_t[] for each field
1315*475e0ac9SMatthew G. Knepley . uOff_x - the offset into u_x[] for each field
1316*475e0ac9SMatthew G. Knepley . u - each field evaluated at the current point
1317*475e0ac9SMatthew G. Knepley . u_t - the time derivative of each field evaluated at the current point
1318*475e0ac9SMatthew G. Knepley . u_x - the gradient of each field evaluated at the current point
1319*475e0ac9SMatthew G. Knepley . aOff - the offset into a[] and a_t[] for each auxiliary field
1320*475e0ac9SMatthew G. Knepley . aOff_x - the offset into a_x[] for each auxiliary field
1321*475e0ac9SMatthew G. Knepley . a - each auxiliary field evaluated at the current point
1322*475e0ac9SMatthew G. Knepley . a_t - the time derivative of each auxiliary field evaluated at the current point
1323*475e0ac9SMatthew G. Knepley . a_x - the gradient of auxiliary each field evaluated at the current point
1324*475e0ac9SMatthew G. Knepley . t - current time
1325*475e0ac9SMatthew G. Knepley . u_tShift - the multiplier a for dF/dU_t
1326*475e0ac9SMatthew G. Knepley . x - coordinates of the current point
1327*475e0ac9SMatthew G. Knepley - g0 - output values at the current point
1328*475e0ac9SMatthew G. Knepley 
1329*475e0ac9SMatthew G. Knepley   Level: intermediate
1330*475e0ac9SMatthew G. Knepley 
1331*475e0ac9SMatthew G. Knepley .seealso: PetscDSSetJacobianPreconditioner(), PetscDSGetJacobian()
1332*475e0ac9SMatthew G. Knepley @*/
1333*475e0ac9SMatthew G. Knepley PetscErrorCode PetscDSGetJacobianPreconditioner(PetscDS prob, PetscInt f, PetscInt g,
1334*475e0ac9SMatthew G. Knepley                                   void (**g0)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1335*475e0ac9SMatthew G. Knepley                                               const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1336*475e0ac9SMatthew G. Knepley                                               const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1337*475e0ac9SMatthew G. Knepley                                               PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscScalar g0[]),
1338*475e0ac9SMatthew G. Knepley                                   void (**g1)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1339*475e0ac9SMatthew G. Knepley                                               const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1340*475e0ac9SMatthew G. Knepley                                               const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1341*475e0ac9SMatthew G. Knepley                                               PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscScalar g1[]),
1342*475e0ac9SMatthew G. Knepley                                   void (**g2)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1343*475e0ac9SMatthew G. Knepley                                               const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1344*475e0ac9SMatthew G. Knepley                                               const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1345*475e0ac9SMatthew G. Knepley                                               PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscScalar g2[]),
1346*475e0ac9SMatthew G. Knepley                                   void (**g3)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1347*475e0ac9SMatthew G. Knepley                                               const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1348*475e0ac9SMatthew G. Knepley                                               const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1349*475e0ac9SMatthew G. Knepley                                               PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscScalar g3[]))
1350*475e0ac9SMatthew G. Knepley {
1351*475e0ac9SMatthew G. Knepley   PetscFunctionBegin;
1352*475e0ac9SMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1353*475e0ac9SMatthew 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);
1354*475e0ac9SMatthew 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);
1355*475e0ac9SMatthew G. Knepley   if (g0) {PetscValidPointer(g0, 4); *g0 = prob->gp[(f*prob->Nf + g)*4+0];}
1356*475e0ac9SMatthew G. Knepley   if (g1) {PetscValidPointer(g1, 5); *g1 = prob->gp[(f*prob->Nf + g)*4+1];}
1357*475e0ac9SMatthew G. Knepley   if (g2) {PetscValidPointer(g2, 6); *g2 = prob->gp[(f*prob->Nf + g)*4+2];}
1358*475e0ac9SMatthew G. Knepley   if (g3) {PetscValidPointer(g3, 7); *g3 = prob->gp[(f*prob->Nf + g)*4+3];}
1359*475e0ac9SMatthew G. Knepley   PetscFunctionReturn(0);
1360*475e0ac9SMatthew G. Knepley }
1361*475e0ac9SMatthew G. Knepley 
1362*475e0ac9SMatthew G. Knepley #undef __FUNCT__
1363*475e0ac9SMatthew G. Knepley #define __FUNCT__ "PetscDSSetJacobianPreconditioner"
1364*475e0ac9SMatthew G. Knepley /*@C
1365*475e0ac9SMatthew 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.
1366*475e0ac9SMatthew G. Knepley 
1367*475e0ac9SMatthew G. Knepley   Not collective
1368*475e0ac9SMatthew G. Knepley 
1369*475e0ac9SMatthew G. Knepley   Input Parameters:
1370*475e0ac9SMatthew G. Knepley + prob - The PetscDS
1371*475e0ac9SMatthew G. Knepley . f    - The test field number
1372*475e0ac9SMatthew G. Knepley . g    - The field number
1373*475e0ac9SMatthew G. Knepley . g0 - integrand for the test and basis function term
1374*475e0ac9SMatthew G. Knepley . g1 - integrand for the test function and basis function gradient term
1375*475e0ac9SMatthew G. Knepley . g2 - integrand for the test function gradient and basis function term
1376*475e0ac9SMatthew G. Knepley - g3 - integrand for the test function gradient and basis function gradient term
1377*475e0ac9SMatthew G. Knepley 
1378*475e0ac9SMatthew G. Knepley   Note: We are using a first order FEM model for the weak form:
1379*475e0ac9SMatthew G. Knepley 
1380*475e0ac9SMatthew 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
1381*475e0ac9SMatthew G. Knepley 
1382*475e0ac9SMatthew G. Knepley The calling sequence for the callbacks g0, g1, g2 and g3 is given by:
1383*475e0ac9SMatthew G. Knepley 
1384*475e0ac9SMatthew G. Knepley $ g0(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1385*475e0ac9SMatthew G. Knepley $    const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1386*475e0ac9SMatthew G. Knepley $    const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1387*475e0ac9SMatthew G. Knepley $    PetscReal t, const PetscReal x[], PetscScalar g0[])
1388*475e0ac9SMatthew G. Knepley 
1389*475e0ac9SMatthew G. Knepley + dim - the spatial dimension
1390*475e0ac9SMatthew G. Knepley . Nf - the number of fields
1391*475e0ac9SMatthew G. Knepley . uOff - the offset into u[] and u_t[] for each field
1392*475e0ac9SMatthew G. Knepley . uOff_x - the offset into u_x[] for each field
1393*475e0ac9SMatthew G. Knepley . u - each field evaluated at the current point
1394*475e0ac9SMatthew G. Knepley . u_t - the time derivative of each field evaluated at the current point
1395*475e0ac9SMatthew G. Knepley . u_x - the gradient of each field evaluated at the current point
1396*475e0ac9SMatthew G. Knepley . aOff - the offset into a[] and a_t[] for each auxiliary field
1397*475e0ac9SMatthew G. Knepley . aOff_x - the offset into a_x[] for each auxiliary field
1398*475e0ac9SMatthew G. Knepley . a - each auxiliary field evaluated at the current point
1399*475e0ac9SMatthew G. Knepley . a_t - the time derivative of each auxiliary field evaluated at the current point
1400*475e0ac9SMatthew G. Knepley . a_x - the gradient of auxiliary each field evaluated at the current point
1401*475e0ac9SMatthew G. Knepley . t - current time
1402*475e0ac9SMatthew G. Knepley . u_tShift - the multiplier a for dF/dU_t
1403*475e0ac9SMatthew G. Knepley . x - coordinates of the current point
1404*475e0ac9SMatthew G. Knepley - g0 - output values at the current point
1405*475e0ac9SMatthew G. Knepley 
1406*475e0ac9SMatthew G. Knepley   Level: intermediate
1407*475e0ac9SMatthew G. Knepley 
1408*475e0ac9SMatthew G. Knepley .seealso: PetscDSGetJacobianPreconditioner(), PetscDSSetJacobian()
1409*475e0ac9SMatthew G. Knepley @*/
1410*475e0ac9SMatthew G. Knepley PetscErrorCode PetscDSSetJacobianPreconditioner(PetscDS prob, PetscInt f, PetscInt g,
1411*475e0ac9SMatthew G. Knepley                                   void (*g0)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1412*475e0ac9SMatthew G. Knepley                                              const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1413*475e0ac9SMatthew G. Knepley                                              const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1414*475e0ac9SMatthew G. Knepley                                              PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscScalar g0[]),
1415*475e0ac9SMatthew G. Knepley                                   void (*g1)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1416*475e0ac9SMatthew G. Knepley                                              const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1417*475e0ac9SMatthew G. Knepley                                              const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1418*475e0ac9SMatthew G. Knepley                                              PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscScalar g1[]),
1419*475e0ac9SMatthew G. Knepley                                   void (*g2)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1420*475e0ac9SMatthew G. Knepley                                              const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1421*475e0ac9SMatthew G. Knepley                                              const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1422*475e0ac9SMatthew G. Knepley                                              PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscScalar g2[]),
1423*475e0ac9SMatthew G. Knepley                                   void (*g3)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1424*475e0ac9SMatthew G. Knepley                                              const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1425*475e0ac9SMatthew G. Knepley                                              const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1426*475e0ac9SMatthew G. Knepley                                              PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscScalar g3[]))
1427*475e0ac9SMatthew G. Knepley {
1428*475e0ac9SMatthew G. Knepley   PetscErrorCode ierr;
1429*475e0ac9SMatthew G. Knepley 
1430*475e0ac9SMatthew G. Knepley   PetscFunctionBegin;
1431*475e0ac9SMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1432*475e0ac9SMatthew G. Knepley   if (g0) PetscValidFunction(g0, 4);
1433*475e0ac9SMatthew G. Knepley   if (g1) PetscValidFunction(g1, 5);
1434*475e0ac9SMatthew G. Knepley   if (g2) PetscValidFunction(g2, 6);
1435*475e0ac9SMatthew G. Knepley   if (g3) PetscValidFunction(g3, 7);
1436*475e0ac9SMatthew G. Knepley   if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f);
1437*475e0ac9SMatthew G. Knepley   if (g < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", g);
1438*475e0ac9SMatthew G. Knepley   ierr = PetscDSEnlarge_Static(prob, PetscMax(f, g)+1);CHKERRQ(ierr);
1439*475e0ac9SMatthew G. Knepley   prob->gp[(f*prob->Nf + g)*4+0] = g0;
1440*475e0ac9SMatthew G. Knepley   prob->gp[(f*prob->Nf + g)*4+1] = g1;
1441*475e0ac9SMatthew G. Knepley   prob->gp[(f*prob->Nf + g)*4+2] = g2;
1442*475e0ac9SMatthew G. Knepley   prob->gp[(f*prob->Nf + g)*4+3] = g3;
1443*475e0ac9SMatthew G. Knepley   PetscFunctionReturn(0);
1444*475e0ac9SMatthew G. Knepley }
1445*475e0ac9SMatthew G. Knepley 
1446*475e0ac9SMatthew G. Knepley #undef __FUNCT__
14470c2f2876SMatthew G. Knepley #define __FUNCT__ "PetscDSGetRiemannSolver"
14480c2f2876SMatthew G. Knepley /*@C
14490c2f2876SMatthew G. Knepley   PetscDSGetRiemannSolver - Returns the Riemann solver for the given field
14500c2f2876SMatthew G. Knepley 
14510c2f2876SMatthew G. Knepley   Not collective
14520c2f2876SMatthew G. Knepley 
14530c2f2876SMatthew G. Knepley   Input Arguments:
14540c2f2876SMatthew G. Knepley + prob - The PetscDS object
14550c2f2876SMatthew G. Knepley - f    - The field number
14560c2f2876SMatthew G. Knepley 
14570c2f2876SMatthew G. Knepley   Output Argument:
14580c2f2876SMatthew G. Knepley . r    - Riemann solver
14590c2f2876SMatthew G. Knepley 
14600c2f2876SMatthew G. Knepley   Calling sequence for r:
14610c2f2876SMatthew G. Knepley 
14625db36cf9SMatthew G. Knepley $ r(PetscInt dim, PetscInt Nf, const PetscReal x[], const PetscReal n[], const PetscScalar uL[], const PetscScalar uR[], PetscScalar flux[], void *ctx)
14630c2f2876SMatthew G. Knepley 
14645db36cf9SMatthew G. Knepley + dim  - The spatial dimension
14655db36cf9SMatthew G. Knepley . Nf   - The number of fields
14665db36cf9SMatthew G. Knepley . x    - The coordinates at a point on the interface
14670c2f2876SMatthew G. Knepley . n    - The normal vector to the interface
14680c2f2876SMatthew G. Knepley . uL   - The state vector to the left of the interface
14690c2f2876SMatthew G. Knepley . uR   - The state vector to the right of the interface
14700c2f2876SMatthew G. Knepley . flux - output array of flux through the interface
14710c2f2876SMatthew G. Knepley - ctx  - optional user context
14720c2f2876SMatthew G. Knepley 
14730c2f2876SMatthew G. Knepley   Level: intermediate
14740c2f2876SMatthew G. Knepley 
14750c2f2876SMatthew G. Knepley .seealso: PetscDSSetRiemannSolver()
14760c2f2876SMatthew G. Knepley @*/
14770c2f2876SMatthew G. Knepley PetscErrorCode PetscDSGetRiemannSolver(PetscDS prob, PetscInt f,
14785db36cf9SMatthew G. Knepley                                        void (**r)(PetscInt dim, PetscInt Nf, const PetscReal x[], const PetscReal n[], const PetscScalar uL[], const PetscScalar uR[], PetscScalar flux[], void *ctx))
14790c2f2876SMatthew G. Knepley {
14800c2f2876SMatthew G. Knepley   PetscFunctionBegin;
14810c2f2876SMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
14820c2f2876SMatthew 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);
14830c2f2876SMatthew G. Knepley   PetscValidPointer(r, 3);
14840c2f2876SMatthew G. Knepley   *r = prob->r[f];
14850c2f2876SMatthew G. Knepley   PetscFunctionReturn(0);
14860c2f2876SMatthew G. Knepley }
14870c2f2876SMatthew G. Knepley 
14880c2f2876SMatthew G. Knepley #undef __FUNCT__
14890c2f2876SMatthew G. Knepley #define __FUNCT__ "PetscDSSetRiemannSolver"
14900c2f2876SMatthew G. Knepley /*@C
14910c2f2876SMatthew G. Knepley   PetscDSSetRiemannSolver - Sets the Riemann solver for the given field
14920c2f2876SMatthew G. Knepley 
14930c2f2876SMatthew G. Knepley   Not collective
14940c2f2876SMatthew G. Knepley 
14950c2f2876SMatthew G. Knepley   Input Arguments:
14960c2f2876SMatthew G. Knepley + prob - The PetscDS object
14970c2f2876SMatthew G. Knepley . f    - The field number
14980c2f2876SMatthew G. Knepley - r    - Riemann solver
14990c2f2876SMatthew G. Knepley 
15000c2f2876SMatthew G. Knepley   Calling sequence for r:
15010c2f2876SMatthew G. Knepley 
15025db36cf9SMatthew G. Knepley $ r(PetscInt dim, PetscInt Nf, const PetscReal x[], const PetscReal n[], const PetscScalar uL[], const PetscScalar uR[], PetscScalar flux[], void *ctx)
15030c2f2876SMatthew G. Knepley 
15045db36cf9SMatthew G. Knepley + dim  - The spatial dimension
15055db36cf9SMatthew G. Knepley . Nf   - The number of fields
15065db36cf9SMatthew G. Knepley . x    - The coordinates at a point on the interface
15070c2f2876SMatthew G. Knepley . n    - The normal vector to the interface
15080c2f2876SMatthew G. Knepley . uL   - The state vector to the left of the interface
15090c2f2876SMatthew G. Knepley . uR   - The state vector to the right of the interface
15100c2f2876SMatthew G. Knepley . flux - output array of flux through the interface
15110c2f2876SMatthew G. Knepley - ctx  - optional user context
15120c2f2876SMatthew G. Knepley 
15130c2f2876SMatthew G. Knepley   Level: intermediate
15140c2f2876SMatthew G. Knepley 
15150c2f2876SMatthew G. Knepley .seealso: PetscDSGetRiemannSolver()
15160c2f2876SMatthew G. Knepley @*/
15170c2f2876SMatthew G. Knepley PetscErrorCode PetscDSSetRiemannSolver(PetscDS prob, PetscInt f,
15185db36cf9SMatthew G. Knepley                                        void (*r)(PetscInt dim, PetscInt Nf, const PetscReal x[], const PetscReal n[], const PetscScalar uL[], const PetscScalar uR[], PetscScalar flux[], void *ctx))
15190c2f2876SMatthew G. Knepley {
15200c2f2876SMatthew G. Knepley   PetscErrorCode ierr;
15210c2f2876SMatthew G. Knepley 
15220c2f2876SMatthew G. Knepley   PetscFunctionBegin;
15230c2f2876SMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
15240c2f2876SMatthew G. Knepley   PetscValidFunction(r, 3);
15250c2f2876SMatthew G. Knepley   if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f);
15260c2f2876SMatthew G. Knepley   ierr = PetscDSEnlarge_Static(prob, f+1);CHKERRQ(ierr);
15270c2f2876SMatthew G. Knepley   prob->r[f] = r;
15280c2f2876SMatthew G. Knepley   PetscFunctionReturn(0);
15290c2f2876SMatthew G. Knepley }
15300c2f2876SMatthew G. Knepley 
15310c2f2876SMatthew G. Knepley #undef __FUNCT__
15320c2f2876SMatthew G. Knepley #define __FUNCT__ "PetscDSGetContext"
15330c2f2876SMatthew G. Knepley PetscErrorCode PetscDSGetContext(PetscDS prob, PetscInt f, void **ctx)
15340c2f2876SMatthew G. Knepley {
15350c2f2876SMatthew G. Knepley   PetscFunctionBegin;
15360c2f2876SMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
15370c2f2876SMatthew 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);
15380c2f2876SMatthew G. Knepley   PetscValidPointer(ctx, 3);
15390c2f2876SMatthew G. Knepley   *ctx = prob->ctx[f];
15400c2f2876SMatthew G. Knepley   PetscFunctionReturn(0);
15410c2f2876SMatthew G. Knepley }
15420c2f2876SMatthew G. Knepley 
15430c2f2876SMatthew G. Knepley #undef __FUNCT__
15440c2f2876SMatthew G. Knepley #define __FUNCT__ "PetscDSSetContext"
15450c2f2876SMatthew G. Knepley PetscErrorCode PetscDSSetContext(PetscDS prob, PetscInt f, void *ctx)
15460c2f2876SMatthew G. Knepley {
15470c2f2876SMatthew G. Knepley   PetscErrorCode ierr;
15480c2f2876SMatthew G. Knepley 
15490c2f2876SMatthew G. Knepley   PetscFunctionBegin;
15500c2f2876SMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
15510c2f2876SMatthew G. Knepley   if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f);
15520c2f2876SMatthew G. Knepley   ierr = PetscDSEnlarge_Static(prob, f+1);CHKERRQ(ierr);
15530c2f2876SMatthew G. Knepley   prob->ctx[f] = ctx;
15540c2f2876SMatthew G. Knepley   PetscFunctionReturn(0);
15550c2f2876SMatthew G. Knepley }
15560c2f2876SMatthew G. Knepley 
15570c2f2876SMatthew G. Knepley #undef __FUNCT__
15582764a2aaSMatthew G. Knepley #define __FUNCT__ "PetscDSGetBdResidual"
1559194d53e6SMatthew G. Knepley /*@C
1560194d53e6SMatthew G. Knepley   PetscDSGetBdResidual - Get the pointwise boundary residual function for a given test field
1561194d53e6SMatthew G. Knepley 
1562194d53e6SMatthew G. Knepley   Not collective
1563194d53e6SMatthew G. Knepley 
1564194d53e6SMatthew G. Knepley   Input Parameters:
1565194d53e6SMatthew G. Knepley + prob - The PetscDS
1566194d53e6SMatthew G. Knepley - f    - The test field number
1567194d53e6SMatthew G. Knepley 
1568194d53e6SMatthew G. Knepley   Output Parameters:
1569194d53e6SMatthew G. Knepley + f0 - boundary integrand for the test function term
1570194d53e6SMatthew G. Knepley - f1 - boundary integrand for the test function gradient term
1571194d53e6SMatthew G. Knepley 
1572194d53e6SMatthew G. Knepley   Note: We are using a first order FEM model for the weak form:
1573194d53e6SMatthew G. Knepley 
1574194d53e6SMatthew 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
1575194d53e6SMatthew G. Knepley 
1576194d53e6SMatthew G. Knepley The calling sequence for the callbacks f0 and f1 is given by:
1577194d53e6SMatthew G. Knepley 
157830b9ff8bSMatthew G. Knepley $ f0(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1579194d53e6SMatthew G. Knepley $    const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1580194d53e6SMatthew G. Knepley $    const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
158130b9ff8bSMatthew G. Knepley $    PetscReal t, const PetscReal x[], const PetscReal n[], PetscScalar f0[])
1582194d53e6SMatthew G. Knepley 
1583194d53e6SMatthew G. Knepley + dim - the spatial dimension
1584194d53e6SMatthew G. Knepley . Nf - the number of fields
1585194d53e6SMatthew G. Knepley . uOff - the offset into u[] and u_t[] for each field
1586194d53e6SMatthew G. Knepley . uOff_x - the offset into u_x[] for each field
1587194d53e6SMatthew G. Knepley . u - each field evaluated at the current point
1588194d53e6SMatthew G. Knepley . u_t - the time derivative of each field evaluated at the current point
1589194d53e6SMatthew G. Knepley . u_x - the gradient of each field evaluated at the current point
1590194d53e6SMatthew G. Knepley . aOff - the offset into a[] and a_t[] for each auxiliary field
1591194d53e6SMatthew G. Knepley . aOff_x - the offset into a_x[] for each auxiliary field
1592194d53e6SMatthew G. Knepley . a - each auxiliary field evaluated at the current point
1593194d53e6SMatthew G. Knepley . a_t - the time derivative of each auxiliary field evaluated at the current point
1594194d53e6SMatthew G. Knepley . a_x - the gradient of auxiliary each field evaluated at the current point
1595194d53e6SMatthew G. Knepley . t - current time
1596194d53e6SMatthew G. Knepley . x - coordinates of the current point
1597194d53e6SMatthew G. Knepley . n - unit normal at the current point
1598194d53e6SMatthew G. Knepley - f0 - output values at the current point
1599194d53e6SMatthew G. Knepley 
1600194d53e6SMatthew G. Knepley   Level: intermediate
1601194d53e6SMatthew G. Knepley 
1602194d53e6SMatthew G. Knepley .seealso: PetscDSSetBdResidual()
1603194d53e6SMatthew G. Knepley @*/
16042764a2aaSMatthew G. Knepley PetscErrorCode PetscDSGetBdResidual(PetscDS prob, PetscInt f,
160530b9ff8bSMatthew G. Knepley                                     void (**f0)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1606194d53e6SMatthew G. Knepley                                                 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1607194d53e6SMatthew G. Knepley                                                 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
160830b9ff8bSMatthew G. Knepley                                                 PetscReal t, const PetscReal x[], const PetscReal n[], PetscScalar f0[]),
160930b9ff8bSMatthew G. Knepley                                     void (**f1)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1610194d53e6SMatthew G. Knepley                                                 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1611194d53e6SMatthew G. Knepley                                                 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
161230b9ff8bSMatthew G. Knepley                                                 PetscReal t, const PetscReal x[], const PetscReal n[], PetscScalar f1[]))
16132764a2aaSMatthew G. Knepley {
16142764a2aaSMatthew G. Knepley   PetscFunctionBegin;
16152764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
16162764a2aaSMatthew 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);
16172764a2aaSMatthew G. Knepley   if (f0) {PetscValidPointer(f0, 3); *f0 = prob->fBd[f*2+0];}
16182764a2aaSMatthew G. Knepley   if (f1) {PetscValidPointer(f1, 4); *f1 = prob->fBd[f*2+1];}
16192764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
16202764a2aaSMatthew G. Knepley }
16212764a2aaSMatthew G. Knepley 
16222764a2aaSMatthew G. Knepley #undef __FUNCT__
16232764a2aaSMatthew G. Knepley #define __FUNCT__ "PetscDSSetBdResidual"
1624194d53e6SMatthew G. Knepley /*@C
1625194d53e6SMatthew G. Knepley   PetscDSSetBdResidual - Get the pointwise boundary residual function for a given test field
1626194d53e6SMatthew G. Knepley 
1627194d53e6SMatthew G. Knepley   Not collective
1628194d53e6SMatthew G. Knepley 
1629194d53e6SMatthew G. Knepley   Input Parameters:
1630194d53e6SMatthew G. Knepley + prob - The PetscDS
1631194d53e6SMatthew G. Knepley . f    - The test field number
1632194d53e6SMatthew G. Knepley . f0 - boundary integrand for the test function term
1633194d53e6SMatthew G. Knepley - f1 - boundary integrand for the test function gradient term
1634194d53e6SMatthew G. Knepley 
1635194d53e6SMatthew G. Knepley   Note: We are using a first order FEM model for the weak form:
1636194d53e6SMatthew G. Knepley 
1637194d53e6SMatthew 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
1638194d53e6SMatthew G. Knepley 
1639194d53e6SMatthew G. Knepley The calling sequence for the callbacks f0 and f1 is given by:
1640194d53e6SMatthew G. Knepley 
164130b9ff8bSMatthew G. Knepley $ f0(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1642194d53e6SMatthew G. Knepley $    const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1643194d53e6SMatthew G. Knepley $    const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
164430b9ff8bSMatthew G. Knepley $    PetscReal t, const PetscReal x[], const PetscReal n[], PetscScalar f0[])
1645194d53e6SMatthew G. Knepley 
1646194d53e6SMatthew G. Knepley + dim - the spatial dimension
1647194d53e6SMatthew G. Knepley . Nf - the number of fields
1648194d53e6SMatthew G. Knepley . uOff - the offset into u[] and u_t[] for each field
1649194d53e6SMatthew G. Knepley . uOff_x - the offset into u_x[] for each field
1650194d53e6SMatthew G. Knepley . u - each field evaluated at the current point
1651194d53e6SMatthew G. Knepley . u_t - the time derivative of each field evaluated at the current point
1652194d53e6SMatthew G. Knepley . u_x - the gradient of each field evaluated at the current point
1653194d53e6SMatthew G. Knepley . aOff - the offset into a[] and a_t[] for each auxiliary field
1654194d53e6SMatthew G. Knepley . aOff_x - the offset into a_x[] for each auxiliary field
1655194d53e6SMatthew G. Knepley . a - each auxiliary field evaluated at the current point
1656194d53e6SMatthew G. Knepley . a_t - the time derivative of each auxiliary field evaluated at the current point
1657194d53e6SMatthew G. Knepley . a_x - the gradient of auxiliary each field evaluated at the current point
1658194d53e6SMatthew G. Knepley . t - current time
1659194d53e6SMatthew G. Knepley . x - coordinates of the current point
1660194d53e6SMatthew G. Knepley . n - unit normal at the current point
1661194d53e6SMatthew G. Knepley - f0 - output values at the current point
1662194d53e6SMatthew G. Knepley 
1663194d53e6SMatthew G. Knepley   Level: intermediate
1664194d53e6SMatthew G. Knepley 
1665194d53e6SMatthew G. Knepley .seealso: PetscDSGetBdResidual()
1666194d53e6SMatthew G. Knepley @*/
16672764a2aaSMatthew G. Knepley PetscErrorCode PetscDSSetBdResidual(PetscDS prob, PetscInt f,
166830b9ff8bSMatthew G. Knepley                                     void (*f0)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1669194d53e6SMatthew G. Knepley                                                const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1670194d53e6SMatthew G. Knepley                                                const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
167130b9ff8bSMatthew G. Knepley                                                PetscReal t, const PetscReal x[], const PetscReal n[], PetscScalar f0[]),
167230b9ff8bSMatthew G. Knepley                                     void (*f1)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1673194d53e6SMatthew G. Knepley                                                const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1674194d53e6SMatthew G. Knepley                                                const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
167530b9ff8bSMatthew G. Knepley                                                PetscReal t, const PetscReal x[], const PetscReal n[], PetscScalar f1[]))
16762764a2aaSMatthew G. Knepley {
16772764a2aaSMatthew G. Knepley   PetscErrorCode ierr;
16782764a2aaSMatthew G. Knepley 
16792764a2aaSMatthew G. Knepley   PetscFunctionBegin;
16802764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
16812764a2aaSMatthew G. Knepley   if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f);
16822764a2aaSMatthew G. Knepley   ierr = PetscDSEnlarge_Static(prob, f+1);CHKERRQ(ierr);
16832764a2aaSMatthew G. Knepley   if (f0) {PetscValidFunction(f0, 3); prob->fBd[f*2+0] = f0;}
16842764a2aaSMatthew G. Knepley   if (f1) {PetscValidFunction(f1, 4); prob->fBd[f*2+1] = f1;}
16852764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
16862764a2aaSMatthew G. Knepley }
16872764a2aaSMatthew G. Knepley 
16882764a2aaSMatthew G. Knepley #undef __FUNCT__
16892764a2aaSMatthew G. Knepley #define __FUNCT__ "PetscDSGetBdJacobian"
1690194d53e6SMatthew G. Knepley /*@C
1691194d53e6SMatthew G. Knepley   PetscDSGetBdJacobian - Get the pointwise boundary Jacobian function for given test and basis field
1692194d53e6SMatthew G. Knepley 
1693194d53e6SMatthew G. Knepley   Not collective
1694194d53e6SMatthew G. Knepley 
1695194d53e6SMatthew G. Knepley   Input Parameters:
1696194d53e6SMatthew G. Knepley + prob - The PetscDS
1697194d53e6SMatthew G. Knepley . f    - The test field number
1698194d53e6SMatthew G. Knepley - g    - The field number
1699194d53e6SMatthew G. Knepley 
1700194d53e6SMatthew G. Knepley   Output Parameters:
1701194d53e6SMatthew G. Knepley + g0 - integrand for the test and basis function term
1702194d53e6SMatthew G. Knepley . g1 - integrand for the test function and basis function gradient term
1703194d53e6SMatthew G. Knepley . g2 - integrand for the test function gradient and basis function term
1704194d53e6SMatthew G. Knepley - g3 - integrand for the test function gradient and basis function gradient term
1705194d53e6SMatthew G. Knepley 
1706194d53e6SMatthew G. Knepley   Note: We are using a first order FEM model for the weak form:
1707194d53e6SMatthew G. Knepley 
1708194d53e6SMatthew 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
1709194d53e6SMatthew G. Knepley 
1710194d53e6SMatthew G. Knepley The calling sequence for the callbacks g0, g1, g2 and g3 is given by:
1711194d53e6SMatthew G. Knepley 
171230b9ff8bSMatthew G. Knepley $ g0(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1713194d53e6SMatthew G. Knepley $    const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1714194d53e6SMatthew G. Knepley $    const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
171530b9ff8bSMatthew G. Knepley $    PetscReal t, const PetscReal x[], const PetscReal n[], PetscScalar g0[])
1716194d53e6SMatthew G. Knepley 
1717194d53e6SMatthew G. Knepley + dim - the spatial dimension
1718194d53e6SMatthew G. Knepley . Nf - the number of fields
1719194d53e6SMatthew G. Knepley . uOff - the offset into u[] and u_t[] for each field
1720194d53e6SMatthew G. Knepley . uOff_x - the offset into u_x[] for each field
1721194d53e6SMatthew G. Knepley . u - each field evaluated at the current point
1722194d53e6SMatthew G. Knepley . u_t - the time derivative of each field evaluated at the current point
1723194d53e6SMatthew G. Knepley . u_x - the gradient of each field evaluated at the current point
1724194d53e6SMatthew G. Knepley . aOff - the offset into a[] and a_t[] for each auxiliary field
1725194d53e6SMatthew G. Knepley . aOff_x - the offset into a_x[] for each auxiliary field
1726194d53e6SMatthew G. Knepley . a - each auxiliary field evaluated at the current point
1727194d53e6SMatthew G. Knepley . a_t - the time derivative of each auxiliary field evaluated at the current point
1728194d53e6SMatthew G. Knepley . a_x - the gradient of auxiliary each field evaluated at the current point
1729194d53e6SMatthew G. Knepley . t - current time
17302aa1fc23SMatthew G. Knepley . u_tShift - the multiplier a for dF/dU_t
1731194d53e6SMatthew G. Knepley . x - coordinates of the current point
1732194d53e6SMatthew G. Knepley . n - normal at the current point
1733194d53e6SMatthew G. Knepley - g0 - output values at the current point
1734194d53e6SMatthew G. Knepley 
1735194d53e6SMatthew G. Knepley   Level: intermediate
1736194d53e6SMatthew G. Knepley 
1737194d53e6SMatthew G. Knepley .seealso: PetscDSSetBdJacobian()
1738194d53e6SMatthew G. Knepley @*/
17392764a2aaSMatthew G. Knepley PetscErrorCode PetscDSGetBdJacobian(PetscDS prob, PetscInt f, PetscInt g,
174030b9ff8bSMatthew G. Knepley                                     void (**g0)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1741194d53e6SMatthew G. Knepley                                                 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1742194d53e6SMatthew G. Knepley                                                 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
17432aa1fc23SMatthew G. Knepley                                                 PetscReal t, PetscReal u_tShift, const PetscReal x[], const PetscReal n[], PetscScalar g0[]),
174430b9ff8bSMatthew G. Knepley                                     void (**g1)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1745194d53e6SMatthew G. Knepley                                                 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1746194d53e6SMatthew G. Knepley                                                 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
17472aa1fc23SMatthew G. Knepley                                                 PetscReal t, PetscReal u_tShift, const PetscReal x[], const PetscReal n[], PetscScalar g1[]),
174830b9ff8bSMatthew G. Knepley                                     void (**g2)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1749194d53e6SMatthew G. Knepley                                                 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1750194d53e6SMatthew G. Knepley                                                 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
17512aa1fc23SMatthew G. Knepley                                                 PetscReal t, PetscReal u_tShift, const PetscReal x[], const PetscReal n[], PetscScalar g2[]),
175230b9ff8bSMatthew G. Knepley                                     void (**g3)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1753194d53e6SMatthew G. Knepley                                                 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1754194d53e6SMatthew G. Knepley                                                 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
17552aa1fc23SMatthew G. Knepley                                                 PetscReal t, PetscReal u_tShift, const PetscReal x[], const PetscReal n[], PetscScalar g3[]))
17562764a2aaSMatthew G. Knepley {
17572764a2aaSMatthew G. Knepley   PetscFunctionBegin;
17582764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
17592764a2aaSMatthew 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);
17602764a2aaSMatthew 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);
17612764a2aaSMatthew G. Knepley   if (g0) {PetscValidPointer(g0, 4); *g0 = prob->gBd[(f*prob->Nf + g)*4+0];}
17622764a2aaSMatthew G. Knepley   if (g1) {PetscValidPointer(g1, 5); *g1 = prob->gBd[(f*prob->Nf + g)*4+1];}
17632764a2aaSMatthew G. Knepley   if (g2) {PetscValidPointer(g2, 6); *g2 = prob->gBd[(f*prob->Nf + g)*4+2];}
17642764a2aaSMatthew G. Knepley   if (g3) {PetscValidPointer(g3, 7); *g3 = prob->gBd[(f*prob->Nf + g)*4+3];}
17652764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
17662764a2aaSMatthew G. Knepley }
17672764a2aaSMatthew G. Knepley 
17682764a2aaSMatthew G. Knepley #undef __FUNCT__
17692764a2aaSMatthew G. Knepley #define __FUNCT__ "PetscDSSetBdJacobian"
1770194d53e6SMatthew G. Knepley /*@C
1771194d53e6SMatthew G. Knepley   PetscDSSetBdJacobian - Set the pointwise boundary Jacobian function for given test and basis field
1772194d53e6SMatthew G. Knepley 
1773194d53e6SMatthew G. Knepley   Not collective
1774194d53e6SMatthew G. Knepley 
1775194d53e6SMatthew G. Knepley   Input Parameters:
1776194d53e6SMatthew G. Knepley + prob - The PetscDS
1777194d53e6SMatthew G. Knepley . f    - The test field number
1778194d53e6SMatthew G. Knepley . g    - The field number
1779194d53e6SMatthew G. Knepley . g0 - integrand for the test and basis function term
1780194d53e6SMatthew G. Knepley . g1 - integrand for the test function and basis function gradient term
1781194d53e6SMatthew G. Knepley . g2 - integrand for the test function gradient and basis function term
1782194d53e6SMatthew G. Knepley - g3 - integrand for the test function gradient and basis function gradient term
1783194d53e6SMatthew G. Knepley 
1784194d53e6SMatthew G. Knepley   Note: We are using a first order FEM model for the weak form:
1785194d53e6SMatthew G. Knepley 
1786194d53e6SMatthew 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
1787194d53e6SMatthew G. Knepley 
1788194d53e6SMatthew G. Knepley The calling sequence for the callbacks g0, g1, g2 and g3 is given by:
1789194d53e6SMatthew G. Knepley 
179030b9ff8bSMatthew G. Knepley $ g0(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1791194d53e6SMatthew G. Knepley $    const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1792194d53e6SMatthew G. Knepley $    const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
179330b9ff8bSMatthew G. Knepley $    PetscReal t, const PetscReal x[], const PetscReal n[], PetscScalar g0[])
1794194d53e6SMatthew G. Knepley 
1795194d53e6SMatthew G. Knepley + dim - the spatial dimension
1796194d53e6SMatthew G. Knepley . Nf - the number of fields
1797194d53e6SMatthew G. Knepley . uOff - the offset into u[] and u_t[] for each field
1798194d53e6SMatthew G. Knepley . uOff_x - the offset into u_x[] for each field
1799194d53e6SMatthew G. Knepley . u - each field evaluated at the current point
1800194d53e6SMatthew G. Knepley . u_t - the time derivative of each field evaluated at the current point
1801194d53e6SMatthew G. Knepley . u_x - the gradient of each field evaluated at the current point
1802194d53e6SMatthew G. Knepley . aOff - the offset into a[] and a_t[] for each auxiliary field
1803194d53e6SMatthew G. Knepley . aOff_x - the offset into a_x[] for each auxiliary field
1804194d53e6SMatthew G. Knepley . a - each auxiliary field evaluated at the current point
1805194d53e6SMatthew G. Knepley . a_t - the time derivative of each auxiliary field evaluated at the current point
1806194d53e6SMatthew G. Knepley . a_x - the gradient of auxiliary each field evaluated at the current point
1807194d53e6SMatthew G. Knepley . t - current time
18082aa1fc23SMatthew G. Knepley . u_tShift - the multiplier a for dF/dU_t
1809194d53e6SMatthew G. Knepley . x - coordinates of the current point
1810194d53e6SMatthew G. Knepley . n - normal at the current point
1811194d53e6SMatthew G. Knepley - g0 - output values at the current point
1812194d53e6SMatthew G. Knepley 
1813194d53e6SMatthew G. Knepley   Level: intermediate
1814194d53e6SMatthew G. Knepley 
1815194d53e6SMatthew G. Knepley .seealso: PetscDSGetBdJacobian()
1816194d53e6SMatthew G. Knepley @*/
18172764a2aaSMatthew G. Knepley PetscErrorCode PetscDSSetBdJacobian(PetscDS prob, PetscInt f, PetscInt g,
181830b9ff8bSMatthew G. Knepley                                     void (*g0)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1819194d53e6SMatthew G. Knepley                                                const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1820194d53e6SMatthew G. Knepley                                                const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
18212aa1fc23SMatthew G. Knepley                                                PetscReal t, PetscReal u_tShift, const PetscReal x[], const PetscReal n[], PetscScalar g0[]),
182230b9ff8bSMatthew G. Knepley                                     void (*g1)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1823194d53e6SMatthew G. Knepley                                                const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1824194d53e6SMatthew G. Knepley                                                const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
18252aa1fc23SMatthew G. Knepley                                                PetscReal t, PetscReal u_tShift, const PetscReal x[], const PetscReal n[], PetscScalar g1[]),
182630b9ff8bSMatthew G. Knepley                                     void (*g2)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1827194d53e6SMatthew G. Knepley                                                const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1828194d53e6SMatthew G. Knepley                                                const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
18292aa1fc23SMatthew G. Knepley                                                PetscReal t, PetscReal u_tShift, const PetscReal x[], const PetscReal n[], PetscScalar g2[]),
183030b9ff8bSMatthew G. Knepley                                     void (*g3)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1831194d53e6SMatthew G. Knepley                                                const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1832194d53e6SMatthew G. Knepley                                                const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
18332aa1fc23SMatthew G. Knepley                                                PetscReal t, PetscReal u_tShift, const PetscReal x[], const PetscReal n[], PetscScalar g3[]))
18342764a2aaSMatthew G. Knepley {
18352764a2aaSMatthew G. Knepley   PetscErrorCode ierr;
18362764a2aaSMatthew G. Knepley 
18372764a2aaSMatthew G. Knepley   PetscFunctionBegin;
18382764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
18392764a2aaSMatthew G. Knepley   if (g0) PetscValidFunction(g0, 4);
18402764a2aaSMatthew G. Knepley   if (g1) PetscValidFunction(g1, 5);
18412764a2aaSMatthew G. Knepley   if (g2) PetscValidFunction(g2, 6);
18422764a2aaSMatthew G. Knepley   if (g3) PetscValidFunction(g3, 7);
18432764a2aaSMatthew G. Knepley   if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f);
18442764a2aaSMatthew G. Knepley   if (g < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", g);
18452764a2aaSMatthew G. Knepley   ierr = PetscDSEnlarge_Static(prob, PetscMax(f, g)+1);CHKERRQ(ierr);
18462764a2aaSMatthew G. Knepley   prob->gBd[(f*prob->Nf + g)*4+0] = g0;
18472764a2aaSMatthew G. Knepley   prob->gBd[(f*prob->Nf + g)*4+1] = g1;
18482764a2aaSMatthew G. Knepley   prob->gBd[(f*prob->Nf + g)*4+2] = g2;
18492764a2aaSMatthew G. Knepley   prob->gBd[(f*prob->Nf + g)*4+3] = g3;
18502764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
18512764a2aaSMatthew G. Knepley }
18522764a2aaSMatthew G. Knepley 
18532764a2aaSMatthew G. Knepley #undef __FUNCT__
18542764a2aaSMatthew G. Knepley #define __FUNCT__ "PetscDSGetFieldOffset"
1855bc4ae4beSMatthew G. Knepley /*@
1856bc4ae4beSMatthew G. Knepley   PetscDSGetFieldOffset - Returns the offset of the given field in the full space basis
1857bc4ae4beSMatthew G. Knepley 
1858bc4ae4beSMatthew G. Knepley   Not collective
1859bc4ae4beSMatthew G. Knepley 
1860bc4ae4beSMatthew G. Knepley   Input Parameters:
1861bc4ae4beSMatthew G. Knepley + prob - The PetscDS object
1862bc4ae4beSMatthew G. Knepley - f - The field number
1863bc4ae4beSMatthew G. Knepley 
1864bc4ae4beSMatthew G. Knepley   Output Parameter:
1865bc4ae4beSMatthew G. Knepley . off - The offset
1866bc4ae4beSMatthew G. Knepley 
1867bc4ae4beSMatthew G. Knepley   Level: beginner
1868bc4ae4beSMatthew G. Knepley 
1869bc4ae4beSMatthew G. Knepley .seealso: PetscDSGetBdFieldOffset(), PetscDSGetNumFields(), PetscDSCreate()
1870bc4ae4beSMatthew G. Knepley @*/
18712764a2aaSMatthew G. Knepley PetscErrorCode PetscDSGetFieldOffset(PetscDS prob, PetscInt f, PetscInt *off)
18722764a2aaSMatthew G. Knepley {
18732764a2aaSMatthew G. Knepley   PetscInt       g;
18742764a2aaSMatthew G. Knepley   PetscErrorCode ierr;
18752764a2aaSMatthew G. Knepley 
18762764a2aaSMatthew G. Knepley   PetscFunctionBegin;
18772764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
18782764a2aaSMatthew G. Knepley   PetscValidPointer(off, 3);
18792764a2aaSMatthew 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);
18802764a2aaSMatthew G. Knepley   *off = 0;
18812764a2aaSMatthew G. Knepley   for (g = 0; g < f; ++g) {
18822764a2aaSMatthew G. Knepley     PetscFE  fe = (PetscFE) prob->disc[g];
18832764a2aaSMatthew G. Knepley     PetscInt Nb, Nc;
18842764a2aaSMatthew G. Knepley 
18852764a2aaSMatthew G. Knepley     ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr);
18862764a2aaSMatthew G. Knepley     ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
18872764a2aaSMatthew G. Knepley     *off += Nb*Nc;
18882764a2aaSMatthew G. Knepley   }
18892764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
18902764a2aaSMatthew G. Knepley }
18912764a2aaSMatthew G. Knepley 
18922764a2aaSMatthew G. Knepley #undef __FUNCT__
18932764a2aaSMatthew G. Knepley #define __FUNCT__ "PetscDSGetBdFieldOffset"
1894bc4ae4beSMatthew G. Knepley /*@
1895c3ac4435SMatthew G. Knepley   PetscDSGetBdFieldOffset - Returns the offset of the given field in the full space boundary basis
1896bc4ae4beSMatthew G. Knepley 
1897bc4ae4beSMatthew G. Knepley   Not collective
1898bc4ae4beSMatthew G. Knepley 
1899bc4ae4beSMatthew G. Knepley   Input Parameters:
1900bc4ae4beSMatthew G. Knepley + prob - The PetscDS object
1901bc4ae4beSMatthew G. Knepley - f - The field number
1902bc4ae4beSMatthew G. Knepley 
1903bc4ae4beSMatthew G. Knepley   Output Parameter:
1904bc4ae4beSMatthew G. Knepley . off - The boundary offset
1905bc4ae4beSMatthew G. Knepley 
1906bc4ae4beSMatthew G. Knepley   Level: beginner
1907bc4ae4beSMatthew G. Knepley 
1908bc4ae4beSMatthew G. Knepley .seealso: PetscDSGetFieldOffset(), PetscDSGetNumFields(), PetscDSCreate()
1909bc4ae4beSMatthew G. Knepley @*/
19102764a2aaSMatthew G. Knepley PetscErrorCode PetscDSGetBdFieldOffset(PetscDS prob, PetscInt f, PetscInt *off)
19112764a2aaSMatthew G. Knepley {
19122764a2aaSMatthew G. Knepley   PetscInt       g;
19132764a2aaSMatthew G. Knepley   PetscErrorCode ierr;
19142764a2aaSMatthew G. Knepley 
19152764a2aaSMatthew G. Knepley   PetscFunctionBegin;
19162764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
19172764a2aaSMatthew G. Knepley   PetscValidPointer(off, 3);
19182764a2aaSMatthew 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);
19192764a2aaSMatthew G. Knepley   *off = 0;
19202764a2aaSMatthew G. Knepley   for (g = 0; g < f; ++g) {
19212764a2aaSMatthew G. Knepley     PetscFE  fe = (PetscFE) prob->discBd[g];
19222764a2aaSMatthew G. Knepley     PetscInt Nb, Nc;
19232764a2aaSMatthew G. Knepley 
19242764a2aaSMatthew G. Knepley     ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr);
19252764a2aaSMatthew G. Knepley     ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
19262764a2aaSMatthew G. Knepley     *off += Nb*Nc;
19276ce16762SMatthew G. Knepley   }
19286ce16762SMatthew G. Knepley   PetscFunctionReturn(0);
19296ce16762SMatthew G. Knepley }
19306ce16762SMatthew G. Knepley 
19316ce16762SMatthew G. Knepley #undef __FUNCT__
19326ce16762SMatthew G. Knepley #define __FUNCT__ "PetscDSGetComponentOffset"
19336ce16762SMatthew G. Knepley /*@
19346ce16762SMatthew G. Knepley   PetscDSGetComponentOffset - Returns the offset of the given field on an evaluation point
19356ce16762SMatthew G. Knepley 
19366ce16762SMatthew G. Knepley   Not collective
19376ce16762SMatthew G. Knepley 
19386ce16762SMatthew G. Knepley   Input Parameters:
19396ce16762SMatthew G. Knepley + prob - The PetscDS object
19406ce16762SMatthew G. Knepley - f - The field number
19416ce16762SMatthew G. Knepley 
19426ce16762SMatthew G. Knepley   Output Parameter:
19436ce16762SMatthew G. Knepley . off - The offset
19446ce16762SMatthew G. Knepley 
19456ce16762SMatthew G. Knepley   Level: beginner
19466ce16762SMatthew G. Knepley 
19476ce16762SMatthew G. Knepley .seealso: PetscDSGetBdFieldOffset(), PetscDSGetNumFields(), PetscDSCreate()
19486ce16762SMatthew G. Knepley @*/
19496ce16762SMatthew G. Knepley PetscErrorCode PetscDSGetComponentOffset(PetscDS prob, PetscInt f, PetscInt *off)
19506ce16762SMatthew G. Knepley {
19516ce16762SMatthew G. Knepley   PetscInt       g;
19526ce16762SMatthew G. Knepley   PetscErrorCode ierr;
19536ce16762SMatthew G. Knepley 
19546ce16762SMatthew G. Knepley   PetscFunctionBegin;
19556ce16762SMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
19566ce16762SMatthew G. Knepley   PetscValidPointer(off, 3);
19576ce16762SMatthew 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);
19586ce16762SMatthew G. Knepley   *off = 0;
19596ce16762SMatthew G. Knepley   for (g = 0; g < f; ++g) {
19606ce16762SMatthew G. Knepley     PetscFE  fe = (PetscFE) prob->disc[g];
19616ce16762SMatthew G. Knepley     PetscInt Nc;
19626ce16762SMatthew G. Knepley 
19636ce16762SMatthew G. Knepley     ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
19646ce16762SMatthew G. Knepley     *off += Nc;
19652764a2aaSMatthew G. Knepley   }
19662764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
19672764a2aaSMatthew G. Knepley }
19682764a2aaSMatthew G. Knepley 
19692764a2aaSMatthew G. Knepley #undef __FUNCT__
1970194d53e6SMatthew G. Knepley #define __FUNCT__ "PetscDSGetComponentOffsets"
1971194d53e6SMatthew G. Knepley /*@
1972194d53e6SMatthew G. Knepley   PetscDSGetComponentOffsets - Returns the offset of each field on an evaluation point
1973194d53e6SMatthew G. Knepley 
1974194d53e6SMatthew G. Knepley   Not collective
1975194d53e6SMatthew G. Knepley 
1976194d53e6SMatthew G. Knepley   Input Parameter:
1977194d53e6SMatthew G. Knepley . prob - The PetscDS object
1978194d53e6SMatthew G. Knepley 
1979194d53e6SMatthew G. Knepley   Output Parameter:
1980194d53e6SMatthew G. Knepley . offsets - The offsets
1981194d53e6SMatthew G. Knepley 
1982194d53e6SMatthew G. Knepley   Level: beginner
1983194d53e6SMatthew G. Knepley 
1984194d53e6SMatthew G. Knepley .seealso: PetscDSGetBdFieldOffset(), PetscDSGetNumFields(), PetscDSCreate()
1985194d53e6SMatthew G. Knepley @*/
1986194d53e6SMatthew G. Knepley PetscErrorCode PetscDSGetComponentOffsets(PetscDS prob, PetscInt *offsets[])
1987194d53e6SMatthew G. Knepley {
1988194d53e6SMatthew G. Knepley   PetscFunctionBegin;
1989194d53e6SMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1990194d53e6SMatthew G. Knepley   PetscValidPointer(offsets, 2);
1991194d53e6SMatthew G. Knepley   *offsets = prob->off;
1992194d53e6SMatthew G. Knepley   PetscFunctionReturn(0);
1993194d53e6SMatthew G. Knepley }
1994194d53e6SMatthew G. Knepley 
1995194d53e6SMatthew G. Knepley #undef __FUNCT__
1996194d53e6SMatthew G. Knepley #define __FUNCT__ "PetscDSGetComponentDerivativeOffsets"
1997194d53e6SMatthew G. Knepley /*@
1998194d53e6SMatthew G. Knepley   PetscDSGetComponentDerivativeOffsets - Returns the offset of each field derivative on an evaluation point
1999194d53e6SMatthew G. Knepley 
2000194d53e6SMatthew G. Knepley   Not collective
2001194d53e6SMatthew G. Knepley 
2002194d53e6SMatthew G. Knepley   Input Parameter:
2003194d53e6SMatthew G. Knepley . prob - The PetscDS object
2004194d53e6SMatthew G. Knepley 
2005194d53e6SMatthew G. Knepley   Output Parameter:
2006194d53e6SMatthew G. Knepley . offsets - The offsets
2007194d53e6SMatthew G. Knepley 
2008194d53e6SMatthew G. Knepley   Level: beginner
2009194d53e6SMatthew G. Knepley 
2010194d53e6SMatthew G. Knepley .seealso: PetscDSGetBdFieldOffset(), PetscDSGetNumFields(), PetscDSCreate()
2011194d53e6SMatthew G. Knepley @*/
2012194d53e6SMatthew G. Knepley PetscErrorCode PetscDSGetComponentDerivativeOffsets(PetscDS prob, PetscInt *offsets[])
2013194d53e6SMatthew G. Knepley {
2014194d53e6SMatthew G. Knepley   PetscFunctionBegin;
2015194d53e6SMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
2016194d53e6SMatthew G. Knepley   PetscValidPointer(offsets, 2);
2017194d53e6SMatthew G. Knepley   *offsets = prob->offDer;
2018194d53e6SMatthew G. Knepley   PetscFunctionReturn(0);
2019194d53e6SMatthew G. Knepley }
2020194d53e6SMatthew G. Knepley 
2021194d53e6SMatthew G. Knepley #undef __FUNCT__
2022194d53e6SMatthew G. Knepley #define __FUNCT__ "PetscDSGetComponentBdOffsets"
2023194d53e6SMatthew G. Knepley /*@
2024194d53e6SMatthew G. Knepley   PetscDSGetComponentBdOffsets - Returns the offset of each field on a boundary evaluation point
2025194d53e6SMatthew G. Knepley 
2026194d53e6SMatthew G. Knepley   Not collective
2027194d53e6SMatthew G. Knepley 
2028194d53e6SMatthew G. Knepley   Input Parameter:
2029194d53e6SMatthew G. Knepley . prob - The PetscDS object
2030194d53e6SMatthew G. Knepley 
2031194d53e6SMatthew G. Knepley   Output Parameter:
2032194d53e6SMatthew G. Knepley . offsets - The offsets
2033194d53e6SMatthew G. Knepley 
2034194d53e6SMatthew G. Knepley   Level: beginner
2035194d53e6SMatthew G. Knepley 
2036194d53e6SMatthew G. Knepley .seealso: PetscDSGetBdFieldOffset(), PetscDSGetNumFields(), PetscDSCreate()
2037194d53e6SMatthew G. Knepley @*/
2038194d53e6SMatthew G. Knepley PetscErrorCode PetscDSGetComponentBdOffsets(PetscDS prob, PetscInt *offsets[])
2039194d53e6SMatthew G. Knepley {
2040194d53e6SMatthew G. Knepley   PetscFunctionBegin;
2041194d53e6SMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
2042194d53e6SMatthew G. Knepley   PetscValidPointer(offsets, 2);
2043194d53e6SMatthew G. Knepley   *offsets = prob->offBd;
2044194d53e6SMatthew G. Knepley   PetscFunctionReturn(0);
2045194d53e6SMatthew G. Knepley }
2046194d53e6SMatthew G. Knepley 
2047194d53e6SMatthew G. Knepley #undef __FUNCT__
2048194d53e6SMatthew G. Knepley #define __FUNCT__ "PetscDSGetComponentBdDerivativeOffsets"
2049194d53e6SMatthew G. Knepley /*@
2050194d53e6SMatthew G. Knepley   PetscDSGetComponentBdDerivativeOffsets - Returns the offset of each field derivative on a boundary evaluation point
2051194d53e6SMatthew G. Knepley 
2052194d53e6SMatthew G. Knepley   Not collective
2053194d53e6SMatthew G. Knepley 
2054194d53e6SMatthew G. Knepley   Input Parameter:
2055194d53e6SMatthew G. Knepley . prob - The PetscDS object
2056194d53e6SMatthew G. Knepley 
2057194d53e6SMatthew G. Knepley   Output Parameter:
2058194d53e6SMatthew G. Knepley . offsets - The offsets
2059194d53e6SMatthew G. Knepley 
2060194d53e6SMatthew G. Knepley   Level: beginner
2061194d53e6SMatthew G. Knepley 
2062194d53e6SMatthew G. Knepley .seealso: PetscDSGetBdFieldOffset(), PetscDSGetNumFields(), PetscDSCreate()
2063194d53e6SMatthew G. Knepley @*/
2064194d53e6SMatthew G. Knepley PetscErrorCode PetscDSGetComponentBdDerivativeOffsets(PetscDS prob, PetscInt *offsets[])
2065194d53e6SMatthew G. Knepley {
2066194d53e6SMatthew G. Knepley   PetscFunctionBegin;
2067194d53e6SMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
2068194d53e6SMatthew G. Knepley   PetscValidPointer(offsets, 2);
2069194d53e6SMatthew G. Knepley   *offsets = prob->offDerBd;
2070194d53e6SMatthew G. Knepley   PetscFunctionReturn(0);
2071194d53e6SMatthew G. Knepley }
2072194d53e6SMatthew G. Knepley 
2073194d53e6SMatthew G. Knepley #undef __FUNCT__
20742764a2aaSMatthew G. Knepley #define __FUNCT__ "PetscDSGetTabulation"
207568c9edb9SMatthew G. Knepley /*@C
207668c9edb9SMatthew G. Knepley   PetscDSGetTabulation - Return the basis tabulation at quadrature points for the volume discretization
207768c9edb9SMatthew G. Knepley 
207868c9edb9SMatthew G. Knepley   Not collective
207968c9edb9SMatthew G. Knepley 
208068c9edb9SMatthew G. Knepley   Input Parameter:
208168c9edb9SMatthew G. Knepley . prob - The PetscDS object
208268c9edb9SMatthew G. Knepley 
208368c9edb9SMatthew G. Knepley   Output Parameters:
208468c9edb9SMatthew G. Knepley + basis - The basis function tabulation at quadrature points
208568c9edb9SMatthew G. Knepley - basisDer - The basis function derivative tabulation at quadrature points
208668c9edb9SMatthew G. Knepley 
208768c9edb9SMatthew G. Knepley   Level: intermediate
208868c9edb9SMatthew G. Knepley 
208968c9edb9SMatthew G. Knepley .seealso: PetscDSGetBdTabulation(), PetscDSCreate()
209068c9edb9SMatthew G. Knepley @*/
20912764a2aaSMatthew G. Knepley PetscErrorCode PetscDSGetTabulation(PetscDS prob, PetscReal ***basis, PetscReal ***basisDer)
20922764a2aaSMatthew G. Knepley {
20932764a2aaSMatthew G. Knepley   PetscErrorCode ierr;
20942764a2aaSMatthew G. Knepley 
20952764a2aaSMatthew G. Knepley   PetscFunctionBegin;
20962764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
20972764a2aaSMatthew G. Knepley   ierr = PetscDSSetUp(prob);CHKERRQ(ierr);
20982764a2aaSMatthew G. Knepley   if (basis)    {PetscValidPointer(basis, 2);    *basis    = prob->basis;}
20992764a2aaSMatthew G. Knepley   if (basisDer) {PetscValidPointer(basisDer, 3); *basisDer = prob->basisDer;}
21002764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
21012764a2aaSMatthew G. Knepley }
21022764a2aaSMatthew G. Knepley 
21032764a2aaSMatthew G. Knepley #undef __FUNCT__
21042764a2aaSMatthew G. Knepley #define __FUNCT__ "PetscDSGetBdTabulation"
210568c9edb9SMatthew G. Knepley /*@C
210668c9edb9SMatthew G. Knepley   PetscDSGetBdTabulation - Return the basis tabulation at quadrature points for the boundary discretization
210768c9edb9SMatthew G. Knepley 
210868c9edb9SMatthew G. Knepley   Not collective
210968c9edb9SMatthew G. Knepley 
211068c9edb9SMatthew G. Knepley   Input Parameter:
211168c9edb9SMatthew G. Knepley . prob - The PetscDS object
211268c9edb9SMatthew G. Knepley 
211368c9edb9SMatthew G. Knepley   Output Parameters:
211468c9edb9SMatthew G. Knepley + basis - The basis function tabulation at quadrature points
211568c9edb9SMatthew G. Knepley - basisDer - The basis function derivative tabulation at quadrature points
211668c9edb9SMatthew G. Knepley 
211768c9edb9SMatthew G. Knepley   Level: intermediate
211868c9edb9SMatthew G. Knepley 
211968c9edb9SMatthew G. Knepley .seealso: PetscDSGetTabulation(), PetscDSCreate()
212068c9edb9SMatthew G. Knepley @*/
21212764a2aaSMatthew G. Knepley PetscErrorCode PetscDSGetBdTabulation(PetscDS prob, PetscReal ***basis, PetscReal ***basisDer)
21222764a2aaSMatthew G. Knepley {
21232764a2aaSMatthew G. Knepley   PetscErrorCode ierr;
21242764a2aaSMatthew G. Knepley 
21252764a2aaSMatthew G. Knepley   PetscFunctionBegin;
21262764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
21272764a2aaSMatthew G. Knepley   ierr = PetscDSSetUp(prob);CHKERRQ(ierr);
21282764a2aaSMatthew G. Knepley   if (basis)    {PetscValidPointer(basis, 2);    *basis    = prob->basisBd;}
21292764a2aaSMatthew G. Knepley   if (basisDer) {PetscValidPointer(basisDer, 3); *basisDer = prob->basisDerBd;}
21302764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
21312764a2aaSMatthew G. Knepley }
21322764a2aaSMatthew G. Knepley 
21332764a2aaSMatthew G. Knepley #undef __FUNCT__
21342764a2aaSMatthew G. Knepley #define __FUNCT__ "PetscDSGetEvaluationArrays"
21352764a2aaSMatthew G. Knepley PetscErrorCode PetscDSGetEvaluationArrays(PetscDS prob, PetscScalar **u, PetscScalar **u_t, PetscScalar **u_x)
21362764a2aaSMatthew G. Knepley {
21372764a2aaSMatthew G. Knepley   PetscErrorCode ierr;
21382764a2aaSMatthew G. Knepley 
21392764a2aaSMatthew G. Knepley   PetscFunctionBegin;
21402764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
21412764a2aaSMatthew G. Knepley   ierr = PetscDSSetUp(prob);CHKERRQ(ierr);
21422764a2aaSMatthew G. Knepley   if (u)   {PetscValidPointer(u, 2);   *u   = prob->u;}
21432764a2aaSMatthew G. Knepley   if (u_t) {PetscValidPointer(u_t, 3); *u_t = prob->u_t;}
21442764a2aaSMatthew G. Knepley   if (u_x) {PetscValidPointer(u_x, 4); *u_x = prob->u_x;}
21452764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
21462764a2aaSMatthew G. Knepley }
21472764a2aaSMatthew G. Knepley 
21482764a2aaSMatthew G. Knepley #undef __FUNCT__
21492764a2aaSMatthew G. Knepley #define __FUNCT__ "PetscDSGetWeakFormArrays"
21502764a2aaSMatthew G. Knepley PetscErrorCode PetscDSGetWeakFormArrays(PetscDS prob, PetscScalar **f0, PetscScalar **f1, PetscScalar **g0, PetscScalar **g1, PetscScalar **g2, PetscScalar **g3)
21512764a2aaSMatthew G. Knepley {
21522764a2aaSMatthew G. Knepley   PetscErrorCode ierr;
21532764a2aaSMatthew G. Knepley 
21542764a2aaSMatthew G. Knepley   PetscFunctionBegin;
21552764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
21562764a2aaSMatthew G. Knepley   ierr = PetscDSSetUp(prob);CHKERRQ(ierr);
21572764a2aaSMatthew G. Knepley   if (f0) {PetscValidPointer(f0, 2); *f0 = prob->f0;}
21582764a2aaSMatthew G. Knepley   if (f1) {PetscValidPointer(f1, 3); *f1 = prob->f1;}
21592764a2aaSMatthew G. Knepley   if (g0) {PetscValidPointer(g0, 4); *g0 = prob->g0;}
21602764a2aaSMatthew G. Knepley   if (g1) {PetscValidPointer(g1, 5); *g1 = prob->g1;}
21612764a2aaSMatthew G. Knepley   if (g2) {PetscValidPointer(g2, 6); *g2 = prob->g2;}
21622764a2aaSMatthew G. Knepley   if (g3) {PetscValidPointer(g3, 7); *g3 = prob->g3;}
21632764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
21642764a2aaSMatthew G. Knepley }
21652764a2aaSMatthew G. Knepley 
21662764a2aaSMatthew G. Knepley #undef __FUNCT__
21672764a2aaSMatthew G. Knepley #define __FUNCT__ "PetscDSGetRefCoordArrays"
21682764a2aaSMatthew G. Knepley PetscErrorCode PetscDSGetRefCoordArrays(PetscDS prob, PetscReal **x, PetscScalar **refSpaceDer)
21692764a2aaSMatthew G. Knepley {
21702764a2aaSMatthew G. Knepley   PetscErrorCode ierr;
21712764a2aaSMatthew G. Knepley 
21722764a2aaSMatthew G. Knepley   PetscFunctionBegin;
21732764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
21742764a2aaSMatthew G. Knepley   ierr = PetscDSSetUp(prob);CHKERRQ(ierr);
21752764a2aaSMatthew G. Knepley   if (x)           {PetscValidPointer(x, 2);           *x           = prob->x;}
21762764a2aaSMatthew G. Knepley   if (refSpaceDer) {PetscValidPointer(refSpaceDer, 3); *refSpaceDer = prob->refSpaceDer;}
21772764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
21782764a2aaSMatthew G. Knepley }
21792764a2aaSMatthew G. Knepley 
21802764a2aaSMatthew G. Knepley #undef __FUNCT__
2181da51fcedSMatthew G. Knepley #define __FUNCT__ "PetscDSCopyEquations"
2182da51fcedSMatthew G. Knepley /*@
2183da51fcedSMatthew G. Knepley   PetscDSCopyEquations - Copy all pointwise function pointers to the new problem
2184da51fcedSMatthew G. Knepley 
2185da51fcedSMatthew G. Knepley   Not collective
2186da51fcedSMatthew G. Knepley 
2187da51fcedSMatthew G. Knepley   Input Parameter:
2188da51fcedSMatthew G. Knepley . prob - The PetscDS object
2189da51fcedSMatthew G. Knepley 
2190da51fcedSMatthew G. Knepley   Output Parameter:
2191da51fcedSMatthew G. Knepley . newprob - The PetscDS copy
2192da51fcedSMatthew G. Knepley 
2193da51fcedSMatthew G. Knepley   Level: intermediate
2194da51fcedSMatthew G. Knepley 
2195da51fcedSMatthew G. Knepley .seealso: PetscDSSetResidual(), PetscDSSetJacobian(), PetscDSSetRiemannSolver(), PetscDSSetBdResidual(), PetscDSSetBdJacobian(), PetscDSCreate()
2196da51fcedSMatthew G. Knepley @*/
2197da51fcedSMatthew G. Knepley PetscErrorCode PetscDSCopyEquations(PetscDS prob, PetscDS newprob)
2198da51fcedSMatthew G. Knepley {
2199da51fcedSMatthew G. Knepley   PetscInt       Nf, Ng, f, g;
2200da51fcedSMatthew G. Knepley   PetscErrorCode ierr;
2201da51fcedSMatthew G. Knepley 
2202da51fcedSMatthew G. Knepley   PetscFunctionBegin;
2203da51fcedSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
2204da51fcedSMatthew G. Knepley   PetscValidHeaderSpecific(newprob, PETSCDS_CLASSID, 2);
2205da51fcedSMatthew G. Knepley   ierr = PetscDSGetNumFields(prob, &Nf);CHKERRQ(ierr);
2206da51fcedSMatthew G. Knepley   ierr = PetscDSGetNumFields(newprob, &Ng);CHKERRQ(ierr);
2207da51fcedSMatthew G. Knepley   if (Nf != Ng) SETERRQ2(PetscObjectComm((PetscObject) prob), PETSC_ERR_ARG_SIZ, "Number of fields must match %D != %D", Nf, Ng);CHKERRQ(ierr);
2208da51fcedSMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
2209da51fcedSMatthew G. Knepley     PetscPointFunc   obj;
2210da51fcedSMatthew G. Knepley     PetscPointFunc   f0, f1;
2211da51fcedSMatthew G. Knepley     PetscPointJac    g0, g1, g2, g3;
2212da51fcedSMatthew G. Knepley     PetscBdPointFunc f0Bd, f1Bd;
2213da51fcedSMatthew G. Knepley     PetscBdPointJac  g0Bd, g1Bd, g2Bd, g3Bd;
2214da51fcedSMatthew G. Knepley     PetscRiemannFunc r;
2215da51fcedSMatthew G. Knepley 
2216da51fcedSMatthew G. Knepley     ierr = PetscDSGetObjective(prob, f, &obj);CHKERRQ(ierr);
2217da51fcedSMatthew G. Knepley     ierr = PetscDSGetResidual(prob, f, &f0, &f1);CHKERRQ(ierr);
2218da51fcedSMatthew G. Knepley     ierr = PetscDSGetBdResidual(prob, f, &f0Bd, &f1Bd);CHKERRQ(ierr);
2219da51fcedSMatthew G. Knepley     ierr = PetscDSGetRiemannSolver(prob, f, &r);CHKERRQ(ierr);
2220da51fcedSMatthew G. Knepley     ierr = PetscDSSetObjective(newprob, f, obj);CHKERRQ(ierr);
2221da51fcedSMatthew G. Knepley     ierr = PetscDSSetResidual(newprob, f, f0, f1);CHKERRQ(ierr);
2222da51fcedSMatthew G. Knepley     ierr = PetscDSSetBdResidual(newprob, f, f0Bd, f1Bd);CHKERRQ(ierr);
2223da51fcedSMatthew G. Knepley     ierr = PetscDSSetRiemannSolver(newprob, f, r);CHKERRQ(ierr);
2224da51fcedSMatthew G. Knepley     for (g = 0; g < Nf; ++g) {
2225da51fcedSMatthew G. Knepley       ierr = PetscDSGetJacobian(prob, f, g, &g0, &g1, &g2, &g3);CHKERRQ(ierr);
2226da51fcedSMatthew G. Knepley       ierr = PetscDSGetBdJacobian(prob, f, g, &g0Bd, &g1Bd, &g2Bd, &g3Bd);CHKERRQ(ierr);
2227da51fcedSMatthew G. Knepley       ierr = PetscDSSetJacobian(newprob, f, g, g0, g1, g2, g3);CHKERRQ(ierr);
2228da51fcedSMatthew G. Knepley       ierr = PetscDSSetBdJacobian(newprob, f, g, g0Bd, g1Bd, g2Bd, g3Bd);CHKERRQ(ierr);
2229da51fcedSMatthew G. Knepley     }
2230da51fcedSMatthew G. Knepley   }
2231da51fcedSMatthew G. Knepley   PetscFunctionReturn(0);
2232da51fcedSMatthew G. Knepley }
2233da51fcedSMatthew G. Knepley 
2234da51fcedSMatthew G. Knepley #undef __FUNCT__
22352764a2aaSMatthew G. Knepley #define __FUNCT__ "PetscDSDestroy_Basic"
2236bc4ae4beSMatthew G. Knepley static PetscErrorCode PetscDSDestroy_Basic(PetscDS prob)
22372764a2aaSMatthew G. Knepley {
22382764a2aaSMatthew G. Knepley   PetscFunctionBegin;
22392764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
22402764a2aaSMatthew G. Knepley }
22412764a2aaSMatthew G. Knepley 
22422764a2aaSMatthew G. Knepley #undef __FUNCT__
22432764a2aaSMatthew G. Knepley #define __FUNCT__ "PetscDSInitialize_Basic"
2244bc4ae4beSMatthew G. Knepley static PetscErrorCode PetscDSInitialize_Basic(PetscDS prob)
22452764a2aaSMatthew G. Knepley {
22462764a2aaSMatthew G. Knepley   PetscFunctionBegin;
22472764a2aaSMatthew G. Knepley   prob->ops->setfromoptions = NULL;
22482764a2aaSMatthew G. Knepley   prob->ops->setup          = NULL;
22492764a2aaSMatthew G. Knepley   prob->ops->view           = NULL;
22502764a2aaSMatthew G. Knepley   prob->ops->destroy        = PetscDSDestroy_Basic;
22512764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
22522764a2aaSMatthew G. Knepley }
22532764a2aaSMatthew G. Knepley 
22542764a2aaSMatthew G. Knepley /*MC
22552764a2aaSMatthew G. Knepley   PETSCDSBASIC = "basic" - A discrete system with pointwise residual and boundary residual functions
22562764a2aaSMatthew G. Knepley 
22572764a2aaSMatthew G. Knepley   Level: intermediate
22582764a2aaSMatthew G. Knepley 
22592764a2aaSMatthew G. Knepley .seealso: PetscDSType, PetscDSCreate(), PetscDSSetType()
22602764a2aaSMatthew G. Knepley M*/
22612764a2aaSMatthew G. Knepley 
22622764a2aaSMatthew G. Knepley #undef __FUNCT__
22632764a2aaSMatthew G. Knepley #define __FUNCT__ "PetscDSCreate_Basic"
22642764a2aaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDSCreate_Basic(PetscDS prob)
22652764a2aaSMatthew G. Knepley {
22662764a2aaSMatthew G. Knepley   PetscDS_Basic *b;
22672764a2aaSMatthew G. Knepley   PetscErrorCode      ierr;
22682764a2aaSMatthew G. Knepley 
22692764a2aaSMatthew G. Knepley   PetscFunctionBegin;
22702764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCSPACE_CLASSID, 1);
22712764a2aaSMatthew G. Knepley   ierr       = PetscNewLog(prob, &b);CHKERRQ(ierr);
22722764a2aaSMatthew G. Knepley   prob->data = b;
22732764a2aaSMatthew G. Knepley 
22742764a2aaSMatthew G. Knepley   ierr = PetscDSInitialize_Basic(prob);CHKERRQ(ierr);
22752764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
22762764a2aaSMatthew G. Knepley }
2277