Actual source code: febasic.c
1: #include <petsc/private/petscfeimpl.h>
2: #include <petscblaslapack.h>
4: static PetscErrorCode PetscFEDestroy_Basic(PetscFE fem)
5: {
6: PetscFE_Basic *b = (PetscFE_Basic *)fem->data;
8: PetscFunctionBegin;
9: PetscCall(PetscFree(b));
10: PetscFunctionReturn(PETSC_SUCCESS);
11: }
13: static PetscErrorCode PetscFEView_Basic_Ascii(PetscFE fe, PetscViewer v)
14: {
15: PetscInt dim, Nc;
16: PetscSpace basis = NULL;
17: PetscDualSpace dual = NULL;
18: PetscQuadrature quad = NULL;
20: PetscFunctionBegin;
21: PetscCall(PetscFEGetSpatialDimension(fe, &dim));
22: PetscCall(PetscFEGetNumComponents(fe, &Nc));
23: PetscCall(PetscFEGetBasisSpace(fe, &basis));
24: PetscCall(PetscFEGetDualSpace(fe, &dual));
25: PetscCall(PetscFEGetQuadrature(fe, &quad));
26: PetscCall(PetscViewerASCIIPushTab(v));
27: PetscCall(PetscViewerASCIIPrintf(v, "Basic Finite Element in %" PetscInt_FMT " dimensions with %" PetscInt_FMT " components\n", dim, Nc));
28: if (basis) PetscCall(PetscSpaceView(basis, v));
29: if (dual) PetscCall(PetscDualSpaceView(dual, v));
30: if (quad) PetscCall(PetscQuadratureView(quad, v));
31: PetscCall(PetscViewerASCIIPopTab(v));
32: PetscFunctionReturn(PETSC_SUCCESS);
33: }
35: static PetscErrorCode PetscFEView_Basic(PetscFE fe, PetscViewer v)
36: {
37: PetscBool isascii;
39: PetscFunctionBegin;
40: PetscCall(PetscObjectTypeCompare((PetscObject)v, PETSCVIEWERASCII, &isascii));
41: if (isascii) PetscCall(PetscFEView_Basic_Ascii(fe, v));
42: PetscFunctionReturn(PETSC_SUCCESS);
43: }
45: /* Construct the change of basis from prime basis to nodal basis */
46: PETSC_INTERN PetscErrorCode PetscFESetUp_Basic(PetscFE fem)
47: {
48: PetscReal *work;
49: PetscBLASInt *pivots;
50: PetscBLASInt n, info;
51: PetscInt pdim;
53: PetscFunctionBegin;
54: PetscCall(PetscDualSpaceGetDimension(fem->dualSpace, &pdim));
55: PetscCall(PetscMalloc1(pdim * pdim, &fem->invV));
56: for (PetscInt j = 0; j < pdim; ++j) {
57: PetscReal *Bf;
58: PetscQuadrature f;
59: const PetscReal *points, *weights;
60: PetscInt Nc, Nq, q, k, c;
62: PetscCall(PetscDualSpaceGetFunctional(fem->dualSpace, j, &f));
63: PetscCall(PetscQuadratureGetData(f, NULL, &Nc, &Nq, &points, &weights));
64: PetscCall(PetscMalloc1(Nc * Nq * pdim, &Bf));
65: PetscCall(PetscSpaceEvaluate(fem->basisSpace, Nq, points, Bf, NULL, NULL));
66: for (k = 0; k < pdim; ++k) {
67: /* V_{jk} = n_j(\phi_k) = \int \phi_k(x) n_j(x) dx */
68: fem->invV[j * pdim + k] = 0.0;
70: for (q = 0; q < Nq; ++q) {
71: for (c = 0; c < Nc; ++c) fem->invV[j * pdim + k] += Bf[(q * pdim + k) * Nc + c] * weights[q * Nc + c];
72: }
73: }
74: PetscCall(PetscFree(Bf));
75: }
77: PetscCall(PetscMalloc2(pdim, &pivots, pdim, &work));
78: PetscCall(PetscBLASIntCast(pdim, &n));
79: PetscCallBLAS("LAPACKgetrf", LAPACKREALgetrf_(&n, &n, fem->invV, &n, pivots, &info));
80: PetscCheck(!info, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error returned from LAPACKgetrf %" PetscBLASInt_FMT, info);
81: PetscCallBLAS("LAPACKgetri", LAPACKREALgetri_(&n, fem->invV, &n, pivots, work, &n, &info));
82: PetscCheck(!info, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error returned from LAPACKgetri %" PetscBLASInt_FMT, info);
83: PetscCall(PetscFree2(pivots, work));
84: PetscFunctionReturn(PETSC_SUCCESS);
85: }
87: PetscErrorCode PetscFEGetDimension_Basic(PetscFE fem, PetscInt *dim)
88: {
89: PetscFunctionBegin;
90: PetscCall(PetscDualSpaceGetDimension(fem->dualSpace, dim));
91: PetscFunctionReturn(PETSC_SUCCESS);
92: }
94: /* Tensor contraction on the middle index,
95: * C[m,n,p] = A[m,k,p] * B[k,n]
96: * where all matrices use C-style ordering.
97: */
98: static PetscErrorCode TensorContract_Private(PetscInt m, PetscInt n, PetscInt p, PetscInt k, const PetscReal *A, const PetscReal *B, PetscReal *C)
99: {
100: PetscFunctionBegin;
101: PetscCheck(n && p, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Empty tensor is not allowed %" PetscInt_FMT " %" PetscInt_FMT, n, p);
102: for (PetscInt i = 0; i < m; i++) {
103: PetscBLASInt n_, p_, k_, lda, ldb, ldc;
104: PetscReal one = 1, zero = 0;
105: /* Taking contiguous submatrices, we wish to comput c[n,p] = a[k,p] * B[k,n]
106: * or, in Fortran ordering, c(p,n) = a(p,k) * B(n,k)
107: */
108: PetscCall(PetscBLASIntCast(n, &n_));
109: PetscCall(PetscBLASIntCast(p, &p_));
110: PetscCall(PetscBLASIntCast(k, &k_));
111: lda = p_;
112: ldb = n_;
113: ldc = p_;
114: PetscCallBLAS("BLASgemm", BLASREALgemm_("N", "T", &p_, &n_, &k_, &one, A + i * k * p, &lda, B, &ldb, &zero, C + i * n * p, &ldc));
115: }
116: PetscCall(PetscLogFlops(2. * m * n * p * k));
117: PetscFunctionReturn(PETSC_SUCCESS);
118: }
120: PETSC_INTERN PetscErrorCode PetscFEComputeTabulation_Basic(PetscFE fem, PetscInt npoints, const PetscReal points[], PetscInt K, PetscTabulation T)
121: {
122: DM dm;
123: PetscInt pdim; /* Dimension of FE space P */
124: PetscInt dim; /* Spatial dimension */
125: PetscInt Nc; /* Field components */
126: PetscReal *B = K >= 0 ? T->T[0] : NULL;
127: PetscReal *D = K >= 1 ? T->T[1] : NULL;
128: PetscReal *H = K >= 2 ? T->T[2] : NULL;
129: PetscReal *tmpB = NULL, *tmpD = NULL, *tmpH = NULL;
131: PetscFunctionBegin;
132: PetscCall(PetscDualSpaceGetDM(fem->dualSpace, &dm));
133: PetscCall(DMGetDimension(dm, &dim));
134: PetscCall(PetscDualSpaceGetDimension(fem->dualSpace, &pdim));
135: PetscCall(PetscFEGetNumComponents(fem, &Nc));
136: /* Evaluate the prime basis functions at all points */
137: if (K >= 0) PetscCall(DMGetWorkArray(dm, npoints * pdim * Nc, MPIU_REAL, &tmpB));
138: if (K >= 1) PetscCall(DMGetWorkArray(dm, npoints * pdim * Nc * dim, MPIU_REAL, &tmpD));
139: if (K >= 2) PetscCall(DMGetWorkArray(dm, npoints * pdim * Nc * dim * dim, MPIU_REAL, &tmpH));
140: PetscCall(PetscSpaceEvaluate(fem->basisSpace, npoints, points, tmpB, tmpD, tmpH));
141: /* Translate from prime to nodal basis */
142: if (B) {
143: /* B[npoints, nodes, Nc] = tmpB[npoints, prime, Nc] * invV[prime, nodes] */
144: PetscCall(TensorContract_Private(npoints, pdim, Nc, pdim, tmpB, fem->invV, B));
145: }
146: if (D && dim) {
147: /* D[npoints, nodes, Nc, dim] = tmpD[npoints, prime, Nc, dim] * invV[prime, nodes] */
148: PetscCall(TensorContract_Private(npoints, pdim, Nc * dim, pdim, tmpD, fem->invV, D));
149: }
150: if (H && dim) {
151: /* H[npoints, nodes, Nc, dim, dim] = tmpH[npoints, prime, Nc, dim, dim] * invV[prime, nodes] */
152: PetscCall(TensorContract_Private(npoints, pdim, Nc * dim * dim, pdim, tmpH, fem->invV, H));
153: }
154: if (K >= 0) PetscCall(DMRestoreWorkArray(dm, npoints * pdim * Nc, MPIU_REAL, &tmpB));
155: if (K >= 1) PetscCall(DMRestoreWorkArray(dm, npoints * pdim * Nc * dim, MPIU_REAL, &tmpD));
156: if (K >= 2) PetscCall(DMRestoreWorkArray(dm, npoints * pdim * Nc * dim * dim, MPIU_REAL, &tmpH));
157: PetscFunctionReturn(PETSC_SUCCESS);
158: }
160: PETSC_INTERN PetscErrorCode PetscFEIntegrate_Basic(PetscDS ds, PetscInt field, PetscInt Ne, PetscFEGeom *cgeom, const PetscScalar coefficients[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscScalar integral[])
161: {
162: const PetscInt debug = ds->printIntegrate;
163: PetscFE fe;
164: PetscPointFn *obj_func;
165: PetscQuadrature quad;
166: PetscTabulation *T, *TAux = NULL;
167: PetscScalar *u, *u_x, *a, *a_x;
168: const PetscScalar *constants;
169: PetscReal *x, cellScale;
170: PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL;
171: PetscInt dim, dE, Np, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, cOffset = 0, cOffsetAux = 0, e;
172: PetscBool isAffine;
173: const PetscReal *quadPoints, *quadWeights;
174: PetscInt qNc, Nq, q;
176: PetscFunctionBegin;
177: PetscCall(PetscDSGetObjective(ds, field, &obj_func));
178: if (!obj_func) PetscFunctionReturn(PETSC_SUCCESS);
179: PetscCall(PetscDSGetDiscretization(ds, field, (PetscObject *)&fe));
180: PetscCall(PetscFEGetSpatialDimension(fe, &dim));
181: cellScale = (PetscReal)PetscPowInt(2, dim);
182: PetscCall(PetscFEGetQuadrature(fe, &quad));
183: PetscCall(PetscDSGetNumFields(ds, &Nf));
184: PetscCall(PetscDSGetTotalDimension(ds, &totDim));
185: PetscCall(PetscDSGetComponentOffsets(ds, &uOff));
186: PetscCall(PetscDSGetComponentDerivativeOffsets(ds, &uOff_x));
187: PetscCall(PetscDSGetTabulation(ds, &T));
188: PetscCall(PetscDSGetEvaluationArrays(ds, &u, NULL, &u_x));
189: PetscCall(PetscDSGetWorkspace(ds, &x, NULL, NULL, NULL, NULL));
190: PetscCall(PetscDSSetIntegrationParameters(ds, field, PETSC_DETERMINE));
191: PetscCall(PetscDSGetConstants(ds, &numConstants, &constants));
192: if (dsAux) {
193: PetscCall(PetscDSGetNumFields(dsAux, &NfAux));
194: PetscCall(PetscDSGetTotalDimension(dsAux, &totDimAux));
195: PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff));
196: PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x));
197: PetscCall(PetscDSGetTabulation(dsAux, &TAux));
198: PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x));
199: PetscCheck(T[0]->Np == TAux[0]->Np, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of tabulation points %" PetscInt_FMT " != %" PetscInt_FMT " number of auxiliary tabulation points", T[0]->Np, TAux[0]->Np);
200: }
201: PetscCall(PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights));
202: PetscCheck(qNc == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %" PetscInt_FMT " components", qNc);
203: Np = cgeom->numPoints;
204: dE = cgeom->dimEmbed;
205: isAffine = cgeom->isAffine;
206: for (e = 0; e < Ne; ++e) {
207: PetscFEGeom fegeom;
209: fegeom.dim = cgeom->dim;
210: fegeom.dimEmbed = cgeom->dimEmbed;
211: fegeom.xi = NULL;
212: if (isAffine) {
213: fegeom.v = x;
214: fegeom.xi = cgeom->xi;
215: fegeom.J = &cgeom->J[e * Np * dE * dE];
216: fegeom.invJ = &cgeom->invJ[e * Np * dE * dE];
217: fegeom.detJ = &cgeom->detJ[e * Np];
218: } else fegeom.xi = NULL;
219: for (q = 0; q < Nq; ++q) {
220: PetscScalar integrand = 0.;
221: PetscReal w;
223: if (isAffine) {
224: CoordinatesRefToReal(dE, dim, fegeom.xi, &cgeom->v[e * Np * dE], fegeom.J, &quadPoints[q * dim], x);
225: } else {
226: fegeom.v = &cgeom->v[(e * Np + q) * dE];
227: fegeom.J = &cgeom->J[(e * Np + q) * dE * dE];
228: fegeom.invJ = &cgeom->invJ[(e * Np + q) * dE * dE];
229: fegeom.detJ = &cgeom->detJ[e * Np + q];
230: }
231: PetscCall(PetscDSSetCellParameters(ds, fegeom.detJ[0] * cellScale));
232: w = fegeom.detJ[0] * quadWeights[q];
233: if (debug > 1 && q < Np) {
234: PetscCall(PetscPrintf(PETSC_COMM_SELF, " detJ: %g\n", (double)fegeom.detJ[0]));
235: #if !defined(PETSC_USE_COMPLEX)
236: PetscCall(DMPrintCellMatrix(e, "invJ", dim, dim, fegeom.invJ));
237: #endif
238: }
239: if (debug) PetscCall(PetscPrintf(PETSC_COMM_SELF, " quad point %" PetscInt_FMT "\n", q));
240: PetscCall(PetscFEEvaluateFieldJets_Internal(ds, Nf, 0, q, T, &fegeom, &coefficients[cOffset], NULL, u, u_x, NULL));
241: if (dsAux) PetscCall(PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, 0, q, TAux, &fegeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL));
242: obj_func(dim, Nf, NfAux, uOff, uOff_x, u, NULL, u_x, aOff, aOff_x, a, NULL, a_x, 0.0, fegeom.v, numConstants, constants, &integrand);
243: integrand *= w;
244: integral[e * Nf + field] += integrand;
245: }
246: if (debug > 1) PetscCall(PetscPrintf(PETSC_COMM_SELF, " Element Field %" PetscInt_FMT " integral: %g\n", Nf, (double)PetscRealPart(integral[e * Nf + field])));
247: cOffset += totDim;
248: cOffsetAux += totDimAux;
249: }
250: PetscFunctionReturn(PETSC_SUCCESS);
251: }
253: PETSC_INTERN PetscErrorCode PetscFEIntegrateBd_Basic(PetscDS ds, PetscInt field, PetscBdPointFn *obj_func, PetscInt Ne, PetscFEGeom *fgeom, const PetscScalar coefficients[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscScalar integral[])
254: {
255: const PetscInt debug = ds->printIntegrate;
256: PetscFE fe;
257: PetscQuadrature quad;
258: PetscTabulation *Tf, *TfAux = NULL;
259: PetscScalar *u, *u_x, *a, *a_x, *basisReal, *basisDerReal;
260: const PetscScalar *constants;
261: PetscReal *x, cellScale;
262: PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL;
263: PetscBool isAffine, auxOnBd;
264: const PetscReal *quadPoints, *quadWeights;
265: PetscInt qNc, Nq, q, Np, dE;
266: PetscInt dim, dimAux, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, cOffset = 0, cOffsetAux = 0, e;
268: PetscFunctionBegin;
269: if (!obj_func) PetscFunctionReturn(PETSC_SUCCESS);
270: PetscCall(PetscDSGetDiscretization(ds, field, (PetscObject *)&fe));
271: PetscCall(PetscFEGetSpatialDimension(fe, &dim));
272: cellScale = (PetscReal)PetscPowInt(2, dim);
273: PetscCall(PetscFEGetFaceQuadrature(fe, &quad));
274: PetscCall(PetscDSGetNumFields(ds, &Nf));
275: PetscCall(PetscDSGetTotalDimension(ds, &totDim));
276: PetscCall(PetscDSGetComponentOffsets(ds, &uOff));
277: PetscCall(PetscDSGetComponentDerivativeOffsets(ds, &uOff_x));
278: PetscCall(PetscDSGetEvaluationArrays(ds, &u, NULL, &u_x));
279: PetscCall(PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, NULL, NULL));
280: PetscCall(PetscDSGetFaceTabulation(ds, &Tf));
281: PetscCall(PetscDSSetIntegrationParameters(ds, field, PETSC_DETERMINE));
282: PetscCall(PetscDSGetConstants(ds, &numConstants, &constants));
283: if (dsAux) {
284: PetscCall(PetscDSGetSpatialDimension(dsAux, &dimAux));
285: PetscCall(PetscDSGetNumFields(dsAux, &NfAux));
286: PetscCall(PetscDSGetTotalDimension(dsAux, &totDimAux));
287: PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff));
288: PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x));
289: PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x));
290: auxOnBd = dimAux < dim ? PETSC_TRUE : PETSC_FALSE;
291: if (auxOnBd) PetscCall(PetscDSGetTabulation(dsAux, &TfAux));
292: else PetscCall(PetscDSGetFaceTabulation(dsAux, &TfAux));
293: PetscCheck(Tf[0]->Np == TfAux[0]->Np, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of tabulation points %" PetscInt_FMT " != %" PetscInt_FMT " number of auxiliary tabulation points", Tf[0]->Np, TfAux[0]->Np);
294: }
295: PetscCall(PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights));
296: PetscCheck(qNc == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %" PetscInt_FMT " components", qNc);
297: if (debug > 1) PetscCall(PetscPrintf(PETSC_COMM_SELF, "Field: %" PetscInt_FMT " Nface: %" PetscInt_FMT " Nq: %" PetscInt_FMT "\n", field, Ne, Nq));
298: Np = fgeom->numPoints;
299: dE = fgeom->dimEmbed;
300: isAffine = fgeom->isAffine;
301: for (e = 0; e < Ne; ++e) {
302: PetscFEGeom fegeom, cgeom;
303: const PetscInt face = fgeom->face[e][0]; /* Local face number in cell */
304: fegeom.n = NULL;
305: fegeom.v = NULL;
306: fegeom.xi = NULL;
307: fegeom.J = NULL;
308: fegeom.invJ = NULL;
309: fegeom.detJ = NULL;
310: fegeom.dim = fgeom->dim;
311: fegeom.dimEmbed = fgeom->dimEmbed;
312: cgeom.dim = fgeom->dim;
313: cgeom.dimEmbed = fgeom->dimEmbed;
314: if (isAffine) {
315: fegeom.v = x;
316: fegeom.xi = fgeom->xi;
317: fegeom.J = &fgeom->J[e * Np * dE * dE];
318: fegeom.invJ = &fgeom->invJ[e * Np * dE * dE];
319: fegeom.detJ = &fgeom->detJ[e * Np];
320: fegeom.n = &fgeom->n[e * Np * dE];
322: cgeom.J = &fgeom->suppJ[0][e * Np * dE * dE];
323: cgeom.invJ = &fgeom->suppInvJ[0][e * Np * dE * dE];
324: cgeom.detJ = &fgeom->suppDetJ[0][e * Np];
325: } else fegeom.xi = NULL;
326: for (q = 0; q < Nq; ++q) {
327: PetscScalar integrand = 0.;
328: PetscReal w;
330: if (isAffine) {
331: CoordinatesRefToReal(dE, dim - 1, fegeom.xi, &fgeom->v[e * Np * dE], fegeom.J, &quadPoints[q * (dim - 1)], x);
332: } else {
333: fegeom.v = &fgeom->v[(e * Np + q) * dE];
334: fegeom.J = &fgeom->J[(e * Np + q) * dE * dE];
335: fegeom.invJ = &fgeom->invJ[(e * Np + q) * dE * dE];
336: fegeom.detJ = &fgeom->detJ[e * Np + q];
337: fegeom.n = &fgeom->n[(e * Np + q) * dE];
339: cgeom.J = &fgeom->suppJ[0][(e * Np + q) * dE * dE];
340: cgeom.invJ = &fgeom->suppInvJ[0][(e * Np + q) * dE * dE];
341: cgeom.detJ = &fgeom->suppDetJ[0][e * Np + q];
342: }
343: PetscCall(PetscDSSetCellParameters(ds, fegeom.detJ[0] * cellScale));
344: w = fegeom.detJ[0] * quadWeights[q];
345: if (debug > 1 && q < Np) {
346: PetscCall(PetscPrintf(PETSC_COMM_SELF, " detJ: %g\n", (double)fegeom.detJ[0]));
347: #if !defined(PETSC_USE_COMPLEX)
348: PetscCall(DMPrintCellMatrix(e, "invJ", dim, dim, fegeom.invJ));
349: #endif
350: }
351: if (debug > 1) PetscCall(PetscPrintf(PETSC_COMM_SELF, " quad point %" PetscInt_FMT "\n", q));
352: if (debug > 3) {
353: PetscCall(PetscPrintf(PETSC_COMM_SELF, " x_q ("));
354: for (PetscInt d = 0; d < dE; ++d) {
355: if (d) PetscCall(PetscPrintf(PETSC_COMM_SELF, ", "));
356: PetscCall(PetscPrintf(PETSC_COMM_SELF, "%g", (double)fegeom.v[d]));
357: }
358: PetscCall(PetscPrintf(PETSC_COMM_SELF, ")\n"));
359: PetscCall(PetscPrintf(PETSC_COMM_SELF, " n_q ("));
360: for (PetscInt d = 0; d < dE; ++d) {
361: if (d) PetscCall(PetscPrintf(PETSC_COMM_SELF, ", "));
362: PetscCall(PetscPrintf(PETSC_COMM_SELF, "%g", (double)fegeom.n[d]));
363: }
364: PetscCall(PetscPrintf(PETSC_COMM_SELF, ")\n"));
365: for (PetscInt f = 0; f < Nf; ++f) {
366: PetscCall(PetscPrintf(PETSC_COMM_SELF, " u_%" PetscInt_FMT " (", f));
367: for (PetscInt c = 0; c < uOff[f + 1] - uOff[f]; ++c) {
368: if (c) PetscCall(PetscPrintf(PETSC_COMM_SELF, ", "));
369: PetscCall(PetscPrintf(PETSC_COMM_SELF, "%g", (double)PetscRealPart(u[uOff[f] + c])));
370: }
371: PetscCall(PetscPrintf(PETSC_COMM_SELF, ")\n"));
372: }
373: }
374: PetscCall(PetscFEEvaluateFieldJets_Internal(ds, Nf, face, q, Tf, &cgeom, &coefficients[cOffset], NULL, u, u_x, NULL));
375: if (dsAux) PetscCall(PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, face, q, TfAux, &cgeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL));
376: obj_func(dim, Nf, NfAux, uOff, uOff_x, u, NULL, u_x, aOff, aOff_x, a, NULL, a_x, 0.0, fegeom.v, fegeom.n, numConstants, constants, &integrand);
377: integrand *= w;
378: integral[e * Nf + field] += integrand;
379: if (debug > 1) PetscCall(PetscPrintf(PETSC_COMM_SELF, " int: %g tot: %g\n", (double)PetscRealPart(integrand), (double)PetscRealPart(integral[e * Nf + field])));
380: }
381: cOffset += totDim;
382: cOffsetAux += totDimAux;
383: }
384: PetscFunctionReturn(PETSC_SUCCESS);
385: }
387: PetscErrorCode PetscFEIntegrateResidual_Basic(PetscDS ds, PetscFormKey key, PetscInt Ne, PetscFEGeom *cgeom, const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscReal t, PetscScalar elemVec[])
388: {
389: const PetscInt debug = ds->printIntegrate;
390: const PetscInt field = key.field;
391: PetscFE fe;
392: PetscWeakForm wf;
393: PetscInt n0, n1, i;
394: PetscPointFn **f0_func, **f1_func;
395: PetscQuadrature quad;
396: PetscTabulation *T, *TAux = NULL;
397: PetscScalar *f0, *f1, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal;
398: const PetscScalar *constants;
399: PetscReal *x, cellScale;
400: PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL;
401: PetscInt dim, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, cOffset = 0, cOffsetAux = 0, fOffset, e;
402: const PetscReal *quadPoints, *quadWeights;
403: PetscInt qdim, qNc, Nq, q, dE;
405: PetscFunctionBegin;
406: PetscCall(PetscDSGetDiscretization(ds, field, (PetscObject *)&fe));
407: PetscCall(PetscFEGetSpatialDimension(fe, &dim));
408: cellScale = (PetscReal)PetscPowInt(2, dim);
409: PetscCall(PetscFEGetQuadrature(fe, &quad));
410: PetscCall(PetscDSGetNumFields(ds, &Nf));
411: PetscCall(PetscDSGetTotalDimension(ds, &totDim));
412: PetscCall(PetscDSGetComponentOffsets(ds, &uOff));
413: PetscCall(PetscDSGetComponentDerivativeOffsets(ds, &uOff_x));
414: PetscCall(PetscDSGetFieldOffset(ds, field, &fOffset));
415: PetscCall(PetscDSGetWeakForm(ds, &wf));
416: PetscCall(PetscWeakFormGetResidual(wf, key.label, key.value, key.field, key.part, &n0, &f0_func, &n1, &f1_func));
417: if (!n0 && !n1) PetscFunctionReturn(PETSC_SUCCESS);
418: PetscCall(PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x));
419: PetscCall(PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, NULL, NULL));
420: PetscCall(PetscDSGetWeakFormArrays(ds, &f0, &f1, NULL, NULL, NULL, NULL));
421: PetscCall(PetscDSGetTabulation(ds, &T));
422: PetscCall(PetscDSSetIntegrationParameters(ds, field, PETSC_DETERMINE));
423: PetscCall(PetscDSGetConstants(ds, &numConstants, &constants));
424: if (dsAux) {
425: PetscCall(PetscDSGetNumFields(dsAux, &NfAux));
426: PetscCall(PetscDSGetTotalDimension(dsAux, &totDimAux));
427: PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff));
428: PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x));
429: PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x));
430: PetscCall(PetscDSGetTabulation(dsAux, &TAux));
431: PetscCheck(T[0]->Np == TAux[0]->Np, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of tabulation points %" PetscInt_FMT " != %" PetscInt_FMT " number of auxiliary tabulation points", T[0]->Np, TAux[0]->Np);
432: }
433: PetscCall(PetscQuadratureGetData(quad, &qdim, &qNc, &Nq, &quadPoints, &quadWeights));
434: PetscCheck(qNc == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %" PetscInt_FMT " components", qNc);
435: dE = cgeom->dimEmbed;
436: PetscCheck(cgeom->dim == qdim, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "FEGeom dim %" PetscInt_FMT " != %" PetscInt_FMT " quadrature dim", cgeom->dim, qdim);
437: for (e = 0; e < Ne; ++e) {
438: PetscFEGeom fegeom;
440: fegeom.v = x; /* workspace */
441: PetscCall(PetscArrayzero(f0, Nq * T[field]->Nc));
442: PetscCall(PetscArrayzero(f1, Nq * T[field]->Nc * dE));
443: for (q = 0; q < Nq; ++q) {
444: PetscReal w;
445: PetscInt c, d;
447: PetscCall(PetscFEGeomGetPoint(cgeom, e, q, &quadPoints[q * cgeom->dim], &fegeom));
448: PetscCall(PetscDSSetCellParameters(ds, fegeom.detJ[0] * cellScale));
449: w = fegeom.detJ[0] * quadWeights[q];
450: if (debug > 1 && q < cgeom->numPoints) {
451: PetscCall(PetscPrintf(PETSC_COMM_SELF, " detJ: %g\n", (double)fegeom.detJ[0]));
452: #if !defined(PETSC_USE_COMPLEX)
453: PetscCall(DMPrintCellMatrix(e, "invJ", dE, dE, fegeom.invJ));
454: #endif
455: }
456: PetscCall(PetscFEEvaluateFieldJets_Internal(ds, Nf, 0, q, T, &fegeom, &coefficients[cOffset], PetscSafePointerPlusOffset(coefficients_t, cOffset), u, u_x, u_t));
457: if (dsAux) PetscCall(PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, 0, q, TAux, &fegeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL));
458: for (i = 0; i < n0; ++i) f0_func[i](dE, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, fegeom.v, numConstants, constants, &f0[q * T[field]->Nc]);
459: for (c = 0; c < T[field]->Nc; ++c) f0[q * T[field]->Nc + c] *= w;
460: for (i = 0; i < n1; ++i) f1_func[i](dE, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, fegeom.v, numConstants, constants, &f1[q * T[field]->Nc * dE]);
461: for (c = 0; c < T[field]->Nc; ++c)
462: for (d = 0; d < dE; ++d) f1[(q * T[field]->Nc + c) * dE + d] *= w;
463: if (debug) {
464: // LCOV_EXCL_START
465: PetscCall(PetscPrintf(PETSC_COMM_SELF, " quad point %" PetscInt_FMT " wt %g x:", q, (double)quadWeights[q]));
466: for (c = 0; c < dE; ++c) PetscCall(PetscPrintf(PETSC_COMM_SELF, " %g", (double)fegeom.v[c]));
467: PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n"));
468: if (debug > 2) {
469: PetscCall(PetscPrintf(PETSC_COMM_SELF, " field %" PetscInt_FMT ":", field));
470: for (c = 0; c < T[field]->Nc; ++c) PetscCall(PetscPrintf(PETSC_COMM_SELF, " %g", (double)PetscRealPart(u[uOff[field] + c])));
471: PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n"));
472: PetscCall(PetscPrintf(PETSC_COMM_SELF, " field der %" PetscInt_FMT ":", field));
473: for (c = 0; c < T[field]->Nc * dE; ++c) PetscCall(PetscPrintf(PETSC_COMM_SELF, " %g", (double)PetscRealPart(u_x[uOff[field] + c])));
474: PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n"));
475: PetscCall(PetscPrintf(PETSC_COMM_SELF, " resid %" PetscInt_FMT ":", field));
476: for (c = 0; c < T[field]->Nc; ++c) PetscCall(PetscPrintf(PETSC_COMM_SELF, " %g", (double)PetscRealPart(f0[q * T[field]->Nc + c])));
477: PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n"));
478: PetscCall(PetscPrintf(PETSC_COMM_SELF, " res der %" PetscInt_FMT ":", field));
479: for (c = 0; c < T[field]->Nc; ++c) {
480: for (d = 0; d < dE; ++d) PetscCall(PetscPrintf(PETSC_COMM_SELF, " %g", (double)PetscRealPart(f1[(q * T[field]->Nc + c) * dE + d])));
481: }
482: PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n"));
483: }
484: // LCOV_EXCL_STOP
485: }
486: }
487: PetscCall(PetscFEUpdateElementVec_Internal(fe, T[field], 0, basisReal, basisDerReal, e, cgeom, f0, f1, &elemVec[cOffset + fOffset]));
488: cOffset += totDim;
489: cOffsetAux += totDimAux;
490: }
491: PetscFunctionReturn(PETSC_SUCCESS);
492: }
494: PetscErrorCode PetscFEIntegrateBdResidual_Basic(PetscDS ds, PetscWeakForm wf, PetscFormKey key, PetscInt Ne, PetscFEGeom *fgeom, const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscReal t, PetscScalar elemVec[])
495: {
496: const PetscInt debug = ds->printIntegrate;
497: const PetscInt field = key.field;
498: PetscFE fe;
499: PetscInt n0, n1, i;
500: PetscBdPointFn **f0_func, **f1_func;
501: PetscQuadrature quad;
502: PetscTabulation *Tf, *TfAux = NULL;
503: PetscScalar *f0, *f1, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal;
504: const PetscScalar *constants;
505: PetscReal *x, cellScale;
506: PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL;
507: PetscInt dim, dimAux, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, cOffset = 0, cOffsetAux = 0, fOffset, e, NcI;
508: PetscBool auxOnBd = PETSC_FALSE;
509: const PetscReal *quadPoints, *quadWeights;
510: PetscInt qdim, qNc, Nq, q, dE;
512: PetscFunctionBegin;
513: PetscCall(PetscDSGetDiscretization(ds, field, (PetscObject *)&fe));
514: PetscCall(PetscFEGetSpatialDimension(fe, &dim));
515: cellScale = (PetscReal)PetscPowInt(2, dim);
516: PetscCall(PetscFEGetFaceQuadrature(fe, &quad));
517: PetscCall(PetscDSGetNumFields(ds, &Nf));
518: PetscCall(PetscDSGetTotalDimension(ds, &totDim));
519: PetscCall(PetscDSGetComponentOffsets(ds, &uOff));
520: PetscCall(PetscDSGetComponentDerivativeOffsets(ds, &uOff_x));
521: PetscCall(PetscDSGetFieldOffset(ds, field, &fOffset));
522: PetscCall(PetscWeakFormGetBdResidual(wf, key.label, key.value, key.field, key.part, &n0, &f0_func, &n1, &f1_func));
523: if (!n0 && !n1) PetscFunctionReturn(PETSC_SUCCESS);
524: PetscCall(PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x));
525: PetscCall(PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, NULL, NULL));
526: PetscCall(PetscDSGetWeakFormArrays(ds, &f0, &f1, NULL, NULL, NULL, NULL));
527: PetscCall(PetscDSGetFaceTabulation(ds, &Tf));
528: PetscCall(PetscDSSetIntegrationParameters(ds, field, PETSC_DETERMINE));
529: PetscCall(PetscDSGetConstants(ds, &numConstants, &constants));
530: if (dsAux) {
531: PetscCall(PetscDSGetSpatialDimension(dsAux, &dimAux));
532: PetscCall(PetscDSGetNumFields(dsAux, &NfAux));
533: PetscCall(PetscDSGetTotalDimension(dsAux, &totDimAux));
534: PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff));
535: PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x));
536: PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x));
537: auxOnBd = dimAux < dim ? PETSC_TRUE : PETSC_FALSE;
538: if (auxOnBd) PetscCall(PetscDSGetTabulation(dsAux, &TfAux));
539: else PetscCall(PetscDSGetFaceTabulation(dsAux, &TfAux));
540: PetscCheck(Tf[0]->Np == TfAux[0]->Np, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of tabulation points %" PetscInt_FMT " != %" PetscInt_FMT " number of auxiliary tabulation points", Tf[0]->Np, TfAux[0]->Np);
541: }
542: NcI = Tf[field]->Nc;
543: PetscCall(PetscQuadratureGetData(quad, &qdim, &qNc, &Nq, &quadPoints, &quadWeights));
544: PetscCheck(qNc == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %" PetscInt_FMT " components", qNc);
545: dE = fgeom->dimEmbed;
546: /* TODO FIX THIS */
547: fgeom->dim = dim - 1;
548: PetscCheck(fgeom->dim == qdim, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "FEGeom dim %" PetscInt_FMT " != %" PetscInt_FMT " quadrature dim", fgeom->dim, qdim);
549: for (e = 0; e < Ne; ++e) {
550: PetscFEGeom fegeom, cgeom;
551: const PetscInt face = fgeom->face[e][0];
552: const PetscInt ornt = fgeom->face[e][1];
554: fegeom.v = x; /* Workspace */
555: PetscCall(PetscArrayzero(f0, Nq * NcI));
556: PetscCall(PetscArrayzero(f1, Nq * NcI * dE));
557: for (q = 0; q < Nq; ++q) {
558: PetscReal w;
559: PetscInt c;
560: const PetscInt qp = ornt < 0 ? (Nq - 1 - q) : q; /* Map physical quadrature index to tabulation index accounting for face orientation */
562: PetscCall(PetscFEGeomGetPoint(fgeom, e, q, &quadPoints[q * fgeom->dim], &fegeom));
563: PetscCall(PetscFEGeomGetCellPoint(fgeom, e, q, &cgeom));
564: PetscCall(PetscDSSetCellParameters(ds, fegeom.detJ[0] * cellScale));
565: w = fegeom.detJ[0] * quadWeights[q];
566: if (debug > 1) {
567: if ((fgeom->isAffine && q == 0) || !fgeom->isAffine) {
568: PetscCall(PetscPrintf(PETSC_COMM_SELF, " detJ: %g\n", (double)fegeom.detJ[0]));
569: #if !defined(PETSC_USE_COMPLEX)
570: PetscCall(DMPrintCellMatrix(e, "invJ", dim, dim, fegeom.invJ));
571: PetscCall(DMPrintCellVector(e, "n", dim, fegeom.n));
572: #endif
573: }
574: }
575: PetscCall(PetscFEEvaluateFieldJets_Internal(ds, Nf, face, qp, Tf, &cgeom, &coefficients[cOffset], PetscSafePointerPlusOffset(coefficients_t, cOffset), u, u_x, u_t));
576: if (dsAux) PetscCall(PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, auxOnBd ? 0 : face, qp, TfAux, &cgeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL));
577: for (i = 0; i < n0; ++i) f0_func[i](dE, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, fegeom.v, fegeom.n, numConstants, constants, &f0[qp * NcI]);
578: for (c = 0; c < NcI; ++c) f0[qp * NcI + c] *= w;
579: for (i = 0; i < n1; ++i) f1_func[i](dE, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, fegeom.v, fegeom.n, numConstants, constants, &f1[qp * NcI * dE]);
580: for (c = 0; c < NcI; ++c)
581: for (PetscInt d = 0; d < dE; ++d) f1[(qp * NcI + c) * dE + d] *= w;
582: if (debug) {
583: PetscCall(PetscPrintf(PETSC_COMM_SELF, " elem %" PetscInt_FMT " quad point %" PetscInt_FMT "\n", e, q));
584: for (c = 0; c < NcI; ++c) {
585: if (n0) PetscCall(PetscPrintf(PETSC_COMM_SELF, " f0[%" PetscInt_FMT "] %g\n", c, (double)PetscRealPart(f0[qp * NcI + c])));
586: if (n1) {
587: for (PetscInt d = 0; d < dim; ++d) PetscCall(PetscPrintf(PETSC_COMM_SELF, " f1[%" PetscInt_FMT ",%" PetscInt_FMT "] %g", c, d, (double)PetscRealPart(f1[(qp * NcI + c) * dim + d])));
588: PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n"));
589: }
590: }
591: }
592: }
593: PetscCall(PetscFEUpdateElementVec_Internal(fe, Tf[field], face, basisReal, basisDerReal, e, fgeom, f0, f1, &elemVec[cOffset + fOffset]));
594: cOffset += totDim;
595: cOffsetAux += totDimAux;
596: }
597: PetscFunctionReturn(PETSC_SUCCESS);
598: }
600: /*
601: BdIntegral: Operates completely in the embedding dimension. The trick is to have special "face quadrature" so we only integrate over the face, but
602: all transforms operate in the full space and are square.
604: HybridIntegral: The discretization is lower dimensional. That means the transforms are non-square.
605: 1) DMPlexGetCellFields() retrieves from the hybrid cell, so it gets fields from both faces
606: 2) We need to assume that the orientation is 0 for both
607: 3) TODO We need to use a non-square Jacobian for the derivative maps, meaning the embedding dimension has to go to EvaluateFieldJets() and UpdateElementVec()
608: */
609: PETSC_INTERN PetscErrorCode PetscFEIntegrateHybridResidual_Basic(PetscDS ds, PetscDS dsIn, PetscFormKey key, PetscInt s, PetscInt Ne, PetscFEGeom *fgeom, PetscFEGeom *nbrgeom, const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscReal t, PetscScalar elemVec[])
610: {
611: const PetscInt debug = ds->printIntegrate;
612: const PetscInt field = key.field;
613: PetscFE fe;
614: PetscWeakForm wf;
615: PetscInt n0, n1, i;
616: PetscBdPointFn **f0_func, **f1_func;
617: PetscQuadrature quad;
618: DMPolytopeType ct;
619: PetscTabulation *Tf, *TfIn, *TfAux = NULL;
620: PetscScalar *f0, *f1, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal;
621: const PetscScalar *constants;
622: PetscReal *x;
623: PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL;
624: PetscInt dim, dimAux, numConstants, Nf, NfAux = 0, totDim, totDimIn, totDimAux = 0, cOffset = 0, cOffsetIn = 0, cOffsetAux = 0, fOffset, e, NcI, NcS;
625: PetscBool isCohesiveField, auxOnBd = PETSC_FALSE;
626: const PetscReal *quadPoints, *quadWeights;
627: PetscInt qdim, qNc, Nq, q, dE;
629: PetscFunctionBegin;
630: /* Hybrid discretization is posed directly on faces */
631: PetscCall(PetscDSGetDiscretization(ds, field, (PetscObject *)&fe));
632: PetscCall(PetscFEGetSpatialDimension(fe, &dim));
633: PetscCall(PetscFEGetQuadrature(fe, &quad));
634: PetscCall(PetscDSGetNumFields(ds, &Nf));
635: PetscCall(PetscDSGetTotalDimension(ds, &totDim));
636: PetscCall(PetscDSGetTotalDimension(dsIn, &totDimIn));
637: PetscCall(PetscDSGetComponentOffsetsCohesive(dsIn, 0, &uOff)); // Change 0 to s for one-sided offsets
638: PetscCall(PetscDSGetComponentDerivativeOffsetsCohesive(dsIn, s, &uOff_x));
639: PetscCall(PetscDSGetFieldOffsetCohesive(ds, field, &fOffset));
640: PetscCall(PetscDSGetWeakForm(ds, &wf));
641: PetscCall(PetscWeakFormGetBdResidual(wf, key.label, key.value, key.field, key.part, &n0, &f0_func, &n1, &f1_func));
642: if (!n0 && !n1) PetscFunctionReturn(PETSC_SUCCESS);
643: PetscCall(PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x));
644: PetscCall(PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, NULL, NULL));
645: PetscCall(PetscDSGetWeakFormArrays(ds, &f0, &f1, NULL, NULL, NULL, NULL));
646: /* NOTE This is a bulk tabulation because the DS is a face discretization */
647: PetscCall(PetscDSGetTabulation(ds, &Tf));
648: PetscCall(PetscDSGetFaceTabulation(dsIn, &TfIn));
649: PetscCall(PetscDSSetIntegrationParameters(ds, field, PETSC_DETERMINE));
650: PetscCall(PetscDSGetConstants(ds, &numConstants, &constants));
651: if (dsAux) {
652: PetscCall(PetscDSGetSpatialDimension(dsAux, &dimAux));
653: PetscCall(PetscDSGetNumFields(dsAux, &NfAux));
654: PetscCall(PetscDSGetTotalDimension(dsAux, &totDimAux));
655: PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff));
656: PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x));
657: PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x));
658: auxOnBd = dimAux == dim ? PETSC_TRUE : PETSC_FALSE;
659: if (auxOnBd) PetscCall(PetscDSGetTabulation(dsAux, &TfAux));
660: else PetscCall(PetscDSGetFaceTabulation(dsAux, &TfAux));
661: PetscCheck(Tf[0]->Np == TfAux[0]->Np, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of tabulation points %" PetscInt_FMT " != %" PetscInt_FMT " number of auxiliary tabulation points", Tf[0]->Np, TfAux[0]->Np);
662: }
663: PetscCall(PetscDSGetCohesive(ds, field, &isCohesiveField));
664: NcI = Tf[field]->Nc;
665: NcS = NcI;
666: if (!isCohesiveField && s == 2) {
667: // If we are integrating over a cohesive cell (s = 2) for a non-cohesive fields, we use both sides
668: NcS *= 2;
669: }
670: PetscCall(PetscQuadratureGetData(quad, &qdim, &qNc, &Nq, &quadPoints, &quadWeights));
671: PetscCall(PetscQuadratureGetCellType(quad, &ct));
672: PetscCheck(qNc == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %" PetscInt_FMT " components", qNc);
673: dE = fgeom->dimEmbed;
674: PetscCheck(fgeom->dim == qdim, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "FEGeom dim %" PetscInt_FMT " != %" PetscInt_FMT " quadrature dim", fgeom->dim, qdim);
675: for (e = 0; e < Ne; ++e) {
676: // In order for the face information to be correct, the support of endcap faces _must_ be correctly oriented
677: PetscFEGeom fegeom, fegeomN[2];
678: const PetscInt face[2] = {fgeom->face[e * 2 + 0][0], fgeom->face[e * 2 + 1][2]};
679: const PetscInt ornt[2] = {fgeom->face[e * 2 + 0][1], fgeom->face[e * 2 + 1][3]};
680: const PetscInt cornt[2] = {fgeom->face[e * 2 + 0][3], fgeom->face[e * 2 + 1][1]};
682: fegeom.v = x; /* Workspace */
683: PetscCall(PetscArrayzero(f0, Nq * NcS));
684: PetscCall(PetscArrayzero(f1, Nq * NcS * dE));
685: if (debug > 2) {
686: PetscCall(PetscPrintf(PETSC_COMM_SELF, "Negative %s face: %" PetscInt_FMT " (%" PetscInt_FMT ") (%" PetscInt_FMT ") perm %" PetscInt_FMT "\n", DMPolytopeTypes[ct], face[0], ornt[0], cornt[0], DMPolytopeTypeComposeOrientationInv(ct, cornt[0], ornt[0])));
687: PetscCall(PetscPrintf(PETSC_COMM_SELF, "Positive %s face: %" PetscInt_FMT " (%" PetscInt_FMT ") (%" PetscInt_FMT ") perm %" PetscInt_FMT "\n", DMPolytopeTypes[ct], face[1], ornt[1], cornt[1], DMPolytopeTypeComposeOrientationInv(ct, cornt[1], ornt[1])));
688: }
689: for (q = 0; q < Nq; ++q) {
690: PetscInt qpt[2];
691: PetscReal w;
692: PetscInt c;
694: PetscCall(PetscDSPermuteQuadPoint(ds, DMPolytopeTypeComposeOrientationInv(ct, cornt[0], ornt[0]), field, q, &qpt[0]));
695: PetscCall(PetscDSPermuteQuadPoint(ds, DMPolytopeTypeComposeOrientationInv(ct, cornt[1], ornt[1]), field, q, &qpt[1]));
696: PetscCall(PetscFEGeomGetPoint(fgeom, e * 2, q, &quadPoints[q * fgeom->dim], &fegeom));
697: PetscCall(PetscFEGeomGetPoint(nbrgeom, e * 2, q, NULL, &fegeomN[0]));
698: PetscCall(PetscFEGeomGetPoint(nbrgeom, e * 2 + 1, q, NULL, &fegeomN[1]));
699: w = fegeom.detJ[0] * quadWeights[q];
700: if (debug > 1 && q < fgeom->numPoints) {
701: PetscCall(PetscPrintf(PETSC_COMM_SELF, " detJ: %g\n", (double)fegeom.detJ[0]));
702: #if !defined(PETSC_USE_COMPLEX)
703: PetscCall(DMPrintCellMatrix(e, "invJ", dim, dE, fegeom.invJ));
704: #endif
705: }
706: if (debug) PetscCall(PetscPrintf(PETSC_COMM_SELF, " quad point %" PetscInt_FMT " weight %g detJ %g\n", q, (double)quadWeights[q], (double)fegeom.detJ[0]));
707: /* TODO Is this cell or face quadrature, meaning should we use 'q' or 'face*Nq+q' */
708: PetscCall(PetscFEEvaluateFieldJets_Hybrid_Internal(dsIn, Nf, 0, q, Tf, face, qpt, TfIn, &fegeom, fegeomN, &coefficients[cOffsetIn], PetscSafePointerPlusOffset(coefficients_t, cOffsetIn), u, u_x, u_t));
709: if (dsAux) PetscCall(PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, auxOnBd ? 0 : face[s], auxOnBd ? q : qpt[s], TfAux, &fegeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL));
710: for (i = 0; i < n0; ++i) f0_func[i](dE, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, fegeom.v, fegeom.n, numConstants, constants, &f0[q * NcS]);
711: for (c = 0; c < NcS; ++c) f0[q * NcS + c] *= w;
712: for (i = 0; i < n1; ++i) f1_func[i](dE, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, fegeom.v, fegeom.n, numConstants, constants, &f1[q * NcS * dE]);
713: for (c = 0; c < NcS; ++c)
714: for (PetscInt d = 0; d < dE; ++d) f1[(q * NcS + c) * dE + d] *= w;
715: if (debug) {
716: PetscCall(PetscPrintf(PETSC_COMM_SELF, " elem %" PetscInt_FMT " quad point %" PetscInt_FMT " field %" PetscInt_FMT " side %" PetscInt_FMT "\n", e, q, field, s));
717: for (PetscInt f = 0; f < Nf; ++f) {
718: PetscCall(PetscPrintf(PETSC_COMM_SELF, " Field %" PetscInt_FMT ":", f));
719: for (PetscInt c = uOff[f]; c < uOff[f + 1]; ++c) PetscCall(PetscPrintf(PETSC_COMM_SELF, " %g", (double)PetscRealPart(u[c])));
720: PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n"));
721: }
722: for (c = 0; c < NcS; ++c) {
723: if (n0) PetscCall(PetscPrintf(PETSC_COMM_SELF, " f0[%" PetscInt_FMT "] %g\n", c, (double)PetscRealPart(f0[q * NcS + c])));
724: if (n1) {
725: for (PetscInt d = 0; d < dE; ++d) PetscCall(PetscPrintf(PETSC_COMM_SELF, " f1[%" PetscInt_FMT ",%" PetscInt_FMT "] %g", c, d, (double)PetscRealPart(f1[(q * NcS + c) * dE + d])));
726: PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n"));
727: }
728: }
729: }
730: }
731: if (isCohesiveField) {
732: PetscCall(PetscFEUpdateElementVec_Internal(fe, Tf[field], 0, basisReal, basisDerReal, e, fgeom, f0, f1, &elemVec[cOffset + fOffset]));
733: } else {
734: PetscCall(PetscFEUpdateElementVec_Hybrid_Internal(fe, Tf[field], 0, s, basisReal, basisDerReal, fgeom, f0, f1, &elemVec[cOffset + fOffset]));
735: }
736: cOffset += totDim;
737: cOffsetIn += totDimIn;
738: cOffsetAux += totDimAux;
739: }
740: PetscFunctionReturn(PETSC_SUCCESS);
741: }
743: PetscErrorCode PetscFEIntegrateJacobian_Basic(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[])
744: {
745: const PetscInt debug = rds->printIntegrate;
746: PetscFE feI, feJ;
747: PetscWeakForm wf;
748: PetscPointJacFn **g0_func, **g1_func, **g2_func, **g3_func;
749: PetscInt n0, n1, n2, n3;
750: PetscInt cOffset = 0; /* Offset into coefficients[] for element e */
751: PetscInt cOffsetAux = 0; /* Offset into coefficientsAux[] for element e */
752: PetscInt eOffset = 0; /* Offset into elemMat[] for element e */
753: PetscInt offsetI = 0; /* Offset into an element vector for fieldI */
754: PetscInt offsetJ = 0; /* Offset into an element vector for fieldJ */
755: PetscQuadrature quad;
756: PetscTabulation *rT, *cT, *TAux = NULL;
757: PetscScalar *g0 = NULL, *g1 = NULL, *g2 = NULL, *g3 = NULL, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal, *testReal, *testDerReal;
758: const PetscScalar *constants;
759: PetscReal *x, cellScale;
760: PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL;
761: PetscInt NcI = 0, NcJ = 0;
762: PetscInt dim, numConstants, Nf, fieldI, fieldJ, NfAux = 0, rtotDim, ctotDim, totDimAux = 0;
763: PetscInt dE, Np;
764: PetscBool isAffine;
765: const PetscReal *quadPoints, *quadWeights;
766: PetscInt qNc, Nq;
768: PetscFunctionBegin;
769: PetscCall(PetscDSGetNumFields(rds, &Nf));
770: fieldI = key.field / Nf;
771: fieldJ = key.field % Nf;
772: PetscCall(PetscDSGetDiscretization(rds, fieldI, (PetscObject *)&feI));
773: PetscCall(PetscDSGetDiscretization(cds, fieldJ, (PetscObject *)&feJ));
774: PetscCall(PetscFEGetSpatialDimension(feI, &dim));
775: cellScale = (PetscReal)PetscPowInt(2, dim);
776: PetscCall(PetscFEGetQuadrature(feI, &quad));
777: PetscCall(PetscDSGetTotalDimension(rds, &rtotDim));
778: PetscCall(PetscDSGetTotalDimension(cds, &ctotDim));
779: PetscCall(PetscDSGetComponentOffsets(rds, &uOff));
780: PetscCall(PetscDSGetComponentDerivativeOffsets(rds, &uOff_x));
781: PetscCall(PetscDSGetWeakForm(rds, &wf));
782: switch (jtype) {
783: case PETSCFE_JACOBIAN_DYN:
784: PetscCall(PetscWeakFormGetDynamicJacobian(wf, key.label, key.value, fieldI, fieldJ, key.part, &n0, &g0_func, &n1, &g1_func, &n2, &g2_func, &n3, &g3_func));
785: break;
786: case PETSCFE_JACOBIAN_PRE:
787: PetscCall(PetscWeakFormGetJacobianPreconditioner(wf, key.label, key.value, fieldI, fieldJ, key.part, &n0, &g0_func, &n1, &g1_func, &n2, &g2_func, &n3, &g3_func));
788: break;
789: case PETSCFE_JACOBIAN:
790: PetscCall(PetscWeakFormGetJacobian(wf, key.label, key.value, fieldI, fieldJ, key.part, &n0, &g0_func, &n1, &g1_func, &n2, &g2_func, &n3, &g3_func));
791: break;
792: }
793: if (!n0 && !n1 && !n2 && !n3) PetscFunctionReturn(PETSC_SUCCESS);
794: PetscCall(PetscDSGetEvaluationArrays(rds, &u, coefficients_t ? &u_t : NULL, &u_x));
795: PetscCall(PetscDSGetWorkspace(rds, &x, &basisReal, &basisDerReal, &testReal, &testDerReal));
796: PetscCall(PetscDSGetWeakFormArrays(rds, NULL, NULL, n0 ? &g0 : NULL, n1 ? &g1 : NULL, n2 ? &g2 : NULL, n3 ? &g3 : NULL));
798: PetscCall(PetscDSGetTabulation(rds, &rT));
799: PetscCall(PetscDSGetTabulation(cds, &cT));
800: PetscCheck(rT[0]->Np == cT[0]->Np, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of row tabulation points %" PetscInt_FMT " != %" PetscInt_FMT " number of col tabulation points", rT[0]->Np, cT[0]->Np);
801: PetscCall(PetscDSGetFieldOffset(rds, fieldI, &offsetI));
802: PetscCall(PetscDSGetFieldOffset(cds, fieldJ, &offsetJ));
803: PetscCall(PetscDSSetIntegrationParameters(rds, fieldI, fieldJ));
804: PetscCall(PetscDSSetIntegrationParameters(cds, fieldI, fieldJ));
805: PetscCall(PetscDSGetConstants(rds, &numConstants, &constants));
806: if (dsAux) {
807: PetscCall(PetscDSGetNumFields(dsAux, &NfAux));
808: PetscCall(PetscDSGetTotalDimension(dsAux, &totDimAux));
809: PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff));
810: PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x));
811: PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x));
812: PetscCall(PetscDSGetTabulation(dsAux, &TAux));
813: PetscCheck(rT[0]->Np == TAux[0]->Np, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of tabulation points %" PetscInt_FMT " != %" PetscInt_FMT " number of auxiliary tabulation points", rT[0]->Np, TAux[0]->Np);
814: }
815: NcI = rT[fieldI]->Nc;
816: NcJ = cT[fieldJ]->Nc;
817: Np = cgeom->numPoints;
818: dE = cgeom->dimEmbed;
819: isAffine = cgeom->isAffine;
820: PetscCall(PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights));
821: PetscCheck(qNc == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %" PetscInt_FMT " components", qNc);
823: for (PetscInt e = 0; e < Ne; ++e) {
824: PetscFEGeom fegeom;
826: fegeom.dim = cgeom->dim;
827: fegeom.dimEmbed = cgeom->dimEmbed;
828: fegeom.xi = NULL;
829: if (isAffine) {
830: fegeom.v = x;
831: fegeom.xi = cgeom->xi;
832: fegeom.J = &cgeom->J[e * Np * dE * dE];
833: fegeom.invJ = &cgeom->invJ[e * Np * dE * dE];
834: fegeom.detJ = &cgeom->detJ[e * Np];
835: } else fegeom.xi = NULL;
836: for (PetscInt q = 0; q < Nq; ++q) {
837: PetscReal w;
839: if (isAffine) {
840: CoordinatesRefToReal(dE, dim, fegeom.xi, &cgeom->v[e * Np * dE], fegeom.J, &quadPoints[q * dim], x);
841: } else {
842: fegeom.v = &cgeom->v[(e * Np + q) * dE];
843: fegeom.J = &cgeom->J[(e * Np + q) * dE * dE];
844: fegeom.invJ = &cgeom->invJ[(e * Np + q) * dE * dE];
845: fegeom.detJ = &cgeom->detJ[e * Np + q];
846: }
847: PetscCall(PetscDSSetCellParameters(rds, fegeom.detJ[0] * cellScale));
848: if (debug) PetscCall(PetscPrintf(PETSC_COMM_SELF, " quad point %" PetscInt_FMT " weight %g detJ %g\n", q, (double)quadWeights[q], (double)fegeom.detJ[0]));
849: w = fegeom.detJ[0] * quadWeights[q];
850: if (coefficients) PetscCall(PetscFEEvaluateFieldJets_Internal(rds, Nf, 0, q, rT, &fegeom, &coefficients[cOffset], PetscSafePointerPlusOffset(coefficients_t, cOffset), u, u_x, u_t));
851: if (dsAux) PetscCall(PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, 0, q, TAux, &fegeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL));
852: if (n0) {
853: PetscCall(PetscArrayzero(g0, NcI * NcJ));
854: for (PetscInt i = 0; i < n0; ++i) g0_func[i](dE, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, u_tshift, fegeom.v, numConstants, constants, g0);
855: for (PetscInt c = 0; c < NcI * NcJ; ++c) g0[c] *= w;
856: }
857: if (n1) {
858: PetscCall(PetscArrayzero(g1, NcI * NcJ * dE));
859: for (PetscInt i = 0; i < n1; ++i) g1_func[i](dE, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, u_tshift, fegeom.v, numConstants, constants, g1);
860: for (PetscInt c = 0; c < NcI * NcJ * dE; ++c) g1[c] *= w;
861: }
862: if (n2) {
863: PetscCall(PetscArrayzero(g2, NcI * NcJ * dE));
864: for (PetscInt i = 0; i < n2; ++i) g2_func[i](dE, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, u_tshift, fegeom.v, numConstants, constants, g2);
865: for (PetscInt c = 0; c < NcI * NcJ * dE; ++c) g2[c] *= w;
866: }
867: if (n3) {
868: PetscCall(PetscArrayzero(g3, NcI * NcJ * dE * dE));
869: for (PetscInt i = 0; i < n3; ++i) g3_func[i](dE, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, u_tshift, fegeom.v, numConstants, constants, g3);
870: for (PetscInt c = 0; c < NcI * NcJ * dE * dE; ++c) g3[c] *= w;
871: }
873: PetscCall(PetscFEUpdateElementMat_Internal(feI, feJ, 0, q, rT[fieldI], basisReal, basisDerReal, cT[fieldJ], testReal, testDerReal, &fegeom, g0, g1, g2, g3, ctotDim, offsetI, offsetJ, elemMat + eOffset));
874: }
875: if (debug > 1) {
876: PetscCall(PetscPrintf(PETSC_COMM_SELF, "Element matrix for fields %" PetscInt_FMT " and %" PetscInt_FMT "\n", fieldI, fieldJ));
877: for (PetscInt f = 0; f < rT[fieldI]->Nb; ++f) {
878: const PetscInt i = offsetI + f;
879: for (PetscInt g = 0; g < cT[fieldJ]->Nb; ++g) {
880: const PetscInt j = offsetJ + g;
881: PetscCall(PetscPrintf(PETSC_COMM_SELF, " elemMat[%" PetscInt_FMT ", %" PetscInt_FMT "]: %g\n", f, g, (double)PetscRealPart(elemMat[eOffset + i * ctotDim + j])));
882: }
883: PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n"));
884: }
885: }
886: cOffset += rtotDim;
887: cOffsetAux += totDimAux;
888: eOffset += rtotDim * ctotDim;
889: }
890: PetscFunctionReturn(PETSC_SUCCESS);
891: }
893: PETSC_INTERN PetscErrorCode PetscFEIntegrateBdJacobian_Basic(PetscDS ds, PetscWeakForm wf, PetscFEJacobianType jtype, PetscFormKey key, PetscInt Ne, PetscFEGeom *fgeom, const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscReal t, PetscReal u_tshift, PetscScalar elemMat[])
894: {
895: const PetscInt debug = ds->printIntegrate;
896: PetscFE feI, feJ;
897: PetscBdPointJacFn **g0_func, **g1_func, **g2_func, **g3_func;
898: PetscInt n0, n1, n2, n3, i;
899: PetscInt cOffset = 0; /* Offset into coefficients[] for element e */
900: PetscInt cOffsetAux = 0; /* Offset into coefficientsAux[] for element e */
901: PetscInt eOffset = 0; /* Offset into elemMat[] for element e */
902: PetscInt offsetI = 0; /* Offset into an element vector for fieldI */
903: PetscInt offsetJ = 0; /* Offset into an element vector for fieldJ */
904: PetscQuadrature quad;
905: PetscTabulation *T, *TAux = NULL;
906: PetscScalar *g0, *g1, *g2, *g3, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal, *testReal, *testDerReal;
907: const PetscScalar *constants;
908: PetscReal *x, cellScale;
909: PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL;
910: PetscInt NcI = 0, NcJ = 0;
911: PetscInt dim, numConstants, Nf, fieldI, fieldJ, NfAux = 0, totDim, totDimAux = 0, e;
912: PetscBool isAffine;
913: const PetscReal *quadPoints, *quadWeights;
914: PetscInt qNc, Nq, q, Np, dE;
916: PetscFunctionBegin;
917: PetscCall(PetscDSGetNumFields(ds, &Nf));
918: fieldI = key.field / Nf;
919: fieldJ = key.field % Nf;
920: PetscCall(PetscDSGetDiscretization(ds, fieldI, (PetscObject *)&feI));
921: PetscCall(PetscDSGetDiscretization(ds, fieldJ, (PetscObject *)&feJ));
922: PetscCall(PetscFEGetSpatialDimension(feI, &dim));
923: cellScale = (PetscReal)PetscPowInt(2, dim);
924: PetscCall(PetscFEGetFaceQuadrature(feI, &quad));
925: PetscCall(PetscDSGetTotalDimension(ds, &totDim));
926: PetscCall(PetscDSGetComponentOffsets(ds, &uOff));
927: PetscCall(PetscDSGetComponentDerivativeOffsets(ds, &uOff_x));
928: PetscCall(PetscDSGetFieldOffset(ds, fieldI, &offsetI));
929: PetscCall(PetscDSGetFieldOffset(ds, fieldJ, &offsetJ));
930: switch (jtype) {
931: case PETSCFE_JACOBIAN_PRE:
932: PetscCall(PetscWeakFormGetBdJacobianPreconditioner(wf, key.label, key.value, fieldI, fieldJ, key.part, &n0, &g0_func, &n1, &g1_func, &n2, &g2_func, &n3, &g3_func));
933: break;
934: case PETSCFE_JACOBIAN:
935: PetscCall(PetscWeakFormGetBdJacobian(wf, key.label, key.value, fieldI, fieldJ, key.part, &n0, &g0_func, &n1, &g1_func, &n2, &g2_func, &n3, &g3_func));
936: break;
937: case PETSCFE_JACOBIAN_DYN:
938: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "PETSCFE_JACOBIAN_DYN is not supported for PetscFEIntegrateBdJacobian()");
939: }
940: if (!n0 && !n1 && !n2 && !n3) PetscFunctionReturn(PETSC_SUCCESS);
941: PetscCall(PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x));
942: PetscCall(PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, &testReal, &testDerReal));
943: PetscCall(PetscDSGetWeakFormArrays(ds, NULL, NULL, &g0, &g1, &g2, &g3));
944: PetscCall(PetscDSGetFaceTabulation(ds, &T));
945: PetscCall(PetscDSSetIntegrationParameters(ds, fieldI, fieldJ));
946: PetscCall(PetscDSGetConstants(ds, &numConstants, &constants));
947: if (dsAux) {
948: PetscCall(PetscDSGetNumFields(dsAux, &NfAux));
949: PetscCall(PetscDSGetTotalDimension(dsAux, &totDimAux));
950: PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff));
951: PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x));
952: PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x));
953: PetscCall(PetscDSGetFaceTabulation(dsAux, &TAux));
954: }
955: NcI = T[fieldI]->Nc, NcJ = T[fieldJ]->Nc;
956: Np = fgeom->numPoints;
957: dE = fgeom->dimEmbed;
958: isAffine = fgeom->isAffine;
959: /* Initialize here in case the function is not defined */
960: PetscCall(PetscArrayzero(g0, NcI * NcJ));
961: PetscCall(PetscArrayzero(g1, NcI * NcJ * dE));
962: PetscCall(PetscArrayzero(g2, NcI * NcJ * dE));
963: PetscCall(PetscArrayzero(g3, NcI * NcJ * dE * dE));
964: PetscCall(PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights));
965: PetscCheck(qNc == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %" PetscInt_FMT " components", qNc);
966: for (e = 0; e < Ne; ++e) {
967: PetscFEGeom fegeom, cgeom;
968: const PetscInt face = fgeom->face[e][0];
969: const PetscInt ornt = fgeom->face[e][1];
970: fegeom.n = NULL;
971: fegeom.v = NULL;
972: fegeom.xi = NULL;
973: fegeom.J = NULL;
974: fegeom.detJ = NULL;
975: fegeom.dim = fgeom->dim;
976: fegeom.dimEmbed = fgeom->dimEmbed;
977: cgeom.dim = fgeom->dim;
978: cgeom.dimEmbed = fgeom->dimEmbed;
979: if (isAffine) {
980: fegeom.v = x;
981: fegeom.xi = fgeom->xi;
982: fegeom.J = &fgeom->J[e * Np * dE * dE];
983: fegeom.invJ = &fgeom->invJ[e * Np * dE * dE];
984: fegeom.detJ = &fgeom->detJ[e * Np];
985: fegeom.n = &fgeom->n[e * Np * dE];
987: cgeom.J = &fgeom->suppJ[0][e * Np * dE * dE];
988: cgeom.invJ = &fgeom->suppInvJ[0][e * Np * dE * dE];
989: cgeom.detJ = &fgeom->suppDetJ[0][e * Np];
990: } else fegeom.xi = NULL;
991: for (q = 0; q < Nq; ++q) {
992: PetscReal w;
993: const PetscInt qp = ornt < 0 ? (Nq - 1 - q) : q; /* Map physical quadrature index to tabulation index accounting for face orientation */
995: if (debug) PetscCall(PetscPrintf(PETSC_COMM_SELF, " quad point %" PetscInt_FMT "\n", q));
996: if (isAffine) {
997: CoordinatesRefToReal(dE, dim - 1, fegeom.xi, &fgeom->v[e * Np * dE], fegeom.J, &quadPoints[q * (dim - 1)], x);
998: } else {
999: fegeom.v = &fgeom->v[(e * Np + q) * dE];
1000: fegeom.J = &fgeom->J[(e * Np + q) * dE * dE];
1001: fegeom.invJ = &fgeom->invJ[(e * Np + q) * dE * dE];
1002: fegeom.detJ = &fgeom->detJ[e * Np + q];
1003: fegeom.n = &fgeom->n[(e * Np + q) * dE];
1005: cgeom.J = &fgeom->suppJ[0][(e * Np + q) * dE * dE];
1006: cgeom.invJ = &fgeom->suppInvJ[0][(e * Np + q) * dE * dE];
1007: cgeom.detJ = &fgeom->suppDetJ[0][e * Np + q];
1008: }
1009: PetscCall(PetscDSSetCellParameters(ds, fegeom.detJ[0] * cellScale));
1010: w = fegeom.detJ[0] * quadWeights[q];
1011: if (coefficients) PetscCall(PetscFEEvaluateFieldJets_Internal(ds, Nf, face, qp, T, &cgeom, &coefficients[cOffset], PetscSafePointerPlusOffset(coefficients_t, cOffset), u, u_x, u_t));
1012: if (dsAux) PetscCall(PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, face, qp, TAux, &cgeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL));
1013: if (n0) {
1014: PetscCall(PetscArrayzero(g0, NcI * NcJ));
1015: for (i = 0; i < n0; ++i) g0_func[i](dE, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, u_tshift, fegeom.v, fegeom.n, numConstants, constants, g0);
1016: for (PetscInt c = 0; c < NcI * NcJ; ++c) g0[c] *= w;
1017: }
1018: if (n1) {
1019: PetscCall(PetscArrayzero(g1, NcI * NcJ * dE));
1020: for (i = 0; i < n1; ++i) g1_func[i](dE, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, u_tshift, fegeom.v, fegeom.n, numConstants, constants, g1);
1021: for (PetscInt c = 0; c < NcI * NcJ * dim; ++c) g1[c] *= w;
1022: }
1023: if (n2) {
1024: PetscCall(PetscArrayzero(g2, NcI * NcJ * dE));
1025: for (i = 0; i < n2; ++i) g2_func[i](dE, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, u_tshift, fegeom.v, fegeom.n, numConstants, constants, g2);
1026: for (PetscInt c = 0; c < NcI * NcJ * dim; ++c) g2[c] *= w;
1027: }
1028: if (n3) {
1029: PetscCall(PetscArrayzero(g3, NcI * NcJ * dE * dE));
1030: for (i = 0; i < n3; ++i) g3_func[i](dE, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, u_tshift, fegeom.v, fegeom.n, numConstants, constants, g3);
1031: for (PetscInt c = 0; c < NcI * NcJ * dim * dim; ++c) g3[c] *= w;
1032: }
1034: PetscCall(PetscFEUpdateElementMat_Internal(feI, feJ, face, qp, T[fieldI], basisReal, basisDerReal, T[fieldJ], testReal, testDerReal, &cgeom, g0, g1, g2, g3, totDim, offsetI, offsetJ, elemMat + eOffset));
1035: }
1036: if (debug > 1) {
1037: PetscInt fc, f, gc, g;
1039: PetscCall(PetscPrintf(PETSC_COMM_SELF, "Element matrix for fields %" PetscInt_FMT " and %" PetscInt_FMT "\n", fieldI, fieldJ));
1040: for (fc = 0; fc < T[fieldI]->Nc; ++fc) {
1041: for (f = 0; f < T[fieldI]->Nb; ++f) {
1042: const PetscInt i = offsetI + f * T[fieldI]->Nc + fc;
1043: for (gc = 0; gc < T[fieldJ]->Nc; ++gc) {
1044: for (g = 0; g < T[fieldJ]->Nb; ++g) {
1045: const PetscInt j = offsetJ + g * T[fieldJ]->Nc + gc;
1046: PetscCall(PetscPrintf(PETSC_COMM_SELF, " elemMat[%" PetscInt_FMT ",%" PetscInt_FMT ",%" PetscInt_FMT ",%" PetscInt_FMT "]: %g\n", f, fc, g, gc, (double)PetscRealPart(elemMat[eOffset + i * totDim + j])));
1047: }
1048: }
1049: PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n"));
1050: }
1051: }
1052: }
1053: cOffset += totDim;
1054: cOffsetAux += totDimAux;
1055: eOffset += PetscSqr(totDim);
1056: }
1057: PetscFunctionReturn(PETSC_SUCCESS);
1058: }
1060: PETSC_INTERN PetscErrorCode PetscFEIntegrateHybridJacobian_Basic(PetscDS ds, PetscDS dsIn, PetscFEJacobianType jtype, PetscFormKey key, PetscInt s, PetscInt Ne, PetscFEGeom *fgeom, PetscFEGeom *nbrgeom, const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscReal t, PetscReal u_tshift, PetscScalar elemMat[])
1061: {
1062: const PetscInt debug = ds->printIntegrate;
1063: PetscFE feI, feJ;
1064: PetscWeakForm wf;
1065: PetscBdPointJacFn **g0_func, **g1_func, **g2_func, **g3_func;
1066: PetscInt n0, n1, n2, n3, i;
1067: PetscInt cOffset = 0; /* Offset into coefficients[] for element e */
1068: PetscInt cOffsetAux = 0; /* Offset into coefficientsAux[] for element e */
1069: PetscInt eOffset = 0; /* Offset into elemMat[] for element e */
1070: PetscInt offsetI = 0; /* Offset into an element vector for fieldI */
1071: PetscInt offsetJ = 0; /* Offset into an element vector for fieldJ */
1072: PetscQuadrature quad;
1073: DMPolytopeType ct;
1074: PetscTabulation *T, *TfIn, *TAux = NULL;
1075: PetscScalar *g0, *g1, *g2, *g3, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal, *testReal, *testDerReal;
1076: const PetscScalar *constants;
1077: PetscReal *x;
1078: PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL;
1079: PetscInt NcI = 0, NcJ = 0, NcS, NcT;
1080: PetscInt dim, dimAux, numConstants, Nf, fieldI, fieldJ, NfAux = 0, totDim, totDimAux = 0, e;
1081: PetscBool isCohesiveFieldI, isCohesiveFieldJ, auxOnBd = PETSC_FALSE;
1082: const PetscReal *quadPoints, *quadWeights;
1083: PetscInt qNc, Nq, q, dE;
1085: PetscFunctionBegin;
1086: PetscCall(PetscDSGetNumFields(ds, &Nf));
1087: fieldI = key.field / Nf;
1088: fieldJ = key.field % Nf;
1089: /* Hybrid discretization is posed directly on faces */
1090: PetscCall(PetscDSGetDiscretization(ds, fieldI, (PetscObject *)&feI));
1091: PetscCall(PetscDSGetDiscretization(ds, fieldJ, (PetscObject *)&feJ));
1092: PetscCall(PetscFEGetSpatialDimension(feI, &dim));
1093: PetscCall(PetscFEGetQuadrature(feI, &quad));
1094: PetscCall(PetscDSGetTotalDimension(ds, &totDim));
1095: PetscCall(PetscDSGetComponentOffsetsCohesive(ds, 0, &uOff)); // Change 0 to s for one-sided offsets
1096: PetscCall(PetscDSGetComponentDerivativeOffsetsCohesive(ds, s, &uOff_x));
1097: PetscCall(PetscDSGetWeakForm(ds, &wf));
1098: switch (jtype) {
1099: case PETSCFE_JACOBIAN_PRE:
1100: PetscCall(PetscWeakFormGetBdJacobianPreconditioner(wf, key.label, key.value, fieldI, fieldJ, key.part, &n0, &g0_func, &n1, &g1_func, &n2, &g2_func, &n3, &g3_func));
1101: break;
1102: case PETSCFE_JACOBIAN:
1103: PetscCall(PetscWeakFormGetBdJacobian(wf, key.label, key.value, fieldI, fieldJ, key.part, &n0, &g0_func, &n1, &g1_func, &n2, &g2_func, &n3, &g3_func));
1104: break;
1105: case PETSCFE_JACOBIAN_DYN:
1106: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "No boundary hybrid Jacobians :)");
1107: }
1108: if (!n0 && !n1 && !n2 && !n3) PetscFunctionReturn(PETSC_SUCCESS);
1109: PetscCall(PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x));
1110: PetscCall(PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, &testReal, &testDerReal));
1111: PetscCall(PetscDSGetWeakFormArrays(ds, NULL, NULL, &g0, &g1, &g2, &g3));
1112: PetscCall(PetscDSGetTabulation(ds, &T));
1113: PetscCall(PetscDSGetFaceTabulation(dsIn, &TfIn));
1114: PetscCall(PetscDSGetFieldOffsetCohesive(ds, fieldI, &offsetI));
1115: PetscCall(PetscDSGetFieldOffsetCohesive(ds, fieldJ, &offsetJ));
1116: PetscCall(PetscDSSetIntegrationParameters(ds, fieldI, fieldJ));
1117: PetscCall(PetscDSGetConstants(ds, &numConstants, &constants));
1118: if (dsAux) {
1119: PetscCall(PetscDSGetSpatialDimension(dsAux, &dimAux));
1120: PetscCall(PetscDSGetNumFields(dsAux, &NfAux));
1121: PetscCall(PetscDSGetTotalDimension(dsAux, &totDimAux));
1122: PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff));
1123: PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x));
1124: PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x));
1125: auxOnBd = dimAux == dim ? PETSC_TRUE : PETSC_FALSE;
1126: if (auxOnBd) PetscCall(PetscDSGetTabulation(dsAux, &TAux));
1127: else PetscCall(PetscDSGetFaceTabulation(dsAux, &TAux));
1128: PetscCheck(T[0]->Np == TAux[0]->Np, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of tabulation points %" PetscInt_FMT " != %" PetscInt_FMT " number of auxiliary tabulation points", T[0]->Np, TAux[0]->Np);
1129: }
1130: PetscCall(PetscDSGetCohesive(ds, fieldI, &isCohesiveFieldI));
1131: PetscCall(PetscDSGetCohesive(ds, fieldJ, &isCohesiveFieldJ));
1132: dE = fgeom->dimEmbed;
1133: NcI = T[fieldI]->Nc;
1134: NcJ = T[fieldJ]->Nc;
1135: NcS = isCohesiveFieldI ? NcI : 2 * NcI;
1136: NcT = isCohesiveFieldJ ? NcJ : 2 * NcJ;
1137: if (!isCohesiveFieldI && s == 2) {
1138: // If we are integrating over a cohesive cell (s = 2) for a non-cohesive fields, we use both sides
1139: NcS *= 2;
1140: }
1141: if (!isCohesiveFieldJ && s == 2) {
1142: // If we are integrating over a cohesive cell (s = 2) for a non-cohesive fields, we use both sides
1143: NcT *= 2;
1144: }
1145: // The derivatives are constrained to be along the cell, so there are dim, not dE, components, even though
1146: // the coordinates are in dE dimensions
1147: PetscCall(PetscArrayzero(g0, NcS * NcT));
1148: PetscCall(PetscArrayzero(g1, NcS * NcT * dim));
1149: PetscCall(PetscArrayzero(g2, NcS * NcT * dim));
1150: PetscCall(PetscArrayzero(g3, NcS * NcT * dim * dim));
1151: PetscCall(PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights));
1152: PetscCall(PetscQuadratureGetCellType(quad, &ct));
1153: PetscCheck(qNc == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %" PetscInt_FMT " components", qNc);
1154: for (e = 0; e < Ne; ++e) {
1155: PetscFEGeom fegeom, fegeomN[2];
1156: const PetscInt face[2] = {fgeom->face[e * 2 + 0][0], fgeom->face[e * 2 + 1][2]};
1157: const PetscInt ornt[2] = {fgeom->face[e * 2 + 0][1], fgeom->face[e * 2 + 1][3]};
1158: const PetscInt cornt[2] = {fgeom->face[e * 2 + 0][3], fgeom->face[e * 2 + 1][1]};
1160: fegeom.v = x; /* Workspace */
1161: for (q = 0; q < Nq; ++q) {
1162: PetscInt qpt[2];
1163: PetscReal w;
1165: PetscCall(PetscDSPermuteQuadPoint(ds, DMPolytopeTypeComposeOrientationInv(ct, cornt[0], ornt[0]), fieldI, q, &qpt[0]));
1166: PetscCall(PetscDSPermuteQuadPoint(ds, DMPolytopeTypeComposeOrientationInv(ct, ornt[1], cornt[1]), fieldI, q, &qpt[1]));
1167: PetscCall(PetscFEGeomGetPoint(fgeom, e * 2, q, &quadPoints[q * fgeom->dim], &fegeom));
1168: PetscCall(PetscFEGeomGetPoint(nbrgeom, e * 2, q, NULL, &fegeomN[0]));
1169: PetscCall(PetscFEGeomGetPoint(nbrgeom, e * 2 + 1, q, NULL, &fegeomN[1]));
1170: w = fegeom.detJ[0] * quadWeights[q];
1171: if (debug > 1 && q < fgeom->numPoints) {
1172: PetscCall(PetscPrintf(PETSC_COMM_SELF, " detJ: %g\n", (double)fegeom.detJ[0]));
1173: #if !defined(PETSC_USE_COMPLEX)
1174: PetscCall(DMPrintCellMatrix(e, "invJ", dim, dim, fegeom.invJ));
1175: #endif
1176: }
1177: if (debug) PetscCall(PetscPrintf(PETSC_COMM_SELF, " quad point %" PetscInt_FMT "\n", q));
1178: if (coefficients) PetscCall(PetscFEEvaluateFieldJets_Hybrid_Internal(dsIn, Nf, 0, q, T, face, qpt, TfIn, &fegeom, fegeomN, &coefficients[cOffset], PetscSafePointerPlusOffset(coefficients_t, cOffset), u, u_x, u_t));
1179: if (dsAux) PetscCall(PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, auxOnBd ? 0 : face[s], auxOnBd ? q : qpt[s], TAux, &fegeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL));
1180: if (n0) {
1181: PetscCall(PetscArrayzero(g0, NcS * NcT));
1182: for (i = 0; i < n0; ++i) g0_func[i](dE, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, u_tshift, fegeom.v, fegeom.n, numConstants, constants, g0);
1183: for (PetscInt c = 0; c < NcS * NcT; ++c) g0[c] *= w;
1184: }
1185: if (n1) {
1186: PetscCall(PetscArrayzero(g1, NcS * NcT * dim));
1187: for (i = 0; i < n1; ++i) g1_func[i](dE, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, u_tshift, fegeom.v, fegeom.n, numConstants, constants, g1);
1188: for (PetscInt c = 0; c < NcS * NcT * dim; ++c) g1[c] *= w;
1189: }
1190: if (n2) {
1191: PetscCall(PetscArrayzero(g2, NcS * NcT * dim));
1192: for (i = 0; i < n2; ++i) g2_func[i](dE, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, u_tshift, fegeom.v, fegeom.n, numConstants, constants, g2);
1193: for (PetscInt c = 0; c < NcS * NcT * dim; ++c) g2[c] *= w;
1194: }
1195: if (n3) {
1196: PetscCall(PetscArrayzero(g3, NcS * NcT * dim * dim));
1197: for (i = 0; i < n3; ++i) g3_func[i](dE, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, u_tshift, fegeom.v, fegeom.n, numConstants, constants, g3);
1198: for (PetscInt c = 0; c < NcS * NcT * dim * dim; ++c) g3[c] *= w;
1199: }
1201: if (isCohesiveFieldI) {
1202: if (isCohesiveFieldJ) {
1203: //PetscCall(PetscFEUpdateElementMat_Internal(feI, feJ, 0, q, T[fieldI], basisReal, basisDerReal, T[fieldJ], testReal, testDerReal, &fegeom, g0, g1, g2, g3, totDim, offsetI, offsetJ, elemMat + eOffset));
1204: PetscCall(PetscFEUpdateElementMat_Hybrid_Internal(feI, isCohesiveFieldI, feJ, isCohesiveFieldJ, 0, 0, 0, q, T[fieldI], basisReal, basisDerReal, T[fieldJ], testReal, testDerReal, &fegeom, g0, g1, g2, g3, eOffset, totDim, offsetI, offsetJ, elemMat));
1205: } else {
1206: PetscCall(PetscFEUpdateElementMat_Hybrid_Internal(feI, isCohesiveFieldI, feJ, isCohesiveFieldJ, 0, 0, 0, q, T[fieldI], basisReal, basisDerReal, T[fieldJ], testReal, testDerReal, &fegeom, g0, g1, g2, g3, eOffset, totDim, offsetI, offsetJ, elemMat));
1207: PetscCall(PetscFEUpdateElementMat_Hybrid_Internal(feI, isCohesiveFieldI, feJ, isCohesiveFieldJ, 0, 1, 1, q, T[fieldI], basisReal, basisDerReal, T[fieldJ], testReal, testDerReal, &fegeom, &g0[NcI * NcJ], &g1[NcI * NcJ * dim], &g2[NcI * NcJ * dim], &g3[NcI * NcJ * dim * dim], eOffset, totDim, offsetI, offsetJ, elemMat));
1208: }
1209: } else {
1210: if (s == 2) {
1211: if (isCohesiveFieldJ) {
1212: PetscCall(PetscFEUpdateElementMat_Hybrid_Internal(feI, isCohesiveFieldI, feJ, isCohesiveFieldJ, 0, 0, 0, q, T[fieldI], basisReal, basisDerReal, T[fieldJ], testReal, testDerReal, &fegeom, g0, g1, g2, g3, eOffset, totDim, offsetI, offsetJ, elemMat));
1213: PetscCall(PetscFEUpdateElementMat_Hybrid_Internal(feI, isCohesiveFieldI, feJ, isCohesiveFieldJ, 0, 1, 1, q, T[fieldI], basisReal, basisDerReal, T[fieldJ], testReal, testDerReal, &fegeom, &g0[NcI * NcJ], &g1[NcI * NcJ * dim], &g2[NcI * NcJ * dim], &g3[NcI * NcJ * dim * dim], eOffset, totDim, offsetI, offsetJ, elemMat));
1214: } else {
1215: PetscCall(PetscFEUpdateElementMat_Hybrid_Internal(feI, isCohesiveFieldI, feJ, isCohesiveFieldJ, 0, 0, 0, q, T[fieldI], basisReal, basisDerReal, T[fieldJ], testReal, testDerReal, &fegeom, g0, g1, g2, g3, eOffset, totDim, offsetI, offsetJ, elemMat));
1216: PetscCall(PetscFEUpdateElementMat_Hybrid_Internal(feI, isCohesiveFieldI, feJ, isCohesiveFieldJ, 0, 0, 1, q, T[fieldI], basisReal, basisDerReal, T[fieldJ], testReal, testDerReal, &fegeom, &g0[NcI * NcJ], &g1[NcI * NcJ * dim], &g2[NcI * NcJ * dim], &g3[NcI * NcJ * dim * dim], eOffset, totDim, offsetI, offsetJ, elemMat));
1217: PetscCall(PetscFEUpdateElementMat_Hybrid_Internal(feI, isCohesiveFieldI, feJ, isCohesiveFieldJ, 0, 1, 0, q, T[fieldI], basisReal, basisDerReal, T[fieldJ], testReal, testDerReal, &fegeom, &g0[NcI * NcJ * 2], &g1[NcI * NcJ * dim * 2], &g2[NcI * NcJ * dim * 2], &g3[NcI * NcJ * dim * dim * 2], eOffset, totDim, offsetI, offsetJ, elemMat));
1218: PetscCall(PetscFEUpdateElementMat_Hybrid_Internal(feI, isCohesiveFieldI, feJ, isCohesiveFieldJ, 0, 1, 1, q, T[fieldI], basisReal, basisDerReal, T[fieldJ], testReal, testDerReal, &fegeom, &g0[NcI * NcJ * 3], &g1[NcI * NcJ * dim * 3], &g2[NcI * NcJ * dim * 3], &g3[NcI * NcJ * dim * dim * 3], eOffset, totDim, offsetI, offsetJ, elemMat));
1219: }
1220: } else
1221: PetscCall(PetscFEUpdateElementMat_Hybrid_Internal(feI, isCohesiveFieldI, feJ, isCohesiveFieldJ, 0, s, s, q, T[fieldI], basisReal, basisDerReal, T[fieldJ], testReal, testDerReal, &fegeom, g0, g1, g2, g3, eOffset, totDim, offsetI, offsetJ, elemMat));
1222: }
1223: }
1224: if (debug > 1) {
1225: const PetscInt fS = 0 + (isCohesiveFieldI ? 0 : (s == 2 ? 0 : s * T[fieldI]->Nb));
1226: const PetscInt fE = T[fieldI]->Nb + (isCohesiveFieldI ? 0 : (s == 2 ? T[fieldI]->Nb : s * T[fieldI]->Nb));
1227: const PetscInt gS = 0 + (isCohesiveFieldJ ? 0 : (s == 2 ? 0 : s * T[fieldJ]->Nb));
1228: const PetscInt gE = T[fieldJ]->Nb + (isCohesiveFieldJ ? 0 : (s == 2 ? T[fieldJ]->Nb : s * T[fieldJ]->Nb));
1230: PetscCall(PetscPrintf(PETSC_COMM_SELF, "Element matrix for fields %" PetscInt_FMT " and %" PetscInt_FMT " s %s totDim %" PetscInt_FMT " offsets (%" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ")\n", fieldI, fieldJ, s ? (s > 1 ? "Coh" : "Pos") : "Neg", totDim, eOffset, offsetI, offsetJ));
1231: for (PetscInt f = fS; f < fE; ++f) {
1232: const PetscInt i = offsetI + f;
1233: for (PetscInt g = gS; g < gE; ++g) {
1234: const PetscInt j = offsetJ + g;
1235: PetscCheck(i < totDim && j < totDim, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Fuck up %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT, f, i, g, j);
1236: PetscCall(PetscPrintf(PETSC_COMM_SELF, " elemMat[%" PetscInt_FMT ",%" PetscInt_FMT ",%" PetscInt_FMT ",%" PetscInt_FMT "]: %g\n", f / NcI, f % NcI, g / NcJ, g % NcJ, (double)PetscRealPart(elemMat[eOffset + i * totDim + j])));
1237: }
1238: PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n"));
1239: }
1240: }
1241: cOffset += totDim;
1242: cOffsetAux += totDimAux;
1243: eOffset += PetscSqr(totDim);
1244: }
1245: PetscFunctionReturn(PETSC_SUCCESS);
1246: }
1248: static PetscErrorCode PetscFEInitialize_Basic(PetscFE fem)
1249: {
1250: PetscFunctionBegin;
1251: fem->ops->setfromoptions = NULL;
1252: fem->ops->setup = PetscFESetUp_Basic;
1253: fem->ops->view = PetscFEView_Basic;
1254: fem->ops->destroy = PetscFEDestroy_Basic;
1255: fem->ops->getdimension = PetscFEGetDimension_Basic;
1256: fem->ops->computetabulation = PetscFEComputeTabulation_Basic;
1257: fem->ops->integrate = PetscFEIntegrate_Basic;
1258: fem->ops->integratebd = PetscFEIntegrateBd_Basic;
1259: fem->ops->integrateresidual = PetscFEIntegrateResidual_Basic;
1260: fem->ops->integratebdresidual = PetscFEIntegrateBdResidual_Basic;
1261: fem->ops->integratehybridresidual = PetscFEIntegrateHybridResidual_Basic;
1262: fem->ops->integratejacobianaction = NULL /* PetscFEIntegrateJacobianAction_Basic */;
1263: fem->ops->integratejacobian = PetscFEIntegrateJacobian_Basic;
1264: fem->ops->integratebdjacobian = PetscFEIntegrateBdJacobian_Basic;
1265: fem->ops->integratehybridjacobian = PetscFEIntegrateHybridJacobian_Basic;
1266: PetscFunctionReturn(PETSC_SUCCESS);
1267: }
1269: /*MC
1270: PETSCFEBASIC = "basic" - A `PetscFE` object that integrates with basic tiling and no vectorization
1272: Level: intermediate
1274: .seealso: `PetscFE`, `PetscFEType`, `PetscFECreate()`, `PetscFESetType()`
1275: M*/
1277: PETSC_EXTERN PetscErrorCode PetscFECreate_Basic(PetscFE fem)
1278: {
1279: PetscFE_Basic *b;
1281: PetscFunctionBegin;
1283: PetscCall(PetscNew(&b));
1284: fem->data = b;
1286: PetscCall(PetscFEInitialize_Basic(fem));
1287: PetscFunctionReturn(PETSC_SUCCESS);
1288: }