Actual source code: dtds.c
petsc-3.6.1 2015-08-06
1: #include <petsc/private/petscdsimpl.h> /*I "petscds.h" I*/
3: PetscClassId PETSCDS_CLASSID = 0;
5: PetscFunctionList PetscDSList = NULL;
6: PetscBool PetscDSRegisterAllCalled = PETSC_FALSE;
10: /*@C
11: PetscDSRegister - Adds a new PetscDS implementation
13: Not Collective
15: Input Parameters:
16: + name - The name of a new user-defined creation routine
17: - create_func - The creation routine itself
19: Notes:
20: PetscDSRegister() may be called multiple times to add several user-defined PetscDSs
22: Sample usage:
23: .vb
24: PetscDSRegister("my_ds", MyPetscDSCreate);
25: .ve
27: Then, your PetscDS type can be chosen with the procedural interface via
28: .vb
29: PetscDSCreate(MPI_Comm, PetscDS *);
30: PetscDSSetType(PetscDS, "my_ds");
31: .ve
32: or at runtime via the option
33: .vb
34: -petscds_type my_ds
35: .ve
37: Level: advanced
39: .keywords: PetscDS, register
40: .seealso: PetscDSRegisterAll(), PetscDSRegisterDestroy()
42: @*/
43: PetscErrorCode PetscDSRegister(const char sname[], PetscErrorCode (*function)(PetscDS))
44: {
48: PetscFunctionListAdd(&PetscDSList, sname, function);
49: return(0);
50: }
54: /*@C
55: PetscDSSetType - Builds a particular PetscDS
57: Collective on PetscDS
59: Input Parameters:
60: + prob - The PetscDS object
61: - name - The kind of system
63: Options Database Key:
64: . -petscds_type <type> - Sets the PetscDS type; use -help for a list of available types
66: Level: intermediate
68: .keywords: PetscDS, set, type
69: .seealso: PetscDSGetType(), PetscDSCreate()
70: @*/
71: PetscErrorCode PetscDSSetType(PetscDS prob, PetscDSType name)
72: {
73: PetscErrorCode (*r)(PetscDS);
74: PetscBool match;
79: PetscObjectTypeCompare((PetscObject) prob, name, &match);
80: if (match) return(0);
82: PetscDSRegisterAll();
83: PetscFunctionListFind(PetscDSList, name, &r);
84: if (!r) SETERRQ1(PetscObjectComm((PetscObject) prob), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown PetscDS type: %s", name);
86: if (prob->ops->destroy) {
87: (*prob->ops->destroy)(prob);
88: prob->ops->destroy = NULL;
89: }
90: (*r)(prob);
91: PetscObjectChangeTypeName((PetscObject) prob, name);
92: return(0);
93: }
97: /*@C
98: PetscDSGetType - Gets the PetscDS type name (as a string) from the object.
100: Not Collective
102: Input Parameter:
103: . prob - The PetscDS
105: Output Parameter:
106: . name - The PetscDS type name
108: Level: intermediate
110: .keywords: PetscDS, get, type, name
111: .seealso: PetscDSSetType(), PetscDSCreate()
112: @*/
113: PetscErrorCode PetscDSGetType(PetscDS prob, PetscDSType *name)
114: {
120: PetscDSRegisterAll();
121: *name = ((PetscObject) prob)->type_name;
122: return(0);
123: }
127: static PetscErrorCode PetscDSView_Ascii(PetscDS prob, PetscViewer viewer)
128: {
129: PetscViewerFormat format;
130: PetscInt f;
131: PetscErrorCode ierr;
134: PetscViewerGetFormat(viewer, &format);
135: PetscViewerASCIIPrintf(viewer, "Discrete System with %d fields\n", prob->Nf);
136: PetscViewerASCIIPushTab(viewer);
137: for (f = 0; f < prob->Nf; ++f) {
138: PetscObject obj;
139: PetscClassId id;
140: const char *name;
141: PetscInt Nc;
143: PetscDSGetDiscretization(prob, f, &obj);
144: PetscObjectGetClassId(obj, &id);
145: PetscObjectGetName(obj, &name);
146: PetscViewerASCIIPrintf(viewer, "Field %s", name ? name : "<unknown>");
147: if (id == PETSCFE_CLASSID) {
148: PetscFEGetNumComponents((PetscFE) obj, &Nc);
149: PetscViewerASCIIPrintf(viewer, " FEM");
150: } else if (id == PETSCFV_CLASSID) {
151: PetscFVGetNumComponents((PetscFV) obj, &Nc);
152: PetscViewerASCIIPrintf(viewer, " FVM");
153: }
154: else SETERRQ1(PetscObjectComm((PetscObject) prob), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", f);
155: if (Nc > 1) {PetscViewerASCIIPrintf(viewer, "%d components", Nc);}
156: else {PetscViewerASCIIPrintf(viewer, "%d component ", Nc);}
157: if (prob->implicit[f]) {PetscViewerASCIIPrintf(viewer, " (implicit)");}
158: else {PetscViewerASCIIPrintf(viewer, " (explicit)");}
159: if (prob->adjacency[f*2+0]) {
160: if (prob->adjacency[f*2+1]) {PetscViewerASCIIPrintf(viewer, " (adj FVM++)");}
161: else {PetscViewerASCIIPrintf(viewer, " (adj FVM)");}
162: } else {
163: if (prob->adjacency[f*2+1]) {PetscViewerASCIIPrintf(viewer, " (adj FEM)");}
164: else {PetscViewerASCIIPrintf(viewer, " (adj FUNKY)");}
165: }
166: PetscViewerASCIIPrintf(viewer, "\n");
167: if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
168: if (id == PETSCFE_CLASSID) {PetscFEView((PetscFE) obj, viewer);}
169: else if (id == PETSCFV_CLASSID) {PetscFVView((PetscFV) obj, viewer);}
170: }
171: }
172: PetscViewerASCIIPopTab(viewer);
173: return(0);
174: }
178: /*@C
179: PetscDSView - Views a PetscDS
181: Collective on PetscDS
183: Input Parameter:
184: + prob - the PetscDS object to view
185: - v - the viewer
187: Level: developer
189: .seealso PetscDSDestroy()
190: @*/
191: PetscErrorCode PetscDSView(PetscDS prob, PetscViewer v)
192: {
193: PetscBool iascii;
198: if (!v) {PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject) prob), &v);}
200: PetscObjectTypeCompare((PetscObject) v, PETSCVIEWERASCII, &iascii);
201: if (iascii) {PetscDSView_Ascii(prob, v);}
202: if (prob->ops->view) {(*prob->ops->view)(prob, v);}
203: return(0);
204: }
208: /*@
209: PetscDSSetFromOptions - sets parameters in a PetscDS from the options database
211: Collective on PetscDS
213: Input Parameter:
214: . prob - the PetscDS object to set options for
216: Options Database:
218: Level: developer
220: .seealso PetscDSView()
221: @*/
222: PetscErrorCode PetscDSSetFromOptions(PetscDS prob)
223: {
224: const char *defaultType;
225: char name[256];
226: PetscBool flg;
231: if (!((PetscObject) prob)->type_name) {
232: defaultType = PETSCDSBASIC;
233: } else {
234: defaultType = ((PetscObject) prob)->type_name;
235: }
236: PetscDSRegisterAll();
238: PetscObjectOptionsBegin((PetscObject) prob);
239: PetscOptionsFList("-petscds_type", "Discrete System", "PetscDSSetType", PetscDSList, defaultType, name, 256, &flg);
240: if (flg) {
241: PetscDSSetType(prob, name);
242: } else if (!((PetscObject) prob)->type_name) {
243: PetscDSSetType(prob, defaultType);
244: }
245: if (prob->ops->setfromoptions) {(*prob->ops->setfromoptions)(prob);}
246: /* process any options handlers added with PetscObjectAddOptionsHandler() */
247: PetscObjectProcessOptionsHandlers((PetscObject) prob);
248: PetscOptionsEnd();
249: PetscDSViewFromOptions(prob, NULL, "-petscds_view");
250: return(0);
251: }
255: /*@C
256: PetscDSSetUp - Construct data structures for the PetscDS
258: Collective on PetscDS
260: Input Parameter:
261: . prob - the PetscDS object to setup
263: Level: developer
265: .seealso PetscDSView(), PetscDSDestroy()
266: @*/
267: PetscErrorCode PetscDSSetUp(PetscDS prob)
268: {
269: const PetscInt Nf = prob->Nf;
270: PetscInt dim, work, NcMax = 0, NqMax = 0, f;
275: if (prob->setup) return(0);
276: /* Calculate sizes */
277: PetscDSGetSpatialDimension(prob, &dim);
278: prob->totDim = prob->totDimBd = prob->totComp = 0;
279: PetscCalloc4(Nf+1,&prob->off,Nf+1,&prob->offDer,Nf+1,&prob->offBd,Nf+1,&prob->offDerBd);
280: PetscMalloc4(Nf,&prob->basis,Nf,&prob->basisDer,Nf,&prob->basisBd,Nf,&prob->basisDerBd);
281: for (f = 0; f < Nf; ++f) {
282: PetscFE feBd = (PetscFE) prob->discBd[f];
283: PetscObject obj;
284: PetscClassId id;
285: PetscQuadrature q;
286: PetscInt Nq = 0, Nb, Nc;
288: PetscDSGetDiscretization(prob, f, &obj);
289: PetscObjectGetClassId(obj, &id);
290: if (id == PETSCFE_CLASSID) {
291: PetscFE fe = (PetscFE) obj;
293: PetscFEGetQuadrature(fe, &q);
294: PetscFEGetDimension(fe, &Nb);
295: PetscFEGetNumComponents(fe, &Nc);
296: PetscFEGetDefaultTabulation(fe, &prob->basis[f], &prob->basisDer[f], NULL);
297: } else if (id == PETSCFV_CLASSID) {
298: PetscFV fv = (PetscFV) obj;
300: PetscFVGetQuadrature(fv, &q);
301: Nb = 1;
302: PetscFVGetNumComponents(fv, &Nc);
303: PetscFVGetDefaultTabulation(fv, &prob->basis[f], &prob->basisDer[f], NULL);
304: } else SETERRQ1(PetscObjectComm((PetscObject) prob), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", f);
305: prob->off[f+1] = Nc + prob->off[f];
306: prob->offDer[f+1] = Nc*dim + prob->offDer[f];
307: if (q) {PetscQuadratureGetData(q, NULL, &Nq, NULL, NULL);}
308: NqMax = PetscMax(NqMax, Nq);
309: NcMax = PetscMax(NcMax, Nc);
310: prob->totDim += Nb*Nc;
311: prob->totComp += Nc;
312: if (feBd) {
313: PetscFEGetDimension(feBd, &Nb);
314: PetscFEGetNumComponents(feBd, &Nc);
315: PetscFEGetDefaultTabulation(feBd, &prob->basisBd[f], &prob->basisDerBd[f], NULL);
316: prob->totDimBd += Nb*Nc;
317: prob->offBd[f+1] = Nc + prob->offBd[f];
318: prob->offDerBd[f+1] = Nc*dim + prob->offDerBd[f];
319: }
320: }
321: work = PetscMax(prob->totComp*dim, PetscSqr(NcMax*dim));
322: /* Allocate works space */
323: PetscMalloc5(prob->totComp,&prob->u,prob->totComp,&prob->u_t,prob->totComp*dim,&prob->u_x,dim,&prob->x,work,&prob->refSpaceDer);
324: 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);
325: if (prob->ops->setup) {(*prob->ops->setup)(prob);}
326: prob->setup = PETSC_TRUE;
327: return(0);
328: }
332: static PetscErrorCode PetscDSDestroyStructs_Static(PetscDS prob)
333: {
337: PetscFree4(prob->off,prob->offDer,prob->offBd,prob->offDerBd);
338: PetscFree4(prob->basis,prob->basisDer,prob->basisBd,prob->basisDerBd);
339: PetscFree5(prob->u,prob->u_t,prob->u_x,prob->x,prob->refSpaceDer);
340: PetscFree6(prob->f0,prob->f1,prob->g0,prob->g1,prob->g2,prob->g3);
341: return(0);
342: }
346: static PetscErrorCode PetscDSEnlarge_Static(PetscDS prob, PetscInt NfNew)
347: {
348: PetscObject *tmpd, *tmpdbd;
349: PetscBool *tmpi, *tmpa;
350: PetscPointFunc *tmpobj, *tmpf;
351: PetscPointJac *tmpg;
352: PetscBdPointFunc *tmpfbd;
353: PetscBdPointJac *tmpgbd;
354: PetscRiemannFunc *tmpr;
355: void **tmpctx;
356: PetscInt Nf = prob->Nf, f, i;
357: PetscErrorCode ierr;
360: if (Nf >= NfNew) return(0);
361: prob->setup = PETSC_FALSE;
362: PetscDSDestroyStructs_Static(prob);
363: PetscMalloc4(NfNew, &tmpd, NfNew, &tmpdbd, NfNew, &tmpi, NfNew*2, &tmpa);
364: for (f = 0; f < Nf; ++f) {tmpd[f] = prob->disc[f]; tmpdbd[f] = prob->discBd[f]; tmpi[f] = prob->implicit[f]; for (i = 0; i < 2; ++i) tmpa[f*2+i] = prob->adjacency[f*2+i];}
365: for (f = Nf; f < NfNew; ++f) {PetscContainerCreate(PetscObjectComm((PetscObject) prob), (PetscContainer *) &tmpd[f]);
366: tmpdbd[f] = NULL; tmpi[f] = PETSC_TRUE; tmpa[f*2+0] = PETSC_FALSE; tmpa[f*2+1] = PETSC_TRUE;}
367: PetscFree4(prob->disc, prob->discBd, prob->implicit, prob->adjacency);
368: prob->Nf = NfNew;
369: prob->disc = tmpd;
370: prob->discBd = tmpdbd;
371: prob->implicit = tmpi;
372: prob->adjacency = tmpa;
373: PetscCalloc5(NfNew, &tmpobj, NfNew*2, &tmpf, NfNew*NfNew*4, &tmpg, NfNew, &tmpr, NfNew, &tmpctx);
374: for (f = 0; f < Nf; ++f) tmpobj[f] = prob->obj[f];
375: for (f = 0; f < Nf*2; ++f) tmpf[f] = prob->f[f];
376: for (f = 0; f < Nf*Nf*4; ++f) tmpg[f] = prob->g[f];
377: for (f = 0; f < Nf; ++f) tmpr[f] = prob->r[f];
378: for (f = 0; f < Nf; ++f) tmpctx[f] = prob->ctx[f];
379: for (f = Nf; f < NfNew; ++f) tmpobj[f] = NULL;
380: for (f = Nf*2; f < NfNew*2; ++f) tmpf[f] = NULL;
381: for (f = Nf*Nf*4; f < NfNew*NfNew*4; ++f) tmpg[f] = NULL;
382: for (f = Nf; f < NfNew; ++f) tmpr[f] = NULL;
383: for (f = Nf; f < NfNew; ++f) tmpctx[f] = NULL;
384: PetscFree5(prob->obj, prob->f, prob->g, prob->r, prob->ctx);
385: prob->obj = tmpobj;
386: prob->f = tmpf;
387: prob->g = tmpg;
388: prob->r = tmpr;
389: prob->ctx = tmpctx;
390: PetscCalloc2(NfNew*2, &tmpfbd, NfNew*NfNew*4, &tmpgbd);
391: for (f = 0; f < Nf*2; ++f) tmpfbd[f] = prob->fBd[f];
392: for (f = 0; f < Nf*Nf*4; ++f) tmpgbd[f] = prob->gBd[f];
393: for (f = Nf*2; f < NfNew*2; ++f) tmpfbd[f] = NULL;
394: for (f = Nf*Nf*4; f < NfNew*NfNew*4; ++f) tmpgbd[f] = NULL;
395: PetscFree2(prob->fBd, prob->gBd);
396: prob->fBd = tmpfbd;
397: prob->gBd = tmpgbd;
398: return(0);
399: }
403: /*@
404: PetscDSDestroy - Destroys a PetscDS object
406: Collective on PetscDS
408: Input Parameter:
409: . prob - the PetscDS object to destroy
411: Level: developer
413: .seealso PetscDSView()
414: @*/
415: PetscErrorCode PetscDSDestroy(PetscDS *prob)
416: {
417: PetscInt f;
421: if (!*prob) return(0);
424: if (--((PetscObject)(*prob))->refct > 0) {*prob = 0; return(0);}
425: ((PetscObject) (*prob))->refct = 0;
426: PetscDSDestroyStructs_Static(*prob);
427: for (f = 0; f < (*prob)->Nf; ++f) {
428: PetscObjectDereference((*prob)->disc[f]);
429: PetscObjectDereference((*prob)->discBd[f]);
430: }
431: PetscFree4((*prob)->disc, (*prob)->discBd, (*prob)->implicit, (*prob)->adjacency);
432: PetscFree5((*prob)->obj,(*prob)->f,(*prob)->g,(*prob)->r,(*prob)->ctx);
433: PetscFree2((*prob)->fBd,(*prob)->gBd);
434: if ((*prob)->ops->destroy) {(*(*prob)->ops->destroy)(*prob);}
435: PetscHeaderDestroy(prob);
436: return(0);
437: }
441: /*@
442: PetscDSCreate - Creates an empty PetscDS object. The type can then be set with PetscDSSetType().
444: Collective on MPI_Comm
446: Input Parameter:
447: . comm - The communicator for the PetscDS object
449: Output Parameter:
450: . prob - The PetscDS object
452: Level: beginner
454: .seealso: PetscDSSetType(), PETSCDSBASIC
455: @*/
456: PetscErrorCode PetscDSCreate(MPI_Comm comm, PetscDS *prob)
457: {
458: PetscDS p;
463: *prob = NULL;
464: PetscDSInitializePackage();
466: PetscHeaderCreate(p, PETSCDS_CLASSID, "PetscDS", "Discrete System", "PetscDS", comm, PetscDSDestroy, PetscDSView);
468: p->Nf = 0;
469: p->setup = PETSC_FALSE;
471: *prob = p;
472: return(0);
473: }
477: /*@
478: PetscDSGetNumFields - Returns the number of fields in the DS
480: Not collective
482: Input Parameter:
483: . prob - The PetscDS object
485: Output Parameter:
486: . Nf - The number of fields
488: Level: beginner
490: .seealso: PetscDSGetSpatialDimension(), PetscDSCreate()
491: @*/
492: PetscErrorCode PetscDSGetNumFields(PetscDS prob, PetscInt *Nf)
493: {
497: *Nf = prob->Nf;
498: return(0);
499: }
503: /*@
504: PetscDSGetSpatialDimension - Returns the spatial dimension of the DS
506: Not collective
508: Input Parameter:
509: . prob - The PetscDS object
511: Output Parameter:
512: . dim - The spatial dimension
514: Level: beginner
516: .seealso: PetscDSGetNumFields(), PetscDSCreate()
517: @*/
518: PetscErrorCode PetscDSGetSpatialDimension(PetscDS prob, PetscInt *dim)
519: {
525: *dim = 0;
526: if (prob->Nf) {
527: PetscObject obj;
528: PetscClassId id;
530: PetscDSGetDiscretization(prob, 0, &obj);
531: PetscObjectGetClassId(obj, &id);
532: if (id == PETSCFE_CLASSID) {PetscFEGetSpatialDimension((PetscFE) obj, dim);}
533: else if (id == PETSCFV_CLASSID) {PetscFVGetSpatialDimension((PetscFV) obj, dim);}
534: else SETERRQ1(PetscObjectComm((PetscObject) prob), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", 0);
535: }
536: return(0);
537: }
541: /*@
542: PetscDSGetTotalDimension - Returns the total size of the approximation space for this system
544: Not collective
546: Input Parameter:
547: . prob - The PetscDS object
549: Output Parameter:
550: . dim - The total problem dimension
552: Level: beginner
554: .seealso: PetscDSGetNumFields(), PetscDSCreate()
555: @*/
556: PetscErrorCode PetscDSGetTotalDimension(PetscDS prob, PetscInt *dim)
557: {
562: PetscDSSetUp(prob);
564: *dim = prob->totDim;
565: return(0);
566: }
570: /*@
571: PetscDSGetTotalBdDimension - Returns the total size of the boundary approximation space for this system
573: Not collective
575: Input Parameter:
576: . prob - The PetscDS object
578: Output Parameter:
579: . dim - The total boundary problem dimension
581: Level: beginner
583: .seealso: PetscDSGetNumFields(), PetscDSCreate()
584: @*/
585: PetscErrorCode PetscDSGetTotalBdDimension(PetscDS prob, PetscInt *dim)
586: {
591: PetscDSSetUp(prob);
593: *dim = prob->totDimBd;
594: return(0);
595: }
599: /*@
600: PetscDSGetTotalComponents - Returns the total number of components in this system
602: Not collective
604: Input Parameter:
605: . prob - The PetscDS object
607: Output Parameter:
608: . dim - The total number of components
610: Level: beginner
612: .seealso: PetscDSGetNumFields(), PetscDSCreate()
613: @*/
614: PetscErrorCode PetscDSGetTotalComponents(PetscDS prob, PetscInt *Nc)
615: {
620: PetscDSSetUp(prob);
622: *Nc = prob->totComp;
623: return(0);
624: }
628: /*@
629: PetscDSGetDiscretization - Returns the discretization object for the given field
631: Not collective
633: Input Parameters:
634: + prob - The PetscDS object
635: - f - The field number
637: Output Parameter:
638: . disc - The discretization object
640: Level: beginner
642: .seealso: PetscDSSetDiscretization(), PetscDSAddDiscretization(), PetscDSGetBdDiscretization(), PetscDSGetNumFields(), PetscDSCreate()
643: @*/
644: PetscErrorCode PetscDSGetDiscretization(PetscDS prob, PetscInt f, PetscObject *disc)
645: {
649: 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);
650: *disc = prob->disc[f];
651: return(0);
652: }
656: /*@
657: PetscDSGetBdDiscretization - Returns the boundary discretization object for the given field
659: Not collective
661: Input Parameters:
662: + prob - The PetscDS object
663: - f - The field number
665: Output Parameter:
666: . disc - The boundary discretization object
668: Level: beginner
670: .seealso: PetscDSSetBdDiscretization(), PetscDSAddBdDiscretization(), PetscDSGetDiscretization(), PetscDSGetNumFields(), PetscDSCreate()
671: @*/
672: PetscErrorCode PetscDSGetBdDiscretization(PetscDS prob, PetscInt f, PetscObject *disc)
673: {
677: 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);
678: *disc = prob->discBd[f];
679: return(0);
680: }
684: /*@
685: PetscDSSetDiscretization - Sets the discretization object for the given field
687: Not collective
689: Input Parameters:
690: + prob - The PetscDS object
691: . f - The field number
692: - disc - The discretization object
694: Level: beginner
696: .seealso: PetscDSGetDiscretization(), PetscDSAddDiscretization(), PetscDSGetNumFields(), PetscDSCreate()
697: @*/
698: PetscErrorCode PetscDSSetDiscretization(PetscDS prob, PetscInt f, PetscObject disc)
699: {
705: if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f);
706: PetscDSEnlarge_Static(prob, f+1);
707: if (prob->disc[f]) {PetscObjectDereference(prob->disc[f]);}
708: prob->disc[f] = disc;
709: PetscObjectReference(disc);
710: {
711: PetscClassId id;
713: PetscObjectGetClassId(disc, &id);
714: if (id == PETSCFV_CLASSID) {
715: prob->implicit[f] = PETSC_FALSE;
716: prob->adjacency[f*2+0] = PETSC_TRUE;
717: prob->adjacency[f*2+1] = PETSC_FALSE;
718: }
719: }
720: return(0);
721: }
725: /*@
726: PetscDSSetBdDiscretization - Sets the boundary discretization object for the given field
728: Not collective
730: Input Parameters:
731: + prob - The PetscDS object
732: . f - The field number
733: - disc - The boundary discretization object
735: Level: beginner
737: .seealso: PetscDSGetBdDiscretization(), PetscDSAddBdDiscretization(), PetscDSGetNumFields(), PetscDSCreate()
738: @*/
739: PetscErrorCode PetscDSSetBdDiscretization(PetscDS prob, PetscInt f, PetscObject disc)
740: {
746: if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f);
747: PetscDSEnlarge_Static(prob, f+1);
748: if (prob->discBd[f]) {PetscObjectDereference(prob->discBd[f]);}
749: prob->discBd[f] = disc;
750: PetscObjectReference(disc);
751: return(0);
752: }
756: /*@
757: PetscDSAddDiscretization - Adds a discretization object
759: Not collective
761: Input Parameters:
762: + prob - The PetscDS object
763: - disc - The boundary discretization object
765: Level: beginner
767: .seealso: PetscDSGetDiscretization(), PetscDSSetDiscretization(), PetscDSGetNumFields(), PetscDSCreate()
768: @*/
769: PetscErrorCode PetscDSAddDiscretization(PetscDS prob, PetscObject disc)
770: {
774: PetscDSSetDiscretization(prob, prob->Nf, disc);
775: return(0);
776: }
780: /*@
781: PetscDSAddBdDiscretization - Adds a boundary discretization object
783: Not collective
785: Input Parameters:
786: + prob - The PetscDS object
787: - disc - The boundary discretization object
789: Level: beginner
791: .seealso: PetscDSGetBdDiscretization(), PetscDSSetBdDiscretization(), PetscDSGetNumFields(), PetscDSCreate()
792: @*/
793: PetscErrorCode PetscDSAddBdDiscretization(PetscDS prob, PetscObject disc)
794: {
798: PetscDSSetBdDiscretization(prob, prob->Nf, disc);
799: return(0);
800: }
804: /*@
805: PetscDSGetImplicit - Returns the flag for implicit solve for this field. This is just a guide for IMEX
807: Not collective
809: Input Parameters:
810: + prob - The PetscDS object
811: - f - The field number
813: Output Parameter:
814: . implicit - The flag indicating what kind of solve to use for this field
816: Level: developer
818: .seealso: PetscDSSetImplicit(), PetscDSSetDiscretization(), PetscDSAddDiscretization(), PetscDSGetBdDiscretization(), PetscDSGetNumFields(), PetscDSCreate()
819: @*/
820: PetscErrorCode PetscDSGetImplicit(PetscDS prob, PetscInt f, PetscBool *implicit)
821: {
825: 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);
826: *implicit = prob->implicit[f];
827: return(0);
828: }
832: /*@
833: PetscDSSetImplicit - Set the flag for implicit solve for this field. This is just a guide for IMEX
835: Not collective
837: Input Parameters:
838: + prob - The PetscDS object
839: . f - The field number
840: - implicit - The flag indicating what kind of solve to use for this field
842: Level: developer
844: .seealso: PetscDSGetImplicit(), PetscDSSetDiscretization(), PetscDSAddDiscretization(), PetscDSGetBdDiscretization(), PetscDSGetNumFields(), PetscDSCreate()
845: @*/
846: PetscErrorCode PetscDSSetImplicit(PetscDS prob, PetscInt f, PetscBool implicit)
847: {
850: 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);
851: prob->implicit[f] = implicit;
852: return(0);
853: }
857: /*@
858: PetscDSGetAdjacency - Returns the flags for determining variable influence
860: Not collective
862: Input Parameters:
863: + prob - The PetscDS object
864: - f - The field number
866: Output Parameter:
867: + useCone - Flag for variable influence starting with the cone operation
868: - useClosure - Flag for variable influence using transitive closure
870: Note: See the discussion in DMPlexGetAdjacencyUseCone() and DMPlexGetAdjacencyUseClosure()
872: Level: developer
874: .seealso: PetscDSSetAdjacency(), DMPlexGetAdjacencyUseCone(), DMPlexGetAdjacencyUseClosure(), PetscDSSetDiscretization(), PetscDSAddDiscretization(), PetscDSGetBdDiscretization(), PetscDSGetNumFields(), PetscDSCreate()
875: @*/
876: PetscErrorCode PetscDSGetAdjacency(PetscDS prob, PetscInt f, PetscBool *useCone, PetscBool *useClosure)
877: {
882: 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);
883: *useCone = prob->adjacency[f*2+0];
884: *useClosure = prob->adjacency[f*2+1];
885: return(0);
886: }
890: /*@
891: PetscDSSetAdjacency - Set the flags for determining variable influence
893: Not collective
895: Input Parameters:
896: + prob - The PetscDS object
897: . f - The field number
898: . useCone - Flag for variable influence starting with the cone operation
899: - useClosure - Flag for variable influence using transitive closure
901: Note: See the discussion in DMPlexGetAdjacencyUseCone() and DMPlexGetAdjacencyUseClosure()
903: Level: developer
905: .seealso: PetscDSGetAdjacency(), DMPlexGetAdjacencyUseCone(), DMPlexGetAdjacencyUseClosure(), PetscDSSetDiscretization(), PetscDSAddDiscretization(), PetscDSGetBdDiscretization(), PetscDSGetNumFields(), PetscDSCreate()
906: @*/
907: PetscErrorCode PetscDSSetAdjacency(PetscDS prob, PetscInt f, PetscBool useCone, PetscBool useClosure)
908: {
911: 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);
912: prob->adjacency[f*2+0] = useCone;
913: prob->adjacency[f*2+1] = useClosure;
914: return(0);
915: }
919: PetscErrorCode PetscDSGetObjective(PetscDS prob, PetscInt f,
920: void (**obj)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
921: const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
922: const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
923: PetscReal t, const PetscReal x[], PetscScalar obj[]))
924: {
928: 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);
929: *obj = prob->obj[f];
930: return(0);
931: }
935: PetscErrorCode PetscDSSetObjective(PetscDS prob, PetscInt f,
936: void (*obj)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
937: const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
938: const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
939: PetscReal t, const PetscReal x[], PetscScalar obj[]))
940: {
946: if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f);
947: PetscDSEnlarge_Static(prob, f+1);
948: prob->obj[f] = obj;
949: return(0);
950: }
954: /*@C
955: PetscDSGetResidual - Get the pointwise residual function for a given test field
957: Not collective
959: Input Parameters:
960: + prob - The PetscDS
961: - f - The test field number
963: Output Parameters:
964: + f0 - integrand for the test function term
965: - f1 - integrand for the test function gradient term
967: Note: We are using a first order FEM model for the weak form:
969: \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)
971: The calling sequence for the callbacks f0 and f1 is given by:
973: $ f0(PetscInt dim, PetscInt Nf, PetscInt NfAux,
974: $ const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
975: $ const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
976: $ PetscReal t, const PetscReal x[], PetscScalar f0[])
978: + dim - the spatial dimension
979: . Nf - the number of fields
980: . uOff - the offset into u[] and u_t[] for each field
981: . uOff_x - the offset into u_x[] for each field
982: . u - each field evaluated at the current point
983: . u_t - the time derivative of each field evaluated at the current point
984: . u_x - the gradient of each field evaluated at the current point
985: . aOff - the offset into a[] and a_t[] for each auxiliary field
986: . aOff_x - the offset into a_x[] for each auxiliary field
987: . a - each auxiliary field evaluated at the current point
988: . a_t - the time derivative of each auxiliary field evaluated at the current point
989: . a_x - the gradient of auxiliary each field evaluated at the current point
990: . t - current time
991: . x - coordinates of the current point
992: - f0 - output values at the current point
994: Level: intermediate
996: .seealso: PetscDSSetResidual()
997: @*/
998: PetscErrorCode PetscDSGetResidual(PetscDS prob, PetscInt f,
999: void (**f0)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1000: const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1001: const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1002: PetscReal t, const PetscReal x[], PetscScalar f0[]),
1003: void (**f1)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1004: const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1005: const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1006: PetscReal t, const PetscReal x[], PetscScalar f1[]))
1007: {
1010: 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);
1013: return(0);
1014: }
1018: /*@C
1019: PetscDSSetResidual - Set the pointwise residual function for a given test field
1021: Not collective
1023: Input Parameters:
1024: + prob - The PetscDS
1025: . f - The test field number
1026: . f0 - integrand for the test function term
1027: - f1 - integrand for the test function gradient term
1029: Note: We are using a first order FEM model for the weak form:
1031: \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)
1033: The calling sequence for the callbacks f0 and f1 is given by:
1035: $ f0(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1036: $ const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1037: $ const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1038: $ PetscReal t, const PetscReal x[], PetscScalar f0[])
1040: + dim - the spatial dimension
1041: . Nf - the number of fields
1042: . uOff - the offset into u[] and u_t[] for each field
1043: . uOff_x - the offset into u_x[] for each field
1044: . u - each field evaluated at the current point
1045: . u_t - the time derivative of each field evaluated at the current point
1046: . u_x - the gradient of each field evaluated at the current point
1047: . aOff - the offset into a[] and a_t[] for each auxiliary field
1048: . aOff_x - the offset into a_x[] for each auxiliary field
1049: . a - each auxiliary field evaluated at the current point
1050: . a_t - the time derivative of each auxiliary field evaluated at the current point
1051: . a_x - the gradient of auxiliary each field evaluated at the current point
1052: . t - current time
1053: . x - coordinates of the current point
1054: - f0 - output values at the current point
1056: Level: intermediate
1058: .seealso: PetscDSGetResidual()
1059: @*/
1060: PetscErrorCode PetscDSSetResidual(PetscDS prob, PetscInt f,
1061: void (*f0)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1062: const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1063: const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1064: PetscReal t, const PetscReal x[], PetscScalar f0[]),
1065: void (*f1)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1066: const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1067: const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1068: PetscReal t, const PetscReal x[], PetscScalar f1[]))
1069: {
1076: if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f);
1077: PetscDSEnlarge_Static(prob, f+1);
1078: prob->f[f*2+0] = f0;
1079: prob->f[f*2+1] = f1;
1080: return(0);
1081: }
1085: /*@C
1086: PetscDSGetJacobian - Get the pointwise Jacobian function for given test and basis field
1088: Not collective
1090: Input Parameters:
1091: + prob - The PetscDS
1092: . f - The test field number
1093: - g - The field number
1095: Output Parameters:
1096: + g0 - integrand for the test and basis function term
1097: . g1 - integrand for the test function and basis function gradient term
1098: . g2 - integrand for the test function gradient and basis function term
1099: - g3 - integrand for the test function gradient and basis function gradient term
1101: Note: We are using a first order FEM model for the weak form:
1103: \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
1105: The calling sequence for the callbacks g0, g1, g2 and g3 is given by:
1107: $ g0(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1108: $ const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1109: $ const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1110: $ PetscReal t, const PetscReal u_tShift, const PetscReal x[], PetscScalar g0[])
1112: + dim - the spatial dimension
1113: . Nf - the number of fields
1114: . uOff - the offset into u[] and u_t[] for each field
1115: . uOff_x - the offset into u_x[] for each field
1116: . u - each field evaluated at the current point
1117: . u_t - the time derivative of each field evaluated at the current point
1118: . u_x - the gradient of each field evaluated at the current point
1119: . aOff - the offset into a[] and a_t[] for each auxiliary field
1120: . aOff_x - the offset into a_x[] for each auxiliary field
1121: . a - each auxiliary field evaluated at the current point
1122: . a_t - the time derivative of each auxiliary field evaluated at the current point
1123: . a_x - the gradient of auxiliary each field evaluated at the current point
1124: . t - current time
1125: . u_tShift - the multiplier a for dF/dU_t
1126: . x - coordinates of the current point
1127: - g0 - output values at the current point
1129: Level: intermediate
1131: .seealso: PetscDSSetJacobian()
1132: @*/
1133: PetscErrorCode PetscDSGetJacobian(PetscDS prob, PetscInt f, PetscInt g,
1134: void (**g0)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1135: const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1136: const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1137: PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscScalar g0[]),
1138: void (**g1)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1139: const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1140: const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1141: PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscScalar g1[]),
1142: void (**g2)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1143: const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1144: const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1145: PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscScalar g2[]),
1146: void (**g3)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1147: const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1148: const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1149: PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscScalar g3[]))
1150: {
1153: 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);
1154: 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);
1159: return(0);
1160: }
1164: /*@C
1165: PetscDSSetJacobian - Set the pointwise Jacobian function for given test and basis fields
1167: Not collective
1169: Input Parameters:
1170: + prob - The PetscDS
1171: . f - The test field number
1172: . g - The field number
1173: . g0 - integrand for the test and basis function term
1174: . g1 - integrand for the test function and basis function gradient term
1175: . g2 - integrand for the test function gradient and basis function term
1176: - g3 - integrand for the test function gradient and basis function gradient term
1178: Note: We are using a first order FEM model for the weak form:
1180: \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
1182: The calling sequence for the callbacks g0, g1, g2 and g3 is given by:
1184: $ g0(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1185: $ const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1186: $ const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1187: $ PetscReal t, const PetscReal x[], PetscScalar g0[])
1189: + dim - the spatial dimension
1190: . Nf - the number of fields
1191: . uOff - the offset into u[] and u_t[] for each field
1192: . uOff_x - the offset into u_x[] for each field
1193: . u - each field evaluated at the current point
1194: . u_t - the time derivative of each field evaluated at the current point
1195: . u_x - the gradient of each field evaluated at the current point
1196: . aOff - the offset into a[] and a_t[] for each auxiliary field
1197: . aOff_x - the offset into a_x[] for each auxiliary field
1198: . a - each auxiliary field evaluated at the current point
1199: . a_t - the time derivative of each auxiliary field evaluated at the current point
1200: . a_x - the gradient of auxiliary each field evaluated at the current point
1201: . t - current time
1202: . u_tShift - the multiplier a for dF/dU_t
1203: . x - coordinates of the current point
1204: - g0 - output values at the current point
1206: Level: intermediate
1208: .seealso: PetscDSGetJacobian()
1209: @*/
1210: PetscErrorCode PetscDSSetJacobian(PetscDS prob, PetscInt f, PetscInt g,
1211: void (*g0)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1212: const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1213: const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1214: PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscScalar g0[]),
1215: void (*g1)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1216: const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1217: const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1218: PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscScalar g1[]),
1219: void (*g2)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1220: const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1221: const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1222: PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscScalar g2[]),
1223: void (*g3)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1224: const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1225: const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1226: PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscScalar g3[]))
1227: {
1236: if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f);
1237: if (g < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", g);
1238: PetscDSEnlarge_Static(prob, PetscMax(f, g)+1);
1239: prob->g[(f*prob->Nf + g)*4+0] = g0;
1240: prob->g[(f*prob->Nf + g)*4+1] = g1;
1241: prob->g[(f*prob->Nf + g)*4+2] = g2;
1242: prob->g[(f*prob->Nf + g)*4+3] = g3;
1243: return(0);
1244: }
1248: /*@C
1249: PetscDSGetRiemannSolver - Returns the Riemann solver for the given field
1251: Not collective
1253: Input Arguments:
1254: + prob - The PetscDS object
1255: - f - The field number
1257: Output Argument:
1258: . r - Riemann solver
1260: Calling sequence for r:
1262: $ r(PetscInt dim, PetscInt Nf, const PetscReal x[], const PetscReal n[], const PetscScalar uL[], const PetscScalar uR[], PetscScalar flux[], void *ctx)
1264: + dim - The spatial dimension
1265: . Nf - The number of fields
1266: . x - The coordinates at a point on the interface
1267: . n - The normal vector to the interface
1268: . uL - The state vector to the left of the interface
1269: . uR - The state vector to the right of the interface
1270: . flux - output array of flux through the interface
1271: - ctx - optional user context
1273: Level: intermediate
1275: .seealso: PetscDSSetRiemannSolver()
1276: @*/
1277: PetscErrorCode PetscDSGetRiemannSolver(PetscDS prob, PetscInt f,
1278: void (**r)(PetscInt dim, PetscInt Nf, const PetscReal x[], const PetscReal n[], const PetscScalar uL[], const PetscScalar uR[], PetscScalar flux[], void *ctx))
1279: {
1282: 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);
1284: *r = prob->r[f];
1285: return(0);
1286: }
1290: /*@C
1291: PetscDSSetRiemannSolver - Sets the Riemann solver for the given field
1293: Not collective
1295: Input Arguments:
1296: + prob - The PetscDS object
1297: . f - The field number
1298: - r - Riemann solver
1300: Calling sequence for r:
1302: $ r(PetscInt dim, PetscInt Nf, const PetscReal x[], const PetscReal n[], const PetscScalar uL[], const PetscScalar uR[], PetscScalar flux[], void *ctx)
1304: + dim - The spatial dimension
1305: . Nf - The number of fields
1306: . x - The coordinates at a point on the interface
1307: . n - The normal vector to the interface
1308: . uL - The state vector to the left of the interface
1309: . uR - The state vector to the right of the interface
1310: . flux - output array of flux through the interface
1311: - ctx - optional user context
1313: Level: intermediate
1315: .seealso: PetscDSGetRiemannSolver()
1316: @*/
1317: PetscErrorCode PetscDSSetRiemannSolver(PetscDS prob, PetscInt f,
1318: void (*r)(PetscInt dim, PetscInt Nf, const PetscReal x[], const PetscReal n[], const PetscScalar uL[], const PetscScalar uR[], PetscScalar flux[], void *ctx))
1319: {
1325: if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f);
1326: PetscDSEnlarge_Static(prob, f+1);
1327: prob->r[f] = r;
1328: return(0);
1329: }
1333: PetscErrorCode PetscDSGetContext(PetscDS prob, PetscInt f, void **ctx)
1334: {
1337: 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);
1339: *ctx = prob->ctx[f];
1340: return(0);
1341: }
1345: PetscErrorCode PetscDSSetContext(PetscDS prob, PetscInt f, void *ctx)
1346: {
1351: if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f);
1352: PetscDSEnlarge_Static(prob, f+1);
1353: prob->ctx[f] = ctx;
1354: return(0);
1355: }
1359: /*@C
1360: PetscDSGetBdResidual - Get the pointwise boundary residual function for a given test field
1362: Not collective
1364: Input Parameters:
1365: + prob - The PetscDS
1366: - f - The test field number
1368: Output Parameters:
1369: + f0 - boundary integrand for the test function term
1370: - f1 - boundary integrand for the test function gradient term
1372: Note: We are using a first order FEM model for the weak form:
1374: \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
1376: The calling sequence for the callbacks f0 and f1 is given by:
1378: $ f0(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1379: $ const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1380: $ const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1381: $ PetscReal t, const PetscReal x[], const PetscReal n[], PetscScalar f0[])
1383: + dim - the spatial dimension
1384: . Nf - the number of fields
1385: . uOff - the offset into u[] and u_t[] for each field
1386: . uOff_x - the offset into u_x[] for each field
1387: . u - each field evaluated at the current point
1388: . u_t - the time derivative of each field evaluated at the current point
1389: . u_x - the gradient of each field evaluated at the current point
1390: . aOff - the offset into a[] and a_t[] for each auxiliary field
1391: . aOff_x - the offset into a_x[] for each auxiliary field
1392: . a - each auxiliary field evaluated at the current point
1393: . a_t - the time derivative of each auxiliary field evaluated at the current point
1394: . a_x - the gradient of auxiliary each field evaluated at the current point
1395: . t - current time
1396: . x - coordinates of the current point
1397: . n - unit normal at the current point
1398: - f0 - output values at the current point
1400: Level: intermediate
1402: .seealso: PetscDSSetBdResidual()
1403: @*/
1404: PetscErrorCode PetscDSGetBdResidual(PetscDS prob, PetscInt f,
1405: void (**f0)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1406: const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1407: const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1408: PetscReal t, const PetscReal x[], const PetscReal n[], PetscScalar f0[]),
1409: void (**f1)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1410: const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1411: const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1412: PetscReal t, const PetscReal x[], const PetscReal n[], PetscScalar f1[]))
1413: {
1416: 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);
1419: return(0);
1420: }
1424: /*@C
1425: PetscDSSetBdResidual - Get the pointwise boundary residual function for a given test field
1427: Not collective
1429: Input Parameters:
1430: + prob - The PetscDS
1431: . f - The test field number
1432: . f0 - boundary integrand for the test function term
1433: - f1 - boundary integrand for the test function gradient term
1435: Note: We are using a first order FEM model for the weak form:
1437: \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
1439: The calling sequence for the callbacks f0 and f1 is given by:
1441: $ f0(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1442: $ const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1443: $ const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1444: $ PetscReal t, const PetscReal x[], const PetscReal n[], PetscScalar f0[])
1446: + dim - the spatial dimension
1447: . Nf - the number of fields
1448: . uOff - the offset into u[] and u_t[] for each field
1449: . uOff_x - the offset into u_x[] for each field
1450: . u - each field evaluated at the current point
1451: . u_t - the time derivative of each field evaluated at the current point
1452: . u_x - the gradient of each field evaluated at the current point
1453: . aOff - the offset into a[] and a_t[] for each auxiliary field
1454: . aOff_x - the offset into a_x[] for each auxiliary field
1455: . a - each auxiliary field evaluated at the current point
1456: . a_t - the time derivative of each auxiliary field evaluated at the current point
1457: . a_x - the gradient of auxiliary each field evaluated at the current point
1458: . t - current time
1459: . x - coordinates of the current point
1460: . n - unit normal at the current point
1461: - f0 - output values at the current point
1463: Level: intermediate
1465: .seealso: PetscDSGetBdResidual()
1466: @*/
1467: PetscErrorCode PetscDSSetBdResidual(PetscDS prob, PetscInt f,
1468: void (*f0)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1469: const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1470: const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1471: PetscReal t, const PetscReal x[], const PetscReal n[], PetscScalar f0[]),
1472: void (*f1)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1473: const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1474: const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1475: PetscReal t, const PetscReal x[], const PetscReal n[], PetscScalar f1[]))
1476: {
1481: if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f);
1482: PetscDSEnlarge_Static(prob, f+1);
1485: return(0);
1486: }
1490: /*@C
1491: PetscDSGetBdJacobian - Get the pointwise boundary Jacobian function for given test and basis field
1493: Not collective
1495: Input Parameters:
1496: + prob - The PetscDS
1497: . f - The test field number
1498: - g - The field number
1500: Output Parameters:
1501: + g0 - integrand for the test and basis function term
1502: . g1 - integrand for the test function and basis function gradient term
1503: . g2 - integrand for the test function gradient and basis function term
1504: - g3 - integrand for the test function gradient and basis function gradient term
1506: Note: We are using a first order FEM model for the weak form:
1508: \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
1510: The calling sequence for the callbacks g0, g1, g2 and g3 is given by:
1512: $ g0(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1513: $ const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1514: $ const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1515: $ PetscReal t, const PetscReal x[], const PetscReal n[], PetscScalar g0[])
1517: + dim - the spatial dimension
1518: . Nf - the number of fields
1519: . uOff - the offset into u[] and u_t[] for each field
1520: . uOff_x - the offset into u_x[] for each field
1521: . u - each field evaluated at the current point
1522: . u_t - the time derivative of each field evaluated at the current point
1523: . u_x - the gradient of each field evaluated at the current point
1524: . aOff - the offset into a[] and a_t[] for each auxiliary field
1525: . aOff_x - the offset into a_x[] for each auxiliary field
1526: . a - each auxiliary field evaluated at the current point
1527: . a_t - the time derivative of each auxiliary field evaluated at the current point
1528: . a_x - the gradient of auxiliary each field evaluated at the current point
1529: . t - current time
1530: . u_tShift - the multiplier a for dF/dU_t
1531: . x - coordinates of the current point
1532: . n - normal at the current point
1533: - g0 - output values at the current point
1535: Level: intermediate
1537: .seealso: PetscDSSetBdJacobian()
1538: @*/
1539: PetscErrorCode PetscDSGetBdJacobian(PetscDS prob, PetscInt f, PetscInt g,
1540: void (**g0)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1541: const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1542: const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1543: PetscReal t, PetscReal u_tShift, const PetscReal x[], const PetscReal n[], PetscScalar g0[]),
1544: void (**g1)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1545: const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1546: const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1547: PetscReal t, PetscReal u_tShift, const PetscReal x[], const PetscReal n[], PetscScalar g1[]),
1548: void (**g2)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1549: const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1550: const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1551: PetscReal t, PetscReal u_tShift, const PetscReal x[], const PetscReal n[], PetscScalar g2[]),
1552: void (**g3)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1553: const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1554: const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1555: PetscReal t, PetscReal u_tShift, const PetscReal x[], const PetscReal n[], PetscScalar g3[]))
1556: {
1559: 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);
1560: 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);
1565: return(0);
1566: }
1570: /*@C
1571: PetscDSSetBdJacobian - Set the pointwise boundary Jacobian function for given test and basis field
1573: Not collective
1575: Input Parameters:
1576: + prob - The PetscDS
1577: . f - The test field number
1578: . g - The field number
1579: . g0 - integrand for the test and basis function term
1580: . g1 - integrand for the test function and basis function gradient term
1581: . g2 - integrand for the test function gradient and basis function term
1582: - g3 - integrand for the test function gradient and basis function gradient term
1584: Note: We are using a first order FEM model for the weak form:
1586: \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
1588: The calling sequence for the callbacks g0, g1, g2 and g3 is given by:
1590: $ g0(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1591: $ const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1592: $ const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1593: $ PetscReal t, const PetscReal x[], const PetscReal n[], PetscScalar g0[])
1595: + dim - the spatial dimension
1596: . Nf - the number of fields
1597: . uOff - the offset into u[] and u_t[] for each field
1598: . uOff_x - the offset into u_x[] for each field
1599: . u - each field evaluated at the current point
1600: . u_t - the time derivative of each field evaluated at the current point
1601: . u_x - the gradient of each field evaluated at the current point
1602: . aOff - the offset into a[] and a_t[] for each auxiliary field
1603: . aOff_x - the offset into a_x[] for each auxiliary field
1604: . a - each auxiliary field evaluated at the current point
1605: . a_t - the time derivative of each auxiliary field evaluated at the current point
1606: . a_x - the gradient of auxiliary each field evaluated at the current point
1607: . t - current time
1608: . u_tShift - the multiplier a for dF/dU_t
1609: . x - coordinates of the current point
1610: . n - normal at the current point
1611: - g0 - output values at the current point
1613: Level: intermediate
1615: .seealso: PetscDSGetBdJacobian()
1616: @*/
1617: PetscErrorCode PetscDSSetBdJacobian(PetscDS prob, PetscInt f, PetscInt g,
1618: void (*g0)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1619: const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1620: const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1621: PetscReal t, PetscReal u_tShift, const PetscReal x[], const PetscReal n[], PetscScalar g0[]),
1622: void (*g1)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1623: const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1624: const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1625: PetscReal t, PetscReal u_tShift, const PetscReal x[], const PetscReal n[], PetscScalar g1[]),
1626: void (*g2)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1627: const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1628: const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1629: PetscReal t, PetscReal u_tShift, const PetscReal x[], const PetscReal n[], PetscScalar g2[]),
1630: void (*g3)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1631: const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1632: const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1633: PetscReal t, PetscReal u_tShift, const PetscReal x[], const PetscReal n[], PetscScalar g3[]))
1634: {
1643: if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f);
1644: if (g < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", g);
1645: PetscDSEnlarge_Static(prob, PetscMax(f, g)+1);
1646: prob->gBd[(f*prob->Nf + g)*4+0] = g0;
1647: prob->gBd[(f*prob->Nf + g)*4+1] = g1;
1648: prob->gBd[(f*prob->Nf + g)*4+2] = g2;
1649: prob->gBd[(f*prob->Nf + g)*4+3] = g3;
1650: return(0);
1651: }
1655: /*@
1656: PetscDSGetFieldOffset - Returns the offset of the given field in the full space basis
1658: Not collective
1660: Input Parameters:
1661: + prob - The PetscDS object
1662: - f - The field number
1664: Output Parameter:
1665: . off - The offset
1667: Level: beginner
1669: .seealso: PetscDSGetBdFieldOffset(), PetscDSGetNumFields(), PetscDSCreate()
1670: @*/
1671: PetscErrorCode PetscDSGetFieldOffset(PetscDS prob, PetscInt f, PetscInt *off)
1672: {
1673: PetscInt g;
1679: 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);
1680: *off = 0;
1681: for (g = 0; g < f; ++g) {
1682: PetscFE fe = (PetscFE) prob->disc[g];
1683: PetscInt Nb, Nc;
1685: PetscFEGetDimension(fe, &Nb);
1686: PetscFEGetNumComponents(fe, &Nc);
1687: *off += Nb*Nc;
1688: }
1689: return(0);
1690: }
1694: /*@
1695: PetscDSGetBdFieldOffset - Returns the offset of the given field in the full space boundary basis
1697: Not collective
1699: Input Parameters:
1700: + prob - The PetscDS object
1701: - f - The field number
1703: Output Parameter:
1704: . off - The boundary offset
1706: Level: beginner
1708: .seealso: PetscDSGetFieldOffset(), PetscDSGetNumFields(), PetscDSCreate()
1709: @*/
1710: PetscErrorCode PetscDSGetBdFieldOffset(PetscDS prob, PetscInt f, PetscInt *off)
1711: {
1712: PetscInt g;
1718: 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);
1719: *off = 0;
1720: for (g = 0; g < f; ++g) {
1721: PetscFE fe = (PetscFE) prob->discBd[g];
1722: PetscInt Nb, Nc;
1724: PetscFEGetDimension(fe, &Nb);
1725: PetscFEGetNumComponents(fe, &Nc);
1726: *off += Nb*Nc;
1727: }
1728: return(0);
1729: }
1733: /*@
1734: PetscDSGetComponentOffset - Returns the offset of the given field on an evaluation point
1736: Not collective
1738: Input Parameters:
1739: + prob - The PetscDS object
1740: - f - The field number
1742: Output Parameter:
1743: . off - The offset
1745: Level: beginner
1747: .seealso: PetscDSGetBdFieldOffset(), PetscDSGetNumFields(), PetscDSCreate()
1748: @*/
1749: PetscErrorCode PetscDSGetComponentOffset(PetscDS prob, PetscInt f, PetscInt *off)
1750: {
1751: PetscInt g;
1757: 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);
1758: *off = 0;
1759: for (g = 0; g < f; ++g) {
1760: PetscFE fe = (PetscFE) prob->disc[g];
1761: PetscInt Nc;
1763: PetscFEGetNumComponents(fe, &Nc);
1764: *off += Nc;
1765: }
1766: return(0);
1767: }
1771: /*@
1772: PetscDSGetComponentOffsets - Returns the offset of each field on an evaluation point
1774: Not collective
1776: Input Parameter:
1777: . prob - The PetscDS object
1779: Output Parameter:
1780: . offsets - The offsets
1782: Level: beginner
1784: .seealso: PetscDSGetBdFieldOffset(), PetscDSGetNumFields(), PetscDSCreate()
1785: @*/
1786: PetscErrorCode PetscDSGetComponentOffsets(PetscDS prob, PetscInt *offsets[])
1787: {
1791: *offsets = prob->off;
1792: return(0);
1793: }
1797: /*@
1798: PetscDSGetComponentDerivativeOffsets - Returns the offset of each field derivative on an evaluation point
1800: Not collective
1802: Input Parameter:
1803: . prob - The PetscDS object
1805: Output Parameter:
1806: . offsets - The offsets
1808: Level: beginner
1810: .seealso: PetscDSGetBdFieldOffset(), PetscDSGetNumFields(), PetscDSCreate()
1811: @*/
1812: PetscErrorCode PetscDSGetComponentDerivativeOffsets(PetscDS prob, PetscInt *offsets[])
1813: {
1817: *offsets = prob->offDer;
1818: return(0);
1819: }
1823: /*@
1824: PetscDSGetComponentBdOffsets - Returns the offset of each field on a boundary evaluation point
1826: Not collective
1828: Input Parameter:
1829: . prob - The PetscDS object
1831: Output Parameter:
1832: . offsets - The offsets
1834: Level: beginner
1836: .seealso: PetscDSGetBdFieldOffset(), PetscDSGetNumFields(), PetscDSCreate()
1837: @*/
1838: PetscErrorCode PetscDSGetComponentBdOffsets(PetscDS prob, PetscInt *offsets[])
1839: {
1843: *offsets = prob->offBd;
1844: return(0);
1845: }
1849: /*@
1850: PetscDSGetComponentBdDerivativeOffsets - Returns the offset of each field derivative on a boundary evaluation point
1852: Not collective
1854: Input Parameter:
1855: . prob - The PetscDS object
1857: Output Parameter:
1858: . offsets - The offsets
1860: Level: beginner
1862: .seealso: PetscDSGetBdFieldOffset(), PetscDSGetNumFields(), PetscDSCreate()
1863: @*/
1864: PetscErrorCode PetscDSGetComponentBdDerivativeOffsets(PetscDS prob, PetscInt *offsets[])
1865: {
1869: *offsets = prob->offDerBd;
1870: return(0);
1871: }
1875: /*@C
1876: PetscDSGetTabulation - Return the basis tabulation at quadrature points for the volume discretization
1878: Not collective
1880: Input Parameter:
1881: . prob - The PetscDS object
1883: Output Parameters:
1884: + basis - The basis function tabulation at quadrature points
1885: - basisDer - The basis function derivative tabulation at quadrature points
1887: Level: intermediate
1889: .seealso: PetscDSGetBdTabulation(), PetscDSCreate()
1890: @*/
1891: PetscErrorCode PetscDSGetTabulation(PetscDS prob, PetscReal ***basis, PetscReal ***basisDer)
1892: {
1897: PetscDSSetUp(prob);
1900: return(0);
1901: }
1905: /*@C
1906: PetscDSGetBdTabulation - Return the basis tabulation at quadrature points for the boundary discretization
1908: Not collective
1910: Input Parameter:
1911: . prob - The PetscDS object
1913: Output Parameters:
1914: + basis - The basis function tabulation at quadrature points
1915: - basisDer - The basis function derivative tabulation at quadrature points
1917: Level: intermediate
1919: .seealso: PetscDSGetTabulation(), PetscDSCreate()
1920: @*/
1921: PetscErrorCode PetscDSGetBdTabulation(PetscDS prob, PetscReal ***basis, PetscReal ***basisDer)
1922: {
1927: PetscDSSetUp(prob);
1930: return(0);
1931: }
1935: PetscErrorCode PetscDSGetEvaluationArrays(PetscDS prob, PetscScalar **u, PetscScalar **u_t, PetscScalar **u_x)
1936: {
1941: PetscDSSetUp(prob);
1945: return(0);
1946: }
1950: PetscErrorCode PetscDSGetWeakFormArrays(PetscDS prob, PetscScalar **f0, PetscScalar **f1, PetscScalar **g0, PetscScalar **g1, PetscScalar **g2, PetscScalar **g3)
1951: {
1956: PetscDSSetUp(prob);
1963: return(0);
1964: }
1968: PetscErrorCode PetscDSGetRefCoordArrays(PetscDS prob, PetscReal **x, PetscScalar **refSpaceDer)
1969: {
1974: PetscDSSetUp(prob);
1977: return(0);
1978: }
1982: static PetscErrorCode PetscDSDestroy_Basic(PetscDS prob)
1983: {
1985: return(0);
1986: }
1990: static PetscErrorCode PetscDSInitialize_Basic(PetscDS prob)
1991: {
1993: prob->ops->setfromoptions = NULL;
1994: prob->ops->setup = NULL;
1995: prob->ops->view = NULL;
1996: prob->ops->destroy = PetscDSDestroy_Basic;
1997: return(0);
1998: }
2000: /*MC
2001: PETSCDSBASIC = "basic" - A discrete system with pointwise residual and boundary residual functions
2003: Level: intermediate
2005: .seealso: PetscDSType, PetscDSCreate(), PetscDSSetType()
2006: M*/
2010: PETSC_EXTERN PetscErrorCode PetscDSCreate_Basic(PetscDS prob)
2011: {
2012: PetscDS_Basic *b;
2013: PetscErrorCode ierr;
2017: PetscNewLog(prob, &b);
2018: prob->data = b;
2020: PetscDSInitialize_Basic(prob);
2021: return(0);
2022: }