Actual source code: dtaltv.c
1: #include <petsc/private/petscimpl.h>
2: #include <petsc/private/dtimpl.h>
4: /*MC
5: PetscDTAltV - An interface for common operations on k-forms, also known as alternating algebraic forms or alternating k-linear maps.
6: The name of the interface comes from the notation "Alt V" for the algebra of all k-forms acting vectors in the space V, also known as the exterior algebra of V*.
8: A recommended reference for this material is Section 2 "Exterior algebra and exterior calculus" in "Finite element
9: exterior calculus, homological techniques, and applications", by Arnold, Falk, & Winther (2006, doi:10.1017/S0962492906210018).
11: A k-form w (k is called the "form degree" of w) is an alternating k-linear map acting on tuples (v_1, ..., v_k) of
12: vectors from a vector space V and producing a real number:
13: - alternating: swapping any two vectors in a tuple reverses the sign of the result, e.g. w(v_1, v_2, ..., v_k) = -w(v_2, v_1, ..., v_k)
14: - k-linear: w acts linear in each vector separately, e.g. w(a*v + b*y, v_2, ..., v_k) = a*w(v,v_2,...,v_k) + b*w(y,v_2,...,v_k)
15: This action is implemented as PetscDTAltVApply.
17: The k-forms on a vector space form a vector space themselves, Alt^k V. The dimension of Alt^k V, if V is N dimensional, is N choose k. (This
18: shows that for an N dimensional space, only 0 <= k <= N are valid form degrees.)
19: The standard basis for Alt^k V, used in PetscDTAltV, has one basis k-form for each ordered subset of k coordinates of the N dimensional space:
20: For example, if the coordinate directions of a four dimensional space are (t, x, y, z), then there are 4 choose 2 = 6 ordered subsets of two coordinates.
21: They are, in lexicographic order, (t, x), (t, y), (t, z), (x, y), (x, z) and (y, z). PetscDTAltV also orders the basis of Alt^k V lexicographically
22: by the associated subsets.
24: The unit basis k-form associated with coordinates (c_1, ..., c_k) acts on a set of k vectors (v_1, ..., v_k) by creating a square matrix V where
25: V[i,j] = v_i[c_j] and taking the determinant of V.
27: If j + k <= N, then a j-form f and a k-form g can be multiplied to create a (j+k)-form using the wedge or exterior product, (f wedge g).
28: This is an anticommutative product, (f wedge g) = -(g wedge f). It is sufficient to describe the wedge product of two basis forms.
29: Let f be the basis j-form associated with coordinates (f_1,...,f_j) and g be the basis k-form associated with coordinates (g_1,...,g_k):
30: - If there is any coordinate in both sets, then (f wedge g) = 0.
31: - Otherwise, (f wedge g) is a multiple of the basis (j+k)-form h associated with (f_1,...,f_j,g_1,...,g_k).
32: - In fact it is equal to either h or -h depending on how (f_1,...,f_j,g_1,...,g_k) compares to the same list of coordinates given in ascending order: if it is an even permutation of that list, then (f wedge g) = h, otherwise (f wedge g) = -h.
33: The wedge product is implemented for either two inputs (f and g) in PetscDTAltVWedge, or for one (just f, giving a
34: matrix to multiply against multiple choices of g) in PetscDTAltVWedgeMatrix.
36: If k > 0, a k-form w and a vector v can combine to make a (k-1)-formm through the interior product, (w int v),
37: defined by (w int v)(v_1,...,v_{k-1}) = w(v,v_1,...,v_{k-1}).
39: The interior product is implemented for either two inputs (w and v) in PetscDTAltVInterior, for one (just v, giving a
40: matrix to multiply against multiple choices of w) in PetscDTAltVInteriorMatrix,
41: or for no inputs (giving the sparsity pattern of PetscDTAltVInteriorMatrix) in PetscDTAltVInteriorPattern.
43: When there is a linear map L: V -> W from an N dimensional vector space to an M dimensional vector space,
44: it induces the linear pullback map L^* : Alt^k W -> Alt^k V, defined by L^* w(v_1,...,v_k) = w(L v_1, ..., L v_k).
45: The pullback is implemented as PetscDTAltVPullback (acting on a known w) and PetscDTAltVPullbackMatrix (creating a matrix that computes the actin of L^*).
47: Alt^k V and Alt^(N-k) V have the same dimension, and the Hodge star operator maps between them. We note that Alt^N V is a one dimensional space, and its
48: basis vector is sometime called vol. The Hodge star operator has the property that (f wedge (star g)) = (f,g) vol, where (f,g) is the simple inner product
49: of the basis coefficients of f and g.
50: Powers of the Hodge star operator can be applied with PetscDTAltVStar
52: Level: intermediate
54: .seealso: PetscDTAltVApply(), PetscDTAltVWedge(), PetscDTAltVInterior(), PetscDTAltVPullback(), PetscDTAltVStar()
55: M*/
57: /*@
58: PetscDTAltVApply - Apply an a k-form (an alternating k-linear map) to a set of k N-dimensional vectors
60: Input Parameters:
61: + N - the dimension of the vector space, N >= 0
62: . k - the degree k of the k-form w, 0 <= k <= N
63: . w - a k-form, size [N choose k] (each degree of freedom of a k-form is associated with a subset of k coordinates of the N-dimensional vectors: the degrees of freedom are ordered lexicographically by their associated subsets)
64: - v - a set of k vectors of size N, size [k x N], each vector stored contiguously
66: Output Parameter:
67: . wv - w(v_1,...,v_k) = \sum_i w_i * det(V_i): the degree of freedom w_i is associated with coordinates [s_{i,1},...,s_{i,k}], and the square matrix V_i has entry (j,k) given by the s_{i,k}'th coordinate of v_j
69: Level: intermediate
71: .seealso: PetscDTAltV, PetscDTAltVPullback(), PetscDTAltVPullbackMatrix()
72: @*/
73: PetscErrorCode PetscDTAltVApply(PetscInt N, PetscInt k, const PetscReal *w, const PetscReal *v, PetscReal *wv)
74: {
77: if (N <= 3) {
78: if (!k) {
79: *wv = w[0];
80: } else {
81: if (N == 1) {*wv = w[0] * v[0];}
82: else if (N == 2) {
83: if (k == 1) {*wv = w[0] * v[0] + w[1] * v[1];}
84: else {*wv = w[0] * (v[0] * v[3] - v[1] * v[2]);}
85: } else {
86: if (k == 1) {*wv = w[0] * v[0] + w[1] * v[1] + w[2] * v[2];}
87: else if (k == 2) {
88: *wv = w[0] * (v[0] * v[4] - v[1] * v[3]) +
89: w[1] * (v[0] * v[5] - v[2] * v[3]) +
90: w[2] * (v[1] * v[5] - v[2] * v[4]);
91: } else {
92: *wv = w[0] * (v[0] * (v[4] * v[8] - v[5] * v[7]) +
93: v[1] * (v[5] * v[6] - v[3] * v[8]) +
94: v[2] * (v[3] * v[7] - v[4] * v[6]));
95: }
96: }
97: }
98: } else {
99: PetscInt Nk, Nf;
100: PetscInt *subset, *perm;
101: PetscInt i, j, l;
102: PetscReal sum = 0.;
104: PetscDTFactorialInt(k, &Nf);
105: PetscDTBinomialInt(N, k, &Nk);
106: PetscMalloc2(k, &subset, k, &perm);
107: for (i = 0; i < Nk; i++) {
108: PetscReal subsum = 0.;
110: PetscDTEnumSubset(N, k, i, subset);
111: for (j = 0; j < Nf; j++) {
112: PetscBool permOdd;
113: PetscReal prod;
115: PetscDTEnumPerm(k, j, perm, &permOdd);
116: prod = permOdd ? -1. : 1.;
117: for (l = 0; l < k; l++) {
118: prod *= v[perm[l] * N + subset[l]];
119: }
120: subsum += prod;
121: }
122: sum += w[i] * subsum;
123: }
124: PetscFree2(subset, perm);
125: *wv = sum;
126: }
127: return 0;
128: }
130: /*@
131: PetscDTAltVWedge - Compute the wedge product of a j-form and a k-form, giving a (j+k) form
133: Input Parameters:
134: + N - the dimension of the vector space, N >= 0
135: . j - the degree j of the j-form a, 0 <= j <= N
136: . k - the degree k of the k-form b, 0 <= k <= N and 0 <= j+k <= N
137: . a - a j-form, size [N choose j]
138: - b - a k-form, size [N choose k]
140: Output Parameter:
141: . awedgeb - the (j+k)-form a wedge b, size [N choose (j+k)]: (a wedge b)(v_1,...,v_{j+k}) = \sum_{s} sign(s) a(v_{s_1},...,v_{s_j}) b(v_{s_{j+1}},...,v_{s_{j+k}}),
142: where the sum is over permutations s such that s_1 < s_2 < ... < s_j and s_{j+1} < s_{j+2} < ... < s_{j+k}.
144: Level: intermediate
146: .seealso: PetscDTAltV, PetscDTAltVWedgeMatrix(), PetscDTAltVPullback(), PetscDTAltVPullbackMatrix()
147: @*/
148: PetscErrorCode PetscDTAltVWedge(PetscInt N, PetscInt j, PetscInt k, const PetscReal *a, const PetscReal *b, PetscReal *awedgeb)
149: {
150: PetscInt i;
155: if (N <= 3) {
156: PetscInt Njk;
158: PetscDTBinomialInt(N, j+k, &Njk);
159: if (!j) {for (i = 0; i < Njk; i++) {awedgeb[i] = a[0] * b[i];}}
160: else if (!k) {for (i = 0; i < Njk; i++) {awedgeb[i] = a[i] * b[0];}}
161: else {
162: if (N == 2) {awedgeb[0] = a[0] * b[1] - a[1] * b[0];}
163: else {
164: if (j+k == 2) {
165: awedgeb[0] = a[0] * b[1] - a[1] * b[0];
166: awedgeb[1] = a[0] * b[2] - a[2] * b[0];
167: awedgeb[2] = a[1] * b[2] - a[2] * b[1];
168: } else {
169: awedgeb[0] = a[0] * b[2] - a[1] * b[1] + a[2] * b[0];
170: }
171: }
172: }
173: } else {
174: PetscInt Njk;
175: PetscInt JKj;
176: PetscInt *subset, *subsetjk, *subsetj, *subsetk;
177: PetscInt i;
179: PetscDTBinomialInt(N, j+k, &Njk);
180: PetscDTBinomialInt(j+k, j, &JKj);
181: PetscMalloc4(j+k, &subset, j+k, &subsetjk, j, &subsetj, k, &subsetk);
182: for (i = 0; i < Njk; i++) {
183: PetscReal sum = 0.;
184: PetscInt l;
186: PetscDTEnumSubset(N, j+k, i, subset);
187: for (l = 0; l < JKj; l++) {
188: PetscBool jkOdd;
189: PetscInt m, jInd, kInd;
191: PetscDTEnumSplit(j+k, j, l, subsetjk, &jkOdd);
192: for (m = 0; m < j; m++) {
193: subsetj[m] = subset[subsetjk[m]];
194: }
195: for (m = 0; m < k; m++) {
196: subsetk[m] = subset[subsetjk[j+m]];
197: }
198: PetscDTSubsetIndex(N, j, subsetj, &jInd);
199: PetscDTSubsetIndex(N, k, subsetk, &kInd);
200: sum += jkOdd ? -(a[jInd] * b[kInd]) : (a[jInd] * b[kInd]);
201: }
202: awedgeb[i] = sum;
203: }
204: PetscFree4(subset, subsetjk, subsetj, subsetk);
205: }
206: return 0;
207: }
209: /*@
210: PetscDTAltVWedgeMatrix - Compute the matrix defined by the wedge product with a given j-form that maps k-forms to (j+k)-forms
212: Input Parameters:
213: + N - the dimension of the vector space, N >= 0
214: . j - the degree j of the j-form a, 0 <= j <= N
215: . k - the degree k of the k-forms that (a wedge) will be applied to, 0 <= k <= N and 0 <= j+k <= N
216: - a - a j-form, size [N choose j]
218: Output Parameter:
219: . awedge - (a wedge), an [(N choose j+k) x (N choose k)] matrix in row-major order, such that (a wedge) * b = a wedge b
221: Level: intermediate
223: .seealso: PetscDTAltV, PetscDTAltVPullback(), PetscDTAltVPullbackMatrix()
224: @*/
225: PetscErrorCode PetscDTAltVWedgeMatrix(PetscInt N, PetscInt j, PetscInt k, const PetscReal *a, PetscReal *awedgeMat)
226: {
227: PetscInt i;
232: if (N <= 3) {
233: PetscInt Njk;
235: PetscDTBinomialInt(N, j+k, &Njk);
236: if (!j) {
237: for (i = 0; i < Njk * Njk; i++) {awedgeMat[i] = 0.;}
238: for (i = 0; i < Njk; i++) {awedgeMat[i * (Njk + 1)] = a[0];}
239: } else if (!k) {
240: for (i = 0; i < Njk; i++) {awedgeMat[i] = a[i];}
241: } else {
242: if (N == 2) {
243: awedgeMat[0] = -a[1]; awedgeMat[1] = a[0];
244: } else {
245: if (j+k == 2) {
246: awedgeMat[0] = -a[1]; awedgeMat[1] = a[0]; awedgeMat[2] = 0.;
247: awedgeMat[3] = -a[2]; awedgeMat[4] = 0.; awedgeMat[5] = a[0];
248: awedgeMat[6] = 0.; awedgeMat[7] = -a[2]; awedgeMat[8] = a[1];
249: } else {
250: awedgeMat[0] = a[2]; awedgeMat[1] = -a[1]; awedgeMat[2] = a[0];
251: }
252: }
253: }
254: } else {
255: PetscInt Njk;
256: PetscInt Nk;
257: PetscInt JKj, i;
258: PetscInt *subset, *subsetjk, *subsetj, *subsetk;
260: PetscDTBinomialInt(N, k, &Nk);
261: PetscDTBinomialInt(N, j+k, &Njk);
262: PetscDTBinomialInt(j+k, j, &JKj);
263: PetscMalloc4(j+k, &subset, j+k, &subsetjk, j, &subsetj, k, &subsetk);
264: for (i = 0; i < Njk * Nk; i++) awedgeMat[i] = 0.;
265: for (i = 0; i < Njk; i++) {
266: PetscInt l;
268: PetscDTEnumSubset(N, j+k, i, subset);
269: for (l = 0; l < JKj; l++) {
270: PetscBool jkOdd;
271: PetscInt m, jInd, kInd;
273: PetscDTEnumSplit(j+k, j, l, subsetjk, &jkOdd);
274: for (m = 0; m < j; m++) {
275: subsetj[m] = subset[subsetjk[m]];
276: }
277: for (m = 0; m < k; m++) {
278: subsetk[m] = subset[subsetjk[j+m]];
279: }
280: PetscDTSubsetIndex(N, j, subsetj, &jInd);
281: PetscDTSubsetIndex(N, k, subsetk, &kInd);
282: awedgeMat[i * Nk + kInd] += jkOdd ? - a[jInd] : a[jInd];
283: }
284: }
285: PetscFree4(subset, subsetjk, subsetj, subsetk);
286: }
287: return 0;
288: }
290: /*@
291: PetscDTAltVPullback - Compute the pullback of a k-form under a linear transformation of the coordinate space
293: Input Parameters:
294: + N - the dimension of the origin vector space of the linear transformation, M >= 0
295: . M - the dimension of the image vector space of the linear transformation, N >= 0
296: . L - a linear transformation, an [M x N] matrix in row-major format
297: . k - the *signed* degree k of the |k|-form w, -(min(M,N)) <= k <= min(M,N). A negative form degree indicates that the pullback should be conjugated by the Hodge star operator (see note).
298: - w - a |k|-form in the image space, size [M choose |k|]
300: Output Parameter:
301: . Lstarw - the pullback of w to a |k|-form in the origin space, size [N choose |k|]: (Lstarw)(v_1,...v_k) = w(L*v_1,...,L*v_k).
303: Level: intermediate
305: Note: negative form degrees accommodate, e.g., H-div conforming vector fields. An H-div conforming vector field stores its degrees of freedom as (dx, dy, dz), like a 1-form,
306: but its normal trace is integrated on faces, like a 2-form. The correct pullback then is to apply the Hodge star transformation from (M-2)-form to 2-form, pullback as a 2-form,
307: then the inverse Hodge star transformation.
309: .seealso: PetscDTAltV, PetscDTAltVPullbackMatrix(), PetscDTAltVStar()
310: @*/
311: PetscErrorCode PetscDTAltVPullback(PetscInt N, PetscInt M, const PetscReal *L, PetscInt k, const PetscReal *w, PetscReal *Lstarw)
312: {
313: PetscInt i, j, Nk, Mk;
317: if (N <= 3 && M <= 3) {
319: PetscDTBinomialInt(M, PetscAbsInt(k), &Mk);
320: PetscDTBinomialInt(N, PetscAbsInt(k), &Nk);
321: if (!k) {
322: Lstarw[0] = w[0];
323: } else if (k == 1) {
324: for (i = 0; i < Nk; i++) {
325: PetscReal sum = 0.;
327: for (j = 0; j < Mk; j++) {sum += L[j * Nk + i] * w[j];}
328: Lstarw[i] = sum;
329: }
330: } else if (k == -1) {
331: PetscReal mult[3] = {1., -1., 1.};
333: for (i = 0; i < Nk; i++) {
334: PetscReal sum = 0.;
336: for (j = 0; j < Mk; j++) {
337: sum += L[(Mk - 1 - j) * Nk + (Nk - 1 - i)] * w[j] * mult[j];
338: }
339: Lstarw[i] = mult[i] * sum;
340: }
341: } else if (k == 2) {
342: PetscInt pairs[3][2] = {{0,1},{0,2},{1,2}};
344: for (i = 0; i < Nk; i++) {
345: PetscReal sum = 0.;
346: for (j = 0; j < Mk; j++) {
347: sum += (L[pairs[j][0] * N + pairs[i][0]] * L[pairs[j][1] * N + pairs[i][1]] -
348: L[pairs[j][1] * N + pairs[i][0]] * L[pairs[j][0] * N + pairs[i][1]]) * w[j];
349: }
350: Lstarw[i] = sum;
351: }
352: } else if (k == -2) {
353: PetscInt pairs[3][2] = {{1,2},{2,0},{0,1}};
354: PetscInt offi = (N == 2) ? 2 : 0;
355: PetscInt offj = (M == 2) ? 2 : 0;
357: for (i = 0; i < Nk; i++) {
358: PetscReal sum = 0.;
360: for (j = 0; j < Mk; j++) {
361: sum += (L[pairs[offj + j][0] * N + pairs[offi + i][0]] *
362: L[pairs[offj + j][1] * N + pairs[offi + i][1]] -
363: L[pairs[offj + j][1] * N + pairs[offi + i][0]] *
364: L[pairs[offj + j][0] * N + pairs[offi + i][1]]) * w[j];
366: }
367: Lstarw[i] = sum;
368: }
369: } else {
370: PetscReal detL = L[0] * (L[4] * L[8] - L[5] * L[7]) +
371: L[1] * (L[5] * L[6] - L[3] * L[8]) +
372: L[2] * (L[3] * L[7] - L[4] * L[6]);
374: for (i = 0; i < Nk; i++) {Lstarw[i] = detL * w[i];}
375: }
376: } else {
377: PetscInt Nf, l, p;
378: PetscReal *Lw, *Lwv;
379: PetscInt *subsetw, *subsetv;
380: PetscInt *perm;
381: PetscReal *walloc = NULL;
382: const PetscReal *ww = NULL;
383: PetscBool negative = PETSC_FALSE;
385: PetscDTBinomialInt(M, PetscAbsInt(k), &Mk);
386: PetscDTBinomialInt(N, PetscAbsInt(k), &Nk);
387: PetscDTFactorialInt(PetscAbsInt(k), &Nf);
388: if (k < 0) {
389: negative = PETSC_TRUE;
390: k = -k;
391: PetscMalloc1(Mk, &walloc);
392: PetscDTAltVStar(M, M - k, 1, w, walloc);
393: ww = walloc;
394: } else {
395: ww = w;
396: }
397: PetscMalloc5(k, &subsetw, k, &subsetv, k, &perm, N * k, &Lw, k * k, &Lwv);
398: for (i = 0; i < Nk; i++) Lstarw[i] = 0.;
399: for (i = 0; i < Mk; i++) {
400: PetscDTEnumSubset(M, k, i, subsetw);
401: for (j = 0; j < Nk; j++) {
402: PetscDTEnumSubset(N, k, j, subsetv);
403: for (p = 0; p < Nf; p++) {
404: PetscReal prod;
405: PetscBool isOdd;
407: PetscDTEnumPerm(k, p, perm, &isOdd);
408: prod = isOdd ? -ww[i] : ww[i];
409: for (l = 0; l < k; l++) {
410: prod *= L[subsetw[perm[l]] * N + subsetv[l]];
411: }
412: Lstarw[j] += prod;
413: }
414: }
415: }
416: if (negative) {
417: PetscReal *sLsw;
419: PetscMalloc1(Nk, &sLsw);
420: PetscDTAltVStar(N, N - k, -1, Lstarw, sLsw);
421: for (i = 0; i < Nk; i++) Lstarw[i] = sLsw[i];
422: PetscFree(sLsw);
423: }
424: PetscFree5(subsetw, subsetv, perm, Lw, Lwv);
425: PetscFree(walloc);
426: }
427: return 0;
428: }
430: /*@
431: PetscDTAltVPullbackMatrix - Compute the pullback matrix for k-forms under a linear transformation
433: Input Parameters:
434: + N - the dimension of the origin vector space of the linear transformation, N >= 0
435: . M - the dimension of the image vector space of the linear transformation, M >= 0
436: . L - a linear transformation, an [M x N] matrix in row-major format
437: - k - the *signed* degree k of the |k|-forms on which Lstar acts, -(min(M,N)) <= k <= min(M,N). A negative form degree indicates that the pullback should be conjugated by the Hodge star operator (see note in PetscDTAltvPullback())
439: Output Parameter:
440: . Lstar - the pullback matrix, an [(N choose |k|) x (M choose |k|)] matrix in row-major format such that Lstar * w = L^* w
442: Level: intermediate
444: .seealso: PetscDTAltV, PetscDTAltVPullback(), PetscDTAltVStar()
445: @*/
446: PetscErrorCode PetscDTAltVPullbackMatrix(PetscInt N, PetscInt M, const PetscReal *L, PetscInt k, PetscReal *Lstar)
447: {
448: PetscInt Nk, Mk, Nf, i, j, l, p;
449: PetscReal *Lw, *Lwv;
450: PetscInt *subsetw, *subsetv;
451: PetscInt *perm;
452: PetscBool negative = PETSC_FALSE;
456: if (N <= 3 && M <= 3) {
457: PetscReal mult[3] = {1., -1., 1.};
459: PetscDTBinomialInt(M, PetscAbsInt(k), &Mk);
460: PetscDTBinomialInt(N, PetscAbsInt(k), &Nk);
461: if (!k) {
462: Lstar[0] = 1.;
463: } else if (k == 1) {
464: for (i = 0; i < Nk; i++) {for (j = 0; j < Mk; j++) {Lstar[i * Mk + j] = L[j * Nk + i];}}
465: } else if (k == -1) {
466: for (i = 0; i < Nk; i++) {
467: for (j = 0; j < Mk; j++) {
468: Lstar[i * Mk + j] = L[(Mk - 1 - j) * Nk + (Nk - 1 - i)] * mult[i] * mult[j];
469: }
470: }
471: } else if (k == 2) {
472: PetscInt pairs[3][2] = {{0,1},{0,2},{1,2}};
474: for (i = 0; i < Nk; i++) {
475: for (j = 0; j < Mk; j++) {
476: Lstar[i * Mk + j] = L[pairs[j][0] * N + pairs[i][0]] *
477: L[pairs[j][1] * N + pairs[i][1]] -
478: L[pairs[j][1] * N + pairs[i][0]] *
479: L[pairs[j][0] * N + pairs[i][1]];
480: }
481: }
482: } else if (k == -2) {
483: PetscInt pairs[3][2] = {{1,2},{2,0},{0,1}};
484: PetscInt offi = (N == 2) ? 2 : 0;
485: PetscInt offj = (M == 2) ? 2 : 0;
487: for (i = 0; i < Nk; i++) {
488: for (j = 0; j < Mk; j++) {
489: Lstar[i * Mk + j] = L[pairs[offj + j][0] * N + pairs[offi + i][0]] *
490: L[pairs[offj + j][1] * N + pairs[offi + i][1]] -
491: L[pairs[offj + j][1] * N + pairs[offi + i][0]] *
492: L[pairs[offj + j][0] * N + pairs[offi + i][1]];
493: }
494: }
495: } else {
496: PetscReal detL = L[0] * (L[4] * L[8] - L[5] * L[7]) +
497: L[1] * (L[5] * L[6] - L[3] * L[8]) +
498: L[2] * (L[3] * L[7] - L[4] * L[6]);
500: for (i = 0; i < Nk; i++) {Lstar[i] = detL;}
501: }
502: } else {
503: if (k < 0) {
504: negative = PETSC_TRUE;
505: k = -k;
506: }
507: PetscDTBinomialInt(M, PetscAbsInt(k), &Mk);
508: PetscDTBinomialInt(N, PetscAbsInt(k), &Nk);
509: PetscDTFactorialInt(PetscAbsInt(k), &Nf);
510: PetscMalloc5(M, &subsetw, N, &subsetv, k, &perm, N * k, &Lw, k * k, &Lwv);
511: for (i = 0; i < Nk * Mk; i++) Lstar[i] = 0.;
512: for (i = 0; i < Mk; i++) {
513: PetscBool iOdd;
514: PetscInt iidx, jidx;
516: PetscDTEnumSplit(M, k, i, subsetw, &iOdd);
517: iidx = negative ? Mk - 1 - i : i;
518: iOdd = negative ? (PetscBool) (iOdd ^ ((k * (M-k)) & 1)) : PETSC_FALSE;
519: for (j = 0; j < Nk; j++) {
520: PetscBool jOdd;
522: PetscDTEnumSplit(N, k, j, subsetv, &jOdd);
523: jidx = negative ? Nk - 1 - j : j;
524: jOdd = negative ? (PetscBool) (iOdd ^ jOdd ^ ((k * (N-k)) & 1)) : PETSC_FALSE;
525: for (p = 0; p < Nf; p++) {
526: PetscReal prod;
527: PetscBool isOdd;
529: PetscDTEnumPerm(k, p, perm, &isOdd);
530: isOdd = (PetscBool) (isOdd ^ jOdd);
531: prod = isOdd ? -1. : 1.;
532: for (l = 0; l < k; l++) {
533: prod *= L[subsetw[perm[l]] * N + subsetv[l]];
534: }
535: Lstar[jidx * Mk + iidx] += prod;
536: }
537: }
538: }
539: PetscFree5(subsetw, subsetv, perm, Lw, Lwv);
540: }
541: return 0;
542: }
544: /*@
545: PetscDTAltVInterior - Compute the interior product of a k-form with a vector
547: Input Parameters:
548: + N - the dimension of the vector space, N >= 0
549: . k - the degree k of the k-form w, 0 <= k <= N
550: . w - a k-form, size [N choose k]
551: - v - an N dimensional vector
553: Output Parameter:
554: . wIntv - the (k-1)-form (w int v), size [N choose (k-1)]: (w int v) is defined by its action on (k-1) vectors {v_1, ..., v_{k-1}} as (w inv v)(v_1, ..., v_{k-1}) = w(v, v_1, ..., v_{k-1}).
556: Level: intermediate
558: .seealso: PetscDTAltV, PetscDTAltVInteriorMatrix(), PetscDTAltVInteriorPattern(), PetscDTAltVPullback(), PetscDTAltVPullbackMatrix()
559: @*/
560: PetscErrorCode PetscDTAltVInterior(PetscInt N, PetscInt k, const PetscReal *w, const PetscReal *v, PetscReal *wIntv)
561: {
562: PetscInt i, Nk, Nkm;
565: PetscDTBinomialInt(N, k, &Nk);
566: PetscDTBinomialInt(N, k-1, &Nkm);
567: if (N <= 3) {
568: if (k == 1) {
569: PetscReal sum = 0.;
571: for (i = 0; i < N; i++) {
572: sum += w[i] * v[i];
573: }
574: wIntv[0] = sum;
575: } else if (k == N) {
576: PetscReal mult[3] = {1., -1., 1.};
578: for (i = 0; i < N; i++) {
579: wIntv[N - 1 - i] = w[0] * v[i] * mult[i];
580: }
581: } else {
582: wIntv[0] = - w[0]*v[1] - w[1]*v[2];
583: wIntv[1] = w[0]*v[0] - w[2]*v[2];
584: wIntv[2] = w[1]*v[0] + w[2]*v[1];
585: }
586: } else {
587: PetscInt *subset, *work;
589: PetscMalloc2(k, &subset, k, &work);
590: for (i = 0; i < Nkm; i++) wIntv[i] = 0.;
591: for (i = 0; i < Nk; i++) {
592: PetscInt j, l, m;
594: PetscDTEnumSubset(N, k, i, subset);
595: for (j = 0; j < k; j++) {
596: PetscInt idx;
597: PetscBool flip = (PetscBool) (j & 1);
599: for (l = 0, m = 0; l < k; l++) {
600: if (l != j) work[m++] = subset[l];
601: }
602: PetscDTSubsetIndex(N, k - 1, work, &idx);
603: wIntv[idx] += flip ? -(w[i] * v[subset[j]]) : (w[i] * v[subset[j]]);
604: }
605: }
606: PetscFree2(subset, work);
607: }
608: return 0;
609: }
611: /*@
612: PetscDTAltVInteriorMatrix - Compute the matrix of the linear transformation induced on a k-form by the interior product with a vector
614: Input Parameters:
615: + N - the dimension of the vector space, N >= 0
616: . k - the degree k of the k-forms on which intvMat acts, 0 <= k <= N
617: - v - an N dimensional vector
619: Output Parameter:
620: . intvMat - an [(N choose (k-1)) x (N choose k)] matrix, row-major: (intvMat) * w = (w int v)
622: Level: intermediate
624: .seealso: PetscDTAltV, PetscDTAltVInterior(), PetscDTAltVInteriorPattern(), PetscDTAltVPullback(), PetscDTAltVPullbackMatrix()
625: @*/
626: PetscErrorCode PetscDTAltVInteriorMatrix(PetscInt N, PetscInt k, const PetscReal *v, PetscReal *intvMat)
627: {
628: PetscInt i, Nk, Nkm;
631: PetscDTBinomialInt(N, k, &Nk);
632: PetscDTBinomialInt(N, k-1, &Nkm);
633: if (N <= 3) {
634: if (k == 1) {
635: for (i = 0; i < N; i++) intvMat[i] = v[i];
636: } else if (k == N) {
637: PetscReal mult[3] = {1., -1., 1.};
639: for (i = 0; i < N; i++) intvMat[N - 1 - i] = v[i] * mult[i];
640: } else {
641: intvMat[0] = -v[1]; intvMat[1] = -v[2]; intvMat[2] = 0.;
642: intvMat[3] = v[0]; intvMat[4] = 0.; intvMat[5] = -v[2];
643: intvMat[6] = 0.; intvMat[7] = v[0]; intvMat[8] = v[1];
644: }
645: } else {
646: PetscInt *subset, *work;
648: PetscMalloc2(k, &subset, k, &work);
649: for (i = 0; i < Nk * Nkm; i++) intvMat[i] = 0.;
650: for (i = 0; i < Nk; i++) {
651: PetscInt j, l, m;
653: PetscDTEnumSubset(N, k, i, subset);
654: for (j = 0; j < k; j++) {
655: PetscInt idx;
656: PetscBool flip = (PetscBool) (j & 1);
658: for (l = 0, m = 0; l < k; l++) {
659: if (l != j) work[m++] = subset[l];
660: }
661: PetscDTSubsetIndex(N, k - 1, work, &idx);
662: intvMat[idx * Nk + i] += flip ? -v[subset[j]] : v[subset[j]];
663: }
664: }
665: PetscFree2(subset, work);
666: }
667: return 0;
668: }
670: /*@
671: PetscDTAltVInteriorPattern - compute the sparsity and sign pattern of the interior product matrix computed in PetscDTAltVInteriorMatrix()
673: Input Parameters:
674: + N - the dimension of the vector space, N >= 0
675: - k - the degree of the k-forms on which intvMat from PetscDTAltVInteriorMatrix() acts, 0 <= k <= N.
677: Output Parameter:
678: . indices - The interior product matrix intvMat has size [(N choose (k-1)) x (N choose k)] and has (N choose k) * k
679: non-zeros. indices[i][0] and indices[i][1] are the row and column of a non-zero, and its value is equal to the vector
680: coordinate v[j] if indices[i][2] = j, or -v[j] if indices[i][2] = -(j+1)
682: Level: intermediate
684: Note: this function is useful when the interior product needs to be computed at multiple locations, as when computing the Koszul differential
686: .seealso: PetscDTAltV, PetscDTAltVInterior(), PetscDTAltVInteriorMatrix(), PetscDTAltVPullback(), PetscDTAltVPullbackMatrix()
687: @*/
688: PetscErrorCode PetscDTAltVInteriorPattern(PetscInt N, PetscInt k, PetscInt (*indices)[3])
689: {
690: PetscInt i, Nk, Nkm;
693: PetscDTBinomialInt(N, k, &Nk);
694: PetscDTBinomialInt(N, k-1, &Nkm);
695: if (N <= 3) {
696: if (k == 1) {
697: for (i = 0; i < N; i++) {
698: indices[i][0] = 0;
699: indices[i][1] = i;
700: indices[i][2] = i;
701: }
702: } else if (k == N) {
703: PetscInt val[3] = {0, -2, 2};
705: for (i = 0; i < N; i++) {
706: indices[i][0] = N - 1 - i;
707: indices[i][1] = 0;
708: indices[i][2] = val[i];
709: }
710: } else {
711: indices[0][0] = 0; indices[0][1] = 0; indices[0][2] = -(1 + 1);
712: indices[1][0] = 0; indices[1][1] = 1; indices[1][2] = -(2 + 1);
713: indices[2][0] = 1; indices[2][1] = 0; indices[2][2] = 0;
714: indices[3][0] = 1; indices[3][1] = 2; indices[3][2] = -(2 + 1);
715: indices[4][0] = 2; indices[4][1] = 1; indices[4][2] = 0;
716: indices[5][0] = 2; indices[5][1] = 2; indices[5][2] = 1;
717: }
718: } else {
719: PetscInt *subset, *work;
721: PetscMalloc2(k, &subset, k, &work);
722: for (i = 0; i < Nk; i++) {
723: PetscInt j, l, m;
725: PetscDTEnumSubset(N, k, i, subset);
726: for (j = 0; j < k; j++) {
727: PetscInt idx;
728: PetscBool flip = (PetscBool) (j & 1);
730: for (l = 0, m = 0; l < k; l++) {
731: if (l != j) work[m++] = subset[l];
732: }
733: PetscDTSubsetIndex(N, k - 1, work, &idx);
734: indices[i * k + j][0] = idx;
735: indices[i * k + j][1] = i;
736: indices[i * k + j][2] = flip ? -(subset[j] + 1) : subset[j];
737: }
738: }
739: PetscFree2(subset, work);
740: }
741: return 0;
742: }
744: /*@
745: PetscDTAltVStar - Apply a power of the Hodge star operator, which maps k-forms to (N-k) forms, to a k-form
747: Input Parameters:
748: + N - the dimension of the vector space, N >= 0
749: . k - the degree k of the k-form w, 0 <= k <= N
750: . pow - the number of times to apply the Hodge star operator: pow < 0 indicates that the inverse of the Hodge star operator should be applied |pow| times.
751: - w - a k-form, size [N choose k]
753: Output Parameter:
754: . starw = (star)^pow w. Each degree of freedom of a k-form is associated with a subset S of k coordinates of the N dimensional vector space: the Hodge start operator (star) maps that degree of freedom to the degree of freedom associated with S', the complement of S, with a sign change if the permutation of coordinates {S[0], ... S[k-1], S'[0], ... S'[N-k- 1]} is an odd permutation. This implies (star)^2 w = (-1)^{k(N-k)} w, and (star)^4 w = w.
756: Level: intermediate
758: .seealso: PetscDTAltV, PetscDTAltVPullback(), PetscDTAltVPullbackMatrix()
759: @*/
760: PetscErrorCode PetscDTAltVStar(PetscInt N, PetscInt k, PetscInt pow, const PetscReal *w, PetscReal *starw)
761: {
762: PetscInt Nk, i;
765: PetscDTBinomialInt(N, k, &Nk);
766: pow = pow % 4;
767: pow = (pow + 4) % 4; /* make non-negative */
768: /* pow is now 0, 1, 2, 3 */
769: if (N <= 3) {
770: if (pow & 1) {
771: PetscReal mult[3] = {1., -1., 1.};
773: for (i = 0; i < Nk; i++) starw[Nk - 1 - i] = w[i] * mult[i];
774: } else {
775: for (i = 0; i < Nk; i++) starw[i] = w[i];
776: }
777: if (pow > 1 && ((k * (N - k)) & 1)) {
778: for (i = 0; i < Nk; i++) starw[i] = -starw[i];
779: }
780: } else {
781: PetscInt *subset;
783: PetscMalloc1(N, &subset);
784: if (pow % 2) {
785: PetscInt l = (pow == 1) ? k : N - k;
786: for (i = 0; i < Nk; i++) {
787: PetscBool sOdd;
788: PetscInt j, idx;
790: PetscDTEnumSplit(N, l, i, subset, &sOdd);
791: PetscDTSubsetIndex(N, l, subset, &idx);
792: PetscDTSubsetIndex(N, N-l, &subset[l], &j);
793: starw[j] = sOdd ? -w[idx] : w[idx];
794: }
795: } else {
796: for (i = 0; i < Nk; i++) starw[i] = w[i];
797: }
798: /* star^2 = -1^(k * (N - k)) */
799: if (pow > 1 && (k * (N - k)) % 2) {
800: for (i = 0; i < Nk; i++) starw[i] = -starw[i];
801: }
802: PetscFree(subset);
803: }
804: return 0;
805: }