xref: /petsc/src/dm/dt/interface/dtds.c (revision c371a6d18df5f2de8bc93ba090321da4edb93f57)
1af0996ceSBarry Smith #include <petsc/private/petscdsimpl.h> /*I "petscds.h" I*/
22764a2aaSMatthew G. Knepley 
32764a2aaSMatthew G. Knepley PetscClassId PETSCDS_CLASSID = 0;
42764a2aaSMatthew G. Knepley 
52764a2aaSMatthew G. Knepley PetscFunctionList PetscDSList              = NULL;
62764a2aaSMatthew G. Knepley PetscBool         PetscDSRegisterAllCalled = PETSC_FALSE;
72764a2aaSMatthew G. Knepley 
82764a2aaSMatthew G. Knepley /*@C
92764a2aaSMatthew G. Knepley   PetscDSRegister - Adds a new PetscDS implementation
102764a2aaSMatthew G. Knepley 
112764a2aaSMatthew G. Knepley   Not Collective
122764a2aaSMatthew G. Knepley 
132764a2aaSMatthew G. Knepley   Input Parameters:
142764a2aaSMatthew G. Knepley + name        - The name of a new user-defined creation routine
152764a2aaSMatthew G. Knepley - create_func - The creation routine itself
162764a2aaSMatthew G. Knepley 
172764a2aaSMatthew G. Knepley   Notes:
182764a2aaSMatthew G. Knepley   PetscDSRegister() may be called multiple times to add several user-defined PetscDSs
192764a2aaSMatthew G. Knepley 
202764a2aaSMatthew G. Knepley   Sample usage:
212764a2aaSMatthew G. Knepley .vb
222764a2aaSMatthew G. Knepley     PetscDSRegister("my_ds", MyPetscDSCreate);
232764a2aaSMatthew G. Knepley .ve
242764a2aaSMatthew G. Knepley 
252764a2aaSMatthew G. Knepley   Then, your PetscDS type can be chosen with the procedural interface via
262764a2aaSMatthew G. Knepley .vb
272764a2aaSMatthew G. Knepley     PetscDSCreate(MPI_Comm, PetscDS *);
282764a2aaSMatthew G. Knepley     PetscDSSetType(PetscDS, "my_ds");
292764a2aaSMatthew G. Knepley .ve
302764a2aaSMatthew G. Knepley    or at runtime via the option
312764a2aaSMatthew G. Knepley .vb
322764a2aaSMatthew G. Knepley     -petscds_type my_ds
332764a2aaSMatthew G. Knepley .ve
342764a2aaSMatthew G. Knepley 
352764a2aaSMatthew G. Knepley   Level: advanced
362764a2aaSMatthew G. Knepley 
372764a2aaSMatthew G. Knepley .keywords: PetscDS, register
382764a2aaSMatthew G. Knepley .seealso: PetscDSRegisterAll(), PetscDSRegisterDestroy()
392764a2aaSMatthew G. Knepley 
402764a2aaSMatthew G. Knepley @*/
412764a2aaSMatthew G. Knepley PetscErrorCode PetscDSRegister(const char sname[], PetscErrorCode (*function)(PetscDS))
422764a2aaSMatthew G. Knepley {
432764a2aaSMatthew G. Knepley   PetscErrorCode ierr;
442764a2aaSMatthew G. Knepley 
452764a2aaSMatthew G. Knepley   PetscFunctionBegin;
462764a2aaSMatthew G. Knepley   ierr = PetscFunctionListAdd(&PetscDSList, sname, function);CHKERRQ(ierr);
472764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
482764a2aaSMatthew G. Knepley }
492764a2aaSMatthew G. Knepley 
502764a2aaSMatthew G. Knepley /*@C
512764a2aaSMatthew G. Knepley   PetscDSSetType - Builds a particular PetscDS
522764a2aaSMatthew G. Knepley 
532764a2aaSMatthew G. Knepley   Collective on PetscDS
542764a2aaSMatthew G. Knepley 
552764a2aaSMatthew G. Knepley   Input Parameters:
562764a2aaSMatthew G. Knepley + prob - The PetscDS object
572764a2aaSMatthew G. Knepley - name - The kind of system
582764a2aaSMatthew G. Knepley 
592764a2aaSMatthew G. Knepley   Options Database Key:
602764a2aaSMatthew G. Knepley . -petscds_type <type> - Sets the PetscDS type; use -help for a list of available types
612764a2aaSMatthew G. Knepley 
622764a2aaSMatthew G. Knepley   Level: intermediate
632764a2aaSMatthew G. Knepley 
642764a2aaSMatthew G. Knepley .keywords: PetscDS, set, type
652764a2aaSMatthew G. Knepley .seealso: PetscDSGetType(), PetscDSCreate()
662764a2aaSMatthew G. Knepley @*/
672764a2aaSMatthew G. Knepley PetscErrorCode PetscDSSetType(PetscDS prob, PetscDSType name)
682764a2aaSMatthew G. Knepley {
692764a2aaSMatthew G. Knepley   PetscErrorCode (*r)(PetscDS);
702764a2aaSMatthew G. Knepley   PetscBool      match;
712764a2aaSMatthew G. Knepley   PetscErrorCode ierr;
722764a2aaSMatthew G. Knepley 
732764a2aaSMatthew G. Knepley   PetscFunctionBegin;
742764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
752764a2aaSMatthew G. Knepley   ierr = PetscObjectTypeCompare((PetscObject) prob, name, &match);CHKERRQ(ierr);
762764a2aaSMatthew G. Knepley   if (match) PetscFunctionReturn(0);
772764a2aaSMatthew G. Knepley 
780f51fdf8SToby Isaac   ierr = PetscDSRegisterAll();CHKERRQ(ierr);
792764a2aaSMatthew G. Knepley   ierr = PetscFunctionListFind(PetscDSList, name, &r);CHKERRQ(ierr);
802764a2aaSMatthew G. Knepley   if (!r) SETERRQ1(PetscObjectComm((PetscObject) prob), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown PetscDS type: %s", name);
812764a2aaSMatthew G. Knepley 
822764a2aaSMatthew G. Knepley   if (prob->ops->destroy) {
832764a2aaSMatthew G. Knepley     ierr             = (*prob->ops->destroy)(prob);CHKERRQ(ierr);
842764a2aaSMatthew G. Knepley     prob->ops->destroy = NULL;
852764a2aaSMatthew G. Knepley   }
862764a2aaSMatthew G. Knepley   ierr = (*r)(prob);CHKERRQ(ierr);
872764a2aaSMatthew G. Knepley   ierr = PetscObjectChangeTypeName((PetscObject) prob, name);CHKERRQ(ierr);
882764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
892764a2aaSMatthew G. Knepley }
902764a2aaSMatthew G. Knepley 
912764a2aaSMatthew G. Knepley /*@C
922764a2aaSMatthew G. Knepley   PetscDSGetType - Gets the PetscDS type name (as a string) from the object.
932764a2aaSMatthew G. Knepley 
942764a2aaSMatthew G. Knepley   Not Collective
952764a2aaSMatthew G. Knepley 
962764a2aaSMatthew G. Knepley   Input Parameter:
972764a2aaSMatthew G. Knepley . prob  - The PetscDS
982764a2aaSMatthew G. Knepley 
992764a2aaSMatthew G. Knepley   Output Parameter:
1002764a2aaSMatthew G. Knepley . name - The PetscDS type name
1012764a2aaSMatthew G. Knepley 
1022764a2aaSMatthew G. Knepley   Level: intermediate
1032764a2aaSMatthew G. Knepley 
1042764a2aaSMatthew G. Knepley .keywords: PetscDS, get, type, name
1052764a2aaSMatthew G. Knepley .seealso: PetscDSSetType(), PetscDSCreate()
1062764a2aaSMatthew G. Knepley @*/
1072764a2aaSMatthew G. Knepley PetscErrorCode PetscDSGetType(PetscDS prob, PetscDSType *name)
1082764a2aaSMatthew G. Knepley {
1092764a2aaSMatthew G. Knepley   PetscErrorCode ierr;
1102764a2aaSMatthew G. Knepley 
1112764a2aaSMatthew G. Knepley   PetscFunctionBegin;
1122764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
113c959eef4SJed Brown   PetscValidPointer(name, 2);
1140f51fdf8SToby Isaac   ierr = PetscDSRegisterAll();CHKERRQ(ierr);
1152764a2aaSMatthew G. Knepley   *name = ((PetscObject) prob)->type_name;
1162764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
1172764a2aaSMatthew G. Knepley }
1182764a2aaSMatthew G. Knepley 
1197d8a60eaSMatthew G. Knepley static PetscErrorCode PetscDSView_Ascii(PetscDS prob, PetscViewer viewer)
1207d8a60eaSMatthew G. Knepley {
1217d8a60eaSMatthew G. Knepley   PetscViewerFormat format;
1227d8a60eaSMatthew G. Knepley   PetscInt          f;
1237d8a60eaSMatthew G. Knepley   PetscErrorCode    ierr;
1247d8a60eaSMatthew G. Knepley 
1257d8a60eaSMatthew G. Knepley   PetscFunctionBegin;
1267d8a60eaSMatthew G. Knepley   ierr = PetscViewerGetFormat(viewer, &format);CHKERRQ(ierr);
1277d8a60eaSMatthew G. Knepley   ierr = PetscViewerASCIIPrintf(viewer, "Discrete System with %d fields\n", prob->Nf);CHKERRQ(ierr);
1287d8a60eaSMatthew G. Knepley   ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
1297d8a60eaSMatthew G. Knepley   for (f = 0; f < prob->Nf; ++f) {
1307d8a60eaSMatthew G. Knepley     PetscObject  obj;
1317d8a60eaSMatthew G. Knepley     PetscClassId id;
1327d8a60eaSMatthew G. Knepley     const char  *name;
1337d8a60eaSMatthew G. Knepley     PetscInt     Nc;
1347d8a60eaSMatthew G. Knepley 
1357d8a60eaSMatthew G. Knepley     ierr = PetscDSGetDiscretization(prob, f, &obj);CHKERRQ(ierr);
1367d8a60eaSMatthew G. Knepley     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
1377d8a60eaSMatthew G. Knepley     ierr = PetscObjectGetName(obj, &name);CHKERRQ(ierr);
1387d8a60eaSMatthew G. Knepley     ierr = PetscViewerASCIIPrintf(viewer, "Field %s", name ? name : "<unknown>");CHKERRQ(ierr);
1397d8a60eaSMatthew G. Knepley     if (id == PETSCFE_CLASSID)      {
1407d8a60eaSMatthew G. Knepley       ierr = PetscFEGetNumComponents((PetscFE) obj, &Nc);CHKERRQ(ierr);
1417d8a60eaSMatthew G. Knepley       ierr = PetscViewerASCIIPrintf(viewer, " FEM");CHKERRQ(ierr);
1427d8a60eaSMatthew G. Knepley     } else if (id == PETSCFV_CLASSID) {
1437d8a60eaSMatthew G. Knepley       ierr = PetscFVGetNumComponents((PetscFV) obj, &Nc);CHKERRQ(ierr);
1447d8a60eaSMatthew G. Knepley       ierr = PetscViewerASCIIPrintf(viewer, " FVM");CHKERRQ(ierr);
1457d8a60eaSMatthew G. Knepley     }
1467d8a60eaSMatthew G. Knepley     else SETERRQ1(PetscObjectComm((PetscObject) prob), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", f);
1477d8a60eaSMatthew G. Knepley     if (Nc > 1) {ierr = PetscViewerASCIIPrintf(viewer, "%d components", Nc);CHKERRQ(ierr);}
1487d8a60eaSMatthew G. Knepley     else        {ierr = PetscViewerASCIIPrintf(viewer, "%d component ", Nc);CHKERRQ(ierr);}
149249df284SMatthew G. Knepley     if (prob->implicit[f]) {ierr = PetscViewerASCIIPrintf(viewer, " (implicit)");CHKERRQ(ierr);}
150249df284SMatthew G. Knepley     else                   {ierr = PetscViewerASCIIPrintf(viewer, " (explicit)");CHKERRQ(ierr);}
151a6cbbb48SMatthew G. Knepley     if (prob->adjacency[f*2+0]) {
152a6cbbb48SMatthew G. Knepley       if (prob->adjacency[f*2+1]) {ierr = PetscViewerASCIIPrintf(viewer, " (adj FVM++)");CHKERRQ(ierr);}
153a6cbbb48SMatthew G. Knepley       else                        {ierr = PetscViewerASCIIPrintf(viewer, " (adj FVM)");CHKERRQ(ierr);}
154a6cbbb48SMatthew G. Knepley     } else {
155a6cbbb48SMatthew G. Knepley       if (prob->adjacency[f*2+1]) {ierr = PetscViewerASCIIPrintf(viewer, " (adj FEM)");CHKERRQ(ierr);}
156a6cbbb48SMatthew G. Knepley       else                        {ierr = PetscViewerASCIIPrintf(viewer, " (adj FUNKY)");CHKERRQ(ierr);}
157a6cbbb48SMatthew G. Knepley     }
1587d8a60eaSMatthew G. Knepley     ierr = PetscViewerASCIIPrintf(viewer, "\n");CHKERRQ(ierr);
1597d8a60eaSMatthew G. Knepley     if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
1607d8a60eaSMatthew G. Knepley       if (id == PETSCFE_CLASSID)      {ierr = PetscFEView((PetscFE) obj, viewer);CHKERRQ(ierr);}
1617d8a60eaSMatthew G. Knepley       else if (id == PETSCFV_CLASSID) {ierr = PetscFVView((PetscFV) obj, viewer);CHKERRQ(ierr);}
1627d8a60eaSMatthew G. Knepley     }
1637d8a60eaSMatthew G. Knepley   }
1647d8a60eaSMatthew G. Knepley   ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1657d8a60eaSMatthew G. Knepley   PetscFunctionReturn(0);
1667d8a60eaSMatthew G. Knepley }
1677d8a60eaSMatthew G. Knepley 
1682764a2aaSMatthew G. Knepley /*@C
1692764a2aaSMatthew G. Knepley   PetscDSView - Views a PetscDS
1702764a2aaSMatthew G. Knepley 
1712764a2aaSMatthew G. Knepley   Collective on PetscDS
1722764a2aaSMatthew G. Knepley 
1732764a2aaSMatthew G. Knepley   Input Parameter:
1742764a2aaSMatthew G. Knepley + prob - the PetscDS object to view
1752764a2aaSMatthew G. Knepley - v  - the viewer
1762764a2aaSMatthew G. Knepley 
1772764a2aaSMatthew G. Knepley   Level: developer
1782764a2aaSMatthew G. Knepley 
1792764a2aaSMatthew G. Knepley .seealso PetscDSDestroy()
1802764a2aaSMatthew G. Knepley @*/
1812764a2aaSMatthew G. Knepley PetscErrorCode PetscDSView(PetscDS prob, PetscViewer v)
1822764a2aaSMatthew G. Knepley {
1837d8a60eaSMatthew G. Knepley   PetscBool      iascii;
1842764a2aaSMatthew G. Knepley   PetscErrorCode ierr;
1852764a2aaSMatthew G. Knepley 
1862764a2aaSMatthew G. Knepley   PetscFunctionBegin;
1872764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1882764a2aaSMatthew G. Knepley   if (!v) {ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject) prob), &v);CHKERRQ(ierr);}
1897d8a60eaSMatthew G. Knepley   else    {PetscValidHeaderSpecific(v, PETSC_VIEWER_CLASSID, 2);}
1907d8a60eaSMatthew G. Knepley   ierr = PetscObjectTypeCompare((PetscObject) v, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
1917d8a60eaSMatthew G. Knepley   if (iascii) {ierr = PetscDSView_Ascii(prob, v);CHKERRQ(ierr);}
1922764a2aaSMatthew G. Knepley   if (prob->ops->view) {ierr = (*prob->ops->view)(prob, v);CHKERRQ(ierr);}
1932764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
1942764a2aaSMatthew G. Knepley }
1952764a2aaSMatthew G. Knepley 
1962764a2aaSMatthew G. Knepley /*@
1972764a2aaSMatthew G. Knepley   PetscDSSetFromOptions - sets parameters in a PetscDS from the options database
1982764a2aaSMatthew G. Knepley 
1992764a2aaSMatthew G. Knepley   Collective on PetscDS
2002764a2aaSMatthew G. Knepley 
2012764a2aaSMatthew G. Knepley   Input Parameter:
2022764a2aaSMatthew G. Knepley . prob - the PetscDS object to set options for
2032764a2aaSMatthew G. Knepley 
2042764a2aaSMatthew G. Knepley   Options Database:
2052764a2aaSMatthew G. Knepley 
2062764a2aaSMatthew G. Knepley   Level: developer
2072764a2aaSMatthew G. Knepley 
2082764a2aaSMatthew G. Knepley .seealso PetscDSView()
2092764a2aaSMatthew G. Knepley @*/
2102764a2aaSMatthew G. Knepley PetscErrorCode PetscDSSetFromOptions(PetscDS prob)
2112764a2aaSMatthew G. Knepley {
212f1fd5e65SToby Isaac   DSBoundary     b;
2132764a2aaSMatthew G. Knepley   const char    *defaultType;
2142764a2aaSMatthew G. Knepley   char           name[256];
2152764a2aaSMatthew G. Knepley   PetscBool      flg;
2162764a2aaSMatthew G. Knepley   PetscErrorCode ierr;
2172764a2aaSMatthew G. Knepley 
2182764a2aaSMatthew G. Knepley   PetscFunctionBegin;
2192764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
2202764a2aaSMatthew G. Knepley   if (!((PetscObject) prob)->type_name) {
2212764a2aaSMatthew G. Knepley     defaultType = PETSCDSBASIC;
2222764a2aaSMatthew G. Knepley   } else {
2232764a2aaSMatthew G. Knepley     defaultType = ((PetscObject) prob)->type_name;
2242764a2aaSMatthew G. Knepley   }
2250f51fdf8SToby Isaac   ierr = PetscDSRegisterAll();CHKERRQ(ierr);
2262764a2aaSMatthew G. Knepley 
2272764a2aaSMatthew G. Knepley   ierr = PetscObjectOptionsBegin((PetscObject) prob);CHKERRQ(ierr);
228f1fd5e65SToby Isaac   for (b = prob->boundary; b; b = b->next) {
229f1fd5e65SToby Isaac     char       optname[1024];
230f1fd5e65SToby Isaac     PetscInt   ids[1024], len = 1024;
231f1fd5e65SToby Isaac     PetscBool  flg;
232f1fd5e65SToby Isaac 
233f1fd5e65SToby Isaac     ierr = PetscSNPrintf(optname, sizeof(optname), "-bc_%s", b->name);CHKERRQ(ierr);
234f1fd5e65SToby Isaac     ierr = PetscMemzero(ids, sizeof(ids));CHKERRQ(ierr);
235f1fd5e65SToby Isaac     ierr = PetscOptionsIntArray(optname, "List of boundary IDs", "", ids, &len, &flg);CHKERRQ(ierr);
236f1fd5e65SToby Isaac     if (flg) {
237f1fd5e65SToby Isaac       b->numids = len;
238f1fd5e65SToby Isaac       ierr = PetscFree(b->ids);CHKERRQ(ierr);
239f1fd5e65SToby Isaac       ierr = PetscMalloc1(len, &b->ids);CHKERRQ(ierr);
240f1fd5e65SToby Isaac       ierr = PetscMemcpy(b->ids, ids, len*sizeof(PetscInt));CHKERRQ(ierr);
241f1fd5e65SToby Isaac     }
242e7b0402cSSander Arens     len = 1024;
243f1fd5e65SToby Isaac     ierr = PetscSNPrintf(optname, sizeof(optname), "-bc_%s_comp", b->name);CHKERRQ(ierr);
244f1fd5e65SToby Isaac     ierr = PetscMemzero(ids, sizeof(ids));CHKERRQ(ierr);
245f1fd5e65SToby Isaac     ierr = PetscOptionsIntArray(optname, "List of boundary field components", "", ids, &len, &flg);CHKERRQ(ierr);
246f1fd5e65SToby Isaac     if (flg) {
247f1fd5e65SToby Isaac       b->numcomps = len;
248f1fd5e65SToby Isaac       ierr = PetscFree(b->comps);CHKERRQ(ierr);
249f1fd5e65SToby Isaac       ierr = PetscMalloc1(len, &b->comps);CHKERRQ(ierr);
250f1fd5e65SToby Isaac       ierr = PetscMemcpy(b->comps, ids, len*sizeof(PetscInt));CHKERRQ(ierr);
251f1fd5e65SToby Isaac     }
252f1fd5e65SToby Isaac   }
2532764a2aaSMatthew G. Knepley   ierr = PetscOptionsFList("-petscds_type", "Discrete System", "PetscDSSetType", PetscDSList, defaultType, name, 256, &flg);CHKERRQ(ierr);
2542764a2aaSMatthew G. Knepley   if (flg) {
2552764a2aaSMatthew G. Knepley     ierr = PetscDSSetType(prob, name);CHKERRQ(ierr);
2562764a2aaSMatthew G. Knepley   } else if (!((PetscObject) prob)->type_name) {
2572764a2aaSMatthew G. Knepley     ierr = PetscDSSetType(prob, defaultType);CHKERRQ(ierr);
2582764a2aaSMatthew G. Knepley   }
2592764a2aaSMatthew G. Knepley   if (prob->ops->setfromoptions) {ierr = (*prob->ops->setfromoptions)(prob);CHKERRQ(ierr);}
2602764a2aaSMatthew G. Knepley   /* process any options handlers added with PetscObjectAddOptionsHandler() */
2610633abcbSJed Brown   ierr = PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject) prob);CHKERRQ(ierr);
2622764a2aaSMatthew G. Knepley   ierr = PetscOptionsEnd();CHKERRQ(ierr);
2632764a2aaSMatthew G. Knepley   ierr = PetscDSViewFromOptions(prob, NULL, "-petscds_view");CHKERRQ(ierr);
2642764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
2652764a2aaSMatthew G. Knepley }
2662764a2aaSMatthew G. Knepley 
2672764a2aaSMatthew G. Knepley /*@C
2682764a2aaSMatthew G. Knepley   PetscDSSetUp - Construct data structures for the PetscDS
2692764a2aaSMatthew G. Knepley 
2702764a2aaSMatthew G. Knepley   Collective on PetscDS
2712764a2aaSMatthew G. Knepley 
2722764a2aaSMatthew G. Knepley   Input Parameter:
2732764a2aaSMatthew G. Knepley . prob - the PetscDS object to setup
2742764a2aaSMatthew G. Knepley 
2752764a2aaSMatthew G. Knepley   Level: developer
2762764a2aaSMatthew G. Knepley 
2772764a2aaSMatthew G. Knepley .seealso PetscDSView(), PetscDSDestroy()
2782764a2aaSMatthew G. Knepley @*/
2792764a2aaSMatthew G. Knepley PetscErrorCode PetscDSSetUp(PetscDS prob)
2802764a2aaSMatthew G. Knepley {
2812764a2aaSMatthew G. Knepley   const PetscInt Nf = prob->Nf;
2822764a2aaSMatthew G. Knepley   PetscInt       dim, work, NcMax = 0, NqMax = 0, f;
2832764a2aaSMatthew G. Knepley   PetscErrorCode ierr;
2842764a2aaSMatthew G. Knepley 
2852764a2aaSMatthew G. Knepley   PetscFunctionBegin;
2862764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
2872764a2aaSMatthew G. Knepley   if (prob->setup) PetscFunctionReturn(0);
2882764a2aaSMatthew G. Knepley   /* Calculate sizes */
2892764a2aaSMatthew G. Knepley   ierr = PetscDSGetSpatialDimension(prob, &dim);CHKERRQ(ierr);
290f744cafaSSander Arens   prob->totDim = prob->totComp = 0;
29147e57110SSander Arens   ierr = PetscMalloc2(Nf,&prob->Nc,Nf,&prob->Nb);CHKERRQ(ierr);
292f744cafaSSander Arens   ierr = PetscCalloc2(Nf+1,&prob->off,Nf+1,&prob->offDer);CHKERRQ(ierr);
293f744cafaSSander Arens   ierr = PetscMalloc4(Nf,&prob->basis,Nf,&prob->basisDer,Nf,&prob->basisFace,Nf,&prob->basisDerFace);CHKERRQ(ierr);
2942764a2aaSMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
2959de99aefSMatthew G. Knepley     PetscObject     obj;
2969de99aefSMatthew G. Knepley     PetscClassId    id;
2972764a2aaSMatthew G. Knepley     PetscQuadrature q;
2989de99aefSMatthew G. Knepley     PetscInt        Nq = 0, Nb, Nc;
2992764a2aaSMatthew G. Knepley 
3009de99aefSMatthew G. Knepley     ierr = PetscDSGetDiscretization(prob, f, &obj);CHKERRQ(ierr);
3019de99aefSMatthew G. Knepley     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
3029de99aefSMatthew G. Knepley     if (id == PETSCFE_CLASSID)      {
3039de99aefSMatthew G. Knepley       PetscFE fe = (PetscFE) obj;
3049de99aefSMatthew G. Knepley 
3052764a2aaSMatthew G. Knepley       ierr = PetscFEGetQuadrature(fe, &q);CHKERRQ(ierr);
3062764a2aaSMatthew G. Knepley       ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr);
3072764a2aaSMatthew G. Knepley       ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
3082764a2aaSMatthew G. Knepley       ierr = PetscFEGetDefaultTabulation(fe, &prob->basis[f], &prob->basisDer[f], NULL);CHKERRQ(ierr);
3094d0b9603SSander Arens       ierr = PetscFEGetFaceTabulation(fe, &prob->basisFace[f], &prob->basisDerFace[f], NULL);CHKERRQ(ierr);
3109de99aefSMatthew G. Knepley     } else if (id == PETSCFV_CLASSID) {
3119de99aefSMatthew G. Knepley       PetscFV fv = (PetscFV) obj;
3129de99aefSMatthew G. Knepley 
3139de99aefSMatthew G. Knepley       ierr = PetscFVGetQuadrature(fv, &q);CHKERRQ(ierr);
3149de99aefSMatthew G. Knepley       ierr = PetscFVGetNumComponents(fv, &Nc);CHKERRQ(ierr);
3159c3cf19fSMatthew G. Knepley       Nb   = Nc;
3166c1a3d01SMatthew G. Knepley       ierr = PetscFVGetDefaultTabulation(fv, &prob->basis[f], &prob->basisDer[f], NULL);CHKERRQ(ierr);
3174d0b9603SSander Arens       /* TODO: should PetscFV also have face tabulation? Otherwise there will be a null pointer in prob->basisFace */
318abac5ca0SMatthew G. Knepley     } else SETERRQ1(PetscObjectComm((PetscObject) prob), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", f);
31947e57110SSander Arens     prob->Nc[f]       = Nc;
32047e57110SSander Arens     prob->Nb[f]       = Nb;
321194d53e6SMatthew G. Knepley     prob->off[f+1]    = Nc     + prob->off[f];
322194d53e6SMatthew G. Knepley     prob->offDer[f+1] = Nc*dim + prob->offDer[f];
323a6b92713SMatthew G. Knepley     if (q) {ierr = PetscQuadratureGetData(q, NULL, NULL, &Nq, NULL, NULL);CHKERRQ(ierr);}
3242764a2aaSMatthew G. Knepley     NqMax          = PetscMax(NqMax, Nq);
3252764a2aaSMatthew G. Knepley     NcMax          = PetscMax(NcMax, Nc);
3269c3cf19fSMatthew G. Knepley     prob->totDim  += Nb;
3272764a2aaSMatthew G. Knepley     prob->totComp += Nc;
3282764a2aaSMatthew G. Knepley   }
3292764a2aaSMatthew G. Knepley   work = PetscMax(prob->totComp*dim, PetscSqr(NcMax*dim));
3302764a2aaSMatthew G. Knepley   /* Allocate works space */
3312764a2aaSMatthew 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);
3322764a2aaSMatthew 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);
3332764a2aaSMatthew G. Knepley   if (prob->ops->setup) {ierr = (*prob->ops->setup)(prob);CHKERRQ(ierr);}
3342764a2aaSMatthew G. Knepley   prob->setup = PETSC_TRUE;
3352764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
3362764a2aaSMatthew G. Knepley }
3372764a2aaSMatthew G. Knepley 
3382764a2aaSMatthew G. Knepley static PetscErrorCode PetscDSDestroyStructs_Static(PetscDS prob)
3392764a2aaSMatthew G. Knepley {
3402764a2aaSMatthew G. Knepley   PetscErrorCode ierr;
3412764a2aaSMatthew G. Knepley 
3422764a2aaSMatthew G. Knepley   PetscFunctionBegin;
34347e57110SSander Arens   ierr = PetscFree2(prob->Nc,prob->Nb);CHKERRQ(ierr);
344f744cafaSSander Arens   ierr = PetscFree2(prob->off,prob->offDer);CHKERRQ(ierr);
345f744cafaSSander Arens   ierr = PetscFree4(prob->basis,prob->basisDer,prob->basisFace,prob->basisDerFace);CHKERRQ(ierr);
3462764a2aaSMatthew G. Knepley   ierr = PetscFree5(prob->u,prob->u_t,prob->u_x,prob->x,prob->refSpaceDer);CHKERRQ(ierr);
3472764a2aaSMatthew G. Knepley   ierr = PetscFree6(prob->f0,prob->f1,prob->g0,prob->g1,prob->g2,prob->g3);CHKERRQ(ierr);
3482764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
3492764a2aaSMatthew G. Knepley }
3502764a2aaSMatthew G. Knepley 
3512764a2aaSMatthew G. Knepley static PetscErrorCode PetscDSEnlarge_Static(PetscDS prob, PetscInt NfNew)
3522764a2aaSMatthew G. Knepley {
353f744cafaSSander Arens   PetscObject      *tmpd;
354a6cbbb48SMatthew G. Knepley   PetscBool        *tmpi, *tmpa;
3552aa1fc23SMatthew G. Knepley   PetscPointFunc   *tmpobj, *tmpf;
356b7e05686SMatthew G. Knepley   PetscPointJac    *tmpg, *tmpgp, *tmpgt;
3572aa1fc23SMatthew G. Knepley   PetscBdPointFunc *tmpfbd;
3582aa1fc23SMatthew G. Knepley   PetscBdPointJac  *tmpgbd;
359194d53e6SMatthew G. Knepley   PetscRiemannFunc *tmpr;
360*c371a6d1SMatthew G. Knepley   PetscSimplePointFunc *tmpexactSol;
3610c2f2876SMatthew G. Knepley   void            **tmpctx;
362a6cbbb48SMatthew G. Knepley   PetscInt          Nf = prob->Nf, f, i;
3632764a2aaSMatthew G. Knepley   PetscErrorCode    ierr;
3642764a2aaSMatthew G. Knepley 
3652764a2aaSMatthew G. Knepley   PetscFunctionBegin;
3662764a2aaSMatthew G. Knepley   if (Nf >= NfNew) PetscFunctionReturn(0);
3672764a2aaSMatthew G. Knepley   prob->setup = PETSC_FALSE;
3682764a2aaSMatthew G. Knepley   ierr = PetscDSDestroyStructs_Static(prob);CHKERRQ(ierr);
369f744cafaSSander Arens   ierr = PetscMalloc3(NfNew, &tmpd, NfNew, &tmpi, NfNew*2, &tmpa);CHKERRQ(ierr);
370f744cafaSSander Arens   for (f = 0; f < Nf; ++f) {tmpd[f] = prob->disc[f]; tmpi[f] = prob->implicit[f]; for (i = 0; i < 2; ++i) tmpa[f*2+i] = prob->adjacency[f*2+i];}
371f744cafaSSander Arens   for (f = Nf; f < NfNew; ++f) {tmpd[f] = NULL; tmpi[f] = PETSC_TRUE; tmpa[f*2+0] = PETSC_FALSE; tmpa[f*2+1] = PETSC_TRUE;}
372f744cafaSSander Arens   ierr = PetscFree3(prob->disc, prob->implicit, prob->adjacency);CHKERRQ(ierr);
3732764a2aaSMatthew G. Knepley   prob->Nf        = NfNew;
3742764a2aaSMatthew G. Knepley   prob->disc      = tmpd;
375249df284SMatthew G. Knepley   prob->implicit  = tmpi;
376a6cbbb48SMatthew G. Knepley   prob->adjacency = tmpa;
377b7e05686SMatthew G. Knepley   ierr = PetscCalloc7(NfNew, &tmpobj, NfNew*2, &tmpf, NfNew*NfNew*4, &tmpg, NfNew*NfNew*4, &tmpgp, NfNew*NfNew*4, &tmpgt, NfNew, &tmpr, NfNew, &tmpctx);CHKERRQ(ierr);
3782764a2aaSMatthew G. Knepley   for (f = 0; f < Nf; ++f) tmpobj[f] = prob->obj[f];
3792764a2aaSMatthew G. Knepley   for (f = 0; f < Nf*2; ++f) tmpf[f] = prob->f[f];
3802764a2aaSMatthew G. Knepley   for (f = 0; f < Nf*Nf*4; ++f) tmpg[f] = prob->g[f];
381475e0ac9SMatthew G. Knepley   for (f = 0; f < Nf*Nf*4; ++f) tmpgp[f] = prob->gp[f];
3820c2f2876SMatthew G. Knepley   for (f = 0; f < Nf; ++f) tmpr[f] = prob->r[f];
3830c2f2876SMatthew G. Knepley   for (f = 0; f < Nf; ++f) tmpctx[f] = prob->ctx[f];
3842764a2aaSMatthew G. Knepley   for (f = Nf; f < NfNew; ++f) tmpobj[f] = NULL;
3852764a2aaSMatthew G. Knepley   for (f = Nf*2; f < NfNew*2; ++f) tmpf[f] = NULL;
3862764a2aaSMatthew G. Knepley   for (f = Nf*Nf*4; f < NfNew*NfNew*4; ++f) tmpg[f] = NULL;
387475e0ac9SMatthew G. Knepley   for (f = Nf*Nf*4; f < NfNew*NfNew*4; ++f) tmpgp[f] = NULL;
388b7e05686SMatthew G. Knepley   for (f = Nf*Nf*4; f < NfNew*NfNew*4; ++f) tmpgt[f] = NULL;
3890c2f2876SMatthew G. Knepley   for (f = Nf; f < NfNew; ++f) tmpr[f] = NULL;
3900c2f2876SMatthew G. Knepley   for (f = Nf; f < NfNew; ++f) tmpctx[f] = NULL;
391b7e05686SMatthew G. Knepley   ierr = PetscFree7(prob->obj, prob->f, prob->g, prob->gp, prob->gt, prob->r, prob->ctx);CHKERRQ(ierr);
3922764a2aaSMatthew G. Knepley   prob->obj = tmpobj;
3932764a2aaSMatthew G. Knepley   prob->f   = tmpf;
3942764a2aaSMatthew G. Knepley   prob->g   = tmpg;
395475e0ac9SMatthew G. Knepley   prob->gp  = tmpgp;
396b7e05686SMatthew G. Knepley   prob->gt  = tmpgt;
3970c2f2876SMatthew G. Knepley   prob->r   = tmpr;
3980c2f2876SMatthew G. Knepley   prob->ctx = tmpctx;
399*c371a6d1SMatthew G. Knepley   ierr = PetscCalloc3(NfNew*2, &tmpfbd, NfNew*NfNew*4, &tmpgbd, NfNew, &tmpexactSol);CHKERRQ(ierr);
4002764a2aaSMatthew G. Knepley   for (f = 0; f < Nf*2; ++f) tmpfbd[f] = prob->fBd[f];
4012764a2aaSMatthew G. Knepley   for (f = 0; f < Nf*Nf*4; ++f) tmpgbd[f] = prob->gBd[f];
402*c371a6d1SMatthew G. Knepley   for (f = 0; f < Nf; ++f) tmpexactSol[f] = prob->exactSol[f];
4032764a2aaSMatthew G. Knepley   for (f = Nf*2; f < NfNew*2; ++f) tmpfbd[f] = NULL;
4042764a2aaSMatthew G. Knepley   for (f = Nf*Nf*4; f < NfNew*NfNew*4; ++f) tmpgbd[f] = NULL;
405*c371a6d1SMatthew G. Knepley   for (f = Nf; f < NfNew; ++f) tmpexactSol[f] = NULL;
406*c371a6d1SMatthew G. Knepley   ierr = PetscFree3(prob->fBd, prob->gBd, prob->exactSol);CHKERRQ(ierr);
4072764a2aaSMatthew G. Knepley   prob->fBd = tmpfbd;
4082764a2aaSMatthew G. Knepley   prob->gBd = tmpgbd;
409*c371a6d1SMatthew G. Knepley   prob->exactSol = tmpexactSol;
4102764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
4112764a2aaSMatthew G. Knepley }
4122764a2aaSMatthew G. Knepley 
4132764a2aaSMatthew G. Knepley /*@
4142764a2aaSMatthew G. Knepley   PetscDSDestroy - Destroys a PetscDS object
4152764a2aaSMatthew G. Knepley 
4162764a2aaSMatthew G. Knepley   Collective on PetscDS
4172764a2aaSMatthew G. Knepley 
4182764a2aaSMatthew G. Knepley   Input Parameter:
4192764a2aaSMatthew G. Knepley . prob - the PetscDS object to destroy
4202764a2aaSMatthew G. Knepley 
4212764a2aaSMatthew G. Knepley   Level: developer
4222764a2aaSMatthew G. Knepley 
4232764a2aaSMatthew G. Knepley .seealso PetscDSView()
4242764a2aaSMatthew G. Knepley @*/
4252764a2aaSMatthew G. Knepley PetscErrorCode PetscDSDestroy(PetscDS *prob)
4262764a2aaSMatthew G. Knepley {
4272764a2aaSMatthew G. Knepley   PetscInt       f;
42858ebd649SToby Isaac   DSBoundary     next;
4292764a2aaSMatthew G. Knepley   PetscErrorCode ierr;
4302764a2aaSMatthew G. Knepley 
4312764a2aaSMatthew G. Knepley   PetscFunctionBegin;
4322764a2aaSMatthew G. Knepley   if (!*prob) PetscFunctionReturn(0);
4332764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific((*prob), PETSCDS_CLASSID, 1);
4342764a2aaSMatthew G. Knepley 
4352764a2aaSMatthew G. Knepley   if (--((PetscObject)(*prob))->refct > 0) {*prob = 0; PetscFunctionReturn(0);}
4362764a2aaSMatthew G. Knepley   ((PetscObject) (*prob))->refct = 0;
4372764a2aaSMatthew G. Knepley   ierr = PetscDSDestroyStructs_Static(*prob);CHKERRQ(ierr);
4382764a2aaSMatthew G. Knepley   for (f = 0; f < (*prob)->Nf; ++f) {
4392764a2aaSMatthew G. Knepley     ierr = PetscObjectDereference((*prob)->disc[f]);CHKERRQ(ierr);
4402764a2aaSMatthew G. Knepley   }
441f744cafaSSander Arens   ierr = PetscFree3((*prob)->disc, (*prob)->implicit, (*prob)->adjacency);CHKERRQ(ierr);
442b7e05686SMatthew G. Knepley   ierr = PetscFree7((*prob)->obj,(*prob)->f,(*prob)->g,(*prob)->gp,(*prob)->gt,(*prob)->r,(*prob)->ctx);CHKERRQ(ierr);
443*c371a6d1SMatthew G. Knepley   ierr = PetscFree3((*prob)->fBd,(*prob)->gBd,(*prob)->exactSol);CHKERRQ(ierr);
4442764a2aaSMatthew G. Knepley   if ((*prob)->ops->destroy) {ierr = (*(*prob)->ops->destroy)(*prob);CHKERRQ(ierr);}
44558ebd649SToby Isaac   next = (*prob)->boundary;
44658ebd649SToby Isaac   while (next) {
44758ebd649SToby Isaac     DSBoundary b = next;
44858ebd649SToby Isaac 
44958ebd649SToby Isaac     next = b->next;
45058ebd649SToby Isaac     ierr = PetscFree(b->comps);CHKERRQ(ierr);
45158ebd649SToby Isaac     ierr = PetscFree(b->ids);CHKERRQ(ierr);
45258ebd649SToby Isaac     ierr = PetscFree(b->name);CHKERRQ(ierr);
45358ebd649SToby Isaac     ierr = PetscFree(b->labelname);CHKERRQ(ierr);
45458ebd649SToby Isaac     ierr = PetscFree(b);CHKERRQ(ierr);
45558ebd649SToby Isaac   }
4562764a2aaSMatthew G. Knepley   ierr = PetscHeaderDestroy(prob);CHKERRQ(ierr);
4572764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
4582764a2aaSMatthew G. Knepley }
4592764a2aaSMatthew G. Knepley 
4602764a2aaSMatthew G. Knepley /*@
4612764a2aaSMatthew G. Knepley   PetscDSCreate - Creates an empty PetscDS object. The type can then be set with PetscDSSetType().
4622764a2aaSMatthew G. Knepley 
4632764a2aaSMatthew G. Knepley   Collective on MPI_Comm
4642764a2aaSMatthew G. Knepley 
4652764a2aaSMatthew G. Knepley   Input Parameter:
4662764a2aaSMatthew G. Knepley . comm - The communicator for the PetscDS object
4672764a2aaSMatthew G. Knepley 
4682764a2aaSMatthew G. Knepley   Output Parameter:
4692764a2aaSMatthew G. Knepley . prob - The PetscDS object
4702764a2aaSMatthew G. Knepley 
4712764a2aaSMatthew G. Knepley   Level: beginner
4722764a2aaSMatthew G. Knepley 
4732764a2aaSMatthew G. Knepley .seealso: PetscDSSetType(), PETSCDSBASIC
4742764a2aaSMatthew G. Knepley @*/
4752764a2aaSMatthew G. Knepley PetscErrorCode PetscDSCreate(MPI_Comm comm, PetscDS *prob)
4762764a2aaSMatthew G. Knepley {
4772764a2aaSMatthew G. Knepley   PetscDS   p;
4782764a2aaSMatthew G. Knepley   PetscErrorCode ierr;
4792764a2aaSMatthew G. Knepley 
4802764a2aaSMatthew G. Knepley   PetscFunctionBegin;
4812764a2aaSMatthew G. Knepley   PetscValidPointer(prob, 2);
4822764a2aaSMatthew G. Knepley   *prob  = NULL;
4832764a2aaSMatthew G. Knepley   ierr = PetscDSInitializePackage();CHKERRQ(ierr);
4842764a2aaSMatthew G. Knepley 
48573107ff1SLisandro Dalcin   ierr = PetscHeaderCreate(p, PETSCDS_CLASSID, "PetscDS", "Discrete System", "PetscDS", comm, PetscDSDestroy, PetscDSView);CHKERRQ(ierr);
4862764a2aaSMatthew G. Knepley 
4872764a2aaSMatthew G. Knepley   p->Nf    = 0;
4882764a2aaSMatthew G. Knepley   p->setup = PETSC_FALSE;
4892764a2aaSMatthew G. Knepley 
4902764a2aaSMatthew G. Knepley   *prob = p;
4912764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
4922764a2aaSMatthew G. Knepley }
4932764a2aaSMatthew G. Knepley 
494bc4ae4beSMatthew G. Knepley /*@
495bc4ae4beSMatthew G. Knepley   PetscDSGetNumFields - Returns the number of fields in the DS
496bc4ae4beSMatthew G. Knepley 
497bc4ae4beSMatthew G. Knepley   Not collective
498bc4ae4beSMatthew G. Knepley 
499bc4ae4beSMatthew G. Knepley   Input Parameter:
500bc4ae4beSMatthew G. Knepley . prob - The PetscDS object
501bc4ae4beSMatthew G. Knepley 
502bc4ae4beSMatthew G. Knepley   Output Parameter:
503bc4ae4beSMatthew G. Knepley . Nf - The number of fields
504bc4ae4beSMatthew G. Knepley 
505bc4ae4beSMatthew G. Knepley   Level: beginner
506bc4ae4beSMatthew G. Knepley 
507bc4ae4beSMatthew G. Knepley .seealso: PetscDSGetSpatialDimension(), PetscDSCreate()
508bc4ae4beSMatthew G. Knepley @*/
5092764a2aaSMatthew G. Knepley PetscErrorCode PetscDSGetNumFields(PetscDS prob, PetscInt *Nf)
5102764a2aaSMatthew G. Knepley {
5112764a2aaSMatthew G. Knepley   PetscFunctionBegin;
5122764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
5132764a2aaSMatthew G. Knepley   PetscValidPointer(Nf, 2);
5142764a2aaSMatthew G. Knepley   *Nf = prob->Nf;
5152764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
5162764a2aaSMatthew G. Knepley }
5172764a2aaSMatthew G. Knepley 
518bc4ae4beSMatthew G. Knepley /*@
519bc4ae4beSMatthew G. Knepley   PetscDSGetSpatialDimension - Returns the spatial dimension of the DS
520bc4ae4beSMatthew G. Knepley 
521bc4ae4beSMatthew G. Knepley   Not collective
522bc4ae4beSMatthew G. Knepley 
523bc4ae4beSMatthew G. Knepley   Input Parameter:
524bc4ae4beSMatthew G. Knepley . prob - The PetscDS object
525bc4ae4beSMatthew G. Knepley 
526bc4ae4beSMatthew G. Knepley   Output Parameter:
527bc4ae4beSMatthew G. Knepley . dim - The spatial dimension
528bc4ae4beSMatthew G. Knepley 
529bc4ae4beSMatthew G. Knepley   Level: beginner
530bc4ae4beSMatthew G. Knepley 
531bc4ae4beSMatthew G. Knepley .seealso: PetscDSGetNumFields(), PetscDSCreate()
532bc4ae4beSMatthew G. Knepley @*/
5332764a2aaSMatthew G. Knepley PetscErrorCode PetscDSGetSpatialDimension(PetscDS prob, PetscInt *dim)
5342764a2aaSMatthew G. Knepley {
5352764a2aaSMatthew G. Knepley   PetscErrorCode ierr;
5362764a2aaSMatthew G. Knepley 
5372764a2aaSMatthew G. Knepley   PetscFunctionBegin;
5382764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
5392764a2aaSMatthew G. Knepley   PetscValidPointer(dim, 2);
5402764a2aaSMatthew G. Knepley   *dim = 0;
5419de99aefSMatthew G. Knepley   if (prob->Nf) {
5429de99aefSMatthew G. Knepley     PetscObject  obj;
5439de99aefSMatthew G. Knepley     PetscClassId id;
5449de99aefSMatthew G. Knepley 
5459de99aefSMatthew G. Knepley     ierr = PetscDSGetDiscretization(prob, 0, &obj);CHKERRQ(ierr);
5469de99aefSMatthew G. Knepley     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
5479de99aefSMatthew G. Knepley     if (id == PETSCFE_CLASSID)      {ierr = PetscFEGetSpatialDimension((PetscFE) obj, dim);CHKERRQ(ierr);}
5489de99aefSMatthew G. Knepley     else if (id == PETSCFV_CLASSID) {ierr = PetscFVGetSpatialDimension((PetscFV) obj, dim);CHKERRQ(ierr);}
5499de99aefSMatthew G. Knepley     else SETERRQ1(PetscObjectComm((PetscObject) prob), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", 0);
5509de99aefSMatthew G. Knepley   }
5512764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
5522764a2aaSMatthew G. Knepley }
5532764a2aaSMatthew G. Knepley 
554bc4ae4beSMatthew G. Knepley /*@
555bc4ae4beSMatthew G. Knepley   PetscDSGetTotalDimension - Returns the total size of the approximation space for this system
556bc4ae4beSMatthew G. Knepley 
557bc4ae4beSMatthew G. Knepley   Not collective
558bc4ae4beSMatthew G. Knepley 
559bc4ae4beSMatthew G. Knepley   Input Parameter:
560bc4ae4beSMatthew G. Knepley . prob - The PetscDS object
561bc4ae4beSMatthew G. Knepley 
562bc4ae4beSMatthew G. Knepley   Output Parameter:
563bc4ae4beSMatthew G. Knepley . dim - The total problem dimension
564bc4ae4beSMatthew G. Knepley 
565bc4ae4beSMatthew G. Knepley   Level: beginner
566bc4ae4beSMatthew G. Knepley 
567bc4ae4beSMatthew G. Knepley .seealso: PetscDSGetNumFields(), PetscDSCreate()
568bc4ae4beSMatthew G. Knepley @*/
5692764a2aaSMatthew G. Knepley PetscErrorCode PetscDSGetTotalDimension(PetscDS prob, PetscInt *dim)
5702764a2aaSMatthew G. Knepley {
5712764a2aaSMatthew G. Knepley   PetscErrorCode ierr;
5722764a2aaSMatthew G. Knepley 
5732764a2aaSMatthew G. Knepley   PetscFunctionBegin;
5742764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
5752764a2aaSMatthew G. Knepley   ierr = PetscDSSetUp(prob);CHKERRQ(ierr);
5762764a2aaSMatthew G. Knepley   PetscValidPointer(dim, 2);
5772764a2aaSMatthew G. Knepley   *dim = prob->totDim;
5782764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
5792764a2aaSMatthew G. Knepley }
5802764a2aaSMatthew G. Knepley 
581bc4ae4beSMatthew G. Knepley /*@
582bc4ae4beSMatthew G. Knepley   PetscDSGetTotalComponents - Returns the total number of components in this system
583bc4ae4beSMatthew G. Knepley 
584bc4ae4beSMatthew G. Knepley   Not collective
585bc4ae4beSMatthew G. Knepley 
586bc4ae4beSMatthew G. Knepley   Input Parameter:
587bc4ae4beSMatthew G. Knepley . prob - The PetscDS object
588bc4ae4beSMatthew G. Knepley 
589bc4ae4beSMatthew G. Knepley   Output Parameter:
590bc4ae4beSMatthew G. Knepley . dim - The total number of components
591bc4ae4beSMatthew G. Knepley 
592bc4ae4beSMatthew G. Knepley   Level: beginner
593bc4ae4beSMatthew G. Knepley 
594bc4ae4beSMatthew G. Knepley .seealso: PetscDSGetNumFields(), PetscDSCreate()
595bc4ae4beSMatthew G. Knepley @*/
5962764a2aaSMatthew G. Knepley PetscErrorCode PetscDSGetTotalComponents(PetscDS prob, PetscInt *Nc)
5972764a2aaSMatthew G. Knepley {
5982764a2aaSMatthew G. Knepley   PetscErrorCode ierr;
5992764a2aaSMatthew G. Knepley 
6002764a2aaSMatthew G. Knepley   PetscFunctionBegin;
6012764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
6022764a2aaSMatthew G. Knepley   ierr = PetscDSSetUp(prob);CHKERRQ(ierr);
6032764a2aaSMatthew G. Knepley   PetscValidPointer(Nc, 2);
6042764a2aaSMatthew G. Knepley   *Nc = prob->totComp;
6052764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
6062764a2aaSMatthew G. Knepley }
6072764a2aaSMatthew G. Knepley 
608bc4ae4beSMatthew G. Knepley /*@
609bc4ae4beSMatthew G. Knepley   PetscDSGetDiscretization - Returns the discretization object for the given field
610bc4ae4beSMatthew G. Knepley 
611bc4ae4beSMatthew G. Knepley   Not collective
612bc4ae4beSMatthew G. Knepley 
613bc4ae4beSMatthew G. Knepley   Input Parameters:
614bc4ae4beSMatthew G. Knepley + prob - The PetscDS object
615bc4ae4beSMatthew G. Knepley - f - The field number
616bc4ae4beSMatthew G. Knepley 
617bc4ae4beSMatthew G. Knepley   Output Parameter:
618bc4ae4beSMatthew G. Knepley . disc - The discretization object
619bc4ae4beSMatthew G. Knepley 
620bc4ae4beSMatthew G. Knepley   Level: beginner
621bc4ae4beSMatthew G. Knepley 
622f744cafaSSander Arens .seealso: PetscDSSetDiscretization(), PetscDSAddDiscretization(), PetscDSGetNumFields(), PetscDSCreate()
623bc4ae4beSMatthew G. Knepley @*/
6242764a2aaSMatthew G. Knepley PetscErrorCode PetscDSGetDiscretization(PetscDS prob, PetscInt f, PetscObject *disc)
6252764a2aaSMatthew G. Knepley {
6262764a2aaSMatthew G. Knepley   PetscFunctionBegin;
6272764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
6282764a2aaSMatthew G. Knepley   PetscValidPointer(disc, 3);
6292764a2aaSMatthew 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);
6302764a2aaSMatthew G. Knepley   *disc = prob->disc[f];
6312764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
6322764a2aaSMatthew G. Knepley }
6332764a2aaSMatthew G. Knepley 
634bc4ae4beSMatthew G. Knepley /*@
635bc4ae4beSMatthew G. Knepley   PetscDSSetDiscretization - Sets the discretization object for the given field
636bc4ae4beSMatthew G. Knepley 
637bc4ae4beSMatthew G. Knepley   Not collective
638bc4ae4beSMatthew G. Knepley 
639bc4ae4beSMatthew G. Knepley   Input Parameters:
640bc4ae4beSMatthew G. Knepley + prob - The PetscDS object
641bc4ae4beSMatthew G. Knepley . f - The field number
642bc4ae4beSMatthew G. Knepley - disc - The discretization object
643bc4ae4beSMatthew G. Knepley 
644bc4ae4beSMatthew G. Knepley   Level: beginner
645bc4ae4beSMatthew G. Knepley 
646bc4ae4beSMatthew G. Knepley .seealso: PetscDSGetDiscretization(), PetscDSAddDiscretization(), PetscDSGetNumFields(), PetscDSCreate()
647bc4ae4beSMatthew G. Knepley @*/
6482764a2aaSMatthew G. Knepley PetscErrorCode PetscDSSetDiscretization(PetscDS prob, PetscInt f, PetscObject disc)
6492764a2aaSMatthew G. Knepley {
6502764a2aaSMatthew G. Knepley   PetscErrorCode ierr;
6512764a2aaSMatthew G. Knepley 
6522764a2aaSMatthew G. Knepley   PetscFunctionBegin;
6532764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
6542764a2aaSMatthew G. Knepley   PetscValidPointer(disc, 3);
6552764a2aaSMatthew G. Knepley   if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f);
6562764a2aaSMatthew G. Knepley   ierr = PetscDSEnlarge_Static(prob, f+1);CHKERRQ(ierr);
6572764a2aaSMatthew G. Knepley   if (prob->disc[f]) {ierr = PetscObjectDereference(prob->disc[f]);CHKERRQ(ierr);}
6582764a2aaSMatthew G. Knepley   prob->disc[f] = disc;
6592764a2aaSMatthew G. Knepley   ierr = PetscObjectReference(disc);CHKERRQ(ierr);
660249df284SMatthew G. Knepley   {
661249df284SMatthew G. Knepley     PetscClassId id;
662249df284SMatthew G. Knepley 
663249df284SMatthew G. Knepley     ierr = PetscObjectGetClassId(disc, &id);CHKERRQ(ierr);
664a6cbbb48SMatthew G. Knepley     if (id == PETSCFV_CLASSID) {
665a6cbbb48SMatthew G. Knepley       prob->implicit[f]      = PETSC_FALSE;
666a6cbbb48SMatthew G. Knepley       prob->adjacency[f*2+0] = PETSC_TRUE;
667a6cbbb48SMatthew G. Knepley       prob->adjacency[f*2+1] = PETSC_FALSE;
668a6cbbb48SMatthew G. Knepley     }
669249df284SMatthew G. Knepley   }
6702764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
6712764a2aaSMatthew G. Knepley }
6722764a2aaSMatthew G. Knepley 
673bc4ae4beSMatthew G. Knepley /*@
674bc4ae4beSMatthew G. Knepley   PetscDSAddDiscretization - Adds a discretization object
675bc4ae4beSMatthew G. Knepley 
676bc4ae4beSMatthew G. Knepley   Not collective
677bc4ae4beSMatthew G. Knepley 
678bc4ae4beSMatthew G. Knepley   Input Parameters:
679bc4ae4beSMatthew G. Knepley + prob - The PetscDS object
680bc4ae4beSMatthew G. Knepley - disc - The boundary discretization object
681bc4ae4beSMatthew G. Knepley 
682bc4ae4beSMatthew G. Knepley   Level: beginner
683bc4ae4beSMatthew G. Knepley 
684bc4ae4beSMatthew G. Knepley .seealso: PetscDSGetDiscretization(), PetscDSSetDiscretization(), PetscDSGetNumFields(), PetscDSCreate()
685bc4ae4beSMatthew G. Knepley @*/
6862764a2aaSMatthew G. Knepley PetscErrorCode PetscDSAddDiscretization(PetscDS prob, PetscObject disc)
6872764a2aaSMatthew G. Knepley {
6882764a2aaSMatthew G. Knepley   PetscErrorCode ierr;
6892764a2aaSMatthew G. Knepley 
6902764a2aaSMatthew G. Knepley   PetscFunctionBegin;
6912764a2aaSMatthew G. Knepley   ierr = PetscDSSetDiscretization(prob, prob->Nf, disc);CHKERRQ(ierr);
6922764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
6932764a2aaSMatthew G. Knepley }
6942764a2aaSMatthew G. Knepley 
695249df284SMatthew G. Knepley /*@
696249df284SMatthew G. Knepley   PetscDSGetImplicit - Returns the flag for implicit solve for this field. This is just a guide for IMEX
697249df284SMatthew G. Knepley 
698249df284SMatthew G. Knepley   Not collective
699249df284SMatthew G. Knepley 
700249df284SMatthew G. Knepley   Input Parameters:
701249df284SMatthew G. Knepley + prob - The PetscDS object
702249df284SMatthew G. Knepley - f - The field number
703249df284SMatthew G. Knepley 
704249df284SMatthew G. Knepley   Output Parameter:
705249df284SMatthew G. Knepley . implicit - The flag indicating what kind of solve to use for this field
706249df284SMatthew G. Knepley 
707249df284SMatthew G. Knepley   Level: developer
708249df284SMatthew G. Knepley 
709f744cafaSSander Arens .seealso: PetscDSSetImplicit(), PetscDSSetDiscretization(), PetscDSAddDiscretization(), PetscDSGetNumFields(), PetscDSCreate()
710249df284SMatthew G. Knepley @*/
711249df284SMatthew G. Knepley PetscErrorCode PetscDSGetImplicit(PetscDS prob, PetscInt f, PetscBool *implicit)
712249df284SMatthew G. Knepley {
713249df284SMatthew G. Knepley   PetscFunctionBegin;
714249df284SMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
715249df284SMatthew G. Knepley   PetscValidPointer(implicit, 3);
716249df284SMatthew 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);
717249df284SMatthew G. Knepley   *implicit = prob->implicit[f];
718249df284SMatthew G. Knepley   PetscFunctionReturn(0);
719249df284SMatthew G. Knepley }
720249df284SMatthew G. Knepley 
721249df284SMatthew G. Knepley /*@
722249df284SMatthew G. Knepley   PetscDSSetImplicit - Set the flag for implicit solve for this field. This is just a guide for IMEX
723249df284SMatthew G. Knepley 
724249df284SMatthew G. Knepley   Not collective
725249df284SMatthew G. Knepley 
726249df284SMatthew G. Knepley   Input Parameters:
727249df284SMatthew G. Knepley + prob - The PetscDS object
728249df284SMatthew G. Knepley . f - The field number
729249df284SMatthew G. Knepley - implicit - The flag indicating what kind of solve to use for this field
730249df284SMatthew G. Knepley 
731249df284SMatthew G. Knepley   Level: developer
732249df284SMatthew G. Knepley 
733f744cafaSSander Arens .seealso: PetscDSGetImplicit(), PetscDSSetDiscretization(), PetscDSAddDiscretization(), PetscDSGetNumFields(), PetscDSCreate()
734249df284SMatthew G. Knepley @*/
735249df284SMatthew G. Knepley PetscErrorCode PetscDSSetImplicit(PetscDS prob, PetscInt f, PetscBool implicit)
736249df284SMatthew G. Knepley {
737249df284SMatthew G. Knepley   PetscFunctionBegin;
738249df284SMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
739249df284SMatthew 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);
740249df284SMatthew G. Knepley   prob->implicit[f] = implicit;
741249df284SMatthew G. Knepley   PetscFunctionReturn(0);
742249df284SMatthew G. Knepley }
743249df284SMatthew G. Knepley 
744a6cbbb48SMatthew G. Knepley /*@
745a6cbbb48SMatthew G. Knepley   PetscDSGetAdjacency - Returns the flags for determining variable influence
746a6cbbb48SMatthew G. Knepley 
747a6cbbb48SMatthew G. Knepley   Not collective
748a6cbbb48SMatthew G. Knepley 
749a6cbbb48SMatthew G. Knepley   Input Parameters:
750a6cbbb48SMatthew G. Knepley + prob - The PetscDS object
751a6cbbb48SMatthew G. Knepley - f - The field number
752a6cbbb48SMatthew G. Knepley 
753a6cbbb48SMatthew G. Knepley   Output Parameter:
754a6cbbb48SMatthew G. Knepley + useCone    - Flag for variable influence starting with the cone operation
755a6cbbb48SMatthew G. Knepley - useClosure - Flag for variable influence using transitive closure
756a6cbbb48SMatthew G. Knepley 
757a6cbbb48SMatthew G. Knepley   Note: See the discussion in DMPlexGetAdjacencyUseCone() and DMPlexGetAdjacencyUseClosure()
758a6cbbb48SMatthew G. Knepley 
759a6cbbb48SMatthew G. Knepley   Level: developer
760a6cbbb48SMatthew G. Knepley 
761f744cafaSSander Arens .seealso: PetscDSSetAdjacency(), DMPlexGetAdjacencyUseCone(), DMPlexGetAdjacencyUseClosure(), PetscDSSetDiscretization(), PetscDSAddDiscretization(), PetscDSGetNumFields(), PetscDSCreate()
762a6cbbb48SMatthew G. Knepley @*/
763a6cbbb48SMatthew G. Knepley PetscErrorCode PetscDSGetAdjacency(PetscDS prob, PetscInt f, PetscBool *useCone, PetscBool *useClosure)
764a6cbbb48SMatthew G. Knepley {
765a6cbbb48SMatthew G. Knepley   PetscFunctionBegin;
766a6cbbb48SMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
767a6cbbb48SMatthew G. Knepley   PetscValidPointer(useCone, 3);
768a6cbbb48SMatthew G. Knepley   PetscValidPointer(useClosure, 4);
769a6cbbb48SMatthew 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);
770a6cbbb48SMatthew G. Knepley   *useCone    = prob->adjacency[f*2+0];
771a6cbbb48SMatthew G. Knepley   *useClosure = prob->adjacency[f*2+1];
772a6cbbb48SMatthew G. Knepley   PetscFunctionReturn(0);
773a6cbbb48SMatthew G. Knepley }
774a6cbbb48SMatthew G. Knepley 
775a6cbbb48SMatthew G. Knepley /*@
776a6cbbb48SMatthew G. Knepley   PetscDSSetAdjacency - Set the flags for determining variable influence
777a6cbbb48SMatthew G. Knepley 
778a6cbbb48SMatthew G. Knepley   Not collective
779a6cbbb48SMatthew G. Knepley 
780a6cbbb48SMatthew G. Knepley   Input Parameters:
781a6cbbb48SMatthew G. Knepley + prob - The PetscDS object
782a6cbbb48SMatthew G. Knepley . f - The field number
783a6cbbb48SMatthew G. Knepley . useCone    - Flag for variable influence starting with the cone operation
784a6cbbb48SMatthew G. Knepley - useClosure - Flag for variable influence using transitive closure
785a6cbbb48SMatthew G. Knepley 
786a6cbbb48SMatthew G. Knepley   Note: See the discussion in DMPlexGetAdjacencyUseCone() and DMPlexGetAdjacencyUseClosure()
787a6cbbb48SMatthew G. Knepley 
788a6cbbb48SMatthew G. Knepley   Level: developer
789a6cbbb48SMatthew G. Knepley 
790f744cafaSSander Arens .seealso: PetscDSGetAdjacency(), DMPlexGetAdjacencyUseCone(), DMPlexGetAdjacencyUseClosure(), PetscDSSetDiscretization(), PetscDSAddDiscretization(), PetscDSGetNumFields(), PetscDSCreate()
791a6cbbb48SMatthew G. Knepley @*/
792a6cbbb48SMatthew G. Knepley PetscErrorCode PetscDSSetAdjacency(PetscDS prob, PetscInt f, PetscBool useCone, PetscBool useClosure)
793a6cbbb48SMatthew G. Knepley {
794a6cbbb48SMatthew G. Knepley   PetscFunctionBegin;
795a6cbbb48SMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
796a6cbbb48SMatthew 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);
797a6cbbb48SMatthew G. Knepley   prob->adjacency[f*2+0] = useCone;
798a6cbbb48SMatthew G. Knepley   prob->adjacency[f*2+1] = useClosure;
799a6cbbb48SMatthew G. Knepley   PetscFunctionReturn(0);
800a6cbbb48SMatthew G. Knepley }
801a6cbbb48SMatthew G. Knepley 
8022764a2aaSMatthew G. Knepley PetscErrorCode PetscDSGetObjective(PetscDS prob, PetscInt f,
80330b9ff8bSMatthew G. Knepley                                    void (**obj)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
804194d53e6SMatthew G. Knepley                                                 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
805194d53e6SMatthew G. Knepley                                                 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
80630b9ff8bSMatthew G. Knepley                                                 PetscReal t, const PetscReal x[], PetscScalar obj[]))
8072764a2aaSMatthew G. Knepley {
8082764a2aaSMatthew G. Knepley   PetscFunctionBegin;
8092764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
8102764a2aaSMatthew G. Knepley   PetscValidPointer(obj, 2);
8112764a2aaSMatthew 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);
8122764a2aaSMatthew G. Knepley   *obj = prob->obj[f];
8132764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
8142764a2aaSMatthew G. Knepley }
8152764a2aaSMatthew G. Knepley 
8162764a2aaSMatthew G. Knepley PetscErrorCode PetscDSSetObjective(PetscDS prob, PetscInt f,
81730b9ff8bSMatthew G. Knepley                                    void (*obj)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
818194d53e6SMatthew G. Knepley                                                const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
819194d53e6SMatthew G. Knepley                                                const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
82030b9ff8bSMatthew G. Knepley                                                PetscReal t, const PetscReal x[], PetscScalar obj[]))
8212764a2aaSMatthew G. Knepley {
8222764a2aaSMatthew G. Knepley   PetscErrorCode ierr;
8232764a2aaSMatthew G. Knepley 
8242764a2aaSMatthew G. Knepley   PetscFunctionBegin;
8252764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
826de716cbcSToby Isaac   if (obj) PetscValidFunction(obj, 2);
8272764a2aaSMatthew G. Knepley   if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f);
8282764a2aaSMatthew G. Knepley   ierr = PetscDSEnlarge_Static(prob, f+1);CHKERRQ(ierr);
8292764a2aaSMatthew G. Knepley   prob->obj[f] = obj;
8302764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
8312764a2aaSMatthew G. Knepley }
8322764a2aaSMatthew G. Knepley 
833194d53e6SMatthew G. Knepley /*@C
834194d53e6SMatthew G. Knepley   PetscDSGetResidual - Get the pointwise residual function for a given test field
835194d53e6SMatthew G. Knepley 
836194d53e6SMatthew G. Knepley   Not collective
837194d53e6SMatthew G. Knepley 
838194d53e6SMatthew G. Knepley   Input Parameters:
839194d53e6SMatthew G. Knepley + prob - The PetscDS
840194d53e6SMatthew G. Knepley - f    - The test field number
841194d53e6SMatthew G. Knepley 
842194d53e6SMatthew G. Knepley   Output Parameters:
843194d53e6SMatthew G. Knepley + f0 - integrand for the test function term
844194d53e6SMatthew G. Knepley - f1 - integrand for the test function gradient term
845194d53e6SMatthew G. Knepley 
846194d53e6SMatthew G. Knepley   Note: We are using a first order FEM model for the weak form:
847194d53e6SMatthew G. Knepley 
848194d53e6SMatthew G. Knepley   \int_\Omega \phi f_0(u, u_t, \nabla u, x, t) + \nabla\phi \cdot {\vec f}_1(u, u_t, \nabla u, x, t)
849194d53e6SMatthew G. Knepley 
850194d53e6SMatthew G. Knepley The calling sequence for the callbacks f0 and f1 is given by:
851194d53e6SMatthew G. Knepley 
85230b9ff8bSMatthew G. Knepley $ f0(PetscInt dim, PetscInt Nf, PetscInt NfAux,
853194d53e6SMatthew G. Knepley $    const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
854194d53e6SMatthew G. Knepley $    const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
85530b9ff8bSMatthew G. Knepley $    PetscReal t, const PetscReal x[], PetscScalar f0[])
856194d53e6SMatthew G. Knepley 
857194d53e6SMatthew G. Knepley + dim - the spatial dimension
858194d53e6SMatthew G. Knepley . Nf - the number of fields
859194d53e6SMatthew G. Knepley . uOff - the offset into u[] and u_t[] for each field
860194d53e6SMatthew G. Knepley . uOff_x - the offset into u_x[] for each field
861194d53e6SMatthew G. Knepley . u - each field evaluated at the current point
862194d53e6SMatthew G. Knepley . u_t - the time derivative of each field evaluated at the current point
863194d53e6SMatthew G. Knepley . u_x - the gradient of each field evaluated at the current point
864194d53e6SMatthew G. Knepley . aOff - the offset into a[] and a_t[] for each auxiliary field
865194d53e6SMatthew G. Knepley . aOff_x - the offset into a_x[] for each auxiliary field
866194d53e6SMatthew G. Knepley . a - each auxiliary field evaluated at the current point
867194d53e6SMatthew G. Knepley . a_t - the time derivative of each auxiliary field evaluated at the current point
868194d53e6SMatthew G. Knepley . a_x - the gradient of auxiliary each field evaluated at the current point
869194d53e6SMatthew G. Knepley . t - current time
870194d53e6SMatthew G. Knepley . x - coordinates of the current point
871194d53e6SMatthew G. Knepley - f0 - output values at the current point
872194d53e6SMatthew G. Knepley 
873194d53e6SMatthew G. Knepley   Level: intermediate
874194d53e6SMatthew G. Knepley 
875194d53e6SMatthew G. Knepley .seealso: PetscDSSetResidual()
876194d53e6SMatthew G. Knepley @*/
8772764a2aaSMatthew G. Knepley PetscErrorCode PetscDSGetResidual(PetscDS prob, PetscInt f,
87830b9ff8bSMatthew G. Knepley                                   void (**f0)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
879194d53e6SMatthew G. Knepley                                               const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
880194d53e6SMatthew G. Knepley                                               const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
88130b9ff8bSMatthew G. Knepley                                               PetscReal t, const PetscReal x[], PetscScalar f0[]),
88230b9ff8bSMatthew G. Knepley                                   void (**f1)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
883194d53e6SMatthew G. Knepley                                               const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
884194d53e6SMatthew G. Knepley                                               const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
88530b9ff8bSMatthew G. Knepley                                               PetscReal t, const PetscReal x[], PetscScalar f1[]))
8862764a2aaSMatthew G. Knepley {
8872764a2aaSMatthew G. Knepley   PetscFunctionBegin;
8882764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
8892764a2aaSMatthew 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);
8902764a2aaSMatthew G. Knepley   if (f0) {PetscValidPointer(f0, 3); *f0 = prob->f[f*2+0];}
8912764a2aaSMatthew G. Knepley   if (f1) {PetscValidPointer(f1, 4); *f1 = prob->f[f*2+1];}
8922764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
8932764a2aaSMatthew G. Knepley }
8942764a2aaSMatthew G. Knepley 
895194d53e6SMatthew G. Knepley /*@C
896194d53e6SMatthew G. Knepley   PetscDSSetResidual - Set the pointwise residual function for a given test field
897194d53e6SMatthew G. Knepley 
898194d53e6SMatthew G. Knepley   Not collective
899194d53e6SMatthew G. Knepley 
900194d53e6SMatthew G. Knepley   Input Parameters:
901194d53e6SMatthew G. Knepley + prob - The PetscDS
902194d53e6SMatthew G. Knepley . f    - The test field number
903194d53e6SMatthew G. Knepley . f0 - integrand for the test function term
904194d53e6SMatthew G. Knepley - f1 - integrand for the test function gradient term
905194d53e6SMatthew G. Knepley 
906194d53e6SMatthew G. Knepley   Note: We are using a first order FEM model for the weak form:
907194d53e6SMatthew G. Knepley 
908194d53e6SMatthew G. Knepley   \int_\Omega \phi f_0(u, u_t, \nabla u, x, t) + \nabla\phi \cdot {\vec f}_1(u, u_t, \nabla u, x, t)
909194d53e6SMatthew G. Knepley 
910194d53e6SMatthew G. Knepley The calling sequence for the callbacks f0 and f1 is given by:
911194d53e6SMatthew G. Knepley 
91230b9ff8bSMatthew G. Knepley $ f0(PetscInt dim, PetscInt Nf, PetscInt NfAux,
913194d53e6SMatthew G. Knepley $    const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
914194d53e6SMatthew G. Knepley $    const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
91530b9ff8bSMatthew G. Knepley $    PetscReal t, const PetscReal x[], PetscScalar f0[])
916194d53e6SMatthew G. Knepley 
917194d53e6SMatthew G. Knepley + dim - the spatial dimension
918194d53e6SMatthew G. Knepley . Nf - the number of fields
919194d53e6SMatthew G. Knepley . uOff - the offset into u[] and u_t[] for each field
920194d53e6SMatthew G. Knepley . uOff_x - the offset into u_x[] for each field
921194d53e6SMatthew G. Knepley . u - each field evaluated at the current point
922194d53e6SMatthew G. Knepley . u_t - the time derivative of each field evaluated at the current point
923194d53e6SMatthew G. Knepley . u_x - the gradient of each field evaluated at the current point
924194d53e6SMatthew G. Knepley . aOff - the offset into a[] and a_t[] for each auxiliary field
925194d53e6SMatthew G. Knepley . aOff_x - the offset into a_x[] for each auxiliary field
926194d53e6SMatthew G. Knepley . a - each auxiliary field evaluated at the current point
927194d53e6SMatthew G. Knepley . a_t - the time derivative of each auxiliary field evaluated at the current point
928194d53e6SMatthew G. Knepley . a_x - the gradient of auxiliary each field evaluated at the current point
929194d53e6SMatthew G. Knepley . t - current time
930194d53e6SMatthew G. Knepley . x - coordinates of the current point
931194d53e6SMatthew G. Knepley - f0 - output values at the current point
932194d53e6SMatthew G. Knepley 
933194d53e6SMatthew G. Knepley   Level: intermediate
934194d53e6SMatthew G. Knepley 
935194d53e6SMatthew G. Knepley .seealso: PetscDSGetResidual()
936194d53e6SMatthew G. Knepley @*/
9372764a2aaSMatthew G. Knepley PetscErrorCode PetscDSSetResidual(PetscDS prob, PetscInt f,
93830b9ff8bSMatthew G. Knepley                                   void (*f0)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
939194d53e6SMatthew G. Knepley                                              const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
940194d53e6SMatthew G. Knepley                                              const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
94130b9ff8bSMatthew G. Knepley                                              PetscReal t, const PetscReal x[], PetscScalar f0[]),
94230b9ff8bSMatthew G. Knepley                                   void (*f1)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
943194d53e6SMatthew G. Knepley                                              const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
944194d53e6SMatthew G. Knepley                                              const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
94530b9ff8bSMatthew G. Knepley                                              PetscReal t, const PetscReal x[], PetscScalar f1[]))
9462764a2aaSMatthew G. Knepley {
9472764a2aaSMatthew G. Knepley   PetscErrorCode ierr;
9482764a2aaSMatthew G. Knepley 
9492764a2aaSMatthew G. Knepley   PetscFunctionBegin;
9502764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
951f866a1d0SMatthew G. Knepley   if (f0) PetscValidFunction(f0, 3);
952f866a1d0SMatthew G. Knepley   if (f1) PetscValidFunction(f1, 4);
9532764a2aaSMatthew G. Knepley   if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f);
9542764a2aaSMatthew G. Knepley   ierr = PetscDSEnlarge_Static(prob, f+1);CHKERRQ(ierr);
9552764a2aaSMatthew G. Knepley   prob->f[f*2+0] = f0;
9562764a2aaSMatthew G. Knepley   prob->f[f*2+1] = f1;
9572764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
9582764a2aaSMatthew G. Knepley }
9592764a2aaSMatthew G. Knepley 
9603e75805dSMatthew G. Knepley /*@C
9613e75805dSMatthew G. Knepley   PetscDSHasJacobian - Signals that Jacobian functions have been set
9623e75805dSMatthew G. Knepley 
9633e75805dSMatthew G. Knepley   Not collective
9643e75805dSMatthew G. Knepley 
9653e75805dSMatthew G. Knepley   Input Parameter:
9663e75805dSMatthew G. Knepley . prob - The PetscDS
9673e75805dSMatthew G. Knepley 
9683e75805dSMatthew G. Knepley   Output Parameter:
9693e75805dSMatthew G. Knepley . hasJac - flag that pointwise function for the Jacobian has been set
9703e75805dSMatthew G. Knepley 
9713e75805dSMatthew G. Knepley   Level: intermediate
9723e75805dSMatthew G. Knepley 
9733e75805dSMatthew G. Knepley .seealso: PetscDSGetJacobianPreconditioner(), PetscDSSetJacobianPreconditioner(), PetscDSGetJacobian()
9743e75805dSMatthew G. Knepley @*/
9753e75805dSMatthew G. Knepley PetscErrorCode PetscDSHasJacobian(PetscDS prob, PetscBool *hasJac)
9763e75805dSMatthew G. Knepley {
9773e75805dSMatthew G. Knepley   PetscInt f, g, h;
9783e75805dSMatthew G. Knepley 
9793e75805dSMatthew G. Knepley   PetscFunctionBegin;
9803e75805dSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
9813e75805dSMatthew G. Knepley   *hasJac = PETSC_FALSE;
9823e75805dSMatthew G. Knepley   for (f = 0; f < prob->Nf; ++f) {
9833e75805dSMatthew G. Knepley     for (g = 0; g < prob->Nf; ++g) {
9843e75805dSMatthew G. Knepley       for (h = 0; h < 4; ++h) {
9853e75805dSMatthew G. Knepley         if (prob->g[(f*prob->Nf + g)*4+h]) *hasJac = PETSC_TRUE;
9863e75805dSMatthew G. Knepley       }
9873e75805dSMatthew G. Knepley     }
9883e75805dSMatthew G. Knepley   }
9893e75805dSMatthew G. Knepley   PetscFunctionReturn(0);
9903e75805dSMatthew G. Knepley }
9913e75805dSMatthew G. Knepley 
992194d53e6SMatthew G. Knepley /*@C
993194d53e6SMatthew G. Knepley   PetscDSGetJacobian - Get the pointwise Jacobian function for given test and basis field
994194d53e6SMatthew G. Knepley 
995194d53e6SMatthew G. Knepley   Not collective
996194d53e6SMatthew G. Knepley 
997194d53e6SMatthew G. Knepley   Input Parameters:
998194d53e6SMatthew G. Knepley + prob - The PetscDS
999194d53e6SMatthew G. Knepley . f    - The test field number
1000194d53e6SMatthew G. Knepley - g    - The field number
1001194d53e6SMatthew G. Knepley 
1002194d53e6SMatthew G. Knepley   Output Parameters:
1003194d53e6SMatthew G. Knepley + g0 - integrand for the test and basis function term
1004194d53e6SMatthew G. Knepley . g1 - integrand for the test function and basis function gradient term
1005194d53e6SMatthew G. Knepley . g2 - integrand for the test function gradient and basis function term
1006194d53e6SMatthew G. Knepley - g3 - integrand for the test function gradient and basis function gradient term
1007194d53e6SMatthew G. Knepley 
1008194d53e6SMatthew G. Knepley   Note: We are using a first order FEM model for the weak form:
1009194d53e6SMatthew G. Knepley 
1010194d53e6SMatthew G. Knepley   \int_\Omega \phi g_0(u, u_t, \nabla u, x, t) \psi + \phi {\vec g}_1(u, u_t, \nabla u, x, t) \nabla \psi + \nabla\phi \cdot {\vec g}_2(u, u_t, \nabla u, x, t) \psi + \nabla\phi \cdot {\overleftrightarrow g}_3(u, u_t, \nabla u, x, t) \cdot \nabla \psi
1011194d53e6SMatthew G. Knepley 
1012194d53e6SMatthew G. Knepley The calling sequence for the callbacks g0, g1, g2 and g3 is given by:
1013194d53e6SMatthew G. Knepley 
101430b9ff8bSMatthew G. Knepley $ g0(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1015194d53e6SMatthew G. Knepley $    const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1016194d53e6SMatthew G. Knepley $    const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
101730b9ff8bSMatthew G. Knepley $    PetscReal t, const PetscReal u_tShift, const PetscReal x[], PetscScalar g0[])
1018194d53e6SMatthew G. Knepley 
1019194d53e6SMatthew G. Knepley + dim - the spatial dimension
1020194d53e6SMatthew G. Knepley . Nf - the number of fields
1021194d53e6SMatthew G. Knepley . uOff - the offset into u[] and u_t[] for each field
1022194d53e6SMatthew G. Knepley . uOff_x - the offset into u_x[] for each field
1023194d53e6SMatthew G. Knepley . u - each field evaluated at the current point
1024194d53e6SMatthew G. Knepley . u_t - the time derivative of each field evaluated at the current point
1025194d53e6SMatthew G. Knepley . u_x - the gradient of each field evaluated at the current point
1026194d53e6SMatthew G. Knepley . aOff - the offset into a[] and a_t[] for each auxiliary field
1027194d53e6SMatthew G. Knepley . aOff_x - the offset into a_x[] for each auxiliary field
1028194d53e6SMatthew G. Knepley . a - each auxiliary field evaluated at the current point
1029194d53e6SMatthew G. Knepley . a_t - the time derivative of each auxiliary field evaluated at the current point
1030194d53e6SMatthew G. Knepley . a_x - the gradient of auxiliary each field evaluated at the current point
1031194d53e6SMatthew G. Knepley . t - current time
10322aa1fc23SMatthew G. Knepley . u_tShift - the multiplier a for dF/dU_t
1033194d53e6SMatthew G. Knepley . x - coordinates of the current point
1034194d53e6SMatthew G. Knepley - g0 - output values at the current point
1035194d53e6SMatthew G. Knepley 
1036194d53e6SMatthew G. Knepley   Level: intermediate
1037194d53e6SMatthew G. Knepley 
1038194d53e6SMatthew G. Knepley .seealso: PetscDSSetJacobian()
1039194d53e6SMatthew G. Knepley @*/
10402764a2aaSMatthew G. Knepley PetscErrorCode PetscDSGetJacobian(PetscDS prob, PetscInt f, PetscInt g,
104130b9ff8bSMatthew G. Knepley                                   void (**g0)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1042194d53e6SMatthew G. Knepley                                               const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1043194d53e6SMatthew G. Knepley                                               const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
10442aa1fc23SMatthew G. Knepley                                               PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscScalar g0[]),
104530b9ff8bSMatthew G. Knepley                                   void (**g1)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1046194d53e6SMatthew G. Knepley                                               const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1047194d53e6SMatthew G. Knepley                                               const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
10482aa1fc23SMatthew G. Knepley                                               PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscScalar g1[]),
104930b9ff8bSMatthew G. Knepley                                   void (**g2)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1050194d53e6SMatthew G. Knepley                                               const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1051194d53e6SMatthew G. Knepley                                               const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
10522aa1fc23SMatthew G. Knepley                                               PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscScalar g2[]),
105330b9ff8bSMatthew G. Knepley                                   void (**g3)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1054194d53e6SMatthew G. Knepley                                               const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1055194d53e6SMatthew G. Knepley                                               const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
10562aa1fc23SMatthew G. Knepley                                               PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscScalar g3[]))
10572764a2aaSMatthew G. Knepley {
10582764a2aaSMatthew G. Knepley   PetscFunctionBegin;
10592764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
10602764a2aaSMatthew 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);
10612764a2aaSMatthew 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);
10622764a2aaSMatthew G. Knepley   if (g0) {PetscValidPointer(g0, 4); *g0 = prob->g[(f*prob->Nf + g)*4+0];}
10632764a2aaSMatthew G. Knepley   if (g1) {PetscValidPointer(g1, 5); *g1 = prob->g[(f*prob->Nf + g)*4+1];}
10642764a2aaSMatthew G. Knepley   if (g2) {PetscValidPointer(g2, 6); *g2 = prob->g[(f*prob->Nf + g)*4+2];}
10652764a2aaSMatthew G. Knepley   if (g3) {PetscValidPointer(g3, 7); *g3 = prob->g[(f*prob->Nf + g)*4+3];}
10662764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
10672764a2aaSMatthew G. Knepley }
10682764a2aaSMatthew G. Knepley 
1069194d53e6SMatthew G. Knepley /*@C
1070194d53e6SMatthew G. Knepley   PetscDSSetJacobian - Set the pointwise Jacobian function for given test and basis fields
1071194d53e6SMatthew G. Knepley 
1072194d53e6SMatthew G. Knepley   Not collective
1073194d53e6SMatthew G. Knepley 
1074194d53e6SMatthew G. Knepley   Input Parameters:
1075194d53e6SMatthew G. Knepley + prob - The PetscDS
1076194d53e6SMatthew G. Knepley . f    - The test field number
1077194d53e6SMatthew G. Knepley . g    - The field number
1078194d53e6SMatthew G. Knepley . g0 - integrand for the test and basis function term
1079194d53e6SMatthew G. Knepley . g1 - integrand for the test function and basis function gradient term
1080194d53e6SMatthew G. Knepley . g2 - integrand for the test function gradient and basis function term
1081194d53e6SMatthew G. Knepley - g3 - integrand for the test function gradient and basis function gradient term
1082194d53e6SMatthew G. Knepley 
1083194d53e6SMatthew G. Knepley   Note: We are using a first order FEM model for the weak form:
1084194d53e6SMatthew G. Knepley 
1085194d53e6SMatthew G. Knepley   \int_\Omega \phi g_0(u, u_t, \nabla u, x, t) \psi + \phi {\vec g}_1(u, u_t, \nabla u, x, t) \nabla \psi + \nabla\phi \cdot {\vec g}_2(u, u_t, \nabla u, x, t) \psi + \nabla\phi \cdot {\overleftrightarrow g}_3(u, u_t, \nabla u, x, t) \cdot \nabla \psi
1086194d53e6SMatthew G. Knepley 
1087194d53e6SMatthew G. Knepley The calling sequence for the callbacks g0, g1, g2 and g3 is given by:
1088194d53e6SMatthew G. Knepley 
108930b9ff8bSMatthew G. Knepley $ g0(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1090194d53e6SMatthew G. Knepley $    const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1091194d53e6SMatthew G. Knepley $    const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
109230b9ff8bSMatthew G. Knepley $    PetscReal t, const PetscReal x[], PetscScalar g0[])
1093194d53e6SMatthew G. Knepley 
1094194d53e6SMatthew G. Knepley + dim - the spatial dimension
1095194d53e6SMatthew G. Knepley . Nf - the number of fields
1096194d53e6SMatthew G. Knepley . uOff - the offset into u[] and u_t[] for each field
1097194d53e6SMatthew G. Knepley . uOff_x - the offset into u_x[] for each field
1098194d53e6SMatthew G. Knepley . u - each field evaluated at the current point
1099194d53e6SMatthew G. Knepley . u_t - the time derivative of each field evaluated at the current point
1100194d53e6SMatthew G. Knepley . u_x - the gradient of each field evaluated at the current point
1101194d53e6SMatthew G. Knepley . aOff - the offset into a[] and a_t[] for each auxiliary field
1102194d53e6SMatthew G. Knepley . aOff_x - the offset into a_x[] for each auxiliary field
1103194d53e6SMatthew G. Knepley . a - each auxiliary field evaluated at the current point
1104194d53e6SMatthew G. Knepley . a_t - the time derivative of each auxiliary field evaluated at the current point
1105194d53e6SMatthew G. Knepley . a_x - the gradient of auxiliary each field evaluated at the current point
1106194d53e6SMatthew G. Knepley . t - current time
11072aa1fc23SMatthew G. Knepley . u_tShift - the multiplier a for dF/dU_t
1108194d53e6SMatthew G. Knepley . x - coordinates of the current point
1109194d53e6SMatthew G. Knepley - g0 - output values at the current point
1110194d53e6SMatthew G. Knepley 
1111194d53e6SMatthew G. Knepley   Level: intermediate
1112194d53e6SMatthew G. Knepley 
1113194d53e6SMatthew G. Knepley .seealso: PetscDSGetJacobian()
1114194d53e6SMatthew G. Knepley @*/
11152764a2aaSMatthew G. Knepley PetscErrorCode PetscDSSetJacobian(PetscDS prob, PetscInt f, PetscInt g,
111630b9ff8bSMatthew G. Knepley                                   void (*g0)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1117194d53e6SMatthew G. Knepley                                              const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1118194d53e6SMatthew G. Knepley                                              const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
111930b9ff8bSMatthew G. Knepley                                              PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscScalar g0[]),
112030b9ff8bSMatthew G. Knepley                                   void (*g1)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1121194d53e6SMatthew G. Knepley                                              const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1122194d53e6SMatthew G. Knepley                                              const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
112330b9ff8bSMatthew G. Knepley                                              PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscScalar g1[]),
112430b9ff8bSMatthew G. Knepley                                   void (*g2)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1125194d53e6SMatthew G. Knepley                                              const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1126194d53e6SMatthew G. Knepley                                              const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
112730b9ff8bSMatthew G. Knepley                                              PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscScalar g2[]),
112830b9ff8bSMatthew G. Knepley                                   void (*g3)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1129194d53e6SMatthew G. Knepley                                              const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1130194d53e6SMatthew G. Knepley                                              const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
113130b9ff8bSMatthew G. Knepley                                              PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscScalar g3[]))
11322764a2aaSMatthew G. Knepley {
11332764a2aaSMatthew G. Knepley   PetscErrorCode ierr;
11342764a2aaSMatthew G. Knepley 
11352764a2aaSMatthew G. Knepley   PetscFunctionBegin;
11362764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
11372764a2aaSMatthew G. Knepley   if (g0) PetscValidFunction(g0, 4);
11382764a2aaSMatthew G. Knepley   if (g1) PetscValidFunction(g1, 5);
11392764a2aaSMatthew G. Knepley   if (g2) PetscValidFunction(g2, 6);
11402764a2aaSMatthew G. Knepley   if (g3) PetscValidFunction(g3, 7);
11412764a2aaSMatthew G. Knepley   if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f);
11422764a2aaSMatthew G. Knepley   if (g < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", g);
11432764a2aaSMatthew G. Knepley   ierr = PetscDSEnlarge_Static(prob, PetscMax(f, g)+1);CHKERRQ(ierr);
11442764a2aaSMatthew G. Knepley   prob->g[(f*prob->Nf + g)*4+0] = g0;
11452764a2aaSMatthew G. Knepley   prob->g[(f*prob->Nf + g)*4+1] = g1;
11462764a2aaSMatthew G. Knepley   prob->g[(f*prob->Nf + g)*4+2] = g2;
11472764a2aaSMatthew G. Knepley   prob->g[(f*prob->Nf + g)*4+3] = g3;
11482764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
11492764a2aaSMatthew G. Knepley }
11502764a2aaSMatthew G. Knepley 
1151475e0ac9SMatthew G. Knepley /*@C
1152475e0ac9SMatthew G. Knepley   PetscDSHasJacobianPreconditioner - Signals that a Jacobian preconditioner matrix has been set
1153475e0ac9SMatthew G. Knepley 
1154475e0ac9SMatthew G. Knepley   Not collective
1155475e0ac9SMatthew G. Knepley 
1156475e0ac9SMatthew G. Knepley   Input Parameter:
1157475e0ac9SMatthew G. Knepley . prob - The PetscDS
1158475e0ac9SMatthew G. Knepley 
1159475e0ac9SMatthew G. Knepley   Output Parameter:
1160475e0ac9SMatthew G. Knepley . hasJacPre - flag that pointwise function for Jacobian preconditioner matrix has been set
1161475e0ac9SMatthew G. Knepley 
1162475e0ac9SMatthew G. Knepley   Level: intermediate
1163475e0ac9SMatthew G. Knepley 
1164475e0ac9SMatthew G. Knepley .seealso: PetscDSGetJacobianPreconditioner(), PetscDSSetJacobianPreconditioner(), PetscDSGetJacobian()
1165475e0ac9SMatthew G. Knepley @*/
1166475e0ac9SMatthew G. Knepley PetscErrorCode PetscDSHasJacobianPreconditioner(PetscDS prob, PetscBool *hasJacPre)
1167475e0ac9SMatthew G. Knepley {
1168475e0ac9SMatthew G. Knepley   PetscInt f, g, h;
1169475e0ac9SMatthew G. Knepley 
1170475e0ac9SMatthew G. Knepley   PetscFunctionBegin;
1171475e0ac9SMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1172475e0ac9SMatthew G. Knepley   *hasJacPre = PETSC_FALSE;
1173475e0ac9SMatthew G. Knepley   for (f = 0; f < prob->Nf; ++f) {
1174475e0ac9SMatthew G. Knepley     for (g = 0; g < prob->Nf; ++g) {
1175475e0ac9SMatthew G. Knepley       for (h = 0; h < 4; ++h) {
1176475e0ac9SMatthew G. Knepley         if (prob->gp[(f*prob->Nf + g)*4+h]) *hasJacPre = PETSC_TRUE;
1177475e0ac9SMatthew G. Knepley       }
1178475e0ac9SMatthew G. Knepley     }
1179475e0ac9SMatthew G. Knepley   }
1180475e0ac9SMatthew G. Knepley   PetscFunctionReturn(0);
1181475e0ac9SMatthew G. Knepley }
1182475e0ac9SMatthew G. Knepley 
1183475e0ac9SMatthew G. Knepley /*@C
1184475e0ac9SMatthew G. Knepley   PetscDSGetJacobianPreconditioner - Get the pointwise Jacobian preconditioner function for given test and basis field. If this is missing, the system matrix is used to build the preconditioner.
1185475e0ac9SMatthew G. Knepley 
1186475e0ac9SMatthew G. Knepley   Not collective
1187475e0ac9SMatthew G. Knepley 
1188475e0ac9SMatthew G. Knepley   Input Parameters:
1189475e0ac9SMatthew G. Knepley + prob - The PetscDS
1190475e0ac9SMatthew G. Knepley . f    - The test field number
1191475e0ac9SMatthew G. Knepley - g    - The field number
1192475e0ac9SMatthew G. Knepley 
1193475e0ac9SMatthew G. Knepley   Output Parameters:
1194475e0ac9SMatthew G. Knepley + g0 - integrand for the test and basis function term
1195475e0ac9SMatthew G. Knepley . g1 - integrand for the test function and basis function gradient term
1196475e0ac9SMatthew G. Knepley . g2 - integrand for the test function gradient and basis function term
1197475e0ac9SMatthew G. Knepley - g3 - integrand for the test function gradient and basis function gradient term
1198475e0ac9SMatthew G. Knepley 
1199475e0ac9SMatthew G. Knepley   Note: We are using a first order FEM model for the weak form:
1200475e0ac9SMatthew G. Knepley 
1201475e0ac9SMatthew G. Knepley   \int_\Omega \phi g_0(u, u_t, \nabla u, x, t) \psi + \phi {\vec g}_1(u, u_t, \nabla u, x, t) \nabla \psi + \nabla\phi \cdot {\vec g}_2(u, u_t, \nabla u, x, t) \psi + \nabla\phi \cdot {\overleftrightarrow g}_3(u, u_t, \nabla u, x, t) \cdot \nabla \psi
1202475e0ac9SMatthew G. Knepley 
1203475e0ac9SMatthew G. Knepley The calling sequence for the callbacks g0, g1, g2 and g3 is given by:
1204475e0ac9SMatthew G. Knepley 
1205475e0ac9SMatthew G. Knepley $ g0(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1206475e0ac9SMatthew G. Knepley $    const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1207475e0ac9SMatthew G. Knepley $    const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1208475e0ac9SMatthew G. Knepley $    PetscReal t, const PetscReal u_tShift, const PetscReal x[], PetscScalar g0[])
1209475e0ac9SMatthew G. Knepley 
1210475e0ac9SMatthew G. Knepley + dim - the spatial dimension
1211475e0ac9SMatthew G. Knepley . Nf - the number of fields
1212475e0ac9SMatthew G. Knepley . uOff - the offset into u[] and u_t[] for each field
1213475e0ac9SMatthew G. Knepley . uOff_x - the offset into u_x[] for each field
1214475e0ac9SMatthew G. Knepley . u - each field evaluated at the current point
1215475e0ac9SMatthew G. Knepley . u_t - the time derivative of each field evaluated at the current point
1216475e0ac9SMatthew G. Knepley . u_x - the gradient of each field evaluated at the current point
1217475e0ac9SMatthew G. Knepley . aOff - the offset into a[] and a_t[] for each auxiliary field
1218475e0ac9SMatthew G. Knepley . aOff_x - the offset into a_x[] for each auxiliary field
1219475e0ac9SMatthew G. Knepley . a - each auxiliary field evaluated at the current point
1220475e0ac9SMatthew G. Knepley . a_t - the time derivative of each auxiliary field evaluated at the current point
1221475e0ac9SMatthew G. Knepley . a_x - the gradient of auxiliary each field evaluated at the current point
1222475e0ac9SMatthew G. Knepley . t - current time
1223475e0ac9SMatthew G. Knepley . u_tShift - the multiplier a for dF/dU_t
1224475e0ac9SMatthew G. Knepley . x - coordinates of the current point
1225475e0ac9SMatthew G. Knepley - g0 - output values at the current point
1226475e0ac9SMatthew G. Knepley 
1227475e0ac9SMatthew G. Knepley   Level: intermediate
1228475e0ac9SMatthew G. Knepley 
1229475e0ac9SMatthew G. Knepley .seealso: PetscDSSetJacobianPreconditioner(), PetscDSGetJacobian()
1230475e0ac9SMatthew G. Knepley @*/
1231475e0ac9SMatthew G. Knepley PetscErrorCode PetscDSGetJacobianPreconditioner(PetscDS prob, PetscInt f, PetscInt g,
1232475e0ac9SMatthew G. Knepley                                   void (**g0)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1233475e0ac9SMatthew G. Knepley                                               const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1234475e0ac9SMatthew G. Knepley                                               const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1235475e0ac9SMatthew G. Knepley                                               PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscScalar g0[]),
1236475e0ac9SMatthew G. Knepley                                   void (**g1)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1237475e0ac9SMatthew G. Knepley                                               const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1238475e0ac9SMatthew G. Knepley                                               const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1239475e0ac9SMatthew G. Knepley                                               PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscScalar g1[]),
1240475e0ac9SMatthew G. Knepley                                   void (**g2)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1241475e0ac9SMatthew G. Knepley                                               const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1242475e0ac9SMatthew G. Knepley                                               const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1243475e0ac9SMatthew G. Knepley                                               PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscScalar g2[]),
1244475e0ac9SMatthew G. Knepley                                   void (**g3)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1245475e0ac9SMatthew G. Knepley                                               const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1246475e0ac9SMatthew G. Knepley                                               const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1247475e0ac9SMatthew G. Knepley                                               PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscScalar g3[]))
1248475e0ac9SMatthew G. Knepley {
1249475e0ac9SMatthew G. Knepley   PetscFunctionBegin;
1250475e0ac9SMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1251475e0ac9SMatthew 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);
1252475e0ac9SMatthew 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);
1253475e0ac9SMatthew G. Knepley   if (g0) {PetscValidPointer(g0, 4); *g0 = prob->gp[(f*prob->Nf + g)*4+0];}
1254475e0ac9SMatthew G. Knepley   if (g1) {PetscValidPointer(g1, 5); *g1 = prob->gp[(f*prob->Nf + g)*4+1];}
1255475e0ac9SMatthew G. Knepley   if (g2) {PetscValidPointer(g2, 6); *g2 = prob->gp[(f*prob->Nf + g)*4+2];}
1256475e0ac9SMatthew G. Knepley   if (g3) {PetscValidPointer(g3, 7); *g3 = prob->gp[(f*prob->Nf + g)*4+3];}
1257475e0ac9SMatthew G. Knepley   PetscFunctionReturn(0);
1258475e0ac9SMatthew G. Knepley }
1259475e0ac9SMatthew G. Knepley 
1260475e0ac9SMatthew G. Knepley /*@C
1261475e0ac9SMatthew G. Knepley   PetscDSSetJacobianPreconditioner - Set the pointwise Jacobian preconditioner function for given test and basis fields. If this is missing, the system matrix is used to build the preconditioner.
1262475e0ac9SMatthew G. Knepley 
1263475e0ac9SMatthew G. Knepley   Not collective
1264475e0ac9SMatthew G. Knepley 
1265475e0ac9SMatthew G. Knepley   Input Parameters:
1266475e0ac9SMatthew G. Knepley + prob - The PetscDS
1267475e0ac9SMatthew G. Knepley . f    - The test field number
1268475e0ac9SMatthew G. Knepley . g    - The field number
1269475e0ac9SMatthew G. Knepley . g0 - integrand for the test and basis function term
1270475e0ac9SMatthew G. Knepley . g1 - integrand for the test function and basis function gradient term
1271475e0ac9SMatthew G. Knepley . g2 - integrand for the test function gradient and basis function term
1272475e0ac9SMatthew G. Knepley - g3 - integrand for the test function gradient and basis function gradient term
1273475e0ac9SMatthew G. Knepley 
1274475e0ac9SMatthew G. Knepley   Note: We are using a first order FEM model for the weak form:
1275475e0ac9SMatthew G. Knepley 
1276475e0ac9SMatthew G. Knepley   \int_\Omega \phi g_0(u, u_t, \nabla u, x, t) \psi + \phi {\vec g}_1(u, u_t, \nabla u, x, t) \nabla \psi + \nabla\phi \cdot {\vec g}_2(u, u_t, \nabla u, x, t) \psi + \nabla\phi \cdot {\overleftrightarrow g}_3(u, u_t, \nabla u, x, t) \cdot \nabla \psi
1277475e0ac9SMatthew G. Knepley 
1278475e0ac9SMatthew G. Knepley The calling sequence for the callbacks g0, g1, g2 and g3 is given by:
1279475e0ac9SMatthew G. Knepley 
1280475e0ac9SMatthew G. Knepley $ g0(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1281475e0ac9SMatthew G. Knepley $    const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1282475e0ac9SMatthew G. Knepley $    const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1283475e0ac9SMatthew G. Knepley $    PetscReal t, const PetscReal x[], PetscScalar g0[])
1284475e0ac9SMatthew G. Knepley 
1285475e0ac9SMatthew G. Knepley + dim - the spatial dimension
1286475e0ac9SMatthew G. Knepley . Nf - the number of fields
1287475e0ac9SMatthew G. Knepley . uOff - the offset into u[] and u_t[] for each field
1288475e0ac9SMatthew G. Knepley . uOff_x - the offset into u_x[] for each field
1289475e0ac9SMatthew G. Knepley . u - each field evaluated at the current point
1290475e0ac9SMatthew G. Knepley . u_t - the time derivative of each field evaluated at the current point
1291475e0ac9SMatthew G. Knepley . u_x - the gradient of each field evaluated at the current point
1292475e0ac9SMatthew G. Knepley . aOff - the offset into a[] and a_t[] for each auxiliary field
1293475e0ac9SMatthew G. Knepley . aOff_x - the offset into a_x[] for each auxiliary field
1294475e0ac9SMatthew G. Knepley . a - each auxiliary field evaluated at the current point
1295475e0ac9SMatthew G. Knepley . a_t - the time derivative of each auxiliary field evaluated at the current point
1296475e0ac9SMatthew G. Knepley . a_x - the gradient of auxiliary each field evaluated at the current point
1297475e0ac9SMatthew G. Knepley . t - current time
1298475e0ac9SMatthew G. Knepley . u_tShift - the multiplier a for dF/dU_t
1299475e0ac9SMatthew G. Knepley . x - coordinates of the current point
1300475e0ac9SMatthew G. Knepley - g0 - output values at the current point
1301475e0ac9SMatthew G. Knepley 
1302475e0ac9SMatthew G. Knepley   Level: intermediate
1303475e0ac9SMatthew G. Knepley 
1304475e0ac9SMatthew G. Knepley .seealso: PetscDSGetJacobianPreconditioner(), PetscDSSetJacobian()
1305475e0ac9SMatthew G. Knepley @*/
1306475e0ac9SMatthew G. Knepley PetscErrorCode PetscDSSetJacobianPreconditioner(PetscDS prob, PetscInt f, PetscInt g,
1307475e0ac9SMatthew G. Knepley                                   void (*g0)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1308475e0ac9SMatthew G. Knepley                                              const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1309475e0ac9SMatthew G. Knepley                                              const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1310475e0ac9SMatthew G. Knepley                                              PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscScalar g0[]),
1311475e0ac9SMatthew G. Knepley                                   void (*g1)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1312475e0ac9SMatthew G. Knepley                                              const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1313475e0ac9SMatthew G. Knepley                                              const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1314475e0ac9SMatthew G. Knepley                                              PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscScalar g1[]),
1315475e0ac9SMatthew G. Knepley                                   void (*g2)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1316475e0ac9SMatthew G. Knepley                                              const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1317475e0ac9SMatthew G. Knepley                                              const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1318475e0ac9SMatthew G. Knepley                                              PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscScalar g2[]),
1319475e0ac9SMatthew G. Knepley                                   void (*g3)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1320475e0ac9SMatthew G. Knepley                                              const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1321475e0ac9SMatthew G. Knepley                                              const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1322475e0ac9SMatthew G. Knepley                                              PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscScalar g3[]))
1323475e0ac9SMatthew G. Knepley {
1324475e0ac9SMatthew G. Knepley   PetscErrorCode ierr;
1325475e0ac9SMatthew G. Knepley 
1326475e0ac9SMatthew G. Knepley   PetscFunctionBegin;
1327475e0ac9SMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1328475e0ac9SMatthew G. Knepley   if (g0) PetscValidFunction(g0, 4);
1329475e0ac9SMatthew G. Knepley   if (g1) PetscValidFunction(g1, 5);
1330475e0ac9SMatthew G. Knepley   if (g2) PetscValidFunction(g2, 6);
1331475e0ac9SMatthew G. Knepley   if (g3) PetscValidFunction(g3, 7);
1332475e0ac9SMatthew G. Knepley   if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f);
1333475e0ac9SMatthew G. Knepley   if (g < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", g);
1334475e0ac9SMatthew G. Knepley   ierr = PetscDSEnlarge_Static(prob, PetscMax(f, g)+1);CHKERRQ(ierr);
1335475e0ac9SMatthew G. Knepley   prob->gp[(f*prob->Nf + g)*4+0] = g0;
1336475e0ac9SMatthew G. Knepley   prob->gp[(f*prob->Nf + g)*4+1] = g1;
1337475e0ac9SMatthew G. Knepley   prob->gp[(f*prob->Nf + g)*4+2] = g2;
1338475e0ac9SMatthew G. Knepley   prob->gp[(f*prob->Nf + g)*4+3] = g3;
1339475e0ac9SMatthew G. Knepley   PetscFunctionReturn(0);
1340475e0ac9SMatthew G. Knepley }
1341475e0ac9SMatthew G. Knepley 
1342b7e05686SMatthew G. Knepley /*@C
1343b7e05686SMatthew G. Knepley   PetscDSHasDynamicJacobian - Signals that a dynamic Jacobian, dF/du_t, has been set
1344b7e05686SMatthew G. Knepley 
1345b7e05686SMatthew G. Knepley   Not collective
1346b7e05686SMatthew G. Knepley 
1347b7e05686SMatthew G. Knepley   Input Parameter:
1348b7e05686SMatthew G. Knepley . prob - The PetscDS
1349b7e05686SMatthew G. Knepley 
1350b7e05686SMatthew G. Knepley   Output Parameter:
1351b7e05686SMatthew G. Knepley . hasDynJac - flag that pointwise function for dynamic Jacobian has been set
1352b7e05686SMatthew G. Knepley 
1353b7e05686SMatthew G. Knepley   Level: intermediate
1354b7e05686SMatthew G. Knepley 
1355b7e05686SMatthew G. Knepley .seealso: PetscDSGetDynamicJacobian(), PetscDSSetDynamicJacobian(), PetscDSGetJacobian()
1356b7e05686SMatthew G. Knepley @*/
1357b7e05686SMatthew G. Knepley PetscErrorCode PetscDSHasDynamicJacobian(PetscDS prob, PetscBool *hasDynJac)
1358b7e05686SMatthew G. Knepley {
1359b7e05686SMatthew G. Knepley   PetscInt f, g, h;
1360b7e05686SMatthew G. Knepley 
1361b7e05686SMatthew G. Knepley   PetscFunctionBegin;
1362b7e05686SMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1363b7e05686SMatthew G. Knepley   *hasDynJac = PETSC_FALSE;
1364b7e05686SMatthew G. Knepley   for (f = 0; f < prob->Nf; ++f) {
1365b7e05686SMatthew G. Knepley     for (g = 0; g < prob->Nf; ++g) {
1366b7e05686SMatthew G. Knepley       for (h = 0; h < 4; ++h) {
1367b7e05686SMatthew G. Knepley         if (prob->gt[(f*prob->Nf + g)*4+h]) *hasDynJac = PETSC_TRUE;
1368b7e05686SMatthew G. Knepley       }
1369b7e05686SMatthew G. Knepley     }
1370b7e05686SMatthew G. Knepley   }
1371b7e05686SMatthew G. Knepley   PetscFunctionReturn(0);
1372b7e05686SMatthew G. Knepley }
1373b7e05686SMatthew G. Knepley 
1374b7e05686SMatthew G. Knepley /*@C
1375b7e05686SMatthew G. Knepley   PetscDSGetDynamicJacobian - Get the pointwise dynamic Jacobian, dF/du_t, function for given test and basis field
1376b7e05686SMatthew G. Knepley 
1377b7e05686SMatthew G. Knepley   Not collective
1378b7e05686SMatthew G. Knepley 
1379b7e05686SMatthew G. Knepley   Input Parameters:
1380b7e05686SMatthew G. Knepley + prob - The PetscDS
1381b7e05686SMatthew G. Knepley . f    - The test field number
1382b7e05686SMatthew G. Knepley - g    - The field number
1383b7e05686SMatthew G. Knepley 
1384b7e05686SMatthew G. Knepley   Output Parameters:
1385b7e05686SMatthew G. Knepley + g0 - integrand for the test and basis function term
1386b7e05686SMatthew G. Knepley . g1 - integrand for the test function and basis function gradient term
1387b7e05686SMatthew G. Knepley . g2 - integrand for the test function gradient and basis function term
1388b7e05686SMatthew G. Knepley - g3 - integrand for the test function gradient and basis function gradient term
1389b7e05686SMatthew G. Knepley 
1390b7e05686SMatthew G. Knepley   Note: We are using a first order FEM model for the weak form:
1391b7e05686SMatthew G. Knepley 
1392b7e05686SMatthew G. Knepley   \int_\Omega \phi g_0(u, u_t, \nabla u, x, t) \psi + \phi {\vec g}_1(u, u_t, \nabla u, x, t) \nabla \psi + \nabla\phi \cdot {\vec g}_2(u, u_t, \nabla u, x, t) \psi + \nabla\phi \cdot {\overleftrightarrow g}_3(u, u_t, \nabla u, x, t) \cdot \nabla \psi
1393b7e05686SMatthew G. Knepley 
1394b7e05686SMatthew G. Knepley The calling sequence for the callbacks g0, g1, g2 and g3 is given by:
1395b7e05686SMatthew G. Knepley 
1396b7e05686SMatthew G. Knepley $ g0(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1397b7e05686SMatthew G. Knepley $    const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1398b7e05686SMatthew G. Knepley $    const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1399b7e05686SMatthew G. Knepley $    PetscReal t, const PetscReal u_tShift, const PetscReal x[], PetscScalar g0[])
1400b7e05686SMatthew G. Knepley 
1401b7e05686SMatthew G. Knepley + dim - the spatial dimension
1402b7e05686SMatthew G. Knepley . Nf - the number of fields
1403b7e05686SMatthew G. Knepley . uOff - the offset into u[] and u_t[] for each field
1404b7e05686SMatthew G. Knepley . uOff_x - the offset into u_x[] for each field
1405b7e05686SMatthew G. Knepley . u - each field evaluated at the current point
1406b7e05686SMatthew G. Knepley . u_t - the time derivative of each field evaluated at the current point
1407b7e05686SMatthew G. Knepley . u_x - the gradient of each field evaluated at the current point
1408b7e05686SMatthew G. Knepley . aOff - the offset into a[] and a_t[] for each auxiliary field
1409b7e05686SMatthew G. Knepley . aOff_x - the offset into a_x[] for each auxiliary field
1410b7e05686SMatthew G. Knepley . a - each auxiliary field evaluated at the current point
1411b7e05686SMatthew G. Knepley . a_t - the time derivative of each auxiliary field evaluated at the current point
1412b7e05686SMatthew G. Knepley . a_x - the gradient of auxiliary each field evaluated at the current point
1413b7e05686SMatthew G. Knepley . t - current time
1414b7e05686SMatthew G. Knepley . u_tShift - the multiplier a for dF/dU_t
1415b7e05686SMatthew G. Knepley . x - coordinates of the current point
1416b7e05686SMatthew G. Knepley - g0 - output values at the current point
1417b7e05686SMatthew G. Knepley 
1418b7e05686SMatthew G. Knepley   Level: intermediate
1419b7e05686SMatthew G. Knepley 
1420b7e05686SMatthew G. Knepley .seealso: PetscDSSetJacobian()
1421b7e05686SMatthew G. Knepley @*/
1422b7e05686SMatthew G. Knepley PetscErrorCode PetscDSGetDynamicJacobian(PetscDS prob, PetscInt f, PetscInt g,
1423b7e05686SMatthew G. Knepley                                          void (**g0)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1424b7e05686SMatthew G. Knepley                                                      const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1425b7e05686SMatthew G. Knepley                                                      const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1426b7e05686SMatthew G. Knepley                                                      PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscScalar g0[]),
1427b7e05686SMatthew G. Knepley                                          void (**g1)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1428b7e05686SMatthew G. Knepley                                                      const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1429b7e05686SMatthew G. Knepley                                                      const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1430b7e05686SMatthew G. Knepley                                                      PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscScalar g1[]),
1431b7e05686SMatthew G. Knepley                                          void (**g2)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1432b7e05686SMatthew G. Knepley                                                      const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1433b7e05686SMatthew G. Knepley                                                      const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1434b7e05686SMatthew G. Knepley                                               PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscScalar g2[]),
1435b7e05686SMatthew G. Knepley                                          void (**g3)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1436b7e05686SMatthew G. Knepley                                                      const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1437b7e05686SMatthew G. Knepley                                                      const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1438b7e05686SMatthew G. Knepley                                                      PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscScalar g3[]))
1439b7e05686SMatthew G. Knepley {
1440b7e05686SMatthew G. Knepley   PetscFunctionBegin;
1441b7e05686SMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1442b7e05686SMatthew 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);
1443b7e05686SMatthew 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);
1444b7e05686SMatthew G. Knepley   if (g0) {PetscValidPointer(g0, 4); *g0 = prob->gt[(f*prob->Nf + g)*4+0];}
1445b7e05686SMatthew G. Knepley   if (g1) {PetscValidPointer(g1, 5); *g1 = prob->gt[(f*prob->Nf + g)*4+1];}
1446b7e05686SMatthew G. Knepley   if (g2) {PetscValidPointer(g2, 6); *g2 = prob->gt[(f*prob->Nf + g)*4+2];}
1447b7e05686SMatthew G. Knepley   if (g3) {PetscValidPointer(g3, 7); *g3 = prob->gt[(f*prob->Nf + g)*4+3];}
1448b7e05686SMatthew G. Knepley   PetscFunctionReturn(0);
1449b7e05686SMatthew G. Knepley }
1450b7e05686SMatthew G. Knepley 
1451b7e05686SMatthew G. Knepley /*@C
1452b7e05686SMatthew G. Knepley   PetscDSSetDynamicJacobian - Set the pointwise dynamic Jacobian, dF/du_t, function for given test and basis fields
1453b7e05686SMatthew G. Knepley 
1454b7e05686SMatthew G. Knepley   Not collective
1455b7e05686SMatthew G. Knepley 
1456b7e05686SMatthew G. Knepley   Input Parameters:
1457b7e05686SMatthew G. Knepley + prob - The PetscDS
1458b7e05686SMatthew G. Knepley . f    - The test field number
1459b7e05686SMatthew G. Knepley . g    - The field number
1460b7e05686SMatthew G. Knepley . g0 - integrand for the test and basis function term
1461b7e05686SMatthew G. Knepley . g1 - integrand for the test function and basis function gradient term
1462b7e05686SMatthew G. Knepley . g2 - integrand for the test function gradient and basis function term
1463b7e05686SMatthew G. Knepley - g3 - integrand for the test function gradient and basis function gradient term
1464b7e05686SMatthew G. Knepley 
1465b7e05686SMatthew G. Knepley   Note: We are using a first order FEM model for the weak form:
1466b7e05686SMatthew G. Knepley 
1467b7e05686SMatthew G. Knepley   \int_\Omega \phi g_0(u, u_t, \nabla u, x, t) \psi + \phi {\vec g}_1(u, u_t, \nabla u, x, t) \nabla \psi + \nabla\phi \cdot {\vec g}_2(u, u_t, \nabla u, x, t) \psi + \nabla\phi \cdot {\overleftrightarrow g}_3(u, u_t, \nabla u, x, t) \cdot \nabla \psi
1468b7e05686SMatthew G. Knepley 
1469b7e05686SMatthew G. Knepley The calling sequence for the callbacks g0, g1, g2 and g3 is given by:
1470b7e05686SMatthew G. Knepley 
1471b7e05686SMatthew G. Knepley $ g0(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1472b7e05686SMatthew G. Knepley $    const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1473b7e05686SMatthew G. Knepley $    const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1474b7e05686SMatthew G. Knepley $    PetscReal t, const PetscReal x[], PetscScalar g0[])
1475b7e05686SMatthew G. Knepley 
1476b7e05686SMatthew G. Knepley + dim - the spatial dimension
1477b7e05686SMatthew G. Knepley . Nf - the number of fields
1478b7e05686SMatthew G. Knepley . uOff - the offset into u[] and u_t[] for each field
1479b7e05686SMatthew G. Knepley . uOff_x - the offset into u_x[] for each field
1480b7e05686SMatthew G. Knepley . u - each field evaluated at the current point
1481b7e05686SMatthew G. Knepley . u_t - the time derivative of each field evaluated at the current point
1482b7e05686SMatthew G. Knepley . u_x - the gradient of each field evaluated at the current point
1483b7e05686SMatthew G. Knepley . aOff - the offset into a[] and a_t[] for each auxiliary field
1484b7e05686SMatthew G. Knepley . aOff_x - the offset into a_x[] for each auxiliary field
1485b7e05686SMatthew G. Knepley . a - each auxiliary field evaluated at the current point
1486b7e05686SMatthew G. Knepley . a_t - the time derivative of each auxiliary field evaluated at the current point
1487b7e05686SMatthew G. Knepley . a_x - the gradient of auxiliary each field evaluated at the current point
1488b7e05686SMatthew G. Knepley . t - current time
1489b7e05686SMatthew G. Knepley . u_tShift - the multiplier a for dF/dU_t
1490b7e05686SMatthew G. Knepley . x - coordinates of the current point
1491b7e05686SMatthew G. Knepley - g0 - output values at the current point
1492b7e05686SMatthew G. Knepley 
1493b7e05686SMatthew G. Knepley   Level: intermediate
1494b7e05686SMatthew G. Knepley 
1495b7e05686SMatthew G. Knepley .seealso: PetscDSGetJacobian()
1496b7e05686SMatthew G. Knepley @*/
1497b7e05686SMatthew G. Knepley PetscErrorCode PetscDSSetDynamicJacobian(PetscDS prob, PetscInt f, PetscInt g,
1498b7e05686SMatthew G. Knepley                                          void (*g0)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1499b7e05686SMatthew G. Knepley                                                     const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1500b7e05686SMatthew G. Knepley                                                     const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1501b7e05686SMatthew G. Knepley                                                     PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscScalar g0[]),
1502b7e05686SMatthew G. Knepley                                          void (*g1)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1503b7e05686SMatthew G. Knepley                                                     const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1504b7e05686SMatthew G. Knepley                                                     const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1505b7e05686SMatthew G. Knepley                                                     PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscScalar g1[]),
1506b7e05686SMatthew G. Knepley                                          void (*g2)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1507b7e05686SMatthew G. Knepley                                                     const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1508b7e05686SMatthew G. Knepley                                                     const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1509b7e05686SMatthew G. Knepley                                                     PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscScalar g2[]),
1510b7e05686SMatthew G. Knepley                                          void (*g3)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1511b7e05686SMatthew G. Knepley                                                     const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1512b7e05686SMatthew G. Knepley                                                     const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1513b7e05686SMatthew G. Knepley                                                     PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscScalar g3[]))
1514b7e05686SMatthew G. Knepley {
1515b7e05686SMatthew G. Knepley   PetscErrorCode ierr;
1516b7e05686SMatthew G. Knepley 
1517b7e05686SMatthew G. Knepley   PetscFunctionBegin;
1518b7e05686SMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1519b7e05686SMatthew G. Knepley   if (g0) PetscValidFunction(g0, 4);
1520b7e05686SMatthew G. Knepley   if (g1) PetscValidFunction(g1, 5);
1521b7e05686SMatthew G. Knepley   if (g2) PetscValidFunction(g2, 6);
1522b7e05686SMatthew G. Knepley   if (g3) PetscValidFunction(g3, 7);
1523b7e05686SMatthew G. Knepley   if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f);
1524b7e05686SMatthew G. Knepley   if (g < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", g);
1525b7e05686SMatthew G. Knepley   ierr = PetscDSEnlarge_Static(prob, PetscMax(f, g)+1);CHKERRQ(ierr);
1526b7e05686SMatthew G. Knepley   prob->gt[(f*prob->Nf + g)*4+0] = g0;
1527b7e05686SMatthew G. Knepley   prob->gt[(f*prob->Nf + g)*4+1] = g1;
1528b7e05686SMatthew G. Knepley   prob->gt[(f*prob->Nf + g)*4+2] = g2;
1529b7e05686SMatthew G. Knepley   prob->gt[(f*prob->Nf + g)*4+3] = g3;
1530b7e05686SMatthew G. Knepley   PetscFunctionReturn(0);
1531b7e05686SMatthew G. Knepley }
1532b7e05686SMatthew G. Knepley 
15330c2f2876SMatthew G. Knepley /*@C
15340c2f2876SMatthew G. Knepley   PetscDSGetRiemannSolver - Returns the Riemann solver for the given field
15350c2f2876SMatthew G. Knepley 
15360c2f2876SMatthew G. Knepley   Not collective
15370c2f2876SMatthew G. Knepley 
15380c2f2876SMatthew G. Knepley   Input Arguments:
15390c2f2876SMatthew G. Knepley + prob - The PetscDS object
15400c2f2876SMatthew G. Knepley - f    - The field number
15410c2f2876SMatthew G. Knepley 
15420c2f2876SMatthew G. Knepley   Output Argument:
15430c2f2876SMatthew G. Knepley . r    - Riemann solver
15440c2f2876SMatthew G. Knepley 
15450c2f2876SMatthew G. Knepley   Calling sequence for r:
15460c2f2876SMatthew G. Knepley 
15475db36cf9SMatthew G. Knepley $ r(PetscInt dim, PetscInt Nf, const PetscReal x[], const PetscReal n[], const PetscScalar uL[], const PetscScalar uR[], PetscScalar flux[], void *ctx)
15480c2f2876SMatthew G. Knepley 
15495db36cf9SMatthew G. Knepley + dim  - The spatial dimension
15505db36cf9SMatthew G. Knepley . Nf   - The number of fields
15515db36cf9SMatthew G. Knepley . x    - The coordinates at a point on the interface
15520c2f2876SMatthew G. Knepley . n    - The normal vector to the interface
15530c2f2876SMatthew G. Knepley . uL   - The state vector to the left of the interface
15540c2f2876SMatthew G. Knepley . uR   - The state vector to the right of the interface
15550c2f2876SMatthew G. Knepley . flux - output array of flux through the interface
15560c2f2876SMatthew G. Knepley - ctx  - optional user context
15570c2f2876SMatthew G. Knepley 
15580c2f2876SMatthew G. Knepley   Level: intermediate
15590c2f2876SMatthew G. Knepley 
15600c2f2876SMatthew G. Knepley .seealso: PetscDSSetRiemannSolver()
15610c2f2876SMatthew G. Knepley @*/
15620c2f2876SMatthew G. Knepley PetscErrorCode PetscDSGetRiemannSolver(PetscDS prob, PetscInt f,
15635db36cf9SMatthew G. Knepley                                        void (**r)(PetscInt dim, PetscInt Nf, const PetscReal x[], const PetscReal n[], const PetscScalar uL[], const PetscScalar uR[], PetscScalar flux[], void *ctx))
15640c2f2876SMatthew G. Knepley {
15650c2f2876SMatthew G. Knepley   PetscFunctionBegin;
15660c2f2876SMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
15670c2f2876SMatthew 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);
15680c2f2876SMatthew G. Knepley   PetscValidPointer(r, 3);
15690c2f2876SMatthew G. Knepley   *r = prob->r[f];
15700c2f2876SMatthew G. Knepley   PetscFunctionReturn(0);
15710c2f2876SMatthew G. Knepley }
15720c2f2876SMatthew G. Knepley 
15730c2f2876SMatthew G. Knepley /*@C
15740c2f2876SMatthew G. Knepley   PetscDSSetRiemannSolver - Sets the Riemann solver for the given field
15750c2f2876SMatthew G. Knepley 
15760c2f2876SMatthew G. Knepley   Not collective
15770c2f2876SMatthew G. Knepley 
15780c2f2876SMatthew G. Knepley   Input Arguments:
15790c2f2876SMatthew G. Knepley + prob - The PetscDS object
15800c2f2876SMatthew G. Knepley . f    - The field number
15810c2f2876SMatthew G. Knepley - r    - Riemann solver
15820c2f2876SMatthew G. Knepley 
15830c2f2876SMatthew G. Knepley   Calling sequence for r:
15840c2f2876SMatthew G. Knepley 
15855db36cf9SMatthew G. Knepley $ r(PetscInt dim, PetscInt Nf, const PetscReal x[], const PetscReal n[], const PetscScalar uL[], const PetscScalar uR[], PetscScalar flux[], void *ctx)
15860c2f2876SMatthew G. Knepley 
15875db36cf9SMatthew G. Knepley + dim  - The spatial dimension
15885db36cf9SMatthew G. Knepley . Nf   - The number of fields
15895db36cf9SMatthew G. Knepley . x    - The coordinates at a point on the interface
15900c2f2876SMatthew G. Knepley . n    - The normal vector to the interface
15910c2f2876SMatthew G. Knepley . uL   - The state vector to the left of the interface
15920c2f2876SMatthew G. Knepley . uR   - The state vector to the right of the interface
15930c2f2876SMatthew G. Knepley . flux - output array of flux through the interface
15940c2f2876SMatthew G. Knepley - ctx  - optional user context
15950c2f2876SMatthew G. Knepley 
15960c2f2876SMatthew G. Knepley   Level: intermediate
15970c2f2876SMatthew G. Knepley 
15980c2f2876SMatthew G. Knepley .seealso: PetscDSGetRiemannSolver()
15990c2f2876SMatthew G. Knepley @*/
16000c2f2876SMatthew G. Knepley PetscErrorCode PetscDSSetRiemannSolver(PetscDS prob, PetscInt f,
16015db36cf9SMatthew G. Knepley                                        void (*r)(PetscInt dim, PetscInt Nf, const PetscReal x[], const PetscReal n[], const PetscScalar uL[], const PetscScalar uR[], PetscScalar flux[], void *ctx))
16020c2f2876SMatthew G. Knepley {
16030c2f2876SMatthew G. Knepley   PetscErrorCode ierr;
16040c2f2876SMatthew G. Knepley 
16050c2f2876SMatthew G. Knepley   PetscFunctionBegin;
16060c2f2876SMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1607de716cbcSToby Isaac   if (r) PetscValidFunction(r, 3);
16080c2f2876SMatthew G. Knepley   if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f);
16090c2f2876SMatthew G. Knepley   ierr = PetscDSEnlarge_Static(prob, f+1);CHKERRQ(ierr);
16100c2f2876SMatthew G. Knepley   prob->r[f] = r;
16110c2f2876SMatthew G. Knepley   PetscFunctionReturn(0);
16120c2f2876SMatthew G. Knepley }
16130c2f2876SMatthew G. Knepley 
16140c2f2876SMatthew G. Knepley PetscErrorCode PetscDSGetContext(PetscDS prob, PetscInt f, void **ctx)
16150c2f2876SMatthew G. Knepley {
16160c2f2876SMatthew G. Knepley   PetscFunctionBegin;
16170c2f2876SMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
16180c2f2876SMatthew 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);
16190c2f2876SMatthew G. Knepley   PetscValidPointer(ctx, 3);
16200c2f2876SMatthew G. Knepley   *ctx = prob->ctx[f];
16210c2f2876SMatthew G. Knepley   PetscFunctionReturn(0);
16220c2f2876SMatthew G. Knepley }
16230c2f2876SMatthew G. Knepley 
16240c2f2876SMatthew G. Knepley PetscErrorCode PetscDSSetContext(PetscDS prob, PetscInt f, void *ctx)
16250c2f2876SMatthew G. Knepley {
16260c2f2876SMatthew G. Knepley   PetscErrorCode ierr;
16270c2f2876SMatthew G. Knepley 
16280c2f2876SMatthew G. Knepley   PetscFunctionBegin;
16290c2f2876SMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
16300c2f2876SMatthew G. Knepley   if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f);
16310c2f2876SMatthew G. Knepley   ierr = PetscDSEnlarge_Static(prob, f+1);CHKERRQ(ierr);
16320c2f2876SMatthew G. Knepley   prob->ctx[f] = ctx;
16330c2f2876SMatthew G. Knepley   PetscFunctionReturn(0);
16340c2f2876SMatthew G. Knepley }
16350c2f2876SMatthew G. Knepley 
1636194d53e6SMatthew G. Knepley /*@C
1637194d53e6SMatthew G. Knepley   PetscDSGetBdResidual - Get the pointwise boundary residual function for a given test field
1638194d53e6SMatthew G. Knepley 
1639194d53e6SMatthew G. Knepley   Not collective
1640194d53e6SMatthew G. Knepley 
1641194d53e6SMatthew G. Knepley   Input Parameters:
1642194d53e6SMatthew G. Knepley + prob - The PetscDS
1643194d53e6SMatthew G. Knepley - f    - The test field number
1644194d53e6SMatthew G. Knepley 
1645194d53e6SMatthew G. Knepley   Output Parameters:
1646194d53e6SMatthew G. Knepley + f0 - boundary integrand for the test function term
1647194d53e6SMatthew G. Knepley - f1 - boundary integrand for the test function gradient term
1648194d53e6SMatthew G. Knepley 
1649194d53e6SMatthew G. Knepley   Note: We are using a first order FEM model for the weak form:
1650194d53e6SMatthew G. Knepley 
1651194d53e6SMatthew G. Knepley   \int_\Gamma \phi {\vec f}_0(u, u_t, \nabla u, x, t) \cdot \hat n + \nabla\phi \cdot {\overleftrightarrow f}_1(u, u_t, \nabla u, x, t) \cdot \hat n
1652194d53e6SMatthew G. Knepley 
1653194d53e6SMatthew G. Knepley The calling sequence for the callbacks f0 and f1 is given by:
1654194d53e6SMatthew G. Knepley 
165530b9ff8bSMatthew G. Knepley $ f0(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1656194d53e6SMatthew G. Knepley $    const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1657194d53e6SMatthew G. Knepley $    const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
165830b9ff8bSMatthew G. Knepley $    PetscReal t, const PetscReal x[], const PetscReal n[], PetscScalar f0[])
1659194d53e6SMatthew G. Knepley 
1660194d53e6SMatthew G. Knepley + dim - the spatial dimension
1661194d53e6SMatthew G. Knepley . Nf - the number of fields
1662194d53e6SMatthew G. Knepley . uOff - the offset into u[] and u_t[] for each field
1663194d53e6SMatthew G. Knepley . uOff_x - the offset into u_x[] for each field
1664194d53e6SMatthew G. Knepley . u - each field evaluated at the current point
1665194d53e6SMatthew G. Knepley . u_t - the time derivative of each field evaluated at the current point
1666194d53e6SMatthew G. Knepley . u_x - the gradient of each field evaluated at the current point
1667194d53e6SMatthew G. Knepley . aOff - the offset into a[] and a_t[] for each auxiliary field
1668194d53e6SMatthew G. Knepley . aOff_x - the offset into a_x[] for each auxiliary field
1669194d53e6SMatthew G. Knepley . a - each auxiliary field evaluated at the current point
1670194d53e6SMatthew G. Knepley . a_t - the time derivative of each auxiliary field evaluated at the current point
1671194d53e6SMatthew G. Knepley . a_x - the gradient of auxiliary each field evaluated at the current point
1672194d53e6SMatthew G. Knepley . t - current time
1673194d53e6SMatthew G. Knepley . x - coordinates of the current point
1674194d53e6SMatthew G. Knepley . n - unit normal at the current point
1675194d53e6SMatthew G. Knepley - f0 - output values at the current point
1676194d53e6SMatthew G. Knepley 
1677194d53e6SMatthew G. Knepley   Level: intermediate
1678194d53e6SMatthew G. Knepley 
1679194d53e6SMatthew G. Knepley .seealso: PetscDSSetBdResidual()
1680194d53e6SMatthew G. Knepley @*/
16812764a2aaSMatthew G. Knepley PetscErrorCode PetscDSGetBdResidual(PetscDS prob, PetscInt f,
168230b9ff8bSMatthew G. Knepley                                     void (**f0)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1683194d53e6SMatthew G. Knepley                                                 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1684194d53e6SMatthew G. Knepley                                                 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
168530b9ff8bSMatthew G. Knepley                                                 PetscReal t, const PetscReal x[], const PetscReal n[], PetscScalar f0[]),
168630b9ff8bSMatthew G. Knepley                                     void (**f1)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1687194d53e6SMatthew G. Knepley                                                 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1688194d53e6SMatthew G. Knepley                                                 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
168930b9ff8bSMatthew G. Knepley                                                 PetscReal t, const PetscReal x[], const PetscReal n[], PetscScalar f1[]))
16902764a2aaSMatthew G. Knepley {
16912764a2aaSMatthew G. Knepley   PetscFunctionBegin;
16922764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
16932764a2aaSMatthew 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);
16942764a2aaSMatthew G. Knepley   if (f0) {PetscValidPointer(f0, 3); *f0 = prob->fBd[f*2+0];}
16952764a2aaSMatthew G. Knepley   if (f1) {PetscValidPointer(f1, 4); *f1 = prob->fBd[f*2+1];}
16962764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
16972764a2aaSMatthew G. Knepley }
16982764a2aaSMatthew G. Knepley 
1699194d53e6SMatthew G. Knepley /*@C
1700194d53e6SMatthew G. Knepley   PetscDSSetBdResidual - Get the pointwise boundary residual function for a given test field
1701194d53e6SMatthew G. Knepley 
1702194d53e6SMatthew G. Knepley   Not collective
1703194d53e6SMatthew G. Knepley 
1704194d53e6SMatthew G. Knepley   Input Parameters:
1705194d53e6SMatthew G. Knepley + prob - The PetscDS
1706194d53e6SMatthew G. Knepley . f    - The test field number
1707194d53e6SMatthew G. Knepley . f0 - boundary integrand for the test function term
1708194d53e6SMatthew G. Knepley - f1 - boundary integrand for the test function gradient term
1709194d53e6SMatthew G. Knepley 
1710194d53e6SMatthew G. Knepley   Note: We are using a first order FEM model for the weak form:
1711194d53e6SMatthew G. Knepley 
1712194d53e6SMatthew G. Knepley   \int_\Gamma \phi {\vec f}_0(u, u_t, \nabla u, x, t) \cdot \hat n + \nabla\phi \cdot {\overleftrightarrow f}_1(u, u_t, \nabla u, x, t) \cdot \hat n
1713194d53e6SMatthew G. Knepley 
1714194d53e6SMatthew G. Knepley The calling sequence for the callbacks f0 and f1 is given by:
1715194d53e6SMatthew G. Knepley 
171630b9ff8bSMatthew G. Knepley $ f0(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1717194d53e6SMatthew G. Knepley $    const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1718194d53e6SMatthew G. Knepley $    const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
171930b9ff8bSMatthew G. Knepley $    PetscReal t, const PetscReal x[], const PetscReal n[], PetscScalar f0[])
1720194d53e6SMatthew G. Knepley 
1721194d53e6SMatthew G. Knepley + dim - the spatial dimension
1722194d53e6SMatthew G. Knepley . Nf - the number of fields
1723194d53e6SMatthew G. Knepley . uOff - the offset into u[] and u_t[] for each field
1724194d53e6SMatthew G. Knepley . uOff_x - the offset into u_x[] for each field
1725194d53e6SMatthew G. Knepley . u - each field evaluated at the current point
1726194d53e6SMatthew G. Knepley . u_t - the time derivative of each field evaluated at the current point
1727194d53e6SMatthew G. Knepley . u_x - the gradient of each field evaluated at the current point
1728194d53e6SMatthew G. Knepley . aOff - the offset into a[] and a_t[] for each auxiliary field
1729194d53e6SMatthew G. Knepley . aOff_x - the offset into a_x[] for each auxiliary field
1730194d53e6SMatthew G. Knepley . a - each auxiliary field evaluated at the current point
1731194d53e6SMatthew G. Knepley . a_t - the time derivative of each auxiliary field evaluated at the current point
1732194d53e6SMatthew G. Knepley . a_x - the gradient of auxiliary each field evaluated at the current point
1733194d53e6SMatthew G. Knepley . t - current time
1734194d53e6SMatthew G. Knepley . x - coordinates of the current point
1735194d53e6SMatthew G. Knepley . n - unit normal at the current point
1736194d53e6SMatthew G. Knepley - f0 - output values at the current point
1737194d53e6SMatthew G. Knepley 
1738194d53e6SMatthew G. Knepley   Level: intermediate
1739194d53e6SMatthew G. Knepley 
1740194d53e6SMatthew G. Knepley .seealso: PetscDSGetBdResidual()
1741194d53e6SMatthew G. Knepley @*/
17422764a2aaSMatthew G. Knepley PetscErrorCode PetscDSSetBdResidual(PetscDS prob, PetscInt f,
174330b9ff8bSMatthew G. Knepley                                     void (*f0)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1744194d53e6SMatthew G. Knepley                                                const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1745194d53e6SMatthew G. Knepley                                                const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
174630b9ff8bSMatthew G. Knepley                                                PetscReal t, const PetscReal x[], const PetscReal n[], PetscScalar f0[]),
174730b9ff8bSMatthew G. Knepley                                     void (*f1)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1748194d53e6SMatthew G. Knepley                                                const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1749194d53e6SMatthew G. Knepley                                                const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
175030b9ff8bSMatthew G. Knepley                                                PetscReal t, const PetscReal x[], const PetscReal n[], PetscScalar f1[]))
17512764a2aaSMatthew G. Knepley {
17522764a2aaSMatthew G. Knepley   PetscErrorCode ierr;
17532764a2aaSMatthew G. Knepley 
17542764a2aaSMatthew G. Knepley   PetscFunctionBegin;
17552764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
17562764a2aaSMatthew G. Knepley   if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f);
17572764a2aaSMatthew G. Knepley   ierr = PetscDSEnlarge_Static(prob, f+1);CHKERRQ(ierr);
17582764a2aaSMatthew G. Knepley   if (f0) {PetscValidFunction(f0, 3); prob->fBd[f*2+0] = f0;}
17592764a2aaSMatthew G. Knepley   if (f1) {PetscValidFunction(f1, 4); prob->fBd[f*2+1] = f1;}
17602764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
17612764a2aaSMatthew G. Knepley }
17622764a2aaSMatthew G. Knepley 
1763194d53e6SMatthew G. Knepley /*@C
1764194d53e6SMatthew G. Knepley   PetscDSGetBdJacobian - Get the pointwise boundary Jacobian function for given test and basis field
1765194d53e6SMatthew G. Knepley 
1766194d53e6SMatthew G. Knepley   Not collective
1767194d53e6SMatthew G. Knepley 
1768194d53e6SMatthew G. Knepley   Input Parameters:
1769194d53e6SMatthew G. Knepley + prob - The PetscDS
1770194d53e6SMatthew G. Knepley . f    - The test field number
1771194d53e6SMatthew G. Knepley - g    - The field number
1772194d53e6SMatthew G. Knepley 
1773194d53e6SMatthew G. Knepley   Output Parameters:
1774194d53e6SMatthew G. Knepley + g0 - integrand for the test and basis function term
1775194d53e6SMatthew G. Knepley . g1 - integrand for the test function and basis function gradient term
1776194d53e6SMatthew G. Knepley . g2 - integrand for the test function gradient and basis function term
1777194d53e6SMatthew G. Knepley - g3 - integrand for the test function gradient and basis function gradient term
1778194d53e6SMatthew G. Knepley 
1779194d53e6SMatthew G. Knepley   Note: We are using a first order FEM model for the weak form:
1780194d53e6SMatthew G. Knepley 
1781194d53e6SMatthew G. Knepley   \int_\Gamma \phi {\vec g}_0(u, u_t, \nabla u, x, t) \cdot \hat n \psi + \phi {\vec g}_1(u, u_t, \nabla u, x, t) \cdot \hat n \nabla \psi + \nabla\phi \cdot {\vec g}_2(u, u_t, \nabla u, x, t) \cdot \hat n \psi + \nabla\phi \cdot {\overleftrightarrow g}_3(u, u_t, \nabla u, x, t) \cdot \hat n \cdot \nabla \psi
1782194d53e6SMatthew G. Knepley 
1783194d53e6SMatthew G. Knepley The calling sequence for the callbacks g0, g1, g2 and g3 is given by:
1784194d53e6SMatthew G. Knepley 
178530b9ff8bSMatthew G. Knepley $ g0(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1786194d53e6SMatthew G. Knepley $    const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1787194d53e6SMatthew G. Knepley $    const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
178830b9ff8bSMatthew G. Knepley $    PetscReal t, const PetscReal x[], const PetscReal n[], PetscScalar g0[])
1789194d53e6SMatthew G. Knepley 
1790194d53e6SMatthew G. Knepley + dim - the spatial dimension
1791194d53e6SMatthew G. Knepley . Nf - the number of fields
1792194d53e6SMatthew G. Knepley . uOff - the offset into u[] and u_t[] for each field
1793194d53e6SMatthew G. Knepley . uOff_x - the offset into u_x[] for each field
1794194d53e6SMatthew G. Knepley . u - each field evaluated at the current point
1795194d53e6SMatthew G. Knepley . u_t - the time derivative of each field evaluated at the current point
1796194d53e6SMatthew G. Knepley . u_x - the gradient of each field evaluated at the current point
1797194d53e6SMatthew G. Knepley . aOff - the offset into a[] and a_t[] for each auxiliary field
1798194d53e6SMatthew G. Knepley . aOff_x - the offset into a_x[] for each auxiliary field
1799194d53e6SMatthew G. Knepley . a - each auxiliary field evaluated at the current point
1800194d53e6SMatthew G. Knepley . a_t - the time derivative of each auxiliary field evaluated at the current point
1801194d53e6SMatthew G. Knepley . a_x - the gradient of auxiliary each field evaluated at the current point
1802194d53e6SMatthew G. Knepley . t - current time
18032aa1fc23SMatthew G. Knepley . u_tShift - the multiplier a for dF/dU_t
1804194d53e6SMatthew G. Knepley . x - coordinates of the current point
1805194d53e6SMatthew G. Knepley . n - normal at the current point
1806194d53e6SMatthew G. Knepley - g0 - output values at the current point
1807194d53e6SMatthew G. Knepley 
1808194d53e6SMatthew G. Knepley   Level: intermediate
1809194d53e6SMatthew G. Knepley 
1810194d53e6SMatthew G. Knepley .seealso: PetscDSSetBdJacobian()
1811194d53e6SMatthew G. Knepley @*/
18122764a2aaSMatthew G. Knepley PetscErrorCode PetscDSGetBdJacobian(PetscDS prob, PetscInt f, PetscInt g,
181330b9ff8bSMatthew G. Knepley                                     void (**g0)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1814194d53e6SMatthew G. Knepley                                                 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1815194d53e6SMatthew G. Knepley                                                 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
18162aa1fc23SMatthew G. Knepley                                                 PetscReal t, PetscReal u_tShift, const PetscReal x[], const PetscReal n[], PetscScalar g0[]),
181730b9ff8bSMatthew G. Knepley                                     void (**g1)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1818194d53e6SMatthew G. Knepley                                                 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1819194d53e6SMatthew G. Knepley                                                 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
18202aa1fc23SMatthew G. Knepley                                                 PetscReal t, PetscReal u_tShift, const PetscReal x[], const PetscReal n[], PetscScalar g1[]),
182130b9ff8bSMatthew G. Knepley                                     void (**g2)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1822194d53e6SMatthew G. Knepley                                                 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1823194d53e6SMatthew G. Knepley                                                 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
18242aa1fc23SMatthew G. Knepley                                                 PetscReal t, PetscReal u_tShift, const PetscReal x[], const PetscReal n[], PetscScalar g2[]),
182530b9ff8bSMatthew G. Knepley                                     void (**g3)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1826194d53e6SMatthew G. Knepley                                                 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1827194d53e6SMatthew G. Knepley                                                 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
18282aa1fc23SMatthew G. Knepley                                                 PetscReal t, PetscReal u_tShift, const PetscReal x[], const PetscReal n[], PetscScalar g3[]))
18292764a2aaSMatthew G. Knepley {
18302764a2aaSMatthew G. Knepley   PetscFunctionBegin;
18312764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
18322764a2aaSMatthew 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);
18332764a2aaSMatthew 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);
18342764a2aaSMatthew G. Knepley   if (g0) {PetscValidPointer(g0, 4); *g0 = prob->gBd[(f*prob->Nf + g)*4+0];}
18352764a2aaSMatthew G. Knepley   if (g1) {PetscValidPointer(g1, 5); *g1 = prob->gBd[(f*prob->Nf + g)*4+1];}
18362764a2aaSMatthew G. Knepley   if (g2) {PetscValidPointer(g2, 6); *g2 = prob->gBd[(f*prob->Nf + g)*4+2];}
18372764a2aaSMatthew G. Knepley   if (g3) {PetscValidPointer(g3, 7); *g3 = prob->gBd[(f*prob->Nf + g)*4+3];}
18382764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
18392764a2aaSMatthew G. Knepley }
18402764a2aaSMatthew G. Knepley 
1841194d53e6SMatthew G. Knepley /*@C
1842194d53e6SMatthew G. Knepley   PetscDSSetBdJacobian - Set the pointwise boundary Jacobian function for given test and basis field
1843194d53e6SMatthew G. Knepley 
1844194d53e6SMatthew G. Knepley   Not collective
1845194d53e6SMatthew G. Knepley 
1846194d53e6SMatthew G. Knepley   Input Parameters:
1847194d53e6SMatthew G. Knepley + prob - The PetscDS
1848194d53e6SMatthew G. Knepley . f    - The test field number
1849194d53e6SMatthew G. Knepley . g    - The field number
1850194d53e6SMatthew G. Knepley . g0 - integrand for the test and basis function term
1851194d53e6SMatthew G. Knepley . g1 - integrand for the test function and basis function gradient term
1852194d53e6SMatthew G. Knepley . g2 - integrand for the test function gradient and basis function term
1853194d53e6SMatthew G. Knepley - g3 - integrand for the test function gradient and basis function gradient term
1854194d53e6SMatthew G. Knepley 
1855194d53e6SMatthew G. Knepley   Note: We are using a first order FEM model for the weak form:
1856194d53e6SMatthew G. Knepley 
1857194d53e6SMatthew G. Knepley   \int_\Gamma \phi {\vec g}_0(u, u_t, \nabla u, x, t) \cdot \hat n \psi + \phi {\vec g}_1(u, u_t, \nabla u, x, t) \cdot \hat n \nabla \psi + \nabla\phi \cdot {\vec g}_2(u, u_t, \nabla u, x, t) \cdot \hat n \psi + \nabla\phi \cdot {\overleftrightarrow g}_3(u, u_t, \nabla u, x, t) \cdot \hat n \cdot \nabla \psi
1858194d53e6SMatthew G. Knepley 
1859194d53e6SMatthew G. Knepley The calling sequence for the callbacks g0, g1, g2 and g3 is given by:
1860194d53e6SMatthew G. Knepley 
186130b9ff8bSMatthew G. Knepley $ g0(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1862194d53e6SMatthew G. Knepley $    const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1863194d53e6SMatthew G. Knepley $    const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
186430b9ff8bSMatthew G. Knepley $    PetscReal t, const PetscReal x[], const PetscReal n[], PetscScalar g0[])
1865194d53e6SMatthew G. Knepley 
1866194d53e6SMatthew G. Knepley + dim - the spatial dimension
1867194d53e6SMatthew G. Knepley . Nf - the number of fields
1868194d53e6SMatthew G. Knepley . uOff - the offset into u[] and u_t[] for each field
1869194d53e6SMatthew G. Knepley . uOff_x - the offset into u_x[] for each field
1870194d53e6SMatthew G. Knepley . u - each field evaluated at the current point
1871194d53e6SMatthew G. Knepley . u_t - the time derivative of each field evaluated at the current point
1872194d53e6SMatthew G. Knepley . u_x - the gradient of each field evaluated at the current point
1873194d53e6SMatthew G. Knepley . aOff - the offset into a[] and a_t[] for each auxiliary field
1874194d53e6SMatthew G. Knepley . aOff_x - the offset into a_x[] for each auxiliary field
1875194d53e6SMatthew G. Knepley . a - each auxiliary field evaluated at the current point
1876194d53e6SMatthew G. Knepley . a_t - the time derivative of each auxiliary field evaluated at the current point
1877194d53e6SMatthew G. Knepley . a_x - the gradient of auxiliary each field evaluated at the current point
1878194d53e6SMatthew G. Knepley . t - current time
18792aa1fc23SMatthew G. Knepley . u_tShift - the multiplier a for dF/dU_t
1880194d53e6SMatthew G. Knepley . x - coordinates of the current point
1881194d53e6SMatthew G. Knepley . n - normal at the current point
1882194d53e6SMatthew G. Knepley - g0 - output values at the current point
1883194d53e6SMatthew G. Knepley 
1884194d53e6SMatthew G. Knepley   Level: intermediate
1885194d53e6SMatthew G. Knepley 
1886194d53e6SMatthew G. Knepley .seealso: PetscDSGetBdJacobian()
1887194d53e6SMatthew G. Knepley @*/
18882764a2aaSMatthew G. Knepley PetscErrorCode PetscDSSetBdJacobian(PetscDS prob, PetscInt f, PetscInt g,
188930b9ff8bSMatthew G. Knepley                                     void (*g0)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1890194d53e6SMatthew G. Knepley                                                const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1891194d53e6SMatthew G. Knepley                                                const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
18922aa1fc23SMatthew G. Knepley                                                PetscReal t, PetscReal u_tShift, const PetscReal x[], const PetscReal n[], PetscScalar g0[]),
189330b9ff8bSMatthew G. Knepley                                     void (*g1)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1894194d53e6SMatthew G. Knepley                                                const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1895194d53e6SMatthew G. Knepley                                                const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
18962aa1fc23SMatthew G. Knepley                                                PetscReal t, PetscReal u_tShift, const PetscReal x[], const PetscReal n[], PetscScalar g1[]),
189730b9ff8bSMatthew G. Knepley                                     void (*g2)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1898194d53e6SMatthew G. Knepley                                                const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1899194d53e6SMatthew G. Knepley                                                const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
19002aa1fc23SMatthew G. Knepley                                                PetscReal t, PetscReal u_tShift, const PetscReal x[], const PetscReal n[], PetscScalar g2[]),
190130b9ff8bSMatthew G. Knepley                                     void (*g3)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1902194d53e6SMatthew G. Knepley                                                const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1903194d53e6SMatthew G. Knepley                                                const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
19042aa1fc23SMatthew G. Knepley                                                PetscReal t, PetscReal u_tShift, const PetscReal x[], const PetscReal n[], PetscScalar g3[]))
19052764a2aaSMatthew G. Knepley {
19062764a2aaSMatthew G. Knepley   PetscErrorCode ierr;
19072764a2aaSMatthew G. Knepley 
19082764a2aaSMatthew G. Knepley   PetscFunctionBegin;
19092764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
19102764a2aaSMatthew G. Knepley   if (g0) PetscValidFunction(g0, 4);
19112764a2aaSMatthew G. Knepley   if (g1) PetscValidFunction(g1, 5);
19122764a2aaSMatthew G. Knepley   if (g2) PetscValidFunction(g2, 6);
19132764a2aaSMatthew G. Knepley   if (g3) PetscValidFunction(g3, 7);
19142764a2aaSMatthew G. Knepley   if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f);
19152764a2aaSMatthew G. Knepley   if (g < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", g);
19162764a2aaSMatthew G. Knepley   ierr = PetscDSEnlarge_Static(prob, PetscMax(f, g)+1);CHKERRQ(ierr);
19172764a2aaSMatthew G. Knepley   prob->gBd[(f*prob->Nf + g)*4+0] = g0;
19182764a2aaSMatthew G. Knepley   prob->gBd[(f*prob->Nf + g)*4+1] = g1;
19192764a2aaSMatthew G. Knepley   prob->gBd[(f*prob->Nf + g)*4+2] = g2;
19202764a2aaSMatthew G. Knepley   prob->gBd[(f*prob->Nf + g)*4+3] = g3;
19212764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
19222764a2aaSMatthew G. Knepley }
19232764a2aaSMatthew G. Knepley 
1924*c371a6d1SMatthew G. Knepley /*@C
1925*c371a6d1SMatthew G. Knepley   PetscDSGetExactSolution - Get the pointwise exact solution function for a given test field
1926*c371a6d1SMatthew G. Knepley 
1927*c371a6d1SMatthew G. Knepley   Not collective
1928*c371a6d1SMatthew G. Knepley 
1929*c371a6d1SMatthew G. Knepley   Input Parameters:
1930*c371a6d1SMatthew G. Knepley + prob - The PetscDS
1931*c371a6d1SMatthew G. Knepley - f    - The test field number
1932*c371a6d1SMatthew G. Knepley 
1933*c371a6d1SMatthew G. Knepley   Output Parameter:
1934*c371a6d1SMatthew G. Knepley . exactSol - exact solution for the test field
1935*c371a6d1SMatthew G. Knepley 
1936*c371a6d1SMatthew G. Knepley   Note: The calling sequence for the solution functions is given by:
1937*c371a6d1SMatthew G. Knepley 
1938*c371a6d1SMatthew G. Knepley $ sol(PetscInt dim, PetscReal t, const PetscReal x[], PetscInt Nc, PetscScalar u[], void *ctx)
1939*c371a6d1SMatthew G. Knepley 
1940*c371a6d1SMatthew G. Knepley + dim - the spatial dimension
1941*c371a6d1SMatthew G. Knepley . t - current time
1942*c371a6d1SMatthew G. Knepley . x - coordinates of the current point
1943*c371a6d1SMatthew G. Knepley . Nc - the number of field components
1944*c371a6d1SMatthew G. Knepley . u - the solution field evaluated at the current point
1945*c371a6d1SMatthew G. Knepley - ctx - a user context
1946*c371a6d1SMatthew G. Knepley 
1947*c371a6d1SMatthew G. Knepley   Level: intermediate
1948*c371a6d1SMatthew G. Knepley 
1949*c371a6d1SMatthew G. Knepley .seealso: PetscDSSetExactSolution()
1950*c371a6d1SMatthew G. Knepley @*/
1951*c371a6d1SMatthew G. Knepley PetscErrorCode PetscDSGetExactSolution(PetscDS prob, PetscInt f, PetscErrorCode (**sol)(PetscInt dim, PetscReal t, const PetscReal x[], PetscInt Nc, PetscScalar u[], void *ctx))
1952*c371a6d1SMatthew G. Knepley {
1953*c371a6d1SMatthew G. Knepley   PetscFunctionBegin;
1954*c371a6d1SMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1955*c371a6d1SMatthew 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);
1956*c371a6d1SMatthew G. Knepley   if (sol) {PetscValidPointer(sol, 3); *sol = prob->exactSol[f];}
1957*c371a6d1SMatthew G. Knepley   PetscFunctionReturn(0);
1958*c371a6d1SMatthew G. Knepley }
1959*c371a6d1SMatthew G. Knepley 
1960*c371a6d1SMatthew G. Knepley /*@C
1961*c371a6d1SMatthew G. Knepley   PetscDSSetExactSolution - Get the pointwise exact solution function for a given test field
1962*c371a6d1SMatthew G. Knepley 
1963*c371a6d1SMatthew G. Knepley   Not collective
1964*c371a6d1SMatthew G. Knepley 
1965*c371a6d1SMatthew G. Knepley   Input Parameters:
1966*c371a6d1SMatthew G. Knepley + prob - The PetscDS
1967*c371a6d1SMatthew G. Knepley . f    - The test field number
1968*c371a6d1SMatthew G. Knepley - sol  - solution function for the test fields
1969*c371a6d1SMatthew G. Knepley 
1970*c371a6d1SMatthew G. Knepley   Note: The calling sequence for solution functions is given by:
1971*c371a6d1SMatthew G. Knepley 
1972*c371a6d1SMatthew G. Knepley $ sol(PetscInt dim, PetscReal t, const PetscReal x[], PetscInt Nc, PetscScalar u[], void *ctx)
1973*c371a6d1SMatthew G. Knepley 
1974*c371a6d1SMatthew G. Knepley + dim - the spatial dimension
1975*c371a6d1SMatthew G. Knepley . t - current time
1976*c371a6d1SMatthew G. Knepley . x - coordinates of the current point
1977*c371a6d1SMatthew G. Knepley . Nc - the number of field components
1978*c371a6d1SMatthew G. Knepley . u - the solution field evaluated at the current point
1979*c371a6d1SMatthew G. Knepley - ctx - a user context
1980*c371a6d1SMatthew G. Knepley 
1981*c371a6d1SMatthew G. Knepley   Level: intermediate
1982*c371a6d1SMatthew G. Knepley 
1983*c371a6d1SMatthew G. Knepley .seealso: PetscDSGetExactSolution()
1984*c371a6d1SMatthew G. Knepley @*/
1985*c371a6d1SMatthew G. Knepley PetscErrorCode PetscDSSetExactSolution(PetscDS prob, PetscInt f, PetscErrorCode (*sol)(PetscInt dim, PetscReal t, const PetscReal x[], PetscInt Nc, PetscScalar u[], void *ctx))
1986*c371a6d1SMatthew G. Knepley {
1987*c371a6d1SMatthew G. Knepley   PetscErrorCode ierr;
1988*c371a6d1SMatthew G. Knepley 
1989*c371a6d1SMatthew G. Knepley   PetscFunctionBegin;
1990*c371a6d1SMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1991*c371a6d1SMatthew G. Knepley   if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f);
1992*c371a6d1SMatthew G. Knepley   ierr = PetscDSEnlarge_Static(prob, f+1);CHKERRQ(ierr);
1993*c371a6d1SMatthew G. Knepley   if (sol) {PetscValidFunction(sol, 3); prob->exactSol[f] = sol;}
1994*c371a6d1SMatthew G. Knepley   PetscFunctionReturn(0);
1995*c371a6d1SMatthew G. Knepley }
1996*c371a6d1SMatthew G. Knepley 
19974cd1e086SMatthew G. Knepley /*@
19984cd1e086SMatthew G. Knepley   PetscDSGetFieldIndex - Returns the index of the given field
19994cd1e086SMatthew G. Knepley 
20004cd1e086SMatthew G. Knepley   Not collective
20014cd1e086SMatthew G. Knepley 
20024cd1e086SMatthew G. Knepley   Input Parameters:
20034cd1e086SMatthew G. Knepley + prob - The PetscDS object
20044cd1e086SMatthew G. Knepley - disc - The discretization object
20054cd1e086SMatthew G. Knepley 
20064cd1e086SMatthew G. Knepley   Output Parameter:
20074cd1e086SMatthew G. Knepley . f - The field number
20084cd1e086SMatthew G. Knepley 
20094cd1e086SMatthew G. Knepley   Level: beginner
20104cd1e086SMatthew G. Knepley 
2011f744cafaSSander Arens .seealso: PetscGetDiscretization(), PetscDSGetNumFields(), PetscDSCreate()
20124cd1e086SMatthew G. Knepley @*/
20134cd1e086SMatthew G. Knepley PetscErrorCode PetscDSGetFieldIndex(PetscDS prob, PetscObject disc, PetscInt *f)
20144cd1e086SMatthew G. Knepley {
20154cd1e086SMatthew G. Knepley   PetscInt g;
20164cd1e086SMatthew G. Knepley 
20174cd1e086SMatthew G. Knepley   PetscFunctionBegin;
20184cd1e086SMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
20194cd1e086SMatthew G. Knepley   PetscValidPointer(f, 3);
20204cd1e086SMatthew G. Knepley   *f = -1;
20214cd1e086SMatthew G. Knepley   for (g = 0; g < prob->Nf; ++g) {if (disc == prob->disc[g]) break;}
20224cd1e086SMatthew G. Knepley   if (g == prob->Nf) SETERRQ(PetscObjectComm((PetscObject) prob), PETSC_ERR_ARG_WRONG, "Field not found in PetscDS.");
20234cd1e086SMatthew G. Knepley   *f = g;
20244cd1e086SMatthew G. Knepley   PetscFunctionReturn(0);
20254cd1e086SMatthew G. Knepley }
20264cd1e086SMatthew G. Knepley 
20274cd1e086SMatthew G. Knepley /*@
20284cd1e086SMatthew G. Knepley   PetscDSGetFieldSize - Returns the size of the given field in the full space basis
20294cd1e086SMatthew G. Knepley 
20304cd1e086SMatthew G. Knepley   Not collective
20314cd1e086SMatthew G. Knepley 
20324cd1e086SMatthew G. Knepley   Input Parameters:
20334cd1e086SMatthew G. Knepley + prob - The PetscDS object
20344cd1e086SMatthew G. Knepley - f - The field number
20354cd1e086SMatthew G. Knepley 
20364cd1e086SMatthew G. Knepley   Output Parameter:
20374cd1e086SMatthew G. Knepley . size - The size
20384cd1e086SMatthew G. Knepley 
20394cd1e086SMatthew G. Knepley   Level: beginner
20404cd1e086SMatthew G. Knepley 
2041f744cafaSSander Arens .seealso: PetscDSGetFieldOffset(), PetscDSGetNumFields(), PetscDSCreate()
20424cd1e086SMatthew G. Knepley @*/
20434cd1e086SMatthew G. Knepley PetscErrorCode PetscDSGetFieldSize(PetscDS prob, PetscInt f, PetscInt *size)
20444cd1e086SMatthew G. Knepley {
20452166fd64SMatthew G. Knepley   PetscErrorCode ierr;
20462166fd64SMatthew G. Knepley 
20474cd1e086SMatthew G. Knepley   PetscFunctionBegin;
20484cd1e086SMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
20494cd1e086SMatthew G. Knepley   PetscValidPointer(size, 3);
20504cd1e086SMatthew 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);
20512166fd64SMatthew G. Knepley   ierr = PetscDSSetUp(prob);CHKERRQ(ierr);
2052d4742ddaSMatthew G. Knepley   *size = prob->Nb[f];
20534cd1e086SMatthew G. Knepley   PetscFunctionReturn(0);
20544cd1e086SMatthew G. Knepley }
20554cd1e086SMatthew G. Knepley 
2056bc4ae4beSMatthew G. Knepley /*@
2057bc4ae4beSMatthew G. Knepley   PetscDSGetFieldOffset - Returns the offset of the given field in the full space basis
2058bc4ae4beSMatthew G. Knepley 
2059bc4ae4beSMatthew G. Knepley   Not collective
2060bc4ae4beSMatthew G. Knepley 
2061bc4ae4beSMatthew G. Knepley   Input Parameters:
2062bc4ae4beSMatthew G. Knepley + prob - The PetscDS object
2063bc4ae4beSMatthew G. Knepley - f - The field number
2064bc4ae4beSMatthew G. Knepley 
2065bc4ae4beSMatthew G. Knepley   Output Parameter:
2066bc4ae4beSMatthew G. Knepley . off - The offset
2067bc4ae4beSMatthew G. Knepley 
2068bc4ae4beSMatthew G. Knepley   Level: beginner
2069bc4ae4beSMatthew G. Knepley 
2070f744cafaSSander Arens .seealso: PetscDSGetFieldSize(), PetscDSGetNumFields(), PetscDSCreate()
2071bc4ae4beSMatthew G. Knepley @*/
20722764a2aaSMatthew G. Knepley PetscErrorCode PetscDSGetFieldOffset(PetscDS prob, PetscInt f, PetscInt *off)
20732764a2aaSMatthew G. Knepley {
20744cd1e086SMatthew G. Knepley   PetscInt       size, g;
20752764a2aaSMatthew G. Knepley   PetscErrorCode ierr;
20762764a2aaSMatthew G. Knepley 
20772764a2aaSMatthew G. Knepley   PetscFunctionBegin;
20782764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
20792764a2aaSMatthew G. Knepley   PetscValidPointer(off, 3);
20802764a2aaSMatthew 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);
20812764a2aaSMatthew G. Knepley   *off = 0;
20822764a2aaSMatthew G. Knepley   for (g = 0; g < f; ++g) {
20834cd1e086SMatthew G. Knepley     ierr = PetscDSGetFieldSize(prob, g, &size);CHKERRQ(ierr);
20844cd1e086SMatthew G. Knepley     *off += size;
20852764a2aaSMatthew G. Knepley   }
20862764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
20872764a2aaSMatthew G. Knepley }
20882764a2aaSMatthew G. Knepley 
2089bc4ae4beSMatthew G. Knepley /*@
209047e57110SSander Arens   PetscDSGetDimensions - Returns the size of the approximation space for each field on an evaluation point
2091bc4ae4beSMatthew G. Knepley 
2092bc4ae4beSMatthew G. Knepley   Not collective
2093bc4ae4beSMatthew G. Knepley 
209447e57110SSander Arens   Input Parameter:
209547e57110SSander Arens . prob - The PetscDS object
2096bc4ae4beSMatthew G. Knepley 
2097bc4ae4beSMatthew G. Knepley   Output Parameter:
209847e57110SSander Arens . dimensions - The number of dimensions
2099bc4ae4beSMatthew G. Knepley 
2100bc4ae4beSMatthew G. Knepley   Level: beginner
2101bc4ae4beSMatthew G. Knepley 
210247e57110SSander Arens .seealso: PetscDSGetComponentOffsets(), PetscDSGetNumFields(), PetscDSCreate()
2103bc4ae4beSMatthew G. Knepley @*/
210447e57110SSander Arens PetscErrorCode PetscDSGetDimensions(PetscDS prob, PetscInt *dimensions[])
21052764a2aaSMatthew G. Knepley {
21062764a2aaSMatthew G. Knepley   PetscErrorCode ierr;
21072764a2aaSMatthew G. Knepley 
21082764a2aaSMatthew G. Knepley   PetscFunctionBegin;
21092764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
211047e57110SSander Arens   ierr = PetscDSSetUp(prob);CHKERRQ(ierr);
211147e57110SSander Arens   PetscValidPointer(dimensions, 2);
211247e57110SSander Arens   *dimensions = prob->Nb;
211347e57110SSander Arens   PetscFunctionReturn(0);
21146ce16762SMatthew G. Knepley }
211547e57110SSander Arens 
211647e57110SSander Arens /*@
211747e57110SSander Arens   PetscDSGetComponents - Returns the number of components for each field on an evaluation point
211847e57110SSander Arens 
211947e57110SSander Arens   Not collective
212047e57110SSander Arens 
212147e57110SSander Arens   Input Parameter:
212247e57110SSander Arens . prob - The PetscDS object
212347e57110SSander Arens 
212447e57110SSander Arens   Output Parameter:
212547e57110SSander Arens . components - The number of components
212647e57110SSander Arens 
212747e57110SSander Arens   Level: beginner
212847e57110SSander Arens 
212947e57110SSander Arens .seealso: PetscDSGetComponentOffsets(), PetscDSGetNumFields(), PetscDSCreate()
213047e57110SSander Arens @*/
213147e57110SSander Arens PetscErrorCode PetscDSGetComponents(PetscDS prob, PetscInt *components[])
213247e57110SSander Arens {
213347e57110SSander Arens   PetscErrorCode ierr;
213447e57110SSander Arens 
213547e57110SSander Arens   PetscFunctionBegin;
213647e57110SSander Arens   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
213747e57110SSander Arens   ierr = PetscDSSetUp(prob);CHKERRQ(ierr);
213847e57110SSander Arens   PetscValidPointer(components, 2);
213947e57110SSander Arens   *components = prob->Nc;
21406ce16762SMatthew G. Knepley   PetscFunctionReturn(0);
21416ce16762SMatthew G. Knepley }
21426ce16762SMatthew G. Knepley 
21436ce16762SMatthew G. Knepley /*@
21446ce16762SMatthew G. Knepley   PetscDSGetComponentOffset - Returns the offset of the given field on an evaluation point
21456ce16762SMatthew G. Knepley 
21466ce16762SMatthew G. Knepley   Not collective
21476ce16762SMatthew G. Knepley 
21486ce16762SMatthew G. Knepley   Input Parameters:
21496ce16762SMatthew G. Knepley + prob - The PetscDS object
21506ce16762SMatthew G. Knepley - f - The field number
21516ce16762SMatthew G. Knepley 
21526ce16762SMatthew G. Knepley   Output Parameter:
21536ce16762SMatthew G. Knepley . off - The offset
21546ce16762SMatthew G. Knepley 
21556ce16762SMatthew G. Knepley   Level: beginner
21566ce16762SMatthew G. Knepley 
2157f744cafaSSander Arens .seealso: PetscDSGetNumFields(), PetscDSCreate()
21586ce16762SMatthew G. Knepley @*/
21596ce16762SMatthew G. Knepley PetscErrorCode PetscDSGetComponentOffset(PetscDS prob, PetscInt f, PetscInt *off)
21606ce16762SMatthew G. Knepley {
21616ce16762SMatthew G. Knepley   PetscFunctionBegin;
21626ce16762SMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
21636ce16762SMatthew G. Knepley   PetscValidPointer(off, 3);
21646ce16762SMatthew 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);
216547e57110SSander Arens   *off = prob->off[f];
21662764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
21672764a2aaSMatthew G. Knepley }
21682764a2aaSMatthew G. Knepley 
2169194d53e6SMatthew G. Knepley /*@
2170194d53e6SMatthew G. Knepley   PetscDSGetComponentOffsets - Returns the offset of each field on an evaluation point
2171194d53e6SMatthew G. Knepley 
2172194d53e6SMatthew G. Knepley   Not collective
2173194d53e6SMatthew G. Knepley 
2174194d53e6SMatthew G. Knepley   Input Parameter:
2175194d53e6SMatthew G. Knepley . prob - The PetscDS object
2176194d53e6SMatthew G. Knepley 
2177194d53e6SMatthew G. Knepley   Output Parameter:
2178194d53e6SMatthew G. Knepley . offsets - The offsets
2179194d53e6SMatthew G. Knepley 
2180194d53e6SMatthew G. Knepley   Level: beginner
2181194d53e6SMatthew G. Knepley 
2182f744cafaSSander Arens .seealso: PetscDSGetNumFields(), PetscDSCreate()
2183194d53e6SMatthew G. Knepley @*/
2184194d53e6SMatthew G. Knepley PetscErrorCode PetscDSGetComponentOffsets(PetscDS prob, PetscInt *offsets[])
2185194d53e6SMatthew G. Knepley {
2186194d53e6SMatthew G. Knepley   PetscFunctionBegin;
2187194d53e6SMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
2188194d53e6SMatthew G. Knepley   PetscValidPointer(offsets, 2);
2189194d53e6SMatthew G. Knepley   *offsets = prob->off;
2190194d53e6SMatthew G. Knepley   PetscFunctionReturn(0);
2191194d53e6SMatthew G. Knepley }
2192194d53e6SMatthew G. Knepley 
2193194d53e6SMatthew G. Knepley /*@
2194194d53e6SMatthew G. Knepley   PetscDSGetComponentDerivativeOffsets - Returns the offset of each field derivative on an evaluation point
2195194d53e6SMatthew G. Knepley 
2196194d53e6SMatthew G. Knepley   Not collective
2197194d53e6SMatthew G. Knepley 
2198194d53e6SMatthew G. Knepley   Input Parameter:
2199194d53e6SMatthew G. Knepley . prob - The PetscDS object
2200194d53e6SMatthew G. Knepley 
2201194d53e6SMatthew G. Knepley   Output Parameter:
2202194d53e6SMatthew G. Knepley . offsets - The offsets
2203194d53e6SMatthew G. Knepley 
2204194d53e6SMatthew G. Knepley   Level: beginner
2205194d53e6SMatthew G. Knepley 
2206f744cafaSSander Arens .seealso: PetscDSGetNumFields(), PetscDSCreate()
2207194d53e6SMatthew G. Knepley @*/
2208194d53e6SMatthew G. Knepley PetscErrorCode PetscDSGetComponentDerivativeOffsets(PetscDS prob, PetscInt *offsets[])
2209194d53e6SMatthew G. Knepley {
2210194d53e6SMatthew G. Knepley   PetscFunctionBegin;
2211194d53e6SMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
2212194d53e6SMatthew G. Knepley   PetscValidPointer(offsets, 2);
2213194d53e6SMatthew G. Knepley   *offsets = prob->offDer;
2214194d53e6SMatthew G. Knepley   PetscFunctionReturn(0);
2215194d53e6SMatthew G. Knepley }
2216194d53e6SMatthew G. Knepley 
221768c9edb9SMatthew G. Knepley /*@C
221868c9edb9SMatthew G. Knepley   PetscDSGetTabulation - Return the basis tabulation at quadrature points for the volume discretization
221968c9edb9SMatthew G. Knepley 
222068c9edb9SMatthew G. Knepley   Not collective
222168c9edb9SMatthew G. Knepley 
222268c9edb9SMatthew G. Knepley   Input Parameter:
222368c9edb9SMatthew G. Knepley . prob - The PetscDS object
222468c9edb9SMatthew G. Knepley 
222568c9edb9SMatthew G. Knepley   Output Parameters:
222668c9edb9SMatthew G. Knepley + basis - The basis function tabulation at quadrature points
222768c9edb9SMatthew G. Knepley - basisDer - The basis function derivative tabulation at quadrature points
222868c9edb9SMatthew G. Knepley 
222968c9edb9SMatthew G. Knepley   Level: intermediate
223068c9edb9SMatthew G. Knepley 
2231f744cafaSSander Arens .seealso: PetscDSCreate()
223268c9edb9SMatthew G. Knepley @*/
22332764a2aaSMatthew G. Knepley PetscErrorCode PetscDSGetTabulation(PetscDS prob, PetscReal ***basis, PetscReal ***basisDer)
22342764a2aaSMatthew G. Knepley {
22352764a2aaSMatthew G. Knepley   PetscErrorCode ierr;
22362764a2aaSMatthew G. Knepley 
22372764a2aaSMatthew G. Knepley   PetscFunctionBegin;
22382764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
22392764a2aaSMatthew G. Knepley   ierr = PetscDSSetUp(prob);CHKERRQ(ierr);
22402764a2aaSMatthew G. Knepley   if (basis)    {PetscValidPointer(basis, 2);    *basis    = prob->basis;}
22412764a2aaSMatthew G. Knepley   if (basisDer) {PetscValidPointer(basisDer, 3); *basisDer = prob->basisDer;}
22422764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
22432764a2aaSMatthew G. Knepley }
22442764a2aaSMatthew G. Knepley 
224568c9edb9SMatthew G. Knepley /*@C
22464d0b9603SSander Arens   PetscDSGetFaceTabulation - Return the basis tabulation at quadrature points on the faces
224768c9edb9SMatthew G. Knepley 
224868c9edb9SMatthew G. Knepley   Not collective
224968c9edb9SMatthew G. Knepley 
225068c9edb9SMatthew G. Knepley   Input Parameter:
225168c9edb9SMatthew G. Knepley . prob - The PetscDS object
225268c9edb9SMatthew G. Knepley 
225368c9edb9SMatthew G. Knepley   Output Parameters:
22544d0b9603SSander Arens + basisFace - The basis function tabulation at quadrature points
22554d0b9603SSander Arens - basisDerFace - The basis function derivative tabulation at quadrature points
225668c9edb9SMatthew G. Knepley 
225768c9edb9SMatthew G. Knepley   Level: intermediate
225868c9edb9SMatthew G. Knepley 
225968c9edb9SMatthew G. Knepley .seealso: PetscDSGetTabulation(), PetscDSCreate()
226068c9edb9SMatthew G. Knepley @*/
22614d0b9603SSander Arens PetscErrorCode PetscDSGetFaceTabulation(PetscDS prob, PetscReal ***basis, PetscReal ***basisDer)
22622764a2aaSMatthew G. Knepley {
22632764a2aaSMatthew G. Knepley   PetscErrorCode ierr;
22642764a2aaSMatthew G. Knepley 
22652764a2aaSMatthew G. Knepley   PetscFunctionBegin;
22662764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
22672764a2aaSMatthew G. Knepley   ierr = PetscDSSetUp(prob);CHKERRQ(ierr);
22684d0b9603SSander Arens   if (basis)    {PetscValidPointer(basis, 2);    *basis    = prob->basisFace;}
22694d0b9603SSander Arens   if (basisDer) {PetscValidPointer(basisDer, 3); *basisDer = prob->basisDerFace;}
22702764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
22712764a2aaSMatthew G. Knepley }
22722764a2aaSMatthew G. Knepley 
22732764a2aaSMatthew G. Knepley PetscErrorCode PetscDSGetEvaluationArrays(PetscDS prob, PetscScalar **u, PetscScalar **u_t, PetscScalar **u_x)
22742764a2aaSMatthew G. Knepley {
22752764a2aaSMatthew G. Knepley   PetscErrorCode ierr;
22762764a2aaSMatthew G. Knepley 
22772764a2aaSMatthew G. Knepley   PetscFunctionBegin;
22782764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
22792764a2aaSMatthew G. Knepley   ierr = PetscDSSetUp(prob);CHKERRQ(ierr);
22802764a2aaSMatthew G. Knepley   if (u)   {PetscValidPointer(u, 2);   *u   = prob->u;}
22812764a2aaSMatthew G. Knepley   if (u_t) {PetscValidPointer(u_t, 3); *u_t = prob->u_t;}
22822764a2aaSMatthew G. Knepley   if (u_x) {PetscValidPointer(u_x, 4); *u_x = prob->u_x;}
22832764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
22842764a2aaSMatthew G. Knepley }
22852764a2aaSMatthew G. Knepley 
22862764a2aaSMatthew G. Knepley PetscErrorCode PetscDSGetWeakFormArrays(PetscDS prob, PetscScalar **f0, PetscScalar **f1, PetscScalar **g0, PetscScalar **g1, PetscScalar **g2, PetscScalar **g3)
22872764a2aaSMatthew G. Knepley {
22882764a2aaSMatthew G. Knepley   PetscErrorCode ierr;
22892764a2aaSMatthew G. Knepley 
22902764a2aaSMatthew G. Knepley   PetscFunctionBegin;
22912764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
22922764a2aaSMatthew G. Knepley   ierr = PetscDSSetUp(prob);CHKERRQ(ierr);
22932764a2aaSMatthew G. Knepley   if (f0) {PetscValidPointer(f0, 2); *f0 = prob->f0;}
22942764a2aaSMatthew G. Knepley   if (f1) {PetscValidPointer(f1, 3); *f1 = prob->f1;}
22952764a2aaSMatthew G. Knepley   if (g0) {PetscValidPointer(g0, 4); *g0 = prob->g0;}
22962764a2aaSMatthew G. Knepley   if (g1) {PetscValidPointer(g1, 5); *g1 = prob->g1;}
22972764a2aaSMatthew G. Knepley   if (g2) {PetscValidPointer(g2, 6); *g2 = prob->g2;}
22982764a2aaSMatthew G. Knepley   if (g3) {PetscValidPointer(g3, 7); *g3 = prob->g3;}
22992764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
23002764a2aaSMatthew G. Knepley }
23012764a2aaSMatthew G. Knepley 
23022764a2aaSMatthew G. Knepley PetscErrorCode PetscDSGetRefCoordArrays(PetscDS prob, PetscReal **x, PetscScalar **refSpaceDer)
23032764a2aaSMatthew G. Knepley {
23042764a2aaSMatthew G. Knepley   PetscErrorCode ierr;
23052764a2aaSMatthew G. Knepley 
23062764a2aaSMatthew G. Knepley   PetscFunctionBegin;
23072764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
23082764a2aaSMatthew G. Knepley   ierr = PetscDSSetUp(prob);CHKERRQ(ierr);
23092764a2aaSMatthew G. Knepley   if (x)           {PetscValidPointer(x, 2);           *x           = prob->x;}
23102764a2aaSMatthew G. Knepley   if (refSpaceDer) {PetscValidPointer(refSpaceDer, 3); *refSpaceDer = prob->refSpaceDer;}
23112764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
23122764a2aaSMatthew G. Knepley }
23132764a2aaSMatthew G. Knepley 
231458ebd649SToby Isaac /*@C
231558ebd649SToby Isaac   PetscDSAddBoundary - Add a boundary condition to the model
231658ebd649SToby Isaac 
231758ebd649SToby Isaac   Input Parameters:
231858ebd649SToby Isaac + ds          - The PetscDS object
23192d47a189SJulian Andrej . type        - The type of condition, e.g. DM_BC_ESSENTIAL/DM_BC_ESSENTIAL_FIELD (Dirichlet), or DM_BC_NATURAL (Neumann)
232058ebd649SToby Isaac . name        - The BC name
232158ebd649SToby Isaac . labelname   - The label defining constrained points
232258ebd649SToby Isaac . field       - The field to constrain
232358ebd649SToby Isaac . numcomps    - The number of constrained field components
232458ebd649SToby Isaac . comps       - An array of constrained component numbers
232558ebd649SToby Isaac . bcFunc      - A pointwise function giving boundary values
232658ebd649SToby Isaac . numids      - The number of DMLabel ids for constrained points
232758ebd649SToby Isaac . ids         - An array of ids for constrained points
232858ebd649SToby Isaac - ctx         - An optional user context for bcFunc
232958ebd649SToby Isaac 
233058ebd649SToby Isaac   Options Database Keys:
233158ebd649SToby Isaac + -bc_<boundary name> <num> - Overrides the boundary ids
233258ebd649SToby Isaac - -bc_<boundary name>_comp <num> - Overrides the boundary components
233358ebd649SToby Isaac 
233458ebd649SToby Isaac   Level: developer
233558ebd649SToby Isaac 
233658ebd649SToby Isaac .seealso: PetscDSGetBoundary()
233758ebd649SToby Isaac @*/
2338a30ec4eaSSatish Balay PetscErrorCode PetscDSAddBoundary(PetscDS ds, DMBoundaryConditionType type, const char name[], const char labelname[], PetscInt field, PetscInt numcomps, const PetscInt *comps, void (*bcFunc)(void), PetscInt numids, const PetscInt *ids, void *ctx)
233958ebd649SToby Isaac {
234058ebd649SToby Isaac   DSBoundary     b;
234158ebd649SToby Isaac   PetscErrorCode ierr;
234258ebd649SToby Isaac 
234358ebd649SToby Isaac   PetscFunctionBegin;
234458ebd649SToby Isaac   PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1);
234558ebd649SToby Isaac   ierr = PetscNew(&b);CHKERRQ(ierr);
234658ebd649SToby Isaac   ierr = PetscStrallocpy(name, (char **) &b->name);CHKERRQ(ierr);
234758ebd649SToby Isaac   ierr = PetscStrallocpy(labelname, (char **) &b->labelname);CHKERRQ(ierr);
234858ebd649SToby Isaac   ierr = PetscMalloc1(numcomps, &b->comps);CHKERRQ(ierr);
234958ebd649SToby Isaac   if (numcomps) {ierr = PetscMemcpy(b->comps, comps, numcomps*sizeof(PetscInt));CHKERRQ(ierr);}
235058ebd649SToby Isaac   ierr = PetscMalloc1(numids, &b->ids);CHKERRQ(ierr);
235158ebd649SToby Isaac   if (numids) {ierr = PetscMemcpy(b->ids, ids, numids*sizeof(PetscInt));CHKERRQ(ierr);}
2352f971fd6bSMatthew G. Knepley   b->type            = type;
235358ebd649SToby Isaac   b->field           = field;
235458ebd649SToby Isaac   b->numcomps        = numcomps;
235558ebd649SToby Isaac   b->func            = bcFunc;
235658ebd649SToby Isaac   b->numids          = numids;
235758ebd649SToby Isaac   b->ctx             = ctx;
235858ebd649SToby Isaac   b->next            = ds->boundary;
235958ebd649SToby Isaac   ds->boundary       = b;
236058ebd649SToby Isaac   PetscFunctionReturn(0);
236158ebd649SToby Isaac }
236258ebd649SToby Isaac 
236358ebd649SToby Isaac /*@
236458ebd649SToby Isaac   PetscDSGetNumBoundary - Get the number of registered BC
236558ebd649SToby Isaac 
236658ebd649SToby Isaac   Input Parameters:
236758ebd649SToby Isaac . ds - The PetscDS object
236858ebd649SToby Isaac 
236958ebd649SToby Isaac   Output Parameters:
237058ebd649SToby Isaac . numBd - The number of BC
237158ebd649SToby Isaac 
237258ebd649SToby Isaac   Level: intermediate
237358ebd649SToby Isaac 
237458ebd649SToby Isaac .seealso: PetscDSAddBoundary(), PetscDSGetBoundary()
237558ebd649SToby Isaac @*/
237658ebd649SToby Isaac PetscErrorCode PetscDSGetNumBoundary(PetscDS ds, PetscInt *numBd)
237758ebd649SToby Isaac {
237858ebd649SToby Isaac   DSBoundary b = ds->boundary;
237958ebd649SToby Isaac 
238058ebd649SToby Isaac   PetscFunctionBegin;
238158ebd649SToby Isaac   PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1);
238258ebd649SToby Isaac   PetscValidPointer(numBd, 2);
238358ebd649SToby Isaac   *numBd = 0;
238458ebd649SToby Isaac   while (b) {++(*numBd); b = b->next;}
238558ebd649SToby Isaac   PetscFunctionReturn(0);
238658ebd649SToby Isaac }
238758ebd649SToby Isaac 
238858ebd649SToby Isaac /*@C
238958ebd649SToby Isaac   PetscDSGetBoundary - Add a boundary condition to the model
239058ebd649SToby Isaac 
239158ebd649SToby Isaac   Input Parameters:
239258ebd649SToby Isaac + ds          - The PetscDS object
239358ebd649SToby Isaac - bd          - The BC number
239458ebd649SToby Isaac 
239558ebd649SToby Isaac   Output Parameters:
23962d47a189SJulian Andrej + type        - The type of condition, e.g. DM_BC_ESSENTIAL/DM_BC_ESSENTIAL_FIELD (Dirichlet), or DM_BC_NATURAL (Neumann)
239758ebd649SToby Isaac . name        - The BC name
239858ebd649SToby Isaac . labelname   - The label defining constrained points
239958ebd649SToby Isaac . field       - The field to constrain
240058ebd649SToby Isaac . numcomps    - The number of constrained field components
240158ebd649SToby Isaac . comps       - An array of constrained component numbers
240258ebd649SToby Isaac . bcFunc      - A pointwise function giving boundary values
240358ebd649SToby Isaac . numids      - The number of DMLabel ids for constrained points
240458ebd649SToby Isaac . ids         - An array of ids for constrained points
240558ebd649SToby Isaac - ctx         - An optional user context for bcFunc
240658ebd649SToby Isaac 
240758ebd649SToby Isaac   Options Database Keys:
240858ebd649SToby Isaac + -bc_<boundary name> <num> - Overrides the boundary ids
240958ebd649SToby Isaac - -bc_<boundary name>_comp <num> - Overrides the boundary components
241058ebd649SToby Isaac 
241158ebd649SToby Isaac   Level: developer
241258ebd649SToby Isaac 
241358ebd649SToby Isaac .seealso: PetscDSAddBoundary()
241458ebd649SToby Isaac @*/
2415a30ec4eaSSatish Balay PetscErrorCode PetscDSGetBoundary(PetscDS ds, PetscInt bd, DMBoundaryConditionType *type, const char **name, const char **labelname, PetscInt *field, PetscInt *numcomps, const PetscInt **comps, void (**func)(void), PetscInt *numids, const PetscInt **ids, void **ctx)
241658ebd649SToby Isaac {
241758ebd649SToby Isaac   DSBoundary b    = ds->boundary;
241858ebd649SToby Isaac   PetscInt   n    = 0;
241958ebd649SToby Isaac 
242058ebd649SToby Isaac   PetscFunctionBegin;
242158ebd649SToby Isaac   PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1);
242258ebd649SToby Isaac   while (b) {
242358ebd649SToby Isaac     if (n == bd) break;
242458ebd649SToby Isaac     b = b->next;
242558ebd649SToby Isaac     ++n;
242658ebd649SToby Isaac   }
242758ebd649SToby Isaac   if (!b) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Boundary %d is not in [0, %d)", bd, n);
2428f971fd6bSMatthew G. Knepley   if (type) {
2429f971fd6bSMatthew G. Knepley     PetscValidPointer(type, 3);
2430f971fd6bSMatthew G. Knepley     *type = b->type;
243158ebd649SToby Isaac   }
243258ebd649SToby Isaac   if (name) {
243358ebd649SToby Isaac     PetscValidPointer(name, 4);
243458ebd649SToby Isaac     *name = b->name;
243558ebd649SToby Isaac   }
243658ebd649SToby Isaac   if (labelname) {
243758ebd649SToby Isaac     PetscValidPointer(labelname, 5);
243858ebd649SToby Isaac     *labelname = b->labelname;
243958ebd649SToby Isaac   }
244058ebd649SToby Isaac   if (field) {
244158ebd649SToby Isaac     PetscValidPointer(field, 6);
244258ebd649SToby Isaac     *field = b->field;
244358ebd649SToby Isaac   }
244458ebd649SToby Isaac   if (numcomps) {
244558ebd649SToby Isaac     PetscValidPointer(numcomps, 7);
244658ebd649SToby Isaac     *numcomps = b->numcomps;
244758ebd649SToby Isaac   }
244858ebd649SToby Isaac   if (comps) {
244958ebd649SToby Isaac     PetscValidPointer(comps, 8);
245058ebd649SToby Isaac     *comps = b->comps;
245158ebd649SToby Isaac   }
245258ebd649SToby Isaac   if (func) {
245358ebd649SToby Isaac     PetscValidPointer(func, 9);
245458ebd649SToby Isaac     *func = b->func;
245558ebd649SToby Isaac   }
245658ebd649SToby Isaac   if (numids) {
245758ebd649SToby Isaac     PetscValidPointer(numids, 10);
245858ebd649SToby Isaac     *numids = b->numids;
245958ebd649SToby Isaac   }
246058ebd649SToby Isaac   if (ids) {
246158ebd649SToby Isaac     PetscValidPointer(ids, 11);
246258ebd649SToby Isaac     *ids = b->ids;
246358ebd649SToby Isaac   }
246458ebd649SToby Isaac   if (ctx) {
246558ebd649SToby Isaac     PetscValidPointer(ctx, 12);
246658ebd649SToby Isaac     *ctx = b->ctx;
246758ebd649SToby Isaac   }
246858ebd649SToby Isaac   PetscFunctionReturn(0);
246958ebd649SToby Isaac }
247058ebd649SToby Isaac 
2471dff059c6SToby Isaac PetscErrorCode PetscDSCopyBoundary(PetscDS probA, PetscDS probB)
2472dff059c6SToby Isaac {
2473dff059c6SToby Isaac   DSBoundary     b, next, *lastnext;
2474dff059c6SToby Isaac   PetscErrorCode ierr;
2475dff059c6SToby Isaac 
2476dff059c6SToby Isaac   PetscFunctionBegin;
2477dff059c6SToby Isaac   PetscValidHeaderSpecific(probA, PETSCDS_CLASSID, 1);
2478dff059c6SToby Isaac   PetscValidHeaderSpecific(probB, PETSCDS_CLASSID, 2);
2479dff059c6SToby Isaac   if (probA == probB) PetscFunctionReturn(0);
2480dff059c6SToby Isaac   next = probB->boundary;
2481dff059c6SToby Isaac   while (next) {
2482dff059c6SToby Isaac     DSBoundary b = next;
2483dff059c6SToby Isaac 
2484dff059c6SToby Isaac     next = b->next;
2485dff059c6SToby Isaac     ierr = PetscFree(b->comps);CHKERRQ(ierr);
2486dff059c6SToby Isaac     ierr = PetscFree(b->ids);CHKERRQ(ierr);
2487dff059c6SToby Isaac     ierr = PetscFree(b->name);CHKERRQ(ierr);
2488dff059c6SToby Isaac     ierr = PetscFree(b->labelname);CHKERRQ(ierr);
2489dff059c6SToby Isaac     ierr = PetscFree(b);CHKERRQ(ierr);
2490dff059c6SToby Isaac   }
2491dff059c6SToby Isaac   lastnext = &(probB->boundary);
2492dff059c6SToby Isaac   for (b = probA->boundary; b; b = b->next) {
2493dff059c6SToby Isaac     DSBoundary bNew;
2494dff059c6SToby Isaac 
2495459726d8SSatish Balay     ierr = PetscNew(&bNew);CHKERRQ(ierr);
2496dff059c6SToby Isaac     bNew->numcomps = b->numcomps;
2497dff059c6SToby Isaac     ierr = PetscMalloc1(bNew->numcomps, &bNew->comps);CHKERRQ(ierr);
2498dff059c6SToby Isaac     ierr = PetscMemcpy(bNew->comps, b->comps, bNew->numcomps*sizeof(PetscInt));CHKERRQ(ierr);
2499dff059c6SToby Isaac     bNew->numids = b->numids;
2500dff059c6SToby Isaac     ierr = PetscMalloc1(bNew->numids, &bNew->ids);CHKERRQ(ierr);
2501dff059c6SToby Isaac     ierr = PetscMemcpy(bNew->ids, b->ids, bNew->numids*sizeof(PetscInt));CHKERRQ(ierr);
2502dff059c6SToby Isaac     ierr = PetscStrallocpy(b->labelname,(char **) &(bNew->labelname));CHKERRQ(ierr);
2503dff059c6SToby Isaac     ierr = PetscStrallocpy(b->name,(char **) &(bNew->name));CHKERRQ(ierr);
2504dff059c6SToby Isaac     bNew->ctx   = b->ctx;
2505f971fd6bSMatthew G. Knepley     bNew->type  = b->type;
2506dff059c6SToby Isaac     bNew->field = b->field;
2507dff059c6SToby Isaac     bNew->func  = b->func;
2508dff059c6SToby Isaac 
2509dff059c6SToby Isaac     *lastnext = bNew;
2510dff059c6SToby Isaac     lastnext = &(bNew->next);
2511dff059c6SToby Isaac   }
2512dff059c6SToby Isaac   PetscFunctionReturn(0);
2513dff059c6SToby Isaac }
2514dff059c6SToby Isaac 
2515da51fcedSMatthew G. Knepley /*@
2516da51fcedSMatthew G. Knepley   PetscDSCopyEquations - Copy all pointwise function pointers to the new problem
2517da51fcedSMatthew G. Knepley 
2518da51fcedSMatthew G. Knepley   Not collective
2519da51fcedSMatthew G. Knepley 
2520da51fcedSMatthew G. Knepley   Input Parameter:
2521da51fcedSMatthew G. Knepley . prob - The PetscDS object
2522da51fcedSMatthew G. Knepley 
2523da51fcedSMatthew G. Knepley   Output Parameter:
2524da51fcedSMatthew G. Knepley . newprob - The PetscDS copy
2525da51fcedSMatthew G. Knepley 
2526da51fcedSMatthew G. Knepley   Level: intermediate
2527da51fcedSMatthew G. Knepley 
2528da51fcedSMatthew G. Knepley .seealso: PetscDSSetResidual(), PetscDSSetJacobian(), PetscDSSetRiemannSolver(), PetscDSSetBdResidual(), PetscDSSetBdJacobian(), PetscDSCreate()
2529da51fcedSMatthew G. Knepley @*/
2530da51fcedSMatthew G. Knepley PetscErrorCode PetscDSCopyEquations(PetscDS prob, PetscDS newprob)
2531da51fcedSMatthew G. Knepley {
2532da51fcedSMatthew G. Knepley   PetscInt       Nf, Ng, f, g;
2533da51fcedSMatthew G. Knepley   PetscErrorCode ierr;
2534da51fcedSMatthew G. Knepley 
2535da51fcedSMatthew G. Knepley   PetscFunctionBegin;
2536da51fcedSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
2537da51fcedSMatthew G. Knepley   PetscValidHeaderSpecific(newprob, PETSCDS_CLASSID, 2);
2538da51fcedSMatthew G. Knepley   ierr = PetscDSGetNumFields(prob, &Nf);CHKERRQ(ierr);
2539da51fcedSMatthew G. Knepley   ierr = PetscDSGetNumFields(newprob, &Ng);CHKERRQ(ierr);
2540da51fcedSMatthew G. Knepley   if (Nf != Ng) SETERRQ2(PetscObjectComm((PetscObject) prob), PETSC_ERR_ARG_SIZ, "Number of fields must match %D != %D", Nf, Ng);CHKERRQ(ierr);
2541da51fcedSMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
2542da51fcedSMatthew G. Knepley     PetscPointFunc   obj;
2543da51fcedSMatthew G. Knepley     PetscPointFunc   f0, f1;
2544da51fcedSMatthew G. Knepley     PetscPointJac    g0, g1, g2, g3;
2545da51fcedSMatthew G. Knepley     PetscBdPointFunc f0Bd, f1Bd;
2546da51fcedSMatthew G. Knepley     PetscBdPointJac  g0Bd, g1Bd, g2Bd, g3Bd;
2547da51fcedSMatthew G. Knepley     PetscRiemannFunc r;
2548da51fcedSMatthew G. Knepley 
2549da51fcedSMatthew G. Knepley     ierr = PetscDSGetObjective(prob, f, &obj);CHKERRQ(ierr);
2550da51fcedSMatthew G. Knepley     ierr = PetscDSGetResidual(prob, f, &f0, &f1);CHKERRQ(ierr);
2551da51fcedSMatthew G. Knepley     ierr = PetscDSGetBdResidual(prob, f, &f0Bd, &f1Bd);CHKERRQ(ierr);
2552da51fcedSMatthew G. Knepley     ierr = PetscDSGetRiemannSolver(prob, f, &r);CHKERRQ(ierr);
2553da51fcedSMatthew G. Knepley     ierr = PetscDSSetObjective(newprob, f, obj);CHKERRQ(ierr);
2554da51fcedSMatthew G. Knepley     ierr = PetscDSSetResidual(newprob, f, f0, f1);CHKERRQ(ierr);
2555da51fcedSMatthew G. Knepley     ierr = PetscDSSetBdResidual(newprob, f, f0Bd, f1Bd);CHKERRQ(ierr);
2556da51fcedSMatthew G. Knepley     ierr = PetscDSSetRiemannSolver(newprob, f, r);CHKERRQ(ierr);
2557da51fcedSMatthew G. Knepley     for (g = 0; g < Nf; ++g) {
2558da51fcedSMatthew G. Knepley       ierr = PetscDSGetJacobian(prob, f, g, &g0, &g1, &g2, &g3);CHKERRQ(ierr);
2559da51fcedSMatthew G. Knepley       ierr = PetscDSGetBdJacobian(prob, f, g, &g0Bd, &g1Bd, &g2Bd, &g3Bd);CHKERRQ(ierr);
2560da51fcedSMatthew G. Knepley       ierr = PetscDSSetJacobian(newprob, f, g, g0, g1, g2, g3);CHKERRQ(ierr);
2561da51fcedSMatthew G. Knepley       ierr = PetscDSSetBdJacobian(newprob, f, g, g0Bd, g1Bd, g2Bd, g3Bd);CHKERRQ(ierr);
2562da51fcedSMatthew G. Knepley     }
2563da51fcedSMatthew G. Knepley   }
2564da51fcedSMatthew G. Knepley   PetscFunctionReturn(0);
2565da51fcedSMatthew G. Knepley }
2566da51fcedSMatthew G. Knepley 
2567bc4ae4beSMatthew G. Knepley static PetscErrorCode PetscDSDestroy_Basic(PetscDS prob)
25682764a2aaSMatthew G. Knepley {
2569931fb3b8SToby Isaac   PetscErrorCode      ierr;
2570931fb3b8SToby Isaac 
25712764a2aaSMatthew G. Knepley   PetscFunctionBegin;
2572931fb3b8SToby Isaac   ierr = PetscFree(prob->data);CHKERRQ(ierr);
25732764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
25742764a2aaSMatthew G. Knepley }
25752764a2aaSMatthew G. Knepley 
2576bc4ae4beSMatthew G. Knepley static PetscErrorCode PetscDSInitialize_Basic(PetscDS prob)
25772764a2aaSMatthew G. Knepley {
25782764a2aaSMatthew G. Knepley   PetscFunctionBegin;
25792764a2aaSMatthew G. Knepley   prob->ops->setfromoptions = NULL;
25802764a2aaSMatthew G. Knepley   prob->ops->setup          = NULL;
25812764a2aaSMatthew G. Knepley   prob->ops->view           = NULL;
25822764a2aaSMatthew G. Knepley   prob->ops->destroy        = PetscDSDestroy_Basic;
25832764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
25842764a2aaSMatthew G. Knepley }
25852764a2aaSMatthew G. Knepley 
25862764a2aaSMatthew G. Knepley /*MC
25872764a2aaSMatthew G. Knepley   PETSCDSBASIC = "basic" - A discrete system with pointwise residual and boundary residual functions
25882764a2aaSMatthew G. Knepley 
25892764a2aaSMatthew G. Knepley   Level: intermediate
25902764a2aaSMatthew G. Knepley 
25912764a2aaSMatthew G. Knepley .seealso: PetscDSType, PetscDSCreate(), PetscDSSetType()
25922764a2aaSMatthew G. Knepley M*/
25932764a2aaSMatthew G. Knepley 
25942764a2aaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDSCreate_Basic(PetscDS prob)
25952764a2aaSMatthew G. Knepley {
25962764a2aaSMatthew G. Knepley   PetscDS_Basic *b;
25972764a2aaSMatthew G. Knepley   PetscErrorCode      ierr;
25982764a2aaSMatthew G. Knepley 
25992764a2aaSMatthew G. Knepley   PetscFunctionBegin;
2600931fb3b8SToby Isaac   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
26012764a2aaSMatthew G. Knepley   ierr       = PetscNewLog(prob, &b);CHKERRQ(ierr);
26022764a2aaSMatthew G. Knepley   prob->data = b;
26032764a2aaSMatthew G. Knepley 
26042764a2aaSMatthew G. Knepley   ierr = PetscDSInitialize_Basic(prob);CHKERRQ(ierr);
26052764a2aaSMatthew G. Knepley   PetscFunctionReturn(0);
26062764a2aaSMatthew G. Knepley }
2607