xref: /petsc/src/dm/dt/interface/dtds.c (revision 73107ff1babc248481eafc1deee82c9703b7a0a0)
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__ "PetscDSViewFromOptions"
2082764a2aaSMatthew G. Knepley /*
2092764a2aaSMatthew G. Knepley   PetscDSViewFromOptions - Processes command line options to determine if/how a PetscDS is to be viewed.
2102764a2aaSMatthew G. Knepley 
2112764a2aaSMatthew G. Knepley   Collective on PetscDS
2122764a2aaSMatthew G. Knepley 
2132764a2aaSMatthew G. Knepley   Input Parameters:
2142764a2aaSMatthew G. Knepley + prob   - the PetscDS
2152764a2aaSMatthew G. Knepley . prefix - prefix to use for viewing, or NULL to use prefix of 'rnd'
2162764a2aaSMatthew G. Knepley - optionname - option to activate viewing
2172764a2aaSMatthew G. Knepley 
2182764a2aaSMatthew G. Knepley   Level: intermediate
2192764a2aaSMatthew G. Knepley 
2202764a2aaSMatthew G. Knepley .keywords: PetscDS, view, options, database
2212764a2aaSMatthew G. Knepley .seealso: VecViewFromOptions(), MatViewFromOptions()
2222764a2aaSMatthew G. Knepley */
2232764a2aaSMatthew G. Knepley PetscErrorCode PetscDSViewFromOptions(PetscDS prob, const char prefix[], const char optionname[])
2242764a2aaSMatthew G. Knepley {
2252764a2aaSMatthew G. Knepley   PetscViewer       viewer;
2262764a2aaSMatthew G. Knepley   PetscViewerFormat format;
2272764a2aaSMatthew G. Knepley   PetscBool         flg;
2282764a2aaSMatthew G. Knepley   PetscErrorCode    ierr;
2292764a2aaSMatthew G. Knepley 
2302764a2aaSMatthew G. Knepley   PetscFunctionBegin;
2312764a2aaSMatthew G. Knepley   if (prefix) {ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject) prob), prefix, optionname, &viewer, &format, &flg);CHKERRQ(ierr);}
2322764a2aaSMatthew G. Knepley   else        {ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject) prob), ((PetscObject) prob)->prefix, optionname, &viewer, &format, &flg);CHKERRQ(ierr);}
2332764a2aaSMatthew G. Knepley   if (flg) {
2342764a2aaSMatthew G. Knepley     ierr = PetscViewerPushFormat(viewer, format);CHKERRQ(ierr);
2352764a2aaSMatthew G. Knepley     ierr = PetscDSView(prob, viewer);CHKERRQ(ierr);
2362764a2aaSMatthew G. Knepley     ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
2372764a2aaSMatthew G. Knepley     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
2382764a2aaSMatthew G. Knepley   }
2392764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
2402764a2aaSMatthew G. Knepley }
2412764a2aaSMatthew G. Knepley 
2422764a2aaSMatthew G. Knepley #undef __FUNCT__
2432764a2aaSMatthew G. Knepley #define __FUNCT__ "PetscDSSetFromOptions"
2442764a2aaSMatthew G. Knepley /*@
2452764a2aaSMatthew G. Knepley   PetscDSSetFromOptions - sets parameters in a PetscDS from the options database
2462764a2aaSMatthew G. Knepley 
2472764a2aaSMatthew G. Knepley   Collective on PetscDS
2482764a2aaSMatthew G. Knepley 
2492764a2aaSMatthew G. Knepley   Input Parameter:
2502764a2aaSMatthew G. Knepley . prob - the PetscDS object to set options for
2512764a2aaSMatthew G. Knepley 
2522764a2aaSMatthew G. Knepley   Options Database:
2532764a2aaSMatthew G. Knepley 
2542764a2aaSMatthew G. Knepley   Level: developer
2552764a2aaSMatthew G. Knepley 
2562764a2aaSMatthew G. Knepley .seealso PetscDSView()
2572764a2aaSMatthew G. Knepley @*/
2582764a2aaSMatthew G. Knepley PetscErrorCode PetscDSSetFromOptions(PetscDS prob)
2592764a2aaSMatthew G. Knepley {
2602764a2aaSMatthew G. Knepley   const char    *defaultType;
2612764a2aaSMatthew G. Knepley   char           name[256];
2622764a2aaSMatthew G. Knepley   PetscBool      flg;
2632764a2aaSMatthew G. Knepley   PetscErrorCode ierr;
2642764a2aaSMatthew G. Knepley 
2652764a2aaSMatthew G. Knepley   PetscFunctionBegin;
2662764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
2672764a2aaSMatthew G. Knepley   if (!((PetscObject) prob)->type_name) {
2682764a2aaSMatthew G. Knepley     defaultType = PETSCDSBASIC;
2692764a2aaSMatthew G. Knepley   } else {
2702764a2aaSMatthew G. Knepley     defaultType = ((PetscObject) prob)->type_name;
2712764a2aaSMatthew G. Knepley   }
2720f51fdf8SToby Isaac   ierr = PetscDSRegisterAll();CHKERRQ(ierr);
2732764a2aaSMatthew G. Knepley 
2742764a2aaSMatthew G. Knepley   ierr = PetscObjectOptionsBegin((PetscObject) prob);CHKERRQ(ierr);
2752764a2aaSMatthew G. Knepley   ierr = PetscOptionsFList("-petscds_type", "Discrete System", "PetscDSSetType", PetscDSList, defaultType, name, 256, &flg);CHKERRQ(ierr);
2762764a2aaSMatthew G. Knepley   if (flg) {
2772764a2aaSMatthew G. Knepley     ierr = PetscDSSetType(prob, name);CHKERRQ(ierr);
2782764a2aaSMatthew G. Knepley   } else if (!((PetscObject) prob)->type_name) {
2792764a2aaSMatthew G. Knepley     ierr = PetscDSSetType(prob, defaultType);CHKERRQ(ierr);
2802764a2aaSMatthew G. Knepley   }
2812764a2aaSMatthew G. Knepley   if (prob->ops->setfromoptions) {ierr = (*prob->ops->setfromoptions)(prob);CHKERRQ(ierr);}
2822764a2aaSMatthew G. Knepley   /* process any options handlers added with PetscObjectAddOptionsHandler() */
2832764a2aaSMatthew G. Knepley   ierr = PetscObjectProcessOptionsHandlers((PetscObject) prob);CHKERRQ(ierr);
2842764a2aaSMatthew G. Knepley   ierr = PetscOptionsEnd();CHKERRQ(ierr);
2852764a2aaSMatthew G. Knepley   ierr = PetscDSViewFromOptions(prob, NULL, "-petscds_view");CHKERRQ(ierr);
2862764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
2872764a2aaSMatthew G. Knepley }
2882764a2aaSMatthew G. Knepley 
2892764a2aaSMatthew G. Knepley #undef __FUNCT__
2902764a2aaSMatthew G. Knepley #define __FUNCT__ "PetscDSSetUp"
2912764a2aaSMatthew G. Knepley /*@C
2922764a2aaSMatthew G. Knepley   PetscDSSetUp - Construct data structures for the PetscDS
2932764a2aaSMatthew G. Knepley 
2942764a2aaSMatthew G. Knepley   Collective on PetscDS
2952764a2aaSMatthew G. Knepley 
2962764a2aaSMatthew G. Knepley   Input Parameter:
2972764a2aaSMatthew G. Knepley . prob - the PetscDS object to setup
2982764a2aaSMatthew G. Knepley 
2992764a2aaSMatthew G. Knepley   Level: developer
3002764a2aaSMatthew G. Knepley 
3012764a2aaSMatthew G. Knepley .seealso PetscDSView(), PetscDSDestroy()
3022764a2aaSMatthew G. Knepley @*/
3032764a2aaSMatthew G. Knepley PetscErrorCode PetscDSSetUp(PetscDS prob)
3042764a2aaSMatthew G. Knepley {
3052764a2aaSMatthew G. Knepley   const PetscInt Nf = prob->Nf;
3062764a2aaSMatthew G. Knepley   PetscInt       dim, work, NcMax = 0, NqMax = 0, f;
3072764a2aaSMatthew G. Knepley   PetscErrorCode ierr;
3082764a2aaSMatthew G. Knepley 
3092764a2aaSMatthew G. Knepley   PetscFunctionBegin;
3102764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
3112764a2aaSMatthew G. Knepley   if (prob->setup) PetscFunctionReturn(0);
3122764a2aaSMatthew G. Knepley   /* Calculate sizes */
3132764a2aaSMatthew G. Knepley   ierr = PetscDSGetSpatialDimension(prob, &dim);CHKERRQ(ierr);
3142764a2aaSMatthew G. Knepley   prob->totDim = prob->totDimBd = prob->totComp = 0;
3152764a2aaSMatthew G. Knepley   ierr = PetscMalloc4(Nf,&prob->basis,Nf,&prob->basisDer,Nf,&prob->basisBd,Nf,&prob->basisDerBd);CHKERRQ(ierr);
3162764a2aaSMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
3172764a2aaSMatthew G. Knepley     PetscFE         feBd = (PetscFE) prob->discBd[f];
3189de99aefSMatthew G. Knepley     PetscObject     obj;
3199de99aefSMatthew G. Knepley     PetscClassId    id;
3202764a2aaSMatthew G. Knepley     PetscQuadrature q;
3219de99aefSMatthew G. Knepley     PetscInt        Nq = 0, Nb, Nc;
3222764a2aaSMatthew G. Knepley 
3239de99aefSMatthew G. Knepley     ierr = PetscDSGetDiscretization(prob, f, &obj);CHKERRQ(ierr);
3249de99aefSMatthew G. Knepley     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
3259de99aefSMatthew G. Knepley     if (id == PETSCFE_CLASSID)      {
3269de99aefSMatthew G. Knepley       PetscFE fe = (PetscFE) obj;
3279de99aefSMatthew G. Knepley 
3282764a2aaSMatthew G. Knepley       ierr = PetscFEGetQuadrature(fe, &q);CHKERRQ(ierr);
3292764a2aaSMatthew G. Knepley       ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr);
3302764a2aaSMatthew G. Knepley       ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
3312764a2aaSMatthew G. Knepley       ierr = PetscFEGetDefaultTabulation(fe, &prob->basis[f], &prob->basisDer[f], NULL);CHKERRQ(ierr);
3329de99aefSMatthew G. Knepley     } else if (id == PETSCFV_CLASSID) {
3339de99aefSMatthew G. Knepley       PetscFV fv = (PetscFV) obj;
3349de99aefSMatthew G. Knepley 
3359de99aefSMatthew G. Knepley       ierr = PetscFVGetQuadrature(fv, &q);CHKERRQ(ierr);
3369de99aefSMatthew G. Knepley       Nb   = 1;
3379de99aefSMatthew G. Knepley       ierr = PetscFVGetNumComponents(fv, &Nc);CHKERRQ(ierr);
3386c1a3d01SMatthew G. Knepley       ierr = PetscFVGetDefaultTabulation(fv, &prob->basis[f], &prob->basisDer[f], NULL);CHKERRQ(ierr);
339abac5ca0SMatthew G. Knepley     } else SETERRQ1(PetscObjectComm((PetscObject) prob), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", f);
3409de99aefSMatthew G. Knepley     if (q) {ierr = PetscQuadratureGetData(q, NULL, &Nq, NULL, NULL);CHKERRQ(ierr);}
3412764a2aaSMatthew G. Knepley     NqMax          = PetscMax(NqMax, Nq);
3422764a2aaSMatthew G. Knepley     NcMax          = PetscMax(NcMax, Nc);
3432764a2aaSMatthew G. Knepley     prob->totDim  += Nb*Nc;
3442764a2aaSMatthew G. Knepley     prob->totComp += Nc;
3452764a2aaSMatthew G. Knepley     if (feBd) {
3462764a2aaSMatthew G. Knepley       ierr = PetscFEGetDimension(feBd, &Nb);CHKERRQ(ierr);
3472764a2aaSMatthew G. Knepley       ierr = PetscFEGetNumComponents(feBd, &Nc);CHKERRQ(ierr);
3482764a2aaSMatthew G. Knepley       ierr = PetscFEGetDefaultTabulation(feBd, &prob->basisBd[f], &prob->basisDerBd[f], NULL);CHKERRQ(ierr);
3492764a2aaSMatthew G. Knepley       prob->totDimBd += Nb*Nc;
3502764a2aaSMatthew G. Knepley     }
3512764a2aaSMatthew G. Knepley   }
3522764a2aaSMatthew G. Knepley   work = PetscMax(prob->totComp*dim, PetscSqr(NcMax*dim));
3532764a2aaSMatthew G. Knepley   /* Allocate works space */
3542764a2aaSMatthew 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);
3552764a2aaSMatthew 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);
3562764a2aaSMatthew G. Knepley   if (prob->ops->setup) {ierr = (*prob->ops->setup)(prob);CHKERRQ(ierr);}
3572764a2aaSMatthew G. Knepley   prob->setup = PETSC_TRUE;
3582764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
3592764a2aaSMatthew G. Knepley }
3602764a2aaSMatthew G. Knepley 
3612764a2aaSMatthew G. Knepley #undef __FUNCT__
3622764a2aaSMatthew G. Knepley #define __FUNCT__ "PetscDSDestroyStructs_Static"
3632764a2aaSMatthew G. Knepley static PetscErrorCode PetscDSDestroyStructs_Static(PetscDS prob)
3642764a2aaSMatthew G. Knepley {
3652764a2aaSMatthew G. Knepley   PetscErrorCode ierr;
3662764a2aaSMatthew G. Knepley 
3672764a2aaSMatthew G. Knepley   PetscFunctionBegin;
3682764a2aaSMatthew G. Knepley   ierr = PetscFree4(prob->basis,prob->basisDer,prob->basisBd,prob->basisDerBd);CHKERRQ(ierr);
3692764a2aaSMatthew G. Knepley   ierr = PetscFree5(prob->u,prob->u_t,prob->u_x,prob->x,prob->refSpaceDer);CHKERRQ(ierr);
3702764a2aaSMatthew G. Knepley   ierr = PetscFree6(prob->f0,prob->f1,prob->g0,prob->g1,prob->g2,prob->g3);CHKERRQ(ierr);
3712764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
3722764a2aaSMatthew G. Knepley }
3732764a2aaSMatthew G. Knepley 
3742764a2aaSMatthew G. Knepley #undef __FUNCT__
3752764a2aaSMatthew G. Knepley #define __FUNCT__ "PetscDSEnlarge_Static"
3762764a2aaSMatthew G. Knepley static PetscErrorCode PetscDSEnlarge_Static(PetscDS prob, PetscInt NfNew)
3772764a2aaSMatthew G. Knepley {
3782764a2aaSMatthew G. Knepley   PetscObject   *tmpd, *tmpdbd;
379a6cbbb48SMatthew G. Knepley   PetscBool     *tmpi, *tmpa;
3802764a2aaSMatthew G. Knepley   PointFunc     *tmpobj, *tmpf, *tmpg;
3812764a2aaSMatthew G. Knepley   BdPointFunc   *tmpfbd, *tmpgbd;
3820c2f2876SMatthew G. Knepley   RiemannFunc   *tmpr;
3830c2f2876SMatthew G. Knepley   void         **tmpctx;
384a6cbbb48SMatthew G. Knepley   PetscInt       Nf = prob->Nf, f, i;
3852764a2aaSMatthew G. Knepley   PetscErrorCode ierr;
3862764a2aaSMatthew G. Knepley 
3872764a2aaSMatthew G. Knepley   PetscFunctionBegin;
3882764a2aaSMatthew G. Knepley   if (Nf >= NfNew) PetscFunctionReturn(0);
3892764a2aaSMatthew G. Knepley   prob->setup = PETSC_FALSE;
3902764a2aaSMatthew G. Knepley   ierr = PetscDSDestroyStructs_Static(prob);CHKERRQ(ierr);
391a6cbbb48SMatthew G. Knepley   ierr = PetscMalloc4(NfNew, &tmpd, NfNew, &tmpdbd, NfNew, &tmpi, NfNew*2, &tmpa);CHKERRQ(ierr);
392a6cbbb48SMatthew 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];}
393a6cbbb48SMatthew G. Knepley   for (f = Nf; f < NfNew; ++f) {ierr = PetscContainerCreate(PetscObjectComm((PetscObject) prob), (PetscContainer *) &tmpd[f]);CHKERRQ(ierr);
394a6cbbb48SMatthew G. Knepley     tmpdbd[f] = NULL; tmpi[f] = PETSC_TRUE; tmpa[f*2+0] = PETSC_FALSE; tmpa[f*2+1] = PETSC_TRUE;}
395a6cbbb48SMatthew G. Knepley   ierr = PetscFree4(prob->disc, prob->discBd, prob->implicit, prob->adjacency);CHKERRQ(ierr);
3962764a2aaSMatthew G. Knepley   prob->Nf        = NfNew;
3972764a2aaSMatthew G. Knepley   prob->disc      = tmpd;
398a6cbbb48SMatthew G. Knepley   prob->discBd    = tmpdbd;
399249df284SMatthew G. Knepley   prob->implicit  = tmpi;
400a6cbbb48SMatthew G. Knepley   prob->adjacency = tmpa;
4010c2f2876SMatthew G. Knepley   ierr = PetscCalloc5(NfNew, &tmpobj, NfNew*2, &tmpf, NfNew*NfNew*4, &tmpg, NfNew, &tmpr, NfNew, &tmpctx);CHKERRQ(ierr);
4022764a2aaSMatthew G. Knepley   for (f = 0; f < Nf; ++f) tmpobj[f] = prob->obj[f];
4032764a2aaSMatthew G. Knepley   for (f = 0; f < Nf*2; ++f) tmpf[f] = prob->f[f];
4042764a2aaSMatthew G. Knepley   for (f = 0; f < Nf*Nf*4; ++f) tmpg[f] = prob->g[f];
4050c2f2876SMatthew G. Knepley   for (f = 0; f < Nf; ++f) tmpr[f] = prob->r[f];
4060c2f2876SMatthew G. Knepley   for (f = 0; f < Nf; ++f) tmpctx[f] = prob->ctx[f];
4072764a2aaSMatthew G. Knepley   for (f = Nf; f < NfNew; ++f) tmpobj[f] = NULL;
4082764a2aaSMatthew G. Knepley   for (f = Nf*2; f < NfNew*2; ++f) tmpf[f] = NULL;
4092764a2aaSMatthew G. Knepley   for (f = Nf*Nf*4; f < NfNew*NfNew*4; ++f) tmpg[f] = NULL;
4100c2f2876SMatthew G. Knepley   for (f = Nf; f < NfNew; ++f) tmpr[f] = NULL;
4110c2f2876SMatthew G. Knepley   for (f = Nf; f < NfNew; ++f) tmpctx[f] = NULL;
4120c2f2876SMatthew G. Knepley   ierr = PetscFree5(prob->obj, prob->f, prob->g, prob->r, prob->ctx);CHKERRQ(ierr);
4132764a2aaSMatthew G. Knepley   prob->obj = tmpobj;
4142764a2aaSMatthew G. Knepley   prob->f   = tmpf;
4152764a2aaSMatthew G. Knepley   prob->g   = tmpg;
4160c2f2876SMatthew G. Knepley   prob->r   = tmpr;
4170c2f2876SMatthew G. Knepley   prob->ctx = tmpctx;
4182764a2aaSMatthew G. Knepley   ierr = PetscCalloc2(NfNew*2, &tmpfbd, NfNew*NfNew*4, &tmpgbd);CHKERRQ(ierr);
4192764a2aaSMatthew G. Knepley   for (f = 0; f < Nf*2; ++f) tmpfbd[f] = prob->fBd[f];
4202764a2aaSMatthew G. Knepley   for (f = 0; f < Nf*Nf*4; ++f) tmpgbd[f] = prob->gBd[f];
4212764a2aaSMatthew G. Knepley   for (f = Nf*2; f < NfNew*2; ++f) tmpfbd[f] = NULL;
4222764a2aaSMatthew G. Knepley   for (f = Nf*Nf*4; f < NfNew*NfNew*4; ++f) tmpgbd[f] = NULL;
4232764a2aaSMatthew G. Knepley   ierr = PetscFree2(prob->fBd, prob->gBd);CHKERRQ(ierr);
4242764a2aaSMatthew G. Knepley   prob->fBd = tmpfbd;
4252764a2aaSMatthew G. Knepley   prob->gBd = tmpgbd;
4262764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
4272764a2aaSMatthew G. Knepley }
4282764a2aaSMatthew G. Knepley 
4292764a2aaSMatthew G. Knepley #undef __FUNCT__
4302764a2aaSMatthew G. Knepley #define __FUNCT__ "PetscDSDestroy"
4312764a2aaSMatthew G. Knepley /*@
4322764a2aaSMatthew G. Knepley   PetscDSDestroy - Destroys a PetscDS object
4332764a2aaSMatthew G. Knepley 
4342764a2aaSMatthew G. Knepley   Collective on PetscDS
4352764a2aaSMatthew G. Knepley 
4362764a2aaSMatthew G. Knepley   Input Parameter:
4372764a2aaSMatthew G. Knepley . prob - the PetscDS object to destroy
4382764a2aaSMatthew G. Knepley 
4392764a2aaSMatthew G. Knepley   Level: developer
4402764a2aaSMatthew G. Knepley 
4412764a2aaSMatthew G. Knepley .seealso PetscDSView()
4422764a2aaSMatthew G. Knepley @*/
4432764a2aaSMatthew G. Knepley PetscErrorCode PetscDSDestroy(PetscDS *prob)
4442764a2aaSMatthew G. Knepley {
4452764a2aaSMatthew G. Knepley   PetscInt       f;
4462764a2aaSMatthew G. Knepley   PetscErrorCode ierr;
4472764a2aaSMatthew G. Knepley 
4482764a2aaSMatthew G. Knepley   PetscFunctionBegin;
4492764a2aaSMatthew G. Knepley   if (!*prob) PetscFunctionReturn(0);
4502764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific((*prob), PETSCDS_CLASSID, 1);
4512764a2aaSMatthew G. Knepley 
4522764a2aaSMatthew G. Knepley   if (--((PetscObject)(*prob))->refct > 0) {*prob = 0; PetscFunctionReturn(0);}
4532764a2aaSMatthew G. Knepley   ((PetscObject) (*prob))->refct = 0;
4542764a2aaSMatthew G. Knepley   ierr = PetscDSDestroyStructs_Static(*prob);CHKERRQ(ierr);
4552764a2aaSMatthew G. Knepley   for (f = 0; f < (*prob)->Nf; ++f) {
4562764a2aaSMatthew G. Knepley     ierr = PetscObjectDereference((*prob)->disc[f]);CHKERRQ(ierr);
4572764a2aaSMatthew G. Knepley     ierr = PetscObjectDereference((*prob)->discBd[f]);CHKERRQ(ierr);
4582764a2aaSMatthew G. Knepley   }
459a6cbbb48SMatthew G. Knepley   ierr = PetscFree4((*prob)->disc, (*prob)->discBd, (*prob)->implicit, (*prob)->adjacency);CHKERRQ(ierr);
4600c2f2876SMatthew G. Knepley   ierr = PetscFree5((*prob)->obj,(*prob)->f,(*prob)->g,(*prob)->r,(*prob)->ctx);CHKERRQ(ierr);
4612764a2aaSMatthew G. Knepley   ierr = PetscFree2((*prob)->fBd,(*prob)->gBd);CHKERRQ(ierr);
4622764a2aaSMatthew G. Knepley   if ((*prob)->ops->destroy) {ierr = (*(*prob)->ops->destroy)(*prob);CHKERRQ(ierr);}
4632764a2aaSMatthew G. Knepley   ierr = PetscHeaderDestroy(prob);CHKERRQ(ierr);
4642764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
4652764a2aaSMatthew G. Knepley }
4662764a2aaSMatthew G. Knepley 
4672764a2aaSMatthew G. Knepley #undef __FUNCT__
4682764a2aaSMatthew G. Knepley #define __FUNCT__ "PetscDSCreate"
4692764a2aaSMatthew G. Knepley /*@
4702764a2aaSMatthew G. Knepley   PetscDSCreate - Creates an empty PetscDS object. The type can then be set with PetscDSSetType().
4712764a2aaSMatthew G. Knepley 
4722764a2aaSMatthew G. Knepley   Collective on MPI_Comm
4732764a2aaSMatthew G. Knepley 
4742764a2aaSMatthew G. Knepley   Input Parameter:
4752764a2aaSMatthew G. Knepley . comm - The communicator for the PetscDS object
4762764a2aaSMatthew G. Knepley 
4772764a2aaSMatthew G. Knepley   Output Parameter:
4782764a2aaSMatthew G. Knepley . prob - The PetscDS object
4792764a2aaSMatthew G. Knepley 
4802764a2aaSMatthew G. Knepley   Level: beginner
4812764a2aaSMatthew G. Knepley 
4822764a2aaSMatthew G. Knepley .seealso: PetscDSSetType(), PETSCDSBASIC
4832764a2aaSMatthew G. Knepley @*/
4842764a2aaSMatthew G. Knepley PetscErrorCode PetscDSCreate(MPI_Comm comm, PetscDS *prob)
4852764a2aaSMatthew G. Knepley {
4862764a2aaSMatthew G. Knepley   PetscDS   p;
4872764a2aaSMatthew G. Knepley   PetscErrorCode ierr;
4882764a2aaSMatthew G. Knepley 
4892764a2aaSMatthew G. Knepley   PetscFunctionBegin;
4902764a2aaSMatthew G. Knepley   PetscValidPointer(prob, 2);
4912764a2aaSMatthew G. Knepley   *prob  = NULL;
4922764a2aaSMatthew G. Knepley   ierr = PetscDSInitializePackage();CHKERRQ(ierr);
4932764a2aaSMatthew G. Knepley 
494*73107ff1SLisandro Dalcin   ierr = PetscHeaderCreate(p, PETSCDS_CLASSID, "PetscDS", "Discrete System", "PetscDS", comm, PetscDSDestroy, PetscDSView);CHKERRQ(ierr);
4952764a2aaSMatthew G. Knepley 
4962764a2aaSMatthew G. Knepley   p->Nf    = 0;
4972764a2aaSMatthew G. Knepley   p->setup = PETSC_FALSE;
4982764a2aaSMatthew G. Knepley 
4992764a2aaSMatthew G. Knepley   *prob = p;
5002764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
5012764a2aaSMatthew G. Knepley }
5022764a2aaSMatthew G. Knepley 
5032764a2aaSMatthew G. Knepley #undef __FUNCT__
5042764a2aaSMatthew G. Knepley #define __FUNCT__ "PetscDSGetNumFields"
505bc4ae4beSMatthew G. Knepley /*@
506bc4ae4beSMatthew G. Knepley   PetscDSGetNumFields - Returns the number of fields in the DS
507bc4ae4beSMatthew G. Knepley 
508bc4ae4beSMatthew G. Knepley   Not collective
509bc4ae4beSMatthew G. Knepley 
510bc4ae4beSMatthew G. Knepley   Input Parameter:
511bc4ae4beSMatthew G. Knepley . prob - The PetscDS object
512bc4ae4beSMatthew G. Knepley 
513bc4ae4beSMatthew G. Knepley   Output Parameter:
514bc4ae4beSMatthew G. Knepley . Nf - The number of fields
515bc4ae4beSMatthew G. Knepley 
516bc4ae4beSMatthew G. Knepley   Level: beginner
517bc4ae4beSMatthew G. Knepley 
518bc4ae4beSMatthew G. Knepley .seealso: PetscDSGetSpatialDimension(), PetscDSCreate()
519bc4ae4beSMatthew G. Knepley @*/
5202764a2aaSMatthew G. Knepley PetscErrorCode PetscDSGetNumFields(PetscDS prob, PetscInt *Nf)
5212764a2aaSMatthew G. Knepley {
5222764a2aaSMatthew G. Knepley   PetscFunctionBegin;
5232764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
5242764a2aaSMatthew G. Knepley   PetscValidPointer(Nf, 2);
5252764a2aaSMatthew G. Knepley   *Nf = prob->Nf;
5262764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
5272764a2aaSMatthew G. Knepley }
5282764a2aaSMatthew G. Knepley 
5292764a2aaSMatthew G. Knepley #undef __FUNCT__
5302764a2aaSMatthew G. Knepley #define __FUNCT__ "PetscDSGetSpatialDimension"
531bc4ae4beSMatthew G. Knepley /*@
532bc4ae4beSMatthew G. Knepley   PetscDSGetSpatialDimension - Returns the spatial dimension of the DS
533bc4ae4beSMatthew G. Knepley 
534bc4ae4beSMatthew G. Knepley   Not collective
535bc4ae4beSMatthew G. Knepley 
536bc4ae4beSMatthew G. Knepley   Input Parameter:
537bc4ae4beSMatthew G. Knepley . prob - The PetscDS object
538bc4ae4beSMatthew G. Knepley 
539bc4ae4beSMatthew G. Knepley   Output Parameter:
540bc4ae4beSMatthew G. Knepley . dim - The spatial dimension
541bc4ae4beSMatthew G. Knepley 
542bc4ae4beSMatthew G. Knepley   Level: beginner
543bc4ae4beSMatthew G. Knepley 
544bc4ae4beSMatthew G. Knepley .seealso: PetscDSGetNumFields(), PetscDSCreate()
545bc4ae4beSMatthew G. Knepley @*/
5462764a2aaSMatthew G. Knepley PetscErrorCode PetscDSGetSpatialDimension(PetscDS prob, PetscInt *dim)
5472764a2aaSMatthew G. Knepley {
5482764a2aaSMatthew G. Knepley   PetscErrorCode ierr;
5492764a2aaSMatthew G. Knepley 
5502764a2aaSMatthew G. Knepley   PetscFunctionBegin;
5512764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
5522764a2aaSMatthew G. Knepley   PetscValidPointer(dim, 2);
5532764a2aaSMatthew G. Knepley   *dim = 0;
5549de99aefSMatthew G. Knepley   if (prob->Nf) {
5559de99aefSMatthew G. Knepley     PetscObject  obj;
5569de99aefSMatthew G. Knepley     PetscClassId id;
5579de99aefSMatthew G. Knepley 
5589de99aefSMatthew G. Knepley     ierr = PetscDSGetDiscretization(prob, 0, &obj);CHKERRQ(ierr);
5599de99aefSMatthew G. Knepley     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
5609de99aefSMatthew G. Knepley     if (id == PETSCFE_CLASSID)      {ierr = PetscFEGetSpatialDimension((PetscFE) obj, dim);CHKERRQ(ierr);}
5619de99aefSMatthew G. Knepley     else if (id == PETSCFV_CLASSID) {ierr = PetscFVGetSpatialDimension((PetscFV) obj, dim);CHKERRQ(ierr);}
5629de99aefSMatthew G. Knepley     else SETERRQ1(PetscObjectComm((PetscObject) prob), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", 0);
5639de99aefSMatthew G. Knepley   }
5642764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
5652764a2aaSMatthew G. Knepley }
5662764a2aaSMatthew G. Knepley 
5672764a2aaSMatthew G. Knepley #undef __FUNCT__
5682764a2aaSMatthew G. Knepley #define __FUNCT__ "PetscDSGetTotalDimension"
569bc4ae4beSMatthew G. Knepley /*@
570bc4ae4beSMatthew G. Knepley   PetscDSGetTotalDimension - Returns the total size of the approximation space for this system
571bc4ae4beSMatthew G. Knepley 
572bc4ae4beSMatthew G. Knepley   Not collective
573bc4ae4beSMatthew G. Knepley 
574bc4ae4beSMatthew G. Knepley   Input Parameter:
575bc4ae4beSMatthew G. Knepley . prob - The PetscDS object
576bc4ae4beSMatthew G. Knepley 
577bc4ae4beSMatthew G. Knepley   Output Parameter:
578bc4ae4beSMatthew G. Knepley . dim - The total problem dimension
579bc4ae4beSMatthew G. Knepley 
580bc4ae4beSMatthew G. Knepley   Level: beginner
581bc4ae4beSMatthew G. Knepley 
582bc4ae4beSMatthew G. Knepley .seealso: PetscDSGetNumFields(), PetscDSCreate()
583bc4ae4beSMatthew G. Knepley @*/
5842764a2aaSMatthew G. Knepley PetscErrorCode PetscDSGetTotalDimension(PetscDS prob, PetscInt *dim)
5852764a2aaSMatthew G. Knepley {
5862764a2aaSMatthew G. Knepley   PetscErrorCode ierr;
5872764a2aaSMatthew G. Knepley 
5882764a2aaSMatthew G. Knepley   PetscFunctionBegin;
5892764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
5902764a2aaSMatthew G. Knepley   ierr = PetscDSSetUp(prob);CHKERRQ(ierr);
5912764a2aaSMatthew G. Knepley   PetscValidPointer(dim, 2);
5922764a2aaSMatthew G. Knepley   *dim = prob->totDim;
5932764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
5942764a2aaSMatthew G. Knepley }
5952764a2aaSMatthew G. Knepley 
5962764a2aaSMatthew G. Knepley #undef __FUNCT__
5972764a2aaSMatthew G. Knepley #define __FUNCT__ "PetscDSGetTotalBdDimension"
598bc4ae4beSMatthew G. Knepley /*@
599c3ac4435SMatthew G. Knepley   PetscDSGetTotalBdDimension - Returns the total size of the boundary approximation space for this system
600bc4ae4beSMatthew G. Knepley 
601bc4ae4beSMatthew G. Knepley   Not collective
602bc4ae4beSMatthew G. Knepley 
603bc4ae4beSMatthew G. Knepley   Input Parameter:
604bc4ae4beSMatthew G. Knepley . prob - The PetscDS object
605bc4ae4beSMatthew G. Knepley 
606bc4ae4beSMatthew G. Knepley   Output Parameter:
607bc4ae4beSMatthew G. Knepley . dim - The total boundary problem dimension
608bc4ae4beSMatthew G. Knepley 
609bc4ae4beSMatthew G. Knepley   Level: beginner
610bc4ae4beSMatthew G. Knepley 
611bc4ae4beSMatthew G. Knepley .seealso: PetscDSGetNumFields(), PetscDSCreate()
612bc4ae4beSMatthew G. Knepley @*/
6132764a2aaSMatthew G. Knepley PetscErrorCode PetscDSGetTotalBdDimension(PetscDS prob, PetscInt *dim)
6142764a2aaSMatthew G. Knepley {
6152764a2aaSMatthew G. Knepley   PetscErrorCode ierr;
6162764a2aaSMatthew G. Knepley 
6172764a2aaSMatthew G. Knepley   PetscFunctionBegin;
6182764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
6192764a2aaSMatthew G. Knepley   ierr = PetscDSSetUp(prob);CHKERRQ(ierr);
6202764a2aaSMatthew G. Knepley   PetscValidPointer(dim, 2);
6212764a2aaSMatthew G. Knepley   *dim = prob->totDimBd;
6222764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
6232764a2aaSMatthew G. Knepley }
6242764a2aaSMatthew G. Knepley 
6252764a2aaSMatthew G. Knepley #undef __FUNCT__
6262764a2aaSMatthew G. Knepley #define __FUNCT__ "PetscDSGetTotalComponents"
627bc4ae4beSMatthew G. Knepley /*@
628bc4ae4beSMatthew G. Knepley   PetscDSGetTotalComponents - Returns the total number of components in this system
629bc4ae4beSMatthew G. Knepley 
630bc4ae4beSMatthew G. Knepley   Not collective
631bc4ae4beSMatthew G. Knepley 
632bc4ae4beSMatthew G. Knepley   Input Parameter:
633bc4ae4beSMatthew G. Knepley . prob - The PetscDS object
634bc4ae4beSMatthew G. Knepley 
635bc4ae4beSMatthew G. Knepley   Output Parameter:
636bc4ae4beSMatthew G. Knepley . dim - The total number of components
637bc4ae4beSMatthew G. Knepley 
638bc4ae4beSMatthew G. Knepley   Level: beginner
639bc4ae4beSMatthew G. Knepley 
640bc4ae4beSMatthew G. Knepley .seealso: PetscDSGetNumFields(), PetscDSCreate()
641bc4ae4beSMatthew G. Knepley @*/
6422764a2aaSMatthew G. Knepley PetscErrorCode PetscDSGetTotalComponents(PetscDS prob, PetscInt *Nc)
6432764a2aaSMatthew G. Knepley {
6442764a2aaSMatthew G. Knepley   PetscErrorCode ierr;
6452764a2aaSMatthew G. Knepley 
6462764a2aaSMatthew G. Knepley   PetscFunctionBegin;
6472764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
6482764a2aaSMatthew G. Knepley   ierr = PetscDSSetUp(prob);CHKERRQ(ierr);
6492764a2aaSMatthew G. Knepley   PetscValidPointer(Nc, 2);
6502764a2aaSMatthew G. Knepley   *Nc = prob->totComp;
6512764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
6522764a2aaSMatthew G. Knepley }
6532764a2aaSMatthew G. Knepley 
6542764a2aaSMatthew G. Knepley #undef __FUNCT__
6552764a2aaSMatthew G. Knepley #define __FUNCT__ "PetscDSGetDiscretization"
656bc4ae4beSMatthew G. Knepley /*@
657bc4ae4beSMatthew G. Knepley   PetscDSGetDiscretization - Returns the discretization object for the given field
658bc4ae4beSMatthew G. Knepley 
659bc4ae4beSMatthew G. Knepley   Not collective
660bc4ae4beSMatthew G. Knepley 
661bc4ae4beSMatthew G. Knepley   Input Parameters:
662bc4ae4beSMatthew G. Knepley + prob - The PetscDS object
663bc4ae4beSMatthew G. Knepley - f - The field number
664bc4ae4beSMatthew G. Knepley 
665bc4ae4beSMatthew G. Knepley   Output Parameter:
666bc4ae4beSMatthew G. Knepley . disc - The discretization object
667bc4ae4beSMatthew G. Knepley 
668bc4ae4beSMatthew G. Knepley   Level: beginner
669bc4ae4beSMatthew G. Knepley 
670bc4ae4beSMatthew G. Knepley .seealso: PetscDSSetDiscretization(), PetscDSAddDiscretization(), PetscDSGetBdDiscretization(), PetscDSGetNumFields(), PetscDSCreate()
671bc4ae4beSMatthew G. Knepley @*/
6722764a2aaSMatthew G. Knepley PetscErrorCode PetscDSGetDiscretization(PetscDS prob, PetscInt f, PetscObject *disc)
6732764a2aaSMatthew G. Knepley {
6742764a2aaSMatthew G. Knepley   PetscFunctionBegin;
6752764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
6762764a2aaSMatthew G. Knepley   PetscValidPointer(disc, 3);
6772764a2aaSMatthew 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);
6782764a2aaSMatthew G. Knepley   *disc = prob->disc[f];
6792764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
6802764a2aaSMatthew G. Knepley }
6812764a2aaSMatthew G. Knepley 
6822764a2aaSMatthew G. Knepley #undef __FUNCT__
6832764a2aaSMatthew G. Knepley #define __FUNCT__ "PetscDSGetBdDiscretization"
684bc4ae4beSMatthew G. Knepley /*@
685bc4ae4beSMatthew G. Knepley   PetscDSGetBdDiscretization - Returns the boundary discretization object for the given field
686bc4ae4beSMatthew G. Knepley 
687bc4ae4beSMatthew G. Knepley   Not collective
688bc4ae4beSMatthew G. Knepley 
689bc4ae4beSMatthew G. Knepley   Input Parameters:
690bc4ae4beSMatthew G. Knepley + prob - The PetscDS object
691bc4ae4beSMatthew G. Knepley - f - The field number
692bc4ae4beSMatthew G. Knepley 
693bc4ae4beSMatthew G. Knepley   Output Parameter:
694bc4ae4beSMatthew G. Knepley . disc - The boundary discretization object
695bc4ae4beSMatthew G. Knepley 
696bc4ae4beSMatthew G. Knepley   Level: beginner
697bc4ae4beSMatthew G. Knepley 
698bc4ae4beSMatthew G. Knepley .seealso: PetscDSSetBdDiscretization(), PetscDSAddBdDiscretization(), PetscDSGetDiscretization(), PetscDSGetNumFields(), PetscDSCreate()
699bc4ae4beSMatthew G. Knepley @*/
7002764a2aaSMatthew G. Knepley PetscErrorCode PetscDSGetBdDiscretization(PetscDS prob, PetscInt f, PetscObject *disc)
7012764a2aaSMatthew G. Knepley {
7022764a2aaSMatthew G. Knepley   PetscFunctionBegin;
7032764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
7042764a2aaSMatthew G. Knepley   PetscValidPointer(disc, 3);
7052764a2aaSMatthew 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);
7062764a2aaSMatthew G. Knepley   *disc = prob->discBd[f];
7072764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
7082764a2aaSMatthew G. Knepley }
7092764a2aaSMatthew G. Knepley 
7102764a2aaSMatthew G. Knepley #undef __FUNCT__
7112764a2aaSMatthew G. Knepley #define __FUNCT__ "PetscDSSetDiscretization"
712bc4ae4beSMatthew G. Knepley /*@
713bc4ae4beSMatthew G. Knepley   PetscDSSetDiscretization - Sets the discretization object for the given field
714bc4ae4beSMatthew G. Knepley 
715bc4ae4beSMatthew G. Knepley   Not collective
716bc4ae4beSMatthew G. Knepley 
717bc4ae4beSMatthew G. Knepley   Input Parameters:
718bc4ae4beSMatthew G. Knepley + prob - The PetscDS object
719bc4ae4beSMatthew G. Knepley . f - The field number
720bc4ae4beSMatthew G. Knepley - disc - The discretization object
721bc4ae4beSMatthew G. Knepley 
722bc4ae4beSMatthew G. Knepley   Level: beginner
723bc4ae4beSMatthew G. Knepley 
724bc4ae4beSMatthew G. Knepley .seealso: PetscDSGetDiscretization(), PetscDSAddDiscretization(), PetscDSGetNumFields(), PetscDSCreate()
725bc4ae4beSMatthew G. Knepley @*/
7262764a2aaSMatthew G. Knepley PetscErrorCode PetscDSSetDiscretization(PetscDS prob, PetscInt f, PetscObject disc)
7272764a2aaSMatthew G. Knepley {
7282764a2aaSMatthew G. Knepley   PetscErrorCode ierr;
7292764a2aaSMatthew G. Knepley 
7302764a2aaSMatthew G. Knepley   PetscFunctionBegin;
7312764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
7322764a2aaSMatthew G. Knepley   PetscValidPointer(disc, 3);
7332764a2aaSMatthew G. Knepley   if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f);
7342764a2aaSMatthew G. Knepley   ierr = PetscDSEnlarge_Static(prob, f+1);CHKERRQ(ierr);
7352764a2aaSMatthew G. Knepley   if (prob->disc[f]) {ierr = PetscObjectDereference(prob->disc[f]);CHKERRQ(ierr);}
7362764a2aaSMatthew G. Knepley   prob->disc[f] = disc;
7372764a2aaSMatthew G. Knepley   ierr = PetscObjectReference(disc);CHKERRQ(ierr);
738249df284SMatthew G. Knepley   {
739249df284SMatthew G. Knepley     PetscClassId id;
740249df284SMatthew G. Knepley 
741249df284SMatthew G. Knepley     ierr = PetscObjectGetClassId(disc, &id);CHKERRQ(ierr);
742a6cbbb48SMatthew G. Knepley     if (id == PETSCFV_CLASSID) {
743a6cbbb48SMatthew G. Knepley       prob->implicit[f]      = PETSC_FALSE;
744a6cbbb48SMatthew G. Knepley       prob->adjacency[f*2+0] = PETSC_TRUE;
745a6cbbb48SMatthew G. Knepley       prob->adjacency[f*2+1] = PETSC_FALSE;
746a6cbbb48SMatthew G. Knepley     }
747249df284SMatthew G. Knepley   }
7482764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
7492764a2aaSMatthew G. Knepley }
7502764a2aaSMatthew G. Knepley 
7512764a2aaSMatthew G. Knepley #undef __FUNCT__
7522764a2aaSMatthew G. Knepley #define __FUNCT__ "PetscDSSetBdDiscretization"
753bc4ae4beSMatthew G. Knepley /*@
754bc4ae4beSMatthew G. Knepley   PetscDSSetBdDiscretization - Sets the boundary discretization object for the given field
755bc4ae4beSMatthew G. Knepley 
756bc4ae4beSMatthew G. Knepley   Not collective
757bc4ae4beSMatthew G. Knepley 
758bc4ae4beSMatthew G. Knepley   Input Parameters:
759bc4ae4beSMatthew G. Knepley + prob - The PetscDS object
760bc4ae4beSMatthew G. Knepley . f - The field number
761bc4ae4beSMatthew G. Knepley - disc - The boundary discretization object
762bc4ae4beSMatthew G. Knepley 
763bc4ae4beSMatthew G. Knepley   Level: beginner
764bc4ae4beSMatthew G. Knepley 
765bc4ae4beSMatthew G. Knepley .seealso: PetscDSGetBdDiscretization(), PetscDSAddBdDiscretization(), PetscDSGetNumFields(), PetscDSCreate()
766bc4ae4beSMatthew G. Knepley @*/
7672764a2aaSMatthew G. Knepley PetscErrorCode PetscDSSetBdDiscretization(PetscDS prob, PetscInt f, PetscObject disc)
7682764a2aaSMatthew G. Knepley {
7692764a2aaSMatthew G. Knepley   PetscErrorCode ierr;
7702764a2aaSMatthew G. Knepley 
7712764a2aaSMatthew G. Knepley   PetscFunctionBegin;
7722764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
7732764a2aaSMatthew G. Knepley   if (disc) PetscValidPointer(disc, 3);
7742764a2aaSMatthew G. Knepley   if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f);
7752764a2aaSMatthew G. Knepley   ierr = PetscDSEnlarge_Static(prob, f+1);CHKERRQ(ierr);
7762764a2aaSMatthew G. Knepley   if (prob->discBd[f]) {ierr = PetscObjectDereference(prob->discBd[f]);CHKERRQ(ierr);}
7772764a2aaSMatthew G. Knepley   prob->discBd[f] = disc;
7782764a2aaSMatthew G. Knepley   ierr = PetscObjectReference(disc);CHKERRQ(ierr);
7792764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
7802764a2aaSMatthew G. Knepley }
7812764a2aaSMatthew G. Knepley 
7822764a2aaSMatthew G. Knepley #undef __FUNCT__
7832764a2aaSMatthew G. Knepley #define __FUNCT__ "PetscDSAddDiscretization"
784bc4ae4beSMatthew G. Knepley /*@
785bc4ae4beSMatthew G. Knepley   PetscDSAddDiscretization - Adds a discretization object
786bc4ae4beSMatthew G. Knepley 
787bc4ae4beSMatthew G. Knepley   Not collective
788bc4ae4beSMatthew G. Knepley 
789bc4ae4beSMatthew G. Knepley   Input Parameters:
790bc4ae4beSMatthew G. Knepley + prob - The PetscDS object
791bc4ae4beSMatthew G. Knepley - disc - The boundary discretization object
792bc4ae4beSMatthew G. Knepley 
793bc4ae4beSMatthew G. Knepley   Level: beginner
794bc4ae4beSMatthew G. Knepley 
795bc4ae4beSMatthew G. Knepley .seealso: PetscDSGetDiscretization(), PetscDSSetDiscretization(), PetscDSGetNumFields(), PetscDSCreate()
796bc4ae4beSMatthew G. Knepley @*/
7972764a2aaSMatthew G. Knepley PetscErrorCode PetscDSAddDiscretization(PetscDS prob, PetscObject disc)
7982764a2aaSMatthew G. Knepley {
7992764a2aaSMatthew G. Knepley   PetscErrorCode ierr;
8002764a2aaSMatthew G. Knepley 
8012764a2aaSMatthew G. Knepley   PetscFunctionBegin;
8022764a2aaSMatthew G. Knepley   ierr = PetscDSSetDiscretization(prob, prob->Nf, disc);CHKERRQ(ierr);
8032764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
8042764a2aaSMatthew G. Knepley }
8052764a2aaSMatthew G. Knepley 
8062764a2aaSMatthew G. Knepley #undef __FUNCT__
8072764a2aaSMatthew G. Knepley #define __FUNCT__ "PetscDSAddBdDiscretization"
808bc4ae4beSMatthew G. Knepley /*@
809bc4ae4beSMatthew G. Knepley   PetscDSAddBdDiscretization - Adds a boundary discretization object
810bc4ae4beSMatthew G. Knepley 
811bc4ae4beSMatthew G. Knepley   Not collective
812bc4ae4beSMatthew G. Knepley 
813bc4ae4beSMatthew G. Knepley   Input Parameters:
814bc4ae4beSMatthew G. Knepley + prob - The PetscDS object
815bc4ae4beSMatthew G. Knepley - disc - The boundary discretization object
816bc4ae4beSMatthew G. Knepley 
817bc4ae4beSMatthew G. Knepley   Level: beginner
818bc4ae4beSMatthew G. Knepley 
819bc4ae4beSMatthew G. Knepley .seealso: PetscDSGetBdDiscretization(), PetscDSSetBdDiscretization(), PetscDSGetNumFields(), PetscDSCreate()
820bc4ae4beSMatthew G. Knepley @*/
8212764a2aaSMatthew G. Knepley PetscErrorCode PetscDSAddBdDiscretization(PetscDS prob, PetscObject disc)
8222764a2aaSMatthew G. Knepley {
8232764a2aaSMatthew G. Knepley   PetscErrorCode ierr;
8242764a2aaSMatthew G. Knepley 
8252764a2aaSMatthew G. Knepley   PetscFunctionBegin;
8262764a2aaSMatthew G. Knepley   ierr = PetscDSSetBdDiscretization(prob, prob->Nf, disc);CHKERRQ(ierr);
8272764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
8282764a2aaSMatthew G. Knepley }
8292764a2aaSMatthew G. Knepley 
8302764a2aaSMatthew G. Knepley #undef __FUNCT__
831249df284SMatthew G. Knepley #define __FUNCT__ "PetscDSGetImplicit"
832249df284SMatthew G. Knepley /*@
833249df284SMatthew G. Knepley   PetscDSGetImplicit - Returns the flag for implicit solve for this field. This is just a guide for IMEX
834249df284SMatthew G. Knepley 
835249df284SMatthew G. Knepley   Not collective
836249df284SMatthew G. Knepley 
837249df284SMatthew G. Knepley   Input Parameters:
838249df284SMatthew G. Knepley + prob - The PetscDS object
839249df284SMatthew G. Knepley - f - The field number
840249df284SMatthew G. Knepley 
841249df284SMatthew G. Knepley   Output Parameter:
842249df284SMatthew G. Knepley . implicit - The flag indicating what kind of solve to use for this field
843249df284SMatthew G. Knepley 
844249df284SMatthew G. Knepley   Level: developer
845249df284SMatthew G. Knepley 
846249df284SMatthew G. Knepley .seealso: PetscDSSetImplicit(), PetscDSSetDiscretization(), PetscDSAddDiscretization(), PetscDSGetBdDiscretization(), PetscDSGetNumFields(), PetscDSCreate()
847249df284SMatthew G. Knepley @*/
848249df284SMatthew G. Knepley PetscErrorCode PetscDSGetImplicit(PetscDS prob, PetscInt f, PetscBool *implicit)
849249df284SMatthew G. Knepley {
850249df284SMatthew G. Knepley   PetscFunctionBegin;
851249df284SMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
852249df284SMatthew G. Knepley   PetscValidPointer(implicit, 3);
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   *implicit = prob->implicit[f];
855249df284SMatthew G. Knepley   PetscFunctionReturn(0);
856249df284SMatthew G. Knepley }
857249df284SMatthew G. Knepley 
858249df284SMatthew G. Knepley #undef __FUNCT__
859249df284SMatthew G. Knepley #define __FUNCT__ "PetscDSSetImplicit"
860249df284SMatthew G. Knepley /*@
861249df284SMatthew G. Knepley   PetscDSSetImplicit - Set the flag for implicit solve for this field. This is just a guide for IMEX
862249df284SMatthew G. Knepley 
863249df284SMatthew G. Knepley   Not collective
864249df284SMatthew G. Knepley 
865249df284SMatthew G. Knepley   Input Parameters:
866249df284SMatthew G. Knepley + prob - The PetscDS object
867249df284SMatthew G. Knepley . f - The field number
868249df284SMatthew G. Knepley - implicit - The flag indicating what kind of solve to use for this field
869249df284SMatthew G. Knepley 
870249df284SMatthew G. Knepley   Level: developer
871249df284SMatthew G. Knepley 
872249df284SMatthew G. Knepley .seealso: PetscDSGetImplicit(), PetscDSSetDiscretization(), PetscDSAddDiscretization(), PetscDSGetBdDiscretization(), PetscDSGetNumFields(), PetscDSCreate()
873249df284SMatthew G. Knepley @*/
874249df284SMatthew G. Knepley PetscErrorCode PetscDSSetImplicit(PetscDS prob, PetscInt f, PetscBool implicit)
875249df284SMatthew G. Knepley {
876249df284SMatthew G. Knepley   PetscFunctionBegin;
877249df284SMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
878249df284SMatthew 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);
879249df284SMatthew G. Knepley   prob->implicit[f] = implicit;
880249df284SMatthew G. Knepley   PetscFunctionReturn(0);
881249df284SMatthew G. Knepley }
882249df284SMatthew G. Knepley 
883249df284SMatthew G. Knepley #undef __FUNCT__
884a6cbbb48SMatthew G. Knepley #define __FUNCT__ "PetscDSGetAdjacency"
885a6cbbb48SMatthew G. Knepley /*@
886a6cbbb48SMatthew G. Knepley   PetscDSGetAdjacency - Returns the flags for determining variable influence
887a6cbbb48SMatthew G. Knepley 
888a6cbbb48SMatthew G. Knepley   Not collective
889a6cbbb48SMatthew G. Knepley 
890a6cbbb48SMatthew G. Knepley   Input Parameters:
891a6cbbb48SMatthew G. Knepley + prob - The PetscDS object
892a6cbbb48SMatthew G. Knepley - f - The field number
893a6cbbb48SMatthew G. Knepley 
894a6cbbb48SMatthew G. Knepley   Output Parameter:
895a6cbbb48SMatthew G. Knepley + useCone    - Flag for variable influence starting with the cone operation
896a6cbbb48SMatthew G. Knepley - useClosure - Flag for variable influence using transitive closure
897a6cbbb48SMatthew G. Knepley 
898a6cbbb48SMatthew G. Knepley   Note: See the discussion in DMPlexGetAdjacencyUseCone() and DMPlexGetAdjacencyUseClosure()
899a6cbbb48SMatthew G. Knepley 
900a6cbbb48SMatthew G. Knepley   Level: developer
901a6cbbb48SMatthew G. Knepley 
902a6cbbb48SMatthew G. Knepley .seealso: PetscDSSetAdjacency(), DMPlexGetAdjacencyUseCone(), DMPlexGetAdjacencyUseClosure(), PetscDSSetDiscretization(), PetscDSAddDiscretization(), PetscDSGetBdDiscretization(), PetscDSGetNumFields(), PetscDSCreate()
903a6cbbb48SMatthew G. Knepley @*/
904a6cbbb48SMatthew G. Knepley PetscErrorCode PetscDSGetAdjacency(PetscDS prob, PetscInt f, PetscBool *useCone, PetscBool *useClosure)
905a6cbbb48SMatthew G. Knepley {
906a6cbbb48SMatthew G. Knepley   PetscFunctionBegin;
907a6cbbb48SMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
908a6cbbb48SMatthew G. Knepley   PetscValidPointer(useCone, 3);
909a6cbbb48SMatthew G. Knepley   PetscValidPointer(useClosure, 4);
910a6cbbb48SMatthew 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);
911a6cbbb48SMatthew G. Knepley   *useCone    = prob->adjacency[f*2+0];
912a6cbbb48SMatthew G. Knepley   *useClosure = prob->adjacency[f*2+1];
913a6cbbb48SMatthew G. Knepley   PetscFunctionReturn(0);
914a6cbbb48SMatthew G. Knepley }
915a6cbbb48SMatthew G. Knepley 
916a6cbbb48SMatthew G. Knepley #undef __FUNCT__
917a6cbbb48SMatthew G. Knepley #define __FUNCT__ "PetscDSSetAdjacency"
918a6cbbb48SMatthew G. Knepley /*@
919a6cbbb48SMatthew G. Knepley   PetscDSSetAdjacency - Set the flags for determining variable influence
920a6cbbb48SMatthew G. Knepley 
921a6cbbb48SMatthew G. Knepley   Not collective
922a6cbbb48SMatthew G. Knepley 
923a6cbbb48SMatthew G. Knepley   Input Parameters:
924a6cbbb48SMatthew G. Knepley + prob - The PetscDS object
925a6cbbb48SMatthew G. Knepley . f - The field number
926a6cbbb48SMatthew G. Knepley . useCone    - Flag for variable influence starting with the cone operation
927a6cbbb48SMatthew G. Knepley - useClosure - Flag for variable influence using transitive closure
928a6cbbb48SMatthew G. Knepley 
929a6cbbb48SMatthew G. Knepley   Note: See the discussion in DMPlexGetAdjacencyUseCone() and DMPlexGetAdjacencyUseClosure()
930a6cbbb48SMatthew G. Knepley 
931a6cbbb48SMatthew G. Knepley   Level: developer
932a6cbbb48SMatthew G. Knepley 
933a6cbbb48SMatthew G. Knepley .seealso: PetscDSGetAdjacency(), DMPlexGetAdjacencyUseCone(), DMPlexGetAdjacencyUseClosure(), PetscDSSetDiscretization(), PetscDSAddDiscretization(), PetscDSGetBdDiscretization(), PetscDSGetNumFields(), PetscDSCreate()
934a6cbbb48SMatthew G. Knepley @*/
935a6cbbb48SMatthew G. Knepley PetscErrorCode PetscDSSetAdjacency(PetscDS prob, PetscInt f, PetscBool useCone, PetscBool useClosure)
936a6cbbb48SMatthew G. Knepley {
937a6cbbb48SMatthew G. Knepley   PetscFunctionBegin;
938a6cbbb48SMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
939a6cbbb48SMatthew 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);
940a6cbbb48SMatthew G. Knepley   prob->adjacency[f*2+0] = useCone;
941a6cbbb48SMatthew G. Knepley   prob->adjacency[f*2+1] = useClosure;
942a6cbbb48SMatthew G. Knepley   PetscFunctionReturn(0);
943a6cbbb48SMatthew G. Knepley }
944a6cbbb48SMatthew G. Knepley 
945a6cbbb48SMatthew G. Knepley #undef __FUNCT__
9462764a2aaSMatthew G. Knepley #define __FUNCT__ "PetscDSGetObjective"
9472764a2aaSMatthew G. Knepley PetscErrorCode PetscDSGetObjective(PetscDS prob, PetscInt f,
9482764a2aaSMatthew G. Knepley                                         void (**obj)(const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], const PetscReal x[], PetscScalar obj[]))
9492764a2aaSMatthew G. Knepley {
9502764a2aaSMatthew G. Knepley   PetscFunctionBegin;
9512764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
9522764a2aaSMatthew G. Knepley   PetscValidPointer(obj, 2);
9532764a2aaSMatthew 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);
9542764a2aaSMatthew G. Knepley   *obj = prob->obj[f];
9552764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
9562764a2aaSMatthew G. Knepley }
9572764a2aaSMatthew G. Knepley 
9582764a2aaSMatthew G. Knepley #undef __FUNCT__
9592764a2aaSMatthew G. Knepley #define __FUNCT__ "PetscDSSetObjective"
9602764a2aaSMatthew G. Knepley PetscErrorCode PetscDSSetObjective(PetscDS prob, PetscInt f,
9612764a2aaSMatthew G. Knepley                                         void (*obj)(const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], const PetscReal x[], PetscScalar obj[]))
9622764a2aaSMatthew G. Knepley {
9632764a2aaSMatthew G. Knepley   PetscErrorCode ierr;
9642764a2aaSMatthew G. Knepley 
9652764a2aaSMatthew G. Knepley   PetscFunctionBegin;
9662764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
9672764a2aaSMatthew G. Knepley   PetscValidFunction(obj, 2);
9682764a2aaSMatthew G. Knepley   if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f);
9692764a2aaSMatthew G. Knepley   ierr = PetscDSEnlarge_Static(prob, f+1);CHKERRQ(ierr);
9702764a2aaSMatthew G. Knepley   prob->obj[f] = obj;
9712764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
9722764a2aaSMatthew G. Knepley }
9732764a2aaSMatthew G. Knepley 
9742764a2aaSMatthew G. Knepley #undef __FUNCT__
9752764a2aaSMatthew G. Knepley #define __FUNCT__ "PetscDSGetResidual"
9762764a2aaSMatthew G. Knepley PetscErrorCode PetscDSGetResidual(PetscDS prob, PetscInt f,
9772764a2aaSMatthew G. Knepley                                        void (**f0)(const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], const PetscReal x[], PetscScalar f0[]),
9782764a2aaSMatthew G. Knepley                                        void (**f1)(const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], const PetscReal x[], PetscScalar f1[]))
9792764a2aaSMatthew G. Knepley {
9802764a2aaSMatthew G. Knepley   PetscFunctionBegin;
9812764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
9822764a2aaSMatthew 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);
9832764a2aaSMatthew G. Knepley   if (f0) {PetscValidPointer(f0, 3); *f0 = prob->f[f*2+0];}
9842764a2aaSMatthew G. Knepley   if (f1) {PetscValidPointer(f1, 4); *f1 = prob->f[f*2+1];}
9852764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
9862764a2aaSMatthew G. Knepley }
9872764a2aaSMatthew G. Knepley 
9882764a2aaSMatthew G. Knepley #undef __FUNCT__
9892764a2aaSMatthew G. Knepley #define __FUNCT__ "PetscDSSetResidual"
9902764a2aaSMatthew G. Knepley PetscErrorCode PetscDSSetResidual(PetscDS prob, PetscInt f,
9912764a2aaSMatthew G. Knepley                                   void (*f0)(const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], const PetscReal x[], PetscScalar f0[]),
9922764a2aaSMatthew G. Knepley                                   void (*f1)(const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], const PetscReal x[], PetscScalar f1[]))
9932764a2aaSMatthew G. Knepley {
9942764a2aaSMatthew G. Knepley   PetscErrorCode ierr;
9952764a2aaSMatthew G. Knepley 
9962764a2aaSMatthew G. Knepley   PetscFunctionBegin;
9972764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
9982764a2aaSMatthew G. Knepley   PetscValidFunction(f0, 3);
9992764a2aaSMatthew G. Knepley   PetscValidFunction(f1, 4);
10002764a2aaSMatthew G. Knepley   if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f);
10012764a2aaSMatthew G. Knepley   ierr = PetscDSEnlarge_Static(prob, f+1);CHKERRQ(ierr);
10022764a2aaSMatthew G. Knepley   prob->f[f*2+0] = f0;
10032764a2aaSMatthew G. Knepley   prob->f[f*2+1] = f1;
10042764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
10052764a2aaSMatthew G. Knepley }
10062764a2aaSMatthew G. Knepley 
10072764a2aaSMatthew G. Knepley #undef __FUNCT__
10082764a2aaSMatthew G. Knepley #define __FUNCT__ "PetscDSGetJacobian"
10092764a2aaSMatthew G. Knepley PetscErrorCode PetscDSGetJacobian(PetscDS prob, PetscInt f, PetscInt g,
10102764a2aaSMatthew G. Knepley                                        void (**g0)(const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], const PetscReal x[], PetscScalar g0[]),
10112764a2aaSMatthew G. Knepley                                        void (**g1)(const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], const PetscReal x[], PetscScalar g1[]),
10122764a2aaSMatthew G. Knepley                                        void (**g2)(const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], const PetscReal x[], PetscScalar g2[]),
10132764a2aaSMatthew G. Knepley                                        void (**g3)(const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], const PetscReal x[], PetscScalar g3[]))
10142764a2aaSMatthew G. Knepley {
10152764a2aaSMatthew G. Knepley   PetscFunctionBegin;
10162764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
10172764a2aaSMatthew 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);
10182764a2aaSMatthew 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);
10192764a2aaSMatthew G. Knepley   if (g0) {PetscValidPointer(g0, 4); *g0 = prob->g[(f*prob->Nf + g)*4+0];}
10202764a2aaSMatthew G. Knepley   if (g1) {PetscValidPointer(g1, 5); *g1 = prob->g[(f*prob->Nf + g)*4+1];}
10212764a2aaSMatthew G. Knepley   if (g2) {PetscValidPointer(g2, 6); *g2 = prob->g[(f*prob->Nf + g)*4+2];}
10222764a2aaSMatthew G. Knepley   if (g3) {PetscValidPointer(g3, 7); *g3 = prob->g[(f*prob->Nf + g)*4+3];}
10232764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
10242764a2aaSMatthew G. Knepley }
10252764a2aaSMatthew G. Knepley 
10262764a2aaSMatthew G. Knepley #undef __FUNCT__
10272764a2aaSMatthew G. Knepley #define __FUNCT__ "PetscDSSetJacobian"
10282764a2aaSMatthew G. Knepley PetscErrorCode PetscDSSetJacobian(PetscDS prob, PetscInt f, PetscInt g,
10292764a2aaSMatthew G. Knepley                                        void (*g0)(const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], const PetscReal x[], PetscScalar g0[]),
10302764a2aaSMatthew G. Knepley                                        void (*g1)(const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], const PetscReal x[], PetscScalar g1[]),
10312764a2aaSMatthew G. Knepley                                        void (*g2)(const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], const PetscReal x[], PetscScalar g2[]),
10322764a2aaSMatthew G. Knepley                                        void (*g3)(const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], const PetscReal x[], PetscScalar g3[]))
10332764a2aaSMatthew G. Knepley {
10342764a2aaSMatthew G. Knepley   PetscErrorCode ierr;
10352764a2aaSMatthew G. Knepley 
10362764a2aaSMatthew G. Knepley   PetscFunctionBegin;
10372764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
10382764a2aaSMatthew G. Knepley   if (g0) PetscValidFunction(g0, 4);
10392764a2aaSMatthew G. Knepley   if (g1) PetscValidFunction(g1, 5);
10402764a2aaSMatthew G. Knepley   if (g2) PetscValidFunction(g2, 6);
10412764a2aaSMatthew G. Knepley   if (g3) PetscValidFunction(g3, 7);
10422764a2aaSMatthew G. Knepley   if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f);
10432764a2aaSMatthew G. Knepley   if (g < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", g);
10442764a2aaSMatthew G. Knepley   ierr = PetscDSEnlarge_Static(prob, PetscMax(f, g)+1);CHKERRQ(ierr);
10452764a2aaSMatthew G. Knepley   prob->g[(f*prob->Nf + g)*4+0] = g0;
10462764a2aaSMatthew G. Knepley   prob->g[(f*prob->Nf + g)*4+1] = g1;
10472764a2aaSMatthew G. Knepley   prob->g[(f*prob->Nf + g)*4+2] = g2;
10482764a2aaSMatthew G. Knepley   prob->g[(f*prob->Nf + g)*4+3] = g3;
10492764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
10502764a2aaSMatthew G. Knepley }
10512764a2aaSMatthew G. Knepley 
10522764a2aaSMatthew G. Knepley #undef __FUNCT__
10530c2f2876SMatthew G. Knepley #define __FUNCT__ "PetscDSGetRiemannSolver"
10540c2f2876SMatthew G. Knepley /*@C
10550c2f2876SMatthew G. Knepley   PetscDSGetRiemannSolver - Returns the Riemann solver for the given field
10560c2f2876SMatthew G. Knepley 
10570c2f2876SMatthew G. Knepley   Not collective
10580c2f2876SMatthew G. Knepley 
10590c2f2876SMatthew G. Knepley   Input Arguments:
10600c2f2876SMatthew G. Knepley + prob - The PetscDS object
10610c2f2876SMatthew G. Knepley - f    - The field number
10620c2f2876SMatthew G. Knepley 
10630c2f2876SMatthew G. Knepley   Output Argument:
10640c2f2876SMatthew G. Knepley . r    - Riemann solver
10650c2f2876SMatthew G. Knepley 
10660c2f2876SMatthew G. Knepley   Calling sequence for r:
10670c2f2876SMatthew G. Knepley 
10680c2f2876SMatthew G. Knepley $ r(const PetscReal x[], const PetscReal n[], const PetscScalar uL[], const PetscScalar uR[], PetscScalar flux[], void *ctx)
10690c2f2876SMatthew G. Knepley 
10700c2f2876SMatthew G. Knepley + x    - The coordinates at a point on the interface
10710c2f2876SMatthew G. Knepley . n    - The normal vector to the interface
10720c2f2876SMatthew G. Knepley . uL   - The state vector to the left of the interface
10730c2f2876SMatthew G. Knepley . uR   - The state vector to the right of the interface
10740c2f2876SMatthew G. Knepley . flux - output array of flux through the interface
10750c2f2876SMatthew G. Knepley - ctx  - optional user context
10760c2f2876SMatthew G. Knepley 
10770c2f2876SMatthew G. Knepley   Level: intermediate
10780c2f2876SMatthew G. Knepley 
10790c2f2876SMatthew G. Knepley .seealso: PetscDSSetRiemannSolver()
10800c2f2876SMatthew G. Knepley @*/
10810c2f2876SMatthew G. Knepley PetscErrorCode PetscDSGetRiemannSolver(PetscDS prob, PetscInt f,
10820c2f2876SMatthew G. Knepley                                        void (**r)(const PetscReal x[], const PetscReal n[], const PetscScalar uL[], const PetscScalar uR[], PetscScalar flux[], void *ctx))
10830c2f2876SMatthew G. Knepley {
10840c2f2876SMatthew G. Knepley   PetscFunctionBegin;
10850c2f2876SMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
10860c2f2876SMatthew 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);
10870c2f2876SMatthew G. Knepley   PetscValidPointer(r, 3);
10880c2f2876SMatthew G. Knepley   *r = prob->r[f];
10890c2f2876SMatthew G. Knepley   PetscFunctionReturn(0);
10900c2f2876SMatthew G. Knepley }
10910c2f2876SMatthew G. Knepley 
10920c2f2876SMatthew G. Knepley #undef __FUNCT__
10930c2f2876SMatthew G. Knepley #define __FUNCT__ "PetscDSSetRiemannSolver"
10940c2f2876SMatthew G. Knepley /*@C
10950c2f2876SMatthew G. Knepley   PetscDSSetRiemannSolver - Sets the Riemann solver for the given field
10960c2f2876SMatthew G. Knepley 
10970c2f2876SMatthew G. Knepley   Not collective
10980c2f2876SMatthew G. Knepley 
10990c2f2876SMatthew G. Knepley   Input Arguments:
11000c2f2876SMatthew G. Knepley + prob - The PetscDS object
11010c2f2876SMatthew G. Knepley . f    - The field number
11020c2f2876SMatthew G. Knepley - r    - Riemann solver
11030c2f2876SMatthew G. Knepley 
11040c2f2876SMatthew G. Knepley   Calling sequence for r:
11050c2f2876SMatthew G. Knepley 
11060c2f2876SMatthew G. Knepley $ r(const PetscReal x[], const PetscReal n[], const PetscScalar uL[], const PetscScalar uR[], PetscScalar flux[], void *ctx)
11070c2f2876SMatthew G. Knepley 
11080c2f2876SMatthew G. Knepley + x    - The coordinates at a point on the interface
11090c2f2876SMatthew G. Knepley . n    - The normal vector to the interface
11100c2f2876SMatthew G. Knepley . uL   - The state vector to the left of the interface
11110c2f2876SMatthew G. Knepley . uR   - The state vector to the right of the interface
11120c2f2876SMatthew G. Knepley . flux - output array of flux through the interface
11130c2f2876SMatthew G. Knepley - ctx  - optional user context
11140c2f2876SMatthew G. Knepley 
11150c2f2876SMatthew G. Knepley   Level: intermediate
11160c2f2876SMatthew G. Knepley 
11170c2f2876SMatthew G. Knepley .seealso: PetscDSGetRiemannSolver()
11180c2f2876SMatthew G. Knepley @*/
11190c2f2876SMatthew G. Knepley PetscErrorCode PetscDSSetRiemannSolver(PetscDS prob, PetscInt f,
11200c2f2876SMatthew G. Knepley                                        void (*r)(const PetscReal x[], const PetscReal n[], const PetscScalar uL[], const PetscScalar uR[], PetscScalar flux[], void *ctx))
11210c2f2876SMatthew G. Knepley {
11220c2f2876SMatthew G. Knepley   PetscErrorCode ierr;
11230c2f2876SMatthew G. Knepley 
11240c2f2876SMatthew G. Knepley   PetscFunctionBegin;
11250c2f2876SMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
11260c2f2876SMatthew G. Knepley   PetscValidFunction(r, 3);
11270c2f2876SMatthew G. Knepley   if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f);
11280c2f2876SMatthew G. Knepley   ierr = PetscDSEnlarge_Static(prob, f+1);CHKERRQ(ierr);
11290c2f2876SMatthew G. Knepley   prob->r[f] = r;
11300c2f2876SMatthew G. Knepley   PetscFunctionReturn(0);
11310c2f2876SMatthew G. Knepley }
11320c2f2876SMatthew G. Knepley 
11330c2f2876SMatthew G. Knepley #undef __FUNCT__
11340c2f2876SMatthew G. Knepley #define __FUNCT__ "PetscDSGetContext"
11350c2f2876SMatthew G. Knepley PetscErrorCode PetscDSGetContext(PetscDS prob, PetscInt f, void **ctx)
11360c2f2876SMatthew G. Knepley {
11370c2f2876SMatthew G. Knepley   PetscFunctionBegin;
11380c2f2876SMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
11390c2f2876SMatthew 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);
11400c2f2876SMatthew G. Knepley   PetscValidPointer(ctx, 3);
11410c2f2876SMatthew G. Knepley   *ctx = prob->ctx[f];
11420c2f2876SMatthew G. Knepley   PetscFunctionReturn(0);
11430c2f2876SMatthew G. Knepley }
11440c2f2876SMatthew G. Knepley 
11450c2f2876SMatthew G. Knepley #undef __FUNCT__
11460c2f2876SMatthew G. Knepley #define __FUNCT__ "PetscDSSetContext"
11470c2f2876SMatthew G. Knepley PetscErrorCode PetscDSSetContext(PetscDS prob, PetscInt f, void *ctx)
11480c2f2876SMatthew G. Knepley {
11490c2f2876SMatthew G. Knepley   PetscErrorCode ierr;
11500c2f2876SMatthew G. Knepley 
11510c2f2876SMatthew G. Knepley   PetscFunctionBegin;
11520c2f2876SMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
11530c2f2876SMatthew G. Knepley   if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f);
11540c2f2876SMatthew G. Knepley   ierr = PetscDSEnlarge_Static(prob, f+1);CHKERRQ(ierr);
11550c2f2876SMatthew G. Knepley   prob->ctx[f] = ctx;
11560c2f2876SMatthew G. Knepley   PetscFunctionReturn(0);
11570c2f2876SMatthew G. Knepley }
11580c2f2876SMatthew G. Knepley 
11590c2f2876SMatthew G. Knepley #undef __FUNCT__
11602764a2aaSMatthew G. Knepley #define __FUNCT__ "PetscDSGetBdResidual"
11612764a2aaSMatthew G. Knepley PetscErrorCode PetscDSGetBdResidual(PetscDS prob, PetscInt f,
11622764a2aaSMatthew G. Knepley                                          void (**f0)(const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], const PetscReal x[], const PetscReal n[], PetscScalar f0[]),
11632764a2aaSMatthew G. Knepley                                          void (**f1)(const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], const PetscReal x[], const PetscReal n[], PetscScalar f1[]))
11642764a2aaSMatthew G. Knepley {
11652764a2aaSMatthew G. Knepley   PetscFunctionBegin;
11662764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
11672764a2aaSMatthew 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);
11682764a2aaSMatthew G. Knepley   if (f0) {PetscValidPointer(f0, 3); *f0 = prob->fBd[f*2+0];}
11692764a2aaSMatthew G. Knepley   if (f1) {PetscValidPointer(f1, 4); *f1 = prob->fBd[f*2+1];}
11702764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
11712764a2aaSMatthew G. Knepley }
11722764a2aaSMatthew G. Knepley 
11732764a2aaSMatthew G. Knepley #undef __FUNCT__
11742764a2aaSMatthew G. Knepley #define __FUNCT__ "PetscDSSetBdResidual"
11752764a2aaSMatthew G. Knepley PetscErrorCode PetscDSSetBdResidual(PetscDS prob, PetscInt f,
11762764a2aaSMatthew G. Knepley                                          void (*f0)(const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], const PetscReal x[], const PetscReal n[], PetscScalar f0[]),
11772764a2aaSMatthew G. Knepley                                          void (*f1)(const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], const PetscReal x[], const PetscReal n[], PetscScalar f1[]))
11782764a2aaSMatthew G. Knepley {
11792764a2aaSMatthew G. Knepley   PetscErrorCode ierr;
11802764a2aaSMatthew G. Knepley 
11812764a2aaSMatthew G. Knepley   PetscFunctionBegin;
11822764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
11832764a2aaSMatthew G. Knepley   if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f);
11842764a2aaSMatthew G. Knepley   ierr = PetscDSEnlarge_Static(prob, f+1);CHKERRQ(ierr);
11852764a2aaSMatthew G. Knepley   if (f0) {PetscValidFunction(f0, 3); prob->fBd[f*2+0] = f0;}
11862764a2aaSMatthew G. Knepley   if (f1) {PetscValidFunction(f1, 4); prob->fBd[f*2+1] = f1;}
11872764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
11882764a2aaSMatthew G. Knepley }
11892764a2aaSMatthew G. Knepley 
11902764a2aaSMatthew G. Knepley #undef __FUNCT__
11912764a2aaSMatthew G. Knepley #define __FUNCT__ "PetscDSGetBdJacobian"
11922764a2aaSMatthew G. Knepley PetscErrorCode PetscDSGetBdJacobian(PetscDS prob, PetscInt f, PetscInt g,
11932764a2aaSMatthew G. Knepley                                          void (**g0)(const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], const PetscReal x[], const PetscReal n[], PetscScalar g0[]),
11942764a2aaSMatthew G. Knepley                                          void (**g1)(const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], const PetscReal x[], const PetscReal n[], PetscScalar g1[]),
11952764a2aaSMatthew G. Knepley                                          void (**g2)(const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], const PetscReal x[], const PetscReal n[], PetscScalar g2[]),
11962764a2aaSMatthew G. Knepley                                          void (**g3)(const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], const PetscReal x[], const PetscReal n[], PetscScalar g3[]))
11972764a2aaSMatthew G. Knepley {
11982764a2aaSMatthew G. Knepley   PetscFunctionBegin;
11992764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
12002764a2aaSMatthew 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);
12012764a2aaSMatthew 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);
12022764a2aaSMatthew G. Knepley   if (g0) {PetscValidPointer(g0, 4); *g0 = prob->gBd[(f*prob->Nf + g)*4+0];}
12032764a2aaSMatthew G. Knepley   if (g1) {PetscValidPointer(g1, 5); *g1 = prob->gBd[(f*prob->Nf + g)*4+1];}
12042764a2aaSMatthew G. Knepley   if (g2) {PetscValidPointer(g2, 6); *g2 = prob->gBd[(f*prob->Nf + g)*4+2];}
12052764a2aaSMatthew G. Knepley   if (g3) {PetscValidPointer(g3, 7); *g3 = prob->gBd[(f*prob->Nf + g)*4+3];}
12062764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
12072764a2aaSMatthew G. Knepley }
12082764a2aaSMatthew G. Knepley 
12092764a2aaSMatthew G. Knepley #undef __FUNCT__
12102764a2aaSMatthew G. Knepley #define __FUNCT__ "PetscDSSetBdJacobian"
12112764a2aaSMatthew G. Knepley PetscErrorCode PetscDSSetBdJacobian(PetscDS prob, PetscInt f, PetscInt g,
12122764a2aaSMatthew G. Knepley                                          void (*g0)(const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], const PetscReal x[], const PetscReal n[], PetscScalar g0[]),
12132764a2aaSMatthew G. Knepley                                          void (*g1)(const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], const PetscReal x[], const PetscReal n[], PetscScalar g1[]),
12142764a2aaSMatthew G. Knepley                                          void (*g2)(const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], const PetscReal x[], const PetscReal n[], PetscScalar g2[]),
12152764a2aaSMatthew G. Knepley                                          void (*g3)(const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], const PetscReal x[], const PetscReal n[], PetscScalar g3[]))
12162764a2aaSMatthew G. Knepley {
12172764a2aaSMatthew G. Knepley   PetscErrorCode ierr;
12182764a2aaSMatthew G. Knepley 
12192764a2aaSMatthew G. Knepley   PetscFunctionBegin;
12202764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
12212764a2aaSMatthew G. Knepley   if (g0) PetscValidFunction(g0, 4);
12222764a2aaSMatthew G. Knepley   if (g1) PetscValidFunction(g1, 5);
12232764a2aaSMatthew G. Knepley   if (g2) PetscValidFunction(g2, 6);
12242764a2aaSMatthew G. Knepley   if (g3) PetscValidFunction(g3, 7);
12252764a2aaSMatthew G. Knepley   if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f);
12262764a2aaSMatthew G. Knepley   if (g < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", g);
12272764a2aaSMatthew G. Knepley   ierr = PetscDSEnlarge_Static(prob, PetscMax(f, g)+1);CHKERRQ(ierr);
12282764a2aaSMatthew G. Knepley   prob->gBd[(f*prob->Nf + g)*4+0] = g0;
12292764a2aaSMatthew G. Knepley   prob->gBd[(f*prob->Nf + g)*4+1] = g1;
12302764a2aaSMatthew G. Knepley   prob->gBd[(f*prob->Nf + g)*4+2] = g2;
12312764a2aaSMatthew G. Knepley   prob->gBd[(f*prob->Nf + g)*4+3] = g3;
12322764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
12332764a2aaSMatthew G. Knepley }
12342764a2aaSMatthew G. Knepley 
12352764a2aaSMatthew G. Knepley #undef __FUNCT__
12362764a2aaSMatthew G. Knepley #define __FUNCT__ "PetscDSGetFieldOffset"
1237bc4ae4beSMatthew G. Knepley /*@
1238bc4ae4beSMatthew G. Knepley   PetscDSGetFieldOffset - Returns the offset of the given field in the full space basis
1239bc4ae4beSMatthew G. Knepley 
1240bc4ae4beSMatthew G. Knepley   Not collective
1241bc4ae4beSMatthew G. Knepley 
1242bc4ae4beSMatthew G. Knepley   Input Parameters:
1243bc4ae4beSMatthew G. Knepley + prob - The PetscDS object
1244bc4ae4beSMatthew G. Knepley - f - The field number
1245bc4ae4beSMatthew G. Knepley 
1246bc4ae4beSMatthew G. Knepley   Output Parameter:
1247bc4ae4beSMatthew G. Knepley . off - The offset
1248bc4ae4beSMatthew G. Knepley 
1249bc4ae4beSMatthew G. Knepley   Level: beginner
1250bc4ae4beSMatthew G. Knepley 
1251bc4ae4beSMatthew G. Knepley .seealso: PetscDSGetBdFieldOffset(), PetscDSGetNumFields(), PetscDSCreate()
1252bc4ae4beSMatthew G. Knepley @*/
12532764a2aaSMatthew G. Knepley PetscErrorCode PetscDSGetFieldOffset(PetscDS prob, PetscInt f, PetscInt *off)
12542764a2aaSMatthew G. Knepley {
12552764a2aaSMatthew G. Knepley   PetscInt       g;
12562764a2aaSMatthew G. Knepley   PetscErrorCode ierr;
12572764a2aaSMatthew G. Knepley 
12582764a2aaSMatthew G. Knepley   PetscFunctionBegin;
12592764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
12602764a2aaSMatthew G. Knepley   PetscValidPointer(off, 3);
12612764a2aaSMatthew 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);
12622764a2aaSMatthew G. Knepley   *off = 0;
12632764a2aaSMatthew G. Knepley   for (g = 0; g < f; ++g) {
12642764a2aaSMatthew G. Knepley     PetscFE  fe = (PetscFE) prob->disc[g];
12652764a2aaSMatthew G. Knepley     PetscInt Nb, Nc;
12662764a2aaSMatthew G. Knepley 
12672764a2aaSMatthew G. Knepley     ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr);
12682764a2aaSMatthew G. Knepley     ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
12692764a2aaSMatthew G. Knepley     *off += Nb*Nc;
12702764a2aaSMatthew G. Knepley   }
12712764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
12722764a2aaSMatthew G. Knepley }
12732764a2aaSMatthew G. Knepley 
12742764a2aaSMatthew G. Knepley #undef __FUNCT__
12752764a2aaSMatthew G. Knepley #define __FUNCT__ "PetscDSGetBdFieldOffset"
1276bc4ae4beSMatthew G. Knepley /*@
1277c3ac4435SMatthew G. Knepley   PetscDSGetBdFieldOffset - Returns the offset of the given field in the full space boundary basis
1278bc4ae4beSMatthew G. Knepley 
1279bc4ae4beSMatthew G. Knepley   Not collective
1280bc4ae4beSMatthew G. Knepley 
1281bc4ae4beSMatthew G. Knepley   Input Parameters:
1282bc4ae4beSMatthew G. Knepley + prob - The PetscDS object
1283bc4ae4beSMatthew G. Knepley - f - The field number
1284bc4ae4beSMatthew G. Knepley 
1285bc4ae4beSMatthew G. Knepley   Output Parameter:
1286bc4ae4beSMatthew G. Knepley . off - The boundary offset
1287bc4ae4beSMatthew G. Knepley 
1288bc4ae4beSMatthew G. Knepley   Level: beginner
1289bc4ae4beSMatthew G. Knepley 
1290bc4ae4beSMatthew G. Knepley .seealso: PetscDSGetFieldOffset(), PetscDSGetNumFields(), PetscDSCreate()
1291bc4ae4beSMatthew G. Knepley @*/
12922764a2aaSMatthew G. Knepley PetscErrorCode PetscDSGetBdFieldOffset(PetscDS prob, PetscInt f, PetscInt *off)
12932764a2aaSMatthew G. Knepley {
12942764a2aaSMatthew G. Knepley   PetscInt       g;
12952764a2aaSMatthew G. Knepley   PetscErrorCode ierr;
12962764a2aaSMatthew G. Knepley 
12972764a2aaSMatthew G. Knepley   PetscFunctionBegin;
12982764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
12992764a2aaSMatthew G. Knepley   PetscValidPointer(off, 3);
13002764a2aaSMatthew 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);
13012764a2aaSMatthew G. Knepley   *off = 0;
13022764a2aaSMatthew G. Knepley   for (g = 0; g < f; ++g) {
13032764a2aaSMatthew G. Knepley     PetscFE  fe = (PetscFE) prob->discBd[g];
13042764a2aaSMatthew G. Knepley     PetscInt Nb, Nc;
13052764a2aaSMatthew G. Knepley 
13062764a2aaSMatthew G. Knepley     ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr);
13072764a2aaSMatthew G. Knepley     ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
13082764a2aaSMatthew G. Knepley     *off += Nb*Nc;
13096ce16762SMatthew G. Knepley   }
13106ce16762SMatthew G. Knepley   PetscFunctionReturn(0);
13116ce16762SMatthew G. Knepley }
13126ce16762SMatthew G. Knepley 
13136ce16762SMatthew G. Knepley #undef __FUNCT__
13146ce16762SMatthew G. Knepley #define __FUNCT__ "PetscDSGetComponentOffset"
13156ce16762SMatthew G. Knepley /*@
13166ce16762SMatthew G. Knepley   PetscDSGetComponentOffset - Returns the offset of the given field on an evaluation point
13176ce16762SMatthew G. Knepley 
13186ce16762SMatthew G. Knepley   Not collective
13196ce16762SMatthew G. Knepley 
13206ce16762SMatthew G. Knepley   Input Parameters:
13216ce16762SMatthew G. Knepley + prob - The PetscDS object
13226ce16762SMatthew G. Knepley - f - The field number
13236ce16762SMatthew G. Knepley 
13246ce16762SMatthew G. Knepley   Output Parameter:
13256ce16762SMatthew G. Knepley . off - The offset
13266ce16762SMatthew G. Knepley 
13276ce16762SMatthew G. Knepley   Level: beginner
13286ce16762SMatthew G. Knepley 
13296ce16762SMatthew G. Knepley .seealso: PetscDSGetBdFieldOffset(), PetscDSGetNumFields(), PetscDSCreate()
13306ce16762SMatthew G. Knepley @*/
13316ce16762SMatthew G. Knepley PetscErrorCode PetscDSGetComponentOffset(PetscDS prob, PetscInt f, PetscInt *off)
13326ce16762SMatthew G. Knepley {
13336ce16762SMatthew G. Knepley   PetscInt       g;
13346ce16762SMatthew G. Knepley   PetscErrorCode ierr;
13356ce16762SMatthew G. Knepley 
13366ce16762SMatthew G. Knepley   PetscFunctionBegin;
13376ce16762SMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
13386ce16762SMatthew G. Knepley   PetscValidPointer(off, 3);
13396ce16762SMatthew 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);
13406ce16762SMatthew G. Knepley   *off = 0;
13416ce16762SMatthew G. Knepley   for (g = 0; g < f; ++g) {
13426ce16762SMatthew G. Knepley     PetscFE  fe = (PetscFE) prob->disc[g];
13436ce16762SMatthew G. Knepley     PetscInt Nc;
13446ce16762SMatthew G. Knepley 
13456ce16762SMatthew G. Knepley     ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
13466ce16762SMatthew G. Knepley     *off += Nc;
13472764a2aaSMatthew G. Knepley   }
13482764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
13492764a2aaSMatthew G. Knepley }
13502764a2aaSMatthew G. Knepley 
13512764a2aaSMatthew G. Knepley #undef __FUNCT__
13522764a2aaSMatthew G. Knepley #define __FUNCT__ "PetscDSGetTabulation"
135368c9edb9SMatthew G. Knepley /*@C
135468c9edb9SMatthew G. Knepley   PetscDSGetTabulation - Return the basis tabulation at quadrature points for the volume discretization
135568c9edb9SMatthew G. Knepley 
135668c9edb9SMatthew G. Knepley   Not collective
135768c9edb9SMatthew G. Knepley 
135868c9edb9SMatthew G. Knepley   Input Parameter:
135968c9edb9SMatthew G. Knepley . prob - The PetscDS object
136068c9edb9SMatthew G. Knepley 
136168c9edb9SMatthew G. Knepley   Output Parameters:
136268c9edb9SMatthew G. Knepley + basis - The basis function tabulation at quadrature points
136368c9edb9SMatthew G. Knepley - basisDer - The basis function derivative tabulation at quadrature points
136468c9edb9SMatthew G. Knepley 
136568c9edb9SMatthew G. Knepley   Level: intermediate
136668c9edb9SMatthew G. Knepley 
136768c9edb9SMatthew G. Knepley .seealso: PetscDSGetBdTabulation(), PetscDSCreate()
136868c9edb9SMatthew G. Knepley @*/
13692764a2aaSMatthew G. Knepley PetscErrorCode PetscDSGetTabulation(PetscDS prob, PetscReal ***basis, PetscReal ***basisDer)
13702764a2aaSMatthew G. Knepley {
13712764a2aaSMatthew G. Knepley   PetscErrorCode ierr;
13722764a2aaSMatthew G. Knepley 
13732764a2aaSMatthew G. Knepley   PetscFunctionBegin;
13742764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
13752764a2aaSMatthew G. Knepley   ierr = PetscDSSetUp(prob);CHKERRQ(ierr);
13762764a2aaSMatthew G. Knepley   if (basis)    {PetscValidPointer(basis, 2);    *basis    = prob->basis;}
13772764a2aaSMatthew G. Knepley   if (basisDer) {PetscValidPointer(basisDer, 3); *basisDer = prob->basisDer;}
13782764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
13792764a2aaSMatthew G. Knepley }
13802764a2aaSMatthew G. Knepley 
13812764a2aaSMatthew G. Knepley #undef __FUNCT__
13822764a2aaSMatthew G. Knepley #define __FUNCT__ "PetscDSGetBdTabulation"
138368c9edb9SMatthew G. Knepley /*@C
138468c9edb9SMatthew G. Knepley   PetscDSGetBdTabulation - Return the basis tabulation at quadrature points for the boundary discretization
138568c9edb9SMatthew G. Knepley 
138668c9edb9SMatthew G. Knepley   Not collective
138768c9edb9SMatthew G. Knepley 
138868c9edb9SMatthew G. Knepley   Input Parameter:
138968c9edb9SMatthew G. Knepley . prob - The PetscDS object
139068c9edb9SMatthew G. Knepley 
139168c9edb9SMatthew G. Knepley   Output Parameters:
139268c9edb9SMatthew G. Knepley + basis - The basis function tabulation at quadrature points
139368c9edb9SMatthew G. Knepley - basisDer - The basis function derivative tabulation at quadrature points
139468c9edb9SMatthew G. Knepley 
139568c9edb9SMatthew G. Knepley   Level: intermediate
139668c9edb9SMatthew G. Knepley 
139768c9edb9SMatthew G. Knepley .seealso: PetscDSGetTabulation(), PetscDSCreate()
139868c9edb9SMatthew G. Knepley @*/
13992764a2aaSMatthew G. Knepley PetscErrorCode PetscDSGetBdTabulation(PetscDS prob, PetscReal ***basis, PetscReal ***basisDer)
14002764a2aaSMatthew G. Knepley {
14012764a2aaSMatthew G. Knepley   PetscErrorCode ierr;
14022764a2aaSMatthew G. Knepley 
14032764a2aaSMatthew G. Knepley   PetscFunctionBegin;
14042764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
14052764a2aaSMatthew G. Knepley   ierr = PetscDSSetUp(prob);CHKERRQ(ierr);
14062764a2aaSMatthew G. Knepley   if (basis)    {PetscValidPointer(basis, 2);    *basis    = prob->basisBd;}
14072764a2aaSMatthew G. Knepley   if (basisDer) {PetscValidPointer(basisDer, 3); *basisDer = prob->basisDerBd;}
14082764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
14092764a2aaSMatthew G. Knepley }
14102764a2aaSMatthew G. Knepley 
14112764a2aaSMatthew G. Knepley #undef __FUNCT__
14122764a2aaSMatthew G. Knepley #define __FUNCT__ "PetscDSGetEvaluationArrays"
14132764a2aaSMatthew G. Knepley PetscErrorCode PetscDSGetEvaluationArrays(PetscDS prob, PetscScalar **u, PetscScalar **u_t, PetscScalar **u_x)
14142764a2aaSMatthew G. Knepley {
14152764a2aaSMatthew G. Knepley   PetscErrorCode ierr;
14162764a2aaSMatthew G. Knepley 
14172764a2aaSMatthew G. Knepley   PetscFunctionBegin;
14182764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
14192764a2aaSMatthew G. Knepley   ierr = PetscDSSetUp(prob);CHKERRQ(ierr);
14202764a2aaSMatthew G. Knepley   if (u)   {PetscValidPointer(u, 2);   *u   = prob->u;}
14212764a2aaSMatthew G. Knepley   if (u_t) {PetscValidPointer(u_t, 3); *u_t = prob->u_t;}
14222764a2aaSMatthew G. Knepley   if (u_x) {PetscValidPointer(u_x, 4); *u_x = prob->u_x;}
14232764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
14242764a2aaSMatthew G. Knepley }
14252764a2aaSMatthew G. Knepley 
14262764a2aaSMatthew G. Knepley #undef __FUNCT__
14272764a2aaSMatthew G. Knepley #define __FUNCT__ "PetscDSGetWeakFormArrays"
14282764a2aaSMatthew G. Knepley PetscErrorCode PetscDSGetWeakFormArrays(PetscDS prob, PetscScalar **f0, PetscScalar **f1, PetscScalar **g0, PetscScalar **g1, PetscScalar **g2, PetscScalar **g3)
14292764a2aaSMatthew G. Knepley {
14302764a2aaSMatthew G. Knepley   PetscErrorCode ierr;
14312764a2aaSMatthew G. Knepley 
14322764a2aaSMatthew G. Knepley   PetscFunctionBegin;
14332764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
14342764a2aaSMatthew G. Knepley   ierr = PetscDSSetUp(prob);CHKERRQ(ierr);
14352764a2aaSMatthew G. Knepley   if (f0) {PetscValidPointer(f0, 2); *f0 = prob->f0;}
14362764a2aaSMatthew G. Knepley   if (f1) {PetscValidPointer(f1, 3); *f1 = prob->f1;}
14372764a2aaSMatthew G. Knepley   if (g0) {PetscValidPointer(g0, 4); *g0 = prob->g0;}
14382764a2aaSMatthew G. Knepley   if (g1) {PetscValidPointer(g1, 5); *g1 = prob->g1;}
14392764a2aaSMatthew G. Knepley   if (g2) {PetscValidPointer(g2, 6); *g2 = prob->g2;}
14402764a2aaSMatthew G. Knepley   if (g3) {PetscValidPointer(g3, 7); *g3 = prob->g3;}
14412764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
14422764a2aaSMatthew G. Knepley }
14432764a2aaSMatthew G. Knepley 
14442764a2aaSMatthew G. Knepley #undef __FUNCT__
14452764a2aaSMatthew G. Knepley #define __FUNCT__ "PetscDSGetRefCoordArrays"
14462764a2aaSMatthew G. Knepley PetscErrorCode PetscDSGetRefCoordArrays(PetscDS prob, PetscReal **x, PetscScalar **refSpaceDer)
14472764a2aaSMatthew G. Knepley {
14482764a2aaSMatthew G. Knepley   PetscErrorCode ierr;
14492764a2aaSMatthew G. Knepley 
14502764a2aaSMatthew G. Knepley   PetscFunctionBegin;
14512764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
14522764a2aaSMatthew G. Knepley   ierr = PetscDSSetUp(prob);CHKERRQ(ierr);
14532764a2aaSMatthew G. Knepley   if (x)           {PetscValidPointer(x, 2);           *x           = prob->x;}
14542764a2aaSMatthew G. Knepley   if (refSpaceDer) {PetscValidPointer(refSpaceDer, 3); *refSpaceDer = prob->refSpaceDer;}
14552764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
14562764a2aaSMatthew G. Knepley }
14572764a2aaSMatthew G. Knepley 
14582764a2aaSMatthew G. Knepley #undef __FUNCT__
14592764a2aaSMatthew G. Knepley #define __FUNCT__ "PetscDSDestroy_Basic"
1460bc4ae4beSMatthew G. Knepley static PetscErrorCode PetscDSDestroy_Basic(PetscDS prob)
14612764a2aaSMatthew G. Knepley {
14622764a2aaSMatthew G. Knepley   PetscFunctionBegin;
14632764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
14642764a2aaSMatthew G. Knepley }
14652764a2aaSMatthew G. Knepley 
14662764a2aaSMatthew G. Knepley #undef __FUNCT__
14672764a2aaSMatthew G. Knepley #define __FUNCT__ "PetscDSInitialize_Basic"
1468bc4ae4beSMatthew G. Knepley static PetscErrorCode PetscDSInitialize_Basic(PetscDS prob)
14692764a2aaSMatthew G. Knepley {
14702764a2aaSMatthew G. Knepley   PetscFunctionBegin;
14712764a2aaSMatthew G. Knepley   prob->ops->setfromoptions = NULL;
14722764a2aaSMatthew G. Knepley   prob->ops->setup          = NULL;
14732764a2aaSMatthew G. Knepley   prob->ops->view           = NULL;
14742764a2aaSMatthew G. Knepley   prob->ops->destroy        = PetscDSDestroy_Basic;
14752764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
14762764a2aaSMatthew G. Knepley }
14772764a2aaSMatthew G. Knepley 
14782764a2aaSMatthew G. Knepley /*MC
14792764a2aaSMatthew G. Knepley   PETSCDSBASIC = "basic" - A discrete system with pointwise residual and boundary residual functions
14802764a2aaSMatthew G. Knepley 
14812764a2aaSMatthew G. Knepley   Level: intermediate
14822764a2aaSMatthew G. Knepley 
14832764a2aaSMatthew G. Knepley .seealso: PetscDSType, PetscDSCreate(), PetscDSSetType()
14842764a2aaSMatthew G. Knepley M*/
14852764a2aaSMatthew G. Knepley 
14862764a2aaSMatthew G. Knepley #undef __FUNCT__
14872764a2aaSMatthew G. Knepley #define __FUNCT__ "PetscDSCreate_Basic"
14882764a2aaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDSCreate_Basic(PetscDS prob)
14892764a2aaSMatthew G. Knepley {
14902764a2aaSMatthew G. Knepley   PetscDS_Basic *b;
14912764a2aaSMatthew G. Knepley   PetscErrorCode      ierr;
14922764a2aaSMatthew G. Knepley 
14932764a2aaSMatthew G. Knepley   PetscFunctionBegin;
14942764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCSPACE_CLASSID, 1);
14952764a2aaSMatthew G. Knepley   ierr       = PetscNewLog(prob, &b);CHKERRQ(ierr);
14962764a2aaSMatthew G. Knepley   prob->data = b;
14972764a2aaSMatthew G. Knepley 
14982764a2aaSMatthew G. Knepley   ierr = PetscDSInitialize_Basic(prob);CHKERRQ(ierr);
14992764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
15002764a2aaSMatthew G. Knepley }
1501