Actual source code: fe.c
1: /* Basis Jet Tabulation
3: We would like to tabulate the nodal basis functions and derivatives at a set of points, usually quadrature points. We
4: follow here the derviation in http://www.math.ttu.edu/~kirby/papers/fiat-toms-2004.pdf. The nodal basis $\psi_i$ can
5: be expressed in terms of a prime basis $\phi_i$ which can be stably evaluated. In PETSc, we will use the Legendre basis
6: as a prime basis.
8: \psi_i = \sum_k \alpha_{ki} \phi_k
10: Our nodal basis is defined in terms of the dual basis $n_j$
12: n_j \cdot \psi_i = \delta_{ji}
14: and we may act on the first equation to obtain
16: n_j \cdot \psi_i = \sum_k \alpha_{ki} n_j \cdot \phi_k
17: \delta_{ji} = \sum_k \alpha_{ki} V_{jk}
18: I = V \alpha
20: so the coefficients of the nodal basis in the prime basis are
22: \alpha = V^{-1}
24: We will define the dual basis vectors $n_j$ using a quadrature rule.
26: Right now, we will just use the polynomial spaces P^k. I know some elements use the space of symmetric polynomials
27: (I think Nedelec), but we will neglect this for now. Constraints in the space, e.g. Arnold-Winther elements, can
28: be implemented exactly as in FIAT using functionals $L_j$.
30: I will have to count the degrees correctly for the Legendre product when we are on simplices.
32: We will have three objects:
33: - Space, P: this just need point evaluation I think
34: - Dual Space, P'+K: This looks like a set of functionals that can act on members of P, each n is defined by a Q
35: - FEM: This keeps {P, P', Q}
36: */
37: #include <petsc/private/petscfeimpl.h>
38: #include <petscdmplex.h>
40: PetscBool FEcite = PETSC_FALSE;
41: const char FECitation[] = "@article{kirby2004,\n"
42: " title = {Algorithm 839: FIAT, a New Paradigm for Computing Finite Element Basis Functions},\n"
43: " journal = {ACM Transactions on Mathematical Software},\n"
44: " author = {Robert C. Kirby},\n"
45: " volume = {30},\n"
46: " number = {4},\n"
47: " pages = {502--516},\n"
48: " doi = {10.1145/1039813.1039820},\n"
49: " year = {2004}\n}\n";
51: PetscClassId PETSCFE_CLASSID = 0;
53: PetscLogEvent PETSCFE_SetUp;
55: PetscFunctionList PetscFEList = NULL;
56: PetscBool PetscFERegisterAllCalled = PETSC_FALSE;
58: /*@C
59: PetscFERegister - Adds a new `PetscFEType`
61: Not Collective, No Fortran Support
63: Input Parameters:
64: + sname - The name of a new user-defined creation routine
65: - function - The creation routine
67: Example Usage:
68: .vb
69: PetscFERegister("my_fe", MyPetscFECreate);
70: .ve
72: Then, your PetscFE type can be chosen with the procedural interface via
73: .vb
74: PetscFECreate(MPI_Comm, PetscFE *);
75: PetscFESetType(PetscFE, "my_fe");
76: .ve
77: or at runtime via the option
78: .vb
79: -petscfe_type my_fe
80: .ve
82: Level: advanced
84: Note:
85: `PetscFERegister()` may be called multiple times to add several user-defined `PetscFE`s
87: .seealso: `PetscFE`, `PetscFEType`, `PetscFERegisterAll()`, `PetscFERegisterDestroy()`
88: @*/
89: PetscErrorCode PetscFERegister(const char sname[], PetscErrorCode (*function)(PetscFE))
90: {
91: PetscFunctionBegin;
92: PetscCall(PetscFunctionListAdd(&PetscFEList, sname, function));
93: PetscFunctionReturn(PETSC_SUCCESS);
94: }
96: /*@
97: PetscFESetType - Builds a particular `PetscFE`
99: Collective
101: Input Parameters:
102: + fem - The `PetscFE` object
103: - name - The kind of FEM space
105: Options Database Key:
106: . -petscfe_type type - Sets the `PetscFE` type; use -help for a list of available types
108: Level: intermediate
110: .seealso: `PetscFEType`, `PetscFE`, `PetscFEGetType()`, `PetscFECreate()`
111: @*/
112: PetscErrorCode PetscFESetType(PetscFE fem, PetscFEType name)
113: {
114: PetscErrorCode (*r)(PetscFE);
115: PetscBool match;
117: PetscFunctionBegin;
119: PetscCall(PetscObjectTypeCompare((PetscObject)fem, name, &match));
120: if (match) PetscFunctionReturn(PETSC_SUCCESS);
122: if (!PetscFERegisterAllCalled) PetscCall(PetscFERegisterAll());
123: PetscCall(PetscFunctionListFind(PetscFEList, name, &r));
124: PetscCheck(r, PetscObjectComm((PetscObject)fem), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown PetscFE type: %s", name);
126: PetscTryTypeMethod(fem, destroy);
127: fem->ops->destroy = NULL;
129: PetscCall((*r)(fem));
130: PetscCall(PetscObjectChangeTypeName((PetscObject)fem, name));
131: PetscFunctionReturn(PETSC_SUCCESS);
132: }
134: /*@
135: PetscFEGetType - Gets the `PetscFEType` (as a string) from the `PetscFE` object.
137: Not Collective
139: Input Parameter:
140: . fem - The `PetscFE`
142: Output Parameter:
143: . name - The `PetscFEType` name
145: Level: intermediate
147: .seealso: `PetscFEType`, `PetscFE`, `PetscFESetType()`, `PetscFECreate()`
148: @*/
149: PetscErrorCode PetscFEGetType(PetscFE fem, PetscFEType *name)
150: {
151: PetscFunctionBegin;
153: PetscAssertPointer(name, 2);
154: if (!PetscFERegisterAllCalled) PetscCall(PetscFERegisterAll());
155: *name = ((PetscObject)fem)->type_name;
156: PetscFunctionReturn(PETSC_SUCCESS);
157: }
159: /*@
160: PetscFEViewFromOptions - View from a `PetscFE` based on values in the options database
162: Collective
164: Input Parameters:
165: + A - the `PetscFE` object
166: . obj - Optional object that provides the options prefix, pass `NULL` to use the options prefix of `A`
167: - name - command line option name
169: Level: intermediate
171: .seealso: `PetscFE`, `PetscFEView()`, `PetscObjectViewFromOptions()`, `PetscFECreate()`
172: @*/
173: PetscErrorCode PetscFEViewFromOptions(PetscFE A, PeOp PetscObject obj, const char name[])
174: {
175: PetscFunctionBegin;
177: PetscCall(PetscObjectViewFromOptions((PetscObject)A, obj, name));
178: PetscFunctionReturn(PETSC_SUCCESS);
179: }
181: /*@
182: PetscFEView - Views a `PetscFE`
184: Collective
186: Input Parameters:
187: + fem - the `PetscFE` object to view
188: - viewer - the viewer
190: Level: beginner
192: .seealso: `PetscFE`, `PetscViewer`, `PetscFEDestroy()`, `PetscFEViewFromOptions()`
193: @*/
194: PetscErrorCode PetscFEView(PetscFE fem, PetscViewer viewer)
195: {
196: PetscBool isascii;
198: PetscFunctionBegin;
201: if (!viewer) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)fem), &viewer));
202: PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)fem, viewer));
203: PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
204: PetscTryTypeMethod(fem, view, viewer);
205: PetscFunctionReturn(PETSC_SUCCESS);
206: }
208: /*@
209: PetscFESetFromOptions - sets parameters in a `PetscFE` from the options database
211: Collective
213: Input Parameter:
214: . fem - the `PetscFE` object to set options for
216: Options Database Keys:
217: + -petscfe_num_blocks - the number of cell blocks to integrate concurrently
218: - -petscfe_num_batches - the number of cell batches to integrate serially
220: Level: intermediate
222: .seealso: `PetscFE`, `PetscFEView()`
223: @*/
224: PetscErrorCode PetscFESetFromOptions(PetscFE fem)
225: {
226: const char *defaultType;
227: char name[256];
228: PetscBool flg;
230: PetscFunctionBegin;
232: if (!((PetscObject)fem)->type_name) defaultType = PETSCFEBASIC;
233: else defaultType = ((PetscObject)fem)->type_name;
234: if (!PetscFERegisterAllCalled) PetscCall(PetscFERegisterAll());
236: PetscObjectOptionsBegin((PetscObject)fem);
237: PetscCall(PetscOptionsFList("-petscfe_type", "Finite element space", "PetscFESetType", PetscFEList, defaultType, name, 256, &flg));
238: if (flg) PetscCall(PetscFESetType(fem, name));
239: else if (!((PetscObject)fem)->type_name) PetscCall(PetscFESetType(fem, defaultType));
240: PetscCall(PetscOptionsBoundedInt("-petscfe_num_blocks", "The number of cell blocks to integrate concurrently", "PetscSpaceSetTileSizes", fem->numBlocks, &fem->numBlocks, NULL, 1));
241: PetscCall(PetscOptionsBoundedInt("-petscfe_num_batches", "The number of cell batches to integrate serially", "PetscSpaceSetTileSizes", fem->numBatches, &fem->numBatches, NULL, 1));
242: PetscTryTypeMethod(fem, setfromoptions, PetscOptionsObject);
243: /* process any options handlers added with PetscObjectAddOptionsHandler() */
244: PetscCall(PetscObjectProcessOptionsHandlers((PetscObject)fem, PetscOptionsObject));
245: PetscOptionsEnd();
246: PetscCall(PetscFEViewFromOptions(fem, NULL, "-petscfe_view"));
247: PetscFunctionReturn(PETSC_SUCCESS);
248: }
250: /*@
251: PetscFESetUp - Construct data structures for the `PetscFE` after the `PetscFEType` has been set
253: Collective
255: Input Parameter:
256: . fem - the `PetscFE` object to setup
258: Level: intermediate
260: .seealso: `PetscFE`, `PetscFEView()`, `PetscFEDestroy()`
261: @*/
262: PetscErrorCode PetscFESetUp(PetscFE fem)
263: {
264: PetscFunctionBegin;
266: if (fem->setupcalled) PetscFunctionReturn(PETSC_SUCCESS);
267: PetscCall(PetscLogEventBegin(PETSCFE_SetUp, fem, 0, 0, 0));
268: fem->setupcalled = PETSC_TRUE;
269: PetscTryTypeMethod(fem, setup);
270: PetscCall(PetscLogEventEnd(PETSCFE_SetUp, fem, 0, 0, 0));
271: PetscFunctionReturn(PETSC_SUCCESS);
272: }
274: /*@
275: PetscFEDestroy - Destroys a `PetscFE` object
277: Collective
279: Input Parameter:
280: . fem - the `PetscFE` object to destroy
282: Level: beginner
284: .seealso: `PetscFE`, `PetscFEView()`
285: @*/
286: PetscErrorCode PetscFEDestroy(PetscFE *fem)
287: {
288: PetscFunctionBegin;
289: if (!*fem) PetscFunctionReturn(PETSC_SUCCESS);
292: if (--((PetscObject)*fem)->refct > 0) {
293: *fem = NULL;
294: PetscFunctionReturn(PETSC_SUCCESS);
295: }
296: ((PetscObject)*fem)->refct = 0;
298: if ((*fem)->subspaces) {
299: PetscInt dim, d;
301: PetscCall(PetscDualSpaceGetDimension((*fem)->dualSpace, &dim));
302: for (d = 0; d < dim; ++d) PetscCall(PetscFEDestroy(&(*fem)->subspaces[d]));
303: }
304: PetscCall(PetscFree((*fem)->subspaces));
305: PetscCall(PetscFree((*fem)->invV));
306: PetscCall(PetscTabulationDestroy(&(*fem)->T));
307: PetscCall(PetscTabulationDestroy(&(*fem)->Tf));
308: PetscCall(PetscTabulationDestroy(&(*fem)->Tc));
309: PetscCall(PetscSpaceDestroy(&(*fem)->basisSpace));
310: PetscCall(PetscDualSpaceDestroy(&(*fem)->dualSpace));
311: PetscCall(PetscQuadratureDestroy(&(*fem)->quadrature));
312: PetscCall(PetscQuadratureDestroy(&(*fem)->faceQuadrature));
313: #ifdef PETSC_HAVE_LIBCEED
314: PetscCallCEED(CeedBasisDestroy(&(*fem)->ceedBasis));
315: PetscCallCEED(CeedDestroy(&(*fem)->ceed));
316: #endif
318: PetscTryTypeMethod(*fem, destroy);
319: PetscCall(PetscHeaderDestroy(fem));
320: PetscFunctionReturn(PETSC_SUCCESS);
321: }
323: /*@
324: PetscFECreate - Creates an empty `PetscFE` object. The type can then be set with `PetscFESetType()`.
326: Collective
328: Input Parameter:
329: . comm - The communicator for the `PetscFE` object
331: Output Parameter:
332: . fem - The `PetscFE` object
334: Level: beginner
336: .seealso: `PetscFE`, `PetscFEType`, `PetscFESetType()`, `PetscFECreateDefault()`, `PETSCFEGALERKIN`
337: @*/
338: PetscErrorCode PetscFECreate(MPI_Comm comm, PetscFE *fem)
339: {
340: PetscFE f;
342: PetscFunctionBegin;
343: PetscAssertPointer(fem, 2);
344: PetscCall(PetscCitationsRegister(FECitation, &FEcite));
345: PetscCall(PetscFEInitializePackage());
347: PetscCall(PetscHeaderCreate(f, PETSCFE_CLASSID, "PetscFE", "Finite Element", "PetscFE", comm, PetscFEDestroy, PetscFEView));
349: f->basisSpace = NULL;
350: f->dualSpace = NULL;
351: f->numComponents = 1;
352: f->subspaces = NULL;
353: f->invV = NULL;
354: f->T = NULL;
355: f->Tf = NULL;
356: f->Tc = NULL;
357: PetscCall(PetscArrayzero(&f->quadrature, 1));
358: PetscCall(PetscArrayzero(&f->faceQuadrature, 1));
359: f->blockSize = 0;
360: f->numBlocks = 1;
361: f->batchSize = 0;
362: f->numBatches = 1;
364: *fem = f;
365: PetscFunctionReturn(PETSC_SUCCESS);
366: }
368: /*@
369: PetscFEGetSpatialDimension - Returns the spatial dimension of the element
371: Not Collective
373: Input Parameter:
374: . fem - The `PetscFE` object
376: Output Parameter:
377: . dim - The spatial dimension
379: Level: intermediate
381: .seealso: `PetscFE`, `PetscFECreate()`
382: @*/
383: PetscErrorCode PetscFEGetSpatialDimension(PetscFE fem, PetscInt *dim)
384: {
385: DM dm;
387: PetscFunctionBegin;
389: PetscAssertPointer(dim, 2);
390: PetscCall(PetscDualSpaceGetDM(fem->dualSpace, &dm));
391: PetscCall(DMGetDimension(dm, dim));
392: PetscFunctionReturn(PETSC_SUCCESS);
393: }
395: /*@
396: PetscFESetNumComponents - Sets the number of field components in the element
398: Not Collective
400: Input Parameters:
401: + fem - The `PetscFE` object
402: - comp - The number of field components
404: Level: intermediate
406: .seealso: `PetscFE`, `PetscFECreate()`, `PetscFEGetSpatialDimension()`, `PetscFEGetNumComponents()`
407: @*/
408: PetscErrorCode PetscFESetNumComponents(PetscFE fem, PetscInt comp)
409: {
410: PetscFunctionBegin;
412: fem->numComponents = comp;
413: PetscFunctionReturn(PETSC_SUCCESS);
414: }
416: /*@
417: PetscFEGetNumComponents - Returns the number of components in the element
419: Not Collective
421: Input Parameter:
422: . fem - The `PetscFE` object
424: Output Parameter:
425: . comp - The number of field components
427: Level: intermediate
429: .seealso: `PetscFE`, `PetscFECreate()`, `PetscFEGetSpatialDimension()`
430: @*/
431: PetscErrorCode PetscFEGetNumComponents(PetscFE fem, PetscInt *comp)
432: {
433: PetscFunctionBegin;
435: PetscAssertPointer(comp, 2);
436: *comp = fem->numComponents;
437: PetscFunctionReturn(PETSC_SUCCESS);
438: }
440: /*@
441: PetscFESetTileSizes - Sets the tile sizes for evaluation
443: Not Collective
445: Input Parameters:
446: + fem - The `PetscFE` object
447: . blockSize - The number of elements in a block
448: . numBlocks - The number of blocks in a batch
449: . batchSize - The number of elements in a batch
450: - numBatches - The number of batches in a chunk
452: Level: intermediate
454: .seealso: `PetscFE`, `PetscFECreate()`, `PetscFEGetTileSizes()`
455: @*/
456: PetscErrorCode PetscFESetTileSizes(PetscFE fem, PetscInt blockSize, PetscInt numBlocks, PetscInt batchSize, PetscInt numBatches)
457: {
458: PetscFunctionBegin;
460: fem->blockSize = blockSize;
461: fem->numBlocks = numBlocks;
462: fem->batchSize = batchSize;
463: fem->numBatches = numBatches;
464: PetscFunctionReturn(PETSC_SUCCESS);
465: }
467: /*@
468: PetscFEGetTileSizes - Returns the tile sizes for evaluation
470: Not Collective
472: Input Parameter:
473: . fem - The `PetscFE` object
475: Output Parameters:
476: + blockSize - The number of elements in a block, pass `NULL` if not needed
477: . numBlocks - The number of blocks in a batch, pass `NULL` if not needed
478: . batchSize - The number of elements in a batch, pass `NULL` if not needed
479: - numBatches - The number of batches in a chunk, pass `NULL` if not needed
481: Level: intermediate
483: .seealso: `PetscFE`, `PetscFECreate()`, `PetscFESetTileSizes()`
484: @*/
485: PetscErrorCode PetscFEGetTileSizes(PetscFE fem, PeOp PetscInt *blockSize, PeOp PetscInt *numBlocks, PeOp PetscInt *batchSize, PeOp PetscInt *numBatches)
486: {
487: PetscFunctionBegin;
489: if (blockSize) PetscAssertPointer(blockSize, 2);
490: if (numBlocks) PetscAssertPointer(numBlocks, 3);
491: if (batchSize) PetscAssertPointer(batchSize, 4);
492: if (numBatches) PetscAssertPointer(numBatches, 5);
493: if (blockSize) *blockSize = fem->blockSize;
494: if (numBlocks) *numBlocks = fem->numBlocks;
495: if (batchSize) *batchSize = fem->batchSize;
496: if (numBatches) *numBatches = fem->numBatches;
497: PetscFunctionReturn(PETSC_SUCCESS);
498: }
500: /*@
501: PetscFEGetBasisSpace - Returns the `PetscSpace` used for the approximation of the solution for the `PetscFE`
503: Not Collective
505: Input Parameter:
506: . fem - The `PetscFE` object
508: Output Parameter:
509: . sp - The `PetscSpace` object
511: Level: intermediate
513: .seealso: `PetscFE`, `PetscSpace`, `PetscFECreate()`
514: @*/
515: PetscErrorCode PetscFEGetBasisSpace(PetscFE fem, PetscSpace *sp)
516: {
517: PetscFunctionBegin;
519: PetscAssertPointer(sp, 2);
520: *sp = fem->basisSpace;
521: PetscFunctionReturn(PETSC_SUCCESS);
522: }
524: /*@
525: PetscFESetBasisSpace - Sets the `PetscSpace` used for the approximation of the solution
527: Not Collective
529: Input Parameters:
530: + fem - The `PetscFE` object
531: - sp - The `PetscSpace` object
533: Level: intermediate
535: Developer Notes:
536: There is `PetscFESetBasisSpace()` but the `PetscFESetDualSpace()`, likely the Basis is unneeded in the function name
538: .seealso: `PetscFE`, `PetscSpace`, `PetscDualSpace`, `PetscFECreate()`, `PetscFESetDualSpace()`
539: @*/
540: PetscErrorCode PetscFESetBasisSpace(PetscFE fem, PetscSpace sp)
541: {
542: PetscFunctionBegin;
545: PetscCall(PetscSpaceDestroy(&fem->basisSpace));
546: fem->basisSpace = sp;
547: PetscCall(PetscObjectReference((PetscObject)fem->basisSpace));
548: PetscFunctionReturn(PETSC_SUCCESS);
549: }
551: /*@
552: PetscFEGetDualSpace - Returns the `PetscDualSpace` used to define the inner product for a `PetscFE`
554: Not Collective
556: Input Parameter:
557: . fem - The `PetscFE` object
559: Output Parameter:
560: . sp - The `PetscDualSpace` object
562: Level: intermediate
564: .seealso: `PetscFE`, `PetscSpace`, `PetscDualSpace`, `PetscFECreate()`
565: @*/
566: PetscErrorCode PetscFEGetDualSpace(PetscFE fem, PetscDualSpace *sp)
567: {
568: PetscFunctionBegin;
570: PetscAssertPointer(sp, 2);
571: *sp = fem->dualSpace;
572: PetscFunctionReturn(PETSC_SUCCESS);
573: }
575: /*@
576: PetscFESetDualSpace - Sets the `PetscDualSpace` used to define the inner product
578: Not Collective
580: Input Parameters:
581: + fem - The `PetscFE` object
582: - sp - The `PetscDualSpace` object
584: Level: intermediate
586: .seealso: `PetscFE`, `PetscSpace`, `PetscDualSpace`, `PetscFECreate()`, `PetscFESetBasisSpace()`
587: @*/
588: PetscErrorCode PetscFESetDualSpace(PetscFE fem, PetscDualSpace sp)
589: {
590: PetscFunctionBegin;
593: PetscCall(PetscDualSpaceDestroy(&fem->dualSpace));
594: fem->dualSpace = sp;
595: PetscCall(PetscObjectReference((PetscObject)fem->dualSpace));
596: PetscFunctionReturn(PETSC_SUCCESS);
597: }
599: /*@
600: PetscFEGetQuadrature - Returns the `PetscQuadrature` used to calculate inner products
602: Not Collective
604: Input Parameter:
605: . fem - The `PetscFE` object
607: Output Parameter:
608: . q - The `PetscQuadrature` object
610: Level: intermediate
612: .seealso: `PetscFE`, `PetscSpace`, `PetscDualSpace`, `PetscQuadrature`, `PetscFECreate()`
613: @*/
614: PetscErrorCode PetscFEGetQuadrature(PetscFE fem, PetscQuadrature *q)
615: {
616: PetscFunctionBegin;
618: PetscAssertPointer(q, 2);
619: *q = fem->quadrature;
620: PetscFunctionReturn(PETSC_SUCCESS);
621: }
623: /*@
624: PetscFESetQuadrature - Sets the `PetscQuadrature` used to calculate inner products
626: Not Collective
628: Input Parameters:
629: + fem - The `PetscFE` object
630: - q - The `PetscQuadrature` object
632: Level: intermediate
634: .seealso: `PetscFE`, `PetscSpace`, `PetscDualSpace`, `PetscQuadrature`, `PetscFECreate()`, `PetscFEGetFaceQuadrature()`
635: @*/
636: PetscErrorCode PetscFESetQuadrature(PetscFE fem, PetscQuadrature q)
637: {
638: PetscInt Nc, qNc;
640: PetscFunctionBegin;
642: if (q == fem->quadrature) PetscFunctionReturn(PETSC_SUCCESS);
643: PetscCall(PetscFEGetNumComponents(fem, &Nc));
644: PetscCall(PetscQuadratureGetNumComponents(q, &qNc));
645: PetscCheck(!(qNc != 1) || !(Nc != qNc), PetscObjectComm((PetscObject)fem), PETSC_ERR_ARG_SIZ, "FE components %" PetscInt_FMT " != Quadrature components %" PetscInt_FMT " and non-scalar quadrature", Nc, qNc);
646: PetscCall(PetscTabulationDestroy(&fem->T));
647: PetscCall(PetscTabulationDestroy(&fem->Tc));
648: PetscCall(PetscObjectReference((PetscObject)q));
649: PetscCall(PetscQuadratureDestroy(&fem->quadrature));
650: fem->quadrature = q;
651: PetscFunctionReturn(PETSC_SUCCESS);
652: }
654: /*@
655: PetscFEGetFaceQuadrature - Returns the `PetscQuadrature` used to calculate inner products on faces
657: Not Collective
659: Input Parameter:
660: . fem - The `PetscFE` object
662: Output Parameter:
663: . q - The `PetscQuadrature` object
665: Level: intermediate
667: Developer Notes:
668: There is a special face quadrature but not edge, likely this API would benefit from a refactorization
670: .seealso: `PetscFE`, `PetscSpace`, `PetscDualSpace`, `PetscQuadrature`, `PetscFECreate()`, `PetscFESetQuadrature()`, `PetscFESetFaceQuadrature()`
671: @*/
672: PetscErrorCode PetscFEGetFaceQuadrature(PetscFE fem, PetscQuadrature *q)
673: {
674: PetscFunctionBegin;
676: PetscAssertPointer(q, 2);
677: *q = fem->faceQuadrature;
678: PetscFunctionReturn(PETSC_SUCCESS);
679: }
681: /*@
682: PetscFESetFaceQuadrature - Sets the `PetscQuadrature` used to calculate inner products on faces
684: Not Collective
686: Input Parameters:
687: + fem - The `PetscFE` object
688: - q - The `PetscQuadrature` object
690: Level: intermediate
692: .seealso: `PetscFE`, `PetscSpace`, `PetscDualSpace`, `PetscQuadrature`, `PetscFECreate()`, `PetscFESetQuadrature()`
693: @*/
694: PetscErrorCode PetscFESetFaceQuadrature(PetscFE fem, PetscQuadrature q)
695: {
696: PetscInt Nc, qNc;
698: PetscFunctionBegin;
700: if (q == fem->faceQuadrature) PetscFunctionReturn(PETSC_SUCCESS);
701: PetscCall(PetscFEGetNumComponents(fem, &Nc));
702: PetscCall(PetscQuadratureGetNumComponents(q, &qNc));
703: PetscCheck(!(qNc != 1) || !(Nc != qNc), PetscObjectComm((PetscObject)fem), PETSC_ERR_ARG_SIZ, "FE components %" PetscInt_FMT " != Quadrature components %" PetscInt_FMT " and non-scalar quadrature", Nc, qNc);
704: PetscCall(PetscTabulationDestroy(&fem->Tf));
705: PetscCall(PetscObjectReference((PetscObject)q));
706: PetscCall(PetscQuadratureDestroy(&fem->faceQuadrature));
707: fem->faceQuadrature = q;
708: PetscFunctionReturn(PETSC_SUCCESS);
709: }
711: /*@
712: PetscFECopyQuadrature - Copy both volumetric and surface quadrature to a new `PetscFE`
714: Not Collective
716: Input Parameters:
717: + sfe - The `PetscFE` source for the quadratures
718: - tfe - The `PetscFE` target for the quadratures
720: Level: intermediate
722: .seealso: `PetscFE`, `PetscSpace`, `PetscDualSpace`, `PetscQuadrature`, `PetscFECreate()`, `PetscFESetQuadrature()`, `PetscFESetFaceQuadrature()`
723: @*/
724: PetscErrorCode PetscFECopyQuadrature(PetscFE sfe, PetscFE tfe)
725: {
726: PetscQuadrature q;
728: PetscFunctionBegin;
731: PetscCall(PetscFEGetQuadrature(sfe, &q));
732: PetscCall(PetscFESetQuadrature(tfe, q));
733: PetscCall(PetscFEGetFaceQuadrature(sfe, &q));
734: PetscCall(PetscFESetFaceQuadrature(tfe, q));
735: PetscFunctionReturn(PETSC_SUCCESS);
736: }
738: /*@C
739: PetscFEGetNumDof - Returns the number of dofs (dual basis vectors) associated to mesh points on the reference cell of a given dimension
741: Not Collective
743: Input Parameter:
744: . fem - The `PetscFE` object
746: Output Parameter:
747: . numDof - Array of length `dim` with the number of dofs in each dimension
749: Level: intermediate
751: .seealso: `PetscFE`, `PetscSpace`, `PetscDualSpace`, `PetscFECreate()`
752: @*/
753: PetscErrorCode PetscFEGetNumDof(PetscFE fem, const PetscInt *numDof[])
754: {
755: PetscFunctionBegin;
757: PetscAssertPointer(numDof, 2);
758: PetscCall(PetscDualSpaceGetNumDof(fem->dualSpace, numDof));
759: PetscFunctionReturn(PETSC_SUCCESS);
760: }
762: /*@C
763: PetscFEGetCellTabulation - Returns the tabulation of the basis functions at the quadrature points on the reference cell
765: Not Collective
767: Input Parameters:
768: + fem - The `PetscFE` object
769: - k - The highest derivative we need to tabulate, very often 1
771: Output Parameter:
772: . T - The basis function values and derivatives at quadrature points
774: Level: intermediate
776: Note:
777: .vb
778: T->T[0] = B[(p*pdim + i)*Nc + c] is the value at point p for basis function i and component c
779: T->T[1] = D[((p*pdim + i)*Nc + c)*dim + d] is the derivative value at point p for basis function i, component c, in direction d
780: T->T[2] = H[(((p*pdim + i)*Nc + c)*dim + d)*dim + e] is the value at point p for basis function i, component c, in directions d and e
781: .ve
783: .seealso: `PetscFE`, `PetscSpace`, `PetscDualSpace`, `PetscTabulation`, `PetscFECreateTabulation()`, `PetscTabulationDestroy()`
784: @*/
785: PetscErrorCode PetscFEGetCellTabulation(PetscFE fem, PetscInt k, PetscTabulation *T)
786: {
787: PetscInt npoints;
788: const PetscReal *points;
790: PetscFunctionBegin;
792: PetscAssertPointer(T, 3);
793: PetscCall(PetscQuadratureGetData(fem->quadrature, NULL, NULL, &npoints, &points, NULL));
794: if (!fem->T) PetscCall(PetscFECreateTabulation(fem, 1, npoints, points, k, &fem->T));
795: PetscCheck(!fem->T || k <= fem->T->K || (!fem->T->cdim && !fem->T->K), PetscObjectComm((PetscObject)fem), PETSC_ERR_ARG_OUTOFRANGE, "Requested %" PetscInt_FMT " derivatives, but only tabulated %" PetscInt_FMT, k, fem->T->K);
796: *T = fem->T;
797: PetscFunctionReturn(PETSC_SUCCESS);
798: }
800: PetscErrorCode PetscFEExpandFaceQuadrature(PetscFE fe, PetscQuadrature fq, PetscQuadrature *efq)
801: {
802: DM dm;
803: PetscDualSpace sp;
804: const PetscInt *faces;
805: const PetscReal *points, *weights;
806: DMPolytopeType ct;
807: PetscReal *facePoints, *faceWeights;
808: PetscInt dim, cStart, Nf, Nc, Np, order;
810: PetscFunctionBegin;
811: PetscCall(PetscFEGetDualSpace(fe, &sp));
812: PetscCall(PetscDualSpaceGetDM(sp, &dm));
813: PetscCall(DMGetDimension(dm, &dim));
814: PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, NULL));
815: PetscCall(DMPlexGetConeSize(dm, cStart, &Nf));
816: PetscCall(DMPlexGetCone(dm, cStart, &faces));
817: PetscCall(PetscQuadratureGetData(fq, NULL, &Nc, &Np, &points, &weights));
818: PetscCall(PetscMalloc1(Nf * Np * dim, &facePoints));
819: PetscCall(PetscMalloc1(Nf * Np * Nc, &faceWeights));
820: for (PetscInt f = 0; f < Nf; ++f) {
821: const PetscReal xi0[3] = {-1., -1., -1.};
822: PetscReal v0[3], J[9], detJ;
824: PetscCall(DMPlexComputeCellGeometryFEM(dm, faces[f], NULL, v0, J, NULL, &detJ));
825: for (PetscInt q = 0; q < Np; ++q) {
826: CoordinatesRefToReal(dim, dim - 1, xi0, v0, J, &points[q * (dim - 1)], &facePoints[(f * Np + q) * dim]);
827: for (PetscInt c = 0; c < Nc; ++c) faceWeights[(f * Np + q) * Nc + c] = weights[q * Nc + c];
828: }
829: }
830: PetscCall(PetscQuadratureCreate(PetscObjectComm((PetscObject)fq), efq));
831: PetscCall(PetscQuadratureGetCellType(fq, &ct));
832: PetscCall(PetscQuadratureSetCellType(*efq, ct));
833: PetscCall(PetscQuadratureGetOrder(fq, &order));
834: PetscCall(PetscQuadratureSetOrder(*efq, order));
835: PetscCall(PetscQuadratureSetData(*efq, dim, Nc, Nf * Np, facePoints, faceWeights));
836: PetscFunctionReturn(PETSC_SUCCESS);
837: }
839: /*@C
840: PetscFEGetFaceTabulation - Returns the tabulation of the basis functions at the face quadrature points for each face of the reference cell
842: Not Collective
844: Input Parameters:
845: + fem - The `PetscFE` object
846: - k - The highest derivative we need to tabulate, very often 1
848: Output Parameter:
849: . Tf - The basis function values and derivatives at face quadrature points
851: Level: intermediate
853: Note:
854: .vb
855: T->T[0] = Bf[((f*Nq + q)*pdim + i)*Nc + c] is the value at point f,q for basis function i and component c
856: T->T[1] = Df[(((f*Nq + q)*pdim + i)*Nc + c)*dim + d] is the derivative value at point f,q for basis function i, component c, in direction d
857: T->T[2] = Hf[((((f*Nq + q)*pdim + i)*Nc + c)*dim + d)*dim + e] is the value at point f,q for basis function i, component c, in directions d and e
858: .ve
860: .seealso: `PetscFE`, `PetscSpace`, `PetscDualSpace`, `PetscTabulation`, `PetscFEGetCellTabulation()`, `PetscFECreateTabulation()`, `PetscTabulationDestroy()`
861: @*/
862: PetscErrorCode PetscFEGetFaceTabulation(PetscFE fem, PetscInt k, PetscTabulation *Tf)
863: {
864: PetscFunctionBegin;
866: PetscAssertPointer(Tf, 3);
867: if (!fem->Tf) {
868: PetscQuadrature fq;
870: PetscCall(PetscFEGetFaceQuadrature(fem, &fq));
871: if (fq) {
872: PetscQuadrature efq;
873: const PetscReal *facePoints;
874: PetscInt Np, eNp;
876: PetscCall(PetscFEExpandFaceQuadrature(fem, fq, &efq));
877: PetscCall(PetscQuadratureGetData(fq, NULL, NULL, &Np, NULL, NULL));
878: PetscCall(PetscQuadratureGetData(efq, NULL, NULL, &eNp, &facePoints, NULL));
879: if (PetscDefined(USE_DEBUG)) {
880: PetscDualSpace sp;
881: DM dm;
882: PetscInt cStart, Nf;
884: PetscCall(PetscFEGetDualSpace(fem, &sp));
885: PetscCall(PetscDualSpaceGetDM(sp, &dm));
886: PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, NULL));
887: PetscCall(DMPlexGetConeSize(dm, cStart, &Nf));
888: PetscCheck(Nf == eNp / Np, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Number of faces %" PetscInt_FMT " != %" PetscInt_FMT " number of quadrature replicas", Nf, eNp / Np);
889: }
890: PetscCall(PetscFECreateTabulation(fem, eNp / Np, Np, facePoints, k, &fem->Tf));
891: PetscCall(PetscQuadratureDestroy(&efq));
892: }
893: }
894: PetscCheck(!fem->Tf || k <= fem->Tf->K, PetscObjectComm((PetscObject)fem), PETSC_ERR_ARG_OUTOFRANGE, "Requested %" PetscInt_FMT " derivatives, but only tabulated %" PetscInt_FMT, k, fem->Tf->K);
895: *Tf = fem->Tf;
896: PetscFunctionReturn(PETSC_SUCCESS);
897: }
899: /*@C
900: PetscFEGetFaceCentroidTabulation - Returns the tabulation of the basis functions at the face centroid points
902: Not Collective
904: Input Parameter:
905: . fem - The `PetscFE` object
907: Output Parameter:
908: . Tc - The basis function values at face centroid points
910: Level: intermediate
912: Note:
913: .vb
914: T->T[0] = Bf[(f*pdim + i)*Nc + c] is the value at point f for basis function i and component c
915: .ve
917: .seealso: `PetscFE`, `PetscSpace`, `PetscDualSpace`, `PetscTabulation`, `PetscFEGetFaceTabulation()`, `PetscFEGetCellTabulation()`, `PetscFECreateTabulation()`, `PetscTabulationDestroy()`
918: @*/
919: PetscErrorCode PetscFEGetFaceCentroidTabulation(PetscFE fem, PetscTabulation *Tc)
920: {
921: PetscFunctionBegin;
923: PetscAssertPointer(Tc, 2);
924: if (!fem->Tc) {
925: PetscDualSpace sp;
926: DM dm;
927: const PetscInt *cone;
928: PetscReal *centroids;
929: PetscInt dim, numFaces, f;
931: PetscCall(PetscFEGetDualSpace(fem, &sp));
932: PetscCall(PetscDualSpaceGetDM(sp, &dm));
933: PetscCall(DMGetDimension(dm, &dim));
934: PetscCall(DMPlexGetConeSize(dm, 0, &numFaces));
935: PetscCall(DMPlexGetCone(dm, 0, &cone));
936: PetscCall(PetscMalloc1(numFaces * dim, ¢roids));
937: for (f = 0; f < numFaces; ++f) PetscCall(DMPlexComputeCellGeometryFVM(dm, cone[f], NULL, ¢roids[f * dim], NULL));
938: PetscCall(PetscFECreateTabulation(fem, 1, numFaces, centroids, 0, &fem->Tc));
939: PetscCall(PetscFree(centroids));
940: }
941: *Tc = fem->Tc;
942: PetscFunctionReturn(PETSC_SUCCESS);
943: }
945: /*@C
946: PetscFECreateTabulation - Tabulates the basis functions, and perhaps derivatives, at the points provided.
948: Not Collective
950: Input Parameters:
951: + fem - The `PetscFE` object
952: . nrepl - The number of replicas
953: . npoints - The number of tabulation points in a replica
954: . points - The tabulation point coordinates
955: - K - The number of derivatives calculated
957: Output Parameter:
958: . T - The basis function values and derivatives at tabulation points
960: Level: intermediate
962: Note:
963: .vb
964: T->T[0] = B[(p*pdim + i)*Nc + c] is the value at point p for basis function i and component c
965: T->T[1] = D[((p*pdim + i)*Nc + c)*dim + d] is the derivative value at point p for basis function i, component c, in direction d
966: T->T[2] = H[(((p*pdim + i)*Nc + c)*dim + d)*dim + e] is the value at point p for basis
967: T->function i, component c, in directions d and e
968: .ve
970: .seealso: `PetscTabulation`, `PetscFEGetCellTabulation()`, `PetscTabulationDestroy()`
971: @*/
972: PetscErrorCode PetscFECreateTabulation(PetscFE fem, PetscInt nrepl, PetscInt npoints, const PetscReal points[], PetscInt K, PetscTabulation *T)
973: {
974: DM dm;
975: PetscDualSpace Q;
976: PetscInt Nb; /* Dimension of FE space P */
977: PetscInt Nc; /* Field components */
978: PetscInt cdim; /* Reference coordinate dimension */
979: PetscInt k;
981: PetscFunctionBegin;
982: if (!npoints || !fem->dualSpace || K < 0) {
983: *T = NULL;
984: PetscFunctionReturn(PETSC_SUCCESS);
985: }
987: PetscAssertPointer(points, 4);
988: PetscAssertPointer(T, 6);
989: PetscCall(PetscFEGetDualSpace(fem, &Q));
990: PetscCall(PetscDualSpaceGetDM(Q, &dm));
991: PetscCall(DMGetDimension(dm, &cdim));
992: PetscCall(PetscDualSpaceGetDimension(Q, &Nb));
993: PetscCall(PetscFEGetNumComponents(fem, &Nc));
994: PetscCall(PetscMalloc1(1, T));
995: (*T)->K = !cdim ? 0 : K;
996: (*T)->Nr = nrepl;
997: (*T)->Np = npoints;
998: (*T)->Nb = Nb;
999: (*T)->Nc = Nc;
1000: (*T)->cdim = cdim;
1001: PetscCall(PetscMalloc1((*T)->K + 1, &(*T)->T));
1002: for (k = 0; k <= (*T)->K; ++k) PetscCall(PetscCalloc1(nrepl * npoints * Nb * Nc * PetscPowInt(cdim, k), &(*T)->T[k]));
1003: PetscUseTypeMethod(fem, computetabulation, nrepl * npoints, points, K, *T);
1004: PetscFunctionReturn(PETSC_SUCCESS);
1005: }
1007: /*@C
1008: PetscFEComputeTabulation - Tabulates the basis functions, and perhaps derivatives, at the points provided.
1010: Not Collective
1012: Input Parameters:
1013: + fem - The `PetscFE` object
1014: . npoints - The number of tabulation points
1015: . points - The tabulation point coordinates
1016: . K - The number of derivatives calculated
1017: - T - An existing tabulation object with enough allocated space
1019: Output Parameter:
1020: . T - The basis function values and derivatives at tabulation points
1022: Level: intermediate
1024: Note:
1025: .vb
1026: T->T[0] = B[(p*pdim + i)*Nc + c] is the value at point p for basis function i and component c
1027: T->T[1] = D[((p*pdim + i)*Nc + c)*dim + d] is the derivative value at point p for basis function i, component c, in direction d
1028: T->T[2] = H[(((p*pdim + i)*Nc + c)*dim + d)*dim + e] is the value at point p for basis function i, component c, in directions d and e
1029: .ve
1031: .seealso: `PetscTabulation`, `PetscFEGetCellTabulation()`, `PetscTabulationDestroy()`
1032: @*/
1033: PetscErrorCode PetscFEComputeTabulation(PetscFE fem, PetscInt npoints, const PetscReal points[], PetscInt K, PetscTabulation T)
1034: {
1035: PetscFunctionBeginHot;
1036: if (!npoints || !fem->dualSpace || K < 0) PetscFunctionReturn(PETSC_SUCCESS);
1038: PetscAssertPointer(points, 3);
1039: PetscAssertPointer(T, 5);
1040: if (PetscDefined(USE_DEBUG)) {
1041: DM dm;
1042: PetscDualSpace Q;
1043: PetscInt Nb; /* Dimension of FE space P */
1044: PetscInt Nc; /* Field components */
1045: PetscInt cdim; /* Reference coordinate dimension */
1047: PetscCall(PetscFEGetDualSpace(fem, &Q));
1048: PetscCall(PetscDualSpaceGetDM(Q, &dm));
1049: PetscCall(DMGetDimension(dm, &cdim));
1050: PetscCall(PetscDualSpaceGetDimension(Q, &Nb));
1051: PetscCall(PetscFEGetNumComponents(fem, &Nc));
1052: PetscCheck(T->K == (!cdim ? 0 : K), PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Tabulation K %" PetscInt_FMT " must match requested K %" PetscInt_FMT, T->K, !cdim ? 0 : K);
1053: PetscCheck(T->Nb == Nb, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Tabulation Nb %" PetscInt_FMT " must match requested Nb %" PetscInt_FMT, T->Nb, Nb);
1054: PetscCheck(T->Nc == Nc, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Tabulation Nc %" PetscInt_FMT " must match requested Nc %" PetscInt_FMT, T->Nc, Nc);
1055: PetscCheck(T->cdim == cdim, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Tabulation cdim %" PetscInt_FMT " must match requested cdim %" PetscInt_FMT, T->cdim, cdim);
1056: }
1057: T->Nr = 1;
1058: T->Np = npoints;
1059: PetscUseTypeMethod(fem, computetabulation, npoints, points, K, T);
1060: PetscFunctionReturn(PETSC_SUCCESS);
1061: }
1063: /*@
1064: PetscTabulationDestroy - Frees memory from the associated tabulation.
1066: Not Collective
1068: Input Parameter:
1069: . T - The tabulation
1071: Level: intermediate
1073: .seealso: `PetscTabulation`, `PetscFECreateTabulation()`, `PetscFEGetCellTabulation()`
1074: @*/
1075: PetscErrorCode PetscTabulationDestroy(PetscTabulation *T)
1076: {
1077: PetscInt k;
1079: PetscFunctionBegin;
1080: PetscAssertPointer(T, 1);
1081: if (!T || !*T) PetscFunctionReturn(PETSC_SUCCESS);
1082: for (k = 0; k <= (*T)->K; ++k) PetscCall(PetscFree((*T)->T[k]));
1083: PetscCall(PetscFree((*T)->T));
1084: PetscCall(PetscFree(*T));
1085: *T = NULL;
1086: PetscFunctionReturn(PETSC_SUCCESS);
1087: }
1089: static PetscErrorCode PetscFECreatePointTraceDefault_Internal(PetscFE fe, PetscInt refPoint, PetscFE *trFE)
1090: {
1091: PetscSpace bsp, bsubsp;
1092: PetscDualSpace dsp, dsubsp;
1093: PetscInt dim, depth, numComp, i, j, coneSize, order;
1094: DM dm;
1095: DMLabel label;
1096: PetscReal *xi, *v, *J, detJ;
1097: const char *name;
1098: PetscQuadrature origin, fullQuad, subQuad;
1100: PetscFunctionBegin;
1101: PetscCall(PetscFEGetBasisSpace(fe, &bsp));
1102: PetscCall(PetscFEGetDualSpace(fe, &dsp));
1103: PetscCall(PetscDualSpaceGetDM(dsp, &dm));
1104: PetscCall(DMGetDimension(dm, &dim));
1105: PetscCall(DMPlexGetDepthLabel(dm, &label));
1106: PetscCall(DMLabelGetValue(label, refPoint, &depth));
1107: PetscCall(PetscCalloc1(depth, &xi));
1108: PetscCall(PetscMalloc1(dim, &v));
1109: PetscCall(PetscMalloc1(dim * dim, &J));
1110: for (i = 0; i < depth; i++) xi[i] = 0.;
1111: PetscCall(PetscQuadratureCreate(PETSC_COMM_SELF, &origin));
1112: PetscCall(PetscQuadratureSetData(origin, depth, 0, 1, xi, NULL));
1113: PetscCall(DMPlexComputeCellGeometryFEM(dm, refPoint, origin, v, J, NULL, &detJ));
1114: /* CellGeometryFEM computes the expanded Jacobian, we want the true jacobian */
1115: for (i = 1; i < dim; i++) {
1116: for (j = 0; j < depth; j++) J[i * depth + j] = J[i * dim + j];
1117: }
1118: PetscCall(PetscQuadratureDestroy(&origin));
1119: PetscCall(PetscDualSpaceGetPointSubspace(dsp, refPoint, &dsubsp));
1120: PetscCall(PetscSpaceCreateSubspace(bsp, dsubsp, v, J, NULL, NULL, PETSC_OWN_POINTER, &bsubsp));
1121: PetscCall(PetscSpaceSetUp(bsubsp));
1122: PetscCall(PetscFECreate(PetscObjectComm((PetscObject)fe), trFE));
1123: PetscCall(PetscFESetType(*trFE, PETSCFEBASIC));
1124: PetscCall(PetscFEGetNumComponents(fe, &numComp));
1125: PetscCall(PetscFESetNumComponents(*trFE, numComp));
1126: PetscCall(PetscFESetBasisSpace(*trFE, bsubsp));
1127: PetscCall(PetscFESetDualSpace(*trFE, dsubsp));
1128: PetscCall(PetscObjectGetName((PetscObject)fe, &name));
1129: if (name) PetscCall(PetscFESetName(*trFE, name));
1130: PetscCall(PetscFEGetQuadrature(fe, &fullQuad));
1131: PetscCall(PetscQuadratureGetOrder(fullQuad, &order));
1132: PetscCall(DMPlexGetConeSize(dm, refPoint, &coneSize));
1133: if (coneSize == 2 * depth) PetscCall(PetscDTGaussTensorQuadrature(depth, 1, (order + 2) / 2, -1., 1., &subQuad));
1134: else PetscCall(PetscDTSimplexQuadrature(depth, order, PETSCDTSIMPLEXQUAD_DEFAULT, &subQuad));
1135: PetscCall(PetscFESetQuadrature(*trFE, subQuad));
1136: PetscCall(PetscFESetUp(*trFE));
1137: PetscCall(PetscQuadratureDestroy(&subQuad));
1138: PetscCall(PetscSpaceDestroy(&bsubsp));
1139: PetscFunctionReturn(PETSC_SUCCESS);
1140: }
1142: PETSC_EXTERN PetscErrorCode PetscFECreatePointTrace(PetscFE fe, PetscInt refPoint, PetscFE *trFE)
1143: {
1144: PetscFunctionBegin;
1146: PetscAssertPointer(trFE, 3);
1147: if (fe->ops->createpointtrace) PetscUseTypeMethod(fe, createpointtrace, refPoint, trFE);
1148: else PetscCall(PetscFECreatePointTraceDefault_Internal(fe, refPoint, trFE));
1149: PetscFunctionReturn(PETSC_SUCCESS);
1150: }
1152: PetscErrorCode PetscFECreateHeightTrace(PetscFE fe, PetscInt height, PetscFE *trFE)
1153: {
1154: PetscInt hStart, hEnd;
1155: PetscDualSpace dsp;
1156: DM dm;
1158: PetscFunctionBegin;
1160: PetscAssertPointer(trFE, 3);
1161: *trFE = NULL;
1162: PetscCall(PetscFEGetDualSpace(fe, &dsp));
1163: PetscCall(PetscDualSpaceGetDM(dsp, &dm));
1164: PetscCall(DMPlexGetHeightStratum(dm, height, &hStart, &hEnd));
1165: if (hEnd <= hStart) PetscFunctionReturn(PETSC_SUCCESS);
1166: PetscCall(PetscFECreatePointTrace(fe, hStart, trFE));
1167: PetscFunctionReturn(PETSC_SUCCESS);
1168: }
1170: /*@
1171: PetscFEGetDimension - Get the dimension of the finite element space on a cell
1173: Not Collective
1175: Input Parameter:
1176: . fem - The `PetscFE`
1178: Output Parameter:
1179: . dim - The dimension
1181: Level: intermediate
1183: .seealso: `PetscFE`, `PetscFECreate()`, `PetscSpaceGetDimension()`, `PetscDualSpaceGetDimension()`
1184: @*/
1185: PetscErrorCode PetscFEGetDimension(PetscFE fem, PetscInt *dim)
1186: {
1187: PetscFunctionBegin;
1189: PetscAssertPointer(dim, 2);
1190: PetscTryTypeMethod(fem, getdimension, dim);
1191: PetscFunctionReturn(PETSC_SUCCESS);
1192: }
1194: /*@
1195: PetscFEPushforward - Map the reference element function to real space
1197: Input Parameters:
1198: + fe - The `PetscFE`
1199: . fegeom - The cell geometry
1200: . Nv - The number of function values
1201: - vals - The function values
1203: Output Parameter:
1204: . vals - The transformed function values
1206: Level: advanced
1208: Notes:
1209: This just forwards the call onto `PetscDualSpacePushforward()`.
1211: It only handles transformations when the embedding dimension of the geometry in fegeom is the same as the reference dimension.
1213: .seealso: `PetscFE`, `PetscFEGeom`, `PetscDualSpace`, `PetscDualSpacePushforward()`
1214: @*/
1215: PetscErrorCode PetscFEPushforward(PetscFE fe, PetscFEGeom *fegeom, PetscInt Nv, PetscScalar vals[])
1216: {
1217: PetscFunctionBeginHot;
1218: PetscCall(PetscDualSpacePushforward(fe->dualSpace, fegeom, Nv, fe->numComponents, vals));
1219: PetscFunctionReturn(PETSC_SUCCESS);
1220: }
1222: /*@
1223: PetscFEPushforwardGradient - Map the reference element function gradient to real space
1225: Input Parameters:
1226: + fe - The `PetscFE`
1227: . fegeom - The cell geometry
1228: . Nv - The number of function gradient values
1229: - vals - The function gradient values
1231: Output Parameter:
1232: . vals - The transformed function gradient values
1234: Level: advanced
1236: Notes:
1237: This just forwards the call onto `PetscDualSpacePushforwardGradient()`.
1239: It only handles transformations when the embedding dimension of the geometry in fegeom is the same as the reference dimension.
1241: .seealso: `PetscFE`, `PetscFEGeom`, `PetscDualSpace`, `PetscFEPushforward()`, `PetscDualSpacePushforwardGradient()`, `PetscDualSpacePushforward()`
1242: @*/
1243: PetscErrorCode PetscFEPushforwardGradient(PetscFE fe, PetscFEGeom *fegeom, PetscInt Nv, PetscScalar vals[])
1244: {
1245: PetscFunctionBeginHot;
1246: PetscCall(PetscDualSpacePushforwardGradient(fe->dualSpace, fegeom, Nv, fe->numComponents, vals));
1247: PetscFunctionReturn(PETSC_SUCCESS);
1248: }
1250: /*@
1251: PetscFEPushforwardHessian - Map the reference element function Hessian to real space
1253: Input Parameters:
1254: + fe - The `PetscFE`
1255: . fegeom - The cell geometry
1256: . Nv - The number of function Hessian values
1257: - vals - The function Hessian values
1259: Output Parameter:
1260: . vals - The transformed function Hessian values
1262: Level: advanced
1264: Notes:
1265: This just forwards the call onto `PetscDualSpacePushforwardHessian()`.
1267: It only handles transformations when the embedding dimension of the geometry in fegeom is the same as the reference dimension.
1269: Developer Notes:
1270: It is unclear why all these one line convenience routines are desirable
1272: .seealso: `PetscFE`, `PetscFEGeom`, `PetscDualSpace`, `PetscFEPushforward()`, `PetscDualSpacePushforwardHessian()`, `PetscDualSpacePushforward()`
1273: @*/
1274: PetscErrorCode PetscFEPushforwardHessian(PetscFE fe, PetscFEGeom *fegeom, PetscInt Nv, PetscScalar vals[])
1275: {
1276: PetscFunctionBeginHot;
1277: PetscCall(PetscDualSpacePushforwardHessian(fe->dualSpace, fegeom, Nv, fe->numComponents, vals));
1278: PetscFunctionReturn(PETSC_SUCCESS);
1279: }
1281: /*
1282: Purpose: Compute element vector for chunk of elements
1284: Input:
1285: Sizes:
1286: Ne: number of elements
1287: Nf: number of fields
1288: PetscFE
1289: dim: spatial dimension
1290: Nb: number of basis functions
1291: Nc: number of field components
1292: PetscQuadrature
1293: Nq: number of quadrature points
1295: Geometry:
1296: PetscFEGeom[Ne] possibly *Nq
1297: PetscReal v0s[dim]
1298: PetscReal n[dim]
1299: PetscReal jacobians[dim*dim]
1300: PetscReal jacobianInverses[dim*dim]
1301: PetscReal jacobianDeterminants
1302: FEM:
1303: PetscFE
1304: PetscQuadrature
1305: PetscReal quadPoints[Nq*dim]
1306: PetscReal quadWeights[Nq]
1307: PetscReal basis[Nq*Nb*Nc]
1308: PetscReal basisDer[Nq*Nb*Nc*dim]
1309: PetscScalar coefficients[Ne*Nb*Nc]
1310: PetscScalar elemVec[Ne*Nb*Nc]
1312: Problem:
1313: PetscInt f: the active field
1314: f0, f1
1316: Work Space:
1317: PetscFE
1318: PetscScalar f0[Nq*dim];
1319: PetscScalar f1[Nq*dim*dim];
1320: PetscScalar u[Nc];
1321: PetscScalar gradU[Nc*dim];
1322: PetscReal x[dim];
1323: PetscScalar realSpaceDer[dim];
1325: Purpose: Compute element vector for N_cb batches of elements
1327: Input:
1328: Sizes:
1329: N_cb: Number of serial cell batches
1331: Geometry:
1332: PetscReal v0s[Ne*dim]
1333: PetscReal jacobians[Ne*dim*dim] possibly *Nq
1334: PetscReal jacobianInverses[Ne*dim*dim] possibly *Nq
1335: PetscReal jacobianDeterminants[Ne] possibly *Nq
1336: FEM:
1337: static PetscReal quadPoints[Nq*dim]
1338: static PetscReal quadWeights[Nq]
1339: static PetscReal basis[Nq*Nb*Nc]
1340: static PetscReal basisDer[Nq*Nb*Nc*dim]
1341: PetscScalar coefficients[Ne*Nb*Nc]
1342: PetscScalar elemVec[Ne*Nb*Nc]
1344: ex62.c:
1345: PetscErrorCode PetscFEIntegrateResidualBatch(PetscInt Ne, PetscInt numFields, PetscInt field, PetscQuadrature quad[], const PetscScalar coefficients[],
1346: const PetscReal v0s[], const PetscReal jacobians[], const PetscReal jacobianInverses[], const PetscReal jacobianDeterminants[],
1347: void (*f0_func)(const PetscScalar u[], const PetscScalar gradU[], const PetscReal x[], PetscScalar f0[]),
1348: void (*f1_func)(const PetscScalar u[], const PetscScalar gradU[], const PetscReal x[], PetscScalar f1[]), PetscScalar elemVec[])
1350: ex52.c:
1351: PetscErrorCode IntegrateLaplacianBatchCPU(PetscInt Ne, PetscInt Nb, const PetscScalar coefficients[], const PetscReal jacobianInverses[], const PetscReal jacobianDeterminants[], PetscInt Nq, const PetscReal quadPoints[], const PetscReal quadWeights[], const PetscReal basisTabulation[], const PetscReal basisDerTabulation[], PetscScalar elemVec[], AppCtx *user)
1352: PetscErrorCode IntegrateElasticityBatchCPU(PetscInt Ne, PetscInt Nb, PetscInt Ncomp, const PetscScalar coefficients[], const PetscReal jacobianInverses[], const PetscReal jacobianDeterminants[], PetscInt Nq, const PetscReal quadPoints[], const PetscReal quadWeights[], const PetscReal basisTabulation[], const PetscReal basisDerTabulation[], PetscScalar elemVec[], AppCtx *user)
1354: ex52_integrateElement.cu
1355: __global__ void integrateElementQuadrature(int N_cb, realType *coefficients, realType *jacobianInverses, realType *jacobianDeterminants, realType *elemVec)
1357: PETSC_EXTERN PetscErrorCode IntegrateElementBatchGPU(PetscInt spatial_dim, PetscInt Ne, PetscInt Ncb, PetscInt Nbc, PetscInt Nbl, const PetscScalar coefficients[],
1358: const PetscReal jacobianInverses[], const PetscReal jacobianDeterminants[], PetscScalar elemVec[],
1359: PetscLogEvent event, PetscInt debug, PetscInt pde_op)
1361: ex52_integrateElementOpenCL.c:
1362: PETSC_EXTERN PetscErrorCode IntegrateElementBatchGPU(PetscInt spatial_dim, PetscInt Ne, PetscInt Ncb, PetscInt Nbc, PetscInt N_bl, const PetscScalar coefficients[],
1363: const PetscReal jacobianInverses[], const PetscReal jacobianDeterminants[], PetscScalar elemVec[],
1364: PetscLogEvent event, PetscInt debug, PetscInt pde_op)
1366: __kernel void integrateElementQuadrature(int N_cb, __global float *coefficients, __global float *jacobianInverses, __global float *jacobianDeterminants, __global float *elemVec)
1367: */
1369: /*@
1370: PetscFEIntegrate - Produce the integral for the given field for a chunk of elements by quadrature integration
1372: Not Collective
1374: Input Parameters:
1375: + prob - The `PetscDS` specifying the discretizations and continuum functions
1376: . field - The field being integrated
1377: . Ne - The number of elements in the chunk
1378: . cgeom - The cell geometry for each cell in the chunk
1379: . coefficients - The array of FEM basis coefficients for the elements
1380: . probAux - The `PetscDS` specifying the auxiliary discretizations
1381: - coefficientsAux - The array of FEM auxiliary basis coefficients for the elements
1383: Output Parameter:
1384: . integral - the integral for this field
1386: Level: intermediate
1388: Developer Notes:
1389: The function name begins with `PetscFE` and yet the first argument is `PetscDS` and it has no `PetscFE` arguments.
1391: .seealso: `PetscFE`, `PetscDS`, `PetscFEIntegrateResidual()`, `PetscFEIntegrateBd()`
1392: @*/
1393: PetscErrorCode PetscFEIntegrate(PetscDS prob, PetscInt field, PetscInt Ne, PetscFEGeom *cgeom, const PetscScalar coefficients[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscScalar integral[])
1394: {
1395: PetscFE fe;
1397: PetscFunctionBegin;
1399: PetscCall(PetscDSGetDiscretization(prob, field, (PetscObject *)&fe));
1400: if (fe->ops->integrate) PetscCall((*fe->ops->integrate)(prob, field, Ne, cgeom, coefficients, probAux, coefficientsAux, integral));
1401: PetscFunctionReturn(PETSC_SUCCESS);
1402: }
1404: /*@C
1405: PetscFEIntegrateBd - Produce the integral for the given field for a chunk of elements by quadrature integration
1407: Not Collective
1409: Input Parameters:
1410: + prob - The `PetscDS` specifying the discretizations and continuum functions
1411: . field - The field being integrated
1412: . obj_func - The function to be integrated
1413: . Ne - The number of elements in the chunk
1414: . geom - The face geometry for each face in the chunk
1415: . coefficients - The array of FEM basis coefficients for the elements
1416: . probAux - The `PetscDS` specifying the auxiliary discretizations
1417: - coefficientsAux - The array of FEM auxiliary basis coefficients for the elements
1419: Output Parameter:
1420: . integral - the integral for this field
1422: Level: intermediate
1424: Developer Notes:
1425: The function name begins with `PetscFE` and yet the first argument is `PetscDS` and it has no `PetscFE` arguments.
1427: .seealso: `PetscFE`, `PetscDS`, `PetscFEIntegrateResidual()`, `PetscFEIntegrate()`
1428: @*/
1429: PetscErrorCode PetscFEIntegrateBd(PetscDS prob, PetscInt field, void (*obj_func)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt Ne, PetscFEGeom *geom, const PetscScalar coefficients[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscScalar integral[])
1430: {
1431: PetscFE fe;
1433: PetscFunctionBegin;
1435: PetscCall(PetscDSGetDiscretization(prob, field, (PetscObject *)&fe));
1436: if (fe->ops->integratebd) PetscCall((*fe->ops->integratebd)(prob, field, obj_func, Ne, geom, coefficients, probAux, coefficientsAux, integral));
1437: PetscFunctionReturn(PETSC_SUCCESS);
1438: }
1440: /*@
1441: PetscFEIntegrateResidual - Produce the element residual vector for a chunk of elements by quadrature integration
1443: Not Collective
1445: Input Parameters:
1446: + ds - The `PetscDS` specifying the discretizations and continuum functions
1447: . key - The (label+value, field) being integrated
1448: . Ne - The number of elements in the chunk
1449: . cgeom - The cell geometry for each cell in the chunk
1450: . coefficients - The array of FEM basis coefficients for the elements
1451: . coefficients_t - The array of FEM basis time derivative coefficients for the elements
1452: . probAux - The `PetscDS` specifying the auxiliary discretizations
1453: . coefficientsAux - The array of FEM auxiliary basis coefficients for the elements
1454: - t - The time
1456: Output Parameter:
1457: . elemVec - the element residual vectors from each element
1459: Level: intermediate
1461: Note:
1462: .vb
1463: Loop over batch of elements (e):
1464: Loop over quadrature points (q):
1465: Make u_q and gradU_q (loops over fields,Nb,Ncomp) and x_q
1466: Call f_0 and f_1
1467: Loop over element vector entries (f,fc --> i):
1468: elemVec[i] += \psi^{fc}_f(q) f0_{fc}(u, \nabla u) + \nabla\psi^{fc}_f(q) \cdot f1_{fc,df}(u, \nabla u)
1469: .ve
1471: .seealso: `PetscFEIntegrateBdResidual()`
1472: @*/
1473: PetscErrorCode PetscFEIntegrateResidual(PetscDS ds, PetscFormKey key, PetscInt Ne, PetscFEGeom *cgeom, const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscReal t, PetscScalar elemVec[])
1474: {
1475: PetscFE fe;
1477: PetscFunctionBeginHot;
1479: PetscCall(PetscDSGetDiscretization(ds, key.field, (PetscObject *)&fe));
1480: if (fe->ops->integrateresidual) PetscCall((*fe->ops->integrateresidual)(ds, key, Ne, cgeom, coefficients, coefficients_t, probAux, coefficientsAux, t, elemVec));
1481: PetscFunctionReturn(PETSC_SUCCESS);
1482: }
1484: /*@
1485: PetscFEIntegrateBdResidual - Produce the element residual vector for a chunk of elements by quadrature integration over a boundary
1487: Not Collective
1489: Input Parameters:
1490: + ds - The `PetscDS` specifying the discretizations and continuum functions
1491: . wf - The PetscWeakForm object holding the pointwise functions
1492: . key - The (label+value, field) being integrated
1493: . Ne - The number of elements in the chunk
1494: . fgeom - The face geometry for each cell in the chunk
1495: . coefficients - The array of FEM basis coefficients for the elements
1496: . coefficients_t - The array of FEM basis time derivative coefficients for the elements
1497: . probAux - The `PetscDS` specifying the auxiliary discretizations
1498: . coefficientsAux - The array of FEM auxiliary basis coefficients for the elements
1499: - t - The time
1501: Output Parameter:
1502: . elemVec - the element residual vectors from each element
1504: Level: intermediate
1506: .seealso: `PetscFEIntegrateResidual()`
1507: @*/
1508: PetscErrorCode PetscFEIntegrateBdResidual(PetscDS ds, PetscWeakForm wf, PetscFormKey key, PetscInt Ne, PetscFEGeom *fgeom, const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscReal t, PetscScalar elemVec[])
1509: {
1510: PetscFE fe;
1512: PetscFunctionBegin;
1514: PetscCall(PetscDSGetDiscretization(ds, key.field, (PetscObject *)&fe));
1515: if (fe->ops->integratebdresidual) PetscCall((*fe->ops->integratebdresidual)(ds, wf, key, Ne, fgeom, coefficients, coefficients_t, probAux, coefficientsAux, t, elemVec));
1516: PetscFunctionReturn(PETSC_SUCCESS);
1517: }
1519: /*@
1520: PetscFEIntegrateHybridResidual - Produce the element residual vector for a chunk of hybrid element faces by quadrature integration
1522: Not Collective
1524: Input Parameters:
1525: + ds - The `PetscDS` specifying the discretizations and continuum functions
1526: . dsIn - The `PetscDS` specifying the discretizations and continuum functions for input
1527: . key - The (label+value, field) being integrated
1528: . s - The side of the cell being integrated, 0 for negative and 1 for positive
1529: . Ne - The number of elements in the chunk
1530: . fgeom - The face geometry for each cell in the chunk
1531: . cgeom - The cell geometry for each neighbor cell in the chunk
1532: . coefficients - The array of FEM basis coefficients for the elements
1533: . coefficients_t - The array of FEM basis time derivative coefficients for the elements
1534: . probAux - The `PetscDS` specifying the auxiliary discretizations
1535: . coefficientsAux - The array of FEM auxiliary basis coefficients for the elements
1536: - t - The time
1538: Output Parameter:
1539: . elemVec - the element residual vectors from each element
1541: Level: developer
1543: .seealso: `PetscFEIntegrateResidual()`
1544: @*/
1545: PetscErrorCode PetscFEIntegrateHybridResidual(PetscDS ds, PetscDS dsIn, PetscFormKey key, PetscInt s, PetscInt Ne, PetscFEGeom *fgeom, PetscFEGeom *cgeom, const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscReal t, PetscScalar elemVec[])
1546: {
1547: PetscFE fe;
1549: PetscFunctionBegin;
1552: PetscCall(PetscDSGetDiscretization(ds, key.field, (PetscObject *)&fe));
1553: if (fe->ops->integratehybridresidual) PetscCall((*fe->ops->integratehybridresidual)(ds, dsIn, key, s, Ne, fgeom, cgeom, coefficients, coefficients_t, probAux, coefficientsAux, t, elemVec));
1554: PetscFunctionReturn(PETSC_SUCCESS);
1555: }
1557: /*@
1558: PetscFEIntegrateJacobian - Produce the element Jacobian for a chunk of elements by quadrature integration
1560: Not Collective
1562: Input Parameters:
1563: + rds - The `PetscDS` specifying the row discretizations and continuum functions
1564: . cds - The `PetscDS` specifying the column discretizations
1565: . jtype - The type of matrix pointwise functions that should be used
1566: . key - The (label+value, fieldI*Nf + fieldJ) being integrated
1567: . Ne - The number of elements in the chunk
1568: . cgeom - The cell geometry for each cell in the chunk
1569: . coefficients - The array of FEM basis coefficients for the elements for the Jacobian evaluation point
1570: . coefficients_t - The array of FEM basis time derivative coefficients for the elements
1571: . dsAux - The `PetscDS` specifying the auxiliary discretizations
1572: . coefficientsAux - The array of FEM auxiliary basis coefficients for the elements
1573: . t - The time
1574: - u_tshift - A multiplier for the dF/du_t term (as opposed to the dF/du term)
1576: Output Parameter:
1577: . elemMat - the element matrices for the Jacobian from each element
1579: Level: intermediate
1581: Note:
1582: .vb
1583: Loop over batch of elements (e):
1584: Loop over element matrix entries (f,fc,g,gc --> i,j):
1585: Loop over quadrature points (q):
1586: Make u_q and gradU_q (loops over fields,Nb,Ncomp)
1587: elemMat[i,j] += \psi^{fc}_f(q) g0_{fc,gc}(u, \nabla u) \phi^{gc}_g(q)
1588: + \psi^{fc}_f(q) \cdot g1_{fc,gc,dg}(u, \nabla u) \nabla\phi^{gc}_g(q)
1589: + \nabla\psi^{fc}_f(q) \cdot g2_{fc,gc,df}(u, \nabla u) \phi^{gc}_g(q)
1590: + \nabla\psi^{fc}_f(q) \cdot g3_{fc,gc,df,dg}(u, \nabla u) \nabla\phi^{gc}_g(q)
1591: .ve
1593: .seealso: `PetscFEIntegrateResidual()`
1594: @*/
1595: PetscErrorCode PetscFEIntegrateJacobian(PetscDS rds, PetscDS cds, PetscFEJacobianType jtype, PetscFormKey key, PetscInt Ne, PetscFEGeom *cgeom, const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscReal t, PetscReal u_tshift, PetscScalar elemMat[])
1596: {
1597: PetscFE fe;
1598: PetscInt Nf;
1600: PetscFunctionBegin;
1603: PetscCall(PetscDSGetNumFields(rds, &Nf));
1604: PetscCall(PetscDSGetDiscretization(rds, key.field / Nf, (PetscObject *)&fe));
1605: if (fe->ops->integratejacobian) PetscCall((*fe->ops->integratejacobian)(rds, cds, jtype, key, Ne, cgeom, coefficients, coefficients_t, dsAux, coefficientsAux, t, u_tshift, elemMat));
1606: PetscFunctionReturn(PETSC_SUCCESS);
1607: }
1609: /*@
1610: PetscFEIntegrateBdJacobian - Produce the boundary element Jacobian for a chunk of elements by quadrature integration
1612: Not Collective
1614: Input Parameters:
1615: + ds - The `PetscDS` specifying the discretizations and continuum functions
1616: . wf - The PetscWeakForm holding the pointwise functions
1617: . jtype - The type of matrix pointwise functions that should be used
1618: . key - The (label+value, fieldI*Nf + fieldJ) being integrated
1619: . Ne - The number of elements in the chunk
1620: . fgeom - The face geometry for each cell in the chunk
1621: . coefficients - The array of FEM basis coefficients for the elements for the Jacobian evaluation point
1622: . coefficients_t - The array of FEM basis time derivative coefficients for the elements
1623: . probAux - The `PetscDS` specifying the auxiliary discretizations
1624: . coefficientsAux - The array of FEM auxiliary basis coefficients for the elements
1625: . t - The time
1626: - u_tshift - A multiplier for the dF/du_t term (as opposed to the dF/du term)
1628: Output Parameter:
1629: . elemMat - the element matrices for the Jacobian from each element
1631: Level: intermediate
1633: Note:
1634: .vb
1635: Loop over batch of elements (e):
1636: Loop over element matrix entries (f,fc,g,gc --> i,j):
1637: Loop over quadrature points (q):
1638: Make u_q and gradU_q (loops over fields,Nb,Ncomp)
1639: elemMat[i,j] += \psi^{fc}_f(q) g0_{fc,gc}(u, \nabla u) \phi^{gc}_g(q)
1640: + \psi^{fc}_f(q) \cdot g1_{fc,gc,dg}(u, \nabla u) \nabla\phi^{gc}_g(q)
1641: + \nabla\psi^{fc}_f(q) \cdot g2_{fc,gc,df}(u, \nabla u) \phi^{gc}_g(q)
1642: + \nabla\psi^{fc}_f(q) \cdot g3_{fc,gc,df,dg}(u, \nabla u) \nabla\phi^{gc}_g(q)
1643: .ve
1645: .seealso: `PetscFEIntegrateJacobian()`, `PetscFEIntegrateResidual()`
1646: @*/
1647: PetscErrorCode PetscFEIntegrateBdJacobian(PetscDS ds, PetscWeakForm wf, PetscFEJacobianType jtype, PetscFormKey key, PetscInt Ne, PetscFEGeom *fgeom, const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscReal t, PetscReal u_tshift, PetscScalar elemMat[])
1648: {
1649: PetscFE fe;
1650: PetscInt Nf;
1652: PetscFunctionBegin;
1654: PetscCall(PetscDSGetNumFields(ds, &Nf));
1655: PetscCall(PetscDSGetDiscretization(ds, key.field / Nf, (PetscObject *)&fe));
1656: if (fe->ops->integratebdjacobian) PetscCall((*fe->ops->integratebdjacobian)(ds, wf, jtype, key, Ne, fgeom, coefficients, coefficients_t, probAux, coefficientsAux, t, u_tshift, elemMat));
1657: PetscFunctionReturn(PETSC_SUCCESS);
1658: }
1660: /*@
1661: PetscFEIntegrateHybridJacobian - Produce the boundary element Jacobian for a chunk of hybrid elements by quadrature integration
1663: Not Collective
1665: Input Parameters:
1666: + ds - The `PetscDS` specifying the discretizations and continuum functions for the output
1667: . dsIn - The `PetscDS` specifying the discretizations and continuum functions for the input
1668: . jtype - The type of matrix pointwise functions that should be used
1669: . key - The (label+value, fieldI*Nf + fieldJ) being integrated
1670: . s - The side of the cell being integrated, 0 for negative and 1 for positive
1671: . Ne - The number of elements in the chunk
1672: . fgeom - The face geometry for each cell in the chunk
1673: . cgeom - The cell geometry for each neighbor cell in the chunk
1674: . coefficients - The array of FEM basis coefficients for the elements for the Jacobian evaluation point
1675: . coefficients_t - The array of FEM basis time derivative coefficients for the elements
1676: . probAux - The `PetscDS` specifying the auxiliary discretizations
1677: . coefficientsAux - The array of FEM auxiliary basis coefficients for the elements
1678: . t - The time
1679: - u_tshift - A multiplier for the dF/du_t term (as opposed to the dF/du term)
1681: Output Parameter:
1682: . elemMat - the element matrices for the Jacobian from each element
1684: Level: developer
1686: Note:
1687: .vb
1688: Loop over batch of elements (e):
1689: Loop over element matrix entries (f,fc,g,gc --> i,j):
1690: Loop over quadrature points (q):
1691: Make u_q and gradU_q (loops over fields,Nb,Ncomp)
1692: elemMat[i,j] += \psi^{fc}_f(q) g0_{fc,gc}(u, \nabla u) \phi^{gc}_g(q)
1693: + \psi^{fc}_f(q) \cdot g1_{fc,gc,dg}(u, \nabla u) \nabla\phi^{gc}_g(q)
1694: + \nabla\psi^{fc}_f(q) \cdot g2_{fc,gc,df}(u, \nabla u) \phi^{gc}_g(q)
1695: + \nabla\psi^{fc}_f(q) \cdot g3_{fc,gc,df,dg}(u, \nabla u) \nabla\phi^{gc}_g(q)
1696: .ve
1698: .seealso: `PetscFEIntegrateJacobian()`, `PetscFEIntegrateResidual()`
1699: @*/
1700: PetscErrorCode PetscFEIntegrateHybridJacobian(PetscDS ds, PetscDS dsIn, PetscFEJacobianType jtype, PetscFormKey key, PetscInt s, PetscInt Ne, PetscFEGeom *fgeom, PetscFEGeom *cgeom, const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscReal t, PetscReal u_tshift, PetscScalar elemMat[])
1701: {
1702: PetscFE fe;
1703: PetscInt Nf;
1705: PetscFunctionBegin;
1707: PetscCall(PetscDSGetNumFields(ds, &Nf));
1708: PetscCall(PetscDSGetDiscretization(ds, key.field / Nf, (PetscObject *)&fe));
1709: if (fe->ops->integratehybridjacobian) PetscCall((*fe->ops->integratehybridjacobian)(ds, dsIn, jtype, key, s, Ne, fgeom, cgeom, coefficients, coefficients_t, probAux, coefficientsAux, t, u_tshift, elemMat));
1710: PetscFunctionReturn(PETSC_SUCCESS);
1711: }
1713: /*@
1714: PetscFEGetHeightSubspace - Get the subspace of this space for a mesh point of a given height
1716: Input Parameters:
1717: + fe - The finite element space
1718: - height - The height of the `DMPLEX` point
1720: Output Parameter:
1721: . subfe - The subspace of this `PetscFE` space
1723: Level: advanced
1725: Note:
1726: For example, if we want the subspace of this space for a face, we would choose height = 1.
1728: .seealso: `PetscFECreateDefault()`
1729: @*/
1730: PetscErrorCode PetscFEGetHeightSubspace(PetscFE fe, PetscInt height, PetscFE *subfe)
1731: {
1732: PetscSpace P, subP;
1733: PetscDualSpace Q, subQ;
1734: PetscQuadrature subq;
1735: PetscInt dim, Nc;
1737: PetscFunctionBegin;
1739: PetscAssertPointer(subfe, 3);
1740: if (height == 0) {
1741: *subfe = fe;
1742: PetscFunctionReturn(PETSC_SUCCESS);
1743: }
1744: PetscCall(PetscFEGetBasisSpace(fe, &P));
1745: PetscCall(PetscFEGetDualSpace(fe, &Q));
1746: PetscCall(PetscFEGetNumComponents(fe, &Nc));
1747: PetscCall(PetscFEGetFaceQuadrature(fe, &subq));
1748: PetscCall(PetscDualSpaceGetDimension(Q, &dim));
1749: PetscCheck(height <= dim && height >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Asked for space at height %" PetscInt_FMT " for dimension %" PetscInt_FMT " space", height, dim);
1750: if (!fe->subspaces) PetscCall(PetscCalloc1(dim, &fe->subspaces));
1751: if (height <= dim) {
1752: if (!fe->subspaces[height - 1]) {
1753: PetscFE sub = NULL;
1754: const char *name;
1756: PetscCall(PetscSpaceGetHeightSubspace(P, height, &subP));
1757: PetscCall(PetscDualSpaceGetHeightSubspace(Q, height, &subQ));
1758: if (subQ) {
1759: PetscCall(PetscObjectReference((PetscObject)subP));
1760: PetscCall(PetscObjectReference((PetscObject)subQ));
1761: PetscCall(PetscObjectReference((PetscObject)subq));
1762: PetscCall(PetscFECreateFromSpaces(subP, subQ, subq, NULL, &sub));
1763: }
1764: if (sub) {
1765: PetscCall(PetscObjectGetName((PetscObject)fe, &name));
1766: if (name) PetscCall(PetscFESetName(sub, name));
1767: }
1768: fe->subspaces[height - 1] = sub;
1769: }
1770: *subfe = fe->subspaces[height - 1];
1771: } else {
1772: *subfe = NULL;
1773: }
1774: PetscFunctionReturn(PETSC_SUCCESS);
1775: }
1777: /*@
1778: PetscFERefine - Create a "refined" `PetscFE` object that refines the reference cell into
1779: smaller copies.
1781: Collective
1783: Input Parameter:
1784: . fe - The initial `PetscFE`
1786: Output Parameter:
1787: . feRef - The refined `PetscFE`
1789: Level: advanced
1791: Notes:
1792: This is typically used to generate a preconditioner for a higher order method from a lower order method on a
1793: refined mesh having the same number of dofs (but more sparsity). It is also used to create an
1794: interpolation between regularly refined meshes.
1796: .seealso: `PetscFEType`, `PetscFECreate()`, `PetscFESetType()`
1797: @*/
1798: PetscErrorCode PetscFERefine(PetscFE fe, PetscFE *feRef)
1799: {
1800: PetscSpace P, Pref;
1801: PetscDualSpace Q, Qref;
1802: DM K, Kref;
1803: PetscQuadrature q, qref;
1804: const PetscReal *v0, *jac;
1805: PetscInt numComp, numSubelements;
1806: PetscInt cStart, cEnd, c;
1807: PetscDualSpace *cellSpaces;
1809: PetscFunctionBegin;
1810: PetscCall(PetscFEGetBasisSpace(fe, &P));
1811: PetscCall(PetscFEGetDualSpace(fe, &Q));
1812: PetscCall(PetscFEGetQuadrature(fe, &q));
1813: PetscCall(PetscDualSpaceGetDM(Q, &K));
1814: /* Create space */
1815: PetscCall(PetscObjectReference((PetscObject)P));
1816: Pref = P;
1817: /* Create dual space */
1818: PetscCall(PetscDualSpaceDuplicate(Q, &Qref));
1819: PetscCall(PetscDualSpaceSetType(Qref, PETSCDUALSPACEREFINED));
1820: PetscCall(DMRefine(K, PetscObjectComm((PetscObject)fe), &Kref));
1821: PetscCall(DMGetCoordinatesLocalSetUp(Kref));
1822: PetscCall(PetscDualSpaceSetDM(Qref, Kref));
1823: PetscCall(DMPlexGetHeightStratum(Kref, 0, &cStart, &cEnd));
1824: PetscCall(PetscMalloc1(cEnd - cStart, &cellSpaces));
1825: /* TODO: fix for non-uniform refinement */
1826: for (c = 0; c < cEnd - cStart; c++) cellSpaces[c] = Q;
1827: PetscCall(PetscDualSpaceRefinedSetCellSpaces(Qref, cellSpaces));
1828: PetscCall(PetscFree(cellSpaces));
1829: PetscCall(DMDestroy(&Kref));
1830: PetscCall(PetscDualSpaceSetUp(Qref));
1831: /* Create element */
1832: PetscCall(PetscFECreate(PetscObjectComm((PetscObject)fe), feRef));
1833: PetscCall(PetscFESetType(*feRef, PETSCFECOMPOSITE));
1834: PetscCall(PetscFESetBasisSpace(*feRef, Pref));
1835: PetscCall(PetscFESetDualSpace(*feRef, Qref));
1836: PetscCall(PetscFEGetNumComponents(fe, &numComp));
1837: PetscCall(PetscFESetNumComponents(*feRef, numComp));
1838: PetscCall(PetscFESetUp(*feRef));
1839: PetscCall(PetscSpaceDestroy(&Pref));
1840: PetscCall(PetscDualSpaceDestroy(&Qref));
1841: /* Create quadrature */
1842: PetscCall(PetscFECompositeGetMapping(*feRef, &numSubelements, &v0, &jac, NULL));
1843: PetscCall(PetscQuadratureExpandComposite(q, numSubelements, v0, jac, &qref));
1844: PetscCall(PetscFESetQuadrature(*feRef, qref));
1845: PetscCall(PetscQuadratureDestroy(&qref));
1846: PetscFunctionReturn(PETSC_SUCCESS);
1847: }
1849: static PetscErrorCode PetscFESetDefaultName_Private(PetscFE fe)
1850: {
1851: PetscSpace P;
1852: PetscDualSpace Q;
1853: DM K;
1854: DMPolytopeType ct;
1855: PetscInt degree;
1856: char name[64];
1858: PetscFunctionBegin;
1859: PetscCall(PetscFEGetBasisSpace(fe, &P));
1860: PetscCall(PetscSpaceGetDegree(P, °ree, NULL));
1861: PetscCall(PetscFEGetDualSpace(fe, &Q));
1862: PetscCall(PetscDualSpaceGetDM(Q, &K));
1863: PetscCall(DMPlexGetCellType(K, 0, &ct));
1864: switch (ct) {
1865: case DM_POLYTOPE_SEGMENT:
1866: case DM_POLYTOPE_POINT_PRISM_TENSOR:
1867: case DM_POLYTOPE_QUADRILATERAL:
1868: case DM_POLYTOPE_SEG_PRISM_TENSOR:
1869: case DM_POLYTOPE_HEXAHEDRON:
1870: case DM_POLYTOPE_QUAD_PRISM_TENSOR:
1871: PetscCall(PetscSNPrintf(name, sizeof(name), "Q%" PetscInt_FMT, degree));
1872: break;
1873: case DM_POLYTOPE_TRIANGLE:
1874: case DM_POLYTOPE_TETRAHEDRON:
1875: PetscCall(PetscSNPrintf(name, sizeof(name), "P%" PetscInt_FMT, degree));
1876: break;
1877: case DM_POLYTOPE_TRI_PRISM:
1878: case DM_POLYTOPE_TRI_PRISM_TENSOR:
1879: PetscCall(PetscSNPrintf(name, sizeof(name), "P%" PetscInt_FMT "xQ%" PetscInt_FMT, degree, degree));
1880: break;
1881: default:
1882: PetscCall(PetscSNPrintf(name, sizeof(name), "FE"));
1883: }
1884: PetscCall(PetscFESetName(fe, name));
1885: PetscFunctionReturn(PETSC_SUCCESS);
1886: }
1888: /*@
1889: PetscFECreateFromSpaces - Create a `PetscFE` from the basis and dual spaces
1891: Collective
1893: Input Parameters:
1894: + P - The basis space
1895: . Q - The dual space
1896: . q - The cell quadrature
1897: - fq - The face quadrature
1899: Output Parameter:
1900: . fem - The `PetscFE` object
1902: Level: beginner
1904: Note:
1905: The `PetscFE` takes ownership of these spaces by calling destroy on each. They should not be used after this call, and for borrowed references from `PetscFEGetSpace()` and the like, the caller must use `PetscObjectReference` before this call.
1907: .seealso: `PetscFE`, `PetscSpace`, `PetscDualSpace`, `PetscQuadrature`,
1908: `PetscFECreateLagrangeByCell()`, `PetscFECreateDefault()`, `PetscFECreateByCell()`, `PetscFECreate()`, `PetscSpaceCreate()`, `PetscDualSpaceCreate()`
1909: @*/
1910: PetscErrorCode PetscFECreateFromSpaces(PetscSpace P, PetscDualSpace Q, PetscQuadrature q, PetscQuadrature fq, PetscFE *fem)
1911: {
1912: PetscInt Nc;
1913: PetscInt p_Ns = -1, p_Nc = -1, q_Ns = -1, q_Nc = -1;
1914: PetscBool p_is_uniform_sum = PETSC_FALSE, p_interleave_basis = PETSC_FALSE, p_interleave_components = PETSC_FALSE;
1915: PetscBool q_is_uniform_sum = PETSC_FALSE, q_interleave_basis = PETSC_FALSE, q_interleave_components = PETSC_FALSE;
1916: const char *prefix;
1918: PetscFunctionBegin;
1919: PetscCall(PetscObjectTypeCompare((PetscObject)P, PETSCSPACESUM, &p_is_uniform_sum));
1920: if (p_is_uniform_sum) {
1921: PetscSpace subsp_0 = NULL;
1922: PetscCall(PetscSpaceSumGetNumSubspaces(P, &p_Ns));
1923: PetscCall(PetscSpaceGetNumComponents(P, &p_Nc));
1924: PetscCall(PetscSpaceSumGetConcatenate(P, &p_is_uniform_sum));
1925: PetscCall(PetscSpaceSumGetInterleave(P, &p_interleave_basis, &p_interleave_components));
1926: for (PetscInt s = 0; s < p_Ns; s++) {
1927: PetscSpace subsp;
1929: PetscCall(PetscSpaceSumGetSubspace(P, s, &subsp));
1930: if (!s) {
1931: subsp_0 = subsp;
1932: } else if (subsp != subsp_0) {
1933: p_is_uniform_sum = PETSC_FALSE;
1934: }
1935: }
1936: }
1937: PetscCall(PetscObjectTypeCompare((PetscObject)Q, PETSCDUALSPACESUM, &q_is_uniform_sum));
1938: if (q_is_uniform_sum) {
1939: PetscDualSpace subsp_0 = NULL;
1940: PetscCall(PetscDualSpaceSumGetNumSubspaces(Q, &q_Ns));
1941: PetscCall(PetscDualSpaceGetNumComponents(Q, &q_Nc));
1942: PetscCall(PetscDualSpaceSumGetConcatenate(Q, &q_is_uniform_sum));
1943: PetscCall(PetscDualSpaceSumGetInterleave(Q, &q_interleave_basis, &q_interleave_components));
1944: for (PetscInt s = 0; s < q_Ns; s++) {
1945: PetscDualSpace subsp;
1947: PetscCall(PetscDualSpaceSumGetSubspace(Q, s, &subsp));
1948: if (!s) {
1949: subsp_0 = subsp;
1950: } else if (subsp != subsp_0) {
1951: q_is_uniform_sum = PETSC_FALSE;
1952: }
1953: }
1954: }
1955: if (p_is_uniform_sum && q_is_uniform_sum && (p_interleave_basis == q_interleave_basis) && (p_interleave_components == q_interleave_components) && (p_Ns == q_Ns) && (p_Nc == q_Nc)) {
1956: PetscSpace scalar_space;
1957: PetscDualSpace scalar_dspace;
1958: PetscFE scalar_fe;
1960: PetscCall(PetscSpaceSumGetSubspace(P, 0, &scalar_space));
1961: PetscCall(PetscDualSpaceSumGetSubspace(Q, 0, &scalar_dspace));
1962: PetscCall(PetscObjectReference((PetscObject)scalar_space));
1963: PetscCall(PetscObjectReference((PetscObject)scalar_dspace));
1964: PetscCall(PetscObjectReference((PetscObject)q));
1965: PetscCall(PetscObjectReference((PetscObject)fq));
1966: PetscCall(PetscFECreateFromSpaces(scalar_space, scalar_dspace, q, fq, &scalar_fe));
1967: PetscCall(PetscFECreateVector(scalar_fe, p_Ns, p_interleave_basis, p_interleave_components, fem));
1968: PetscCall(PetscFEDestroy(&scalar_fe));
1969: } else {
1970: PetscCall(PetscFECreate(PetscObjectComm((PetscObject)P), fem));
1971: PetscCall(PetscFESetType(*fem, PETSCFEBASIC));
1972: }
1973: PetscCall(PetscSpaceGetNumComponents(P, &Nc));
1974: PetscCall(PetscFESetNumComponents(*fem, Nc));
1975: PetscCall(PetscFESetBasisSpace(*fem, P));
1976: PetscCall(PetscFESetDualSpace(*fem, Q));
1977: PetscCall(PetscObjectGetOptionsPrefix((PetscObject)P, &prefix));
1978: PetscCall(PetscObjectSetOptionsPrefix((PetscObject)*fem, prefix));
1979: PetscCall(PetscFESetUp(*fem));
1980: PetscCall(PetscSpaceDestroy(&P));
1981: PetscCall(PetscDualSpaceDestroy(&Q));
1982: PetscCall(PetscFESetQuadrature(*fem, q));
1983: PetscCall(PetscFESetFaceQuadrature(*fem, fq));
1984: PetscCall(PetscQuadratureDestroy(&q));
1985: PetscCall(PetscQuadratureDestroy(&fq));
1986: PetscCall(PetscFESetDefaultName_Private(*fem));
1987: PetscFunctionReturn(PETSC_SUCCESS);
1988: }
1990: static PetscErrorCode PetscFECreate_Internal(MPI_Comm comm, PetscInt dim, PetscInt Nc, DMPolytopeType ct, const char prefix[], PetscInt degree, PetscInt qorder, PetscBool setFromOptions, PetscFE *fem)
1991: {
1992: DM K;
1993: PetscSpace P;
1994: PetscDualSpace Q;
1995: PetscQuadrature q, fq;
1996: PetscBool tensor;
1998: PetscFunctionBegin;
1999: if (prefix) PetscAssertPointer(prefix, 5);
2000: PetscAssertPointer(fem, 9);
2001: switch (ct) {
2002: case DM_POLYTOPE_SEGMENT:
2003: case DM_POLYTOPE_POINT_PRISM_TENSOR:
2004: case DM_POLYTOPE_QUADRILATERAL:
2005: case DM_POLYTOPE_SEG_PRISM_TENSOR:
2006: case DM_POLYTOPE_HEXAHEDRON:
2007: case DM_POLYTOPE_QUAD_PRISM_TENSOR:
2008: tensor = PETSC_TRUE;
2009: break;
2010: default:
2011: tensor = PETSC_FALSE;
2012: }
2013: /* Create space */
2014: PetscCall(PetscSpaceCreate(comm, &P));
2015: PetscCall(PetscSpaceSetType(P, PETSCSPACEPOLYNOMIAL));
2016: PetscCall(PetscObjectSetOptionsPrefix((PetscObject)P, prefix));
2017: PetscCall(PetscSpacePolynomialSetTensor(P, tensor));
2018: PetscCall(PetscSpaceSetNumComponents(P, Nc));
2019: PetscCall(PetscSpaceSetNumVariables(P, dim));
2020: if (degree >= 0) {
2021: PetscCall(PetscSpaceSetDegree(P, degree, PETSC_DETERMINE));
2022: if (ct == DM_POLYTOPE_TRI_PRISM || ct == DM_POLYTOPE_TRI_PRISM_TENSOR) {
2023: PetscSpace Pend, Pside;
2025: PetscCall(PetscSpaceSetNumComponents(P, 1));
2026: PetscCall(PetscSpaceCreate(comm, &Pend));
2027: PetscCall(PetscSpaceSetType(Pend, PETSCSPACEPOLYNOMIAL));
2028: PetscCall(PetscSpacePolynomialSetTensor(Pend, PETSC_FALSE));
2029: PetscCall(PetscSpaceSetNumComponents(Pend, 1));
2030: PetscCall(PetscSpaceSetNumVariables(Pend, dim - 1));
2031: PetscCall(PetscSpaceSetDegree(Pend, degree, PETSC_DETERMINE));
2032: PetscCall(PetscSpaceCreate(comm, &Pside));
2033: PetscCall(PetscSpaceSetType(Pside, PETSCSPACEPOLYNOMIAL));
2034: PetscCall(PetscSpacePolynomialSetTensor(Pside, PETSC_FALSE));
2035: PetscCall(PetscSpaceSetNumComponents(Pside, 1));
2036: PetscCall(PetscSpaceSetNumVariables(Pside, 1));
2037: PetscCall(PetscSpaceSetDegree(Pside, degree, PETSC_DETERMINE));
2038: PetscCall(PetscSpaceSetType(P, PETSCSPACETENSOR));
2039: PetscCall(PetscSpaceTensorSetNumSubspaces(P, 2));
2040: PetscCall(PetscSpaceTensorSetSubspace(P, 0, Pend));
2041: PetscCall(PetscSpaceTensorSetSubspace(P, 1, Pside));
2042: PetscCall(PetscSpaceDestroy(&Pend));
2043: PetscCall(PetscSpaceDestroy(&Pside));
2045: if (Nc > 1) {
2046: PetscSpace scalar_P = P;
2048: PetscCall(PetscSpaceCreate(comm, &P));
2049: PetscCall(PetscSpaceSetNumVariables(P, dim));
2050: PetscCall(PetscSpaceSetNumComponents(P, Nc));
2051: PetscCall(PetscSpaceSetType(P, PETSCSPACESUM));
2052: PetscCall(PetscSpaceSumSetNumSubspaces(P, Nc));
2053: PetscCall(PetscSpaceSumSetConcatenate(P, PETSC_TRUE));
2054: PetscCall(PetscSpaceSumSetInterleave(P, PETSC_TRUE, PETSC_FALSE));
2055: for (PetscInt i = 0; i < Nc; i++) PetscCall(PetscSpaceSumSetSubspace(P, i, scalar_P));
2056: PetscCall(PetscSpaceDestroy(&scalar_P));
2057: }
2058: }
2059: }
2060: if (setFromOptions) PetscCall(PetscSpaceSetFromOptions(P));
2061: PetscCall(PetscSpaceSetUp(P));
2062: PetscCall(PetscSpaceGetDegree(P, °ree, NULL));
2063: PetscCall(PetscSpacePolynomialGetTensor(P, &tensor));
2064: PetscCall(PetscSpaceGetNumComponents(P, &Nc));
2065: /* Create dual space */
2066: PetscCall(PetscDualSpaceCreate(comm, &Q));
2067: PetscCall(PetscDualSpaceSetType(Q, PETSCDUALSPACELAGRANGE));
2068: PetscCall(PetscObjectSetOptionsPrefix((PetscObject)Q, prefix));
2069: PetscCall(DMPlexCreateReferenceCell(PETSC_COMM_SELF, ct, &K));
2070: PetscCall(PetscDualSpaceSetDM(Q, K));
2071: PetscCall(DMDestroy(&K));
2072: PetscCall(PetscDualSpaceSetNumComponents(Q, Nc));
2073: PetscCall(PetscDualSpaceSetOrder(Q, degree));
2074: PetscCall(PetscDualSpaceLagrangeSetTensor(Q, (tensor || (ct == DM_POLYTOPE_TRI_PRISM)) ? PETSC_TRUE : PETSC_FALSE));
2075: if (setFromOptions) PetscCall(PetscDualSpaceSetFromOptions(Q));
2076: PetscCall(PetscDualSpaceSetUp(Q));
2077: /* Create quadrature */
2078: PetscDTSimplexQuadratureType qtype = PETSCDTSIMPLEXQUAD_DEFAULT;
2080: qorder = qorder >= 0 ? qorder : degree;
2081: if (setFromOptions) {
2082: PetscObjectOptionsBegin((PetscObject)P);
2083: PetscCall(PetscOptionsBoundedInt("-petscfe_default_quadrature_order", "Quadrature order is one less than quadrature points per edge", "PetscFECreateDefault", qorder, &qorder, NULL, 0));
2084: PetscCall(PetscOptionsEnum("-petscfe_default_quadrature_type", "Simplex quadrature type", "PetscDTSimplexQuadratureType", PetscDTSimplexQuadratureTypes, (PetscEnum)qtype, (PetscEnum *)&qtype, NULL));
2085: PetscOptionsEnd();
2086: }
2087: PetscCall(PetscDTCreateQuadratureByCell(ct, qorder, qtype, &q, &fq));
2088: /* Create finite element */
2089: PetscCall(PetscFECreateFromSpaces(P, Q, q, fq, fem));
2090: if (setFromOptions) PetscCall(PetscFESetFromOptions(*fem));
2091: PetscFunctionReturn(PETSC_SUCCESS);
2092: }
2094: /*@
2095: PetscFECreateDefault - Create a `PetscFE` for basic FEM computation
2097: Collective
2099: Input Parameters:
2100: + comm - The MPI comm
2101: . dim - The spatial dimension
2102: . Nc - The number of components
2103: . isSimplex - Flag for simplex reference cell, otherwise its a tensor product
2104: . prefix - The options prefix, or `NULL`
2105: - qorder - The quadrature order or `PETSC_DETERMINE` to use `PetscSpace` polynomial degree
2107: Output Parameter:
2108: . fem - The `PetscFE` object
2110: Level: beginner
2112: Note:
2113: Each subobject is SetFromOption() during creation, so that the object may be customized from the command line, using the prefix specified above. See the links below for the particular options available.
2115: .seealso: `PetscFECreateLagrange()`, `PetscFECreateByCell()`, `PetscSpaceSetFromOptions()`, `PetscDualSpaceSetFromOptions()`, `PetscFESetFromOptions()`, `PetscFECreate()`, `PetscSpaceCreate()`, `PetscDualSpaceCreate()`
2116: @*/
2117: PetscErrorCode PetscFECreateDefault(MPI_Comm comm, PetscInt dim, PetscInt Nc, PetscBool isSimplex, const char prefix[], PetscInt qorder, PetscFE *fem)
2118: {
2119: PetscFunctionBegin;
2120: PetscCall(PetscFECreate_Internal(comm, dim, Nc, DMPolytopeTypeSimpleShape(dim, isSimplex), prefix, PETSC_DECIDE, qorder, PETSC_TRUE, fem));
2121: PetscFunctionReturn(PETSC_SUCCESS);
2122: }
2124: /*@
2125: PetscFECreateByCell - Create a `PetscFE` for basic FEM computation
2127: Collective
2129: Input Parameters:
2130: + comm - The MPI comm
2131: . dim - The spatial dimension
2132: . Nc - The number of components
2133: . ct - The celltype of the reference cell
2134: . prefix - The options prefix, or `NULL`
2135: - qorder - The quadrature order or `PETSC_DETERMINE` to use `PetscSpace` polynomial degree
2137: Output Parameter:
2138: . fem - The `PetscFE` object
2140: Level: beginner
2142: Note:
2143: Each subobject is SetFromOption() during creation, so that the object may be customized from the command line, using the prefix specified above. See the links below for the particular options available.
2145: .seealso: `PetscFECreateDefault()`, `PetscFECreateLagrange()`, `PetscSpaceSetFromOptions()`, `PetscDualSpaceSetFromOptions()`, `PetscFESetFromOptions()`, `PetscFECreate()`, `PetscSpaceCreate()`, `PetscDualSpaceCreate()`
2146: @*/
2147: PetscErrorCode PetscFECreateByCell(MPI_Comm comm, PetscInt dim, PetscInt Nc, DMPolytopeType ct, const char prefix[], PetscInt qorder, PetscFE *fem)
2148: {
2149: PetscFunctionBegin;
2150: PetscCall(PetscFECreate_Internal(comm, dim, Nc, ct, prefix, PETSC_DECIDE, qorder, PETSC_TRUE, fem));
2151: PetscFunctionReturn(PETSC_SUCCESS);
2152: }
2154: /*@
2155: PetscFECreateLagrange - Create a `PetscFE` for the basic Lagrange space of degree k
2157: Collective
2159: Input Parameters:
2160: + comm - The MPI comm
2161: . dim - The spatial dimension
2162: . Nc - The number of components
2163: . isSimplex - Flag for simplex reference cell, otherwise its a tensor product
2164: . k - The degree k of the space
2165: - qorder - The quadrature order or `PETSC_DETERMINE` to use `PetscSpace` polynomial degree
2167: Output Parameter:
2168: . fem - The `PetscFE` object
2170: Level: beginner
2172: Note:
2173: For simplices, this element is the space of maximum polynomial degree k, otherwise it is a tensor product of 1D polynomials, each with maximal degree k.
2175: .seealso: `PetscFECreateLagrangeByCell()`, `PetscFECreateDefault()`, `PetscFECreateByCell()`, `PetscFECreate()`, `PetscSpaceCreate()`, `PetscDualSpaceCreate()`
2176: @*/
2177: PetscErrorCode PetscFECreateLagrange(MPI_Comm comm, PetscInt dim, PetscInt Nc, PetscBool isSimplex, PetscInt k, PetscInt qorder, PetscFE *fem)
2178: {
2179: PetscFunctionBegin;
2180: PetscCall(PetscFECreate_Internal(comm, dim, Nc, DMPolytopeTypeSimpleShape(dim, isSimplex), NULL, k, qorder, PETSC_FALSE, fem));
2181: PetscFunctionReturn(PETSC_SUCCESS);
2182: }
2184: /*@
2185: PetscFECreateLagrangeByCell - Create a `PetscFE` for the basic Lagrange space of degree k
2187: Collective
2189: Input Parameters:
2190: + comm - The MPI comm
2191: . dim - The spatial dimension
2192: . Nc - The number of components
2193: . ct - The celltype of the reference cell
2194: . k - The degree k of the space
2195: - qorder - The quadrature order or `PETSC_DETERMINE` to use `PetscSpace` polynomial degree
2197: Output Parameter:
2198: . fem - The `PetscFE` object
2200: Level: beginner
2202: Note:
2203: For simplices, this element is the space of maximum polynomial degree k, otherwise it is a tensor product of 1D polynomials, each with maximal degree k.
2205: .seealso: `PetscFECreateLagrange()`, `PetscFECreateDefault()`, `PetscFECreateByCell()`, `PetscFECreate()`, `PetscSpaceCreate()`, `PetscDualSpaceCreate()`
2206: @*/
2207: PetscErrorCode PetscFECreateLagrangeByCell(MPI_Comm comm, PetscInt dim, PetscInt Nc, DMPolytopeType ct, PetscInt k, PetscInt qorder, PetscFE *fem)
2208: {
2209: PetscFunctionBegin;
2210: PetscCall(PetscFECreate_Internal(comm, dim, Nc, ct, NULL, k, qorder, PETSC_FALSE, fem));
2211: PetscFunctionReturn(PETSC_SUCCESS);
2212: }
2214: /*@
2215: PetscFELimitDegree - Copy a `PetscFE` but limit the degree to be in the given range
2217: Collective
2219: Input Parameters:
2220: + fe - The `PetscFE`
2221: . minDegree - The minimum degree, or `PETSC_DETERMINE` for no limit
2222: - maxDegree - The maximum degree, or `PETSC_DETERMINE` for no limit
2224: Output Parameter:
2225: . newfe - The `PetscFE` object
2227: Level: advanced
2229: Note:
2230: This currently only works for Lagrange elements.
2232: .seealso: `PetscFECreateLagrange()`, `PetscFECreateDefault()`, `PetscFECreateByCell()`, `PetscFECreate()`, `PetscSpaceCreate()`, `PetscDualSpaceCreate()`
2233: @*/
2234: PetscErrorCode PetscFELimitDegree(PetscFE fe, PetscInt minDegree, PetscInt maxDegree, PetscFE *newfe)
2235: {
2236: PetscDualSpace Q;
2237: PetscBool islag, issum;
2238: PetscInt oldk = 0, k;
2240: PetscFunctionBegin;
2241: PetscCall(PetscFEGetDualSpace(fe, &Q));
2242: PetscCall(PetscObjectTypeCompare((PetscObject)Q, PETSCDUALSPACELAGRANGE, &islag));
2243: PetscCall(PetscObjectTypeCompare((PetscObject)Q, PETSCDUALSPACESUM, &issum));
2244: if (islag) {
2245: PetscCall(PetscDualSpaceGetOrder(Q, &oldk));
2246: } else if (issum) {
2247: PetscDualSpace subQ;
2249: PetscCall(PetscDualSpaceSumGetSubspace(Q, 0, &subQ));
2250: PetscCall(PetscDualSpaceGetOrder(subQ, &oldk));
2251: } else {
2252: PetscCall(PetscObjectReference((PetscObject)fe));
2253: *newfe = fe;
2254: PetscFunctionReturn(PETSC_SUCCESS);
2255: }
2256: k = oldk;
2257: if (minDegree >= 0) k = PetscMax(k, minDegree);
2258: if (maxDegree >= 0) k = PetscMin(k, maxDegree);
2259: if (k != oldk) {
2260: DM K;
2261: PetscSpace P;
2262: PetscQuadrature q;
2263: DMPolytopeType ct;
2264: PetscInt dim, Nc;
2266: PetscCall(PetscFEGetBasisSpace(fe, &P));
2267: PetscCall(PetscSpaceGetNumVariables(P, &dim));
2268: PetscCall(PetscSpaceGetNumComponents(P, &Nc));
2269: PetscCall(PetscDualSpaceGetDM(Q, &K));
2270: PetscCall(DMPlexGetCellType(K, 0, &ct));
2271: PetscCall(PetscFECreateLagrangeByCell(PetscObjectComm((PetscObject)fe), dim, Nc, ct, k, PETSC_DETERMINE, newfe));
2272: PetscCall(PetscFEGetQuadrature(fe, &q));
2273: PetscCall(PetscFESetQuadrature(*newfe, q));
2274: } else {
2275: PetscCall(PetscObjectReference((PetscObject)fe));
2276: *newfe = fe;
2277: }
2278: PetscFunctionReturn(PETSC_SUCCESS);
2279: }
2281: /*@
2282: PetscFECreateBrokenElement - Create a discontinuous version of the input `PetscFE`
2284: Collective
2286: Input Parameters:
2287: . cgfe - The continuous `PetscFE` object
2289: Output Parameter:
2290: . dgfe - The discontinuous `PetscFE` object
2292: Level: advanced
2294: Note:
2295: This only works for Lagrange elements.
2297: .seealso: `PetscFECreate()`, `PetscSpaceCreate()`, `PetscDualSpaceCreate()`, `PetscFECreateLagrange()`, `PetscFECreateLagrangeByCell()`, `PetscDualSpaceLagrangeSetContinuity()`
2298: @*/
2299: PetscErrorCode PetscFECreateBrokenElement(PetscFE cgfe, PetscFE *dgfe)
2300: {
2301: PetscSpace P;
2302: PetscDualSpace Q, dgQ;
2303: PetscQuadrature q, fq;
2304: PetscBool is_lagrange, is_sum;
2306: PetscFunctionBegin;
2307: PetscCall(PetscFEGetBasisSpace(cgfe, &P));
2308: PetscCall(PetscObjectReference((PetscObject)P));
2309: PetscCall(PetscFEGetDualSpace(cgfe, &Q));
2310: PetscCall(PetscObjectTypeCompare((PetscObject)Q, PETSCDUALSPACELAGRANGE, &is_lagrange));
2311: PetscCall(PetscObjectTypeCompare((PetscObject)Q, PETSCDUALSPACESUM, &is_sum));
2312: PetscCheck(is_lagrange || is_sum, PETSC_COMM_SELF, PETSC_ERR_SUP, "Can only create broken elements of Lagrange elements");
2313: PetscCall(PetscDualSpaceDuplicate(Q, &dgQ));
2314: PetscCall(PetscDualSpaceLagrangeSetContinuity(dgQ, PETSC_FALSE));
2315: PetscCall(PetscDualSpaceSetUp(dgQ));
2316: PetscCall(PetscFEGetQuadrature(cgfe, &q));
2317: PetscCall(PetscObjectReference((PetscObject)q));
2318: PetscCall(PetscFEGetFaceQuadrature(cgfe, &fq));
2319: PetscCall(PetscObjectReference((PetscObject)fq));
2320: PetscCall(PetscFECreateFromSpaces(P, dgQ, q, fq, dgfe));
2321: PetscFunctionReturn(PETSC_SUCCESS);
2322: }
2324: /*@
2325: PetscFESetName - Names the `PetscFE` and its subobjects
2327: Not Collective
2329: Input Parameters:
2330: + fe - The `PetscFE`
2331: - name - The name
2333: Level: intermediate
2335: .seealso: `PetscFECreate()`, `PetscSpaceCreate()`, `PetscDualSpaceCreate()`
2336: @*/
2337: PetscErrorCode PetscFESetName(PetscFE fe, const char name[])
2338: {
2339: PetscSpace P;
2340: PetscDualSpace Q;
2342: PetscFunctionBegin;
2343: PetscCall(PetscFEGetBasisSpace(fe, &P));
2344: PetscCall(PetscFEGetDualSpace(fe, &Q));
2345: PetscCall(PetscObjectSetName((PetscObject)fe, name));
2346: PetscCall(PetscObjectSetName((PetscObject)P, name));
2347: PetscCall(PetscObjectSetName((PetscObject)Q, name));
2348: PetscFunctionReturn(PETSC_SUCCESS);
2349: }
2351: PetscErrorCode PetscFEEvaluateFieldJets_Internal(PetscDS ds, PetscInt Nf, PetscInt r, PetscInt q, PetscTabulation T[], PetscFEGeom *fegeom, const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscScalar u[], PetscScalar u_x[], PetscScalar u_t[])
2352: {
2353: PetscInt dOffset = 0, fOffset = 0, f, g;
2355: for (f = 0; f < Nf; ++f) {
2356: PetscCheck(r < T[f]->Nr, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Replica number %" PetscInt_FMT " should be in [0, %" PetscInt_FMT ")", r, T[f]->Nr);
2357: PetscCheck(q < T[f]->Np, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Point number %" PetscInt_FMT " should be in [0, %" PetscInt_FMT ")", q, T[f]->Np);
2358: PetscFE fe;
2359: const PetscInt k = ds->jetDegree[f];
2360: const PetscInt cdim = T[f]->cdim;
2361: const PetscInt dE = fegeom->dimEmbed;
2362: const PetscInt Nq = T[f]->Np;
2363: const PetscInt Nbf = T[f]->Nb;
2364: const PetscInt Ncf = T[f]->Nc;
2365: const PetscReal *Bq = &T[f]->T[0][(r * Nq + q) * Nbf * Ncf];
2366: const PetscReal *Dq = &T[f]->T[1][(r * Nq + q) * Nbf * Ncf * cdim];
2367: const PetscReal *Hq = k > 1 ? &T[f]->T[2][(r * Nq + q) * Nbf * Ncf * cdim * cdim] : NULL;
2368: PetscInt hOffset = 0, b, c, d;
2370: PetscCall(PetscDSGetDiscretization(ds, f, (PetscObject *)&fe));
2371: for (c = 0; c < Ncf; ++c) u[fOffset + c] = 0.0;
2372: for (d = 0; d < dE * Ncf; ++d) u_x[fOffset * dE + d] = 0.0;
2373: for (b = 0; b < Nbf; ++b) {
2374: for (c = 0; c < Ncf; ++c) {
2375: const PetscInt cidx = b * Ncf + c;
2377: u[fOffset + c] += Bq[cidx] * coefficients[dOffset + b];
2378: for (d = 0; d < cdim; ++d) u_x[(fOffset + c) * dE + d] += Dq[cidx * cdim + d] * coefficients[dOffset + b];
2379: }
2380: }
2381: if (k > 1) {
2382: for (g = 0; g < Nf; ++g) hOffset += T[g]->Nc * dE;
2383: for (d = 0; d < dE * dE * Ncf; ++d) u_x[hOffset + fOffset * dE * dE + d] = 0.0;
2384: for (b = 0; b < Nbf; ++b) {
2385: for (c = 0; c < Ncf; ++c) {
2386: const PetscInt cidx = b * Ncf + c;
2388: for (d = 0; d < cdim * cdim; ++d) u_x[hOffset + (fOffset + c) * dE * dE + d] += Hq[cidx * cdim * cdim + d] * coefficients[dOffset + b];
2389: }
2390: }
2391: PetscCall(PetscFEPushforwardHessian(fe, fegeom, 1, &u_x[hOffset + fOffset * dE * dE]));
2392: }
2393: PetscCall(PetscFEPushforward(fe, fegeom, 1, &u[fOffset]));
2394: PetscCall(PetscFEPushforwardGradient(fe, fegeom, 1, &u_x[fOffset * dE]));
2395: if (u_t) {
2396: for (c = 0; c < Ncf; ++c) u_t[fOffset + c] = 0.0;
2397: for (b = 0; b < Nbf; ++b) {
2398: for (c = 0; c < Ncf; ++c) {
2399: const PetscInt cidx = b * Ncf + c;
2401: u_t[fOffset + c] += Bq[cidx] * coefficients_t[dOffset + b];
2402: }
2403: }
2404: PetscCall(PetscFEPushforward(fe, fegeom, 1, &u_t[fOffset]));
2405: }
2406: fOffset += Ncf;
2407: dOffset += Nbf;
2408: }
2409: return PETSC_SUCCESS;
2410: }
2412: PetscErrorCode PetscFEEvaluateFieldJets_Hybrid_Internal(PetscDS ds, PetscInt Nf, PetscInt rc, PetscInt qc, PetscTabulation Tab[], const PetscInt rf[], const PetscInt qf[], PetscTabulation Tabf[], PetscFEGeom *fegeom, PetscFEGeom *fegeomNbr, const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscScalar u[], PetscScalar u_x[], PetscScalar u_t[])
2413: {
2414: PetscInt dOffset = 0, fOffset = 0, f;
2416: /* f is the field number in the DS */
2417: for (f = 0; f < Nf; ++f) {
2418: PetscBool isCohesive;
2419: PetscInt Ns, s;
2421: if (!Tab[f]) continue;
2422: PetscCall(PetscDSGetCohesive(ds, f, &isCohesive));
2423: Ns = isCohesive ? 1 : 2;
2424: {
2425: PetscTabulation T = isCohesive ? Tab[f] : Tabf[f];
2426: PetscFE fe = (PetscFE)ds->disc[f];
2427: const PetscInt dEt = T->cdim;
2428: const PetscInt dE = fegeom->dimEmbed;
2429: const PetscInt Nq = T->Np;
2430: const PetscInt Nbf = T->Nb;
2431: const PetscInt Ncf = T->Nc;
2433: for (s = 0; s < Ns; ++s) {
2434: const PetscInt r = isCohesive ? rc : rf[s];
2435: const PetscInt q = isCohesive ? qc : qf[s];
2436: const PetscReal *Bq = &T->T[0][(r * Nq + q) * Nbf * Ncf];
2437: const PetscReal *Dq = &T->T[1][(r * Nq + q) * Nbf * Ncf * dEt];
2438: PetscInt b, c, d;
2440: PetscCheck(r < T->Nr, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field %" PetscInt_FMT " Side %" PetscInt_FMT " Replica number %" PetscInt_FMT " should be in [0, %" PetscInt_FMT ")", f, s, r, T->Nr);
2441: PetscCheck(q < T->Np, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field %" PetscInt_FMT " Side %" PetscInt_FMT " Point number %" PetscInt_FMT " should be in [0, %" PetscInt_FMT ")", f, s, q, T->Np);
2442: for (c = 0; c < Ncf; ++c) u[fOffset + c] = 0.0;
2443: for (d = 0; d < dE * Ncf; ++d) u_x[fOffset * dE + d] = 0.0;
2444: for (b = 0; b < Nbf; ++b) {
2445: for (c = 0; c < Ncf; ++c) {
2446: const PetscInt cidx = b * Ncf + c;
2448: u[fOffset + c] += Bq[cidx] * coefficients[dOffset + b];
2449: for (d = 0; d < dEt; ++d) u_x[(fOffset + c) * dE + d] += Dq[cidx * dEt + d] * coefficients[dOffset + b];
2450: }
2451: }
2452: PetscCall(PetscFEPushforward(fe, isCohesive ? fegeom : &fegeomNbr[s], 1, &u[fOffset]));
2453: PetscCall(PetscFEPushforwardGradient(fe, isCohesive ? fegeom : &fegeomNbr[s], 1, &u_x[fOffset * dE]));
2454: if (u_t) {
2455: for (c = 0; c < Ncf; ++c) u_t[fOffset + c] = 0.0;
2456: for (b = 0; b < Nbf; ++b) {
2457: for (c = 0; c < Ncf; ++c) {
2458: const PetscInt cidx = b * Ncf + c;
2460: u_t[fOffset + c] += Bq[cidx] * coefficients_t[dOffset + b];
2461: }
2462: }
2463: PetscCall(PetscFEPushforward(fe, fegeom, 1, &u_t[fOffset]));
2464: }
2465: fOffset += Ncf;
2466: dOffset += Nbf;
2467: }
2468: }
2469: }
2470: return PETSC_SUCCESS;
2471: }
2473: PetscErrorCode PetscFEEvaluateFaceFields_Internal(PetscDS prob, PetscInt field, PetscInt faceLoc, const PetscScalar coefficients[], PetscScalar u[])
2474: {
2475: PetscFE fe;
2476: PetscTabulation Tc;
2477: PetscInt b, c;
2479: if (!prob) return PETSC_SUCCESS;
2480: PetscCall(PetscDSGetDiscretization(prob, field, (PetscObject *)&fe));
2481: PetscCall(PetscFEGetFaceCentroidTabulation(fe, &Tc));
2482: {
2483: const PetscReal *faceBasis = Tc->T[0];
2484: const PetscInt Nb = Tc->Nb;
2485: const PetscInt Nc = Tc->Nc;
2487: for (c = 0; c < Nc; ++c) u[c] = 0.0;
2488: for (b = 0; b < Nb; ++b) {
2489: for (c = 0; c < Nc; ++c) u[c] += coefficients[b] * faceBasis[(faceLoc * Nb + b) * Nc + c];
2490: }
2491: }
2492: return PETSC_SUCCESS;
2493: }
2495: PetscErrorCode PetscFEUpdateElementVec_Internal(PetscFE fe, PetscTabulation T, PetscInt r, PetscScalar tmpBasis[], PetscScalar tmpBasisDer[], PetscInt e, PetscFEGeom *fegeom, PetscScalar f0[], PetscScalar f1[], PetscScalar elemVec[])
2496: {
2497: PetscFEGeom pgeom;
2498: const PetscInt dEt = T->cdim;
2499: const PetscInt dE = fegeom->dimEmbed;
2500: const PetscInt Nq = T->Np;
2501: const PetscInt Nb = T->Nb;
2502: const PetscInt Nc = T->Nc;
2503: const PetscReal *basis = &T->T[0][r * Nq * Nb * Nc];
2504: const PetscReal *basisDer = &T->T[1][r * Nq * Nb * Nc * dEt];
2505: PetscInt q, b, c, d;
2507: for (q = 0; q < Nq; ++q) {
2508: for (b = 0; b < Nb; ++b) {
2509: for (c = 0; c < Nc; ++c) {
2510: const PetscInt bcidx = b * Nc + c;
2512: tmpBasis[bcidx] = basis[q * Nb * Nc + bcidx];
2513: for (d = 0; d < dEt; ++d) tmpBasisDer[bcidx * dE + d] = basisDer[q * Nb * Nc * dEt + bcidx * dEt + d];
2514: for (d = dEt; d < dE; ++d) tmpBasisDer[bcidx * dE + d] = 0.0;
2515: }
2516: }
2517: PetscCall(PetscFEGeomGetCellPoint(fegeom, e, q, &pgeom));
2518: PetscCall(PetscFEPushforward(fe, &pgeom, Nb, tmpBasis));
2519: PetscCall(PetscFEPushforwardGradient(fe, &pgeom, Nb, tmpBasisDer));
2520: for (b = 0; b < Nb; ++b) {
2521: for (c = 0; c < Nc; ++c) {
2522: const PetscInt bcidx = b * Nc + c;
2523: const PetscInt qcidx = q * Nc + c;
2525: elemVec[b] += tmpBasis[bcidx] * f0[qcidx];
2526: for (d = 0; d < dE; ++d) elemVec[b] += tmpBasisDer[bcidx * dE + d] * f1[qcidx * dE + d];
2527: }
2528: }
2529: }
2530: return PETSC_SUCCESS;
2531: }
2533: PetscErrorCode PetscFEUpdateElementVec_Hybrid_Internal(PetscFE fe, PetscTabulation T, PetscInt r, PetscInt side, PetscScalar tmpBasis[], PetscScalar tmpBasisDer[], PetscFEGeom *fegeom, PetscScalar f0[], PetscScalar f1[], PetscScalar elemVec[])
2534: {
2535: const PetscInt dE = T->cdim;
2536: const PetscInt Nq = T->Np;
2537: const PetscInt Nb = T->Nb;
2538: const PetscInt Nc = T->Nc;
2539: const PetscReal *basis = &T->T[0][r * Nq * Nb * Nc];
2540: const PetscReal *basisDer = &T->T[1][r * Nq * Nb * Nc * dE];
2542: for (PetscInt q = 0; q < Nq; ++q) {
2543: for (PetscInt b = 0; b < Nb; ++b) {
2544: for (PetscInt c = 0; c < Nc; ++c) {
2545: const PetscInt bcidx = b * Nc + c;
2547: tmpBasis[bcidx] = basis[q * Nb * Nc + bcidx];
2548: for (PetscInt d = 0; d < dE; ++d) tmpBasisDer[bcidx * dE + d] = basisDer[q * Nb * Nc * dE + bcidx * dE + d];
2549: }
2550: }
2551: PetscCall(PetscFEPushforward(fe, fegeom, Nb, tmpBasis));
2552: // TODO This is currently broken since we do not pull the geometry down to the lower dimension
2553: // PetscCall(PetscFEPushforwardGradient(fe, fegeom, Nb, tmpBasisDer));
2554: if (side == 2) {
2555: // Integrating over whole cohesive cell, so insert for both sides
2556: for (PetscInt s = 0; s < 2; ++s) {
2557: for (PetscInt b = 0; b < Nb; ++b) {
2558: for (PetscInt c = 0; c < Nc; ++c) {
2559: const PetscInt bcidx = b * Nc + c;
2560: const PetscInt qcidx = (q * 2 + s) * Nc + c;
2562: elemVec[Nb * s + b] += tmpBasis[bcidx] * f0[qcidx];
2563: for (PetscInt d = 0; d < dE; ++d) elemVec[Nb * s + b] += tmpBasisDer[bcidx * dE + d] * f1[qcidx * dE + d];
2564: }
2565: }
2566: }
2567: } else {
2568: // Integrating over endcaps of cohesive cell, so insert for correct side
2569: for (PetscInt b = 0; b < Nb; ++b) {
2570: for (PetscInt c = 0; c < Nc; ++c) {
2571: const PetscInt bcidx = b * Nc + c;
2572: const PetscInt qcidx = q * Nc + c;
2574: elemVec[Nb * side + b] += tmpBasis[bcidx] * f0[qcidx];
2575: for (PetscInt d = 0; d < dE; ++d) elemVec[Nb * side + b] += tmpBasisDer[bcidx * dE + d] * f1[qcidx * dE + d];
2576: }
2577: }
2578: }
2579: }
2580: return PETSC_SUCCESS;
2581: }
2583: #define petsc_elemmat_kernel_g1(_NbI, _NcI, _NbJ, _NcJ, _dE) \
2584: do { \
2585: for (PetscInt fc = 0; fc < (_NcI); ++fc) { \
2586: for (PetscInt gc = 0; gc < (_NcJ); ++gc) { \
2587: const PetscScalar *G = g1 + (fc * (_NcJ) + gc) * _dE; \
2588: for (PetscInt f = 0; f < (_NbI); ++f) { \
2589: const PetscScalar tBIv = tmpBasisI[f * (_NcI) + fc]; \
2590: for (PetscInt g = 0; g < (_NbJ); ++g) { \
2591: const PetscScalar *tBDJ = tmpBasisDerJ + (g * (_NcJ) + gc) * (_dE); \
2592: PetscScalar s = 0.0; \
2593: for (PetscInt df = 0; df < _dE; ++df) s += G[df] * tBDJ[df]; \
2594: elemMat[(offsetI + f) * totDim + (offsetJ + g)] += s * tBIv; \
2595: } \
2596: } \
2597: } \
2598: } \
2599: } while (0)
2601: #define petsc_elemmat_kernel_g2(_NbI, _NcI, _NbJ, _NcJ, _dE) \
2602: do { \
2603: for (PetscInt gc = 0; gc < (_NcJ); ++gc) { \
2604: for (PetscInt fc = 0; fc < (_NcI); ++fc) { \
2605: const PetscScalar *G = g2 + (fc * (_NcJ) + gc) * _dE; \
2606: for (PetscInt g = 0; g < (_NbJ); ++g) { \
2607: const PetscScalar tBJv = tmpBasisJ[g * (_NcJ) + gc]; \
2608: for (PetscInt f = 0; f < (_NbI); ++f) { \
2609: const PetscScalar *tBDI = tmpBasisDerI + (f * (_NcI) + fc) * (_dE); \
2610: PetscScalar s = 0.0; \
2611: for (PetscInt df = 0; df < _dE; ++df) s += tBDI[df] * G[df]; \
2612: elemMat[(offsetI + f) * totDim + (offsetJ + g)] += s * tBJv; \
2613: } \
2614: } \
2615: } \
2616: } \
2617: } while (0)
2619: #define petsc_elemmat_kernel_g3(_NbI, _NcI, _NbJ, _NcJ, _dE) \
2620: do { \
2621: for (PetscInt fc = 0; fc < (_NcI); ++fc) { \
2622: for (PetscInt gc = 0; gc < (_NcJ); ++gc) { \
2623: const PetscScalar *G = g3 + (fc * (_NcJ) + gc) * (_dE) * (_dE); \
2624: for (PetscInt f = 0; f < (_NbI); ++f) { \
2625: const PetscScalar *tBDI = tmpBasisDerI + (f * (_NcI) + fc) * (_dE); \
2626: for (PetscInt g = 0; g < (_NbJ); ++g) { \
2627: PetscScalar s = 0.0; \
2628: const PetscScalar *tBDJ = tmpBasisDerJ + (g * (_NcJ) + gc) * (_dE); \
2629: for (PetscInt df = 0; df < (_dE); ++df) { \
2630: for (PetscInt dg = 0; dg < (_dE); ++dg) s += tBDI[df] * G[df * (_dE) + dg] * tBDJ[dg]; \
2631: } \
2632: elemMat[(offsetI + f) * totDim + (offsetJ + g)] += s; \
2633: } \
2634: } \
2635: } \
2636: } \
2637: } while (0)
2639: PetscErrorCode PetscFEUpdateElementMat_Internal(PetscFE feI, PetscFE feJ, PetscInt r, PetscInt q, PetscTabulation TI, PetscScalar tmpBasisI[], PetscScalar tmpBasisDerI[], PetscTabulation TJ, PetscScalar tmpBasisJ[], PetscScalar tmpBasisDerJ[], PetscFEGeom *fegeom, const PetscScalar g0[], const PetscScalar g1[], const PetscScalar g2[], const PetscScalar g3[], PetscInt totDim, PetscInt offsetI, PetscInt offsetJ, PetscScalar elemMat[])
2640: {
2641: const PetscInt cdim = TI->cdim;
2642: const PetscInt dE = fegeom->dimEmbed;
2643: const PetscInt NqI = TI->Np;
2644: const PetscInt NbI = TI->Nb;
2645: const PetscInt NcI = TI->Nc;
2646: const PetscReal *basisI = &TI->T[0][(r * NqI + q) * NbI * NcI];
2647: const PetscReal *basisDerI = &TI->T[1][(r * NqI + q) * NbI * NcI * cdim];
2648: const PetscInt NqJ = TJ->Np;
2649: const PetscInt NbJ = TJ->Nb;
2650: const PetscInt NcJ = TJ->Nc;
2651: const PetscReal *basisJ = &TJ->T[0][(r * NqJ + q) * NbJ * NcJ];
2652: const PetscReal *basisDerJ = &TJ->T[1][(r * NqJ + q) * NbJ * NcJ * cdim];
2654: for (PetscInt f = 0; f < NbI; ++f) {
2655: for (PetscInt fc = 0; fc < NcI; ++fc) {
2656: const PetscInt fidx = f * NcI + fc; /* Test function basis index */
2658: tmpBasisI[fidx] = basisI[fidx];
2659: for (PetscInt df = 0; df < cdim; ++df) tmpBasisDerI[fidx * dE + df] = basisDerI[fidx * cdim + df];
2660: }
2661: }
2662: PetscCall(PetscFEPushforward(feI, fegeom, NbI, tmpBasisI));
2663: PetscCall(PetscFEPushforwardGradient(feI, fegeom, NbI, tmpBasisDerI));
2664: if (feI != feJ) {
2665: for (PetscInt g = 0; g < NbJ; ++g) {
2666: for (PetscInt gc = 0; gc < NcJ; ++gc) {
2667: const PetscInt gidx = g * NcJ + gc; /* Trial function basis index */
2669: tmpBasisJ[gidx] = basisJ[gidx];
2670: for (PetscInt dg = 0; dg < cdim; ++dg) tmpBasisDerJ[gidx * dE + dg] = basisDerJ[gidx * cdim + dg];
2671: }
2672: }
2673: PetscCall(PetscFEPushforward(feJ, fegeom, NbJ, tmpBasisJ));
2674: PetscCall(PetscFEPushforwardGradient(feJ, fegeom, NbJ, tmpBasisDerJ));
2675: } else {
2676: tmpBasisJ = tmpBasisI;
2677: tmpBasisDerJ = tmpBasisDerI;
2678: }
2679: if (PetscUnlikely(g0)) {
2680: for (PetscInt f = 0; f < NbI; ++f) {
2681: const PetscInt i = offsetI + f; /* Element matrix row */
2683: for (PetscInt fc = 0; fc < NcI; ++fc) {
2684: const PetscScalar bI = tmpBasisI[f * NcI + fc]; /* Test function basis value */
2686: for (PetscInt g = 0; g < NbJ; ++g) {
2687: const PetscInt j = offsetJ + g; /* Element matrix column */
2688: const PetscInt fOff = i * totDim + j;
2690: for (PetscInt gc = 0; gc < NcJ; ++gc) elemMat[fOff] += bI * g0[fc * NcJ + gc] * tmpBasisJ[g * NcJ + gc];
2691: }
2692: }
2693: }
2694: }
2695: if (PetscUnlikely(g1)) {
2696: #if 1
2697: if (dE == 2) {
2698: petsc_elemmat_kernel_g1(NbI, NcI, NbJ, NcJ, 2);
2699: } else if (dE == 3) {
2700: petsc_elemmat_kernel_g1(NbI, NcI, NbJ, NcJ, 3);
2701: } else {
2702: petsc_elemmat_kernel_g1(NbI, NcI, NbJ, NcJ, dE);
2703: }
2704: #else
2705: for (PetscInt f = 0; f < NbI; ++f) {
2706: const PetscInt i = offsetI + f; /* Element matrix row */
2708: for (PetscInt fc = 0; fc < NcI; ++fc) {
2709: const PetscScalar bI = tmpBasisI[f * NcI + fc]; /* Test function basis value */
2711: for (PetscInt g = 0; g < NbJ; ++g) {
2712: const PetscInt j = offsetJ + g; /* Element matrix column */
2713: const PetscInt fOff = i * totDim + j;
2715: for (PetscInt gc = 0; gc < NcJ; ++gc) {
2716: const PetscInt gidx = g * NcJ + gc; /* Trial function basis index */
2718: for (PetscInt df = 0; df < dE; ++df) elemMat[fOff] += bI * g1[(fc * NcJ + gc) * dE + df] * tmpBasisDerJ[gidx * dE + df];
2719: }
2720: }
2721: }
2722: }
2723: #endif
2724: }
2725: if (PetscUnlikely(g2)) {
2726: #if 1
2727: if (dE == 2) {
2728: petsc_elemmat_kernel_g2(NbI, NcI, NbJ, NcJ, 2);
2729: } else if (dE == 3) {
2730: petsc_elemmat_kernel_g2(NbI, NcI, NbJ, NcJ, 3);
2731: } else {
2732: petsc_elemmat_kernel_g2(NbI, NcI, NbJ, NcJ, dE);
2733: }
2734: #else
2735: for (PetscInt g = 0; g < NbJ; ++g) {
2736: const PetscInt j = offsetJ + g; /* Element matrix column */
2738: for (PetscInt gc = 0; gc < NcJ; ++gc) {
2739: const PetscScalar bJ = tmpBasisJ[g * NcJ + gc]; /* Trial function basis value */
2741: for (PetscInt f = 0; f < NbI; ++f) {
2742: const PetscInt i = offsetI + f; /* Element matrix row */
2743: const PetscInt fOff = i * totDim + j;
2745: for (PetscInt fc = 0; fc < NcI; ++fc) {
2746: const PetscInt fidx = f * NcI + fc; /* Test function basis index */
2748: for (PetscInt df = 0; df < dE; ++df) elemMat[fOff] += tmpBasisDerI[fidx * dE + df] * g2[(fc * NcJ + gc) * dE + df] * bJ;
2749: }
2750: }
2751: }
2752: }
2753: #endif
2754: }
2755: if (PetscUnlikely(g3)) {
2756: #if 1
2757: if (dE == 2) {
2758: petsc_elemmat_kernel_g3(NbI, NcI, NbJ, NcJ, 2);
2759: } else if (dE == 3) {
2760: petsc_elemmat_kernel_g3(NbI, NcI, NbJ, NcJ, 3);
2761: } else {
2762: petsc_elemmat_kernel_g3(NbI, NcI, NbJ, NcJ, dE);
2763: }
2764: #else
2765: for (PetscInt f = 0; f < NbI; ++f) {
2766: const PetscInt i = offsetI + f; /* Element matrix row */
2768: for (PetscInt fc = 0; fc < NcI; ++fc) {
2769: const PetscInt fidx = f * NcI + fc; /* Test function basis index */
2771: for (PetscInt g = 0; g < NbJ; ++g) {
2772: const PetscInt j = offsetJ + g; /* Element matrix column */
2773: const PetscInt fOff = i * totDim + j;
2775: for (PetscInt gc = 0; gc < NcJ; ++gc) {
2776: const PetscInt gidx = g * NcJ + gc; /* Trial function basis index */
2778: for (PetscInt df = 0; df < dE; ++df) {
2779: for (PetscInt dg = 0; dg < dE; ++dg) elemMat[fOff] += tmpBasisDerI[fidx * dE + df] * g3[((fc * NcJ + gc) * dE + df) * dE + dg] * tmpBasisDerJ[gidx * dE + dg];
2780: }
2781: }
2782: }
2783: }
2784: }
2785: #endif
2786: }
2787: return PETSC_SUCCESS;
2788: }
2790: #undef petsc_elemmat_kernel_g1
2791: #undef petsc_elemmat_kernel_g2
2792: #undef petsc_elemmat_kernel_g3
2794: PetscErrorCode PetscFEUpdateElementMat_Hybrid_Internal(PetscFE feI, PetscBool isHybridI, PetscFE feJ, PetscBool isHybridJ, PetscInt r, PetscInt s, PetscInt t, PetscInt q, PetscTabulation TI, PetscScalar tmpBasisI[], PetscScalar tmpBasisDerI[], PetscTabulation TJ, PetscScalar tmpBasisJ[], PetscScalar tmpBasisDerJ[], PetscFEGeom *fegeom, const PetscScalar g0[], const PetscScalar g1[], const PetscScalar g2[], const PetscScalar g3[], PetscInt eOffset, PetscInt totDim, PetscInt offsetI, PetscInt offsetJ, PetscScalar elemMat[])
2795: {
2796: const PetscInt dE = TI->cdim;
2797: const PetscInt NqI = TI->Np;
2798: const PetscInt NbI = TI->Nb;
2799: const PetscInt NcI = TI->Nc;
2800: const PetscReal *basisI = &TI->T[0][(r * NqI + q) * NbI * NcI];
2801: const PetscReal *basisDerI = &TI->T[1][(r * NqI + q) * NbI * NcI * dE];
2802: const PetscInt NqJ = TJ->Np;
2803: const PetscInt NbJ = TJ->Nb;
2804: const PetscInt NcJ = TJ->Nc;
2805: const PetscReal *basisJ = &TJ->T[0][(r * NqJ + q) * NbJ * NcJ];
2806: const PetscReal *basisDerJ = &TJ->T[1][(r * NqJ + q) * NbJ * NcJ * dE];
2807: const PetscInt so = isHybridI ? 0 : s;
2808: const PetscInt to = isHybridJ ? 0 : t;
2809: PetscInt f, fc, g, gc, df, dg;
2811: for (f = 0; f < NbI; ++f) {
2812: for (fc = 0; fc < NcI; ++fc) {
2813: const PetscInt fidx = f * NcI + fc; /* Test function basis index */
2815: tmpBasisI[fidx] = basisI[fidx];
2816: for (df = 0; df < dE; ++df) tmpBasisDerI[fidx * dE + df] = basisDerI[fidx * dE + df];
2817: }
2818: }
2819: PetscCall(PetscFEPushforward(feI, fegeom, NbI, tmpBasisI));
2820: PetscCall(PetscFEPushforwardGradient(feI, fegeom, NbI, tmpBasisDerI));
2821: for (g = 0; g < NbJ; ++g) {
2822: for (gc = 0; gc < NcJ; ++gc) {
2823: const PetscInt gidx = g * NcJ + gc; /* Trial function basis index */
2825: tmpBasisJ[gidx] = basisJ[gidx];
2826: for (dg = 0; dg < dE; ++dg) tmpBasisDerJ[gidx * dE + dg] = basisDerJ[gidx * dE + dg];
2827: }
2828: }
2829: PetscCall(PetscFEPushforward(feJ, fegeom, NbJ, tmpBasisJ));
2830: // TODO This is currently broken since we do not pull the geometry down to the lower dimension
2831: // PetscCall(PetscFEPushforwardGradient(feJ, fegeom, NbJ, tmpBasisDerJ));
2832: for (f = 0; f < NbI; ++f) {
2833: for (fc = 0; fc < NcI; ++fc) {
2834: const PetscInt fidx = f * NcI + fc; /* Test function basis index */
2835: const PetscInt i = offsetI + NbI * so + f; /* Element matrix row */
2836: for (g = 0; g < NbJ; ++g) {
2837: for (gc = 0; gc < NcJ; ++gc) {
2838: const PetscInt gidx = g * NcJ + gc; /* Trial function basis index */
2839: const PetscInt j = offsetJ + NbJ * to + g; /* Element matrix column */
2840: const PetscInt fOff = eOffset + i * totDim + j;
2842: elemMat[fOff] += tmpBasisI[fidx] * g0[fc * NcJ + gc] * tmpBasisJ[gidx];
2843: for (df = 0; df < dE; ++df) {
2844: elemMat[fOff] += tmpBasisI[fidx] * g1[(fc * NcJ + gc) * dE + df] * tmpBasisDerJ[gidx * dE + df];
2845: elemMat[fOff] += tmpBasisDerI[fidx * dE + df] * g2[(fc * NcJ + gc) * dE + df] * tmpBasisJ[gidx];
2846: for (dg = 0; dg < dE; ++dg) elemMat[fOff] += tmpBasisDerI[fidx * dE + df] * g3[((fc * NcJ + gc) * dE + df) * dE + dg] * tmpBasisDerJ[gidx * dE + dg];
2847: }
2848: }
2849: }
2850: }
2851: }
2852: return PETSC_SUCCESS;
2853: }
2855: PetscErrorCode PetscFECreateCellGeometry(PetscFE fe, PetscQuadrature quad, PetscFEGeom *cgeom)
2856: {
2857: PetscDualSpace dsp;
2858: DM dm;
2859: PetscQuadrature quadDef;
2860: PetscInt dim, cdim, Nq;
2862: PetscFunctionBegin;
2863: PetscCall(PetscFEGetDualSpace(fe, &dsp));
2864: PetscCall(PetscDualSpaceGetDM(dsp, &dm));
2865: PetscCall(DMGetDimension(dm, &dim));
2866: PetscCall(DMGetCoordinateDim(dm, &cdim));
2867: PetscCall(PetscFEGetQuadrature(fe, &quadDef));
2868: quad = quad ? quad : quadDef;
2869: PetscCall(PetscQuadratureGetData(quad, NULL, NULL, &Nq, NULL, NULL));
2870: PetscCall(PetscMalloc1(Nq * cdim, &cgeom->v));
2871: PetscCall(PetscMalloc1(Nq * cdim * cdim, &cgeom->J));
2872: PetscCall(PetscMalloc1(Nq * cdim * cdim, &cgeom->invJ));
2873: PetscCall(PetscMalloc1(Nq, &cgeom->detJ));
2874: cgeom->dim = dim;
2875: cgeom->dimEmbed = cdim;
2876: cgeom->numCells = 1;
2877: cgeom->numPoints = Nq;
2878: PetscCall(DMPlexComputeCellGeometryFEM(dm, 0, quad, cgeom->v, cgeom->J, cgeom->invJ, cgeom->detJ));
2879: PetscFunctionReturn(PETSC_SUCCESS);
2880: }
2882: PetscErrorCode PetscFEDestroyCellGeometry(PetscFE fe, PetscFEGeom *cgeom)
2883: {
2884: PetscFunctionBegin;
2885: PetscCall(PetscFree(cgeom->v));
2886: PetscCall(PetscFree(cgeom->J));
2887: PetscCall(PetscFree(cgeom->invJ));
2888: PetscCall(PetscFree(cgeom->detJ));
2889: PetscFunctionReturn(PETSC_SUCCESS);
2890: }
2892: #if 0
2893: PetscErrorCode PetscFEUpdateElementMat_Internal_SparseIndices(PetscTabulation TI, PetscTabulation TJ, PetscInt dimEmbed, const PetscInt g0[], const PetscInt g1[], const PetscInt g2[], const PetscInt g3[], PetscInt totDim, PetscInt offsetI, PetscInt offsetJ, PetscInt *n_g0, PetscInt **g0_idxs_out, PetscInt *n_g1, PetscInt **g1_idxs_out, PetscInt *n_g2, PetscInt **g2_idxs_out, PetscInt *n_g3, PetscInt **g3_idxs_out)
2894: {
2895: const PetscInt dE = dimEmbed;
2896: const PetscInt NbI = TI->Nb;
2897: const PetscInt NcI = TI->Nc;
2898: const PetscInt NbJ = TJ->Nb;
2899: const PetscInt NcJ = TJ->Nc;
2900: PetscBool has_g0 = g0 ? PETSC_TRUE : PETSC_FALSE;
2901: PetscBool has_g1 = g1 ? PETSC_TRUE : PETSC_FALSE;
2902: PetscBool has_g2 = g2 ? PETSC_TRUE : PETSC_FALSE;
2903: PetscBool has_g3 = g3 ? PETSC_TRUE : PETSC_FALSE;
2904: PetscInt *g0_idxs = NULL, *g1_idxs = NULL, *g2_idxs = NULL, *g3_idxs = NULL;
2905: PetscInt g0_i, g1_i, g2_i, g3_i;
2907: PetscFunctionBegin;
2908: g0_i = g1_i = g2_i = g3_i = 0;
2909: if (has_g0)
2910: for (PetscInt i = 0; i < NcI * NcJ; i++)
2911: if (g0[i]) g0_i += NbI * NbJ;
2912: if (has_g1)
2913: for (PetscInt i = 0; i < NcI * NcJ * dE; i++)
2914: if (g1[i]) g1_i += NbI * NbJ;
2915: if (has_g2)
2916: for (PetscInt i = 0; i < NcI * NcJ * dE; i++)
2917: if (g2[i]) g2_i += NbI * NbJ;
2918: if (has_g3)
2919: for (PetscInt i = 0; i < NcI * NcJ * dE * dE; i++)
2920: if (g3[i]) g3_i += NbI * NbJ;
2921: if (g0_i == NbI * NbJ * NcI * NcJ) g0_i = 0;
2922: if (g1_i == NbI * NbJ * NcI * NcJ * dE) g1_i = 0;
2923: if (g2_i == NbI * NbJ * NcI * NcJ * dE) g2_i = 0;
2924: if (g3_i == NbI * NbJ * NcI * NcJ * dE * dE) g3_i = 0;
2925: has_g0 = g0_i ? PETSC_TRUE : PETSC_FALSE;
2926: has_g1 = g1_i ? PETSC_TRUE : PETSC_FALSE;
2927: has_g2 = g2_i ? PETSC_TRUE : PETSC_FALSE;
2928: has_g3 = g3_i ? PETSC_TRUE : PETSC_FALSE;
2929: if (has_g0) PetscCall(PetscMalloc1(4 * g0_i, &g0_idxs));
2930: if (has_g1) PetscCall(PetscMalloc1(4 * g1_i, &g1_idxs));
2931: if (has_g2) PetscCall(PetscMalloc1(4 * g2_i, &g2_idxs));
2932: if (has_g3) PetscCall(PetscMalloc1(4 * g3_i, &g3_idxs));
2933: g0_i = g1_i = g2_i = g3_i = 0;
2935: for (PetscInt f = 0; f < NbI; ++f) {
2936: const PetscInt i = offsetI + f; /* Element matrix row */
2937: for (PetscInt fc = 0; fc < NcI; ++fc) {
2938: const PetscInt fidx = f * NcI + fc; /* Test function basis index */
2940: for (PetscInt g = 0; g < NbJ; ++g) {
2941: const PetscInt j = offsetJ + g; /* Element matrix column */
2942: const PetscInt fOff = i * totDim + j;
2943: for (PetscInt gc = 0; gc < NcJ; ++gc) {
2944: const PetscInt gidx = g * NcJ + gc; /* Trial function basis index */
2946: if (has_g0) {
2947: if (g0[fc * NcJ + gc]) {
2948: g0_idxs[4 * g0_i + 0] = fidx;
2949: g0_idxs[4 * g0_i + 1] = fc * NcJ + gc;
2950: g0_idxs[4 * g0_i + 2] = gidx;
2951: g0_idxs[4 * g0_i + 3] = fOff;
2952: g0_i++;
2953: }
2954: }
2956: for (PetscInt df = 0; df < dE; ++df) {
2957: if (has_g1) {
2958: if (g1[(fc * NcJ + gc) * dE + df]) {
2959: g1_idxs[4 * g1_i + 0] = fidx;
2960: g1_idxs[4 * g1_i + 1] = (fc * NcJ + gc) * dE + df;
2961: g1_idxs[4 * g1_i + 2] = gidx * dE + df;
2962: g1_idxs[4 * g1_i + 3] = fOff;
2963: g1_i++;
2964: }
2965: }
2966: if (has_g2) {
2967: if (g2[(fc * NcJ + gc) * dE + df]) {
2968: g2_idxs[4 * g2_i + 0] = fidx * dE + df;
2969: g2_idxs[4 * g2_i + 1] = (fc * NcJ + gc) * dE + df;
2970: g2_idxs[4 * g2_i + 2] = gidx;
2971: g2_idxs[4 * g2_i + 3] = fOff;
2972: g2_i++;
2973: }
2974: }
2975: if (has_g3) {
2976: for (PetscInt dg = 0; dg < dE; ++dg) {
2977: if (g3[((fc * NcJ + gc) * dE + df) * dE + dg]) {
2978: g3_idxs[4 * g3_i + 0] = fidx * dE + df;
2979: g3_idxs[4 * g3_i + 1] = ((fc * NcJ + gc) * dE + df) * dE + dg;
2980: g3_idxs[4 * g3_i + 2] = gidx * dE + dg;
2981: g3_idxs[4 * g3_i + 3] = fOff;
2982: g3_i++;
2983: }
2984: }
2985: }
2986: }
2987: }
2988: }
2989: }
2990: }
2991: *n_g0 = g0_i;
2992: *n_g1 = g1_i;
2993: *n_g2 = g2_i;
2994: *n_g3 = g3_i;
2996: *g0_idxs_out = g0_idxs;
2997: *g1_idxs_out = g1_idxs;
2998: *g2_idxs_out = g2_idxs;
2999: *g3_idxs_out = g3_idxs;
3000: PetscFunctionReturn(PETSC_SUCCESS);
3001: }
3003: //example HOW TO USE
3004: for (PetscInt i = 0; i < g0_sparse_n; i++) {
3005: PetscInt bM = g0_sparse_idxs[4 * i + 0];
3006: PetscInt bN = g0_sparse_idxs[4 * i + 1];
3007: PetscInt bK = g0_sparse_idxs[4 * i + 2];
3008: PetscInt bO = g0_sparse_idxs[4 * i + 3];
3009: elemMat[bO] += tmpBasisI[bM] * g0[bN] * tmpBasisJ[bK];
3010: }
3011: #endif